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