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