xref: /petsc/src/snes/impls/composite/snescomposite.c (revision eed5f15bc2a3dc26d54d7c2989531114c489d821)
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 
7 const char *const        SNESCompositeTypes[]   = {"ADDITIVE","MULTIPLICATIVE","SNESCompositeType","SNES_COMPOSITE",0};
8 
9 typedef struct _SNES_CompositeLink *SNES_CompositeLink;
10 struct _SNES_CompositeLink {
11   SNES               snes;
12   SNES_CompositeLink next;
13   SNES_CompositeLink previous;
14 };
15 
16 typedef struct {
17   SNES_CompositeLink head;
18   SNESCompositeType  type;
19   Vec                Xorig;
20 } SNES_Composite;
21 
22 #undef __FUNCT__
23 #define __FUNCT__ "SNESCompositeApply_Multiplicative"
24 static PetscErrorCode SNESCompositeApply_Multiplicative(SNES snes,Vec X,Vec B,Vec F,PetscReal *fnorm)
25 {
26   PetscErrorCode     ierr;
27   SNES_Composite     *jac = (SNES_Composite*)snes->data;
28   SNES_CompositeLink next = jac->head;
29   Vec                FSub;
30   PetscReal          fsubnorm;
31 
32   PetscFunctionBegin;
33   if (!next) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"No composite SNESes supplied via SNESCompositeAddSNES() or -snes_composite_sneses");
34   if (snes->normtype == SNES_NORM_FUNCTION) {
35     ierr = SNESSetInitialFunction(next->snes,F);CHKERRQ(ierr);
36     if (fnorm) {ierr = SNESSetInitialFunctionNorm(next->snes,*fnorm);CHKERRQ(ierr);}
37   }
38   ierr = SNESSolve(next->snes,B,X);CHKERRQ(ierr);
39 
40   while (next->next) {
41     /* only copy the function over in the case where the functions correspond */
42     if (next->snes->pcside == PC_RIGHT && next->snes->normtype != SNES_NORM_NONE) {
43       ierr = SNESGetFunction(next->snes,&FSub,NULL,NULL);CHKERRQ(ierr);
44       ierr = SNESGetFunctionNorm(next->snes,&fsubnorm);CHKERRQ(ierr);
45       next = next->next;
46       ierr = SNESSetInitialFunction(next->snes,FSub);CHKERRQ(ierr);
47       ierr = SNESSetInitialFunctionNorm(next->snes,fsubnorm);CHKERRQ(ierr);
48     } else {
49       next = next->next;
50     }
51     ierr = SNESSolve(next->snes,B,X);CHKERRQ(ierr);
52   }
53   if (next->snes->pcside == PC_RIGHT) {
54     ierr = SNESGetFunction(next->snes,&FSub,NULL,NULL);CHKERRQ(ierr);
55     ierr = VecCopy(FSub,F);CHKERRQ(ierr);
56     if (fnorm) {ierr = SNESGetFunctionNorm(next->snes,fnorm);CHKERRQ(ierr);}
57   } else if (snes->normtype == SNES_NORM_FUNCTION) {
58     SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
59     if (fnorm) {ierr = VecNorm(F,NORM_2,fnorm);CHKERRQ(ierr);}
60   }
61   PetscFunctionReturn(0);
62 }
63 
64 #undef __FUNCT__
65 #define __FUNCT__ "SNESCompositeApply_Additive"
66 static PetscErrorCode SNESCompositeApply_Additive(SNES snes,Vec X,Vec B,Vec F,PetscReal *fnorm)
67 {
68   PetscErrorCode     ierr;
69   SNES_Composite     *jac = (SNES_Composite*)snes->data;
70   SNES_CompositeLink next = jac->head;
71   Vec                Y,Xorig;
72 
73   PetscFunctionBegin;
74   Y = snes->vec_sol_update;
75   if (!jac->Xorig) {ierr = VecDuplicate(X,&jac->Xorig);CHKERRQ(ierr);}
76   Xorig = jac->Xorig;
77   ierr = VecCopy(X,Xorig);
78   if (!next) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"No composite SNESes supplied via SNESCompositeAddSNES() or -snes_composite_sneses");
79   ierr = SNESSolve(next->snes,B,X);CHKERRQ(ierr);
80   if (snes->normtype == SNES_NORM_FUNCTION) {
81     ierr = SNESSetInitialFunction(next->snes,F);CHKERRQ(ierr);
82     if (fnorm) {ierr = SNESSetInitialFunctionNorm(next->snes,*fnorm);CHKERRQ(ierr);}
83     while (next->next) {
84       next = next->next;
85       ierr = SNESSetInitialFunction(next->snes,F);CHKERRQ(ierr);
86       if (fnorm) {ierr = SNESSetInitialFunctionNorm(next->snes,*fnorm);CHKERRQ(ierr);}
87     }
88   }
89   next = jac->head;
90   while (next->next) {
91     next = next->next;
92     ierr = VecCopy(X,Y);CHKERRQ(ierr);
93     ierr = SNESSolve(next->snes,B,Y);CHKERRQ(ierr);
94     ierr = VecAXPY(Y,-1.0,Xorig);CHKERRQ(ierr);
95     ierr = VecAXPY(Y,1.0,X);CHKERRQ(ierr);
96   }
97   if (snes->normtype == SNES_NORM_FUNCTION) {
98     ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
99     if (fnorm) {ierr = VecNorm(F,NORM_2,fnorm);CHKERRQ(ierr);}
100   }
101   PetscFunctionReturn(0);
102 }
103 #undef __FUNCT__
104 #define __FUNCT__ "SNESSetUp_Composite"
105 static PetscErrorCode SNESSetUp_Composite(SNES snes)
106 {
107   PetscErrorCode   ierr;
108   SNES_Composite     *jac = (SNES_Composite*)snes->data;
109   SNES_CompositeLink next = jac->head;
110 
111   PetscFunctionBegin;
112   while (next) {
113     ierr = SNESSetFromOptions(next->snes);CHKERRQ(ierr);
114     next = next->next;
115   }
116   PetscFunctionReturn(0);
117 }
118 
119 #undef __FUNCT__
120 #define __FUNCT__ "SNESReset_Composite"
121 static PetscErrorCode SNESReset_Composite(SNES snes)
122 {
123   SNES_Composite     *jac = (SNES_Composite*)snes->data;
124   PetscErrorCode   ierr;
125   SNES_CompositeLink next = jac->head;
126 
127   PetscFunctionBegin;
128   while (next) {
129     ierr = SNESReset(next->snes);CHKERRQ(ierr);
130     next = next->next;
131   }
132   ierr = VecDestroy(&jac->Xorig);CHKERRQ(ierr);
133   PetscFunctionReturn(0);
134 }
135 
136 #undef __FUNCT__
137 #define __FUNCT__ "SNESDestroy_Composite"
138 static PetscErrorCode SNESDestroy_Composite(SNES snes)
139 {
140   SNES_Composite     *jac = (SNES_Composite*)snes->data;
141   PetscErrorCode   ierr;
142   SNES_CompositeLink next = jac->head,next_tmp;
143 
144   PetscFunctionBegin;
145   ierr = SNESReset_Composite(snes);CHKERRQ(ierr);
146   while (next) {
147     ierr     = SNESDestroy(&next->snes);CHKERRQ(ierr);
148     next_tmp = next;
149     next     = next->next;
150     ierr     = PetscFree(next_tmp);CHKERRQ(ierr);
151   }
152   ierr = PetscFree(snes->data);CHKERRQ(ierr);
153   PetscFunctionReturn(0);
154 }
155 
156 #undef __FUNCT__
157 #define __FUNCT__ "SNESSetFromOptions_Composite"
158 static PetscErrorCode SNESSetFromOptions_Composite(SNES snes)
159 {
160   SNES_Composite     *jac = (SNES_Composite*)snes->data;
161   PetscErrorCode     ierr;
162   PetscInt           nmax = 8,i;
163   SNES_CompositeLink next;
164   char               *sneses[8];
165   PetscBool          flg;
166 
167   PetscFunctionBegin;
168   ierr = PetscOptionsHead("Composite preconditioner options");CHKERRQ(ierr);
169   ierr = PetscOptionsEnum("-snes_composite_type","Type of composition","SNESCompositeSetType",SNESCompositeTypes,(PetscEnum)jac->type,(PetscEnum*)&jac->type,&flg);CHKERRQ(ierr);
170   if (flg) {
171     ierr = SNESCompositeSetType(snes,jac->type);CHKERRQ(ierr);
172   }
173   ierr = PetscOptionsStringArray("-snes_composite_sneses","List of composite solvers","SNESCompositeAddSNES",sneses,&nmax,&flg);CHKERRQ(ierr);
174   if (flg) {
175     for (i=0; i<nmax; i++) {
176       ierr = SNESCompositeAddSNES(snes,sneses[i]);CHKERRQ(ierr);
177       ierr = PetscFree(sneses[i]);CHKERRQ(ierr);   /* deallocate string sneses[i], which is allocated in PetscOptionsStringArray() */
178     }
179   }
180   ierr = PetscOptionsTail();CHKERRQ(ierr);
181 
182   next = jac->head;
183   while (next) {
184     ierr = SNESSetFromOptions(next->snes);CHKERRQ(ierr);
185     next = next->next;
186   }
187   PetscFunctionReturn(0);
188 }
189 
190 #undef __FUNCT__
191 #define __FUNCT__ "SNESView_Composite"
192 static PetscErrorCode SNESView_Composite(SNES snes,PetscViewer viewer)
193 {
194   SNES_Composite     *jac = (SNES_Composite*)snes->data;
195   PetscErrorCode   ierr;
196   SNES_CompositeLink next = jac->head;
197   PetscBool        iascii;
198 
199   PetscFunctionBegin;
200   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
201   if (iascii) {
202     ierr = PetscViewerASCIIPrintf(viewer,"Composite SNES type - %s\n",SNESCompositeTypes[jac->type]);CHKERRQ(ierr);
203     ierr = PetscViewerASCIIPrintf(viewer,"SNESes on composite preconditioner follow\n");CHKERRQ(ierr);
204     ierr = PetscViewerASCIIPrintf(viewer,"---------------------------------\n");CHKERRQ(ierr);
205   }
206   if (iascii) {
207     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
208   }
209   while (next) {
210     ierr = SNESView(next->snes,viewer);CHKERRQ(ierr);
211     next = next->next;
212   }
213   if (iascii) {
214     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
215     ierr = PetscViewerASCIIPrintf(viewer,"---------------------------------\n");CHKERRQ(ierr);
216   }
217   PetscFunctionReturn(0);
218 }
219 
220 /* ------------------------------------------------------------------------------*/
221 
222 #undef __FUNCT__
223 #define __FUNCT__ "SNESCompositeSetType_Composite"
224 static PetscErrorCode  SNESCompositeSetType_Composite(SNES snes,SNESCompositeType type)
225 {
226   SNES_Composite *jac = (SNES_Composite*)snes->data;
227 
228   PetscFunctionBegin;
229   jac->type = type;
230   PetscFunctionReturn(0);
231 }
232 
233 #undef __FUNCT__
234 #define __FUNCT__ "SNESCompositeAddSNES_Composite"
235 static PetscErrorCode  SNESCompositeAddSNES_Composite(SNES snes,SNESType type)
236 {
237   SNES_Composite     *jac;
238   SNES_CompositeLink next,ilink;
239   PetscErrorCode   ierr;
240   PetscInt         cnt = 0;
241   const char       *prefix;
242   char             newprefix[8];
243   DM               dm;
244 
245   PetscFunctionBegin;
246   ierr        = PetscNewLog(snes,struct _SNES_CompositeLink,&ilink);CHKERRQ(ierr);
247   ilink->next = 0;
248   ierr        = SNESCreate(PetscObjectComm((PetscObject)snes),&ilink->snes);CHKERRQ(ierr);
249   ierr        = PetscLogObjectParent((PetscObject)snes,(PetscObject)ilink->snes);CHKERRQ(ierr);
250   ierr        = SNESGetDM(snes,&dm);CHKERRQ(ierr);
251   ierr        = SNESSetDM(ilink->snes,dm);CHKERRQ(ierr);
252 
253   jac  = (SNES_Composite*)snes->data;
254   next = jac->head;
255   if (!next) {
256     jac->head       = ilink;
257     ilink->previous = NULL;
258   } else {
259     cnt++;
260     while (next->next) {
261       next = next->next;
262       cnt++;
263     }
264     next->next      = ilink;
265     ilink->previous = next;
266   }
267   ierr = SNESGetOptionsPrefix(snes,&prefix);CHKERRQ(ierr);
268   ierr = SNESSetOptionsPrefix(ilink->snes,prefix);CHKERRQ(ierr);
269   sprintf(newprefix,"sub_%d_",(int)cnt);
270   ierr = SNESAppendOptionsPrefix(ilink->snes,newprefix);CHKERRQ(ierr);
271   ierr = SNESSetType(ilink->snes,type);CHKERRQ(ierr);
272   PetscFunctionReturn(0);
273 }
274 
275 #undef __FUNCT__
276 #define __FUNCT__ "SNESCompositeGetSNES_Composite"
277 static PetscErrorCode  SNESCompositeGetSNES_Composite(SNES snes,PetscInt n,SNES *subsnes)
278 {
279   SNES_Composite     *jac;
280   SNES_CompositeLink next;
281   PetscInt         i;
282 
283   PetscFunctionBegin;
284   jac  = (SNES_Composite*)snes->data;
285   next = jac->head;
286   for (i=0; i<n; i++) {
287     if (!next->next) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_INCOMP,"Not enough SNESes in composite preconditioner");
288     next = next->next;
289   }
290   *subsnes = next->snes;
291   PetscFunctionReturn(0);
292 }
293 
294 /* -------------------------------------------------------------------------------- */
295 #undef __FUNCT__
296 #define __FUNCT__ "SNESCompositeSetType"
297 /*@
298    SNESCompositeSetType - Sets the type of composite preconditioner.
299 
300    Logically Collective on SNES
301 
302    Input Parameter:
303 +  snes - the preconditioner context
304 -  type - SNES_COMPOSITE_ADDITIVE (default), SNES_COMPOSITE_MULTIPLICATIVE
305 
306    Options Database Key:
307 .  -snes_composite_type <type: one of multiplicative, additive, special> - Sets composite preconditioner type
308 
309    Level: Developer
310 
311 .keywords: SNES, set, type, composite preconditioner, additive, multiplicative
312 @*/
313 PetscErrorCode  SNESCompositeSetType(SNES snes,SNESCompositeType type)
314 {
315   PetscErrorCode ierr;
316 
317   PetscFunctionBegin;
318   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
319   PetscValidLogicalCollectiveEnum(snes,type,2);
320   ierr = PetscTryMethod(snes,"SNESCompositeSetType_C",(SNES,SNESCompositeType),(snes,type));CHKERRQ(ierr);
321   PetscFunctionReturn(0);
322 }
323 
324 #undef __FUNCT__
325 #define __FUNCT__ "SNESCompositeAddSNES"
326 /*@C
327    SNESCompositeAddSNES - Adds another SNES to the composite SNES.
328 
329    Collective on SNES
330 
331    Input Parameters:
332 +  snes - the preconditioner context
333 -  type - the type of the new preconditioner
334 
335    Level: Developer
336 
337 .keywords: SNES, composite preconditioner, add
338 @*/
339 PetscErrorCode  SNESCompositeAddSNES(SNES snes,SNESType type)
340 {
341   PetscErrorCode ierr;
342 
343   PetscFunctionBegin;
344   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
345   ierr = PetscTryMethod(snes,"SNESCompositeAddSNES_C",(SNES,SNESType),(snes,type));CHKERRQ(ierr);
346   PetscFunctionReturn(0);
347 }
348 
349 #undef __FUNCT__
350 #define __FUNCT__ "SNESCompositeGetSNES"
351 /*@
352    SNESCompositeGetSNES - Gets one of the SNES objects in the composite SNES.
353 
354    Not Collective
355 
356    Input Parameter:
357 +  snes - the preconditioner context
358 -  n - the number of the snes requested
359 
360    Output Parameters:
361 .  subsnes - the SNES requested
362 
363    Level: Developer
364 
365 .keywords: SNES, get, composite preconditioner, sub preconditioner
366 
367 .seealso: SNESCompositeAddSNES()
368 @*/
369 PetscErrorCode  SNESCompositeGetSNES(SNES snes,PetscInt n,SNES *subsnes)
370 {
371   PetscErrorCode ierr;
372 
373   PetscFunctionBegin;
374   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
375   PetscValidPointer(subsnes,3);
376   ierr = PetscUseMethod(snes,"SNESCompositeGetSNES_C",(SNES,PetscInt,SNES*),(snes,n,subsnes));CHKERRQ(ierr);
377   PetscFunctionReturn(0);
378 }
379 
380 #undef __FUNCT__
381 #define __FUNCT__ "SNESSolve_Composite"
382 PetscErrorCode SNESSolve_Composite(SNES snes)
383 {
384   Vec            F;
385   Vec            X;
386   Vec            B;
387   PetscInt       i;
388   PetscReal      fnorm = 0.0;
389   PetscErrorCode ierr;
390   SNESNormType   normtype;
391   SNES_Composite *comp = (SNES_Composite*)snes->data;
392 
393   PetscFunctionBegin;
394   X = snes->vec_sol;
395   F = snes->vec_func;
396   B = snes->vec_rhs;
397 
398   ierr         = PetscObjectAMSTakeAccess((PetscObject)snes);CHKERRQ(ierr);
399   snes->iter   = 0;
400   snes->norm   = 0.;
401   ierr         = PetscObjectAMSGrantAccess((PetscObject)snes);CHKERRQ(ierr);
402   snes->reason = SNES_CONVERGED_ITERATING;
403   ierr         = SNESGetNormType(snes, &normtype);CHKERRQ(ierr);
404   if (normtype == SNES_NORM_FUNCTION || normtype == SNES_NORM_INITIAL_ONLY || normtype == SNES_NORM_INITIAL_FINAL_ONLY) {
405     if (!snes->vec_func_init_set) {
406       ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
407       if (snes->domainerror) {
408         snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
409         PetscFunctionReturn(0);
410       }
411     } else snes->vec_func_init_set = PETSC_FALSE;
412 
413     /* convergence test */
414     if (!snes->norm_init_set) {
415       ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F||  */
416       if (PetscIsInfOrNanReal(fnorm)) {
417         snes->reason = SNES_DIVERGED_FNORM_NAN;
418         PetscFunctionReturn(0);
419       }
420     } else {
421       fnorm               = snes->norm_init;
422       snes->norm_init_set = PETSC_FALSE;
423     }
424     ierr       = PetscObjectAMSTakeAccess((PetscObject)snes);CHKERRQ(ierr);
425     snes->iter = 0;
426     snes->norm = fnorm;
427     ierr       = PetscObjectAMSGrantAccess((PetscObject)snes);CHKERRQ(ierr);
428     ierr       = SNESLogConvergenceHistory(snes,snes->norm,0);CHKERRQ(ierr);
429     ierr       = SNESMonitor(snes,0,snes->norm);CHKERRQ(ierr);
430 
431     /* set parameter for default relative tolerance convergence test */
432     snes->ttol = fnorm*snes->rtol;
433 
434     /* test convergence */
435     ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
436     if (snes->reason) PetscFunctionReturn(0);
437   } else {
438     ierr = PetscObjectAMSGrantAccess((PetscObject)snes);CHKERRQ(ierr);
439     ierr = SNESLogConvergenceHistory(snes,snes->norm,0);CHKERRQ(ierr);
440     ierr = SNESMonitor(snes,0,snes->norm);CHKERRQ(ierr);
441   }
442 
443   /* Call general purpose update function */
444   if (snes->ops->update) {
445     ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
446   }
447 
448   for (i = 0; i < snes->max_its; i++) {
449     if (comp->type == SNES_COMPOSITE_ADDITIVE) {
450       ierr = SNESCompositeApply_Additive(snes,X,B,F,&fnorm);CHKERRQ(ierr);
451     } else if (comp->type == SNES_COMPOSITE_MULTIPLICATIVE) {
452       ierr = SNESCompositeApply_Multiplicative(snes,X,B,F,&fnorm);CHKERRQ(ierr);
453     } else {
454     }
455     if ((i == snes->max_its - 1) && (normtype == SNES_NORM_INITIAL_FINAL_ONLY || normtype == SNES_NORM_FINAL_ONLY)) {
456       ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
457       if (snes->domainerror) {
458         snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
459         break;
460       }
461       ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F||  */
462       if (PetscIsInfOrNanReal(fnorm)) {
463         snes->reason = SNES_DIVERGED_FNORM_NAN;
464         break;
465       }
466     }
467     /* Monitor convergence */
468     ierr       = PetscObjectAMSTakeAccess((PetscObject)snes);CHKERRQ(ierr);
469     snes->iter = i+1;
470     snes->norm = fnorm;
471     ierr       = PetscObjectAMSGrantAccess((PetscObject)snes);CHKERRQ(ierr);
472     ierr       = SNESLogConvergenceHistory(snes,snes->norm,0);CHKERRQ(ierr);
473     ierr       = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr);
474     /* Test for convergence */
475     if (normtype == SNES_NORM_FUNCTION) {ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);}
476     if (snes->reason) break;
477     /* Call general purpose update function */
478     if (snes->ops->update) {ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);}
479   }
480   if (normtype == SNES_NORM_FUNCTION) {
481     if (i == snes->max_its) {
482       ierr = PetscInfo1(snes,"Maximum number of iterations has been reached: %D\n",snes->max_its);CHKERRQ(ierr);
483       if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
484     }
485   } else if (!snes->reason) snes->reason = SNES_CONVERGED_ITS;
486   PetscFunctionReturn(0);
487 }
488 
489 /* -------------------------------------------------------------------------------------------*/
490 
491 /*MC
492      SNESCOMPOSITE - Build a preconditioner by composing together several nonlinear solvers
493 
494    Options Database Keys:
495 +  -snes_composite_type <type: one of multiplicative, additive, symmetric_multiplicative, special> - Sets composite preconditioner type
496 -  -snes_composite_sneses - <snes0,snes1,...> list of SNESes to compose
497 
498    Level: intermediate
499 
500    Concepts: composing solvers
501 
502 .seealso:  SNESCreate(), SNESSetType(), SNESType (for list of available types), SNES,
503            SNESSHELL, SNESCompositeSetType(), SNESCompositeSpecialSetAlpha(), SNESCompositeAddSNES(),
504            SNESCompositeGetSNES()
505 
506 M*/
507 
508 #undef __FUNCT__
509 #define __FUNCT__ "SNESCreate_Composite"
510 PETSC_EXTERN PetscErrorCode SNESCreate_Composite(SNES snes)
511 {
512   PetscErrorCode ierr;
513   SNES_Composite   *jac;
514 
515   PetscFunctionBegin;
516   ierr = PetscNewLog(snes,SNES_Composite,&jac);CHKERRQ(ierr);
517 
518   snes->ops->solve           = SNESSolve_Composite;
519   snes->ops->setup           = SNESSetUp_Composite;
520   snes->ops->reset           = SNESReset_Composite;
521   snes->ops->destroy         = SNESDestroy_Composite;
522   snes->ops->setfromoptions  = SNESSetFromOptions_Composite;
523   snes->ops->view            = SNESView_Composite;
524 
525   snes->data = (void*)jac;
526   jac->type  = SNES_COMPOSITE_ADDITIVE;
527   jac->head  = 0;
528 
529   ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESCompositeSetType_C","SNESCompositeSetType_Composite",SNESCompositeSetType_Composite);CHKERRQ(ierr);
530   ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESCompositeAddSNES_C","SNESCompositeAddSNES_Composite",SNESCompositeAddSNES_Composite);CHKERRQ(ierr);
531   ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESCompositeGetSNES_C","SNESCompositeGetSNES_Composite",SNESCompositeGetSNES_Composite);CHKERRQ(ierr);
532   PetscFunctionReturn(0);
533 }
534 
535