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_INTERN 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
 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 (*matmatmultsymbolic)(Mat, Mat, Mat, PetscReal, Mat);
114:   /*64*/
115:   PetscErrorCode (*matmatmultnumeric)(Mat, Mat, Mat, Mat);
116:   PetscErrorCode (*setlocaltoglobalmapping)(Mat, ISLocalToGlobalMapping, ISLocalToGlobalMapping);
117:   PetscErrorCode (*setvalueslocal)(Mat, PetscInt, const PetscInt[], PetscInt, const PetscInt[], const PetscScalar[], InsertMode);
118:   PetscErrorCode (*zerorowslocal)(Mat, PetscInt, const PetscInt[], PetscScalar, Vec, Vec);
119:   PetscErrorCode (*getrowmaxabs)(Mat, Vec, PetscInt[]);
120:   /*69*/
121:   PetscErrorCode (*getrowminabs)(Mat, Vec, PetscInt[]);
122:   PetscErrorCode (*convert)(Mat, MatType, MatReuse, Mat *);
123:   PetscErrorCode (*hasoperation)(Mat, MatOperation, PetscBool *);
124:   PetscErrorCode (*fdcoloringapply)(Mat, MatFDColoring, Vec, void *);
125:   PetscErrorCode (*setfromoptions)(Mat, PetscOptionItems);
126:   /*74*/
127:   PetscErrorCode (*findzerodiagonals)(Mat, IS *);
128:   PetscErrorCode (*mults)(Mat, Vecs, Vecs);
129:   PetscErrorCode (*solves)(Mat, Vecs, Vecs);
130:   PetscErrorCode (*getinertia)(Mat, PetscInt *, PetscInt *, PetscInt *);
131:   PetscErrorCode (*load)(Mat, PetscViewer);
132:   /*79*/
133:   PetscErrorCode (*issymmetric)(Mat, PetscReal, PetscBool *);
134:   PetscErrorCode (*ishermitian)(Mat, PetscReal, PetscBool *);
135:   PetscErrorCode (*isstructurallysymmetric)(Mat, PetscBool *);
136:   PetscErrorCode (*setvaluesblockedlocal)(Mat, PetscInt, const PetscInt[], PetscInt, const PetscInt[], const PetscScalar[], InsertMode);
137:   PetscErrorCode (*getvecs)(Mat, Vec *, Vec *);
138:   /*84*/
139:   PetscErrorCode (*matmultsymbolic)(Mat, Mat, PetscReal, Mat);
140:   PetscErrorCode (*matmultnumeric)(Mat, Mat, Mat);
141:   PetscErrorCode (*ptapnumeric)(Mat, Mat, Mat); /* double dispatch wrapper routine */
142:   PetscErrorCode (*mattransposemultsymbolic)(Mat, Mat, PetscReal, Mat);
143:   PetscErrorCode (*mattransposemultnumeric)(Mat, Mat, Mat);
144:   /*89*/
145:   PetscErrorCode (*bindtocpu)(Mat, PetscBool);
146:   PetscErrorCode (*productsetfromoptions)(Mat);
147:   PetscErrorCode (*productsymbolic)(Mat);
148:   PetscErrorCode (*productnumeric)(Mat);
149:   PetscErrorCode (*conjugate)(Mat); /* complex conjugate */
150:   /*94*/
151:   PetscErrorCode (*viewnative)(Mat, PetscViewer);
152:   PetscErrorCode (*setvaluesrow)(Mat, PetscInt, const PetscScalar[]);
153:   PetscErrorCode (*realpart)(Mat);
154:   PetscErrorCode (*imaginarypart)(Mat);
155:   PetscErrorCode (*getrowuppertriangular)(Mat);
156:   /*99*/
157:   PetscErrorCode (*restorerowuppertriangular)(Mat);
158:   PetscErrorCode (*matsolve)(Mat, Mat, Mat);
159:   PetscErrorCode (*matsolvetranspose)(Mat, Mat, Mat);
160:   PetscErrorCode (*getrowmin)(Mat, Vec, PetscInt[]);
161:   PetscErrorCode (*getcolumnvector)(Mat, Vec, PetscInt);
162:   /*104*/
163:   PetscErrorCode (*missingdiagonal)(Mat, PetscBool *, PetscInt *);
164:   PetscErrorCode (*getseqnonzerostructure)(Mat, Mat *);
165:   PetscErrorCode (*create)(Mat);
166:   PetscErrorCode (*getghosts)(Mat, PetscInt *, const PetscInt *[]);
167:   PetscErrorCode (*getlocalsubmatrix)(Mat, IS, IS, Mat *);
168:   /*109*/
169:   PetscErrorCode (*restorelocalsubmatrix)(Mat, IS, IS, Mat *);
170:   PetscErrorCode (*multdiagonalblock)(Mat, Vec, Vec);
171:   PetscErrorCode (*hermitiantranspose)(Mat, MatReuse, Mat *);
172:   PetscErrorCode (*multhermitiantranspose)(Mat, Vec, Vec);
173:   PetscErrorCode (*multhermitiantransposeadd)(Mat, Vec, Vec, Vec);
174:   /*114*/
175:   PetscErrorCode (*getmultiprocblock)(Mat, MPI_Comm, MatReuse, Mat *);
176:   PetscErrorCode (*findnonzerorows)(Mat, IS *);
177:   PetscErrorCode (*getcolumnreductions)(Mat, PetscInt, PetscReal *);
178:   PetscErrorCode (*invertblockdiagonal)(Mat, const PetscScalar **);
179:   PetscErrorCode (*invertvariableblockdiagonal)(Mat, PetscInt, const PetscInt *, PetscScalar *);
180:   /*119*/
181:   PetscErrorCode (*createsubmatricesmpi)(Mat, PetscInt, const IS[], const IS[], MatReuse, Mat **);
182:   PetscErrorCode (*setvaluesbatch)(Mat, PetscInt, PetscInt, PetscInt *, const PetscScalar *);
183:   PetscErrorCode (*transposematmultsymbolic)(Mat, Mat, PetscReal, Mat);
184:   PetscErrorCode (*transposematmultnumeric)(Mat, Mat, Mat);
185:   PetscErrorCode (*transposecoloringcreate)(Mat, ISColoring, MatTransposeColoring);
186:   /*124*/
187:   PetscErrorCode (*transcoloringapplysptoden)(MatTransposeColoring, Mat, Mat);
188:   PetscErrorCode (*transcoloringapplydentosp)(MatTransposeColoring, Mat, Mat);
189:   PetscErrorCode (*rartnumeric)(Mat, Mat, Mat); /* double dispatch wrapper routine */
190:   PetscErrorCode (*setblocksizes)(Mat, PetscInt, PetscInt);
191:   PetscErrorCode (*residual)(Mat, Vec, Vec, Vec);
192:   /*129*/
193:   PetscErrorCode (*fdcoloringsetup)(Mat, ISColoring, MatFDColoring);
194:   PetscErrorCode (*findoffblockdiagonalentries)(Mat, IS *);
195:   PetscErrorCode (*creatempimatconcatenateseqmat)(MPI_Comm, Mat, PetscInt, MatReuse, Mat *);
196:   PetscErrorCode (*destroysubmatrices)(PetscInt, Mat *[]);
197:   PetscErrorCode (*mattransposesolve)(Mat, Mat, Mat);
198:   /*134*/
199:   PetscErrorCode (*getvalueslocal)(Mat, PetscInt, const PetscInt[], PetscInt, const PetscInt[], PetscScalar[]);
200:   PetscErrorCode (*creategraph)(Mat, PetscBool, PetscBool, PetscReal, PetscInt, PetscInt[], Mat *);
201:   PetscErrorCode (*transposesymbolic)(Mat, Mat *);
202:   PetscErrorCode (*eliminatezeros)(Mat, PetscBool);
203:   PetscErrorCode (*getrowsumabs)(Mat, Vec);
204:   /*139*/
205:   PetscErrorCode (*getfactor)(Mat, MatSolverType, MatFactorType, Mat *);
206:   PetscErrorCode (*getblockdiagonal)(Mat, Mat *);  // NOTE: the caller of get{block, vblock}diagonal owns the returned matrix;
207:   PetscErrorCode (*getvblockdiagonal)(Mat, Mat *); // they must destroy it after use
208:   PetscErrorCode (*copyhashtoxaij)(Mat, Mat);
209:   PetscErrorCode (*getcurrentmemtype)(Mat, PetscMemType *);
210:   /*144*/
211:   PetscErrorCode (*zerorowscolumnslocal)(Mat, PetscInt, const PetscInt[], PetscScalar, Vec, Vec);
212: };
213: /*
214:     If you add MatOps entries above also add them to the MATOP enum
215:     in include/petscmat.h
216: */

218: #include <petscsys.h>

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

226: PETSC_EXTERN MatRootName MatRootNameList;

228: /*
229:    Utility private matrix routines used outside Mat
230: */
231: PETSC_SINGLE_LIBRARY_INTERN PetscErrorCode MatFindNonzeroRowsOrCols_Basic(Mat, PetscBool, PetscReal, IS *);
232: PETSC_EXTERN PetscErrorCode                MatShellGetScalingShifts(Mat, PetscScalar *, PetscScalar *, Vec *, Vec *, Vec *, Mat *, IS *, IS *);

234: #define MAT_SHELL_NOT_ALLOWED (void *)-1

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

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

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

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

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

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

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

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

308: typedef struct _MatStashSpace *PetscMatStashSpace;

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

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

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

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

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

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

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

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

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

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

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

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

435: typedef struct { /* used by MatProduct() */
436:   MatProductType type;
437:   char          *alg;
438:   Mat            A, B, C, Dwork;
439:   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, .. */
440:   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 .. */
441:   PetscBool      symbolic_used_the_fact_C_is_symmetric; /* MatMatMult(A,B,MAT_REUSE_MATRIX,..&C) is still legitimate), we need to redo symbolic! */
442:   PetscObjectParameterDeclare(PetscReal, fill);
443:   PetscBool api_user; /* used to distinguish command line options and to indicate the matrix values are ready to be consumed at symbolic phase if needed */
444:   PetscBool setfromoptionscalled;

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

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

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

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

518: PETSC_INTERN PetscErrorCode MatSetUp_Default(Mat);

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

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

531: /*
532:     Object for partitioning graphs
533: */

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

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

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

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

570: #define MAT_COARSEN_STRENGTH_INDEX_SIZE 3
571: struct _p_MatCoarsen {
572:   PETSCHEADER(struct _MatCoarsenOps);
573:   Mat   graph;
574:   void *subctx;
575:   /* */
576:   PetscBool         strict_aggs;
577:   IS                perm;
578:   PetscCoarsenData *agg_lists;
579:   PetscInt          max_it;    /* number of iterations in HEM */
580:   PetscReal         threshold; /* HEM can filter interim graphs */
581:   PetscInt          strength_index_size;
582:   PetscInt          strength_index[MAT_COARSEN_STRENGTH_INDEX_SIZE];
583: };

585: PETSC_EXTERN PetscErrorCode MatCoarsenMISKSetDistance(MatCoarsen, PetscInt);
586: PETSC_EXTERN PetscErrorCode MatCoarsenMISKGetDistance(MatCoarsen, PetscInt *);

588: /*
589:     Used in aijdevice.h
590: */
591: typedef struct {
592:   PetscInt    *i;
593:   PetscInt    *j;
594:   PetscScalar *a;
595:   PetscInt     n;
596:   PetscInt     ignorezeroentries;
597: } PetscCSRDataStructure;

599: /*
600:     MatFDColoring is used to compute Jacobian matrices efficiently
601:   via coloring. The data structure is explained below in an example.

603:    Color =   0    1     0    2   |   2      3       0
604:    ---------------------------------------------------
605:             00   01              |          05
606:             10   11              |   14     15               Processor  0
607:                        22    23  |          25
608:                        32    33  |
609:    ===================================================
610:                                  |   44     45     46
611:             50                   |          55               Processor 1
612:                                  |   64            66
613:    ---------------------------------------------------

615:     ncolors = 4;

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

624:     ncolumns      = {1,0,1,1}
625:     columns       = {{6},{},{4},{5}}
626:     nrows         = {3,0,2,2}
627:     rows          = {{0,1,2},{},{1,2},{1,2}}
628:     vwscale       = {dx(4),dx(5),dx(6)}              MPI Vec
629:     vscale        = {dx(0),dx(4),dx(5),dx(6)}        Seq Vec

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

634: */
635: typedef struct {
636:   PetscInt     row;
637:   PetscInt     col;
638:   PetscScalar *valaddr; /* address of value */
639: } MatEntry;

641: typedef struct {
642:   PetscInt     row;
643:   PetscScalar *valaddr; /* address of value */
644: } MatEntry2;

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

676: typedef struct _MatColoringOps *MatColoringOps;
677: struct _MatColoringOps {
678:   PetscErrorCode (*destroy)(MatColoring);
679:   PetscErrorCode (*setfromoptions)(MatColoring, PetscOptionItems);
680:   PetscErrorCode (*view)(MatColoring, PetscViewer);
681:   PetscErrorCode (*apply)(MatColoring, ISColoring *);
682:   PetscErrorCode (*weights)(MatColoring, PetscReal **, PetscInt **);
683: };

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

698: struct _p_MatTransposeColoring {
699:   PETSCHEADER(int);
700:   PetscInt       M, N, m;      /* total rows, columns; local rows */
701:   PetscInt       rstart;       /* first row owned by local processor */
702:   PetscInt       ncolors;      /* number of colors */
703:   PetscInt      *ncolumns;     /* number of local columns for a color */
704:   PetscInt      *nrows;        /* number of local rows for each color */
705:   PetscInt       currentcolor; /* color for which function evaluation is being done now */
706:   ISColoringType ctype;        /* IS_COLORING_GLOBAL or IS_COLORING_LOCAL */

708:   PetscInt *colorforrow, *colorforcol; /* pointer to rows and columns */
709:   PetscInt *rows;                      /* lists the local rows for each color (using the local row numbering) */
710:   PetscInt *den2sp;                    /* maps (row,color) in the dense matrix to index of sparse matrix array a->a */
711:   PetscInt *columns;                   /* lists the local columns of each color (using global column numbering) */
712:   PetscInt  brows;                     /* number of rows for efficient implementation of MatTransColoringApplyDenToSp() */
713:   PetscInt *lstart;                    /* array used for loop over row blocks of Csparse */
714: };

716: /*
717:    Null space context for preconditioner/operators
718: */
719: struct _p_MatNullSpace {
720:   PETSCHEADER(int);
721:   PetscBool             has_cnst;
722:   PetscInt              n;
723:   Vec                  *vecs;
724:   PetscScalar          *alpha;  /* for projections */
725:   MatNullSpaceRemoveFn *remove; /* for user provided removal function */
726:   void                 *rmctx;  /* context for remove() function */
727: };

729: /*
730:    Internal data structure for MATMPIDENSE
731: */
732: typedef struct {
733:   Mat A; /* local submatrix */

735:   /* The following variables are used for matrix assembly */
736:   PetscBool    donotstash;        /* Flag indicating if values should be stashed */
737:   MPI_Request *send_waits;        /* array of send requests */
738:   MPI_Request *recv_waits;        /* array of receive requests */
739:   PetscInt     nsends, nrecvs;    /* numbers of sends and receives */
740:   PetscScalar *svalues, *rvalues; /* sending and receiving data */
741:   PetscInt     rmax;              /* maximum message length */

743:   /* The following variables are used for matrix-vector products */
744:   Vec       lvec;        /* local vector */
745:   PetscSF   Mvctx;       /* for mat-mult communications */
746:   PetscBool roworiented; /* if true, row-oriented input (default) */

748:   /* Support for MatDenseGetColumnVec and MatDenseGetSubMatrix */
749:   Mat                cmat;     /* matrix representation of a given subset of columns */
750:   Vec                cvec;     /* vector representation of a given column */
751:   const PetscScalar *ptrinuse; /* holds array to be restored (just a placeholder) */
752:   PetscInt           vecinuse; /* if cvec is in use (col = vecinuse-1) */
753:   PetscInt           matinuse; /* if cmat is in use (cbegin = matinuse-1) */
754:   /* if this is from MatDenseGetSubMatrix, which columns and rows does it correspond to? */
755:   PetscInt sub_rbegin;
756:   PetscInt sub_rend;
757:   PetscInt sub_cbegin;
758:   PetscInt sub_cend;
759: } Mat_MPIDense;

761: /*
762:    Checking zero pivot for LU, ILU preconditioners.
763: */
764: typedef struct {
765:   PetscInt    nshift, nshift_max;
766:   PetscReal   shift_amount, shift_lo, shift_hi, shift_top, shift_fraction;
767:   PetscBool   newshift;
768:   PetscReal   rs; /* active row sum of abs(off-diagonals) */
769:   PetscScalar pv; /* pivot of the active row */
770: } FactorShiftCtx;

772: PETSC_SINGLE_LIBRARY_INTERN PetscErrorCode MatTransposeCheckNonzeroState_Private(Mat, Mat);

774: /*
775:  Used by MatTranspose() and potentially other functions to track the matrix used in the generation of another matrix
776: */
777: typedef struct {
778:   PetscObjectId    id;
779:   PetscObjectState state;
780:   PetscObjectState nonzerostate;
781: } MatParentState;

783: PETSC_EXTERN PetscErrorCode MatFactorDumpMatrix(Mat);
784: PETSC_INTERN PetscErrorCode MatSetBlockSizes_Default(Mat, PetscInt, PetscInt);

786: PETSC_SINGLE_LIBRARY_INTERN PetscErrorCode MatShift_Basic(Mat, PetscScalar);

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

793:   PetscFunctionBegin;
794:   if (PetscAbsScalar(sctx->pv) <= _zero && !PetscIsNanScalar(sctx->pv)) {
795:     /* force |diag| > zeropivot*rs */
796:     if (!sctx->nshift) sctx->shift_amount = info->shiftamount;
797:     else sctx->shift_amount *= 2.0;
798:     sctx->newshift = PETSC_TRUE;
799:     (sctx->nshift)++;
800:   } else {
801:     sctx->newshift = PETSC_FALSE;
802:   }
803:   PetscFunctionReturn(PETSC_SUCCESS);
804: }

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

811:   PetscFunctionBegin;
812:   if (PetscRealPart(sctx->pv) <= _zero && !PetscIsNanScalar(sctx->pv)) {
813:     /* force matfactor to be diagonally dominant */
814:     if (sctx->nshift == sctx->nshift_max) {
815:       sctx->shift_fraction = sctx->shift_hi;
816:     } else {
817:       sctx->shift_lo       = sctx->shift_fraction;
818:       sctx->shift_fraction = (sctx->shift_hi + sctx->shift_lo) / (PetscReal)2.;
819:     }
820:     sctx->shift_amount = sctx->shift_fraction * sctx->shift_top;
821:     sctx->nshift++;
822:     sctx->newshift = PETSC_TRUE;
823:   } else {
824:     sctx->newshift = PETSC_FALSE;
825:   }
826:   PetscFunctionReturn(PETSC_SUCCESS);
827: }

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

833:   PetscFunctionBegin;
834:   if (PetscAbsScalar(sctx->pv) <= _zero && !PetscIsNanScalar(sctx->pv)) {
835:     sctx->pv += info->shiftamount;
836:     sctx->shift_amount = 0.0;
837:     sctx->nshift++;
838:   }
839:   sctx->newshift = PETSC_FALSE;
840:   PetscFunctionReturn(PETSC_SUCCESS);
841: }

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

847:   PetscFunctionBegin;
848:   sctx->newshift = PETSC_FALSE;
849:   if (PetscAbsScalar(sctx->pv) <= _zero && !PetscIsNanScalar(sctx->pv)) {
850:     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);
851:     PetscCall(PetscInfo(mat, "Detected zero pivot in factorization in row %" PetscInt_FMT " value %g tolerance %g\n", row, (double)PetscAbsScalar(sctx->pv), (double)_zero));
852:     fact->factorerrortype             = MAT_FACTOR_NUMERIC_ZEROPIVOT;
853:     fact->factorerror_zeropivot_value = PetscAbsScalar(sctx->pv);
854:     fact->factorerror_zeropivot_row   = row;
855:   }
856:   PetscFunctionReturn(PETSC_SUCCESS);
857: }

859: static inline PetscErrorCode MatPivotCheck(Mat fact, Mat mat, const MatFactorInfo *info, FactorShiftCtx *sctx, PetscInt row)
860: {
861:   PetscFunctionBegin;
862:   if (info->shifttype == (PetscReal)MAT_SHIFT_NONZERO) PetscCall(MatPivotCheck_nz(mat, info, sctx, row));
863:   else if (info->shifttype == (PetscReal)MAT_SHIFT_POSITIVE_DEFINITE) PetscCall(MatPivotCheck_pd(mat, info, sctx, row));
864:   else if (info->shifttype == (PetscReal)MAT_SHIFT_INBLOCKS) PetscCall(MatPivotCheck_inblocks(mat, info, sctx, row));
865:   else PetscCall(MatPivotCheck_none(fact, mat, info, sctx, row));
866:   PetscFunctionReturn(PETSC_SUCCESS);
867: }

869: #include <petscbt.h>
870: /*
871:   Create and initialize a linked list
872:   Input Parameters:
873:     idx_start - starting index of the list
874:     lnk_max   - max value of lnk indicating the end of the list
875:     nlnk      - max length of the list
876:   Output Parameters:
877:     lnk       - list initialized
878:     bt        - PetscBT (bitarray) with all bits set to false
879:     lnk_empty - flg indicating the list is empty
880: */
881: #define PetscLLCreate(idx_start, lnk_max, nlnk, lnk, bt) ((PetscErrorCode)(PetscMalloc1(nlnk, &lnk) || PetscBTCreate(nlnk, &(bt)) || (lnk[idx_start] = lnk_max, PETSC_SUCCESS)))

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

885: 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)
886: {
887:   PetscInt location;

889:   PetscFunctionBegin;
890:   /* start from the beginning if entry < previous entry */
891:   if (!assume_sorted && k && entry < *lnkdata) *lnkdata = idx_start;
892:   /* search for insertion location */
893:   do {
894:     location = *lnkdata;
895:     *lnkdata = lnk[location];
896:   } while (entry > *lnkdata);
897:   /* insertion location is found, add entry into lnk */
898:   lnk[location] = entry;
899:   lnk[entry]    = *lnkdata;
900:   ++(*nlnk);
901:   *lnkdata = entry; /* next search starts from here if next_entry > entry */
902:   PetscFunctionReturn(PETSC_SUCCESS);
903: }

905: 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)
906: {
907:   PetscFunctionBegin;
908:   *nlnk = 0;
909:   for (PetscInt k = 0, lnkdata = idx_start; k < nidx; ++k) {
910:     const PetscInt entry = indices[k];

912:     if (!PetscBTLookupSet(bt, entry)) PetscCall(PetscLLInsertLocation_Private(assume_sorted, k, idx_start, entry, nlnk, &lnkdata, lnk));
913:   }
914:   PetscFunctionReturn(PETSC_SUCCESS);
915: }

917: /*
918:   Add an index set into a sorted linked list
919:   Input Parameters:
920:     nidx      - number of input indices
921:     indices   - 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 PetscLLAdd(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_FALSE));
934:   PetscFunctionReturn(PETSC_SUCCESS);
935: }

937: /*
938:   Add a SORTED ascending index set into a sorted linked list - same as PetscLLAdd() bus skip 'if (_k && _entry < _lnkdata) _lnkdata  = idx_start;'
939:   Input Parameters:
940:     nidx      - number of input indices
941:     indices   - sorted integer array
942:     idx_start - starting index of the list
943:     lnk       - linked list(an integer array) that is created
944:     bt        - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
945:   output Parameters:
946:     nlnk      - number of newly added indices
947:     lnk       - the sorted(increasing order) linked list containing new and non-redundate entries from indices
948:     bt        - updated PetscBT (bitarray)
949: */
950: static inline PetscErrorCode PetscLLAddSorted(PetscInt nidx, const PetscInt *PETSC_RESTRICT indices, PetscInt idx_start, PetscInt *PETSC_RESTRICT nlnk, PetscInt *PETSC_RESTRICT lnk, PetscBT bt)
951: {
952:   PetscFunctionBegin;
953:   PetscCall(PetscLLAdd_Private(nidx, indices, idx_start, nlnk, lnk, bt, PETSC_TRUE));
954:   PetscFunctionReturn(PETSC_SUCCESS);
955: }

957: /*
958:   Add a permuted index set into a sorted linked list
959:   Input Parameters:
960:     nidx      - number of input indices
961:     indices   - integer array
962:     perm      - permutation of indices
963:     idx_start - starting index of the list
964:     lnk       - linked list(an integer array) that is created
965:     bt        - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
966:   output Parameters:
967:     nlnk      - number of newly added indices
968:     lnk       - the sorted(increasing order) linked list containing new and non-redundate entries from indices
969:     bt        - updated PetscBT (bitarray)
970: */
971: 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)
972: {
973:   PetscFunctionBegin;
974:   *nlnk = 0;
975:   for (PetscInt k = 0, lnkdata = idx_start; k < nidx; ++k) {
976:     const PetscInt entry = perm[indices[k]];

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

983: #if 0
984: /* this appears to be unused? */
985: static inline PetscErrorCode PetscLLAddSorted_new(PetscInt nidx, PetscInt *indices, PetscInt idx_start, PetscBool *lnk_empty, PetscInt *nlnk, PetscInt *lnk, PetscBT bt)
986: {
987:   PetscInt lnkdata = idx_start;

989:   PetscFunctionBegin;
990:   if (*lnk_empty) {
991:     for (PetscInt k = 0; k < nidx; ++k) {
992:       const PetscInt entry = indices[k], location = lnkdata;

994:       PetscCall(PetscBTSet(bt,entry)); /* mark the new entry */
995:       lnkdata       = lnk[location];
996:       /* insertion location is found, add entry into lnk */
997:       lnk[location] = entry;
998:       lnk[entry]    = lnkdata;
999:       lnkdata       = entry; /* next search starts from here */
1000:     }
1001:     /* lnk[indices[nidx-1]] = lnk[idx_start];
1002:        lnk[idx_start]       = indices[0];
1003:        PetscCall(PetscBTSet(bt,indices[0]));
1004:        for (_k=1; _k<nidx; _k++) {
1005:        PetscCall(PetscBTSet(bt,indices[_k]));
1006:        lnk[indices[_k-1]] = indices[_k];
1007:        }
1008:     */
1009:     *nlnk      = nidx;
1010:     *lnk_empty = PETSC_FALSE;
1011:   } else {
1012:     *nlnk = 0;
1013:     for (PetscInt k = 0; k < nidx; ++k) {
1014:       const PetscInt entry = indices[k];

1016:       if (!PetscBTLookupSet(bt,entry)) PetscCall(PetscLLInsertLocation_Private(PETSC_TRUE,k,idx_start,entry,nlnk,&lnkdata,lnk));
1017:     }
1018:   }
1019:   PetscFunctionReturn(PETSC_SUCCESS);
1020: }
1021: #endif

1023: /*
1024:   Add a SORTED index set into a sorted linked list used for LUFactorSymbolic()
1025:   Same as PetscLLAddSorted() with an additional operation:
1026:        count the number of input indices that are no larger than 'diag'
1027:   Input Parameters:
1028:     indices   - sorted integer array
1029:     idx_start - starting index of the list, index of pivot row
1030:     lnk       - linked list(an integer array) that is created
1031:     bt        - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
1032:     diag      - index of the active row in LUFactorSymbolic
1033:     nzbd      - number of input indices with indices <= idx_start
1034:     im        - im[idx_start] is initialized as num of nonzero entries in row=idx_start
1035:   output Parameters:
1036:     nlnk      - number of newly added indices
1037:     lnk       - the sorted(increasing order) linked list containing new and non-redundate entries from indices
1038:     bt        - updated PetscBT (bitarray)
1039:     im        - im[idx_start]: unchanged if diag is not an entry
1040:                              : num of entries with indices <= diag if diag is an entry
1041: */
1042: 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)
1043: {
1044:   const PetscInt nidx = im[idx_start] - nzbd; /* num of entries with idx_start < index <= diag */

1046:   PetscFunctionBegin;
1047:   *nlnk = 0;
1048:   for (PetscInt k = 0, lnkdata = idx_start; k < nidx; ++k) {
1049:     const PetscInt entry = indices[k];

1051:     ++nzbd;
1052:     if (entry == diag) im[idx_start] = nzbd;
1053:     if (!PetscBTLookupSet(bt, entry)) PetscCall(PetscLLInsertLocation_Private(PETSC_TRUE, k, idx_start, entry, nlnk, &lnkdata, lnk));
1054:   }
1055:   PetscFunctionReturn(PETSC_SUCCESS);
1056: }

1058: /*
1059:   Copy data on the list into an array, then initialize the list
1060:   Input Parameters:
1061:     idx_start - starting index of the list
1062:     lnk_max   - max value of lnk indicating the end of the list
1063:     nlnk      - number of data on the list to be copied
1064:     lnk       - linked list
1065:     bt        - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
1066:   output Parameters:
1067:     indices   - array that contains the copied data
1068:     lnk       - linked list that is cleaned and initialize
1069:     bt        - PetscBT (bitarray) with all bits set to false
1070: */
1071: static inline PetscErrorCode PetscLLClean(PetscInt idx_start, PetscInt lnk_max, PetscInt nlnk, PetscInt *PETSC_RESTRICT lnk, PetscInt *PETSC_RESTRICT indices, PetscBT bt)
1072: {
1073:   PetscFunctionBegin;
1074:   for (PetscInt j = 0, idx = idx_start; j < nlnk; ++j) {
1075:     idx        = lnk[idx];
1076:     indices[j] = idx;
1077:     PetscCall(PetscBTClear(bt, idx));
1078:   }
1079:   lnk[idx_start] = lnk_max;
1080:   PetscFunctionReturn(PETSC_SUCCESS);
1081: }

1083: /*
1084:   Free memories used by the list
1085: */
1086: #define PetscLLDestroy(lnk, bt) ((PetscErrorCode)(PetscFree(lnk) || PetscBTDestroy(&(bt))))

1088: /* Routines below are used for incomplete matrix factorization */
1089: /*
1090:   Create and initialize a linked list and its levels
1091:   Input Parameters:
1092:     idx_start - starting index of the list
1093:     lnk_max   - max value of lnk indicating the end of the list
1094:     nlnk      - max length of the list
1095:   Output Parameters:
1096:     lnk       - list initialized
1097:     lnk_lvl   - array of size nlnk for storing levels of lnk
1098:     bt        - PetscBT (bitarray) with all bits set to false
1099: */
1100: #define PetscIncompleteLLCreate(idx_start, lnk_max, nlnk, lnk, lnk_lvl, bt) \
1101:   ((PetscErrorCode)(PetscIntMultError(2, nlnk, NULL) || PetscMalloc1(2 * nlnk, &lnk) || PetscBTCreate(nlnk, &(bt)) || (lnk[idx_start] = lnk_max, lnk_lvl = lnk + nlnk, PETSC_SUCCESS)))

1103: 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)
1104: {
1105:   PetscFunctionBegin;
1106:   PetscCall(PetscLLInsertLocation_Private(assume_sorted, k, idx_start, entry, nlnk, lnkdata, lnk));
1107:   lnklvl[entry] = newval;
1108:   PetscFunctionReturn(PETSC_SUCCESS);
1109: }

1111: /*
1112:   Initialize a sorted linked list used for ILU and ICC
1113:   Input Parameters:
1114:     nidx      - number of input idx
1115:     idx       - integer array used for storing column indices
1116:     idx_start - starting index of the list
1117:     perm      - indices of an IS
1118:     lnk       - linked list(an integer array) that is created
1119:     lnklvl    - levels of lnk
1120:     bt        - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
1121:   output Parameters:
1122:     nlnk     - number of newly added idx
1123:     lnk      - the sorted(increasing order) linked list containing new and non-redundate entries from idx
1124:     lnklvl   - levels of lnk
1125:     bt       - updated PetscBT (bitarray)
1126: */
1127: 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)
1128: {
1129:   PetscFunctionBegin;
1130:   *nlnk = 0;
1131:   for (PetscInt k = 0, lnkdata = idx_start; k < nidx; ++k) {
1132:     const PetscInt entry = perm[idx[k]];

1134:     if (!PetscBTLookupSet(bt, entry)) PetscCall(PetscIncompleteLLInsertLocation_Private(PETSC_FALSE, k, idx_start, entry, nlnk, &lnkdata, lnk, lnklvl, 0));
1135:   }
1136:   PetscFunctionReturn(PETSC_SUCCESS);
1137: }

1139: 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)
1140: {
1141:   PetscFunctionBegin;
1142:   *nlnk = 0;
1143:   for (PetscInt k = 0, lnkdata = idx_start; k < nidx; ++k) {
1144:     const PetscInt incrlev = idxlvl[k] + prow_offset + 1;

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

1149:       if (!PetscBTLookupSet(bt, entry)) PetscCall(PetscIncompleteLLInsertLocation_Private(assume_sorted, k, idx_start, entry, nlnk, &lnkdata, lnk, lnklvl, incrlev));
1150:       else if (lnklvl[entry] > incrlev) lnklvl[entry] = incrlev; /* existing entry */
1151:     }
1152:   }
1153:   PetscFunctionReturn(PETSC_SUCCESS);
1154: }

1156: /*
1157:   Add a SORTED index set into a sorted linked list for ICC
1158:   Input Parameters:
1159:     nidx      - number of input indices
1160:     idx       - sorted integer array used for storing column indices
1161:     level     - level of fill, e.g., ICC(level)
1162:     idxlvl    - level of idx
1163:     idx_start - starting index of the list
1164:     lnk       - linked list(an integer array) that is created
1165:     lnklvl    - levels of lnk
1166:     bt        - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
1167:     idxlvl_prow - idxlvl[prow], where prow is the row number of the idx
1168:   output Parameters:
1169:     nlnk   - number of newly added indices
1170:     lnk    - the sorted(increasing order) linked list containing new and non-redundate entries from idx
1171:     lnklvl - levels of lnk
1172:     bt     - updated PetscBT (bitarray)
1173:   Note: the level of U(i,j) is set as lvl(i,j) = min{ lvl(i,j), lvl(prow,i)+lvl(prow,j)+1)
1174:         where idx = non-zero columns of U(prow,prow+1:n-1), prow<i
1175: */
1176: 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)
1177: {
1178:   PetscFunctionBegin;
1179:   PetscCall(PetscIncompleteLLAdd_Private(nidx, idx, level, idxlvl, idx_start, nlnk, lnk, lnklvl, bt, idxlvl_prow, PETSC_TRUE));
1180:   PetscFunctionReturn(PETSC_SUCCESS);
1181: }

1183: /*
1184:   Add a SORTED index set into a sorted linked list for ILU
1185:   Input Parameters:
1186:     nidx      - number of input indices
1187:     idx       - sorted integer array used for storing column indices
1188:     level     - level of fill, e.g., ICC(level)
1189:     idxlvl    - level of idx
1190:     idx_start - starting index of the list
1191:     lnk       - linked list(an integer array) that is created
1192:     lnklvl    - levels of lnk
1193:     bt        - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
1194:     prow      - the row number of idx
1195:   output Parameters:
1196:     nlnk     - number of newly added idx
1197:     lnk      - the sorted(increasing order) linked list containing new and non-redundate entries from idx
1198:     lnklvl   - levels of lnk
1199:     bt       - updated PetscBT (bitarray)

1201:   Note: the level of factor(i,j) is set as lvl(i,j) = min{ lvl(i,j), lvl(i,prow)+lvl(prow,j)+1)
1202:         where idx = non-zero columns of U(prow,prow+1:n-1), prow<i
1203: */
1204: 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)
1205: {
1206:   PetscFunctionBegin;
1207:   PetscCall(PetscIncompleteLLAdd_Private(nidx, idx, level, idxlvl, idx_start, nlnk, lnk, lnklvl, bt, lnklvl[prow], PETSC_TRUE));
1208:   PetscFunctionReturn(PETSC_SUCCESS);
1209: }

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

1235: /*
1236:   Add a SORTED index set into a sorted linked list
1237:   Input Parameters:
1238:     nidx      - number of input indices
1239:     idx   - sorted integer array used for storing column indices
1240:     level     - level of fill, e.g., ICC(level)
1241:     idxlvl - level of idx
1242:     idx_start - starting index of the list
1243:     lnk       - linked list(an integer array) that is created
1244:     lnklvl    - levels of lnk
1245:     bt        - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
1246:   output Parameters:
1247:     nlnk      - number of newly added idx
1248:     lnk       - the sorted(increasing order) linked list containing new and non-redundate entries from idx
1249:     lnklvl    - levels of lnk
1250:     bt        - updated PetscBT (bitarray)
1251: */
1252: 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)
1253: {
1254:   PetscFunctionBegin;
1255:   PetscCall(PetscIncompleteLLAdd_Private(nidx, idx, level, idxlvl, idx_start, nlnk, lnk, lnklvl, bt, 0, PETSC_TRUE));
1256:   PetscFunctionReturn(PETSC_SUCCESS);
1257: }

1259: /*
1260:   Copy data on the list into an array, then initialize the list
1261:   Input Parameters:
1262:     idx_start - starting index of the list
1263:     lnk_max   - max value of lnk indicating the end of the list
1264:     nlnk      - number of data on the list to be copied
1265:     lnk       - linked list
1266:     lnklvl    - level of lnk
1267:     bt        - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
1268:   output Parameters:
1269:     indices - array that contains the copied data
1270:     lnk     - linked list that is cleaned and initialize
1271:     lnklvl  - level of lnk that is reinitialized
1272:     bt      - PetscBT (bitarray) with all bits set to false
1273: */
1274: 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)
1275: {
1276:   PetscFunctionBegin;
1277:   for (PetscInt j = 0, idx = idx_start; j < nlnk; ++j) {
1278:     idx           = lnk[idx];
1279:     indices[j]    = idx;
1280:     indiceslvl[j] = lnklvl[idx];
1281:     lnklvl[idx]   = -1;
1282:     PetscCall(PetscBTClear(bt, idx));
1283:   }
1284:   lnk[idx_start] = lnk_max;
1285:   PetscFunctionReturn(PETSC_SUCCESS);
1286: }

1288: /*
1289:   Free memories used by the list
1290: */
1291: #define PetscIncompleteLLDestroy(lnk, bt) ((PetscErrorCode)(PetscFree(lnk) || PetscBTDestroy(&(bt))))

1293: #if !defined(PETSC_CLANG_STATIC_ANALYZER)
1294:   #define MatCheckSameLocalSize(A, ar1, B, ar2) \
1295:     do { \
1296:       PetscCheckSameComm(A, ar1, B, ar2); \
1297:       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, \
1298:                  (A)->rmap->n, (A)->cmap->n, ar2, (B)->rmap->n, (B)->cmap->n); \
1299:     } while (0)
1300:   #define MatCheckSameSize(A, ar1, B, ar2) \
1301:     do { \
1302:       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, \
1303:                  (A)->rmap->N, (A)->cmap->N, ar2, (B)->rmap->N, (B)->cmap->N); \
1304:       MatCheckSameLocalSize(A, ar1, B, ar2); \
1305:     } while (0)
1306: #else
1307: template <typename Tm>
1308: extern void MatCheckSameLocalSize(Tm, int, Tm, int);
1309: template <typename Tm>
1310: extern void MatCheckSameSize(Tm, int, Tm, int);
1311: #endif

1313: #define VecCheckMatCompatible(M, x, ar1, b, ar2) \
1314:   do { \
1315:     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, \
1316:                (M)->cmap->N); \
1317:     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, \
1318:                (M)->rmap->N); \
1319:   } while (0)

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

1328:       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
1329:       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
1330:       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
1331:       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.
1332:       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
1333:       to each other so memory access is much better than using the big array.

1335:   Example:
1336:      nlnk_max=5, lnk_max=36:
1337:      Initial list: [0, 0 | 36, 2 | 0, 0 | 0, 0 | 0, 0 | 0, 0 | 0, 0]
1338:      here, head_node has index 2 with value lnk[2]=lnk_max=36,
1339:            0-th entry is used to store the number of entries in the list,
1340:      The initial lnk represents head -> tail(marked by 36) with number of entries = lnk[0]=0.

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

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

1350:   Input Parameters:
1351:     nlnk_max  - max length of the list
1352:     lnk_max   - max value of the entries
1353:   Output Parameters:
1354:     lnk       - list created and initialized
1355:     bt        - PetscBT (bitarray) with all bits set to false. Note: bt has size lnk_max, not nln_max!
1356: */
1357: static inline PetscErrorCode PetscLLCondensedCreate(PetscInt nlnk_max, PetscInt lnk_max, PetscInt **lnk, PetscBT *bt)
1358: {
1359:   PetscInt *llnk, lsize = 0;

1361:   PetscFunctionBegin;
1362:   PetscCall(PetscIntMultError(2, nlnk_max + 2, &lsize));
1363:   PetscCall(PetscMalloc1(lsize, lnk));
1364:   PetscCall(PetscBTCreate(lnk_max, bt));
1365:   llnk    = *lnk;
1366:   llnk[0] = 0;       /* number of entries on the list */
1367:   llnk[2] = lnk_max; /* value in the head node */
1368:   llnk[3] = 2;       /* next for the head node */
1369:   PetscFunctionReturn(PETSC_SUCCESS);
1370: }

1372: /*
1373:   Add a SORTED ascending index set into a sorted linked list. See PetscLLCondensedCreate() for detailed description.
1374:   Input Parameters:
1375:     nidx      - number of input indices
1376:     indices   - sorted integer array
1377:     lnk       - condensed linked list(an integer array) that is created
1378:     bt        - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
1379:   output Parameters:
1380:     lnk       - the sorted(increasing order) linked list containing previous and newly added non-redundate indices
1381:     bt        - updated PetscBT (bitarray)
1382: */
1383: static inline PetscErrorCode PetscLLCondensedAddSorted(PetscInt nidx, const PetscInt indices[], PetscInt lnk[], PetscBT bt)
1384: {
1385:   PetscInt location = 2;      /* head */
1386:   PetscInt nlnk     = lnk[0]; /* num of entries on the input lnk */

1388:   PetscFunctionBegin;
1389:   for (PetscInt k = 0; k < nidx; k++) {
1390:     const PetscInt entry = indices[k];
1391:     if (!PetscBTLookupSet(bt, entry)) { /* new entry */
1392:       PetscInt next, lnkdata;

1394:       /* search for insertion location */
1395:       do {
1396:         next     = location + 1;  /* link from previous node to next node */
1397:         location = lnk[next];     /* idx of next node */
1398:         lnkdata  = lnk[location]; /* value of next node */
1399:       } while (entry > lnkdata);
1400:       /* insertion location is found, add entry into lnk */
1401:       const PetscInt newnode = 2 * (nlnk + 2); /* index for this new node */
1402:       lnk[next]              = newnode;        /* connect previous node to the new node */
1403:       lnk[newnode]           = entry;          /* set value of the new node */
1404:       lnk[newnode + 1]       = location;       /* connect new node to next node */
1405:       location               = newnode;        /* next search starts from the new node */
1406:       nlnk++;
1407:     }
1408:   }
1409:   lnk[0] = nlnk; /* number of entries in the list */
1410:   PetscFunctionReturn(PETSC_SUCCESS);
1411: }

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

1418:   PetscFunctionBegin;
1419:   for (PetscInt k = 0; k < nlnk; k++) {
1420:     indices[k] = lnk[next];
1421:     next       = lnk[next + 1];
1422:     PetscCall(PetscBTClear(bt, indices[k]));
1423:   }
1424:   lnk[0] = 0;       /* num of entries on the list */
1425:   lnk[2] = lnk_max; /* initialize head node */
1426:   lnk[3] = 2;       /* head node */
1427:   PetscFunctionReturn(PETSC_SUCCESS);
1428: }

1430: static inline PetscErrorCode PetscLLCondensedView(PetscInt *lnk)
1431: {
1432:   PetscFunctionBegin;
1433:   PetscCall(PetscPrintf(PETSC_COMM_SELF, "LLCondensed of size %" PetscInt_FMT ", (val,  next)\n", lnk[0]));
1434:   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]));
1435:   PetscFunctionReturn(PETSC_SUCCESS);
1436: }

1438: /*
1439:   Free memories used by the list
1440: */
1441: static inline PetscErrorCode PetscLLCondensedDestroy(PetscInt *lnk, PetscBT bt)
1442: {
1443:   PetscFunctionBegin;
1444:   PetscCall(PetscFree(lnk));
1445:   PetscCall(PetscBTDestroy(&bt));
1446:   PetscFunctionReturn(PETSC_SUCCESS);
1447: }

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

1460:   PetscFunctionBegin;
1461:   PetscCall(PetscIntMultError(2, nlnk_max + 2, &lsize));
1462:   PetscCall(PetscMalloc1(lsize, lnk));
1463:   llnk    = *lnk;
1464:   llnk[0] = 0;             /* number of entries on the list */
1465:   llnk[2] = PETSC_INT_MAX; /* value in the head node */
1466:   llnk[3] = 2;             /* next for the head node */
1467:   PetscFunctionReturn(PETSC_SUCCESS);
1468: }

1470: static inline PetscErrorCode PetscLLCondensedExpand_Scalable(PetscInt nlnk_max, PetscInt **lnk)
1471: {
1472:   PetscInt lsize = 0;

1474:   PetscFunctionBegin;
1475:   PetscCall(PetscIntMultError(2, nlnk_max + 2, &lsize));
1476:   PetscCall(PetscRealloc(lsize * sizeof(PetscInt), lnk));
1477:   PetscFunctionReturn(PETSC_SUCCESS);
1478: }

1480: static inline PetscErrorCode PetscLLCondensedAddSorted_Scalable(PetscInt nidx, const PetscInt indices[], PetscInt lnk[])
1481: {
1482:   PetscInt location = 2;      /* head */
1483:   PetscInt nlnk     = lnk[0]; /* num of entries on the input lnk */

1485:   for (PetscInt k = 0; k < nidx; k++) {
1486:     const PetscInt entry = indices[k];
1487:     PetscInt       next, lnkdata;

1489:     /* search for insertion location */
1490:     do {
1491:       next     = location + 1;  /* link from previous node to next node */
1492:       location = lnk[next];     /* idx of next node */
1493:       lnkdata  = lnk[location]; /* value of next node */
1494:     } while (entry > lnkdata);
1495:     if (entry < lnkdata) {
1496:       /* insertion location is found, add entry into lnk */
1497:       const PetscInt newnode = 2 * (nlnk + 2); /* index for this new node */
1498:       lnk[next]              = newnode;        /* connect previous node to the new node */
1499:       lnk[newnode]           = entry;          /* set value of the new node */
1500:       lnk[newnode + 1]       = location;       /* connect new node to next node */
1501:       location               = newnode;        /* next search starts from the new node */
1502:       nlnk++;
1503:     }
1504:   }
1505:   lnk[0] = nlnk; /* number of entries in the list */
1506:   return PETSC_SUCCESS;
1507: }

1509: static inline PetscErrorCode PetscLLCondensedClean_Scalable(PETSC_UNUSED PetscInt nidx, PetscInt *indices, PetscInt *lnk)
1510: {
1511:   const PetscInt nlnk = lnk[0];
1512:   PetscInt       next = lnk[3]; /* head node */

1514:   for (PetscInt k = 0; k < nlnk; k++) {
1515:     indices[k] = lnk[next];
1516:     next       = lnk[next + 1];
1517:   }
1518:   lnk[0] = 0; /* num of entries on the list */
1519:   lnk[3] = 2; /* head node */
1520:   return PETSC_SUCCESS;
1521: }

1523: static inline PetscErrorCode PetscLLCondensedDestroy_Scalable(PetscInt *lnk)
1524: {
1525:   return PetscFree(lnk);
1526: }

1528: /*
1529:       lnk[0]   number of links
1530:       lnk[1]   number of entries
1531:       lnk[3n]  value
1532:       lnk[3n+1] len
1533:       lnk[3n+2] link to next value

1535:       The next three are always the first link

1537:       lnk[3]    PETSC_INT_MIN+1
1538:       lnk[4]    1
1539:       lnk[5]    link to first real entry

1541:       The next three are always the last link

1543:       lnk[6]    PETSC_INT_MAX - 1
1544:       lnk[7]    1
1545:       lnk[8]    next valid link (this is the same as lnk[0] but without the decreases)
1546: */

1548: static inline PetscErrorCode PetscLLCondensedCreate_fast(PetscInt nlnk_max, PetscInt **lnk)
1549: {
1550:   PetscInt *llnk;
1551:   PetscInt  lsize = 0;

1553:   PetscFunctionBegin;
1554:   PetscCall(PetscIntMultError(3, nlnk_max + 3, &lsize));
1555:   PetscCall(PetscMalloc1(lsize, lnk));
1556:   llnk    = *lnk;
1557:   llnk[0] = 0;                 /* nlnk: number of entries on the list */
1558:   llnk[1] = 0;                 /* number of integer entries represented in list */
1559:   llnk[3] = PETSC_INT_MIN + 1; /* value in the first node */
1560:   llnk[4] = 1;                 /* count for the first node */
1561:   llnk[5] = 6;                 /* next for the first node */
1562:   llnk[6] = PETSC_INT_MAX - 1; /* value in the last node */
1563:   llnk[7] = 1;                 /* count for the last node */
1564:   llnk[8] = 0;                 /* next valid node to be used */
1565:   PetscFunctionReturn(PETSC_SUCCESS);
1566: }

1568: static inline PetscErrorCode PetscLLCondensedAddSorted_fast(PetscInt nidx, const PetscInt indices[], PetscInt lnk[])
1569: {
1570:   for (PetscInt k = 0, prev = 3 /* first value */; k < nidx; k++) {
1571:     const PetscInt entry = indices[k];
1572:     PetscInt       next  = lnk[prev + 2];

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

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

1617:   for (PetscInt k = 0, cnt = 0; k < nlnk; k++) {
1618:     for (PetscInt j = 0; j < lnk[next + 1]; j++) indices[cnt++] = lnk[next] + j;
1619:     next = lnk[next + 2];
1620:   }
1621:   lnk[0] = 0;                 /* nlnk: number of links */
1622:   lnk[1] = 0;                 /* number of integer entries represented in list */
1623:   lnk[3] = PETSC_INT_MIN + 1; /* value in the first node */
1624:   lnk[4] = 1;                 /* count for the first node */
1625:   lnk[5] = 6;                 /* next for the first node */
1626:   lnk[6] = PETSC_INT_MAX - 1; /* value in the last node */
1627:   lnk[7] = 1;                 /* count for the last node */
1628:   lnk[8] = 0;                 /* next valid location to make link */
1629:   return PETSC_SUCCESS;
1630: }

1632: static inline PetscErrorCode PetscLLCondensedView_fast(const PetscInt *lnk)
1633: {
1634:   const PetscInt nlnk = lnk[0];
1635:   PetscInt       next = lnk[5]; /* first node */

1637:   for (PetscInt k = 0; k < nlnk; k++) {
1638: #if 0 /* Debugging code */
1639:     printf("%d value %d len %d next %d\n", next, lnk[next], lnk[next + 1], lnk[next + 2]);
1640: #endif
1641:     next = lnk[next + 2];
1642:   }
1643:   return PETSC_SUCCESS;
1644: }

1646: static inline PetscErrorCode PetscLLCondensedDestroy_fast(PetscInt *lnk)
1647: {
1648:   return PetscFree(lnk);
1649: }

1651: PETSC_EXTERN PetscErrorCode PetscCDCreate(PetscInt, PetscCoarsenData **);
1652: PETSC_EXTERN PetscErrorCode PetscCDDestroy(PetscCoarsenData *);
1653: PETSC_EXTERN PetscErrorCode PetscCDIntNdSetID(PetscCDIntNd *, PetscInt);
1654: PETSC_EXTERN PetscErrorCode PetscCDIntNdGetID(const PetscCDIntNd *, PetscInt *);
1655: PETSC_EXTERN PetscErrorCode PetscCDAppendID(PetscCoarsenData *, PetscInt, PetscInt);
1656: PETSC_EXTERN PetscErrorCode PetscCDMoveAppend(PetscCoarsenData *, PetscInt, PetscInt);
1657: PETSC_EXTERN PetscErrorCode PetscCDAppendNode(PetscCoarsenData *, PetscInt, PetscCDIntNd *);
1658: PETSC_EXTERN PetscErrorCode PetscCDRemoveNextNode(PetscCoarsenData *, PetscInt, PetscCDIntNd *);
1659: PETSC_EXTERN PetscErrorCode PetscCDCountAt(const PetscCoarsenData *, PetscInt, PetscInt *);
1660: PETSC_EXTERN PetscErrorCode PetscCDIsEmptyAt(const PetscCoarsenData *, PetscInt, PetscBool *);
1661: PETSC_EXTERN PetscErrorCode PetscCDSetChunkSize(PetscCoarsenData *, PetscInt);
1662: PETSC_EXTERN PetscErrorCode PetscCDPrint(const PetscCoarsenData *, PetscInt, MPI_Comm);
1663: PETSC_EXTERN PetscErrorCode PetscCDGetNonemptyIS(PetscCoarsenData *, IS *);
1664: PETSC_EXTERN PetscErrorCode PetscCDGetMat(PetscCoarsenData *, Mat *);
1665: PETSC_EXTERN PetscErrorCode PetscCDSetMat(PetscCoarsenData *, Mat);
1666: PETSC_EXTERN PetscErrorCode PetscCDClearMat(PetscCoarsenData *);
1667: PETSC_EXTERN PetscErrorCode PetscCDRemoveAllAt(PetscCoarsenData *, PetscInt);
1668: PETSC_EXTERN PetscErrorCode PetscCDCount(const PetscCoarsenData *, PetscInt *_sz);

1670: PETSC_EXTERN PetscErrorCode PetscCDGetHeadPos(const PetscCoarsenData *, PetscInt, PetscCDIntNd **);
1671: PETSC_EXTERN PetscErrorCode PetscCDGetNextPos(const PetscCoarsenData *, PetscInt, PetscCDIntNd **);
1672: PETSC_EXTERN PetscErrorCode PetscCDGetASMBlocks(const PetscCoarsenData *, const PetscInt, PetscInt *, IS **);

1674: PETSC_SINGLE_LIBRARY_VISIBILITY_INTERNAL PetscErrorCode MatFDColoringApply_AIJ(Mat, MatFDColoring, Vec, void *);

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