Actual source code: matimpl.h
petsc-3.11.1 2019-04-12
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 (*freeintermediatedatastructures)(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_MatRootName* MatRootName;
221: struct _p_MatRootName {
222: char *rname,*sname,*mname;
223: MatRootName next;
224: };
226: PETSC_EXTERN MatRootName MatRootNameList;
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_CUDA)
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);
420: PETSC_INTERN PetscErrorCode MatAXPY_Basic_Preallocate(Mat,Mat,Mat*);
422: /*
423: Utility for MatFactor (Schur complement)
424: */
425: PETSC_INTERN PetscErrorCode MatFactorFactorizeSchurComplement_Private(Mat);
426: PETSC_INTERN PetscErrorCode MatFactorInvertSchurComplement_Private(Mat);
427: PETSC_INTERN PetscErrorCode MatFactorUpdateSchurStatus_Private(Mat);
428: PETSC_INTERN PetscErrorCode MatFactorSetUpInPlaceSchur_Private(Mat);
430: /*
431: Utility for MatZeroRows
432: */
433: PETSC_INTERN PetscErrorCode MatZeroRowsMapLocal_Private(Mat,PetscInt,const PetscInt*,PetscInt*,PetscInt**);
435: /*
436: Object for partitioning graphs
437: */
439: typedef struct _MatPartitioningOps *MatPartitioningOps;
440: struct _MatPartitioningOps {
441: PetscErrorCode (*apply)(MatPartitioning,IS*);
442: PetscErrorCode (*applynd)(MatPartitioning,IS*);
443: PetscErrorCode (*setfromoptions)(PetscOptionItems*,MatPartitioning);
444: PetscErrorCode (*destroy)(MatPartitioning);
445: PetscErrorCode (*view)(MatPartitioning,PetscViewer);
446: PetscErrorCode (*improve)(MatPartitioning,IS*);
447: };
449: struct _p_MatPartitioning {
450: PETSCHEADER(struct _MatPartitioningOps);
451: Mat adj;
452: PetscInt *vertex_weights;
453: PetscReal *part_weights;
454: PetscInt n; /* number of partitions */
455: void *data;
456: PetscInt setupcalled;
457: };
459: /* needed for parallel nested dissection by ParMetis and PTSCOTCH */
460: PETSC_INTERN PetscErrorCode MatPartitioningSizesToSep_Private(PetscInt,PetscInt[],PetscInt[],PetscInt[]);
462: /*
463: Object for coarsen graphs
464: */
465: typedef struct _MatCoarsenOps *MatCoarsenOps;
466: struct _MatCoarsenOps {
467: PetscErrorCode (*apply)(MatCoarsen);
468: PetscErrorCode (*setfromoptions)(PetscOptionItems*,MatCoarsen);
469: PetscErrorCode (*destroy)(MatCoarsen);
470: PetscErrorCode (*view)(MatCoarsen,PetscViewer);
471: };
473: struct _p_MatCoarsen {
474: PETSCHEADER(struct _MatCoarsenOps);
475: Mat graph;
476: PetscInt setupcalled;
477: void *subctx;
478: /* */
479: PetscBool strict_aggs;
480: IS perm;
481: PetscCoarsenData *agg_lists;
482: };
484: /*
485: MatFDColoring is used to compute Jacobian matrices efficiently
486: via coloring. The data structure is explained below in an example.
488: Color = 0 1 0 2 | 2 3 0
489: ---------------------------------------------------
490: 00 01 | 05
491: 10 11 | 14 15 Processor 0
492: 22 23 | 25
493: 32 33 |
494: ===================================================
495: | 44 45 46
496: 50 | 55 Processor 1
497: | 64 66
498: ---------------------------------------------------
500: ncolors = 4;
502: ncolumns = {2,1,1,0}
503: columns = {{0,2},{1},{3},{}}
504: nrows = {4,2,3,3}
505: rows = {{0,1,2,3},{0,1},{1,2,3},{0,1,2}}
506: vwscale = {dx(0),dx(1),dx(2),dx(3)} MPI Vec
507: vscale = {dx(0),dx(1),dx(2),dx(3),dx(4),dx(5)} Seq Vec
509: ncolumns = {1,0,1,1}
510: columns = {{6},{},{4},{5}}
511: nrows = {3,0,2,2}
512: rows = {{0,1,2},{},{1,2},{1,2}}
513: vwscale = {dx(4),dx(5),dx(6)} MPI Vec
514: vscale = {dx(0),dx(4),dx(5),dx(6)} Seq Vec
516: See the routine MatFDColoringApply() for how this data is used
517: to compute the Jacobian.
519: */
520: typedef struct {
521: PetscInt row;
522: PetscInt col;
523: PetscScalar *valaddr; /* address of value */
524: } MatEntry;
526: typedef struct {
527: PetscInt row;
528: PetscScalar *valaddr; /* address of value */
529: } MatEntry2;
531: struct _p_MatFDColoring{
532: PETSCHEADER(int);
533: PetscInt M,N,m; /* total rows, columns; local rows */
534: PetscInt rstart; /* first row owned by local processor */
535: PetscInt ncolors; /* number of colors */
536: PetscInt *ncolumns; /* number of local columns for a color */
537: PetscInt **columns; /* lists the local columns of each color (using global column numbering) */
538: PetscInt *nrows; /* number of local rows for each color */
539: MatEntry *matentry; /* holds (row, column, address of value) for Jacobian matrix entry */
540: MatEntry2 *matentry2; /* holds (row, address of value) for Jacobian matrix entry */
541: PetscScalar *dy; /* store a block of F(x+dx)-F(x) when J is in BAIJ format */
542: PetscReal error_rel; /* square root of relative error in computing function */
543: PetscReal umin; /* minimum allowable u'dx value */
544: Vec w1,w2,w3; /* work vectors used in computing Jacobian */
545: PetscBool fset; /* indicates that the initial function value F(X) is set */
546: PetscErrorCode (*f)(void); /* function that defines Jacobian */
547: void *fctx; /* optional user-defined context for use by the function f */
548: Vec vscale; /* holds FD scaling, i.e. 1/dx for each perturbed column */
549: PetscInt currentcolor; /* color for which function evaluation is being done now */
550: const char *htype; /* "wp" or "ds" */
551: ISColoringType ctype; /* IS_COLORING_GLOBAL or IS_COLORING_LOCAL */
552: PetscInt brows,bcols; /* number of block rows or columns for speedup inserting the dense matrix into sparse Jacobian */
553: PetscBool setupcalled; /* true if setup has been called */
554: PetscBool viewed; /* true if the -mat_fd_coloring_view has been triggered already */
555: void (*ftn_func_pointer)(void),*ftn_func_cntx; /* serve the same purpose as *fortran_func_pointers in PETSc objects */
556: };
558: typedef struct _MatColoringOps *MatColoringOps;
559: struct _MatColoringOps {
560: PetscErrorCode (*destroy)(MatColoring);
561: PetscErrorCode (*setfromoptions)(PetscOptionItems*,MatColoring);
562: PetscErrorCode (*view)(MatColoring,PetscViewer);
563: PetscErrorCode (*apply)(MatColoring,ISColoring*);
564: PetscErrorCode (*weights)(MatColoring,PetscReal**,PetscInt**);
565: };
567: struct _p_MatColoring {
568: PETSCHEADER(struct _MatColoringOps);
569: Mat mat;
570: PetscInt dist; /* distance of the coloring */
571: PetscInt maxcolors; /* the maximum number of colors returned, maxcolors=1 for MIS */
572: void *data; /* inner context */
573: PetscBool valid; /* check to see if what is produced is a valid coloring */
574: MatColoringWeightType weight_type; /* type of weight computation to be performed */
575: PetscReal *user_weights; /* custom weights and permutation */
576: PetscInt *user_lperm;
577: PetscBool valid_iscoloring; /* check to see if matcoloring is produced a valid iscoloring */
578: };
580: struct _p_MatTransposeColoring{
581: PETSCHEADER(int);
582: PetscInt M,N,m; /* total rows, columns; local rows */
583: PetscInt rstart; /* first row owned by local processor */
584: PetscInt ncolors; /* number of colors */
585: PetscInt *ncolumns; /* number of local columns for a color */
586: PetscInt *nrows; /* number of local rows for each color */
587: PetscInt currentcolor; /* color for which function evaluation is being done now */
588: ISColoringType ctype; /* IS_COLORING_GLOBAL or IS_COLORING_LOCAL */
590: PetscInt *colorforrow,*colorforcol; /* pointer to rows and columns */
591: PetscInt *rows; /* lists the local rows for each color (using the local row numbering) */
592: PetscInt *den2sp; /* maps (row,color) in the dense matrix to index of sparse matrix array a->a */
593: PetscInt *columns; /* lists the local columns of each color (using global column numbering) */
594: PetscInt brows; /* number of rows for efficient implementation of MatTransColoringApplyDenToSp() */
595: PetscInt *lstart; /* array used for loop over row blocks of Csparse */
596: };
598: /*
599: Null space context for preconditioner/operators
600: */
601: struct _p_MatNullSpace {
602: PETSCHEADER(int);
603: PetscBool has_cnst;
604: PetscInt n;
605: Vec* vecs;
606: PetscScalar* alpha; /* for projections */
607: PetscErrorCode (*remove)(MatNullSpace,Vec,void*); /* for user provided removal function */
608: void* rmctx; /* context for remove() function */
609: };
611: /*
612: Checking zero pivot for LU, ILU preconditioners.
613: */
614: typedef struct {
615: PetscInt nshift,nshift_max;
616: PetscReal shift_amount,shift_lo,shift_hi,shift_top,shift_fraction;
617: PetscBool newshift;
618: PetscReal rs; /* active row sum of abs(offdiagonals) */
619: PetscScalar pv; /* pivot of the active row */
620: } FactorShiftCtx;
622: /*
623: Used by MatCreateSubMatrices_MPIXAIJ_Local()
624: */
625: #include <petscctable.h>
626: typedef struct { /* used by MatCreateSubMatrices_MPIAIJ_SingleIS_Local() and MatCreateSubMatrices_MPIAIJ_Local */
627: PetscInt id; /* index of submats, only submats[0] is responsible for deleting some arrays below */
628: PetscInt nrqs,nrqr;
629: PetscInt **rbuf1,**rbuf2,**rbuf3,**sbuf1,**sbuf2;
630: PetscInt **ptr;
631: PetscInt *tmp;
632: PetscInt *ctr;
633: PetscInt *pa; /* proc array */
634: PetscInt *req_size,*req_source1,*req_source2;
635: PetscBool allcolumns,allrows;
636: PetscBool singleis;
637: PetscInt *row2proc; /* row to proc map */
638: PetscInt nstages;
639: #if defined(PETSC_USE_CTABLE)
640: PetscTable cmap,rmap;
641: PetscInt *cmap_loc,*rmap_loc;
642: #else
643: PetscInt *cmap,*rmap;
644: #endif
646: PetscErrorCode (*destroy)(Mat);
647: } Mat_SubSppt;
649: PETSC_EXTERN PetscErrorCode MatFactorDumpMatrix(Mat);
650: PETSC_INTERN PetscErrorCode MatShift_Basic(Mat,PetscScalar);
651: PETSC_INTERN PetscErrorCode MatSetBlockSizes_Default(Mat,PetscInt,PetscInt);
653: PETSC_STATIC_INLINE PetscErrorCode MatPivotCheck_nz(Mat mat,const MatFactorInfo *info,FactorShiftCtx *sctx,PetscInt row)
654: {
655: PetscReal _rs = sctx->rs;
656: PetscReal _zero = info->zeropivot*_rs;
659: if (PetscAbsScalar(sctx->pv) <= _zero && !PetscIsNanScalar(sctx->pv)){
660: /* force |diag| > zeropivot*rs */
661: if (!sctx->nshift) sctx->shift_amount = info->shiftamount;
662: else sctx->shift_amount *= 2.0;
663: sctx->newshift = PETSC_TRUE;
664: (sctx->nshift)++;
665: } else {
666: sctx->newshift = PETSC_FALSE;
667: }
668: return(0);
669: }
671: PETSC_STATIC_INLINE PetscErrorCode MatPivotCheck_pd(Mat mat,const MatFactorInfo *info,FactorShiftCtx *sctx,PetscInt row)
672: {
673: PetscReal _rs = sctx->rs;
674: PetscReal _zero = info->zeropivot*_rs;
677: if (PetscRealPart(sctx->pv) <= _zero && !PetscIsNanScalar(sctx->pv)){
678: /* force matfactor to be diagonally dominant */
679: if (sctx->nshift == sctx->nshift_max) {
680: sctx->shift_fraction = sctx->shift_hi;
681: } else {
682: sctx->shift_lo = sctx->shift_fraction;
683: sctx->shift_fraction = (sctx->shift_hi+sctx->shift_lo)/2.;
684: }
685: sctx->shift_amount = sctx->shift_fraction * sctx->shift_top;
686: sctx->nshift++;
687: sctx->newshift = PETSC_TRUE;
688: } else {
689: sctx->newshift = PETSC_FALSE;
690: }
691: return(0);
692: }
694: PETSC_STATIC_INLINE PetscErrorCode MatPivotCheck_inblocks(Mat mat,const MatFactorInfo *info,FactorShiftCtx *sctx,PetscInt row)
695: {
696: PetscReal _zero = info->zeropivot;
699: if (PetscAbsScalar(sctx->pv) <= _zero && !PetscIsNanScalar(sctx->pv)){
700: sctx->pv += info->shiftamount;
701: sctx->shift_amount = 0.0;
702: sctx->nshift++;
703: }
704: sctx->newshift = PETSC_FALSE;
705: return(0);
706: }
708: PETSC_STATIC_INLINE PetscErrorCode MatPivotCheck_none(Mat fact,Mat mat,const MatFactorInfo *info,FactorShiftCtx *sctx,PetscInt row)
709: {
710: PetscReal _zero = info->zeropivot;
714: sctx->newshift = PETSC_FALSE;
715: if (PetscAbsScalar(sctx->pv) <= _zero && !PetscIsNanScalar(sctx->pv)) {
716: if (!mat->erroriffailure) {
717: PetscInfo3(mat,"Detected zero pivot in factorization in row %D value %g tolerance %g\n",row,(double)PetscAbsScalar(sctx->pv),(double)_zero);
718: fact->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
719: fact->factorerror_zeropivot_value = PetscAbsScalar(sctx->pv);
720: fact->factorerror_zeropivot_row = row;
721: } 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);
722: }
723: return(0);
724: }
726: PETSC_STATIC_INLINE PetscErrorCode MatPivotCheck(Mat fact,Mat mat,const MatFactorInfo *info,FactorShiftCtx *sctx,PetscInt row)
727: {
731: if (info->shifttype == (PetscReal) MAT_SHIFT_NONZERO){
732: MatPivotCheck_nz(mat,info,sctx,row);
733: } else if (info->shifttype == (PetscReal) MAT_SHIFT_POSITIVE_DEFINITE){
734: MatPivotCheck_pd(mat,info,sctx,row);
735: } else if (info->shifttype == (PetscReal) MAT_SHIFT_INBLOCKS){
736: MatPivotCheck_inblocks(mat,info,sctx,row);
737: } else {
738: MatPivotCheck_none(fact,mat,info,sctx,row);
739: }
740: return(0);
741: }
743: /*
744: Create and initialize a linked list
745: Input Parameters:
746: idx_start - starting index of the list
747: lnk_max - max value of lnk indicating the end of the list
748: nlnk - max length of the list
749: Output Parameters:
750: lnk - list initialized
751: bt - PetscBT (bitarray) with all bits set to false
752: lnk_empty - flg indicating the list is empty
753: */
754: #define PetscLLCreate(idx_start,lnk_max,nlnk,lnk,bt) \
755: (PetscMalloc1(nlnk,&lnk) || PetscBTCreate(nlnk,&(bt)) || (lnk[idx_start] = lnk_max,0))
757: #define PetscLLCreate_new(idx_start,lnk_max,nlnk,lnk,bt,lnk_empty)\
758: (PetscMalloc1(nlnk,&lnk) || PetscBTCreate(nlnk,&(bt)) || (lnk_empty = PETSC_TRUE,0) ||(lnk[idx_start] = lnk_max,0))
760: /*
761: Add an index set into a sorted linked list
762: Input Parameters:
763: nidx - number of input indices
764: indices - integer array
765: idx_start - starting index of the list
766: lnk - linked list(an integer array) that is created
767: bt - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
768: output Parameters:
769: nlnk - number of newly added indices
770: lnk - the sorted(increasing order) linked list containing new and non-redundate entries from indices
771: bt - updated PetscBT (bitarray)
772: */
773: #define PetscLLAdd(nidx,indices,idx_start,nlnk,lnk,bt) 0;\
774: {\
775: PetscInt _k,_entry,_location,_lnkdata;\
776: nlnk = 0;\
777: _lnkdata = idx_start;\
778: for (_k=0; _k<nidx; _k++){\
779: _entry = indices[_k];\
780: if (!PetscBTLookupSet(bt,_entry)){ /* new entry */\
781: /* search for insertion location */\
782: /* start from the beginning if _entry < previous _entry */\
783: if (_k && _entry < _lnkdata) _lnkdata = idx_start;\
784: do {\
785: _location = _lnkdata;\
786: _lnkdata = lnk[_location];\
787: } while (_entry > _lnkdata);\
788: /* insertion location is found, add entry into lnk */\
789: lnk[_location] = _entry;\
790: lnk[_entry] = _lnkdata;\
791: nlnk++;\
792: _lnkdata = _entry; /* next search starts from here if next_entry > _entry */\
793: }\
794: }\
795: }
797: /*
798: Add a permuted index set into a sorted linked list
799: Input Parameters:
800: nidx - number of input indices
801: indices - integer array
802: perm - permutation of indices
803: idx_start - starting index of the list
804: lnk - linked list(an integer array) that is created
805: bt - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
806: output Parameters:
807: nlnk - number of newly added indices
808: lnk - the sorted(increasing order) linked list containing new and non-redundate entries from indices
809: bt - updated PetscBT (bitarray)
810: */
811: #define PetscLLAddPerm(nidx,indices,perm,idx_start,nlnk,lnk,bt) 0;\
812: {\
813: PetscInt _k,_entry,_location,_lnkdata;\
814: nlnk = 0;\
815: _lnkdata = idx_start;\
816: for (_k=0; _k<nidx; _k++){\
817: _entry = perm[indices[_k]];\
818: if (!PetscBTLookupSet(bt,_entry)){ /* new entry */\
819: /* search for insertion location */\
820: /* start from the beginning if _entry < previous _entry */\
821: if (_k && _entry < _lnkdata) _lnkdata = idx_start;\
822: do {\
823: _location = _lnkdata;\
824: _lnkdata = lnk[_location];\
825: } while (_entry > _lnkdata);\
826: /* insertion location is found, add entry into lnk */\
827: lnk[_location] = _entry;\
828: lnk[_entry] = _lnkdata;\
829: nlnk++;\
830: _lnkdata = _entry; /* next search starts from here if next_entry > _entry */\
831: }\
832: }\
833: }
835: /*
836: Add a SORTED ascending index set into a sorted linked list - same as PetscLLAdd() bus skip 'if (_k && _entry < _lnkdata) _lnkdata = idx_start;'
837: Input Parameters:
838: nidx - number of input indices
839: indices - sorted integer array
840: idx_start - starting index of the list
841: lnk - linked list(an integer array) that is created
842: bt - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
843: output Parameters:
844: nlnk - number of newly added indices
845: lnk - the sorted(increasing order) linked list containing new and non-redundate entries from indices
846: bt - updated PetscBT (bitarray)
847: */
848: #define PetscLLAddSorted(nidx,indices,idx_start,nlnk,lnk,bt) 0;\
849: {\
850: PetscInt _k,_entry,_location,_lnkdata;\
851: nlnk = 0;\
852: _lnkdata = idx_start;\
853: for (_k=0; _k<nidx; _k++){\
854: _entry = indices[_k];\
855: if (!PetscBTLookupSet(bt,_entry)){ /* new entry */\
856: /* search for insertion location */\
857: do {\
858: _location = _lnkdata;\
859: _lnkdata = lnk[_location];\
860: } while (_entry > _lnkdata);\
861: /* insertion location is found, add entry into lnk */\
862: lnk[_location] = _entry;\
863: lnk[_entry] = _lnkdata;\
864: nlnk++;\
865: _lnkdata = _entry; /* next search starts from here */\
866: }\
867: }\
868: }
870: #define PetscLLAddSorted_new(nidx,indices,idx_start,lnk_empty,nlnk,lnk,bt) 0; \
871: {\
872: PetscInt _k,_entry,_location,_lnkdata;\
873: if (lnk_empty){\
874: _lnkdata = idx_start; \
875: for (_k=0; _k<nidx; _k++){ \
876: _entry = indices[_k]; \
877: PetscBTSet(bt,_entry); /* mark the new entry */ \
878: _location = _lnkdata; \
879: _lnkdata = lnk[_location]; \
880: /* insertion location is found, add entry into lnk */ \
881: lnk[_location] = _entry; \
882: lnk[_entry] = _lnkdata; \
883: _lnkdata = _entry; /* next search starts from here */ \
884: } \
885: /*\
886: lnk[indices[nidx-1]] = lnk[idx_start];\
887: lnk[idx_start] = indices[0];\
888: PetscBTSet(bt,indices[0]); \
889: for (_k=1; _k<nidx; _k++){ \
890: PetscBTSet(bt,indices[_k]); \
891: lnk[indices[_k-1]] = indices[_k]; \
892: } \
893: */\
894: nlnk = nidx;\
895: lnk_empty = PETSC_FALSE;\
896: } else {\
897: nlnk = 0; \
898: _lnkdata = idx_start; \
899: for (_k=0; _k<nidx; _k++){ \
900: _entry = indices[_k]; \
901: if (!PetscBTLookupSet(bt,_entry)){ /* new entry */ \
902: /* search for insertion location */ \
903: do { \
904: _location = _lnkdata; \
905: _lnkdata = lnk[_location]; \
906: } while (_entry > _lnkdata); \
907: /* insertion location is found, add entry into lnk */ \
908: lnk[_location] = _entry; \
909: lnk[_entry] = _lnkdata; \
910: nlnk++; \
911: _lnkdata = _entry; /* next search starts from here */ \
912: } \
913: } \
914: } \
915: }
917: /*
918: Add a SORTED index set into a sorted linked list used for LUFactorSymbolic()
919: Same as PetscLLAddSorted() with an additional operation:
920: count the number of input indices that are no larger than 'diag'
921: Input Parameters:
922: indices - sorted integer array
923: idx_start - starting index of the list, index of pivot row
924: lnk - linked list(an integer array) that is created
925: bt - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
926: diag - index of the active row in LUFactorSymbolic
927: nzbd - number of input indices with indices <= idx_start
928: im - im[idx_start] is initialized as num of nonzero entries in row=idx_start
929: output Parameters:
930: nlnk - number of newly added indices
931: lnk - the sorted(increasing order) linked list containing new and non-redundate entries from indices
932: bt - updated PetscBT (bitarray)
933: im - im[idx_start]: unchanged if diag is not an entry
934: : num of entries with indices <= diag if diag is an entry
935: */
936: #define PetscLLAddSortedLU(indices,idx_start,nlnk,lnk,bt,diag,nzbd,im) 0;\
937: {\
938: PetscInt _k,_entry,_location,_lnkdata,_nidx;\
939: nlnk = 0;\
940: _lnkdata = idx_start;\
941: _nidx = im[idx_start] - nzbd; /* num of entries with idx_start < index <= diag */\
942: for (_k=0; _k<_nidx; _k++){\
943: _entry = indices[_k];\
944: nzbd++;\
945: if ( _entry== diag) im[idx_start] = nzbd;\
946: if (!PetscBTLookupSet(bt,_entry)){ /* new entry */\
947: /* search for insertion location */\
948: do {\
949: _location = _lnkdata;\
950: _lnkdata = lnk[_location];\
951: } while (_entry > _lnkdata);\
952: /* insertion location is found, add entry into lnk */\
953: lnk[_location] = _entry;\
954: lnk[_entry] = _lnkdata;\
955: nlnk++;\
956: _lnkdata = _entry; /* next search starts from here */\
957: }\
958: }\
959: }
961: /*
962: Copy data on the list into an array, then initialize the list
963: Input Parameters:
964: idx_start - starting index of the list
965: lnk_max - max value of lnk indicating the end of the list
966: nlnk - number of data on the list to be copied
967: lnk - linked list
968: bt - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
969: output Parameters:
970: indices - array that contains the copied data
971: lnk - linked list that is cleaned and initialize
972: bt - PetscBT (bitarray) with all bits set to false
973: */
974: #define PetscLLClean(idx_start,lnk_max,nlnk,lnk,indices,bt) 0;\
975: {\
976: PetscInt _j,_idx=idx_start;\
977: for (_j=0; _j<nlnk; _j++){\
978: _idx = lnk[_idx];\
979: indices[_j] = _idx;\
980: PetscBTClear(bt,_idx);\
981: }\
982: lnk[idx_start] = lnk_max;\
983: }
984: /*
985: Free memories used by the list
986: */
987: #define PetscLLDestroy(lnk,bt) (PetscFree(lnk) || PetscBTDestroy(&(bt)))
989: /* Routines below are used for incomplete matrix factorization */
990: /*
991: Create and initialize a linked list and its levels
992: Input Parameters:
993: idx_start - starting index of the list
994: lnk_max - max value of lnk indicating the end of the list
995: nlnk - max length of the list
996: Output Parameters:
997: lnk - list initialized
998: lnk_lvl - array of size nlnk for storing levels of lnk
999: bt - PetscBT (bitarray) with all bits set to false
1000: */
1001: #define PetscIncompleteLLCreate(idx_start,lnk_max,nlnk,lnk,lnk_lvl,bt)\
1002: (PetscIntMultError(2,nlnk,NULL) || PetscMalloc1(2*nlnk,&lnk) || PetscBTCreate(nlnk,&(bt)) || (lnk[idx_start] = lnk_max,lnk_lvl = lnk + nlnk,0))
1004: /*
1005: Initialize a sorted linked list used for ILU and ICC
1006: Input Parameters:
1007: nidx - number of input idx
1008: idx - integer array used for storing column indices
1009: idx_start - starting index of the list
1010: perm - indices of an IS
1011: lnk - linked list(an integer array) that is created
1012: lnklvl - levels of lnk
1013: bt - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
1014: output Parameters:
1015: nlnk - number of newly added idx
1016: lnk - the sorted(increasing order) linked list containing new and non-redundate entries from idx
1017: lnklvl - levels of lnk
1018: bt - updated PetscBT (bitarray)
1019: */
1020: #define PetscIncompleteLLInit(nidx,idx,idx_start,perm,nlnk,lnk,lnklvl,bt) 0;\
1021: {\
1022: PetscInt _k,_entry,_location,_lnkdata;\
1023: nlnk = 0;\
1024: _lnkdata = idx_start;\
1025: for (_k=0; _k<nidx; _k++){\
1026: _entry = perm[idx[_k]];\
1027: if (!PetscBTLookupSet(bt,_entry)){ /* new entry */\
1028: /* search for insertion location */\
1029: if (_k && _entry < _lnkdata) _lnkdata = idx_start;\
1030: do {\
1031: _location = _lnkdata;\
1032: _lnkdata = lnk[_location];\
1033: } while (_entry > _lnkdata);\
1034: /* insertion location is found, add entry into lnk */\
1035: lnk[_location] = _entry;\
1036: lnk[_entry] = _lnkdata;\
1037: lnklvl[_entry] = 0;\
1038: nlnk++;\
1039: _lnkdata = _entry; /* next search starts from here if next_entry > _entry */\
1040: }\
1041: }\
1042: }
1044: /*
1045: Add a SORTED index set into a sorted linked list for ILU
1046: Input Parameters:
1047: nidx - number of input indices
1048: idx - sorted integer array used for storing column indices
1049: level - level of fill, e.g., ICC(level)
1050: idxlvl - level of idx
1051: idx_start - starting index of the list
1052: lnk - linked list(an integer array) that is created
1053: lnklvl - levels of lnk
1054: bt - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
1055: prow - the row number of idx
1056: output Parameters:
1057: nlnk - number of newly added idx
1058: lnk - the sorted(increasing order) linked list containing new and non-redundate entries from idx
1059: lnklvl - levels of lnk
1060: bt - updated PetscBT (bitarray)
1062: Note: the level of factor(i,j) is set as lvl(i,j) = min{ lvl(i,j), lvl(i,prow)+lvl(prow,j)+1)
1063: where idx = non-zero columns of U(prow,prow+1:n-1), prow<i
1064: */
1065: #define PetscILULLAddSorted(nidx,idx,level,idxlvl,idx_start,nlnk,lnk,lnklvl,bt,lnklvl_prow) 0;\
1066: {\
1067: PetscInt _k,_entry,_location,_lnkdata,_incrlev,_lnklvl_prow=lnklvl[prow];\
1068: nlnk = 0;\
1069: _lnkdata = idx_start;\
1070: for (_k=0; _k<nidx; _k++){\
1071: _incrlev = idxlvl[_k] + _lnklvl_prow + 1;\
1072: if (_incrlev > level) continue;\
1073: _entry = idx[_k];\
1074: if (!PetscBTLookupSet(bt,_entry)){ /* new entry */\
1075: /* search for insertion location */\
1076: do {\
1077: _location = _lnkdata;\
1078: _lnkdata = lnk[_location];\
1079: } while (_entry > _lnkdata);\
1080: /* insertion location is found, add entry into lnk */\
1081: lnk[_location] = _entry;\
1082: lnk[_entry] = _lnkdata;\
1083: lnklvl[_entry] = _incrlev;\
1084: nlnk++;\
1085: _lnkdata = _entry; /* next search starts from here if next_entry > _entry */\
1086: } else { /* existing entry: update lnklvl */\
1087: if (lnklvl[_entry] > _incrlev) lnklvl[_entry] = _incrlev;\
1088: }\
1089: }\
1090: }
1092: /*
1093: Add a index set into a sorted linked list
1094: Input Parameters:
1095: nidx - number of input idx
1096: idx - integer array used for storing column indices
1097: level - level of fill, e.g., ICC(level)
1098: idxlvl - level of idx
1099: idx_start - starting index of the list
1100: lnk - linked list(an integer array) that is created
1101: lnklvl - levels of lnk
1102: bt - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
1103: output Parameters:
1104: nlnk - number of newly added idx
1105: lnk - the sorted(increasing order) linked list containing new and non-redundate entries from idx
1106: lnklvl - levels of lnk
1107: bt - updated PetscBT (bitarray)
1108: */
1109: #define PetscIncompleteLLAdd(nidx,idx,level,idxlvl,idx_start,nlnk,lnk,lnklvl,bt) 0;\
1110: {\
1111: PetscInt _k,_entry,_location,_lnkdata,_incrlev;\
1112: nlnk = 0;\
1113: _lnkdata = idx_start;\
1114: for (_k=0; _k<nidx; _k++){\
1115: _incrlev = idxlvl[_k] + 1;\
1116: if (_incrlev > level) continue;\
1117: _entry = idx[_k];\
1118: if (!PetscBTLookupSet(bt,_entry)){ /* new entry */\
1119: /* search for insertion location */\
1120: if (_k && _entry < _lnkdata) _lnkdata = idx_start;\
1121: do {\
1122: _location = _lnkdata;\
1123: _lnkdata = lnk[_location];\
1124: } while (_entry > _lnkdata);\
1125: /* insertion location is found, add entry into lnk */\
1126: lnk[_location] = _entry;\
1127: lnk[_entry] = _lnkdata;\
1128: lnklvl[_entry] = _incrlev;\
1129: nlnk++;\
1130: _lnkdata = _entry; /* next search starts from here if next_entry > _entry */\
1131: } else { /* existing entry: update lnklvl */\
1132: if (lnklvl[_entry] > _incrlev) lnklvl[_entry] = _incrlev;\
1133: }\
1134: }\
1135: }
1137: /*
1138: Add a SORTED index set into a sorted linked list
1139: Input Parameters:
1140: nidx - number of input indices
1141: idx - sorted integer array used for storing column indices
1142: level - level of fill, e.g., ICC(level)
1143: idxlvl - level of idx
1144: idx_start - starting index of the list
1145: lnk - linked list(an integer array) that is created
1146: lnklvl - levels of lnk
1147: bt - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
1148: output Parameters:
1149: nlnk - number of newly added idx
1150: lnk - the sorted(increasing order) linked list containing new and non-redundate entries from idx
1151: lnklvl - levels of lnk
1152: bt - updated PetscBT (bitarray)
1153: */
1154: #define PetscIncompleteLLAddSorted(nidx,idx,level,idxlvl,idx_start,nlnk,lnk,lnklvl,bt) 0;\
1155: {\
1156: PetscInt _k,_entry,_location,_lnkdata,_incrlev;\
1157: nlnk = 0;\
1158: _lnkdata = idx_start;\
1159: for (_k=0; _k<nidx; _k++){\
1160: _incrlev = idxlvl[_k] + 1;\
1161: if (_incrlev > level) continue;\
1162: _entry = idx[_k];\
1163: if (!PetscBTLookupSet(bt,_entry)){ /* new entry */\
1164: /* search for insertion location */\
1165: do {\
1166: _location = _lnkdata;\
1167: _lnkdata = lnk[_location];\
1168: } while (_entry > _lnkdata);\
1169: /* insertion location is found, add entry into lnk */\
1170: lnk[_location] = _entry;\
1171: lnk[_entry] = _lnkdata;\
1172: lnklvl[_entry] = _incrlev;\
1173: nlnk++;\
1174: _lnkdata = _entry; /* next search starts from here */\
1175: } else { /* existing entry: update lnklvl */\
1176: if (lnklvl[_entry] > _incrlev) lnklvl[_entry] = _incrlev;\
1177: }\
1178: }\
1179: }
1181: /*
1182: Add a SORTED index set into a sorted linked list for ICC
1183: Input Parameters:
1184: nidx - number of input indices
1185: idx - sorted integer array used for storing column indices
1186: level - level of fill, e.g., ICC(level)
1187: idxlvl - level of idx
1188: idx_start - starting index of the list
1189: lnk - linked list(an integer array) that is created
1190: lnklvl - levels of lnk
1191: bt - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
1192: idxlvl_prow - idxlvl[prow], where prow is the row number of the idx
1193: output Parameters:
1194: nlnk - number of newly added indices
1195: lnk - the sorted(increasing order) linked list containing new and non-redundate entries from idx
1196: lnklvl - levels of lnk
1197: bt - updated PetscBT (bitarray)
1198: Note: the level of U(i,j) is set as lvl(i,j) = min{ lvl(i,j), lvl(prow,i)+lvl(prow,j)+1)
1199: where idx = non-zero columns of U(prow,prow+1:n-1), prow<i
1200: */
1201: #define PetscICCLLAddSorted(nidx,idx,level,idxlvl,idx_start,nlnk,lnk,lnklvl,bt,idxlvl_prow) 0;\
1202: {\
1203: PetscInt _k,_entry,_location,_lnkdata,_incrlev;\
1204: nlnk = 0;\
1205: _lnkdata = idx_start;\
1206: for (_k=0; _k<nidx; _k++){\
1207: _incrlev = idxlvl[_k] + idxlvl_prow + 1;\
1208: if (_incrlev > level) continue;\
1209: _entry = idx[_k];\
1210: if (!PetscBTLookupSet(bt,_entry)){ /* new entry */\
1211: /* search for insertion location */\
1212: do {\
1213: _location = _lnkdata;\
1214: _lnkdata = lnk[_location];\
1215: } while (_entry > _lnkdata);\
1216: /* insertion location is found, add entry into lnk */\
1217: lnk[_location] = _entry;\
1218: lnk[_entry] = _lnkdata;\
1219: lnklvl[_entry] = _incrlev;\
1220: nlnk++;\
1221: _lnkdata = _entry; /* next search starts from here */\
1222: } else { /* existing entry: update lnklvl */\
1223: if (lnklvl[_entry] > _incrlev) lnklvl[_entry] = _incrlev;\
1224: }\
1225: }\
1226: }
1228: /*
1229: Copy data on the list into an array, then initialize the list
1230: Input Parameters:
1231: idx_start - starting index of the list
1232: lnk_max - max value of lnk indicating the end of the list
1233: nlnk - number of data on the list to be copied
1234: lnk - linked list
1235: lnklvl - level of lnk
1236: bt - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
1237: output Parameters:
1238: indices - array that contains the copied data
1239: lnk - linked list that is cleaned and initialize
1240: lnklvl - level of lnk that is reinitialized
1241: bt - PetscBT (bitarray) with all bits set to false
1242: */
1243: #define PetscIncompleteLLClean(idx_start,lnk_max,nlnk,lnk,lnklvl,indices,indiceslvl,bt) 0;\
1244: {\
1245: PetscInt _j,_idx=idx_start;\
1246: for (_j=0; _j<nlnk; _j++){\
1247: _idx = lnk[_idx];\
1248: *(indices+_j) = _idx;\
1249: *(indiceslvl+_j) = lnklvl[_idx];\
1250: lnklvl[_idx] = -1;\
1251: PetscBTClear(bt,_idx);\
1252: }\
1253: lnk[idx_start] = lnk_max;\
1254: }
1255: /*
1256: Free memories used by the list
1257: */
1258: #define PetscIncompleteLLDestroy(lnk,bt) (PetscFree(lnk) || PetscBTDestroy(&(bt)))
1260: #define MatCheckSameLocalSize(A,ar1,B,ar2) \
1262: 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);
1263:
1264: #define MatCheckSameSize(A,ar1,B,ar2) \
1265: 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);\
1266: MatCheckSameLocalSize(A,ar1,B,ar2);
1267:
1268: #define VecCheckMatCompatible(M,x,ar1,b,ar2) \
1269: 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);\
1270: 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);
1272: /* -------------------------------------------------------------------------------------------------------*/
1273: #include <petscbt.h>
1274: /*
1275: Create and initialize a condensed linked list -
1276: same as PetscLLCreate(), but uses a scalable array 'lnk' with size of max number of entries, not O(N).
1277: Barry suggested this approach (Dec. 6, 2011):
1278: I've thought of an alternative way of representing a linked list that is efficient but doesn't have the O(N) scaling issue
1279: (it may be faster than the O(N) even sequentially due to less crazy memory access).
1281: 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
1282: 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
1283: 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
1284: 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.
1285: 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
1286: to each other so memory access is much better than using the big array.
1288: Example:
1289: nlnk_max=5, lnk_max=36:
1290: Initial list: [0, 0 | 36, 2 | 0, 0 | 0, 0 | 0, 0 | 0, 0 | 0, 0]
1291: here, head_node has index 2 with value lnk[2]=lnk_max=36,
1292: 0-th entry is used to store the number of entries in the list,
1293: The initial lnk represents head -> tail(marked by 36) with number of entries = lnk[0]=0.
1295: Now adding a sorted set {2,4}, the list becomes
1296: [2, 0 | 36, 4 |2, 6 | 4, 2 | 0, 0 | 0, 0 | 0, 0 ]
1297: represents head -> 2 -> 4 -> tail with number of entries = lnk[0]=2.
1299: Then adding a sorted set {0,3,35}, the list
1300: [5, 0 | 36, 8 | 2, 10 | 4, 12 | 0, 4 | 3, 6 | 35, 2 ]
1301: represents head -> 0 -> 2 -> 3 -> 4 -> 35 -> tail with number of entries = lnk[0]=5.
1303: Input Parameters:
1304: nlnk_max - max length of the list
1305: lnk_max - max value of the entries
1306: Output Parameters:
1307: lnk - list created and initialized
1308: bt - PetscBT (bitarray) with all bits set to false. Note: bt has size lnk_max, not nln_max!
1309: */
1310: PETSC_STATIC_INLINE PetscErrorCode PetscLLCondensedCreate(PetscInt nlnk_max,PetscInt lnk_max,PetscInt **lnk,PetscBT *bt)
1311: {
1313: PetscInt *llnk,lsize = 0;
1316: PetscIntMultError(2,nlnk_max+2,&lsize);
1317: PetscMalloc1(lsize,lnk);
1318: PetscBTCreate(lnk_max,bt);
1319: llnk = *lnk;
1320: llnk[0] = 0; /* number of entries on the list */
1321: llnk[2] = lnk_max; /* value in the head node */
1322: llnk[3] = 2; /* next for the head node */
1323: return(0);
1324: }
1326: /*
1327: Add a SORTED ascending index set into a sorted linked list. See PetscLLCondensedCreate() for detailed description.
1328: Input Parameters:
1329: nidx - number of input indices
1330: indices - sorted integer array
1331: lnk - condensed linked list(an integer array) that is created
1332: bt - PetscBT (bitarray), bt[idx]=true marks idx is in lnk
1333: output Parameters:
1334: lnk - the sorted(increasing order) linked list containing previous and newly added non-redundate indices
1335: bt - updated PetscBT (bitarray)
1336: */
1337: PETSC_STATIC_INLINE PetscErrorCode PetscLLCondensedAddSorted(PetscInt nidx,const PetscInt indices[],PetscInt lnk[],PetscBT bt)
1338: {
1339: PetscInt _k,_entry,_location,_next,_lnkdata,_nlnk,_newnode;
1342: _nlnk = lnk[0]; /* num of entries on the input lnk */
1343: _location = 2; /* head */
1344: for (_k=0; _k<nidx; _k++){
1345: _entry = indices[_k];
1346: if (!PetscBTLookupSet(bt,_entry)){ /* new entry */
1347: /* search for insertion location */
1348: do {
1349: _next = _location + 1; /* link from previous node to next node */
1350: _location = lnk[_next]; /* idx of next node */
1351: _lnkdata = lnk[_location];/* value of next node */
1352: } while (_entry > _lnkdata);
1353: /* insertion location is found, add entry into lnk */
1354: _newnode = 2*(_nlnk+2); /* index for this new node */
1355: lnk[_next] = _newnode; /* connect previous node to the new node */
1356: lnk[_newnode] = _entry; /* set value of the new node */
1357: lnk[_newnode+1] = _location; /* connect new node to next node */
1358: _location = _newnode; /* next search starts from the new node */
1359: _nlnk++;
1360: } \
1361: }\
1362: lnk[0] = _nlnk; /* number of entries in the list */
1363: return(0);
1364: }
1366: PETSC_STATIC_INLINE PetscErrorCode PetscLLCondensedClean(PetscInt lnk_max,PetscInt nidx,PetscInt *indices,PetscInt lnk[],PetscBT bt)
1367: {
1369: PetscInt _k,_next,_nlnk;
1372: _next = lnk[3]; /* head node */
1373: _nlnk = lnk[0]; /* num of entries on the list */
1374: for (_k=0; _k<_nlnk; _k++){
1375: indices[_k] = lnk[_next];
1376: _next = lnk[_next + 1];
1377: PetscBTClear(bt,indices[_k]);
1378: }
1379: lnk[0] = 0; /* num of entries on the list */
1380: lnk[2] = lnk_max; /* initialize head node */
1381: lnk[3] = 2; /* head node */
1382: return(0);
1383: }
1385: PETSC_STATIC_INLINE PetscErrorCode PetscLLCondensedView(PetscInt *lnk)
1386: {
1388: PetscInt k;
1391: PetscPrintf(PETSC_COMM_SELF,"LLCondensed of size %D, (val, next)\n",lnk[0]);
1392: for (k=2; k< lnk[0]+2; k++){
1393: PetscPrintf(PETSC_COMM_SELF," %D: (%D, %D)\n",2*k,lnk[2*k],lnk[2*k+1]);
1394: }
1395: return(0);
1396: }
1398: /*
1399: Free memories used by the list
1400: */
1401: PETSC_STATIC_INLINE PetscErrorCode PetscLLCondensedDestroy(PetscInt *lnk,PetscBT bt)
1402: {
1406: PetscFree(lnk);
1407: PetscBTDestroy(&bt);
1408: return(0);
1409: }
1411: /* -------------------------------------------------------------------------------------------------------*/
1412: /*
1413: Same as PetscLLCondensedCreate(), but does not use non-scalable O(lnk_max) bitarray
1414: Input Parameters:
1415: nlnk_max - max length of the list
1416: Output Parameters:
1417: lnk - list created and initialized
1418: */
1419: PETSC_STATIC_INLINE PetscErrorCode PetscLLCondensedCreate_Scalable(PetscInt nlnk_max,PetscInt **lnk)
1420: {
1422: PetscInt *llnk,lsize = 0;
1425: PetscIntMultError(2,nlnk_max+2,&lsize);
1426: PetscMalloc1(lsize,lnk);
1427: llnk = *lnk;
1428: llnk[0] = 0; /* number of entries on the list */
1429: llnk[2] = PETSC_MAX_INT; /* value in the head node */
1430: llnk[3] = 2; /* next for the head node */
1431: return(0);
1432: }
1434: PETSC_STATIC_INLINE PetscErrorCode PetscLLCondensedExpand_Scalable(PetscInt nlnk_max,PetscInt **lnk)
1435: {
1437: PetscInt lsize = 0;
1440: PetscIntMultError(2,nlnk_max+2,&lsize);
1441: PetscRealloc(lsize*sizeof(PetscInt),lnk);
1442: return(0);
1443: }
1445: PETSC_STATIC_INLINE PetscErrorCode PetscLLCondensedAddSorted_Scalable(PetscInt nidx,const PetscInt indices[],PetscInt lnk[])
1446: {
1447: PetscInt _k,_entry,_location,_next,_lnkdata,_nlnk,_newnode;
1448: _nlnk = lnk[0]; /* num of entries on the input lnk */
1449: _location = 2; /* head */ \
1450: for (_k=0; _k<nidx; _k++){
1451: _entry = indices[_k];
1452: /* search for insertion location */
1453: do {
1454: _next = _location + 1; /* link from previous node to next node */
1455: _location = lnk[_next]; /* idx of next node */
1456: _lnkdata = lnk[_location];/* value of next node */
1457: } while (_entry > _lnkdata);
1458: if (_entry < _lnkdata) {
1459: /* insertion location is found, add entry into lnk */
1460: _newnode = 2*(_nlnk+2); /* index for this new node */
1461: lnk[_next] = _newnode; /* connect previous node to the new node */
1462: lnk[_newnode] = _entry; /* set value of the new node */
1463: lnk[_newnode+1] = _location; /* connect new node to next node */
1464: _location = _newnode; /* next search starts from the new node */
1465: _nlnk++;
1466: }
1467: }
1468: lnk[0] = _nlnk; /* number of entries in the list */
1469: return 0;
1470: }
1472: PETSC_STATIC_INLINE PetscErrorCode PetscLLCondensedClean_Scalable(PetscInt nidx,PetscInt *indices,PetscInt *lnk)
1473: {
1474: PetscInt _k,_next,_nlnk;
1475: _next = lnk[3]; /* head node */
1476: _nlnk = lnk[0];
1477: for (_k=0; _k<_nlnk; _k++){
1478: indices[_k] = lnk[_next];
1479: _next = lnk[_next + 1];
1480: }
1481: lnk[0] = 0; /* num of entries on the list */
1482: lnk[3] = 2; /* head node */
1483: return 0;
1484: }
1486: PETSC_STATIC_INLINE PetscErrorCode PetscLLCondensedDestroy_Scalable(PetscInt *lnk)
1487: {
1488: return PetscFree(lnk);
1489: }
1491: /* -------------------------------------------------------------------------------------------------------*/
1492: /*
1493: lnk[0] number of links
1494: lnk[1] number of entries
1495: lnk[3n] value
1496: lnk[3n+1] len
1497: lnk[3n+2] link to next value
1499: The next three are always the first link
1501: lnk[3] PETSC_MIN_INT+1
1502: lnk[4] 1
1503: lnk[5] link to first real entry
1505: The next three are always the last link
1507: lnk[6] PETSC_MAX_INT - 1
1508: lnk[7] 1
1509: lnk[8] next valid link (this is the same as lnk[0] but without the decreases)
1510: */
1512: PETSC_STATIC_INLINE PetscErrorCode PetscLLCondensedCreate_fast(PetscInt nlnk_max,PetscInt **lnk)
1513: {
1515: PetscInt *llnk,lsize = 0;
1518: PetscIntMultError(3,nlnk_max+3,&lsize);
1519: PetscMalloc1(lsize,lnk);
1520: llnk = *lnk;
1521: llnk[0] = 0; /* nlnk: number of entries on the list */
1522: llnk[1] = 0; /* number of integer entries represented in list */
1523: llnk[3] = PETSC_MIN_INT+1; /* value in the first node */
1524: llnk[4] = 1; /* count for the first node */
1525: llnk[5] = 6; /* next for the first node */
1526: llnk[6] = PETSC_MAX_INT-1; /* value in the last node */
1527: llnk[7] = 1; /* count for the last node */
1528: llnk[8] = 0; /* next valid node to be used */
1529: return(0);
1530: }
1532: PETSC_STATIC_INLINE PetscErrorCode PetscLLCondensedAddSorted_fast(PetscInt nidx,const PetscInt indices[],PetscInt lnk[])
1533: {
1534: PetscInt k,entry,prev,next;
1535: prev = 3; /* first value */
1536: next = lnk[prev+2];
1537: for (k=0; k<nidx; k++){
1538: entry = indices[k];
1539: /* search for insertion location */
1540: while (entry >= lnk[next]) {
1541: prev = next;
1542: next = lnk[next+2];
1543: }
1544: /* entry is in range of previous list */
1545: if (entry < lnk[prev]+lnk[prev+1]) continue;
1546: lnk[1]++;
1547: /* entry is right after previous list */
1548: if (entry == lnk[prev]+lnk[prev+1]) {
1549: lnk[prev+1]++;
1550: if (lnk[next] == entry+1) { /* combine two contiguous strings */
1551: lnk[prev+1] += lnk[next+1];
1552: lnk[prev+2] = lnk[next+2];
1553: next = lnk[next+2];
1554: lnk[0]--;
1555: }
1556: continue;
1557: }
1558: /* entry is right before next list */
1559: if (entry == lnk[next]-1) {
1560: lnk[next]--;
1561: lnk[next+1]++;
1562: prev = next;
1563: next = lnk[prev+2];
1564: continue;
1565: }
1566: /* add entry into lnk */
1567: lnk[prev+2] = 3*((lnk[8]++)+3); /* connect previous node to the new node */
1568: prev = lnk[prev+2];
1569: lnk[prev] = entry; /* set value of the new node */
1570: lnk[prev+1] = 1; /* number of values in contiguous string is one to start */
1571: lnk[prev+2] = next; /* connect new node to next node */
1572: lnk[0]++;
1573: }
1574: return 0;
1575: }
1577: PETSC_STATIC_INLINE PetscErrorCode PetscLLCondensedClean_fast(PetscInt nidx,PetscInt *indices,PetscInt *lnk)
1578: {
1579: PetscInt _k,_next,_nlnk,cnt,j;
1580: _next = lnk[5]; /* first node */
1581: _nlnk = lnk[0];
1582: cnt = 0;
1583: for (_k=0; _k<_nlnk; _k++){
1584: for (j=0; j<lnk[_next+1]; j++) {
1585: indices[cnt++] = lnk[_next] + j;
1586: }
1587: _next = lnk[_next + 2];
1588: }
1589: lnk[0] = 0; /* nlnk: number of links */
1590: lnk[1] = 0; /* number of integer entries represented in list */
1591: lnk[3] = PETSC_MIN_INT+1; /* value in the first node */
1592: lnk[4] = 1; /* count for the first node */
1593: lnk[5] = 6; /* next for the first node */
1594: lnk[6] = PETSC_MAX_INT-1; /* value in the last node */
1595: lnk[7] = 1; /* count for the last node */
1596: lnk[8] = 0; /* next valid location to make link */
1597: return 0;
1598: }
1600: PETSC_STATIC_INLINE PetscErrorCode PetscLLCondensedView_fast(PetscInt *lnk)
1601: {
1602: PetscInt k,next,nlnk;
1603: next = lnk[5]; /* first node */
1604: nlnk = lnk[0];
1605: for (k=0; k<nlnk; k++){
1606: #if 0 /* Debugging code */
1607: printf("%d value %d len %d next %d\n",next,lnk[next],lnk[next+1],lnk[next+2]);
1608: #endif
1609: next = lnk[next + 2];
1610: }
1611: return 0;
1612: }
1614: PETSC_STATIC_INLINE PetscErrorCode PetscLLCondensedDestroy_fast(PetscInt *lnk)
1615: {
1616: return PetscFree(lnk);
1617: }
1619: /* this is extern because it is used in MatFDColoringUseDM() which is in the DM library */
1620: PETSC_EXTERN PetscErrorCode MatFDColoringApply_AIJ(Mat,MatFDColoring,Vec,void*);
1622: PETSC_EXTERN PetscLogEvent MAT_Mult;
1623: PETSC_EXTERN PetscLogEvent MAT_MultMatrixFree;
1624: PETSC_EXTERN PetscLogEvent MAT_Mults;
1625: PETSC_EXTERN PetscLogEvent MAT_MultConstrained;
1626: PETSC_EXTERN PetscLogEvent MAT_MultAdd;
1627: PETSC_EXTERN PetscLogEvent MAT_MultTranspose;
1628: PETSC_EXTERN PetscLogEvent MAT_MultTransposeConstrained;
1629: PETSC_EXTERN PetscLogEvent MAT_MultTransposeAdd;
1630: PETSC_EXTERN PetscLogEvent MAT_Solve;
1631: PETSC_EXTERN PetscLogEvent MAT_Solves;
1632: PETSC_EXTERN PetscLogEvent MAT_SolveAdd;
1633: PETSC_EXTERN PetscLogEvent MAT_SolveTranspose;
1634: PETSC_EXTERN PetscLogEvent MAT_SolveTransposeAdd;
1635: PETSC_EXTERN PetscLogEvent MAT_SOR;
1636: PETSC_EXTERN PetscLogEvent MAT_ForwardSolve;
1637: PETSC_EXTERN PetscLogEvent MAT_BackwardSolve;
1638: PETSC_EXTERN PetscLogEvent MAT_LUFactor;
1639: PETSC_EXTERN PetscLogEvent MAT_LUFactorSymbolic;
1640: PETSC_EXTERN PetscLogEvent MAT_LUFactorNumeric;
1641: PETSC_EXTERN PetscLogEvent MAT_CholeskyFactor;
1642: PETSC_EXTERN PetscLogEvent MAT_CholeskyFactorSymbolic;
1643: PETSC_EXTERN PetscLogEvent MAT_CholeskyFactorNumeric;
1644: PETSC_EXTERN PetscLogEvent MAT_ILUFactor;
1645: PETSC_EXTERN PetscLogEvent MAT_ILUFactorSymbolic;
1646: PETSC_EXTERN PetscLogEvent MAT_ICCFactorSymbolic;
1647: PETSC_EXTERN PetscLogEvent MAT_Copy;
1648: PETSC_EXTERN PetscLogEvent MAT_Convert;
1649: PETSC_EXTERN PetscLogEvent MAT_Scale;
1650: PETSC_EXTERN PetscLogEvent MAT_AssemblyBegin;
1651: PETSC_EXTERN PetscLogEvent MAT_AssemblyEnd;
1652: PETSC_EXTERN PetscLogEvent MAT_SetValues;
1653: PETSC_EXTERN PetscLogEvent MAT_GetValues;
1654: PETSC_EXTERN PetscLogEvent MAT_GetRow;
1655: PETSC_EXTERN PetscLogEvent MAT_GetRowIJ;
1656: PETSC_EXTERN PetscLogEvent MAT_CreateSubMats;
1657: PETSC_EXTERN PetscLogEvent MAT_GetColoring;
1658: PETSC_EXTERN PetscLogEvent MAT_GetOrdering;
1659: PETSC_EXTERN PetscLogEvent MAT_RedundantMat;
1660: PETSC_EXTERN PetscLogEvent MAT_IncreaseOverlap;
1661: PETSC_EXTERN PetscLogEvent MAT_Partitioning;
1662: PETSC_EXTERN PetscLogEvent MAT_PartitioningND;
1663: PETSC_EXTERN PetscLogEvent MAT_Coarsen;
1664: PETSC_EXTERN PetscLogEvent MAT_ZeroEntries;
1665: PETSC_EXTERN PetscLogEvent MAT_Load;
1666: PETSC_EXTERN PetscLogEvent MAT_View;
1667: PETSC_EXTERN PetscLogEvent MAT_AXPY;
1668: PETSC_EXTERN PetscLogEvent MAT_FDColoringCreate;
1669: PETSC_EXTERN PetscLogEvent MAT_TransposeColoringCreate;
1670: PETSC_EXTERN PetscLogEvent MAT_FDColoringSetUp;
1671: PETSC_EXTERN PetscLogEvent MAT_FDColoringApply;
1672: PETSC_EXTERN PetscLogEvent MAT_Transpose;
1673: PETSC_EXTERN PetscLogEvent MAT_FDColoringFunction;
1674: PETSC_EXTERN PetscLogEvent MAT_CreateSubMat;
1675: PETSC_EXTERN PetscLogEvent MAT_MatMult;
1676: PETSC_EXTERN PetscLogEvent MAT_MatSolve;
1677: PETSC_EXTERN PetscLogEvent MAT_MatTrSolve;
1678: PETSC_EXTERN PetscLogEvent MAT_MatMultSymbolic;
1679: PETSC_EXTERN PetscLogEvent MAT_MatMultNumeric;
1680: PETSC_EXTERN PetscLogEvent MAT_Getlocalmatcondensed;
1681: PETSC_EXTERN PetscLogEvent MAT_GetBrowsOfAcols;
1682: PETSC_EXTERN PetscLogEvent MAT_GetBrowsOfAocols;
1683: PETSC_EXTERN PetscLogEvent MAT_PtAP;
1684: PETSC_EXTERN PetscLogEvent MAT_PtAPSymbolic;
1685: PETSC_EXTERN PetscLogEvent MAT_PtAPNumeric;
1686: PETSC_EXTERN PetscLogEvent MAT_Seqstompinum;
1687: PETSC_EXTERN PetscLogEvent MAT_Seqstompisym;
1688: PETSC_EXTERN PetscLogEvent MAT_Seqstompi;
1689: PETSC_EXTERN PetscLogEvent MAT_Getlocalmat;
1690: PETSC_EXTERN PetscLogEvent MAT_RARt;
1691: PETSC_EXTERN PetscLogEvent MAT_RARtSymbolic;
1692: PETSC_EXTERN PetscLogEvent MAT_RARtNumeric;
1693: PETSC_EXTERN PetscLogEvent MAT_MatTransposeMult;
1694: PETSC_EXTERN PetscLogEvent MAT_MatTransposeMultSymbolic;
1695: PETSC_EXTERN PetscLogEvent MAT_MatTransposeMultNumeric;
1696: PETSC_EXTERN PetscLogEvent MAT_TransposeMatMult;
1697: PETSC_EXTERN PetscLogEvent MAT_TransposeMatMultSymbolic;
1698: PETSC_EXTERN PetscLogEvent MAT_TransposeMatMultNumeric;
1699: PETSC_EXTERN PetscLogEvent MAT_MatMatMult;
1700: PETSC_EXTERN PetscLogEvent MAT_MatMatMultSymbolic;
1701: PETSC_EXTERN PetscLogEvent MAT_MatMatMultNumeric;
1702: PETSC_EXTERN PetscLogEvent MAT_Applypapt;
1703: PETSC_EXTERN PetscLogEvent MAT_Applypapt_symbolic;
1704: PETSC_EXTERN PetscLogEvent MAT_Applypapt_numeric;
1705: PETSC_EXTERN PetscLogEvent MAT_Getsymtranspose;
1706: PETSC_EXTERN PetscLogEvent MAT_Transpose_SeqAIJ;
1707: PETSC_EXTERN PetscLogEvent MAT_Getsymtransreduced;
1708: PETSC_EXTERN PetscLogEvent MAT_GetSequentialNonzeroStructure;
1709: PETSC_EXTERN PetscLogEvent MATMFFD_Mult;
1710: PETSC_EXTERN PetscLogEvent MAT_GetMultiProcBlock;
1711: PETSC_EXTERN PetscLogEvent MAT_CUSPARSECopyToGPU;
1712: PETSC_EXTERN PetscLogEvent MAT_SetValuesBatch;
1713: PETSC_EXTERN PetscLogEvent MAT_ViennaCLCopyToGPU;
1714: PETSC_EXTERN PetscLogEvent MAT_Merge;
1715: PETSC_EXTERN PetscLogEvent MAT_Residual;
1716: PETSC_EXTERN PetscLogEvent MAT_SetRandom;
1717: PETSC_EXTERN PetscLogEvent MATCOLORING_Apply;
1718: PETSC_EXTERN PetscLogEvent MATCOLORING_Comm;
1719: PETSC_EXTERN PetscLogEvent MATCOLORING_Local;
1720: PETSC_EXTERN PetscLogEvent MATCOLORING_ISCreate;
1721: PETSC_EXTERN PetscLogEvent MATCOLORING_SetUp;
1722: PETSC_EXTERN PetscLogEvent MATCOLORING_Weights;
1724: #endif