xref: /petsc/src/mat/impls/baij/seq/baijfact.c (revision e2ee2a47cae8248a2ea96b9d6f35ff3f37ecfc00)
1 #define PETSCMAT_DLL
2 
3 /*
4     Factorization code for BAIJ format.
5 */
6 #include "src/mat/impls/baij/seq/baij.h"
7 #include "src/inline/ilu.h"
8 
9 /* ------------------------------------------------------------*/
10 /*
11       Version for when blocks are 2 by 2
12 */
13 #undef __FUNCT__
14 #define __FUNCT__ "MatLUFactorNumeric_SeqBAIJ_2"
15 PetscErrorCode MatLUFactorNumeric_SeqBAIJ_2(Mat A,MatFactorInfo *info,Mat *B)
16 {
17   Mat            C = *B;
18   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ *)C->data;
19   IS             isrow = b->row,isicol = b->icol;
20   PetscErrorCode ierr;
21   PetscInt       *r,*ic,i,j,n = a->mbs,*bi = b->i,*bj = b->j;
22   PetscInt       *ajtmpold,*ajtmp,nz,row;
23   PetscInt       *diag_offset=b->diag,idx,*ai=a->i,*aj=a->j,*pj;
24   MatScalar      *pv,*v,*rtmp,m1,m2,m3,m4,*pc,*w,*x,x1,x2,x3,x4;
25   MatScalar      p1,p2,p3,p4;
26   MatScalar      *ba = b->a,*aa = a->a;
27 
28   PetscFunctionBegin;
29   ierr  = ISGetIndices(isrow,&r);CHKERRQ(ierr);
30   ierr  = ISGetIndices(isicol,&ic);CHKERRQ(ierr);
31   ierr  = PetscMalloc(4*(n+1)*sizeof(MatScalar),&rtmp);CHKERRQ(ierr);
32 
33   for (i=0; i<n; i++) {
34     nz    = bi[i+1] - bi[i];
35     ajtmp = bj + bi[i];
36     for  (j=0; j<nz; j++) {
37       x = rtmp+4*ajtmp[j]; x[0] = x[1] = x[2] = x[3] = 0.0;
38     }
39     /* load in initial (unfactored row) */
40     idx      = r[i];
41     nz       = ai[idx+1] - ai[idx];
42     ajtmpold = aj + ai[idx];
43     v        = aa + 4*ai[idx];
44     for (j=0; j<nz; j++) {
45       x    = rtmp+4*ic[ajtmpold[j]];
46       x[0] = v[0]; x[1] = v[1]; x[2] = v[2]; x[3] = v[3];
47       v    += 4;
48     }
49     row = *ajtmp++;
50     while (row < i) {
51       pc = rtmp + 4*row;
52       p1 = pc[0]; p2 = pc[1]; p3 = pc[2]; p4 = pc[3];
53       if (p1 != 0.0 || p2 != 0.0 || p3 != 0.0 || p4 != 0.0) {
54         pv = ba + 4*diag_offset[row];
55         pj = bj + diag_offset[row] + 1;
56         x1 = pv[0]; x2 = pv[1]; x3 = pv[2]; x4 = pv[3];
57         pc[0] = m1 = p1*x1 + p3*x2;
58         pc[1] = m2 = p2*x1 + p4*x2;
59         pc[2] = m3 = p1*x3 + p3*x4;
60         pc[3] = m4 = p2*x3 + p4*x4;
61         nz = bi[row+1] - diag_offset[row] - 1;
62         pv += 4;
63         for (j=0; j<nz; j++) {
64           x1   = pv[0]; x2 = pv[1]; x3 = pv[2]; x4 = pv[3];
65           x    = rtmp + 4*pj[j];
66           x[0] -= m1*x1 + m3*x2;
67           x[1] -= m2*x1 + m4*x2;
68           x[2] -= m1*x3 + m3*x4;
69           x[3] -= m2*x3 + m4*x4;
70           pv   += 4;
71         }
72         ierr = PetscLogFlops(16*nz+12);CHKERRQ(ierr);
73       }
74       row = *ajtmp++;
75     }
76     /* finished row so stick it into b->a */
77     pv = ba + 4*bi[i];
78     pj = bj + bi[i];
79     nz = bi[i+1] - bi[i];
80     for (j=0; j<nz; j++) {
81       x     = rtmp+4*pj[j];
82       pv[0] = x[0]; pv[1] = x[1]; pv[2] = x[2]; pv[3] = x[3];
83       pv   += 4;
84     }
85     /* invert diagonal block */
86     w = ba + 4*diag_offset[i];
87     ierr = Kernel_A_gets_inverse_A_2(w);CHKERRQ(ierr);
88   }
89 
90   ierr = PetscFree(rtmp);CHKERRQ(ierr);
91   ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr);
92   ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr);
93   C->factor = FACTOR_LU;
94   C->assembled = PETSC_TRUE;
95   ierr = PetscLogFlops(1.3333*8*b->mbs);CHKERRQ(ierr); /* from inverting diagonal blocks */
96   PetscFunctionReturn(0);
97 }
98 /*
99       Version for when blocks are 2 by 2 Using natural ordering
100 */
101 #undef __FUNCT__
102 #define __FUNCT__ "MatLUFactorNumeric_SeqBAIJ_2_NaturalOrdering"
103 PetscErrorCode MatLUFactorNumeric_SeqBAIJ_2_NaturalOrdering(Mat A,MatFactorInfo *info,Mat *B)
104 {
105   Mat            C = *B;
106   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ *)C->data;
107   PetscErrorCode ierr;
108   PetscInt       i,j,n = a->mbs,*bi = b->i,*bj = b->j;
109   PetscInt       *ajtmpold,*ajtmp,nz,row;
110   PetscInt       *diag_offset = b->diag,*ai=a->i,*aj=a->j,*pj;
111   MatScalar      *pv,*v,*rtmp,*pc,*w,*x;
112   MatScalar      p1,p2,p3,p4,m1,m2,m3,m4,x1,x2,x3,x4;
113   MatScalar      *ba = b->a,*aa = a->a;
114 
115   PetscFunctionBegin;
116   ierr = PetscMalloc(4*(n+1)*sizeof(MatScalar),&rtmp);CHKERRQ(ierr);
117 
118   for (i=0; i<n; i++) {
119     nz    = bi[i+1] - bi[i];
120     ajtmp = bj + bi[i];
121     for  (j=0; j<nz; j++) {
122       x = rtmp+4*ajtmp[j];
123       x[0]  = x[1]  = x[2]  = x[3]  = 0.0;
124     }
125     /* load in initial (unfactored row) */
126     nz       = ai[i+1] - ai[i];
127     ajtmpold = aj + ai[i];
128     v        = aa + 4*ai[i];
129     for (j=0; j<nz; j++) {
130       x    = rtmp+4*ajtmpold[j];
131       x[0]  = v[0];  x[1]  = v[1];  x[2]  = v[2];  x[3]  = v[3];
132       v    += 4;
133     }
134     row = *ajtmp++;
135     while (row < i) {
136       pc  = rtmp + 4*row;
137       p1  = pc[0];  p2  = pc[1];  p3  = pc[2];  p4  = pc[3];
138       if (p1 != 0.0 || p2 != 0.0 || p3 != 0.0 || p4 != 0.0) {
139         pv = ba + 4*diag_offset[row];
140         pj = bj + diag_offset[row] + 1;
141         x1  = pv[0];  x2  = pv[1];  x3  = pv[2];  x4  = pv[3];
142         pc[0] = m1 = p1*x1 + p3*x2;
143         pc[1] = m2 = p2*x1 + p4*x2;
144         pc[2] = m3 = p1*x3 + p3*x4;
145         pc[3] = m4 = p2*x3 + p4*x4;
146         nz = bi[row+1] - diag_offset[row] - 1;
147         pv += 4;
148         for (j=0; j<nz; j++) {
149           x1   = pv[0];  x2  = pv[1];   x3 = pv[2];  x4  = pv[3];
150           x    = rtmp + 4*pj[j];
151           x[0] -= m1*x1 + m3*x2;
152           x[1] -= m2*x1 + m4*x2;
153           x[2] -= m1*x3 + m3*x4;
154           x[3] -= m2*x3 + m4*x4;
155           pv   += 4;
156         }
157         ierr = PetscLogFlops(16*nz+12);CHKERRQ(ierr);
158       }
159       row = *ajtmp++;
160     }
161     /* finished row so stick it into b->a */
162     pv = ba + 4*bi[i];
163     pj = bj + bi[i];
164     nz = bi[i+1] - bi[i];
165     for (j=0; j<nz; j++) {
166       x      = rtmp+4*pj[j];
167       pv[0]  = x[0];  pv[1]  = x[1];  pv[2]  = x[2];  pv[3]  = x[3];
168       pv   += 4;
169     }
170     /* invert diagonal block */
171     w = ba + 4*diag_offset[i];
172     ierr = Kernel_A_gets_inverse_A_2(w);CHKERRQ(ierr);
173     /*Kernel_A_gets_inverse_A(bs,w,v_pivots,v_work);*/
174   }
175 
176   ierr = PetscFree(rtmp);CHKERRQ(ierr);
177   C->factor    = FACTOR_LU;
178   C->assembled = PETSC_TRUE;
179   ierr = PetscLogFlops(1.3333*8*b->mbs);CHKERRQ(ierr); /* from inverting diagonal blocks */
180   PetscFunctionReturn(0);
181 }
182 
183 /* ----------------------------------------------------------- */
184 /*
185      Version for when blocks are 1 by 1.
186 */
187 #undef __FUNCT__
188 #define __FUNCT__ "MatLUFactorNumeric_SeqBAIJ_1"
189 PetscErrorCode MatLUFactorNumeric_SeqBAIJ_1(Mat A,MatFactorInfo *info,Mat *B)
190 {
191   Mat            C = *B;
192   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ *)C->data;
193   IS             isrow = b->row,isicol = b->icol;
194   PetscErrorCode ierr;
195   PetscInt       *r,*ic,i,j,n = a->mbs,*bi = b->i,*bj = b->j;
196   PetscInt       *ajtmpold,*ajtmp,nz,row,*ai = a->i,*aj = a->j;
197   PetscInt       *diag_offset = b->diag,diag,*pj;
198   MatScalar      *pv,*v,*rtmp,multiplier,*pc;
199   MatScalar      *ba = b->a,*aa = a->a;
200 
201   PetscFunctionBegin;
202   ierr  = ISGetIndices(isrow,&r);CHKERRQ(ierr);
203   ierr  = ISGetIndices(isicol,&ic);CHKERRQ(ierr);
204   ierr  = PetscMalloc((n+1)*sizeof(MatScalar),&rtmp);CHKERRQ(ierr);
205 
206   for (i=0; i<n; i++) {
207     nz    = bi[i+1] - bi[i];
208     ajtmp = bj + bi[i];
209     for  (j=0; j<nz; j++) rtmp[ajtmp[j]] = 0.0;
210 
211     /* load in initial (unfactored row) */
212     nz       = ai[r[i]+1] - ai[r[i]];
213     ajtmpold = aj + ai[r[i]];
214     v        = aa + ai[r[i]];
215     for (j=0; j<nz; j++) rtmp[ic[ajtmpold[j]]] =  v[j];
216 
217     row = *ajtmp++;
218     while (row < i) {
219       pc = rtmp + row;
220       if (*pc != 0.0) {
221         pv         = ba + diag_offset[row];
222         pj         = bj + diag_offset[row] + 1;
223         multiplier = *pc * *pv++;
224         *pc        = multiplier;
225         nz         = bi[row+1] - diag_offset[row] - 1;
226         for (j=0; j<nz; j++) rtmp[pj[j]] -= multiplier * pv[j];
227         ierr = PetscLogFlops(1+2*nz);CHKERRQ(ierr);
228       }
229       row = *ajtmp++;
230     }
231     /* finished row so stick it into b->a */
232     pv = ba + bi[i];
233     pj = bj + bi[i];
234     nz = bi[i+1] - bi[i];
235     for (j=0; j<nz; j++) {pv[j] = rtmp[pj[j]];}
236     diag = diag_offset[i] - bi[i];
237     /* check pivot entry for current row */
238     if (pv[diag] == 0.0) {
239       SETERRQ(PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot");
240     }
241     pv[diag] = 1.0/pv[diag];
242   }
243 
244   ierr = PetscFree(rtmp);CHKERRQ(ierr);
245   ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr);
246   ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr);
247   C->factor    = FACTOR_LU;
248   C->assembled = PETSC_TRUE;
249   ierr = PetscLogFlops(C->cmap.n);CHKERRQ(ierr);
250   PetscFunctionReturn(0);
251 }
252 
253 
254 /* ----------------------------------------------------------- */
255 #undef __FUNCT__
256 #define __FUNCT__ "MatLUFactor_SeqBAIJ"
257 PetscErrorCode MatLUFactor_SeqBAIJ(Mat A,IS row,IS col,MatFactorInfo *info)
258 {
259   PetscErrorCode ierr;
260   Mat            C;
261 
262   PetscFunctionBegin;
263   ierr = MatLUFactorSymbolic(A,row,col,info,&C);CHKERRQ(ierr);
264   ierr = MatLUFactorNumeric(A,info,&C);CHKERRQ(ierr);
265   ierr = MatHeaderCopy(A,C);CHKERRQ(ierr);
266   ierr = PetscLogObjectParent(A,((Mat_SeqBAIJ*)(A->data))->icol);CHKERRQ(ierr);
267   PetscFunctionReturn(0);
268 }
269 
270 #include "src/mat/impls/sbaij/seq/sbaij.h"
271 #undef __FUNCT__
272 #define __FUNCT__ "MatCholeskyFactorNumeric_SeqBAIJ_N"
273 PetscErrorCode MatCholeskyFactorNumeric_SeqBAIJ_N(Mat A,MatFactorInfo *info,Mat *B)
274 {
275   PetscErrorCode ierr;
276   Mat            C = *B;
277   Mat_SeqBAIJ    *a=(Mat_SeqBAIJ*)A->data;
278   Mat_SeqSBAIJ   *b=(Mat_SeqSBAIJ*)C->data;
279   IS             ip=b->row;
280   PetscInt       *rip,i,j,mbs=a->mbs,bs=A->rmap.bs,*bi=b->i,*bj=b->j,*bcol;
281   PetscInt       *ai=a->i,*aj=a->j;
282   PetscInt       k,jmin,jmax,*jl,*il,col,nexti,ili,nz;
283   MatScalar      *rtmp,*ba=b->a,*bval,*aa=a->a,dk,uikdi;
284   PetscReal      zeropivot,rs,shiftnz;
285   PetscReal      shiftpd;
286   ChShift_Ctx    sctx;
287   PetscInt       newshift;
288 
289   PetscFunctionBegin;
290   if (bs > 1) {
291     if (!a->sbaijMat){
292       ierr = MatConvert(A,MATSEQSBAIJ,MAT_INITIAL_MATRIX,&a->sbaijMat);CHKERRQ(ierr);
293     }
294     ierr = (a->sbaijMat)->ops->choleskyfactornumeric(a->sbaijMat,info,B);CHKERRQ(ierr);
295     ierr = MatDestroy(a->sbaijMat);CHKERRQ(ierr);
296     a->sbaijMat = PETSC_NULL;
297     PetscFunctionReturn(0);
298   }
299 
300   /* initialization */
301   shiftnz   = info->shiftnz;
302   shiftpd   = info->shiftpd;
303   zeropivot = info->zeropivot;
304 
305   ierr  = ISGetIndices(ip,&rip);CHKERRQ(ierr);
306   nz   = (2*mbs+1)*sizeof(PetscInt)+mbs*sizeof(MatScalar);
307   ierr = PetscMalloc(nz,&il);CHKERRQ(ierr);
308   jl   = il + mbs;
309   rtmp = (MatScalar*)(jl + mbs);
310 
311   sctx.shift_amount = 0;
312   sctx.nshift       = 0;
313   do {
314     sctx.chshift = PETSC_FALSE;
315     for (i=0; i<mbs; i++) {
316       rtmp[i] = 0.0; jl[i] = mbs; il[0] = 0;
317     }
318 
319     for (k = 0; k<mbs; k++){
320       bval = ba + bi[k];
321       /* initialize k-th row by the perm[k]-th row of A */
322       jmin = ai[rip[k]]; jmax = ai[rip[k]+1];
323       for (j = jmin; j < jmax; j++){
324         col = rip[aj[j]];
325         if (col >= k){ /* only take upper triangular entry */
326           rtmp[col] = aa[j];
327           *bval++  = 0.0; /* for in-place factorization */
328         }
329       }
330 
331       /* shift the diagonal of the matrix */
332       if (sctx.nshift) rtmp[k] += sctx.shift_amount;
333 
334       /* modify k-th row by adding in those rows i with U(i,k)!=0 */
335       dk = rtmp[k];
336       i = jl[k]; /* first row to be added to k_th row  */
337 
338       while (i < k){
339         nexti = jl[i]; /* next row to be added to k_th row */
340 
341         /* compute multiplier, update diag(k) and U(i,k) */
342         ili = il[i];  /* index of first nonzero element in U(i,k:bms-1) */
343         uikdi = - ba[ili]*ba[bi[i]];  /* diagonal(k) */
344         dk += uikdi*ba[ili];
345         ba[ili] = uikdi; /* -U(i,k) */
346 
347         /* add multiple of row i to k-th row */
348         jmin = ili + 1; jmax = bi[i+1];
349         if (jmin < jmax){
350           for (j=jmin; j<jmax; j++) rtmp[bj[j]] += uikdi*ba[j];
351           /* update il and jl for row i */
352           il[i] = jmin;
353           j = bj[jmin]; jl[i] = jl[j]; jl[j] = i;
354         }
355         i = nexti;
356       }
357 
358       /* shift the diagonals when zero pivot is detected */
359       /* compute rs=sum of abs(off-diagonal) */
360       rs   = 0.0;
361       jmin = bi[k]+1;
362       nz   = bi[k+1] - jmin;
363       if (nz){
364         bcol = bj + jmin;
365         while (nz--){
366           rs += PetscAbsScalar(rtmp[*bcol]);
367           bcol++;
368         }
369       }
370 
371       sctx.rs = rs;
372       sctx.pv = dk;
373       ierr = MatCholeskyCheckShift_inline(info,sctx,k,newshift);CHKERRQ(ierr);
374       if (newshift == 1) break;
375 
376       /* copy data into U(k,:) */
377       ba[bi[k]] = 1.0/dk; /* U(k,k) */
378       jmin = bi[k]+1; jmax = bi[k+1];
379       if (jmin < jmax) {
380         for (j=jmin; j<jmax; j++){
381           col = bj[j]; ba[j] = rtmp[col]; rtmp[col] = 0.0;
382         }
383         /* add the k-th row into il and jl */
384         il[k] = jmin;
385         i = bj[jmin]; jl[k] = jl[i]; jl[i] = k;
386       }
387     }
388   } while (sctx.chshift);
389   ierr = PetscFree(il);CHKERRQ(ierr);
390 
391   ierr = ISRestoreIndices(ip,&rip);CHKERRQ(ierr);
392   C->factor       = FACTOR_CHOLESKY;
393   C->assembled    = PETSC_TRUE;
394   C->preallocated = PETSC_TRUE;
395   ierr = PetscLogFlops(C->rmap.N);CHKERRQ(ierr);
396   if (sctx.nshift){
397     if (shiftnz) {
398       ierr = PetscInfo2(0,"number of shiftnz tries %D, shift_amount %G\n",sctx.nshift,sctx.shift_amount);CHKERRQ(ierr);
399     } else if (shiftpd) {
400       ierr = PetscInfo2(0,"number of shiftpd tries %D, shift_amount %G\n",sctx.nshift,sctx.shift_amount);CHKERRQ(ierr);
401     }
402   }
403   PetscFunctionReturn(0);
404 }
405 
406 #undef __FUNCT__
407 #define __FUNCT__ "MatCholeskyFactorNumeric_SeqBAIJ_N_NaturalOrdering"
408 PetscErrorCode MatCholeskyFactorNumeric_SeqBAIJ_N_NaturalOrdering(Mat A,MatFactorInfo *info,Mat *fact)
409 {
410   Mat            C = *fact;
411   Mat_SeqBAIJ    *a=(Mat_SeqBAIJ*)A->data;
412   Mat_SeqSBAIJ   *b=(Mat_SeqSBAIJ*)C->data;
413   PetscErrorCode ierr;
414   PetscInt       i,j,am=a->mbs;
415   PetscInt       *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j;
416   PetscInt       k,jmin,*jl,*il,nexti,ili,*acol,*bcol,nz;
417   MatScalar      *rtmp,*ba=b->a,*aa=a->a,dk,uikdi,*aval,*bval;
418   PetscReal      zeropivot,rs,shiftnz;
419   PetscReal      shiftpd;
420   ChShift_Ctx    sctx;
421   PetscInt       newshift;
422 
423   PetscFunctionBegin;
424   /* initialization */
425   shiftnz   = info->shiftnz;
426   shiftpd   = info->shiftpd;
427   zeropivot = info->zeropivot;
428 
429   nz   = (2*am+1)*sizeof(PetscInt)+am*sizeof(MatScalar);
430   ierr = PetscMalloc(nz,&il);CHKERRQ(ierr);
431   jl   = il + am;
432   rtmp = (MatScalar*)(jl + am);
433 
434   sctx.shift_amount = 0;
435   sctx.nshift       = 0;
436   do {
437     sctx.chshift = PETSC_FALSE;
438     for (i=0; i<am; i++) {
439       rtmp[i] = 0.0; jl[i] = am; il[0] = 0;
440     }
441 
442     for (k = 0; k<am; k++){
443     /* initialize k-th row with elements nonzero in row perm(k) of A */
444       nz   = ai[k+1] - ai[k];
445       acol = aj + ai[k];
446       aval = aa + ai[k];
447       bval = ba + bi[k];
448       while (nz -- ){
449         if (*acol < k) { /* skip lower triangular entries */
450           acol++; aval++;
451         } else {
452           rtmp[*acol++] = *aval++;
453           *bval++       = 0.0; /* for in-place factorization */
454         }
455       }
456 
457       /* shift the diagonal of the matrix */
458       if (sctx.nshift) rtmp[k] += sctx.shift_amount;
459 
460       /* modify k-th row by adding in those rows i with U(i,k)!=0 */
461       dk = rtmp[k];
462       i  = jl[k]; /* first row to be added to k_th row  */
463 
464       while (i < k){
465         nexti = jl[i]; /* next row to be added to k_th row */
466         /* compute multiplier, update D(k) and U(i,k) */
467         ili   = il[i];  /* index of first nonzero element in U(i,k:bms-1) */
468         uikdi = - ba[ili]*ba[bi[i]];
469         dk   += uikdi*ba[ili];
470         ba[ili] = uikdi; /* -U(i,k) */
471 
472         /* add multiple of row i to k-th row ... */
473         jmin = ili + 1;
474         nz   = bi[i+1] - jmin;
475         if (nz > 0){
476           bcol = bj + jmin;
477           bval = ba + jmin;
478           while (nz --) rtmp[*bcol++] += uikdi*(*bval++);
479           /* update il and jl for i-th row */
480           il[i] = jmin;
481           j = bj[jmin]; jl[i] = jl[j]; jl[j] = i;
482         }
483         i = nexti;
484       }
485 
486       /* shift the diagonals when zero pivot is detected */
487       /* compute rs=sum of abs(off-diagonal) */
488       rs   = 0.0;
489       jmin = bi[k]+1;
490       nz   = bi[k+1] - jmin;
491       if (nz){
492         bcol = bj + jmin;
493         while (nz--){
494           rs += PetscAbsScalar(rtmp[*bcol]);
495           bcol++;
496         }
497       }
498 
499       sctx.rs = rs;
500       sctx.pv = dk;
501       ierr = MatCholeskyCheckShift_inline(info,sctx,k,newshift);CHKERRQ(ierr);
502       if (newshift == 1) break;    /* sctx.shift_amount is updated */
503 
504       /* copy data into U(k,:) */
505       ba[bi[k]] = 1.0/dk;
506       jmin      = bi[k]+1;
507       nz        = bi[k+1] - jmin;
508       if (nz){
509         bcol = bj + jmin;
510         bval = ba + jmin;
511         while (nz--){
512           *bval++       = rtmp[*bcol];
513           rtmp[*bcol++] = 0.0;
514         }
515         /* add k-th row into il and jl */
516         il[k] = jmin;
517         i = bj[jmin]; jl[k] = jl[i]; jl[i] = k;
518       }
519     }
520   } while (sctx.chshift);
521   ierr = PetscFree(il);CHKERRQ(ierr);
522 
523   C->factor       = FACTOR_CHOLESKY;
524   C->assembled    = PETSC_TRUE;
525   C->preallocated = PETSC_TRUE;
526   ierr = PetscLogFlops(C->rmap.N);CHKERRQ(ierr);
527     if (sctx.nshift){
528     if (shiftnz) {
529       ierr = PetscInfo2(0,"number of shiftnz tries %D, shift_amount %G\n",sctx.nshift,sctx.shift_amount);CHKERRQ(ierr);
530     } else if (shiftpd) {
531       ierr = PetscInfo2(0,"number of shiftpd tries %D, shift_amount %G\n",sctx.nshift,sctx.shift_amount);CHKERRQ(ierr);
532     }
533   }
534   PetscFunctionReturn(0);
535 }
536 
537 #include "petscbt.h"
538 #include "src/mat/utils/freespace.h"
539 #undef __FUNCT__
540 #define __FUNCT__ "MatICCFactorSymbolic_SeqBAIJ"
541 PetscErrorCode MatICCFactorSymbolic_SeqBAIJ(Mat A,IS perm,MatFactorInfo *info,Mat *fact)
542 {
543   Mat_SeqBAIJ        *a = (Mat_SeqBAIJ*)A->data;
544   Mat_SeqSBAIJ       *b;
545   Mat                B;
546   PetscErrorCode     ierr;
547   PetscTruth         perm_identity;
548   PetscInt           reallocs=0,*rip,i,*ai=a->i,*aj=a->j,am=a->mbs,bs=A->rmap.bs,*ui;
549   PetscInt           jmin,jmax,nzk,k,j,*jl,prow,*il,nextprow;
550   PetscInt           nlnk,*lnk,*lnk_lvl=PETSC_NULL,ncols,ncols_upper,*cols,*cols_lvl,*uj,**uj_ptr,**uj_lvl_ptr;
551   PetscReal          fill=info->fill,levels=info->levels;
552   PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL;
553   PetscFreeSpaceList free_space_lvl=PETSC_NULL,current_space_lvl=PETSC_NULL;
554   PetscBT            lnkbt;
555 
556   PetscFunctionBegin;
557   if (bs > 1){
558     if (!a->sbaijMat){
559       ierr = MatConvert(A,MATSEQSBAIJ,MAT_INITIAL_MATRIX,&a->sbaijMat);CHKERRQ(ierr);
560     }
561     ierr = MatICCFactorSymbolic(a->sbaijMat,perm,info,fact);CHKERRQ(ierr);
562     B = *fact;
563     B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqBAIJ_N;
564     PetscFunctionReturn(0);
565   }
566 
567   ierr = ISIdentity(perm,&perm_identity);CHKERRQ(ierr);
568   ierr = ISGetIndices(perm,&rip);CHKERRQ(ierr);
569 
570   /* special case that simply copies fill pattern */
571   if (!levels && perm_identity) {
572     ierr = MatMarkDiagonal_SeqBAIJ(A);CHKERRQ(ierr);
573     ierr = PetscMalloc((am+1)*sizeof(PetscInt),&ui);CHKERRQ(ierr);
574     for (i=0; i<am; i++) {
575       ui[i] = ai[i+1] - a->diag[i]; /* ui: rowlengths - changes when !perm_identity */
576     }
577     ierr = MatCreate(PETSC_COMM_SELF,fact);CHKERRQ(ierr);
578     ierr = MatSetSizes(*fact,am,am,am,am);CHKERRQ(ierr);
579     B = *fact;
580     ierr = MatSetType(B,MATSEQSBAIJ);CHKERRQ(ierr);
581     ierr = MatSeqSBAIJSetPreallocation(B,1,0,ui);CHKERRQ(ierr);
582 
583     b  = (Mat_SeqSBAIJ*)B->data;
584     uj = b->j;
585     for (i=0; i<am; i++) {
586       aj = a->j + a->diag[i];
587       for (j=0; j<ui[i]; j++){
588         *uj++ = *aj++;
589       }
590       b->ilen[i] = ui[i];
591     }
592     ierr = PetscFree(ui);CHKERRQ(ierr);
593     ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
594     ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
595 
596     B->ops->solve                 = MatSolve_SeqSBAIJ_1_NaturalOrdering;
597     B->ops->solvetranspose        = MatSolve_SeqSBAIJ_1_NaturalOrdering;
598     B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqBAIJ_N_NaturalOrdering;
599     PetscFunctionReturn(0);
600   }
601 
602   /* initialization */
603   ierr  = PetscMalloc((am+1)*sizeof(PetscInt),&ui);CHKERRQ(ierr);
604   ui[0] = 0;
605   ierr  = PetscMalloc((2*am+1)*sizeof(PetscInt),&cols_lvl);CHKERRQ(ierr);
606 
607   /* jl: linked list for storing indices of the pivot rows
608      il: il[i] points to the 1st nonzero entry of U(i,k:am-1) */
609   ierr = PetscMalloc((2*am+1)*sizeof(PetscInt)+2*am*sizeof(PetscInt*),&jl);CHKERRQ(ierr);
610   il         = jl + am;
611   uj_ptr     = (PetscInt**)(il + am);
612   uj_lvl_ptr = (PetscInt**)(uj_ptr + am);
613   for (i=0; i<am; i++){
614     jl[i] = am; il[i] = 0;
615   }
616 
617   /* create and initialize a linked list for storing column indices of the active row k */
618   nlnk = am + 1;
619   ierr = PetscIncompleteLLCreate(am,am,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr);
620 
621   /* initial FreeSpace size is fill*(ai[am]+1) */
622   ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+1)),&free_space);CHKERRQ(ierr);
623   current_space = free_space;
624   ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+1)),&free_space_lvl);CHKERRQ(ierr);
625   current_space_lvl = free_space_lvl;
626 
627   for (k=0; k<am; k++){  /* for each active row k */
628     /* initialize lnk by the column indices of row rip[k] of A */
629     nzk   = 0;
630     ncols = ai[rip[k]+1] - ai[rip[k]];
631     ncols_upper = 0;
632     cols        = cols_lvl + am;
633     for (j=0; j<ncols; j++){
634       i = rip[*(aj + ai[rip[k]] + j)];
635       if (i >= k){ /* only take upper triangular entry */
636         cols[ncols_upper] = i;
637         cols_lvl[ncols_upper] = -1;  /* initialize level for nonzero entries */
638         ncols_upper++;
639       }
640     }
641     ierr = PetscIncompleteLLAdd(ncols_upper,cols,levels,cols_lvl,am,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr);
642     nzk += nlnk;
643 
644     /* update lnk by computing fill-in for each pivot row to be merged in */
645     prow = jl[k]; /* 1st pivot row */
646 
647     while (prow < k){
648       nextprow = jl[prow];
649 
650       /* merge prow into k-th row */
651       jmin = il[prow] + 1;  /* index of the 2nd nzero entry in U(prow,k:am-1) */
652       jmax = ui[prow+1];
653       ncols = jmax-jmin;
654       i     = jmin - ui[prow];
655       cols = uj_ptr[prow] + i; /* points to the 2nd nzero entry in U(prow,k:am-1) */
656       for (j=0; j<ncols; j++) cols_lvl[j] = *(uj_lvl_ptr[prow] + i + j);
657       ierr = PetscIncompleteLLAddSorted(ncols,cols,levels,cols_lvl,am,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr);
658       nzk += nlnk;
659 
660       /* update il and jl for prow */
661       if (jmin < jmax){
662         il[prow] = jmin;
663         j = *cols; jl[prow] = jl[j]; jl[j] = prow;
664       }
665       prow = nextprow;
666     }
667 
668     /* if free space is not available, make more free space */
669     if (current_space->local_remaining<nzk) {
670       i = am - k + 1; /* num of unfactored rows */
671       i = PetscMin(i*nzk, i*(i-1)); /* i*nzk, i*(i-1): estimated and max additional space needed */
672       ierr = PetscFreeSpaceGet(i,&current_space);CHKERRQ(ierr);
673       ierr = PetscFreeSpaceGet(i,&current_space_lvl);CHKERRQ(ierr);
674       reallocs++;
675     }
676 
677     /* copy data into free_space and free_space_lvl, then initialize lnk */
678     ierr = PetscIncompleteLLClean(am,am,nzk,lnk,lnk_lvl,current_space->array,current_space_lvl->array,lnkbt);CHKERRQ(ierr);
679 
680     /* add the k-th row into il and jl */
681     if (nzk-1 > 0){
682       i = current_space->array[1]; /* col value of the first nonzero element in U(k, k+1:am-1) */
683       jl[k] = jl[i]; jl[i] = k;
684       il[k] = ui[k] + 1;
685     }
686     uj_ptr[k]     = current_space->array;
687     uj_lvl_ptr[k] = current_space_lvl->array;
688 
689     current_space->array           += nzk;
690     current_space->local_used      += nzk;
691     current_space->local_remaining -= nzk;
692 
693     current_space_lvl->array           += nzk;
694     current_space_lvl->local_used      += nzk;
695     current_space_lvl->local_remaining -= nzk;
696 
697     ui[k+1] = ui[k] + nzk;
698   }
699 
700 #if defined(PETSC_USE_INFO)
701   if (ai[am] != 0) {
702     PetscReal af = ((PetscReal)(2*ui[am]-am))/((PetscReal)ai[am]);
703     ierr = PetscInfo3(A,"Reallocs %D Fill ratio:given %G needed %G\n",reallocs,fill,af);CHKERRQ(ierr);
704     ierr = PetscInfo1(A,"Run with -pc_factor_fill %G or use \n",af);CHKERRQ(ierr);
705     ierr = PetscInfo1(A,"PCFactorSetFill(pc,%G) for best performance.\n",af);CHKERRQ(ierr);
706   } else {
707     ierr = PetscInfo(A,"Empty matrix.\n");CHKERRQ(ierr);
708   }
709 #endif
710 
711   ierr = ISRestoreIndices(perm,&rip);CHKERRQ(ierr);
712   ierr = PetscFree(jl);CHKERRQ(ierr);
713   ierr = PetscFree(cols_lvl);CHKERRQ(ierr);
714 
715   /* destroy list of free space and other temporary array(s) */
716   ierr = PetscMalloc((ui[am]+1)*sizeof(PetscInt),&uj);CHKERRQ(ierr);
717   ierr = PetscFreeSpaceContiguous(&free_space,uj);CHKERRQ(ierr);
718   ierr = PetscIncompleteLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
719   ierr = PetscFreeSpaceDestroy(free_space_lvl);CHKERRQ(ierr);
720 
721   /* put together the new matrix in MATSEQSBAIJ format */
722   ierr = MatCreate(PETSC_COMM_SELF,fact);CHKERRQ(ierr);
723   ierr = MatSetSizes(*fact,am,am,am,am);CHKERRQ(ierr);
724   B = *fact;
725   ierr = MatSetType(B,MATSEQSBAIJ);CHKERRQ(ierr);
726   ierr = MatSeqSBAIJSetPreallocation(B,1,MAT_SKIP_ALLOCATION,PETSC_NULL);CHKERRQ(ierr);
727 
728   b = (Mat_SeqSBAIJ*)B->data;
729   b->singlemalloc = PETSC_FALSE;
730   b->freedata     = PETSC_TRUE;
731   ierr = PetscMalloc((ui[am]+1)*sizeof(MatScalar),&b->a);CHKERRQ(ierr);
732   b->j    = uj;
733   b->i    = ui;
734   b->diag = 0;
735   b->ilen = 0;
736   b->imax = 0;
737   b->row  = perm;
738   b->pivotinblocks = PETSC_FALSE; /* need to get from MatFactorInfo */
739   ierr    = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
740   b->icol = perm;
741   ierr    = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
742   ierr    = PetscMalloc((am+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
743   ierr    = PetscLogObjectMemory(B,(ui[am]-am)*(sizeof(PetscInt)+sizeof(MatScalar)));CHKERRQ(ierr);
744   b->maxnz = b->nz = ui[am];
745 
746   B->factor                 = FACTOR_CHOLESKY;
747   B->info.factor_mallocs    = reallocs;
748   B->info.fill_ratio_given  = fill;
749   if (ai[am] != 0) {
750     B->info.fill_ratio_needed = ((PetscReal)ui[am])/((PetscReal)ai[am]);
751   } else {
752     B->info.fill_ratio_needed = 0.0;
753   }
754   if (perm_identity){
755     B->ops->solve           = MatSolve_SeqSBAIJ_1_NaturalOrdering;
756     B->ops->solvetranspose  = MatSolve_SeqSBAIJ_1_NaturalOrdering;
757     B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqBAIJ_N_NaturalOrdering;
758   } else {
759     (*fact)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqBAIJ_N;
760   }
761   PetscFunctionReturn(0);
762 }
763 
764 #undef __FUNCT__
765 #define __FUNCT__ "MatCholeskyFactorSymbolic_SeqBAIJ"
766 PetscErrorCode MatCholeskyFactorSymbolic_SeqBAIJ(Mat A,IS perm,MatFactorInfo *info,Mat *fact)
767 {
768   Mat_SeqBAIJ        *a = (Mat_SeqBAIJ*)A->data;
769   Mat_SeqSBAIJ       *b;
770   Mat                B;
771   PetscErrorCode     ierr;
772   PetscTruth         perm_identity;
773   PetscReal          fill = info->fill;
774   PetscInt           *rip,*riip,i,mbs=a->mbs,bs=A->rmap.bs,*ai=a->i,*aj=a->j,reallocs=0,prow;
775   PetscInt           *jl,jmin,jmax,nzk,*ui,k,j,*il,nextprow;
776   PetscInt           nlnk,*lnk,ncols,ncols_upper,*cols,*uj,**ui_ptr,*uj_ptr;
777   PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL;
778   PetscBT            lnkbt;
779   IS                 iperm;
780 
781   PetscFunctionBegin;
782   if (bs > 1) { /* convert to seqsbaij */
783     if (!a->sbaijMat){
784       ierr = MatConvert(A,MATSEQSBAIJ,MAT_INITIAL_MATRIX,&a->sbaijMat);CHKERRQ(ierr);
785     }
786     ierr = MatCholeskyFactorSymbolic(a->sbaijMat,perm,info,fact);CHKERRQ(ierr);
787     B    = *fact;
788     B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqBAIJ_N;
789     PetscFunctionReturn(0);
790   }
791 
792   /* check whether perm is the identity mapping */
793   ierr = ISIdentity(perm,&perm_identity);CHKERRQ(ierr);
794   ierr = ISGetIndices(perm,&rip);CHKERRQ(ierr);
795 
796   if (!perm_identity){
797     /* check if perm is symmetric! */
798     ierr = ISInvertPermutation(perm,PETSC_DECIDE,&iperm);CHKERRQ(ierr);
799     ierr = ISGetIndices(iperm,&riip);CHKERRQ(ierr);
800     for (i=0; i<mbs; i++) {
801       if (rip[i] != riip[i]) SETERRQ(PETSC_ERR_ARG_INCOMP,"Non-symmetric permutation, must use symmetric permutation");
802     }
803     ierr = ISRestoreIndices(iperm,&riip);CHKERRQ(ierr);
804     ierr = ISDestroy(iperm);CHKERRQ(ierr);
805   }
806 
807   /* initialization */
808   ierr  = PetscMalloc((mbs+1)*sizeof(PetscInt),&ui);CHKERRQ(ierr);
809   ui[0] = 0;
810 
811   /* jl: linked list for storing indices of the pivot rows
812      il: il[i] points to the 1st nonzero entry of U(i,k:mbs-1) */
813   ierr = PetscMalloc((3*mbs+1)*sizeof(PetscInt)+mbs*sizeof(PetscInt*),&jl);CHKERRQ(ierr);
814   il     = jl + mbs;
815   cols   = il + mbs;
816   ui_ptr = (PetscInt**)(cols + mbs);
817   for (i=0; i<mbs; i++){
818     jl[i] = mbs; il[i] = 0;
819   }
820 
821   /* create and initialize a linked list for storing column indices of the active row k */
822   nlnk = mbs + 1;
823   ierr = PetscLLCreate(mbs,mbs,nlnk,lnk,lnkbt);CHKERRQ(ierr);
824 
825   /* initial FreeSpace size is fill*(ai[mbs]+1) */
826   ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[mbs]+1)),&free_space);CHKERRQ(ierr);
827   current_space = free_space;
828 
829   for (k=0; k<mbs; k++){  /* for each active row k */
830     /* initialize lnk by the column indices of row rip[k] of A */
831     nzk   = 0;
832     ncols = ai[rip[k]+1] - ai[rip[k]];
833     ncols_upper = 0;
834     for (j=0; j<ncols; j++){
835       i = rip[*(aj + ai[rip[k]] + j)];
836       if (i >= k){ /* only take upper triangular entry */
837         cols[ncols_upper] = i;
838         ncols_upper++;
839       }
840     }
841     ierr = PetscLLAdd(ncols_upper,cols,mbs,nlnk,lnk,lnkbt);CHKERRQ(ierr);
842     nzk += nlnk;
843 
844     /* update lnk by computing fill-in for each pivot row to be merged in */
845     prow = jl[k]; /* 1st pivot row */
846 
847     while (prow < k){
848       nextprow = jl[prow];
849       /* merge prow into k-th row */
850       jmin = il[prow] + 1;  /* index of the 2nd nzero entry in U(prow,k:mbs-1) */
851       jmax = ui[prow+1];
852       ncols = jmax-jmin;
853       uj_ptr = ui_ptr[prow] + jmin - ui[prow]; /* points to the 2nd nzero entry in U(prow,k:mbs-1) */
854       ierr = PetscLLAddSorted(ncols,uj_ptr,mbs,nlnk,lnk,lnkbt);CHKERRQ(ierr);
855       nzk += nlnk;
856 
857       /* update il and jl for prow */
858       if (jmin < jmax){
859         il[prow] = jmin;
860         j = *uj_ptr; jl[prow] = jl[j]; jl[j] = prow;
861       }
862       prow = nextprow;
863     }
864 
865     /* if free space is not available, make more free space */
866     if (current_space->local_remaining<nzk) {
867       i = mbs - k + 1; /* num of unfactored rows */
868       i = PetscMin(i*nzk, i*(i-1)); /* i*nzk, i*(i-1): estimated and max additional space needed */
869       ierr = PetscFreeSpaceGet(i,&current_space);CHKERRQ(ierr);
870       reallocs++;
871     }
872 
873     /* copy data into free space, then initialize lnk */
874     ierr = PetscLLClean(mbs,mbs,nzk,lnk,current_space->array,lnkbt);CHKERRQ(ierr);
875 
876     /* add the k-th row into il and jl */
877     if (nzk-1 > 0){
878       i = current_space->array[1]; /* col value of the first nonzero element in U(k, k+1:mbs-1) */
879       jl[k] = jl[i]; jl[i] = k;
880       il[k] = ui[k] + 1;
881     }
882     ui_ptr[k] = current_space->array;
883     current_space->array           += nzk;
884     current_space->local_used      += nzk;
885     current_space->local_remaining -= nzk;
886 
887     ui[k+1] = ui[k] + nzk;
888   }
889 
890 #if defined(PETSC_USE_INFO)
891   if (ai[mbs] != 0) {
892     PetscReal af = ((PetscReal)ui[mbs])/((PetscReal)ai[mbs]);
893     ierr = PetscInfo3(A,"Reallocs %D Fill ratio:given %G needed %G\n",reallocs,fill,af);CHKERRQ(ierr);
894     ierr = PetscInfo1(A,"Run with -pc_factor_fill %G or use \n",af);CHKERRQ(ierr);
895     ierr = PetscInfo1(A,"PCFactorSetFill(pc,%G) for best performance.\n",af);CHKERRQ(ierr);
896   } else {
897     ierr = PetscInfo(A,"Empty matrix.\n");CHKERRQ(ierr);
898   }
899 #endif
900 
901   ierr = ISRestoreIndices(perm,&rip);CHKERRQ(ierr);
902   ierr = PetscFree(jl);CHKERRQ(ierr);
903 
904   /* destroy list of free space and other temporary array(s) */
905   ierr = PetscMalloc((ui[mbs]+1)*sizeof(PetscInt),&uj);CHKERRQ(ierr);
906   ierr = PetscFreeSpaceContiguous(&free_space,uj);CHKERRQ(ierr);
907   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
908 
909   /* put together the new matrix in MATSEQSBAIJ format */
910   ierr = MatCreate(PETSC_COMM_SELF,fact);CHKERRQ(ierr);
911   ierr = MatSetSizes(*fact,mbs,mbs,mbs,mbs);CHKERRQ(ierr);
912   B    = *fact;
913   ierr = MatSetType(B,MATSEQSBAIJ);CHKERRQ(ierr);
914   ierr = MatSeqSBAIJSetPreallocation(B,bs,MAT_SKIP_ALLOCATION,PETSC_NULL);CHKERRQ(ierr);
915 
916   b = (Mat_SeqSBAIJ*)B->data;
917   b->singlemalloc = PETSC_FALSE;
918   b->freedata     = PETSC_TRUE;
919   ierr = PetscMalloc((ui[mbs]+1)*sizeof(MatScalar),&b->a);CHKERRQ(ierr);
920   b->j    = uj;
921   b->i    = ui;
922   b->diag = 0;
923   b->ilen = 0;
924   b->imax = 0;
925   b->row  = perm;
926   b->pivotinblocks = PETSC_FALSE; /* need to get from MatFactorInfo */
927   ierr    = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
928   b->icol = perm;
929   ierr    = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
930   ierr    = PetscMalloc((mbs+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
931   ierr    = PetscLogObjectMemory(B,(ui[mbs]-mbs)*(sizeof(PetscInt)+sizeof(MatScalar)));CHKERRQ(ierr);
932   b->maxnz = b->nz = ui[mbs];
933 
934   B->factor                 = FACTOR_CHOLESKY;
935   B->info.factor_mallocs    = reallocs;
936   B->info.fill_ratio_given  = fill;
937   if (ai[mbs] != 0) {
938     B->info.fill_ratio_needed = ((PetscReal)ui[mbs])/((PetscReal)ai[mbs]);
939   } else {
940     B->info.fill_ratio_needed = 0.0;
941   }
942   if (perm_identity){
943     B->ops->solve           = MatSolve_SeqSBAIJ_1_NaturalOrdering;
944     B->ops->solvetranspose  = MatSolve_SeqSBAIJ_1_NaturalOrdering;
945     B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqBAIJ_N_NaturalOrdering;
946   } else {
947     B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqBAIJ_N;
948   }
949   PetscFunctionReturn(0);
950 }
951