Actual source code: matimpl.h
2: #ifndef __MATIMPL_H
5: #include <petscmat.h>
6: #include <petscmatcoarsen.h>
7: #include <petsc/private/petscimpl.h>
9: PETSC_EXTERN PetscBool MatRegisterAllCalled;
10: PETSC_EXTERN PetscBool MatSeqAIJRegisterAllCalled;
11: PETSC_EXTERN PetscBool MatOrderingRegisterAllCalled;
12: PETSC_EXTERN PetscBool MatColoringRegisterAllCalled;
13: PETSC_EXTERN PetscBool MatPartitioningRegisterAllCalled;
14: PETSC_EXTERN PetscBool MatCoarsenRegisterAllCalled;
15: PETSC_EXTERN PetscErrorCode MatRegisterAll(void);
16: PETSC_EXTERN PetscErrorCode MatOrderingRegisterAll(void);
17: PETSC_EXTERN PetscErrorCode MatColoringRegisterAll(void);
18: PETSC_EXTERN PetscErrorCode MatPartitioningRegisterAll(void);
19: PETSC_EXTERN PetscErrorCode MatCoarsenRegisterAll(void);
20: PETSC_EXTERN PetscErrorCode MatSeqAIJRegisterAll(void);
22: /* Gets the root type of the input matrix's type (e.g., MATAIJ for MATSEQAIJ) */
23: PETSC_EXTERN PetscErrorCode MatGetRootType_Private(Mat, MatType *);
25: /* Gets the MPI type corresponding to the input matrix's type (e.g., MATMPIAIJ for MATSEQAIJ) */
26: PETSC_EXTERN PetscErrorCode MatGetMPIMatType_Private(Mat, MatType *);
28: /*
29: This file defines the parts of the matrix data structure that are
30: shared by all matrix types.
31: */
33: /*
34: If you add entries here also add them to the MATOP enum
35: in include/petscmat.h and src/mat/f90-mod/petscmat.h
36: */
37: typedef struct _MatOps *MatOps;
38: struct _MatOps {
39: /* 0*/
40: PetscErrorCode (*setvalues)(Mat, PetscInt, const PetscInt[], PetscInt, const PetscInt[], const PetscScalar[], InsertMode);
41: PetscErrorCode (*getrow)(Mat, PetscInt, PetscInt *, PetscInt *[], PetscScalar *[]);
42: PetscErrorCode (*restorerow)(Mat, PetscInt, PetscInt *, PetscInt *[], PetscScalar *[]);
43: PetscErrorCode (*mult)(Mat, Vec, Vec);
44: PetscErrorCode (*multadd)(Mat, Vec, Vec, Vec);
45: /* 5*/
46: PetscErrorCode (*multtranspose)(Mat, Vec, Vec);
47: PetscErrorCode (*multtransposeadd)(Mat, Vec, Vec, Vec);
48: PetscErrorCode (*solve)(Mat, Vec, Vec);
49: PetscErrorCode (*solveadd)(Mat, Vec, Vec, Vec);
50: PetscErrorCode (*solvetranspose)(Mat, Vec, Vec);
51: /*10*/
52: PetscErrorCode (*solvetransposeadd)(Mat, Vec, Vec, Vec);
53: PetscErrorCode (*lufactor)(Mat, IS, IS, const MatFactorInfo *);
54: PetscErrorCode (*choleskyfactor)(Mat, IS, const MatFactorInfo *);
55: PetscErrorCode (*sor)(Mat, Vec, PetscReal, MatSORType, PetscReal, PetscInt, PetscInt, Vec);
56: PetscErrorCode (*transpose)(Mat, MatReuse, Mat *);
57: /*15*/
58: PetscErrorCode (*getinfo)(Mat, MatInfoType, MatInfo *);
59: PetscErrorCode (*equal)(Mat, Mat, PetscBool *);
60: PetscErrorCode (*getdiagonal)(Mat, Vec);
61: PetscErrorCode (*diagonalscale)(Mat, Vec, Vec);
62: PetscErrorCode (*norm)(Mat, NormType, PetscReal *);
63: /*20*/
64: PetscErrorCode (*assemblybegin)(Mat, MatAssemblyType);
65: PetscErrorCode (*assemblyend)(Mat, MatAssemblyType);
66: PetscErrorCode (*setoption)(Mat, MatOption, PetscBool);
67: PetscErrorCode (*zeroentries)(Mat);
68: /*24*/
69: PetscErrorCode (*zerorows)(Mat, PetscInt, const PetscInt[], PetscScalar, Vec, Vec);
70: PetscErrorCode (*lufactorsymbolic)(Mat, Mat, IS, IS, const MatFactorInfo *);
71: PetscErrorCode (*lufactornumeric)(Mat, Mat, const MatFactorInfo *);
72: PetscErrorCode (*choleskyfactorsymbolic)(Mat, Mat, IS, const MatFactorInfo *);
73: PetscErrorCode (*choleskyfactornumeric)(Mat, Mat, const MatFactorInfo *);
74: /*29*/
75: PetscErrorCode (*setup)(Mat);
76: PetscErrorCode (*ilufactorsymbolic)(Mat, Mat, IS, IS, const MatFactorInfo *);
77: PetscErrorCode (*iccfactorsymbolic)(Mat, Mat, IS, const MatFactorInfo *);
78: PetscErrorCode (*getdiagonalblock)(Mat, Mat *);
79: PetscErrorCode (*setinf)(Mat);
80: /*34*/
81: PetscErrorCode (*duplicate)(Mat, MatDuplicateOption, Mat *);
82: PetscErrorCode (*forwardsolve)(Mat, Vec, Vec);
83: PetscErrorCode (*backwardsolve)(Mat, Vec, Vec);
84: PetscErrorCode (*ilufactor)(Mat, IS, IS, const MatFactorInfo *);
85: PetscErrorCode (*iccfactor)(Mat, IS, const MatFactorInfo *);
86: /*39*/
87: PetscErrorCode (*axpy)(Mat, PetscScalar, Mat, MatStructure);
88: PetscErrorCode (*createsubmatrices)(Mat, PetscInt, const IS[], const IS[], MatReuse, Mat *[]);
89: PetscErrorCode (*increaseoverlap)(Mat, PetscInt, IS[], PetscInt);
90: PetscErrorCode (*getvalues)(Mat, PetscInt, const PetscInt[], PetscInt, const PetscInt[], PetscScalar[]);
91: PetscErrorCode (*copy)(Mat, Mat, MatStructure);
92: /*44*/
93: PetscErrorCode (*getrowmax)(Mat, Vec, PetscInt[]);
94: PetscErrorCode (*scale)(Mat, PetscScalar);
95: PetscErrorCode (*shift)(Mat, PetscScalar);
96: PetscErrorCode (*diagonalset)(Mat, Vec, InsertMode);
97: PetscErrorCode (*zerorowscolumns)(Mat, PetscInt, const PetscInt[], PetscScalar, Vec, Vec);
98: /*49*/
99: PetscErrorCode (*setrandom)(Mat, PetscRandom);
100: PetscErrorCode (*getrowij)(Mat, PetscInt, PetscBool, PetscBool, PetscInt *, const PetscInt *[], const PetscInt *[], PetscBool *);
101: PetscErrorCode (*restorerowij)(Mat, PetscInt, PetscBool, PetscBool, PetscInt *, const PetscInt *[], const PetscInt *[], PetscBool *);
102: PetscErrorCode (*getcolumnij)(Mat, PetscInt, PetscBool, PetscBool, PetscInt *, const PetscInt *[], const PetscInt *[], PetscBool *);
103: PetscErrorCode (*restorecolumnij)(Mat, PetscInt, PetscBool, PetscBool, PetscInt *, const PetscInt *[], const PetscInt *[], PetscBool *);
104: /*54*/
105: PetscErrorCode (*fdcoloringcreate)(Mat, ISColoring, MatFDColoring);
106: PetscErrorCode (*coloringpatch)(Mat, PetscInt, PetscInt, ISColoringValue[], ISColoring *);
107: PetscErrorCode (*setunfactored)(Mat);
108: PetscErrorCode (*permute)(Mat, IS, IS, Mat *);
109: PetscErrorCode (*setvaluesblocked)(Mat, PetscInt, const PetscInt[], PetscInt, const PetscInt[], const PetscScalar[], InsertMode);
110: /*59*/
111: PetscErrorCode (*createsubmatrix)(Mat, IS, IS, MatReuse, Mat *);
112: PetscErrorCode (*destroy)(Mat);
113: PetscErrorCode (*view)(Mat, PetscViewer);
114: PetscErrorCode (*convertfrom)(Mat, MatType, MatReuse, Mat *);
115: PetscErrorCode (*placeholder_63)(void);
116: /*64*/
117: PetscErrorCode (*matmatmultsymbolic)(Mat, Mat, Mat, PetscReal, Mat);
118: PetscErrorCode (*matmatmultnumeric)(Mat, Mat, Mat, Mat);
119: PetscErrorCode (*setlocaltoglobalmapping)(Mat, ISLocalToGlobalMapping, ISLocalToGlobalMapping);
120: PetscErrorCode (*setvalueslocal)(Mat, PetscInt, const PetscInt[], PetscInt, const PetscInt[], const PetscScalar[], InsertMode);
121: PetscErrorCode (*zerorowslocal)(Mat, PetscInt, const PetscInt[], PetscScalar, Vec, Vec);
122: /*69*/
123: PetscErrorCode (*getrowmaxabs)(Mat, Vec, PetscInt[]);
124: PetscErrorCode (*getrowminabs)(Mat, Vec, PetscInt[]);
125: PetscErrorCode (*convert)(Mat, MatType, MatReuse, Mat *);
126: PetscErrorCode (*hasoperation)(Mat, MatOperation, PetscBool *);
127: PetscErrorCode (*placeholder_73)(void);
128: /*74*/
129: PetscErrorCode (*setvaluesadifor)(Mat, PetscInt, void *);
130: PetscErrorCode (*fdcoloringapply)(Mat, MatFDColoring, Vec, void *);
131: PetscErrorCode (*setfromoptions)(Mat, PetscOptionItems *);
132: PetscErrorCode (*placeholder_77)(void);
133: PetscErrorCode (*placeholder_78)(void);
134: /*79*/
135: PetscErrorCode (*findzerodiagonals)(Mat, IS *);
136: PetscErrorCode (*mults)(Mat, Vecs, Vecs);
137: PetscErrorCode (*solves)(Mat, Vecs, Vecs);
138: PetscErrorCode (*getinertia)(Mat, PetscInt *, PetscInt *, PetscInt *);
139: PetscErrorCode (*load)(Mat, PetscViewer);
140: /*84*/
141: PetscErrorCode (*issymmetric)(Mat, PetscReal, PetscBool *);
142: PetscErrorCode (*ishermitian)(Mat, PetscReal, PetscBool *);
143: PetscErrorCode (*isstructurallysymmetric)(Mat, PetscBool *);
144: PetscErrorCode (*setvaluesblockedlocal)(Mat, PetscInt, const PetscInt[], PetscInt, const PetscInt[], const PetscScalar[], InsertMode);
145: PetscErrorCode (*getvecs)(Mat, Vec *, Vec *);
146: /*89*/
147: PetscErrorCode (*placeholder_89)(void);
148: PetscErrorCode (*matmultsymbolic)(Mat, Mat, PetscReal, Mat);
149: PetscErrorCode (*matmultnumeric)(Mat, Mat, Mat);
150: PetscErrorCode (*placeholder_92)(void);
151: PetscErrorCode (*ptapsymbolic)(Mat, Mat, PetscReal, Mat); /* double dispatch wrapper routine */
152: /*94*/
153: PetscErrorCode (*ptapnumeric)(Mat, Mat, Mat); /* double dispatch wrapper routine */
154: PetscErrorCode (*placeholder_95)(void);
155: PetscErrorCode (*mattransposemultsymbolic)(Mat, Mat, PetscReal, Mat);
156: PetscErrorCode (*mattransposemultnumeric)(Mat, Mat, Mat);
157: PetscErrorCode (*bindtocpu)(Mat, PetscBool);
158: /*99*/
159: PetscErrorCode (*productsetfromoptions)(Mat);
160: PetscErrorCode (*productsymbolic)(Mat);
161: PetscErrorCode (*productnumeric)(Mat);
162: PetscErrorCode (*conjugate)(Mat); /* complex conjugate */
163: PetscErrorCode (*viewnative)(Mat, PetscViewer);
164: /*104*/
165: PetscErrorCode (*setvaluesrow)(Mat, PetscInt, const PetscScalar[]);
166: PetscErrorCode (*realpart)(Mat);
167: PetscErrorCode (*imaginarypart)(Mat);
168: PetscErrorCode (*getrowuppertriangular)(Mat);
169: PetscErrorCode (*restorerowuppertriangular)(Mat);
170: /*109*/
171: PetscErrorCode (*matsolve)(Mat, Mat, Mat);
172: PetscErrorCode (*matsolvetranspose)(Mat, Mat, Mat);
173: PetscErrorCode (*getrowmin)(Mat, Vec, PetscInt[]);
174: PetscErrorCode (*getcolumnvector)(Mat, Vec, PetscInt);
175: PetscErrorCode (*missingdiagonal)(Mat, PetscBool *, PetscInt *);
176: /*114*/
177: PetscErrorCode (*getseqnonzerostructure)(Mat, Mat *);
178: PetscErrorCode (*create)(Mat);
179: PetscErrorCode (*getghosts)(Mat, PetscInt *, const PetscInt *[]);
180: PetscErrorCode (*getlocalsubmatrix)(Mat, IS, IS, Mat *);
181: PetscErrorCode (*restorelocalsubmatrix)(Mat, IS, IS, Mat *);
182: /*119*/
183: PetscErrorCode (*multdiagonalblock)(Mat, Vec, Vec);
184: PetscErrorCode (*hermitiantranspose)(Mat, MatReuse, Mat *);
185: PetscErrorCode (*multhermitiantranspose)(Mat, Vec, Vec);
186: PetscErrorCode (*multhermitiantransposeadd)(Mat, Vec, Vec, Vec);
187: PetscErrorCode (*getmultiprocblock)(Mat, MPI_Comm, MatReuse, Mat *);
188: /*124*/
189: PetscErrorCode (*findnonzerorows)(Mat, IS *);
190: PetscErrorCode (*getcolumnreductions)(Mat, PetscInt, PetscReal *);
191: PetscErrorCode (*invertblockdiagonal)(Mat, const PetscScalar **);
192: PetscErrorCode (*invertvariableblockdiagonal)(Mat, PetscInt, const PetscInt *, PetscScalar *);
193: PetscErrorCode (*createsubmatricesmpi)(Mat, PetscInt, const IS[], const IS[], MatReuse, Mat **);
194: /*129*/
195: PetscErrorCode (*setvaluesbatch)(Mat, PetscInt, PetscInt, PetscInt *, const PetscScalar *);
196: PetscErrorCode (*placeholder_130)(void);
197: PetscErrorCode (*transposematmultsymbolic)(Mat, Mat, PetscReal, Mat);
198: PetscErrorCode (*transposematmultnumeric)(Mat, Mat, Mat);
199: PetscErrorCode (*transposecoloringcreate)(Mat, ISColoring, MatTransposeColoring);
200: /*134*/
201: PetscErrorCode (*transcoloringapplysptoden)(MatTransposeColoring, Mat, Mat);
202: PetscErrorCode (*transcoloringapplydentosp)(MatTransposeColoring, Mat, Mat);
203: PetscErrorCode (*placeholder_136)(void);
204: PetscErrorCode (*rartsymbolic)(Mat, Mat, PetscReal, Mat); /* double dispatch wrapper routine */
205: PetscErrorCode (*rartnumeric)(Mat, Mat, Mat); /* double dispatch wrapper routine */
206: /*139*/
207: PetscErrorCode (*setblocksizes)(Mat, PetscInt, PetscInt);
208: PetscErrorCode (*aypx)(Mat, PetscScalar, Mat, MatStructure);
209: PetscErrorCode (*residual)(Mat, Vec, Vec, Vec);
210: PetscErrorCode (*fdcoloringsetup)(Mat, ISColoring, MatFDColoring);
211: PetscErrorCode (*findoffblockdiagonalentries)(Mat, IS *);
212: PetscErrorCode (*creatempimatconcatenateseqmat)(MPI_Comm, Mat, PetscInt, MatReuse, Mat *);
213: /*145*/
214: PetscErrorCode (*destroysubmatrices)(PetscInt, Mat *[]);
215: PetscErrorCode (*mattransposesolve)(Mat, Mat, Mat);
216: PetscErrorCode (*getvalueslocal)(Mat, PetscInt, const PetscInt[], PetscInt, const PetscInt[], PetscScalar[]);
217: PetscErrorCode (*creategraph)(Mat, PetscBool, PetscBool, PetscReal, Mat *);
218: PetscErrorCode (*dummy)(Mat);
219: /*150*/
220: PetscErrorCode (*transposesymbolic)(Mat, Mat *);
221: PetscErrorCode (*eliminatezeros)(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>
229: PETSC_EXTERN PetscErrorCode MatRegisterOp(MPI_Comm, const char[], PetscVoidFunction, const char[], PetscInt, ...);
230: PETSC_EXTERN PetscErrorCode MatQueryOp(MPI_Comm, PetscVoidFunction *, const char[], PetscInt, ...);
232: typedef struct _p_MatRootName *MatRootName;
233: struct _p_MatRootName {
234: char *rname, *sname, *mname;
235: MatRootName next;
236: };
238: PETSC_EXTERN MatRootName MatRootNameList;
240: /*
241: Utility private matrix routines
242: */
243: PETSC_INTERN PetscErrorCode MatFindNonzeroRowsOrCols_Basic(Mat, PetscBool, PetscReal, IS *);
244: PETSC_INTERN PetscErrorCode MatConvert_Basic(Mat, MatType, MatReuse, Mat *);
245: PETSC_INTERN PetscErrorCode MatConvert_Shell(Mat, MatType, MatReuse, Mat *);
246: PETSC_INTERN PetscErrorCode MatConvertFrom_Shell(Mat, MatType, MatReuse, Mat *);
247: PETSC_INTERN PetscErrorCode MatCopy_Basic(Mat, Mat, MatStructure);
248: PETSC_INTERN PetscErrorCode MatDiagonalSet_Default(Mat, Vec, InsertMode);
249: #if defined(PETSC_HAVE_SCALAPACK)
250: PETSC_INTERN PetscErrorCode MatConvert_Dense_ScaLAPACK(Mat, MatType, MatReuse, Mat *);
251: #endif
252: PETSC_INTERN PetscErrorCode MatSetPreallocationCOO_Basic(Mat, PetscCount, const PetscInt[], const PetscInt[]);
253: PETSC_INTERN PetscErrorCode MatSetValuesCOO_Basic(Mat, const PetscScalar[], InsertMode);
255: /* these callbacks rely on the old matrix function pointers for
256: matmat operations. They are unsafe, and should be removed.
257: However, the amount of work needed to clean up all the
258: implementations is not negligible */
259: PETSC_INTERN PetscErrorCode MatProductSymbolic_AB(Mat);
260: PETSC_INTERN PetscErrorCode MatProductNumeric_AB(Mat);
261: PETSC_INTERN PetscErrorCode MatProductSymbolic_AtB(Mat);
262: PETSC_INTERN PetscErrorCode MatProductNumeric_AtB(Mat);
263: PETSC_INTERN PetscErrorCode MatProductSymbolic_ABt(Mat);
264: PETSC_INTERN PetscErrorCode MatProductNumeric_ABt(Mat);
265: PETSC_INTERN PetscErrorCode MatProductNumeric_PtAP(Mat);
266: PETSC_INTERN PetscErrorCode MatProductNumeric_RARt(Mat);
267: PETSC_INTERN PetscErrorCode MatProductSymbolic_ABC(Mat);
268: PETSC_INTERN PetscErrorCode MatProductNumeric_ABC(Mat);
270: PETSC_INTERN PetscErrorCode MatProductCreate_Private(Mat, Mat, Mat, Mat);
271: /* this callback handles all the different triple products and
272: does not rely on the function pointers; used by cuSPARSE/hipSPARSE and KOKKOS-KERNELS */
273: PETSC_INTERN PetscErrorCode MatProductSymbolic_ABC_Basic(Mat);
275: /* CreateGraph is common to AIJ seq and mpi */
276: PETSC_INTERN PetscErrorCode MatCreateGraph_Simple_AIJ(Mat, PetscBool, PetscBool, PetscReal, Mat *);
278: #if defined(PETSC_CLANG_STATIC_ANALYZER)
279: template <typename Tm>
280: void MatCheckPreallocated(Tm, int);
281: template <typename Tm>
282: void MatCheckProduct(Tm, int);
283: #else /* PETSC_CLANG_STATIC_ANALYZER */
284: #define MatCheckPreallocated(A, arg) \
285: do { \
286: if (!(A)->preallocated) PetscCall(MatSetUp(A)); \
287: } while (0)
289: #if defined(PETSC_USE_DEBUG)
290: #define MatCheckProduct(A, arg) \
291: do { \
292: PetscCheck((A)->product, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Argument %d \"%s\" is not a matrix obtained from MatProductCreate()", (arg), #A); \
293: } while (0)
294: #else
295: #define MatCheckProduct(A, arg) \
296: do { \
297: } while (0)
298: #endif
299: #endif /* PETSC_CLANG_STATIC_ANALYZER */
301: /*
302: The stash is used to temporarily store inserted matrix values that
303: belong to another processor. During the assembly phase the stashed
304: values are moved to the correct processor and
305: */
307: typedef struct _MatStashSpace *PetscMatStashSpace;
309: struct _MatStashSpace {
310: PetscMatStashSpace next;
311: PetscScalar *space_head, *val;
312: PetscInt *idx, *idy;
313: PetscInt total_space_size;
314: PetscInt local_used;
315: PetscInt local_remaining;
316: };
318: PETSC_EXTERN PetscErrorCode PetscMatStashSpaceGet(PetscInt, PetscInt, PetscMatStashSpace *);
319: PETSC_EXTERN PetscErrorCode PetscMatStashSpaceContiguous(PetscInt, PetscMatStashSpace *, PetscScalar *, PetscInt *, PetscInt *);
320: PETSC_EXTERN PetscErrorCode PetscMatStashSpaceDestroy(PetscMatStashSpace *);
322: typedef struct {
323: PetscInt count;
324: } MatStashHeader;
326: typedef struct {
327: void *buffer; /* Of type blocktype, dynamically constructed */
328: PetscInt count;
329: char pending;
330: } MatStashFrame;
332: typedef struct _MatStash MatStash;
333: struct _MatStash {
334: PetscInt nmax; /* maximum stash size */
335: PetscInt umax; /* user specified max-size */
336: PetscInt oldnmax; /* the nmax value used previously */
337: PetscInt n; /* stash size */
338: PetscInt bs; /* block size of the stash */
339: PetscInt reallocs; /* preserve the no of mallocs invoked */
340: PetscMatStashSpace space_head, space; /* linked list to hold stashed global row/column numbers and matrix values */
342: PetscErrorCode (*ScatterBegin)(Mat, MatStash *, PetscInt *);
343: PetscErrorCode (*ScatterGetMesg)(MatStash *, PetscMPIInt *, PetscInt **, PetscInt **, PetscScalar **, PetscInt *);
344: PetscErrorCode (*ScatterEnd)(MatStash *);
345: PetscErrorCode (*ScatterDestroy)(MatStash *);
347: /* The following variables are used for communication */
348: MPI_Comm comm;
349: PetscMPIInt size, rank;
350: PetscMPIInt tag1, tag2;
351: MPI_Request *send_waits; /* array of send requests */
352: MPI_Request *recv_waits; /* array of receive requests */
353: MPI_Status *send_status; /* array of send status */
354: PetscInt nsends, nrecvs; /* numbers of sends and receives */
355: PetscScalar *svalues; /* sending data */
356: PetscInt *sindices;
357: PetscScalar **rvalues; /* receiving data (values) */
358: PetscInt **rindices; /* receiving data (indices) */
359: PetscInt nprocessed; /* number of messages already processed */
360: PetscMPIInt *flg_v; /* indicates what messages have arrived so far and from whom */
361: PetscBool reproduce;
362: PetscInt reproduce_count;
364: /* The following variables are used for BTS communication */
365: PetscBool first_assembly_done; /* Is the first time matrix assembly done? */
366: PetscBool use_status; /* Use MPI_Status to determine number of items in each message */
367: PetscMPIInt nsendranks;
368: PetscMPIInt nrecvranks;
369: PetscMPIInt *sendranks;
370: PetscMPIInt *recvranks;
371: MatStashHeader *sendhdr, *recvhdr;
372: MatStashFrame *sendframes; /* pointers to the main messages */
373: MatStashFrame *recvframes;
374: MatStashFrame *recvframe_active;
375: PetscInt recvframe_i; /* index of block within active frame */
376: PetscMPIInt recvframe_count; /* Count actually sent for current frame */
377: PetscInt recvcount; /* Number of receives processed so far */
378: PetscMPIInt *some_indices; /* From last call to MPI_Waitsome */
379: MPI_Status *some_statuses; /* Statuses from last call to MPI_Waitsome */
380: PetscMPIInt some_count; /* Number of requests completed in last call to MPI_Waitsome */
381: PetscMPIInt some_i; /* Index of request currently being processed */
382: MPI_Request *sendreqs;
383: MPI_Request *recvreqs;
384: PetscSegBuffer segsendblocks;
385: PetscSegBuffer segrecvframe;
386: PetscSegBuffer segrecvblocks;
387: MPI_Datatype blocktype;
388: size_t blocktype_size;
389: InsertMode *insertmode; /* Pointer to check mat->insertmode and set upon message arrival in case no local values have been set. */
390: };
392: #if !defined(PETSC_HAVE_MPIUNI)
393: PETSC_INTERN PetscErrorCode MatStashScatterDestroy_BTS(MatStash *);
394: #endif
395: PETSC_INTERN PetscErrorCode MatStashCreate_Private(MPI_Comm, PetscInt, MatStash *);
396: PETSC_INTERN PetscErrorCode MatStashDestroy_Private(MatStash *);
397: PETSC_INTERN PetscErrorCode MatStashScatterEnd_Private(MatStash *);
398: PETSC_INTERN PetscErrorCode MatStashSetInitialSize_Private(MatStash *, PetscInt);
399: PETSC_INTERN PetscErrorCode MatStashGetInfo_Private(MatStash *, PetscInt *, PetscInt *);
400: PETSC_INTERN PetscErrorCode MatStashValuesRow_Private(MatStash *, PetscInt, PetscInt, const PetscInt[], const PetscScalar[], PetscBool);
401: PETSC_INTERN PetscErrorCode MatStashValuesCol_Private(MatStash *, PetscInt, PetscInt, const PetscInt[], const PetscScalar[], PetscInt, PetscBool);
402: PETSC_INTERN PetscErrorCode MatStashValuesRowBlocked_Private(MatStash *, PetscInt, PetscInt, const PetscInt[], const PetscScalar[], PetscInt, PetscInt, PetscInt);
403: PETSC_INTERN PetscErrorCode MatStashValuesColBlocked_Private(MatStash *, PetscInt, PetscInt, const PetscInt[], const PetscScalar[], PetscInt, PetscInt, PetscInt);
404: PETSC_INTERN PetscErrorCode MatStashScatterBegin_Private(Mat, MatStash *, PetscInt *);
405: PETSC_INTERN PetscErrorCode MatStashScatterGetMesg_Private(MatStash *, PetscMPIInt *, PetscInt **, PetscInt **, PetscScalar **, PetscInt *);
406: PETSC_INTERN PetscErrorCode MatGetInfo_External(Mat, MatInfoType, MatInfo *);
408: typedef struct {
409: PetscInt dim;
410: PetscInt dims[4];
411: PetscInt starts[4];
412: PetscBool noc; /* this is a single component problem, hence user will not set MatStencil.c */
413: } MatStencilInfo;
415: /* Info about using compressed row format */
416: typedef struct {
417: PetscBool use; /* indicates compressed rows have been checked and will be used */
418: PetscInt nrows; /* number of non-zero rows */
419: PetscInt *i; /* compressed row pointer */
420: PetscInt *rindex; /* compressed row index */
421: } Mat_CompressedRow;
422: PETSC_EXTERN PetscErrorCode MatCheckCompressedRow(Mat, PetscInt, Mat_CompressedRow *, PetscInt *, PetscInt, PetscReal);
424: typedef struct { /* used by MatCreateRedundantMatrix() for reusing matredundant */
425: PetscInt nzlocal, nsends, nrecvs;
426: PetscMPIInt *send_rank, *recv_rank;
427: PetscInt *sbuf_nz, *rbuf_nz, *sbuf_j, **rbuf_j;
428: PetscScalar *sbuf_a, **rbuf_a;
429: MPI_Comm subcomm; /* when user does not provide a subcomm */
430: IS isrow, iscol;
431: Mat *matseq;
432: } Mat_Redundant;
434: typedef struct { /* used by MatProduct() */
435: MatProductType type;
436: char *alg;
437: Mat A, B, C, Dwork;
438: 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, .. */
439: 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 .. */
440: PetscBool symbolic_used_the_fact_C_is_symmetric; /* MatMatMult(A,B,MAT_REUSE_MATRIX,..&C) is still legitimate), we need to redo symbolic! */
441: PetscReal fill;
442: PetscBool api_user; /* used to distinguish command line options and to indicate the matrix values are ready to be consumed at symbolic phase if needed */
444: /* Some products may display the information on the algorithm used */
445: PetscErrorCode (*view)(Mat, PetscViewer);
447: /* many products have intermediate data structures, each specific to Mat types and product type */
448: PetscBool clear; /* whether or not to clear the data structures after MatProductNumeric has been called */
449: void *data; /* where to stash those structures */
450: PetscErrorCode (*destroy)(void *); /* destroy routine */
451: } Mat_Product;
453: struct _p_Mat {
454: PETSCHEADER(struct _MatOps);
455: PetscLayout rmap, cmap;
456: void *data; /* implementation-specific data */
457: MatFactorType factortype; /* MAT_FACTOR_LU, ILU, CHOLESKY or ICC */
458: PetscBool trivialsymbolic; /* indicates the symbolic factorization doesn't actually do a symbolic factorization, it is delayed to the numeric factorization */
459: PetscBool canuseordering; /* factorization can use ordering provide to routine (most PETSc implementations) */
460: MatOrderingType preferredordering[MAT_FACTOR_NUM_TYPES]; /* what is the preferred (or default) ordering for the matrix solver type */
461: PetscBool assembled; /* is the matrix assembled? */
462: PetscBool was_assembled; /* new values inserted into assembled mat */
463: PetscInt num_ass; /* number of times matrix has been assembled */
464: PetscObjectState nonzerostate; /* each time new nonzeros locations are introduced into the matrix this is updated */
465: PetscObjectState ass_nonzerostate; /* nonzero state at last assembly */
466: MatInfo info; /* matrix information */
467: InsertMode insertmode; /* have values been inserted in matrix or added? */
468: MatStash stash, bstash; /* used for assembling off-proc mat emements */
469: MatNullSpace nullsp; /* null space (operator is singular) */
470: MatNullSpace transnullsp; /* null space of transpose of operator */
471: MatNullSpace nearnullsp; /* near null space to be used by multigrid methods */
472: PetscInt congruentlayouts; /* are the rows and columns layouts congruent? */
473: PetscBool preallocated;
474: MatStencilInfo stencil; /* information for structured grid */
475: PetscBool3 symmetric, hermitian, structurally_symmetric, spd;
476: PetscBool symmetry_eternal, structural_symmetry_eternal, spd_eternal;
477: PetscBool nooffprocentries, nooffproczerorows;
478: PetscBool assembly_subset; /* set by MAT_SUBSET_OFF_PROC_ENTRIES */
479: PetscBool submat_singleis; /* for efficient PCSetUp_ASM() */
480: PetscBool structure_only;
481: PetscBool sortedfull; /* full, sorted rows are inserted */
482: PetscBool force_diagonals; /* set by MAT_FORCE_DIAGONAL_ENTRIES */
483: #if defined(PETSC_HAVE_DEVICE)
484: PetscOffloadMask offloadmask; /* a mask which indicates where the valid matrix data is (GPU, CPU or both) */
485: PetscBool boundtocpu;
486: PetscBool bindingpropagates;
487: #endif
488: char *defaultrandtype;
489: void *spptr; /* pointer for special library like SuperLU */
490: char *solvertype;
491: PetscBool checksymmetryonassembly, checknullspaceonassembly;
492: PetscReal checksymmetrytol;
493: Mat schur; /* Schur complement matrix */
494: MatFactorSchurStatus schur_status; /* status of the Schur complement matrix */
495: Mat_Redundant *redundant; /* used by MatCreateRedundantMatrix() */
496: PetscBool erroriffailure; /* Generate an error if detected (for example a zero pivot) instead of returning */
497: MatFactorError factorerrortype; /* type of error in factorization */
498: PetscReal factorerror_zeropivot_value; /* If numerical zero pivot was detected this is the computed value */
499: PetscInt factorerror_zeropivot_row; /* Row where zero pivot was detected */
500: PetscInt nblocks, *bsizes; /* support for MatSetVariableBlockSizes() */
501: PetscInt p_cstart, p_rank, p_cend, n_rank; /* Information from parallel MatComputeVariableBlockEnvelope() */
502: PetscBool p_parallel;
503: char *defaultvectype;
504: Mat_Product *product;
505: PetscBool form_explicit_transpose; /* hint to generate an explicit mat tranpsose for operations like MatMultTranspose() */
506: PetscBool transupdated; /* whether or not the explicitly generated transpose is up-to-date */
507: char *factorprefix; /* the prefix to use with factored matrix that is created */
508: PetscBool hash_active; /* indicates MatSetValues() is being handled by hashing */
509: };
511: PETSC_INTERN PetscErrorCode MatAXPY_Basic(Mat, PetscScalar, Mat, MatStructure);
512: PETSC_INTERN PetscErrorCode MatAXPY_BasicWithPreallocation(Mat, Mat, PetscScalar, Mat, MatStructure);
513: PETSC_INTERN PetscErrorCode MatAXPY_Basic_Preallocate(Mat, Mat, Mat *);
514: PETSC_INTERN PetscErrorCode MatAXPY_Dense_Nest(Mat, PetscScalar, Mat);
516: PETSC_INTERN PetscErrorCode MatSetUp_Default(Mat);
518: /*
519: Utility for MatZeroRows
520: */
521: PETSC_INTERN PetscErrorCode MatZeroRowsMapLocal_Private(Mat, PetscInt, const PetscInt *, PetscInt *, PetscInt **);
523: /*
524: Utility for MatView/MatLoad
525: */
526: PETSC_INTERN PetscErrorCode MatView_Binary_BlockSizes(Mat, PetscViewer);
527: PETSC_INTERN PetscErrorCode MatLoad_Binary_BlockSizes(Mat, PetscViewer);
529: /*
530: Object for partitioning graphs
531: */
533: typedef struct _MatPartitioningOps *MatPartitioningOps;
534: struct _MatPartitioningOps {
535: PetscErrorCode (*apply)(MatPartitioning, IS *);
536: PetscErrorCode (*applynd)(MatPartitioning, IS *);
537: PetscErrorCode (*setfromoptions)(MatPartitioning, PetscOptionItems *);
538: PetscErrorCode (*destroy)(MatPartitioning);
539: PetscErrorCode (*view)(MatPartitioning, PetscViewer);
540: PetscErrorCode (*improve)(MatPartitioning, IS *);
541: };
543: struct _p_MatPartitioning {
544: PETSCHEADER(struct _MatPartitioningOps);
545: Mat adj;
546: PetscInt *vertex_weights;
547: PetscReal *part_weights;
548: PetscInt n; /* number of partitions */
549: PetscInt ncon; /* number of vertex weights per vertex */
550: void *data;
551: PetscInt setupcalled;
552: PetscBool use_edge_weights; /* A flag indicates whether or not to use edge weights */
553: };
555: /* needed for parallel nested dissection by ParMetis and PTSCOTCH */
556: PETSC_INTERN PetscErrorCode MatPartitioningSizesToSep_Private(PetscInt, PetscInt[], PetscInt[], PetscInt[]);
558: /*
559: Object for coarsen graphs
560: */
561: typedef struct _MatCoarsenOps *MatCoarsenOps;
562: struct _MatCoarsenOps {
563: PetscErrorCode (*apply)(MatCoarsen);
564: PetscErrorCode (*setfromoptions)(MatCoarsen, PetscOptionItems *);
565: PetscErrorCode (*destroy)(MatCoarsen);
566: PetscErrorCode (*view)(MatCoarsen, PetscViewer);
567: };
569: struct _p_MatCoarsen {
570: PETSCHEADER(struct _MatCoarsenOps);
571: Mat graph;
572: void *subctx;
573: /* */
574: PetscBool strict_aggs;
575: IS perm;
576: PetscCoarsenData *agg_lists;
577: };
579: PETSC_EXTERN PetscErrorCode MatCoarsenMISKSetDistance(MatCoarsen, PetscInt);
580: PETSC_EXTERN PetscErrorCode MatCoarsenMISKGetDistance(MatCoarsen, PetscInt *);
582: /*
583: Used in aijdevice.h
584: */
585: typedef struct {
586: PetscInt *i;
587: PetscInt *j;
588: PetscScalar *a;
589: PetscInt n;
590: PetscInt ignorezeroentries;
591: } PetscCSRDataStructure;
593: struct _n_SplitCSRMat {
594: PetscInt cstart, cend, rstart, rend;
595: PetscCSRDataStructure diag, offdiag;
596: PetscInt *colmap;
597: PetscInt M; // number of columns for out of bounds check
598: PetscMPIInt rank;
599: PetscBool allocated_indices;
600: };
602: /*
603: MatFDColoring is used to compute Jacobian matrices efficiently
604: via coloring. The data structure is explained below in an example.
606: Color = 0 1 0 2 | 2 3 0
607: ---------------------------------------------------
608: 00 01 | 05
609: 10 11 | 14 15 Processor 0
610: 22 23 | 25
611: 32 33 |
612: ===================================================
613: | 44 45 46
614: 50 | 55 Processor 1
615: | 64 66
616: ---------------------------------------------------
618: ncolors = 4;
620: ncolumns = {2,1,1,0}
621: columns = {{0,2},{1},{3},{}}
622: nrows = {4,2,3,3}
623: rows = {{0,1,2,3},{0,1},{1,2,3},{0,1,2}}
624: vwscale = {dx(0),dx(1),dx(2),dx(3)} MPI Vec
625: vscale = {dx(0),dx(1),dx(2),dx(3),dx(4),dx(5)} Seq Vec
627: ncolumns = {1,0,1,1}
628: columns = {{6},{},{4},{5}}
629: nrows = {3,0,2,2}
630: rows = {{0,1,2},{},{1,2},{1,2}}
631: vwscale = {dx(4),dx(5),dx(6)} MPI Vec
632: vscale = {dx(0),dx(4),dx(5),dx(6)} Seq Vec
634: See the routine MatFDColoringApply() for how this data is used
635: to compute the Jacobian.
637: */
638: typedef struct {
639: PetscInt row;
640: PetscInt col;
641: PetscScalar *valaddr; /* address of value */
642: } MatEntry;
644: typedef struct {
645: PetscInt row;
646: PetscScalar *valaddr; /* address of value */
647: } MatEntry2;
649: struct _p_MatFDColoring {
650: PETSCHEADER(int);
651: PetscInt M, N, m; /* total rows, columns; local rows */
652: PetscInt rstart; /* first row owned by local processor */
653: PetscInt ncolors; /* number of colors */
654: PetscInt *ncolumns; /* number of local columns for a color */
655: PetscInt **columns; /* lists the local columns of each color (using global column numbering) */
656: IS *isa; /* these are the IS that contain the column values given in columns */
657: PetscInt *nrows; /* number of local rows for each color */
658: MatEntry *matentry; /* holds (row, column, address of value) for Jacobian matrix entry */
659: MatEntry2 *matentry2; /* holds (row, address of value) for Jacobian matrix entry */
660: PetscScalar *dy; /* store a block of F(x+dx)-F(x) when J is in BAIJ format */
661: PetscReal error_rel; /* square root of relative error in computing function */
662: PetscReal umin; /* minimum allowable u'dx value */
663: Vec w1, w2, w3; /* work vectors used in computing Jacobian */
664: PetscBool fset; /* indicates that the initial function value F(X) is set */
665: PetscErrorCode (*f)(void); /* function that defines Jacobian */
666: void *fctx; /* optional user-defined context for use by the function f */
667: Vec vscale; /* holds FD scaling, i.e. 1/dx for each perturbed column */
668: PetscInt currentcolor; /* color for which function evaluation is being done now */
669: const char *htype; /* "wp" or "ds" */
670: ISColoringType ctype; /* IS_COLORING_GLOBAL or IS_COLORING_LOCAL */
671: PetscInt brows, bcols; /* number of block rows or columns for speedup inserting the dense matrix into sparse Jacobian */
672: PetscBool setupcalled; /* true if setup has been called */
673: PetscBool viewed; /* true if the -mat_fd_coloring_view has been triggered already */
674: void (*ftn_func_pointer)(void), *ftn_func_cntx; /* serve the same purpose as *fortran_func_pointers in PETSc objects */
675: PetscObjectId matid; /* matrix this object was created with, must always be the same */
676: };
678: typedef struct _MatColoringOps *MatColoringOps;
679: struct _MatColoringOps {
680: PetscErrorCode (*destroy)(MatColoring);
681: PetscErrorCode (*setfromoptions)(MatColoring, PetscOptionItems *);
682: PetscErrorCode (*view)(MatColoring, PetscViewer);
683: PetscErrorCode (*apply)(MatColoring, ISColoring *);
684: PetscErrorCode (*weights)(MatColoring, PetscReal **, PetscInt **);
685: };
687: struct _p_MatColoring {
688: PETSCHEADER(struct _MatColoringOps);
689: Mat mat;
690: PetscInt dist; /* distance of the coloring */
691: PetscInt maxcolors; /* the maximum number of colors returned, maxcolors=1 for MIS */
692: void *data; /* inner context */
693: PetscBool valid; /* check to see if what is produced is a valid coloring */
694: MatColoringWeightType weight_type; /* type of weight computation to be performed */
695: PetscReal *user_weights; /* custom weights and permutation */
696: PetscInt *user_lperm;
697: PetscBool valid_iscoloring; /* check to see if matcoloring is produced a valid iscoloring */
698: };
700: struct _p_MatTransposeColoring {
701: PETSCHEADER(int);
702: PetscInt M, N, m; /* total rows, columns; local rows */
703: PetscInt rstart; /* first row owned by local processor */
704: PetscInt ncolors; /* number of colors */
705: PetscInt *ncolumns; /* number of local columns for a color */
706: PetscInt *nrows; /* number of local rows for each color */
707: PetscInt currentcolor; /* color for which function evaluation is being done now */
708: ISColoringType ctype; /* IS_COLORING_GLOBAL or IS_COLORING_LOCAL */
710: PetscInt *colorforrow, *colorforcol; /* pointer to rows and columns */
711: PetscInt *rows; /* lists the local rows for each color (using the local row numbering) */
712: PetscInt *den2sp; /* maps (row,color) in the dense matrix to index of sparse matrix array a->a */
713: PetscInt *columns; /* lists the local columns of each color (using global column numbering) */
714: PetscInt brows; /* number of rows for efficient implementation of MatTransColoringApplyDenToSp() */
715: PetscInt *lstart; /* array used for loop over row blocks of Csparse */
716: };
718: /*
719: Null space context for preconditioner/operators
720: */
721: struct _p_MatNullSpace {
722: PETSCHEADER(int);
723: PetscBool has_cnst;
724: PetscInt n;
725: Vec *vecs;
726: PetscScalar *alpha; /* for projections */
727: PetscErrorCode (*remove)(MatNullSpace, Vec, void *); /* for user provided removal function */
728: void *rmctx; /* context for remove() function */
729: };
731: /*
732: Checking zero pivot for LU, ILU preconditioners.
733: */
734: typedef struct {
735: PetscInt nshift, nshift_max;
736: PetscReal shift_amount, shift_lo, shift_hi, shift_top, shift_fraction;
737: PetscBool newshift;
738: PetscReal rs; /* active row sum of abs(offdiagonals) */
739: PetscScalar pv; /* pivot of the active row */
740: } FactorShiftCtx;
742: PETSC_EXTERN PetscErrorCode MatTransposeCheckNonzeroState_Private(Mat, Mat);
744: /*
745: Used by MatTranspose() and potentially other functions to track the matrix used in the generation of another matrix
746: */
747: typedef struct {
748: PetscObjectId id;
749: PetscObjectState state;
750: PetscObjectState nonzerostate;
751: } MatParentState;
753: PETSC_EXTERN PetscErrorCode MatFactorDumpMatrix(Mat);
754: PETSC_INTERN PetscErrorCode MatShift_Basic(Mat, PetscScalar);
755: PETSC_INTERN PetscErrorCode MatSetBlockSizes_Default(Mat, PetscInt, PetscInt);
757: static inline PetscErrorCode MatPivotCheck_nz(PETSC_UNUSED Mat mat, const MatFactorInfo *info, FactorShiftCtx *sctx, PETSC_UNUSED PetscInt row)
758: {
759: PetscReal _rs = sctx->rs;
760: PetscReal _zero = info->zeropivot * _rs;
762: PetscFunctionBegin;
763: if (PetscAbsScalar(sctx->pv) <= _zero && !PetscIsNanScalar(sctx->pv)) {
764: /* force |diag| > zeropivot*rs */
765: if (!sctx->nshift) sctx->shift_amount = info->shiftamount;
766: else sctx->shift_amount *= 2.0;
767: sctx->newshift = PETSC_TRUE;
768: (sctx->nshift)++;
769: } else {
770: sctx->newshift = PETSC_FALSE;
771: }
772: PetscFunctionReturn(PETSC_SUCCESS);
773: }
775: static inline PetscErrorCode MatPivotCheck_pd(PETSC_UNUSED Mat mat, const MatFactorInfo *info, FactorShiftCtx *sctx, PETSC_UNUSED PetscInt row)
776: {
777: PetscReal _rs = sctx->rs;
778: PetscReal _zero = info->zeropivot * _rs;
780: PetscFunctionBegin;
781: if (PetscRealPart(sctx->pv) <= _zero && !PetscIsNanScalar(sctx->pv)) {
782: /* force matfactor to be diagonally dominant */
783: if (sctx->nshift == sctx->nshift_max) {
784: sctx->shift_fraction = sctx->shift_hi;
785: } else {
786: sctx->shift_lo = sctx->shift_fraction;
787: sctx->shift_fraction = (sctx->shift_hi + sctx->shift_lo) / 2.;
788: }
789: sctx->shift_amount = sctx->shift_fraction * sctx->shift_top;
790: sctx->nshift++;
791: sctx->newshift = PETSC_TRUE;
792: } else {
793: sctx->newshift = PETSC_FALSE;
794: }
795: PetscFunctionReturn(PETSC_SUCCESS);
796: }
798: static inline PetscErrorCode MatPivotCheck_inblocks(PETSC_UNUSED Mat mat, const MatFactorInfo *info, FactorShiftCtx *sctx, PETSC_UNUSED PetscInt row)
799: {
800: PetscReal _zero = info->zeropivot;
802: PetscFunctionBegin;
803: if (PetscAbsScalar(sctx->pv) <= _zero && !PetscIsNanScalar(sctx->pv)) {
804: sctx->pv += info->shiftamount;
805: sctx->shift_amount = 0.0;
806: sctx->nshift++;
807: }
808: sctx->newshift = PETSC_FALSE;
809: PetscFunctionReturn(PETSC_SUCCESS);
810: }
812: static inline PetscErrorCode MatPivotCheck_none(Mat fact, Mat mat, const MatFactorInfo *info, FactorShiftCtx *sctx, PetscInt row)
813: {
814: PetscReal _zero = info->zeropivot;
816: PetscFunctionBegin;
817: sctx->newshift = PETSC_FALSE;
818: if (PetscAbsScalar(sctx->pv) <= _zero && !PetscIsNanScalar(sctx->pv)) {
819: 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);
820: PetscCall(PetscInfo(mat, "Detected zero pivot in factorization in row %" PetscInt_FMT " value %g tolerance %g\n", row, (double)PetscAbsScalar(sctx->pv), (double)_zero));
821: fact->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
822: fact->factorerror_zeropivot_value = PetscAbsScalar(sctx->pv);
823: fact->factorerror_zeropivot_row = row;
824: }
825: PetscFunctionReturn(PETSC_SUCCESS);
826: }
828: static inline PetscErrorCode MatPivotCheck(Mat fact, Mat mat, const MatFactorInfo *info, FactorShiftCtx *sctx, PetscInt row)
829: {
830: PetscFunctionBegin;
831: if (info->shifttype == (PetscReal)MAT_SHIFT_NONZERO) PetscCall(MatPivotCheck_nz(mat, info, sctx, row));
832: else if (info->shifttype == (PetscReal)MAT_SHIFT_POSITIVE_DEFINITE) PetscCall(MatPivotCheck_pd(mat, info, sctx, row));
833: else if (info->shifttype == (PetscReal)MAT_SHIFT_INBLOCKS) PetscCall(MatPivotCheck_inblocks(mat, info, sctx, row));
834: else PetscCall(MatPivotCheck_none(fact, mat, info, sctx, row));
835: PetscFunctionReturn(PETSC_SUCCESS);
836: }
838: #include <petscbt.h>
839: /*
840: Create and initialize a linked list
841: Input Parameters:
842: idx_start - starting index of the list
843: lnk_max - max value of lnk indicating the end of the list
844: nlnk - max length of the list
845: Output Parameters:
846: lnk - list initialized
847: bt - PetscBT (bitarray) with all bits set to false
848: lnk_empty - flg indicating the list is empty
849: */
850: #define PetscLLCreate(idx_start, lnk_max, nlnk, lnk, bt) ((PetscErrorCode)(PetscMalloc1(nlnk, &lnk) || PetscBTCreate(nlnk, &(bt)) || (lnk[idx_start] = lnk_max, PETSC_SUCCESS)))
852: #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)))
854: 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)
855: {
856: PetscInt location;
858: PetscFunctionBegin;
859: /* start from the beginning if entry < previous entry */
860: if (!assume_sorted && k && entry < *lnkdata) *lnkdata = idx_start;
861: /* search for insertion location */
862: do {
863: location = *lnkdata;
864: *lnkdata = lnk[location];
865: } while (entry > *lnkdata);
866: /* insertion location is found, add entry into lnk */
867: lnk[location] = entry;
868: lnk[entry] = *lnkdata;
869: ++(*nlnk);
870: *lnkdata = entry; /* next search starts from here if next_entry > entry */
871: PetscFunctionReturn(PETSC_SUCCESS);
872: }
874: 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)
875: {
876: PetscFunctionBegin;
877: *nlnk = 0;
878: for (PetscInt k = 0, lnkdata = idx_start; k < nidx; ++k) {
879: const PetscInt entry = indices[k];
881: if (!PetscBTLookupSet(bt, entry)) PetscCall(PetscLLInsertLocation_Private(assume_sorted, k, idx_start, entry, nlnk, &lnkdata, lnk));
882: }
883: PetscFunctionReturn(PETSC_SUCCESS);
884: }
886: /*
887: Add an index set into a sorted linked list
888: Input Parameters:
889: nidx - number of input indices
890: indices - integer array
891: idx_start - starting index of the list
892: lnk - linked list(an integer array) that is created
893: bt - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
894: output Parameters:
895: nlnk - number of newly added indices
896: lnk - the sorted(increasing order) linked list containing new and non-redundate entries from indices
897: bt - updated PetscBT (bitarray)
898: */
899: static inline PetscErrorCode PetscLLAdd(PetscInt nidx, const PetscInt *PETSC_RESTRICT indices, PetscInt idx_start, PetscInt *PETSC_RESTRICT nlnk, PetscInt *PETSC_RESTRICT lnk, PetscBT bt)
900: {
901: PetscFunctionBegin;
902: PetscCall(PetscLLAdd_Private(nidx, indices, idx_start, nlnk, lnk, bt, PETSC_FALSE));
903: PetscFunctionReturn(PETSC_SUCCESS);
904: }
906: /*
907: Add a SORTED ascending index set into a sorted linked list - same as PetscLLAdd() bus skip 'if (_k && _entry < _lnkdata) _lnkdata = idx_start;'
908: Input Parameters:
909: nidx - number of input indices
910: indices - sorted integer array
911: idx_start - starting index of the list
912: lnk - linked list(an integer array) that is created
913: bt - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
914: output Parameters:
915: nlnk - number of newly added indices
916: lnk - the sorted(increasing order) linked list containing new and non-redundate entries from indices
917: bt - updated PetscBT (bitarray)
918: */
919: static inline PetscErrorCode PetscLLAddSorted(PetscInt nidx, const PetscInt *PETSC_RESTRICT indices, PetscInt idx_start, PetscInt *PETSC_RESTRICT nlnk, PetscInt *PETSC_RESTRICT lnk, PetscBT bt)
920: {
921: PetscFunctionBegin;
922: PetscCall(PetscLLAdd_Private(nidx, indices, idx_start, nlnk, lnk, bt, PETSC_TRUE));
923: PetscFunctionReturn(PETSC_SUCCESS);
924: }
926: /*
927: Add a permuted index set into a sorted linked list
928: Input Parameters:
929: nidx - number of input indices
930: indices - integer array
931: perm - permutation of indices
932: idx_start - starting index of the list
933: lnk - linked list(an integer array) that is created
934: bt - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
935: output Parameters:
936: nlnk - number of newly added indices
937: lnk - the sorted(increasing order) linked list containing new and non-redundate entries from indices
938: bt - updated PetscBT (bitarray)
939: */
940: 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)
941: {
942: PetscFunctionBegin;
943: *nlnk = 0;
944: for (PetscInt k = 0, lnkdata = idx_start; k < nidx; ++k) {
945: const PetscInt entry = perm[indices[k]];
947: if (!PetscBTLookupSet(bt, entry)) PetscCall(PetscLLInsertLocation_Private(PETSC_FALSE, k, idx_start, entry, nlnk, &lnkdata, lnk));
948: }
949: PetscFunctionReturn(PETSC_SUCCESS);
950: }
952: #if 0
953: /* this appears to be unused? */
954: static inline PetscErrorCode PetscLLAddSorted_new(PetscInt nidx, PetscInt *indices, PetscInt idx_start, PetscBool *lnk_empty, PetscInt *nlnk, PetscInt *lnk, PetscBT bt)
955: {
956: PetscInt lnkdata = idx_start;
958: PetscFunctionBegin;
959: if (*lnk_empty) {
960: for (PetscInt k = 0; k < nidx; ++k) {
961: const PetscInt entry = indices[k], location = lnkdata;
963: PetscCall(PetscBTSet(bt,entry)); /* mark the new entry */
964: lnkdata = lnk[location];
965: /* insertion location is found, add entry into lnk */
966: lnk[location] = entry;
967: lnk[entry] = lnkdata;
968: lnkdata = entry; /* next search starts from here */
969: }
970: /* lnk[indices[nidx-1]] = lnk[idx_start];
971: lnk[idx_start] = indices[0];
972: PetscCall(PetscBTSet(bt,indices[0]));
973: for (_k=1; _k<nidx; _k++) {
974: PetscCall(PetscBTSet(bt,indices[_k]));
975: lnk[indices[_k-1]] = indices[_k];
976: }
977: */
978: *nlnk = nidx;
979: *lnk_empty = PETSC_FALSE;
980: } else {
981: *nlnk = 0;
982: for (PetscInt k = 0; k < nidx; ++k) {
983: const PetscInt entry = indices[k];
985: if (!PetscBTLookupSet(bt,entry)) PetscCall(PetscLLInsertLocation_Private(PETSC_TRUE,k,idx_start,entry,nlnk,&lnkdata,lnk));
986: }
987: }
988: PetscFunctionReturn(PETSC_SUCCESS);
989: }
990: #endif
992: /*
993: Add a SORTED index set into a sorted linked list used for LUFactorSymbolic()
994: Same as PetscLLAddSorted() with an additional operation:
995: count the number of input indices that are no larger than 'diag'
996: Input Parameters:
997: indices - sorted integer array
998: idx_start - starting index of the list, index of pivot row
999: lnk - linked list(an integer array) that is created
1000: bt - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
1001: diag - index of the active row in LUFactorSymbolic
1002: nzbd - number of input indices with indices <= idx_start
1003: im - im[idx_start] is initialized as num of nonzero entries in row=idx_start
1004: output Parameters:
1005: nlnk - number of newly added indices
1006: lnk - the sorted(increasing order) linked list containing new and non-redundate entries from indices
1007: bt - updated PetscBT (bitarray)
1008: im - im[idx_start]: unchanged if diag is not an entry
1009: : num of entries with indices <= diag if diag is an entry
1010: */
1011: 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)
1012: {
1013: const PetscInt nidx = im[idx_start] - nzbd; /* num of entries with idx_start < index <= diag */
1015: PetscFunctionBegin;
1016: *nlnk = 0;
1017: for (PetscInt k = 0, lnkdata = idx_start; k < nidx; ++k) {
1018: const PetscInt entry = indices[k];
1020: ++nzbd;
1021: if (entry == diag) im[idx_start] = nzbd;
1022: if (!PetscBTLookupSet(bt, entry)) PetscCall(PetscLLInsertLocation_Private(PETSC_TRUE, k, idx_start, entry, nlnk, &lnkdata, lnk));
1023: }
1024: PetscFunctionReturn(PETSC_SUCCESS);
1025: }
1027: /*
1028: Copy data on the list into an array, then initialize the list
1029: Input Parameters:
1030: idx_start - starting index of the list
1031: lnk_max - max value of lnk indicating the end of the list
1032: nlnk - number of data on the list to be copied
1033: lnk - linked list
1034: bt - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
1035: output Parameters:
1036: indices - array that contains the copied data
1037: lnk - linked list that is cleaned and initialize
1038: bt - PetscBT (bitarray) with all bits set to false
1039: */
1040: static inline PetscErrorCode PetscLLClean(PetscInt idx_start, PetscInt lnk_max, PetscInt nlnk, PetscInt *PETSC_RESTRICT lnk, PetscInt *PETSC_RESTRICT indices, PetscBT bt)
1041: {
1042: PetscFunctionBegin;
1043: for (PetscInt j = 0, idx = idx_start; j < nlnk; ++j) {
1044: idx = lnk[idx];
1045: indices[j] = idx;
1046: PetscCall(PetscBTClear(bt, idx));
1047: }
1048: lnk[idx_start] = lnk_max;
1049: PetscFunctionReturn(PETSC_SUCCESS);
1050: }
1052: /*
1053: Free memories used by the list
1054: */
1055: #define PetscLLDestroy(lnk, bt) ((PetscErrorCode)(PetscFree(lnk) || PetscBTDestroy(&(bt))))
1057: /* Routines below are used for incomplete matrix factorization */
1058: /*
1059: Create and initialize a linked list and its levels
1060: Input Parameters:
1061: idx_start - starting index of the list
1062: lnk_max - max value of lnk indicating the end of the list
1063: nlnk - max length of the list
1064: Output Parameters:
1065: lnk - list initialized
1066: lnk_lvl - array of size nlnk for storing levels of lnk
1067: bt - PetscBT (bitarray) with all bits set to false
1068: */
1069: #define PetscIncompleteLLCreate(idx_start, lnk_max, nlnk, lnk, lnk_lvl, bt) \
1070: ((PetscErrorCode)(PetscIntMultError(2, nlnk, NULL) || PetscMalloc1(2 * nlnk, &lnk) || PetscBTCreate(nlnk, &(bt)) || (lnk[idx_start] = lnk_max, lnk_lvl = lnk + nlnk, PETSC_SUCCESS)))
1072: 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)
1073: {
1074: PetscFunctionBegin;
1075: PetscCall(PetscLLInsertLocation_Private(assume_sorted, k, idx_start, entry, nlnk, lnkdata, lnk));
1076: lnklvl[entry] = newval;
1077: PetscFunctionReturn(PETSC_SUCCESS);
1078: }
1080: /*
1081: Initialize a sorted linked list used for ILU and ICC
1082: Input Parameters:
1083: nidx - number of input idx
1084: idx - integer array used for storing column indices
1085: idx_start - starting index of the list
1086: perm - indices of an IS
1087: lnk - linked list(an integer array) that is created
1088: lnklvl - levels of lnk
1089: bt - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
1090: output Parameters:
1091: nlnk - number of newly added idx
1092: lnk - the sorted(increasing order) linked list containing new and non-redundate entries from idx
1093: lnklvl - levels of lnk
1094: bt - updated PetscBT (bitarray)
1095: */
1096: 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)
1097: {
1098: PetscFunctionBegin;
1099: *nlnk = 0;
1100: for (PetscInt k = 0, lnkdata = idx_start; k < nidx; ++k) {
1101: const PetscInt entry = perm[idx[k]];
1103: if (!PetscBTLookupSet(bt, entry)) PetscCall(PetscIncompleteLLInsertLocation_Private(PETSC_FALSE, k, idx_start, entry, nlnk, &lnkdata, lnk, lnklvl, 0));
1104: }
1105: PetscFunctionReturn(PETSC_SUCCESS);
1106: }
1108: 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)
1109: {
1110: PetscFunctionBegin;
1111: *nlnk = 0;
1112: for (PetscInt k = 0, lnkdata = idx_start; k < nidx; ++k) {
1113: const PetscInt incrlev = idxlvl[k] + prow_offset + 1;
1115: if (incrlev <= level) {
1116: const PetscInt entry = idx[k];
1118: if (!PetscBTLookupSet(bt, entry)) PetscCall(PetscIncompleteLLInsertLocation_Private(assume_sorted, k, idx_start, entry, nlnk, &lnkdata, lnk, lnklvl, incrlev));
1119: else if (lnklvl[entry] > incrlev) lnklvl[entry] = incrlev; /* existing entry */
1120: }
1121: }
1122: PetscFunctionReturn(PETSC_SUCCESS);
1123: }
1125: /*
1126: Add a SORTED index set into a sorted linked list for ICC
1127: Input Parameters:
1128: nidx - number of input indices
1129: idx - sorted integer array used for storing column indices
1130: level - level of fill, e.g., ICC(level)
1131: idxlvl - level of idx
1132: idx_start - starting index of the list
1133: lnk - linked list(an integer array) that is created
1134: lnklvl - levels of lnk
1135: bt - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
1136: idxlvl_prow - idxlvl[prow], where prow is the row number of the idx
1137: output Parameters:
1138: nlnk - number of newly added indices
1139: lnk - the sorted(increasing order) linked list containing new and non-redundate entries from idx
1140: lnklvl - levels of lnk
1141: bt - updated PetscBT (bitarray)
1142: Note: the level of U(i,j) is set as lvl(i,j) = min{ lvl(i,j), lvl(prow,i)+lvl(prow,j)+1)
1143: where idx = non-zero columns of U(prow,prow+1:n-1), prow<i
1144: */
1145: 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)
1146: {
1147: PetscFunctionBegin;
1148: PetscCall(PetscIncompleteLLAdd_Private(nidx, idx, level, idxlvl, idx_start, nlnk, lnk, lnklvl, bt, idxlvl_prow, PETSC_TRUE));
1149: PetscFunctionReturn(PETSC_SUCCESS);
1150: }
1152: /*
1153: Add a SORTED index set into a sorted linked list for ILU
1154: Input Parameters:
1155: nidx - number of input indices
1156: idx - sorted integer array used for storing column indices
1157: level - level of fill, e.g., ICC(level)
1158: idxlvl - level of idx
1159: idx_start - starting index of the list
1160: lnk - linked list(an integer array) that is created
1161: lnklvl - levels of lnk
1162: bt - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
1163: prow - the row number of idx
1164: output Parameters:
1165: nlnk - number of newly added idx
1166: lnk - the sorted(increasing order) linked list containing new and non-redundate entries from idx
1167: lnklvl - levels of lnk
1168: bt - updated PetscBT (bitarray)
1170: Note: the level of factor(i,j) is set as lvl(i,j) = min{ lvl(i,j), lvl(i,prow)+lvl(prow,j)+1)
1171: where idx = non-zero columns of U(prow,prow+1:n-1), prow<i
1172: */
1173: 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)
1174: {
1175: PetscFunctionBegin;
1176: PetscCall(PetscIncompleteLLAdd_Private(nidx, idx, level, idxlvl, idx_start, nlnk, lnk, lnklvl, bt, lnklvl[prow], PETSC_TRUE));
1177: PetscFunctionReturn(PETSC_SUCCESS);
1178: }
1180: /*
1181: Add a index set into a sorted linked list
1182: Input Parameters:
1183: nidx - number of input idx
1184: idx - integer array used for storing column indices
1185: level - level of fill, e.g., ICC(level)
1186: idxlvl - level of idx
1187: idx_start - starting index of the list
1188: lnk - linked list(an integer array) that is created
1189: lnklvl - levels of lnk
1190: bt - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
1191: output Parameters:
1192: nlnk - number of newly added idx
1193: lnk - the sorted(increasing order) linked list containing new and non-redundate entries from idx
1194: lnklvl - levels of lnk
1195: bt - updated PetscBT (bitarray)
1196: */
1197: 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)
1198: {
1199: PetscFunctionBegin;
1200: PetscCall(PetscIncompleteLLAdd_Private(nidx, idx, level, idxlvl, idx_start, nlnk, lnk, lnklvl, bt, 0, PETSC_FALSE));
1201: PetscFunctionReturn(PETSC_SUCCESS);
1202: }
1204: /*
1205: Add a SORTED index set into a sorted linked list
1206: Input Parameters:
1207: nidx - number of input indices
1208: idx - sorted integer array used for storing column indices
1209: level - level of fill, e.g., ICC(level)
1210: idxlvl - level of idx
1211: idx_start - starting index of the list
1212: lnk - linked list(an integer array) that is created
1213: lnklvl - levels of lnk
1214: bt - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
1215: output Parameters:
1216: nlnk - number of newly added idx
1217: lnk - the sorted(increasing order) linked list containing new and non-redundate entries from idx
1218: lnklvl - levels of lnk
1219: bt - updated PetscBT (bitarray)
1220: */
1221: 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)
1222: {
1223: PetscFunctionBegin;
1224: PetscCall(PetscIncompleteLLAdd_Private(nidx, idx, level, idxlvl, idx_start, nlnk, lnk, lnklvl, bt, 0, PETSC_TRUE));
1225: PetscFunctionReturn(PETSC_SUCCESS);
1226: }
1228: /*
1229: Copy data on the list into an array, then initialize the list
1230: Input Parameters:
1231: idx_start - starting index of the list
1232: lnk_max - max value of lnk indicating the end of the list
1233: nlnk - number of data on the list to be copied
1234: lnk - linked list
1235: lnklvl - level of lnk
1236: bt - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
1237: output Parameters:
1238: indices - array that contains the copied data
1239: lnk - linked list that is cleaned and initialize
1240: lnklvl - level of lnk that is reinitialized
1241: bt - PetscBT (bitarray) with all bits set to false
1242: */
1243: 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)
1244: {
1245: PetscFunctionBegin;
1246: for (PetscInt j = 0, idx = idx_start; j < nlnk; ++j) {
1247: idx = lnk[idx];
1248: indices[j] = idx;
1249: indiceslvl[j] = lnklvl[idx];
1250: lnklvl[idx] = -1;
1251: PetscCall(PetscBTClear(bt, idx));
1252: }
1253: lnk[idx_start] = lnk_max;
1254: PetscFunctionReturn(PETSC_SUCCESS);
1255: }
1257: /*
1258: Free memories used by the list
1259: */
1260: #define PetscIncompleteLLDestroy(lnk, bt) ((PetscErrorCode)(PetscFree(lnk) || PetscBTDestroy(&(bt))))
1262: #if !defined(PETSC_CLANG_STATIC_ANALYZER)
1263: #define MatCheckSameLocalSize(A, ar1, B, ar2) \
1264: do { \
1265: PetscCheckSameComm(A, ar1, B, ar2); \
1266: 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, \
1267: (A)->rmap->n, (A)->cmap->n, ar2, (B)->rmap->n, (B)->cmap->n); \
1268: } while (0)
1269: #define MatCheckSameSize(A, ar1, B, ar2) \
1270: do { \
1271: 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, \
1272: (A)->rmap->N, (A)->cmap->N, ar2, (B)->rmap->N, (B)->cmap->N); \
1273: MatCheckSameLocalSize(A, ar1, B, ar2); \
1274: } while (0)
1275: #else
1276: template <typename Tm>
1277: void MatCheckSameLocalSize(Tm, int, Tm, int);
1278: template <typename Tm>
1279: void MatCheckSameSize(Tm, int, Tm, int);
1280: #endif
1282: #define VecCheckMatCompatible(M, x, ar1, b, ar2) \
1283: do { \
1284: 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, \
1285: (M)->cmap->N); \
1286: 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, \
1287: (M)->rmap->N); \
1288: } while (0)
1290: /* -------------------------------------------------------------------------------------------------------*/
1291: /*
1292: Create and initialize a condensed linked list -
1293: same as PetscLLCreate(), but uses a scalable array 'lnk' with size of max number of entries, not O(N).
1294: Barry suggested this approach (Dec. 6, 2011):
1295: I've thought of an alternative way of representing a linked list that is efficient but doesn't have the O(N) scaling issue
1296: (it may be faster than the O(N) even sequentially due to less crazy memory access).
1298: 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
1299: 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
1300: 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
1301: 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.
1302: 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
1303: to each other so memory access is much better than using the big array.
1305: Example:
1306: nlnk_max=5, lnk_max=36:
1307: Initial list: [0, 0 | 36, 2 | 0, 0 | 0, 0 | 0, 0 | 0, 0 | 0, 0]
1308: here, head_node has index 2 with value lnk[2]=lnk_max=36,
1309: 0-th entry is used to store the number of entries in the list,
1310: The initial lnk represents head -> tail(marked by 36) with number of entries = lnk[0]=0.
1312: Now adding a sorted set {2,4}, the list becomes
1313: [2, 0 | 36, 4 |2, 6 | 4, 2 | 0, 0 | 0, 0 | 0, 0 ]
1314: represents head -> 2 -> 4 -> tail with number of entries = lnk[0]=2.
1316: Then adding a sorted set {0,3,35}, the list
1317: [5, 0 | 36, 8 | 2, 10 | 4, 12 | 0, 4 | 3, 6 | 35, 2 ]
1318: represents head -> 0 -> 2 -> 3 -> 4 -> 35 -> tail with number of entries = lnk[0]=5.
1320: Input Parameters:
1321: nlnk_max - max length of the list
1322: lnk_max - max value of the entries
1323: Output Parameters:
1324: lnk - list created and initialized
1325: bt - PetscBT (bitarray) with all bits set to false. Note: bt has size lnk_max, not nln_max!
1326: */
1327: static inline PetscErrorCode PetscLLCondensedCreate(PetscInt nlnk_max, PetscInt lnk_max, PetscInt **lnk, PetscBT *bt)
1328: {
1329: PetscInt *llnk, lsize = 0;
1331: PetscFunctionBegin;
1332: PetscCall(PetscIntMultError(2, nlnk_max + 2, &lsize));
1333: PetscCall(PetscMalloc1(lsize, lnk));
1334: PetscCall(PetscBTCreate(lnk_max, bt));
1335: llnk = *lnk;
1336: llnk[0] = 0; /* number of entries on the list */
1337: llnk[2] = lnk_max; /* value in the head node */
1338: llnk[3] = 2; /* next for the head node */
1339: PetscFunctionReturn(PETSC_SUCCESS);
1340: }
1342: /*
1343: Add a SORTED ascending index set into a sorted linked list. See PetscLLCondensedCreate() for detailed description.
1344: Input Parameters:
1345: nidx - number of input indices
1346: indices - sorted integer array
1347: lnk - condensed linked list(an integer array) that is created
1348: bt - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
1349: output Parameters:
1350: lnk - the sorted(increasing order) linked list containing previous and newly added non-redundate indices
1351: bt - updated PetscBT (bitarray)
1352: */
1353: static inline PetscErrorCode PetscLLCondensedAddSorted(PetscInt nidx, const PetscInt indices[], PetscInt lnk[], PetscBT bt)
1354: {
1355: PetscInt location = 2; /* head */
1356: PetscInt nlnk = lnk[0]; /* num of entries on the input lnk */
1358: PetscFunctionBegin;
1359: for (PetscInt k = 0; k < nidx; k++) {
1360: const PetscInt entry = indices[k];
1361: if (!PetscBTLookupSet(bt, entry)) { /* new entry */
1362: PetscInt next, lnkdata;
1364: /* search for insertion location */
1365: do {
1366: next = location + 1; /* link from previous node to next node */
1367: location = lnk[next]; /* idx of next node */
1368: lnkdata = lnk[location]; /* value of next node */
1369: } while (entry > lnkdata);
1370: /* insertion location is found, add entry into lnk */
1371: const PetscInt newnode = 2 * (nlnk + 2); /* index for this new node */
1372: lnk[next] = newnode; /* connect previous node to the new node */
1373: lnk[newnode] = entry; /* set value of the new node */
1374: lnk[newnode + 1] = location; /* connect new node to next node */
1375: location = newnode; /* next search starts from the new node */
1376: nlnk++;
1377: }
1378: }
1379: lnk[0] = nlnk; /* number of entries in the list */
1380: PetscFunctionReturn(PETSC_SUCCESS);
1381: }
1383: static inline PetscErrorCode PetscLLCondensedClean(PetscInt lnk_max, PETSC_UNUSED PetscInt nidx, PetscInt *indices, PetscInt lnk[], PetscBT bt)
1384: {
1385: const PetscInt nlnk = lnk[0]; /* num of entries on the list */
1386: PetscInt next = lnk[3]; /* head node */
1388: PetscFunctionBegin;
1389: for (PetscInt k = 0; k < nlnk; k++) {
1390: indices[k] = lnk[next];
1391: next = lnk[next + 1];
1392: PetscCall(PetscBTClear(bt, indices[k]));
1393: }
1394: lnk[0] = 0; /* num of entries on the list */
1395: lnk[2] = lnk_max; /* initialize head node */
1396: lnk[3] = 2; /* head node */
1397: PetscFunctionReturn(PETSC_SUCCESS);
1398: }
1400: static inline PetscErrorCode PetscLLCondensedView(PetscInt *lnk)
1401: {
1402: PetscFunctionBegin;
1403: PetscCall(PetscPrintf(PETSC_COMM_SELF, "LLCondensed of size %" PetscInt_FMT ", (val, next)\n", lnk[0]));
1404: 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]));
1405: PetscFunctionReturn(PETSC_SUCCESS);
1406: }
1408: /*
1409: Free memories used by the list
1410: */
1411: static inline PetscErrorCode PetscLLCondensedDestroy(PetscInt *lnk, PetscBT bt)
1412: {
1413: PetscFunctionBegin;
1414: PetscCall(PetscFree(lnk));
1415: PetscCall(PetscBTDestroy(&bt));
1416: PetscFunctionReturn(PETSC_SUCCESS);
1417: }
1419: /* -------------------------------------------------------------------------------------------------------*/
1420: /*
1421: Same as PetscLLCondensedCreate(), but does not use non-scalable O(lnk_max) bitarray
1422: Input Parameters:
1423: nlnk_max - max length of the list
1424: Output Parameters:
1425: lnk - list created and initialized
1426: */
1427: static inline PetscErrorCode PetscLLCondensedCreate_Scalable(PetscInt nlnk_max, PetscInt **lnk)
1428: {
1429: PetscInt *llnk, lsize = 0;
1431: PetscFunctionBegin;
1432: PetscCall(PetscIntMultError(2, nlnk_max + 2, &lsize));
1433: PetscCall(PetscMalloc1(lsize, lnk));
1434: llnk = *lnk;
1435: llnk[0] = 0; /* number of entries on the list */
1436: llnk[2] = PETSC_MAX_INT; /* value in the head node */
1437: llnk[3] = 2; /* next for the head node */
1438: PetscFunctionReturn(PETSC_SUCCESS);
1439: }
1441: static inline PetscErrorCode PetscLLCondensedExpand_Scalable(PetscInt nlnk_max, PetscInt **lnk)
1442: {
1443: PetscInt lsize = 0;
1445: PetscFunctionBegin;
1446: PetscCall(PetscIntMultError(2, nlnk_max + 2, &lsize));
1447: PetscCall(PetscRealloc(lsize * sizeof(PetscInt), lnk));
1448: PetscFunctionReturn(PETSC_SUCCESS);
1449: }
1451: static inline PetscErrorCode PetscLLCondensedAddSorted_Scalable(PetscInt nidx, const PetscInt indices[], PetscInt lnk[])
1452: {
1453: PetscInt location = 2; /* head */
1454: PetscInt nlnk = lnk[0]; /* num of entries on the input lnk */
1456: for (PetscInt k = 0; k < nidx; k++) {
1457: const PetscInt entry = indices[k];
1458: PetscInt next, lnkdata;
1460: /* search for insertion location */
1461: do {
1462: next = location + 1; /* link from previous node to next node */
1463: location = lnk[next]; /* idx of next node */
1464: lnkdata = lnk[location]; /* value of next node */
1465: } while (entry > lnkdata);
1466: if (entry < lnkdata) {
1467: /* insertion location is found, add entry into lnk */
1468: const PetscInt newnode = 2 * (nlnk + 2); /* index for this new node */
1469: lnk[next] = newnode; /* connect previous node to the new node */
1470: lnk[newnode] = entry; /* set value of the new node */
1471: lnk[newnode + 1] = location; /* connect new node to next node */
1472: location = newnode; /* next search starts from the new node */
1473: nlnk++;
1474: }
1475: }
1476: lnk[0] = nlnk; /* number of entries in the list */
1477: return PETSC_SUCCESS;
1478: }
1480: static inline PetscErrorCode PetscLLCondensedClean_Scalable(PETSC_UNUSED PetscInt nidx, PetscInt *indices, PetscInt *lnk)
1481: {
1482: const PetscInt nlnk = lnk[0];
1483: PetscInt next = lnk[3]; /* head node */
1485: for (PetscInt k = 0; k < nlnk; k++) {
1486: indices[k] = lnk[next];
1487: next = lnk[next + 1];
1488: }
1489: lnk[0] = 0; /* num of entries on the list */
1490: lnk[3] = 2; /* head node */
1491: return PETSC_SUCCESS;
1492: }
1494: static inline PetscErrorCode PetscLLCondensedDestroy_Scalable(PetscInt *lnk)
1495: {
1496: return PetscFree(lnk);
1497: }
1499: /* -------------------------------------------------------------------------------------------------------*/
1500: /*
1501: lnk[0] number of links
1502: lnk[1] number of entries
1503: lnk[3n] value
1504: lnk[3n+1] len
1505: lnk[3n+2] link to next value
1507: The next three are always the first link
1509: lnk[3] PETSC_MIN_INT+1
1510: lnk[4] 1
1511: lnk[5] link to first real entry
1513: The next three are always the last link
1515: lnk[6] PETSC_MAX_INT - 1
1516: lnk[7] 1
1517: lnk[8] next valid link (this is the same as lnk[0] but without the decreases)
1518: */
1520: static inline PetscErrorCode PetscLLCondensedCreate_fast(PetscInt nlnk_max, PetscInt **lnk)
1521: {
1522: PetscInt *llnk;
1523: PetscInt lsize = 0;
1525: PetscFunctionBegin;
1526: PetscCall(PetscIntMultError(3, nlnk_max + 3, &lsize));
1527: PetscCall(PetscMalloc1(lsize, lnk));
1528: llnk = *lnk;
1529: llnk[0] = 0; /* nlnk: number of entries on the list */
1530: llnk[1] = 0; /* number of integer entries represented in list */
1531: llnk[3] = PETSC_MIN_INT + 1; /* value in the first node */
1532: llnk[4] = 1; /* count for the first node */
1533: llnk[5] = 6; /* next for the first node */
1534: llnk[6] = PETSC_MAX_INT - 1; /* value in the last node */
1535: llnk[7] = 1; /* count for the last node */
1536: llnk[8] = 0; /* next valid node to be used */
1537: PetscFunctionReturn(PETSC_SUCCESS);
1538: }
1540: static inline PetscErrorCode PetscLLCondensedAddSorted_fast(PetscInt nidx, const PetscInt indices[], PetscInt lnk[])
1541: {
1542: for (PetscInt k = 0, prev = 3 /* first value */; k < nidx; k++) {
1543: const PetscInt entry = indices[k];
1544: PetscInt next = lnk[prev + 2];
1546: /* search for insertion location */
1547: while (entry >= lnk[next]) {
1548: prev = next;
1549: next = lnk[next + 2];
1550: }
1551: /* entry is in range of previous list */
1552: if (entry < lnk[prev] + lnk[prev + 1]) continue;
1553: lnk[1]++;
1554: /* entry is right after previous list */
1555: if (entry == lnk[prev] + lnk[prev + 1]) {
1556: lnk[prev + 1]++;
1557: if (lnk[next] == entry + 1) { /* combine two contiguous strings */
1558: lnk[prev + 1] += lnk[next + 1];
1559: lnk[prev + 2] = lnk[next + 2];
1560: next = lnk[next + 2];
1561: lnk[0]--;
1562: }
1563: continue;
1564: }
1565: /* entry is right before next list */
1566: if (entry == lnk[next] - 1) {
1567: lnk[next]--;
1568: lnk[next + 1]++;
1569: prev = next;
1570: next = lnk[prev + 2];
1571: continue;
1572: }
1573: /* add entry into lnk */
1574: lnk[prev + 2] = 3 * ((lnk[8]++) + 3); /* connect previous node to the new node */
1575: prev = lnk[prev + 2];
1576: lnk[prev] = entry; /* set value of the new node */
1577: lnk[prev + 1] = 1; /* number of values in contiguous string is one to start */
1578: lnk[prev + 2] = next; /* connect new node to next node */
1579: lnk[0]++;
1580: }
1581: return PETSC_SUCCESS;
1582: }
1584: static inline PetscErrorCode PetscLLCondensedClean_fast(PETSC_UNUSED PetscInt nidx, PetscInt *indices, PetscInt *lnk)
1585: {
1586: const PetscInt nlnk = lnk[0];
1587: PetscInt next = lnk[5]; /* first node */
1589: for (PetscInt k = 0, cnt = 0; k < nlnk; k++) {
1590: for (PetscInt j = 0; j < lnk[next + 1]; j++) indices[cnt++] = lnk[next] + j;
1591: next = lnk[next + 2];
1592: }
1593: lnk[0] = 0; /* nlnk: number of links */
1594: lnk[1] = 0; /* number of integer entries represented in list */
1595: lnk[3] = PETSC_MIN_INT + 1; /* value in the first node */
1596: lnk[4] = 1; /* count for the first node */
1597: lnk[5] = 6; /* next for the first node */
1598: lnk[6] = PETSC_MAX_INT - 1; /* value in the last node */
1599: lnk[7] = 1; /* count for the last node */
1600: lnk[8] = 0; /* next valid location to make link */
1601: return PETSC_SUCCESS;
1602: }
1604: static inline PetscErrorCode PetscLLCondensedView_fast(const PetscInt *lnk)
1605: {
1606: const PetscInt nlnk = lnk[0];
1607: PetscInt next = lnk[5]; /* first node */
1609: for (PetscInt k = 0; k < nlnk; k++) {
1610: #if 0 /* Debugging code */
1611: printf("%d value %d len %d next %d\n", next, lnk[next], lnk[next + 1], lnk[next + 2]);
1612: #endif
1613: next = lnk[next + 2];
1614: }
1615: return PETSC_SUCCESS;
1616: }
1618: static inline PetscErrorCode PetscLLCondensedDestroy_fast(PetscInt *lnk)
1619: {
1620: return PetscFree(lnk);
1621: }
1623: /* this is extern because it is used in MatFDColoringUseDM() which is in the DM library */
1624: PETSC_EXTERN PetscErrorCode MatFDColoringApply_AIJ(Mat, MatFDColoring, Vec, void *);
1626: #if defined(PETSC_HAVE_KOKKOS_KERNELS)
1627: PETSC_INTERN PetscErrorCode MatSeqAIJMoveDiagonalValuesFront_SeqAIJKokkos(Mat, const PetscInt *);
1628: #endif
1630: PETSC_EXTERN PetscLogEvent MAT_Mult;
1631: PETSC_EXTERN PetscLogEvent MAT_MultAdd;
1632: PETSC_EXTERN PetscLogEvent MAT_MultTranspose;
1633: PETSC_EXTERN PetscLogEvent MAT_MultHermitianTranspose;
1634: PETSC_EXTERN PetscLogEvent MAT_MultTransposeAdd;
1635: PETSC_EXTERN PetscLogEvent MAT_MultHermitianTransposeAdd;
1636: PETSC_EXTERN PetscLogEvent MAT_Solve;
1637: PETSC_EXTERN PetscLogEvent MAT_Solves;
1638: PETSC_EXTERN PetscLogEvent MAT_SolveAdd;
1639: PETSC_EXTERN PetscLogEvent MAT_SolveTranspose;
1640: PETSC_EXTERN PetscLogEvent MAT_SolveTransposeAdd;
1641: PETSC_EXTERN PetscLogEvent MAT_SOR;
1642: PETSC_EXTERN PetscLogEvent MAT_ForwardSolve;
1643: PETSC_EXTERN PetscLogEvent MAT_BackwardSolve;
1644: PETSC_EXTERN PetscLogEvent MAT_LUFactor;
1645: PETSC_EXTERN PetscLogEvent MAT_LUFactorSymbolic;
1646: PETSC_EXTERN PetscLogEvent MAT_LUFactorNumeric;
1647: PETSC_EXTERN PetscLogEvent MAT_QRFactor;
1648: PETSC_EXTERN PetscLogEvent MAT_QRFactorSymbolic;
1649: PETSC_EXTERN PetscLogEvent MAT_QRFactorNumeric;
1650: PETSC_EXTERN PetscLogEvent MAT_CholeskyFactor;
1651: PETSC_EXTERN PetscLogEvent MAT_CholeskyFactorSymbolic;
1652: PETSC_EXTERN PetscLogEvent MAT_CholeskyFactorNumeric;
1653: PETSC_EXTERN PetscLogEvent MAT_ILUFactor;
1654: PETSC_EXTERN PetscLogEvent MAT_ILUFactorSymbolic;
1655: PETSC_EXTERN PetscLogEvent MAT_ICCFactorSymbolic;
1656: PETSC_EXTERN PetscLogEvent MAT_Copy;
1657: PETSC_EXTERN PetscLogEvent MAT_Convert;
1658: PETSC_EXTERN PetscLogEvent MAT_Scale;
1659: PETSC_EXTERN PetscLogEvent MAT_AssemblyBegin;
1660: PETSC_EXTERN PetscLogEvent MAT_AssemblyEnd;
1661: PETSC_EXTERN PetscLogEvent MAT_SetValues;
1662: PETSC_EXTERN PetscLogEvent MAT_GetValues;
1663: PETSC_EXTERN PetscLogEvent MAT_GetRow;
1664: PETSC_EXTERN PetscLogEvent MAT_GetRowIJ;
1665: PETSC_EXTERN PetscLogEvent MAT_CreateSubMats;
1666: PETSC_EXTERN PetscLogEvent MAT_GetOrdering;
1667: PETSC_EXTERN PetscLogEvent MAT_RedundantMat;
1668: PETSC_EXTERN PetscLogEvent MAT_IncreaseOverlap;
1669: PETSC_EXTERN PetscLogEvent MAT_Partitioning;
1670: PETSC_EXTERN PetscLogEvent MAT_PartitioningND;
1671: PETSC_EXTERN PetscLogEvent MAT_Coarsen;
1672: PETSC_EXTERN PetscLogEvent MAT_ZeroEntries;
1673: PETSC_EXTERN PetscLogEvent MAT_Load;
1674: PETSC_EXTERN PetscLogEvent MAT_View;
1675: PETSC_EXTERN PetscLogEvent MAT_AXPY;
1676: PETSC_EXTERN PetscLogEvent MAT_FDColoringCreate;
1677: PETSC_EXTERN PetscLogEvent MAT_TransposeColoringCreate;
1678: PETSC_EXTERN PetscLogEvent MAT_FDColoringSetUp;
1679: PETSC_EXTERN PetscLogEvent MAT_FDColoringApply;
1680: PETSC_EXTERN PetscLogEvent MAT_Transpose;
1681: PETSC_EXTERN PetscLogEvent MAT_FDColoringFunction;
1682: PETSC_EXTERN PetscLogEvent MAT_CreateSubMat;
1683: PETSC_EXTERN PetscLogEvent MAT_MatSolve;
1684: PETSC_EXTERN PetscLogEvent MAT_MatTrSolve;
1685: PETSC_EXTERN PetscLogEvent MAT_MatMultSymbolic;
1686: PETSC_EXTERN PetscLogEvent MAT_MatMultNumeric;
1687: PETSC_EXTERN PetscLogEvent MAT_Getlocalmatcondensed;
1688: PETSC_EXTERN PetscLogEvent MAT_GetBrowsOfAcols;
1689: PETSC_EXTERN PetscLogEvent MAT_GetBrowsOfAocols;
1690: PETSC_EXTERN PetscLogEvent MAT_PtAPSymbolic;
1691: PETSC_EXTERN PetscLogEvent MAT_PtAPNumeric;
1692: PETSC_EXTERN PetscLogEvent MAT_Seqstompinum;
1693: PETSC_EXTERN PetscLogEvent MAT_Seqstompisym;
1694: PETSC_EXTERN PetscLogEvent MAT_Seqstompi;
1695: PETSC_EXTERN PetscLogEvent MAT_Getlocalmat;
1696: PETSC_EXTERN PetscLogEvent MAT_RARtSymbolic;
1697: PETSC_EXTERN PetscLogEvent MAT_RARtNumeric;
1698: PETSC_EXTERN PetscLogEvent MAT_MatTransposeMultSymbolic;
1699: PETSC_EXTERN PetscLogEvent MAT_MatTransposeMultNumeric;
1700: PETSC_EXTERN PetscLogEvent MAT_TransposeMatMultSymbolic;
1701: PETSC_EXTERN PetscLogEvent MAT_TransposeMatMultNumeric;
1702: PETSC_EXTERN PetscLogEvent MAT_MatMatMultSymbolic;
1703: PETSC_EXTERN PetscLogEvent MAT_MatMatMultNumeric;
1704: PETSC_EXTERN PetscLogEvent MAT_Getsymtransreduced;
1705: PETSC_EXTERN PetscLogEvent MAT_GetSeqNonzeroStructure;
1706: PETSC_EXTERN PetscLogEvent MATMFFD_Mult;
1707: PETSC_EXTERN PetscLogEvent MAT_GetMultiProcBlock;
1708: PETSC_EXTERN PetscLogEvent MAT_CUSPARSECopyToGPU;
1709: PETSC_EXTERN PetscLogEvent MAT_CUSPARSECopyFromGPU;
1710: PETSC_EXTERN PetscLogEvent MAT_CUSPARSEGenerateTranspose;
1711: PETSC_EXTERN PetscLogEvent MAT_CUSPARSESolveAnalysis;
1712: PETSC_EXTERN PetscLogEvent MAT_HIPSPARSECopyToGPU;
1713: PETSC_EXTERN PetscLogEvent MAT_HIPSPARSECopyFromGPU;
1714: PETSC_EXTERN PetscLogEvent MAT_HIPSPARSEGenerateTranspose;
1715: PETSC_EXTERN PetscLogEvent MAT_HIPSPARSESolveAnalysis;
1716: PETSC_EXTERN PetscLogEvent MAT_SetValuesBatch;
1717: PETSC_EXTERN PetscLogEvent MAT_ViennaCLCopyToGPU;
1718: PETSC_EXTERN PetscLogEvent MAT_DenseCopyToGPU;
1719: PETSC_EXTERN PetscLogEvent MAT_DenseCopyFromGPU;
1720: PETSC_EXTERN PetscLogEvent MAT_Merge;
1721: PETSC_EXTERN PetscLogEvent MAT_Residual;
1722: PETSC_EXTERN PetscLogEvent MAT_SetRandom;
1723: PETSC_EXTERN PetscLogEvent MAT_FactorFactS;
1724: PETSC_EXTERN PetscLogEvent MAT_FactorInvS;
1725: PETSC_EXTERN PetscLogEvent MAT_PreallCOO;
1726: PETSC_EXTERN PetscLogEvent MAT_SetVCOO;
1727: PETSC_EXTERN PetscLogEvent MATCOLORING_Apply;
1728: PETSC_EXTERN PetscLogEvent MATCOLORING_Comm;
1729: PETSC_EXTERN PetscLogEvent MATCOLORING_Local;
1730: PETSC_EXTERN PetscLogEvent MATCOLORING_ISCreate;
1731: PETSC_EXTERN PetscLogEvent MATCOLORING_SetUp;
1732: PETSC_EXTERN PetscLogEvent MATCOLORING_Weights;
1733: PETSC_EXTERN PetscLogEvent MAT_H2Opus_Build;
1734: PETSC_EXTERN PetscLogEvent MAT_H2Opus_Compress;
1735: PETSC_EXTERN PetscLogEvent MAT_H2Opus_Orthog;
1736: PETSC_EXTERN PetscLogEvent MAT_H2Opus_LR;
1738: #endif