Actual source code: matimpl.h
1: #pragma once
3: #include <petscmat.h>
4: #include <petscmatcoarsen.h>
5: #include <petsc/private/petscimpl.h>
7: PETSC_EXTERN PetscBool MatRegisterAllCalled;
8: PETSC_EXTERN PetscBool MatSeqAIJRegisterAllCalled;
9: PETSC_EXTERN PetscBool MatOrderingRegisterAllCalled;
10: PETSC_EXTERN PetscBool MatColoringRegisterAllCalled;
11: PETSC_EXTERN PetscBool MatPartitioningRegisterAllCalled;
12: PETSC_EXTERN PetscBool MatCoarsenRegisterAllCalled;
13: PETSC_EXTERN PetscErrorCode MatRegisterAll(void);
14: PETSC_EXTERN PetscErrorCode MatOrderingRegisterAll(void);
15: PETSC_EXTERN PetscErrorCode MatColoringRegisterAll(void);
16: PETSC_EXTERN PetscErrorCode MatPartitioningRegisterAll(void);
17: PETSC_EXTERN PetscErrorCode MatCoarsenRegisterAll(void);
18: PETSC_EXTERN PetscErrorCode MatSeqAIJRegisterAll(void);
20: /* Gets the root type of the input matrix's type (e.g., MATAIJ for MATSEQAIJ) */
21: PETSC_EXTERN PetscErrorCode MatGetRootType_Private(Mat, MatType *);
23: /* Gets the MPI type corresponding to the input matrix's type (e.g., MATMPIAIJ for MATSEQAIJ) */
24: PETSC_INTERN PetscErrorCode MatGetMPIMatType_Private(Mat, MatType *);
26: /*
27: This file defines the parts of the matrix data structure that are
28: shared by all matrix types.
29: */
31: /*
32: If you add entries here also add them to the MATOP enum
33: in include/petscmat.h
34: */
35: typedef struct _MatOps *MatOps;
36: struct _MatOps {
37: /* 0*/
38: PetscErrorCode (*setvalues)(Mat, PetscInt, const PetscInt[], PetscInt, const PetscInt[], const PetscScalar[], InsertMode);
39: PetscErrorCode (*getrow)(Mat, PetscInt, PetscInt *, PetscInt *[], PetscScalar *[]);
40: PetscErrorCode (*restorerow)(Mat, PetscInt, PetscInt *, PetscInt *[], PetscScalar *[]);
41: PetscErrorCode (*mult)(Mat, Vec, Vec);
42: PetscErrorCode (*multadd)(Mat, Vec, Vec, Vec);
43: /* 5*/
44: PetscErrorCode (*multtranspose)(Mat, Vec, Vec);
45: PetscErrorCode (*multtransposeadd)(Mat, Vec, Vec, Vec);
46: PetscErrorCode (*solve)(Mat, Vec, Vec);
47: PetscErrorCode (*solveadd)(Mat, Vec, Vec, Vec);
48: PetscErrorCode (*solvetranspose)(Mat, Vec, Vec);
49: /*10*/
50: PetscErrorCode (*solvetransposeadd)(Mat, Vec, Vec, Vec);
51: PetscErrorCode (*lufactor)(Mat, IS, IS, const MatFactorInfo *);
52: PetscErrorCode (*choleskyfactor)(Mat, IS, const MatFactorInfo *);
53: PetscErrorCode (*sor)(Mat, Vec, PetscReal, MatSORType, PetscReal, PetscInt, PetscInt, Vec);
54: PetscErrorCode (*transpose)(Mat, MatReuse, Mat *);
55: /*15*/
56: PetscErrorCode (*getinfo)(Mat, MatInfoType, MatInfo *);
57: PetscErrorCode (*equal)(Mat, Mat, PetscBool *);
58: PetscErrorCode (*getdiagonal)(Mat, Vec);
59: PetscErrorCode (*diagonalscale)(Mat, Vec, Vec);
60: PetscErrorCode (*norm)(Mat, NormType, PetscReal *);
61: /*20*/
62: PetscErrorCode (*assemblybegin)(Mat, MatAssemblyType);
63: PetscErrorCode (*assemblyend)(Mat, MatAssemblyType);
64: PetscErrorCode (*setoption)(Mat, MatOption, PetscBool);
65: PetscErrorCode (*zeroentries)(Mat);
66: /*24*/
67: PetscErrorCode (*zerorows)(Mat, PetscInt, const PetscInt[], PetscScalar, Vec, Vec);
68: PetscErrorCode (*lufactorsymbolic)(Mat, Mat, IS, IS, const MatFactorInfo *);
69: PetscErrorCode (*lufactornumeric)(Mat, Mat, const MatFactorInfo *);
70: PetscErrorCode (*choleskyfactorsymbolic)(Mat, Mat, IS, const MatFactorInfo *);
71: PetscErrorCode (*choleskyfactornumeric)(Mat, Mat, const MatFactorInfo *);
72: /*29*/
73: PetscErrorCode (*setup)(Mat);
74: PetscErrorCode (*ilufactorsymbolic)(Mat, Mat, IS, IS, const MatFactorInfo *);
75: PetscErrorCode (*iccfactorsymbolic)(Mat, Mat, IS, const MatFactorInfo *);
76: PetscErrorCode (*getdiagonalblock)(Mat, Mat *);
77: PetscErrorCode (*setinf)(Mat);
78: /*34*/
79: PetscErrorCode (*duplicate)(Mat, MatDuplicateOption, Mat *);
80: PetscErrorCode (*forwardsolve)(Mat, Vec, Vec);
81: PetscErrorCode (*backwardsolve)(Mat, Vec, Vec);
82: PetscErrorCode (*ilufactor)(Mat, IS, IS, const MatFactorInfo *);
83: PetscErrorCode (*iccfactor)(Mat, IS, const MatFactorInfo *);
84: /*39*/
85: PetscErrorCode (*axpy)(Mat, PetscScalar, Mat, MatStructure);
86: PetscErrorCode (*createsubmatrices)(Mat, PetscInt, const IS[], const IS[], MatReuse, Mat *[]);
87: PetscErrorCode (*increaseoverlap)(Mat, PetscInt, IS[], PetscInt);
88: PetscErrorCode (*getvalues)(Mat, PetscInt, const PetscInt[], PetscInt, const PetscInt[], PetscScalar[]);
89: PetscErrorCode (*copy)(Mat, Mat, MatStructure);
90: /*44*/
91: PetscErrorCode (*getrowmax)(Mat, Vec, PetscInt[]);
92: PetscErrorCode (*scale)(Mat, PetscScalar);
93: PetscErrorCode (*shift)(Mat, PetscScalar);
94: PetscErrorCode (*diagonalset)(Mat, Vec, InsertMode);
95: PetscErrorCode (*zerorowscolumns)(Mat, PetscInt, const PetscInt[], PetscScalar, Vec, Vec);
96: /*49*/
97: PetscErrorCode (*setrandom)(Mat, PetscRandom);
98: PetscErrorCode (*getrowij)(Mat, PetscInt, PetscBool, PetscBool, PetscInt *, const PetscInt *[], const PetscInt *[], PetscBool *);
99: PetscErrorCode (*restorerowij)(Mat, PetscInt, PetscBool, PetscBool, PetscInt *, const PetscInt *[], const PetscInt *[], PetscBool *);
100: PetscErrorCode (*getcolumnij)(Mat, PetscInt, PetscBool, PetscBool, PetscInt *, const PetscInt *[], const PetscInt *[], PetscBool *);
101: PetscErrorCode (*restorecolumnij)(Mat, PetscInt, PetscBool, PetscBool, PetscInt *, const PetscInt *[], const PetscInt *[], PetscBool *);
102: /*54*/
103: PetscErrorCode (*fdcoloringcreate)(Mat, ISColoring, MatFDColoring);
104: PetscErrorCode (*coloringpatch)(Mat, PetscInt, PetscInt, ISColoringValue[], ISColoring *);
105: PetscErrorCode (*setunfactored)(Mat);
106: PetscErrorCode (*permute)(Mat, IS, IS, Mat *);
107: PetscErrorCode (*setvaluesblocked)(Mat, PetscInt, const PetscInt[], PetscInt, const PetscInt[], const PetscScalar[], InsertMode);
108: /*59*/
109: PetscErrorCode (*createsubmatrix)(Mat, IS, IS, MatReuse, Mat *);
110: PetscErrorCode (*destroy)(Mat);
111: PetscErrorCode (*view)(Mat, PetscViewer);
112: PetscErrorCode (*convertfrom)(Mat, MatType, MatReuse, Mat *);
113: PetscErrorCode (*matmatmultsymbolic)(Mat, Mat, Mat, PetscReal, Mat);
114: /*64*/
115: PetscErrorCode (*matmatmultnumeric)(Mat, Mat, Mat, Mat);
116: PetscErrorCode (*setlocaltoglobalmapping)(Mat, ISLocalToGlobalMapping, ISLocalToGlobalMapping);
117: PetscErrorCode (*setvalueslocal)(Mat, PetscInt, const PetscInt[], PetscInt, const PetscInt[], const PetscScalar[], InsertMode);
118: PetscErrorCode (*zerorowslocal)(Mat, PetscInt, const PetscInt[], PetscScalar, Vec, Vec);
119: PetscErrorCode (*getrowmaxabs)(Mat, Vec, PetscInt[]);
120: /*69*/
121: PetscErrorCode (*getrowminabs)(Mat, Vec, PetscInt[]);
122: PetscErrorCode (*convert)(Mat, MatType, MatReuse, Mat *);
123: PetscErrorCode (*hasoperation)(Mat, MatOperation, PetscBool *);
124: PetscErrorCode (*fdcoloringapply)(Mat, MatFDColoring, Vec, void *);
125: PetscErrorCode (*setfromoptions)(Mat, PetscOptionItems);
126: /*74*/
127: PetscErrorCode (*findzerodiagonals)(Mat, IS *);
128: PetscErrorCode (*mults)(Mat, Vecs, Vecs);
129: PetscErrorCode (*solves)(Mat, Vecs, Vecs);
130: PetscErrorCode (*getinertia)(Mat, PetscInt *, PetscInt *, PetscInt *);
131: PetscErrorCode (*load)(Mat, PetscViewer);
132: /*79*/
133: PetscErrorCode (*issymmetric)(Mat, PetscReal, PetscBool *);
134: PetscErrorCode (*ishermitian)(Mat, PetscReal, PetscBool *);
135: PetscErrorCode (*isstructurallysymmetric)(Mat, PetscBool *);
136: PetscErrorCode (*setvaluesblockedlocal)(Mat, PetscInt, const PetscInt[], PetscInt, const PetscInt[], const PetscScalar[], InsertMode);
137: PetscErrorCode (*getvecs)(Mat, Vec *, Vec *);
138: /*84*/
139: PetscErrorCode (*matmultsymbolic)(Mat, Mat, PetscReal, Mat);
140: PetscErrorCode (*matmultnumeric)(Mat, Mat, Mat);
141: PetscErrorCode (*ptapnumeric)(Mat, Mat, Mat); /* double dispatch wrapper routine */
142: PetscErrorCode (*mattransposemultsymbolic)(Mat, Mat, PetscReal, Mat);
143: PetscErrorCode (*mattransposemultnumeric)(Mat, Mat, Mat);
144: /*89*/
145: PetscErrorCode (*bindtocpu)(Mat, PetscBool);
146: PetscErrorCode (*productsetfromoptions)(Mat);
147: PetscErrorCode (*productsymbolic)(Mat);
148: PetscErrorCode (*productnumeric)(Mat);
149: PetscErrorCode (*conjugate)(Mat); /* complex conjugate */
150: /*94*/
151: PetscErrorCode (*viewnative)(Mat, PetscViewer);
152: PetscErrorCode (*setvaluesrow)(Mat, PetscInt, const PetscScalar[]);
153: PetscErrorCode (*realpart)(Mat);
154: PetscErrorCode (*imaginarypart)(Mat);
155: PetscErrorCode (*getrowuppertriangular)(Mat);
156: /*99*/
157: PetscErrorCode (*restorerowuppertriangular)(Mat);
158: PetscErrorCode (*matsolve)(Mat, Mat, Mat);
159: PetscErrorCode (*matsolvetranspose)(Mat, Mat, Mat);
160: PetscErrorCode (*getrowmin)(Mat, Vec, PetscInt[]);
161: PetscErrorCode (*getcolumnvector)(Mat, Vec, PetscInt);
162: /*104*/
163: PetscErrorCode (*missingdiagonal)(Mat, PetscBool *, PetscInt *);
164: PetscErrorCode (*getseqnonzerostructure)(Mat, Mat *);
165: PetscErrorCode (*create)(Mat);
166: PetscErrorCode (*getghosts)(Mat, PetscInt *, const PetscInt *[]);
167: PetscErrorCode (*getlocalsubmatrix)(Mat, IS, IS, Mat *);
168: /*109*/
169: PetscErrorCode (*restorelocalsubmatrix)(Mat, IS, IS, Mat *);
170: PetscErrorCode (*multdiagonalblock)(Mat, Vec, Vec);
171: PetscErrorCode (*hermitiantranspose)(Mat, MatReuse, Mat *);
172: PetscErrorCode (*multhermitiantranspose)(Mat, Vec, Vec);
173: PetscErrorCode (*multhermitiantransposeadd)(Mat, Vec, Vec, Vec);
174: /*114*/
175: PetscErrorCode (*getmultiprocblock)(Mat, MPI_Comm, MatReuse, Mat *);
176: PetscErrorCode (*findnonzerorows)(Mat, IS *);
177: PetscErrorCode (*getcolumnreductions)(Mat, PetscInt, PetscReal *);
178: PetscErrorCode (*invertblockdiagonal)(Mat, const PetscScalar **);
179: PetscErrorCode (*invertvariableblockdiagonal)(Mat, PetscInt, const PetscInt *, PetscScalar *);
180: /*119*/
181: PetscErrorCode (*createsubmatricesmpi)(Mat, PetscInt, const IS[], const IS[], MatReuse, Mat **);
182: PetscErrorCode (*setvaluesbatch)(Mat, PetscInt, PetscInt, PetscInt *, const PetscScalar *);
183: PetscErrorCode (*transposematmultsymbolic)(Mat, Mat, PetscReal, Mat);
184: PetscErrorCode (*transposematmultnumeric)(Mat, Mat, Mat);
185: PetscErrorCode (*transposecoloringcreate)(Mat, ISColoring, MatTransposeColoring);
186: /*124*/
187: PetscErrorCode (*transcoloringapplysptoden)(MatTransposeColoring, Mat, Mat);
188: PetscErrorCode (*transcoloringapplydentosp)(MatTransposeColoring, Mat, Mat);
189: PetscErrorCode (*rartnumeric)(Mat, Mat, Mat); /* double dispatch wrapper routine */
190: PetscErrorCode (*setblocksizes)(Mat, PetscInt, PetscInt);
191: PetscErrorCode (*residual)(Mat, Vec, Vec, Vec);
192: /*129*/
193: PetscErrorCode (*fdcoloringsetup)(Mat, ISColoring, MatFDColoring);
194: PetscErrorCode (*findoffblockdiagonalentries)(Mat, IS *);
195: PetscErrorCode (*creatempimatconcatenateseqmat)(MPI_Comm, Mat, PetscInt, MatReuse, Mat *);
196: PetscErrorCode (*destroysubmatrices)(PetscInt, Mat *[]);
197: PetscErrorCode (*mattransposesolve)(Mat, Mat, Mat);
198: /*134*/
199: PetscErrorCode (*getvalueslocal)(Mat, PetscInt, const PetscInt[], PetscInt, const PetscInt[], PetscScalar[]);
200: PetscErrorCode (*creategraph)(Mat, PetscBool, PetscBool, PetscReal, PetscInt, PetscInt[], Mat *);
201: PetscErrorCode (*transposesymbolic)(Mat, Mat *);
202: PetscErrorCode (*eliminatezeros)(Mat, PetscBool);
203: PetscErrorCode (*getrowsumabs)(Mat, Vec);
204: /*139*/
205: PetscErrorCode (*getfactor)(Mat, MatSolverType, MatFactorType, Mat *);
206: PetscErrorCode (*getblockdiagonal)(Mat, Mat *); // NOTE: the caller of get{block, vblock}diagonal owns the returned matrix;
207: PetscErrorCode (*getvblockdiagonal)(Mat, Mat *); // they must destroy it after use
208: PetscErrorCode (*copyhashtoxaij)(Mat, Mat);
209: PetscErrorCode (*getcurrentmemtype)(Mat, PetscMemType *);
210: };
211: /*
212: If you add MatOps entries above also add them to the MATOP enum
213: in include/petscmat.h
214: */
216: #include <petscsys.h>
218: typedef struct _p_MatRootName *MatRootName;
219: struct _p_MatRootName {
220: char *rname, *sname, *mname;
221: MatRootName next;
222: };
224: PETSC_EXTERN MatRootName MatRootNameList;
226: /*
227: Utility private matrix routines used outside Mat
228: */
229: PETSC_SINGLE_LIBRARY_INTERN PetscErrorCode MatFindNonzeroRowsOrCols_Basic(Mat, PetscBool, PetscReal, IS *);
230: PETSC_EXTERN PetscErrorCode MatShellGetScalingShifts(Mat, PetscScalar *, PetscScalar *, Vec *, Vec *, Vec *, Mat *, IS *, IS *);
232: #define MAT_SHELL_NOT_ALLOWED (void *)-1
234: /*
235: Utility private matrix routines
236: */
237: PETSC_INTERN PetscErrorCode MatConvert_Basic(Mat, MatType, MatReuse, Mat *);
238: PETSC_INTERN PetscErrorCode MatConvert_Shell(Mat, MatType, MatReuse, Mat *);
239: PETSC_INTERN PetscErrorCode MatConvertFrom_Shell(Mat, MatType, MatReuse, Mat *);
240: PETSC_INTERN PetscErrorCode MatShellSetContext_Immutable(Mat, void *);
241: PETSC_INTERN PetscErrorCode MatShellSetContextDestroy_Immutable(Mat, PetscCtxDestroyFn *);
242: PETSC_INTERN PetscErrorCode MatShellSetManageScalingShifts_Immutable(Mat);
243: PETSC_INTERN PetscErrorCode MatCopy_Basic(Mat, Mat, MatStructure);
244: PETSC_INTERN PetscErrorCode MatDiagonalSet_Default(Mat, Vec, InsertMode);
245: #if defined(PETSC_HAVE_SCALAPACK)
246: PETSC_INTERN PetscErrorCode MatConvert_Dense_ScaLAPACK(Mat, MatType, MatReuse, Mat *);
247: #endif
248: PETSC_INTERN PetscErrorCode MatSetPreallocationCOO_Basic(Mat, PetscCount, PetscInt[], PetscInt[]);
249: PETSC_INTERN PetscErrorCode MatSetValuesCOO_Basic(Mat, const PetscScalar[], InsertMode);
251: /* This can be moved to the public header after implementing some missing MatProducts */
252: PETSC_INTERN PetscErrorCode MatCreateFromISLocalToGlobalMapping(ISLocalToGlobalMapping, Mat, PetscBool, PetscBool, MatType, Mat *);
254: /* these callbacks rely on the old matrix function pointers for
255: matmat operations. They are unsafe, and should be removed.
256: However, the amount of work needed to clean up all the
257: implementations is not negligible */
258: PETSC_INTERN PetscErrorCode MatProductSymbolic_AB(Mat);
259: PETSC_INTERN PetscErrorCode MatProductNumeric_AB(Mat);
260: PETSC_INTERN PetscErrorCode MatProductSymbolic_AtB(Mat);
261: PETSC_INTERN PetscErrorCode MatProductNumeric_AtB(Mat);
262: PETSC_INTERN PetscErrorCode MatProductSymbolic_ABt(Mat);
263: PETSC_INTERN PetscErrorCode MatProductNumeric_ABt(Mat);
264: PETSC_INTERN PetscErrorCode MatProductNumeric_PtAP(Mat);
265: PETSC_INTERN PetscErrorCode MatProductNumeric_RARt(Mat);
266: PETSC_INTERN PetscErrorCode MatProductSymbolic_ABC(Mat);
267: PETSC_INTERN PetscErrorCode MatProductNumeric_ABC(Mat);
269: PETSC_INTERN PetscErrorCode MatProductCreate_Private(Mat, Mat, Mat, Mat);
270: /* this callback handles all the different triple products and
271: does not rely on the function pointers; used by cuSPARSE/hipSPARSE and KOKKOS-KERNELS */
272: PETSC_INTERN PetscErrorCode MatProductSymbolic_ABC_Basic(Mat);
274: /* CreateGraph is common to AIJ seq and mpi */
275: PETSC_INTERN PetscErrorCode MatCreateGraph_Simple_AIJ(Mat, PetscBool, PetscBool, PetscReal, PetscInt, PetscInt[], Mat *);
277: #if defined(PETSC_CLANG_STATIC_ANALYZER)
278: template <typename Tm>
279: extern void MatCheckPreallocated(Tm, int);
280: template <typename Tm>
281: extern void MatCheckProduct(Tm, int);
282: #else /* PETSC_CLANG_STATIC_ANALYZER */
283: #define MatCheckPreallocated(A, arg) \
284: do { \
285: if (!(A)->preallocated) PetscCall(MatSetUp(A)); \
286: } while (0)
288: #if defined(PETSC_USE_DEBUG)
289: #define MatCheckProduct(A, arg) \
290: do { \
291: PetscCheck((A)->product, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Argument %d \"%s\" is not a matrix obtained from MatProductCreate()", (arg), #A); \
292: } while (0)
293: #else
294: #define MatCheckProduct(A, arg) \
295: do { \
296: } while (0)
297: #endif
298: #endif /* PETSC_CLANG_STATIC_ANALYZER */
300: /*
301: The stash is used to temporarily store inserted matrix values that
302: belong to another processor. During the assembly phase the stashed
303: values are moved to the correct processor and
304: */
306: typedef struct _MatStashSpace *PetscMatStashSpace;
308: struct _MatStashSpace {
309: PetscMatStashSpace next;
310: PetscScalar *space_head, *val;
311: PetscInt *idx, *idy;
312: PetscInt total_space_size;
313: PetscInt local_used;
314: PetscInt local_remaining;
315: };
317: PETSC_EXTERN PetscErrorCode PetscMatStashSpaceGet(PetscInt, PetscInt, PetscMatStashSpace *);
318: PETSC_EXTERN PetscErrorCode PetscMatStashSpaceContiguous(PetscInt, PetscMatStashSpace *, PetscScalar *, PetscInt *, PetscInt *);
319: PETSC_EXTERN PetscErrorCode PetscMatStashSpaceDestroy(PetscMatStashSpace *);
321: typedef struct {
322: PetscInt count;
323: } MatStashHeader;
325: typedef struct {
326: void *buffer; /* Of type blocktype, dynamically constructed */
327: PetscInt count;
328: char pending;
329: } MatStashFrame;
331: typedef struct _MatStash MatStash;
332: struct _MatStash {
333: PetscInt nmax; /* maximum stash size */
334: PetscInt umax; /* user specified max-size */
335: PetscInt oldnmax; /* the nmax value used previously */
336: PetscInt n; /* stash size */
337: PetscInt bs; /* block size of the stash */
338: PetscInt reallocs; /* preserve the no of mallocs invoked */
339: PetscMatStashSpace space_head, space; /* linked list to hold stashed global row/column numbers and matrix values */
341: PetscErrorCode (*ScatterBegin)(Mat, MatStash *, PetscInt *);
342: PetscErrorCode (*ScatterGetMesg)(MatStash *, PetscMPIInt *, PetscInt **, PetscInt **, PetscScalar **, PetscInt *);
343: PetscErrorCode (*ScatterEnd)(MatStash *);
344: PetscErrorCode (*ScatterDestroy)(MatStash *);
346: /* The following variables are used for communication */
347: MPI_Comm comm;
348: PetscMPIInt size, rank;
349: PetscMPIInt tag1, tag2;
350: MPI_Request *send_waits; /* array of send requests */
351: MPI_Request *recv_waits; /* array of receive requests */
352: MPI_Status *send_status; /* array of send status */
353: PetscMPIInt nsends, nrecvs; /* numbers of sends and receives */
354: PetscScalar *svalues; /* sending data */
355: PetscInt *sindices;
356: PetscScalar **rvalues; /* receiving data (values) */
357: PetscInt **rindices; /* receiving data (indices) */
358: PetscMPIInt nprocessed; /* number of messages already processed */
359: PetscMPIInt *flg_v; /* indicates what messages have arrived so far and from whom */
360: PetscBool reproduce;
361: PetscMPIInt reproduce_count;
363: /* The following variables are used for BTS communication */
364: PetscBool first_assembly_done; /* Is the first time matrix assembly done? */
365: PetscBool use_status; /* Use MPI_Status to determine number of items in each message */
366: PetscMPIInt nsendranks;
367: PetscMPIInt nrecvranks;
368: PetscMPIInt *sendranks;
369: PetscMPIInt *recvranks;
370: MatStashHeader *sendhdr, *recvhdr;
371: MatStashFrame *sendframes; /* pointers to the main messages */
372: MatStashFrame *recvframes;
373: MatStashFrame *recvframe_active;
374: PetscInt recvframe_i; /* index of block within active frame */
375: PetscInt recvframe_count; /* Count actually sent for current frame */
376: PetscMPIInt recvcount; /* Number of receives processed so far */
377: PetscMPIInt *some_indices; /* From last call to MPI_Waitsome */
378: MPI_Status *some_statuses; /* Statuses from last call to MPI_Waitsome */
379: PetscMPIInt some_count; /* Number of requests completed in last call to MPI_Waitsome */
380: PetscMPIInt some_i; /* Index of request currently being processed */
381: MPI_Request *sendreqs;
382: MPI_Request *recvreqs;
383: PetscSegBuffer segsendblocks;
384: PetscSegBuffer segrecvframe;
385: PetscSegBuffer segrecvblocks;
386: MPI_Datatype blocktype;
387: size_t blocktype_size;
388: InsertMode *insertmode; /* Pointer to check mat->insertmode and set upon message arrival in case no local values have been set. */
389: };
391: #if !defined(PETSC_HAVE_MPIUNI)
392: PETSC_INTERN PetscErrorCode MatStashScatterDestroy_BTS(MatStash *);
393: #endif
394: PETSC_INTERN PetscErrorCode MatStashCreate_Private(MPI_Comm, PetscInt, MatStash *);
395: PETSC_INTERN PetscErrorCode MatStashDestroy_Private(MatStash *);
396: PETSC_INTERN PetscErrorCode MatStashScatterEnd_Private(MatStash *);
397: PETSC_INTERN PetscErrorCode MatStashSetInitialSize_Private(MatStash *, PetscInt);
398: PETSC_INTERN PetscErrorCode MatStashGetInfo_Private(MatStash *, PetscInt *, PetscInt *);
399: PETSC_INTERN PetscErrorCode MatStashValuesRow_Private(MatStash *, PetscInt, PetscInt, const PetscInt[], const PetscScalar[], PetscBool);
400: PETSC_INTERN PetscErrorCode MatStashValuesCol_Private(MatStash *, PetscInt, PetscInt, const PetscInt[], const PetscScalar[], PetscInt, PetscBool);
401: PETSC_INTERN PetscErrorCode MatStashValuesRowBlocked_Private(MatStash *, PetscInt, PetscInt, const PetscInt[], const PetscScalar[], PetscInt, PetscInt, PetscInt);
402: PETSC_INTERN PetscErrorCode MatStashValuesColBlocked_Private(MatStash *, PetscInt, PetscInt, const PetscInt[], const PetscScalar[], PetscInt, PetscInt, PetscInt);
403: PETSC_INTERN PetscErrorCode MatStashScatterBegin_Private(Mat, MatStash *, PetscInt *);
404: PETSC_INTERN PetscErrorCode MatStashScatterGetMesg_Private(MatStash *, PetscMPIInt *, PetscInt **, PetscInt **, PetscScalar **, PetscInt *);
405: PETSC_INTERN PetscErrorCode MatGetInfo_External(Mat, MatInfoType, MatInfo *);
407: typedef struct {
408: PetscInt dim;
409: PetscInt dims[4];
410: PetscInt starts[4];
411: PetscBool noc; /* this is a single component problem, hence user will not set MatStencil.c */
412: } MatStencilInfo;
414: /* Info about using compressed row format */
415: typedef struct {
416: PetscBool use; /* indicates compressed rows have been checked and will be used */
417: PetscInt nrows; /* number of non-zero rows */
418: PetscInt *i; /* compressed row pointer */
419: PetscInt *rindex; /* compressed row index */
420: } Mat_CompressedRow;
421: PETSC_EXTERN PetscErrorCode MatCheckCompressedRow(Mat, PetscInt, Mat_CompressedRow *, PetscInt *, PetscInt, PetscReal);
423: typedef struct { /* used by MatCreateRedundantMatrix() for reusing matredundant */
424: PetscInt nzlocal, nsends, nrecvs;
425: PetscMPIInt *send_rank, *recv_rank;
426: PetscInt *sbuf_nz, *rbuf_nz, *sbuf_j, **rbuf_j;
427: PetscScalar *sbuf_a, **rbuf_a;
428: MPI_Comm subcomm; /* when user does not provide a subcomm */
429: IS isrow, iscol;
430: Mat *matseq;
431: } Mat_Redundant;
433: typedef struct { /* used by MatProduct() */
434: MatProductType type;
435: char *alg;
436: Mat A, B, C, Dwork;
437: 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, .. */
438: 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 .. */
439: PetscBool symbolic_used_the_fact_C_is_symmetric; /* MatMatMult(A,B,MAT_REUSE_MATRIX,..&C) is still legitimate), we need to redo symbolic! */
440: PetscObjectParameterDeclare(PetscReal, fill);
441: 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 */
442: PetscBool setfromoptionscalled;
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: PetscBool use_edge_weights; /* A flag indicates whether or not to use edge weights */
552: };
554: /* needed for parallel nested dissection by ParMetis and PTSCOTCH */
555: PETSC_INTERN PetscErrorCode MatPartitioningSizesToSep_Private(PetscInt, PetscInt[], PetscInt[], PetscInt[]);
557: /*
558: Object for coarsen graphs
559: */
560: typedef struct _MatCoarsenOps *MatCoarsenOps;
561: struct _MatCoarsenOps {
562: PetscErrorCode (*apply)(MatCoarsen);
563: PetscErrorCode (*setfromoptions)(MatCoarsen, PetscOptionItems);
564: PetscErrorCode (*destroy)(MatCoarsen);
565: PetscErrorCode (*view)(MatCoarsen, PetscViewer);
566: };
568: #define MAT_COARSEN_STRENGTH_INDEX_SIZE 3
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: PetscInt max_it; /* number of iterations in HEM */
578: PetscReal threshold; /* HEM can filter interim graphs */
579: PetscInt strength_index_size;
580: PetscInt strength_index[MAT_COARSEN_STRENGTH_INDEX_SIZE];
581: };
583: PETSC_EXTERN PetscErrorCode MatCoarsenMISKSetDistance(MatCoarsen, PetscInt);
584: PETSC_EXTERN PetscErrorCode MatCoarsenMISKGetDistance(MatCoarsen, PetscInt *);
586: /*
587: Used in aijdevice.h
588: */
589: typedef struct {
590: PetscInt *i;
591: PetscInt *j;
592: PetscScalar *a;
593: PetscInt n;
594: PetscInt ignorezeroentries;
595: } PetscCSRDataStructure;
597: /*
598: MatFDColoring is used to compute Jacobian matrices efficiently
599: via coloring. The data structure is explained below in an example.
601: Color = 0 1 0 2 | 2 3 0
602: ---------------------------------------------------
603: 00 01 | 05
604: 10 11 | 14 15 Processor 0
605: 22 23 | 25
606: 32 33 |
607: ===================================================
608: | 44 45 46
609: 50 | 55 Processor 1
610: | 64 66
611: ---------------------------------------------------
613: ncolors = 4;
615: ncolumns = {2,1,1,0}
616: columns = {{0,2},{1},{3},{}}
617: nrows = {4,2,3,3}
618: rows = {{0,1,2,3},{0,1},{1,2,3},{0,1,2}}
619: vwscale = {dx(0),dx(1),dx(2),dx(3)} MPI Vec
620: vscale = {dx(0),dx(1),dx(2),dx(3),dx(4),dx(5)} Seq Vec
622: ncolumns = {1,0,1,1}
623: columns = {{6},{},{4},{5}}
624: nrows = {3,0,2,2}
625: rows = {{0,1,2},{},{1,2},{1,2}}
626: vwscale = {dx(4),dx(5),dx(6)} MPI Vec
627: vscale = {dx(0),dx(4),dx(5),dx(6)} Seq Vec
629: See the routine MatFDColoringApply() for how this data is used
630: to compute the Jacobian.
632: */
633: typedef struct {
634: PetscInt row;
635: PetscInt col;
636: PetscScalar *valaddr; /* address of value */
637: } MatEntry;
639: typedef struct {
640: PetscInt row;
641: PetscScalar *valaddr; /* address of value */
642: } MatEntry2;
644: struct _p_MatFDColoring {
645: PETSCHEADER(int);
646: PetscInt M, N, m; /* total rows, columns; local rows */
647: PetscInt rstart; /* first row owned by local processor */
648: PetscInt ncolors; /* number of colors */
649: PetscInt *ncolumns; /* number of local columns for a color */
650: PetscInt **columns; /* lists the local columns of each color (using global column numbering) */
651: IS *isa; /* these are the IS that contain the column values given in columns */
652: PetscInt *nrows; /* number of local rows for each color */
653: MatEntry *matentry; /* holds (row, column, address of value) for Jacobian matrix entry */
654: MatEntry2 *matentry2; /* holds (row, address of value) for Jacobian matrix entry */
655: PetscScalar *dy; /* store a block of F(x+dx)-F(x) when J is in BAIJ format */
656: PetscReal error_rel; /* square root of relative error in computing function */
657: PetscReal umin; /* minimum allowable u'dx value */
658: Vec w1, w2, w3; /* work vectors used in computing Jacobian */
659: PetscBool fset; /* indicates that the initial function value F(X) is set */
660: MatFDColoringFn *f; /* function that defines Jacobian */
661: void *fctx; /* optional user-defined context for use by the function f */
662: Vec vscale; /* holds FD scaling, i.e. 1/dx for each perturbed column */
663: PetscInt currentcolor; /* color for which function evaluation is being done now */
664: const char *htype; /* "wp" or "ds" */
665: ISColoringType ctype; /* IS_COLORING_GLOBAL or IS_COLORING_LOCAL */
666: PetscInt brows, bcols; /* number of block rows or columns for speedup inserting the dense matrix into sparse Jacobian */
667: PetscBool setupcalled; /* true if setup has been called */
668: PetscBool viewed; /* true if the -mat_fd_coloring_view has been triggered already */
669: void (*ftn_func_pointer)(void), *ftn_func_cntx; /* serve the same purpose as *fortran_func_pointers in PETSc objects */
670: PetscObjectId matid; /* matrix this object was created with, must always be the same */
671: };
673: typedef struct _MatColoringOps *MatColoringOps;
674: struct _MatColoringOps {
675: PetscErrorCode (*destroy)(MatColoring);
676: PetscErrorCode (*setfromoptions)(MatColoring, PetscOptionItems);
677: PetscErrorCode (*view)(MatColoring, PetscViewer);
678: PetscErrorCode (*apply)(MatColoring, ISColoring *);
679: PetscErrorCode (*weights)(MatColoring, PetscReal **, PetscInt **);
680: };
682: struct _p_MatColoring {
683: PETSCHEADER(struct _MatColoringOps);
684: Mat mat;
685: PetscInt dist; /* distance of the coloring */
686: PetscInt maxcolors; /* the maximum number of colors returned, maxcolors=1 for MIS */
687: void *data; /* inner context */
688: PetscBool valid; /* check to see if what is produced is a valid coloring */
689: MatColoringWeightType weight_type; /* type of weight computation to be performed */
690: PetscReal *user_weights; /* custom weights and permutation */
691: PetscInt *user_lperm;
692: PetscBool valid_iscoloring; /* check to see if matcoloring is produced a valid iscoloring */
693: };
695: struct _p_MatTransposeColoring {
696: PETSCHEADER(int);
697: PetscInt M, N, m; /* total rows, columns; local rows */
698: PetscInt rstart; /* first row owned by local processor */
699: PetscInt ncolors; /* number of colors */
700: PetscInt *ncolumns; /* number of local columns for a color */
701: PetscInt *nrows; /* number of local rows for each color */
702: PetscInt currentcolor; /* color for which function evaluation is being done now */
703: ISColoringType ctype; /* IS_COLORING_GLOBAL or IS_COLORING_LOCAL */
705: PetscInt *colorforrow, *colorforcol; /* pointer to rows and columns */
706: PetscInt *rows; /* lists the local rows for each color (using the local row numbering) */
707: PetscInt *den2sp; /* maps (row,color) in the dense matrix to index of sparse matrix array a->a */
708: PetscInt *columns; /* lists the local columns of each color (using global column numbering) */
709: PetscInt brows; /* number of rows for efficient implementation of MatTransColoringApplyDenToSp() */
710: PetscInt *lstart; /* array used for loop over row blocks of Csparse */
711: };
713: /*
714: Null space context for preconditioner/operators
715: */
716: struct _p_MatNullSpace {
717: PETSCHEADER(int);
718: PetscBool has_cnst;
719: PetscInt n;
720: Vec *vecs;
721: PetscScalar *alpha; /* for projections */
722: MatNullSpaceRemoveFn *remove; /* for user provided removal function */
723: void *rmctx; /* context for remove() function */
724: };
726: /*
727: Internal data structure for MATMPIDENSE
728: */
729: typedef struct {
730: Mat A; /* local submatrix */
732: /* The following variables are used for matrix assembly */
733: PetscBool donotstash; /* Flag indicating if values should be stashed */
734: MPI_Request *send_waits; /* array of send requests */
735: MPI_Request *recv_waits; /* array of receive requests */
736: PetscInt nsends, nrecvs; /* numbers of sends and receives */
737: PetscScalar *svalues, *rvalues; /* sending and receiving data */
738: PetscInt rmax; /* maximum message length */
740: /* The following variables are used for matrix-vector products */
741: Vec lvec; /* local vector */
742: PetscSF Mvctx; /* for mat-mult communications */
743: PetscBool roworiented; /* if true, row-oriented input (default) */
745: /* Support for MatDenseGetColumnVec and MatDenseGetSubMatrix */
746: Mat cmat; /* matrix representation of a given subset of columns */
747: Vec cvec; /* vector representation of a given column */
748: const PetscScalar *ptrinuse; /* holds array to be restored (just a placeholder) */
749: PetscInt vecinuse; /* if cvec is in use (col = vecinuse-1) */
750: PetscInt matinuse; /* if cmat is in use (cbegin = matinuse-1) */
751: /* if this is from MatDenseGetSubMatrix, which columns and rows does it correspond to? */
752: PetscInt sub_rbegin;
753: PetscInt sub_rend;
754: PetscInt sub_cbegin;
755: PetscInt sub_cend;
756: } Mat_MPIDense;
758: /*
759: Checking zero pivot for LU, ILU preconditioners.
760: */
761: typedef struct {
762: PetscInt nshift, nshift_max;
763: PetscReal shift_amount, shift_lo, shift_hi, shift_top, shift_fraction;
764: PetscBool newshift;
765: PetscReal rs; /* active row sum of abs(off-diagonals) */
766: PetscScalar pv; /* pivot of the active row */
767: } FactorShiftCtx;
769: PETSC_SINGLE_LIBRARY_INTERN PetscErrorCode MatTransposeCheckNonzeroState_Private(Mat, Mat);
771: /*
772: Used by MatTranspose() and potentially other functions to track the matrix used in the generation of another matrix
773: */
774: typedef struct {
775: PetscObjectId id;
776: PetscObjectState state;
777: PetscObjectState nonzerostate;
778: } MatParentState;
780: PETSC_EXTERN PetscErrorCode MatFactorDumpMatrix(Mat);
781: PETSC_INTERN PetscErrorCode MatSetBlockSizes_Default(Mat, PetscInt, PetscInt);
783: PETSC_SINGLE_LIBRARY_INTERN PetscErrorCode MatShift_Basic(Mat, PetscScalar);
785: static inline PetscErrorCode MatPivotCheck_nz(PETSC_UNUSED Mat mat, const MatFactorInfo *info, FactorShiftCtx *sctx, PETSC_UNUSED PetscInt row)
786: {
787: PetscReal _rs = sctx->rs;
788: PetscReal _zero = info->zeropivot * _rs;
790: PetscFunctionBegin;
791: if (PetscAbsScalar(sctx->pv) <= _zero && !PetscIsNanScalar(sctx->pv)) {
792: /* force |diag| > zeropivot*rs */
793: if (!sctx->nshift) sctx->shift_amount = info->shiftamount;
794: else sctx->shift_amount *= 2.0;
795: sctx->newshift = PETSC_TRUE;
796: (sctx->nshift)++;
797: } else {
798: sctx->newshift = PETSC_FALSE;
799: }
800: PetscFunctionReturn(PETSC_SUCCESS);
801: }
803: static inline PetscErrorCode MatPivotCheck_pd(PETSC_UNUSED Mat mat, const MatFactorInfo *info, FactorShiftCtx *sctx, PETSC_UNUSED PetscInt row)
804: {
805: PetscReal _rs = sctx->rs;
806: PetscReal _zero = info->zeropivot * _rs;
808: PetscFunctionBegin;
809: if (PetscRealPart(sctx->pv) <= _zero && !PetscIsNanScalar(sctx->pv)) {
810: /* force matfactor to be diagonally dominant */
811: if (sctx->nshift == sctx->nshift_max) {
812: sctx->shift_fraction = sctx->shift_hi;
813: } else {
814: sctx->shift_lo = sctx->shift_fraction;
815: sctx->shift_fraction = (sctx->shift_hi + sctx->shift_lo) / (PetscReal)2.;
816: }
817: sctx->shift_amount = sctx->shift_fraction * sctx->shift_top;
818: sctx->nshift++;
819: sctx->newshift = PETSC_TRUE;
820: } else {
821: sctx->newshift = PETSC_FALSE;
822: }
823: PetscFunctionReturn(PETSC_SUCCESS);
824: }
826: static inline PetscErrorCode MatPivotCheck_inblocks(PETSC_UNUSED Mat mat, const MatFactorInfo *info, FactorShiftCtx *sctx, PETSC_UNUSED PetscInt row)
827: {
828: PetscReal _zero = info->zeropivot;
830: PetscFunctionBegin;
831: if (PetscAbsScalar(sctx->pv) <= _zero && !PetscIsNanScalar(sctx->pv)) {
832: sctx->pv += info->shiftamount;
833: sctx->shift_amount = 0.0;
834: sctx->nshift++;
835: }
836: sctx->newshift = PETSC_FALSE;
837: PetscFunctionReturn(PETSC_SUCCESS);
838: }
840: static inline PetscErrorCode MatPivotCheck_none(Mat fact, Mat mat, const MatFactorInfo *info, FactorShiftCtx *sctx, PetscInt row)
841: {
842: PetscReal _zero = info->zeropivot;
844: PetscFunctionBegin;
845: sctx->newshift = PETSC_FALSE;
846: if (PetscAbsScalar(sctx->pv) <= _zero && !PetscIsNanScalar(sctx->pv)) {
847: 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);
848: PetscCall(PetscInfo(mat, "Detected zero pivot in factorization in row %" PetscInt_FMT " value %g tolerance %g\n", row, (double)PetscAbsScalar(sctx->pv), (double)_zero));
849: fact->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
850: fact->factorerror_zeropivot_value = PetscAbsScalar(sctx->pv);
851: fact->factorerror_zeropivot_row = row;
852: }
853: PetscFunctionReturn(PETSC_SUCCESS);
854: }
856: static inline PetscErrorCode MatPivotCheck(Mat fact, Mat mat, const MatFactorInfo *info, FactorShiftCtx *sctx, PetscInt row)
857: {
858: PetscFunctionBegin;
859: if (info->shifttype == (PetscReal)MAT_SHIFT_NONZERO) PetscCall(MatPivotCheck_nz(mat, info, sctx, row));
860: else if (info->shifttype == (PetscReal)MAT_SHIFT_POSITIVE_DEFINITE) PetscCall(MatPivotCheck_pd(mat, info, sctx, row));
861: else if (info->shifttype == (PetscReal)MAT_SHIFT_INBLOCKS) PetscCall(MatPivotCheck_inblocks(mat, info, sctx, row));
862: else PetscCall(MatPivotCheck_none(fact, mat, info, sctx, row));
863: PetscFunctionReturn(PETSC_SUCCESS);
864: }
866: #include <petscbt.h>
867: /*
868: Create and initialize a linked list
869: Input Parameters:
870: idx_start - starting index of the list
871: lnk_max - max value of lnk indicating the end of the list
872: nlnk - max length of the list
873: Output Parameters:
874: lnk - list initialized
875: bt - PetscBT (bitarray) with all bits set to false
876: lnk_empty - flg indicating the list is empty
877: */
878: #define PetscLLCreate(idx_start, lnk_max, nlnk, lnk, bt) ((PetscErrorCode)(PetscMalloc1(nlnk, &lnk) || PetscBTCreate(nlnk, &(bt)) || (lnk[idx_start] = lnk_max, PETSC_SUCCESS)))
880: #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)))
882: 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)
883: {
884: PetscInt location;
886: PetscFunctionBegin;
887: /* start from the beginning if entry < previous entry */
888: if (!assume_sorted && k && entry < *lnkdata) *lnkdata = idx_start;
889: /* search for insertion location */
890: do {
891: location = *lnkdata;
892: *lnkdata = lnk[location];
893: } while (entry > *lnkdata);
894: /* insertion location is found, add entry into lnk */
895: lnk[location] = entry;
896: lnk[entry] = *lnkdata;
897: ++(*nlnk);
898: *lnkdata = entry; /* next search starts from here if next_entry > entry */
899: PetscFunctionReturn(PETSC_SUCCESS);
900: }
902: 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)
903: {
904: PetscFunctionBegin;
905: *nlnk = 0;
906: for (PetscInt k = 0, lnkdata = idx_start; k < nidx; ++k) {
907: const PetscInt entry = indices[k];
909: if (!PetscBTLookupSet(bt, entry)) PetscCall(PetscLLInsertLocation_Private(assume_sorted, k, idx_start, entry, nlnk, &lnkdata, lnk));
910: }
911: PetscFunctionReturn(PETSC_SUCCESS);
912: }
914: /*
915: Add an index set into a sorted linked list
916: Input Parameters:
917: nidx - number of input indices
918: indices - integer array
919: idx_start - starting index of the list
920: lnk - linked list(an integer array) that is created
921: bt - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
922: output Parameters:
923: nlnk - number of newly added indices
924: lnk - the sorted(increasing order) linked list containing new and non-redundate entries from indices
925: bt - updated PetscBT (bitarray)
926: */
927: static inline PetscErrorCode PetscLLAdd(PetscInt nidx, const PetscInt *PETSC_RESTRICT indices, PetscInt idx_start, PetscInt *PETSC_RESTRICT nlnk, PetscInt *PETSC_RESTRICT lnk, PetscBT bt)
928: {
929: PetscFunctionBegin;
930: PetscCall(PetscLLAdd_Private(nidx, indices, idx_start, nlnk, lnk, bt, PETSC_FALSE));
931: PetscFunctionReturn(PETSC_SUCCESS);
932: }
934: /*
935: Add a SORTED ascending index set into a sorted linked list - same as PetscLLAdd() bus skip 'if (_k && _entry < _lnkdata) _lnkdata = idx_start;'
936: Input Parameters:
937: nidx - number of input indices
938: indices - sorted integer array
939: idx_start - starting index of the list
940: lnk - linked list(an integer array) that is created
941: bt - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
942: output Parameters:
943: nlnk - number of newly added indices
944: lnk - the sorted(increasing order) linked list containing new and non-redundate entries from indices
945: bt - updated PetscBT (bitarray)
946: */
947: static inline PetscErrorCode PetscLLAddSorted(PetscInt nidx, const PetscInt *PETSC_RESTRICT indices, PetscInt idx_start, PetscInt *PETSC_RESTRICT nlnk, PetscInt *PETSC_RESTRICT lnk, PetscBT bt)
948: {
949: PetscFunctionBegin;
950: PetscCall(PetscLLAdd_Private(nidx, indices, idx_start, nlnk, lnk, bt, PETSC_TRUE));
951: PetscFunctionReturn(PETSC_SUCCESS);
952: }
954: /*
955: Add a permuted index set into a sorted linked list
956: Input Parameters:
957: nidx - number of input indices
958: indices - integer array
959: perm - permutation of indices
960: idx_start - starting index of the list
961: lnk - linked list(an integer array) that is created
962: bt - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
963: output Parameters:
964: nlnk - number of newly added indices
965: lnk - the sorted(increasing order) linked list containing new and non-redundate entries from indices
966: bt - updated PetscBT (bitarray)
967: */
968: 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)
969: {
970: PetscFunctionBegin;
971: *nlnk = 0;
972: for (PetscInt k = 0, lnkdata = idx_start; k < nidx; ++k) {
973: const PetscInt entry = perm[indices[k]];
975: if (!PetscBTLookupSet(bt, entry)) PetscCall(PetscLLInsertLocation_Private(PETSC_FALSE, k, idx_start, entry, nlnk, &lnkdata, lnk));
976: }
977: PetscFunctionReturn(PETSC_SUCCESS);
978: }
980: #if 0
981: /* this appears to be unused? */
982: static inline PetscErrorCode PetscLLAddSorted_new(PetscInt nidx, PetscInt *indices, PetscInt idx_start, PetscBool *lnk_empty, PetscInt *nlnk, PetscInt *lnk, PetscBT bt)
983: {
984: PetscInt lnkdata = idx_start;
986: PetscFunctionBegin;
987: if (*lnk_empty) {
988: for (PetscInt k = 0; k < nidx; ++k) {
989: const PetscInt entry = indices[k], location = lnkdata;
991: PetscCall(PetscBTSet(bt,entry)); /* mark the new entry */
992: lnkdata = lnk[location];
993: /* insertion location is found, add entry into lnk */
994: lnk[location] = entry;
995: lnk[entry] = lnkdata;
996: lnkdata = entry; /* next search starts from here */
997: }
998: /* lnk[indices[nidx-1]] = lnk[idx_start];
999: lnk[idx_start] = indices[0];
1000: PetscCall(PetscBTSet(bt,indices[0]));
1001: for (_k=1; _k<nidx; _k++) {
1002: PetscCall(PetscBTSet(bt,indices[_k]));
1003: lnk[indices[_k-1]] = indices[_k];
1004: }
1005: */
1006: *nlnk = nidx;
1007: *lnk_empty = PETSC_FALSE;
1008: } else {
1009: *nlnk = 0;
1010: for (PetscInt k = 0; k < nidx; ++k) {
1011: const PetscInt entry = indices[k];
1013: if (!PetscBTLookupSet(bt,entry)) PetscCall(PetscLLInsertLocation_Private(PETSC_TRUE,k,idx_start,entry,nlnk,&lnkdata,lnk));
1014: }
1015: }
1016: PetscFunctionReturn(PETSC_SUCCESS);
1017: }
1018: #endif
1020: /*
1021: Add a SORTED index set into a sorted linked list used for LUFactorSymbolic()
1022: Same as PetscLLAddSorted() with an additional operation:
1023: count the number of input indices that are no larger than 'diag'
1024: Input Parameters:
1025: indices - sorted integer array
1026: idx_start - starting index of the list, index of pivot row
1027: lnk - linked list(an integer array) that is created
1028: bt - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
1029: diag - index of the active row in LUFactorSymbolic
1030: nzbd - number of input indices with indices <= idx_start
1031: im - im[idx_start] is initialized as num of nonzero entries in row=idx_start
1032: output Parameters:
1033: nlnk - number of newly added indices
1034: lnk - the sorted(increasing order) linked list containing new and non-redundate entries from indices
1035: bt - updated PetscBT (bitarray)
1036: im - im[idx_start]: unchanged if diag is not an entry
1037: : num of entries with indices <= diag if diag is an entry
1038: */
1039: 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)
1040: {
1041: const PetscInt nidx = im[idx_start] - nzbd; /* num of entries with idx_start < index <= diag */
1043: PetscFunctionBegin;
1044: *nlnk = 0;
1045: for (PetscInt k = 0, lnkdata = idx_start; k < nidx; ++k) {
1046: const PetscInt entry = indices[k];
1048: ++nzbd;
1049: if (entry == diag) im[idx_start] = nzbd;
1050: if (!PetscBTLookupSet(bt, entry)) PetscCall(PetscLLInsertLocation_Private(PETSC_TRUE, k, idx_start, entry, nlnk, &lnkdata, lnk));
1051: }
1052: PetscFunctionReturn(PETSC_SUCCESS);
1053: }
1055: /*
1056: Copy data on the list into an array, then initialize the list
1057: Input Parameters:
1058: idx_start - starting index of the list
1059: lnk_max - max value of lnk indicating the end of the list
1060: nlnk - number of data on the list to be copied
1061: lnk - linked list
1062: bt - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
1063: output Parameters:
1064: indices - array that contains the copied data
1065: lnk - linked list that is cleaned and initialize
1066: bt - PetscBT (bitarray) with all bits set to false
1067: */
1068: static inline PetscErrorCode PetscLLClean(PetscInt idx_start, PetscInt lnk_max, PetscInt nlnk, PetscInt *PETSC_RESTRICT lnk, PetscInt *PETSC_RESTRICT indices, PetscBT bt)
1069: {
1070: PetscFunctionBegin;
1071: for (PetscInt j = 0, idx = idx_start; j < nlnk; ++j) {
1072: idx = lnk[idx];
1073: indices[j] = idx;
1074: PetscCall(PetscBTClear(bt, idx));
1075: }
1076: lnk[idx_start] = lnk_max;
1077: PetscFunctionReturn(PETSC_SUCCESS);
1078: }
1080: /*
1081: Free memories used by the list
1082: */
1083: #define PetscLLDestroy(lnk, bt) ((PetscErrorCode)(PetscFree(lnk) || PetscBTDestroy(&(bt))))
1085: /* Routines below are used for incomplete matrix factorization */
1086: /*
1087: Create and initialize a linked list and its levels
1088: Input Parameters:
1089: idx_start - starting index of the list
1090: lnk_max - max value of lnk indicating the end of the list
1091: nlnk - max length of the list
1092: Output Parameters:
1093: lnk - list initialized
1094: lnk_lvl - array of size nlnk for storing levels of lnk
1095: bt - PetscBT (bitarray) with all bits set to false
1096: */
1097: #define PetscIncompleteLLCreate(idx_start, lnk_max, nlnk, lnk, lnk_lvl, bt) \
1098: ((PetscErrorCode)(PetscIntMultError(2, nlnk, NULL) || PetscMalloc1(2 * nlnk, &lnk) || PetscBTCreate(nlnk, &(bt)) || (lnk[idx_start] = lnk_max, lnk_lvl = lnk + nlnk, PETSC_SUCCESS)))
1100: 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)
1101: {
1102: PetscFunctionBegin;
1103: PetscCall(PetscLLInsertLocation_Private(assume_sorted, k, idx_start, entry, nlnk, lnkdata, lnk));
1104: lnklvl[entry] = newval;
1105: PetscFunctionReturn(PETSC_SUCCESS);
1106: }
1108: /*
1109: Initialize a sorted linked list used for ILU and ICC
1110: Input Parameters:
1111: nidx - number of input idx
1112: idx - integer array used for storing column indices
1113: idx_start - starting index of the list
1114: perm - indices of an IS
1115: lnk - linked list(an integer array) that is created
1116: lnklvl - levels of lnk
1117: bt - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
1118: output Parameters:
1119: nlnk - number of newly added idx
1120: lnk - the sorted(increasing order) linked list containing new and non-redundate entries from idx
1121: lnklvl - levels of lnk
1122: bt - updated PetscBT (bitarray)
1123: */
1124: 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)
1125: {
1126: PetscFunctionBegin;
1127: *nlnk = 0;
1128: for (PetscInt k = 0, lnkdata = idx_start; k < nidx; ++k) {
1129: const PetscInt entry = perm[idx[k]];
1131: if (!PetscBTLookupSet(bt, entry)) PetscCall(PetscIncompleteLLInsertLocation_Private(PETSC_FALSE, k, idx_start, entry, nlnk, &lnkdata, lnk, lnklvl, 0));
1132: }
1133: PetscFunctionReturn(PETSC_SUCCESS);
1134: }
1136: 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)
1137: {
1138: PetscFunctionBegin;
1139: *nlnk = 0;
1140: for (PetscInt k = 0, lnkdata = idx_start; k < nidx; ++k) {
1141: const PetscInt incrlev = idxlvl[k] + prow_offset + 1;
1143: if (incrlev <= level) {
1144: const PetscInt entry = idx[k];
1146: if (!PetscBTLookupSet(bt, entry)) PetscCall(PetscIncompleteLLInsertLocation_Private(assume_sorted, k, idx_start, entry, nlnk, &lnkdata, lnk, lnklvl, incrlev));
1147: else if (lnklvl[entry] > incrlev) lnklvl[entry] = incrlev; /* existing entry */
1148: }
1149: }
1150: PetscFunctionReturn(PETSC_SUCCESS);
1151: }
1153: /*
1154: Add a SORTED index set into a sorted linked list for ICC
1155: Input Parameters:
1156: nidx - number of input indices
1157: idx - sorted integer array used for storing column indices
1158: level - level of fill, e.g., ICC(level)
1159: idxlvl - level of idx
1160: idx_start - starting index of the list
1161: lnk - linked list(an integer array) that is created
1162: lnklvl - levels of lnk
1163: bt - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
1164: idxlvl_prow - idxlvl[prow], where prow is the row number of the idx
1165: output Parameters:
1166: nlnk - number of newly added indices
1167: lnk - the sorted(increasing order) linked list containing new and non-redundate entries from idx
1168: lnklvl - levels of lnk
1169: bt - updated PetscBT (bitarray)
1170: Note: the level of U(i,j) is set as lvl(i,j) = min{ lvl(i,j), lvl(prow,i)+lvl(prow,j)+1)
1171: where idx = non-zero columns of U(prow,prow+1:n-1), prow<i
1172: */
1173: 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)
1174: {
1175: PetscFunctionBegin;
1176: PetscCall(PetscIncompleteLLAdd_Private(nidx, idx, level, idxlvl, idx_start, nlnk, lnk, lnklvl, bt, idxlvl_prow, PETSC_TRUE));
1177: PetscFunctionReturn(PETSC_SUCCESS);
1178: }
1180: /*
1181: Add a SORTED index set into a sorted linked list for ILU
1182: Input Parameters:
1183: nidx - number of input indices
1184: idx - sorted 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: prow - the row number of idx
1192: output Parameters:
1193: nlnk - number of newly added idx
1194: lnk - the sorted(increasing order) linked list containing new and non-redundate entries from idx
1195: lnklvl - levels of lnk
1196: bt - updated PetscBT (bitarray)
1198: Note: the level of factor(i,j) is set as lvl(i,j) = min{ lvl(i,j), lvl(i,prow)+lvl(prow,j)+1)
1199: where idx = non-zero columns of U(prow,prow+1:n-1), prow<i
1200: */
1201: 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)
1202: {
1203: PetscFunctionBegin;
1204: PetscCall(PetscIncompleteLLAdd_Private(nidx, idx, level, idxlvl, idx_start, nlnk, lnk, lnklvl, bt, lnklvl[prow], PETSC_TRUE));
1205: PetscFunctionReturn(PETSC_SUCCESS);
1206: }
1208: /*
1209: Add a index set into a sorted linked list
1210: Input Parameters:
1211: nidx - number of input idx
1212: idx - integer array used for storing column indices
1213: level - level of fill, e.g., ICC(level)
1214: idxlvl - level of idx
1215: idx_start - starting index of the list
1216: lnk - linked list(an integer array) that is created
1217: lnklvl - levels of lnk
1218: bt - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
1219: output Parameters:
1220: nlnk - number of newly added idx
1221: lnk - the sorted(increasing order) linked list containing new and non-redundate entries from idx
1222: lnklvl - levels of lnk
1223: bt - updated PetscBT (bitarray)
1224: */
1225: 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)
1226: {
1227: PetscFunctionBegin;
1228: PetscCall(PetscIncompleteLLAdd_Private(nidx, idx, level, idxlvl, idx_start, nlnk, lnk, lnklvl, bt, 0, PETSC_FALSE));
1229: PetscFunctionReturn(PETSC_SUCCESS);
1230: }
1232: /*
1233: Add a SORTED index set into a sorted linked list
1234: Input Parameters:
1235: nidx - number of input indices
1236: idx - sorted integer array used for storing column indices
1237: level - level of fill, e.g., ICC(level)
1238: idxlvl - level of idx
1239: idx_start - starting index of the list
1240: lnk - linked list(an integer array) that is created
1241: lnklvl - levels of lnk
1242: bt - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
1243: output Parameters:
1244: nlnk - number of newly added idx
1245: lnk - the sorted(increasing order) linked list containing new and non-redundate entries from idx
1246: lnklvl - levels of lnk
1247: bt - updated PetscBT (bitarray)
1248: */
1249: 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)
1250: {
1251: PetscFunctionBegin;
1252: PetscCall(PetscIncompleteLLAdd_Private(nidx, idx, level, idxlvl, idx_start, nlnk, lnk, lnklvl, bt, 0, PETSC_TRUE));
1253: PetscFunctionReturn(PETSC_SUCCESS);
1254: }
1256: /*
1257: Copy data on the list into an array, then initialize the list
1258: Input Parameters:
1259: idx_start - starting index of the list
1260: lnk_max - max value of lnk indicating the end of the list
1261: nlnk - number of data on the list to be copied
1262: lnk - linked list
1263: lnklvl - level of lnk
1264: bt - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
1265: output Parameters:
1266: indices - array that contains the copied data
1267: lnk - linked list that is cleaned and initialize
1268: lnklvl - level of lnk that is reinitialized
1269: bt - PetscBT (bitarray) with all bits set to false
1270: */
1271: 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)
1272: {
1273: PetscFunctionBegin;
1274: for (PetscInt j = 0, idx = idx_start; j < nlnk; ++j) {
1275: idx = lnk[idx];
1276: indices[j] = idx;
1277: indiceslvl[j] = lnklvl[idx];
1278: lnklvl[idx] = -1;
1279: PetscCall(PetscBTClear(bt, idx));
1280: }
1281: lnk[idx_start] = lnk_max;
1282: PetscFunctionReturn(PETSC_SUCCESS);
1283: }
1285: /*
1286: Free memories used by the list
1287: */
1288: #define PetscIncompleteLLDestroy(lnk, bt) ((PetscErrorCode)(PetscFree(lnk) || PetscBTDestroy(&(bt))))
1290: #if !defined(PETSC_CLANG_STATIC_ANALYZER)
1291: #define MatCheckSameLocalSize(A, ar1, B, ar2) \
1292: do { \
1293: PetscCheckSameComm(A, ar1, B, ar2); \
1294: 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, \
1295: (A)->rmap->n, (A)->cmap->n, ar2, (B)->rmap->n, (B)->cmap->n); \
1296: } while (0)
1297: #define MatCheckSameSize(A, ar1, B, ar2) \
1298: do { \
1299: 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, \
1300: (A)->rmap->N, (A)->cmap->N, ar2, (B)->rmap->N, (B)->cmap->N); \
1301: MatCheckSameLocalSize(A, ar1, B, ar2); \
1302: } while (0)
1303: #else
1304: template <typename Tm>
1305: extern void MatCheckSameLocalSize(Tm, int, Tm, int);
1306: template <typename Tm>
1307: extern void MatCheckSameSize(Tm, int, Tm, int);
1308: #endif
1310: #define VecCheckMatCompatible(M, x, ar1, b, ar2) \
1311: do { \
1312: 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, \
1313: (M)->cmap->N); \
1314: 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, \
1315: (M)->rmap->N); \
1316: } while (0)
1318: /*
1319: Create and initialize a condensed linked list -
1320: same as PetscLLCreate(), but uses a scalable array 'lnk' with size of max number of entries, not O(N).
1321: Barry suggested this approach (Dec. 6, 2011):
1322: I've thought of an alternative way of representing a linked list that is efficient but doesn't have the O(N) scaling issue
1323: (it may be faster than the O(N) even sequentially due to less crazy memory access).
1325: 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
1326: 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
1327: 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
1328: 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.
1329: 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
1330: to each other so memory access is much better than using the big array.
1332: Example:
1333: nlnk_max=5, lnk_max=36:
1334: Initial list: [0, 0 | 36, 2 | 0, 0 | 0, 0 | 0, 0 | 0, 0 | 0, 0]
1335: here, head_node has index 2 with value lnk[2]=lnk_max=36,
1336: 0-th entry is used to store the number of entries in the list,
1337: The initial lnk represents head -> tail(marked by 36) with number of entries = lnk[0]=0.
1339: Now adding a sorted set {2,4}, the list becomes
1340: [2, 0 | 36, 4 |2, 6 | 4, 2 | 0, 0 | 0, 0 | 0, 0 ]
1341: represents head -> 2 -> 4 -> tail with number of entries = lnk[0]=2.
1343: Then adding a sorted set {0,3,35}, the list
1344: [5, 0 | 36, 8 | 2, 10 | 4, 12 | 0, 4 | 3, 6 | 35, 2 ]
1345: represents head -> 0 -> 2 -> 3 -> 4 -> 35 -> tail with number of entries = lnk[0]=5.
1347: Input Parameters:
1348: nlnk_max - max length of the list
1349: lnk_max - max value of the entries
1350: Output Parameters:
1351: lnk - list created and initialized
1352: bt - PetscBT (bitarray) with all bits set to false. Note: bt has size lnk_max, not nln_max!
1353: */
1354: static inline PetscErrorCode PetscLLCondensedCreate(PetscInt nlnk_max, PetscInt lnk_max, PetscInt **lnk, PetscBT *bt)
1355: {
1356: PetscInt *llnk, lsize = 0;
1358: PetscFunctionBegin;
1359: PetscCall(PetscIntMultError(2, nlnk_max + 2, &lsize));
1360: PetscCall(PetscMalloc1(lsize, lnk));
1361: PetscCall(PetscBTCreate(lnk_max, bt));
1362: llnk = *lnk;
1363: llnk[0] = 0; /* number of entries on the list */
1364: llnk[2] = lnk_max; /* value in the head node */
1365: llnk[3] = 2; /* next for the head node */
1366: PetscFunctionReturn(PETSC_SUCCESS);
1367: }
1369: /*
1370: Add a SORTED ascending index set into a sorted linked list. See PetscLLCondensedCreate() for detailed description.
1371: Input Parameters:
1372: nidx - number of input indices
1373: indices - sorted integer array
1374: lnk - condensed linked list(an integer array) that is created
1375: bt - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
1376: output Parameters:
1377: lnk - the sorted(increasing order) linked list containing previous and newly added non-redundate indices
1378: bt - updated PetscBT (bitarray)
1379: */
1380: static inline PetscErrorCode PetscLLCondensedAddSorted(PetscInt nidx, const PetscInt indices[], PetscInt lnk[], PetscBT bt)
1381: {
1382: PetscInt location = 2; /* head */
1383: PetscInt nlnk = lnk[0]; /* num of entries on the input lnk */
1385: PetscFunctionBegin;
1386: for (PetscInt k = 0; k < nidx; k++) {
1387: const PetscInt entry = indices[k];
1388: if (!PetscBTLookupSet(bt, entry)) { /* new entry */
1389: PetscInt next, lnkdata;
1391: /* search for insertion location */
1392: do {
1393: next = location + 1; /* link from previous node to next node */
1394: location = lnk[next]; /* idx of next node */
1395: lnkdata = lnk[location]; /* value of next node */
1396: } while (entry > lnkdata);
1397: /* insertion location is found, add entry into lnk */
1398: const PetscInt newnode = 2 * (nlnk + 2); /* index for this new node */
1399: lnk[next] = newnode; /* connect previous node to the new node */
1400: lnk[newnode] = entry; /* set value of the new node */
1401: lnk[newnode + 1] = location; /* connect new node to next node */
1402: location = newnode; /* next search starts from the new node */
1403: nlnk++;
1404: }
1405: }
1406: lnk[0] = nlnk; /* number of entries in the list */
1407: PetscFunctionReturn(PETSC_SUCCESS);
1408: }
1410: static inline PetscErrorCode PetscLLCondensedClean(PetscInt lnk_max, PETSC_UNUSED PetscInt nidx, PetscInt *indices, PetscInt lnk[], PetscBT bt)
1411: {
1412: const PetscInt nlnk = lnk[0]; /* num of entries on the list */
1413: PetscInt next = lnk[3]; /* head node */
1415: PetscFunctionBegin;
1416: for (PetscInt k = 0; k < nlnk; k++) {
1417: indices[k] = lnk[next];
1418: next = lnk[next + 1];
1419: PetscCall(PetscBTClear(bt, indices[k]));
1420: }
1421: lnk[0] = 0; /* num of entries on the list */
1422: lnk[2] = lnk_max; /* initialize head node */
1423: lnk[3] = 2; /* head node */
1424: PetscFunctionReturn(PETSC_SUCCESS);
1425: }
1427: static inline PetscErrorCode PetscLLCondensedView(PetscInt *lnk)
1428: {
1429: PetscFunctionBegin;
1430: PetscCall(PetscPrintf(PETSC_COMM_SELF, "LLCondensed of size %" PetscInt_FMT ", (val, next)\n", lnk[0]));
1431: 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]));
1432: PetscFunctionReturn(PETSC_SUCCESS);
1433: }
1435: /*
1436: Free memories used by the list
1437: */
1438: static inline PetscErrorCode PetscLLCondensedDestroy(PetscInt *lnk, PetscBT bt)
1439: {
1440: PetscFunctionBegin;
1441: PetscCall(PetscFree(lnk));
1442: PetscCall(PetscBTDestroy(&bt));
1443: PetscFunctionReturn(PETSC_SUCCESS);
1444: }
1446: /*
1447: Same as PetscLLCondensedCreate(), but does not use non-scalable O(lnk_max) bitarray
1448: Input Parameters:
1449: nlnk_max - max length of the list
1450: Output Parameters:
1451: lnk - list created and initialized
1452: */
1453: static inline PetscErrorCode PetscLLCondensedCreate_Scalable(PetscInt nlnk_max, PetscInt **lnk)
1454: {
1455: PetscInt *llnk, lsize = 0;
1457: PetscFunctionBegin;
1458: PetscCall(PetscIntMultError(2, nlnk_max + 2, &lsize));
1459: PetscCall(PetscMalloc1(lsize, lnk));
1460: llnk = *lnk;
1461: llnk[0] = 0; /* number of entries on the list */
1462: llnk[2] = PETSC_INT_MAX; /* value in the head node */
1463: llnk[3] = 2; /* next for the head node */
1464: PetscFunctionReturn(PETSC_SUCCESS);
1465: }
1467: static inline PetscErrorCode PetscLLCondensedExpand_Scalable(PetscInt nlnk_max, PetscInt **lnk)
1468: {
1469: PetscInt lsize = 0;
1471: PetscFunctionBegin;
1472: PetscCall(PetscIntMultError(2, nlnk_max + 2, &lsize));
1473: PetscCall(PetscRealloc(lsize * sizeof(PetscInt), lnk));
1474: PetscFunctionReturn(PETSC_SUCCESS);
1475: }
1477: static inline PetscErrorCode PetscLLCondensedAddSorted_Scalable(PetscInt nidx, const PetscInt indices[], PetscInt lnk[])
1478: {
1479: PetscInt location = 2; /* head */
1480: PetscInt nlnk = lnk[0]; /* num of entries on the input lnk */
1482: for (PetscInt k = 0; k < nidx; k++) {
1483: const PetscInt entry = indices[k];
1484: PetscInt next, lnkdata;
1486: /* search for insertion location */
1487: do {
1488: next = location + 1; /* link from previous node to next node */
1489: location = lnk[next]; /* idx of next node */
1490: lnkdata = lnk[location]; /* value of next node */
1491: } while (entry > lnkdata);
1492: if (entry < lnkdata) {
1493: /* insertion location is found, add entry into lnk */
1494: const PetscInt newnode = 2 * (nlnk + 2); /* index for this new node */
1495: lnk[next] = newnode; /* connect previous node to the new node */
1496: lnk[newnode] = entry; /* set value of the new node */
1497: lnk[newnode + 1] = location; /* connect new node to next node */
1498: location = newnode; /* next search starts from the new node */
1499: nlnk++;
1500: }
1501: }
1502: lnk[0] = nlnk; /* number of entries in the list */
1503: return PETSC_SUCCESS;
1504: }
1506: static inline PetscErrorCode PetscLLCondensedClean_Scalable(PETSC_UNUSED PetscInt nidx, PetscInt *indices, PetscInt *lnk)
1507: {
1508: const PetscInt nlnk = lnk[0];
1509: PetscInt next = lnk[3]; /* head node */
1511: for (PetscInt k = 0; k < nlnk; k++) {
1512: indices[k] = lnk[next];
1513: next = lnk[next + 1];
1514: }
1515: lnk[0] = 0; /* num of entries on the list */
1516: lnk[3] = 2; /* head node */
1517: return PETSC_SUCCESS;
1518: }
1520: static inline PetscErrorCode PetscLLCondensedDestroy_Scalable(PetscInt *lnk)
1521: {
1522: return PetscFree(lnk);
1523: }
1525: /*
1526: lnk[0] number of links
1527: lnk[1] number of entries
1528: lnk[3n] value
1529: lnk[3n+1] len
1530: lnk[3n+2] link to next value
1532: The next three are always the first link
1534: lnk[3] PETSC_INT_MIN+1
1535: lnk[4] 1
1536: lnk[5] link to first real entry
1538: The next three are always the last link
1540: lnk[6] PETSC_INT_MAX - 1
1541: lnk[7] 1
1542: lnk[8] next valid link (this is the same as lnk[0] but without the decreases)
1543: */
1545: static inline PetscErrorCode PetscLLCondensedCreate_fast(PetscInt nlnk_max, PetscInt **lnk)
1546: {
1547: PetscInt *llnk;
1548: PetscInt lsize = 0;
1550: PetscFunctionBegin;
1551: PetscCall(PetscIntMultError(3, nlnk_max + 3, &lsize));
1552: PetscCall(PetscMalloc1(lsize, lnk));
1553: llnk = *lnk;
1554: llnk[0] = 0; /* nlnk: number of entries on the list */
1555: llnk[1] = 0; /* number of integer entries represented in list */
1556: llnk[3] = PETSC_INT_MIN + 1; /* value in the first node */
1557: llnk[4] = 1; /* count for the first node */
1558: llnk[5] = 6; /* next for the first node */
1559: llnk[6] = PETSC_INT_MAX - 1; /* value in the last node */
1560: llnk[7] = 1; /* count for the last node */
1561: llnk[8] = 0; /* next valid node to be used */
1562: PetscFunctionReturn(PETSC_SUCCESS);
1563: }
1565: static inline PetscErrorCode PetscLLCondensedAddSorted_fast(PetscInt nidx, const PetscInt indices[], PetscInt lnk[])
1566: {
1567: for (PetscInt k = 0, prev = 3 /* first value */; k < nidx; k++) {
1568: const PetscInt entry = indices[k];
1569: PetscInt next = lnk[prev + 2];
1571: /* search for insertion location */
1572: while (entry >= lnk[next]) {
1573: prev = next;
1574: next = lnk[next + 2];
1575: }
1576: /* entry is in range of previous list */
1577: if (entry < lnk[prev] + lnk[prev + 1]) continue;
1578: lnk[1]++;
1579: /* entry is right after previous list */
1580: if (entry == lnk[prev] + lnk[prev + 1]) {
1581: lnk[prev + 1]++;
1582: if (lnk[next] == entry + 1) { /* combine two contiguous strings */
1583: lnk[prev + 1] += lnk[next + 1];
1584: lnk[prev + 2] = lnk[next + 2];
1585: next = lnk[next + 2];
1586: lnk[0]--;
1587: }
1588: continue;
1589: }
1590: /* entry is right before next list */
1591: if (entry == lnk[next] - 1) {
1592: lnk[next]--;
1593: lnk[next + 1]++;
1594: prev = next;
1595: next = lnk[prev + 2];
1596: continue;
1597: }
1598: /* add entry into lnk */
1599: lnk[prev + 2] = 3 * ((lnk[8]++) + 3); /* connect previous node to the new node */
1600: prev = lnk[prev + 2];
1601: lnk[prev] = entry; /* set value of the new node */
1602: lnk[prev + 1] = 1; /* number of values in contiguous string is one to start */
1603: lnk[prev + 2] = next; /* connect new node to next node */
1604: lnk[0]++;
1605: }
1606: return PETSC_SUCCESS;
1607: }
1609: static inline PetscErrorCode PetscLLCondensedClean_fast(PETSC_UNUSED PetscInt nidx, PetscInt *indices, PetscInt *lnk)
1610: {
1611: const PetscInt nlnk = lnk[0];
1612: PetscInt next = lnk[5]; /* first node */
1614: for (PetscInt k = 0, cnt = 0; k < nlnk; k++) {
1615: for (PetscInt j = 0; j < lnk[next + 1]; j++) indices[cnt++] = lnk[next] + j;
1616: next = lnk[next + 2];
1617: }
1618: lnk[0] = 0; /* nlnk: number of links */
1619: lnk[1] = 0; /* number of integer entries represented in list */
1620: lnk[3] = PETSC_INT_MIN + 1; /* value in the first node */
1621: lnk[4] = 1; /* count for the first node */
1622: lnk[5] = 6; /* next for the first node */
1623: lnk[6] = PETSC_INT_MAX - 1; /* value in the last node */
1624: lnk[7] = 1; /* count for the last node */
1625: lnk[8] = 0; /* next valid location to make link */
1626: return PETSC_SUCCESS;
1627: }
1629: static inline PetscErrorCode PetscLLCondensedView_fast(const PetscInt *lnk)
1630: {
1631: const PetscInt nlnk = lnk[0];
1632: PetscInt next = lnk[5]; /* first node */
1634: for (PetscInt k = 0; k < nlnk; k++) {
1635: #if 0 /* Debugging code */
1636: printf("%d value %d len %d next %d\n", next, lnk[next], lnk[next + 1], lnk[next + 2]);
1637: #endif
1638: next = lnk[next + 2];
1639: }
1640: return PETSC_SUCCESS;
1641: }
1643: static inline PetscErrorCode PetscLLCondensedDestroy_fast(PetscInt *lnk)
1644: {
1645: return PetscFree(lnk);
1646: }
1648: PETSC_EXTERN PetscErrorCode PetscCDCreate(PetscInt, PetscCoarsenData **);
1649: PETSC_EXTERN PetscErrorCode PetscCDDestroy(PetscCoarsenData *);
1650: PETSC_EXTERN PetscErrorCode PetscCDIntNdSetID(PetscCDIntNd *, PetscInt);
1651: PETSC_EXTERN PetscErrorCode PetscCDIntNdGetID(const PetscCDIntNd *, PetscInt *);
1652: PETSC_EXTERN PetscErrorCode PetscCDAppendID(PetscCoarsenData *, PetscInt, PetscInt);
1653: PETSC_EXTERN PetscErrorCode PetscCDMoveAppend(PetscCoarsenData *, PetscInt, PetscInt);
1654: PETSC_EXTERN PetscErrorCode PetscCDAppendNode(PetscCoarsenData *, PetscInt, PetscCDIntNd *);
1655: PETSC_EXTERN PetscErrorCode PetscCDRemoveNextNode(PetscCoarsenData *, PetscInt, PetscCDIntNd *);
1656: PETSC_EXTERN PetscErrorCode PetscCDCountAt(const PetscCoarsenData *, PetscInt, PetscInt *);
1657: PETSC_EXTERN PetscErrorCode PetscCDIsEmptyAt(const PetscCoarsenData *, PetscInt, PetscBool *);
1658: PETSC_EXTERN PetscErrorCode PetscCDSetChunkSize(PetscCoarsenData *, PetscInt);
1659: PETSC_EXTERN PetscErrorCode PetscCDPrint(const PetscCoarsenData *, PetscInt, MPI_Comm);
1660: PETSC_EXTERN PetscErrorCode PetscCDGetNonemptyIS(PetscCoarsenData *, IS *);
1661: PETSC_EXTERN PetscErrorCode PetscCDGetMat(PetscCoarsenData *, Mat *);
1662: PETSC_EXTERN PetscErrorCode PetscCDSetMat(PetscCoarsenData *, Mat);
1663: PETSC_EXTERN PetscErrorCode PetscCDClearMat(PetscCoarsenData *);
1664: PETSC_EXTERN PetscErrorCode PetscCDRemoveAllAt(PetscCoarsenData *, PetscInt);
1665: PETSC_EXTERN PetscErrorCode PetscCDCount(const PetscCoarsenData *, PetscInt *_sz);
1667: PETSC_EXTERN PetscErrorCode PetscCDGetHeadPos(const PetscCoarsenData *, PetscInt, PetscCDIntNd **);
1668: PETSC_EXTERN PetscErrorCode PetscCDGetNextPos(const PetscCoarsenData *, PetscInt, PetscCDIntNd **);
1669: PETSC_EXTERN PetscErrorCode PetscCDGetASMBlocks(const PetscCoarsenData *, const PetscInt, PetscInt *, IS **);
1671: PETSC_SINGLE_LIBRARY_VISIBILITY_INTERNAL PetscErrorCode MatFDColoringApply_AIJ(Mat, MatFDColoring, Vec, void *);
1673: PETSC_EXTERN PetscLogEvent MAT_Mult;
1674: PETSC_EXTERN PetscLogEvent MAT_MultAdd;
1675: PETSC_EXTERN PetscLogEvent MAT_MultTranspose;
1676: PETSC_EXTERN PetscLogEvent MAT_MultHermitianTranspose;
1677: PETSC_EXTERN PetscLogEvent MAT_MultTransposeAdd;
1678: PETSC_EXTERN PetscLogEvent MAT_MultHermitianTransposeAdd;
1679: PETSC_EXTERN PetscLogEvent MAT_Solve;
1680: PETSC_EXTERN PetscLogEvent MAT_Solves;
1681: PETSC_EXTERN PetscLogEvent MAT_SolveAdd;
1682: PETSC_EXTERN PetscLogEvent MAT_SolveTranspose;
1683: PETSC_EXTERN PetscLogEvent MAT_SolveTransposeAdd;
1684: PETSC_EXTERN PetscLogEvent MAT_SOR;
1685: PETSC_EXTERN PetscLogEvent MAT_ForwardSolve;
1686: PETSC_EXTERN PetscLogEvent MAT_BackwardSolve;
1687: PETSC_EXTERN PetscLogEvent MAT_LUFactor;
1688: PETSC_EXTERN PetscLogEvent MAT_LUFactorSymbolic;
1689: PETSC_EXTERN PetscLogEvent MAT_LUFactorNumeric;
1690: PETSC_EXTERN PetscLogEvent MAT_QRFactor;
1691: PETSC_EXTERN PetscLogEvent MAT_QRFactorSymbolic;
1692: PETSC_EXTERN PetscLogEvent MAT_QRFactorNumeric;
1693: PETSC_EXTERN PetscLogEvent MAT_CholeskyFactor;
1694: PETSC_EXTERN PetscLogEvent MAT_CholeskyFactorSymbolic;
1695: PETSC_EXTERN PetscLogEvent MAT_CholeskyFactorNumeric;
1696: PETSC_EXTERN PetscLogEvent MAT_ILUFactor;
1697: PETSC_EXTERN PetscLogEvent MAT_ILUFactorSymbolic;
1698: PETSC_EXTERN PetscLogEvent MAT_ICCFactorSymbolic;
1699: PETSC_EXTERN PetscLogEvent MAT_Copy;
1700: PETSC_EXTERN PetscLogEvent MAT_Convert;
1701: PETSC_EXTERN PetscLogEvent MAT_Scale;
1702: PETSC_EXTERN PetscLogEvent MAT_AssemblyBegin;
1703: PETSC_EXTERN PetscLogEvent MAT_AssemblyEnd;
1704: PETSC_EXTERN PetscLogEvent MAT_SetValues;
1705: PETSC_EXTERN PetscLogEvent MAT_GetValues;
1706: PETSC_EXTERN PetscLogEvent MAT_GetRow;
1707: PETSC_EXTERN PetscLogEvent MAT_GetRowIJ;
1708: PETSC_EXTERN PetscLogEvent MAT_CreateSubMats;
1709: PETSC_EXTERN PetscLogEvent MAT_GetOrdering;
1710: PETSC_EXTERN PetscLogEvent MAT_RedundantMat;
1711: PETSC_EXTERN PetscLogEvent MAT_IncreaseOverlap;
1712: PETSC_EXTERN PetscLogEvent MAT_Partitioning;
1713: PETSC_EXTERN PetscLogEvent MAT_PartitioningND;
1714: PETSC_EXTERN PetscLogEvent MAT_Coarsen;
1715: PETSC_EXTERN PetscLogEvent MAT_ZeroEntries;
1716: PETSC_EXTERN PetscLogEvent MAT_Load;
1717: PETSC_EXTERN PetscLogEvent MAT_View;
1718: PETSC_EXTERN PetscLogEvent MAT_AXPY;
1719: PETSC_EXTERN PetscLogEvent MAT_FDColoringCreate;
1720: PETSC_EXTERN PetscLogEvent MAT_TransposeColoringCreate;
1721: PETSC_EXTERN PetscLogEvent MAT_FDColoringSetUp;
1722: PETSC_EXTERN PetscLogEvent MAT_FDColoringApply;
1723: PETSC_EXTERN PetscLogEvent MAT_Transpose;
1724: PETSC_EXTERN PetscLogEvent MAT_FDColoringFunction;
1725: PETSC_EXTERN PetscLogEvent MAT_CreateSubMat;
1726: PETSC_EXTERN PetscLogEvent MAT_MatSolve;
1727: PETSC_EXTERN PetscLogEvent MAT_MatTrSolve;
1728: PETSC_EXTERN PetscLogEvent MAT_MatMultSymbolic;
1729: PETSC_EXTERN PetscLogEvent MAT_MatMultNumeric;
1730: PETSC_EXTERN PetscLogEvent MAT_Getlocalmatcondensed;
1731: PETSC_EXTERN PetscLogEvent MAT_GetBrowsOfAcols;
1732: PETSC_EXTERN PetscLogEvent MAT_GetBrowsOfAocols;
1733: PETSC_EXTERN PetscLogEvent MAT_PtAPSymbolic;
1734: PETSC_EXTERN PetscLogEvent MAT_PtAPNumeric;
1735: PETSC_EXTERN PetscLogEvent MAT_Seqstompinum;
1736: PETSC_EXTERN PetscLogEvent MAT_Seqstompisym;
1737: PETSC_EXTERN PetscLogEvent MAT_Seqstompi;
1738: PETSC_EXTERN PetscLogEvent MAT_Getlocalmat;
1739: PETSC_EXTERN PetscLogEvent MAT_RARtSymbolic;
1740: PETSC_EXTERN PetscLogEvent MAT_RARtNumeric;
1741: PETSC_EXTERN PetscLogEvent MAT_MatTransposeMultSymbolic;
1742: PETSC_EXTERN PetscLogEvent MAT_MatTransposeMultNumeric;
1743: PETSC_EXTERN PetscLogEvent MAT_TransposeMatMultSymbolic;
1744: PETSC_EXTERN PetscLogEvent MAT_TransposeMatMultNumeric;
1745: PETSC_EXTERN PetscLogEvent MAT_MatMatMultSymbolic;
1746: PETSC_EXTERN PetscLogEvent MAT_MatMatMultNumeric;
1747: PETSC_EXTERN PetscLogEvent MAT_Getsymtransreduced;
1748: PETSC_EXTERN PetscLogEvent MAT_GetSeqNonzeroStructure;
1749: PETSC_EXTERN PetscLogEvent MATMFFD_Mult;
1750: PETSC_EXTERN PetscLogEvent MAT_GetMultiProcBlock;
1751: PETSC_EXTERN PetscLogEvent MAT_CUSPARSECopyToGPU;
1752: PETSC_EXTERN PetscLogEvent MAT_CUSPARSECopyFromGPU;
1753: PETSC_EXTERN PetscLogEvent MAT_CUSPARSEGenerateTranspose;
1754: PETSC_EXTERN PetscLogEvent MAT_CUSPARSESolveAnalysis;
1755: PETSC_EXTERN PetscLogEvent MAT_HIPSPARSECopyToGPU;
1756: PETSC_EXTERN PetscLogEvent MAT_HIPSPARSECopyFromGPU;
1757: PETSC_EXTERN PetscLogEvent MAT_HIPSPARSEGenerateTranspose;
1758: PETSC_EXTERN PetscLogEvent MAT_HIPSPARSESolveAnalysis;
1759: PETSC_EXTERN PetscLogEvent MAT_SetValuesBatch;
1760: PETSC_EXTERN PetscLogEvent MAT_CreateGraph;
1761: PETSC_EXTERN PetscLogEvent MAT_ViennaCLCopyToGPU;
1762: PETSC_EXTERN PetscLogEvent MAT_DenseCopyToGPU;
1763: PETSC_EXTERN PetscLogEvent MAT_DenseCopyFromGPU;
1764: PETSC_EXTERN PetscLogEvent MAT_Merge;
1765: PETSC_EXTERN PetscLogEvent MAT_Residual;
1766: PETSC_EXTERN PetscLogEvent MAT_SetRandom;
1767: PETSC_EXTERN PetscLogEvent MAT_FactorFactS;
1768: PETSC_EXTERN PetscLogEvent MAT_FactorInvS;
1769: PETSC_EXTERN PetscLogEvent MAT_PreallCOO;
1770: PETSC_EXTERN PetscLogEvent MAT_SetVCOO;
1771: PETSC_EXTERN PetscLogEvent MATCOLORING_Apply;
1772: PETSC_EXTERN PetscLogEvent MATCOLORING_Comm;
1773: PETSC_EXTERN PetscLogEvent MATCOLORING_Local;
1774: PETSC_EXTERN PetscLogEvent MATCOLORING_ISCreate;
1775: PETSC_EXTERN PetscLogEvent MATCOLORING_SetUp;
1776: PETSC_EXTERN PetscLogEvent MATCOLORING_Weights;
1777: PETSC_EXTERN PetscLogEvent MAT_H2Opus_Build;
1778: PETSC_EXTERN PetscLogEvent MAT_H2Opus_Compress;
1779: PETSC_EXTERN PetscLogEvent MAT_H2Opus_Orthog;
1780: PETSC_EXTERN PetscLogEvent MAT_H2Opus_LR;
1781: PETSC_EXTERN PetscLogEvent MAT_CUDACopyToGPU;
1782: PETSC_EXTERN PetscLogEvent MAT_HIPCopyToGPU;