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