Actual source code: sbaijfact.c

petsc-3.6.4 2016-04-12
Report Typos and Errors
  2: #include <../src/mat/impls/baij/seq/baij.h>
  3: #include <../src/mat/impls/sbaij/seq/sbaij.h>
  4: #include <petsc/private/kernels/blockinvert.h>
  5: #include <petscis.h>

  7: /*
  8:   input:
  9:    F -- numeric factor
 10:   output:
 11:    nneg, nzero, npos: matrix inertia
 12: */

 16: PetscErrorCode MatGetInertia_SeqSBAIJ(Mat F,PetscInt *nneig,PetscInt *nzero,PetscInt *npos)
 17: {
 18:   Mat_SeqSBAIJ *fact_ptr = (Mat_SeqSBAIJ*)F->data;
 19:   MatScalar    *dd       = fact_ptr->a;
 20:   PetscInt     mbs       =fact_ptr->mbs,bs=F->rmap->bs,i,nneig_tmp,npos_tmp,*fi = fact_ptr->diag;

 23:   if (bs != 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support for bs: %D >1 yet",bs);
 24:   nneig_tmp = 0; npos_tmp = 0;
 25:   for (i=0; i<mbs; i++) {
 26:     if (PetscRealPart(dd[*fi]) > 0.0) npos_tmp++;
 27:     else if (PetscRealPart(dd[*fi]) < 0.0) nneig_tmp++;
 28:     fi++;
 29:   }
 30:   if (nneig) *nneig = nneig_tmp;
 31:   if (npos)  *npos  = npos_tmp;
 32:   if (nzero) *nzero = mbs - nneig_tmp - npos_tmp;
 33:   return(0);
 34: }

 36: /*
 37:   Symbolic U^T*D*U factorization for SBAIJ format. Modified from SSF of YSMP.
 38:   Use Modified Sparse Row (MSR) storage for u and ju. See page 85, "Iterative Methods ..." by Saad.
 39: */
 42: PetscErrorCode MatCholeskyFactorSymbolic_SeqSBAIJ_MSR(Mat F,Mat A,IS perm,const MatFactorInfo *info)
 43: {
 44:   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)A->data,*b;
 46:   const PetscInt *rip,*ai,*aj;
 47:   PetscInt       i,mbs = a->mbs,*jutmp,bs = A->rmap->bs,bs2=a->bs2;
 48:   PetscInt       m,reallocs = 0,prow;
 49:   PetscInt       *jl,*q,jmin,jmax,juidx,nzk,qm,*iu,*ju,k,j,vj,umax,maxadd;
 50:   PetscReal      f = info->fill;
 51:   PetscBool      perm_identity;

 54:   /* check whether perm is the identity mapping */
 55:   ISIdentity(perm,&perm_identity);
 56:   ISGetIndices(perm,&rip);

 58:   if (perm_identity) { /* without permutation */
 59:     a->permute = PETSC_FALSE;

 61:     ai = a->i; aj = a->j;
 62:   } else {            /* non-trivial permutation */
 63:     a->permute = PETSC_TRUE;

 65:     MatReorderingSeqSBAIJ(A,perm);

 67:     ai = a->inew; aj = a->jnew;
 68:   }

 70:   /* initialization */
 71:   PetscMalloc1(mbs+1,&iu);
 72:   umax  = (PetscInt)(f*ai[mbs] + 1); umax += mbs + 1;
 73:   PetscMalloc1(umax,&ju);
 74:   iu[0] = mbs+1;
 75:   juidx = mbs + 1; /* index for ju */
 76:   /* jl linked list for pivot row -- linked list for col index */
 77:   PetscMalloc2(mbs,&jl,mbs,&q);
 78:   for (i=0; i<mbs; i++) {
 79:     jl[i] = mbs;
 80:     q[i]  = 0;
 81:   }

 83:   /* for each row k */
 84:   for (k=0; k<mbs; k++) {
 85:     for (i=0; i<mbs; i++) q[i] = 0;  /* to be removed! */
 86:     nzk  = 0; /* num. of nz blocks in k-th block row with diagonal block excluded */
 87:     q[k] = mbs;
 88:     /* initialize nonzero structure of k-th row to row rip[k] of A */
 89:     jmin = ai[rip[k]] +1; /* exclude diag[k] */
 90:     jmax = ai[rip[k]+1];
 91:     for (j=jmin; j<jmax; j++) {
 92:       vj = rip[aj[j]]; /* col. value */
 93:       if (vj > k) {
 94:         qm = k;
 95:         do {
 96:           m = qm; qm = q[m];
 97:         } while (qm < vj);
 98:         if (qm == vj) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Duplicate entry in A\n");
 99:         nzk++;
100:         q[m]  = vj;
101:         q[vj] = qm;
102:       } /* if (vj > k) */
103:     } /* for (j=jmin; j<jmax; j++) */

105:     /* modify nonzero structure of k-th row by computing fill-in
106:        for each row i to be merged in */
107:     prow = k;
108:     prow = jl[prow]; /* next pivot row (== mbs for symbolic factorization) */

110:     while (prow < k) {
111:       /* merge row prow into k-th row */
112:       jmin = iu[prow] + 1; jmax = iu[prow+1];
113:       qm   = k;
114:       for (j=jmin; j<jmax; j++) {
115:         vj = ju[j];
116:         do {
117:           m = qm; qm = q[m];
118:         } while (qm < vj);
119:         if (qm != vj) {
120:           nzk++; q[m] = vj; q[vj] = qm; qm = vj;
121:         }
122:       }
123:       prow = jl[prow]; /* next pivot row */
124:     }

126:     /* add k to row list for first nonzero element in k-th row */
127:     if (nzk > 0) {
128:       i     = q[k]; /* col value of first nonzero element in U(k, k+1:mbs-1) */
129:       jl[k] = jl[i]; jl[i] = k;
130:     }
131:     iu[k+1] = iu[k] + nzk;

133:     /* allocate more space to ju if needed */
134:     if (iu[k+1] > umax) {
135:       /* estimate how much additional space we will need */
136:       /* use the strategy suggested by David Hysom <hysom@perch-t.icase.edu> */
137:       /* just double the memory each time */
138:       maxadd = umax;
139:       if (maxadd < nzk) maxadd = (mbs-k)*(nzk+1)/2;
140:       umax += maxadd;

142:       /* allocate a longer ju */
143:       PetscMalloc1(umax,&jutmp);
144:       PetscMemcpy(jutmp,ju,iu[k]*sizeof(PetscInt));
145:       PetscFree(ju);
146:       ju   = jutmp;
147:       reallocs++; /* count how many times we realloc */
148:     }

150:     /* save nonzero structure of k-th row in ju */
151:     i=k;
152:     while (nzk--) {
153:       i           = q[i];
154:       ju[juidx++] = i;
155:     }
156:   }

158: #if defined(PETSC_USE_INFO)
159:   if (ai[mbs] != 0) {
160:     PetscReal af = ((PetscReal)iu[mbs])/((PetscReal)ai[mbs]);
161:     PetscInfo3(A,"Reallocs %D Fill ratio:given %g needed %g\n",reallocs,(double)f,(double)af);
162:     PetscInfo1(A,"Run with -pc_factor_fill %g or use \n",(double)af);
163:     PetscInfo1(A,"PCFactorSetFill(pc,%g);\n",(double)af);
164:     PetscInfo(A,"for best performance.\n");
165:   } else {
166:     PetscInfo(A,"Empty matrix.\n");
167:   }
168: #endif

170:   ISRestoreIndices(perm,&rip);
171:   PetscFree2(jl,q);

173:   /* put together the new matrix */
174:   MatSeqSBAIJSetPreallocation_SeqSBAIJ(F,bs,MAT_SKIP_ALLOCATION,NULL);

176:   /* PetscLogObjectParent((PetscObject)B,(PetscObject)iperm); */
177:   b                = (Mat_SeqSBAIJ*)(F)->data;
178:   b->singlemalloc  = PETSC_FALSE;
179:   b->free_a        = PETSC_TRUE;
180:   b->free_ij       = PETSC_TRUE;

182:   PetscMalloc1((iu[mbs]+1)*bs2,&b->a);
183:   b->j    = ju;
184:   b->i    = iu;
185:   b->diag = 0;
186:   b->ilen = 0;
187:   b->imax = 0;
188:   b->row  = perm;

190:   b->pivotinblocks = PETSC_FALSE; /* need to get from MatFactorInfo */

192:   PetscObjectReference((PetscObject)perm);

194:   b->icol = perm;
195:   PetscObjectReference((PetscObject)perm);
196:   PetscMalloc1(bs*mbs+bs,&b->solve_work);
197:   /* In b structure:  Free imax, ilen, old a, old j.
198:      Allocate idnew, solve_work, new a, new j */
199:   PetscLogObjectMemory((PetscObject)F,(iu[mbs]-mbs)*(sizeof(PetscInt)+sizeof(MatScalar)));
200:   b->maxnz = b->nz = iu[mbs];

202:   (F)->info.factor_mallocs   = reallocs;
203:   (F)->info.fill_ratio_given = f;
204:   if (ai[mbs] != 0) {
205:     (F)->info.fill_ratio_needed = ((PetscReal)iu[mbs])/((PetscReal)ai[mbs]);
206:   } else {
207:     (F)->info.fill_ratio_needed = 0.0;
208:   }
209:   MatSeqSBAIJSetNumericFactorization_inplace(F,perm_identity);
210:   return(0);
211: }
212: /*
213:     Symbolic U^T*D*U factorization for SBAIJ format.
214:     See MatICCFactorSymbolic_SeqAIJ() for description of its data structure.
215: */
216: #include <petscbt.h>
217: #include <../src/mat/utils/freespace.h>
220: PetscErrorCode MatCholeskyFactorSymbolic_SeqSBAIJ(Mat fact,Mat A,IS perm,const MatFactorInfo *info)
221: {
222:   Mat_SeqSBAIJ       *a = (Mat_SeqSBAIJ*)A->data;
223:   Mat_SeqSBAIJ       *b;
224:   PetscErrorCode     ierr;
225:   PetscBool          perm_identity,missing;
226:   PetscReal          fill = info->fill;
227:   const PetscInt     *rip,*ai=a->i,*aj=a->j;
228:   PetscInt           i,mbs=a->mbs,bs=A->rmap->bs,reallocs=0,prow;
229:   PetscInt           *jl,jmin,jmax,nzk,*ui,k,j,*il,nextprow;
230:   PetscInt           nlnk,*lnk,ncols,*cols,*uj,**ui_ptr,*uj_ptr,*udiag;
231:   PetscFreeSpaceList free_space=NULL,current_space=NULL;
232:   PetscBT            lnkbt;

235:   if (A->rmap->n != A->cmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Must be square matrix, rows %D columns %D",A->rmap->n,A->cmap->n);
236:   MatMissingDiagonal(A,&missing,&i);
237:   if (missing) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Matrix is missing diagonal entry %D",i);
238:   if (bs > 1) {
239:     MatCholeskyFactorSymbolic_SeqSBAIJ_inplace(fact,A,perm,info);
240:     return(0);
241:   }

243:   /* check whether perm is the identity mapping */
244:   ISIdentity(perm,&perm_identity);
245:   if (!perm_identity) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Matrix reordering is not supported for sbaij matrix. Use aij format");
246:   a->permute = PETSC_FALSE;
247:   ISGetIndices(perm,&rip);

249:   /* initialization */
250:   PetscMalloc1(mbs+1,&ui);
251:   PetscMalloc1(mbs+1,&udiag);
252:   ui[0] = 0;

254:   /* jl: linked list for storing indices of the pivot rows
255:      il: il[i] points to the 1st nonzero entry of U(i,k:mbs-1) */
256:   PetscMalloc4(mbs,&ui_ptr,mbs,&il,mbs,&jl,mbs,&cols);
257:   for (i=0; i<mbs; i++) {
258:     jl[i] = mbs; il[i] = 0;
259:   }

261:   /* create and initialize a linked list for storing column indices of the active row k */
262:   nlnk = mbs + 1;
263:   PetscLLCreate(mbs,mbs,nlnk,lnk,lnkbt);

265:   /* initial FreeSpace size is fill*(ai[mbs]+1) */
266:   PetscFreeSpaceGet((PetscInt)(fill*(ai[mbs]+1)),&free_space);
267:   current_space = free_space;

269:   for (k=0; k<mbs; k++) {  /* for each active row k */
270:     /* initialize lnk by the column indices of row rip[k] of A */
271:     nzk   = 0;
272:     ncols = ai[k+1] - ai[k];
273:     if (!ncols) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_MAT_CH_ZRPVT,"Empty row %D in matrix ",k);
274:     for (j=0; j<ncols; j++) {
275:       i       = *(aj + ai[k] + j);
276:       cols[j] = i;
277:     }
278:     PetscLLAdd(ncols,cols,mbs,nlnk,lnk,lnkbt);
279:     nzk += nlnk;

281:     /* update lnk by computing fill-in for each pivot row to be merged in */
282:     prow = jl[k]; /* 1st pivot row */

284:     while (prow < k) {
285:       nextprow = jl[prow];
286:       /* merge prow into k-th row */
287:       jmin   = il[prow] + 1; /* index of the 2nd nzero entry in U(prow,k:mbs-1) */
288:       jmax   = ui[prow+1];
289:       ncols  = jmax-jmin;
290:       uj_ptr = ui_ptr[prow] + jmin - ui[prow]; /* points to the 2nd nzero entry in U(prow,k:mbs-1) */
291:       PetscLLAddSorted(ncols,uj_ptr,mbs,nlnk,lnk,lnkbt);
292:       nzk   += nlnk;

294:       /* update il and jl for prow */
295:       if (jmin < jmax) {
296:         il[prow] = jmin;
297:         j        = *uj_ptr; jl[prow] = jl[j]; jl[j] = prow;
298:       }
299:       prow = nextprow;
300:     }

302:     /* if free space is not available, make more free space */
303:     if (current_space->local_remaining<nzk) {
304:       i    = mbs - k + 1; /* num of unfactored rows */
305:       i   *= PetscMin(nzk, i-1); /* i*nzk, i*(i-1): estimated and max additional space needed */
306:       PetscFreeSpaceGet(i,&current_space);
307:       reallocs++;
308:     }

310:     /* copy data into free space, then initialize lnk */
311:     PetscLLClean(mbs,mbs,nzk,lnk,current_space->array,lnkbt);

313:     /* add the k-th row into il and jl */
314:     if (nzk > 1) {
315:       i     = current_space->array[1]; /* col value of the first nonzero element in U(k, k+1:mbs-1) */
316:       jl[k] = jl[i]; jl[i] = k;
317:       il[k] = ui[k] + 1;
318:     }
319:     ui_ptr[k] = current_space->array;

321:     current_space->array           += nzk;
322:     current_space->local_used      += nzk;
323:     current_space->local_remaining -= nzk;

325:     ui[k+1] = ui[k] + nzk;
326:   }

328:   ISRestoreIndices(perm,&rip);
329:   PetscFree4(ui_ptr,il,jl,cols);

331:   /* destroy list of free space and other temporary array(s) */
332:   PetscMalloc1(ui[mbs]+1,&uj);
333:   PetscFreeSpaceContiguous_Cholesky(&free_space,uj,mbs,ui,udiag); /* store matrix factor */
334:   PetscLLDestroy(lnk,lnkbt);

336:   /* put together the new matrix in MATSEQSBAIJ format */
337:   MatSeqSBAIJSetPreallocation_SeqSBAIJ(fact,bs,MAT_SKIP_ALLOCATION,NULL);

339:   b               = (Mat_SeqSBAIJ*)fact->data;
340:   b->singlemalloc = PETSC_FALSE;
341:   b->free_a       = PETSC_TRUE;
342:   b->free_ij      = PETSC_TRUE;

344:   PetscMalloc1(ui[mbs]+1,&b->a);

346:   b->j         = uj;
347:   b->i         = ui;
348:   b->diag      = udiag;
349:   b->free_diag = PETSC_TRUE;
350:   b->ilen      = 0;
351:   b->imax      = 0;
352:   b->row       = perm;
353:   b->icol      = perm;

355:   PetscObjectReference((PetscObject)perm);
356:   PetscObjectReference((PetscObject)perm);

358:   b->pivotinblocks = PETSC_FALSE; /* need to get from MatFactorInfo */

360:   PetscMalloc1(mbs+1,&b->solve_work);
361:   PetscLogObjectMemory((PetscObject)fact,ui[mbs]*(sizeof(PetscInt)+sizeof(MatScalar)));

363:   b->maxnz = b->nz = ui[mbs];

365:   fact->info.factor_mallocs   = reallocs;
366:   fact->info.fill_ratio_given = fill;
367:   if (ai[mbs] != 0) {
368:     fact->info.fill_ratio_needed = ((PetscReal)ui[mbs])/ai[mbs];
369:   } else {
370:     fact->info.fill_ratio_needed = 0.0;
371:   }
372: #if defined(PETSC_USE_INFO)
373:   if (ai[mbs] != 0) {
374:     PetscReal af = fact->info.fill_ratio_needed;
375:     PetscInfo3(A,"Reallocs %D Fill ratio:given %g needed %g\n",reallocs,(double)fill,(double)af);
376:     PetscInfo1(A,"Run with -pc_factor_fill %g or use \n",(double)af);
377:     PetscInfo1(A,"PCFactorSetFill(pc,%g) for best performance.\n",(double)af);
378:   } else {
379:     PetscInfo(A,"Empty matrix.\n");
380:   }
381: #endif
382:   fact->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_1_NaturalOrdering;
383:   return(0);
384: }

388: PetscErrorCode MatCholeskyFactorSymbolic_SeqSBAIJ_inplace(Mat fact,Mat A,IS perm,const MatFactorInfo *info)
389: {
390:   Mat_SeqSBAIJ       *a = (Mat_SeqSBAIJ*)A->data;
391:   Mat_SeqSBAIJ       *b;
392:   PetscErrorCode     ierr;
393:   PetscBool          perm_identity,missing;
394:   PetscReal          fill = info->fill;
395:   const PetscInt     *rip,*ai,*aj;
396:   PetscInt           i,mbs=a->mbs,bs=A->rmap->bs,reallocs=0,prow,d;
397:   PetscInt           *jl,jmin,jmax,nzk,*ui,k,j,*il,nextprow;
398:   PetscInt           nlnk,*lnk,ncols,*cols,*uj,**ui_ptr,*uj_ptr;
399:   PetscFreeSpaceList free_space=NULL,current_space=NULL;
400:   PetscBT            lnkbt;

403:   MatMissingDiagonal(A,&missing,&d);
404:   if (missing) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Matrix is missing diagonal entry %D",d);

406:   /*
407:    This code originally uses Modified Sparse Row (MSR) storage
408:    (see page 85, "Iterative Methods ..." by Saad) for the output matrix B - bad choise!
409:    Then it is rewritten so the factor B takes seqsbaij format. However the associated
410:    MatCholeskyFactorNumeric_() have not been modified for the cases of bs>1 or !perm_identity,
411:    thus the original code in MSR format is still used for these cases.
412:    The code below should replace MatCholeskyFactorSymbolic_SeqSBAIJ_MSR() whenever
413:    MatCholeskyFactorNumeric_() is modified for using sbaij symbolic factor.
414:   */
415:   if (bs > 1) {
416:     MatCholeskyFactorSymbolic_SeqSBAIJ_MSR(fact,A,perm,info);
417:     return(0);
418:   }

420:   /* check whether perm is the identity mapping */
421:   ISIdentity(perm,&perm_identity);

423:   if (perm_identity) {
424:     a->permute = PETSC_FALSE;

426:     ai = a->i; aj = a->j;
427:   } else {
428:     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Matrix reordering is not supported for sbaij matrix. Use aij format");
429: #if 0
430:     /* There are bugs for reordeing. Needs further work.
431:        MatReordering for sbaij cannot be efficient. User should use aij formt! */
432:     a->permute = PETSC_TRUE;

434:     MatReorderingSeqSBAIJ(A,perm);
435:     ai   = a->inew; aj = a->jnew;
436: #endif
437:   }
438:   ISGetIndices(perm,&rip);

440:   /* initialization */
441:   PetscMalloc1(mbs+1,&ui);
442:   ui[0] = 0;

444:   /* jl: linked list for storing indices of the pivot rows
445:      il: il[i] points to the 1st nonzero entry of U(i,k:mbs-1) */
446:   PetscMalloc4(mbs,&ui_ptr,mbs,&il,mbs,&jl,mbs,&cols);
447:   for (i=0; i<mbs; i++) {
448:     jl[i] = mbs; il[i] = 0;
449:   }

451:   /* create and initialize a linked list for storing column indices of the active row k */
452:   nlnk = mbs + 1;
453:   PetscLLCreate(mbs,mbs,nlnk,lnk,lnkbt);

455:   /* initial FreeSpace size is fill*(ai[mbs]+1) */
456:   PetscFreeSpaceGet((PetscInt)(fill*(ai[mbs]+1)),&free_space);
457:   current_space = free_space;

459:   for (k=0; k<mbs; k++) {  /* for each active row k */
460:     /* initialize lnk by the column indices of row rip[k] of A */
461:     nzk   = 0;
462:     ncols = ai[rip[k]+1] - ai[rip[k]];
463:     for (j=0; j<ncols; j++) {
464:       i       = *(aj + ai[rip[k]] + j);
465:       cols[j] = rip[i];
466:     }
467:     PetscLLAdd(ncols,cols,mbs,nlnk,lnk,lnkbt);
468:     nzk += nlnk;

470:     /* update lnk by computing fill-in for each pivot row to be merged in */
471:     prow = jl[k]; /* 1st pivot row */

473:     while (prow < k) {
474:       nextprow = jl[prow];
475:       /* merge prow into k-th row */
476:       jmin   = il[prow] + 1; /* index of the 2nd nzero entry in U(prow,k:mbs-1) */
477:       jmax   = ui[prow+1];
478:       ncols  = jmax-jmin;
479:       uj_ptr = ui_ptr[prow] + jmin - ui[prow]; /* points to the 2nd nzero entry in U(prow,k:mbs-1) */
480:       PetscLLAddSorted(ncols,uj_ptr,mbs,nlnk,lnk,lnkbt);
481:       nzk   += nlnk;

483:       /* update il and jl for prow */
484:       if (jmin < jmax) {
485:         il[prow] = jmin;

487:         j = *uj_ptr; jl[prow] = jl[j]; jl[j] = prow;
488:       }
489:       prow = nextprow;
490:     }

492:     /* if free space is not available, make more free space */
493:     if (current_space->local_remaining<nzk) {
494:       i    = mbs - k + 1; /* num of unfactored rows */
495:       i    = PetscMin(i*nzk, i*(i-1)); /* i*nzk, i*(i-1): estimated and max additional space needed */
496:       PetscFreeSpaceGet(i,&current_space);
497:       reallocs++;
498:     }

500:     /* copy data into free space, then initialize lnk */
501:     PetscLLClean(mbs,mbs,nzk,lnk,current_space->array,lnkbt);

503:     /* add the k-th row into il and jl */
504:     if (nzk-1 > 0) {
505:       i     = current_space->array[1]; /* col value of the first nonzero element in U(k, k+1:mbs-1) */
506:       jl[k] = jl[i]; jl[i] = k;
507:       il[k] = ui[k] + 1;
508:     }
509:     ui_ptr[k] = current_space->array;

511:     current_space->array           += nzk;
512:     current_space->local_used      += nzk;
513:     current_space->local_remaining -= nzk;

515:     ui[k+1] = ui[k] + nzk;
516:   }

518:   ISRestoreIndices(perm,&rip);
519:   PetscFree4(ui_ptr,il,jl,cols);

521:   /* destroy list of free space and other temporary array(s) */
522:   PetscMalloc1(ui[mbs]+1,&uj);
523:   PetscFreeSpaceContiguous(&free_space,uj);
524:   PetscLLDestroy(lnk,lnkbt);

526:   /* put together the new matrix in MATSEQSBAIJ format */
527:   MatSeqSBAIJSetPreallocation_SeqSBAIJ(fact,bs,MAT_SKIP_ALLOCATION,NULL);

529:   b               = (Mat_SeqSBAIJ*)fact->data;
530:   b->singlemalloc = PETSC_FALSE;
531:   b->free_a       = PETSC_TRUE;
532:   b->free_ij      = PETSC_TRUE;

534:   PetscMalloc1(ui[mbs]+1,&b->a);

536:   b->j    = uj;
537:   b->i    = ui;
538:   b->diag = 0;
539:   b->ilen = 0;
540:   b->imax = 0;
541:   b->row  = perm;

543:   b->pivotinblocks = PETSC_FALSE; /* need to get from MatFactorInfo */

545:   PetscObjectReference((PetscObject)perm);
546:   b->icol  = perm;
547:   PetscObjectReference((PetscObject)perm);
548:   PetscMalloc1(mbs+1,&b->solve_work);
549:   PetscLogObjectMemory((PetscObject)fact,(ui[mbs]-mbs)*(sizeof(PetscInt)+sizeof(MatScalar)));
550:   b->maxnz = b->nz = ui[mbs];

552:   fact->info.factor_mallocs   = reallocs;
553:   fact->info.fill_ratio_given = fill;
554:   if (ai[mbs] != 0) {
555:     fact->info.fill_ratio_needed = ((PetscReal)ui[mbs])/ai[mbs];
556:   } else {
557:     fact->info.fill_ratio_needed = 0.0;
558:   }
559: #if defined(PETSC_USE_INFO)
560:   if (ai[mbs] != 0) {
561:     PetscReal af = fact->info.fill_ratio_needed;
562:     PetscInfo3(A,"Reallocs %D Fill ratio:given %g needed %g\n",reallocs,(double)fill,(double)af);
563:     PetscInfo1(A,"Run with -pc_factor_fill %g or use \n",(double)af);
564:     PetscInfo1(A,"PCFactorSetFill(pc,%g) for best performance.\n",(double)af);
565:   } else {
566:     PetscInfo(A,"Empty matrix.\n");
567:   }
568: #endif
569:   MatSeqSBAIJSetNumericFactorization_inplace(fact,perm_identity);
570:   return(0);
571: }

575: PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_N(Mat C,Mat A,const MatFactorInfo *info)
576: {
577:   Mat_SeqSBAIJ   *a   = (Mat_SeqSBAIJ*)A->data,*b = (Mat_SeqSBAIJ*)C->data;
578:   IS             perm = b->row;
580:   const PetscInt *ai,*aj,*perm_ptr,mbs=a->mbs,*bi=b->i,*bj=b->j;
581:   PetscInt       i,j;
582:   PetscInt       *a2anew,k,k1,jmin,jmax,*jl,*il,vj,nexti,ili;
583:   PetscInt       bs  =A->rmap->bs,bs2 = a->bs2,bslog = 0;
584:   MatScalar      *ba = b->a,*aa,*ap,*dk,*uik;
585:   MatScalar      *u,*diag,*rtmp,*rtmp_ptr;
586:   MatScalar      *work;
587:   PetscInt       *pivots;

590:   /* initialization */
591:   PetscCalloc1(bs2*mbs,&rtmp);
592:   PetscMalloc2(mbs,&il,mbs,&jl);
593:   for (i=0; i<mbs; i++) {
594:     jl[i] = mbs; il[0] = 0;
595:   }
596:   PetscMalloc3(bs2,&dk,bs2,&uik,bs,&work);
597:   PetscMalloc1(bs,&pivots);

599:   ISGetIndices(perm,&perm_ptr);

601:   /* check permutation */
602:   if (!a->permute) {
603:     ai = a->i; aj = a->j; aa = a->a;
604:   } else {
605:     ai   = a->inew; aj = a->jnew;
606:     PetscMalloc1(bs2*ai[mbs],&aa);
607:     PetscMemcpy(aa,a->a,bs2*ai[mbs]*sizeof(MatScalar));
608:     PetscMalloc1(ai[mbs],&a2anew);
609:     PetscMemcpy(a2anew,a->a2anew,(ai[mbs])*sizeof(PetscInt));

611:     /* flops in while loop */
612:     bslog = 2*bs*bs2;

614:     for (i=0; i<mbs; i++) {
615:       jmin = ai[i]; jmax = ai[i+1];
616:       for (j=jmin; j<jmax; j++) {
617:         while (a2anew[j] != j) {
618:           k = a2anew[j]; a2anew[j] = a2anew[k]; a2anew[k] = k;
619:           for (k1=0; k1<bs2; k1++) {
620:             dk[k1]       = aa[k*bs2+k1];
621:             aa[k*bs2+k1] = aa[j*bs2+k1];
622:             aa[j*bs2+k1] = dk[k1];
623:           }
624:         }
625:         /* transform columnoriented blocks that lie in the lower triangle to roworiented blocks */
626:         if (i > aj[j]) {
627:           /* printf("change orientation, row: %d, col: %d\n",i,aj[j]); */
628:           ap = aa + j*bs2;                     /* ptr to the beginning of j-th block of aa */
629:           for (k=0; k<bs2; k++) dk[k] = ap[k]; /* dk <- j-th block of aa */
630:           for (k=0; k<bs; k++) {               /* j-th block of aa <- dk^T */
631:             for (k1=0; k1<bs; k1++) *ap++ = dk[k + bs*k1];
632:           }
633:         }
634:       }
635:     }
636:     PetscFree(a2anew);
637:   }

639:   /* for each row k */
640:   for (k = 0; k<mbs; k++) {

642:     /*initialize k-th row with elements nonzero in row perm(k) of A */
643:     jmin = ai[perm_ptr[k]]; jmax = ai[perm_ptr[k]+1];

645:     ap = aa + jmin*bs2;
646:     for (j = jmin; j < jmax; j++) {
647:       vj       = perm_ptr[aj[j]];   /* block col. index */
648:       rtmp_ptr = rtmp + vj*bs2;
649:       for (i=0; i<bs2; i++) *rtmp_ptr++ = *ap++;
650:     }

652:     /* modify k-th row by adding in those rows i with U(i,k) != 0 */
653:     PetscMemcpy(dk,rtmp+k*bs2,bs2*sizeof(MatScalar));
654:     i    = jl[k]; /* first row to be added to k_th row  */

656:     while (i < k) {
657:       nexti = jl[i]; /* next row to be added to k_th row */

659:       /* compute multiplier */
660:       ili = il[i];  /* index of first nonzero element in U(i,k:bms-1) */

662:       /* uik = -inv(Di)*U_bar(i,k) */
663:       diag = ba + i*bs2;
664:       u    = ba + ili*bs2;
665:       PetscMemzero(uik,bs2*sizeof(MatScalar));
666:       PetscKernel_A_gets_A_minus_B_times_C(bs,uik,diag,u);

668:       /* update D(k) += -U(i,k)^T * U_bar(i,k) */
669:       PetscKernel_A_gets_A_plus_Btranspose_times_C(bs,dk,uik,u);
670:       PetscLogFlops(bslog*2.0);

672:       /* update -U(i,k) */
673:       PetscMemcpy(ba+ili*bs2,uik,bs2*sizeof(MatScalar));

675:       /* add multiple of row i to k-th row ... */
676:       jmin = ili + 1; jmax = bi[i+1];
677:       if (jmin < jmax) {
678:         for (j=jmin; j<jmax; j++) {
679:           /* rtmp += -U(i,k)^T * U_bar(i,j) */
680:           rtmp_ptr = rtmp + bj[j]*bs2;
681:           u        = ba + j*bs2;
682:           PetscKernel_A_gets_A_plus_Btranspose_times_C(bs,rtmp_ptr,uik,u);
683:         }
684:         PetscLogFlops(bslog*(jmax-jmin));

686:         /* ... add i to row list for next nonzero entry */
687:         il[i] = jmin;             /* update il(i) in column k+1, ... mbs-1 */
688:         j     = bj[jmin];
689:         jl[i] = jl[j]; jl[j] = i; /* update jl */
690:       }
691:       i = nexti;
692:     }

694:     /* save nonzero entries in k-th row of U ... */

696:     /* invert diagonal block */
697:     diag = ba+k*bs2;
698:     PetscMemcpy(diag,dk,bs2*sizeof(MatScalar));
699:     PetscKernel_A_gets_inverse_A(bs,diag,pivots,work);

701:     jmin = bi[k]; jmax = bi[k+1];
702:     if (jmin < jmax) {
703:       for (j=jmin; j<jmax; j++) {
704:         vj       = bj[j];      /* block col. index of U */
705:         u        = ba + j*bs2;
706:         rtmp_ptr = rtmp + vj*bs2;
707:         for (k1=0; k1<bs2; k1++) {
708:           *u++        = *rtmp_ptr;
709:           *rtmp_ptr++ = 0.0;
710:         }
711:       }

713:       /* ... add k to row list for first nonzero entry in k-th row */
714:       il[k] = jmin;
715:       i     = bj[jmin];
716:       jl[k] = jl[i]; jl[i] = k;
717:     }
718:   }

720:   PetscFree(rtmp);
721:   PetscFree2(il,jl);
722:   PetscFree3(dk,uik,work);
723:   PetscFree(pivots);
724:   if (a->permute) {
725:     PetscFree(aa);
726:   }

728:   ISRestoreIndices(perm,&perm_ptr);

730:   C->ops->solve          = MatSolve_SeqSBAIJ_N_inplace;
731:   C->ops->solvetranspose = MatSolve_SeqSBAIJ_N_inplace;
732:   C->ops->forwardsolve   = MatForwardSolve_SeqSBAIJ_N_inplace;
733:   C->ops->backwardsolve  = MatBackwardSolve_SeqSBAIJ_N_inplace;

735:   C->assembled    = PETSC_TRUE;
736:   C->preallocated = PETSC_TRUE;

738:   PetscLogFlops(1.3333*bs*bs2*b->mbs); /* from inverting diagonal blocks */
739:   return(0);
740: }

744: PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_N_NaturalOrdering(Mat C,Mat A,const MatFactorInfo *info)
745: {
746:   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)A->data,*b = (Mat_SeqSBAIJ*)C->data;
748:   PetscInt       i,j,mbs=a->mbs,*bi=b->i,*bj=b->j;
749:   PetscInt       *ai,*aj,k,k1,jmin,jmax,*jl,*il,vj,nexti,ili;
750:   PetscInt       bs  =A->rmap->bs,bs2 = a->bs2,bslog;
751:   MatScalar      *ba = b->a,*aa,*ap,*dk,*uik;
752:   MatScalar      *u,*diag,*rtmp,*rtmp_ptr;
753:   MatScalar      *work;
754:   PetscInt       *pivots;

757:   PetscCalloc1(bs2*mbs,&rtmp);
758:   PetscMalloc2(mbs,&il,mbs,&jl);
759:   for (i=0; i<mbs; i++) {
760:     jl[i] = mbs; il[0] = 0;
761:   }
762:   PetscMalloc3(bs2,&dk,bs2,&uik,bs,&work);
763:   PetscMalloc1(bs,&pivots);

765:   ai = a->i; aj = a->j; aa = a->a;

767:   /* flops in while loop */
768:   bslog = 2*bs*bs2;

770:   /* for each row k */
771:   for (k = 0; k<mbs; k++) {

773:     /*initialize k-th row with elements nonzero in row k of A */
774:     jmin = ai[k]; jmax = ai[k+1];
775:     ap   = aa + jmin*bs2;
776:     for (j = jmin; j < jmax; j++) {
777:       vj       = aj[j];   /* block col. index */
778:       rtmp_ptr = rtmp + vj*bs2;
779:       for (i=0; i<bs2; i++) *rtmp_ptr++ = *ap++;
780:     }

782:     /* modify k-th row by adding in those rows i with U(i,k) != 0 */
783:     PetscMemcpy(dk,rtmp+k*bs2,bs2*sizeof(MatScalar));
784:     i    = jl[k]; /* first row to be added to k_th row  */

786:     while (i < k) {
787:       nexti = jl[i]; /* next row to be added to k_th row */

789:       /* compute multiplier */
790:       ili = il[i];  /* index of first nonzero element in U(i,k:bms-1) */

792:       /* uik = -inv(Di)*U_bar(i,k) */
793:       diag = ba + i*bs2;
794:       u    = ba + ili*bs2;
795:       PetscMemzero(uik,bs2*sizeof(MatScalar));
796:       PetscKernel_A_gets_A_minus_B_times_C(bs,uik,diag,u);

798:       /* update D(k) += -U(i,k)^T * U_bar(i,k) */
799:       PetscKernel_A_gets_A_plus_Btranspose_times_C(bs,dk,uik,u);
800:       PetscLogFlops(bslog*2.0);

802:       /* update -U(i,k) */
803:       PetscMemcpy(ba+ili*bs2,uik,bs2*sizeof(MatScalar));

805:       /* add multiple of row i to k-th row ... */
806:       jmin = ili + 1; jmax = bi[i+1];
807:       if (jmin < jmax) {
808:         for (j=jmin; j<jmax; j++) {
809:           /* rtmp += -U(i,k)^T * U_bar(i,j) */
810:           rtmp_ptr = rtmp + bj[j]*bs2;
811:           u        = ba + j*bs2;
812:           PetscKernel_A_gets_A_plus_Btranspose_times_C(bs,rtmp_ptr,uik,u);
813:         }
814:         PetscLogFlops(bslog*(jmax-jmin));

816:         /* ... add i to row list for next nonzero entry */
817:         il[i] = jmin;             /* update il(i) in column k+1, ... mbs-1 */
818:         j     = bj[jmin];
819:         jl[i] = jl[j]; jl[j] = i; /* update jl */
820:       }
821:       i = nexti;
822:     }

824:     /* save nonzero entries in k-th row of U ... */

826:     /* invert diagonal block */
827:     diag = ba+k*bs2;
828:     PetscMemcpy(diag,dk,bs2*sizeof(MatScalar));
829:     PetscKernel_A_gets_inverse_A(bs,diag,pivots,work);

831:     jmin = bi[k]; jmax = bi[k+1];
832:     if (jmin < jmax) {
833:       for (j=jmin; j<jmax; j++) {
834:         vj       = bj[j];      /* block col. index of U */
835:         u        = ba + j*bs2;
836:         rtmp_ptr = rtmp + vj*bs2;
837:         for (k1=0; k1<bs2; k1++) {
838:           *u++        = *rtmp_ptr;
839:           *rtmp_ptr++ = 0.0;
840:         }
841:       }

843:       /* ... add k to row list for first nonzero entry in k-th row */
844:       il[k] = jmin;
845:       i     = bj[jmin];
846:       jl[k] = jl[i]; jl[i] = k;
847:     }
848:   }

850:   PetscFree(rtmp);
851:   PetscFree2(il,jl);
852:   PetscFree3(dk,uik,work);
853:   PetscFree(pivots);

855:   C->ops->solve          = MatSolve_SeqSBAIJ_N_NaturalOrdering_inplace;
856:   C->ops->solvetranspose = MatSolve_SeqSBAIJ_N_NaturalOrdering_inplace;
857:   C->ops->forwardsolve   = MatForwardSolve_SeqSBAIJ_N_NaturalOrdering_inplace;
858:   C->ops->backwardsolve  = MatBackwardSolve_SeqSBAIJ_N_NaturalOrdering_inplace;
859:   C->assembled           = PETSC_TRUE;
860:   C->preallocated        = PETSC_TRUE;

862:   PetscLogFlops(1.3333*bs*bs2*b->mbs); /* from inverting diagonal blocks */
863:   return(0);
864: }

866: /*
867:     Numeric U^T*D*U factorization for SBAIJ format. Modified from SNF of YSMP.
868:     Version for blocks 2 by 2.
869: */
872: PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_2(Mat C,Mat A,const MatFactorInfo *info)
873: {
874:   Mat_SeqSBAIJ   *a   = (Mat_SeqSBAIJ*)A->data,*b = (Mat_SeqSBAIJ*)C->data;
875:   IS             perm = b->row;
877:   const PetscInt *ai,*aj,*perm_ptr;
878:   PetscInt       i,j,mbs=a->mbs,*bi=b->i,*bj=b->j;
879:   PetscInt       *a2anew,k,k1,jmin,jmax,*jl,*il,vj,nexti,ili;
880:   MatScalar      *ba = b->a,*aa,*ap;
881:   MatScalar      *u,*diag,*rtmp,*rtmp_ptr,dk[4],uik[4];
882:   PetscReal      shift = info->shiftamount;

885:   /* initialization */
886:   /* il and jl record the first nonzero element in each row of the accessing
887:      window U(0:k, k:mbs-1).
888:      jl:    list of rows to be added to uneliminated rows
889:             i>= k: jl(i) is the first row to be added to row i
890:             i<  k: jl(i) is the row following row i in some list of rows
891:             jl(i) = mbs indicates the end of a list
892:      il(i): points to the first nonzero element in columns k,...,mbs-1 of
893:             row i of U */
894:   PetscCalloc1(4*mbs,&rtmp);
895:   PetscMalloc2(mbs,&il,mbs,&jl);
896:   for (i=0; i<mbs; i++) {
897:     jl[i] = mbs; il[0] = 0;
898:   }
899:   ISGetIndices(perm,&perm_ptr);

901:   /* check permutation */
902:   if (!a->permute) {
903:     ai = a->i; aj = a->j; aa = a->a;
904:   } else {
905:     ai   = a->inew; aj = a->jnew;
906:     PetscMalloc1(4*ai[mbs],&aa);
907:     PetscMemcpy(aa,a->a,4*ai[mbs]*sizeof(MatScalar));
908:     PetscMalloc1(ai[mbs],&a2anew);
909:     PetscMemcpy(a2anew,a->a2anew,(ai[mbs])*sizeof(PetscInt));

911:     for (i=0; i<mbs; i++) {
912:       jmin = ai[i]; jmax = ai[i+1];
913:       for (j=jmin; j<jmax; j++) {
914:         while (a2anew[j] != j) {
915:           k = a2anew[j]; a2anew[j] = a2anew[k]; a2anew[k] = k;
916:           for (k1=0; k1<4; k1++) {
917:             dk[k1]     = aa[k*4+k1];
918:             aa[k*4+k1] = aa[j*4+k1];
919:             aa[j*4+k1] = dk[k1];
920:           }
921:         }
922:         /* transform columnoriented blocks that lie in the lower triangle to roworiented blocks */
923:         if (i > aj[j]) {
924:           /* printf("change orientation, row: %d, col: %d\n",i,aj[j]); */
925:           ap    = aa + j*4;  /* ptr to the beginning of the block */
926:           dk[1] = ap[1];     /* swap ap[1] and ap[2] */
927:           ap[1] = ap[2];
928:           ap[2] = dk[1];
929:         }
930:       }
931:     }
932:     PetscFree(a2anew);
933:   }

935:   /* for each row k */
936:   for (k = 0; k<mbs; k++) {

938:     /*initialize k-th row with elements nonzero in row perm(k) of A */
939:     jmin = ai[perm_ptr[k]]; jmax = ai[perm_ptr[k]+1];
940:     ap   = aa + jmin*4;
941:     for (j = jmin; j < jmax; j++) {
942:       vj       = perm_ptr[aj[j]];   /* block col. index */
943:       rtmp_ptr = rtmp + vj*4;
944:       for (i=0; i<4; i++) *rtmp_ptr++ = *ap++;
945:     }

947:     /* modify k-th row by adding in those rows i with U(i,k) != 0 */
948:     PetscMemcpy(dk,rtmp+k*4,4*sizeof(MatScalar));
949:     i    = jl[k]; /* first row to be added to k_th row  */

951:     while (i < k) {
952:       nexti = jl[i]; /* next row to be added to k_th row */

954:       /* compute multiplier */
955:       ili = il[i];  /* index of first nonzero element in U(i,k:bms-1) */

957:       /* uik = -inv(Di)*U_bar(i,k): - ba[ili]*ba[i] */
958:       diag   = ba + i*4;
959:       u      = ba + ili*4;
960:       uik[0] = -(diag[0]*u[0] + diag[2]*u[1]);
961:       uik[1] = -(diag[1]*u[0] + diag[3]*u[1]);
962:       uik[2] = -(diag[0]*u[2] + diag[2]*u[3]);
963:       uik[3] = -(diag[1]*u[2] + diag[3]*u[3]);

965:       /* update D(k) += -U(i,k)^T * U_bar(i,k): dk += uik*ba[ili] */
966:       dk[0] += uik[0]*u[0] + uik[1]*u[1];
967:       dk[1] += uik[2]*u[0] + uik[3]*u[1];
968:       dk[2] += uik[0]*u[2] + uik[1]*u[3];
969:       dk[3] += uik[2]*u[2] + uik[3]*u[3];

971:       PetscLogFlops(16.0*2.0);

973:       /* update -U(i,k): ba[ili] = uik */
974:       PetscMemcpy(ba+ili*4,uik,4*sizeof(MatScalar));

976:       /* add multiple of row i to k-th row ... */
977:       jmin = ili + 1; jmax = bi[i+1];
978:       if (jmin < jmax) {
979:         for (j=jmin; j<jmax; j++) {
980:           /* rtmp += -U(i,k)^T * U_bar(i,j): rtmp[bj[j]] += uik*ba[j]; */
981:           rtmp_ptr     = rtmp + bj[j]*4;
982:           u            = ba + j*4;
983:           rtmp_ptr[0] += uik[0]*u[0] + uik[1]*u[1];
984:           rtmp_ptr[1] += uik[2]*u[0] + uik[3]*u[1];
985:           rtmp_ptr[2] += uik[0]*u[2] + uik[1]*u[3];
986:           rtmp_ptr[3] += uik[2]*u[2] + uik[3]*u[3];
987:         }
988:         PetscLogFlops(16.0*(jmax-jmin));

990:         /* ... add i to row list for next nonzero entry */
991:         il[i] = jmin;             /* update il(i) in column k+1, ... mbs-1 */
992:         j     = bj[jmin];
993:         jl[i] = jl[j]; jl[j] = i; /* update jl */
994:       }
995:       i = nexti;
996:     }

998:     /* save nonzero entries in k-th row of U ... */

1000:     /* invert diagonal block */
1001:     diag = ba+k*4;
1002:     PetscMemcpy(diag,dk,4*sizeof(MatScalar));
1003:     PetscKernel_A_gets_inverse_A_2(diag,shift);

1005:     jmin = bi[k]; jmax = bi[k+1];
1006:     if (jmin < jmax) {
1007:       for (j=jmin; j<jmax; j++) {
1008:         vj       = bj[j];      /* block col. index of U */
1009:         u        = ba + j*4;
1010:         rtmp_ptr = rtmp + vj*4;
1011:         for (k1=0; k1<4; k1++) {
1012:           *u++        = *rtmp_ptr;
1013:           *rtmp_ptr++ = 0.0;
1014:         }
1015:       }

1017:       /* ... add k to row list for first nonzero entry in k-th row */
1018:       il[k] = jmin;
1019:       i     = bj[jmin];
1020:       jl[k] = jl[i]; jl[i] = k;
1021:     }
1022:   }

1024:   PetscFree(rtmp);
1025:   PetscFree2(il,jl);
1026:   if (a->permute) {
1027:     PetscFree(aa);
1028:   }
1029:   ISRestoreIndices(perm,&perm_ptr);

1031:   C->ops->solve          = MatSolve_SeqSBAIJ_2_inplace;
1032:   C->ops->solvetranspose = MatSolve_SeqSBAIJ_2_inplace;
1033:   C->assembled           = PETSC_TRUE;
1034:   C->preallocated        = PETSC_TRUE;

1036:   PetscLogFlops(1.3333*8*b->mbs); /* from inverting diagonal blocks */
1037:   return(0);
1038: }

1040: /*
1041:       Version for when blocks are 2 by 2 Using natural ordering
1042: */
1045: PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_2_NaturalOrdering(Mat C,Mat A,const MatFactorInfo *info)
1046: {
1047:   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)A->data,*b = (Mat_SeqSBAIJ*)C->data;
1049:   PetscInt       i,j,mbs=a->mbs,*bi=b->i,*bj=b->j;
1050:   PetscInt       *ai,*aj,k,k1,jmin,jmax,*jl,*il,vj,nexti,ili;
1051:   MatScalar      *ba = b->a,*aa,*ap,dk[8],uik[8];
1052:   MatScalar      *u,*diag,*rtmp,*rtmp_ptr;
1053:   PetscReal      shift = info->shiftamount;

1056:   /* initialization */
1057:   /* il and jl record the first nonzero element in each row of the accessing
1058:      window U(0:k, k:mbs-1).
1059:      jl:    list of rows to be added to uneliminated rows
1060:             i>= k: jl(i) is the first row to be added to row i
1061:             i<  k: jl(i) is the row following row i in some list of rows
1062:             jl(i) = mbs indicates the end of a list
1063:      il(i): points to the first nonzero element in columns k,...,mbs-1 of
1064:             row i of U */
1065:   PetscCalloc1(4*mbs,&rtmp);
1066:   PetscMalloc2(mbs,&il,mbs,&jl);
1067:   for (i=0; i<mbs; i++) {
1068:     jl[i] = mbs; il[0] = 0;
1069:   }
1070:   ai = a->i; aj = a->j; aa = a->a;

1072:   /* for each row k */
1073:   for (k = 0; k<mbs; k++) {

1075:     /*initialize k-th row with elements nonzero in row k of A */
1076:     jmin = ai[k]; jmax = ai[k+1];
1077:     ap   = aa + jmin*4;
1078:     for (j = jmin; j < jmax; j++) {
1079:       vj       = aj[j];   /* block col. index */
1080:       rtmp_ptr = rtmp + vj*4;
1081:       for (i=0; i<4; i++) *rtmp_ptr++ = *ap++;
1082:     }

1084:     /* modify k-th row by adding in those rows i with U(i,k) != 0 */
1085:     PetscMemcpy(dk,rtmp+k*4,4*sizeof(MatScalar));
1086:     i    = jl[k]; /* first row to be added to k_th row  */

1088:     while (i < k) {
1089:       nexti = jl[i]; /* next row to be added to k_th row */

1091:       /* compute multiplier */
1092:       ili = il[i];  /* index of first nonzero element in U(i,k:bms-1) */

1094:       /* uik = -inv(Di)*U_bar(i,k): - ba[ili]*ba[i] */
1095:       diag   = ba + i*4;
1096:       u      = ba + ili*4;
1097:       uik[0] = -(diag[0]*u[0] + diag[2]*u[1]);
1098:       uik[1] = -(diag[1]*u[0] + diag[3]*u[1]);
1099:       uik[2] = -(diag[0]*u[2] + diag[2]*u[3]);
1100:       uik[3] = -(diag[1]*u[2] + diag[3]*u[3]);

1102:       /* update D(k) += -U(i,k)^T * U_bar(i,k): dk += uik*ba[ili] */
1103:       dk[0] += uik[0]*u[0] + uik[1]*u[1];
1104:       dk[1] += uik[2]*u[0] + uik[3]*u[1];
1105:       dk[2] += uik[0]*u[2] + uik[1]*u[3];
1106:       dk[3] += uik[2]*u[2] + uik[3]*u[3];

1108:       PetscLogFlops(16.0*2.0);

1110:       /* update -U(i,k): ba[ili] = uik */
1111:       PetscMemcpy(ba+ili*4,uik,4*sizeof(MatScalar));

1113:       /* add multiple of row i to k-th row ... */
1114:       jmin = ili + 1; jmax = bi[i+1];
1115:       if (jmin < jmax) {
1116:         for (j=jmin; j<jmax; j++) {
1117:           /* rtmp += -U(i,k)^T * U_bar(i,j): rtmp[bj[j]] += uik*ba[j]; */
1118:           rtmp_ptr     = rtmp + bj[j]*4;
1119:           u            = ba + j*4;
1120:           rtmp_ptr[0] += uik[0]*u[0] + uik[1]*u[1];
1121:           rtmp_ptr[1] += uik[2]*u[0] + uik[3]*u[1];
1122:           rtmp_ptr[2] += uik[0]*u[2] + uik[1]*u[3];
1123:           rtmp_ptr[3] += uik[2]*u[2] + uik[3]*u[3];
1124:         }
1125:         PetscLogFlops(16.0*(jmax-jmin));

1127:         /* ... add i to row list for next nonzero entry */
1128:         il[i] = jmin;             /* update il(i) in column k+1, ... mbs-1 */
1129:         j     = bj[jmin];
1130:         jl[i] = jl[j]; jl[j] = i; /* update jl */
1131:       }
1132:       i = nexti;
1133:     }

1135:     /* save nonzero entries in k-th row of U ... */

1137:     /* invert diagonal block */
1138:     diag = ba+k*4;
1139:     PetscMemcpy(diag,dk,4*sizeof(MatScalar));
1140:     PetscKernel_A_gets_inverse_A_2(diag,shift);

1142:     jmin = bi[k]; jmax = bi[k+1];
1143:     if (jmin < jmax) {
1144:       for (j=jmin; j<jmax; j++) {
1145:         vj       = bj[j];      /* block col. index of U */
1146:         u        = ba + j*4;
1147:         rtmp_ptr = rtmp + vj*4;
1148:         for (k1=0; k1<4; k1++) {
1149:           *u++        = *rtmp_ptr;
1150:           *rtmp_ptr++ = 0.0;
1151:         }
1152:       }

1154:       /* ... add k to row list for first nonzero entry in k-th row */
1155:       il[k] = jmin;
1156:       i     = bj[jmin];
1157:       jl[k] = jl[i]; jl[i] = k;
1158:     }
1159:   }

1161:   PetscFree(rtmp);
1162:   PetscFree2(il,jl);

1164:   C->ops->solve          = MatSolve_SeqSBAIJ_2_NaturalOrdering_inplace;
1165:   C->ops->solvetranspose = MatSolve_SeqSBAIJ_2_NaturalOrdering_inplace;
1166:   C->ops->forwardsolve   = MatForwardSolve_SeqSBAIJ_2_NaturalOrdering_inplace;
1167:   C->ops->backwardsolve  = MatBackwardSolve_SeqSBAIJ_2_NaturalOrdering_inplace;
1168:   C->assembled           = PETSC_TRUE;
1169:   C->preallocated        = PETSC_TRUE;

1171:   PetscLogFlops(1.3333*8*b->mbs); /* from inverting diagonal blocks */
1172:   return(0);
1173: }

1175: /*
1176:     Numeric U^T*D*U factorization for SBAIJ format.
1177:     Version for blocks are 1 by 1.
1178: */
1181: PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_1_inplace(Mat C,Mat A,const MatFactorInfo *info)
1182: {
1183:   Mat_SeqSBAIJ   *a=(Mat_SeqSBAIJ*)A->data,*b=(Mat_SeqSBAIJ*)C->data;
1184:   IS             ip=b->row;
1186:   const PetscInt *ai,*aj,*rip;
1187:   PetscInt       *a2anew,i,j,mbs=a->mbs,*bi=b->i,*bj=b->j,*bcol;
1188:   PetscInt       k,jmin,jmax,*jl,*il,col,nexti,ili,nz;
1189:   MatScalar      *rtmp,*ba=b->a,*bval,*aa,dk,uikdi;
1190:   PetscReal      rs;
1191:   FactorShiftCtx sctx;

1194:   /* MatPivotSetUp(): initialize shift context sctx */
1195:   PetscMemzero(&sctx,sizeof(FactorShiftCtx));

1197:   ISGetIndices(ip,&rip);
1198:   if (!a->permute) {
1199:     ai = a->i; aj = a->j; aa = a->a;
1200:   } else {
1201:     ai     = a->inew; aj = a->jnew;
1202:     nz     = ai[mbs];
1203:     PetscMalloc1(nz,&aa);
1204:     a2anew = a->a2anew;
1205:     bval   = a->a;
1206:     for (j=0; j<nz; j++) {
1207:       aa[a2anew[j]] = *(bval++);
1208:     }
1209:   }

1211:   /* initialization */
1212:   /* il and jl record the first nonzero element in each row of the accessing
1213:      window U(0:k, k:mbs-1).
1214:      jl:    list of rows to be added to uneliminated rows
1215:             i>= k: jl(i) is the first row to be added to row i
1216:             i<  k: jl(i) is the row following row i in some list of rows
1217:             jl(i) = mbs indicates the end of a list
1218:      il(i): points to the first nonzero element in columns k,...,mbs-1 of
1219:             row i of U */
1220:   PetscMalloc3(mbs,&rtmp,mbs,&il,mbs,&jl);

1222:   do {
1223:     sctx.newshift = PETSC_FALSE;
1224:     for (i=0; i<mbs; i++) {
1225:       rtmp[i] = 0.0; jl[i] = mbs; il[0] = 0;
1226:     }

1228:     for (k = 0; k<mbs; k++) {
1229:       /*initialize k-th row by the perm[k]-th row of A */
1230:       jmin = ai[rip[k]]; jmax = ai[rip[k]+1];
1231:       bval = ba + bi[k];
1232:       for (j = jmin; j < jmax; j++) {
1233:         col       = rip[aj[j]];
1234:         rtmp[col] = aa[j];
1235:         *bval++   = 0.0; /* for in-place factorization */
1236:       }

1238:       /* shift the diagonal of the matrix */
1239:       if (sctx.nshift) rtmp[k] += sctx.shift_amount;

1241:       /* modify k-th row by adding in those rows i with U(i,k)!=0 */
1242:       dk = rtmp[k];
1243:       i  = jl[k]; /* first row to be added to k_th row  */

1245:       while (i < k) {
1246:         nexti = jl[i]; /* next row to be added to k_th row */

1248:         /* compute multiplier, update diag(k) and U(i,k) */
1249:         ili     = il[i]; /* index of first nonzero element in U(i,k:bms-1) */
1250:         uikdi   = -ba[ili]*ba[bi[i]]; /* diagonal(k) */
1251:         dk     += uikdi*ba[ili];
1252:         ba[ili] = uikdi; /* -U(i,k) */

1254:         /* add multiple of row i to k-th row */
1255:         jmin = ili + 1; jmax = bi[i+1];
1256:         if (jmin < jmax) {
1257:           for (j=jmin; j<jmax; j++) rtmp[bj[j]] += uikdi*ba[j];
1258:           PetscLogFlops(2.0*(jmax-jmin));

1260:           /* update il and jl for row i */
1261:           il[i] = jmin;
1262:           j     = bj[jmin]; jl[i] = jl[j]; jl[j] = i;
1263:         }
1264:         i = nexti;
1265:       }

1267:       /* shift the diagonals when zero pivot is detected */
1268:       /* compute rs=sum of abs(off-diagonal) */
1269:       rs   = 0.0;
1270:       jmin = bi[k]+1;
1271:       nz   = bi[k+1] - jmin;
1272:       if (nz) {
1273:         bcol = bj + jmin;
1274:         while (nz--) {
1275:           rs += PetscAbsScalar(rtmp[*bcol]);
1276:           bcol++;
1277:         }
1278:       }

1280:       sctx.rs = rs;
1281:       sctx.pv = dk;
1282:       MatPivotCheck(A,info,&sctx,k);
1283:       if (sctx.newshift) break;    /* sctx.shift_amount is updated */
1284:       dk = sctx.pv;

1286:       /* copy data into U(k,:) */
1287:       ba[bi[k]] = 1.0/dk; /* U(k,k) */
1288:       jmin      = bi[k]+1; jmax = bi[k+1];
1289:       if (jmin < jmax) {
1290:         for (j=jmin; j<jmax; j++) {
1291:           col = bj[j]; ba[j] = rtmp[col]; rtmp[col] = 0.0;
1292:         }
1293:         /* add the k-th row into il and jl */
1294:         il[k] = jmin;
1295:         i     = bj[jmin]; jl[k] = jl[i]; jl[i] = k;
1296:       }
1297:     }
1298:   } while (sctx.newshift);
1299:   PetscFree3(rtmp,il,jl);
1300:   if (a->permute) {PetscFree(aa);}

1302:   ISRestoreIndices(ip,&rip);

1304:   C->ops->solve          = MatSolve_SeqSBAIJ_1_inplace;
1305:   C->ops->solves         = MatSolves_SeqSBAIJ_1_inplace;
1306:   C->ops->solvetranspose = MatSolve_SeqSBAIJ_1_inplace;
1307:   C->ops->forwardsolve   = MatForwardSolve_SeqSBAIJ_1_inplace;
1308:   C->ops->backwardsolve  = MatBackwardSolve_SeqSBAIJ_1_inplace;
1309:   C->assembled           = PETSC_TRUE;
1310:   C->preallocated        = PETSC_TRUE;

1312:   PetscLogFlops(C->rmap->N);
1313:   if (sctx.nshift) {
1314:     if (info->shifttype == (PetscReal)MAT_SHIFT_NONZERO) {
1315:       PetscInfo2(A,"number of shiftnz tries %D, shift_amount %g\n",sctx.nshift,(double)sctx.shift_amount);
1316:     } else if (info->shifttype == (PetscReal)MAT_SHIFT_POSITIVE_DEFINITE) {
1317:       PetscInfo2(A,"number of shiftpd tries %D, shift_amount %g\n",sctx.nshift,(double)sctx.shift_amount);
1318:     }
1319:   }
1320:   return(0);
1321: }

1323: /*
1324:   Version for when blocks are 1 by 1 Using natural ordering under new datastructure
1325:   Modified from MatCholeskyFactorNumeric_SeqAIJ()
1326: */
1329: PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_1_NaturalOrdering(Mat B,Mat A,const MatFactorInfo *info)
1330: {
1331:   Mat_SeqSBAIJ   *a=(Mat_SeqSBAIJ*)A->data;
1332:   Mat_SeqSBAIJ   *b=(Mat_SeqSBAIJ*)B->data;
1334:   PetscInt       i,j,mbs=A->rmap->n,*bi=b->i,*bj=b->j,*bdiag=b->diag,*bjtmp;
1335:   PetscInt       *ai=a->i,*aj=a->j,*ajtmp;
1336:   PetscInt       k,jmin,jmax,*c2r,*il,col,nexti,ili,nz;
1337:   MatScalar      *rtmp,*ba=b->a,*bval,*aa=a->a,dk,uikdi;
1338:   FactorShiftCtx sctx;
1339:   PetscReal      rs;
1340:   MatScalar      d,*v;

1343:   PetscMalloc3(mbs,&rtmp,mbs,&il,mbs,&c2r);

1345:   /* MatPivotSetUp(): initialize shift context sctx */
1346:   PetscMemzero(&sctx,sizeof(FactorShiftCtx));

1348:   if (info->shifttype == (PetscReal)MAT_SHIFT_POSITIVE_DEFINITE) { /* set sctx.shift_top=max{rs} */
1349:     sctx.shift_top = info->zeropivot;

1351:     PetscMemzero(rtmp,mbs*sizeof(MatScalar));

1353:     for (i=0; i<mbs; i++) {
1354:       /* calculate sum(|aij|)-RealPart(aii), amt of shift needed for this row */
1355:       d        = (aa)[a->diag[i]];
1356:       rtmp[i] += -PetscRealPart(d);  /* diagonal entry */
1357:       ajtmp    = aj + ai[i] + 1;     /* exclude diagonal */
1358:       v        = aa + ai[i] + 1;
1359:       nz       = ai[i+1] - ai[i] - 1;
1360:       for (j=0; j<nz; j++) {
1361:         rtmp[i]        += PetscAbsScalar(v[j]);
1362:         rtmp[ajtmp[j]] += PetscAbsScalar(v[j]);
1363:       }
1364:       if (PetscRealPart(rtmp[i]) > sctx.shift_top) sctx.shift_top = PetscRealPart(rtmp[i]);
1365:     }
1366:     sctx.shift_top *= 1.1;
1367:     sctx.nshift_max = 5;
1368:     sctx.shift_lo   = 0.;
1369:     sctx.shift_hi   = 1.;
1370:   }

1372:   /* allocate working arrays
1373:      c2r: linked list, keep track of pivot rows for a given column. c2r[col]: head of the list for a given col
1374:      il:  for active k row, il[i] gives the index of the 1st nonzero entry in U[i,k:n-1] in bj and ba arrays
1375:   */
1376:   do {
1377:     sctx.newshift = PETSC_FALSE;

1379:     for (i=0; i<mbs; i++) c2r[i] = mbs;
1380:     if (mbs) il[0] = 0;

1382:     for (k = 0; k<mbs; k++) {
1383:       /* zero rtmp */
1384:       nz    = bi[k+1] - bi[k];
1385:       bjtmp = bj + bi[k];
1386:       for (j=0; j<nz; j++) rtmp[bjtmp[j]] = 0.0;

1388:       /* load in initial unfactored row */
1389:       bval = ba + bi[k];
1390:       jmin = ai[k]; jmax = ai[k+1];
1391:       for (j = jmin; j < jmax; j++) {
1392:         col       = aj[j];
1393:         rtmp[col] = aa[j];
1394:         *bval++   = 0.0; /* for in-place factorization */
1395:       }
1396:       /* shift the diagonal of the matrix: ZeropivotApply() */
1397:       rtmp[k] += sctx.shift_amount;  /* shift the diagonal of the matrix */

1399:       /* modify k-th row by adding in those rows i with U(i,k)!=0 */
1400:       dk = rtmp[k];
1401:       i  = c2r[k]; /* first row to be added to k_th row  */

1403:       while (i < k) {
1404:         nexti = c2r[i]; /* next row to be added to k_th row */

1406:         /* compute multiplier, update diag(k) and U(i,k) */
1407:         ili     = il[i]; /* index of first nonzero element in U(i,k:bms-1) */
1408:         uikdi   = -ba[ili]*ba[bdiag[i]]; /* diagonal(k) */
1409:         dk     += uikdi*ba[ili]; /* update diag[k] */
1410:         ba[ili] = uikdi; /* -U(i,k) */

1412:         /* add multiple of row i to k-th row */
1413:         jmin = ili + 1; jmax = bi[i+1];
1414:         if (jmin < jmax) {
1415:           for (j=jmin; j<jmax; j++) rtmp[bj[j]] += uikdi*ba[j];
1416:           /* update il and c2r for row i */
1417:           il[i] = jmin;
1418:           j     = bj[jmin]; c2r[i] = c2r[j]; c2r[j] = i;
1419:         }
1420:         i = nexti;
1421:       }

1423:       /* copy data into U(k,:) */
1424:       rs   = 0.0;
1425:       jmin = bi[k]; jmax = bi[k+1]-1;
1426:       if (jmin < jmax) {
1427:         for (j=jmin; j<jmax; j++) {
1428:           col = bj[j]; ba[j] = rtmp[col]; rs += PetscAbsScalar(ba[j]);
1429:         }
1430:         /* add the k-th row into il and c2r */
1431:         il[k] = jmin;
1432:         i     = bj[jmin]; c2r[k] = c2r[i]; c2r[i] = k;
1433:       }

1435:       sctx.rs = rs;
1436:       sctx.pv = dk;
1437:       MatPivotCheck(A,info,&sctx,k);
1438:       if (sctx.newshift) break;
1439:       dk = sctx.pv;

1441:       ba[bdiag[k]] = 1.0/dk; /* U(k,k) */
1442:     }
1443:   } while (sctx.newshift);

1445:   PetscFree3(rtmp,il,c2r);

1447:   B->ops->solve          = MatSolve_SeqSBAIJ_1_NaturalOrdering;
1448:   B->ops->solves         = MatSolves_SeqSBAIJ_1;
1449:   B->ops->solvetranspose = MatSolve_SeqSBAIJ_1_NaturalOrdering;
1450:   B->ops->forwardsolve   = MatForwardSolve_SeqSBAIJ_1_NaturalOrdering;
1451:   B->ops->backwardsolve  = MatBackwardSolve_SeqSBAIJ_1_NaturalOrdering;

1453:   B->assembled    = PETSC_TRUE;
1454:   B->preallocated = PETSC_TRUE;

1456:   PetscLogFlops(B->rmap->n);

1458:   /* MatPivotView() */
1459:   if (sctx.nshift) {
1460:     if (info->shifttype == (PetscReal)MAT_SHIFT_POSITIVE_DEFINITE) {
1461:       PetscInfo4(A,"number of shift_pd tries %D, shift_amount %g, diagonal shifted up by %e fraction top_value %e\n",sctx.nshift,(double)sctx.shift_amount,(double)sctx.shift_fraction,(double)sctx.shift_top);
1462:     } else if (info->shifttype == (PetscReal)MAT_SHIFT_NONZERO) {
1463:       PetscInfo2(A,"number of shift_nz tries %D, shift_amount %g\n",sctx.nshift,(double)sctx.shift_amount);
1464:     } else if (info->shifttype == (PetscReal)MAT_SHIFT_INBLOCKS) {
1465:       PetscInfo2(A,"number of shift_inblocks applied %D, each shift_amount %g\n",sctx.nshift,(double)info->shiftamount);
1466:     }
1467:   }
1468:   return(0);
1469: }

1473: PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_1_NaturalOrdering_inplace(Mat C,Mat A,const MatFactorInfo *info)
1474: {
1475:   Mat_SeqSBAIJ   *a=(Mat_SeqSBAIJ*)A->data,*b=(Mat_SeqSBAIJ*)C->data;
1477:   PetscInt       i,j,mbs = a->mbs;
1478:   PetscInt       *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j;
1479:   PetscInt       k,jmin,*jl,*il,nexti,ili,*acol,*bcol,nz;
1480:   MatScalar      *rtmp,*ba=b->a,*aa=a->a,dk,uikdi,*aval,*bval;
1481:   PetscReal      rs;
1482:   FactorShiftCtx sctx;

1485:   /* MatPivotSetUp(): initialize shift context sctx */
1486:   PetscMemzero(&sctx,sizeof(FactorShiftCtx));

1488:   /* initialization */
1489:   /* il and jl record the first nonzero element in each row of the accessing
1490:      window U(0:k, k:mbs-1).
1491:      jl:    list of rows to be added to uneliminated rows
1492:             i>= k: jl(i) is the first row to be added to row i
1493:             i<  k: jl(i) is the row following row i in some list of rows
1494:             jl(i) = mbs indicates the end of a list
1495:      il(i): points to the first nonzero element in U(i,k:mbs-1)
1496:   */
1497:   PetscMalloc1(mbs,&rtmp);
1498:   PetscMalloc2(mbs,&il,mbs,&jl);

1500:   do {
1501:     sctx.newshift = PETSC_FALSE;
1502:     for (i=0; i<mbs; i++) {
1503:       rtmp[i] = 0.0; jl[i] = mbs; il[0] = 0;
1504:     }

1506:     for (k = 0; k<mbs; k++) {
1507:       /*initialize k-th row with elements nonzero in row perm(k) of A */
1508:       nz   = ai[k+1] - ai[k];
1509:       acol = aj + ai[k];
1510:       aval = aa + ai[k];
1511:       bval = ba + bi[k];
1512:       while (nz--) {
1513:         rtmp[*acol++] = *aval++;
1514:         *bval++       = 0.0; /* for in-place factorization */
1515:       }

1517:       /* shift the diagonal of the matrix */
1518:       if (sctx.nshift) rtmp[k] += sctx.shift_amount;

1520:       /* modify k-th row by adding in those rows i with U(i,k)!=0 */
1521:       dk = rtmp[k];
1522:       i  = jl[k]; /* first row to be added to k_th row  */

1524:       while (i < k) {
1525:         nexti = jl[i]; /* next row to be added to k_th row */
1526:         /* compute multiplier, update D(k) and U(i,k) */
1527:         ili     = il[i]; /* index of first nonzero element in U(i,k:bms-1) */
1528:         uikdi   = -ba[ili]*ba[bi[i]];
1529:         dk     += uikdi*ba[ili];
1530:         ba[ili] = uikdi; /* -U(i,k) */

1532:         /* add multiple of row i to k-th row ... */
1533:         jmin = ili + 1;
1534:         nz   = bi[i+1] - jmin;
1535:         if (nz > 0) {
1536:           bcol = bj + jmin;
1537:           bval = ba + jmin;
1538:           PetscLogFlops(2.0*nz);
1539:           while (nz--) rtmp[*bcol++] += uikdi*(*bval++);

1541:           /* update il and jl for i-th row */
1542:           il[i] = jmin;
1543:           j     = bj[jmin]; jl[i] = jl[j]; jl[j] = i;
1544:         }
1545:         i = nexti;
1546:       }

1548:       /* shift the diagonals when zero pivot is detected */
1549:       /* compute rs=sum of abs(off-diagonal) */
1550:       rs   = 0.0;
1551:       jmin = bi[k]+1;
1552:       nz   = bi[k+1] - jmin;
1553:       if (nz) {
1554:         bcol = bj + jmin;
1555:         while (nz--) {
1556:           rs += PetscAbsScalar(rtmp[*bcol]);
1557:           bcol++;
1558:         }
1559:       }

1561:       sctx.rs = rs;
1562:       sctx.pv = dk;
1563:       MatPivotCheck(A,info,&sctx,k);
1564:       if (sctx.newshift) break;    /* sctx.shift_amount is updated */
1565:       dk = sctx.pv;

1567:       /* copy data into U(k,:) */
1568:       ba[bi[k]] = 1.0/dk;
1569:       jmin      = bi[k]+1;
1570:       nz        = bi[k+1] - jmin;
1571:       if (nz) {
1572:         bcol = bj + jmin;
1573:         bval = ba + jmin;
1574:         while (nz--) {
1575:           *bval++       = rtmp[*bcol];
1576:           rtmp[*bcol++] = 0.0;
1577:         }
1578:         /* add k-th row into il and jl */
1579:         il[k] = jmin;
1580:         i     = bj[jmin]; jl[k] = jl[i]; jl[i] = k;
1581:       }
1582:     } /* end of for (k = 0; k<mbs; k++) */
1583:   } while (sctx.newshift);
1584:   PetscFree(rtmp);
1585:   PetscFree2(il,jl);

1587:   C->ops->solve          = MatSolve_SeqSBAIJ_1_NaturalOrdering_inplace;
1588:   C->ops->solves         = MatSolves_SeqSBAIJ_1_inplace;
1589:   C->ops->solvetranspose = MatSolve_SeqSBAIJ_1_NaturalOrdering_inplace;
1590:   C->ops->forwardsolve   = MatForwardSolve_SeqSBAIJ_1_NaturalOrdering_inplace;
1591:   C->ops->backwardsolve  = MatBackwardSolve_SeqSBAIJ_1_NaturalOrdering_inplace;

1593:   C->assembled    = PETSC_TRUE;
1594:   C->preallocated = PETSC_TRUE;

1596:   PetscLogFlops(C->rmap->N);
1597:   if (sctx.nshift) {
1598:     if (info->shifttype == (PetscReal)MAT_SHIFT_NONZERO) {
1599:       PetscInfo2(A,"number of shiftnz tries %D, shift_amount %g\n",sctx.nshift,(double)sctx.shift_amount);
1600:     } else if (info->shifttype == (PetscReal)MAT_SHIFT_POSITIVE_DEFINITE) {
1601:       PetscInfo2(A,"number of shiftpd tries %D, shift_amount %g\n",sctx.nshift,(double)sctx.shift_amount);
1602:     }
1603:   }
1604:   return(0);
1605: }

1609: PetscErrorCode MatCholeskyFactor_SeqSBAIJ(Mat A,IS perm,const MatFactorInfo *info)
1610: {
1612:   Mat            C;

1615:   MatGetFactor(A,"petsc",MAT_FACTOR_CHOLESKY,&C);
1616:   MatCholeskyFactorSymbolic(C,A,perm,info);
1617:   MatCholeskyFactorNumeric(C,A,info);

1619:   A->ops->solve          = C->ops->solve;
1620:   A->ops->solvetranspose = C->ops->solvetranspose;

1622:   MatHeaderMerge(A,C);
1623:   return(0);
1624: }