Actual source code: nasm.c

petsc-3.6.4 2016-04-12
Report Typos and Errors
  1: #include <petsc/private/snesimpl.h>             /*I   "petscsnes.h"   I*/
  2: #include <petscdm.h>

  4: typedef struct {
  5:   PetscInt   n;                   /* local subdomains */
  6:   SNES       *subsnes;            /* nonlinear solvers for each subdomain */
  7:   Vec        *x;                  /* solution vectors */
  8:   Vec        *xl;                 /* solution local vectors */
  9:   Vec        *y;                  /* step vectors */
 10:   Vec        *b;                  /* rhs vectors */
 11:   VecScatter *oscatter;           /* scatter from global space to the subdomain global space */
 12:   VecScatter *iscatter;           /* scatter from global space to the nonoverlapping subdomain space */
 13:   VecScatter *gscatter;           /* scatter from global space to the subdomain local space */
 14:   PCASMType  type;                /* ASM type */
 15:   PetscBool  usesdm;              /* use the DM for setting up the subproblems */
 16:   PetscBool  finaljacobian;       /* compute the jacobian of the converged solution */
 17:   PetscReal  damping;             /* damping parameter for updates from the blocks */
 18:   PetscBool  same_local_solves;   /* flag to determine if the solvers have been individually modified */

 20:   /* logging events */
 21:   PetscLogEvent eventrestrictinterp;
 22:   PetscLogEvent eventsubsolve;

 24:   PetscInt      fjtype;            /* type of computed jacobian */
 25:   Vec           xinit;             /* initial solution in case the final jacobian type is computed as first */
 26: } SNES_NASM;

 28: const char *const SNESNASMTypes[] = {"NONE","RESTRICT","INTERPOLATE","BASIC","PCASMType","PC_ASM_",0};
 29: const char *const SNESNASMFJTypes[] = {"FINALOUTER","FINALINNER","INITIAL"};

 33: PetscErrorCode SNESReset_NASM(SNES snes)
 34: {
 35:   SNES_NASM      *nasm = (SNES_NASM*)snes->data;
 37:   PetscInt       i;

 40:   for (i=0; i<nasm->n; i++) {
 41:     if (nasm->xl) { VecDestroy(&nasm->xl[i]); }
 42:     if (nasm->x) { VecDestroy(&nasm->x[i]); }
 43:     if (nasm->y) { VecDestroy(&nasm->y[i]); }
 44:     if (nasm->b) { VecDestroy(&nasm->b[i]); }

 46:     if (nasm->subsnes) { SNESDestroy(&nasm->subsnes[i]); }
 47:     if (nasm->oscatter) { VecScatterDestroy(&nasm->oscatter[i]); }
 48:     if (nasm->iscatter) { VecScatterDestroy(&nasm->iscatter[i]); }
 49:     if (nasm->gscatter) { VecScatterDestroy(&nasm->gscatter[i]); }
 50:   }

 52:   if (nasm->x) {PetscFree(nasm->x);}
 53:   if (nasm->xl) {PetscFree(nasm->xl);}
 54:   if (nasm->y) {PetscFree(nasm->y);}
 55:   if (nasm->b) {PetscFree(nasm->b);}

 57:   if (nasm->xinit) {VecDestroy(&nasm->xinit);}

 59:   if (nasm->subsnes) {PetscFree(nasm->subsnes);}
 60:   if (nasm->oscatter) {PetscFree(nasm->oscatter);}
 61:   if (nasm->iscatter) {PetscFree(nasm->iscatter);}
 62:   if (nasm->gscatter) {PetscFree(nasm->gscatter);}

 64:   nasm->eventrestrictinterp = 0;
 65:   nasm->eventsubsolve = 0;
 66:   return(0);
 67: }

 71: PetscErrorCode SNESDestroy_NASM(SNES snes)
 72: {

 76:   SNESReset_NASM(snes);
 77:   PetscFree(snes->data);
 78:   return(0);
 79: }

 83: PetscErrorCode DMGlobalToLocalSubDomainDirichletHook_Private(DM dm,Vec g,InsertMode mode,Vec l,void *ctx)
 84: {
 86:   Vec            bcs = (Vec)ctx;

 89:   VecCopy(bcs,l);
 90:   return(0);
 91: }

 95: PetscErrorCode SNESSetUp_NASM(SNES snes)
 96: {
 97:   SNES_NASM      *nasm = (SNES_NASM*)snes->data;
 99:   DM             dm,subdm;
100:   DM             *subdms;
101:   PetscInt       i;
102:   const char     *optionsprefix;
103:   Vec            F;
104:   PetscMPIInt    size;
105:   KSP            ksp;
106:   PC             pc;

109:   if (!nasm->subsnes) {
110:     SNESGetDM(snes,&dm);
111:     if (dm) {
112:       nasm->usesdm = PETSC_TRUE;
113:       DMCreateDomainDecomposition(dm,&nasm->n,NULL,NULL,NULL,&subdms);
114:       if (!subdms) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"DM has no default decomposition defined.  Set subsolves manually with SNESNASMSetSubdomains().");
115:       DMCreateDomainDecompositionScatters(dm,nasm->n,subdms,&nasm->iscatter,&nasm->oscatter,&nasm->gscatter);

117:       SNESGetOptionsPrefix(snes, &optionsprefix);
118:       PetscMalloc1(nasm->n,&nasm->subsnes);
119:       for (i=0; i<nasm->n; i++) {
120:         SNESCreate(PETSC_COMM_SELF,&nasm->subsnes[i]);
121:         SNESAppendOptionsPrefix(nasm->subsnes[i],optionsprefix);
122:         SNESAppendOptionsPrefix(nasm->subsnes[i],"sub_");
123:         SNESSetDM(nasm->subsnes[i],subdms[i]);
124:         MPI_Comm_size(PetscObjectComm((PetscObject)nasm->subsnes[i]),&size);
125:         if (size == 1) {
126:           SNESGetKSP(nasm->subsnes[i],&ksp);
127:           KSPGetPC(ksp,&pc);
128:           KSPSetType(ksp,KSPPREONLY);
129:           PCSetType(pc,PCLU);
130:         }
131:         SNESSetFromOptions(nasm->subsnes[i]);
132:         DMDestroy(&subdms[i]);
133:       }
134:       PetscFree(subdms);
135:     } else SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"Cannot construct local problems automatically without a DM!");
136:   } else SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"Must set subproblems manually if there is no DM!");
137:   /* allocate the global vectors */
138:   if (!nasm->x) {
139:     PetscCalloc1(nasm->n,&nasm->x);
140:   }
141:   if (!nasm->xl) {
142:     PetscCalloc1(nasm->n,&nasm->xl);
143:   }
144:   if (!nasm->y) {
145:     PetscCalloc1(nasm->n,&nasm->y);
146:   }
147:   if (!nasm->b) {
148:     PetscCalloc1(nasm->n,&nasm->b);
149:   }

151:   for (i=0; i<nasm->n; i++) {
152:     SNESGetFunction(nasm->subsnes[i],&F,NULL,NULL);
153:     if (!nasm->x[i]) {VecDuplicate(F,&nasm->x[i]);}
154:     if (!nasm->y[i]) {VecDuplicate(F,&nasm->y[i]);}
155:     if (!nasm->b[i]) {VecDuplicate(F,&nasm->b[i]);}
156:     if (!nasm->xl[i]) {
157:       SNESGetDM(nasm->subsnes[i],&subdm);
158:       DMCreateLocalVector(subdm,&nasm->xl[i]);
159:     }
160:     DMGlobalToLocalHookAdd(subdm,DMGlobalToLocalSubDomainDirichletHook_Private,NULL,nasm->xl[i]);
161:   }
162:   if (nasm->finaljacobian) {
163:     SNESSetUpMatrices(snes);
164:     if (nasm->fjtype == 2) {
165:       VecDuplicate(snes->vec_sol,&nasm->xinit);
166:     }
167:     for (i=0; i<nasm->n;i++) {
168:       SNESSetUpMatrices(nasm->subsnes[i]);
169:     }
170:   }
171:   return(0);
172: }

176: PetscErrorCode SNESSetFromOptions_NASM(PetscOptions *PetscOptionsObject,SNES snes)
177: {
178:   PetscErrorCode    ierr;
179:   PCASMType         asmtype;
180:   PetscBool         flg,monflg,subviewflg;
181:   SNES_NASM         *nasm = (SNES_NASM*)snes->data;

184:   PetscOptionsHead(PetscOptionsObject,"Nonlinear Additive Schwartz options");
185:   PetscOptionsEnum("-snes_nasm_type","Type of restriction/extension","",SNESNASMTypes,(PetscEnum)nasm->type,(PetscEnum*)&asmtype,&flg);
186:   if (flg) {SNESNASMSetType(snes,asmtype);}
187:   flg    = PETSC_FALSE;
188:   monflg = PETSC_TRUE;
189:   PetscOptionsReal("-snes_nasm_damping","Log times for subSNES solves and restriction","SNESNASMSetDamping",nasm->damping,&nasm->damping,&flg);
190:   if (flg) {SNESNASMSetDamping(snes,nasm->damping);}
191:   subviewflg = PETSC_FALSE;
192:   PetscOptionsBool("-snes_nasm_sub_view","Print detailed information for every processor when using -snes_view","",subviewflg,&subviewflg,&flg);
193:   if (flg) {
194:     nasm->same_local_solves = PETSC_FALSE;
195:     if (!subviewflg) {
196:       nasm->same_local_solves = PETSC_TRUE;
197:     }
198:   }
199:   PetscOptionsBool("-snes_nasm_finaljacobian","Compute the global jacobian of the final iterate (for ASPIN)","",nasm->finaljacobian,&nasm->finaljacobian,NULL);
200:   PetscOptionsEList("-snes_nasm_finaljacobian_type","The type of the final jacobian computed.","",SNESNASMFJTypes,3,SNESNASMFJTypes[0],&nasm->fjtype,NULL);
201:   PetscOptionsBool("-snes_nasm_log","Log times for subSNES solves and restriction","",monflg,&monflg,&flg);
202:   if (flg) {
203:     PetscLogEventRegister("SNESNASMSubSolve",((PetscObject)snes)->classid,&nasm->eventsubsolve);
204:     PetscLogEventRegister("SNESNASMRestrict",((PetscObject)snes)->classid,&nasm->eventrestrictinterp);
205:   }
206:   PetscOptionsTail();
207:   return(0);
208: }

212: PetscErrorCode SNESView_NASM(SNES snes, PetscViewer viewer)
213: {
214:   SNES_NASM      *nasm = (SNES_NASM*)snes->data;
216:   PetscMPIInt    rank,size;
217:   PetscInt       i,N,bsz;
218:   PetscBool      iascii,isstring;
219:   PetscViewer    sviewer;
220:   MPI_Comm       comm;

223:   PetscObjectGetComm((PetscObject)snes,&comm);
224:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
225:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring);
226:   MPI_Comm_rank(comm,&rank);
227:   MPI_Comm_size(comm,&size);
228:   MPI_Allreduce(&nasm->n,&N,1,MPIU_INT,MPI_SUM,comm);
229:   if (iascii) {
230:     PetscViewerASCIIPrintf(viewer, "  Nonlinear Additive Schwarz: total subdomain blocks = %D\n",N);
231:     if (nasm->same_local_solves) {
232:       if (nasm->subsnes) {
233:         PetscViewerASCIIPrintf(viewer,"  Local solve is the same for all blocks:\n");
234:         PetscViewerASCIIPushTab(viewer);
235:         PetscViewerGetSingleton(viewer,&sviewer);
236:         if (!rank) {
237:           PetscViewerASCIIPushTab(viewer);
238:           SNESView(nasm->subsnes[0],sviewer);
239:           PetscViewerASCIIPopTab(viewer);
240:         }
241:         PetscViewerRestoreSingleton(viewer,&sviewer);
242:         PetscViewerASCIIPopTab(viewer);
243:       }
244:     } else {
245:       /* print the solver on each block */
246:       PetscViewerASCIISynchronizedAllow(viewer,PETSC_TRUE);
247:       PetscViewerASCIISynchronizedPrintf(viewer,"  [%d] number of local blocks = %D\n",(int)rank,nasm->n);
248:       PetscViewerFlush(viewer);
249:       PetscViewerASCIISynchronizedAllow(viewer,PETSC_FALSE);
250:       PetscViewerASCIIPrintf(viewer,"  Local solve info for each block is in the following SNES objects:\n");
251:       PetscViewerASCIIPushTab(viewer);
252:       PetscViewerASCIIPrintf(viewer,"- - - - - - - - - - - - - - - - - -\n");
253:       PetscViewerGetSingleton(viewer,&sviewer);
254:       for (i=0; i<nasm->n; i++) {
255:         VecGetLocalSize(nasm->x[i],&bsz);
256:         PetscViewerASCIIPrintf(sviewer,"[%d] local block number %D, size = %D\n",(int)rank,i,bsz);
257:         SNESView(nasm->subsnes[i],sviewer);
258:         PetscViewerASCIIPrintf(sviewer,"- - - - - - - - - - - - - - - - - -\n");
259:       }
260:       PetscViewerRestoreSingleton(viewer,&sviewer);
261:       PetscViewerFlush(viewer);
262:       PetscViewerASCIIPopTab(viewer);
263:     }
264:   } else if (isstring) {
265:     PetscViewerStringSPrintf(viewer," blocks=%D,type=%s",N,SNESNASMTypes[nasm->type]);
266:     PetscViewerGetSingleton(viewer,&sviewer);
267:     if (nasm->subsnes && !rank) {SNESView(nasm->subsnes[0],sviewer);}
268:     PetscViewerRestoreSingleton(viewer,&sviewer);
269:   }
270:   return(0);
271: }

275: /*@
276:    SNESNASMSetType - Set the type of subdomain update used

278:    Logically Collective on SNES

280:    Input Parameters:
281: +  SNES - the SNES context
282: -  type - the type of update, PC_ASM_BASIC or PC_ASM_RESTRICT

284:    Level: intermediate

286: .keywords: SNES, NASM

288: .seealso: SNESNASM, SNESNASMGetType(), PCASMSetType()
289: @*/
290: PetscErrorCode SNESNASMSetType(SNES snes,PCASMType type)
291: {
293:   PetscErrorCode (*f)(SNES,PCASMType);

296:   PetscObjectQueryFunction((PetscObject)snes,"SNESNASMSetType_C",&f);
297:   if (f) {(f)(snes,type);}
298:   return(0);
299: }

303: PetscErrorCode SNESNASMSetType_NASM(SNES snes,PCASMType type)
304: {
305:   SNES_NASM      *nasm = (SNES_NASM*)snes->data;

308:   if (type != PC_ASM_BASIC && type != PC_ASM_RESTRICT) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_OUTOFRANGE,"SNESNASM only supports basic and restrict types");
309:   nasm->type = type;
310:   return(0);
311: }

315: /*@
316:    SNESNASMGetType - Get the type of subdomain update used

318:    Logically Collective on SNES

320:    Input Parameters:
321: .  SNES - the SNES context

323:    Output Parameters:
324: .  type - the type of update

326:    Level: intermediate

328: .keywords: SNES, NASM

330: .seealso: SNESNASM, SNESNASMSetType(), PCASMGetType()
331: @*/
332: PetscErrorCode SNESNASMGetType(SNES snes,PCASMType *type)
333: {
335:   PetscErrorCode (*f)(SNES,PCASMType*);

338:   PetscObjectQueryFunction((PetscObject)snes,"SNESNASMGetType_C",&f);
339:   if (f) {(f)(snes,type);}
340:   return(0);
341: }

345: PetscErrorCode SNESNASMGetType_NASM(SNES snes,PCASMType *type)
346: {
347:   SNES_NASM      *nasm = (SNES_NASM*)snes->data;

350:   *type = nasm->type;
351:   return(0);
352: }

356: /*@
357:    SNESNASMSetSubdomains - Manually Set the context required to restrict and solve subdomain problems.

359:    Not Collective

361:    Input Parameters:
362: +  SNES - the SNES context
363: .  n - the number of local subdomains
364: .  subsnes - solvers defined on the local subdomains
365: .  iscatter - scatters into the nonoverlapping portions of the local subdomains
366: .  oscatter - scatters into the overlapping portions of the local subdomains
367: -  gscatter - scatters into the (ghosted) local vector of the local subdomain

369:    Level: intermediate

371: .keywords: SNES, NASM

373: .seealso: SNESNASM, SNESNASMGetSubdomains()
374: @*/
375: PetscErrorCode SNESNASMSetSubdomains(SNES snes,PetscInt n,SNES subsnes[],VecScatter iscatter[],VecScatter oscatter[],VecScatter gscatter[])
376: {
378:   PetscErrorCode (*f)(SNES,PetscInt,SNES*,VecScatter*,VecScatter*,VecScatter*);

381:   PetscObjectQueryFunction((PetscObject)snes,"SNESNASMSetSubdomains_C",&f);
382:   if (f) {(f)(snes,n,subsnes,iscatter,oscatter,gscatter);}
383:   return(0);
384: }

388: PetscErrorCode SNESNASMSetSubdomains_NASM(SNES snes,PetscInt n,SNES subsnes[],VecScatter iscatter[],VecScatter oscatter[],VecScatter gscatter[])
389: {
390:   PetscInt       i;
392:   SNES_NASM      *nasm = (SNES_NASM*)snes->data;

395:   if (snes->setupcalled) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"SNESNASMSetSubdomains() should be called before calling SNESSetUp().");

397:   /* tear down the previously set things */
398:   SNESReset(snes);

400:   nasm->n = n;
401:   if (oscatter) {
402:     for (i=0; i<n; i++) {PetscObjectReference((PetscObject)oscatter[i]);}
403:   }
404:   if (iscatter) {
405:     for (i=0; i<n; i++) {PetscObjectReference((PetscObject)iscatter[i]);}
406:   }
407:   if (gscatter) {
408:     for (i=0; i<n; i++) {PetscObjectReference((PetscObject)gscatter[i]);}
409:   }
410:   if (oscatter) {
411:     PetscMalloc1(n,&nasm->oscatter);
412:     for (i=0; i<n; i++) {
413:       nasm->oscatter[i] = oscatter[i];
414:     }
415:   }
416:   if (iscatter) {
417:     PetscMalloc1(n,&nasm->iscatter);
418:     for (i=0; i<n; i++) {
419:       nasm->iscatter[i] = iscatter[i];
420:     }
421:   }
422:   if (gscatter) {
423:     PetscMalloc1(n,&nasm->gscatter);
424:     for (i=0; i<n; i++) {
425:       nasm->gscatter[i] = gscatter[i];
426:     }
427:   }

429:   if (subsnes) {
430:     PetscMalloc1(n,&nasm->subsnes);
431:     for (i=0; i<n; i++) {
432:       nasm->subsnes[i] = subsnes[i];
433:     }
434:     nasm->same_local_solves = PETSC_FALSE;
435:   }
436:   return(0);
437: }

441: /*@
442:    SNESNASMGetSubdomains - Get the local subdomain context.

444:    Not Collective

446:    Input Parameters:
447: .  SNES - the SNES context

449:    Output Parameters:
450: +  n - the number of local subdomains
451: .  subsnes - solvers defined on the local subdomains
452: .  iscatter - scatters into the nonoverlapping portions of the local subdomains
453: .  oscatter - scatters into the overlapping portions of the local subdomains
454: -  gscatter - scatters into the (ghosted) local vector of the local subdomain

456:    Level: intermediate

458: .keywords: SNES, NASM

460: .seealso: SNESNASM, SNESNASMSetSubdomains()
461: @*/
462: PetscErrorCode SNESNASMGetSubdomains(SNES snes,PetscInt *n,SNES *subsnes[],VecScatter *iscatter[],VecScatter *oscatter[],VecScatter *gscatter[])
463: {
465:   PetscErrorCode (*f)(SNES,PetscInt*,SNES**,VecScatter**,VecScatter**,VecScatter**);

468:   PetscObjectQueryFunction((PetscObject)snes,"SNESNASMGetSubdomains_C",&f);
469:   if (f) {(f)(snes,n,subsnes,iscatter,oscatter,gscatter);}
470:   return(0);
471: }

475: PetscErrorCode SNESNASMGetSubdomains_NASM(SNES snes,PetscInt *n,SNES *subsnes[],VecScatter *iscatter[],VecScatter *oscatter[],VecScatter *gscatter[])
476: {
477:   SNES_NASM      *nasm = (SNES_NASM*)snes->data;

480:   if (n) *n = nasm->n;
481:   if (oscatter) *oscatter = nasm->oscatter;
482:   if (iscatter) *iscatter = nasm->iscatter;
483:   if (gscatter) *gscatter = nasm->gscatter;
484:   if (subsnes)  {
485:     *subsnes  = nasm->subsnes;
486:     nasm->same_local_solves = PETSC_FALSE;
487:   }
488:   return(0);
489: }

493: /*@
494:    SNESNASMGetSubdomainVecs - Get the processor-local subdomain vectors

496:    Not Collective

498:    Input Parameters:
499: .  SNES - the SNES context

501:    Output Parameters:
502: +  n - the number of local subdomains
503: .  x - The subdomain solution vector
504: .  y - The subdomain step vector
505: .  b - The subdomain RHS vector
506: -  xl - The subdomain local vectors (ghosted)

508:    Level: developer

510: .keywords: SNES, NASM

512: .seealso: SNESNASM, SNESNASMGetSubdomains()
513: @*/
514: PetscErrorCode SNESNASMGetSubdomainVecs(SNES snes,PetscInt *n,Vec **x,Vec **y,Vec **b, Vec **xl)
515: {
517:   PetscErrorCode (*f)(SNES,PetscInt*,Vec**,Vec**,Vec**,Vec**);

520:   PetscObjectQueryFunction((PetscObject)snes,"SNESNASMGetSubdomainVecs_C",&f);
521:   if (f) {(f)(snes,n,x,y,b,xl);}
522:   return(0);
523: }

527: PetscErrorCode SNESNASMGetSubdomainVecs_NASM(SNES snes,PetscInt *n,Vec **x,Vec **y,Vec **b,Vec **xl)
528: {
529:   SNES_NASM      *nasm = (SNES_NASM*)snes->data;

532:   if (n)  *n  = nasm->n;
533:   if (x)  *x  = nasm->x;
534:   if (y)  *y  = nasm->y;
535:   if (b)  *b  = nasm->b;
536:   if (xl) *xl = nasm->xl;
537:   return(0);
538: }

542: /*@
543:    SNESNASMSetComputeFinalJacobian - Schedules the computation of the global and subdomain jacobians upon convergence

545:    Collective on SNES

547:    Input Parameters:
548: +  SNES - the SNES context
549: -  flg - indication of whether to compute the jacobians or not

551:    Level: developer

553:    Notes: This is used almost exclusively in the implementation of ASPIN, where the converged subdomain and global jacobian
554:    is needed at each linear iteration.

556: .keywords: SNES, NASM, ASPIN

558: .seealso: SNESNASM, SNESNASMGetSubdomains()
559: @*/
560: PetscErrorCode SNESNASMSetComputeFinalJacobian(SNES snes,PetscBool flg)
561: {
562:   PetscErrorCode (*f)(SNES,PetscBool);

566:   PetscObjectQueryFunction((PetscObject)snes,"SNESNASMSetComputeFinalJacobian_C",&f);
567:   if (f) {(f)(snes,flg);}
568:   return(0);
569: }

573: PetscErrorCode SNESNASMSetComputeFinalJacobian_NASM(SNES snes,PetscBool flg)
574: {
575:   SNES_NASM      *nasm = (SNES_NASM*)snes->data;

578:   nasm->finaljacobian = flg;
579:   if (flg) snes->usesksp = PETSC_TRUE;
580:   return(0);
581: }

585: /*@
586:    SNESNASMSetDamping - Sets the update damping for NASM

588:    Logically collective on SNES

590:    Input Parameters:
591: +  SNES - the SNES context
592: -  dmp - damping

594:    Level: intermediate

596: .keywords: SNES, NASM, damping

598: .seealso: SNESNASM, SNESNASMGetDamping()
599: @*/
600: PetscErrorCode SNESNASMSetDamping(SNES snes,PetscReal dmp)
601: {
602:   PetscErrorCode (*f)(SNES,PetscReal);

606:   PetscObjectQueryFunction((PetscObject)snes,"SNESNASMSetDamping_C",(void (**)(void))&f);
607:   if (f) {(f)(snes,dmp);}
608:   return(0);
609: }

613: PetscErrorCode SNESNASMSetDamping_NASM(SNES snes,PetscReal dmp)
614: {
615:   SNES_NASM      *nasm = (SNES_NASM*)snes->data;

618:   nasm->damping = dmp;
619:   return(0);
620: }

624: /*@
625:    SNESNASMGetDamping - Gets the update damping for NASM

627:    Not Collective

629:    Input Parameters:
630: +  SNES - the SNES context
631: -  dmp - damping

633:    Level: intermediate

635: .keywords: SNES, NASM, damping

637: .seealso: SNESNASM, SNESNASMSetDamping()
638: @*/
639: PetscErrorCode SNESNASMGetDamping(SNES snes,PetscReal *dmp)
640: {
641:   PetscErrorCode (*f)(SNES,PetscReal*);

645:   PetscObjectQueryFunction((PetscObject)snes,"SNESNASMGetDamping_C",(void (**)(void))&f);
646:   if (f) {(f)(snes,dmp);}
647:   return(0);
648: }

652: PetscErrorCode SNESNASMGetDamping_NASM(SNES snes,PetscReal *dmp)
653: {
654:   SNES_NASM      *nasm = (SNES_NASM*)snes->data;

657:   *dmp = nasm->damping;
658:   return(0);
659: }


664: /*
665:   Input Parameters:
666: + snes - The solver
667: . B - The RHS vector
668: - X - The initial guess

670:   Output Parameters:
671: . Y - The solution update

673:   TODO: All scatters should be packed into one
674: */
675: PetscErrorCode SNESNASMSolveLocal_Private(SNES snes,Vec B,Vec Y,Vec X)
676: {
677:   SNES_NASM      *nasm = (SNES_NASM*)snes->data;
678:   SNES           subsnes;
679:   PetscInt       i;
680:   PetscReal      dmp;
682:   Vec            Xlloc,Xl,Bl,Yl;
683:   VecScatter     iscat,oscat,gscat;
684:   DM             dm,subdm;
685:   PCASMType      type;

688:   SNESNASMGetType(snes,&type);
689:   SNESGetDM(snes,&dm);
690:   SNESNASMGetDamping(snes,&dmp);
691:   VecSet(Y,0);
692:   if (nasm->eventrestrictinterp) {PetscLogEventBegin(nasm->eventrestrictinterp,snes,0,0,0);}
693:   for (i=0; i<nasm->n; i++) {
694:     /* scatter the solution to the local solution */
695:     Xlloc = nasm->xl[i];
696:     gscat   = nasm->gscatter[i];
697:     oscat   = nasm->oscatter[i];
698:     VecScatterBegin(gscat,X,Xlloc,INSERT_VALUES,SCATTER_FORWARD);
699:     if (B) {
700:       /* scatter the RHS to the local RHS */
701:       Bl   = nasm->b[i];
702:       VecScatterBegin(oscat,B,Bl,INSERT_VALUES,SCATTER_FORWARD);
703:     }
704:   }
705:   if (nasm->eventrestrictinterp) {PetscLogEventEnd(nasm->eventrestrictinterp,snes,0,0,0);}


708:   if (nasm->eventsubsolve) {PetscLogEventBegin(nasm->eventsubsolve,snes,0,0,0);}
709:   for (i=0; i<nasm->n; i++) {
710:     Xl    = nasm->x[i];
711:     Xlloc = nasm->xl[i];
712:     Yl    = nasm->y[i];
713:     subsnes = nasm->subsnes[i];
714:     SNESGetDM(subsnes,&subdm);
715:     iscat   = nasm->iscatter[i];
716:     oscat   = nasm->oscatter[i];
717:     gscat   = nasm->gscatter[i];
718:     VecScatterEnd(gscat,X,Xlloc,INSERT_VALUES,SCATTER_FORWARD);
719:     if (B) {
720:       Bl   = nasm->b[i];
721:       VecScatterEnd(oscat,B,Bl,INSERT_VALUES,SCATTER_FORWARD);
722:     } else Bl = NULL;
723:     DMSubDomainRestrict(dm,oscat,gscat,subdm);
724:     /* Could scatter directly from X */
725:     DMLocalToGlobalBegin(subdm,Xlloc,INSERT_VALUES,Xl);
726:     DMLocalToGlobalEnd(subdm,Xlloc,INSERT_VALUES,Xl);
727:     VecCopy(Xl,Yl);
728:     SNESSolve(subsnes,Bl,Xl);
729:     VecAYPX(Yl,-1.0,Xl);
730:     if (type == PC_ASM_BASIC) {
731:       VecScatterBegin(oscat,Yl,Y,ADD_VALUES,SCATTER_REVERSE);
732:     } else if (type == PC_ASM_RESTRICT) {
733:       VecScatterBegin(iscat,Yl,Y,ADD_VALUES,SCATTER_REVERSE);
734:     } else SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"Only basic and restrict types are supported for SNESNASM");
735:   }
736:   if (nasm->eventsubsolve) {PetscLogEventEnd(nasm->eventsubsolve,snes,0,0,0);}
737:   if (nasm->eventrestrictinterp) {PetscLogEventBegin(nasm->eventrestrictinterp,snes,0,0,0);}
738:   for (i=0; i<nasm->n; i++) {
739:     Yl    = nasm->y[i];
740:     iscat   = nasm->iscatter[i];
741:     oscat   = nasm->oscatter[i];
742:     if (type == PC_ASM_BASIC) {
743:       VecScatterEnd(oscat,Yl,Y,ADD_VALUES,SCATTER_REVERSE);
744:     } else if (type == PC_ASM_RESTRICT) {
745:       VecScatterEnd(iscat,Yl,Y,ADD_VALUES,SCATTER_REVERSE);
746:     } else SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"Only basic and restrict types are supported for SNESNASM");
747:   }
748:   if (nasm->eventrestrictinterp) {PetscLogEventEnd(nasm->eventrestrictinterp,snes,0,0,0);}
749:   VecAXPY(X,dmp,Y);
750:   return(0);
751: }

755: PetscErrorCode SNESNASMComputeFinalJacobian_Private(SNES snes, Vec Xfinal)
756: {
757:   Vec            X = Xfinal;
758:   SNES_NASM      *nasm = (SNES_NASM*)snes->data;
759:   SNES           subsnes;
760:   PetscInt       i,lag = 1;
762:   Vec            Xlloc,Xl,Fl,F;
763:   VecScatter     oscat,gscat;
764:   DM             dm,subdm;

767:   if (nasm->fjtype == 2) X = nasm->xinit;
768:   F = snes->vec_func;
769:   if (snes->normschedule == SNES_NORM_NONE) {SNESComputeFunction(snes,X,F);}
770:   SNESComputeJacobian(snes,X,snes->jacobian,snes->jacobian_pre);
771:   SNESGetDM(snes,&dm);
772:   if (nasm->eventrestrictinterp) {PetscLogEventBegin(nasm->eventrestrictinterp,snes,0,0,0);}
773:   if (nasm->fjtype != 1) {
774:     for (i=0; i<nasm->n; i++) {
775:       Xlloc = nasm->xl[i];
776:       gscat = nasm->gscatter[i];
777:       oscat = nasm->oscatter[i];
778:       VecScatterBegin(gscat,X,Xlloc,INSERT_VALUES,SCATTER_FORWARD);
779:     }
780:   }
781:   if (nasm->eventrestrictinterp) {PetscLogEventEnd(nasm->eventrestrictinterp,snes,0,0,0);}
782:   for (i=0; i<nasm->n; i++) {
783:     Fl      = nasm->subsnes[i]->vec_func;
784:     Xl      = nasm->x[i];
785:     Xlloc   = nasm->xl[i];
786:     subsnes = nasm->subsnes[i];
787:     oscat   = nasm->oscatter[i];
788:     gscat   = nasm->gscatter[i];
789:     if (nasm->fjtype != 1) {VecScatterEnd(gscat,X,Xlloc,INSERT_VALUES,SCATTER_FORWARD);}
790:     SNESGetDM(subsnes,&subdm);
791:     DMSubDomainRestrict(dm,oscat,gscat,subdm);
792:     if (nasm->fjtype != 1) {
793:       DMLocalToGlobalBegin(subdm,Xlloc,INSERT_VALUES,Xl);
794:       DMLocalToGlobalEnd(subdm,Xlloc,INSERT_VALUES,Xl);
795:     }
796:     if (subsnes->lagjacobian == -1)    subsnes->lagjacobian = -2;
797:     else if (subsnes->lagjacobian > 1) lag = subsnes->lagjacobian;
798:     SNESComputeFunction(subsnes,Xl,Fl);
799:     SNESComputeJacobian(subsnes,Xl,subsnes->jacobian,subsnes->jacobian_pre);
800:     if (lag > 1) subsnes->lagjacobian = lag;
801:   }
802:   return(0);
803: }

807: PetscErrorCode SNESSolve_NASM(SNES snes)
808: {
809:   Vec              F;
810:   Vec              X;
811:   Vec              B;
812:   Vec              Y;
813:   PetscInt         i;
814:   PetscReal        fnorm = 0.0;
815:   PetscErrorCode   ierr;
816:   SNESNormSchedule normschedule;
817:   SNES_NASM        *nasm = (SNES_NASM*)snes->data;


821:   if (snes->xl || snes->xu || snes->ops->computevariablebounds) {
822:     SETERRQ1(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE, "SNES solver %s does not support bounds", ((PetscObject)snes)->type_name);
823:   }

825:   PetscCitationsRegister(SNESCitation,&SNEScite);
826:   X = snes->vec_sol;
827:   Y = snes->vec_sol_update;
828:   F = snes->vec_func;
829:   B = snes->vec_rhs;

831:   PetscObjectSAWsTakeAccess((PetscObject)snes);
832:   snes->iter   = 0;
833:   snes->norm   = 0.;
834:   PetscObjectSAWsGrantAccess((PetscObject)snes);
835:   snes->reason = SNES_CONVERGED_ITERATING;
836:   SNESGetNormSchedule(snes, &normschedule);
837:   if (normschedule == SNES_NORM_ALWAYS || normschedule == SNES_NORM_INITIAL_ONLY || normschedule == SNES_NORM_INITIAL_FINAL_ONLY) {
838:     /* compute the initial function and preconditioned update delX */
839:     if (!snes->vec_func_init_set) {
840:       SNESComputeFunction(snes,X,F);
841:     } else snes->vec_func_init_set = PETSC_FALSE;

843:     VecNorm(F, NORM_2, &fnorm); /* fnorm <- ||F||  */
844:     SNESCheckFunctionNorm(snes,fnorm);
845:     PetscObjectSAWsTakeAccess((PetscObject)snes);
846:     snes->iter = 0;
847:     snes->norm = fnorm;
848:     PetscObjectSAWsGrantAccess((PetscObject)snes);
849:     SNESLogConvergenceHistory(snes,snes->norm,0);
850:     SNESMonitor(snes,0,snes->norm);

852:     /* test convergence */
853:     (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);
854:     if (snes->reason) return(0);
855:   } else {
856:     PetscObjectSAWsGrantAccess((PetscObject)snes);
857:     SNESLogConvergenceHistory(snes,snes->norm,0);
858:     SNESMonitor(snes,0,snes->norm);
859:   }

861:   /* Call general purpose update function */
862:   if (snes->ops->update) {
863:     (*snes->ops->update)(snes, snes->iter);
864:   }
865:   /* copy the initial solution over for later */
866:   if (nasm->fjtype == 2) {VecCopy(X,nasm->xinit);}

868:   for (i = 0; i < snes->max_its; i++) {
869:     SNESNASMSolveLocal_Private(snes,B,Y,X);
870:     if (normschedule == SNES_NORM_ALWAYS || ((i == snes->max_its - 1) && (normschedule == SNES_NORM_INITIAL_FINAL_ONLY || normschedule == SNES_NORM_FINAL_ONLY))) {
871:       SNESComputeFunction(snes,X,F);
872:       VecNorm(F, NORM_2, &fnorm); /* fnorm <- ||F||  */
873:       SNESCheckFunctionNorm(snes,fnorm);
874:     }
875:     /* Monitor convergence */
876:     PetscObjectSAWsTakeAccess((PetscObject)snes);
877:     snes->iter = i+1;
878:     snes->norm = fnorm;
879:     PetscObjectSAWsGrantAccess((PetscObject)snes);
880:     SNESLogConvergenceHistory(snes,snes->norm,0);
881:     SNESMonitor(snes,snes->iter,snes->norm);
882:     /* Test for convergence */
883:     if (normschedule == SNES_NORM_ALWAYS) {(*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP);}
884:     if (snes->reason) break;
885:     /* Call general purpose update function */
886:     if (snes->ops->update) {(*snes->ops->update)(snes, snes->iter);}
887:   }
888:   if (nasm->finaljacobian) {SNESNASMComputeFinalJacobian_Private(snes,X);}
889:   if (normschedule == SNES_NORM_ALWAYS) {
890:     if (i == snes->max_its) {
891:       PetscInfo1(snes,"Maximum number of iterations has been reached: %D\n",snes->max_its);
892:       if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
893:     }
894:   } else if (!snes->reason) snes->reason = SNES_CONVERGED_ITS; /* NASM is meant to be used as a preconditioner */
895:   return(0);
896: }

898: /*MC
899:   SNESNASM - Nonlinear Additive Schwartz

901:    Options Database:
902: +  -snes_nasm_log - enable logging events for the communication and solve stages
903: .  -snes_nasm_type <basic,restrict> - type of subdomain update used
904: .  -snes_nasm_finaljacobian - compute the local and global jacobians of the final iterate
905: .  -snes_nasm_finaljacobian_type <finalinner,finalouter,initial> pick state the jacobian is calculated at
906: .  -sub_snes_ - options prefix of the subdomain nonlinear solves
907: .  -sub_ksp_ - options prefix of the subdomain Krylov solver
908: -  -sub_pc_ - options prefix of the subdomain preconditioner

910:    Level: advanced

912: .seealso: SNESCreate(), SNES, SNESSetType(), SNESType (for list of available types)
913: M*/

917: PETSC_EXTERN PetscErrorCode SNESCreate_NASM(SNES snes)
918: {
919:   SNES_NASM      *nasm;

923:   PetscNewLog(snes,&nasm);
924:   snes->data = (void*)nasm;

926:   nasm->n        = PETSC_DECIDE;
927:   nasm->subsnes  = 0;
928:   nasm->x        = 0;
929:   nasm->xl       = 0;
930:   nasm->y        = 0;
931:   nasm->b        = 0;
932:   nasm->oscatter = 0;
933:   nasm->iscatter = 0;
934:   nasm->gscatter = 0;
935:   nasm->damping  = 1.;

937:   nasm->type = PC_ASM_BASIC;
938:   nasm->finaljacobian = PETSC_FALSE;
939:   nasm->same_local_solves = PETSC_TRUE;

941:   snes->ops->destroy        = SNESDestroy_NASM;
942:   snes->ops->setup          = SNESSetUp_NASM;
943:   snes->ops->setfromoptions = SNESSetFromOptions_NASM;
944:   snes->ops->view           = SNESView_NASM;
945:   snes->ops->solve          = SNESSolve_NASM;
946:   snes->ops->reset          = SNESReset_NASM;

948:   snes->usesksp = PETSC_FALSE;
949:   snes->usespc  = PETSC_FALSE;

951:   nasm->fjtype              = 0;
952:   nasm->xinit               = NULL;
953:   nasm->eventrestrictinterp = 0;
954:   nasm->eventsubsolve       = 0;

956:   if (!snes->tolerancesset) {
957:     snes->max_its   = 10000;
958:     snes->max_funcs = 10000;
959:   }

961:   PetscObjectComposeFunction((PetscObject)snes,"SNESNASMSetType_C",SNESNASMSetType_NASM);
962:   PetscObjectComposeFunction((PetscObject)snes,"SNESNASMGetType_C",SNESNASMGetType_NASM);
963:   PetscObjectComposeFunction((PetscObject)snes,"SNESNASMSetSubdomains_C",SNESNASMSetSubdomains_NASM);
964:   PetscObjectComposeFunction((PetscObject)snes,"SNESNASMGetSubdomains_C",SNESNASMGetSubdomains_NASM);
965:   PetscObjectComposeFunction((PetscObject)snes,"SNESNASMSetDamping_C",SNESNASMSetDamping_NASM);
966:   PetscObjectComposeFunction((PetscObject)snes,"SNESNASMGetDamping_C",SNESNASMGetDamping_NASM);
967:   PetscObjectComposeFunction((PetscObject)snes,"SNESNASMGetSubdomainVecs_C",SNESNASMGetSubdomainVecs_NASM);
968:   PetscObjectComposeFunction((PetscObject)snes,"SNESNASMSetComputeFinalJacobian_C",SNESNASMSetComputeFinalJacobian_NASM);
969:   return(0);
970: }