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: };
222: /*
223:     If you add MatOps entries above also add them to the MATOP enum
224:     in include/petscmat.h and src/mat/f90-mod/petscmat.h
225: */

227: #include <petscsys.h>

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

235: PETSC_EXTERN MatRootName MatRootNameList;

237: /*
238:    Utility private matrix routines used outside Mat
239: */
240: PETSC_SINGLE_LIBRARY_INTERN PetscErrorCode MatFindNonzeroRowsOrCols_Basic(Mat, PetscBool, PetscReal, IS *);
241: PETSC_SINGLE_LIBRARY_INTERN PetscErrorCode MatShellGetScalingShifts(Mat, PetscScalar *, PetscScalar *, Vec *, Vec *, Vec *, Mat *, IS *, IS *);

243: #define MAT_SHELL_NOT_ALLOWED (void *)-1

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

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

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

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

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

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

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

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

317: typedef struct _MatStashSpace *PetscMatStashSpace;

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

328: PETSC_EXTERN PetscErrorCode PetscMatStashSpaceGet(PetscInt, PetscInt, PetscMatStashSpace *);
329: PETSC_EXTERN PetscErrorCode PetscMatStashSpaceContiguous(PetscInt, PetscMatStashSpace *, PetscScalar *, PetscInt *, PetscInt *);
330: PETSC_EXTERN PetscErrorCode PetscMatStashSpaceDestroy(PetscMatStashSpace *);

332: typedef struct {
333:   PetscInt count;
334: } MatStashHeader;

336: typedef struct {
337:   void    *buffer; /* Of type blocktype, dynamically constructed  */
338:   PetscInt count;
339:   char     pending;
340: } MatStashFrame;

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

352:   PetscErrorCode (*ScatterBegin)(Mat, MatStash *, PetscInt *);
353:   PetscErrorCode (*ScatterGetMesg)(MatStash *, PetscMPIInt *, PetscInt **, PetscInt **, PetscScalar **, PetscInt *);
354:   PetscErrorCode (*ScatterEnd)(MatStash *);
355:   PetscErrorCode (*ScatterDestroy)(MatStash *);

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

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

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

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

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

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

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

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

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

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

522: PETSC_INTERN PetscErrorCode MatAXPY_Basic(Mat, PetscScalar, Mat, MatStructure);
523: PETSC_INTERN PetscErrorCode MatAXPY_BasicWithPreallocation(Mat, Mat, PetscScalar, Mat, MatStructure);
524: PETSC_INTERN PetscErrorCode MatAXPY_Basic_Preallocate(Mat, Mat, Mat *);
525: PETSC_INTERN PetscErrorCode MatAXPY_Dense_Nest(Mat, PetscScalar, Mat);

527: PETSC_INTERN PetscErrorCode MatSetUp_Default(Mat);

529: /*
530:     Utility for MatZeroRows
531: */
532: PETSC_INTERN PetscErrorCode MatZeroRowsMapLocal_Private(Mat, PetscInt, const PetscInt *, PetscInt *, PetscInt **);

534: /*
535:     Utility for MatView/MatLoad
536: */
537: PETSC_INTERN PetscErrorCode MatView_Binary_BlockSizes(Mat, PetscViewer);
538: PETSC_INTERN PetscErrorCode MatLoad_Binary_BlockSizes(Mat, PetscViewer);

540: /*
541:     Object for partitioning graphs
542: */

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

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

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

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

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

595: PETSC_EXTERN PetscErrorCode MatCoarsenMISKSetDistance(MatCoarsen, PetscInt);
596: PETSC_EXTERN PetscErrorCode MatCoarsenMISKGetDistance(MatCoarsen, PetscInt *);

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

609: /*
610:     MatFDColoring is used to compute Jacobian matrices efficiently
611:   via coloring. The data structure is explained below in an example.

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

625:     ncolors = 4;

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

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

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

644: */
645: typedef struct {
646:   PetscInt     row;
647:   PetscInt     col;
648:   PetscScalar *valaddr; /* address of value */
649: } MatEntry;

651: typedef struct {
652:   PetscInt     row;
653:   PetscScalar *valaddr; /* address of value */
654: } MatEntry2;

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

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

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

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

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

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

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

749: PETSC_EXTERN PetscErrorCode MatTransposeCheckNonzeroState_Private(Mat, Mat);

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

760: PETSC_EXTERN PetscErrorCode MatFactorDumpMatrix(Mat);
761: PETSC_INTERN PetscErrorCode MatSetBlockSizes_Default(Mat, PetscInt, PetscInt);

763: PETSC_SINGLE_LIBRARY_INTERN PetscErrorCode MatShift_Basic(Mat, PetscScalar);

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

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

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

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

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

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

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

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

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

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

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

862: 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)
863: {
864:   PetscInt location;

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

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

889:     if (!PetscBTLookupSet(bt, entry)) PetscCall(PetscLLInsertLocation_Private(assume_sorted, k, idx_start, entry, nlnk, &lnkdata, lnk));
890:   }
891:   PetscFunctionReturn(PETSC_SUCCESS);
892: }

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

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

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

955:     if (!PetscBTLookupSet(bt, entry)) PetscCall(PetscLLInsertLocation_Private(PETSC_FALSE, k, idx_start, entry, nlnk, &lnkdata, lnk));
956:   }
957:   PetscFunctionReturn(PETSC_SUCCESS);
958: }

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

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

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

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

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

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

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

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

1060: /*
1061:   Free memories used by the list
1062: */
1063: #define PetscLLDestroy(lnk, bt) ((PetscErrorCode)(PetscFree(lnk) || PetscBTDestroy(&(bt))))

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

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

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

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

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

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

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

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

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

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

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

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

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

1265: /*
1266:   Free memories used by the list
1267: */
1268: #define PetscIncompleteLLDestroy(lnk, bt) ((PetscErrorCode)(PetscFree(lnk) || PetscBTDestroy(&(bt))))

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

1449: static inline PetscErrorCode PetscLLCondensedExpand_Scalable(PetscInt nlnk_max, PetscInt **lnk)
1450: {
1451:   PetscInt lsize = 0;

1453:   PetscFunctionBegin;
1454:   PetscCall(PetscIntMultError(2, nlnk_max + 2, &lsize));
1455:   PetscCall(PetscRealloc(lsize * sizeof(PetscInt), lnk));
1456:   PetscFunctionReturn(PETSC_SUCCESS);
1457: }

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

1464:   for (PetscInt k = 0; k < nidx; k++) {
1465:     const PetscInt entry = indices[k];
1466:     PetscInt       next, lnkdata;

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

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

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

1502: static inline PetscErrorCode PetscLLCondensedDestroy_Scalable(PetscInt *lnk)
1503: {
1504:   return PetscFree(lnk);
1505: }

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

1515:       The next three are always the first link

1517:       lnk[3]    PETSC_MIN_INT+1
1518:       lnk[4]    1
1519:       lnk[5]    link to first real entry

1521:       The next three are always the last link

1523:       lnk[6]    PETSC_MAX_INT - 1
1524:       lnk[7]    1
1525:       lnk[8]    next valid link (this is the same as lnk[0] but without the decreases)
1526: */

1528: static inline PetscErrorCode PetscLLCondensedCreate_fast(PetscInt nlnk_max, PetscInt **lnk)
1529: {
1530:   PetscInt *llnk;
1531:   PetscInt  lsize = 0;

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

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

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

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

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

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

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

1626: static inline PetscErrorCode PetscLLCondensedDestroy_fast(PetscInt *lnk)
1627: {
1628:   return PetscFree(lnk);
1629: }

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

1650: PETSC_EXTERN PetscErrorCode PetscCDGetHeadPos(const PetscCoarsenData *, PetscInt, PetscCDIntNd **);
1651: PETSC_EXTERN PetscErrorCode PetscCDGetNextPos(const PetscCoarsenData *, PetscInt, PetscCDIntNd **);
1652: PETSC_EXTERN PetscErrorCode PetscCDGetASMBlocks(const PetscCoarsenData *, const PetscInt, PetscInt *, IS **);

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

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