Actual source code: matimpl.h


  2: #ifndef __MATIMPL_H

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

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

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

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

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

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

223: #include <petscsys.h>
224: PETSC_EXTERN PetscErrorCode MatRegisterOp(MPI_Comm, const char[], PetscVoidFunction, const char[], PetscInt, ...);
225: PETSC_EXTERN PetscErrorCode MatQueryOp(MPI_Comm, PetscVoidFunction*, const char[], PetscInt, ...);

227: typedef struct _p_MatRootName* MatRootName;
228: struct _p_MatRootName {
229:   char        *rname,*sname,*mname;
230:   MatRootName next;
231: };

233: PETSC_EXTERN MatRootName MatRootNameList;

235: /*
236:    Utility private matrix routines
237: */
238: PETSC_INTERN PetscErrorCode MatFindNonzeroRowsOrCols_Basic(Mat,PetscBool,PetscReal,IS*);
239: PETSC_INTERN PetscErrorCode MatConvert_Basic(Mat,MatType,MatReuse,Mat*);
240: PETSC_INTERN PetscErrorCode MatConvert_Shell(Mat,MatType,MatReuse,Mat*);
241: PETSC_INTERN PetscErrorCode MatConvertFrom_Shell(Mat,MatType,MatReuse,Mat*);
242: PETSC_INTERN PetscErrorCode MatCopy_Basic(Mat,Mat,MatStructure);
243: PETSC_INTERN PetscErrorCode MatDiagonalSet_Default(Mat,Vec,InsertMode);
244: #if defined(PETSC_HAVE_SCALAPACK)
245: PETSC_INTERN PetscErrorCode MatConvert_Dense_ScaLAPACK(Mat,MatType,MatReuse,Mat*);
246: #endif
247: PETSC_INTERN PetscErrorCode MatSetPreallocationCOO_Basic(Mat,PetscCount,const PetscInt[],const PetscInt[]);
248: PETSC_INTERN PetscErrorCode MatSetValuesCOO_Basic(Mat,const PetscScalar[],InsertMode);

250: /* these callbacks rely on the old matrix function pointers for
251:    matmat operations. They are unsafe, and should be removed.
252:    However, the amount of work needed to clean up all the
253:    implementations is not negligible */
254: PETSC_INTERN PetscErrorCode MatProductSymbolic_AB(Mat);
255: PETSC_INTERN PetscErrorCode MatProductNumeric_AB(Mat);
256: PETSC_INTERN PetscErrorCode MatProductSymbolic_AtB(Mat);
257: PETSC_INTERN PetscErrorCode MatProductNumeric_AtB(Mat);
258: PETSC_INTERN PetscErrorCode MatProductSymbolic_ABt(Mat);
259: PETSC_INTERN PetscErrorCode MatProductNumeric_ABt(Mat);
260: PETSC_INTERN PetscErrorCode MatProductNumeric_PtAP(Mat);
261: PETSC_INTERN PetscErrorCode MatProductNumeric_RARt(Mat);
262: PETSC_INTERN PetscErrorCode MatProductSymbolic_ABC(Mat);
263: PETSC_INTERN PetscErrorCode MatProductNumeric_ABC(Mat);

265: PETSC_INTERN PetscErrorCode MatProductCreate_Private(Mat,Mat,Mat,Mat);
266: /* this callback handles all the different triple products and
267:    does not rely on the function pointers; used by cuSPARSE and KOKKOS-KERNELS */
268: PETSC_INTERN PetscErrorCode MatProductSymbolic_ABC_Basic(Mat);

270: #if defined(PETSC_CLANG_STATIC_ANALYZER)
271: template <typename Tm> void MatCheckPreallocated(Tm,int);
272: template <typename Tm> void MatCheckProduct(Tm,int);
273: #else/* PETSC_CLANG_STATIC_ANALYZER */
274: #if defined(PETSC_USE_DEBUG)
275: #  define MatCheckPreallocated(A,arg) do {                              \
277:   } while (0)
278: #else
279: #  define MatCheckPreallocated(A,arg) do {} while (0)
280: #endif

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

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

297: typedef struct _MatStashSpace *PetscMatStashSpace;

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

308: PETSC_EXTERN PetscErrorCode PetscMatStashSpaceGet(PetscInt,PetscInt,PetscMatStashSpace *);
309: PETSC_EXTERN PetscErrorCode PetscMatStashSpaceContiguous(PetscInt,PetscMatStashSpace *,PetscScalar *,PetscInt *,PetscInt *);
310: PETSC_EXTERN PetscErrorCode PetscMatStashSpaceDestroy(PetscMatStashSpace*);

312: typedef struct {
313:   PetscInt    count;
314: } MatStashHeader;

316: typedef struct {
317:   void        *buffer;          /* Of type blocktype, dynamically constructed  */
318:   PetscInt    count;
319:   char        pending;
320: } MatStashFrame;

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

332:   PetscErrorCode (*ScatterBegin)(Mat,MatStash*,PetscInt*);
333:   PetscErrorCode (*ScatterGetMesg)(MatStash*,PetscMPIInt*,PetscInt**,PetscInt**,PetscScalar**,PetscInt*);
334:   PetscErrorCode (*ScatterEnd)(MatStash*);
335:   PetscErrorCode (*ScatterDestroy)(MatStash*);

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

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

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

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

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

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

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

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

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

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

497: PETSC_INTERN PetscErrorCode MatAXPY_Basic(Mat,PetscScalar,Mat,MatStructure);
498: PETSC_INTERN PetscErrorCode MatAXPY_BasicWithPreallocation(Mat,Mat,PetscScalar,Mat,MatStructure);
499: PETSC_INTERN PetscErrorCode MatAXPY_Basic_Preallocate(Mat,Mat,Mat*);
500: PETSC_INTERN PetscErrorCode MatAXPY_Dense_Nest(Mat,PetscScalar,Mat);

502: /*
503:     Utility for MatFactor (Schur complement)
504: */
505: PETSC_INTERN PetscErrorCode MatFactorFactorizeSchurComplement_Private(Mat);
506: PETSC_INTERN PetscErrorCode MatFactorInvertSchurComplement_Private(Mat);
507: PETSC_INTERN PetscErrorCode MatFactorUpdateSchurStatus_Private(Mat);
508: PETSC_INTERN PetscErrorCode MatFactorSetUpInPlaceSchur_Private(Mat);

510: /*
511:     Utility for MatZeroRows
512: */
513: PETSC_INTERN PetscErrorCode MatZeroRowsMapLocal_Private(Mat,PetscInt,const PetscInt*,PetscInt*,PetscInt**);

515: /*
516:     Utility for MatView/MatLoad
517: */
518: PETSC_INTERN PetscErrorCode MatView_Binary_BlockSizes(Mat,PetscViewer);
519: PETSC_INTERN PetscErrorCode MatLoad_Binary_BlockSizes(Mat,PetscViewer);

521: /*
522:     Object for partitioning graphs
523: */

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

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

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

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

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

570: /*
571:     MatFDColoring is used to compute Jacobian matrices efficiently
572:   via coloring. The data structure is explained below in an example.

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

586:     ncolors = 4;

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

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

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

605: */
606: typedef struct {
607:   PetscInt     row;
608:   PetscInt     col;
609:   PetscScalar  *valaddr;   /* address of value */
610: } MatEntry;

612: typedef struct {
613:   PetscInt     row;
614:   PetscScalar  *valaddr;   /* address of value */
615: } MatEntry2;

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

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

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

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

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

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

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

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

734:   PetscErrorCode (*destroy)(Mat);
735: } Mat_SubSppt;

737: PETSC_EXTERN PetscErrorCode MatFactorDumpMatrix(Mat);
738: PETSC_INTERN PetscErrorCode MatShift_Basic(Mat,PetscScalar);
739: PETSC_INTERN PetscErrorCode MatSetBlockSizes_Default(Mat,PetscInt,PetscInt);

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

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

758: static inline PetscErrorCode MatPivotCheck_pd(Mat mat,const MatFactorInfo *info,FactorShiftCtx *sctx,PetscInt row)
759: {
760:   PetscReal _rs   = sctx->rs;
761:   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: static inline PetscErrorCode MatPivotCheck_inblocks(Mat mat,const MatFactorInfo *info,FactorShiftCtx *sctx,PetscInt row)
781: {
782:   PetscReal _zero = info->zeropivot;

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

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

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

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

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

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

840: 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)
841: {
842:   PetscInt location;

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

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

865:     if (!PetscBTLookupSet(bt,entry)) PetscLLInsertLocation_Private(assume_sorted,k,idx_start,entry,nlnk,&lnkdata,lnk);
866:   }
867:   return 0;
868: }

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

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

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

928:     if (!PetscBTLookupSet(bt,entry)) PetscLLInsertLocation_Private(PETSC_FALSE,k,idx_start,entry,nlnk,&lnkdata,lnk);
929:   }
930:   return 0;
931: }

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

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

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

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

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

995:   *nlnk = 0;
996:   for (PetscInt k = 0, lnkdata = idx_start; k < nidx; ++k) {
997:     const PetscInt entry = indices[k];

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

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

1030: /*
1031:   Free memories used by the list
1032: */
1033: #define PetscLLDestroy(lnk,bt) (PetscFree(lnk) || PetscBTDestroy(&(bt)))

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

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

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

1079:     if (!PetscBTLookupSet(bt,entry)) PetscIncompleteLLInsertLocation_Private(PETSC_FALSE,k,idx_start,entry,nlnk,&lnkdata,lnk,lnklvl,0);
1080:   }
1081:   return 0;
1082: }

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

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

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

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

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

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

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

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

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

1227: /*
1228:   Free memories used by the list
1229: */
1230: #define PetscIncompleteLLDestroy(lnk,bt) (PetscFree(lnk) || PetscBTDestroy(&(bt)))

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

1399: static inline PetscErrorCode PetscLLCondensedExpand_Scalable(PetscInt nlnk_max,PetscInt **lnk)
1400: {
1401:   PetscInt lsize = 0;

1403:   PetscIntMultError(2,nlnk_max+2,&lsize);
1404:   PetscRealloc(lsize*sizeof(PetscInt),lnk);
1405:   return 0;
1406: }

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

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

1449: static inline PetscErrorCode PetscLLCondensedDestroy_Scalable(PetscInt *lnk)
1450: {
1451:   return PetscFree(lnk);
1452: }

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

1462:       The next three are always the first link

1464:       lnk[3]    PETSC_MIN_INT+1
1465:       lnk[4]    1
1466:       lnk[5]    link to first real entry

1468:       The next three are always the last link

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

1475: static inline PetscErrorCode PetscLLCondensedCreate_fast(PetscInt nlnk_max,PetscInt **lnk)
1476: {
1477:   PetscInt *llnk,lsize = 0;

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

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

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

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

1575: static inline PetscErrorCode PetscLLCondensedDestroy_fast(PetscInt *lnk)
1576: {
1577:   return PetscFree(lnk);
1578: }

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

1583: #if defined(PETSC_HAVE_KOKKOS_KERNELS)
1584: PETSC_INTERN PetscErrorCode MatSeqAIJMoveDiagonalValuesFront_SeqAIJKokkos(Mat,const PetscInt*);
1585: #endif

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

1696: #endif