xref: /petsc/src/snes/impls/composite/snescomposite.c (revision b00a91154f763f12aa55f3d53a3f2776f15f49e3)
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;
55d5a53f06SPeter Brune   SNESConvergedReason reason;
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   }
62eed5f15bSPeter Brune   ierr = SNESSolve(next->snes,B,X);CHKERRQ(ierr);
63d5a53f06SPeter Brune   ierr = SNESGetConvergedReason(next->snes,&reason);CHKERRQ(ierr);
64d5a53f06SPeter Brune   if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
65d5a53f06SPeter Brune     snes->reason = SNES_DIVERGED_INNER;
66d5a53f06SPeter Brune     PetscFunctionReturn(0);
67d5a53f06SPeter Brune   }
68eed5f15bSPeter Brune 
69eed5f15bSPeter Brune   while (next->next) {
70eed5f15bSPeter Brune     /* only copy the function over in the case where the functions correspond */
7172edecb9SPeter Brune     if (next->snes->pcside == PC_RIGHT && next->snes->normschedule != SNES_NORM_NONE) {
72eed5f15bSPeter Brune       ierr = SNESGetFunction(next->snes,&FSub,NULL,NULL);CHKERRQ(ierr);
73eed5f15bSPeter Brune       next = next->next;
74eed5f15bSPeter Brune       ierr = SNESSetInitialFunction(next->snes,FSub);CHKERRQ(ierr);
75eed5f15bSPeter Brune     } else {
76eed5f15bSPeter Brune       next = next->next;
77eed5f15bSPeter Brune     }
78eed5f15bSPeter Brune     ierr = SNESSolve(next->snes,B,X);CHKERRQ(ierr);
79d5a53f06SPeter Brune     ierr = SNESGetConvergedReason(next->snes,&reason);CHKERRQ(ierr);
80d5a53f06SPeter Brune     if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
81d5a53f06SPeter Brune       snes->reason = SNES_DIVERGED_INNER;
82d5a53f06SPeter Brune       PetscFunctionReturn(0);
83d5a53f06SPeter Brune     }
84eed5f15bSPeter Brune   }
85eed5f15bSPeter Brune   if (next->snes->pcside == PC_RIGHT) {
86eed5f15bSPeter Brune     ierr = SNESGetFunction(next->snes,&FSub,NULL,NULL);CHKERRQ(ierr);
87eed5f15bSPeter Brune     ierr = VecCopy(FSub,F);CHKERRQ(ierr);
88d5a53f06SPeter Brune     if (fnorm) {
89d5a53f06SPeter Brune       ierr = SNESGetFunctionNorm(next->snes,fnorm);CHKERRQ(ierr);
90d5a53f06SPeter Brune       if (PetscIsInfOrNanReal(*fnorm)) {
91d5a53f06SPeter Brune         snes->reason = SNES_DIVERGED_FNORM_NAN;
92d5a53f06SPeter Brune         PetscFunctionReturn(0);
93d5a53f06SPeter Brune       }
94d5a53f06SPeter Brune     }
9572edecb9SPeter Brune   } else if (snes->normschedule == SNES_NORM_ALWAYS) {
96eed5f15bSPeter Brune     SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
97d5a53f06SPeter Brune     if (snes->domainerror) {
98d5a53f06SPeter Brune       snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
99d5a53f06SPeter Brune       PetscFunctionReturn(0);
100d5a53f06SPeter Brune     }
101d5a53f06SPeter Brune     if (fnorm) {
102d5a53f06SPeter Brune       ierr = VecNorm(F,NORM_2,fnorm);CHKERRQ(ierr);
103d5a53f06SPeter Brune       if (PetscIsInfOrNanReal(*fnorm)) {
104d5a53f06SPeter Brune         snes->reason = SNES_DIVERGED_FNORM_NAN;
105d5a53f06SPeter Brune         PetscFunctionReturn(0);
106d5a53f06SPeter Brune       }
107d5a53f06SPeter Brune     }
108eed5f15bSPeter Brune   }
109eed5f15bSPeter Brune   PetscFunctionReturn(0);
110eed5f15bSPeter Brune }
111eed5f15bSPeter Brune 
112eed5f15bSPeter Brune #undef __FUNCT__
113eed5f15bSPeter Brune #define __FUNCT__ "SNESCompositeApply_Additive"
114eed5f15bSPeter Brune static PetscErrorCode SNESCompositeApply_Additive(SNES snes,Vec X,Vec B,Vec F,PetscReal *fnorm)
115eed5f15bSPeter Brune {
116eed5f15bSPeter Brune   PetscErrorCode      ierr;
117eed5f15bSPeter Brune   SNES_Composite      *jac = (SNES_Composite*)snes->data;
118eed5f15bSPeter Brune   SNES_CompositeLink  next = jac->head;
119eed5f15bSPeter Brune   Vec                 Y,Xorig;
120d5a53f06SPeter Brune   SNESConvergedReason reason;
121eed5f15bSPeter Brune 
122eed5f15bSPeter Brune   PetscFunctionBegin;
123eed5f15bSPeter Brune   Y = snes->vec_sol_update;
124eed5f15bSPeter Brune   if (!jac->Xorig) {ierr = VecDuplicate(X,&jac->Xorig);CHKERRQ(ierr);}
125eed5f15bSPeter Brune   Xorig = jac->Xorig;
126eed5f15bSPeter Brune   ierr = VecCopy(X,Xorig);
127eed5f15bSPeter Brune   if (!next) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"No composite SNESes supplied via SNESCompositeAddSNES() or -snes_composite_sneses");
12872edecb9SPeter Brune   if (snes->normschedule == SNES_NORM_ALWAYS) {
129eed5f15bSPeter Brune     ierr = SNESSetInitialFunction(next->snes,F);CHKERRQ(ierr);
130eed5f15bSPeter Brune     while (next->next) {
131eed5f15bSPeter Brune       next = next->next;
132eed5f15bSPeter Brune       ierr = SNESSetInitialFunction(next->snes,F);CHKERRQ(ierr);
133eed5f15bSPeter Brune     }
134eed5f15bSPeter Brune   }
135eed5f15bSPeter Brune   next = jac->head;
1368f626970SPeter Brune   ierr = VecCopy(Xorig,Y);CHKERRQ(ierr);
137eed5f15bSPeter Brune   ierr = SNESSolve(next->snes,B,Y);CHKERRQ(ierr);
138d5a53f06SPeter Brune   ierr = SNESGetConvergedReason(next->snes,&reason);CHKERRQ(ierr);
139d5a53f06SPeter Brune   if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
140d5a53f06SPeter Brune     snes->reason = SNES_DIVERGED_INNER;
141d5a53f06SPeter Brune     PetscFunctionReturn(0);
142d5a53f06SPeter Brune   }
143eed5f15bSPeter Brune   ierr = VecAXPY(Y,-1.0,Xorig);CHKERRQ(ierr);
1448f626970SPeter Brune   ierr = VecAXPY(X,next->dmp,Y);CHKERRQ(ierr);
1458f626970SPeter Brune   while (next->next) {
1468f626970SPeter Brune     next = next->next;
1478f626970SPeter Brune     ierr = VecCopy(Xorig,Y);CHKERRQ(ierr);
1488f626970SPeter Brune     ierr = SNESSolve(next->snes,B,Y);CHKERRQ(ierr);
149d5a53f06SPeter Brune     ierr = SNESGetConvergedReason(next->snes,&reason);CHKERRQ(ierr);
150d5a53f06SPeter Brune     if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
151d5a53f06SPeter Brune       snes->reason = SNES_DIVERGED_INNER;
152d5a53f06SPeter Brune       PetscFunctionReturn(0);
153d5a53f06SPeter Brune     }
1548f626970SPeter Brune     ierr = VecAXPY(Y,-1.0,Xorig);CHKERRQ(ierr);
1558f626970SPeter Brune     ierr = VecAXPY(X,next->dmp,Y);CHKERRQ(ierr);
156eed5f15bSPeter Brune   }
15772edecb9SPeter Brune   if (snes->normschedule == SNES_NORM_ALWAYS) {
158eed5f15bSPeter Brune     ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
159eed5f15bSPeter Brune     if (fnorm) {ierr = VecNorm(F,NORM_2,fnorm);CHKERRQ(ierr);}
160eed5f15bSPeter Brune   }
161eed5f15bSPeter Brune   PetscFunctionReturn(0);
162eed5f15bSPeter Brune }
16390a8ba9bSPeter Brune 
16490a8ba9bSPeter Brune #undef __FUNCT__
16590a8ba9bSPeter Brune #define __FUNCT__ "SNESCompositeApply_AdditiveOptimal"
16690a8ba9bSPeter Brune /* approximately solve the overdetermined system:
16790a8ba9bSPeter Brune 
16890a8ba9bSPeter Brune  2*F(x_i)\cdot F(\x_j)\alpha_i = 0
16990a8ba9bSPeter Brune  \alpha_i                      = 1
17090a8ba9bSPeter Brune 
17190a8ba9bSPeter Brune  Which minimizes the L2 norm of the linearization of:
17290a8ba9bSPeter Brune  ||F(\sum_i \alpha_i*x_i)||^2
17390a8ba9bSPeter Brune 
17490a8ba9bSPeter Brune  With the constraint that \sum_i\alpha_i = 1
17590a8ba9bSPeter Brune  Where x_i is the solution from the ith subsolver.
17690a8ba9bSPeter Brune  */
17790a8ba9bSPeter Brune static PetscErrorCode SNESCompositeApply_AdditiveOptimal(SNES snes,Vec X,Vec B,Vec F,PetscReal *fnorm)
17890a8ba9bSPeter Brune {
17990a8ba9bSPeter Brune   PetscErrorCode      ierr;
18090a8ba9bSPeter Brune   SNES_Composite      *jac = (SNES_Composite*)snes->data;
18190a8ba9bSPeter Brune   SNES_CompositeLink  next = jac->head;
18290a8ba9bSPeter Brune   Vec                 *Xes = jac->Xes,*Fes = jac->Fes;
18390a8ba9bSPeter Brune   PetscInt            i,j;
1845e806d2eSPeter Brune   PetscScalar         tot,total,ftf;
1855e806d2eSPeter Brune   PetscReal           min_fnorm;
1865e806d2eSPeter Brune   PetscInt            min_i;
187d5a53f06SPeter Brune   SNESConvergedReason reason;
18890a8ba9bSPeter Brune 
18990a8ba9bSPeter Brune   PetscFunctionBegin;
19090a8ba9bSPeter Brune   if (!next) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"No composite SNESes supplied via SNESCompositeAddSNES() or -snes_composite_sneses");
19158bdfa14SPeter Brune 
19258bdfa14SPeter Brune   if (snes->normschedule == SNES_NORM_ALWAYS) {
19358bdfa14SPeter Brune     next = jac->head;
19458bdfa14SPeter Brune     ierr = SNESSetInitialFunction(next->snes,F);CHKERRQ(ierr);
19558bdfa14SPeter Brune     while (next->next) {
19658bdfa14SPeter Brune       next = next->next;
19758bdfa14SPeter Brune       ierr = SNESSetInitialFunction(next->snes,F);CHKERRQ(ierr);
19858bdfa14SPeter Brune     }
19958bdfa14SPeter Brune   }
20058bdfa14SPeter Brune 
20190a8ba9bSPeter Brune   next = jac->head;
20290a8ba9bSPeter Brune   i = 0;
20390a8ba9bSPeter Brune   ierr = VecCopy(X,Xes[i]);CHKERRQ(ierr);
20490a8ba9bSPeter Brune   ierr = SNESSolve(next->snes,B,Xes[i]);CHKERRQ(ierr);
205d5a53f06SPeter Brune   ierr = SNESGetConvergedReason(next->snes,&reason);CHKERRQ(ierr);
206d5a53f06SPeter Brune   if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
207d5a53f06SPeter Brune     snes->reason = SNES_DIVERGED_INNER;
208d5a53f06SPeter Brune     PetscFunctionReturn(0);
209d5a53f06SPeter Brune   }
21090a8ba9bSPeter Brune   while (next->next) {
21190a8ba9bSPeter Brune     i++;
21290a8ba9bSPeter Brune     next = next->next;
21390a8ba9bSPeter Brune     ierr = VecCopy(X,Xes[i]);CHKERRQ(ierr);
21490a8ba9bSPeter Brune     ierr = SNESSolve(next->snes,B,Xes[i]);CHKERRQ(ierr);
215d5a53f06SPeter Brune     ierr = SNESGetConvergedReason(next->snes,&reason);CHKERRQ(ierr);
216d5a53f06SPeter Brune     if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
217d5a53f06SPeter Brune       snes->reason = SNES_DIVERGED_INNER;
218d5a53f06SPeter Brune       PetscFunctionReturn(0);
219d5a53f06SPeter Brune     }
22090a8ba9bSPeter Brune   }
22190a8ba9bSPeter Brune 
22290a8ba9bSPeter Brune   /* all the solutions are collected; combine optimally */
22390a8ba9bSPeter Brune   for (i=0;i<jac->n;i++) {
22490a8ba9bSPeter Brune     for (j=0;j<i+1;j++) {
2259c2f3473SPeter Brune       ierr = VecDotBegin(Fes[i],Fes[j],&jac->h[i + j*jac->n]);CHKERRQ(ierr);
22690a8ba9bSPeter Brune     }
2279c2f3473SPeter Brune     ierr = VecDotBegin(Fes[i],F,&jac->g[i]);CHKERRQ(ierr);
22890a8ba9bSPeter Brune   }
2295e806d2eSPeter Brune 
23090a8ba9bSPeter Brune   for (i=0;i<jac->n;i++) {
23190a8ba9bSPeter Brune     for (j=0;j<i+1;j++) {
2329c2f3473SPeter Brune       ierr = VecDotEnd(Fes[i],Fes[j],&jac->h[i + j*jac->n]);CHKERRQ(ierr);
2339c2f3473SPeter Brune       if (i == j) jac->fnorms[i] = PetscSqrtReal(PetscRealPart(jac->h[i + j*jac->n]));
23490a8ba9bSPeter Brune     }
2359c2f3473SPeter Brune     ierr = VecDotEnd(Fes[i],F,&jac->g[i]);CHKERRQ(ierr);
23690a8ba9bSPeter Brune   }
23790a8ba9bSPeter Brune 
2389c2f3473SPeter Brune   ftf = (*fnorm)*(*fnorm);
23990a8ba9bSPeter Brune 
24090a8ba9bSPeter Brune   for (i=0; i<jac->n; i++) {
24190a8ba9bSPeter Brune     for (j=i+1;j<jac->n;j++) {
2429c2f3473SPeter Brune       jac->h[i + j*jac->n] = jac->h[j + i*jac->n];
24390a8ba9bSPeter Brune     }
24490a8ba9bSPeter Brune   }
24590a8ba9bSPeter Brune 
24690a8ba9bSPeter Brune   for (i=0; i<jac->n; i++) {
2479c2f3473SPeter Brune     for (j=0; j<jac->n; j++) {
2489c2f3473SPeter Brune       jac->h[i + j*jac->n] = jac->h[i + j*jac->n] - jac->g[j] - jac->g[i] + ftf;
24990a8ba9bSPeter Brune     }
2509c2f3473SPeter Brune     jac->beta[i] = ftf - jac->g[i];
2519c2f3473SPeter Brune   }
25290a8ba9bSPeter Brune 
25390a8ba9bSPeter Brune #if defined(PETSC_MISSING_LAPACK_GELSS)
25490a8ba9bSPeter Brune   SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_SUP,"SNESCOMPOSITE with ADDITIVEOPTIMAL requires the LAPACK GELSS routine.");
25590a8ba9bSPeter Brune #else
25690a8ba9bSPeter Brune   jac->info  = 0;
25790a8ba9bSPeter Brune   jac->rcond = -1.;
25890a8ba9bSPeter Brune   ierr          = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr);
25990a8ba9bSPeter Brune #if defined(PETSC_USE_COMPLEX)
2609c2f3473SPeter 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));
26190a8ba9bSPeter Brune #else
2629c2f3473SPeter 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));
26390a8ba9bSPeter Brune #endif
26490a8ba9bSPeter Brune   ierr = PetscFPTrapPop();CHKERRQ(ierr);
26590a8ba9bSPeter Brune   if (jac->info < 0) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_LIB,"Bad argument to GELSS");
26690a8ba9bSPeter Brune   if (jac->info > 0) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_LIB,"SVD failed to converge");
26790a8ba9bSPeter Brune #endif
2689c2f3473SPeter Brune   tot = 0.;
2695e806d2eSPeter Brune   total = 0.;
27090a8ba9bSPeter Brune   for (i=0; i<jac->n; i++) {
27190a8ba9bSPeter Brune     if (PetscIsInfOrNanScalar(jac->beta[i])) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_LIB,"SVD generated inconsistent output");
272b3000dc2SPeter Brune     ierr = PetscInfo2(snes,"%d: %f\n",i,PetscRealPart(jac->beta[i]));CHKERRQ(ierr);
2739c2f3473SPeter Brune     tot += jac->beta[i];
2745e806d2eSPeter Brune     total += PetscAbsScalar(jac->beta[i]);
27590a8ba9bSPeter Brune   }
2769c2f3473SPeter Brune   ierr = VecScale(X,(1. - tot));CHKERRQ(ierr);
27790a8ba9bSPeter Brune   ierr = VecMAXPY(X,jac->n,jac->beta,Xes);CHKERRQ(ierr);
27890a8ba9bSPeter Brune   ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
2799c2f3473SPeter Brune   ierr = VecNorm(F,NORM_2,fnorm);CHKERRQ(ierr);
28090a8ba9bSPeter Brune 
2815e806d2eSPeter Brune   /* take the minimum-normed candidate if it beats the combination by a factor of rtol or the combination has stagnated */
2825e806d2eSPeter Brune   min_fnorm = jac->fnorms[0];
2835e806d2eSPeter Brune   min_i     = 0;
2849c2f3473SPeter Brune   for (i=0; i<jac->n; i++) {
2855e806d2eSPeter Brune     if (jac->fnorms[i] < min_fnorm) {
2865e806d2eSPeter Brune       min_fnorm = jac->fnorms[i];
2875e806d2eSPeter Brune       min_i     = i;
2889c2f3473SPeter Brune     }
2899c2f3473SPeter Brune   }
2905e806d2eSPeter Brune 
2915e806d2eSPeter Brune   /* stagnation or divergence restart to the solution of the solver that failed the least */
2925e806d2eSPeter Brune   if (PetscRealPart(total) < jac->stol || min_fnorm*jac->rtol < *fnorm) {
29358bdfa14SPeter Brune     ierr = VecCopy(jac->Xes[min_i],X);CHKERRQ(ierr);
29458bdfa14SPeter Brune     ierr = VecCopy(jac->Fes[min_i],F);CHKERRQ(ierr);
2955e806d2eSPeter Brune     *fnorm = min_fnorm;
2965e806d2eSPeter Brune   }
29790a8ba9bSPeter Brune   PetscFunctionReturn(0);
29890a8ba9bSPeter Brune }
29990a8ba9bSPeter Brune 
300eed5f15bSPeter Brune #undef __FUNCT__
301eed5f15bSPeter Brune #define __FUNCT__ "SNESSetUp_Composite"
302eed5f15bSPeter Brune static PetscErrorCode SNESSetUp_Composite(SNES snes)
303eed5f15bSPeter Brune {
304eed5f15bSPeter Brune   PetscErrorCode     ierr;
305dd515a93SPeter Brune   DM                 dm;
306eed5f15bSPeter Brune   SNES_Composite     *jac = (SNES_Composite*)snes->data;
307eed5f15bSPeter Brune   SNES_CompositeLink next = jac->head;
30890a8ba9bSPeter Brune   PetscInt           n=0,i;
30990a8ba9bSPeter Brune   Vec                F;
310eed5f15bSPeter Brune 
311eed5f15bSPeter Brune   PetscFunctionBegin;
312dd515a93SPeter Brune   ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr);
313eed5f15bSPeter Brune   while (next) {
31490a8ba9bSPeter Brune     n++;
315dd515a93SPeter Brune     ierr = SNESSetDM(next->snes,dm);CHKERRQ(ierr);
316eed5f15bSPeter Brune     ierr = SNESSetFromOptions(next->snes);CHKERRQ(ierr);
317eed5f15bSPeter Brune     next = next->next;
318eed5f15bSPeter Brune   }
31990a8ba9bSPeter Brune   jac->nsnes = n;
32090a8ba9bSPeter Brune   ierr = SNESGetFunction(snes,&F,NULL,NULL);CHKERRQ(ierr);
32190a8ba9bSPeter Brune   if (jac->type == SNES_COMPOSITE_ADDITIVEOPTIMAL) {
32290a8ba9bSPeter Brune     ierr = VecDuplicateVecs(F,jac->nsnes,&jac->Xes);CHKERRQ(ierr);
32390a8ba9bSPeter Brune     ierr = PetscMalloc(sizeof(Vec)*n,&jac->Fes);CHKERRQ(ierr);
3249c2f3473SPeter Brune     ierr = PetscMalloc(sizeof(PetscReal)*n,&jac->fnorms);CHKERRQ(ierr);
32590a8ba9bSPeter Brune     next = jac->head;
32690a8ba9bSPeter Brune     i = 0;
32790a8ba9bSPeter Brune     while (next) {
32890a8ba9bSPeter Brune       ierr = SNESGetFunction(next->snes,&F,NULL,NULL);CHKERRQ(ierr);
32990a8ba9bSPeter Brune       jac->Fes[i] = F;
33090a8ba9bSPeter Brune       ierr = PetscObjectReference((PetscObject)F);CHKERRQ(ierr);
33190a8ba9bSPeter Brune       next = next->next;
33290a8ba9bSPeter Brune       i++;
33390a8ba9bSPeter Brune     }
33490a8ba9bSPeter Brune     /* allocate the subspace direct solve area */
33590a8ba9bSPeter Brune     jac->nrhs  = 1;
3369c2f3473SPeter Brune     jac->lda   = jac->nsnes;
33790a8ba9bSPeter Brune     jac->ldb   = jac->nsnes;
33890a8ba9bSPeter Brune     jac->n     = jac->nsnes;
33990a8ba9bSPeter Brune 
340785e854fSJed Brown     ierr = PetscMalloc1(jac->n*jac->n,&jac->h);CHKERRQ(ierr);
341785e854fSJed Brown     ierr = PetscMalloc1(jac->n,&jac->beta);CHKERRQ(ierr);
342785e854fSJed Brown     ierr = PetscMalloc1(jac->n,&jac->s);CHKERRQ(ierr);
343785e854fSJed Brown     ierr = PetscMalloc1(jac->n,&jac->g);CHKERRQ(ierr);
3449c2f3473SPeter Brune     jac->lwork = 12*jac->n;
34590a8ba9bSPeter Brune #if PETSC_USE_COMPLEX
34690a8ba9bSPeter Brune     ierr = PetscMalloc(sizeof(PetscReal)*jac->lwork,&jac->rwork);CHKERRQ(ierr);
34790a8ba9bSPeter Brune #endif
34890a8ba9bSPeter Brune     ierr = PetscMalloc(sizeof(PetscScalar)*jac->lwork,&jac->work);CHKERRQ(ierr);
34990a8ba9bSPeter Brune   }
35090a8ba9bSPeter Brune 
351eed5f15bSPeter Brune   PetscFunctionReturn(0);
352eed5f15bSPeter Brune }
353eed5f15bSPeter Brune 
354eed5f15bSPeter Brune #undef __FUNCT__
355eed5f15bSPeter Brune #define __FUNCT__ "SNESReset_Composite"
356eed5f15bSPeter Brune static PetscErrorCode SNESReset_Composite(SNES snes)
357eed5f15bSPeter Brune {
358eed5f15bSPeter Brune   SNES_Composite     *jac = (SNES_Composite*)snes->data;
359eed5f15bSPeter Brune   PetscErrorCode   ierr;
360eed5f15bSPeter Brune   SNES_CompositeLink next = jac->head;
361eed5f15bSPeter Brune 
362eed5f15bSPeter Brune   PetscFunctionBegin;
363eed5f15bSPeter Brune   while (next) {
364eed5f15bSPeter Brune     ierr = SNESReset(next->snes);CHKERRQ(ierr);
365eed5f15bSPeter Brune     next = next->next;
366eed5f15bSPeter Brune   }
367eed5f15bSPeter Brune   ierr = VecDestroy(&jac->Xorig);CHKERRQ(ierr);
36890a8ba9bSPeter Brune   if (jac->Xes) {ierr = VecDestroyVecs(jac->nsnes,&jac->Xes);CHKERRQ(ierr);}
36990a8ba9bSPeter Brune   if (jac->Fes) {ierr = VecDestroyVecs(jac->nsnes,&jac->Fes);CHKERRQ(ierr);}
3709c2f3473SPeter Brune   ierr = PetscFree(jac->fnorms);CHKERRQ(ierr);
37190a8ba9bSPeter Brune   ierr = PetscFree(jac->h);CHKERRQ(ierr);
37290a8ba9bSPeter Brune   ierr = PetscFree(jac->s);CHKERRQ(ierr);
3739c2f3473SPeter Brune   ierr = PetscFree(jac->g);CHKERRQ(ierr);
37490a8ba9bSPeter Brune   ierr = PetscFree(jac->beta);CHKERRQ(ierr);
37590a8ba9bSPeter Brune   ierr = PetscFree(jac->work);CHKERRQ(ierr);
37690a8ba9bSPeter Brune   ierr = PetscFree(jac->rwork);CHKERRQ(ierr);
37790a8ba9bSPeter Brune 
378eed5f15bSPeter Brune   PetscFunctionReturn(0);
379eed5f15bSPeter Brune }
380eed5f15bSPeter Brune 
381eed5f15bSPeter Brune #undef __FUNCT__
382eed5f15bSPeter Brune #define __FUNCT__ "SNESDestroy_Composite"
383eed5f15bSPeter Brune static PetscErrorCode SNESDestroy_Composite(SNES snes)
384eed5f15bSPeter Brune {
385eed5f15bSPeter Brune   SNES_Composite     *jac = (SNES_Composite*)snes->data;
386eed5f15bSPeter Brune   PetscErrorCode   ierr;
387eed5f15bSPeter Brune   SNES_CompositeLink next = jac->head,next_tmp;
388eed5f15bSPeter Brune 
389eed5f15bSPeter Brune   PetscFunctionBegin;
390eed5f15bSPeter Brune   ierr = SNESReset_Composite(snes);CHKERRQ(ierr);
391eed5f15bSPeter Brune   while (next) {
392eed5f15bSPeter Brune     ierr     = SNESDestroy(&next->snes);CHKERRQ(ierr);
393eed5f15bSPeter Brune     next_tmp = next;
394eed5f15bSPeter Brune     next     = next->next;
395eed5f15bSPeter Brune     ierr     = PetscFree(next_tmp);CHKERRQ(ierr);
396eed5f15bSPeter Brune   }
397eed5f15bSPeter Brune   ierr = PetscFree(snes->data);CHKERRQ(ierr);
398eed5f15bSPeter Brune   PetscFunctionReturn(0);
399eed5f15bSPeter Brune }
400eed5f15bSPeter Brune 
401eed5f15bSPeter Brune #undef __FUNCT__
402eed5f15bSPeter Brune #define __FUNCT__ "SNESSetFromOptions_Composite"
403eed5f15bSPeter Brune static PetscErrorCode SNESSetFromOptions_Composite(SNES snes)
404eed5f15bSPeter Brune {
405eed5f15bSPeter Brune   SNES_Composite     *jac = (SNES_Composite*)snes->data;
406eed5f15bSPeter Brune   PetscErrorCode     ierr;
407eed5f15bSPeter Brune   PetscInt           nmax = 8,i;
408eed5f15bSPeter Brune   SNES_CompositeLink next;
409eed5f15bSPeter Brune   char               *sneses[8];
4108f626970SPeter Brune   PetscReal          dmps[8];
411eed5f15bSPeter Brune   PetscBool          flg;
412eed5f15bSPeter Brune 
413eed5f15bSPeter Brune   PetscFunctionBegin;
414eed5f15bSPeter Brune   ierr = PetscOptionsHead("Composite preconditioner options");CHKERRQ(ierr);
415eed5f15bSPeter Brune   ierr = PetscOptionsEnum("-snes_composite_type","Type of composition","SNESCompositeSetType",SNESCompositeTypes,(PetscEnum)jac->type,(PetscEnum*)&jac->type,&flg);CHKERRQ(ierr);
416eed5f15bSPeter Brune   if (flg) {
417eed5f15bSPeter Brune     ierr = SNESCompositeSetType(snes,jac->type);CHKERRQ(ierr);
418eed5f15bSPeter Brune   }
419eed5f15bSPeter Brune   ierr = PetscOptionsStringArray("-snes_composite_sneses","List of composite solvers","SNESCompositeAddSNES",sneses,&nmax,&flg);CHKERRQ(ierr);
420eed5f15bSPeter Brune   if (flg) {
421eed5f15bSPeter Brune     for (i=0; i<nmax; i++) {
422eed5f15bSPeter Brune       ierr = SNESCompositeAddSNES(snes,sneses[i]);CHKERRQ(ierr);
423eed5f15bSPeter Brune       ierr = PetscFree(sneses[i]);CHKERRQ(ierr);   /* deallocate string sneses[i], which is allocated in PetscOptionsStringArray() */
424eed5f15bSPeter Brune     }
425eed5f15bSPeter Brune   }
4268f626970SPeter Brune   ierr = PetscOptionsRealArray("-snes_composite_damping","Damping of the additive composite solvers","SNESCompositeSetDamping",dmps,&nmax,&flg);CHKERRQ(ierr);
4278f626970SPeter Brune   if (flg) {
4288f626970SPeter Brune     for (i=0; i<nmax; i++) {
4298f626970SPeter Brune       ierr = SNESCompositeSetDamping(snes,i,dmps[i]);CHKERRQ(ierr);
4308f626970SPeter Brune     }
4318f626970SPeter Brune   }
4325e806d2eSPeter Brune   ierr = PetscOptionsReal("-snes_composite_stol","Step tolerance for restart on the additive composite solvers","",jac->stol,&jac->stol,NULL);CHKERRQ(ierr);
4335e806d2eSPeter Brune   ierr = PetscOptionsReal("-snes_composite_rtol","Residual tolerance for the additive composite solvers","",jac->rtol,&jac->rtol,NULL);CHKERRQ(ierr);
434eed5f15bSPeter Brune   ierr = PetscOptionsTail();CHKERRQ(ierr);
435eed5f15bSPeter Brune 
436eed5f15bSPeter Brune   next = jac->head;
437eed5f15bSPeter Brune   while (next) {
438eed5f15bSPeter Brune     ierr = SNESSetFromOptions(next->snes);CHKERRQ(ierr);
439eed5f15bSPeter Brune     next = next->next;
440eed5f15bSPeter Brune   }
441eed5f15bSPeter Brune   PetscFunctionReturn(0);
442eed5f15bSPeter Brune }
443eed5f15bSPeter Brune 
444eed5f15bSPeter Brune #undef __FUNCT__
445eed5f15bSPeter Brune #define __FUNCT__ "SNESView_Composite"
446eed5f15bSPeter Brune static PetscErrorCode SNESView_Composite(SNES snes,PetscViewer viewer)
447eed5f15bSPeter Brune {
448eed5f15bSPeter Brune   SNES_Composite     *jac = (SNES_Composite*)snes->data;
449eed5f15bSPeter Brune   PetscErrorCode   ierr;
450eed5f15bSPeter Brune   SNES_CompositeLink next = jac->head;
451eed5f15bSPeter Brune   PetscBool        iascii;
452eed5f15bSPeter Brune 
453eed5f15bSPeter Brune   PetscFunctionBegin;
454eed5f15bSPeter Brune   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
455eed5f15bSPeter Brune   if (iascii) {
456eed5f15bSPeter Brune     ierr = PetscViewerASCIIPrintf(viewer,"Composite SNES type - %s\n",SNESCompositeTypes[jac->type]);CHKERRQ(ierr);
457eed5f15bSPeter Brune     ierr = PetscViewerASCIIPrintf(viewer,"SNESes on composite preconditioner follow\n");CHKERRQ(ierr);
458eed5f15bSPeter Brune     ierr = PetscViewerASCIIPrintf(viewer,"---------------------------------\n");CHKERRQ(ierr);
459eed5f15bSPeter Brune   }
460eed5f15bSPeter Brune   if (iascii) {
461eed5f15bSPeter Brune     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
462eed5f15bSPeter Brune   }
463eed5f15bSPeter Brune   while (next) {
464eed5f15bSPeter Brune     ierr = SNESView(next->snes,viewer);CHKERRQ(ierr);
465eed5f15bSPeter Brune     next = next->next;
466eed5f15bSPeter Brune   }
467eed5f15bSPeter Brune   if (iascii) {
468eed5f15bSPeter Brune     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
469eed5f15bSPeter Brune     ierr = PetscViewerASCIIPrintf(viewer,"---------------------------------\n");CHKERRQ(ierr);
470eed5f15bSPeter Brune   }
471eed5f15bSPeter Brune   PetscFunctionReturn(0);
472eed5f15bSPeter Brune }
473eed5f15bSPeter Brune 
474eed5f15bSPeter Brune /* ------------------------------------------------------------------------------*/
475eed5f15bSPeter Brune 
476eed5f15bSPeter Brune #undef __FUNCT__
477eed5f15bSPeter Brune #define __FUNCT__ "SNESCompositeSetType_Composite"
478eed5f15bSPeter Brune static PetscErrorCode  SNESCompositeSetType_Composite(SNES snes,SNESCompositeType type)
479eed5f15bSPeter Brune {
480eed5f15bSPeter Brune   SNES_Composite *jac = (SNES_Composite*)snes->data;
481eed5f15bSPeter Brune 
482eed5f15bSPeter Brune   PetscFunctionBegin;
483eed5f15bSPeter Brune   jac->type = type;
484eed5f15bSPeter Brune   PetscFunctionReturn(0);
485eed5f15bSPeter Brune }
486eed5f15bSPeter Brune 
487eed5f15bSPeter Brune #undef __FUNCT__
488eed5f15bSPeter Brune #define __FUNCT__ "SNESCompositeAddSNES_Composite"
489eed5f15bSPeter Brune static PetscErrorCode  SNESCompositeAddSNES_Composite(SNES snes,SNESType type)
490eed5f15bSPeter Brune {
491eed5f15bSPeter Brune   SNES_Composite     *jac;
492eed5f15bSPeter Brune   SNES_CompositeLink next,ilink;
493eed5f15bSPeter Brune   PetscErrorCode   ierr;
494eed5f15bSPeter Brune   PetscInt         cnt = 0;
495eed5f15bSPeter Brune   const char       *prefix;
496eed5f15bSPeter Brune   char             newprefix[8];
497eed5f15bSPeter Brune   DM               dm;
498eed5f15bSPeter Brune 
499eed5f15bSPeter Brune   PetscFunctionBegin;
500*b00a9115SJed Brown   ierr        = PetscNewLog(snes,&ilink);CHKERRQ(ierr);
501eed5f15bSPeter Brune   ilink->next = 0;
502eed5f15bSPeter Brune   ierr        = SNESCreate(PetscObjectComm((PetscObject)snes),&ilink->snes);CHKERRQ(ierr);
503eed5f15bSPeter Brune   ierr        = PetscLogObjectParent((PetscObject)snes,(PetscObject)ilink->snes);CHKERRQ(ierr);
504eed5f15bSPeter Brune   ierr        = SNESGetDM(snes,&dm);CHKERRQ(ierr);
505eed5f15bSPeter Brune   ierr        = SNESSetDM(ilink->snes,dm);CHKERRQ(ierr);
506cf5b3eb5SPeter Brune   ierr        = SNESSetTolerances(ilink->snes,ilink->snes->abstol,ilink->snes->rtol,ilink->snes->stol,1,ilink->snes->max_funcs);CHKERRQ(ierr);
507eed5f15bSPeter Brune   jac  = (SNES_Composite*)snes->data;
508eed5f15bSPeter Brune   next = jac->head;
509eed5f15bSPeter Brune   if (!next) {
510eed5f15bSPeter Brune     jac->head       = ilink;
511eed5f15bSPeter Brune     ilink->previous = NULL;
512eed5f15bSPeter Brune   } else {
513eed5f15bSPeter Brune     cnt++;
514eed5f15bSPeter Brune     while (next->next) {
515eed5f15bSPeter Brune       next = next->next;
516eed5f15bSPeter Brune       cnt++;
517eed5f15bSPeter Brune     }
518eed5f15bSPeter Brune     next->next      = ilink;
519eed5f15bSPeter Brune     ilink->previous = next;
520eed5f15bSPeter Brune   }
521eed5f15bSPeter Brune   ierr = SNESGetOptionsPrefix(snes,&prefix);CHKERRQ(ierr);
522eed5f15bSPeter Brune   ierr = SNESSetOptionsPrefix(ilink->snes,prefix);CHKERRQ(ierr);
523eed5f15bSPeter Brune   sprintf(newprefix,"sub_%d_",(int)cnt);
524eed5f15bSPeter Brune   ierr = SNESAppendOptionsPrefix(ilink->snes,newprefix);CHKERRQ(ierr);
5258f626970SPeter Brune   ierr = PetscObjectIncrementTabLevel((PetscObject)ilink->snes,(PetscObject)snes,1);CHKERRQ(ierr);
526eed5f15bSPeter Brune   ierr = SNESSetType(ilink->snes,type);CHKERRQ(ierr);
52772edecb9SPeter Brune   ierr = SNESSetNormSchedule(ilink->snes, SNES_NORM_FINAL_ONLY);CHKERRQ(ierr);
5288f626970SPeter Brune   ilink->dmp = 1.0;
52990a8ba9bSPeter Brune   jac->nsnes++;
530eed5f15bSPeter Brune   PetscFunctionReturn(0);
531eed5f15bSPeter Brune }
532eed5f15bSPeter Brune 
533eed5f15bSPeter Brune #undef __FUNCT__
534eed5f15bSPeter Brune #define __FUNCT__ "SNESCompositeGetSNES_Composite"
535eed5f15bSPeter Brune static PetscErrorCode  SNESCompositeGetSNES_Composite(SNES snes,PetscInt n,SNES *subsnes)
536eed5f15bSPeter Brune {
537eed5f15bSPeter Brune   SNES_Composite     *jac;
538eed5f15bSPeter Brune   SNES_CompositeLink next;
539eed5f15bSPeter Brune   PetscInt         i;
540eed5f15bSPeter Brune 
541eed5f15bSPeter Brune   PetscFunctionBegin;
542eed5f15bSPeter Brune   jac  = (SNES_Composite*)snes->data;
543eed5f15bSPeter Brune   next = jac->head;
544eed5f15bSPeter Brune   for (i=0; i<n; i++) {
545eed5f15bSPeter Brune     if (!next->next) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_INCOMP,"Not enough SNESes in composite preconditioner");
546eed5f15bSPeter Brune     next = next->next;
547eed5f15bSPeter Brune   }
548eed5f15bSPeter Brune   *subsnes = next->snes;
549eed5f15bSPeter Brune   PetscFunctionReturn(0);
550eed5f15bSPeter Brune }
551eed5f15bSPeter Brune 
552eed5f15bSPeter Brune /* -------------------------------------------------------------------------------- */
553eed5f15bSPeter Brune #undef __FUNCT__
554eed5f15bSPeter Brune #define __FUNCT__ "SNESCompositeSetType"
555e31ab4f9SPeter Brune /*@C
556eed5f15bSPeter Brune    SNESCompositeSetType - Sets the type of composite preconditioner.
557eed5f15bSPeter Brune 
558eed5f15bSPeter Brune    Logically Collective on SNES
559eed5f15bSPeter Brune 
560eed5f15bSPeter Brune    Input Parameter:
561eed5f15bSPeter Brune +  snes - the preconditioner context
562eed5f15bSPeter Brune -  type - SNES_COMPOSITE_ADDITIVE (default), SNES_COMPOSITE_MULTIPLICATIVE
563eed5f15bSPeter Brune 
564eed5f15bSPeter Brune    Options Database Key:
565eed5f15bSPeter Brune .  -snes_composite_type <type: one of multiplicative, additive, special> - Sets composite preconditioner type
566eed5f15bSPeter Brune 
567eed5f15bSPeter Brune    Level: Developer
568eed5f15bSPeter Brune 
569eed5f15bSPeter Brune .keywords: SNES, set, type, composite preconditioner, additive, multiplicative
570eed5f15bSPeter Brune @*/
571eed5f15bSPeter Brune PetscErrorCode  SNESCompositeSetType(SNES snes,SNESCompositeType type)
572eed5f15bSPeter Brune {
573eed5f15bSPeter Brune   PetscErrorCode ierr;
574eed5f15bSPeter Brune 
575eed5f15bSPeter Brune   PetscFunctionBegin;
576eed5f15bSPeter Brune   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
577eed5f15bSPeter Brune   PetscValidLogicalCollectiveEnum(snes,type,2);
578eed5f15bSPeter Brune   ierr = PetscTryMethod(snes,"SNESCompositeSetType_C",(SNES,SNESCompositeType),(snes,type));CHKERRQ(ierr);
579eed5f15bSPeter Brune   PetscFunctionReturn(0);
580eed5f15bSPeter Brune }
581eed5f15bSPeter Brune 
582eed5f15bSPeter Brune #undef __FUNCT__
583eed5f15bSPeter Brune #define __FUNCT__ "SNESCompositeAddSNES"
584eed5f15bSPeter Brune /*@C
585eed5f15bSPeter Brune    SNESCompositeAddSNES - Adds another SNES to the composite SNES.
586eed5f15bSPeter Brune 
587eed5f15bSPeter Brune    Collective on SNES
588eed5f15bSPeter Brune 
589eed5f15bSPeter Brune    Input Parameters:
590eed5f15bSPeter Brune +  snes - the preconditioner context
591eed5f15bSPeter Brune -  type - the type of the new preconditioner
592eed5f15bSPeter Brune 
593eed5f15bSPeter Brune    Level: Developer
594eed5f15bSPeter Brune 
595eed5f15bSPeter Brune .keywords: SNES, composite preconditioner, add
596eed5f15bSPeter Brune @*/
597eed5f15bSPeter Brune PetscErrorCode  SNESCompositeAddSNES(SNES snes,SNESType type)
598eed5f15bSPeter Brune {
599eed5f15bSPeter Brune   PetscErrorCode ierr;
600eed5f15bSPeter Brune 
601eed5f15bSPeter Brune   PetscFunctionBegin;
602eed5f15bSPeter Brune   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
603eed5f15bSPeter Brune   ierr = PetscTryMethod(snes,"SNESCompositeAddSNES_C",(SNES,SNESType),(snes,type));CHKERRQ(ierr);
604eed5f15bSPeter Brune   PetscFunctionReturn(0);
605eed5f15bSPeter Brune }
606eed5f15bSPeter Brune #undef __FUNCT__
607eed5f15bSPeter Brune #define __FUNCT__ "SNESCompositeGetSNES"
608eed5f15bSPeter Brune /*@
609eed5f15bSPeter Brune    SNESCompositeGetSNES - Gets one of the SNES objects in the composite SNES.
610eed5f15bSPeter Brune 
611eed5f15bSPeter Brune    Not Collective
612eed5f15bSPeter Brune 
613eed5f15bSPeter Brune    Input Parameter:
614eed5f15bSPeter Brune +  snes - the preconditioner context
615eed5f15bSPeter Brune -  n - the number of the snes requested
616eed5f15bSPeter Brune 
617eed5f15bSPeter Brune    Output Parameters:
618eed5f15bSPeter Brune .  subsnes - the SNES requested
619eed5f15bSPeter Brune 
620eed5f15bSPeter Brune    Level: Developer
621eed5f15bSPeter Brune 
622eed5f15bSPeter Brune .keywords: SNES, get, composite preconditioner, sub preconditioner
623eed5f15bSPeter Brune 
624eed5f15bSPeter Brune .seealso: SNESCompositeAddSNES()
625eed5f15bSPeter Brune @*/
626eed5f15bSPeter Brune PetscErrorCode  SNESCompositeGetSNES(SNES snes,PetscInt n,SNES *subsnes)
627eed5f15bSPeter Brune {
628eed5f15bSPeter Brune   PetscErrorCode ierr;
629eed5f15bSPeter Brune 
630eed5f15bSPeter Brune   PetscFunctionBegin;
631eed5f15bSPeter Brune   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
632eed5f15bSPeter Brune   PetscValidPointer(subsnes,3);
633eed5f15bSPeter Brune   ierr = PetscUseMethod(snes,"SNESCompositeGetSNES_C",(SNES,PetscInt,SNES*),(snes,n,subsnes));CHKERRQ(ierr);
634eed5f15bSPeter Brune   PetscFunctionReturn(0);
635eed5f15bSPeter Brune }
636eed5f15bSPeter Brune 
637eed5f15bSPeter Brune #undef __FUNCT__
6388f626970SPeter Brune #define __FUNCT__ "SNESCompositeSetDamping_Composite"
6398f626970SPeter Brune static PetscErrorCode  SNESCompositeSetDamping_Composite(SNES snes,PetscInt n,PetscReal dmp)
6408f626970SPeter Brune {
6418f626970SPeter Brune   SNES_Composite     *jac;
6428f626970SPeter Brune   SNES_CompositeLink next;
6438f626970SPeter Brune   PetscInt         i;
6448f626970SPeter Brune 
6458f626970SPeter Brune   PetscFunctionBegin;
6468f626970SPeter Brune   jac  = (SNES_Composite*)snes->data;
6478f626970SPeter Brune   next = jac->head;
6488f626970SPeter Brune   for (i=0; i<n; i++) {
6498f626970SPeter Brune     if (!next->next) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_INCOMP,"Not enough SNESes in composite preconditioner");
6508f626970SPeter Brune     next = next->next;
6518f626970SPeter Brune   }
6528f626970SPeter Brune   next->dmp = dmp;
6538f626970SPeter Brune   PetscFunctionReturn(0);
6548f626970SPeter Brune }
6558f626970SPeter Brune 
6568f626970SPeter Brune #undef __FUNCT__
6578f626970SPeter Brune #define __FUNCT__ "SNESCompositeSetDamping"
6588f626970SPeter Brune /*@
6598f626970SPeter Brune    SNESCompositeSetDamping - Sets the damping of a subsolver when using additive composite SNES.
6608f626970SPeter Brune 
6618f626970SPeter Brune    Not Collective
6628f626970SPeter Brune 
6638f626970SPeter Brune    Input Parameter:
6648f626970SPeter Brune +  snes - the preconditioner context
6658f626970SPeter Brune .  n - the number of the snes requested
6668f626970SPeter Brune -  dmp - the damping
6678f626970SPeter Brune 
6688f626970SPeter Brune    Level: Developer
6698f626970SPeter Brune 
6708f626970SPeter Brune .keywords: SNES, get, composite preconditioner, sub preconditioner
6718f626970SPeter Brune 
6728f626970SPeter Brune .seealso: SNESCompositeAddSNES()
6738f626970SPeter Brune @*/
6748f626970SPeter Brune PetscErrorCode  SNESCompositeSetDamping(SNES snes,PetscInt n,PetscReal dmp)
6758f626970SPeter Brune {
6768f626970SPeter Brune   PetscErrorCode ierr;
6778f626970SPeter Brune 
6788f626970SPeter Brune   PetscFunctionBegin;
6798f626970SPeter Brune   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
6808f626970SPeter Brune   ierr = PetscUseMethod(snes,"SNESCompositeSetDamping_C",(SNES,PetscInt,PetscReal),(snes,n,dmp));CHKERRQ(ierr);
6818f626970SPeter Brune   PetscFunctionReturn(0);
6828f626970SPeter Brune }
6838f626970SPeter Brune 
6848f626970SPeter Brune #undef __FUNCT__
685eed5f15bSPeter Brune #define __FUNCT__ "SNESSolve_Composite"
686eed5f15bSPeter Brune PetscErrorCode SNESSolve_Composite(SNES snes)
687eed5f15bSPeter Brune {
688eed5f15bSPeter Brune   Vec            F;
689eed5f15bSPeter Brune   Vec            X;
690eed5f15bSPeter Brune   Vec            B;
691eed5f15bSPeter Brune   PetscInt       i;
692eed5f15bSPeter Brune   PetscReal      fnorm = 0.0;
693eed5f15bSPeter Brune   PetscErrorCode ierr;
69472edecb9SPeter Brune   SNESNormSchedule normtype;
695eed5f15bSPeter Brune   SNES_Composite *comp = (SNES_Composite*)snes->data;
696eed5f15bSPeter Brune 
697eed5f15bSPeter Brune   PetscFunctionBegin;
698eed5f15bSPeter Brune   X = snes->vec_sol;
699eed5f15bSPeter Brune   F = snes->vec_func;
700eed5f15bSPeter Brune   B = snes->vec_rhs;
701eed5f15bSPeter Brune 
702e04113cfSBarry Smith   ierr         = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr);
703eed5f15bSPeter Brune   snes->iter   = 0;
704eed5f15bSPeter Brune   snes->norm   = 0.;
705e04113cfSBarry Smith   ierr         = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr);
706eed5f15bSPeter Brune   snes->reason = SNES_CONVERGED_ITERATING;
70772edecb9SPeter Brune   ierr         = SNESGetNormSchedule(snes, &normtype);CHKERRQ(ierr);
70872edecb9SPeter Brune   if (normtype == SNES_NORM_ALWAYS || normtype == SNES_NORM_INITIAL_ONLY || normtype == SNES_NORM_INITIAL_FINAL_ONLY) {
709eed5f15bSPeter Brune     if (!snes->vec_func_init_set) {
710eed5f15bSPeter Brune       ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
711eed5f15bSPeter Brune       if (snes->domainerror) {
712eed5f15bSPeter Brune         snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
713eed5f15bSPeter Brune         PetscFunctionReturn(0);
714eed5f15bSPeter Brune       }
715eed5f15bSPeter Brune     } else snes->vec_func_init_set = PETSC_FALSE;
716eed5f15bSPeter Brune 
717eed5f15bSPeter Brune     ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F||  */
718eed5f15bSPeter Brune     if (PetscIsInfOrNanReal(fnorm)) {
719eed5f15bSPeter Brune       snes->reason = SNES_DIVERGED_FNORM_NAN;
720eed5f15bSPeter Brune       PetscFunctionReturn(0);
721eed5f15bSPeter Brune     }
722e04113cfSBarry Smith     ierr       = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr);
723eed5f15bSPeter Brune     snes->iter = 0;
724eed5f15bSPeter Brune     snes->norm = fnorm;
725e04113cfSBarry Smith     ierr       = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr);
726eed5f15bSPeter Brune     ierr       = SNESLogConvergenceHistory(snes,snes->norm,0);CHKERRQ(ierr);
727eed5f15bSPeter Brune     ierr       = SNESMonitor(snes,0,snes->norm);CHKERRQ(ierr);
728eed5f15bSPeter Brune 
729eed5f15bSPeter Brune     /* test convergence */
730eed5f15bSPeter Brune     ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
731eed5f15bSPeter Brune     if (snes->reason) PetscFunctionReturn(0);
732eed5f15bSPeter Brune   } else {
733e04113cfSBarry Smith     ierr = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr);
734eed5f15bSPeter Brune     ierr = SNESLogConvergenceHistory(snes,snes->norm,0);CHKERRQ(ierr);
735eed5f15bSPeter Brune     ierr = SNESMonitor(snes,0,snes->norm);CHKERRQ(ierr);
736eed5f15bSPeter Brune   }
737eed5f15bSPeter Brune 
738eed5f15bSPeter Brune   /* Call general purpose update function */
739eed5f15bSPeter Brune   if (snes->ops->update) {
740eed5f15bSPeter Brune     ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
741eed5f15bSPeter Brune   }
742eed5f15bSPeter Brune 
743eed5f15bSPeter Brune   for (i = 0; i < snes->max_its; i++) {
744eed5f15bSPeter Brune     if (comp->type == SNES_COMPOSITE_ADDITIVE) {
745eed5f15bSPeter Brune       ierr = SNESCompositeApply_Additive(snes,X,B,F,&fnorm);CHKERRQ(ierr);
746eed5f15bSPeter Brune     } else if (comp->type == SNES_COMPOSITE_MULTIPLICATIVE) {
747eed5f15bSPeter Brune       ierr = SNESCompositeApply_Multiplicative(snes,X,B,F,&fnorm);CHKERRQ(ierr);
74890a8ba9bSPeter Brune     } else if (comp->type == SNES_COMPOSITE_ADDITIVEOPTIMAL) {
74990a8ba9bSPeter Brune       ierr = SNESCompositeApply_AdditiveOptimal(snes,X,B,F,&fnorm);CHKERRQ(ierr);
750eed5f15bSPeter Brune     } else {
75190a8ba9bSPeter Brune       SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"Unsupported SNESComposite type");
752eed5f15bSPeter Brune     }
753d5a53f06SPeter Brune     if (snes->reason < 0) break;
754d5a53f06SPeter Brune 
755eed5f15bSPeter Brune     if ((i == snes->max_its - 1) && (normtype == SNES_NORM_INITIAL_FINAL_ONLY || normtype == SNES_NORM_FINAL_ONLY)) {
756eed5f15bSPeter Brune       ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
757eed5f15bSPeter Brune       if (snes->domainerror) {
758eed5f15bSPeter Brune         snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
759eed5f15bSPeter Brune         break;
760eed5f15bSPeter Brune       }
761d5a53f06SPeter Brune       ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr);
762eed5f15bSPeter Brune       if (PetscIsInfOrNanReal(fnorm)) {
763eed5f15bSPeter Brune         snes->reason = SNES_DIVERGED_FNORM_NAN;
764eed5f15bSPeter Brune         break;
765eed5f15bSPeter Brune       }
766eed5f15bSPeter Brune     }
767eed5f15bSPeter Brune     /* Monitor convergence */
768e04113cfSBarry Smith     ierr       = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr);
769eed5f15bSPeter Brune     snes->iter = i+1;
770eed5f15bSPeter Brune     snes->norm = fnorm;
771e04113cfSBarry Smith     ierr       = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr);
772eed5f15bSPeter Brune     ierr       = SNESLogConvergenceHistory(snes,snes->norm,0);CHKERRQ(ierr);
773eed5f15bSPeter Brune     ierr       = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr);
774eed5f15bSPeter Brune     /* Test for convergence */
77572edecb9SPeter Brune     if (normtype == SNES_NORM_ALWAYS) {ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);}
776eed5f15bSPeter Brune     if (snes->reason) break;
777eed5f15bSPeter Brune     /* Call general purpose update function */
778eed5f15bSPeter Brune     if (snes->ops->update) {ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);}
779eed5f15bSPeter Brune   }
78072edecb9SPeter Brune   if (normtype == SNES_NORM_ALWAYS) {
781eed5f15bSPeter Brune     if (i == snes->max_its) {
782eed5f15bSPeter Brune       ierr = PetscInfo1(snes,"Maximum number of iterations has been reached: %D\n",snes->max_its);CHKERRQ(ierr);
783eed5f15bSPeter Brune       if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
784eed5f15bSPeter Brune     }
785eed5f15bSPeter Brune   } else if (!snes->reason) snes->reason = SNES_CONVERGED_ITS;
786eed5f15bSPeter Brune   PetscFunctionReturn(0);
787eed5f15bSPeter Brune }
788eed5f15bSPeter Brune 
789eed5f15bSPeter Brune /* -------------------------------------------------------------------------------------------*/
790eed5f15bSPeter Brune 
791eed5f15bSPeter Brune /*MC
792eed5f15bSPeter Brune      SNESCOMPOSITE - Build a preconditioner by composing together several nonlinear solvers
793eed5f15bSPeter Brune 
794eed5f15bSPeter Brune    Options Database Keys:
795eed5f15bSPeter Brune +  -snes_composite_type <type: one of multiplicative, additive, symmetric_multiplicative, special> - Sets composite preconditioner type
796eed5f15bSPeter Brune -  -snes_composite_sneses - <snes0,snes1,...> list of SNESes to compose
797eed5f15bSPeter Brune 
798eed5f15bSPeter Brune    Level: intermediate
799eed5f15bSPeter Brune 
800eed5f15bSPeter Brune    Concepts: composing solvers
801eed5f15bSPeter Brune 
802eed5f15bSPeter Brune .seealso:  SNESCreate(), SNESSetType(), SNESType (for list of available types), SNES,
803eed5f15bSPeter Brune            SNESSHELL, SNESCompositeSetType(), SNESCompositeSpecialSetAlpha(), SNESCompositeAddSNES(),
804eed5f15bSPeter Brune            SNESCompositeGetSNES()
805eed5f15bSPeter Brune 
806eed5f15bSPeter Brune M*/
807eed5f15bSPeter Brune 
808eed5f15bSPeter Brune #undef __FUNCT__
809eed5f15bSPeter Brune #define __FUNCT__ "SNESCreate_Composite"
810eed5f15bSPeter Brune PETSC_EXTERN PetscErrorCode SNESCreate_Composite(SNES snes)
811eed5f15bSPeter Brune {
812eed5f15bSPeter Brune   PetscErrorCode ierr;
813eed5f15bSPeter Brune   SNES_Composite   *jac;
814eed5f15bSPeter Brune 
815eed5f15bSPeter Brune   PetscFunctionBegin;
816*b00a9115SJed Brown   ierr = PetscNewLog(snes,&jac);CHKERRQ(ierr);
817eed5f15bSPeter Brune 
818eed5f15bSPeter Brune   snes->ops->solve           = SNESSolve_Composite;
819eed5f15bSPeter Brune   snes->ops->setup           = SNESSetUp_Composite;
820eed5f15bSPeter Brune   snes->ops->reset           = SNESReset_Composite;
821eed5f15bSPeter Brune   snes->ops->destroy         = SNESDestroy_Composite;
822eed5f15bSPeter Brune   snes->ops->setfromoptions  = SNESSetFromOptions_Composite;
823eed5f15bSPeter Brune   snes->ops->view            = SNESView_Composite;
824eed5f15bSPeter Brune 
825eed5f15bSPeter Brune   snes->data = (void*)jac;
82690a8ba9bSPeter Brune   jac->type  = SNES_COMPOSITE_ADDITIVEOPTIMAL;
82790a8ba9bSPeter Brune   jac->Fes   = NULL;
82890a8ba9bSPeter Brune   jac->Xes   = NULL;
8299c2f3473SPeter Brune   jac->fnorms = NULL;
83090a8ba9bSPeter Brune   jac->nsnes = 0;
831eed5f15bSPeter Brune   jac->head  = 0;
8325e806d2eSPeter Brune   jac->stol  = 0.1;
8335e806d2eSPeter Brune   jac->rtol  = 1.1;
834eed5f15bSPeter Brune 
83590a8ba9bSPeter Brune   jac->h     = NULL;
83690a8ba9bSPeter Brune   jac->s     = NULL;
83790a8ba9bSPeter Brune   jac->beta  = NULL;
83890a8ba9bSPeter Brune   jac->work  = NULL;
83990a8ba9bSPeter Brune   jac->rwork = NULL;
84090a8ba9bSPeter Brune 
8418f626970SPeter Brune   ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESCompositeSetType_C",SNESCompositeSetType_Composite);CHKERRQ(ierr);
8428f626970SPeter Brune   ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESCompositeAddSNES_C",SNESCompositeAddSNES_Composite);CHKERRQ(ierr);
8438f626970SPeter Brune   ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESCompositeGetSNES_C",SNESCompositeGetSNES_Composite);CHKERRQ(ierr);
8448f626970SPeter Brune   ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESCompositeSetDamping_C",SNESCompositeSetDamping_Composite);CHKERRQ(ierr);
845eed5f15bSPeter Brune   PetscFunctionReturn(0);
846eed5f15bSPeter Brune }
847eed5f15bSPeter Brune 
848