Actual source code: matimpl.h


  2: #ifndef __MATIMPL_H

  5: #include <petscmat.h>
  6: #include <petscmatcoarsen.h>
  7: #include <petsc/private/petscimpl.h>

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

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

 25: /*
 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 (*placeholder_77)(void);
130:   PetscErrorCode (*placeholder_78)(void);
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,PetscCount,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: template <typename Tm> void MatCheckPreallocated(Tm,int);
269: template <typename Tm> void MatCheckProduct(Tm,int);
270: #else/* PETSC_CLANG_STATIC_ANALYZER */
271: #if defined(PETSC_USE_DEBUG)
272: #  define MatCheckPreallocated(A,arg) do {                              \
274:   } while (0)
275: #else
276: #  define MatCheckPreallocated(A,arg) do {} while (0)
277: #endif

279: #if defined(PETSC_USE_DEBUG)
280: #  define MatCheckProduct(A,arg) do {                              \
282:   } while (0)
283: #else
284: #  define MatCheckProduct(A,arg) do {} while (0)
285: #endif
286: #endif /* PETSC_CLANG_STATIC_ANALYZER */

288: /*
289:   The stash is used to temporarily store inserted matrix values that
290:   belong to another processor. During the assembly phase the stashed
291:   values are moved to the correct processor and
292: */

294: typedef struct _MatStashSpace *PetscMatStashSpace;

296: struct _MatStashSpace {
297:   PetscMatStashSpace next;
298:   PetscScalar        *space_head,*val;
299:   PetscInt           *idx,*idy;
300:   PetscInt           total_space_size;
301:   PetscInt           local_used;
302:   PetscInt           local_remaining;
303: };

305: PETSC_EXTERN PetscErrorCode PetscMatStashSpaceGet(PetscInt,PetscInt,PetscMatStashSpace *);
306: PETSC_EXTERN PetscErrorCode PetscMatStashSpaceContiguous(PetscInt,PetscMatStashSpace *,PetscScalar *,PetscInt *,PetscInt *);
307: PETSC_EXTERN PetscErrorCode PetscMatStashSpaceDestroy(PetscMatStashSpace*);

309: typedef struct {
310:   PetscInt    count;
311: } MatStashHeader;

313: typedef struct {
314:   void        *buffer;          /* Of type blocktype, dynamically constructed  */
315:   PetscInt    count;
316:   char        pending;
317: } MatStashFrame;

319: typedef struct _MatStash MatStash;
320: struct _MatStash {
321:   PetscInt      nmax;                   /* maximum stash size */
322:   PetscInt      umax;                   /* user specified max-size */
323:   PetscInt      oldnmax;                /* the nmax value used previously */
324:   PetscInt      n;                      /* stash size */
325:   PetscInt      bs;                     /* block size of the stash */
326:   PetscInt      reallocs;               /* preserve the no of mallocs invoked */
327:   PetscMatStashSpace space_head,space;  /* linked list to hold stashed global row/column numbers and matrix values */

329:   PetscErrorCode (*ScatterBegin)(Mat,MatStash*,PetscInt*);
330:   PetscErrorCode (*ScatterGetMesg)(MatStash*,PetscMPIInt*,PetscInt**,PetscInt**,PetscScalar**,PetscInt*);
331:   PetscErrorCode (*ScatterEnd)(MatStash*);
332:   PetscErrorCode (*ScatterDestroy)(MatStash*);

334:   /* The following variables are used for communication */
335:   MPI_Comm      comm;
336:   PetscMPIInt   size,rank;
337:   PetscMPIInt   tag1,tag2;
338:   MPI_Request   *send_waits;            /* array of send requests */
339:   MPI_Request   *recv_waits;            /* array of receive requests */
340:   MPI_Status    *send_status;           /* array of send status */
341:   PetscInt      nsends,nrecvs;          /* numbers of sends and receives */
342:   PetscScalar   *svalues;               /* sending data */
343:   PetscInt      *sindices;
344:   PetscScalar   **rvalues;              /* receiving data (values) */
345:   PetscInt      **rindices;             /* receiving data (indices) */
346:   PetscInt      nprocessed;             /* number of messages already processed */
347:   PetscMPIInt   *flg_v;                 /* indicates what messages have arrived so far and from whom */
348:   PetscBool     reproduce;
349:   PetscInt      reproduce_count;

351:   /* The following variables are used for BTS communication */
352:   PetscBool      first_assembly_done;   /* Is the first time matrix assembly done? */
353:   PetscBool      use_status;            /* Use MPI_Status to determine number of items in each message */
354:   PetscMPIInt    nsendranks;
355:   PetscMPIInt    nrecvranks;
356:   PetscMPIInt    *sendranks;
357:   PetscMPIInt    *recvranks;
358:   MatStashHeader *sendhdr,*recvhdr;
359:   MatStashFrame  *sendframes;   /* pointers to the main messages */
360:   MatStashFrame  *recvframes;
361:   MatStashFrame  *recvframe_active;
362:   PetscInt       recvframe_i;     /* index of block within active frame */
363:   PetscMPIInt    recvframe_count; /* Count actually sent for current frame */
364:   PetscInt       recvcount;       /* Number of receives processed so far */
365:   PetscMPIInt    *some_indices;   /* From last call to MPI_Waitsome */
366:   MPI_Status     *some_statuses;  /* Statuses from last call to MPI_Waitsome */
367:   PetscMPIInt    some_count;      /* Number of requests completed in last call to MPI_Waitsome */
368:   PetscMPIInt    some_i;          /* Index of request currently being processed */
369:   MPI_Request    *sendreqs;
370:   MPI_Request    *recvreqs;
371:   PetscSegBuffer segsendblocks;
372:   PetscSegBuffer segrecvframe;
373:   PetscSegBuffer segrecvblocks;
374:   MPI_Datatype   blocktype;
375:   size_t         blocktype_size;
376:   InsertMode     *insertmode;   /* Pointer to check mat->insertmode and set upon message arrival in case no local values have been set. */
377: };

379: #if !defined(PETSC_HAVE_MPIUNI)
380: PETSC_INTERN PetscErrorCode MatStashScatterDestroy_BTS(MatStash*);
381: #endif
382: PETSC_INTERN PetscErrorCode MatStashCreate_Private(MPI_Comm,PetscInt,MatStash*);
383: PETSC_INTERN PetscErrorCode MatStashDestroy_Private(MatStash*);
384: PETSC_INTERN PetscErrorCode MatStashScatterEnd_Private(MatStash*);
385: PETSC_INTERN PetscErrorCode MatStashSetInitialSize_Private(MatStash*,PetscInt);
386: PETSC_INTERN PetscErrorCode MatStashGetInfo_Private(MatStash*,PetscInt*,PetscInt*);
387: PETSC_INTERN PetscErrorCode MatStashValuesRow_Private(MatStash*,PetscInt,PetscInt,const PetscInt[],const PetscScalar[],PetscBool);
388: PETSC_INTERN PetscErrorCode MatStashValuesCol_Private(MatStash*,PetscInt,PetscInt,const PetscInt[],const PetscScalar[],PetscInt,PetscBool);
389: PETSC_INTERN PetscErrorCode MatStashValuesRowBlocked_Private(MatStash*,PetscInt,PetscInt,const PetscInt[],const PetscScalar[],PetscInt,PetscInt,PetscInt);
390: PETSC_INTERN PetscErrorCode MatStashValuesColBlocked_Private(MatStash*,PetscInt,PetscInt,const PetscInt[],const PetscScalar[],PetscInt,PetscInt,PetscInt);
391: PETSC_INTERN PetscErrorCode MatStashScatterBegin_Private(Mat,MatStash*,PetscInt*);
392: PETSC_INTERN PetscErrorCode MatStashScatterGetMesg_Private(MatStash*,PetscMPIInt*,PetscInt**,PetscInt**,PetscScalar**,PetscInt*);
393: PETSC_INTERN PetscErrorCode MatGetInfo_External(Mat,MatInfoType,MatInfo*);

395: typedef struct {
396:   PetscInt   dim;
397:   PetscInt   dims[4];
398:   PetscInt   starts[4];
399:   PetscBool  noc;        /* this is a single component problem, hence user will not set MatStencil.c */
400: } MatStencilInfo;

402: /* Info about using compressed row format */
403: typedef struct {
404:   PetscBool  use;                           /* indicates compressed rows have been checked and will be used */
405:   PetscInt   nrows;                         /* number of non-zero rows */
406:   PetscInt   *i;                            /* compressed row pointer  */
407:   PetscInt   *rindex;                       /* compressed row index               */
408: } Mat_CompressedRow;
409: PETSC_EXTERN PetscErrorCode MatCheckCompressedRow(Mat,PetscInt,Mat_CompressedRow*,PetscInt*,PetscInt,PetscReal);

411: typedef struct { /* used by MatCreateRedundantMatrix() for reusing matredundant */
412:   PetscInt     nzlocal,nsends,nrecvs;
413:   PetscMPIInt  *send_rank,*recv_rank;
414:   PetscInt     *sbuf_nz,*rbuf_nz,*sbuf_j,**rbuf_j;
415:   PetscScalar  *sbuf_a,**rbuf_a;
416:   MPI_Comm     subcomm;   /* when user does not provide a subcomm */
417:   IS           isrow,iscol;
418:   Mat          *matseq;
419: } Mat_Redundant;

421: typedef struct { /* used by MatProduct() */
422:   MatProductType type;
423:   char           *alg;
424:   Mat            A,B,C,Dwork;
425:   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, .. */
426:   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 .. */
427:   PetscBool      symbolic_used_the_fact_C_is_symmetric; /* MatMatMult(A,B,MAT_REUSE_MATRIX,..&C) is still legitimate), we need to redo symbolic! */
428:   PetscReal      fill;
429:   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 */

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

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

440: struct _p_Mat {
441:   PETSCHEADER(struct _MatOps);
442:   PetscLayout            rmap,cmap;
443:   void                   *data;            /* implementation-specific data */
444:   MatFactorType          factortype;       /* MAT_FACTOR_LU, ILU, CHOLESKY or ICC */
445:   PetscBool              trivialsymbolic;  /* indicates the symbolic factorization doesn't actually do a symbolic factorization, it is delayed to the numeric factorization */
446:   PetscBool              canuseordering;   /* factorization can use ordering provide to routine (most PETSc implementations) */
447:   MatOrderingType        preferredordering[MAT_FACTOR_NUM_TYPES] ;/* what is the preferred (or default) ordering for the matrix solver type */
448:   PetscBool              assembled;        /* is the matrix assembled? */
449:   PetscBool              was_assembled;    /* new values inserted into assembled mat */
450:   PetscInt               num_ass;          /* number of times matrix has been assembled */
451:   PetscObjectState       nonzerostate;     /* each time new nonzeros locations are introduced into the matrix this is updated */
452:   PetscObjectState       ass_nonzerostate; /* nonzero state at last assembly */
453:   MatInfo                info;             /* matrix information */
454:   InsertMode             insertmode;       /* have values been inserted in matrix or added? */
455:   MatStash               stash,bstash;     /* used for assembling off-proc mat emements */
456:   MatNullSpace           nullsp;           /* null space (operator is singular) */
457:   MatNullSpace           transnullsp;      /* null space of transpose of operator */
458:   MatNullSpace           nearnullsp;       /* near null space to be used by multigrid methods */
459:   PetscInt               congruentlayouts; /* are the rows and columns layouts congruent? */
460:   PetscBool              preallocated;
461:   MatStencilInfo         stencil;          /* information for structured grid */
462:   PetscBool              symmetric,hermitian,structurally_symmetric,spd;
463:   PetscBool              symmetric_set,hermitian_set,structurally_symmetric_set,spd_set; /* if true, then corresponding flag is correct*/
464:   PetscBool              symmetric_eternal;
465:   PetscBool              nooffprocentries,nooffproczerorows;
466:   PetscBool              assembly_subset;  /* set by MAT_SUBSET_OFF_PROC_ENTRIES */
467:   PetscBool              submat_singleis;  /* for efficient PCSetUp_ASM() */
468:   PetscBool              structure_only;
469:   PetscBool              sortedfull;       /* full, sorted rows are inserted */
470:   PetscBool              force_diagonals;  /* set by MAT_FORCE_DIAGONAL_ENTRIES */
471: #if defined(PETSC_HAVE_DEVICE)
472:   PetscOffloadMask       offloadmask;      /* a mask which indicates where the valid matrix data is (GPU, CPU or both) */
473:   PetscBool              boundtocpu;
474:   PetscBool              bindingpropagates;
475: #endif
476:   void                   *spptr;          /* pointer for special library like SuperLU */
477:   char                   *solvertype;
478:   PetscBool              checksymmetryonassembly,checknullspaceonassembly;
479:   PetscReal              checksymmetrytol;
480:   Mat                    schur;             /* Schur complement matrix */
481:   MatFactorSchurStatus   schur_status;      /* status of the Schur complement matrix */
482:   Mat_Redundant          *redundant;        /* used by MatCreateRedundantMatrix() */
483:   PetscBool              erroriffailure;    /* Generate an error if detected (for example a zero pivot) instead of returning */
484:   MatFactorError         factorerrortype;               /* type of error in factorization */
485:   PetscReal              factorerror_zeropivot_value;   /* If numerical zero pivot was detected this is the computed value */
486:   PetscInt               factorerror_zeropivot_row;     /* Row where zero pivot was detected */
487:   PetscInt               nblocks,*bsizes;   /* support for MatSetVariableBlockSizes() */
488:   char                   *defaultvectype;
489:   Mat_Product            *product;
490:   PetscBool              form_explicit_transpose; /* hint to generate an explicit mat tranpsose for operations like MatMultTranspose() */
491:   PetscBool              transupdated;            /* whether or not the explicitly generated transpose is up-to-date */
492: };

494: PETSC_INTERN PetscErrorCode MatAXPY_Basic(Mat,PetscScalar,Mat,MatStructure);
495: PETSC_INTERN PetscErrorCode MatAXPY_BasicWithPreallocation(Mat,Mat,PetscScalar,Mat,MatStructure);
496: PETSC_INTERN PetscErrorCode MatAXPY_Basic_Preallocate(Mat,Mat,Mat*);
497: PETSC_INTERN PetscErrorCode MatAXPY_Dense_Nest(Mat,PetscScalar,Mat);

499: /*
500:     Utility for MatFactor (Schur complement)
501: */
502: PETSC_INTERN PetscErrorCode MatFactorFactorizeSchurComplement_Private(Mat);
503: PETSC_INTERN PetscErrorCode MatFactorInvertSchurComplement_Private(Mat);
504: PETSC_INTERN PetscErrorCode MatFactorUpdateSchurStatus_Private(Mat);
505: PETSC_INTERN PetscErrorCode MatFactorSetUpInPlaceSchur_Private(Mat);

507: /*
508:     Utility for MatZeroRows
509: */
510: PETSC_INTERN PetscErrorCode MatZeroRowsMapLocal_Private(Mat,PetscInt,const PetscInt*,PetscInt*,PetscInt**);

512: /*
513:     Utility for MatView/MatLoad
514: */
515: PETSC_INTERN PetscErrorCode MatView_Binary_BlockSizes(Mat,PetscViewer);
516: PETSC_INTERN PetscErrorCode MatLoad_Binary_BlockSizes(Mat,PetscViewer);

518: /*
519:     Object for partitioning graphs
520: */

522: typedef struct _MatPartitioningOps *MatPartitioningOps;
523: struct _MatPartitioningOps {
524:   PetscErrorCode (*apply)(MatPartitioning,IS*);
525:   PetscErrorCode (*applynd)(MatPartitioning,IS*);
526:   PetscErrorCode (*setfromoptions)(PetscOptionItems*,MatPartitioning);
527:   PetscErrorCode (*destroy)(MatPartitioning);
528:   PetscErrorCode (*view)(MatPartitioning,PetscViewer);
529:   PetscErrorCode (*improve)(MatPartitioning,IS*);
530: };

532: struct _p_MatPartitioning {
533:   PETSCHEADER(struct _MatPartitioningOps);
534:   Mat         adj;
535:   PetscInt    *vertex_weights;
536:   PetscReal   *part_weights;
537:   PetscInt    n;                                 /* number of partitions */
538:   void        *data;
539:   PetscInt    setupcalled;
540:   PetscBool   use_edge_weights;  /* A flag indicates whether or not to use edge weights */
541: };

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

546: /*
547:     Object for coarsen graphs
548: */
549: typedef struct _MatCoarsenOps *MatCoarsenOps;
550: struct _MatCoarsenOps {
551:   PetscErrorCode (*apply)(MatCoarsen);
552:   PetscErrorCode (*setfromoptions)(PetscOptionItems*,MatCoarsen);
553:   PetscErrorCode (*destroy)(MatCoarsen);
554:   PetscErrorCode (*view)(MatCoarsen,PetscViewer);
555: };

557: struct _p_MatCoarsen {
558:   PETSCHEADER(struct _MatCoarsenOps);
559:   Mat              graph;
560:   void             *subctx;
561:   /* */
562:   PetscBool        strict_aggs;
563:   IS               perm;
564:   PetscCoarsenData *agg_lists;
565: };

567: /*
568:     MatFDColoring is used to compute Jacobian matrices efficiently
569:   via coloring. The data structure is explained below in an example.

571:    Color =   0    1     0    2   |   2      3       0
572:    ---------------------------------------------------
573:             00   01              |          05
574:             10   11              |   14     15               Processor  0
575:                        22    23  |          25
576:                        32    33  |
577:    ===================================================
578:                                  |   44     45     46
579:             50                   |          55               Processor 1
580:                                  |   64            66
581:    ---------------------------------------------------

583:     ncolors = 4;

585:     ncolumns      = {2,1,1,0}
586:     columns       = {{0,2},{1},{3},{}}
587:     nrows         = {4,2,3,3}
588:     rows          = {{0,1,2,3},{0,1},{1,2,3},{0,1,2}}
589:     vwscale       = {dx(0),dx(1),dx(2),dx(3)}               MPI Vec
590:     vscale        = {dx(0),dx(1),dx(2),dx(3),dx(4),dx(5)}   Seq Vec

592:     ncolumns      = {1,0,1,1}
593:     columns       = {{6},{},{4},{5}}
594:     nrows         = {3,0,2,2}
595:     rows          = {{0,1,2},{},{1,2},{1,2}}
596:     vwscale       = {dx(4),dx(5),dx(6)}              MPI Vec
597:     vscale        = {dx(0),dx(4),dx(5),dx(6)}        Seq Vec

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

602: */
603: typedef struct {
604:   PetscInt     row;
605:   PetscInt     col;
606:   PetscScalar  *valaddr;   /* address of value */
607: } MatEntry;

609: typedef struct {
610:   PetscInt     row;
611:   PetscScalar  *valaddr;   /* address of value */
612: } MatEntry2;

614: struct  _p_MatFDColoring{
615:   PETSCHEADER(int);
616:   PetscInt       M,N,m;            /* total rows, columns; local rows */
617:   PetscInt       rstart;           /* first row owned by local processor */
618:   PetscInt       ncolors;          /* number of colors */
619:   PetscInt       *ncolumns;        /* number of local columns for a color */
620:   PetscInt       **columns;        /* lists the local columns of each color (using global column numbering) */
621:   IS             *isa;             /* these are the IS that contain the column values given in columns */
622:   PetscInt       *nrows;           /* number of local rows for each color */
623:   MatEntry       *matentry;        /* holds (row, column, address of value) for Jacobian matrix entry */
624:   MatEntry2      *matentry2;       /* holds (row, address of value) for Jacobian matrix entry */
625:   PetscScalar    *dy;              /* store a block of F(x+dx)-F(x) when J is in BAIJ format */
626:   PetscReal      error_rel;        /* square root of relative error in computing function */
627:   PetscReal      umin;             /* minimum allowable u'dx value */
628:   Vec            w1,w2,w3;         /* work vectors used in computing Jacobian */
629:   PetscBool      fset;             /* indicates that the initial function value F(X) is set */
630:   PetscErrorCode (*f)(void);       /* function that defines Jacobian */
631:   void           *fctx;            /* optional user-defined context for use by the function f */
632:   Vec            vscale;           /* holds FD scaling, i.e. 1/dx for each perturbed column */
633:   PetscInt       currentcolor;     /* color for which function evaluation is being done now */
634:   const char     *htype;           /* "wp" or "ds" */
635:   ISColoringType ctype;            /* IS_COLORING_GLOBAL or IS_COLORING_LOCAL */
636:   PetscInt       brows,bcols;      /* number of block rows or columns for speedup inserting the dense matrix into sparse Jacobian */
637:   PetscBool      setupcalled;      /* true if setup has been called */
638:   PetscBool      viewed;           /* true if the -mat_fd_coloring_view has been triggered already */
639:   void           (*ftn_func_pointer)(void),*ftn_func_cntx; /* serve the same purpose as *fortran_func_pointers in PETSc objects */
640:   PetscObjectId  matid;            /* matrix this object was created with, must always be the same */
641: };

643: typedef struct _MatColoringOps *MatColoringOps;
644: struct _MatColoringOps {
645:   PetscErrorCode (*destroy)(MatColoring);
646:   PetscErrorCode (*setfromoptions)(PetscOptionItems*,MatColoring);
647:   PetscErrorCode (*view)(MatColoring,PetscViewer);
648:   PetscErrorCode (*apply)(MatColoring,ISColoring*);
649:   PetscErrorCode (*weights)(MatColoring,PetscReal**,PetscInt**);
650: };

652: struct _p_MatColoring {
653:   PETSCHEADER(struct _MatColoringOps);
654:   Mat                   mat;
655:   PetscInt              dist;             /* distance of the coloring */
656:   PetscInt              maxcolors;        /* the maximum number of colors returned, maxcolors=1 for MIS */
657:   void                  *data;            /* inner context */
658:   PetscBool             valid;            /* check to see if what is produced is a valid coloring */
659:   MatColoringWeightType weight_type;      /* type of weight computation to be performed */
660:   PetscReal             *user_weights;    /* custom weights and permutation */
661:   PetscInt              *user_lperm;
662:   PetscBool             valid_iscoloring; /* check to see if matcoloring is produced a valid iscoloring */
663: };

665: struct  _p_MatTransposeColoring{
666:   PETSCHEADER(int);
667:   PetscInt       M,N,m;            /* total rows, columns; local rows */
668:   PetscInt       rstart;           /* first row owned by local processor */
669:   PetscInt       ncolors;          /* number of colors */
670:   PetscInt       *ncolumns;        /* number of local columns for a color */
671:   PetscInt       *nrows;           /* number of local rows for each color */
672:   PetscInt       currentcolor;     /* color for which function evaluation is being done now */
673:   ISColoringType ctype;            /* IS_COLORING_GLOBAL or IS_COLORING_LOCAL */

675:   PetscInt       *colorforrow,*colorforcol;  /* pointer to rows and columns */
676:   PetscInt       *rows;                      /* lists the local rows for each color (using the local row numbering) */
677:   PetscInt       *den2sp;                    /* maps (row,color) in the dense matrix to index of sparse matrix array a->a */
678:   PetscInt       *columns;                   /* lists the local columns of each color (using global column numbering) */
679:   PetscInt       brows;                      /* number of rows for efficient implementation of MatTransColoringApplyDenToSp() */
680:   PetscInt       *lstart;                    /* array used for loop over row blocks of Csparse */
681: };

683: /*
684:    Null space context for preconditioner/operators
685: */
686: struct _p_MatNullSpace {
687:   PETSCHEADER(int);
688:   PetscBool      has_cnst;
689:   PetscInt       n;
690:   Vec*           vecs;
691:   PetscScalar*   alpha;                 /* for projections */
692:   PetscErrorCode (*remove)(MatNullSpace,Vec,void*);  /* for user provided removal function */
693:   void*          rmctx;                 /* context for remove() function */
694: };

696: /*
697:    Checking zero pivot for LU, ILU preconditioners.
698: */
699: typedef struct {
700:   PetscInt       nshift,nshift_max;
701:   PetscReal      shift_amount,shift_lo,shift_hi,shift_top,shift_fraction;
702:   PetscBool      newshift;
703:   PetscReal      rs;  /* active row sum of abs(offdiagonals) */
704:   PetscScalar    pv;  /* pivot of the active row */
705: } FactorShiftCtx;

707: /*
708:  Used by MatCreateSubMatrices_MPIXAIJ_Local()
709: */
710: #include <petscctable.h>
711: typedef struct { /* used by MatCreateSubMatrices_MPIAIJ_SingleIS_Local() and MatCreateSubMatrices_MPIAIJ_Local */
712:   PetscInt   id;   /* index of submats, only submats[0] is responsible for deleting some arrays below */
713:   PetscInt   nrqs,nrqr;
714:   PetscInt   **rbuf1,**rbuf2,**rbuf3,**sbuf1,**sbuf2;
715:   PetscInt   **ptr;
716:   PetscInt   *tmp;
717:   PetscInt   *ctr;
718:   PetscInt   *pa; /* proc array */
719:   PetscInt   *req_size,*req_source1,*req_source2;
720:   PetscBool  allcolumns,allrows;
721:   PetscBool  singleis;
722:   PetscInt   *row2proc; /* row to proc map */
723:   PetscInt   nstages;
724: #if defined(PETSC_USE_CTABLE)
725:   PetscTable cmap,rmap;
726:   PetscInt   *cmap_loc,*rmap_loc;
727: #else
728:   PetscInt   *cmap,*rmap;
729: #endif

731:   PetscErrorCode (*destroy)(Mat);
732: } Mat_SubSppt;

734: PETSC_EXTERN PetscErrorCode MatFactorDumpMatrix(Mat);
735: PETSC_INTERN PetscErrorCode MatShift_Basic(Mat,PetscScalar);
736: PETSC_INTERN PetscErrorCode MatSetBlockSizes_Default(Mat,PetscInt,PetscInt);

738: static inline PetscErrorCode MatPivotCheck_nz(Mat mat,const MatFactorInfo *info,FactorShiftCtx *sctx,PetscInt row)
739: {
740:   PetscReal _rs   = sctx->rs;
741:   PetscReal _zero = info->zeropivot*_rs;

743:   if (PetscAbsScalar(sctx->pv) <= _zero && !PetscIsNanScalar(sctx->pv)) {
744:     /* force |diag| > zeropivot*rs */
745:     if (!sctx->nshift) sctx->shift_amount = info->shiftamount;
746:     else sctx->shift_amount *= 2.0;
747:     sctx->newshift = PETSC_TRUE;
748:     (sctx->nshift)++;
749:   } else {
750:     sctx->newshift = PETSC_FALSE;
751:   }
752:   return 0;
753: }

755: static inline PetscErrorCode MatPivotCheck_pd(Mat mat,const MatFactorInfo *info,FactorShiftCtx *sctx,PetscInt row)
756: {
757:   PetscReal _rs   = sctx->rs;
758:   PetscReal _zero = info->zeropivot*_rs;

760:   if (PetscRealPart(sctx->pv) <= _zero && !PetscIsNanScalar(sctx->pv)) {
761:     /* force matfactor to be diagonally dominant */
762:     if (sctx->nshift == sctx->nshift_max) {
763:       sctx->shift_fraction = sctx->shift_hi;
764:     } else {
765:       sctx->shift_lo = sctx->shift_fraction;
766:       sctx->shift_fraction = (sctx->shift_hi+sctx->shift_lo)/2.;
767:     }
768:     sctx->shift_amount = sctx->shift_fraction * sctx->shift_top;
769:     sctx->nshift++;
770:     sctx->newshift = PETSC_TRUE;
771:   } else {
772:     sctx->newshift = PETSC_FALSE;
773:   }
774:   return 0;
775: }

777: static inline PetscErrorCode MatPivotCheck_inblocks(Mat mat,const MatFactorInfo *info,FactorShiftCtx *sctx,PetscInt row)
778: {
779:   PetscReal _zero = info->zeropivot;

781:   if (PetscAbsScalar(sctx->pv) <= _zero && !PetscIsNanScalar(sctx->pv)) {
782:     sctx->pv          += info->shiftamount;
783:     sctx->shift_amount = 0.0;
784:     sctx->nshift++;
785:   }
786:   sctx->newshift = PETSC_FALSE;
787:   return 0;
788: }

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

794:   sctx->newshift = PETSC_FALSE;
795:   if (PetscAbsScalar(sctx->pv) <= _zero && !PetscIsNanScalar(sctx->pv)) {
797:     PetscInfo(mat,"Detected zero pivot in factorization in row %" PetscInt_FMT " value %g tolerance %g\n",row,(double)PetscAbsScalar(sctx->pv),(double)_zero);
798:     fact->factorerrortype             = MAT_FACTOR_NUMERIC_ZEROPIVOT;
799:     fact->factorerror_zeropivot_value = PetscAbsScalar(sctx->pv);
800:     fact->factorerror_zeropivot_row   = row;
801:   }
802:   return 0;
803: }

805: static inline PetscErrorCode MatPivotCheck(Mat fact,Mat mat,const MatFactorInfo *info,FactorShiftCtx *sctx,PetscInt row)
806: {
807:   if (info->shifttype == (PetscReal) MAT_SHIFT_NONZERO) {
808:     MatPivotCheck_nz(mat,info,sctx,row);
809:   } else if (info->shifttype == (PetscReal) MAT_SHIFT_POSITIVE_DEFINITE) {
810:     MatPivotCheck_pd(mat,info,sctx,row);
811:   } else if (info->shifttype == (PetscReal) MAT_SHIFT_INBLOCKS) {
812:     MatPivotCheck_inblocks(mat,info,sctx,row);
813:   } else {
814:     MatPivotCheck_none(fact,mat,info,sctx,row);
815:   }
816:   return 0;
817: }

819: #include <petscbt.h>
820: /*
821:   Create and initialize a linked list
822:   Input Parameters:
823:     idx_start - starting index of the list
824:     lnk_max   - max value of lnk indicating the end of the list
825:     nlnk      - max length of the list
826:   Output Parameters:
827:     lnk       - list initialized
828:     bt        - PetscBT (bitarray) with all bits set to false
829:     lnk_empty - flg indicating the list is empty
830: */
831: #define PetscLLCreate(idx_start,lnk_max,nlnk,lnk,bt) \
832:   (PetscMalloc1(nlnk,&lnk) || PetscBTCreate(nlnk,&(bt)) || (lnk[idx_start] = lnk_max,0))

834: #define PetscLLCreate_new(idx_start,lnk_max,nlnk,lnk,bt,lnk_empty)\
835:   (PetscMalloc1(nlnk,&lnk) || PetscBTCreate(nlnk,&(bt)) || (lnk_empty = PETSC_TRUE,0) ||(lnk[idx_start] = lnk_max,0))

837: 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)
838: {
839:   PetscInt location;

841:   /* start from the beginning if entry < previous entry */
842:   if (!assume_sorted && k && entry < *lnkdata) *lnkdata = idx_start;
843:   /* search for insertion location */
844:   do {
845:     location = *lnkdata;
846:     *lnkdata = lnk[location];
847:   } while (entry > *lnkdata);
848:   /* insertion location is found, add entry into lnk */
849:   lnk[location] = entry;
850:   lnk[entry]    = *lnkdata;
851:   ++(*nlnk);
852:   *lnkdata      = entry; /* next search starts from here if next_entry > entry */
853:   return 0;
854: }

856: 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)
857: {
858:   *nlnk = 0;
859:   for (PetscInt k = 0, lnkdata = idx_start; k < nidx; ++k) {
860:     const PetscInt entry = indices[k];

862:     if (!PetscBTLookupSet(bt,entry)) PetscLLInsertLocation_Private(assume_sorted,k,idx_start,entry,nlnk,&lnkdata,lnk);
863:   }
864:   return 0;
865: }

867: /*
868:   Add an index set into a sorted linked list
869:   Input Parameters:
870:     nidx      - number of input indices
871:     indices   - integer array
872:     idx_start - starting index of the list
873:     lnk       - linked list(an integer array) that is created
874:     bt        - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
875:   output Parameters:
876:     nlnk      - number of newly added indices
877:     lnk       - the sorted(increasing order) linked list containing new and non-redundate entries from indices
878:     bt        - updated PetscBT (bitarray)
879: */
880: static inline PetscErrorCode PetscLLAdd(PetscInt nidx, const PetscInt *PETSC_RESTRICT indices, PetscInt idx_start, PetscInt *PETSC_RESTRICT nlnk, PetscInt *PETSC_RESTRICT lnk, PetscBT bt)
881: {
882:   PetscLLAdd_Private(nidx,indices,idx_start,nlnk,lnk,bt,PETSC_FALSE);
883:   return 0;
884: }

886: /*
887:   Add a SORTED ascending index set into a sorted linked list - same as PetscLLAdd() bus skip 'if (_k && _entry < _lnkdata) _lnkdata  = idx_start;'
888:   Input Parameters:
889:     nidx      - number of input indices
890:     indices   - sorted integer array
891:     idx_start - starting index of the list
892:     lnk       - linked list(an integer array) that is created
893:     bt        - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
894:   output Parameters:
895:     nlnk      - number of newly added indices
896:     lnk       - the sorted(increasing order) linked list containing new and non-redundate entries from indices
897:     bt        - updated PetscBT (bitarray)
898: */
899: static inline PetscErrorCode PetscLLAddSorted(PetscInt nidx, const PetscInt *PETSC_RESTRICT indices, PetscInt idx_start, PetscInt *PETSC_RESTRICT nlnk, PetscInt *PETSC_RESTRICT lnk, PetscBT bt)
900: {
901:   PetscLLAdd_Private(nidx,indices,idx_start,nlnk,lnk,bt,PETSC_TRUE);
902:   return 0;
903: }

905: /*
906:   Add a permuted index set into a sorted linked list
907:   Input Parameters:
908:     nidx      - number of input indices
909:     indices   - integer array
910:     perm      - permutation of indices
911:     idx_start - starting index of the list
912:     lnk       - linked list(an integer array) that is created
913:     bt        - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
914:   output Parameters:
915:     nlnk      - number of newly added indices
916:     lnk       - the sorted(increasing order) linked list containing new and non-redundate entries from indices
917:     bt        - updated PetscBT (bitarray)
918: */
919: static inline PetscErrorCode 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)
920: {
921:   *nlnk = 0;
922:   for (PetscInt k = 0, lnkdata = idx_start; k < nidx; ++k) {
923:     const PetscInt entry = perm[indices[k]];

925:     if (!PetscBTLookupSet(bt,entry)) PetscLLInsertLocation_Private(PETSC_FALSE,k,idx_start,entry,nlnk,&lnkdata,lnk);
926:   }
927:   return 0;
928: }

930: #if 0
931: /* this appears to be unused? */
932: static inline PetscErrorCode PetscLLAddSorted_new(PetscInt nidx, PetscInt *indices, PetscInt idx_start, PetscBool *lnk_empty, PetscInt *nlnk, PetscInt *lnk, PetscBT bt)
933: {
934:   PetscInt lnkdata = idx_start;

936:   if (*lnk_empty) {
937:     for (PetscInt k = 0; k < nidx; ++k) {
938:       const PetscInt entry = indices[k], location = lnkdata;

940:       PetscBTSet(bt,entry); /* mark the new entry */
941:       lnkdata       = lnk[location];
942:       /* insertion location is found, add entry into lnk */
943:       lnk[location] = entry;
944:       lnk[entry]    = lnkdata;
945:       lnkdata       = entry; /* next search starts from here */
946:     }
947:     /* lnk[indices[nidx-1]] = lnk[idx_start];
948:        lnk[idx_start]       = indices[0];
949:        PetscBTSet(bt,indices[0]);
950:        for (_k=1; _k<nidx; _k++) {
951:        PetscBTSet(bt,indices[_k]);
952:        lnk[indices[_k-1]] = indices[_k];
953:        }
954:     */
955:     *nlnk      = nidx;
956:     *lnk_empty = PETSC_FALSE;
957:   } else {
958:     *nlnk = 0;
959:     for (PetscInt k = 0; k < nidx; ++k) {
960:       const PetscInt entry = indices[k];

962:       if (!PetscBTLookupSet(bt,entry)) PetscLLInsertLocation_Private(PETSC_TRUE,k,idx_start,entry,nlnk,&lnkdata,lnk);
963:     }
964:   }
965:   return 0;
966: }
967: #endif

969: /*
970:   Add a SORTED index set into a sorted linked list used for LUFactorSymbolic()
971:   Same as PetscLLAddSorted() with an additional operation:
972:        count the number of input indices that are no larger than 'diag'
973:   Input Parameters:
974:     indices   - sorted integer array
975:     idx_start - starting index of the list, index of pivot row
976:     lnk       - linked list(an integer array) that is created
977:     bt        - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
978:     diag      - index of the active row in LUFactorSymbolic
979:     nzbd      - number of input indices with indices <= idx_start
980:     im        - im[idx_start] is initialized as num of nonzero entries in row=idx_start
981:   output Parameters:
982:     nlnk      - number of newly added indices
983:     lnk       - the sorted(increasing order) linked list containing new and non-redundate entries from indices
984:     bt        - updated PetscBT (bitarray)
985:     im        - im[idx_start]: unchanged if diag is not an entry
986:                              : num of entries with indices <= diag if diag is an entry
987: */
988: 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)
989: {
990:   const PetscInt nidx = im[idx_start]-nzbd; /* num of entries with idx_start < index <= diag */

992:   *nlnk = 0;
993:   for (PetscInt k = 0, lnkdata = idx_start; k < nidx; ++k) {
994:     const PetscInt entry = indices[k];

996:     ++nzbd;
997:     if (entry == diag) im[idx_start] = nzbd;
998:     if (!PetscBTLookupSet(bt,entry)) PetscLLInsertLocation_Private(PETSC_TRUE,k,idx_start,entry,nlnk,&lnkdata,lnk);
999:   }
1000:   return 0;
1001: }

1003: /*
1004:   Copy data on the list into an array, then initialize the list
1005:   Input Parameters:
1006:     idx_start - starting index of the list
1007:     lnk_max   - max value of lnk indicating the end of the list
1008:     nlnk      - number of data on the list to be copied
1009:     lnk       - linked list
1010:     bt        - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
1011:   output Parameters:
1012:     indices   - array that contains the copied data
1013:     lnk       - linked list that is cleaned and initialize
1014:     bt        - PetscBT (bitarray) with all bits set to false
1015: */
1016: static inline PetscErrorCode PetscLLClean(PetscInt idx_start, PetscInt lnk_max, PetscInt nlnk, PetscInt *PETSC_RESTRICT lnk, PetscInt *PETSC_RESTRICT indices, PetscBT bt)
1017: {
1018:   for (PetscInt j = 0, idx = idx_start; j < nlnk; ++j) {
1019:     idx        = lnk[idx];
1020:     indices[j] = idx;
1021:     PetscBTClear(bt,idx);
1022:   }
1023:   lnk[idx_start] = lnk_max;
1024:   return 0;
1025: }

1027: /*
1028:   Free memories used by the list
1029: */
1030: #define PetscLLDestroy(lnk,bt) (PetscFree(lnk) || PetscBTDestroy(&(bt)))

1032: /* Routines below are used for incomplete matrix factorization */
1033: /*
1034:   Create and initialize a linked list and its levels
1035:   Input Parameters:
1036:     idx_start - starting index of the list
1037:     lnk_max   - max value of lnk indicating the end of the list
1038:     nlnk      - max length of the list
1039:   Output Parameters:
1040:     lnk       - list initialized
1041:     lnk_lvl   - array of size nlnk for storing levels of lnk
1042:     bt        - PetscBT (bitarray) with all bits set to false
1043: */
1044: #define PetscIncompleteLLCreate(idx_start,lnk_max,nlnk,lnk,lnk_lvl,bt)\
1045:   (PetscIntMultError(2,nlnk,NULL) || PetscMalloc1(2*nlnk,&lnk) || PetscBTCreate(nlnk,&(bt)) || (lnk[idx_start] = lnk_max,lnk_lvl = lnk + nlnk,0))

1047: 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)
1048: {
1049:   PetscLLInsertLocation_Private(assume_sorted,k,idx_start,entry,nlnk,lnkdata,lnk);
1050:   lnklvl[entry] = newval;
1051:   return 0;
1052: }

1054: /*
1055:   Initialize a sorted linked list used for ILU and ICC
1056:   Input Parameters:
1057:     nidx      - number of input idx
1058:     idx       - integer array used for storing column indices
1059:     idx_start - starting index of the list
1060:     perm      - indices of an IS
1061:     lnk       - linked list(an integer array) that is created
1062:     lnklvl    - levels of lnk
1063:     bt        - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
1064:   output Parameters:
1065:     nlnk     - number of newly added idx
1066:     lnk      - the sorted(increasing order) linked list containing new and non-redundate entries from idx
1067:     lnklvl   - levels of lnk
1068:     bt       - updated PetscBT (bitarray)
1069: */
1070: 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)
1071: {
1072:   *nlnk = 0;
1073:   for (PetscInt k = 0, lnkdata = idx_start; k < nidx; ++k) {
1074:     const PetscInt entry = perm[idx[k]];

1076:     if (!PetscBTLookupSet(bt,entry)) PetscIncompleteLLInsertLocation_Private(PETSC_FALSE,k,idx_start,entry,nlnk,&lnkdata,lnk,lnklvl,0);
1077:   }
1078:   return 0;
1079: }

1081: 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)
1082: {
1083:   *nlnk = 0;
1084:   for (PetscInt k = 0, lnkdata = idx_start; k < nidx; ++k) {
1085:     const PetscInt incrlev = idxlvl[k]+prow_offset+1;

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

1090:       if (!PetscBTLookupSet(bt,entry)) PetscIncompleteLLInsertLocation_Private(assume_sorted,k,idx_start,entry,nlnk,&lnkdata,lnk,lnklvl,incrlev);
1091:       else if (lnklvl[entry] > incrlev) lnklvl[entry] = incrlev;  /* existing entry */
1092:     }
1093:   }
1094:   return 0;
1095: }

1097: /*
1098:   Add a SORTED index set into a sorted linked list for ICC
1099:   Input Parameters:
1100:     nidx      - number of input indices
1101:     idx       - sorted integer array used for storing column indices
1102:     level     - level of fill, e.g., ICC(level)
1103:     idxlvl    - level of idx
1104:     idx_start - starting index of the list
1105:     lnk       - linked list(an integer array) that is created
1106:     lnklvl    - levels of lnk
1107:     bt        - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
1108:     idxlvl_prow - idxlvl[prow], where prow is the row number of the idx
1109:   output Parameters:
1110:     nlnk   - number of newly added indices
1111:     lnk    - the sorted(increasing order) linked list containing new and non-redundate entries from idx
1112:     lnklvl - levels of lnk
1113:     bt     - updated PetscBT (bitarray)
1114:   Note: the level of U(i,j) is set as lvl(i,j) = min{ lvl(i,j), lvl(prow,i)+lvl(prow,j)+1)
1115:         where idx = non-zero columns of U(prow,prow+1:n-1), prow<i
1116: */
1117: 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)
1118: {
1119:   PetscIncompleteLLAdd_Private(nidx,idx,level,idxlvl,idx_start,nlnk,lnk,lnklvl,bt,idxlvl_prow,PETSC_TRUE);
1120:   return 0;
1121: }

1123: /*
1124:   Add a SORTED index set into a sorted linked list for ILU
1125:   Input Parameters:
1126:     nidx      - number of input indices
1127:     idx       - sorted integer array used for storing column indices
1128:     level     - level of fill, e.g., ICC(level)
1129:     idxlvl    - level of idx
1130:     idx_start - starting index of the list
1131:     lnk       - linked list(an integer array) that is created
1132:     lnklvl    - levels of lnk
1133:     bt        - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
1134:     prow      - the row number of idx
1135:   output Parameters:
1136:     nlnk     - number of newly added idx
1137:     lnk      - the sorted(increasing order) linked list containing new and non-redundate entries from idx
1138:     lnklvl   - levels of lnk
1139:     bt       - updated PetscBT (bitarray)

1141:   Note: the level of factor(i,j) is set as lvl(i,j) = min{ lvl(i,j), lvl(i,prow)+lvl(prow,j)+1)
1142:         where idx = non-zero columns of U(prow,prow+1:n-1), prow<i
1143: */
1144: 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)
1145: {
1146:   PetscIncompleteLLAdd_Private(nidx,idx,level,idxlvl,idx_start,nlnk,lnk,lnklvl,bt,lnklvl[prow],PETSC_TRUE);
1147:   return 0;
1148: }

1150: /*
1151:   Add a index set into a sorted linked list
1152:   Input Parameters:
1153:     nidx      - number of input idx
1154:     idx   - integer array used for storing column indices
1155:     level     - level of fill, e.g., ICC(level)
1156:     idxlvl - level of idx
1157:     idx_start - starting index of the list
1158:     lnk       - linked list(an integer array) that is created
1159:     lnklvl   - levels of lnk
1160:     bt        - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
1161:   output Parameters:
1162:     nlnk      - number of newly added idx
1163:     lnk       - the sorted(increasing order) linked list containing new and non-redundate entries from idx
1164:     lnklvl   - levels of lnk
1165:     bt        - updated PetscBT (bitarray)
1166: */
1167: 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)
1168: {
1169:   PetscIncompleteLLAdd_Private(nidx,idx,level,idxlvl,idx_start,nlnk,lnk,lnklvl,bt,0,PETSC_FALSE);
1170:   return 0;
1171: }

1173: /*
1174:   Add a SORTED index set into a sorted linked list
1175:   Input Parameters:
1176:     nidx      - number of input indices
1177:     idx   - sorted integer array used for storing column indices
1178:     level     - level of fill, e.g., ICC(level)
1179:     idxlvl - level of idx
1180:     idx_start - starting index of the list
1181:     lnk       - linked list(an integer array) that is created
1182:     lnklvl    - levels of lnk
1183:     bt        - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
1184:   output Parameters:
1185:     nlnk      - number of newly added idx
1186:     lnk       - the sorted(increasing order) linked list containing new and non-redundate entries from idx
1187:     lnklvl    - levels of lnk
1188:     bt        - updated PetscBT (bitarray)
1189: */
1190: 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)
1191: {
1192:   PetscIncompleteLLAdd_Private(nidx,idx,level,idxlvl,idx_start,nlnk,lnk,lnklvl,bt,0,PETSC_TRUE);
1193:   return 0;
1194: }

1196: /*
1197:   Copy data on the list into an array, then initialize the list
1198:   Input Parameters:
1199:     idx_start - starting index of the list
1200:     lnk_max   - max value of lnk indicating the end of the list
1201:     nlnk      - number of data on the list to be copied
1202:     lnk       - linked list
1203:     lnklvl    - level of lnk
1204:     bt        - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
1205:   output Parameters:
1206:     indices - array that contains the copied data
1207:     lnk     - linked list that is cleaned and initialize
1208:     lnklvl  - level of lnk that is reinitialized
1209:     bt      - PetscBT (bitarray) with all bits set to false
1210: */
1211: 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)
1212: {
1213:   for (PetscInt j = 0, idx = idx_start; j < nlnk; ++j) {
1214:     idx           = lnk[idx];
1215:     indices[j]    = idx;
1216:     indiceslvl[j] = lnklvl[idx];
1217:     lnklvl[idx]   = -1;
1218:     PetscBTClear(bt,idx);
1219:   }
1220:   lnk[idx_start] = lnk_max;
1221:   return 0;
1222: }

1224: /*
1225:   Free memories used by the list
1226: */
1227: #define PetscIncompleteLLDestroy(lnk,bt) (PetscFree(lnk) || PetscBTDestroy(&(bt)))

1229: #if !defined(PETSC_CLANG_STATIC_ANALYZER)
1230: #define MatCheckSameLocalSize(A,ar1,B,ar2) do {                         \
1233:   } while (0)
1234: #define MatCheckSameSize(A,ar1,B,ar2) do {                              \
1236:     MatCheckSameLocalSize(A,ar1,B,ar2);                                 \
1237:   } while (0)
1238: #else
1239: template <typename Tm>
1240: void MatCheckSameLocalSize(Tm,int,Tm,int);
1241: template <typename Tm>
1242: void MatCheckSameSize(Tm,int,Tm,int);
1243: #endif

1245: #define VecCheckMatCompatible(M,x,ar1,b,ar2) do { \
1248:   } while (0)

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

1258:       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
1259:       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
1260:       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
1261:       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.
1262:       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
1263:       to each other so memory access is much better than using the big array.

1265:   Example:
1266:      nlnk_max=5, lnk_max=36:
1267:      Initial list: [0, 0 | 36, 2 | 0, 0 | 0, 0 | 0, 0 | 0, 0 | 0, 0]
1268:      here, head_node has index 2 with value lnk[2]=lnk_max=36,
1269:            0-th entry is used to store the number of entries in the list,
1270:      The initial lnk represents head -> tail(marked by 36) with number of entries = lnk[0]=0.

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

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

1280:   Input Parameters:
1281:     nlnk_max  - max length of the list
1282:     lnk_max   - max value of the entries
1283:   Output Parameters:
1284:     lnk       - list created and initialized
1285:     bt        - PetscBT (bitarray) with all bits set to false. Note: bt has size lnk_max, not nln_max!
1286: */
1287: static inline PetscErrorCode PetscLLCondensedCreate(PetscInt nlnk_max,PetscInt lnk_max,PetscInt **lnk,PetscBT *bt)
1288: {
1289:   PetscInt *llnk,lsize = 0;

1291:   PetscIntMultError(2,nlnk_max+2,&lsize);
1292:   PetscMalloc1(lsize,lnk);
1293:   PetscBTCreate(lnk_max,bt);
1294:   llnk = *lnk;
1295:   llnk[0] = 0;         /* number of entries on the list */
1296:   llnk[2] = lnk_max;   /* value in the head node */
1297:   llnk[3] = 2;         /* next for the head node */
1298:   return 0;
1299: }

1301: /*
1302:   Add a SORTED ascending index set into a sorted linked list. See PetscLLCondensedCreate() for detailed description.
1303:   Input Parameters:
1304:     nidx      - number of input indices
1305:     indices   - sorted integer array
1306:     lnk       - condensed linked list(an integer array) that is created
1307:     bt        - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
1308:   output Parameters:
1309:     lnk       - the sorted(increasing order) linked list containing previous and newly added non-redundate indices
1310:     bt        - updated PetscBT (bitarray)
1311: */
1312: static inline PetscErrorCode PetscLLCondensedAddSorted(PetscInt nidx,const PetscInt indices[],PetscInt lnk[],PetscBT bt)
1313: {
1314:   PetscInt _k,_entry,_location,_next,_lnkdata,_nlnk,_newnode;

1316:   _nlnk     = lnk[0]; /* num of entries on the input lnk */
1317:   _location = 2; /* head */
1318:     for (_k=0; _k<nidx; _k++) {
1319:       _entry = indices[_k];
1320:       if (!PetscBTLookupSet(bt,_entry)) {  /* new entry */
1321:         /* search for insertion location */
1322:         do {
1323:           _next     = _location + 1; /* link from previous node to next node */
1324:           _location = lnk[_next];    /* idx of next node */
1325:           _lnkdata  = lnk[_location];/* value of next node */
1326:         } while (_entry > _lnkdata);
1327:         /* insertion location is found, add entry into lnk */
1328:         _newnode        = 2*(_nlnk+2);   /* index for this new node */
1329:         lnk[_next]      = _newnode;      /* connect previous node to the new node */
1330:         lnk[_newnode]   = _entry;        /* set value of the new node */
1331:         lnk[_newnode+1] = _location;     /* connect new node to next node */
1332:         _location       = _newnode;      /* next search starts from the new node */
1333:         _nlnk++;
1334:       }   \
1335:     }\
1336:   lnk[0]   = _nlnk;   /* number of entries in the list */
1337:   return 0;
1338: }

1340: static inline PetscErrorCode PetscLLCondensedClean(PetscInt lnk_max,PetscInt nidx,PetscInt *indices,PetscInt lnk[],PetscBT bt)
1341: {
1342:   PetscInt _next = lnk[3]; /* head node */
1343:   PetscInt _nlnk = lnk[0]; /* num of entries on the list */

1345:   for (PetscInt _k=0; _k<_nlnk; _k++) {
1346:     indices[_k] = lnk[_next];
1347:     _next       = lnk[_next + 1];
1348:     PetscBTClear(bt,indices[_k]);
1349:   }
1350:   lnk[0] = 0;          /* num of entries on the list */
1351:   lnk[2] = lnk_max;    /* initialize head node */
1352:   lnk[3] = 2;          /* head node */
1353:   return 0;
1354: }

1356: static inline PetscErrorCode PetscLLCondensedView(PetscInt *lnk)
1357: {
1358:   PetscPrintf(PETSC_COMM_SELF,"LLCondensed of size %" PetscInt_FMT ", (val,  next)\n",lnk[0]);
1359:   for (PetscInt k = 2; k < lnk[0]+2; ++k) {
1360:     PetscPrintf(PETSC_COMM_SELF," %" PetscInt_FMT ": (%" PetscInt_FMT ", %" PetscInt_FMT")\n",2*k,lnk[2*k],lnk[2*k+1]);
1361:   }
1362:   return 0;
1363: }

1365: /*
1366:   Free memories used by the list
1367: */
1368: static inline PetscErrorCode PetscLLCondensedDestroy(PetscInt *lnk,PetscBT bt)
1369: {
1370:   PetscFree(lnk);
1371:   PetscBTDestroy(&bt);
1372:   return 0;
1373: }

1375: /* -------------------------------------------------------------------------------------------------------*/
1376: /*
1377:  Same as PetscLLCondensedCreate(), but does not use non-scalable O(lnk_max) bitarray
1378:   Input Parameters:
1379:     nlnk_max  - max length of the list
1380:   Output Parameters:
1381:     lnk       - list created and initialized
1382: */
1383: static inline PetscErrorCode PetscLLCondensedCreate_Scalable(PetscInt nlnk_max,PetscInt **lnk)
1384: {
1385:   PetscInt *llnk,lsize = 0;

1387:   PetscIntMultError(2,nlnk_max+2,&lsize);
1388:   PetscMalloc1(lsize,lnk);
1389:   llnk = *lnk;
1390:   llnk[0] = 0;               /* number of entries on the list */
1391:   llnk[2] = PETSC_MAX_INT;   /* value in the head node */
1392:   llnk[3] = 2;               /* next for the head node */
1393:   return 0;
1394: }

1396: static inline PetscErrorCode PetscLLCondensedExpand_Scalable(PetscInt nlnk_max,PetscInt **lnk)
1397: {
1398:   PetscInt lsize = 0;

1400:   PetscIntMultError(2,nlnk_max+2,&lsize);
1401:   PetscRealloc(lsize*sizeof(PetscInt),lnk);
1402:   return 0;
1403: }

1405: static inline PetscErrorCode PetscLLCondensedAddSorted_Scalable(PetscInt nidx,const PetscInt indices[],PetscInt lnk[])
1406: {
1407:   PetscInt _k,_entry,_location,_next,_lnkdata,_nlnk,_newnode;
1408:   _nlnk     = lnk[0]; /* num of entries on the input lnk */
1409:   _location = 2; /* head */ \
1410:     for (_k=0; _k<nidx; _k++) {
1411:       _entry = indices[_k];
1412:       /* search for insertion location */
1413:       do {
1414:         _next     = _location + 1; /* link from previous node to next node */
1415:         _location = lnk[_next];    /* idx of next node */
1416:         _lnkdata  = lnk[_location];/* value of next node */
1417:       } while (_entry > _lnkdata);
1418:       if (_entry < _lnkdata) {
1419:         /* insertion location is found, add entry into lnk */
1420:         _newnode        = 2*(_nlnk+2);   /* index for this new node */
1421:         lnk[_next]      = _newnode;      /* connect previous node to the new node */
1422:         lnk[_newnode]   = _entry;        /* set value of the new node */
1423:         lnk[_newnode+1] = _location;     /* connect new node to next node */
1424:         _location       = _newnode;      /* next search starts from the new node */
1425:         _nlnk++;
1426:       }
1427:     }
1428:   lnk[0]   = _nlnk;   /* number of entries in the list */
1429:   return 0;
1430: }

1432: static inline PetscErrorCode PetscLLCondensedClean_Scalable(PetscInt nidx,PetscInt *indices,PetscInt *lnk)
1433: {
1434:   PetscInt _k,_next,_nlnk;
1435:   _next = lnk[3];       /* head node */
1436:   _nlnk = lnk[0];
1437:   for (_k=0; _k<_nlnk; _k++) {
1438:     indices[_k] = lnk[_next];
1439:     _next       = lnk[_next + 1];
1440:   }
1441:   lnk[0] = 0;          /* num of entries on the list */
1442:   lnk[3] = 2;          /* head node */
1443:   return 0;
1444: }

1446: static inline PetscErrorCode PetscLLCondensedDestroy_Scalable(PetscInt *lnk)
1447: {
1448:   return PetscFree(lnk);
1449: }

1451: /* -------------------------------------------------------------------------------------------------------*/
1452: /*
1453:       lnk[0]   number of links
1454:       lnk[1]   number of entries
1455:       lnk[3n]  value
1456:       lnk[3n+1] len
1457:       lnk[3n+2] link to next value

1459:       The next three are always the first link

1461:       lnk[3]    PETSC_MIN_INT+1
1462:       lnk[4]    1
1463:       lnk[5]    link to first real entry

1465:       The next three are always the last link

1467:       lnk[6]    PETSC_MAX_INT - 1
1468:       lnk[7]    1
1469:       lnk[8]    next valid link (this is the same as lnk[0] but without the decreases)
1470: */

1472: static inline PetscErrorCode PetscLLCondensedCreate_fast(PetscInt nlnk_max,PetscInt **lnk)
1473: {
1474:   PetscInt *llnk,lsize = 0;

1476:   PetscIntMultError(3,nlnk_max+3,&lsize);
1477:   PetscMalloc1(lsize,lnk);
1478:   llnk = *lnk;
1479:   llnk[0] = 0;   /* nlnk: number of entries on the list */
1480:   llnk[1] = 0;          /* number of integer entries represented in list */
1481:   llnk[3] = PETSC_MIN_INT+1;   /* value in the first node */
1482:   llnk[4] = 1;           /* count for the first node */
1483:   llnk[5] = 6;         /* next for the first node */
1484:   llnk[6] = PETSC_MAX_INT-1;   /* value in the last node */
1485:   llnk[7] = 1;           /* count for the last node */
1486:   llnk[8] = 0;         /* next valid node to be used */
1487:   return 0;
1488: }

1490: static inline PetscErrorCode PetscLLCondensedAddSorted_fast(PetscInt nidx,const PetscInt indices[],PetscInt lnk[])
1491: {
1492:   PetscInt k,entry,prev,next;
1493:   prev      = 3;      /* first value */
1494:   next      = lnk[prev+2];
1495:   for (k=0; k<nidx; k++) {
1496:     entry = indices[k];
1497:     /* search for insertion location */
1498:     while (entry >= lnk[next]) {
1499:       prev = next;
1500:       next = lnk[next+2];
1501:     }
1502:     /* entry is in range of previous list */
1503:     if (entry < lnk[prev]+lnk[prev+1]) continue;
1504:     lnk[1]++;
1505:     /* entry is right after previous list */
1506:     if (entry == lnk[prev]+lnk[prev+1]) {
1507:       lnk[prev+1]++;
1508:       if (lnk[next] == entry+1) { /* combine two contiguous strings */
1509:         lnk[prev+1] += lnk[next+1];
1510:         lnk[prev+2]  = lnk[next+2];
1511:         next         = lnk[next+2];
1512:         lnk[0]--;
1513:       }
1514:       continue;
1515:     }
1516:     /* entry is right before next list */
1517:     if (entry == lnk[next]-1) {
1518:       lnk[next]--;
1519:       lnk[next+1]++;
1520:       prev = next;
1521:       next = lnk[prev+2];
1522:       continue;
1523:     }
1524:     /*  add entry into lnk */
1525:     lnk[prev+2]    = 3*((lnk[8]++)+3);      /* connect previous node to the new node */
1526:     prev           = lnk[prev+2];
1527:     lnk[prev]      = entry;        /* set value of the new node */
1528:     lnk[prev+1]    = 1;             /* number of values in contiguous string is one to start */
1529:     lnk[prev+2]    = next;          /* connect new node to next node */
1530:     lnk[0]++;
1531:   }
1532:   return 0;
1533: }

1535: static inline PetscErrorCode PetscLLCondensedClean_fast(PetscInt nidx,PetscInt *indices,PetscInt *lnk)
1536: {
1537:   PetscInt _k,_next,_nlnk,cnt,j;
1538:   _next = lnk[5];       /* first node */
1539:   _nlnk = lnk[0];
1540:   cnt   = 0;
1541:   for (_k=0; _k<_nlnk; _k++) {
1542:     for (j=0; j<lnk[_next+1]; j++) {
1543:       indices[cnt++] = lnk[_next] + j;
1544:     }
1545:     _next       = lnk[_next + 2];
1546:   }
1547:   lnk[0] = 0;   /* nlnk: number of links */
1548:   lnk[1] = 0;          /* number of integer entries represented in list */
1549:   lnk[3] = PETSC_MIN_INT+1;   /* value in the first node */
1550:   lnk[4] = 1;           /* count for the first node */
1551:   lnk[5] = 6;         /* next for the first node */
1552:   lnk[6] = PETSC_MAX_INT-1;   /* value in the last node */
1553:   lnk[7] = 1;           /* count for the last node */
1554:   lnk[8] = 0;         /* next valid location to make link */
1555:   return 0;
1556: }

1558: static inline PetscErrorCode PetscLLCondensedView_fast(PetscInt *lnk)
1559: {
1560:   PetscInt k,next,nlnk;
1561:   next = lnk[5];       /* first node */
1562:   nlnk = lnk[0];
1563:   for (k=0; k<nlnk; k++) {
1564: #if 0                           /* Debugging code */
1565:     printf("%d value %d len %d next %d\n",next,lnk[next],lnk[next+1],lnk[next+2]);
1566: #endif
1567:     next = lnk[next + 2];
1568:   }
1569:   return 0;
1570: }

1572: static inline PetscErrorCode PetscLLCondensedDestroy_fast(PetscInt *lnk)
1573: {
1574:   return PetscFree(lnk);
1575: }

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

1580: #if defined(PETSC_HAVE_KOKKOS_KERNELS)
1581: PETSC_INTERN PetscErrorCode MatSeqAIJMoveDiagonalValuesFront_SeqAIJKokkos(Mat,const PetscInt*);
1582: #endif

1584: PETSC_EXTERN PetscLogEvent MAT_Mult;
1585: PETSC_EXTERN PetscLogEvent MAT_MultMatrixFree;
1586: PETSC_EXTERN PetscLogEvent MAT_Mults;
1587: PETSC_EXTERN PetscLogEvent MAT_MultAdd;
1588: PETSC_EXTERN PetscLogEvent MAT_MultTranspose;
1589: PETSC_EXTERN PetscLogEvent MAT_MultTransposeAdd;
1590: PETSC_EXTERN PetscLogEvent MAT_Solve;
1591: PETSC_EXTERN PetscLogEvent MAT_Solves;
1592: PETSC_EXTERN PetscLogEvent MAT_SolveAdd;
1593: PETSC_EXTERN PetscLogEvent MAT_SolveTranspose;
1594: PETSC_EXTERN PetscLogEvent MAT_SolveTransposeAdd;
1595: PETSC_EXTERN PetscLogEvent MAT_SOR;
1596: PETSC_EXTERN PetscLogEvent MAT_ForwardSolve;
1597: PETSC_EXTERN PetscLogEvent MAT_BackwardSolve;
1598: PETSC_EXTERN PetscLogEvent MAT_LUFactor;
1599: PETSC_EXTERN PetscLogEvent MAT_LUFactorSymbolic;
1600: PETSC_EXTERN PetscLogEvent MAT_LUFactorNumeric;
1601: PETSC_EXTERN PetscLogEvent MAT_QRFactor;
1602: PETSC_EXTERN PetscLogEvent MAT_QRFactorSymbolic;
1603: PETSC_EXTERN PetscLogEvent MAT_QRFactorNumeric;
1604: PETSC_EXTERN PetscLogEvent MAT_CholeskyFactor;
1605: PETSC_EXTERN PetscLogEvent MAT_CholeskyFactorSymbolic;
1606: PETSC_EXTERN PetscLogEvent MAT_CholeskyFactorNumeric;
1607: PETSC_EXTERN PetscLogEvent MAT_ILUFactor;
1608: PETSC_EXTERN PetscLogEvent MAT_ILUFactorSymbolic;
1609: PETSC_EXTERN PetscLogEvent MAT_ICCFactorSymbolic;
1610: PETSC_EXTERN PetscLogEvent MAT_Copy;
1611: PETSC_EXTERN PetscLogEvent MAT_Convert;
1612: PETSC_EXTERN PetscLogEvent MAT_Scale;
1613: PETSC_EXTERN PetscLogEvent MAT_AssemblyBegin;
1614: PETSC_EXTERN PetscLogEvent MAT_AssemblyEnd;
1615: PETSC_EXTERN PetscLogEvent MAT_SetValues;
1616: PETSC_EXTERN PetscLogEvent MAT_GetValues;
1617: PETSC_EXTERN PetscLogEvent MAT_GetRow;
1618: PETSC_EXTERN PetscLogEvent MAT_GetRowIJ;
1619: PETSC_EXTERN PetscLogEvent MAT_CreateSubMats;
1620: PETSC_EXTERN PetscLogEvent MAT_GetColoring;
1621: PETSC_EXTERN PetscLogEvent MAT_GetOrdering;
1622: PETSC_EXTERN PetscLogEvent MAT_RedundantMat;
1623: PETSC_EXTERN PetscLogEvent MAT_IncreaseOverlap;
1624: PETSC_EXTERN PetscLogEvent MAT_Partitioning;
1625: PETSC_EXTERN PetscLogEvent MAT_PartitioningND;
1626: PETSC_EXTERN PetscLogEvent MAT_Coarsen;
1627: PETSC_EXTERN PetscLogEvent MAT_ZeroEntries;
1628: PETSC_EXTERN PetscLogEvent MAT_Load;
1629: PETSC_EXTERN PetscLogEvent MAT_View;
1630: PETSC_EXTERN PetscLogEvent MAT_AXPY;
1631: PETSC_EXTERN PetscLogEvent MAT_FDColoringCreate;
1632: PETSC_EXTERN PetscLogEvent MAT_TransposeColoringCreate;
1633: PETSC_EXTERN PetscLogEvent MAT_FDColoringSetUp;
1634: PETSC_EXTERN PetscLogEvent MAT_FDColoringApply;
1635: PETSC_EXTERN PetscLogEvent MAT_Transpose;
1636: PETSC_EXTERN PetscLogEvent MAT_FDColoringFunction;
1637: PETSC_EXTERN PetscLogEvent MAT_CreateSubMat;
1638: PETSC_EXTERN PetscLogEvent MAT_MatSolve;
1639: PETSC_EXTERN PetscLogEvent MAT_MatTrSolve;
1640: PETSC_EXTERN PetscLogEvent MAT_MatMultSymbolic;
1641: PETSC_EXTERN PetscLogEvent MAT_MatMultNumeric;
1642: PETSC_EXTERN PetscLogEvent MAT_Getlocalmatcondensed;
1643: PETSC_EXTERN PetscLogEvent MAT_GetBrowsOfAcols;
1644: PETSC_EXTERN PetscLogEvent MAT_GetBrowsOfAocols;
1645: PETSC_EXTERN PetscLogEvent MAT_PtAPSymbolic;
1646: PETSC_EXTERN PetscLogEvent MAT_PtAPNumeric;
1647: PETSC_EXTERN PetscLogEvent MAT_Seqstompinum;
1648: PETSC_EXTERN PetscLogEvent MAT_Seqstompisym;
1649: PETSC_EXTERN PetscLogEvent MAT_Seqstompi;
1650: PETSC_EXTERN PetscLogEvent MAT_Getlocalmat;
1651: PETSC_EXTERN PetscLogEvent MAT_RARtSymbolic;
1652: PETSC_EXTERN PetscLogEvent MAT_RARtNumeric;
1653: PETSC_EXTERN PetscLogEvent MAT_MatTransposeMultSymbolic;
1654: PETSC_EXTERN PetscLogEvent MAT_MatTransposeMultNumeric;
1655: PETSC_EXTERN PetscLogEvent MAT_TransposeMatMultSymbolic;
1656: PETSC_EXTERN PetscLogEvent MAT_TransposeMatMultNumeric;
1657: PETSC_EXTERN PetscLogEvent MAT_MatMatMultSymbolic;
1658: PETSC_EXTERN PetscLogEvent MAT_MatMatMultNumeric;
1659: PETSC_EXTERN PetscLogEvent MAT_Applypapt;
1660: PETSC_EXTERN PetscLogEvent MAT_Applypapt_symbolic;
1661: PETSC_EXTERN PetscLogEvent MAT_Applypapt_numeric;
1662: PETSC_EXTERN PetscLogEvent MAT_Getsymtranspose;
1663: PETSC_EXTERN PetscLogEvent MAT_Getsymtransreduced;
1664: PETSC_EXTERN PetscLogEvent MAT_GetSequentialNonzeroStructure;
1665: PETSC_EXTERN PetscLogEvent MATMFFD_Mult;
1666: PETSC_EXTERN PetscLogEvent MAT_GetMultiProcBlock;
1667: PETSC_EXTERN PetscLogEvent MAT_CUSPARSECopyToGPU;
1668: PETSC_EXTERN PetscLogEvent MAT_CUSPARSECopyFromGPU;
1669: PETSC_EXTERN PetscLogEvent MAT_CUSPARSEGenerateTranspose;
1670: PETSC_EXTERN PetscLogEvent MAT_CUSPARSESolveAnalysis;
1671: PETSC_EXTERN PetscLogEvent MAT_SetValuesBatch;
1672: PETSC_EXTERN PetscLogEvent MAT_ViennaCLCopyToGPU;
1673: PETSC_EXTERN PetscLogEvent MAT_DenseCopyToGPU;
1674: PETSC_EXTERN PetscLogEvent MAT_DenseCopyFromGPU;
1675: PETSC_EXTERN PetscLogEvent MAT_Merge;
1676: PETSC_EXTERN PetscLogEvent MAT_Residual;
1677: PETSC_EXTERN PetscLogEvent MAT_SetRandom;
1678: PETSC_EXTERN PetscLogEvent MAT_FactorFactS;
1679: PETSC_EXTERN PetscLogEvent MAT_FactorInvS;
1680: PETSC_EXTERN PetscLogEvent MAT_PreallCOO;
1681: PETSC_EXTERN PetscLogEvent MAT_SetVCOO;
1682: PETSC_EXTERN PetscLogEvent MATCOLORING_Apply;
1683: PETSC_EXTERN PetscLogEvent MATCOLORING_Comm;
1684: PETSC_EXTERN PetscLogEvent MATCOLORING_Local;
1685: PETSC_EXTERN PetscLogEvent MATCOLORING_ISCreate;
1686: PETSC_EXTERN PetscLogEvent MATCOLORING_SetUp;
1687: PETSC_EXTERN PetscLogEvent MATCOLORING_Weights;
1688: PETSC_EXTERN PetscLogEvent MAT_H2Opus_Build;
1689: PETSC_EXTERN PetscLogEvent MAT_H2Opus_Compress;
1690: PETSC_EXTERN PetscLogEvent MAT_H2Opus_Orthog;
1691: PETSC_EXTERN PetscLogEvent MAT_H2Opus_LR;

1693: #endif