Actual source code: matimpl.h

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

467: PETSC_INTERN PetscErrorCode MatAXPY_Basic(Mat,PetscScalar,Mat,MatStructure);
468: PETSC_INTERN PetscErrorCode MatAXPY_BasicWithPreallocation(Mat,Mat,PetscScalar,Mat,MatStructure);
469: PETSC_INTERN PetscErrorCode MatAXPY_Basic_Preallocate(Mat,Mat,Mat*);

471: /*
472:     Utility for MatFactor (Schur complement)
473: */
474: PETSC_INTERN PetscErrorCode MatFactorFactorizeSchurComplement_Private(Mat);
475: PETSC_INTERN PetscErrorCode MatFactorInvertSchurComplement_Private(Mat);
476: PETSC_INTERN PetscErrorCode MatFactorUpdateSchurStatus_Private(Mat);
477: PETSC_INTERN PetscErrorCode MatFactorSetUpInPlaceSchur_Private(Mat);

479: /*
480:     Utility for MatZeroRows
481: */
482: PETSC_INTERN PetscErrorCode MatZeroRowsMapLocal_Private(Mat,PetscInt,const PetscInt*,PetscInt*,PetscInt**);

484: /*
485:     Utility for MatView/MatLoad
486: */
487: PETSC_INTERN PetscErrorCode MatView_Binary_BlockSizes(Mat,PetscViewer);
488: PETSC_INTERN PetscErrorCode MatLoad_Binary_BlockSizes(Mat,PetscViewer);


491: /*
492:     Object for partitioning graphs
493: */

495: typedef struct _MatPartitioningOps *MatPartitioningOps;
496: struct _MatPartitioningOps {
497:   PetscErrorCode (*apply)(MatPartitioning,IS*);
498:   PetscErrorCode (*applynd)(MatPartitioning,IS*);
499:   PetscErrorCode (*setfromoptions)(PetscOptionItems*,MatPartitioning);
500:   PetscErrorCode (*destroy)(MatPartitioning);
501:   PetscErrorCode (*view)(MatPartitioning,PetscViewer);
502:   PetscErrorCode (*improve)(MatPartitioning,IS*);
503: };

505: struct _p_MatPartitioning {
506:   PETSCHEADER(struct _MatPartitioningOps);
507:   Mat         adj;
508:   PetscInt    *vertex_weights;
509:   PetscReal   *part_weights;
510:   PetscInt    n;                                 /* number of partitions */
511:   void        *data;
512:   PetscInt    setupcalled;
513:   PetscBool   use_edge_weights;  /* A flag indicates whether or not to use edge weights */
514: };

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

519: /*
520:     Object for coarsen graphs
521: */
522: typedef struct _MatCoarsenOps *MatCoarsenOps;
523: struct _MatCoarsenOps {
524:   PetscErrorCode (*apply)(MatCoarsen);
525:   PetscErrorCode (*setfromoptions)(PetscOptionItems*,MatCoarsen);
526:   PetscErrorCode (*destroy)(MatCoarsen);
527:   PetscErrorCode (*view)(MatCoarsen,PetscViewer);
528: };

530: struct _p_MatCoarsen {
531:   PETSCHEADER(struct _MatCoarsenOps);
532:   Mat              graph;
533:   PetscInt         setupcalled;
534:   void             *subctx;
535:   /* */
536:   PetscBool        strict_aggs;
537:   IS               perm;
538:   PetscCoarsenData *agg_lists;
539: };

541: /*
542:     MatFDColoring is used to compute Jacobian matrices efficiently
543:   via coloring. The data structure is explained below in an example.

545:    Color =   0    1     0    2   |   2      3       0
546:    ---------------------------------------------------
547:             00   01              |          05
548:             10   11              |   14     15               Processor  0
549:                        22    23  |          25
550:                        32    33  |
551:    ===================================================
552:                                  |   44     45     46
553:             50                   |          55               Processor 1
554:                                  |   64            66
555:    ---------------------------------------------------

557:     ncolors = 4;

559:     ncolumns      = {2,1,1,0}
560:     columns       = {{0,2},{1},{3},{}}
561:     nrows         = {4,2,3,3}
562:     rows          = {{0,1,2,3},{0,1},{1,2,3},{0,1,2}}
563:     vwscale       = {dx(0),dx(1),dx(2),dx(3)}               MPI Vec
564:     vscale        = {dx(0),dx(1),dx(2),dx(3),dx(4),dx(5)}   Seq Vec

566:     ncolumns      = {1,0,1,1}
567:     columns       = {{6},{},{4},{5}}
568:     nrows         = {3,0,2,2}
569:     rows          = {{0,1,2},{},{1,2},{1,2}}
570:     vwscale       = {dx(4),dx(5),dx(6)}              MPI Vec
571:     vscale        = {dx(0),dx(4),dx(5),dx(6)}        Seq Vec

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

576: */
577: typedef struct {
578:   PetscInt     row;
579:   PetscInt     col;
580:   PetscScalar  *valaddr;   /* address of value */
581: } MatEntry;

583: typedef struct {
584:   PetscInt     row;
585:   PetscScalar  *valaddr;   /* address of value */
586: } MatEntry2;

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

617: typedef struct _MatColoringOps *MatColoringOps;
618: struct _MatColoringOps {
619:   PetscErrorCode (*destroy)(MatColoring);
620:   PetscErrorCode (*setfromoptions)(PetscOptionItems*,MatColoring);
621:   PetscErrorCode (*view)(MatColoring,PetscViewer);
622:   PetscErrorCode (*apply)(MatColoring,ISColoring*);
623:   PetscErrorCode (*weights)(MatColoring,PetscReal**,PetscInt**);
624: };

626: struct _p_MatColoring {
627:   PETSCHEADER(struct _MatColoringOps);
628:   Mat                   mat;
629:   PetscInt              dist;             /* distance of the coloring */
630:   PetscInt              maxcolors;        /* the maximum number of colors returned, maxcolors=1 for MIS */
631:   void                  *data;            /* inner context */
632:   PetscBool             valid;            /* check to see if what is produced is a valid coloring */
633:   MatColoringWeightType weight_type;      /* type of weight computation to be performed */
634:   PetscReal             *user_weights;    /* custom weights and permutation */
635:   PetscInt              *user_lperm;
636:   PetscBool             valid_iscoloring; /* check to see if matcoloring is produced a valid iscoloring */
637: };

639: struct  _p_MatTransposeColoring{
640:   PETSCHEADER(int);
641:   PetscInt       M,N,m;            /* total rows, columns; local rows */
642:   PetscInt       rstart;           /* first row owned by local processor */
643:   PetscInt       ncolors;          /* number of colors */
644:   PetscInt       *ncolumns;        /* number of local columns for a color */
645:   PetscInt       *nrows;           /* number of local rows for each color */
646:   PetscInt       currentcolor;     /* color for which function evaluation is being done now */
647:   ISColoringType ctype;            /* IS_COLORING_GLOBAL or IS_COLORING_LOCAL */

649:   PetscInt       *colorforrow,*colorforcol;  /* pointer to rows and columns */
650:   PetscInt       *rows;                      /* lists the local rows for each color (using the local row numbering) */
651:   PetscInt       *den2sp;                    /* maps (row,color) in the dense matrix to index of sparse matrix array a->a */
652:   PetscInt       *columns;                   /* lists the local columns of each color (using global column numbering) */
653:   PetscInt       brows;                      /* number of rows for efficient implementation of MatTransColoringApplyDenToSp() */
654:   PetscInt       *lstart;                    /* array used for loop over row blocks of Csparse */
655: };

657: /*
658:    Null space context for preconditioner/operators
659: */
660: struct _p_MatNullSpace {
661:   PETSCHEADER(int);
662:   PetscBool      has_cnst;
663:   PetscInt       n;
664:   Vec*           vecs;
665:   PetscScalar*   alpha;                 /* for projections */
666:   PetscErrorCode (*remove)(MatNullSpace,Vec,void*);  /* for user provided removal function */
667:   void*          rmctx;                 /* context for remove() function */
668: };

670: /*
671:    Checking zero pivot for LU, ILU preconditioners.
672: */
673: typedef struct {
674:   PetscInt       nshift,nshift_max;
675:   PetscReal      shift_amount,shift_lo,shift_hi,shift_top,shift_fraction;
676:   PetscBool      newshift;
677:   PetscReal      rs;  /* active row sum of abs(offdiagonals) */
678:   PetscScalar    pv;  /* pivot of the active row */
679: } FactorShiftCtx;

681: /*
682:  Used by MatCreateSubMatrices_MPIXAIJ_Local()
683: */
684: #include <petscctable.h>
685: typedef struct { /* used by MatCreateSubMatrices_MPIAIJ_SingleIS_Local() and MatCreateSubMatrices_MPIAIJ_Local */
686:   PetscInt   id;   /* index of submats, only submats[0] is responsible for deleting some arrays below */
687:   PetscInt   nrqs,nrqr;
688:   PetscInt   **rbuf1,**rbuf2,**rbuf3,**sbuf1,**sbuf2;
689:   PetscInt   **ptr;
690:   PetscInt   *tmp;
691:   PetscInt   *ctr;
692:   PetscInt   *pa; /* proc array */
693:   PetscInt   *req_size,*req_source1,*req_source2;
694:   PetscBool  allcolumns,allrows;
695:   PetscBool  singleis;
696:   PetscInt   *row2proc; /* row to proc map */
697:   PetscInt   nstages;
698: #if defined(PETSC_USE_CTABLE)
699:   PetscTable cmap,rmap;
700:   PetscInt   *cmap_loc,*rmap_loc;
701: #else
702:   PetscInt   *cmap,*rmap;
703: #endif

705:   PetscErrorCode (*destroy)(Mat);
706: } Mat_SubSppt;

708: PETSC_EXTERN PetscErrorCode MatFactorDumpMatrix(Mat);
709: PETSC_INTERN PetscErrorCode MatShift_Basic(Mat,PetscScalar);
710: PETSC_INTERN PetscErrorCode MatSetBlockSizes_Default(Mat,PetscInt,PetscInt);

712: PETSC_STATIC_INLINE PetscErrorCode MatPivotCheck_nz(Mat mat,const MatFactorInfo *info,FactorShiftCtx *sctx,PetscInt row)
713: {
714:   PetscReal _rs   = sctx->rs;
715:   PetscReal _zero = info->zeropivot*_rs;

718:   if (PetscAbsScalar(sctx->pv) <= _zero && !PetscIsNanScalar(sctx->pv)){
719:     /* force |diag| > zeropivot*rs */
720:     if (!sctx->nshift) sctx->shift_amount = info->shiftamount;
721:     else sctx->shift_amount *= 2.0;
722:     sctx->newshift = PETSC_TRUE;
723:     (sctx->nshift)++;
724:   } else {
725:     sctx->newshift = PETSC_FALSE;
726:   }
727:   return(0);
728: }

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

736:   if (PetscRealPart(sctx->pv) <= _zero && !PetscIsNanScalar(sctx->pv)){
737:     /* force matfactor to be diagonally dominant */
738:     if (sctx->nshift == sctx->nshift_max) {
739:       sctx->shift_fraction = sctx->shift_hi;
740:     } else {
741:       sctx->shift_lo = sctx->shift_fraction;
742:       sctx->shift_fraction = (sctx->shift_hi+sctx->shift_lo)/2.;
743:     }
744:     sctx->shift_amount = sctx->shift_fraction * sctx->shift_top;
745:     sctx->nshift++;
746:     sctx->newshift = PETSC_TRUE;
747:   } else {
748:     sctx->newshift = PETSC_FALSE;
749:   }
750:   return(0);
751: }

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

758:   if (PetscAbsScalar(sctx->pv) <= _zero && !PetscIsNanScalar(sctx->pv)){
759:     sctx->pv          += info->shiftamount;
760:     sctx->shift_amount = 0.0;
761:     sctx->nshift++;
762:   }
763:   sctx->newshift = PETSC_FALSE;
764:   return(0);
765: }

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

773:   sctx->newshift = PETSC_FALSE;
774:   if (PetscAbsScalar(sctx->pv) <= _zero && !PetscIsNanScalar(sctx->pv)) {
775:     if (!mat->erroriffailure) {
776:       PetscInfo3(mat,"Detected zero pivot in factorization in row %D value %g tolerance %g\n",row,(double)PetscAbsScalar(sctx->pv),(double)_zero);
777:       fact->factorerrortype             = MAT_FACTOR_NUMERIC_ZEROPIVOT;
778:       fact->factorerror_zeropivot_value = PetscAbsScalar(sctx->pv);
779:       fact->factorerror_zeropivot_row   = row;
780:     } 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);
781:   }
782:   return(0);
783: }

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

790:   if (info->shifttype == (PetscReal) MAT_SHIFT_NONZERO){
791:     MatPivotCheck_nz(mat,info,sctx,row);
792:   } else if (info->shifttype == (PetscReal) MAT_SHIFT_POSITIVE_DEFINITE){
793:     MatPivotCheck_pd(mat,info,sctx,row);
794:   } else if (info->shifttype == (PetscReal) MAT_SHIFT_INBLOCKS){
795:     MatPivotCheck_inblocks(mat,info,sctx,row);
796:   } else {
797:     MatPivotCheck_none(fact,mat,info,sctx,row);
798:   }
799:   return(0);
800: }

802: /*
803:   Create and initialize a linked list
804:   Input Parameters:
805:     idx_start - starting index of the list
806:     lnk_max   - max value of lnk indicating the end of the list
807:     nlnk      - max length of the list
808:   Output Parameters:
809:     lnk       - list initialized
810:     bt        - PetscBT (bitarray) with all bits set to false
811:     lnk_empty - flg indicating the list is empty
812: */
813: #define PetscLLCreate(idx_start,lnk_max,nlnk,lnk,bt) \
814:   (PetscMalloc1(nlnk,&lnk) || PetscBTCreate(nlnk,&(bt)) || (lnk[idx_start] = lnk_max,0))

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

819: /*
820:   Add an index set into a sorted linked list
821:   Input Parameters:
822:     nidx      - number of input indices
823:     indices   - integer array
824:     idx_start - starting index of the list
825:     lnk       - linked list(an integer array) that is created
826:     bt        - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
827:   output Parameters:
828:     nlnk      - number of newly added indices
829:     lnk       - the sorted(increasing order) linked list containing new and non-redundate entries from indices
830:     bt        - updated PetscBT (bitarray)
831: */
832: #define PetscLLAdd(nidx,indices,idx_start,nlnk,lnk,bt) 0;\
833: {\
834:   PetscInt _k,_entry,_location,_lnkdata;\
835:   nlnk     = 0;\
836:   _lnkdata = idx_start;\
837:   for (_k=0; _k<nidx; _k++){\
838:     _entry = indices[_k];\
839:     if (!PetscBTLookupSet(bt,_entry)){  /* new entry */\
840:       /* search for insertion location */\
841:       /* start from the beginning if _entry < previous _entry */\
842:       if (_k && _entry < _lnkdata) _lnkdata  = idx_start;\
843:       do {\
844:         _location = _lnkdata;\
845:         _lnkdata  = lnk[_location];\
846:       } while (_entry > _lnkdata);\
847:       /* insertion location is found, add entry into lnk */\
848:       lnk[_location] = _entry;\
849:       lnk[_entry]    = _lnkdata;\
850:       nlnk++;\
851:       _lnkdata = _entry; /* next search starts from here if next_entry > _entry */\
852:     }\
853:   }\
854: }

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

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

929: #define PetscLLAddSorted_new(nidx,indices,idx_start,lnk_empty,nlnk,lnk,bt) 0; \
930: {\
931:   PetscInt _k,_entry,_location,_lnkdata;\
932:   if (lnk_empty){\
933:     _lnkdata  = idx_start;                      \
934:     for (_k=0; _k<nidx; _k++){                  \
935:       _entry = indices[_k];                             \
936:       PetscBTSet(bt,_entry);  /* mark the new entry */          \
937:           _location = _lnkdata;                                 \
938:           _lnkdata  = lnk[_location];                           \
939:         /* insertion location is found, add entry into lnk */   \
940:         lnk[_location] = _entry;                                \
941:         lnk[_entry]    = _lnkdata;                              \
942:         _lnkdata = _entry; /* next search starts from here */   \
943:     }                                                           \
944:     /*\
945:     lnk[indices[nidx-1]] = lnk[idx_start];\
946:     lnk[idx_start]       = indices[0];\
947:     PetscBTSet(bt,indices[0]);  \
948:     for (_k=1; _k<nidx; _k++){                  \
949:       PetscBTSet(bt,indices[_k]);                                          \
950:       lnk[indices[_k-1]] = indices[_k];                                  \
951:     }                                                           \
952:      */\
953:     nlnk      = nidx;\
954:     lnk_empty = PETSC_FALSE;\
955:   } else {\
956:     nlnk      = 0;                              \
957:     _lnkdata  = idx_start;                      \
958:     for (_k=0; _k<nidx; _k++){                  \
959:       _entry = indices[_k];                             \
960:       if (!PetscBTLookupSet(bt,_entry)){  /* new entry */       \
961:         /* search for insertion location */                     \
962:         do {                                                    \
963:           _location = _lnkdata;                                 \
964:           _lnkdata  = lnk[_location];                           \
965:         } while (_entry > _lnkdata);                            \
966:         /* insertion location is found, add entry into lnk */   \
967:         lnk[_location] = _entry;                                \
968:         lnk[_entry]    = _lnkdata;                              \
969:         nlnk++;                                                 \
970:         _lnkdata = _entry; /* next search starts from here */   \
971:       }                                                         \
972:     }                                                           \
973:   }                                                             \
974: }

976: /*
977:   Add a SORTED index set into a sorted linked list used for LUFactorSymbolic()
978:   Same as PetscLLAddSorted() with an additional operation:
979:        count the number of input indices that are no larger than 'diag'
980:   Input Parameters:
981:     indices   - sorted integer array
982:     idx_start - starting index of the list, index of pivot row
983:     lnk       - linked list(an integer array) that is created
984:     bt        - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
985:     diag      - index of the active row in LUFactorSymbolic
986:     nzbd      - number of input indices with indices <= idx_start
987:     im        - im[idx_start] is initialized as num of nonzero entries in row=idx_start
988:   output Parameters:
989:     nlnk      - number of newly added indices
990:     lnk       - the sorted(increasing order) linked list containing new and non-redundate entries from indices
991:     bt        - updated PetscBT (bitarray)
992:     im        - im[idx_start]: unchanged if diag is not an entry
993:                              : num of entries with indices <= diag if diag is an entry
994: */
995: #define PetscLLAddSortedLU(indices,idx_start,nlnk,lnk,bt,diag,nzbd,im) 0;\
996: {\
997:   PetscInt _k,_entry,_location,_lnkdata,_nidx;\
998:   nlnk     = 0;\
999:   _lnkdata = idx_start;\
1000:   _nidx = im[idx_start] - nzbd; /* num of entries with idx_start < index <= diag */\
1001:   for (_k=0; _k<_nidx; _k++){\
1002:     _entry = indices[_k];\
1003:     nzbd++;\
1004:     if (_entry== diag) im[idx_start] = nzbd;\
1005:     if (!PetscBTLookupSet(bt,_entry)){  /* new entry */\
1006:       /* search for insertion location */\
1007:       do {\
1008:         _location = _lnkdata;\
1009:         _lnkdata  = lnk[_location];\
1010:       } while (_entry > _lnkdata);\
1011:       /* insertion location is found, add entry into lnk */\
1012:       lnk[_location] = _entry;\
1013:       lnk[_entry]    = _lnkdata;\
1014:       nlnk++;\
1015:       _lnkdata = _entry; /* next search starts from here */\
1016:     }\
1017:   }\
1018: }

1020: /*
1021:   Copy data on the list into an array, then initialize the list
1022:   Input Parameters:
1023:     idx_start - starting index of the list
1024:     lnk_max   - max value of lnk indicating the end of the list
1025:     nlnk      - number of data on the list to be copied
1026:     lnk       - linked list
1027:     bt        - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
1028:   output Parameters:
1029:     indices   - array that contains the copied data
1030:     lnk       - linked list that is cleaned and initialize
1031:     bt        - PetscBT (bitarray) with all bits set to false
1032: */
1033: #define PetscLLClean(idx_start,lnk_max,nlnk,lnk,indices,bt) 0;\
1034: {\
1035:   PetscInt _j,_idx=idx_start;\
1036:   for (_j=0; _j<nlnk; _j++){\
1037:     _idx = lnk[_idx];\
1038:     indices[_j] = _idx;\
1039:     PetscBTClear(bt,_idx);\
1040:   }\
1041:   lnk[idx_start] = lnk_max;\
1042: }
1043: /*
1044:   Free memories used by the list
1045: */
1046: #define PetscLLDestroy(lnk,bt) (PetscFree(lnk) || PetscBTDestroy(&(bt)))

1048: /* Routines below are used for incomplete matrix factorization */
1049: /*
1050:   Create and initialize a linked list and its levels
1051:   Input Parameters:
1052:     idx_start - starting index of the list
1053:     lnk_max   - max value of lnk indicating the end of the list
1054:     nlnk      - max length of the list
1055:   Output Parameters:
1056:     lnk       - list initialized
1057:     lnk_lvl   - array of size nlnk for storing levels of lnk
1058:     bt        - PetscBT (bitarray) with all bits set to false
1059: */
1060: #define PetscIncompleteLLCreate(idx_start,lnk_max,nlnk,lnk,lnk_lvl,bt)\
1061:   (PetscIntMultError(2,nlnk,NULL) || PetscMalloc1(2*nlnk,&lnk) || PetscBTCreate(nlnk,&(bt)) || (lnk[idx_start] = lnk_max,lnk_lvl = lnk + nlnk,0))

1063: /*
1064:   Initialize a sorted linked list used for ILU and ICC
1065:   Input Parameters:
1066:     nidx      - number of input idx
1067:     idx       - integer array used for storing column indices
1068:     idx_start - starting index of the list
1069:     perm      - indices of an IS
1070:     lnk       - linked list(an integer array) that is created
1071:     lnklvl    - levels of lnk
1072:     bt        - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
1073:   output Parameters:
1074:     nlnk     - number of newly added idx
1075:     lnk      - the sorted(increasing order) linked list containing new and non-redundate entries from idx
1076:     lnklvl   - levels of lnk
1077:     bt       - updated PetscBT (bitarray)
1078: */
1079: #define PetscIncompleteLLInit(nidx,idx,idx_start,perm,nlnk,lnk,lnklvl,bt) 0;\
1080: {\
1081:   PetscInt _k,_entry,_location,_lnkdata;\
1082:   nlnk     = 0;\
1083:   _lnkdata = idx_start;\
1084:   for (_k=0; _k<nidx; _k++){\
1085:     _entry = perm[idx[_k]];\
1086:     if (!PetscBTLookupSet(bt,_entry)){  /* new entry */\
1087:       /* search for insertion location */\
1088:       if (_k && _entry < _lnkdata) _lnkdata  = idx_start;\
1089:       do {\
1090:         _location = _lnkdata;\
1091:         _lnkdata  = lnk[_location];\
1092:       } while (_entry > _lnkdata);\
1093:       /* insertion location is found, add entry into lnk */\
1094:       lnk[_location]  = _entry;\
1095:       lnk[_entry]     = _lnkdata;\
1096:       lnklvl[_entry] = 0;\
1097:       nlnk++;\
1098:       _lnkdata = _entry; /* next search starts from here if next_entry > _entry */\
1099:     }\
1100:   }\
1101: }

1103: /*
1104:   Add a SORTED index set into a sorted linked list for ILU
1105:   Input Parameters:
1106:     nidx      - number of input indices
1107:     idx       - sorted integer array used for storing column indices
1108:     level     - level of fill, e.g., ICC(level)
1109:     idxlvl    - level of idx
1110:     idx_start - starting index of the list
1111:     lnk       - linked list(an integer array) that is created
1112:     lnklvl    - levels of lnk
1113:     bt        - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
1114:     prow      - the row number of idx
1115:   output Parameters:
1116:     nlnk     - number of newly added idx
1117:     lnk      - the sorted(increasing order) linked list containing new and non-redundate entries from idx
1118:     lnklvl   - levels of lnk
1119:     bt       - updated PetscBT (bitarray)

1121:   Note: the level of factor(i,j) is set as lvl(i,j) = min{ lvl(i,j), lvl(i,prow)+lvl(prow,j)+1)
1122:         where idx = non-zero columns of U(prow,prow+1:n-1), prow<i
1123: */
1124: #define PetscILULLAddSorted(nidx,idx,level,idxlvl,idx_start,nlnk,lnk,lnklvl,bt,lnklvl_prow) 0;\
1125: {\
1126:   PetscInt _k,_entry,_location,_lnkdata,_incrlev,_lnklvl_prow=lnklvl[prow];\
1127:   nlnk     = 0;\
1128:   _lnkdata = idx_start;\
1129:   for (_k=0; _k<nidx; _k++){\
1130:     _incrlev = idxlvl[_k] + _lnklvl_prow + 1;\
1131:     if (_incrlev > level) continue;\
1132:     _entry = idx[_k];\
1133:     if (!PetscBTLookupSet(bt,_entry)){  /* new entry */\
1134:       /* search for insertion location */\
1135:       do {\
1136:         _location = _lnkdata;\
1137:         _lnkdata  = lnk[_location];\
1138:       } while (_entry > _lnkdata);\
1139:       /* insertion location is found, add entry into lnk */\
1140:       lnk[_location]  = _entry;\
1141:       lnk[_entry]     = _lnkdata;\
1142:       lnklvl[_entry] = _incrlev;\
1143:       nlnk++;\
1144:       _lnkdata = _entry; /* next search starts from here if next_entry > _entry */\
1145:     } else { /* existing entry: update lnklvl */\
1146:       if (lnklvl[_entry] > _incrlev) lnklvl[_entry] = _incrlev;\
1147:     }\
1148:   }\
1149: }

1151: /*
1152:   Add a index set into a sorted linked list
1153:   Input Parameters:
1154:     nidx      - number of input idx
1155:     idx   - integer array used for storing column indices
1156:     level     - level of fill, e.g., ICC(level)
1157:     idxlvl - level of idx
1158:     idx_start - starting index of the list
1159:     lnk       - linked list(an integer array) that is created
1160:     lnklvl   - levels of lnk
1161:     bt        - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
1162:   output Parameters:
1163:     nlnk      - number of newly added idx
1164:     lnk       - the sorted(increasing order) linked list containing new and non-redundate entries from idx
1165:     lnklvl   - levels of lnk
1166:     bt        - updated PetscBT (bitarray)
1167: */
1168: #define PetscIncompleteLLAdd(nidx,idx,level,idxlvl,idx_start,nlnk,lnk,lnklvl,bt) 0;\
1169: {\
1170:   PetscInt _k,_entry,_location,_lnkdata,_incrlev;\
1171:   nlnk     = 0;\
1172:   _lnkdata = idx_start;\
1173:   for (_k=0; _k<nidx; _k++){\
1174:     _incrlev = idxlvl[_k] + 1;\
1175:     if (_incrlev > level) continue;\
1176:     _entry = idx[_k];\
1177:     if (!PetscBTLookupSet(bt,_entry)){  /* new entry */\
1178:       /* search for insertion location */\
1179:       if (_k && _entry < _lnkdata) _lnkdata  = idx_start;\
1180:       do {\
1181:         _location = _lnkdata;\
1182:         _lnkdata  = lnk[_location];\
1183:       } while (_entry > _lnkdata);\
1184:       /* insertion location is found, add entry into lnk */\
1185:       lnk[_location]  = _entry;\
1186:       lnk[_entry]     = _lnkdata;\
1187:       lnklvl[_entry] = _incrlev;\
1188:       nlnk++;\
1189:       _lnkdata = _entry; /* next search starts from here if next_entry > _entry */\
1190:     } else { /* existing entry: update lnklvl */\
1191:       if (lnklvl[_entry] > _incrlev) lnklvl[_entry] = _incrlev;\
1192:     }\
1193:   }\
1194: }

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

1240: /*
1241:   Add a SORTED index set into a sorted linked list for ICC
1242:   Input Parameters:
1243:     nidx      - number of input indices
1244:     idx       - sorted integer array used for storing column indices
1245:     level     - level of fill, e.g., ICC(level)
1246:     idxlvl    - level of idx
1247:     idx_start - starting index of the list
1248:     lnk       - linked list(an integer array) that is created
1249:     lnklvl    - levels of lnk
1250:     bt        - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
1251:     idxlvl_prow - idxlvl[prow], where prow is the row number of the idx
1252:   output Parameters:
1253:     nlnk   - number of newly added indices
1254:     lnk    - the sorted(increasing order) linked list containing new and non-redundate entries from idx
1255:     lnklvl - levels of lnk
1256:     bt     - updated PetscBT (bitarray)
1257:   Note: the level of U(i,j) is set as lvl(i,j) = min{ lvl(i,j), lvl(prow,i)+lvl(prow,j)+1)
1258:         where idx = non-zero columns of U(prow,prow+1:n-1), prow<i
1259: */
1260: #define PetscICCLLAddSorted(nidx,idx,level,idxlvl,idx_start,nlnk,lnk,lnklvl,bt,idxlvl_prow) 0;\
1261: {\
1262:   PetscInt _k,_entry,_location,_lnkdata,_incrlev;\
1263:   nlnk = 0;\
1264:   _lnkdata = idx_start;\
1265:   for (_k=0; _k<nidx; _k++){\
1266:     _incrlev = idxlvl[_k] + idxlvl_prow + 1;\
1267:     if (_incrlev > level) continue;\
1268:     _entry = idx[_k];\
1269:     if (!PetscBTLookupSet(bt,_entry)){  /* new entry */\
1270:       /* search for insertion location */\
1271:       do {\
1272:         _location = _lnkdata;\
1273:         _lnkdata  = lnk[_location];\
1274:       } while (_entry > _lnkdata);\
1275:       /* insertion location is found, add entry into lnk */\
1276:       lnk[_location] = _entry;\
1277:       lnk[_entry]    = _lnkdata;\
1278:       lnklvl[_entry] = _incrlev;\
1279:       nlnk++;\
1280:       _lnkdata = _entry; /* next search starts from here */\
1281:     } else { /* existing entry: update lnklvl */\
1282:       if (lnklvl[_entry] > _incrlev) lnklvl[_entry] = _incrlev;\
1283:     }\
1284:   }\
1285: }

1287: /*
1288:   Copy data on the list into an array, then initialize the list
1289:   Input Parameters:
1290:     idx_start - starting index of the list
1291:     lnk_max   - max value of lnk indicating the end of the list
1292:     nlnk      - number of data on the list to be copied
1293:     lnk       - linked list
1294:     lnklvl    - level of lnk
1295:     bt        - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
1296:   output Parameters:
1297:     indices - array that contains the copied data
1298:     lnk     - linked list that is cleaned and initialize
1299:     lnklvl  - level of lnk that is reinitialized
1300:     bt      - PetscBT (bitarray) with all bits set to false
1301: */
1302: #define PetscIncompleteLLClean(idx_start,lnk_max,nlnk,lnk,lnklvl,indices,indiceslvl,bt) 0;\
1303: do {\
1304:   PetscInt _j,_idx=idx_start;\
1305:   for (_j=0; _j<nlnk; _j++){\
1306:     _idx = lnk[_idx];\
1307:     *(indices+_j) = _idx;\
1308:     *(indiceslvl+_j) = lnklvl[_idx];\
1309:     lnklvl[_idx] = -1;\
1310:     PetscBTClear(bt,_idx);\
1311:   }\
1312:   lnk[idx_start] = lnk_max;\
1313: } while (0)
1314: /*
1315:   Free memories used by the list
1316: */
1317: #define PetscIncompleteLLDestroy(lnk,bt) (PetscFree(lnk) || PetscBTDestroy(&(bt)))

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

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

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

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

1340:       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
1341:       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
1342:       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
1343:       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.
1344:       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
1345:       to each other so memory access is much better than using the big array.

1347:   Example:
1348:      nlnk_max=5, lnk_max=36:
1349:      Initial list: [0, 0 | 36, 2 | 0, 0 | 0, 0 | 0, 0 | 0, 0 | 0, 0]
1350:      here, head_node has index 2 with value lnk[2]=lnk_max=36,
1351:            0-th entry is used to store the number of entries in the list,
1352:      The initial lnk represents head -> tail(marked by 36) with number of entries = lnk[0]=0.

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

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

1362:   Input Parameters:
1363:     nlnk_max  - max length of the list
1364:     lnk_max   - max value of the entries
1365:   Output Parameters:
1366:     lnk       - list created and initialized
1367:     bt        - PetscBT (bitarray) with all bits set to false. Note: bt has size lnk_max, not nln_max!
1368: */
1369: PETSC_STATIC_INLINE PetscErrorCode PetscLLCondensedCreate(PetscInt nlnk_max,PetscInt lnk_max,PetscInt **lnk,PetscBT *bt)
1370: {
1372:   PetscInt       *llnk,lsize = 0;

1375:   PetscIntMultError(2,nlnk_max+2,&lsize);
1376:   PetscMalloc1(lsize,lnk);
1377:   PetscBTCreate(lnk_max,bt);
1378:   llnk = *lnk;
1379:   llnk[0] = 0;         /* number of entries on the list */
1380:   llnk[2] = lnk_max;   /* value in the head node */
1381:   llnk[3] = 2;         /* next for the head node */
1382:   return(0);
1383: }

1385: /*
1386:   Add a SORTED ascending index set into a sorted linked list. See PetscLLCondensedCreate() for detailed description.
1387:   Input Parameters:
1388:     nidx      - number of input indices
1389:     indices   - sorted integer array
1390:     lnk       - condensed linked list(an integer array) that is created
1391:     bt        - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
1392:   output Parameters:
1393:     lnk       - the sorted(increasing order) linked list containing previous and newly added non-redundate indices
1394:     bt        - updated PetscBT (bitarray)
1395: */
1396: PETSC_STATIC_INLINE PetscErrorCode PetscLLCondensedAddSorted(PetscInt nidx,const PetscInt indices[],PetscInt lnk[],PetscBT bt)
1397: {
1398:   PetscInt _k,_entry,_location,_next,_lnkdata,_nlnk,_newnode;

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

1425: PETSC_STATIC_INLINE PetscErrorCode PetscLLCondensedClean(PetscInt lnk_max,PetscInt nidx,PetscInt *indices,PetscInt lnk[],PetscBT bt)
1426: {
1428:   PetscInt       _k,_next,_nlnk;

1431:   _next = lnk[3];       /* head node */
1432:   _nlnk = lnk[0];       /* num of entries on the list */
1433:   for (_k=0; _k<_nlnk; _k++){
1434:     indices[_k] = lnk[_next];
1435:     _next       = lnk[_next + 1];
1436:     PetscBTClear(bt,indices[_k]);
1437:   }
1438:   lnk[0] = 0;          /* num of entries on the list */
1439:   lnk[2] = lnk_max;    /* initialize head node */
1440:   lnk[3] = 2;          /* head node */
1441:   return(0);
1442: }

1444: PETSC_STATIC_INLINE PetscErrorCode PetscLLCondensedView(PetscInt *lnk)
1445: {
1447:   PetscInt       k;

1450:   PetscPrintf(PETSC_COMM_SELF,"LLCondensed of size %D, (val,  next)\n",lnk[0]);
1451:   for (k=2; k< lnk[0]+2; k++){
1452:     PetscPrintf(PETSC_COMM_SELF," %D: (%D, %D)\n",2*k,lnk[2*k],lnk[2*k+1]);
1453:   }
1454:   return(0);
1455: }

1457: /*
1458:   Free memories used by the list
1459: */
1460: PETSC_STATIC_INLINE PetscErrorCode PetscLLCondensedDestroy(PetscInt *lnk,PetscBT bt)
1461: {

1465:   PetscFree(lnk);
1466:   PetscBTDestroy(&bt);
1467:   return(0);
1468: }

1470: /* -------------------------------------------------------------------------------------------------------*/
1471: /*
1472:  Same as PetscLLCondensedCreate(), but does not use non-scalable O(lnk_max) bitarray
1473:   Input Parameters:
1474:     nlnk_max  - max length of the list
1475:   Output Parameters:
1476:     lnk       - list created and initialized
1477: */
1478: PETSC_STATIC_INLINE PetscErrorCode PetscLLCondensedCreate_Scalable(PetscInt nlnk_max,PetscInt **lnk)
1479: {
1481:   PetscInt       *llnk,lsize = 0;

1484:   PetscIntMultError(2,nlnk_max+2,&lsize);
1485:   PetscMalloc1(lsize,lnk);
1486:   llnk = *lnk;
1487:   llnk[0] = 0;               /* number of entries on the list */
1488:   llnk[2] = PETSC_MAX_INT;   /* value in the head node */
1489:   llnk[3] = 2;               /* next for the head node */
1490:   return(0);
1491: }

1493: PETSC_STATIC_INLINE PetscErrorCode PetscLLCondensedExpand_Scalable(PetscInt nlnk_max,PetscInt **lnk)
1494: {
1496:   PetscInt       lsize = 0;

1499:   PetscIntMultError(2,nlnk_max+2,&lsize);
1500:   PetscRealloc(lsize*sizeof(PetscInt),lnk);
1501:   return(0);
1502: }

1504: PETSC_STATIC_INLINE PetscErrorCode PetscLLCondensedAddSorted_Scalable(PetscInt nidx,const PetscInt indices[],PetscInt lnk[])
1505: {
1506:   PetscInt _k,_entry,_location,_next,_lnkdata,_nlnk,_newnode;
1507:   _nlnk     = lnk[0]; /* num of entries on the input lnk */
1508:   _location = 2; /* head */ \
1509:     for (_k=0; _k<nidx; _k++){
1510:       _entry = indices[_k];
1511:       /* search for insertion location */
1512:       do {
1513:         _next     = _location + 1; /* link from previous node to next node */
1514:         _location = lnk[_next];    /* idx of next node */
1515:         _lnkdata  = lnk[_location];/* value of next node */
1516:       } while (_entry > _lnkdata);
1517:       if (_entry < _lnkdata) {
1518:         /* insertion location is found, add entry into lnk */
1519:         _newnode        = 2*(_nlnk+2);   /* index for this new node */
1520:         lnk[_next]      = _newnode;      /* connect previous node to the new node */
1521:         lnk[_newnode]   = _entry;        /* set value of the new node */
1522:         lnk[_newnode+1] = _location;     /* connect new node to next node */
1523:         _location       = _newnode;      /* next search starts from the new node */
1524:         _nlnk++;
1525:       }
1526:     }
1527:   lnk[0]   = _nlnk;   /* number of entries in the list */
1528:   return 0;
1529: }

1531: PETSC_STATIC_INLINE PetscErrorCode PetscLLCondensedClean_Scalable(PetscInt nidx,PetscInt *indices,PetscInt *lnk)
1532: {
1533:   PetscInt _k,_next,_nlnk;
1534:   _next = lnk[3];       /* head node */
1535:   _nlnk = lnk[0];
1536:   for (_k=0; _k<_nlnk; _k++){
1537:     indices[_k] = lnk[_next];
1538:     _next       = lnk[_next + 1];
1539:   }
1540:   lnk[0] = 0;          /* num of entries on the list */
1541:   lnk[3] = 2;          /* head node */
1542:   return 0;
1543: }

1545: PETSC_STATIC_INLINE PetscErrorCode PetscLLCondensedDestroy_Scalable(PetscInt *lnk)
1546: {
1547:   return PetscFree(lnk);
1548: }

1550: /* -------------------------------------------------------------------------------------------------------*/
1551: /*
1552:       lnk[0]   number of links
1553:       lnk[1]   number of entries
1554:       lnk[3n]  value
1555:       lnk[3n+1] len
1556:       lnk[3n+2] link to next value

1558:       The next three are always the first link

1560:       lnk[3]    PETSC_MIN_INT+1
1561:       lnk[4]    1
1562:       lnk[5]    link to first real entry

1564:       The next three are always the last link

1566:       lnk[6]    PETSC_MAX_INT - 1
1567:       lnk[7]    1
1568:       lnk[8]    next valid link (this is the same as lnk[0] but without the decreases)
1569: */

1571: PETSC_STATIC_INLINE PetscErrorCode PetscLLCondensedCreate_fast(PetscInt nlnk_max,PetscInt **lnk)
1572: {
1574:   PetscInt       *llnk,lsize = 0;

1577:   PetscIntMultError(3,nlnk_max+3,&lsize);
1578:   PetscMalloc1(lsize,lnk);
1579:   llnk = *lnk;
1580:   llnk[0] = 0;   /* nlnk: number of entries on the list */
1581:   llnk[1] = 0;          /* number of integer entries represented in list */
1582:   llnk[3] = PETSC_MIN_INT+1;   /* value in the first node */
1583:   llnk[4] = 1;           /* count for the first node */
1584:   llnk[5] = 6;         /* next for the first node */
1585:   llnk[6] = PETSC_MAX_INT-1;   /* value in the last node */
1586:   llnk[7] = 1;           /* count for the last node */
1587:   llnk[8] = 0;         /* next valid node to be used */
1588:   return(0);
1589: }

1591: PETSC_STATIC_INLINE PetscErrorCode PetscLLCondensedAddSorted_fast(PetscInt nidx,const PetscInt indices[],PetscInt lnk[])
1592: {
1593:   PetscInt k,entry,prev,next;
1594:   prev      = 3;      /* first value */
1595:   next      = lnk[prev+2];
1596:   for (k=0; k<nidx; k++){
1597:     entry = indices[k];
1598:     /* search for insertion location */
1599:     while (entry >= lnk[next]) {
1600:       prev = next;
1601:       next = lnk[next+2];
1602:     }
1603:     /* entry is in range of previous list */
1604:     if (entry < lnk[prev]+lnk[prev+1]) continue;
1605:     lnk[1]++;
1606:     /* entry is right after previous list */
1607:     if (entry == lnk[prev]+lnk[prev+1]) {
1608:       lnk[prev+1]++;
1609:       if (lnk[next] == entry+1) { /* combine two contiguous strings */
1610:         lnk[prev+1] += lnk[next+1];
1611:         lnk[prev+2]  = lnk[next+2];
1612:         next         = lnk[next+2];
1613:         lnk[0]--;
1614:       }
1615:       continue;
1616:     }
1617:     /* entry is right before next list */
1618:     if (entry == lnk[next]-1) {
1619:       lnk[next]--;
1620:       lnk[next+1]++;
1621:       prev = next;
1622:       next = lnk[prev+2];
1623:       continue;
1624:     }
1625:     /*  add entry into lnk */
1626:     lnk[prev+2]    = 3*((lnk[8]++)+3);      /* connect previous node to the new node */
1627:     prev           = lnk[prev+2];
1628:     lnk[prev]      = entry;        /* set value of the new node */
1629:     lnk[prev+1]    = 1;             /* number of values in contiguous string is one to start */
1630:     lnk[prev+2]    = next;          /* connect new node to next node */
1631:     lnk[0]++;
1632:   }
1633:   return 0;
1634: }

1636: PETSC_STATIC_INLINE PetscErrorCode PetscLLCondensedClean_fast(PetscInt nidx,PetscInt *indices,PetscInt *lnk)
1637: {
1638:   PetscInt _k,_next,_nlnk,cnt,j;
1639:   _next = lnk[5];       /* first node */
1640:   _nlnk = lnk[0];
1641:   cnt   = 0;
1642:   for (_k=0; _k<_nlnk; _k++){
1643:     for (j=0; j<lnk[_next+1]; j++) {
1644:       indices[cnt++] = lnk[_next] + j;
1645:     }
1646:     _next       = lnk[_next + 2];
1647:   }
1648:   lnk[0] = 0;   /* nlnk: number of links */
1649:   lnk[1] = 0;          /* number of integer entries represented in list */
1650:   lnk[3] = PETSC_MIN_INT+1;   /* value in the first node */
1651:   lnk[4] = 1;           /* count for the first node */
1652:   lnk[5] = 6;         /* next for the first node */
1653:   lnk[6] = PETSC_MAX_INT-1;   /* value in the last node */
1654:   lnk[7] = 1;           /* count for the last node */
1655:   lnk[8] = 0;         /* next valid location to make link */
1656:   return 0;
1657: }

1659: PETSC_STATIC_INLINE PetscErrorCode PetscLLCondensedView_fast(PetscInt *lnk)
1660: {
1661:   PetscInt k,next,nlnk;
1662:   next = lnk[5];       /* first node */
1663:   nlnk = lnk[0];
1664:   for (k=0; k<nlnk; k++){
1665: #if 0                           /* Debugging code */
1666:     printf("%d value %d len %d next %d\n",next,lnk[next],lnk[next+1],lnk[next+2]);
1667: #endif
1668:     next = lnk[next + 2];
1669:   }
1670:   return 0;
1671: }

1673: PETSC_STATIC_INLINE PetscErrorCode PetscLLCondensedDestroy_fast(PetscInt *lnk)
1674: {
1675:   return PetscFree(lnk);
1676: }

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

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

1781: #endif