Actual source code: matimpl.h

petsc-3.11.4 2019-09-28
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 (*placeholder_98)(Mat);
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 MatCopy_Basic(Mat,Mat,MatStructure);
235: PETSC_INTERN PetscErrorCode MatDiagonalSet_Default(Mat,Vec,InsertMode);

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

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

251: typedef struct _MatStashSpace *PetscMatStashSpace;

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

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

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

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

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

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

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

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

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

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

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

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

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

421: PETSC_INTERN PetscErrorCode MatAXPY_Basic(Mat,PetscScalar,Mat,MatStructure);
422: PETSC_INTERN PetscErrorCode MatAXPY_BasicWithPreallocation(Mat,Mat,PetscScalar,Mat,MatStructure);
423: PETSC_INTERN PetscErrorCode MatAXPY_Basic_Preallocate(Mat,Mat,Mat*);

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

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

438: /*
439:     Object for partitioning graphs
440: */

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

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

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

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

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

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

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

503:     ncolors = 4;

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

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

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

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

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

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

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

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

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

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

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

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

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

649:   PetscErrorCode (*destroy)(Mat);
650: } Mat_SubSppt;

652: PETSC_EXTERN PetscErrorCode MatFactorDumpMatrix(Mat);
653: PETSC_INTERN PetscErrorCode MatShift_Basic(Mat,PetscScalar);
654: PETSC_INTERN PetscErrorCode MatSetBlockSizes_Default(Mat,PetscInt,PetscInt);

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

1263: #define MatCheckSameLocalSize(A,ar1,B,ar2) \
1265:   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);
1266: 
1267: #define MatCheckSameSize(A,ar1,B,ar2) \
1268:   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);\
1269:   MatCheckSameLocalSize(A,ar1,B,ar2);
1270: 
1271: #define VecCheckMatCompatible(M,x,ar1,b,ar2)                               \
1272:   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);\
1273:   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);

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

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

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

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

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

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

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

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

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

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

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

1388: PETSC_STATIC_INLINE PetscErrorCode PetscLLCondensedView(PetscInt *lnk)
1389: {
1391:   PetscInt       k;

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

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

1409:   PetscFree(lnk);
1410:   PetscBTDestroy(&bt);
1411:   return(0);
1412: }

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

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

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

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

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

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

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

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

1502:       The next three are always the first link

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

1508:       The next three are always the last link

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

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

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

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

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

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

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

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

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

1727: #endif