Actual source code: matimpl.h

  1: #pragma once

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

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

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

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

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

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

228: #include <petscsys.h>

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

236: PETSC_EXTERN MatRootName MatRootNameList;

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

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

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

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

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

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

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

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

306: /*
307:   The stash is used to temporarily store inserted matrix values that
308:   belong to another processor. During the assembly phase the stashed
309:   values are moved to the correct processor and
310: */

312: typedef struct _MatStashSpace *PetscMatStashSpace;

314: struct _MatStashSpace {
315:   PetscMatStashSpace next;
316:   PetscScalar       *space_head, *val;
317:   PetscInt          *idx, *idy;
318:   PetscInt           total_space_size;
319:   PetscInt           local_used;
320:   PetscInt           local_remaining;
321: };

323: PETSC_EXTERN PetscErrorCode PetscMatStashSpaceGet(PetscInt, PetscInt, PetscMatStashSpace *);
324: PETSC_EXTERN PetscErrorCode PetscMatStashSpaceContiguous(PetscInt, PetscMatStashSpace *, PetscScalar *, PetscInt *, PetscInt *);
325: PETSC_EXTERN PetscErrorCode PetscMatStashSpaceDestroy(PetscMatStashSpace *);

327: typedef struct {
328:   PetscInt count;
329: } MatStashHeader;

331: typedef struct {
332:   void    *buffer; /* Of type blocktype, dynamically constructed  */
333:   PetscInt count;
334:   char     pending;
335: } MatStashFrame;

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

347:   PetscErrorCode (*ScatterBegin)(Mat, MatStash *, PetscInt *);
348:   PetscErrorCode (*ScatterGetMesg)(MatStash *, PetscMPIInt *, PetscInt **, PetscInt **, PetscScalar **, PetscInt *);
349:   PetscErrorCode (*ScatterEnd)(MatStash *);
350:   PetscErrorCode (*ScatterDestroy)(MatStash *);

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

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

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

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

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

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

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

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

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

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

517: PETSC_INTERN PetscErrorCode MatAXPY_Basic(Mat, PetscScalar, Mat, MatStructure);
518: PETSC_INTERN PetscErrorCode MatAXPY_BasicWithPreallocation(Mat, Mat, PetscScalar, Mat, MatStructure);
519: PETSC_INTERN PetscErrorCode MatAXPY_Basic_Preallocate(Mat, Mat, Mat *);
520: PETSC_INTERN PetscErrorCode MatAXPY_Dense_Nest(Mat, PetscScalar, Mat);

522: PETSC_INTERN PetscErrorCode MatSetUp_Default(Mat);

524: /*
525:     Utility for MatZeroRows
526: */
527: PETSC_INTERN PetscErrorCode MatZeroRowsMapLocal_Private(Mat, PetscInt, const PetscInt *, PetscInt *, PetscInt **);

529: /*
530:     Utility for MatView/MatLoad
531: */
532: PETSC_INTERN PetscErrorCode MatView_Binary_BlockSizes(Mat, PetscViewer);
533: PETSC_INTERN PetscErrorCode MatLoad_Binary_BlockSizes(Mat, PetscViewer);

535: /*
536:     Object for partitioning graphs
537: */

539: typedef struct _MatPartitioningOps *MatPartitioningOps;
540: struct _MatPartitioningOps {
541:   PetscErrorCode (*apply)(MatPartitioning, IS *);
542:   PetscErrorCode (*applynd)(MatPartitioning, IS *);
543:   PetscErrorCode (*setfromoptions)(MatPartitioning, PetscOptionItems *);
544:   PetscErrorCode (*destroy)(MatPartitioning);
545:   PetscErrorCode (*view)(MatPartitioning, PetscViewer);
546:   PetscErrorCode (*improve)(MatPartitioning, IS *);
547: };

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

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

564: /*
565:     Object for coarsen graphs
566: */
567: typedef struct _MatCoarsenOps *MatCoarsenOps;
568: struct _MatCoarsenOps {
569:   PetscErrorCode (*apply)(MatCoarsen);
570:   PetscErrorCode (*setfromoptions)(MatCoarsen, PetscOptionItems *);
571:   PetscErrorCode (*destroy)(MatCoarsen);
572:   PetscErrorCode (*view)(MatCoarsen, PetscViewer);
573: };

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

590: PETSC_EXTERN PetscErrorCode MatCoarsenMISKSetDistance(MatCoarsen, PetscInt);
591: PETSC_EXTERN PetscErrorCode MatCoarsenMISKGetDistance(MatCoarsen, PetscInt *);

593: /*
594:     Used in aijdevice.h
595: */
596: typedef struct {
597:   PetscInt    *i;
598:   PetscInt    *j;
599:   PetscScalar *a;
600:   PetscInt     n;
601:   PetscInt     ignorezeroentries;
602: } PetscCSRDataStructure;

604: /*
605:     MatFDColoring is used to compute Jacobian matrices efficiently
606:   via coloring. The data structure is explained below in an example.

608:    Color =   0    1     0    2   |   2      3       0
609:    ---------------------------------------------------
610:             00   01              |          05
611:             10   11              |   14     15               Processor  0
612:                        22    23  |          25
613:                        32    33  |
614:    ===================================================
615:                                  |   44     45     46
616:             50                   |          55               Processor 1
617:                                  |   64            66
618:    ---------------------------------------------------

620:     ncolors = 4;

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

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

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

639: */
640: typedef struct {
641:   PetscInt     row;
642:   PetscInt     col;
643:   PetscScalar *valaddr; /* address of value */
644: } MatEntry;

646: typedef struct {
647:   PetscInt     row;
648:   PetscScalar *valaddr; /* address of value */
649: } MatEntry2;

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

680: typedef struct _MatColoringOps *MatColoringOps;
681: struct _MatColoringOps {
682:   PetscErrorCode (*destroy)(MatColoring);
683:   PetscErrorCode (*setfromoptions)(MatColoring, PetscOptionItems *);
684:   PetscErrorCode (*view)(MatColoring, PetscViewer);
685:   PetscErrorCode (*apply)(MatColoring, ISColoring *);
686:   PetscErrorCode (*weights)(MatColoring, PetscReal **, PetscInt **);
687: };

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

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

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

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

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

744: PETSC_EXTERN PetscErrorCode MatTransposeCheckNonzeroState_Private(Mat, Mat);

746: /*
747:  Used by MatTranspose() and potentially other functions to track the matrix used in the generation of another matrix
748: */
749: typedef struct {
750:   PetscObjectId    id;
751:   PetscObjectState state;
752:   PetscObjectState nonzerostate;
753: } MatParentState;

755: PETSC_EXTERN PetscErrorCode MatFactorDumpMatrix(Mat);
756: PETSC_INTERN PetscErrorCode MatSetBlockSizes_Default(Mat, PetscInt, PetscInt);

758: PETSC_SINGLE_LIBRARY_INTERN PetscErrorCode MatShift_Basic(Mat, PetscScalar);

760: static inline PetscErrorCode MatPivotCheck_nz(PETSC_UNUSED Mat mat, const MatFactorInfo *info, FactorShiftCtx *sctx, PETSC_UNUSED PetscInt row)
761: {
762:   PetscReal _rs   = sctx->rs;
763:   PetscReal _zero = info->zeropivot * _rs;

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

778: static inline PetscErrorCode MatPivotCheck_pd(PETSC_UNUSED Mat mat, const MatFactorInfo *info, FactorShiftCtx *sctx, PETSC_UNUSED PetscInt row)
779: {
780:   PetscReal _rs   = sctx->rs;
781:   PetscReal _zero = info->zeropivot * _rs;

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

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

805:   PetscFunctionBegin;
806:   if (PetscAbsScalar(sctx->pv) <= _zero && !PetscIsNanScalar(sctx->pv)) {
807:     sctx->pv += info->shiftamount;
808:     sctx->shift_amount = 0.0;
809:     sctx->nshift++;
810:   }
811:   sctx->newshift = PETSC_FALSE;
812:   PetscFunctionReturn(PETSC_SUCCESS);
813: }

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

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

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

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

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

857: 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)
858: {
859:   PetscInt location;

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

877: 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)
878: {
879:   PetscFunctionBegin;
880:   *nlnk = 0;
881:   for (PetscInt k = 0, lnkdata = idx_start; k < nidx; ++k) {
882:     const PetscInt entry = indices[k];

884:     if (!PetscBTLookupSet(bt, entry)) PetscCall(PetscLLInsertLocation_Private(assume_sorted, k, idx_start, entry, nlnk, &lnkdata, lnk));
885:   }
886:   PetscFunctionReturn(PETSC_SUCCESS);
887: }

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

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

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

950:     if (!PetscBTLookupSet(bt, entry)) PetscCall(PetscLLInsertLocation_Private(PETSC_FALSE, k, idx_start, entry, nlnk, &lnkdata, lnk));
951:   }
952:   PetscFunctionReturn(PETSC_SUCCESS);
953: }

955: #if 0
956: /* this appears to be unused? */
957: static inline PetscErrorCode PetscLLAddSorted_new(PetscInt nidx, PetscInt *indices, PetscInt idx_start, PetscBool *lnk_empty, PetscInt *nlnk, PetscInt *lnk, PetscBT bt)
958: {
959:   PetscInt lnkdata = idx_start;

961:   PetscFunctionBegin;
962:   if (*lnk_empty) {
963:     for (PetscInt k = 0; k < nidx; ++k) {
964:       const PetscInt entry = indices[k], location = lnkdata;

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

988:       if (!PetscBTLookupSet(bt,entry)) PetscCall(PetscLLInsertLocation_Private(PETSC_TRUE,k,idx_start,entry,nlnk,&lnkdata,lnk));
989:     }
990:   }
991:   PetscFunctionReturn(PETSC_SUCCESS);
992: }
993: #endif

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

1018:   PetscFunctionBegin;
1019:   *nlnk = 0;
1020:   for (PetscInt k = 0, lnkdata = idx_start; k < nidx; ++k) {
1021:     const PetscInt entry = indices[k];

1023:     ++nzbd;
1024:     if (entry == diag) im[idx_start] = nzbd;
1025:     if (!PetscBTLookupSet(bt, entry)) PetscCall(PetscLLInsertLocation_Private(PETSC_TRUE, k, idx_start, entry, nlnk, &lnkdata, lnk));
1026:   }
1027:   PetscFunctionReturn(PETSC_SUCCESS);
1028: }

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

1055: /*
1056:   Free memories used by the list
1057: */
1058: #define PetscLLDestroy(lnk, bt) ((PetscErrorCode)(PetscFree(lnk) || PetscBTDestroy(&(bt))))

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

1075: 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)
1076: {
1077:   PetscFunctionBegin;
1078:   PetscCall(PetscLLInsertLocation_Private(assume_sorted, k, idx_start, entry, nlnk, lnkdata, lnk));
1079:   lnklvl[entry] = newval;
1080:   PetscFunctionReturn(PETSC_SUCCESS);
1081: }

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

1106:     if (!PetscBTLookupSet(bt, entry)) PetscCall(PetscIncompleteLLInsertLocation_Private(PETSC_FALSE, k, idx_start, entry, nlnk, &lnkdata, lnk, lnklvl, 0));
1107:   }
1108:   PetscFunctionReturn(PETSC_SUCCESS);
1109: }

1111: 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)
1112: {
1113:   PetscFunctionBegin;
1114:   *nlnk = 0;
1115:   for (PetscInt k = 0, lnkdata = idx_start; k < nidx; ++k) {
1116:     const PetscInt incrlev = idxlvl[k] + prow_offset + 1;

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

1121:       if (!PetscBTLookupSet(bt, entry)) PetscCall(PetscIncompleteLLInsertLocation_Private(assume_sorted, k, idx_start, entry, nlnk, &lnkdata, lnk, lnklvl, incrlev));
1122:       else if (lnklvl[entry] > incrlev) lnklvl[entry] = incrlev; /* existing entry */
1123:     }
1124:   }
1125:   PetscFunctionReturn(PETSC_SUCCESS);
1126: }

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

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

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

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

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

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

1260: /*
1261:   Free memories used by the list
1262: */
1263: #define PetscIncompleteLLDestroy(lnk, bt) ((PetscErrorCode)(PetscFree(lnk) || PetscBTDestroy(&(bt))))

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

1285: #define VecCheckMatCompatible(M, x, ar1, b, ar2) \
1286:   do { \
1287:     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, \
1288:                (M)->cmap->N); \
1289:     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, \
1290:                (M)->rmap->N); \
1291:   } while (0)

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

1301:       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
1302:       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
1303:       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
1304:       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.
1305:       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
1306:       to each other so memory access is much better than using the big array.

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

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

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

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

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

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

1361:   PetscFunctionBegin;
1362:   for (PetscInt k = 0; k < nidx; k++) {
1363:     const PetscInt entry = indices[k];
1364:     if (!PetscBTLookupSet(bt, entry)) { /* new entry */
1365:       PetscInt next, lnkdata;

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

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

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

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

1411: /*
1412:   Free memories used by the list
1413: */
1414: static inline PetscErrorCode PetscLLCondensedDestroy(PetscInt *lnk, PetscBT bt)
1415: {
1416:   PetscFunctionBegin;
1417:   PetscCall(PetscFree(lnk));
1418:   PetscCall(PetscBTDestroy(&bt));
1419:   PetscFunctionReturn(PETSC_SUCCESS);
1420: }

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

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

1444: static inline PetscErrorCode PetscLLCondensedExpand_Scalable(PetscInt nlnk_max, PetscInt **lnk)
1445: {
1446:   PetscInt lsize = 0;

1448:   PetscFunctionBegin;
1449:   PetscCall(PetscIntMultError(2, nlnk_max + 2, &lsize));
1450:   PetscCall(PetscRealloc(lsize * sizeof(PetscInt), lnk));
1451:   PetscFunctionReturn(PETSC_SUCCESS);
1452: }

1454: static inline PetscErrorCode PetscLLCondensedAddSorted_Scalable(PetscInt nidx, const PetscInt indices[], PetscInt lnk[])
1455: {
1456:   PetscInt location = 2;      /* head */
1457:   PetscInt nlnk     = lnk[0]; /* num of entries on the input lnk */

1459:   for (PetscInt k = 0; k < nidx; k++) {
1460:     const PetscInt entry = indices[k];
1461:     PetscInt       next, lnkdata;

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

1483: static inline PetscErrorCode PetscLLCondensedClean_Scalable(PETSC_UNUSED PetscInt nidx, PetscInt *indices, PetscInt *lnk)
1484: {
1485:   const PetscInt nlnk = lnk[0];
1486:   PetscInt       next = lnk[3]; /* head node */

1488:   for (PetscInt k = 0; k < nlnk; k++) {
1489:     indices[k] = lnk[next];
1490:     next       = lnk[next + 1];
1491:   }
1492:   lnk[0] = 0; /* num of entries on the list */
1493:   lnk[3] = 2; /* head node */
1494:   return PETSC_SUCCESS;
1495: }

1497: static inline PetscErrorCode PetscLLCondensedDestroy_Scalable(PetscInt *lnk)
1498: {
1499:   return PetscFree(lnk);
1500: }

1502: /* -------------------------------------------------------------------------------------------------------*/
1503: /*
1504:       lnk[0]   number of links
1505:       lnk[1]   number of entries
1506:       lnk[3n]  value
1507:       lnk[3n+1] len
1508:       lnk[3n+2] link to next value

1510:       The next three are always the first link

1512:       lnk[3]    PETSC_MIN_INT+1
1513:       lnk[4]    1
1514:       lnk[5]    link to first real entry

1516:       The next three are always the last link

1518:       lnk[6]    PETSC_MAX_INT - 1
1519:       lnk[7]    1
1520:       lnk[8]    next valid link (this is the same as lnk[0] but without the decreases)
1521: */

1523: static inline PetscErrorCode PetscLLCondensedCreate_fast(PetscInt nlnk_max, PetscInt **lnk)
1524: {
1525:   PetscInt *llnk;
1526:   PetscInt  lsize = 0;

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

1543: static inline PetscErrorCode PetscLLCondensedAddSorted_fast(PetscInt nidx, const PetscInt indices[], PetscInt lnk[])
1544: {
1545:   for (PetscInt k = 0, prev = 3 /* first value */; k < nidx; k++) {
1546:     const PetscInt entry = indices[k];
1547:     PetscInt       next  = lnk[prev + 2];

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

1587: static inline PetscErrorCode PetscLLCondensedClean_fast(PETSC_UNUSED PetscInt nidx, PetscInt *indices, PetscInt *lnk)
1588: {
1589:   const PetscInt nlnk = lnk[0];
1590:   PetscInt       next = lnk[5]; /* first node */

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

1607: static inline PetscErrorCode PetscLLCondensedView_fast(const PetscInt *lnk)
1608: {
1609:   const PetscInt nlnk = lnk[0];
1610:   PetscInt       next = lnk[5]; /* first node */

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

1621: static inline PetscErrorCode PetscLLCondensedDestroy_fast(PetscInt *lnk)
1622: {
1623:   return PetscFree(lnk);
1624: }

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

1645: PETSC_EXTERN PetscErrorCode PetscCDGetHeadPos(const PetscCoarsenData *, PetscInt, PetscCDIntNd **);
1646: PETSC_EXTERN PetscErrorCode PetscCDGetNextPos(const PetscCoarsenData *, PetscInt, PetscCDIntNd **);
1647: PETSC_EXTERN PetscErrorCode PetscCDGetASMBlocks(const PetscCoarsenData *, const PetscInt, PetscInt *, IS **);

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

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