xref: /petsc/src/snes/impls/composite/snescomposite.c (revision 58bdfa14ad7d3f9f16e20f54039c675862e9487f)
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 */
279c2f3473SPeter Brune   PetscReal          *fnorms;        /* norms of the solutions */
2890a8ba9bSPeter Brune   PetscScalar        *h;             /* the matrix formed as q_ij = (rdot_i, rdot_j) */
299c2f3473SPeter Brune   PetscScalar        *g;             /* the dotproducts of the previous function with the candidate functions */
309c2f3473SPeter 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 */
359c2f3473SPeter 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 
435e806d2eSPeter Brune   PetscReal          rtol;           /* restart tolerance for accepting the combination */
445e806d2eSPeter Brune   PetscReal          stol;           /* restart tolerance for the combination */
45eed5f15bSPeter Brune } SNES_Composite;
46eed5f15bSPeter Brune 
47eed5f15bSPeter Brune #undef __FUNCT__
48eed5f15bSPeter Brune #define __FUNCT__ "SNESCompositeApply_Multiplicative"
49eed5f15bSPeter Brune static PetscErrorCode SNESCompositeApply_Multiplicative(SNES snes,Vec X,Vec B,Vec F,PetscReal *fnorm)
50eed5f15bSPeter Brune {
51eed5f15bSPeter Brune   PetscErrorCode     ierr;
52eed5f15bSPeter Brune   SNES_Composite     *jac = (SNES_Composite*)snes->data;
53eed5f15bSPeter Brune   SNES_CompositeLink next = jac->head;
54eed5f15bSPeter Brune   Vec                FSub;
55eed5f15bSPeter Brune   PetscReal          fsubnorm;
56eed5f15bSPeter Brune 
57eed5f15bSPeter Brune   PetscFunctionBegin;
58eed5f15bSPeter Brune   if (!next) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"No composite SNESes supplied via SNESCompositeAddSNES() or -snes_composite_sneses");
5972edecb9SPeter Brune   if (snes->normschedule == SNES_NORM_ALWAYS) {
60eed5f15bSPeter Brune     ierr = SNESSetInitialFunction(next->snes,F);CHKERRQ(ierr);
61eed5f15bSPeter Brune     if (fnorm) {ierr = SNESSetInitialFunctionNorm(next->snes,*fnorm);CHKERRQ(ierr);}
62eed5f15bSPeter Brune   }
63eed5f15bSPeter Brune   ierr = SNESSolve(next->snes,B,X);CHKERRQ(ierr);
64eed5f15bSPeter Brune 
65eed5f15bSPeter Brune   while (next->next) {
66eed5f15bSPeter Brune     /* only copy the function over in the case where the functions correspond */
6772edecb9SPeter Brune     if (next->snes->pcside == PC_RIGHT && next->snes->normschedule != SNES_NORM_NONE) {
68eed5f15bSPeter Brune       ierr = SNESGetFunction(next->snes,&FSub,NULL,NULL);CHKERRQ(ierr);
69eed5f15bSPeter Brune       ierr = SNESGetFunctionNorm(next->snes,&fsubnorm);CHKERRQ(ierr);
70eed5f15bSPeter Brune       next = next->next;
71eed5f15bSPeter Brune       ierr = SNESSetInitialFunction(next->snes,FSub);CHKERRQ(ierr);
72eed5f15bSPeter Brune       ierr = SNESSetInitialFunctionNorm(next->snes,fsubnorm);CHKERRQ(ierr);
73eed5f15bSPeter Brune     } else {
74eed5f15bSPeter Brune       next = next->next;
75eed5f15bSPeter Brune     }
76eed5f15bSPeter Brune     ierr = SNESSolve(next->snes,B,X);CHKERRQ(ierr);
77eed5f15bSPeter Brune   }
78eed5f15bSPeter Brune   if (next->snes->pcside == PC_RIGHT) {
79eed5f15bSPeter Brune     ierr = SNESGetFunction(next->snes,&FSub,NULL,NULL);CHKERRQ(ierr);
80eed5f15bSPeter Brune     ierr = VecCopy(FSub,F);CHKERRQ(ierr);
81eed5f15bSPeter Brune     if (fnorm) {ierr = SNESGetFunctionNorm(next->snes,fnorm);CHKERRQ(ierr);}
8272edecb9SPeter Brune   } else if (snes->normschedule == SNES_NORM_ALWAYS) {
83eed5f15bSPeter Brune     SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
84eed5f15bSPeter Brune     if (fnorm) {ierr = VecNorm(F,NORM_2,fnorm);CHKERRQ(ierr);}
85eed5f15bSPeter Brune   }
86eed5f15bSPeter Brune   PetscFunctionReturn(0);
87eed5f15bSPeter Brune }
88eed5f15bSPeter Brune 
89eed5f15bSPeter Brune #undef __FUNCT__
90eed5f15bSPeter Brune #define __FUNCT__ "SNESCompositeApply_Additive"
91eed5f15bSPeter Brune static PetscErrorCode SNESCompositeApply_Additive(SNES snes,Vec X,Vec B,Vec F,PetscReal *fnorm)
92eed5f15bSPeter Brune {
93eed5f15bSPeter Brune   PetscErrorCode     ierr;
94eed5f15bSPeter Brune   SNES_Composite     *jac = (SNES_Composite*)snes->data;
95eed5f15bSPeter Brune   SNES_CompositeLink next = jac->head;
96eed5f15bSPeter Brune   Vec                Y,Xorig;
97eed5f15bSPeter Brune 
98eed5f15bSPeter Brune   PetscFunctionBegin;
99eed5f15bSPeter Brune   Y = snes->vec_sol_update;
100eed5f15bSPeter Brune   if (!jac->Xorig) {ierr = VecDuplicate(X,&jac->Xorig);CHKERRQ(ierr);}
101eed5f15bSPeter Brune   Xorig = jac->Xorig;
102eed5f15bSPeter Brune   ierr = VecCopy(X,Xorig);
103eed5f15bSPeter Brune   if (!next) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"No composite SNESes supplied via SNESCompositeAddSNES() or -snes_composite_sneses");
10472edecb9SPeter Brune   if (snes->normschedule == SNES_NORM_ALWAYS) {
105eed5f15bSPeter Brune     ierr = SNESSetInitialFunction(next->snes,F);CHKERRQ(ierr);
106eed5f15bSPeter Brune     if (fnorm) {ierr = SNESSetInitialFunctionNorm(next->snes,*fnorm);CHKERRQ(ierr);}
107eed5f15bSPeter Brune     while (next->next) {
108eed5f15bSPeter Brune       next = next->next;
109eed5f15bSPeter Brune       ierr = SNESSetInitialFunction(next->snes,F);CHKERRQ(ierr);
110eed5f15bSPeter Brune       if (fnorm) {ierr = SNESSetInitialFunctionNorm(next->snes,*fnorm);CHKERRQ(ierr);}
111eed5f15bSPeter Brune     }
112eed5f15bSPeter Brune   }
113eed5f15bSPeter Brune   next = jac->head;
1148f626970SPeter Brune   ierr = VecCopy(Xorig,Y);CHKERRQ(ierr);
115eed5f15bSPeter Brune   ierr = SNESSolve(next->snes,B,Y);CHKERRQ(ierr);
116eed5f15bSPeter Brune   ierr = VecAXPY(Y,-1.0,Xorig);CHKERRQ(ierr);
1178f626970SPeter Brune   ierr = VecAXPY(X,next->dmp,Y);CHKERRQ(ierr);
1188f626970SPeter Brune   while (next->next) {
1198f626970SPeter Brune     next = next->next;
1208f626970SPeter Brune     ierr = VecCopy(Xorig,Y);CHKERRQ(ierr);
1218f626970SPeter Brune     ierr = SNESSolve(next->snes,B,Y);CHKERRQ(ierr);
1228f626970SPeter Brune     ierr = VecAXPY(Y,-1.0,Xorig);CHKERRQ(ierr);
1238f626970SPeter Brune     ierr = VecAXPY(X,next->dmp,Y);CHKERRQ(ierr);
124eed5f15bSPeter Brune   }
12572edecb9SPeter Brune   if (snes->normschedule == SNES_NORM_ALWAYS) {
126eed5f15bSPeter Brune     ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
127eed5f15bSPeter Brune     if (fnorm) {ierr = VecNorm(F,NORM_2,fnorm);CHKERRQ(ierr);}
128eed5f15bSPeter Brune   }
129eed5f15bSPeter Brune   PetscFunctionReturn(0);
130eed5f15bSPeter Brune }
13190a8ba9bSPeter Brune 
13290a8ba9bSPeter Brune #undef __FUNCT__
13390a8ba9bSPeter Brune #define __FUNCT__ "SNESCompositeApply_AdditiveOptimal"
13490a8ba9bSPeter Brune /* approximately solve the overdetermined system:
13590a8ba9bSPeter Brune 
13690a8ba9bSPeter Brune  2*F(x_i)\cdot F(\x_j)\alpha_i = 0
13790a8ba9bSPeter Brune  \alpha_i                      = 1
13890a8ba9bSPeter Brune 
13990a8ba9bSPeter Brune  Which minimizes the L2 norm of the linearization of:
14090a8ba9bSPeter Brune  ||F(\sum_i \alpha_i*x_i)||^2
14190a8ba9bSPeter Brune 
14290a8ba9bSPeter Brune  With the constraint that \sum_i\alpha_i = 1
14390a8ba9bSPeter Brune  Where x_i is the solution from the ith subsolver.
14490a8ba9bSPeter Brune  */
14590a8ba9bSPeter Brune static PetscErrorCode SNESCompositeApply_AdditiveOptimal(SNES snes,Vec X,Vec B,Vec F,PetscReal *fnorm)
14690a8ba9bSPeter Brune {
14790a8ba9bSPeter Brune   PetscErrorCode     ierr;
14890a8ba9bSPeter Brune   SNES_Composite     *jac = (SNES_Composite*)snes->data;
14990a8ba9bSPeter Brune   SNES_CompositeLink next = jac->head;
15090a8ba9bSPeter Brune   Vec                *Xes = jac->Xes,*Fes = jac->Fes;
15190a8ba9bSPeter Brune   PetscInt           i,j;
1525e806d2eSPeter Brune   PetscScalar        tot,total,ftf;
1535e806d2eSPeter Brune   PetscReal          min_fnorm;
1545e806d2eSPeter Brune   PetscInt           min_i;
15590a8ba9bSPeter Brune 
15690a8ba9bSPeter Brune   PetscFunctionBegin;
15790a8ba9bSPeter Brune   if (!next) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"No composite SNESes supplied via SNESCompositeAddSNES() or -snes_composite_sneses");
158*58bdfa14SPeter Brune 
159*58bdfa14SPeter Brune   if (snes->normschedule == SNES_NORM_ALWAYS) {
160*58bdfa14SPeter Brune     next = jac->head;
161*58bdfa14SPeter Brune     ierr = SNESSetInitialFunction(next->snes,F);CHKERRQ(ierr);
162*58bdfa14SPeter Brune     if (fnorm) {ierr = SNESSetInitialFunctionNorm(next->snes,*fnorm);CHKERRQ(ierr);}
163*58bdfa14SPeter Brune     while (next->next) {
164*58bdfa14SPeter Brune       next = next->next;
165*58bdfa14SPeter Brune       ierr = SNESSetInitialFunction(next->snes,F);CHKERRQ(ierr);
166*58bdfa14SPeter Brune       if (fnorm) {ierr = SNESSetInitialFunctionNorm(next->snes,*fnorm);CHKERRQ(ierr);}
167*58bdfa14SPeter Brune     }
168*58bdfa14SPeter Brune   }
169*58bdfa14SPeter Brune 
17090a8ba9bSPeter Brune   next = jac->head;
17190a8ba9bSPeter Brune   i = 0;
17290a8ba9bSPeter Brune   ierr = VecCopy(X,Xes[i]);CHKERRQ(ierr);
17390a8ba9bSPeter Brune   ierr = SNESSolve(next->snes,B,Xes[i]);CHKERRQ(ierr);
17490a8ba9bSPeter Brune   while (next->next) {
17590a8ba9bSPeter Brune     i++;
17690a8ba9bSPeter Brune     next = next->next;
17790a8ba9bSPeter Brune     ierr = VecCopy(X,Xes[i]);CHKERRQ(ierr);
17890a8ba9bSPeter Brune     ierr = SNESSolve(next->snes,B,Xes[i]);CHKERRQ(ierr);
17990a8ba9bSPeter Brune   }
18090a8ba9bSPeter Brune 
18190a8ba9bSPeter Brune   /* all the solutions are collected; combine optimally */
18290a8ba9bSPeter Brune   for (i=0;i<jac->n;i++) {
18390a8ba9bSPeter Brune     for (j=0;j<i+1;j++) {
1849c2f3473SPeter Brune       ierr = VecDotBegin(Fes[i],Fes[j],&jac->h[i + j*jac->n]);CHKERRQ(ierr);
18590a8ba9bSPeter Brune     }
1869c2f3473SPeter Brune     ierr = VecDotBegin(Fes[i],F,&jac->g[i]);CHKERRQ(ierr);
18790a8ba9bSPeter Brune   }
1885e806d2eSPeter Brune 
18990a8ba9bSPeter Brune   for (i=0;i<jac->n;i++) {
19090a8ba9bSPeter Brune     for (j=0;j<i+1;j++) {
1919c2f3473SPeter Brune       ierr = VecDotEnd(Fes[i],Fes[j],&jac->h[i + j*jac->n]);CHKERRQ(ierr);
1929c2f3473SPeter Brune       if (i == j) jac->fnorms[i] = PetscSqrtReal(PetscRealPart(jac->h[i + j*jac->n]));
19390a8ba9bSPeter Brune     }
1949c2f3473SPeter Brune     ierr = VecDotEnd(Fes[i],F,&jac->g[i]);CHKERRQ(ierr);
19590a8ba9bSPeter Brune   }
19690a8ba9bSPeter Brune 
1979c2f3473SPeter Brune   ftf = (*fnorm)*(*fnorm);
19890a8ba9bSPeter Brune 
19990a8ba9bSPeter Brune   for (i=0; i<jac->n; i++) {
20090a8ba9bSPeter Brune     for (j=i+1;j<jac->n;j++) {
2019c2f3473SPeter Brune       jac->h[i + j*jac->n] = jac->h[j + i*jac->n];
20290a8ba9bSPeter Brune     }
20390a8ba9bSPeter Brune   }
20490a8ba9bSPeter Brune 
20590a8ba9bSPeter Brune   for (i=0; i<jac->n; i++) {
2069c2f3473SPeter Brune     for (j=0; j<jac->n; j++) {
2079c2f3473SPeter Brune       jac->h[i + j*jac->n] = jac->h[i + j*jac->n] - jac->g[j] - jac->g[i] + ftf;
20890a8ba9bSPeter Brune     }
2099c2f3473SPeter Brune     jac->beta[i] = ftf - jac->g[i];
2109c2f3473SPeter Brune   }
21190a8ba9bSPeter Brune 
21290a8ba9bSPeter Brune #if defined(PETSC_MISSING_LAPACK_GELSS)
21390a8ba9bSPeter Brune   SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_SUP,"SNESCOMPOSITE with ADDITIVEOPTIMAL requires the LAPACK GELSS routine.");
21490a8ba9bSPeter Brune #else
21590a8ba9bSPeter Brune   jac->info  = 0;
21690a8ba9bSPeter Brune   jac->rcond = -1.;
21790a8ba9bSPeter Brune   ierr          = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
21890a8ba9bSPeter Brune #if defined(PETSC_USE_COMPLEX)
2199c2f3473SPeter 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));
22090a8ba9bSPeter Brune #else
2219c2f3473SPeter 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));
22290a8ba9bSPeter Brune #endif
22390a8ba9bSPeter Brune   ierr = PetscFPTrapPop();CHKERRQ(ierr);
22490a8ba9bSPeter Brune   if (jac->info < 0) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_LIB,"Bad argument to GELSS");
22590a8ba9bSPeter Brune   if (jac->info > 0) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_LIB,"SVD failed to converge");
22690a8ba9bSPeter Brune #endif
2279c2f3473SPeter Brune   tot = 0.;
2285e806d2eSPeter Brune   total = 0.;
22990a8ba9bSPeter Brune   for (i=0; i<jac->n; i++) {
23090a8ba9bSPeter Brune     if (PetscIsInfOrNanScalar(jac->beta[i])) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_LIB,"SVD generated inconsistent output");
231b3000dc2SPeter Brune     ierr = PetscInfo2(snes,"%d: %f\n",i,PetscRealPart(jac->beta[i]));CHKERRQ(ierr);
2329c2f3473SPeter Brune     tot += jac->beta[i];
2335e806d2eSPeter Brune     total += PetscAbsScalar(jac->beta[i]);
23490a8ba9bSPeter Brune   }
2359c2f3473SPeter Brune   ierr = VecScale(X,(1. - tot));CHKERRQ(ierr);
23690a8ba9bSPeter Brune   ierr = VecMAXPY(X,jac->n,jac->beta,Xes);CHKERRQ(ierr);
23790a8ba9bSPeter Brune   ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
2389c2f3473SPeter Brune   ierr = VecNorm(F,NORM_2,fnorm);CHKERRQ(ierr);
23990a8ba9bSPeter Brune 
2405e806d2eSPeter Brune   /* take the minimum-normed candidate if it beats the combination by a factor of rtol or the combination has stagnated */
2415e806d2eSPeter Brune   min_fnorm = jac->fnorms[0];
2425e806d2eSPeter Brune   min_i     = 0;
2439c2f3473SPeter Brune   for (i=0; i<jac->n; i++) {
2445e806d2eSPeter Brune     if (jac->fnorms[i] < min_fnorm) {
2455e806d2eSPeter Brune       min_fnorm = jac->fnorms[i];
2465e806d2eSPeter Brune       min_i     = i;
2479c2f3473SPeter Brune     }
2489c2f3473SPeter Brune   }
2495e806d2eSPeter Brune 
2505e806d2eSPeter Brune   /* stagnation or divergence restart to the solution of the solver that failed the least */
2515e806d2eSPeter Brune   if (PetscRealPart(total) < jac->stol || min_fnorm*jac->rtol < *fnorm) {
252*58bdfa14SPeter Brune     ierr = VecCopy(jac->Xes[min_i],X);CHKERRQ(ierr);
253*58bdfa14SPeter Brune     ierr = VecCopy(jac->Fes[min_i],F);CHKERRQ(ierr);
2545e806d2eSPeter Brune     *fnorm = min_fnorm;
2555e806d2eSPeter Brune   }
25690a8ba9bSPeter Brune   PetscFunctionReturn(0);
25790a8ba9bSPeter Brune }
25890a8ba9bSPeter Brune 
259eed5f15bSPeter Brune #undef __FUNCT__
260eed5f15bSPeter Brune #define __FUNCT__ "SNESSetUp_Composite"
261eed5f15bSPeter Brune static PetscErrorCode SNESSetUp_Composite(SNES snes)
262eed5f15bSPeter Brune {
263eed5f15bSPeter Brune   PetscErrorCode     ierr;
264dd515a93SPeter Brune   DM                 dm;
265eed5f15bSPeter Brune   SNES_Composite     *jac = (SNES_Composite*)snes->data;
266eed5f15bSPeter Brune   SNES_CompositeLink next = jac->head;
26790a8ba9bSPeter Brune   PetscInt           n=0,i;
26890a8ba9bSPeter Brune   Vec                F;
269eed5f15bSPeter Brune 
270eed5f15bSPeter Brune   PetscFunctionBegin;
271dd515a93SPeter Brune   ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr);
272eed5f15bSPeter Brune   while (next) {
27390a8ba9bSPeter Brune     n++;
274dd515a93SPeter Brune     ierr = SNESSetDM(next->snes,dm);CHKERRQ(ierr);
275eed5f15bSPeter Brune     ierr = SNESSetFromOptions(next->snes);CHKERRQ(ierr);
276eed5f15bSPeter Brune     next = next->next;
277eed5f15bSPeter Brune   }
27890a8ba9bSPeter Brune   jac->nsnes = n;
27990a8ba9bSPeter Brune   ierr = SNESGetFunction(snes,&F,NULL,NULL);CHKERRQ(ierr);
28090a8ba9bSPeter Brune   if (jac->type == SNES_COMPOSITE_ADDITIVEOPTIMAL) {
28190a8ba9bSPeter Brune     ierr = VecDuplicateVecs(F,jac->nsnes,&jac->Xes);CHKERRQ(ierr);
28290a8ba9bSPeter Brune     ierr = PetscMalloc(sizeof(Vec)*n,&jac->Fes);CHKERRQ(ierr);
2839c2f3473SPeter Brune     ierr = PetscMalloc(sizeof(PetscReal)*n,&jac->fnorms);CHKERRQ(ierr);
28490a8ba9bSPeter Brune     next = jac->head;
28590a8ba9bSPeter Brune     i = 0;
28690a8ba9bSPeter Brune     while (next) {
28790a8ba9bSPeter Brune       ierr = SNESGetFunction(next->snes,&F,NULL,NULL);CHKERRQ(ierr);
28890a8ba9bSPeter Brune       jac->Fes[i] = F;
28990a8ba9bSPeter Brune       ierr = PetscObjectReference((PetscObject)F);CHKERRQ(ierr);
29090a8ba9bSPeter Brune       next = next->next;
29190a8ba9bSPeter Brune       i++;
29290a8ba9bSPeter Brune     }
29390a8ba9bSPeter Brune     /* allocate the subspace direct solve area */
29490a8ba9bSPeter Brune     jac->nrhs  = 1;
2959c2f3473SPeter Brune     jac->lda   = jac->nsnes;
29690a8ba9bSPeter Brune     jac->ldb   = jac->nsnes;
29790a8ba9bSPeter Brune     jac->n     = jac->nsnes;
29890a8ba9bSPeter Brune 
2999c2f3473SPeter Brune     ierr = PetscMalloc(jac->n*jac->n*sizeof(PetscScalar),&jac->h);CHKERRQ(ierr);
3009c2f3473SPeter Brune     ierr = PetscMalloc(jac->n*sizeof(PetscScalar),&jac->beta);CHKERRQ(ierr);
3019c2f3473SPeter Brune     ierr = PetscMalloc(jac->n*sizeof(PetscScalar),&jac->s);CHKERRQ(ierr);
3029c2f3473SPeter Brune     ierr = PetscMalloc(jac->n*sizeof(PetscScalar),&jac->g);CHKERRQ(ierr);
3039c2f3473SPeter Brune     jac->lwork = 12*jac->n;
30490a8ba9bSPeter Brune #if PETSC_USE_COMPLEX
30590a8ba9bSPeter Brune     ierr = PetscMalloc(sizeof(PetscReal)*jac->lwork,&jac->rwork);CHKERRQ(ierr);
30690a8ba9bSPeter Brune #endif
30790a8ba9bSPeter Brune     ierr = PetscMalloc(sizeof(PetscScalar)*jac->lwork,&jac->work);CHKERRQ(ierr);
30890a8ba9bSPeter Brune   }
30990a8ba9bSPeter Brune 
310eed5f15bSPeter Brune   PetscFunctionReturn(0);
311eed5f15bSPeter Brune }
312eed5f15bSPeter Brune 
313eed5f15bSPeter Brune #undef __FUNCT__
314eed5f15bSPeter Brune #define __FUNCT__ "SNESReset_Composite"
315eed5f15bSPeter Brune static PetscErrorCode SNESReset_Composite(SNES snes)
316eed5f15bSPeter Brune {
317eed5f15bSPeter Brune   SNES_Composite     *jac = (SNES_Composite*)snes->data;
318eed5f15bSPeter Brune   PetscErrorCode   ierr;
319eed5f15bSPeter Brune   SNES_CompositeLink next = jac->head;
320eed5f15bSPeter Brune 
321eed5f15bSPeter Brune   PetscFunctionBegin;
322eed5f15bSPeter Brune   while (next) {
323eed5f15bSPeter Brune     ierr = SNESReset(next->snes);CHKERRQ(ierr);
324eed5f15bSPeter Brune     next = next->next;
325eed5f15bSPeter Brune   }
326eed5f15bSPeter Brune   ierr = VecDestroy(&jac->Xorig);CHKERRQ(ierr);
32790a8ba9bSPeter Brune   if (jac->Xes) {ierr = VecDestroyVecs(jac->nsnes,&jac->Xes);CHKERRQ(ierr);}
32890a8ba9bSPeter Brune   if (jac->Fes) {ierr = VecDestroyVecs(jac->nsnes,&jac->Fes);CHKERRQ(ierr);}
3299c2f3473SPeter Brune   ierr = PetscFree(jac->fnorms);CHKERRQ(ierr);
33090a8ba9bSPeter Brune   ierr = PetscFree(jac->h);CHKERRQ(ierr);
33190a8ba9bSPeter Brune   ierr = PetscFree(jac->s);CHKERRQ(ierr);
3329c2f3473SPeter Brune   ierr = PetscFree(jac->g);CHKERRQ(ierr);
33390a8ba9bSPeter Brune   ierr = PetscFree(jac->beta);CHKERRQ(ierr);
33490a8ba9bSPeter Brune   ierr = PetscFree(jac->work);CHKERRQ(ierr);
33590a8ba9bSPeter Brune   ierr = PetscFree(jac->rwork);CHKERRQ(ierr);
33690a8ba9bSPeter Brune 
337eed5f15bSPeter Brune   PetscFunctionReturn(0);
338eed5f15bSPeter Brune }
339eed5f15bSPeter Brune 
340eed5f15bSPeter Brune #undef __FUNCT__
341eed5f15bSPeter Brune #define __FUNCT__ "SNESDestroy_Composite"
342eed5f15bSPeter Brune static PetscErrorCode SNESDestroy_Composite(SNES snes)
343eed5f15bSPeter Brune {
344eed5f15bSPeter Brune   SNES_Composite     *jac = (SNES_Composite*)snes->data;
345eed5f15bSPeter Brune   PetscErrorCode   ierr;
346eed5f15bSPeter Brune   SNES_CompositeLink next = jac->head,next_tmp;
347eed5f15bSPeter Brune 
348eed5f15bSPeter Brune   PetscFunctionBegin;
349eed5f15bSPeter Brune   ierr = SNESReset_Composite(snes);CHKERRQ(ierr);
350eed5f15bSPeter Brune   while (next) {
351eed5f15bSPeter Brune     ierr     = SNESDestroy(&next->snes);CHKERRQ(ierr);
352eed5f15bSPeter Brune     next_tmp = next;
353eed5f15bSPeter Brune     next     = next->next;
354eed5f15bSPeter Brune     ierr     = PetscFree(next_tmp);CHKERRQ(ierr);
355eed5f15bSPeter Brune   }
356eed5f15bSPeter Brune   ierr = PetscFree(snes->data);CHKERRQ(ierr);
357eed5f15bSPeter Brune   PetscFunctionReturn(0);
358eed5f15bSPeter Brune }
359eed5f15bSPeter Brune 
360eed5f15bSPeter Brune #undef __FUNCT__
361eed5f15bSPeter Brune #define __FUNCT__ "SNESSetFromOptions_Composite"
362eed5f15bSPeter Brune static PetscErrorCode SNESSetFromOptions_Composite(SNES snes)
363eed5f15bSPeter Brune {
364eed5f15bSPeter Brune   SNES_Composite     *jac = (SNES_Composite*)snes->data;
365eed5f15bSPeter Brune   PetscErrorCode     ierr;
366eed5f15bSPeter Brune   PetscInt           nmax = 8,i;
367eed5f15bSPeter Brune   SNES_CompositeLink next;
368eed5f15bSPeter Brune   char               *sneses[8];
3698f626970SPeter Brune   PetscReal          dmps[8];
370eed5f15bSPeter Brune   PetscBool          flg;
371eed5f15bSPeter Brune 
372eed5f15bSPeter Brune   PetscFunctionBegin;
373eed5f15bSPeter Brune   ierr = PetscOptionsHead("Composite preconditioner options");CHKERRQ(ierr);
374eed5f15bSPeter Brune   ierr = PetscOptionsEnum("-snes_composite_type","Type of composition","SNESCompositeSetType",SNESCompositeTypes,(PetscEnum)jac->type,(PetscEnum*)&jac->type,&flg);CHKERRQ(ierr);
375eed5f15bSPeter Brune   if (flg) {
376eed5f15bSPeter Brune     ierr = SNESCompositeSetType(snes,jac->type);CHKERRQ(ierr);
377eed5f15bSPeter Brune   }
378eed5f15bSPeter Brune   ierr = PetscOptionsStringArray("-snes_composite_sneses","List of composite solvers","SNESCompositeAddSNES",sneses,&nmax,&flg);CHKERRQ(ierr);
379eed5f15bSPeter Brune   if (flg) {
380eed5f15bSPeter Brune     for (i=0; i<nmax; i++) {
381eed5f15bSPeter Brune       ierr = SNESCompositeAddSNES(snes,sneses[i]);CHKERRQ(ierr);
382eed5f15bSPeter Brune       ierr = PetscFree(sneses[i]);CHKERRQ(ierr);   /* deallocate string sneses[i], which is allocated in PetscOptionsStringArray() */
383eed5f15bSPeter Brune     }
384eed5f15bSPeter Brune   }
3858f626970SPeter Brune   ierr = PetscOptionsRealArray("-snes_composite_damping","Damping of the additive composite solvers","SNESCompositeSetDamping",dmps,&nmax,&flg);CHKERRQ(ierr);
3868f626970SPeter Brune   if (flg) {
3878f626970SPeter Brune     for (i=0; i<nmax; i++) {
3888f626970SPeter Brune       ierr = SNESCompositeSetDamping(snes,i,dmps[i]);CHKERRQ(ierr);
3898f626970SPeter Brune     }
3908f626970SPeter Brune   }
3915e806d2eSPeter Brune   ierr = PetscOptionsReal("-snes_composite_stol","Step tolerance for restart on the additive composite solvers","",jac->stol,&jac->stol,NULL);CHKERRQ(ierr);
3925e806d2eSPeter Brune   ierr = PetscOptionsReal("-snes_composite_rtol","Residual tolerance for the additive composite solvers","",jac->rtol,&jac->rtol,NULL);CHKERRQ(ierr);
393eed5f15bSPeter Brune   ierr = PetscOptionsTail();CHKERRQ(ierr);
394eed5f15bSPeter Brune 
395eed5f15bSPeter Brune   next = jac->head;
396eed5f15bSPeter Brune   while (next) {
397eed5f15bSPeter Brune     ierr = SNESSetFromOptions(next->snes);CHKERRQ(ierr);
398eed5f15bSPeter Brune     next = next->next;
399eed5f15bSPeter Brune   }
400eed5f15bSPeter Brune   PetscFunctionReturn(0);
401eed5f15bSPeter Brune }
402eed5f15bSPeter Brune 
403eed5f15bSPeter Brune #undef __FUNCT__
404eed5f15bSPeter Brune #define __FUNCT__ "SNESView_Composite"
405eed5f15bSPeter Brune static PetscErrorCode SNESView_Composite(SNES snes,PetscViewer viewer)
406eed5f15bSPeter Brune {
407eed5f15bSPeter Brune   SNES_Composite     *jac = (SNES_Composite*)snes->data;
408eed5f15bSPeter Brune   PetscErrorCode   ierr;
409eed5f15bSPeter Brune   SNES_CompositeLink next = jac->head;
410eed5f15bSPeter Brune   PetscBool        iascii;
411eed5f15bSPeter Brune 
412eed5f15bSPeter Brune   PetscFunctionBegin;
413eed5f15bSPeter Brune   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
414eed5f15bSPeter Brune   if (iascii) {
415eed5f15bSPeter Brune     ierr = PetscViewerASCIIPrintf(viewer,"Composite SNES type - %s\n",SNESCompositeTypes[jac->type]);CHKERRQ(ierr);
416eed5f15bSPeter Brune     ierr = PetscViewerASCIIPrintf(viewer,"SNESes on composite preconditioner follow\n");CHKERRQ(ierr);
417eed5f15bSPeter Brune     ierr = PetscViewerASCIIPrintf(viewer,"---------------------------------\n");CHKERRQ(ierr);
418eed5f15bSPeter Brune   }
419eed5f15bSPeter Brune   if (iascii) {
420eed5f15bSPeter Brune     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
421eed5f15bSPeter Brune   }
422eed5f15bSPeter Brune   while (next) {
423eed5f15bSPeter Brune     ierr = SNESView(next->snes,viewer);CHKERRQ(ierr);
424eed5f15bSPeter Brune     next = next->next;
425eed5f15bSPeter Brune   }
426eed5f15bSPeter Brune   if (iascii) {
427eed5f15bSPeter Brune     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
428eed5f15bSPeter Brune     ierr = PetscViewerASCIIPrintf(viewer,"---------------------------------\n");CHKERRQ(ierr);
429eed5f15bSPeter Brune   }
430eed5f15bSPeter Brune   PetscFunctionReturn(0);
431eed5f15bSPeter Brune }
432eed5f15bSPeter Brune 
433eed5f15bSPeter Brune /* ------------------------------------------------------------------------------*/
434eed5f15bSPeter Brune 
435eed5f15bSPeter Brune #undef __FUNCT__
436eed5f15bSPeter Brune #define __FUNCT__ "SNESCompositeSetType_Composite"
437eed5f15bSPeter Brune static PetscErrorCode  SNESCompositeSetType_Composite(SNES snes,SNESCompositeType type)
438eed5f15bSPeter Brune {
439eed5f15bSPeter Brune   SNES_Composite *jac = (SNES_Composite*)snes->data;
440eed5f15bSPeter Brune 
441eed5f15bSPeter Brune   PetscFunctionBegin;
442eed5f15bSPeter Brune   jac->type = type;
443eed5f15bSPeter Brune   PetscFunctionReturn(0);
444eed5f15bSPeter Brune }
445eed5f15bSPeter Brune 
446eed5f15bSPeter Brune #undef __FUNCT__
447eed5f15bSPeter Brune #define __FUNCT__ "SNESCompositeAddSNES_Composite"
448eed5f15bSPeter Brune static PetscErrorCode  SNESCompositeAddSNES_Composite(SNES snes,SNESType type)
449eed5f15bSPeter Brune {
450eed5f15bSPeter Brune   SNES_Composite     *jac;
451eed5f15bSPeter Brune   SNES_CompositeLink next,ilink;
452eed5f15bSPeter Brune   PetscErrorCode   ierr;
453eed5f15bSPeter Brune   PetscInt         cnt = 0;
454eed5f15bSPeter Brune   const char       *prefix;
455eed5f15bSPeter Brune   char             newprefix[8];
456eed5f15bSPeter Brune   DM               dm;
457eed5f15bSPeter Brune 
458eed5f15bSPeter Brune   PetscFunctionBegin;
459eed5f15bSPeter Brune   ierr        = PetscNewLog(snes,struct _SNES_CompositeLink,&ilink);CHKERRQ(ierr);
460eed5f15bSPeter Brune   ilink->next = 0;
461eed5f15bSPeter Brune   ierr        = SNESCreate(PetscObjectComm((PetscObject)snes),&ilink->snes);CHKERRQ(ierr);
462eed5f15bSPeter Brune   ierr        = PetscLogObjectParent((PetscObject)snes,(PetscObject)ilink->snes);CHKERRQ(ierr);
463eed5f15bSPeter Brune   ierr        = SNESGetDM(snes,&dm);CHKERRQ(ierr);
464eed5f15bSPeter Brune   ierr        = SNESSetDM(ilink->snes,dm);CHKERRQ(ierr);
465cf5b3eb5SPeter Brune   ierr        = SNESSetTolerances(ilink->snes,ilink->snes->abstol,ilink->snes->rtol,ilink->snes->stol,1,ilink->snes->max_funcs);CHKERRQ(ierr);
466eed5f15bSPeter Brune   jac  = (SNES_Composite*)snes->data;
467eed5f15bSPeter Brune   next = jac->head;
468eed5f15bSPeter Brune   if (!next) {
469eed5f15bSPeter Brune     jac->head       = ilink;
470eed5f15bSPeter Brune     ilink->previous = NULL;
471eed5f15bSPeter Brune   } else {
472eed5f15bSPeter Brune     cnt++;
473eed5f15bSPeter Brune     while (next->next) {
474eed5f15bSPeter Brune       next = next->next;
475eed5f15bSPeter Brune       cnt++;
476eed5f15bSPeter Brune     }
477eed5f15bSPeter Brune     next->next      = ilink;
478eed5f15bSPeter Brune     ilink->previous = next;
479eed5f15bSPeter Brune   }
480eed5f15bSPeter Brune   ierr = SNESGetOptionsPrefix(snes,&prefix);CHKERRQ(ierr);
481eed5f15bSPeter Brune   ierr = SNESSetOptionsPrefix(ilink->snes,prefix);CHKERRQ(ierr);
482eed5f15bSPeter Brune   sprintf(newprefix,"sub_%d_",(int)cnt);
483eed5f15bSPeter Brune   ierr = SNESAppendOptionsPrefix(ilink->snes,newprefix);CHKERRQ(ierr);
4848f626970SPeter Brune   ierr = PetscObjectIncrementTabLevel((PetscObject)ilink->snes,(PetscObject)snes,1);CHKERRQ(ierr);
485eed5f15bSPeter Brune   ierr = SNESSetType(ilink->snes,type);CHKERRQ(ierr);
48672edecb9SPeter Brune   ierr = SNESSetNormSchedule(ilink->snes, SNES_NORM_FINAL_ONLY);CHKERRQ(ierr);
4878f626970SPeter Brune   ilink->dmp = 1.0;
48890a8ba9bSPeter Brune   jac->nsnes++;
489eed5f15bSPeter Brune   PetscFunctionReturn(0);
490eed5f15bSPeter Brune }
491eed5f15bSPeter Brune 
492eed5f15bSPeter Brune #undef __FUNCT__
493eed5f15bSPeter Brune #define __FUNCT__ "SNESCompositeGetSNES_Composite"
494eed5f15bSPeter Brune static PetscErrorCode  SNESCompositeGetSNES_Composite(SNES snes,PetscInt n,SNES *subsnes)
495eed5f15bSPeter Brune {
496eed5f15bSPeter Brune   SNES_Composite     *jac;
497eed5f15bSPeter Brune   SNES_CompositeLink next;
498eed5f15bSPeter Brune   PetscInt         i;
499eed5f15bSPeter Brune 
500eed5f15bSPeter Brune   PetscFunctionBegin;
501eed5f15bSPeter Brune   jac  = (SNES_Composite*)snes->data;
502eed5f15bSPeter Brune   next = jac->head;
503eed5f15bSPeter Brune   for (i=0; i<n; i++) {
504eed5f15bSPeter Brune     if (!next->next) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_INCOMP,"Not enough SNESes in composite preconditioner");
505eed5f15bSPeter Brune     next = next->next;
506eed5f15bSPeter Brune   }
507eed5f15bSPeter Brune   *subsnes = next->snes;
508eed5f15bSPeter Brune   PetscFunctionReturn(0);
509eed5f15bSPeter Brune }
510eed5f15bSPeter Brune 
511eed5f15bSPeter Brune /* -------------------------------------------------------------------------------- */
512eed5f15bSPeter Brune #undef __FUNCT__
513eed5f15bSPeter Brune #define __FUNCT__ "SNESCompositeSetType"
514e31ab4f9SPeter Brune /*@C
515eed5f15bSPeter Brune    SNESCompositeSetType - Sets the type of composite preconditioner.
516eed5f15bSPeter Brune 
517eed5f15bSPeter Brune    Logically Collective on SNES
518eed5f15bSPeter Brune 
519eed5f15bSPeter Brune    Input Parameter:
520eed5f15bSPeter Brune +  snes - the preconditioner context
521eed5f15bSPeter Brune -  type - SNES_COMPOSITE_ADDITIVE (default), SNES_COMPOSITE_MULTIPLICATIVE
522eed5f15bSPeter Brune 
523eed5f15bSPeter Brune    Options Database Key:
524eed5f15bSPeter Brune .  -snes_composite_type <type: one of multiplicative, additive, special> - Sets composite preconditioner type
525eed5f15bSPeter Brune 
526eed5f15bSPeter Brune    Level: Developer
527eed5f15bSPeter Brune 
528eed5f15bSPeter Brune .keywords: SNES, set, type, composite preconditioner, additive, multiplicative
529eed5f15bSPeter Brune @*/
530eed5f15bSPeter Brune PetscErrorCode  SNESCompositeSetType(SNES snes,SNESCompositeType type)
531eed5f15bSPeter Brune {
532eed5f15bSPeter Brune   PetscErrorCode ierr;
533eed5f15bSPeter Brune 
534eed5f15bSPeter Brune   PetscFunctionBegin;
535eed5f15bSPeter Brune   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
536eed5f15bSPeter Brune   PetscValidLogicalCollectiveEnum(snes,type,2);
537eed5f15bSPeter Brune   ierr = PetscTryMethod(snes,"SNESCompositeSetType_C",(SNES,SNESCompositeType),(snes,type));CHKERRQ(ierr);
538eed5f15bSPeter Brune   PetscFunctionReturn(0);
539eed5f15bSPeter Brune }
540eed5f15bSPeter Brune 
541eed5f15bSPeter Brune #undef __FUNCT__
542eed5f15bSPeter Brune #define __FUNCT__ "SNESCompositeAddSNES"
543eed5f15bSPeter Brune /*@C
544eed5f15bSPeter Brune    SNESCompositeAddSNES - Adds another SNES to the composite SNES.
545eed5f15bSPeter Brune 
546eed5f15bSPeter Brune    Collective on SNES
547eed5f15bSPeter Brune 
548eed5f15bSPeter Brune    Input Parameters:
549eed5f15bSPeter Brune +  snes - the preconditioner context
550eed5f15bSPeter Brune -  type - the type of the new preconditioner
551eed5f15bSPeter Brune 
552eed5f15bSPeter Brune    Level: Developer
553eed5f15bSPeter Brune 
554eed5f15bSPeter Brune .keywords: SNES, composite preconditioner, add
555eed5f15bSPeter Brune @*/
556eed5f15bSPeter Brune PetscErrorCode  SNESCompositeAddSNES(SNES snes,SNESType type)
557eed5f15bSPeter Brune {
558eed5f15bSPeter Brune   PetscErrorCode ierr;
559eed5f15bSPeter Brune 
560eed5f15bSPeter Brune   PetscFunctionBegin;
561eed5f15bSPeter Brune   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
562eed5f15bSPeter Brune   ierr = PetscTryMethod(snes,"SNESCompositeAddSNES_C",(SNES,SNESType),(snes,type));CHKERRQ(ierr);
563eed5f15bSPeter Brune   PetscFunctionReturn(0);
564eed5f15bSPeter Brune }
565eed5f15bSPeter Brune #undef __FUNCT__
566eed5f15bSPeter Brune #define __FUNCT__ "SNESCompositeGetSNES"
567eed5f15bSPeter Brune /*@
568eed5f15bSPeter Brune    SNESCompositeGetSNES - Gets one of the SNES objects in the composite SNES.
569eed5f15bSPeter Brune 
570eed5f15bSPeter Brune    Not Collective
571eed5f15bSPeter Brune 
572eed5f15bSPeter Brune    Input Parameter:
573eed5f15bSPeter Brune +  snes - the preconditioner context
574eed5f15bSPeter Brune -  n - the number of the snes requested
575eed5f15bSPeter Brune 
576eed5f15bSPeter Brune    Output Parameters:
577eed5f15bSPeter Brune .  subsnes - the SNES requested
578eed5f15bSPeter Brune 
579eed5f15bSPeter Brune    Level: Developer
580eed5f15bSPeter Brune 
581eed5f15bSPeter Brune .keywords: SNES, get, composite preconditioner, sub preconditioner
582eed5f15bSPeter Brune 
583eed5f15bSPeter Brune .seealso: SNESCompositeAddSNES()
584eed5f15bSPeter Brune @*/
585eed5f15bSPeter Brune PetscErrorCode  SNESCompositeGetSNES(SNES snes,PetscInt n,SNES *subsnes)
586eed5f15bSPeter Brune {
587eed5f15bSPeter Brune   PetscErrorCode ierr;
588eed5f15bSPeter Brune 
589eed5f15bSPeter Brune   PetscFunctionBegin;
590eed5f15bSPeter Brune   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
591eed5f15bSPeter Brune   PetscValidPointer(subsnes,3);
592eed5f15bSPeter Brune   ierr = PetscUseMethod(snes,"SNESCompositeGetSNES_C",(SNES,PetscInt,SNES*),(snes,n,subsnes));CHKERRQ(ierr);
593eed5f15bSPeter Brune   PetscFunctionReturn(0);
594eed5f15bSPeter Brune }
595eed5f15bSPeter Brune 
596eed5f15bSPeter Brune #undef __FUNCT__
5978f626970SPeter Brune #define __FUNCT__ "SNESCompositeSetDamping_Composite"
5988f626970SPeter Brune static PetscErrorCode  SNESCompositeSetDamping_Composite(SNES snes,PetscInt n,PetscReal dmp)
5998f626970SPeter Brune {
6008f626970SPeter Brune   SNES_Composite     *jac;
6018f626970SPeter Brune   SNES_CompositeLink next;
6028f626970SPeter Brune   PetscInt         i;
6038f626970SPeter Brune 
6048f626970SPeter Brune   PetscFunctionBegin;
6058f626970SPeter Brune   jac  = (SNES_Composite*)snes->data;
6068f626970SPeter Brune   next = jac->head;
6078f626970SPeter Brune   for (i=0; i<n; i++) {
6088f626970SPeter Brune     if (!next->next) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_INCOMP,"Not enough SNESes in composite preconditioner");
6098f626970SPeter Brune     next = next->next;
6108f626970SPeter Brune   }
6118f626970SPeter Brune   next->dmp = dmp;
6128f626970SPeter Brune   PetscFunctionReturn(0);
6138f626970SPeter Brune }
6148f626970SPeter Brune 
6158f626970SPeter Brune #undef __FUNCT__
6168f626970SPeter Brune #define __FUNCT__ "SNESCompositeSetDamping"
6178f626970SPeter Brune /*@
6188f626970SPeter Brune    SNESCompositeSetDamping - Sets the damping of a subsolver when using additive composite SNES.
6198f626970SPeter Brune 
6208f626970SPeter Brune    Not Collective
6218f626970SPeter Brune 
6228f626970SPeter Brune    Input Parameter:
6238f626970SPeter Brune +  snes - the preconditioner context
6248f626970SPeter Brune .  n - the number of the snes requested
6258f626970SPeter Brune -  dmp - the damping
6268f626970SPeter Brune 
6278f626970SPeter Brune    Level: Developer
6288f626970SPeter Brune 
6298f626970SPeter Brune .keywords: SNES, get, composite preconditioner, sub preconditioner
6308f626970SPeter Brune 
6318f626970SPeter Brune .seealso: SNESCompositeAddSNES()
6328f626970SPeter Brune @*/
6338f626970SPeter Brune PetscErrorCode  SNESCompositeSetDamping(SNES snes,PetscInt n,PetscReal dmp)
6348f626970SPeter Brune {
6358f626970SPeter Brune   PetscErrorCode ierr;
6368f626970SPeter Brune 
6378f626970SPeter Brune   PetscFunctionBegin;
6388f626970SPeter Brune   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
6398f626970SPeter Brune   ierr = PetscUseMethod(snes,"SNESCompositeSetDamping_C",(SNES,PetscInt,PetscReal),(snes,n,dmp));CHKERRQ(ierr);
6408f626970SPeter Brune   PetscFunctionReturn(0);
6418f626970SPeter Brune }
6428f626970SPeter Brune 
6438f626970SPeter Brune #undef __FUNCT__
644eed5f15bSPeter Brune #define __FUNCT__ "SNESSolve_Composite"
645eed5f15bSPeter Brune PetscErrorCode SNESSolve_Composite(SNES snes)
646eed5f15bSPeter Brune {
647eed5f15bSPeter Brune   Vec            F;
648eed5f15bSPeter Brune   Vec            X;
649eed5f15bSPeter Brune   Vec            B;
650eed5f15bSPeter Brune   PetscInt       i;
651eed5f15bSPeter Brune   PetscReal      fnorm = 0.0;
652eed5f15bSPeter Brune   PetscErrorCode ierr;
65372edecb9SPeter Brune   SNESNormSchedule normtype;
654eed5f15bSPeter Brune   SNES_Composite *comp = (SNES_Composite*)snes->data;
655eed5f15bSPeter Brune 
656eed5f15bSPeter Brune   PetscFunctionBegin;
657eed5f15bSPeter Brune   X = snes->vec_sol;
658eed5f15bSPeter Brune   F = snes->vec_func;
659eed5f15bSPeter Brune   B = snes->vec_rhs;
660eed5f15bSPeter Brune 
661eed5f15bSPeter Brune   ierr         = PetscObjectAMSTakeAccess((PetscObject)snes);CHKERRQ(ierr);
662eed5f15bSPeter Brune   snes->iter   = 0;
663eed5f15bSPeter Brune   snes->norm   = 0.;
664eed5f15bSPeter Brune   ierr         = PetscObjectAMSGrantAccess((PetscObject)snes);CHKERRQ(ierr);
665eed5f15bSPeter Brune   snes->reason = SNES_CONVERGED_ITERATING;
66672edecb9SPeter Brune   ierr         = SNESGetNormSchedule(snes, &normtype);CHKERRQ(ierr);
66772edecb9SPeter Brune   if (normtype == SNES_NORM_ALWAYS || normtype == SNES_NORM_INITIAL_ONLY || normtype == SNES_NORM_INITIAL_FINAL_ONLY) {
668eed5f15bSPeter Brune     if (!snes->vec_func_init_set) {
669eed5f15bSPeter Brune       ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
670eed5f15bSPeter Brune       if (snes->domainerror) {
671eed5f15bSPeter Brune         snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
672eed5f15bSPeter Brune         PetscFunctionReturn(0);
673eed5f15bSPeter Brune       }
674eed5f15bSPeter Brune     } else snes->vec_func_init_set = PETSC_FALSE;
675eed5f15bSPeter Brune 
676eed5f15bSPeter Brune     /* convergence test */
677eed5f15bSPeter Brune     if (!snes->norm_init_set) {
678eed5f15bSPeter Brune       ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F||  */
679eed5f15bSPeter Brune       if (PetscIsInfOrNanReal(fnorm)) {
680eed5f15bSPeter Brune         snes->reason = SNES_DIVERGED_FNORM_NAN;
681eed5f15bSPeter Brune         PetscFunctionReturn(0);
682eed5f15bSPeter Brune       }
683eed5f15bSPeter Brune     } else {
684eed5f15bSPeter Brune       fnorm               = snes->norm_init;
685eed5f15bSPeter Brune       snes->norm_init_set = PETSC_FALSE;
686eed5f15bSPeter Brune     }
687eed5f15bSPeter Brune     ierr       = PetscObjectAMSTakeAccess((PetscObject)snes);CHKERRQ(ierr);
688eed5f15bSPeter Brune     snes->iter = 0;
689eed5f15bSPeter Brune     snes->norm = fnorm;
690eed5f15bSPeter Brune     ierr       = PetscObjectAMSGrantAccess((PetscObject)snes);CHKERRQ(ierr);
691eed5f15bSPeter Brune     ierr       = SNESLogConvergenceHistory(snes,snes->norm,0);CHKERRQ(ierr);
692eed5f15bSPeter Brune     ierr       = SNESMonitor(snes,0,snes->norm);CHKERRQ(ierr);
693eed5f15bSPeter Brune 
694eed5f15bSPeter Brune     /* set parameter for default relative tolerance convergence test */
695eed5f15bSPeter Brune     snes->ttol = fnorm*snes->rtol;
696eed5f15bSPeter Brune 
697eed5f15bSPeter Brune     /* test convergence */
698eed5f15bSPeter Brune     ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
699eed5f15bSPeter Brune     if (snes->reason) PetscFunctionReturn(0);
700eed5f15bSPeter Brune   } else {
701eed5f15bSPeter Brune     ierr = PetscObjectAMSGrantAccess((PetscObject)snes);CHKERRQ(ierr);
702eed5f15bSPeter Brune     ierr = SNESLogConvergenceHistory(snes,snes->norm,0);CHKERRQ(ierr);
703eed5f15bSPeter Brune     ierr = SNESMonitor(snes,0,snes->norm);CHKERRQ(ierr);
704eed5f15bSPeter Brune   }
705eed5f15bSPeter Brune 
706eed5f15bSPeter Brune   /* Call general purpose update function */
707eed5f15bSPeter Brune   if (snes->ops->update) {
708eed5f15bSPeter Brune     ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
709eed5f15bSPeter Brune   }
710eed5f15bSPeter Brune 
711eed5f15bSPeter Brune   for (i = 0; i < snes->max_its; i++) {
712eed5f15bSPeter Brune     if (comp->type == SNES_COMPOSITE_ADDITIVE) {
713eed5f15bSPeter Brune       ierr = SNESCompositeApply_Additive(snes,X,B,F,&fnorm);CHKERRQ(ierr);
714eed5f15bSPeter Brune     } else if (comp->type == SNES_COMPOSITE_MULTIPLICATIVE) {
715eed5f15bSPeter Brune       ierr = SNESCompositeApply_Multiplicative(snes,X,B,F,&fnorm);CHKERRQ(ierr);
71690a8ba9bSPeter Brune     } else if (comp->type == SNES_COMPOSITE_ADDITIVEOPTIMAL) {
71790a8ba9bSPeter Brune       ierr = SNESCompositeApply_AdditiveOptimal(snes,X,B,F,&fnorm);CHKERRQ(ierr);
718eed5f15bSPeter Brune     } else {
71990a8ba9bSPeter Brune       SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"Unsupported SNESComposite type");
720eed5f15bSPeter Brune     }
721eed5f15bSPeter Brune     if ((i == snes->max_its - 1) && (normtype == SNES_NORM_INITIAL_FINAL_ONLY || normtype == SNES_NORM_FINAL_ONLY)) {
722eed5f15bSPeter Brune       ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
723eed5f15bSPeter Brune       if (snes->domainerror) {
724eed5f15bSPeter Brune         snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
725eed5f15bSPeter Brune         break;
726eed5f15bSPeter Brune       }
727eed5f15bSPeter Brune       ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F||  */
728eed5f15bSPeter Brune       if (PetscIsInfOrNanReal(fnorm)) {
729eed5f15bSPeter Brune         snes->reason = SNES_DIVERGED_FNORM_NAN;
730eed5f15bSPeter Brune         break;
731eed5f15bSPeter Brune       }
732eed5f15bSPeter Brune     }
733eed5f15bSPeter Brune     /* Monitor convergence */
734eed5f15bSPeter Brune     ierr       = PetscObjectAMSTakeAccess((PetscObject)snes);CHKERRQ(ierr);
735eed5f15bSPeter Brune     snes->iter = i+1;
736eed5f15bSPeter Brune     snes->norm = fnorm;
737eed5f15bSPeter Brune     ierr       = PetscObjectAMSGrantAccess((PetscObject)snes);CHKERRQ(ierr);
738eed5f15bSPeter Brune     ierr       = SNESLogConvergenceHistory(snes,snes->norm,0);CHKERRQ(ierr);
739eed5f15bSPeter Brune     ierr       = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr);
740eed5f15bSPeter Brune     /* Test for convergence */
74172edecb9SPeter Brune     if (normtype == SNES_NORM_ALWAYS) {ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);}
742eed5f15bSPeter Brune     if (snes->reason) break;
743eed5f15bSPeter Brune     /* Call general purpose update function */
744eed5f15bSPeter Brune     if (snes->ops->update) {ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);}
745eed5f15bSPeter Brune   }
74672edecb9SPeter Brune   if (normtype == SNES_NORM_ALWAYS) {
747eed5f15bSPeter Brune     if (i == snes->max_its) {
748eed5f15bSPeter Brune       ierr = PetscInfo1(snes,"Maximum number of iterations has been reached: %D\n",snes->max_its);CHKERRQ(ierr);
749eed5f15bSPeter Brune       if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
750eed5f15bSPeter Brune     }
751eed5f15bSPeter Brune   } else if (!snes->reason) snes->reason = SNES_CONVERGED_ITS;
752eed5f15bSPeter Brune   PetscFunctionReturn(0);
753eed5f15bSPeter Brune }
754eed5f15bSPeter Brune 
755eed5f15bSPeter Brune /* -------------------------------------------------------------------------------------------*/
756eed5f15bSPeter Brune 
757eed5f15bSPeter Brune /*MC
758eed5f15bSPeter Brune      SNESCOMPOSITE - Build a preconditioner by composing together several nonlinear solvers
759eed5f15bSPeter Brune 
760eed5f15bSPeter Brune    Options Database Keys:
761eed5f15bSPeter Brune +  -snes_composite_type <type: one of multiplicative, additive, symmetric_multiplicative, special> - Sets composite preconditioner type
762eed5f15bSPeter Brune -  -snes_composite_sneses - <snes0,snes1,...> list of SNESes to compose
763eed5f15bSPeter Brune 
764eed5f15bSPeter Brune    Level: intermediate
765eed5f15bSPeter Brune 
766eed5f15bSPeter Brune    Concepts: composing solvers
767eed5f15bSPeter Brune 
768eed5f15bSPeter Brune .seealso:  SNESCreate(), SNESSetType(), SNESType (for list of available types), SNES,
769eed5f15bSPeter Brune            SNESSHELL, SNESCompositeSetType(), SNESCompositeSpecialSetAlpha(), SNESCompositeAddSNES(),
770eed5f15bSPeter Brune            SNESCompositeGetSNES()
771eed5f15bSPeter Brune 
772eed5f15bSPeter Brune M*/
773eed5f15bSPeter Brune 
774eed5f15bSPeter Brune #undef __FUNCT__
775eed5f15bSPeter Brune #define __FUNCT__ "SNESCreate_Composite"
776eed5f15bSPeter Brune PETSC_EXTERN PetscErrorCode SNESCreate_Composite(SNES snes)
777eed5f15bSPeter Brune {
778eed5f15bSPeter Brune   PetscErrorCode ierr;
779eed5f15bSPeter Brune   SNES_Composite   *jac;
780eed5f15bSPeter Brune 
781eed5f15bSPeter Brune   PetscFunctionBegin;
782eed5f15bSPeter Brune   ierr = PetscNewLog(snes,SNES_Composite,&jac);CHKERRQ(ierr);
783eed5f15bSPeter Brune 
784eed5f15bSPeter Brune   snes->ops->solve           = SNESSolve_Composite;
785eed5f15bSPeter Brune   snes->ops->setup           = SNESSetUp_Composite;
786eed5f15bSPeter Brune   snes->ops->reset           = SNESReset_Composite;
787eed5f15bSPeter Brune   snes->ops->destroy         = SNESDestroy_Composite;
788eed5f15bSPeter Brune   snes->ops->setfromoptions  = SNESSetFromOptions_Composite;
789eed5f15bSPeter Brune   snes->ops->view            = SNESView_Composite;
790eed5f15bSPeter Brune 
791eed5f15bSPeter Brune   snes->data = (void*)jac;
79290a8ba9bSPeter Brune   jac->type  = SNES_COMPOSITE_ADDITIVEOPTIMAL;
79390a8ba9bSPeter Brune   jac->Fes   = NULL;
79490a8ba9bSPeter Brune   jac->Xes   = NULL;
7959c2f3473SPeter Brune   jac->fnorms = NULL;
79690a8ba9bSPeter Brune   jac->nsnes = 0;
797eed5f15bSPeter Brune   jac->head  = 0;
7985e806d2eSPeter Brune   jac->stol  = 0.1;
7995e806d2eSPeter Brune   jac->rtol  = 1.1;
800eed5f15bSPeter Brune 
80190a8ba9bSPeter Brune   jac->h     = NULL;
80290a8ba9bSPeter Brune   jac->s     = NULL;
80390a8ba9bSPeter Brune   jac->beta  = NULL;
80490a8ba9bSPeter Brune   jac->work  = NULL;
80590a8ba9bSPeter Brune   jac->rwork = NULL;
80690a8ba9bSPeter Brune 
8078f626970SPeter Brune   ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESCompositeSetType_C",SNESCompositeSetType_Composite);CHKERRQ(ierr);
8088f626970SPeter Brune   ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESCompositeAddSNES_C",SNESCompositeAddSNES_Composite);CHKERRQ(ierr);
8098f626970SPeter Brune   ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESCompositeGetSNES_C",SNESCompositeGetSNES_Composite);CHKERRQ(ierr);
8108f626970SPeter Brune   ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESCompositeSetDamping_C",SNESCompositeSetDamping_Composite);CHKERRQ(ierr);
811eed5f15bSPeter Brune   PetscFunctionReturn(0);
812eed5f15bSPeter Brune }
813eed5f15bSPeter Brune 
814