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