Actual source code: parms.c
petsc-3.8.4 2018-03-24
1: #define PETSCKSP_DLL
3: /*
4: Provides an interface to pARMS.
5: Requires pARMS 3.2 or later.
6: */
8: #include <petsc/private/pcimpl.h>
10: #if defined(PETSC_USE_COMPLEX)
11: #define DBL_CMPLX
12: #else
13: #define DBL
14: #endif
15: #define USE_MPI
16: #define REAL double
17: #define HAS_BLAS
18: #define FORTRAN_UNDERSCORE
19: #include "parms_sys.h"
20: #undef FLOAT
21: #define FLOAT PetscScalar
22: #include <parms.h>
24: /*
25: Private context (data structure) for the preconditioner.
26: */
27: typedef struct {
28: parms_Map map;
29: parms_Mat A;
30: parms_PC pc;
31: PCPARMSGlobalType global;
32: PCPARMSLocalType local;
33: PetscInt levels, blocksize, maxdim, maxits, lfil[7];
34: PetscBool nonsymperm, meth[8];
35: PetscReal solvetol, indtol, droptol[7];
36: PetscScalar *lvec0, *lvec1;
37: } PC_PARMS;
40: static PetscErrorCode PCSetUp_PARMS(PC pc)
41: {
42: Mat pmat;
43: PC_PARMS *parms = (PC_PARMS*)pc->data;
44: const PetscInt *mapptr0;
45: PetscInt n, lsize, low, high, i, pos, ncols, length;
46: int *maptmp, *mapptr, *ia, *ja, *ja1, *im;
47: PetscScalar *aa, *aa1;
48: const PetscInt *cols;
49: PetscInt meth[8];
50: const PetscScalar *values;
51: PetscErrorCode ierr;
52: MatInfo matinfo;
53: PetscMPIInt rank, npro;
56: /* Get preconditioner matrix from PETSc and setup pARMS structs */
57: PCGetOperators(pc,NULL,&pmat);
58: MPI_Comm_size(PetscObjectComm((PetscObject)pmat),&npro);
59: MPI_Comm_rank(PetscObjectComm((PetscObject)pmat),&rank);
61: MatGetSize(pmat,&n,NULL);
62: PetscMalloc1(npro+1,&mapptr);
63: PetscMalloc1(n,&maptmp);
64: MatGetOwnershipRanges(pmat,&mapptr0);
65: low = mapptr0[rank];
66: high = mapptr0[rank+1];
67: lsize = high - low;
69: for (i=0; i<npro+1; i++) mapptr[i] = mapptr0[i]+1;
70: for (i = 0; i<n; i++) maptmp[i] = i+1;
72: /* if created, destroy the previous map */
73: if (parms->map) {
74: parms_MapFree(&parms->map);
75: parms->map = NULL;
76: }
78: /* create pARMS map object */
79: parms_MapCreateFromPtr(&parms->map,(int)n,maptmp,mapptr,PetscObjectComm((PetscObject)pmat),1,NONINTERLACED);
81: /* if created, destroy the previous pARMS matrix */
82: if (parms->A) {
83: parms_MatFree(&parms->A);
84: parms->A = NULL;
85: }
87: /* create pARMS mat object */
88: parms_MatCreate(&parms->A,parms->map);
90: /* setup and copy csr data structure for pARMS */
91: PetscMalloc1(lsize+1,&ia);
92: ia[0] = 1;
93: MatGetInfo(pmat,MAT_LOCAL,&matinfo);
94: length = matinfo.nz_used;
95: PetscMalloc1(length,&ja);
96: PetscMalloc1(length,&aa);
98: for (i = low; i<high; i++) {
99: pos = ia[i-low]-1;
100: MatGetRow(pmat,i,&ncols,&cols,&values);
101: ia[i-low+1] = ia[i-low] + ncols;
103: if (ia[i-low+1] >= length) {
104: length += ncols;
105: PetscMalloc1(length,&ja1);
106: PetscMemcpy(ja1,ja,(ia[i-low]-1)*sizeof(int));
107: PetscFree(ja);
108: ja = ja1;
109: PetscMalloc1(length,&aa1);
110: PetscMemcpy(aa1,aa,(ia[i-low]-1)*sizeof(PetscScalar));
111: PetscFree(aa);
112: aa = aa1;
113: }
114: PetscMemcpy(&ja[pos],cols,ncols*sizeof(int));
115: PetscMemcpy(&aa[pos],values,ncols*sizeof(PetscScalar));
116: MatRestoreRow(pmat,i,&ncols,&cols,&values);
117: }
119: /* csr info is for local matrix so initialize im[] locally */
120: PetscMalloc1(lsize,&im);
121: PetscMemcpy(im,&maptmp[mapptr[rank]-1],lsize*sizeof(int));
123: /* 1-based indexing */
124: for (i=0; i<ia[lsize]-1; i++) ja[i] = ja[i]+1;
126: /* Now copy csr matrix to parms_mat object */
127: parms_MatSetValues(parms->A,(int)lsize,im,ia,ja,aa,INSERT);
129: /* free memory */
130: PetscFree(maptmp);
131: PetscFree(mapptr);
132: PetscFree(aa);
133: PetscFree(ja);
134: PetscFree(ia);
135: PetscFree(im);
137: /* setup parms matrix */
138: parms_MatSetup(parms->A);
140: /* if created, destroy the previous pARMS pc */
141: if (parms->pc) {
142: parms_PCFree(&parms->pc);
143: parms->pc = NULL;
144: }
146: /* Now create pARMS preconditioner object based on A */
147: parms_PCCreate(&parms->pc,parms->A);
149: /* Transfer options from PC to pARMS */
150: switch (parms->global) {
151: case 0: parms_PCSetType(parms->pc, PCRAS); break;
152: case 1: parms_PCSetType(parms->pc, PCSCHUR); break;
153: case 2: parms_PCSetType(parms->pc, PCBJ); break;
154: }
155: switch (parms->local) {
156: case 0: parms_PCSetILUType(parms->pc, PCILU0); break;
157: case 1: parms_PCSetILUType(parms->pc, PCILUK); break;
158: case 2: parms_PCSetILUType(parms->pc, PCILUT); break;
159: case 3: parms_PCSetILUType(parms->pc, PCARMS); break;
160: }
161: parms_PCSetInnerEps(parms->pc, parms->solvetol);
162: parms_PCSetNlevels(parms->pc, parms->levels);
163: parms_PCSetPermType(parms->pc, parms->nonsymperm ? 1 : 0);
164: parms_PCSetBsize(parms->pc, parms->blocksize);
165: parms_PCSetTolInd(parms->pc, parms->indtol);
166: parms_PCSetInnerKSize(parms->pc, parms->maxdim);
167: parms_PCSetInnerMaxits(parms->pc, parms->maxits);
168: for (i=0; i<8; i++) meth[i] = parms->meth[i] ? 1 : 0;
169: parms_PCSetPermScalOptions(parms->pc, &meth[0], 1);
170: parms_PCSetPermScalOptions(parms->pc, &meth[4], 0);
171: parms_PCSetFill(parms->pc, parms->lfil);
172: parms_PCSetTol(parms->pc, parms->droptol);
174: parms_PCSetup(parms->pc);
176: /* Allocate two auxiliary vector of length lsize */
177: if (parms->lvec0) { PetscFree(parms->lvec0); }
178: PetscMalloc1(lsize, &parms->lvec0);
179: if (parms->lvec1) { PetscFree(parms->lvec1); }
180: PetscMalloc1(lsize, &parms->lvec1);
181: return(0);
182: }
184: static PetscErrorCode PCView_PARMS(PC pc,PetscViewer viewer)
185: {
187: PetscBool iascii;
188: PC_PARMS *parms = (PC_PARMS*)pc->data;
189: char *str;
190: double fill_fact;
193: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
194: if (iascii) {
195: parms_PCGetName(parms->pc,&str);
196: PetscViewerASCIIPrintf(viewer," global preconditioner: %s\n",str);
197: parms_PCILUGetName(parms->pc,&str);
198: PetscViewerASCIIPrintf(viewer," local preconditioner: %s\n",str);
199: parms_PCGetRatio(parms->pc,&fill_fact);
200: PetscViewerASCIIPrintf(viewer," non-zero elements/original non-zero entries: %-4.2f\n",fill_fact);
201: PetscViewerASCIIPrintf(viewer," Tolerance for local solve: %g\n",parms->solvetol);
202: PetscViewerASCIIPrintf(viewer," Number of levels: %d\n",parms->levels);
203: if (parms->nonsymperm) {
204: PetscViewerASCIIPrintf(viewer," Using nonsymmetric permutation\n");
205: }
206: PetscViewerASCIIPrintf(viewer," Block size: %d\n",parms->blocksize);
207: PetscViewerASCIIPrintf(viewer," Tolerance for independent sets: %g\n",parms->indtol);
208: PetscViewerASCIIPrintf(viewer," Inner Krylov dimension: %d\n",parms->maxdim);
209: PetscViewerASCIIPrintf(viewer," Maximum number of inner iterations: %d\n",parms->maxits);
210: if (parms->meth[0]) {
211: PetscViewerASCIIPrintf(viewer," Using nonsymmetric permutation for interlevel blocks\n");
212: }
213: if (parms->meth[1]) {
214: PetscViewerASCIIPrintf(viewer," Using column permutation for interlevel blocks\n");
215: }
216: if (parms->meth[2]) {
217: PetscViewerASCIIPrintf(viewer," Using row scaling for interlevel blocks\n");
218: }
219: if (parms->meth[3]) {
220: PetscViewerASCIIPrintf(viewer," Using column scaling for interlevel blocks\n");
221: }
222: if (parms->meth[4]) {
223: PetscViewerASCIIPrintf(viewer," Using nonsymmetric permutation for last level blocks\n");
224: }
225: if (parms->meth[5]) {
226: PetscViewerASCIIPrintf(viewer," Using column permutation for last level blocks\n");
227: }
228: if (parms->meth[6]) {
229: PetscViewerASCIIPrintf(viewer," Using row scaling for last level blocks\n");
230: }
231: if (parms->meth[7]) {
232: PetscViewerASCIIPrintf(viewer," Using column scaling for last level blocks\n");
233: }
234: PetscViewerASCIIPrintf(viewer," amount of fill-in for ilut, iluk and arms: %d\n",parms->lfil[0]);
235: PetscViewerASCIIPrintf(viewer," amount of fill-in for schur: %d\n",parms->lfil[4]);
236: PetscViewerASCIIPrintf(viewer," amount of fill-in for ILUT L and U: %d\n",parms->lfil[5]);
237: PetscViewerASCIIPrintf(viewer," drop tolerance for L, U, L^{-1}F and EU^{-1}: %g\n",parms->droptol[0]);
238: PetscViewerASCIIPrintf(viewer," drop tolerance for schur complement at each level: %g\n",parms->droptol[4]);
239: PetscViewerASCIIPrintf(viewer," drop tolerance for ILUT in last level schur complement: %g\n",parms->droptol[5]);
240: }
241: return(0);
242: }
244: static PetscErrorCode PCDestroy_PARMS(PC pc)
245: {
246: PC_PARMS *parms = (PC_PARMS*)pc->data;
250: if (parms->map) parms_MapFree(&parms->map);
251: if (parms->A) parms_MatFree(&parms->A);
252: if (parms->pc) parms_PCFree(&parms->pc);
253: if (parms->lvec0) {
254: PetscFree(parms->lvec0);
255: }
256: if (parms->lvec1) {
257: PetscFree(parms->lvec1);
258: }
259: PetscFree(pc->data);
261: PetscObjectChangeTypeName((PetscObject)pc,0);
262: PetscObjectComposeFunction((PetscObject)pc,"PCPARMSSetGlobal_C",NULL);
263: PetscObjectComposeFunction((PetscObject)pc,"PCPARMSSetLocal_C",NULL);
264: PetscObjectComposeFunction((PetscObject)pc,"PCPARMSSetSolveTolerances_C",NULL);
265: PetscObjectComposeFunction((PetscObject)pc,"PCPARMSSetSolveRestart_C",NULL);
266: PetscObjectComposeFunction((PetscObject)pc,"PCPARMSSetNonsymPerm_C",NULL);
267: PetscObjectComposeFunction((PetscObject)pc,"PCPARMSSetFill_C",NULL);
268: return(0);
269: }
271: static PetscErrorCode PCSetFromOptions_PARMS(PetscOptionItems *PetscOptionsObject,PC pc)
272: {
273: PC_PARMS *parms = (PC_PARMS*)pc->data;
274: PetscBool flag;
275: PCPARMSGlobalType global;
276: PCPARMSLocalType local;
277: PetscErrorCode ierr;
280: PetscOptionsHead(PetscOptionsObject,"PARMS Options");
281: PetscOptionsEnum("-pc_parms_global","Global preconditioner","PCPARMSSetGlobal",PCPARMSGlobalTypes,(PetscEnum)parms->global,(PetscEnum*)&global,&flag);
282: if (flag) {PCPARMSSetGlobal(pc,global);}
283: PetscOptionsEnum("-pc_parms_local","Local preconditioner","PCPARMSSetLocal",PCPARMSLocalTypes,(PetscEnum)parms->local,(PetscEnum*)&local,&flag);
284: if (flag) {PCPARMSSetLocal(pc,local);}
285: PetscOptionsReal("-pc_parms_solve_tol","Tolerance for local solve","PCPARMSSetSolveTolerances",parms->solvetol,&parms->solvetol,NULL);
286: PetscOptionsInt("-pc_parms_levels","Number of levels","None",parms->levels,&parms->levels,NULL);
287: PetscOptionsBool("-pc_parms_nonsymmetric_perm","Use nonsymmetric permutation","PCPARMSSetNonsymPerm",parms->nonsymperm,&parms->nonsymperm,NULL);
288: PetscOptionsInt("-pc_parms_blocksize","Block size","None",parms->blocksize,&parms->blocksize,NULL);
289: PetscOptionsReal("-pc_parms_ind_tol","Tolerance for independent sets","None",parms->indtol,&parms->indtol,NULL);
290: PetscOptionsInt("-pc_parms_max_dim","Inner Krylov dimension","PCPARMSSetSolveRestart",parms->maxdim,&parms->maxdim,NULL);
291: PetscOptionsInt("-pc_parms_max_it","Maximum number of inner iterations","PCPARMSSetSolveTolerances",parms->maxits,&parms->maxits,NULL);
292: PetscOptionsBool("-pc_parms_inter_nonsymmetric_perm","nonsymmetric permutation for interlevel blocks","None",parms->meth[0],&parms->meth[0],NULL);
293: PetscOptionsBool("-pc_parms_inter_column_perm","column permutation for interlevel blocks","None",parms->meth[1],&parms->meth[1],NULL);
294: PetscOptionsBool("-pc_parms_inter_row_scaling","row scaling for interlevel blocks","None",parms->meth[2],&parms->meth[2],NULL);
295: PetscOptionsBool("-pc_parms_inter_column_scaling","column scaling for interlevel blocks","None",parms->meth[3],&parms->meth[3],NULL);
296: PetscOptionsBool("-pc_parms_last_nonsymmetric_perm","nonsymmetric permutation for last level blocks","None",parms->meth[4],&parms->meth[4],NULL);
297: PetscOptionsBool("-pc_parms_last_column_perm","column permutation for last level blocks","None",parms->meth[5],&parms->meth[5],NULL);
298: PetscOptionsBool("-pc_parms_last_row_scaling","row scaling for last level blocks","None",parms->meth[6],&parms->meth[6],NULL);
299: PetscOptionsBool("-pc_parms_last_column_scaling","column scaling for last level blocks","None",parms->meth[7],&parms->meth[7],NULL);
300: PetscOptionsInt("-pc_parms_lfil_ilu_arms","amount of fill-in for ilut, iluk and arms","PCPARMSSetFill",parms->lfil[0],&parms->lfil[0],&flag);
301: if (flag) parms->lfil[1] = parms->lfil[2] = parms->lfil[3] = parms->lfil[0];
302: PetscOptionsInt("-pc_parms_lfil_schur","amount of fill-in for schur","PCPARMSSetFill",parms->lfil[4],&parms->lfil[4],NULL);
303: PetscOptionsInt("-pc_parms_lfil_ilut_L_U","amount of fill-in for ILUT L and U","PCPARMSSetFill",parms->lfil[5],&parms->lfil[5],&flag);
304: if (flag) parms->lfil[6] = parms->lfil[5];
305: PetscOptionsReal("-pc_parms_droptol_factors","drop tolerance for L, U, L^{-1}F and EU^{-1}","None",parms->droptol[0],&parms->droptol[0],NULL);
306: PetscOptionsReal("-pc_parms_droptol_schur_compl","drop tolerance for schur complement at each level","None",parms->droptol[4],&parms->droptol[4],&flag);
307: if (flag) parms->droptol[1] = parms->droptol[2] = parms->droptol[3] = parms->droptol[0];
308: PetscOptionsReal("-pc_parms_droptol_last_schur","drop tolerance for ILUT in last level schur complement","None",parms->droptol[5],&parms->droptol[5],&flag);
309: if (flag) parms->droptol[6] = parms->droptol[5];
310: PetscOptionsTail();
311: return(0);
312: }
314: static PetscErrorCode PCApply_PARMS(PC pc,Vec b,Vec x)
315: {
316: PetscErrorCode ierr;
317: PC_PARMS *parms = (PC_PARMS*)pc->data;
318: const PetscScalar *b1;
319: PetscScalar *x1;
322: VecGetArrayRead(b,&b1);
323: VecGetArray(x,&x1);
324: parms_VecPermAux((PetscScalar*)b1,parms->lvec0,parms->map);
325: parms_PCApply(parms->pc,parms->lvec0,parms->lvec1);
326: parms_VecInvPermAux(parms->lvec1,x1,parms->map);
327: VecRestoreArrayRead(b,&b1);
328: VecRestoreArray(x,&x1);
329: return(0);
330: }
332: static PetscErrorCode PCPARMSSetGlobal_PARMS(PC pc,PCPARMSGlobalType type)
333: {
334: PC_PARMS *parms = (PC_PARMS*)pc->data;
337: if (type != parms->global) {
338: parms->global = type;
339: pc->setupcalled = 0;
340: }
341: return(0);
342: }
344: /*@
345: PCPARMSSetGlobal - Sets the global preconditioner to be used in PARMS.
347: Collective on PC
349: Input Parameters:
350: + pc - the preconditioner context
351: - type - the global preconditioner type, one of
352: .vb
353: PC_PARMS_GLOBAL_RAS - Restricted additive Schwarz
354: PC_PARMS_GLOBAL_SCHUR - Schur complement
355: PC_PARMS_GLOBAL_BJ - Block Jacobi
356: .ve
358: Options Database Keys:
359: -pc_parms_global [ras,schur,bj] - Sets global preconditioner
361: Level: intermediate
363: Notes:
364: See the pARMS function parms_PCSetType for more information.
366: .seealso: PCPARMS, PCPARMSSetLocal()
367: @*/
368: PetscErrorCode PCPARMSSetGlobal(PC pc,PCPARMSGlobalType type)
369: {
375: PetscTryMethod(pc,"PCPARMSSetGlobal_C",(PC,PCPARMSGlobalType),(pc,type));
376: return(0);
377: }
379: static PetscErrorCode PCPARMSSetLocal_PARMS(PC pc,PCPARMSLocalType type)
380: {
381: PC_PARMS *parms = (PC_PARMS*)pc->data;
384: if (type != parms->local) {
385: parms->local = type;
386: pc->setupcalled = 0;
387: }
388: return(0);
389: }
391: /*@
392: PCPARMSSetLocal - Sets the local preconditioner to be used in PARMS.
394: Collective on PC
396: Input Parameters:
397: + pc - the preconditioner context
398: - type - the local preconditioner type, one of
399: .vb
400: PC_PARMS_LOCAL_ILU0 - ILU0 preconditioner
401: PC_PARMS_LOCAL_ILUK - ILU(k) preconditioner
402: PC_PARMS_LOCAL_ILUT - ILUT preconditioner
403: PC_PARMS_LOCAL_ARMS - ARMS preconditioner
404: .ve
406: Options Database Keys:
407: -pc_parms_local [ilu0,iluk,ilut,arms] - Sets local preconditioner
409: Level: intermediate
411: Notes:
412: For the ARMS preconditioner, one can use either the symmetric ARMS or the non-symmetric
413: variant (ARMS-ddPQ) by setting the permutation type with PCPARMSSetNonsymPerm().
415: See the pARMS function parms_PCILUSetType for more information.
417: .seealso: PCPARMS, PCPARMSSetGlobal(), PCPARMSSetNonsymPerm()
419: @*/
420: PetscErrorCode PCPARMSSetLocal(PC pc,PCPARMSLocalType type)
421: {
427: PetscTryMethod(pc,"PCPARMSSetLocal_C",(PC,PCPARMSLocalType),(pc,type));
428: return(0);
429: }
431: static PetscErrorCode PCPARMSSetSolveTolerances_PARMS(PC pc,PetscReal tol,PetscInt maxits)
432: {
433: PC_PARMS *parms = (PC_PARMS*)pc->data;
436: if (tol != parms->solvetol) {
437: parms->solvetol = tol;
438: pc->setupcalled = 0;
439: }
440: if (maxits != parms->maxits) {
441: parms->maxits = maxits;
442: pc->setupcalled = 0;
443: }
444: return(0);
445: }
447: /*@
448: PCPARMSSetSolveTolerances - Sets the convergence tolerance and the maximum iterations for the
449: inner GMRES solver, when the Schur global preconditioner is used.
451: Collective on PC
453: Input Parameters:
454: + pc - the preconditioner context
455: . tol - the convergence tolerance
456: - maxits - the maximum number of iterations to use
458: Options Database Keys:
459: + -pc_parms_solve_tol - set the tolerance for local solve
460: - -pc_parms_max_it - set the maximum number of inner iterations
462: Level: intermediate
464: Notes:
465: See the pARMS functions parms_PCSetInnerEps and parms_PCSetInnerMaxits for more information.
467: .seealso: PCPARMS, PCPARMSSetSolveRestart()
468: @*/
469: PetscErrorCode PCPARMSSetSolveTolerances(PC pc,PetscReal tol,PetscInt maxits)
470: {
475: PetscTryMethod(pc,"PCPARMSSetSolveTolerances_C",(PC,PetscReal,PetscInt),(pc,tol,maxits));
476: return(0);
477: }
479: static PetscErrorCode PCPARMSSetSolveRestart_PARMS(PC pc,PetscInt restart)
480: {
481: PC_PARMS *parms = (PC_PARMS*)pc->data;
484: if (restart != parms->maxdim) {
485: parms->maxdim = restart;
486: pc->setupcalled = 0;
487: }
488: return(0);
489: }
491: /*@
492: PCPARMSSetSolveRestart - Sets the number of iterations at which the
493: inner GMRES solver restarts.
495: Collective on PC
497: Input Parameters:
498: + pc - the preconditioner context
499: - restart - maximum dimension of the Krylov subspace
501: Options Database Keys:
502: . -pc_parms_max_dim - sets the inner Krylov dimension
504: Level: intermediate
506: Notes:
507: See the pARMS function parms_PCSetInnerKSize for more information.
509: .seealso: PCPARMS, PCPARMSSetSolveTolerances()
510: @*/
511: PetscErrorCode PCPARMSSetSolveRestart(PC pc,PetscInt restart)
512: {
517: PetscTryMethod(pc,"PCPARMSSetSolveRestart_C",(PC,PetscInt),(pc,restart));
518: return(0);
519: }
521: static PetscErrorCode PCPARMSSetNonsymPerm_PARMS(PC pc,PetscBool nonsym)
522: {
523: PC_PARMS *parms = (PC_PARMS*)pc->data;
526: if ((nonsym && !parms->nonsymperm) || (!nonsym && parms->nonsymperm)) {
527: parms->nonsymperm = nonsym;
528: pc->setupcalled = 0;
529: }
530: return(0);
531: }
533: /*@
534: PCPARMSSetNonsymPerm - Sets the type of permutation for the ARMS preconditioner: the standard
535: symmetric ARMS or the non-symmetric ARMS (ARMS-ddPQ).
537: Collective on PC
539: Input Parameters:
540: + pc - the preconditioner context
541: - nonsym - PETSC_TRUE indicates the non-symmetric ARMS is used;
542: PETSC_FALSE indicates the symmetric ARMS is used
544: Options Database Keys:
545: . -pc_parms_nonsymmetric_perm - sets the use of nonsymmetric permutation
547: Level: intermediate
549: Notes:
550: See the pARMS function parms_PCSetPermType for more information.
552: .seealso: PCPARMS
553: @*/
554: PetscErrorCode PCPARMSSetNonsymPerm(PC pc,PetscBool nonsym)
555: {
560: PetscTryMethod(pc,"PCPARMSSetNonsymPerm_C",(PC,PetscBool),(pc,nonsym));
561: return(0);
562: }
564: static PetscErrorCode PCPARMSSetFill_PARMS(PC pc,PetscInt lfil0,PetscInt lfil1,PetscInt lfil2)
565: {
566: PC_PARMS *parms = (PC_PARMS*)pc->data;
569: if (lfil0 != parms->lfil[0] || lfil0 != parms->lfil[1] || lfil0 != parms->lfil[2] || lfil0 != parms->lfil[3]) {
570: parms->lfil[1] = parms->lfil[2] = parms->lfil[3] = parms->lfil[0] = lfil0;
571: pc->setupcalled = 0;
572: }
573: if (lfil1 != parms->lfil[4]) {
574: parms->lfil[4] = lfil1;
575: pc->setupcalled = 0;
576: }
577: if (lfil2 != parms->lfil[5] || lfil2 != parms->lfil[6]) {
578: parms->lfil[5] = parms->lfil[6] = lfil2;
579: pc->setupcalled = 0;
580: }
581: return(0);
582: }
584: /*@
585: PCPARMSSetFill - Sets the fill-in parameters for ILUT, ILUK and ARMS preconditioners.
586: Consider the original matrix A = [B F; E C] and the approximate version
587: M = [LB 0; E/UB I]*[UB LB\F; 0 S].
589: Collective on PC
591: Input Parameters:
592: + pc - the preconditioner context
593: . fil0 - the level of fill-in kept in LB, UB, E/UB and LB\F
594: . fil1 - the level of fill-in kept in S
595: - fil2 - the level of fill-in kept in the L and U parts of the LU factorization of S
597: Options Database Keys:
598: + -pc_parms_lfil_ilu_arms - set the amount of fill-in for ilut, iluk and arms
599: . -pc_parms_lfil_schur - set the amount of fill-in for schur
600: - -pc_parms_lfil_ilut_L_U - set the amount of fill-in for ILUT L and U
602: Level: intermediate
604: Notes:
605: See the pARMS function parms_PCSetFill for more information.
607: .seealso: PCPARMS
608: @*/
609: PetscErrorCode PCPARMSSetFill(PC pc,PetscInt lfil0,PetscInt lfil1,PetscInt lfil2)
610: {
615: PetscTryMethod(pc,"PCPARMSSetFill_C",(PC,PetscInt,PetscInt,PetscInt),(pc,lfil0,lfil1,lfil2));
616: return(0);
617: }
619: /*MC
620: PCPARMS - Allows the use of the parallel Algebraic Recursive Multilevel Solvers
621: available in the package pARMS
623: Options Database Keys:
624: + -pc_parms_global - one of ras, schur, bj
625: . -pc_parms_local - one of ilu0, iluk, ilut, arms
626: . -pc_parms_solve_tol - set the tolerance for local solve
627: . -pc_parms_levels - set the number of levels
628: . -pc_parms_nonsymmetric_perm - set the use of nonsymmetric permutation
629: . -pc_parms_blocksize - set the block size
630: . -pc_parms_ind_tol - set the tolerance for independent sets
631: . -pc_parms_max_dim - set the inner krylov dimension
632: . -pc_parms_max_it - set the maximum number of inner iterations
633: . -pc_parms_inter_nonsymmetric_perm - set the use of nonsymmetric permutation for interlevel blocks
634: . -pc_parms_inter_column_perm - set the use of column permutation for interlevel blocks
635: . -pc_parms_inter_row_scaling - set the use of row scaling for interlevel blocks
636: . -pc_parms_inter_column_scaling - set the use of column scaling for interlevel blocks
637: . -pc_parms_last_nonsymmetric_perm - set the use of nonsymmetric permutation for last level blocks
638: . -pc_parms_last_column_perm - set the use of column permutation for last level blocks
639: . -pc_parms_last_row_scaling - set the use of row scaling for last level blocks
640: . -pc_parms_last_column_scaling - set the use of column scaling for last level blocks
641: . -pc_parms_lfil_ilu_arms - set the amount of fill-in for ilut, iluk and arms
642: . -pc_parms_lfil_schur - set the amount of fill-in for schur
643: . -pc_parms_lfil_ilut_L_U - set the amount of fill-in for ILUT L and U
644: . -pc_parms_droptol_factors - set the drop tolerance for L, U, L^{-1}F and EU^{-1}
645: . -pc_parms_droptol_schur_compl - set the drop tolerance for schur complement at each level
646: - -pc_parms_droptol_last_schur - set the drop tolerance for ILUT in last level schur complement
648: IMPORTANT:
649: Unless configured appropriately, this preconditioner performs an inexact solve
650: as part of the preconditioner application. Therefore, it must be used in combination
651: with flexible variants of iterative solvers, such as KSPFGMRES or KSPCGR.
653: Level: intermediate
655: .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC
656: M*/
658: PETSC_EXTERN PetscErrorCode PCCreate_PARMS(PC pc)
659: {
660: PC_PARMS *parms;
664: PetscNewLog(pc,&parms);
666: parms->map = 0;
667: parms->A = 0;
668: parms->pc = 0;
669: parms->global = PC_PARMS_GLOBAL_RAS;
670: parms->local = PC_PARMS_LOCAL_ARMS;
671: parms->levels = 10;
672: parms->nonsymperm = PETSC_TRUE;
673: parms->blocksize = 250;
674: parms->maxdim = 0;
675: parms->maxits = 0;
676: parms->meth[0] = PETSC_FALSE;
677: parms->meth[1] = PETSC_FALSE;
678: parms->meth[2] = PETSC_FALSE;
679: parms->meth[3] = PETSC_FALSE;
680: parms->meth[4] = PETSC_FALSE;
681: parms->meth[5] = PETSC_FALSE;
682: parms->meth[6] = PETSC_FALSE;
683: parms->meth[7] = PETSC_FALSE;
684: parms->solvetol = 0.01;
685: parms->indtol = 0.4;
686: parms->lfil[0] = parms->lfil[1] = parms->lfil[2] = parms->lfil[3] = 20;
687: parms->lfil[4] = parms->lfil[5] = parms->lfil[6] = 20;
688: parms->droptol[0] = parms->droptol[1] = parms->droptol[2] = parms->droptol[3] = 0.00001;
689: parms->droptol[4] = 0.001;
690: parms->droptol[5] = parms->droptol[6] = 0.001;
692: pc->data = parms;
693: pc->ops->destroy = PCDestroy_PARMS;
694: pc->ops->setfromoptions = PCSetFromOptions_PARMS;
695: pc->ops->setup = PCSetUp_PARMS;
696: pc->ops->apply = PCApply_PARMS;
697: pc->ops->view = PCView_PARMS;
699: PetscObjectComposeFunction((PetscObject)pc,"PCPARMSSetGlobal_C",PCPARMSSetGlobal_PARMS);
700: PetscObjectComposeFunction((PetscObject)pc,"PCPARMSSetLocal_C",PCPARMSSetLocal_PARMS);
701: PetscObjectComposeFunction((PetscObject)pc,"PCPARMSSetSolveTolerances_C",PCPARMSSetSolveTolerances_PARMS);
702: PetscObjectComposeFunction((PetscObject)pc,"PCPARMSSetSolveRestart_C",PCPARMSSetSolveRestart_PARMS);
703: PetscObjectComposeFunction((PetscObject)pc,"PCPARMSSetNonsymPerm_C",PCPARMSSetNonsymPerm_PARMS);
704: PetscObjectComposeFunction((PetscObject)pc,"PCPARMSSetFill_C",PCPARMSSetFill_PARMS);
705: return(0);
706: }