Actual source code: matimpl.h


 5:  #include petscmat.h

  7: /*
  8:   This file defines the parts of the matrix data structure that are 
  9:   shared by all matrix types.
 10: */

 12: /*
 13:     If you add entries here also add them to the MATOP enum
 14:     in include/petscmat.h and include/finclude/petscmat.h
 15: */
 16: typedef struct _MatOps *MatOps;
 17: struct _MatOps {
 18:   /* 0*/
 19:   PetscErrorCode (*setvalues)(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const MatScalar[],InsertMode);
 20:   PetscErrorCode (*getrow)(Mat,PetscInt,PetscInt *,PetscInt*[],PetscScalar*[]);
 21:   PetscErrorCode (*restorerow)(Mat,PetscInt,PetscInt *,PetscInt *[],PetscScalar *[]);
 22:   PetscErrorCode (*mult)(Mat,Vec,Vec);
 23:   PetscErrorCode (*multadd)(Mat,Vec,Vec,Vec);
 24:   /* 5*/
 25:   PetscErrorCode (*multtranspose)(Mat,Vec,Vec);
 26:   PetscErrorCode (*multtransposeadd)(Mat,Vec,Vec,Vec);
 27:   PetscErrorCode (*solve)(Mat,Vec,Vec);
 28:   PetscErrorCode (*solveadd)(Mat,Vec,Vec,Vec);
 29:   PetscErrorCode (*solvetranspose)(Mat,Vec,Vec);
 30:   /*10*/
 31:   PetscErrorCode (*solvetransposeadd)(Mat,Vec,Vec,Vec);
 32:   PetscErrorCode (*lufactor)(Mat,IS,IS,MatFactorInfo*);
 33:   PetscErrorCode (*choleskyfactor)(Mat,IS,MatFactorInfo*);
 34:   PetscErrorCode (*relax)(Mat,Vec,PetscReal,MatSORType,PetscReal,PetscInt,PetscInt,Vec);
 35:   PetscErrorCode (*transpose)(Mat,Mat *);
 36:   /*15*/
 37:   PetscErrorCode (*getinfo)(Mat,MatInfoType,MatInfo*);
 38:   PetscErrorCode (*equal)(Mat,Mat,PetscTruth *);
 39:   PetscErrorCode (*getdiagonal)(Mat,Vec);
 40:   PetscErrorCode (*diagonalscale)(Mat,Vec,Vec);
 41:   PetscErrorCode (*norm)(Mat,NormType,PetscReal*);
 42:   /*20*/
 43:   PetscErrorCode (*assemblybegin)(Mat,MatAssemblyType);
 44:   PetscErrorCode (*assemblyend)(Mat,MatAssemblyType);
 45:   PetscErrorCode (*compress)(Mat);
 46:   PetscErrorCode (*setoption)(Mat,MatOption);
 47:   PetscErrorCode (*zeroentries)(Mat);
 48:   /*25*/
 49:   PetscErrorCode (*zerorows)(Mat,PetscInt,const PetscInt[],PetscScalar);
 50:   PetscErrorCode (*lufactorsymbolic)(Mat,IS,IS,MatFactorInfo*,Mat*);
 51:   PetscErrorCode (*lufactornumeric)(Mat,MatFactorInfo*,Mat*);
 52:   PetscErrorCode (*choleskyfactorsymbolic)(Mat,IS,MatFactorInfo*,Mat*);
 53:   PetscErrorCode (*choleskyfactornumeric)(Mat,MatFactorInfo*,Mat*);
 54:   /*30*/
 55:   PetscErrorCode (*setuppreallocation)(Mat);
 56:   PetscErrorCode (*ilufactorsymbolic)(Mat,IS,IS,MatFactorInfo*,Mat*);
 57:   PetscErrorCode (*iccfactorsymbolic)(Mat,IS,MatFactorInfo*,Mat*);
 58:   PetscErrorCode (*getarray)(Mat,PetscScalar**);
 59:   PetscErrorCode (*restorearray)(Mat,PetscScalar**);
 60:   /*35*/
 61:   PetscErrorCode (*duplicate)(Mat,MatDuplicateOption,Mat*);
 62:   PetscErrorCode (*forwardsolve)(Mat,Vec,Vec);
 63:   PetscErrorCode (*backwardsolve)(Mat,Vec,Vec);
 64:   PetscErrorCode (*ilufactor)(Mat,IS,IS,MatFactorInfo*);
 65:   PetscErrorCode (*iccfactor)(Mat,IS,MatFactorInfo*);
 66:   /*40*/
 67:   PetscErrorCode (*axpy)(Mat,PetscScalar,Mat,MatStructure);
 68:   PetscErrorCode (*getsubmatrices)(Mat,PetscInt,const IS[],const IS[],MatReuse,Mat *[]);
 69:   PetscErrorCode (*increaseoverlap)(Mat,PetscInt,IS[],PetscInt);
 70:   PetscErrorCode (*getvalues)(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],PetscScalar []);
 71:   PetscErrorCode (*copy)(Mat,Mat,MatStructure);
 72:   /*45*/
 73:   PetscErrorCode (*dummy)(void);
 74:   PetscErrorCode (*scale)(Mat,PetscScalar);
 75:   PetscErrorCode (*shift)(Mat,PetscScalar);
 76:   PetscErrorCode (*diagonalset)(Mat,Vec,InsertMode);
 77:   PetscErrorCode (*iludtfactor)(Mat,IS,IS,MatFactorInfo*,Mat *);
 78:   /*50*/
 79:   PetscErrorCode (*setblocksize)(Mat,PetscInt);
 80:   PetscErrorCode (*getrowij)(Mat,PetscInt,PetscTruth,PetscInt*,PetscInt *[],PetscInt *[],PetscTruth *);
 81:   PetscErrorCode (*restorerowij)(Mat,PetscInt,PetscTruth,PetscInt *,PetscInt *[],PetscInt *[],PetscTruth *);
 82:   PetscErrorCode (*getcolumnij)(Mat,PetscInt,PetscTruth,PetscInt*,PetscInt *[],PetscInt *[],PetscTruth *);
 83:   PetscErrorCode (*restorecolumnij)(Mat,PetscInt,PetscTruth,PetscInt*,PetscInt *[],PetscInt *[],PetscTruth *);
 84:   /*55*/
 85:   PetscErrorCode (*fdcoloringcreate)(Mat,ISColoring,MatFDColoring);
 86:   PetscErrorCode (*coloringpatch)(Mat,PetscInt,PetscInt,ISColoringValue[],ISColoring*);
 87:   PetscErrorCode (*setunfactored)(Mat);
 88:   PetscErrorCode (*permute)(Mat,IS,IS,Mat*);
 89:   PetscErrorCode (*setvaluesblocked)(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
 90:   /*60*/
 91:   PetscErrorCode (*getsubmatrix)(Mat,IS,IS,PetscInt,MatReuse,Mat*);
 92:   PetscErrorCode (*destroy)(Mat);
 93:   PetscErrorCode (*view)(Mat,PetscViewer);
 94:   PetscErrorCode (*convertfrom)(Mat, MatType,MatReuse,Mat*);
 95:   PetscErrorCode (*usescaledform)(Mat,PetscTruth);
 96:   /*65*/
 97:   PetscErrorCode (*scalesystem)(Mat,Vec,Vec);
 98:   PetscErrorCode (*unscalesystem)(Mat,Vec,Vec);
 99:   PetscErrorCode (*setlocaltoglobalmapping)(Mat,ISLocalToGlobalMapping);
100:   PetscErrorCode (*setvalueslocal)(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
101:   PetscErrorCode (*zerorowslocal)(Mat,PetscInt,const PetscInt[],PetscScalar);
102:   /*70*/
103:   PetscErrorCode (*getrowmax)(Mat,Vec);
104:   PetscErrorCode (*convert)(Mat, MatType,MatReuse,Mat*);
105:   PetscErrorCode (*setcoloring)(Mat,ISColoring);
106:   PetscErrorCode (*setvaluesadic)(Mat,void*);
107:   PetscErrorCode (*setvaluesadifor)(Mat,PetscInt,void*);
108:   /*75*/
109:   PetscErrorCode (*fdcoloringapply)(Mat,MatFDColoring,Vec,MatStructure*,void*);
110:   PetscErrorCode (*setfromoptions)(Mat);
111:   PetscErrorCode (*multconstrained)(Mat,Vec,Vec);
112:   PetscErrorCode (*multtransposeconstrained)(Mat,Vec,Vec);
113:   PetscErrorCode (*ilufactorsymbolicconstrained)(Mat,IS,IS,double,PetscInt,PetscInt,Mat *);
114:   /*80*/
115:   PetscErrorCode (*permutesparsify)(Mat, PetscInt, double, double, IS, IS, Mat *);
116:   PetscErrorCode (*mults)(Mat, Vecs, Vecs);
117:   PetscErrorCode (*solves)(Mat, Vecs, Vecs);
118:   PetscErrorCode (*getinertia)(Mat,PetscInt*,PetscInt*,PetscInt*);
119:   PetscErrorCode (*load)(PetscViewer, MatType,Mat*);
120:   /*85*/
121:   PetscErrorCode (*issymmetric)(Mat,PetscReal,PetscTruth*);
122:   PetscErrorCode (*ishermitian)(Mat,PetscTruth*);
123:   PetscErrorCode (*isstructurallysymmetric)(Mat,PetscTruth*);
124:   PetscErrorCode (*pbrelax)(Mat,Vec,PetscReal,MatSORType,PetscReal,PetscInt,PetscInt,Vec);
125:   PetscErrorCode (*getvecs)(Mat,Vec*,Vec*);
126:   /*90*/
127:   PetscErrorCode (*matmult)(Mat,Mat,MatReuse,PetscReal,Mat*);
128:   PetscErrorCode (*matmultsymbolic)(Mat,Mat,PetscReal,Mat*);
129:   PetscErrorCode (*matmultnumeric)(Mat,Mat,Mat);
130:   PetscErrorCode (*ptap)(Mat,Mat,MatReuse,PetscReal,Mat*);
131:   PetscErrorCode (*ptapsymbolic)(Mat,Mat,PetscReal,Mat*); /* double dispatch wrapper routine */
132:   /*95*/
133:   PetscErrorCode (*ptapnumeric)(Mat,Mat,Mat);             /* double dispatch wrapper routine */
134:   PetscErrorCode (*matmulttranspose)(Mat,Mat,MatReuse,PetscReal,Mat*);
135:   PetscErrorCode (*matmulttransposesymbolic)(Mat,Mat,PetscReal,Mat*);
136:   PetscErrorCode (*matmulttransposenumeric)(Mat,Mat,Mat);
137:   PetscErrorCode (*ptapsymbolic_seqaij)(Mat,Mat,PetscReal,Mat*); /* actual implememtation, A=seqaij */
138:   /*100*/
139:   PetscErrorCode (*ptapnumeric_seqaij)(Mat,Mat,Mat);             /* actual implememtation, A=seqaij */
140:   PetscErrorCode (*ptapsymbolic_mpiaij)(Mat,Mat,PetscReal,Mat*); /* actual implememtation, A=mpiaij */
141:   PetscErrorCode (*ptapnumeric_mpiaij)(Mat,Mat,Mat);             /* actual implememtation, A=mpiaij */
142:   PetscErrorCode (*conjugate)(Mat);                              /* complex conjugate */
143:   PetscErrorCode (*setsizes)(Mat,PetscInt,PetscInt,PetscInt,PetscInt);
144:   /*105*/
145:   PetscErrorCode (*setvaluesrow)(Mat,PetscInt,const MatScalar[]);
146:   PetscErrorCode (*realpart)(Mat);
147:   PetscErrorCode (*imaginarypart)(Mat);
148:   PetscErrorCode (*getrowuppertriangular)(Mat);
149:   PetscErrorCode (*restorerowuppertriangular)(Mat);
150:   /*110*/
151:   PetscErrorCode (*matsolve)(Mat,Mat,Mat);
152: };
153: /*
154:     If you add MatOps entries above also add them to the MATOP enum
155:     in include/petscmat.h and include/finclude/petscmat.h
156: */

158: /*
159:    Utility private matrix routines
160: */
161: EXTERN PetscErrorCode MatConvert_Basic(Mat, MatType,MatReuse,Mat*);
162: EXTERN PetscErrorCode MatCopy_Basic(Mat,Mat,MatStructure);
163: EXTERN PetscErrorCode MatView_Private(Mat);

165: EXTERN PetscErrorCode MatHeaderCopy(Mat,Mat);
166: EXTERN PetscErrorCode MatHeaderReplace(Mat,Mat);
167: EXTERN PetscErrorCode MatAXPYGetxtoy_Private(PetscInt,PetscInt*,PetscInt*,PetscInt*, PetscInt*,PetscInt*,PetscInt*, PetscInt**);
168: EXTERN PetscErrorCode MatPtAP_Basic(Mat,Mat,MatReuse,PetscReal,Mat*);
169: EXTERN PetscErrorCode MatDiagonalSet_Default(Mat,Vec,InsertMode);

171: /* 
172:   The stash is used to temporarily store inserted matrix values that 
173:   belong to another processor. During the assembly phase the stashed 
174:   values are moved to the correct processor and 
175: */
176:  #include src/mat/utils/matstashspace.h
177: typedef struct {
178:   PetscInt      nmax;                   /* maximum stash size */
179:   PetscInt      umax;                   /* user specified max-size */
180:   PetscInt      oldnmax;                /* the nmax value used previously */
181:   PetscInt      n;                      /* stash size */
182:   PetscInt      bs;                     /* block size of the stash */
183:   PetscInt      reallocs;               /* preserve the no of mallocs invoked */
184:   PetscMatStashSpace space_head,space;  /* linked list to hold stashed global row/column numbers and matrix values */
185:   /* The following variables are used for communication */
186:   MPI_Comm      comm;
187:   PetscMPIInt   size,rank;
188:   PetscMPIInt   tag1,tag2;
189:   MPI_Request   *send_waits;            /* array of send requests */
190:   MPI_Request   *recv_waits;            /* array of receive requests */
191:   MPI_Status    *send_status;           /* array of send status */
192:   PetscInt      nsends,nrecvs;          /* numbers of sends and receives */
193:   MatScalar     *svalues;               /* sending data */
194:   MatScalar     **rvalues;              /* receiving data (values) */
195:   PetscInt      **rindices;             /* receiving data (indices) */
196:   PetscMPIInt   *nprocs;                /* tmp data used both during scatterbegin and end */
197:   PetscInt      nprocessed;             /* number of messages already processed */
198: } MatStash;

200: EXTERN PetscErrorCode MatStashCreate_Private(MPI_Comm,PetscInt,MatStash*);
201: EXTERN PetscErrorCode MatStashDestroy_Private(MatStash*);
202: EXTERN PetscErrorCode MatStashScatterEnd_Private(MatStash*);
203: EXTERN PetscErrorCode MatStashSetInitialSize_Private(MatStash*,PetscInt);
204: EXTERN PetscErrorCode MatStashGetInfo_Private(MatStash*,PetscInt*,PetscInt*);
205: EXTERN PetscErrorCode MatStashValuesRow_Private(MatStash*,PetscInt,PetscInt,const PetscInt[],const MatScalar[]);
206: EXTERN PetscErrorCode MatStashValuesCol_Private(MatStash*,PetscInt,PetscInt,const PetscInt[],const MatScalar[],PetscInt);
207: EXTERN PetscErrorCode MatStashValuesRowBlocked_Private(MatStash*,PetscInt,PetscInt,const PetscInt[],const MatScalar[],PetscInt,PetscInt,PetscInt);
208: EXTERN PetscErrorCode MatStashValuesColBlocked_Private(MatStash*,PetscInt,PetscInt,const PetscInt[],const MatScalar[],PetscInt,PetscInt,PetscInt);
209: EXTERN PetscErrorCode MatStashScatterBegin_Private(MatStash*,PetscInt*);
210: EXTERN PetscErrorCode MatStashScatterGetMesg_Private(MatStash*,PetscMPIInt*,PetscInt**,PetscInt**,MatScalar**,PetscInt*);

212: #define FACTOR_LU       1
213: #define FACTOR_CHOLESKY 2

215: typedef struct {
216:   PetscInt   dim;
217:   PetscInt   dims[4];
218:   PetscInt   starts[4];
219:   PetscTruth noc;        /* this is a single component problem, hence user will not set MatStencil.c */
220: } MatStencilInfo;

222: /* Info about using compressed row format */
223: typedef struct {
224:   PetscTruth use;
225:   PetscInt   nrows;                         /* number of non-zero rows */
226:   PetscInt   *i;                            /* compressed row pointer  */
227:   PetscInt   *rindex;                       /* compressed row index               */
228:   PetscTruth checked;                       /* if compressed row format have been checked for */
229: } Mat_CompressedRow;
230: EXTERN PetscErrorCode Mat_CheckCompressedRow(Mat,Mat_CompressedRow*,PetscInt*,PetscInt,PetscReal);

232: struct _p_Mat {
233:   PETSCHEADER(struct _MatOps);
234:   PetscMap               rmap,cmap;
235:   void                   *data;            /* implementation-specific data */
236:   PetscInt               factor;           /* 0, FACTOR_LU, or FACTOR_CHOLESKY */
237:   PetscTruth             assembled;        /* is the matrix assembled? */
238:   PetscTruth             was_assembled;    /* new values inserted into assembled mat */
239:   PetscInt               num_ass;          /* number of times matrix has been assembled */
240:   PetscTruth             same_nonzero;     /* matrix has same nonzero pattern as previous */
241:   MatInfo                info;             /* matrix information */
242:   ISLocalToGlobalMapping mapping;          /* mapping used in MatSetValuesLocal() */
243:   ISLocalToGlobalMapping bmapping;         /* mapping used in MatSetValuesBlockedLocal() */
244:   InsertMode             insertmode;       /* have values been inserted in matrix or added? */
245:   MatStash               stash,bstash;     /* used for assembling off-proc mat emements */
246:   MatNullSpace           nullsp;
247:   PetscTruth             preallocated;
248:   MatStencilInfo         stencil;          /* information for structured grid */
249:   PetscTruth             symmetric,hermitian,structurally_symmetric;
250:   PetscTruth             symmetric_set,hermitian_set,structurally_symmetric_set; /* if true, then corresponding flag is correct*/
251:   PetscTruth             symmetric_eternal;
252:   void                   *spptr;          /* pointer for special library like SuperLU */
253: };

255: #define MatPreallocated(A)  ((!(A)->preallocated) ? MatSetUpPreallocation(A) : 0)

258: /*
259:     Object for partitioning graphs
260: */

262: typedef struct _MatPartitioningOps *MatPartitioningOps;
263: struct _MatPartitioningOps {
264:   PetscErrorCode (*apply)(MatPartitioning,IS*);
265:   PetscErrorCode (*setfromoptions)(MatPartitioning);
266:   PetscErrorCode (*destroy)(MatPartitioning);
267:   PetscErrorCode (*view)(MatPartitioning,PetscViewer);
268: };

270: struct _p_MatPartitioning {
271:   PETSCHEADER(struct _MatPartitioningOps);
272:   Mat         adj;
273:   PetscInt    *vertex_weights;
274:   PetscReal   *part_weights;
275:   PetscInt    n;                                 /* number of partitions */
276:   void        *data;
277:   PetscInt    setupcalled;
278: };

280: /*
281:     MatFDColoring is used to compute Jacobian matrices efficiently
282:   via coloring. The data structure is explained below in an example.

284:    Color =   0    1     0    2   |   2      3       0 
285:    ---------------------------------------------------
286:             00   01              |          05
287:             10   11              |   14     15               Processor  0
288:                        22    23  |          25
289:                        32    33  | 
290:    ===================================================
291:                                  |   44     45     46
292:             50                   |          55               Processor 1
293:                                  |   64            66
294:    ---------------------------------------------------

296:     ncolors = 4;

298:     ncolumns      = {2,1,1,0}
299:     columns       = {{0,2},{1},{3},{}}
300:     nrows         = {4,2,3,3}
301:     rows          = {{0,1,2,3},{0,1},{1,2,3},{0,1,2}}
302:     columnsforrow = {{0,0,2,2},{1,1},{4,3,3},{5,5,5}}
303:     vscaleforrow  = {{,,,},{,},{,,},{,,}}
304:     vwscale       = {dx(0),dx(1),dx(2),dx(3)}               MPI Vec
305:     vscale        = {dx(0),dx(1),dx(2),dx(3),dx(4),dx(5)}   Seq Vec

307:     ncolumns      = {1,0,1,1}
308:     columns       = {{6},{},{4},{5}}
309:     nrows         = {3,0,2,2}
310:     rows          = {{0,1,2},{},{1,2},{1,2}}
311:     columnsforrow = {{6,0,6},{},{4,4},{5,5}}
312:     vscaleforrow =  {{,,},{},{,},{,}}
313:     vwscale       = {dx(4),dx(5),dx(6)}              MPI Vec
314:     vscale        = {dx(0),dx(4),dx(5),dx(6)}        Seq Vec

316:     See the routine MatFDColoringApply() for how this data is used
317:     to compute the Jacobian.

319: */

321: struct  _p_MatFDColoring{
322:   PETSCHEADER(int);
323:   PetscInt       M,N,m;            /* total rows, columns; local rows */
324:   PetscInt       rstart;           /* first row owned by local processor */
325:   PetscInt       ncolors;          /* number of colors */
326:   PetscInt       *ncolumns;        /* number of local columns for a color */
327:   PetscInt       **columns;        /* lists the local columns of each color (using global column numbering) */
328:   PetscInt       *nrows;           /* number of local rows for each color */
329:   PetscInt       **rows;           /* lists the local rows for each color (using the local row numbering) */
330:   PetscInt       **columnsforrow;  /* lists the corresponding columns for those rows (using the global column) */
331:   PetscReal      error_rel;        /* square root of relative error in computing function */
332:   PetscReal      umin;             /* minimum allowable u'dx value */
333:   PetscInt       freq;             /* frequency at which new Jacobian is computed */
334:   Vec            w1,w2,w3;         /* work vectors used in computing Jacobian */
335:   PetscErrorCode (*f)(void);       /* function that defines Jacobian */
336:   void           *fctx;            /* optional user-defined context for use by the function f */
337:   PetscInt       **vscaleforrow;   /* location in vscale for each columnsforrow[] entry */
338:   Vec            vscale;           /* holds FD scaling, i.e. 1/dx for each perturbed column */
339:   PetscTruth     usersetsrecompute;/* user determines when Jacobian is recomputed, via MatFDColoringSetRecompute() */
340:   PetscTruth     recompute;        /* used with usersetrecompute to determine if Jacobian should be recomputed */
341:   Vec            F;                /* current value of user provided function; can set with MatFDColoringSetF() */
342:   PetscInt       currentcolor;     /* color for which function evaluation is being done now */
343:   const char     *htype;            /* "wp" or "ds" */
344:   ISColoringType ctype;            /* IS_COLORING_GLOBAL or IS_COLORING_GHOSTED */
345: };

347: /*
348:    Null space context for preconditioner/operators
349: */
350: struct _p_MatNullSpace {
351:   PETSCHEADER(int);
352:   PetscTruth     has_cnst;
353:   PetscInt       n;
354:   Vec*           vecs;
355:   Vec            vec;                   /* for out of place removals */
356:   PetscErrorCode (*remove)(Vec,void*);  /* for user provided removal function */
357:   void*          rmctx;                 /* context for remove() function */
358: };

360: /* 
361:    Checking zero pivot for LU, ILU preconditioners.
362: */
363: typedef struct {
364:   PetscInt       nshift,nshift_max;
365:   PetscReal      shift_amount,shift_lo,shift_hi,shift_top;
366:   PetscTruth     lushift;
367:   PetscReal      rs;  /* active row sum of abs(offdiagonals) */
368:   PetscScalar    pv;  /* pivot of the active row */
369: } LUShift_Ctx;

371: EXTERN PetscErrorCode MatFactorDumpMatrix(Mat);

375: /*@C
376:    MatLUCheckShift_inline - shift the diagonals when zero pivot is detected on LU factor

378:    Collective on Mat

380:    Input Parameters:
381: +  info - information about the matrix factorization 
382: .  sctx - pointer to the struct LUShift_Ctx
383: -  row  - active row index

385:    Output  Parameter:
386: +  newshift - 0: shift is unchanged; 1: shft is updated; -1: zeropivot  

388:    Level: developer
389: @*/
390: #define MatLUCheckShift_inline(info,sctx,row,newshift) 0;\
391: {\
392:   PetscInt  _newshift;\
393:   PetscReal _rs   = sctx.rs;\
394:   PetscReal _zero = info->zeropivot*_rs;\
395:   if (info->shiftnz && PetscAbsScalar(sctx.pv) <= _zero){\
396:     /* force |diag| > zeropivot*rs */\
397:     if (!sctx.nshift){\
398:       sctx.shift_amount = info->shiftnz;\
399:     } else {\
400:       sctx.shift_amount *= 2.0;\
401:     }\
402:     sctx.lushift = PETSC_TRUE;\
403:     (sctx.nshift)++;\
404:     _newshift = 1;\
405:   } else if (info->shiftpd && PetscRealPart(sctx.pv) <= _zero){\
406:     /* force matfactor to be diagonally dominant */\
407:     if (sctx.nshift > sctx.nshift_max) {\
408:       MatFactorDumpMatrix(A);\
409:       SETERRQ1(PETSC_ERR_CONV_FAILED,"Unable to determine shift to enforce positive definite preconditioner after %d tries",sctx.nshift);\
410:     } else if (sctx.nshift == sctx.nshift_max) {\
411:       info->shift_fraction = sctx.shift_hi;\
412:       sctx.lushift        = PETSC_TRUE;\
413:     } else {\
414:       sctx.shift_lo = info->shift_fraction;\
415:       info->shift_fraction = (sctx.shift_hi+sctx.shift_lo)/2.;\
416:       sctx.lushift  = PETSC_TRUE;\
417:     }\
418:     sctx.shift_amount = info->shift_fraction * sctx.shift_top;\
419:     sctx.nshift++;\
420:     _newshift = 1;\
421:   } else if (PetscAbsScalar(sctx.pv) <= _zero){\
422:     MatFactorDumpMatrix(A);\
423:     SETERRQ4(PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot row %D value %G tolerance %G * rowsum %G",row,PetscAbsScalar(sctx.pv),_zero,_rs); \
424:   } else {\
425:     _newshift = 0;\
426:   }\
427:   newshift = _newshift;\
428: }

430: /* 
431:    Checking zero pivot for Cholesky, ICC preconditioners.
432: */
433: typedef struct {
434:   PetscInt       nshift;
435:   PetscReal      shift_amount;
436:   PetscTruth     chshift;
437:   PetscReal      rs;  /* active row sum of abs(offdiagonals) */
438:   PetscScalar    pv;  /* pivot of the active row */
439: } ChShift_Ctx;

443: /*@C
444:    MatCholeskyCheckShift_inline -  shift the diagonals when zero pivot is detected on Cholesky factor

446:    Collective on Mat

448:    Input Parameters:
449: +  info - information about the matrix factorization 
450: .  sctx - pointer to the struct CholeskyShift_Ctx
451: .  row  - pivot row
452: -  newshift - 0: shift is unchanged; 1: shft is updated; -1: zeropivot  

454:    Level: developer
455:    Note: Unlike in the ILU case there is no exit condition on nshift:
456:        we increase the shift until it converges. There is no guarantee that
457:        this algorithm converges faster or slower, or is better or worse
458:        than the ILU algorithm. 
459: @*/
460: #define MatCholeskyCheckShift_inline(info,sctx,row,newshift) 0;        \
461: {\
462:   PetscInt  _newshift;\
463:   PetscReal _rs   = sctx.rs;\
464:   PetscReal _zero = info->zeropivot*_rs;\
465:   if (info->shiftnz && PetscAbsScalar(sctx.pv) <= _zero){\
466:     /* force |diag| > zeropivot*sctx.rs */\
467:     if (!sctx.nshift){\
468:       sctx.shift_amount = info->shiftnz;\
469:     } else {\
470:       sctx.shift_amount *= 2.0;\
471:     }\
472:     sctx.chshift = PETSC_TRUE;\
473:     sctx.nshift++;\
474:     _newshift = 1;\
475:   } else if (info->shiftpd && PetscRealPart(sctx.pv) <= _zero){\
476:     /* calculate a shift that would make this row diagonally dominant */\
477:     sctx.shift_amount = PetscMax(_rs+PetscAbs(PetscRealPart(sctx.pv)),1.1*sctx.shift_amount);\
478:     sctx.chshift      = PETSC_TRUE;\
479:     sctx.nshift++;\
480:     _newshift = 1;\
481:   } else if (PetscAbsScalar(sctx.pv) <= _zero){\
482:     SETERRQ4(PETSC_ERR_MAT_CH_ZRPVT,"Zero pivot row %D value %G tolerance %G * rowsum %G",row,PetscAbsScalar(sctx.pv),_zero,_rs); \
483:   } else {\
484:     _newshift = 0; \
485:   }\
486:   newshift = _newshift;\
487: }

489: /* 
490:   Create and initialize a linked list 
491:   Input Parameters:
492:     idx_start - starting index of the list
493:     lnk_max   - max value of lnk indicating the end of the list
494:     nlnk      - max length of the list
495:   Output Parameters:
496:     lnk       - list initialized
497:     bt        - PetscBT (bitarray) with all bits set to false
498: */
499: #define PetscLLCreate(idx_start,lnk_max,nlnk,lnk,bt) \
500:   (PetscMalloc(nlnk*sizeof(PetscInt),&lnk) || PetscBTCreate(nlnk,bt) || PetscBTMemzero(nlnk,bt) || (lnk[idx_start] = lnk_max,0))

502: /*
503:   Add an index set into a sorted linked list
504:   Input Parameters:
505:     nidx      - number of input indices
506:     indices   - interger array
507:     idx_start - starting index of the list
508:     lnk       - linked list(an integer array) that is created
509:     bt        - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
510:   output Parameters:
511:     nlnk      - number of newly added indices
512:     lnk       - the sorted(increasing order) linked list containing new and non-redundate entries from indices
513:     bt        - updated PetscBT (bitarray) 
514: */
515: #define PetscLLAdd(nidx,indices,idx_start,nlnk,lnk,bt) 0;\
516: {\
517:   PetscInt _k,_entry,_location,_lnkdata;\
518:   nlnk     = 0;\
519:   _lnkdata = idx_start;\
520:   for (_k=0; _k<nidx; _k++){\
521:     _entry = indices[_k];\
522:     if (!PetscBTLookupSet(bt,_entry)){  /* new entry */\
523:       /* search for insertion location */\
524:       /* start from the beginning if _entry < previous _entry */\
525:       if (_k && _entry < _lnkdata) _lnkdata  = idx_start;\
526:       do {\
527:         _location = _lnkdata;\
528:         _lnkdata  = lnk[_location];\
529:       } while (_entry > _lnkdata);\
530:       /* insertion location is found, add entry into lnk */\
531:       lnk[_location] = _entry;\
532:       lnk[_entry]    = _lnkdata;\
533:       nlnk++;\
534:       _lnkdata = _entry; /* next search starts from here if next_entry > _entry */\
535:     }\
536:   }\
537: }

539: /*
540:   Add a permuted index set into a sorted linked list
541:   Input Parameters:
542:     nidx      - number of input indices
543:     indices   - interger array
544:     perm      - permutation of indices
545:     idx_start - starting index of the list
546:     lnk       - linked list(an integer array) that is created
547:     bt        - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
548:   output Parameters:
549:     nlnk      - number of newly added indices
550:     lnk       - the sorted(increasing order) linked list containing new and non-redundate entries from indices
551:     bt        - updated PetscBT (bitarray) 
552: */
553: #define PetscLLAddPerm(nidx,indices,perm,idx_start,nlnk,lnk,bt) 0;\
554: {\
555:   PetscInt _k,_entry,_location,_lnkdata;\
556:   nlnk     = 0;\
557:   _lnkdata = idx_start;\
558:   for (_k=0; _k<nidx; _k++){\
559:     _entry = perm[indices[_k]];\
560:     if (!PetscBTLookupSet(bt,_entry)){  /* new entry */\
561:       /* search for insertion location */\
562:       /* start from the beginning if _entry < previous _entry */\
563:       if (_k && _entry < _lnkdata) _lnkdata  = idx_start;\
564:       do {\
565:         _location = _lnkdata;\
566:         _lnkdata  = lnk[_location];\
567:       } while (_entry > _lnkdata);\
568:       /* insertion location is found, add entry into lnk */\
569:       lnk[_location] = _entry;\
570:       lnk[_entry]    = _lnkdata;\
571:       nlnk++;\
572:       _lnkdata = _entry; /* next search starts from here if next_entry > _entry */\
573:     }\
574:   }\
575: }

577: /*
578:   Add a SORTED index set into a sorted linked list
579:   Input Parameters:
580:     nidx      - number of input indices
581:     indices   - sorted interger array 
582:     idx_start - starting index of the list
583:     lnk       - linked list(an integer array) that is created
584:     bt        - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
585:   output Parameters:
586:     nlnk      - number of newly added indices
587:     lnk       - the sorted(increasing order) linked list containing new and non-redundate entries from indices
588:     bt        - updated PetscBT (bitarray) 
589: */
590: #define PetscLLAddSorted(nidx,indices,idx_start,nlnk,lnk,bt) 0;\
591: {\
592:   PetscInt _k,_entry,_location,_lnkdata;\
593:   nlnk      = 0;\
594:   _lnkdata  = idx_start;\
595:   for (_k=0; _k<nidx; _k++){\
596:     _entry = indices[_k];\
597:     if (!PetscBTLookupSet(bt,_entry)){  /* new entry */\
598:       /* search for insertion location */\
599:       do {\
600:         _location = _lnkdata;\
601:         _lnkdata  = lnk[_location];\
602:       } while (_entry > _lnkdata);\
603:       /* insertion location is found, add entry into lnk */\
604:       lnk[_location] = _entry;\
605:       lnk[_entry]    = _lnkdata;\
606:       nlnk++;\
607:       _lnkdata = _entry; /* next search starts from here */\
608:     }\
609:   }\
610: }

612: /*
613:   Add a SORTED index set into a sorted linked list used for LUFactorSymbolic()
614:   Same as PetscLLAddSorted() with an additional operation:
615:        count the number of input indices that are no larger than 'diag'
616:   Input Parameters:
617:     indices   - sorted interger array 
618:     idx_start - starting index of the list
619:     lnk       - linked list(an integer array) that is created
620:     bt        - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
621:     diag      - index of the active row in LUFactorSymbolic
622:     nzbd      - number of input indices with indices <= idx_start
623:   output Parameters:
624:     nlnk      - number of newly added indices
625:     lnk       - the sorted(increasing order) linked list containing new and non-redundate entries from indices
626:     bt        - updated PetscBT (bitarray) 
627:     im        - im[idx_start] =  num of entries with indices <= diag
628: */
629: #define PetscLLAddSortedLU(indices,idx_start,nlnk,lnk,bt,diag,nzbd,im) 0;\
630: {\
631:   PetscInt _k,_entry,_location,_lnkdata,_nidx;\
632:   nlnk     = 0;\
633:   _lnkdata = idx_start;\
634:   _nidx = im[idx_start] - nzbd; /* num of entries with idx_start < index <= diag */\
635:   for (_k=0; _k<_nidx; _k++){\
636:     _entry = indices[_k];\
637:     nzbd++;\
638:     if ( _entry== diag) im[idx_start] = nzbd;\
639:     if (!PetscBTLookupSet(bt,_entry)){  /* new entry */\
640:       /* search for insertion location */\
641:       do {\
642:         _location = _lnkdata;\
643:         _lnkdata  = lnk[_location];\
644:       } while (_entry > _lnkdata);\
645:       /* insertion location is found, add entry into lnk */\
646:       lnk[_location] = _entry;\
647:       lnk[_entry]    = _lnkdata;\
648:       nlnk++;\
649:       _lnkdata = _entry; /* next search starts from here */\
650:     }\
651:   }\
652: }

654: /*
655:   Copy data on the list into an array, then initialize the list 
656:   Input Parameters:
657:     idx_start - starting index of the list 
658:     lnk_max   - max value of lnk indicating the end of the list 
659:     nlnk      - number of data on the list to be copied
660:     lnk       - linked list
661:     bt        - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
662:   output Parameters:
663:     indices   - array that contains the copied data
664:     lnk       - linked list that is cleaned and initialize
665:     bt        - PetscBT (bitarray) with all bits set to false
666: */
667: #define PetscLLClean(idx_start,lnk_max,nlnk,lnk,indices,bt) 0;\
668: {\
669:   PetscInt _j,_idx=idx_start;\
670:   for (_j=0; _j<nlnk; _j++){\
671:     _idx = lnk[_idx];\
672:     *(indices+_j) = _idx;\
673:     PetscBTClear(bt,_idx);\
674:   }\
675:   lnk[idx_start] = lnk_max;\
676: }
677: /*
678:   Free memories used by the list
679: */
680: #define PetscLLDestroy(lnk,bt) (PetscFree(lnk) || PetscBTDestroy(bt))

682: /* Routines below are used for incomplete matrix factorization */
683: /* 
684:   Create and initialize a linked list and its levels
685:   Input Parameters:
686:     idx_start - starting index of the list
687:     lnk_max   - max value of lnk indicating the end of the list
688:     nlnk      - max length of the list
689:   Output Parameters:
690:     lnk       - list initialized
691:     lnk_lvl   - array of size nlnk for storing levels of lnk
692:     bt        - PetscBT (bitarray) with all bits set to false
693: */
694: #define PetscIncompleteLLCreate(idx_start,lnk_max,nlnk,lnk,lnk_lvl,bt)\
695:   (PetscMalloc(2*nlnk*sizeof(PetscInt),&lnk) || PetscBTCreate(nlnk,bt) || PetscBTMemzero(nlnk,bt) || (lnk[idx_start] = lnk_max,lnk_lvl = lnk + nlnk,0))

697: /*
698:   Initialize a sorted linked list used for ILU and ICC
699:   Input Parameters:
700:     nidx      - number of input idx
701:     idx       - interger array used for storing column indices
702:     idx_start - starting index of the list
703:     perm      - indices of an IS
704:     lnk       - linked list(an integer array) that is created
705:     lnklvl    - levels of lnk
706:     bt        - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
707:   output Parameters:
708:     nlnk     - number of newly added idx
709:     lnk      - the sorted(increasing order) linked list containing new and non-redundate entries from idx
710:     lnklvl   - levels of lnk
711:     bt       - updated PetscBT (bitarray) 
712: */
713: #define PetscIncompleteLLInit(nidx,idx,idx_start,perm,nlnk,lnk,lnklvl,bt) 0;\
714: {\
715:   PetscInt _k,_entry,_location,_lnkdata;\
716:   nlnk     = 0;\
717:   _lnkdata = idx_start;\
718:   for (_k=0; _k<nidx; _k++){\
719:     _entry = perm[idx[_k]];\
720:     if (!PetscBTLookupSet(bt,_entry)){  /* new entry */\
721:       /* search for insertion location */\
722:       if (_k && _entry < _lnkdata) _lnkdata  = idx_start;\
723:       do {\
724:         _location = _lnkdata;\
725:         _lnkdata  = lnk[_location];\
726:       } while (_entry > _lnkdata);\
727:       /* insertion location is found, add entry into lnk */\
728:       lnk[_location]  = _entry;\
729:       lnk[_entry]     = _lnkdata;\
730:       lnklvl[_entry] = 0;\
731:       nlnk++;\
732:       _lnkdata = _entry; /* next search starts from here if next_entry > _entry */\
733:     }\
734:   }\
735: }

737: /*
738:   Add a SORTED index set into a sorted linked list for ILU
739:   Input Parameters:
740:     nidx      - number of input indices
741:     idx       - sorted interger array used for storing column indices
742:     level     - level of fill, e.g., ICC(level)
743:     idxlvl    - level of idx 
744:     idx_start - starting index of the list
745:     lnk       - linked list(an integer array) that is created
746:     lnklvl    - levels of lnk
747:     bt        - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
748:     prow      - the row number of idx
749:   output Parameters:
750:     nlnk     - number of newly added idx
751:     lnk      - the sorted(increasing order) linked list containing new and non-redundate entries from idx
752:     lnklvl   - levels of lnk
753:     bt       - updated PetscBT (bitarray) 

755:   Note: the level of factor(i,j) is set as lvl(i,j) = min{ lvl(i,j), lvl(i,prow)+lvl(prow,j)+1)
756:         where idx = non-zero columns of U(prow,prow+1:n-1), prow<i
757: */
758: #define PetscILULLAddSorted(nidx,idx,level,idxlvl,idx_start,nlnk,lnk,lnklvl,bt,lnklvl_prow) 0;\
759: {\
760:   PetscInt _k,_entry,_location,_lnkdata,_incrlev,_lnklvl_prow=lnklvl[prow];\
761:   nlnk     = 0;\
762:   _lnkdata = idx_start;\
763:   for (_k=0; _k<nidx; _k++){\
764:     _incrlev = idxlvl[_k] + _lnklvl_prow + 1;\
765:     if (_incrlev > level) continue;\
766:     _entry = idx[_k];\
767:     if (!PetscBTLookupSet(bt,_entry)){  /* new entry */\
768:       /* search for insertion location */\
769:       do {\
770:         _location = _lnkdata;\
771:         _lnkdata  = lnk[_location];\
772:       } while (_entry > _lnkdata);\
773:       /* insertion location is found, add entry into lnk */\
774:       lnk[_location]  = _entry;\
775:       lnk[_entry]     = _lnkdata;\
776:       lnklvl[_entry] = _incrlev;\
777:       nlnk++;\
778:       _lnkdata = _entry; /* next search starts from here if next_entry > _entry */\
779:     } else { /* existing entry: update lnklvl */\
780:       if (lnklvl[_entry] > _incrlev) lnklvl[_entry] = _incrlev;\
781:     }\
782:   }\
783: }

785: /*
786:   Add a index set into a sorted linked list
787:   Input Parameters:
788:     nidx      - number of input idx
789:     idx   - interger array used for storing column indices
790:     level     - level of fill, e.g., ICC(level)
791:     idxlvl - level of idx 
792:     idx_start - starting index of the list
793:     lnk       - linked list(an integer array) that is created
794:     lnklvl   - levels of lnk
795:     bt        - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
796:   output Parameters:
797:     nlnk      - number of newly added idx
798:     lnk       - the sorted(increasing order) linked list containing new and non-redundate entries from idx
799:     lnklvl   - levels of lnk
800:     bt        - updated PetscBT (bitarray) 
801: */
802: #define PetscIncompleteLLAdd(nidx,idx,level,idxlvl,idx_start,nlnk,lnk,lnklvl,bt) 0;\
803: {\
804:   PetscInt _k,_entry,_location,_lnkdata,_incrlev;\
805:   nlnk     = 0;\
806:   _lnkdata = idx_start;\
807:   for (_k=0; _k<nidx; _k++){\
808:     _incrlev = idxlvl[_k] + 1;\
809:     if (_incrlev > level) continue;\
810:     _entry = idx[_k];\
811:     if (!PetscBTLookupSet(bt,_entry)){  /* new entry */\
812:       /* search for insertion location */\
813:       if (_k && _entry < _lnkdata) _lnkdata  = idx_start;\
814:       do {\
815:         _location = _lnkdata;\
816:         _lnkdata  = lnk[_location];\
817:       } while (_entry > _lnkdata);\
818:       /* insertion location is found, add entry into lnk */\
819:       lnk[_location]  = _entry;\
820:       lnk[_entry]     = _lnkdata;\
821:       lnklvl[_entry] = _incrlev;\
822:       nlnk++;\
823:       _lnkdata = _entry; /* next search starts from here if next_entry > _entry */\
824:     } else { /* existing entry: update lnklvl */\
825:       if (lnklvl[_entry] > _incrlev) lnklvl[_entry] = _incrlev;\
826:     }\
827:   }\
828: }

830: /*
831:   Add a SORTED index set into a sorted linked list
832:   Input Parameters:
833:     nidx      - number of input indices
834:     idx   - sorted interger array used for storing column indices
835:     level     - level of fill, e.g., ICC(level)
836:     idxlvl - level of idx 
837:     idx_start - starting index of the list
838:     lnk       - linked list(an integer array) that is created
839:     lnklvl    - levels of lnk
840:     bt        - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
841:   output Parameters:
842:     nlnk      - number of newly added idx
843:     lnk       - the sorted(increasing order) linked list containing new and non-redundate entries from idx
844:     lnklvl    - levels of lnk
845:     bt        - updated PetscBT (bitarray) 
846: */
847: #define PetscIncompleteLLAddSorted(nidx,idx,level,idxlvl,idx_start,nlnk,lnk,lnklvl,bt) 0;\
848: {\
849:   PetscInt _k,_entry,_location,_lnkdata,_incrlev;\
850:   nlnk = 0;\
851:   _lnkdata = idx_start;\
852:   for (_k=0; _k<nidx; _k++){\
853:     _incrlev = idxlvl[_k] + 1;\
854:     if (_incrlev > level) continue;\
855:     _entry = idx[_k];\
856:     if (!PetscBTLookupSet(bt,_entry)){  /* new entry */\
857:       /* search for insertion location */\
858:       do {\
859:         _location = _lnkdata;\
860:         _lnkdata  = lnk[_location];\
861:       } while (_entry > _lnkdata);\
862:       /* insertion location is found, add entry into lnk */\
863:       lnk[_location] = _entry;\
864:       lnk[_entry]    = _lnkdata;\
865:       lnklvl[_entry] = _incrlev;\
866:       nlnk++;\
867:       _lnkdata = _entry; /* next search starts from here */\
868:     } else { /* existing entry: update lnklvl */\
869:       if (lnklvl[_entry] > _incrlev) lnklvl[_entry] = _incrlev;\
870:     }\
871:   }\
872: }

874: /*
875:   Add a SORTED index set into a sorted linked list for ICC
876:   Input Parameters:
877:     nidx      - number of input indices
878:     idx       - sorted interger array used for storing column indices
879:     level     - level of fill, e.g., ICC(level)
880:     idxlvl    - level of idx 
881:     idx_start - starting index of the list
882:     lnk       - linked list(an integer array) that is created
883:     lnklvl    - levels of lnk
884:     bt        - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
885:     idxlvl_prow - idxlvl[prow], where prow is the row number of the idx
886:   output Parameters:
887:     nlnk   - number of newly added indices
888:     lnk    - the sorted(increasing order) linked list containing new and non-redundate entries from idx
889:     lnklvl - levels of lnk
890:     bt     - updated PetscBT (bitarray) 
891:   Note: the level of U(i,j) is set as lvl(i,j) = min{ lvl(i,j), lvl(prow,i)+lvl(prow,j)+1)
892:         where idx = non-zero columns of U(prow,prow+1:n-1), prow<i
893: */
894: #define PetscICCLLAddSorted(nidx,idx,level,idxlvl,idx_start,nlnk,lnk,lnklvl,bt,idxlvl_prow) 0;\
895: {\
896:   PetscInt _k,_entry,_location,_lnkdata,_incrlev;\
897:   nlnk = 0;\
898:   _lnkdata = idx_start;\
899:   for (_k=0; _k<nidx; _k++){\
900:     _incrlev = idxlvl[_k] + idxlvl_prow + 1;\
901:     if (_incrlev > level) continue;\
902:     _entry = idx[_k];\
903:     if (!PetscBTLookupSet(bt,_entry)){  /* new entry */\
904:       /* search for insertion location */\
905:       do {\
906:         _location = _lnkdata;\
907:         _lnkdata  = lnk[_location];\
908:       } while (_entry > _lnkdata);\
909:       /* insertion location is found, add entry into lnk */\
910:       lnk[_location] = _entry;\
911:       lnk[_entry]    = _lnkdata;\
912:       lnklvl[_entry] = _incrlev;\
913:       nlnk++;\
914:       _lnkdata = _entry; /* next search starts from here */\
915:     } else { /* existing entry: update lnklvl */\
916:       if (lnklvl[_entry] > _incrlev) lnklvl[_entry] = _incrlev;\
917:     }\
918:   }\
919: }

921: /*
922:   Copy data on the list into an array, then initialize the list 
923:   Input Parameters:
924:     idx_start - starting index of the list 
925:     lnk_max   - max value of lnk indicating the end of the list 
926:     nlnk      - number of data on the list to be copied
927:     lnk       - linked list
928:     lnklvl    - level of lnk
929:     bt        - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
930:   output Parameters:
931:     indices - array that contains the copied data
932:     lnk     - linked list that is cleaned and initialize
933:     lnklvl  - level of lnk that is reinitialized 
934:     bt      - PetscBT (bitarray) with all bits set to false
935: */
936: #define PetscIncompleteLLClean(idx_start,lnk_max,nlnk,lnk,lnklvl,indices,indiceslvl,bt) 0;\
937: {\
938:   PetscInt _j,_idx=idx_start;\
939:   for (_j=0; _j<nlnk; _j++){\
940:     _idx = lnk[_idx];\
941:     *(indices+_j) = _idx;\
942:     *(indiceslvl+_j) = lnklvl[_idx];\
943:     lnklvl[_idx] = -1;\
944:     PetscBTClear(bt,_idx);\
945:   }\
946:   lnk[idx_start] = lnk_max;\
947: }
948: /*
949:   Free memories used by the list
950: */
951: #define PetscIncompleteLLDestroy(lnk,bt) (PetscFree(lnk) || PetscBTDestroy(bt))


965: #endif