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