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