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