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