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