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