Actual source code: matimpl.h

petsc-3.14.0 2020-09-29
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

242: PETSC_INTERN PetscErrorCode MatProductSymbolic_AB(Mat);
243: PETSC_INTERN PetscErrorCode MatProductNumeric_AB(Mat);
244: PETSC_INTERN PetscErrorCode MatProductSymbolic_AtB(Mat);
245: PETSC_INTERN PetscErrorCode MatProductNumeric_AtB(Mat);
246: PETSC_INTERN PetscErrorCode MatProductSymbolic_ABt(Mat);
247: PETSC_INTERN PetscErrorCode MatProductNumeric_ABt(Mat);
248: PETSC_INTERN PetscErrorCode MatProductNumeric_PtAP(Mat);
249: PETSC_INTERN PetscErrorCode MatProductNumeric_RARt(Mat);
250: PETSC_INTERN PetscErrorCode MatProductSymbolic_ABC(Mat);
251: PETSC_INTERN PetscErrorCode MatProductNumeric_ABC(Mat);
252: PETSC_INTERN PetscErrorCode MatProductCreate_Private(Mat,Mat,Mat,Mat);

254: #if defined(PETSC_USE_DEBUG)
255: #  define MatCheckPreallocated(A,arg) do {                              \
256:     if (PetscUnlikely(!(A)->preallocated)) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call MatXXXSetPreallocation() or MatSetUp() on argument %D \"%s\" before %s()",(arg),#A,PETSC_FUNCTION_NAME); \
257:   } while (0)
258: #else
259: #  define MatCheckPreallocated(A,arg) do {} while (0)
260: #endif

262: #if defined(PETSC_USE_DEBUG)
263: #  define MatCheckProduct(A,arg) do {                              \
264:     if (PetscUnlikely(!(A)->product)) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Argument %D \"%s\" is not a matrix obtained from MatProductCreate()",(arg),#A); \
265:   } while (0)
266: #else
267: #  define MatCheckProduct(A,arg) do {} while (0)
268: #endif

270: /*
271:   The stash is used to temporarily store inserted matrix values that
272:   belong to another processor. During the assembly phase the stashed
273:   values are moved to the correct processor and
274: */

276: typedef struct _MatStashSpace *PetscMatStashSpace;

278: struct _MatStashSpace {
279:   PetscMatStashSpace next;
280:   PetscScalar        *space_head,*val;
281:   PetscInt           *idx,*idy;
282:   PetscInt           total_space_size;
283:   PetscInt           local_used;
284:   PetscInt           local_remaining;
285: };

287: PETSC_EXTERN PetscErrorCode PetscMatStashSpaceGet(PetscInt,PetscInt,PetscMatStashSpace *);
288: PETSC_EXTERN PetscErrorCode PetscMatStashSpaceContiguous(PetscInt,PetscMatStashSpace *,PetscScalar *,PetscInt *,PetscInt *);
289: PETSC_EXTERN PetscErrorCode PetscMatStashSpaceDestroy(PetscMatStashSpace*);

291: typedef struct {
292:   PetscInt    count;
293: } MatStashHeader;

295: typedef struct {
296:   void        *buffer;          /* Of type blocktype, dynamically constructed  */
297:   PetscInt    count;
298:   char        pending;
299: } MatStashFrame;

301: typedef struct _MatStash MatStash;
302: struct _MatStash {
303:   PetscInt      nmax;                   /* maximum stash size */
304:   PetscInt      umax;                   /* user specified max-size */
305:   PetscInt      oldnmax;                /* the nmax value used previously */
306:   PetscInt      n;                      /* stash size */
307:   PetscInt      bs;                     /* block size of the stash */
308:   PetscInt      reallocs;               /* preserve the no of mallocs invoked */
309:   PetscMatStashSpace space_head,space;  /* linked list to hold stashed global row/column numbers and matrix values */

311:   PetscErrorCode (*ScatterBegin)(Mat,MatStash*,PetscInt*);
312:   PetscErrorCode (*ScatterGetMesg)(MatStash*,PetscMPIInt*,PetscInt**,PetscInt**,PetscScalar**,PetscInt*);
313:   PetscErrorCode (*ScatterEnd)(MatStash*);
314:   PetscErrorCode (*ScatterDestroy)(MatStash*);

316:   /* The following variables are used for communication */
317:   MPI_Comm      comm;
318:   PetscMPIInt   size,rank;
319:   PetscMPIInt   tag1,tag2;
320:   MPI_Request   *send_waits;            /* array of send requests */
321:   MPI_Request   *recv_waits;            /* array of receive requests */
322:   MPI_Status    *send_status;           /* array of send status */
323:   PetscInt      nsends,nrecvs;          /* numbers of sends and receives */
324:   PetscScalar   *svalues;               /* sending data */
325:   PetscInt      *sindices;
326:   PetscScalar   **rvalues;              /* receiving data (values) */
327:   PetscInt      **rindices;             /* receiving data (indices) */
328:   PetscInt      nprocessed;             /* number of messages already processed */
329:   PetscMPIInt   *flg_v;                 /* indicates what messages have arrived so far and from whom */
330:   PetscBool     reproduce;
331:   PetscInt      reproduce_count;

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

361: #if !defined(PETSC_HAVE_MPIUNI)
362: PETSC_INTERN PetscErrorCode MatStashScatterDestroy_BTS(MatStash*);
363: #endif
364: PETSC_INTERN PetscErrorCode MatStashCreate_Private(MPI_Comm,PetscInt,MatStash*);
365: PETSC_INTERN PetscErrorCode MatStashDestroy_Private(MatStash*);
366: PETSC_INTERN PetscErrorCode MatStashScatterEnd_Private(MatStash*);
367: PETSC_INTERN PetscErrorCode MatStashSetInitialSize_Private(MatStash*,PetscInt);
368: PETSC_INTERN PetscErrorCode MatStashGetInfo_Private(MatStash*,PetscInt*,PetscInt*);
369: PETSC_INTERN PetscErrorCode MatStashValuesRow_Private(MatStash*,PetscInt,PetscInt,const PetscInt[],const PetscScalar[],PetscBool);
370: PETSC_INTERN PetscErrorCode MatStashValuesCol_Private(MatStash*,PetscInt,PetscInt,const PetscInt[],const PetscScalar[],PetscInt,PetscBool);
371: PETSC_INTERN PetscErrorCode MatStashValuesRowBlocked_Private(MatStash*,PetscInt,PetscInt,const PetscInt[],const PetscScalar[],PetscInt,PetscInt,PetscInt);
372: PETSC_INTERN PetscErrorCode MatStashValuesColBlocked_Private(MatStash*,PetscInt,PetscInt,const PetscInt[],const PetscScalar[],PetscInt,PetscInt,PetscInt);
373: PETSC_INTERN PetscErrorCode MatStashScatterBegin_Private(Mat,MatStash*,PetscInt*);
374: PETSC_INTERN PetscErrorCode MatStashScatterGetMesg_Private(MatStash*,PetscMPIInt*,PetscInt**,PetscInt**,PetscScalar**,PetscInt*);
375: PETSC_INTERN PetscErrorCode MatGetInfo_External(Mat,MatInfoType,MatInfo*);

377: typedef struct {
378:   PetscInt   dim;
379:   PetscInt   dims[4];
380:   PetscInt   starts[4];
381:   PetscBool  noc;        /* this is a single component problem, hence user will not set MatStencil.c */
382: } MatStencilInfo;

384: /* Info about using compressed row format */
385: typedef struct {
386:   PetscBool  use;                           /* indicates compressed rows have been checked and will be used */
387:   PetscInt   nrows;                         /* number of non-zero rows */
388:   PetscInt   *i;                            /* compressed row pointer  */
389:   PetscInt   *rindex;                       /* compressed row index               */
390: } Mat_CompressedRow;
391: PETSC_EXTERN PetscErrorCode MatCheckCompressedRow(Mat,PetscInt,Mat_CompressedRow*,PetscInt*,PetscInt,PetscReal);

393: typedef struct { /* used by MatCreateRedundantMatrix() for reusing matredundant */
394:   PetscInt     nzlocal,nsends,nrecvs;
395:   PetscMPIInt  *send_rank,*recv_rank;
396:   PetscInt     *sbuf_nz,*rbuf_nz,*sbuf_j,**rbuf_j;
397:   PetscScalar  *sbuf_a,**rbuf_a;
398:   MPI_Comm     subcomm;   /* when user does not provide a subcomm */
399:   IS           isrow,iscol;
400:   Mat          *matseq;
401: } Mat_Redundant;

403: typedef struct { /* used by MatProduct() */
404:   MatProductType type;
405:   char           *alg;
406:   Mat            A,B,C,Dwork;
407:   PetscReal      fill;
408:   PetscBool      api_user; /* used by MatProductSetFromOptions_xxx() to distinguish command line options */

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

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

419: #define CSRDataStructure(datatype)  \
420:   int         *i; \
421:   int         *j; \
422:   datatype    *a;\
423:   PetscInt    n;\
424:   PetscInt    ignorezeroentries;

426: typedef struct {
427:   CSRDataStructure(PetscScalar)
428: } PetscCSRDataStructure;

430: struct _p_SplitCSRMat {
431:   PetscInt              cstart,cend,rstart,rend;
432:   PetscCSRDataStructure diag,offdiag;
433:   PetscInt              *colmap;
434:   PetscBool             seq;
435:   PetscMPIInt           rank;
436:   PetscInt              nonzerostate;
437: };

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

487: PETSC_INTERN PetscErrorCode MatAXPY_Basic(Mat,PetscScalar,Mat,MatStructure);
488: PETSC_INTERN PetscErrorCode MatAXPY_BasicWithPreallocation(Mat,Mat,PetscScalar,Mat,MatStructure);
489: PETSC_INTERN PetscErrorCode MatAXPY_Basic_Preallocate(Mat,Mat,Mat*);

491: /*
492:     Utility for MatFactor (Schur complement)
493: */
494: PETSC_INTERN PetscErrorCode MatFactorFactorizeSchurComplement_Private(Mat);
495: PETSC_INTERN PetscErrorCode MatFactorInvertSchurComplement_Private(Mat);
496: PETSC_INTERN PetscErrorCode MatFactorUpdateSchurStatus_Private(Mat);
497: PETSC_INTERN PetscErrorCode MatFactorSetUpInPlaceSchur_Private(Mat);

499: /*
500:     Utility for MatZeroRows
501: */
502: PETSC_INTERN PetscErrorCode MatZeroRowsMapLocal_Private(Mat,PetscInt,const PetscInt*,PetscInt*,PetscInt**);

504: /*
505:     Utility for MatView/MatLoad
506: */
507: PETSC_INTERN PetscErrorCode MatView_Binary_BlockSizes(Mat,PetscViewer);
508: PETSC_INTERN PetscErrorCode MatLoad_Binary_BlockSizes(Mat,PetscViewer);


511: /*
512:     Object for partitioning graphs
513: */

515: typedef struct _MatPartitioningOps *MatPartitioningOps;
516: struct _MatPartitioningOps {
517:   PetscErrorCode (*apply)(MatPartitioning,IS*);
518:   PetscErrorCode (*applynd)(MatPartitioning,IS*);
519:   PetscErrorCode (*setfromoptions)(PetscOptionItems*,MatPartitioning);
520:   PetscErrorCode (*destroy)(MatPartitioning);
521:   PetscErrorCode (*view)(MatPartitioning,PetscViewer);
522:   PetscErrorCode (*improve)(MatPartitioning,IS*);
523: };

525: struct _p_MatPartitioning {
526:   PETSCHEADER(struct _MatPartitioningOps);
527:   Mat         adj;
528:   PetscInt    *vertex_weights;
529:   PetscReal   *part_weights;
530:   PetscInt    n;                                 /* number of partitions */
531:   void        *data;
532:   PetscInt    setupcalled;
533:   PetscBool   use_edge_weights;  /* A flag indicates whether or not to use edge weights */
534: };

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

539: /*
540:     Object for coarsen graphs
541: */
542: typedef struct _MatCoarsenOps *MatCoarsenOps;
543: struct _MatCoarsenOps {
544:   PetscErrorCode (*apply)(MatCoarsen);
545:   PetscErrorCode (*setfromoptions)(PetscOptionItems*,MatCoarsen);
546:   PetscErrorCode (*destroy)(MatCoarsen);
547:   PetscErrorCode (*view)(MatCoarsen,PetscViewer);
548: };

550: struct _p_MatCoarsen {
551:   PETSCHEADER(struct _MatCoarsenOps);
552:   Mat              graph;
553:   PetscInt         setupcalled;
554:   void             *subctx;
555:   /* */
556:   PetscBool        strict_aggs;
557:   IS               perm;
558:   PetscCoarsenData *agg_lists;
559: };

561: /*
562:     MatFDColoring is used to compute Jacobian matrices efficiently
563:   via coloring. The data structure is explained below in an example.

565:    Color =   0    1     0    2   |   2      3       0
566:    ---------------------------------------------------
567:             00   01              |          05
568:             10   11              |   14     15               Processor  0
569:                        22    23  |          25
570:                        32    33  |
571:    ===================================================
572:                                  |   44     45     46
573:             50                   |          55               Processor 1
574:                                  |   64            66
575:    ---------------------------------------------------

577:     ncolors = 4;

579:     ncolumns      = {2,1,1,0}
580:     columns       = {{0,2},{1},{3},{}}
581:     nrows         = {4,2,3,3}
582:     rows          = {{0,1,2,3},{0,1},{1,2,3},{0,1,2}}
583:     vwscale       = {dx(0),dx(1),dx(2),dx(3)}               MPI Vec
584:     vscale        = {dx(0),dx(1),dx(2),dx(3),dx(4),dx(5)}   Seq Vec

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

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

596: */
597: typedef struct {
598:   PetscInt     row;
599:   PetscInt     col;
600:   PetscScalar  *valaddr;   /* address of value */
601: } MatEntry;

603: typedef struct {
604:   PetscInt     row;
605:   PetscScalar  *valaddr;   /* address of value */
606: } MatEntry2;

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

637: typedef struct _MatColoringOps *MatColoringOps;
638: struct _MatColoringOps {
639:   PetscErrorCode (*destroy)(MatColoring);
640:   PetscErrorCode (*setfromoptions)(PetscOptionItems*,MatColoring);
641:   PetscErrorCode (*view)(MatColoring,PetscViewer);
642:   PetscErrorCode (*apply)(MatColoring,ISColoring*);
643:   PetscErrorCode (*weights)(MatColoring,PetscReal**,PetscInt**);
644: };

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

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

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

677: /*
678:    Null space context for preconditioner/operators
679: */
680: struct _p_MatNullSpace {
681:   PETSCHEADER(int);
682:   PetscBool      has_cnst;
683:   PetscInt       n;
684:   Vec*           vecs;
685:   PetscScalar*   alpha;                 /* for projections */
686:   PetscErrorCode (*remove)(MatNullSpace,Vec,void*);  /* for user provided removal function */
687:   void*          rmctx;                 /* context for remove() function */
688: };

690: /*
691:    Checking zero pivot for LU, ILU preconditioners.
692: */
693: typedef struct {
694:   PetscInt       nshift,nshift_max;
695:   PetscReal      shift_amount,shift_lo,shift_hi,shift_top,shift_fraction;
696:   PetscBool      newshift;
697:   PetscReal      rs;  /* active row sum of abs(offdiagonals) */
698:   PetscScalar    pv;  /* pivot of the active row */
699: } FactorShiftCtx;

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

725:   PetscErrorCode (*destroy)(Mat);
726: } Mat_SubSppt;

728: PETSC_EXTERN PetscErrorCode MatFactorDumpMatrix(Mat);
729: PETSC_INTERN PetscErrorCode MatShift_Basic(Mat,PetscScalar);
730: PETSC_INTERN PetscErrorCode MatSetBlockSizes_Default(Mat,PetscInt,PetscInt);

732: PETSC_STATIC_INLINE PetscErrorCode MatPivotCheck_nz(Mat mat,const MatFactorInfo *info,FactorShiftCtx *sctx,PetscInt row)
733: {
734:   PetscReal _rs   = sctx->rs;
735:   PetscReal _zero = info->zeropivot*_rs;

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

750: PETSC_STATIC_INLINE PetscErrorCode MatPivotCheck_pd(Mat mat,const MatFactorInfo *info,FactorShiftCtx *sctx,PetscInt row)
751: {
752:   PetscReal _rs   = sctx->rs;
753:   PetscReal _zero = info->zeropivot*_rs;

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

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

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

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

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

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

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

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

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

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

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

949: #define PetscLLAddSorted_new(nidx,indices,idx_start,lnk_empty,nlnk,lnk,bt) 0; \
950: {\
951:   PetscInt _k,_entry,_location,_lnkdata;\
952:   if (lnk_empty){\
953:     _lnkdata  = idx_start;                      \
954:     for (_k=0; _k<nidx; _k++){                  \
955:       _entry = indices[_k];                             \
956:       PetscBTSet(bt,_entry);  /* mark the new entry */          \
957:           _location = _lnkdata;                                 \
958:           _lnkdata  = lnk[_location];                           \
959:         /* insertion location is found, add entry into lnk */   \
960:         lnk[_location] = _entry;                                \
961:         lnk[_entry]    = _lnkdata;                              \
962:         _lnkdata = _entry; /* next search starts from here */   \
963:     }                                                           \
964:     /*\
965:     lnk[indices[nidx-1]] = lnk[idx_start];\
966:     lnk[idx_start]       = indices[0];\
967:     PetscBTSet(bt,indices[0]);  \
968:     for (_k=1; _k<nidx; _k++){                  \
969:       PetscBTSet(bt,indices[_k]);                                          \
970:       lnk[indices[_k-1]] = indices[_k];                                  \
971:     }                                                           \
972:      */\
973:     nlnk      = nidx;\
974:     lnk_empty = PETSC_FALSE;\
975:   } else {\
976:     nlnk      = 0;                              \
977:     _lnkdata  = idx_start;                      \
978:     for (_k=0; _k<nidx; _k++){                  \
979:       _entry = indices[_k];                             \
980:       if (!PetscBTLookupSet(bt,_entry)){  /* new entry */       \
981:         /* search for insertion location */                     \
982:         do {                                                    \
983:           _location = _lnkdata;                                 \
984:           _lnkdata  = lnk[_location];                           \
985:         } while (_entry > _lnkdata);                            \
986:         /* insertion location is found, add entry into lnk */   \
987:         lnk[_location] = _entry;                                \
988:         lnk[_entry]    = _lnkdata;                              \
989:         nlnk++;                                                 \
990:         _lnkdata = _entry; /* next search starts from here */   \
991:       }                                                         \
992:     }                                                           \
993:   }                                                             \
994: }

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

1040: /*
1041:   Copy data on the list into an array, then initialize the list
1042:   Input Parameters:
1043:     idx_start - starting index of the list
1044:     lnk_max   - max value of lnk indicating the end of the list
1045:     nlnk      - number of data on the list to be copied
1046:     lnk       - linked list
1047:     bt        - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
1048:   output Parameters:
1049:     indices   - array that contains the copied data
1050:     lnk       - linked list that is cleaned and initialize
1051:     bt        - PetscBT (bitarray) with all bits set to false
1052: */
1053: #define PetscLLClean(idx_start,lnk_max,nlnk,lnk,indices,bt) 0;\
1054: {\
1055:   PetscInt _j,_idx=idx_start;\
1056:   for (_j=0; _j<nlnk; _j++){\
1057:     _idx = lnk[_idx];\
1058:     indices[_j] = _idx;\
1059:     PetscBTClear(bt,_idx);\
1060:   }\
1061:   lnk[idx_start] = lnk_max;\
1062: }
1063: /*
1064:   Free memories used by the list
1065: */
1066: #define PetscLLDestroy(lnk,bt) (PetscFree(lnk) || PetscBTDestroy(&(bt)))

1068: /* Routines below are used for incomplete matrix factorization */
1069: /*
1070:   Create and initialize a linked list and its levels
1071:   Input Parameters:
1072:     idx_start - starting index of the list
1073:     lnk_max   - max value of lnk indicating the end of the list
1074:     nlnk      - max length of the list
1075:   Output Parameters:
1076:     lnk       - list initialized
1077:     lnk_lvl   - array of size nlnk for storing levels of lnk
1078:     bt        - PetscBT (bitarray) with all bits set to false
1079: */
1080: #define PetscIncompleteLLCreate(idx_start,lnk_max,nlnk,lnk,lnk_lvl,bt)\
1081:   (PetscIntMultError(2,nlnk,NULL) || PetscMalloc1(2*nlnk,&lnk) || PetscBTCreate(nlnk,&(bt)) || (lnk[idx_start] = lnk_max,lnk_lvl = lnk + nlnk,0))

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

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

1141:   Note: the level of factor(i,j) is set as lvl(i,j) = min{ lvl(i,j), lvl(i,prow)+lvl(prow,j)+1)
1142:         where idx = non-zero columns of U(prow,prow+1:n-1), prow<i
1143: */
1144: #define PetscILULLAddSorted(nidx,idx,level,idxlvl,idx_start,nlnk,lnk,lnklvl,bt,lnklvl_prow) 0;\
1145: {\
1146:   PetscInt _k,_entry,_location,_lnkdata,_incrlev,_lnklvl_prow=lnklvl[prow];\
1147:   nlnk     = 0;\
1148:   _lnkdata = idx_start;\
1149:   for (_k=0; _k<nidx; _k++){\
1150:     _incrlev = idxlvl[_k] + _lnklvl_prow + 1;\
1151:     if (_incrlev > level) continue;\
1152:     _entry = idx[_k];\
1153:     if (!PetscBTLookupSet(bt,_entry)){  /* new entry */\
1154:       /* search for insertion location */\
1155:       do {\
1156:         _location = _lnkdata;\
1157:         _lnkdata  = lnk[_location];\
1158:       } while (_entry > _lnkdata);\
1159:       /* insertion location is found, add entry into lnk */\
1160:       lnk[_location]  = _entry;\
1161:       lnk[_entry]     = _lnkdata;\
1162:       lnklvl[_entry] = _incrlev;\
1163:       nlnk++;\
1164:       _lnkdata = _entry; /* next search starts from here if next_entry > _entry */\
1165:     } else { /* existing entry: update lnklvl */\
1166:       if (lnklvl[_entry] > _incrlev) lnklvl[_entry] = _incrlev;\
1167:     }\
1168:   }\
1169: }

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

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

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

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

1339: #define MatCheckSameLocalSize(A,ar1,B,ar2) do { \
1341:   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)

1343: #define MatCheckSameSize(A,ar1,B,ar2) do { \
1344:   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);\
1345:   MatCheckSameLocalSize(A,ar1,B,ar2);} while (0)

1347: #define VecCheckMatCompatible(M,x,ar1,b,ar2) do { \
1348:   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); \
1349:   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)

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

1360:       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
1361:       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
1362:       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
1363:       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.
1364:       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
1365:       to each other so memory access is much better than using the big array.

1367:   Example:
1368:      nlnk_max=5, lnk_max=36:
1369:      Initial list: [0, 0 | 36, 2 | 0, 0 | 0, 0 | 0, 0 | 0, 0 | 0, 0]
1370:      here, head_node has index 2 with value lnk[2]=lnk_max=36,
1371:            0-th entry is used to store the number of entries in the list,
1372:      The initial lnk represents head -> tail(marked by 36) with number of entries = lnk[0]=0.

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

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

1382:   Input Parameters:
1383:     nlnk_max  - max length of the list
1384:     lnk_max   - max value of the entries
1385:   Output Parameters:
1386:     lnk       - list created and initialized
1387:     bt        - PetscBT (bitarray) with all bits set to false. Note: bt has size lnk_max, not nln_max!
1388: */
1389: PETSC_STATIC_INLINE PetscErrorCode PetscLLCondensedCreate(PetscInt nlnk_max,PetscInt lnk_max,PetscInt **lnk,PetscBT *bt)
1390: {
1392:   PetscInt       *llnk,lsize = 0;

1395:   PetscIntMultError(2,nlnk_max+2,&lsize);
1396:   PetscMalloc1(lsize,lnk);
1397:   PetscBTCreate(lnk_max,bt);
1398:   llnk = *lnk;
1399:   llnk[0] = 0;         /* number of entries on the list */
1400:   llnk[2] = lnk_max;   /* value in the head node */
1401:   llnk[3] = 2;         /* next for the head node */
1402:   return(0);
1403: }

1405: /*
1406:   Add a SORTED ascending index set into a sorted linked list. See PetscLLCondensedCreate() for detailed description.
1407:   Input Parameters:
1408:     nidx      - number of input indices
1409:     indices   - sorted integer array
1410:     lnk       - condensed linked list(an integer array) that is created
1411:     bt        - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
1412:   output Parameters:
1413:     lnk       - the sorted(increasing order) linked list containing previous and newly added non-redundate indices
1414:     bt        - updated PetscBT (bitarray)
1415: */
1416: PETSC_STATIC_INLINE PetscErrorCode PetscLLCondensedAddSorted(PetscInt nidx,const PetscInt indices[],PetscInt lnk[],PetscBT bt)
1417: {
1418:   PetscInt _k,_entry,_location,_next,_lnkdata,_nlnk,_newnode;

1421:   _nlnk     = lnk[0]; /* num of entries on the input lnk */
1422:   _location = 2; /* head */
1423:     for (_k=0; _k<nidx; _k++){
1424:       _entry = indices[_k];
1425:       if (!PetscBTLookupSet(bt,_entry)){  /* new entry */
1426:         /* search for insertion location */
1427:         do {
1428:           _next     = _location + 1; /* link from previous node to next node */
1429:           _location = lnk[_next];    /* idx of next node */
1430:           _lnkdata  = lnk[_location];/* value of next node */
1431:         } while (_entry > _lnkdata);
1432:         /* insertion location is found, add entry into lnk */
1433:         _newnode        = 2*(_nlnk+2);   /* index for this new node */
1434:         lnk[_next]      = _newnode;      /* connect previous node to the new node */
1435:         lnk[_newnode]   = _entry;        /* set value of the new node */
1436:         lnk[_newnode+1] = _location;     /* connect new node to next node */
1437:         _location       = _newnode;      /* next search starts from the new node */
1438:         _nlnk++;
1439:       }   \
1440:     }\
1441:   lnk[0]   = _nlnk;   /* number of entries in the list */
1442:   return(0);
1443: }

1445: PETSC_STATIC_INLINE PetscErrorCode PetscLLCondensedClean(PetscInt lnk_max,PetscInt nidx,PetscInt *indices,PetscInt lnk[],PetscBT bt)
1446: {
1448:   PetscInt       _k,_next,_nlnk;

1451:   _next = lnk[3];       /* head node */
1452:   _nlnk = lnk[0];       /* num of entries on the list */
1453:   for (_k=0; _k<_nlnk; _k++){
1454:     indices[_k] = lnk[_next];
1455:     _next       = lnk[_next + 1];
1456:     PetscBTClear(bt,indices[_k]);
1457:   }
1458:   lnk[0] = 0;          /* num of entries on the list */
1459:   lnk[2] = lnk_max;    /* initialize head node */
1460:   lnk[3] = 2;          /* head node */
1461:   return(0);
1462: }

1464: PETSC_STATIC_INLINE PetscErrorCode PetscLLCondensedView(PetscInt *lnk)
1465: {
1467:   PetscInt       k;

1470:   PetscPrintf(PETSC_COMM_SELF,"LLCondensed of size %D, (val,  next)\n",lnk[0]);
1471:   for (k=2; k< lnk[0]+2; k++){
1472:     PetscPrintf(PETSC_COMM_SELF," %D: (%D, %D)\n",2*k,lnk[2*k],lnk[2*k+1]);
1473:   }
1474:   return(0);
1475: }

1477: /*
1478:   Free memories used by the list
1479: */
1480: PETSC_STATIC_INLINE PetscErrorCode PetscLLCondensedDestroy(PetscInt *lnk,PetscBT bt)
1481: {

1485:   PetscFree(lnk);
1486:   PetscBTDestroy(&bt);
1487:   return(0);
1488: }

1490: /* -------------------------------------------------------------------------------------------------------*/
1491: /*
1492:  Same as PetscLLCondensedCreate(), but does not use non-scalable O(lnk_max) bitarray
1493:   Input Parameters:
1494:     nlnk_max  - max length of the list
1495:   Output Parameters:
1496:     lnk       - list created and initialized
1497: */
1498: PETSC_STATIC_INLINE PetscErrorCode PetscLLCondensedCreate_Scalable(PetscInt nlnk_max,PetscInt **lnk)
1499: {
1501:   PetscInt       *llnk,lsize = 0;

1504:   PetscIntMultError(2,nlnk_max+2,&lsize);
1505:   PetscMalloc1(lsize,lnk);
1506:   llnk = *lnk;
1507:   llnk[0] = 0;               /* number of entries on the list */
1508:   llnk[2] = PETSC_MAX_INT;   /* value in the head node */
1509:   llnk[3] = 2;               /* next for the head node */
1510:   return(0);
1511: }

1513: PETSC_STATIC_INLINE PetscErrorCode PetscLLCondensedExpand_Scalable(PetscInt nlnk_max,PetscInt **lnk)
1514: {
1516:   PetscInt       lsize = 0;

1519:   PetscIntMultError(2,nlnk_max+2,&lsize);
1520:   PetscRealloc(lsize*sizeof(PetscInt),lnk);
1521:   return(0);
1522: }

1524: PETSC_STATIC_INLINE PetscErrorCode PetscLLCondensedAddSorted_Scalable(PetscInt nidx,const PetscInt indices[],PetscInt lnk[])
1525: {
1526:   PetscInt _k,_entry,_location,_next,_lnkdata,_nlnk,_newnode;
1527:   _nlnk     = lnk[0]; /* num of entries on the input lnk */
1528:   _location = 2; /* head */ \
1529:     for (_k=0; _k<nidx; _k++){
1530:       _entry = indices[_k];
1531:       /* search for insertion location */
1532:       do {
1533:         _next     = _location + 1; /* link from previous node to next node */
1534:         _location = lnk[_next];    /* idx of next node */
1535:         _lnkdata  = lnk[_location];/* value of next node */
1536:       } while (_entry > _lnkdata);
1537:       if (_entry < _lnkdata) {
1538:         /* insertion location is found, add entry into lnk */
1539:         _newnode        = 2*(_nlnk+2);   /* index for this new node */
1540:         lnk[_next]      = _newnode;      /* connect previous node to the new node */
1541:         lnk[_newnode]   = _entry;        /* set value of the new node */
1542:         lnk[_newnode+1] = _location;     /* connect new node to next node */
1543:         _location       = _newnode;      /* next search starts from the new node */
1544:         _nlnk++;
1545:       }
1546:     }
1547:   lnk[0]   = _nlnk;   /* number of entries in the list */
1548:   return 0;
1549: }

1551: PETSC_STATIC_INLINE PetscErrorCode PetscLLCondensedClean_Scalable(PetscInt nidx,PetscInt *indices,PetscInt *lnk)
1552: {
1553:   PetscInt _k,_next,_nlnk;
1554:   _next = lnk[3];       /* head node */
1555:   _nlnk = lnk[0];
1556:   for (_k=0; _k<_nlnk; _k++){
1557:     indices[_k] = lnk[_next];
1558:     _next       = lnk[_next + 1];
1559:   }
1560:   lnk[0] = 0;          /* num of entries on the list */
1561:   lnk[3] = 2;          /* head node */
1562:   return 0;
1563: }

1565: PETSC_STATIC_INLINE PetscErrorCode PetscLLCondensedDestroy_Scalable(PetscInt *lnk)
1566: {
1567:   return PetscFree(lnk);
1568: }

1570: /* -------------------------------------------------------------------------------------------------------*/
1571: /*
1572:       lnk[0]   number of links
1573:       lnk[1]   number of entries
1574:       lnk[3n]  value
1575:       lnk[3n+1] len
1576:       lnk[3n+2] link to next value

1578:       The next three are always the first link

1580:       lnk[3]    PETSC_MIN_INT+1
1581:       lnk[4]    1
1582:       lnk[5]    link to first real entry

1584:       The next three are always the last link

1586:       lnk[6]    PETSC_MAX_INT - 1
1587:       lnk[7]    1
1588:       lnk[8]    next valid link (this is the same as lnk[0] but without the decreases)
1589: */

1591: PETSC_STATIC_INLINE PetscErrorCode PetscLLCondensedCreate_fast(PetscInt nlnk_max,PetscInt **lnk)
1592: {
1594:   PetscInt       *llnk,lsize = 0;

1597:   PetscIntMultError(3,nlnk_max+3,&lsize);
1598:   PetscMalloc1(lsize,lnk);
1599:   llnk = *lnk;
1600:   llnk[0] = 0;   /* nlnk: number of entries on the list */
1601:   llnk[1] = 0;          /* number of integer entries represented in list */
1602:   llnk[3] = PETSC_MIN_INT+1;   /* value in the first node */
1603:   llnk[4] = 1;           /* count for the first node */
1604:   llnk[5] = 6;         /* next for the first node */
1605:   llnk[6] = PETSC_MAX_INT-1;   /* value in the last node */
1606:   llnk[7] = 1;           /* count for the last node */
1607:   llnk[8] = 0;         /* next valid node to be used */
1608:   return(0);
1609: }

1611: PETSC_STATIC_INLINE PetscErrorCode PetscLLCondensedAddSorted_fast(PetscInt nidx,const PetscInt indices[],PetscInt lnk[])
1612: {
1613:   PetscInt k,entry,prev,next;
1614:   prev      = 3;      /* first value */
1615:   next      = lnk[prev+2];
1616:   for (k=0; k<nidx; k++){
1617:     entry = indices[k];
1618:     /* search for insertion location */
1619:     while (entry >= lnk[next]) {
1620:       prev = next;
1621:       next = lnk[next+2];
1622:     }
1623:     /* entry is in range of previous list */
1624:     if (entry < lnk[prev]+lnk[prev+1]) continue;
1625:     lnk[1]++;
1626:     /* entry is right after previous list */
1627:     if (entry == lnk[prev]+lnk[prev+1]) {
1628:       lnk[prev+1]++;
1629:       if (lnk[next] == entry+1) { /* combine two contiguous strings */
1630:         lnk[prev+1] += lnk[next+1];
1631:         lnk[prev+2]  = lnk[next+2];
1632:         next         = lnk[next+2];
1633:         lnk[0]--;
1634:       }
1635:       continue;
1636:     }
1637:     /* entry is right before next list */
1638:     if (entry == lnk[next]-1) {
1639:       lnk[next]--;
1640:       lnk[next+1]++;
1641:       prev = next;
1642:       next = lnk[prev+2];
1643:       continue;
1644:     }
1645:     /*  add entry into lnk */
1646:     lnk[prev+2]    = 3*((lnk[8]++)+3);      /* connect previous node to the new node */
1647:     prev           = lnk[prev+2];
1648:     lnk[prev]      = entry;        /* set value of the new node */
1649:     lnk[prev+1]    = 1;             /* number of values in contiguous string is one to start */
1650:     lnk[prev+2]    = next;          /* connect new node to next node */
1651:     lnk[0]++;
1652:   }
1653:   return 0;
1654: }

1656: PETSC_STATIC_INLINE PetscErrorCode PetscLLCondensedClean_fast(PetscInt nidx,PetscInt *indices,PetscInt *lnk)
1657: {
1658:   PetscInt _k,_next,_nlnk,cnt,j;
1659:   _next = lnk[5];       /* first node */
1660:   _nlnk = lnk[0];
1661:   cnt   = 0;
1662:   for (_k=0; _k<_nlnk; _k++){
1663:     for (j=0; j<lnk[_next+1]; j++) {
1664:       indices[cnt++] = lnk[_next] + j;
1665:     }
1666:     _next       = lnk[_next + 2];
1667:   }
1668:   lnk[0] = 0;   /* nlnk: number of links */
1669:   lnk[1] = 0;          /* number of integer entries represented in list */
1670:   lnk[3] = PETSC_MIN_INT+1;   /* value in the first node */
1671:   lnk[4] = 1;           /* count for the first node */
1672:   lnk[5] = 6;         /* next for the first node */
1673:   lnk[6] = PETSC_MAX_INT-1;   /* value in the last node */
1674:   lnk[7] = 1;           /* count for the last node */
1675:   lnk[8] = 0;         /* next valid location to make link */
1676:   return 0;
1677: }

1679: PETSC_STATIC_INLINE PetscErrorCode PetscLLCondensedView_fast(PetscInt *lnk)
1680: {
1681:   PetscInt k,next,nlnk;
1682:   next = lnk[5];       /* first node */
1683:   nlnk = lnk[0];
1684:   for (k=0; k<nlnk; k++){
1685: #if 0                           /* Debugging code */
1686:     printf("%d value %d len %d next %d\n",next,lnk[next],lnk[next+1],lnk[next+2]);
1687: #endif
1688:     next = lnk[next + 2];
1689:   }
1690:   return 0;
1691: }

1693: PETSC_STATIC_INLINE PetscErrorCode PetscLLCondensedDestroy_fast(PetscInt *lnk)
1694: {
1695:   return PetscFree(lnk);
1696: }

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

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

1801: #endif