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 used outside Mat
238: */
239: PETSC_SINGLE_LIBRARY_INTERN PetscErrorCode MatFindNonzeroRowsOrCols_Basic(Mat, PetscBool, PetscReal, IS *);

241: /*
242:    Utility private matrix routines
243: */
244: PETSC_INTERN PetscErrorCode MatConvert_Basic(Mat, MatType, MatReuse, Mat *);
245: PETSC_INTERN PetscErrorCode MatConvert_Shell(Mat, MatType, MatReuse, Mat *);
246: PETSC_INTERN PetscErrorCode MatConvertFrom_Shell(Mat, MatType, MatReuse, Mat *);
247: PETSC_INTERN PetscErrorCode MatCopy_Basic(Mat, Mat, MatStructure);
248: PETSC_INTERN PetscErrorCode MatDiagonalSet_Default(Mat, Vec, InsertMode);
249: #if defined(PETSC_HAVE_SCALAPACK)
250: PETSC_INTERN PetscErrorCode MatConvert_Dense_ScaLAPACK(Mat, MatType, MatReuse, Mat *);
251: #endif
252: PETSC_INTERN PetscErrorCode MatSetPreallocationCOO_Basic(Mat, PetscCount, PetscInt[], PetscInt[]);
253: PETSC_INTERN PetscErrorCode MatSetValuesCOO_Basic(Mat, const PetscScalar[], InsertMode);

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

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

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

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

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

301: /*
302:   The stash is used to temporarily store inserted matrix values that
303:   belong to another processor. During the assembly phase the stashed
304:   values are moved to the correct processor and
305: */

307: typedef struct _MatStashSpace *PetscMatStashSpace;

309: struct _MatStashSpace {
310:   PetscMatStashSpace next;
311:   PetscScalar       *space_head, *val;
312:   PetscInt          *idx, *idy;
313:   PetscInt           total_space_size;
314:   PetscInt           local_used;
315:   PetscInt           local_remaining;
316: };

318: PETSC_EXTERN PetscErrorCode PetscMatStashSpaceGet(PetscInt, PetscInt, PetscMatStashSpace *);
319: PETSC_EXTERN PetscErrorCode PetscMatStashSpaceContiguous(PetscInt, PetscMatStashSpace *, PetscScalar *, PetscInt *, PetscInt *);
320: PETSC_EXTERN PetscErrorCode PetscMatStashSpaceDestroy(PetscMatStashSpace *);

322: typedef struct {
323:   PetscInt count;
324: } MatStashHeader;

326: typedef struct {
327:   void    *buffer; /* Of type blocktype, dynamically constructed  */
328:   PetscInt count;
329:   char     pending;
330: } MatStashFrame;

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

342:   PetscErrorCode (*ScatterBegin)(Mat, MatStash *, PetscInt *);
343:   PetscErrorCode (*ScatterGetMesg)(MatStash *, PetscMPIInt *, PetscInt **, PetscInt **, PetscScalar **, PetscInt *);
344:   PetscErrorCode (*ScatterEnd)(MatStash *);
345:   PetscErrorCode (*ScatterDestroy)(MatStash *);

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

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

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

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

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

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

434: typedef struct { /* used by MatProduct() */
435:   MatProductType type;
436:   char          *alg;
437:   Mat            A, B, C, Dwork;
438:   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, .. */
439:   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 .. */
440:   PetscBool      symbolic_used_the_fact_C_is_symmetric; /* MatMatMult(A,B,MAT_REUSE_MATRIX,..&C) is still legitimate), we need to redo symbolic! */
441:   PetscReal      fill;
442:   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 */

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

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

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

511: PETSC_INTERN PetscErrorCode MatAXPY_Basic(Mat, PetscScalar, Mat, MatStructure);
512: PETSC_INTERN PetscErrorCode MatAXPY_BasicWithPreallocation(Mat, Mat, PetscScalar, Mat, MatStructure);
513: PETSC_INTERN PetscErrorCode MatAXPY_Basic_Preallocate(Mat, Mat, Mat *);
514: PETSC_INTERN PetscErrorCode MatAXPY_Dense_Nest(Mat, PetscScalar, Mat);

516: PETSC_INTERN PetscErrorCode MatSetUp_Default(Mat);

518: /*
519:     Utility for MatZeroRows
520: */
521: PETSC_INTERN PetscErrorCode MatZeroRowsMapLocal_Private(Mat, PetscInt, const PetscInt *, PetscInt *, PetscInt **);

523: /*
524:     Utility for MatView/MatLoad
525: */
526: PETSC_INTERN PetscErrorCode MatView_Binary_BlockSizes(Mat, PetscViewer);
527: PETSC_INTERN PetscErrorCode MatLoad_Binary_BlockSizes(Mat, PetscViewer);

529: /*
530:     Object for partitioning graphs
531: */

533: typedef struct _MatPartitioningOps *MatPartitioningOps;
534: struct _MatPartitioningOps {
535:   PetscErrorCode (*apply)(MatPartitioning, IS *);
536:   PetscErrorCode (*applynd)(MatPartitioning, IS *);
537:   PetscErrorCode (*setfromoptions)(MatPartitioning, PetscOptionItems *);
538:   PetscErrorCode (*destroy)(MatPartitioning);
539:   PetscErrorCode (*view)(MatPartitioning, PetscViewer);
540:   PetscErrorCode (*improve)(MatPartitioning, IS *);
541: };

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

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

558: /*
559:     Object for coarsen graphs
560: */
561: typedef struct _MatCoarsenOps *MatCoarsenOps;
562: struct _MatCoarsenOps {
563:   PetscErrorCode (*apply)(MatCoarsen);
564:   PetscErrorCode (*setfromoptions)(MatCoarsen, PetscOptionItems *);
565:   PetscErrorCode (*destroy)(MatCoarsen);
566:   PetscErrorCode (*view)(MatCoarsen, PetscViewer);
567: };

569: struct _p_MatCoarsen {
570:   PETSCHEADER(struct _MatCoarsenOps);
571:   Mat   graph;
572:   void *subctx;
573:   /* */
574:   PetscBool         strict_aggs;
575:   IS                perm;
576:   PetscCoarsenData *agg_lists;
577: };

579: PETSC_EXTERN PetscErrorCode MatCoarsenMISKSetDistance(MatCoarsen, PetscInt);
580: PETSC_EXTERN PetscErrorCode MatCoarsenMISKGetDistance(MatCoarsen, PetscInt *);

582: /*
583:     Used in aijdevice.h
584: */
585: typedef struct {
586:   PetscInt    *i;
587:   PetscInt    *j;
588:   PetscScalar *a;
589:   PetscInt     n;
590:   PetscInt     ignorezeroentries;
591: } PetscCSRDataStructure;

593: /*
594:     MatFDColoring is used to compute Jacobian matrices efficiently
595:   via coloring. The data structure is explained below in an example.

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

609:     ncolors = 4;

611:     ncolumns      = {2,1,1,0}
612:     columns       = {{0,2},{1},{3},{}}
613:     nrows         = {4,2,3,3}
614:     rows          = {{0,1,2,3},{0,1},{1,2,3},{0,1,2}}
615:     vwscale       = {dx(0),dx(1),dx(2),dx(3)}               MPI Vec
616:     vscale        = {dx(0),dx(1),dx(2),dx(3),dx(4),dx(5)}   Seq Vec

618:     ncolumns      = {1,0,1,1}
619:     columns       = {{6},{},{4},{5}}
620:     nrows         = {3,0,2,2}
621:     rows          = {{0,1,2},{},{1,2},{1,2}}
622:     vwscale       = {dx(4),dx(5),dx(6)}              MPI Vec
623:     vscale        = {dx(0),dx(4),dx(5),dx(6)}        Seq Vec

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

628: */
629: typedef struct {
630:   PetscInt     row;
631:   PetscInt     col;
632:   PetscScalar *valaddr; /* address of value */
633: } MatEntry;

635: typedef struct {
636:   PetscInt     row;
637:   PetscScalar *valaddr; /* address of value */
638: } MatEntry2;

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

669: typedef struct _MatColoringOps *MatColoringOps;
670: struct _MatColoringOps {
671:   PetscErrorCode (*destroy)(MatColoring);
672:   PetscErrorCode (*setfromoptions)(MatColoring, PetscOptionItems *);
673:   PetscErrorCode (*view)(MatColoring, PetscViewer);
674:   PetscErrorCode (*apply)(MatColoring, ISColoring *);
675:   PetscErrorCode (*weights)(MatColoring, PetscReal **, PetscInt **);
676: };

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

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

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

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

722: /*
723:    Checking zero pivot for LU, ILU preconditioners.
724: */
725: typedef struct {
726:   PetscInt    nshift, nshift_max;
727:   PetscReal   shift_amount, shift_lo, shift_hi, shift_top, shift_fraction;
728:   PetscBool   newshift;
729:   PetscReal   rs; /* active row sum of abs(off-diagonals) */
730:   PetscScalar pv; /* pivot of the active row */
731: } FactorShiftCtx;

733: PETSC_EXTERN PetscErrorCode MatTransposeCheckNonzeroState_Private(Mat, Mat);

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

744: PETSC_EXTERN PetscErrorCode MatFactorDumpMatrix(Mat);
745: PETSC_INTERN PetscErrorCode MatShift_Basic(Mat, PetscScalar);
746: PETSC_INTERN PetscErrorCode MatSetBlockSizes_Default(Mat, PetscInt, PetscInt);

748: static inline PetscErrorCode MatPivotCheck_nz(PETSC_UNUSED Mat mat, const MatFactorInfo *info, FactorShiftCtx *sctx, PETSC_UNUSED PetscInt row)
749: {
750:   PetscReal _rs   = sctx->rs;
751:   PetscReal _zero = info->zeropivot * _rs;

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

766: static inline PetscErrorCode MatPivotCheck_pd(PETSC_UNUSED Mat mat, const MatFactorInfo *info, FactorShiftCtx *sctx, PETSC_UNUSED PetscInt row)
767: {
768:   PetscReal _rs   = sctx->rs;
769:   PetscReal _zero = info->zeropivot * _rs;

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

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

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

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

807:   PetscFunctionBegin;
808:   sctx->newshift = PETSC_FALSE;
809:   if (PetscAbsScalar(sctx->pv) <= _zero && !PetscIsNanScalar(sctx->pv)) {
810:     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);
811:     PetscCall(PetscInfo(mat, "Detected zero pivot in factorization in row %" PetscInt_FMT " value %g tolerance %g\n", row, (double)PetscAbsScalar(sctx->pv), (double)_zero));
812:     fact->factorerrortype             = MAT_FACTOR_NUMERIC_ZEROPIVOT;
813:     fact->factorerror_zeropivot_value = PetscAbsScalar(sctx->pv);
814:     fact->factorerror_zeropivot_row   = row;
815:   }
816:   PetscFunctionReturn(PETSC_SUCCESS);
817: }

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

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

843: #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)))

845: 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)
846: {
847:   PetscInt location;

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

865: 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)
866: {
867:   PetscFunctionBegin;
868:   *nlnk = 0;
869:   for (PetscInt k = 0, lnkdata = idx_start; k < nidx; ++k) {
870:     const PetscInt entry = indices[k];

872:     if (!PetscBTLookupSet(bt, entry)) PetscCall(PetscLLInsertLocation_Private(assume_sorted, k, idx_start, entry, nlnk, &lnkdata, lnk));
873:   }
874:   PetscFunctionReturn(PETSC_SUCCESS);
875: }

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

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

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

938:     if (!PetscBTLookupSet(bt, entry)) PetscCall(PetscLLInsertLocation_Private(PETSC_FALSE, k, idx_start, entry, nlnk, &lnkdata, lnk));
939:   }
940:   PetscFunctionReturn(PETSC_SUCCESS);
941: }

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

949:   PetscFunctionBegin;
950:   if (*lnk_empty) {
951:     for (PetscInt k = 0; k < nidx; ++k) {
952:       const PetscInt entry = indices[k], location = lnkdata;

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

976:       if (!PetscBTLookupSet(bt,entry)) PetscCall(PetscLLInsertLocation_Private(PETSC_TRUE,k,idx_start,entry,nlnk,&lnkdata,lnk));
977:     }
978:   }
979:   PetscFunctionReturn(PETSC_SUCCESS);
980: }
981: #endif

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

1006:   PetscFunctionBegin;
1007:   *nlnk = 0;
1008:   for (PetscInt k = 0, lnkdata = idx_start; k < nidx; ++k) {
1009:     const PetscInt entry = indices[k];

1011:     ++nzbd;
1012:     if (entry == diag) im[idx_start] = nzbd;
1013:     if (!PetscBTLookupSet(bt, entry)) PetscCall(PetscLLInsertLocation_Private(PETSC_TRUE, k, idx_start, entry, nlnk, &lnkdata, lnk));
1014:   }
1015:   PetscFunctionReturn(PETSC_SUCCESS);
1016: }

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

1043: /*
1044:   Free memories used by the list
1045: */
1046: #define PetscLLDestroy(lnk, bt) ((PetscErrorCode)(PetscFree(lnk) || PetscBTDestroy(&(bt))))

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

1063: 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)
1064: {
1065:   PetscFunctionBegin;
1066:   PetscCall(PetscLLInsertLocation_Private(assume_sorted, k, idx_start, entry, nlnk, lnkdata, lnk));
1067:   lnklvl[entry] = newval;
1068:   PetscFunctionReturn(PETSC_SUCCESS);
1069: }

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

1094:     if (!PetscBTLookupSet(bt, entry)) PetscCall(PetscIncompleteLLInsertLocation_Private(PETSC_FALSE, k, idx_start, entry, nlnk, &lnkdata, lnk, lnklvl, 0));
1095:   }
1096:   PetscFunctionReturn(PETSC_SUCCESS);
1097: }

1099: 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)
1100: {
1101:   PetscFunctionBegin;
1102:   *nlnk = 0;
1103:   for (PetscInt k = 0, lnkdata = idx_start; k < nidx; ++k) {
1104:     const PetscInt incrlev = idxlvl[k] + prow_offset + 1;

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

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

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

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

1161:   Note: the level of factor(i,j) is set as lvl(i,j) = min{ lvl(i,j), lvl(i,prow)+lvl(prow,j)+1)
1162:         where idx = non-zero columns of U(prow,prow+1:n-1), prow<i
1163: */
1164: 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)
1165: {
1166:   PetscFunctionBegin;
1167:   PetscCall(PetscIncompleteLLAdd_Private(nidx, idx, level, idxlvl, idx_start, nlnk, lnk, lnklvl, bt, lnklvl[prow], PETSC_TRUE));
1168:   PetscFunctionReturn(PETSC_SUCCESS);
1169: }

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

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

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

1248: /*
1249:   Free memories used by the list
1250: */
1251: #define PetscIncompleteLLDestroy(lnk, bt) ((PetscErrorCode)(PetscFree(lnk) || PetscBTDestroy(&(bt))))

1253: #if !defined(PETSC_CLANG_STATIC_ANALYZER)
1254:   #define MatCheckSameLocalSize(A, ar1, B, ar2) \
1255:     do { \
1256:       PetscCheckSameComm(A, ar1, B, ar2); \
1257:       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, \
1258:                  (A)->rmap->n, (A)->cmap->n, ar2, (B)->rmap->n, (B)->cmap->n); \
1259:     } while (0)
1260:   #define MatCheckSameSize(A, ar1, B, ar2) \
1261:     do { \
1262:       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, \
1263:                  (A)->rmap->N, (A)->cmap->N, ar2, (B)->rmap->N, (B)->cmap->N); \
1264:       MatCheckSameLocalSize(A, ar1, B, ar2); \
1265:     } while (0)
1266: #else
1267: template <typename Tm>
1268: extern void MatCheckSameLocalSize(Tm, int, Tm, int);
1269: template <typename Tm>
1270: extern void MatCheckSameSize(Tm, int, Tm, int);
1271: #endif

1273: #define VecCheckMatCompatible(M, x, ar1, b, ar2) \
1274:   do { \
1275:     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, \
1276:                (M)->cmap->N); \
1277:     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, \
1278:                (M)->rmap->N); \
1279:   } while (0)

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

1289:       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
1290:       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
1291:       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
1292:       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.
1293:       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
1294:       to each other so memory access is much better than using the big array.

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

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

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

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

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

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

1349:   PetscFunctionBegin;
1350:   for (PetscInt k = 0; k < nidx; k++) {
1351:     const PetscInt entry = indices[k];
1352:     if (!PetscBTLookupSet(bt, entry)) { /* new entry */
1353:       PetscInt next, lnkdata;

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

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

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

1391: static inline PetscErrorCode PetscLLCondensedView(PetscInt *lnk)
1392: {
1393:   PetscFunctionBegin;
1394:   PetscCall(PetscPrintf(PETSC_COMM_SELF, "LLCondensed of size %" PetscInt_FMT ", (val,  next)\n", lnk[0]));
1395:   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]));
1396:   PetscFunctionReturn(PETSC_SUCCESS);
1397: }

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

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

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

1432: static inline PetscErrorCode PetscLLCondensedExpand_Scalable(PetscInt nlnk_max, PetscInt **lnk)
1433: {
1434:   PetscInt lsize = 0;

1436:   PetscFunctionBegin;
1437:   PetscCall(PetscIntMultError(2, nlnk_max + 2, &lsize));
1438:   PetscCall(PetscRealloc(lsize * sizeof(PetscInt), lnk));
1439:   PetscFunctionReturn(PETSC_SUCCESS);
1440: }

1442: static inline PetscErrorCode PetscLLCondensedAddSorted_Scalable(PetscInt nidx, const PetscInt indices[], PetscInt lnk[])
1443: {
1444:   PetscInt location = 2;      /* head */
1445:   PetscInt nlnk     = lnk[0]; /* num of entries on the input lnk */

1447:   for (PetscInt k = 0; k < nidx; k++) {
1448:     const PetscInt entry = indices[k];
1449:     PetscInt       next, lnkdata;

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

1471: static inline PetscErrorCode PetscLLCondensedClean_Scalable(PETSC_UNUSED PetscInt nidx, PetscInt *indices, PetscInt *lnk)
1472: {
1473:   const PetscInt nlnk = lnk[0];
1474:   PetscInt       next = lnk[3]; /* head node */

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

1485: static inline PetscErrorCode PetscLLCondensedDestroy_Scalable(PetscInt *lnk)
1486: {
1487:   return PetscFree(lnk);
1488: }

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

1498:       The next three are always the first link

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

1504:       The next three are always the last link

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

1511: static inline PetscErrorCode PetscLLCondensedCreate_fast(PetscInt nlnk_max, PetscInt **lnk)
1512: {
1513:   PetscInt *llnk;
1514:   PetscInt  lsize = 0;

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

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

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

1575: static inline PetscErrorCode PetscLLCondensedClean_fast(PETSC_UNUSED PetscInt nidx, PetscInt *indices, PetscInt *lnk)
1576: {
1577:   const PetscInt nlnk = lnk[0];
1578:   PetscInt       next = lnk[5]; /* first node */

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

1595: static inline PetscErrorCode PetscLLCondensedView_fast(const PetscInt *lnk)
1596: {
1597:   const PetscInt nlnk = lnk[0];
1598:   PetscInt       next = lnk[5]; /* first node */

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

1609: static inline PetscErrorCode PetscLLCondensedDestroy_fast(PetscInt *lnk)
1610: {
1611:   return PetscFree(lnk);
1612: }

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

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