Actual source code: matimpl.h

petsc-3.10.3 2018-12-18
Report Typos and Errors

  2: #ifndef __MATIMPL_H

  5:  #include <petscmat.h>
  6:  #include <petscmatcoarsen.h>
  7:  #include <petsc/private/petscimpl.h>

  9: PETSC_EXTERN PetscBool MatRegisterAllCalled;
 10: PETSC_EXTERN PetscBool MatSeqAIJRegisterAllCalled;
 11: PETSC_EXTERN PetscBool MatOrderingRegisterAllCalled;
 12: PETSC_EXTERN PetscBool MatColoringRegisterAllCalled;
 13: PETSC_EXTERN PetscBool MatPartitioningRegisterAllCalled;
 14: PETSC_EXTERN PetscBool MatCoarsenRegisterAllCalled;
 15: PETSC_EXTERN PetscErrorCode MatRegisterAll(void);
 16: PETSC_EXTERN PetscErrorCode MatOrderingRegisterAll(void);
 17: PETSC_EXTERN PetscErrorCode MatColoringRegisterAll(void);
 18: PETSC_EXTERN PetscErrorCode MatPartitioningRegisterAll(void);
 19: PETSC_EXTERN PetscErrorCode MatCoarsenRegisterAll(void);
 20: PETSC_EXTERN PetscErrorCode MatSeqAIJRegisterAll(void);

 22: /*
 23:   This file defines the parts of the matrix data structure that are
 24:   shared by all matrix types.
 25: */

 27: /*
 28:     If you add entries here also add them to the MATOP enum
 29:     in include/petscmat.h and include/petsc/finclude/petscmat.h
 30: */
 31: typedef struct _MatOps *MatOps;
 32: struct _MatOps {
 33:   /* 0*/
 34:   PetscErrorCode (*setvalues)(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
 35:   PetscErrorCode (*getrow)(Mat,PetscInt,PetscInt *,PetscInt*[],PetscScalar*[]);
 36:   PetscErrorCode (*restorerow)(Mat,PetscInt,PetscInt *,PetscInt *[],PetscScalar *[]);
 37:   PetscErrorCode (*mult)(Mat,Vec,Vec);
 38:   PetscErrorCode (*multadd)(Mat,Vec,Vec,Vec);
 39:   /* 5*/
 40:   PetscErrorCode (*multtranspose)(Mat,Vec,Vec);
 41:   PetscErrorCode (*multtransposeadd)(Mat,Vec,Vec,Vec);
 42:   PetscErrorCode (*solve)(Mat,Vec,Vec);
 43:   PetscErrorCode (*solveadd)(Mat,Vec,Vec,Vec);
 44:   PetscErrorCode (*solvetranspose)(Mat,Vec,Vec);
 45:   /*10*/
 46:   PetscErrorCode (*solvetransposeadd)(Mat,Vec,Vec,Vec);
 47:   PetscErrorCode (*lufactor)(Mat,IS,IS,const MatFactorInfo*);
 48:   PetscErrorCode (*choleskyfactor)(Mat,IS,const MatFactorInfo*);
 49:   PetscErrorCode (*sor)(Mat,Vec,PetscReal,MatSORType,PetscReal,PetscInt,PetscInt,Vec);
 50:   PetscErrorCode (*transpose)(Mat,MatReuse,Mat *);
 51:   /*15*/
 52:   PetscErrorCode (*getinfo)(Mat,MatInfoType,MatInfo*);
 53:   PetscErrorCode (*equal)(Mat,Mat,PetscBool  *);
 54:   PetscErrorCode (*getdiagonal)(Mat,Vec);
 55:   PetscErrorCode (*diagonalscale)(Mat,Vec,Vec);
 56:   PetscErrorCode (*norm)(Mat,NormType,PetscReal*);
 57:   /*20*/
 58:   PetscErrorCode (*assemblybegin)(Mat,MatAssemblyType);
 59:   PetscErrorCode (*assemblyend)(Mat,MatAssemblyType);
 60:   PetscErrorCode (*setoption)(Mat,MatOption,PetscBool );
 61:   PetscErrorCode (*zeroentries)(Mat);
 62:   /*24*/
 63:   PetscErrorCode (*zerorows)(Mat,PetscInt,const PetscInt[],PetscScalar,Vec,Vec);
 64:   PetscErrorCode (*lufactorsymbolic)(Mat,Mat,IS,IS,const MatFactorInfo*);
 65:   PetscErrorCode (*lufactornumeric)(Mat,Mat,const MatFactorInfo*);
 66:   PetscErrorCode (*choleskyfactorsymbolic)(Mat,Mat,IS,const MatFactorInfo*);
 67:   PetscErrorCode (*choleskyfactornumeric)(Mat,Mat,const MatFactorInfo*);
 68:   /*29*/
 69:   PetscErrorCode (*setup)(Mat);
 70:   PetscErrorCode (*ilufactorsymbolic)(Mat,Mat,IS,IS,const MatFactorInfo*);
 71:   PetscErrorCode (*iccfactorsymbolic)(Mat,Mat,IS,const MatFactorInfo*);
 72:   PetscErrorCode (*getdiagonalblock)(Mat,Mat*);
 73:   PetscErrorCode (*placeholder_33)(Mat);
 74:   /*34*/
 75:   PetscErrorCode (*duplicate)(Mat,MatDuplicateOption,Mat*);
 76:   PetscErrorCode (*forwardsolve)(Mat,Vec,Vec);
 77:   PetscErrorCode (*backwardsolve)(Mat,Vec,Vec);
 78:   PetscErrorCode (*ilufactor)(Mat,IS,IS,const MatFactorInfo*);
 79:   PetscErrorCode (*iccfactor)(Mat,IS,const MatFactorInfo*);
 80:   /*39*/
 81:   PetscErrorCode (*axpy)(Mat,PetscScalar,Mat,MatStructure);
 82:   PetscErrorCode (*createsubmatrices)(Mat,PetscInt,const IS[],const IS[],MatReuse,Mat *[]);
 83:   PetscErrorCode (*increaseoverlap)(Mat,PetscInt,IS[],PetscInt);
 84:   PetscErrorCode (*getvalues)(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],PetscScalar []);
 85:   PetscErrorCode (*copy)(Mat,Mat,MatStructure);
 86:   /*44*/
 87:   PetscErrorCode (*getrowmax)(Mat,Vec,PetscInt[]);
 88:   PetscErrorCode (*scale)(Mat,PetscScalar);
 89:   PetscErrorCode (*shift)(Mat,PetscScalar);
 90:   PetscErrorCode (*diagonalset)(Mat,Vec,InsertMode);
 91:   PetscErrorCode (*zerorowscolumns)(Mat,PetscInt,const PetscInt[],PetscScalar,Vec,Vec);
 92:   /*49*/
 93:   PetscErrorCode (*setrandom)(Mat,PetscRandom);
 94:   PetscErrorCode (*getrowij)(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt*,const PetscInt *[],const PetscInt *[],PetscBool  *);
 95:   PetscErrorCode (*restorerowij)(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt *,const PetscInt *[],const PetscInt *[],PetscBool  *);
 96:   PetscErrorCode (*getcolumnij)(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt*,const PetscInt *[],const PetscInt *[],PetscBool  *);
 97:   PetscErrorCode (*restorecolumnij)(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt*,const PetscInt *[],const PetscInt *[],PetscBool  *);
 98:   /*54*/
 99:   PetscErrorCode (*fdcoloringcreate)(Mat,ISColoring,MatFDColoring);
100:   PetscErrorCode (*coloringpatch)(Mat,PetscInt,PetscInt,ISColoringValue[],ISColoring*);
101:   PetscErrorCode (*setunfactored)(Mat);
102:   PetscErrorCode (*permute)(Mat,IS,IS,Mat*);
103:   PetscErrorCode (*setvaluesblocked)(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
104:   /*59*/
105:   PetscErrorCode (*createsubmatrix)(Mat,IS,IS,MatReuse,Mat*);
106:   PetscErrorCode (*destroy)(Mat);
107:   PetscErrorCode (*view)(Mat,PetscViewer);
108:   PetscErrorCode (*convertfrom)(Mat, MatType,MatReuse,Mat*);
109:   PetscErrorCode (*matmatmult)(Mat,Mat,Mat,MatReuse,PetscReal,Mat*);
110:   /*64*/
111:   PetscErrorCode (*matmatmultsymbolic)(Mat,Mat,Mat,PetscReal,Mat*);
112:   PetscErrorCode (*matmatmultnumeric)(Mat,Mat,Mat,Mat);
113:   PetscErrorCode (*setlocaltoglobalmapping)(Mat,ISLocalToGlobalMapping,ISLocalToGlobalMapping);
114:   PetscErrorCode (*setvalueslocal)(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
115:   PetscErrorCode (*zerorowslocal)(Mat,PetscInt,const PetscInt[],PetscScalar,Vec,Vec);
116:   /*69*/
117:   PetscErrorCode (*getrowmaxabs)(Mat,Vec,PetscInt[]);
118:   PetscErrorCode (*getrowminabs)(Mat,Vec,PetscInt[]);
119:   PetscErrorCode (*convert)(Mat, MatType,MatReuse,Mat*);
120:   PetscErrorCode (*hasoperation)(Mat,MatOperation,PetscBool*);
121:   PetscErrorCode (*placeholder_73)(Mat,void*);
122:   /*74*/
123:   PetscErrorCode (*setvaluesadifor)(Mat,PetscInt,void*);
124:   PetscErrorCode (*fdcoloringapply)(Mat,MatFDColoring,Vec,void*);
125:   PetscErrorCode (*setfromoptions)(PetscOptionItems*,Mat);
126:   PetscErrorCode (*multconstrained)(Mat,Vec,Vec);
127:   PetscErrorCode (*multtransposeconstrained)(Mat,Vec,Vec);
128:   /*79*/
129:   PetscErrorCode (*findzerodiagonals)(Mat,IS*);
130:   PetscErrorCode (*mults)(Mat, Vecs, Vecs);
131:   PetscErrorCode (*solves)(Mat, Vecs, Vecs);
132:   PetscErrorCode (*getinertia)(Mat,PetscInt*,PetscInt*,PetscInt*);
133:   PetscErrorCode (*load)(Mat, PetscViewer);
134:   /*84*/
135:   PetscErrorCode (*issymmetric)(Mat,PetscReal,PetscBool *);
136:   PetscErrorCode (*ishermitian)(Mat,PetscReal,PetscBool *);
137:   PetscErrorCode (*isstructurallysymmetric)(Mat,PetscBool *);
138:   PetscErrorCode (*setvaluesblockedlocal)(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
139:   PetscErrorCode (*getvecs)(Mat,Vec*,Vec*);
140:   /*89*/
141:   PetscErrorCode (*matmult)(Mat,Mat,MatReuse,PetscReal,Mat*);
142:   PetscErrorCode (*matmultsymbolic)(Mat,Mat,PetscReal,Mat*);
143:   PetscErrorCode (*matmultnumeric)(Mat,Mat,Mat);
144:   PetscErrorCode (*ptap)(Mat,Mat,MatReuse,PetscReal,Mat*);
145:   PetscErrorCode (*ptapsymbolic)(Mat,Mat,PetscReal,Mat*); /* double dispatch wrapper routine */
146:   /*94*/
147:   PetscErrorCode (*ptapnumeric)(Mat,Mat,Mat);             /* double dispatch wrapper routine */
148:   PetscErrorCode (*mattransposemult)(Mat,Mat,MatReuse,PetscReal,Mat*);
149:   PetscErrorCode (*mattransposemultsymbolic)(Mat,Mat,PetscReal,Mat*);
150:   PetscErrorCode (*mattransposemultnumeric)(Mat,Mat,Mat);
151:   PetscErrorCode (*placeholder_98)(Mat);
152:   /*99*/
153:   PetscErrorCode (*placeholder_99)(Mat);
154:   PetscErrorCode (*placeholder_100)(Mat);
155:   PetscErrorCode (*placeholder_101)(Mat);
156:   PetscErrorCode (*conjugate)(Mat);                              /* complex conjugate */
157:   PetscErrorCode (*viewnative)(Mat,PetscViewer);
158:   /*104*/
159:   PetscErrorCode (*setvaluesrow)(Mat,PetscInt,const PetscScalar[]);
160:   PetscErrorCode (*realpart)(Mat);
161:   PetscErrorCode (*imaginarypart)(Mat);
162:   PetscErrorCode (*getrowuppertriangular)(Mat);
163:   PetscErrorCode (*restorerowuppertriangular)(Mat);
164:   /*109*/
165:   PetscErrorCode (*matsolve)(Mat,Mat,Mat);
166:   PetscErrorCode (*matsolvetranspose)(Mat,Mat,Mat);
167:   PetscErrorCode (*getrowmin)(Mat,Vec,PetscInt[]);
168:   PetscErrorCode (*getcolumnvector)(Mat,Vec,PetscInt);
169:   PetscErrorCode (*missingdiagonal)(Mat,PetscBool *,PetscInt*);
170:   /*114*/
171:   PetscErrorCode (*getseqnonzerostructure)(Mat,Mat *);
172:   PetscErrorCode (*create)(Mat);
173:   PetscErrorCode (*getghosts)(Mat,PetscInt*,const PetscInt *[]);
174:   PetscErrorCode (*getlocalsubmatrix)(Mat,IS,IS,Mat*);
175:   PetscErrorCode (*restorelocalsubmatrix)(Mat,IS,IS,Mat*);
176:   /*119*/
177:   PetscErrorCode (*multdiagonalblock)(Mat,Vec,Vec);
178:   PetscErrorCode (*hermitiantranspose)(Mat,MatReuse,Mat*);
179:   PetscErrorCode (*multhermitiantranspose)(Mat,Vec,Vec);
180:   PetscErrorCode (*multhermitiantransposeadd)(Mat,Vec,Vec,Vec);
181:   PetscErrorCode (*getmultiprocblock)(Mat,MPI_Comm,MatReuse,Mat*);
182:   /*124*/
183:   PetscErrorCode (*findnonzerorows)(Mat,IS*);
184:   PetscErrorCode (*getcolumnnorms)(Mat,NormType,PetscReal*);
185:   PetscErrorCode (*invertblockdiagonal)(Mat,const PetscScalar**);
186:   PetscErrorCode (*invertvariableblockdiagonal)(Mat,PetscInt,const PetscInt*,PetscScalar*);
187:   PetscErrorCode (*createsubmatricesmpi)(Mat,PetscInt,const IS[], const IS[], MatReuse, Mat**);
188:   /*129*/
189:   PetscErrorCode (*setvaluesbatch)(Mat,PetscInt,PetscInt,PetscInt*,const PetscScalar*);
190:   PetscErrorCode (*transposematmult)(Mat,Mat,MatReuse,PetscReal,Mat*);
191:   PetscErrorCode (*transposematmultsymbolic)(Mat,Mat,PetscReal,Mat*);
192:   PetscErrorCode (*transposematmultnumeric)(Mat,Mat,Mat);
193:   PetscErrorCode (*transposecoloringcreate)(Mat,ISColoring,MatTransposeColoring);
194:   /*134*/
195:   PetscErrorCode (*transcoloringapplysptoden)(MatTransposeColoring,Mat,Mat);
196:   PetscErrorCode (*transcoloringapplydentosp)(MatTransposeColoring,Mat,Mat);
197:   PetscErrorCode (*rart)(Mat,Mat,MatReuse,PetscReal,Mat*);
198:   PetscErrorCode (*rartsymbolic)(Mat,Mat,PetscReal,Mat*); /* double dispatch wrapper routine */
199:   PetscErrorCode (*rartnumeric)(Mat,Mat,Mat);             /* double dispatch wrapper routine */
200:   /*139*/
201:   PetscErrorCode (*setblocksizes)(Mat,PetscInt,PetscInt);
202:   PetscErrorCode (*aypx)(Mat,PetscScalar,Mat,MatStructure);
203:   PetscErrorCode (*residual)(Mat,Vec,Vec,Vec);
204:   PetscErrorCode (*fdcoloringsetup)(Mat,ISColoring,MatFDColoring);
205:   PetscErrorCode (*findoffblockdiagonalentries)(Mat,IS*);
206:   /*144*/
207:   PetscErrorCode (*creatempimatconcatenateseqmat)(MPI_Comm,Mat,PetscInt,MatReuse,Mat*);
208:   PetscErrorCode (*destroysubmatrices)(PetscInt,Mat*[]);
209:   PetscErrorCode (*mattransposesolve)(Mat,Mat,Mat);
210: };
211: /*
212:     If you add MatOps entries above also add them to the MATOP enum
213:     in include/petscmat.h and include/petsc/finclude/petscmat.h
214: */

216:  #include <petscsys.h>
217: PETSC_EXTERN PetscErrorCode MatRegisterOp(MPI_Comm, const char[], PetscVoidFunction, const char[], PetscInt, ...);
218: PETSC_EXTERN PetscErrorCode MatQueryOp(MPI_Comm, PetscVoidFunction*, const char[], PetscInt, ...);

220: typedef struct _p_MatBaseName* MatBaseName;
221: struct _p_MatBaseName {
222:   char        *bname,*sname,*mname;
223:   MatBaseName next;
224: };

226: PETSC_EXTERN MatBaseName MatBaseNameList;

228: /*
229:    Utility private matrix routines
230: */
231: PETSC_INTERN PetscErrorCode MatFindNonzeroRowsOrCols_Basic(Mat,PetscBool,PetscReal,IS*);
232: PETSC_INTERN PetscErrorCode MatConvert_Basic(Mat,MatType,MatReuse,Mat*);
233: PETSC_INTERN PetscErrorCode MatCopy_Basic(Mat,Mat,MatStructure);
234: PETSC_INTERN PetscErrorCode MatDiagonalSet_Default(Mat,Vec,InsertMode);

236: #if defined(PETSC_USE_DEBUG)
237: #  define MatCheckPreallocated(A,arg) do {                              \
238:     if (PetscUnlikely(!(A)->preallocated)) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call MatXXXSetPreallocation() or MatSetUp() on argument %D \"%s\" before %s()",(arg),#A,PETSC_FUNCTION_NAME); \
239:   } while (0)
240: #else
241: #  define MatCheckPreallocated(A,arg) do {} while (0)
242: #endif

244: /*
245:   The stash is used to temporarily store inserted matrix values that
246:   belong to another processor. During the assembly phase the stashed
247:   values are moved to the correct processor and
248: */

250: typedef struct _MatStashSpace *PetscMatStashSpace;

252: struct _MatStashSpace {
253:   PetscMatStashSpace next;
254:   PetscScalar        *space_head,*val;
255:   PetscInt           *idx,*idy;
256:   PetscInt           total_space_size;
257:   PetscInt           local_used;
258:   PetscInt           local_remaining;
259: };

261: PETSC_EXTERN PetscErrorCode PetscMatStashSpaceGet(PetscInt,PetscInt,PetscMatStashSpace *);
262: PETSC_EXTERN PetscErrorCode PetscMatStashSpaceContiguous(PetscInt,PetscMatStashSpace *,PetscScalar *,PetscInt *,PetscInt *);
263: PETSC_EXTERN PetscErrorCode PetscMatStashSpaceDestroy(PetscMatStashSpace*);

265: typedef struct {
266:   PetscInt    count;
267: } MatStashHeader;

269: typedef struct {
270:   void        *buffer;          /* Of type blocktype, dynamically constructed  */
271:   PetscInt    count;
272:   char        pending;
273: } MatStashFrame;

275: typedef struct _MatStash MatStash;
276: struct _MatStash {
277:   PetscInt      nmax;                   /* maximum stash size */
278:   PetscInt      umax;                   /* user specified max-size */
279:   PetscInt      oldnmax;                /* the nmax value used previously */
280:   PetscInt      n;                      /* stash size */
281:   PetscInt      bs;                     /* block size of the stash */
282:   PetscInt      reallocs;               /* preserve the no of mallocs invoked */
283:   PetscMatStashSpace space_head,space;  /* linked list to hold stashed global row/column numbers and matrix values */

285:   PetscErrorCode (*ScatterBegin)(Mat,MatStash*,PetscInt*);
286:   PetscErrorCode (*ScatterGetMesg)(MatStash*,PetscMPIInt*,PetscInt**,PetscInt**,PetscScalar**,PetscInt*);
287:   PetscErrorCode (*ScatterEnd)(MatStash*);
288:   PetscErrorCode (*ScatterDestroy)(MatStash*);

290:   /* The following variables are used for communication */
291:   MPI_Comm      comm;
292:   PetscMPIInt   size,rank;
293:   PetscMPIInt   tag1,tag2;
294:   MPI_Request   *send_waits;            /* array of send requests */
295:   MPI_Request   *recv_waits;            /* array of receive requests */
296:   MPI_Status    *send_status;           /* array of send status */
297:   PetscInt      nsends,nrecvs;          /* numbers of sends and receives */
298:   PetscScalar   *svalues;               /* sending data */
299:   PetscInt      *sindices;
300:   PetscScalar   **rvalues;              /* receiving data (values) */
301:   PetscInt      **rindices;             /* receiving data (indices) */
302:   PetscInt      nprocessed;             /* number of messages already processed */
303:   PetscMPIInt   *flg_v;                 /* indicates what messages have arrived so far and from whom */
304:   PetscBool     reproduce;
305:   PetscInt      reproduce_count;

307:   /* The following variables are used for BTS communication */
308:   PetscBool      subset_off_proc; /* Subsequent assemblies will set a subset (perhaps equal) of off-process entries set on first assembly */
309:   PetscBool      use_status;      /* Use MPI_Status to determine number of items in each message */
310:   PetscMPIInt    nsendranks;
311:   PetscMPIInt    nrecvranks;
312:   PetscMPIInt    *sendranks;
313:   PetscMPIInt    *recvranks;
314:   MatStashHeader *sendhdr,*recvhdr;
315:   MatStashFrame  *sendframes;   /* pointers to the main messages */
316:   MatStashFrame  *recvframes;
317:   MatStashFrame  *recvframe_active;
318:   PetscInt       recvframe_i;     /* index of block within active frame */
319:   PetscMPIInt    recvframe_count; /* Count actually sent for current frame */
320:   PetscInt       recvcount;       /* Number of receives processed so far */
321:   PetscMPIInt    *some_indices;   /* From last call to MPI_Waitsome */
322:   MPI_Status     *some_statuses;  /* Statuses from last call to MPI_Waitsome */
323:   PetscMPIInt    some_count;      /* Number of requests completed in last call to MPI_Waitsome */
324:   PetscMPIInt    some_i;          /* Index of request currently being processed */
325:   MPI_Request    *sendreqs;
326:   MPI_Request    *recvreqs;
327:   PetscSegBuffer segsendblocks;
328:   PetscSegBuffer segrecvframe;
329:   PetscSegBuffer segrecvblocks;
330:   MPI_Datatype   blocktype;
331:   size_t         blocktype_size;
332:   InsertMode     *insertmode;   /* Pointer to check mat->insertmode and set upon message arrival in case no local values have been set. */
333: };

335: PETSC_INTERN PetscErrorCode MatStashCreate_Private(MPI_Comm,PetscInt,MatStash*);
336: PETSC_INTERN PetscErrorCode MatStashDestroy_Private(MatStash*);
337: PETSC_INTERN PetscErrorCode MatStashScatterEnd_Private(MatStash*);
338: PETSC_INTERN PetscErrorCode MatStashSetInitialSize_Private(MatStash*,PetscInt);
339: PETSC_INTERN PetscErrorCode MatStashGetInfo_Private(MatStash*,PetscInt*,PetscInt*);
340: PETSC_INTERN PetscErrorCode MatStashValuesRow_Private(MatStash*,PetscInt,PetscInt,const PetscInt[],const PetscScalar[],PetscBool );
341: PETSC_INTERN PetscErrorCode MatStashValuesCol_Private(MatStash*,PetscInt,PetscInt,const PetscInt[],const PetscScalar[],PetscInt,PetscBool );
342: PETSC_INTERN PetscErrorCode MatStashValuesRowBlocked_Private(MatStash*,PetscInt,PetscInt,const PetscInt[],const PetscScalar[],PetscInt,PetscInt,PetscInt);
343: PETSC_INTERN PetscErrorCode MatStashValuesColBlocked_Private(MatStash*,PetscInt,PetscInt,const PetscInt[],const PetscScalar[],PetscInt,PetscInt,PetscInt);
344: PETSC_INTERN PetscErrorCode MatStashScatterBegin_Private(Mat,MatStash*,PetscInt*);
345: PETSC_INTERN PetscErrorCode MatStashScatterGetMesg_Private(MatStash*,PetscMPIInt*,PetscInt**,PetscInt**,PetscScalar**,PetscInt*);
346: PETSC_INTERN PetscErrorCode MatGetInfo_External(Mat,MatInfoType,MatInfo*);

348: typedef struct {
349:   PetscInt   dim;
350:   PetscInt   dims[4];
351:   PetscInt   starts[4];
352:   PetscBool  noc;        /* this is a single component problem, hence user will not set MatStencil.c */
353: } MatStencilInfo;

355: /* Info about using compressed row format */
356: typedef struct {
357:   PetscBool  use;                           /* indicates compressed rows have been checked and will be used */
358:   PetscInt   nrows;                         /* number of non-zero rows */
359:   PetscInt   *i;                            /* compressed row pointer  */
360:   PetscInt   *rindex;                       /* compressed row index               */
361: } Mat_CompressedRow;
362: PETSC_EXTERN PetscErrorCode MatCheckCompressedRow(Mat,PetscInt,Mat_CompressedRow*,PetscInt*,PetscInt,PetscReal);

364: typedef struct { /* used by MatCreateRedundantMatrix() for reusing matredundant */
365:   PetscInt     nzlocal,nsends,nrecvs;
366:   PetscMPIInt  *send_rank,*recv_rank;
367:   PetscInt     *sbuf_nz,*rbuf_nz,*sbuf_j,**rbuf_j;
368:   PetscScalar  *sbuf_a,**rbuf_a;
369:   MPI_Comm     subcomm;   /* when user does not provide a subcomm */
370:   IS           isrow,iscol;
371:   Mat          *matseq;
372: } Mat_Redundant;

374: struct _p_Mat {
375:   PETSCHEADER(struct _MatOps);
376:   PetscLayout            rmap,cmap;
377:   void                   *data;            /* implementation-specific data */
378:   MatFactorType          factortype;       /* MAT_FACTOR_LU, ILU, CHOLESKY or ICC */
379:   PetscBool              assembled;        /* is the matrix assembled? */
380:   PetscBool              was_assembled;    /* new values inserted into assembled mat */
381:   PetscInt               num_ass;          /* number of times matrix has been assembled */
382:   PetscObjectState       nonzerostate;     /* each time new nonzeros locations are introduced into the matrix this is updated */
383:   MatInfo                info;             /* matrix information */
384:   InsertMode             insertmode;       /* have values been inserted in matrix or added? */
385:   MatStash               stash,bstash;     /* used for assembling off-proc mat emements */
386:   MatNullSpace           nullsp;           /* null space (operator is singular) */
387:   MatNullSpace           transnullsp;      /* null space of transpose of operator */
388:   MatNullSpace           nearnullsp;       /* near null space to be used by multigrid methods */
389:   PetscInt               congruentlayouts; /* are the rows and columns layouts congruent? */
390:   PetscBool              preallocated;
391:   MatStencilInfo         stencil;          /* information for structured grid */
392:   PetscBool              symmetric,hermitian,structurally_symmetric,spd;
393:   PetscBool              symmetric_set,hermitian_set,structurally_symmetric_set,spd_set; /* if true, then corresponding flag is correct*/
394:   PetscBool              symmetric_eternal;
395:   PetscBool              nooffprocentries,nooffproczerorows;
396:   PetscBool              subsetoffprocentries;
397:   PetscBool              submat_singleis; /* for efficient PCSetUP_ASM() */
398:   PetscBool              structure_only;
399: #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_VECCUDA)
400:   PetscOffloadFlag       valid_GPU_matrix; /* flag pointing to the matrix on the gpu*/
401: #endif
402:   void                   *spptr;          /* pointer for special library like SuperLU */
403:   char                   *solvertype;
404:   PetscBool              checksymmetryonassembly,checknullspaceonassembly;
405:   PetscReal              checksymmetrytol;
406:   Mat                    schur;             /* Schur complement matrix */
407:   MatFactorSchurStatus   schur_status;      /* status of the Schur complement matrix */
408:   Mat_Redundant          *redundant;        /* used by MatCreateRedundantMatrix() */
409:   PetscBool              erroriffailure;    /* Generate an error if detected (for example a zero pivot) instead of returning */
410:   MatFactorError         factorerrortype;               /* type of error in factorization */
411:   PetscReal              factorerror_zeropivot_value;   /* If numerical zero pivot was detected this is the computed value */
412:   PetscInt               factorerror_zeropivot_row;     /* Row where zero pivot was detected */
413:   PetscInt               nblocks,*bsizes;   /* support for MatSetVariableBlockSizes() */
414:   char                   *defaultvectype;
415: };

417: PETSC_INTERN PetscErrorCode MatAXPY_Basic(Mat,PetscScalar,Mat,MatStructure);
418: PETSC_INTERN PetscErrorCode MatAXPY_BasicWithPreallocation(Mat,Mat,PetscScalar,Mat,MatStructure);

420: /*
421:     Utility for MatFactor (Schur complement)
422: */
423: PETSC_INTERN PetscErrorCode MatFactorFactorizeSchurComplement_Private(Mat);
424: PETSC_INTERN PetscErrorCode MatFactorInvertSchurComplement_Private(Mat);
425: PETSC_INTERN PetscErrorCode MatFactorUpdateSchurStatus_Private(Mat);
426: PETSC_INTERN PetscErrorCode MatFactorSetUpInPlaceSchur_Private(Mat);

428: /*
429:     Utility for MatZeroRows
430: */
431: PETSC_INTERN PetscErrorCode MatZeroRowsMapLocal_Private(Mat,PetscInt,const PetscInt*,PetscInt*,PetscInt**);

433: /*
434:     Object for partitioning graphs
435: */

437: typedef struct _MatPartitioningOps *MatPartitioningOps;
438: struct _MatPartitioningOps {
439:   PetscErrorCode (*apply)(MatPartitioning,IS*);
440:   PetscErrorCode (*applynd)(MatPartitioning,IS*);
441:   PetscErrorCode (*setfromoptions)(PetscOptionItems*,MatPartitioning);
442:   PetscErrorCode (*destroy)(MatPartitioning);
443:   PetscErrorCode (*view)(MatPartitioning,PetscViewer);
444: };

446: struct _p_MatPartitioning {
447:   PETSCHEADER(struct _MatPartitioningOps);
448:   Mat         adj;
449:   PetscInt    *vertex_weights;
450:   PetscReal   *part_weights;
451:   PetscInt    n;                                 /* number of partitions */
452:   void        *data;
453:   PetscInt    setupcalled;
454: };

456: /* needed for parallel nested dissection by ParMetis and PTSCOTCH */
457: PETSC_INTERN PetscErrorCode MatPartitioningSizesToSep_Private(PetscInt,PetscInt[],PetscInt[],PetscInt[]);

459: /*
460:     Object for coarsen graphs
461: */
462: typedef struct _MatCoarsenOps *MatCoarsenOps;
463: struct _MatCoarsenOps {
464:   PetscErrorCode (*apply)(MatCoarsen);
465:   PetscErrorCode (*setfromoptions)(PetscOptionItems*,MatCoarsen);
466:   PetscErrorCode (*destroy)(MatCoarsen);
467:   PetscErrorCode (*view)(MatCoarsen,PetscViewer);
468: };

470: struct _p_MatCoarsen {
471:   PETSCHEADER(struct _MatCoarsenOps);
472:   Mat              graph;
473:   PetscInt         setupcalled;
474:   void             *subctx;
475:   /* */
476:   PetscBool        strict_aggs;
477:   IS               perm;
478:   PetscCoarsenData *agg_lists;
479: };

481: /*
482:     MatFDColoring is used to compute Jacobian matrices efficiently
483:   via coloring. The data structure is explained below in an example.

485:    Color =   0    1     0    2   |   2      3       0
486:    ---------------------------------------------------
487:             00   01              |          05
488:             10   11              |   14     15               Processor  0
489:                        22    23  |          25
490:                        32    33  |
491:    ===================================================
492:                                  |   44     45     46
493:             50                   |          55               Processor 1
494:                                  |   64            66
495:    ---------------------------------------------------

497:     ncolors = 4;

499:     ncolumns      = {2,1,1,0}
500:     columns       = {{0,2},{1},{3},{}}
501:     nrows         = {4,2,3,3}
502:     rows          = {{0,1,2,3},{0,1},{1,2,3},{0,1,2}}
503:     vwscale       = {dx(0),dx(1),dx(2),dx(3)}               MPI Vec
504:     vscale        = {dx(0),dx(1),dx(2),dx(3),dx(4),dx(5)}   Seq Vec

506:     ncolumns      = {1,0,1,1}
507:     columns       = {{6},{},{4},{5}}
508:     nrows         = {3,0,2,2}
509:     rows          = {{0,1,2},{},{1,2},{1,2}}
510:     vwscale       = {dx(4),dx(5),dx(6)}              MPI Vec
511:     vscale        = {dx(0),dx(4),dx(5),dx(6)}        Seq Vec

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

516: */
517: typedef struct {
518:   PetscInt     row;
519:   PetscInt     col;
520:   PetscScalar  *valaddr;   /* address of value */
521: } MatEntry;

523: typedef struct {
524:   PetscInt     row;
525:   PetscScalar  *valaddr;   /* address of value */
526: } MatEntry2;

528: struct  _p_MatFDColoring{
529:   PETSCHEADER(int);
530:   PetscInt       M,N,m;            /* total rows, columns; local rows */
531:   PetscInt       rstart;           /* first row owned by local processor */
532:   PetscInt       ncolors;          /* number of colors */
533:   PetscInt       *ncolumns;        /* number of local columns for a color */
534:   PetscInt       **columns;        /* lists the local columns of each color (using global column numbering) */
535:   PetscInt       *nrows;           /* number of local rows for each color */
536:   MatEntry       *matentry;        /* holds (row, column, address of value) for Jacobian matrix entry */
537:   MatEntry2      *matentry2;       /* holds (row, address of value) for Jacobian matrix entry */
538:   PetscScalar    *dy;              /* store a block of F(x+dx)-F(x) when J is in BAIJ format */
539:   PetscReal      error_rel;        /* square root of relative error in computing function */
540:   PetscReal      umin;             /* minimum allowable u'dx value */
541:   Vec            w1,w2,w3;         /* work vectors used in computing Jacobian */
542:   PetscBool      fset;             /* indicates that the initial function value F(X) is set */
543:   PetscErrorCode (*f)(void);       /* function that defines Jacobian */
544:   void           *fctx;            /* optional user-defined context for use by the function f */
545:   Vec            vscale;           /* holds FD scaling, i.e. 1/dx for each perturbed column */
546:   PetscInt       currentcolor;     /* color for which function evaluation is being done now */
547:   const char     *htype;           /* "wp" or "ds" */
548:   ISColoringType ctype;            /* IS_COLORING_GLOBAL or IS_COLORING_LOCAL */
549:   PetscInt       brows,bcols;      /* number of block rows or columns for speedup inserting the dense matrix into sparse Jacobian */
550:   PetscBool      setupcalled;      /* true if setup has been called */
551:   PetscBool      viewed;           /* true if the -mat_fd_coloring_view has been triggered already */
552:   void           (*ftn_func_pointer)(void),*ftn_func_cntx; /* serve the same purpose as *fortran_func_pointers in PETSc objects */
553: };

555: typedef struct _MatColoringOps *MatColoringOps;
556: struct _MatColoringOps {
557:   PetscErrorCode (*destroy)(MatColoring);
558:   PetscErrorCode (*setfromoptions)(PetscOptionItems*,MatColoring);
559:   PetscErrorCode (*view)(MatColoring,PetscViewer);
560:   PetscErrorCode (*apply)(MatColoring,ISColoring*);
561:   PetscErrorCode (*weights)(MatColoring,PetscReal**,PetscInt**);
562: };

564: struct _p_MatColoring {
565:   PETSCHEADER(struct _MatColoringOps);
566:   Mat                   mat;
567:   PetscInt              dist;             /* distance of the coloring */
568:   PetscInt              maxcolors;        /* the maximum number of colors returned, maxcolors=1 for MIS */
569:   void                  *data;            /* inner context */
570:   PetscBool             valid;            /* check to see if what is produced is a valid coloring */
571:   MatColoringWeightType weight_type;      /* type of weight computation to be performed */
572:   PetscReal             *user_weights;    /* custom weights and permutation */
573:   PetscInt              *user_lperm;
574:   PetscBool             valid_iscoloring; /* check to see if matcoloring is produced a valid iscoloring */
575: };

577: struct  _p_MatTransposeColoring{
578:   PETSCHEADER(int);
579:   PetscInt       M,N,m;            /* total rows, columns; local rows */
580:   PetscInt       rstart;           /* first row owned by local processor */
581:   PetscInt       ncolors;          /* number of colors */
582:   PetscInt       *ncolumns;        /* number of local columns for a color */
583:   PetscInt       *nrows;           /* number of local rows for each color */
584:   PetscInt       currentcolor;     /* color for which function evaluation is being done now */
585:   ISColoringType ctype;            /* IS_COLORING_GLOBAL or IS_COLORING_LOCAL */

587:   PetscInt       *colorforrow,*colorforcol;  /* pointer to rows and columns */
588:   PetscInt       *rows;                      /* lists the local rows for each color (using the local row numbering) */
589:   PetscInt       *den2sp;                    /* maps (row,color) in the dense matrix to index of sparse matrix array a->a */
590:   PetscInt       *columns;                   /* lists the local columns of each color (using global column numbering) */
591:   PetscInt       brows;                      /* number of rows for efficient implementation of MatTransColoringApplyDenToSp() */
592:   PetscInt       *lstart;                    /* array used for loop over row blocks of Csparse */
593: };

595: /*
596:    Null space context for preconditioner/operators
597: */
598: struct _p_MatNullSpace {
599:   PETSCHEADER(int);
600:   PetscBool      has_cnst;
601:   PetscInt       n;
602:   Vec*           vecs;
603:   PetscScalar*   alpha;                 /* for projections */
604:   PetscErrorCode (*remove)(MatNullSpace,Vec,void*);  /* for user provided removal function */
605:   void*          rmctx;                 /* context for remove() function */
606: };

608: /*
609:    Checking zero pivot for LU, ILU preconditioners.
610: */
611: typedef struct {
612:   PetscInt       nshift,nshift_max;
613:   PetscReal      shift_amount,shift_lo,shift_hi,shift_top,shift_fraction;
614:   PetscBool      newshift;
615:   PetscReal      rs;  /* active row sum of abs(offdiagonals) */
616:   PetscScalar    pv;  /* pivot of the active row */
617: } FactorShiftCtx;

619: /*
620:  Used by MatCreateSubMatrices_MPIXAIJ_Local()
621: */
622:  #include <petscctable.h>
623: typedef struct { /* used by MatCreateSubMatrices_MPIAIJ_SingleIS_Local() and MatCreateSubMatrices_MPIAIJ_Local */
624:   PetscInt   id;   /* index of submats, only submats[0] is responsible for deleting some arrays below */
625:   PetscInt   nrqs,nrqr;
626:   PetscInt   **rbuf1,**rbuf2,**rbuf3,**sbuf1,**sbuf2;
627:   PetscInt   **ptr;
628:   PetscInt   *tmp;
629:   PetscInt   *ctr;
630:   PetscInt   *pa; /* proc array */
631:   PetscInt   *req_size,*req_source1,*req_source2;
632:   PetscBool  allcolumns,allrows;
633:   PetscBool  singleis;
634:   PetscInt   *row2proc; /* row to proc map */
635:   PetscInt   nstages;
636: #if defined(PETSC_USE_CTABLE)
637:   PetscTable cmap,rmap;
638:   PetscInt   *cmap_loc,*rmap_loc;
639: #else
640:   PetscInt   *cmap,*rmap;
641: #endif

643:   PetscErrorCode (*destroy)(Mat);
644: } Mat_SubSppt;

646: PETSC_EXTERN PetscErrorCode MatFactorDumpMatrix(Mat);
647: PETSC_INTERN PetscErrorCode MatShift_Basic(Mat,PetscScalar);
648: PETSC_INTERN PetscErrorCode MatSetBlockSizes_Default(Mat,PetscInt,PetscInt);

650: PETSC_STATIC_INLINE PetscErrorCode MatPivotCheck_nz(Mat mat,const MatFactorInfo *info,FactorShiftCtx *sctx,PetscInt row)
651: {
652:   PetscReal _rs   = sctx->rs;
653:   PetscReal _zero = info->zeropivot*_rs;

656:   if (PetscAbsScalar(sctx->pv) <= _zero && !PetscIsNanScalar(sctx->pv)){
657:     /* force |diag| > zeropivot*rs */
658:     if (!sctx->nshift) sctx->shift_amount = info->shiftamount;
659:     else sctx->shift_amount *= 2.0;
660:     sctx->newshift = PETSC_TRUE;
661:     (sctx->nshift)++;
662:   } else {
663:     sctx->newshift = PETSC_FALSE;
664:   }
665:   return(0);
666: }

668: PETSC_STATIC_INLINE PetscErrorCode MatPivotCheck_pd(Mat mat,const MatFactorInfo *info,FactorShiftCtx *sctx,PetscInt row)
669: {
670:   PetscReal _rs   = sctx->rs;
671:   PetscReal _zero = info->zeropivot*_rs;

674:   if (PetscRealPart(sctx->pv) <= _zero && !PetscIsNanScalar(sctx->pv)){
675:     /* force matfactor to be diagonally dominant */
676:     if (sctx->nshift == sctx->nshift_max) {
677:       sctx->shift_fraction = sctx->shift_hi;
678:     } else {
679:       sctx->shift_lo = sctx->shift_fraction;
680:       sctx->shift_fraction = (sctx->shift_hi+sctx->shift_lo)/2.;
681:     }
682:     sctx->shift_amount = sctx->shift_fraction * sctx->shift_top;
683:     sctx->nshift++;
684:     sctx->newshift = PETSC_TRUE;
685:   } else {
686:     sctx->newshift = PETSC_FALSE;
687:   }
688:   return(0);
689: }

691: PETSC_STATIC_INLINE PetscErrorCode MatPivotCheck_inblocks(Mat mat,const MatFactorInfo *info,FactorShiftCtx *sctx,PetscInt row)
692: {
693:   PetscReal _zero = info->zeropivot;

696:   if (PetscAbsScalar(sctx->pv) <= _zero && !PetscIsNanScalar(sctx->pv)){
697:     sctx->pv          += info->shiftamount;
698:     sctx->shift_amount = 0.0;
699:     sctx->nshift++;
700:   }
701:   sctx->newshift = PETSC_FALSE;
702:   return(0);
703: }

705: PETSC_STATIC_INLINE PetscErrorCode MatPivotCheck_none(Mat fact,Mat mat,const MatFactorInfo *info,FactorShiftCtx *sctx,PetscInt row)
706: {
707:   PetscReal      _zero = info->zeropivot;

711:   sctx->newshift = PETSC_FALSE;
712:   if (PetscAbsScalar(sctx->pv) <= _zero && !PetscIsNanScalar(sctx->pv)) {
713:     if (!mat->erroriffailure) {
714:       PetscInfo3(mat,"Detected zero pivot in factorization in row %D value %g tolerance %g\n",row,(double)PetscAbsScalar(sctx->pv),(double)_zero);
715:       fact->factorerrortype             = MAT_FACTOR_NUMERIC_ZEROPIVOT;
716:       fact->factorerror_zeropivot_value = PetscAbsScalar(sctx->pv);
717:       fact->factorerror_zeropivot_row   = row;
718:     } else SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot row %D value %g tolerance %g\n",row,(double)PetscAbsScalar(sctx->pv),(double)_zero);
719:   }
720:   return(0);
721: }

723: PETSC_STATIC_INLINE PetscErrorCode MatPivotCheck(Mat fact,Mat mat,const MatFactorInfo *info,FactorShiftCtx *sctx,PetscInt row)
724: {

728:   if (info->shifttype == (PetscReal) MAT_SHIFT_NONZERO){
729:     MatPivotCheck_nz(mat,info,sctx,row);
730:   } else if (info->shifttype == (PetscReal) MAT_SHIFT_POSITIVE_DEFINITE){
731:     MatPivotCheck_pd(mat,info,sctx,row);
732:   } else if (info->shifttype == (PetscReal) MAT_SHIFT_INBLOCKS){
733:     MatPivotCheck_inblocks(mat,info,sctx,row);
734:   } else {
735:     MatPivotCheck_none(fact,mat,info,sctx,row);
736:   }
737:   return(0);
738: }

740: /*
741:   Create and initialize a linked list
742:   Input Parameters:
743:     idx_start - starting index of the list
744:     lnk_max   - max value of lnk indicating the end of the list
745:     nlnk      - max length of the list
746:   Output Parameters:
747:     lnk       - list initialized
748:     bt        - PetscBT (bitarray) with all bits set to false
749:     lnk_empty - flg indicating the list is empty
750: */
751: #define PetscLLCreate(idx_start,lnk_max,nlnk,lnk,bt) \
752:   (PetscMalloc1(nlnk,&lnk) || PetscBTCreate(nlnk,&(bt)) || (lnk[idx_start] = lnk_max,0))

754: #define PetscLLCreate_new(idx_start,lnk_max,nlnk,lnk,bt,lnk_empty)\
755:   (PetscMalloc1(nlnk,&lnk) || PetscBTCreate(nlnk,&(bt)) || (lnk_empty = PETSC_TRUE,0) ||(lnk[idx_start] = lnk_max,0))

757: /*
758:   Add an index set into a sorted linked list
759:   Input Parameters:
760:     nidx      - number of input indices
761:     indices   - integer array
762:     idx_start - starting index of the list
763:     lnk       - linked list(an integer array) that is created
764:     bt        - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
765:   output Parameters:
766:     nlnk      - number of newly added indices
767:     lnk       - the sorted(increasing order) linked list containing new and non-redundate entries from indices
768:     bt        - updated PetscBT (bitarray)
769: */
770: #define PetscLLAdd(nidx,indices,idx_start,nlnk,lnk,bt) 0;\
771: {\
772:   PetscInt _k,_entry,_location,_lnkdata;\
773:   nlnk     = 0;\
774:   _lnkdata = idx_start;\
775:   for (_k=0; _k<nidx; _k++){\
776:     _entry = indices[_k];\
777:     if (!PetscBTLookupSet(bt,_entry)){  /* new entry */\
778:       /* search for insertion location */\
779:       /* start from the beginning if _entry < previous _entry */\
780:       if (_k && _entry < _lnkdata) _lnkdata  = idx_start;\
781:       do {\
782:         _location = _lnkdata;\
783:         _lnkdata  = lnk[_location];\
784:       } while (_entry > _lnkdata);\
785:       /* insertion location is found, add entry into lnk */\
786:       lnk[_location] = _entry;\
787:       lnk[_entry]    = _lnkdata;\
788:       nlnk++;\
789:       _lnkdata = _entry; /* next search starts from here if next_entry > _entry */\
790:     }\
791:   }\
792: }

794: /*
795:   Add a permuted index set into a sorted linked list
796:   Input Parameters:
797:     nidx      - number of input indices
798:     indices   - integer array
799:     perm      - permutation of indices
800:     idx_start - starting index of the list
801:     lnk       - linked list(an integer array) that is created
802:     bt        - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
803:   output Parameters:
804:     nlnk      - number of newly added indices
805:     lnk       - the sorted(increasing order) linked list containing new and non-redundate entries from indices
806:     bt        - updated PetscBT (bitarray)
807: */
808: #define PetscLLAddPerm(nidx,indices,perm,idx_start,nlnk,lnk,bt) 0;\
809: {\
810:   PetscInt _k,_entry,_location,_lnkdata;\
811:   nlnk     = 0;\
812:   _lnkdata = idx_start;\
813:   for (_k=0; _k<nidx; _k++){\
814:     _entry = perm[indices[_k]];\
815:     if (!PetscBTLookupSet(bt,_entry)){  /* new entry */\
816:       /* search for insertion location */\
817:       /* start from the beginning if _entry < previous _entry */\
818:       if (_k && _entry < _lnkdata) _lnkdata  = idx_start;\
819:       do {\
820:         _location = _lnkdata;\
821:         _lnkdata  = lnk[_location];\
822:       } while (_entry > _lnkdata);\
823:       /* insertion location is found, add entry into lnk */\
824:       lnk[_location] = _entry;\
825:       lnk[_entry]    = _lnkdata;\
826:       nlnk++;\
827:       _lnkdata = _entry; /* next search starts from here if next_entry > _entry */\
828:     }\
829:   }\
830: }

832: /*
833:   Add a SORTED ascending index set into a sorted linked list - same as PetscLLAdd() bus skip 'if (_k && _entry < _lnkdata) _lnkdata  = idx_start;'
834:   Input Parameters:
835:     nidx      - number of input indices
836:     indices   - sorted integer array
837:     idx_start - starting index of the list
838:     lnk       - linked list(an integer array) that is created
839:     bt        - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
840:   output Parameters:
841:     nlnk      - number of newly added indices
842:     lnk       - the sorted(increasing order) linked list containing new and non-redundate entries from indices
843:     bt        - updated PetscBT (bitarray)
844: */
845: #define PetscLLAddSorted(nidx,indices,idx_start,nlnk,lnk,bt) 0;\
846: {\
847:   PetscInt _k,_entry,_location,_lnkdata;\
848:   nlnk      = 0;\
849:   _lnkdata  = idx_start;\
850:   for (_k=0; _k<nidx; _k++){\
851:     _entry = indices[_k];\
852:     if (!PetscBTLookupSet(bt,_entry)){  /* new entry */\
853:       /* search for insertion location */\
854:       do {\
855:         _location = _lnkdata;\
856:         _lnkdata  = lnk[_location];\
857:       } while (_entry > _lnkdata);\
858:       /* insertion location is found, add entry into lnk */\
859:       lnk[_location] = _entry;\
860:       lnk[_entry]    = _lnkdata;\
861:       nlnk++;\
862:       _lnkdata = _entry; /* next search starts from here */\
863:     }\
864:   }\
865: }

867: #define PetscLLAddSorted_new(nidx,indices,idx_start,lnk_empty,nlnk,lnk,bt) 0; \
868: {\
869:   PetscInt _k,_entry,_location,_lnkdata;\
870:   if (lnk_empty){\
871:     _lnkdata  = idx_start;                      \
872:     for (_k=0; _k<nidx; _k++){                  \
873:       _entry = indices[_k];                             \
874:       PetscBTSet(bt,_entry);  /* mark the new entry */          \
875:           _location = _lnkdata;                                 \
876:           _lnkdata  = lnk[_location];                           \
877:         /* insertion location is found, add entry into lnk */   \
878:         lnk[_location] = _entry;                                \
879:         lnk[_entry]    = _lnkdata;                              \
880:         _lnkdata = _entry; /* next search starts from here */   \
881:     }                                                           \
882:     /*\
883:     lnk[indices[nidx-1]] = lnk[idx_start];\
884:     lnk[idx_start]       = indices[0];\
885:     PetscBTSet(bt,indices[0]);  \
886:     for (_k=1; _k<nidx; _k++){                  \
887:       PetscBTSet(bt,indices[_k]);                                          \
888:       lnk[indices[_k-1]] = indices[_k];                                  \
889:     }                                                           \
890:      */\
891:     nlnk      = nidx;\
892:     lnk_empty = PETSC_FALSE;\
893:   } else {\
894:     nlnk      = 0;                              \
895:     _lnkdata  = idx_start;                      \
896:     for (_k=0; _k<nidx; _k++){                  \
897:       _entry = indices[_k];                             \
898:       if (!PetscBTLookupSet(bt,_entry)){  /* new entry */       \
899:         /* search for insertion location */                     \
900:         do {                                                    \
901:           _location = _lnkdata;                                 \
902:           _lnkdata  = lnk[_location];                           \
903:         } while (_entry > _lnkdata);                            \
904:         /* insertion location is found, add entry into lnk */   \
905:         lnk[_location] = _entry;                                \
906:         lnk[_entry]    = _lnkdata;                              \
907:         nlnk++;                                                 \
908:         _lnkdata = _entry; /* next search starts from here */   \
909:       }                                                         \
910:     }                                                           \
911:   }                                                             \
912: }

914: /*
915:   Add a SORTED index set into a sorted linked list used for LUFactorSymbolic()
916:   Same as PetscLLAddSorted() with an additional operation:
917:        count the number of input indices that are no larger than 'diag'
918:   Input Parameters:
919:     indices   - sorted integer array
920:     idx_start - starting index of the list, index of pivot row
921:     lnk       - linked list(an integer array) that is created
922:     bt        - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
923:     diag      - index of the active row in LUFactorSymbolic
924:     nzbd      - number of input indices with indices <= idx_start
925:     im        - im[idx_start] is initialized as num of nonzero entries in row=idx_start
926:   output Parameters:
927:     nlnk      - number of newly added indices
928:     lnk       - the sorted(increasing order) linked list containing new and non-redundate entries from indices
929:     bt        - updated PetscBT (bitarray)
930:     im        - im[idx_start]: unchanged if diag is not an entry
931:                              : num of entries with indices <= diag if diag is an entry
932: */
933: #define PetscLLAddSortedLU(indices,idx_start,nlnk,lnk,bt,diag,nzbd,im) 0;\
934: {\
935:   PetscInt _k,_entry,_location,_lnkdata,_nidx;\
936:   nlnk     = 0;\
937:   _lnkdata = idx_start;\
938:   _nidx = im[idx_start] - nzbd; /* num of entries with idx_start < index <= diag */\
939:   for (_k=0; _k<_nidx; _k++){\
940:     _entry = indices[_k];\
941:     nzbd++;\
942:     if ( _entry== diag) im[idx_start] = nzbd;\
943:     if (!PetscBTLookupSet(bt,_entry)){  /* new entry */\
944:       /* search for insertion location */\
945:       do {\
946:         _location = _lnkdata;\
947:         _lnkdata  = lnk[_location];\
948:       } while (_entry > _lnkdata);\
949:       /* insertion location is found, add entry into lnk */\
950:       lnk[_location] = _entry;\
951:       lnk[_entry]    = _lnkdata;\
952:       nlnk++;\
953:       _lnkdata = _entry; /* next search starts from here */\
954:     }\
955:   }\
956: }

958: /*
959:   Copy data on the list into an array, then initialize the list
960:   Input Parameters:
961:     idx_start - starting index of the list
962:     lnk_max   - max value of lnk indicating the end of the list
963:     nlnk      - number of data on the list to be copied
964:     lnk       - linked list
965:     bt        - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
966:   output Parameters:
967:     indices   - array that contains the copied data
968:     lnk       - linked list that is cleaned and initialize
969:     bt        - PetscBT (bitarray) with all bits set to false
970: */
971: #define PetscLLClean(idx_start,lnk_max,nlnk,lnk,indices,bt) 0;\
972: {\
973:   PetscInt _j,_idx=idx_start;\
974:   for (_j=0; _j<nlnk; _j++){\
975:     _idx = lnk[_idx];\
976:     indices[_j] = _idx;\
977:     PetscBTClear(bt,_idx);\
978:   }\
979:   lnk[idx_start] = lnk_max;\
980: }
981: /*
982:   Free memories used by the list
983: */
984: #define PetscLLDestroy(lnk,bt) (PetscFree(lnk) || PetscBTDestroy(&(bt)))

986: /* Routines below are used for incomplete matrix factorization */
987: /*
988:   Create and initialize a linked list and its levels
989:   Input Parameters:
990:     idx_start - starting index of the list
991:     lnk_max   - max value of lnk indicating the end of the list
992:     nlnk      - max length of the list
993:   Output Parameters:
994:     lnk       - list initialized
995:     lnk_lvl   - array of size nlnk for storing levels of lnk
996:     bt        - PetscBT (bitarray) with all bits set to false
997: */
998: #define PetscIncompleteLLCreate(idx_start,lnk_max,nlnk,lnk,lnk_lvl,bt)\
999:   (PetscIntMultError(2,nlnk,NULL) || PetscMalloc1(2*nlnk,&lnk) || PetscBTCreate(nlnk,&(bt)) || (lnk[idx_start] = lnk_max,lnk_lvl = lnk + nlnk,0))

1001: /*
1002:   Initialize a sorted linked list used for ILU and ICC
1003:   Input Parameters:
1004:     nidx      - number of input idx
1005:     idx       - integer array used for storing column indices
1006:     idx_start - starting index of the list
1007:     perm      - indices of an IS
1008:     lnk       - linked list(an integer array) that is created
1009:     lnklvl    - levels of lnk
1010:     bt        - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
1011:   output Parameters:
1012:     nlnk     - number of newly added idx
1013:     lnk      - the sorted(increasing order) linked list containing new and non-redundate entries from idx
1014:     lnklvl   - levels of lnk
1015:     bt       - updated PetscBT (bitarray)
1016: */
1017: #define PetscIncompleteLLInit(nidx,idx,idx_start,perm,nlnk,lnk,lnklvl,bt) 0;\
1018: {\
1019:   PetscInt _k,_entry,_location,_lnkdata;\
1020:   nlnk     = 0;\
1021:   _lnkdata = idx_start;\
1022:   for (_k=0; _k<nidx; _k++){\
1023:     _entry = perm[idx[_k]];\
1024:     if (!PetscBTLookupSet(bt,_entry)){  /* new entry */\
1025:       /* search for insertion location */\
1026:       if (_k && _entry < _lnkdata) _lnkdata  = idx_start;\
1027:       do {\
1028:         _location = _lnkdata;\
1029:         _lnkdata  = lnk[_location];\
1030:       } while (_entry > _lnkdata);\
1031:       /* insertion location is found, add entry into lnk */\
1032:       lnk[_location]  = _entry;\
1033:       lnk[_entry]     = _lnkdata;\
1034:       lnklvl[_entry] = 0;\
1035:       nlnk++;\
1036:       _lnkdata = _entry; /* next search starts from here if next_entry > _entry */\
1037:     }\
1038:   }\
1039: }

1041: /*
1042:   Add a SORTED index set into a sorted linked list for ILU
1043:   Input Parameters:
1044:     nidx      - number of input indices
1045:     idx       - sorted integer array used for storing column indices
1046:     level     - level of fill, e.g., ICC(level)
1047:     idxlvl    - level of idx
1048:     idx_start - starting index of the list
1049:     lnk       - linked list(an integer array) that is created
1050:     lnklvl    - levels of lnk
1051:     bt        - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
1052:     prow      - the row number of idx
1053:   output Parameters:
1054:     nlnk     - number of newly added idx
1055:     lnk      - the sorted(increasing order) linked list containing new and non-redundate entries from idx
1056:     lnklvl   - levels of lnk
1057:     bt       - updated PetscBT (bitarray)

1059:   Note: the level of factor(i,j) is set as lvl(i,j) = min{ lvl(i,j), lvl(i,prow)+lvl(prow,j)+1)
1060:         where idx = non-zero columns of U(prow,prow+1:n-1), prow<i
1061: */
1062: #define PetscILULLAddSorted(nidx,idx,level,idxlvl,idx_start,nlnk,lnk,lnklvl,bt,lnklvl_prow) 0;\
1063: {\
1064:   PetscInt _k,_entry,_location,_lnkdata,_incrlev,_lnklvl_prow=lnklvl[prow];\
1065:   nlnk     = 0;\
1066:   _lnkdata = idx_start;\
1067:   for (_k=0; _k<nidx; _k++){\
1068:     _incrlev = idxlvl[_k] + _lnklvl_prow + 1;\
1069:     if (_incrlev > level) continue;\
1070:     _entry = idx[_k];\
1071:     if (!PetscBTLookupSet(bt,_entry)){  /* new entry */\
1072:       /* search for insertion location */\
1073:       do {\
1074:         _location = _lnkdata;\
1075:         _lnkdata  = lnk[_location];\
1076:       } while (_entry > _lnkdata);\
1077:       /* insertion location is found, add entry into lnk */\
1078:       lnk[_location]  = _entry;\
1079:       lnk[_entry]     = _lnkdata;\
1080:       lnklvl[_entry] = _incrlev;\
1081:       nlnk++;\
1082:       _lnkdata = _entry; /* next search starts from here if next_entry > _entry */\
1083:     } else { /* existing entry: update lnklvl */\
1084:       if (lnklvl[_entry] > _incrlev) lnklvl[_entry] = _incrlev;\
1085:     }\
1086:   }\
1087: }

1089: /*
1090:   Add a index set into a sorted linked list
1091:   Input Parameters:
1092:     nidx      - number of input idx
1093:     idx   - integer array used for storing column indices
1094:     level     - level of fill, e.g., ICC(level)
1095:     idxlvl - level of idx
1096:     idx_start - starting index of the list
1097:     lnk       - linked list(an integer array) that is created
1098:     lnklvl   - levels of lnk
1099:     bt        - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
1100:   output Parameters:
1101:     nlnk      - number of newly added idx
1102:     lnk       - the sorted(increasing order) linked list containing new and non-redundate entries from idx
1103:     lnklvl   - levels of lnk
1104:     bt        - updated PetscBT (bitarray)
1105: */
1106: #define PetscIncompleteLLAdd(nidx,idx,level,idxlvl,idx_start,nlnk,lnk,lnklvl,bt) 0;\
1107: {\
1108:   PetscInt _k,_entry,_location,_lnkdata,_incrlev;\
1109:   nlnk     = 0;\
1110:   _lnkdata = idx_start;\
1111:   for (_k=0; _k<nidx; _k++){\
1112:     _incrlev = idxlvl[_k] + 1;\
1113:     if (_incrlev > level) continue;\
1114:     _entry = idx[_k];\
1115:     if (!PetscBTLookupSet(bt,_entry)){  /* new entry */\
1116:       /* search for insertion location */\
1117:       if (_k && _entry < _lnkdata) _lnkdata  = idx_start;\
1118:       do {\
1119:         _location = _lnkdata;\
1120:         _lnkdata  = lnk[_location];\
1121:       } while (_entry > _lnkdata);\
1122:       /* insertion location is found, add entry into lnk */\
1123:       lnk[_location]  = _entry;\
1124:       lnk[_entry]     = _lnkdata;\
1125:       lnklvl[_entry] = _incrlev;\
1126:       nlnk++;\
1127:       _lnkdata = _entry; /* next search starts from here if next_entry > _entry */\
1128:     } else { /* existing entry: update lnklvl */\
1129:       if (lnklvl[_entry] > _incrlev) lnklvl[_entry] = _incrlev;\
1130:     }\
1131:   }\
1132: }

1134: /*
1135:   Add a SORTED index set into a sorted linked list
1136:   Input Parameters:
1137:     nidx      - number of input indices
1138:     idx   - sorted integer array used for storing column indices
1139:     level     - level of fill, e.g., ICC(level)
1140:     idxlvl - level of idx
1141:     idx_start - starting index of the list
1142:     lnk       - linked list(an integer array) that is created
1143:     lnklvl    - levels of lnk
1144:     bt        - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
1145:   output Parameters:
1146:     nlnk      - number of newly added idx
1147:     lnk       - the sorted(increasing order) linked list containing new and non-redundate entries from idx
1148:     lnklvl    - levels of lnk
1149:     bt        - updated PetscBT (bitarray)
1150: */
1151: #define PetscIncompleteLLAddSorted(nidx,idx,level,idxlvl,idx_start,nlnk,lnk,lnklvl,bt) 0;\
1152: {\
1153:   PetscInt _k,_entry,_location,_lnkdata,_incrlev;\
1154:   nlnk = 0;\
1155:   _lnkdata = idx_start;\
1156:   for (_k=0; _k<nidx; _k++){\
1157:     _incrlev = idxlvl[_k] + 1;\
1158:     if (_incrlev > level) continue;\
1159:     _entry = idx[_k];\
1160:     if (!PetscBTLookupSet(bt,_entry)){  /* new entry */\
1161:       /* search for insertion location */\
1162:       do {\
1163:         _location = _lnkdata;\
1164:         _lnkdata  = lnk[_location];\
1165:       } while (_entry > _lnkdata);\
1166:       /* insertion location is found, add entry into lnk */\
1167:       lnk[_location] = _entry;\
1168:       lnk[_entry]    = _lnkdata;\
1169:       lnklvl[_entry] = _incrlev;\
1170:       nlnk++;\
1171:       _lnkdata = _entry; /* next search starts from here */\
1172:     } else { /* existing entry: update lnklvl */\
1173:       if (lnklvl[_entry] > _incrlev) lnklvl[_entry] = _incrlev;\
1174:     }\
1175:   }\
1176: }

1178: /*
1179:   Add a SORTED index set into a sorted linked list for ICC
1180:   Input Parameters:
1181:     nidx      - number of input indices
1182:     idx       - sorted integer array used for storing column indices
1183:     level     - level of fill, e.g., ICC(level)
1184:     idxlvl    - level of idx
1185:     idx_start - starting index of the list
1186:     lnk       - linked list(an integer array) that is created
1187:     lnklvl    - levels of lnk
1188:     bt        - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
1189:     idxlvl_prow - idxlvl[prow], where prow is the row number of the idx
1190:   output Parameters:
1191:     nlnk   - number of newly added indices
1192:     lnk    - the sorted(increasing order) linked list containing new and non-redundate entries from idx
1193:     lnklvl - levels of lnk
1194:     bt     - updated PetscBT (bitarray)
1195:   Note: the level of U(i,j) is set as lvl(i,j) = min{ lvl(i,j), lvl(prow,i)+lvl(prow,j)+1)
1196:         where idx = non-zero columns of U(prow,prow+1:n-1), prow<i
1197: */
1198: #define PetscICCLLAddSorted(nidx,idx,level,idxlvl,idx_start,nlnk,lnk,lnklvl,bt,idxlvl_prow) 0;\
1199: {\
1200:   PetscInt _k,_entry,_location,_lnkdata,_incrlev;\
1201:   nlnk = 0;\
1202:   _lnkdata = idx_start;\
1203:   for (_k=0; _k<nidx; _k++){\
1204:     _incrlev = idxlvl[_k] + idxlvl_prow + 1;\
1205:     if (_incrlev > level) continue;\
1206:     _entry = idx[_k];\
1207:     if (!PetscBTLookupSet(bt,_entry)){  /* new entry */\
1208:       /* search for insertion location */\
1209:       do {\
1210:         _location = _lnkdata;\
1211:         _lnkdata  = lnk[_location];\
1212:       } while (_entry > _lnkdata);\
1213:       /* insertion location is found, add entry into lnk */\
1214:       lnk[_location] = _entry;\
1215:       lnk[_entry]    = _lnkdata;\
1216:       lnklvl[_entry] = _incrlev;\
1217:       nlnk++;\
1218:       _lnkdata = _entry; /* next search starts from here */\
1219:     } else { /* existing entry: update lnklvl */\
1220:       if (lnklvl[_entry] > _incrlev) lnklvl[_entry] = _incrlev;\
1221:     }\
1222:   }\
1223: }

1225: /*
1226:   Copy data on the list into an array, then initialize the list
1227:   Input Parameters:
1228:     idx_start - starting index of the list
1229:     lnk_max   - max value of lnk indicating the end of the list
1230:     nlnk      - number of data on the list to be copied
1231:     lnk       - linked list
1232:     lnklvl    - level of lnk
1233:     bt        - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
1234:   output Parameters:
1235:     indices - array that contains the copied data
1236:     lnk     - linked list that is cleaned and initialize
1237:     lnklvl  - level of lnk that is reinitialized
1238:     bt      - PetscBT (bitarray) with all bits set to false
1239: */
1240: #define PetscIncompleteLLClean(idx_start,lnk_max,nlnk,lnk,lnklvl,indices,indiceslvl,bt) 0;\
1241: {\
1242:   PetscInt _j,_idx=idx_start;\
1243:   for (_j=0; _j<nlnk; _j++){\
1244:     _idx = lnk[_idx];\
1245:     *(indices+_j) = _idx;\
1246:     *(indiceslvl+_j) = lnklvl[_idx];\
1247:     lnklvl[_idx] = -1;\
1248:     PetscBTClear(bt,_idx);\
1249:   }\
1250:   lnk[idx_start] = lnk_max;\
1251: }
1252: /*
1253:   Free memories used by the list
1254: */
1255: #define PetscIncompleteLLDestroy(lnk,bt) (PetscFree(lnk) || PetscBTDestroy(&(bt)))

1257: #define MatCheckSameLocalSize(A,ar1,B,ar2) \
1259:   if ((A->rmap->n != B->rmap->n) || (A->cmap->n != B->cmap->n)) SETERRQ6(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Incompatible matrix local sizes: parameter # %d (%D x %D) != parameter # %d (%D x %D)",ar1,A->rmap->n,A->cmap->n,ar2,B->rmap->n,B->cmap->n);
1260: 
1261: #define MatCheckSameSize(A,ar1,B,ar2) \
1262:   if ((A->rmap->N != B->rmap->N) || (A->cmap->N != B->cmap->N)) SETERRQ6(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_INCOMP,"Incompatible matrix global sizes: parameter # %d (%D x %D) != parameter # %d (%D x %D)",ar1,A->rmap->N,A->cmap->N,ar2,B->rmap->N,B->cmap->N);\
1263:   MatCheckSameLocalSize(A,ar1,B,ar2);
1264: 
1265: #define VecCheckMatCompatible(M,x,ar1,b,ar2)                               \
1266:   if (M->cmap->N != x->map->N) SETERRQ3(PetscObjectComm((PetscObject)M),PETSC_ERR_ARG_SIZ,"Vector global length incompatible with matrix: parameter # %d global size %D != matrix column global size %D",ar1,x->map->N,M->cmap->N);\
1267:   if (M->rmap->N != b->map->N) SETERRQ3(PetscObjectComm((PetscObject)M),PETSC_ERR_ARG_SIZ,"Vector global length incompatible with matrix: parameter # %d global size %D != matrix row global size %D",ar2,b->map->N,M->rmap->N);

1269: /* -------------------------------------------------------------------------------------------------------*/
1270:  #include <petscbt.h>
1271: /*
1272:   Create and initialize a condensed linked list -
1273:     same as PetscLLCreate(), but uses a scalable array 'lnk' with size of max number of entries, not O(N).
1274:     Barry suggested this approach (Dec. 6, 2011):
1275:       I've thought of an alternative way of representing a linked list that is efficient but doesn't have the O(N) scaling issue
1276:       (it may be faster than the O(N) even sequentially due to less crazy memory access).

1278:       Instead of having some like  a  2  -> 4 -> 11 ->  22  list that uses slot 2  4 11 and 22 in a big array use a small array with two slots
1279:       for each entry for example  [ 2 1 | 4 3 | 22 -1 | 11 2]   so the first number (of the pair) is the value while the second tells you where
1280:       in the list the next entry is. Inserting a new link means just append another pair at the end. For example say we want to insert 13 into the
1281:       list it would then become [2 1 | 4 3 | 22 -1 | 11 4 | 13 2 ] you just add a pair at the end and fix the point for the one that points to it.
1282:       That is 11 use to point to the 2 slot, after the change 11 points to the 4th slot which has the value 13. Note that values are always next
1283:       to each other so memory access is much better than using the big array.

1285:   Example:
1286:      nlnk_max=5, lnk_max=36:
1287:      Initial list: [0, 0 | 36, 2 | 0, 0 | 0, 0 | 0, 0 | 0, 0 | 0, 0]
1288:      here, head_node has index 2 with value lnk[2]=lnk_max=36,
1289:            0-th entry is used to store the number of entries in the list,
1290:      The initial lnk represents head -> tail(marked by 36) with number of entries = lnk[0]=0.

1292:      Now adding a sorted set {2,4}, the list becomes
1293:      [2, 0 | 36, 4 |2, 6 | 4, 2 | 0, 0 | 0, 0 | 0, 0 ]
1294:      represents head -> 2 -> 4 -> tail with number of entries = lnk[0]=2.

1296:      Then adding a sorted set {0,3,35}, the list
1297:      [5, 0 | 36, 8 | 2, 10 | 4, 12 | 0, 4 | 3, 6 | 35, 2 ]
1298:      represents head -> 0 -> 2 -> 3 -> 4 -> 35 -> tail with number of entries = lnk[0]=5.

1300:   Input Parameters:
1301:     nlnk_max  - max length of the list
1302:     lnk_max   - max value of the entries
1303:   Output Parameters:
1304:     lnk       - list created and initialized
1305:     bt        - PetscBT (bitarray) with all bits set to false. Note: bt has size lnk_max, not nln_max!
1306: */
1307: PETSC_STATIC_INLINE PetscErrorCode PetscLLCondensedCreate(PetscInt nlnk_max,PetscInt lnk_max,PetscInt **lnk,PetscBT *bt)
1308: {
1310:   PetscInt       *llnk,lsize = 0;

1313:   PetscIntMultError(2,nlnk_max+2,&lsize);
1314:   PetscMalloc1(lsize,lnk);
1315:   PetscBTCreate(lnk_max,bt);
1316:   llnk = *lnk;
1317:   llnk[0] = 0;         /* number of entries on the list */
1318:   llnk[2] = lnk_max;   /* value in the head node */
1319:   llnk[3] = 2;         /* next for the head node */
1320:   return(0);
1321: }

1323: /*
1324:   Add a SORTED ascending index set into a sorted linked list. See PetscLLCondensedCreate() for detailed description.
1325:   Input Parameters:
1326:     nidx      - number of input indices
1327:     indices   - sorted integer array
1328:     lnk       - condensed linked list(an integer array) that is created
1329:     bt        - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
1330:   output Parameters:
1331:     lnk       - the sorted(increasing order) linked list containing previous and newly added non-redundate indices
1332:     bt        - updated PetscBT (bitarray)
1333: */
1334: PETSC_STATIC_INLINE PetscErrorCode PetscLLCondensedAddSorted(PetscInt nidx,const PetscInt indices[],PetscInt lnk[],PetscBT bt)
1335: {
1336:   PetscInt _k,_entry,_location,_next,_lnkdata,_nlnk,_newnode;

1339:   _nlnk     = lnk[0]; /* num of entries on the input lnk */
1340:   _location = 2; /* head */
1341:     for (_k=0; _k<nidx; _k++){
1342:       _entry = indices[_k];
1343:       if (!PetscBTLookupSet(bt,_entry)){  /* new entry */
1344:         /* search for insertion location */
1345:         do {
1346:           _next     = _location + 1; /* link from previous node to next node */
1347:           _location = lnk[_next];    /* idx of next node */
1348:           _lnkdata  = lnk[_location];/* value of next node */
1349:         } while (_entry > _lnkdata);
1350:         /* insertion location is found, add entry into lnk */
1351:         _newnode        = 2*(_nlnk+2);   /* index for this new node */
1352:         lnk[_next]      = _newnode;      /* connect previous node to the new node */
1353:         lnk[_newnode]   = _entry;        /* set value of the new node */
1354:         lnk[_newnode+1] = _location;     /* connect new node to next node */
1355:         _location       = _newnode;      /* next search starts from the new node */
1356:         _nlnk++;
1357:       }   \
1358:     }\
1359:   lnk[0]   = _nlnk;   /* number of entries in the list */
1360:   return(0);
1361: }

1363: PETSC_STATIC_INLINE PetscErrorCode PetscLLCondensedClean(PetscInt lnk_max,PetscInt nidx,PetscInt *indices,PetscInt lnk[],PetscBT bt)
1364: {
1366:   PetscInt       _k,_next,_nlnk;

1369:   _next = lnk[3];       /* head node */
1370:   _nlnk = lnk[0];       /* num of entries on the list */
1371:   for (_k=0; _k<_nlnk; _k++){
1372:     indices[_k] = lnk[_next];
1373:     _next       = lnk[_next + 1];
1374:     PetscBTClear(bt,indices[_k]);
1375:   }
1376:   lnk[0] = 0;          /* num of entries on the list */
1377:   lnk[2] = lnk_max;    /* initialize head node */
1378:   lnk[3] = 2;          /* head node */
1379:   return(0);
1380: }

1382: PETSC_STATIC_INLINE PetscErrorCode PetscLLCondensedView(PetscInt *lnk)
1383: {
1385:   PetscInt       k;

1388:   PetscPrintf(PETSC_COMM_SELF,"LLCondensed of size %D, (val,  next)\n",lnk[0]);
1389:   for (k=2; k< lnk[0]+2; k++){
1390:     PetscPrintf(PETSC_COMM_SELF," %D: (%D, %D)\n",2*k,lnk[2*k],lnk[2*k+1]);
1391:   }
1392:   return(0);
1393: }

1395: /*
1396:   Free memories used by the list
1397: */
1398: PETSC_STATIC_INLINE PetscErrorCode PetscLLCondensedDestroy(PetscInt *lnk,PetscBT bt)
1399: {

1403:   PetscFree(lnk);
1404:   PetscBTDestroy(&bt);
1405:   return(0);
1406: }

1408: /* -------------------------------------------------------------------------------------------------------*/
1409: /*
1410:  Same as PetscLLCondensedCreate(), but does not use non-scalable O(lnk_max) bitarray
1411:   Input Parameters:
1412:     nlnk_max  - max length of the list
1413:   Output Parameters:
1414:     lnk       - list created and initialized
1415: */
1416: PETSC_STATIC_INLINE PetscErrorCode PetscLLCondensedCreate_Scalable(PetscInt nlnk_max,PetscInt **lnk)
1417: {
1419:   PetscInt       *llnk,lsize = 0;

1422:   PetscIntMultError(2,nlnk_max+2,&lsize);
1423:   PetscMalloc1(lsize,lnk);
1424:   llnk = *lnk;
1425:   llnk[0] = 0;               /* number of entries on the list */
1426:   llnk[2] = PETSC_MAX_INT;   /* value in the head node */
1427:   llnk[3] = 2;               /* next for the head node */
1428:   return(0);
1429: }

1431: PETSC_STATIC_INLINE PetscErrorCode PetscLLCondensedExpand_Scalable(PetscInt nlnk_max,PetscInt **lnk)
1432: {
1434:   PetscInt       lsize = 0;

1437:   PetscIntMultError(2,nlnk_max+2,&lsize);
1438:   PetscRealloc(lsize*sizeof(PetscInt),lnk);
1439:   return(0);
1440: }

1442: PETSC_STATIC_INLINE PetscErrorCode PetscLLCondensedAddSorted_Scalable(PetscInt nidx,const PetscInt indices[],PetscInt lnk[])
1443: {
1444:   PetscInt _k,_entry,_location,_next,_lnkdata,_nlnk,_newnode;
1445:   _nlnk     = lnk[0]; /* num of entries on the input lnk */
1446:   _location = 2; /* head */ \
1447:     for (_k=0; _k<nidx; _k++){
1448:       _entry = indices[_k];
1449:       /* search for insertion location */
1450:       do {
1451:         _next     = _location + 1; /* link from previous node to next node */
1452:         _location = lnk[_next];    /* idx of next node */
1453:         _lnkdata  = lnk[_location];/* value of next node */
1454:       } while (_entry > _lnkdata);
1455:       if (_entry < _lnkdata) {
1456:         /* insertion location is found, add entry into lnk */
1457:         _newnode        = 2*(_nlnk+2);   /* index for this new node */
1458:         lnk[_next]      = _newnode;      /* connect previous node to the new node */
1459:         lnk[_newnode]   = _entry;        /* set value of the new node */
1460:         lnk[_newnode+1] = _location;     /* connect new node to next node */
1461:         _location       = _newnode;      /* next search starts from the new node */
1462:         _nlnk++;
1463:       }
1464:     }
1465:   lnk[0]   = _nlnk;   /* number of entries in the list */
1466:   return 0;
1467: }

1469: PETSC_STATIC_INLINE PetscErrorCode PetscLLCondensedClean_Scalable(PetscInt nidx,PetscInt *indices,PetscInt *lnk)
1470: {
1471:   PetscInt _k,_next,_nlnk;
1472:   _next = lnk[3];       /* head node */
1473:   _nlnk = lnk[0];
1474:   for (_k=0; _k<_nlnk; _k++){
1475:     indices[_k] = lnk[_next];
1476:     _next       = lnk[_next + 1];
1477:   }
1478:   lnk[0] = 0;          /* num of entries on the list */
1479:   lnk[3] = 2;          /* head node */
1480:   return 0;
1481: }

1483: PETSC_STATIC_INLINE PetscErrorCode PetscLLCondensedDestroy_Scalable(PetscInt *lnk)
1484: {
1485:   return PetscFree(lnk);
1486: }

1488: /* -------------------------------------------------------------------------------------------------------*/
1489: /*
1490:       lnk[0]   number of links
1491:       lnk[1]   number of entries
1492:       lnk[3n]  value
1493:       lnk[3n+1] len
1494:       lnk[3n+2] link to next value

1496:       The next three are always the first link

1498:       lnk[3]    PETSC_MIN_INT+1
1499:       lnk[4]    1
1500:       lnk[5]    link to first real entry

1502:       The next three are always the last link

1504:       lnk[6]    PETSC_MAX_INT - 1
1505:       lnk[7]    1
1506:       lnk[8]    next valid link (this is the same as lnk[0] but without the decreases)
1507: */

1509: PETSC_STATIC_INLINE PetscErrorCode PetscLLCondensedCreate_fast(PetscInt nlnk_max,PetscInt **lnk)
1510: {
1512:   PetscInt       *llnk,lsize = 0;

1515:   PetscIntMultError(3,nlnk_max+3,&lsize);
1516:   PetscMalloc1(lsize,lnk);
1517:   llnk = *lnk;
1518:   llnk[0] = 0;   /* nlnk: number of entries on the list */
1519:   llnk[1] = 0;          /* number of integer entries represented in list */
1520:   llnk[3] = PETSC_MIN_INT+1;   /* value in the first node */
1521:   llnk[4] = 1;           /* count for the first node */
1522:   llnk[5] = 6;         /* next for the first node */
1523:   llnk[6] = PETSC_MAX_INT-1;   /* value in the last node */
1524:   llnk[7] = 1;           /* count for the last node */
1525:   llnk[8] = 0;         /* next valid node to be used */
1526:   return(0);
1527: }

1529: PETSC_STATIC_INLINE PetscErrorCode PetscLLCondensedAddSorted_fast(PetscInt nidx,const PetscInt indices[],PetscInt lnk[])
1530: {
1531:   PetscInt k,entry,prev,next;
1532:   prev      = 3;      /* first value */
1533:   next      = lnk[prev+2];
1534:   for (k=0; k<nidx; k++){
1535:     entry = indices[k];
1536:     /* search for insertion location */
1537:     while (entry >= lnk[next]) {
1538:       prev = next;
1539:       next = lnk[next+2];
1540:     }
1541:     /* entry is in range of previous list */
1542:     if (entry < lnk[prev]+lnk[prev+1]) continue;
1543:     lnk[1]++;
1544:     /* entry is right after previous list */
1545:     if (entry == lnk[prev]+lnk[prev+1]) {
1546:       lnk[prev+1]++;
1547:       if (lnk[next] == entry+1) { /* combine two contiguous strings */
1548:         lnk[prev+1] += lnk[next+1];
1549:         lnk[prev+2]  = lnk[next+2];
1550:         next         = lnk[next+2];
1551:         lnk[0]--;
1552:       }
1553:       continue;
1554:     }
1555:     /* entry is right before next list */
1556:     if (entry == lnk[next]-1) {
1557:       lnk[next]--;
1558:       lnk[next+1]++;
1559:       prev = next;
1560:       next = lnk[prev+2];
1561:       continue;
1562:     }
1563:     /*  add entry into lnk */
1564:     lnk[prev+2]    = 3*((lnk[8]++)+3);      /* connect previous node to the new node */
1565:     prev           = lnk[prev+2];
1566:     lnk[prev]      = entry;        /* set value of the new node */
1567:     lnk[prev+1]    = 1;             /* number of values in contiguous string is one to start */
1568:     lnk[prev+2]    = next;          /* connect new node to next node */
1569:     lnk[0]++;
1570:   }
1571:   return 0;
1572: }

1574: PETSC_STATIC_INLINE PetscErrorCode PetscLLCondensedClean_fast(PetscInt nidx,PetscInt *indices,PetscInt *lnk)
1575: {
1576:   PetscInt _k,_next,_nlnk,cnt,j;
1577:   _next = lnk[5];       /* first node */
1578:   _nlnk = lnk[0];
1579:   cnt   = 0;
1580:   for (_k=0; _k<_nlnk; _k++){
1581:     for (j=0; j<lnk[_next+1]; j++) {
1582:       indices[cnt++] = lnk[_next] + j;
1583:     }
1584:     _next       = lnk[_next + 2];
1585:   }
1586:   lnk[0] = 0;   /* nlnk: number of links */
1587:   lnk[1] = 0;          /* number of integer entries represented in list */
1588:   lnk[3] = PETSC_MIN_INT+1;   /* value in the first node */
1589:   lnk[4] = 1;           /* count for the first node */
1590:   lnk[5] = 6;         /* next for the first node */
1591:   lnk[6] = PETSC_MAX_INT-1;   /* value in the last node */
1592:   lnk[7] = 1;           /* count for the last node */
1593:   lnk[8] = 0;         /* next valid location to make link */
1594:   return 0;
1595: }

1597: PETSC_STATIC_INLINE PetscErrorCode PetscLLCondensedView_fast(PetscInt *lnk)
1598: {
1599:   PetscInt k,next,nlnk;
1600:   next = lnk[5];       /* first node */
1601:   nlnk = lnk[0];
1602:   for (k=0; k<nlnk; k++){
1603: #if 0                           /* Debugging code */
1604:     printf("%d value %d len %d next %d\n",next,lnk[next],lnk[next+1],lnk[next+2]);
1605: #endif
1606:     next = lnk[next + 2];
1607:   }
1608:   return 0;
1609: }

1611: PETSC_STATIC_INLINE PetscErrorCode PetscLLCondensedDestroy_fast(PetscInt *lnk)
1612: {
1613:   return PetscFree(lnk);
1614: }

1616: /* this is extern because it is used in MatFDColoringUseDM() which is in the DM library */
1617: PETSC_EXTERN PetscErrorCode MatFDColoringApply_AIJ(Mat,MatFDColoring,Vec,void*);

1619: PETSC_EXTERN PetscLogEvent MAT_Mult;
1620: PETSC_EXTERN PetscLogEvent MAT_MultMatrixFree;
1621: PETSC_EXTERN PetscLogEvent MAT_Mults;
1622: PETSC_EXTERN PetscLogEvent MAT_MultConstrained;
1623: PETSC_EXTERN PetscLogEvent MAT_MultAdd;
1624: PETSC_EXTERN PetscLogEvent MAT_MultTranspose;
1625: PETSC_EXTERN PetscLogEvent MAT_MultTransposeConstrained;
1626: PETSC_EXTERN PetscLogEvent MAT_MultTransposeAdd;
1627: PETSC_EXTERN PetscLogEvent MAT_Solve;
1628: PETSC_EXTERN PetscLogEvent MAT_Solves;
1629: PETSC_EXTERN PetscLogEvent MAT_SolveAdd;
1630: PETSC_EXTERN PetscLogEvent MAT_SolveTranspose;
1631: PETSC_EXTERN PetscLogEvent MAT_SolveTransposeAdd;
1632: PETSC_EXTERN PetscLogEvent MAT_SOR;
1633: PETSC_EXTERN PetscLogEvent MAT_ForwardSolve;
1634: PETSC_EXTERN PetscLogEvent MAT_BackwardSolve;
1635: PETSC_EXTERN PetscLogEvent MAT_LUFactor;
1636: PETSC_EXTERN PetscLogEvent MAT_LUFactorSymbolic;
1637: PETSC_EXTERN PetscLogEvent MAT_LUFactorNumeric;
1638: PETSC_EXTERN PetscLogEvent MAT_CholeskyFactor;
1639: PETSC_EXTERN PetscLogEvent MAT_CholeskyFactorSymbolic;
1640: PETSC_EXTERN PetscLogEvent MAT_CholeskyFactorNumeric;
1641: PETSC_EXTERN PetscLogEvent MAT_ILUFactor;
1642: PETSC_EXTERN PetscLogEvent MAT_ILUFactorSymbolic;
1643: PETSC_EXTERN PetscLogEvent MAT_ICCFactorSymbolic;
1644: PETSC_EXTERN PetscLogEvent MAT_Copy;
1645: PETSC_EXTERN PetscLogEvent MAT_Convert;
1646: PETSC_EXTERN PetscLogEvent MAT_Scale;
1647: PETSC_EXTERN PetscLogEvent MAT_AssemblyBegin;
1648: PETSC_EXTERN PetscLogEvent MAT_AssemblyEnd;
1649: PETSC_EXTERN PetscLogEvent MAT_SetValues;
1650: PETSC_EXTERN PetscLogEvent MAT_GetValues;
1651: PETSC_EXTERN PetscLogEvent MAT_GetRow;
1652: PETSC_EXTERN PetscLogEvent MAT_GetRowIJ;
1653: PETSC_EXTERN PetscLogEvent MAT_CreateSubMats;
1654: PETSC_EXTERN PetscLogEvent MAT_GetColoring;
1655: PETSC_EXTERN PetscLogEvent MAT_GetOrdering;
1656: PETSC_EXTERN PetscLogEvent MAT_RedundantMat;
1657: PETSC_EXTERN PetscLogEvent MAT_IncreaseOverlap;
1658: PETSC_EXTERN PetscLogEvent MAT_Partitioning;
1659: PETSC_EXTERN PetscLogEvent MAT_PartitioningND;
1660: PETSC_EXTERN PetscLogEvent MAT_Coarsen;
1661: PETSC_EXTERN PetscLogEvent MAT_ZeroEntries;
1662: PETSC_EXTERN PetscLogEvent MAT_Load;
1663: PETSC_EXTERN PetscLogEvent MAT_View;
1664: PETSC_EXTERN PetscLogEvent MAT_AXPY;
1665: PETSC_EXTERN PetscLogEvent MAT_FDColoringCreate;
1666: PETSC_EXTERN PetscLogEvent MAT_TransposeColoringCreate;
1667: PETSC_EXTERN PetscLogEvent MAT_FDColoringSetUp;
1668: PETSC_EXTERN PetscLogEvent MAT_FDColoringApply;
1669: PETSC_EXTERN PetscLogEvent MAT_Transpose;
1670: PETSC_EXTERN PetscLogEvent MAT_FDColoringFunction;
1671: PETSC_EXTERN PetscLogEvent MAT_CreateSubMat;
1672: PETSC_EXTERN PetscLogEvent MAT_MatMult;
1673: PETSC_EXTERN PetscLogEvent MAT_MatSolve;
1674: PETSC_EXTERN PetscLogEvent MAT_MatTrSolve;
1675: PETSC_EXTERN PetscLogEvent MAT_MatMultSymbolic;
1676: PETSC_EXTERN PetscLogEvent MAT_MatMultNumeric;
1677: PETSC_EXTERN PetscLogEvent MAT_Getlocalmatcondensed;
1678: PETSC_EXTERN PetscLogEvent MAT_GetBrowsOfAcols;
1679: PETSC_EXTERN PetscLogEvent MAT_GetBrowsOfAocols;
1680: PETSC_EXTERN PetscLogEvent MAT_PtAP;
1681: PETSC_EXTERN PetscLogEvent MAT_PtAPSymbolic;
1682: PETSC_EXTERN PetscLogEvent MAT_PtAPNumeric;
1683: PETSC_EXTERN PetscLogEvent MAT_Seqstompinum;
1684: PETSC_EXTERN PetscLogEvent MAT_Seqstompisym;
1685: PETSC_EXTERN PetscLogEvent MAT_Seqstompi;
1686: PETSC_EXTERN PetscLogEvent MAT_Getlocalmat;
1687: PETSC_EXTERN PetscLogEvent MAT_RARt;
1688: PETSC_EXTERN PetscLogEvent MAT_RARtSymbolic;
1689: PETSC_EXTERN PetscLogEvent MAT_RARtNumeric;
1690: PETSC_EXTERN PetscLogEvent MAT_MatTransposeMult;
1691: PETSC_EXTERN PetscLogEvent MAT_MatTransposeMultSymbolic;
1692: PETSC_EXTERN PetscLogEvent MAT_MatTransposeMultNumeric;
1693: PETSC_EXTERN PetscLogEvent MAT_TransposeMatMult;
1694: PETSC_EXTERN PetscLogEvent MAT_TransposeMatMultSymbolic;
1695: PETSC_EXTERN PetscLogEvent MAT_TransposeMatMultNumeric;
1696: PETSC_EXTERN PetscLogEvent MAT_MatMatMult;
1697: PETSC_EXTERN PetscLogEvent MAT_MatMatMultSymbolic;
1698: PETSC_EXTERN PetscLogEvent MAT_MatMatMultNumeric;
1699: PETSC_EXTERN PetscLogEvent MAT_Applypapt;
1700: PETSC_EXTERN PetscLogEvent MAT_Applypapt_symbolic;
1701: PETSC_EXTERN PetscLogEvent MAT_Applypapt_numeric;
1702: PETSC_EXTERN PetscLogEvent MAT_Getsymtranspose;
1703: PETSC_EXTERN PetscLogEvent MAT_Transpose_SeqAIJ;
1704: PETSC_EXTERN PetscLogEvent MAT_Getsymtransreduced;
1705: PETSC_EXTERN PetscLogEvent MAT_GetSequentialNonzeroStructure;
1706: PETSC_EXTERN PetscLogEvent MATMFFD_Mult;
1707: PETSC_EXTERN PetscLogEvent MAT_GetMultiProcBlock;
1708: PETSC_EXTERN PetscLogEvent MAT_CUSPARSECopyToGPU;
1709: PETSC_EXTERN PetscLogEvent MAT_SetValuesBatch;
1710: PETSC_EXTERN PetscLogEvent MAT_ViennaCLCopyToGPU;
1711: PETSC_EXTERN PetscLogEvent MAT_Merge;
1712: PETSC_EXTERN PetscLogEvent MAT_Residual;
1713: PETSC_EXTERN PetscLogEvent MAT_SetRandom;
1714: PETSC_EXTERN PetscLogEvent MATCOLORING_Apply;
1715: PETSC_EXTERN PetscLogEvent MATCOLORING_Comm;
1716: PETSC_EXTERN PetscLogEvent MATCOLORING_Local;
1717: PETSC_EXTERN PetscLogEvent MATCOLORING_ISCreate;
1718: PETSC_EXTERN PetscLogEvent MATCOLORING_SetUp;
1719: PETSC_EXTERN PetscLogEvent MATCOLORING_Weights;

1721: #endif