Actual source code: matimpl.h

  1: #pragma once

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

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

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

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

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

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

232: #include <petscsys.h>

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

240: PETSC_EXTERN MatRootName MatRootNameList;

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

248: #define MAT_SHELL_NOT_ALLOWED (void *)-1

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

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

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

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

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

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

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

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

322: typedef struct _MatStashSpace *PetscMatStashSpace;

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

333: PETSC_EXTERN PetscErrorCode PetscMatStashSpaceGet(PetscInt, PetscInt, PetscMatStashSpace *);
334: PETSC_EXTERN PetscErrorCode PetscMatStashSpaceContiguous(PetscInt, PetscMatStashSpace *, PetscScalar *, PetscInt *, PetscInt *);
335: PETSC_EXTERN PetscErrorCode PetscMatStashSpaceDestroy(PetscMatStashSpace *);

337: typedef struct {
338:   PetscInt count;
339: } MatStashHeader;

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

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

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

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

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

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

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

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

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

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

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

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

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

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

532: PETSC_INTERN PetscErrorCode MatSetUp_Default(Mat);

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

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

545: /*
546:     Object for partitioning graphs
547: */

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

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

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

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

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

600: PETSC_EXTERN PetscErrorCode MatCoarsenMISKSetDistance(MatCoarsen, PetscInt);
601: PETSC_EXTERN PetscErrorCode MatCoarsenMISKGetDistance(MatCoarsen, PetscInt *);

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

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

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

630:     ncolors = 4;

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

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

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

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

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

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

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

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

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

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

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

743: /*
744:    Internal data structure for MATMPIDENSE
745: */
746: typedef struct {
747:   Mat A; /* local submatrix */

749:   /* The following variables are used for matrix assembly */
750:   PetscBool    donotstash;        /* Flag indicating if values should be stashed */
751:   MPI_Request *send_waits;        /* array of send requests */
752:   MPI_Request *recv_waits;        /* array of receive requests */
753:   PetscInt     nsends, nrecvs;    /* numbers of sends and receives */
754:   PetscScalar *svalues, *rvalues; /* sending and receiving data */
755:   PetscInt     rmax;              /* maximum message length */

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

762:   /* Support for MatDenseGetColumnVec and MatDenseGetSubMatrix */
763:   Mat                cmat;     /* matrix representation of a given subset of columns */
764:   Vec                cvec;     /* vector representation of a given column */
765:   const PetscScalar *ptrinuse; /* holds array to be restored (just a placeholder) */
766:   PetscInt           vecinuse; /* if cvec is in use (col = vecinuse-1) */
767:   PetscInt           matinuse; /* if cmat is in use (cbegin = matinuse-1) */
768: } Mat_MPIDense;

770: /*
771:    Checking zero pivot for LU, ILU preconditioners.
772: */
773: typedef struct {
774:   PetscInt    nshift, nshift_max;
775:   PetscReal   shift_amount, shift_lo, shift_hi, shift_top, shift_fraction;
776:   PetscBool   newshift;
777:   PetscReal   rs; /* active row sum of abs(off-diagonals) */
778:   PetscScalar pv; /* pivot of the active row */
779: } FactorShiftCtx;

781: PETSC_EXTERN PetscErrorCode MatTransposeCheckNonzeroState_Private(Mat, Mat);

783: /*
784:  Used by MatTranspose() and potentially other functions to track the matrix used in the generation of another matrix
785: */
786: typedef struct {
787:   PetscObjectId    id;
788:   PetscObjectState state;
789:   PetscObjectState nonzerostate;
790: } MatParentState;

792: PETSC_EXTERN PetscErrorCode MatFactorDumpMatrix(Mat);
793: PETSC_INTERN PetscErrorCode MatSetBlockSizes_Default(Mat, PetscInt, PetscInt);

795: PETSC_SINGLE_LIBRARY_INTERN PetscErrorCode MatShift_Basic(Mat, PetscScalar);

797: static inline PetscErrorCode MatPivotCheck_nz(PETSC_UNUSED Mat mat, const MatFactorInfo *info, FactorShiftCtx *sctx, PETSC_UNUSED PetscInt row)
798: {
799:   PetscReal _rs   = sctx->rs;
800:   PetscReal _zero = info->zeropivot * _rs;

802:   PetscFunctionBegin;
803:   if (PetscAbsScalar(sctx->pv) <= _zero && !PetscIsNanScalar(sctx->pv)) {
804:     /* force |diag| > zeropivot*rs */
805:     if (!sctx->nshift) sctx->shift_amount = info->shiftamount;
806:     else sctx->shift_amount *= 2.0;
807:     sctx->newshift = PETSC_TRUE;
808:     (sctx->nshift)++;
809:   } else {
810:     sctx->newshift = PETSC_FALSE;
811:   }
812:   PetscFunctionReturn(PETSC_SUCCESS);
813: }

815: static inline PetscErrorCode MatPivotCheck_pd(PETSC_UNUSED Mat mat, const MatFactorInfo *info, FactorShiftCtx *sctx, PETSC_UNUSED PetscInt row)
816: {
817:   PetscReal _rs   = sctx->rs;
818:   PetscReal _zero = info->zeropivot * _rs;

820:   PetscFunctionBegin;
821:   if (PetscRealPart(sctx->pv) <= _zero && !PetscIsNanScalar(sctx->pv)) {
822:     /* force matfactor to be diagonally dominant */
823:     if (sctx->nshift == sctx->nshift_max) {
824:       sctx->shift_fraction = sctx->shift_hi;
825:     } else {
826:       sctx->shift_lo       = sctx->shift_fraction;
827:       sctx->shift_fraction = (sctx->shift_hi + sctx->shift_lo) / (PetscReal)2.;
828:     }
829:     sctx->shift_amount = sctx->shift_fraction * sctx->shift_top;
830:     sctx->nshift++;
831:     sctx->newshift = PETSC_TRUE;
832:   } else {
833:     sctx->newshift = PETSC_FALSE;
834:   }
835:   PetscFunctionReturn(PETSC_SUCCESS);
836: }

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

842:   PetscFunctionBegin;
843:   if (PetscAbsScalar(sctx->pv) <= _zero && !PetscIsNanScalar(sctx->pv)) {
844:     sctx->pv += info->shiftamount;
845:     sctx->shift_amount = 0.0;
846:     sctx->nshift++;
847:   }
848:   sctx->newshift = PETSC_FALSE;
849:   PetscFunctionReturn(PETSC_SUCCESS);
850: }

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

856:   PetscFunctionBegin;
857:   sctx->newshift = PETSC_FALSE;
858:   if (PetscAbsScalar(sctx->pv) <= _zero && !PetscIsNanScalar(sctx->pv)) {
859:     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);
860:     PetscCall(PetscInfo(mat, "Detected zero pivot in factorization in row %" PetscInt_FMT " value %g tolerance %g\n", row, (double)PetscAbsScalar(sctx->pv), (double)_zero));
861:     fact->factorerrortype             = MAT_FACTOR_NUMERIC_ZEROPIVOT;
862:     fact->factorerror_zeropivot_value = PetscAbsScalar(sctx->pv);
863:     fact->factorerror_zeropivot_row   = row;
864:   }
865:   PetscFunctionReturn(PETSC_SUCCESS);
866: }

868: static inline PetscErrorCode MatPivotCheck(Mat fact, Mat mat, const MatFactorInfo *info, FactorShiftCtx *sctx, PetscInt row)
869: {
870:   PetscFunctionBegin;
871:   if (info->shifttype == (PetscReal)MAT_SHIFT_NONZERO) PetscCall(MatPivotCheck_nz(mat, info, sctx, row));
872:   else if (info->shifttype == (PetscReal)MAT_SHIFT_POSITIVE_DEFINITE) PetscCall(MatPivotCheck_pd(mat, info, sctx, row));
873:   else if (info->shifttype == (PetscReal)MAT_SHIFT_INBLOCKS) PetscCall(MatPivotCheck_inblocks(mat, info, sctx, row));
874:   else PetscCall(MatPivotCheck_none(fact, mat, info, sctx, row));
875:   PetscFunctionReturn(PETSC_SUCCESS);
876: }

878: #include <petscbt.h>
879: /*
880:   Create and initialize a linked list
881:   Input Parameters:
882:     idx_start - starting index of the list
883:     lnk_max   - max value of lnk indicating the end of the list
884:     nlnk      - max length of the list
885:   Output Parameters:
886:     lnk       - list initialized
887:     bt        - PetscBT (bitarray) with all bits set to false
888:     lnk_empty - flg indicating the list is empty
889: */
890: #define PetscLLCreate(idx_start, lnk_max, nlnk, lnk, bt) ((PetscErrorCode)(PetscMalloc1(nlnk, &lnk) || PetscBTCreate(nlnk, &(bt)) || (lnk[idx_start] = lnk_max, PETSC_SUCCESS)))

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

894: 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)
895: {
896:   PetscInt location;

898:   PetscFunctionBegin;
899:   /* start from the beginning if entry < previous entry */
900:   if (!assume_sorted && k && entry < *lnkdata) *lnkdata = idx_start;
901:   /* search for insertion location */
902:   do {
903:     location = *lnkdata;
904:     *lnkdata = lnk[location];
905:   } while (entry > *lnkdata);
906:   /* insertion location is found, add entry into lnk */
907:   lnk[location] = entry;
908:   lnk[entry]    = *lnkdata;
909:   ++(*nlnk);
910:   *lnkdata = entry; /* next search starts from here if next_entry > entry */
911:   PetscFunctionReturn(PETSC_SUCCESS);
912: }

914: 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)
915: {
916:   PetscFunctionBegin;
917:   *nlnk = 0;
918:   for (PetscInt k = 0, lnkdata = idx_start; k < nidx; ++k) {
919:     const PetscInt entry = indices[k];

921:     if (!PetscBTLookupSet(bt, entry)) PetscCall(PetscLLInsertLocation_Private(assume_sorted, k, idx_start, entry, nlnk, &lnkdata, lnk));
922:   }
923:   PetscFunctionReturn(PETSC_SUCCESS);
924: }

926: /*
927:   Add an index set into a sorted linked list
928:   Input Parameters:
929:     nidx      - number of input indices
930:     indices   - integer array
931:     idx_start - starting index of the list
932:     lnk       - linked list(an integer array) that is created
933:     bt        - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
934:   output Parameters:
935:     nlnk      - number of newly added indices
936:     lnk       - the sorted(increasing order) linked list containing new and non-redundate entries from indices
937:     bt        - updated PetscBT (bitarray)
938: */
939: static inline PetscErrorCode PetscLLAdd(PetscInt nidx, const PetscInt *PETSC_RESTRICT indices, PetscInt idx_start, PetscInt *PETSC_RESTRICT nlnk, PetscInt *PETSC_RESTRICT lnk, PetscBT bt)
940: {
941:   PetscFunctionBegin;
942:   PetscCall(PetscLLAdd_Private(nidx, indices, idx_start, nlnk, lnk, bt, PETSC_FALSE));
943:   PetscFunctionReturn(PETSC_SUCCESS);
944: }

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

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

987:     if (!PetscBTLookupSet(bt, entry)) PetscCall(PetscLLInsertLocation_Private(PETSC_FALSE, k, idx_start, entry, nlnk, &lnkdata, lnk));
988:   }
989:   PetscFunctionReturn(PETSC_SUCCESS);
990: }

992: #if 0
993: /* this appears to be unused? */
994: static inline PetscErrorCode PetscLLAddSorted_new(PetscInt nidx, PetscInt *indices, PetscInt idx_start, PetscBool *lnk_empty, PetscInt *nlnk, PetscInt *lnk, PetscBT bt)
995: {
996:   PetscInt lnkdata = idx_start;

998:   PetscFunctionBegin;
999:   if (*lnk_empty) {
1000:     for (PetscInt k = 0; k < nidx; ++k) {
1001:       const PetscInt entry = indices[k], location = lnkdata;

1003:       PetscCall(PetscBTSet(bt,entry)); /* mark the new entry */
1004:       lnkdata       = lnk[location];
1005:       /* insertion location is found, add entry into lnk */
1006:       lnk[location] = entry;
1007:       lnk[entry]    = lnkdata;
1008:       lnkdata       = entry; /* next search starts from here */
1009:     }
1010:     /* lnk[indices[nidx-1]] = lnk[idx_start];
1011:        lnk[idx_start]       = indices[0];
1012:        PetscCall(PetscBTSet(bt,indices[0]));
1013:        for (_k=1; _k<nidx; _k++) {
1014:        PetscCall(PetscBTSet(bt,indices[_k]));
1015:        lnk[indices[_k-1]] = indices[_k];
1016:        }
1017:     */
1018:     *nlnk      = nidx;
1019:     *lnk_empty = PETSC_FALSE;
1020:   } else {
1021:     *nlnk = 0;
1022:     for (PetscInt k = 0; k < nidx; ++k) {
1023:       const PetscInt entry = indices[k];

1025:       if (!PetscBTLookupSet(bt,entry)) PetscCall(PetscLLInsertLocation_Private(PETSC_TRUE,k,idx_start,entry,nlnk,&lnkdata,lnk));
1026:     }
1027:   }
1028:   PetscFunctionReturn(PETSC_SUCCESS);
1029: }
1030: #endif

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

1055:   PetscFunctionBegin;
1056:   *nlnk = 0;
1057:   for (PetscInt k = 0, lnkdata = idx_start; k < nidx; ++k) {
1058:     const PetscInt entry = indices[k];

1060:     ++nzbd;
1061:     if (entry == diag) im[idx_start] = nzbd;
1062:     if (!PetscBTLookupSet(bt, entry)) PetscCall(PetscLLInsertLocation_Private(PETSC_TRUE, k, idx_start, entry, nlnk, &lnkdata, lnk));
1063:   }
1064:   PetscFunctionReturn(PETSC_SUCCESS);
1065: }

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

1092: /*
1093:   Free memories used by the list
1094: */
1095: #define PetscLLDestroy(lnk, bt) ((PetscErrorCode)(PetscFree(lnk) || PetscBTDestroy(&(bt))))

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

1112: 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)
1113: {
1114:   PetscFunctionBegin;
1115:   PetscCall(PetscLLInsertLocation_Private(assume_sorted, k, idx_start, entry, nlnk, lnkdata, lnk));
1116:   lnklvl[entry] = newval;
1117:   PetscFunctionReturn(PETSC_SUCCESS);
1118: }

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

1143:     if (!PetscBTLookupSet(bt, entry)) PetscCall(PetscIncompleteLLInsertLocation_Private(PETSC_FALSE, k, idx_start, entry, nlnk, &lnkdata, lnk, lnklvl, 0));
1144:   }
1145:   PetscFunctionReturn(PETSC_SUCCESS);
1146: }

1148: 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)
1149: {
1150:   PetscFunctionBegin;
1151:   *nlnk = 0;
1152:   for (PetscInt k = 0, lnkdata = idx_start; k < nidx; ++k) {
1153:     const PetscInt incrlev = idxlvl[k] + prow_offset + 1;

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

1158:       if (!PetscBTLookupSet(bt, entry)) PetscCall(PetscIncompleteLLInsertLocation_Private(assume_sorted, k, idx_start, entry, nlnk, &lnkdata, lnk, lnklvl, incrlev));
1159:       else if (lnklvl[entry] > incrlev) lnklvl[entry] = incrlev; /* existing entry */
1160:     }
1161:   }
1162:   PetscFunctionReturn(PETSC_SUCCESS);
1163: }

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

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

1210:   Note: the level of factor(i,j) is set as lvl(i,j) = min{ lvl(i,j), lvl(i,prow)+lvl(prow,j)+1)
1211:         where idx = non-zero columns of U(prow,prow+1:n-1), prow<i
1212: */
1213: 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)
1214: {
1215:   PetscFunctionBegin;
1216:   PetscCall(PetscIncompleteLLAdd_Private(nidx, idx, level, idxlvl, idx_start, nlnk, lnk, lnklvl, bt, lnklvl[prow], PETSC_TRUE));
1217:   PetscFunctionReturn(PETSC_SUCCESS);
1218: }

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

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

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

1297: /*
1298:   Free memories used by the list
1299: */
1300: #define PetscIncompleteLLDestroy(lnk, bt) ((PetscErrorCode)(PetscFree(lnk) || PetscBTDestroy(&(bt))))

1302: #if !defined(PETSC_CLANG_STATIC_ANALYZER)
1303:   #define MatCheckSameLocalSize(A, ar1, B, ar2) \
1304:     do { \
1305:       PetscCheckSameComm(A, ar1, B, ar2); \
1306:       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, \
1307:                  (A)->rmap->n, (A)->cmap->n, ar2, (B)->rmap->n, (B)->cmap->n); \
1308:     } while (0)
1309:   #define MatCheckSameSize(A, ar1, B, ar2) \
1310:     do { \
1311:       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, \
1312:                  (A)->rmap->N, (A)->cmap->N, ar2, (B)->rmap->N, (B)->cmap->N); \
1313:       MatCheckSameLocalSize(A, ar1, B, ar2); \
1314:     } while (0)
1315: #else
1316: template <typename Tm>
1317: extern void MatCheckSameLocalSize(Tm, int, Tm, int);
1318: template <typename Tm>
1319: extern void MatCheckSameSize(Tm, int, Tm, int);
1320: #endif

1322: #define VecCheckMatCompatible(M, x, ar1, b, ar2) \
1323:   do { \
1324:     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, \
1325:                (M)->cmap->N); \
1326:     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, \
1327:                (M)->rmap->N); \
1328:   } while (0)

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

1337:       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
1338:       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
1339:       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
1340:       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.
1341:       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
1342:       to each other so memory access is much better than using the big array.

1344:   Example:
1345:      nlnk_max=5, lnk_max=36:
1346:      Initial list: [0, 0 | 36, 2 | 0, 0 | 0, 0 | 0, 0 | 0, 0 | 0, 0]
1347:      here, head_node has index 2 with value lnk[2]=lnk_max=36,
1348:            0-th entry is used to store the number of entries in the list,
1349:      The initial lnk represents head -> tail(marked by 36) with number of entries = lnk[0]=0.

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

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

1359:   Input Parameters:
1360:     nlnk_max  - max length of the list
1361:     lnk_max   - max value of the entries
1362:   Output Parameters:
1363:     lnk       - list created and initialized
1364:     bt        - PetscBT (bitarray) with all bits set to false. Note: bt has size lnk_max, not nln_max!
1365: */
1366: static inline PetscErrorCode PetscLLCondensedCreate(PetscInt nlnk_max, PetscInt lnk_max, PetscInt **lnk, PetscBT *bt)
1367: {
1368:   PetscInt *llnk, lsize = 0;

1370:   PetscFunctionBegin;
1371:   PetscCall(PetscIntMultError(2, nlnk_max + 2, &lsize));
1372:   PetscCall(PetscMalloc1(lsize, lnk));
1373:   PetscCall(PetscBTCreate(lnk_max, bt));
1374:   llnk    = *lnk;
1375:   llnk[0] = 0;       /* number of entries on the list */
1376:   llnk[2] = lnk_max; /* value in the head node */
1377:   llnk[3] = 2;       /* next for the head node */
1378:   PetscFunctionReturn(PETSC_SUCCESS);
1379: }

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

1397:   PetscFunctionBegin;
1398:   for (PetscInt k = 0; k < nidx; k++) {
1399:     const PetscInt entry = indices[k];
1400:     if (!PetscBTLookupSet(bt, entry)) { /* new entry */
1401:       PetscInt next, lnkdata;

1403:       /* search for insertion location */
1404:       do {
1405:         next     = location + 1;  /* link from previous node to next node */
1406:         location = lnk[next];     /* idx of next node */
1407:         lnkdata  = lnk[location]; /* value of next node */
1408:       } while (entry > lnkdata);
1409:       /* insertion location is found, add entry into lnk */
1410:       const PetscInt newnode = 2 * (nlnk + 2); /* index for this new node */
1411:       lnk[next]              = newnode;        /* connect previous node to the new node */
1412:       lnk[newnode]           = entry;          /* set value of the new node */
1413:       lnk[newnode + 1]       = location;       /* connect new node to next node */
1414:       location               = newnode;        /* next search starts from the new node */
1415:       nlnk++;
1416:     }
1417:   }
1418:   lnk[0] = nlnk; /* number of entries in the list */
1419:   PetscFunctionReturn(PETSC_SUCCESS);
1420: }

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

1427:   PetscFunctionBegin;
1428:   for (PetscInt k = 0; k < nlnk; k++) {
1429:     indices[k] = lnk[next];
1430:     next       = lnk[next + 1];
1431:     PetscCall(PetscBTClear(bt, indices[k]));
1432:   }
1433:   lnk[0] = 0;       /* num of entries on the list */
1434:   lnk[2] = lnk_max; /* initialize head node */
1435:   lnk[3] = 2;       /* head node */
1436:   PetscFunctionReturn(PETSC_SUCCESS);
1437: }

1439: static inline PetscErrorCode PetscLLCondensedView(PetscInt *lnk)
1440: {
1441:   PetscFunctionBegin;
1442:   PetscCall(PetscPrintf(PETSC_COMM_SELF, "LLCondensed of size %" PetscInt_FMT ", (val,  next)\n", lnk[0]));
1443:   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]));
1444:   PetscFunctionReturn(PETSC_SUCCESS);
1445: }

1447: /*
1448:   Free memories used by the list
1449: */
1450: static inline PetscErrorCode PetscLLCondensedDestroy(PetscInt *lnk, PetscBT bt)
1451: {
1452:   PetscFunctionBegin;
1453:   PetscCall(PetscFree(lnk));
1454:   PetscCall(PetscBTDestroy(&bt));
1455:   PetscFunctionReturn(PETSC_SUCCESS);
1456: }

1458: /*
1459:  Same as PetscLLCondensedCreate(), but does not use non-scalable O(lnk_max) bitarray
1460:   Input Parameters:
1461:     nlnk_max  - max length of the list
1462:   Output Parameters:
1463:     lnk       - list created and initialized
1464: */
1465: static inline PetscErrorCode PetscLLCondensedCreate_Scalable(PetscInt nlnk_max, PetscInt **lnk)
1466: {
1467:   PetscInt *llnk, lsize = 0;

1469:   PetscFunctionBegin;
1470:   PetscCall(PetscIntMultError(2, nlnk_max + 2, &lsize));
1471:   PetscCall(PetscMalloc1(lsize, lnk));
1472:   llnk    = *lnk;
1473:   llnk[0] = 0;             /* number of entries on the list */
1474:   llnk[2] = PETSC_INT_MAX; /* value in the head node */
1475:   llnk[3] = 2;             /* next for the head node */
1476:   PetscFunctionReturn(PETSC_SUCCESS);
1477: }

1479: static inline PetscErrorCode PetscLLCondensedExpand_Scalable(PetscInt nlnk_max, PetscInt **lnk)
1480: {
1481:   PetscInt lsize = 0;

1483:   PetscFunctionBegin;
1484:   PetscCall(PetscIntMultError(2, nlnk_max + 2, &lsize));
1485:   PetscCall(PetscRealloc(lsize * sizeof(PetscInt), lnk));
1486:   PetscFunctionReturn(PETSC_SUCCESS);
1487: }

1489: static inline PetscErrorCode PetscLLCondensedAddSorted_Scalable(PetscInt nidx, const PetscInt indices[], PetscInt lnk[])
1490: {
1491:   PetscInt location = 2;      /* head */
1492:   PetscInt nlnk     = lnk[0]; /* num of entries on the input lnk */

1494:   for (PetscInt k = 0; k < nidx; k++) {
1495:     const PetscInt entry = indices[k];
1496:     PetscInt       next, lnkdata;

1498:     /* search for insertion location */
1499:     do {
1500:       next     = location + 1;  /* link from previous node to next node */
1501:       location = lnk[next];     /* idx of next node */
1502:       lnkdata  = lnk[location]; /* value of next node */
1503:     } while (entry > lnkdata);
1504:     if (entry < lnkdata) {
1505:       /* insertion location is found, add entry into lnk */
1506:       const PetscInt newnode = 2 * (nlnk + 2); /* index for this new node */
1507:       lnk[next]              = newnode;        /* connect previous node to the new node */
1508:       lnk[newnode]           = entry;          /* set value of the new node */
1509:       lnk[newnode + 1]       = location;       /* connect new node to next node */
1510:       location               = newnode;        /* next search starts from the new node */
1511:       nlnk++;
1512:     }
1513:   }
1514:   lnk[0] = nlnk; /* number of entries in the list */
1515:   return PETSC_SUCCESS;
1516: }

1518: static inline PetscErrorCode PetscLLCondensedClean_Scalable(PETSC_UNUSED PetscInt nidx, PetscInt *indices, PetscInt *lnk)
1519: {
1520:   const PetscInt nlnk = lnk[0];
1521:   PetscInt       next = lnk[3]; /* head node */

1523:   for (PetscInt k = 0; k < nlnk; k++) {
1524:     indices[k] = lnk[next];
1525:     next       = lnk[next + 1];
1526:   }
1527:   lnk[0] = 0; /* num of entries on the list */
1528:   lnk[3] = 2; /* head node */
1529:   return PETSC_SUCCESS;
1530: }

1532: static inline PetscErrorCode PetscLLCondensedDestroy_Scalable(PetscInt *lnk)
1533: {
1534:   return PetscFree(lnk);
1535: }

1537: /*
1538:       lnk[0]   number of links
1539:       lnk[1]   number of entries
1540:       lnk[3n]  value
1541:       lnk[3n+1] len
1542:       lnk[3n+2] link to next value

1544:       The next three are always the first link

1546:       lnk[3]    PETSC_INT_MIN+1
1547:       lnk[4]    1
1548:       lnk[5]    link to first real entry

1550:       The next three are always the last link

1552:       lnk[6]    PETSC_INT_MAX - 1
1553:       lnk[7]    1
1554:       lnk[8]    next valid link (this is the same as lnk[0] but without the decreases)
1555: */

1557: static inline PetscErrorCode PetscLLCondensedCreate_fast(PetscInt nlnk_max, PetscInt **lnk)
1558: {
1559:   PetscInt *llnk;
1560:   PetscInt  lsize = 0;

1562:   PetscFunctionBegin;
1563:   PetscCall(PetscIntMultError(3, nlnk_max + 3, &lsize));
1564:   PetscCall(PetscMalloc1(lsize, lnk));
1565:   llnk    = *lnk;
1566:   llnk[0] = 0;                 /* nlnk: number of entries on the list */
1567:   llnk[1] = 0;                 /* number of integer entries represented in list */
1568:   llnk[3] = PETSC_INT_MIN + 1; /* value in the first node */
1569:   llnk[4] = 1;                 /* count for the first node */
1570:   llnk[5] = 6;                 /* next for the first node */
1571:   llnk[6] = PETSC_INT_MAX - 1; /* value in the last node */
1572:   llnk[7] = 1;                 /* count for the last node */
1573:   llnk[8] = 0;                 /* next valid node to be used */
1574:   PetscFunctionReturn(PETSC_SUCCESS);
1575: }

1577: static inline PetscErrorCode PetscLLCondensedAddSorted_fast(PetscInt nidx, const PetscInt indices[], PetscInt lnk[])
1578: {
1579:   for (PetscInt k = 0, prev = 3 /* first value */; k < nidx; k++) {
1580:     const PetscInt entry = indices[k];
1581:     PetscInt       next  = lnk[prev + 2];

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

1621: static inline PetscErrorCode PetscLLCondensedClean_fast(PETSC_UNUSED PetscInt nidx, PetscInt *indices, PetscInt *lnk)
1622: {
1623:   const PetscInt nlnk = lnk[0];
1624:   PetscInt       next = lnk[5]; /* first node */

1626:   for (PetscInt k = 0, cnt = 0; k < nlnk; k++) {
1627:     for (PetscInt j = 0; j < lnk[next + 1]; j++) indices[cnt++] = lnk[next] + j;
1628:     next = lnk[next + 2];
1629:   }
1630:   lnk[0] = 0;                 /* nlnk: number of links */
1631:   lnk[1] = 0;                 /* number of integer entries represented in list */
1632:   lnk[3] = PETSC_INT_MIN + 1; /* value in the first node */
1633:   lnk[4] = 1;                 /* count for the first node */
1634:   lnk[5] = 6;                 /* next for the first node */
1635:   lnk[6] = PETSC_INT_MAX - 1; /* value in the last node */
1636:   lnk[7] = 1;                 /* count for the last node */
1637:   lnk[8] = 0;                 /* next valid location to make link */
1638:   return PETSC_SUCCESS;
1639: }

1641: static inline PetscErrorCode PetscLLCondensedView_fast(const PetscInt *lnk)
1642: {
1643:   const PetscInt nlnk = lnk[0];
1644:   PetscInt       next = lnk[5]; /* first node */

1646:   for (PetscInt k = 0; k < nlnk; k++) {
1647: #if 0 /* Debugging code */
1648:     printf("%d value %d len %d next %d\n", next, lnk[next], lnk[next + 1], lnk[next + 2]);
1649: #endif
1650:     next = lnk[next + 2];
1651:   }
1652:   return PETSC_SUCCESS;
1653: }

1655: static inline PetscErrorCode PetscLLCondensedDestroy_fast(PetscInt *lnk)
1656: {
1657:   return PetscFree(lnk);
1658: }

1660: PETSC_EXTERN PetscErrorCode PetscCDCreate(PetscInt, PetscCoarsenData **);
1661: PETSC_EXTERN PetscErrorCode PetscCDDestroy(PetscCoarsenData *);
1662: PETSC_EXTERN PetscErrorCode PetscCDIntNdSetID(PetscCDIntNd *, PetscInt);
1663: PETSC_EXTERN PetscErrorCode PetscCDIntNdGetID(const PetscCDIntNd *, PetscInt *);
1664: PETSC_EXTERN PetscErrorCode PetscCDAppendID(PetscCoarsenData *, PetscInt, PetscInt);
1665: PETSC_EXTERN PetscErrorCode PetscCDMoveAppend(PetscCoarsenData *, PetscInt, PetscInt);
1666: PETSC_EXTERN PetscErrorCode PetscCDAppendNode(PetscCoarsenData *, PetscInt, PetscCDIntNd *);
1667: PETSC_EXTERN PetscErrorCode PetscCDRemoveNextNode(PetscCoarsenData *, PetscInt, PetscCDIntNd *);
1668: PETSC_EXTERN PetscErrorCode PetscCDCountAt(const PetscCoarsenData *, PetscInt, PetscInt *);
1669: PETSC_EXTERN PetscErrorCode PetscCDIsEmptyAt(const PetscCoarsenData *, PetscInt, PetscBool *);
1670: PETSC_EXTERN PetscErrorCode PetscCDSetChunkSize(PetscCoarsenData *, PetscInt);
1671: PETSC_EXTERN PetscErrorCode PetscCDPrint(const PetscCoarsenData *, PetscInt, MPI_Comm);
1672: PETSC_EXTERN PetscErrorCode PetscCDGetNonemptyIS(PetscCoarsenData *, IS *);
1673: PETSC_EXTERN PetscErrorCode PetscCDGetMat(PetscCoarsenData *, Mat *);
1674: PETSC_EXTERN PetscErrorCode PetscCDSetMat(PetscCoarsenData *, Mat);
1675: PETSC_EXTERN PetscErrorCode PetscCDClearMat(PetscCoarsenData *);
1676: PETSC_EXTERN PetscErrorCode PetscCDRemoveAllAt(PetscCoarsenData *, PetscInt);
1677: PETSC_EXTERN PetscErrorCode PetscCDCount(const PetscCoarsenData *, PetscInt *_sz);

1679: PETSC_EXTERN PetscErrorCode PetscCDGetHeadPos(const PetscCoarsenData *, PetscInt, PetscCDIntNd **);
1680: PETSC_EXTERN PetscErrorCode PetscCDGetNextPos(const PetscCoarsenData *, PetscInt, PetscCDIntNd **);
1681: PETSC_EXTERN PetscErrorCode PetscCDGetASMBlocks(const PetscCoarsenData *, const PetscInt, PetscInt *, IS **);

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

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