Actual source code: parms.c
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;
39: static PetscErrorCode PCSetUp_PARMS(PC pc)
40: {
41: Mat pmat;
42: PC_PARMS *parms = (PC_PARMS*)pc->data;
43: const PetscInt *mapptr0;
44: PetscInt n, lsize, low, high, i, pos, ncols, length;
45: int *maptmp, *mapptr, *ia, *ja, *ja1, *im;
46: PetscScalar *aa, *aa1;
47: const PetscInt *cols;
48: PetscInt meth[8];
49: const PetscScalar *values;
50: PetscErrorCode ierr;
51: MatInfo matinfo;
52: PetscMPIInt rank, npro;
55: /* Get preconditioner matrix from PETSc and setup pARMS structs */
56: PCGetOperators(pc,NULL,&pmat);
57: MPI_Comm_size(PetscObjectComm((PetscObject)pmat),&npro);
58: MPI_Comm_rank(PetscObjectComm((PetscObject)pmat),&rank);
60: MatGetSize(pmat,&n,NULL);
61: PetscMalloc1(npro+1,&mapptr);
62: PetscMalloc1(n,&maptmp);
63: MatGetOwnershipRanges(pmat,&mapptr0);
64: low = mapptr0[rank];
65: high = mapptr0[rank+1];
66: lsize = high - low;
68: for (i=0; i<npro+1; i++) mapptr[i] = mapptr0[i]+1;
69: for (i = 0; i<n; i++) maptmp[i] = i+1;
71: /* if created, destroy the previous map */
72: if (parms->map) {
73: parms_MapFree(&parms->map);
74: parms->map = NULL;
75: }
77: /* create pARMS map object */
78: parms_MapCreateFromPtr(&parms->map,(int)n,maptmp,mapptr,PetscObjectComm((PetscObject)pmat),1,NONINTERLACED);
80: /* if created, destroy the previous pARMS matrix */
81: if (parms->A) {
82: parms_MatFree(&parms->A);
83: parms->A = NULL;
84: }
86: /* create pARMS mat object */
87: parms_MatCreate(&parms->A,parms->map);
89: /* setup and copy csr data structure for pARMS */
90: PetscMalloc1(lsize+1,&ia);
91: ia[0] = 1;
92: MatGetInfo(pmat,MAT_LOCAL,&matinfo);
93: length = matinfo.nz_used;
94: PetscMalloc1(length,&ja);
95: PetscMalloc1(length,&aa);
97: for (i = low; i<high; i++) {
98: pos = ia[i-low]-1;
99: MatGetRow(pmat,i,&ncols,&cols,&values);
100: ia[i-low+1] = ia[i-low] + ncols;
102: if (ia[i-low+1] >= length) {
103: length += ncols;
104: PetscMalloc1(length,&ja1);
105: PetscArraycpy(ja1,ja,ia[i-low]-1);
106: PetscFree(ja);
107: ja = ja1;
108: PetscMalloc1(length,&aa1);
109: PetscArraycpy(aa1,aa,ia[i-low]-1);
110: PetscFree(aa);
111: aa = aa1;
112: }
113: PetscArraycpy(&ja[pos],cols,ncols);
114: PetscArraycpy(&aa[pos],values,ncols);
115: MatRestoreRow(pmat,i,&ncols,&cols,&values);
116: }
118: /* csr info is for local matrix so initialize im[] locally */
119: PetscMalloc1(lsize,&im);
120: PetscArraycpy(im,&maptmp[mapptr[rank]-1],lsize);
122: /* 1-based indexing */
123: for (i=0; i<ia[lsize]-1; i++) ja[i] = ja[i]+1;
125: /* Now copy csr matrix to parms_mat object */
126: parms_MatSetValues(parms->A,(int)lsize,im,ia,ja,aa,INSERT);
128: /* free memory */
129: PetscFree(maptmp);
130: PetscFree(mapptr);
131: PetscFree(aa);
132: PetscFree(ja);
133: PetscFree(ia);
134: PetscFree(im);
136: /* setup parms matrix */
137: parms_MatSetup(parms->A);
139: /* if created, destroy the previous pARMS pc */
140: if (parms->pc) {
141: parms_PCFree(&parms->pc);
142: parms->pc = NULL;
143: }
145: /* Now create pARMS preconditioner object based on A */
146: parms_PCCreate(&parms->pc,parms->A);
148: /* Transfer options from PC to pARMS */
149: switch (parms->global) {
150: case 0: parms_PCSetType(parms->pc, PCRAS); break;
151: case 1: parms_PCSetType(parms->pc, PCSCHUR); break;
152: case 2: parms_PCSetType(parms->pc, PCBJ); break;
153: }
154: switch (parms->local) {
155: case 0: parms_PCSetILUType(parms->pc, PCILU0); break;
156: case 1: parms_PCSetILUType(parms->pc, PCILUK); break;
157: case 2: parms_PCSetILUType(parms->pc, PCILUT); break;
158: case 3: parms_PCSetILUType(parms->pc, PCARMS); break;
159: }
160: parms_PCSetInnerEps(parms->pc, parms->solvetol);
161: parms_PCSetNlevels(parms->pc, parms->levels);
162: parms_PCSetPermType(parms->pc, parms->nonsymperm ? 1 : 0);
163: parms_PCSetBsize(parms->pc, parms->blocksize);
164: parms_PCSetTolInd(parms->pc, parms->indtol);
165: parms_PCSetInnerKSize(parms->pc, parms->maxdim);
166: parms_PCSetInnerMaxits(parms->pc, parms->maxits);
167: for (i=0; i<8; i++) meth[i] = parms->meth[i] ? 1 : 0;
168: parms_PCSetPermScalOptions(parms->pc, &meth[0], 1);
169: parms_PCSetPermScalOptions(parms->pc, &meth[4], 0);
170: parms_PCSetFill(parms->pc, parms->lfil);
171: parms_PCSetTol(parms->pc, parms->droptol);
173: parms_PCSetup(parms->pc);
175: /* Allocate two auxiliary vector of length lsize */
176: if (parms->lvec0) { PetscFree(parms->lvec0); }
177: PetscMalloc1(lsize, &parms->lvec0);
178: if (parms->lvec1) { PetscFree(parms->lvec1); }
179: PetscMalloc1(lsize, &parms->lvec1);
180: return(0);
181: }
183: static PetscErrorCode PCView_PARMS(PC pc,PetscViewer viewer)
184: {
186: PetscBool iascii;
187: PC_PARMS *parms = (PC_PARMS*)pc->data;
188: char *str;
189: double fill_fact;
192: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
193: if (iascii) {
194: parms_PCGetName(parms->pc,&str);
195: PetscViewerASCIIPrintf(viewer," global preconditioner: %s\n",str);
196: parms_PCILUGetName(parms->pc,&str);
197: PetscViewerASCIIPrintf(viewer," local preconditioner: %s\n",str);
198: parms_PCGetRatio(parms->pc,&fill_fact);
199: PetscViewerASCIIPrintf(viewer," non-zero elements/original non-zero entries: %-4.2f\n",fill_fact);
200: PetscViewerASCIIPrintf(viewer," Tolerance for local solve: %g\n",parms->solvetol);
201: PetscViewerASCIIPrintf(viewer," Number of levels: %d\n",parms->levels);
202: if (parms->nonsymperm) {
203: PetscViewerASCIIPrintf(viewer," Using nonsymmetric permutation\n");
204: }
205: PetscViewerASCIIPrintf(viewer," Block size: %d\n",parms->blocksize);
206: PetscViewerASCIIPrintf(viewer," Tolerance for independent sets: %g\n",parms->indtol);
207: PetscViewerASCIIPrintf(viewer," Inner Krylov dimension: %d\n",parms->maxdim);
208: PetscViewerASCIIPrintf(viewer," Maximum number of inner iterations: %d\n",parms->maxits);
209: if (parms->meth[0]) {
210: PetscViewerASCIIPrintf(viewer," Using nonsymmetric permutation for interlevel blocks\n");
211: }
212: if (parms->meth[1]) {
213: PetscViewerASCIIPrintf(viewer," Using column permutation for interlevel blocks\n");
214: }
215: if (parms->meth[2]) {
216: PetscViewerASCIIPrintf(viewer," Using row scaling for interlevel blocks\n");
217: }
218: if (parms->meth[3]) {
219: PetscViewerASCIIPrintf(viewer," Using column scaling for interlevel blocks\n");
220: }
221: if (parms->meth[4]) {
222: PetscViewerASCIIPrintf(viewer," Using nonsymmetric permutation for last level blocks\n");
223: }
224: if (parms->meth[5]) {
225: PetscViewerASCIIPrintf(viewer," Using column permutation for last level blocks\n");
226: }
227: if (parms->meth[6]) {
228: PetscViewerASCIIPrintf(viewer," Using row scaling for last level blocks\n");
229: }
230: if (parms->meth[7]) {
231: PetscViewerASCIIPrintf(viewer," Using column scaling for last level blocks\n");
232: }
233: PetscViewerASCIIPrintf(viewer," amount of fill-in for ilut, iluk and arms: %d\n",parms->lfil[0]);
234: PetscViewerASCIIPrintf(viewer," amount of fill-in for schur: %d\n",parms->lfil[4]);
235: PetscViewerASCIIPrintf(viewer," amount of fill-in for ILUT L and U: %d\n",parms->lfil[5]);
236: PetscViewerASCIIPrintf(viewer," drop tolerance for L, U, L^{-1}F and EU^{-1}: %g\n",parms->droptol[0]);
237: PetscViewerASCIIPrintf(viewer," drop tolerance for schur complement at each level: %g\n",parms->droptol[4]);
238: PetscViewerASCIIPrintf(viewer," drop tolerance for ILUT in last level schur complement: %g\n",parms->droptol[5]);
239: }
240: return(0);
241: }
243: static PetscErrorCode PCDestroy_PARMS(PC pc)
244: {
245: PC_PARMS *parms = (PC_PARMS*)pc->data;
249: if (parms->map) parms_MapFree(&parms->map);
250: if (parms->A) parms_MatFree(&parms->A);
251: if (parms->pc) parms_PCFree(&parms->pc);
252: if (parms->lvec0) {
253: PetscFree(parms->lvec0);
254: }
255: if (parms->lvec1) {
256: PetscFree(parms->lvec1);
257: }
258: PetscFree(pc->data);
260: PetscObjectChangeTypeName((PetscObject)pc,0);
261: PetscObjectComposeFunction((PetscObject)pc,"PCPARMSSetGlobal_C",NULL);
262: PetscObjectComposeFunction((PetscObject)pc,"PCPARMSSetLocal_C",NULL);
263: PetscObjectComposeFunction((PetscObject)pc,"PCPARMSSetSolveTolerances_C",NULL);
264: PetscObjectComposeFunction((PetscObject)pc,"PCPARMSSetSolveRestart_C",NULL);
265: PetscObjectComposeFunction((PetscObject)pc,"PCPARMSSetNonsymPerm_C",NULL);
266: PetscObjectComposeFunction((PetscObject)pc,"PCPARMSSetFill_C",NULL);
267: return(0);
268: }
270: static PetscErrorCode PCSetFromOptions_PARMS(PetscOptionItems *PetscOptionsObject,PC pc)
271: {
272: PC_PARMS *parms = (PC_PARMS*)pc->data;
273: PetscBool flag;
274: PCPARMSGlobalType global;
275: PCPARMSLocalType local;
276: PetscErrorCode ierr;
279: PetscOptionsHead(PetscOptionsObject,"PARMS Options");
280: PetscOptionsEnum("-pc_parms_global","Global preconditioner","PCPARMSSetGlobal",PCPARMSGlobalTypes,(PetscEnum)parms->global,(PetscEnum*)&global,&flag);
281: if (flag) {PCPARMSSetGlobal(pc,global);}
282: PetscOptionsEnum("-pc_parms_local","Local preconditioner","PCPARMSSetLocal",PCPARMSLocalTypes,(PetscEnum)parms->local,(PetscEnum*)&local,&flag);
283: if (flag) {PCPARMSSetLocal(pc,local);}
284: PetscOptionsReal("-pc_parms_solve_tol","Tolerance for local solve","PCPARMSSetSolveTolerances",parms->solvetol,&parms->solvetol,NULL);
285: PetscOptionsInt("-pc_parms_levels","Number of levels","None",parms->levels,&parms->levels,NULL);
286: PetscOptionsBool("-pc_parms_nonsymmetric_perm","Use nonsymmetric permutation","PCPARMSSetNonsymPerm",parms->nonsymperm,&parms->nonsymperm,NULL);
287: PetscOptionsInt("-pc_parms_blocksize","Block size","None",parms->blocksize,&parms->blocksize,NULL);
288: PetscOptionsReal("-pc_parms_ind_tol","Tolerance for independent sets","None",parms->indtol,&parms->indtol,NULL);
289: PetscOptionsInt("-pc_parms_max_dim","Inner Krylov dimension","PCPARMSSetSolveRestart",parms->maxdim,&parms->maxdim,NULL);
290: PetscOptionsInt("-pc_parms_max_it","Maximum number of inner iterations","PCPARMSSetSolveTolerances",parms->maxits,&parms->maxits,NULL);
291: PetscOptionsBool("-pc_parms_inter_nonsymmetric_perm","nonsymmetric permutation for interlevel blocks","None",parms->meth[0],&parms->meth[0],NULL);
292: PetscOptionsBool("-pc_parms_inter_column_perm","column permutation for interlevel blocks","None",parms->meth[1],&parms->meth[1],NULL);
293: PetscOptionsBool("-pc_parms_inter_row_scaling","row scaling for interlevel blocks","None",parms->meth[2],&parms->meth[2],NULL);
294: PetscOptionsBool("-pc_parms_inter_column_scaling","column scaling for interlevel blocks","None",parms->meth[3],&parms->meth[3],NULL);
295: PetscOptionsBool("-pc_parms_last_nonsymmetric_perm","nonsymmetric permutation for last level blocks","None",parms->meth[4],&parms->meth[4],NULL);
296: PetscOptionsBool("-pc_parms_last_column_perm","column permutation for last level blocks","None",parms->meth[5],&parms->meth[5],NULL);
297: PetscOptionsBool("-pc_parms_last_row_scaling","row scaling for last level blocks","None",parms->meth[6],&parms->meth[6],NULL);
298: PetscOptionsBool("-pc_parms_last_column_scaling","column scaling for last level blocks","None",parms->meth[7],&parms->meth[7],NULL);
299: PetscOptionsInt("-pc_parms_lfil_ilu_arms","amount of fill-in for ilut, iluk and arms","PCPARMSSetFill",parms->lfil[0],&parms->lfil[0],&flag);
300: if (flag) parms->lfil[1] = parms->lfil[2] = parms->lfil[3] = parms->lfil[0];
301: PetscOptionsInt("-pc_parms_lfil_schur","amount of fill-in for schur","PCPARMSSetFill",parms->lfil[4],&parms->lfil[4],NULL);
302: PetscOptionsInt("-pc_parms_lfil_ilut_L_U","amount of fill-in for ILUT L and U","PCPARMSSetFill",parms->lfil[5],&parms->lfil[5],&flag);
303: if (flag) parms->lfil[6] = parms->lfil[5];
304: PetscOptionsReal("-pc_parms_droptol_factors","drop tolerance for L, U, L^{-1}F and EU^{-1}","None",parms->droptol[0],&parms->droptol[0],NULL);
305: PetscOptionsReal("-pc_parms_droptol_schur_compl","drop tolerance for schur complement at each level","None",parms->droptol[4],&parms->droptol[4],&flag);
306: if (flag) parms->droptol[1] = parms->droptol[2] = parms->droptol[3] = parms->droptol[0];
307: PetscOptionsReal("-pc_parms_droptol_last_schur","drop tolerance for ILUT in last level schur complement","None",parms->droptol[5],&parms->droptol[5],&flag);
308: if (flag) parms->droptol[6] = parms->droptol[5];
309: PetscOptionsTail();
310: return(0);
311: }
313: static PetscErrorCode PCApply_PARMS(PC pc,Vec b,Vec x)
314: {
315: PetscErrorCode ierr;
316: PC_PARMS *parms = (PC_PARMS*)pc->data;
317: const PetscScalar *b1;
318: PetscScalar *x1;
321: VecGetArrayRead(b,&b1);
322: VecGetArray(x,&x1);
323: parms_VecPermAux((PetscScalar*)b1,parms->lvec0,parms->map);
324: parms_PCApply(parms->pc,parms->lvec0,parms->lvec1);
325: parms_VecInvPermAux(parms->lvec1,x1,parms->map);
326: VecRestoreArrayRead(b,&b1);
327: VecRestoreArray(x,&x1);
328: return(0);
329: }
331: static PetscErrorCode PCPARMSSetGlobal_PARMS(PC pc,PCPARMSGlobalType type)
332: {
333: PC_PARMS *parms = (PC_PARMS*)pc->data;
336: if (type != parms->global) {
337: parms->global = type;
338: pc->setupcalled = 0;
339: }
340: return(0);
341: }
343: /*@
344: PCPARMSSetGlobal - Sets the global preconditioner to be used in PARMS.
346: Collective on PC
348: Input Parameters:
349: + pc - the preconditioner context
350: - type - the global preconditioner type, one of
351: .vb
352: PC_PARMS_GLOBAL_RAS - Restricted additive Schwarz
353: PC_PARMS_GLOBAL_SCHUR - Schur complement
354: PC_PARMS_GLOBAL_BJ - Block Jacobi
355: .ve
357: Options Database Keys:
358: -pc_parms_global [ras,schur,bj] - Sets global preconditioner
360: Level: intermediate
362: Notes:
363: See the pARMS function parms_PCSetType for more information.
365: .seealso: PCPARMS, PCPARMSSetLocal()
366: @*/
367: PetscErrorCode PCPARMSSetGlobal(PC pc,PCPARMSGlobalType type)
368: {
374: PetscTryMethod(pc,"PCPARMSSetGlobal_C",(PC,PCPARMSGlobalType),(pc,type));
375: return(0);
376: }
378: static PetscErrorCode PCPARMSSetLocal_PARMS(PC pc,PCPARMSLocalType type)
379: {
380: PC_PARMS *parms = (PC_PARMS*)pc->data;
383: if (type != parms->local) {
384: parms->local = type;
385: pc->setupcalled = 0;
386: }
387: return(0);
388: }
390: /*@
391: PCPARMSSetLocal - Sets the local preconditioner to be used in PARMS.
393: Collective on PC
395: Input Parameters:
396: + pc - the preconditioner context
397: - type - the local preconditioner type, one of
398: .vb
399: PC_PARMS_LOCAL_ILU0 - ILU0 preconditioner
400: PC_PARMS_LOCAL_ILUK - ILU(k) preconditioner
401: PC_PARMS_LOCAL_ILUT - ILUT preconditioner
402: PC_PARMS_LOCAL_ARMS - ARMS preconditioner
403: .ve
405: Options Database Keys:
406: -pc_parms_local [ilu0,iluk,ilut,arms] - Sets local preconditioner
408: Level: intermediate
410: Notes:
411: For the ARMS preconditioner, one can use either the symmetric ARMS or the non-symmetric
412: variant (ARMS-ddPQ) by setting the permutation type with PCPARMSSetNonsymPerm().
414: See the pARMS function parms_PCILUSetType for more information.
416: .seealso: PCPARMS, PCPARMSSetGlobal(), PCPARMSSetNonsymPerm()
418: @*/
419: PetscErrorCode PCPARMSSetLocal(PC pc,PCPARMSLocalType type)
420: {
426: PetscTryMethod(pc,"PCPARMSSetLocal_C",(PC,PCPARMSLocalType),(pc,type));
427: return(0);
428: }
430: static PetscErrorCode PCPARMSSetSolveTolerances_PARMS(PC pc,PetscReal tol,PetscInt maxits)
431: {
432: PC_PARMS *parms = (PC_PARMS*)pc->data;
435: if (tol != parms->solvetol) {
436: parms->solvetol = tol;
437: pc->setupcalled = 0;
438: }
439: if (maxits != parms->maxits) {
440: parms->maxits = maxits;
441: pc->setupcalled = 0;
442: }
443: return(0);
444: }
446: /*@
447: PCPARMSSetSolveTolerances - Sets the convergence tolerance and the maximum iterations for the
448: inner GMRES solver, when the Schur global preconditioner is used.
450: Collective on PC
452: Input Parameters:
453: + pc - the preconditioner context
454: . tol - the convergence tolerance
455: - maxits - the maximum number of iterations to use
457: Options Database Keys:
458: + -pc_parms_solve_tol - set the tolerance for local solve
459: - -pc_parms_max_it - set the maximum number of inner iterations
461: Level: intermediate
463: Notes:
464: See the pARMS functions parms_PCSetInnerEps and parms_PCSetInnerMaxits for more information.
466: .seealso: PCPARMS, PCPARMSSetSolveRestart()
467: @*/
468: PetscErrorCode PCPARMSSetSolveTolerances(PC pc,PetscReal tol,PetscInt maxits)
469: {
474: PetscTryMethod(pc,"PCPARMSSetSolveTolerances_C",(PC,PetscReal,PetscInt),(pc,tol,maxits));
475: return(0);
476: }
478: static PetscErrorCode PCPARMSSetSolveRestart_PARMS(PC pc,PetscInt restart)
479: {
480: PC_PARMS *parms = (PC_PARMS*)pc->data;
483: if (restart != parms->maxdim) {
484: parms->maxdim = restart;
485: pc->setupcalled = 0;
486: }
487: return(0);
488: }
490: /*@
491: PCPARMSSetSolveRestart - Sets the number of iterations at which the
492: inner GMRES solver restarts.
494: Collective on PC
496: Input Parameters:
497: + pc - the preconditioner context
498: - restart - maximum dimension of the Krylov subspace
500: Options Database Keys:
501: . -pc_parms_max_dim - sets the inner Krylov dimension
503: Level: intermediate
505: Notes:
506: See the pARMS function parms_PCSetInnerKSize for more information.
508: .seealso: PCPARMS, PCPARMSSetSolveTolerances()
509: @*/
510: PetscErrorCode PCPARMSSetSolveRestart(PC pc,PetscInt restart)
511: {
516: PetscTryMethod(pc,"PCPARMSSetSolveRestart_C",(PC,PetscInt),(pc,restart));
517: return(0);
518: }
520: static PetscErrorCode PCPARMSSetNonsymPerm_PARMS(PC pc,PetscBool nonsym)
521: {
522: PC_PARMS *parms = (PC_PARMS*)pc->data;
525: if ((nonsym && !parms->nonsymperm) || (!nonsym && parms->nonsymperm)) {
526: parms->nonsymperm = nonsym;
527: pc->setupcalled = 0;
528: }
529: return(0);
530: }
532: /*@
533: PCPARMSSetNonsymPerm - Sets the type of permutation for the ARMS preconditioner: the standard
534: symmetric ARMS or the non-symmetric ARMS (ARMS-ddPQ).
536: Collective on PC
538: Input Parameters:
539: + pc - the preconditioner context
540: - nonsym - PETSC_TRUE indicates the non-symmetric ARMS is used;
541: PETSC_FALSE indicates the symmetric ARMS is used
543: Options Database Keys:
544: . -pc_parms_nonsymmetric_perm - sets the use of nonsymmetric permutation
546: Level: intermediate
548: Notes:
549: See the pARMS function parms_PCSetPermType for more information.
551: .seealso: PCPARMS
552: @*/
553: PetscErrorCode PCPARMSSetNonsymPerm(PC pc,PetscBool nonsym)
554: {
559: PetscTryMethod(pc,"PCPARMSSetNonsymPerm_C",(PC,PetscBool),(pc,nonsym));
560: return(0);
561: }
563: static PetscErrorCode PCPARMSSetFill_PARMS(PC pc,PetscInt lfil0,PetscInt lfil1,PetscInt lfil2)
564: {
565: PC_PARMS *parms = (PC_PARMS*)pc->data;
568: if (lfil0 != parms->lfil[0] || lfil0 != parms->lfil[1] || lfil0 != parms->lfil[2] || lfil0 != parms->lfil[3]) {
569: parms->lfil[1] = parms->lfil[2] = parms->lfil[3] = parms->lfil[0] = lfil0;
570: pc->setupcalled = 0;
571: }
572: if (lfil1 != parms->lfil[4]) {
573: parms->lfil[4] = lfil1;
574: pc->setupcalled = 0;
575: }
576: if (lfil2 != parms->lfil[5] || lfil2 != parms->lfil[6]) {
577: parms->lfil[5] = parms->lfil[6] = lfil2;
578: pc->setupcalled = 0;
579: }
580: return(0);
581: }
583: /*@
584: PCPARMSSetFill - Sets the fill-in parameters for ILUT, ILUK and ARMS preconditioners.
585: Consider the original matrix A = [B F; E C] and the approximate version
586: M = [LB 0; E/UB I]*[UB LB\F; 0 S].
588: Collective on PC
590: Input Parameters:
591: + pc - the preconditioner context
592: . fil0 - the level of fill-in kept in LB, UB, E/UB and LB\F
593: . fil1 - the level of fill-in kept in S
594: - fil2 - the level of fill-in kept in the L and U parts of the LU factorization of S
596: Options Database Keys:
597: + -pc_parms_lfil_ilu_arms - set the amount of fill-in for ilut, iluk and arms
598: . -pc_parms_lfil_schur - set the amount of fill-in for schur
599: - -pc_parms_lfil_ilut_L_U - set the amount of fill-in for ILUT L and U
601: Level: intermediate
603: Notes:
604: See the pARMS function parms_PCSetFill for more information.
606: .seealso: PCPARMS
607: @*/
608: PetscErrorCode PCPARMSSetFill(PC pc,PetscInt lfil0,PetscInt lfil1,PetscInt lfil2)
609: {
614: PetscTryMethod(pc,"PCPARMSSetFill_C",(PC,PetscInt,PetscInt,PetscInt),(pc,lfil0,lfil1,lfil2));
615: return(0);
616: }
618: /*MC
619: PCPARMS - Allows the use of the parallel Algebraic Recursive Multilevel Solvers
620: available in the package pARMS
622: Options Database Keys:
623: + -pc_parms_global - one of ras, schur, bj
624: . -pc_parms_local - one of ilu0, iluk, ilut, arms
625: . -pc_parms_solve_tol - set the tolerance for local solve
626: . -pc_parms_levels - set the number of levels
627: . -pc_parms_nonsymmetric_perm - set the use of nonsymmetric permutation
628: . -pc_parms_blocksize - set the block size
629: . -pc_parms_ind_tol - set the tolerance for independent sets
630: . -pc_parms_max_dim - set the inner krylov dimension
631: . -pc_parms_max_it - set the maximum number of inner iterations
632: . -pc_parms_inter_nonsymmetric_perm - set the use of nonsymmetric permutation for interlevel blocks
633: . -pc_parms_inter_column_perm - set the use of column permutation for interlevel blocks
634: . -pc_parms_inter_row_scaling - set the use of row scaling for interlevel blocks
635: . -pc_parms_inter_column_scaling - set the use of column scaling for interlevel blocks
636: . -pc_parms_last_nonsymmetric_perm - set the use of nonsymmetric permutation for last level blocks
637: . -pc_parms_last_column_perm - set the use of column permutation for last level blocks
638: . -pc_parms_last_row_scaling - set the use of row scaling for last level blocks
639: . -pc_parms_last_column_scaling - set the use of column scaling for last level blocks
640: . -pc_parms_lfil_ilu_arms - set the amount of fill-in for ilut, iluk and arms
641: . -pc_parms_lfil_schur - set the amount of fill-in for schur
642: . -pc_parms_lfil_ilut_L_U - set the amount of fill-in for ILUT L and U
643: . -pc_parms_droptol_factors - set the drop tolerance for L, U, L^{-1}F and EU^{-1}
644: . -pc_parms_droptol_schur_compl - set the drop tolerance for schur complement at each level
645: - -pc_parms_droptol_last_schur - set the drop tolerance for ILUT in last level schur complement
647: IMPORTANT:
648: Unless configured appropriately, this preconditioner performs an inexact solve
649: as part of the preconditioner application. Therefore, it must be used in combination
650: with flexible variants of iterative solvers, such as KSPFGMRES or KSPGCR.
652: Level: intermediate
654: .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC
655: M*/
657: PETSC_EXTERN PetscErrorCode PCCreate_PARMS(PC pc)
658: {
659: PC_PARMS *parms;
663: PetscNewLog(pc,&parms);
665: parms->map = 0;
666: parms->A = 0;
667: parms->pc = 0;
668: parms->global = PC_PARMS_GLOBAL_RAS;
669: parms->local = PC_PARMS_LOCAL_ARMS;
670: parms->levels = 10;
671: parms->nonsymperm = PETSC_TRUE;
672: parms->blocksize = 250;
673: parms->maxdim = 0;
674: parms->maxits = 0;
675: parms->meth[0] = PETSC_FALSE;
676: parms->meth[1] = PETSC_FALSE;
677: parms->meth[2] = PETSC_FALSE;
678: parms->meth[3] = PETSC_FALSE;
679: parms->meth[4] = PETSC_FALSE;
680: parms->meth[5] = PETSC_FALSE;
681: parms->meth[6] = PETSC_FALSE;
682: parms->meth[7] = PETSC_FALSE;
683: parms->solvetol = 0.01;
684: parms->indtol = 0.4;
685: parms->lfil[0] = parms->lfil[1] = parms->lfil[2] = parms->lfil[3] = 20;
686: parms->lfil[4] = parms->lfil[5] = parms->lfil[6] = 20;
687: parms->droptol[0] = parms->droptol[1] = parms->droptol[2] = parms->droptol[3] = 0.00001;
688: parms->droptol[4] = 0.001;
689: parms->droptol[5] = parms->droptol[6] = 0.001;
691: pc->data = parms;
692: pc->ops->destroy = PCDestroy_PARMS;
693: pc->ops->setfromoptions = PCSetFromOptions_PARMS;
694: pc->ops->setup = PCSetUp_PARMS;
695: pc->ops->apply = PCApply_PARMS;
696: pc->ops->view = PCView_PARMS;
698: PetscObjectComposeFunction((PetscObject)pc,"PCPARMSSetGlobal_C",PCPARMSSetGlobal_PARMS);
699: PetscObjectComposeFunction((PetscObject)pc,"PCPARMSSetLocal_C",PCPARMSSetLocal_PARMS);
700: PetscObjectComposeFunction((PetscObject)pc,"PCPARMSSetSolveTolerances_C",PCPARMSSetSolveTolerances_PARMS);
701: PetscObjectComposeFunction((PetscObject)pc,"PCPARMSSetSolveRestart_C",PCPARMSSetSolveRestart_PARMS);
702: PetscObjectComposeFunction((PetscObject)pc,"PCPARMSSetNonsymPerm_C",PCPARMSSetNonsymPerm_PARMS);
703: PetscObjectComposeFunction((PetscObject)pc,"PCPARMSSetFill_C",PCPARMSSetFill_PARMS);
704: return(0);
705: }