Actual source code: matimpl.h

petsc-master 2019-06-26
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 include/petsc/finclude/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 (*freeintermediatedatastructures)(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 (*matmatmult)(Mat,Mat,Mat,MatReuse,PetscReal,Mat*);
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)(Mat,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 (*matmult)(Mat,Mat,MatReuse,PetscReal,Mat*);
142:   PetscErrorCode (*matmultsymbolic)(Mat,Mat,PetscReal,Mat*);
143:   PetscErrorCode (*matmultnumeric)(Mat,Mat,Mat);
144:   PetscErrorCode (*ptap)(Mat,Mat,MatReuse,PetscReal,Mat*);
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 (*mattransposemult)(Mat,Mat,MatReuse,PetscReal,Mat*);
149:   PetscErrorCode (*mattransposemultsymbolic)(Mat,Mat,PetscReal,Mat*);
150:   PetscErrorCode (*mattransposemultnumeric)(Mat,Mat,Mat);
151:   PetscErrorCode (*pintocpu)(Mat,PetscBool);
152:   /*99*/
153:   PetscErrorCode (*placeholder_99)(Mat);
154:   PetscErrorCode (*placeholder_100)(Mat);
155:   PetscErrorCode (*placeholder_101)(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 (*transposematmult)(Mat,Mat,MatReuse,PetscReal,Mat*);
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 (*rart)(Mat,Mat,MatReuse,PetscReal,Mat*);
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:   /*144*/
207:   PetscErrorCode (*creatempimatconcatenateseqmat)(MPI_Comm,Mat,PetscInt,MatReuse,Mat*);
208:   PetscErrorCode (*destroysubmatrices)(PetscInt,Mat*[]);
209:   PetscErrorCode (*mattransposesolve)(Mat,Mat,Mat);
210: };
211: /*
212:     If you add MatOps entries above also add them to the MATOP enum
213:     in include/petscmat.h and include/petsc/finclude/petscmat.h
214: */

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

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

226: PETSC_EXTERN MatRootName MatRootNameList;

228: /*
229:    Utility private matrix routines
230: */
231: PETSC_INTERN PetscErrorCode MatFindNonzeroRowsOrCols_Basic(Mat,PetscBool,PetscReal,IS*);
232: PETSC_INTERN PetscErrorCode MatConvert_Basic(Mat,MatType,MatReuse,Mat*);
233: PETSC_INTERN PetscErrorCode MatConvert_Shell(Mat,MatType,MatReuse,Mat*);
234: PETSC_INTERN PetscErrorCode MatConvertFrom_Shell(Mat,MatType,MatReuse,Mat*);
235: PETSC_INTERN PetscErrorCode MatCopy_Basic(Mat,Mat,MatStructure);
236: PETSC_INTERN PetscErrorCode MatDiagonalSet_Default(Mat,Vec,InsertMode);

238: #if defined(PETSC_USE_DEBUG)
239: #  define MatCheckPreallocated(A,arg) do {                              \
240:     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); \
241:   } while (0)
242: #else
243: #  define MatCheckPreallocated(A,arg) do {} while (0)
244: #endif

246: /*
247:   The stash is used to temporarily store inserted matrix values that
248:   belong to another processor. During the assembly phase the stashed
249:   values are moved to the correct processor and
250: */

252: typedef struct _MatStashSpace *PetscMatStashSpace;

254: struct _MatStashSpace {
255:   PetscMatStashSpace next;
256:   PetscScalar        *space_head,*val;
257:   PetscInt           *idx,*idy;
258:   PetscInt           total_space_size;
259:   PetscInt           local_used;
260:   PetscInt           local_remaining;
261: };

263: PETSC_EXTERN PetscErrorCode PetscMatStashSpaceGet(PetscInt,PetscInt,PetscMatStashSpace *);
264: PETSC_EXTERN PetscErrorCode PetscMatStashSpaceContiguous(PetscInt,PetscMatStashSpace *,PetscScalar *,PetscInt *,PetscInt *);
265: PETSC_EXTERN PetscErrorCode PetscMatStashSpaceDestroy(PetscMatStashSpace*);

267: typedef struct {
268:   PetscInt    count;
269: } MatStashHeader;

271: typedef struct {
272:   void        *buffer;          /* Of type blocktype, dynamically constructed  */
273:   PetscInt    count;
274:   char        pending;
275: } MatStashFrame;

277: typedef struct _MatStash MatStash;
278: struct _MatStash {
279:   PetscInt      nmax;                   /* maximum stash size */
280:   PetscInt      umax;                   /* user specified max-size */
281:   PetscInt      oldnmax;                /* the nmax value used previously */
282:   PetscInt      n;                      /* stash size */
283:   PetscInt      bs;                     /* block size of the stash */
284:   PetscInt      reallocs;               /* preserve the no of mallocs invoked */
285:   PetscMatStashSpace space_head,space;  /* linked list to hold stashed global row/column numbers and matrix values */

287:   PetscErrorCode (*ScatterBegin)(Mat,MatStash*,PetscInt*);
288:   PetscErrorCode (*ScatterGetMesg)(MatStash*,PetscMPIInt*,PetscInt**,PetscInt**,PetscScalar**,PetscInt*);
289:   PetscErrorCode (*ScatterEnd)(MatStash*);
290:   PetscErrorCode (*ScatterDestroy)(MatStash*);

292:   /* The following variables are used for communication */
293:   MPI_Comm      comm;
294:   PetscMPIInt   size,rank;
295:   PetscMPIInt   tag1,tag2;
296:   MPI_Request   *send_waits;            /* array of send requests */
297:   MPI_Request   *recv_waits;            /* array of receive requests */
298:   MPI_Status    *send_status;           /* array of send status */
299:   PetscInt      nsends,nrecvs;          /* numbers of sends and receives */
300:   PetscScalar   *svalues;               /* sending data */
301:   PetscInt      *sindices;
302:   PetscScalar   **rvalues;              /* receiving data (values) */
303:   PetscInt      **rindices;             /* receiving data (indices) */
304:   PetscInt      nprocessed;             /* number of messages already processed */
305:   PetscMPIInt   *flg_v;                 /* indicates what messages have arrived so far and from whom */
306:   PetscBool     reproduce;
307:   PetscInt      reproduce_count;

309:   /* The following variables are used for BTS communication */
310:   PetscBool      first_assembly_done;   /* Is the first time matrix assembly done? */
311:   PetscBool      use_status;            /* Use MPI_Status to determine number of items in each message */
312:   PetscMPIInt    nsendranks;
313:   PetscMPIInt    nrecvranks;
314:   PetscMPIInt    *sendranks;
315:   PetscMPIInt    *recvranks;
316:   MatStashHeader *sendhdr,*recvhdr;
317:   MatStashFrame  *sendframes;   /* pointers to the main messages */
318:   MatStashFrame  *recvframes;
319:   MatStashFrame  *recvframe_active;
320:   PetscInt       recvframe_i;     /* index of block within active frame */
321:   PetscMPIInt    recvframe_count; /* Count actually sent for current frame */
322:   PetscInt       recvcount;       /* Number of receives processed so far */
323:   PetscMPIInt    *some_indices;   /* From last call to MPI_Waitsome */
324:   MPI_Status     *some_statuses;  /* Statuses from last call to MPI_Waitsome */
325:   PetscMPIInt    some_count;      /* Number of requests completed in last call to MPI_Waitsome */
326:   PetscMPIInt    some_i;          /* Index of request currently being processed */
327:   MPI_Request    *sendreqs;
328:   MPI_Request    *recvreqs;
329:   PetscSegBuffer segsendblocks;
330:   PetscSegBuffer segrecvframe;
331:   PetscSegBuffer segrecvblocks;
332:   MPI_Datatype   blocktype;
333:   size_t         blocktype_size;
334:   InsertMode     *insertmode;   /* Pointer to check mat->insertmode and set upon message arrival in case no local values have been set. */
335: };

337: #if !defined(PETSC_HAVE_MPIUNI)
338: PETSC_INTERN PetscErrorCode MatStashScatterDestroy_BTS(MatStash*);
339: #endif
340: PETSC_INTERN PetscErrorCode MatStashCreate_Private(MPI_Comm,PetscInt,MatStash*);
341: PETSC_INTERN PetscErrorCode MatStashDestroy_Private(MatStash*);
342: PETSC_INTERN PetscErrorCode MatStashScatterEnd_Private(MatStash*);
343: PETSC_INTERN PetscErrorCode MatStashSetInitialSize_Private(MatStash*,PetscInt);
344: PETSC_INTERN PetscErrorCode MatStashGetInfo_Private(MatStash*,PetscInt*,PetscInt*);
345: PETSC_INTERN PetscErrorCode MatStashValuesRow_Private(MatStash*,PetscInt,PetscInt,const PetscInt[],const PetscScalar[],PetscBool );
346: PETSC_INTERN PetscErrorCode MatStashValuesCol_Private(MatStash*,PetscInt,PetscInt,const PetscInt[],const PetscScalar[],PetscInt,PetscBool );
347: PETSC_INTERN PetscErrorCode MatStashValuesRowBlocked_Private(MatStash*,PetscInt,PetscInt,const PetscInt[],const PetscScalar[],PetscInt,PetscInt,PetscInt);
348: PETSC_INTERN PetscErrorCode MatStashValuesColBlocked_Private(MatStash*,PetscInt,PetscInt,const PetscInt[],const PetscScalar[],PetscInt,PetscInt,PetscInt);
349: PETSC_INTERN PetscErrorCode MatStashScatterBegin_Private(Mat,MatStash*,PetscInt*);
350: PETSC_INTERN PetscErrorCode MatStashScatterGetMesg_Private(MatStash*,PetscMPIInt*,PetscInt**,PetscInt**,PetscScalar**,PetscInt*);
351: PETSC_INTERN PetscErrorCode MatGetInfo_External(Mat,MatInfoType,MatInfo*);

353: typedef struct {
354:   PetscInt   dim;
355:   PetscInt   dims[4];
356:   PetscInt   starts[4];
357:   PetscBool  noc;        /* this is a single component problem, hence user will not set MatStencil.c */
358: } MatStencilInfo;

360: /* Info about using compressed row format */
361: typedef struct {
362:   PetscBool  use;                           /* indicates compressed rows have been checked and will be used */
363:   PetscInt   nrows;                         /* number of non-zero rows */
364:   PetscInt   *i;                            /* compressed row pointer  */
365:   PetscInt   *rindex;                       /* compressed row index               */
366: } Mat_CompressedRow;
367: PETSC_EXTERN PetscErrorCode MatCheckCompressedRow(Mat,PetscInt,Mat_CompressedRow*,PetscInt*,PetscInt,PetscReal);

369: typedef struct { /* used by MatCreateRedundantMatrix() for reusing matredundant */
370:   PetscInt     nzlocal,nsends,nrecvs;
371:   PetscMPIInt  *send_rank,*recv_rank;
372:   PetscInt     *sbuf_nz,*rbuf_nz,*sbuf_j,**rbuf_j;
373:   PetscScalar  *sbuf_a,**rbuf_a;
374:   MPI_Comm     subcomm;   /* when user does not provide a subcomm */
375:   IS           isrow,iscol;
376:   Mat          *matseq;
377: } Mat_Redundant;

379: struct _p_Mat {
380:   PETSCHEADER(struct _MatOps);
381:   PetscLayout            rmap,cmap;
382:   void                   *data;            /* implementation-specific data */
383:   MatFactorType          factortype;       /* MAT_FACTOR_LU, ILU, CHOLESKY or ICC */
384:   PetscBool              assembled;        /* is the matrix assembled? */
385:   PetscBool              was_assembled;    /* new values inserted into assembled mat */
386:   PetscInt               num_ass;          /* number of times matrix has been assembled */
387:   PetscObjectState       nonzerostate;     /* each time new nonzeros locations are introduced into the matrix this is updated */
388:   PetscObjectState       ass_nonzerostate; /* nonzero state at last assembly */
389:   MatInfo                info;             /* matrix information */
390:   InsertMode             insertmode;       /* have values been inserted in matrix or added? */
391:   MatStash               stash,bstash;     /* used for assembling off-proc mat emements */
392:   MatNullSpace           nullsp;           /* null space (operator is singular) */
393:   MatNullSpace           transnullsp;      /* null space of transpose of operator */
394:   MatNullSpace           nearnullsp;       /* near null space to be used by multigrid methods */
395:   PetscInt               congruentlayouts; /* are the rows and columns layouts congruent? */
396:   PetscBool              preallocated;
397:   MatStencilInfo         stencil;          /* information for structured grid */
398:   PetscBool              symmetric,hermitian,structurally_symmetric,spd;
399:   PetscBool              symmetric_set,hermitian_set,structurally_symmetric_set,spd_set; /* if true, then corresponding flag is correct*/
400:   PetscBool              symmetric_eternal;
401:   PetscBool              nooffprocentries,nooffproczerorows;
402:   PetscBool              assembly_subset;  /* set by MAT_SUBSET_OFF_PROC_ENTRIES */
403:   PetscBool              submat_singleis;  /* for efficient PCSetUP_ASM() */
404:   PetscBool              structure_only;
405:   PetscBool              sortedfull;       /* full, sorted rows are inserted */
406: #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_CUDA)
407:   PetscOffloadFlag       valid_GPU_matrix; /* flag pointing to the matrix on the gpu*/
408:   PetscBool              pinnedtocpu;
409: #endif
410:   void                   *spptr;          /* pointer for special library like SuperLU */
411:   char                   *solvertype;
412:   PetscBool              checksymmetryonassembly,checknullspaceonassembly;
413:   PetscReal              checksymmetrytol;
414:   Mat                    schur;             /* Schur complement matrix */
415:   MatFactorSchurStatus   schur_status;      /* status of the Schur complement matrix */
416:   Mat_Redundant          *redundant;        /* used by MatCreateRedundantMatrix() */
417:   PetscBool              erroriffailure;    /* Generate an error if detected (for example a zero pivot) instead of returning */
418:   MatFactorError         factorerrortype;               /* type of error in factorization */
419:   PetscReal              factorerror_zeropivot_value;   /* If numerical zero pivot was detected this is the computed value */
420:   PetscInt               factorerror_zeropivot_row;     /* Row where zero pivot was detected */
421:   PetscInt               nblocks,*bsizes;   /* support for MatSetVariableBlockSizes() */
422:   char                   *defaultvectype;
423: };

425: PETSC_INTERN PetscErrorCode MatAXPY_Basic(Mat,PetscScalar,Mat,MatStructure);
426: PETSC_INTERN PetscErrorCode MatAXPY_BasicWithPreallocation(Mat,Mat,PetscScalar,Mat,MatStructure);
427: PETSC_INTERN PetscErrorCode MatAXPY_Basic_Preallocate(Mat,Mat,Mat*);

429: /*
430:     Utility for MatFactor (Schur complement)
431: */
432: PETSC_INTERN PetscErrorCode MatFactorFactorizeSchurComplement_Private(Mat);
433: PETSC_INTERN PetscErrorCode MatFactorInvertSchurComplement_Private(Mat);
434: PETSC_INTERN PetscErrorCode MatFactorUpdateSchurStatus_Private(Mat);
435: PETSC_INTERN PetscErrorCode MatFactorSetUpInPlaceSchur_Private(Mat);

437: /*
438:     Utility for MatZeroRows
439: */
440: PETSC_INTERN PetscErrorCode MatZeroRowsMapLocal_Private(Mat,PetscInt,const PetscInt*,PetscInt*,PetscInt**);

442: /*
443:     Object for partitioning graphs
444: */

446: typedef struct _MatPartitioningOps *MatPartitioningOps;
447: struct _MatPartitioningOps {
448:   PetscErrorCode (*apply)(MatPartitioning,IS*);
449:   PetscErrorCode (*applynd)(MatPartitioning,IS*);
450:   PetscErrorCode (*setfromoptions)(PetscOptionItems*,MatPartitioning);
451:   PetscErrorCode (*destroy)(MatPartitioning);
452:   PetscErrorCode (*view)(MatPartitioning,PetscViewer);
453:   PetscErrorCode (*improve)(MatPartitioning,IS*);
454: };

456: struct _p_MatPartitioning {
457:   PETSCHEADER(struct _MatPartitioningOps);
458:   Mat         adj;
459:   PetscInt    *vertex_weights;
460:   PetscReal   *part_weights;
461:   PetscInt    n;                                 /* number of partitions */
462:   void        *data;
463:   PetscInt    setupcalled;
464: };

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

469: /*
470:     Object for coarsen graphs
471: */
472: typedef struct _MatCoarsenOps *MatCoarsenOps;
473: struct _MatCoarsenOps {
474:   PetscErrorCode (*apply)(MatCoarsen);
475:   PetscErrorCode (*setfromoptions)(PetscOptionItems*,MatCoarsen);
476:   PetscErrorCode (*destroy)(MatCoarsen);
477:   PetscErrorCode (*view)(MatCoarsen,PetscViewer);
478: };

480: struct _p_MatCoarsen {
481:   PETSCHEADER(struct _MatCoarsenOps);
482:   Mat              graph;
483:   PetscInt         setupcalled;
484:   void             *subctx;
485:   /* */
486:   PetscBool        strict_aggs;
487:   IS               perm;
488:   PetscCoarsenData *agg_lists;
489: };

491: /*
492:     MatFDColoring is used to compute Jacobian matrices efficiently
493:   via coloring. The data structure is explained below in an example.

495:    Color =   0    1     0    2   |   2      3       0
496:    ---------------------------------------------------
497:             00   01              |          05
498:             10   11              |   14     15               Processor  0
499:                        22    23  |          25
500:                        32    33  |
501:    ===================================================
502:                                  |   44     45     46
503:             50                   |          55               Processor 1
504:                                  |   64            66
505:    ---------------------------------------------------

507:     ncolors = 4;

509:     ncolumns      = {2,1,1,0}
510:     columns       = {{0,2},{1},{3},{}}
511:     nrows         = {4,2,3,3}
512:     rows          = {{0,1,2,3},{0,1},{1,2,3},{0,1,2}}
513:     vwscale       = {dx(0),dx(1),dx(2),dx(3)}               MPI Vec
514:     vscale        = {dx(0),dx(1),dx(2),dx(3),dx(4),dx(5)}   Seq Vec

516:     ncolumns      = {1,0,1,1}
517:     columns       = {{6},{},{4},{5}}
518:     nrows         = {3,0,2,2}
519:     rows          = {{0,1,2},{},{1,2},{1,2}}
520:     vwscale       = {dx(4),dx(5),dx(6)}              MPI Vec
521:     vscale        = {dx(0),dx(4),dx(5),dx(6)}        Seq Vec

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

526: */
527: typedef struct {
528:   PetscInt     row;
529:   PetscInt     col;
530:   PetscScalar  *valaddr;   /* address of value */
531: } MatEntry;

533: typedef struct {
534:   PetscInt     row;
535:   PetscScalar  *valaddr;   /* address of value */
536: } MatEntry2;

538: struct  _p_MatFDColoring{
539:   PETSCHEADER(int);
540:   PetscInt       M,N,m;            /* total rows, columns; local rows */
541:   PetscInt       rstart;           /* first row owned by local processor */
542:   PetscInt       ncolors;          /* number of colors */
543:   PetscInt       *ncolumns;        /* number of local columns for a color */
544:   PetscInt       **columns;        /* lists the local columns of each color (using global column numbering) */
545:   IS             *isa;             /* these are the IS that contain the column values given in columns */
546:   PetscInt       *nrows;           /* number of local rows for each color */
547:   MatEntry       *matentry;        /* holds (row, column, address of value) for Jacobian matrix entry */
548:   MatEntry2      *matentry2;       /* holds (row, address of value) for Jacobian matrix entry */
549:   PetscScalar    *dy;              /* store a block of F(x+dx)-F(x) when J is in BAIJ format */
550:   PetscReal      error_rel;        /* square root of relative error in computing function */
551:   PetscReal      umin;             /* minimum allowable u'dx value */
552:   Vec            w1,w2,w3;         /* work vectors used in computing Jacobian */
553:   PetscBool      fset;             /* indicates that the initial function value F(X) is set */
554:   PetscErrorCode (*f)(void);       /* function that defines Jacobian */
555:   void           *fctx;            /* optional user-defined context for use by the function f */
556:   Vec            vscale;           /* holds FD scaling, i.e. 1/dx for each perturbed column */
557:   PetscInt       currentcolor;     /* color for which function evaluation is being done now */
558:   const char     *htype;           /* "wp" or "ds" */
559:   ISColoringType ctype;            /* IS_COLORING_GLOBAL or IS_COLORING_LOCAL */
560:   PetscInt       brows,bcols;      /* number of block rows or columns for speedup inserting the dense matrix into sparse Jacobian */
561:   PetscBool      setupcalled;      /* true if setup has been called */
562:   PetscBool      viewed;           /* true if the -mat_fd_coloring_view has been triggered already */
563:   void           (*ftn_func_pointer)(void),*ftn_func_cntx; /* serve the same purpose as *fortran_func_pointers in PETSc objects */
564: };

566: typedef struct _MatColoringOps *MatColoringOps;
567: struct _MatColoringOps {
568:   PetscErrorCode (*destroy)(MatColoring);
569:   PetscErrorCode (*setfromoptions)(PetscOptionItems*,MatColoring);
570:   PetscErrorCode (*view)(MatColoring,PetscViewer);
571:   PetscErrorCode (*apply)(MatColoring,ISColoring*);
572:   PetscErrorCode (*weights)(MatColoring,PetscReal**,PetscInt**);
573: };

575: struct _p_MatColoring {
576:   PETSCHEADER(struct _MatColoringOps);
577:   Mat                   mat;
578:   PetscInt              dist;             /* distance of the coloring */
579:   PetscInt              maxcolors;        /* the maximum number of colors returned, maxcolors=1 for MIS */
580:   void                  *data;            /* inner context */
581:   PetscBool             valid;            /* check to see if what is produced is a valid coloring */
582:   MatColoringWeightType weight_type;      /* type of weight computation to be performed */
583:   PetscReal             *user_weights;    /* custom weights and permutation */
584:   PetscInt              *user_lperm;
585:   PetscBool             valid_iscoloring; /* check to see if matcoloring is produced a valid iscoloring */
586: };

588: struct  _p_MatTransposeColoring{
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       *nrows;           /* number of local rows for each color */
595:   PetscInt       currentcolor;     /* color for which function evaluation is being done now */
596:   ISColoringType ctype;            /* IS_COLORING_GLOBAL or IS_COLORING_LOCAL */

598:   PetscInt       *colorforrow,*colorforcol;  /* pointer to rows and columns */
599:   PetscInt       *rows;                      /* lists the local rows for each color (using the local row numbering) */
600:   PetscInt       *den2sp;                    /* maps (row,color) in the dense matrix to index of sparse matrix array a->a */
601:   PetscInt       *columns;                   /* lists the local columns of each color (using global column numbering) */
602:   PetscInt       brows;                      /* number of rows for efficient implementation of MatTransColoringApplyDenToSp() */
603:   PetscInt       *lstart;                    /* array used for loop over row blocks of Csparse */
604: };

606: /*
607:    Null space context for preconditioner/operators
608: */
609: struct _p_MatNullSpace {
610:   PETSCHEADER(int);
611:   PetscBool      has_cnst;
612:   PetscInt       n;
613:   Vec*           vecs;
614:   PetscScalar*   alpha;                 /* for projections */
615:   PetscErrorCode (*remove)(MatNullSpace,Vec,void*);  /* for user provided removal function */
616:   void*          rmctx;                 /* context for remove() function */
617: };

619: /*
620:    Checking zero pivot for LU, ILU preconditioners.
621: */
622: typedef struct {
623:   PetscInt       nshift,nshift_max;
624:   PetscReal      shift_amount,shift_lo,shift_hi,shift_top,shift_fraction;
625:   PetscBool      newshift;
626:   PetscReal      rs;  /* active row sum of abs(offdiagonals) */
627:   PetscScalar    pv;  /* pivot of the active row */
628: } FactorShiftCtx;

630: /*
631:  Used by MatCreateSubMatrices_MPIXAIJ_Local()
632: */
633:  #include <petscctable.h>
634: typedef struct { /* used by MatCreateSubMatrices_MPIAIJ_SingleIS_Local() and MatCreateSubMatrices_MPIAIJ_Local */
635:   PetscInt   id;   /* index of submats, only submats[0] is responsible for deleting some arrays below */
636:   PetscInt   nrqs,nrqr;
637:   PetscInt   **rbuf1,**rbuf2,**rbuf3,**sbuf1,**sbuf2;
638:   PetscInt   **ptr;
639:   PetscInt   *tmp;
640:   PetscInt   *ctr;
641:   PetscInt   *pa; /* proc array */
642:   PetscInt   *req_size,*req_source1,*req_source2;
643:   PetscBool  allcolumns,allrows;
644:   PetscBool  singleis;
645:   PetscInt   *row2proc; /* row to proc map */
646:   PetscInt   nstages;
647: #if defined(PETSC_USE_CTABLE)
648:   PetscTable cmap,rmap;
649:   PetscInt   *cmap_loc,*rmap_loc;
650: #else
651:   PetscInt   *cmap,*rmap;
652: #endif

654:   PetscErrorCode (*destroy)(Mat);
655: } Mat_SubSppt;

657: PETSC_EXTERN PetscErrorCode MatFactorDumpMatrix(Mat);
658: PETSC_INTERN PetscErrorCode MatShift_Basic(Mat,PetscScalar);
659: PETSC_INTERN PetscErrorCode MatSetBlockSizes_Default(Mat,PetscInt,PetscInt);

661: PETSC_STATIC_INLINE PetscErrorCode MatPivotCheck_nz(Mat mat,const MatFactorInfo *info,FactorShiftCtx *sctx,PetscInt row)
662: {
663:   PetscReal _rs   = sctx->rs;
664:   PetscReal _zero = info->zeropivot*_rs;

667:   if (PetscAbsScalar(sctx->pv) <= _zero && !PetscIsNanScalar(sctx->pv)){
668:     /* force |diag| > zeropivot*rs */
669:     if (!sctx->nshift) sctx->shift_amount = info->shiftamount;
670:     else sctx->shift_amount *= 2.0;
671:     sctx->newshift = PETSC_TRUE;
672:     (sctx->nshift)++;
673:   } else {
674:     sctx->newshift = PETSC_FALSE;
675:   }
676:   return(0);
677: }

679: PETSC_STATIC_INLINE PetscErrorCode MatPivotCheck_pd(Mat mat,const MatFactorInfo *info,FactorShiftCtx *sctx,PetscInt row)
680: {
681:   PetscReal _rs   = sctx->rs;
682:   PetscReal _zero = info->zeropivot*_rs;

685:   if (PetscRealPart(sctx->pv) <= _zero && !PetscIsNanScalar(sctx->pv)){
686:     /* force matfactor to be diagonally dominant */
687:     if (sctx->nshift == sctx->nshift_max) {
688:       sctx->shift_fraction = sctx->shift_hi;
689:     } else {
690:       sctx->shift_lo = sctx->shift_fraction;
691:       sctx->shift_fraction = (sctx->shift_hi+sctx->shift_lo)/2.;
692:     }
693:     sctx->shift_amount = sctx->shift_fraction * sctx->shift_top;
694:     sctx->nshift++;
695:     sctx->newshift = PETSC_TRUE;
696:   } else {
697:     sctx->newshift = PETSC_FALSE;
698:   }
699:   return(0);
700: }

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

707:   if (PetscAbsScalar(sctx->pv) <= _zero && !PetscIsNanScalar(sctx->pv)){
708:     sctx->pv          += info->shiftamount;
709:     sctx->shift_amount = 0.0;
710:     sctx->nshift++;
711:   }
712:   sctx->newshift = PETSC_FALSE;
713:   return(0);
714: }

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

722:   sctx->newshift = PETSC_FALSE;
723:   if (PetscAbsScalar(sctx->pv) <= _zero && !PetscIsNanScalar(sctx->pv)) {
724:     if (!mat->erroriffailure) {
725:       PetscInfo3(mat,"Detected zero pivot in factorization in row %D value %g tolerance %g\n",row,(double)PetscAbsScalar(sctx->pv),(double)_zero);
726:       fact->factorerrortype             = MAT_FACTOR_NUMERIC_ZEROPIVOT;
727:       fact->factorerror_zeropivot_value = PetscAbsScalar(sctx->pv);
728:       fact->factorerror_zeropivot_row   = row;
729:     } 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);
730:   }
731:   return(0);
732: }

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

739:   if (info->shifttype == (PetscReal) MAT_SHIFT_NONZERO){
740:     MatPivotCheck_nz(mat,info,sctx,row);
741:   } else if (info->shifttype == (PetscReal) MAT_SHIFT_POSITIVE_DEFINITE){
742:     MatPivotCheck_pd(mat,info,sctx,row);
743:   } else if (info->shifttype == (PetscReal) MAT_SHIFT_INBLOCKS){
744:     MatPivotCheck_inblocks(mat,info,sctx,row);
745:   } else {
746:     MatPivotCheck_none(fact,mat,info,sctx,row);
747:   }
748:   return(0);
749: }

751: /*
752:   Create and initialize a linked list
753:   Input Parameters:
754:     idx_start - starting index of the list
755:     lnk_max   - max value of lnk indicating the end of the list
756:     nlnk      - max length of the list
757:   Output Parameters:
758:     lnk       - list initialized
759:     bt        - PetscBT (bitarray) with all bits set to false
760:     lnk_empty - flg indicating the list is empty
761: */
762: #define PetscLLCreate(idx_start,lnk_max,nlnk,lnk,bt) \
763:   (PetscMalloc1(nlnk,&lnk) || PetscBTCreate(nlnk,&(bt)) || (lnk[idx_start] = lnk_max,0))

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

768: /*
769:   Add an index set into a sorted linked list
770:   Input Parameters:
771:     nidx      - number of input indices
772:     indices   - integer array
773:     idx_start - starting index of the list
774:     lnk       - linked list(an integer array) that is created
775:     bt        - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
776:   output Parameters:
777:     nlnk      - number of newly added indices
778:     lnk       - the sorted(increasing order) linked list containing new and non-redundate entries from indices
779:     bt        - updated PetscBT (bitarray)
780: */
781: #define PetscLLAdd(nidx,indices,idx_start,nlnk,lnk,bt) 0;\
782: {\
783:   PetscInt _k,_entry,_location,_lnkdata;\
784:   nlnk     = 0;\
785:   _lnkdata = idx_start;\
786:   for (_k=0; _k<nidx; _k++){\
787:     _entry = indices[_k];\
788:     if (!PetscBTLookupSet(bt,_entry)){  /* new entry */\
789:       /* search for insertion location */\
790:       /* start from the beginning if _entry < previous _entry */\
791:       if (_k && _entry < _lnkdata) _lnkdata  = idx_start;\
792:       do {\
793:         _location = _lnkdata;\
794:         _lnkdata  = lnk[_location];\
795:       } while (_entry > _lnkdata);\
796:       /* insertion location is found, add entry into lnk */\
797:       lnk[_location] = _entry;\
798:       lnk[_entry]    = _lnkdata;\
799:       nlnk++;\
800:       _lnkdata = _entry; /* next search starts from here if next_entry > _entry */\
801:     }\
802:   }\
803: }

805: /*
806:   Add a permuted index set into a sorted linked list
807:   Input Parameters:
808:     nidx      - number of input indices
809:     indices   - integer array
810:     perm      - permutation of indices
811:     idx_start - starting index of the list
812:     lnk       - linked list(an integer array) that is created
813:     bt        - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
814:   output Parameters:
815:     nlnk      - number of newly added indices
816:     lnk       - the sorted(increasing order) linked list containing new and non-redundate entries from indices
817:     bt        - updated PetscBT (bitarray)
818: */
819: #define PetscLLAddPerm(nidx,indices,perm,idx_start,nlnk,lnk,bt) 0;\
820: {\
821:   PetscInt _k,_entry,_location,_lnkdata;\
822:   nlnk     = 0;\
823:   _lnkdata = idx_start;\
824:   for (_k=0; _k<nidx; _k++){\
825:     _entry = perm[indices[_k]];\
826:     if (!PetscBTLookupSet(bt,_entry)){  /* new entry */\
827:       /* search for insertion location */\
828:       /* start from the beginning if _entry < previous _entry */\
829:       if (_k && _entry < _lnkdata) _lnkdata  = idx_start;\
830:       do {\
831:         _location = _lnkdata;\
832:         _lnkdata  = lnk[_location];\
833:       } while (_entry > _lnkdata);\
834:       /* insertion location is found, add entry into lnk */\
835:       lnk[_location] = _entry;\
836:       lnk[_entry]    = _lnkdata;\
837:       nlnk++;\
838:       _lnkdata = _entry; /* next search starts from here if next_entry > _entry */\
839:     }\
840:   }\
841: }

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

878: #define PetscLLAddSorted_new(nidx,indices,idx_start,lnk_empty,nlnk,lnk,bt) 0; \
879: {\
880:   PetscInt _k,_entry,_location,_lnkdata;\
881:   if (lnk_empty){\
882:     _lnkdata  = idx_start;                      \
883:     for (_k=0; _k<nidx; _k++){                  \
884:       _entry = indices[_k];                             \
885:       PetscBTSet(bt,_entry);  /* mark the new entry */          \
886:           _location = _lnkdata;                                 \
887:           _lnkdata  = lnk[_location];                           \
888:         /* insertion location is found, add entry into lnk */   \
889:         lnk[_location] = _entry;                                \
890:         lnk[_entry]    = _lnkdata;                              \
891:         _lnkdata = _entry; /* next search starts from here */   \
892:     }                                                           \
893:     /*\
894:     lnk[indices[nidx-1]] = lnk[idx_start];\
895:     lnk[idx_start]       = indices[0];\
896:     PetscBTSet(bt,indices[0]);  \
897:     for (_k=1; _k<nidx; _k++){                  \
898:       PetscBTSet(bt,indices[_k]);                                          \
899:       lnk[indices[_k-1]] = indices[_k];                                  \
900:     }                                                           \
901:      */\
902:     nlnk      = nidx;\
903:     lnk_empty = PETSC_FALSE;\
904:   } else {\
905:     nlnk      = 0;                              \
906:     _lnkdata  = idx_start;                      \
907:     for (_k=0; _k<nidx; _k++){                  \
908:       _entry = indices[_k];                             \
909:       if (!PetscBTLookupSet(bt,_entry)){  /* new entry */       \
910:         /* search for insertion location */                     \
911:         do {                                                    \
912:           _location = _lnkdata;                                 \
913:           _lnkdata  = lnk[_location];                           \
914:         } while (_entry > _lnkdata);                            \
915:         /* insertion location is found, add entry into lnk */   \
916:         lnk[_location] = _entry;                                \
917:         lnk[_entry]    = _lnkdata;                              \
918:         nlnk++;                                                 \
919:         _lnkdata = _entry; /* next search starts from here */   \
920:       }                                                         \
921:     }                                                           \
922:   }                                                             \
923: }

925: /*
926:   Add a SORTED index set into a sorted linked list used for LUFactorSymbolic()
927:   Same as PetscLLAddSorted() with an additional operation:
928:        count the number of input indices that are no larger than 'diag'
929:   Input Parameters:
930:     indices   - sorted integer array
931:     idx_start - starting index of the list, index of pivot row
932:     lnk       - linked list(an integer array) that is created
933:     bt        - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
934:     diag      - index of the active row in LUFactorSymbolic
935:     nzbd      - number of input indices with indices <= idx_start
936:     im        - im[idx_start] is initialized as num of nonzero entries in row=idx_start
937:   output Parameters:
938:     nlnk      - number of newly added indices
939:     lnk       - the sorted(increasing order) linked list containing new and non-redundate entries from indices
940:     bt        - updated PetscBT (bitarray)
941:     im        - im[idx_start]: unchanged if diag is not an entry
942:                              : num of entries with indices <= diag if diag is an entry
943: */
944: #define PetscLLAddSortedLU(indices,idx_start,nlnk,lnk,bt,diag,nzbd,im) 0;\
945: {\
946:   PetscInt _k,_entry,_location,_lnkdata,_nidx;\
947:   nlnk     = 0;\
948:   _lnkdata = idx_start;\
949:   _nidx = im[idx_start] - nzbd; /* num of entries with idx_start < index <= diag */\
950:   for (_k=0; _k<_nidx; _k++){\
951:     _entry = indices[_k];\
952:     nzbd++;\
953:     if ( _entry== diag) im[idx_start] = nzbd;\
954:     if (!PetscBTLookupSet(bt,_entry)){  /* new entry */\
955:       /* search for insertion location */\
956:       do {\
957:         _location = _lnkdata;\
958:         _lnkdata  = lnk[_location];\
959:       } while (_entry > _lnkdata);\
960:       /* insertion location is found, add entry into lnk */\
961:       lnk[_location] = _entry;\
962:       lnk[_entry]    = _lnkdata;\
963:       nlnk++;\
964:       _lnkdata = _entry; /* next search starts from here */\
965:     }\
966:   }\
967: }

969: /*
970:   Copy data on the list into an array, then initialize the list
971:   Input Parameters:
972:     idx_start - starting index of the list
973:     lnk_max   - max value of lnk indicating the end of the list
974:     nlnk      - number of data on the list to be copied
975:     lnk       - linked list
976:     bt        - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
977:   output Parameters:
978:     indices   - array that contains the copied data
979:     lnk       - linked list that is cleaned and initialize
980:     bt        - PetscBT (bitarray) with all bits set to false
981: */
982: #define PetscLLClean(idx_start,lnk_max,nlnk,lnk,indices,bt) 0;\
983: {\
984:   PetscInt _j,_idx=idx_start;\
985:   for (_j=0; _j<nlnk; _j++){\
986:     _idx = lnk[_idx];\
987:     indices[_j] = _idx;\
988:     PetscBTClear(bt,_idx);\
989:   }\
990:   lnk[idx_start] = lnk_max;\
991: }
992: /*
993:   Free memories used by the list
994: */
995: #define PetscLLDestroy(lnk,bt) (PetscFree(lnk) || PetscBTDestroy(&(bt)))

997: /* Routines below are used for incomplete matrix factorization */
998: /*
999:   Create and initialize a linked list and its levels
1000:   Input Parameters:
1001:     idx_start - starting index of the list
1002:     lnk_max   - max value of lnk indicating the end of the list
1003:     nlnk      - max length of the list
1004:   Output Parameters:
1005:     lnk       - list initialized
1006:     lnk_lvl   - array of size nlnk for storing levels of lnk
1007:     bt        - PetscBT (bitarray) with all bits set to false
1008: */
1009: #define PetscIncompleteLLCreate(idx_start,lnk_max,nlnk,lnk,lnk_lvl,bt)\
1010:   (PetscIntMultError(2,nlnk,NULL) || PetscMalloc1(2*nlnk,&lnk) || PetscBTCreate(nlnk,&(bt)) || (lnk[idx_start] = lnk_max,lnk_lvl = lnk + nlnk,0))

1012: /*
1013:   Initialize a sorted linked list used for ILU and ICC
1014:   Input Parameters:
1015:     nidx      - number of input idx
1016:     idx       - integer array used for storing column indices
1017:     idx_start - starting index of the list
1018:     perm      - indices of an IS
1019:     lnk       - linked list(an integer array) that is created
1020:     lnklvl    - levels of lnk
1021:     bt        - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
1022:   output Parameters:
1023:     nlnk     - number of newly added idx
1024:     lnk      - the sorted(increasing order) linked list containing new and non-redundate entries from idx
1025:     lnklvl   - levels of lnk
1026:     bt       - updated PetscBT (bitarray)
1027: */
1028: #define PetscIncompleteLLInit(nidx,idx,idx_start,perm,nlnk,lnk,lnklvl,bt) 0;\
1029: {\
1030:   PetscInt _k,_entry,_location,_lnkdata;\
1031:   nlnk     = 0;\
1032:   _lnkdata = idx_start;\
1033:   for (_k=0; _k<nidx; _k++){\
1034:     _entry = perm[idx[_k]];\
1035:     if (!PetscBTLookupSet(bt,_entry)){  /* new entry */\
1036:       /* search for insertion location */\
1037:       if (_k && _entry < _lnkdata) _lnkdata  = idx_start;\
1038:       do {\
1039:         _location = _lnkdata;\
1040:         _lnkdata  = lnk[_location];\
1041:       } while (_entry > _lnkdata);\
1042:       /* insertion location is found, add entry into lnk */\
1043:       lnk[_location]  = _entry;\
1044:       lnk[_entry]     = _lnkdata;\
1045:       lnklvl[_entry] = 0;\
1046:       nlnk++;\
1047:       _lnkdata = _entry; /* next search starts from here if next_entry > _entry */\
1048:     }\
1049:   }\
1050: }

1052: /*
1053:   Add a SORTED index set into a sorted linked list for ILU
1054:   Input Parameters:
1055:     nidx      - number of input indices
1056:     idx       - sorted integer array used for storing column indices
1057:     level     - level of fill, e.g., ICC(level)
1058:     idxlvl    - level of idx
1059:     idx_start - starting index of the list
1060:     lnk       - linked list(an integer array) that is created
1061:     lnklvl    - levels of lnk
1062:     bt        - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
1063:     prow      - the row number of idx
1064:   output Parameters:
1065:     nlnk     - number of newly added idx
1066:     lnk      - the sorted(increasing order) linked list containing new and non-redundate entries from idx
1067:     lnklvl   - levels of lnk
1068:     bt       - updated PetscBT (bitarray)

1070:   Note: the level of factor(i,j) is set as lvl(i,j) = min{ lvl(i,j), lvl(i,prow)+lvl(prow,j)+1)
1071:         where idx = non-zero columns of U(prow,prow+1:n-1), prow<i
1072: */
1073: #define PetscILULLAddSorted(nidx,idx,level,idxlvl,idx_start,nlnk,lnk,lnklvl,bt,lnklvl_prow) 0;\
1074: {\
1075:   PetscInt _k,_entry,_location,_lnkdata,_incrlev,_lnklvl_prow=lnklvl[prow];\
1076:   nlnk     = 0;\
1077:   _lnkdata = idx_start;\
1078:   for (_k=0; _k<nidx; _k++){\
1079:     _incrlev = idxlvl[_k] + _lnklvl_prow + 1;\
1080:     if (_incrlev > level) continue;\
1081:     _entry = idx[_k];\
1082:     if (!PetscBTLookupSet(bt,_entry)){  /* new entry */\
1083:       /* search for insertion location */\
1084:       do {\
1085:         _location = _lnkdata;\
1086:         _lnkdata  = lnk[_location];\
1087:       } while (_entry > _lnkdata);\
1088:       /* insertion location is found, add entry into lnk */\
1089:       lnk[_location]  = _entry;\
1090:       lnk[_entry]     = _lnkdata;\
1091:       lnklvl[_entry] = _incrlev;\
1092:       nlnk++;\
1093:       _lnkdata = _entry; /* next search starts from here if next_entry > _entry */\
1094:     } else { /* existing entry: update lnklvl */\
1095:       if (lnklvl[_entry] > _incrlev) lnklvl[_entry] = _incrlev;\
1096:     }\
1097:   }\
1098: }

1100: /*
1101:   Add a index set into a sorted linked list
1102:   Input Parameters:
1103:     nidx      - number of input idx
1104:     idx   - integer array used for storing column indices
1105:     level     - level of fill, e.g., ICC(level)
1106:     idxlvl - level of idx
1107:     idx_start - starting index of the list
1108:     lnk       - linked list(an integer array) that is created
1109:     lnklvl   - levels of lnk
1110:     bt        - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
1111:   output Parameters:
1112:     nlnk      - number of newly added idx
1113:     lnk       - the sorted(increasing order) linked list containing new and non-redundate entries from idx
1114:     lnklvl   - levels of lnk
1115:     bt        - updated PetscBT (bitarray)
1116: */
1117: #define PetscIncompleteLLAdd(nidx,idx,level,idxlvl,idx_start,nlnk,lnk,lnklvl,bt) 0;\
1118: {\
1119:   PetscInt _k,_entry,_location,_lnkdata,_incrlev;\
1120:   nlnk     = 0;\
1121:   _lnkdata = idx_start;\
1122:   for (_k=0; _k<nidx; _k++){\
1123:     _incrlev = idxlvl[_k] + 1;\
1124:     if (_incrlev > level) continue;\
1125:     _entry = idx[_k];\
1126:     if (!PetscBTLookupSet(bt,_entry)){  /* new entry */\
1127:       /* search for insertion location */\
1128:       if (_k && _entry < _lnkdata) _lnkdata  = idx_start;\
1129:       do {\
1130:         _location = _lnkdata;\
1131:         _lnkdata  = lnk[_location];\
1132:       } while (_entry > _lnkdata);\
1133:       /* insertion location is found, add entry into lnk */\
1134:       lnk[_location]  = _entry;\
1135:       lnk[_entry]     = _lnkdata;\
1136:       lnklvl[_entry] = _incrlev;\
1137:       nlnk++;\
1138:       _lnkdata = _entry; /* next search starts from here if next_entry > _entry */\
1139:     } else { /* existing entry: update lnklvl */\
1140:       if (lnklvl[_entry] > _incrlev) lnklvl[_entry] = _incrlev;\
1141:     }\
1142:   }\
1143: }

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

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

1236: /*
1237:   Copy data on the list into an array, then initialize the list
1238:   Input Parameters:
1239:     idx_start - starting index of the list
1240:     lnk_max   - max value of lnk indicating the end of the list
1241:     nlnk      - number of data on the list to be copied
1242:     lnk       - linked list
1243:     lnklvl    - level of lnk
1244:     bt        - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
1245:   output Parameters:
1246:     indices - array that contains the copied data
1247:     lnk     - linked list that is cleaned and initialize
1248:     lnklvl  - level of lnk that is reinitialized
1249:     bt      - PetscBT (bitarray) with all bits set to false
1250: */
1251: #define PetscIncompleteLLClean(idx_start,lnk_max,nlnk,lnk,lnklvl,indices,indiceslvl,bt) 0;\
1252: do {\
1253:   PetscInt _j,_idx=idx_start;\
1254:   for (_j=0; _j<nlnk; _j++){\
1255:     _idx = lnk[_idx];\
1256:     *(indices+_j) = _idx;\
1257:     *(indiceslvl+_j) = lnklvl[_idx];\
1258:     lnklvl[_idx] = -1;\
1259:     PetscBTClear(bt,_idx);\
1260:   }\
1261:   lnk[idx_start] = lnk_max;\
1262: } while (0)
1263: /*
1264:   Free memories used by the list
1265: */
1266: #define PetscIncompleteLLDestroy(lnk,bt) (PetscFree(lnk) || PetscBTDestroy(&(bt)))

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

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

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

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

1289:       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
1290:       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
1291:       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
1292:       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.
1293:       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
1294:       to each other so memory access is much better than using the big array.

1296:   Example:
1297:      nlnk_max=5, lnk_max=36:
1298:      Initial list: [0, 0 | 36, 2 | 0, 0 | 0, 0 | 0, 0 | 0, 0 | 0, 0]
1299:      here, head_node has index 2 with value lnk[2]=lnk_max=36,
1300:            0-th entry is used to store the number of entries in the list,
1301:      The initial lnk represents head -> tail(marked by 36) with number of entries = lnk[0]=0.

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

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

1311:   Input Parameters:
1312:     nlnk_max  - max length of the list
1313:     lnk_max   - max value of the entries
1314:   Output Parameters:
1315:     lnk       - list created and initialized
1316:     bt        - PetscBT (bitarray) with all bits set to false. Note: bt has size lnk_max, not nln_max!
1317: */
1318: PETSC_STATIC_INLINE PetscErrorCode PetscLLCondensedCreate(PetscInt nlnk_max,PetscInt lnk_max,PetscInt **lnk,PetscBT *bt)
1319: {
1321:   PetscInt       *llnk,lsize = 0;

1324:   PetscIntMultError(2,nlnk_max+2,&lsize);
1325:   PetscMalloc1(lsize,lnk);
1326:   PetscBTCreate(lnk_max,bt);
1327:   llnk = *lnk;
1328:   llnk[0] = 0;         /* number of entries on the list */
1329:   llnk[2] = lnk_max;   /* value in the head node */
1330:   llnk[3] = 2;         /* next for the head node */
1331:   return(0);
1332: }

1334: /*
1335:   Add a SORTED ascending index set into a sorted linked list. See PetscLLCondensedCreate() for detailed description.
1336:   Input Parameters:
1337:     nidx      - number of input indices
1338:     indices   - sorted integer array
1339:     lnk       - condensed linked list(an integer array) that is created
1340:     bt        - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
1341:   output Parameters:
1342:     lnk       - the sorted(increasing order) linked list containing previous and newly added non-redundate indices
1343:     bt        - updated PetscBT (bitarray)
1344: */
1345: PETSC_STATIC_INLINE PetscErrorCode PetscLLCondensedAddSorted(PetscInt nidx,const PetscInt indices[],PetscInt lnk[],PetscBT bt)
1346: {
1347:   PetscInt _k,_entry,_location,_next,_lnkdata,_nlnk,_newnode;

1350:   _nlnk     = lnk[0]; /* num of entries on the input lnk */
1351:   _location = 2; /* head */
1352:     for (_k=0; _k<nidx; _k++){
1353:       _entry = indices[_k];
1354:       if (!PetscBTLookupSet(bt,_entry)){  /* new entry */
1355:         /* search for insertion location */
1356:         do {
1357:           _next     = _location + 1; /* link from previous node to next node */
1358:           _location = lnk[_next];    /* idx of next node */
1359:           _lnkdata  = lnk[_location];/* value of next node */
1360:         } while (_entry > _lnkdata);
1361:         /* insertion location is found, add entry into lnk */
1362:         _newnode        = 2*(_nlnk+2);   /* index for this new node */
1363:         lnk[_next]      = _newnode;      /* connect previous node to the new node */
1364:         lnk[_newnode]   = _entry;        /* set value of the new node */
1365:         lnk[_newnode+1] = _location;     /* connect new node to next node */
1366:         _location       = _newnode;      /* next search starts from the new node */
1367:         _nlnk++;
1368:       }   \
1369:     }\
1370:   lnk[0]   = _nlnk;   /* number of entries in the list */
1371:   return(0);
1372: }

1374: PETSC_STATIC_INLINE PetscErrorCode PetscLLCondensedClean(PetscInt lnk_max,PetscInt nidx,PetscInt *indices,PetscInt lnk[],PetscBT bt)
1375: {
1377:   PetscInt       _k,_next,_nlnk;

1380:   _next = lnk[3];       /* head node */
1381:   _nlnk = lnk[0];       /* num of entries on the list */
1382:   for (_k=0; _k<_nlnk; _k++){
1383:     indices[_k] = lnk[_next];
1384:     _next       = lnk[_next + 1];
1385:     PetscBTClear(bt,indices[_k]);
1386:   }
1387:   lnk[0] = 0;          /* num of entries on the list */
1388:   lnk[2] = lnk_max;    /* initialize head node */
1389:   lnk[3] = 2;          /* head node */
1390:   return(0);
1391: }

1393: PETSC_STATIC_INLINE PetscErrorCode PetscLLCondensedView(PetscInt *lnk)
1394: {
1396:   PetscInt       k;

1399:   PetscPrintf(PETSC_COMM_SELF,"LLCondensed of size %D, (val,  next)\n",lnk[0]);
1400:   for (k=2; k< lnk[0]+2; k++){
1401:     PetscPrintf(PETSC_COMM_SELF," %D: (%D, %D)\n",2*k,lnk[2*k],lnk[2*k+1]);
1402:   }
1403:   return(0);
1404: }

1406: /*
1407:   Free memories used by the list
1408: */
1409: PETSC_STATIC_INLINE PetscErrorCode PetscLLCondensedDestroy(PetscInt *lnk,PetscBT bt)
1410: {

1414:   PetscFree(lnk);
1415:   PetscBTDestroy(&bt);
1416:   return(0);
1417: }

1419: /* -------------------------------------------------------------------------------------------------------*/
1420: /*
1421:  Same as PetscLLCondensedCreate(), but does not use non-scalable O(lnk_max) bitarray
1422:   Input Parameters:
1423:     nlnk_max  - max length of the list
1424:   Output Parameters:
1425:     lnk       - list created and initialized
1426: */
1427: PETSC_STATIC_INLINE PetscErrorCode PetscLLCondensedCreate_Scalable(PetscInt nlnk_max,PetscInt **lnk)
1428: {
1430:   PetscInt       *llnk,lsize = 0;

1433:   PetscIntMultError(2,nlnk_max+2,&lsize);
1434:   PetscMalloc1(lsize,lnk);
1435:   llnk = *lnk;
1436:   llnk[0] = 0;               /* number of entries on the list */
1437:   llnk[2] = PETSC_MAX_INT;   /* value in the head node */
1438:   llnk[3] = 2;               /* next for the head node */
1439:   return(0);
1440: }

1442: PETSC_STATIC_INLINE PetscErrorCode PetscLLCondensedExpand_Scalable(PetscInt nlnk_max,PetscInt **lnk)
1443: {
1445:   PetscInt       lsize = 0;

1448:   PetscIntMultError(2,nlnk_max+2,&lsize);
1449:   PetscRealloc(lsize*sizeof(PetscInt),lnk);
1450:   return(0);
1451: }

1453: PETSC_STATIC_INLINE PetscErrorCode PetscLLCondensedAddSorted_Scalable(PetscInt nidx,const PetscInt indices[],PetscInt lnk[])
1454: {
1455:   PetscInt _k,_entry,_location,_next,_lnkdata,_nlnk,_newnode;
1456:   _nlnk     = lnk[0]; /* num of entries on the input lnk */
1457:   _location = 2; /* head */ \
1458:     for (_k=0; _k<nidx; _k++){
1459:       _entry = indices[_k];
1460:       /* search for insertion location */
1461:       do {
1462:         _next     = _location + 1; /* link from previous node to next node */
1463:         _location = lnk[_next];    /* idx of next node */
1464:         _lnkdata  = lnk[_location];/* value of next node */
1465:       } while (_entry > _lnkdata);
1466:       if (_entry < _lnkdata) {
1467:         /* insertion location is found, add entry into lnk */
1468:         _newnode        = 2*(_nlnk+2);   /* index for this new node */
1469:         lnk[_next]      = _newnode;      /* connect previous node to the new node */
1470:         lnk[_newnode]   = _entry;        /* set value of the new node */
1471:         lnk[_newnode+1] = _location;     /* connect new node to next node */
1472:         _location       = _newnode;      /* next search starts from the new node */
1473:         _nlnk++;
1474:       }
1475:     }
1476:   lnk[0]   = _nlnk;   /* number of entries in the list */
1477:   return 0;
1478: }

1480: PETSC_STATIC_INLINE PetscErrorCode PetscLLCondensedClean_Scalable(PetscInt nidx,PetscInt *indices,PetscInt *lnk)
1481: {
1482:   PetscInt _k,_next,_nlnk;
1483:   _next = lnk[3];       /* head node */
1484:   _nlnk = lnk[0];
1485:   for (_k=0; _k<_nlnk; _k++){
1486:     indices[_k] = lnk[_next];
1487:     _next       = lnk[_next + 1];
1488:   }
1489:   lnk[0] = 0;          /* num of entries on the list */
1490:   lnk[3] = 2;          /* head node */
1491:   return 0;
1492: }

1494: PETSC_STATIC_INLINE PetscErrorCode PetscLLCondensedDestroy_Scalable(PetscInt *lnk)
1495: {
1496:   return PetscFree(lnk);
1497: }

1499: /* -------------------------------------------------------------------------------------------------------*/
1500: /*
1501:       lnk[0]   number of links
1502:       lnk[1]   number of entries
1503:       lnk[3n]  value
1504:       lnk[3n+1] len
1505:       lnk[3n+2] link to next value

1507:       The next three are always the first link

1509:       lnk[3]    PETSC_MIN_INT+1
1510:       lnk[4]    1
1511:       lnk[5]    link to first real entry

1513:       The next three are always the last link

1515:       lnk[6]    PETSC_MAX_INT - 1
1516:       lnk[7]    1
1517:       lnk[8]    next valid link (this is the same as lnk[0] but without the decreases)
1518: */

1520: PETSC_STATIC_INLINE PetscErrorCode PetscLLCondensedCreate_fast(PetscInt nlnk_max,PetscInt **lnk)
1521: {
1523:   PetscInt       *llnk,lsize = 0;

1526:   PetscIntMultError(3,nlnk_max+3,&lsize);
1527:   PetscMalloc1(lsize,lnk);
1528:   llnk = *lnk;
1529:   llnk[0] = 0;   /* nlnk: number of entries on the list */
1530:   llnk[1] = 0;          /* number of integer entries represented in list */
1531:   llnk[3] = PETSC_MIN_INT+1;   /* value in the first node */
1532:   llnk[4] = 1;           /* count for the first node */
1533:   llnk[5] = 6;         /* next for the first node */
1534:   llnk[6] = PETSC_MAX_INT-1;   /* value in the last node */
1535:   llnk[7] = 1;           /* count for the last node */
1536:   llnk[8] = 0;         /* next valid node to be used */
1537:   return(0);
1538: }

1540: PETSC_STATIC_INLINE PetscErrorCode PetscLLCondensedAddSorted_fast(PetscInt nidx,const PetscInt indices[],PetscInt lnk[])
1541: {
1542:   PetscInt k,entry,prev,next;
1543:   prev      = 3;      /* first value */
1544:   next      = lnk[prev+2];
1545:   for (k=0; k<nidx; k++){
1546:     entry = indices[k];
1547:     /* search for insertion location */
1548:     while (entry >= lnk[next]) {
1549:       prev = next;
1550:       next = lnk[next+2];
1551:     }
1552:     /* entry is in range of previous list */
1553:     if (entry < lnk[prev]+lnk[prev+1]) continue;
1554:     lnk[1]++;
1555:     /* entry is right after previous list */
1556:     if (entry == lnk[prev]+lnk[prev+1]) {
1557:       lnk[prev+1]++;
1558:       if (lnk[next] == entry+1) { /* combine two contiguous strings */
1559:         lnk[prev+1] += lnk[next+1];
1560:         lnk[prev+2]  = lnk[next+2];
1561:         next         = lnk[next+2];
1562:         lnk[0]--;
1563:       }
1564:       continue;
1565:     }
1566:     /* entry is right before next list */
1567:     if (entry == lnk[next]-1) {
1568:       lnk[next]--;
1569:       lnk[next+1]++;
1570:       prev = next;
1571:       next = lnk[prev+2];
1572:       continue;
1573:     }
1574:     /*  add entry into lnk */
1575:     lnk[prev+2]    = 3*((lnk[8]++)+3);      /* connect previous node to the new node */
1576:     prev           = lnk[prev+2];
1577:     lnk[prev]      = entry;        /* set value of the new node */
1578:     lnk[prev+1]    = 1;             /* number of values in contiguous string is one to start */
1579:     lnk[prev+2]    = next;          /* connect new node to next node */
1580:     lnk[0]++;
1581:   }
1582:   return 0;
1583: }

1585: PETSC_STATIC_INLINE PetscErrorCode PetscLLCondensedClean_fast(PetscInt nidx,PetscInt *indices,PetscInt *lnk)
1586: {
1587:   PetscInt _k,_next,_nlnk,cnt,j;
1588:   _next = lnk[5];       /* first node */
1589:   _nlnk = lnk[0];
1590:   cnt   = 0;
1591:   for (_k=0; _k<_nlnk; _k++){
1592:     for (j=0; j<lnk[_next+1]; j++) {
1593:       indices[cnt++] = lnk[_next] + j;
1594:     }
1595:     _next       = lnk[_next + 2];
1596:   }
1597:   lnk[0] = 0;   /* nlnk: number of links */
1598:   lnk[1] = 0;          /* number of integer entries represented in list */
1599:   lnk[3] = PETSC_MIN_INT+1;   /* value in the first node */
1600:   lnk[4] = 1;           /* count for the first node */
1601:   lnk[5] = 6;         /* next for the first node */
1602:   lnk[6] = PETSC_MAX_INT-1;   /* value in the last node */
1603:   lnk[7] = 1;           /* count for the last node */
1604:   lnk[8] = 0;         /* next valid location to make link */
1605:   return 0;
1606: }

1608: PETSC_STATIC_INLINE PetscErrorCode PetscLLCondensedView_fast(PetscInt *lnk)
1609: {
1610:   PetscInt k,next,nlnk;
1611:   next = lnk[5];       /* first node */
1612:   nlnk = lnk[0];
1613:   for (k=0; k<nlnk; k++){
1614: #if 0                           /* Debugging code */
1615:     printf("%d value %d len %d next %d\n",next,lnk[next],lnk[next+1],lnk[next+2]);
1616: #endif
1617:     next = lnk[next + 2];
1618:   }
1619:   return 0;
1620: }

1622: PETSC_STATIC_INLINE PetscErrorCode PetscLLCondensedDestroy_fast(PetscInt *lnk)
1623: {
1624:   return PetscFree(lnk);
1625: }

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

1630: PETSC_EXTERN PetscLogEvent MAT_Mult;
1631: PETSC_EXTERN PetscLogEvent MAT_MultMatrixFree;
1632: PETSC_EXTERN PetscLogEvent MAT_Mults;
1633: PETSC_EXTERN PetscLogEvent MAT_MultConstrained;
1634: PETSC_EXTERN PetscLogEvent MAT_MultAdd;
1635: PETSC_EXTERN PetscLogEvent MAT_MultTranspose;
1636: PETSC_EXTERN PetscLogEvent MAT_MultTransposeConstrained;
1637: PETSC_EXTERN PetscLogEvent MAT_MultTransposeAdd;
1638: PETSC_EXTERN PetscLogEvent MAT_Solve;
1639: PETSC_EXTERN PetscLogEvent MAT_Solves;
1640: PETSC_EXTERN PetscLogEvent MAT_SolveAdd;
1641: PETSC_EXTERN PetscLogEvent MAT_SolveTranspose;
1642: PETSC_EXTERN PetscLogEvent MAT_SolveTransposeAdd;
1643: PETSC_EXTERN PetscLogEvent MAT_SOR;
1644: PETSC_EXTERN PetscLogEvent MAT_ForwardSolve;
1645: PETSC_EXTERN PetscLogEvent MAT_BackwardSolve;
1646: PETSC_EXTERN PetscLogEvent MAT_LUFactor;
1647: PETSC_EXTERN PetscLogEvent MAT_LUFactorSymbolic;
1648: PETSC_EXTERN PetscLogEvent MAT_LUFactorNumeric;
1649: PETSC_EXTERN PetscLogEvent MAT_CholeskyFactor;
1650: PETSC_EXTERN PetscLogEvent MAT_CholeskyFactorSymbolic;
1651: PETSC_EXTERN PetscLogEvent MAT_CholeskyFactorNumeric;
1652: PETSC_EXTERN PetscLogEvent MAT_ILUFactor;
1653: PETSC_EXTERN PetscLogEvent MAT_ILUFactorSymbolic;
1654: PETSC_EXTERN PetscLogEvent MAT_ICCFactorSymbolic;
1655: PETSC_EXTERN PetscLogEvent MAT_Copy;
1656: PETSC_EXTERN PetscLogEvent MAT_Convert;
1657: PETSC_EXTERN PetscLogEvent MAT_Scale;
1658: PETSC_EXTERN PetscLogEvent MAT_AssemblyBegin;
1659: PETSC_EXTERN PetscLogEvent MAT_AssemblyEnd;
1660: PETSC_EXTERN PetscLogEvent MAT_SetValues;
1661: PETSC_EXTERN PetscLogEvent MAT_GetValues;
1662: PETSC_EXTERN PetscLogEvent MAT_GetRow;
1663: PETSC_EXTERN PetscLogEvent MAT_GetRowIJ;
1664: PETSC_EXTERN PetscLogEvent MAT_CreateSubMats;
1665: PETSC_EXTERN PetscLogEvent MAT_GetColoring;
1666: PETSC_EXTERN PetscLogEvent MAT_GetOrdering;
1667: PETSC_EXTERN PetscLogEvent MAT_RedundantMat;
1668: PETSC_EXTERN PetscLogEvent MAT_IncreaseOverlap;
1669: PETSC_EXTERN PetscLogEvent MAT_Partitioning;
1670: PETSC_EXTERN PetscLogEvent MAT_PartitioningND;
1671: PETSC_EXTERN PetscLogEvent MAT_Coarsen;
1672: PETSC_EXTERN PetscLogEvent MAT_ZeroEntries;
1673: PETSC_EXTERN PetscLogEvent MAT_Load;
1674: PETSC_EXTERN PetscLogEvent MAT_View;
1675: PETSC_EXTERN PetscLogEvent MAT_AXPY;
1676: PETSC_EXTERN PetscLogEvent MAT_FDColoringCreate;
1677: PETSC_EXTERN PetscLogEvent MAT_TransposeColoringCreate;
1678: PETSC_EXTERN PetscLogEvent MAT_FDColoringSetUp;
1679: PETSC_EXTERN PetscLogEvent MAT_FDColoringApply;
1680: PETSC_EXTERN PetscLogEvent MAT_Transpose;
1681: PETSC_EXTERN PetscLogEvent MAT_FDColoringFunction;
1682: PETSC_EXTERN PetscLogEvent MAT_CreateSubMat;
1683: PETSC_EXTERN PetscLogEvent MAT_MatMult;
1684: PETSC_EXTERN PetscLogEvent MAT_MatSolve;
1685: PETSC_EXTERN PetscLogEvent MAT_MatTrSolve;
1686: PETSC_EXTERN PetscLogEvent MAT_MatMultSymbolic;
1687: PETSC_EXTERN PetscLogEvent MAT_MatMultNumeric;
1688: PETSC_EXTERN PetscLogEvent MAT_Getlocalmatcondensed;
1689: PETSC_EXTERN PetscLogEvent MAT_GetBrowsOfAcols;
1690: PETSC_EXTERN PetscLogEvent MAT_GetBrowsOfAocols;
1691: PETSC_EXTERN PetscLogEvent MAT_PtAP;
1692: PETSC_EXTERN PetscLogEvent MAT_PtAPSymbolic;
1693: PETSC_EXTERN PetscLogEvent MAT_PtAPNumeric;
1694: PETSC_EXTERN PetscLogEvent MAT_Seqstompinum;
1695: PETSC_EXTERN PetscLogEvent MAT_Seqstompisym;
1696: PETSC_EXTERN PetscLogEvent MAT_Seqstompi;
1697: PETSC_EXTERN PetscLogEvent MAT_Getlocalmat;
1698: PETSC_EXTERN PetscLogEvent MAT_RARt;
1699: PETSC_EXTERN PetscLogEvent MAT_RARtSymbolic;
1700: PETSC_EXTERN PetscLogEvent MAT_RARtNumeric;
1701: PETSC_EXTERN PetscLogEvent MAT_MatTransposeMult;
1702: PETSC_EXTERN PetscLogEvent MAT_MatTransposeMultSymbolic;
1703: PETSC_EXTERN PetscLogEvent MAT_MatTransposeMultNumeric;
1704: PETSC_EXTERN PetscLogEvent MAT_TransposeMatMult;
1705: PETSC_EXTERN PetscLogEvent MAT_TransposeMatMultSymbolic;
1706: PETSC_EXTERN PetscLogEvent MAT_TransposeMatMultNumeric;
1707: PETSC_EXTERN PetscLogEvent MAT_MatMatMult;
1708: PETSC_EXTERN PetscLogEvent MAT_MatMatMultSymbolic;
1709: PETSC_EXTERN PetscLogEvent MAT_MatMatMultNumeric;
1710: PETSC_EXTERN PetscLogEvent MAT_Applypapt;
1711: PETSC_EXTERN PetscLogEvent MAT_Applypapt_symbolic;
1712: PETSC_EXTERN PetscLogEvent MAT_Applypapt_numeric;
1713: PETSC_EXTERN PetscLogEvent MAT_Getsymtranspose;
1714: PETSC_EXTERN PetscLogEvent MAT_Getsymtransreduced;
1715: PETSC_EXTERN PetscLogEvent MAT_GetSequentialNonzeroStructure;
1716: PETSC_EXTERN PetscLogEvent MATMFFD_Mult;
1717: PETSC_EXTERN PetscLogEvent MAT_GetMultiProcBlock;
1718: PETSC_EXTERN PetscLogEvent MAT_CUSPARSECopyToGPU;
1719: PETSC_EXTERN PetscLogEvent MAT_SetValuesBatch;
1720: PETSC_EXTERN PetscLogEvent MAT_ViennaCLCopyToGPU;
1721: PETSC_EXTERN PetscLogEvent MAT_Merge;
1722: PETSC_EXTERN PetscLogEvent MAT_Residual;
1723: PETSC_EXTERN PetscLogEvent MAT_SetRandom;
1724: PETSC_EXTERN PetscLogEvent MATCOLORING_Apply;
1725: PETSC_EXTERN PetscLogEvent MATCOLORING_Comm;
1726: PETSC_EXTERN PetscLogEvent MATCOLORING_Local;
1727: PETSC_EXTERN PetscLogEvent MATCOLORING_ISCreate;
1728: PETSC_EXTERN PetscLogEvent MATCOLORING_SetUp;
1729: PETSC_EXTERN PetscLogEvent MATCOLORING_Weights;

1731: #endif