xref: /petsc/src/mat/impls/sbaij/seq/sbaijfact.c (revision 6024bd2c273c6b4ac03c2bd49d1968cdcd0f9df3)
1 #define PETSCMAT_DLL
2 
3 #include "src/mat/impls/baij/seq/baij.h"
4 #include "src/mat/impls/sbaij/seq/sbaij.h"
5 #include "src/inline/ilu.h"
6 #include "include/petscis.h"
7 
8 #if !defined(PETSC_USE_COMPLEX)
9 /*
10   input:
11    F -- numeric factor
12   output:
13    nneg, nzero, npos: matrix inertia
14 */
15 
16 #undef __FUNCT__
17 #define __FUNCT__ "MatGetInertia_SeqSBAIJ"
18 PetscErrorCode MatGetInertia_SeqSBAIJ(Mat F,PetscInt *nneig,PetscInt *nzero,PetscInt *npos)
19 {
20   Mat_SeqSBAIJ *fact_ptr = (Mat_SeqSBAIJ*)F->data;
21   PetscScalar  *dd = fact_ptr->a;
22   PetscInt     mbs=fact_ptr->mbs,bs=F->bs,i,nneig_tmp,npos_tmp,*fi = fact_ptr->i;
23 
24   PetscFunctionBegin;
25   if (bs != 1) SETERRQ1(PETSC_ERR_SUP,"No support for bs: %D >1 yet",bs);
26   nneig_tmp = 0; npos_tmp = 0;
27   for (i=0; i<mbs; i++){
28     if (PetscRealPart(dd[*fi]) > 0.0){
29       npos_tmp++;
30     } else if (PetscRealPart(dd[*fi]) < 0.0){
31       nneig_tmp++;
32     }
33     fi++;
34   }
35   if (nneig) *nneig = nneig_tmp;
36   if (npos)  *npos  = npos_tmp;
37   if (nzero) *nzero = mbs - nneig_tmp - npos_tmp;
38 
39   PetscFunctionReturn(0);
40 }
41 #endif /* !defined(PETSC_USE_COMPLEX) */
42 
43 /*
44   Symbolic U^T*D*U factorization for SBAIJ format. Modified from SSF of YSMP.
45   Use Modified Sparse Row (MSR) storage for u and ju. See page 85, "Iterative Methods ..." by Saad.
46 */
47 #undef __FUNCT__
48 #define __FUNCT__ "MatCholeskyFactorSymbolic_SeqSBAIJ_MSR"
49 PetscErrorCode MatCholeskyFactorSymbolic_SeqSBAIJ_MSR(Mat A,IS perm,MatFactorInfo *info,Mat *B)
50 {
51   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)A->data,*b;
52   PetscErrorCode ierr;
53   PetscInt       *rip,i,mbs = a->mbs,*ai,*aj;
54   PetscInt       *jutmp,bs = A->bs,bs2=a->bs2;
55   PetscInt       m,reallocs = 0,prow;
56   PetscInt       *jl,*q,jmin,jmax,juidx,nzk,qm,*iu,*ju,k,j,vj,umax,maxadd;
57   PetscReal      f = info->fill;
58   PetscTruth     perm_identity;
59 
60   PetscFunctionBegin;
61   /* check whether perm is the identity mapping */
62   ierr = ISIdentity(perm,&perm_identity);CHKERRQ(ierr);
63   ierr = ISGetIndices(perm,&rip);CHKERRQ(ierr);
64 
65   if (perm_identity){ /* without permutation */
66     a->permute = PETSC_FALSE;
67     ai = a->i; aj = a->j;
68   } else {            /* non-trivial permutation */
69     a->permute = PETSC_TRUE;
70     ierr = MatReorderingSeqSBAIJ(A,perm);CHKERRQ(ierr);
71     ai = a->inew; aj = a->jnew;
72   }
73 
74   /* initialization */
75   ierr  = PetscMalloc((mbs+1)*sizeof(PetscInt),&iu);CHKERRQ(ierr);
76   umax  = (PetscInt)(f*ai[mbs] + 1); umax += mbs + 1;
77   ierr  = PetscMalloc(umax*sizeof(PetscInt),&ju);CHKERRQ(ierr);
78   iu[0] = mbs+1;
79   juidx = mbs + 1; /* index for ju */
80   ierr  = PetscMalloc(2*mbs*sizeof(PetscInt),&jl);CHKERRQ(ierr); /* linked list for pivot row */
81   q     = jl + mbs;   /* linked list for col index */
82   for (i=0; i<mbs; i++){
83     jl[i] = mbs;
84     q[i] = 0;
85   }
86 
87   /* for each row k */
88   for (k=0; k<mbs; k++){
89     for (i=0; i<mbs; i++) q[i] = 0;  /* to be removed! */
90     nzk  = 0; /* num. of nz blocks in k-th block row with diagonal block excluded */
91     q[k] = mbs;
92     /* initialize nonzero structure of k-th row to row rip[k] of A */
93     jmin = ai[rip[k]] +1; /* exclude diag[k] */
94     jmax = ai[rip[k]+1];
95     for (j=jmin; j<jmax; j++){
96       vj = rip[aj[j]]; /* col. value */
97       if(vj > k){
98         qm = k;
99         do {
100           m  = qm; qm = q[m];
101         } while(qm < vj);
102         if (qm == vj) {
103           SETERRQ(PETSC_ERR_PLIB,"Duplicate entry in A\n");
104         }
105         nzk++;
106         q[m]  = vj;
107         q[vj] = qm;
108       } /* if(vj > k) */
109     } /* for (j=jmin; j<jmax; j++) */
110 
111     /* modify nonzero structure of k-th row by computing fill-in
112        for each row i to be merged in */
113     prow = k;
114     prow = jl[prow]; /* next pivot row (== mbs for symbolic factorization) */
115 
116     while (prow < k){
117       /* merge row prow into k-th row */
118       jmin = iu[prow] + 1; jmax = iu[prow+1];
119       qm = k;
120       for (j=jmin; j<jmax; j++){
121         vj = ju[j];
122         do {
123           m = qm; qm = q[m];
124         } while (qm < vj);
125         if (qm != vj){
126          nzk++; q[m] = vj; q[vj] = qm; qm = vj;
127         }
128       }
129       prow = jl[prow]; /* next pivot row */
130     }
131 
132     /* add k to row list for first nonzero element in k-th row */
133     if (nzk > 0){
134       i = q[k]; /* col value of first nonzero element in U(k, k+1:mbs-1) */
135       jl[k] = jl[i]; jl[i] = k;
136     }
137     iu[k+1] = iu[k] + nzk;
138 
139     /* allocate more space to ju if needed */
140     if (iu[k+1] > umax) {
141       /* estimate how much additional space we will need */
142       /* use the strategy suggested by David Hysom <hysom@perch-t.icase.edu> */
143       /* just double the memory each time */
144       maxadd = umax;
145       if (maxadd < nzk) maxadd = (mbs-k)*(nzk+1)/2;
146       umax += maxadd;
147 
148       /* allocate a longer ju */
149       ierr = PetscMalloc(umax*sizeof(PetscInt),&jutmp);CHKERRQ(ierr);
150       ierr = PetscMemcpy(jutmp,ju,iu[k]*sizeof(PetscInt));CHKERRQ(ierr);
151       ierr = PetscFree(ju);CHKERRQ(ierr);
152       ju   = jutmp;
153       reallocs++; /* count how many times we realloc */
154     }
155 
156     /* save nonzero structure of k-th row in ju */
157     i=k;
158     while (nzk --) {
159       i           = q[i];
160       ju[juidx++] = i;
161     }
162   }
163 
164 #if defined(PETSC_USE_DEBUG)
165   if (ai[mbs] != 0) {
166     PetscReal af = ((PetscReal)iu[mbs])/((PetscReal)ai[mbs]);
167     ierr = PetscLogInfo((A,"MatCholeskyFactorSymbolic_SeqSBAIJ:Reallocs %D Fill ratio:given %g needed %g\n",reallocs,f,af));CHKERRQ(ierr);
168     ierr = PetscLogInfo((A,"MatCholeskyFactorSymbolic_SeqSBAIJ:Run with -pc_cholesky_fill %g or use \n",af));CHKERRQ(ierr);
169     ierr = PetscLogInfo((A,"MatCholeskyFactorSymbolic_SeqSBAIJ:PCCholeskySetFill(pc,%g);\n",af));CHKERRQ(ierr);
170     ierr = PetscLogInfo((A,"MatCholeskyFactorSymbolic_SeqSBAIJ:for best performance.\n"));CHKERRQ(ierr);
171   } else {
172     ierr = PetscLogInfo((A,"MatCholeskyFactorSymbolic_SeqSBAIJ:Empty matrix.\n"));CHKERRQ(ierr);
173   }
174 #endif
175 
176   ierr = ISRestoreIndices(perm,&rip);CHKERRQ(ierr);
177   ierr = PetscFree(jl);CHKERRQ(ierr);
178 
179   /* put together the new matrix */
180   ierr = MatCreate(A->comm,bs*mbs,bs*mbs,bs*mbs,bs*mbs,B);CHKERRQ(ierr);
181   ierr = MatSetType(*B,A->type_name);CHKERRQ(ierr);
182   ierr = MatSeqSBAIJSetPreallocation(*B,bs,0,PETSC_NULL);CHKERRQ(ierr);
183 
184   /* ierr = PetscLogObjectParent(*B,iperm);CHKERRQ(ierr); */
185   b = (Mat_SeqSBAIJ*)(*B)->data;
186   ierr = PetscFree(b->imax);CHKERRQ(ierr);
187   b->singlemalloc = PETSC_FALSE;
188   /* the next line frees the default space generated by the Create() */
189   ierr = PetscFree(b->a);CHKERRQ(ierr);
190   ierr = PetscFree(b->ilen);CHKERRQ(ierr);
191   ierr = PetscMalloc((iu[mbs]+1)*sizeof(MatScalar)*bs2,&b->a);CHKERRQ(ierr);
192   b->j    = ju;
193   b->i    = iu;
194   b->diag = 0;
195   b->ilen = 0;
196   b->imax = 0;
197   b->row  = perm;
198   b->pivotinblocks = PETSC_FALSE; /* need to get from MatFactorInfo */
199   ierr    = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
200   b->icol = perm;
201   ierr    = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
202   ierr    = PetscMalloc((bs*mbs+bs)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
203   /* In b structure:  Free imax, ilen, old a, old j.
204      Allocate idnew, solve_work, new a, new j */
205   ierr = PetscLogObjectMemory(*B,(iu[mbs]-mbs)*(sizeof(PetscInt)+sizeof(MatScalar)));CHKERRQ(ierr);
206   b->maxnz = b->nz = iu[mbs];
207 
208   (*B)->factor                 = FACTOR_CHOLESKY;
209   (*B)->info.factor_mallocs    = reallocs;
210   (*B)->info.fill_ratio_given  = f;
211   if (ai[mbs] != 0) {
212     (*B)->info.fill_ratio_needed = ((PetscReal)iu[mbs])/((PetscReal)ai[mbs]);
213   } else {
214     (*B)->info.fill_ratio_needed = 0.0;
215   }
216 
217   if (perm_identity){
218     switch (bs) {
219       case 1:
220         (*B)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_1_NaturalOrdering;
221         (*B)->ops->solve           = MatSolve_SeqSBAIJ_1_NaturalOrdering;
222         ierr = PetscLogInfo((A,"MatICCFactorSymbolic_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=1\n"));CHKERRQ(ierr);
223         break;
224       case 2:
225         (*B)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_2_NaturalOrdering;
226         (*B)->ops->solve           = MatSolve_SeqSBAIJ_2_NaturalOrdering;
227         ierr = PetscLogInfo((A,"MatICCFactorSymbolic_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=2\n"));CHKERRQ(ierr);
228         break;
229       case 3:
230         (*B)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_3_NaturalOrdering;
231         (*B)->ops->solve           = MatSolve_SeqSBAIJ_3_NaturalOrdering;
232         ierr = PetscLogInfo((A,"MatICCFactorSymbolic_SeqSBAIJ:sing special in-place natural ordering factor and solve BS=3\n"));CHKERRQ(ierr);
233         break;
234       case 4:
235         (*B)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_4_NaturalOrdering;
236         (*B)->ops->solve           = MatSolve_SeqSBAIJ_4_NaturalOrdering;
237         ierr = PetscLogInfo((A,"MatICCFactorSymbolic_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=4\n"));CHKERRQ(ierr);
238         break;
239       case 5:
240         (*B)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_5_NaturalOrdering;
241         (*B)->ops->solve           = MatSolve_SeqSBAIJ_5_NaturalOrdering;
242         ierr = PetscLogInfo((A,"MatICCFactorSymbolic_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=5\n"));CHKERRQ(ierr);
243         break;
244       case 6:
245         (*B)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_6_NaturalOrdering;
246         (*B)->ops->solve           = MatSolve_SeqSBAIJ_6_NaturalOrdering;
247         ierr = PetscLogInfo((A,"MatICCFactorSymbolic_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=6\n"));CHKERRQ(ierr);
248         break;
249       case 7:
250         (*B)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_7_NaturalOrdering;
251         (*B)->ops->solve           = MatSolve_SeqSBAIJ_7_NaturalOrdering;
252         ierr = PetscLogInfo((A,"MatICCFactorSymbolic_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=7\n"));CHKERRQ(ierr);
253       break;
254       default:
255         (*B)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_N_NaturalOrdering;
256         (*B)->ops->solve           = MatSolve_SeqSBAIJ_N_NaturalOrdering;
257         ierr = PetscLogInfo((A,"MatICCFactorSymbolic_SeqSBAIJ:Using special in-place natural ordering factor and solve BS>7\n"));CHKERRQ(ierr);
258       break;
259     }
260   }
261   PetscFunctionReturn(0);
262 }
263 /*
264     Symbolic U^T*D*U factorization for SBAIJ format.
265 */
266 #include "petscbt.h"
267 #include "src/mat/utils/freespace.h"
268 #undef __FUNCT__
269 #define __FUNCT__ "MatCholeskyFactorSymbolic_SeqSBAIJ"
270 PetscErrorCode MatCholeskyFactorSymbolic_SeqSBAIJ(Mat A,IS perm,MatFactorInfo *info,Mat *fact)
271 {
272   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)A->data;
273   Mat_SeqSBAIJ   *b;
274   Mat            B;
275   PetscErrorCode ierr;
276   PetscTruth     perm_identity;
277   PetscReal      fill = info->fill;
278   PetscInt       *rip,i,mbs=a->mbs,bs=A->bs,*ai,*aj,reallocs=0,prow;
279   PetscInt       *jl,jmin,jmax,nzk,*ui,k,j,*il,nextprow;
280   PetscInt       nlnk,*lnk,ncols,*cols,*uj,**ui_ptr,*uj_ptr;
281   FreeSpaceList  free_space=PETSC_NULL,current_space=PETSC_NULL;
282   PetscBT        lnkbt;
283 
284   PetscFunctionBegin;
285   /*
286    This code originally uses Modified Sparse Row (MSR) storage
287    (see page 85, "Iterative Methods ..." by Saad) for the output matrix B - bad choise!
288    Then it is rewritten so the factor B takes seqsbaij format. However the associated
289    MatCholeskyFactorNumeric_() have not been modified for the cases of bs>1 or !perm_identity,
290    thus the original code in MSR format is still used for these cases.
291    The code below should replace MatCholeskyFactorSymbolic_SeqSBAIJ_MSR() whenever
292    MatCholeskyFactorNumeric_() is modified for using sbaij symbolic factor.
293   */
294   if (bs > 1){
295     ierr = MatCholeskyFactorSymbolic_SeqSBAIJ_MSR(A,perm,info,fact);CHKERRQ(ierr);
296     PetscFunctionReturn(0);
297   }
298 
299   /* check whether perm is the identity mapping */
300   ierr = ISIdentity(perm,&perm_identity);CHKERRQ(ierr);
301 
302   if (perm_identity){
303     a->permute = PETSC_FALSE;
304     ai = a->i; aj = a->j;
305   } else {
306     a->permute = PETSC_TRUE;
307     ierr = MatReorderingSeqSBAIJ(A,perm);CHKERRQ(ierr);
308     ai = a->inew; aj = a->jnew;
309   }
310   ierr = ISGetIndices(perm,&rip);CHKERRQ(ierr);
311 
312   /* initialization */
313   ierr  = PetscMalloc((mbs+1)*sizeof(PetscInt),&ui);CHKERRQ(ierr);
314   ui[0] = 0;
315 
316   /* jl: linked list for storing indices of the pivot rows
317      il: il[i] points to the 1st nonzero entry of U(i,k:mbs-1) */
318   ierr = PetscMalloc((3*mbs+1)*sizeof(PetscInt)+mbs*sizeof(PetscInt*),&jl);CHKERRQ(ierr);
319   il     = jl + mbs;
320   cols   = il + mbs;
321   ui_ptr = (PetscInt**)(cols + mbs);
322 
323   for (i=0; i<mbs; i++){
324     jl[i] = mbs; il[i] = 0;
325   }
326 
327   /* create and initialize a linked list for storing column indices of the active row k */
328   nlnk = mbs + 1;
329   ierr = PetscLLCreate(mbs,mbs,nlnk,lnk,lnkbt);CHKERRQ(ierr);
330 
331   /* initial FreeSpace size is fill*(ai[mbs]+1) */
332   ierr = GetMoreSpace((PetscInt)(fill*(ai[mbs]+1)),&free_space);CHKERRQ(ierr);
333   current_space = free_space;
334 
335   for (k=0; k<mbs; k++){  /* for each active row k */
336     /* initialize lnk by the column indices of row rip[k] of A */
337     nzk   = 0;
338     ncols = ai[rip[k]+1] - ai[rip[k]];
339     for (j=0; j<ncols; j++){
340       i = *(aj + ai[rip[k]] + j);
341       cols[j] = rip[i];
342     }
343     ierr = PetscLLAdd(ncols,cols,mbs,nlnk,lnk,lnkbt);CHKERRQ(ierr);
344     nzk += nlnk;
345 
346     /* update lnk by computing fill-in for each pivot row to be merged in */
347     prow = jl[k]; /* 1st pivot row */
348 
349     while (prow < k){
350       nextprow = jl[prow];
351       /* merge prow into k-th row */
352       jmin = il[prow] + 1;  /* index of the 2nd nzero entry in U(prow,k:mbs-1) */
353       jmax = ui[prow+1];
354       ncols = jmax-jmin;
355       uj_ptr = ui_ptr[prow] + jmin - ui[prow]; /* points to the 2nd nzero entry in U(prow,k:mbs-1) */
356       ierr = PetscLLAddSorted(ncols,uj_ptr,mbs,nlnk,lnk,lnkbt);CHKERRQ(ierr);
357       nzk += nlnk;
358 
359       /* update il and jl for prow */
360       if (jmin < jmax){
361         il[prow] = jmin;
362         j = *uj_ptr; jl[prow] = jl[j]; jl[j] = prow;
363       }
364       prow = nextprow;
365     }
366 
367     /* if free space is not available, make more free space */
368     if (current_space->local_remaining<nzk) {
369       i = mbs - k + 1; /* num of unfactored rows */
370       i = PetscMin(i*nzk, i*(i-1)); /* i*nzk, i*(i-1): estimated and max additional space needed */
371       ierr = GetMoreSpace(i,&current_space);CHKERRQ(ierr);
372       reallocs++;
373     }
374 
375     /* copy data into free space, then initialize lnk */
376     ierr = PetscLLClean(mbs,mbs,nzk,lnk,current_space->array,lnkbt);CHKERRQ(ierr);
377 
378     /* add the k-th row into il and jl */
379     if (nzk-1 > 0){
380       i = current_space->array[1]; /* col value of the first nonzero element in U(k, k+1:mbs-1) */
381       jl[k] = jl[i]; jl[i] = k;
382       il[k] = ui[k] + 1;
383     }
384     ui_ptr[k] = current_space->array;
385     current_space->array           += nzk;
386     current_space->local_used      += nzk;
387     current_space->local_remaining -= nzk;
388 
389     ui[k+1] = ui[k] + nzk;
390   }
391 
392 #if defined(PETSC_USE_DEBUG)
393   if (ai[mbs] != 0) {
394     PetscReal af = ((PetscReal)ui[mbs])/((PetscReal)ai[mbs]);
395     ierr = PetscLogInfo((A,"MatCholeskyFactorSymbolic_SeqSBAIJ:Reallocs %D Fill ratio:given %g needed %g\n",reallocs,fill,af));CHKERRQ(ierr);
396     ierr = PetscLogInfo((A,"MatCholeskyFactorSymbolic_SeqSBAIJ:Run with -pc_cholesky_fill %g or use \n",af));CHKERRQ(ierr);
397     ierr = PetscLogInfo((A,"MatCholeskyFactorSymbolic_SeqSBAIJ:PCCholeskySetFill(pc,%g) for best performance.\n",af));CHKERRQ(ierr);
398   } else {
399     ierr = PetscLogInfo((A,"MatCholeskyFactorSymbolic_SeqSBAIJ:Empty matrix.\n"));CHKERRQ(ierr);
400   }
401 #endif
402 
403   ierr = ISRestoreIndices(perm,&rip);CHKERRQ(ierr);
404   ierr = PetscFree(jl);CHKERRQ(ierr);
405 
406   /* destroy list of free space and other temporary array(s) */
407   ierr = PetscMalloc((ui[mbs]+1)*sizeof(PetscInt),&uj);CHKERRQ(ierr);
408   ierr = MakeSpaceContiguous(&free_space,uj);CHKERRQ(ierr);
409   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
410 
411   /* put together the new matrix in MATSEQSBAIJ format */
412   ierr = MatCreate(PETSC_COMM_SELF,mbs,mbs,mbs,mbs,fact);CHKERRQ(ierr);
413   B    = *fact;
414   ierr = MatSetType(B,MATSEQSBAIJ);CHKERRQ(ierr);
415   ierr = MatSeqSBAIJSetPreallocation(B,bs,0,PETSC_NULL);CHKERRQ(ierr);
416 
417   b = (Mat_SeqSBAIJ*)B->data;
418   ierr = PetscFree(b->imax);CHKERRQ(ierr);
419   b->singlemalloc = PETSC_FALSE;
420   /* the next line frees the default space generated by the Create() */
421   ierr = PetscFree(b->a);CHKERRQ(ierr);
422   ierr = PetscFree(b->ilen);CHKERRQ(ierr);
423   ierr = PetscMalloc((ui[mbs]+1)*sizeof(MatScalar),&b->a);CHKERRQ(ierr);
424   b->j    = uj;
425   b->i    = ui;
426   b->diag = 0;
427   b->ilen = 0;
428   b->imax = 0;
429   b->row  = perm;
430   b->pivotinblocks = PETSC_FALSE; /* need to get from MatFactorInfo */
431   ierr    = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
432   b->icol = perm;
433   ierr    = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr);
434   ierr    = PetscMalloc((mbs+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr);
435   ierr    = PetscLogObjectMemory(B,(ui[mbs]-mbs)*(sizeof(PetscInt)+sizeof(MatScalar)));CHKERRQ(ierr);
436   b->maxnz = b->nz = ui[mbs];
437 
438   B->factor                 = FACTOR_CHOLESKY;
439   B->info.factor_mallocs    = reallocs;
440   B->info.fill_ratio_given  = fill;
441   if (ai[mbs] != 0) {
442     B->info.fill_ratio_needed = ((PetscReal)ui[mbs])/((PetscReal)ai[mbs]);
443   } else {
444     B->info.fill_ratio_needed = 0.0;
445   }
446 
447   if (perm_identity){
448     switch (bs) {
449       case 1:
450         B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_1_NaturalOrdering;
451         B->ops->solve           = MatSolve_SeqSBAIJ_1_NaturalOrdering;
452         ierr = PetscLogInfo((A,"MatCholeskyFactorSymbolic_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=1\n"));CHKERRQ(ierr);
453         break;
454       case 2:
455         B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_2_NaturalOrdering;
456         B->ops->solve           = MatSolve_SeqSBAIJ_2_NaturalOrdering;
457         ierr = PetscLogInfo((A,"MatCholeskyFactorSymbolic_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=2\n"));CHKERRQ(ierr);
458         break;
459       case 3:
460         B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_3_NaturalOrdering;
461         B->ops->solve           = MatSolve_SeqSBAIJ_3_NaturalOrdering;
462         ierr = PetscLogInfo((A,"MatCholeskyFactorSymbolic_SeqSBAIJ:sing special in-place natural ordering factor and solve BS=3\n"));CHKERRQ(ierr);
463         break;
464       case 4:
465         B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_4_NaturalOrdering;
466         B->ops->solve           = MatSolve_SeqSBAIJ_4_NaturalOrdering;
467         ierr = PetscLogInfo((A,"MatCholeskyFactorSymbolic_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=4\n"));CHKERRQ(ierr);
468         break;
469       case 5:
470         B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_5_NaturalOrdering;
471         B->ops->solve           = MatSolve_SeqSBAIJ_5_NaturalOrdering;
472         ierr = PetscLogInfo((A,"MatCholeskyFactorSymbolic_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=5\n"));CHKERRQ(ierr);
473         break;
474       case 6:
475         B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_6_NaturalOrdering;
476         B->ops->solve           = MatSolve_SeqSBAIJ_6_NaturalOrdering;
477         ierr = PetscLogInfo((A,"MatCholeskyFactorSymbolic_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=6\n"));CHKERRQ(ierr);
478         break;
479       case 7:
480         B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_7_NaturalOrdering;
481         B->ops->solve           = MatSolve_SeqSBAIJ_7_NaturalOrdering;
482         ierr = PetscLogInfo((A,"MatCholeskyFactorSymbolic_SeqSBAIJ:Using special in-place natural ordering factor and solve BS=7\n"));CHKERRQ(ierr);
483       break;
484       default:
485         B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_N_NaturalOrdering;
486         B->ops->solve           = MatSolve_SeqSBAIJ_N_NaturalOrdering;
487         ierr = PetscLogInfo((A,"MatCholeskyFactorSymbolic_SeqSBAIJ:Using special in-place natural ordering factor and solve BS>7\n"));CHKERRQ(ierr);
488       break;
489     }
490   }
491   PetscFunctionReturn(0);
492 }
493 #undef __FUNCT__
494 #define __FUNCT__ "MatCholeskyFactorNumeric_SeqSBAIJ_N"
495 PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_N(Mat A,MatFactorInfo *info,Mat *B)
496 {
497   Mat            C = *B;
498   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)A->data,*b = (Mat_SeqSBAIJ *)C->data;
499   IS             perm = b->row;
500   PetscErrorCode ierr;
501   PetscInt       *perm_ptr,i,j,mbs=a->mbs,*bi=b->i,*bj=b->j;
502   PetscInt       *ai,*aj,*a2anew,k,k1,jmin,jmax,*jl,*il,vj,nexti,ili;
503   PetscInt       bs=A->bs,bs2 = a->bs2;
504   MatScalar      *ba = b->a,*aa,*ap,*dk,*uik;
505   MatScalar      *u,*diag,*rtmp,*rtmp_ptr;
506   MatScalar      *work;
507   PetscInt       *pivots;
508 
509   PetscFunctionBegin;
510   /* initialization */
511   ierr = PetscMalloc(bs2*mbs*sizeof(MatScalar),&rtmp);CHKERRQ(ierr);
512   ierr = PetscMemzero(rtmp,bs2*mbs*sizeof(MatScalar));CHKERRQ(ierr);
513   ierr = PetscMalloc(2*mbs*sizeof(PetscInt),&il);CHKERRQ(ierr);
514   jl   = il + mbs;
515   for (i=0; i<mbs; i++) {
516     jl[i] = mbs; il[0] = 0;
517   }
518   ierr = PetscMalloc((2*bs2+bs)*sizeof(MatScalar),&dk);CHKERRQ(ierr);
519   uik  = dk + bs2;
520   work = uik + bs2;
521   ierr = PetscMalloc(bs*sizeof(PetscInt),&pivots);CHKERRQ(ierr);
522 
523   ierr  = ISGetIndices(perm,&perm_ptr);CHKERRQ(ierr);
524 
525   /* check permutation */
526   if (!a->permute){
527     ai = a->i; aj = a->j; aa = a->a;
528   } else {
529     ai   = a->inew; aj = a->jnew;
530     ierr = PetscMalloc(bs2*ai[mbs]*sizeof(MatScalar),&aa);CHKERRQ(ierr);
531     ierr = PetscMemcpy(aa,a->a,bs2*ai[mbs]*sizeof(MatScalar));CHKERRQ(ierr);
532     ierr = PetscMalloc(ai[mbs]*sizeof(PetscInt),&a2anew);CHKERRQ(ierr);
533     ierr = PetscMemcpy(a2anew,a->a2anew,(ai[mbs])*sizeof(PetscInt));CHKERRQ(ierr);
534 
535     for (i=0; i<mbs; i++){
536       jmin = ai[i]; jmax = ai[i+1];
537       for (j=jmin; j<jmax; j++){
538         while (a2anew[j] != j){
539           k = a2anew[j]; a2anew[j] = a2anew[k]; a2anew[k] = k;
540           for (k1=0; k1<bs2; k1++){
541             dk[k1]       = aa[k*bs2+k1];
542             aa[k*bs2+k1] = aa[j*bs2+k1];
543             aa[j*bs2+k1] = dk[k1];
544           }
545         }
546         /* transform columnoriented blocks that lie in the lower triangle to roworiented blocks */
547         if (i > aj[j]){
548           /* printf("change orientation, row: %d, col: %d\n",i,aj[j]); */
549           ap = aa + j*bs2;                     /* ptr to the beginning of j-th block of aa */
550           for (k=0; k<bs2; k++) dk[k] = ap[k]; /* dk <- j-th block of aa */
551           for (k=0; k<bs; k++){               /* j-th block of aa <- dk^T */
552             for (k1=0; k1<bs; k1++) *ap++ = dk[k + bs*k1];
553           }
554         }
555       }
556     }
557     ierr = PetscFree(a2anew);CHKERRQ(ierr);
558   }
559 
560   /* for each row k */
561   for (k = 0; k<mbs; k++){
562 
563     /*initialize k-th row with elements nonzero in row perm(k) of A */
564     jmin = ai[perm_ptr[k]]; jmax = ai[perm_ptr[k]+1];
565 
566     ap = aa + jmin*bs2;
567     for (j = jmin; j < jmax; j++){
568       vj = perm_ptr[aj[j]];         /* block col. index */
569       rtmp_ptr = rtmp + vj*bs2;
570       for (i=0; i<bs2; i++) *rtmp_ptr++ = *ap++;
571     }
572 
573     /* modify k-th row by adding in those rows i with U(i,k) != 0 */
574     ierr = PetscMemcpy(dk,rtmp+k*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr);
575     i = jl[k]; /* first row to be added to k_th row  */
576 
577     while (i < k){
578       nexti = jl[i]; /* next row to be added to k_th row */
579 
580       /* compute multiplier */
581       ili = il[i];  /* index of first nonzero element in U(i,k:bms-1) */
582 
583       /* uik = -inv(Di)*U_bar(i,k) */
584       diag = ba + i*bs2;
585       u    = ba + ili*bs2;
586       ierr = PetscMemzero(uik,bs2*sizeof(MatScalar));CHKERRQ(ierr);
587       Kernel_A_gets_A_minus_B_times_C(bs,uik,diag,u);
588 
589       /* update D(k) += -U(i,k)^T * U_bar(i,k) */
590       Kernel_A_gets_A_plus_Btranspose_times_C(bs,dk,uik,u);
591 
592       /* update -U(i,k) */
593       ierr = PetscMemcpy(ba+ili*bs2,uik,bs2*sizeof(MatScalar));CHKERRQ(ierr);
594 
595       /* add multiple of row i to k-th row ... */
596       jmin = ili + 1; jmax = bi[i+1];
597       if (jmin < jmax){
598         for (j=jmin; j<jmax; j++) {
599           /* rtmp += -U(i,k)^T * U_bar(i,j) */
600           rtmp_ptr = rtmp + bj[j]*bs2;
601           u = ba + j*bs2;
602           Kernel_A_gets_A_plus_Btranspose_times_C(bs,rtmp_ptr,uik,u);
603         }
604 
605         /* ... add i to row list for next nonzero entry */
606         il[i] = jmin;             /* update il(i) in column k+1, ... mbs-1 */
607         j     = bj[jmin];
608         jl[i] = jl[j]; jl[j] = i; /* update jl */
609       }
610       i = nexti;
611     }
612 
613     /* save nonzero entries in k-th row of U ... */
614 
615     /* invert diagonal block */
616     diag = ba+k*bs2;
617     ierr = PetscMemcpy(diag,dk,bs2*sizeof(MatScalar));CHKERRQ(ierr);
618     ierr = Kernel_A_gets_inverse_A(bs,diag,pivots,work);CHKERRQ(ierr);
619 
620     jmin = bi[k]; jmax = bi[k+1];
621     if (jmin < jmax) {
622       for (j=jmin; j<jmax; j++){
623          vj = bj[j];           /* block col. index of U */
624          u   = ba + j*bs2;
625          rtmp_ptr = rtmp + vj*bs2;
626          for (k1=0; k1<bs2; k1++){
627            *u++        = *rtmp_ptr;
628            *rtmp_ptr++ = 0.0;
629          }
630       }
631 
632       /* ... add k to row list for first nonzero entry in k-th row */
633       il[k] = jmin;
634       i     = bj[jmin];
635       jl[k] = jl[i]; jl[i] = k;
636     }
637   }
638 
639   ierr = PetscFree(rtmp);CHKERRQ(ierr);
640   ierr = PetscFree(il);CHKERRQ(ierr);
641   ierr = PetscFree(dk);CHKERRQ(ierr);
642   ierr = PetscFree(pivots);CHKERRQ(ierr);
643   if (a->permute){
644     ierr = PetscFree(aa);CHKERRQ(ierr);
645   }
646 
647   ierr = ISRestoreIndices(perm,&perm_ptr);CHKERRQ(ierr);
648   C->factor       = FACTOR_CHOLESKY;
649   C->assembled    = PETSC_TRUE;
650   C->preallocated = PETSC_TRUE;
651   ierr = PetscLogFlops(1.3333*bs*bs2*b->mbs);CHKERRQ(ierr); /* from inverting diagonal blocks */
652   PetscFunctionReturn(0);
653 }
654 
655 #undef __FUNCT__
656 #define __FUNCT__ "MatCholeskyFactorNumeric_SeqSBAIJ_N_NaturalOrdering"
657 PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_N_NaturalOrdering(Mat A,MatFactorInfo *info,Mat *B)
658 {
659   Mat            C = *B;
660   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)A->data,*b = (Mat_SeqSBAIJ *)C->data;
661   PetscErrorCode ierr;
662   PetscInt       i,j,mbs=a->mbs,*bi=b->i,*bj=b->j;
663   PetscInt       *ai,*aj,k,k1,jmin,jmax,*jl,*il,vj,nexti,ili;
664   PetscInt       bs=A->bs,bs2 = a->bs2;
665   MatScalar      *ba = b->a,*aa,*ap,*dk,*uik;
666   MatScalar      *u,*diag,*rtmp,*rtmp_ptr;
667   MatScalar      *work;
668   PetscInt       *pivots;
669 
670   PetscFunctionBegin;
671   ierr = PetscMalloc(bs2*mbs*sizeof(MatScalar),&rtmp);CHKERRQ(ierr);
672   ierr = PetscMemzero(rtmp,bs2*mbs*sizeof(MatScalar));CHKERRQ(ierr);
673   ierr = PetscMalloc(2*mbs*sizeof(PetscInt),&il);CHKERRQ(ierr);
674   jl   = il + mbs;
675   for (i=0; i<mbs; i++) {
676     jl[i] = mbs; il[0] = 0;
677   }
678   ierr = PetscMalloc((2*bs2+bs)*sizeof(MatScalar),&dk);CHKERRQ(ierr);
679   uik  = dk + bs2;
680   work = uik + bs2;
681   ierr = PetscMalloc(bs*sizeof(PetscInt),&pivots);CHKERRQ(ierr);
682 
683   ai = a->i; aj = a->j; aa = a->a;
684 
685   /* for each row k */
686   for (k = 0; k<mbs; k++){
687 
688     /*initialize k-th row with elements nonzero in row k of A */
689     jmin = ai[k]; jmax = ai[k+1];
690     ap = aa + jmin*bs2;
691     for (j = jmin; j < jmax; j++){
692       vj = aj[j];         /* block col. index */
693       rtmp_ptr = rtmp + vj*bs2;
694       for (i=0; i<bs2; i++) *rtmp_ptr++ = *ap++;
695     }
696 
697     /* modify k-th row by adding in those rows i with U(i,k) != 0 */
698     ierr = PetscMemcpy(dk,rtmp+k*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr);
699     i = jl[k]; /* first row to be added to k_th row  */
700 
701     while (i < k){
702       nexti = jl[i]; /* next row to be added to k_th row */
703 
704       /* compute multiplier */
705       ili = il[i];  /* index of first nonzero element in U(i,k:bms-1) */
706 
707       /* uik = -inv(Di)*U_bar(i,k) */
708       diag = ba + i*bs2;
709       u    = ba + ili*bs2;
710       ierr = PetscMemzero(uik,bs2*sizeof(MatScalar));CHKERRQ(ierr);
711       Kernel_A_gets_A_minus_B_times_C(bs,uik,diag,u);
712 
713       /* update D(k) += -U(i,k)^T * U_bar(i,k) */
714       Kernel_A_gets_A_plus_Btranspose_times_C(bs,dk,uik,u);
715 
716       /* update -U(i,k) */
717       ierr = PetscMemcpy(ba+ili*bs2,uik,bs2*sizeof(MatScalar));CHKERRQ(ierr);
718 
719       /* add multiple of row i to k-th row ... */
720       jmin = ili + 1; jmax = bi[i+1];
721       if (jmin < jmax){
722         for (j=jmin; j<jmax; j++) {
723           /* rtmp += -U(i,k)^T * U_bar(i,j) */
724           rtmp_ptr = rtmp + bj[j]*bs2;
725           u = ba + j*bs2;
726           Kernel_A_gets_A_plus_Btranspose_times_C(bs,rtmp_ptr,uik,u);
727         }
728 
729         /* ... add i to row list for next nonzero entry */
730         il[i] = jmin;             /* update il(i) in column k+1, ... mbs-1 */
731         j     = bj[jmin];
732         jl[i] = jl[j]; jl[j] = i; /* update jl */
733       }
734       i = nexti;
735     }
736 
737     /* save nonzero entries in k-th row of U ... */
738 
739     /* invert diagonal block */
740     diag = ba+k*bs2;
741     ierr = PetscMemcpy(diag,dk,bs2*sizeof(MatScalar));CHKERRQ(ierr);
742     ierr = Kernel_A_gets_inverse_A(bs,diag,pivots,work);CHKERRQ(ierr);
743 
744     jmin = bi[k]; jmax = bi[k+1];
745     if (jmin < jmax) {
746       for (j=jmin; j<jmax; j++){
747          vj = bj[j];           /* block col. index of U */
748          u   = ba + j*bs2;
749          rtmp_ptr = rtmp + vj*bs2;
750          for (k1=0; k1<bs2; k1++){
751            *u++        = *rtmp_ptr;
752            *rtmp_ptr++ = 0.0;
753          }
754       }
755 
756       /* ... add k to row list for first nonzero entry in k-th row */
757       il[k] = jmin;
758       i     = bj[jmin];
759       jl[k] = jl[i]; jl[i] = k;
760     }
761   }
762 
763   ierr = PetscFree(rtmp);CHKERRQ(ierr);
764   ierr = PetscFree(il);CHKERRQ(ierr);
765   ierr = PetscFree(dk);CHKERRQ(ierr);
766   ierr = PetscFree(pivots);CHKERRQ(ierr);
767 
768   C->factor    = FACTOR_CHOLESKY;
769   C->assembled = PETSC_TRUE;
770   C->preallocated = PETSC_TRUE;
771   ierr = PetscLogFlops(1.3333*bs*bs2*b->mbs);CHKERRQ(ierr); /* from inverting diagonal blocks */
772   PetscFunctionReturn(0);
773 }
774 
775 /*
776     Numeric U^T*D*U factorization for SBAIJ format. Modified from SNF of YSMP.
777     Version for blocks 2 by 2.
778 */
779 #undef __FUNCT__
780 #define __FUNCT__ "MatCholeskyFactorNumeric_SeqSBAIJ_2"
781 PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_2(Mat A,MatFactorInfo *info,Mat *B)
782 {
783   Mat            C = *B;
784   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)A->data,*b = (Mat_SeqSBAIJ *)C->data;
785   IS             perm = b->row;
786   PetscErrorCode ierr;
787   PetscInt       *perm_ptr,i,j,mbs=a->mbs,*bi=b->i,*bj=b->j;
788   PetscInt       *ai,*aj,*a2anew,k,k1,jmin,jmax,*jl,*il,vj,nexti,ili;
789   MatScalar      *ba = b->a,*aa,*ap,*dk,*uik;
790   MatScalar      *u,*diag,*rtmp,*rtmp_ptr;
791 
792   PetscFunctionBegin;
793   /* initialization */
794   /* il and jl record the first nonzero element in each row of the accessing
795      window U(0:k, k:mbs-1).
796      jl:    list of rows to be added to uneliminated rows
797             i>= k: jl(i) is the first row to be added to row i
798             i<  k: jl(i) is the row following row i in some list of rows
799             jl(i) = mbs indicates the end of a list
800      il(i): points to the first nonzero element in columns k,...,mbs-1 of
801             row i of U */
802   ierr = PetscMalloc(4*mbs*sizeof(MatScalar),&rtmp);CHKERRQ(ierr);
803   ierr = PetscMemzero(rtmp,4*mbs*sizeof(MatScalar));CHKERRQ(ierr);
804   ierr = PetscMalloc(2*mbs*sizeof(PetscInt),&il);CHKERRQ(ierr);
805   jl   = il + mbs;
806   for (i=0; i<mbs; i++) {
807     jl[i] = mbs; il[0] = 0;
808   }
809   ierr = PetscMalloc(8*sizeof(MatScalar),&dk);CHKERRQ(ierr);
810   uik  = dk + 4;
811   ierr = ISGetIndices(perm,&perm_ptr);CHKERRQ(ierr);
812 
813   /* check permutation */
814   if (!a->permute){
815     ai = a->i; aj = a->j; aa = a->a;
816   } else {
817     ai   = a->inew; aj = a->jnew;
818     ierr = PetscMalloc(4*ai[mbs]*sizeof(MatScalar),&aa);CHKERRQ(ierr);
819     ierr = PetscMemcpy(aa,a->a,4*ai[mbs]*sizeof(MatScalar));CHKERRQ(ierr);
820     ierr = PetscMalloc(ai[mbs]*sizeof(PetscInt),&a2anew);CHKERRQ(ierr);
821     ierr = PetscMemcpy(a2anew,a->a2anew,(ai[mbs])*sizeof(PetscInt));CHKERRQ(ierr);
822 
823     for (i=0; i<mbs; i++){
824       jmin = ai[i]; jmax = ai[i+1];
825       for (j=jmin; j<jmax; j++){
826         while (a2anew[j] != j){
827           k = a2anew[j]; a2anew[j] = a2anew[k]; a2anew[k] = k;
828           for (k1=0; k1<4; k1++){
829             dk[k1]       = aa[k*4+k1];
830             aa[k*4+k1] = aa[j*4+k1];
831             aa[j*4+k1] = dk[k1];
832           }
833         }
834         /* transform columnoriented blocks that lie in the lower triangle to roworiented blocks */
835         if (i > aj[j]){
836           /* printf("change orientation, row: %d, col: %d\n",i,aj[j]); */
837           ap = aa + j*4;     /* ptr to the beginning of the block */
838           dk[1] = ap[1];     /* swap ap[1] and ap[2] */
839           ap[1] = ap[2];
840           ap[2] = dk[1];
841         }
842       }
843     }
844     ierr = PetscFree(a2anew);CHKERRQ(ierr);
845   }
846 
847   /* for each row k */
848   for (k = 0; k<mbs; k++){
849 
850     /*initialize k-th row with elements nonzero in row perm(k) of A */
851     jmin = ai[perm_ptr[k]]; jmax = ai[perm_ptr[k]+1];
852     ap = aa + jmin*4;
853     for (j = jmin; j < jmax; j++){
854       vj = perm_ptr[aj[j]];         /* block col. index */
855       rtmp_ptr = rtmp + vj*4;
856       for (i=0; i<4; i++) *rtmp_ptr++ = *ap++;
857     }
858 
859     /* modify k-th row by adding in those rows i with U(i,k) != 0 */
860     ierr = PetscMemcpy(dk,rtmp+k*4,4*sizeof(MatScalar));CHKERRQ(ierr);
861     i = jl[k]; /* first row to be added to k_th row  */
862 
863     while (i < k){
864       nexti = jl[i]; /* next row to be added to k_th row */
865 
866       /* compute multiplier */
867       ili = il[i];  /* index of first nonzero element in U(i,k:bms-1) */
868 
869       /* uik = -inv(Di)*U_bar(i,k): - ba[ili]*ba[i] */
870       diag = ba + i*4;
871       u    = ba + ili*4;
872       uik[0] = -(diag[0]*u[0] + diag[2]*u[1]);
873       uik[1] = -(diag[1]*u[0] + diag[3]*u[1]);
874       uik[2] = -(diag[0]*u[2] + diag[2]*u[3]);
875       uik[3] = -(diag[1]*u[2] + diag[3]*u[3]);
876 
877       /* update D(k) += -U(i,k)^T * U_bar(i,k): dk += uik*ba[ili] */
878       dk[0] += uik[0]*u[0] + uik[1]*u[1];
879       dk[1] += uik[2]*u[0] + uik[3]*u[1];
880       dk[2] += uik[0]*u[2] + uik[1]*u[3];
881       dk[3] += uik[2]*u[2] + uik[3]*u[3];
882 
883       /* update -U(i,k): ba[ili] = uik */
884       ierr = PetscMemcpy(ba+ili*4,uik,4*sizeof(MatScalar));CHKERRQ(ierr);
885 
886       /* add multiple of row i to k-th row ... */
887       jmin = ili + 1; jmax = bi[i+1];
888       if (jmin < jmax){
889         for (j=jmin; j<jmax; j++) {
890           /* rtmp += -U(i,k)^T * U_bar(i,j): rtmp[bj[j]] += uik*ba[j]; */
891           rtmp_ptr = rtmp + bj[j]*4;
892           u = ba + j*4;
893           rtmp_ptr[0] += uik[0]*u[0] + uik[1]*u[1];
894           rtmp_ptr[1] += uik[2]*u[0] + uik[3]*u[1];
895           rtmp_ptr[2] += uik[0]*u[2] + uik[1]*u[3];
896           rtmp_ptr[3] += uik[2]*u[2] + uik[3]*u[3];
897         }
898 
899         /* ... add i to row list for next nonzero entry */
900         il[i] = jmin;             /* update il(i) in column k+1, ... mbs-1 */
901         j     = bj[jmin];
902         jl[i] = jl[j]; jl[j] = i; /* update jl */
903       }
904       i = nexti;
905     }
906 
907     /* save nonzero entries in k-th row of U ... */
908 
909     /* invert diagonal block */
910     diag = ba+k*4;
911     ierr = PetscMemcpy(diag,dk,4*sizeof(MatScalar));CHKERRQ(ierr);
912     ierr = Kernel_A_gets_inverse_A_2(diag);CHKERRQ(ierr);
913 
914     jmin = bi[k]; jmax = bi[k+1];
915     if (jmin < jmax) {
916       for (j=jmin; j<jmax; j++){
917          vj = bj[j];           /* block col. index of U */
918          u   = ba + j*4;
919          rtmp_ptr = rtmp + vj*4;
920          for (k1=0; k1<4; k1++){
921            *u++        = *rtmp_ptr;
922            *rtmp_ptr++ = 0.0;
923          }
924       }
925 
926       /* ... add k to row list for first nonzero entry in k-th row */
927       il[k] = jmin;
928       i     = bj[jmin];
929       jl[k] = jl[i]; jl[i] = k;
930     }
931   }
932 
933   ierr = PetscFree(rtmp);CHKERRQ(ierr);
934   ierr = PetscFree(il);CHKERRQ(ierr);
935   ierr = PetscFree(dk);CHKERRQ(ierr);
936   if (a->permute) {
937     ierr = PetscFree(aa);CHKERRQ(ierr);
938   }
939   ierr = ISRestoreIndices(perm,&perm_ptr);CHKERRQ(ierr);
940   C->factor    = FACTOR_CHOLESKY;
941   C->assembled = PETSC_TRUE;
942   C->preallocated = PETSC_TRUE;
943   ierr = PetscLogFlops(1.3333*8*b->mbs);CHKERRQ(ierr); /* from inverting diagonal blocks */
944   PetscFunctionReturn(0);
945 }
946 
947 /*
948       Version for when blocks are 2 by 2 Using natural ordering
949 */
950 #undef __FUNCT__
951 #define __FUNCT__ "MatCholeskyFactorNumeric_SeqSBAIJ_2_NaturalOrdering"
952 PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_2_NaturalOrdering(Mat A,MatFactorInfo *info,Mat *B)
953 {
954   Mat            C = *B;
955   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)A->data,*b = (Mat_SeqSBAIJ *)C->data;
956   PetscErrorCode ierr;
957   PetscInt       i,j,mbs=a->mbs,*bi=b->i,*bj=b->j;
958   PetscInt       *ai,*aj,k,k1,jmin,jmax,*jl,*il,vj,nexti,ili;
959   MatScalar      *ba = b->a,*aa,*ap,*dk,*uik;
960   MatScalar      *u,*diag,*rtmp,*rtmp_ptr;
961 
962   PetscFunctionBegin;
963   /* initialization */
964   /* il and jl record the first nonzero element in each row of the accessing
965      window U(0:k, k:mbs-1).
966      jl:    list of rows to be added to uneliminated rows
967             i>= k: jl(i) is the first row to be added to row i
968             i<  k: jl(i) is the row following row i in some list of rows
969             jl(i) = mbs indicates the end of a list
970      il(i): points to the first nonzero element in columns k,...,mbs-1 of
971             row i of U */
972   ierr = PetscMalloc(4*mbs*sizeof(MatScalar),&rtmp);CHKERRQ(ierr);
973   ierr = PetscMemzero(rtmp,4*mbs*sizeof(MatScalar));CHKERRQ(ierr);
974   ierr = PetscMalloc(2*mbs*sizeof(PetscInt),&il);CHKERRQ(ierr);
975   jl   = il + mbs;
976   for (i=0; i<mbs; i++) {
977     jl[i] = mbs; il[0] = 0;
978   }
979   ierr = PetscMalloc(8*sizeof(MatScalar),&dk);CHKERRQ(ierr);
980   uik  = dk + 4;
981 
982   ai = a->i; aj = a->j; aa = a->a;
983 
984   /* for each row k */
985   for (k = 0; k<mbs; k++){
986 
987     /*initialize k-th row with elements nonzero in row k of A */
988     jmin = ai[k]; jmax = ai[k+1];
989     ap = aa + jmin*4;
990     for (j = jmin; j < jmax; j++){
991       vj = aj[j];         /* block col. index */
992       rtmp_ptr = rtmp + vj*4;
993       for (i=0; i<4; i++) *rtmp_ptr++ = *ap++;
994     }
995 
996     /* modify k-th row by adding in those rows i with U(i,k) != 0 */
997     ierr = PetscMemcpy(dk,rtmp+k*4,4*sizeof(MatScalar));CHKERRQ(ierr);
998     i = jl[k]; /* first row to be added to k_th row  */
999 
1000     while (i < k){
1001       nexti = jl[i]; /* next row to be added to k_th row */
1002 
1003       /* compute multiplier */
1004       ili = il[i];  /* index of first nonzero element in U(i,k:bms-1) */
1005 
1006       /* uik = -inv(Di)*U_bar(i,k): - ba[ili]*ba[i] */
1007       diag = ba + i*4;
1008       u    = ba + ili*4;
1009       uik[0] = -(diag[0]*u[0] + diag[2]*u[1]);
1010       uik[1] = -(diag[1]*u[0] + diag[3]*u[1]);
1011       uik[2] = -(diag[0]*u[2] + diag[2]*u[3]);
1012       uik[3] = -(diag[1]*u[2] + diag[3]*u[3]);
1013 
1014       /* update D(k) += -U(i,k)^T * U_bar(i,k): dk += uik*ba[ili] */
1015       dk[0] += uik[0]*u[0] + uik[1]*u[1];
1016       dk[1] += uik[2]*u[0] + uik[3]*u[1];
1017       dk[2] += uik[0]*u[2] + uik[1]*u[3];
1018       dk[3] += uik[2]*u[2] + uik[3]*u[3];
1019 
1020       /* update -U(i,k): ba[ili] = uik */
1021       ierr = PetscMemcpy(ba+ili*4,uik,4*sizeof(MatScalar));CHKERRQ(ierr);
1022 
1023       /* add multiple of row i to k-th row ... */
1024       jmin = ili + 1; jmax = bi[i+1];
1025       if (jmin < jmax){
1026         for (j=jmin; j<jmax; j++) {
1027           /* rtmp += -U(i,k)^T * U_bar(i,j): rtmp[bj[j]] += uik*ba[j]; */
1028           rtmp_ptr = rtmp + bj[j]*4;
1029           u = ba + j*4;
1030           rtmp_ptr[0] += uik[0]*u[0] + uik[1]*u[1];
1031           rtmp_ptr[1] += uik[2]*u[0] + uik[3]*u[1];
1032           rtmp_ptr[2] += uik[0]*u[2] + uik[1]*u[3];
1033           rtmp_ptr[3] += uik[2]*u[2] + uik[3]*u[3];
1034         }
1035 
1036         /* ... add i to row list for next nonzero entry */
1037         il[i] = jmin;             /* update il(i) in column k+1, ... mbs-1 */
1038         j     = bj[jmin];
1039         jl[i] = jl[j]; jl[j] = i; /* update jl */
1040       }
1041       i = nexti;
1042     }
1043 
1044     /* save nonzero entries in k-th row of U ... */
1045 
1046     /* invert diagonal block */
1047     diag = ba+k*4;
1048     ierr = PetscMemcpy(diag,dk,4*sizeof(MatScalar));CHKERRQ(ierr);
1049     ierr = Kernel_A_gets_inverse_A_2(diag);CHKERRQ(ierr);
1050 
1051     jmin = bi[k]; jmax = bi[k+1];
1052     if (jmin < jmax) {
1053       for (j=jmin; j<jmax; j++){
1054          vj = bj[j];           /* block col. index of U */
1055          u   = ba + j*4;
1056          rtmp_ptr = rtmp + vj*4;
1057          for (k1=0; k1<4; k1++){
1058            *u++        = *rtmp_ptr;
1059            *rtmp_ptr++ = 0.0;
1060          }
1061       }
1062 
1063       /* ... add k to row list for first nonzero entry in k-th row */
1064       il[k] = jmin;
1065       i     = bj[jmin];
1066       jl[k] = jl[i]; jl[i] = k;
1067     }
1068   }
1069 
1070   ierr = PetscFree(rtmp);CHKERRQ(ierr);
1071   ierr = PetscFree(il);CHKERRQ(ierr);
1072   ierr = PetscFree(dk);CHKERRQ(ierr);
1073 
1074   C->factor    = FACTOR_CHOLESKY;
1075   C->assembled = PETSC_TRUE;
1076   C->preallocated = PETSC_TRUE;
1077   ierr = PetscLogFlops(1.3333*8*b->mbs);CHKERRQ(ierr); /* from inverting diagonal blocks */
1078   PetscFunctionReturn(0);
1079 }
1080 
1081 /*
1082     Numeric U^T*D*U factorization for SBAIJ format.
1083     Version for blocks are 1 by 1.
1084 */
1085 #undef __FUNCT__
1086 #define __FUNCT__ "MatCholeskyFactorNumeric_SeqSBAIJ_1"
1087 PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_1(Mat A,MatFactorInfo *info,Mat *B)
1088 {
1089   Mat            C = *B;
1090   Mat_SeqSBAIJ   *a=(Mat_SeqSBAIJ*)A->data,*b=(Mat_SeqSBAIJ *)C->data;
1091   IS             ip=b->row;
1092   PetscErrorCode ierr;
1093   PetscInt       *rip,i,j,mbs=a->mbs,*bi=b->i,*bj=b->j,*bcol;
1094   PetscInt       *ai,*aj,*a2anew;
1095   PetscInt       k,jmin,jmax,*jl,*il,col,nexti,ili,nz;
1096   MatScalar      *rtmp,*ba=b->a,*bval,*aa,dk,uikdi;
1097   PetscReal      zeropivot,rs,shiftnz;
1098   PetscTruth     shiftpd;
1099   ChShift_Ctx    sctx;
1100   PetscInt       newshift;
1101 
1102   PetscFunctionBegin;
1103   /* initialization */
1104   shiftnz   = info->shiftnz;
1105   shiftpd   = info->shiftpd;
1106   zeropivot = info->zeropivot;
1107 
1108   ierr  = ISGetIndices(ip,&rip);CHKERRQ(ierr);
1109   if (!a->permute){
1110     ai = a->i; aj = a->j; aa = a->a;
1111   } else {
1112     ai = a->inew; aj = a->jnew;
1113     nz = ai[mbs];
1114     ierr = PetscMalloc(nz*sizeof(MatScalar),&aa);CHKERRQ(ierr);
1115     a2anew = a->a2anew;
1116     bval   = a->a;
1117     for (j=0; j<nz; j++){
1118       aa[a2anew[j]] = *(bval++);
1119     }
1120   }
1121 
1122   /* initialization */
1123   /* il and jl record the first nonzero element in each row of the accessing
1124      window U(0:k, k:mbs-1).
1125      jl:    list of rows to be added to uneliminated rows
1126             i>= k: jl(i) is the first row to be added to row i
1127             i<  k: jl(i) is the row following row i in some list of rows
1128             jl(i) = mbs indicates the end of a list
1129      il(i): points to the first nonzero element in columns k,...,mbs-1 of
1130             row i of U */
1131   nz   = (2*mbs+1)*sizeof(PetscInt)+mbs*sizeof(MatScalar);
1132   ierr = PetscMalloc(nz,&il);CHKERRQ(ierr);
1133   jl   = il + mbs;
1134   rtmp = (MatScalar*)(jl + mbs);
1135 
1136   sctx.shift_amount = 0;
1137   sctx.nshift       = 0;
1138   do {
1139     sctx.chshift = PETSC_FALSE;
1140     for (i=0; i<mbs; i++) {
1141       rtmp[i] = 0.0; jl[i] = mbs; il[0] = 0;
1142     }
1143 
1144     for (k = 0; k<mbs; k++){
1145       /*initialize k-th row by the perm[k]-th row of A */
1146       jmin = ai[rip[k]]; jmax = ai[rip[k]+1];
1147       bval = ba + bi[k];
1148       for (j = jmin; j < jmax; j++){
1149         col = rip[aj[j]];
1150         rtmp[col] = aa[j];
1151         *bval++  = 0.0; /* for in-place factorization */
1152       }
1153 
1154       /* shift the diagonal of the matrix */
1155       if (sctx.nshift) rtmp[k] += sctx.shift_amount;
1156 
1157       /* modify k-th row by adding in those rows i with U(i,k)!=0 */
1158       dk = rtmp[k];
1159       i = jl[k]; /* first row to be added to k_th row  */
1160 
1161       while (i < k){
1162         nexti = jl[i]; /* next row to be added to k_th row */
1163 
1164         /* compute multiplier, update diag(k) and U(i,k) */
1165         ili = il[i];  /* index of first nonzero element in U(i,k:bms-1) */
1166         uikdi = - ba[ili]*ba[bi[i]];  /* diagonal(k) */
1167         dk += uikdi*ba[ili];
1168         ba[ili] = uikdi; /* -U(i,k) */
1169 
1170         /* add multiple of row i to k-th row */
1171         jmin = ili + 1; jmax = bi[i+1];
1172         if (jmin < jmax){
1173           for (j=jmin; j<jmax; j++) rtmp[bj[j]] += uikdi*ba[j];
1174           /* update il and jl for row i */
1175           il[i] = jmin;
1176           j = bj[jmin]; jl[i] = jl[j]; jl[j] = i;
1177         }
1178         i = nexti;
1179       }
1180 
1181       /* shift the diagonals when zero pivot is detected */
1182       /* compute rs=sum of abs(off-diagonal) */
1183       rs   = 0.0;
1184       jmin = bi[k]+1;
1185       nz   = bi[k+1] - jmin;
1186       if (nz){
1187         bcol = bj + jmin;
1188         while (nz--){
1189           rs += PetscAbsScalar(rtmp[*bcol]);
1190           bcol++;
1191         }
1192       }
1193 
1194       sctx.rs = rs;
1195       sctx.pv = dk;
1196       ierr = MatCholeskyCheckShift_inline(info,sctx,newshift);CHKERRQ(ierr);
1197       if (newshift == 1){
1198         break;    /* sctx.shift_amount is updated */
1199       } else if (newshift == -1){
1200         SETERRQ4(PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot row %D value %g tolerance %g * rs %g",k,PetscAbsScalar(dk),zeropivot,rs);
1201       }
1202 
1203       /* copy data into U(k,:) */
1204       ba[bi[k]] = 1.0/dk; /* U(k,k) */
1205       jmin = bi[k]+1; jmax = bi[k+1];
1206       if (jmin < jmax) {
1207         for (j=jmin; j<jmax; j++){
1208           col = bj[j]; ba[j] = rtmp[col]; rtmp[col] = 0.0;
1209         }
1210         /* add the k-th row into il and jl */
1211         il[k] = jmin;
1212         i = bj[jmin]; jl[k] = jl[i]; jl[i] = k;
1213       }
1214     }
1215   } while (sctx.chshift);
1216   ierr = PetscFree(il);CHKERRQ(ierr);
1217   if (a->permute){ierr = PetscFree(aa);CHKERRQ(ierr);}
1218 
1219   ierr = ISRestoreIndices(ip,&rip);CHKERRQ(ierr);
1220   C->factor       = FACTOR_CHOLESKY;
1221   C->assembled    = PETSC_TRUE;
1222   C->preallocated = PETSC_TRUE;
1223   ierr = PetscLogFlops(C->m);CHKERRQ(ierr);
1224     if (sctx.nshift){
1225     if (shiftnz) {
1226       ierr = PetscLogInfo((0,"MatCholeskyFactorNumeric_SeqSBAIJ_1: number of shiftnz tries %D, shift_amount %g\n",sctx.nshift,sctx.shift_amount));CHKERRQ(ierr);
1227     } else if (shiftpd) {
1228       ierr = PetscLogInfo((0,"MatCholeskyFactorNumeric_SeqSBAIJ_1: number of shiftpd tries %D, shift_amount %g\n",sctx.nshift,sctx.shift_amount));CHKERRQ(ierr);
1229     }
1230   }
1231   PetscFunctionReturn(0);
1232 }
1233 
1234 /*
1235   Version for when blocks are 1 by 1 Using natural ordering
1236 */
1237 #undef __FUNCT__
1238 #define __FUNCT__ "MatCholeskyFactorNumeric_SeqSBAIJ_1_NaturalOrdering"
1239 PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_1_NaturalOrdering(Mat A,MatFactorInfo *info,Mat *B)
1240 {
1241   Mat            C = *B;
1242   Mat_SeqSBAIJ   *a=(Mat_SeqSBAIJ*)A->data,*b=(Mat_SeqSBAIJ *)C->data;
1243   PetscErrorCode ierr;
1244   PetscInt       i,j,mbs = a->mbs;
1245   PetscInt       *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j;
1246   PetscInt       k,jmin,*jl,*il,nexti,ili,*acol,*bcol,nz;
1247   MatScalar      *rtmp,*ba=b->a,*aa=a->a,dk,uikdi,*aval,*bval;
1248   PetscReal      zeropivot,rs,shiftnz;
1249   PetscTruth     shiftpd;
1250   ChShift_Ctx    sctx;
1251   PetscInt       newshift;
1252 
1253   PetscFunctionBegin;
1254   /* initialization */
1255   shiftnz   = info->shiftnz;
1256   shiftpd   = info->shiftpd;
1257   zeropivot = info->zeropivot;
1258 
1259   /* il and jl record the first nonzero element in each row of the accessing
1260      window U(0:k, k:mbs-1).
1261      jl:    list of rows to be added to uneliminated rows
1262             i>= k: jl(i) is the first row to be added to row i
1263             i<  k: jl(i) is the row following row i in some list of rows
1264             jl(i) = mbs indicates the end of a list
1265      il(i): points to the first nonzero element in U(i,k:mbs-1)
1266   */
1267   ierr = PetscMalloc(mbs*sizeof(MatScalar),&rtmp);CHKERRQ(ierr);
1268   ierr = PetscMalloc(2*mbs*sizeof(PetscInt),&il);CHKERRQ(ierr);
1269   jl   = il + mbs;
1270 
1271   sctx.shift_amount = 0;
1272   sctx.nshift       = 0;
1273   do {
1274     sctx.chshift = PETSC_FALSE;
1275     for (i=0; i<mbs; i++) {
1276       rtmp[i] = 0.0; jl[i] = mbs; il[0] = 0;
1277     }
1278 
1279     for (k = 0; k<mbs; k++){
1280       /*initialize k-th row with elements nonzero in row perm(k) of A */
1281       nz   = ai[k+1] - ai[k];
1282       acol = aj + ai[k];
1283       aval = aa + ai[k];
1284       bval = ba + bi[k];
1285       while (nz -- ){
1286         rtmp[*acol++] = *aval++;
1287         *bval++       = 0.0; /* for in-place factorization */
1288       }
1289 
1290       /* shift the diagonal of the matrix */
1291       if (sctx.nshift) rtmp[k] += sctx.shift_amount;
1292 
1293       /* modify k-th row by adding in those rows i with U(i,k)!=0 */
1294       dk = rtmp[k];
1295       i  = jl[k]; /* first row to be added to k_th row  */
1296 
1297       while (i < k){
1298         nexti = jl[i]; /* next row to be added to k_th row */
1299         /* compute multiplier, update D(k) and U(i,k) */
1300         ili   = il[i];  /* index of first nonzero element in U(i,k:bms-1) */
1301         uikdi = - ba[ili]*ba[bi[i]];
1302         dk   += uikdi*ba[ili];
1303         ba[ili] = uikdi; /* -U(i,k) */
1304 
1305         /* add multiple of row i to k-th row ... */
1306         jmin = ili + 1;
1307         nz   = bi[i+1] - jmin;
1308         if (nz > 0){
1309           bcol = bj + jmin;
1310           bval = ba + jmin;
1311           while (nz --) rtmp[*bcol++] += uikdi*(*bval++);
1312           /* update il and jl for i-th row */
1313           il[i] = jmin;
1314           j = bj[jmin]; jl[i] = jl[j]; jl[j] = i;
1315         }
1316         i = nexti;
1317       }
1318 
1319       /* shift the diagonals when zero pivot is detected */
1320       /* compute rs=sum of abs(off-diagonal) */
1321       rs   = 0.0;
1322       jmin = bi[k]+1;
1323       nz   = bi[k+1] - jmin;
1324       if (nz){
1325         bcol = bj + jmin;
1326         while (nz--){
1327           rs += PetscAbsScalar(rtmp[*bcol]);
1328           bcol++;
1329         }
1330       }
1331 
1332       sctx.rs = rs;
1333       sctx.pv = dk;
1334       ierr = MatCholeskyCheckShift_inline(info,sctx,newshift);CHKERRQ(ierr);
1335       if (newshift == 1){
1336         break;    /* sctx.shift_amount is updated */
1337       } else if (newshift == -1){
1338         SETERRQ4(PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot row %D value %g tolerance %g * rs %g",k,PetscAbsScalar(dk),zeropivot,rs);
1339       }
1340 
1341       /* copy data into U(k,:) */
1342       ba[bi[k]] = 1.0/dk;
1343       jmin      = bi[k]+1;
1344       nz        = bi[k+1] - jmin;
1345       if (nz){
1346         bcol = bj + jmin;
1347         bval = ba + jmin;
1348         while (nz--){
1349           *bval++       = rtmp[*bcol];
1350           rtmp[*bcol++] = 0.0;
1351         }
1352         /* add k-th row into il and jl */
1353         il[k] = jmin;
1354         i = bj[jmin]; jl[k] = jl[i]; jl[i] = k;
1355       }
1356     } /* end of for (k = 0; k<mbs; k++) */
1357   } while (sctx.chshift);
1358   ierr = PetscFree(rtmp);CHKERRQ(ierr);
1359   ierr = PetscFree(il);CHKERRQ(ierr);
1360 
1361   C->factor       = FACTOR_CHOLESKY;
1362   C->assembled    = PETSC_TRUE;
1363   C->preallocated = PETSC_TRUE;
1364   ierr = PetscLogFlops(C->m);CHKERRQ(ierr);
1365   if (sctx.nshift){
1366     if (shiftnz) {
1367       ierr = PetscLogInfo((0,"MatCholeskyFactorNumeric_SeqSBAIJ_1_NaturalOrdering: number of shiftnz tries %D, shift_amount %g\n",sctx.nshift,sctx.shift_amount));CHKERRQ(ierr);
1368     } else if (shiftpd) {
1369       ierr = PetscLogInfo((0,"MatCholeskyFactorNumeric_SeqSBAIJ_1_NaturalOrdering: number of shiftpd tries %D, shift_amount %g\n",sctx.nshift,sctx.shift_amount));CHKERRQ(ierr);
1370     }
1371   }
1372   PetscFunctionReturn(0);
1373 }
1374 
1375 #undef __FUNCT__
1376 #define __FUNCT__ "MatCholeskyFactor_SeqSBAIJ"
1377 PetscErrorCode MatCholeskyFactor_SeqSBAIJ(Mat A,IS perm,MatFactorInfo *info)
1378 {
1379   PetscErrorCode ierr;
1380   Mat            C;
1381 
1382   PetscFunctionBegin;
1383   ierr = MatCholeskyFactorSymbolic(A,perm,info,&C);CHKERRQ(ierr);
1384   ierr = MatCholeskyFactorNumeric(A,info,&C);CHKERRQ(ierr);
1385   ierr = MatHeaderCopy(A,C);CHKERRQ(ierr);
1386   PetscFunctionReturn(0);
1387 }
1388 
1389 
1390