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