Actual source code: matimpl.h

petsc-3.13.6 2020-09-29
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 src/mat/f90-mod/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 (*freeintermediatedatastructures)(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 (*placeholder_63)(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 (*placeholder_89)(Mat,void*);
142:   PetscErrorCode (*matmultsymbolic)(Mat,Mat,PetscReal,Mat);
143:   PetscErrorCode (*matmultnumeric)(Mat,Mat,Mat);
144:   PetscErrorCode (*placeholder_92)(Mat,void*);
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 (*placeholder_95)(Mat,void*);
149:   PetscErrorCode (*mattransposemultsymbolic)(Mat,Mat,PetscReal,Mat);
150:   PetscErrorCode (*mattransposemultnumeric)(Mat,Mat,Mat);
151:   PetscErrorCode (*bindtocpu)(Mat,PetscBool);
152:   /*99*/
153:   PetscErrorCode (*productsetfromoptions)(Mat);
154:   PetscErrorCode (*productsymbolic)(Mat);
155:   PetscErrorCode (*productnumeric)(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 (*placeholder_130)(Mat,void*);
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 (*placeholder_136)(Mat,void*);
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:   PetscErrorCode (*creatempimatconcatenateseqmat)(MPI_Comm,Mat,PetscInt,MatReuse,Mat*);
207:   /*145*/
208:   PetscErrorCode (*destroysubmatrices)(PetscInt,Mat*[]);
209:   PetscErrorCode (*mattransposesolve)(Mat,Mat,Mat);
210:   PetscErrorCode (*getvalueslocal)(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],PetscScalar[]);
211: };
212: /*
213:     If you add MatOps entries above also add them to the MATOP enum
214:     in include/petscmat.h and src/mat/f90-mod/petscmat.h
215: */

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

221: typedef struct _p_MatRootName* MatRootName;
222: struct _p_MatRootName {
223:   char        *rname,*sname,*mname;
224:   MatRootName next;
225: };

227: PETSC_EXTERN MatRootName MatRootNameList;

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

239: PETSC_INTERN PetscErrorCode MatProductSymbolic_Basic(Mat);
240: PETSC_EXTERN PetscErrorCode MatProductSymbolic_AB(Mat);
241: PETSC_EXTERN PetscErrorCode MatProductNumeric_AB(Mat);
242: PETSC_EXTERN PetscErrorCode MatProductSymbolic_AtB(Mat);
243: PETSC_EXTERN PetscErrorCode MatProductNumeric_AtB(Mat);
244: PETSC_EXTERN PetscErrorCode MatProductSymbolic_ABt(Mat);
245: PETSC_EXTERN PetscErrorCode MatProductNumeric_ABt(Mat);
246: PETSC_EXTERN PetscErrorCode MatProductNumeric_PtAP(Mat);
247: PETSC_EXTERN PetscErrorCode MatProductNumeric_RARt(Mat);
248: PETSC_EXTERN PetscErrorCode MatProductSymbolic_ABC(Mat);
249: PETSC_EXTERN PetscErrorCode MatProductNumeric_ABC(Mat);
250: PETSC_EXTERN PetscErrorCode MatProductCreate_Private(Mat,Mat,Mat,Mat);

252: #if defined(PETSC_USE_DEBUG)
253: #  define MatCheckPreallocated(A,arg) do {                              \
254:     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); \
255:   } while (0)
256: #else
257: #  define MatCheckPreallocated(A,arg) do {} while (0)
258: #endif

260: /*
261:   The stash is used to temporarily store inserted matrix values that
262:   belong to another processor. During the assembly phase the stashed
263:   values are moved to the correct processor and
264: */

266: typedef struct _MatStashSpace *PetscMatStashSpace;

268: struct _MatStashSpace {
269:   PetscMatStashSpace next;
270:   PetscScalar        *space_head,*val;
271:   PetscInt           *idx,*idy;
272:   PetscInt           total_space_size;
273:   PetscInt           local_used;
274:   PetscInt           local_remaining;
275: };

277: PETSC_EXTERN PetscErrorCode PetscMatStashSpaceGet(PetscInt,PetscInt,PetscMatStashSpace *);
278: PETSC_EXTERN PetscErrorCode PetscMatStashSpaceContiguous(PetscInt,PetscMatStashSpace *,PetscScalar *,PetscInt *,PetscInt *);
279: PETSC_EXTERN PetscErrorCode PetscMatStashSpaceDestroy(PetscMatStashSpace*);

281: typedef struct {
282:   PetscInt    count;
283: } MatStashHeader;

285: typedef struct {
286:   void        *buffer;          /* Of type blocktype, dynamically constructed  */
287:   PetscInt    count;
288:   char        pending;
289: } MatStashFrame;

291: typedef struct _MatStash MatStash;
292: struct _MatStash {
293:   PetscInt      nmax;                   /* maximum stash size */
294:   PetscInt      umax;                   /* user specified max-size */
295:   PetscInt      oldnmax;                /* the nmax value used previously */
296:   PetscInt      n;                      /* stash size */
297:   PetscInt      bs;                     /* block size of the stash */
298:   PetscInt      reallocs;               /* preserve the no of mallocs invoked */
299:   PetscMatStashSpace space_head,space;  /* linked list to hold stashed global row/column numbers and matrix values */

301:   PetscErrorCode (*ScatterBegin)(Mat,MatStash*,PetscInt*);
302:   PetscErrorCode (*ScatterGetMesg)(MatStash*,PetscMPIInt*,PetscInt**,PetscInt**,PetscScalar**,PetscInt*);
303:   PetscErrorCode (*ScatterEnd)(MatStash*);
304:   PetscErrorCode (*ScatterDestroy)(MatStash*);

306:   /* The following variables are used for communication */
307:   MPI_Comm      comm;
308:   PetscMPIInt   size,rank;
309:   PetscMPIInt   tag1,tag2;
310:   MPI_Request   *send_waits;            /* array of send requests */
311:   MPI_Request   *recv_waits;            /* array of receive requests */
312:   MPI_Status    *send_status;           /* array of send status */
313:   PetscInt      nsends,nrecvs;          /* numbers of sends and receives */
314:   PetscScalar   *svalues;               /* sending data */
315:   PetscInt      *sindices;
316:   PetscScalar   **rvalues;              /* receiving data (values) */
317:   PetscInt      **rindices;             /* receiving data (indices) */
318:   PetscInt      nprocessed;             /* number of messages already processed */
319:   PetscMPIInt   *flg_v;                 /* indicates what messages have arrived so far and from whom */
320:   PetscBool     reproduce;
321:   PetscInt      reproduce_count;

323:   /* The following variables are used for BTS communication */
324:   PetscBool      first_assembly_done;   /* Is the first time matrix assembly done? */
325:   PetscBool      use_status;            /* Use MPI_Status to determine number of items in each message */
326:   PetscMPIInt    nsendranks;
327:   PetscMPIInt    nrecvranks;
328:   PetscMPIInt    *sendranks;
329:   PetscMPIInt    *recvranks;
330:   MatStashHeader *sendhdr,*recvhdr;
331:   MatStashFrame  *sendframes;   /* pointers to the main messages */
332:   MatStashFrame  *recvframes;
333:   MatStashFrame  *recvframe_active;
334:   PetscInt       recvframe_i;     /* index of block within active frame */
335:   PetscMPIInt    recvframe_count; /* Count actually sent for current frame */
336:   PetscInt       recvcount;       /* Number of receives processed so far */
337:   PetscMPIInt    *some_indices;   /* From last call to MPI_Waitsome */
338:   MPI_Status     *some_statuses;  /* Statuses from last call to MPI_Waitsome */
339:   PetscMPIInt    some_count;      /* Number of requests completed in last call to MPI_Waitsome */
340:   PetscMPIInt    some_i;          /* Index of request currently being processed */
341:   MPI_Request    *sendreqs;
342:   MPI_Request    *recvreqs;
343:   PetscSegBuffer segsendblocks;
344:   PetscSegBuffer segrecvframe;
345:   PetscSegBuffer segrecvblocks;
346:   MPI_Datatype   blocktype;
347:   size_t         blocktype_size;
348:   InsertMode     *insertmode;   /* Pointer to check mat->insertmode and set upon message arrival in case no local values have been set. */
349: };

351: #if !defined(PETSC_HAVE_MPIUNI)
352: PETSC_INTERN PetscErrorCode MatStashScatterDestroy_BTS(MatStash*);
353: #endif
354: PETSC_INTERN PetscErrorCode MatStashCreate_Private(MPI_Comm,PetscInt,MatStash*);
355: PETSC_INTERN PetscErrorCode MatStashDestroy_Private(MatStash*);
356: PETSC_INTERN PetscErrorCode MatStashScatterEnd_Private(MatStash*);
357: PETSC_INTERN PetscErrorCode MatStashSetInitialSize_Private(MatStash*,PetscInt);
358: PETSC_INTERN PetscErrorCode MatStashGetInfo_Private(MatStash*,PetscInt*,PetscInt*);
359: PETSC_INTERN PetscErrorCode MatStashValuesRow_Private(MatStash*,PetscInt,PetscInt,const PetscInt[],const PetscScalar[],PetscBool );
360: PETSC_INTERN PetscErrorCode MatStashValuesCol_Private(MatStash*,PetscInt,PetscInt,const PetscInt[],const PetscScalar[],PetscInt,PetscBool );
361: PETSC_INTERN PetscErrorCode MatStashValuesRowBlocked_Private(MatStash*,PetscInt,PetscInt,const PetscInt[],const PetscScalar[],PetscInt,PetscInt,PetscInt);
362: PETSC_INTERN PetscErrorCode MatStashValuesColBlocked_Private(MatStash*,PetscInt,PetscInt,const PetscInt[],const PetscScalar[],PetscInt,PetscInt,PetscInt);
363: PETSC_INTERN PetscErrorCode MatStashScatterBegin_Private(Mat,MatStash*,PetscInt*);
364: PETSC_INTERN PetscErrorCode MatStashScatterGetMesg_Private(MatStash*,PetscMPIInt*,PetscInt**,PetscInt**,PetscScalar**,PetscInt*);
365: PETSC_INTERN PetscErrorCode MatGetInfo_External(Mat,MatInfoType,MatInfo*);

367: typedef struct {
368:   PetscInt   dim;
369:   PetscInt   dims[4];
370:   PetscInt   starts[4];
371:   PetscBool  noc;        /* this is a single component problem, hence user will not set MatStencil.c */
372: } MatStencilInfo;

374: /* Info about using compressed row format */
375: typedef struct {
376:   PetscBool  use;                           /* indicates compressed rows have been checked and will be used */
377:   PetscInt   nrows;                         /* number of non-zero rows */
378:   PetscInt   *i;                            /* compressed row pointer  */
379:   PetscInt   *rindex;                       /* compressed row index               */
380: } Mat_CompressedRow;
381: PETSC_EXTERN PetscErrorCode MatCheckCompressedRow(Mat,PetscInt,Mat_CompressedRow*,PetscInt*,PetscInt,PetscReal);

383: typedef struct { /* used by MatCreateRedundantMatrix() for reusing matredundant */
384:   PetscInt     nzlocal,nsends,nrecvs;
385:   PetscMPIInt  *send_rank,*recv_rank;
386:   PetscInt     *sbuf_nz,*rbuf_nz,*sbuf_j,**rbuf_j;
387:   PetscScalar  *sbuf_a,**rbuf_a;
388:   MPI_Comm     subcomm;   /* when user does not provide a subcomm */
389:   IS           isrow,iscol;
390:   Mat          *matseq;
391: } Mat_Redundant;

393: typedef struct { /* used by MatProduct() */
394:   MatProductType       type;
395:   MatProductAlgorithm  alg;
396:   Mat                  A,B,C,Dwork;
397:   PetscReal            fill;
398:   PetscBool            Areplaced,Breplaced; /* if an internal implementation changes user's input A or B, these matrices cannot be called by MatProductReplaceMats(). */
399:   PetscBool            api_user; /* used by MatProductSetFromOptions_xxx() */
400: } Mat_Product;

402: struct _p_Mat {
403:   PETSCHEADER(struct _MatOps);
404:   PetscLayout            rmap,cmap;
405:   void                   *data;            /* implementation-specific data */
406:   MatFactorType          factortype;       /* MAT_FACTOR_LU, ILU, CHOLESKY or ICC */
407:   PetscBool              assembled;        /* is the matrix assembled? */
408:   PetscBool              was_assembled;    /* new values inserted into assembled mat */
409:   PetscInt               num_ass;          /* number of times matrix has been assembled */
410:   PetscObjectState       nonzerostate;     /* each time new nonzeros locations are introduced into the matrix this is updated */
411:   PetscObjectState       ass_nonzerostate; /* nonzero state at last assembly */
412:   MatInfo                info;             /* matrix information */
413:   InsertMode             insertmode;       /* have values been inserted in matrix or added? */
414:   MatStash               stash,bstash;     /* used for assembling off-proc mat emements */
415:   MatNullSpace           nullsp;           /* null space (operator is singular) */
416:   MatNullSpace           transnullsp;      /* null space of transpose of operator */
417:   MatNullSpace           nearnullsp;       /* near null space to be used by multigrid methods */
418:   PetscInt               congruentlayouts; /* are the rows and columns layouts congruent? */
419:   PetscBool              preallocated;
420:   MatStencilInfo         stencil;          /* information for structured grid */
421:   PetscBool              symmetric,hermitian,structurally_symmetric,spd;
422:   PetscBool              symmetric_set,hermitian_set,structurally_symmetric_set,spd_set; /* if true, then corresponding flag is correct*/
423:   PetscBool              symmetric_eternal;
424:   PetscBool              nooffprocentries,nooffproczerorows;
425:   PetscBool              assembly_subset;  /* set by MAT_SUBSET_OFF_PROC_ENTRIES */
426:   PetscBool              submat_singleis;  /* for efficient PCSetUp_ASM() */
427:   PetscBool              structure_only;
428:   PetscBool              sortedfull;       /* full, sorted rows are inserted */
429: #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_CUDA)
430:   PetscOffloadMask       offloadmask;      /* a mask which indicates where the valid matrix data is (GPU, CPU or both) */
431:   PetscBool              boundtocpu;
432: #endif
433:   void                   *spptr;          /* pointer for special library like SuperLU */
434:   char                   *solvertype;
435:   PetscBool              checksymmetryonassembly,checknullspaceonassembly;
436:   PetscReal              checksymmetrytol;
437:   Mat                    schur;             /* Schur complement matrix */
438:   MatFactorSchurStatus   schur_status;      /* status of the Schur complement matrix */
439:   Mat_Redundant          *redundant;        /* used by MatCreateRedundantMatrix() */
440:   PetscBool              erroriffailure;    /* Generate an error if detected (for example a zero pivot) instead of returning */
441:   MatFactorError         factorerrortype;               /* type of error in factorization */
442:   PetscReal              factorerror_zeropivot_value;   /* If numerical zero pivot was detected this is the computed value */
443:   PetscInt               factorerror_zeropivot_row;     /* Row where zero pivot was detected */
444:   PetscInt               nblocks,*bsizes;   /* support for MatSetVariableBlockSizes() */
445:   char                   *defaultvectype;
446:   Mat_Product            *product;
447: };

449: PETSC_INTERN PetscErrorCode MatAXPY_Basic(Mat,PetscScalar,Mat,MatStructure);
450: PETSC_INTERN PetscErrorCode MatAXPY_BasicWithPreallocation(Mat,Mat,PetscScalar,Mat,MatStructure);
451: PETSC_INTERN PetscErrorCode MatAXPY_Basic_Preallocate(Mat,Mat,Mat*);

453: /*
454:     Utility for MatFactor (Schur complement)
455: */
456: PETSC_INTERN PetscErrorCode MatFactorFactorizeSchurComplement_Private(Mat);
457: PETSC_INTERN PetscErrorCode MatFactorInvertSchurComplement_Private(Mat);
458: PETSC_INTERN PetscErrorCode MatFactorUpdateSchurStatus_Private(Mat);
459: PETSC_INTERN PetscErrorCode MatFactorSetUpInPlaceSchur_Private(Mat);

461: /*
462:     Utility for MatZeroRows
463: */
464: PETSC_INTERN PetscErrorCode MatZeroRowsMapLocal_Private(Mat,PetscInt,const PetscInt*,PetscInt*,PetscInt**);

466: /*
467:     Utility for MatView/MatLoad
468: */
469: PETSC_INTERN PetscErrorCode MatView_Binary_BlockSizes(Mat,PetscViewer);
470: PETSC_INTERN PetscErrorCode MatLoad_Binary_BlockSizes(Mat,PetscViewer);


473: /*
474:     Object for partitioning graphs
475: */

477: typedef struct _MatPartitioningOps *MatPartitioningOps;
478: struct _MatPartitioningOps {
479:   PetscErrorCode (*apply)(MatPartitioning,IS*);
480:   PetscErrorCode (*applynd)(MatPartitioning,IS*);
481:   PetscErrorCode (*setfromoptions)(PetscOptionItems*,MatPartitioning);
482:   PetscErrorCode (*destroy)(MatPartitioning);
483:   PetscErrorCode (*view)(MatPartitioning,PetscViewer);
484:   PetscErrorCode (*improve)(MatPartitioning,IS*);
485: };

487: struct _p_MatPartitioning {
488:   PETSCHEADER(struct _MatPartitioningOps);
489:   Mat         adj;
490:   PetscInt    *vertex_weights;
491:   PetscReal   *part_weights;
492:   PetscInt    n;                                 /* number of partitions */
493:   void        *data;
494:   PetscInt    setupcalled;
495:   PetscBool   use_edge_weights;  /* A flag indicates whether or not to use edge weights */
496: };

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

501: /*
502:     Object for coarsen graphs
503: */
504: typedef struct _MatCoarsenOps *MatCoarsenOps;
505: struct _MatCoarsenOps {
506:   PetscErrorCode (*apply)(MatCoarsen);
507:   PetscErrorCode (*setfromoptions)(PetscOptionItems*,MatCoarsen);
508:   PetscErrorCode (*destroy)(MatCoarsen);
509:   PetscErrorCode (*view)(MatCoarsen,PetscViewer);
510: };

512: struct _p_MatCoarsen {
513:   PETSCHEADER(struct _MatCoarsenOps);
514:   Mat              graph;
515:   PetscInt         setupcalled;
516:   void             *subctx;
517:   /* */
518:   PetscBool        strict_aggs;
519:   IS               perm;
520:   PetscCoarsenData *agg_lists;
521: };

523: /*
524:     MatFDColoring is used to compute Jacobian matrices efficiently
525:   via coloring. The data structure is explained below in an example.

527:    Color =   0    1     0    2   |   2      3       0
528:    ---------------------------------------------------
529:             00   01              |          05
530:             10   11              |   14     15               Processor  0
531:                        22    23  |          25
532:                        32    33  |
533:    ===================================================
534:                                  |   44     45     46
535:             50                   |          55               Processor 1
536:                                  |   64            66
537:    ---------------------------------------------------

539:     ncolors = 4;

541:     ncolumns      = {2,1,1,0}
542:     columns       = {{0,2},{1},{3},{}}
543:     nrows         = {4,2,3,3}
544:     rows          = {{0,1,2,3},{0,1},{1,2,3},{0,1,2}}
545:     vwscale       = {dx(0),dx(1),dx(2),dx(3)}               MPI Vec
546:     vscale        = {dx(0),dx(1),dx(2),dx(3),dx(4),dx(5)}   Seq Vec

548:     ncolumns      = {1,0,1,1}
549:     columns       = {{6},{},{4},{5}}
550:     nrows         = {3,0,2,2}
551:     rows          = {{0,1,2},{},{1,2},{1,2}}
552:     vwscale       = {dx(4),dx(5),dx(6)}              MPI Vec
553:     vscale        = {dx(0),dx(4),dx(5),dx(6)}        Seq Vec

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

558: */
559: typedef struct {
560:   PetscInt     row;
561:   PetscInt     col;
562:   PetscScalar  *valaddr;   /* address of value */
563: } MatEntry;

565: typedef struct {
566:   PetscInt     row;
567:   PetscScalar  *valaddr;   /* address of value */
568: } MatEntry2;

570: struct  _p_MatFDColoring{
571:   PETSCHEADER(int);
572:   PetscInt       M,N,m;            /* total rows, columns; local rows */
573:   PetscInt       rstart;           /* first row owned by local processor */
574:   PetscInt       ncolors;          /* number of colors */
575:   PetscInt       *ncolumns;        /* number of local columns for a color */
576:   PetscInt       **columns;        /* lists the local columns of each color (using global column numbering) */
577:   IS             *isa;             /* these are the IS that contain the column values given in columns */
578:   PetscInt       *nrows;           /* number of local rows for each color */
579:   MatEntry       *matentry;        /* holds (row, column, address of value) for Jacobian matrix entry */
580:   MatEntry2      *matentry2;       /* holds (row, address of value) for Jacobian matrix entry */
581:   PetscScalar    *dy;              /* store a block of F(x+dx)-F(x) when J is in BAIJ format */
582:   PetscReal      error_rel;        /* square root of relative error in computing function */
583:   PetscReal      umin;             /* minimum allowable u'dx value */
584:   Vec            w1,w2,w3;         /* work vectors used in computing Jacobian */
585:   PetscBool      fset;             /* indicates that the initial function value F(X) is set */
586:   PetscErrorCode (*f)(void);       /* function that defines Jacobian */
587:   void           *fctx;            /* optional user-defined context for use by the function f */
588:   Vec            vscale;           /* holds FD scaling, i.e. 1/dx for each perturbed column */
589:   PetscInt       currentcolor;     /* color for which function evaluation is being done now */
590:   const char     *htype;           /* "wp" or "ds" */
591:   ISColoringType ctype;            /* IS_COLORING_GLOBAL or IS_COLORING_LOCAL */
592:   PetscInt       brows,bcols;      /* number of block rows or columns for speedup inserting the dense matrix into sparse Jacobian */
593:   PetscBool      setupcalled;      /* true if setup has been called */
594:   PetscBool      viewed;           /* true if the -mat_fd_coloring_view has been triggered already */
595:   void           (*ftn_func_pointer)(void),*ftn_func_cntx; /* serve the same purpose as *fortran_func_pointers in PETSc objects */
596:   PetscObjectId  matid;            /* matrix this object was created with, must always be the same */
597: };

599: typedef struct _MatColoringOps *MatColoringOps;
600: struct _MatColoringOps {
601:   PetscErrorCode (*destroy)(MatColoring);
602:   PetscErrorCode (*setfromoptions)(PetscOptionItems*,MatColoring);
603:   PetscErrorCode (*view)(MatColoring,PetscViewer);
604:   PetscErrorCode (*apply)(MatColoring,ISColoring*);
605:   PetscErrorCode (*weights)(MatColoring,PetscReal**,PetscInt**);
606: };

608: struct _p_MatColoring {
609:   PETSCHEADER(struct _MatColoringOps);
610:   Mat                   mat;
611:   PetscInt              dist;             /* distance of the coloring */
612:   PetscInt              maxcolors;        /* the maximum number of colors returned, maxcolors=1 for MIS */
613:   void                  *data;            /* inner context */
614:   PetscBool             valid;            /* check to see if what is produced is a valid coloring */
615:   MatColoringWeightType weight_type;      /* type of weight computation to be performed */
616:   PetscReal             *user_weights;    /* custom weights and permutation */
617:   PetscInt              *user_lperm;
618:   PetscBool             valid_iscoloring; /* check to see if matcoloring is produced a valid iscoloring */
619: };

621: struct  _p_MatTransposeColoring{
622:   PETSCHEADER(int);
623:   PetscInt       M,N,m;            /* total rows, columns; local rows */
624:   PetscInt       rstart;           /* first row owned by local processor */
625:   PetscInt       ncolors;          /* number of colors */
626:   PetscInt       *ncolumns;        /* number of local columns for a color */
627:   PetscInt       *nrows;           /* number of local rows for each color */
628:   PetscInt       currentcolor;     /* color for which function evaluation is being done now */
629:   ISColoringType ctype;            /* IS_COLORING_GLOBAL or IS_COLORING_LOCAL */

631:   PetscInt       *colorforrow,*colorforcol;  /* pointer to rows and columns */
632:   PetscInt       *rows;                      /* lists the local rows for each color (using the local row numbering) */
633:   PetscInt       *den2sp;                    /* maps (row,color) in the dense matrix to index of sparse matrix array a->a */
634:   PetscInt       *columns;                   /* lists the local columns of each color (using global column numbering) */
635:   PetscInt       brows;                      /* number of rows for efficient implementation of MatTransColoringApplyDenToSp() */
636:   PetscInt       *lstart;                    /* array used for loop over row blocks of Csparse */
637: };

639: /*
640:    Null space context for preconditioner/operators
641: */
642: struct _p_MatNullSpace {
643:   PETSCHEADER(int);
644:   PetscBool      has_cnst;
645:   PetscInt       n;
646:   Vec*           vecs;
647:   PetscScalar*   alpha;                 /* for projections */
648:   PetscErrorCode (*remove)(MatNullSpace,Vec,void*);  /* for user provided removal function */
649:   void*          rmctx;                 /* context for remove() function */
650: };

652: /*
653:    Checking zero pivot for LU, ILU preconditioners.
654: */
655: typedef struct {
656:   PetscInt       nshift,nshift_max;
657:   PetscReal      shift_amount,shift_lo,shift_hi,shift_top,shift_fraction;
658:   PetscBool      newshift;
659:   PetscReal      rs;  /* active row sum of abs(offdiagonals) */
660:   PetscScalar    pv;  /* pivot of the active row */
661: } FactorShiftCtx;

663: /*
664:  Used by MatCreateSubMatrices_MPIXAIJ_Local()
665: */
666:  #include <petscctable.h>
667: typedef struct { /* used by MatCreateSubMatrices_MPIAIJ_SingleIS_Local() and MatCreateSubMatrices_MPIAIJ_Local */
668:   PetscInt   id;   /* index of submats, only submats[0] is responsible for deleting some arrays below */
669:   PetscInt   nrqs,nrqr;
670:   PetscInt   **rbuf1,**rbuf2,**rbuf3,**sbuf1,**sbuf2;
671:   PetscInt   **ptr;
672:   PetscInt   *tmp;
673:   PetscInt   *ctr;
674:   PetscInt   *pa; /* proc array */
675:   PetscInt   *req_size,*req_source1,*req_source2;
676:   PetscBool  allcolumns,allrows;
677:   PetscBool  singleis;
678:   PetscInt   *row2proc; /* row to proc map */
679:   PetscInt   nstages;
680: #if defined(PETSC_USE_CTABLE)
681:   PetscTable cmap,rmap;
682:   PetscInt   *cmap_loc,*rmap_loc;
683: #else
684:   PetscInt   *cmap,*rmap;
685: #endif

687:   PetscErrorCode (*destroy)(Mat);
688: } Mat_SubSppt;

690: PETSC_EXTERN PetscErrorCode MatFactorDumpMatrix(Mat);
691: PETSC_INTERN PetscErrorCode MatShift_Basic(Mat,PetscScalar);
692: PETSC_INTERN PetscErrorCode MatSetBlockSizes_Default(Mat,PetscInt,PetscInt);

694: PETSC_STATIC_INLINE PetscErrorCode MatPivotCheck_nz(Mat mat,const MatFactorInfo *info,FactorShiftCtx *sctx,PetscInt row)
695: {
696:   PetscReal _rs   = sctx->rs;
697:   PetscReal _zero = info->zeropivot*_rs;

700:   if (PetscAbsScalar(sctx->pv) <= _zero && !PetscIsNanScalar(sctx->pv)){
701:     /* force |diag| > zeropivot*rs */
702:     if (!sctx->nshift) sctx->shift_amount = info->shiftamount;
703:     else sctx->shift_amount *= 2.0;
704:     sctx->newshift = PETSC_TRUE;
705:     (sctx->nshift)++;
706:   } else {
707:     sctx->newshift = PETSC_FALSE;
708:   }
709:   return(0);
710: }

712: PETSC_STATIC_INLINE PetscErrorCode MatPivotCheck_pd(Mat mat,const MatFactorInfo *info,FactorShiftCtx *sctx,PetscInt row)
713: {
714:   PetscReal _rs   = sctx->rs;
715:   PetscReal _zero = info->zeropivot*_rs;

718:   if (PetscRealPart(sctx->pv) <= _zero && !PetscIsNanScalar(sctx->pv)){
719:     /* force matfactor to be diagonally dominant */
720:     if (sctx->nshift == sctx->nshift_max) {
721:       sctx->shift_fraction = sctx->shift_hi;
722:     } else {
723:       sctx->shift_lo = sctx->shift_fraction;
724:       sctx->shift_fraction = (sctx->shift_hi+sctx->shift_lo)/2.;
725:     }
726:     sctx->shift_amount = sctx->shift_fraction * sctx->shift_top;
727:     sctx->nshift++;
728:     sctx->newshift = PETSC_TRUE;
729:   } else {
730:     sctx->newshift = PETSC_FALSE;
731:   }
732:   return(0);
733: }

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

740:   if (PetscAbsScalar(sctx->pv) <= _zero && !PetscIsNanScalar(sctx->pv)){
741:     sctx->pv          += info->shiftamount;
742:     sctx->shift_amount = 0.0;
743:     sctx->nshift++;
744:   }
745:   sctx->newshift = PETSC_FALSE;
746:   return(0);
747: }

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

755:   sctx->newshift = PETSC_FALSE;
756:   if (PetscAbsScalar(sctx->pv) <= _zero && !PetscIsNanScalar(sctx->pv)) {
757:     if (!mat->erroriffailure) {
758:       PetscInfo3(mat,"Detected zero pivot in factorization in row %D value %g tolerance %g\n",row,(double)PetscAbsScalar(sctx->pv),(double)_zero);
759:       fact->factorerrortype             = MAT_FACTOR_NUMERIC_ZEROPIVOT;
760:       fact->factorerror_zeropivot_value = PetscAbsScalar(sctx->pv);
761:       fact->factorerror_zeropivot_row   = row;
762:     } 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);
763:   }
764:   return(0);
765: }

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

772:   if (info->shifttype == (PetscReal) MAT_SHIFT_NONZERO){
773:     MatPivotCheck_nz(mat,info,sctx,row);
774:   } else if (info->shifttype == (PetscReal) MAT_SHIFT_POSITIVE_DEFINITE){
775:     MatPivotCheck_pd(mat,info,sctx,row);
776:   } else if (info->shifttype == (PetscReal) MAT_SHIFT_INBLOCKS){
777:     MatPivotCheck_inblocks(mat,info,sctx,row);
778:   } else {
779:     MatPivotCheck_none(fact,mat,info,sctx,row);
780:   }
781:   return(0);
782: }

784: /*
785:   Create and initialize a linked list
786:   Input Parameters:
787:     idx_start - starting index of the list
788:     lnk_max   - max value of lnk indicating the end of the list
789:     nlnk      - max length of the list
790:   Output Parameters:
791:     lnk       - list initialized
792:     bt        - PetscBT (bitarray) with all bits set to false
793:     lnk_empty - flg indicating the list is empty
794: */
795: #define PetscLLCreate(idx_start,lnk_max,nlnk,lnk,bt) \
796:   (PetscMalloc1(nlnk,&lnk) || PetscBTCreate(nlnk,&(bt)) || (lnk[idx_start] = lnk_max,0))

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

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

838: /*
839:   Add a permuted index set into a sorted linked list
840:   Input Parameters:
841:     nidx      - number of input indices
842:     indices   - integer array
843:     perm      - permutation of indices
844:     idx_start - starting index of the list
845:     lnk       - linked list(an integer array) that is created
846:     bt        - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
847:   output Parameters:
848:     nlnk      - number of newly added indices
849:     lnk       - the sorted(increasing order) linked list containing new and non-redundate entries from indices
850:     bt        - updated PetscBT (bitarray)
851: */
852: #define PetscLLAddPerm(nidx,indices,perm,idx_start,nlnk,lnk,bt) 0;\
853: {\
854:   PetscInt _k,_entry,_location,_lnkdata;\
855:   nlnk     = 0;\
856:   _lnkdata = idx_start;\
857:   for (_k=0; _k<nidx; _k++){\
858:     _entry = perm[indices[_k]];\
859:     if (!PetscBTLookupSet(bt,_entry)){  /* new entry */\
860:       /* search for insertion location */\
861:       /* start from the beginning if _entry < previous _entry */\
862:       if (_k && _entry < _lnkdata) _lnkdata  = idx_start;\
863:       do {\
864:         _location = _lnkdata;\
865:         _lnkdata  = lnk[_location];\
866:       } while (_entry > _lnkdata);\
867:       /* insertion location is found, add entry into lnk */\
868:       lnk[_location] = _entry;\
869:       lnk[_entry]    = _lnkdata;\
870:       nlnk++;\
871:       _lnkdata = _entry; /* next search starts from here if next_entry > _entry */\
872:     }\
873:   }\
874: }

876: /*
877:   Add a SORTED ascending index set into a sorted linked list - same as PetscLLAdd() bus skip 'if (_k && _entry < _lnkdata) _lnkdata  = idx_start;'
878:   Input Parameters:
879:     nidx      - number of input indices
880:     indices   - sorted integer array
881:     idx_start - starting index of the list
882:     lnk       - linked list(an integer array) that is created
883:     bt        - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
884:   output Parameters:
885:     nlnk      - number of newly added indices
886:     lnk       - the sorted(increasing order) linked list containing new and non-redundate entries from indices
887:     bt        - updated PetscBT (bitarray)
888: */
889: #define PetscLLAddSorted(nidx,indices,idx_start,nlnk,lnk,bt) 0;\
890: {\
891:   PetscInt _k,_entry,_location,_lnkdata;\
892:   nlnk      = 0;\
893:   _lnkdata  = idx_start;\
894:   for (_k=0; _k<nidx; _k++){\
895:     _entry = indices[_k];\
896:     if (!PetscBTLookupSet(bt,_entry)){  /* new entry */\
897:       /* search for insertion location */\
898:       do {\
899:         _location = _lnkdata;\
900:         _lnkdata  = lnk[_location];\
901:       } while (_entry > _lnkdata);\
902:       /* insertion location is found, add entry into lnk */\
903:       lnk[_location] = _entry;\
904:       lnk[_entry]    = _lnkdata;\
905:       nlnk++;\
906:       _lnkdata = _entry; /* next search starts from here */\
907:     }\
908:   }\
909: }

911: #define PetscLLAddSorted_new(nidx,indices,idx_start,lnk_empty,nlnk,lnk,bt) 0; \
912: {\
913:   PetscInt _k,_entry,_location,_lnkdata;\
914:   if (lnk_empty){\
915:     _lnkdata  = idx_start;                      \
916:     for (_k=0; _k<nidx; _k++){                  \
917:       _entry = indices[_k];                             \
918:       PetscBTSet(bt,_entry);  /* mark the new entry */          \
919:           _location = _lnkdata;                                 \
920:           _lnkdata  = lnk[_location];                           \
921:         /* insertion location is found, add entry into lnk */   \
922:         lnk[_location] = _entry;                                \
923:         lnk[_entry]    = _lnkdata;                              \
924:         _lnkdata = _entry; /* next search starts from here */   \
925:     }                                                           \
926:     /*\
927:     lnk[indices[nidx-1]] = lnk[idx_start];\
928:     lnk[idx_start]       = indices[0];\
929:     PetscBTSet(bt,indices[0]);  \
930:     for (_k=1; _k<nidx; _k++){                  \
931:       PetscBTSet(bt,indices[_k]);                                          \
932:       lnk[indices[_k-1]] = indices[_k];                                  \
933:     }                                                           \
934:      */\
935:     nlnk      = nidx;\
936:     lnk_empty = PETSC_FALSE;\
937:   } else {\
938:     nlnk      = 0;                              \
939:     _lnkdata  = idx_start;                      \
940:     for (_k=0; _k<nidx; _k++){                  \
941:       _entry = indices[_k];                             \
942:       if (!PetscBTLookupSet(bt,_entry)){  /* new entry */       \
943:         /* search for insertion location */                     \
944:         do {                                                    \
945:           _location = _lnkdata;                                 \
946:           _lnkdata  = lnk[_location];                           \
947:         } while (_entry > _lnkdata);                            \
948:         /* insertion location is found, add entry into lnk */   \
949:         lnk[_location] = _entry;                                \
950:         lnk[_entry]    = _lnkdata;                              \
951:         nlnk++;                                                 \
952:         _lnkdata = _entry; /* next search starts from here */   \
953:       }                                                         \
954:     }                                                           \
955:   }                                                             \
956: }

958: /*
959:   Add a SORTED index set into a sorted linked list used for LUFactorSymbolic()
960:   Same as PetscLLAddSorted() with an additional operation:
961:        count the number of input indices that are no larger than 'diag'
962:   Input Parameters:
963:     indices   - sorted integer array
964:     idx_start - starting index of the list, index of pivot row
965:     lnk       - linked list(an integer array) that is created
966:     bt        - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
967:     diag      - index of the active row in LUFactorSymbolic
968:     nzbd      - number of input indices with indices <= idx_start
969:     im        - im[idx_start] is initialized as num of nonzero entries in row=idx_start
970:   output Parameters:
971:     nlnk      - number of newly added indices
972:     lnk       - the sorted(increasing order) linked list containing new and non-redundate entries from indices
973:     bt        - updated PetscBT (bitarray)
974:     im        - im[idx_start]: unchanged if diag is not an entry
975:                              : num of entries with indices <= diag if diag is an entry
976: */
977: #define PetscLLAddSortedLU(indices,idx_start,nlnk,lnk,bt,diag,nzbd,im) 0;\
978: {\
979:   PetscInt _k,_entry,_location,_lnkdata,_nidx;\
980:   nlnk     = 0;\
981:   _lnkdata = idx_start;\
982:   _nidx = im[idx_start] - nzbd; /* num of entries with idx_start < index <= diag */\
983:   for (_k=0; _k<_nidx; _k++){\
984:     _entry = indices[_k];\
985:     nzbd++;\
986:     if ( _entry== diag) im[idx_start] = nzbd;\
987:     if (!PetscBTLookupSet(bt,_entry)){  /* new entry */\
988:       /* search for insertion location */\
989:       do {\
990:         _location = _lnkdata;\
991:         _lnkdata  = lnk[_location];\
992:       } while (_entry > _lnkdata);\
993:       /* insertion location is found, add entry into lnk */\
994:       lnk[_location] = _entry;\
995:       lnk[_entry]    = _lnkdata;\
996:       nlnk++;\
997:       _lnkdata = _entry; /* next search starts from here */\
998:     }\
999:   }\
1000: }

1002: /*
1003:   Copy data on the list into an array, then initialize the list
1004:   Input Parameters:
1005:     idx_start - starting index of the list
1006:     lnk_max   - max value of lnk indicating the end of the list
1007:     nlnk      - number of data on the list to be copied
1008:     lnk       - linked list
1009:     bt        - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
1010:   output Parameters:
1011:     indices   - array that contains the copied data
1012:     lnk       - linked list that is cleaned and initialize
1013:     bt        - PetscBT (bitarray) with all bits set to false
1014: */
1015: #define PetscLLClean(idx_start,lnk_max,nlnk,lnk,indices,bt) 0;\
1016: {\
1017:   PetscInt _j,_idx=idx_start;\
1018:   for (_j=0; _j<nlnk; _j++){\
1019:     _idx = lnk[_idx];\
1020:     indices[_j] = _idx;\
1021:     PetscBTClear(bt,_idx);\
1022:   }\
1023:   lnk[idx_start] = lnk_max;\
1024: }
1025: /*
1026:   Free memories used by the list
1027: */
1028: #define PetscLLDestroy(lnk,bt) (PetscFree(lnk) || PetscBTDestroy(&(bt)))

1030: /* Routines below are used for incomplete matrix factorization */
1031: /*
1032:   Create and initialize a linked list and its levels
1033:   Input Parameters:
1034:     idx_start - starting index of the list
1035:     lnk_max   - max value of lnk indicating the end of the list
1036:     nlnk      - max length of the list
1037:   Output Parameters:
1038:     lnk       - list initialized
1039:     lnk_lvl   - array of size nlnk for storing levels of lnk
1040:     bt        - PetscBT (bitarray) with all bits set to false
1041: */
1042: #define PetscIncompleteLLCreate(idx_start,lnk_max,nlnk,lnk,lnk_lvl,bt)\
1043:   (PetscIntMultError(2,nlnk,NULL) || PetscMalloc1(2*nlnk,&lnk) || PetscBTCreate(nlnk,&(bt)) || (lnk[idx_start] = lnk_max,lnk_lvl = lnk + nlnk,0))

1045: /*
1046:   Initialize a sorted linked list used for ILU and ICC
1047:   Input Parameters:
1048:     nidx      - number of input idx
1049:     idx       - integer array used for storing column indices
1050:     idx_start - starting index of the list
1051:     perm      - indices of an IS
1052:     lnk       - linked list(an integer array) that is created
1053:     lnklvl    - levels of lnk
1054:     bt        - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
1055:   output Parameters:
1056:     nlnk     - number of newly added idx
1057:     lnk      - the sorted(increasing order) linked list containing new and non-redundate entries from idx
1058:     lnklvl   - levels of lnk
1059:     bt       - updated PetscBT (bitarray)
1060: */
1061: #define PetscIncompleteLLInit(nidx,idx,idx_start,perm,nlnk,lnk,lnklvl,bt) 0;\
1062: {\
1063:   PetscInt _k,_entry,_location,_lnkdata;\
1064:   nlnk     = 0;\
1065:   _lnkdata = idx_start;\
1066:   for (_k=0; _k<nidx; _k++){\
1067:     _entry = perm[idx[_k]];\
1068:     if (!PetscBTLookupSet(bt,_entry)){  /* new entry */\
1069:       /* search for insertion location */\
1070:       if (_k && _entry < _lnkdata) _lnkdata  = idx_start;\
1071:       do {\
1072:         _location = _lnkdata;\
1073:         _lnkdata  = lnk[_location];\
1074:       } while (_entry > _lnkdata);\
1075:       /* insertion location is found, add entry into lnk */\
1076:       lnk[_location]  = _entry;\
1077:       lnk[_entry]     = _lnkdata;\
1078:       lnklvl[_entry] = 0;\
1079:       nlnk++;\
1080:       _lnkdata = _entry; /* next search starts from here if next_entry > _entry */\
1081:     }\
1082:   }\
1083: }

1085: /*
1086:   Add a SORTED index set into a sorted linked list for ILU
1087:   Input Parameters:
1088:     nidx      - number of input indices
1089:     idx       - sorted integer array used for storing column indices
1090:     level     - level of fill, e.g., ICC(level)
1091:     idxlvl    - level of idx
1092:     idx_start - starting index of the list
1093:     lnk       - linked list(an integer array) that is created
1094:     lnklvl    - levels of lnk
1095:     bt        - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
1096:     prow      - the row number of idx
1097:   output Parameters:
1098:     nlnk     - number of newly added idx
1099:     lnk      - the sorted(increasing order) linked list containing new and non-redundate entries from idx
1100:     lnklvl   - levels of lnk
1101:     bt       - updated PetscBT (bitarray)

1103:   Note: the level of factor(i,j) is set as lvl(i,j) = min{ lvl(i,j), lvl(i,prow)+lvl(prow,j)+1)
1104:         where idx = non-zero columns of U(prow,prow+1:n-1), prow<i
1105: */
1106: #define PetscILULLAddSorted(nidx,idx,level,idxlvl,idx_start,nlnk,lnk,lnklvl,bt,lnklvl_prow) 0;\
1107: {\
1108:   PetscInt _k,_entry,_location,_lnkdata,_incrlev,_lnklvl_prow=lnklvl[prow];\
1109:   nlnk     = 0;\
1110:   _lnkdata = idx_start;\
1111:   for (_k=0; _k<nidx; _k++){\
1112:     _incrlev = idxlvl[_k] + _lnklvl_prow + 1;\
1113:     if (_incrlev > level) continue;\
1114:     _entry = idx[_k];\
1115:     if (!PetscBTLookupSet(bt,_entry)){  /* new entry */\
1116:       /* search for insertion location */\
1117:       do {\
1118:         _location = _lnkdata;\
1119:         _lnkdata  = lnk[_location];\
1120:       } while (_entry > _lnkdata);\
1121:       /* insertion location is found, add entry into lnk */\
1122:       lnk[_location]  = _entry;\
1123:       lnk[_entry]     = _lnkdata;\
1124:       lnklvl[_entry] = _incrlev;\
1125:       nlnk++;\
1126:       _lnkdata = _entry; /* next search starts from here if next_entry > _entry */\
1127:     } else { /* existing entry: update lnklvl */\
1128:       if (lnklvl[_entry] > _incrlev) lnklvl[_entry] = _incrlev;\
1129:     }\
1130:   }\
1131: }

1133: /*
1134:   Add a index set into a sorted linked list
1135:   Input Parameters:
1136:     nidx      - number of input idx
1137:     idx   - integer array used for storing column indices
1138:     level     - level of fill, e.g., ICC(level)
1139:     idxlvl - level of idx
1140:     idx_start - starting index of the list
1141:     lnk       - linked list(an integer array) that is created
1142:     lnklvl   - levels of lnk
1143:     bt        - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
1144:   output Parameters:
1145:     nlnk      - number of newly added idx
1146:     lnk       - the sorted(increasing order) linked list containing new and non-redundate entries from idx
1147:     lnklvl   - levels of lnk
1148:     bt        - updated PetscBT (bitarray)
1149: */
1150: #define PetscIncompleteLLAdd(nidx,idx,level,idxlvl,idx_start,nlnk,lnk,lnklvl,bt) 0;\
1151: {\
1152:   PetscInt _k,_entry,_location,_lnkdata,_incrlev;\
1153:   nlnk     = 0;\
1154:   _lnkdata = idx_start;\
1155:   for (_k=0; _k<nidx; _k++){\
1156:     _incrlev = idxlvl[_k] + 1;\
1157:     if (_incrlev > level) continue;\
1158:     _entry = idx[_k];\
1159:     if (!PetscBTLookupSet(bt,_entry)){  /* new entry */\
1160:       /* search for insertion location */\
1161:       if (_k && _entry < _lnkdata) _lnkdata  = idx_start;\
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 if next_entry > _entry */\
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
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:   output Parameters:
1190:     nlnk      - number of newly added idx
1191:     lnk       - the sorted(increasing order) linked list containing new and non-redundate entries from idx
1192:     lnklvl    - levels of lnk
1193:     bt        - updated PetscBT (bitarray)
1194: */
1195: #define PetscIncompleteLLAddSorted(nidx,idx,level,idxlvl,idx_start,nlnk,lnk,lnklvl,bt) 0;\
1196: {\
1197:   PetscInt _k,_entry,_location,_lnkdata,_incrlev;\
1198:   nlnk = 0;\
1199:   _lnkdata = idx_start;\
1200:   for (_k=0; _k<nidx; _k++){\
1201:     _incrlev = idxlvl[_k] + 1;\
1202:     if (_incrlev > level) continue;\
1203:     _entry = idx[_k];\
1204:     if (!PetscBTLookupSet(bt,_entry)){  /* new entry */\
1205:       /* search for insertion location */\
1206:       do {\
1207:         _location = _lnkdata;\
1208:         _lnkdata  = lnk[_location];\
1209:       } while (_entry > _lnkdata);\
1210:       /* insertion location is found, add entry into lnk */\
1211:       lnk[_location] = _entry;\
1212:       lnk[_entry]    = _lnkdata;\
1213:       lnklvl[_entry] = _incrlev;\
1214:       nlnk++;\
1215:       _lnkdata = _entry; /* next search starts from here */\
1216:     } else { /* existing entry: update lnklvl */\
1217:       if (lnklvl[_entry] > _incrlev) lnklvl[_entry] = _incrlev;\
1218:     }\
1219:   }\
1220: }

1222: /*
1223:   Add a SORTED index set into a sorted linked list for ICC
1224:   Input Parameters:
1225:     nidx      - number of input indices
1226:     idx       - sorted integer array used for storing column indices
1227:     level     - level of fill, e.g., ICC(level)
1228:     idxlvl    - level of idx
1229:     idx_start - starting index of the list
1230:     lnk       - linked list(an integer array) that is created
1231:     lnklvl    - levels of lnk
1232:     bt        - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
1233:     idxlvl_prow - idxlvl[prow], where prow is the row number of the idx
1234:   output Parameters:
1235:     nlnk   - number of newly added indices
1236:     lnk    - the sorted(increasing order) linked list containing new and non-redundate entries from idx
1237:     lnklvl - levels of lnk
1238:     bt     - updated PetscBT (bitarray)
1239:   Note: the level of U(i,j) is set as lvl(i,j) = min{ lvl(i,j), lvl(prow,i)+lvl(prow,j)+1)
1240:         where idx = non-zero columns of U(prow,prow+1:n-1), prow<i
1241: */
1242: #define PetscICCLLAddSorted(nidx,idx,level,idxlvl,idx_start,nlnk,lnk,lnklvl,bt,idxlvl_prow) 0;\
1243: {\
1244:   PetscInt _k,_entry,_location,_lnkdata,_incrlev;\
1245:   nlnk = 0;\
1246:   _lnkdata = idx_start;\
1247:   for (_k=0; _k<nidx; _k++){\
1248:     _incrlev = idxlvl[_k] + idxlvl_prow + 1;\
1249:     if (_incrlev > level) continue;\
1250:     _entry = idx[_k];\
1251:     if (!PetscBTLookupSet(bt,_entry)){  /* new entry */\
1252:       /* search for insertion location */\
1253:       do {\
1254:         _location = _lnkdata;\
1255:         _lnkdata  = lnk[_location];\
1256:       } while (_entry > _lnkdata);\
1257:       /* insertion location is found, add entry into lnk */\
1258:       lnk[_location] = _entry;\
1259:       lnk[_entry]    = _lnkdata;\
1260:       lnklvl[_entry] = _incrlev;\
1261:       nlnk++;\
1262:       _lnkdata = _entry; /* next search starts from here */\
1263:     } else { /* existing entry: update lnklvl */\
1264:       if (lnklvl[_entry] > _incrlev) lnklvl[_entry] = _incrlev;\
1265:     }\
1266:   }\
1267: }

1269: /*
1270:   Copy data on the list into an array, then initialize the list
1271:   Input Parameters:
1272:     idx_start - starting index of the list
1273:     lnk_max   - max value of lnk indicating the end of the list
1274:     nlnk      - number of data on the list to be copied
1275:     lnk       - linked list
1276:     lnklvl    - level of lnk
1277:     bt        - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
1278:   output Parameters:
1279:     indices - array that contains the copied data
1280:     lnk     - linked list that is cleaned and initialize
1281:     lnklvl  - level of lnk that is reinitialized
1282:     bt      - PetscBT (bitarray) with all bits set to false
1283: */
1284: #define PetscIncompleteLLClean(idx_start,lnk_max,nlnk,lnk,lnklvl,indices,indiceslvl,bt) 0;\
1285: do {\
1286:   PetscInt _j,_idx=idx_start;\
1287:   for (_j=0; _j<nlnk; _j++){\
1288:     _idx = lnk[_idx];\
1289:     *(indices+_j) = _idx;\
1290:     *(indiceslvl+_j) = lnklvl[_idx];\
1291:     lnklvl[_idx] = -1;\
1292:     PetscBTClear(bt,_idx);\
1293:   }\
1294:   lnk[idx_start] = lnk_max;\
1295: } while (0)
1296: /*
1297:   Free memories used by the list
1298: */
1299: #define PetscIncompleteLLDestroy(lnk,bt) (PetscFree(lnk) || PetscBTDestroy(&(bt)))

1301: #define MatCheckSameLocalSize(A,ar1,B,ar2) do { \
1303:   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);} while (0)

1305: #define MatCheckSameSize(A,ar1,B,ar2) do { \
1306:   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);\
1307:   MatCheckSameLocalSize(A,ar1,B,ar2);} while (0)

1309: #define VecCheckMatCompatible(M,x,ar1,b,ar2) do { \
1310:   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); \
1311:   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);} while (0)

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

1322:       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
1323:       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
1324:       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
1325:       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.
1326:       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
1327:       to each other so memory access is much better than using the big array.

1329:   Example:
1330:      nlnk_max=5, lnk_max=36:
1331:      Initial list: [0, 0 | 36, 2 | 0, 0 | 0, 0 | 0, 0 | 0, 0 | 0, 0]
1332:      here, head_node has index 2 with value lnk[2]=lnk_max=36,
1333:            0-th entry is used to store the number of entries in the list,
1334:      The initial lnk represents head -> tail(marked by 36) with number of entries = lnk[0]=0.

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

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

1344:   Input Parameters:
1345:     nlnk_max  - max length of the list
1346:     lnk_max   - max value of the entries
1347:   Output Parameters:
1348:     lnk       - list created and initialized
1349:     bt        - PetscBT (bitarray) with all bits set to false. Note: bt has size lnk_max, not nln_max!
1350: */
1351: PETSC_STATIC_INLINE PetscErrorCode PetscLLCondensedCreate(PetscInt nlnk_max,PetscInt lnk_max,PetscInt **lnk,PetscBT *bt)
1352: {
1354:   PetscInt       *llnk,lsize = 0;

1357:   PetscIntMultError(2,nlnk_max+2,&lsize);
1358:   PetscMalloc1(lsize,lnk);
1359:   PetscBTCreate(lnk_max,bt);
1360:   llnk = *lnk;
1361:   llnk[0] = 0;         /* number of entries on the list */
1362:   llnk[2] = lnk_max;   /* value in the head node */
1363:   llnk[3] = 2;         /* next for the head node */
1364:   return(0);
1365: }

1367: /*
1368:   Add a SORTED ascending index set into a sorted linked list. See PetscLLCondensedCreate() for detailed description.
1369:   Input Parameters:
1370:     nidx      - number of input indices
1371:     indices   - sorted integer array
1372:     lnk       - condensed linked list(an integer array) that is created
1373:     bt        - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
1374:   output Parameters:
1375:     lnk       - the sorted(increasing order) linked list containing previous and newly added non-redundate indices
1376:     bt        - updated PetscBT (bitarray)
1377: */
1378: PETSC_STATIC_INLINE PetscErrorCode PetscLLCondensedAddSorted(PetscInt nidx,const PetscInt indices[],PetscInt lnk[],PetscBT bt)
1379: {
1380:   PetscInt _k,_entry,_location,_next,_lnkdata,_nlnk,_newnode;

1383:   _nlnk     = lnk[0]; /* num of entries on the input lnk */
1384:   _location = 2; /* head */
1385:     for (_k=0; _k<nidx; _k++){
1386:       _entry = indices[_k];
1387:       if (!PetscBTLookupSet(bt,_entry)){  /* new entry */
1388:         /* search for insertion location */
1389:         do {
1390:           _next     = _location + 1; /* link from previous node to next node */
1391:           _location = lnk[_next];    /* idx of next node */
1392:           _lnkdata  = lnk[_location];/* value of next node */
1393:         } while (_entry > _lnkdata);
1394:         /* insertion location is found, add entry into lnk */
1395:         _newnode        = 2*(_nlnk+2);   /* index for this new node */
1396:         lnk[_next]      = _newnode;      /* connect previous node to the new node */
1397:         lnk[_newnode]   = _entry;        /* set value of the new node */
1398:         lnk[_newnode+1] = _location;     /* connect new node to next node */
1399:         _location       = _newnode;      /* next search starts from the new node */
1400:         _nlnk++;
1401:       }   \
1402:     }\
1403:   lnk[0]   = _nlnk;   /* number of entries in the list */
1404:   return(0);
1405: }

1407: PETSC_STATIC_INLINE PetscErrorCode PetscLLCondensedClean(PetscInt lnk_max,PetscInt nidx,PetscInt *indices,PetscInt lnk[],PetscBT bt)
1408: {
1410:   PetscInt       _k,_next,_nlnk;

1413:   _next = lnk[3];       /* head node */
1414:   _nlnk = lnk[0];       /* num of entries on the list */
1415:   for (_k=0; _k<_nlnk; _k++){
1416:     indices[_k] = lnk[_next];
1417:     _next       = lnk[_next + 1];
1418:     PetscBTClear(bt,indices[_k]);
1419:   }
1420:   lnk[0] = 0;          /* num of entries on the list */
1421:   lnk[2] = lnk_max;    /* initialize head node */
1422:   lnk[3] = 2;          /* head node */
1423:   return(0);
1424: }

1426: PETSC_STATIC_INLINE PetscErrorCode PetscLLCondensedView(PetscInt *lnk)
1427: {
1429:   PetscInt       k;

1432:   PetscPrintf(PETSC_COMM_SELF,"LLCondensed of size %D, (val,  next)\n",lnk[0]);
1433:   for (k=2; k< lnk[0]+2; k++){
1434:     PetscPrintf(PETSC_COMM_SELF," %D: (%D, %D)\n",2*k,lnk[2*k],lnk[2*k+1]);
1435:   }
1436:   return(0);
1437: }

1439: /*
1440:   Free memories used by the list
1441: */
1442: PETSC_STATIC_INLINE PetscErrorCode PetscLLCondensedDestroy(PetscInt *lnk,PetscBT bt)
1443: {

1447:   PetscFree(lnk);
1448:   PetscBTDestroy(&bt);
1449:   return(0);
1450: }

1452: /* -------------------------------------------------------------------------------------------------------*/
1453: /*
1454:  Same as PetscLLCondensedCreate(), but does not use non-scalable O(lnk_max) bitarray
1455:   Input Parameters:
1456:     nlnk_max  - max length of the list
1457:   Output Parameters:
1458:     lnk       - list created and initialized
1459: */
1460: PETSC_STATIC_INLINE PetscErrorCode PetscLLCondensedCreate_Scalable(PetscInt nlnk_max,PetscInt **lnk)
1461: {
1463:   PetscInt       *llnk,lsize = 0;

1466:   PetscIntMultError(2,nlnk_max+2,&lsize);
1467:   PetscMalloc1(lsize,lnk);
1468:   llnk = *lnk;
1469:   llnk[0] = 0;               /* number of entries on the list */
1470:   llnk[2] = PETSC_MAX_INT;   /* value in the head node */
1471:   llnk[3] = 2;               /* next for the head node */
1472:   return(0);
1473: }

1475: PETSC_STATIC_INLINE PetscErrorCode PetscLLCondensedExpand_Scalable(PetscInt nlnk_max,PetscInt **lnk)
1476: {
1478:   PetscInt       lsize = 0;

1481:   PetscIntMultError(2,nlnk_max+2,&lsize);
1482:   PetscRealloc(lsize*sizeof(PetscInt),lnk);
1483:   return(0);
1484: }

1486: PETSC_STATIC_INLINE PetscErrorCode PetscLLCondensedAddSorted_Scalable(PetscInt nidx,const PetscInt indices[],PetscInt lnk[])
1487: {
1488:   PetscInt _k,_entry,_location,_next,_lnkdata,_nlnk,_newnode;
1489:   _nlnk     = lnk[0]; /* num of entries on the input lnk */
1490:   _location = 2; /* head */ \
1491:     for (_k=0; _k<nidx; _k++){
1492:       _entry = indices[_k];
1493:       /* search for insertion location */
1494:       do {
1495:         _next     = _location + 1; /* link from previous node to next node */
1496:         _location = lnk[_next];    /* idx of next node */
1497:         _lnkdata  = lnk[_location];/* value of next node */
1498:       } while (_entry > _lnkdata);
1499:       if (_entry < _lnkdata) {
1500:         /* insertion location is found, add entry into lnk */
1501:         _newnode        = 2*(_nlnk+2);   /* index for this new node */
1502:         lnk[_next]      = _newnode;      /* connect previous node to the new node */
1503:         lnk[_newnode]   = _entry;        /* set value of the new node */
1504:         lnk[_newnode+1] = _location;     /* connect new node to next node */
1505:         _location       = _newnode;      /* next search starts from the new node */
1506:         _nlnk++;
1507:       }
1508:     }
1509:   lnk[0]   = _nlnk;   /* number of entries in the list */
1510:   return 0;
1511: }

1513: PETSC_STATIC_INLINE PetscErrorCode PetscLLCondensedClean_Scalable(PetscInt nidx,PetscInt *indices,PetscInt *lnk)
1514: {
1515:   PetscInt _k,_next,_nlnk;
1516:   _next = lnk[3];       /* head node */
1517:   _nlnk = lnk[0];
1518:   for (_k=0; _k<_nlnk; _k++){
1519:     indices[_k] = lnk[_next];
1520:     _next       = lnk[_next + 1];
1521:   }
1522:   lnk[0] = 0;          /* num of entries on the list */
1523:   lnk[3] = 2;          /* head node */
1524:   return 0;
1525: }

1527: PETSC_STATIC_INLINE PetscErrorCode PetscLLCondensedDestroy_Scalable(PetscInt *lnk)
1528: {
1529:   return PetscFree(lnk);
1530: }

1532: /* -------------------------------------------------------------------------------------------------------*/
1533: /*
1534:       lnk[0]   number of links
1535:       lnk[1]   number of entries
1536:       lnk[3n]  value
1537:       lnk[3n+1] len
1538:       lnk[3n+2] link to next value

1540:       The next three are always the first link

1542:       lnk[3]    PETSC_MIN_INT+1
1543:       lnk[4]    1
1544:       lnk[5]    link to first real entry

1546:       The next three are always the last link

1548:       lnk[6]    PETSC_MAX_INT - 1
1549:       lnk[7]    1
1550:       lnk[8]    next valid link (this is the same as lnk[0] but without the decreases)
1551: */

1553: PETSC_STATIC_INLINE PetscErrorCode PetscLLCondensedCreate_fast(PetscInt nlnk_max,PetscInt **lnk)
1554: {
1556:   PetscInt       *llnk,lsize = 0;

1559:   PetscIntMultError(3,nlnk_max+3,&lsize);
1560:   PetscMalloc1(lsize,lnk);
1561:   llnk = *lnk;
1562:   llnk[0] = 0;   /* nlnk: number of entries on the list */
1563:   llnk[1] = 0;          /* number of integer entries represented in list */
1564:   llnk[3] = PETSC_MIN_INT+1;   /* value in the first node */
1565:   llnk[4] = 1;           /* count for the first node */
1566:   llnk[5] = 6;         /* next for the first node */
1567:   llnk[6] = PETSC_MAX_INT-1;   /* value in the last node */
1568:   llnk[7] = 1;           /* count for the last node */
1569:   llnk[8] = 0;         /* next valid node to be used */
1570:   return(0);
1571: }

1573: PETSC_STATIC_INLINE PetscErrorCode PetscLLCondensedAddSorted_fast(PetscInt nidx,const PetscInt indices[],PetscInt lnk[])
1574: {
1575:   PetscInt k,entry,prev,next;
1576:   prev      = 3;      /* first value */
1577:   next      = lnk[prev+2];
1578:   for (k=0; k<nidx; k++){
1579:     entry = indices[k];
1580:     /* search for insertion location */
1581:     while (entry >= lnk[next]) {
1582:       prev = next;
1583:       next = lnk[next+2];
1584:     }
1585:     /* entry is in range of previous list */
1586:     if (entry < lnk[prev]+lnk[prev+1]) continue;
1587:     lnk[1]++;
1588:     /* entry is right after previous list */
1589:     if (entry == lnk[prev]+lnk[prev+1]) {
1590:       lnk[prev+1]++;
1591:       if (lnk[next] == entry+1) { /* combine two contiguous strings */
1592:         lnk[prev+1] += lnk[next+1];
1593:         lnk[prev+2]  = lnk[next+2];
1594:         next         = lnk[next+2];
1595:         lnk[0]--;
1596:       }
1597:       continue;
1598:     }
1599:     /* entry is right before next list */
1600:     if (entry == lnk[next]-1) {
1601:       lnk[next]--;
1602:       lnk[next+1]++;
1603:       prev = next;
1604:       next = lnk[prev+2];
1605:       continue;
1606:     }
1607:     /*  add entry into lnk */
1608:     lnk[prev+2]    = 3*((lnk[8]++)+3);      /* connect previous node to the new node */
1609:     prev           = lnk[prev+2];
1610:     lnk[prev]      = entry;        /* set value of the new node */
1611:     lnk[prev+1]    = 1;             /* number of values in contiguous string is one to start */
1612:     lnk[prev+2]    = next;          /* connect new node to next node */
1613:     lnk[0]++;
1614:   }
1615:   return 0;
1616: }

1618: PETSC_STATIC_INLINE PetscErrorCode PetscLLCondensedClean_fast(PetscInt nidx,PetscInt *indices,PetscInt *lnk)
1619: {
1620:   PetscInt _k,_next,_nlnk,cnt,j;
1621:   _next = lnk[5];       /* first node */
1622:   _nlnk = lnk[0];
1623:   cnt   = 0;
1624:   for (_k=0; _k<_nlnk; _k++){
1625:     for (j=0; j<lnk[_next+1]; j++) {
1626:       indices[cnt++] = lnk[_next] + j;
1627:     }
1628:     _next       = lnk[_next + 2];
1629:   }
1630:   lnk[0] = 0;   /* nlnk: number of links */
1631:   lnk[1] = 0;          /* number of integer entries represented in list */
1632:   lnk[3] = PETSC_MIN_INT+1;   /* value in the first node */
1633:   lnk[4] = 1;           /* count for the first node */
1634:   lnk[5] = 6;         /* next for the first node */
1635:   lnk[6] = PETSC_MAX_INT-1;   /* value in the last node */
1636:   lnk[7] = 1;           /* count for the last node */
1637:   lnk[8] = 0;         /* next valid location to make link */
1638:   return 0;
1639: }

1641: PETSC_STATIC_INLINE PetscErrorCode PetscLLCondensedView_fast(PetscInt *lnk)
1642: {
1643:   PetscInt k,next,nlnk;
1644:   next = lnk[5];       /* first node */
1645:   nlnk = lnk[0];
1646:   for (k=0; k<nlnk; k++){
1647: #if 0                           /* Debugging code */
1648:     printf("%d value %d len %d next %d\n",next,lnk[next],lnk[next+1],lnk[next+2]);
1649: #endif
1650:     next = lnk[next + 2];
1651:   }
1652:   return 0;
1653: }

1655: PETSC_STATIC_INLINE PetscErrorCode PetscLLCondensedDestroy_fast(PetscInt *lnk)
1656: {
1657:   return PetscFree(lnk);
1658: }

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

1663: PETSC_EXTERN PetscLogEvent MAT_Mult;
1664: PETSC_EXTERN PetscLogEvent MAT_MultMatrixFree;
1665: PETSC_EXTERN PetscLogEvent MAT_Mults;
1666: PETSC_EXTERN PetscLogEvent MAT_MultConstrained;
1667: PETSC_EXTERN PetscLogEvent MAT_MultAdd;
1668: PETSC_EXTERN PetscLogEvent MAT_MultTranspose;
1669: PETSC_EXTERN PetscLogEvent MAT_MultTransposeConstrained;
1670: PETSC_EXTERN PetscLogEvent MAT_MultTransposeAdd;
1671: PETSC_EXTERN PetscLogEvent MAT_Solve;
1672: PETSC_EXTERN PetscLogEvent MAT_Solves;
1673: PETSC_EXTERN PetscLogEvent MAT_SolveAdd;
1674: PETSC_EXTERN PetscLogEvent MAT_SolveTranspose;
1675: PETSC_EXTERN PetscLogEvent MAT_SolveTransposeAdd;
1676: PETSC_EXTERN PetscLogEvent MAT_SOR;
1677: PETSC_EXTERN PetscLogEvent MAT_ForwardSolve;
1678: PETSC_EXTERN PetscLogEvent MAT_BackwardSolve;
1679: PETSC_EXTERN PetscLogEvent MAT_LUFactor;
1680: PETSC_EXTERN PetscLogEvent MAT_LUFactorSymbolic;
1681: PETSC_EXTERN PetscLogEvent MAT_LUFactorNumeric;
1682: PETSC_EXTERN PetscLogEvent MAT_CholeskyFactor;
1683: PETSC_EXTERN PetscLogEvent MAT_CholeskyFactorSymbolic;
1684: PETSC_EXTERN PetscLogEvent MAT_CholeskyFactorNumeric;
1685: PETSC_EXTERN PetscLogEvent MAT_ILUFactor;
1686: PETSC_EXTERN PetscLogEvent MAT_ILUFactorSymbolic;
1687: PETSC_EXTERN PetscLogEvent MAT_ICCFactorSymbolic;
1688: PETSC_EXTERN PetscLogEvent MAT_Copy;
1689: PETSC_EXTERN PetscLogEvent MAT_Convert;
1690: PETSC_EXTERN PetscLogEvent MAT_Scale;
1691: PETSC_EXTERN PetscLogEvent MAT_AssemblyBegin;
1692: PETSC_EXTERN PetscLogEvent MAT_AssemblyEnd;
1693: PETSC_EXTERN PetscLogEvent MAT_SetValues;
1694: PETSC_EXTERN PetscLogEvent MAT_GetValues;
1695: PETSC_EXTERN PetscLogEvent MAT_GetRow;
1696: PETSC_EXTERN PetscLogEvent MAT_GetRowIJ;
1697: PETSC_EXTERN PetscLogEvent MAT_CreateSubMats;
1698: PETSC_EXTERN PetscLogEvent MAT_GetColoring;
1699: PETSC_EXTERN PetscLogEvent MAT_GetOrdering;
1700: PETSC_EXTERN PetscLogEvent MAT_RedundantMat;
1701: PETSC_EXTERN PetscLogEvent MAT_IncreaseOverlap;
1702: PETSC_EXTERN PetscLogEvent MAT_Partitioning;
1703: PETSC_EXTERN PetscLogEvent MAT_PartitioningND;
1704: PETSC_EXTERN PetscLogEvent MAT_Coarsen;
1705: PETSC_EXTERN PetscLogEvent MAT_ZeroEntries;
1706: PETSC_EXTERN PetscLogEvent MAT_Load;
1707: PETSC_EXTERN PetscLogEvent MAT_View;
1708: PETSC_EXTERN PetscLogEvent MAT_AXPY;
1709: PETSC_EXTERN PetscLogEvent MAT_FDColoringCreate;
1710: PETSC_EXTERN PetscLogEvent MAT_TransposeColoringCreate;
1711: PETSC_EXTERN PetscLogEvent MAT_FDColoringSetUp;
1712: PETSC_EXTERN PetscLogEvent MAT_FDColoringApply;
1713: PETSC_EXTERN PetscLogEvent MAT_Transpose;
1714: PETSC_EXTERN PetscLogEvent MAT_FDColoringFunction;
1715: PETSC_EXTERN PetscLogEvent MAT_CreateSubMat;
1716: PETSC_EXTERN PetscLogEvent MAT_MatSolve;
1717: PETSC_EXTERN PetscLogEvent MAT_MatTrSolve;
1718: PETSC_EXTERN PetscLogEvent MAT_MatMultSymbolic;
1719: PETSC_EXTERN PetscLogEvent MAT_MatMultNumeric;
1720: PETSC_EXTERN PetscLogEvent MAT_Getlocalmatcondensed;
1721: PETSC_EXTERN PetscLogEvent MAT_GetBrowsOfAcols;
1722: PETSC_EXTERN PetscLogEvent MAT_GetBrowsOfAocols;
1723: PETSC_EXTERN PetscLogEvent MAT_PtAPSymbolic;
1724: PETSC_EXTERN PetscLogEvent MAT_PtAPNumeric;
1725: PETSC_EXTERN PetscLogEvent MAT_Seqstompinum;
1726: PETSC_EXTERN PetscLogEvent MAT_Seqstompisym;
1727: PETSC_EXTERN PetscLogEvent MAT_Seqstompi;
1728: PETSC_EXTERN PetscLogEvent MAT_Getlocalmat;
1729: PETSC_EXTERN PetscLogEvent MAT_RARtSymbolic;
1730: PETSC_EXTERN PetscLogEvent MAT_RARtNumeric;
1731: PETSC_EXTERN PetscLogEvent MAT_MatTransposeMultSymbolic;
1732: PETSC_EXTERN PetscLogEvent MAT_MatTransposeMultNumeric;
1733: PETSC_EXTERN PetscLogEvent MAT_TransposeMatMultSymbolic;
1734: PETSC_EXTERN PetscLogEvent MAT_TransposeMatMultNumeric;
1735: PETSC_EXTERN PetscLogEvent MAT_MatMatMultSymbolic;
1736: PETSC_EXTERN PetscLogEvent MAT_MatMatMultNumeric;
1737: PETSC_EXTERN PetscLogEvent MAT_Applypapt;
1738: PETSC_EXTERN PetscLogEvent MAT_Applypapt_symbolic;
1739: PETSC_EXTERN PetscLogEvent MAT_Applypapt_numeric;
1740: PETSC_EXTERN PetscLogEvent MAT_Getsymtranspose;
1741: PETSC_EXTERN PetscLogEvent MAT_Getsymtransreduced;
1742: PETSC_EXTERN PetscLogEvent MAT_GetSequentialNonzeroStructure;
1743: PETSC_EXTERN PetscLogEvent MATMFFD_Mult;
1744: PETSC_EXTERN PetscLogEvent MAT_GetMultiProcBlock;
1745: PETSC_EXTERN PetscLogEvent MAT_CUSPARSECopyToGPU;
1746: PETSC_EXTERN PetscLogEvent MAT_SetValuesBatch;
1747: PETSC_EXTERN PetscLogEvent MAT_ViennaCLCopyToGPU;
1748: PETSC_EXTERN PetscLogEvent MAT_DenseCopyToGPU;
1749: PETSC_EXTERN PetscLogEvent MAT_DenseCopyFromGPU;
1750: PETSC_EXTERN PetscLogEvent MAT_Merge;
1751: PETSC_EXTERN PetscLogEvent MAT_Residual;
1752: PETSC_EXTERN PetscLogEvent MAT_SetRandom;
1753: PETSC_EXTERN PetscLogEvent MAT_FactorFactS;
1754: PETSC_EXTERN PetscLogEvent MAT_FactorInvS;
1755: PETSC_EXTERN PetscLogEvent MATCOLORING_Apply;
1756: PETSC_EXTERN PetscLogEvent MATCOLORING_Comm;
1757: PETSC_EXTERN PetscLogEvent MATCOLORING_Local;
1758: PETSC_EXTERN PetscLogEvent MATCOLORING_ISCreate;
1759: PETSC_EXTERN PetscLogEvent MATCOLORING_SetUp;
1760: PETSC_EXTERN PetscLogEvent MATCOLORING_Weights;

1762: #endif