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