Actual source code: petsc-matimpl.h

petsc-3.4.2 2013-07-02
  2: #ifndef __MATIMPL_H

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

  8: /*
  9:   This file defines the parts of the matrix data structure that are
 10:   shared by all matrix types.
 11: */

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

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

199: typedef struct _p_MatBaseName* MatBaseName;
200: struct _p_MatBaseName {
201:   char        *bname,*sname,*mname;
202:   MatBaseName next;
203: };

205: PETSC_EXTERN MatBaseName MatBaseNameList;

207: /*
208:    Utility private matrix routines
209: */
210: PETSC_INTERN PetscErrorCode MatConvert_Basic(Mat, MatType,MatReuse,Mat*);
211: PETSC_INTERN PetscErrorCode MatCopy_Basic(Mat,Mat,MatStructure);
212: PETSC_INTERN PetscErrorCode MatHeaderMerge(Mat,Mat);
213: PETSC_EXTERN PetscErrorCode MatHeaderReplace(Mat,Mat);
214: PETSC_INTERN PetscErrorCode MatAXPYGetxtoy_Private(PetscInt,PetscInt*,PetscInt*,PetscInt*, PetscInt*,PetscInt*,PetscInt*, PetscInt**);
215: PETSC_INTERN PetscErrorCode MatDiagonalSet_Default(Mat,Vec,InsertMode);

217: #if defined(PETSC_USE_DEBUG)
218: #  define MatCheckPreallocated(A,arg) do {                              \
219:     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); \
220:   } while (0)
221: #else
222: #  define MatCheckPreallocated(A,arg) do {} while (0)
223: #endif

225: /*
226:   The stash is used to temporarily store inserted matrix values that
227:   belong to another processor. During the assembly phase the stashed
228:   values are moved to the correct processor and
229: */

231: typedef struct _MatStashSpace *PetscMatStashSpace;

233: struct _MatStashSpace {
234:   PetscMatStashSpace next;
235:   PetscScalar        *space_head,*val;
236:   PetscInt           *idx,*idy;
237:   PetscInt           total_space_size;
238:   PetscInt           local_used;
239:   PetscInt           local_remaining;
240: };

242: PETSC_EXTERN PetscErrorCode PetscMatStashSpaceGet(PetscInt,PetscInt,PetscMatStashSpace *);
243: PETSC_EXTERN PetscErrorCode PetscMatStashSpaceContiguous(PetscInt,PetscMatStashSpace *,PetscScalar *,PetscInt *,PetscInt *);
244: PETSC_EXTERN PetscErrorCode PetscMatStashSpaceDestroy(PetscMatStashSpace*);

246: typedef struct {
247:   PetscInt      nmax;                   /* maximum stash size */
248:   PetscInt      umax;                   /* user specified max-size */
249:   PetscInt      oldnmax;                /* the nmax value used previously */
250:   PetscInt      n;                      /* stash size */
251:   PetscInt      bs;                     /* block size of the stash */
252:   PetscInt      reallocs;               /* preserve the no of mallocs invoked */
253:   PetscMatStashSpace space_head,space;  /* linked list to hold stashed global row/column numbers and matrix values */
254:   /* The following variables are used for communication */
255:   MPI_Comm      comm;
256:   PetscMPIInt   size,rank;
257:   PetscMPIInt   tag1,tag2;
258:   MPI_Request   *send_waits;            /* array of send requests */
259:   MPI_Request   *recv_waits;            /* array of receive requests */
260:   MPI_Status    *send_status;           /* array of send status */
261:   PetscInt      nsends,nrecvs;          /* numbers of sends and receives */
262:   PetscScalar   *svalues;               /* sending data */
263:   PetscInt      *sindices;
264:   PetscScalar   **rvalues;              /* receiving data (values) */
265:   PetscInt      **rindices;             /* receiving data (indices) */
266:   PetscInt      nprocessed;             /* number of messages already processed */
267:   PetscMPIInt   *flg_v;                 /* indicates what messages have arrived so far and from whom */
268:   PetscBool     reproduce;
269:   PetscInt      reproduce_count;
270: } MatStash;

272: PETSC_INTERN PetscErrorCode MatStashCreate_Private(MPI_Comm,PetscInt,MatStash*);
273: PETSC_INTERN PetscErrorCode MatStashDestroy_Private(MatStash*);
274: PETSC_INTERN PetscErrorCode MatStashScatterEnd_Private(MatStash*);
275: PETSC_INTERN PetscErrorCode MatStashSetInitialSize_Private(MatStash*,PetscInt);
276: PETSC_INTERN PetscErrorCode MatStashGetInfo_Private(MatStash*,PetscInt*,PetscInt*);
277: PETSC_INTERN PetscErrorCode MatStashValuesRow_Private(MatStash*,PetscInt,PetscInt,const PetscInt[],const PetscScalar[],PetscBool );
278: PETSC_INTERN PetscErrorCode MatStashValuesCol_Private(MatStash*,PetscInt,PetscInt,const PetscInt[],const PetscScalar[],PetscInt,PetscBool );
279: PETSC_INTERN PetscErrorCode MatStashValuesRowBlocked_Private(MatStash*,PetscInt,PetscInt,const PetscInt[],const PetscScalar[],PetscInt,PetscInt,PetscInt);
280: PETSC_INTERN PetscErrorCode MatStashValuesColBlocked_Private(MatStash*,PetscInt,PetscInt,const PetscInt[],const PetscScalar[],PetscInt,PetscInt,PetscInt);
281: PETSC_INTERN PetscErrorCode MatStashScatterBegin_Private(Mat,MatStash*,PetscInt*);
282: PETSC_INTERN PetscErrorCode MatStashScatterGetMesg_Private(MatStash*,PetscMPIInt*,PetscInt**,PetscInt**,PetscScalar**,PetscInt*);

284: typedef struct {
285:   PetscInt   dim;
286:   PetscInt   dims[4];
287:   PetscInt   starts[4];
288:   PetscBool  noc;        /* this is a single component problem, hence user will not set MatStencil.c */
289: } MatStencilInfo;

291: /* Info about using compressed row format */
292: typedef struct {
293:   PetscBool  check;                         /* indicates that at MatAssembly() it should check if compressed rows will be efficient */
294:   PetscBool  use;                           /* indicates compressed rows have been checked and will be used */
295:   PetscInt   nrows;                         /* number of non-zero rows */
296:   PetscInt   *i;                            /* compressed row pointer  */
297:   PetscInt   *rindex;                       /* compressed row index               */
298: } Mat_CompressedRow;
299: PETSC_EXTERN PetscErrorCode MatCheckCompressedRow(Mat,Mat_CompressedRow*,PetscInt*,PetscInt,PetscReal);

301: struct _p_Mat {
302:   PETSCHEADER(struct _MatOps);
303:   PetscLayout            rmap,cmap;
304:   void                   *data;            /* implementation-specific data */
305:   MatFactorType          factortype;       /* MAT_FACTOR_LU, ILU, CHOLESKY or ICC */
306:   PetscBool              assembled;        /* is the matrix assembled? */
307:   PetscBool              was_assembled;    /* new values inserted into assembled mat */
308:   PetscInt               num_ass;          /* number of times matrix has been assembled */
309:   PetscBool              same_nonzero;     /* matrix has same nonzero pattern as previous */
310:   MatInfo                info;             /* matrix information */
311:   InsertMode             insertmode;       /* have values been inserted in matrix or added? */
312:   MatStash               stash,bstash;     /* used for assembling off-proc mat emements */
313:   MatNullSpace           nullsp;           /* null space (operator is singular) */
314:   MatNullSpace           nearnullsp;       /* near null space to be used by multigrid methods */
315:   PetscBool              preallocated;
316:   MatStencilInfo         stencil;          /* information for structured grid */
317:   PetscBool              symmetric,hermitian,structurally_symmetric,spd;
318:   PetscBool              symmetric_set,hermitian_set,structurally_symmetric_set,spd_set; /* if true, then corresponding flag is correct*/
319:   PetscBool              symmetric_eternal;
320:   PetscBool              nooffprocentries,nooffproczerorows;
321: #if defined(PETSC_HAVE_CUSP)
322:   PetscCUSPFlag          valid_GPU_matrix; /* flag pointing to the matrix on the gpu*/
323: #endif
324:   void                   *spptr;          /* pointer for special library like SuperLU */
325:   MatSolverPackage       solvertype;
326:   PetscViewer            viewonassembly;         /* the following are set in MatSetFromOptions() and used in MatAssemblyEnd() */
327:   PetscViewerFormat      viewformatonassembly;
328:   PetscBool              checksymmetryonassembly,checknullspaceonassembly;
329:   PetscReal              checksymmetrytol;
330:   };

332: PETSC_INTERN PetscErrorCode MatAXPY_Basic(Mat,PetscScalar,Mat,MatStructure);
333: PETSC_INTERN PetscErrorCode MatAXPY_BasicWithPreallocation(Mat,Mat,PetscScalar,Mat,MatStructure);
334: /*
335:     Object for partitioning graphs
336: */

338: typedef struct _MatPartitioningOps *MatPartitioningOps;
339: struct _MatPartitioningOps {
340:   PetscErrorCode (*apply)(MatPartitioning,IS*);
341:   PetscErrorCode (*setfromoptions)(MatPartitioning);
342:   PetscErrorCode (*destroy)(MatPartitioning);
343:   PetscErrorCode (*view)(MatPartitioning,PetscViewer);
344: };

346: struct _p_MatPartitioning {
347:   PETSCHEADER(struct _MatPartitioningOps);
348:   Mat         adj;
349:   PetscInt    *vertex_weights;
350:   PetscReal   *part_weights;
351:   PetscInt    n;                                 /* number of partitions */
352:   void        *data;
353:   PetscInt    setupcalled;
354: };

356: /*
357:     Object for coarsen graphs
358: */
359: typedef struct _MatCoarsenOps *MatCoarsenOps;
360: struct _MatCoarsenOps {
361:   PetscErrorCode (*apply)(MatCoarsen);
362:   PetscErrorCode (*setfromoptions)(MatCoarsen);
363:   PetscErrorCode (*destroy)(MatCoarsen);
364:   PetscErrorCode (*view)(MatCoarsen,PetscViewer);
365: };

367: struct _p_MatCoarsen {
368:   PETSCHEADER(struct _MatCoarsenOps);
369:   Mat         graph;
370:   PetscInt    verbose;
371:   PetscInt    setupcalled;
372:   void        *subctx;
373:   /* */
374:   PetscBool strict_aggs;
375:   IS perm;
376:   PetscCoarsenData *agg_lists;
377: };

379: PETSC_EXTERN PetscErrorCode PetscCDCreate(PetscInt,PetscCoarsenData**);
380: PETSC_EXTERN PetscErrorCode PetscCDDestroy(PetscCoarsenData*);
381: PETSC_EXTERN PetscErrorCode PetscLLNSetID(PetscCDIntNd*,PetscInt);
382: PETSC_EXTERN PetscErrorCode PetscLLNGetID(const PetscCDIntNd*,PetscInt*);
383: PETSC_EXTERN PetscErrorCode PetscCDAppendID(PetscCoarsenData*,PetscInt,PetscInt);
384: PETSC_EXTERN PetscErrorCode PetscCDAppendRemove(PetscCoarsenData*,PetscInt,PetscInt);
385: PETSC_EXTERN PetscErrorCode PetscCDAppendNode(PetscCoarsenData*,PetscInt,PetscCDIntNd*);
386: PETSC_EXTERN PetscErrorCode PetscCDRemoveNextNode(PetscCoarsenData*,PetscInt,PetscCDIntNd*);
387: PETSC_EXTERN PetscErrorCode PetscCDRemoveAllAt(PetscCoarsenData*,PetscInt);
388: PETSC_EXTERN PetscErrorCode PetscCDSizeAt(const PetscCoarsenData*,PetscInt,PetscInt*);
389: PETSC_EXTERN PetscErrorCode PetscCDEmptyAt(const PetscCoarsenData*,PetscInt,PetscBool*);
390: PETSC_EXTERN PetscErrorCode PetscCDSetChuckSize(PetscCoarsenData*,PetscInt);
391: PETSC_EXTERN PetscErrorCode PetscCDPrint(const PetscCoarsenData*,MPI_Comm);
392: PETSC_EXTERN PetscErrorCode PetscCDGetMIS(PetscCoarsenData*,IS*);
393: PETSC_EXTERN PetscErrorCode PetscCDGetMat(const PetscCoarsenData*,Mat*);
394: PETSC_EXTERN PetscErrorCode PetscCDSetMat(PetscCoarsenData*,Mat);

396: typedef PetscCDIntNd *PetscCDPos;
397: PETSC_EXTERN PetscErrorCode PetscCDGetHeadPos(const PetscCoarsenData*,PetscInt,PetscCDPos*);
398: PETSC_EXTERN PetscErrorCode PetscCDGetNextPos(const PetscCoarsenData*,PetscInt,PetscCDPos*);
399: PETSC_EXTERN PetscErrorCode PetscCDGetASMBlocks(const PetscCoarsenData*,const PetscInt,PetscInt*,IS**);
400: /* PetscErrorCode PetscCDSetRemovedIS( PetscCoarsenData *ail, MPI_Comm, const PetscInt, PetscInt[] ); */
401: /* PetscErrorCode PetscCDGetRemovedIS( PetscCoarsenData *ail, IS * ); */

403: /*
404:     MatFDColoring is used to compute Jacobian matrices efficiently
405:   via coloring. The data structure is explained below in an example.

407:    Color =   0    1     0    2   |   2      3       0
408:    ---------------------------------------------------
409:             00   01              |          05
410:             10   11              |   14     15               Processor  0
411:                        22    23  |          25
412:                        32    33  |
413:    ===================================================
414:                                  |   44     45     46
415:             50                   |          55               Processor 1
416:                                  |   64            66
417:    ---------------------------------------------------

419:     ncolors = 4;

421:     ncolumns      = {2,1,1,0}
422:     columns       = {{0,2},{1},{3},{}}
423:     nrows         = {4,2,3,3}
424:     rows          = {{0,1,2,3},{0,1},{1,2,3},{0,1,2}}
425:     columnsforrow = {{0,0,2,2},{1,1},{4,3,3},{5,5,5}}
426:     vscaleforrow  = {{,,,},{,},{,,},{,,}}
427:     vwscale       = {dx(0),dx(1),dx(2),dx(3)}               MPI Vec
428:     vscale        = {dx(0),dx(1),dx(2),dx(3),dx(4),dx(5)}   Seq Vec

430:     ncolumns      = {1,0,1,1}
431:     columns       = {{6},{},{4},{5}}
432:     nrows         = {3,0,2,2}
433:     rows          = {{0,1,2},{},{1,2},{1,2}}
434:     columnsforrow = {{6,0,6},{},{4,4},{5,5}}
435:     vscaleforrow =  {{,,},{},{,},{,}}
436:     vwscale       = {dx(4),dx(5),dx(6)}              MPI Vec
437:     vscale        = {dx(0),dx(4),dx(5),dx(6)}        Seq Vec

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

442: */

444: struct  _p_MatFDColoring{
445:   PETSCHEADER(int);
446:   PetscInt       M,N,m;            /* total rows, columns; local rows */
447:   PetscInt       rstart;           /* first row owned by local processor */
448:   PetscInt       ncolors;          /* number of colors */
449:   PetscInt       *ncolumns;        /* number of local columns for a color */
450:   PetscInt       **columns;        /* lists the local columns of each color (using global column numbering) */
451:   PetscInt       *nrows;           /* number of local rows for each color */
452:   PetscInt       **rows;           /* lists the local rows for each color (using the local row numbering) */
453:   PetscInt       **columnsforrow;  /* lists the corresponding columns for those rows (using the global column) */
454:   PetscReal      error_rel;        /* square root of relative error in computing function */
455:   PetscReal      umin;             /* minimum allowable u'dx value */
456:   Vec            w1,w2,w3;         /* work vectors used in computing Jacobian */
457:   PetscBool      fset;             /* indicates that the initial function value F(X) is set */
458:   PetscErrorCode (*f)(void);       /* function that defines Jacobian */
459:   void           *fctx;            /* optional user-defined context for use by the function f */
460:   PetscInt       **vscaleforrow;   /* location in vscale for each columnsforrow[] entry */
461:   Vec            vscale;           /* holds FD scaling, i.e. 1/dx for each perturbed column */
462:   PetscInt       currentcolor;     /* color for which function evaluation is being done now */
463:   const char     *htype;            /* "wp" or "ds" */
464:   ISColoringType ctype;            /* IS_COLORING_GLOBAL or IS_COLORING_GHOSTED */

466:   void           *ftn_func_pointer,*ftn_func_cntx; /* serve the same purpose as *fortran_func_pointers in PETSc objects */
467: };

469: struct  _p_MatTransposeColoring{
470:   PETSCHEADER(int);
471:   PetscInt       M,N,m;            /* total rows, columns; local rows */
472:   PetscInt       rstart;           /* first row owned by local processor */
473:   PetscInt       ncolors;          /* number of colors */
474:   PetscInt       *ncolumns;        /* number of local columns for a color */
475:   PetscInt       *nrows;           /* number of local rows for each color */
476:   PetscInt       currentcolor;     /* color for which function evaluation is being done now */
477:   ISColoringType ctype;            /* IS_COLORING_GLOBAL or IS_COLORING_GHOSTED */

479:   PetscInt       *colorforrow,*colorforcol;  /* pointer to rows and columns */
480:   PetscInt       *rows;                  /* lists the local rows for each color (using the local row numbering) */
481:   PetscInt       *columnsforspidx;       /* maps (row,color) in the dense matrix to index of sparse matrix arrays a->j and a->a */
482:   PetscInt       *columns;               /* lists the local columns of each color (using global column numbering) */
483: };

485: /*
486:    Null space context for preconditioner/operators
487: */
488: struct _p_MatNullSpace {
489:   PETSCHEADER(int);
490:   PetscBool      has_cnst;
491:   PetscInt       n;
492:   Vec*           vecs;
493:   PetscScalar*   alpha;                 /* for projections */
494:   Vec            vec;                   /* for out of place removals */
495:   PetscErrorCode (*remove)(MatNullSpace,Vec,void*);  /* for user provided removal function */
496:   void*          rmctx;                 /* context for remove() function */
497: };

499: /*
500:    Checking zero pivot for LU, ILU preconditioners.
501: */
502: typedef struct {
503:   PetscInt       nshift,nshift_max;
504:   PetscReal      shift_amount,shift_lo,shift_hi,shift_top,shift_fraction;
505:   PetscBool      newshift;
506:   PetscReal      rs;  /* active row sum of abs(offdiagonals) */
507:   PetscScalar    pv;  /* pivot of the active row */
508: } FactorShiftCtx;

510: PETSC_EXTERN PetscErrorCode MatFactorDumpMatrix(Mat);

514: PETSC_STATIC_INLINE PetscErrorCode MatPivotCheck_nz(Mat mat,const MatFactorInfo *info,FactorShiftCtx *sctx,PetscInt row)
515: {
516:   PetscReal _rs   = sctx->rs;
517:   PetscReal _zero = info->zeropivot*_rs;

520:   if (PetscAbsScalar(sctx->pv) <= _zero){
521:     /* force |diag| > zeropivot*rs */
522:     if (!sctx->nshift) sctx->shift_amount = info->shiftamount;
523:     else sctx->shift_amount *= 2.0;
524:     sctx->newshift = PETSC_TRUE;
525:     (sctx->nshift)++;
526:   } else {
527:     sctx->newshift = PETSC_FALSE;
528:   }
529:   return(0);
530: }

534: PETSC_STATIC_INLINE PetscErrorCode MatPivotCheck_pd(Mat mat,const MatFactorInfo *info,FactorShiftCtx *sctx,PetscInt row)
535: {
536:   PetscReal _rs   = sctx->rs;
537:   PetscReal _zero = info->zeropivot*_rs;

540:   if (PetscRealPart(sctx->pv) <= _zero){
541:     /* force matfactor to be diagonally dominant */
542:     if (sctx->nshift == sctx->nshift_max) {
543:       sctx->shift_fraction = sctx->shift_hi;
544:     } else {
545:       sctx->shift_lo = sctx->shift_fraction;
546:       sctx->shift_fraction = (sctx->shift_hi+sctx->shift_lo)/2.;
547:     }
548:     sctx->shift_amount = sctx->shift_fraction * sctx->shift_top;
549:     sctx->nshift++;
550:     sctx->newshift = PETSC_TRUE;
551:   } else {
552:     sctx->newshift = PETSC_FALSE;
553:   }
554:   return(0);
555: }

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

564:   if (PetscAbsScalar(sctx->pv) <= _zero){
565:     sctx->pv          += info->shiftamount;
566:     sctx->shift_amount = 0.0;
567:     sctx->nshift++;
568:   }
569:   sctx->newshift = PETSC_FALSE;
570:   return(0);
571: }

575: PETSC_STATIC_INLINE PetscErrorCode MatPivotCheck_none(Mat mat,const MatFactorInfo *info,FactorShiftCtx *sctx,PetscInt row)
576: {
577:   PetscReal _zero = info->zeropivot;

580:   sctx->newshift = PETSC_FALSE;
581:   if (PetscAbsScalar(sctx->pv) <= _zero) {
583:     PetscBool      flg = PETSC_FALSE;

585:     PetscOptionsGetBool(NULL,"-mat_dump",&flg,NULL);
586:     if (flg) {
587:       MatView(mat,PETSC_VIEWER_BINARY_(PetscObjectComm((PetscObject)mat)));
588:     }
589:     SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot row %D value %G tolerance %G",row,PetscAbsScalar(sctx->pv),_zero);
590:   }
591:   return(0);
592: }

596: PETSC_STATIC_INLINE PetscErrorCode MatPivotCheck(Mat mat,const MatFactorInfo *info,FactorShiftCtx *sctx,PetscInt row)
597: {

601:   if (info->shifttype == (PetscReal) MAT_SHIFT_NONZERO){
602:     MatPivotCheck_nz(mat,info,sctx,row);
603:   } else if (info->shifttype == (PetscReal) MAT_SHIFT_POSITIVE_DEFINITE){
604:     MatPivotCheck_pd(mat,info,sctx,row);
605:   } else if (info->shifttype == (PetscReal) MAT_SHIFT_INBLOCKS){
606:     MatPivotCheck_inblocks(mat,info,sctx,row);
607:   } else {
608:     MatPivotCheck_none(mat,info,sctx,row);
609:   }
610:   return(0);
611: }

613: /*
614:   Create and initialize a linked list
615:   Input Parameters:
616:     idx_start - starting index of the list
617:     lnk_max   - max value of lnk indicating the end of the list
618:     nlnk      - max length of the list
619:   Output Parameters:
620:     lnk       - list initialized
621:     bt        - PetscBT (bitarray) with all bits set to false
622:     lnk_empty - flg indicating the list is empty
623: */
624: #define PetscLLCreate(idx_start,lnk_max,nlnk,lnk,bt) \
625:   (PetscMalloc(nlnk*sizeof(PetscInt),&lnk) || PetscBTCreate(nlnk,&(bt)) || (lnk[idx_start] = lnk_max,0))

627: #define PetscLLCreate_new(idx_start,lnk_max,nlnk,lnk,bt,lnk_empty)\
628:   (PetscMalloc(nlnk*sizeof(PetscInt),&lnk) || PetscBTCreate(nlnk,&(bt)) || (lnk_empty = PETSC_TRUE,0) ||(lnk[idx_start] = lnk_max,0))

630: /*
631:   Add an index set into a sorted linked list
632:   Input Parameters:
633:     nidx      - number of input indices
634:     indices   - interger array
635:     idx_start - starting index of the list
636:     lnk       - linked list(an integer array) that is created
637:     bt        - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
638:   output Parameters:
639:     nlnk      - number of newly added indices
640:     lnk       - the sorted(increasing order) linked list containing new and non-redundate entries from indices
641:     bt        - updated PetscBT (bitarray)
642: */
643: #define PetscLLAdd(nidx,indices,idx_start,nlnk,lnk,bt) 0;\
644: {\
645:   PetscInt _k,_entry,_location,_lnkdata;\
646:   nlnk     = 0;\
647:   _lnkdata = idx_start;\
648:   for (_k=0; _k<nidx; _k++){\
649:     _entry = indices[_k];\
650:     if (!PetscBTLookupSet(bt,_entry)){  /* new entry */\
651:       /* search for insertion location */\
652:       /* start from the beginning if _entry < previous _entry */\
653:       if (_k && _entry < _lnkdata) _lnkdata  = idx_start;\
654:       do {\
655:         _location = _lnkdata;\
656:         _lnkdata  = lnk[_location];\
657:       } while (_entry > _lnkdata);\
658:       /* insertion location is found, add entry into lnk */\
659:       lnk[_location] = _entry;\
660:       lnk[_entry]    = _lnkdata;\
661:       nlnk++;\
662:       _lnkdata = _entry; /* next search starts from here if next_entry > _entry */\
663:     }\
664:   }\
665: }

667: /*
668:   Add a permuted index set into a sorted linked list
669:   Input Parameters:
670:     nidx      - number of input indices
671:     indices   - interger array
672:     perm      - permutation of indices
673:     idx_start - starting index of the list
674:     lnk       - linked list(an integer array) that is created
675:     bt        - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
676:   output Parameters:
677:     nlnk      - number of newly added indices
678:     lnk       - the sorted(increasing order) linked list containing new and non-redundate entries from indices
679:     bt        - updated PetscBT (bitarray)
680: */
681: #define PetscLLAddPerm(nidx,indices,perm,idx_start,nlnk,lnk,bt) 0;\
682: {\
683:   PetscInt _k,_entry,_location,_lnkdata;\
684:   nlnk     = 0;\
685:   _lnkdata = idx_start;\
686:   for (_k=0; _k<nidx; _k++){\
687:     _entry = perm[indices[_k]];\
688:     if (!PetscBTLookupSet(bt,_entry)){  /* new entry */\
689:       /* search for insertion location */\
690:       /* start from the beginning if _entry < previous _entry */\
691:       if (_k && _entry < _lnkdata) _lnkdata  = idx_start;\
692:       do {\
693:         _location = _lnkdata;\
694:         _lnkdata  = lnk[_location];\
695:       } while (_entry > _lnkdata);\
696:       /* insertion location is found, add entry into lnk */\
697:       lnk[_location] = _entry;\
698:       lnk[_entry]    = _lnkdata;\
699:       nlnk++;\
700:       _lnkdata = _entry; /* next search starts from here if next_entry > _entry */\
701:     }\
702:   }\
703: }

705: /*
706:   Add a SORTED ascending index set into a sorted linked list - same as PetscLLAdd() bus skip 'if (_k && _entry < _lnkdata) _lnkdata  = idx_start;'
707:   Input Parameters:
708:     nidx      - number of input indices
709:     indices   - sorted interger array
710:     idx_start - starting index of the list
711:     lnk       - linked list(an integer array) that is created
712:     bt        - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
713:   output Parameters:
714:     nlnk      - number of newly added indices
715:     lnk       - the sorted(increasing order) linked list containing new and non-redundate entries from indices
716:     bt        - updated PetscBT (bitarray)
717: */
718: #define PetscLLAddSorted(nidx,indices,idx_start,nlnk,lnk,bt) 0;\
719: {\
720:   PetscInt _k,_entry,_location,_lnkdata;\
721:   nlnk      = 0;\
722:   _lnkdata  = idx_start;\
723:   for (_k=0; _k<nidx; _k++){\
724:     _entry = indices[_k];\
725:     if (!PetscBTLookupSet(bt,_entry)){  /* new entry */\
726:       /* search for insertion location */\
727:       do {\
728:         _location = _lnkdata;\
729:         _lnkdata  = lnk[_location];\
730:       } while (_entry > _lnkdata);\
731:       /* insertion location is found, add entry into lnk */\
732:       lnk[_location] = _entry;\
733:       lnk[_entry]    = _lnkdata;\
734:       nlnk++;\
735:       _lnkdata = _entry; /* next search starts from here */\
736:     }\
737:   }\
738: }

740: #define PetscLLAddSorted_new(nidx,indices,idx_start,lnk_empty,nlnk,lnk,bt) 0; \
741: {\
742:   PetscInt _k,_entry,_location,_lnkdata;\
743:   if (lnk_empty){\
744:     _lnkdata  = idx_start;                      \
745:     for (_k=0; _k<nidx; _k++){                  \
746:       _entry = indices[_k];                             \
747:       PetscBTSet(bt,_entry);  /* mark the new entry */          \
748:           _location = _lnkdata;                                 \
749:           _lnkdata  = lnk[_location];                           \
750:         /* insertion location is found, add entry into lnk */   \
751:         lnk[_location] = _entry;                                \
752:         lnk[_entry]    = _lnkdata;                              \
753:         _lnkdata = _entry; /* next search starts from here */   \
754:     }                                                           \
755:     /*\
756:     lnk[indices[nidx-1]] = lnk[idx_start];\
757:     lnk[idx_start]       = indices[0];\
758:     PetscBTSet(bt,indices[0]);  \
759:     for (_k=1; _k<nidx; _k++){                  \
760:       PetscBTSet(bt,indices[_k]);                                          \
761:       lnk[indices[_k-1]] = indices[_k];                                  \
762:     }                                                           \
763:      */\
764:     nlnk      = nidx;\
765:     lnk_empty = PETSC_FALSE;\
766:   } else {\
767:     nlnk      = 0;                              \
768:     _lnkdata  = idx_start;                      \
769:     for (_k=0; _k<nidx; _k++){                  \
770:       _entry = indices[_k];                             \
771:       if (!PetscBTLookupSet(bt,_entry)){  /* new entry */       \
772:         /* search for insertion location */                     \
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 */   \
782:       }                                                         \
783:     }                                                           \
784:   }                                                             \
785: }

787: /*
788:   Add a SORTED index set into a sorted linked list used for LUFactorSymbolic()
789:   Same as PetscLLAddSorted() with an additional operation:
790:        count the number of input indices that are no larger than 'diag'
791:   Input Parameters:
792:     indices   - sorted interger array
793:     idx_start - starting index of the list, index of pivot row
794:     lnk       - linked list(an integer array) that is created
795:     bt        - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
796:     diag      - index of the active row in LUFactorSymbolic
797:     nzbd      - number of input indices with indices <= idx_start
798:     im        - im[idx_start] is initialized as num of nonzero entries in row=idx_start
799:   output Parameters:
800:     nlnk      - number of newly added indices
801:     lnk       - the sorted(increasing order) linked list containing new and non-redundate entries from indices
802:     bt        - updated PetscBT (bitarray)
803:     im        - im[idx_start]: unchanged if diag is not an entry
804:                              : num of entries with indices <= diag if diag is an entry
805: */
806: #define PetscLLAddSortedLU(indices,idx_start,nlnk,lnk,bt,diag,nzbd,im) 0;\
807: {\
808:   PetscInt _k,_entry,_location,_lnkdata,_nidx;\
809:   nlnk     = 0;\
810:   _lnkdata = idx_start;\
811:   _nidx = im[idx_start] - nzbd; /* num of entries with idx_start < index <= diag */\
812:   for (_k=0; _k<_nidx; _k++){\
813:     _entry = indices[_k];\
814:     nzbd++;\
815:     if ( _entry== diag) im[idx_start] = nzbd;\
816:     if (!PetscBTLookupSet(bt,_entry)){  /* new entry */\
817:       /* search for insertion location */\
818:       do {\
819:         _location = _lnkdata;\
820:         _lnkdata  = lnk[_location];\
821:       } while (_entry > _lnkdata);\
822:       /* insertion location is found, add entry into lnk */\
823:       lnk[_location] = _entry;\
824:       lnk[_entry]    = _lnkdata;\
825:       nlnk++;\
826:       _lnkdata = _entry; /* next search starts from here */\
827:     }\
828:   }\
829: }

831: /*
832:   Copy data on the list into an array, then initialize the list
833:   Input Parameters:
834:     idx_start - starting index of the list
835:     lnk_max   - max value of lnk indicating the end of the list
836:     nlnk      - number of data on the list to be copied
837:     lnk       - linked list
838:     bt        - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
839:   output Parameters:
840:     indices   - array that contains the copied data
841:     lnk       - linked list that is cleaned and initialize
842:     bt        - PetscBT (bitarray) with all bits set to false
843: */
844: #define PetscLLClean(idx_start,lnk_max,nlnk,lnk,indices,bt) 0;\
845: {\
846:   PetscInt _j,_idx=idx_start;\
847:   for (_j=0; _j<nlnk; _j++){\
848:     _idx = lnk[_idx];\
849:     indices[_j] = _idx;\
850:     PetscBTClear(bt,_idx);\
851:   }\
852:   lnk[idx_start] = lnk_max;\
853: }
854: /*
855:   Free memories used by the list
856: */
857: #define PetscLLDestroy(lnk,bt) (PetscFree(lnk) || PetscBTDestroy(&(bt)))

859: /* Routines below are used for incomplete matrix factorization */
860: /*
861:   Create and initialize a linked list and its levels
862:   Input Parameters:
863:     idx_start - starting index of the list
864:     lnk_max   - max value of lnk indicating the end of the list
865:     nlnk      - max length of the list
866:   Output Parameters:
867:     lnk       - list initialized
868:     lnk_lvl   - array of size nlnk for storing levels of lnk
869:     bt        - PetscBT (bitarray) with all bits set to false
870: */
871: #define PetscIncompleteLLCreate(idx_start,lnk_max,nlnk,lnk,lnk_lvl,bt)\
872:   (PetscMalloc(2*nlnk*sizeof(PetscInt),&lnk) || PetscBTCreate(nlnk,&(bt)) || (lnk[idx_start] = lnk_max,lnk_lvl = lnk + nlnk,0))

874: /*
875:   Initialize a sorted linked list used for ILU and ICC
876:   Input Parameters:
877:     nidx      - number of input idx
878:     idx       - interger array used for storing column indices
879:     idx_start - starting index of the list
880:     perm      - indices of an IS
881:     lnk       - linked list(an integer array) that is created
882:     lnklvl    - levels of lnk
883:     bt        - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
884:   output Parameters:
885:     nlnk     - number of newly added idx
886:     lnk      - the sorted(increasing order) linked list containing new and non-redundate entries from idx
887:     lnklvl   - levels of lnk
888:     bt       - updated PetscBT (bitarray)
889: */
890: #define PetscIncompleteLLInit(nidx,idx,idx_start,perm,nlnk,lnk,lnklvl,bt) 0;\
891: {\
892:   PetscInt _k,_entry,_location,_lnkdata;\
893:   nlnk     = 0;\
894:   _lnkdata = idx_start;\
895:   for (_k=0; _k<nidx; _k++){\
896:     _entry = perm[idx[_k]];\
897:     if (!PetscBTLookupSet(bt,_entry)){  /* new entry */\
898:       /* search for insertion location */\
899:       if (_k && _entry < _lnkdata) _lnkdata  = idx_start;\
900:       do {\
901:         _location = _lnkdata;\
902:         _lnkdata  = lnk[_location];\
903:       } while (_entry > _lnkdata);\
904:       /* insertion location is found, add entry into lnk */\
905:       lnk[_location]  = _entry;\
906:       lnk[_entry]     = _lnkdata;\
907:       lnklvl[_entry] = 0;\
908:       nlnk++;\
909:       _lnkdata = _entry; /* next search starts from here if next_entry > _entry */\
910:     }\
911:   }\
912: }

914: /*
915:   Add a SORTED index set into a sorted linked list for ILU
916:   Input Parameters:
917:     nidx      - number of input indices
918:     idx       - sorted interger array used for storing column indices
919:     level     - level of fill, e.g., ICC(level)
920:     idxlvl    - level of idx
921:     idx_start - starting index of the list
922:     lnk       - linked list(an integer array) that is created
923:     lnklvl    - levels of lnk
924:     bt        - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
925:     prow      - the row number of idx
926:   output Parameters:
927:     nlnk     - number of newly added idx
928:     lnk      - the sorted(increasing order) linked list containing new and non-redundate entries from idx
929:     lnklvl   - levels of lnk
930:     bt       - updated PetscBT (bitarray)

932:   Note: the level of factor(i,j) is set as lvl(i,j) = min{ lvl(i,j), lvl(i,prow)+lvl(prow,j)+1)
933:         where idx = non-zero columns of U(prow,prow+1:n-1), prow<i
934: */
935: #define PetscILULLAddSorted(nidx,idx,level,idxlvl,idx_start,nlnk,lnk,lnklvl,bt,lnklvl_prow) 0;\
936: {\
937:   PetscInt _k,_entry,_location,_lnkdata,_incrlev,_lnklvl_prow=lnklvl[prow];\
938:   nlnk     = 0;\
939:   _lnkdata = idx_start;\
940:   for (_k=0; _k<nidx; _k++){\
941:     _incrlev = idxlvl[_k] + _lnklvl_prow + 1;\
942:     if (_incrlev > level) continue;\
943:     _entry = idx[_k];\
944:     if (!PetscBTLookupSet(bt,_entry)){  /* new entry */\
945:       /* search for insertion location */\
946:       do {\
947:         _location = _lnkdata;\
948:         _lnkdata  = lnk[_location];\
949:       } while (_entry > _lnkdata);\
950:       /* insertion location is found, add entry into lnk */\
951:       lnk[_location]  = _entry;\
952:       lnk[_entry]     = _lnkdata;\
953:       lnklvl[_entry] = _incrlev;\
954:       nlnk++;\
955:       _lnkdata = _entry; /* next search starts from here if next_entry > _entry */\
956:     } else { /* existing entry: update lnklvl */\
957:       if (lnklvl[_entry] > _incrlev) lnklvl[_entry] = _incrlev;\
958:     }\
959:   }\
960: }

962: /*
963:   Add a index set into a sorted linked list
964:   Input Parameters:
965:     nidx      - number of input idx
966:     idx   - interger array used for storing column indices
967:     level     - level of fill, e.g., ICC(level)
968:     idxlvl - level of idx
969:     idx_start - starting index of the list
970:     lnk       - linked list(an integer array) that is created
971:     lnklvl   - levels of lnk
972:     bt        - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
973:   output Parameters:
974:     nlnk      - number of newly added idx
975:     lnk       - the sorted(increasing order) linked list containing new and non-redundate entries from idx
976:     lnklvl   - levels of lnk
977:     bt        - updated PetscBT (bitarray)
978: */
979: #define PetscIncompleteLLAdd(nidx,idx,level,idxlvl,idx_start,nlnk,lnk,lnklvl,bt) 0;\
980: {\
981:   PetscInt _k,_entry,_location,_lnkdata,_incrlev;\
982:   nlnk     = 0;\
983:   _lnkdata = idx_start;\
984:   for (_k=0; _k<nidx; _k++){\
985:     _incrlev = idxlvl[_k] + 1;\
986:     if (_incrlev > level) continue;\
987:     _entry = idx[_k];\
988:     if (!PetscBTLookupSet(bt,_entry)){  /* new entry */\
989:       /* search for insertion location */\
990:       if (_k && _entry < _lnkdata) _lnkdata  = idx_start;\
991:       do {\
992:         _location = _lnkdata;\
993:         _lnkdata  = lnk[_location];\
994:       } while (_entry > _lnkdata);\
995:       /* insertion location is found, add entry into lnk */\
996:       lnk[_location]  = _entry;\
997:       lnk[_entry]     = _lnkdata;\
998:       lnklvl[_entry] = _incrlev;\
999:       nlnk++;\
1000:       _lnkdata = _entry; /* next search starts from here if next_entry > _entry */\
1001:     } else { /* existing entry: update lnklvl */\
1002:       if (lnklvl[_entry] > _incrlev) lnklvl[_entry] = _incrlev;\
1003:     }\
1004:   }\
1005: }

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

1051: /*
1052:   Add a SORTED index set into a sorted linked list for ICC
1053:   Input Parameters:
1054:     nidx      - number of input indices
1055:     idx       - sorted interger array used for storing column indices
1056:     level     - level of fill, e.g., ICC(level)
1057:     idxlvl    - level of idx
1058:     idx_start - starting index of the list
1059:     lnk       - linked list(an integer array) that is created
1060:     lnklvl    - levels of lnk
1061:     bt        - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
1062:     idxlvl_prow - idxlvl[prow], where prow is the row number of the idx
1063:   output Parameters:
1064:     nlnk   - number of newly added indices
1065:     lnk    - the sorted(increasing order) linked list containing new and non-redundate entries from idx
1066:     lnklvl - levels of lnk
1067:     bt     - updated PetscBT (bitarray)
1068:   Note: the level of U(i,j) is set as lvl(i,j) = min{ lvl(i,j), lvl(prow,i)+lvl(prow,j)+1)
1069:         where idx = non-zero columns of U(prow,prow+1:n-1), prow<i
1070: */
1071: #define PetscICCLLAddSorted(nidx,idx,level,idxlvl,idx_start,nlnk,lnk,lnklvl,bt,idxlvl_prow) 0;\
1072: {\
1073:   PetscInt _k,_entry,_location,_lnkdata,_incrlev;\
1074:   nlnk = 0;\
1075:   _lnkdata = idx_start;\
1076:   for (_k=0; _k<nidx; _k++){\
1077:     _incrlev = idxlvl[_k] + idxlvl_prow + 1;\
1078:     if (_incrlev > level) continue;\
1079:     _entry = idx[_k];\
1080:     if (!PetscBTLookupSet(bt,_entry)){  /* new entry */\
1081:       /* search for insertion location */\
1082:       do {\
1083:         _location = _lnkdata;\
1084:         _lnkdata  = lnk[_location];\
1085:       } while (_entry > _lnkdata);\
1086:       /* insertion location is found, add entry into lnk */\
1087:       lnk[_location] = _entry;\
1088:       lnk[_entry]    = _lnkdata;\
1089:       lnklvl[_entry] = _incrlev;\
1090:       nlnk++;\
1091:       _lnkdata = _entry; /* next search starts from here */\
1092:     } else { /* existing entry: update lnklvl */\
1093:       if (lnklvl[_entry] > _incrlev) lnklvl[_entry] = _incrlev;\
1094:     }\
1095:   }\
1096: }

1098: /*
1099:   Copy data on the list into an array, then initialize the list
1100:   Input Parameters:
1101:     idx_start - starting index of the list
1102:     lnk_max   - max value of lnk indicating the end of the list
1103:     nlnk      - number of data on the list to be copied
1104:     lnk       - linked list
1105:     lnklvl    - level of lnk
1106:     bt        - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
1107:   output Parameters:
1108:     indices - array that contains the copied data
1109:     lnk     - linked list that is cleaned and initialize
1110:     lnklvl  - level of lnk that is reinitialized
1111:     bt      - PetscBT (bitarray) with all bits set to false
1112: */
1113: #define PetscIncompleteLLClean(idx_start,lnk_max,nlnk,lnk,lnklvl,indices,indiceslvl,bt) 0;\
1114: {\
1115:   PetscInt _j,_idx=idx_start;\
1116:   for (_j=0; _j<nlnk; _j++){\
1117:     _idx = lnk[_idx];\
1118:     *(indices+_j) = _idx;\
1119:     *(indiceslvl+_j) = lnklvl[_idx];\
1120:     lnklvl[_idx] = -1;\
1121:     PetscBTClear(bt,_idx);\
1122:   }\
1123:   lnk[idx_start] = lnk_max;\
1124: }
1125: /*
1126:   Free memories used by the list
1127: */
1128: #define PetscIncompleteLLDestroy(lnk,bt) (PetscFree(lnk) || PetscBTDestroy(&(bt)))

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

1141:       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
1142:       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
1143:       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
1144:       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.
1145:       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
1146:       to each other so memory access is much better than using the big array.

1148:   Example:
1149:      nlnk_max=5, lnk_max=36:
1150:      Initial list: [0, 0 | 36, 2 | 0, 0 | 0, 0 | 0, 0 | 0, 0 | 0, 0]
1151:      here, head_node has index 2 with value lnk[2]=lnk_max=36,
1152:            0-th entry is used to store the number of entries in the list,
1153:      The initial lnk represents head -> tail(marked by 36) with number of entries = lnk[0]=0.

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

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

1163:   Input Parameters:
1164:     nlnk_max  - max length of the list
1165:     lnk_max   - max value of the entries
1166:   Output Parameters:
1167:     lnk       - list created and initialized
1168:     bt        - PetscBT (bitarray) with all bits set to false. Note: bt has size lnk_max, not nln_max!
1169: */
1170: PETSC_STATIC_INLINE PetscErrorCode PetscLLCondensedCreate(PetscInt nlnk_max,PetscInt lnk_max,PetscInt **lnk,PetscBT *bt)
1171: {
1173:   PetscInt       *llnk;

1176:   PetscMalloc(2*(nlnk_max+2)*sizeof(PetscInt),lnk);
1177:   PetscBTCreate(lnk_max,bt);
1178:   llnk = *lnk;
1179:   llnk[0] = 0;         /* number of entries on the list */
1180:   llnk[2] = lnk_max;   /* value in the head node */
1181:   llnk[3] = 2;         /* next for the head node */
1182:   return(0);
1183: }

1187: /*
1188:   Add a SORTED ascending index set into a sorted linked list. See PetscLLCondensedCreate() for detailed description.
1189:   Input Parameters:
1190:     nidx      - number of input indices
1191:     indices   - sorted interger array
1192:     lnk       - condensed linked list(an integer array) that is created
1193:     bt        - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
1194:   output Parameters:
1195:     lnk       - the sorted(increasing order) linked list containing previous and newly added non-redundate indices
1196:     bt        - updated PetscBT (bitarray)
1197: */
1198: PETSC_STATIC_INLINE PetscErrorCode PetscLLCondensedAddSorted(PetscInt nidx,const PetscInt indices[],PetscInt lnk[],PetscBT bt)
1199: {
1200:   PetscInt _k,_entry,_location,_next,_lnkdata,_nlnk,_newnode;

1203:   _nlnk     = lnk[0]; /* num of entries on the input lnk */
1204:   _location = 2; /* head */
1205:     for (_k=0; _k<nidx; _k++){
1206:       _entry = indices[_k];
1207:       if (!PetscBTLookupSet(bt,_entry)){  /* new entry */
1208:         /* search for insertion location */
1209:         do {
1210:           _next     = _location + 1; /* link from previous node to next node */
1211:           _location = lnk[_next];    /* idx of next node */
1212:           _lnkdata  = lnk[_location];/* value of next node */
1213:         } while (_entry > _lnkdata);
1214:         /* insertion location is found, add entry into lnk */
1215:         _newnode        = 2*(_nlnk+2);   /* index for this new node */
1216:         lnk[_next]      = _newnode;      /* connect previous node to the new node */
1217:         lnk[_newnode]   = _entry;        /* set value of the new node */
1218:         lnk[_newnode+1] = _location;     /* connect new node to next node */
1219:         _location       = _newnode;      /* next search starts from the new node */
1220:         _nlnk++;
1221:       }   \
1222:     }\
1223:   lnk[0]   = _nlnk;   /* number of entries in the list */
1224:   return(0);
1225: }

1229: PETSC_STATIC_INLINE PetscErrorCode PetscLLCondensedClean(PetscInt lnk_max,PetscInt nidx,PetscInt *indices,PetscInt lnk[],PetscBT bt)
1230: {
1232:   PetscInt       _k,_next,_nlnk;

1235:   _next = lnk[3];       /* head node */
1236:   _nlnk = lnk[0];       /* num of entries on the list */
1237:   for (_k=0; _k<_nlnk; _k++){
1238:     indices[_k] = lnk[_next];
1239:     _next       = lnk[_next + 1];
1240:     PetscBTClear(bt,indices[_k]);
1241:   }
1242:   lnk[0] = 0;          /* num of entries on the list */
1243:   lnk[2] = lnk_max;    /* initialize head node */
1244:   lnk[3] = 2;          /* head node */
1245:   return(0);
1246: }

1250: PETSC_STATIC_INLINE PetscErrorCode PetscLLCondensedView(PetscInt *lnk)
1251: {
1253:   PetscInt       k;

1256:   PetscPrintf(PETSC_COMM_SELF,"LLCondensed of size %d, (val,  next)\n",lnk[0]);
1257:   for (k=2; k< lnk[0]+2; k++){
1258:     PetscPrintf(PETSC_COMM_SELF," %D: (%D, %D)\n",2*k,lnk[2*k],lnk[2*k+1]);
1259:   }
1260:   return(0);
1261: }

1265: /*
1266:   Free memories used by the list
1267: */
1268: PETSC_STATIC_INLINE PetscErrorCode PetscLLCondensedDestroy(PetscInt *lnk,PetscBT bt)
1269: {

1273:   PetscFree(lnk);
1274:   PetscBTDestroy(&bt);
1275:   return(0);
1276: }

1278: /* -------------------------------------------------------------------------------------------------------*/
1281: /*
1282:  Same as PetscLLCondensedCreate(), but does not use non-scalable O(lnk_max) bitarray
1283:   Input Parameters:
1284:     nlnk_max  - max length of the list
1285:   Output Parameters:
1286:     lnk       - list created and initialized
1287: */
1288: PETSC_STATIC_INLINE PetscErrorCode PetscLLCondensedCreate_Scalable(PetscInt nlnk_max,PetscInt **lnk)
1289: {
1291:   PetscInt       *llnk;

1294:   PetscMalloc(2*(nlnk_max+2)*sizeof(PetscInt),lnk);
1295:   llnk = *lnk;
1296:   llnk[0] = 0;               /* number of entries on the list */
1297:   llnk[2] = PETSC_MAX_INT;   /* value in the head node */
1298:   llnk[3] = 2;               /* next for the head node */
1299:   return(0);
1300: }

1304: PETSC_STATIC_INLINE PetscErrorCode PetscLLCondensedAddSorted_Scalable(PetscInt nidx,const PetscInt indices[],PetscInt lnk[])
1305: {
1306:   PetscInt _k,_entry,_location,_next,_lnkdata,_nlnk,_newnode;
1307:   _nlnk     = lnk[0]; /* num of entries on the input lnk */
1308:   _location = 2; /* head */ \
1309:     for (_k=0; _k<nidx; _k++){
1310:       _entry = indices[_k];
1311:       /* search for insertion location */
1312:       do {
1313:         _next     = _location + 1; /* link from previous node to next node */
1314:         _location = lnk[_next];    /* idx of next node */
1315:         _lnkdata  = lnk[_location];/* value of next node */
1316:       } while (_entry > _lnkdata);
1317:       if (_entry < _lnkdata) {
1318:         /* insertion location is found, add entry into lnk */
1319:         _newnode        = 2*(_nlnk+2);   /* index for this new node */
1320:         lnk[_next]      = _newnode;      /* connect previous node to the new node */
1321:         lnk[_newnode]   = _entry;        /* set value of the new node */
1322:         lnk[_newnode+1] = _location;     /* connect new node to next node */
1323:         _location       = _newnode;      /* next search starts from the new node */
1324:         _nlnk++;
1325:       }
1326:     }
1327:   lnk[0]   = _nlnk;   /* number of entries in the list */
1328:   return 0;
1329: }

1333: PETSC_STATIC_INLINE PetscErrorCode PetscLLCondensedClean_Scalable(PetscInt nidx,PetscInt *indices,PetscInt *lnk)
1334: {
1335:   PetscInt _k,_next,_nlnk;
1336:   _next = lnk[3];       /* head node */
1337:   _nlnk = lnk[0];
1338:   for (_k=0; _k<_nlnk; _k++){
1339:     indices[_k] = lnk[_next];
1340:     _next       = lnk[_next + 1];
1341:   }
1342:   lnk[0] = 0;          /* num of entries on the list */
1343:   lnk[3] = 2;          /* head node */
1344:   return 0;
1345: }

1349: PETSC_STATIC_INLINE PetscErrorCode PetscLLCondensedDestroy_Scalable(PetscInt *lnk)
1350: {
1351:   return PetscFree(lnk);
1352: }

1354: /* -------------------------------------------------------------------------------------------------------*/
1355: /*
1356:       lnk[0]   number of links
1357:       lnk[1]   number of entries
1358:       lnk[3n]  value
1359:       lnk[3n+1] len
1360:       lnk[3n+2] link to next value

1362:       The next three are always the first link

1364:       lnk[3]    PETSC_MIN_INT+1
1365:       lnk[4]    1
1366:       lnk[5]    link to first real entry

1368:       The next three are always the last link

1370:       lnk[6]    PETSC_MAX_INT - 1
1371:       lnk[7]    1
1372:       lnk[8]    next valid link (this is the same as lnk[0] but without the decreases)
1373: */

1377: PETSC_STATIC_INLINE PetscErrorCode PetscLLCondensedCreate_fast(PetscInt nlnk_max,PetscInt **lnk)
1378: {
1380:   PetscInt       *llnk;

1383:   PetscMalloc(3*(nlnk_max+3)*sizeof(PetscInt),lnk);
1384:   llnk = *lnk;
1385:   llnk[0] = 0;   /* nlnk: number of entries on the list */
1386:   llnk[1] = 0;          /* number of integer entries represented in list */
1387:   llnk[3] = PETSC_MIN_INT+1;   /* value in the first node */
1388:   llnk[4] = 1;           /* count for the first node */
1389:   llnk[5] = 6;         /* next for the first node */
1390:   llnk[6] = PETSC_MAX_INT-1;   /* value in the last node */
1391:   llnk[7] = 1;           /* count for the last node */
1392:   llnk[8] = 0;         /* next valid node to be used */
1393:   return(0);
1394: }

1396: PETSC_STATIC_INLINE PetscErrorCode PetscLLCondensedAddSorted_fast(PetscInt nidx,const PetscInt indices[],PetscInt lnk[])
1397: {
1398:   PetscInt k,entry,prev,next;
1399:   prev      = 3;      /* first value */
1400:   next      = lnk[prev+2];
1401:   for (k=0; k<nidx; k++){
1402:     entry = indices[k];
1403:     /* search for insertion location */
1404:     while (entry >= lnk[next]) {
1405:       prev = next;
1406:       next = lnk[next+2];
1407:     }
1408:     /* entry is in range of previous list */
1409:     if (entry < lnk[prev]+lnk[prev+1]) continue;
1410:     lnk[1]++;
1411:     /* entry is right after previous list */
1412:     if (entry == lnk[prev]+lnk[prev+1]) {
1413:       lnk[prev+1]++;
1414:       if (lnk[next] == entry+1) { /* combine two contiquous strings */
1415:         lnk[prev+1] += lnk[next+1];
1416:         lnk[prev+2]  = lnk[next+2];
1417:         next         = lnk[next+2];
1418:         lnk[0]--;
1419:       }
1420:       continue;
1421:     }
1422:     /* entry is right before next list */
1423:     if (entry == lnk[next]-1) {
1424:       lnk[next]--;
1425:       lnk[next+1]++;
1426:       prev = next;
1427:       next = lnk[prev+2];
1428:       continue;
1429:     }
1430:     /*  add entry into lnk */
1431:     lnk[prev+2]    = 3*((lnk[8]++)+3);      /* connect previous node to the new node */
1432:     prev           = lnk[prev+2];
1433:     lnk[prev]      = entry;        /* set value of the new node */
1434:     lnk[prev+1]    = 1;             /* number of values in contiquous string is one to start */
1435:     lnk[prev+2]    = next;          /* connect new node to next node */
1436:     lnk[0]++;
1437:   }
1438:   return 0;
1439: }

1441: PETSC_STATIC_INLINE PetscErrorCode PetscLLCondensedClean_fast(PetscInt nidx,PetscInt *indices,PetscInt *lnk)
1442: {
1443:   PetscInt _k,_next,_nlnk,cnt,j;
1444:   _next = lnk[5];       /* first node */
1445:   _nlnk = lnk[0];
1446:   cnt   = 0;
1447:   for (_k=0; _k<_nlnk; _k++){
1448:     for (j=0; j<lnk[_next+1]; j++) {
1449:       indices[cnt++] = lnk[_next] + j;
1450:     }
1451:     _next       = lnk[_next + 2];
1452:   }
1453:   lnk[0] = 0;   /* nlnk: number of links */
1454:   lnk[1] = 0;          /* number of integer entries represented in list */
1455:   lnk[3] = PETSC_MIN_INT+1;   /* value in the first node */
1456:   lnk[4] = 1;           /* count for the first node */
1457:   lnk[5] = 6;         /* next for the first node */
1458:   lnk[6] = PETSC_MAX_INT-1;   /* value in the last node */
1459:   lnk[7] = 1;           /* count for the last node */
1460:   lnk[8] = 0;         /* next valid location to make link */
1461:   return 0;
1462: }

1464: PETSC_STATIC_INLINE PetscErrorCode PetscLLCondensedView_fast(PetscInt *lnk)
1465: {
1466:   PetscInt k,next,nlnk;
1467:   next = lnk[5];       /* first node */
1468:   nlnk = lnk[0];
1469:   for (k=0; k<nlnk; k++){
1470: #if 0                           /* Debugging code */
1471:     printf("%d value %d len %d next %d\n",next,lnk[next],lnk[next+1],lnk[next+2]);
1472: #endif
1473:     next = lnk[next + 2];
1474:   }
1475:   return 0;
1476: }

1478: PETSC_STATIC_INLINE PetscErrorCode PetscLLCondensedDestroy_fast(PetscInt *lnk)
1479: {
1480:   return PetscFree(lnk);
1481: }

1483: PETSC_EXTERN PetscLogEvent MAT_Mult, MAT_MultMatrixFree, MAT_Mults, MAT_MultConstrained, MAT_MultAdd, MAT_MultTranspose;
1484: PETSC_EXTERN PetscLogEvent MAT_MultTransposeConstrained, MAT_MultTransposeAdd, MAT_Solve, MAT_Solves, MAT_SolveAdd, MAT_SolveTranspose;
1485: PETSC_EXTERN PetscLogEvent MAT_SolveTransposeAdd, MAT_SOR, MAT_ForwardSolve, MAT_BackwardSolve, MAT_LUFactor, MAT_LUFactorSymbolic;
1486: PETSC_EXTERN PetscLogEvent MAT_LUFactorNumeric, MAT_CholeskyFactor, MAT_CholeskyFactorSymbolic, MAT_CholeskyFactorNumeric, MAT_ILUFactor;
1487: PETSC_EXTERN PetscLogEvent MAT_ILUFactorSymbolic, MAT_ICCFactorSymbolic, MAT_Copy, MAT_Convert, MAT_Scale, MAT_AssemblyBegin;
1488: PETSC_EXTERN PetscLogEvent MAT_AssemblyEnd, MAT_SetValues, MAT_GetValues, MAT_GetRow, MAT_GetRowIJ, MAT_GetSubMatrices, MAT_GetColoring, MAT_GetOrdering, MAT_GetRedundantMatrix;
1489: PETSC_EXTERN PetscLogEvent MAT_IncreaseOverlap, MAT_Partitioning, MAT_Coarsen, MAT_ZeroEntries, MAT_Load, MAT_View, MAT_AXPY, MAT_FDColoringCreate, MAT_TransposeColoringCreate;
1490: PETSC_EXTERN PetscLogEvent MAT_FDColoringApply, MAT_Transpose, MAT_FDColoringFunction;
1491: PETSC_EXTERN PetscLogEvent MAT_MatMult, MAT_MatSolve,MAT_MatMultSymbolic, MAT_MatMultNumeric,MAT_Getlocalmatcondensed,MAT_GetBrowsOfAcols,MAT_GetBrowsOfAocols;
1492: PETSC_EXTERN PetscLogEvent MAT_PtAP, MAT_PtAPSymbolic, MAT_PtAPNumeric,MAT_Seqstompinum,MAT_Seqstompisym,MAT_Seqstompi,MAT_Getlocalmat;
1493: PETSC_EXTERN PetscLogEvent MAT_RARt, MAT_RARtSymbolic, MAT_RARtNumeric;
1494: PETSC_EXTERN PetscLogEvent MAT_MatTransposeMult, MAT_MatTransposeMultSymbolic, MAT_MatTransposeMultNumeric;
1495: PETSC_EXTERN PetscLogEvent MAT_TransposeMatMult, MAT_TransposeMatMultSymbolic, MAT_TransposeMatMultNumeric;
1496: PETSC_EXTERN PetscLogEvent MAT_MatMatMult, MAT_MatMatMultSymbolic, MAT_MatMatMultNumeric;
1497: PETSC_EXTERN PetscLogEvent MAT_Applypapt, MAT_Applypapt_symbolic, MAT_Applypapt_numeric;
1498: PETSC_EXTERN PetscLogEvent MAT_Getsymtranspose, MAT_Transpose_SeqAIJ, MAT_Getsymtransreduced,MAT_GetSequentialNonzeroStructure;

1500: PETSC_EXTERN PetscLogEvent MATMFFD_Mult;
1501: PETSC_EXTERN PetscLogEvent MAT_GetMultiProcBlock;
1502: PETSC_EXTERN PetscLogEvent MAT_CUSPCopyToGPU, MAT_CUSPARSECopyToGPU, MAT_SetValuesBatch, MAT_SetValuesBatchI, MAT_SetValuesBatchII, MAT_SetValuesBatchIII, MAT_SetValuesBatchIV;
1503: PETSC_EXTERN PetscLogEvent MAT_Merge;

1505: #endif