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