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_LOCAL 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
384: -  idiag - index of diagonals in array aval

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

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

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

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

447:    Collective on Mat

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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


966: #endif