Actual source code: matimpl.h

petsc-master 2019-05-24
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:   MatInfo                info;             /* matrix information */
389:   InsertMode             insertmode;       /* have values been inserted in matrix or added? */
390:   MatStash               stash,bstash;     /* used for assembling off-proc mat emements */
391:   MatNullSpace           nullsp;           /* null space (operator is singular) */
392:   MatNullSpace           transnullsp;      /* null space of transpose of operator */
393:   MatNullSpace           nearnullsp;       /* near null space to be used by multigrid methods */
394:   PetscInt               congruentlayouts; /* are the rows and columns layouts congruent? */
395:   PetscBool              preallocated;
396:   MatStencilInfo         stencil;          /* information for structured grid */
397:   PetscBool              symmetric,hermitian,structurally_symmetric,spd;
398:   PetscBool              symmetric_set,hermitian_set,structurally_symmetric_set,spd_set; /* if true, then corresponding flag is correct*/
399:   PetscBool              symmetric_eternal;
400:   PetscBool              nooffprocentries,nooffproczerorows;
401:   PetscBool              assembly_subset; /* set by MAT_SUBSET_OFF_PROC_ENTRIES */
402:   PetscBool              submat_singleis; /* for efficient PCSetUP_ASM() */
403:   PetscBool              structure_only;
404: #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_CUDA)
405:   PetscOffloadFlag       valid_GPU_matrix; /* flag pointing to the matrix on the gpu*/
406:   PetscBool              pinnedtocpu;
407: #endif
408:   void                   *spptr;          /* pointer for special library like SuperLU */
409:   char                   *solvertype;
410:   PetscBool              checksymmetryonassembly,checknullspaceonassembly;
411:   PetscReal              checksymmetrytol;
412:   Mat                    schur;             /* Schur complement matrix */
413:   MatFactorSchurStatus   schur_status;      /* status of the Schur complement matrix */
414:   Mat_Redundant          *redundant;        /* used by MatCreateRedundantMatrix() */
415:   PetscBool              erroriffailure;    /* Generate an error if detected (for example a zero pivot) instead of returning */
416:   MatFactorError         factorerrortype;               /* type of error in factorization */
417:   PetscReal              factorerror_zeropivot_value;   /* If numerical zero pivot was detected this is the computed value */
418:   PetscInt               factorerror_zeropivot_row;     /* Row where zero pivot was detected */
419:   PetscInt               nblocks,*bsizes;   /* support for MatSetVariableBlockSizes() */
420:   char                   *defaultvectype;
421: };

423: PETSC_INTERN PetscErrorCode MatAXPY_Basic(Mat,PetscScalar,Mat,MatStructure);
424: PETSC_INTERN PetscErrorCode MatAXPY_BasicWithPreallocation(Mat,Mat,PetscScalar,Mat,MatStructure);
425: PETSC_INTERN PetscErrorCode MatAXPY_Basic_Preallocate(Mat,Mat,Mat*);

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

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

440: /*
441:     Object for partitioning graphs
442: */

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

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

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

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

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

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

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

505:     ncolors = 4;

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

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

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

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

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

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

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

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

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

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

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

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

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

651:   PetscErrorCode (*destroy)(Mat);
652: } Mat_SubSppt;

654: PETSC_EXTERN PetscErrorCode MatFactorDumpMatrix(Mat);
655: PETSC_INTERN PetscErrorCode MatShift_Basic(Mat,PetscScalar);
656: PETSC_INTERN PetscErrorCode MatSetBlockSizes_Default(Mat,PetscInt,PetscInt);

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

1390: PETSC_STATIC_INLINE PetscErrorCode PetscLLCondensedView(PetscInt *lnk)
1391: {
1393:   PetscInt       k;

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

1403: /*
1404:   Free memories used by the list
1405: */
1406: PETSC_STATIC_INLINE PetscErrorCode PetscLLCondensedDestroy(PetscInt *lnk,PetscBT bt)
1407: {

1411:   PetscFree(lnk);
1412:   PetscBTDestroy(&bt);
1413:   return(0);
1414: }

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

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

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

1445:   PetscIntMultError(2,nlnk_max+2,&lsize);
1446:   PetscRealloc(lsize*sizeof(PetscInt),lnk);
1447:   return(0);
1448: }

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

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

1491: PETSC_STATIC_INLINE PetscErrorCode PetscLLCondensedDestroy_Scalable(PetscInt *lnk)
1492: {
1493:   return PetscFree(lnk);
1494: }

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

1504:       The next three are always the first link

1506:       lnk[3]    PETSC_MIN_INT+1
1507:       lnk[4]    1
1508:       lnk[5]    link to first real entry

1510:       The next three are always the last link

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

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

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

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

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

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

1619: PETSC_STATIC_INLINE PetscErrorCode PetscLLCondensedDestroy_fast(PetscInt *lnk)
1620: {
1621:   return PetscFree(lnk);
1622: }

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

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

1728: #endif