Actual source code: matimpl.h

  1: #pragma once

  3: #include <petscmat.h>
  4: #include <petscmatcoarsen.h>
  5: #include <petsc/private/petscimpl.h>

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

 20: /* Gets the root type of the input matrix's type (e.g., MATAIJ for MATSEQAIJ) */
 21: PETSC_EXTERN PetscErrorCode MatGetRootType_Private(Mat, MatType *);

 23: /* Gets the MPI type corresponding to the input matrix's type (e.g., MATMPIAIJ for MATSEQAIJ) */
 24: PETSC_EXTERN PetscErrorCode MatGetMPIMatType_Private(Mat, MatType *);

 26: /*
 27:   This file defines the parts of the matrix data structure that are
 28:   shared by all matrix types.
 29: */

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

226: #include <petscsys.h>

228: typedef struct _p_MatRootName *MatRootName;
229: struct _p_MatRootName {
230:   char       *rname, *sname, *mname;
231:   MatRootName next;
232: };

234: PETSC_EXTERN MatRootName MatRootNameList;

236: /*
237:    Utility private matrix routines
238: */
239: PETSC_INTERN PetscErrorCode MatFindNonzeroRowsOrCols_Basic(Mat, PetscBool, PetscReal, IS *);
240: PETSC_INTERN PetscErrorCode MatConvert_Basic(Mat, MatType, MatReuse, Mat *);
241: PETSC_INTERN PetscErrorCode MatConvert_Shell(Mat, MatType, MatReuse, Mat *);
242: PETSC_INTERN PetscErrorCode MatConvertFrom_Shell(Mat, MatType, MatReuse, Mat *);
243: PETSC_INTERN PetscErrorCode MatCopy_Basic(Mat, Mat, MatStructure);
244: PETSC_INTERN PetscErrorCode MatDiagonalSet_Default(Mat, Vec, InsertMode);
245: #if defined(PETSC_HAVE_SCALAPACK)
246: PETSC_INTERN PetscErrorCode MatConvert_Dense_ScaLAPACK(Mat, MatType, MatReuse, Mat *);
247: #endif
248: PETSC_INTERN PetscErrorCode MatSetPreallocationCOO_Basic(Mat, PetscCount, PetscInt[], PetscInt[]);
249: PETSC_INTERN PetscErrorCode MatSetValuesCOO_Basic(Mat, const PetscScalar[], InsertMode);

251: /* these callbacks rely on the old matrix function pointers for
252:    matmat operations. They are unsafe, and should be removed.
253:    However, the amount of work needed to clean up all the
254:    implementations is not negligible */
255: PETSC_INTERN PetscErrorCode MatProductSymbolic_AB(Mat);
256: PETSC_INTERN PetscErrorCode MatProductNumeric_AB(Mat);
257: PETSC_INTERN PetscErrorCode MatProductSymbolic_AtB(Mat);
258: PETSC_INTERN PetscErrorCode MatProductNumeric_AtB(Mat);
259: PETSC_INTERN PetscErrorCode MatProductSymbolic_ABt(Mat);
260: PETSC_INTERN PetscErrorCode MatProductNumeric_ABt(Mat);
261: PETSC_INTERN PetscErrorCode MatProductNumeric_PtAP(Mat);
262: PETSC_INTERN PetscErrorCode MatProductNumeric_RARt(Mat);
263: PETSC_INTERN PetscErrorCode MatProductSymbolic_ABC(Mat);
264: PETSC_INTERN PetscErrorCode MatProductNumeric_ABC(Mat);

266: PETSC_INTERN PetscErrorCode MatProductCreate_Private(Mat, Mat, Mat, Mat);
267: /* this callback handles all the different triple products and
268:    does not rely on the function pointers; used by cuSPARSE/hipSPARSE and KOKKOS-KERNELS */
269: PETSC_INTERN PetscErrorCode MatProductSymbolic_ABC_Basic(Mat);

271: /* CreateGraph is common to AIJ seq and mpi */
272: PETSC_INTERN PetscErrorCode MatCreateGraph_Simple_AIJ(Mat, PetscBool, PetscBool, PetscReal, Mat *);

274: #if defined(PETSC_CLANG_STATIC_ANALYZER)
275: template <typename Tm>
276: extern void MatCheckPreallocated(Tm, int);
277: template <typename Tm>
278: extern void MatCheckProduct(Tm, int);
279: #else /* PETSC_CLANG_STATIC_ANALYZER */
280:   #define MatCheckPreallocated(A, arg) \
281:     do { \
282:       if (!(A)->preallocated) PetscCall(MatSetUp(A)); \
283:     } while (0)

285:   #if defined(PETSC_USE_DEBUG)
286:     #define MatCheckProduct(A, arg) \
287:       do { \
288:         PetscCheck((A)->product, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Argument %d \"%s\" is not a matrix obtained from MatProductCreate()", (arg), #A); \
289:       } while (0)
290:   #else
291:     #define MatCheckProduct(A, arg) \
292:       do { \
293:       } while (0)
294:   #endif
295: #endif /* PETSC_CLANG_STATIC_ANALYZER */

297: /*
298:   The stash is used to temporarily store inserted matrix values that
299:   belong to another processor. During the assembly phase the stashed
300:   values are moved to the correct processor and
301: */

303: typedef struct _MatStashSpace *PetscMatStashSpace;

305: struct _MatStashSpace {
306:   PetscMatStashSpace next;
307:   PetscScalar       *space_head, *val;
308:   PetscInt          *idx, *idy;
309:   PetscInt           total_space_size;
310:   PetscInt           local_used;
311:   PetscInt           local_remaining;
312: };

314: PETSC_EXTERN PetscErrorCode PetscMatStashSpaceGet(PetscInt, PetscInt, PetscMatStashSpace *);
315: PETSC_EXTERN PetscErrorCode PetscMatStashSpaceContiguous(PetscInt, PetscMatStashSpace *, PetscScalar *, PetscInt *, PetscInt *);
316: PETSC_EXTERN PetscErrorCode PetscMatStashSpaceDestroy(PetscMatStashSpace *);

318: typedef struct {
319:   PetscInt count;
320: } MatStashHeader;

322: typedef struct {
323:   void    *buffer; /* Of type blocktype, dynamically constructed  */
324:   PetscInt count;
325:   char     pending;
326: } MatStashFrame;

328: typedef struct _MatStash MatStash;
329: struct _MatStash {
330:   PetscInt           nmax;              /* maximum stash size */
331:   PetscInt           umax;              /* user specified max-size */
332:   PetscInt           oldnmax;           /* the nmax value used previously */
333:   PetscInt           n;                 /* stash size */
334:   PetscInt           bs;                /* block size of the stash */
335:   PetscInt           reallocs;          /* preserve the no of mallocs invoked */
336:   PetscMatStashSpace space_head, space; /* linked list to hold stashed global row/column numbers and matrix values */

338:   PetscErrorCode (*ScatterBegin)(Mat, MatStash *, PetscInt *);
339:   PetscErrorCode (*ScatterGetMesg)(MatStash *, PetscMPIInt *, PetscInt **, PetscInt **, PetscScalar **, PetscInt *);
340:   PetscErrorCode (*ScatterEnd)(MatStash *);
341:   PetscErrorCode (*ScatterDestroy)(MatStash *);

343:   /* The following variables are used for communication */
344:   MPI_Comm      comm;
345:   PetscMPIInt   size, rank;
346:   PetscMPIInt   tag1, tag2;
347:   MPI_Request  *send_waits;     /* array of send requests */
348:   MPI_Request  *recv_waits;     /* array of receive requests */
349:   MPI_Status   *send_status;    /* array of send status */
350:   PetscInt      nsends, nrecvs; /* numbers of sends and receives */
351:   PetscScalar  *svalues;        /* sending data */
352:   PetscInt     *sindices;
353:   PetscScalar **rvalues;    /* receiving data (values) */
354:   PetscInt    **rindices;   /* receiving data (indices) */
355:   PetscInt      nprocessed; /* number of messages already processed */
356:   PetscMPIInt  *flg_v;      /* indicates what messages have arrived so far and from whom */
357:   PetscBool     reproduce;
358:   PetscInt      reproduce_count;

360:   /* The following variables are used for BTS communication */
361:   PetscBool       first_assembly_done; /* Is the first time matrix assembly done? */
362:   PetscBool       use_status;          /* Use MPI_Status to determine number of items in each message */
363:   PetscMPIInt     nsendranks;
364:   PetscMPIInt     nrecvranks;
365:   PetscMPIInt    *sendranks;
366:   PetscMPIInt    *recvranks;
367:   MatStashHeader *sendhdr, *recvhdr;
368:   MatStashFrame  *sendframes; /* pointers to the main messages */
369:   MatStashFrame  *recvframes;
370:   MatStashFrame  *recvframe_active;
371:   PetscInt        recvframe_i;     /* index of block within active frame */
372:   PetscMPIInt     recvframe_count; /* Count actually sent for current frame */
373:   PetscInt        recvcount;       /* Number of receives processed so far */
374:   PetscMPIInt    *some_indices;    /* From last call to MPI_Waitsome */
375:   MPI_Status     *some_statuses;   /* Statuses from last call to MPI_Waitsome */
376:   PetscMPIInt     some_count;      /* Number of requests completed in last call to MPI_Waitsome */
377:   PetscMPIInt     some_i;          /* Index of request currently being processed */
378:   MPI_Request    *sendreqs;
379:   MPI_Request    *recvreqs;
380:   PetscSegBuffer  segsendblocks;
381:   PetscSegBuffer  segrecvframe;
382:   PetscSegBuffer  segrecvblocks;
383:   MPI_Datatype    blocktype;
384:   size_t          blocktype_size;
385:   InsertMode     *insertmode; /* Pointer to check mat->insertmode and set upon message arrival in case no local values have been set. */
386: };

388: #if !defined(PETSC_HAVE_MPIUNI)
389: PETSC_INTERN PetscErrorCode MatStashScatterDestroy_BTS(MatStash *);
390: #endif
391: PETSC_INTERN PetscErrorCode MatStashCreate_Private(MPI_Comm, PetscInt, MatStash *);
392: PETSC_INTERN PetscErrorCode MatStashDestroy_Private(MatStash *);
393: PETSC_INTERN PetscErrorCode MatStashScatterEnd_Private(MatStash *);
394: PETSC_INTERN PetscErrorCode MatStashSetInitialSize_Private(MatStash *, PetscInt);
395: PETSC_INTERN PetscErrorCode MatStashGetInfo_Private(MatStash *, PetscInt *, PetscInt *);
396: PETSC_INTERN PetscErrorCode MatStashValuesRow_Private(MatStash *, PetscInt, PetscInt, const PetscInt[], const PetscScalar[], PetscBool);
397: PETSC_INTERN PetscErrorCode MatStashValuesCol_Private(MatStash *, PetscInt, PetscInt, const PetscInt[], const PetscScalar[], PetscInt, PetscBool);
398: PETSC_INTERN PetscErrorCode MatStashValuesRowBlocked_Private(MatStash *, PetscInt, PetscInt, const PetscInt[], const PetscScalar[], PetscInt, PetscInt, PetscInt);
399: PETSC_INTERN PetscErrorCode MatStashValuesColBlocked_Private(MatStash *, PetscInt, PetscInt, const PetscInt[], const PetscScalar[], PetscInt, PetscInt, PetscInt);
400: PETSC_INTERN PetscErrorCode MatStashScatterBegin_Private(Mat, MatStash *, PetscInt *);
401: PETSC_INTERN PetscErrorCode MatStashScatterGetMesg_Private(MatStash *, PetscMPIInt *, PetscInt **, PetscInt **, PetscScalar **, PetscInt *);
402: PETSC_INTERN PetscErrorCode MatGetInfo_External(Mat, MatInfoType, MatInfo *);

404: typedef struct {
405:   PetscInt  dim;
406:   PetscInt  dims[4];
407:   PetscInt  starts[4];
408:   PetscBool noc; /* this is a single component problem, hence user will not set MatStencil.c */
409: } MatStencilInfo;

411: /* Info about using compressed row format */
412: typedef struct {
413:   PetscBool use;    /* indicates compressed rows have been checked and will be used */
414:   PetscInt  nrows;  /* number of non-zero rows */
415:   PetscInt *i;      /* compressed row pointer  */
416:   PetscInt *rindex; /* compressed row index               */
417: } Mat_CompressedRow;
418: PETSC_EXTERN PetscErrorCode MatCheckCompressedRow(Mat, PetscInt, Mat_CompressedRow *, PetscInt *, PetscInt, PetscReal);

420: typedef struct { /* used by MatCreateRedundantMatrix() for reusing matredundant */
421:   PetscInt     nzlocal, nsends, nrecvs;
422:   PetscMPIInt *send_rank, *recv_rank;
423:   PetscInt    *sbuf_nz, *rbuf_nz, *sbuf_j, **rbuf_j;
424:   PetscScalar *sbuf_a, **rbuf_a;
425:   MPI_Comm     subcomm; /* when user does not provide a subcomm */
426:   IS           isrow, iscol;
427:   Mat         *matseq;
428: } Mat_Redundant;

430: typedef struct { /* used by MatProduct() */
431:   MatProductType type;
432:   char          *alg;
433:   Mat            A, B, C, Dwork;
434:   PetscBool      symbolic_used_the_fact_A_is_symmetric; /* Symbolic phase took advantage of the fact that A is symmetric, and optimized e.g. AtB as AB. Then, .. */
435:   PetscBool      symbolic_used_the_fact_B_is_symmetric; /* .. in the numeric phase, if a new A is not symmetric (but has the same sparsity as the old A therefore .. */
436:   PetscBool      symbolic_used_the_fact_C_is_symmetric; /* MatMatMult(A,B,MAT_REUSE_MATRIX,..&C) is still legitimate), we need to redo symbolic! */
437:   PetscReal      fill;
438:   PetscBool      api_user; /* used to distinguish command line options and to indicate the matrix values are ready to be consumed at symbolic phase if needed */

440:   /* Some products may display the information on the algorithm used */
441:   PetscErrorCode (*view)(Mat, PetscViewer);

443:   /* many products have intermediate data structures, each specific to Mat types and product type */
444:   PetscBool clear;                   /* whether or not to clear the data structures after MatProductNumeric has been called */
445:   void     *data;                    /* where to stash those structures */
446:   PetscErrorCode (*destroy)(void *); /* destroy routine */
447: } Mat_Product;

449: struct _p_Mat {
450:   PETSCHEADER(struct _MatOps);
451:   PetscLayout      rmap, cmap;
452:   void            *data;                                    /* implementation-specific data */
453:   MatFactorType    factortype;                              /* MAT_FACTOR_LU, ILU, CHOLESKY or ICC */
454:   PetscBool        trivialsymbolic;                         /* indicates the symbolic factorization doesn't actually do a symbolic factorization, it is delayed to the numeric factorization */
455:   PetscBool        canuseordering;                          /* factorization can use ordering provide to routine (most PETSc implementations) */
456:   MatOrderingType  preferredordering[MAT_FACTOR_NUM_TYPES]; /* what is the preferred (or default) ordering for the matrix solver type */
457:   PetscBool        assembled;                               /* is the matrix assembled? */
458:   PetscBool        was_assembled;                           /* new values inserted into assembled mat */
459:   PetscInt         num_ass;                                 /* number of times matrix has been assembled */
460:   PetscObjectState nonzerostate;                            /* each time new nonzeros locations are introduced into the matrix this is updated */
461:   PetscObjectState ass_nonzerostate;                        /* nonzero state at last assembly */
462:   MatInfo          info;                                    /* matrix information */
463:   InsertMode       insertmode;                              /* have values been inserted in matrix or added? */
464:   MatStash         stash, bstash;                           /* used for assembling off-proc mat emements */
465:   MatNullSpace     nullsp;                                  /* null space (operator is singular) */
466:   MatNullSpace     transnullsp;                             /* null space of transpose of operator */
467:   MatNullSpace     nearnullsp;                              /* near null space to be used by multigrid methods */
468:   PetscInt         congruentlayouts;                        /* are the rows and columns layouts congruent? */
469:   PetscBool        preallocated;
470:   MatStencilInfo   stencil; /* information for structured grid */
471:   PetscBool3       symmetric, hermitian, structurally_symmetric, spd;
472:   PetscBool        symmetry_eternal, structural_symmetry_eternal, spd_eternal;
473:   PetscBool        nooffprocentries, nooffproczerorows;
474:   PetscBool        assembly_subset; /* set by MAT_SUBSET_OFF_PROC_ENTRIES */
475:   PetscBool        submat_singleis; /* for efficient PCSetUp_ASM() */
476:   PetscBool        structure_only;
477:   PetscBool        sortedfull;      /* full, sorted rows are inserted */
478:   PetscBool        force_diagonals; /* set by MAT_FORCE_DIAGONAL_ENTRIES */
479: #if defined(PETSC_HAVE_DEVICE)
480:   PetscOffloadMask offloadmask; /* a mask which indicates where the valid matrix data is (GPU, CPU or both) */
481:   PetscBool        boundtocpu;
482:   PetscBool        bindingpropagates;
483: #endif
484:   char                *defaultrandtype;
485:   void                *spptr; /* pointer for special library like SuperLU */
486:   char                *solvertype;
487:   PetscBool            checksymmetryonassembly, checknullspaceonassembly;
488:   PetscReal            checksymmetrytol;
489:   Mat                  schur;                            /* Schur complement matrix */
490:   MatFactorSchurStatus schur_status;                     /* status of the Schur complement matrix */
491:   Mat_Redundant       *redundant;                        /* used by MatCreateRedundantMatrix() */
492:   PetscBool            erroriffailure;                   /* Generate an error if detected (for example a zero pivot) instead of returning */
493:   MatFactorError       factorerrortype;                  /* type of error in factorization */
494:   PetscReal            factorerror_zeropivot_value;      /* If numerical zero pivot was detected this is the computed value */
495:   PetscInt             factorerror_zeropivot_row;        /* Row where zero pivot was detected */
496:   PetscInt             nblocks, *bsizes;                 /* support for MatSetVariableBlockSizes() */
497:   PetscInt             p_cstart, p_rank, p_cend, n_rank; /* Information from parallel MatComputeVariableBlockEnvelope() */
498:   PetscBool            p_parallel;
499:   char                *defaultvectype;
500:   Mat_Product         *product;
501:   PetscBool            form_explicit_transpose; /* hint to generate an explicit mat tranpsose for operations like MatMultTranspose() */
502:   PetscBool            transupdated;            /* whether or not the explicitly generated transpose is up-to-date */
503:   char                *factorprefix;            /* the prefix to use with factored matrix that is created */
504:   PetscBool            hash_active;             /* indicates MatSetValues() is being handled by hashing */
505: };

507: PETSC_INTERN PetscErrorCode MatAXPY_Basic(Mat, PetscScalar, Mat, MatStructure);
508: PETSC_INTERN PetscErrorCode MatAXPY_BasicWithPreallocation(Mat, Mat, PetscScalar, Mat, MatStructure);
509: PETSC_INTERN PetscErrorCode MatAXPY_Basic_Preallocate(Mat, Mat, Mat *);
510: PETSC_INTERN PetscErrorCode MatAXPY_Dense_Nest(Mat, PetscScalar, Mat);

512: PETSC_INTERN PetscErrorCode MatSetUp_Default(Mat);

514: /*
515:     Utility for MatZeroRows
516: */
517: PETSC_INTERN PetscErrorCode MatZeroRowsMapLocal_Private(Mat, PetscInt, const PetscInt *, PetscInt *, PetscInt **);

519: /*
520:     Utility for MatView/MatLoad
521: */
522: PETSC_INTERN PetscErrorCode MatView_Binary_BlockSizes(Mat, PetscViewer);
523: PETSC_INTERN PetscErrorCode MatLoad_Binary_BlockSizes(Mat, PetscViewer);

525: /*
526:     Object for partitioning graphs
527: */

529: typedef struct _MatPartitioningOps *MatPartitioningOps;
530: struct _MatPartitioningOps {
531:   PetscErrorCode (*apply)(MatPartitioning, IS *);
532:   PetscErrorCode (*applynd)(MatPartitioning, IS *);
533:   PetscErrorCode (*setfromoptions)(MatPartitioning, PetscOptionItems *);
534:   PetscErrorCode (*destroy)(MatPartitioning);
535:   PetscErrorCode (*view)(MatPartitioning, PetscViewer);
536:   PetscErrorCode (*improve)(MatPartitioning, IS *);
537: };

539: struct _p_MatPartitioning {
540:   PETSCHEADER(struct _MatPartitioningOps);
541:   Mat        adj;
542:   PetscInt  *vertex_weights;
543:   PetscReal *part_weights;
544:   PetscInt   n;    /* number of partitions */
545:   PetscInt   ncon; /* number of vertex weights per vertex */
546:   void      *data;
547:   PetscInt   setupcalled;
548:   PetscBool  use_edge_weights; /* A flag indicates whether or not to use edge weights */
549: };

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

554: /*
555:     Object for coarsen graphs
556: */
557: typedef struct _MatCoarsenOps *MatCoarsenOps;
558: struct _MatCoarsenOps {
559:   PetscErrorCode (*apply)(MatCoarsen);
560:   PetscErrorCode (*setfromoptions)(MatCoarsen, PetscOptionItems *);
561:   PetscErrorCode (*destroy)(MatCoarsen);
562:   PetscErrorCode (*view)(MatCoarsen, PetscViewer);
563: };

565: struct _p_MatCoarsen {
566:   PETSCHEADER(struct _MatCoarsenOps);
567:   Mat   graph;
568:   void *subctx;
569:   /* */
570:   PetscBool         strict_aggs;
571:   IS                perm;
572:   PetscCoarsenData *agg_lists;
573: };

575: PETSC_EXTERN PetscErrorCode MatCoarsenMISKSetDistance(MatCoarsen, PetscInt);
576: PETSC_EXTERN PetscErrorCode MatCoarsenMISKGetDistance(MatCoarsen, PetscInt *);

578: /*
579:     Used in aijdevice.h
580: */
581: typedef struct {
582:   PetscInt    *i;
583:   PetscInt    *j;
584:   PetscScalar *a;
585:   PetscInt     n;
586:   PetscInt     ignorezeroentries;
587: } PetscCSRDataStructure;

589: /*
590:     MatFDColoring is used to compute Jacobian matrices efficiently
591:   via coloring. The data structure is explained below in an example.

593:    Color =   0    1     0    2   |   2      3       0
594:    ---------------------------------------------------
595:             00   01              |          05
596:             10   11              |   14     15               Processor  0
597:                        22    23  |          25
598:                        32    33  |
599:    ===================================================
600:                                  |   44     45     46
601:             50                   |          55               Processor 1
602:                                  |   64            66
603:    ---------------------------------------------------

605:     ncolors = 4;

607:     ncolumns      = {2,1,1,0}
608:     columns       = {{0,2},{1},{3},{}}
609:     nrows         = {4,2,3,3}
610:     rows          = {{0,1,2,3},{0,1},{1,2,3},{0,1,2}}
611:     vwscale       = {dx(0),dx(1),dx(2),dx(3)}               MPI Vec
612:     vscale        = {dx(0),dx(1),dx(2),dx(3),dx(4),dx(5)}   Seq Vec

614:     ncolumns      = {1,0,1,1}
615:     columns       = {{6},{},{4},{5}}
616:     nrows         = {3,0,2,2}
617:     rows          = {{0,1,2},{},{1,2},{1,2}}
618:     vwscale       = {dx(4),dx(5),dx(6)}              MPI Vec
619:     vscale        = {dx(0),dx(4),dx(5),dx(6)}        Seq Vec

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

624: */
625: typedef struct {
626:   PetscInt     row;
627:   PetscInt     col;
628:   PetscScalar *valaddr; /* address of value */
629: } MatEntry;

631: typedef struct {
632:   PetscInt     row;
633:   PetscScalar *valaddr; /* address of value */
634: } MatEntry2;

636: struct _p_MatFDColoring {
637:   PETSCHEADER(int);
638:   PetscInt     M, N, m;                           /* total rows, columns; local rows */
639:   PetscInt     rstart;                            /* first row owned by local processor */
640:   PetscInt     ncolors;                           /* number of colors */
641:   PetscInt    *ncolumns;                          /* number of local columns for a color */
642:   PetscInt   **columns;                           /* lists the local columns of each color (using global column numbering) */
643:   IS          *isa;                               /* these are the IS that contain the column values given in columns */
644:   PetscInt    *nrows;                             /* number of local rows for each color */
645:   MatEntry    *matentry;                          /* holds (row, column, address of value) for Jacobian matrix entry */
646:   MatEntry2   *matentry2;                         /* holds (row, address of value) for Jacobian matrix entry */
647:   PetscScalar *dy;                                /* store a block of F(x+dx)-F(x) when J is in BAIJ format */
648:   PetscReal    error_rel;                         /* square root of relative error in computing function */
649:   PetscReal    umin;                              /* minimum allowable u'dx value */
650:   Vec          w1, w2, w3;                        /* work vectors used in computing Jacobian */
651:   PetscBool    fset;                              /* indicates that the initial function value F(X) is set */
652:   PetscErrorCode (*f)(void);                      /* function that defines Jacobian */
653:   void          *fctx;                            /* optional user-defined context for use by the function f */
654:   Vec            vscale;                          /* holds FD scaling, i.e. 1/dx for each perturbed column */
655:   PetscInt       currentcolor;                    /* color for which function evaluation is being done now */
656:   const char    *htype;                           /* "wp" or "ds" */
657:   ISColoringType ctype;                           /* IS_COLORING_GLOBAL or IS_COLORING_LOCAL */
658:   PetscInt       brows, bcols;                    /* number of block rows or columns for speedup inserting the dense matrix into sparse Jacobian */
659:   PetscBool      setupcalled;                     /* true if setup has been called */
660:   PetscBool      viewed;                          /* true if the -mat_fd_coloring_view has been triggered already */
661:   void (*ftn_func_pointer)(void), *ftn_func_cntx; /* serve the same purpose as *fortran_func_pointers in PETSc objects */
662:   PetscObjectId matid;                            /* matrix this object was created with, must always be the same */
663: };

665: typedef struct _MatColoringOps *MatColoringOps;
666: struct _MatColoringOps {
667:   PetscErrorCode (*destroy)(MatColoring);
668:   PetscErrorCode (*setfromoptions)(MatColoring, PetscOptionItems *);
669:   PetscErrorCode (*view)(MatColoring, PetscViewer);
670:   PetscErrorCode (*apply)(MatColoring, ISColoring *);
671:   PetscErrorCode (*weights)(MatColoring, PetscReal **, PetscInt **);
672: };

674: struct _p_MatColoring {
675:   PETSCHEADER(struct _MatColoringOps);
676:   Mat                   mat;
677:   PetscInt              dist;         /* distance of the coloring */
678:   PetscInt              maxcolors;    /* the maximum number of colors returned, maxcolors=1 for MIS */
679:   void                 *data;         /* inner context */
680:   PetscBool             valid;        /* check to see if what is produced is a valid coloring */
681:   MatColoringWeightType weight_type;  /* type of weight computation to be performed */
682:   PetscReal            *user_weights; /* custom weights and permutation */
683:   PetscInt             *user_lperm;
684:   PetscBool             valid_iscoloring; /* check to see if matcoloring is produced a valid iscoloring */
685: };

687: struct _p_MatTransposeColoring {
688:   PETSCHEADER(int);
689:   PetscInt       M, N, m;      /* total rows, columns; local rows */
690:   PetscInt       rstart;       /* first row owned by local processor */
691:   PetscInt       ncolors;      /* number of colors */
692:   PetscInt      *ncolumns;     /* number of local columns for a color */
693:   PetscInt      *nrows;        /* number of local rows for each color */
694:   PetscInt       currentcolor; /* color for which function evaluation is being done now */
695:   ISColoringType ctype;        /* IS_COLORING_GLOBAL or IS_COLORING_LOCAL */

697:   PetscInt *colorforrow, *colorforcol; /* pointer to rows and columns */
698:   PetscInt *rows;                      /* lists the local rows for each color (using the local row numbering) */
699:   PetscInt *den2sp;                    /* maps (row,color) in the dense matrix to index of sparse matrix array a->a */
700:   PetscInt *columns;                   /* lists the local columns of each color (using global column numbering) */
701:   PetscInt  brows;                     /* number of rows for efficient implementation of MatTransColoringApplyDenToSp() */
702:   PetscInt *lstart;                    /* array used for loop over row blocks of Csparse */
703: };

705: /*
706:    Null space context for preconditioner/operators
707: */
708: struct _p_MatNullSpace {
709:   PETSCHEADER(int);
710:   PetscBool    has_cnst;
711:   PetscInt     n;
712:   Vec         *vecs;
713:   PetscScalar *alpha;                                  /* for projections */
714:   PetscErrorCode (*remove)(MatNullSpace, Vec, void *); /* for user provided removal function */
715:   void *rmctx;                                         /* context for remove() function */
716: };

718: /*
719:    Checking zero pivot for LU, ILU preconditioners.
720: */
721: typedef struct {
722:   PetscInt    nshift, nshift_max;
723:   PetscReal   shift_amount, shift_lo, shift_hi, shift_top, shift_fraction;
724:   PetscBool   newshift;
725:   PetscReal   rs; /* active row sum of abs(offdiagonals) */
726:   PetscScalar pv; /* pivot of the active row */
727: } FactorShiftCtx;

729: PETSC_EXTERN PetscErrorCode MatTransposeCheckNonzeroState_Private(Mat, Mat);

731: /*
732:  Used by MatTranspose() and potentially other functions to track the matrix used in the generation of another matrix
733: */
734: typedef struct {
735:   PetscObjectId    id;
736:   PetscObjectState state;
737:   PetscObjectState nonzerostate;
738: } MatParentState;

740: PETSC_EXTERN PetscErrorCode MatFactorDumpMatrix(Mat);
741: PETSC_INTERN PetscErrorCode MatShift_Basic(Mat, PetscScalar);
742: PETSC_INTERN PetscErrorCode MatSetBlockSizes_Default(Mat, PetscInt, PetscInt);

744: static inline PetscErrorCode MatPivotCheck_nz(PETSC_UNUSED Mat mat, const MatFactorInfo *info, FactorShiftCtx *sctx, PETSC_UNUSED PetscInt row)
745: {
746:   PetscReal _rs   = sctx->rs;
747:   PetscReal _zero = info->zeropivot * _rs;

749:   PetscFunctionBegin;
750:   if (PetscAbsScalar(sctx->pv) <= _zero && !PetscIsNanScalar(sctx->pv)) {
751:     /* force |diag| > zeropivot*rs */
752:     if (!sctx->nshift) sctx->shift_amount = info->shiftamount;
753:     else sctx->shift_amount *= 2.0;
754:     sctx->newshift = PETSC_TRUE;
755:     (sctx->nshift)++;
756:   } else {
757:     sctx->newshift = PETSC_FALSE;
758:   }
759:   PetscFunctionReturn(PETSC_SUCCESS);
760: }

762: static inline PetscErrorCode MatPivotCheck_pd(PETSC_UNUSED Mat mat, const MatFactorInfo *info, FactorShiftCtx *sctx, PETSC_UNUSED PetscInt row)
763: {
764:   PetscReal _rs   = sctx->rs;
765:   PetscReal _zero = info->zeropivot * _rs;

767:   PetscFunctionBegin;
768:   if (PetscRealPart(sctx->pv) <= _zero && !PetscIsNanScalar(sctx->pv)) {
769:     /* force matfactor to be diagonally dominant */
770:     if (sctx->nshift == sctx->nshift_max) {
771:       sctx->shift_fraction = sctx->shift_hi;
772:     } else {
773:       sctx->shift_lo       = sctx->shift_fraction;
774:       sctx->shift_fraction = (sctx->shift_hi + sctx->shift_lo) / (PetscReal)2.;
775:     }
776:     sctx->shift_amount = sctx->shift_fraction * sctx->shift_top;
777:     sctx->nshift++;
778:     sctx->newshift = PETSC_TRUE;
779:   } else {
780:     sctx->newshift = PETSC_FALSE;
781:   }
782:   PetscFunctionReturn(PETSC_SUCCESS);
783: }

785: static inline PetscErrorCode MatPivotCheck_inblocks(PETSC_UNUSED Mat mat, const MatFactorInfo *info, FactorShiftCtx *sctx, PETSC_UNUSED PetscInt row)
786: {
787:   PetscReal _zero = info->zeropivot;

789:   PetscFunctionBegin;
790:   if (PetscAbsScalar(sctx->pv) <= _zero && !PetscIsNanScalar(sctx->pv)) {
791:     sctx->pv += info->shiftamount;
792:     sctx->shift_amount = 0.0;
793:     sctx->nshift++;
794:   }
795:   sctx->newshift = PETSC_FALSE;
796:   PetscFunctionReturn(PETSC_SUCCESS);
797: }

799: static inline PetscErrorCode MatPivotCheck_none(Mat fact, Mat mat, const MatFactorInfo *info, FactorShiftCtx *sctx, PetscInt row)
800: {
801:   PetscReal _zero = info->zeropivot;

803:   PetscFunctionBegin;
804:   sctx->newshift = PETSC_FALSE;
805:   if (PetscAbsScalar(sctx->pv) <= _zero && !PetscIsNanScalar(sctx->pv)) {
806:     PetscCheck(!mat->erroriffailure, PETSC_COMM_SELF, PETSC_ERR_MAT_LU_ZRPVT, "Zero pivot row %" PetscInt_FMT " value %g tolerance %g", row, (double)PetscAbsScalar(sctx->pv), (double)_zero);
807:     PetscCall(PetscInfo(mat, "Detected zero pivot in factorization in row %" PetscInt_FMT " value %g tolerance %g\n", row, (double)PetscAbsScalar(sctx->pv), (double)_zero));
808:     fact->factorerrortype             = MAT_FACTOR_NUMERIC_ZEROPIVOT;
809:     fact->factorerror_zeropivot_value = PetscAbsScalar(sctx->pv);
810:     fact->factorerror_zeropivot_row   = row;
811:   }
812:   PetscFunctionReturn(PETSC_SUCCESS);
813: }

815: static inline PetscErrorCode MatPivotCheck(Mat fact, Mat mat, const MatFactorInfo *info, FactorShiftCtx *sctx, PetscInt row)
816: {
817:   PetscFunctionBegin;
818:   if (info->shifttype == (PetscReal)MAT_SHIFT_NONZERO) PetscCall(MatPivotCheck_nz(mat, info, sctx, row));
819:   else if (info->shifttype == (PetscReal)MAT_SHIFT_POSITIVE_DEFINITE) PetscCall(MatPivotCheck_pd(mat, info, sctx, row));
820:   else if (info->shifttype == (PetscReal)MAT_SHIFT_INBLOCKS) PetscCall(MatPivotCheck_inblocks(mat, info, sctx, row));
821:   else PetscCall(MatPivotCheck_none(fact, mat, info, sctx, row));
822:   PetscFunctionReturn(PETSC_SUCCESS);
823: }

825: #include <petscbt.h>
826: /*
827:   Create and initialize a linked list
828:   Input Parameters:
829:     idx_start - starting index of the list
830:     lnk_max   - max value of lnk indicating the end of the list
831:     nlnk      - max length of the list
832:   Output Parameters:
833:     lnk       - list initialized
834:     bt        - PetscBT (bitarray) with all bits set to false
835:     lnk_empty - flg indicating the list is empty
836: */
837: #define PetscLLCreate(idx_start, lnk_max, nlnk, lnk, bt) ((PetscErrorCode)(PetscMalloc1(nlnk, &lnk) || PetscBTCreate(nlnk, &(bt)) || (lnk[idx_start] = lnk_max, PETSC_SUCCESS)))

839: #define PetscLLCreate_new(idx_start, lnk_max, nlnk, lnk, bt, lnk_empty) ((PetscErrorCode)(PetscMalloc1(nlnk, &lnk) || PetscBTCreate(nlnk, &(bt)) || (lnk_empty = PETSC_TRUE, 0) || (lnk[idx_start] = lnk_max, PETSC_SUCCESS)))

841: static inline PetscErrorCode PetscLLInsertLocation_Private(PetscBool assume_sorted, PetscInt k, PetscInt idx_start, PetscInt entry, PetscInt *PETSC_RESTRICT nlnk, PetscInt *PETSC_RESTRICT lnkdata, PetscInt *PETSC_RESTRICT lnk)
842: {
843:   PetscInt location;

845:   PetscFunctionBegin;
846:   /* start from the beginning if entry < previous entry */
847:   if (!assume_sorted && k && entry < *lnkdata) *lnkdata = idx_start;
848:   /* search for insertion location */
849:   do {
850:     location = *lnkdata;
851:     *lnkdata = lnk[location];
852:   } while (entry > *lnkdata);
853:   /* insertion location is found, add entry into lnk */
854:   lnk[location] = entry;
855:   lnk[entry]    = *lnkdata;
856:   ++(*nlnk);
857:   *lnkdata = entry; /* next search starts from here if next_entry > entry */
858:   PetscFunctionReturn(PETSC_SUCCESS);
859: }

861: static inline PetscErrorCode PetscLLAdd_Private(PetscInt nidx, const PetscInt *PETSC_RESTRICT indices, PetscInt idx_start, PetscInt *PETSC_RESTRICT nlnk, PetscInt *PETSC_RESTRICT lnk, PetscBT bt, PetscBool assume_sorted)
862: {
863:   PetscFunctionBegin;
864:   *nlnk = 0;
865:   for (PetscInt k = 0, lnkdata = idx_start; k < nidx; ++k) {
866:     const PetscInt entry = indices[k];

868:     if (!PetscBTLookupSet(bt, entry)) PetscCall(PetscLLInsertLocation_Private(assume_sorted, k, idx_start, entry, nlnk, &lnkdata, lnk));
869:   }
870:   PetscFunctionReturn(PETSC_SUCCESS);
871: }

873: /*
874:   Add an index set into a sorted linked list
875:   Input Parameters:
876:     nidx      - number of input indices
877:     indices   - integer array
878:     idx_start - starting index of the list
879:     lnk       - linked list(an integer array) that is created
880:     bt        - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
881:   output Parameters:
882:     nlnk      - number of newly added indices
883:     lnk       - the sorted(increasing order) linked list containing new and non-redundate entries from indices
884:     bt        - updated PetscBT (bitarray)
885: */
886: static inline PetscErrorCode PetscLLAdd(PetscInt nidx, const PetscInt *PETSC_RESTRICT indices, PetscInt idx_start, PetscInt *PETSC_RESTRICT nlnk, PetscInt *PETSC_RESTRICT lnk, PetscBT bt)
887: {
888:   PetscFunctionBegin;
889:   PetscCall(PetscLLAdd_Private(nidx, indices, idx_start, nlnk, lnk, bt, PETSC_FALSE));
890:   PetscFunctionReturn(PETSC_SUCCESS);
891: }

893: /*
894:   Add a SORTED ascending index set into a sorted linked list - same as PetscLLAdd() bus skip 'if (_k && _entry < _lnkdata) _lnkdata  = idx_start;'
895:   Input Parameters:
896:     nidx      - number of input indices
897:     indices   - sorted integer array
898:     idx_start - starting index of the list
899:     lnk       - linked list(an integer array) that is created
900:     bt        - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
901:   output Parameters:
902:     nlnk      - number of newly added indices
903:     lnk       - the sorted(increasing order) linked list containing new and non-redundate entries from indices
904:     bt        - updated PetscBT (bitarray)
905: */
906: static inline PetscErrorCode PetscLLAddSorted(PetscInt nidx, const PetscInt *PETSC_RESTRICT indices, PetscInt idx_start, PetscInt *PETSC_RESTRICT nlnk, PetscInt *PETSC_RESTRICT lnk, PetscBT bt)
907: {
908:   PetscFunctionBegin;
909:   PetscCall(PetscLLAdd_Private(nidx, indices, idx_start, nlnk, lnk, bt, PETSC_TRUE));
910:   PetscFunctionReturn(PETSC_SUCCESS);
911: }

913: /*
914:   Add a permuted index set into a sorted linked list
915:   Input Parameters:
916:     nidx      - number of input indices
917:     indices   - integer array
918:     perm      - permutation of indices
919:     idx_start - starting index of the list
920:     lnk       - linked list(an integer array) that is created
921:     bt        - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
922:   output Parameters:
923:     nlnk      - number of newly added indices
924:     lnk       - the sorted(increasing order) linked list containing new and non-redundate entries from indices
925:     bt        - updated PetscBT (bitarray)
926: */
927: static inline PetscErrorCode PetscLLAddPerm(PetscInt nidx, const PetscInt *PETSC_RESTRICT indices, const PetscInt *PETSC_RESTRICT perm, PetscInt idx_start, PetscInt *PETSC_RESTRICT nlnk, PetscInt *PETSC_RESTRICT lnk, PetscBT bt)
928: {
929:   PetscFunctionBegin;
930:   *nlnk = 0;
931:   for (PetscInt k = 0, lnkdata = idx_start; k < nidx; ++k) {
932:     const PetscInt entry = perm[indices[k]];

934:     if (!PetscBTLookupSet(bt, entry)) PetscCall(PetscLLInsertLocation_Private(PETSC_FALSE, k, idx_start, entry, nlnk, &lnkdata, lnk));
935:   }
936:   PetscFunctionReturn(PETSC_SUCCESS);
937: }

939: #if 0
940: /* this appears to be unused? */
941: static inline PetscErrorCode PetscLLAddSorted_new(PetscInt nidx, PetscInt *indices, PetscInt idx_start, PetscBool *lnk_empty, PetscInt *nlnk, PetscInt *lnk, PetscBT bt)
942: {
943:   PetscInt lnkdata = idx_start;

945:   PetscFunctionBegin;
946:   if (*lnk_empty) {
947:     for (PetscInt k = 0; k < nidx; ++k) {
948:       const PetscInt entry = indices[k], location = lnkdata;

950:       PetscCall(PetscBTSet(bt,entry)); /* mark the new entry */
951:       lnkdata       = lnk[location];
952:       /* insertion location is found, add entry into lnk */
953:       lnk[location] = entry;
954:       lnk[entry]    = lnkdata;
955:       lnkdata       = entry; /* next search starts from here */
956:     }
957:     /* lnk[indices[nidx-1]] = lnk[idx_start];
958:        lnk[idx_start]       = indices[0];
959:        PetscCall(PetscBTSet(bt,indices[0]));
960:        for (_k=1; _k<nidx; _k++) {
961:        PetscCall(PetscBTSet(bt,indices[_k]));
962:        lnk[indices[_k-1]] = indices[_k];
963:        }
964:     */
965:     *nlnk      = nidx;
966:     *lnk_empty = PETSC_FALSE;
967:   } else {
968:     *nlnk = 0;
969:     for (PetscInt k = 0; k < nidx; ++k) {
970:       const PetscInt entry = indices[k];

972:       if (!PetscBTLookupSet(bt,entry)) PetscCall(PetscLLInsertLocation_Private(PETSC_TRUE,k,idx_start,entry,nlnk,&lnkdata,lnk));
973:     }
974:   }
975:   PetscFunctionReturn(PETSC_SUCCESS);
976: }
977: #endif

979: /*
980:   Add a SORTED index set into a sorted linked list used for LUFactorSymbolic()
981:   Same as PetscLLAddSorted() with an additional operation:
982:        count the number of input indices that are no larger than 'diag'
983:   Input Parameters:
984:     indices   - sorted integer array
985:     idx_start - starting index of the list, index of pivot row
986:     lnk       - linked list(an integer array) that is created
987:     bt        - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
988:     diag      - index of the active row in LUFactorSymbolic
989:     nzbd      - number of input indices with indices <= idx_start
990:     im        - im[idx_start] is initialized as num of nonzero entries in row=idx_start
991:   output Parameters:
992:     nlnk      - number of newly added indices
993:     lnk       - the sorted(increasing order) linked list containing new and non-redundate entries from indices
994:     bt        - updated PetscBT (bitarray)
995:     im        - im[idx_start]: unchanged if diag is not an entry
996:                              : num of entries with indices <= diag if diag is an entry
997: */
998: static inline PetscErrorCode PetscLLAddSortedLU(const PetscInt *PETSC_RESTRICT indices, PetscInt idx_start, PetscInt *PETSC_RESTRICT nlnk, PetscInt *PETSC_RESTRICT lnk, PetscBT bt, PetscInt diag, PetscInt nzbd, PetscInt *PETSC_RESTRICT im)
999: {
1000:   const PetscInt nidx = im[idx_start] - nzbd; /* num of entries with idx_start < index <= diag */

1002:   PetscFunctionBegin;
1003:   *nlnk = 0;
1004:   for (PetscInt k = 0, lnkdata = idx_start; k < nidx; ++k) {
1005:     const PetscInt entry = indices[k];

1007:     ++nzbd;
1008:     if (entry == diag) im[idx_start] = nzbd;
1009:     if (!PetscBTLookupSet(bt, entry)) PetscCall(PetscLLInsertLocation_Private(PETSC_TRUE, k, idx_start, entry, nlnk, &lnkdata, lnk));
1010:   }
1011:   PetscFunctionReturn(PETSC_SUCCESS);
1012: }

1014: /*
1015:   Copy data on the list into an array, then initialize the list
1016:   Input Parameters:
1017:     idx_start - starting index of the list
1018:     lnk_max   - max value of lnk indicating the end of the list
1019:     nlnk      - number of data on the list to be copied
1020:     lnk       - linked list
1021:     bt        - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
1022:   output Parameters:
1023:     indices   - array that contains the copied data
1024:     lnk       - linked list that is cleaned and initialize
1025:     bt        - PetscBT (bitarray) with all bits set to false
1026: */
1027: static inline PetscErrorCode PetscLLClean(PetscInt idx_start, PetscInt lnk_max, PetscInt nlnk, PetscInt *PETSC_RESTRICT lnk, PetscInt *PETSC_RESTRICT indices, PetscBT bt)
1028: {
1029:   PetscFunctionBegin;
1030:   for (PetscInt j = 0, idx = idx_start; j < nlnk; ++j) {
1031:     idx        = lnk[idx];
1032:     indices[j] = idx;
1033:     PetscCall(PetscBTClear(bt, idx));
1034:   }
1035:   lnk[idx_start] = lnk_max;
1036:   PetscFunctionReturn(PETSC_SUCCESS);
1037: }

1039: /*
1040:   Free memories used by the list
1041: */
1042: #define PetscLLDestroy(lnk, bt) ((PetscErrorCode)(PetscFree(lnk) || PetscBTDestroy(&(bt))))

1044: /* Routines below are used for incomplete matrix factorization */
1045: /*
1046:   Create and initialize a linked list and its levels
1047:   Input Parameters:
1048:     idx_start - starting index of the list
1049:     lnk_max   - max value of lnk indicating the end of the list
1050:     nlnk      - max length of the list
1051:   Output Parameters:
1052:     lnk       - list initialized
1053:     lnk_lvl   - array of size nlnk for storing levels of lnk
1054:     bt        - PetscBT (bitarray) with all bits set to false
1055: */
1056: #define PetscIncompleteLLCreate(idx_start, lnk_max, nlnk, lnk, lnk_lvl, bt) \
1057:   ((PetscErrorCode)(PetscIntMultError(2, nlnk, NULL) || PetscMalloc1(2 * nlnk, &lnk) || PetscBTCreate(nlnk, &(bt)) || (lnk[idx_start] = lnk_max, lnk_lvl = lnk + nlnk, PETSC_SUCCESS)))

1059: static inline PetscErrorCode PetscIncompleteLLInsertLocation_Private(PetscBool assume_sorted, PetscInt k, PetscInt idx_start, PetscInt entry, PetscInt *PETSC_RESTRICT nlnk, PetscInt *PETSC_RESTRICT lnkdata, PetscInt *PETSC_RESTRICT lnk, PetscInt *PETSC_RESTRICT lnklvl, PetscInt newval)
1060: {
1061:   PetscFunctionBegin;
1062:   PetscCall(PetscLLInsertLocation_Private(assume_sorted, k, idx_start, entry, nlnk, lnkdata, lnk));
1063:   lnklvl[entry] = newval;
1064:   PetscFunctionReturn(PETSC_SUCCESS);
1065: }

1067: /*
1068:   Initialize a sorted linked list used for ILU and ICC
1069:   Input Parameters:
1070:     nidx      - number of input idx
1071:     idx       - integer array used for storing column indices
1072:     idx_start - starting index of the list
1073:     perm      - indices of an IS
1074:     lnk       - linked list(an integer array) that is created
1075:     lnklvl    - levels of lnk
1076:     bt        - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
1077:   output Parameters:
1078:     nlnk     - number of newly added idx
1079:     lnk      - the sorted(increasing order) linked list containing new and non-redundate entries from idx
1080:     lnklvl   - levels of lnk
1081:     bt       - updated PetscBT (bitarray)
1082: */
1083: static inline PetscErrorCode PetscIncompleteLLInit(PetscInt nidx, const PetscInt *PETSC_RESTRICT idx, PetscInt idx_start, const PetscInt *PETSC_RESTRICT perm, PetscInt *PETSC_RESTRICT nlnk, PetscInt *PETSC_RESTRICT lnk, PetscInt *PETSC_RESTRICT lnklvl, PetscBT bt)
1084: {
1085:   PetscFunctionBegin;
1086:   *nlnk = 0;
1087:   for (PetscInt k = 0, lnkdata = idx_start; k < nidx; ++k) {
1088:     const PetscInt entry = perm[idx[k]];

1090:     if (!PetscBTLookupSet(bt, entry)) PetscCall(PetscIncompleteLLInsertLocation_Private(PETSC_FALSE, k, idx_start, entry, nlnk, &lnkdata, lnk, lnklvl, 0));
1091:   }
1092:   PetscFunctionReturn(PETSC_SUCCESS);
1093: }

1095: static inline PetscErrorCode PetscIncompleteLLAdd_Private(PetscInt nidx, const PetscInt *PETSC_RESTRICT idx, PetscReal level, const PetscInt *PETSC_RESTRICT idxlvl, PetscInt idx_start, PetscInt *PETSC_RESTRICT nlnk, PetscInt *PETSC_RESTRICT lnk, PetscInt *PETSC_RESTRICT lnklvl, PetscBT bt, PetscInt prow_offset, PetscBool assume_sorted)
1096: {
1097:   PetscFunctionBegin;
1098:   *nlnk = 0;
1099:   for (PetscInt k = 0, lnkdata = idx_start; k < nidx; ++k) {
1100:     const PetscInt incrlev = idxlvl[k] + prow_offset + 1;

1102:     if (incrlev <= level) {
1103:       const PetscInt entry = idx[k];

1105:       if (!PetscBTLookupSet(bt, entry)) PetscCall(PetscIncompleteLLInsertLocation_Private(assume_sorted, k, idx_start, entry, nlnk, &lnkdata, lnk, lnklvl, incrlev));
1106:       else if (lnklvl[entry] > incrlev) lnklvl[entry] = incrlev; /* existing entry */
1107:     }
1108:   }
1109:   PetscFunctionReturn(PETSC_SUCCESS);
1110: }

1112: /*
1113:   Add a SORTED index set into a sorted linked list for ICC
1114:   Input Parameters:
1115:     nidx      - number of input indices
1116:     idx       - sorted integer array used for storing column indices
1117:     level     - level of fill, e.g., ICC(level)
1118:     idxlvl    - level of idx
1119:     idx_start - starting index of the list
1120:     lnk       - linked list(an integer array) that is created
1121:     lnklvl    - levels of lnk
1122:     bt        - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
1123:     idxlvl_prow - idxlvl[prow], where prow is the row number of the idx
1124:   output Parameters:
1125:     nlnk   - number of newly added indices
1126:     lnk    - the sorted(increasing order) linked list containing new and non-redundate entries from idx
1127:     lnklvl - levels of lnk
1128:     bt     - updated PetscBT (bitarray)
1129:   Note: the level of U(i,j) is set as lvl(i,j) = min{ lvl(i,j), lvl(prow,i)+lvl(prow,j)+1)
1130:         where idx = non-zero columns of U(prow,prow+1:n-1), prow<i
1131: */
1132: static inline PetscErrorCode PetscICCLLAddSorted(PetscInt nidx, const PetscInt *PETSC_RESTRICT idx, PetscReal level, const PetscInt *PETSC_RESTRICT idxlvl, PetscInt idx_start, PetscInt *PETSC_RESTRICT nlnk, PetscInt *PETSC_RESTRICT lnk, PetscInt *PETSC_RESTRICT lnklvl, PetscBT bt, PetscInt idxlvl_prow)
1133: {
1134:   PetscFunctionBegin;
1135:   PetscCall(PetscIncompleteLLAdd_Private(nidx, idx, level, idxlvl, idx_start, nlnk, lnk, lnklvl, bt, idxlvl_prow, PETSC_TRUE));
1136:   PetscFunctionReturn(PETSC_SUCCESS);
1137: }

1139: /*
1140:   Add a SORTED index set into a sorted linked list for ILU
1141:   Input Parameters:
1142:     nidx      - number of input indices
1143:     idx       - sorted integer array used for storing column indices
1144:     level     - level of fill, e.g., ICC(level)
1145:     idxlvl    - level of idx
1146:     idx_start - starting index of the list
1147:     lnk       - linked list(an integer array) that is created
1148:     lnklvl    - levels of lnk
1149:     bt        - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
1150:     prow      - the row number of idx
1151:   output Parameters:
1152:     nlnk     - number of newly added idx
1153:     lnk      - the sorted(increasing order) linked list containing new and non-redundate entries from idx
1154:     lnklvl   - levels of lnk
1155:     bt       - updated PetscBT (bitarray)

1157:   Note: the level of factor(i,j) is set as lvl(i,j) = min{ lvl(i,j), lvl(i,prow)+lvl(prow,j)+1)
1158:         where idx = non-zero columns of U(prow,prow+1:n-1), prow<i
1159: */
1160: static inline PetscErrorCode PetscILULLAddSorted(PetscInt nidx, const PetscInt *PETSC_RESTRICT idx, PetscInt level, const PetscInt *PETSC_RESTRICT idxlvl, PetscInt idx_start, PetscInt *PETSC_RESTRICT nlnk, PetscInt *PETSC_RESTRICT lnk, PetscInt *PETSC_RESTRICT lnklvl, PetscBT bt, PetscInt prow)
1161: {
1162:   PetscFunctionBegin;
1163:   PetscCall(PetscIncompleteLLAdd_Private(nidx, idx, level, idxlvl, idx_start, nlnk, lnk, lnklvl, bt, lnklvl[prow], PETSC_TRUE));
1164:   PetscFunctionReturn(PETSC_SUCCESS);
1165: }

1167: /*
1168:   Add a index set into a sorted linked list
1169:   Input Parameters:
1170:     nidx      - number of input idx
1171:     idx   - integer array used for storing column indices
1172:     level     - level of fill, e.g., ICC(level)
1173:     idxlvl - level of idx
1174:     idx_start - starting index of the list
1175:     lnk       - linked list(an integer array) that is created
1176:     lnklvl   - levels of lnk
1177:     bt        - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
1178:   output Parameters:
1179:     nlnk      - number of newly added idx
1180:     lnk       - the sorted(increasing order) linked list containing new and non-redundate entries from idx
1181:     lnklvl   - levels of lnk
1182:     bt        - updated PetscBT (bitarray)
1183: */
1184: static inline PetscErrorCode PetscIncompleteLLAdd(PetscInt nidx, const PetscInt *PETSC_RESTRICT idx, PetscReal level, const PetscInt *PETSC_RESTRICT idxlvl, PetscInt idx_start, PetscInt *PETSC_RESTRICT nlnk, PetscInt *PETSC_RESTRICT lnk, PetscInt *PETSC_RESTRICT lnklvl, PetscBT bt)
1185: {
1186:   PetscFunctionBegin;
1187:   PetscCall(PetscIncompleteLLAdd_Private(nidx, idx, level, idxlvl, idx_start, nlnk, lnk, lnklvl, bt, 0, PETSC_FALSE));
1188:   PetscFunctionReturn(PETSC_SUCCESS);
1189: }

1191: /*
1192:   Add a SORTED index set into a sorted linked list
1193:   Input Parameters:
1194:     nidx      - number of input indices
1195:     idx   - sorted integer array used for storing column indices
1196:     level     - level of fill, e.g., ICC(level)
1197:     idxlvl - level of idx
1198:     idx_start - starting index of the list
1199:     lnk       - linked list(an integer array) that is created
1200:     lnklvl    - levels of lnk
1201:     bt        - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
1202:   output Parameters:
1203:     nlnk      - number of newly added idx
1204:     lnk       - the sorted(increasing order) linked list containing new and non-redundate entries from idx
1205:     lnklvl    - levels of lnk
1206:     bt        - updated PetscBT (bitarray)
1207: */
1208: static inline PetscErrorCode PetscIncompleteLLAddSorted(PetscInt nidx, const PetscInt *PETSC_RESTRICT idx, PetscReal level, const PetscInt *PETSC_RESTRICT idxlvl, PetscInt idx_start, PetscInt *PETSC_RESTRICT nlnk, PetscInt *PETSC_RESTRICT lnk, PetscInt *PETSC_RESTRICT lnklvl, PetscBT bt)
1209: {
1210:   PetscFunctionBegin;
1211:   PetscCall(PetscIncompleteLLAdd_Private(nidx, idx, level, idxlvl, idx_start, nlnk, lnk, lnklvl, bt, 0, PETSC_TRUE));
1212:   PetscFunctionReturn(PETSC_SUCCESS);
1213: }

1215: /*
1216:   Copy data on the list into an array, then initialize the list
1217:   Input Parameters:
1218:     idx_start - starting index of the list
1219:     lnk_max   - max value of lnk indicating the end of the list
1220:     nlnk      - number of data on the list to be copied
1221:     lnk       - linked list
1222:     lnklvl    - level of lnk
1223:     bt        - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
1224:   output Parameters:
1225:     indices - array that contains the copied data
1226:     lnk     - linked list that is cleaned and initialize
1227:     lnklvl  - level of lnk that is reinitialized
1228:     bt      - PetscBT (bitarray) with all bits set to false
1229: */
1230: static inline PetscErrorCode PetscIncompleteLLClean(PetscInt idx_start, PetscInt lnk_max, PetscInt nlnk, PetscInt *PETSC_RESTRICT lnk, PetscInt *PETSC_RESTRICT lnklvl, PetscInt *PETSC_RESTRICT indices, PetscInt *PETSC_RESTRICT indiceslvl, PetscBT bt)
1231: {
1232:   PetscFunctionBegin;
1233:   for (PetscInt j = 0, idx = idx_start; j < nlnk; ++j) {
1234:     idx           = lnk[idx];
1235:     indices[j]    = idx;
1236:     indiceslvl[j] = lnklvl[idx];
1237:     lnklvl[idx]   = -1;
1238:     PetscCall(PetscBTClear(bt, idx));
1239:   }
1240:   lnk[idx_start] = lnk_max;
1241:   PetscFunctionReturn(PETSC_SUCCESS);
1242: }

1244: /*
1245:   Free memories used by the list
1246: */
1247: #define PetscIncompleteLLDestroy(lnk, bt) ((PetscErrorCode)(PetscFree(lnk) || PetscBTDestroy(&(bt))))

1249: #if !defined(PETSC_CLANG_STATIC_ANALYZER)
1250:   #define MatCheckSameLocalSize(A, ar1, B, ar2) \
1251:     do { \
1252:       PetscCheckSameComm(A, ar1, B, ar2); \
1253:       PetscCheck(((A)->rmap->n == (B)->rmap->n) && ((A)->cmap->n == (B)->cmap->n), PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Incompatible matrix local sizes: parameter # %d (%" PetscInt_FMT " x %" PetscInt_FMT ") != parameter # %d (%" PetscInt_FMT " x %" PetscInt_FMT ")", ar1, \
1254:                  (A)->rmap->n, (A)->cmap->n, ar2, (B)->rmap->n, (B)->cmap->n); \
1255:     } while (0)
1256:   #define MatCheckSameSize(A, ar1, B, ar2) \
1257:     do { \
1258:       PetscCheck(((A)->rmap->N == (B)->rmap->N) && ((A)->cmap->N == (B)->cmap->N), PetscObjectComm((PetscObject)(A)), PETSC_ERR_ARG_INCOMP, "Incompatible matrix global sizes: parameter # %d (%" PetscInt_FMT " x %" PetscInt_FMT ") != parameter # %d (%" PetscInt_FMT " x %" PetscInt_FMT ")", ar1, \
1259:                  (A)->rmap->N, (A)->cmap->N, ar2, (B)->rmap->N, (B)->cmap->N); \
1260:       MatCheckSameLocalSize(A, ar1, B, ar2); \
1261:     } while (0)
1262: #else
1263: template <typename Tm>
1264: extern void MatCheckSameLocalSize(Tm, int, Tm, int);
1265: template <typename Tm>
1266: extern void MatCheckSameSize(Tm, int, Tm, int);
1267: #endif

1269: #define VecCheckMatCompatible(M, x, ar1, b, ar2) \
1270:   do { \
1271:     PetscCheck((M)->cmap->N == (x)->map->N, PetscObjectComm((PetscObject)(M)), PETSC_ERR_ARG_SIZ, "Vector global length incompatible with matrix: parameter # %d global size %" PetscInt_FMT " != matrix column global size %" PetscInt_FMT, ar1, (x)->map->N, \
1272:                (M)->cmap->N); \
1273:     PetscCheck((M)->rmap->N == (b)->map->N, PetscObjectComm((PetscObject)(M)), PETSC_ERR_ARG_SIZ, "Vector global length incompatible with matrix: parameter # %d global size %" PetscInt_FMT " != matrix row global size %" PetscInt_FMT, ar2, (b)->map->N, \
1274:                (M)->rmap->N); \
1275:   } while (0)

1277: /* -------------------------------------------------------------------------------------------------------*/
1278: /*
1279:   Create and initialize a condensed linked list -
1280:     same as PetscLLCreate(), but uses a scalable array 'lnk' with size of max number of entries, not O(N).
1281:     Barry suggested this approach (Dec. 6, 2011):
1282:       I've thought of an alternative way of representing a linked list that is efficient but doesn't have the O(N) scaling issue
1283:       (it may be faster than the O(N) even sequentially due to less crazy memory access).

1285:       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
1286:       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
1287:       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
1288:       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.
1289:       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
1290:       to each other so memory access is much better than using the big array.

1292:   Example:
1293:      nlnk_max=5, lnk_max=36:
1294:      Initial list: [0, 0 | 36, 2 | 0, 0 | 0, 0 | 0, 0 | 0, 0 | 0, 0]
1295:      here, head_node has index 2 with value lnk[2]=lnk_max=36,
1296:            0-th entry is used to store the number of entries in the list,
1297:      The initial lnk represents head -> tail(marked by 36) with number of entries = lnk[0]=0.

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

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

1307:   Input Parameters:
1308:     nlnk_max  - max length of the list
1309:     lnk_max   - max value of the entries
1310:   Output Parameters:
1311:     lnk       - list created and initialized
1312:     bt        - PetscBT (bitarray) with all bits set to false. Note: bt has size lnk_max, not nln_max!
1313: */
1314: static inline PetscErrorCode PetscLLCondensedCreate(PetscInt nlnk_max, PetscInt lnk_max, PetscInt **lnk, PetscBT *bt)
1315: {
1316:   PetscInt *llnk, lsize = 0;

1318:   PetscFunctionBegin;
1319:   PetscCall(PetscIntMultError(2, nlnk_max + 2, &lsize));
1320:   PetscCall(PetscMalloc1(lsize, lnk));
1321:   PetscCall(PetscBTCreate(lnk_max, bt));
1322:   llnk    = *lnk;
1323:   llnk[0] = 0;       /* number of entries on the list */
1324:   llnk[2] = lnk_max; /* value in the head node */
1325:   llnk[3] = 2;       /* next for the head node */
1326:   PetscFunctionReturn(PETSC_SUCCESS);
1327: }

1329: /*
1330:   Add a SORTED ascending index set into a sorted linked list. See PetscLLCondensedCreate() for detailed description.
1331:   Input Parameters:
1332:     nidx      - number of input indices
1333:     indices   - sorted integer array
1334:     lnk       - condensed linked list(an integer array) that is created
1335:     bt        - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
1336:   output Parameters:
1337:     lnk       - the sorted(increasing order) linked list containing previous and newly added non-redundate indices
1338:     bt        - updated PetscBT (bitarray)
1339: */
1340: static inline PetscErrorCode PetscLLCondensedAddSorted(PetscInt nidx, const PetscInt indices[], PetscInt lnk[], PetscBT bt)
1341: {
1342:   PetscInt location = 2;      /* head */
1343:   PetscInt nlnk     = lnk[0]; /* num of entries on the input lnk */

1345:   PetscFunctionBegin;
1346:   for (PetscInt k = 0; k < nidx; k++) {
1347:     const PetscInt entry = indices[k];
1348:     if (!PetscBTLookupSet(bt, entry)) { /* new entry */
1349:       PetscInt next, lnkdata;

1351:       /* search for insertion location */
1352:       do {
1353:         next     = location + 1;  /* link from previous node to next node */
1354:         location = lnk[next];     /* idx of next node */
1355:         lnkdata  = lnk[location]; /* value of next node */
1356:       } while (entry > lnkdata);
1357:       /* insertion location is found, add entry into lnk */
1358:       const PetscInt newnode = 2 * (nlnk + 2); /* index for this new node */
1359:       lnk[next]              = newnode;        /* connect previous node to the new node */
1360:       lnk[newnode]           = entry;          /* set value of the new node */
1361:       lnk[newnode + 1]       = location;       /* connect new node to next node */
1362:       location               = newnode;        /* next search starts from the new node */
1363:       nlnk++;
1364:     }
1365:   }
1366:   lnk[0] = nlnk; /* number of entries in the list */
1367:   PetscFunctionReturn(PETSC_SUCCESS);
1368: }

1370: static inline PetscErrorCode PetscLLCondensedClean(PetscInt lnk_max, PETSC_UNUSED PetscInt nidx, PetscInt *indices, PetscInt lnk[], PetscBT bt)
1371: {
1372:   const PetscInt nlnk = lnk[0]; /* num of entries on the list */
1373:   PetscInt       next = lnk[3]; /* head node */

1375:   PetscFunctionBegin;
1376:   for (PetscInt k = 0; k < nlnk; k++) {
1377:     indices[k] = lnk[next];
1378:     next       = lnk[next + 1];
1379:     PetscCall(PetscBTClear(bt, indices[k]));
1380:   }
1381:   lnk[0] = 0;       /* num of entries on the list */
1382:   lnk[2] = lnk_max; /* initialize head node */
1383:   lnk[3] = 2;       /* head node */
1384:   PetscFunctionReturn(PETSC_SUCCESS);
1385: }

1387: static inline PetscErrorCode PetscLLCondensedView(PetscInt *lnk)
1388: {
1389:   PetscFunctionBegin;
1390:   PetscCall(PetscPrintf(PETSC_COMM_SELF, "LLCondensed of size %" PetscInt_FMT ", (val,  next)\n", lnk[0]));
1391:   for (PetscInt k = 2; k < lnk[0] + 2; ++k) PetscCall(PetscPrintf(PETSC_COMM_SELF, " %" PetscInt_FMT ": (%" PetscInt_FMT ", %" PetscInt_FMT ")\n", 2 * k, lnk[2 * k], lnk[2 * k + 1]));
1392:   PetscFunctionReturn(PETSC_SUCCESS);
1393: }

1395: /*
1396:   Free memories used by the list
1397: */
1398: static inline PetscErrorCode PetscLLCondensedDestroy(PetscInt *lnk, PetscBT bt)
1399: {
1400:   PetscFunctionBegin;
1401:   PetscCall(PetscFree(lnk));
1402:   PetscCall(PetscBTDestroy(&bt));
1403:   PetscFunctionReturn(PETSC_SUCCESS);
1404: }

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

1418:   PetscFunctionBegin;
1419:   PetscCall(PetscIntMultError(2, nlnk_max + 2, &lsize));
1420:   PetscCall(PetscMalloc1(lsize, lnk));
1421:   llnk    = *lnk;
1422:   llnk[0] = 0;             /* number of entries on the list */
1423:   llnk[2] = PETSC_MAX_INT; /* value in the head node */
1424:   llnk[3] = 2;             /* next for the head node */
1425:   PetscFunctionReturn(PETSC_SUCCESS);
1426: }

1428: static inline PetscErrorCode PetscLLCondensedExpand_Scalable(PetscInt nlnk_max, PetscInt **lnk)
1429: {
1430:   PetscInt lsize = 0;

1432:   PetscFunctionBegin;
1433:   PetscCall(PetscIntMultError(2, nlnk_max + 2, &lsize));
1434:   PetscCall(PetscRealloc(lsize * sizeof(PetscInt), lnk));
1435:   PetscFunctionReturn(PETSC_SUCCESS);
1436: }

1438: static inline PetscErrorCode PetscLLCondensedAddSorted_Scalable(PetscInt nidx, const PetscInt indices[], PetscInt lnk[])
1439: {
1440:   PetscInt location = 2;      /* head */
1441:   PetscInt nlnk     = lnk[0]; /* num of entries on the input lnk */

1443:   for (PetscInt k = 0; k < nidx; k++) {
1444:     const PetscInt entry = indices[k];
1445:     PetscInt       next, lnkdata;

1447:     /* search for insertion location */
1448:     do {
1449:       next     = location + 1;  /* link from previous node to next node */
1450:       location = lnk[next];     /* idx of next node */
1451:       lnkdata  = lnk[location]; /* value of next node */
1452:     } while (entry > lnkdata);
1453:     if (entry < lnkdata) {
1454:       /* insertion location is found, add entry into lnk */
1455:       const PetscInt newnode = 2 * (nlnk + 2); /* index for this new node */
1456:       lnk[next]              = newnode;        /* connect previous node to the new node */
1457:       lnk[newnode]           = entry;          /* set value of the new node */
1458:       lnk[newnode + 1]       = location;       /* connect new node to next node */
1459:       location               = newnode;        /* next search starts from the new node */
1460:       nlnk++;
1461:     }
1462:   }
1463:   lnk[0] = nlnk; /* number of entries in the list */
1464:   return PETSC_SUCCESS;
1465: }

1467: static inline PetscErrorCode PetscLLCondensedClean_Scalable(PETSC_UNUSED PetscInt nidx, PetscInt *indices, PetscInt *lnk)
1468: {
1469:   const PetscInt nlnk = lnk[0];
1470:   PetscInt       next = lnk[3]; /* head node */

1472:   for (PetscInt k = 0; k < nlnk; k++) {
1473:     indices[k] = lnk[next];
1474:     next       = lnk[next + 1];
1475:   }
1476:   lnk[0] = 0; /* num of entries on the list */
1477:   lnk[3] = 2; /* head node */
1478:   return PETSC_SUCCESS;
1479: }

1481: static inline PetscErrorCode PetscLLCondensedDestroy_Scalable(PetscInt *lnk)
1482: {
1483:   return PetscFree(lnk);
1484: }

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

1494:       The next three are always the first link

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

1500:       The next three are always the last link

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

1507: static inline PetscErrorCode PetscLLCondensedCreate_fast(PetscInt nlnk_max, PetscInt **lnk)
1508: {
1509:   PetscInt *llnk;
1510:   PetscInt  lsize = 0;

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

1527: static inline PetscErrorCode PetscLLCondensedAddSorted_fast(PetscInt nidx, const PetscInt indices[], PetscInt lnk[])
1528: {
1529:   for (PetscInt k = 0, prev = 3 /* first value */; k < nidx; k++) {
1530:     const PetscInt entry = indices[k];
1531:     PetscInt       next  = lnk[prev + 2];

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

1571: static inline PetscErrorCode PetscLLCondensedClean_fast(PETSC_UNUSED PetscInt nidx, PetscInt *indices, PetscInt *lnk)
1572: {
1573:   const PetscInt nlnk = lnk[0];
1574:   PetscInt       next = lnk[5]; /* first node */

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

1591: static inline PetscErrorCode PetscLLCondensedView_fast(const PetscInt *lnk)
1592: {
1593:   const PetscInt nlnk = lnk[0];
1594:   PetscInt       next = lnk[5]; /* first node */

1596:   for (PetscInt k = 0; k < nlnk; k++) {
1597: #if 0 /* Debugging code */
1598:     printf("%d value %d len %d next %d\n", next, lnk[next], lnk[next + 1], lnk[next + 2]);
1599: #endif
1600:     next = lnk[next + 2];
1601:   }
1602:   return PETSC_SUCCESS;
1603: }

1605: static inline PetscErrorCode PetscLLCondensedDestroy_fast(PetscInt *lnk)
1606: {
1607:   return PetscFree(lnk);
1608: }

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

1613: PETSC_EXTERN PetscLogEvent MAT_Mult;
1614: PETSC_EXTERN PetscLogEvent MAT_MultAdd;
1615: PETSC_EXTERN PetscLogEvent MAT_MultTranspose;
1616: PETSC_EXTERN PetscLogEvent MAT_MultHermitianTranspose;
1617: PETSC_EXTERN PetscLogEvent MAT_MultTransposeAdd;
1618: PETSC_EXTERN PetscLogEvent MAT_MultHermitianTransposeAdd;
1619: PETSC_EXTERN PetscLogEvent MAT_Solve;
1620: PETSC_EXTERN PetscLogEvent MAT_Solves;
1621: PETSC_EXTERN PetscLogEvent MAT_SolveAdd;
1622: PETSC_EXTERN PetscLogEvent MAT_SolveTranspose;
1623: PETSC_EXTERN PetscLogEvent MAT_SolveTransposeAdd;
1624: PETSC_EXTERN PetscLogEvent MAT_SOR;
1625: PETSC_EXTERN PetscLogEvent MAT_ForwardSolve;
1626: PETSC_EXTERN PetscLogEvent MAT_BackwardSolve;
1627: PETSC_EXTERN PetscLogEvent MAT_LUFactor;
1628: PETSC_EXTERN PetscLogEvent MAT_LUFactorSymbolic;
1629: PETSC_EXTERN PetscLogEvent MAT_LUFactorNumeric;
1630: PETSC_EXTERN PetscLogEvent MAT_QRFactor;
1631: PETSC_EXTERN PetscLogEvent MAT_QRFactorSymbolic;
1632: PETSC_EXTERN PetscLogEvent MAT_QRFactorNumeric;
1633: PETSC_EXTERN PetscLogEvent MAT_CholeskyFactor;
1634: PETSC_EXTERN PetscLogEvent MAT_CholeskyFactorSymbolic;
1635: PETSC_EXTERN PetscLogEvent MAT_CholeskyFactorNumeric;
1636: PETSC_EXTERN PetscLogEvent MAT_ILUFactor;
1637: PETSC_EXTERN PetscLogEvent MAT_ILUFactorSymbolic;
1638: PETSC_EXTERN PetscLogEvent MAT_ICCFactorSymbolic;
1639: PETSC_EXTERN PetscLogEvent MAT_Copy;
1640: PETSC_EXTERN PetscLogEvent MAT_Convert;
1641: PETSC_EXTERN PetscLogEvent MAT_Scale;
1642: PETSC_EXTERN PetscLogEvent MAT_AssemblyBegin;
1643: PETSC_EXTERN PetscLogEvent MAT_AssemblyEnd;
1644: PETSC_EXTERN PetscLogEvent MAT_SetValues;
1645: PETSC_EXTERN PetscLogEvent MAT_GetValues;
1646: PETSC_EXTERN PetscLogEvent MAT_GetRow;
1647: PETSC_EXTERN PetscLogEvent MAT_GetRowIJ;
1648: PETSC_EXTERN PetscLogEvent MAT_CreateSubMats;
1649: PETSC_EXTERN PetscLogEvent MAT_GetOrdering;
1650: PETSC_EXTERN PetscLogEvent MAT_RedundantMat;
1651: PETSC_EXTERN PetscLogEvent MAT_IncreaseOverlap;
1652: PETSC_EXTERN PetscLogEvent MAT_Partitioning;
1653: PETSC_EXTERN PetscLogEvent MAT_PartitioningND;
1654: PETSC_EXTERN PetscLogEvent MAT_Coarsen;
1655: PETSC_EXTERN PetscLogEvent MAT_ZeroEntries;
1656: PETSC_EXTERN PetscLogEvent MAT_Load;
1657: PETSC_EXTERN PetscLogEvent MAT_View;
1658: PETSC_EXTERN PetscLogEvent MAT_AXPY;
1659: PETSC_EXTERN PetscLogEvent MAT_FDColoringCreate;
1660: PETSC_EXTERN PetscLogEvent MAT_TransposeColoringCreate;
1661: PETSC_EXTERN PetscLogEvent MAT_FDColoringSetUp;
1662: PETSC_EXTERN PetscLogEvent MAT_FDColoringApply;
1663: PETSC_EXTERN PetscLogEvent MAT_Transpose;
1664: PETSC_EXTERN PetscLogEvent MAT_FDColoringFunction;
1665: PETSC_EXTERN PetscLogEvent MAT_CreateSubMat;
1666: PETSC_EXTERN PetscLogEvent MAT_MatSolve;
1667: PETSC_EXTERN PetscLogEvent MAT_MatTrSolve;
1668: PETSC_EXTERN PetscLogEvent MAT_MatMultSymbolic;
1669: PETSC_EXTERN PetscLogEvent MAT_MatMultNumeric;
1670: PETSC_EXTERN PetscLogEvent MAT_Getlocalmatcondensed;
1671: PETSC_EXTERN PetscLogEvent MAT_GetBrowsOfAcols;
1672: PETSC_EXTERN PetscLogEvent MAT_GetBrowsOfAocols;
1673: PETSC_EXTERN PetscLogEvent MAT_PtAPSymbolic;
1674: PETSC_EXTERN PetscLogEvent MAT_PtAPNumeric;
1675: PETSC_EXTERN PetscLogEvent MAT_Seqstompinum;
1676: PETSC_EXTERN PetscLogEvent MAT_Seqstompisym;
1677: PETSC_EXTERN PetscLogEvent MAT_Seqstompi;
1678: PETSC_EXTERN PetscLogEvent MAT_Getlocalmat;
1679: PETSC_EXTERN PetscLogEvent MAT_RARtSymbolic;
1680: PETSC_EXTERN PetscLogEvent MAT_RARtNumeric;
1681: PETSC_EXTERN PetscLogEvent MAT_MatTransposeMultSymbolic;
1682: PETSC_EXTERN PetscLogEvent MAT_MatTransposeMultNumeric;
1683: PETSC_EXTERN PetscLogEvent MAT_TransposeMatMultSymbolic;
1684: PETSC_EXTERN PetscLogEvent MAT_TransposeMatMultNumeric;
1685: PETSC_EXTERN PetscLogEvent MAT_MatMatMultSymbolic;
1686: PETSC_EXTERN PetscLogEvent MAT_MatMatMultNumeric;
1687: PETSC_EXTERN PetscLogEvent MAT_Getsymtransreduced;
1688: PETSC_EXTERN PetscLogEvent MAT_GetSeqNonzeroStructure;
1689: PETSC_EXTERN PetscLogEvent MATMFFD_Mult;
1690: PETSC_EXTERN PetscLogEvent MAT_GetMultiProcBlock;
1691: PETSC_EXTERN PetscLogEvent MAT_CUSPARSECopyToGPU;
1692: PETSC_EXTERN PetscLogEvent MAT_CUSPARSECopyFromGPU;
1693: PETSC_EXTERN PetscLogEvent MAT_CUSPARSEGenerateTranspose;
1694: PETSC_EXTERN PetscLogEvent MAT_CUSPARSESolveAnalysis;
1695: PETSC_EXTERN PetscLogEvent MAT_HIPSPARSECopyToGPU;
1696: PETSC_EXTERN PetscLogEvent MAT_HIPSPARSECopyFromGPU;
1697: PETSC_EXTERN PetscLogEvent MAT_HIPSPARSEGenerateTranspose;
1698: PETSC_EXTERN PetscLogEvent MAT_HIPSPARSESolveAnalysis;
1699: PETSC_EXTERN PetscLogEvent MAT_SetValuesBatch;
1700: PETSC_EXTERN PetscLogEvent MAT_ViennaCLCopyToGPU;
1701: PETSC_EXTERN PetscLogEvent MAT_DenseCopyToGPU;
1702: PETSC_EXTERN PetscLogEvent MAT_DenseCopyFromGPU;
1703: PETSC_EXTERN PetscLogEvent MAT_Merge;
1704: PETSC_EXTERN PetscLogEvent MAT_Residual;
1705: PETSC_EXTERN PetscLogEvent MAT_SetRandom;
1706: PETSC_EXTERN PetscLogEvent MAT_FactorFactS;
1707: PETSC_EXTERN PetscLogEvent MAT_FactorInvS;
1708: PETSC_EXTERN PetscLogEvent MAT_PreallCOO;
1709: PETSC_EXTERN PetscLogEvent MAT_SetVCOO;
1710: PETSC_EXTERN PetscLogEvent MATCOLORING_Apply;
1711: PETSC_EXTERN PetscLogEvent MATCOLORING_Comm;
1712: PETSC_EXTERN PetscLogEvent MATCOLORING_Local;
1713: PETSC_EXTERN PetscLogEvent MATCOLORING_ISCreate;
1714: PETSC_EXTERN PetscLogEvent MATCOLORING_SetUp;
1715: PETSC_EXTERN PetscLogEvent MATCOLORING_Weights;
1716: PETSC_EXTERN PetscLogEvent MAT_H2Opus_Build;
1717: PETSC_EXTERN PetscLogEvent MAT_H2Opus_Compress;
1718: PETSC_EXTERN PetscLogEvent MAT_H2Opus_Orthog;
1719: PETSC_EXTERN PetscLogEvent MAT_H2Opus_LR;
1720: PETSC_EXTERN PetscLogEvent MAT_CUDACopyToGPU;