Actual source code: matimpl.h
 
   petsc-3.12.4 2020-02-04
   
  2: #ifndef __MATIMPL_H
  5:  #include <petscmat.h>
  6:  #include <petscmatcoarsen.h>
  7:  #include <petsc/private/petscimpl.h>
  9: PETSC_EXTERN PetscBool MatRegisterAllCalled;
 10: PETSC_EXTERN PetscBool MatSeqAIJRegisterAllCalled;
 11: PETSC_EXTERN PetscBool MatOrderingRegisterAllCalled;
 12: PETSC_EXTERN PetscBool MatColoringRegisterAllCalled;
 13: PETSC_EXTERN PetscBool MatPartitioningRegisterAllCalled;
 14: PETSC_EXTERN PetscBool MatCoarsenRegisterAllCalled;
 15: PETSC_EXTERN PetscErrorCode MatRegisterAll(void);
 16: PETSC_EXTERN PetscErrorCode MatOrderingRegisterAll(void);
 17: PETSC_EXTERN PetscErrorCode MatColoringRegisterAll(void);
 18: PETSC_EXTERN PetscErrorCode MatPartitioningRegisterAll(void);
 19: PETSC_EXTERN PetscErrorCode MatCoarsenRegisterAll(void);
 20: PETSC_EXTERN PetscErrorCode MatSeqAIJRegisterAll(void);
 22: /*
 23:   This file defines the parts of the matrix data structure that are
 24:   shared by all matrix types.
 25: */
 27: /*
 28:     If you add entries here also add them to the MATOP enum
 29:     in include/petscmat.h and include/petsc/finclude/petscmat.h
 30: */
 31: typedef struct _MatOps *MatOps;
 32: struct _MatOps {
 33:   /* 0*/
 34:   PetscErrorCode (*setvalues)(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
 35:   PetscErrorCode (*getrow)(Mat,PetscInt,PetscInt *,PetscInt*[],PetscScalar*[]);
 36:   PetscErrorCode (*restorerow)(Mat,PetscInt,PetscInt *,PetscInt *[],PetscScalar *[]);
 37:   PetscErrorCode (*mult)(Mat,Vec,Vec);
 38:   PetscErrorCode (*multadd)(Mat,Vec,Vec,Vec);
 39:   /* 5*/
 40:   PetscErrorCode (*multtranspose)(Mat,Vec,Vec);
 41:   PetscErrorCode (*multtransposeadd)(Mat,Vec,Vec,Vec);
 42:   PetscErrorCode (*solve)(Mat,Vec,Vec);
 43:   PetscErrorCode (*solveadd)(Mat,Vec,Vec,Vec);
 44:   PetscErrorCode (*solvetranspose)(Mat,Vec,Vec);
 45:   /*10*/
 46:   PetscErrorCode (*solvetransposeadd)(Mat,Vec,Vec,Vec);
 47:   PetscErrorCode (*lufactor)(Mat,IS,IS,const MatFactorInfo*);
 48:   PetscErrorCode (*choleskyfactor)(Mat,IS,const MatFactorInfo*);
 49:   PetscErrorCode (*sor)(Mat,Vec,PetscReal,MatSORType,PetscReal,PetscInt,PetscInt,Vec);
 50:   PetscErrorCode (*transpose)(Mat,MatReuse,Mat *);
 51:   /*15*/
 52:   PetscErrorCode (*getinfo)(Mat,MatInfoType,MatInfo*);
 53:   PetscErrorCode (*equal)(Mat,Mat,PetscBool  *);
 54:   PetscErrorCode (*getdiagonal)(Mat,Vec);
 55:   PetscErrorCode (*diagonalscale)(Mat,Vec,Vec);
 56:   PetscErrorCode (*norm)(Mat,NormType,PetscReal*);
 57:   /*20*/
 58:   PetscErrorCode (*assemblybegin)(Mat,MatAssemblyType);
 59:   PetscErrorCode (*assemblyend)(Mat,MatAssemblyType);
 60:   PetscErrorCode (*setoption)(Mat,MatOption,PetscBool );
 61:   PetscErrorCode (*zeroentries)(Mat);
 62:   /*24*/
 63:   PetscErrorCode (*zerorows)(Mat,PetscInt,const PetscInt[],PetscScalar,Vec,Vec);
 64:   PetscErrorCode (*lufactorsymbolic)(Mat,Mat,IS,IS,const MatFactorInfo*);
 65:   PetscErrorCode (*lufactornumeric)(Mat,Mat,const MatFactorInfo*);
 66:   PetscErrorCode (*choleskyfactorsymbolic)(Mat,Mat,IS,const MatFactorInfo*);
 67:   PetscErrorCode (*choleskyfactornumeric)(Mat,Mat,const MatFactorInfo*);
 68:   /*29*/
 69:   PetscErrorCode (*setup)(Mat);
 70:   PetscErrorCode (*ilufactorsymbolic)(Mat,Mat,IS,IS,const MatFactorInfo*);
 71:   PetscErrorCode (*iccfactorsymbolic)(Mat,Mat,IS,const MatFactorInfo*);
 72:   PetscErrorCode (*getdiagonalblock)(Mat,Mat*);
 73:   PetscErrorCode (*freeintermediatedatastructures)(Mat);
 74:   /*34*/
 75:   PetscErrorCode (*duplicate)(Mat,MatDuplicateOption,Mat*);
 76:   PetscErrorCode (*forwardsolve)(Mat,Vec,Vec);
 77:   PetscErrorCode (*backwardsolve)(Mat,Vec,Vec);
 78:   PetscErrorCode (*ilufactor)(Mat,IS,IS,const MatFactorInfo*);
 79:   PetscErrorCode (*iccfactor)(Mat,IS,const MatFactorInfo*);
 80:   /*39*/
 81:   PetscErrorCode (*axpy)(Mat,PetscScalar,Mat,MatStructure);
 82:   PetscErrorCode (*createsubmatrices)(Mat,PetscInt,const IS[],const IS[],MatReuse,Mat *[]);
 83:   PetscErrorCode (*increaseoverlap)(Mat,PetscInt,IS[],PetscInt);
 84:   PetscErrorCode (*getvalues)(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],PetscScalar []);
 85:   PetscErrorCode (*copy)(Mat,Mat,MatStructure);
 86:   /*44*/
 87:   PetscErrorCode (*getrowmax)(Mat,Vec,PetscInt[]);
 88:   PetscErrorCode (*scale)(Mat,PetscScalar);
 89:   PetscErrorCode (*shift)(Mat,PetscScalar);
 90:   PetscErrorCode (*diagonalset)(Mat,Vec,InsertMode);
 91:   PetscErrorCode (*zerorowscolumns)(Mat,PetscInt,const PetscInt[],PetscScalar,Vec,Vec);
 92:   /*49*/
 93:   PetscErrorCode (*setrandom)(Mat,PetscRandom);
 94:   PetscErrorCode (*getrowij)(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt*,const PetscInt *[],const PetscInt *[],PetscBool  *);
 95:   PetscErrorCode (*restorerowij)(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt *,const PetscInt *[],const PetscInt *[],PetscBool  *);
 96:   PetscErrorCode (*getcolumnij)(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt*,const PetscInt *[],const PetscInt *[],PetscBool  *);
 97:   PetscErrorCode (*restorecolumnij)(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt*,const PetscInt *[],const PetscInt *[],PetscBool  *);
 98:   /*54*/
 99:   PetscErrorCode (*fdcoloringcreate)(Mat,ISColoring,MatFDColoring);
100:   PetscErrorCode (*coloringpatch)(Mat,PetscInt,PetscInt,ISColoringValue[],ISColoring*);
101:   PetscErrorCode (*setunfactored)(Mat);
102:   PetscErrorCode (*permute)(Mat,IS,IS,Mat*);
103:   PetscErrorCode (*setvaluesblocked)(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
104:   /*59*/
105:   PetscErrorCode (*createsubmatrix)(Mat,IS,IS,MatReuse,Mat*);
106:   PetscErrorCode (*destroy)(Mat);
107:   PetscErrorCode (*view)(Mat,PetscViewer);
108:   PetscErrorCode (*convertfrom)(Mat,MatType,MatReuse,Mat*);
109:   PetscErrorCode (*matmatmult)(Mat,Mat,Mat,MatReuse,PetscReal,Mat*);
110:   /*64*/
111:   PetscErrorCode (*matmatmultsymbolic)(Mat,Mat,Mat,PetscReal,Mat*);
112:   PetscErrorCode (*matmatmultnumeric)(Mat,Mat,Mat,Mat);
113:   PetscErrorCode (*setlocaltoglobalmapping)(Mat,ISLocalToGlobalMapping,ISLocalToGlobalMapping);
114:   PetscErrorCode (*setvalueslocal)(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
115:   PetscErrorCode (*zerorowslocal)(Mat,PetscInt,const PetscInt[],PetscScalar,Vec,Vec);
116:   /*69*/
117:   PetscErrorCode (*getrowmaxabs)(Mat,Vec,PetscInt[]);
118:   PetscErrorCode (*getrowminabs)(Mat,Vec,PetscInt[]);
119:   PetscErrorCode (*convert)(Mat, MatType,MatReuse,Mat*);
120:   PetscErrorCode (*hasoperation)(Mat,MatOperation,PetscBool*);
121:   PetscErrorCode (*placeholder_73)(Mat,void*);
122:   /*74*/
123:   PetscErrorCode (*setvaluesadifor)(Mat,PetscInt,void*);
124:   PetscErrorCode (*fdcoloringapply)(Mat,MatFDColoring,Vec,void*);
125:   PetscErrorCode (*setfromoptions)(PetscOptionItems*,Mat);
126:   PetscErrorCode (*multconstrained)(Mat,Vec,Vec);
127:   PetscErrorCode (*multtransposeconstrained)(Mat,Vec,Vec);
128:   /*79*/
129:   PetscErrorCode (*findzerodiagonals)(Mat,IS*);
130:   PetscErrorCode (*mults)(Mat, Vecs, Vecs);
131:   PetscErrorCode (*solves)(Mat, Vecs, Vecs);
132:   PetscErrorCode (*getinertia)(Mat,PetscInt*,PetscInt*,PetscInt*);
133:   PetscErrorCode (*load)(Mat, PetscViewer);
134:   /*84*/
135:   PetscErrorCode (*issymmetric)(Mat,PetscReal,PetscBool *);
136:   PetscErrorCode (*ishermitian)(Mat,PetscReal,PetscBool *);
137:   PetscErrorCode (*isstructurallysymmetric)(Mat,PetscBool *);
138:   PetscErrorCode (*setvaluesblockedlocal)(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
139:   PetscErrorCode (*getvecs)(Mat,Vec*,Vec*);
140:   /*89*/
141:   PetscErrorCode (*matmult)(Mat,Mat,MatReuse,PetscReal,Mat*);
142:   PetscErrorCode (*matmultsymbolic)(Mat,Mat,PetscReal,Mat*);
143:   PetscErrorCode (*matmultnumeric)(Mat,Mat,Mat);
144:   PetscErrorCode (*ptap)(Mat,Mat,MatReuse,PetscReal,Mat*);
145:   PetscErrorCode (*ptapsymbolic)(Mat,Mat,PetscReal,Mat*); /* double dispatch wrapper routine */
146:   /*94*/
147:   PetscErrorCode (*ptapnumeric)(Mat,Mat,Mat);             /* double dispatch wrapper routine */
148:   PetscErrorCode (*mattransposemult)(Mat,Mat,MatReuse,PetscReal,Mat*);
149:   PetscErrorCode (*mattransposemultsymbolic)(Mat,Mat,PetscReal,Mat*);
150:   PetscErrorCode (*mattransposemultnumeric)(Mat,Mat,Mat);
151:   PetscErrorCode (*pintocpu)(Mat,PetscBool);
152:   /*99*/
153:   PetscErrorCode (*placeholder_99)(Mat);
154:   PetscErrorCode (*placeholder_100)(Mat);
155:   PetscErrorCode (*placeholder_101)(Mat);
156:   PetscErrorCode (*conjugate)(Mat);                              /* complex conjugate */
157:   PetscErrorCode (*viewnative)(Mat,PetscViewer);
158:   /*104*/
159:   PetscErrorCode (*setvaluesrow)(Mat,PetscInt,const PetscScalar[]);
160:   PetscErrorCode (*realpart)(Mat);
161:   PetscErrorCode (*imaginarypart)(Mat);
162:   PetscErrorCode (*getrowuppertriangular)(Mat);
163:   PetscErrorCode (*restorerowuppertriangular)(Mat);
164:   /*109*/
165:   PetscErrorCode (*matsolve)(Mat,Mat,Mat);
166:   PetscErrorCode (*matsolvetranspose)(Mat,Mat,Mat);
167:   PetscErrorCode (*getrowmin)(Mat,Vec,PetscInt[]);
168:   PetscErrorCode (*getcolumnvector)(Mat,Vec,PetscInt);
169:   PetscErrorCode (*missingdiagonal)(Mat,PetscBool *,PetscInt*);
170:   /*114*/
171:   PetscErrorCode (*getseqnonzerostructure)(Mat,Mat *);
172:   PetscErrorCode (*create)(Mat);
173:   PetscErrorCode (*getghosts)(Mat,PetscInt*,const PetscInt *[]);
174:   PetscErrorCode (*getlocalsubmatrix)(Mat,IS,IS,Mat*);
175:   PetscErrorCode (*restorelocalsubmatrix)(Mat,IS,IS,Mat*);
176:   /*119*/
177:   PetscErrorCode (*multdiagonalblock)(Mat,Vec,Vec);
178:   PetscErrorCode (*hermitiantranspose)(Mat,MatReuse,Mat*);
179:   PetscErrorCode (*multhermitiantranspose)(Mat,Vec,Vec);
180:   PetscErrorCode (*multhermitiantransposeadd)(Mat,Vec,Vec,Vec);
181:   PetscErrorCode (*getmultiprocblock)(Mat,MPI_Comm,MatReuse,Mat*);
182:   /*124*/
183:   PetscErrorCode (*findnonzerorows)(Mat,IS*);
184:   PetscErrorCode (*getcolumnnorms)(Mat,NormType,PetscReal*);
185:   PetscErrorCode (*invertblockdiagonal)(Mat,const PetscScalar**);
186:   PetscErrorCode (*invertvariableblockdiagonal)(Mat,PetscInt,const PetscInt*,PetscScalar*);
187:   PetscErrorCode (*createsubmatricesmpi)(Mat,PetscInt,const IS[], const IS[], MatReuse, Mat**);
188:   /*129*/
189:   PetscErrorCode (*setvaluesbatch)(Mat,PetscInt,PetscInt,PetscInt*,const PetscScalar*);
190:   PetscErrorCode (*transposematmult)(Mat,Mat,MatReuse,PetscReal,Mat*);
191:   PetscErrorCode (*transposematmultsymbolic)(Mat,Mat,PetscReal,Mat*);
192:   PetscErrorCode (*transposematmultnumeric)(Mat,Mat,Mat);
193:   PetscErrorCode (*transposecoloringcreate)(Mat,ISColoring,MatTransposeColoring);
194:   /*134*/
195:   PetscErrorCode (*transcoloringapplysptoden)(MatTransposeColoring,Mat,Mat);
196:   PetscErrorCode (*transcoloringapplydentosp)(MatTransposeColoring,Mat,Mat);
197:   PetscErrorCode (*rart)(Mat,Mat,MatReuse,PetscReal,Mat*);
198:   PetscErrorCode (*rartsymbolic)(Mat,Mat,PetscReal,Mat*); /* double dispatch wrapper routine */
199:   PetscErrorCode (*rartnumeric)(Mat,Mat,Mat);             /* double dispatch wrapper routine */
200:   /*139*/
201:   PetscErrorCode (*setblocksizes)(Mat,PetscInt,PetscInt);
202:   PetscErrorCode (*aypx)(Mat,PetscScalar,Mat,MatStructure);
203:   PetscErrorCode (*residual)(Mat,Vec,Vec,Vec);
204:   PetscErrorCode (*fdcoloringsetup)(Mat,ISColoring,MatFDColoring);
205:   PetscErrorCode (*findoffblockdiagonalentries)(Mat,IS*);
206:   /*144*/
207:   PetscErrorCode (*creatempimatconcatenateseqmat)(MPI_Comm,Mat,PetscInt,MatReuse,Mat*);
208:   PetscErrorCode (*destroysubmatrices)(PetscInt,Mat*[]);
209:   PetscErrorCode (*mattransposesolve)(Mat,Mat,Mat);
210: };
211: /*
212:     If you add MatOps entries above also add them to the MATOP enum
213:     in include/petscmat.h and include/petsc/finclude/petscmat.h
214: */
216:  #include <petscsys.h>
217: PETSC_EXTERN PetscErrorCode MatRegisterOp(MPI_Comm, const char[], PetscVoidFunction, const char[], PetscInt, ...);
218: PETSC_EXTERN PetscErrorCode MatQueryOp(MPI_Comm, PetscVoidFunction*, const char[], PetscInt, ...);
220: typedef struct _p_MatRootName* MatRootName;
221: struct _p_MatRootName {
222:   char        *rname,*sname,*mname;
223:   MatRootName next;
224: };
226: PETSC_EXTERN MatRootName MatRootNameList;
228: /*
229:    Utility private matrix routines
230: */
231: PETSC_INTERN PetscErrorCode MatFindNonzeroRowsOrCols_Basic(Mat,PetscBool,PetscReal,IS*);
232: PETSC_INTERN PetscErrorCode MatConvert_Basic(Mat,MatType,MatReuse,Mat*);
233: PETSC_INTERN PetscErrorCode MatConvert_Shell(Mat,MatType,MatReuse,Mat*);
234: PETSC_INTERN PetscErrorCode MatConvertFrom_Shell(Mat,MatType,MatReuse,Mat*);
235: PETSC_INTERN PetscErrorCode MatCopy_Basic(Mat,Mat,MatStructure);
236: PETSC_INTERN PetscErrorCode MatDiagonalSet_Default(Mat,Vec,InsertMode);
238: #if defined(PETSC_USE_DEBUG)
239: #  define MatCheckPreallocated(A,arg) do {                              \
240:     if (PetscUnlikely(!(A)->preallocated)) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call MatXXXSetPreallocation() or MatSetUp() on argument %D \"%s\" before %s()",(arg),#A,PETSC_FUNCTION_NAME); \
241:   } while (0)
242: #else
243: #  define MatCheckPreallocated(A,arg) do {} while (0)
244: #endif
246: /*
247:   The stash is used to temporarily store inserted matrix values that
248:   belong to another processor. During the assembly phase the stashed
249:   values are moved to the correct processor and
250: */
252: typedef struct _MatStashSpace *PetscMatStashSpace;
254: struct _MatStashSpace {
255:   PetscMatStashSpace next;
256:   PetscScalar        *space_head,*val;
257:   PetscInt           *idx,*idy;
258:   PetscInt           total_space_size;
259:   PetscInt           local_used;
260:   PetscInt           local_remaining;
261: };
263: PETSC_EXTERN PetscErrorCode PetscMatStashSpaceGet(PetscInt,PetscInt,PetscMatStashSpace *);
264: PETSC_EXTERN PetscErrorCode PetscMatStashSpaceContiguous(PetscInt,PetscMatStashSpace *,PetscScalar *,PetscInt *,PetscInt *);
265: PETSC_EXTERN PetscErrorCode PetscMatStashSpaceDestroy(PetscMatStashSpace*);
267: typedef struct {
268:   PetscInt    count;
269: } MatStashHeader;
271: typedef struct {
272:   void        *buffer;          /* Of type blocktype, dynamically constructed  */
273:   PetscInt    count;
274:   char        pending;
275: } MatStashFrame;
277: typedef struct _MatStash MatStash;
278: struct _MatStash {
279:   PetscInt      nmax;                   /* maximum stash size */
280:   PetscInt      umax;                   /* user specified max-size */
281:   PetscInt      oldnmax;                /* the nmax value used previously */
282:   PetscInt      n;                      /* stash size */
283:   PetscInt      bs;                     /* block size of the stash */
284:   PetscInt      reallocs;               /* preserve the no of mallocs invoked */
285:   PetscMatStashSpace space_head,space;  /* linked list to hold stashed global row/column numbers and matrix values */
287:   PetscErrorCode (*ScatterBegin)(Mat,MatStash*,PetscInt*);
288:   PetscErrorCode (*ScatterGetMesg)(MatStash*,PetscMPIInt*,PetscInt**,PetscInt**,PetscScalar**,PetscInt*);
289:   PetscErrorCode (*ScatterEnd)(MatStash*);
290:   PetscErrorCode (*ScatterDestroy)(MatStash*);
292:   /* The following variables are used for communication */
293:   MPI_Comm      comm;
294:   PetscMPIInt   size,rank;
295:   PetscMPIInt   tag1,tag2;
296:   MPI_Request   *send_waits;            /* array of send requests */
297:   MPI_Request   *recv_waits;            /* array of receive requests */
298:   MPI_Status    *send_status;           /* array of send status */
299:   PetscInt      nsends,nrecvs;          /* numbers of sends and receives */
300:   PetscScalar   *svalues;               /* sending data */
301:   PetscInt      *sindices;
302:   PetscScalar   **rvalues;              /* receiving data (values) */
303:   PetscInt      **rindices;             /* receiving data (indices) */
304:   PetscInt      nprocessed;             /* number of messages already processed */
305:   PetscMPIInt   *flg_v;                 /* indicates what messages have arrived so far and from whom */
306:   PetscBool     reproduce;
307:   PetscInt      reproduce_count;
309:   /* The following variables are used for BTS communication */
310:   PetscBool      first_assembly_done;   /* Is the first time matrix assembly done? */
311:   PetscBool      use_status;            /* Use MPI_Status to determine number of items in each message */
312:   PetscMPIInt    nsendranks;
313:   PetscMPIInt    nrecvranks;
314:   PetscMPIInt    *sendranks;
315:   PetscMPIInt    *recvranks;
316:   MatStashHeader *sendhdr,*recvhdr;
317:   MatStashFrame  *sendframes;   /* pointers to the main messages */
318:   MatStashFrame  *recvframes;
319:   MatStashFrame  *recvframe_active;
320:   PetscInt       recvframe_i;     /* index of block within active frame */
321:   PetscMPIInt    recvframe_count; /* Count actually sent for current frame */
322:   PetscInt       recvcount;       /* Number of receives processed so far */
323:   PetscMPIInt    *some_indices;   /* From last call to MPI_Waitsome */
324:   MPI_Status     *some_statuses;  /* Statuses from last call to MPI_Waitsome */
325:   PetscMPIInt    some_count;      /* Number of requests completed in last call to MPI_Waitsome */
326:   PetscMPIInt    some_i;          /* Index of request currently being processed */
327:   MPI_Request    *sendreqs;
328:   MPI_Request    *recvreqs;
329:   PetscSegBuffer segsendblocks;
330:   PetscSegBuffer segrecvframe;
331:   PetscSegBuffer segrecvblocks;
332:   MPI_Datatype   blocktype;
333:   size_t         blocktype_size;
334:   InsertMode     *insertmode;   /* Pointer to check mat->insertmode and set upon message arrival in case no local values have been set. */
335: };
337: #if !defined(PETSC_HAVE_MPIUNI)
338: PETSC_INTERN PetscErrorCode MatStashScatterDestroy_BTS(MatStash*);
339: #endif
340: PETSC_INTERN PetscErrorCode MatStashCreate_Private(MPI_Comm,PetscInt,MatStash*);
341: PETSC_INTERN PetscErrorCode MatStashDestroy_Private(MatStash*);
342: PETSC_INTERN PetscErrorCode MatStashScatterEnd_Private(MatStash*);
343: PETSC_INTERN PetscErrorCode MatStashSetInitialSize_Private(MatStash*,PetscInt);
344: PETSC_INTERN PetscErrorCode MatStashGetInfo_Private(MatStash*,PetscInt*,PetscInt*);
345: PETSC_INTERN PetscErrorCode MatStashValuesRow_Private(MatStash*,PetscInt,PetscInt,const PetscInt[],const PetscScalar[],PetscBool );
346: PETSC_INTERN PetscErrorCode MatStashValuesCol_Private(MatStash*,PetscInt,PetscInt,const PetscInt[],const PetscScalar[],PetscInt,PetscBool );
347: PETSC_INTERN PetscErrorCode MatStashValuesRowBlocked_Private(MatStash*,PetscInt,PetscInt,const PetscInt[],const PetscScalar[],PetscInt,PetscInt,PetscInt);
348: PETSC_INTERN PetscErrorCode MatStashValuesColBlocked_Private(MatStash*,PetscInt,PetscInt,const PetscInt[],const PetscScalar[],PetscInt,PetscInt,PetscInt);
349: PETSC_INTERN PetscErrorCode MatStashScatterBegin_Private(Mat,MatStash*,PetscInt*);
350: PETSC_INTERN PetscErrorCode MatStashScatterGetMesg_Private(MatStash*,PetscMPIInt*,PetscInt**,PetscInt**,PetscScalar**,PetscInt*);
351: PETSC_INTERN PetscErrorCode MatGetInfo_External(Mat,MatInfoType,MatInfo*);
353: typedef struct {
354:   PetscInt   dim;
355:   PetscInt   dims[4];
356:   PetscInt   starts[4];
357:   PetscBool  noc;        /* this is a single component problem, hence user will not set MatStencil.c */
358: } MatStencilInfo;
360: /* Info about using compressed row format */
361: typedef struct {
362:   PetscBool  use;                           /* indicates compressed rows have been checked and will be used */
363:   PetscInt   nrows;                         /* number of non-zero rows */
364:   PetscInt   *i;                            /* compressed row pointer  */
365:   PetscInt   *rindex;                       /* compressed row index               */
366: } Mat_CompressedRow;
367: PETSC_EXTERN PetscErrorCode MatCheckCompressedRow(Mat,PetscInt,Mat_CompressedRow*,PetscInt*,PetscInt,PetscReal);
369: typedef struct { /* used by MatCreateRedundantMatrix() for reusing matredundant */
370:   PetscInt     nzlocal,nsends,nrecvs;
371:   PetscMPIInt  *send_rank,*recv_rank;
372:   PetscInt     *sbuf_nz,*rbuf_nz,*sbuf_j,**rbuf_j;
373:   PetscScalar  *sbuf_a,**rbuf_a;
374:   MPI_Comm     subcomm;   /* when user does not provide a subcomm */
375:   IS           isrow,iscol;
376:   Mat          *matseq;
377: } Mat_Redundant;
379: struct _p_Mat {
380:   PETSCHEADER(struct _MatOps);
381:   PetscLayout            rmap,cmap;
382:   void                   *data;            /* implementation-specific data */
383:   MatFactorType          factortype;       /* MAT_FACTOR_LU, ILU, CHOLESKY or ICC */
384:   PetscBool              assembled;        /* is the matrix assembled? */
385:   PetscBool              was_assembled;    /* new values inserted into assembled mat */
386:   PetscInt               num_ass;          /* number of times matrix has been assembled */
387:   PetscObjectState       nonzerostate;     /* each time new nonzeros locations are introduced into the matrix this is updated */
388:   PetscObjectState       ass_nonzerostate; /* nonzero state at last assembly */
389:   MatInfo                info;             /* matrix information */
390:   InsertMode             insertmode;       /* have values been inserted in matrix or added? */
391:   MatStash               stash,bstash;     /* used for assembling off-proc mat emements */
392:   MatNullSpace           nullsp;           /* null space (operator is singular) */
393:   MatNullSpace           transnullsp;      /* null space of transpose of operator */
394:   MatNullSpace           nearnullsp;       /* near null space to be used by multigrid methods */
395:   PetscInt               congruentlayouts; /* are the rows and columns layouts congruent? */
396:   PetscBool              preallocated;
397:   MatStencilInfo         stencil;          /* information for structured grid */
398:   PetscBool              symmetric,hermitian,structurally_symmetric,spd;
399:   PetscBool              symmetric_set,hermitian_set,structurally_symmetric_set,spd_set; /* if true, then corresponding flag is correct*/
400:   PetscBool              symmetric_eternal;
401:   PetscBool              nooffprocentries,nooffproczerorows;
402:   PetscBool              assembly_subset;  /* set by MAT_SUBSET_OFF_PROC_ENTRIES */
403:   PetscBool              submat_singleis;  /* for efficient PCSetUp_ASM() */
404:   PetscBool              structure_only;
405:   PetscBool              sortedfull;       /* full, sorted rows are inserted */
406: #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_CUDA)
407:   PetscOffloadMask       offloadmask;      /* a mask which indicates where the valid matrix data is (GPU, CPU or both) */
408:   PetscBool              pinnedtocpu;
409: #endif
410:   void                   *spptr;          /* pointer for special library like SuperLU */
411:   char                   *solvertype;
412:   PetscBool              checksymmetryonassembly,checknullspaceonassembly;
413:   PetscReal              checksymmetrytol;
414:   Mat                    schur;             /* Schur complement matrix */
415:   MatFactorSchurStatus   schur_status;      /* status of the Schur complement matrix */
416:   Mat_Redundant          *redundant;        /* used by MatCreateRedundantMatrix() */
417:   PetscBool              erroriffailure;    /* Generate an error if detected (for example a zero pivot) instead of returning */
418:   MatFactorError         factorerrortype;               /* type of error in factorization */
419:   PetscReal              factorerror_zeropivot_value;   /* If numerical zero pivot was detected this is the computed value */
420:   PetscInt               factorerror_zeropivot_row;     /* Row where zero pivot was detected */
421:   PetscInt               nblocks,*bsizes;   /* support for MatSetVariableBlockSizes() */
422:   char                   *defaultvectype;
423: };
425: PETSC_INTERN PetscErrorCode MatPtAP_Basic(Mat,Mat,MatReuse,PetscReal,Mat*);
426: PETSC_INTERN PetscErrorCode MatAXPY_Basic(Mat,PetscScalar,Mat,MatStructure);
427: PETSC_INTERN PetscErrorCode MatAXPY_BasicWithPreallocation(Mat,Mat,PetscScalar,Mat,MatStructure);
428: PETSC_INTERN PetscErrorCode MatAXPY_Basic_Preallocate(Mat,Mat,Mat*);
430: /*
431:     Utility for MatFactor (Schur complement)
432: */
433: PETSC_INTERN PetscErrorCode MatFactorFactorizeSchurComplement_Private(Mat);
434: PETSC_INTERN PetscErrorCode MatFactorInvertSchurComplement_Private(Mat);
435: PETSC_INTERN PetscErrorCode MatFactorUpdateSchurStatus_Private(Mat);
436: PETSC_INTERN PetscErrorCode MatFactorSetUpInPlaceSchur_Private(Mat);
438: /*
439:     Utility for MatZeroRows
440: */
441: PETSC_INTERN PetscErrorCode MatZeroRowsMapLocal_Private(Mat,PetscInt,const PetscInt*,PetscInt*,PetscInt**);
443: /*
444:     Object for partitioning graphs
445: */
447: typedef struct _MatPartitioningOps *MatPartitioningOps;
448: struct _MatPartitioningOps {
449:   PetscErrorCode (*apply)(MatPartitioning,IS*);
450:   PetscErrorCode (*applynd)(MatPartitioning,IS*);
451:   PetscErrorCode (*setfromoptions)(PetscOptionItems*,MatPartitioning);
452:   PetscErrorCode (*destroy)(MatPartitioning);
453:   PetscErrorCode (*view)(MatPartitioning,PetscViewer);
454:   PetscErrorCode (*improve)(MatPartitioning,IS*);
455: };
457: struct _p_MatPartitioning {
458:   PETSCHEADER(struct _MatPartitioningOps);
459:   Mat         adj;
460:   PetscInt    *vertex_weights;
461:   PetscReal   *part_weights;
462:   PetscInt    n;                                 /* number of partitions */
463:   void        *data;
464:   PetscInt    setupcalled;
465:   PetscBool   use_edge_weights;  /* A flag indicates whether or not to use edge weights */
466: };
468: /* needed for parallel nested dissection by ParMetis and PTSCOTCH */
469: PETSC_INTERN PetscErrorCode MatPartitioningSizesToSep_Private(PetscInt,PetscInt[],PetscInt[],PetscInt[]);
471: /*
472:     Object for coarsen graphs
473: */
474: typedef struct _MatCoarsenOps *MatCoarsenOps;
475: struct _MatCoarsenOps {
476:   PetscErrorCode (*apply)(MatCoarsen);
477:   PetscErrorCode (*setfromoptions)(PetscOptionItems*,MatCoarsen);
478:   PetscErrorCode (*destroy)(MatCoarsen);
479:   PetscErrorCode (*view)(MatCoarsen,PetscViewer);
480: };
482: struct _p_MatCoarsen {
483:   PETSCHEADER(struct _MatCoarsenOps);
484:   Mat              graph;
485:   PetscInt         setupcalled;
486:   void             *subctx;
487:   /* */
488:   PetscBool        strict_aggs;
489:   IS               perm;
490:   PetscCoarsenData *agg_lists;
491: };
493: /*
494:     MatFDColoring is used to compute Jacobian matrices efficiently
495:   via coloring. The data structure is explained below in an example.
497:    Color =   0    1     0    2   |   2      3       0
498:    ---------------------------------------------------
499:             00   01              |          05
500:             10   11              |   14     15               Processor  0
501:                        22    23  |          25
502:                        32    33  |
503:    ===================================================
504:                                  |   44     45     46
505:             50                   |          55               Processor 1
506:                                  |   64            66
507:    ---------------------------------------------------
509:     ncolors = 4;
511:     ncolumns      = {2,1,1,0}
512:     columns       = {{0,2},{1},{3},{}}
513:     nrows         = {4,2,3,3}
514:     rows          = {{0,1,2,3},{0,1},{1,2,3},{0,1,2}}
515:     vwscale       = {dx(0),dx(1),dx(2),dx(3)}               MPI Vec
516:     vscale        = {dx(0),dx(1),dx(2),dx(3),dx(4),dx(5)}   Seq Vec
518:     ncolumns      = {1,0,1,1}
519:     columns       = {{6},{},{4},{5}}
520:     nrows         = {3,0,2,2}
521:     rows          = {{0,1,2},{},{1,2},{1,2}}
522:     vwscale       = {dx(4),dx(5),dx(6)}              MPI Vec
523:     vscale        = {dx(0),dx(4),dx(5),dx(6)}        Seq Vec
525:     See the routine MatFDColoringApply() for how this data is used
526:     to compute the Jacobian.
528: */
529: typedef struct {
530:   PetscInt     row;
531:   PetscInt     col;
532:   PetscScalar  *valaddr;   /* address of value */
533: } MatEntry;
535: typedef struct {
536:   PetscInt     row;
537:   PetscScalar  *valaddr;   /* address of value */
538: } MatEntry2;
540: struct  _p_MatFDColoring{
541:   PETSCHEADER(int);
542:   PetscInt       M,N,m;            /* total rows, columns; local rows */
543:   PetscInt       rstart;           /* first row owned by local processor */
544:   PetscInt       ncolors;          /* number of colors */
545:   PetscInt       *ncolumns;        /* number of local columns for a color */
546:   PetscInt       **columns;        /* lists the local columns of each color (using global column numbering) */
547:   IS             *isa;             /* these are the IS that contain the column values given in columns */
548:   PetscInt       *nrows;           /* number of local rows for each color */
549:   MatEntry       *matentry;        /* holds (row, column, address of value) for Jacobian matrix entry */
550:   MatEntry2      *matentry2;       /* holds (row, address of value) for Jacobian matrix entry */
551:   PetscScalar    *dy;              /* store a block of F(x+dx)-F(x) when J is in BAIJ format */
552:   PetscReal      error_rel;        /* square root of relative error in computing function */
553:   PetscReal      umin;             /* minimum allowable u'dx value */
554:   Vec            w1,w2,w3;         /* work vectors used in computing Jacobian */
555:   PetscBool      fset;             /* indicates that the initial function value F(X) is set */
556:   PetscErrorCode (*f)(void);       /* function that defines Jacobian */
557:   void           *fctx;            /* optional user-defined context for use by the function f */
558:   Vec            vscale;           /* holds FD scaling, i.e. 1/dx for each perturbed column */
559:   PetscInt       currentcolor;     /* color for which function evaluation is being done now */
560:   const char     *htype;           /* "wp" or "ds" */
561:   ISColoringType ctype;            /* IS_COLORING_GLOBAL or IS_COLORING_LOCAL */
562:   PetscInt       brows,bcols;      /* number of block rows or columns for speedup inserting the dense matrix into sparse Jacobian */
563:   PetscBool      setupcalled;      /* true if setup has been called */
564:   PetscBool      viewed;           /* true if the -mat_fd_coloring_view has been triggered already */
565:   void           (*ftn_func_pointer)(void),*ftn_func_cntx; /* serve the same purpose as *fortran_func_pointers in PETSc objects */
566:   PetscObjectId  matid;            /* matrix this object was created with, must always be the same */
567: };
569: typedef struct _MatColoringOps *MatColoringOps;
570: struct _MatColoringOps {
571:   PetscErrorCode (*destroy)(MatColoring);
572:   PetscErrorCode (*setfromoptions)(PetscOptionItems*,MatColoring);
573:   PetscErrorCode (*view)(MatColoring,PetscViewer);
574:   PetscErrorCode (*apply)(MatColoring,ISColoring*);
575:   PetscErrorCode (*weights)(MatColoring,PetscReal**,PetscInt**);
576: };
578: struct _p_MatColoring {
579:   PETSCHEADER(struct _MatColoringOps);
580:   Mat                   mat;
581:   PetscInt              dist;             /* distance of the coloring */
582:   PetscInt              maxcolors;        /* the maximum number of colors returned, maxcolors=1 for MIS */
583:   void                  *data;            /* inner context */
584:   PetscBool             valid;            /* check to see if what is produced is a valid coloring */
585:   MatColoringWeightType weight_type;      /* type of weight computation to be performed */
586:   PetscReal             *user_weights;    /* custom weights and permutation */
587:   PetscInt              *user_lperm;
588:   PetscBool             valid_iscoloring; /* check to see if matcoloring is produced a valid iscoloring */
589: };
591: struct  _p_MatTransposeColoring{
592:   PETSCHEADER(int);
593:   PetscInt       M,N,m;            /* total rows, columns; local rows */
594:   PetscInt       rstart;           /* first row owned by local processor */
595:   PetscInt       ncolors;          /* number of colors */
596:   PetscInt       *ncolumns;        /* number of local columns for a color */
597:   PetscInt       *nrows;           /* number of local rows for each color */
598:   PetscInt       currentcolor;     /* color for which function evaluation is being done now */
599:   ISColoringType ctype;            /* IS_COLORING_GLOBAL or IS_COLORING_LOCAL */
601:   PetscInt       *colorforrow,*colorforcol;  /* pointer to rows and columns */
602:   PetscInt       *rows;                      /* lists the local rows for each color (using the local row numbering) */
603:   PetscInt       *den2sp;                    /* maps (row,color) in the dense matrix to index of sparse matrix array a->a */
604:   PetscInt       *columns;                   /* lists the local columns of each color (using global column numbering) */
605:   PetscInt       brows;                      /* number of rows for efficient implementation of MatTransColoringApplyDenToSp() */
606:   PetscInt       *lstart;                    /* array used for loop over row blocks of Csparse */
607: };
609: /*
610:    Null space context for preconditioner/operators
611: */
612: struct _p_MatNullSpace {
613:   PETSCHEADER(int);
614:   PetscBool      has_cnst;
615:   PetscInt       n;
616:   Vec*           vecs;
617:   PetscScalar*   alpha;                 /* for projections */
618:   PetscErrorCode (*remove)(MatNullSpace,Vec,void*);  /* for user provided removal function */
619:   void*          rmctx;                 /* context for remove() function */
620: };
622: /*
623:    Checking zero pivot for LU, ILU preconditioners.
624: */
625: typedef struct {
626:   PetscInt       nshift,nshift_max;
627:   PetscReal      shift_amount,shift_lo,shift_hi,shift_top,shift_fraction;
628:   PetscBool      newshift;
629:   PetscReal      rs;  /* active row sum of abs(offdiagonals) */
630:   PetscScalar    pv;  /* pivot of the active row */
631: } FactorShiftCtx;
633: /*
634:  Used by MatCreateSubMatrices_MPIXAIJ_Local()
635: */
636:  #include <petscctable.h>
637: typedef struct { /* used by MatCreateSubMatrices_MPIAIJ_SingleIS_Local() and MatCreateSubMatrices_MPIAIJ_Local */
638:   PetscInt   id;   /* index of submats, only submats[0] is responsible for deleting some arrays below */
639:   PetscInt   nrqs,nrqr;
640:   PetscInt   **rbuf1,**rbuf2,**rbuf3,**sbuf1,**sbuf2;
641:   PetscInt   **ptr;
642:   PetscInt   *tmp;
643:   PetscInt   *ctr;
644:   PetscInt   *pa; /* proc array */
645:   PetscInt   *req_size,*req_source1,*req_source2;
646:   PetscBool  allcolumns,allrows;
647:   PetscBool  singleis;
648:   PetscInt   *row2proc; /* row to proc map */
649:   PetscInt   nstages;
650: #if defined(PETSC_USE_CTABLE)
651:   PetscTable cmap,rmap;
652:   PetscInt   *cmap_loc,*rmap_loc;
653: #else
654:   PetscInt   *cmap,*rmap;
655: #endif
657:   PetscErrorCode (*destroy)(Mat);
658: } Mat_SubSppt;
660: PETSC_EXTERN PetscErrorCode MatFactorDumpMatrix(Mat);
661: PETSC_INTERN PetscErrorCode MatShift_Basic(Mat,PetscScalar);
662: PETSC_INTERN PetscErrorCode MatSetBlockSizes_Default(Mat,PetscInt,PetscInt);
664: PETSC_STATIC_INLINE PetscErrorCode MatPivotCheck_nz(Mat mat,const MatFactorInfo *info,FactorShiftCtx *sctx,PetscInt row)
665: {
666:   PetscReal _rs   = sctx->rs;
667:   PetscReal _zero = info->zeropivot*_rs;
670:   if (PetscAbsScalar(sctx->pv) <= _zero && !PetscIsNanScalar(sctx->pv)){
671:     /* force |diag| > zeropivot*rs */
672:     if (!sctx->nshift) sctx->shift_amount = info->shiftamount;
673:     else sctx->shift_amount *= 2.0;
674:     sctx->newshift = PETSC_TRUE;
675:     (sctx->nshift)++;
676:   } else {
677:     sctx->newshift = PETSC_FALSE;
678:   }
679:   return(0);
680: }
682: PETSC_STATIC_INLINE PetscErrorCode MatPivotCheck_pd(Mat mat,const MatFactorInfo *info,FactorShiftCtx *sctx,PetscInt row)
683: {
684:   PetscReal _rs   = sctx->rs;
685:   PetscReal _zero = info->zeropivot*_rs;
688:   if (PetscRealPart(sctx->pv) <= _zero && !PetscIsNanScalar(sctx->pv)){
689:     /* force matfactor to be diagonally dominant */
690:     if (sctx->nshift == sctx->nshift_max) {
691:       sctx->shift_fraction = sctx->shift_hi;
692:     } else {
693:       sctx->shift_lo = sctx->shift_fraction;
694:       sctx->shift_fraction = (sctx->shift_hi+sctx->shift_lo)/2.;
695:     }
696:     sctx->shift_amount = sctx->shift_fraction * sctx->shift_top;
697:     sctx->nshift++;
698:     sctx->newshift = PETSC_TRUE;
699:   } else {
700:     sctx->newshift = PETSC_FALSE;
701:   }
702:   return(0);
703: }
705: PETSC_STATIC_INLINE PetscErrorCode MatPivotCheck_inblocks(Mat mat,const MatFactorInfo *info,FactorShiftCtx *sctx,PetscInt row)
706: {
707:   PetscReal _zero = info->zeropivot;
710:   if (PetscAbsScalar(sctx->pv) <= _zero && !PetscIsNanScalar(sctx->pv)){
711:     sctx->pv          += info->shiftamount;
712:     sctx->shift_amount = 0.0;
713:     sctx->nshift++;
714:   }
715:   sctx->newshift = PETSC_FALSE;
716:   return(0);
717: }
719: PETSC_STATIC_INLINE PetscErrorCode MatPivotCheck_none(Mat fact,Mat mat,const MatFactorInfo *info,FactorShiftCtx *sctx,PetscInt row)
720: {
721:   PetscReal      _zero = info->zeropivot;
725:   sctx->newshift = PETSC_FALSE;
726:   if (PetscAbsScalar(sctx->pv) <= _zero && !PetscIsNanScalar(sctx->pv)) {
727:     if (!mat->erroriffailure) {
728:       PetscInfo3(mat,"Detected zero pivot in factorization in row %D value %g tolerance %g\n",row,(double)PetscAbsScalar(sctx->pv),(double)_zero);
729:       fact->factorerrortype             = MAT_FACTOR_NUMERIC_ZEROPIVOT;
730:       fact->factorerror_zeropivot_value = PetscAbsScalar(sctx->pv);
731:       fact->factorerror_zeropivot_row   = row;
732:     } 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);
733:   }
734:   return(0);
735: }
737: PETSC_STATIC_INLINE PetscErrorCode MatPivotCheck(Mat fact,Mat mat,const MatFactorInfo *info,FactorShiftCtx *sctx,PetscInt row)
738: {
742:   if (info->shifttype == (PetscReal) MAT_SHIFT_NONZERO){
743:     MatPivotCheck_nz(mat,info,sctx,row);
744:   } else if (info->shifttype == (PetscReal) MAT_SHIFT_POSITIVE_DEFINITE){
745:     MatPivotCheck_pd(mat,info,sctx,row);
746:   } else if (info->shifttype == (PetscReal) MAT_SHIFT_INBLOCKS){
747:     MatPivotCheck_inblocks(mat,info,sctx,row);
748:   } else {
749:     MatPivotCheck_none(fact,mat,info,sctx,row);
750:   }
751:   return(0);
752: }
754: /*
755:   Create and initialize a linked list
756:   Input Parameters:
757:     idx_start - starting index of the list
758:     lnk_max   - max value of lnk indicating the end of the list
759:     nlnk      - max length of the list
760:   Output Parameters:
761:     lnk       - list initialized
762:     bt        - PetscBT (bitarray) with all bits set to false
763:     lnk_empty - flg indicating the list is empty
764: */
765: #define PetscLLCreate(idx_start,lnk_max,nlnk,lnk,bt) \
766:   (PetscMalloc1(nlnk,&lnk) || PetscBTCreate(nlnk,&(bt)) || (lnk[idx_start] = lnk_max,0))
768: #define PetscLLCreate_new(idx_start,lnk_max,nlnk,lnk,bt,lnk_empty)\
769:   (PetscMalloc1(nlnk,&lnk) || PetscBTCreate(nlnk,&(bt)) || (lnk_empty = PETSC_TRUE,0) ||(lnk[idx_start] = lnk_max,0))
771: /*
772:   Add an index set into a sorted linked list
773:   Input Parameters:
774:     nidx      - number of input indices
775:     indices   - integer array
776:     idx_start - starting index of the list
777:     lnk       - linked list(an integer array) that is created
778:     bt        - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
779:   output Parameters:
780:     nlnk      - number of newly added indices
781:     lnk       - the sorted(increasing order) linked list containing new and non-redundate entries from indices
782:     bt        - updated PetscBT (bitarray)
783: */
784: #define PetscLLAdd(nidx,indices,idx_start,nlnk,lnk,bt) 0;\
785: {\
786:   PetscInt _k,_entry,_location,_lnkdata;\
787:   nlnk     = 0;\
788:   _lnkdata = idx_start;\
789:   for (_k=0; _k<nidx; _k++){\
790:     _entry = indices[_k];\
791:     if (!PetscBTLookupSet(bt,_entry)){  /* new entry */\
792:       /* search for insertion location */\
793:       /* start from the beginning if _entry < previous _entry */\
794:       if (_k && _entry < _lnkdata) _lnkdata  = idx_start;\
795:       do {\
796:         _location = _lnkdata;\
797:         _lnkdata  = lnk[_location];\
798:       } while (_entry > _lnkdata);\
799:       /* insertion location is found, add entry into lnk */\
800:       lnk[_location] = _entry;\
801:       lnk[_entry]    = _lnkdata;\
802:       nlnk++;\
803:       _lnkdata = _entry; /* next search starts from here if next_entry > _entry */\
804:     }\
805:   }\
806: }
808: /*
809:   Add a permuted index set into a sorted linked list
810:   Input Parameters:
811:     nidx      - number of input indices
812:     indices   - integer array
813:     perm      - permutation of indices
814:     idx_start - starting index of the list
815:     lnk       - linked list(an integer array) that is created
816:     bt        - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
817:   output Parameters:
818:     nlnk      - number of newly added indices
819:     lnk       - the sorted(increasing order) linked list containing new and non-redundate entries from indices
820:     bt        - updated PetscBT (bitarray)
821: */
822: #define PetscLLAddPerm(nidx,indices,perm,idx_start,nlnk,lnk,bt) 0;\
823: {\
824:   PetscInt _k,_entry,_location,_lnkdata;\
825:   nlnk     = 0;\
826:   _lnkdata = idx_start;\
827:   for (_k=0; _k<nidx; _k++){\
828:     _entry = perm[indices[_k]];\
829:     if (!PetscBTLookupSet(bt,_entry)){  /* new entry */\
830:       /* search for insertion location */\
831:       /* start from the beginning if _entry < previous _entry */\
832:       if (_k && _entry < _lnkdata) _lnkdata  = idx_start;\
833:       do {\
834:         _location = _lnkdata;\
835:         _lnkdata  = lnk[_location];\
836:       } while (_entry > _lnkdata);\
837:       /* insertion location is found, add entry into lnk */\
838:       lnk[_location] = _entry;\
839:       lnk[_entry]    = _lnkdata;\
840:       nlnk++;\
841:       _lnkdata = _entry; /* next search starts from here if next_entry > _entry */\
842:     }\
843:   }\
844: }
846: /*
847:   Add a SORTED ascending index set into a sorted linked list - same as PetscLLAdd() bus skip 'if (_k && _entry < _lnkdata) _lnkdata  = idx_start;'
848:   Input Parameters:
849:     nidx      - number of input indices
850:     indices   - sorted integer array
851:     idx_start - starting index of the list
852:     lnk       - linked list(an integer array) that is created
853:     bt        - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
854:   output Parameters:
855:     nlnk      - number of newly added indices
856:     lnk       - the sorted(increasing order) linked list containing new and non-redundate entries from indices
857:     bt        - updated PetscBT (bitarray)
858: */
859: #define PetscLLAddSorted(nidx,indices,idx_start,nlnk,lnk,bt) 0;\
860: {\
861:   PetscInt _k,_entry,_location,_lnkdata;\
862:   nlnk      = 0;\
863:   _lnkdata  = idx_start;\
864:   for (_k=0; _k<nidx; _k++){\
865:     _entry = indices[_k];\
866:     if (!PetscBTLookupSet(bt,_entry)){  /* new entry */\
867:       /* search for insertion location */\
868:       do {\
869:         _location = _lnkdata;\
870:         _lnkdata  = lnk[_location];\
871:       } while (_entry > _lnkdata);\
872:       /* insertion location is found, add entry into lnk */\
873:       lnk[_location] = _entry;\
874:       lnk[_entry]    = _lnkdata;\
875:       nlnk++;\
876:       _lnkdata = _entry; /* next search starts from here */\
877:     }\
878:   }\
879: }
881: #define PetscLLAddSorted_new(nidx,indices,idx_start,lnk_empty,nlnk,lnk,bt) 0; \
882: {\
883:   PetscInt _k,_entry,_location,_lnkdata;\
884:   if (lnk_empty){\
885:     _lnkdata  = idx_start;                      \
886:     for (_k=0; _k<nidx; _k++){                  \
887:       _entry = indices[_k];                             \
888:       PetscBTSet(bt,_entry);  /* mark the new entry */          \
889:           _location = _lnkdata;                                 \
890:           _lnkdata  = lnk[_location];                           \
891:         /* insertion location is found, add entry into lnk */   \
892:         lnk[_location] = _entry;                                \
893:         lnk[_entry]    = _lnkdata;                              \
894:         _lnkdata = _entry; /* next search starts from here */   \
895:     }                                                           \
896:     /*\
897:     lnk[indices[nidx-1]] = lnk[idx_start];\
898:     lnk[idx_start]       = indices[0];\
899:     PetscBTSet(bt,indices[0]);  \
900:     for (_k=1; _k<nidx; _k++){                  \
901:       PetscBTSet(bt,indices[_k]);                                          \
902:       lnk[indices[_k-1]] = indices[_k];                                  \
903:     }                                                           \
904:      */\
905:     nlnk      = nidx;\
906:     lnk_empty = PETSC_FALSE;\
907:   } else {\
908:     nlnk      = 0;                              \
909:     _lnkdata  = idx_start;                      \
910:     for (_k=0; _k<nidx; _k++){                  \
911:       _entry = indices[_k];                             \
912:       if (!PetscBTLookupSet(bt,_entry)){  /* new entry */       \
913:         /* search for insertion location */                     \
914:         do {                                                    \
915:           _location = _lnkdata;                                 \
916:           _lnkdata  = lnk[_location];                           \
917:         } while (_entry > _lnkdata);                            \
918:         /* insertion location is found, add entry into lnk */   \
919:         lnk[_location] = _entry;                                \
920:         lnk[_entry]    = _lnkdata;                              \
921:         nlnk++;                                                 \
922:         _lnkdata = _entry; /* next search starts from here */   \
923:       }                                                         \
924:     }                                                           \
925:   }                                                             \
926: }
928: /*
929:   Add a SORTED index set into a sorted linked list used for LUFactorSymbolic()
930:   Same as PetscLLAddSorted() with an additional operation:
931:        count the number of input indices that are no larger than 'diag'
932:   Input Parameters:
933:     indices   - sorted integer array
934:     idx_start - starting index of the list, index of pivot row
935:     lnk       - linked list(an integer array) that is created
936:     bt        - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
937:     diag      - index of the active row in LUFactorSymbolic
938:     nzbd      - number of input indices with indices <= idx_start
939:     im        - im[idx_start] is initialized as num of nonzero entries in row=idx_start
940:   output Parameters:
941:     nlnk      - number of newly added indices
942:     lnk       - the sorted(increasing order) linked list containing new and non-redundate entries from indices
943:     bt        - updated PetscBT (bitarray)
944:     im        - im[idx_start]: unchanged if diag is not an entry
945:                              : num of entries with indices <= diag if diag is an entry
946: */
947: #define PetscLLAddSortedLU(indices,idx_start,nlnk,lnk,bt,diag,nzbd,im) 0;\
948: {\
949:   PetscInt _k,_entry,_location,_lnkdata,_nidx;\
950:   nlnk     = 0;\
951:   _lnkdata = idx_start;\
952:   _nidx = im[idx_start] - nzbd; /* num of entries with idx_start < index <= diag */\
953:   for (_k=0; _k<_nidx; _k++){\
954:     _entry = indices[_k];\
955:     nzbd++;\
956:     if ( _entry== diag) im[idx_start] = nzbd;\
957:     if (!PetscBTLookupSet(bt,_entry)){  /* new entry */\
958:       /* search for insertion location */\
959:       do {\
960:         _location = _lnkdata;\
961:         _lnkdata  = lnk[_location];\
962:       } while (_entry > _lnkdata);\
963:       /* insertion location is found, add entry into lnk */\
964:       lnk[_location] = _entry;\
965:       lnk[_entry]    = _lnkdata;\
966:       nlnk++;\
967:       _lnkdata = _entry; /* next search starts from here */\
968:     }\
969:   }\
970: }
972: /*
973:   Copy data on the list into an array, then initialize the list
974:   Input Parameters:
975:     idx_start - starting index of the list
976:     lnk_max   - max value of lnk indicating the end of the list
977:     nlnk      - number of data on the list to be copied
978:     lnk       - linked list
979:     bt        - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
980:   output Parameters:
981:     indices   - array that contains the copied data
982:     lnk       - linked list that is cleaned and initialize
983:     bt        - PetscBT (bitarray) with all bits set to false
984: */
985: #define PetscLLClean(idx_start,lnk_max,nlnk,lnk,indices,bt) 0;\
986: {\
987:   PetscInt _j,_idx=idx_start;\
988:   for (_j=0; _j<nlnk; _j++){\
989:     _idx = lnk[_idx];\
990:     indices[_j] = _idx;\
991:     PetscBTClear(bt,_idx);\
992:   }\
993:   lnk[idx_start] = lnk_max;\
994: }
995: /*
996:   Free memories used by the list
997: */
998: #define PetscLLDestroy(lnk,bt) (PetscFree(lnk) || PetscBTDestroy(&(bt)))
1000: /* Routines below are used for incomplete matrix factorization */
1001: /*
1002:   Create and initialize a linked list and its levels
1003:   Input Parameters:
1004:     idx_start - starting index of the list
1005:     lnk_max   - max value of lnk indicating the end of the list
1006:     nlnk      - max length of the list
1007:   Output Parameters:
1008:     lnk       - list initialized
1009:     lnk_lvl   - array of size nlnk for storing levels of lnk
1010:     bt        - PetscBT (bitarray) with all bits set to false
1011: */
1012: #define PetscIncompleteLLCreate(idx_start,lnk_max,nlnk,lnk,lnk_lvl,bt)\
1013:   (PetscIntMultError(2,nlnk,NULL) || PetscMalloc1(2*nlnk,&lnk) || PetscBTCreate(nlnk,&(bt)) || (lnk[idx_start] = lnk_max,lnk_lvl = lnk + nlnk,0))
1015: /*
1016:   Initialize a sorted linked list used for ILU and ICC
1017:   Input Parameters:
1018:     nidx      - number of input idx
1019:     idx       - integer array used for storing column indices
1020:     idx_start - starting index of the list
1021:     perm      - indices of an IS
1022:     lnk       - linked list(an integer array) that is created
1023:     lnklvl    - levels of lnk
1024:     bt        - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
1025:   output Parameters:
1026:     nlnk     - number of newly added idx
1027:     lnk      - the sorted(increasing order) linked list containing new and non-redundate entries from idx
1028:     lnklvl   - levels of lnk
1029:     bt       - updated PetscBT (bitarray)
1030: */
1031: #define PetscIncompleteLLInit(nidx,idx,idx_start,perm,nlnk,lnk,lnklvl,bt) 0;\
1032: {\
1033:   PetscInt _k,_entry,_location,_lnkdata;\
1034:   nlnk     = 0;\
1035:   _lnkdata = idx_start;\
1036:   for (_k=0; _k<nidx; _k++){\
1037:     _entry = perm[idx[_k]];\
1038:     if (!PetscBTLookupSet(bt,_entry)){  /* new entry */\
1039:       /* search for insertion location */\
1040:       if (_k && _entry < _lnkdata) _lnkdata  = idx_start;\
1041:       do {\
1042:         _location = _lnkdata;\
1043:         _lnkdata  = lnk[_location];\
1044:       } while (_entry > _lnkdata);\
1045:       /* insertion location is found, add entry into lnk */\
1046:       lnk[_location]  = _entry;\
1047:       lnk[_entry]     = _lnkdata;\
1048:       lnklvl[_entry] = 0;\
1049:       nlnk++;\
1050:       _lnkdata = _entry; /* next search starts from here if next_entry > _entry */\
1051:     }\
1052:   }\
1053: }
1055: /*
1056:   Add a SORTED index set into a sorted linked list for ILU
1057:   Input Parameters:
1058:     nidx      - number of input indices
1059:     idx       - sorted integer array used for storing column indices
1060:     level     - level of fill, e.g., ICC(level)
1061:     idxlvl    - level of idx
1062:     idx_start - starting index of the list
1063:     lnk       - linked list(an integer array) that is created
1064:     lnklvl    - levels of lnk
1065:     bt        - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
1066:     prow      - the row number of idx
1067:   output Parameters:
1068:     nlnk     - number of newly added idx
1069:     lnk      - the sorted(increasing order) linked list containing new and non-redundate entries from idx
1070:     lnklvl   - levels of lnk
1071:     bt       - updated PetscBT (bitarray)
1073:   Note: the level of factor(i,j) is set as lvl(i,j) = min{ lvl(i,j), lvl(i,prow)+lvl(prow,j)+1)
1074:         where idx = non-zero columns of U(prow,prow+1:n-1), prow<i
1075: */
1076: #define PetscILULLAddSorted(nidx,idx,level,idxlvl,idx_start,nlnk,lnk,lnklvl,bt,lnklvl_prow) 0;\
1077: {\
1078:   PetscInt _k,_entry,_location,_lnkdata,_incrlev,_lnklvl_prow=lnklvl[prow];\
1079:   nlnk     = 0;\
1080:   _lnkdata = idx_start;\
1081:   for (_k=0; _k<nidx; _k++){\
1082:     _incrlev = idxlvl[_k] + _lnklvl_prow + 1;\
1083:     if (_incrlev > level) continue;\
1084:     _entry = idx[_k];\
1085:     if (!PetscBTLookupSet(bt,_entry)){  /* new entry */\
1086:       /* search for insertion location */\
1087:       do {\
1088:         _location = _lnkdata;\
1089:         _lnkdata  = lnk[_location];\
1090:       } while (_entry > _lnkdata);\
1091:       /* insertion location is found, add entry into lnk */\
1092:       lnk[_location]  = _entry;\
1093:       lnk[_entry]     = _lnkdata;\
1094:       lnklvl[_entry] = _incrlev;\
1095:       nlnk++;\
1096:       _lnkdata = _entry; /* next search starts from here if next_entry > _entry */\
1097:     } else { /* existing entry: update lnklvl */\
1098:       if (lnklvl[_entry] > _incrlev) lnklvl[_entry] = _incrlev;\
1099:     }\
1100:   }\
1101: }
1103: /*
1104:   Add a index set into a sorted linked list
1105:   Input Parameters:
1106:     nidx      - number of input idx
1107:     idx   - integer array used for storing column indices
1108:     level     - level of fill, e.g., ICC(level)
1109:     idxlvl - level of idx
1110:     idx_start - starting index of the list
1111:     lnk       - linked list(an integer array) that is created
1112:     lnklvl   - levels of lnk
1113:     bt        - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
1114:   output Parameters:
1115:     nlnk      - number of newly added idx
1116:     lnk       - the sorted(increasing order) linked list containing new and non-redundate entries from idx
1117:     lnklvl   - levels of lnk
1118:     bt        - updated PetscBT (bitarray)
1119: */
1120: #define PetscIncompleteLLAdd(nidx,idx,level,idxlvl,idx_start,nlnk,lnk,lnklvl,bt) 0;\
1121: {\
1122:   PetscInt _k,_entry,_location,_lnkdata,_incrlev;\
1123:   nlnk     = 0;\
1124:   _lnkdata = idx_start;\
1125:   for (_k=0; _k<nidx; _k++){\
1126:     _incrlev = idxlvl[_k] + 1;\
1127:     if (_incrlev > level) continue;\
1128:     _entry = idx[_k];\
1129:     if (!PetscBTLookupSet(bt,_entry)){  /* new entry */\
1130:       /* search for insertion location */\
1131:       if (_k && _entry < _lnkdata) _lnkdata  = idx_start;\
1132:       do {\
1133:         _location = _lnkdata;\
1134:         _lnkdata  = lnk[_location];\
1135:       } while (_entry > _lnkdata);\
1136:       /* insertion location is found, add entry into lnk */\
1137:       lnk[_location]  = _entry;\
1138:       lnk[_entry]     = _lnkdata;\
1139:       lnklvl[_entry] = _incrlev;\
1140:       nlnk++;\
1141:       _lnkdata = _entry; /* next search starts from here if next_entry > _entry */\
1142:     } else { /* existing entry: update lnklvl */\
1143:       if (lnklvl[_entry] > _incrlev) lnklvl[_entry] = _incrlev;\
1144:     }\
1145:   }\
1146: }
1148: /*
1149:   Add a SORTED index set into a sorted linked list
1150:   Input Parameters:
1151:     nidx      - number of input indices
1152:     idx   - sorted integer array used for storing column indices
1153:     level     - level of fill, e.g., ICC(level)
1154:     idxlvl - level of idx
1155:     idx_start - starting index of the list
1156:     lnk       - linked list(an integer array) that is created
1157:     lnklvl    - levels of lnk
1158:     bt        - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
1159:   output Parameters:
1160:     nlnk      - number of newly added idx
1161:     lnk       - the sorted(increasing order) linked list containing new and non-redundate entries from idx
1162:     lnklvl    - levels of lnk
1163:     bt        - updated PetscBT (bitarray)
1164: */
1165: #define PetscIncompleteLLAddSorted(nidx,idx,level,idxlvl,idx_start,nlnk,lnk,lnklvl,bt) 0;\
1166: {\
1167:   PetscInt _k,_entry,_location,_lnkdata,_incrlev;\
1168:   nlnk = 0;\
1169:   _lnkdata = idx_start;\
1170:   for (_k=0; _k<nidx; _k++){\
1171:     _incrlev = idxlvl[_k] + 1;\
1172:     if (_incrlev > level) continue;\
1173:     _entry = idx[_k];\
1174:     if (!PetscBTLookupSet(bt,_entry)){  /* new entry */\
1175:       /* search for insertion location */\
1176:       do {\
1177:         _location = _lnkdata;\
1178:         _lnkdata  = lnk[_location];\
1179:       } while (_entry > _lnkdata);\
1180:       /* insertion location is found, add entry into lnk */\
1181:       lnk[_location] = _entry;\
1182:       lnk[_entry]    = _lnkdata;\
1183:       lnklvl[_entry] = _incrlev;\
1184:       nlnk++;\
1185:       _lnkdata = _entry; /* next search starts from here */\
1186:     } else { /* existing entry: update lnklvl */\
1187:       if (lnklvl[_entry] > _incrlev) lnklvl[_entry] = _incrlev;\
1188:     }\
1189:   }\
1190: }
1192: /*
1193:   Add a SORTED index set into a sorted linked list for ICC
1194:   Input Parameters:
1195:     nidx      - number of input indices
1196:     idx       - sorted integer array used for storing column indices
1197:     level     - level of fill, e.g., ICC(level)
1198:     idxlvl    - level of idx
1199:     idx_start - starting index of the list
1200:     lnk       - linked list(an integer array) that is created
1201:     lnklvl    - levels of lnk
1202:     bt        - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
1203:     idxlvl_prow - idxlvl[prow], where prow is the row number of the idx
1204:   output Parameters:
1205:     nlnk   - number of newly added indices
1206:     lnk    - the sorted(increasing order) linked list containing new and non-redundate entries from idx
1207:     lnklvl - levels of lnk
1208:     bt     - updated PetscBT (bitarray)
1209:   Note: the level of U(i,j) is set as lvl(i,j) = min{ lvl(i,j), lvl(prow,i)+lvl(prow,j)+1)
1210:         where idx = non-zero columns of U(prow,prow+1:n-1), prow<i
1211: */
1212: #define PetscICCLLAddSorted(nidx,idx,level,idxlvl,idx_start,nlnk,lnk,lnklvl,bt,idxlvl_prow) 0;\
1213: {\
1214:   PetscInt _k,_entry,_location,_lnkdata,_incrlev;\
1215:   nlnk = 0;\
1216:   _lnkdata = idx_start;\
1217:   for (_k=0; _k<nidx; _k++){\
1218:     _incrlev = idxlvl[_k] + idxlvl_prow + 1;\
1219:     if (_incrlev > level) continue;\
1220:     _entry = idx[_k];\
1221:     if (!PetscBTLookupSet(bt,_entry)){  /* new entry */\
1222:       /* search for insertion location */\
1223:       do {\
1224:         _location = _lnkdata;\
1225:         _lnkdata  = lnk[_location];\
1226:       } while (_entry > _lnkdata);\
1227:       /* insertion location is found, add entry into lnk */\
1228:       lnk[_location] = _entry;\
1229:       lnk[_entry]    = _lnkdata;\
1230:       lnklvl[_entry] = _incrlev;\
1231:       nlnk++;\
1232:       _lnkdata = _entry; /* next search starts from here */\
1233:     } else { /* existing entry: update lnklvl */\
1234:       if (lnklvl[_entry] > _incrlev) lnklvl[_entry] = _incrlev;\
1235:     }\
1236:   }\
1237: }
1239: /*
1240:   Copy data on the list into an array, then initialize the list
1241:   Input Parameters:
1242:     idx_start - starting index of the list
1243:     lnk_max   - max value of lnk indicating the end of the list
1244:     nlnk      - number of data on the list to be copied
1245:     lnk       - linked list
1246:     lnklvl    - level of lnk
1247:     bt        - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
1248:   output Parameters:
1249:     indices - array that contains the copied data
1250:     lnk     - linked list that is cleaned and initialize
1251:     lnklvl  - level of lnk that is reinitialized
1252:     bt      - PetscBT (bitarray) with all bits set to false
1253: */
1254: #define PetscIncompleteLLClean(idx_start,lnk_max,nlnk,lnk,lnklvl,indices,indiceslvl,bt) 0;\
1255: do {\
1256:   PetscInt _j,_idx=idx_start;\
1257:   for (_j=0; _j<nlnk; _j++){\
1258:     _idx = lnk[_idx];\
1259:     *(indices+_j) = _idx;\
1260:     *(indiceslvl+_j) = lnklvl[_idx];\
1261:     lnklvl[_idx] = -1;\
1262:     PetscBTClear(bt,_idx);\
1263:   }\
1264:   lnk[idx_start] = lnk_max;\
1265: } while (0)
1266: /*
1267:   Free memories used by the list
1268: */
1269: #define PetscIncompleteLLDestroy(lnk,bt) (PetscFree(lnk) || PetscBTDestroy(&(bt)))
1271: #define MatCheckSameLocalSize(A,ar1,B,ar2) do { \
1273:   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)
1275: #define MatCheckSameSize(A,ar1,B,ar2) do { \
1276:   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);\
1277:   MatCheckSameLocalSize(A,ar1,B,ar2);} while (0)
1279: #define VecCheckMatCompatible(M,x,ar1,b,ar2) do { \
1280:   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); \
1281:   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)
1283: /* -------------------------------------------------------------------------------------------------------*/
1284:  #include <petscbt.h>
1285: /*
1286:   Create and initialize a condensed linked list -
1287:     same as PetscLLCreate(), but uses a scalable array 'lnk' with size of max number of entries, not O(N).
1288:     Barry suggested this approach (Dec. 6, 2011):
1289:       I've thought of an alternative way of representing a linked list that is efficient but doesn't have the O(N) scaling issue
1290:       (it may be faster than the O(N) even sequentially due to less crazy memory access).
1292:       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
1293:       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
1294:       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
1295:       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.
1296:       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
1297:       to each other so memory access is much better than using the big array.
1299:   Example:
1300:      nlnk_max=5, lnk_max=36:
1301:      Initial list: [0, 0 | 36, 2 | 0, 0 | 0, 0 | 0, 0 | 0, 0 | 0, 0]
1302:      here, head_node has index 2 with value lnk[2]=lnk_max=36,
1303:            0-th entry is used to store the number of entries in the list,
1304:      The initial lnk represents head -> tail(marked by 36) with number of entries = lnk[0]=0.
1306:      Now adding a sorted set {2,4}, the list becomes
1307:      [2, 0 | 36, 4 |2, 6 | 4, 2 | 0, 0 | 0, 0 | 0, 0 ]
1308:      represents head -> 2 -> 4 -> tail with number of entries = lnk[0]=2.
1310:      Then adding a sorted set {0,3,35}, the list
1311:      [5, 0 | 36, 8 | 2, 10 | 4, 12 | 0, 4 | 3, 6 | 35, 2 ]
1312:      represents head -> 0 -> 2 -> 3 -> 4 -> 35 -> tail with number of entries = lnk[0]=5.
1314:   Input Parameters:
1315:     nlnk_max  - max length of the list
1316:     lnk_max   - max value of the entries
1317:   Output Parameters:
1318:     lnk       - list created and initialized
1319:     bt        - PetscBT (bitarray) with all bits set to false. Note: bt has size lnk_max, not nln_max!
1320: */
1321: PETSC_STATIC_INLINE PetscErrorCode PetscLLCondensedCreate(PetscInt nlnk_max,PetscInt lnk_max,PetscInt **lnk,PetscBT *bt)
1322: {
1324:   PetscInt       *llnk,lsize = 0;
1327:   PetscIntMultError(2,nlnk_max+2,&lsize);
1328:   PetscMalloc1(lsize,lnk);
1329:   PetscBTCreate(lnk_max,bt);
1330:   llnk = *lnk;
1331:   llnk[0] = 0;         /* number of entries on the list */
1332:   llnk[2] = lnk_max;   /* value in the head node */
1333:   llnk[3] = 2;         /* next for the head node */
1334:   return(0);
1335: }
1337: /*
1338:   Add a SORTED ascending index set into a sorted linked list. See PetscLLCondensedCreate() for detailed description.
1339:   Input Parameters:
1340:     nidx      - number of input indices
1341:     indices   - sorted integer array
1342:     lnk       - condensed linked list(an integer array) that is created
1343:     bt        - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
1344:   output Parameters:
1345:     lnk       - the sorted(increasing order) linked list containing previous and newly added non-redundate indices
1346:     bt        - updated PetscBT (bitarray)
1347: */
1348: PETSC_STATIC_INLINE PetscErrorCode PetscLLCondensedAddSorted(PetscInt nidx,const PetscInt indices[],PetscInt lnk[],PetscBT bt)
1349: {
1350:   PetscInt _k,_entry,_location,_next,_lnkdata,_nlnk,_newnode;
1353:   _nlnk     = lnk[0]; /* num of entries on the input lnk */
1354:   _location = 2; /* head */
1355:     for (_k=0; _k<nidx; _k++){
1356:       _entry = indices[_k];
1357:       if (!PetscBTLookupSet(bt,_entry)){  /* new entry */
1358:         /* search for insertion location */
1359:         do {
1360:           _next     = _location + 1; /* link from previous node to next node */
1361:           _location = lnk[_next];    /* idx of next node */
1362:           _lnkdata  = lnk[_location];/* value of next node */
1363:         } while (_entry > _lnkdata);
1364:         /* insertion location is found, add entry into lnk */
1365:         _newnode        = 2*(_nlnk+2);   /* index for this new node */
1366:         lnk[_next]      = _newnode;      /* connect previous node to the new node */
1367:         lnk[_newnode]   = _entry;        /* set value of the new node */
1368:         lnk[_newnode+1] = _location;     /* connect new node to next node */
1369:         _location       = _newnode;      /* next search starts from the new node */
1370:         _nlnk++;
1371:       }   \
1372:     }\
1373:   lnk[0]   = _nlnk;   /* number of entries in the list */
1374:   return(0);
1375: }
1377: PETSC_STATIC_INLINE PetscErrorCode PetscLLCondensedClean(PetscInt lnk_max,PetscInt nidx,PetscInt *indices,PetscInt lnk[],PetscBT bt)
1378: {
1380:   PetscInt       _k,_next,_nlnk;
1383:   _next = lnk[3];       /* head node */
1384:   _nlnk = lnk[0];       /* num of entries on the list */
1385:   for (_k=0; _k<_nlnk; _k++){
1386:     indices[_k] = lnk[_next];
1387:     _next       = lnk[_next + 1];
1388:     PetscBTClear(bt,indices[_k]);
1389:   }
1390:   lnk[0] = 0;          /* num of entries on the list */
1391:   lnk[2] = lnk_max;    /* initialize head node */
1392:   lnk[3] = 2;          /* head node */
1393:   return(0);
1394: }
1396: PETSC_STATIC_INLINE PetscErrorCode PetscLLCondensedView(PetscInt *lnk)
1397: {
1399:   PetscInt       k;
1402:   PetscPrintf(PETSC_COMM_SELF,"LLCondensed of size %D, (val,  next)\n",lnk[0]);
1403:   for (k=2; k< lnk[0]+2; k++){
1404:     PetscPrintf(PETSC_COMM_SELF," %D: (%D, %D)\n",2*k,lnk[2*k],lnk[2*k+1]);
1405:   }
1406:   return(0);
1407: }
1409: /*
1410:   Free memories used by the list
1411: */
1412: PETSC_STATIC_INLINE PetscErrorCode PetscLLCondensedDestroy(PetscInt *lnk,PetscBT bt)
1413: {
1417:   PetscFree(lnk);
1418:   PetscBTDestroy(&bt);
1419:   return(0);
1420: }
1422: /* -------------------------------------------------------------------------------------------------------*/
1423: /*
1424:  Same as PetscLLCondensedCreate(), but does not use non-scalable O(lnk_max) bitarray
1425:   Input Parameters:
1426:     nlnk_max  - max length of the list
1427:   Output Parameters:
1428:     lnk       - list created and initialized
1429: */
1430: PETSC_STATIC_INLINE PetscErrorCode PetscLLCondensedCreate_Scalable(PetscInt nlnk_max,PetscInt **lnk)
1431: {
1433:   PetscInt       *llnk,lsize = 0;
1436:   PetscIntMultError(2,nlnk_max+2,&lsize);
1437:   PetscMalloc1(lsize,lnk);
1438:   llnk = *lnk;
1439:   llnk[0] = 0;               /* number of entries on the list */
1440:   llnk[2] = PETSC_MAX_INT;   /* value in the head node */
1441:   llnk[3] = 2;               /* next for the head node */
1442:   return(0);
1443: }
1445: PETSC_STATIC_INLINE PetscErrorCode PetscLLCondensedExpand_Scalable(PetscInt nlnk_max,PetscInt **lnk)
1446: {
1448:   PetscInt       lsize = 0;
1451:   PetscIntMultError(2,nlnk_max+2,&lsize);
1452:   PetscRealloc(lsize*sizeof(PetscInt),lnk);
1453:   return(0);
1454: }
1456: PETSC_STATIC_INLINE PetscErrorCode PetscLLCondensedAddSorted_Scalable(PetscInt nidx,const PetscInt indices[],PetscInt lnk[])
1457: {
1458:   PetscInt _k,_entry,_location,_next,_lnkdata,_nlnk,_newnode;
1459:   _nlnk     = lnk[0]; /* num of entries on the input lnk */
1460:   _location = 2; /* head */ \
1461:     for (_k=0; _k<nidx; _k++){
1462:       _entry = indices[_k];
1463:       /* search for insertion location */
1464:       do {
1465:         _next     = _location + 1; /* link from previous node to next node */
1466:         _location = lnk[_next];    /* idx of next node */
1467:         _lnkdata  = lnk[_location];/* value of next node */
1468:       } while (_entry > _lnkdata);
1469:       if (_entry < _lnkdata) {
1470:         /* insertion location is found, add entry into lnk */
1471:         _newnode        = 2*(_nlnk+2);   /* index for this new node */
1472:         lnk[_next]      = _newnode;      /* connect previous node to the new node */
1473:         lnk[_newnode]   = _entry;        /* set value of the new node */
1474:         lnk[_newnode+1] = _location;     /* connect new node to next node */
1475:         _location       = _newnode;      /* next search starts from the new node */
1476:         _nlnk++;
1477:       }
1478:     }
1479:   lnk[0]   = _nlnk;   /* number of entries in the list */
1480:   return 0;
1481: }
1483: PETSC_STATIC_INLINE PetscErrorCode PetscLLCondensedClean_Scalable(PetscInt nidx,PetscInt *indices,PetscInt *lnk)
1484: {
1485:   PetscInt _k,_next,_nlnk;
1486:   _next = lnk[3];       /* head node */
1487:   _nlnk = lnk[0];
1488:   for (_k=0; _k<_nlnk; _k++){
1489:     indices[_k] = lnk[_next];
1490:     _next       = lnk[_next + 1];
1491:   }
1492:   lnk[0] = 0;          /* num of entries on the list */
1493:   lnk[3] = 2;          /* head node */
1494:   return 0;
1495: }
1497: PETSC_STATIC_INLINE PetscErrorCode PetscLLCondensedDestroy_Scalable(PetscInt *lnk)
1498: {
1499:   return PetscFree(lnk);
1500: }
1502: /* -------------------------------------------------------------------------------------------------------*/
1503: /*
1504:       lnk[0]   number of links
1505:       lnk[1]   number of entries
1506:       lnk[3n]  value
1507:       lnk[3n+1] len
1508:       lnk[3n+2] link to next value
1510:       The next three are always the first link
1512:       lnk[3]    PETSC_MIN_INT+1
1513:       lnk[4]    1
1514:       lnk[5]    link to first real entry
1516:       The next three are always the last link
1518:       lnk[6]    PETSC_MAX_INT - 1
1519:       lnk[7]    1
1520:       lnk[8]    next valid link (this is the same as lnk[0] but without the decreases)
1521: */
1523: PETSC_STATIC_INLINE PetscErrorCode PetscLLCondensedCreate_fast(PetscInt nlnk_max,PetscInt **lnk)
1524: {
1526:   PetscInt       *llnk,lsize = 0;
1529:   PetscIntMultError(3,nlnk_max+3,&lsize);
1530:   PetscMalloc1(lsize,lnk);
1531:   llnk = *lnk;
1532:   llnk[0] = 0;   /* nlnk: number of entries on the list */
1533:   llnk[1] = 0;          /* number of integer entries represented in list */
1534:   llnk[3] = PETSC_MIN_INT+1;   /* value in the first node */
1535:   llnk[4] = 1;           /* count for the first node */
1536:   llnk[5] = 6;         /* next for the first node */
1537:   llnk[6] = PETSC_MAX_INT-1;   /* value in the last node */
1538:   llnk[7] = 1;           /* count for the last node */
1539:   llnk[8] = 0;         /* next valid node to be used */
1540:   return(0);
1541: }
1543: PETSC_STATIC_INLINE PetscErrorCode PetscLLCondensedAddSorted_fast(PetscInt nidx,const PetscInt indices[],PetscInt lnk[])
1544: {
1545:   PetscInt k,entry,prev,next;
1546:   prev      = 3;      /* first value */
1547:   next      = lnk[prev+2];
1548:   for (k=0; k<nidx; k++){
1549:     entry = indices[k];
1550:     /* search for insertion location */
1551:     while (entry >= lnk[next]) {
1552:       prev = next;
1553:       next = lnk[next+2];
1554:     }
1555:     /* entry is in range of previous list */
1556:     if (entry < lnk[prev]+lnk[prev+1]) continue;
1557:     lnk[1]++;
1558:     /* entry is right after previous list */
1559:     if (entry == lnk[prev]+lnk[prev+1]) {
1560:       lnk[prev+1]++;
1561:       if (lnk[next] == entry+1) { /* combine two contiguous strings */
1562:         lnk[prev+1] += lnk[next+1];
1563:         lnk[prev+2]  = lnk[next+2];
1564:         next         = lnk[next+2];
1565:         lnk[0]--;
1566:       }
1567:       continue;
1568:     }
1569:     /* entry is right before next list */
1570:     if (entry == lnk[next]-1) {
1571:       lnk[next]--;
1572:       lnk[next+1]++;
1573:       prev = next;
1574:       next = lnk[prev+2];
1575:       continue;
1576:     }
1577:     /*  add entry into lnk */
1578:     lnk[prev+2]    = 3*((lnk[8]++)+3);      /* connect previous node to the new node */
1579:     prev           = lnk[prev+2];
1580:     lnk[prev]      = entry;        /* set value of the new node */
1581:     lnk[prev+1]    = 1;             /* number of values in contiguous string is one to start */
1582:     lnk[prev+2]    = next;          /* connect new node to next node */
1583:     lnk[0]++;
1584:   }
1585:   return 0;
1586: }
1588: PETSC_STATIC_INLINE PetscErrorCode PetscLLCondensedClean_fast(PetscInt nidx,PetscInt *indices,PetscInt *lnk)
1589: {
1590:   PetscInt _k,_next,_nlnk,cnt,j;
1591:   _next = lnk[5];       /* first node */
1592:   _nlnk = lnk[0];
1593:   cnt   = 0;
1594:   for (_k=0; _k<_nlnk; _k++){
1595:     for (j=0; j<lnk[_next+1]; j++) {
1596:       indices[cnt++] = lnk[_next] + j;
1597:     }
1598:     _next       = lnk[_next + 2];
1599:   }
1600:   lnk[0] = 0;   /* nlnk: number of links */
1601:   lnk[1] = 0;          /* number of integer entries represented in list */
1602:   lnk[3] = PETSC_MIN_INT+1;   /* value in the first node */
1603:   lnk[4] = 1;           /* count for the first node */
1604:   lnk[5] = 6;         /* next for the first node */
1605:   lnk[6] = PETSC_MAX_INT-1;   /* value in the last node */
1606:   lnk[7] = 1;           /* count for the last node */
1607:   lnk[8] = 0;         /* next valid location to make link */
1608:   return 0;
1609: }
1611: PETSC_STATIC_INLINE PetscErrorCode PetscLLCondensedView_fast(PetscInt *lnk)
1612: {
1613:   PetscInt k,next,nlnk;
1614:   next = lnk[5];       /* first node */
1615:   nlnk = lnk[0];
1616:   for (k=0; k<nlnk; k++){
1617: #if 0                           /* Debugging code */
1618:     printf("%d value %d len %d next %d\n",next,lnk[next],lnk[next+1],lnk[next+2]);
1619: #endif
1620:     next = lnk[next + 2];
1621:   }
1622:   return 0;
1623: }
1625: PETSC_STATIC_INLINE PetscErrorCode PetscLLCondensedDestroy_fast(PetscInt *lnk)
1626: {
1627:   return PetscFree(lnk);
1628: }
1630: /* this is extern because it is used in MatFDColoringUseDM() which is in the DM library */
1631: PETSC_EXTERN PetscErrorCode MatFDColoringApply_AIJ(Mat,MatFDColoring,Vec,void*);
1633: PETSC_EXTERN PetscLogEvent MAT_Mult;
1634: PETSC_EXTERN PetscLogEvent MAT_MultMatrixFree;
1635: PETSC_EXTERN PetscLogEvent MAT_Mults;
1636: PETSC_EXTERN PetscLogEvent MAT_MultConstrained;
1637: PETSC_EXTERN PetscLogEvent MAT_MultAdd;
1638: PETSC_EXTERN PetscLogEvent MAT_MultTranspose;
1639: PETSC_EXTERN PetscLogEvent MAT_MultTransposeConstrained;
1640: PETSC_EXTERN PetscLogEvent MAT_MultTransposeAdd;
1641: PETSC_EXTERN PetscLogEvent MAT_Solve;
1642: PETSC_EXTERN PetscLogEvent MAT_Solves;
1643: PETSC_EXTERN PetscLogEvent MAT_SolveAdd;
1644: PETSC_EXTERN PetscLogEvent MAT_SolveTranspose;
1645: PETSC_EXTERN PetscLogEvent MAT_SolveTransposeAdd;
1646: PETSC_EXTERN PetscLogEvent MAT_SOR;
1647: PETSC_EXTERN PetscLogEvent MAT_ForwardSolve;
1648: PETSC_EXTERN PetscLogEvent MAT_BackwardSolve;
1649: PETSC_EXTERN PetscLogEvent MAT_LUFactor;
1650: PETSC_EXTERN PetscLogEvent MAT_LUFactorSymbolic;
1651: PETSC_EXTERN PetscLogEvent MAT_LUFactorNumeric;
1652: PETSC_EXTERN PetscLogEvent MAT_CholeskyFactor;
1653: PETSC_EXTERN PetscLogEvent MAT_CholeskyFactorSymbolic;
1654: PETSC_EXTERN PetscLogEvent MAT_CholeskyFactorNumeric;
1655: PETSC_EXTERN PetscLogEvent MAT_ILUFactor;
1656: PETSC_EXTERN PetscLogEvent MAT_ILUFactorSymbolic;
1657: PETSC_EXTERN PetscLogEvent MAT_ICCFactorSymbolic;
1658: PETSC_EXTERN PetscLogEvent MAT_Copy;
1659: PETSC_EXTERN PetscLogEvent MAT_Convert;
1660: PETSC_EXTERN PetscLogEvent MAT_Scale;
1661: PETSC_EXTERN PetscLogEvent MAT_AssemblyBegin;
1662: PETSC_EXTERN PetscLogEvent MAT_AssemblyEnd;
1663: PETSC_EXTERN PetscLogEvent MAT_SetValues;
1664: PETSC_EXTERN PetscLogEvent MAT_GetValues;
1665: PETSC_EXTERN PetscLogEvent MAT_GetRow;
1666: PETSC_EXTERN PetscLogEvent MAT_GetRowIJ;
1667: PETSC_EXTERN PetscLogEvent MAT_CreateSubMats;
1668: PETSC_EXTERN PetscLogEvent MAT_GetColoring;
1669: PETSC_EXTERN PetscLogEvent MAT_GetOrdering;
1670: PETSC_EXTERN PetscLogEvent MAT_RedundantMat;
1671: PETSC_EXTERN PetscLogEvent MAT_IncreaseOverlap;
1672: PETSC_EXTERN PetscLogEvent MAT_Partitioning;
1673: PETSC_EXTERN PetscLogEvent MAT_PartitioningND;
1674: PETSC_EXTERN PetscLogEvent MAT_Coarsen;
1675: PETSC_EXTERN PetscLogEvent MAT_ZeroEntries;
1676: PETSC_EXTERN PetscLogEvent MAT_Load;
1677: PETSC_EXTERN PetscLogEvent MAT_View;
1678: PETSC_EXTERN PetscLogEvent MAT_AXPY;
1679: PETSC_EXTERN PetscLogEvent MAT_FDColoringCreate;
1680: PETSC_EXTERN PetscLogEvent MAT_TransposeColoringCreate;
1681: PETSC_EXTERN PetscLogEvent MAT_FDColoringSetUp;
1682: PETSC_EXTERN PetscLogEvent MAT_FDColoringApply;
1683: PETSC_EXTERN PetscLogEvent MAT_Transpose;
1684: PETSC_EXTERN PetscLogEvent MAT_FDColoringFunction;
1685: PETSC_EXTERN PetscLogEvent MAT_CreateSubMat;
1686: PETSC_EXTERN PetscLogEvent MAT_MatMult;
1687: PETSC_EXTERN PetscLogEvent MAT_MatSolve;
1688: PETSC_EXTERN PetscLogEvent MAT_MatTrSolve;
1689: PETSC_EXTERN PetscLogEvent MAT_MatMultSymbolic;
1690: PETSC_EXTERN PetscLogEvent MAT_MatMultNumeric;
1691: PETSC_EXTERN PetscLogEvent MAT_Getlocalmatcondensed;
1692: PETSC_EXTERN PetscLogEvent MAT_GetBrowsOfAcols;
1693: PETSC_EXTERN PetscLogEvent MAT_GetBrowsOfAocols;
1694: PETSC_EXTERN PetscLogEvent MAT_PtAP;
1695: PETSC_EXTERN PetscLogEvent MAT_PtAPSymbolic;
1696: PETSC_EXTERN PetscLogEvent MAT_PtAPNumeric;
1697: PETSC_EXTERN PetscLogEvent MAT_Seqstompinum;
1698: PETSC_EXTERN PetscLogEvent MAT_Seqstompisym;
1699: PETSC_EXTERN PetscLogEvent MAT_Seqstompi;
1700: PETSC_EXTERN PetscLogEvent MAT_Getlocalmat;
1701: PETSC_EXTERN PetscLogEvent MAT_RARt;
1702: PETSC_EXTERN PetscLogEvent MAT_RARtSymbolic;
1703: PETSC_EXTERN PetscLogEvent MAT_RARtNumeric;
1704: PETSC_EXTERN PetscLogEvent MAT_MatTransposeMult;
1705: PETSC_EXTERN PetscLogEvent MAT_MatTransposeMultSymbolic;
1706: PETSC_EXTERN PetscLogEvent MAT_MatTransposeMultNumeric;
1707: PETSC_EXTERN PetscLogEvent MAT_TransposeMatMult;
1708: PETSC_EXTERN PetscLogEvent MAT_TransposeMatMultSymbolic;
1709: PETSC_EXTERN PetscLogEvent MAT_TransposeMatMultNumeric;
1710: PETSC_EXTERN PetscLogEvent MAT_MatMatMult;
1711: PETSC_EXTERN PetscLogEvent MAT_MatMatMultSymbolic;
1712: PETSC_EXTERN PetscLogEvent MAT_MatMatMultNumeric;
1713: PETSC_EXTERN PetscLogEvent MAT_Applypapt;
1714: PETSC_EXTERN PetscLogEvent MAT_Applypapt_symbolic;
1715: PETSC_EXTERN PetscLogEvent MAT_Applypapt_numeric;
1716: PETSC_EXTERN PetscLogEvent MAT_Getsymtranspose;
1717: PETSC_EXTERN PetscLogEvent MAT_Getsymtransreduced;
1718: PETSC_EXTERN PetscLogEvent MAT_GetSequentialNonzeroStructure;
1719: PETSC_EXTERN PetscLogEvent MATMFFD_Mult;
1720: PETSC_EXTERN PetscLogEvent MAT_GetMultiProcBlock;
1721: PETSC_EXTERN PetscLogEvent MAT_CUSPARSECopyToGPU;
1722: PETSC_EXTERN PetscLogEvent MAT_SetValuesBatch;
1723: PETSC_EXTERN PetscLogEvent MAT_ViennaCLCopyToGPU;
1724: PETSC_EXTERN PetscLogEvent MAT_DenseCopyToGPU;
1725: PETSC_EXTERN PetscLogEvent MAT_DenseCopyFromGPU;
1726: PETSC_EXTERN PetscLogEvent MAT_Merge;
1727: PETSC_EXTERN PetscLogEvent MAT_Residual;
1728: PETSC_EXTERN PetscLogEvent MAT_SetRandom;
1729: PETSC_EXTERN PetscLogEvent MAT_FactorFactS;
1730: PETSC_EXTERN PetscLogEvent MAT_FactorInvS;
1731: PETSC_EXTERN PetscLogEvent MATCOLORING_Apply;
1732: PETSC_EXTERN PetscLogEvent MATCOLORING_Comm;
1733: PETSC_EXTERN PetscLogEvent MATCOLORING_Local;
1734: PETSC_EXTERN PetscLogEvent MATCOLORING_ISCreate;
1735: PETSC_EXTERN PetscLogEvent MATCOLORING_SetUp;
1736: PETSC_EXTERN PetscLogEvent MATCOLORING_Weights;
1738: #endif