xref: /petsc/src/snes/impls/composite/snescomposite.c (revision a13bc4e38d03c664bfce6f01234033d0aae2418e)
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   SNES_Composite     *jac = (SNES_Composite*)snes->data;
253   SNES_CompositeLink next = jac->head;
254   PetscInt           n=0,i;
255   Vec                F;
256 
257   PetscFunctionBegin;
258   while (next) {
259     n++;
260     ierr = SNESSetFromOptions(next->snes);CHKERRQ(ierr);
261     next = next->next;
262   }
263   jac->nsnes = n;
264   ierr = SNESGetFunction(snes,&F,NULL,NULL);CHKERRQ(ierr);
265   if (jac->type == SNES_COMPOSITE_ADDITIVEOPTIMAL) {
266     ierr = VecDuplicateVecs(F,jac->nsnes,&jac->Xes);CHKERRQ(ierr);
267     ierr = PetscMalloc(sizeof(Vec)*n,&jac->Fes);CHKERRQ(ierr);
268     ierr = PetscMalloc(sizeof(PetscReal)*n,&jac->fnorms);CHKERRQ(ierr);
269     next = jac->head;
270     i = 0;
271     while (next) {
272       ierr = SNESGetFunction(next->snes,&F,NULL,NULL);CHKERRQ(ierr);
273       jac->Fes[i] = F;
274       ierr = PetscObjectReference((PetscObject)F);CHKERRQ(ierr);
275       next = next->next;
276       i++;
277     }
278     /* allocate the subspace direct solve area */
279     jac->nrhs  = 1;
280     jac->lda   = jac->nsnes;
281     jac->ldb   = jac->nsnes;
282     jac->n     = jac->nsnes;
283 
284     ierr = PetscMalloc(jac->n*jac->n*sizeof(PetscScalar),&jac->h);CHKERRQ(ierr);
285     ierr = PetscMalloc(jac->n*sizeof(PetscScalar),&jac->beta);CHKERRQ(ierr);
286     ierr = PetscMalloc(jac->n*sizeof(PetscScalar),&jac->s);CHKERRQ(ierr);
287     ierr = PetscMalloc(jac->n*sizeof(PetscScalar),&jac->g);CHKERRQ(ierr);
288     jac->lwork = 12*jac->n;
289 #if PETSC_USE_COMPLEX
290     ierr = PetscMalloc(sizeof(PetscReal)*jac->lwork,&jac->rwork);CHKERRQ(ierr);
291 #endif
292     ierr = PetscMalloc(sizeof(PetscScalar)*jac->lwork,&jac->work);CHKERRQ(ierr);
293   }
294 
295   PetscFunctionReturn(0);
296 }
297 
298 #undef __FUNCT__
299 #define __FUNCT__ "SNESReset_Composite"
300 static PetscErrorCode SNESReset_Composite(SNES snes)
301 {
302   SNES_Composite     *jac = (SNES_Composite*)snes->data;
303   PetscErrorCode   ierr;
304   SNES_CompositeLink next = jac->head;
305 
306   PetscFunctionBegin;
307   while (next) {
308     ierr = SNESReset(next->snes);CHKERRQ(ierr);
309     next = next->next;
310   }
311   ierr = VecDestroy(&jac->Xorig);CHKERRQ(ierr);
312   if (jac->Xes) {ierr = VecDestroyVecs(jac->nsnes,&jac->Xes);CHKERRQ(ierr);}
313   if (jac->Fes) {ierr = VecDestroyVecs(jac->nsnes,&jac->Fes);CHKERRQ(ierr);}
314   ierr = PetscFree(jac->fnorms);CHKERRQ(ierr);
315   ierr = PetscFree(jac->h);CHKERRQ(ierr);
316   ierr = PetscFree(jac->s);CHKERRQ(ierr);
317   ierr = PetscFree(jac->g);CHKERRQ(ierr);
318   ierr = PetscFree(jac->beta);CHKERRQ(ierr);
319   ierr = PetscFree(jac->work);CHKERRQ(ierr);
320   ierr = PetscFree(jac->rwork);CHKERRQ(ierr);
321 
322   PetscFunctionReturn(0);
323 }
324 
325 #undef __FUNCT__
326 #define __FUNCT__ "SNESDestroy_Composite"
327 static PetscErrorCode SNESDestroy_Composite(SNES snes)
328 {
329   SNES_Composite     *jac = (SNES_Composite*)snes->data;
330   PetscErrorCode   ierr;
331   SNES_CompositeLink next = jac->head,next_tmp;
332 
333   PetscFunctionBegin;
334   ierr = SNESReset_Composite(snes);CHKERRQ(ierr);
335   while (next) {
336     ierr     = SNESDestroy(&next->snes);CHKERRQ(ierr);
337     next_tmp = next;
338     next     = next->next;
339     ierr     = PetscFree(next_tmp);CHKERRQ(ierr);
340   }
341   ierr = PetscFree(snes->data);CHKERRQ(ierr);
342   PetscFunctionReturn(0);
343 }
344 
345 #undef __FUNCT__
346 #define __FUNCT__ "SNESSetFromOptions_Composite"
347 static PetscErrorCode SNESSetFromOptions_Composite(SNES snes)
348 {
349   SNES_Composite     *jac = (SNES_Composite*)snes->data;
350   PetscErrorCode     ierr;
351   PetscInt           nmax = 8,i;
352   SNES_CompositeLink next;
353   char               *sneses[8];
354   PetscReal          dmps[8];
355   PetscBool          flg;
356 
357   PetscFunctionBegin;
358   ierr = PetscOptionsHead("Composite preconditioner options");CHKERRQ(ierr);
359   ierr = PetscOptionsEnum("-snes_composite_type","Type of composition","SNESCompositeSetType",SNESCompositeTypes,(PetscEnum)jac->type,(PetscEnum*)&jac->type,&flg);CHKERRQ(ierr);
360   if (flg) {
361     ierr = SNESCompositeSetType(snes,jac->type);CHKERRQ(ierr);
362   }
363   ierr = PetscOptionsStringArray("-snes_composite_sneses","List of composite solvers","SNESCompositeAddSNES",sneses,&nmax,&flg);CHKERRQ(ierr);
364   if (flg) {
365     for (i=0; i<nmax; i++) {
366       ierr = SNESCompositeAddSNES(snes,sneses[i]);CHKERRQ(ierr);
367       ierr = PetscFree(sneses[i]);CHKERRQ(ierr);   /* deallocate string sneses[i], which is allocated in PetscOptionsStringArray() */
368     }
369   }
370   ierr = PetscOptionsRealArray("-snes_composite_damping","Damping of the additive composite solvers","SNESCompositeSetDamping",dmps,&nmax,&flg);CHKERRQ(ierr);
371   if (flg) {
372     for (i=0; i<nmax; i++) {
373       ierr = SNESCompositeSetDamping(snes,i,dmps[i]);CHKERRQ(ierr);
374     }
375   }
376   ierr = PetscOptionsReal("-snes_composite_stol","Step tolerance for restart on the additive composite solvers","",jac->stol,&jac->stol,NULL);CHKERRQ(ierr);
377   ierr = PetscOptionsReal("-snes_composite_rtol","Residual tolerance for the additive composite solvers","",jac->rtol,&jac->rtol,NULL);CHKERRQ(ierr);
378   ierr = PetscOptionsTail();CHKERRQ(ierr);
379 
380   next = jac->head;
381   while (next) {
382     ierr = SNESSetFromOptions(next->snes);CHKERRQ(ierr);
383     next = next->next;
384   }
385   PetscFunctionReturn(0);
386 }
387 
388 #undef __FUNCT__
389 #define __FUNCT__ "SNESView_Composite"
390 static PetscErrorCode SNESView_Composite(SNES snes,PetscViewer viewer)
391 {
392   SNES_Composite     *jac = (SNES_Composite*)snes->data;
393   PetscErrorCode   ierr;
394   SNES_CompositeLink next = jac->head;
395   PetscBool        iascii;
396 
397   PetscFunctionBegin;
398   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
399   if (iascii) {
400     ierr = PetscViewerASCIIPrintf(viewer,"Composite SNES type - %s\n",SNESCompositeTypes[jac->type]);CHKERRQ(ierr);
401     ierr = PetscViewerASCIIPrintf(viewer,"SNESes on composite preconditioner follow\n");CHKERRQ(ierr);
402     ierr = PetscViewerASCIIPrintf(viewer,"---------------------------------\n");CHKERRQ(ierr);
403   }
404   if (iascii) {
405     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
406   }
407   while (next) {
408     ierr = SNESView(next->snes,viewer);CHKERRQ(ierr);
409     next = next->next;
410   }
411   if (iascii) {
412     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
413     ierr = PetscViewerASCIIPrintf(viewer,"---------------------------------\n");CHKERRQ(ierr);
414   }
415   PetscFunctionReturn(0);
416 }
417 
418 /* ------------------------------------------------------------------------------*/
419 
420 #undef __FUNCT__
421 #define __FUNCT__ "SNESCompositeSetType_Composite"
422 static PetscErrorCode  SNESCompositeSetType_Composite(SNES snes,SNESCompositeType type)
423 {
424   SNES_Composite *jac = (SNES_Composite*)snes->data;
425 
426   PetscFunctionBegin;
427   jac->type = type;
428   PetscFunctionReturn(0);
429 }
430 
431 #undef __FUNCT__
432 #define __FUNCT__ "SNESCompositeAddSNES_Composite"
433 static PetscErrorCode  SNESCompositeAddSNES_Composite(SNES snes,SNESType type)
434 {
435   SNES_Composite     *jac;
436   SNES_CompositeLink next,ilink;
437   PetscErrorCode   ierr;
438   PetscInt         cnt = 0;
439   const char       *prefix;
440   char             newprefix[8];
441   DM               dm;
442 
443   PetscFunctionBegin;
444   ierr        = PetscNewLog(snes,struct _SNES_CompositeLink,&ilink);CHKERRQ(ierr);
445   ilink->next = 0;
446   ierr        = SNESCreate(PetscObjectComm((PetscObject)snes),&ilink->snes);CHKERRQ(ierr);
447   ierr        = PetscLogObjectParent((PetscObject)snes,(PetscObject)ilink->snes);CHKERRQ(ierr);
448   ierr        = SNESGetDM(snes,&dm);CHKERRQ(ierr);
449   ierr        = SNESSetDM(ilink->snes,dm);CHKERRQ(ierr);
450   ierr        = SNESSetTolerances(ilink->snes,ilink->snes->abstol,ilink->snes->rtol,ilink->snes->stol,1,ilink->snes->max_funcs);CHKERRQ(ierr);
451   jac  = (SNES_Composite*)snes->data;
452   next = jac->head;
453   if (!next) {
454     jac->head       = ilink;
455     ilink->previous = NULL;
456   } else {
457     cnt++;
458     while (next->next) {
459       next = next->next;
460       cnt++;
461     }
462     next->next      = ilink;
463     ilink->previous = next;
464   }
465   ierr = SNESGetOptionsPrefix(snes,&prefix);CHKERRQ(ierr);
466   ierr = SNESSetOptionsPrefix(ilink->snes,prefix);CHKERRQ(ierr);
467   sprintf(newprefix,"sub_%d_",(int)cnt);
468   ierr = SNESAppendOptionsPrefix(ilink->snes,newprefix);CHKERRQ(ierr);
469   ierr = PetscObjectIncrementTabLevel((PetscObject)ilink->snes,(PetscObject)snes,1);CHKERRQ(ierr);
470   ierr = SNESSetType(ilink->snes,type);CHKERRQ(ierr);
471   ierr = SNESSetNormSchedule(ilink->snes, SNES_NORM_FINAL_ONLY);CHKERRQ(ierr);
472   ilink->dmp = 1.0;
473   jac->nsnes++;
474   PetscFunctionReturn(0);
475 }
476 
477 #undef __FUNCT__
478 #define __FUNCT__ "SNESCompositeGetSNES_Composite"
479 static PetscErrorCode  SNESCompositeGetSNES_Composite(SNES snes,PetscInt n,SNES *subsnes)
480 {
481   SNES_Composite     *jac;
482   SNES_CompositeLink next;
483   PetscInt         i;
484 
485   PetscFunctionBegin;
486   jac  = (SNES_Composite*)snes->data;
487   next = jac->head;
488   for (i=0; i<n; i++) {
489     if (!next->next) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_INCOMP,"Not enough SNESes in composite preconditioner");
490     next = next->next;
491   }
492   *subsnes = next->snes;
493   PetscFunctionReturn(0);
494 }
495 
496 /* -------------------------------------------------------------------------------- */
497 #undef __FUNCT__
498 #define __FUNCT__ "SNESCompositeSetType"
499 /*@C
500    SNESCompositeSetType - Sets the type of composite preconditioner.
501 
502    Logically Collective on SNES
503 
504    Input Parameter:
505 +  snes - the preconditioner context
506 -  type - SNES_COMPOSITE_ADDITIVE (default), SNES_COMPOSITE_MULTIPLICATIVE
507 
508    Options Database Key:
509 .  -snes_composite_type <type: one of multiplicative, additive, special> - Sets composite preconditioner type
510 
511    Level: Developer
512 
513 .keywords: SNES, set, type, composite preconditioner, additive, multiplicative
514 @*/
515 PetscErrorCode  SNESCompositeSetType(SNES snes,SNESCompositeType type)
516 {
517   PetscErrorCode ierr;
518 
519   PetscFunctionBegin;
520   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
521   PetscValidLogicalCollectiveEnum(snes,type,2);
522   ierr = PetscTryMethod(snes,"SNESCompositeSetType_C",(SNES,SNESCompositeType),(snes,type));CHKERRQ(ierr);
523   PetscFunctionReturn(0);
524 }
525 
526 #undef __FUNCT__
527 #define __FUNCT__ "SNESCompositeAddSNES"
528 /*@C
529    SNESCompositeAddSNES - Adds another SNES to the composite SNES.
530 
531    Collective on SNES
532 
533    Input Parameters:
534 +  snes - the preconditioner context
535 -  type - the type of the new preconditioner
536 
537    Level: Developer
538 
539 .keywords: SNES, composite preconditioner, add
540 @*/
541 PetscErrorCode  SNESCompositeAddSNES(SNES snes,SNESType type)
542 {
543   PetscErrorCode ierr;
544 
545   PetscFunctionBegin;
546   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
547   ierr = PetscTryMethod(snes,"SNESCompositeAddSNES_C",(SNES,SNESType),(snes,type));CHKERRQ(ierr);
548   PetscFunctionReturn(0);
549 }
550 #undef __FUNCT__
551 #define __FUNCT__ "SNESCompositeGetSNES"
552 /*@
553    SNESCompositeGetSNES - Gets one of the SNES objects in the composite SNES.
554 
555    Not Collective
556 
557    Input Parameter:
558 +  snes - the preconditioner context
559 -  n - the number of the snes requested
560 
561    Output Parameters:
562 .  subsnes - the SNES requested
563 
564    Level: Developer
565 
566 .keywords: SNES, get, composite preconditioner, sub preconditioner
567 
568 .seealso: SNESCompositeAddSNES()
569 @*/
570 PetscErrorCode  SNESCompositeGetSNES(SNES snes,PetscInt n,SNES *subsnes)
571 {
572   PetscErrorCode ierr;
573 
574   PetscFunctionBegin;
575   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
576   PetscValidPointer(subsnes,3);
577   ierr = PetscUseMethod(snes,"SNESCompositeGetSNES_C",(SNES,PetscInt,SNES*),(snes,n,subsnes));CHKERRQ(ierr);
578   PetscFunctionReturn(0);
579 }
580 
581 #undef __FUNCT__
582 #define __FUNCT__ "SNESCompositeSetDamping_Composite"
583 static PetscErrorCode  SNESCompositeSetDamping_Composite(SNES snes,PetscInt n,PetscReal dmp)
584 {
585   SNES_Composite     *jac;
586   SNES_CompositeLink next;
587   PetscInt         i;
588 
589   PetscFunctionBegin;
590   jac  = (SNES_Composite*)snes->data;
591   next = jac->head;
592   for (i=0; i<n; i++) {
593     if (!next->next) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_INCOMP,"Not enough SNESes in composite preconditioner");
594     next = next->next;
595   }
596   next->dmp = dmp;
597   PetscFunctionReturn(0);
598 }
599 
600 #undef __FUNCT__
601 #define __FUNCT__ "SNESCompositeSetDamping"
602 /*@
603    SNESCompositeSetDamping - Sets the damping of a subsolver when using additive composite SNES.
604 
605    Not Collective
606 
607    Input Parameter:
608 +  snes - the preconditioner context
609 .  n - the number of the snes requested
610 -  dmp - the damping
611 
612    Level: Developer
613 
614 .keywords: SNES, get, composite preconditioner, sub preconditioner
615 
616 .seealso: SNESCompositeAddSNES()
617 @*/
618 PetscErrorCode  SNESCompositeSetDamping(SNES snes,PetscInt n,PetscReal dmp)
619 {
620   PetscErrorCode ierr;
621 
622   PetscFunctionBegin;
623   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
624   ierr = PetscUseMethod(snes,"SNESCompositeSetDamping_C",(SNES,PetscInt,PetscReal),(snes,n,dmp));CHKERRQ(ierr);
625   PetscFunctionReturn(0);
626 }
627 
628 #undef __FUNCT__
629 #define __FUNCT__ "SNESSolve_Composite"
630 PetscErrorCode SNESSolve_Composite(SNES snes)
631 {
632   Vec            F;
633   Vec            X;
634   Vec            B;
635   PetscInt       i;
636   PetscReal      fnorm = 0.0;
637   PetscErrorCode ierr;
638   SNESNormSchedule normtype;
639   SNES_Composite *comp = (SNES_Composite*)snes->data;
640 
641   PetscFunctionBegin;
642   X = snes->vec_sol;
643   F = snes->vec_func;
644   B = snes->vec_rhs;
645 
646   ierr         = PetscObjectAMSTakeAccess((PetscObject)snes);CHKERRQ(ierr);
647   snes->iter   = 0;
648   snes->norm   = 0.;
649   ierr         = PetscObjectAMSGrantAccess((PetscObject)snes);CHKERRQ(ierr);
650   snes->reason = SNES_CONVERGED_ITERATING;
651   ierr         = SNESGetNormSchedule(snes, &normtype);CHKERRQ(ierr);
652   if (normtype == SNES_NORM_ALWAYS || normtype == SNES_NORM_INITIAL_ONLY || normtype == SNES_NORM_INITIAL_FINAL_ONLY) {
653     if (!snes->vec_func_init_set) {
654       ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
655       if (snes->domainerror) {
656         snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
657         PetscFunctionReturn(0);
658       }
659     } else snes->vec_func_init_set = PETSC_FALSE;
660 
661     /* convergence test */
662     if (!snes->norm_init_set) {
663       ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F||  */
664       if (PetscIsInfOrNanReal(fnorm)) {
665         snes->reason = SNES_DIVERGED_FNORM_NAN;
666         PetscFunctionReturn(0);
667       }
668     } else {
669       fnorm               = snes->norm_init;
670       snes->norm_init_set = PETSC_FALSE;
671     }
672     ierr       = PetscObjectAMSTakeAccess((PetscObject)snes);CHKERRQ(ierr);
673     snes->iter = 0;
674     snes->norm = fnorm;
675     ierr       = PetscObjectAMSGrantAccess((PetscObject)snes);CHKERRQ(ierr);
676     ierr       = SNESLogConvergenceHistory(snes,snes->norm,0);CHKERRQ(ierr);
677     ierr       = SNESMonitor(snes,0,snes->norm);CHKERRQ(ierr);
678 
679     /* set parameter for default relative tolerance convergence test */
680     snes->ttol = fnorm*snes->rtol;
681 
682     /* test convergence */
683     ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
684     if (snes->reason) PetscFunctionReturn(0);
685   } else {
686     ierr = PetscObjectAMSGrantAccess((PetscObject)snes);CHKERRQ(ierr);
687     ierr = SNESLogConvergenceHistory(snes,snes->norm,0);CHKERRQ(ierr);
688     ierr = SNESMonitor(snes,0,snes->norm);CHKERRQ(ierr);
689   }
690 
691   /* Call general purpose update function */
692   if (snes->ops->update) {
693     ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
694   }
695 
696   for (i = 0; i < snes->max_its; i++) {
697     if (comp->type == SNES_COMPOSITE_ADDITIVE) {
698       ierr = SNESCompositeApply_Additive(snes,X,B,F,&fnorm);CHKERRQ(ierr);
699     } else if (comp->type == SNES_COMPOSITE_MULTIPLICATIVE) {
700       ierr = SNESCompositeApply_Multiplicative(snes,X,B,F,&fnorm);CHKERRQ(ierr);
701     } else if (comp->type == SNES_COMPOSITE_ADDITIVEOPTIMAL) {
702       ierr = SNESCompositeApply_AdditiveOptimal(snes,X,B,F,&fnorm);CHKERRQ(ierr);
703     } else {
704       SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"Unsupported SNESComposite type");
705     }
706     if ((i == snes->max_its - 1) && (normtype == SNES_NORM_INITIAL_FINAL_ONLY || normtype == SNES_NORM_FINAL_ONLY)) {
707       ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
708       if (snes->domainerror) {
709         snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
710         break;
711       }
712       ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F||  */
713       if (PetscIsInfOrNanReal(fnorm)) {
714         snes->reason = SNES_DIVERGED_FNORM_NAN;
715         break;
716       }
717     }
718     /* Monitor convergence */
719     ierr       = PetscObjectAMSTakeAccess((PetscObject)snes);CHKERRQ(ierr);
720     snes->iter = i+1;
721     snes->norm = fnorm;
722     ierr       = PetscObjectAMSGrantAccess((PetscObject)snes);CHKERRQ(ierr);
723     ierr       = SNESLogConvergenceHistory(snes,snes->norm,0);CHKERRQ(ierr);
724     ierr       = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr);
725     /* Test for convergence */
726     if (normtype == SNES_NORM_ALWAYS) {ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);}
727     if (snes->reason) break;
728     /* Call general purpose update function */
729     if (snes->ops->update) {ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);}
730   }
731   if (normtype == SNES_NORM_ALWAYS) {
732     if (i == snes->max_its) {
733       ierr = PetscInfo1(snes,"Maximum number of iterations has been reached: %D\n",snes->max_its);CHKERRQ(ierr);
734       if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
735     }
736   } else if (!snes->reason) snes->reason = SNES_CONVERGED_ITS;
737   PetscFunctionReturn(0);
738 }
739 
740 /* -------------------------------------------------------------------------------------------*/
741 
742 /*MC
743      SNESCOMPOSITE - Build a preconditioner by composing together several nonlinear solvers
744 
745    Options Database Keys:
746 +  -snes_composite_type <type: one of multiplicative, additive, symmetric_multiplicative, special> - Sets composite preconditioner type
747 -  -snes_composite_sneses - <snes0,snes1,...> list of SNESes to compose
748 
749    Level: intermediate
750 
751    Concepts: composing solvers
752 
753 .seealso:  SNESCreate(), SNESSetType(), SNESType (for list of available types), SNES,
754            SNESSHELL, SNESCompositeSetType(), SNESCompositeSpecialSetAlpha(), SNESCompositeAddSNES(),
755            SNESCompositeGetSNES()
756 
757 M*/
758 
759 #undef __FUNCT__
760 #define __FUNCT__ "SNESCreate_Composite"
761 PETSC_EXTERN PetscErrorCode SNESCreate_Composite(SNES snes)
762 {
763   PetscErrorCode ierr;
764   SNES_Composite   *jac;
765 
766   PetscFunctionBegin;
767   ierr = PetscNewLog(snes,SNES_Composite,&jac);CHKERRQ(ierr);
768 
769   snes->ops->solve           = SNESSolve_Composite;
770   snes->ops->setup           = SNESSetUp_Composite;
771   snes->ops->reset           = SNESReset_Composite;
772   snes->ops->destroy         = SNESDestroy_Composite;
773   snes->ops->setfromoptions  = SNESSetFromOptions_Composite;
774   snes->ops->view            = SNESView_Composite;
775 
776   snes->data = (void*)jac;
777   jac->type  = SNES_COMPOSITE_ADDITIVEOPTIMAL;
778   jac->Fes   = NULL;
779   jac->Xes   = NULL;
780   jac->fnorms = NULL;
781   jac->nsnes = 0;
782   jac->head  = 0;
783   jac->stol  = 0.1;
784   jac->rtol  = 1.1;
785 
786   jac->h     = NULL;
787   jac->s     = NULL;
788   jac->beta  = NULL;
789   jac->work  = NULL;
790   jac->rwork = NULL;
791 
792   ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESCompositeSetType_C",SNESCompositeSetType_Composite);CHKERRQ(ierr);
793   ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESCompositeAddSNES_C",SNESCompositeAddSNES_Composite);CHKERRQ(ierr);
794   ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESCompositeGetSNES_C",SNESCompositeGetSNES_Composite);CHKERRQ(ierr);
795   ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESCompositeSetDamping_C",SNESCompositeSetDamping_Composite);CHKERRQ(ierr);
796   PetscFunctionReturn(0);
797 }
798 
799