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