Actual source code: matimpl.h
petsc-3.10.5 2019-03-28
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 (*invertvariableblockdiagonal)(Mat,PetscInt,const PetscInt*,PetscScalar*);
187: PetscErrorCode (*createsubmatricesmpi)(Mat,PetscInt,const IS[], const IS[], MatReuse, Mat**);
188: /*129*/
189: PetscErrorCode (*setvaluesbatch)(Mat,PetscInt,PetscInt,PetscInt*,const PetscScalar*);
190: PetscErrorCode (*transposematmult)(Mat,Mat,MatReuse,PetscReal,Mat*);
191: PetscErrorCode (*transposematmultsymbolic)(Mat,Mat,PetscReal,Mat*);
192: PetscErrorCode (*transposematmultnumeric)(Mat,Mat,Mat);
193: PetscErrorCode (*transposecoloringcreate)(Mat,ISColoring,MatTransposeColoring);
194: /*134*/
195: PetscErrorCode (*transcoloringapplysptoden)(MatTransposeColoring,Mat,Mat);
196: PetscErrorCode (*transcoloringapplydentosp)(MatTransposeColoring,Mat,Mat);
197: PetscErrorCode (*rart)(Mat,Mat,MatReuse,PetscReal,Mat*);
198: PetscErrorCode (*rartsymbolic)(Mat,Mat,PetscReal,Mat*); /* double dispatch wrapper routine */
199: PetscErrorCode (*rartnumeric)(Mat,Mat,Mat); /* double dispatch wrapper routine */
200: /*139*/
201: PetscErrorCode (*setblocksizes)(Mat,PetscInt,PetscInt);
202: PetscErrorCode (*aypx)(Mat,PetscScalar,Mat,MatStructure);
203: PetscErrorCode (*residual)(Mat,Vec,Vec,Vec);
204: PetscErrorCode (*fdcoloringsetup)(Mat,ISColoring,MatFDColoring);
205: PetscErrorCode (*findoffblockdiagonalentries)(Mat,IS*);
206: /*144*/
207: PetscErrorCode (*creatempimatconcatenateseqmat)(MPI_Comm,Mat,PetscInt,MatReuse,Mat*);
208: PetscErrorCode (*destroysubmatrices)(PetscInt,Mat*[]);
209: PetscErrorCode (*mattransposesolve)(Mat,Mat,Mat);
210: };
211: /*
212: If you add MatOps entries above also add them to the MATOP enum
213: in include/petscmat.h and include/petsc/finclude/petscmat.h
214: */
216: #include <petscsys.h>
217: PETSC_EXTERN PetscErrorCode MatRegisterOp(MPI_Comm, const char[], PetscVoidFunction, const char[], PetscInt, ...);
218: PETSC_EXTERN PetscErrorCode MatQueryOp(MPI_Comm, PetscVoidFunction*, const char[], PetscInt, ...);
220: typedef struct _p_MatBaseName* MatBaseName;
221: struct _p_MatBaseName {
222: char *bname,*sname,*mname;
223: MatBaseName next;
224: };
226: PETSC_EXTERN MatBaseName MatBaseNameList;
228: /*
229: Utility private matrix routines
230: */
231: PETSC_INTERN PetscErrorCode MatFindNonzeroRowsOrCols_Basic(Mat,PetscBool,PetscReal,IS*);
232: PETSC_INTERN PetscErrorCode MatConvert_Basic(Mat,MatType,MatReuse,Mat*);
233: PETSC_INTERN PetscErrorCode MatConvert_Shell(Mat, MatType,MatReuse,Mat*);
234: PETSC_INTERN PetscErrorCode MatCopy_Basic(Mat,Mat,MatStructure);
235: PETSC_INTERN PetscErrorCode MatDiagonalSet_Default(Mat,Vec,InsertMode);
237: #if defined(PETSC_USE_DEBUG)
238: # define MatCheckPreallocated(A,arg) do { \
239: 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); \
240: } while (0)
241: #else
242: # define MatCheckPreallocated(A,arg) do {} while (0)
243: #endif
245: /*
246: The stash is used to temporarily store inserted matrix values that
247: belong to another processor. During the assembly phase the stashed
248: values are moved to the correct processor and
249: */
251: typedef struct _MatStashSpace *PetscMatStashSpace;
253: struct _MatStashSpace {
254: PetscMatStashSpace next;
255: PetscScalar *space_head,*val;
256: PetscInt *idx,*idy;
257: PetscInt total_space_size;
258: PetscInt local_used;
259: PetscInt local_remaining;
260: };
262: PETSC_EXTERN PetscErrorCode PetscMatStashSpaceGet(PetscInt,PetscInt,PetscMatStashSpace *);
263: PETSC_EXTERN PetscErrorCode PetscMatStashSpaceContiguous(PetscInt,PetscMatStashSpace *,PetscScalar *,PetscInt *,PetscInt *);
264: PETSC_EXTERN PetscErrorCode PetscMatStashSpaceDestroy(PetscMatStashSpace*);
266: typedef struct {
267: PetscInt count;
268: } MatStashHeader;
270: typedef struct {
271: void *buffer; /* Of type blocktype, dynamically constructed */
272: PetscInt count;
273: char pending;
274: } MatStashFrame;
276: typedef struct _MatStash MatStash;
277: struct _MatStash {
278: PetscInt nmax; /* maximum stash size */
279: PetscInt umax; /* user specified max-size */
280: PetscInt oldnmax; /* the nmax value used previously */
281: PetscInt n; /* stash size */
282: PetscInt bs; /* block size of the stash */
283: PetscInt reallocs; /* preserve the no of mallocs invoked */
284: PetscMatStashSpace space_head,space; /* linked list to hold stashed global row/column numbers and matrix values */
286: PetscErrorCode (*ScatterBegin)(Mat,MatStash*,PetscInt*);
287: PetscErrorCode (*ScatterGetMesg)(MatStash*,PetscMPIInt*,PetscInt**,PetscInt**,PetscScalar**,PetscInt*);
288: PetscErrorCode (*ScatterEnd)(MatStash*);
289: PetscErrorCode (*ScatterDestroy)(MatStash*);
291: /* The following variables are used for communication */
292: MPI_Comm comm;
293: PetscMPIInt size,rank;
294: PetscMPIInt tag1,tag2;
295: MPI_Request *send_waits; /* array of send requests */
296: MPI_Request *recv_waits; /* array of receive requests */
297: MPI_Status *send_status; /* array of send status */
298: PetscInt nsends,nrecvs; /* numbers of sends and receives */
299: PetscScalar *svalues; /* sending data */
300: PetscInt *sindices;
301: PetscScalar **rvalues; /* receiving data (values) */
302: PetscInt **rindices; /* receiving data (indices) */
303: PetscInt nprocessed; /* number of messages already processed */
304: PetscMPIInt *flg_v; /* indicates what messages have arrived so far and from whom */
305: PetscBool reproduce;
306: PetscInt reproduce_count;
308: /* The following variables are used for BTS communication */
309: PetscBool subset_off_proc; /* Subsequent assemblies will set a subset (perhaps equal) of off-process entries set on first assembly */
310: PetscBool use_status; /* Use MPI_Status to determine number of items in each message */
311: PetscMPIInt nsendranks;
312: PetscMPIInt nrecvranks;
313: PetscMPIInt *sendranks;
314: PetscMPIInt *recvranks;
315: MatStashHeader *sendhdr,*recvhdr;
316: MatStashFrame *sendframes; /* pointers to the main messages */
317: MatStashFrame *recvframes;
318: MatStashFrame *recvframe_active;
319: PetscInt recvframe_i; /* index of block within active frame */
320: PetscMPIInt recvframe_count; /* Count actually sent for current frame */
321: PetscInt recvcount; /* Number of receives processed so far */
322: PetscMPIInt *some_indices; /* From last call to MPI_Waitsome */
323: MPI_Status *some_statuses; /* Statuses from last call to MPI_Waitsome */
324: PetscMPIInt some_count; /* Number of requests completed in last call to MPI_Waitsome */
325: PetscMPIInt some_i; /* Index of request currently being processed */
326: MPI_Request *sendreqs;
327: MPI_Request *recvreqs;
328: PetscSegBuffer segsendblocks;
329: PetscSegBuffer segrecvframe;
330: PetscSegBuffer segrecvblocks;
331: MPI_Datatype blocktype;
332: size_t blocktype_size;
333: InsertMode *insertmode; /* Pointer to check mat->insertmode and set upon message arrival in case no local values have been set. */
334: };
336: PETSC_INTERN PetscErrorCode MatStashCreate_Private(MPI_Comm,PetscInt,MatStash*);
337: PETSC_INTERN PetscErrorCode MatStashDestroy_Private(MatStash*);
338: PETSC_INTERN PetscErrorCode MatStashScatterEnd_Private(MatStash*);
339: PETSC_INTERN PetscErrorCode MatStashSetInitialSize_Private(MatStash*,PetscInt);
340: PETSC_INTERN PetscErrorCode MatStashGetInfo_Private(MatStash*,PetscInt*,PetscInt*);
341: PETSC_INTERN PetscErrorCode MatStashValuesRow_Private(MatStash*,PetscInt,PetscInt,const PetscInt[],const PetscScalar[],PetscBool );
342: PETSC_INTERN PetscErrorCode MatStashValuesCol_Private(MatStash*,PetscInt,PetscInt,const PetscInt[],const PetscScalar[],PetscInt,PetscBool );
343: PETSC_INTERN PetscErrorCode MatStashValuesRowBlocked_Private(MatStash*,PetscInt,PetscInt,const PetscInt[],const PetscScalar[],PetscInt,PetscInt,PetscInt);
344: PETSC_INTERN PetscErrorCode MatStashValuesColBlocked_Private(MatStash*,PetscInt,PetscInt,const PetscInt[],const PetscScalar[],PetscInt,PetscInt,PetscInt);
345: PETSC_INTERN PetscErrorCode MatStashScatterBegin_Private(Mat,MatStash*,PetscInt*);
346: PETSC_INTERN PetscErrorCode MatStashScatterGetMesg_Private(MatStash*,PetscMPIInt*,PetscInt**,PetscInt**,PetscScalar**,PetscInt*);
347: PETSC_INTERN PetscErrorCode MatGetInfo_External(Mat,MatInfoType,MatInfo*);
349: typedef struct {
350: PetscInt dim;
351: PetscInt dims[4];
352: PetscInt starts[4];
353: PetscBool noc; /* this is a single component problem, hence user will not set MatStencil.c */
354: } MatStencilInfo;
356: /* Info about using compressed row format */
357: typedef struct {
358: PetscBool use; /* indicates compressed rows have been checked and will be used */
359: PetscInt nrows; /* number of non-zero rows */
360: PetscInt *i; /* compressed row pointer */
361: PetscInt *rindex; /* compressed row index */
362: } Mat_CompressedRow;
363: PETSC_EXTERN PetscErrorCode MatCheckCompressedRow(Mat,PetscInt,Mat_CompressedRow*,PetscInt*,PetscInt,PetscReal);
365: typedef struct { /* used by MatCreateRedundantMatrix() for reusing matredundant */
366: PetscInt nzlocal,nsends,nrecvs;
367: PetscMPIInt *send_rank,*recv_rank;
368: PetscInt *sbuf_nz,*rbuf_nz,*sbuf_j,**rbuf_j;
369: PetscScalar *sbuf_a,**rbuf_a;
370: MPI_Comm subcomm; /* when user does not provide a subcomm */
371: IS isrow,iscol;
372: Mat *matseq;
373: } Mat_Redundant;
375: struct _p_Mat {
376: PETSCHEADER(struct _MatOps);
377: PetscLayout rmap,cmap;
378: void *data; /* implementation-specific data */
379: MatFactorType factortype; /* MAT_FACTOR_LU, ILU, CHOLESKY or ICC */
380: PetscBool assembled; /* is the matrix assembled? */
381: PetscBool was_assembled; /* new values inserted into assembled mat */
382: PetscInt num_ass; /* number of times matrix has been assembled */
383: PetscObjectState nonzerostate; /* each time new nonzeros locations are introduced into the matrix this is updated */
384: MatInfo info; /* matrix information */
385: InsertMode insertmode; /* have values been inserted in matrix or added? */
386: MatStash stash,bstash; /* used for assembling off-proc mat emements */
387: MatNullSpace nullsp; /* null space (operator is singular) */
388: MatNullSpace transnullsp; /* null space of transpose of operator */
389: MatNullSpace nearnullsp; /* near null space to be used by multigrid methods */
390: PetscInt congruentlayouts; /* are the rows and columns layouts congruent? */
391: PetscBool preallocated;
392: MatStencilInfo stencil; /* information for structured grid */
393: PetscBool symmetric,hermitian,structurally_symmetric,spd;
394: PetscBool symmetric_set,hermitian_set,structurally_symmetric_set,spd_set; /* if true, then corresponding flag is correct*/
395: PetscBool symmetric_eternal;
396: PetscBool nooffprocentries,nooffproczerorows;
397: PetscBool subsetoffprocentries;
398: PetscBool submat_singleis; /* for efficient PCSetUP_ASM() */
399: PetscBool structure_only;
400: #if defined(PETSC_HAVE_VIENNACL) || defined(PETSC_HAVE_VECCUDA)
401: PetscOffloadFlag valid_GPU_matrix; /* flag pointing to the matrix on the gpu*/
402: #endif
403: void *spptr; /* pointer for special library like SuperLU */
404: char *solvertype;
405: PetscBool checksymmetryonassembly,checknullspaceonassembly;
406: PetscReal checksymmetrytol;
407: Mat schur; /* Schur complement matrix */
408: MatFactorSchurStatus schur_status; /* status of the Schur complement matrix */
409: Mat_Redundant *redundant; /* used by MatCreateRedundantMatrix() */
410: PetscBool erroriffailure; /* Generate an error if detected (for example a zero pivot) instead of returning */
411: MatFactorError factorerrortype; /* type of error in factorization */
412: PetscReal factorerror_zeropivot_value; /* If numerical zero pivot was detected this is the computed value */
413: PetscInt factorerror_zeropivot_row; /* Row where zero pivot was detected */
414: PetscInt nblocks,*bsizes; /* support for MatSetVariableBlockSizes() */
415: char *defaultvectype;
416: };
418: PETSC_INTERN PetscErrorCode MatAXPY_Basic(Mat,PetscScalar,Mat,MatStructure);
419: PETSC_INTERN PetscErrorCode MatAXPY_BasicWithPreallocation(Mat,Mat,PetscScalar,Mat,MatStructure);
421: /*
422: Utility for MatFactor (Schur complement)
423: */
424: PETSC_INTERN PetscErrorCode MatFactorFactorizeSchurComplement_Private(Mat);
425: PETSC_INTERN PetscErrorCode MatFactorInvertSchurComplement_Private(Mat);
426: PETSC_INTERN PetscErrorCode MatFactorUpdateSchurStatus_Private(Mat);
427: PETSC_INTERN PetscErrorCode MatFactorSetUpInPlaceSchur_Private(Mat);
429: /*
430: Utility for MatZeroRows
431: */
432: PETSC_INTERN PetscErrorCode MatZeroRowsMapLocal_Private(Mat,PetscInt,const PetscInt*,PetscInt*,PetscInt**);
434: /*
435: Object for partitioning graphs
436: */
438: typedef struct _MatPartitioningOps *MatPartitioningOps;
439: struct _MatPartitioningOps {
440: PetscErrorCode (*apply)(MatPartitioning,IS*);
441: PetscErrorCode (*applynd)(MatPartitioning,IS*);
442: PetscErrorCode (*setfromoptions)(PetscOptionItems*,MatPartitioning);
443: PetscErrorCode (*destroy)(MatPartitioning);
444: PetscErrorCode (*view)(MatPartitioning,PetscViewer);
445: };
447: struct _p_MatPartitioning {
448: PETSCHEADER(struct _MatPartitioningOps);
449: Mat adj;
450: PetscInt *vertex_weights;
451: PetscReal *part_weights;
452: PetscInt n; /* number of partitions */
453: void *data;
454: PetscInt setupcalled;
455: };
457: /* needed for parallel nested dissection by ParMetis and PTSCOTCH */
458: PETSC_INTERN PetscErrorCode MatPartitioningSizesToSep_Private(PetscInt,PetscInt[],PetscInt[],PetscInt[]);
460: /*
461: Object for coarsen graphs
462: */
463: typedef struct _MatCoarsenOps *MatCoarsenOps;
464: struct _MatCoarsenOps {
465: PetscErrorCode (*apply)(MatCoarsen);
466: PetscErrorCode (*setfromoptions)(PetscOptionItems*,MatCoarsen);
467: PetscErrorCode (*destroy)(MatCoarsen);
468: PetscErrorCode (*view)(MatCoarsen,PetscViewer);
469: };
471: struct _p_MatCoarsen {
472: PETSCHEADER(struct _MatCoarsenOps);
473: Mat graph;
474: PetscInt setupcalled;
475: void *subctx;
476: /* */
477: PetscBool strict_aggs;
478: IS perm;
479: PetscCoarsenData *agg_lists;
480: };
482: /*
483: MatFDColoring is used to compute Jacobian matrices efficiently
484: via coloring. The data structure is explained below in an example.
486: Color = 0 1 0 2 | 2 3 0
487: ---------------------------------------------------
488: 00 01 | 05
489: 10 11 | 14 15 Processor 0
490: 22 23 | 25
491: 32 33 |
492: ===================================================
493: | 44 45 46
494: 50 | 55 Processor 1
495: | 64 66
496: ---------------------------------------------------
498: ncolors = 4;
500: ncolumns = {2,1,1,0}
501: columns = {{0,2},{1},{3},{}}
502: nrows = {4,2,3,3}
503: rows = {{0,1,2,3},{0,1},{1,2,3},{0,1,2}}
504: vwscale = {dx(0),dx(1),dx(2),dx(3)} MPI Vec
505: vscale = {dx(0),dx(1),dx(2),dx(3),dx(4),dx(5)} Seq Vec
507: ncolumns = {1,0,1,1}
508: columns = {{6},{},{4},{5}}
509: nrows = {3,0,2,2}
510: rows = {{0,1,2},{},{1,2},{1,2}}
511: vwscale = {dx(4),dx(5),dx(6)} MPI Vec
512: vscale = {dx(0),dx(4),dx(5),dx(6)} Seq Vec
514: See the routine MatFDColoringApply() for how this data is used
515: to compute the Jacobian.
517: */
518: typedef struct {
519: PetscInt row;
520: PetscInt col;
521: PetscScalar *valaddr; /* address of value */
522: } MatEntry;
524: typedef struct {
525: PetscInt row;
526: PetscScalar *valaddr; /* address of value */
527: } MatEntry2;
529: struct _p_MatFDColoring{
530: PETSCHEADER(int);
531: PetscInt M,N,m; /* total rows, columns; local rows */
532: PetscInt rstart; /* first row owned by local processor */
533: PetscInt ncolors; /* number of colors */
534: PetscInt *ncolumns; /* number of local columns for a color */
535: PetscInt **columns; /* lists the local columns of each color (using global column numbering) */
536: PetscInt *nrows; /* number of local rows for each color */
537: MatEntry *matentry; /* holds (row, column, address of value) for Jacobian matrix entry */
538: MatEntry2 *matentry2; /* holds (row, address of value) for Jacobian matrix entry */
539: PetscScalar *dy; /* store a block of F(x+dx)-F(x) when J is in BAIJ format */
540: PetscReal error_rel; /* square root of relative error in computing function */
541: PetscReal umin; /* minimum allowable u'dx value */
542: Vec w1,w2,w3; /* work vectors used in computing Jacobian */
543: PetscBool fset; /* indicates that the initial function value F(X) is set */
544: PetscErrorCode (*f)(void); /* function that defines Jacobian */
545: void *fctx; /* optional user-defined context for use by the function f */
546: Vec vscale; /* holds FD scaling, i.e. 1/dx for each perturbed column */
547: PetscInt currentcolor; /* color for which function evaluation is being done now */
548: const char *htype; /* "wp" or "ds" */
549: ISColoringType ctype; /* IS_COLORING_GLOBAL or IS_COLORING_LOCAL */
550: PetscInt brows,bcols; /* number of block rows or columns for speedup inserting the dense matrix into sparse Jacobian */
551: PetscBool setupcalled; /* true if setup has been called */
552: PetscBool viewed; /* true if the -mat_fd_coloring_view has been triggered already */
553: void (*ftn_func_pointer)(void),*ftn_func_cntx; /* serve the same purpose as *fortran_func_pointers in PETSc objects */
554: };
556: typedef struct _MatColoringOps *MatColoringOps;
557: struct _MatColoringOps {
558: PetscErrorCode (*destroy)(MatColoring);
559: PetscErrorCode (*setfromoptions)(PetscOptionItems*,MatColoring);
560: PetscErrorCode (*view)(MatColoring,PetscViewer);
561: PetscErrorCode (*apply)(MatColoring,ISColoring*);
562: PetscErrorCode (*weights)(MatColoring,PetscReal**,PetscInt**);
563: };
565: struct _p_MatColoring {
566: PETSCHEADER(struct _MatColoringOps);
567: Mat mat;
568: PetscInt dist; /* distance of the coloring */
569: PetscInt maxcolors; /* the maximum number of colors returned, maxcolors=1 for MIS */
570: void *data; /* inner context */
571: PetscBool valid; /* check to see if what is produced is a valid coloring */
572: MatColoringWeightType weight_type; /* type of weight computation to be performed */
573: PetscReal *user_weights; /* custom weights and permutation */
574: PetscInt *user_lperm;
575: PetscBool valid_iscoloring; /* check to see if matcoloring is produced a valid iscoloring */
576: };
578: struct _p_MatTransposeColoring{
579: PETSCHEADER(int);
580: PetscInt M,N,m; /* total rows, columns; local rows */
581: PetscInt rstart; /* first row owned by local processor */
582: PetscInt ncolors; /* number of colors */
583: PetscInt *ncolumns; /* number of local columns for a color */
584: PetscInt *nrows; /* number of local rows for each color */
585: PetscInt currentcolor; /* color for which function evaluation is being done now */
586: ISColoringType ctype; /* IS_COLORING_GLOBAL or IS_COLORING_LOCAL */
588: PetscInt *colorforrow,*colorforcol; /* pointer to rows and columns */
589: PetscInt *rows; /* lists the local rows for each color (using the local row numbering) */
590: PetscInt *den2sp; /* maps (row,color) in the dense matrix to index of sparse matrix array a->a */
591: PetscInt *columns; /* lists the local columns of each color (using global column numbering) */
592: PetscInt brows; /* number of rows for efficient implementation of MatTransColoringApplyDenToSp() */
593: PetscInt *lstart; /* array used for loop over row blocks of Csparse */
594: };
596: /*
597: Null space context for preconditioner/operators
598: */
599: struct _p_MatNullSpace {
600: PETSCHEADER(int);
601: PetscBool has_cnst;
602: PetscInt n;
603: Vec* vecs;
604: PetscScalar* alpha; /* for projections */
605: PetscErrorCode (*remove)(MatNullSpace,Vec,void*); /* for user provided removal function */
606: void* rmctx; /* context for remove() function */
607: };
609: /*
610: Checking zero pivot for LU, ILU preconditioners.
611: */
612: typedef struct {
613: PetscInt nshift,nshift_max;
614: PetscReal shift_amount,shift_lo,shift_hi,shift_top,shift_fraction;
615: PetscBool newshift;
616: PetscReal rs; /* active row sum of abs(offdiagonals) */
617: PetscScalar pv; /* pivot of the active row */
618: } FactorShiftCtx;
620: /*
621: Used by MatCreateSubMatrices_MPIXAIJ_Local()
622: */
623: #include <petscctable.h>
624: typedef struct { /* used by MatCreateSubMatrices_MPIAIJ_SingleIS_Local() and MatCreateSubMatrices_MPIAIJ_Local */
625: PetscInt id; /* index of submats, only submats[0] is responsible for deleting some arrays below */
626: PetscInt nrqs,nrqr;
627: PetscInt **rbuf1,**rbuf2,**rbuf3,**sbuf1,**sbuf2;
628: PetscInt **ptr;
629: PetscInt *tmp;
630: PetscInt *ctr;
631: PetscInt *pa; /* proc array */
632: PetscInt *req_size,*req_source1,*req_source2;
633: PetscBool allcolumns,allrows;
634: PetscBool singleis;
635: PetscInt *row2proc; /* row to proc map */
636: PetscInt nstages;
637: #if defined(PETSC_USE_CTABLE)
638: PetscTable cmap,rmap;
639: PetscInt *cmap_loc,*rmap_loc;
640: #else
641: PetscInt *cmap,*rmap;
642: #endif
644: PetscErrorCode (*destroy)(Mat);
645: } Mat_SubSppt;
647: PETSC_EXTERN PetscErrorCode MatFactorDumpMatrix(Mat);
648: PETSC_INTERN PetscErrorCode MatShift_Basic(Mat,PetscScalar);
649: PETSC_INTERN PetscErrorCode MatSetBlockSizes_Default(Mat,PetscInt,PetscInt);
651: PETSC_STATIC_INLINE PetscErrorCode MatPivotCheck_nz(Mat mat,const MatFactorInfo *info,FactorShiftCtx *sctx,PetscInt row)
652: {
653: PetscReal _rs = sctx->rs;
654: PetscReal _zero = info->zeropivot*_rs;
657: if (PetscAbsScalar(sctx->pv) <= _zero && !PetscIsNanScalar(sctx->pv)){
658: /* force |diag| > zeropivot*rs */
659: if (!sctx->nshift) sctx->shift_amount = info->shiftamount;
660: else sctx->shift_amount *= 2.0;
661: sctx->newshift = PETSC_TRUE;
662: (sctx->nshift)++;
663: } else {
664: sctx->newshift = PETSC_FALSE;
665: }
666: return(0);
667: }
669: PETSC_STATIC_INLINE PetscErrorCode MatPivotCheck_pd(Mat mat,const MatFactorInfo *info,FactorShiftCtx *sctx,PetscInt row)
670: {
671: PetscReal _rs = sctx->rs;
672: PetscReal _zero = info->zeropivot*_rs;
675: if (PetscRealPart(sctx->pv) <= _zero && !PetscIsNanScalar(sctx->pv)){
676: /* force matfactor to be diagonally dominant */
677: if (sctx->nshift == sctx->nshift_max) {
678: sctx->shift_fraction = sctx->shift_hi;
679: } else {
680: sctx->shift_lo = sctx->shift_fraction;
681: sctx->shift_fraction = (sctx->shift_hi+sctx->shift_lo)/2.;
682: }
683: sctx->shift_amount = sctx->shift_fraction * sctx->shift_top;
684: sctx->nshift++;
685: sctx->newshift = PETSC_TRUE;
686: } else {
687: sctx->newshift = PETSC_FALSE;
688: }
689: return(0);
690: }
692: PETSC_STATIC_INLINE PetscErrorCode MatPivotCheck_inblocks(Mat mat,const MatFactorInfo *info,FactorShiftCtx *sctx,PetscInt row)
693: {
694: PetscReal _zero = info->zeropivot;
697: if (PetscAbsScalar(sctx->pv) <= _zero && !PetscIsNanScalar(sctx->pv)){
698: sctx->pv += info->shiftamount;
699: sctx->shift_amount = 0.0;
700: sctx->nshift++;
701: }
702: sctx->newshift = PETSC_FALSE;
703: return(0);
704: }
706: PETSC_STATIC_INLINE PetscErrorCode MatPivotCheck_none(Mat fact,Mat mat,const MatFactorInfo *info,FactorShiftCtx *sctx,PetscInt row)
707: {
708: PetscReal _zero = info->zeropivot;
712: sctx->newshift = PETSC_FALSE;
713: if (PetscAbsScalar(sctx->pv) <= _zero && !PetscIsNanScalar(sctx->pv)) {
714: if (!mat->erroriffailure) {
715: PetscInfo3(mat,"Detected zero pivot in factorization in row %D value %g tolerance %g\n",row,(double)PetscAbsScalar(sctx->pv),(double)_zero);
716: fact->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
717: fact->factorerror_zeropivot_value = PetscAbsScalar(sctx->pv);
718: fact->factorerror_zeropivot_row = row;
719: } else SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot row %D value %g tolerance %g\n",row,(double)PetscAbsScalar(sctx->pv),(double)_zero);
720: }
721: return(0);
722: }
724: PETSC_STATIC_INLINE PetscErrorCode MatPivotCheck(Mat fact,Mat mat,const MatFactorInfo *info,FactorShiftCtx *sctx,PetscInt row)
725: {
729: if (info->shifttype == (PetscReal) MAT_SHIFT_NONZERO){
730: MatPivotCheck_nz(mat,info,sctx,row);
731: } else if (info->shifttype == (PetscReal) MAT_SHIFT_POSITIVE_DEFINITE){
732: MatPivotCheck_pd(mat,info,sctx,row);
733: } else if (info->shifttype == (PetscReal) MAT_SHIFT_INBLOCKS){
734: MatPivotCheck_inblocks(mat,info,sctx,row);
735: } else {
736: MatPivotCheck_none(fact,mat,info,sctx,row);
737: }
738: return(0);
739: }
741: /*
742: Create and initialize a linked list
743: Input Parameters:
744: idx_start - starting index of the list
745: lnk_max - max value of lnk indicating the end of the list
746: nlnk - max length of the list
747: Output Parameters:
748: lnk - list initialized
749: bt - PetscBT (bitarray) with all bits set to false
750: lnk_empty - flg indicating the list is empty
751: */
752: #define PetscLLCreate(idx_start,lnk_max,nlnk,lnk,bt) \
753: (PetscMalloc1(nlnk,&lnk) || PetscBTCreate(nlnk,&(bt)) || (lnk[idx_start] = lnk_max,0))
755: #define PetscLLCreate_new(idx_start,lnk_max,nlnk,lnk,bt,lnk_empty)\
756: (PetscMalloc1(nlnk,&lnk) || PetscBTCreate(nlnk,&(bt)) || (lnk_empty = PETSC_TRUE,0) ||(lnk[idx_start] = lnk_max,0))
758: /*
759: Add an index set into a sorted linked list
760: Input Parameters:
761: nidx - number of input indices
762: indices - integer array
763: idx_start - starting index of the list
764: lnk - linked list(an integer array) that is created
765: bt - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
766: output Parameters:
767: nlnk - number of newly added indices
768: lnk - the sorted(increasing order) linked list containing new and non-redundate entries from indices
769: bt - updated PetscBT (bitarray)
770: */
771: #define PetscLLAdd(nidx,indices,idx_start,nlnk,lnk,bt) 0;\
772: {\
773: PetscInt _k,_entry,_location,_lnkdata;\
774: nlnk = 0;\
775: _lnkdata = idx_start;\
776: for (_k=0; _k<nidx; _k++){\
777: _entry = indices[_k];\
778: if (!PetscBTLookupSet(bt,_entry)){ /* new entry */\
779: /* search for insertion location */\
780: /* start from the beginning if _entry < previous _entry */\
781: if (_k && _entry < _lnkdata) _lnkdata = idx_start;\
782: do {\
783: _location = _lnkdata;\
784: _lnkdata = lnk[_location];\
785: } while (_entry > _lnkdata);\
786: /* insertion location is found, add entry into lnk */\
787: lnk[_location] = _entry;\
788: lnk[_entry] = _lnkdata;\
789: nlnk++;\
790: _lnkdata = _entry; /* next search starts from here if next_entry > _entry */\
791: }\
792: }\
793: }
795: /*
796: Add a permuted index set into a sorted linked list
797: Input Parameters:
798: nidx - number of input indices
799: indices - integer array
800: perm - permutation of indices
801: idx_start - starting index of the list
802: lnk - linked list(an integer array) that is created
803: bt - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
804: output Parameters:
805: nlnk - number of newly added indices
806: lnk - the sorted(increasing order) linked list containing new and non-redundate entries from indices
807: bt - updated PetscBT (bitarray)
808: */
809: #define PetscLLAddPerm(nidx,indices,perm,idx_start,nlnk,lnk,bt) 0;\
810: {\
811: PetscInt _k,_entry,_location,_lnkdata;\
812: nlnk = 0;\
813: _lnkdata = idx_start;\
814: for (_k=0; _k<nidx; _k++){\
815: _entry = perm[indices[_k]];\
816: if (!PetscBTLookupSet(bt,_entry)){ /* new entry */\
817: /* search for insertion location */\
818: /* start from the beginning if _entry < previous _entry */\
819: if (_k && _entry < _lnkdata) _lnkdata = idx_start;\
820: do {\
821: _location = _lnkdata;\
822: _lnkdata = lnk[_location];\
823: } while (_entry > _lnkdata);\
824: /* insertion location is found, add entry into lnk */\
825: lnk[_location] = _entry;\
826: lnk[_entry] = _lnkdata;\
827: nlnk++;\
828: _lnkdata = _entry; /* next search starts from here if next_entry > _entry */\
829: }\
830: }\
831: }
833: /*
834: Add a SORTED ascending index set into a sorted linked list - same as PetscLLAdd() bus skip 'if (_k && _entry < _lnkdata) _lnkdata = idx_start;'
835: Input Parameters:
836: nidx - number of input indices
837: indices - sorted integer array
838: idx_start - starting index of the list
839: lnk - linked list(an integer array) that is created
840: bt - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
841: output Parameters:
842: nlnk - number of newly added indices
843: lnk - the sorted(increasing order) linked list containing new and non-redundate entries from indices
844: bt - updated PetscBT (bitarray)
845: */
846: #define PetscLLAddSorted(nidx,indices,idx_start,nlnk,lnk,bt) 0;\
847: {\
848: PetscInt _k,_entry,_location,_lnkdata;\
849: nlnk = 0;\
850: _lnkdata = idx_start;\
851: for (_k=0; _k<nidx; _k++){\
852: _entry = indices[_k];\
853: if (!PetscBTLookupSet(bt,_entry)){ /* new entry */\
854: /* search for insertion location */\
855: do {\
856: _location = _lnkdata;\
857: _lnkdata = lnk[_location];\
858: } while (_entry > _lnkdata);\
859: /* insertion location is found, add entry into lnk */\
860: lnk[_location] = _entry;\
861: lnk[_entry] = _lnkdata;\
862: nlnk++;\
863: _lnkdata = _entry; /* next search starts from here */\
864: }\
865: }\
866: }
868: #define PetscLLAddSorted_new(nidx,indices,idx_start,lnk_empty,nlnk,lnk,bt) 0; \
869: {\
870: PetscInt _k,_entry,_location,_lnkdata;\
871: if (lnk_empty){\
872: _lnkdata = idx_start; \
873: for (_k=0; _k<nidx; _k++){ \
874: _entry = indices[_k]; \
875: PetscBTSet(bt,_entry); /* mark the new entry */ \
876: _location = _lnkdata; \
877: _lnkdata = lnk[_location]; \
878: /* insertion location is found, add entry into lnk */ \
879: lnk[_location] = _entry; \
880: lnk[_entry] = _lnkdata; \
881: _lnkdata = _entry; /* next search starts from here */ \
882: } \
883: /*\
884: lnk[indices[nidx-1]] = lnk[idx_start];\
885: lnk[idx_start] = indices[0];\
886: PetscBTSet(bt,indices[0]); \
887: for (_k=1; _k<nidx; _k++){ \
888: PetscBTSet(bt,indices[_k]); \
889: lnk[indices[_k-1]] = indices[_k]; \
890: } \
891: */\
892: nlnk = nidx;\
893: lnk_empty = PETSC_FALSE;\
894: } else {\
895: nlnk = 0; \
896: _lnkdata = idx_start; \
897: for (_k=0; _k<nidx; _k++){ \
898: _entry = indices[_k]; \
899: if (!PetscBTLookupSet(bt,_entry)){ /* new entry */ \
900: /* search for insertion location */ \
901: do { \
902: _location = _lnkdata; \
903: _lnkdata = lnk[_location]; \
904: } while (_entry > _lnkdata); \
905: /* insertion location is found, add entry into lnk */ \
906: lnk[_location] = _entry; \
907: lnk[_entry] = _lnkdata; \
908: nlnk++; \
909: _lnkdata = _entry; /* next search starts from here */ \
910: } \
911: } \
912: } \
913: }
915: /*
916: Add a SORTED index set into a sorted linked list used for LUFactorSymbolic()
917: Same as PetscLLAddSorted() with an additional operation:
918: count the number of input indices that are no larger than 'diag'
919: Input Parameters:
920: indices - sorted integer array
921: idx_start - starting index of the list, index of pivot row
922: lnk - linked list(an integer array) that is created
923: bt - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
924: diag - index of the active row in LUFactorSymbolic
925: nzbd - number of input indices with indices <= idx_start
926: im - im[idx_start] is initialized as num of nonzero entries in row=idx_start
927: output Parameters:
928: nlnk - number of newly added indices
929: lnk - the sorted(increasing order) linked list containing new and non-redundate entries from indices
930: bt - updated PetscBT (bitarray)
931: im - im[idx_start]: unchanged if diag is not an entry
932: : num of entries with indices <= diag if diag is an entry
933: */
934: #define PetscLLAddSortedLU(indices,idx_start,nlnk,lnk,bt,diag,nzbd,im) 0;\
935: {\
936: PetscInt _k,_entry,_location,_lnkdata,_nidx;\
937: nlnk = 0;\
938: _lnkdata = idx_start;\
939: _nidx = im[idx_start] - nzbd; /* num of entries with idx_start < index <= diag */\
940: for (_k=0; _k<_nidx; _k++){\
941: _entry = indices[_k];\
942: nzbd++;\
943: if ( _entry== diag) im[idx_start] = nzbd;\
944: if (!PetscBTLookupSet(bt,_entry)){ /* new entry */\
945: /* search for insertion location */\
946: do {\
947: _location = _lnkdata;\
948: _lnkdata = lnk[_location];\
949: } while (_entry > _lnkdata);\
950: /* insertion location is found, add entry into lnk */\
951: lnk[_location] = _entry;\
952: lnk[_entry] = _lnkdata;\
953: nlnk++;\
954: _lnkdata = _entry; /* next search starts from here */\
955: }\
956: }\
957: }
959: /*
960: Copy data on the list into an array, then initialize the list
961: Input Parameters:
962: idx_start - starting index of the list
963: lnk_max - max value of lnk indicating the end of the list
964: nlnk - number of data on the list to be copied
965: lnk - linked list
966: bt - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
967: output Parameters:
968: indices - array that contains the copied data
969: lnk - linked list that is cleaned and initialize
970: bt - PetscBT (bitarray) with all bits set to false
971: */
972: #define PetscLLClean(idx_start,lnk_max,nlnk,lnk,indices,bt) 0;\
973: {\
974: PetscInt _j,_idx=idx_start;\
975: for (_j=0; _j<nlnk; _j++){\
976: _idx = lnk[_idx];\
977: indices[_j] = _idx;\
978: PetscBTClear(bt,_idx);\
979: }\
980: lnk[idx_start] = lnk_max;\
981: }
982: /*
983: Free memories used by the list
984: */
985: #define PetscLLDestroy(lnk,bt) (PetscFree(lnk) || PetscBTDestroy(&(bt)))
987: /* Routines below are used for incomplete matrix factorization */
988: /*
989: Create and initialize a linked list and its levels
990: Input Parameters:
991: idx_start - starting index of the list
992: lnk_max - max value of lnk indicating the end of the list
993: nlnk - max length of the list
994: Output Parameters:
995: lnk - list initialized
996: lnk_lvl - array of size nlnk for storing levels of lnk
997: bt - PetscBT (bitarray) with all bits set to false
998: */
999: #define PetscIncompleteLLCreate(idx_start,lnk_max,nlnk,lnk,lnk_lvl,bt)\
1000: (PetscIntMultError(2,nlnk,NULL) || PetscMalloc1(2*nlnk,&lnk) || PetscBTCreate(nlnk,&(bt)) || (lnk[idx_start] = lnk_max,lnk_lvl = lnk + nlnk,0))
1002: /*
1003: Initialize a sorted linked list used for ILU and ICC
1004: Input Parameters:
1005: nidx - number of input idx
1006: idx - integer array used for storing column indices
1007: idx_start - starting index of the list
1008: perm - indices of an IS
1009: lnk - linked list(an integer array) that is created
1010: lnklvl - levels of lnk
1011: bt - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
1012: output Parameters:
1013: nlnk - number of newly added idx
1014: lnk - the sorted(increasing order) linked list containing new and non-redundate entries from idx
1015: lnklvl - levels of lnk
1016: bt - updated PetscBT (bitarray)
1017: */
1018: #define PetscIncompleteLLInit(nidx,idx,idx_start,perm,nlnk,lnk,lnklvl,bt) 0;\
1019: {\
1020: PetscInt _k,_entry,_location,_lnkdata;\
1021: nlnk = 0;\
1022: _lnkdata = idx_start;\
1023: for (_k=0; _k<nidx; _k++){\
1024: _entry = perm[idx[_k]];\
1025: if (!PetscBTLookupSet(bt,_entry)){ /* new entry */\
1026: /* search for insertion location */\
1027: if (_k && _entry < _lnkdata) _lnkdata = idx_start;\
1028: do {\
1029: _location = _lnkdata;\
1030: _lnkdata = lnk[_location];\
1031: } while (_entry > _lnkdata);\
1032: /* insertion location is found, add entry into lnk */\
1033: lnk[_location] = _entry;\
1034: lnk[_entry] = _lnkdata;\
1035: lnklvl[_entry] = 0;\
1036: nlnk++;\
1037: _lnkdata = _entry; /* next search starts from here if next_entry > _entry */\
1038: }\
1039: }\
1040: }
1042: /*
1043: Add a SORTED index set into a sorted linked list for ILU
1044: Input Parameters:
1045: nidx - number of input indices
1046: idx - sorted integer array used for storing column indices
1047: level - level of fill, e.g., ICC(level)
1048: idxlvl - level of idx
1049: idx_start - starting index of the list
1050: lnk - linked list(an integer array) that is created
1051: lnklvl - levels of lnk
1052: bt - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
1053: prow - the row number of idx
1054: output Parameters:
1055: nlnk - number of newly added idx
1056: lnk - the sorted(increasing order) linked list containing new and non-redundate entries from idx
1057: lnklvl - levels of lnk
1058: bt - updated PetscBT (bitarray)
1060: Note: the level of factor(i,j) is set as lvl(i,j) = min{ lvl(i,j), lvl(i,prow)+lvl(prow,j)+1)
1061: where idx = non-zero columns of U(prow,prow+1:n-1), prow<i
1062: */
1063: #define PetscILULLAddSorted(nidx,idx,level,idxlvl,idx_start,nlnk,lnk,lnklvl,bt,lnklvl_prow) 0;\
1064: {\
1065: PetscInt _k,_entry,_location,_lnkdata,_incrlev,_lnklvl_prow=lnklvl[prow];\
1066: nlnk = 0;\
1067: _lnkdata = idx_start;\
1068: for (_k=0; _k<nidx; _k++){\
1069: _incrlev = idxlvl[_k] + _lnklvl_prow + 1;\
1070: if (_incrlev > level) continue;\
1071: _entry = idx[_k];\
1072: if (!PetscBTLookupSet(bt,_entry)){ /* new entry */\
1073: /* search for insertion location */\
1074: do {\
1075: _location = _lnkdata;\
1076: _lnkdata = lnk[_location];\
1077: } while (_entry > _lnkdata);\
1078: /* insertion location is found, add entry into lnk */\
1079: lnk[_location] = _entry;\
1080: lnk[_entry] = _lnkdata;\
1081: lnklvl[_entry] = _incrlev;\
1082: nlnk++;\
1083: _lnkdata = _entry; /* next search starts from here if next_entry > _entry */\
1084: } else { /* existing entry: update lnklvl */\
1085: if (lnklvl[_entry] > _incrlev) lnklvl[_entry] = _incrlev;\
1086: }\
1087: }\
1088: }
1090: /*
1091: Add a index set into a sorted linked list
1092: Input Parameters:
1093: nidx - number of input idx
1094: idx - integer array used for storing column indices
1095: level - level of fill, e.g., ICC(level)
1096: idxlvl - level of idx
1097: idx_start - starting index of the list
1098: lnk - linked list(an integer array) that is created
1099: lnklvl - levels of lnk
1100: bt - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
1101: output Parameters:
1102: nlnk - number of newly added idx
1103: lnk - the sorted(increasing order) linked list containing new and non-redundate entries from idx
1104: lnklvl - levels of lnk
1105: bt - updated PetscBT (bitarray)
1106: */
1107: #define PetscIncompleteLLAdd(nidx,idx,level,idxlvl,idx_start,nlnk,lnk,lnklvl,bt) 0;\
1108: {\
1109: PetscInt _k,_entry,_location,_lnkdata,_incrlev;\
1110: nlnk = 0;\
1111: _lnkdata = idx_start;\
1112: for (_k=0; _k<nidx; _k++){\
1113: _incrlev = idxlvl[_k] + 1;\
1114: if (_incrlev > level) continue;\
1115: _entry = idx[_k];\
1116: if (!PetscBTLookupSet(bt,_entry)){ /* new entry */\
1117: /* search for insertion location */\
1118: if (_k && _entry < _lnkdata) _lnkdata = idx_start;\
1119: do {\
1120: _location = _lnkdata;\
1121: _lnkdata = lnk[_location];\
1122: } while (_entry > _lnkdata);\
1123: /* insertion location is found, add entry into lnk */\
1124: lnk[_location] = _entry;\
1125: lnk[_entry] = _lnkdata;\
1126: lnklvl[_entry] = _incrlev;\
1127: nlnk++;\
1128: _lnkdata = _entry; /* next search starts from here if next_entry > _entry */\
1129: } else { /* existing entry: update lnklvl */\
1130: if (lnklvl[_entry] > _incrlev) lnklvl[_entry] = _incrlev;\
1131: }\
1132: }\
1133: }
1135: /*
1136: Add a SORTED index set into a sorted linked list
1137: Input Parameters:
1138: nidx - number of input indices
1139: idx - sorted integer array used for storing column indices
1140: level - level of fill, e.g., ICC(level)
1141: idxlvl - level of idx
1142: idx_start - starting index of the list
1143: lnk - linked list(an integer array) that is created
1144: lnklvl - levels of lnk
1145: bt - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
1146: output Parameters:
1147: nlnk - number of newly added idx
1148: lnk - the sorted(increasing order) linked list containing new and non-redundate entries from idx
1149: lnklvl - levels of lnk
1150: bt - updated PetscBT (bitarray)
1151: */
1152: #define PetscIncompleteLLAddSorted(nidx,idx,level,idxlvl,idx_start,nlnk,lnk,lnklvl,bt) 0;\
1153: {\
1154: PetscInt _k,_entry,_location,_lnkdata,_incrlev;\
1155: nlnk = 0;\
1156: _lnkdata = idx_start;\
1157: for (_k=0; _k<nidx; _k++){\
1158: _incrlev = idxlvl[_k] + 1;\
1159: if (_incrlev > level) continue;\
1160: _entry = idx[_k];\
1161: if (!PetscBTLookupSet(bt,_entry)){ /* new entry */\
1162: /* search for insertion location */\
1163: do {\
1164: _location = _lnkdata;\
1165: _lnkdata = lnk[_location];\
1166: } while (_entry > _lnkdata);\
1167: /* insertion location is found, add entry into lnk */\
1168: lnk[_location] = _entry;\
1169: lnk[_entry] = _lnkdata;\
1170: lnklvl[_entry] = _incrlev;\
1171: nlnk++;\
1172: _lnkdata = _entry; /* next search starts from here */\
1173: } else { /* existing entry: update lnklvl */\
1174: if (lnklvl[_entry] > _incrlev) lnklvl[_entry] = _incrlev;\
1175: }\
1176: }\
1177: }
1179: /*
1180: Add a SORTED index set into a sorted linked list for ICC
1181: Input Parameters:
1182: nidx - number of input indices
1183: idx - sorted integer array used for storing column indices
1184: level - level of fill, e.g., ICC(level)
1185: idxlvl - level of idx
1186: idx_start - starting index of the list
1187: lnk - linked list(an integer array) that is created
1188: lnklvl - levels of lnk
1189: bt - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
1190: idxlvl_prow - idxlvl[prow], where prow is the row number of the idx
1191: output Parameters:
1192: nlnk - number of newly added indices
1193: lnk - the sorted(increasing order) linked list containing new and non-redundate entries from idx
1194: lnklvl - levels of lnk
1195: bt - updated PetscBT (bitarray)
1196: Note: the level of U(i,j) is set as lvl(i,j) = min{ lvl(i,j), lvl(prow,i)+lvl(prow,j)+1)
1197: where idx = non-zero columns of U(prow,prow+1:n-1), prow<i
1198: */
1199: #define PetscICCLLAddSorted(nidx,idx,level,idxlvl,idx_start,nlnk,lnk,lnklvl,bt,idxlvl_prow) 0;\
1200: {\
1201: PetscInt _k,_entry,_location,_lnkdata,_incrlev;\
1202: nlnk = 0;\
1203: _lnkdata = idx_start;\
1204: for (_k=0; _k<nidx; _k++){\
1205: _incrlev = idxlvl[_k] + idxlvl_prow + 1;\
1206: if (_incrlev > level) continue;\
1207: _entry = idx[_k];\
1208: if (!PetscBTLookupSet(bt,_entry)){ /* new entry */\
1209: /* search for insertion location */\
1210: do {\
1211: _location = _lnkdata;\
1212: _lnkdata = lnk[_location];\
1213: } while (_entry > _lnkdata);\
1214: /* insertion location is found, add entry into lnk */\
1215: lnk[_location] = _entry;\
1216: lnk[_entry] = _lnkdata;\
1217: lnklvl[_entry] = _incrlev;\
1218: nlnk++;\
1219: _lnkdata = _entry; /* next search starts from here */\
1220: } else { /* existing entry: update lnklvl */\
1221: if (lnklvl[_entry] > _incrlev) lnklvl[_entry] = _incrlev;\
1222: }\
1223: }\
1224: }
1226: /*
1227: Copy data on the list into an array, then initialize the list
1228: Input Parameters:
1229: idx_start - starting index of the list
1230: lnk_max - max value of lnk indicating the end of the list
1231: nlnk - number of data on the list to be copied
1232: lnk - linked list
1233: lnklvl - level of lnk
1234: bt - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
1235: output Parameters:
1236: indices - array that contains the copied data
1237: lnk - linked list that is cleaned and initialize
1238: lnklvl - level of lnk that is reinitialized
1239: bt - PetscBT (bitarray) with all bits set to false
1240: */
1241: #define PetscIncompleteLLClean(idx_start,lnk_max,nlnk,lnk,lnklvl,indices,indiceslvl,bt) 0;\
1242: {\
1243: PetscInt _j,_idx=idx_start;\
1244: for (_j=0; _j<nlnk; _j++){\
1245: _idx = lnk[_idx];\
1246: *(indices+_j) = _idx;\
1247: *(indiceslvl+_j) = lnklvl[_idx];\
1248: lnklvl[_idx] = -1;\
1249: PetscBTClear(bt,_idx);\
1250: }\
1251: lnk[idx_start] = lnk_max;\
1252: }
1253: /*
1254: Free memories used by the list
1255: */
1256: #define PetscIncompleteLLDestroy(lnk,bt) (PetscFree(lnk) || PetscBTDestroy(&(bt)))
1258: #define MatCheckSameLocalSize(A,ar1,B,ar2) \
1260: if ((A->rmap->n != B->rmap->n) || (A->cmap->n != B->cmap->n)) SETERRQ6(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Incompatible matrix local sizes: parameter # %d (%D x %D) != parameter # %d (%D x %D)",ar1,A->rmap->n,A->cmap->n,ar2,B->rmap->n,B->cmap->n);
1261:
1262: #define MatCheckSameSize(A,ar1,B,ar2) \
1263: if ((A->rmap->N != B->rmap->N) || (A->cmap->N != B->cmap->N)) SETERRQ6(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_INCOMP,"Incompatible matrix global sizes: parameter # %d (%D x %D) != parameter # %d (%D x %D)",ar1,A->rmap->N,A->cmap->N,ar2,B->rmap->N,B->cmap->N);\
1264: MatCheckSameLocalSize(A,ar1,B,ar2);
1265:
1266: #define VecCheckMatCompatible(M,x,ar1,b,ar2) \
1267: if (M->cmap->N != x->map->N) SETERRQ3(PetscObjectComm((PetscObject)M),PETSC_ERR_ARG_SIZ,"Vector global length incompatible with matrix: parameter # %d global size %D != matrix column global size %D",ar1,x->map->N,M->cmap->N);\
1268: if (M->rmap->N != b->map->N) SETERRQ3(PetscObjectComm((PetscObject)M),PETSC_ERR_ARG_SIZ,"Vector global length incompatible with matrix: parameter # %d global size %D != matrix row global size %D",ar2,b->map->N,M->rmap->N);
1270: /* -------------------------------------------------------------------------------------------------------*/
1271: #include <petscbt.h>
1272: /*
1273: Create and initialize a condensed linked list -
1274: same as PetscLLCreate(), but uses a scalable array 'lnk' with size of max number of entries, not O(N).
1275: Barry suggested this approach (Dec. 6, 2011):
1276: I've thought of an alternative way of representing a linked list that is efficient but doesn't have the O(N) scaling issue
1277: (it may be faster than the O(N) even sequentially due to less crazy memory access).
1279: 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
1280: 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
1281: 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
1282: 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.
1283: 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
1284: to each other so memory access is much better than using the big array.
1286: Example:
1287: nlnk_max=5, lnk_max=36:
1288: Initial list: [0, 0 | 36, 2 | 0, 0 | 0, 0 | 0, 0 | 0, 0 | 0, 0]
1289: here, head_node has index 2 with value lnk[2]=lnk_max=36,
1290: 0-th entry is used to store the number of entries in the list,
1291: The initial lnk represents head -> tail(marked by 36) with number of entries = lnk[0]=0.
1293: Now adding a sorted set {2,4}, the list becomes
1294: [2, 0 | 36, 4 |2, 6 | 4, 2 | 0, 0 | 0, 0 | 0, 0 ]
1295: represents head -> 2 -> 4 -> tail with number of entries = lnk[0]=2.
1297: Then adding a sorted set {0,3,35}, the list
1298: [5, 0 | 36, 8 | 2, 10 | 4, 12 | 0, 4 | 3, 6 | 35, 2 ]
1299: represents head -> 0 -> 2 -> 3 -> 4 -> 35 -> tail with number of entries = lnk[0]=5.
1301: Input Parameters:
1302: nlnk_max - max length of the list
1303: lnk_max - max value of the entries
1304: Output Parameters:
1305: lnk - list created and initialized
1306: bt - PetscBT (bitarray) with all bits set to false. Note: bt has size lnk_max, not nln_max!
1307: */
1308: PETSC_STATIC_INLINE PetscErrorCode PetscLLCondensedCreate(PetscInt nlnk_max,PetscInt lnk_max,PetscInt **lnk,PetscBT *bt)
1309: {
1311: PetscInt *llnk,lsize = 0;
1314: PetscIntMultError(2,nlnk_max+2,&lsize);
1315: PetscMalloc1(lsize,lnk);
1316: PetscBTCreate(lnk_max,bt);
1317: llnk = *lnk;
1318: llnk[0] = 0; /* number of entries on the list */
1319: llnk[2] = lnk_max; /* value in the head node */
1320: llnk[3] = 2; /* next for the head node */
1321: return(0);
1322: }
1324: /*
1325: Add a SORTED ascending index set into a sorted linked list. See PetscLLCondensedCreate() for detailed description.
1326: Input Parameters:
1327: nidx - number of input indices
1328: indices - sorted integer array
1329: lnk - condensed linked list(an integer array) that is created
1330: bt - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
1331: output Parameters:
1332: lnk - the sorted(increasing order) linked list containing previous and newly added non-redundate indices
1333: bt - updated PetscBT (bitarray)
1334: */
1335: PETSC_STATIC_INLINE PetscErrorCode PetscLLCondensedAddSorted(PetscInt nidx,const PetscInt indices[],PetscInt lnk[],PetscBT bt)
1336: {
1337: PetscInt _k,_entry,_location,_next,_lnkdata,_nlnk,_newnode;
1340: _nlnk = lnk[0]; /* num of entries on the input lnk */
1341: _location = 2; /* head */
1342: for (_k=0; _k<nidx; _k++){
1343: _entry = indices[_k];
1344: if (!PetscBTLookupSet(bt,_entry)){ /* new entry */
1345: /* search for insertion location */
1346: do {
1347: _next = _location + 1; /* link from previous node to next node */
1348: _location = lnk[_next]; /* idx of next node */
1349: _lnkdata = lnk[_location];/* value of next node */
1350: } while (_entry > _lnkdata);
1351: /* insertion location is found, add entry into lnk */
1352: _newnode = 2*(_nlnk+2); /* index for this new node */
1353: lnk[_next] = _newnode; /* connect previous node to the new node */
1354: lnk[_newnode] = _entry; /* set value of the new node */
1355: lnk[_newnode+1] = _location; /* connect new node to next node */
1356: _location = _newnode; /* next search starts from the new node */
1357: _nlnk++;
1358: } \
1359: }\
1360: lnk[0] = _nlnk; /* number of entries in the list */
1361: return(0);
1362: }
1364: PETSC_STATIC_INLINE PetscErrorCode PetscLLCondensedClean(PetscInt lnk_max,PetscInt nidx,PetscInt *indices,PetscInt lnk[],PetscBT bt)
1365: {
1367: PetscInt _k,_next,_nlnk;
1370: _next = lnk[3]; /* head node */
1371: _nlnk = lnk[0]; /* num of entries on the list */
1372: for (_k=0; _k<_nlnk; _k++){
1373: indices[_k] = lnk[_next];
1374: _next = lnk[_next + 1];
1375: PetscBTClear(bt,indices[_k]);
1376: }
1377: lnk[0] = 0; /* num of entries on the list */
1378: lnk[2] = lnk_max; /* initialize head node */
1379: lnk[3] = 2; /* head node */
1380: return(0);
1381: }
1383: PETSC_STATIC_INLINE PetscErrorCode PetscLLCondensedView(PetscInt *lnk)
1384: {
1386: PetscInt k;
1389: PetscPrintf(PETSC_COMM_SELF,"LLCondensed of size %D, (val, next)\n",lnk[0]);
1390: for (k=2; k< lnk[0]+2; k++){
1391: PetscPrintf(PETSC_COMM_SELF," %D: (%D, %D)\n",2*k,lnk[2*k],lnk[2*k+1]);
1392: }
1393: return(0);
1394: }
1396: /*
1397: Free memories used by the list
1398: */
1399: PETSC_STATIC_INLINE PetscErrorCode PetscLLCondensedDestroy(PetscInt *lnk,PetscBT bt)
1400: {
1404: PetscFree(lnk);
1405: PetscBTDestroy(&bt);
1406: return(0);
1407: }
1409: /* -------------------------------------------------------------------------------------------------------*/
1410: /*
1411: Same as PetscLLCondensedCreate(), but does not use non-scalable O(lnk_max) bitarray
1412: Input Parameters:
1413: nlnk_max - max length of the list
1414: Output Parameters:
1415: lnk - list created and initialized
1416: */
1417: PETSC_STATIC_INLINE PetscErrorCode PetscLLCondensedCreate_Scalable(PetscInt nlnk_max,PetscInt **lnk)
1418: {
1420: PetscInt *llnk,lsize = 0;
1423: PetscIntMultError(2,nlnk_max+2,&lsize);
1424: PetscMalloc1(lsize,lnk);
1425: llnk = *lnk;
1426: llnk[0] = 0; /* number of entries on the list */
1427: llnk[2] = PETSC_MAX_INT; /* value in the head node */
1428: llnk[3] = 2; /* next for the head node */
1429: return(0);
1430: }
1432: PETSC_STATIC_INLINE PetscErrorCode PetscLLCondensedExpand_Scalable(PetscInt nlnk_max,PetscInt **lnk)
1433: {
1435: PetscInt lsize = 0;
1438: PetscIntMultError(2,nlnk_max+2,&lsize);
1439: PetscRealloc(lsize*sizeof(PetscInt),lnk);
1440: return(0);
1441: }
1443: PETSC_STATIC_INLINE PetscErrorCode PetscLLCondensedAddSorted_Scalable(PetscInt nidx,const PetscInt indices[],PetscInt lnk[])
1444: {
1445: PetscInt _k,_entry,_location,_next,_lnkdata,_nlnk,_newnode;
1446: _nlnk = lnk[0]; /* num of entries on the input lnk */
1447: _location = 2; /* head */ \
1448: for (_k=0; _k<nidx; _k++){
1449: _entry = indices[_k];
1450: /* search for insertion location */
1451: do {
1452: _next = _location + 1; /* link from previous node to next node */
1453: _location = lnk[_next]; /* idx of next node */
1454: _lnkdata = lnk[_location];/* value of next node */
1455: } while (_entry > _lnkdata);
1456: if (_entry < _lnkdata) {
1457: /* insertion location is found, add entry into lnk */
1458: _newnode = 2*(_nlnk+2); /* index for this new node */
1459: lnk[_next] = _newnode; /* connect previous node to the new node */
1460: lnk[_newnode] = _entry; /* set value of the new node */
1461: lnk[_newnode+1] = _location; /* connect new node to next node */
1462: _location = _newnode; /* next search starts from the new node */
1463: _nlnk++;
1464: }
1465: }
1466: lnk[0] = _nlnk; /* number of entries in the list */
1467: return 0;
1468: }
1470: PETSC_STATIC_INLINE PetscErrorCode PetscLLCondensedClean_Scalable(PetscInt nidx,PetscInt *indices,PetscInt *lnk)
1471: {
1472: PetscInt _k,_next,_nlnk;
1473: _next = lnk[3]; /* head node */
1474: _nlnk = lnk[0];
1475: for (_k=0; _k<_nlnk; _k++){
1476: indices[_k] = lnk[_next];
1477: _next = lnk[_next + 1];
1478: }
1479: lnk[0] = 0; /* num of entries on the list */
1480: lnk[3] = 2; /* head node */
1481: return 0;
1482: }
1484: PETSC_STATIC_INLINE PetscErrorCode PetscLLCondensedDestroy_Scalable(PetscInt *lnk)
1485: {
1486: return PetscFree(lnk);
1487: }
1489: /* -------------------------------------------------------------------------------------------------------*/
1490: /*
1491: lnk[0] number of links
1492: lnk[1] number of entries
1493: lnk[3n] value
1494: lnk[3n+1] len
1495: lnk[3n+2] link to next value
1497: The next three are always the first link
1499: lnk[3] PETSC_MIN_INT+1
1500: lnk[4] 1
1501: lnk[5] link to first real entry
1503: The next three are always the last link
1505: lnk[6] PETSC_MAX_INT - 1
1506: lnk[7] 1
1507: lnk[8] next valid link (this is the same as lnk[0] but without the decreases)
1508: */
1510: PETSC_STATIC_INLINE PetscErrorCode PetscLLCondensedCreate_fast(PetscInt nlnk_max,PetscInt **lnk)
1511: {
1513: PetscInt *llnk,lsize = 0;
1516: PetscIntMultError(3,nlnk_max+3,&lsize);
1517: PetscMalloc1(lsize,lnk);
1518: llnk = *lnk;
1519: llnk[0] = 0; /* nlnk: number of entries on the list */
1520: llnk[1] = 0; /* number of integer entries represented in list */
1521: llnk[3] = PETSC_MIN_INT+1; /* value in the first node */
1522: llnk[4] = 1; /* count for the first node */
1523: llnk[5] = 6; /* next for the first node */
1524: llnk[6] = PETSC_MAX_INT-1; /* value in the last node */
1525: llnk[7] = 1; /* count for the last node */
1526: llnk[8] = 0; /* next valid node to be used */
1527: return(0);
1528: }
1530: PETSC_STATIC_INLINE PetscErrorCode PetscLLCondensedAddSorted_fast(PetscInt nidx,const PetscInt indices[],PetscInt lnk[])
1531: {
1532: PetscInt k,entry,prev,next;
1533: prev = 3; /* first value */
1534: next = lnk[prev+2];
1535: for (k=0; k<nidx; k++){
1536: entry = indices[k];
1537: /* search for insertion location */
1538: while (entry >= lnk[next]) {
1539: prev = next;
1540: next = lnk[next+2];
1541: }
1542: /* entry is in range of previous list */
1543: if (entry < lnk[prev]+lnk[prev+1]) continue;
1544: lnk[1]++;
1545: /* entry is right after previous list */
1546: if (entry == lnk[prev]+lnk[prev+1]) {
1547: lnk[prev+1]++;
1548: if (lnk[next] == entry+1) { /* combine two contiguous strings */
1549: lnk[prev+1] += lnk[next+1];
1550: lnk[prev+2] = lnk[next+2];
1551: next = lnk[next+2];
1552: lnk[0]--;
1553: }
1554: continue;
1555: }
1556: /* entry is right before next list */
1557: if (entry == lnk[next]-1) {
1558: lnk[next]--;
1559: lnk[next+1]++;
1560: prev = next;
1561: next = lnk[prev+2];
1562: continue;
1563: }
1564: /* add entry into lnk */
1565: lnk[prev+2] = 3*((lnk[8]++)+3); /* connect previous node to the new node */
1566: prev = lnk[prev+2];
1567: lnk[prev] = entry; /* set value of the new node */
1568: lnk[prev+1] = 1; /* number of values in contiguous string is one to start */
1569: lnk[prev+2] = next; /* connect new node to next node */
1570: lnk[0]++;
1571: }
1572: return 0;
1573: }
1575: PETSC_STATIC_INLINE PetscErrorCode PetscLLCondensedClean_fast(PetscInt nidx,PetscInt *indices,PetscInt *lnk)
1576: {
1577: PetscInt _k,_next,_nlnk,cnt,j;
1578: _next = lnk[5]; /* first node */
1579: _nlnk = lnk[0];
1580: cnt = 0;
1581: for (_k=0; _k<_nlnk; _k++){
1582: for (j=0; j<lnk[_next+1]; j++) {
1583: indices[cnt++] = lnk[_next] + j;
1584: }
1585: _next = lnk[_next + 2];
1586: }
1587: lnk[0] = 0; /* nlnk: number of links */
1588: lnk[1] = 0; /* number of integer entries represented in list */
1589: lnk[3] = PETSC_MIN_INT+1; /* value in the first node */
1590: lnk[4] = 1; /* count for the first node */
1591: lnk[5] = 6; /* next for the first node */
1592: lnk[6] = PETSC_MAX_INT-1; /* value in the last node */
1593: lnk[7] = 1; /* count for the last node */
1594: lnk[8] = 0; /* next valid location to make link */
1595: return 0;
1596: }
1598: PETSC_STATIC_INLINE PetscErrorCode PetscLLCondensedView_fast(PetscInt *lnk)
1599: {
1600: PetscInt k,next,nlnk;
1601: next = lnk[5]; /* first node */
1602: nlnk = lnk[0];
1603: for (k=0; k<nlnk; k++){
1604: #if 0 /* Debugging code */
1605: printf("%d value %d len %d next %d\n",next,lnk[next],lnk[next+1],lnk[next+2]);
1606: #endif
1607: next = lnk[next + 2];
1608: }
1609: return 0;
1610: }
1612: PETSC_STATIC_INLINE PetscErrorCode PetscLLCondensedDestroy_fast(PetscInt *lnk)
1613: {
1614: return PetscFree(lnk);
1615: }
1617: /* this is extern because it is used in MatFDColoringUseDM() which is in the DM library */
1618: PETSC_EXTERN PetscErrorCode MatFDColoringApply_AIJ(Mat,MatFDColoring,Vec,void*);
1620: PETSC_EXTERN PetscLogEvent MAT_Mult;
1621: PETSC_EXTERN PetscLogEvent MAT_MultMatrixFree;
1622: PETSC_EXTERN PetscLogEvent MAT_Mults;
1623: PETSC_EXTERN PetscLogEvent MAT_MultConstrained;
1624: PETSC_EXTERN PetscLogEvent MAT_MultAdd;
1625: PETSC_EXTERN PetscLogEvent MAT_MultTranspose;
1626: PETSC_EXTERN PetscLogEvent MAT_MultTransposeConstrained;
1627: PETSC_EXTERN PetscLogEvent MAT_MultTransposeAdd;
1628: PETSC_EXTERN PetscLogEvent MAT_Solve;
1629: PETSC_EXTERN PetscLogEvent MAT_Solves;
1630: PETSC_EXTERN PetscLogEvent MAT_SolveAdd;
1631: PETSC_EXTERN PetscLogEvent MAT_SolveTranspose;
1632: PETSC_EXTERN PetscLogEvent MAT_SolveTransposeAdd;
1633: PETSC_EXTERN PetscLogEvent MAT_SOR;
1634: PETSC_EXTERN PetscLogEvent MAT_ForwardSolve;
1635: PETSC_EXTERN PetscLogEvent MAT_BackwardSolve;
1636: PETSC_EXTERN PetscLogEvent MAT_LUFactor;
1637: PETSC_EXTERN PetscLogEvent MAT_LUFactorSymbolic;
1638: PETSC_EXTERN PetscLogEvent MAT_LUFactorNumeric;
1639: PETSC_EXTERN PetscLogEvent MAT_CholeskyFactor;
1640: PETSC_EXTERN PetscLogEvent MAT_CholeskyFactorSymbolic;
1641: PETSC_EXTERN PetscLogEvent MAT_CholeskyFactorNumeric;
1642: PETSC_EXTERN PetscLogEvent MAT_ILUFactor;
1643: PETSC_EXTERN PetscLogEvent MAT_ILUFactorSymbolic;
1644: PETSC_EXTERN PetscLogEvent MAT_ICCFactorSymbolic;
1645: PETSC_EXTERN PetscLogEvent MAT_Copy;
1646: PETSC_EXTERN PetscLogEvent MAT_Convert;
1647: PETSC_EXTERN PetscLogEvent MAT_Scale;
1648: PETSC_EXTERN PetscLogEvent MAT_AssemblyBegin;
1649: PETSC_EXTERN PetscLogEvent MAT_AssemblyEnd;
1650: PETSC_EXTERN PetscLogEvent MAT_SetValues;
1651: PETSC_EXTERN PetscLogEvent MAT_GetValues;
1652: PETSC_EXTERN PetscLogEvent MAT_GetRow;
1653: PETSC_EXTERN PetscLogEvent MAT_GetRowIJ;
1654: PETSC_EXTERN PetscLogEvent MAT_CreateSubMats;
1655: PETSC_EXTERN PetscLogEvent MAT_GetColoring;
1656: PETSC_EXTERN PetscLogEvent MAT_GetOrdering;
1657: PETSC_EXTERN PetscLogEvent MAT_RedundantMat;
1658: PETSC_EXTERN PetscLogEvent MAT_IncreaseOverlap;
1659: PETSC_EXTERN PetscLogEvent MAT_Partitioning;
1660: PETSC_EXTERN PetscLogEvent MAT_PartitioningND;
1661: PETSC_EXTERN PetscLogEvent MAT_Coarsen;
1662: PETSC_EXTERN PetscLogEvent MAT_ZeroEntries;
1663: PETSC_EXTERN PetscLogEvent MAT_Load;
1664: PETSC_EXTERN PetscLogEvent MAT_View;
1665: PETSC_EXTERN PetscLogEvent MAT_AXPY;
1666: PETSC_EXTERN PetscLogEvent MAT_FDColoringCreate;
1667: PETSC_EXTERN PetscLogEvent MAT_TransposeColoringCreate;
1668: PETSC_EXTERN PetscLogEvent MAT_FDColoringSetUp;
1669: PETSC_EXTERN PetscLogEvent MAT_FDColoringApply;
1670: PETSC_EXTERN PetscLogEvent MAT_Transpose;
1671: PETSC_EXTERN PetscLogEvent MAT_FDColoringFunction;
1672: PETSC_EXTERN PetscLogEvent MAT_CreateSubMat;
1673: PETSC_EXTERN PetscLogEvent MAT_MatMult;
1674: PETSC_EXTERN PetscLogEvent MAT_MatSolve;
1675: PETSC_EXTERN PetscLogEvent MAT_MatTrSolve;
1676: PETSC_EXTERN PetscLogEvent MAT_MatMultSymbolic;
1677: PETSC_EXTERN PetscLogEvent MAT_MatMultNumeric;
1678: PETSC_EXTERN PetscLogEvent MAT_Getlocalmatcondensed;
1679: PETSC_EXTERN PetscLogEvent MAT_GetBrowsOfAcols;
1680: PETSC_EXTERN PetscLogEvent MAT_GetBrowsOfAocols;
1681: PETSC_EXTERN PetscLogEvent MAT_PtAP;
1682: PETSC_EXTERN PetscLogEvent MAT_PtAPSymbolic;
1683: PETSC_EXTERN PetscLogEvent MAT_PtAPNumeric;
1684: PETSC_EXTERN PetscLogEvent MAT_Seqstompinum;
1685: PETSC_EXTERN PetscLogEvent MAT_Seqstompisym;
1686: PETSC_EXTERN PetscLogEvent MAT_Seqstompi;
1687: PETSC_EXTERN PetscLogEvent MAT_Getlocalmat;
1688: PETSC_EXTERN PetscLogEvent MAT_RARt;
1689: PETSC_EXTERN PetscLogEvent MAT_RARtSymbolic;
1690: PETSC_EXTERN PetscLogEvent MAT_RARtNumeric;
1691: PETSC_EXTERN PetscLogEvent MAT_MatTransposeMult;
1692: PETSC_EXTERN PetscLogEvent MAT_MatTransposeMultSymbolic;
1693: PETSC_EXTERN PetscLogEvent MAT_MatTransposeMultNumeric;
1694: PETSC_EXTERN PetscLogEvent MAT_TransposeMatMult;
1695: PETSC_EXTERN PetscLogEvent MAT_TransposeMatMultSymbolic;
1696: PETSC_EXTERN PetscLogEvent MAT_TransposeMatMultNumeric;
1697: PETSC_EXTERN PetscLogEvent MAT_MatMatMult;
1698: PETSC_EXTERN PetscLogEvent MAT_MatMatMultSymbolic;
1699: PETSC_EXTERN PetscLogEvent MAT_MatMatMultNumeric;
1700: PETSC_EXTERN PetscLogEvent MAT_Applypapt;
1701: PETSC_EXTERN PetscLogEvent MAT_Applypapt_symbolic;
1702: PETSC_EXTERN PetscLogEvent MAT_Applypapt_numeric;
1703: PETSC_EXTERN PetscLogEvent MAT_Getsymtranspose;
1704: PETSC_EXTERN PetscLogEvent MAT_Transpose_SeqAIJ;
1705: PETSC_EXTERN PetscLogEvent MAT_Getsymtransreduced;
1706: PETSC_EXTERN PetscLogEvent MAT_GetSequentialNonzeroStructure;
1707: PETSC_EXTERN PetscLogEvent MATMFFD_Mult;
1708: PETSC_EXTERN PetscLogEvent MAT_GetMultiProcBlock;
1709: PETSC_EXTERN PetscLogEvent MAT_CUSPARSECopyToGPU;
1710: PETSC_EXTERN PetscLogEvent MAT_SetValuesBatch;
1711: PETSC_EXTERN PetscLogEvent MAT_ViennaCLCopyToGPU;
1712: PETSC_EXTERN PetscLogEvent MAT_Merge;
1713: PETSC_EXTERN PetscLogEvent MAT_Residual;
1714: PETSC_EXTERN PetscLogEvent MAT_SetRandom;
1715: PETSC_EXTERN PetscLogEvent MATCOLORING_Apply;
1716: PETSC_EXTERN PetscLogEvent MATCOLORING_Comm;
1717: PETSC_EXTERN PetscLogEvent MATCOLORING_Local;
1718: PETSC_EXTERN PetscLogEvent MATCOLORING_ISCreate;
1719: PETSC_EXTERN PetscLogEvent MATCOLORING_SetUp;
1720: PETSC_EXTERN PetscLogEvent MATCOLORING_Weights;
1722: #endif