Actual source code: matimpl.h

petsc-main 2021-04-20
Report Typos and Errors

  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: /*
 23:   This file defines the parts of the matrix data structure that are
 24:   shared by all matrix types.
 25: */

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

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

221: typedef struct _p_MatRootName* MatRootName;
222: struct _p_MatRootName {
223:   char        *rname,*sname,*mname;
224:   MatRootName next;
225: };

227: PETSC_EXTERN MatRootName MatRootNameList;

229: /*
230:    Utility private matrix routines
231: */
232: PETSC_INTERN PetscErrorCode MatFindNonzeroRowsOrCols_Basic(Mat,PetscBool,PetscReal,IS*);
233: PETSC_INTERN PetscErrorCode MatConvert_Basic(Mat,MatType,MatReuse,Mat*);
234: PETSC_INTERN PetscErrorCode MatConvert_Shell(Mat,MatType,MatReuse,Mat*);
235: PETSC_INTERN PetscErrorCode MatConvertFrom_Shell(Mat,MatType,MatReuse,Mat*);
236: PETSC_INTERN PetscErrorCode MatCopy_Basic(Mat,Mat,MatStructure);
237: PETSC_INTERN PetscErrorCode MatDiagonalSet_Default(Mat,Vec,InsertMode);
238: #if defined(PETSC_HAVE_SCALAPACK)
239: PETSC_INTERN PetscErrorCode MatConvert_Dense_ScaLAPACK(Mat,MatType,MatReuse,Mat*);
240: #endif
241: PETSC_INTERN PetscErrorCode MatSetPreallocationCOO_Basic(Mat,PetscInt,const PetscInt[],const PetscInt[]);
242: PETSC_INTERN PetscErrorCode MatSetValuesCOO_Basic(Mat,const PetscScalar[],InsertMode);

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

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


265: #if defined(PETSC_USE_DEBUG)
266: #  define MatCheckPreallocated(A,arg) do {                              \
267:     if (PetscUnlikely(!(A)->preallocated)) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call MatXXXSetPreallocation(), MatSetUp() or the matrix has not yet been factored on argument %D \"%s\" before %s()",(arg),#A,PETSC_FUNCTION_NAME); \
268:   } while (0)
269: #else
270: #  define MatCheckPreallocated(A,arg) do {} while (0)
271: #endif

273: #if defined(PETSC_USE_DEBUG)
274: #  define MatCheckProduct(A,arg) do {                              \
275:     if (PetscUnlikely(!(A)->product)) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Argument %D \"%s\" is not a matrix obtained from MatProductCreate()",(arg),#A); \
276:   } while (0)
277: #else
278: #  define MatCheckProduct(A,arg) do {} while (0)
279: #endif

281: /*
282:   The stash is used to temporarily store inserted matrix values that
283:   belong to another processor. During the assembly phase the stashed
284:   values are moved to the correct processor and
285: */

287: typedef struct _MatStashSpace *PetscMatStashSpace;

289: struct _MatStashSpace {
290:   PetscMatStashSpace next;
291:   PetscScalar        *space_head,*val;
292:   PetscInt           *idx,*idy;
293:   PetscInt           total_space_size;
294:   PetscInt           local_used;
295:   PetscInt           local_remaining;
296: };

298: PETSC_EXTERN PetscErrorCode PetscMatStashSpaceGet(PetscInt,PetscInt,PetscMatStashSpace *);
299: PETSC_EXTERN PetscErrorCode PetscMatStashSpaceContiguous(PetscInt,PetscMatStashSpace *,PetscScalar *,PetscInt *,PetscInt *);
300: PETSC_EXTERN PetscErrorCode PetscMatStashSpaceDestroy(PetscMatStashSpace*);

302: typedef struct {
303:   PetscInt    count;
304: } MatStashHeader;

306: typedef struct {
307:   void        *buffer;          /* Of type blocktype, dynamically constructed  */
308:   PetscInt    count;
309:   char        pending;
310: } MatStashFrame;

312: typedef struct _MatStash MatStash;
313: struct _MatStash {
314:   PetscInt      nmax;                   /* maximum stash size */
315:   PetscInt      umax;                   /* user specified max-size */
316:   PetscInt      oldnmax;                /* the nmax value used previously */
317:   PetscInt      n;                      /* stash size */
318:   PetscInt      bs;                     /* block size of the stash */
319:   PetscInt      reallocs;               /* preserve the no of mallocs invoked */
320:   PetscMatStashSpace space_head,space;  /* linked list to hold stashed global row/column numbers and matrix values */

322:   PetscErrorCode (*ScatterBegin)(Mat,MatStash*,PetscInt*);
323:   PetscErrorCode (*ScatterGetMesg)(MatStash*,PetscMPIInt*,PetscInt**,PetscInt**,PetscScalar**,PetscInt*);
324:   PetscErrorCode (*ScatterEnd)(MatStash*);
325:   PetscErrorCode (*ScatterDestroy)(MatStash*);

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

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

372: #if !defined(PETSC_HAVE_MPIUNI)
373: PETSC_INTERN PetscErrorCode MatStashScatterDestroy_BTS(MatStash*);
374: #endif
375: PETSC_INTERN PetscErrorCode MatStashCreate_Private(MPI_Comm,PetscInt,MatStash*);
376: PETSC_INTERN PetscErrorCode MatStashDestroy_Private(MatStash*);
377: PETSC_INTERN PetscErrorCode MatStashScatterEnd_Private(MatStash*);
378: PETSC_INTERN PetscErrorCode MatStashSetInitialSize_Private(MatStash*,PetscInt);
379: PETSC_INTERN PetscErrorCode MatStashGetInfo_Private(MatStash*,PetscInt*,PetscInt*);
380: PETSC_INTERN PetscErrorCode MatStashValuesRow_Private(MatStash*,PetscInt,PetscInt,const PetscInt[],const PetscScalar[],PetscBool);
381: PETSC_INTERN PetscErrorCode MatStashValuesCol_Private(MatStash*,PetscInt,PetscInt,const PetscInt[],const PetscScalar[],PetscInt,PetscBool);
382: PETSC_INTERN PetscErrorCode MatStashValuesRowBlocked_Private(MatStash*,PetscInt,PetscInt,const PetscInt[],const PetscScalar[],PetscInt,PetscInt,PetscInt);
383: PETSC_INTERN PetscErrorCode MatStashValuesColBlocked_Private(MatStash*,PetscInt,PetscInt,const PetscInt[],const PetscScalar[],PetscInt,PetscInt,PetscInt);
384: PETSC_INTERN PetscErrorCode MatStashScatterBegin_Private(Mat,MatStash*,PetscInt*);
385: PETSC_INTERN PetscErrorCode MatStashScatterGetMesg_Private(MatStash*,PetscMPIInt*,PetscInt**,PetscInt**,PetscScalar**,PetscInt*);
386: PETSC_INTERN PetscErrorCode MatGetInfo_External(Mat,MatInfoType,MatInfo*);

388: typedef struct {
389:   PetscInt   dim;
390:   PetscInt   dims[4];
391:   PetscInt   starts[4];
392:   PetscBool  noc;        /* this is a single component problem, hence user will not set MatStencil.c */
393: } MatStencilInfo;

395: /* Info about using compressed row format */
396: typedef struct {
397:   PetscBool  use;                           /* indicates compressed rows have been checked and will be used */
398:   PetscInt   nrows;                         /* number of non-zero rows */
399:   PetscInt   *i;                            /* compressed row pointer  */
400:   PetscInt   *rindex;                       /* compressed row index               */
401: } Mat_CompressedRow;
402: PETSC_EXTERN PetscErrorCode MatCheckCompressedRow(Mat,PetscInt,Mat_CompressedRow*,PetscInt*,PetscInt,PetscReal);

404: typedef struct { /* used by MatCreateRedundantMatrix() for reusing matredundant */
405:   PetscInt     nzlocal,nsends,nrecvs;
406:   PetscMPIInt  *send_rank,*recv_rank;
407:   PetscInt     *sbuf_nz,*rbuf_nz,*sbuf_j,**rbuf_j;
408:   PetscScalar  *sbuf_a,**rbuf_a;
409:   MPI_Comm     subcomm;   /* when user does not provide a subcomm */
410:   IS           isrow,iscol;
411:   Mat          *matseq;
412: } Mat_Redundant;

414: typedef struct { /* used by MatProduct() */
415:   MatProductType type;
416:   char           *alg;
417:   Mat            A,B,C,Dwork;
418:   PetscReal      fill;
419:   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 */

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

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

430: #define CSRDataStructure(datatype)  \
431:   int         *i; \
432:   int         *j; \
433:   datatype    *a;\
434:   PetscInt    n;\
435:   PetscInt    ignorezeroentries;

437: typedef struct {
438:   CSRDataStructure(PetscScalar)
439: } PetscCSRDataStructure;

441: struct _p_SplitCSRMat {
442:   PetscInt              cstart,cend,rstart,rend;
443:   PetscCSRDataStructure diag,offdiag;
444:   PetscInt              *colmap;
445:   PetscMPIInt           rank;
446: };

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

500: PETSC_INTERN PetscErrorCode MatAXPY_Basic(Mat,PetscScalar,Mat,MatStructure);
501: PETSC_INTERN PetscErrorCode MatAXPY_BasicWithPreallocation(Mat,Mat,PetscScalar,Mat,MatStructure);
502: PETSC_INTERN PetscErrorCode MatAXPY_Basic_Preallocate(Mat,Mat,Mat*);
503: PETSC_INTERN PetscErrorCode MatAXPY_Dense_Nest(Mat,PetscScalar,Mat);

505: /*
506:     Utility for MatFactor (Schur complement)
507: */
508: PETSC_INTERN PetscErrorCode MatFactorFactorizeSchurComplement_Private(Mat);
509: PETSC_INTERN PetscErrorCode MatFactorInvertSchurComplement_Private(Mat);
510: PETSC_INTERN PetscErrorCode MatFactorUpdateSchurStatus_Private(Mat);
511: PETSC_INTERN PetscErrorCode MatFactorSetUpInPlaceSchur_Private(Mat);

513: /*
514:     Utility for MatZeroRows
515: */
516: PETSC_INTERN PetscErrorCode MatZeroRowsMapLocal_Private(Mat,PetscInt,const PetscInt*,PetscInt*,PetscInt**);

518: /*
519:     Utility for MatView/MatLoad
520: */
521: PETSC_INTERN PetscErrorCode MatView_Binary_BlockSizes(Mat,PetscViewer);
522: PETSC_INTERN PetscErrorCode MatLoad_Binary_BlockSizes(Mat,PetscViewer);


525: /*
526:     Object for partitioning graphs
527: */

529: typedef struct _MatPartitioningOps *MatPartitioningOps;
530: struct _MatPartitioningOps {
531:   PetscErrorCode (*apply)(MatPartitioning,IS*);
532:   PetscErrorCode (*applynd)(MatPartitioning,IS*);
533:   PetscErrorCode (*setfromoptions)(PetscOptionItems*,MatPartitioning);
534:   PetscErrorCode (*destroy)(MatPartitioning);
535:   PetscErrorCode (*view)(MatPartitioning,PetscViewer);
536:   PetscErrorCode (*improve)(MatPartitioning,IS*);
537: };

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

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

553: /*
554:     Object for coarsen graphs
555: */
556: typedef struct _MatCoarsenOps *MatCoarsenOps;
557: struct _MatCoarsenOps {
558:   PetscErrorCode (*apply)(MatCoarsen);
559:   PetscErrorCode (*setfromoptions)(PetscOptionItems*,MatCoarsen);
560:   PetscErrorCode (*destroy)(MatCoarsen);
561:   PetscErrorCode (*view)(MatCoarsen,PetscViewer);
562: };

564: struct _p_MatCoarsen {
565:   PETSCHEADER(struct _MatCoarsenOps);
566:   Mat              graph;
567:   void             *subctx;
568:   /* */
569:   PetscBool        strict_aggs;
570:   IS               perm;
571:   PetscCoarsenData *agg_lists;
572: };

574: /*
575:     MatFDColoring is used to compute Jacobian matrices efficiently
576:   via coloring. The data structure is explained below in an example.

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

590:     ncolors = 4;

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

599:     ncolumns      = {1,0,1,1}
600:     columns       = {{6},{},{4},{5}}
601:     nrows         = {3,0,2,2}
602:     rows          = {{0,1,2},{},{1,2},{1,2}}
603:     vwscale       = {dx(4),dx(5),dx(6)}              MPI Vec
604:     vscale        = {dx(0),dx(4),dx(5),dx(6)}        Seq Vec

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

609: */
610: typedef struct {
611:   PetscInt     row;
612:   PetscInt     col;
613:   PetscScalar  *valaddr;   /* address of value */
614: } MatEntry;

616: typedef struct {
617:   PetscInt     row;
618:   PetscScalar  *valaddr;   /* address of value */
619: } MatEntry2;

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

650: typedef struct _MatColoringOps *MatColoringOps;
651: struct _MatColoringOps {
652:   PetscErrorCode (*destroy)(MatColoring);
653:   PetscErrorCode (*setfromoptions)(PetscOptionItems*,MatColoring);
654:   PetscErrorCode (*view)(MatColoring,PetscViewer);
655:   PetscErrorCode (*apply)(MatColoring,ISColoring*);
656:   PetscErrorCode (*weights)(MatColoring,PetscReal**,PetscInt**);
657: };

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

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

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

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

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

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

738:   PetscErrorCode (*destroy)(Mat);
739: } Mat_SubSppt;

741: PETSC_EXTERN PetscErrorCode MatFactorDumpMatrix(Mat);
742: PETSC_INTERN PetscErrorCode MatShift_Basic(Mat,PetscScalar);
743: PETSC_INTERN PetscErrorCode MatSetBlockSizes_Default(Mat,PetscInt,PetscInt);

745: PETSC_STATIC_INLINE PetscErrorCode MatPivotCheck_nz(Mat mat,const MatFactorInfo *info,FactorShiftCtx *sctx,PetscInt row)
746: {
747:   PetscReal _rs   = sctx->rs;
748:   PetscReal _zero = info->zeropivot*_rs;

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

763: PETSC_STATIC_INLINE PetscErrorCode MatPivotCheck_pd(Mat mat,const MatFactorInfo *info,FactorShiftCtx *sctx,PetscInt row)
764: {
765:   PetscReal _rs   = sctx->rs;
766:   PetscReal _zero = info->zeropivot*_rs;

769:   if (PetscRealPart(sctx->pv) <= _zero && !PetscIsNanScalar(sctx->pv)){
770:     /* force matfactor to be diagonally dominant */
771:     if (sctx->nshift == sctx->nshift_max) {
772:       sctx->shift_fraction = sctx->shift_hi;
773:     } else {
774:       sctx->shift_lo = sctx->shift_fraction;
775:       sctx->shift_fraction = (sctx->shift_hi+sctx->shift_lo)/2.;
776:     }
777:     sctx->shift_amount = sctx->shift_fraction * sctx->shift_top;
778:     sctx->nshift++;
779:     sctx->newshift = PETSC_TRUE;
780:   } else {
781:     sctx->newshift = PETSC_FALSE;
782:   }
783:   return(0);
784: }

786: PETSC_STATIC_INLINE PetscErrorCode MatPivotCheck_inblocks(Mat mat,const MatFactorInfo *info,FactorShiftCtx *sctx,PetscInt row)
787: {
788:   PetscReal _zero = info->zeropivot;

791:   if (PetscAbsScalar(sctx->pv) <= _zero && !PetscIsNanScalar(sctx->pv)){
792:     sctx->pv          += info->shiftamount;
793:     sctx->shift_amount = 0.0;
794:     sctx->nshift++;
795:   }
796:   sctx->newshift = PETSC_FALSE;
797:   return(0);
798: }

800: PETSC_STATIC_INLINE PetscErrorCode MatPivotCheck_none(Mat fact,Mat mat,const MatFactorInfo *info,FactorShiftCtx *sctx,PetscInt row)
801: {
802:   PetscReal      _zero = info->zeropivot;

806:   sctx->newshift = PETSC_FALSE;
807:   if (PetscAbsScalar(sctx->pv) <= _zero && !PetscIsNanScalar(sctx->pv)) {
808:     if (!mat->erroriffailure) {
809:       PetscInfo3(mat,"Detected zero pivot in factorization in row %D value %g tolerance %g\n",row,(double)PetscAbsScalar(sctx->pv),(double)_zero);
810:       fact->factorerrortype             = MAT_FACTOR_NUMERIC_ZEROPIVOT;
811:       fact->factorerror_zeropivot_value = PetscAbsScalar(sctx->pv);
812:       fact->factorerror_zeropivot_row   = row;
813:     } else SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot row %D value %g tolerance %g\n",row,(double)PetscAbsScalar(sctx->pv),(double)_zero);
814:   }
815:   return(0);
816: }

818: PETSC_STATIC_INLINE PetscErrorCode MatPivotCheck(Mat fact,Mat mat,const MatFactorInfo *info,FactorShiftCtx *sctx,PetscInt row)
819: {

823:   if (info->shifttype == (PetscReal) MAT_SHIFT_NONZERO){
824:     MatPivotCheck_nz(mat,info,sctx,row);
825:   } else if (info->shifttype == (PetscReal) MAT_SHIFT_POSITIVE_DEFINITE){
826:     MatPivotCheck_pd(mat,info,sctx,row);
827:   } else if (info->shifttype == (PetscReal) MAT_SHIFT_INBLOCKS){
828:     MatPivotCheck_inblocks(mat,info,sctx,row);
829:   } else {
830:     MatPivotCheck_none(fact,mat,info,sctx,row);
831:   }
832:   return(0);
833: }

835: /*
836:   Create and initialize a linked list
837:   Input Parameters:
838:     idx_start - starting index of the list
839:     lnk_max   - max value of lnk indicating the end of the list
840:     nlnk      - max length of the list
841:   Output Parameters:
842:     lnk       - list initialized
843:     bt        - PetscBT (bitarray) with all bits set to false
844:     lnk_empty - flg indicating the list is empty
845: */
846: #define PetscLLCreate(idx_start,lnk_max,nlnk,lnk,bt) \
847:   (PetscMalloc1(nlnk,&lnk) || PetscBTCreate(nlnk,&(bt)) || (lnk[idx_start] = lnk_max,0))

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

852: /*
853:   Add an index set into a sorted linked list
854:   Input Parameters:
855:     nidx      - number of input indices
856:     indices   - integer array
857:     idx_start - starting index of the list
858:     lnk       - linked list(an integer array) that is created
859:     bt        - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
860:   output Parameters:
861:     nlnk      - number of newly added indices
862:     lnk       - the sorted(increasing order) linked list containing new and non-redundate entries from indices
863:     bt        - updated PetscBT (bitarray)
864: */
865: #define PetscLLAdd(nidx,indices,idx_start,nlnk,lnk,bt) 0;\
866: {\
867:   PetscInt _k,_entry,_location,_lnkdata;\
868:   nlnk     = 0;\
869:   _lnkdata = idx_start;\
870:   for (_k=0; _k<nidx; _k++){\
871:     _entry = indices[_k];\
872:     if (!PetscBTLookupSet(bt,_entry)){  /* new entry */\
873:       /* search for insertion location */\
874:       /* start from the beginning if _entry < previous _entry */\
875:       if (_k && _entry < _lnkdata) _lnkdata  = idx_start;\
876:       do {\
877:         _location = _lnkdata;\
878:         _lnkdata  = lnk[_location];\
879:       } while (_entry > _lnkdata);\
880:       /* insertion location is found, add entry into lnk */\
881:       lnk[_location] = _entry;\
882:       lnk[_entry]    = _lnkdata;\
883:       nlnk++;\
884:       _lnkdata = _entry; /* next search starts from here if next_entry > _entry */\
885:     }\
886:   }\
887: }

889: /*
890:   Add a permuted index set into a sorted linked list
891:   Input Parameters:
892:     nidx      - number of input indices
893:     indices   - integer array
894:     perm      - permutation of indices
895:     idx_start - starting index of the list
896:     lnk       - linked list(an integer array) that is created
897:     bt        - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
898:   output Parameters:
899:     nlnk      - number of newly added indices
900:     lnk       - the sorted(increasing order) linked list containing new and non-redundate entries from indices
901:     bt        - updated PetscBT (bitarray)
902: */
903: #define PetscLLAddPerm(nidx,indices,perm,idx_start,nlnk,lnk,bt) 0;\
904: {\
905:   PetscInt _k,_entry,_location,_lnkdata;\
906:   nlnk     = 0;\
907:   _lnkdata = idx_start;\
908:   for (_k=0; _k<nidx; _k++){\
909:     _entry = perm[indices[_k]];\
910:     if (!PetscBTLookupSet(bt,_entry)){  /* new entry */\
911:       /* search for insertion location */\
912:       /* start from the beginning if _entry < previous _entry */\
913:       if (_k && _entry < _lnkdata) _lnkdata  = idx_start;\
914:       do {\
915:         _location = _lnkdata;\
916:         _lnkdata  = lnk[_location];\
917:       } while (_entry > _lnkdata);\
918:       /* insertion location is found, add entry into lnk */\
919:       lnk[_location] = _entry;\
920:       lnk[_entry]    = _lnkdata;\
921:       nlnk++;\
922:       _lnkdata = _entry; /* next search starts from here if next_entry > _entry */\
923:     }\
924:   }\
925: }

927: /*
928:   Add a SORTED ascending index set into a sorted linked list - same as PetscLLAdd() bus skip 'if (_k && _entry < _lnkdata) _lnkdata  = idx_start;'
929:   Input Parameters:
930:     nidx      - number of input indices
931:     indices   - sorted integer array
932:     idx_start - starting index of the list
933:     lnk       - linked list(an integer array) that is created
934:     bt        - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
935:   output Parameters:
936:     nlnk      - number of newly added indices
937:     lnk       - the sorted(increasing order) linked list containing new and non-redundate entries from indices
938:     bt        - updated PetscBT (bitarray)
939: */
940: #define PetscLLAddSorted(nidx,indices,idx_start,nlnk,lnk,bt) 0;\
941: {\
942:   PetscInt _k,_entry,_location,_lnkdata;\
943:   nlnk      = 0;\
944:   _lnkdata  = idx_start;\
945:   for (_k=0; _k<nidx; _k++){\
946:     _entry = indices[_k];\
947:     if (!PetscBTLookupSet(bt,_entry)){  /* new entry */\
948:       /* search for insertion location */\
949:       do {\
950:         _location = _lnkdata;\
951:         _lnkdata  = lnk[_location];\
952:       } while (_entry > _lnkdata);\
953:       /* insertion location is found, add entry into lnk */\
954:       lnk[_location] = _entry;\
955:       lnk[_entry]    = _lnkdata;\
956:       nlnk++;\
957:       _lnkdata = _entry; /* next search starts from here */\
958:     }\
959:   }\
960: }

962: #define PetscLLAddSorted_new(nidx,indices,idx_start,lnk_empty,nlnk,lnk,bt) 0; \
963: {\
964:   PetscInt _k,_entry,_location,_lnkdata;\
965:   if (lnk_empty){\
966:     _lnkdata  = idx_start;                      \
967:     for (_k=0; _k<nidx; _k++){                  \
968:       _entry = indices[_k];                             \
969:       PetscBTSet(bt,_entry);  /* mark the new entry */          \
970:           _location = _lnkdata;                                 \
971:           _lnkdata  = lnk[_location];                           \
972:         /* insertion location is found, add entry into lnk */   \
973:         lnk[_location] = _entry;                                \
974:         lnk[_entry]    = _lnkdata;                              \
975:         _lnkdata = _entry; /* next search starts from here */   \
976:     }                                                           \
977:     /*\
978:     lnk[indices[nidx-1]] = lnk[idx_start];\
979:     lnk[idx_start]       = indices[0];\
980:     PetscBTSet(bt,indices[0]);  \
981:     for (_k=1; _k<nidx; _k++){                  \
982:       PetscBTSet(bt,indices[_k]);                                          \
983:       lnk[indices[_k-1]] = indices[_k];                                  \
984:     }                                                           \
985:      */\
986:     nlnk      = nidx;\
987:     lnk_empty = PETSC_FALSE;\
988:   } else {\
989:     nlnk      = 0;                              \
990:     _lnkdata  = idx_start;                      \
991:     for (_k=0; _k<nidx; _k++){                  \
992:       _entry = indices[_k];                             \
993:       if (!PetscBTLookupSet(bt,_entry)){  /* new entry */       \
994:         /* search for insertion location */                     \
995:         do {                                                    \
996:           _location = _lnkdata;                                 \
997:           _lnkdata  = lnk[_location];                           \
998:         } while (_entry > _lnkdata);                            \
999:         /* insertion location is found, add entry into lnk */   \
1000:         lnk[_location] = _entry;                                \
1001:         lnk[_entry]    = _lnkdata;                              \
1002:         nlnk++;                                                 \
1003:         _lnkdata = _entry; /* next search starts from here */   \
1004:       }                                                         \
1005:     }                                                           \
1006:   }                                                             \
1007: }

1009: /*
1010:   Add a SORTED index set into a sorted linked list used for LUFactorSymbolic()
1011:   Same as PetscLLAddSorted() with an additional operation:
1012:        count the number of input indices that are no larger than 'diag'
1013:   Input Parameters:
1014:     indices   - sorted integer array
1015:     idx_start - starting index of the list, index of pivot row
1016:     lnk       - linked list(an integer array) that is created
1017:     bt        - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
1018:     diag      - index of the active row in LUFactorSymbolic
1019:     nzbd      - number of input indices with indices <= idx_start
1020:     im        - im[idx_start] is initialized as num of nonzero entries in row=idx_start
1021:   output Parameters:
1022:     nlnk      - number of newly added indices
1023:     lnk       - the sorted(increasing order) linked list containing new and non-redundate entries from indices
1024:     bt        - updated PetscBT (bitarray)
1025:     im        - im[idx_start]: unchanged if diag is not an entry
1026:                              : num of entries with indices <= diag if diag is an entry
1027: */
1028: #define PetscLLAddSortedLU(indices,idx_start,nlnk,lnk,bt,diag,nzbd,im) 0;\
1029: {\
1030:   PetscInt _k,_entry,_location,_lnkdata,_nidx;\
1031:   nlnk     = 0;\
1032:   _lnkdata = idx_start;\
1033:   _nidx = im[idx_start] - nzbd; /* num of entries with idx_start < index <= diag */\
1034:   for (_k=0; _k<_nidx; _k++){\
1035:     _entry = indices[_k];\
1036:     nzbd++;\
1037:     if (_entry== diag) im[idx_start] = nzbd;\
1038:     if (!PetscBTLookupSet(bt,_entry)){  /* new entry */\
1039:       /* search for insertion location */\
1040:       do {\
1041:         _location = _lnkdata;\
1042:         _lnkdata  = lnk[_location];\
1043:       } while (_entry > _lnkdata);\
1044:       /* insertion location is found, add entry into lnk */\
1045:       lnk[_location] = _entry;\
1046:       lnk[_entry]    = _lnkdata;\
1047:       nlnk++;\
1048:       _lnkdata = _entry; /* next search starts from here */\
1049:     }\
1050:   }\
1051: }

1053: /*
1054:   Copy data on the list into an array, then initialize the list
1055:   Input Parameters:
1056:     idx_start - starting index of the list
1057:     lnk_max   - max value of lnk indicating the end of the list
1058:     nlnk      - number of data on the list to be copied
1059:     lnk       - linked list
1060:     bt        - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
1061:   output Parameters:
1062:     indices   - array that contains the copied data
1063:     lnk       - linked list that is cleaned and initialize
1064:     bt        - PetscBT (bitarray) with all bits set to false
1065: */
1066: #define PetscLLClean(idx_start,lnk_max,nlnk,lnk,indices,bt) 0;\
1067: {\
1068:   PetscInt _j,_idx=idx_start;\
1069:   for (_j=0; _j<nlnk; _j++){\
1070:     _idx = lnk[_idx];\
1071:     indices[_j] = _idx;\
1072:     PetscBTClear(bt,_idx);\
1073:   }\
1074:   lnk[idx_start] = lnk_max;\
1075: }
1076: /*
1077:   Free memories used by the list
1078: */
1079: #define PetscLLDestroy(lnk,bt) (PetscFree(lnk) || PetscBTDestroy(&(bt)))

1081: /* Routines below are used for incomplete matrix factorization */
1082: /*
1083:   Create and initialize a linked list and its levels
1084:   Input Parameters:
1085:     idx_start - starting index of the list
1086:     lnk_max   - max value of lnk indicating the end of the list
1087:     nlnk      - max length of the list
1088:   Output Parameters:
1089:     lnk       - list initialized
1090:     lnk_lvl   - array of size nlnk for storing levels of lnk
1091:     bt        - PetscBT (bitarray) with all bits set to false
1092: */
1093: #define PetscIncompleteLLCreate(idx_start,lnk_max,nlnk,lnk,lnk_lvl,bt)\
1094:   (PetscIntMultError(2,nlnk,NULL) || PetscMalloc1(2*nlnk,&lnk) || PetscBTCreate(nlnk,&(bt)) || (lnk[idx_start] = lnk_max,lnk_lvl = lnk + nlnk,0))

1096: /*
1097:   Initialize a sorted linked list used for ILU and ICC
1098:   Input Parameters:
1099:     nidx      - number of input idx
1100:     idx       - integer array used for storing column indices
1101:     idx_start - starting index of the list
1102:     perm      - indices of an IS
1103:     lnk       - linked list(an integer array) that is created
1104:     lnklvl    - levels of lnk
1105:     bt        - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
1106:   output Parameters:
1107:     nlnk     - number of newly added idx
1108:     lnk      - the sorted(increasing order) linked list containing new and non-redundate entries from idx
1109:     lnklvl   - levels of lnk
1110:     bt       - updated PetscBT (bitarray)
1111: */
1112: #define PetscIncompleteLLInit(nidx,idx,idx_start,perm,nlnk,lnk,lnklvl,bt) 0;\
1113: {\
1114:   PetscInt _k,_entry,_location,_lnkdata;\
1115:   nlnk     = 0;\
1116:   _lnkdata = idx_start;\
1117:   for (_k=0; _k<nidx; _k++){\
1118:     _entry = perm[idx[_k]];\
1119:     if (!PetscBTLookupSet(bt,_entry)){  /* new entry */\
1120:       /* search for insertion location */\
1121:       if (_k && _entry < _lnkdata) _lnkdata  = idx_start;\
1122:       do {\
1123:         _location = _lnkdata;\
1124:         _lnkdata  = lnk[_location];\
1125:       } while (_entry > _lnkdata);\
1126:       /* insertion location is found, add entry into lnk */\
1127:       lnk[_location]  = _entry;\
1128:       lnk[_entry]     = _lnkdata;\
1129:       lnklvl[_entry] = 0;\
1130:       nlnk++;\
1131:       _lnkdata = _entry; /* next search starts from here if next_entry > _entry */\
1132:     }\
1133:   }\
1134: }

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

1154:   Note: the level of factor(i,j) is set as lvl(i,j) = min{ lvl(i,j), lvl(i,prow)+lvl(prow,j)+1)
1155:         where idx = non-zero columns of U(prow,prow+1:n-1), prow<i
1156: */
1157: #define PetscILULLAddSorted(nidx,idx,level,idxlvl,idx_start,nlnk,lnk,lnklvl,bt,lnklvl_prow) 0;\
1158: {\
1159:   PetscInt _k,_entry,_location,_lnkdata,_incrlev,_lnklvl_prow=lnklvl[prow];\
1160:   nlnk     = 0;\
1161:   _lnkdata = idx_start;\
1162:   for (_k=0; _k<nidx; _k++){\
1163:     _incrlev = idxlvl[_k] + _lnklvl_prow + 1;\
1164:     if (_incrlev > level) continue;\
1165:     _entry = idx[_k];\
1166:     if (!PetscBTLookupSet(bt,_entry)){  /* new entry */\
1167:       /* search for insertion location */\
1168:       do {\
1169:         _location = _lnkdata;\
1170:         _lnkdata  = lnk[_location];\
1171:       } while (_entry > _lnkdata);\
1172:       /* insertion location is found, add entry into lnk */\
1173:       lnk[_location]  = _entry;\
1174:       lnk[_entry]     = _lnkdata;\
1175:       lnklvl[_entry] = _incrlev;\
1176:       nlnk++;\
1177:       _lnkdata = _entry; /* next search starts from here if next_entry > _entry */\
1178:     } else { /* existing entry: update lnklvl */\
1179:       if (lnklvl[_entry] > _incrlev) lnklvl[_entry] = _incrlev;\
1180:     }\
1181:   }\
1182: }

1184: /*
1185:   Add a index set into a sorted linked list
1186:   Input Parameters:
1187:     nidx      - number of input idx
1188:     idx   - integer array used for storing column indices
1189:     level     - level of fill, e.g., ICC(level)
1190:     idxlvl - level of idx
1191:     idx_start - starting index of the list
1192:     lnk       - linked list(an integer array) that is created
1193:     lnklvl   - levels of lnk
1194:     bt        - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
1195:   output Parameters:
1196:     nlnk      - number of newly added idx
1197:     lnk       - the sorted(increasing order) linked list containing new and non-redundate entries from idx
1198:     lnklvl   - levels of lnk
1199:     bt        - updated PetscBT (bitarray)
1200: */
1201: #define PetscIncompleteLLAdd(nidx,idx,level,idxlvl,idx_start,nlnk,lnk,lnklvl,bt) 0;\
1202: {\
1203:   PetscInt _k,_entry,_location,_lnkdata,_incrlev;\
1204:   nlnk     = 0;\
1205:   _lnkdata = idx_start;\
1206:   for (_k=0; _k<nidx; _k++){\
1207:     _incrlev = idxlvl[_k] + 1;\
1208:     if (_incrlev > level) continue;\
1209:     _entry = idx[_k];\
1210:     if (!PetscBTLookupSet(bt,_entry)){  /* new entry */\
1211:       /* search for insertion location */\
1212:       if (_k && _entry < _lnkdata) _lnkdata  = idx_start;\
1213:       do {\
1214:         _location = _lnkdata;\
1215:         _lnkdata  = lnk[_location];\
1216:       } while (_entry > _lnkdata);\
1217:       /* insertion location is found, add entry into lnk */\
1218:       lnk[_location]  = _entry;\
1219:       lnk[_entry]     = _lnkdata;\
1220:       lnklvl[_entry] = _incrlev;\
1221:       nlnk++;\
1222:       _lnkdata = _entry; /* next search starts from here if next_entry > _entry */\
1223:     } else { /* existing entry: update lnklvl */\
1224:       if (lnklvl[_entry] > _incrlev) lnklvl[_entry] = _incrlev;\
1225:     }\
1226:   }\
1227: }

1229: /*
1230:   Add a SORTED index set into a sorted linked list
1231:   Input Parameters:
1232:     nidx      - number of input indices
1233:     idx   - sorted integer array used for storing column indices
1234:     level     - level of fill, e.g., ICC(level)
1235:     idxlvl - level of idx
1236:     idx_start - starting index of the list
1237:     lnk       - linked list(an integer array) that is created
1238:     lnklvl    - levels of lnk
1239:     bt        - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
1240:   output Parameters:
1241:     nlnk      - number of newly added idx
1242:     lnk       - the sorted(increasing order) linked list containing new and non-redundate entries from idx
1243:     lnklvl    - levels of lnk
1244:     bt        - updated PetscBT (bitarray)
1245: */
1246: #define PetscIncompleteLLAddSorted(nidx,idx,level,idxlvl,idx_start,nlnk,lnk,lnklvl,bt) 0;\
1247: {\
1248:   PetscInt _k,_entry,_location,_lnkdata,_incrlev;\
1249:   nlnk = 0;\
1250:   _lnkdata = idx_start;\
1251:   for (_k=0; _k<nidx; _k++){\
1252:     _incrlev = idxlvl[_k] + 1;\
1253:     if (_incrlev > level) continue;\
1254:     _entry = idx[_k];\
1255:     if (!PetscBTLookupSet(bt,_entry)){  /* new entry */\
1256:       /* search for insertion location */\
1257:       do {\
1258:         _location = _lnkdata;\
1259:         _lnkdata  = lnk[_location];\
1260:       } while (_entry > _lnkdata);\
1261:       /* insertion location is found, add entry into lnk */\
1262:       lnk[_location] = _entry;\
1263:       lnk[_entry]    = _lnkdata;\
1264:       lnklvl[_entry] = _incrlev;\
1265:       nlnk++;\
1266:       _lnkdata = _entry; /* next search starts from here */\
1267:     } else { /* existing entry: update lnklvl */\
1268:       if (lnklvl[_entry] > _incrlev) lnklvl[_entry] = _incrlev;\
1269:     }\
1270:   }\
1271: }

1273: /*
1274:   Add a SORTED index set into a sorted linked list for ICC
1275:   Input Parameters:
1276:     nidx      - number of input indices
1277:     idx       - sorted integer array used for storing column indices
1278:     level     - level of fill, e.g., ICC(level)
1279:     idxlvl    - level of idx
1280:     idx_start - starting index of the list
1281:     lnk       - linked list(an integer array) that is created
1282:     lnklvl    - levels of lnk
1283:     bt        - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
1284:     idxlvl_prow - idxlvl[prow], where prow is the row number of the idx
1285:   output Parameters:
1286:     nlnk   - number of newly added indices
1287:     lnk    - the sorted(increasing order) linked list containing new and non-redundate entries from idx
1288:     lnklvl - levels of lnk
1289:     bt     - updated PetscBT (bitarray)
1290:   Note: the level of U(i,j) is set as lvl(i,j) = min{ lvl(i,j), lvl(prow,i)+lvl(prow,j)+1)
1291:         where idx = non-zero columns of U(prow,prow+1:n-1), prow<i
1292: */
1293: #define PetscICCLLAddSorted(nidx,idx,level,idxlvl,idx_start,nlnk,lnk,lnklvl,bt,idxlvl_prow) 0;\
1294: {\
1295:   PetscInt _k,_entry,_location,_lnkdata,_incrlev;\
1296:   nlnk = 0;\
1297:   _lnkdata = idx_start;\
1298:   for (_k=0; _k<nidx; _k++){\
1299:     _incrlev = idxlvl[_k] + idxlvl_prow + 1;\
1300:     if (_incrlev > level) continue;\
1301:     _entry = idx[_k];\
1302:     if (!PetscBTLookupSet(bt,_entry)){  /* new entry */\
1303:       /* search for insertion location */\
1304:       do {\
1305:         _location = _lnkdata;\
1306:         _lnkdata  = lnk[_location];\
1307:       } while (_entry > _lnkdata);\
1308:       /* insertion location is found, add entry into lnk */\
1309:       lnk[_location] = _entry;\
1310:       lnk[_entry]    = _lnkdata;\
1311:       lnklvl[_entry] = _incrlev;\
1312:       nlnk++;\
1313:       _lnkdata = _entry; /* next search starts from here */\
1314:     } else { /* existing entry: update lnklvl */\
1315:       if (lnklvl[_entry] > _incrlev) lnklvl[_entry] = _incrlev;\
1316:     }\
1317:   }\
1318: }

1320: /*
1321:   Copy data on the list into an array, then initialize the list
1322:   Input Parameters:
1323:     idx_start - starting index of the list
1324:     lnk_max   - max value of lnk indicating the end of the list
1325:     nlnk      - number of data on the list to be copied
1326:     lnk       - linked list
1327:     lnklvl    - level of lnk
1328:     bt        - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
1329:   output Parameters:
1330:     indices - array that contains the copied data
1331:     lnk     - linked list that is cleaned and initialize
1332:     lnklvl  - level of lnk that is reinitialized
1333:     bt      - PetscBT (bitarray) with all bits set to false
1334: */
1335: #define PetscIncompleteLLClean(idx_start,lnk_max,nlnk,lnk,lnklvl,indices,indiceslvl,bt) 0;\
1336: do {\
1337:   PetscInt _j,_idx=idx_start;\
1338:   for (_j=0; _j<nlnk; _j++){\
1339:     _idx = lnk[_idx];\
1340:     *(indices+_j) = _idx;\
1341:     *(indiceslvl+_j) = lnklvl[_idx];\
1342:     lnklvl[_idx] = -1;\
1343:     PetscBTClear(bt,_idx);\
1344:   }\
1345:   lnk[idx_start] = lnk_max;\
1346: } while (0)
1347: /*
1348:   Free memories used by the list
1349: */
1350: #define PetscIncompleteLLDestroy(lnk,bt) (PetscFree(lnk) || PetscBTDestroy(&(bt)))

1352: #define MatCheckSameLocalSize(A,ar1,B,ar2) do { \
1354:   if ((A->rmap->n != B->rmap->n) || (A->cmap->n != B->cmap->n)) SETERRQ6(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Incompatible matrix local sizes: parameter # %d (%D x %D) != parameter # %d (%D x %D)",ar1,A->rmap->n,A->cmap->n,ar2,B->rmap->n,B->cmap->n);} while (0)

1356: #define MatCheckSameSize(A,ar1,B,ar2) do { \
1357:   if ((A->rmap->N != B->rmap->N) || (A->cmap->N != B->cmap->N)) SETERRQ6(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_INCOMP,"Incompatible matrix global sizes: parameter # %d (%D x %D) != parameter # %d (%D x %D)",ar1,A->rmap->N,A->cmap->N,ar2,B->rmap->N,B->cmap->N);\
1358:   MatCheckSameLocalSize(A,ar1,B,ar2);} while (0)

1360: #define VecCheckMatCompatible(M,x,ar1,b,ar2) do { \
1361:   if (M->cmap->N != x->map->N) SETERRQ3(PetscObjectComm((PetscObject)M),PETSC_ERR_ARG_SIZ,"Vector global length incompatible with matrix: parameter # %d global size %D != matrix column global size %D",ar1,x->map->N,M->cmap->N); \
1362:   if (M->rmap->N != b->map->N) SETERRQ3(PetscObjectComm((PetscObject)M),PETSC_ERR_ARG_SIZ,"Vector global length incompatible with matrix: parameter # %d global size %D != matrix row global size %D",ar2,b->map->N,M->rmap->N);} while (0)

1364: /* -------------------------------------------------------------------------------------------------------*/
1365: #include <petscbt.h>
1366: /*
1367:   Create and initialize a condensed linked list -
1368:     same as PetscLLCreate(), but uses a scalable array 'lnk' with size of max number of entries, not O(N).
1369:     Barry suggested this approach (Dec. 6, 2011):
1370:       I've thought of an alternative way of representing a linked list that is efficient but doesn't have the O(N) scaling issue
1371:       (it may be faster than the O(N) even sequentially due to less crazy memory access).

1373:       Instead of having some like  a  2  -> 4 -> 11 ->  22  list that uses slot 2  4 11 and 22 in a big array use a small array with two slots
1374:       for each entry for example  [ 2 1 | 4 3 | 22 -1 | 11 2]   so the first number (of the pair) is the value while the second tells you where
1375:       in the list the next entry is. Inserting a new link means just append another pair at the end. For example say we want to insert 13 into the
1376:       list it would then become [2 1 | 4 3 | 22 -1 | 11 4 | 13 2 ] you just add a pair at the end and fix the point for the one that points to it.
1377:       That is 11 use to point to the 2 slot, after the change 11 points to the 4th slot which has the value 13. Note that values are always next
1378:       to each other so memory access is much better than using the big array.

1380:   Example:
1381:      nlnk_max=5, lnk_max=36:
1382:      Initial list: [0, 0 | 36, 2 | 0, 0 | 0, 0 | 0, 0 | 0, 0 | 0, 0]
1383:      here, head_node has index 2 with value lnk[2]=lnk_max=36,
1384:            0-th entry is used to store the number of entries in the list,
1385:      The initial lnk represents head -> tail(marked by 36) with number of entries = lnk[0]=0.

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

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

1395:   Input Parameters:
1396:     nlnk_max  - max length of the list
1397:     lnk_max   - max value of the entries
1398:   Output Parameters:
1399:     lnk       - list created and initialized
1400:     bt        - PetscBT (bitarray) with all bits set to false. Note: bt has size lnk_max, not nln_max!
1401: */
1402: PETSC_STATIC_INLINE PetscErrorCode PetscLLCondensedCreate(PetscInt nlnk_max,PetscInt lnk_max,PetscInt **lnk,PetscBT *bt)
1403: {
1405:   PetscInt       *llnk,lsize = 0;

1408:   PetscIntMultError(2,nlnk_max+2,&lsize);
1409:   PetscMalloc1(lsize,lnk);
1410:   PetscBTCreate(lnk_max,bt);
1411:   llnk = *lnk;
1412:   llnk[0] = 0;         /* number of entries on the list */
1413:   llnk[2] = lnk_max;   /* value in the head node */
1414:   llnk[3] = 2;         /* next for the head node */
1415:   return(0);
1416: }

1418: /*
1419:   Add a SORTED ascending index set into a sorted linked list. See PetscLLCondensedCreate() for detailed description.
1420:   Input Parameters:
1421:     nidx      - number of input indices
1422:     indices   - sorted integer array
1423:     lnk       - condensed linked list(an integer array) that is created
1424:     bt        - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
1425:   output Parameters:
1426:     lnk       - the sorted(increasing order) linked list containing previous and newly added non-redundate indices
1427:     bt        - updated PetscBT (bitarray)
1428: */
1429: PETSC_STATIC_INLINE PetscErrorCode PetscLLCondensedAddSorted(PetscInt nidx,const PetscInt indices[],PetscInt lnk[],PetscBT bt)
1430: {
1431:   PetscInt _k,_entry,_location,_next,_lnkdata,_nlnk,_newnode;

1434:   _nlnk     = lnk[0]; /* num of entries on the input lnk */
1435:   _location = 2; /* head */
1436:     for (_k=0; _k<nidx; _k++){
1437:       _entry = indices[_k];
1438:       if (!PetscBTLookupSet(bt,_entry)){  /* new entry */
1439:         /* search for insertion location */
1440:         do {
1441:           _next     = _location + 1; /* link from previous node to next node */
1442:           _location = lnk[_next];    /* idx of next node */
1443:           _lnkdata  = lnk[_location];/* value of next node */
1444:         } while (_entry > _lnkdata);
1445:         /* insertion location is found, add entry into lnk */
1446:         _newnode        = 2*(_nlnk+2);   /* index for this new node */
1447:         lnk[_next]      = _newnode;      /* connect previous node to the new node */
1448:         lnk[_newnode]   = _entry;        /* set value of the new node */
1449:         lnk[_newnode+1] = _location;     /* connect new node to next node */
1450:         _location       = _newnode;      /* next search starts from the new node */
1451:         _nlnk++;
1452:       }   \
1453:     }\
1454:   lnk[0]   = _nlnk;   /* number of entries in the list */
1455:   return(0);
1456: }

1458: PETSC_STATIC_INLINE PetscErrorCode PetscLLCondensedClean(PetscInt lnk_max,PetscInt nidx,PetscInt *indices,PetscInt lnk[],PetscBT bt)
1459: {
1461:   PetscInt       _k,_next,_nlnk;

1464:   _next = lnk[3];       /* head node */
1465:   _nlnk = lnk[0];       /* num of entries on the list */
1466:   for (_k=0; _k<_nlnk; _k++){
1467:     indices[_k] = lnk[_next];
1468:     _next       = lnk[_next + 1];
1469:     PetscBTClear(bt,indices[_k]);
1470:   }
1471:   lnk[0] = 0;          /* num of entries on the list */
1472:   lnk[2] = lnk_max;    /* initialize head node */
1473:   lnk[3] = 2;          /* head node */
1474:   return(0);
1475: }

1477: PETSC_STATIC_INLINE PetscErrorCode PetscLLCondensedView(PetscInt *lnk)
1478: {
1480:   PetscInt       k;

1483:   PetscPrintf(PETSC_COMM_SELF,"LLCondensed of size %D, (val,  next)\n",lnk[0]);
1484:   for (k=2; k< lnk[0]+2; k++){
1485:     PetscPrintf(PETSC_COMM_SELF," %D: (%D, %D)\n",2*k,lnk[2*k],lnk[2*k+1]);
1486:   }
1487:   return(0);
1488: }

1490: /*
1491:   Free memories used by the list
1492: */
1493: PETSC_STATIC_INLINE PetscErrorCode PetscLLCondensedDestroy(PetscInt *lnk,PetscBT bt)
1494: {

1498:   PetscFree(lnk);
1499:   PetscBTDestroy(&bt);
1500:   return(0);
1501: }

1503: /* -------------------------------------------------------------------------------------------------------*/
1504: /*
1505:  Same as PetscLLCondensedCreate(), but does not use non-scalable O(lnk_max) bitarray
1506:   Input Parameters:
1507:     nlnk_max  - max length of the list
1508:   Output Parameters:
1509:     lnk       - list created and initialized
1510: */
1511: PETSC_STATIC_INLINE PetscErrorCode PetscLLCondensedCreate_Scalable(PetscInt nlnk_max,PetscInt **lnk)
1512: {
1514:   PetscInt       *llnk,lsize = 0;

1517:   PetscIntMultError(2,nlnk_max+2,&lsize);
1518:   PetscMalloc1(lsize,lnk);
1519:   llnk = *lnk;
1520:   llnk[0] = 0;               /* number of entries on the list */
1521:   llnk[2] = PETSC_MAX_INT;   /* value in the head node */
1522:   llnk[3] = 2;               /* next for the head node */
1523:   return(0);
1524: }

1526: PETSC_STATIC_INLINE PetscErrorCode PetscLLCondensedExpand_Scalable(PetscInt nlnk_max,PetscInt **lnk)
1527: {
1529:   PetscInt       lsize = 0;

1532:   PetscIntMultError(2,nlnk_max+2,&lsize);
1533:   PetscRealloc(lsize*sizeof(PetscInt),lnk);
1534:   return(0);
1535: }

1537: PETSC_STATIC_INLINE PetscErrorCode PetscLLCondensedAddSorted_Scalable(PetscInt nidx,const PetscInt indices[],PetscInt lnk[])
1538: {
1539:   PetscInt _k,_entry,_location,_next,_lnkdata,_nlnk,_newnode;
1540:   _nlnk     = lnk[0]; /* num of entries on the input lnk */
1541:   _location = 2; /* head */ \
1542:     for (_k=0; _k<nidx; _k++){
1543:       _entry = indices[_k];
1544:       /* search for insertion location */
1545:       do {
1546:         _next     = _location + 1; /* link from previous node to next node */
1547:         _location = lnk[_next];    /* idx of next node */
1548:         _lnkdata  = lnk[_location];/* value of next node */
1549:       } while (_entry > _lnkdata);
1550:       if (_entry < _lnkdata) {
1551:         /* insertion location is found, add entry into lnk */
1552:         _newnode        = 2*(_nlnk+2);   /* index for this new node */
1553:         lnk[_next]      = _newnode;      /* connect previous node to the new node */
1554:         lnk[_newnode]   = _entry;        /* set value of the new node */
1555:         lnk[_newnode+1] = _location;     /* connect new node to next node */
1556:         _location       = _newnode;      /* next search starts from the new node */
1557:         _nlnk++;
1558:       }
1559:     }
1560:   lnk[0]   = _nlnk;   /* number of entries in the list */
1561:   return 0;
1562: }

1564: PETSC_STATIC_INLINE PetscErrorCode PetscLLCondensedClean_Scalable(PetscInt nidx,PetscInt *indices,PetscInt *lnk)
1565: {
1566:   PetscInt _k,_next,_nlnk;
1567:   _next = lnk[3];       /* head node */
1568:   _nlnk = lnk[0];
1569:   for (_k=0; _k<_nlnk; _k++){
1570:     indices[_k] = lnk[_next];
1571:     _next       = lnk[_next + 1];
1572:   }
1573:   lnk[0] = 0;          /* num of entries on the list */
1574:   lnk[3] = 2;          /* head node */
1575:   return 0;
1576: }

1578: PETSC_STATIC_INLINE PetscErrorCode PetscLLCondensedDestroy_Scalable(PetscInt *lnk)
1579: {
1580:   return PetscFree(lnk);
1581: }

1583: /* -------------------------------------------------------------------------------------------------------*/
1584: /*
1585:       lnk[0]   number of links
1586:       lnk[1]   number of entries
1587:       lnk[3n]  value
1588:       lnk[3n+1] len
1589:       lnk[3n+2] link to next value

1591:       The next three are always the first link

1593:       lnk[3]    PETSC_MIN_INT+1
1594:       lnk[4]    1
1595:       lnk[5]    link to first real entry

1597:       The next three are always the last link

1599:       lnk[6]    PETSC_MAX_INT - 1
1600:       lnk[7]    1
1601:       lnk[8]    next valid link (this is the same as lnk[0] but without the decreases)
1602: */

1604: PETSC_STATIC_INLINE PetscErrorCode PetscLLCondensedCreate_fast(PetscInt nlnk_max,PetscInt **lnk)
1605: {
1607:   PetscInt       *llnk,lsize = 0;

1610:   PetscIntMultError(3,nlnk_max+3,&lsize);
1611:   PetscMalloc1(lsize,lnk);
1612:   llnk = *lnk;
1613:   llnk[0] = 0;   /* nlnk: number of entries on the list */
1614:   llnk[1] = 0;          /* number of integer entries represented in list */
1615:   llnk[3] = PETSC_MIN_INT+1;   /* value in the first node */
1616:   llnk[4] = 1;           /* count for the first node */
1617:   llnk[5] = 6;         /* next for the first node */
1618:   llnk[6] = PETSC_MAX_INT-1;   /* value in the last node */
1619:   llnk[7] = 1;           /* count for the last node */
1620:   llnk[8] = 0;         /* next valid node to be used */
1621:   return(0);
1622: }

1624: PETSC_STATIC_INLINE PetscErrorCode PetscLLCondensedAddSorted_fast(PetscInt nidx,const PetscInt indices[],PetscInt lnk[])
1625: {
1626:   PetscInt k,entry,prev,next;
1627:   prev      = 3;      /* first value */
1628:   next      = lnk[prev+2];
1629:   for (k=0; k<nidx; k++){
1630:     entry = indices[k];
1631:     /* search for insertion location */
1632:     while (entry >= lnk[next]) {
1633:       prev = next;
1634:       next = lnk[next+2];
1635:     }
1636:     /* entry is in range of previous list */
1637:     if (entry < lnk[prev]+lnk[prev+1]) continue;
1638:     lnk[1]++;
1639:     /* entry is right after previous list */
1640:     if (entry == lnk[prev]+lnk[prev+1]) {
1641:       lnk[prev+1]++;
1642:       if (lnk[next] == entry+1) { /* combine two contiguous strings */
1643:         lnk[prev+1] += lnk[next+1];
1644:         lnk[prev+2]  = lnk[next+2];
1645:         next         = lnk[next+2];
1646:         lnk[0]--;
1647:       }
1648:       continue;
1649:     }
1650:     /* entry is right before next list */
1651:     if (entry == lnk[next]-1) {
1652:       lnk[next]--;
1653:       lnk[next+1]++;
1654:       prev = next;
1655:       next = lnk[prev+2];
1656:       continue;
1657:     }
1658:     /*  add entry into lnk */
1659:     lnk[prev+2]    = 3*((lnk[8]++)+3);      /* connect previous node to the new node */
1660:     prev           = lnk[prev+2];
1661:     lnk[prev]      = entry;        /* set value of the new node */
1662:     lnk[prev+1]    = 1;             /* number of values in contiguous string is one to start */
1663:     lnk[prev+2]    = next;          /* connect new node to next node */
1664:     lnk[0]++;
1665:   }
1666:   return 0;
1667: }

1669: PETSC_STATIC_INLINE PetscErrorCode PetscLLCondensedClean_fast(PetscInt nidx,PetscInt *indices,PetscInt *lnk)
1670: {
1671:   PetscInt _k,_next,_nlnk,cnt,j;
1672:   _next = lnk[5];       /* first node */
1673:   _nlnk = lnk[0];
1674:   cnt   = 0;
1675:   for (_k=0; _k<_nlnk; _k++){
1676:     for (j=0; j<lnk[_next+1]; j++) {
1677:       indices[cnt++] = lnk[_next] + j;
1678:     }
1679:     _next       = lnk[_next + 2];
1680:   }
1681:   lnk[0] = 0;   /* nlnk: number of links */
1682:   lnk[1] = 0;          /* number of integer entries represented in list */
1683:   lnk[3] = PETSC_MIN_INT+1;   /* value in the first node */
1684:   lnk[4] = 1;           /* count for the first node */
1685:   lnk[5] = 6;         /* next for the first node */
1686:   lnk[6] = PETSC_MAX_INT-1;   /* value in the last node */
1687:   lnk[7] = 1;           /* count for the last node */
1688:   lnk[8] = 0;         /* next valid location to make link */
1689:   return 0;
1690: }

1692: PETSC_STATIC_INLINE PetscErrorCode PetscLLCondensedView_fast(PetscInt *lnk)
1693: {
1694:   PetscInt k,next,nlnk;
1695:   next = lnk[5];       /* first node */
1696:   nlnk = lnk[0];
1697:   for (k=0; k<nlnk; k++){
1698: #if 0                           /* Debugging code */
1699:     printf("%d value %d len %d next %d\n",next,lnk[next],lnk[next+1],lnk[next+2]);
1700: #endif
1701:     next = lnk[next + 2];
1702:   }
1703:   return 0;
1704: }

1706: PETSC_STATIC_INLINE PetscErrorCode PetscLLCondensedDestroy_fast(PetscInt *lnk)
1707: {
1708:   return PetscFree(lnk);
1709: }

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

1714: PETSC_EXTERN PetscLogEvent MAT_Mult;
1715: PETSC_EXTERN PetscLogEvent MAT_MultMatrixFree;
1716: PETSC_EXTERN PetscLogEvent MAT_Mults;
1717: PETSC_EXTERN PetscLogEvent MAT_MultConstrained;
1718: PETSC_EXTERN PetscLogEvent MAT_MultAdd;
1719: PETSC_EXTERN PetscLogEvent MAT_MultTranspose;
1720: PETSC_EXTERN PetscLogEvent MAT_MultTransposeConstrained;
1721: PETSC_EXTERN PetscLogEvent MAT_MultTransposeAdd;
1722: PETSC_EXTERN PetscLogEvent MAT_Solve;
1723: PETSC_EXTERN PetscLogEvent MAT_Solves;
1724: PETSC_EXTERN PetscLogEvent MAT_SolveAdd;
1725: PETSC_EXTERN PetscLogEvent MAT_SolveTranspose;
1726: PETSC_EXTERN PetscLogEvent MAT_SolveTransposeAdd;
1727: PETSC_EXTERN PetscLogEvent MAT_SOR;
1728: PETSC_EXTERN PetscLogEvent MAT_ForwardSolve;
1729: PETSC_EXTERN PetscLogEvent MAT_BackwardSolve;
1730: PETSC_EXTERN PetscLogEvent MAT_LUFactor;
1731: PETSC_EXTERN PetscLogEvent MAT_LUFactorSymbolic;
1732: PETSC_EXTERN PetscLogEvent MAT_LUFactorNumeric;
1733: PETSC_EXTERN PetscLogEvent MAT_CholeskyFactor;
1734: PETSC_EXTERN PetscLogEvent MAT_CholeskyFactorSymbolic;
1735: PETSC_EXTERN PetscLogEvent MAT_CholeskyFactorNumeric;
1736: PETSC_EXTERN PetscLogEvent MAT_ILUFactor;
1737: PETSC_EXTERN PetscLogEvent MAT_ILUFactorSymbolic;
1738: PETSC_EXTERN PetscLogEvent MAT_ICCFactorSymbolic;
1739: PETSC_EXTERN PetscLogEvent MAT_Copy;
1740: PETSC_EXTERN PetscLogEvent MAT_Convert;
1741: PETSC_EXTERN PetscLogEvent MAT_Scale;
1742: PETSC_EXTERN PetscLogEvent MAT_AssemblyBegin;
1743: PETSC_EXTERN PetscLogEvent MAT_AssemblyEnd;
1744: PETSC_EXTERN PetscLogEvent MAT_SetValues;
1745: PETSC_EXTERN PetscLogEvent MAT_GetValues;
1746: PETSC_EXTERN PetscLogEvent MAT_GetRow;
1747: PETSC_EXTERN PetscLogEvent MAT_GetRowIJ;
1748: PETSC_EXTERN PetscLogEvent MAT_CreateSubMats;
1749: PETSC_EXTERN PetscLogEvent MAT_GetColoring;
1750: PETSC_EXTERN PetscLogEvent MAT_GetOrdering;
1751: PETSC_EXTERN PetscLogEvent MAT_RedundantMat;
1752: PETSC_EXTERN PetscLogEvent MAT_IncreaseOverlap;
1753: PETSC_EXTERN PetscLogEvent MAT_Partitioning;
1754: PETSC_EXTERN PetscLogEvent MAT_PartitioningND;
1755: PETSC_EXTERN PetscLogEvent MAT_Coarsen;
1756: PETSC_EXTERN PetscLogEvent MAT_ZeroEntries;
1757: PETSC_EXTERN PetscLogEvent MAT_Load;
1758: PETSC_EXTERN PetscLogEvent MAT_View;
1759: PETSC_EXTERN PetscLogEvent MAT_AXPY;
1760: PETSC_EXTERN PetscLogEvent MAT_FDColoringCreate;
1761: PETSC_EXTERN PetscLogEvent MAT_TransposeColoringCreate;
1762: PETSC_EXTERN PetscLogEvent MAT_FDColoringSetUp;
1763: PETSC_EXTERN PetscLogEvent MAT_FDColoringApply;
1764: PETSC_EXTERN PetscLogEvent MAT_Transpose;
1765: PETSC_EXTERN PetscLogEvent MAT_FDColoringFunction;
1766: PETSC_EXTERN PetscLogEvent MAT_CreateSubMat;
1767: PETSC_EXTERN PetscLogEvent MAT_MatSolve;
1768: PETSC_EXTERN PetscLogEvent MAT_MatTrSolve;
1769: PETSC_EXTERN PetscLogEvent MAT_MatMultSymbolic;
1770: PETSC_EXTERN PetscLogEvent MAT_MatMultNumeric;
1771: PETSC_EXTERN PetscLogEvent MAT_Getlocalmatcondensed;
1772: PETSC_EXTERN PetscLogEvent MAT_GetBrowsOfAcols;
1773: PETSC_EXTERN PetscLogEvent MAT_GetBrowsOfAocols;
1774: PETSC_EXTERN PetscLogEvent MAT_PtAPSymbolic;
1775: PETSC_EXTERN PetscLogEvent MAT_PtAPNumeric;
1776: PETSC_EXTERN PetscLogEvent MAT_Seqstompinum;
1777: PETSC_EXTERN PetscLogEvent MAT_Seqstompisym;
1778: PETSC_EXTERN PetscLogEvent MAT_Seqstompi;
1779: PETSC_EXTERN PetscLogEvent MAT_Getlocalmat;
1780: PETSC_EXTERN PetscLogEvent MAT_RARtSymbolic;
1781: PETSC_EXTERN PetscLogEvent MAT_RARtNumeric;
1782: PETSC_EXTERN PetscLogEvent MAT_MatTransposeMultSymbolic;
1783: PETSC_EXTERN PetscLogEvent MAT_MatTransposeMultNumeric;
1784: PETSC_EXTERN PetscLogEvent MAT_TransposeMatMultSymbolic;
1785: PETSC_EXTERN PetscLogEvent MAT_TransposeMatMultNumeric;
1786: PETSC_EXTERN PetscLogEvent MAT_MatMatMultSymbolic;
1787: PETSC_EXTERN PetscLogEvent MAT_MatMatMultNumeric;
1788: PETSC_EXTERN PetscLogEvent MAT_Applypapt;
1789: PETSC_EXTERN PetscLogEvent MAT_Applypapt_symbolic;
1790: PETSC_EXTERN PetscLogEvent MAT_Applypapt_numeric;
1791: PETSC_EXTERN PetscLogEvent MAT_Getsymtranspose;
1792: PETSC_EXTERN PetscLogEvent MAT_Getsymtransreduced;
1793: PETSC_EXTERN PetscLogEvent MAT_GetSequentialNonzeroStructure;
1794: PETSC_EXTERN PetscLogEvent MATMFFD_Mult;
1795: PETSC_EXTERN PetscLogEvent MAT_GetMultiProcBlock;
1796: PETSC_EXTERN PetscLogEvent MAT_CUSPARSECopyToGPU;
1797: PETSC_EXTERN PetscLogEvent MAT_CUSPARSECopyFromGPU;
1798: PETSC_EXTERN PetscLogEvent MAT_CUSPARSEGenerateTranspose;
1799: PETSC_EXTERN PetscLogEvent MAT_CUSPARSESolveAnalysis;
1800: PETSC_EXTERN PetscLogEvent MAT_SetValuesBatch;
1801: PETSC_EXTERN PetscLogEvent MAT_ViennaCLCopyToGPU;
1802: PETSC_EXTERN PetscLogEvent MAT_DenseCopyToGPU;
1803: PETSC_EXTERN PetscLogEvent MAT_DenseCopyFromGPU;
1804: PETSC_EXTERN PetscLogEvent MAT_Merge;
1805: PETSC_EXTERN PetscLogEvent MAT_Residual;
1806: PETSC_EXTERN PetscLogEvent MAT_SetRandom;
1807: PETSC_EXTERN PetscLogEvent MAT_FactorFactS;
1808: PETSC_EXTERN PetscLogEvent MAT_FactorInvS;
1809: PETSC_EXTERN PetscLogEvent MAT_PreallCOO;
1810: PETSC_EXTERN PetscLogEvent MAT_SetVCOO;
1811: PETSC_EXTERN PetscLogEvent MATCOLORING_Apply;
1812: PETSC_EXTERN PetscLogEvent MATCOLORING_Comm;
1813: PETSC_EXTERN PetscLogEvent MATCOLORING_Local;
1814: PETSC_EXTERN PetscLogEvent MATCOLORING_ISCreate;
1815: PETSC_EXTERN PetscLogEvent MATCOLORING_SetUp;
1816: PETSC_EXTERN PetscLogEvent MATCOLORING_Weights;

1818: #endif