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, PetscInt, PetscInt[], Mat *);
216:   PetscErrorCode (*dummy)(Mat);
217:   /*150*/
218:   PetscErrorCode (*transposesymbolic)(Mat, Mat *);
219:   PetscErrorCode (*eliminatezeros)(Mat, PetscBool);
220:   PetscErrorCode (*getrowsumabs)(Mat, Vec);
221:   PetscErrorCode (*getfactor)(Mat, MatSolverType, MatFactorType, Mat *);
222:   PetscErrorCode (*getblockdiagonal)(Mat, Mat *);  // NOTE: the caller of get{block, vblock}diagonal owns the returned matrix;
223:   PetscErrorCode (*getvblockdiagonal)(Mat, Mat *); // they must destroy it after use
224: };
225: /*
226:     If you add MatOps entries above also add them to the MATOP enum
227:     in include/petscmat.h and src/mat/f90-mod/petscmat.h
228: */

230: #include <petscsys.h>

232: typedef struct _p_MatRootName *MatRootName;
233: struct _p_MatRootName {
234:   char       *rname, *sname, *mname;
235:   MatRootName next;
236: };

238: PETSC_EXTERN MatRootName MatRootNameList;

240: /*
241:    Utility private matrix routines used outside Mat
242: */
243: PETSC_SINGLE_LIBRARY_INTERN PetscErrorCode MatFindNonzeroRowsOrCols_Basic(Mat, PetscBool, PetscReal, IS *);
244: PETSC_EXTERN PetscErrorCode                MatShellGetScalingShifts(Mat, PetscScalar *, PetscScalar *, Vec *, Vec *, Vec *, Mat *, IS *, IS *);

246: #define MAT_SHELL_NOT_ALLOWED (void *)-1

248: /*
249:    Utility private matrix routines
250: */
251: PETSC_INTERN PetscErrorCode MatConvert_Basic(Mat, MatType, MatReuse, Mat *);
252: PETSC_INTERN PetscErrorCode MatConvert_Shell(Mat, MatType, MatReuse, Mat *);
253: PETSC_INTERN PetscErrorCode MatConvertFrom_Shell(Mat, MatType, MatReuse, Mat *);
254: PETSC_INTERN PetscErrorCode MatShellSetContext_Immutable(Mat X, void *ctx);
255: PETSC_INTERN PetscErrorCode MatShellSetContextDestroy_Immutable(Mat X, PetscErrorCode (*f)(void *));
256: PETSC_INTERN PetscErrorCode MatShellSetManageScalingShifts_Immutable(Mat X);
257: PETSC_INTERN PetscErrorCode MatCopy_Basic(Mat, Mat, MatStructure);
258: PETSC_INTERN PetscErrorCode MatDiagonalSet_Default(Mat, Vec, InsertMode);
259: #if defined(PETSC_HAVE_SCALAPACK)
260: PETSC_INTERN PetscErrorCode MatConvert_Dense_ScaLAPACK(Mat, MatType, MatReuse, Mat *);
261: #endif
262: PETSC_INTERN PetscErrorCode MatSetPreallocationCOO_Basic(Mat, PetscCount, PetscInt[], PetscInt[]);
263: PETSC_INTERN PetscErrorCode MatSetValuesCOO_Basic(Mat, const PetscScalar[], InsertMode);

265: /* This can be moved to the public header after implementing some missing MatProducts */
266: PETSC_INTERN PetscErrorCode MatCreateFromISLocalToGlobalMapping(ISLocalToGlobalMapping, Mat, PetscBool, PetscBool, MatType, Mat *);

268: /* these callbacks rely on the old matrix function pointers for
269:    matmat operations. They are unsafe, and should be removed.
270:    However, the amount of work needed to clean up all the
271:    implementations is not negligible */
272: PETSC_INTERN PetscErrorCode MatProductSymbolic_AB(Mat);
273: PETSC_INTERN PetscErrorCode MatProductNumeric_AB(Mat);
274: PETSC_INTERN PetscErrorCode MatProductSymbolic_AtB(Mat);
275: PETSC_INTERN PetscErrorCode MatProductNumeric_AtB(Mat);
276: PETSC_INTERN PetscErrorCode MatProductSymbolic_ABt(Mat);
277: PETSC_INTERN PetscErrorCode MatProductNumeric_ABt(Mat);
278: PETSC_INTERN PetscErrorCode MatProductNumeric_PtAP(Mat);
279: PETSC_INTERN PetscErrorCode MatProductNumeric_RARt(Mat);
280: PETSC_INTERN PetscErrorCode MatProductSymbolic_ABC(Mat);
281: PETSC_INTERN PetscErrorCode MatProductNumeric_ABC(Mat);

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

288: /* CreateGraph is common to AIJ seq and mpi */
289: PETSC_INTERN PetscErrorCode MatCreateGraph_Simple_AIJ(Mat, PetscBool, PetscBool, PetscReal, PetscInt, PetscInt[], Mat *);

291: #if defined(PETSC_CLANG_STATIC_ANALYZER)
292: template <typename Tm>
293: extern void MatCheckPreallocated(Tm, int);
294: template <typename Tm>
295: extern void MatCheckProduct(Tm, int);
296: #else /* PETSC_CLANG_STATIC_ANALYZER */
297:   #define MatCheckPreallocated(A, arg) \
298:     do { \
299:       if (!(A)->preallocated) PetscCall(MatSetUp(A)); \
300:     } while (0)

302:   #if defined(PETSC_USE_DEBUG)
303:     #define MatCheckProduct(A, arg) \
304:       do { \
305:         PetscCheck((A)->product, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Argument %d \"%s\" is not a matrix obtained from MatProductCreate()", (arg), #A); \
306:       } while (0)
307:   #else
308:     #define MatCheckProduct(A, arg) \
309:       do { \
310:       } while (0)
311:   #endif
312: #endif /* PETSC_CLANG_STATIC_ANALYZER */

314: /*
315:   The stash is used to temporarily store inserted matrix values that
316:   belong to another processor. During the assembly phase the stashed
317:   values are moved to the correct processor and
318: */

320: typedef struct _MatStashSpace *PetscMatStashSpace;

322: struct _MatStashSpace {
323:   PetscMatStashSpace next;
324:   PetscScalar       *space_head, *val;
325:   PetscInt          *idx, *idy;
326:   PetscInt           total_space_size;
327:   PetscInt           local_used;
328:   PetscInt           local_remaining;
329: };

331: PETSC_EXTERN PetscErrorCode PetscMatStashSpaceGet(PetscInt, PetscInt, PetscMatStashSpace *);
332: PETSC_EXTERN PetscErrorCode PetscMatStashSpaceContiguous(PetscInt, PetscMatStashSpace *, PetscScalar *, PetscInt *, PetscInt *);
333: PETSC_EXTERN PetscErrorCode PetscMatStashSpaceDestroy(PetscMatStashSpace *);

335: typedef struct {
336:   PetscInt count;
337: } MatStashHeader;

339: typedef struct {
340:   void    *buffer; /* Of type blocktype, dynamically constructed  */
341:   PetscInt count;
342:   char     pending;
343: } MatStashFrame;

345: typedef struct _MatStash MatStash;
346: struct _MatStash {
347:   PetscInt           nmax;              /* maximum stash size */
348:   PetscInt           umax;              /* user specified max-size */
349:   PetscInt           oldnmax;           /* the nmax value used previously */
350:   PetscInt           n;                 /* stash size */
351:   PetscInt           bs;                /* block size of the stash */
352:   PetscInt           reallocs;          /* preserve the no of mallocs invoked */
353:   PetscMatStashSpace space_head, space; /* linked list to hold stashed global row/column numbers and matrix values */

355:   PetscErrorCode (*ScatterBegin)(Mat, MatStash *, PetscInt *);
356:   PetscErrorCode (*ScatterGetMesg)(MatStash *, PetscMPIInt *, PetscInt **, PetscInt **, PetscScalar **, PetscInt *);
357:   PetscErrorCode (*ScatterEnd)(MatStash *);
358:   PetscErrorCode (*ScatterDestroy)(MatStash *);

360:   /* The following variables are used for communication */
361:   MPI_Comm      comm;
362:   PetscMPIInt   size, rank;
363:   PetscMPIInt   tag1, tag2;
364:   MPI_Request  *send_waits;     /* array of send requests */
365:   MPI_Request  *recv_waits;     /* array of receive requests */
366:   MPI_Status   *send_status;    /* array of send status */
367:   PetscMPIInt   nsends, nrecvs; /* numbers of sends and receives */
368:   PetscScalar  *svalues;        /* sending data */
369:   PetscInt     *sindices;
370:   PetscScalar **rvalues;    /* receiving data (values) */
371:   PetscInt    **rindices;   /* receiving data (indices) */
372:   PetscMPIInt   nprocessed; /* number of messages already processed */
373:   PetscMPIInt  *flg_v;      /* indicates what messages have arrived so far and from whom */
374:   PetscBool     reproduce;
375:   PetscMPIInt   reproduce_count;

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

405: #if !defined(PETSC_HAVE_MPIUNI)
406: PETSC_INTERN PetscErrorCode MatStashScatterDestroy_BTS(MatStash *);
407: #endif
408: PETSC_INTERN PetscErrorCode MatStashCreate_Private(MPI_Comm, PetscInt, MatStash *);
409: PETSC_INTERN PetscErrorCode MatStashDestroy_Private(MatStash *);
410: PETSC_INTERN PetscErrorCode MatStashScatterEnd_Private(MatStash *);
411: PETSC_INTERN PetscErrorCode MatStashSetInitialSize_Private(MatStash *, PetscInt);
412: PETSC_INTERN PetscErrorCode MatStashGetInfo_Private(MatStash *, PetscInt *, PetscInt *);
413: PETSC_INTERN PetscErrorCode MatStashValuesRow_Private(MatStash *, PetscInt, PetscInt, const PetscInt[], const PetscScalar[], PetscBool);
414: PETSC_INTERN PetscErrorCode MatStashValuesCol_Private(MatStash *, PetscInt, PetscInt, const PetscInt[], const PetscScalar[], PetscInt, PetscBool);
415: PETSC_INTERN PetscErrorCode MatStashValuesRowBlocked_Private(MatStash *, PetscInt, PetscInt, const PetscInt[], const PetscScalar[], PetscInt, PetscInt, PetscInt);
416: PETSC_INTERN PetscErrorCode MatStashValuesColBlocked_Private(MatStash *, PetscInt, PetscInt, const PetscInt[], const PetscScalar[], PetscInt, PetscInt, PetscInt);
417: PETSC_INTERN PetscErrorCode MatStashScatterBegin_Private(Mat, MatStash *, PetscInt *);
418: PETSC_INTERN PetscErrorCode MatStashScatterGetMesg_Private(MatStash *, PetscMPIInt *, PetscInt **, PetscInt **, PetscScalar **, PetscInt *);
419: PETSC_INTERN PetscErrorCode MatGetInfo_External(Mat, MatInfoType, MatInfo *);

421: typedef struct {
422:   PetscInt  dim;
423:   PetscInt  dims[4];
424:   PetscInt  starts[4];
425:   PetscBool noc; /* this is a single component problem, hence user will not set MatStencil.c */
426: } MatStencilInfo;

428: /* Info about using compressed row format */
429: typedef struct {
430:   PetscBool use;    /* indicates compressed rows have been checked and will be used */
431:   PetscInt  nrows;  /* number of non-zero rows */
432:   PetscInt *i;      /* compressed row pointer  */
433:   PetscInt *rindex; /* compressed row index               */
434: } Mat_CompressedRow;
435: PETSC_EXTERN PetscErrorCode MatCheckCompressedRow(Mat, PetscInt, Mat_CompressedRow *, PetscInt *, PetscInt, PetscReal);

437: typedef struct { /* used by MatCreateRedundantMatrix() for reusing matredundant */
438:   PetscInt     nzlocal, nsends, nrecvs;
439:   PetscMPIInt *send_rank, *recv_rank;
440:   PetscInt    *sbuf_nz, *rbuf_nz, *sbuf_j, **rbuf_j;
441:   PetscScalar *sbuf_a, **rbuf_a;
442:   MPI_Comm     subcomm; /* when user does not provide a subcomm */
443:   IS           isrow, iscol;
444:   Mat         *matseq;
445: } Mat_Redundant;

447: typedef struct { /* used by MatProduct() */
448:   MatProductType type;
449:   char          *alg;
450:   Mat            A, B, C, Dwork;
451:   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, .. */
452:   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 .. */
453:   PetscBool      symbolic_used_the_fact_C_is_symmetric; /* MatMatMult(A,B,MAT_REUSE_MATRIX,..&C) is still legitimate), we need to redo symbolic! */
454:   PetscObjectParameterDeclare(PetscReal, fill);
455:   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 */
456:   PetscBool setfromoptionscalled;

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

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

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

525: PETSC_INTERN PetscErrorCode MatAXPY_Basic(Mat, PetscScalar, Mat, MatStructure);
526: PETSC_INTERN PetscErrorCode MatAXPY_BasicWithPreallocation(Mat, Mat, PetscScalar, Mat, MatStructure);
527: PETSC_INTERN PetscErrorCode MatAXPY_Basic_Preallocate(Mat, Mat, Mat *);
528: PETSC_INTERN PetscErrorCode MatAXPY_Dense_Nest(Mat, PetscScalar, Mat);

530: PETSC_INTERN PetscErrorCode MatSetUp_Default(Mat);

532: /*
533:     Utility for MatZeroRows
534: */
535: PETSC_INTERN PetscErrorCode MatZeroRowsMapLocal_Private(Mat, PetscInt, const PetscInt *, PetscInt *, PetscInt **);

537: /*
538:     Utility for MatView/MatLoad
539: */
540: PETSC_INTERN PetscErrorCode MatView_Binary_BlockSizes(Mat, PetscViewer);
541: PETSC_INTERN PetscErrorCode MatLoad_Binary_BlockSizes(Mat, PetscViewer);

543: /*
544:     Object for partitioning graphs
545: */

547: typedef struct _MatPartitioningOps *MatPartitioningOps;
548: struct _MatPartitioningOps {
549:   PetscErrorCode (*apply)(MatPartitioning, IS *);
550:   PetscErrorCode (*applynd)(MatPartitioning, IS *);
551:   PetscErrorCode (*setfromoptions)(MatPartitioning, PetscOptionItems *);
552:   PetscErrorCode (*destroy)(MatPartitioning);
553:   PetscErrorCode (*view)(MatPartitioning, PetscViewer);
554:   PetscErrorCode (*improve)(MatPartitioning, IS *);
555: };

557: struct _p_MatPartitioning {
558:   PETSCHEADER(struct _MatPartitioningOps);
559:   Mat        adj;
560:   PetscInt  *vertex_weights;
561:   PetscReal *part_weights;
562:   PetscInt   n;    /* number of partitions */
563:   PetscInt   ncon; /* number of vertex weights per vertex */
564:   void      *data;
565:   PetscInt   setupcalled;
566:   PetscBool  use_edge_weights; /* A flag indicates whether or not to use edge weights */
567: };

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

572: /*
573:     Object for coarsen graphs
574: */
575: typedef struct _MatCoarsenOps *MatCoarsenOps;
576: struct _MatCoarsenOps {
577:   PetscErrorCode (*apply)(MatCoarsen);
578:   PetscErrorCode (*setfromoptions)(MatCoarsen, PetscOptionItems *);
579:   PetscErrorCode (*destroy)(MatCoarsen);
580:   PetscErrorCode (*view)(MatCoarsen, PetscViewer);
581: };

583: #define MAT_COARSEN_STRENGTH_INDEX_SIZE 3
584: struct _p_MatCoarsen {
585:   PETSCHEADER(struct _MatCoarsenOps);
586:   Mat   graph;
587:   void *subctx;
588:   /* */
589:   PetscBool         strict_aggs;
590:   IS                perm;
591:   PetscCoarsenData *agg_lists;
592:   PetscInt          max_it;    /* number of iterations in HEM */
593:   PetscReal         threshold; /* HEM can filter interim graphs */
594:   PetscInt          strength_index_size;
595:   PetscInt          strength_index[MAT_COARSEN_STRENGTH_INDEX_SIZE];
596: };

598: PETSC_EXTERN PetscErrorCode MatCoarsenMISKSetDistance(MatCoarsen, PetscInt);
599: PETSC_EXTERN PetscErrorCode MatCoarsenMISKGetDistance(MatCoarsen, PetscInt *);

601: /*
602:     Used in aijdevice.h
603: */
604: typedef struct {
605:   PetscInt    *i;
606:   PetscInt    *j;
607:   PetscScalar *a;
608:   PetscInt     n;
609:   PetscInt     ignorezeroentries;
610: } PetscCSRDataStructure;

612: /*
613:     MatFDColoring is used to compute Jacobian matrices efficiently
614:   via coloring. The data structure is explained below in an example.

616:    Color =   0    1     0    2   |   2      3       0
617:    ---------------------------------------------------
618:             00   01              |          05
619:             10   11              |   14     15               Processor  0
620:                        22    23  |          25
621:                        32    33  |
622:    ===================================================
623:                                  |   44     45     46
624:             50                   |          55               Processor 1
625:                                  |   64            66
626:    ---------------------------------------------------

628:     ncolors = 4;

630:     ncolumns      = {2,1,1,0}
631:     columns       = {{0,2},{1},{3},{}}
632:     nrows         = {4,2,3,3}
633:     rows          = {{0,1,2,3},{0,1},{1,2,3},{0,1,2}}
634:     vwscale       = {dx(0),dx(1),dx(2),dx(3)}               MPI Vec
635:     vscale        = {dx(0),dx(1),dx(2),dx(3),dx(4),dx(5)}   Seq Vec

637:     ncolumns      = {1,0,1,1}
638:     columns       = {{6},{},{4},{5}}
639:     nrows         = {3,0,2,2}
640:     rows          = {{0,1,2},{},{1,2},{1,2}}
641:     vwscale       = {dx(4),dx(5),dx(6)}              MPI Vec
642:     vscale        = {dx(0),dx(4),dx(5),dx(6)}        Seq Vec

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

647: */
648: typedef struct {
649:   PetscInt     row;
650:   PetscInt     col;
651:   PetscScalar *valaddr; /* address of value */
652: } MatEntry;

654: typedef struct {
655:   PetscInt     row;
656:   PetscScalar *valaddr; /* address of value */
657: } MatEntry2;

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

688: typedef struct _MatColoringOps *MatColoringOps;
689: struct _MatColoringOps {
690:   PetscErrorCode (*destroy)(MatColoring);
691:   PetscErrorCode (*setfromoptions)(MatColoring, PetscOptionItems *);
692:   PetscErrorCode (*view)(MatColoring, PetscViewer);
693:   PetscErrorCode (*apply)(MatColoring, ISColoring *);
694:   PetscErrorCode (*weights)(MatColoring, PetscReal **, PetscInt **);
695: };

697: struct _p_MatColoring {
698:   PETSCHEADER(struct _MatColoringOps);
699:   Mat                   mat;
700:   PetscInt              dist;         /* distance of the coloring */
701:   PetscInt              maxcolors;    /* the maximum number of colors returned, maxcolors=1 for MIS */
702:   void                 *data;         /* inner context */
703:   PetscBool             valid;        /* check to see if what is produced is a valid coloring */
704:   MatColoringWeightType weight_type;  /* type of weight computation to be performed */
705:   PetscReal            *user_weights; /* custom weights and permutation */
706:   PetscInt             *user_lperm;
707:   PetscBool             valid_iscoloring; /* check to see if matcoloring is produced a valid iscoloring */
708: };

710: struct _p_MatTransposeColoring {
711:   PETSCHEADER(int);
712:   PetscInt       M, N, m;      /* total rows, columns; local rows */
713:   PetscInt       rstart;       /* first row owned by local processor */
714:   PetscInt       ncolors;      /* number of colors */
715:   PetscInt      *ncolumns;     /* number of local columns for a color */
716:   PetscInt      *nrows;        /* number of local rows for each color */
717:   PetscInt       currentcolor; /* color for which function evaluation is being done now */
718:   ISColoringType ctype;        /* IS_COLORING_GLOBAL or IS_COLORING_LOCAL */

720:   PetscInt *colorforrow, *colorforcol; /* pointer to rows and columns */
721:   PetscInt *rows;                      /* lists the local rows for each color (using the local row numbering) */
722:   PetscInt *den2sp;                    /* maps (row,color) in the dense matrix to index of sparse matrix array a->a */
723:   PetscInt *columns;                   /* lists the local columns of each color (using global column numbering) */
724:   PetscInt  brows;                     /* number of rows for efficient implementation of MatTransColoringApplyDenToSp() */
725:   PetscInt *lstart;                    /* array used for loop over row blocks of Csparse */
726: };

728: /*
729:    Null space context for preconditioner/operators
730: */
731: struct _p_MatNullSpace {
732:   PETSCHEADER(int);
733:   PetscBool    has_cnst;
734:   PetscInt     n;
735:   Vec         *vecs;
736:   PetscScalar *alpha;                                  /* for projections */
737:   PetscErrorCode (*remove)(MatNullSpace, Vec, void *); /* for user provided removal function */
738:   void *rmctx;                                         /* context for remove() function */
739: };

741: /*
742:    Checking zero pivot for LU, ILU preconditioners.
743: */
744: typedef struct {
745:   PetscInt    nshift, nshift_max;
746:   PetscReal   shift_amount, shift_lo, shift_hi, shift_top, shift_fraction;
747:   PetscBool   newshift;
748:   PetscReal   rs; /* active row sum of abs(off-diagonals) */
749:   PetscScalar pv; /* pivot of the active row */
750: } FactorShiftCtx;

752: PETSC_EXTERN PetscErrorCode MatTransposeCheckNonzeroState_Private(Mat, Mat);

754: /*
755:  Used by MatTranspose() and potentially other functions to track the matrix used in the generation of another matrix
756: */
757: typedef struct {
758:   PetscObjectId    id;
759:   PetscObjectState state;
760:   PetscObjectState nonzerostate;
761: } MatParentState;

763: PETSC_EXTERN PetscErrorCode MatFactorDumpMatrix(Mat);
764: PETSC_INTERN PetscErrorCode MatSetBlockSizes_Default(Mat, PetscInt, PetscInt);

766: PETSC_SINGLE_LIBRARY_INTERN PetscErrorCode MatShift_Basic(Mat, PetscScalar);

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

773:   PetscFunctionBegin;
774:   if (PetscAbsScalar(sctx->pv) <= _zero && !PetscIsNanScalar(sctx->pv)) {
775:     /* force |diag| > zeropivot*rs */
776:     if (!sctx->nshift) sctx->shift_amount = info->shiftamount;
777:     else sctx->shift_amount *= 2.0;
778:     sctx->newshift = PETSC_TRUE;
779:     (sctx->nshift)++;
780:   } else {
781:     sctx->newshift = PETSC_FALSE;
782:   }
783:   PetscFunctionReturn(PETSC_SUCCESS);
784: }

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

791:   PetscFunctionBegin;
792:   if (PetscRealPart(sctx->pv) <= _zero && !PetscIsNanScalar(sctx->pv)) {
793:     /* force matfactor to be diagonally dominant */
794:     if (sctx->nshift == sctx->nshift_max) {
795:       sctx->shift_fraction = sctx->shift_hi;
796:     } else {
797:       sctx->shift_lo       = sctx->shift_fraction;
798:       sctx->shift_fraction = (sctx->shift_hi + sctx->shift_lo) / (PetscReal)2.;
799:     }
800:     sctx->shift_amount = sctx->shift_fraction * sctx->shift_top;
801:     sctx->nshift++;
802:     sctx->newshift = PETSC_TRUE;
803:   } else {
804:     sctx->newshift = PETSC_FALSE;
805:   }
806:   PetscFunctionReturn(PETSC_SUCCESS);
807: }

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

813:   PetscFunctionBegin;
814:   if (PetscAbsScalar(sctx->pv) <= _zero && !PetscIsNanScalar(sctx->pv)) {
815:     sctx->pv += info->shiftamount;
816:     sctx->shift_amount = 0.0;
817:     sctx->nshift++;
818:   }
819:   sctx->newshift = PETSC_FALSE;
820:   PetscFunctionReturn(PETSC_SUCCESS);
821: }

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

827:   PetscFunctionBegin;
828:   sctx->newshift = PETSC_FALSE;
829:   if (PetscAbsScalar(sctx->pv) <= _zero && !PetscIsNanScalar(sctx->pv)) {
830:     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);
831:     PetscCall(PetscInfo(mat, "Detected zero pivot in factorization in row %" PetscInt_FMT " value %g tolerance %g\n", row, (double)PetscAbsScalar(sctx->pv), (double)_zero));
832:     fact->factorerrortype             = MAT_FACTOR_NUMERIC_ZEROPIVOT;
833:     fact->factorerror_zeropivot_value = PetscAbsScalar(sctx->pv);
834:     fact->factorerror_zeropivot_row   = row;
835:   }
836:   PetscFunctionReturn(PETSC_SUCCESS);
837: }

839: static inline PetscErrorCode MatPivotCheck(Mat fact, Mat mat, const MatFactorInfo *info, FactorShiftCtx *sctx, PetscInt row)
840: {
841:   PetscFunctionBegin;
842:   if (info->shifttype == (PetscReal)MAT_SHIFT_NONZERO) PetscCall(MatPivotCheck_nz(mat, info, sctx, row));
843:   else if (info->shifttype == (PetscReal)MAT_SHIFT_POSITIVE_DEFINITE) PetscCall(MatPivotCheck_pd(mat, info, sctx, row));
844:   else if (info->shifttype == (PetscReal)MAT_SHIFT_INBLOCKS) PetscCall(MatPivotCheck_inblocks(mat, info, sctx, row));
845:   else PetscCall(MatPivotCheck_none(fact, mat, info, sctx, row));
846:   PetscFunctionReturn(PETSC_SUCCESS);
847: }

849: #include <petscbt.h>
850: /*
851:   Create and initialize a linked list
852:   Input Parameters:
853:     idx_start - starting index of the list
854:     lnk_max   - max value of lnk indicating the end of the list
855:     nlnk      - max length of the list
856:   Output Parameters:
857:     lnk       - list initialized
858:     bt        - PetscBT (bitarray) with all bits set to false
859:     lnk_empty - flg indicating the list is empty
860: */
861: #define PetscLLCreate(idx_start, lnk_max, nlnk, lnk, bt) ((PetscErrorCode)(PetscMalloc1(nlnk, &lnk) || PetscBTCreate(nlnk, &(bt)) || (lnk[idx_start] = lnk_max, PETSC_SUCCESS)))

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

865: 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)
866: {
867:   PetscInt location;

869:   PetscFunctionBegin;
870:   /* start from the beginning if entry < previous entry */
871:   if (!assume_sorted && k && entry < *lnkdata) *lnkdata = idx_start;
872:   /* search for insertion location */
873:   do {
874:     location = *lnkdata;
875:     *lnkdata = lnk[location];
876:   } while (entry > *lnkdata);
877:   /* insertion location is found, add entry into lnk */
878:   lnk[location] = entry;
879:   lnk[entry]    = *lnkdata;
880:   ++(*nlnk);
881:   *lnkdata = entry; /* next search starts from here if next_entry > entry */
882:   PetscFunctionReturn(PETSC_SUCCESS);
883: }

885: 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)
886: {
887:   PetscFunctionBegin;
888:   *nlnk = 0;
889:   for (PetscInt k = 0, lnkdata = idx_start; k < nidx; ++k) {
890:     const PetscInt entry = indices[k];

892:     if (!PetscBTLookupSet(bt, entry)) PetscCall(PetscLLInsertLocation_Private(assume_sorted, k, idx_start, entry, nlnk, &lnkdata, lnk));
893:   }
894:   PetscFunctionReturn(PETSC_SUCCESS);
895: }

897: /*
898:   Add an index set into a sorted linked list
899:   Input Parameters:
900:     nidx      - number of input indices
901:     indices   - 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 PetscLLAdd(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_FALSE));
914:   PetscFunctionReturn(PETSC_SUCCESS);
915: }

917: /*
918:   Add a SORTED ascending index set into a sorted linked list - same as PetscLLAdd() bus skip 'if (_k && _entry < _lnkdata) _lnkdata  = idx_start;'
919:   Input Parameters:
920:     nidx      - number of input indices
921:     indices   - sorted integer array
922:     idx_start - starting index of the list
923:     lnk       - linked list(an integer array) that is created
924:     bt        - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
925:   output Parameters:
926:     nlnk      - number of newly added indices
927:     lnk       - the sorted(increasing order) linked list containing new and non-redundate entries from indices
928:     bt        - updated PetscBT (bitarray)
929: */
930: static inline PetscErrorCode PetscLLAddSorted(PetscInt nidx, const PetscInt *PETSC_RESTRICT indices, PetscInt idx_start, PetscInt *PETSC_RESTRICT nlnk, PetscInt *PETSC_RESTRICT lnk, PetscBT bt)
931: {
932:   PetscFunctionBegin;
933:   PetscCall(PetscLLAdd_Private(nidx, indices, idx_start, nlnk, lnk, bt, PETSC_TRUE));
934:   PetscFunctionReturn(PETSC_SUCCESS);
935: }

937: /*
938:   Add a permuted index set into a sorted linked list
939:   Input Parameters:
940:     nidx      - number of input indices
941:     indices   - integer array
942:     perm      - permutation of indices
943:     idx_start - starting index of the list
944:     lnk       - linked list(an integer array) that is created
945:     bt        - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
946:   output Parameters:
947:     nlnk      - number of newly added indices
948:     lnk       - the sorted(increasing order) linked list containing new and non-redundate entries from indices
949:     bt        - updated PetscBT (bitarray)
950: */
951: 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)
952: {
953:   PetscFunctionBegin;
954:   *nlnk = 0;
955:   for (PetscInt k = 0, lnkdata = idx_start; k < nidx; ++k) {
956:     const PetscInt entry = perm[indices[k]];

958:     if (!PetscBTLookupSet(bt, entry)) PetscCall(PetscLLInsertLocation_Private(PETSC_FALSE, k, idx_start, entry, nlnk, &lnkdata, lnk));
959:   }
960:   PetscFunctionReturn(PETSC_SUCCESS);
961: }

963: #if 0
964: /* this appears to be unused? */
965: static inline PetscErrorCode PetscLLAddSorted_new(PetscInt nidx, PetscInt *indices, PetscInt idx_start, PetscBool *lnk_empty, PetscInt *nlnk, PetscInt *lnk, PetscBT bt)
966: {
967:   PetscInt lnkdata = idx_start;

969:   PetscFunctionBegin;
970:   if (*lnk_empty) {
971:     for (PetscInt k = 0; k < nidx; ++k) {
972:       const PetscInt entry = indices[k], location = lnkdata;

974:       PetscCall(PetscBTSet(bt,entry)); /* mark the new entry */
975:       lnkdata       = lnk[location];
976:       /* insertion location is found, add entry into lnk */
977:       lnk[location] = entry;
978:       lnk[entry]    = lnkdata;
979:       lnkdata       = entry; /* next search starts from here */
980:     }
981:     /* lnk[indices[nidx-1]] = lnk[idx_start];
982:        lnk[idx_start]       = indices[0];
983:        PetscCall(PetscBTSet(bt,indices[0]));
984:        for (_k=1; _k<nidx; _k++) {
985:        PetscCall(PetscBTSet(bt,indices[_k]));
986:        lnk[indices[_k-1]] = indices[_k];
987:        }
988:     */
989:     *nlnk      = nidx;
990:     *lnk_empty = PETSC_FALSE;
991:   } else {
992:     *nlnk = 0;
993:     for (PetscInt k = 0; k < nidx; ++k) {
994:       const PetscInt entry = indices[k];

996:       if (!PetscBTLookupSet(bt,entry)) PetscCall(PetscLLInsertLocation_Private(PETSC_TRUE,k,idx_start,entry,nlnk,&lnkdata,lnk));
997:     }
998:   }
999:   PetscFunctionReturn(PETSC_SUCCESS);
1000: }
1001: #endif

1003: /*
1004:   Add a SORTED index set into a sorted linked list used for LUFactorSymbolic()
1005:   Same as PetscLLAddSorted() with an additional operation:
1006:        count the number of input indices that are no larger than 'diag'
1007:   Input Parameters:
1008:     indices   - sorted integer array
1009:     idx_start - starting index of the list, index of pivot row
1010:     lnk       - linked list(an integer array) that is created
1011:     bt        - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
1012:     diag      - index of the active row in LUFactorSymbolic
1013:     nzbd      - number of input indices with indices <= idx_start
1014:     im        - im[idx_start] is initialized as num of nonzero entries in row=idx_start
1015:   output Parameters:
1016:     nlnk      - number of newly added indices
1017:     lnk       - the sorted(increasing order) linked list containing new and non-redundate entries from indices
1018:     bt        - updated PetscBT (bitarray)
1019:     im        - im[idx_start]: unchanged if diag is not an entry
1020:                              : num of entries with indices <= diag if diag is an entry
1021: */
1022: 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)
1023: {
1024:   const PetscInt nidx = im[idx_start] - nzbd; /* num of entries with idx_start < index <= diag */

1026:   PetscFunctionBegin;
1027:   *nlnk = 0;
1028:   for (PetscInt k = 0, lnkdata = idx_start; k < nidx; ++k) {
1029:     const PetscInt entry = indices[k];

1031:     ++nzbd;
1032:     if (entry == diag) im[idx_start] = nzbd;
1033:     if (!PetscBTLookupSet(bt, entry)) PetscCall(PetscLLInsertLocation_Private(PETSC_TRUE, k, idx_start, entry, nlnk, &lnkdata, lnk));
1034:   }
1035:   PetscFunctionReturn(PETSC_SUCCESS);
1036: }

1038: /*
1039:   Copy data on the list into an array, then initialize the list
1040:   Input Parameters:
1041:     idx_start - starting index of the list
1042:     lnk_max   - max value of lnk indicating the end of the list
1043:     nlnk      - number of data on the list to be copied
1044:     lnk       - linked list
1045:     bt        - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
1046:   output Parameters:
1047:     indices   - array that contains the copied data
1048:     lnk       - linked list that is cleaned and initialize
1049:     bt        - PetscBT (bitarray) with all bits set to false
1050: */
1051: static inline PetscErrorCode PetscLLClean(PetscInt idx_start, PetscInt lnk_max, PetscInt nlnk, PetscInt *PETSC_RESTRICT lnk, PetscInt *PETSC_RESTRICT indices, PetscBT bt)
1052: {
1053:   PetscFunctionBegin;
1054:   for (PetscInt j = 0, idx = idx_start; j < nlnk; ++j) {
1055:     idx        = lnk[idx];
1056:     indices[j] = idx;
1057:     PetscCall(PetscBTClear(bt, idx));
1058:   }
1059:   lnk[idx_start] = lnk_max;
1060:   PetscFunctionReturn(PETSC_SUCCESS);
1061: }

1063: /*
1064:   Free memories used by the list
1065: */
1066: #define PetscLLDestroy(lnk, bt) ((PetscErrorCode)(PetscFree(lnk) || PetscBTDestroy(&(bt))))

1068: /* Routines below are used for incomplete matrix factorization */
1069: /*
1070:   Create and initialize a linked list and its levels
1071:   Input Parameters:
1072:     idx_start - starting index of the list
1073:     lnk_max   - max value of lnk indicating the end of the list
1074:     nlnk      - max length of the list
1075:   Output Parameters:
1076:     lnk       - list initialized
1077:     lnk_lvl   - array of size nlnk for storing levels of lnk
1078:     bt        - PetscBT (bitarray) with all bits set to false
1079: */
1080: #define PetscIncompleteLLCreate(idx_start, lnk_max, nlnk, lnk, lnk_lvl, bt) \
1081:   ((PetscErrorCode)(PetscIntMultError(2, nlnk, NULL) || PetscMalloc1(2 * nlnk, &lnk) || PetscBTCreate(nlnk, &(bt)) || (lnk[idx_start] = lnk_max, lnk_lvl = lnk + nlnk, PETSC_SUCCESS)))

1083: 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)
1084: {
1085:   PetscFunctionBegin;
1086:   PetscCall(PetscLLInsertLocation_Private(assume_sorted, k, idx_start, entry, nlnk, lnkdata, lnk));
1087:   lnklvl[entry] = newval;
1088:   PetscFunctionReturn(PETSC_SUCCESS);
1089: }

1091: /*
1092:   Initialize a sorted linked list used for ILU and ICC
1093:   Input Parameters:
1094:     nidx      - number of input idx
1095:     idx       - integer array used for storing column indices
1096:     idx_start - starting index of the list
1097:     perm      - indices of an IS
1098:     lnk       - linked list(an integer array) that is created
1099:     lnklvl    - levels of lnk
1100:     bt        - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
1101:   output Parameters:
1102:     nlnk     - number of newly added idx
1103:     lnk      - the sorted(increasing order) linked list containing new and non-redundate entries from idx
1104:     lnklvl   - levels of lnk
1105:     bt       - updated PetscBT (bitarray)
1106: */
1107: 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)
1108: {
1109:   PetscFunctionBegin;
1110:   *nlnk = 0;
1111:   for (PetscInt k = 0, lnkdata = idx_start; k < nidx; ++k) {
1112:     const PetscInt entry = perm[idx[k]];

1114:     if (!PetscBTLookupSet(bt, entry)) PetscCall(PetscIncompleteLLInsertLocation_Private(PETSC_FALSE, k, idx_start, entry, nlnk, &lnkdata, lnk, lnklvl, 0));
1115:   }
1116:   PetscFunctionReturn(PETSC_SUCCESS);
1117: }

1119: 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)
1120: {
1121:   PetscFunctionBegin;
1122:   *nlnk = 0;
1123:   for (PetscInt k = 0, lnkdata = idx_start; k < nidx; ++k) {
1124:     const PetscInt incrlev = idxlvl[k] + prow_offset + 1;

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

1129:       if (!PetscBTLookupSet(bt, entry)) PetscCall(PetscIncompleteLLInsertLocation_Private(assume_sorted, k, idx_start, entry, nlnk, &lnkdata, lnk, lnklvl, incrlev));
1130:       else if (lnklvl[entry] > incrlev) lnklvl[entry] = incrlev; /* existing entry */
1131:     }
1132:   }
1133:   PetscFunctionReturn(PETSC_SUCCESS);
1134: }

1136: /*
1137:   Add a SORTED index set into a sorted linked list for ICC
1138:   Input Parameters:
1139:     nidx      - number of input indices
1140:     idx       - sorted integer array used for storing column indices
1141:     level     - level of fill, e.g., ICC(level)
1142:     idxlvl    - level of idx
1143:     idx_start - starting index of the list
1144:     lnk       - linked list(an integer array) that is created
1145:     lnklvl    - levels of lnk
1146:     bt        - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
1147:     idxlvl_prow - idxlvl[prow], where prow is the row number of the idx
1148:   output Parameters:
1149:     nlnk   - number of newly added indices
1150:     lnk    - the sorted(increasing order) linked list containing new and non-redundate entries from idx
1151:     lnklvl - levels of lnk
1152:     bt     - updated PetscBT (bitarray)
1153:   Note: the level of U(i,j) is set as lvl(i,j) = min{ lvl(i,j), lvl(prow,i)+lvl(prow,j)+1)
1154:         where idx = non-zero columns of U(prow,prow+1:n-1), prow<i
1155: */
1156: 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)
1157: {
1158:   PetscFunctionBegin;
1159:   PetscCall(PetscIncompleteLLAdd_Private(nidx, idx, level, idxlvl, idx_start, nlnk, lnk, lnklvl, bt, idxlvl_prow, PETSC_TRUE));
1160:   PetscFunctionReturn(PETSC_SUCCESS);
1161: }

1163: /*
1164:   Add a SORTED index set into a sorted linked list for ILU
1165:   Input Parameters:
1166:     nidx      - number of input indices
1167:     idx       - sorted integer array used for storing column indices
1168:     level     - level of fill, e.g., ICC(level)
1169:     idxlvl    - level of idx
1170:     idx_start - starting index of the list
1171:     lnk       - linked list(an integer array) that is created
1172:     lnklvl    - levels of lnk
1173:     bt        - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
1174:     prow      - the row number of idx
1175:   output Parameters:
1176:     nlnk     - number of newly added idx
1177:     lnk      - the sorted(increasing order) linked list containing new and non-redundate entries from idx
1178:     lnklvl   - levels of lnk
1179:     bt       - updated PetscBT (bitarray)

1181:   Note: the level of factor(i,j) is set as lvl(i,j) = min{ lvl(i,j), lvl(i,prow)+lvl(prow,j)+1)
1182:         where idx = non-zero columns of U(prow,prow+1:n-1), prow<i
1183: */
1184: 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)
1185: {
1186:   PetscFunctionBegin;
1187:   PetscCall(PetscIncompleteLLAdd_Private(nidx, idx, level, idxlvl, idx_start, nlnk, lnk, lnklvl, bt, lnklvl[prow], PETSC_TRUE));
1188:   PetscFunctionReturn(PETSC_SUCCESS);
1189: }

1191: /*
1192:   Add a index set into a sorted linked list
1193:   Input Parameters:
1194:     nidx      - number of input idx
1195:     idx   - 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 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)
1209: {
1210:   PetscFunctionBegin;
1211:   PetscCall(PetscIncompleteLLAdd_Private(nidx, idx, level, idxlvl, idx_start, nlnk, lnk, lnklvl, bt, 0, PETSC_FALSE));
1212:   PetscFunctionReturn(PETSC_SUCCESS);
1213: }

1215: /*
1216:   Add a SORTED index set into a sorted linked list
1217:   Input Parameters:
1218:     nidx      - number of input indices
1219:     idx   - sorted integer array used for storing column indices
1220:     level     - level of fill, e.g., ICC(level)
1221:     idxlvl - level of idx
1222:     idx_start - starting index of the list
1223:     lnk       - linked list(an integer array) that is created
1224:     lnklvl    - levels of lnk
1225:     bt        - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
1226:   output Parameters:
1227:     nlnk      - number of newly added idx
1228:     lnk       - the sorted(increasing order) linked list containing new and non-redundate entries from idx
1229:     lnklvl    - levels of lnk
1230:     bt        - updated PetscBT (bitarray)
1231: */
1232: 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)
1233: {
1234:   PetscFunctionBegin;
1235:   PetscCall(PetscIncompleteLLAdd_Private(nidx, idx, level, idxlvl, idx_start, nlnk, lnk, lnklvl, bt, 0, PETSC_TRUE));
1236:   PetscFunctionReturn(PETSC_SUCCESS);
1237: }

1239: /*
1240:   Copy data on the list into an array, then initialize the list
1241:   Input Parameters:
1242:     idx_start - starting index of the list
1243:     lnk_max   - max value of lnk indicating the end of the list
1244:     nlnk      - number of data on the list to be copied
1245:     lnk       - linked list
1246:     lnklvl    - level of lnk
1247:     bt        - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
1248:   output Parameters:
1249:     indices - array that contains the copied data
1250:     lnk     - linked list that is cleaned and initialize
1251:     lnklvl  - level of lnk that is reinitialized
1252:     bt      - PetscBT (bitarray) with all bits set to false
1253: */
1254: 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)
1255: {
1256:   PetscFunctionBegin;
1257:   for (PetscInt j = 0, idx = idx_start; j < nlnk; ++j) {
1258:     idx           = lnk[idx];
1259:     indices[j]    = idx;
1260:     indiceslvl[j] = lnklvl[idx];
1261:     lnklvl[idx]   = -1;
1262:     PetscCall(PetscBTClear(bt, idx));
1263:   }
1264:   lnk[idx_start] = lnk_max;
1265:   PetscFunctionReturn(PETSC_SUCCESS);
1266: }

1268: /*
1269:   Free memories used by the list
1270: */
1271: #define PetscIncompleteLLDestroy(lnk, bt) ((PetscErrorCode)(PetscFree(lnk) || PetscBTDestroy(&(bt))))

1273: #if !defined(PETSC_CLANG_STATIC_ANALYZER)
1274:   #define MatCheckSameLocalSize(A, ar1, B, ar2) \
1275:     do { \
1276:       PetscCheckSameComm(A, ar1, B, ar2); \
1277:       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, \
1278:                  (A)->rmap->n, (A)->cmap->n, ar2, (B)->rmap->n, (B)->cmap->n); \
1279:     } while (0)
1280:   #define MatCheckSameSize(A, ar1, B, ar2) \
1281:     do { \
1282:       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, \
1283:                  (A)->rmap->N, (A)->cmap->N, ar2, (B)->rmap->N, (B)->cmap->N); \
1284:       MatCheckSameLocalSize(A, ar1, B, ar2); \
1285:     } while (0)
1286: #else
1287: template <typename Tm>
1288: extern void MatCheckSameLocalSize(Tm, int, Tm, int);
1289: template <typename Tm>
1290: extern void MatCheckSameSize(Tm, int, Tm, int);
1291: #endif

1293: #define VecCheckMatCompatible(M, x, ar1, b, ar2) \
1294:   do { \
1295:     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, \
1296:                (M)->cmap->N); \
1297:     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, \
1298:                (M)->rmap->N); \
1299:   } while (0)

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

1309:       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
1310:       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
1311:       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
1312:       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.
1313:       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
1314:       to each other so memory access is much better than using the big array.

1316:   Example:
1317:      nlnk_max=5, lnk_max=36:
1318:      Initial list: [0, 0 | 36, 2 | 0, 0 | 0, 0 | 0, 0 | 0, 0 | 0, 0]
1319:      here, head_node has index 2 with value lnk[2]=lnk_max=36,
1320:            0-th entry is used to store the number of entries in the list,
1321:      The initial lnk represents head -> tail(marked by 36) with number of entries = lnk[0]=0.

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

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

1331:   Input Parameters:
1332:     nlnk_max  - max length of the list
1333:     lnk_max   - max value of the entries
1334:   Output Parameters:
1335:     lnk       - list created and initialized
1336:     bt        - PetscBT (bitarray) with all bits set to false. Note: bt has size lnk_max, not nln_max!
1337: */
1338: static inline PetscErrorCode PetscLLCondensedCreate(PetscInt nlnk_max, PetscInt lnk_max, PetscInt **lnk, PetscBT *bt)
1339: {
1340:   PetscInt *llnk, lsize = 0;

1342:   PetscFunctionBegin;
1343:   PetscCall(PetscIntMultError(2, nlnk_max + 2, &lsize));
1344:   PetscCall(PetscMalloc1(lsize, lnk));
1345:   PetscCall(PetscBTCreate(lnk_max, bt));
1346:   llnk    = *lnk;
1347:   llnk[0] = 0;       /* number of entries on the list */
1348:   llnk[2] = lnk_max; /* value in the head node */
1349:   llnk[3] = 2;       /* next for the head node */
1350:   PetscFunctionReturn(PETSC_SUCCESS);
1351: }

1353: /*
1354:   Add a SORTED ascending index set into a sorted linked list. See PetscLLCondensedCreate() for detailed description.
1355:   Input Parameters:
1356:     nidx      - number of input indices
1357:     indices   - sorted integer array
1358:     lnk       - condensed linked list(an integer array) that is created
1359:     bt        - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
1360:   output Parameters:
1361:     lnk       - the sorted(increasing order) linked list containing previous and newly added non-redundate indices
1362:     bt        - updated PetscBT (bitarray)
1363: */
1364: static inline PetscErrorCode PetscLLCondensedAddSorted(PetscInt nidx, const PetscInt indices[], PetscInt lnk[], PetscBT bt)
1365: {
1366:   PetscInt location = 2;      /* head */
1367:   PetscInt nlnk     = lnk[0]; /* num of entries on the input lnk */

1369:   PetscFunctionBegin;
1370:   for (PetscInt k = 0; k < nidx; k++) {
1371:     const PetscInt entry = indices[k];
1372:     if (!PetscBTLookupSet(bt, entry)) { /* new entry */
1373:       PetscInt next, lnkdata;

1375:       /* search for insertion location */
1376:       do {
1377:         next     = location + 1;  /* link from previous node to next node */
1378:         location = lnk[next];     /* idx of next node */
1379:         lnkdata  = lnk[location]; /* value of next node */
1380:       } while (entry > lnkdata);
1381:       /* insertion location is found, add entry into lnk */
1382:       const PetscInt newnode = 2 * (nlnk + 2); /* index for this new node */
1383:       lnk[next]              = newnode;        /* connect previous node to the new node */
1384:       lnk[newnode]           = entry;          /* set value of the new node */
1385:       lnk[newnode + 1]       = location;       /* connect new node to next node */
1386:       location               = newnode;        /* next search starts from the new node */
1387:       nlnk++;
1388:     }
1389:   }
1390:   lnk[0] = nlnk; /* number of entries in the list */
1391:   PetscFunctionReturn(PETSC_SUCCESS);
1392: }

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

1399:   PetscFunctionBegin;
1400:   for (PetscInt k = 0; k < nlnk; k++) {
1401:     indices[k] = lnk[next];
1402:     next       = lnk[next + 1];
1403:     PetscCall(PetscBTClear(bt, indices[k]));
1404:   }
1405:   lnk[0] = 0;       /* num of entries on the list */
1406:   lnk[2] = lnk_max; /* initialize head node */
1407:   lnk[3] = 2;       /* head node */
1408:   PetscFunctionReturn(PETSC_SUCCESS);
1409: }

1411: static inline PetscErrorCode PetscLLCondensedView(PetscInt *lnk)
1412: {
1413:   PetscFunctionBegin;
1414:   PetscCall(PetscPrintf(PETSC_COMM_SELF, "LLCondensed of size %" PetscInt_FMT ", (val,  next)\n", lnk[0]));
1415:   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]));
1416:   PetscFunctionReturn(PETSC_SUCCESS);
1417: }

1419: /*
1420:   Free memories used by the list
1421: */
1422: static inline PetscErrorCode PetscLLCondensedDestroy(PetscInt *lnk, PetscBT bt)
1423: {
1424:   PetscFunctionBegin;
1425:   PetscCall(PetscFree(lnk));
1426:   PetscCall(PetscBTDestroy(&bt));
1427:   PetscFunctionReturn(PETSC_SUCCESS);
1428: }

1430: /* -------------------------------------------------------------------------------------------------------*/
1431: /*
1432:  Same as PetscLLCondensedCreate(), but does not use non-scalable O(lnk_max) bitarray
1433:   Input Parameters:
1434:     nlnk_max  - max length of the list
1435:   Output Parameters:
1436:     lnk       - list created and initialized
1437: */
1438: static inline PetscErrorCode PetscLLCondensedCreate_Scalable(PetscInt nlnk_max, PetscInt **lnk)
1439: {
1440:   PetscInt *llnk, lsize = 0;

1442:   PetscFunctionBegin;
1443:   PetscCall(PetscIntMultError(2, nlnk_max + 2, &lsize));
1444:   PetscCall(PetscMalloc1(lsize, lnk));
1445:   llnk    = *lnk;
1446:   llnk[0] = 0;             /* number of entries on the list */
1447:   llnk[2] = PETSC_INT_MAX; /* value in the head node */
1448:   llnk[3] = 2;             /* next for the head node */
1449:   PetscFunctionReturn(PETSC_SUCCESS);
1450: }

1452: static inline PetscErrorCode PetscLLCondensedExpand_Scalable(PetscInt nlnk_max, PetscInt **lnk)
1453: {
1454:   PetscInt lsize = 0;

1456:   PetscFunctionBegin;
1457:   PetscCall(PetscIntMultError(2, nlnk_max + 2, &lsize));
1458:   PetscCall(PetscRealloc(lsize * sizeof(PetscInt), lnk));
1459:   PetscFunctionReturn(PETSC_SUCCESS);
1460: }

1462: static inline PetscErrorCode PetscLLCondensedAddSorted_Scalable(PetscInt nidx, const PetscInt indices[], PetscInt lnk[])
1463: {
1464:   PetscInt location = 2;      /* head */
1465:   PetscInt nlnk     = lnk[0]; /* num of entries on the input lnk */

1467:   for (PetscInt k = 0; k < nidx; k++) {
1468:     const PetscInt entry = indices[k];
1469:     PetscInt       next, lnkdata;

1471:     /* search for insertion location */
1472:     do {
1473:       next     = location + 1;  /* link from previous node to next node */
1474:       location = lnk[next];     /* idx of next node */
1475:       lnkdata  = lnk[location]; /* value of next node */
1476:     } while (entry > lnkdata);
1477:     if (entry < lnkdata) {
1478:       /* insertion location is found, add entry into lnk */
1479:       const PetscInt newnode = 2 * (nlnk + 2); /* index for this new node */
1480:       lnk[next]              = newnode;        /* connect previous node to the new node */
1481:       lnk[newnode]           = entry;          /* set value of the new node */
1482:       lnk[newnode + 1]       = location;       /* connect new node to next node */
1483:       location               = newnode;        /* next search starts from the new node */
1484:       nlnk++;
1485:     }
1486:   }
1487:   lnk[0] = nlnk; /* number of entries in the list */
1488:   return PETSC_SUCCESS;
1489: }

1491: static inline PetscErrorCode PetscLLCondensedClean_Scalable(PETSC_UNUSED PetscInt nidx, PetscInt *indices, PetscInt *lnk)
1492: {
1493:   const PetscInt nlnk = lnk[0];
1494:   PetscInt       next = lnk[3]; /* head node */

1496:   for (PetscInt k = 0; k < nlnk; k++) {
1497:     indices[k] = lnk[next];
1498:     next       = lnk[next + 1];
1499:   }
1500:   lnk[0] = 0; /* num of entries on the list */
1501:   lnk[3] = 2; /* head node */
1502:   return PETSC_SUCCESS;
1503: }

1505: static inline PetscErrorCode PetscLLCondensedDestroy_Scalable(PetscInt *lnk)
1506: {
1507:   return PetscFree(lnk);
1508: }

1510: /* -------------------------------------------------------------------------------------------------------*/
1511: /*
1512:       lnk[0]   number of links
1513:       lnk[1]   number of entries
1514:       lnk[3n]  value
1515:       lnk[3n+1] len
1516:       lnk[3n+2] link to next value

1518:       The next three are always the first link

1520:       lnk[3]    PETSC_INT_MIN+1
1521:       lnk[4]    1
1522:       lnk[5]    link to first real entry

1524:       The next three are always the last link

1526:       lnk[6]    PETSC_INT_MAX - 1
1527:       lnk[7]    1
1528:       lnk[8]    next valid link (this is the same as lnk[0] but without the decreases)
1529: */

1531: static inline PetscErrorCode PetscLLCondensedCreate_fast(PetscInt nlnk_max, PetscInt **lnk)
1532: {
1533:   PetscInt *llnk;
1534:   PetscInt  lsize = 0;

1536:   PetscFunctionBegin;
1537:   PetscCall(PetscIntMultError(3, nlnk_max + 3, &lsize));
1538:   PetscCall(PetscMalloc1(lsize, lnk));
1539:   llnk    = *lnk;
1540:   llnk[0] = 0;                 /* nlnk: number of entries on the list */
1541:   llnk[1] = 0;                 /* number of integer entries represented in list */
1542:   llnk[3] = PETSC_INT_MIN + 1; /* value in the first node */
1543:   llnk[4] = 1;                 /* count for the first node */
1544:   llnk[5] = 6;                 /* next for the first node */
1545:   llnk[6] = PETSC_INT_MAX - 1; /* value in the last node */
1546:   llnk[7] = 1;                 /* count for the last node */
1547:   llnk[8] = 0;                 /* next valid node to be used */
1548:   PetscFunctionReturn(PETSC_SUCCESS);
1549: }

1551: static inline PetscErrorCode PetscLLCondensedAddSorted_fast(PetscInt nidx, const PetscInt indices[], PetscInt lnk[])
1552: {
1553:   for (PetscInt k = 0, prev = 3 /* first value */; k < nidx; k++) {
1554:     const PetscInt entry = indices[k];
1555:     PetscInt       next  = lnk[prev + 2];

1557:     /* search for insertion location */
1558:     while (entry >= lnk[next]) {
1559:       prev = next;
1560:       next = lnk[next + 2];
1561:     }
1562:     /* entry is in range of previous list */
1563:     if (entry < lnk[prev] + lnk[prev + 1]) continue;
1564:     lnk[1]++;
1565:     /* entry is right after previous list */
1566:     if (entry == lnk[prev] + lnk[prev + 1]) {
1567:       lnk[prev + 1]++;
1568:       if (lnk[next] == entry + 1) { /* combine two contiguous strings */
1569:         lnk[prev + 1] += lnk[next + 1];
1570:         lnk[prev + 2] = lnk[next + 2];
1571:         next          = lnk[next + 2];
1572:         lnk[0]--;
1573:       }
1574:       continue;
1575:     }
1576:     /* entry is right before next list */
1577:     if (entry == lnk[next] - 1) {
1578:       lnk[next]--;
1579:       lnk[next + 1]++;
1580:       prev = next;
1581:       next = lnk[prev + 2];
1582:       continue;
1583:     }
1584:     /*  add entry into lnk */
1585:     lnk[prev + 2] = 3 * ((lnk[8]++) + 3); /* connect previous node to the new node */
1586:     prev          = lnk[prev + 2];
1587:     lnk[prev]     = entry; /* set value of the new node */
1588:     lnk[prev + 1] = 1;     /* number of values in contiguous string is one to start */
1589:     lnk[prev + 2] = next;  /* connect new node to next node */
1590:     lnk[0]++;
1591:   }
1592:   return PETSC_SUCCESS;
1593: }

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

1600:   for (PetscInt k = 0, cnt = 0; k < nlnk; k++) {
1601:     for (PetscInt j = 0; j < lnk[next + 1]; j++) indices[cnt++] = lnk[next] + j;
1602:     next = lnk[next + 2];
1603:   }
1604:   lnk[0] = 0;                 /* nlnk: number of links */
1605:   lnk[1] = 0;                 /* number of integer entries represented in list */
1606:   lnk[3] = PETSC_INT_MIN + 1; /* value in the first node */
1607:   lnk[4] = 1;                 /* count for the first node */
1608:   lnk[5] = 6;                 /* next for the first node */
1609:   lnk[6] = PETSC_INT_MAX - 1; /* value in the last node */
1610:   lnk[7] = 1;                 /* count for the last node */
1611:   lnk[8] = 0;                 /* next valid location to make link */
1612:   return PETSC_SUCCESS;
1613: }

1615: static inline PetscErrorCode PetscLLCondensedView_fast(const PetscInt *lnk)
1616: {
1617:   const PetscInt nlnk = lnk[0];
1618:   PetscInt       next = lnk[5]; /* first node */

1620:   for (PetscInt k = 0; k < nlnk; k++) {
1621: #if 0 /* Debugging code */
1622:     printf("%d value %d len %d next %d\n", next, lnk[next], lnk[next + 1], lnk[next + 2]);
1623: #endif
1624:     next = lnk[next + 2];
1625:   }
1626:   return PETSC_SUCCESS;
1627: }

1629: static inline PetscErrorCode PetscLLCondensedDestroy_fast(PetscInt *lnk)
1630: {
1631:   return PetscFree(lnk);
1632: }

1634: PETSC_EXTERN PetscErrorCode PetscCDCreate(PetscInt, PetscCoarsenData **);
1635: PETSC_EXTERN PetscErrorCode PetscCDDestroy(PetscCoarsenData *);
1636: PETSC_EXTERN PetscErrorCode PetscCDIntNdSetID(PetscCDIntNd *, PetscInt);
1637: PETSC_EXTERN PetscErrorCode PetscCDIntNdGetID(const PetscCDIntNd *, PetscInt *);
1638: PETSC_EXTERN PetscErrorCode PetscCDAppendID(PetscCoarsenData *, PetscInt, PetscInt);
1639: PETSC_EXTERN PetscErrorCode PetscCDMoveAppend(PetscCoarsenData *, PetscInt, PetscInt);
1640: PETSC_EXTERN PetscErrorCode PetscCDAppendNode(PetscCoarsenData *, PetscInt, PetscCDIntNd *);
1641: PETSC_EXTERN PetscErrorCode PetscCDRemoveNextNode(PetscCoarsenData *, PetscInt, PetscCDIntNd *);
1642: PETSC_EXTERN PetscErrorCode PetscCDCountAt(const PetscCoarsenData *, PetscInt, PetscInt *);
1643: PETSC_EXTERN PetscErrorCode PetscCDIsEmptyAt(const PetscCoarsenData *, PetscInt, PetscBool *);
1644: PETSC_EXTERN PetscErrorCode PetscCDSetChunkSize(PetscCoarsenData *, PetscInt);
1645: PETSC_EXTERN PetscErrorCode PetscCDPrint(const PetscCoarsenData *, PetscInt, MPI_Comm);
1646: PETSC_EXTERN PetscErrorCode PetscCDGetNonemptyIS(PetscCoarsenData *, IS *);
1647: PETSC_EXTERN PetscErrorCode PetscCDGetMat(PetscCoarsenData *, Mat *);
1648: PETSC_EXTERN PetscErrorCode PetscCDSetMat(PetscCoarsenData *, Mat);
1649: PETSC_EXTERN PetscErrorCode PetscCDClearMat(PetscCoarsenData *);
1650: PETSC_EXTERN PetscErrorCode PetscCDRemoveAllAt(PetscCoarsenData *, PetscInt);
1651: PETSC_EXTERN PetscErrorCode PetscCDCount(const PetscCoarsenData *, PetscInt *_sz);

1653: PETSC_EXTERN PetscErrorCode PetscCDGetHeadPos(const PetscCoarsenData *, PetscInt, PetscCDIntNd **);
1654: PETSC_EXTERN PetscErrorCode PetscCDGetNextPos(const PetscCoarsenData *, PetscInt, PetscCDIntNd **);
1655: PETSC_EXTERN PetscErrorCode PetscCDGetASMBlocks(const PetscCoarsenData *, const PetscInt, PetscInt *, IS **);

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

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