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