Actual source code: matimpl.h

petsc-3.9.4 2018-09-11
Report Typos and Errors

  2: #ifndef __MATIMPL_H

  5:  #include <petscmat.h>
  6:  #include <petscmatcoarsen.h>
  7:  #include <petsc/private/petscimpl.h>

  9: PETSC_EXTERN PetscBool MatRegisterAllCalled;
 10: PETSC_EXTERN PetscBool MatSeqAIJRegisterAllCalled;
 11: PETSC_EXTERN PetscBool MatOrderingRegisterAllCalled;
 12: PETSC_EXTERN PetscBool MatColoringRegisterAllCalled;
 13: PETSC_EXTERN PetscBool MatPartitioningRegisterAllCalled;
 14: PETSC_EXTERN PetscBool MatCoarsenRegisterAllCalled;
 15: PETSC_EXTERN PetscErrorCode MatRegisterAll(void);
 16: PETSC_EXTERN PetscErrorCode MatOrderingRegisterAll(void);
 17: PETSC_EXTERN PetscErrorCode MatColoringRegisterAll(void);
 18: PETSC_EXTERN PetscErrorCode MatPartitioningRegisterAll(void);
 19: PETSC_EXTERN PetscErrorCode MatCoarsenRegisterAll(void);
 20: PETSC_EXTERN PetscErrorCode MatSeqAIJRegisterAll(void);

 22: /*
 23:   This file defines the parts of the matrix data structure that are
 24:   shared by all matrix types.
 25: */

 27: /*
 28:     If you add entries here also add them to the MATOP enum
 29:     in include/petscmat.h and include/petsc/finclude/petscmat.h
 30: */
 31: typedef struct _MatOps *MatOps;
 32: struct _MatOps {
 33:   /* 0*/
 34:   PetscErrorCode (*setvalues)(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
 35:   PetscErrorCode (*getrow)(Mat,PetscInt,PetscInt *,PetscInt*[],PetscScalar*[]);
 36:   PetscErrorCode (*restorerow)(Mat,PetscInt,PetscInt *,PetscInt *[],PetscScalar *[]);
 37:   PetscErrorCode (*mult)(Mat,Vec,Vec);
 38:   PetscErrorCode (*multadd)(Mat,Vec,Vec,Vec);
 39:   /* 5*/
 40:   PetscErrorCode (*multtranspose)(Mat,Vec,Vec);
 41:   PetscErrorCode (*multtransposeadd)(Mat,Vec,Vec,Vec);
 42:   PetscErrorCode (*solve)(Mat,Vec,Vec);
 43:   PetscErrorCode (*solveadd)(Mat,Vec,Vec,Vec);
 44:   PetscErrorCode (*solvetranspose)(Mat,Vec,Vec);
 45:   /*10*/
 46:   PetscErrorCode (*solvetransposeadd)(Mat,Vec,Vec,Vec);
 47:   PetscErrorCode (*lufactor)(Mat,IS,IS,const MatFactorInfo*);
 48:   PetscErrorCode (*choleskyfactor)(Mat,IS,const MatFactorInfo*);
 49:   PetscErrorCode (*sor)(Mat,Vec,PetscReal,MatSORType,PetscReal,PetscInt,PetscInt,Vec);
 50:   PetscErrorCode (*transpose)(Mat,MatReuse,Mat *);
 51:   /*15*/
 52:   PetscErrorCode (*getinfo)(Mat,MatInfoType,MatInfo*);
 53:   PetscErrorCode (*equal)(Mat,Mat,PetscBool  *);
 54:   PetscErrorCode (*getdiagonal)(Mat,Vec);
 55:   PetscErrorCode (*diagonalscale)(Mat,Vec,Vec);
 56:   PetscErrorCode (*norm)(Mat,NormType,PetscReal*);
 57:   /*20*/
 58:   PetscErrorCode (*assemblybegin)(Mat,MatAssemblyType);
 59:   PetscErrorCode (*assemblyend)(Mat,MatAssemblyType);
 60:   PetscErrorCode (*setoption)(Mat,MatOption,PetscBool );
 61:   PetscErrorCode (*zeroentries)(Mat);
 62:   /*24*/
 63:   PetscErrorCode (*zerorows)(Mat,PetscInt,const PetscInt[],PetscScalar,Vec,Vec);
 64:   PetscErrorCode (*lufactorsymbolic)(Mat,Mat,IS,IS,const MatFactorInfo*);
 65:   PetscErrorCode (*lufactornumeric)(Mat,Mat,const MatFactorInfo*);
 66:   PetscErrorCode (*choleskyfactorsymbolic)(Mat,Mat,IS,const MatFactorInfo*);
 67:   PetscErrorCode (*choleskyfactornumeric)(Mat,Mat,const MatFactorInfo*);
 68:   /*29*/
 69:   PetscErrorCode (*setup)(Mat);
 70:   PetscErrorCode (*ilufactorsymbolic)(Mat,Mat,IS,IS,const MatFactorInfo*);
 71:   PetscErrorCode (*iccfactorsymbolic)(Mat,Mat,IS,const MatFactorInfo*);
 72:   PetscErrorCode (*getdiagonalblock)(Mat,Mat*);
 73:   PetscErrorCode (*placeholder_33)(Mat);
 74:   /*34*/
 75:   PetscErrorCode (*duplicate)(Mat,MatDuplicateOption,Mat*);
 76:   PetscErrorCode (*forwardsolve)(Mat,Vec,Vec);
 77:   PetscErrorCode (*backwardsolve)(Mat,Vec,Vec);
 78:   PetscErrorCode (*ilufactor)(Mat,IS,IS,const MatFactorInfo*);
 79:   PetscErrorCode (*iccfactor)(Mat,IS,const MatFactorInfo*);
 80:   /*39*/
 81:   PetscErrorCode (*axpy)(Mat,PetscScalar,Mat,MatStructure);
 82:   PetscErrorCode (*createsubmatrices)(Mat,PetscInt,const IS[],const IS[],MatReuse,Mat *[]);
 83:   PetscErrorCode (*increaseoverlap)(Mat,PetscInt,IS[],PetscInt);
 84:   PetscErrorCode (*getvalues)(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],PetscScalar []);
 85:   PetscErrorCode (*copy)(Mat,Mat,MatStructure);
 86:   /*44*/
 87:   PetscErrorCode (*getrowmax)(Mat,Vec,PetscInt[]);
 88:   PetscErrorCode (*scale)(Mat,PetscScalar);
 89:   PetscErrorCode (*shift)(Mat,PetscScalar);
 90:   PetscErrorCode (*diagonalset)(Mat,Vec,InsertMode);
 91:   PetscErrorCode (*zerorowscolumns)(Mat,PetscInt,const PetscInt[],PetscScalar,Vec,Vec);
 92:   /*49*/
 93:   PetscErrorCode (*setrandom)(Mat,PetscRandom);
 94:   PetscErrorCode (*getrowij)(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt*,const PetscInt *[],const PetscInt *[],PetscBool  *);
 95:   PetscErrorCode (*restorerowij)(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt *,const PetscInt *[],const PetscInt *[],PetscBool  *);
 96:   PetscErrorCode (*getcolumnij)(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt*,const PetscInt *[],const PetscInt *[],PetscBool  *);
 97:   PetscErrorCode (*restorecolumnij)(Mat,PetscInt,PetscBool ,PetscBool ,PetscInt*,const PetscInt *[],const PetscInt *[],PetscBool  *);
 98:   /*54*/
 99:   PetscErrorCode (*fdcoloringcreate)(Mat,ISColoring,MatFDColoring);
100:   PetscErrorCode (*coloringpatch)(Mat,PetscInt,PetscInt,ISColoringValue[],ISColoring*);
101:   PetscErrorCode (*setunfactored)(Mat);
102:   PetscErrorCode (*permute)(Mat,IS,IS,Mat*);
103:   PetscErrorCode (*setvaluesblocked)(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
104:   /*59*/
105:   PetscErrorCode (*createsubmatrix)(Mat,IS,IS,MatReuse,Mat*);
106:   PetscErrorCode (*destroy)(Mat);
107:   PetscErrorCode (*view)(Mat,PetscViewer);
108:   PetscErrorCode (*convertfrom)(Mat, MatType,MatReuse,Mat*);
109:   PetscErrorCode (*matmatmult)(Mat,Mat,Mat,MatReuse,PetscReal,Mat*);
110:   /*64*/
111:   PetscErrorCode (*matmatmultsymbolic)(Mat,Mat,Mat,PetscReal,Mat*);
112:   PetscErrorCode (*matmatmultnumeric)(Mat,Mat,Mat,Mat);
113:   PetscErrorCode (*setlocaltoglobalmapping)(Mat,ISLocalToGlobalMapping,ISLocalToGlobalMapping);
114:   PetscErrorCode (*setvalueslocal)(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
115:   PetscErrorCode (*zerorowslocal)(Mat,PetscInt,const PetscInt[],PetscScalar,Vec,Vec);
116:   /*69*/
117:   PetscErrorCode (*getrowmaxabs)(Mat,Vec,PetscInt[]);
118:   PetscErrorCode (*getrowminabs)(Mat,Vec,PetscInt[]);
119:   PetscErrorCode (*convert)(Mat, MatType,MatReuse,Mat*);
120:   PetscErrorCode (*hasoperation)(Mat,MatOperation,PetscBool*);
121:   PetscErrorCode (*placeholder_73)(Mat,void*);
122:   /*74*/
123:   PetscErrorCode (*setvaluesadifor)(Mat,PetscInt,void*);
124:   PetscErrorCode (*fdcoloringapply)(Mat,MatFDColoring,Vec,void*);
125:   PetscErrorCode (*setfromoptions)(PetscOptionItems*,Mat);
126:   PetscErrorCode (*multconstrained)(Mat,Vec,Vec);
127:   PetscErrorCode (*multtransposeconstrained)(Mat,Vec,Vec);
128:   /*79*/
129:   PetscErrorCode (*findzerodiagonals)(Mat,IS*);
130:   PetscErrorCode (*mults)(Mat, Vecs, Vecs);
131:   PetscErrorCode (*solves)(Mat, Vecs, Vecs);
132:   PetscErrorCode (*getinertia)(Mat,PetscInt*,PetscInt*,PetscInt*);
133:   PetscErrorCode (*load)(Mat, PetscViewer);
134:   /*84*/
135:   PetscErrorCode (*issymmetric)(Mat,PetscReal,PetscBool *);
136:   PetscErrorCode (*ishermitian)(Mat,PetscReal,PetscBool *);
137:   PetscErrorCode (*isstructurallysymmetric)(Mat,PetscBool *);
138:   PetscErrorCode (*setvaluesblockedlocal)(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode);
139:   PetscErrorCode (*getvecs)(Mat,Vec*,Vec*);
140:   /*89*/
141:   PetscErrorCode (*matmult)(Mat,Mat,MatReuse,PetscReal,Mat*);
142:   PetscErrorCode (*matmultsymbolic)(Mat,Mat,PetscReal,Mat*);
143:   PetscErrorCode (*matmultnumeric)(Mat,Mat,Mat);
144:   PetscErrorCode (*ptap)(Mat,Mat,MatReuse,PetscReal,Mat*);
145:   PetscErrorCode (*ptapsymbolic)(Mat,Mat,PetscReal,Mat*); /* double dispatch wrapper routine */
146:   /*94*/
147:   PetscErrorCode (*ptapnumeric)(Mat,Mat,Mat);             /* double dispatch wrapper routine */
148:   PetscErrorCode (*mattransposemult)(Mat,Mat,MatReuse,PetscReal,Mat*);
149:   PetscErrorCode (*mattransposemultsymbolic)(Mat,Mat,PetscReal,Mat*);
150:   PetscErrorCode (*mattransposemultnumeric)(Mat,Mat,Mat);
151:   PetscErrorCode (*placeholder_98)(Mat);
152:   /*99*/
153:   PetscErrorCode (*placeholder_99)(Mat);
154:   PetscErrorCode (*placeholder_100)(Mat);
155:   PetscErrorCode (*placeholder_101)(Mat);
156:   PetscErrorCode (*conjugate)(Mat);                              /* complex conjugate */
157:   PetscErrorCode (*viewnative)(Mat,PetscViewer);
158:   /*104*/
159:   PetscErrorCode (*setvaluesrow)(Mat,PetscInt,const PetscScalar[]);
160:   PetscErrorCode (*realpart)(Mat);
161:   PetscErrorCode (*imaginarypart)(Mat);
162:   PetscErrorCode (*getrowuppertriangular)(Mat);
163:   PetscErrorCode (*restorerowuppertriangular)(Mat);
164:   /*109*/
165:   PetscErrorCode (*matsolve)(Mat,Mat,Mat);
166:   PetscErrorCode (*matsolvetranspose)(Mat,Mat,Mat);
167:   PetscErrorCode (*getrowmin)(Mat,Vec,PetscInt[]);
168:   PetscErrorCode (*getcolumnvector)(Mat,Vec,PetscInt);
169:   PetscErrorCode (*missingdiagonal)(Mat,PetscBool *,PetscInt*);
170:   /*114*/
171:   PetscErrorCode (*getseqnonzerostructure)(Mat,Mat *);
172:   PetscErrorCode (*create)(Mat);
173:   PetscErrorCode (*getghosts)(Mat,PetscInt*,const PetscInt *[]);
174:   PetscErrorCode (*getlocalsubmatrix)(Mat,IS,IS,Mat*);
175:   PetscErrorCode (*restorelocalsubmatrix)(Mat,IS,IS,Mat*);
176:   /*119*/
177:   PetscErrorCode (*multdiagonalblock)(Mat,Vec,Vec);
178:   PetscErrorCode (*hermitiantranspose)(Mat,MatReuse,Mat*);
179:   PetscErrorCode (*multhermitiantranspose)(Mat,Vec,Vec);
180:   PetscErrorCode (*multhermitiantransposeadd)(Mat,Vec,Vec,Vec);
181:   PetscErrorCode (*getmultiprocblock)(Mat,MPI_Comm,MatReuse,Mat*);
182:   /*124*/
183:   PetscErrorCode (*findnonzerorows)(Mat,IS*);
184:   PetscErrorCode (*getcolumnnorms)(Mat,NormType,PetscReal*);
185:   PetscErrorCode (*invertblockdiagonal)(Mat,const PetscScalar**);
186:   PetscErrorCode (*placeholder_127)(Mat,Vec,Vec,Vec);
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: };
210: /*
211:     If you add MatOps entries above also add them to the MATOP enum
212:     in include/petscmat.h and include/petsc/finclude/petscmat.h
213: */

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

219: typedef struct _p_MatBaseName* MatBaseName;
220: struct _p_MatBaseName {
221:   char        *bname,*sname,*mname;
222:   MatBaseName next;
223: };

225: PETSC_EXTERN MatBaseName MatBaseNameList;

227: /*
228:    Utility private matrix routines
229: */
230: PETSC_INTERN PetscErrorCode MatConvert_Basic(Mat, MatType,MatReuse,Mat*);
231: PETSC_INTERN PetscErrorCode MatCopy_Basic(Mat,Mat,MatStructure);
232: PETSC_INTERN PetscErrorCode MatDiagonalSet_Default(Mat,Vec,InsertMode);

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

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

248: typedef struct _MatStashSpace *PetscMatStashSpace;

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

259: PETSC_EXTERN PetscErrorCode PetscMatStashSpaceGet(PetscInt,PetscInt,PetscMatStashSpace *);
260: PETSC_EXTERN PetscErrorCode PetscMatStashSpaceContiguous(PetscInt,PetscMatStashSpace *,PetscScalar *,PetscInt *,PetscInt *);
261: PETSC_EXTERN PetscErrorCode PetscMatStashSpaceDestroy(PetscMatStashSpace*);

263: typedef struct {
264:   PetscInt    count;
265: } MatStashHeader;

267: typedef struct {
268:   void        *buffer;          /* Of type blocktype, dynamically constructed  */
269:   PetscInt    count;
270:   char        pending;
271: } MatStashFrame;

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

283:   PetscErrorCode (*ScatterBegin)(Mat,MatStash*,PetscInt*);
284:   PetscErrorCode (*ScatterGetMesg)(MatStash*,PetscMPIInt*,PetscInt**,PetscInt**,PetscScalar**,PetscInt*);
285:   PetscErrorCode (*ScatterEnd)(MatStash*);
286:   PetscErrorCode (*ScatterDestroy)(MatStash*);

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

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

333: PETSC_INTERN PetscErrorCode MatStashCreate_Private(MPI_Comm,PetscInt,MatStash*);
334: PETSC_INTERN PetscErrorCode MatStashDestroy_Private(MatStash*);
335: PETSC_INTERN PetscErrorCode MatStashScatterEnd_Private(MatStash*);
336: PETSC_INTERN PetscErrorCode MatStashSetInitialSize_Private(MatStash*,PetscInt);
337: PETSC_INTERN PetscErrorCode MatStashGetInfo_Private(MatStash*,PetscInt*,PetscInt*);
338: PETSC_INTERN PetscErrorCode MatStashValuesRow_Private(MatStash*,PetscInt,PetscInt,const PetscInt[],const PetscScalar[],PetscBool );
339: PETSC_INTERN PetscErrorCode MatStashValuesCol_Private(MatStash*,PetscInt,PetscInt,const PetscInt[],const PetscScalar[],PetscInt,PetscBool );
340: PETSC_INTERN PetscErrorCode MatStashValuesRowBlocked_Private(MatStash*,PetscInt,PetscInt,const PetscInt[],const PetscScalar[],PetscInt,PetscInt,PetscInt);
341: PETSC_INTERN PetscErrorCode MatStashValuesColBlocked_Private(MatStash*,PetscInt,PetscInt,const PetscInt[],const PetscScalar[],PetscInt,PetscInt,PetscInt);
342: PETSC_INTERN PetscErrorCode MatStashScatterBegin_Private(Mat,MatStash*,PetscInt*);
343: PETSC_INTERN PetscErrorCode MatStashScatterGetMesg_Private(MatStash*,PetscMPIInt*,PetscInt**,PetscInt**,PetscScalar**,PetscInt*);
344: PETSC_INTERN PetscErrorCode MatGetInfo_External(Mat,MatInfoType,MatInfo*);

346: typedef struct {
347:   PetscInt   dim;
348:   PetscInt   dims[4];
349:   PetscInt   starts[4];
350:   PetscBool  noc;        /* this is a single component problem, hence user will not set MatStencil.c */
351: } MatStencilInfo;

353: /* Info about using compressed row format */
354: typedef struct {
355:   PetscBool  use;                           /* indicates compressed rows have been checked and will be used */
356:   PetscInt   nrows;                         /* number of non-zero rows */
357:   PetscInt   *i;                            /* compressed row pointer  */
358:   PetscInt   *rindex;                       /* compressed row index               */
359: } Mat_CompressedRow;
360: PETSC_EXTERN PetscErrorCode MatCheckCompressedRow(Mat,PetscInt,Mat_CompressedRow*,PetscInt*,PetscInt,PetscReal);

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

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

413: PETSC_INTERN PetscErrorCode MatAXPY_Basic(Mat,PetscScalar,Mat,MatStructure);
414: PETSC_INTERN PetscErrorCode MatAXPY_BasicWithPreallocation(Mat,Mat,PetscScalar,Mat,MatStructure);

416: /*
417:     Utility for MatFactor (Schur complement)
418: */
419: PETSC_INTERN PetscErrorCode MatFactorFactorizeSchurComplement_Private(Mat);
420: PETSC_INTERN PetscErrorCode MatFactorInvertSchurComplement_Private(Mat);
421: PETSC_INTERN PetscErrorCode MatFactorUpdateSchurStatus_Private(Mat);
422: PETSC_INTERN PetscErrorCode MatFactorSetUpInPlaceSchur_Private(Mat);

424: /*
425:     Utility for MatZeroRows
426: */
427: PETSC_INTERN PetscErrorCode MatZeroRowsMapLocal_Private(Mat,PetscInt,const PetscInt*,PetscInt*,PetscInt**);

429: /*
430:     Object for partitioning graphs
431: */

433: typedef struct _MatPartitioningOps *MatPartitioningOps;
434: struct _MatPartitioningOps {
435:   PetscErrorCode (*apply)(MatPartitioning,IS*);
436:   PetscErrorCode (*setfromoptions)(PetscOptionItems*,MatPartitioning);
437:   PetscErrorCode (*destroy)(MatPartitioning);
438:   PetscErrorCode (*view)(MatPartitioning,PetscViewer);
439: };

441: struct _p_MatPartitioning {
442:   PETSCHEADER(struct _MatPartitioningOps);
443:   Mat         adj;
444:   PetscInt    *vertex_weights;
445:   PetscReal   *part_weights;
446:   PetscInt    n;                                 /* number of partitions */
447:   void        *data;
448:   PetscInt    setupcalled;
449: };

451: /*
452:     Object for coarsen graphs
453: */
454: typedef struct _MatCoarsenOps *MatCoarsenOps;
455: struct _MatCoarsenOps {
456:   PetscErrorCode (*apply)(MatCoarsen);
457:   PetscErrorCode (*setfromoptions)(PetscOptionItems*,MatCoarsen);
458:   PetscErrorCode (*destroy)(MatCoarsen);
459:   PetscErrorCode (*view)(MatCoarsen,PetscViewer);
460: };

462: struct _p_MatCoarsen {
463:   PETSCHEADER(struct _MatCoarsenOps);
464:   Mat              graph;
465:   PetscInt         setupcalled;
466:   void             *subctx;
467:   /* */
468:   PetscBool        strict_aggs;
469:   IS               perm;
470:   PetscCoarsenData *agg_lists;
471: };

473: /*
474:     MatFDColoring is used to compute Jacobian matrices efficiently
475:   via coloring. The data structure is explained below in an example.

477:    Color =   0    1     0    2   |   2      3       0
478:    ---------------------------------------------------
479:             00   01              |          05
480:             10   11              |   14     15               Processor  0
481:                        22    23  |          25
482:                        32    33  |
483:    ===================================================
484:                                  |   44     45     46
485:             50                   |          55               Processor 1
486:                                  |   64            66
487:    ---------------------------------------------------

489:     ncolors = 4;

491:     ncolumns      = {2,1,1,0}
492:     columns       = {{0,2},{1},{3},{}}
493:     nrows         = {4,2,3,3}
494:     rows          = {{0,1,2,3},{0,1},{1,2,3},{0,1,2}}
495:     vwscale       = {dx(0),dx(1),dx(2),dx(3)}               MPI Vec
496:     vscale        = {dx(0),dx(1),dx(2),dx(3),dx(4),dx(5)}   Seq Vec

498:     ncolumns      = {1,0,1,1}
499:     columns       = {{6},{},{4},{5}}
500:     nrows         = {3,0,2,2}
501:     rows          = {{0,1,2},{},{1,2},{1,2}}
502:     vwscale       = {dx(4),dx(5),dx(6)}              MPI Vec
503:     vscale        = {dx(0),dx(4),dx(5),dx(6)}        Seq Vec

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

508: */
509: typedef struct {
510:   PetscInt     row;
511:   PetscInt     col;
512:   PetscScalar  *valaddr;   /* address of value */
513: } MatEntry;

515: typedef struct {
516:   PetscInt     row;
517:   PetscScalar  *valaddr;   /* address of value */
518: } MatEntry2;

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

547: typedef struct _MatColoringOps *MatColoringOps;
548: struct _MatColoringOps {
549:   PetscErrorCode (*destroy)(MatColoring);
550:   PetscErrorCode (*setfromoptions)(PetscOptionItems*,MatColoring);
551:   PetscErrorCode (*view)(MatColoring,PetscViewer);
552:   PetscErrorCode (*apply)(MatColoring,ISColoring*);
553:   PetscErrorCode (*weights)(MatColoring,PetscReal**,PetscInt**);
554: };

556: struct _p_MatColoring {
557:   PETSCHEADER(struct _MatColoringOps);
558:   Mat                   mat;
559:   PetscInt              dist;             /* distance of the coloring */
560:   PetscInt              maxcolors;        /* the maximum number of colors returned, maxcolors=1 for MIS */
561:   void                  *data;            /* inner context */
562:   PetscBool             valid;            /* check to see if what is produced is a valid coloring */
563:   MatColoringWeightType weight_type;      /* type of weight computation to be performed */
564:   PetscReal             *user_weights;    /* custom weights and permutation */
565:   PetscInt              *user_lperm;
566:   PetscBool             valid_iscoloring; /* check to see if matcoloring is produced a valid iscoloring */
567: };

569: struct  _p_MatTransposeColoring{
570:   PETSCHEADER(int);
571:   PetscInt       M,N,m;            /* total rows, columns; local rows */
572:   PetscInt       rstart;           /* first row owned by local processor */
573:   PetscInt       ncolors;          /* number of colors */
574:   PetscInt       *ncolumns;        /* number of local columns for a color */
575:   PetscInt       *nrows;           /* number of local rows for each color */
576:   PetscInt       currentcolor;     /* color for which function evaluation is being done now */
577:   ISColoringType ctype;            /* IS_COLORING_GLOBAL or IS_COLORING_LOCAL */

579:   PetscInt       *colorforrow,*colorforcol;  /* pointer to rows and columns */
580:   PetscInt       *rows;                      /* lists the local rows for each color (using the local row numbering) */
581:   PetscInt       *den2sp;                    /* maps (row,color) in the dense matrix to index of sparse matrix array a->a */
582:   PetscInt       *columns;                   /* lists the local columns of each color (using global column numbering) */
583:   PetscInt       brows;                      /* number of rows for efficient implementation of MatTransColoringApplyDenToSp() */
584:   PetscInt       *lstart;                    /* array used for loop over row blocks of Csparse */
585: };

587: /*
588:    Null space context for preconditioner/operators
589: */
590: struct _p_MatNullSpace {
591:   PETSCHEADER(int);
592:   PetscBool      has_cnst;
593:   PetscInt       n;
594:   Vec*           vecs;
595:   PetscScalar*   alpha;                 /* for projections */
596:   PetscErrorCode (*remove)(MatNullSpace,Vec,void*);  /* for user provided removal function */
597:   void*          rmctx;                 /* context for remove() function */
598: };

600: /*
601:    Checking zero pivot for LU, ILU preconditioners.
602: */
603: typedef struct {
604:   PetscInt       nshift,nshift_max;
605:   PetscReal      shift_amount,shift_lo,shift_hi,shift_top,shift_fraction;
606:   PetscBool      newshift;
607:   PetscReal      rs;  /* active row sum of abs(offdiagonals) */
608:   PetscScalar    pv;  /* pivot of the active row */
609: } FactorShiftCtx;

611: /*
612:  Used by MatCreateSubMatrices_MPIXAIJ_Local()
613: */
614:  #include <petscctable.h>
615: typedef struct { /* used by MatCreateSubMatrices_MPIAIJ_SingleIS_Local() and MatCreateSubMatrices_MPIAIJ_Local */
616:   PetscInt   id;   /* index of submats, only submats[0] is responsible for deleting some arrays below */
617:   PetscInt   nrqs,nrqr;
618:   PetscInt   **rbuf1,**rbuf2,**rbuf3,**sbuf1,**sbuf2;
619:   PetscInt   **ptr;
620:   PetscInt   *tmp;
621:   PetscInt   *ctr;
622:   PetscInt   *pa; /* proc array */
623:   PetscInt   *req_size,*req_source1,*req_source2;
624:   PetscBool  allcolumns,allrows;
625:   PetscBool  singleis;
626:   PetscInt   *row2proc; /* row to proc map */
627:   PetscInt   nstages;
628: #if defined(PETSC_USE_CTABLE)
629:   PetscTable cmap,rmap;
630:   PetscInt   *cmap_loc,*rmap_loc;
631: #else
632:   PetscInt   *cmap,*rmap;
633: #endif

635:   PetscErrorCode (*destroy)(Mat);
636: } Mat_SubSppt;

638: PETSC_EXTERN PetscErrorCode MatFactorDumpMatrix(Mat);
639: PETSC_INTERN PetscErrorCode MatShift_Basic(Mat,PetscScalar);
640: PETSC_INTERN PetscErrorCode MatSetBlockSizes_Default(Mat,PetscInt,PetscInt);

642: PETSC_STATIC_INLINE PetscErrorCode MatPivotCheck_nz(Mat mat,const MatFactorInfo *info,FactorShiftCtx *sctx,PetscInt row)
643: {
644:   PetscReal _rs   = sctx->rs;
645:   PetscReal _zero = info->zeropivot*_rs;

648:   if (PetscAbsScalar(sctx->pv) <= _zero && !PetscIsNanScalar(sctx->pv)){
649:     /* force |diag| > zeropivot*rs */
650:     if (!sctx->nshift) sctx->shift_amount = info->shiftamount;
651:     else sctx->shift_amount *= 2.0;
652:     sctx->newshift = PETSC_TRUE;
653:     (sctx->nshift)++;
654:   } else {
655:     sctx->newshift = PETSC_FALSE;
656:   }
657:   return(0);
658: }

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

666:   if (PetscRealPart(sctx->pv) <= _zero && !PetscIsNanScalar(sctx->pv)){
667:     /* force matfactor to be diagonally dominant */
668:     if (sctx->nshift == sctx->nshift_max) {
669:       sctx->shift_fraction = sctx->shift_hi;
670:     } else {
671:       sctx->shift_lo = sctx->shift_fraction;
672:       sctx->shift_fraction = (sctx->shift_hi+sctx->shift_lo)/2.;
673:     }
674:     sctx->shift_amount = sctx->shift_fraction * sctx->shift_top;
675:     sctx->nshift++;
676:     sctx->newshift = PETSC_TRUE;
677:   } else {
678:     sctx->newshift = PETSC_FALSE;
679:   }
680:   return(0);
681: }

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

688:   if (PetscAbsScalar(sctx->pv) <= _zero && !PetscIsNanScalar(sctx->pv)){
689:     sctx->pv          += info->shiftamount;
690:     sctx->shift_amount = 0.0;
691:     sctx->nshift++;
692:   }
693:   sctx->newshift = PETSC_FALSE;
694:   return(0);
695: }

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

703:   sctx->newshift = PETSC_FALSE;
704:   if (PetscAbsScalar(sctx->pv) <= _zero && !PetscIsNanScalar(sctx->pv)) {
705:     if (!mat->erroriffailure) {
706:       PetscInfo3(mat,"Detected zero pivot in factorization in row %D value %g tolerance %g\n",row,(double)PetscAbsScalar(sctx->pv),(double)_zero);
707:       fact->factorerrortype             = MAT_FACTOR_NUMERIC_ZEROPIVOT;
708:       fact->factorerror_zeropivot_value = PetscAbsScalar(sctx->pv);
709:       fact->factorerror_zeropivot_row   = row;
710:     } 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);
711:   }
712:   return(0);
713: }

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

720:   if (info->shifttype == (PetscReal) MAT_SHIFT_NONZERO){
721:     MatPivotCheck_nz(mat,info,sctx,row);
722:   } else if (info->shifttype == (PetscReal) MAT_SHIFT_POSITIVE_DEFINITE){
723:     MatPivotCheck_pd(mat,info,sctx,row);
724:   } else if (info->shifttype == (PetscReal) MAT_SHIFT_INBLOCKS){
725:     MatPivotCheck_inblocks(mat,info,sctx,row);
726:   } else {
727:     MatPivotCheck_none(fact,mat,info,sctx,row);
728:   }
729:   return(0);
730: }

732: /*
733:   Create and initialize a linked list
734:   Input Parameters:
735:     idx_start - starting index of the list
736:     lnk_max   - max value of lnk indicating the end of the list
737:     nlnk      - max length of the list
738:   Output Parameters:
739:     lnk       - list initialized
740:     bt        - PetscBT (bitarray) with all bits set to false
741:     lnk_empty - flg indicating the list is empty
742: */
743: #define PetscLLCreate(idx_start,lnk_max,nlnk,lnk,bt) \
744:   (PetscMalloc1(nlnk,&lnk) || PetscBTCreate(nlnk,&(bt)) || (lnk[idx_start] = lnk_max,0))

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

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

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

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

859: #define PetscLLAddSorted_new(nidx,indices,idx_start,lnk_empty,nlnk,lnk,bt) 0; \
860: {\
861:   PetscInt _k,_entry,_location,_lnkdata;\
862:   if (lnk_empty){\
863:     _lnkdata  = idx_start;                      \
864:     for (_k=0; _k<nidx; _k++){                  \
865:       _entry = indices[_k];                             \
866:       PetscBTSet(bt,_entry);  /* mark the new entry */          \
867:           _location = _lnkdata;                                 \
868:           _lnkdata  = lnk[_location];                           \
869:         /* insertion location is found, add entry into lnk */   \
870:         lnk[_location] = _entry;                                \
871:         lnk[_entry]    = _lnkdata;                              \
872:         _lnkdata = _entry; /* next search starts from here */   \
873:     }                                                           \
874:     /*\
875:     lnk[indices[nidx-1]] = lnk[idx_start];\
876:     lnk[idx_start]       = indices[0];\
877:     PetscBTSet(bt,indices[0]);  \
878:     for (_k=1; _k<nidx; _k++){                  \
879:       PetscBTSet(bt,indices[_k]);                                          \
880:       lnk[indices[_k-1]] = indices[_k];                                  \
881:     }                                                           \
882:      */\
883:     nlnk      = nidx;\
884:     lnk_empty = PETSC_FALSE;\
885:   } else {\
886:     nlnk      = 0;                              \
887:     _lnkdata  = idx_start;                      \
888:     for (_k=0; _k<nidx; _k++){                  \
889:       _entry = indices[_k];                             \
890:       if (!PetscBTLookupSet(bt,_entry)){  /* new entry */       \
891:         /* search for insertion location */                     \
892:         do {                                                    \
893:           _location = _lnkdata;                                 \
894:           _lnkdata  = lnk[_location];                           \
895:         } while (_entry > _lnkdata);                            \
896:         /* insertion location is found, add entry into lnk */   \
897:         lnk[_location] = _entry;                                \
898:         lnk[_entry]    = _lnkdata;                              \
899:         nlnk++;                                                 \
900:         _lnkdata = _entry; /* next search starts from here */   \
901:       }                                                         \
902:     }                                                           \
903:   }                                                             \
904: }

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

950: /*
951:   Copy data on the list into an array, then initialize the list
952:   Input Parameters:
953:     idx_start - starting index of the list
954:     lnk_max   - max value of lnk indicating the end of the list
955:     nlnk      - number of data on the list to be copied
956:     lnk       - linked list
957:     bt        - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
958:   output Parameters:
959:     indices   - array that contains the copied data
960:     lnk       - linked list that is cleaned and initialize
961:     bt        - PetscBT (bitarray) with all bits set to false
962: */
963: #define PetscLLClean(idx_start,lnk_max,nlnk,lnk,indices,bt) 0;\
964: {\
965:   PetscInt _j,_idx=idx_start;\
966:   for (_j=0; _j<nlnk; _j++){\
967:     _idx = lnk[_idx];\
968:     indices[_j] = _idx;\
969:     PetscBTClear(bt,_idx);\
970:   }\
971:   lnk[idx_start] = lnk_max;\
972: }
973: /*
974:   Free memories used by the list
975: */
976: #define PetscLLDestroy(lnk,bt) (PetscFree(lnk) || PetscBTDestroy(&(bt)))

978: /* Routines below are used for incomplete matrix factorization */
979: /*
980:   Create and initialize a linked list and its levels
981:   Input Parameters:
982:     idx_start - starting index of the list
983:     lnk_max   - max value of lnk indicating the end of the list
984:     nlnk      - max length of the list
985:   Output Parameters:
986:     lnk       - list initialized
987:     lnk_lvl   - array of size nlnk for storing levels of lnk
988:     bt        - PetscBT (bitarray) with all bits set to false
989: */
990: #define PetscIncompleteLLCreate(idx_start,lnk_max,nlnk,lnk,lnk_lvl,bt)\
991:   (PetscIntMultError(2,nlnk,NULL) || PetscMalloc1(2*nlnk,&lnk) || PetscBTCreate(nlnk,&(bt)) || (lnk[idx_start] = lnk_max,lnk_lvl = lnk + nlnk,0))

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

1033: /*
1034:   Add a SORTED index set into a sorted linked list for ILU
1035:   Input Parameters:
1036:     nidx      - number of input indices
1037:     idx       - sorted integer array used for storing column indices
1038:     level     - level of fill, e.g., ICC(level)
1039:     idxlvl    - level of idx
1040:     idx_start - starting index of the list
1041:     lnk       - linked list(an integer array) that is created
1042:     lnklvl    - levels of lnk
1043:     bt        - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
1044:     prow      - the row number of idx
1045:   output Parameters:
1046:     nlnk     - number of newly added idx
1047:     lnk      - the sorted(increasing order) linked list containing new and non-redundate entries from idx
1048:     lnklvl   - levels of lnk
1049:     bt       - updated PetscBT (bitarray)

1051:   Note: the level of factor(i,j) is set as lvl(i,j) = min{ lvl(i,j), lvl(i,prow)+lvl(prow,j)+1)
1052:         where idx = non-zero columns of U(prow,prow+1:n-1), prow<i
1053: */
1054: #define PetscILULLAddSorted(nidx,idx,level,idxlvl,idx_start,nlnk,lnk,lnklvl,bt,lnklvl_prow) 0;\
1055: {\
1056:   PetscInt _k,_entry,_location,_lnkdata,_incrlev,_lnklvl_prow=lnklvl[prow];\
1057:   nlnk     = 0;\
1058:   _lnkdata = idx_start;\
1059:   for (_k=0; _k<nidx; _k++){\
1060:     _incrlev = idxlvl[_k] + _lnklvl_prow + 1;\
1061:     if (_incrlev > level) continue;\
1062:     _entry = idx[_k];\
1063:     if (!PetscBTLookupSet(bt,_entry)){  /* new entry */\
1064:       /* search for insertion location */\
1065:       do {\
1066:         _location = _lnkdata;\
1067:         _lnkdata  = lnk[_location];\
1068:       } while (_entry > _lnkdata);\
1069:       /* insertion location is found, add entry into lnk */\
1070:       lnk[_location]  = _entry;\
1071:       lnk[_entry]     = _lnkdata;\
1072:       lnklvl[_entry] = _incrlev;\
1073:       nlnk++;\
1074:       _lnkdata = _entry; /* next search starts from here if next_entry > _entry */\
1075:     } else { /* existing entry: update lnklvl */\
1076:       if (lnklvl[_entry] > _incrlev) lnklvl[_entry] = _incrlev;\
1077:     }\
1078:   }\
1079: }

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

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

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

1217: /*
1218:   Copy data on the list into an array, then initialize the list
1219:   Input Parameters:
1220:     idx_start - starting index of the list
1221:     lnk_max   - max value of lnk indicating the end of the list
1222:     nlnk      - number of data on the list to be copied
1223:     lnk       - linked list
1224:     lnklvl    - level of lnk
1225:     bt        - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
1226:   output Parameters:
1227:     indices - array that contains the copied data
1228:     lnk     - linked list that is cleaned and initialize
1229:     lnklvl  - level of lnk that is reinitialized
1230:     bt      - PetscBT (bitarray) with all bits set to false
1231: */
1232: #define PetscIncompleteLLClean(idx_start,lnk_max,nlnk,lnk,lnklvl,indices,indiceslvl,bt) 0;\
1233: {\
1234:   PetscInt _j,_idx=idx_start;\
1235:   for (_j=0; _j<nlnk; _j++){\
1236:     _idx = lnk[_idx];\
1237:     *(indices+_j) = _idx;\
1238:     *(indiceslvl+_j) = lnklvl[_idx];\
1239:     lnklvl[_idx] = -1;\
1240:     PetscBTClear(bt,_idx);\
1241:   }\
1242:   lnk[idx_start] = lnk_max;\
1243: }
1244: /*
1245:   Free memories used by the list
1246: */
1247: #define PetscIncompleteLLDestroy(lnk,bt) (PetscFree(lnk) || PetscBTDestroy(&(bt)))

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

1258:       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
1259:       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
1260:       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
1261:       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.
1262:       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
1263:       to each other so memory access is much better than using the big array.

1265:   Example:
1266:      nlnk_max=5, lnk_max=36:
1267:      Initial list: [0, 0 | 36, 2 | 0, 0 | 0, 0 | 0, 0 | 0, 0 | 0, 0]
1268:      here, head_node has index 2 with value lnk[2]=lnk_max=36,
1269:            0-th entry is used to store the number of entries in the list,
1270:      The initial lnk represents head -> tail(marked by 36) with number of entries = lnk[0]=0.

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

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

1280:   Input Parameters:
1281:     nlnk_max  - max length of the list
1282:     lnk_max   - max value of the entries
1283:   Output Parameters:
1284:     lnk       - list created and initialized
1285:     bt        - PetscBT (bitarray) with all bits set to false. Note: bt has size lnk_max, not nln_max!
1286: */
1287: PETSC_STATIC_INLINE PetscErrorCode PetscLLCondensedCreate(PetscInt nlnk_max,PetscInt lnk_max,PetscInt **lnk,PetscBT *bt)
1288: {
1290:   PetscInt       *llnk,lsize = 0;

1293:   PetscIntMultError(2,nlnk_max+2,&lsize);
1294:   PetscMalloc1(lsize,lnk);
1295:   PetscBTCreate(lnk_max,bt);
1296:   llnk = *lnk;
1297:   llnk[0] = 0;         /* number of entries on the list */
1298:   llnk[2] = lnk_max;   /* value in the head node */
1299:   llnk[3] = 2;         /* next for the head node */
1300:   return(0);
1301: }

1303: /*
1304:   Add a SORTED ascending index set into a sorted linked list. See PetscLLCondensedCreate() for detailed description.
1305:   Input Parameters:
1306:     nidx      - number of input indices
1307:     indices   - sorted integer array
1308:     lnk       - condensed linked list(an integer array) that is created
1309:     bt        - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
1310:   output Parameters:
1311:     lnk       - the sorted(increasing order) linked list containing previous and newly added non-redundate indices
1312:     bt        - updated PetscBT (bitarray)
1313: */
1314: PETSC_STATIC_INLINE PetscErrorCode PetscLLCondensedAddSorted(PetscInt nidx,const PetscInt indices[],PetscInt lnk[],PetscBT bt)
1315: {
1316:   PetscInt _k,_entry,_location,_next,_lnkdata,_nlnk,_newnode;

1319:   _nlnk     = lnk[0]; /* num of entries on the input lnk */
1320:   _location = 2; /* head */
1321:     for (_k=0; _k<nidx; _k++){
1322:       _entry = indices[_k];
1323:       if (!PetscBTLookupSet(bt,_entry)){  /* new entry */
1324:         /* search for insertion location */
1325:         do {
1326:           _next     = _location + 1; /* link from previous node to next node */
1327:           _location = lnk[_next];    /* idx of next node */
1328:           _lnkdata  = lnk[_location];/* value of next node */
1329:         } while (_entry > _lnkdata);
1330:         /* insertion location is found, add entry into lnk */
1331:         _newnode        = 2*(_nlnk+2);   /* index for this new node */
1332:         lnk[_next]      = _newnode;      /* connect previous node to the new node */
1333:         lnk[_newnode]   = _entry;        /* set value of the new node */
1334:         lnk[_newnode+1] = _location;     /* connect new node to next node */
1335:         _location       = _newnode;      /* next search starts from the new node */
1336:         _nlnk++;
1337:       }   \
1338:     }\
1339:   lnk[0]   = _nlnk;   /* number of entries in the list */
1340:   return(0);
1341: }

1343: PETSC_STATIC_INLINE PetscErrorCode PetscLLCondensedClean(PetscInt lnk_max,PetscInt nidx,PetscInt *indices,PetscInt lnk[],PetscBT bt)
1344: {
1346:   PetscInt       _k,_next,_nlnk;

1349:   _next = lnk[3];       /* head node */
1350:   _nlnk = lnk[0];       /* num of entries on the list */
1351:   for (_k=0; _k<_nlnk; _k++){
1352:     indices[_k] = lnk[_next];
1353:     _next       = lnk[_next + 1];
1354:     PetscBTClear(bt,indices[_k]);
1355:   }
1356:   lnk[0] = 0;          /* num of entries on the list */
1357:   lnk[2] = lnk_max;    /* initialize head node */
1358:   lnk[3] = 2;          /* head node */
1359:   return(0);
1360: }

1362: PETSC_STATIC_INLINE PetscErrorCode PetscLLCondensedView(PetscInt *lnk)
1363: {
1365:   PetscInt       k;

1368:   PetscPrintf(PETSC_COMM_SELF,"LLCondensed of size %D, (val,  next)\n",lnk[0]);
1369:   for (k=2; k< lnk[0]+2; k++){
1370:     PetscPrintf(PETSC_COMM_SELF," %D: (%D, %D)\n",2*k,lnk[2*k],lnk[2*k+1]);
1371:   }
1372:   return(0);
1373: }

1375: /*
1376:   Free memories used by the list
1377: */
1378: PETSC_STATIC_INLINE PetscErrorCode PetscLLCondensedDestroy(PetscInt *lnk,PetscBT bt)
1379: {

1383:   PetscFree(lnk);
1384:   PetscBTDestroy(&bt);
1385:   return(0);
1386: }

1388: /* -------------------------------------------------------------------------------------------------------*/
1389: /*
1390:  Same as PetscLLCondensedCreate(), but does not use non-scalable O(lnk_max) bitarray
1391:   Input Parameters:
1392:     nlnk_max  - max length of the list
1393:   Output Parameters:
1394:     lnk       - list created and initialized
1395: */
1396: PETSC_STATIC_INLINE PetscErrorCode PetscLLCondensedCreate_Scalable(PetscInt nlnk_max,PetscInt **lnk)
1397: {
1399:   PetscInt       *llnk,lsize = 0;

1402:   PetscIntMultError(2,nlnk_max+2,&lsize);
1403:   PetscMalloc1(lsize,lnk);
1404:   llnk = *lnk;
1405:   llnk[0] = 0;               /* number of entries on the list */
1406:   llnk[2] = PETSC_MAX_INT;   /* value in the head node */
1407:   llnk[3] = 2;               /* next for the head node */
1408:   return(0);
1409: }

1411: PETSC_STATIC_INLINE PetscErrorCode PetscLLCondensedAddSorted_Scalable(PetscInt nidx,const PetscInt indices[],PetscInt lnk[])
1412: {
1413:   PetscInt _k,_entry,_location,_next,_lnkdata,_nlnk,_newnode;
1414:   _nlnk     = lnk[0]; /* num of entries on the input lnk */
1415:   _location = 2; /* head */ \
1416:     for (_k=0; _k<nidx; _k++){
1417:       _entry = indices[_k];
1418:       /* search for insertion location */
1419:       do {
1420:         _next     = _location + 1; /* link from previous node to next node */
1421:         _location = lnk[_next];    /* idx of next node */
1422:         _lnkdata  = lnk[_location];/* value of next node */
1423:       } while (_entry > _lnkdata);
1424:       if (_entry < _lnkdata) {
1425:         /* insertion location is found, add entry into lnk */
1426:         _newnode        = 2*(_nlnk+2);   /* index for this new node */
1427:         lnk[_next]      = _newnode;      /* connect previous node to the new node */
1428:         lnk[_newnode]   = _entry;        /* set value of the new node */
1429:         lnk[_newnode+1] = _location;     /* connect new node to next node */
1430:         _location       = _newnode;      /* next search starts from the new node */
1431:         _nlnk++;
1432:       }
1433:     }
1434:   lnk[0]   = _nlnk;   /* number of entries in the list */
1435:   return 0;
1436: }

1438: PETSC_STATIC_INLINE PetscErrorCode PetscLLCondensedClean_Scalable(PetscInt nidx,PetscInt *indices,PetscInt *lnk)
1439: {
1440:   PetscInt _k,_next,_nlnk;
1441:   _next = lnk[3];       /* head node */
1442:   _nlnk = lnk[0];
1443:   for (_k=0; _k<_nlnk; _k++){
1444:     indices[_k] = lnk[_next];
1445:     _next       = lnk[_next + 1];
1446:   }
1447:   lnk[0] = 0;          /* num of entries on the list */
1448:   lnk[3] = 2;          /* head node */
1449:   return 0;
1450: }

1452: PETSC_STATIC_INLINE PetscErrorCode PetscLLCondensedDestroy_Scalable(PetscInt *lnk)
1453: {
1454:   return PetscFree(lnk);
1455: }

1457: /* -------------------------------------------------------------------------------------------------------*/
1458: /*
1459:       lnk[0]   number of links
1460:       lnk[1]   number of entries
1461:       lnk[3n]  value
1462:       lnk[3n+1] len
1463:       lnk[3n+2] link to next value

1465:       The next three are always the first link

1467:       lnk[3]    PETSC_MIN_INT+1
1468:       lnk[4]    1
1469:       lnk[5]    link to first real entry

1471:       The next three are always the last link

1473:       lnk[6]    PETSC_MAX_INT - 1
1474:       lnk[7]    1
1475:       lnk[8]    next valid link (this is the same as lnk[0] but without the decreases)
1476: */

1478: PETSC_STATIC_INLINE PetscErrorCode PetscLLCondensedCreate_fast(PetscInt nlnk_max,PetscInt **lnk)
1479: {
1481:   PetscInt       *llnk,lsize = 0;

1484:   PetscIntMultError(3,nlnk_max+3,&lsize);
1485:   PetscMalloc1(lsize,lnk);
1486:   llnk = *lnk;
1487:   llnk[0] = 0;   /* nlnk: number of entries on the list */
1488:   llnk[1] = 0;          /* number of integer entries represented in list */
1489:   llnk[3] = PETSC_MIN_INT+1;   /* value in the first node */
1490:   llnk[4] = 1;           /* count for the first node */
1491:   llnk[5] = 6;         /* next for the first node */
1492:   llnk[6] = PETSC_MAX_INT-1;   /* value in the last node */
1493:   llnk[7] = 1;           /* count for the last node */
1494:   llnk[8] = 0;         /* next valid node to be used */
1495:   return(0);
1496: }

1498: PETSC_STATIC_INLINE PetscErrorCode PetscLLCondensedAddSorted_fast(PetscInt nidx,const PetscInt indices[],PetscInt lnk[])
1499: {
1500:   PetscInt k,entry,prev,next;
1501:   prev      = 3;      /* first value */
1502:   next      = lnk[prev+2];
1503:   for (k=0; k<nidx; k++){
1504:     entry = indices[k];
1505:     /* search for insertion location */
1506:     while (entry >= lnk[next]) {
1507:       prev = next;
1508:       next = lnk[next+2];
1509:     }
1510:     /* entry is in range of previous list */
1511:     if (entry < lnk[prev]+lnk[prev+1]) continue;
1512:     lnk[1]++;
1513:     /* entry is right after previous list */
1514:     if (entry == lnk[prev]+lnk[prev+1]) {
1515:       lnk[prev+1]++;
1516:       if (lnk[next] == entry+1) { /* combine two contiguous strings */
1517:         lnk[prev+1] += lnk[next+1];
1518:         lnk[prev+2]  = lnk[next+2];
1519:         next         = lnk[next+2];
1520:         lnk[0]--;
1521:       }
1522:       continue;
1523:     }
1524:     /* entry is right before next list */
1525:     if (entry == lnk[next]-1) {
1526:       lnk[next]--;
1527:       lnk[next+1]++;
1528:       prev = next;
1529:       next = lnk[prev+2];
1530:       continue;
1531:     }
1532:     /*  add entry into lnk */
1533:     lnk[prev+2]    = 3*((lnk[8]++)+3);      /* connect previous node to the new node */
1534:     prev           = lnk[prev+2];
1535:     lnk[prev]      = entry;        /* set value of the new node */
1536:     lnk[prev+1]    = 1;             /* number of values in contiguous string is one to start */
1537:     lnk[prev+2]    = next;          /* connect new node to next node */
1538:     lnk[0]++;
1539:   }
1540:   return 0;
1541: }

1543: PETSC_STATIC_INLINE PetscErrorCode PetscLLCondensedClean_fast(PetscInt nidx,PetscInt *indices,PetscInt *lnk)
1544: {
1545:   PetscInt _k,_next,_nlnk,cnt,j;
1546:   _next = lnk[5];       /* first node */
1547:   _nlnk = lnk[0];
1548:   cnt   = 0;
1549:   for (_k=0; _k<_nlnk; _k++){
1550:     for (j=0; j<lnk[_next+1]; j++) {
1551:       indices[cnt++] = lnk[_next] + j;
1552:     }
1553:     _next       = lnk[_next + 2];
1554:   }
1555:   lnk[0] = 0;   /* nlnk: number of links */
1556:   lnk[1] = 0;          /* number of integer entries represented in list */
1557:   lnk[3] = PETSC_MIN_INT+1;   /* value in the first node */
1558:   lnk[4] = 1;           /* count for the first node */
1559:   lnk[5] = 6;         /* next for the first node */
1560:   lnk[6] = PETSC_MAX_INT-1;   /* value in the last node */
1561:   lnk[7] = 1;           /* count for the last node */
1562:   lnk[8] = 0;         /* next valid location to make link */
1563:   return 0;
1564: }

1566: PETSC_STATIC_INLINE PetscErrorCode PetscLLCondensedView_fast(PetscInt *lnk)
1567: {
1568:   PetscInt k,next,nlnk;
1569:   next = lnk[5];       /* first node */
1570:   nlnk = lnk[0];
1571:   for (k=0; k<nlnk; k++){
1572: #if 0                           /* Debugging code */
1573:     printf("%d value %d len %d next %d\n",next,lnk[next],lnk[next+1],lnk[next+2]);
1574: #endif
1575:     next = lnk[next + 2];
1576:   }
1577:   return 0;
1578: }

1580: PETSC_STATIC_INLINE PetscErrorCode PetscLLCondensedDestroy_fast(PetscInt *lnk)
1581: {
1582:   return PetscFree(lnk);
1583: }

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

1588: PETSC_EXTERN PetscLogEvent MAT_Mult;
1589: PETSC_EXTERN PetscLogEvent MAT_MultMatrixFree;
1590: PETSC_EXTERN PetscLogEvent MAT_Mults;
1591: PETSC_EXTERN PetscLogEvent MAT_MultConstrained;
1592: PETSC_EXTERN PetscLogEvent MAT_MultAdd;
1593: PETSC_EXTERN PetscLogEvent MAT_MultTranspose;
1594: PETSC_EXTERN PetscLogEvent MAT_MultTransposeConstrained;
1595: PETSC_EXTERN PetscLogEvent MAT_MultTransposeAdd;
1596: PETSC_EXTERN PetscLogEvent MAT_Solve;
1597: PETSC_EXTERN PetscLogEvent MAT_Solves;
1598: PETSC_EXTERN PetscLogEvent MAT_SolveAdd;
1599: PETSC_EXTERN PetscLogEvent MAT_SolveTranspose;
1600: PETSC_EXTERN PetscLogEvent MAT_SolveTransposeAdd;
1601: PETSC_EXTERN PetscLogEvent MAT_SOR;
1602: PETSC_EXTERN PetscLogEvent MAT_ForwardSolve;
1603: PETSC_EXTERN PetscLogEvent MAT_BackwardSolve;
1604: PETSC_EXTERN PetscLogEvent MAT_LUFactor;
1605: PETSC_EXTERN PetscLogEvent MAT_LUFactorSymbolic;
1606: PETSC_EXTERN PetscLogEvent MAT_LUFactorNumeric;
1607: PETSC_EXTERN PetscLogEvent MAT_CholeskyFactor;
1608: PETSC_EXTERN PetscLogEvent MAT_CholeskyFactorSymbolic;
1609: PETSC_EXTERN PetscLogEvent MAT_CholeskyFactorNumeric;
1610: PETSC_EXTERN PetscLogEvent MAT_ILUFactor;
1611: PETSC_EXTERN PetscLogEvent MAT_ILUFactorSymbolic;
1612: PETSC_EXTERN PetscLogEvent MAT_ICCFactorSymbolic;
1613: PETSC_EXTERN PetscLogEvent MAT_Copy;
1614: PETSC_EXTERN PetscLogEvent MAT_Convert;
1615: PETSC_EXTERN PetscLogEvent MAT_Scale;
1616: PETSC_EXTERN PetscLogEvent MAT_AssemblyBegin;
1617: PETSC_EXTERN PetscLogEvent MAT_AssemblyEnd;
1618: PETSC_EXTERN PetscLogEvent MAT_SetValues;
1619: PETSC_EXTERN PetscLogEvent MAT_GetValues;
1620: PETSC_EXTERN PetscLogEvent MAT_GetRow;
1621: PETSC_EXTERN PetscLogEvent MAT_GetRowIJ;
1622: PETSC_EXTERN PetscLogEvent MAT_CreateSubMats;
1623: PETSC_EXTERN PetscLogEvent MAT_GetColoring;
1624: PETSC_EXTERN PetscLogEvent MAT_GetOrdering;
1625: PETSC_EXTERN PetscLogEvent MAT_RedundantMat;
1626: PETSC_EXTERN PetscLogEvent MAT_IncreaseOverlap;
1627: PETSC_EXTERN PetscLogEvent MAT_Partitioning;
1628: PETSC_EXTERN PetscLogEvent MAT_Coarsen;
1629: PETSC_EXTERN PetscLogEvent MAT_ZeroEntries;
1630: PETSC_EXTERN PetscLogEvent MAT_Load;
1631: PETSC_EXTERN PetscLogEvent MAT_View;
1632: PETSC_EXTERN PetscLogEvent MAT_AXPY;
1633: PETSC_EXTERN PetscLogEvent MAT_FDColoringCreate;
1634: PETSC_EXTERN PetscLogEvent MAT_TransposeColoringCreate;
1635: PETSC_EXTERN PetscLogEvent MAT_FDColoringSetUp;
1636: PETSC_EXTERN PetscLogEvent MAT_FDColoringApply;
1637: PETSC_EXTERN PetscLogEvent MAT_Transpose;
1638: PETSC_EXTERN PetscLogEvent MAT_FDColoringFunction;
1639: PETSC_EXTERN PetscLogEvent MAT_CreateSubMat;
1640: PETSC_EXTERN PetscLogEvent MAT_MatMult;
1641: PETSC_EXTERN PetscLogEvent MAT_MatSolve;
1642: PETSC_EXTERN PetscLogEvent MAT_MatMultSymbolic;
1643: PETSC_EXTERN PetscLogEvent MAT_MatMultNumeric;
1644: PETSC_EXTERN PetscLogEvent MAT_Getlocalmatcondensed;
1645: PETSC_EXTERN PetscLogEvent MAT_GetBrowsOfAcols;
1646: PETSC_EXTERN PetscLogEvent MAT_GetBrowsOfAocols;
1647: PETSC_EXTERN PetscLogEvent MAT_PtAP;
1648: PETSC_EXTERN PetscLogEvent MAT_PtAPSymbolic;
1649: PETSC_EXTERN PetscLogEvent MAT_PtAPNumeric;
1650: PETSC_EXTERN PetscLogEvent MAT_Seqstompinum;
1651: PETSC_EXTERN PetscLogEvent MAT_Seqstompisym;
1652: PETSC_EXTERN PetscLogEvent MAT_Seqstompi;
1653: PETSC_EXTERN PetscLogEvent MAT_Getlocalmat;
1654: PETSC_EXTERN PetscLogEvent MAT_RARt;
1655: PETSC_EXTERN PetscLogEvent MAT_RARtSymbolic;
1656: PETSC_EXTERN PetscLogEvent MAT_RARtNumeric;
1657: PETSC_EXTERN PetscLogEvent MAT_MatTransposeMult;
1658: PETSC_EXTERN PetscLogEvent MAT_MatTransposeMultSymbolic;
1659: PETSC_EXTERN PetscLogEvent MAT_MatTransposeMultNumeric;
1660: PETSC_EXTERN PetscLogEvent MAT_TransposeMatMult;
1661: PETSC_EXTERN PetscLogEvent MAT_TransposeMatMultSymbolic;
1662: PETSC_EXTERN PetscLogEvent MAT_TransposeMatMultNumeric;
1663: PETSC_EXTERN PetscLogEvent MAT_MatMatMult;
1664: PETSC_EXTERN PetscLogEvent MAT_MatMatMultSymbolic;
1665: PETSC_EXTERN PetscLogEvent MAT_MatMatMultNumeric;
1666: PETSC_EXTERN PetscLogEvent MAT_Applypapt;
1667: PETSC_EXTERN PetscLogEvent MAT_Applypapt_symbolic;
1668: PETSC_EXTERN PetscLogEvent MAT_Applypapt_numeric;
1669: PETSC_EXTERN PetscLogEvent MAT_Getsymtranspose;
1670: PETSC_EXTERN PetscLogEvent MAT_Transpose_SeqAIJ;
1671: PETSC_EXTERN PetscLogEvent MAT_Getsymtransreduced;
1672: PETSC_EXTERN PetscLogEvent MAT_GetSequentialNonzeroStructure;
1673: PETSC_EXTERN PetscLogEvent MATMFFD_Mult;
1674: PETSC_EXTERN PetscLogEvent MAT_GetMultiProcBlock;
1675: PETSC_EXTERN PetscLogEvent MAT_CUSPARSECopyToGPU;
1676: PETSC_EXTERN PetscLogEvent MAT_SetValuesBatch;
1677: PETSC_EXTERN PetscLogEvent MAT_SetValuesBatchI;
1678: PETSC_EXTERN PetscLogEvent MAT_SetValuesBatchII;
1679: PETSC_EXTERN PetscLogEvent MAT_SetValuesBatchIII;
1680: PETSC_EXTERN PetscLogEvent MAT_SetValuesBatchIV;
1681: PETSC_EXTERN PetscLogEvent MAT_ViennaCLCopyToGPU;
1682: PETSC_EXTERN PetscLogEvent MAT_Merge;
1683: PETSC_EXTERN PetscLogEvent MAT_Residual;
1684: PETSC_EXTERN PetscLogEvent MAT_SetRandom;
1685: PETSC_EXTERN PetscLogEvent MATCOLORING_Apply;
1686: PETSC_EXTERN PetscLogEvent MATCOLORING_Comm;
1687: PETSC_EXTERN PetscLogEvent MATCOLORING_Local;
1688: PETSC_EXTERN PetscLogEvent MATCOLORING_ISCreate;
1689: PETSC_EXTERN PetscLogEvent MATCOLORING_SetUp;
1690: PETSC_EXTERN PetscLogEvent MATCOLORING_Weights;

1692: #endif