xref: /petsc/src/snes/impls/nasm/nasm.c (revision 1795a4d16c893ec2fc06bbbc6c5ce592a2de75d4)
1 #include <petsc-private/snesimpl.h>             /*I   "petscsnes.h"   I*/
2 #include <petscdm.h>
3 
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 */
19 
20   /* logging events */
21   PetscLogEvent eventrestrictinterp;
22   PetscLogEvent eventsubsolve;
23 
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;
27 
28 const char *const SNESNASMTypes[] = {"NONE","RESTRICT","INTERPOLATE","BASIC","PCASMType","PC_ASM_",0};
29 const char *const SNESNASMFJTypes[] = {"FINALOUTER","FINALINNER","INITIAL"};
30 
31 #undef __FUNCT__
32 #define __FUNCT__ "SNESReset_NASM"
33 PetscErrorCode SNESReset_NASM(SNES snes)
34 {
35   SNES_NASM      *nasm = (SNES_NASM*)snes->data;
36   PetscErrorCode ierr;
37   PetscInt       i;
38 
39   PetscFunctionBegin;
40   for (i=0; i<nasm->n; i++) {
41     if (nasm->xl) { ierr = VecDestroy(&nasm->xl[i]);CHKERRQ(ierr); }
42     if (nasm->x) { ierr = VecDestroy(&nasm->x[i]);CHKERRQ(ierr); }
43     if (nasm->y) { ierr = VecDestroy(&nasm->y[i]);CHKERRQ(ierr); }
44     if (nasm->b) { ierr = VecDestroy(&nasm->b[i]);CHKERRQ(ierr); }
45 
46     if (nasm->subsnes) { ierr = SNESDestroy(&nasm->subsnes[i]);CHKERRQ(ierr); }
47     if (nasm->oscatter) { ierr = VecScatterDestroy(&nasm->oscatter[i]);CHKERRQ(ierr); }
48     if (nasm->iscatter) { ierr = VecScatterDestroy(&nasm->iscatter[i]);CHKERRQ(ierr); }
49     if (nasm->gscatter) { ierr = VecScatterDestroy(&nasm->gscatter[i]);CHKERRQ(ierr); }
50   }
51 
52   if (nasm->x) {ierr = PetscFree(nasm->x);CHKERRQ(ierr);}
53   if (nasm->xl) {ierr = PetscFree(nasm->xl);CHKERRQ(ierr);}
54   if (nasm->y) {ierr = PetscFree(nasm->y);CHKERRQ(ierr);}
55   if (nasm->b) {ierr = PetscFree(nasm->b);CHKERRQ(ierr);}
56 
57   if (nasm->xinit) {ierr = VecDestroy(&nasm->xinit);CHKERRQ(ierr);}
58 
59   if (nasm->subsnes) {ierr = PetscFree(nasm->subsnes);CHKERRQ(ierr);}
60   if (nasm->oscatter) {ierr = PetscFree(nasm->oscatter);CHKERRQ(ierr);}
61   if (nasm->iscatter) {ierr = PetscFree(nasm->iscatter);CHKERRQ(ierr);}
62   if (nasm->gscatter) {ierr = PetscFree(nasm->gscatter);CHKERRQ(ierr);}
63 
64   nasm->eventrestrictinterp = 0;
65   nasm->eventsubsolve = 0;
66   PetscFunctionReturn(0);
67 }
68 
69 #undef __FUNCT__
70 #define __FUNCT__ "SNESDestroy_NASM"
71 PetscErrorCode SNESDestroy_NASM(SNES snes)
72 {
73   PetscErrorCode ierr;
74 
75   PetscFunctionBegin;
76   ierr = SNESReset_NASM(snes);CHKERRQ(ierr);
77   ierr = PetscFree(snes->data);CHKERRQ(ierr);
78   PetscFunctionReturn(0);
79 }
80 
81 #undef __FUNCT__
82 #define __FUNCT__ "DMGlobalToLocalSubDomainDirichletHook_Private"
83 PetscErrorCode DMGlobalToLocalSubDomainDirichletHook_Private(DM dm,Vec g,InsertMode mode,Vec l,void *ctx)
84 {
85   PetscErrorCode ierr;
86   Vec            bcs = (Vec)ctx;
87 
88   PetscFunctionBegin;
89   ierr = VecCopy(bcs,l);CHKERRQ(ierr);
90   PetscFunctionReturn(0);
91 }
92 
93 #undef __FUNCT__
94 #define __FUNCT__ "SNESSetUp_NASM"
95 PetscErrorCode SNESSetUp_NASM(SNES snes)
96 {
97   SNES_NASM      *nasm = (SNES_NASM*)snes->data;
98   PetscErrorCode ierr;
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;
107 
108   PetscFunctionBegin;
109   if (!nasm->subsnes) {
110     ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr);
111     if (dm) {
112       nasm->usesdm = PETSC_TRUE;
113       ierr         = DMCreateDomainDecomposition(dm,&nasm->n,NULL,NULL,NULL,&subdms);CHKERRQ(ierr);
114       if (!subdms) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"DM has no default decomposition defined.  Set subsolves manually with SNESNASMSetSubdomains().");
115       ierr = DMCreateDomainDecompositionScatters(dm,nasm->n,subdms,&nasm->iscatter,&nasm->oscatter,&nasm->gscatter);CHKERRQ(ierr);
116 
117       ierr = SNESGetOptionsPrefix(snes, &optionsprefix);CHKERRQ(ierr);
118       ierr = PetscMalloc1(nasm->n,&nasm->subsnes);CHKERRQ(ierr);
119       for (i=0; i<nasm->n; i++) {
120         ierr = SNESCreate(PETSC_COMM_SELF,&nasm->subsnes[i]);CHKERRQ(ierr);
121         ierr = SNESAppendOptionsPrefix(nasm->subsnes[i],optionsprefix);CHKERRQ(ierr);
122         ierr = SNESAppendOptionsPrefix(nasm->subsnes[i],"sub_");CHKERRQ(ierr);
123         ierr = SNESSetDM(nasm->subsnes[i],subdms[i]);CHKERRQ(ierr);
124         ierr = MPI_Comm_size(PetscObjectComm((PetscObject)nasm->subsnes[i]),&size);CHKERRQ(ierr);
125         if (size == 1) {
126           ierr = SNESGetKSP(nasm->subsnes[i],&ksp);CHKERRQ(ierr);
127           ierr = KSPGetPC(ksp,&pc);CHKERRQ(ierr);
128           ierr = KSPSetType(ksp,KSPPREONLY);CHKERRQ(ierr);
129           ierr = PCSetType(pc,PCLU);CHKERRQ(ierr);
130         }
131         ierr = SNESSetFromOptions(nasm->subsnes[i]);CHKERRQ(ierr);
132         ierr = DMDestroy(&subdms[i]);CHKERRQ(ierr);
133       }
134       ierr = PetscFree(subdms);CHKERRQ(ierr);
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     ierr = PetscCalloc1(nasm->n,&nasm->x);CHKERRQ(ierr);
140   }
141   if (!nasm->xl) {
142     ierr = PetscCalloc1(nasm->n,&nasm->xl);CHKERRQ(ierr);
143   }
144   if (!nasm->y) {
145     ierr = PetscCalloc1(nasm->n,&nasm->y);CHKERRQ(ierr);
146   }
147   if (!nasm->b) {
148     ierr = PetscCalloc1(nasm->n,&nasm->b);CHKERRQ(ierr);
149   }
150 
151   for (i=0; i<nasm->n; i++) {
152     ierr = SNESGetFunction(nasm->subsnes[i],&F,NULL,NULL);CHKERRQ(ierr);
153     if (!nasm->x[i]) {ierr = VecDuplicate(F,&nasm->x[i]);CHKERRQ(ierr);}
154     if (!nasm->y[i]) {ierr = VecDuplicate(F,&nasm->y[i]);CHKERRQ(ierr);}
155     if (!nasm->b[i]) {ierr = VecDuplicate(F,&nasm->b[i]);CHKERRQ(ierr);}
156     if (!nasm->xl[i]) {
157       ierr = SNESGetDM(nasm->subsnes[i],&subdm);CHKERRQ(ierr);
158       ierr = DMCreateLocalVector(subdm,&nasm->xl[i]);CHKERRQ(ierr);
159     }
160     ierr = DMGlobalToLocalHookAdd(subdm,DMGlobalToLocalSubDomainDirichletHook_Private,NULL,nasm->xl[i]);CHKERRQ(ierr);
161   }
162   if (nasm->finaljacobian) {
163     ierr = SNESSetUpMatrices(snes);CHKERRQ(ierr);
164     if (nasm->fjtype == 2) {
165       ierr = VecDuplicate(snes->vec_sol,&nasm->xinit);CHKERRQ(ierr);
166     }
167     for (i=0; i<nasm->n;i++) {
168       ierr = SNESSetUpMatrices(nasm->subsnes[i]);CHKERRQ(ierr);
169     }
170   }
171   PetscFunctionReturn(0);
172 }
173 
174 #undef __FUNCT__
175 #define __FUNCT__ "SNESSetFromOptions_NASM"
176 PetscErrorCode SNESSetFromOptions_NASM(SNES snes)
177 {
178   PetscErrorCode    ierr;
179   PCASMType         asmtype;
180   PetscBool         flg,monflg,subviewflg;
181   SNES_NASM         *nasm = (SNES_NASM*)snes->data;
182 
183   PetscFunctionBegin;
184   ierr = PetscOptionsHead("Nonlinear Additive Schwartz options");CHKERRQ(ierr);
185   ierr = PetscOptionsEnum("-snes_nasm_type","Type of restriction/extension","",SNESNASMTypes,(PetscEnum)nasm->type,(PetscEnum*)&asmtype,&flg);CHKERRQ(ierr);
186   if (flg) nasm->type = asmtype;
187   flg    = PETSC_FALSE;
188   monflg = PETSC_TRUE;
189   ierr   = PetscOptionsReal("-snes_nasm_damping","Log times for subSNES solves and restriction","SNESNASMSetDamping",nasm->damping,&nasm->damping,&flg);CHKERRQ(ierr);
190   if (flg) {ierr = SNESNASMSetDamping(snes,nasm->damping);CHKERRQ(ierr);}
191   subviewflg = PETSC_FALSE;
192   ierr   = PetscOptionsBool("-snes_nasm_sub_view","Print detailed information for every processor when using -snes_view","",subviewflg,&subviewflg,&flg);CHKERRQ(ierr);
193   if (flg) {
194     nasm->same_local_solves = PETSC_FALSE;
195     if (!subviewflg) {
196       nasm->same_local_solves = PETSC_TRUE;
197     }
198   }
199   ierr   = PetscOptionsBool("-snes_nasm_finaljacobian","Compute the global jacobian of the final iterate (for ASPIN)","",nasm->finaljacobian,&nasm->finaljacobian,NULL);CHKERRQ(ierr);
200   ierr   = PetscOptionsEList("-snes_nasm_finaljacobian_type","The type of the final jacobian computed.","",SNESNASMFJTypes,3,SNESNASMFJTypes[0],&nasm->fjtype,NULL);CHKERRQ(ierr);
201   ierr   = PetscOptionsBool("-snes_nasm_log","Log times for subSNES solves and restriction","",monflg,&monflg,&flg);CHKERRQ(ierr);
202   if (flg) {
203     ierr = PetscLogEventRegister("SNESNASMSubSolve",((PetscObject)snes)->classid,&nasm->eventsubsolve);CHKERRQ(ierr);
204     ierr = PetscLogEventRegister("SNESNASMRestrict",((PetscObject)snes)->classid,&nasm->eventrestrictinterp);CHKERRQ(ierr);
205   }
206   ierr = PetscOptionsTail();CHKERRQ(ierr);
207   PetscFunctionReturn(0);
208 }
209 
210 #undef __FUNCT__
211 #define __FUNCT__ "SNESView_NASM"
212 PetscErrorCode SNESView_NASM(SNES snes, PetscViewer viewer)
213 {
214   SNES_NASM      *nasm = (SNES_NASM*)snes->data;
215   PetscErrorCode ierr;
216   PetscMPIInt    rank,size;
217   PetscInt       i,j,N,bsz;
218   PetscBool      iascii,isstring;
219   PetscViewer    sviewer;
220   MPI_Comm       comm;
221 
222   PetscFunctionBegin;
223   ierr = PetscObjectGetComm((PetscObject)snes,&comm);CHKERRQ(ierr);
224   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
225   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring);CHKERRQ(ierr);
226   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
227   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
228   ierr = MPI_Allreduce(&nasm->n,&N,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(ierr);
229   if (iascii) {
230     ierr = PetscViewerASCIIPrintf(viewer, "  Nonlinear Additive Schwarz: total subdomain blocks = %D\n",N);CHKERRQ(ierr);
231     if (nasm->same_local_solves) {
232       if (nasm->subsnes) {
233         ierr = PetscViewerASCIIPrintf(viewer,"  Local solve is the same for all blocks:\n");CHKERRQ(ierr);
234         ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
235         ierr = PetscViewerGetSingleton(viewer,&sviewer);CHKERRQ(ierr);
236         if (!rank) {
237           ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
238           ierr = SNESView(nasm->subsnes[0],sviewer);CHKERRQ(ierr);
239           ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
240         }
241         ierr = PetscViewerRestoreSingleton(viewer,&sviewer);CHKERRQ(ierr);
242         ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
243       }
244     } else {
245       /* print the solver on each block */
246       ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_TRUE);CHKERRQ(ierr);
247       ierr = PetscViewerASCIISynchronizedPrintf(viewer,"  [%d] number of local blocks = %D\n",(int)rank,nasm->n);CHKERRQ(ierr);
248       ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
249       ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_FALSE);CHKERRQ(ierr);
250       ierr = PetscViewerASCIIPrintf(viewer,"  Local solve info for each block is in the following SNES objects:\n");CHKERRQ(ierr);
251       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
252       ierr = PetscViewerASCIIPrintf(viewer,"- - - - - - - - - - - - - - - - - -\n");CHKERRQ(ierr);
253       for (j=0; j<size; j++) {
254         ierr = PetscViewerGetSingleton(viewer,&sviewer);CHKERRQ(ierr);
255         if (rank == j) {
256           for (i=0; i<nasm->n; i++) {
257             ierr = VecGetLocalSize(nasm->x[i],&bsz);CHKERRQ(ierr);
258             ierr = PetscViewerASCIIPrintf(sviewer,"[%d] local block number %D, size = %D\n",(int)rank,i,bsz);CHKERRQ(ierr);
259             ierr = SNESView(nasm->subsnes[i],sviewer);CHKERRQ(ierr);
260             ierr = PetscViewerASCIIPrintf(sviewer,"- - - - - - - - - - - - - - - - - -\n");CHKERRQ(ierr);
261           }
262         }
263         ierr = PetscViewerRestoreSingleton(viewer,&sviewer);CHKERRQ(ierr);
264         ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
265       }
266       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
267     }
268   } else if (isstring) {
269     ierr = PetscViewerStringSPrintf(viewer," blocks=%D,type=%s",N,SNESNASMTypes[nasm->type]);CHKERRQ(ierr);
270     ierr = PetscViewerGetSingleton(viewer,&sviewer);CHKERRQ(ierr);
271     if (nasm->subsnes && !rank) {ierr = SNESView(nasm->subsnes[0],sviewer);CHKERRQ(ierr);}
272     ierr = PetscViewerRestoreSingleton(viewer,&sviewer);CHKERRQ(ierr);
273   }
274   PetscFunctionReturn(0);
275 }
276 
277 #undef __FUNCT__
278 #define __FUNCT__ "SNESNASMSetSubdomains"
279 /*@
280    SNESNASMSetSubdomains - Manually Set the context required to restrict and solve subdomain problems.
281 
282    Not Collective
283 
284    Input Parameters:
285 +  SNES - the SNES context
286 .  n - the number of local subdomains
287 .  subsnes - solvers defined on the local subdomains
288 .  iscatter - scatters into the nonoverlapping portions of the local subdomains
289 .  oscatter - scatters into the overlapping portions of the local subdomains
290 -  gscatter - scatters into the (ghosted) local vector of the local subdomain
291 
292    Level: intermediate
293 
294 .keywords: SNES, NASM
295 
296 .seealso: SNESNASM, SNESNASMGetSubdomains()
297 @*/
298 PetscErrorCode SNESNASMSetSubdomains(SNES snes,PetscInt n,SNES subsnes[],VecScatter iscatter[],VecScatter oscatter[],VecScatter gscatter[])
299 {
300   PetscErrorCode ierr;
301   PetscErrorCode (*f)(SNES,PetscInt,SNES*,VecScatter*,VecScatter*,VecScatter*);
302 
303   PetscFunctionBegin;
304   ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESNASMSetSubdomains_C",&f);CHKERRQ(ierr);
305   if (f) {ierr = (f)(snes,n,subsnes,iscatter,oscatter,gscatter);CHKERRQ(ierr);}
306   PetscFunctionReturn(0);
307 }
308 
309 #undef __FUNCT__
310 #define __FUNCT__ "SNESNASMSetSubdomains_NASM"
311 PetscErrorCode SNESNASMSetSubdomains_NASM(SNES snes,PetscInt n,SNES subsnes[],VecScatter iscatter[],VecScatter oscatter[],VecScatter gscatter[])
312 {
313   PetscInt       i;
314   PetscErrorCode ierr;
315   SNES_NASM      *nasm = (SNES_NASM*)snes->data;
316 
317   PetscFunctionBegin;
318   if (snes->setupcalled) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"SNESNASMSetSubdomains() should be called before calling SNESSetUp().");
319 
320   /* tear down the previously set things */
321   ierr = SNESReset(snes);CHKERRQ(ierr);
322 
323   nasm->n = n;
324   if (oscatter) {
325     for (i=0; i<n; i++) {ierr = PetscObjectReference((PetscObject)oscatter[i]);CHKERRQ(ierr);}
326   }
327   if (iscatter) {
328     for (i=0; i<n; i++) {ierr = PetscObjectReference((PetscObject)iscatter[i]);CHKERRQ(ierr);}
329   }
330   if (gscatter) {
331     for (i=0; i<n; i++) {ierr = PetscObjectReference((PetscObject)gscatter[i]);CHKERRQ(ierr);}
332   }
333   if (oscatter) {
334     ierr = PetscMalloc1(n,&nasm->oscatter);CHKERRQ(ierr);
335     for (i=0; i<n; i++) {
336       nasm->oscatter[i] = oscatter[i];
337     }
338   }
339   if (iscatter) {
340     ierr = PetscMalloc1(n,&nasm->iscatter);CHKERRQ(ierr);
341     for (i=0; i<n; i++) {
342       nasm->iscatter[i] = iscatter[i];
343     }
344   }
345   if (gscatter) {
346     ierr = PetscMalloc1(n,&nasm->gscatter);CHKERRQ(ierr);
347     for (i=0; i<n; i++) {
348       nasm->gscatter[i] = gscatter[i];
349     }
350   }
351 
352   if (subsnes) {
353     ierr = PetscMalloc1(n,&nasm->subsnes);CHKERRQ(ierr);
354     for (i=0; i<n; i++) {
355       nasm->subsnes[i] = subsnes[i];
356     }
357     nasm->same_local_solves = PETSC_FALSE;
358   }
359   PetscFunctionReturn(0);
360 }
361 
362 #undef __FUNCT__
363 #define __FUNCT__ "SNESNASMGetSubdomains"
364 /*@
365    SNESNASMGetSubdomains - Get the local subdomain context.
366 
367    Not Collective
368 
369    Input Parameters:
370 .  SNES - the SNES context
371 
372    Output Parameters:
373 +  n - the number of local subdomains
374 .  subsnes - solvers defined on the local subdomains
375 .  iscatter - scatters into the nonoverlapping portions of the local subdomains
376 .  oscatter - scatters into the overlapping portions of the local subdomains
377 -  gscatter - scatters into the (ghosted) local vector of the local subdomain
378 
379    Level: intermediate
380 
381 .keywords: SNES, NASM
382 
383 .seealso: SNESNASM, SNESNASMSetSubdomains()
384 @*/
385 PetscErrorCode SNESNASMGetSubdomains(SNES snes,PetscInt *n,SNES *subsnes[],VecScatter *iscatter[],VecScatter *oscatter[],VecScatter *gscatter[])
386 {
387   PetscErrorCode ierr;
388   PetscErrorCode (*f)(SNES,PetscInt*,SNES**,VecScatter**,VecScatter**,VecScatter**);
389 
390   PetscFunctionBegin;
391   ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESNASMGetSubdomains_C",&f);CHKERRQ(ierr);
392   if (f) {ierr = (f)(snes,n,subsnes,iscatter,oscatter,gscatter);CHKERRQ(ierr);}
393   PetscFunctionReturn(0);
394 }
395 
396 #undef __FUNCT__
397 #define __FUNCT__ "SNESNASMGetSubdomains_NASM"
398 PetscErrorCode SNESNASMGetSubdomains_NASM(SNES snes,PetscInt *n,SNES *subsnes[],VecScatter *iscatter[],VecScatter *oscatter[],VecScatter *gscatter[])
399 {
400   SNES_NASM      *nasm = (SNES_NASM*)snes->data;
401 
402   PetscFunctionBegin;
403   if (n) *n = nasm->n;
404   if (oscatter) *oscatter = nasm->oscatter;
405   if (iscatter) *iscatter = nasm->iscatter;
406   if (gscatter) *gscatter = nasm->gscatter;
407   if (subsnes)  {
408     *subsnes  = nasm->subsnes;
409     nasm->same_local_solves = PETSC_FALSE;
410   }
411   PetscFunctionReturn(0);
412 }
413 
414 #undef __FUNCT__
415 #define __FUNCT__ "SNESNASMGetSubdomainVecs"
416 /*@
417    SNESNASMGetSubdomainVecs - Get the processor-local subdomain vectors
418 
419    Not Collective
420 
421    Input Parameters:
422 .  SNES - the SNES context
423 
424    Output Parameters:
425 +  n - the number of local subdomains
426 .  x - The subdomain solution vector
427 .  y - The subdomain step vector
428 .  b - The subdomain RHS vector
429 -  xl - The subdomain local vectors (ghosted)
430 
431    Level: developer
432 
433 .keywords: SNES, NASM
434 
435 .seealso: SNESNASM, SNESNASMGetSubdomains()
436 @*/
437 PetscErrorCode SNESNASMGetSubdomainVecs(SNES snes,PetscInt *n,Vec **x,Vec **y,Vec **b, Vec **xl)
438 {
439   PetscErrorCode ierr;
440   PetscErrorCode (*f)(SNES,PetscInt*,Vec**,Vec**,Vec**,Vec**);
441 
442   PetscFunctionBegin;
443   ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESNASMGetSubdomainVecs_C",&f);CHKERRQ(ierr);
444   if (f) {ierr = (f)(snes,n,x,y,b,xl);CHKERRQ(ierr);}
445   PetscFunctionReturn(0);
446 }
447 
448 #undef __FUNCT__
449 #define __FUNCT__ "SNESNASMGetSubdomainVecs_NASM"
450 PetscErrorCode SNESNASMGetSubdomainVecs_NASM(SNES snes,PetscInt *n,Vec **x,Vec **y,Vec **b,Vec **xl)
451 {
452   SNES_NASM      *nasm = (SNES_NASM*)snes->data;
453 
454   PetscFunctionBegin;
455   if (n)  *n  = nasm->n;
456   if (x)  *x  = nasm->x;
457   if (y)  *y  = nasm->y;
458   if (b)  *b  = nasm->b;
459   if (xl) *xl = nasm->xl;
460   PetscFunctionReturn(0);
461 }
462 
463 #undef __FUNCT__
464 #define __FUNCT__ "SNESNASMSetComputeFinalJacobian"
465 /*@
466    SNESNASMSetComputeFinalJacobian - Schedules the computation of the global and subdomain jacobians upon convergence
467 
468    Collective on SNES
469 
470    Input Parameters:
471 +  SNES - the SNES context
472 -  flg - indication of whether to compute the jacobians or not
473 
474    Level: developer
475 
476    Notes: This is used almost exclusively in the implementation of ASPIN, where the converged subdomain and global jacobian
477    is needed at each linear iteration.
478 
479 .keywords: SNES, NASM, ASPIN
480 
481 .seealso: SNESNASM, SNESNASMGetSubdomains()
482 @*/
483 PetscErrorCode SNESNASMSetComputeFinalJacobian(SNES snes,PetscBool flg)
484 {
485   PetscErrorCode (*f)(SNES,PetscBool);
486   PetscErrorCode ierr;
487 
488   PetscFunctionBegin;
489   ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESNASMSetComputeFinalJacobian_C",&f);CHKERRQ(ierr);
490   if (f) {ierr = (f)(snes,flg);CHKERRQ(ierr);}
491   PetscFunctionReturn(0);
492 }
493 
494 #undef __FUNCT__
495 #define __FUNCT__ "SNESNASMSetComputeFinalJacobian_NASM"
496 PetscErrorCode SNESNASMSetComputeFinalJacobian_NASM(SNES snes,PetscBool flg)
497 {
498   SNES_NASM      *nasm = (SNES_NASM*)snes->data;
499 
500   PetscFunctionBegin;
501   nasm->finaljacobian = flg;
502   if (flg) snes->usesksp = PETSC_TRUE;
503   PetscFunctionReturn(0);
504 }
505 
506 #undef __FUNCT__
507 #define __FUNCT__ "SNESNASMSetDamping"
508 /*@
509    SNESNASMSetDamping - Sets the update damping for NASM
510 
511    Logically collective on SNES
512 
513    Input Parameters:
514 +  SNES - the SNES context
515 -  dmp - damping
516 
517    Level: intermediate
518 
519 .keywords: SNES, NASM, damping
520 
521 .seealso: SNESNASM, SNESNASMGetDamping()
522 @*/
523 PetscErrorCode SNESNASMSetDamping(SNES snes,PetscReal dmp)
524 {
525   PetscErrorCode (*f)(SNES,PetscReal);
526   PetscErrorCode ierr;
527 
528   PetscFunctionBegin;
529   ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESNASMSetDamping_C",(void (**)(void))&f);CHKERRQ(ierr);
530   if (f) {ierr = (f)(snes,dmp);CHKERRQ(ierr);}
531   PetscFunctionReturn(0);
532 }
533 
534 #undef __FUNCT__
535 #define __FUNCT__ "SNESNASMSetDamping_NASM"
536 PetscErrorCode SNESNASMSetDamping_NASM(SNES snes,PetscReal dmp)
537 {
538   SNES_NASM      *nasm = (SNES_NASM*)snes->data;
539 
540   PetscFunctionBegin;
541   nasm->damping = dmp;
542   PetscFunctionReturn(0);
543 }
544 
545 #undef __FUNCT__
546 #define __FUNCT__ "SNESNASMGetDamping"
547 /*@
548    SNESNASMGetDamping - Gets the update damping for NASM
549 
550    Not Collective
551 
552    Input Parameters:
553 +  SNES - the SNES context
554 -  dmp - damping
555 
556    Level: intermediate
557 
558 .keywords: SNES, NASM, damping
559 
560 .seealso: SNESNASM, SNESNASMSetDamping()
561 @*/
562 PetscErrorCode SNESNASMGetDamping(SNES snes,PetscReal *dmp)
563 {
564   PetscErrorCode (*f)(SNES,PetscReal*);
565   PetscErrorCode ierr;
566 
567   PetscFunctionBegin;
568   ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESNASMGetDamping_C",(void (**)(void))&f);CHKERRQ(ierr);
569   if (f) {ierr = (f)(snes,dmp);CHKERRQ(ierr);}
570   PetscFunctionReturn(0);
571 }
572 
573 #undef __FUNCT__
574 #define __FUNCT__ "SNESNASMGetDamping_NASM"
575 PetscErrorCode SNESNASMGetDamping_NASM(SNES snes,PetscReal *dmp)
576 {
577   SNES_NASM      *nasm = (SNES_NASM*)snes->data;
578 
579   PetscFunctionBegin;
580   *dmp = nasm->damping;
581   PetscFunctionReturn(0);
582 }
583 
584 
585 #undef __FUNCT__
586 #define __FUNCT__ "SNESNASMSolveLocal_Private"
587 PetscErrorCode SNESNASMSolveLocal_Private(SNES snes,Vec B,Vec Y,Vec X)
588 {
589   SNES_NASM      *nasm = (SNES_NASM*)snes->data;
590   SNES           subsnes;
591   PetscInt       i;
592   PetscReal      dmp;
593   PetscErrorCode ierr;
594   Vec            Xlloc,Xl,Bl,Yl;
595   VecScatter     iscat,oscat,gscat;
596   DM             dm,subdm;
597 
598   PetscFunctionBegin;
599   ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr);
600   ierr = SNESNASMGetDamping(snes,&dmp);CHKERRQ(ierr);
601   ierr = VecSet(Y,0);CHKERRQ(ierr);
602   if (nasm->eventrestrictinterp) {ierr = PetscLogEventBegin(nasm->eventrestrictinterp,snes,0,0,0);CHKERRQ(ierr);}
603   for (i=0; i<nasm->n; i++) {
604     /* scatter the solution to the local solution */
605     Xlloc = nasm->xl[i];
606     gscat   = nasm->gscatter[i];
607     oscat   = nasm->oscatter[i];
608     ierr = VecScatterBegin(gscat,X,Xlloc,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
609     if (B) {
610       /* scatter the RHS to the local RHS */
611       Bl   = nasm->b[i];
612       ierr = VecScatterBegin(oscat,B,Bl,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
613     }
614   }
615   if (nasm->eventrestrictinterp) {ierr = PetscLogEventEnd(nasm->eventrestrictinterp,snes,0,0,0);CHKERRQ(ierr);}
616 
617 
618   if (nasm->eventsubsolve) {ierr = PetscLogEventBegin(nasm->eventsubsolve,snes,0,0,0);CHKERRQ(ierr);}
619   for (i=0; i<nasm->n; i++) {
620     Xl    = nasm->x[i];
621     Xlloc = nasm->xl[i];
622     Yl    = nasm->y[i];
623     subsnes = nasm->subsnes[i];
624     ierr    = SNESGetDM(subsnes,&subdm);CHKERRQ(ierr);
625     iscat   = nasm->iscatter[i];
626     oscat   = nasm->oscatter[i];
627     gscat   = nasm->gscatter[i];
628     ierr = VecScatterEnd(gscat,X,Xlloc,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
629     if (B) {
630       Bl   = nasm->b[i];
631       ierr = VecScatterEnd(oscat,B,Bl,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
632     } else Bl = NULL;
633     ierr = DMSubDomainRestrict(dm,oscat,gscat,subdm);CHKERRQ(ierr);
634     ierr = DMLocalToGlobalBegin(subdm,Xlloc,INSERT_VALUES,Xl);CHKERRQ(ierr);
635     ierr = DMLocalToGlobalEnd(subdm,Xlloc,INSERT_VALUES,Xl);CHKERRQ(ierr);
636     ierr = VecCopy(Xl,Yl);CHKERRQ(ierr);
637     ierr = SNESSolve(subsnes,Bl,Xl);CHKERRQ(ierr);
638     ierr = VecAYPX(Yl,-1.0,Xl);CHKERRQ(ierr);
639     if (nasm->type == PC_ASM_BASIC) {
640       ierr = VecScatterBegin(oscat,Yl,Y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
641     } else if (nasm->type == PC_ASM_RESTRICT) {
642       ierr = VecScatterBegin(iscat,Yl,Y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
643     } else SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"Only basic and restrict types are supported for SNESNASM");
644   }
645   if (nasm->eventsubsolve) {ierr = PetscLogEventEnd(nasm->eventsubsolve,snes,0,0,0);CHKERRQ(ierr);}
646   if (nasm->eventrestrictinterp) {ierr = PetscLogEventBegin(nasm->eventrestrictinterp,snes,0,0,0);CHKERRQ(ierr);}
647   for (i=0; i<nasm->n; i++) {
648     Yl    = nasm->y[i];
649     iscat   = nasm->iscatter[i];
650     oscat   = nasm->oscatter[i];
651     if (nasm->type == PC_ASM_BASIC) {
652       ierr = VecScatterEnd(oscat,Yl,Y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
653     } else if (nasm->type == PC_ASM_RESTRICT) {
654       ierr = VecScatterEnd(iscat,Yl,Y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
655     } else SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"Only basic and restrict types are supported for SNESNASM");
656   }
657   if (nasm->eventrestrictinterp) {ierr = PetscLogEventEnd(nasm->eventrestrictinterp,snes,0,0,0);CHKERRQ(ierr);}
658   ierr = VecAXPY(X,dmp,Y);CHKERRQ(ierr);
659   PetscFunctionReturn(0);
660 }
661 
662 #undef __FUNCT__
663 #define __FUNCT__ "SNESNASMComputeFinalJacobian_Private"
664 PetscErrorCode SNESNASMComputeFinalJacobian_Private(SNES snes, Vec Xfinal)
665 {
666   Vec            X = Xfinal;
667   SNES_NASM      *nasm = (SNES_NASM*)snes->data;
668   SNES           subsnes;
669   PetscInt       i,lag = 1;
670   PetscErrorCode ierr;
671   Vec            Xlloc,Xl,Fl,F;
672   VecScatter     oscat,gscat;
673   DM             dm,subdm;
674   MatStructure   flg = DIFFERENT_NONZERO_PATTERN;
675   PetscFunctionBegin;
676   if (nasm->fjtype == 2) X = nasm->xinit;
677   F = snes->vec_func;
678   if (snes->normschedule == SNES_NORM_NONE) {ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);}
679   ierr = SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg);CHKERRQ(ierr);
680   ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr);
681   if (nasm->eventrestrictinterp) {ierr = PetscLogEventBegin(nasm->eventrestrictinterp,snes,0,0,0);CHKERRQ(ierr);}
682   if (nasm->fjtype != 1) {
683     for (i=0; i<nasm->n; i++) {
684       Xlloc = nasm->xl[i];
685       gscat = nasm->gscatter[i];
686       oscat = nasm->oscatter[i];
687       ierr = VecScatterBegin(gscat,X,Xlloc,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
688     }
689   }
690   if (nasm->eventrestrictinterp) {ierr = PetscLogEventEnd(nasm->eventrestrictinterp,snes,0,0,0);CHKERRQ(ierr);}
691   for (i=0; i<nasm->n; i++) {
692     Fl      = nasm->subsnes[i]->vec_func;
693     Xl      = nasm->x[i];
694     Xlloc   = nasm->xl[i];
695     subsnes = nasm->subsnes[i];
696     oscat   = nasm->oscatter[i];
697     gscat   = nasm->gscatter[i];
698     if (nasm->fjtype != 1) {ierr = VecScatterEnd(gscat,X,Xlloc,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);}
699     ierr = SNESGetDM(subsnes,&subdm);CHKERRQ(ierr);
700     ierr = DMSubDomainRestrict(dm,oscat,gscat,subdm);CHKERRQ(ierr);
701     if (nasm->fjtype != 1) {
702       ierr = DMLocalToGlobalBegin(subdm,Xlloc,INSERT_VALUES,Xl);CHKERRQ(ierr);
703       ierr = DMLocalToGlobalEnd(subdm,Xlloc,INSERT_VALUES,Xl);CHKERRQ(ierr);
704     }
705     if (subsnes->lagjacobian == -1)    subsnes->lagjacobian = -2;
706     else if (subsnes->lagjacobian > 1) lag = subsnes->lagjacobian;
707     ierr = SNESComputeFunction(subsnes,Xl,Fl);CHKERRQ(ierr);
708     ierr = SNESComputeJacobian(subsnes,Xl,&subsnes->jacobian,&subsnes->jacobian_pre,&flg);CHKERRQ(ierr);
709     if (lag > 1) subsnes->lagjacobian = lag;
710     ierr = KSPSetOperators(subsnes->ksp,subsnes->jacobian,subsnes->jacobian_pre,flg);CHKERRQ(ierr);
711   }
712   PetscFunctionReturn(0);
713 }
714 
715 #undef __FUNCT__
716 #define __FUNCT__ "SNESSolve_NASM"
717 PetscErrorCode SNESSolve_NASM(SNES snes)
718 {
719   Vec              F;
720   Vec              X;
721   Vec              B;
722   Vec              Y;
723   PetscInt         i;
724   PetscReal        fnorm = 0.0;
725   PetscErrorCode   ierr;
726   SNESNormSchedule normschedule;
727   SNES_NASM        *nasm = (SNES_NASM*)snes->data;
728 
729   PetscFunctionBegin;
730   X = snes->vec_sol;
731   Y = snes->vec_sol_update;
732   F = snes->vec_func;
733   B = snes->vec_rhs;
734 
735   ierr         = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr);
736   snes->iter   = 0;
737   snes->norm   = 0.;
738   ierr         = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr);
739   snes->reason = SNES_CONVERGED_ITERATING;
740   ierr         = SNESGetNormSchedule(snes, &normschedule);CHKERRQ(ierr);
741   if (normschedule == SNES_NORM_ALWAYS || normschedule == SNES_NORM_INITIAL_ONLY || normschedule == SNES_NORM_INITIAL_FINAL_ONLY) {
742     /* compute the initial function and preconditioned update delX */
743     if (!snes->vec_func_init_set) {
744       ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
745       if (snes->domainerror) {
746         snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
747         PetscFunctionReturn(0);
748       }
749     } else snes->vec_func_init_set = PETSC_FALSE;
750 
751     ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F||  */
752     if (PetscIsInfOrNanReal(fnorm)) {
753       snes->reason = SNES_DIVERGED_FNORM_NAN;
754       PetscFunctionReturn(0);
755     }
756     ierr       = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr);
757     snes->iter = 0;
758     snes->norm = fnorm;
759     ierr       = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr);
760     ierr       = SNESLogConvergenceHistory(snes,snes->norm,0);CHKERRQ(ierr);
761     ierr       = SNESMonitor(snes,0,snes->norm);CHKERRQ(ierr);
762 
763     /* test convergence */
764     ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
765     if (snes->reason) PetscFunctionReturn(0);
766   } else {
767     ierr = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr);
768     ierr = SNESLogConvergenceHistory(snes,snes->norm,0);CHKERRQ(ierr);
769     ierr = SNESMonitor(snes,0,snes->norm);CHKERRQ(ierr);
770   }
771 
772   /* Call general purpose update function */
773   if (snes->ops->update) {
774     ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
775   }
776   /* copy the initial solution over for later */
777   if (nasm->fjtype == 2) {ierr = VecCopy(X,nasm->xinit);CHKERRQ(ierr);}
778 
779   for (i = 0; i < snes->max_its; i++) {
780     ierr = SNESNASMSolveLocal_Private(snes,B,Y,X);CHKERRQ(ierr);
781     if (normschedule == SNES_NORM_ALWAYS || ((i == snes->max_its - 1) && (normschedule == SNES_NORM_INITIAL_FINAL_ONLY || normschedule == SNES_NORM_FINAL_ONLY))) {
782       ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
783       if (snes->domainerror) {
784         snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
785         break;
786       }
787       ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F||  */
788       if (PetscIsInfOrNanReal(fnorm)) {
789         snes->reason = SNES_DIVERGED_FNORM_NAN;
790         break;
791       }
792     }
793     /* Monitor convergence */
794     ierr       = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr);
795     snes->iter = i+1;
796     snes->norm = fnorm;
797     ierr       = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr);
798     ierr       = SNESLogConvergenceHistory(snes,snes->norm,0);CHKERRQ(ierr);
799     ierr       = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr);
800     /* Test for convergence */
801     if (normschedule == SNES_NORM_ALWAYS) {ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);}
802     if (snes->reason) break;
803     /* Call general purpose update function */
804     if (snes->ops->update) {ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);}
805   }
806   if (nasm->finaljacobian) {ierr = SNESNASMComputeFinalJacobian_Private(snes,X);CHKERRQ(ierr);}
807   if (normschedule == SNES_NORM_ALWAYS) {
808     if (i == snes->max_its) {
809       ierr = PetscInfo1(snes,"Maximum number of iterations has been reached: %D\n",snes->max_its);CHKERRQ(ierr);
810       if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
811     }
812   } else if (!snes->reason) snes->reason = SNES_CONVERGED_ITS; /* NASM is meant to be used as a preconditioner */
813   PetscFunctionReturn(0);
814 }
815 
816 /*MC
817   SNESNASM - Nonlinear Additive Schwartz
818 
819    Options Database:
820 +  -snes_nasm_log - enable logging events for the communication and solve stages
821 .  -snes_nasm_type <basic,restrict> - type of subdomain update used
822 .  -snes_nasm_finaljacobian - compute the local and global jacobians of the final iterate
823 .  -snes_nasm_finaljacobian_type <finalinner,finalouter,initial> pick state the jacobian is calculated at
824 .  -sub_snes_ - options prefix of the subdomain nonlinear solves
825 .  -sub_ksp_ - options prefix of the subdomain Krylov solver
826 -  -sub_pc_ - options prefix of the subdomain preconditioner
827 
828    Level: advanced
829 
830 .seealso: SNESCreate(), SNES, SNESSetType(), SNESType (for list of available types)
831 M*/
832 
833 #undef __FUNCT__
834 #define __FUNCT__ "SNESCreate_NASM"
835 PETSC_EXTERN PetscErrorCode SNESCreate_NASM(SNES snes)
836 {
837   SNES_NASM      *nasm;
838   PetscErrorCode ierr;
839 
840   PetscFunctionBegin;
841   ierr       = PetscNewLog(snes, SNES_NASM, &nasm);CHKERRQ(ierr);
842   snes->data = (void*)nasm;
843 
844   nasm->n        = PETSC_DECIDE;
845   nasm->subsnes  = 0;
846   nasm->x        = 0;
847   nasm->xl       = 0;
848   nasm->y        = 0;
849   nasm->b        = 0;
850   nasm->oscatter = 0;
851   nasm->iscatter = 0;
852   nasm->gscatter = 0;
853   nasm->damping  = 1.;
854 
855   nasm->type = PC_ASM_BASIC;
856   nasm->finaljacobian = PETSC_FALSE;
857   nasm->same_local_solves = PETSC_TRUE;
858 
859   snes->ops->destroy        = SNESDestroy_NASM;
860   snes->ops->setup          = SNESSetUp_NASM;
861   snes->ops->setfromoptions = SNESSetFromOptions_NASM;
862   snes->ops->view           = SNESView_NASM;
863   snes->ops->solve          = SNESSolve_NASM;
864   snes->ops->reset          = SNESReset_NASM;
865 
866   snes->usesksp = PETSC_FALSE;
867   snes->usespc  = PETSC_FALSE;
868 
869   nasm->fjtype              = 0;
870   nasm->xinit               = NULL;
871   nasm->eventrestrictinterp = 0;
872   nasm->eventsubsolve       = 0;
873 
874   if (!snes->tolerancesset) {
875     snes->max_its   = 10000;
876     snes->max_funcs = 10000;
877   }
878 
879   ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESNASMSetSubdomains_C",SNESNASMSetSubdomains_NASM);CHKERRQ(ierr);
880   ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESNASMGetSubdomains_C",SNESNASMGetSubdomains_NASM);CHKERRQ(ierr);
881   ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESNASMSetDamping_C",SNESNASMSetDamping_NASM);CHKERRQ(ierr);
882   ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESNASMGetDamping_C",SNESNASMGetDamping_NASM);CHKERRQ(ierr);
883   ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESNASMGetSubdomainVecs_C",SNESNASMGetSubdomainVecs_NASM);CHKERRQ(ierr);
884   ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESNASMSetComputeFinalJacobian_C",SNESNASMSetComputeFinalJacobian_NASM);CHKERRQ(ierr);
885   PetscFunctionReturn(0);
886 }
887 
888