xref: /petsc/src/snes/impls/composite/snescomposite.c (revision 9c2f34738574ff127efffa739c2c7ed1973d0335)
1eed5f15bSPeter Brune 
2eed5f15bSPeter Brune /*
3eed5f15bSPeter Brune       Defines a SNES that can consist of a collection of SNESes
4eed5f15bSPeter Brune */
5eed5f15bSPeter Brune #include <petsc-private/snesimpl.h> /*I "petscsnes.h" I*/
690a8ba9bSPeter Brune #include <petscblaslapack.h>
7eed5f15bSPeter Brune 
890a8ba9bSPeter Brune const char *const        SNESCompositeTypes[]   = {"ADDITIVE","MULTIPLICATIVE","ADDITIVEOPTIMAL","SNESCompositeType","SNES_COMPOSITE",0};
9eed5f15bSPeter Brune 
10eed5f15bSPeter Brune typedef struct _SNES_CompositeLink *SNES_CompositeLink;
11eed5f15bSPeter Brune struct _SNES_CompositeLink {
12eed5f15bSPeter Brune   SNES               snes;
138f626970SPeter Brune   PetscReal          dmp;
1490a8ba9bSPeter Brune   Vec                X;
15eed5f15bSPeter Brune   SNES_CompositeLink next;
16eed5f15bSPeter Brune   SNES_CompositeLink previous;
17eed5f15bSPeter Brune };
18eed5f15bSPeter Brune 
19eed5f15bSPeter Brune typedef struct {
20eed5f15bSPeter Brune   SNES_CompositeLink head;
2190a8ba9bSPeter Brune   PetscInt           nsnes;
22eed5f15bSPeter Brune   SNESCompositeType  type;
23eed5f15bSPeter Brune   Vec                Xorig;
2490a8ba9bSPeter Brune 
2590a8ba9bSPeter Brune   /* context for ADDITIVEOPTIMAL */
2690a8ba9bSPeter Brune   Vec                *Xes,*Fes;      /* solution and residual vectors for the subsolvers */
27*9c2f3473SPeter Brune   PetscReal          *fnorms;        /* norms of the solutions */
2890a8ba9bSPeter Brune   PetscScalar        *h;             /* the matrix formed as q_ij = (rdot_i, rdot_j) */
29*9c2f3473SPeter Brune   PetscScalar        *g;             /* the dotproducts of the previous function with the candidate functions */
30*9c2f3473SPeter Brune   PetscBLASInt       n;              /* matrix dimension -- nsnes */
3190a8ba9bSPeter Brune   PetscBLASInt       nrhs;           /* the number of right hand sides */
3290a8ba9bSPeter Brune   PetscBLASInt       lda;            /* the padded matrix dimension */
3390a8ba9bSPeter Brune   PetscBLASInt       ldb;            /* the padded vector dimension */
3490a8ba9bSPeter Brune   PetscReal          *s;             /* the singular values */
35*9c2f3473SPeter Brune   PetscScalar        *beta;          /* the RHS and combination */
3690a8ba9bSPeter Brune   PetscReal          rcond;          /* the exit condition */
3790a8ba9bSPeter Brune   PetscBLASInt       rank;           /* the effective rank */
3890a8ba9bSPeter Brune   PetscScalar        *work;          /* the work vector */
3990a8ba9bSPeter Brune   PetscReal          *rwork;         /* the real work vector used for complex */
4090a8ba9bSPeter Brune   PetscBLASInt       lwork;          /* the size of the work vector */
4190a8ba9bSPeter Brune   PetscBLASInt       info;           /* the output condition */
4290a8ba9bSPeter Brune 
43eed5f15bSPeter Brune } SNES_Composite;
44eed5f15bSPeter Brune 
45eed5f15bSPeter Brune #undef __FUNCT__
46eed5f15bSPeter Brune #define __FUNCT__ "SNESCompositeApply_Multiplicative"
47eed5f15bSPeter Brune static PetscErrorCode SNESCompositeApply_Multiplicative(SNES snes,Vec X,Vec B,Vec F,PetscReal *fnorm)
48eed5f15bSPeter Brune {
49eed5f15bSPeter Brune   PetscErrorCode     ierr;
50eed5f15bSPeter Brune   SNES_Composite     *jac = (SNES_Composite*)snes->data;
51eed5f15bSPeter Brune   SNES_CompositeLink next = jac->head;
52eed5f15bSPeter Brune   Vec                FSub;
53eed5f15bSPeter Brune   PetscReal          fsubnorm;
54eed5f15bSPeter Brune 
55eed5f15bSPeter Brune   PetscFunctionBegin;
56eed5f15bSPeter Brune   if (!next) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"No composite SNESes supplied via SNESCompositeAddSNES() or -snes_composite_sneses");
57eed5f15bSPeter Brune   if (snes->normtype == SNES_NORM_FUNCTION) {
58eed5f15bSPeter Brune     ierr = SNESSetInitialFunction(next->snes,F);CHKERRQ(ierr);
59eed5f15bSPeter Brune     if (fnorm) {ierr = SNESSetInitialFunctionNorm(next->snes,*fnorm);CHKERRQ(ierr);}
60eed5f15bSPeter Brune   }
61eed5f15bSPeter Brune   ierr = SNESSolve(next->snes,B,X);CHKERRQ(ierr);
62eed5f15bSPeter Brune 
63eed5f15bSPeter Brune   while (next->next) {
64eed5f15bSPeter Brune     /* only copy the function over in the case where the functions correspond */
65eed5f15bSPeter Brune     if (next->snes->pcside == PC_RIGHT && next->snes->normtype != SNES_NORM_NONE) {
66eed5f15bSPeter Brune       ierr = SNESGetFunction(next->snes,&FSub,NULL,NULL);CHKERRQ(ierr);
67eed5f15bSPeter Brune       ierr = SNESGetFunctionNorm(next->snes,&fsubnorm);CHKERRQ(ierr);
68eed5f15bSPeter Brune       next = next->next;
69eed5f15bSPeter Brune       ierr = SNESSetInitialFunction(next->snes,FSub);CHKERRQ(ierr);
70eed5f15bSPeter Brune       ierr = SNESSetInitialFunctionNorm(next->snes,fsubnorm);CHKERRQ(ierr);
71eed5f15bSPeter Brune     } else {
72eed5f15bSPeter Brune       next = next->next;
73eed5f15bSPeter Brune     }
74eed5f15bSPeter Brune     ierr = SNESSolve(next->snes,B,X);CHKERRQ(ierr);
75eed5f15bSPeter Brune   }
76eed5f15bSPeter Brune   if (next->snes->pcside == PC_RIGHT) {
77eed5f15bSPeter Brune     ierr = SNESGetFunction(next->snes,&FSub,NULL,NULL);CHKERRQ(ierr);
78eed5f15bSPeter Brune     ierr = VecCopy(FSub,F);CHKERRQ(ierr);
79eed5f15bSPeter Brune     if (fnorm) {ierr = SNESGetFunctionNorm(next->snes,fnorm);CHKERRQ(ierr);}
80eed5f15bSPeter Brune   } else if (snes->normtype == SNES_NORM_FUNCTION) {
81eed5f15bSPeter Brune     SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
82eed5f15bSPeter Brune     if (fnorm) {ierr = VecNorm(F,NORM_2,fnorm);CHKERRQ(ierr);}
83eed5f15bSPeter Brune   }
84eed5f15bSPeter Brune   PetscFunctionReturn(0);
85eed5f15bSPeter Brune }
86eed5f15bSPeter Brune 
87eed5f15bSPeter Brune #undef __FUNCT__
88eed5f15bSPeter Brune #define __FUNCT__ "SNESCompositeApply_Additive"
89eed5f15bSPeter Brune static PetscErrorCode SNESCompositeApply_Additive(SNES snes,Vec X,Vec B,Vec F,PetscReal *fnorm)
90eed5f15bSPeter Brune {
91eed5f15bSPeter Brune   PetscErrorCode     ierr;
92eed5f15bSPeter Brune   SNES_Composite     *jac = (SNES_Composite*)snes->data;
93eed5f15bSPeter Brune   SNES_CompositeLink next = jac->head;
94eed5f15bSPeter Brune   Vec                Y,Xorig;
95eed5f15bSPeter Brune 
96eed5f15bSPeter Brune   PetscFunctionBegin;
97eed5f15bSPeter Brune   Y = snes->vec_sol_update;
98eed5f15bSPeter Brune   if (!jac->Xorig) {ierr = VecDuplicate(X,&jac->Xorig);CHKERRQ(ierr);}
99eed5f15bSPeter Brune   Xorig = jac->Xorig;
100eed5f15bSPeter Brune   ierr = VecCopy(X,Xorig);
101eed5f15bSPeter Brune   if (!next) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"No composite SNESes supplied via SNESCompositeAddSNES() or -snes_composite_sneses");
102eed5f15bSPeter Brune   if (snes->normtype == SNES_NORM_FUNCTION) {
103eed5f15bSPeter Brune     ierr = SNESSetInitialFunction(next->snes,F);CHKERRQ(ierr);
104eed5f15bSPeter Brune     if (fnorm) {ierr = SNESSetInitialFunctionNorm(next->snes,*fnorm);CHKERRQ(ierr);}
105eed5f15bSPeter Brune     while (next->next) {
106eed5f15bSPeter Brune       next = next->next;
107eed5f15bSPeter Brune       ierr = SNESSetInitialFunction(next->snes,F);CHKERRQ(ierr);
108eed5f15bSPeter Brune       if (fnorm) {ierr = SNESSetInitialFunctionNorm(next->snes,*fnorm);CHKERRQ(ierr);}
109eed5f15bSPeter Brune     }
110eed5f15bSPeter Brune   }
111eed5f15bSPeter Brune   next = jac->head;
1128f626970SPeter Brune   ierr = VecCopy(Xorig,Y);CHKERRQ(ierr);
113eed5f15bSPeter Brune   ierr = SNESSolve(next->snes,B,Y);CHKERRQ(ierr);
114eed5f15bSPeter Brune   ierr = VecAXPY(Y,-1.0,Xorig);CHKERRQ(ierr);
1158f626970SPeter Brune   ierr = VecAXPY(X,next->dmp,Y);CHKERRQ(ierr);
1168f626970SPeter Brune   while (next->next) {
1178f626970SPeter Brune     next = next->next;
1188f626970SPeter Brune     ierr = VecCopy(Xorig,Y);CHKERRQ(ierr);
1198f626970SPeter Brune     ierr = SNESSolve(next->snes,B,Y);CHKERRQ(ierr);
1208f626970SPeter Brune     ierr = VecAXPY(Y,-1.0,Xorig);CHKERRQ(ierr);
1218f626970SPeter Brune     ierr = VecAXPY(X,next->dmp,Y);CHKERRQ(ierr);
122eed5f15bSPeter Brune   }
123eed5f15bSPeter Brune   if (snes->normtype == SNES_NORM_FUNCTION) {
124eed5f15bSPeter Brune     ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
125eed5f15bSPeter Brune     if (fnorm) {ierr = VecNorm(F,NORM_2,fnorm);CHKERRQ(ierr);}
126eed5f15bSPeter Brune   }
127eed5f15bSPeter Brune   PetscFunctionReturn(0);
128eed5f15bSPeter Brune }
12990a8ba9bSPeter Brune 
13090a8ba9bSPeter Brune #undef __FUNCT__
13190a8ba9bSPeter Brune #define __FUNCT__ "SNESCompositeApply_AdditiveOptimal"
13290a8ba9bSPeter Brune /* approximately solve the overdetermined system:
13390a8ba9bSPeter Brune 
13490a8ba9bSPeter Brune  2*F(x_i)\cdot F(\x_j)\alpha_i = 0
13590a8ba9bSPeter Brune  \alpha_i                      = 1
13690a8ba9bSPeter Brune 
13790a8ba9bSPeter Brune  Which minimizes the L2 norm of the linearization of:
13890a8ba9bSPeter Brune  ||F(\sum_i \alpha_i*x_i)||^2
13990a8ba9bSPeter Brune 
14090a8ba9bSPeter Brune  With the constraint that \sum_i\alpha_i = 1
14190a8ba9bSPeter Brune  Where x_i is the solution from the ith subsolver.
14290a8ba9bSPeter Brune  */
14390a8ba9bSPeter Brune static PetscErrorCode SNESCompositeApply_AdditiveOptimal(SNES snes,Vec X,Vec B,Vec F,PetscReal *fnorm)
14490a8ba9bSPeter Brune {
14590a8ba9bSPeter Brune   PetscErrorCode     ierr;
14690a8ba9bSPeter Brune   SNES_Composite     *jac = (SNES_Composite*)snes->data;
14790a8ba9bSPeter Brune   SNES_CompositeLink next = jac->head;
14890a8ba9bSPeter Brune   Vec                *Xes = jac->Xes,*Fes = jac->Fes;
14990a8ba9bSPeter Brune   PetscInt           i,j;
150*9c2f3473SPeter Brune   PetscScalar        tot,ftf;
15190a8ba9bSPeter Brune 
15290a8ba9bSPeter Brune   PetscFunctionBegin;
15390a8ba9bSPeter Brune   if (!next) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"No composite SNESes supplied via SNESCompositeAddSNES() or -snes_composite_sneses");
15490a8ba9bSPeter Brune   next = jac->head;
15590a8ba9bSPeter Brune   i = 0;
15690a8ba9bSPeter Brune   ierr = VecCopy(X,Xes[i]);CHKERRQ(ierr);
15790a8ba9bSPeter Brune   ierr = SNESSolve(next->snes,B,Xes[i]);CHKERRQ(ierr);
15890a8ba9bSPeter Brune   while (next->next) {
15990a8ba9bSPeter Brune     i++;
16090a8ba9bSPeter Brune     next = next->next;
16190a8ba9bSPeter Brune     ierr = VecCopy(X,Xes[i]);CHKERRQ(ierr);
16290a8ba9bSPeter Brune     ierr = SNESSolve(next->snes,B,Xes[i]);CHKERRQ(ierr);
16390a8ba9bSPeter Brune   }
16490a8ba9bSPeter Brune 
16590a8ba9bSPeter Brune   /* all the solutions are collected; combine optimally */
16690a8ba9bSPeter Brune   for (i=0;i<jac->n;i++) {
16790a8ba9bSPeter Brune     for (j=0;j<i+1;j++) {
168*9c2f3473SPeter Brune       ierr = VecDotBegin(Fes[i],Fes[j],&jac->h[i + j*jac->n]);CHKERRQ(ierr);
16990a8ba9bSPeter Brune     }
170*9c2f3473SPeter Brune     ierr = VecDotBegin(Fes[i],F,&jac->g[i]);CHKERRQ(ierr);
17190a8ba9bSPeter Brune   }
17290a8ba9bSPeter Brune   for (i=0;i<jac->n;i++) {
17390a8ba9bSPeter Brune     for (j=0;j<i+1;j++) {
174*9c2f3473SPeter Brune       ierr = VecDotEnd(Fes[i],Fes[j],&jac->h[i + j*jac->n]);CHKERRQ(ierr);
175*9c2f3473SPeter Brune       if (i == j) jac->fnorms[i] = PetscSqrtReal(PetscRealPart(jac->h[i + j*jac->n]));
17690a8ba9bSPeter Brune     }
177*9c2f3473SPeter Brune     ierr = VecDotEnd(Fes[i],F,&jac->g[i]);CHKERRQ(ierr);
17890a8ba9bSPeter Brune   }
17990a8ba9bSPeter Brune 
180*9c2f3473SPeter Brune   ftf = (*fnorm)*(*fnorm);
18190a8ba9bSPeter Brune 
18290a8ba9bSPeter Brune   for (i=0; i<jac->n; i++) {
18390a8ba9bSPeter Brune     for (j=i+1;j<jac->n;j++) {
184*9c2f3473SPeter Brune       jac->h[i + j*jac->n] = jac->h[j + i*jac->n];
18590a8ba9bSPeter Brune     }
18690a8ba9bSPeter Brune   }
18790a8ba9bSPeter Brune 
18890a8ba9bSPeter Brune   for (i=0; i<jac->n; i++) {
189*9c2f3473SPeter Brune     for (j=0; j<jac->n; j++) {
190*9c2f3473SPeter Brune       jac->h[i + j*jac->n] = jac->h[i + j*jac->n] - jac->g[j] - jac->g[i] + ftf;
19190a8ba9bSPeter Brune     }
192*9c2f3473SPeter Brune     jac->beta[i] = ftf - jac->g[i];
193*9c2f3473SPeter Brune   }
19490a8ba9bSPeter Brune 
19590a8ba9bSPeter Brune #if defined(PETSC_MISSING_LAPACK_GELSS)
19690a8ba9bSPeter Brune   SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_SUP,"SNESCOMPOSITE with ADDITIVEOPTIMAL requires the LAPACK GELSS routine.");
19790a8ba9bSPeter Brune #else
19890a8ba9bSPeter Brune   jac->info  = 0;
19990a8ba9bSPeter Brune   jac->rcond = -1.;
20090a8ba9bSPeter Brune   ierr          = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
20190a8ba9bSPeter Brune #if defined(PETSC_USE_COMPLEX)
202*9c2f3473SPeter Brune   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));
20390a8ba9bSPeter Brune #else
204*9c2f3473SPeter Brune   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));
20590a8ba9bSPeter Brune #endif
20690a8ba9bSPeter Brune   ierr = PetscFPTrapPop();CHKERRQ(ierr);
20790a8ba9bSPeter Brune   if (jac->info < 0) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_LIB,"Bad argument to GELSS");
20890a8ba9bSPeter Brune   if (jac->info > 0) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_LIB,"SVD failed to converge");
20990a8ba9bSPeter Brune #endif
210*9c2f3473SPeter Brune   tot = 0.;
21190a8ba9bSPeter Brune   for (i=0; i<jac->n; i++) {
21290a8ba9bSPeter Brune     if (PetscIsInfOrNanScalar(jac->beta[i])) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_LIB,"SVD generated inconsistent output");
213*9c2f3473SPeter Brune     ierr = PetscInfo2(snes,"%d: %f\n",i,jac->beta[i]);CHKERRQ(ierr);
214*9c2f3473SPeter Brune     tot += jac->beta[i];
21590a8ba9bSPeter Brune   }
216*9c2f3473SPeter Brune   ierr = VecScale(X,(1. - tot));CHKERRQ(ierr);
21790a8ba9bSPeter Brune   ierr = VecMAXPY(X,jac->n,jac->beta,Xes);CHKERRQ(ierr);
21890a8ba9bSPeter Brune   ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
219*9c2f3473SPeter Brune   ierr = VecNorm(F,NORM_2,fnorm);CHKERRQ(ierr);
22090a8ba9bSPeter Brune 
221*9c2f3473SPeter Brune   /* take the minimum-normed candidate if it beats the combination */
222*9c2f3473SPeter Brune   for (i=0; i<jac->n; i++) {
223*9c2f3473SPeter Brune     if (jac->fnorms[i] < *fnorm) {
224*9c2f3473SPeter Brune       ierr = VecCopy(Xes[i],X);CHKERRQ(ierr);
225*9c2f3473SPeter Brune       ierr = VecCopy(Fes[i],F);CHKERRQ(ierr);
226*9c2f3473SPeter Brune       *fnorm = jac->fnorms[i];
227*9c2f3473SPeter Brune     }
228*9c2f3473SPeter Brune   }
22990a8ba9bSPeter Brune   PetscFunctionReturn(0);
23090a8ba9bSPeter Brune }
23190a8ba9bSPeter Brune 
232eed5f15bSPeter Brune #undef __FUNCT__
233eed5f15bSPeter Brune #define __FUNCT__ "SNESSetUp_Composite"
234eed5f15bSPeter Brune static PetscErrorCode SNESSetUp_Composite(SNES snes)
235eed5f15bSPeter Brune {
236eed5f15bSPeter Brune   PetscErrorCode     ierr;
237eed5f15bSPeter Brune   SNES_Composite     *jac = (SNES_Composite*)snes->data;
238eed5f15bSPeter Brune   SNES_CompositeLink next = jac->head;
23990a8ba9bSPeter Brune   PetscInt           n=0,i;
24090a8ba9bSPeter Brune   Vec                F;
241eed5f15bSPeter Brune 
242eed5f15bSPeter Brune   PetscFunctionBegin;
243eed5f15bSPeter Brune   while (next) {
24490a8ba9bSPeter Brune     n++;
245eed5f15bSPeter Brune     ierr = SNESSetFromOptions(next->snes);CHKERRQ(ierr);
246eed5f15bSPeter Brune     next = next->next;
247eed5f15bSPeter Brune   }
24890a8ba9bSPeter Brune   jac->nsnes = n;
24990a8ba9bSPeter Brune   ierr = SNESGetFunction(snes,&F,NULL,NULL);CHKERRQ(ierr);
25090a8ba9bSPeter Brune   if (jac->type == SNES_COMPOSITE_ADDITIVEOPTIMAL) {
25190a8ba9bSPeter Brune     ierr = VecDuplicateVecs(F,jac->nsnes,&jac->Xes);CHKERRQ(ierr);
25290a8ba9bSPeter Brune     ierr = PetscMalloc(sizeof(Vec)*n,&jac->Fes);CHKERRQ(ierr);
253*9c2f3473SPeter Brune     ierr = PetscMalloc(sizeof(PetscReal)*n,&jac->fnorms);CHKERRQ(ierr);
25490a8ba9bSPeter Brune     next = jac->head;
25590a8ba9bSPeter Brune     i = 0;
25690a8ba9bSPeter Brune     while (next) {
25790a8ba9bSPeter Brune       ierr = SNESGetFunction(next->snes,&F,NULL,NULL);CHKERRQ(ierr);
25890a8ba9bSPeter Brune       jac->Fes[i] = F;
25990a8ba9bSPeter Brune       ierr = PetscObjectReference((PetscObject)F);CHKERRQ(ierr);
26090a8ba9bSPeter Brune       next = next->next;
26190a8ba9bSPeter Brune       i++;
26290a8ba9bSPeter Brune     }
26390a8ba9bSPeter Brune     /* allocate the subspace direct solve area */
26490a8ba9bSPeter Brune     jac->nrhs  = 1;
265*9c2f3473SPeter Brune     jac->lda   = jac->nsnes;
26690a8ba9bSPeter Brune     jac->ldb   = jac->nsnes;
26790a8ba9bSPeter Brune     jac->n     = jac->nsnes;
26890a8ba9bSPeter Brune 
269*9c2f3473SPeter Brune     ierr = PetscMalloc(jac->n*jac->n*sizeof(PetscScalar),&jac->h);CHKERRQ(ierr);
270*9c2f3473SPeter Brune     ierr = PetscMalloc(jac->n*sizeof(PetscScalar),&jac->beta);CHKERRQ(ierr);
271*9c2f3473SPeter Brune     ierr = PetscMalloc(jac->n*sizeof(PetscScalar),&jac->s);CHKERRQ(ierr);
272*9c2f3473SPeter Brune     ierr = PetscMalloc(jac->n*sizeof(PetscScalar),&jac->g);CHKERRQ(ierr);
273*9c2f3473SPeter Brune     jac->lwork = 12*jac->n;
27490a8ba9bSPeter Brune #if PETSC_USE_COMPLEX
27590a8ba9bSPeter Brune     ierr = PetscMalloc(sizeof(PetscReal)*jac->lwork,&jac->rwork);CHKERRQ(ierr);
27690a8ba9bSPeter Brune #endif
27790a8ba9bSPeter Brune     ierr = PetscMalloc(sizeof(PetscScalar)*jac->lwork,&jac->work);CHKERRQ(ierr);
27890a8ba9bSPeter Brune   }
27990a8ba9bSPeter Brune 
280eed5f15bSPeter Brune   PetscFunctionReturn(0);
281eed5f15bSPeter Brune }
282eed5f15bSPeter Brune 
283eed5f15bSPeter Brune #undef __FUNCT__
284eed5f15bSPeter Brune #define __FUNCT__ "SNESReset_Composite"
285eed5f15bSPeter Brune static PetscErrorCode SNESReset_Composite(SNES snes)
286eed5f15bSPeter Brune {
287eed5f15bSPeter Brune   SNES_Composite     *jac = (SNES_Composite*)snes->data;
288eed5f15bSPeter Brune   PetscErrorCode   ierr;
289eed5f15bSPeter Brune   SNES_CompositeLink next = jac->head;
290eed5f15bSPeter Brune 
291eed5f15bSPeter Brune   PetscFunctionBegin;
292eed5f15bSPeter Brune   while (next) {
293eed5f15bSPeter Brune     ierr = SNESReset(next->snes);CHKERRQ(ierr);
294eed5f15bSPeter Brune     next = next->next;
295eed5f15bSPeter Brune   }
296eed5f15bSPeter Brune   ierr = VecDestroy(&jac->Xorig);CHKERRQ(ierr);
29790a8ba9bSPeter Brune   if (jac->Xes) {ierr = VecDestroyVecs(jac->nsnes,&jac->Xes);CHKERRQ(ierr);}
29890a8ba9bSPeter Brune   if (jac->Fes) {ierr = VecDestroyVecs(jac->nsnes,&jac->Fes);CHKERRQ(ierr);}
299*9c2f3473SPeter Brune   ierr = PetscFree(jac->fnorms);CHKERRQ(ierr);
30090a8ba9bSPeter Brune   ierr = PetscFree(jac->h);CHKERRQ(ierr);
30190a8ba9bSPeter Brune   ierr = PetscFree(jac->s);CHKERRQ(ierr);
302*9c2f3473SPeter Brune   ierr = PetscFree(jac->g);CHKERRQ(ierr);
30390a8ba9bSPeter Brune   ierr = PetscFree(jac->beta);CHKERRQ(ierr);
30490a8ba9bSPeter Brune   ierr = PetscFree(jac->work);CHKERRQ(ierr);
30590a8ba9bSPeter Brune   ierr = PetscFree(jac->rwork);CHKERRQ(ierr);
30690a8ba9bSPeter Brune 
307eed5f15bSPeter Brune   PetscFunctionReturn(0);
308eed5f15bSPeter Brune }
309eed5f15bSPeter Brune 
310eed5f15bSPeter Brune #undef __FUNCT__
311eed5f15bSPeter Brune #define __FUNCT__ "SNESDestroy_Composite"
312eed5f15bSPeter Brune static PetscErrorCode SNESDestroy_Composite(SNES snes)
313eed5f15bSPeter Brune {
314eed5f15bSPeter Brune   SNES_Composite     *jac = (SNES_Composite*)snes->data;
315eed5f15bSPeter Brune   PetscErrorCode   ierr;
316eed5f15bSPeter Brune   SNES_CompositeLink next = jac->head,next_tmp;
317eed5f15bSPeter Brune 
318eed5f15bSPeter Brune   PetscFunctionBegin;
319eed5f15bSPeter Brune   ierr = SNESReset_Composite(snes);CHKERRQ(ierr);
320eed5f15bSPeter Brune   while (next) {
321eed5f15bSPeter Brune     ierr     = SNESDestroy(&next->snes);CHKERRQ(ierr);
322eed5f15bSPeter Brune     next_tmp = next;
323eed5f15bSPeter Brune     next     = next->next;
324eed5f15bSPeter Brune     ierr     = PetscFree(next_tmp);CHKERRQ(ierr);
325eed5f15bSPeter Brune   }
326eed5f15bSPeter Brune   ierr = PetscFree(snes->data);CHKERRQ(ierr);
327eed5f15bSPeter Brune   PetscFunctionReturn(0);
328eed5f15bSPeter Brune }
329eed5f15bSPeter Brune 
330eed5f15bSPeter Brune #undef __FUNCT__
331eed5f15bSPeter Brune #define __FUNCT__ "SNESSetFromOptions_Composite"
332eed5f15bSPeter Brune static PetscErrorCode SNESSetFromOptions_Composite(SNES snes)
333eed5f15bSPeter Brune {
334eed5f15bSPeter Brune   SNES_Composite     *jac = (SNES_Composite*)snes->data;
335eed5f15bSPeter Brune   PetscErrorCode     ierr;
336eed5f15bSPeter Brune   PetscInt           nmax = 8,i;
337eed5f15bSPeter Brune   SNES_CompositeLink next;
338eed5f15bSPeter Brune   char               *sneses[8];
3398f626970SPeter Brune   PetscReal          dmps[8];
340eed5f15bSPeter Brune   PetscBool          flg;
341eed5f15bSPeter Brune 
342eed5f15bSPeter Brune   PetscFunctionBegin;
343eed5f15bSPeter Brune   ierr = PetscOptionsHead("Composite preconditioner options");CHKERRQ(ierr);
344eed5f15bSPeter Brune   ierr = PetscOptionsEnum("-snes_composite_type","Type of composition","SNESCompositeSetType",SNESCompositeTypes,(PetscEnum)jac->type,(PetscEnum*)&jac->type,&flg);CHKERRQ(ierr);
345eed5f15bSPeter Brune   if (flg) {
346eed5f15bSPeter Brune     ierr = SNESCompositeSetType(snes,jac->type);CHKERRQ(ierr);
347eed5f15bSPeter Brune   }
348eed5f15bSPeter Brune   ierr = PetscOptionsStringArray("-snes_composite_sneses","List of composite solvers","SNESCompositeAddSNES",sneses,&nmax,&flg);CHKERRQ(ierr);
349eed5f15bSPeter Brune   if (flg) {
350eed5f15bSPeter Brune     for (i=0; i<nmax; i++) {
351eed5f15bSPeter Brune       ierr = SNESCompositeAddSNES(snes,sneses[i]);CHKERRQ(ierr);
352eed5f15bSPeter Brune       ierr = PetscFree(sneses[i]);CHKERRQ(ierr);   /* deallocate string sneses[i], which is allocated in PetscOptionsStringArray() */
353eed5f15bSPeter Brune     }
354eed5f15bSPeter Brune   }
3558f626970SPeter Brune   ierr = PetscOptionsRealArray("-snes_composite_damping","Damping of the additive composite solvers","SNESCompositeSetDamping",dmps,&nmax,&flg);CHKERRQ(ierr);
3568f626970SPeter Brune   if (flg) {
3578f626970SPeter Brune     for (i=0; i<nmax; i++) {
3588f626970SPeter Brune       ierr = SNESCompositeSetDamping(snes,i,dmps[i]);CHKERRQ(ierr);
3598f626970SPeter Brune     }
3608f626970SPeter Brune   }
361eed5f15bSPeter Brune   ierr = PetscOptionsTail();CHKERRQ(ierr);
362eed5f15bSPeter Brune 
363eed5f15bSPeter Brune   next = jac->head;
364eed5f15bSPeter Brune   while (next) {
365eed5f15bSPeter Brune     ierr = SNESSetFromOptions(next->snes);CHKERRQ(ierr);
366eed5f15bSPeter Brune     next = next->next;
367eed5f15bSPeter Brune   }
368eed5f15bSPeter Brune   PetscFunctionReturn(0);
369eed5f15bSPeter Brune }
370eed5f15bSPeter Brune 
371eed5f15bSPeter Brune #undef __FUNCT__
372eed5f15bSPeter Brune #define __FUNCT__ "SNESView_Composite"
373eed5f15bSPeter Brune static PetscErrorCode SNESView_Composite(SNES snes,PetscViewer viewer)
374eed5f15bSPeter Brune {
375eed5f15bSPeter Brune   SNES_Composite     *jac = (SNES_Composite*)snes->data;
376eed5f15bSPeter Brune   PetscErrorCode   ierr;
377eed5f15bSPeter Brune   SNES_CompositeLink next = jac->head;
378eed5f15bSPeter Brune   PetscBool        iascii;
379eed5f15bSPeter Brune 
380eed5f15bSPeter Brune   PetscFunctionBegin;
381eed5f15bSPeter Brune   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
382eed5f15bSPeter Brune   if (iascii) {
383eed5f15bSPeter Brune     ierr = PetscViewerASCIIPrintf(viewer,"Composite SNES type - %s\n",SNESCompositeTypes[jac->type]);CHKERRQ(ierr);
384eed5f15bSPeter Brune     ierr = PetscViewerASCIIPrintf(viewer,"SNESes on composite preconditioner follow\n");CHKERRQ(ierr);
385eed5f15bSPeter Brune     ierr = PetscViewerASCIIPrintf(viewer,"---------------------------------\n");CHKERRQ(ierr);
386eed5f15bSPeter Brune   }
387eed5f15bSPeter Brune   if (iascii) {
388eed5f15bSPeter Brune     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
389eed5f15bSPeter Brune   }
390eed5f15bSPeter Brune   while (next) {
391eed5f15bSPeter Brune     ierr = SNESView(next->snes,viewer);CHKERRQ(ierr);
392eed5f15bSPeter Brune     next = next->next;
393eed5f15bSPeter Brune   }
394eed5f15bSPeter Brune   if (iascii) {
395eed5f15bSPeter Brune     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
396eed5f15bSPeter Brune     ierr = PetscViewerASCIIPrintf(viewer,"---------------------------------\n");CHKERRQ(ierr);
397eed5f15bSPeter Brune   }
398eed5f15bSPeter Brune   PetscFunctionReturn(0);
399eed5f15bSPeter Brune }
400eed5f15bSPeter Brune 
401eed5f15bSPeter Brune /* ------------------------------------------------------------------------------*/
402eed5f15bSPeter Brune 
403eed5f15bSPeter Brune #undef __FUNCT__
404eed5f15bSPeter Brune #define __FUNCT__ "SNESCompositeSetType_Composite"
405eed5f15bSPeter Brune static PetscErrorCode  SNESCompositeSetType_Composite(SNES snes,SNESCompositeType type)
406eed5f15bSPeter Brune {
407eed5f15bSPeter Brune   SNES_Composite *jac = (SNES_Composite*)snes->data;
408eed5f15bSPeter Brune 
409eed5f15bSPeter Brune   PetscFunctionBegin;
410eed5f15bSPeter Brune   jac->type = type;
411eed5f15bSPeter Brune   PetscFunctionReturn(0);
412eed5f15bSPeter Brune }
413eed5f15bSPeter Brune 
414eed5f15bSPeter Brune #undef __FUNCT__
415eed5f15bSPeter Brune #define __FUNCT__ "SNESCompositeAddSNES_Composite"
416eed5f15bSPeter Brune static PetscErrorCode  SNESCompositeAddSNES_Composite(SNES snes,SNESType type)
417eed5f15bSPeter Brune {
418eed5f15bSPeter Brune   SNES_Composite     *jac;
419eed5f15bSPeter Brune   SNES_CompositeLink next,ilink;
420eed5f15bSPeter Brune   PetscErrorCode   ierr;
421eed5f15bSPeter Brune   PetscInt         cnt = 0;
422eed5f15bSPeter Brune   const char       *prefix;
423eed5f15bSPeter Brune   char             newprefix[8];
424eed5f15bSPeter Brune   DM               dm;
425eed5f15bSPeter Brune 
426eed5f15bSPeter Brune   PetscFunctionBegin;
427eed5f15bSPeter Brune   ierr        = PetscNewLog(snes,struct _SNES_CompositeLink,&ilink);CHKERRQ(ierr);
428eed5f15bSPeter Brune   ilink->next = 0;
429eed5f15bSPeter Brune   ierr        = SNESCreate(PetscObjectComm((PetscObject)snes),&ilink->snes);CHKERRQ(ierr);
430eed5f15bSPeter Brune   ierr        = PetscLogObjectParent((PetscObject)snes,(PetscObject)ilink->snes);CHKERRQ(ierr);
431eed5f15bSPeter Brune   ierr        = SNESGetDM(snes,&dm);CHKERRQ(ierr);
432eed5f15bSPeter Brune   ierr        = SNESSetDM(ilink->snes,dm);CHKERRQ(ierr);
433eed5f15bSPeter Brune 
434eed5f15bSPeter Brune   jac  = (SNES_Composite*)snes->data;
435eed5f15bSPeter Brune   next = jac->head;
436eed5f15bSPeter Brune   if (!next) {
437eed5f15bSPeter Brune     jac->head       = ilink;
438eed5f15bSPeter Brune     ilink->previous = NULL;
439eed5f15bSPeter Brune   } else {
440eed5f15bSPeter Brune     cnt++;
441eed5f15bSPeter Brune     while (next->next) {
442eed5f15bSPeter Brune       next = next->next;
443eed5f15bSPeter Brune       cnt++;
444eed5f15bSPeter Brune     }
445eed5f15bSPeter Brune     next->next      = ilink;
446eed5f15bSPeter Brune     ilink->previous = next;
447eed5f15bSPeter Brune   }
448eed5f15bSPeter Brune   ierr = SNESGetOptionsPrefix(snes,&prefix);CHKERRQ(ierr);
449eed5f15bSPeter Brune   ierr = SNESSetOptionsPrefix(ilink->snes,prefix);CHKERRQ(ierr);
450eed5f15bSPeter Brune   sprintf(newprefix,"sub_%d_",(int)cnt);
451eed5f15bSPeter Brune   ierr = SNESAppendOptionsPrefix(ilink->snes,newprefix);CHKERRQ(ierr);
4528f626970SPeter Brune   ierr = PetscObjectIncrementTabLevel((PetscObject)ilink->snes,(PetscObject)snes,1);CHKERRQ(ierr);
453eed5f15bSPeter Brune   ierr = SNESSetType(ilink->snes,type);CHKERRQ(ierr);
454*9c2f3473SPeter Brune   ierr = SNESSetNormType(ilink->snes, SNES_NORM_FINAL_ONLY);CHKERRQ(ierr);
4558f626970SPeter Brune   ilink->dmp = 1.0;
45690a8ba9bSPeter Brune   jac->nsnes++;
457eed5f15bSPeter Brune   PetscFunctionReturn(0);
458eed5f15bSPeter Brune }
459eed5f15bSPeter Brune 
460eed5f15bSPeter Brune #undef __FUNCT__
461eed5f15bSPeter Brune #define __FUNCT__ "SNESCompositeGetSNES_Composite"
462eed5f15bSPeter Brune static PetscErrorCode  SNESCompositeGetSNES_Composite(SNES snes,PetscInt n,SNES *subsnes)
463eed5f15bSPeter Brune {
464eed5f15bSPeter Brune   SNES_Composite     *jac;
465eed5f15bSPeter Brune   SNES_CompositeLink next;
466eed5f15bSPeter Brune   PetscInt         i;
467eed5f15bSPeter Brune 
468eed5f15bSPeter Brune   PetscFunctionBegin;
469eed5f15bSPeter Brune   jac  = (SNES_Composite*)snes->data;
470eed5f15bSPeter Brune   next = jac->head;
471eed5f15bSPeter Brune   for (i=0; i<n; i++) {
472eed5f15bSPeter Brune     if (!next->next) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_INCOMP,"Not enough SNESes in composite preconditioner");
473eed5f15bSPeter Brune     next = next->next;
474eed5f15bSPeter Brune   }
475eed5f15bSPeter Brune   *subsnes = next->snes;
476eed5f15bSPeter Brune   PetscFunctionReturn(0);
477eed5f15bSPeter Brune }
478eed5f15bSPeter Brune 
479eed5f15bSPeter Brune /* -------------------------------------------------------------------------------- */
480eed5f15bSPeter Brune #undef __FUNCT__
481eed5f15bSPeter Brune #define __FUNCT__ "SNESCompositeSetType"
482e31ab4f9SPeter Brune /*@C
483eed5f15bSPeter Brune    SNESCompositeSetType - Sets the type of composite preconditioner.
484eed5f15bSPeter Brune 
485eed5f15bSPeter Brune    Logically Collective on SNES
486eed5f15bSPeter Brune 
487eed5f15bSPeter Brune    Input Parameter:
488eed5f15bSPeter Brune +  snes - the preconditioner context
489eed5f15bSPeter Brune -  type - SNES_COMPOSITE_ADDITIVE (default), SNES_COMPOSITE_MULTIPLICATIVE
490eed5f15bSPeter Brune 
491eed5f15bSPeter Brune    Options Database Key:
492eed5f15bSPeter Brune .  -snes_composite_type <type: one of multiplicative, additive, special> - Sets composite preconditioner type
493eed5f15bSPeter Brune 
494eed5f15bSPeter Brune    Level: Developer
495eed5f15bSPeter Brune 
496eed5f15bSPeter Brune .keywords: SNES, set, type, composite preconditioner, additive, multiplicative
497eed5f15bSPeter Brune @*/
498eed5f15bSPeter Brune PetscErrorCode  SNESCompositeSetType(SNES snes,SNESCompositeType type)
499eed5f15bSPeter Brune {
500eed5f15bSPeter Brune   PetscErrorCode ierr;
501eed5f15bSPeter Brune 
502eed5f15bSPeter Brune   PetscFunctionBegin;
503eed5f15bSPeter Brune   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
504eed5f15bSPeter Brune   PetscValidLogicalCollectiveEnum(snes,type,2);
505eed5f15bSPeter Brune   ierr = PetscTryMethod(snes,"SNESCompositeSetType_C",(SNES,SNESCompositeType),(snes,type));CHKERRQ(ierr);
506eed5f15bSPeter Brune   PetscFunctionReturn(0);
507eed5f15bSPeter Brune }
508eed5f15bSPeter Brune 
509eed5f15bSPeter Brune #undef __FUNCT__
510eed5f15bSPeter Brune #define __FUNCT__ "SNESCompositeAddSNES"
511eed5f15bSPeter Brune /*@C
512eed5f15bSPeter Brune    SNESCompositeAddSNES - Adds another SNES to the composite SNES.
513eed5f15bSPeter Brune 
514eed5f15bSPeter Brune    Collective on SNES
515eed5f15bSPeter Brune 
516eed5f15bSPeter Brune    Input Parameters:
517eed5f15bSPeter Brune +  snes - the preconditioner context
518eed5f15bSPeter Brune -  type - the type of the new preconditioner
519eed5f15bSPeter Brune 
520eed5f15bSPeter Brune    Level: Developer
521eed5f15bSPeter Brune 
522eed5f15bSPeter Brune .keywords: SNES, composite preconditioner, add
523eed5f15bSPeter Brune @*/
524eed5f15bSPeter Brune PetscErrorCode  SNESCompositeAddSNES(SNES snes,SNESType type)
525eed5f15bSPeter Brune {
526eed5f15bSPeter Brune   PetscErrorCode ierr;
527eed5f15bSPeter Brune 
528eed5f15bSPeter Brune   PetscFunctionBegin;
529eed5f15bSPeter Brune   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
530eed5f15bSPeter Brune   ierr = PetscTryMethod(snes,"SNESCompositeAddSNES_C",(SNES,SNESType),(snes,type));CHKERRQ(ierr);
531eed5f15bSPeter Brune   PetscFunctionReturn(0);
532eed5f15bSPeter Brune }
533eed5f15bSPeter Brune #undef __FUNCT__
534eed5f15bSPeter Brune #define __FUNCT__ "SNESCompositeGetSNES"
535eed5f15bSPeter Brune /*@
536eed5f15bSPeter Brune    SNESCompositeGetSNES - Gets one of the SNES objects in the composite SNES.
537eed5f15bSPeter Brune 
538eed5f15bSPeter Brune    Not Collective
539eed5f15bSPeter Brune 
540eed5f15bSPeter Brune    Input Parameter:
541eed5f15bSPeter Brune +  snes - the preconditioner context
542eed5f15bSPeter Brune -  n - the number of the snes requested
543eed5f15bSPeter Brune 
544eed5f15bSPeter Brune    Output Parameters:
545eed5f15bSPeter Brune .  subsnes - the SNES requested
546eed5f15bSPeter Brune 
547eed5f15bSPeter Brune    Level: Developer
548eed5f15bSPeter Brune 
549eed5f15bSPeter Brune .keywords: SNES, get, composite preconditioner, sub preconditioner
550eed5f15bSPeter Brune 
551eed5f15bSPeter Brune .seealso: SNESCompositeAddSNES()
552eed5f15bSPeter Brune @*/
553eed5f15bSPeter Brune PetscErrorCode  SNESCompositeGetSNES(SNES snes,PetscInt n,SNES *subsnes)
554eed5f15bSPeter Brune {
555eed5f15bSPeter Brune   PetscErrorCode ierr;
556eed5f15bSPeter Brune 
557eed5f15bSPeter Brune   PetscFunctionBegin;
558eed5f15bSPeter Brune   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
559eed5f15bSPeter Brune   PetscValidPointer(subsnes,3);
560eed5f15bSPeter Brune   ierr = PetscUseMethod(snes,"SNESCompositeGetSNES_C",(SNES,PetscInt,SNES*),(snes,n,subsnes));CHKERRQ(ierr);
561eed5f15bSPeter Brune   PetscFunctionReturn(0);
562eed5f15bSPeter Brune }
563eed5f15bSPeter Brune 
564eed5f15bSPeter Brune #undef __FUNCT__
5658f626970SPeter Brune #define __FUNCT__ "SNESCompositeSetDamping_Composite"
5668f626970SPeter Brune static PetscErrorCode  SNESCompositeSetDamping_Composite(SNES snes,PetscInt n,PetscReal dmp)
5678f626970SPeter Brune {
5688f626970SPeter Brune   SNES_Composite     *jac;
5698f626970SPeter Brune   SNES_CompositeLink next;
5708f626970SPeter Brune   PetscInt         i;
5718f626970SPeter Brune 
5728f626970SPeter Brune   PetscFunctionBegin;
5738f626970SPeter Brune   jac  = (SNES_Composite*)snes->data;
5748f626970SPeter Brune   next = jac->head;
5758f626970SPeter Brune   for (i=0; i<n; i++) {
5768f626970SPeter Brune     if (!next->next) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_INCOMP,"Not enough SNESes in composite preconditioner");
5778f626970SPeter Brune     next = next->next;
5788f626970SPeter Brune   }
5798f626970SPeter Brune   next->dmp = dmp;
5808f626970SPeter Brune   PetscFunctionReturn(0);
5818f626970SPeter Brune }
5828f626970SPeter Brune 
5838f626970SPeter Brune #undef __FUNCT__
5848f626970SPeter Brune #define __FUNCT__ "SNESCompositeSetDamping"
5858f626970SPeter Brune /*@
5868f626970SPeter Brune    SNESCompositeSetDamping - Sets the damping of a subsolver when using additive composite SNES.
5878f626970SPeter Brune 
5888f626970SPeter Brune    Not Collective
5898f626970SPeter Brune 
5908f626970SPeter Brune    Input Parameter:
5918f626970SPeter Brune +  snes - the preconditioner context
5928f626970SPeter Brune .  n - the number of the snes requested
5938f626970SPeter Brune -  dmp - the damping
5948f626970SPeter Brune 
5958f626970SPeter Brune    Level: Developer
5968f626970SPeter Brune 
5978f626970SPeter Brune .keywords: SNES, get, composite preconditioner, sub preconditioner
5988f626970SPeter Brune 
5998f626970SPeter Brune .seealso: SNESCompositeAddSNES()
6008f626970SPeter Brune @*/
6018f626970SPeter Brune PetscErrorCode  SNESCompositeSetDamping(SNES snes,PetscInt n,PetscReal dmp)
6028f626970SPeter Brune {
6038f626970SPeter Brune   PetscErrorCode ierr;
6048f626970SPeter Brune 
6058f626970SPeter Brune   PetscFunctionBegin;
6068f626970SPeter Brune   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
6078f626970SPeter Brune   ierr = PetscUseMethod(snes,"SNESCompositeSetDamping_C",(SNES,PetscInt,PetscReal),(snes,n,dmp));CHKERRQ(ierr);
6088f626970SPeter Brune   PetscFunctionReturn(0);
6098f626970SPeter Brune }
6108f626970SPeter Brune 
6118f626970SPeter Brune #undef __FUNCT__
612eed5f15bSPeter Brune #define __FUNCT__ "SNESSolve_Composite"
613eed5f15bSPeter Brune PetscErrorCode SNESSolve_Composite(SNES snes)
614eed5f15bSPeter Brune {
615eed5f15bSPeter Brune   Vec            F;
616eed5f15bSPeter Brune   Vec            X;
617eed5f15bSPeter Brune   Vec            B;
618eed5f15bSPeter Brune   PetscInt       i;
619eed5f15bSPeter Brune   PetscReal      fnorm = 0.0;
620eed5f15bSPeter Brune   PetscErrorCode ierr;
621eed5f15bSPeter Brune   SNESNormType   normtype;
622eed5f15bSPeter Brune   SNES_Composite *comp = (SNES_Composite*)snes->data;
623eed5f15bSPeter Brune 
624eed5f15bSPeter Brune   PetscFunctionBegin;
625eed5f15bSPeter Brune   X = snes->vec_sol;
626eed5f15bSPeter Brune   F = snes->vec_func;
627eed5f15bSPeter Brune   B = snes->vec_rhs;
628eed5f15bSPeter Brune 
629eed5f15bSPeter Brune   ierr         = PetscObjectAMSTakeAccess((PetscObject)snes);CHKERRQ(ierr);
630eed5f15bSPeter Brune   snes->iter   = 0;
631eed5f15bSPeter Brune   snes->norm   = 0.;
632eed5f15bSPeter Brune   ierr         = PetscObjectAMSGrantAccess((PetscObject)snes);CHKERRQ(ierr);
633eed5f15bSPeter Brune   snes->reason = SNES_CONVERGED_ITERATING;
634eed5f15bSPeter Brune   ierr         = SNESGetNormType(snes, &normtype);CHKERRQ(ierr);
635eed5f15bSPeter Brune   if (normtype == SNES_NORM_FUNCTION || normtype == SNES_NORM_INITIAL_ONLY || normtype == SNES_NORM_INITIAL_FINAL_ONLY) {
636eed5f15bSPeter Brune     if (!snes->vec_func_init_set) {
637eed5f15bSPeter Brune       ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
638eed5f15bSPeter Brune       if (snes->domainerror) {
639eed5f15bSPeter Brune         snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
640eed5f15bSPeter Brune         PetscFunctionReturn(0);
641eed5f15bSPeter Brune       }
642eed5f15bSPeter Brune     } else snes->vec_func_init_set = PETSC_FALSE;
643eed5f15bSPeter Brune 
644eed5f15bSPeter Brune     /* convergence test */
645eed5f15bSPeter Brune     if (!snes->norm_init_set) {
646eed5f15bSPeter Brune       ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F||  */
647eed5f15bSPeter Brune       if (PetscIsInfOrNanReal(fnorm)) {
648eed5f15bSPeter Brune         snes->reason = SNES_DIVERGED_FNORM_NAN;
649eed5f15bSPeter Brune         PetscFunctionReturn(0);
650eed5f15bSPeter Brune       }
651eed5f15bSPeter Brune     } else {
652eed5f15bSPeter Brune       fnorm               = snes->norm_init;
653eed5f15bSPeter Brune       snes->norm_init_set = PETSC_FALSE;
654eed5f15bSPeter Brune     }
655eed5f15bSPeter Brune     ierr       = PetscObjectAMSTakeAccess((PetscObject)snes);CHKERRQ(ierr);
656eed5f15bSPeter Brune     snes->iter = 0;
657eed5f15bSPeter Brune     snes->norm = fnorm;
658eed5f15bSPeter Brune     ierr       = PetscObjectAMSGrantAccess((PetscObject)snes);CHKERRQ(ierr);
659eed5f15bSPeter Brune     ierr       = SNESLogConvergenceHistory(snes,snes->norm,0);CHKERRQ(ierr);
660eed5f15bSPeter Brune     ierr       = SNESMonitor(snes,0,snes->norm);CHKERRQ(ierr);
661eed5f15bSPeter Brune 
662eed5f15bSPeter Brune     /* set parameter for default relative tolerance convergence test */
663eed5f15bSPeter Brune     snes->ttol = fnorm*snes->rtol;
664eed5f15bSPeter Brune 
665eed5f15bSPeter Brune     /* test convergence */
666eed5f15bSPeter Brune     ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
667eed5f15bSPeter Brune     if (snes->reason) PetscFunctionReturn(0);
668eed5f15bSPeter Brune   } else {
669eed5f15bSPeter Brune     ierr = PetscObjectAMSGrantAccess((PetscObject)snes);CHKERRQ(ierr);
670eed5f15bSPeter Brune     ierr = SNESLogConvergenceHistory(snes,snes->norm,0);CHKERRQ(ierr);
671eed5f15bSPeter Brune     ierr = SNESMonitor(snes,0,snes->norm);CHKERRQ(ierr);
672eed5f15bSPeter Brune   }
673eed5f15bSPeter Brune 
674eed5f15bSPeter Brune   /* Call general purpose update function */
675eed5f15bSPeter Brune   if (snes->ops->update) {
676eed5f15bSPeter Brune     ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
677eed5f15bSPeter Brune   }
678eed5f15bSPeter Brune 
679eed5f15bSPeter Brune   for (i = 0; i < snes->max_its; i++) {
680eed5f15bSPeter Brune     if (comp->type == SNES_COMPOSITE_ADDITIVE) {
681eed5f15bSPeter Brune       ierr = SNESCompositeApply_Additive(snes,X,B,F,&fnorm);CHKERRQ(ierr);
682eed5f15bSPeter Brune     } else if (comp->type == SNES_COMPOSITE_MULTIPLICATIVE) {
683eed5f15bSPeter Brune       ierr = SNESCompositeApply_Multiplicative(snes,X,B,F,&fnorm);CHKERRQ(ierr);
68490a8ba9bSPeter Brune     } else if (comp->type == SNES_COMPOSITE_ADDITIVEOPTIMAL) {
68590a8ba9bSPeter Brune       ierr = SNESCompositeApply_AdditiveOptimal(snes,X,B,F,&fnorm);CHKERRQ(ierr);
686eed5f15bSPeter Brune     } else {
68790a8ba9bSPeter Brune       SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"Unsupported SNESComposite type");
688eed5f15bSPeter Brune     }
689eed5f15bSPeter Brune     if ((i == snes->max_its - 1) && (normtype == SNES_NORM_INITIAL_FINAL_ONLY || normtype == SNES_NORM_FINAL_ONLY)) {
690eed5f15bSPeter Brune       ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
691eed5f15bSPeter Brune       if (snes->domainerror) {
692eed5f15bSPeter Brune         snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
693eed5f15bSPeter Brune         break;
694eed5f15bSPeter Brune       }
695eed5f15bSPeter Brune       ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F||  */
696eed5f15bSPeter Brune       if (PetscIsInfOrNanReal(fnorm)) {
697eed5f15bSPeter Brune         snes->reason = SNES_DIVERGED_FNORM_NAN;
698eed5f15bSPeter Brune         break;
699eed5f15bSPeter Brune       }
700eed5f15bSPeter Brune     }
701eed5f15bSPeter Brune     /* Monitor convergence */
702eed5f15bSPeter Brune     ierr       = PetscObjectAMSTakeAccess((PetscObject)snes);CHKERRQ(ierr);
703eed5f15bSPeter Brune     snes->iter = i+1;
704eed5f15bSPeter Brune     snes->norm = fnorm;
705eed5f15bSPeter Brune     ierr       = PetscObjectAMSGrantAccess((PetscObject)snes);CHKERRQ(ierr);
706eed5f15bSPeter Brune     ierr       = SNESLogConvergenceHistory(snes,snes->norm,0);CHKERRQ(ierr);
707eed5f15bSPeter Brune     ierr       = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr);
708eed5f15bSPeter Brune     /* Test for convergence */
709eed5f15bSPeter Brune     if (normtype == SNES_NORM_FUNCTION) {ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);}
710eed5f15bSPeter Brune     if (snes->reason) break;
711eed5f15bSPeter Brune     /* Call general purpose update function */
712eed5f15bSPeter Brune     if (snes->ops->update) {ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);}
713eed5f15bSPeter Brune   }
714eed5f15bSPeter Brune   if (normtype == SNES_NORM_FUNCTION) {
715eed5f15bSPeter Brune     if (i == snes->max_its) {
716eed5f15bSPeter Brune       ierr = PetscInfo1(snes,"Maximum number of iterations has been reached: %D\n",snes->max_its);CHKERRQ(ierr);
717eed5f15bSPeter Brune       if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
718eed5f15bSPeter Brune     }
719eed5f15bSPeter Brune   } else if (!snes->reason) snes->reason = SNES_CONVERGED_ITS;
720eed5f15bSPeter Brune   PetscFunctionReturn(0);
721eed5f15bSPeter Brune }
722eed5f15bSPeter Brune 
723eed5f15bSPeter Brune /* -------------------------------------------------------------------------------------------*/
724eed5f15bSPeter Brune 
725eed5f15bSPeter Brune /*MC
726eed5f15bSPeter Brune      SNESCOMPOSITE - Build a preconditioner by composing together several nonlinear solvers
727eed5f15bSPeter Brune 
728eed5f15bSPeter Brune    Options Database Keys:
729eed5f15bSPeter Brune +  -snes_composite_type <type: one of multiplicative, additive, symmetric_multiplicative, special> - Sets composite preconditioner type
730eed5f15bSPeter Brune -  -snes_composite_sneses - <snes0,snes1,...> list of SNESes to compose
731eed5f15bSPeter Brune 
732eed5f15bSPeter Brune    Level: intermediate
733eed5f15bSPeter Brune 
734eed5f15bSPeter Brune    Concepts: composing solvers
735eed5f15bSPeter Brune 
736eed5f15bSPeter Brune .seealso:  SNESCreate(), SNESSetType(), SNESType (for list of available types), SNES,
737eed5f15bSPeter Brune            SNESSHELL, SNESCompositeSetType(), SNESCompositeSpecialSetAlpha(), SNESCompositeAddSNES(),
738eed5f15bSPeter Brune            SNESCompositeGetSNES()
739eed5f15bSPeter Brune 
740eed5f15bSPeter Brune M*/
741eed5f15bSPeter Brune 
742eed5f15bSPeter Brune #undef __FUNCT__
743eed5f15bSPeter Brune #define __FUNCT__ "SNESCreate_Composite"
744eed5f15bSPeter Brune PETSC_EXTERN PetscErrorCode SNESCreate_Composite(SNES snes)
745eed5f15bSPeter Brune {
746eed5f15bSPeter Brune   PetscErrorCode ierr;
747eed5f15bSPeter Brune   SNES_Composite   *jac;
748eed5f15bSPeter Brune 
749eed5f15bSPeter Brune   PetscFunctionBegin;
750eed5f15bSPeter Brune   ierr = PetscNewLog(snes,SNES_Composite,&jac);CHKERRQ(ierr);
751eed5f15bSPeter Brune 
752eed5f15bSPeter Brune   snes->ops->solve           = SNESSolve_Composite;
753eed5f15bSPeter Brune   snes->ops->setup           = SNESSetUp_Composite;
754eed5f15bSPeter Brune   snes->ops->reset           = SNESReset_Composite;
755eed5f15bSPeter Brune   snes->ops->destroy         = SNESDestroy_Composite;
756eed5f15bSPeter Brune   snes->ops->setfromoptions  = SNESSetFromOptions_Composite;
757eed5f15bSPeter Brune   snes->ops->view            = SNESView_Composite;
758eed5f15bSPeter Brune 
759eed5f15bSPeter Brune   snes->data = (void*)jac;
76090a8ba9bSPeter Brune   jac->type  = SNES_COMPOSITE_ADDITIVEOPTIMAL;
76190a8ba9bSPeter Brune   jac->Fes   = NULL;
76290a8ba9bSPeter Brune   jac->Xes   = NULL;
763*9c2f3473SPeter Brune   jac->fnorms = NULL;
76490a8ba9bSPeter Brune   jac->nsnes = 0;
765eed5f15bSPeter Brune   jac->head  = 0;
766eed5f15bSPeter Brune 
76790a8ba9bSPeter Brune   jac->h     = NULL;
76890a8ba9bSPeter Brune   jac->s     = NULL;
76990a8ba9bSPeter Brune   jac->beta  = NULL;
77090a8ba9bSPeter Brune   jac->work  = NULL;
77190a8ba9bSPeter Brune   jac->rwork = NULL;
77290a8ba9bSPeter Brune 
7738f626970SPeter Brune   ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESCompositeSetType_C",SNESCompositeSetType_Composite);CHKERRQ(ierr);
7748f626970SPeter Brune   ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESCompositeAddSNES_C",SNESCompositeAddSNES_Composite);CHKERRQ(ierr);
7758f626970SPeter Brune   ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESCompositeGetSNES_C",SNESCompositeGetSNES_Composite);CHKERRQ(ierr);
7768f626970SPeter Brune   ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESCompositeSetDamping_C",SNESCompositeSetDamping_Composite);CHKERRQ(ierr);
777eed5f15bSPeter Brune   PetscFunctionReturn(0);
778eed5f15bSPeter Brune }
779eed5f15bSPeter Brune 
780