Actual source code: parms.c

petsc-3.6.1 2015-08-06
Report Typos and Errors
  1: #define PETSCKSP_DLL

  3: /*
  4:    Provides an interface to pARMS.
  5:    Requires pARMS 3.2 or later.
  6: */

  8: #include <petsc/private/pcimpl.h>          /*I "petscpc.h" I*/

 10: #if defined(PETSC_USE_COMPLEX)
 11: #define DBL_CMPLX
 12: #else
 13: #define DBL
 14: #endif
 15: #define USE_MPI
 16: #define REAL double
 17: #define HAS_BLAS
 18: #define FORTRAN_UNDERSCORE
 19: #include "parms_sys.h"
 20: #undef FLOAT
 21: #define FLOAT PetscScalar
 22: #include <parms.h>

 24: /*
 25:    Private context (data structure) for the  preconditioner.
 26: */
 27: typedef struct {
 28:   parms_Map         map;
 29:   parms_Mat         A;
 30:   parms_PC          pc;
 31:   PCPARMSGlobalType global;
 32:   PCPARMSLocalType  local;
 33:   PetscInt          levels, blocksize, maxdim, maxits, lfil[7];
 34:   PetscBool         nonsymperm, meth[8];
 35:   PetscReal         solvetol, indtol, droptol[7];
 36:   PetscScalar       *lvec0, *lvec1;
 37: } PC_PARMS;


 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;

 58:   /* Get preconditioner matrix from PETSc and setup pARMS structs */
 59:   PCGetOperators(pc,NULL,&pmat);
 60:   MPI_Comm_size(PetscObjectComm((PetscObject)pmat),&npro);
 61:   MPI_Comm_rank(PetscObjectComm((PetscObject)pmat),&rank);

 63:   MatGetSize(pmat,&n,NULL);
 64:   PetscMalloc1(npro+1,&mapptr);
 65:   PetscMalloc1(n,&maptmp);
 66:   MatGetOwnershipRanges(pmat,&mapptr0);
 67:   low   = mapptr0[rank];
 68:   high  = mapptr0[rank+1];
 69:   lsize = high - low;

 71:   for (i=0; i<npro+1; i++) mapptr[i] = mapptr0[i]+1;
 72:   for (i = 0; i<n; i++) maptmp[i] = i+1;

 74:   /* if created, destroy the previous map */
 75:   if (parms->map) {
 76:     parms_MapFree(&parms->map);
 77:     parms->map = NULL;
 78:   }

 80:   /* create pARMS map object */
 81:   parms_MapCreateFromPtr(&parms->map,(int)n,maptmp,mapptr,PetscObjectComm((PetscObject)pmat),1,NONINTERLACED);

 83:   /* if created, destroy the previous pARMS matrix */
 84:   if (parms->A) {
 85:     parms_MatFree(&parms->A);
 86:     parms->A = NULL;
 87:   }

 89:   /* create pARMS mat object */
 90:   parms_MatCreate(&parms->A,parms->map);

 92:   /* setup and copy csr data structure for pARMS */
 93:   PetscMalloc1(lsize+1,&ia);
 94:   ia[0]  = 1;
 95:   MatGetInfo(pmat,MAT_LOCAL,&matinfo);
 96:   length = matinfo.nz_used;
 97:   PetscMalloc1(length,&ja);
 98:   PetscMalloc1(length,&aa);

100:   for (i = low; i<high; i++) {
101:     pos         = ia[i-low]-1;
102:     MatGetRow(pmat,i,&ncols,&cols,&values);
103:     ia[i-low+1] = ia[i-low] + ncols;

105:     if (ia[i-low+1] >= length) {
106:       length += ncols;
107:       PetscMalloc1(length,&ja1);
108:       PetscMemcpy(ja1,ja,(ia[i-low]-1)*sizeof(int));
109:       PetscFree(ja);
110:       ja      = ja1;
111:       PetscMalloc1(length,&aa1);
112:       PetscMemcpy(aa1,aa,(ia[i-low]-1)*sizeof(PetscScalar));
113:       PetscFree(aa);
114:       aa      = aa1;
115:     }
116:     PetscMemcpy(&ja[pos],cols,ncols*sizeof(int));
117:     PetscMemcpy(&aa[pos],values,ncols*sizeof(PetscScalar));
118:     MatRestoreRow(pmat,i,&ncols,&cols,&values);
119:   }

121:   /* csr info is for local matrix so initialize im[] locally */
122:   PetscMalloc1(lsize,&im);
123:   PetscMemcpy(im,&maptmp[mapptr[rank]-1],lsize*sizeof(int));

125:   /* 1-based indexing */
126:   for (i=0; i<ia[lsize]-1; i++) ja[i] = ja[i]+1;

128:   /* Now copy csr matrix to parms_mat object */
129:   parms_MatSetValues(parms->A,(int)lsize,im,ia,ja,aa,INSERT);

131:   /* free memory */
132:   PetscFree(maptmp);
133:   PetscFree(mapptr);
134:   PetscFree(aa);
135:   PetscFree(ja);
136:   PetscFree(ia);
137:   PetscFree(im);

139:   /* setup parms matrix */
140:   parms_MatSetup(parms->A);

142:   /* if created, destroy the previous pARMS pc */
143:   if (parms->pc) {
144:     parms_PCFree(&parms->pc);
145:     parms->pc = NULL;
146:   }

148:   /* Now create pARMS preconditioner object based on A */
149:   parms_PCCreate(&parms->pc,parms->A);

151:   /* Transfer options from PC to pARMS */
152:   switch (parms->global) {
153:   case 0: parms_PCSetType(parms->pc, PCRAS); break;
154:   case 1: parms_PCSetType(parms->pc, PCSCHUR); break;
155:   case 2: parms_PCSetType(parms->pc, PCBJ); break;
156:   }
157:   switch (parms->local) {
158:   case 0: parms_PCSetILUType(parms->pc, PCILU0); break;
159:   case 1: parms_PCSetILUType(parms->pc, PCILUK); break;
160:   case 2: parms_PCSetILUType(parms->pc, PCILUT); break;
161:   case 3: parms_PCSetILUType(parms->pc, PCARMS); break;
162:   }
163:   parms_PCSetInnerEps(parms->pc, parms->solvetol);
164:   parms_PCSetNlevels(parms->pc, parms->levels);
165:   parms_PCSetPermType(parms->pc, parms->nonsymperm ? 1 : 0);
166:   parms_PCSetBsize(parms->pc, parms->blocksize);
167:   parms_PCSetTolInd(parms->pc, parms->indtol);
168:   parms_PCSetInnerKSize(parms->pc, parms->maxdim);
169:   parms_PCSetInnerMaxits(parms->pc, parms->maxits);
170:   for (i=0; i<8; i++) meth[i] = parms->meth[i] ? 1 : 0;
171:   parms_PCSetPermScalOptions(parms->pc, &meth[0], 1);
172:   parms_PCSetPermScalOptions(parms->pc, &meth[4], 0);
173:   parms_PCSetFill(parms->pc, parms->lfil);
174:   parms_PCSetTol(parms->pc, parms->droptol);

176:   parms_PCSetup(parms->pc);

178:   /* Allocate two auxiliary vector of length lsize */
179:   if (parms->lvec0) { PetscFree(parms->lvec0); }
180:   PetscMalloc1(lsize, &parms->lvec0);
181:   if (parms->lvec1) { PetscFree(parms->lvec1); }
182:   PetscMalloc1(lsize, &parms->lvec1);
183:   return(0);
184: }

188: static PetscErrorCode PCView_PARMS(PC pc,PetscViewer viewer)
189: {
191:   PetscBool      iascii;
192:   PC_PARMS       *parms = (PC_PARMS*)pc->data;
193:   char           *str;
194:   double         fill_fact;

197:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
198:   if (iascii) {
199:     parms_PCGetName(parms->pc,&str);
200:     PetscViewerASCIIPrintf(viewer,"  global preconditioner: %s\n",str);
201:     parms_PCILUGetName(parms->pc,&str);
202:     PetscViewerASCIIPrintf(viewer,"  local preconditioner: %s\n",str);
203:     parms_PCGetRatio(parms->pc,&fill_fact);
204:     PetscViewerASCIIPrintf(viewer,"  non-zero elements/original non-zero entries: %-4.2f\n",fill_fact);
205:     PetscViewerASCIIPrintf(viewer,"  Tolerance for local solve: %g\n",parms->solvetol);
206:     PetscViewerASCIIPrintf(viewer,"  Number of levels: %d\n",parms->levels);
207:     if (parms->nonsymperm) {
208:       PetscViewerASCIIPrintf(viewer,"  Using nonsymmetric permutation\n");
209:     }
210:     PetscViewerASCIIPrintf(viewer,"  Block size: %d\n",parms->blocksize);
211:     PetscViewerASCIIPrintf(viewer,"  Tolerance for independent sets: %g\n",parms->indtol);
212:     PetscViewerASCIIPrintf(viewer,"  Inner Krylov dimension: %d\n",parms->maxdim);
213:     PetscViewerASCIIPrintf(viewer,"  Maximum number of inner iterations: %d\n",parms->maxits);
214:     if (parms->meth[0]) {
215:       PetscViewerASCIIPrintf(viewer,"  Using nonsymmetric permutation for interlevel blocks\n");
216:     }
217:     if (parms->meth[1]) {
218:       PetscViewerASCIIPrintf(viewer,"  Using column permutation for interlevel blocks\n");
219:     }
220:     if (parms->meth[2]) {
221:       PetscViewerASCIIPrintf(viewer,"  Using row scaling for interlevel blocks\n");
222:     }
223:     if (parms->meth[3]) {
224:       PetscViewerASCIIPrintf(viewer,"  Using column scaling for interlevel blocks\n");
225:     }
226:     if (parms->meth[4]) {
227:       PetscViewerASCIIPrintf(viewer,"  Using nonsymmetric permutation for last level blocks\n");
228:     }
229:     if (parms->meth[5]) {
230:       PetscViewerASCIIPrintf(viewer,"  Using column permutation for last level blocks\n");
231:     }
232:     if (parms->meth[6]) {
233:       PetscViewerASCIIPrintf(viewer,"  Using row scaling for last level blocks\n");
234:     }
235:     if (parms->meth[7]) {
236:       PetscViewerASCIIPrintf(viewer,"  Using column scaling for last level blocks\n");
237:     }
238:     PetscViewerASCIIPrintf(viewer,"  amount of fill-in for ilut, iluk and arms: %d\n",parms->lfil[0]);
239:     PetscViewerASCIIPrintf(viewer,"  amount of fill-in for schur: %d\n",parms->lfil[4]);
240:     PetscViewerASCIIPrintf(viewer,"  amount of fill-in for ILUT L and U: %d\n",parms->lfil[5]);
241:     PetscViewerASCIIPrintf(viewer,"  drop tolerance for L, U, L^{-1}F and EU^{-1}: %g\n",parms->droptol[0]);
242:     PetscViewerASCIIPrintf(viewer,"  drop tolerance for schur complement at each level: %g\n",parms->droptol[4]);
243:     PetscViewerASCIIPrintf(viewer,"  drop tolerance for ILUT in last level schur complement: %g\n",parms->droptol[5]);
244:   }
245:   return(0);
246: }

250: static PetscErrorCode PCDestroy_PARMS(PC pc)
251: {
252:   PC_PARMS       *parms = (PC_PARMS*)pc->data;

256:   if (parms->map) parms_MapFree(&parms->map);
257:   if (parms->A) parms_MatFree(&parms->A);
258:   if (parms->pc) parms_PCFree(&parms->pc);
259:   if (parms->lvec0) {
260:     PetscFree(parms->lvec0);
261:   }
262:   if (parms->lvec1) {
263:     PetscFree(parms->lvec1);
264:   }
265:   PetscFree(pc->data);

267:   PetscObjectChangeTypeName((PetscObject)pc,0);
268:   PetscObjectComposeFunction((PetscObject)pc,"PCPARMSSetGlobal_C",NULL);
269:   PetscObjectComposeFunction((PetscObject)pc,"PCPARMSSetLocal_C",NULL);
270:   PetscObjectComposeFunction((PetscObject)pc,"PCPARMSSetSolveTolerances_C",NULL);
271:   PetscObjectComposeFunction((PetscObject)pc,"PCPARMSSetSolveRestart_C",NULL);
272:   PetscObjectComposeFunction((PetscObject)pc,"PCPARMSSetNonsymPerm_C",NULL);
273:   PetscObjectComposeFunction((PetscObject)pc,"PCPARMSSetFill_C",NULL);
274:   return(0);
275: }

279: static PetscErrorCode PCSetFromOptions_PARMS(PetscOptions *PetscOptionsObject,PC pc)
280: {
281:   PC_PARMS          *parms = (PC_PARMS*)pc->data;
282:   PetscBool         flag;
283:   PCPARMSGlobalType global;
284:   PCPARMSLocalType  local;
285:   PetscErrorCode    ierr;

288:   PetscOptionsHead(PetscOptionsObject,"PARMS Options");
289:   PetscOptionsEnum("-pc_parms_global","Global preconditioner","PCPARMSSetGlobal",PCPARMSGlobalTypes,(PetscEnum)parms->global,(PetscEnum*)&global,&flag);
290:   if (flag) {PCPARMSSetGlobal(pc,global);}
291:   PetscOptionsEnum("-pc_parms_local","Local preconditioner","PCPARMSSetLocal",PCPARMSLocalTypes,(PetscEnum)parms->local,(PetscEnum*)&local,&flag);
292:   if (flag) {PCPARMSSetLocal(pc,local);}
293:   PetscOptionsReal("-pc_parms_solve_tol","Tolerance for local solve","PCPARMSSetSolveTolerances",parms->solvetol,&parms->solvetol,NULL);
294:   PetscOptionsInt("-pc_parms_levels","Number of levels","None",parms->levels,&parms->levels,NULL);
295:   PetscOptionsBool("-pc_parms_nonsymmetric_perm","Use nonsymmetric permutation","PCPARMSSetNonsymPerm",parms->nonsymperm,&parms->nonsymperm,NULL);
296:   PetscOptionsInt("-pc_parms_blocksize","Block size","None",parms->blocksize,&parms->blocksize,NULL);
297:   PetscOptionsReal("-pc_parms_ind_tol","Tolerance for independent sets","None",parms->indtol,&parms->indtol,NULL);
298:   PetscOptionsInt("-pc_parms_max_dim","Inner Krylov dimension","PCPARMSSetSolveRestart",parms->maxdim,&parms->maxdim,NULL);
299:   PetscOptionsInt("-pc_parms_max_it","Maximum number of inner iterations","PCPARMSSetSolveTolerances",parms->maxits,&parms->maxits,NULL);
300:   PetscOptionsBool("-pc_parms_inter_nonsymmetric_perm","nonsymmetric permutation for interlevel blocks","None",parms->meth[0],&parms->meth[0],NULL);
301:   PetscOptionsBool("-pc_parms_inter_column_perm","column permutation for interlevel blocks","None",parms->meth[1],&parms->meth[1],NULL);
302:   PetscOptionsBool("-pc_parms_inter_row_scaling","row scaling for interlevel blocks","None",parms->meth[2],&parms->meth[2],NULL);
303:   PetscOptionsBool("-pc_parms_inter_column_scaling","column scaling for interlevel blocks","None",parms->meth[3],&parms->meth[3],NULL);
304:   PetscOptionsBool("-pc_parms_last_nonsymmetric_perm","nonsymmetric permutation for last level blocks","None",parms->meth[4],&parms->meth[4],NULL);
305:   PetscOptionsBool("-pc_parms_last_column_perm","column permutation for last level blocks","None",parms->meth[5],&parms->meth[5],NULL);
306:   PetscOptionsBool("-pc_parms_last_row_scaling","row scaling for last level blocks","None",parms->meth[6],&parms->meth[6],NULL);
307:   PetscOptionsBool("-pc_parms_last_column_scaling","column scaling for last level blocks","None",parms->meth[7],&parms->meth[7],NULL);
308:   PetscOptionsInt("-pc_parms_lfil_ilu_arms","amount of fill-in for ilut, iluk and arms","PCPARMSSetFill",parms->lfil[0],&parms->lfil[0],&flag);
309:   if (flag) parms->lfil[1] = parms->lfil[2] = parms->lfil[3] = parms->lfil[0];
310:   PetscOptionsInt("-pc_parms_lfil_schur","amount of fill-in for schur","PCPARMSSetFill",parms->lfil[4],&parms->lfil[4],NULL);
311:   PetscOptionsInt("-pc_parms_lfil_ilut_L_U","amount of fill-in for ILUT L and U","PCPARMSSetFill",parms->lfil[5],&parms->lfil[5],&flag);
312:   if (flag) parms->lfil[6] = parms->lfil[5];
313:   PetscOptionsReal("-pc_parms_droptol_factors","drop tolerance for L, U, L^{-1}F and EU^{-1}","None",parms->droptol[0],&parms->droptol[0],NULL);
314:   PetscOptionsReal("-pc_parms_droptol_schur_compl","drop tolerance for schur complement at each level","None",parms->droptol[4],&parms->droptol[4],&flag);
315:   if (flag) parms->droptol[1] = parms->droptol[2] = parms->droptol[3] = parms->droptol[0];
316:   PetscOptionsReal("-pc_parms_droptol_last_schur","drop tolerance for ILUT in last level schur complement","None",parms->droptol[5],&parms->droptol[5],&flag);
317:   if (flag) parms->droptol[6] = parms->droptol[5];
318:   PetscOptionsTail();
319:   return(0);
320: }

324: static PetscErrorCode PCApply_PARMS(PC pc,Vec b,Vec x)
325: {
326:   PetscErrorCode    ierr;
327:   PC_PARMS          *parms = (PC_PARMS*)pc->data;
328:   const PetscScalar *b1;
329:   PetscScalar       *x1;

332:   VecGetArrayRead(b,&b1);
333:   VecGetArray(x,&x1);
334:   parms_VecPermAux((PetscScalar*)b1,parms->lvec0,parms->map);
335:   parms_PCApply(parms->pc,parms->lvec0,parms->lvec1);
336:   parms_VecInvPermAux(parms->lvec1,x1,parms->map);
337:   VecRestoreArrayRead(b,&b1);
338:   VecRestoreArray(x,&x1);
339:   return(0);
340: }

344: static PetscErrorCode PCPARMSSetGlobal_PARMS(PC pc,PCPARMSGlobalType type)
345: {
346:   PC_PARMS *parms = (PC_PARMS*)pc->data;

349:   if (type != parms->global) {
350:     parms->global   = type;
351:     pc->setupcalled = 0;
352:   }
353:   return(0);
354: }

358: /*@
359:    PCPARMSSetGlobal - Sets the global preconditioner to be used in PARMS.

361:    Collective on PC

363:    Input Parameters:
364: +  pc - the preconditioner context
365: -  type - the global preconditioner type, one of
366: .vb
367:      PC_PARMS_GLOBAL_RAS   - Restricted additive Schwarz
368:      PC_PARMS_GLOBAL_SCHUR - Schur complement
369:      PC_PARMS_GLOBAL_BJ    - Block Jacobi
370: .ve

372:    Options Database Keys:
373:    -pc_parms_global [ras,schur,bj] - Sets global preconditioner

375:    Level: intermediate

377:    Notes:
378:    See the pARMS function parms_PCSetType for more information.

380: .seealso: PCPARMS, PCPARMSSetLocal()
381: @*/
382: PetscErrorCode PCPARMSSetGlobal(PC pc,PCPARMSGlobalType type)
383: {

389:   PetscTryMethod(pc,"PCPARMSSetGlobal_C",(PC,PCPARMSGlobalType),(pc,type));
390:   return(0);
391: }

395: static PetscErrorCode PCPARMSSetLocal_PARMS(PC pc,PCPARMSLocalType type)
396: {
397:   PC_PARMS *parms = (PC_PARMS*)pc->data;

400:   if (type != parms->local) {
401:     parms->local    = type;
402:     pc->setupcalled = 0;
403:   }
404:   return(0);
405: }

409: /*@
410:    PCPARMSSetLocal - Sets the local preconditioner to be used in PARMS.

412:    Collective on PC

414:    Input Parameters:
415: +  pc - the preconditioner context
416: -  type - the local preconditioner type, one of
417: .vb
418:      PC_PARMS_LOCAL_ILU0   - ILU0 preconditioner
419:      PC_PARMS_LOCAL_ILUK   - ILU(k) preconditioner
420:      PC_PARMS_LOCAL_ILUT   - ILUT preconditioner
421:      PC_PARMS_LOCAL_ARMS   - ARMS preconditioner
422: .ve

424:    Options Database Keys:
425:    -pc_parms_local [ilu0,iluk,ilut,arms] - Sets local preconditioner

427:    Level: intermediate

429:    Notes:
430:    For the ARMS preconditioner, one can use either the symmetric ARMS or the non-symmetric
431:    variant (ARMS-ddPQ) by setting the permutation type with PCPARMSSetNonsymPerm().

433:    See the pARMS function parms_PCILUSetType for more information.

435: .seealso: PCPARMS, PCPARMSSetGlobal(), PCPARMSSetNonsymPerm()

437: @*/
438: PetscErrorCode PCPARMSSetLocal(PC pc,PCPARMSLocalType type)
439: {

445:   PetscTryMethod(pc,"PCPARMSSetLocal_C",(PC,PCPARMSLocalType),(pc,type));
446:   return(0);
447: }

451: static PetscErrorCode PCPARMSSetSolveTolerances_PARMS(PC pc,PetscReal tol,PetscInt maxits)
452: {
453:   PC_PARMS *parms = (PC_PARMS*)pc->data;

456:   if (tol != parms->solvetol) {
457:     parms->solvetol = tol;
458:     pc->setupcalled = 0;
459:   }
460:   if (maxits != parms->maxits) {
461:     parms->maxits   = maxits;
462:     pc->setupcalled = 0;
463:   }
464:   return(0);
465: }

469: /*@
470:    PCPARMSSetSolveTolerances - Sets the convergence tolerance and the maximum iterations for the
471:    inner GMRES solver, when the Schur global preconditioner is used.

473:    Collective on PC

475:    Input Parameters:
476: +  pc - the preconditioner context
477: .  tol - the convergence tolerance
478: -  maxits - the maximum number of iterations to use

480:    Options Database Keys:
481: +  -pc_parms_solve_tol - set the tolerance for local solve
482: -  -pc_parms_max_it - set the maximum number of inner iterations

484:    Level: intermediate

486:    Notes:
487:    See the pARMS functions parms_PCSetInnerEps and parms_PCSetInnerMaxits for more information.

489: .seealso: PCPARMS, PCPARMSSetSolveRestart()
490: @*/
491: PetscErrorCode PCPARMSSetSolveTolerances(PC pc,PetscReal tol,PetscInt maxits)
492: {

497:   PetscTryMethod(pc,"PCPARMSSetSolveTolerances_C",(PC,PetscReal,PetscInt),(pc,tol,maxits));
498:   return(0);
499: }

503: static PetscErrorCode PCPARMSSetSolveRestart_PARMS(PC pc,PetscInt restart)
504: {
505:   PC_PARMS *parms = (PC_PARMS*)pc->data;

508:   if (restart != parms->maxdim) {
509:     parms->maxdim   = restart;
510:     pc->setupcalled = 0;
511:   }
512:   return(0);
513: }

517: /*@
518:    PCPARMSSetSolveRestart - Sets the number of iterations at which the
519:    inner GMRES solver restarts.

521:    Collective on PC

523:    Input Parameters:
524: +  pc - the preconditioner context
525: -  restart - maximum dimension of the Krylov subspace

527:    Options Database Keys:
528: .  -pc_parms_max_dim - sets the inner Krylov dimension

530:    Level: intermediate

532:    Notes:
533:    See the pARMS function parms_PCSetInnerKSize for more information.

535: .seealso: PCPARMS, PCPARMSSetSolveTolerances()
536: @*/
537: PetscErrorCode PCPARMSSetSolveRestart(PC pc,PetscInt restart)
538: {

543:   PetscTryMethod(pc,"PCPARMSSetSolveRestart_C",(PC,PetscInt),(pc,restart));
544:   return(0);
545: }

549: static PetscErrorCode PCPARMSSetNonsymPerm_PARMS(PC pc,PetscBool nonsym)
550: {
551:   PC_PARMS *parms = (PC_PARMS*)pc->data;

554:   if ((nonsym && !parms->nonsymperm) || (!nonsym && parms->nonsymperm)) {
555:     parms->nonsymperm = nonsym;
556:     pc->setupcalled   = 0;
557:   }
558:   return(0);
559: }

563: /*@
564:    PCPARMSSetNonsymPerm - Sets the type of permutation for the ARMS preconditioner: the standard
565:    symmetric ARMS or the non-symmetric ARMS (ARMS-ddPQ).

567:    Collective on PC

569:    Input Parameters:
570: +  pc - the preconditioner context
571: -  nonsym - PETSC_TRUE indicates the non-symmetric ARMS is used;
572:             PETSC_FALSE indicates the symmetric ARMS is used

574:    Options Database Keys:
575: .  -pc_parms_nonsymmetric_perm - sets the use of nonsymmetric permutation

577:    Level: intermediate

579:    Notes:
580:    See the pARMS function parms_PCSetPermType for more information.

582: .seealso: PCPARMS
583: @*/
584: PetscErrorCode PCPARMSSetNonsymPerm(PC pc,PetscBool nonsym)
585: {

590:   PetscTryMethod(pc,"PCPARMSSetNonsymPerm_C",(PC,PetscBool),(pc,nonsym));
591:   return(0);
592: }

596: static PetscErrorCode PCPARMSSetFill_PARMS(PC pc,PetscInt lfil0,PetscInt lfil1,PetscInt lfil2)
597: {
598:   PC_PARMS *parms = (PC_PARMS*)pc->data;

601:   if (lfil0 != parms->lfil[0] || lfil0 != parms->lfil[1] || lfil0 != parms->lfil[2] || lfil0 != parms->lfil[3]) {
602:     parms->lfil[1]  = parms->lfil[2] = parms->lfil[3] = parms->lfil[0] = lfil0;
603:     pc->setupcalled = 0;
604:   }
605:   if (lfil1 != parms->lfil[4]) {
606:     parms->lfil[4]  = lfil1;
607:     pc->setupcalled = 0;
608:   }
609:   if (lfil2 != parms->lfil[5] || lfil2 != parms->lfil[6]) {
610:     parms->lfil[5]  = parms->lfil[6] = lfil2;
611:     pc->setupcalled = 0;
612:   }
613:   return(0);
614: }

618: /*@
619:    PCPARMSSetFill - Sets the fill-in parameters for ILUT, ILUK and ARMS preconditioners.
620:    Consider the original matrix A = [B F; E C] and the approximate version
621:    M = [LB 0; E/UB I]*[UB LB\F; 0 S].

623:    Collective on PC

625:    Input Parameters:
626: +  pc - the preconditioner context
627: .  fil0 - the level of fill-in kept in LB, UB, E/UB and LB\F
628: .  fil1 - the level of fill-in kept in S
629: -  fil2 - the level of fill-in kept in the L and U parts of the LU factorization of S

631:    Options Database Keys:
632: +  -pc_parms_lfil_ilu_arms - set the amount of fill-in for ilut, iluk and arms
633: .  -pc_parms_lfil_schur - set the amount of fill-in for schur
634: -  -pc_parms_lfil_ilut_L_U - set the amount of fill-in for ILUT L and U

636:    Level: intermediate

638:    Notes:
639:    See the pARMS function parms_PCSetFill for more information.

641: .seealso: PCPARMS
642: @*/
643: PetscErrorCode PCPARMSSetFill(PC pc,PetscInt lfil0,PetscInt lfil1,PetscInt lfil2)
644: {

649:   PetscTryMethod(pc,"PCPARMSSetFill_C",(PC,PetscInt,PetscInt,PetscInt),(pc,lfil0,lfil1,lfil2));
650:   return(0);
651: }

653: /*MC
654:    PCPARMS - Allows the use of the parallel Algebraic Recursive Multilevel Solvers
655:       available in the package pARMS

657:    Options Database Keys:
658: +  -pc_parms_global - one of ras, schur, bj
659: .  -pc_parms_local - one of ilu0, iluk, ilut, arms
660: .  -pc_parms_solve_tol - set the tolerance for local solve
661: .  -pc_parms_levels - set the number of levels
662: .  -pc_parms_nonsymmetric_perm - set the use of nonsymmetric permutation
663: .  -pc_parms_blocksize - set the block size
664: .  -pc_parms_ind_tol - set the tolerance for independent sets
665: .  -pc_parms_max_dim - set the inner krylov dimension
666: .  -pc_parms_max_it - set the maximum number of inner iterations
667: .  -pc_parms_inter_nonsymmetric_perm - set the use of nonsymmetric permutation for interlevel blocks
668: .  -pc_parms_inter_column_perm - set the use of column permutation for interlevel blocks
669: .  -pc_parms_inter_row_scaling - set the use of row scaling for interlevel blocks
670: .  -pc_parms_inter_column_scaling - set the use of column scaling for interlevel blocks
671: .  -pc_parms_last_nonsymmetric_perm - set the use of nonsymmetric permutation for last level blocks
672: .  -pc_parms_last_column_perm - set the use of column permutation for last level blocks
673: .  -pc_parms_last_row_scaling - set the use of row scaling for last level blocks
674: .  -pc_parms_last_column_scaling - set the use of column scaling for last level blocks
675: .  -pc_parms_lfil_ilu_arms - set the amount of fill-in for ilut, iluk and arms
676: .  -pc_parms_lfil_schur - set the amount of fill-in for schur
677: .  -pc_parms_lfil_ilut_L_U - set the amount of fill-in for ILUT L and U
678: .  -pc_parms_droptol_factors - set the drop tolerance for L, U, L^{-1}F and EU^{-1}
679: .  -pc_parms_droptol_schur_compl - set the drop tolerance for schur complement at each level
680: -  -pc_parms_droptol_last_schur - set the drop tolerance for ILUT in last level schur complement

682:    IMPORTANT:
683:    Unless configured appropriately, this preconditioner performs an inexact solve
684:    as part of the preconditioner application. Therefore, it must be used in combination
685:    with flexible variants of iterative solvers, such as KSPFGMRES or KSPCGR.

687:    Level: intermediate

689: .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC
690: M*/

694: PETSC_EXTERN PetscErrorCode PCCreate_PARMS(PC pc)
695: {
696:   PC_PARMS       *parms;

700:   PetscNewLog(pc,&parms);

702:   parms->map        = 0;
703:   parms->A          = 0;
704:   parms->pc         = 0;
705:   parms->global     = PC_PARMS_GLOBAL_RAS;
706:   parms->local      = PC_PARMS_LOCAL_ARMS;
707:   parms->levels     = 10;
708:   parms->nonsymperm = PETSC_TRUE;
709:   parms->blocksize  = 250;
710:   parms->maxdim     = 0;
711:   parms->maxits     = 0;
712:   parms->meth[0]    = PETSC_FALSE;
713:   parms->meth[1]    = PETSC_FALSE;
714:   parms->meth[2]    = PETSC_FALSE;
715:   parms->meth[3]    = PETSC_FALSE;
716:   parms->meth[4]    = PETSC_FALSE;
717:   parms->meth[5]    = PETSC_FALSE;
718:   parms->meth[6]    = PETSC_FALSE;
719:   parms->meth[7]    = PETSC_FALSE;
720:   parms->solvetol   = 0.01;
721:   parms->indtol     = 0.4;
722:   parms->lfil[0]    = parms->lfil[1] = parms->lfil[2] = parms->lfil[3] = 20;
723:   parms->lfil[4]    = parms->lfil[5] = parms->lfil[6] = 20;
724:   parms->droptol[0] = parms->droptol[1] = parms->droptol[2] = parms->droptol[3] = 0.00001;
725:   parms->droptol[4] = 0.001;
726:   parms->droptol[5] = parms->droptol[6] = 0.001;

728:   pc->data                = parms;
729:   pc->ops->destroy        = PCDestroy_PARMS;
730:   pc->ops->setfromoptions = PCSetFromOptions_PARMS;
731:   pc->ops->setup          = PCSetUp_PARMS;
732:   pc->ops->apply          = PCApply_PARMS;
733:   pc->ops->view           = PCView_PARMS;

735:   PetscObjectComposeFunction((PetscObject)pc,"PCPARMSSetGlobal_C",PCPARMSSetGlobal_PARMS);
736:   PetscObjectComposeFunction((PetscObject)pc,"PCPARMSSetLocal_C",PCPARMSSetLocal_PARMS);
737:   PetscObjectComposeFunction((PetscObject)pc,"PCPARMSSetSolveTolerances_C",PCPARMSSetSolveTolerances_PARMS);
738:   PetscObjectComposeFunction((PetscObject)pc,"PCPARMSSetSolveRestart_C",PCPARMSSetSolveRestart_PARMS);
739:   PetscObjectComposeFunction((PetscObject)pc,"PCPARMSSetNonsymPerm_C",PCPARMSSetNonsymPerm_PARMS);
740:   PetscObjectComposeFunction((PetscObject)pc,"PCPARMSSetFill_C",PCPARMSSetFill_PARMS);
741:   return(0);
742: }