Actual source code: matimpl.h

petsc-master 2020-01-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 src/mat/f90-mod/petscmat.h
 30: */
 31: typedef struct _MatOps *MatOps;
 32: struct _MatOps {
 33:   /* 0*/
 34:   PetscErrorCode (*setvalues)(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
 35:   PetscErrorCode (*getrow)(Mat,PetscInt,PetscInt *,PetscInt*[],PetscScalar*[]);
 36:   PetscErrorCode (*restorerow)(Mat,PetscInt,PetscInt *,PetscInt *[],PetscScalar *[]);
 37:   PetscErrorCode (*mult)(Mat,Vec,Vec);
 38:   PetscErrorCode (*multadd)(Mat,Vec,Vec,Vec);
 39:   /* 5*/
 40:   PetscErrorCode (*multtranspose)(Mat,Vec,Vec);
 41:   PetscErrorCode (*multtransposeadd)(Mat,Vec,Vec,Vec);
 42:   PetscErrorCode (*solve)(Mat,Vec,Vec);
 43:   PetscErrorCode (*solveadd)(Mat,Vec,Vec,Vec);
 44:   PetscErrorCode (*solvetranspose)(Mat,Vec,Vec);
 45:   /*10*/
 46:   PetscErrorCode (*solvetransposeadd)(Mat,Vec,Vec,Vec);
 47:   PetscErrorCode (*lufactor)(Mat,IS,IS,const MatFactorInfo*);
 48:   PetscErrorCode (*choleskyfactor)(Mat,IS,const MatFactorInfo*);
 49:   PetscErrorCode (*sor)(Mat,Vec,PetscReal,MatSORType,PetscReal,PetscInt,PetscInt,Vec);
 50:   PetscErrorCode (*transpose)(Mat,MatReuse,Mat *);
 51:   /*15*/
 52:   PetscErrorCode (*getinfo)(Mat,MatInfoType,MatInfo*);
 53:   PetscErrorCode (*equal)(Mat,Mat,PetscBool  *);
 54:   PetscErrorCode (*getdiagonal)(Mat,Vec);
 55:   PetscErrorCode (*diagonalscale)(Mat,Vec,Vec);
 56:   PetscErrorCode (*norm)(Mat,NormType,PetscReal*);
 57:   /*20*/
 58:   PetscErrorCode (*assemblybegin)(Mat,MatAssemblyType);
 59:   PetscErrorCode (*assemblyend)(Mat,MatAssemblyType);
 60:   PetscErrorCode (*setoption)(Mat,MatOption,PetscBool );
 61:   PetscErrorCode (*zeroentries)(Mat);
 62:   /*24*/
 63:   PetscErrorCode (*zerorows)(Mat,PetscInt,const PetscInt[],PetscScalar,Vec,Vec);
 64:   PetscErrorCode (*lufactorsymbolic)(Mat,Mat,IS,IS,const MatFactorInfo*);
 65:   PetscErrorCode (*lufactornumeric)(Mat,Mat,const MatFactorInfo*);
 66:   PetscErrorCode (*choleskyfactorsymbolic)(Mat,Mat,IS,const MatFactorInfo*);
 67:   PetscErrorCode (*choleskyfactornumeric)(Mat,Mat,const MatFactorInfo*);
 68:   /*29*/
 69:   PetscErrorCode (*setup)(Mat);
 70:   PetscErrorCode (*ilufactorsymbolic)(Mat,Mat,IS,IS,const MatFactorInfo*);
 71:   PetscErrorCode (*iccfactorsymbolic)(Mat,Mat,IS,const MatFactorInfo*);
 72:   PetscErrorCode (*getdiagonalblock)(Mat,Mat*);
 73:   PetscErrorCode (*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:   PetscErrorCode (*getvalueslocal)(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],PetscScalar[]);
211: };
212: /*
213:     If you add MatOps entries above also add them to the MATOP enum
214:     in include/petscmat.h and src/mat/f90-mod/petscmat.h
215: */

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

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

227: PETSC_EXTERN MatRootName MatRootNameList;

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

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

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

253: typedef struct _MatStashSpace *PetscMatStashSpace;

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

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

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

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

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

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

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

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

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

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

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

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

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

426: PETSC_INTERN PetscErrorCode MatPtAP_Basic(Mat,Mat,MatReuse,PetscReal,Mat*);
427: PETSC_INTERN PetscErrorCode MatAXPY_Basic(Mat,PetscScalar,Mat,MatStructure);
428: PETSC_INTERN PetscErrorCode MatAXPY_BasicWithPreallocation(Mat,Mat,PetscScalar,Mat,MatStructure);
429: PETSC_INTERN PetscErrorCode MatAXPY_Basic_Preallocate(Mat,Mat,Mat*);

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

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

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

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

458: struct _p_MatPartitioning {
459:   PETSCHEADER(struct _MatPartitioningOps);
460:   Mat         adj;
461:   PetscInt    *vertex_weights;
462:   PetscReal   *part_weights;
463:   PetscInt    n;                                 /* number of partitions */
464:   void        *data;
465:   PetscInt    setupcalled;
466:   PetscBool   use_edge_weights;  /* A flag indicates whether or not to use edge weights */
467: };

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

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

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

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

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

510:     ncolors = 4;

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

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

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

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

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

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

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

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

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

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

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

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

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

658:   PetscErrorCode (*destroy)(Mat);
659: } Mat_SubSppt;

661: PETSC_EXTERN PetscErrorCode MatFactorDumpMatrix(Mat);
662: PETSC_INTERN PetscErrorCode MatShift_Basic(Mat,PetscScalar);
663: PETSC_INTERN PetscErrorCode MatSetBlockSizes_Default(Mat,PetscInt,PetscInt);

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

1397: PETSC_STATIC_INLINE PetscErrorCode PetscLLCondensedView(PetscInt *lnk)
1398: {
1400:   PetscInt       k;

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

1410: /*
1411:   Free memories used by the list
1412: */
1413: PETSC_STATIC_INLINE PetscErrorCode PetscLLCondensedDestroy(PetscInt *lnk,PetscBT bt)
1414: {

1418:   PetscFree(lnk);
1419:   PetscBTDestroy(&bt);
1420:   return(0);
1421: }

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

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

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

1452:   PetscIntMultError(2,nlnk_max+2,&lsize);
1453:   PetscRealloc(lsize*sizeof(PetscInt),lnk);
1454:   return(0);
1455: }

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

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

1498: PETSC_STATIC_INLINE PetscErrorCode PetscLLCondensedDestroy_Scalable(PetscInt *lnk)
1499: {
1500:   return PetscFree(lnk);
1501: }

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

1511:       The next three are always the first link

1513:       lnk[3]    PETSC_MIN_INT+1
1514:       lnk[4]    1
1515:       lnk[5]    link to first real entry

1517:       The next three are always the last link

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

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

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

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

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

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

1626: PETSC_STATIC_INLINE PetscErrorCode PetscLLCondensedDestroy_fast(PetscInt *lnk)
1627: {
1628:   return PetscFree(lnk);
1629: }

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

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

1739: #endif