Actual source code: aijfact.c

petsc-3.11.4 2019-09-28
Report Typos and Errors

  2:  #include <../src/mat/impls/aij/seq/aij.h>
  3:  #include <../src/mat/impls/sbaij/seq/sbaij.h>
  4:  #include <petscbt.h>
  5:  #include <../src/mat/utils/freespace.h>

  7: /*
  8:       Computes an ordering to get most of the large numerical values in the lower triangular part of the matrix

 10:       This code does not work and is not called anywhere. It would be registered with MatOrderingRegisterAll()
 11: */
 12: PetscErrorCode MatGetOrdering_Flow_SeqAIJ(Mat mat,MatOrderingType type,IS *irow,IS *icol)
 13: {
 14:   Mat_SeqAIJ        *a = (Mat_SeqAIJ*)mat->data;
 15:   PetscErrorCode    ierr;
 16:   PetscInt          i,j,jj,k, kk,n = mat->rmap->n, current = 0, newcurrent = 0,*order;
 17:   const PetscInt    *ai = a->i, *aj = a->j;
 18:   const PetscScalar *aa = a->a;
 19:   PetscBool         *done;
 20:   PetscReal         best,past = 0,future;

 23:   /* pick initial row */
 24:   best = -1;
 25:   for (i=0; i<n; i++) {
 26:     future = 0.0;
 27:     for (j=ai[i]; j<ai[i+1]; j++) {
 28:       if (aj[j] != i) future += PetscAbsScalar(aa[j]);
 29:       else              past  = PetscAbsScalar(aa[j]);
 30:     }
 31:     if (!future) future = 1.e-10; /* if there is zero in the upper diagonal part want to rank this row high */
 32:     if (past/future > best) {
 33:       best    = past/future;
 34:       current = i;
 35:     }
 36:   }

 38:   PetscMalloc1(n,&done);
 39:   PetscMemzero(done,n*sizeof(PetscBool));
 40:   PetscMalloc1(n,&order);
 41:   order[0] = current;
 42:   for (i=0; i<n-1; i++) {
 43:     done[current] = PETSC_TRUE;
 44:     best          = -1;
 45:     /* loop over all neighbors of current pivot */
 46:     for (j=ai[current]; j<ai[current+1]; j++) {
 47:       jj = aj[j];
 48:       if (done[jj]) continue;
 49:       /* loop over columns of potential next row computing weights for below and above diagonal */
 50:       past = future = 0.0;
 51:       for (k=ai[jj]; k<ai[jj+1]; k++) {
 52:         kk = aj[k];
 53:         if (done[kk]) past += PetscAbsScalar(aa[k]);
 54:         else if (kk != jj) future += PetscAbsScalar(aa[k]);
 55:       }
 56:       if (!future) future = 1.e-10; /* if there is zero in the upper diagonal part want to rank this row high */
 57:       if (past/future > best) {
 58:         best       = past/future;
 59:         newcurrent = jj;
 60:       }
 61:     }
 62:     if (best == -1) { /* no neighbors to select from so select best of all that remain */
 63:       best = -1;
 64:       for (k=0; k<n; k++) {
 65:         if (done[k]) continue;
 66:         future = 0.0;
 67:         past   = 0.0;
 68:         for (j=ai[k]; j<ai[k+1]; j++) {
 69:           kk = aj[j];
 70:           if (done[kk])       past += PetscAbsScalar(aa[j]);
 71:           else if (kk != k) future += PetscAbsScalar(aa[j]);
 72:         }
 73:         if (!future) future = 1.e-10; /* if there is zero in the upper diagonal part want to rank this row high */
 74:         if (past/future > best) {
 75:           best       = past/future;
 76:           newcurrent = k;
 77:         }
 78:       }
 79:     }
 80:     if (current == newcurrent) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"newcurrent cannot be current");
 81:     current    = newcurrent;
 82:     order[i+1] = current;
 83:   }
 84:   ISCreateGeneral(PETSC_COMM_SELF,n,order,PETSC_COPY_VALUES,irow);
 85:   *icol = *irow;
 86:   PetscObjectReference((PetscObject)*irow);
 87:   PetscFree(done);
 88:   PetscFree(order);
 89:   return(0);
 90: }

 92: PETSC_INTERN PetscErrorCode MatGetFactor_seqaij_petsc(Mat A,MatFactorType ftype,Mat *B)
 93: {
 94:   PetscInt       n = A->rmap->n;

 98: #if defined(PETSC_USE_COMPLEX)
 99:   if (A->hermitian && (ftype == MAT_FACTOR_CHOLESKY || ftype == MAT_FACTOR_ICC)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Hermitian Factor is not supported");
100: #endif
101:   MatCreate(PetscObjectComm((PetscObject)A),B);
102:   MatSetSizes(*B,n,n,n,n);
103:   if (ftype == MAT_FACTOR_LU || ftype == MAT_FACTOR_ILU || ftype == MAT_FACTOR_ILUDT) {
104:     MatSetType(*B,MATSEQAIJ);

106:     (*B)->ops->ilufactorsymbolic = MatILUFactorSymbolic_SeqAIJ;
107:     (*B)->ops->lufactorsymbolic  = MatLUFactorSymbolic_SeqAIJ;

109:     MatSetBlockSizesFromMats(*B,A,A);
110:   } else if (ftype == MAT_FACTOR_CHOLESKY || ftype == MAT_FACTOR_ICC) {
111:     MatSetType(*B,MATSEQSBAIJ);
112:     MatSeqSBAIJSetPreallocation(*B,1,MAT_SKIP_ALLOCATION,NULL);

114:     (*B)->ops->iccfactorsymbolic      = MatICCFactorSymbolic_SeqAIJ;
115:     (*B)->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_SeqAIJ;
116:   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Factor type not supported");
117:   (*B)->factortype = ftype;

119:   PetscFree((*B)->solvertype);
120:   PetscStrallocpy(MATSOLVERPETSC,&(*B)->solvertype);
121:   return(0);
122: }

124: PetscErrorCode MatLUFactorSymbolic_SeqAIJ_inplace(Mat B,Mat A,IS isrow,IS iscol,const MatFactorInfo *info)
125: {
126:   Mat_SeqAIJ         *a = (Mat_SeqAIJ*)A->data,*b;
127:   IS                 isicol;
128:   PetscErrorCode     ierr;
129:   const PetscInt     *r,*ic;
130:   PetscInt           i,n=A->rmap->n,*ai=a->i,*aj=a->j;
131:   PetscInt           *bi,*bj,*ajtmp;
132:   PetscInt           *bdiag,row,nnz,nzi,reallocs=0,nzbd,*im;
133:   PetscReal          f;
134:   PetscInt           nlnk,*lnk,k,**bi_ptr;
135:   PetscFreeSpaceList free_space=NULL,current_space=NULL;
136:   PetscBT            lnkbt;
137:   PetscBool          missing;

140:   if (A->rmap->N != A->cmap->N) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"matrix must be square");
141:   MatMissingDiagonal(A,&missing,&i);
142:   if (missing) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Matrix is missing diagonal entry %D",i);

144:   ISInvertPermutation(iscol,PETSC_DECIDE,&isicol);
145:   ISGetIndices(isrow,&r);
146:   ISGetIndices(isicol,&ic);

148:   /* get new row pointers */
149:   PetscMalloc1(n+1,&bi);
150:   bi[0] = 0;

152:   /* bdiag is location of diagonal in factor */
153:   PetscMalloc1(n+1,&bdiag);
154:   bdiag[0] = 0;

156:   /* linked list for storing column indices of the active row */
157:   nlnk = n + 1;
158:   PetscLLCreate(n,n,nlnk,lnk,lnkbt);

160:   PetscMalloc2(n+1,&bi_ptr,n+1,&im);

162:   /* initial FreeSpace size is f*(ai[n]+1) */
163:   f             = info->fill;
164:   PetscFreeSpaceGet(PetscRealIntMultTruncate(f,ai[n]+1),&free_space);
165:   current_space = free_space;

167:   for (i=0; i<n; i++) {
168:     /* copy previous fill into linked list */
169:     nzi = 0;
170:     nnz = ai[r[i]+1] - ai[r[i]];
171:     ajtmp = aj + ai[r[i]];
172:     PetscLLAddPerm(nnz,ajtmp,ic,n,nlnk,lnk,lnkbt);
173:     nzi  += nlnk;

175:     /* add pivot rows into linked list */
176:     row = lnk[n];
177:     while (row < i) {
178:       nzbd  = bdiag[row] - bi[row] + 1;   /* num of entries in the row with column index <= row */
179:       ajtmp = bi_ptr[row] + nzbd;   /* points to the entry next to the diagonal */
180:       PetscLLAddSortedLU(ajtmp,row,nlnk,lnk,lnkbt,i,nzbd,im);
181:       nzi  += nlnk;
182:       row   = lnk[row];
183:     }
184:     bi[i+1] = bi[i] + nzi;
185:     im[i]   = nzi;

187:     /* mark bdiag */
188:     nzbd = 0;
189:     nnz  = nzi;
190:     k    = lnk[n];
191:     while (nnz-- && k < i) {
192:       nzbd++;
193:       k = lnk[k];
194:     }
195:     bdiag[i] = bi[i] + nzbd;

197:     /* if free space is not available, make more free space */
198:     if (current_space->local_remaining<nzi) {
199:       nnz  = PetscIntMultTruncate(n - i,nzi); /* estimated and max additional space needed */
200:       PetscFreeSpaceGet(nnz,&current_space);
201:       reallocs++;
202:     }

204:     /* copy data into free space, then initialize lnk */
205:     PetscLLClean(n,n,nzi,lnk,current_space->array,lnkbt);

207:     bi_ptr[i]                       = current_space->array;
208:     current_space->array           += nzi;
209:     current_space->local_used      += nzi;
210:     current_space->local_remaining -= nzi;
211:   }
212: #if defined(PETSC_USE_INFO)
213:   if (ai[n] != 0) {
214:     PetscReal af = ((PetscReal)bi[n])/((PetscReal)ai[n]);
215:     PetscInfo3(A,"Reallocs %D Fill ratio:given %g needed %g\n",reallocs,(double)f,(double)af);
216:     PetscInfo1(A,"Run with -pc_factor_fill %g or use \n",(double)af);
217:     PetscInfo1(A,"PCFactorSetFill(pc,%g);\n",(double)af);
218:     PetscInfo(A,"for best performance.\n");
219:   } else {
220:     PetscInfo(A,"Empty matrix\n");
221:   }
222: #endif

224:   ISRestoreIndices(isrow,&r);
225:   ISRestoreIndices(isicol,&ic);

227:   /* destroy list of free space and other temporary array(s) */
228:   PetscMalloc1(bi[n]+1,&bj);
229:   PetscFreeSpaceContiguous(&free_space,bj);
230:   PetscLLDestroy(lnk,lnkbt);
231:   PetscFree2(bi_ptr,im);

233:   /* put together the new matrix */
234:   MatSeqAIJSetPreallocation_SeqAIJ(B,MAT_SKIP_ALLOCATION,NULL);
235:   PetscLogObjectParent((PetscObject)B,(PetscObject)isicol);
236:   b    = (Mat_SeqAIJ*)(B)->data;

238:   b->free_a       = PETSC_TRUE;
239:   b->free_ij      = PETSC_TRUE;
240:   b->singlemalloc = PETSC_FALSE;

242:   PetscMalloc1(bi[n]+1,&b->a);
243:   b->j    = bj;
244:   b->i    = bi;
245:   b->diag = bdiag;
246:   b->ilen = 0;
247:   b->imax = 0;
248:   b->row  = isrow;
249:   b->col  = iscol;
250:   PetscObjectReference((PetscObject)isrow);
251:   PetscObjectReference((PetscObject)iscol);
252:   b->icol = isicol;
253:   PetscMalloc1(n+1,&b->solve_work);

255:   /* In b structure:  Free imax, ilen, old a, old j.  Allocate solve_work, new a, new j */
256:   PetscLogObjectMemory((PetscObject)B,(bi[n]-n)*(sizeof(PetscInt)+sizeof(PetscScalar)));
257:   b->maxnz = b->nz = bi[n];

259:   (B)->factortype            = MAT_FACTOR_LU;
260:   (B)->info.factor_mallocs   = reallocs;
261:   (B)->info.fill_ratio_given = f;

263:   if (ai[n]) {
264:     (B)->info.fill_ratio_needed = ((PetscReal)bi[n])/((PetscReal)ai[n]);
265:   } else {
266:     (B)->info.fill_ratio_needed = 0.0;
267:   }
268:   (B)->ops->lufactornumeric = MatLUFactorNumeric_SeqAIJ_inplace;
269:   if (a->inode.size) {
270:     (B)->ops->lufactornumeric = MatLUFactorNumeric_SeqAIJ_Inode_inplace;
271:   }
272:   return(0);
273: }

275: PetscErrorCode MatLUFactorSymbolic_SeqAIJ(Mat B,Mat A,IS isrow,IS iscol,const MatFactorInfo *info)
276: {
277:   Mat_SeqAIJ         *a = (Mat_SeqAIJ*)A->data,*b;
278:   IS                 isicol;
279:   PetscErrorCode     ierr;
280:   const PetscInt     *r,*ic,*ai=a->i,*aj=a->j,*ajtmp;
281:   PetscInt           i,n=A->rmap->n;
282:   PetscInt           *bi,*bj;
283:   PetscInt           *bdiag,row,nnz,nzi,reallocs=0,nzbd,*im;
284:   PetscReal          f;
285:   PetscInt           nlnk,*lnk,k,**bi_ptr;
286:   PetscFreeSpaceList free_space=NULL,current_space=NULL;
287:   PetscBT            lnkbt;
288:   PetscBool          missing;

291:   if (A->rmap->N != A->cmap->N) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"matrix must be square");
292:   MatMissingDiagonal(A,&missing,&i);
293:   if (missing) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Matrix is missing diagonal entry %D",i);
294: 
295:   ISInvertPermutation(iscol,PETSC_DECIDE,&isicol);
296:   ISGetIndices(isrow,&r);
297:   ISGetIndices(isicol,&ic);

299:   /* get new row and diagonal pointers, must be allocated separately because they will be given to the Mat_SeqAIJ and freed separately */
300:   PetscMalloc1(n+1,&bi);
301:   PetscMalloc1(n+1,&bdiag);
302:   bi[0] = bdiag[0] = 0;

304:   /* linked list for storing column indices of the active row */
305:   nlnk = n + 1;
306:   PetscLLCreate(n,n,nlnk,lnk,lnkbt);

308:   PetscMalloc2(n+1,&bi_ptr,n+1,&im);

310:   /* initial FreeSpace size is f*(ai[n]+1) */
311:   f             = info->fill;
312:   PetscFreeSpaceGet(PetscRealIntMultTruncate(f,ai[n]+1),&free_space);
313:   current_space = free_space;

315:   for (i=0; i<n; i++) {
316:     /* copy previous fill into linked list */
317:     nzi = 0;
318:     nnz = ai[r[i]+1] - ai[r[i]];
319:     ajtmp = aj + ai[r[i]];
320:     PetscLLAddPerm(nnz,ajtmp,ic,n,nlnk,lnk,lnkbt);
321:     nzi  += nlnk;

323:     /* add pivot rows into linked list */
324:     row = lnk[n];
325:     while (row < i) {
326:       nzbd  = bdiag[row] + 1; /* num of entries in the row with column index <= row */
327:       ajtmp = bi_ptr[row] + nzbd; /* points to the entry next to the diagonal */
328:       PetscLLAddSortedLU(ajtmp,row,nlnk,lnk,lnkbt,i,nzbd,im);
329:       nzi  += nlnk;
330:       row   = lnk[row];
331:     }
332:     bi[i+1] = bi[i] + nzi;
333:     im[i]   = nzi;

335:     /* mark bdiag */
336:     nzbd = 0;
337:     nnz  = nzi;
338:     k    = lnk[n];
339:     while (nnz-- && k < i) {
340:       nzbd++;
341:       k = lnk[k];
342:     }
343:     bdiag[i] = nzbd; /* note: bdiag[i] = nnzL as input for PetscFreeSpaceContiguous_LU() */

345:     /* if free space is not available, make more free space */
346:     if (current_space->local_remaining<nzi) {
347:       /* estimated additional space needed */
348:       nnz  = PetscIntMultTruncate(2,PetscIntMultTruncate(n-1,nzi));
349:       PetscFreeSpaceGet(nnz,&current_space);
350:       reallocs++;
351:     }

353:     /* copy data into free space, then initialize lnk */
354:     PetscLLClean(n,n,nzi,lnk,current_space->array,lnkbt);

356:     bi_ptr[i]                       = current_space->array;
357:     current_space->array           += nzi;
358:     current_space->local_used      += nzi;
359:     current_space->local_remaining -= nzi;
360:   }

362:   ISRestoreIndices(isrow,&r);
363:   ISRestoreIndices(isicol,&ic);

365:   /*   copy free_space into bj and free free_space; set bi, bj, bdiag in new datastructure; */
366:   PetscMalloc1(bi[n]+1,&bj);
367:   PetscFreeSpaceContiguous_LU(&free_space,bj,n,bi,bdiag);
368:   PetscLLDestroy(lnk,lnkbt);
369:   PetscFree2(bi_ptr,im);

371:   /* put together the new matrix */
372:   MatSeqAIJSetPreallocation_SeqAIJ(B,MAT_SKIP_ALLOCATION,NULL);
373:   PetscLogObjectParent((PetscObject)B,(PetscObject)isicol);
374:   b    = (Mat_SeqAIJ*)(B)->data;

376:   b->free_a       = PETSC_TRUE;
377:   b->free_ij      = PETSC_TRUE;
378:   b->singlemalloc = PETSC_FALSE;

380:   PetscMalloc1(bdiag[0]+1,&b->a);

382:   b->j    = bj;
383:   b->i    = bi;
384:   b->diag = bdiag;
385:   b->ilen = 0;
386:   b->imax = 0;
387:   b->row  = isrow;
388:   b->col  = iscol;
389:   PetscObjectReference((PetscObject)isrow);
390:   PetscObjectReference((PetscObject)iscol);
391:   b->icol = isicol;
392:   PetscMalloc1(n+1,&b->solve_work);

394:   /* In b structure:  Free imax, ilen, old a, old j.  Allocate solve_work, new a, new j */
395:   PetscLogObjectMemory((PetscObject)B,(bdiag[0]+1)*(sizeof(PetscInt)+sizeof(PetscScalar)));
396:   b->maxnz = b->nz = bdiag[0]+1;

398:   B->factortype            = MAT_FACTOR_LU;
399:   B->info.factor_mallocs   = reallocs;
400:   B->info.fill_ratio_given = f;

402:   if (ai[n]) {
403:     B->info.fill_ratio_needed = ((PetscReal)(bdiag[0]+1))/((PetscReal)ai[n]);
404:   } else {
405:     B->info.fill_ratio_needed = 0.0;
406:   }
407: #if defined(PETSC_USE_INFO)
408:   if (ai[n] != 0) {
409:     PetscReal af = B->info.fill_ratio_needed;
410:     PetscInfo3(A,"Reallocs %D Fill ratio:given %g needed %g\n",reallocs,(double)f,(double)af);
411:     PetscInfo1(A,"Run with -pc_factor_fill %g or use \n",(double)af);
412:     PetscInfo1(A,"PCFactorSetFill(pc,%g);\n",(double)af);
413:     PetscInfo(A,"for best performance.\n");
414:   } else {
415:     PetscInfo(A,"Empty matrix\n");
416:   }
417: #endif
418:   B->ops->lufactornumeric = MatLUFactorNumeric_SeqAIJ;
419:   if (a->inode.size) {
420:     B->ops->lufactornumeric = MatLUFactorNumeric_SeqAIJ_Inode;
421:   }
422:   MatSeqAIJCheckInode_FactorLU(B);
423:   return(0);
424: }

426: /*
427:     Trouble in factorization, should we dump the original matrix?
428: */
429: PetscErrorCode MatFactorDumpMatrix(Mat A)
430: {
432:   PetscBool      flg = PETSC_FALSE;

435:   PetscOptionsGetBool(((PetscObject)A)->options,NULL,"-mat_factor_dump_on_error",&flg,NULL);
436:   if (flg) {
437:     PetscViewer viewer;
438:     char        filename[PETSC_MAX_PATH_LEN];

440:     PetscSNPrintf(filename,PETSC_MAX_PATH_LEN,"matrix_factor_error.%d",PetscGlobalRank);
441:     PetscViewerBinaryOpen(PetscObjectComm((PetscObject)A),filename,FILE_MODE_WRITE,&viewer);
442:     MatView(A,viewer);
443:     PetscViewerDestroy(&viewer);
444:   }
445:   return(0);
446: }

448: PetscErrorCode MatLUFactorNumeric_SeqAIJ(Mat B,Mat A,const MatFactorInfo *info)
449: {
450:   Mat             C     =B;
451:   Mat_SeqAIJ      *a    =(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)C->data;
452:   IS              isrow = b->row,isicol = b->icol;
453:   PetscErrorCode  ierr;
454:   const PetscInt  *r,*ic,*ics;
455:   const PetscInt  n=A->rmap->n,*ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j,*bdiag=b->diag;
456:   PetscInt        i,j,k,nz,nzL,row,*pj;
457:   const PetscInt  *ajtmp,*bjtmp;
458:   MatScalar       *rtmp,*pc,multiplier,*pv;
459:   const MatScalar *aa=a->a,*v;
460:   PetscBool       row_identity,col_identity;
461:   FactorShiftCtx  sctx;
462:   const PetscInt  *ddiag;
463:   PetscReal       rs;
464:   MatScalar       d;

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

470:   if (info->shifttype == (PetscReal) MAT_SHIFT_POSITIVE_DEFINITE) { /* set sctx.shift_top=max{rs} */
471:     ddiag          = a->diag;
472:     sctx.shift_top = info->zeropivot;
473:     for (i=0; i<n; i++) {
474:       /* calculate sum(|aij|)-RealPart(aii), amt of shift needed for this row */
475:       d  = (aa)[ddiag[i]];
476:       rs = -PetscAbsScalar(d) - PetscRealPart(d);
477:       v  = aa+ai[i];
478:       nz = ai[i+1] - ai[i];
479:       for (j=0; j<nz; j++) rs += PetscAbsScalar(v[j]);
480:       if (rs>sctx.shift_top) sctx.shift_top = rs;
481:     }
482:     sctx.shift_top *= 1.1;
483:     sctx.nshift_max = 5;
484:     sctx.shift_lo   = 0.;
485:     sctx.shift_hi   = 1.;
486:   }

488:   ISGetIndices(isrow,&r);
489:   ISGetIndices(isicol,&ic);
490:   PetscMalloc1(n+1,&rtmp);
491:   ics  = ic;

493:   do {
494:     sctx.newshift = PETSC_FALSE;
495:     for (i=0; i<n; i++) {
496:       /* zero rtmp */
497:       /* L part */
498:       nz    = bi[i+1] - bi[i];
499:       bjtmp = bj + bi[i];
500:       for  (j=0; j<nz; j++) rtmp[bjtmp[j]] = 0.0;

502:       /* U part */
503:       nz    = bdiag[i]-bdiag[i+1];
504:       bjtmp = bj + bdiag[i+1]+1;
505:       for  (j=0; j<nz; j++) rtmp[bjtmp[j]] = 0.0;

507:       /* load in initial (unfactored row) */
508:       nz    = ai[r[i]+1] - ai[r[i]];
509:       ajtmp = aj + ai[r[i]];
510:       v     = aa + ai[r[i]];
511:       for (j=0; j<nz; j++) {
512:         rtmp[ics[ajtmp[j]]] = v[j];
513:       }
514:       /* ZeropivotApply() */
515:       rtmp[i] += sctx.shift_amount;  /* shift the diagonal of the matrix */

517:       /* elimination */
518:       bjtmp = bj + bi[i];
519:       row   = *bjtmp++;
520:       nzL   = bi[i+1] - bi[i];
521:       for (k=0; k < nzL; k++) {
522:         pc = rtmp + row;
523:         if (*pc != 0.0) {
524:           pv         = b->a + bdiag[row];
525:           multiplier = *pc * (*pv);
526:           *pc        = multiplier;

528:           pj = b->j + bdiag[row+1]+1; /* beginning of U(row,:) */
529:           pv = b->a + bdiag[row+1]+1;
530:           nz = bdiag[row]-bdiag[row+1]-1; /* num of entries in U(row,:) excluding diag */

532:           for (j=0; j<nz; j++) rtmp[pj[j]] -= multiplier * pv[j];
533:           PetscLogFlops(1+2*nz);
534:         }
535:         row = *bjtmp++;
536:       }

538:       /* finished row so stick it into b->a */
539:       rs = 0.0;
540:       /* L part */
541:       pv = b->a + bi[i];
542:       pj = b->j + bi[i];
543:       nz = bi[i+1] - bi[i];
544:       for (j=0; j<nz; j++) {
545:         pv[j] = rtmp[pj[j]]; rs += PetscAbsScalar(pv[j]);
546:       }

548:       /* U part */
549:       pv = b->a + bdiag[i+1]+1;
550:       pj = b->j + bdiag[i+1]+1;
551:       nz = bdiag[i] - bdiag[i+1]-1;
552:       for (j=0; j<nz; j++) {
553:         pv[j] = rtmp[pj[j]]; rs += PetscAbsScalar(pv[j]);
554:       }

556:       sctx.rs = rs;
557:       sctx.pv = rtmp[i];
558:       MatPivotCheck(B,A,info,&sctx,i);
559:       if (sctx.newshift) break; /* break for-loop */
560:       rtmp[i] = sctx.pv; /* sctx.pv might be updated in the case of MAT_SHIFT_INBLOCKS */

562:       /* Mark diagonal and invert diagonal for simplier triangular solves */
563:       pv  = b->a + bdiag[i];
564:       *pv = 1.0/rtmp[i];

566:     } /* endof for (i=0; i<n; i++) { */

568:     /* MatPivotRefine() */
569:     if (info->shifttype == (PetscReal)MAT_SHIFT_POSITIVE_DEFINITE && !sctx.newshift && sctx.shift_fraction>0 && sctx.nshift<sctx.nshift_max) {
570:       /*
571:        * if no shift in this attempt & shifting & started shifting & can refine,
572:        * then try lower shift
573:        */
574:       sctx.shift_hi       = sctx.shift_fraction;
575:       sctx.shift_fraction = (sctx.shift_hi+sctx.shift_lo)/2.;
576:       sctx.shift_amount   = sctx.shift_fraction * sctx.shift_top;
577:       sctx.newshift       = PETSC_TRUE;
578:       sctx.nshift++;
579:     }
580:   } while (sctx.newshift);

582:   PetscFree(rtmp);
583:   ISRestoreIndices(isicol,&ic);
584:   ISRestoreIndices(isrow,&r);

586:   ISIdentity(isrow,&row_identity);
587:   ISIdentity(isicol,&col_identity);
588:   if (b->inode.size) {
589:     C->ops->solve = MatSolve_SeqAIJ_Inode;
590:   } else if (row_identity && col_identity) {
591:     C->ops->solve = MatSolve_SeqAIJ_NaturalOrdering;
592:   } else {
593:     C->ops->solve = MatSolve_SeqAIJ;
594:   }
595:   C->ops->solveadd          = MatSolveAdd_SeqAIJ;
596:   C->ops->solvetranspose    = MatSolveTranspose_SeqAIJ;
597:   C->ops->solvetransposeadd = MatSolveTransposeAdd_SeqAIJ;
598:   C->ops->matsolve          = MatMatSolve_SeqAIJ;
599:   C->assembled              = PETSC_TRUE;
600:   C->preallocated           = PETSC_TRUE;

602:   PetscLogFlops(C->cmap->n);

604:   /* MatShiftView(A,info,&sctx) */
605:   if (sctx.nshift) {
606:     if (info->shifttype == (PetscReal)MAT_SHIFT_POSITIVE_DEFINITE) {
607:       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);
608:     } else if (info->shifttype == (PetscReal)MAT_SHIFT_NONZERO) {
609:       PetscInfo2(A,"number of shift_nz tries %D, shift_amount %g\n",sctx.nshift,(double)sctx.shift_amount);
610:     } else if (info->shifttype == (PetscReal)MAT_SHIFT_INBLOCKS) {
611:       PetscInfo2(A,"number of shift_inblocks applied %D, each shift_amount %g\n",sctx.nshift,(double)info->shiftamount);
612:     }
613:   }
614:   return(0);
615: }

617: PetscErrorCode MatLUFactorNumeric_SeqAIJ_inplace(Mat B,Mat A,const MatFactorInfo *info)
618: {
619:   Mat             C     =B;
620:   Mat_SeqAIJ      *a    =(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)C->data;
621:   IS              isrow = b->row,isicol = b->icol;
622:   PetscErrorCode  ierr;
623:   const PetscInt  *r,*ic,*ics;
624:   PetscInt        nz,row,i,j,n=A->rmap->n,diag;
625:   const PetscInt  *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j;
626:   const PetscInt  *ajtmp,*bjtmp,*diag_offset = b->diag,*pj;
627:   MatScalar       *pv,*rtmp,*pc,multiplier,d;
628:   const MatScalar *v,*aa=a->a;
629:   PetscReal       rs=0.0;
630:   FactorShiftCtx  sctx;
631:   const PetscInt  *ddiag;
632:   PetscBool       row_identity, col_identity;

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

638:   if (info->shifttype == (PetscReal) MAT_SHIFT_POSITIVE_DEFINITE) { /* set sctx.shift_top=max{rs} */
639:     ddiag          = a->diag;
640:     sctx.shift_top = info->zeropivot;
641:     for (i=0; i<n; i++) {
642:       /* calculate sum(|aij|)-RealPart(aii), amt of shift needed for this row */
643:       d  = (aa)[ddiag[i]];
644:       rs = -PetscAbsScalar(d) - PetscRealPart(d);
645:       v  = aa+ai[i];
646:       nz = ai[i+1] - ai[i];
647:       for (j=0; j<nz; j++) rs += PetscAbsScalar(v[j]);
648:       if (rs>sctx.shift_top) sctx.shift_top = rs;
649:     }
650:     sctx.shift_top *= 1.1;
651:     sctx.nshift_max = 5;
652:     sctx.shift_lo   = 0.;
653:     sctx.shift_hi   = 1.;
654:   }

656:   ISGetIndices(isrow,&r);
657:   ISGetIndices(isicol,&ic);
658:   PetscMalloc1(n+1,&rtmp);
659:   ics  = ic;

661:   do {
662:     sctx.newshift = PETSC_FALSE;
663:     for (i=0; i<n; i++) {
664:       nz    = bi[i+1] - bi[i];
665:       bjtmp = bj + bi[i];
666:       for  (j=0; j<nz; j++) rtmp[bjtmp[j]] = 0.0;

668:       /* load in initial (unfactored row) */
669:       nz    = ai[r[i]+1] - ai[r[i]];
670:       ajtmp = aj + ai[r[i]];
671:       v     = aa + ai[r[i]];
672:       for (j=0; j<nz; j++) {
673:         rtmp[ics[ajtmp[j]]] = v[j];
674:       }
675:       rtmp[ics[r[i]]] += sctx.shift_amount; /* shift the diagonal of the matrix */

677:       row = *bjtmp++;
678:       while  (row < i) {
679:         pc = rtmp + row;
680:         if (*pc != 0.0) {
681:           pv         = b->a + diag_offset[row];
682:           pj         = b->j + diag_offset[row] + 1;
683:           multiplier = *pc / *pv++;
684:           *pc        = multiplier;
685:           nz         = bi[row+1] - diag_offset[row] - 1;
686:           for (j=0; j<nz; j++) rtmp[pj[j]] -= multiplier * pv[j];
687:           PetscLogFlops(1+2*nz);
688:         }
689:         row = *bjtmp++;
690:       }
691:       /* finished row so stick it into b->a */
692:       pv   = b->a + bi[i];
693:       pj   = b->j + bi[i];
694:       nz   = bi[i+1] - bi[i];
695:       diag = diag_offset[i] - bi[i];
696:       rs   = 0.0;
697:       for (j=0; j<nz; j++) {
698:         pv[j] = rtmp[pj[j]];
699:         rs   += PetscAbsScalar(pv[j]);
700:       }
701:       rs -= PetscAbsScalar(pv[diag]);

703:       sctx.rs = rs;
704:       sctx.pv = pv[diag];
705:       MatPivotCheck(B,A,info,&sctx,i);
706:       if (sctx.newshift) break;
707:       pv[diag] = sctx.pv;
708:     }

710:     if (info->shifttype == (PetscReal)MAT_SHIFT_POSITIVE_DEFINITE && !sctx.newshift && sctx.shift_fraction>0 && sctx.nshift<sctx.nshift_max) {
711:       /*
712:        * if no shift in this attempt & shifting & started shifting & can refine,
713:        * then try lower shift
714:        */
715:       sctx.shift_hi       = sctx.shift_fraction;
716:       sctx.shift_fraction = (sctx.shift_hi+sctx.shift_lo)/2.;
717:       sctx.shift_amount   = sctx.shift_fraction * sctx.shift_top;
718:       sctx.newshift       = PETSC_TRUE;
719:       sctx.nshift++;
720:     }
721:   } while (sctx.newshift);

723:   /* invert diagonal entries for simplier triangular solves */
724:   for (i=0; i<n; i++) {
725:     b->a[diag_offset[i]] = 1.0/b->a[diag_offset[i]];
726:   }
727:   PetscFree(rtmp);
728:   ISRestoreIndices(isicol,&ic);
729:   ISRestoreIndices(isrow,&r);

731:   ISIdentity(isrow,&row_identity);
732:   ISIdentity(isicol,&col_identity);
733:   if (row_identity && col_identity) {
734:     C->ops->solve = MatSolve_SeqAIJ_NaturalOrdering_inplace;
735:   } else {
736:     C->ops->solve = MatSolve_SeqAIJ_inplace;
737:   }
738:   C->ops->solveadd          = MatSolveAdd_SeqAIJ_inplace;
739:   C->ops->solvetranspose    = MatSolveTranspose_SeqAIJ_inplace;
740:   C->ops->solvetransposeadd = MatSolveTransposeAdd_SeqAIJ_inplace;
741:   C->ops->matsolve          = MatMatSolve_SeqAIJ_inplace;

743:   C->assembled    = PETSC_TRUE;
744:   C->preallocated = PETSC_TRUE;

746:   PetscLogFlops(C->cmap->n);
747:   if (sctx.nshift) {
748:     if (info->shifttype == (PetscReal)MAT_SHIFT_POSITIVE_DEFINITE) {
749:       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);
750:     } else if (info->shifttype == (PetscReal)MAT_SHIFT_NONZERO) {
751:       PetscInfo2(A,"number of shift_nz tries %D, shift_amount %g\n",sctx.nshift,(double)sctx.shift_amount);
752:     }
753:   }
754:   (C)->ops->solve          = MatSolve_SeqAIJ_inplace;
755:   (C)->ops->solvetranspose = MatSolveTranspose_SeqAIJ_inplace;

757:   MatSeqAIJCheckInode(C);
758:   return(0);
759: }

761: /*
762:    This routine implements inplace ILU(0) with row or/and column permutations.
763:    Input:
764:      A - original matrix
765:    Output;
766:      A - a->i (rowptr) is same as original rowptr, but factored i-the row is stored in rowperm[i]
767:          a->j (col index) is permuted by the inverse of colperm, then sorted
768:          a->a reordered accordingly with a->j
769:          a->diag (ptr to diagonal elements) is updated.
770: */
771: PetscErrorCode MatLUFactorNumeric_SeqAIJ_InplaceWithPerm(Mat B,Mat A,const MatFactorInfo *info)
772: {
773:   Mat_SeqAIJ      *a    =(Mat_SeqAIJ*)A->data;
774:   IS              isrow = a->row,isicol = a->icol;
775:   PetscErrorCode  ierr;
776:   const PetscInt  *r,*ic,*ics;
777:   PetscInt        i,j,n=A->rmap->n,*ai=a->i,*aj=a->j;
778:   PetscInt        *ajtmp,nz,row;
779:   PetscInt        *diag = a->diag,nbdiag,*pj;
780:   PetscScalar     *rtmp,*pc,multiplier,d;
781:   MatScalar       *pv,*v;
782:   PetscReal       rs;
783:   FactorShiftCtx  sctx;
784:   const MatScalar *aa=a->a,*vtmp;

787:   if (A != B) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"input and output matrix must have same address");

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

792:   if (info->shifttype == (PetscReal) MAT_SHIFT_POSITIVE_DEFINITE) { /* set sctx.shift_top=max{rs} */
793:     const PetscInt *ddiag = a->diag;
794:     sctx.shift_top = info->zeropivot;
795:     for (i=0; i<n; i++) {
796:       /* calculate sum(|aij|)-RealPart(aii), amt of shift needed for this row */
797:       d    = (aa)[ddiag[i]];
798:       rs   = -PetscAbsScalar(d) - PetscRealPart(d);
799:       vtmp = aa+ai[i];
800:       nz   = ai[i+1] - ai[i];
801:       for (j=0; j<nz; j++) rs += PetscAbsScalar(vtmp[j]);
802:       if (rs>sctx.shift_top) sctx.shift_top = rs;
803:     }
804:     sctx.shift_top *= 1.1;
805:     sctx.nshift_max = 5;
806:     sctx.shift_lo   = 0.;
807:     sctx.shift_hi   = 1.;
808:   }

810:   ISGetIndices(isrow,&r);
811:   ISGetIndices(isicol,&ic);
812:   PetscMalloc1(n+1,&rtmp);
813:   PetscMemzero(rtmp,(n+1)*sizeof(PetscScalar));
814:   ics  = ic;

816: #if defined(MV)
817:   sctx.shift_top      = 0.;
818:   sctx.nshift_max     = 0;
819:   sctx.shift_lo       = 0.;
820:   sctx.shift_hi       = 0.;
821:   sctx.shift_fraction = 0.;

823:   if (info->shifttype == (PetscReal)MAT_SHIFT_POSITIVE_DEFINITE) { /* set sctx.shift_top=max{rs} */
824:     sctx.shift_top = 0.;
825:     for (i=0; i<n; i++) {
826:       /* calculate sum(|aij|)-RealPart(aii), amt of shift needed for this row */
827:       d  = (a->a)[diag[i]];
828:       rs = -PetscAbsScalar(d) - PetscRealPart(d);
829:       v  = a->a+ai[i];
830:       nz = ai[i+1] - ai[i];
831:       for (j=0; j<nz; j++) rs += PetscAbsScalar(v[j]);
832:       if (rs>sctx.shift_top) sctx.shift_top = rs;
833:     }
834:     if (sctx.shift_top < info->zeropivot) sctx.shift_top = info->zeropivot;
835:     sctx.shift_top *= 1.1;
836:     sctx.nshift_max = 5;
837:     sctx.shift_lo   = 0.;
838:     sctx.shift_hi   = 1.;
839:   }

841:   sctx.shift_amount = 0.;
842:   sctx.nshift       = 0;
843: #endif

845:   do {
846:     sctx.newshift = PETSC_FALSE;
847:     for (i=0; i<n; i++) {
848:       /* load in initial unfactored row */
849:       nz    = ai[r[i]+1] - ai[r[i]];
850:       ajtmp = aj + ai[r[i]];
851:       v     = a->a + ai[r[i]];
852:       /* sort permuted ajtmp and values v accordingly */
853:       for (j=0; j<nz; j++) ajtmp[j] = ics[ajtmp[j]];
854:       PetscSortIntWithScalarArray(nz,ajtmp,v);

856:       diag[r[i]] = ai[r[i]];
857:       for (j=0; j<nz; j++) {
858:         rtmp[ajtmp[j]] = v[j];
859:         if (ajtmp[j] < i) diag[r[i]]++; /* update a->diag */
860:       }
861:       rtmp[r[i]] += sctx.shift_amount; /* shift the diagonal of the matrix */

863:       row = *ajtmp++;
864:       while  (row < i) {
865:         pc = rtmp + row;
866:         if (*pc != 0.0) {
867:           pv = a->a + diag[r[row]];
868:           pj = aj + diag[r[row]] + 1;

870:           multiplier = *pc / *pv++;
871:           *pc        = multiplier;
872:           nz         = ai[r[row]+1] - diag[r[row]] - 1;
873:           for (j=0; j<nz; j++) rtmp[pj[j]] -= multiplier * pv[j];
874:           PetscLogFlops(1+2*nz);
875:         }
876:         row = *ajtmp++;
877:       }
878:       /* finished row so overwrite it onto a->a */
879:       pv     = a->a + ai[r[i]];
880:       pj     = aj + ai[r[i]];
881:       nz     = ai[r[i]+1] - ai[r[i]];
882:       nbdiag = diag[r[i]] - ai[r[i]]; /* num of entries before the diagonal */

884:       rs = 0.0;
885:       for (j=0; j<nz; j++) {
886:         pv[j] = rtmp[pj[j]];
887:         if (j != nbdiag) rs += PetscAbsScalar(pv[j]);
888:       }

890:       sctx.rs = rs;
891:       sctx.pv = pv[nbdiag];
892:       MatPivotCheck(B,A,info,&sctx,i);
893:       if (sctx.newshift) break;
894:       pv[nbdiag] = sctx.pv;
895:     }

897:     if (info->shifttype == (PetscReal)MAT_SHIFT_POSITIVE_DEFINITE && !sctx.newshift && sctx.shift_fraction>0 && sctx.nshift<sctx.nshift_max) {
898:       /*
899:        * if no shift in this attempt & shifting & started shifting & can refine,
900:        * then try lower shift
901:        */
902:       sctx.shift_hi       = sctx.shift_fraction;
903:       sctx.shift_fraction = (sctx.shift_hi+sctx.shift_lo)/2.;
904:       sctx.shift_amount   = sctx.shift_fraction * sctx.shift_top;
905:       sctx.newshift       = PETSC_TRUE;
906:       sctx.nshift++;
907:     }
908:   } while (sctx.newshift);

910:   /* invert diagonal entries for simplier triangular solves */
911:   for (i=0; i<n; i++) {
912:     a->a[diag[r[i]]] = 1.0/a->a[diag[r[i]]];
913:   }

915:   PetscFree(rtmp);
916:   ISRestoreIndices(isicol,&ic);
917:   ISRestoreIndices(isrow,&r);

919:   A->ops->solve             = MatSolve_SeqAIJ_InplaceWithPerm;
920:   A->ops->solveadd          = MatSolveAdd_SeqAIJ_inplace;
921:   A->ops->solvetranspose    = MatSolveTranspose_SeqAIJ_inplace;
922:   A->ops->solvetransposeadd = MatSolveTransposeAdd_SeqAIJ_inplace;

924:   A->assembled    = PETSC_TRUE;
925:   A->preallocated = PETSC_TRUE;

927:   PetscLogFlops(A->cmap->n);
928:   if (sctx.nshift) {
929:     if (info->shifttype == (PetscReal)MAT_SHIFT_POSITIVE_DEFINITE) {
930:       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);
931:     } else if (info->shifttype == (PetscReal)MAT_SHIFT_NONZERO) {
932:       PetscInfo2(A,"number of shift_nz tries %D, shift_amount %g\n",sctx.nshift,(double)sctx.shift_amount);
933:     }
934:   }
935:   return(0);
936: }

938: /* ----------------------------------------------------------- */
939: PetscErrorCode MatLUFactor_SeqAIJ(Mat A,IS row,IS col,const MatFactorInfo *info)
940: {
942:   Mat            C;

945:   MatGetFactor(A,MATSOLVERPETSC,MAT_FACTOR_LU,&C);
946:   MatLUFactorSymbolic(C,A,row,col,info);
947:   MatLUFactorNumeric(C,A,info);

949:   A->ops->solve          = C->ops->solve;
950:   A->ops->solvetranspose = C->ops->solvetranspose;

952:   MatHeaderMerge(A,&C);
953:   PetscLogObjectParent((PetscObject)A,(PetscObject)((Mat_SeqAIJ*)(A->data))->icol);
954:   return(0);
955: }
956: /* ----------------------------------------------------------- */


959: PetscErrorCode MatSolve_SeqAIJ_inplace(Mat A,Vec bb,Vec xx)
960: {
961:   Mat_SeqAIJ        *a    = (Mat_SeqAIJ*)A->data;
962:   IS                iscol = a->col,isrow = a->row;
963:   PetscErrorCode    ierr;
964:   PetscInt          i, n = A->rmap->n,*vi,*ai = a->i,*aj = a->j;
965:   PetscInt          nz;
966:   const PetscInt    *rout,*cout,*r,*c;
967:   PetscScalar       *x,*tmp,*tmps,sum;
968:   const PetscScalar *b;
969:   const MatScalar   *aa = a->a,*v;

972:   if (!n) return(0);

974:   VecGetArrayRead(bb,&b);
975:   VecGetArray(xx,&x);
976:   tmp  = a->solve_work;

978:   ISGetIndices(isrow,&rout); r = rout;
979:   ISGetIndices(iscol,&cout); c = cout + (n-1);

981:   /* forward solve the lower triangular */
982:   tmp[0] = b[*r++];
983:   tmps   = tmp;
984:   for (i=1; i<n; i++) {
985:     v   = aa + ai[i];
986:     vi  = aj + ai[i];
987:     nz  = a->diag[i] - ai[i];
988:     sum = b[*r++];
989:     PetscSparseDenseMinusDot(sum,tmps,v,vi,nz);
990:     tmp[i] = sum;
991:   }

993:   /* backward solve the upper triangular */
994:   for (i=n-1; i>=0; i--) {
995:     v   = aa + a->diag[i] + 1;
996:     vi  = aj + a->diag[i] + 1;
997:     nz  = ai[i+1] - a->diag[i] - 1;
998:     sum = tmp[i];
999:     PetscSparseDenseMinusDot(sum,tmps,v,vi,nz);
1000:     x[*c--] = tmp[i] = sum*aa[a->diag[i]];
1001:   }

1003:   ISRestoreIndices(isrow,&rout);
1004:   ISRestoreIndices(iscol,&cout);
1005:   VecRestoreArrayRead(bb,&b);
1006:   VecRestoreArray(xx,&x);
1007:   PetscLogFlops(2.0*a->nz - A->cmap->n);
1008:   return(0);
1009: }

1011: PetscErrorCode MatMatSolve_SeqAIJ_inplace(Mat A,Mat B,Mat X)
1012: {
1013:   Mat_SeqAIJ      *a    = (Mat_SeqAIJ*)A->data;
1014:   IS              iscol = a->col,isrow = a->row;
1015:   PetscErrorCode  ierr;
1016:   PetscInt        i, n = A->rmap->n,*vi,*ai = a->i,*aj = a->j;
1017:   PetscInt        nz,neq;
1018:   const PetscInt  *rout,*cout,*r,*c;
1019:   PetscScalar     *x,*b,*tmp,*tmps,sum;
1020:   const MatScalar *aa = a->a,*v;
1021:   PetscBool       bisdense,xisdense;

1024:   if (!n) return(0);

1026:   PetscObjectTypeCompare((PetscObject)B,MATSEQDENSE,&bisdense);
1027:   if (!bisdense) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"B matrix must be a SeqDense matrix");
1028:   PetscObjectTypeCompare((PetscObject)X,MATSEQDENSE,&xisdense);
1029:   if (!xisdense) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"X matrix must be a SeqDense matrix");

1031:   MatDenseGetArray(B,&b);
1032:   MatDenseGetArray(X,&x);

1034:   tmp  = a->solve_work;
1035:   ISGetIndices(isrow,&rout); r = rout;
1036:   ISGetIndices(iscol,&cout); c = cout;

1038:   for (neq=0; neq<B->cmap->n; neq++) {
1039:     /* forward solve the lower triangular */
1040:     tmp[0] = b[r[0]];
1041:     tmps   = tmp;
1042:     for (i=1; i<n; i++) {
1043:       v   = aa + ai[i];
1044:       vi  = aj + ai[i];
1045:       nz  = a->diag[i] - ai[i];
1046:       sum = b[r[i]];
1047:       PetscSparseDenseMinusDot(sum,tmps,v,vi,nz);
1048:       tmp[i] = sum;
1049:     }
1050:     /* backward solve the upper triangular */
1051:     for (i=n-1; i>=0; i--) {
1052:       v   = aa + a->diag[i] + 1;
1053:       vi  = aj + a->diag[i] + 1;
1054:       nz  = ai[i+1] - a->diag[i] - 1;
1055:       sum = tmp[i];
1056:       PetscSparseDenseMinusDot(sum,tmps,v,vi,nz);
1057:       x[c[i]] = tmp[i] = sum*aa[a->diag[i]];
1058:     }

1060:     b += n;
1061:     x += n;
1062:   }
1063:   ISRestoreIndices(isrow,&rout);
1064:   ISRestoreIndices(iscol,&cout);
1065:   MatDenseRestoreArray(B,&b);
1066:   MatDenseRestoreArray(X,&x);
1067:   PetscLogFlops(B->cmap->n*(2.0*a->nz - n));
1068:   return(0);
1069: }

1071: PetscErrorCode MatMatSolve_SeqAIJ(Mat A,Mat B,Mat X)
1072: {
1073:   Mat_SeqAIJ      *a    = (Mat_SeqAIJ*)A->data;
1074:   IS              iscol = a->col,isrow = a->row;
1075:   PetscErrorCode  ierr;
1076:   PetscInt        i, n = A->rmap->n,*vi,*ai = a->i,*aj = a->j,*adiag = a->diag;
1077:   PetscInt        nz,neq;
1078:   const PetscInt  *rout,*cout,*r,*c;
1079:   PetscScalar     *x,*b,*tmp,sum;
1080:   const MatScalar *aa = a->a,*v;
1081:   PetscBool       bisdense,xisdense;

1084:   if (!n) return(0);

1086:   PetscObjectTypeCompare((PetscObject)B,MATSEQDENSE,&bisdense);
1087:   if (!bisdense) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"B matrix must be a SeqDense matrix");
1088:   PetscObjectTypeCompare((PetscObject)X,MATSEQDENSE,&xisdense);
1089:   if (!xisdense) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"X matrix must be a SeqDense matrix");

1091:   MatDenseGetArray(B,&b);
1092:   MatDenseGetArray(X,&x);

1094:   tmp  = a->solve_work;
1095:   ISGetIndices(isrow,&rout); r = rout;
1096:   ISGetIndices(iscol,&cout); c = cout;

1098:   for (neq=0; neq<B->cmap->n; neq++) {
1099:     /* forward solve the lower triangular */
1100:     tmp[0] = b[r[0]];
1101:     v      = aa;
1102:     vi     = aj;
1103:     for (i=1; i<n; i++) {
1104:       nz  = ai[i+1] - ai[i];
1105:       sum = b[r[i]];
1106:       PetscSparseDenseMinusDot(sum,tmp,v,vi,nz);
1107:       tmp[i] = sum;
1108:       v     += nz; vi += nz;
1109:     }

1111:     /* backward solve the upper triangular */
1112:     for (i=n-1; i>=0; i--) {
1113:       v   = aa + adiag[i+1]+1;
1114:       vi  = aj + adiag[i+1]+1;
1115:       nz  = adiag[i]-adiag[i+1]-1;
1116:       sum = tmp[i];
1117:       PetscSparseDenseMinusDot(sum,tmp,v,vi,nz);
1118:       x[c[i]] = tmp[i] = sum*v[nz]; /* v[nz] = aa[adiag[i]] */
1119:     }

1121:     b += n;
1122:     x += n;
1123:   }
1124:   ISRestoreIndices(isrow,&rout);
1125:   ISRestoreIndices(iscol,&cout);
1126:   MatDenseRestoreArray(B,&b);
1127:   MatDenseRestoreArray(X,&x);
1128:   PetscLogFlops(B->cmap->n*(2.0*a->nz - n));
1129:   return(0);
1130: }

1132: PetscErrorCode MatSolve_SeqAIJ_InplaceWithPerm(Mat A,Vec bb,Vec xx)
1133: {
1134:   Mat_SeqAIJ      *a    = (Mat_SeqAIJ*)A->data;
1135:   IS              iscol = a->col,isrow = a->row;
1136:   PetscErrorCode  ierr;
1137:   const PetscInt  *r,*c,*rout,*cout;
1138:   PetscInt        i, n = A->rmap->n,*vi,*ai = a->i,*aj = a->j;
1139:   PetscInt        nz,row;
1140:   PetscScalar     *x,*b,*tmp,*tmps,sum;
1141:   const MatScalar *aa = a->a,*v;

1144:   if (!n) return(0);

1146:   VecGetArray(bb,&b);
1147:   VecGetArray(xx,&x);
1148:   tmp  = a->solve_work;

1150:   ISGetIndices(isrow,&rout); r = rout;
1151:   ISGetIndices(iscol,&cout); c = cout + (n-1);

1153:   /* forward solve the lower triangular */
1154:   tmp[0] = b[*r++];
1155:   tmps   = tmp;
1156:   for (row=1; row<n; row++) {
1157:     i   = rout[row]; /* permuted row */
1158:     v   = aa + ai[i];
1159:     vi  = aj + ai[i];
1160:     nz  = a->diag[i] - ai[i];
1161:     sum = b[*r++];
1162:     PetscSparseDenseMinusDot(sum,tmps,v,vi,nz);
1163:     tmp[row] = sum;
1164:   }

1166:   /* backward solve the upper triangular */
1167:   for (row=n-1; row>=0; row--) {
1168:     i   = rout[row]; /* permuted row */
1169:     v   = aa + a->diag[i] + 1;
1170:     vi  = aj + a->diag[i] + 1;
1171:     nz  = ai[i+1] - a->diag[i] - 1;
1172:     sum = tmp[row];
1173:     PetscSparseDenseMinusDot(sum,tmps,v,vi,nz);
1174:     x[*c--] = tmp[row] = sum*aa[a->diag[i]];
1175:   }

1177:   ISRestoreIndices(isrow,&rout);
1178:   ISRestoreIndices(iscol,&cout);
1179:   VecRestoreArray(bb,&b);
1180:   VecRestoreArray(xx,&x);
1181:   PetscLogFlops(2.0*a->nz - A->cmap->n);
1182:   return(0);
1183: }

1185: /* ----------------------------------------------------------- */
1186: #include <../src/mat/impls/aij/seq/ftn-kernels/fsolve.h>
1187: PetscErrorCode MatSolve_SeqAIJ_NaturalOrdering_inplace(Mat A,Vec bb,Vec xx)
1188: {
1189:   Mat_SeqAIJ        *a = (Mat_SeqAIJ*)A->data;
1190:   PetscErrorCode    ierr;
1191:   PetscInt          n   = A->rmap->n;
1192:   const PetscInt    *ai = a->i,*aj = a->j,*adiag = a->diag;
1193:   PetscScalar       *x;
1194:   const PetscScalar *b;
1195:   const MatScalar   *aa = a->a;
1196: #if !defined(PETSC_USE_FORTRAN_KERNEL_SOLVEAIJ)
1197:   PetscInt        adiag_i,i,nz,ai_i;
1198:   const PetscInt  *vi;
1199:   const MatScalar *v;
1200:   PetscScalar     sum;
1201: #endif

1204:   if (!n) return(0);

1206:   VecGetArrayRead(bb,&b);
1207:   VecGetArray(xx,&x);

1209: #if defined(PETSC_USE_FORTRAN_KERNEL_SOLVEAIJ)
1210:   fortransolveaij_(&n,x,ai,aj,adiag,aa,b);
1211: #else
1212:   /* forward solve the lower triangular */
1213:   x[0] = b[0];
1214:   for (i=1; i<n; i++) {
1215:     ai_i = ai[i];
1216:     v    = aa + ai_i;
1217:     vi   = aj + ai_i;
1218:     nz   = adiag[i] - ai_i;
1219:     sum  = b[i];
1220:     PetscSparseDenseMinusDot(sum,x,v,vi,nz);
1221:     x[i] = sum;
1222:   }

1224:   /* backward solve the upper triangular */
1225:   for (i=n-1; i>=0; i--) {
1226:     adiag_i = adiag[i];
1227:     v       = aa + adiag_i + 1;
1228:     vi      = aj + adiag_i + 1;
1229:     nz      = ai[i+1] - adiag_i - 1;
1230:     sum     = x[i];
1231:     PetscSparseDenseMinusDot(sum,x,v,vi,nz);
1232:     x[i] = sum*aa[adiag_i];
1233:   }
1234: #endif
1235:   PetscLogFlops(2.0*a->nz - A->cmap->n);
1236:   VecRestoreArrayRead(bb,&b);
1237:   VecRestoreArray(xx,&x);
1238:   return(0);
1239: }

1241: PetscErrorCode MatSolveAdd_SeqAIJ_inplace(Mat A,Vec bb,Vec yy,Vec xx)
1242: {
1243:   Mat_SeqAIJ        *a    = (Mat_SeqAIJ*)A->data;
1244:   IS                iscol = a->col,isrow = a->row;
1245:   PetscErrorCode    ierr;
1246:   PetscInt          i, n = A->rmap->n,j;
1247:   PetscInt          nz;
1248:   const PetscInt    *rout,*cout,*r,*c,*vi,*ai = a->i,*aj = a->j;
1249:   PetscScalar       *x,*tmp,sum;
1250:   const PetscScalar *b;
1251:   const MatScalar   *aa = a->a,*v;

1254:   if (yy != xx) {VecCopy(yy,xx);}

1256:   VecGetArrayRead(bb,&b);
1257:   VecGetArray(xx,&x);
1258:   tmp  = a->solve_work;

1260:   ISGetIndices(isrow,&rout); r = rout;
1261:   ISGetIndices(iscol,&cout); c = cout + (n-1);

1263:   /* forward solve the lower triangular */
1264:   tmp[0] = b[*r++];
1265:   for (i=1; i<n; i++) {
1266:     v   = aa + ai[i];
1267:     vi  = aj + ai[i];
1268:     nz  = a->diag[i] - ai[i];
1269:     sum = b[*r++];
1270:     for (j=0; j<nz; j++) sum -= v[j]*tmp[vi[j]];
1271:     tmp[i] = sum;
1272:   }

1274:   /* backward solve the upper triangular */
1275:   for (i=n-1; i>=0; i--) {
1276:     v   = aa + a->diag[i] + 1;
1277:     vi  = aj + a->diag[i] + 1;
1278:     nz  = ai[i+1] - a->diag[i] - 1;
1279:     sum = tmp[i];
1280:     for (j=0; j<nz; j++) sum -= v[j]*tmp[vi[j]];
1281:     tmp[i]   = sum*aa[a->diag[i]];
1282:     x[*c--] += tmp[i];
1283:   }

1285:   ISRestoreIndices(isrow,&rout);
1286:   ISRestoreIndices(iscol,&cout);
1287:   VecRestoreArrayRead(bb,&b);
1288:   VecRestoreArray(xx,&x);
1289:   PetscLogFlops(2.0*a->nz);
1290:   return(0);
1291: }

1293: PetscErrorCode MatSolveAdd_SeqAIJ(Mat A,Vec bb,Vec yy,Vec xx)
1294: {
1295:   Mat_SeqAIJ        *a    = (Mat_SeqAIJ*)A->data;
1296:   IS                iscol = a->col,isrow = a->row;
1297:   PetscErrorCode    ierr;
1298:   PetscInt          i, n = A->rmap->n,j;
1299:   PetscInt          nz;
1300:   const PetscInt    *rout,*cout,*r,*c,*vi,*ai = a->i,*aj = a->j,*adiag = a->diag;
1301:   PetscScalar       *x,*tmp,sum;
1302:   const PetscScalar *b;
1303:   const MatScalar   *aa = a->a,*v;

1306:   if (yy != xx) {VecCopy(yy,xx);}

1308:   VecGetArrayRead(bb,&b);
1309:   VecGetArray(xx,&x);
1310:   tmp  = a->solve_work;

1312:   ISGetIndices(isrow,&rout); r = rout;
1313:   ISGetIndices(iscol,&cout); c = cout;

1315:   /* forward solve the lower triangular */
1316:   tmp[0] = b[r[0]];
1317:   v      = aa;
1318:   vi     = aj;
1319:   for (i=1; i<n; i++) {
1320:     nz  = ai[i+1] - ai[i];
1321:     sum = b[r[i]];
1322:     for (j=0; j<nz; j++) sum -= v[j]*tmp[vi[j]];
1323:     tmp[i] = sum;
1324:     v     += nz;
1325:     vi    += nz;
1326:   }

1328:   /* backward solve the upper triangular */
1329:   v  = aa + adiag[n-1];
1330:   vi = aj + adiag[n-1];
1331:   for (i=n-1; i>=0; i--) {
1332:     nz  = adiag[i] - adiag[i+1] - 1;
1333:     sum = tmp[i];
1334:     for (j=0; j<nz; j++) sum -= v[j]*tmp[vi[j]];
1335:     tmp[i]   = sum*v[nz];
1336:     x[c[i]] += tmp[i];
1337:     v       += nz+1; vi += nz+1;
1338:   }

1340:   ISRestoreIndices(isrow,&rout);
1341:   ISRestoreIndices(iscol,&cout);
1342:   VecRestoreArrayRead(bb,&b);
1343:   VecRestoreArray(xx,&x);
1344:   PetscLogFlops(2.0*a->nz);
1345:   return(0);
1346: }

1348: PetscErrorCode MatSolveTranspose_SeqAIJ_inplace(Mat A,Vec bb,Vec xx)
1349: {
1350:   Mat_SeqAIJ        *a    = (Mat_SeqAIJ*)A->data;
1351:   IS                iscol = a->col,isrow = a->row;
1352:   PetscErrorCode    ierr;
1353:   const PetscInt    *rout,*cout,*r,*c,*diag = a->diag,*ai = a->i,*aj = a->j,*vi;
1354:   PetscInt          i,n = A->rmap->n,j;
1355:   PetscInt          nz;
1356:   PetscScalar       *x,*tmp,s1;
1357:   const MatScalar   *aa = a->a,*v;
1358:   const PetscScalar *b;

1361:   VecGetArrayRead(bb,&b);
1362:   VecGetArray(xx,&x);
1363:   tmp  = a->solve_work;

1365:   ISGetIndices(isrow,&rout); r = rout;
1366:   ISGetIndices(iscol,&cout); c = cout;

1368:   /* copy the b into temp work space according to permutation */
1369:   for (i=0; i<n; i++) tmp[i] = b[c[i]];

1371:   /* forward solve the U^T */
1372:   for (i=0; i<n; i++) {
1373:     v   = aa + diag[i];
1374:     vi  = aj + diag[i] + 1;
1375:     nz  = ai[i+1] - diag[i] - 1;
1376:     s1  = tmp[i];
1377:     s1 *= (*v++);  /* multiply by inverse of diagonal entry */
1378:     for (j=0; j<nz; j++) tmp[vi[j]] -= s1*v[j];
1379:     tmp[i] = s1;
1380:   }

1382:   /* backward solve the L^T */
1383:   for (i=n-1; i>=0; i--) {
1384:     v  = aa + diag[i] - 1;
1385:     vi = aj + diag[i] - 1;
1386:     nz = diag[i] - ai[i];
1387:     s1 = tmp[i];
1388:     for (j=0; j>-nz; j--) tmp[vi[j]] -= s1*v[j];
1389:   }

1391:   /* copy tmp into x according to permutation */
1392:   for (i=0; i<n; i++) x[r[i]] = tmp[i];

1394:   ISRestoreIndices(isrow,&rout);
1395:   ISRestoreIndices(iscol,&cout);
1396:   VecRestoreArrayRead(bb,&b);
1397:   VecRestoreArray(xx,&x);

1399:   PetscLogFlops(2.0*a->nz-A->cmap->n);
1400:   return(0);
1401: }

1403: PetscErrorCode MatSolveTranspose_SeqAIJ(Mat A,Vec bb,Vec xx)
1404: {
1405:   Mat_SeqAIJ        *a    = (Mat_SeqAIJ*)A->data;
1406:   IS                iscol = a->col,isrow = a->row;
1407:   PetscErrorCode    ierr;
1408:   const PetscInt    *rout,*cout,*r,*c,*adiag = a->diag,*ai = a->i,*aj = a->j,*vi;
1409:   PetscInt          i,n = A->rmap->n,j;
1410:   PetscInt          nz;
1411:   PetscScalar       *x,*tmp,s1;
1412:   const MatScalar   *aa = a->a,*v;
1413:   const PetscScalar *b;

1416:   VecGetArrayRead(bb,&b);
1417:   VecGetArray(xx,&x);
1418:   tmp  = a->solve_work;

1420:   ISGetIndices(isrow,&rout); r = rout;
1421:   ISGetIndices(iscol,&cout); c = cout;

1423:   /* copy the b into temp work space according to permutation */
1424:   for (i=0; i<n; i++) tmp[i] = b[c[i]];

1426:   /* forward solve the U^T */
1427:   for (i=0; i<n; i++) {
1428:     v   = aa + adiag[i+1] + 1;
1429:     vi  = aj + adiag[i+1] + 1;
1430:     nz  = adiag[i] - adiag[i+1] - 1;
1431:     s1  = tmp[i];
1432:     s1 *= v[nz];  /* multiply by inverse of diagonal entry */
1433:     for (j=0; j<nz; j++) tmp[vi[j]] -= s1*v[j];
1434:     tmp[i] = s1;
1435:   }

1437:   /* backward solve the L^T */
1438:   for (i=n-1; i>=0; i--) {
1439:     v  = aa + ai[i];
1440:     vi = aj + ai[i];
1441:     nz = ai[i+1] - ai[i];
1442:     s1 = tmp[i];
1443:     for (j=0; j<nz; j++) tmp[vi[j]] -= s1*v[j];
1444:   }

1446:   /* copy tmp into x according to permutation */
1447:   for (i=0; i<n; i++) x[r[i]] = tmp[i];

1449:   ISRestoreIndices(isrow,&rout);
1450:   ISRestoreIndices(iscol,&cout);
1451:   VecRestoreArrayRead(bb,&b);
1452:   VecRestoreArray(xx,&x);

1454:   PetscLogFlops(2.0*a->nz-A->cmap->n);
1455:   return(0);
1456: }

1458: PetscErrorCode MatSolveTransposeAdd_SeqAIJ_inplace(Mat A,Vec bb,Vec zz,Vec xx)
1459: {
1460:   Mat_SeqAIJ        *a    = (Mat_SeqAIJ*)A->data;
1461:   IS                iscol = a->col,isrow = a->row;
1462:   PetscErrorCode    ierr;
1463:   const PetscInt    *rout,*cout,*r,*c,*diag = a->diag,*ai = a->i,*aj = a->j,*vi;
1464:   PetscInt          i,n = A->rmap->n,j;
1465:   PetscInt          nz;
1466:   PetscScalar       *x,*tmp,s1;
1467:   const MatScalar   *aa = a->a,*v;
1468:   const PetscScalar *b;

1471:   if (zz != xx) {VecCopy(zz,xx);}
1472:   VecGetArrayRead(bb,&b);
1473:   VecGetArray(xx,&x);
1474:   tmp  = a->solve_work;

1476:   ISGetIndices(isrow,&rout); r = rout;
1477:   ISGetIndices(iscol,&cout); c = cout;

1479:   /* copy the b into temp work space according to permutation */
1480:   for (i=0; i<n; i++) tmp[i] = b[c[i]];

1482:   /* forward solve the U^T */
1483:   for (i=0; i<n; i++) {
1484:     v   = aa + diag[i];
1485:     vi  = aj + diag[i] + 1;
1486:     nz  = ai[i+1] - diag[i] - 1;
1487:     s1  = tmp[i];
1488:     s1 *= (*v++);  /* multiply by inverse of diagonal entry */
1489:     for (j=0; j<nz; j++) tmp[vi[j]] -= s1*v[j];
1490:     tmp[i] = s1;
1491:   }

1493:   /* backward solve the L^T */
1494:   for (i=n-1; i>=0; i--) {
1495:     v  = aa + diag[i] - 1;
1496:     vi = aj + diag[i] - 1;
1497:     nz = diag[i] - ai[i];
1498:     s1 = tmp[i];
1499:     for (j=0; j>-nz; j--) tmp[vi[j]] -= s1*v[j];
1500:   }

1502:   /* copy tmp into x according to permutation */
1503:   for (i=0; i<n; i++) x[r[i]] += tmp[i];

1505:   ISRestoreIndices(isrow,&rout);
1506:   ISRestoreIndices(iscol,&cout);
1507:   VecRestoreArrayRead(bb,&b);
1508:   VecRestoreArray(xx,&x);

1510:   PetscLogFlops(2.0*a->nz-A->cmap->n);
1511:   return(0);
1512: }

1514: PetscErrorCode MatSolveTransposeAdd_SeqAIJ(Mat A,Vec bb,Vec zz,Vec xx)
1515: {
1516:   Mat_SeqAIJ        *a    = (Mat_SeqAIJ*)A->data;
1517:   IS                iscol = a->col,isrow = a->row;
1518:   PetscErrorCode    ierr;
1519:   const PetscInt    *rout,*cout,*r,*c,*adiag = a->diag,*ai = a->i,*aj = a->j,*vi;
1520:   PetscInt          i,n = A->rmap->n,j;
1521:   PetscInt          nz;
1522:   PetscScalar       *x,*tmp,s1;
1523:   const MatScalar   *aa = a->a,*v;
1524:   const PetscScalar *b;

1527:   if (zz != xx) {VecCopy(zz,xx);}
1528:   VecGetArrayRead(bb,&b);
1529:   VecGetArray(xx,&x);
1530:   tmp  = a->solve_work;

1532:   ISGetIndices(isrow,&rout); r = rout;
1533:   ISGetIndices(iscol,&cout); c = cout;

1535:   /* copy the b into temp work space according to permutation */
1536:   for (i=0; i<n; i++) tmp[i] = b[c[i]];

1538:   /* forward solve the U^T */
1539:   for (i=0; i<n; i++) {
1540:     v   = aa + adiag[i+1] + 1;
1541:     vi  = aj + adiag[i+1] + 1;
1542:     nz  = adiag[i] - adiag[i+1] - 1;
1543:     s1  = tmp[i];
1544:     s1 *= v[nz];  /* multiply by inverse of diagonal entry */
1545:     for (j=0; j<nz; j++) tmp[vi[j]] -= s1*v[j];
1546:     tmp[i] = s1;
1547:   }


1550:   /* backward solve the L^T */
1551:   for (i=n-1; i>=0; i--) {
1552:     v  = aa + ai[i];
1553:     vi = aj + ai[i];
1554:     nz = ai[i+1] - ai[i];
1555:     s1 = tmp[i];
1556:     for (j=0; j<nz; j++) tmp[vi[j]] -= s1*v[j];
1557:   }

1559:   /* copy tmp into x according to permutation */
1560:   for (i=0; i<n; i++) x[r[i]] += tmp[i];

1562:   ISRestoreIndices(isrow,&rout);
1563:   ISRestoreIndices(iscol,&cout);
1564:   VecRestoreArrayRead(bb,&b);
1565:   VecRestoreArray(xx,&x);

1567:   PetscLogFlops(2.0*a->nz-A->cmap->n);
1568:   return(0);
1569: }

1571: /* ----------------------------------------------------------------*/

1573: /*
1574:    ilu() under revised new data structure.
1575:    Factored arrays bj and ba are stored as
1576:      L(0,:), L(1,:), ...,L(n-1,:),  U(n-1,:),...,U(i,:),U(i-1,:),...,U(0,:)

1578:    bi=fact->i is an array of size n+1, in which
1579:    bi+
1580:      bi[i]:  points to 1st entry of L(i,:),i=0,...,n-1
1581:      bi[n]:  points to L(n-1,n-1)+1

1583:   bdiag=fact->diag is an array of size n+1,in which
1584:      bdiag[i]: points to diagonal of U(i,:), i=0,...,n-1
1585:      bdiag[n]: points to entry of U(n-1,0)-1

1587:    U(i,:) contains bdiag[i] as its last entry, i.e.,
1588:     U(i,:) = (u[i,i+1],...,u[i,n-1],diag[i])
1589: */
1590: PetscErrorCode MatILUFactorSymbolic_SeqAIJ_ilu0(Mat fact,Mat A,IS isrow,IS iscol,const MatFactorInfo *info)
1591: {
1592:   Mat_SeqAIJ     *a = (Mat_SeqAIJ*)A->data,*b;
1594:   const PetscInt n=A->rmap->n,*ai=a->i,*aj,*adiag=a->diag;
1595:   PetscInt       i,j,k=0,nz,*bi,*bj,*bdiag;
1596:   IS             isicol;

1599:   ISInvertPermutation(iscol,PETSC_DECIDE,&isicol);
1600:   MatDuplicateNoCreate_SeqAIJ(fact,A,MAT_DO_NOT_COPY_VALUES,PETSC_FALSE);
1601:   b    = (Mat_SeqAIJ*)(fact)->data;

1603:   /* allocate matrix arrays for new data structure */
1604:   PetscMalloc3(ai[n]+1,&b->a,ai[n]+1,&b->j,n+1,&b->i);
1605:   PetscLogObjectMemory((PetscObject)fact,ai[n]*(sizeof(PetscScalar)+sizeof(PetscInt))+(n+1)*sizeof(PetscInt));

1607:   b->singlemalloc = PETSC_TRUE;
1608:   if (!b->diag) {
1609:     PetscMalloc1(n+1,&b->diag);
1610:     PetscLogObjectMemory((PetscObject)fact,(n+1)*sizeof(PetscInt));
1611:   }
1612:   bdiag = b->diag;

1614:   if (n > 0) {
1615:     PetscMemzero(b->a,(ai[n])*sizeof(MatScalar));
1616:   }

1618:   /* set bi and bj with new data structure */
1619:   bi = b->i;
1620:   bj = b->j;

1622:   /* L part */
1623:   bi[0] = 0;
1624:   for (i=0; i<n; i++) {
1625:     nz      = adiag[i] - ai[i];
1626:     bi[i+1] = bi[i] + nz;
1627:     aj      = a->j + ai[i];
1628:     for (j=0; j<nz; j++) {
1629:       /*   *bj = aj[j]; bj++; */
1630:       bj[k++] = aj[j];
1631:     }
1632:   }

1634:   /* U part */
1635:   bdiag[n] = bi[n]-1;
1636:   for (i=n-1; i>=0; i--) {
1637:     nz = ai[i+1] - adiag[i] - 1;
1638:     aj = a->j + adiag[i] + 1;
1639:     for (j=0; j<nz; j++) {
1640:       /*      *bj = aj[j]; bj++; */
1641:       bj[k++] = aj[j];
1642:     }
1643:     /* diag[i] */
1644:     /*    *bj = i; bj++; */
1645:     bj[k++]  = i;
1646:     bdiag[i] = bdiag[i+1] + nz + 1;
1647:   }

1649:   fact->factortype             = MAT_FACTOR_ILU;
1650:   fact->info.factor_mallocs    = 0;
1651:   fact->info.fill_ratio_given  = info->fill;
1652:   fact->info.fill_ratio_needed = 1.0;
1653:   fact->ops->lufactornumeric   = MatLUFactorNumeric_SeqAIJ;
1654:   MatSeqAIJCheckInode_FactorLU(fact);

1656:   b       = (Mat_SeqAIJ*)(fact)->data;
1657:   b->row  = isrow;
1658:   b->col  = iscol;
1659:   b->icol = isicol;
1660:   PetscMalloc1(fact->rmap->n+1,&b->solve_work);
1661:   PetscObjectReference((PetscObject)isrow);
1662:   PetscObjectReference((PetscObject)iscol);
1663:   return(0);
1664: }

1666: PetscErrorCode MatILUFactorSymbolic_SeqAIJ(Mat fact,Mat A,IS isrow,IS iscol,const MatFactorInfo *info)
1667: {
1668:   Mat_SeqAIJ         *a = (Mat_SeqAIJ*)A->data,*b;
1669:   IS                 isicol;
1670:   PetscErrorCode     ierr;
1671:   const PetscInt     *r,*ic;
1672:   PetscInt           n=A->rmap->n,*ai=a->i,*aj=a->j;
1673:   PetscInt           *bi,*cols,nnz,*cols_lvl;
1674:   PetscInt           *bdiag,prow,fm,nzbd,reallocs=0,dcount=0;
1675:   PetscInt           i,levels,diagonal_fill;
1676:   PetscBool          col_identity,row_identity,missing;
1677:   PetscReal          f;
1678:   PetscInt           nlnk,*lnk,*lnk_lvl=NULL;
1679:   PetscBT            lnkbt;
1680:   PetscInt           nzi,*bj,**bj_ptr,**bjlvl_ptr;
1681:   PetscFreeSpaceList free_space    =NULL,current_space=NULL;
1682:   PetscFreeSpaceList free_space_lvl=NULL,current_space_lvl=NULL;

1685:   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);
1686:   MatMissingDiagonal(A,&missing,&i);
1687:   if (missing) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Matrix is missing diagonal entry %D",i);
1688: 
1689:   levels = (PetscInt)info->levels;
1690:   ISIdentity(isrow,&row_identity);
1691:   ISIdentity(iscol,&col_identity);
1692:   if (!levels && row_identity && col_identity) {
1693:     /* special case: ilu(0) with natural ordering */
1694:     MatILUFactorSymbolic_SeqAIJ_ilu0(fact,A,isrow,iscol,info);
1695:     if (a->inode.size) {
1696:       fact->ops->lufactornumeric = MatLUFactorNumeric_SeqAIJ_Inode;
1697:     }
1698:     return(0);
1699:   }

1701:   ISInvertPermutation(iscol,PETSC_DECIDE,&isicol);
1702:   ISGetIndices(isrow,&r);
1703:   ISGetIndices(isicol,&ic);

1705:   /* get new row and diagonal pointers, must be allocated separately because they will be given to the Mat_SeqAIJ and freed separately */
1706:   PetscMalloc1(n+1,&bi);
1707:   PetscMalloc1(n+1,&bdiag);
1708:   bi[0] = bdiag[0] = 0;
1709:   PetscMalloc2(n,&bj_ptr,n,&bjlvl_ptr);

1711:   /* create a linked list for storing column indices of the active row */
1712:   nlnk = n + 1;
1713:   PetscIncompleteLLCreate(n,n,nlnk,lnk,lnk_lvl,lnkbt);

1715:   /* initial FreeSpace size is f*(ai[n]+1) */
1716:   f                 = info->fill;
1717:   diagonal_fill     = (PetscInt)info->diagonal_fill;
1718:   PetscFreeSpaceGet(PetscRealIntMultTruncate(f,ai[n]+1),&free_space);
1719:   current_space     = free_space;
1720:   PetscFreeSpaceGet(PetscRealIntMultTruncate(f,ai[n]+1),&free_space_lvl);
1721:   current_space_lvl = free_space_lvl;
1722:   for (i=0; i<n; i++) {
1723:     nzi = 0;
1724:     /* copy current row into linked list */
1725:     nnz = ai[r[i]+1] - ai[r[i]];
1726:     if (!nnz) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_MAT_LU_ZRPVT,"Empty row in matrix: row in original ordering %D in permuted ordering %D",r[i],i);
1727:     cols   = aj + ai[r[i]];
1728:     lnk[i] = -1; /* marker to indicate if diagonal exists */
1729:     PetscIncompleteLLInit(nnz,cols,n,ic,nlnk,lnk,lnk_lvl,lnkbt);
1730:     nzi   += nlnk;

1732:     /* make sure diagonal entry is included */
1733:     if (diagonal_fill && lnk[i] == -1) {
1734:       fm = n;
1735:       while (lnk[fm] < i) fm = lnk[fm];
1736:       lnk[i]     = lnk[fm]; /* insert diagonal into linked list */
1737:       lnk[fm]    = i;
1738:       lnk_lvl[i] = 0;
1739:       nzi++; dcount++;
1740:     }

1742:     /* add pivot rows into the active row */
1743:     nzbd = 0;
1744:     prow = lnk[n];
1745:     while (prow < i) {
1746:       nnz      = bdiag[prow];
1747:       cols     = bj_ptr[prow] + nnz + 1;
1748:       cols_lvl = bjlvl_ptr[prow] + nnz + 1;
1749:       nnz      = bi[prow+1] - bi[prow] - nnz - 1;
1750:       PetscILULLAddSorted(nnz,cols,levels,cols_lvl,prow,nlnk,lnk,lnk_lvl,lnkbt,prow);
1751:       nzi     += nlnk;
1752:       prow     = lnk[prow];
1753:       nzbd++;
1754:     }
1755:     bdiag[i] = nzbd;
1756:     bi[i+1]  = bi[i] + nzi;
1757:     /* if free space is not available, make more free space */
1758:     if (current_space->local_remaining<nzi) {
1759:       nnz  = PetscIntMultTruncate(2,PetscIntMultTruncate(nzi,n - i)); /* estimated and max additional space needed */
1760:       PetscFreeSpaceGet(nnz,&current_space);
1761:       PetscFreeSpaceGet(nnz,&current_space_lvl);
1762:       reallocs++;
1763:     }

1765:     /* copy data into free_space and free_space_lvl, then initialize lnk */
1766:     PetscIncompleteLLClean(n,n,nzi,lnk,lnk_lvl,current_space->array,current_space_lvl->array,lnkbt);
1767:     bj_ptr[i]    = current_space->array;
1768:     bjlvl_ptr[i] = current_space_lvl->array;

1770:     /* make sure the active row i has diagonal entry */
1771:     if (*(bj_ptr[i]+bdiag[i]) != i) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_MAT_LU_ZRPVT,"Row %D has missing diagonal in factored matrix\ntry running with -pc_factor_nonzeros_along_diagonal or -pc_factor_diagonal_fill",i);

1773:     current_space->array               += nzi;
1774:     current_space->local_used          += nzi;
1775:     current_space->local_remaining     -= nzi;
1776:     current_space_lvl->array           += nzi;
1777:     current_space_lvl->local_used      += nzi;
1778:     current_space_lvl->local_remaining -= nzi;
1779:   }

1781:   ISRestoreIndices(isrow,&r);
1782:   ISRestoreIndices(isicol,&ic);
1783:   /* copy free_space into bj and free free_space; set bi, bj, bdiag in new datastructure; */
1784:   PetscMalloc1(bi[n]+1,&bj);
1785:   PetscFreeSpaceContiguous_LU(&free_space,bj,n,bi,bdiag);

1787:   PetscIncompleteLLDestroy(lnk,lnkbt);
1788:   PetscFreeSpaceDestroy(free_space_lvl);
1789:   PetscFree2(bj_ptr,bjlvl_ptr);

1791: #if defined(PETSC_USE_INFO)
1792:   {
1793:     PetscReal af = ((PetscReal)(bdiag[0]+1))/((PetscReal)ai[n]);
1794:     PetscInfo3(A,"Reallocs %D Fill ratio:given %g needed %g\n",reallocs,(double)f,(double)af);
1795:     PetscInfo1(A,"Run with -[sub_]pc_factor_fill %g or use \n",(double)af);
1796:     PetscInfo1(A,"PCFactorSetFill([sub]pc,%g);\n",(double)af);
1797:     PetscInfo(A,"for best performance.\n");
1798:     if (diagonal_fill) {
1799:       PetscInfo1(A,"Detected and replaced %D missing diagonals\n",dcount);
1800:     }
1801:   }
1802: #endif
1803:   /* put together the new matrix */
1804:   MatSeqAIJSetPreallocation_SeqAIJ(fact,MAT_SKIP_ALLOCATION,NULL);
1805:   PetscLogObjectParent((PetscObject)fact,(PetscObject)isicol);
1806:   b    = (Mat_SeqAIJ*)(fact)->data;

1808:   b->free_a       = PETSC_TRUE;
1809:   b->free_ij      = PETSC_TRUE;
1810:   b->singlemalloc = PETSC_FALSE;

1812:   PetscMalloc1(bdiag[0]+1,&b->a);

1814:   b->j    = bj;
1815:   b->i    = bi;
1816:   b->diag = bdiag;
1817:   b->ilen = 0;
1818:   b->imax = 0;
1819:   b->row  = isrow;
1820:   b->col  = iscol;
1821:   PetscObjectReference((PetscObject)isrow);
1822:   PetscObjectReference((PetscObject)iscol);
1823:   b->icol = isicol;

1825:   PetscMalloc1(n+1,&b->solve_work);
1826:   /* In b structure:  Free imax, ilen, old a, old j.
1827:      Allocate bdiag, solve_work, new a, new j */
1828:   PetscLogObjectMemory((PetscObject)fact,(bdiag[0]+1)*(sizeof(PetscInt)+sizeof(PetscScalar)));
1829:   b->maxnz = b->nz = bdiag[0]+1;

1831:   (fact)->info.factor_mallocs    = reallocs;
1832:   (fact)->info.fill_ratio_given  = f;
1833:   (fact)->info.fill_ratio_needed = ((PetscReal)(bdiag[0]+1))/((PetscReal)ai[n]);
1834:   (fact)->ops->lufactornumeric   = MatLUFactorNumeric_SeqAIJ;
1835:   if (a->inode.size) {
1836:     (fact)->ops->lufactornumeric = MatLUFactorNumeric_SeqAIJ_Inode;
1837:   }
1838:   MatSeqAIJCheckInode_FactorLU(fact);
1839:   return(0);
1840: }

1842: PetscErrorCode MatILUFactorSymbolic_SeqAIJ_inplace(Mat fact,Mat A,IS isrow,IS iscol,const MatFactorInfo *info)
1843: {
1844:   Mat_SeqAIJ         *a = (Mat_SeqAIJ*)A->data,*b;
1845:   IS                 isicol;
1846:   PetscErrorCode     ierr;
1847:   const PetscInt     *r,*ic;
1848:   PetscInt           n=A->rmap->n,*ai=a->i,*aj=a->j;
1849:   PetscInt           *bi,*cols,nnz,*cols_lvl;
1850:   PetscInt           *bdiag,prow,fm,nzbd,reallocs=0,dcount=0;
1851:   PetscInt           i,levels,diagonal_fill;
1852:   PetscBool          col_identity,row_identity;
1853:   PetscReal          f;
1854:   PetscInt           nlnk,*lnk,*lnk_lvl=NULL;
1855:   PetscBT            lnkbt;
1856:   PetscInt           nzi,*bj,**bj_ptr,**bjlvl_ptr;
1857:   PetscFreeSpaceList free_space    =NULL,current_space=NULL;
1858:   PetscFreeSpaceList free_space_lvl=NULL,current_space_lvl=NULL;
1859:   PetscBool          missing;

1862:   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);
1863:   MatMissingDiagonal(A,&missing,&i);
1864:   if (missing) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Matrix is missing diagonal entry %D",i);

1866:   f             = info->fill;
1867:   levels        = (PetscInt)info->levels;
1868:   diagonal_fill = (PetscInt)info->diagonal_fill;

1870:   ISInvertPermutation(iscol,PETSC_DECIDE,&isicol);

1872:   ISIdentity(isrow,&row_identity);
1873:   ISIdentity(iscol,&col_identity);
1874:   if (!levels && row_identity && col_identity) { /* special case: ilu(0) with natural ordering */
1875:     MatDuplicateNoCreate_SeqAIJ(fact,A,MAT_DO_NOT_COPY_VALUES,PETSC_TRUE);

1877:     (fact)->ops->lufactornumeric =  MatLUFactorNumeric_SeqAIJ_inplace;
1878:     if (a->inode.size) {
1879:       (fact)->ops->lufactornumeric = MatLUFactorNumeric_SeqAIJ_Inode_inplace;
1880:     }
1881:     fact->factortype               = MAT_FACTOR_ILU;
1882:     (fact)->info.factor_mallocs    = 0;
1883:     (fact)->info.fill_ratio_given  = info->fill;
1884:     (fact)->info.fill_ratio_needed = 1.0;

1886:     b    = (Mat_SeqAIJ*)(fact)->data;
1887:     b->row  = isrow;
1888:     b->col  = iscol;
1889:     b->icol = isicol;
1890:     PetscMalloc1((fact)->rmap->n+1,&b->solve_work);
1891:     PetscObjectReference((PetscObject)isrow);
1892:     PetscObjectReference((PetscObject)iscol);
1893:     return(0);
1894:   }

1896:   ISGetIndices(isrow,&r);
1897:   ISGetIndices(isicol,&ic);

1899:   /* get new row and diagonal pointers, must be allocated separately because they will be given to the Mat_SeqAIJ and freed separately */
1900:   PetscMalloc1(n+1,&bi);
1901:   PetscMalloc1(n+1,&bdiag);
1902:   bi[0] = bdiag[0] = 0;

1904:   PetscMalloc2(n,&bj_ptr,n,&bjlvl_ptr);

1906:   /* create a linked list for storing column indices of the active row */
1907:   nlnk = n + 1;
1908:   PetscIncompleteLLCreate(n,n,nlnk,lnk,lnk_lvl,lnkbt);

1910:   /* initial FreeSpace size is f*(ai[n]+1) */
1911:   PetscFreeSpaceGet(PetscRealIntMultTruncate(f,ai[n]+1),&free_space);
1912:   current_space     = free_space;
1913:   PetscFreeSpaceGet(PetscRealIntMultTruncate(f,ai[n]+1),&free_space_lvl);
1914:   current_space_lvl = free_space_lvl;

1916:   for (i=0; i<n; i++) {
1917:     nzi = 0;
1918:     /* copy current row into linked list */
1919:     nnz = ai[r[i]+1] - ai[r[i]];
1920:     if (!nnz) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_MAT_LU_ZRPVT,"Empty row in matrix: row in original ordering %D in permuted ordering %D",r[i],i);
1921:     cols   = aj + ai[r[i]];
1922:     lnk[i] = -1; /* marker to indicate if diagonal exists */
1923:     PetscIncompleteLLInit(nnz,cols,n,ic,nlnk,lnk,lnk_lvl,lnkbt);
1924:     nzi   += nlnk;

1926:     /* make sure diagonal entry is included */
1927:     if (diagonal_fill && lnk[i] == -1) {
1928:       fm = n;
1929:       while (lnk[fm] < i) fm = lnk[fm];
1930:       lnk[i]     = lnk[fm]; /* insert diagonal into linked list */
1931:       lnk[fm]    = i;
1932:       lnk_lvl[i] = 0;
1933:       nzi++; dcount++;
1934:     }

1936:     /* add pivot rows into the active row */
1937:     nzbd = 0;
1938:     prow = lnk[n];
1939:     while (prow < i) {
1940:       nnz      = bdiag[prow];
1941:       cols     = bj_ptr[prow] + nnz + 1;
1942:       cols_lvl = bjlvl_ptr[prow] + nnz + 1;
1943:       nnz      = bi[prow+1] - bi[prow] - nnz - 1;
1944:       PetscILULLAddSorted(nnz,cols,levels,cols_lvl,prow,nlnk,lnk,lnk_lvl,lnkbt,prow);
1945:       nzi     += nlnk;
1946:       prow     = lnk[prow];
1947:       nzbd++;
1948:     }
1949:     bdiag[i] = nzbd;
1950:     bi[i+1]  = bi[i] + nzi;

1952:     /* if free space is not available, make more free space */
1953:     if (current_space->local_remaining<nzi) {
1954:       nnz  = PetscIntMultTruncate(nzi,n - i); /* estimated and max additional space needed */
1955:       PetscFreeSpaceGet(nnz,&current_space);
1956:       PetscFreeSpaceGet(nnz,&current_space_lvl);
1957:       reallocs++;
1958:     }

1960:     /* copy data into free_space and free_space_lvl, then initialize lnk */
1961:     PetscIncompleteLLClean(n,n,nzi,lnk,lnk_lvl,current_space->array,current_space_lvl->array,lnkbt);
1962:     bj_ptr[i]    = current_space->array;
1963:     bjlvl_ptr[i] = current_space_lvl->array;

1965:     /* make sure the active row i has diagonal entry */
1966:     if (*(bj_ptr[i]+bdiag[i]) != i) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_MAT_LU_ZRPVT,"Row %D has missing diagonal in factored matrix\ntry running with -pc_factor_nonzeros_along_diagonal or -pc_factor_diagonal_fill",i);

1968:     current_space->array               += nzi;
1969:     current_space->local_used          += nzi;
1970:     current_space->local_remaining     -= nzi;
1971:     current_space_lvl->array           += nzi;
1972:     current_space_lvl->local_used      += nzi;
1973:     current_space_lvl->local_remaining -= nzi;
1974:   }

1976:   ISRestoreIndices(isrow,&r);
1977:   ISRestoreIndices(isicol,&ic);

1979:   /* destroy list of free space and other temporary arrays */
1980:   PetscMalloc1(bi[n]+1,&bj);
1981:   PetscFreeSpaceContiguous(&free_space,bj); /* copy free_space -> bj */
1982:   PetscIncompleteLLDestroy(lnk,lnkbt);
1983:   PetscFreeSpaceDestroy(free_space_lvl);
1984:   PetscFree2(bj_ptr,bjlvl_ptr);

1986: #if defined(PETSC_USE_INFO)
1987:   {
1988:     PetscReal af = ((PetscReal)bi[n])/((PetscReal)ai[n]);
1989:     PetscInfo3(A,"Reallocs %D Fill ratio:given %g needed %g\n",reallocs,(double)f,(double)af);
1990:     PetscInfo1(A,"Run with -[sub_]pc_factor_fill %g or use \n",(double)af);
1991:     PetscInfo1(A,"PCFactorSetFill([sub]pc,%g);\n",(double)af);
1992:     PetscInfo(A,"for best performance.\n");
1993:     if (diagonal_fill) {
1994:       PetscInfo1(A,"Detected and replaced %D missing diagonals\n",dcount);
1995:     }
1996:   }
1997: #endif

1999:   /* put together the new matrix */
2000:   MatSeqAIJSetPreallocation_SeqAIJ(fact,MAT_SKIP_ALLOCATION,NULL);
2001:   PetscLogObjectParent((PetscObject)fact,(PetscObject)isicol);
2002:   b    = (Mat_SeqAIJ*)(fact)->data;

2004:   b->free_a       = PETSC_TRUE;
2005:   b->free_ij      = PETSC_TRUE;
2006:   b->singlemalloc = PETSC_FALSE;

2008:   PetscMalloc1(bi[n],&b->a);
2009:   b->j = bj;
2010:   b->i = bi;
2011:   for (i=0; i<n; i++) bdiag[i] += bi[i];
2012:   b->diag = bdiag;
2013:   b->ilen = 0;
2014:   b->imax = 0;
2015:   b->row  = isrow;
2016:   b->col  = iscol;
2017:   PetscObjectReference((PetscObject)isrow);
2018:   PetscObjectReference((PetscObject)iscol);
2019:   b->icol = isicol;
2020:   PetscMalloc1(n+1,&b->solve_work);
2021:   /* In b structure:  Free imax, ilen, old a, old j.
2022:      Allocate bdiag, solve_work, new a, new j */
2023:   PetscLogObjectMemory((PetscObject)fact,(bi[n]-n) * (sizeof(PetscInt)+sizeof(PetscScalar)));
2024:   b->maxnz = b->nz = bi[n];

2026:   (fact)->info.factor_mallocs    = reallocs;
2027:   (fact)->info.fill_ratio_given  = f;
2028:   (fact)->info.fill_ratio_needed = ((PetscReal)bi[n])/((PetscReal)ai[n]);
2029:   (fact)->ops->lufactornumeric   =  MatLUFactorNumeric_SeqAIJ_inplace;
2030:   if (a->inode.size) {
2031:     (fact)->ops->lufactornumeric = MatLUFactorNumeric_SeqAIJ_Inode_inplace;
2032:   }
2033:   return(0);
2034: }

2036: PetscErrorCode MatCholeskyFactorNumeric_SeqAIJ(Mat B,Mat A,const MatFactorInfo *info)
2037: {
2038:   Mat            C = B;
2039:   Mat_SeqAIJ     *a=(Mat_SeqAIJ*)A->data;
2040:   Mat_SeqSBAIJ   *b=(Mat_SeqSBAIJ*)C->data;
2041:   IS             ip=b->row,iip = b->icol;
2043:   const PetscInt *rip,*riip;
2044:   PetscInt       i,j,mbs=A->rmap->n,*bi=b->i,*bj=b->j,*bdiag=b->diag,*bjtmp;
2045:   PetscInt       *ai=a->i,*aj=a->j;
2046:   PetscInt       k,jmin,jmax,*c2r,*il,col,nexti,ili,nz;
2047:   MatScalar      *rtmp,*ba=b->a,*bval,*aa=a->a,dk,uikdi;
2048:   PetscBool      perm_identity;
2049:   FactorShiftCtx sctx;
2050:   PetscReal      rs;
2051:   MatScalar      d,*v;

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

2057:   if (info->shifttype == (PetscReal)MAT_SHIFT_POSITIVE_DEFINITE) { /* set sctx.shift_top=max{rs} */
2058:     sctx.shift_top = info->zeropivot;
2059:     for (i=0; i<mbs; i++) {
2060:       /* calculate sum(|aij|)-RealPart(aii), amt of shift needed for this row */
2061:       d  = (aa)[a->diag[i]];
2062:       rs = -PetscAbsScalar(d) - PetscRealPart(d);
2063:       v  = aa+ai[i];
2064:       nz = ai[i+1] - ai[i];
2065:       for (j=0; j<nz; j++) rs += PetscAbsScalar(v[j]);
2066:       if (rs>sctx.shift_top) sctx.shift_top = rs;
2067:     }
2068:     sctx.shift_top *= 1.1;
2069:     sctx.nshift_max = 5;
2070:     sctx.shift_lo   = 0.;
2071:     sctx.shift_hi   = 1.;
2072:   }

2074:   ISGetIndices(ip,&rip);
2075:   ISGetIndices(iip,&riip);

2077:   /* allocate working arrays
2078:      c2r: linked list, keep track of pivot rows for a given column. c2r[col]: head of the list for a given col
2079:      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
2080:   */
2081:   PetscMalloc3(mbs,&rtmp,mbs,&il,mbs,&c2r);

2083:   do {
2084:     sctx.newshift = PETSC_FALSE;

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

2089:     for (k = 0; k<mbs; k++) {
2090:       /* zero rtmp */
2091:       nz    = bi[k+1] - bi[k];
2092:       bjtmp = bj + bi[k];
2093:       for (j=0; j<nz; j++) rtmp[bjtmp[j]] = 0.0;

2095:       /* load in initial unfactored row */
2096:       bval = ba + bi[k];
2097:       jmin = ai[rip[k]]; jmax = ai[rip[k]+1];
2098:       for (j = jmin; j < jmax; j++) {
2099:         col = riip[aj[j]];
2100:         if (col >= k) { /* only take upper triangular entry */
2101:           rtmp[col] = aa[j];
2102:           *bval++   = 0.0; /* for in-place factorization */
2103:         }
2104:       }
2105:       /* shift the diagonal of the matrix: ZeropivotApply() */
2106:       rtmp[k] += sctx.shift_amount;  /* shift the diagonal of the matrix */

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

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

2115:         /* compute multiplier, update diag(k) and U(i,k) */
2116:         ili     = il[i]; /* index of first nonzero element in U(i,k:bms-1) */
2117:         uikdi   = -ba[ili]*ba[bdiag[i]]; /* diagonal(k) */
2118:         dk     += uikdi*ba[ili]; /* update diag[k] */
2119:         ba[ili] = uikdi; /* -U(i,k) */

2121:         /* add multiple of row i to k-th row */
2122:         jmin = ili + 1; jmax = bi[i+1];
2123:         if (jmin < jmax) {
2124:           for (j=jmin; j<jmax; j++) rtmp[bj[j]] += uikdi*ba[j];
2125:           /* update il and c2r for row i */
2126:           il[i] = jmin;
2127:           j     = bj[jmin]; c2r[i] = c2r[j]; c2r[j] = i;
2128:         }
2129:         i = nexti;
2130:       }

2132:       /* copy data into U(k,:) */
2133:       rs   = 0.0;
2134:       jmin = bi[k]; jmax = bi[k+1]-1;
2135:       if (jmin < jmax) {
2136:         for (j=jmin; j<jmax; j++) {
2137:           col = bj[j]; ba[j] = rtmp[col]; rs += PetscAbsScalar(ba[j]);
2138:         }
2139:         /* add the k-th row into il and c2r */
2140:         il[k] = jmin;
2141:         i     = bj[jmin]; c2r[k] = c2r[i]; c2r[i] = k;
2142:       }

2144:       /* MatPivotCheck() */
2145:       sctx.rs = rs;
2146:       sctx.pv = dk;
2147:       MatPivotCheck(B,A,info,&sctx,i);
2148:       if (sctx.newshift) break;
2149:       dk = sctx.pv;

2151:       ba[bdiag[k]] = 1.0/dk; /* U(k,k) */
2152:     }
2153:   } while (sctx.newshift);

2155:   PetscFree3(rtmp,il,c2r);
2156:   ISRestoreIndices(ip,&rip);
2157:   ISRestoreIndices(iip,&riip);

2159:   ISIdentity(ip,&perm_identity);
2160:   if (perm_identity) {
2161:     B->ops->solve          = MatSolve_SeqSBAIJ_1_NaturalOrdering;
2162:     B->ops->solvetranspose = MatSolve_SeqSBAIJ_1_NaturalOrdering;
2163:     B->ops->forwardsolve   = MatForwardSolve_SeqSBAIJ_1_NaturalOrdering;
2164:     B->ops->backwardsolve  = MatBackwardSolve_SeqSBAIJ_1_NaturalOrdering;
2165:   } else {
2166:     B->ops->solve          = MatSolve_SeqSBAIJ_1;
2167:     B->ops->solvetranspose = MatSolve_SeqSBAIJ_1;
2168:     B->ops->forwardsolve   = MatForwardSolve_SeqSBAIJ_1;
2169:     B->ops->backwardsolve  = MatBackwardSolve_SeqSBAIJ_1;
2170:   }

2172:   C->assembled    = PETSC_TRUE;
2173:   C->preallocated = PETSC_TRUE;

2175:   PetscLogFlops(C->rmap->n);

2177:   /* MatPivotView() */
2178:   if (sctx.nshift) {
2179:     if (info->shifttype == (PetscReal)MAT_SHIFT_POSITIVE_DEFINITE) {
2180:       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);
2181:     } else if (info->shifttype == (PetscReal)MAT_SHIFT_NONZERO) {
2182:       PetscInfo2(A,"number of shift_nz tries %D, shift_amount %g\n",sctx.nshift,(double)sctx.shift_amount);
2183:     } else if (info->shifttype == (PetscReal)MAT_SHIFT_INBLOCKS) {
2184:       PetscInfo2(A,"number of shift_inblocks applied %D, each shift_amount %g\n",sctx.nshift,(double)info->shiftamount);
2185:     }
2186:   }
2187:   return(0);
2188: }

2190: PetscErrorCode MatCholeskyFactorNumeric_SeqAIJ_inplace(Mat B,Mat A,const MatFactorInfo *info)
2191: {
2192:   Mat            C = B;
2193:   Mat_SeqAIJ     *a=(Mat_SeqAIJ*)A->data;
2194:   Mat_SeqSBAIJ   *b=(Mat_SeqSBAIJ*)C->data;
2195:   IS             ip=b->row,iip = b->icol;
2197:   const PetscInt *rip,*riip;
2198:   PetscInt       i,j,mbs=A->rmap->n,*bi=b->i,*bj=b->j,*bcol,*bjtmp;
2199:   PetscInt       *ai=a->i,*aj=a->j;
2200:   PetscInt       k,jmin,jmax,*jl,*il,col,nexti,ili,nz;
2201:   MatScalar      *rtmp,*ba=b->a,*bval,*aa=a->a,dk,uikdi;
2202:   PetscBool      perm_identity;
2203:   FactorShiftCtx sctx;
2204:   PetscReal      rs;
2205:   MatScalar      d,*v;

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

2211:   if (info->shifttype == (PetscReal)MAT_SHIFT_POSITIVE_DEFINITE) { /* set sctx.shift_top=max{rs} */
2212:     sctx.shift_top = info->zeropivot;
2213:     for (i=0; i<mbs; i++) {
2214:       /* calculate sum(|aij|)-RealPart(aii), amt of shift needed for this row */
2215:       d  = (aa)[a->diag[i]];
2216:       rs = -PetscAbsScalar(d) - PetscRealPart(d);
2217:       v  = aa+ai[i];
2218:       nz = ai[i+1] - ai[i];
2219:       for (j=0; j<nz; j++) rs += PetscAbsScalar(v[j]);
2220:       if (rs>sctx.shift_top) sctx.shift_top = rs;
2221:     }
2222:     sctx.shift_top *= 1.1;
2223:     sctx.nshift_max = 5;
2224:     sctx.shift_lo   = 0.;
2225:     sctx.shift_hi   = 1.;
2226:   }

2228:   ISGetIndices(ip,&rip);
2229:   ISGetIndices(iip,&riip);

2231:   /* initialization */
2232:   PetscMalloc3(mbs,&rtmp,mbs,&il,mbs,&jl);

2234:   do {
2235:     sctx.newshift = PETSC_FALSE;

2237:     for (i=0; i<mbs; i++) jl[i] = mbs;
2238:     il[0] = 0;

2240:     for (k = 0; k<mbs; k++) {
2241:       /* zero rtmp */
2242:       nz    = bi[k+1] - bi[k];
2243:       bjtmp = bj + bi[k];
2244:       for (j=0; j<nz; j++) rtmp[bjtmp[j]] = 0.0;

2246:       bval = ba + bi[k];
2247:       /* initialize k-th row by the perm[k]-th row of A */
2248:       jmin = ai[rip[k]]; jmax = ai[rip[k]+1];
2249:       for (j = jmin; j < jmax; j++) {
2250:         col = riip[aj[j]];
2251:         if (col >= k) { /* only take upper triangular entry */
2252:           rtmp[col] = aa[j];
2253:           *bval++   = 0.0; /* for in-place factorization */
2254:         }
2255:       }
2256:       /* shift the diagonal of the matrix */
2257:       if (sctx.nshift) rtmp[k] += sctx.shift_amount;

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

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

2266:         /* compute multiplier, update diag(k) and U(i,k) */
2267:         ili     = il[i]; /* index of first nonzero element in U(i,k:bms-1) */
2268:         uikdi   = -ba[ili]*ba[bi[i]]; /* diagonal(k) */
2269:         dk     += uikdi*ba[ili];
2270:         ba[ili] = uikdi; /* -U(i,k) */

2272:         /* add multiple of row i to k-th row */
2273:         jmin = ili + 1; jmax = bi[i+1];
2274:         if (jmin < jmax) {
2275:           for (j=jmin; j<jmax; j++) rtmp[bj[j]] += uikdi*ba[j];
2276:           /* update il and jl for row i */
2277:           il[i] = jmin;
2278:           j     = bj[jmin]; jl[i] = jl[j]; jl[j] = i;
2279:         }
2280:         i = nexti;
2281:       }

2283:       /* shift the diagonals when zero pivot is detected */
2284:       /* compute rs=sum of abs(off-diagonal) */
2285:       rs   = 0.0;
2286:       jmin = bi[k]+1;
2287:       nz   = bi[k+1] - jmin;
2288:       bcol = bj + jmin;
2289:       for (j=0; j<nz; j++) {
2290:         rs += PetscAbsScalar(rtmp[bcol[j]]);
2291:       }

2293:       sctx.rs = rs;
2294:       sctx.pv = dk;
2295:       MatPivotCheck(B,A,info,&sctx,k);
2296:       if (sctx.newshift) break;
2297:       dk = sctx.pv;

2299:       /* copy data into U(k,:) */
2300:       ba[bi[k]] = 1.0/dk; /* U(k,k) */
2301:       jmin      = bi[k]+1; jmax = bi[k+1];
2302:       if (jmin < jmax) {
2303:         for (j=jmin; j<jmax; j++) {
2304:           col = bj[j]; ba[j] = rtmp[col];
2305:         }
2306:         /* add the k-th row into il and jl */
2307:         il[k] = jmin;
2308:         i     = bj[jmin]; jl[k] = jl[i]; jl[i] = k;
2309:       }
2310:     }
2311:   } while (sctx.newshift);

2313:   PetscFree3(rtmp,il,jl);
2314:   ISRestoreIndices(ip,&rip);
2315:   ISRestoreIndices(iip,&riip);

2317:   ISIdentity(ip,&perm_identity);
2318:   if (perm_identity) {
2319:     B->ops->solve          = MatSolve_SeqSBAIJ_1_NaturalOrdering_inplace;
2320:     B->ops->solvetranspose = MatSolve_SeqSBAIJ_1_NaturalOrdering_inplace;
2321:     B->ops->forwardsolve   = MatForwardSolve_SeqSBAIJ_1_NaturalOrdering_inplace;
2322:     B->ops->backwardsolve  = MatBackwardSolve_SeqSBAIJ_1_NaturalOrdering_inplace;
2323:   } else {
2324:     B->ops->solve          = MatSolve_SeqSBAIJ_1_inplace;
2325:     B->ops->solvetranspose = MatSolve_SeqSBAIJ_1_inplace;
2326:     B->ops->forwardsolve   = MatForwardSolve_SeqSBAIJ_1_inplace;
2327:     B->ops->backwardsolve  = MatBackwardSolve_SeqSBAIJ_1_inplace;
2328:   }

2330:   C->assembled    = PETSC_TRUE;
2331:   C->preallocated = PETSC_TRUE;

2333:   PetscLogFlops(C->rmap->n);
2334:   if (sctx.nshift) {
2335:     if (info->shifttype == (PetscReal)MAT_SHIFT_NONZERO) {
2336:       PetscInfo2(A,"number of shiftnz tries %D, shift_amount %g\n",sctx.nshift,(double)sctx.shift_amount);
2337:     } else if (info->shifttype == (PetscReal)MAT_SHIFT_POSITIVE_DEFINITE) {
2338:       PetscInfo2(A,"number of shiftpd tries %D, shift_amount %g\n",sctx.nshift,(double)sctx.shift_amount);
2339:     }
2340:   }
2341:   return(0);
2342: }

2344: /*
2345:    icc() under revised new data structure.
2346:    Factored arrays bj and ba are stored as
2347:      U(0,:),...,U(i,:),U(n-1,:)

2349:    ui=fact->i is an array of size n+1, in which
2350:    ui+
2351:      ui[i]:  points to 1st entry of U(i,:),i=0,...,n-1
2352:      ui[n]:  points to U(n-1,n-1)+1

2354:   udiag=fact->diag is an array of size n,in which
2355:      udiag[i]: points to diagonal of U(i,:), i=0,...,n-1

2357:    U(i,:) contains udiag[i] as its last entry, i.e.,
2358:     U(i,:) = (u[i,i+1],...,u[i,n-1],diag[i])
2359: */

2361: PetscErrorCode MatICCFactorSymbolic_SeqAIJ(Mat fact,Mat A,IS perm,const MatFactorInfo *info)
2362: {
2363:   Mat_SeqAIJ         *a = (Mat_SeqAIJ*)A->data;
2364:   Mat_SeqSBAIJ       *b;
2365:   PetscErrorCode     ierr;
2366:   PetscBool          perm_identity,missing;
2367:   PetscInt           reallocs=0,i,*ai=a->i,*aj=a->j,am=A->rmap->n,*ui,*udiag;
2368:   const PetscInt     *rip,*riip;
2369:   PetscInt           jmin,jmax,nzk,k,j,*jl,prow,*il,nextprow;
2370:   PetscInt           nlnk,*lnk,*lnk_lvl=NULL,d;
2371:   PetscInt           ncols,ncols_upper,*cols,*ajtmp,*uj,**uj_ptr,**uj_lvl_ptr;
2372:   PetscReal          fill          =info->fill,levels=info->levels;
2373:   PetscFreeSpaceList free_space    =NULL,current_space=NULL;
2374:   PetscFreeSpaceList free_space_lvl=NULL,current_space_lvl=NULL;
2375:   PetscBT            lnkbt;
2376:   IS                 iperm;

2379:   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);
2380:   MatMissingDiagonal(A,&missing,&d);
2381:   if (missing) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Matrix is missing diagonal entry %D",d);
2382:   ISIdentity(perm,&perm_identity);
2383:   ISInvertPermutation(perm,PETSC_DECIDE,&iperm);

2385:   PetscMalloc1(am+1,&ui);
2386:   PetscMalloc1(am+1,&udiag);
2387:   ui[0] = 0;

2389:   /* ICC(0) without matrix ordering: simply rearrange column indices */
2390:   if (!levels && perm_identity) {
2391:     for (i=0; i<am; i++) {
2392:       ncols    = ai[i+1] - a->diag[i];
2393:       ui[i+1]  = ui[i] + ncols;
2394:       udiag[i] = ui[i+1] - 1; /* points to the last entry of U(i,:) */
2395:     }
2396:     PetscMalloc1(ui[am]+1,&uj);
2397:     cols = uj;
2398:     for (i=0; i<am; i++) {
2399:       aj    = a->j + a->diag[i] + 1; /* 1st entry of U(i,:) without diagonal */
2400:       ncols = ai[i+1] - a->diag[i] -1;
2401:       for (j=0; j<ncols; j++) *cols++ = aj[j];
2402:       *cols++ = i; /* diagoanl is located as the last entry of U(i,:) */
2403:     }
2404:   } else { /* case: levels>0 || (levels=0 && !perm_identity) */
2405:     ISGetIndices(iperm,&riip);
2406:     ISGetIndices(perm,&rip);

2408:     /* initialization */
2409:     PetscMalloc1(am+1,&ajtmp);

2411:     /* jl: linked list for storing indices of the pivot rows
2412:        il: il[i] points to the 1st nonzero entry of U(i,k:am-1) */
2413:     PetscMalloc4(am,&uj_ptr,am,&uj_lvl_ptr,am,&jl,am,&il);
2414:     for (i=0; i<am; i++) {
2415:       jl[i] = am; il[i] = 0;
2416:     }

2418:     /* create and initialize a linked list for storing column indices of the active row k */
2419:     nlnk = am + 1;
2420:     PetscIncompleteLLCreate(am,am,nlnk,lnk,lnk_lvl,lnkbt);

2422:     /* initial FreeSpace size is fill*(ai[am]+am)/2 */
2423:     PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,(ai[am]+am)/2),&free_space);
2424:     current_space     = free_space;
2425:     PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,(ai[am]+am)/2),&free_space_lvl);
2426:     current_space_lvl = free_space_lvl;

2428:     for (k=0; k<am; k++) {  /* for each active row k */
2429:       /* initialize lnk by the column indices of row rip[k] of A */
2430:       nzk   = 0;
2431:       ncols = ai[rip[k]+1] - ai[rip[k]];
2432:       if (!ncols) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_MAT_CH_ZRPVT,"Empty row in matrix: row in original ordering %D in permuted ordering %D",rip[k],k);
2433:       ncols_upper = 0;
2434:       for (j=0; j<ncols; j++) {
2435:         i = *(aj + ai[rip[k]] + j); /* unpermuted column index */
2436:         if (riip[i] >= k) { /* only take upper triangular entry */
2437:           ajtmp[ncols_upper] = i;
2438:           ncols_upper++;
2439:         }
2440:       }
2441:       PetscIncompleteLLInit(ncols_upper,ajtmp,am,riip,nlnk,lnk,lnk_lvl,lnkbt);
2442:       nzk += nlnk;

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

2447:       while (prow < k) {
2448:         nextprow = jl[prow];

2450:         /* merge prow into k-th row */
2451:         jmin  = il[prow] + 1; /* index of the 2nd nzero entry in U(prow,k:am-1) */
2452:         jmax  = ui[prow+1];
2453:         ncols = jmax-jmin;
2454:         i     = jmin - ui[prow];
2455:         cols  = uj_ptr[prow] + i; /* points to the 2nd nzero entry in U(prow,k:am-1) */
2456:         uj    = uj_lvl_ptr[prow] + i; /* levels of cols */
2457:         j     = *(uj - 1);
2458:         PetscICCLLAddSorted(ncols,cols,levels,uj,am,nlnk,lnk,lnk_lvl,lnkbt,j);
2459:         nzk  += nlnk;

2461:         /* update il and jl for prow */
2462:         if (jmin < jmax) {
2463:           il[prow] = jmin;
2464:           j        = *cols; jl[prow] = jl[j]; jl[j] = prow;
2465:         }
2466:         prow = nextprow;
2467:       }

2469:       /* if free space is not available, make more free space */
2470:       if (current_space->local_remaining<nzk) {
2471:         i    = am - k + 1; /* num of unfactored rows */
2472:         i    = PetscIntMultTruncate(i,PetscMin(nzk, i-1)); /* i*nzk, i*(i-1): estimated and max additional space needed */
2473:         PetscFreeSpaceGet(i,&current_space);
2474:         PetscFreeSpaceGet(i,&current_space_lvl);
2475:         reallocs++;
2476:       }

2478:       /* copy data into free_space and free_space_lvl, then initialize lnk */
2479:       if (nzk == 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Empty row %D in ICC matrix factor",k);
2480:       PetscIncompleteLLClean(am,am,nzk,lnk,lnk_lvl,current_space->array,current_space_lvl->array,lnkbt);

2482:       /* add the k-th row into il and jl */
2483:       if (nzk > 1) {
2484:         i     = current_space->array[1]; /* col value of the first nonzero element in U(k, k+1:am-1) */
2485:         jl[k] = jl[i]; jl[i] = k;
2486:         il[k] = ui[k] + 1;
2487:       }
2488:       uj_ptr[k]     = current_space->array;
2489:       uj_lvl_ptr[k] = current_space_lvl->array;

2491:       current_space->array           += nzk;
2492:       current_space->local_used      += nzk;
2493:       current_space->local_remaining -= nzk;

2495:       current_space_lvl->array           += nzk;
2496:       current_space_lvl->local_used      += nzk;
2497:       current_space_lvl->local_remaining -= nzk;

2499:       ui[k+1] = ui[k] + nzk;
2500:     }

2502:     ISRestoreIndices(perm,&rip);
2503:     ISRestoreIndices(iperm,&riip);
2504:     PetscFree4(uj_ptr,uj_lvl_ptr,jl,il);
2505:     PetscFree(ajtmp);

2507:     /* copy free_space into uj and free free_space; set ui, uj, udiag in new datastructure; */
2508:     PetscMalloc1(ui[am]+1,&uj);
2509:     PetscFreeSpaceContiguous_Cholesky(&free_space,uj,am,ui,udiag); /* store matrix factor  */
2510:     PetscIncompleteLLDestroy(lnk,lnkbt);
2511:     PetscFreeSpaceDestroy(free_space_lvl);

2513:   } /* end of case: levels>0 || (levels=0 && !perm_identity) */

2515:   /* put together the new matrix in MATSEQSBAIJ format */
2516:   b               = (Mat_SeqSBAIJ*)(fact)->data;
2517:   b->singlemalloc = PETSC_FALSE;

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

2521:   b->j             = uj;
2522:   b->i             = ui;
2523:   b->diag          = udiag;
2524:   b->free_diag     = PETSC_TRUE;
2525:   b->ilen          = 0;
2526:   b->imax          = 0;
2527:   b->row           = perm;
2528:   b->col           = perm;
2529:   PetscObjectReference((PetscObject)perm);
2530:   PetscObjectReference((PetscObject)perm);
2531:   b->icol          = iperm;
2532:   b->pivotinblocks = PETSC_FALSE; /* need to get from MatFactorInfo */

2534:   PetscMalloc1(am+1,&b->solve_work);
2535:   PetscLogObjectMemory((PetscObject)fact,ui[am]*(sizeof(PetscInt)+sizeof(MatScalar)));

2537:   b->maxnz   = b->nz = ui[am];
2538:   b->free_a  = PETSC_TRUE;
2539:   b->free_ij = PETSC_TRUE;

2541:   fact->info.factor_mallocs   = reallocs;
2542:   fact->info.fill_ratio_given = fill;
2543:   if (ai[am] != 0) {
2544:     /* nonzeros in lower triangular part of A (including diagonals) = (ai[am]+am)/2 */
2545:     fact->info.fill_ratio_needed = ((PetscReal)2*ui[am])/(ai[am]+am);
2546:   } else {
2547:     fact->info.fill_ratio_needed = 0.0;
2548:   }
2549: #if defined(PETSC_USE_INFO)
2550:   if (ai[am] != 0) {
2551:     PetscReal af = fact->info.fill_ratio_needed;
2552:     PetscInfo3(A,"Reallocs %D Fill ratio:given %g needed %g\n",reallocs,(double)fill,(double)af);
2553:     PetscInfo1(A,"Run with -pc_factor_fill %g or use \n",(double)af);
2554:     PetscInfo1(A,"PCFactorSetFill(pc,%g) for best performance.\n",(double)af);
2555:   } else {
2556:     PetscInfo(A,"Empty matrix.\n");
2557:   }
2558: #endif
2559:   fact->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqAIJ;
2560:   return(0);
2561: }

2563: PetscErrorCode MatICCFactorSymbolic_SeqAIJ_inplace(Mat fact,Mat A,IS perm,const MatFactorInfo *info)
2564: {
2565:   Mat_SeqAIJ         *a = (Mat_SeqAIJ*)A->data;
2566:   Mat_SeqSBAIJ       *b;
2567:   PetscErrorCode     ierr;
2568:   PetscBool          perm_identity,missing;
2569:   PetscInt           reallocs=0,i,*ai=a->i,*aj=a->j,am=A->rmap->n,*ui,*udiag;
2570:   const PetscInt     *rip,*riip;
2571:   PetscInt           jmin,jmax,nzk,k,j,*jl,prow,*il,nextprow;
2572:   PetscInt           nlnk,*lnk,*lnk_lvl=NULL,d;
2573:   PetscInt           ncols,ncols_upper,*cols,*ajtmp,*uj,**uj_ptr,**uj_lvl_ptr;
2574:   PetscReal          fill          =info->fill,levels=info->levels;
2575:   PetscFreeSpaceList free_space    =NULL,current_space=NULL;
2576:   PetscFreeSpaceList free_space_lvl=NULL,current_space_lvl=NULL;
2577:   PetscBT            lnkbt;
2578:   IS                 iperm;

2581:   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);
2582:   MatMissingDiagonal(A,&missing,&d);
2583:   if (missing) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Matrix is missing diagonal entry %D",d);
2584:   ISIdentity(perm,&perm_identity);
2585:   ISInvertPermutation(perm,PETSC_DECIDE,&iperm);

2587:   PetscMalloc1(am+1,&ui);
2588:   PetscMalloc1(am+1,&udiag);
2589:   ui[0] = 0;

2591:   /* ICC(0) without matrix ordering: simply copies fill pattern */
2592:   if (!levels && perm_identity) {

2594:     for (i=0; i<am; i++) {
2595:       ui[i+1]  = ui[i] + ai[i+1] - a->diag[i];
2596:       udiag[i] = ui[i];
2597:     }
2598:     PetscMalloc1(ui[am]+1,&uj);
2599:     cols = uj;
2600:     for (i=0; i<am; i++) {
2601:       aj    = a->j + a->diag[i];
2602:       ncols = ui[i+1] - ui[i];
2603:       for (j=0; j<ncols; j++) *cols++ = *aj++;
2604:     }
2605:   } else { /* case: levels>0 || (levels=0 && !perm_identity) */
2606:     ISGetIndices(iperm,&riip);
2607:     ISGetIndices(perm,&rip);

2609:     /* initialization */
2610:     PetscMalloc1(am+1,&ajtmp);

2612:     /* jl: linked list for storing indices of the pivot rows
2613:        il: il[i] points to the 1st nonzero entry of U(i,k:am-1) */
2614:     PetscMalloc4(am,&uj_ptr,am,&uj_lvl_ptr,am,&jl,am,&il);
2615:     for (i=0; i<am; i++) {
2616:       jl[i] = am; il[i] = 0;
2617:     }

2619:     /* create and initialize a linked list for storing column indices of the active row k */
2620:     nlnk = am + 1;
2621:     PetscIncompleteLLCreate(am,am,nlnk,lnk,lnk_lvl,lnkbt);

2623:     /* initial FreeSpace size is fill*(ai[am]+1) */
2624:     PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,ai[am]+1),&free_space);
2625:     current_space     = free_space;
2626:     PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,ai[am]+1),&free_space_lvl);
2627:     current_space_lvl = free_space_lvl;

2629:     for (k=0; k<am; k++) {  /* for each active row k */
2630:       /* initialize lnk by the column indices of row rip[k] of A */
2631:       nzk   = 0;
2632:       ncols = ai[rip[k]+1] - ai[rip[k]];
2633:       if (!ncols) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_MAT_CH_ZRPVT,"Empty row in matrix: row in original ordering %D in permuted ordering %D",rip[k],k);
2634:       ncols_upper = 0;
2635:       for (j=0; j<ncols; j++) {
2636:         i = *(aj + ai[rip[k]] + j); /* unpermuted column index */
2637:         if (riip[i] >= k) { /* only take upper triangular entry */
2638:           ajtmp[ncols_upper] = i;
2639:           ncols_upper++;
2640:         }
2641:       }
2642:       PetscIncompleteLLInit(ncols_upper,ajtmp,am,riip,nlnk,lnk,lnk_lvl,lnkbt);
2643:       nzk += nlnk;

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

2648:       while (prow < k) {
2649:         nextprow = jl[prow];

2651:         /* merge prow into k-th row */
2652:         jmin  = il[prow] + 1; /* index of the 2nd nzero entry in U(prow,k:am-1) */
2653:         jmax  = ui[prow+1];
2654:         ncols = jmax-jmin;
2655:         i     = jmin - ui[prow];
2656:         cols  = uj_ptr[prow] + i; /* points to the 2nd nzero entry in U(prow,k:am-1) */
2657:         uj    = uj_lvl_ptr[prow] + i; /* levels of cols */
2658:         j     = *(uj - 1);
2659:         PetscICCLLAddSorted(ncols,cols,levels,uj,am,nlnk,lnk,lnk_lvl,lnkbt,j);
2660:         nzk  += nlnk;

2662:         /* update il and jl for prow */
2663:         if (jmin < jmax) {
2664:           il[prow] = jmin;
2665:           j        = *cols; jl[prow] = jl[j]; jl[j] = prow;
2666:         }
2667:         prow = nextprow;
2668:       }

2670:       /* if free space is not available, make more free space */
2671:       if (current_space->local_remaining<nzk) {
2672:         i    = am - k + 1; /* num of unfactored rows */
2673:         i    = PetscIntMultTruncate(i,PetscMin(nzk, i-1)); /* i*nzk, i*(i-1): estimated and max additional space needed */
2674:         PetscFreeSpaceGet(i,&current_space);
2675:         PetscFreeSpaceGet(i,&current_space_lvl);
2676:         reallocs++;
2677:       }

2679:       /* copy data into free_space and free_space_lvl, then initialize lnk */
2680:       if (!nzk) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Empty row %D in ICC matrix factor",k);
2681:       PetscIncompleteLLClean(am,am,nzk,lnk,lnk_lvl,current_space->array,current_space_lvl->array,lnkbt);

2683:       /* add the k-th row into il and jl */
2684:       if (nzk > 1) {
2685:         i     = current_space->array[1]; /* col value of the first nonzero element in U(k, k+1:am-1) */
2686:         jl[k] = jl[i]; jl[i] = k;
2687:         il[k] = ui[k] + 1;
2688:       }
2689:       uj_ptr[k]     = current_space->array;
2690:       uj_lvl_ptr[k] = current_space_lvl->array;

2692:       current_space->array           += nzk;
2693:       current_space->local_used      += nzk;
2694:       current_space->local_remaining -= nzk;

2696:       current_space_lvl->array           += nzk;
2697:       current_space_lvl->local_used      += nzk;
2698:       current_space_lvl->local_remaining -= nzk;

2700:       ui[k+1] = ui[k] + nzk;
2701:     }

2703: #if defined(PETSC_USE_INFO)
2704:     if (ai[am] != 0) {
2705:       PetscReal af = (PetscReal)ui[am]/((PetscReal)ai[am]);
2706:       PetscInfo3(A,"Reallocs %D Fill ratio:given %g needed %g\n",reallocs,(double)fill,(double)af);
2707:       PetscInfo1(A,"Run with -pc_factor_fill %g or use \n",(double)af);
2708:       PetscInfo1(A,"PCFactorSetFill(pc,%g) for best performance.\n",(double)af);
2709:     } else {
2710:       PetscInfo(A,"Empty matrix.\n");
2711:     }
2712: #endif

2714:     ISRestoreIndices(perm,&rip);
2715:     ISRestoreIndices(iperm,&riip);
2716:     PetscFree4(uj_ptr,uj_lvl_ptr,jl,il);
2717:     PetscFree(ajtmp);

2719:     /* destroy list of free space and other temporary array(s) */
2720:     PetscMalloc1(ui[am]+1,&uj);
2721:     PetscFreeSpaceContiguous(&free_space,uj);
2722:     PetscIncompleteLLDestroy(lnk,lnkbt);
2723:     PetscFreeSpaceDestroy(free_space_lvl);

2725:   } /* end of case: levels>0 || (levels=0 && !perm_identity) */

2727:   /* put together the new matrix in MATSEQSBAIJ format */

2729:   b               = (Mat_SeqSBAIJ*)fact->data;
2730:   b->singlemalloc = PETSC_FALSE;

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

2734:   b->j         = uj;
2735:   b->i         = ui;
2736:   b->diag      = udiag;
2737:   b->free_diag = PETSC_TRUE;
2738:   b->ilen      = 0;
2739:   b->imax      = 0;
2740:   b->row       = perm;
2741:   b->col       = perm;

2743:   PetscObjectReference((PetscObject)perm);
2744:   PetscObjectReference((PetscObject)perm);

2746:   b->icol          = iperm;
2747:   b->pivotinblocks = PETSC_FALSE; /* need to get from MatFactorInfo */
2748:   PetscMalloc1(am+1,&b->solve_work);
2749:   PetscLogObjectMemory((PetscObject)fact,(ui[am]-am)*(sizeof(PetscInt)+sizeof(MatScalar)));
2750:   b->maxnz         = b->nz = ui[am];
2751:   b->free_a        = PETSC_TRUE;
2752:   b->free_ij       = PETSC_TRUE;

2754:   fact->info.factor_mallocs   = reallocs;
2755:   fact->info.fill_ratio_given = fill;
2756:   if (ai[am] != 0) {
2757:     fact->info.fill_ratio_needed = ((PetscReal)ui[am])/((PetscReal)ai[am]);
2758:   } else {
2759:     fact->info.fill_ratio_needed = 0.0;
2760:   }
2761:   fact->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqAIJ_inplace;
2762:   return(0);
2763: }

2765: PetscErrorCode MatCholeskyFactorSymbolic_SeqAIJ(Mat fact,Mat A,IS perm,const MatFactorInfo *info)
2766: {
2767:   Mat_SeqAIJ         *a = (Mat_SeqAIJ*)A->data;
2768:   Mat_SeqSBAIJ       *b;
2769:   PetscErrorCode     ierr;
2770:   PetscBool          perm_identity,missing;
2771:   PetscReal          fill = info->fill;
2772:   const PetscInt     *rip,*riip;
2773:   PetscInt           i,am=A->rmap->n,*ai=a->i,*aj=a->j,reallocs=0,prow;
2774:   PetscInt           *jl,jmin,jmax,nzk,*ui,k,j,*il,nextprow;
2775:   PetscInt           nlnk,*lnk,ncols,ncols_upper,*cols,*uj,**ui_ptr,*uj_ptr,*udiag;
2776:   PetscFreeSpaceList free_space=NULL,current_space=NULL;
2777:   PetscBT            lnkbt;
2778:   IS                 iperm;

2781:   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);
2782:   MatMissingDiagonal(A,&missing,&i);
2783:   if (missing) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Matrix is missing diagonal entry %D",i);

2785:   /* check whether perm is the identity mapping */
2786:   ISIdentity(perm,&perm_identity);
2787:   ISInvertPermutation(perm,PETSC_DECIDE,&iperm);
2788:   ISGetIndices(iperm,&riip);
2789:   ISGetIndices(perm,&rip);

2791:   /* initialization */
2792:   PetscMalloc1(am+1,&ui);
2793:   PetscMalloc1(am+1,&udiag);
2794:   ui[0] = 0;

2796:   /* jl: linked list for storing indices of the pivot rows
2797:      il: il[i] points to the 1st nonzero entry of U(i,k:am-1) */
2798:   PetscMalloc4(am,&ui_ptr,am,&jl,am,&il,am,&cols);
2799:   for (i=0; i<am; i++) {
2800:     jl[i] = am; il[i] = 0;
2801:   }

2803:   /* create and initialize a linked list for storing column indices of the active row k */
2804:   nlnk = am + 1;
2805:   PetscLLCreate(am,am,nlnk,lnk,lnkbt);

2807:   /* initial FreeSpace size is fill*(ai[am]+am)/2 */
2808:   PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,(ai[am]+am)/2),&free_space);
2809:   current_space = free_space;

2811:   for (k=0; k<am; k++) {  /* for each active row k */
2812:     /* initialize lnk by the column indices of row rip[k] of A */
2813:     nzk   = 0;
2814:     ncols = ai[rip[k]+1] - ai[rip[k]];
2815:     if (!ncols) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_MAT_CH_ZRPVT,"Empty row in matrix: row in original ordering %D in permuted ordering %D",rip[k],k);
2816:     ncols_upper = 0;
2817:     for (j=0; j<ncols; j++) {
2818:       i = riip[*(aj + ai[rip[k]] + j)];
2819:       if (i >= k) { /* only take upper triangular entry */
2820:         cols[ncols_upper] = i;
2821:         ncols_upper++;
2822:       }
2823:     }
2824:     PetscLLAdd(ncols_upper,cols,am,nlnk,lnk,lnkbt);
2825:     nzk += nlnk;

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

2830:     while (prow < k) {
2831:       nextprow = jl[prow];
2832:       /* merge prow into k-th row */
2833:       jmin   = il[prow] + 1; /* index of the 2nd nzero entry in U(prow,k:am-1) */
2834:       jmax   = ui[prow+1];
2835:       ncols  = jmax-jmin;
2836:       uj_ptr = ui_ptr[prow] + jmin - ui[prow]; /* points to the 2nd nzero entry in U(prow,k:am-1) */
2837:       PetscLLAddSorted(ncols,uj_ptr,am,nlnk,lnk,lnkbt);
2838:       nzk   += nlnk;

2840:       /* update il and jl for prow */
2841:       if (jmin < jmax) {
2842:         il[prow] = jmin;
2843:         j        = *uj_ptr;
2844:         jl[prow] = jl[j];
2845:         jl[j]    = prow;
2846:       }
2847:       prow = nextprow;
2848:     }

2850:     /* if free space is not available, make more free space */
2851:     if (current_space->local_remaining<nzk) {
2852:       i    = am - k + 1; /* num of unfactored rows */
2853:       i    = PetscIntMultTruncate(i,PetscMin(nzk,i-1)); /* i*nzk, i*(i-1): estimated and max additional space needed */
2854:       PetscFreeSpaceGet(i,&current_space);
2855:       reallocs++;
2856:     }

2858:     /* copy data into free space, then initialize lnk */
2859:     PetscLLClean(am,am,nzk,lnk,current_space->array,lnkbt);

2861:     /* add the k-th row into il and jl */
2862:     if (nzk > 1) {
2863:       i     = current_space->array[1]; /* col value of the first nonzero element in U(k, k+1:am-1) */
2864:       jl[k] = jl[i]; jl[i] = k;
2865:       il[k] = ui[k] + 1;
2866:     }
2867:     ui_ptr[k] = current_space->array;

2869:     current_space->array           += nzk;
2870:     current_space->local_used      += nzk;
2871:     current_space->local_remaining -= nzk;

2873:     ui[k+1] = ui[k] + nzk;
2874:   }

2876:   ISRestoreIndices(perm,&rip);
2877:   ISRestoreIndices(iperm,&riip);
2878:   PetscFree4(ui_ptr,jl,il,cols);

2880:   /* copy free_space into uj and free free_space; set ui, uj, udiag in new datastructure; */
2881:   PetscMalloc1(ui[am]+1,&uj);
2882:   PetscFreeSpaceContiguous_Cholesky(&free_space,uj,am,ui,udiag); /* store matrix factor */
2883:   PetscLLDestroy(lnk,lnkbt);

2885:   /* put together the new matrix in MATSEQSBAIJ format */

2887:   b               = (Mat_SeqSBAIJ*)fact->data;
2888:   b->singlemalloc = PETSC_FALSE;
2889:   b->free_a       = PETSC_TRUE;
2890:   b->free_ij      = PETSC_TRUE;

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

2894:   b->j         = uj;
2895:   b->i         = ui;
2896:   b->diag      = udiag;
2897:   b->free_diag = PETSC_TRUE;
2898:   b->ilen      = 0;
2899:   b->imax      = 0;
2900:   b->row       = perm;
2901:   b->col       = perm;

2903:   PetscObjectReference((PetscObject)perm);
2904:   PetscObjectReference((PetscObject)perm);

2906:   b->icol          = iperm;
2907:   b->pivotinblocks = PETSC_FALSE; /* need to get from MatFactorInfo */

2909:   PetscMalloc1(am+1,&b->solve_work);
2910:   PetscLogObjectMemory((PetscObject)fact,ui[am]*(sizeof(PetscInt)+sizeof(MatScalar)));

2912:   b->maxnz = b->nz = ui[am];

2914:   fact->info.factor_mallocs   = reallocs;
2915:   fact->info.fill_ratio_given = fill;
2916:   if (ai[am] != 0) {
2917:     /* nonzeros in lower triangular part of A (including diagonals) = (ai[am]+am)/2 */
2918:     fact->info.fill_ratio_needed = ((PetscReal)2*ui[am])/(ai[am]+am);
2919:   } else {
2920:     fact->info.fill_ratio_needed = 0.0;
2921:   }
2922: #if defined(PETSC_USE_INFO)
2923:   if (ai[am] != 0) {
2924:     PetscReal af = fact->info.fill_ratio_needed;
2925:     PetscInfo3(A,"Reallocs %D Fill ratio:given %g needed %g\n",reallocs,(double)fill,(double)af);
2926:     PetscInfo1(A,"Run with -pc_factor_fill %g or use \n",(double)af);
2927:     PetscInfo1(A,"PCFactorSetFill(pc,%g) for best performance.\n",(double)af);
2928:   } else {
2929:     PetscInfo(A,"Empty matrix.\n");
2930:   }
2931: #endif
2932:   fact->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqAIJ;
2933:   return(0);
2934: }

2936: PetscErrorCode MatCholeskyFactorSymbolic_SeqAIJ_inplace(Mat fact,Mat A,IS perm,const MatFactorInfo *info)
2937: {
2938:   Mat_SeqAIJ         *a = (Mat_SeqAIJ*)A->data;
2939:   Mat_SeqSBAIJ       *b;
2940:   PetscErrorCode     ierr;
2941:   PetscBool          perm_identity,missing;
2942:   PetscReal          fill = info->fill;
2943:   const PetscInt     *rip,*riip;
2944:   PetscInt           i,am=A->rmap->n,*ai=a->i,*aj=a->j,reallocs=0,prow;
2945:   PetscInt           *jl,jmin,jmax,nzk,*ui,k,j,*il,nextprow;
2946:   PetscInt           nlnk,*lnk,ncols,ncols_upper,*cols,*uj,**ui_ptr,*uj_ptr;
2947:   PetscFreeSpaceList free_space=NULL,current_space=NULL;
2948:   PetscBT            lnkbt;
2949:   IS                 iperm;

2952:   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);
2953:   MatMissingDiagonal(A,&missing,&i);
2954:   if (missing) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Matrix is missing diagonal entry %D",i);

2956:   /* check whether perm is the identity mapping */
2957:   ISIdentity(perm,&perm_identity);
2958:   ISInvertPermutation(perm,PETSC_DECIDE,&iperm);
2959:   ISGetIndices(iperm,&riip);
2960:   ISGetIndices(perm,&rip);

2962:   /* initialization */
2963:   PetscMalloc1(am+1,&ui);
2964:   ui[0] = 0;

2966:   /* jl: linked list for storing indices of the pivot rows
2967:      il: il[i] points to the 1st nonzero entry of U(i,k:am-1) */
2968:   PetscMalloc4(am,&ui_ptr,am,&jl,am,&il,am,&cols);
2969:   for (i=0; i<am; i++) {
2970:     jl[i] = am; il[i] = 0;
2971:   }

2973:   /* create and initialize a linked list for storing column indices of the active row k */
2974:   nlnk = am + 1;
2975:   PetscLLCreate(am,am,nlnk,lnk,lnkbt);

2977:   /* initial FreeSpace size is fill*(ai[am]+1) */
2978:   PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,ai[am]+1),&free_space);
2979:   current_space = free_space;

2981:   for (k=0; k<am; k++) {  /* for each active row k */
2982:     /* initialize lnk by the column indices of row rip[k] of A */
2983:     nzk   = 0;
2984:     ncols = ai[rip[k]+1] - ai[rip[k]];
2985:     if (!ncols) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_MAT_CH_ZRPVT,"Empty row in matrix: row in original ordering %D in permuted ordering %D",rip[k],k);
2986:     ncols_upper = 0;
2987:     for (j=0; j<ncols; j++) {
2988:       i = riip[*(aj + ai[rip[k]] + j)];
2989:       if (i >= k) { /* only take upper triangular entry */
2990:         cols[ncols_upper] = i;
2991:         ncols_upper++;
2992:       }
2993:     }
2994:     PetscLLAdd(ncols_upper,cols,am,nlnk,lnk,lnkbt);
2995:     nzk += nlnk;

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

3000:     while (prow < k) {
3001:       nextprow = jl[prow];
3002:       /* merge prow into k-th row */
3003:       jmin   = il[prow] + 1; /* index of the 2nd nzero entry in U(prow,k:am-1) */
3004:       jmax   = ui[prow+1];
3005:       ncols  = jmax-jmin;
3006:       uj_ptr = ui_ptr[prow] + jmin - ui[prow]; /* points to the 2nd nzero entry in U(prow,k:am-1) */
3007:       PetscLLAddSorted(ncols,uj_ptr,am,nlnk,lnk,lnkbt);
3008:       nzk   += nlnk;

3010:       /* update il and jl for prow */
3011:       if (jmin < jmax) {
3012:         il[prow] = jmin;
3013:         j        = *uj_ptr; jl[prow] = jl[j]; jl[j] = prow;
3014:       }
3015:       prow = nextprow;
3016:     }

3018:     /* if free space is not available, make more free space */
3019:     if (current_space->local_remaining<nzk) {
3020:       i    = am - k + 1; /* num of unfactored rows */
3021:       i    = PetscMin(i*nzk, i*(i-1)); /* i*nzk, i*(i-1): estimated and max additional space needed */
3022:       PetscFreeSpaceGet(i,&current_space);
3023:       reallocs++;
3024:     }

3026:     /* copy data into free space, then initialize lnk */
3027:     PetscLLClean(am,am,nzk,lnk,current_space->array,lnkbt);

3029:     /* add the k-th row into il and jl */
3030:     if (nzk-1 > 0) {
3031:       i     = current_space->array[1]; /* col value of the first nonzero element in U(k, k+1:am-1) */
3032:       jl[k] = jl[i]; jl[i] = k;
3033:       il[k] = ui[k] + 1;
3034:     }
3035:     ui_ptr[k] = current_space->array;

3037:     current_space->array           += nzk;
3038:     current_space->local_used      += nzk;
3039:     current_space->local_remaining -= nzk;

3041:     ui[k+1] = ui[k] + nzk;
3042:   }

3044: #if defined(PETSC_USE_INFO)
3045:   if (ai[am] != 0) {
3046:     PetscReal af = (PetscReal)(ui[am])/((PetscReal)ai[am]);
3047:     PetscInfo3(A,"Reallocs %D Fill ratio:given %g needed %g\n",reallocs,(double)fill,(double)af);
3048:     PetscInfo1(A,"Run with -pc_factor_fill %g or use \n",(double)af);
3049:     PetscInfo1(A,"PCFactorSetFill(pc,%g) for best performance.\n",(double)af);
3050:   } else {
3051:     PetscInfo(A,"Empty matrix.\n");
3052:   }
3053: #endif

3055:   ISRestoreIndices(perm,&rip);
3056:   ISRestoreIndices(iperm,&riip);
3057:   PetscFree4(ui_ptr,jl,il,cols);

3059:   /* destroy list of free space and other temporary array(s) */
3060:   PetscMalloc1(ui[am]+1,&uj);
3061:   PetscFreeSpaceContiguous(&free_space,uj);
3062:   PetscLLDestroy(lnk,lnkbt);

3064:   /* put together the new matrix in MATSEQSBAIJ format */

3066:   b               = (Mat_SeqSBAIJ*)fact->data;
3067:   b->singlemalloc = PETSC_FALSE;
3068:   b->free_a       = PETSC_TRUE;
3069:   b->free_ij      = PETSC_TRUE;

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

3073:   b->j    = uj;
3074:   b->i    = ui;
3075:   b->diag = 0;
3076:   b->ilen = 0;
3077:   b->imax = 0;
3078:   b->row  = perm;
3079:   b->col  = perm;

3081:   PetscObjectReference((PetscObject)perm);
3082:   PetscObjectReference((PetscObject)perm);

3084:   b->icol          = iperm;
3085:   b->pivotinblocks = PETSC_FALSE; /* need to get from MatFactorInfo */

3087:   PetscMalloc1(am+1,&b->solve_work);
3088:   PetscLogObjectMemory((PetscObject)fact,(ui[am]-am)*(sizeof(PetscInt)+sizeof(MatScalar)));
3089:   b->maxnz = b->nz = ui[am];

3091:   fact->info.factor_mallocs   = reallocs;
3092:   fact->info.fill_ratio_given = fill;
3093:   if (ai[am] != 0) {
3094:     fact->info.fill_ratio_needed = ((PetscReal)ui[am])/((PetscReal)ai[am]);
3095:   } else {
3096:     fact->info.fill_ratio_needed = 0.0;
3097:   }
3098:   fact->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqAIJ_inplace;
3099:   return(0);
3100: }

3102: PetscErrorCode MatSolve_SeqAIJ_NaturalOrdering(Mat A,Vec bb,Vec xx)
3103: {
3104:   Mat_SeqAIJ        *a = (Mat_SeqAIJ*)A->data;
3105:   PetscErrorCode    ierr;
3106:   PetscInt          n   = A->rmap->n;
3107:   const PetscInt    *ai = a->i,*aj = a->j,*adiag = a->diag,*vi;
3108:   PetscScalar       *x,sum;
3109:   const PetscScalar *b;
3110:   const MatScalar   *aa = a->a,*v;
3111:   PetscInt          i,nz;

3114:   if (!n) return(0);

3116:   VecGetArrayRead(bb,&b);
3117:   VecGetArray(xx,&x);

3119:   /* forward solve the lower triangular */
3120:   x[0] = b[0];
3121:   v    = aa;
3122:   vi   = aj;
3123:   for (i=1; i<n; i++) {
3124:     nz  = ai[i+1] - ai[i];
3125:     sum = b[i];
3126:     PetscSparseDenseMinusDot(sum,x,v,vi,nz);
3127:     v   += nz;
3128:     vi  += nz;
3129:     x[i] = sum;
3130:   }

3132:   /* backward solve the upper triangular */
3133:   for (i=n-1; i>=0; i--) {
3134:     v   = aa + adiag[i+1] + 1;
3135:     vi  = aj + adiag[i+1] + 1;
3136:     nz  = adiag[i] - adiag[i+1]-1;
3137:     sum = x[i];
3138:     PetscSparseDenseMinusDot(sum,x,v,vi,nz);
3139:     x[i] = sum*v[nz]; /* x[i]=aa[adiag[i]]*sum; v++; */
3140:   }

3142:   PetscLogFlops(2.0*a->nz - A->cmap->n);
3143:   VecRestoreArrayRead(bb,&b);
3144:   VecRestoreArray(xx,&x);
3145:   return(0);
3146: }

3148: PetscErrorCode MatSolve_SeqAIJ(Mat A,Vec bb,Vec xx)
3149: {
3150:   Mat_SeqAIJ        *a    = (Mat_SeqAIJ*)A->data;
3151:   IS                iscol = a->col,isrow = a->row;
3152:   PetscErrorCode    ierr;
3153:   PetscInt          i,n=A->rmap->n,*vi,*ai=a->i,*aj=a->j,*adiag = a->diag,nz;
3154:   const PetscInt    *rout,*cout,*r,*c;
3155:   PetscScalar       *x,*tmp,sum;
3156:   const PetscScalar *b;
3157:   const MatScalar   *aa = a->a,*v;

3160:   if (!n) return(0);

3162:   VecGetArrayRead(bb,&b);
3163:   VecGetArray(xx,&x);
3164:   tmp  = a->solve_work;

3166:   ISGetIndices(isrow,&rout); r = rout;
3167:   ISGetIndices(iscol,&cout); c = cout;

3169:   /* forward solve the lower triangular */
3170:   tmp[0] = b[r[0]];
3171:   v      = aa;
3172:   vi     = aj;
3173:   for (i=1; i<n; i++) {
3174:     nz  = ai[i+1] - ai[i];
3175:     sum = b[r[i]];
3176:     PetscSparseDenseMinusDot(sum,tmp,v,vi,nz);
3177:     tmp[i] = sum;
3178:     v     += nz; vi += nz;
3179:   }

3181:   /* backward solve the upper triangular */
3182:   for (i=n-1; i>=0; i--) {
3183:     v   = aa + adiag[i+1]+1;
3184:     vi  = aj + adiag[i+1]+1;
3185:     nz  = adiag[i]-adiag[i+1]-1;
3186:     sum = tmp[i];
3187:     PetscSparseDenseMinusDot(sum,tmp,v,vi,nz);
3188:     x[c[i]] = tmp[i] = sum*v[nz]; /* v[nz] = aa[adiag[i]] */
3189:   }

3191:   ISRestoreIndices(isrow,&rout);
3192:   ISRestoreIndices(iscol,&cout);
3193:   VecRestoreArrayRead(bb,&b);
3194:   VecRestoreArray(xx,&x);
3195:   PetscLogFlops(2*a->nz - A->cmap->n);
3196:   return(0);
3197: }

3199: /*
3200:     This will get a new name and become a varient of MatILUFactor_SeqAIJ() there is no longer separate functions in the matrix function table for dt factors
3201: */
3202: PetscErrorCode MatILUDTFactor_SeqAIJ(Mat A,IS isrow,IS iscol,const MatFactorInfo *info,Mat *fact)
3203: {
3204:   Mat            B = *fact;
3205:   Mat_SeqAIJ     *a=(Mat_SeqAIJ*)A->data,*b;
3206:   IS             isicol;
3208:   const PetscInt *r,*ic;
3209:   PetscInt       i,n=A->rmap->n,*ai=a->i,*aj=a->j,*ajtmp,*adiag;
3210:   PetscInt       *bi,*bj,*bdiag,*bdiag_rev;
3211:   PetscInt       row,nzi,nzi_bl,nzi_bu,*im,nzi_al,nzi_au;
3212:   PetscInt       nlnk,*lnk;
3213:   PetscBT        lnkbt;
3214:   PetscBool      row_identity,icol_identity;
3215:   MatScalar      *aatmp,*pv,*batmp,*ba,*rtmp,*pc,multiplier,*vtmp,diag_tmp;
3216:   const PetscInt *ics;
3217:   PetscInt       j,nz,*pj,*bjtmp,k,ncut,*jtmp;
3218:   PetscReal      dt     =info->dt,shift=info->shiftamount;
3219:   PetscInt       dtcount=(PetscInt)info->dtcount,nnz_max;
3220:   PetscBool      missing;

3223:   if (dt      == PETSC_DEFAULT) dt = 0.005;
3224:   if (dtcount == PETSC_DEFAULT) dtcount = (PetscInt)(1.5*a->rmax);

3226:   /* ------- symbolic factorization, can be reused ---------*/
3227:   MatMissingDiagonal(A,&missing,&i);
3228:   if (missing) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Matrix is missing diagonal entry %D",i);
3229:   adiag=a->diag;

3231:   ISInvertPermutation(iscol,PETSC_DECIDE,&isicol);

3233:   /* bdiag is location of diagonal in factor */
3234:   PetscMalloc1(n+1,&bdiag);     /* becomes b->diag */
3235:   PetscMalloc1(n+1,&bdiag_rev); /* temporary */

3237:   /* allocate row pointers bi */
3238:   PetscMalloc1(2*n+2,&bi);

3240:   /* allocate bj and ba; max num of nonzero entries is (ai[n]+2*n*dtcount+2) */
3241:   if (dtcount > n-1) dtcount = n-1; /* diagonal is excluded */
3242:   nnz_max = ai[n]+2*n*dtcount+2;

3244:   PetscMalloc1(nnz_max+1,&bj);
3245:   PetscMalloc1(nnz_max+1,&ba);

3247:   /* put together the new matrix */
3248:   MatSeqAIJSetPreallocation_SeqAIJ(B,MAT_SKIP_ALLOCATION,NULL);
3249:   PetscLogObjectParent((PetscObject)B,(PetscObject)isicol);
3250:   b    = (Mat_SeqAIJ*)B->data;

3252:   b->free_a       = PETSC_TRUE;
3253:   b->free_ij      = PETSC_TRUE;
3254:   b->singlemalloc = PETSC_FALSE;

3256:   b->a    = ba;
3257:   b->j    = bj;
3258:   b->i    = bi;
3259:   b->diag = bdiag;
3260:   b->ilen = 0;
3261:   b->imax = 0;
3262:   b->row  = isrow;
3263:   b->col  = iscol;
3264:   PetscObjectReference((PetscObject)isrow);
3265:   PetscObjectReference((PetscObject)iscol);
3266:   b->icol = isicol;

3268:   PetscMalloc1(n+1,&b->solve_work);
3269:   PetscLogObjectMemory((PetscObject)B,nnz_max*(sizeof(PetscInt)+sizeof(MatScalar)));
3270:   b->maxnz = nnz_max;

3272:   B->factortype            = MAT_FACTOR_ILUDT;
3273:   B->info.factor_mallocs   = 0;
3274:   B->info.fill_ratio_given = ((PetscReal)nnz_max)/((PetscReal)ai[n]);
3275:   /* ------- end of symbolic factorization ---------*/

3277:   ISGetIndices(isrow,&r);
3278:   ISGetIndices(isicol,&ic);
3279:   ics  = ic;

3281:   /* linked list for storing column indices of the active row */
3282:   nlnk = n + 1;
3283:   PetscLLCreate(n,n,nlnk,lnk,lnkbt);

3285:   /* im: used by PetscLLAddSortedLU(); jtmp: working array for column indices of active row */
3286:   PetscMalloc2(n,&im,n,&jtmp);
3287:   /* rtmp, vtmp: working arrays for sparse and contiguous row entries of active row */
3288:   PetscMalloc2(n,&rtmp,n,&vtmp);
3289:   PetscMemzero(rtmp,n*sizeof(MatScalar));

3291:   bi[0]        = 0;
3292:   bdiag[0]     = nnz_max-1; /* location of diag[0] in factor B */
3293:   bdiag_rev[n] = bdiag[0];
3294:   bi[2*n+1]    = bdiag[0]+1; /* endof bj and ba array */
3295:   for (i=0; i<n; i++) {
3296:     /* copy initial fill into linked list */
3297:     nzi = ai[r[i]+1] - ai[r[i]];
3298:     if (!nzi) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_MAT_LU_ZRPVT,"Empty row in matrix: row in original ordering %D in permuted ordering %D",r[i],i);
3299:     nzi_al = adiag[r[i]] - ai[r[i]];
3300:     nzi_au = ai[r[i]+1] - adiag[r[i]] -1;
3301:     ajtmp  = aj + ai[r[i]];
3302:     PetscLLAddPerm(nzi,ajtmp,ic,n,nlnk,lnk,lnkbt);

3304:     /* load in initial (unfactored row) */
3305:     aatmp = a->a + ai[r[i]];
3306:     for (j=0; j<nzi; j++) {
3307:       rtmp[ics[*ajtmp++]] = *aatmp++;
3308:     }

3310:     /* add pivot rows into linked list */
3311:     row = lnk[n];
3312:     while (row < i) {
3313:       nzi_bl = bi[row+1] - bi[row] + 1;
3314:       bjtmp  = bj + bdiag[row+1]+1; /* points to 1st column next to the diagonal in U */
3315:       PetscLLAddSortedLU(bjtmp,row,nlnk,lnk,lnkbt,i,nzi_bl,im);
3316:       nzi   += nlnk;
3317:       row    = lnk[row];
3318:     }

3320:     /* copy data from lnk into jtmp, then initialize lnk */
3321:     PetscLLClean(n,n,nzi,lnk,jtmp,lnkbt);

3323:     /* numerical factorization */
3324:     bjtmp = jtmp;
3325:     row   = *bjtmp++; /* 1st pivot row */
3326:     while (row < i) {
3327:       pc         = rtmp + row;
3328:       pv         = ba + bdiag[row]; /* 1./(diag of the pivot row) */
3329:       multiplier = (*pc) * (*pv);
3330:       *pc        = multiplier;
3331:       if (PetscAbsScalar(*pc) > dt) { /* apply tolerance dropping rule */
3332:         pj = bj + bdiag[row+1] + 1;         /* point to 1st entry of U(row,:) */
3333:         pv = ba + bdiag[row+1] + 1;
3334:         /* if (multiplier < -1.0 or multiplier >1.0) printf("row/prow %d, %d, multiplier %g\n",i,row,multiplier); */
3335:         nz = bdiag[row] - bdiag[row+1] - 1;         /* num of entries in U(row,:), excluding diagonal */
3336:         for (j=0; j<nz; j++) rtmp[*pj++] -= multiplier * (*pv++);
3337:         PetscLogFlops(1+2*nz);
3338:       }
3339:       row = *bjtmp++;
3340:     }

3342:     /* copy sparse rtmp into contiguous vtmp; separate L and U part */
3343:     diag_tmp = rtmp[i];  /* save diagonal value - may not needed?? */
3344:     nzi_bl   = 0; j = 0;
3345:     while (jtmp[j] < i) { /* Note: jtmp is sorted */
3346:       vtmp[j] = rtmp[jtmp[j]]; rtmp[jtmp[j]]=0.0;
3347:       nzi_bl++; j++;
3348:     }
3349:     nzi_bu = nzi - nzi_bl -1;
3350:     while (j < nzi) {
3351:       vtmp[j] = rtmp[jtmp[j]]; rtmp[jtmp[j]]=0.0;
3352:       j++;
3353:     }

3355:     bjtmp = bj + bi[i];
3356:     batmp = ba + bi[i];
3357:     /* apply level dropping rule to L part */
3358:     ncut = nzi_al + dtcount;
3359:     if (ncut < nzi_bl) {
3360:       PetscSortSplit(ncut,nzi_bl,vtmp,jtmp);
3361:       PetscSortIntWithScalarArray(ncut,jtmp,vtmp);
3362:     } else {
3363:       ncut = nzi_bl;
3364:     }
3365:     for (j=0; j<ncut; j++) {
3366:       bjtmp[j] = jtmp[j];
3367:       batmp[j] = vtmp[j];
3368:       /* printf(" (%d,%g),",bjtmp[j],batmp[j]); */
3369:     }
3370:     bi[i+1] = bi[i] + ncut;
3371:     nzi     = ncut + 1;

3373:     /* apply level dropping rule to U part */
3374:     ncut = nzi_au + dtcount;
3375:     if (ncut < nzi_bu) {
3376:       PetscSortSplit(ncut,nzi_bu,vtmp+nzi_bl+1,jtmp+nzi_bl+1);
3377:       PetscSortIntWithScalarArray(ncut,jtmp+nzi_bl+1,vtmp+nzi_bl+1);
3378:     } else {
3379:       ncut = nzi_bu;
3380:     }
3381:     nzi += ncut;

3383:     /* mark bdiagonal */
3384:     bdiag[i+1]       = bdiag[i] - (ncut + 1);
3385:     bdiag_rev[n-i-1] = bdiag[i+1];
3386:     bi[2*n - i]      = bi[2*n - i +1] - (ncut + 1);
3387:     bjtmp            = bj + bdiag[i];
3388:     batmp            = ba + bdiag[i];
3389:     *bjtmp           = i;
3390:     *batmp           = diag_tmp; /* rtmp[i]; */
3391:     if (*batmp == 0.0) {
3392:       *batmp = dt+shift;
3393:       /* printf(" row %d add shift %g\n",i,shift); */
3394:     }
3395:     *batmp = 1.0/(*batmp); /* invert diagonal entries for simplier triangular solves */
3396:     /* printf(" (%d,%g),",*bjtmp,*batmp); */

3398:     bjtmp = bj + bdiag[i+1]+1;
3399:     batmp = ba + bdiag[i+1]+1;
3400:     for (k=0; k<ncut; k++) {
3401:       bjtmp[k] = jtmp[nzi_bl+1+k];
3402:       batmp[k] = vtmp[nzi_bl+1+k];
3403:       /* printf(" (%d,%g),",bjtmp[k],batmp[k]); */
3404:     }
3405:     /* printf("\n"); */

3407:     im[i] = nzi;   /* used by PetscLLAddSortedLU() */
3408:     /*
3409:     printf("row %d: bi %d, bdiag %d\n",i,bi[i],bdiag[i]);
3410:     printf(" ----------------------------\n");
3411:     */
3412:   } /* for (i=0; i<n; i++) */
3413:     /* printf("end of L %d, beginning of U %d\n",bi[n],bdiag[n]); */
3414:   if (bi[n] >= bdiag[n]) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"end of L array %d cannot >= the beginning of U array %d",bi[n],bdiag[n]);

3416:   ISRestoreIndices(isrow,&r);
3417:   ISRestoreIndices(isicol,&ic);

3419:   PetscLLDestroy(lnk,lnkbt);
3420:   PetscFree2(im,jtmp);
3421:   PetscFree2(rtmp,vtmp);
3422:   PetscFree(bdiag_rev);

3424:   PetscLogFlops(B->cmap->n);
3425:   b->maxnz = b->nz = bi[n] + bdiag[0] - bdiag[n];

3427:   ISIdentity(isrow,&row_identity);
3428:   ISIdentity(isicol,&icol_identity);
3429:   if (row_identity && icol_identity) {
3430:     B->ops->solve = MatSolve_SeqAIJ_NaturalOrdering;
3431:   } else {
3432:     B->ops->solve = MatSolve_SeqAIJ;
3433:   }

3435:   B->ops->solveadd          = 0;
3436:   B->ops->solvetranspose    = 0;
3437:   B->ops->solvetransposeadd = 0;
3438:   B->ops->matsolve          = 0;
3439:   B->assembled              = PETSC_TRUE;
3440:   B->preallocated           = PETSC_TRUE;
3441:   return(0);
3442: }

3444: /* a wraper of MatILUDTFactor_SeqAIJ() */
3445: /*
3446:     This will get a new name and become a varient of MatILUFactor_SeqAIJ() there is no longer separate functions in the matrix function table for dt factors
3447: */

3449: PetscErrorCode  MatILUDTFactorSymbolic_SeqAIJ(Mat fact,Mat A,IS row,IS col,const MatFactorInfo *info)
3450: {

3454:   MatILUDTFactor_SeqAIJ(A,row,col,info,&fact);
3455:   return(0);
3456: }

3458: /*
3459:    same as MatLUFactorNumeric_SeqAIJ(), except using contiguous array matrix factors
3460:    - intend to replace existing MatLUFactorNumeric_SeqAIJ()
3461: */
3462: /*
3463:     This will get a new name and become a varient of MatILUFactor_SeqAIJ() there is no longer separate functions in the matrix function table for dt factors
3464: */

3466: PetscErrorCode  MatILUDTFactorNumeric_SeqAIJ(Mat fact,Mat A,const MatFactorInfo *info)
3467: {
3468:   Mat            C     =fact;
3469:   Mat_SeqAIJ     *a    =(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)C->data;
3470:   IS             isrow = b->row,isicol = b->icol;
3472:   const PetscInt *r,*ic,*ics;
3473:   PetscInt       i,j,k,n=A->rmap->n,*ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j;
3474:   PetscInt       *ajtmp,*bjtmp,nz,nzl,nzu,row,*bdiag = b->diag,*pj;
3475:   MatScalar      *rtmp,*pc,multiplier,*v,*pv,*aa=a->a;
3476:   PetscReal      dt=info->dt,shift=info->shiftamount;
3477:   PetscBool      row_identity, col_identity;

3480:   ISGetIndices(isrow,&r);
3481:   ISGetIndices(isicol,&ic);
3482:   PetscMalloc1(n+1,&rtmp);
3483:   ics  = ic;

3485:   for (i=0; i<n; i++) {
3486:     /* initialize rtmp array */
3487:     nzl   = bi[i+1] - bi[i];       /* num of nozeros in L(i,:) */
3488:     bjtmp = bj + bi[i];
3489:     for  (j=0; j<nzl; j++) rtmp[*bjtmp++] = 0.0;
3490:     rtmp[i] = 0.0;
3491:     nzu     = bdiag[i] - bdiag[i+1]; /* num of nozeros in U(i,:) */
3492:     bjtmp   = bj + bdiag[i+1] + 1;
3493:     for  (j=0; j<nzu; j++) rtmp[*bjtmp++] = 0.0;

3495:     /* load in initial unfactored row of A */
3496:     /* printf("row %d\n",i); */
3497:     nz    = ai[r[i]+1] - ai[r[i]];
3498:     ajtmp = aj + ai[r[i]];
3499:     v     = aa + ai[r[i]];
3500:     for (j=0; j<nz; j++) {
3501:       rtmp[ics[*ajtmp++]] = v[j];
3502:       /* printf(" (%d,%g),",ics[ajtmp[j]],rtmp[ics[ajtmp[j]]]); */
3503:     }
3504:     /* printf("\n"); */

3506:     /* numerical factorization */
3507:     bjtmp = bj + bi[i]; /* point to 1st entry of L(i,:) */
3508:     nzl   = bi[i+1] - bi[i]; /* num of entries in L(i,:) */
3509:     k     = 0;
3510:     while (k < nzl) {
3511:       row = *bjtmp++;
3512:       /* printf("  prow %d\n",row); */
3513:       pc         = rtmp + row;
3514:       pv         = b->a + bdiag[row]; /* 1./(diag of the pivot row) */
3515:       multiplier = (*pc) * (*pv);
3516:       *pc        = multiplier;
3517:       if (PetscAbsScalar(multiplier) > dt) {
3518:         pj = bj + bdiag[row+1] + 1;         /* point to 1st entry of U(row,:) */
3519:         pv = b->a + bdiag[row+1] + 1;
3520:         nz = bdiag[row] - bdiag[row+1] - 1;         /* num of entries in U(row,:), excluding diagonal */
3521:         for (j=0; j<nz; j++) rtmp[*pj++] -= multiplier * (*pv++);
3522:         PetscLogFlops(1+2*nz);
3523:       }
3524:       k++;
3525:     }

3527:     /* finished row so stick it into b->a */
3528:     /* L-part */
3529:     pv  = b->a + bi[i];
3530:     pj  = bj + bi[i];
3531:     nzl = bi[i+1] - bi[i];
3532:     for (j=0; j<nzl; j++) {
3533:       pv[j] = rtmp[pj[j]];
3534:       /* printf(" (%d,%g),",pj[j],pv[j]); */
3535:     }

3537:     /* diagonal: invert diagonal entries for simplier triangular solves */
3538:     if (rtmp[i] == 0.0) rtmp[i] = dt+shift;
3539:     b->a[bdiag[i]] = 1.0/rtmp[i];
3540:     /* printf(" (%d,%g),",i,b->a[bdiag[i]]); */

3542:     /* U-part */
3543:     pv  = b->a + bdiag[i+1] + 1;
3544:     pj  = bj + bdiag[i+1] + 1;
3545:     nzu = bdiag[i] - bdiag[i+1] - 1;
3546:     for (j=0; j<nzu; j++) {
3547:       pv[j] = rtmp[pj[j]];
3548:       /* printf(" (%d,%g),",pj[j],pv[j]); */
3549:     }
3550:     /* printf("\n"); */
3551:   }

3553:   PetscFree(rtmp);
3554:   ISRestoreIndices(isicol,&ic);
3555:   ISRestoreIndices(isrow,&r);

3557:   ISIdentity(isrow,&row_identity);
3558:   ISIdentity(isicol,&col_identity);
3559:   if (row_identity && col_identity) {
3560:     C->ops->solve = MatSolve_SeqAIJ_NaturalOrdering;
3561:   } else {
3562:     C->ops->solve = MatSolve_SeqAIJ;
3563:   }
3564:   C->ops->solveadd          = 0;
3565:   C->ops->solvetranspose    = 0;
3566:   C->ops->solvetransposeadd = 0;
3567:   C->ops->matsolve          = 0;
3568:   C->assembled              = PETSC_TRUE;
3569:   C->preallocated           = PETSC_TRUE;

3571:   PetscLogFlops(C->cmap->n);
3572:   return(0);
3573: }