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