Actual source code: matimpl.h

petsc-master 2019-08-22
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 MatPtAP_Basic(Mat,Mat,MatReuse,PetscReal,Mat*);
426: PETSC_INTERN PetscErrorCode MatAXPY_Basic(Mat,PetscScalar,Mat,MatStructure);
427: PETSC_INTERN PetscErrorCode MatAXPY_BasicWithPreallocation(Mat,Mat,PetscScalar,Mat,MatStructure);
428: PETSC_INTERN PetscErrorCode MatAXPY_Basic_Preallocate(Mat,Mat,Mat*);

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

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

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

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

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

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

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

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

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

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

508:     ncolors = 4;

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

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

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

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

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

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

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

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

589: struct  _p_MatTransposeColoring{
590:   PETSCHEADER(int);
591:   PetscInt       M,N,m;            /* total rows, columns; local rows */
592:   PetscInt       rstart;           /* first row owned by local processor */
593:   PetscInt       ncolors;          /* number of colors */
594:   PetscInt       *ncolumns;        /* number of local columns for a color */
595:   PetscInt       *nrows;           /* number of local rows for each color */
596:   PetscInt       currentcolor;     /* color for which function evaluation is being done now */
597:   ISColoringType ctype;            /* IS_COLORING_GLOBAL or IS_COLORING_LOCAL */

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

1394: PETSC_STATIC_INLINE PetscErrorCode PetscLLCondensedView(PetscInt *lnk)
1395: {
1397:   PetscInt       k;

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

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

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

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

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

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

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

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

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

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

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

1508:       The next three are always the first link

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

1514:       The next three are always the last link

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

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

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

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

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

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

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

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

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

1734: #endif