Actual source code: petsc-matimpl.h
petsc-dev 2014-02-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 (*placeholder_32)(Mat);
59: PetscErrorCode (*placeholder_33)(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 (*placeholder_73)(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 (*placeholder_98)(Mat);
138: /*99*/
139: PetscErrorCode (*placeholder_99)(Mat);
140: PetscErrorCode (*placeholder_100)(Mat);
141: PetscErrorCode (*placeholder_101)(Mat);
142: PetscErrorCode (*conjugate)(Mat); /* complex conjugate */
143: PetscErrorCode (*placeholder_103)(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,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 (*placeholder_127)(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: PetscErrorCode (*residual)(Mat,Vec,Vec,Vec);
190: PetscErrorCode (*fdcoloringsetup)(Mat,ISColoring,MatFDColoring);
191: };
192: /*
193: If you add MatOps entries above also add them to the MATOP enum
194: in include/petscmat.h and include/finclude/petscmat.h
195: */
197: #include <petscsys.h>
198: PETSC_EXTERN PetscErrorCode MatRegisterOp(MPI_Comm, const char[], PetscVoidFunction, const char[], PetscInt, ...);
199: PETSC_EXTERN PetscErrorCode MatQueryOp(MPI_Comm, PetscVoidFunction*, const char[], PetscInt, ...);
201: typedef struct _p_MatBaseName* MatBaseName;
202: struct _p_MatBaseName {
203: char *bname,*sname,*mname;
204: MatBaseName next;
205: };
207: PETSC_EXTERN MatBaseName MatBaseNameList;
209: /*
210: Utility private matrix routines
211: */
212: PETSC_INTERN PetscErrorCode MatConvert_Basic(Mat, MatType,MatReuse,Mat*);
213: PETSC_INTERN PetscErrorCode MatCopy_Basic(Mat,Mat,MatStructure);
214: PETSC_INTERN PetscErrorCode MatHeaderMerge(Mat,Mat);
215: PETSC_EXTERN PetscErrorCode MatHeaderReplace(Mat,Mat);
216: PETSC_INTERN PetscErrorCode MatAXPYGetxtoy_Private(PetscInt,PetscInt*,PetscInt*,PetscInt*, PetscInt*,PetscInt*,PetscInt*, PetscInt**);
217: PETSC_INTERN PetscErrorCode MatDiagonalSet_Default(Mat,Vec,InsertMode);
219: #if defined(PETSC_USE_DEBUG)
220: # define MatCheckPreallocated(A,arg) do { \
221: 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); \
222: } while (0)
223: #else
224: # define MatCheckPreallocated(A,arg) do {} while (0)
225: #endif
227: /*
228: The stash is used to temporarily store inserted matrix values that
229: belong to another processor. During the assembly phase the stashed
230: values are moved to the correct processor and
231: */
233: typedef struct _MatStashSpace *PetscMatStashSpace;
235: struct _MatStashSpace {
236: PetscMatStashSpace next;
237: PetscScalar *space_head,*val;
238: PetscInt *idx,*idy;
239: PetscInt total_space_size;
240: PetscInt local_used;
241: PetscInt local_remaining;
242: };
244: PETSC_EXTERN PetscErrorCode PetscMatStashSpaceGet(PetscInt,PetscInt,PetscMatStashSpace *);
245: PETSC_EXTERN PetscErrorCode PetscMatStashSpaceContiguous(PetscInt,PetscMatStashSpace *,PetscScalar *,PetscInt *,PetscInt *);
246: PETSC_EXTERN PetscErrorCode PetscMatStashSpaceDestroy(PetscMatStashSpace*);
248: typedef struct {
249: PetscInt nmax; /* maximum stash size */
250: PetscInt umax; /* user specified max-size */
251: PetscInt oldnmax; /* the nmax value used previously */
252: PetscInt n; /* stash size */
253: PetscInt bs; /* block size of the stash */
254: PetscInt reallocs; /* preserve the no of mallocs invoked */
255: PetscMatStashSpace space_head,space; /* linked list to hold stashed global row/column numbers and matrix values */
256: /* The following variables are used for communication */
257: MPI_Comm comm;
258: PetscMPIInt size,rank;
259: PetscMPIInt tag1,tag2;
260: MPI_Request *send_waits; /* array of send requests */
261: MPI_Request *recv_waits; /* array of receive requests */
262: MPI_Status *send_status; /* array of send status */
263: PetscInt nsends,nrecvs; /* numbers of sends and receives */
264: PetscScalar *svalues; /* sending data */
265: PetscInt *sindices;
266: PetscScalar **rvalues; /* receiving data (values) */
267: PetscInt **rindices; /* receiving data (indices) */
268: PetscInt nprocessed; /* number of messages already processed */
269: PetscMPIInt *flg_v; /* indicates what messages have arrived so far and from whom */
270: PetscBool reproduce;
271: PetscInt reproduce_count;
272: } MatStash;
274: PETSC_INTERN PetscErrorCode MatStashCreate_Private(MPI_Comm,PetscInt,MatStash*);
275: PETSC_INTERN PetscErrorCode MatStashDestroy_Private(MatStash*);
276: PETSC_INTERN PetscErrorCode MatStashScatterEnd_Private(MatStash*);
277: PETSC_INTERN PetscErrorCode MatStashSetInitialSize_Private(MatStash*,PetscInt);
278: PETSC_INTERN PetscErrorCode MatStashGetInfo_Private(MatStash*,PetscInt*,PetscInt*);
279: PETSC_INTERN PetscErrorCode MatStashValuesRow_Private(MatStash*,PetscInt,PetscInt,const PetscInt[],const PetscScalar[],PetscBool );
280: PETSC_INTERN PetscErrorCode MatStashValuesCol_Private(MatStash*,PetscInt,PetscInt,const PetscInt[],const PetscScalar[],PetscInt,PetscBool );
281: PETSC_INTERN PetscErrorCode MatStashValuesRowBlocked_Private(MatStash*,PetscInt,PetscInt,const PetscInt[],const PetscScalar[],PetscInt,PetscInt,PetscInt);
282: PETSC_INTERN PetscErrorCode MatStashValuesColBlocked_Private(MatStash*,PetscInt,PetscInt,const PetscInt[],const PetscScalar[],PetscInt,PetscInt,PetscInt);
283: PETSC_INTERN PetscErrorCode MatStashScatterBegin_Private(Mat,MatStash*,PetscInt*);
284: PETSC_INTERN PetscErrorCode MatStashScatterGetMesg_Private(MatStash*,PetscMPIInt*,PetscInt**,PetscInt**,PetscScalar**,PetscInt*);
286: typedef struct {
287: PetscInt dim;
288: PetscInt dims[4];
289: PetscInt starts[4];
290: PetscBool noc; /* this is a single component problem, hence user will not set MatStencil.c */
291: } MatStencilInfo;
293: /* Info about using compressed row format */
294: typedef struct {
295: PetscBool use; /* indicates compressed rows have been checked and will be used */
296: PetscInt nrows; /* number of non-zero rows */
297: PetscInt *i; /* compressed row pointer */
298: PetscInt *rindex; /* compressed row index */
299: } Mat_CompressedRow;
300: PETSC_EXTERN PetscErrorCode MatCheckCompressedRow(Mat,PetscInt,Mat_CompressedRow*,PetscInt*,PetscInt,PetscReal);
302: struct _p_Mat {
303: PETSCHEADER(struct _MatOps);
304: PetscLayout rmap,cmap;
305: void *data; /* implementation-specific data */
306: MatFactorType factortype; /* MAT_FACTOR_LU, ILU, CHOLESKY or ICC */
307: PetscBool assembled; /* is the matrix assembled? */
308: PetscBool was_assembled; /* new values inserted into assembled mat */
309: PetscInt num_ass; /* number of times matrix has been assembled */
310: PetscBool same_nonzero; /* matrix has same nonzero pattern as previous */
311: MatInfo info; /* matrix information */
312: InsertMode insertmode; /* have values been inserted in matrix or added? */
313: MatStash stash,bstash; /* used for assembling off-proc mat emements */
314: MatNullSpace nullsp; /* null space (operator is singular) */
315: MatNullSpace nearnullsp; /* near null space to be used by multigrid methods */
316: PetscBool preallocated;
317: MatStencilInfo stencil; /* information for structured grid */
318: PetscBool symmetric,hermitian,structurally_symmetric,spd;
319: PetscBool symmetric_set,hermitian_set,structurally_symmetric_set,spd_set; /* if true, then corresponding flag is correct*/
320: PetscBool symmetric_eternal;
321: PetscBool nooffprocentries,nooffproczerorows;
322: #if defined(PETSC_HAVE_CUSP)
323: PetscCUSPFlag valid_GPU_matrix; /* flag pointing to the matrix on the gpu*/
324: #endif
325: #if defined(PETSC_HAVE_VIENNACL)
326: PetscViennaCLFlag valid_GPU_matrix; /* flag pointing to the matrix on the gpu*/
327: #endif
328: void *spptr; /* pointer for special library like SuperLU */
329: MatSolverPackage solvertype;
330: PetscBool checksymmetryonassembly,checknullspaceonassembly;
331: PetscReal checksymmetrytol;
332: };
334: PETSC_INTERN PetscErrorCode MatAXPY_Basic(Mat,PetscScalar,Mat,MatStructure);
335: PETSC_INTERN PetscErrorCode MatAXPY_BasicWithPreallocation(Mat,Mat,PetscScalar,Mat,MatStructure);
336: /*
337: Object for partitioning graphs
338: */
340: typedef struct _MatPartitioningOps *MatPartitioningOps;
341: struct _MatPartitioningOps {
342: PetscErrorCode (*apply)(MatPartitioning,IS*);
343: PetscErrorCode (*setfromoptions)(MatPartitioning);
344: PetscErrorCode (*destroy)(MatPartitioning);
345: PetscErrorCode (*view)(MatPartitioning,PetscViewer);
346: };
348: struct _p_MatPartitioning {
349: PETSCHEADER(struct _MatPartitioningOps);
350: Mat adj;
351: PetscInt *vertex_weights;
352: PetscReal *part_weights;
353: PetscInt n; /* number of partitions */
354: void *data;
355: PetscInt setupcalled;
356: };
358: /*
359: Object for coarsen graphs
360: */
361: typedef struct _MatCoarsenOps *MatCoarsenOps;
362: struct _MatCoarsenOps {
363: PetscErrorCode (*apply)(MatCoarsen);
364: PetscErrorCode (*setfromoptions)(MatCoarsen);
365: PetscErrorCode (*destroy)(MatCoarsen);
366: PetscErrorCode (*view)(MatCoarsen,PetscViewer);
367: };
369: struct _p_MatCoarsen {
370: PETSCHEADER(struct _MatCoarsenOps);
371: Mat graph;
372: PetscInt verbose;
373: PetscInt setupcalled;
374: void *subctx;
375: /* */
376: PetscBool strict_aggs;
377: IS perm;
378: PetscCoarsenData *agg_lists;
379: };
381: PETSC_EXTERN PetscErrorCode PetscCDCreate(PetscInt,PetscCoarsenData**);
382: PETSC_EXTERN PetscErrorCode PetscCDDestroy(PetscCoarsenData*);
383: PETSC_EXTERN PetscErrorCode PetscLLNSetID(PetscCDIntNd*,PetscInt);
384: PETSC_EXTERN PetscErrorCode PetscLLNGetID(const PetscCDIntNd*,PetscInt*);
385: PETSC_EXTERN PetscErrorCode PetscCDAppendID(PetscCoarsenData*,PetscInt,PetscInt);
386: PETSC_EXTERN PetscErrorCode PetscCDAppendRemove(PetscCoarsenData*,PetscInt,PetscInt);
387: PETSC_EXTERN PetscErrorCode PetscCDAppendNode(PetscCoarsenData*,PetscInt,PetscCDIntNd*);
388: PETSC_EXTERN PetscErrorCode PetscCDRemoveNextNode(PetscCoarsenData*,PetscInt,PetscCDIntNd*);
389: PETSC_EXTERN PetscErrorCode PetscCDRemoveAllAt(PetscCoarsenData*,PetscInt);
390: PETSC_EXTERN PetscErrorCode PetscCDSizeAt(const PetscCoarsenData*,PetscInt,PetscInt*);
391: PETSC_EXTERN PetscErrorCode PetscCDEmptyAt(const PetscCoarsenData*,PetscInt,PetscBool*);
392: PETSC_EXTERN PetscErrorCode PetscCDSetChuckSize(PetscCoarsenData*,PetscInt);
393: PETSC_EXTERN PetscErrorCode PetscCDPrint(const PetscCoarsenData*,MPI_Comm);
394: PETSC_EXTERN PetscErrorCode PetscCDGetMIS(PetscCoarsenData*,IS*);
395: PETSC_EXTERN PetscErrorCode PetscCDGetMat(const PetscCoarsenData*,Mat*);
396: PETSC_EXTERN PetscErrorCode PetscCDSetMat(PetscCoarsenData*,Mat);
398: typedef PetscCDIntNd *PetscCDPos;
399: PETSC_EXTERN PetscErrorCode PetscCDGetHeadPos(const PetscCoarsenData*,PetscInt,PetscCDPos*);
400: PETSC_EXTERN PetscErrorCode PetscCDGetNextPos(const PetscCoarsenData*,PetscInt,PetscCDPos*);
401: PETSC_EXTERN PetscErrorCode PetscCDGetASMBlocks(const PetscCoarsenData*,const PetscInt,PetscInt*,IS**);
402: /* PetscErrorCode PetscCDSetRemovedIS( PetscCoarsenData *ail, MPI_Comm, const PetscInt, PetscInt[] ); */
403: /* PetscErrorCode PetscCDGetRemovedIS( PetscCoarsenData *ail, IS * ); */
405: /*
406: MatFDColoring is used to compute Jacobian matrices efficiently
407: via coloring. The data structure is explained below in an example.
409: Color = 0 1 0 2 | 2 3 0
410: ---------------------------------------------------
411: 00 01 | 05
412: 10 11 | 14 15 Processor 0
413: 22 23 | 25
414: 32 33 |
415: ===================================================
416: | 44 45 46
417: 50 | 55 Processor 1
418: | 64 66
419: ---------------------------------------------------
421: ncolors = 4;
423: ncolumns = {2,1,1,0}
424: columns = {{0,2},{1},{3},{}}
425: nrows = {4,2,3,3}
426: rows = {{0,1,2,3},{0,1},{1,2,3},{0,1,2}}
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: vwscale = {dx(4),dx(5),dx(6)} MPI Vec
435: vscale = {dx(0),dx(4),dx(5),dx(6)} Seq Vec
437: See the routine MatFDColoringApply() for how this data is used
438: to compute the Jacobian.
440: */
441: typedef struct {
442: PetscInt row;
443: PetscInt col;
444: PetscScalar *valaddr; /* address of value */
445: } MatEntry;
447: typedef struct {
448: PetscInt row;
449: PetscScalar *valaddr; /* address of value */
450: } MatEntry2;
452: struct _p_MatFDColoring{
453: PETSCHEADER(int);
454: PetscInt M,N,m; /* total rows, columns; local rows */
455: PetscInt rstart; /* first row owned by local processor */
456: PetscInt ncolors; /* number of colors */
457: PetscInt *ncolumns; /* number of local columns for a color */
458: PetscInt **columns; /* lists the local columns of each color (using global column numbering) */
459: PetscInt *nrows; /* number of local rows for each color */
460: MatEntry *matentry; /* holds (row, column, address of value) for Jacobian matrix entry */
461: MatEntry2 *matentry2; /* holds (row, address of value) for Jacobian matrix entry */
462: PetscScalar *dy; /* store a block of F(x+dx)-F(x) when J is in BAIJ format */
463: PetscReal error_rel; /* square root of relative error in computing function */
464: PetscReal umin; /* minimum allowable u'dx value */
465: Vec w1,w2,w3; /* work vectors used in computing Jacobian */
466: PetscBool fset; /* indicates that the initial function value F(X) is set */
467: PetscErrorCode (*f)(void); /* function that defines Jacobian */
468: void *fctx; /* optional user-defined context for use by the function f */
469: Vec vscale; /* holds FD scaling, i.e. 1/dx for each perturbed column */
470: PetscInt currentcolor; /* color for which function evaluation is being done now */
471: const char *htype; /* "wp" or "ds" */
472: ISColoringType ctype; /* IS_COLORING_GLOBAL or IS_COLORING_GHOSTED */
473: PetscInt brows,bcols; /* number of block rows or columns for speedup inserting the dense matrix into sparse Jacobian */
474: PetscBool setupcalled; /* true if setup has been called */
475: void (*ftn_func_pointer)(void),*ftn_func_cntx; /* serve the same purpose as *fortran_func_pointers in PETSc objects */
476: };
478: typedef struct _MatColoringOps *MatColoringOps;
479: struct _MatColoringOps {
480: PetscErrorCode (*destroy)(MatColoring);
481: PetscErrorCode (*setfromoptions)(MatColoring);
482: PetscErrorCode (*view)(MatColoring,PetscViewer);
483: PetscErrorCode (*apply)(MatColoring,ISColoring*);
484: };
486: struct _p_MatColoring {
487: PETSCHEADER(struct _MatColoringOps);
488: Mat mat;
489: PetscInt dist; /* distance of the coloring */
490: PetscInt maxcolors; /* the maximum number of colors returned, maxcolors=1 for MIS */
491: void *data; /* inner context */
492: PetscBool valid; /* check to see if what is produced is a valid coloring */
493: };
495: struct _p_MatTransposeColoring{
496: PETSCHEADER(int);
497: PetscInt M,N,m; /* total rows, columns; local rows */
498: PetscInt rstart; /* first row owned by local processor */
499: PetscInt ncolors; /* number of colors */
500: PetscInt *ncolumns; /* number of local columns for a color */
501: PetscInt *nrows; /* number of local rows for each color */
502: PetscInt currentcolor; /* color for which function evaluation is being done now */
503: ISColoringType ctype; /* IS_COLORING_GLOBAL or IS_COLORING_GHOSTED */
505: PetscInt *colorforrow,*colorforcol; /* pointer to rows and columns */
506: PetscInt *rows; /* lists the local rows for each color (using the local row numbering) */
507: PetscInt *den2sp; /* maps (row,color) in the dense matrix to index of sparse matrix array a->a */
508: PetscInt *columns; /* lists the local columns of each color (using global column numbering) */
509: PetscInt brows; /* number of rows for efficient implementation of MatTransColoringApplyDenToSp() */
510: PetscInt *lstart; /* array used for loop over row blocks of Csparse */
511: };
513: /*
514: Null space context for preconditioner/operators
515: */
516: struct _p_MatNullSpace {
517: PETSCHEADER(int);
518: PetscBool has_cnst;
519: PetscInt n;
520: Vec* vecs;
521: PetscScalar* alpha; /* for projections */
522: PetscErrorCode (*remove)(MatNullSpace,Vec,void*); /* for user provided removal function */
523: void* rmctx; /* context for remove() function */
524: };
526: /*
527: Checking zero pivot for LU, ILU preconditioners.
528: */
529: typedef struct {
530: PetscInt nshift,nshift_max;
531: PetscReal shift_amount,shift_lo,shift_hi,shift_top,shift_fraction;
532: PetscBool newshift;
533: PetscReal rs; /* active row sum of abs(offdiagonals) */
534: PetscScalar pv; /* pivot of the active row */
535: } FactorShiftCtx;
537: PETSC_EXTERN PetscErrorCode MatFactorDumpMatrix(Mat);
541: PETSC_STATIC_INLINE PetscErrorCode MatPivotCheck_nz(Mat mat,const MatFactorInfo *info,FactorShiftCtx *sctx,PetscInt row)
542: {
543: PetscReal _rs = sctx->rs;
544: PetscReal _zero = info->zeropivot*_rs;
547: if (PetscAbsScalar(sctx->pv) <= _zero){
548: /* force |diag| > zeropivot*rs */
549: if (!sctx->nshift) sctx->shift_amount = info->shiftamount;
550: else sctx->shift_amount *= 2.0;
551: sctx->newshift = PETSC_TRUE;
552: (sctx->nshift)++;
553: } else {
554: sctx->newshift = PETSC_FALSE;
555: }
556: return(0);
557: }
561: PETSC_STATIC_INLINE PetscErrorCode MatPivotCheck_pd(Mat mat,const MatFactorInfo *info,FactorShiftCtx *sctx,PetscInt row)
562: {
563: PetscReal _rs = sctx->rs;
564: PetscReal _zero = info->zeropivot*_rs;
567: if (PetscRealPart(sctx->pv) <= _zero){
568: /* force matfactor to be diagonally dominant */
569: if (sctx->nshift == sctx->nshift_max) {
570: sctx->shift_fraction = sctx->shift_hi;
571: } else {
572: sctx->shift_lo = sctx->shift_fraction;
573: sctx->shift_fraction = (sctx->shift_hi+sctx->shift_lo)/2.;
574: }
575: sctx->shift_amount = sctx->shift_fraction * sctx->shift_top;
576: sctx->nshift++;
577: sctx->newshift = PETSC_TRUE;
578: } else {
579: sctx->newshift = PETSC_FALSE;
580: }
581: return(0);
582: }
586: PETSC_STATIC_INLINE PetscErrorCode MatPivotCheck_inblocks(Mat mat,const MatFactorInfo *info,FactorShiftCtx *sctx,PetscInt row)
587: {
588: PetscReal _zero = info->zeropivot;
591: if (PetscAbsScalar(sctx->pv) <= _zero){
592: sctx->pv += info->shiftamount;
593: sctx->shift_amount = 0.0;
594: sctx->nshift++;
595: }
596: sctx->newshift = PETSC_FALSE;
597: return(0);
598: }
602: PETSC_STATIC_INLINE PetscErrorCode MatPivotCheck_none(Mat mat,const MatFactorInfo *info,FactorShiftCtx *sctx,PetscInt row)
603: {
604: PetscReal _zero = info->zeropivot;
607: sctx->newshift = PETSC_FALSE;
608: if (PetscAbsScalar(sctx->pv) <= _zero) {
610: PetscBool flg = PETSC_FALSE;
612: PetscOptionsGetBool(NULL,"-mat_dump",&flg,NULL);
613: if (flg) {
614: MatView(mat,PETSC_VIEWER_BINARY_(PetscObjectComm((PetscObject)mat)));
615: }
616: SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot row %D value %g tolerance %g",row,(double)PetscAbsScalar(sctx->pv),(double)_zero);
617: }
618: return(0);
619: }
623: PETSC_STATIC_INLINE PetscErrorCode MatPivotCheck(Mat mat,const MatFactorInfo *info,FactorShiftCtx *sctx,PetscInt row)
624: {
628: if (info->shifttype == (PetscReal) MAT_SHIFT_NONZERO){
629: MatPivotCheck_nz(mat,info,sctx,row);
630: } else if (info->shifttype == (PetscReal) MAT_SHIFT_POSITIVE_DEFINITE){
631: MatPivotCheck_pd(mat,info,sctx,row);
632: } else if (info->shifttype == (PetscReal) MAT_SHIFT_INBLOCKS){
633: MatPivotCheck_inblocks(mat,info,sctx,row);
634: } else {
635: MatPivotCheck_none(mat,info,sctx,row);
636: }
637: return(0);
638: }
640: /*
641: Create and initialize a linked list
642: Input Parameters:
643: idx_start - starting index of the list
644: lnk_max - max value of lnk indicating the end of the list
645: nlnk - max length of the list
646: Output Parameters:
647: lnk - list initialized
648: bt - PetscBT (bitarray) with all bits set to false
649: lnk_empty - flg indicating the list is empty
650: */
651: #define PetscLLCreate(idx_start,lnk_max,nlnk,lnk,bt) \
652: (PetscMalloc1(nlnk,&lnk) || PetscBTCreate(nlnk,&(bt)) || (lnk[idx_start] = lnk_max,0))
654: #define PetscLLCreate_new(idx_start,lnk_max,nlnk,lnk,bt,lnk_empty)\
655: (PetscMalloc1(nlnk,&lnk) || PetscBTCreate(nlnk,&(bt)) || (lnk_empty = PETSC_TRUE,0) ||(lnk[idx_start] = lnk_max,0))
657: /*
658: Add an index set into a sorted linked list
659: Input Parameters:
660: nidx - number of input indices
661: indices - interger array
662: idx_start - starting index of the list
663: lnk - linked list(an integer array) that is created
664: bt - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
665: output Parameters:
666: nlnk - number of newly added indices
667: lnk - the sorted(increasing order) linked list containing new and non-redundate entries from indices
668: bt - updated PetscBT (bitarray)
669: */
670: #define PetscLLAdd(nidx,indices,idx_start,nlnk,lnk,bt) 0;\
671: {\
672: PetscInt _k,_entry,_location,_lnkdata;\
673: nlnk = 0;\
674: _lnkdata = idx_start;\
675: for (_k=0; _k<nidx; _k++){\
676: _entry = indices[_k];\
677: if (!PetscBTLookupSet(bt,_entry)){ /* new entry */\
678: /* search for insertion location */\
679: /* start from the beginning if _entry < previous _entry */\
680: if (_k && _entry < _lnkdata) _lnkdata = idx_start;\
681: do {\
682: _location = _lnkdata;\
683: _lnkdata = lnk[_location];\
684: } while (_entry > _lnkdata);\
685: /* insertion location is found, add entry into lnk */\
686: lnk[_location] = _entry;\
687: lnk[_entry] = _lnkdata;\
688: nlnk++;\
689: _lnkdata = _entry; /* next search starts from here if next_entry > _entry */\
690: }\
691: }\
692: }
694: /*
695: Add a permuted index set into a sorted linked list
696: Input Parameters:
697: nidx - number of input indices
698: indices - interger array
699: perm - permutation of indices
700: idx_start - starting index of the list
701: lnk - linked list(an integer array) that is created
702: bt - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
703: output Parameters:
704: nlnk - number of newly added indices
705: lnk - the sorted(increasing order) linked list containing new and non-redundate entries from indices
706: bt - updated PetscBT (bitarray)
707: */
708: #define PetscLLAddPerm(nidx,indices,perm,idx_start,nlnk,lnk,bt) 0;\
709: {\
710: PetscInt _k,_entry,_location,_lnkdata;\
711: nlnk = 0;\
712: _lnkdata = idx_start;\
713: for (_k=0; _k<nidx; _k++){\
714: _entry = perm[indices[_k]];\
715: if (!PetscBTLookupSet(bt,_entry)){ /* new entry */\
716: /* search for insertion location */\
717: /* start from the beginning if _entry < previous _entry */\
718: if (_k && _entry < _lnkdata) _lnkdata = idx_start;\
719: do {\
720: _location = _lnkdata;\
721: _lnkdata = lnk[_location];\
722: } while (_entry > _lnkdata);\
723: /* insertion location is found, add entry into lnk */\
724: lnk[_location] = _entry;\
725: lnk[_entry] = _lnkdata;\
726: nlnk++;\
727: _lnkdata = _entry; /* next search starts from here if next_entry > _entry */\
728: }\
729: }\
730: }
732: /*
733: Add a SORTED ascending index set into a sorted linked list - same as PetscLLAdd() bus skip 'if (_k && _entry < _lnkdata) _lnkdata = idx_start;'
734: Input Parameters:
735: nidx - number of input indices
736: indices - sorted interger array
737: idx_start - starting index of the list
738: lnk - linked list(an integer array) that is created
739: bt - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
740: output Parameters:
741: nlnk - number of newly added indices
742: lnk - the sorted(increasing order) linked list containing new and non-redundate entries from indices
743: bt - updated PetscBT (bitarray)
744: */
745: #define PetscLLAddSorted(nidx,indices,idx_start,nlnk,lnk,bt) 0;\
746: {\
747: PetscInt _k,_entry,_location,_lnkdata;\
748: nlnk = 0;\
749: _lnkdata = idx_start;\
750: for (_k=0; _k<nidx; _k++){\
751: _entry = indices[_k];\
752: if (!PetscBTLookupSet(bt,_entry)){ /* new entry */\
753: /* search for insertion location */\
754: do {\
755: _location = _lnkdata;\
756: _lnkdata = lnk[_location];\
757: } while (_entry > _lnkdata);\
758: /* insertion location is found, add entry into lnk */\
759: lnk[_location] = _entry;\
760: lnk[_entry] = _lnkdata;\
761: nlnk++;\
762: _lnkdata = _entry; /* next search starts from here */\
763: }\
764: }\
765: }
767: #define PetscLLAddSorted_new(nidx,indices,idx_start,lnk_empty,nlnk,lnk,bt) 0; \
768: {\
769: PetscInt _k,_entry,_location,_lnkdata;\
770: if (lnk_empty){\
771: _lnkdata = idx_start; \
772: for (_k=0; _k<nidx; _k++){ \
773: _entry = indices[_k]; \
774: PetscBTSet(bt,_entry); /* mark the new entry */ \
775: _location = _lnkdata; \
776: _lnkdata = lnk[_location]; \
777: /* insertion location is found, add entry into lnk */ \
778: lnk[_location] = _entry; \
779: lnk[_entry] = _lnkdata; \
780: _lnkdata = _entry; /* next search starts from here */ \
781: } \
782: /*\
783: lnk[indices[nidx-1]] = lnk[idx_start];\
784: lnk[idx_start] = indices[0];\
785: PetscBTSet(bt,indices[0]); \
786: for (_k=1; _k<nidx; _k++){ \
787: PetscBTSet(bt,indices[_k]); \
788: lnk[indices[_k-1]] = indices[_k]; \
789: } \
790: */\
791: nlnk = nidx;\
792: lnk_empty = PETSC_FALSE;\
793: } else {\
794: nlnk = 0; \
795: _lnkdata = idx_start; \
796: for (_k=0; _k<nidx; _k++){ \
797: _entry = indices[_k]; \
798: if (!PetscBTLookupSet(bt,_entry)){ /* new entry */ \
799: /* search for insertion location */ \
800: do { \
801: _location = _lnkdata; \
802: _lnkdata = lnk[_location]; \
803: } while (_entry > _lnkdata); \
804: /* insertion location is found, add entry into lnk */ \
805: lnk[_location] = _entry; \
806: lnk[_entry] = _lnkdata; \
807: nlnk++; \
808: _lnkdata = _entry; /* next search starts from here */ \
809: } \
810: } \
811: } \
812: }
814: /*
815: Add a SORTED index set into a sorted linked list used for LUFactorSymbolic()
816: Same as PetscLLAddSorted() with an additional operation:
817: count the number of input indices that are no larger than 'diag'
818: Input Parameters:
819: indices - sorted interger array
820: idx_start - starting index of the list, index of pivot row
821: lnk - linked list(an integer array) that is created
822: bt - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
823: diag - index of the active row in LUFactorSymbolic
824: nzbd - number of input indices with indices <= idx_start
825: im - im[idx_start] is initialized as num of nonzero entries in row=idx_start
826: output Parameters:
827: nlnk - number of newly added indices
828: lnk - the sorted(increasing order) linked list containing new and non-redundate entries from indices
829: bt - updated PetscBT (bitarray)
830: im - im[idx_start]: unchanged if diag is not an entry
831: : num of entries with indices <= diag if diag is an entry
832: */
833: #define PetscLLAddSortedLU(indices,idx_start,nlnk,lnk,bt,diag,nzbd,im) 0;\
834: {\
835: PetscInt _k,_entry,_location,_lnkdata,_nidx;\
836: nlnk = 0;\
837: _lnkdata = idx_start;\
838: _nidx = im[idx_start] - nzbd; /* num of entries with idx_start < index <= diag */\
839: for (_k=0; _k<_nidx; _k++){\
840: _entry = indices[_k];\
841: nzbd++;\
842: if ( _entry== diag) im[idx_start] = nzbd;\
843: if (!PetscBTLookupSet(bt,_entry)){ /* new entry */\
844: /* search for insertion location */\
845: do {\
846: _location = _lnkdata;\
847: _lnkdata = lnk[_location];\
848: } while (_entry > _lnkdata);\
849: /* insertion location is found, add entry into lnk */\
850: lnk[_location] = _entry;\
851: lnk[_entry] = _lnkdata;\
852: nlnk++;\
853: _lnkdata = _entry; /* next search starts from here */\
854: }\
855: }\
856: }
858: /*
859: Copy data on the list into an array, then initialize the list
860: Input Parameters:
861: idx_start - starting index of the list
862: lnk_max - max value of lnk indicating the end of the list
863: nlnk - number of data on the list to be copied
864: lnk - linked list
865: bt - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
866: output Parameters:
867: indices - array that contains the copied data
868: lnk - linked list that is cleaned and initialize
869: bt - PetscBT (bitarray) with all bits set to false
870: */
871: #define PetscLLClean(idx_start,lnk_max,nlnk,lnk,indices,bt) 0;\
872: {\
873: PetscInt _j,_idx=idx_start;\
874: for (_j=0; _j<nlnk; _j++){\
875: _idx = lnk[_idx];\
876: indices[_j] = _idx;\
877: PetscBTClear(bt,_idx);\
878: }\
879: lnk[idx_start] = lnk_max;\
880: }
881: /*
882: Free memories used by the list
883: */
884: #define PetscLLDestroy(lnk,bt) (PetscFree(lnk) || PetscBTDestroy(&(bt)))
886: /* Routines below are used for incomplete matrix factorization */
887: /*
888: Create and initialize a linked list and its levels
889: Input Parameters:
890: idx_start - starting index of the list
891: lnk_max - max value of lnk indicating the end of the list
892: nlnk - max length of the list
893: Output Parameters:
894: lnk - list initialized
895: lnk_lvl - array of size nlnk for storing levels of lnk
896: bt - PetscBT (bitarray) with all bits set to false
897: */
898: #define PetscIncompleteLLCreate(idx_start,lnk_max,nlnk,lnk,lnk_lvl,bt)\
899: (PetscMalloc1(2*nlnk,&lnk) || PetscBTCreate(nlnk,&(bt)) || (lnk[idx_start] = lnk_max,lnk_lvl = lnk + nlnk,0))
901: /*
902: Initialize a sorted linked list used for ILU and ICC
903: Input Parameters:
904: nidx - number of input idx
905: idx - interger array used for storing column indices
906: idx_start - starting index of the list
907: perm - indices of an IS
908: lnk - linked list(an integer array) that is created
909: lnklvl - levels of lnk
910: bt - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
911: output Parameters:
912: nlnk - number of newly added idx
913: lnk - the sorted(increasing order) linked list containing new and non-redundate entries from idx
914: lnklvl - levels of lnk
915: bt - updated PetscBT (bitarray)
916: */
917: #define PetscIncompleteLLInit(nidx,idx,idx_start,perm,nlnk,lnk,lnklvl,bt) 0;\
918: {\
919: PetscInt _k,_entry,_location,_lnkdata;\
920: nlnk = 0;\
921: _lnkdata = idx_start;\
922: for (_k=0; _k<nidx; _k++){\
923: _entry = perm[idx[_k]];\
924: if (!PetscBTLookupSet(bt,_entry)){ /* new entry */\
925: /* search for insertion location */\
926: if (_k && _entry < _lnkdata) _lnkdata = idx_start;\
927: do {\
928: _location = _lnkdata;\
929: _lnkdata = lnk[_location];\
930: } while (_entry > _lnkdata);\
931: /* insertion location is found, add entry into lnk */\
932: lnk[_location] = _entry;\
933: lnk[_entry] = _lnkdata;\
934: lnklvl[_entry] = 0;\
935: nlnk++;\
936: _lnkdata = _entry; /* next search starts from here if next_entry > _entry */\
937: }\
938: }\
939: }
941: /*
942: Add a SORTED index set into a sorted linked list for ILU
943: Input Parameters:
944: nidx - number of input indices
945: idx - sorted interger array used for storing column indices
946: level - level of fill, e.g., ICC(level)
947: idxlvl - level of idx
948: idx_start - starting index of the list
949: lnk - linked list(an integer array) that is created
950: lnklvl - levels of lnk
951: bt - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
952: prow - the row number of idx
953: output Parameters:
954: nlnk - number of newly added idx
955: lnk - the sorted(increasing order) linked list containing new and non-redundate entries from idx
956: lnklvl - levels of lnk
957: bt - updated PetscBT (bitarray)
959: Note: the level of factor(i,j) is set as lvl(i,j) = min{ lvl(i,j), lvl(i,prow)+lvl(prow,j)+1)
960: where idx = non-zero columns of U(prow,prow+1:n-1), prow<i
961: */
962: #define PetscILULLAddSorted(nidx,idx,level,idxlvl,idx_start,nlnk,lnk,lnklvl,bt,lnklvl_prow) 0;\
963: {\
964: PetscInt _k,_entry,_location,_lnkdata,_incrlev,_lnklvl_prow=lnklvl[prow];\
965: nlnk = 0;\
966: _lnkdata = idx_start;\
967: for (_k=0; _k<nidx; _k++){\
968: _incrlev = idxlvl[_k] + _lnklvl_prow + 1;\
969: if (_incrlev > level) continue;\
970: _entry = idx[_k];\
971: if (!PetscBTLookupSet(bt,_entry)){ /* new entry */\
972: /* search for insertion location */\
973: do {\
974: _location = _lnkdata;\
975: _lnkdata = lnk[_location];\
976: } while (_entry > _lnkdata);\
977: /* insertion location is found, add entry into lnk */\
978: lnk[_location] = _entry;\
979: lnk[_entry] = _lnkdata;\
980: lnklvl[_entry] = _incrlev;\
981: nlnk++;\
982: _lnkdata = _entry; /* next search starts from here if next_entry > _entry */\
983: } else { /* existing entry: update lnklvl */\
984: if (lnklvl[_entry] > _incrlev) lnklvl[_entry] = _incrlev;\
985: }\
986: }\
987: }
989: /*
990: Add a index set into a sorted linked list
991: Input Parameters:
992: nidx - number of input idx
993: idx - interger array used for storing column indices
994: level - level of fill, e.g., ICC(level)
995: idxlvl - level of idx
996: idx_start - starting index of the list
997: lnk - linked list(an integer array) that is created
998: lnklvl - levels of lnk
999: bt - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
1000: output Parameters:
1001: nlnk - number of newly added idx
1002: lnk - the sorted(increasing order) linked list containing new and non-redundate entries from idx
1003: lnklvl - levels of lnk
1004: bt - updated PetscBT (bitarray)
1005: */
1006: #define PetscIncompleteLLAdd(nidx,idx,level,idxlvl,idx_start,nlnk,lnk,lnklvl,bt) 0;\
1007: {\
1008: PetscInt _k,_entry,_location,_lnkdata,_incrlev;\
1009: nlnk = 0;\
1010: _lnkdata = idx_start;\
1011: for (_k=0; _k<nidx; _k++){\
1012: _incrlev = idxlvl[_k] + 1;\
1013: if (_incrlev > level) continue;\
1014: _entry = idx[_k];\
1015: if (!PetscBTLookupSet(bt,_entry)){ /* new entry */\
1016: /* search for insertion location */\
1017: if (_k && _entry < _lnkdata) _lnkdata = idx_start;\
1018: do {\
1019: _location = _lnkdata;\
1020: _lnkdata = lnk[_location];\
1021: } while (_entry > _lnkdata);\
1022: /* insertion location is found, add entry into lnk */\
1023: lnk[_location] = _entry;\
1024: lnk[_entry] = _lnkdata;\
1025: lnklvl[_entry] = _incrlev;\
1026: nlnk++;\
1027: _lnkdata = _entry; /* next search starts from here if next_entry > _entry */\
1028: } else { /* existing entry: update lnklvl */\
1029: if (lnklvl[_entry] > _incrlev) lnklvl[_entry] = _incrlev;\
1030: }\
1031: }\
1032: }
1034: /*
1035: Add a SORTED index set into a sorted linked list
1036: Input Parameters:
1037: nidx - number of input indices
1038: idx - sorted interger array used for storing column indices
1039: level - level of fill, e.g., ICC(level)
1040: idxlvl - level of idx
1041: idx_start - starting index of the list
1042: lnk - linked list(an integer array) that is created
1043: lnklvl - levels of lnk
1044: bt - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
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)
1050: */
1051: #define PetscIncompleteLLAddSorted(nidx,idx,level,idxlvl,idx_start,nlnk,lnk,lnklvl,bt) 0;\
1052: {\
1053: PetscInt _k,_entry,_location,_lnkdata,_incrlev;\
1054: nlnk = 0;\
1055: _lnkdata = idx_start;\
1056: for (_k=0; _k<nidx; _k++){\
1057: _incrlev = idxlvl[_k] + 1;\
1058: if (_incrlev > level) continue;\
1059: _entry = idx[_k];\
1060: if (!PetscBTLookupSet(bt,_entry)){ /* new entry */\
1061: /* search for insertion location */\
1062: do {\
1063: _location = _lnkdata;\
1064: _lnkdata = lnk[_location];\
1065: } while (_entry > _lnkdata);\
1066: /* insertion location is found, add entry into lnk */\
1067: lnk[_location] = _entry;\
1068: lnk[_entry] = _lnkdata;\
1069: lnklvl[_entry] = _incrlev;\
1070: nlnk++;\
1071: _lnkdata = _entry; /* next search starts from here */\
1072: } else { /* existing entry: update lnklvl */\
1073: if (lnklvl[_entry] > _incrlev) lnklvl[_entry] = _incrlev;\
1074: }\
1075: }\
1076: }
1078: /*
1079: Add a SORTED index set into a sorted linked list for ICC
1080: Input Parameters:
1081: nidx - number of input indices
1082: idx - sorted interger array used for storing column indices
1083: level - level of fill, e.g., ICC(level)
1084: idxlvl - level of idx
1085: idx_start - starting index of the list
1086: lnk - linked list(an integer array) that is created
1087: lnklvl - levels of lnk
1088: bt - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
1089: idxlvl_prow - idxlvl[prow], where prow is the row number of the idx
1090: output Parameters:
1091: nlnk - number of newly added indices
1092: lnk - the sorted(increasing order) linked list containing new and non-redundate entries from idx
1093: lnklvl - levels of lnk
1094: bt - updated PetscBT (bitarray)
1095: Note: the level of U(i,j) is set as lvl(i,j) = min{ lvl(i,j), lvl(prow,i)+lvl(prow,j)+1)
1096: where idx = non-zero columns of U(prow,prow+1:n-1), prow<i
1097: */
1098: #define PetscICCLLAddSorted(nidx,idx,level,idxlvl,idx_start,nlnk,lnk,lnklvl,bt,idxlvl_prow) 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] + idxlvl_prow + 1;\
1105: if (_incrlev > level) continue;\
1106: _entry = idx[_k];\
1107: if (!PetscBTLookupSet(bt,_entry)){ /* new entry */\
1108: /* search for insertion location */\
1109: do {\
1110: _location = _lnkdata;\
1111: _lnkdata = lnk[_location];\
1112: } while (_entry > _lnkdata);\
1113: /* insertion location is found, add entry into lnk */\
1114: lnk[_location] = _entry;\
1115: lnk[_entry] = _lnkdata;\
1116: lnklvl[_entry] = _incrlev;\
1117: nlnk++;\
1118: _lnkdata = _entry; /* next search starts from here */\
1119: } else { /* existing entry: update lnklvl */\
1120: if (lnklvl[_entry] > _incrlev) lnklvl[_entry] = _incrlev;\
1121: }\
1122: }\
1123: }
1125: /*
1126: Copy data on the list into an array, then initialize the list
1127: Input Parameters:
1128: idx_start - starting index of the list
1129: lnk_max - max value of lnk indicating the end of the list
1130: nlnk - number of data on the list to be copied
1131: lnk - linked list
1132: lnklvl - level of lnk
1133: bt - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
1134: output Parameters:
1135: indices - array that contains the copied data
1136: lnk - linked list that is cleaned and initialize
1137: lnklvl - level of lnk that is reinitialized
1138: bt - PetscBT (bitarray) with all bits set to false
1139: */
1140: #define PetscIncompleteLLClean(idx_start,lnk_max,nlnk,lnk,lnklvl,indices,indiceslvl,bt) 0;\
1141: {\
1142: PetscInt _j,_idx=idx_start;\
1143: for (_j=0; _j<nlnk; _j++){\
1144: _idx = lnk[_idx];\
1145: *(indices+_j) = _idx;\
1146: *(indiceslvl+_j) = lnklvl[_idx];\
1147: lnklvl[_idx] = -1;\
1148: PetscBTClear(bt,_idx);\
1149: }\
1150: lnk[idx_start] = lnk_max;\
1151: }
1152: /*
1153: Free memories used by the list
1154: */
1155: #define PetscIncompleteLLDestroy(lnk,bt) (PetscFree(lnk) || PetscBTDestroy(&(bt)))
1157: /* -------------------------------------------------------------------------------------------------------*/
1158: #include <petscbt.h>
1161: /*
1162: Create and initialize a condensed linked list -
1163: same as PetscLLCreate(), but uses a scalable array 'lnk' with size of max number of entries, not O(N).
1164: Barry suggested this approach (Dec. 6, 2011):
1165: I've thought of an alternative way of representing a linked list that is efficient but doesn't have the O(N) scaling issue
1166: (it may be faster than the O(N) even sequentially due to less crazy memory access).
1168: 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
1169: 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
1170: 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
1171: 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.
1172: 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
1173: to each other so memory access is much better than using the big array.
1175: Example:
1176: nlnk_max=5, lnk_max=36:
1177: Initial list: [0, 0 | 36, 2 | 0, 0 | 0, 0 | 0, 0 | 0, 0 | 0, 0]
1178: here, head_node has index 2 with value lnk[2]=lnk_max=36,
1179: 0-th entry is used to store the number of entries in the list,
1180: The initial lnk represents head -> tail(marked by 36) with number of entries = lnk[0]=0.
1182: Now adding a sorted set {2,4}, the list becomes
1183: [2, 0 | 36, 4 |2, 6 | 4, 2 | 0, 0 | 0, 0 | 0, 0 ]
1184: represents head -> 2 -> 4 -> tail with number of entries = lnk[0]=2.
1186: Then adding a sorted set {0,3,35}, the list
1187: [5, 0 | 36, 8 | 2, 10 | 4, 12 | 0, 4 | 3, 6 | 35, 2 ]
1188: represents head -> 0 -> 2 -> 3 -> 4 -> 35 -> tail with number of entries = lnk[0]=5.
1190: Input Parameters:
1191: nlnk_max - max length of the list
1192: lnk_max - max value of the entries
1193: Output Parameters:
1194: lnk - list created and initialized
1195: bt - PetscBT (bitarray) with all bits set to false. Note: bt has size lnk_max, not nln_max!
1196: */
1197: PETSC_STATIC_INLINE PetscErrorCode PetscLLCondensedCreate(PetscInt nlnk_max,PetscInt lnk_max,PetscInt **lnk,PetscBT *bt)
1198: {
1200: PetscInt *llnk;
1203: PetscMalloc1(2*(nlnk_max+2),lnk);
1204: PetscBTCreate(lnk_max,bt);
1205: llnk = *lnk;
1206: llnk[0] = 0; /* number of entries on the list */
1207: llnk[2] = lnk_max; /* value in the head node */
1208: llnk[3] = 2; /* next for the head node */
1209: return(0);
1210: }
1214: /*
1215: Add a SORTED ascending index set into a sorted linked list. See PetscLLCondensedCreate() for detailed description.
1216: Input Parameters:
1217: nidx - number of input indices
1218: indices - sorted interger array
1219: lnk - condensed linked list(an integer array) that is created
1220: bt - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
1221: output Parameters:
1222: lnk - the sorted(increasing order) linked list containing previous and newly added non-redundate indices
1223: bt - updated PetscBT (bitarray)
1224: */
1225: PETSC_STATIC_INLINE PetscErrorCode PetscLLCondensedAddSorted(PetscInt nidx,const PetscInt indices[],PetscInt lnk[],PetscBT bt)
1226: {
1227: PetscInt _k,_entry,_location,_next,_lnkdata,_nlnk,_newnode;
1230: _nlnk = lnk[0]; /* num of entries on the input lnk */
1231: _location = 2; /* head */
1232: for (_k=0; _k<nidx; _k++){
1233: _entry = indices[_k];
1234: if (!PetscBTLookupSet(bt,_entry)){ /* new entry */
1235: /* search for insertion location */
1236: do {
1237: _next = _location + 1; /* link from previous node to next node */
1238: _location = lnk[_next]; /* idx of next node */
1239: _lnkdata = lnk[_location];/* value of next node */
1240: } while (_entry > _lnkdata);
1241: /* insertion location is found, add entry into lnk */
1242: _newnode = 2*(_nlnk+2); /* index for this new node */
1243: lnk[_next] = _newnode; /* connect previous node to the new node */
1244: lnk[_newnode] = _entry; /* set value of the new node */
1245: lnk[_newnode+1] = _location; /* connect new node to next node */
1246: _location = _newnode; /* next search starts from the new node */
1247: _nlnk++;
1248: } \
1249: }\
1250: lnk[0] = _nlnk; /* number of entries in the list */
1251: return(0);
1252: }
1256: PETSC_STATIC_INLINE PetscErrorCode PetscLLCondensedClean(PetscInt lnk_max,PetscInt nidx,PetscInt *indices,PetscInt lnk[],PetscBT bt)
1257: {
1259: PetscInt _k,_next,_nlnk;
1262: _next = lnk[3]; /* head node */
1263: _nlnk = lnk[0]; /* num of entries on the list */
1264: for (_k=0; _k<_nlnk; _k++){
1265: indices[_k] = lnk[_next];
1266: _next = lnk[_next + 1];
1267: PetscBTClear(bt,indices[_k]);
1268: }
1269: lnk[0] = 0; /* num of entries on the list */
1270: lnk[2] = lnk_max; /* initialize head node */
1271: lnk[3] = 2; /* head node */
1272: return(0);
1273: }
1277: PETSC_STATIC_INLINE PetscErrorCode PetscLLCondensedView(PetscInt *lnk)
1278: {
1280: PetscInt k;
1283: PetscPrintf(PETSC_COMM_SELF,"LLCondensed of size %d, (val, next)\n",lnk[0]);
1284: for (k=2; k< lnk[0]+2; k++){
1285: PetscPrintf(PETSC_COMM_SELF," %D: (%D, %D)\n",2*k,lnk[2*k],lnk[2*k+1]);
1286: }
1287: return(0);
1288: }
1292: /*
1293: Free memories used by the list
1294: */
1295: PETSC_STATIC_INLINE PetscErrorCode PetscLLCondensedDestroy(PetscInt *lnk,PetscBT bt)
1296: {
1300: PetscFree(lnk);
1301: PetscBTDestroy(&bt);
1302: return(0);
1303: }
1305: /* -------------------------------------------------------------------------------------------------------*/
1308: /*
1309: Same as PetscLLCondensedCreate(), but does not use non-scalable O(lnk_max) bitarray
1310: Input Parameters:
1311: nlnk_max - max length of the list
1312: Output Parameters:
1313: lnk - list created and initialized
1314: */
1315: PETSC_STATIC_INLINE PetscErrorCode PetscLLCondensedCreate_Scalable(PetscInt nlnk_max,PetscInt **lnk)
1316: {
1318: PetscInt *llnk;
1321: PetscMalloc1(2*(nlnk_max+2),lnk);
1322: llnk = *lnk;
1323: llnk[0] = 0; /* number of entries on the list */
1324: llnk[2] = PETSC_MAX_INT; /* value in the head node */
1325: llnk[3] = 2; /* next for the head node */
1326: return(0);
1327: }
1331: PETSC_STATIC_INLINE PetscErrorCode PetscLLCondensedAddSorted_Scalable(PetscInt nidx,const PetscInt indices[],PetscInt lnk[])
1332: {
1333: PetscInt _k,_entry,_location,_next,_lnkdata,_nlnk,_newnode;
1334: _nlnk = lnk[0]; /* num of entries on the input lnk */
1335: _location = 2; /* head */ \
1336: for (_k=0; _k<nidx; _k++){
1337: _entry = indices[_k];
1338: /* search for insertion location */
1339: do {
1340: _next = _location + 1; /* link from previous node to next node */
1341: _location = lnk[_next]; /* idx of next node */
1342: _lnkdata = lnk[_location];/* value of next node */
1343: } while (_entry > _lnkdata);
1344: if (_entry < _lnkdata) {
1345: /* insertion location is found, add entry into lnk */
1346: _newnode = 2*(_nlnk+2); /* index for this new node */
1347: lnk[_next] = _newnode; /* connect previous node to the new node */
1348: lnk[_newnode] = _entry; /* set value of the new node */
1349: lnk[_newnode+1] = _location; /* connect new node to next node */
1350: _location = _newnode; /* next search starts from the new node */
1351: _nlnk++;
1352: }
1353: }
1354: lnk[0] = _nlnk; /* number of entries in the list */
1355: return 0;
1356: }
1360: PETSC_STATIC_INLINE PetscErrorCode PetscLLCondensedClean_Scalable(PetscInt nidx,PetscInt *indices,PetscInt *lnk)
1361: {
1362: PetscInt _k,_next,_nlnk;
1363: _next = lnk[3]; /* head node */
1364: _nlnk = lnk[0];
1365: for (_k=0; _k<_nlnk; _k++){
1366: indices[_k] = lnk[_next];
1367: _next = lnk[_next + 1];
1368: }
1369: lnk[0] = 0; /* num of entries on the list */
1370: lnk[3] = 2; /* head node */
1371: return 0;
1372: }
1376: PETSC_STATIC_INLINE PetscErrorCode PetscLLCondensedDestroy_Scalable(PetscInt *lnk)
1377: {
1378: return PetscFree(lnk);
1379: }
1381: /* -------------------------------------------------------------------------------------------------------*/
1382: /*
1383: lnk[0] number of links
1384: lnk[1] number of entries
1385: lnk[3n] value
1386: lnk[3n+1] len
1387: lnk[3n+2] link to next value
1389: The next three are always the first link
1391: lnk[3] PETSC_MIN_INT+1
1392: lnk[4] 1
1393: lnk[5] link to first real entry
1395: The next three are always the last link
1397: lnk[6] PETSC_MAX_INT - 1
1398: lnk[7] 1
1399: lnk[8] next valid link (this is the same as lnk[0] but without the decreases)
1400: */
1404: PETSC_STATIC_INLINE PetscErrorCode PetscLLCondensedCreate_fast(PetscInt nlnk_max,PetscInt **lnk)
1405: {
1407: PetscInt *llnk;
1410: PetscMalloc1(3*(nlnk_max+3),lnk);
1411: llnk = *lnk;
1412: llnk[0] = 0; /* nlnk: number of entries on the list */
1413: llnk[1] = 0; /* number of integer entries represented in list */
1414: llnk[3] = PETSC_MIN_INT+1; /* value in the first node */
1415: llnk[4] = 1; /* count for the first node */
1416: llnk[5] = 6; /* next for the first node */
1417: llnk[6] = PETSC_MAX_INT-1; /* value in the last node */
1418: llnk[7] = 1; /* count for the last node */
1419: llnk[8] = 0; /* next valid node to be used */
1420: return(0);
1421: }
1423: PETSC_STATIC_INLINE PetscErrorCode PetscLLCondensedAddSorted_fast(PetscInt nidx,const PetscInt indices[],PetscInt lnk[])
1424: {
1425: PetscInt k,entry,prev,next;
1426: prev = 3; /* first value */
1427: next = lnk[prev+2];
1428: for (k=0; k<nidx; k++){
1429: entry = indices[k];
1430: /* search for insertion location */
1431: while (entry >= lnk[next]) {
1432: prev = next;
1433: next = lnk[next+2];
1434: }
1435: /* entry is in range of previous list */
1436: if (entry < lnk[prev]+lnk[prev+1]) continue;
1437: lnk[1]++;
1438: /* entry is right after previous list */
1439: if (entry == lnk[prev]+lnk[prev+1]) {
1440: lnk[prev+1]++;
1441: if (lnk[next] == entry+1) { /* combine two contiquous strings */
1442: lnk[prev+1] += lnk[next+1];
1443: lnk[prev+2] = lnk[next+2];
1444: next = lnk[next+2];
1445: lnk[0]--;
1446: }
1447: continue;
1448: }
1449: /* entry is right before next list */
1450: if (entry == lnk[next]-1) {
1451: lnk[next]--;
1452: lnk[next+1]++;
1453: prev = next;
1454: next = lnk[prev+2];
1455: continue;
1456: }
1457: /* add entry into lnk */
1458: lnk[prev+2] = 3*((lnk[8]++)+3); /* connect previous node to the new node */
1459: prev = lnk[prev+2];
1460: lnk[prev] = entry; /* set value of the new node */
1461: lnk[prev+1] = 1; /* number of values in contiquous string is one to start */
1462: lnk[prev+2] = next; /* connect new node to next node */
1463: lnk[0]++;
1464: }
1465: return 0;
1466: }
1468: PETSC_STATIC_INLINE PetscErrorCode PetscLLCondensedClean_fast(PetscInt nidx,PetscInt *indices,PetscInt *lnk)
1469: {
1470: PetscInt _k,_next,_nlnk,cnt,j;
1471: _next = lnk[5]; /* first node */
1472: _nlnk = lnk[0];
1473: cnt = 0;
1474: for (_k=0; _k<_nlnk; _k++){
1475: for (j=0; j<lnk[_next+1]; j++) {
1476: indices[cnt++] = lnk[_next] + j;
1477: }
1478: _next = lnk[_next + 2];
1479: }
1480: lnk[0] = 0; /* nlnk: number of links */
1481: lnk[1] = 0; /* number of integer entries represented in list */
1482: lnk[3] = PETSC_MIN_INT+1; /* value in the first node */
1483: lnk[4] = 1; /* count for the first node */
1484: lnk[5] = 6; /* next for the first node */
1485: lnk[6] = PETSC_MAX_INT-1; /* value in the last node */
1486: lnk[7] = 1; /* count for the last node */
1487: lnk[8] = 0; /* next valid location to make link */
1488: return 0;
1489: }
1491: PETSC_STATIC_INLINE PetscErrorCode PetscLLCondensedView_fast(PetscInt *lnk)
1492: {
1493: PetscInt k,next,nlnk;
1494: next = lnk[5]; /* first node */
1495: nlnk = lnk[0];
1496: for (k=0; k<nlnk; k++){
1497: #if 0 /* Debugging code */
1498: printf("%d value %d len %d next %d\n",next,lnk[next],lnk[next+1],lnk[next+2]);
1499: #endif
1500: next = lnk[next + 2];
1501: }
1502: return 0;
1503: }
1505: PETSC_STATIC_INLINE PetscErrorCode PetscLLCondensedDestroy_fast(PetscInt *lnk)
1506: {
1507: return PetscFree(lnk);
1508: }
1510: PETSC_EXTERN PetscLogEvent MAT_Mult, MAT_MultMatrixFree, MAT_Mults, MAT_MultConstrained, MAT_MultAdd, MAT_MultTranspose;
1511: PETSC_EXTERN PetscLogEvent MAT_MultTransposeConstrained, MAT_MultTransposeAdd, MAT_Solve, MAT_Solves, MAT_SolveAdd, MAT_SolveTranspose;
1512: PETSC_EXTERN PetscLogEvent MAT_SolveTransposeAdd, MAT_SOR, MAT_ForwardSolve, MAT_BackwardSolve, MAT_LUFactor, MAT_LUFactorSymbolic;
1513: PETSC_EXTERN PetscLogEvent MAT_LUFactorNumeric, MAT_CholeskyFactor, MAT_CholeskyFactorSymbolic, MAT_CholeskyFactorNumeric, MAT_ILUFactor;
1514: PETSC_EXTERN PetscLogEvent MAT_ILUFactorSymbolic, MAT_ICCFactorSymbolic, MAT_Copy, MAT_Convert, MAT_Scale, MAT_AssemblyBegin;
1515: PETSC_EXTERN PetscLogEvent MAT_AssemblyEnd, MAT_SetValues, MAT_GetValues, MAT_GetRow, MAT_GetRowIJ, MAT_GetSubMatrices, MAT_GetColoring, MAT_GetOrdering, MAT_GetRedundantMatrix;
1516: PETSC_EXTERN PetscLogEvent MAT_IncreaseOverlap, MAT_Partitioning, MAT_Coarsen, MAT_ZeroEntries, MAT_Load, MAT_View, MAT_AXPY, MAT_FDColoringCreate, MAT_TransposeColoringCreate;
1517: PETSC_EXTERN PetscLogEvent MAT_FDColoringSetUp, MAT_FDColoringApply, MAT_Transpose, MAT_FDColoringFunction;
1518: PETSC_EXTERN PetscLogEvent MAT_MatMult, MAT_MatSolve,MAT_MatMultSymbolic, MAT_MatMultNumeric,MAT_Getlocalmatcondensed,MAT_GetBrowsOfAcols,MAT_GetBrowsOfAocols;
1519: PETSC_EXTERN PetscLogEvent MAT_PtAP, MAT_PtAPSymbolic, MAT_PtAPNumeric,MAT_Seqstompinum,MAT_Seqstompisym,MAT_Seqstompi,MAT_Getlocalmat;
1520: PETSC_EXTERN PetscLogEvent MAT_RARt, MAT_RARtSymbolic, MAT_RARtNumeric;
1521: PETSC_EXTERN PetscLogEvent MAT_MatTransposeMult, MAT_MatTransposeMultSymbolic, MAT_MatTransposeMultNumeric;
1522: PETSC_EXTERN PetscLogEvent MAT_TransposeMatMult, MAT_TransposeMatMultSymbolic, MAT_TransposeMatMultNumeric;
1523: PETSC_EXTERN PetscLogEvent MAT_MatMatMult, MAT_MatMatMultSymbolic, MAT_MatMatMultNumeric;
1524: PETSC_EXTERN PetscLogEvent MAT_Applypapt, MAT_Applypapt_symbolic, MAT_Applypapt_numeric;
1525: PETSC_EXTERN PetscLogEvent MAT_Getsymtranspose, MAT_Transpose_SeqAIJ, MAT_Getsymtransreduced,MAT_GetSequentialNonzeroStructure;
1527: PETSC_EXTERN PetscLogEvent MATMFFD_Mult;
1528: PETSC_EXTERN PetscLogEvent MAT_GetMultiProcBlock;
1529: PETSC_EXTERN PetscLogEvent MAT_CUSPCopyToGPU, MAT_CUSPARSECopyToGPU, MAT_SetValuesBatch, MAT_SetValuesBatchI, MAT_SetValuesBatchII, MAT_SetValuesBatchIII, MAT_SetValuesBatchIV;
1530: PETSC_EXTERN PetscLogEvent MAT_ViennaCLCopyToGPU;
1531: PETSC_EXTERN PetscLogEvent MAT_Merge,MAT_Residual;
1532: PETSC_EXTERN PetscLogEvent Mat_Coloring_Apply,Mat_Coloring_Comm,Mat_Coloring_Local,Mat_Coloring_ISCreate,Mat_Coloring_SetUp;
1534: #endif