xref: /petsc/src/snes/impls/composite/snescomposite.c (revision dfd676b1a855b7f967ece75a22ee7f6626d10f89)
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 #include <petscblaslapack.h>
7 
8 const char *const        SNESCompositeTypes[]   = {"ADDITIVE","MULTIPLICATIVE","ADDITIVEOPTIMAL","SNESCompositeType","SNES_COMPOSITE",NULL};
9 
10 typedef struct _SNES_CompositeLink *SNES_CompositeLink;
11 struct _SNES_CompositeLink {
12   SNES               snes;
13   PetscReal          dmp;
14   Vec                X;
15   SNES_CompositeLink next;
16   SNES_CompositeLink previous;
17 };
18 
19 typedef struct {
20   SNES_CompositeLink head;
21   PetscInt           nsnes;
22   SNESCompositeType  type;
23   Vec                Xorig;
24   PetscInt           innerFailures; /* the number of inner failures we've seen */
25 
26   /* context for ADDITIVEOPTIMAL */
27   Vec                *Xes,*Fes;      /* solution and residual vectors for the subsolvers */
28   PetscReal          *fnorms;        /* norms of the solutions */
29   PetscScalar        *h;             /* the matrix formed as q_ij = (rdot_i, rdot_j) */
30   PetscScalar        *g;             /* the dotproducts of the previous function with the candidate functions */
31   PetscBLASInt       n;              /* matrix dimension -- nsnes */
32   PetscBLASInt       nrhs;           /* the number of right hand sides */
33   PetscBLASInt       lda;            /* the padded matrix dimension */
34   PetscBLASInt       ldb;            /* the padded vector dimension */
35   PetscReal          *s;             /* the singular values */
36   PetscScalar        *beta;          /* the RHS and combination */
37   PetscReal          rcond;          /* the exit condition */
38   PetscBLASInt       rank;           /* the effective rank */
39   PetscScalar        *work;          /* the work vector */
40   PetscReal          *rwork;         /* the real work vector used for complex */
41   PetscBLASInt       lwork;          /* the size of the work vector */
42   PetscBLASInt       info;           /* the output condition */
43 
44   PetscReal          rtol;           /* restart tolerance for accepting the combination */
45   PetscReal          stol;           /* restart tolerance for the combination */
46 } SNES_Composite;
47 
48 static PetscErrorCode SNESCompositeApply_Multiplicative(SNES snes,Vec X,Vec B,Vec F,PetscReal *fnorm)
49 {
50   SNES_Composite      *jac = (SNES_Composite*)snes->data;
51   SNES_CompositeLink  next = jac->head;
52   Vec                 FSub;
53   SNESConvergedReason reason;
54 
55   PetscFunctionBegin;
56   PetscCheck(next,PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"No composite SNESes supplied via SNESCompositeAddSNES() or -snes_composite_sneses");
57   if (snes->normschedule == SNES_NORM_ALWAYS) {
58     PetscCall(SNESSetInitialFunction(next->snes,F));
59   }
60   PetscCall(SNESSolve(next->snes,B,X));
61   PetscCall(SNESGetConvergedReason(next->snes,&reason));
62   if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
63     jac->innerFailures++;
64     if (jac->innerFailures >= snes->maxFailures) {
65       snes->reason = SNES_DIVERGED_INNER;
66       PetscFunctionReturn(0);
67     }
68   }
69 
70   while (next->next) {
71     /* only copy the function over in the case where the functions correspond */
72     if (next->snes->npcside== PC_RIGHT && next->snes->normschedule != SNES_NORM_NONE) {
73       PetscCall(SNESGetFunction(next->snes,&FSub,NULL,NULL));
74       next = next->next;
75       PetscCall(SNESSetInitialFunction(next->snes,FSub));
76     } else {
77       next = next->next;
78     }
79     PetscCall(SNESSolve(next->snes,B,X));
80     PetscCall(SNESGetConvergedReason(next->snes,&reason));
81     if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
82       jac->innerFailures++;
83       if (jac->innerFailures >= snes->maxFailures) {
84         snes->reason = SNES_DIVERGED_INNER;
85         PetscFunctionReturn(0);
86       }
87     }
88   }
89   if (next->snes->npcside== PC_RIGHT) {
90     PetscCall(SNESGetFunction(next->snes,&FSub,NULL,NULL));
91     PetscCall(VecCopy(FSub,F));
92     if (fnorm) {
93       if (snes->xl && snes->xu) {
94         PetscCall(SNESVIComputeInactiveSetFnorm(snes, F, X, fnorm));
95       } else {
96         PetscCall(VecNorm(F, NORM_2, fnorm));
97       }
98       SNESCheckFunctionNorm(snes,*fnorm);
99     }
100   } else if (snes->normschedule == SNES_NORM_ALWAYS) {
101     PetscCall(SNESComputeFunction(snes,X,F));
102     if (fnorm) {
103       if (snes->xl && snes->xu) {
104         PetscCall(SNESVIComputeInactiveSetFnorm(snes, F, X, fnorm));
105       } else {
106         PetscCall(VecNorm(F, NORM_2, fnorm));
107       }
108       SNESCheckFunctionNorm(snes,*fnorm);
109     }
110   }
111   PetscFunctionReturn(0);
112 }
113 
114 static PetscErrorCode SNESCompositeApply_Additive(SNES snes,Vec X,Vec B,Vec F,PetscReal *fnorm)
115 {
116   SNES_Composite      *jac = (SNES_Composite*)snes->data;
117   SNES_CompositeLink  next = jac->head;
118   Vec                 Y,Xorig;
119   SNESConvergedReason reason;
120 
121   PetscFunctionBegin;
122   Y = snes->vec_sol_update;
123   if (!jac->Xorig) PetscCall(VecDuplicate(X,&jac->Xorig));
124   Xorig = jac->Xorig;
125   PetscCall(VecCopy(X,Xorig));
126   PetscCheck(next,PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"No composite SNESes supplied via SNESCompositeAddSNES() or -snes_composite_sneses");
127   if (snes->normschedule == SNES_NORM_ALWAYS) {
128     PetscCall(SNESSetInitialFunction(next->snes,F));
129     while (next->next) {
130       next = next->next;
131       PetscCall(SNESSetInitialFunction(next->snes,F));
132     }
133   }
134   next = jac->head;
135   PetscCall(VecCopy(Xorig,Y));
136   PetscCall(SNESSolve(next->snes,B,Y));
137   PetscCall(SNESGetConvergedReason(next->snes,&reason));
138   if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
139     jac->innerFailures++;
140     if (jac->innerFailures >= snes->maxFailures) {
141       snes->reason = SNES_DIVERGED_INNER;
142       PetscFunctionReturn(0);
143     }
144   }
145   PetscCall(VecAXPY(Y,-1.0,Xorig));
146   PetscCall(VecAXPY(X,next->dmp,Y));
147   while (next->next) {
148     next = next->next;
149     PetscCall(VecCopy(Xorig,Y));
150     PetscCall(SNESSolve(next->snes,B,Y));
151     PetscCall(SNESGetConvergedReason(next->snes,&reason));
152     if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
153       jac->innerFailures++;
154       if (jac->innerFailures >= snes->maxFailures) {
155         snes->reason = SNES_DIVERGED_INNER;
156         PetscFunctionReturn(0);
157       }
158     }
159     PetscCall(VecAXPY(Y,-1.0,Xorig));
160     PetscCall(VecAXPY(X,next->dmp,Y));
161   }
162   if (snes->normschedule == SNES_NORM_ALWAYS) {
163     PetscCall(SNESComputeFunction(snes,X,F));
164     if (fnorm) {
165       if (snes->xl && snes->xu) {
166         PetscCall(SNESVIComputeInactiveSetFnorm(snes, F, X, fnorm));
167       } else {
168         PetscCall(VecNorm(F, NORM_2, fnorm));
169       }
170       SNESCheckFunctionNorm(snes,*fnorm);
171     }
172   }
173   PetscFunctionReturn(0);
174 }
175 
176 /* approximately solve the overdetermined system:
177 
178  2*F(x_i)\cdot F(\x_j)\alpha_i = 0
179  \alpha_i                      = 1
180 
181  Which minimizes the L2 norm of the linearization of:
182  ||F(\sum_i \alpha_i*x_i)||^2
183 
184  With the constraint that \sum_i\alpha_i = 1
185  Where x_i is the solution from the ith subsolver.
186  */
187 static PetscErrorCode SNESCompositeApply_AdditiveOptimal(SNES snes,Vec X,Vec B,Vec F,PetscReal *fnorm)
188 {
189   SNES_Composite      *jac = (SNES_Composite*)snes->data;
190   SNES_CompositeLink  next = jac->head;
191   Vec                 *Xes = jac->Xes,*Fes = jac->Fes;
192   PetscInt            i,j;
193   PetscScalar         tot,total,ftf;
194   PetscReal           min_fnorm;
195   PetscInt            min_i;
196   SNESConvergedReason reason;
197 
198   PetscFunctionBegin;
199   PetscCheck(next,PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"No composite SNESes supplied via SNESCompositeAddSNES() or -snes_composite_sneses");
200 
201   if (snes->normschedule == SNES_NORM_ALWAYS) {
202     next = jac->head;
203     PetscCall(SNESSetInitialFunction(next->snes,F));
204     while (next->next) {
205       next = next->next;
206       PetscCall(SNESSetInitialFunction(next->snes,F));
207     }
208   }
209 
210   next = jac->head;
211   i = 0;
212   PetscCall(VecCopy(X,Xes[i]));
213   PetscCall(SNESSolve(next->snes,B,Xes[i]));
214   PetscCall(SNESGetConvergedReason(next->snes,&reason));
215   if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
216     jac->innerFailures++;
217     if (jac->innerFailures >= snes->maxFailures) {
218       snes->reason = SNES_DIVERGED_INNER;
219       PetscFunctionReturn(0);
220     }
221   }
222   while (next->next) {
223     i++;
224     next = next->next;
225     PetscCall(VecCopy(X,Xes[i]));
226     PetscCall(SNESSolve(next->snes,B,Xes[i]));
227     PetscCall(SNESGetConvergedReason(next->snes,&reason));
228     if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
229       jac->innerFailures++;
230       if (jac->innerFailures >= snes->maxFailures) {
231         snes->reason = SNES_DIVERGED_INNER;
232         PetscFunctionReturn(0);
233       }
234     }
235   }
236 
237   /* all the solutions are collected; combine optimally */
238   for (i=0;i<jac->n;i++) {
239     for (j=0;j<i+1;j++) {
240       PetscCall(VecDotBegin(Fes[i],Fes[j],&jac->h[i + j*jac->n]));
241     }
242     PetscCall(VecDotBegin(Fes[i],F,&jac->g[i]));
243   }
244 
245   for (i=0;i<jac->n;i++) {
246     for (j=0;j<i+1;j++) {
247       PetscCall(VecDotEnd(Fes[i],Fes[j],&jac->h[i + j*jac->n]));
248       if (i == j) jac->fnorms[i] = PetscSqrtReal(PetscRealPart(jac->h[i + j*jac->n]));
249     }
250     PetscCall(VecDotEnd(Fes[i],F,&jac->g[i]));
251   }
252 
253   ftf = (*fnorm)*(*fnorm);
254 
255   for (i=0; i<jac->n; i++) {
256     for (j=i+1;j<jac->n;j++) {
257       jac->h[i + j*jac->n] = jac->h[j + i*jac->n];
258     }
259   }
260 
261   for (i=0; i<jac->n; i++) {
262     for (j=0; j<jac->n; j++) {
263       jac->h[i + j*jac->n] = jac->h[i + j*jac->n] - jac->g[j] - jac->g[i] + ftf;
264     }
265     jac->beta[i] = ftf - jac->g[i];
266   }
267 
268   jac->info  = 0;
269   jac->rcond = -1.;
270   PetscCall(PetscFPTrapPush(PETSC_FP_TRAP_OFF));
271 #if defined(PETSC_USE_COMPLEX)
272   PetscStackCallBLAS("LAPACKgelss",LAPACKgelss_(&jac->n,&jac->n,&jac->nrhs,jac->h,&jac->lda,jac->beta,&jac->lda,jac->s,&jac->rcond,&jac->rank,jac->work,&jac->lwork,jac->rwork,&jac->info));
273 #else
274   PetscStackCallBLAS("LAPACKgelss",LAPACKgelss_(&jac->n,&jac->n,&jac->nrhs,jac->h,&jac->lda,jac->beta,&jac->lda,jac->s,&jac->rcond,&jac->rank,jac->work,&jac->lwork,&jac->info));
275 #endif
276   PetscCall(PetscFPTrapPop());
277   PetscCheck(jac->info >= 0,PetscObjectComm((PetscObject)snes),PETSC_ERR_LIB,"Bad argument to GELSS");
278   PetscCheck(jac->info <= 0,PetscObjectComm((PetscObject)snes),PETSC_ERR_LIB,"SVD failed to converge");
279   tot = 0.;
280   total = 0.;
281   for (i=0; i<jac->n; i++) {
282     PetscCheck(!snes->errorifnotconverged || !PetscIsInfOrNanScalar(jac->beta[i]),PetscObjectComm((PetscObject)snes),PETSC_ERR_LIB,"SVD generated inconsistent output");
283     PetscCall(PetscInfo(snes,"%" PetscInt_FMT ": %g\n",i,(double)PetscRealPart(jac->beta[i])));
284     tot += jac->beta[i];
285     total += PetscAbsScalar(jac->beta[i]);
286   }
287   PetscCall(VecScale(X,(1. - tot)));
288   PetscCall(VecMAXPY(X,jac->n,jac->beta,Xes));
289   PetscCall(SNESComputeFunction(snes,X,F));
290 
291   if (snes->xl && snes->xu) {
292     PetscCall(SNESVIComputeInactiveSetFnorm(snes, F, X, fnorm));
293   } else {
294     PetscCall(VecNorm(F, NORM_2, fnorm));
295   }
296 
297   /* take the minimum-normed candidate if it beats the combination by a factor of rtol or the combination has stagnated */
298   min_fnorm = jac->fnorms[0];
299   min_i     = 0;
300   for (i=0; i<jac->n; i++) {
301     if (jac->fnorms[i] < min_fnorm) {
302       min_fnorm = jac->fnorms[i];
303       min_i     = i;
304     }
305   }
306 
307   /* stagnation or divergence restart to the solution of the solver that failed the least */
308   if (PetscRealPart(total) < jac->stol || min_fnorm*jac->rtol < *fnorm) {
309     PetscCall(VecCopy(jac->Xes[min_i],X));
310     PetscCall(VecCopy(jac->Fes[min_i],F));
311     *fnorm = min_fnorm;
312   }
313   PetscFunctionReturn(0);
314 }
315 
316 static PetscErrorCode SNESSetUp_Composite(SNES snes)
317 {
318   DM                 dm;
319   SNES_Composite     *jac = (SNES_Composite*)snes->data;
320   SNES_CompositeLink next = jac->head;
321   PetscInt           n=0,i;
322   Vec                F;
323 
324   PetscFunctionBegin;
325   PetscCall(SNESGetDM(snes,&dm));
326 
327   if (snes->ops->computevariablebounds) {
328     /* SNESVI only ever calls computevariablebounds once, so calling it once here is justified */
329     if (!snes->xl) PetscCall(VecDuplicate(snes->vec_sol,&snes->xl));
330     if (!snes->xu) PetscCall(VecDuplicate(snes->vec_sol,&snes->xu));
331     PetscCall((*snes->ops->computevariablebounds)(snes,snes->xl,snes->xu));
332   }
333 
334   while (next) {
335     n++;
336     PetscCall(SNESSetDM(next->snes,dm));
337     PetscCall(SNESSetApplicationContext(next->snes, snes->user));
338     if (snes->xl && snes->xu) {
339       if (snes->ops->computevariablebounds) {
340         PetscCall(SNESVISetComputeVariableBounds(next->snes, snes->ops->computevariablebounds));
341       } else {
342         PetscCall(SNESVISetVariableBounds(next->snes,snes->xl,snes->xu));
343       }
344     }
345 
346     next = next->next;
347   }
348   jac->nsnes = n;
349   PetscCall(SNESGetFunction(snes,&F,NULL,NULL));
350   if (jac->type == SNES_COMPOSITE_ADDITIVEOPTIMAL) {
351     PetscCall(VecDuplicateVecs(F,jac->nsnes,&jac->Xes));
352     PetscCall(PetscMalloc1(n,&jac->Fes));
353     PetscCall(PetscMalloc1(n,&jac->fnorms));
354     next = jac->head;
355     i = 0;
356     while (next) {
357       PetscCall(SNESGetFunction(next->snes,&F,NULL,NULL));
358       jac->Fes[i] = F;
359       PetscCall(PetscObjectReference((PetscObject)F));
360       next = next->next;
361       i++;
362     }
363     /* allocate the subspace direct solve area */
364     jac->nrhs  = 1;
365     jac->lda   = jac->nsnes;
366     jac->ldb   = jac->nsnes;
367     jac->n     = jac->nsnes;
368 
369     PetscCall(PetscMalloc1(jac->n*jac->n,&jac->h));
370     PetscCall(PetscMalloc1(jac->n,&jac->beta));
371     PetscCall(PetscMalloc1(jac->n,&jac->s));
372     PetscCall(PetscMalloc1(jac->n,&jac->g));
373     jac->lwork = 12*jac->n;
374 #if defined(PETSC_USE_COMPLEX)
375     PetscCall(PetscMalloc1(jac->lwork,&jac->rwork));
376 #endif
377     PetscCall(PetscMalloc1(jac->lwork,&jac->work));
378   }
379 
380   PetscFunctionReturn(0);
381 }
382 
383 static PetscErrorCode SNESReset_Composite(SNES snes)
384 {
385   SNES_Composite     *jac = (SNES_Composite*)snes->data;
386   SNES_CompositeLink next = jac->head;
387 
388   PetscFunctionBegin;
389   while (next) {
390     PetscCall(SNESReset(next->snes));
391     next = next->next;
392   }
393   PetscCall(VecDestroy(&jac->Xorig));
394   if (jac->Xes) PetscCall(VecDestroyVecs(jac->nsnes,&jac->Xes));
395   if (jac->Fes) PetscCall(VecDestroyVecs(jac->nsnes,&jac->Fes));
396   PetscCall(PetscFree(jac->fnorms));
397   PetscCall(PetscFree(jac->h));
398   PetscCall(PetscFree(jac->s));
399   PetscCall(PetscFree(jac->g));
400   PetscCall(PetscFree(jac->beta));
401   PetscCall(PetscFree(jac->work));
402   PetscCall(PetscFree(jac->rwork));
403   PetscFunctionReturn(0);
404 }
405 
406 static PetscErrorCode SNESDestroy_Composite(SNES snes)
407 {
408   SNES_Composite     *jac = (SNES_Composite*)snes->data;
409   SNES_CompositeLink next = jac->head,next_tmp;
410 
411   PetscFunctionBegin;
412   PetscCall(SNESReset_Composite(snes));
413   while (next) {
414     PetscCall(SNESDestroy(&next->snes));
415     next_tmp = next;
416     next     = next->next;
417     PetscCall(PetscFree(next_tmp));
418   }
419   PetscCall(PetscFree(snes->data));
420   PetscFunctionReturn(0);
421 }
422 
423 static PetscErrorCode SNESSetFromOptions_Composite(PetscOptionItems *PetscOptionsObject,SNES snes)
424 {
425   SNES_Composite     *jac = (SNES_Composite*)snes->data;
426   PetscInt           nmax = 8,i;
427   SNES_CompositeLink next;
428   char               *sneses[8];
429   PetscReal          dmps[8];
430   PetscBool          flg;
431 
432   PetscFunctionBegin;
433   PetscOptionsHeadBegin(PetscOptionsObject,"Composite preconditioner options");
434   PetscCall(PetscOptionsEnum("-snes_composite_type","Type of composition","SNESCompositeSetType",SNESCompositeTypes,(PetscEnum)jac->type,(PetscEnum*)&jac->type,&flg));
435   if (flg) {
436     PetscCall(SNESCompositeSetType(snes,jac->type));
437   }
438   PetscCall(PetscOptionsStringArray("-snes_composite_sneses","List of composite solvers","SNESCompositeAddSNES",sneses,&nmax,&flg));
439   if (flg) {
440     for (i=0; i<nmax; i++) {
441       PetscCall(SNESCompositeAddSNES(snes,sneses[i]));
442       PetscCall(PetscFree(sneses[i]));   /* deallocate string sneses[i], which is allocated in PetscOptionsStringArray() */
443     }
444   }
445   PetscCall(PetscOptionsRealArray("-snes_composite_damping","Damping of the additive composite solvers","SNESCompositeSetDamping",dmps,&nmax,&flg));
446   if (flg) {
447     for (i=0; i<nmax; i++) {
448       PetscCall(SNESCompositeSetDamping(snes,i,dmps[i]));
449     }
450   }
451   PetscCall(PetscOptionsReal("-snes_composite_stol","Step tolerance for restart on the additive composite solvers","",jac->stol,&jac->stol,NULL));
452   PetscCall(PetscOptionsReal("-snes_composite_rtol","Residual tolerance for the additive composite solvers","",jac->rtol,&jac->rtol,NULL));
453   PetscOptionsHeadEnd();
454 
455   next = jac->head;
456   while (next) {
457     PetscCall(SNESSetFromOptions(next->snes));
458     next = next->next;
459   }
460   PetscFunctionReturn(0);
461 }
462 
463 static PetscErrorCode SNESView_Composite(SNES snes,PetscViewer viewer)
464 {
465   SNES_Composite     *jac = (SNES_Composite*)snes->data;
466   SNES_CompositeLink next = jac->head;
467   PetscBool          iascii;
468 
469   PetscFunctionBegin;
470   PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii));
471   if (iascii) {
472     PetscCall(PetscViewerASCIIPrintf(viewer,"  type - %s\n",SNESCompositeTypes[jac->type]));
473     PetscCall(PetscViewerASCIIPrintf(viewer,"  SNESes on composite preconditioner follow\n"));
474     PetscCall(PetscViewerASCIIPrintf(viewer,"  ---------------------------------\n"));
475   }
476   if (iascii) {
477     PetscCall(PetscViewerASCIIPushTab(viewer));
478   }
479   while (next) {
480     PetscCall(SNESView(next->snes,viewer));
481     next = next->next;
482   }
483   if (iascii) {
484     PetscCall(PetscViewerASCIIPopTab(viewer));
485     PetscCall(PetscViewerASCIIPrintf(viewer,"  ---------------------------------\n"));
486   }
487   PetscFunctionReturn(0);
488 }
489 
490 /* ------------------------------------------------------------------------------*/
491 
492 static PetscErrorCode  SNESCompositeSetType_Composite(SNES snes,SNESCompositeType type)
493 {
494   SNES_Composite *jac = (SNES_Composite*)snes->data;
495 
496   PetscFunctionBegin;
497   jac->type = type;
498   PetscFunctionReturn(0);
499 }
500 
501 static PetscErrorCode  SNESCompositeAddSNES_Composite(SNES snes,SNESType type)
502 {
503   SNES_Composite     *jac;
504   SNES_CompositeLink next,ilink;
505   PetscInt           cnt = 0;
506   const char         *prefix;
507   char               newprefix[20];
508   DM                 dm;
509 
510   PetscFunctionBegin;
511   PetscCall(PetscNewLog(snes,&ilink));
512   ilink->next = NULL;
513   PetscCall(SNESCreate(PetscObjectComm((PetscObject)snes),&ilink->snes));
514   PetscCall(PetscLogObjectParent((PetscObject)snes,(PetscObject)ilink->snes));
515   PetscCall(SNESGetDM(snes,&dm));
516   PetscCall(SNESSetDM(ilink->snes,dm));
517   PetscCall(SNESSetTolerances(ilink->snes,snes->abstol,snes->rtol,snes->stol,1,snes->max_funcs));
518   PetscCall(PetscObjectCopyFortranFunctionPointers((PetscObject)snes,(PetscObject)ilink->snes));
519   jac  = (SNES_Composite*)snes->data;
520   next = jac->head;
521   if (!next) {
522     jac->head       = ilink;
523     ilink->previous = NULL;
524   } else {
525     cnt++;
526     while (next->next) {
527       next = next->next;
528       cnt++;
529     }
530     next->next      = ilink;
531     ilink->previous = next;
532   }
533   PetscCall(SNESGetOptionsPrefix(snes,&prefix));
534   PetscCall(SNESSetOptionsPrefix(ilink->snes,prefix));
535   PetscCall(PetscSNPrintf(newprefix,sizeof(newprefix),"sub_%d_",(int)cnt));
536   PetscCall(SNESAppendOptionsPrefix(ilink->snes,newprefix));
537   PetscCall(PetscObjectIncrementTabLevel((PetscObject)ilink->snes,(PetscObject)snes,1));
538   PetscCall(SNESSetType(ilink->snes,type));
539   PetscCall(SNESSetNormSchedule(ilink->snes, SNES_NORM_FINAL_ONLY));
540 
541   ilink->dmp = 1.0;
542   jac->nsnes++;
543   PetscFunctionReturn(0);
544 }
545 
546 static PetscErrorCode  SNESCompositeGetSNES_Composite(SNES snes,PetscInt n,SNES *subsnes)
547 {
548   SNES_Composite     *jac;
549   SNES_CompositeLink next;
550   PetscInt           i;
551 
552   PetscFunctionBegin;
553   jac  = (SNES_Composite*)snes->data;
554   next = jac->head;
555   for (i=0; i<n; i++) {
556     PetscCheck(next->next,PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_INCOMP,"Not enough SNESes in composite preconditioner");
557     next = next->next;
558   }
559   *subsnes = next->snes;
560   PetscFunctionReturn(0);
561 }
562 
563 /* -------------------------------------------------------------------------------- */
564 /*@C
565    SNESCompositeSetType - Sets the type of composite preconditioner.
566 
567    Logically Collective on SNES
568 
569    Input Parameters:
570 +  snes - the preconditioner context
571 -  type - SNES_COMPOSITE_ADDITIVE (default), SNES_COMPOSITE_MULTIPLICATIVE
572 
573    Options Database Key:
574 .  -snes_composite_type <type: one of multiplicative, additive, special> - Sets composite preconditioner type
575 
576    Level: Developer
577 
578 @*/
579 PetscErrorCode  SNESCompositeSetType(SNES snes,SNESCompositeType type)
580 {
581   PetscFunctionBegin;
582   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
583   PetscValidLogicalCollectiveEnum(snes,type,2);
584   PetscTryMethod(snes,"SNESCompositeSetType_C",(SNES,SNESCompositeType),(snes,type));
585   PetscFunctionReturn(0);
586 }
587 
588 /*@C
589    SNESCompositeAddSNES - Adds another SNES to the composite SNES.
590 
591    Collective on SNES
592 
593    Input Parameters:
594 +  snes - the preconditioner context
595 -  type - the type of the new preconditioner
596 
597    Level: Developer
598 
599 @*/
600 PetscErrorCode  SNESCompositeAddSNES(SNES snes,SNESType type)
601 {
602   PetscFunctionBegin;
603   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
604   PetscTryMethod(snes,"SNESCompositeAddSNES_C",(SNES,SNESType),(snes,type));
605   PetscFunctionReturn(0);
606 }
607 
608 /*@
609    SNESCompositeGetSNES - Gets one of the SNES objects in the composite SNES.
610 
611    Not Collective
612 
613    Input Parameters:
614 +  snes - the preconditioner context
615 -  n - the number of the snes requested
616 
617    Output Parameters:
618 .  subsnes - the SNES requested
619 
620    Level: Developer
621 
622 .seealso: SNESCompositeAddSNES()
623 @*/
624 PetscErrorCode  SNESCompositeGetSNES(SNES snes,PetscInt n,SNES *subsnes)
625 {
626   PetscFunctionBegin;
627   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
628   PetscValidPointer(subsnes,3);
629   PetscUseMethod(snes,"SNESCompositeGetSNES_C",(SNES,PetscInt,SNES*),(snes,n,subsnes));
630   PetscFunctionReturn(0);
631 }
632 
633 /*@
634    SNESCompositeGetNumber - Get the number of subsolvers in the composite SNES.
635 
636    Logically Collective on SNES
637 
638    Input Parameter:
639    snes - the preconditioner context
640 
641    Output Parameter:
642    n - the number of subsolvers
643 
644    Level: Developer
645 
646 @*/
647 PetscErrorCode  SNESCompositeGetNumber(SNES snes,PetscInt *n)
648 {
649   SNES_Composite     *jac;
650   SNES_CompositeLink next;
651 
652   PetscFunctionBegin;
653   jac  = (SNES_Composite*)snes->data;
654   next = jac->head;
655 
656   *n = 0;
657   while (next) {
658     *n = *n + 1;
659     next = next->next;
660   }
661   PetscFunctionReturn(0);
662 }
663 
664 static PetscErrorCode  SNESCompositeSetDamping_Composite(SNES snes,PetscInt n,PetscReal dmp)
665 {
666   SNES_Composite     *jac;
667   SNES_CompositeLink next;
668   PetscInt           i;
669 
670   PetscFunctionBegin;
671   jac  = (SNES_Composite*)snes->data;
672   next = jac->head;
673   for (i=0; i<n; i++) {
674     PetscCheck(next->next,PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_INCOMP,"Not enough SNESes in composite preconditioner");
675     next = next->next;
676   }
677   next->dmp = dmp;
678   PetscFunctionReturn(0);
679 }
680 
681 /*@
682    SNESCompositeSetDamping - Sets the damping of a subsolver when using additive composite SNES.
683 
684    Not Collective
685 
686    Input Parameters:
687 +  snes - the preconditioner context
688 .  n - the number of the snes requested
689 -  dmp - the damping
690 
691    Level: Developer
692 
693 .seealso: SNESCompositeAddSNES()
694 @*/
695 PetscErrorCode  SNESCompositeSetDamping(SNES snes,PetscInt n,PetscReal dmp)
696 {
697   PetscFunctionBegin;
698   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
699   PetscUseMethod(snes,"SNESCompositeSetDamping_C",(SNES,PetscInt,PetscReal),(snes,n,dmp));
700   PetscFunctionReturn(0);
701 }
702 
703 static PetscErrorCode SNESSolve_Composite(SNES snes)
704 {
705   Vec              F,X,B,Y;
706   PetscInt         i;
707   PetscReal        fnorm = 0.0, xnorm = 0.0, snorm = 0.0;
708   SNESNormSchedule normtype;
709   SNES_Composite   *comp = (SNES_Composite*)snes->data;
710 
711   PetscFunctionBegin;
712   X = snes->vec_sol;
713   F = snes->vec_func;
714   B = snes->vec_rhs;
715   Y = snes->vec_sol_update;
716 
717   PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
718   snes->iter   = 0;
719   snes->norm   = 0.;
720   comp->innerFailures = 0;
721   PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
722   snes->reason = SNES_CONVERGED_ITERATING;
723   PetscCall(SNESGetNormSchedule(snes, &normtype));
724   if (normtype == SNES_NORM_ALWAYS || normtype == SNES_NORM_INITIAL_ONLY || normtype == SNES_NORM_INITIAL_FINAL_ONLY) {
725     if (!snes->vec_func_init_set) {
726       PetscCall(SNESComputeFunction(snes,X,F));
727     } else snes->vec_func_init_set = PETSC_FALSE;
728 
729     if (snes->xl && snes->xu) {
730       PetscCall(SNESVIComputeInactiveSetFnorm(snes, F, X, &fnorm));
731     } else {
732       PetscCall(VecNorm(F, NORM_2, &fnorm)); /* fnorm <- ||F||  */
733     }
734     SNESCheckFunctionNorm(snes,fnorm);
735     PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
736     snes->iter = 0;
737     snes->norm = fnorm;
738     PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
739     PetscCall(SNESLogConvergenceHistory(snes,snes->norm,0));
740     PetscCall(SNESMonitor(snes,0,snes->norm));
741 
742     /* test convergence */
743     PetscCall((*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP));
744     if (snes->reason) PetscFunctionReturn(0);
745   } else {
746     PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
747     PetscCall(SNESLogConvergenceHistory(snes,snes->norm,0));
748     PetscCall(SNESMonitor(snes,0,snes->norm));
749   }
750 
751   for (i = 0; i < snes->max_its; i++) {
752     /* Call general purpose update function */
753     if (snes->ops->update) {
754       PetscCall((*snes->ops->update)(snes, snes->iter));
755     }
756 
757     /* Copy the state before modification by application of the composite solver;
758        we will subtract the new state after application */
759     PetscCall(VecCopy(X, Y));
760 
761     if (comp->type == SNES_COMPOSITE_ADDITIVE) {
762       PetscCall(SNESCompositeApply_Additive(snes,X,B,F,&fnorm));
763     } else if (comp->type == SNES_COMPOSITE_MULTIPLICATIVE) {
764       PetscCall(SNESCompositeApply_Multiplicative(snes,X,B,F,&fnorm));
765     } else if (comp->type == SNES_COMPOSITE_ADDITIVEOPTIMAL) {
766       PetscCall(SNESCompositeApply_AdditiveOptimal(snes,X,B,F,&fnorm));
767     } else SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"Unsupported SNESComposite type");
768     if (snes->reason < 0) break;
769 
770     /* Compute the solution update for convergence testing */
771     PetscCall(VecAYPX(Y, -1.0, X));
772 
773     if ((i == snes->max_its - 1) && (normtype == SNES_NORM_INITIAL_FINAL_ONLY || normtype == SNES_NORM_FINAL_ONLY)) {
774       PetscCall(SNESComputeFunction(snes,X,F));
775 
776       if (snes->xl && snes->xu) {
777         PetscCall(VecNormBegin(X, NORM_2, &xnorm));
778         PetscCall(VecNormBegin(Y, NORM_2, &snorm));
779         PetscCall(SNESVIComputeInactiveSetFnorm(snes, F, X, &fnorm));
780         PetscCall(VecNormEnd(X, NORM_2, &xnorm));
781         PetscCall(VecNormEnd(Y, NORM_2, &snorm));
782       } else {
783         PetscCall(VecNormBegin(F, NORM_2, &fnorm));
784         PetscCall(VecNormBegin(X, NORM_2, &xnorm));
785         PetscCall(VecNormBegin(Y, NORM_2, &snorm));
786 
787         PetscCall(VecNormEnd(F, NORM_2, &fnorm));
788         PetscCall(VecNormEnd(X, NORM_2, &xnorm));
789         PetscCall(VecNormEnd(Y, NORM_2, &snorm));
790       }
791       SNESCheckFunctionNorm(snes,fnorm);
792     } else if (normtype == SNES_NORM_ALWAYS) {
793       PetscCall(VecNormBegin(X, NORM_2, &xnorm));
794       PetscCall(VecNormBegin(Y, NORM_2, &snorm));
795       PetscCall(VecNormEnd(X, NORM_2, &xnorm));
796       PetscCall(VecNormEnd(Y, NORM_2, &snorm));
797     }
798     /* Monitor convergence */
799     PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
800     snes->iter = i+1;
801     snes->norm = fnorm;
802     snes->xnorm = xnorm;
803     snes->ynorm = snorm;
804     PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
805     PetscCall(SNESLogConvergenceHistory(snes,snes->norm,0));
806     PetscCall(SNESMonitor(snes,snes->iter,snes->norm));
807     /* Test for convergence */
808     if (normtype == SNES_NORM_ALWAYS) PetscCall((*snes->ops->converged)(snes,snes->iter,xnorm,snorm,fnorm,&snes->reason,snes->cnvP));
809     if (snes->reason) break;
810   }
811   if (normtype == SNES_NORM_ALWAYS) {
812     if (i == snes->max_its) {
813       PetscCall(PetscInfo(snes,"Maximum number of iterations has been reached: %" PetscInt_FMT "\n",snes->max_its));
814       if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
815     }
816   } else if (!snes->reason) snes->reason = SNES_CONVERGED_ITS;
817   PetscFunctionReturn(0);
818 }
819 
820 /* -------------------------------------------------------------------------------------------*/
821 
822 /*MC
823      SNESCOMPOSITE - Build a preconditioner by composing together several nonlinear solvers
824 
825    Options Database Keys:
826 +  -snes_composite_type <type: one of multiplicative, additive, symmetric_multiplicative, special> - Sets composite preconditioner type
827 -  -snes_composite_sneses - <snes0,snes1,...> list of SNESes to compose
828 
829    Level: intermediate
830 
831 .seealso:  SNESCreate(), SNESSetType(), SNESType (for list of available types), SNES,
832            SNESSHELL, SNESCompositeSetType(), SNESCompositeSpecialSetAlpha(), SNESCompositeAddSNES(),
833            SNESCompositeGetSNES()
834 
835    References:
836 .  * - Peter R. Brune, Matthew G. Knepley, Barry F. Smith, and Xuemin Tu, "Composing Scalable Nonlinear Algebraic Solvers",
837    SIAM Review, 57(4), 2015
838 
839 M*/
840 
841 PETSC_EXTERN PetscErrorCode SNESCreate_Composite(SNES snes)
842 {
843   SNES_Composite *jac;
844 
845   PetscFunctionBegin;
846   PetscCall(PetscNewLog(snes,&jac));
847 
848   snes->ops->solve           = SNESSolve_Composite;
849   snes->ops->setup           = SNESSetUp_Composite;
850   snes->ops->reset           = SNESReset_Composite;
851   snes->ops->destroy         = SNESDestroy_Composite;
852   snes->ops->setfromoptions  = SNESSetFromOptions_Composite;
853   snes->ops->view            = SNESView_Composite;
854 
855   snes->usesksp        = PETSC_FALSE;
856 
857   snes->alwayscomputesfinalresidual = PETSC_FALSE;
858 
859   snes->data = (void*)jac;
860   jac->type  = SNES_COMPOSITE_ADDITIVEOPTIMAL;
861   jac->Fes   = NULL;
862   jac->Xes   = NULL;
863   jac->fnorms = NULL;
864   jac->nsnes = 0;
865   jac->head  = NULL;
866   jac->stol  = 0.1;
867   jac->rtol  = 1.1;
868 
869   jac->h     = NULL;
870   jac->s     = NULL;
871   jac->beta  = NULL;
872   jac->work  = NULL;
873   jac->rwork = NULL;
874 
875   PetscCall(PetscObjectComposeFunction((PetscObject)snes,"SNESCompositeSetType_C",SNESCompositeSetType_Composite));
876   PetscCall(PetscObjectComposeFunction((PetscObject)snes,"SNESCompositeAddSNES_C",SNESCompositeAddSNES_Composite));
877   PetscCall(PetscObjectComposeFunction((PetscObject)snes,"SNESCompositeGetSNES_C",SNESCompositeGetSNES_Composite));
878   PetscCall(PetscObjectComposeFunction((PetscObject)snes,"SNESCompositeSetDamping_C",SNESCompositeSetDamping_Composite));
879   PetscFunctionReturn(0);
880 }
881