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