Actual source code: matimpl.h

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

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

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

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

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


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

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

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

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

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

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

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

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

556:     ncolors = 4;

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

1443: PETSC_STATIC_INLINE PetscErrorCode PetscLLCondensedView(PetscInt *lnk)
1444: {
1446:   PetscInt       k;

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

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

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

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

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

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

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

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

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

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

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

1557:       The next three are always the first link

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

1563:       The next three are always the last link

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

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

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

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

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

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

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

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

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

1780: #endif