xref: /petsc/src/snes/impls/composite/snescomposite.c (revision e31ab4f968abf01cfa239171b74d7046a5872add)
1 
2 /*
3       Defines a SNES that can consist of a collection of SNESes
4 */
5 #include <petsc-private/snesimpl.h> /*I "petscsnes.h" I*/
6 
7 const char *const        SNESCompositeTypes[]   = {"ADDITIVE","MULTIPLICATIVE","SNESCompositeType","SNES_COMPOSITE",0};
8 
9 typedef struct _SNES_CompositeLink *SNES_CompositeLink;
10 struct _SNES_CompositeLink {
11   SNES               snes;
12   PetscReal          dmp;
13   SNES_CompositeLink next;
14   SNES_CompositeLink previous;
15 };
16 
17 typedef struct {
18   SNES_CompositeLink head;
19   SNESCompositeType  type;
20   Vec                Xorig;
21 } SNES_Composite;
22 
23 #undef __FUNCT__
24 #define __FUNCT__ "SNESCompositeApply_Multiplicative"
25 static PetscErrorCode SNESCompositeApply_Multiplicative(SNES snes,Vec X,Vec B,Vec F,PetscReal *fnorm)
26 {
27   PetscErrorCode     ierr;
28   SNES_Composite     *jac = (SNES_Composite*)snes->data;
29   SNES_CompositeLink next = jac->head;
30   Vec                FSub;
31   PetscReal          fsubnorm;
32 
33   PetscFunctionBegin;
34   if (!next) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"No composite SNESes supplied via SNESCompositeAddSNES() or -snes_composite_sneses");
35   if (snes->normtype == SNES_NORM_FUNCTION) {
36     ierr = SNESSetInitialFunction(next->snes,F);CHKERRQ(ierr);
37     if (fnorm) {ierr = SNESSetInitialFunctionNorm(next->snes,*fnorm);CHKERRQ(ierr);}
38   }
39   ierr = SNESSolve(next->snes,B,X);CHKERRQ(ierr);
40 
41   while (next->next) {
42     /* only copy the function over in the case where the functions correspond */
43     if (next->snes->pcside == PC_RIGHT && next->snes->normtype != SNES_NORM_NONE) {
44       ierr = SNESGetFunction(next->snes,&FSub,NULL,NULL);CHKERRQ(ierr);
45       ierr = SNESGetFunctionNorm(next->snes,&fsubnorm);CHKERRQ(ierr);
46       next = next->next;
47       ierr = SNESSetInitialFunction(next->snes,FSub);CHKERRQ(ierr);
48       ierr = SNESSetInitialFunctionNorm(next->snes,fsubnorm);CHKERRQ(ierr);
49     } else {
50       next = next->next;
51     }
52     ierr = SNESSolve(next->snes,B,X);CHKERRQ(ierr);
53   }
54   if (next->snes->pcside == PC_RIGHT) {
55     ierr = SNESGetFunction(next->snes,&FSub,NULL,NULL);CHKERRQ(ierr);
56     ierr = VecCopy(FSub,F);CHKERRQ(ierr);
57     if (fnorm) {ierr = SNESGetFunctionNorm(next->snes,fnorm);CHKERRQ(ierr);}
58   } else if (snes->normtype == SNES_NORM_FUNCTION) {
59     SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
60     if (fnorm) {ierr = VecNorm(F,NORM_2,fnorm);CHKERRQ(ierr);}
61   }
62   PetscFunctionReturn(0);
63 }
64 
65 #undef __FUNCT__
66 #define __FUNCT__ "SNESCompositeApply_Additive"
67 static PetscErrorCode SNESCompositeApply_Additive(SNES snes,Vec X,Vec B,Vec F,PetscReal *fnorm)
68 {
69   PetscErrorCode     ierr;
70   SNES_Composite     *jac = (SNES_Composite*)snes->data;
71   SNES_CompositeLink next = jac->head;
72   Vec                Y,Xorig;
73 
74   PetscFunctionBegin;
75   Y = snes->vec_sol_update;
76   if (!jac->Xorig) {ierr = VecDuplicate(X,&jac->Xorig);CHKERRQ(ierr);}
77   Xorig = jac->Xorig;
78   ierr = VecCopy(X,Xorig);
79   if (!next) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"No composite SNESes supplied via SNESCompositeAddSNES() or -snes_composite_sneses");
80   if (snes->normtype == SNES_NORM_FUNCTION) {
81     ierr = SNESSetInitialFunction(next->snes,F);CHKERRQ(ierr);
82     if (fnorm) {ierr = SNESSetInitialFunctionNorm(next->snes,*fnorm);CHKERRQ(ierr);}
83     while (next->next) {
84       next = next->next;
85       ierr = SNESSetInitialFunction(next->snes,F);CHKERRQ(ierr);
86       if (fnorm) {ierr = SNESSetInitialFunctionNorm(next->snes,*fnorm);CHKERRQ(ierr);}
87     }
88   }
89   next = jac->head;
90   ierr = VecCopy(Xorig,Y);CHKERRQ(ierr);
91   ierr = SNESSolve(next->snes,B,Y);CHKERRQ(ierr);
92   ierr = VecAXPY(Y,-1.0,Xorig);CHKERRQ(ierr);
93   ierr = VecAXPY(X,next->dmp,Y);CHKERRQ(ierr);
94   while (next->next) {
95     next = next->next;
96     ierr = VecCopy(Xorig,Y);CHKERRQ(ierr);
97     ierr = SNESSolve(next->snes,B,Y);CHKERRQ(ierr);
98     ierr = VecAXPY(Y,-1.0,Xorig);CHKERRQ(ierr);
99     ierr = VecAXPY(X,next->dmp,Y);CHKERRQ(ierr);
100   }
101   if (snes->normtype == SNES_NORM_FUNCTION) {
102     ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
103     if (fnorm) {ierr = VecNorm(F,NORM_2,fnorm);CHKERRQ(ierr);}
104   }
105   PetscFunctionReturn(0);
106 }
107 #undef __FUNCT__
108 #define __FUNCT__ "SNESSetUp_Composite"
109 static PetscErrorCode SNESSetUp_Composite(SNES snes)
110 {
111   PetscErrorCode   ierr;
112   SNES_Composite     *jac = (SNES_Composite*)snes->data;
113   SNES_CompositeLink next = jac->head;
114 
115   PetscFunctionBegin;
116   while (next) {
117     ierr = SNESSetFromOptions(next->snes);CHKERRQ(ierr);
118     next = next->next;
119   }
120   PetscFunctionReturn(0);
121 }
122 
123 #undef __FUNCT__
124 #define __FUNCT__ "SNESReset_Composite"
125 static PetscErrorCode SNESReset_Composite(SNES snes)
126 {
127   SNES_Composite     *jac = (SNES_Composite*)snes->data;
128   PetscErrorCode   ierr;
129   SNES_CompositeLink next = jac->head;
130 
131   PetscFunctionBegin;
132   while (next) {
133     ierr = SNESReset(next->snes);CHKERRQ(ierr);
134     next = next->next;
135   }
136   ierr = VecDestroy(&jac->Xorig);CHKERRQ(ierr);
137   PetscFunctionReturn(0);
138 }
139 
140 #undef __FUNCT__
141 #define __FUNCT__ "SNESDestroy_Composite"
142 static PetscErrorCode SNESDestroy_Composite(SNES snes)
143 {
144   SNES_Composite     *jac = (SNES_Composite*)snes->data;
145   PetscErrorCode   ierr;
146   SNES_CompositeLink next = jac->head,next_tmp;
147 
148   PetscFunctionBegin;
149   ierr = SNESReset_Composite(snes);CHKERRQ(ierr);
150   while (next) {
151     ierr     = SNESDestroy(&next->snes);CHKERRQ(ierr);
152     next_tmp = next;
153     next     = next->next;
154     ierr     = PetscFree(next_tmp);CHKERRQ(ierr);
155   }
156   ierr = PetscFree(snes->data);CHKERRQ(ierr);
157   PetscFunctionReturn(0);
158 }
159 
160 #undef __FUNCT__
161 #define __FUNCT__ "SNESSetFromOptions_Composite"
162 static PetscErrorCode SNESSetFromOptions_Composite(SNES snes)
163 {
164   SNES_Composite     *jac = (SNES_Composite*)snes->data;
165   PetscErrorCode     ierr;
166   PetscInt           nmax = 8,i;
167   SNES_CompositeLink next;
168   char               *sneses[8];
169   PetscReal          dmps[8];
170   PetscBool          flg;
171 
172   PetscFunctionBegin;
173   ierr = PetscOptionsHead("Composite preconditioner options");CHKERRQ(ierr);
174   ierr = PetscOptionsEnum("-snes_composite_type","Type of composition","SNESCompositeSetType",SNESCompositeTypes,(PetscEnum)jac->type,(PetscEnum*)&jac->type,&flg);CHKERRQ(ierr);
175   if (flg) {
176     ierr = SNESCompositeSetType(snes,jac->type);CHKERRQ(ierr);
177   }
178   ierr = PetscOptionsStringArray("-snes_composite_sneses","List of composite solvers","SNESCompositeAddSNES",sneses,&nmax,&flg);CHKERRQ(ierr);
179   if (flg) {
180     for (i=0; i<nmax; i++) {
181       ierr = SNESCompositeAddSNES(snes,sneses[i]);CHKERRQ(ierr);
182       ierr = PetscFree(sneses[i]);CHKERRQ(ierr);   /* deallocate string sneses[i], which is allocated in PetscOptionsStringArray() */
183     }
184   }
185   ierr = PetscOptionsRealArray("-snes_composite_damping","Damping of the additive composite solvers","SNESCompositeSetDamping",dmps,&nmax,&flg);CHKERRQ(ierr);
186   if (flg) {
187     for (i=0; i<nmax; i++) {
188       ierr = SNESCompositeSetDamping(snes,i,dmps[i]);CHKERRQ(ierr);
189     }
190   }
191   ierr = PetscOptionsTail();CHKERRQ(ierr);
192 
193   next = jac->head;
194   while (next) {
195     ierr = SNESSetFromOptions(next->snes);CHKERRQ(ierr);
196     next = next->next;
197   }
198   PetscFunctionReturn(0);
199 }
200 
201 #undef __FUNCT__
202 #define __FUNCT__ "SNESView_Composite"
203 static PetscErrorCode SNESView_Composite(SNES snes,PetscViewer viewer)
204 {
205   SNES_Composite     *jac = (SNES_Composite*)snes->data;
206   PetscErrorCode   ierr;
207   SNES_CompositeLink next = jac->head;
208   PetscBool        iascii;
209 
210   PetscFunctionBegin;
211   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
212   if (iascii) {
213     ierr = PetscViewerASCIIPrintf(viewer,"Composite SNES type - %s\n",SNESCompositeTypes[jac->type]);CHKERRQ(ierr);
214     ierr = PetscViewerASCIIPrintf(viewer,"SNESes on composite preconditioner follow\n");CHKERRQ(ierr);
215     ierr = PetscViewerASCIIPrintf(viewer,"---------------------------------\n");CHKERRQ(ierr);
216   }
217   if (iascii) {
218     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
219   }
220   while (next) {
221     ierr = SNESView(next->snes,viewer);CHKERRQ(ierr);
222     next = next->next;
223   }
224   if (iascii) {
225     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
226     ierr = PetscViewerASCIIPrintf(viewer,"---------------------------------\n");CHKERRQ(ierr);
227   }
228   PetscFunctionReturn(0);
229 }
230 
231 /* ------------------------------------------------------------------------------*/
232 
233 #undef __FUNCT__
234 #define __FUNCT__ "SNESCompositeSetType_Composite"
235 static PetscErrorCode  SNESCompositeSetType_Composite(SNES snes,SNESCompositeType type)
236 {
237   SNES_Composite *jac = (SNES_Composite*)snes->data;
238 
239   PetscFunctionBegin;
240   jac->type = type;
241   PetscFunctionReturn(0);
242 }
243 
244 #undef __FUNCT__
245 #define __FUNCT__ "SNESCompositeAddSNES_Composite"
246 static PetscErrorCode  SNESCompositeAddSNES_Composite(SNES snes,SNESType type)
247 {
248   SNES_Composite     *jac;
249   SNES_CompositeLink next,ilink;
250   PetscErrorCode   ierr;
251   PetscInt         cnt = 0;
252   const char       *prefix;
253   char             newprefix[8];
254   DM               dm;
255 
256   PetscFunctionBegin;
257   ierr        = PetscNewLog(snes,struct _SNES_CompositeLink,&ilink);CHKERRQ(ierr);
258   ilink->next = 0;
259   ierr        = SNESCreate(PetscObjectComm((PetscObject)snes),&ilink->snes);CHKERRQ(ierr);
260   ierr        = PetscLogObjectParent((PetscObject)snes,(PetscObject)ilink->snes);CHKERRQ(ierr);
261   ierr        = SNESGetDM(snes,&dm);CHKERRQ(ierr);
262   ierr        = SNESSetDM(ilink->snes,dm);CHKERRQ(ierr);
263 
264   jac  = (SNES_Composite*)snes->data;
265   next = jac->head;
266   if (!next) {
267     jac->head       = ilink;
268     ilink->previous = NULL;
269   } else {
270     cnt++;
271     while (next->next) {
272       next = next->next;
273       cnt++;
274     }
275     next->next      = ilink;
276     ilink->previous = next;
277   }
278   ierr = SNESGetOptionsPrefix(snes,&prefix);CHKERRQ(ierr);
279   ierr = SNESSetOptionsPrefix(ilink->snes,prefix);CHKERRQ(ierr);
280   sprintf(newprefix,"sub_%d_",(int)cnt);
281   ierr = SNESAppendOptionsPrefix(ilink->snes,newprefix);CHKERRQ(ierr);
282   ierr = PetscObjectIncrementTabLevel((PetscObject)ilink->snes,(PetscObject)snes,1);CHKERRQ(ierr);
283   ierr = SNESSetType(ilink->snes,type);CHKERRQ(ierr);
284   ilink->dmp = 1.0;
285   PetscFunctionReturn(0);
286 }
287 
288 #undef __FUNCT__
289 #define __FUNCT__ "SNESCompositeGetSNES_Composite"
290 static PetscErrorCode  SNESCompositeGetSNES_Composite(SNES snes,PetscInt n,SNES *subsnes)
291 {
292   SNES_Composite     *jac;
293   SNES_CompositeLink next;
294   PetscInt         i;
295 
296   PetscFunctionBegin;
297   jac  = (SNES_Composite*)snes->data;
298   next = jac->head;
299   for (i=0; i<n; i++) {
300     if (!next->next) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_INCOMP,"Not enough SNESes in composite preconditioner");
301     next = next->next;
302   }
303   *subsnes = next->snes;
304   PetscFunctionReturn(0);
305 }
306 
307 /* -------------------------------------------------------------------------------- */
308 #undef __FUNCT__
309 #define __FUNCT__ "SNESCompositeSetType"
310 /*@C
311    SNESCompositeSetType - Sets the type of composite preconditioner.
312 
313    Logically Collective on SNES
314 
315    Input Parameter:
316 +  snes - the preconditioner context
317 -  type - SNES_COMPOSITE_ADDITIVE (default), SNES_COMPOSITE_MULTIPLICATIVE
318 
319    Options Database Key:
320 .  -snes_composite_type <type: one of multiplicative, additive, special> - Sets composite preconditioner type
321 
322    Level: Developer
323 
324 .keywords: SNES, set, type, composite preconditioner, additive, multiplicative
325 @*/
326 PetscErrorCode  SNESCompositeSetType(SNES snes,SNESCompositeType type)
327 {
328   PetscErrorCode ierr;
329 
330   PetscFunctionBegin;
331   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
332   PetscValidLogicalCollectiveEnum(snes,type,2);
333   ierr = PetscTryMethod(snes,"SNESCompositeSetType_C",(SNES,SNESCompositeType),(snes,type));CHKERRQ(ierr);
334   PetscFunctionReturn(0);
335 }
336 
337 #undef __FUNCT__
338 #define __FUNCT__ "SNESCompositeAddSNES"
339 /*@C
340    SNESCompositeAddSNES - Adds another SNES to the composite SNES.
341 
342    Collective on SNES
343 
344    Input Parameters:
345 +  snes - the preconditioner context
346 -  type - the type of the new preconditioner
347 
348    Level: Developer
349 
350 .keywords: SNES, composite preconditioner, add
351 @*/
352 PetscErrorCode  SNESCompositeAddSNES(SNES snes,SNESType type)
353 {
354   PetscErrorCode ierr;
355 
356   PetscFunctionBegin;
357   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
358   ierr = PetscTryMethod(snes,"SNESCompositeAddSNES_C",(SNES,SNESType),(snes,type));CHKERRQ(ierr);
359   PetscFunctionReturn(0);
360 }
361 #undef __FUNCT__
362 #define __FUNCT__ "SNESCompositeGetSNES"
363 /*@
364    SNESCompositeGetSNES - Gets one of the SNES objects in the composite SNES.
365 
366    Not Collective
367 
368    Input Parameter:
369 +  snes - the preconditioner context
370 -  n - the number of the snes requested
371 
372    Output Parameters:
373 .  subsnes - the SNES requested
374 
375    Level: Developer
376 
377 .keywords: SNES, get, composite preconditioner, sub preconditioner
378 
379 .seealso: SNESCompositeAddSNES()
380 @*/
381 PetscErrorCode  SNESCompositeGetSNES(SNES snes,PetscInt n,SNES *subsnes)
382 {
383   PetscErrorCode ierr;
384 
385   PetscFunctionBegin;
386   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
387   PetscValidPointer(subsnes,3);
388   ierr = PetscUseMethod(snes,"SNESCompositeGetSNES_C",(SNES,PetscInt,SNES*),(snes,n,subsnes));CHKERRQ(ierr);
389   PetscFunctionReturn(0);
390 }
391 
392 #undef __FUNCT__
393 #define __FUNCT__ "SNESCompositeSetDamping_Composite"
394 static PetscErrorCode  SNESCompositeSetDamping_Composite(SNES snes,PetscInt n,PetscReal dmp)
395 {
396   SNES_Composite     *jac;
397   SNES_CompositeLink next;
398   PetscInt         i;
399 
400   PetscFunctionBegin;
401   jac  = (SNES_Composite*)snes->data;
402   next = jac->head;
403   for (i=0; i<n; i++) {
404     if (!next->next) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_INCOMP,"Not enough SNESes in composite preconditioner");
405     next = next->next;
406   }
407   next->dmp = dmp;
408   PetscFunctionReturn(0);
409 }
410 
411 #undef __FUNCT__
412 #define __FUNCT__ "SNESCompositeSetDamping"
413 /*@
414    SNESCompositeSetDamping - Sets the damping of a subsolver when using additive composite SNES.
415 
416    Not Collective
417 
418    Input Parameter:
419 +  snes - the preconditioner context
420 .  n - the number of the snes requested
421 -  dmp - the damping
422 
423    Level: Developer
424 
425 .keywords: SNES, get, composite preconditioner, sub preconditioner
426 
427 .seealso: SNESCompositeAddSNES()
428 @*/
429 PetscErrorCode  SNESCompositeSetDamping(SNES snes,PetscInt n,PetscReal dmp)
430 {
431   PetscErrorCode ierr;
432 
433   PetscFunctionBegin;
434   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
435   ierr = PetscUseMethod(snes,"SNESCompositeSetDamping_C",(SNES,PetscInt,PetscReal),(snes,n,dmp));CHKERRQ(ierr);
436   PetscFunctionReturn(0);
437 }
438 
439 #undef __FUNCT__
440 #define __FUNCT__ "SNESSolve_Composite"
441 PetscErrorCode SNESSolve_Composite(SNES snes)
442 {
443   Vec            F;
444   Vec            X;
445   Vec            B;
446   PetscInt       i;
447   PetscReal      fnorm = 0.0;
448   PetscErrorCode ierr;
449   SNESNormType   normtype;
450   SNES_Composite *comp = (SNES_Composite*)snes->data;
451 
452   PetscFunctionBegin;
453   X = snes->vec_sol;
454   F = snes->vec_func;
455   B = snes->vec_rhs;
456 
457   ierr         = PetscObjectAMSTakeAccess((PetscObject)snes);CHKERRQ(ierr);
458   snes->iter   = 0;
459   snes->norm   = 0.;
460   ierr         = PetscObjectAMSGrantAccess((PetscObject)snes);CHKERRQ(ierr);
461   snes->reason = SNES_CONVERGED_ITERATING;
462   ierr         = SNESGetNormType(snes, &normtype);CHKERRQ(ierr);
463   if (normtype == SNES_NORM_FUNCTION || normtype == SNES_NORM_INITIAL_ONLY || normtype == SNES_NORM_INITIAL_FINAL_ONLY) {
464     if (!snes->vec_func_init_set) {
465       ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
466       if (snes->domainerror) {
467         snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
468         PetscFunctionReturn(0);
469       }
470     } else snes->vec_func_init_set = PETSC_FALSE;
471 
472     /* convergence test */
473     if (!snes->norm_init_set) {
474       ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F||  */
475       if (PetscIsInfOrNanReal(fnorm)) {
476         snes->reason = SNES_DIVERGED_FNORM_NAN;
477         PetscFunctionReturn(0);
478       }
479     } else {
480       fnorm               = snes->norm_init;
481       snes->norm_init_set = PETSC_FALSE;
482     }
483     ierr       = PetscObjectAMSTakeAccess((PetscObject)snes);CHKERRQ(ierr);
484     snes->iter = 0;
485     snes->norm = fnorm;
486     ierr       = PetscObjectAMSGrantAccess((PetscObject)snes);CHKERRQ(ierr);
487     ierr       = SNESLogConvergenceHistory(snes,snes->norm,0);CHKERRQ(ierr);
488     ierr       = SNESMonitor(snes,0,snes->norm);CHKERRQ(ierr);
489 
490     /* set parameter for default relative tolerance convergence test */
491     snes->ttol = fnorm*snes->rtol;
492 
493     /* test convergence */
494     ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
495     if (snes->reason) PetscFunctionReturn(0);
496   } else {
497     ierr = PetscObjectAMSGrantAccess((PetscObject)snes);CHKERRQ(ierr);
498     ierr = SNESLogConvergenceHistory(snes,snes->norm,0);CHKERRQ(ierr);
499     ierr = SNESMonitor(snes,0,snes->norm);CHKERRQ(ierr);
500   }
501 
502   /* Call general purpose update function */
503   if (snes->ops->update) {
504     ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
505   }
506 
507   for (i = 0; i < snes->max_its; i++) {
508     if (comp->type == SNES_COMPOSITE_ADDITIVE) {
509       ierr = SNESCompositeApply_Additive(snes,X,B,F,&fnorm);CHKERRQ(ierr);
510     } else if (comp->type == SNES_COMPOSITE_MULTIPLICATIVE) {
511       ierr = SNESCompositeApply_Multiplicative(snes,X,B,F,&fnorm);CHKERRQ(ierr);
512     } else {
513     }
514     if ((i == snes->max_its - 1) && (normtype == SNES_NORM_INITIAL_FINAL_ONLY || normtype == SNES_NORM_FINAL_ONLY)) {
515       ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
516       if (snes->domainerror) {
517         snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
518         break;
519       }
520       ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F||  */
521       if (PetscIsInfOrNanReal(fnorm)) {
522         snes->reason = SNES_DIVERGED_FNORM_NAN;
523         break;
524       }
525     }
526     /* Monitor convergence */
527     ierr       = PetscObjectAMSTakeAccess((PetscObject)snes);CHKERRQ(ierr);
528     snes->iter = i+1;
529     snes->norm = fnorm;
530     ierr       = PetscObjectAMSGrantAccess((PetscObject)snes);CHKERRQ(ierr);
531     ierr       = SNESLogConvergenceHistory(snes,snes->norm,0);CHKERRQ(ierr);
532     ierr       = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr);
533     /* Test for convergence */
534     if (normtype == SNES_NORM_FUNCTION) {ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);}
535     if (snes->reason) break;
536     /* Call general purpose update function */
537     if (snes->ops->update) {ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);}
538   }
539   if (normtype == SNES_NORM_FUNCTION) {
540     if (i == snes->max_its) {
541       ierr = PetscInfo1(snes,"Maximum number of iterations has been reached: %D\n",snes->max_its);CHKERRQ(ierr);
542       if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
543     }
544   } else if (!snes->reason) snes->reason = SNES_CONVERGED_ITS;
545   PetscFunctionReturn(0);
546 }
547 
548 /* -------------------------------------------------------------------------------------------*/
549 
550 /*MC
551      SNESCOMPOSITE - Build a preconditioner by composing together several nonlinear solvers
552 
553    Options Database Keys:
554 +  -snes_composite_type <type: one of multiplicative, additive, symmetric_multiplicative, special> - Sets composite preconditioner type
555 -  -snes_composite_sneses - <snes0,snes1,...> list of SNESes to compose
556 
557    Level: intermediate
558 
559    Concepts: composing solvers
560 
561 .seealso:  SNESCreate(), SNESSetType(), SNESType (for list of available types), SNES,
562            SNESSHELL, SNESCompositeSetType(), SNESCompositeSpecialSetAlpha(), SNESCompositeAddSNES(),
563            SNESCompositeGetSNES()
564 
565 M*/
566 
567 #undef __FUNCT__
568 #define __FUNCT__ "SNESCreate_Composite"
569 PETSC_EXTERN PetscErrorCode SNESCreate_Composite(SNES snes)
570 {
571   PetscErrorCode ierr;
572   SNES_Composite   *jac;
573 
574   PetscFunctionBegin;
575   ierr = PetscNewLog(snes,SNES_Composite,&jac);CHKERRQ(ierr);
576 
577   snes->ops->solve           = SNESSolve_Composite;
578   snes->ops->setup           = SNESSetUp_Composite;
579   snes->ops->reset           = SNESReset_Composite;
580   snes->ops->destroy         = SNESDestroy_Composite;
581   snes->ops->setfromoptions  = SNESSetFromOptions_Composite;
582   snes->ops->view            = SNESView_Composite;
583 
584   snes->data = (void*)jac;
585   jac->type  = SNES_COMPOSITE_ADDITIVE;
586   jac->head  = 0;
587 
588   ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESCompositeSetType_C",SNESCompositeSetType_Composite);CHKERRQ(ierr);
589   ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESCompositeAddSNES_C",SNESCompositeAddSNES_Composite);CHKERRQ(ierr);
590   ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESCompositeGetSNES_C",SNESCompositeGetSNES_Composite);CHKERRQ(ierr);
591   ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESCompositeSetDamping_C",SNESCompositeSetDamping_Composite);CHKERRQ(ierr);
592   PetscFunctionReturn(0);
593 }
594 
595