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; 56*d5a53f06SPeter Brune SNESConvergedReason reason; 57eed5f15bSPeter Brune 58eed5f15bSPeter Brune PetscFunctionBegin; 59eed5f15bSPeter Brune if (!next) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"No composite SNESes supplied via SNESCompositeAddSNES() or -snes_composite_sneses"); 6072edecb9SPeter Brune if (snes->normschedule == SNES_NORM_ALWAYS) { 61eed5f15bSPeter Brune ierr = SNESSetInitialFunction(next->snes,F);CHKERRQ(ierr); 62eed5f15bSPeter Brune if (fnorm) {ierr = SNESSetInitialFunctionNorm(next->snes,*fnorm);CHKERRQ(ierr);} 63eed5f15bSPeter Brune } 64eed5f15bSPeter Brune ierr = SNESSolve(next->snes,B,X);CHKERRQ(ierr); 65*d5a53f06SPeter Brune ierr = SNESGetConvergedReason(next->snes,&reason);CHKERRQ(ierr); 66*d5a53f06SPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 67*d5a53f06SPeter Brune snes->reason = SNES_DIVERGED_INNER; 68*d5a53f06SPeter Brune PetscFunctionReturn(0); 69*d5a53f06SPeter Brune } 70eed5f15bSPeter Brune 71eed5f15bSPeter Brune while (next->next) { 72eed5f15bSPeter Brune /* only copy the function over in the case where the functions correspond */ 7372edecb9SPeter Brune if (next->snes->pcside == PC_RIGHT && next->snes->normschedule != SNES_NORM_NONE) { 74eed5f15bSPeter Brune ierr = SNESGetFunction(next->snes,&FSub,NULL,NULL);CHKERRQ(ierr); 75eed5f15bSPeter Brune ierr = SNESGetFunctionNorm(next->snes,&fsubnorm);CHKERRQ(ierr); 76eed5f15bSPeter Brune next = next->next; 77eed5f15bSPeter Brune ierr = SNESSetInitialFunction(next->snes,FSub);CHKERRQ(ierr); 78eed5f15bSPeter Brune ierr = SNESSetInitialFunctionNorm(next->snes,fsubnorm);CHKERRQ(ierr); 79eed5f15bSPeter Brune } else { 80eed5f15bSPeter Brune next = next->next; 81eed5f15bSPeter Brune } 82eed5f15bSPeter Brune ierr = SNESSolve(next->snes,B,X);CHKERRQ(ierr); 83*d5a53f06SPeter Brune ierr = SNESGetConvergedReason(next->snes,&reason);CHKERRQ(ierr); 84*d5a53f06SPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 85*d5a53f06SPeter Brune snes->reason = SNES_DIVERGED_INNER; 86*d5a53f06SPeter Brune PetscFunctionReturn(0); 87*d5a53f06SPeter Brune } 88eed5f15bSPeter Brune } 89eed5f15bSPeter Brune if (next->snes->pcside == PC_RIGHT) { 90eed5f15bSPeter Brune ierr = SNESGetFunction(next->snes,&FSub,NULL,NULL);CHKERRQ(ierr); 91eed5f15bSPeter Brune ierr = VecCopy(FSub,F);CHKERRQ(ierr); 92*d5a53f06SPeter Brune if (fnorm) { 93*d5a53f06SPeter Brune ierr = SNESGetFunctionNorm(next->snes,fnorm);CHKERRQ(ierr); 94*d5a53f06SPeter Brune if (PetscIsInfOrNanReal(*fnorm)) { 95*d5a53f06SPeter Brune snes->reason = SNES_DIVERGED_FNORM_NAN; 96*d5a53f06SPeter Brune PetscFunctionReturn(0); 97*d5a53f06SPeter Brune } 98*d5a53f06SPeter Brune } 9972edecb9SPeter Brune } else if (snes->normschedule == SNES_NORM_ALWAYS) { 100eed5f15bSPeter Brune SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 101*d5a53f06SPeter Brune if (snes->domainerror) { 102*d5a53f06SPeter Brune snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 103*d5a53f06SPeter Brune PetscFunctionReturn(0); 104*d5a53f06SPeter Brune } 105*d5a53f06SPeter Brune if (fnorm) { 106*d5a53f06SPeter Brune ierr = VecNorm(F,NORM_2,fnorm);CHKERRQ(ierr); 107*d5a53f06SPeter Brune if (PetscIsInfOrNanReal(*fnorm)) { 108*d5a53f06SPeter Brune snes->reason = SNES_DIVERGED_FNORM_NAN; 109*d5a53f06SPeter Brune PetscFunctionReturn(0); 110*d5a53f06SPeter Brune } 111*d5a53f06SPeter Brune } 112eed5f15bSPeter Brune } 113eed5f15bSPeter Brune PetscFunctionReturn(0); 114eed5f15bSPeter Brune } 115eed5f15bSPeter Brune 116eed5f15bSPeter Brune #undef __FUNCT__ 117eed5f15bSPeter Brune #define __FUNCT__ "SNESCompositeApply_Additive" 118eed5f15bSPeter Brune static PetscErrorCode SNESCompositeApply_Additive(SNES snes,Vec X,Vec B,Vec F,PetscReal *fnorm) 119eed5f15bSPeter Brune { 120eed5f15bSPeter Brune PetscErrorCode ierr; 121eed5f15bSPeter Brune SNES_Composite *jac = (SNES_Composite*)snes->data; 122eed5f15bSPeter Brune SNES_CompositeLink next = jac->head; 123eed5f15bSPeter Brune Vec Y,Xorig; 124*d5a53f06SPeter Brune SNESConvergedReason reason; 125eed5f15bSPeter Brune 126eed5f15bSPeter Brune PetscFunctionBegin; 127eed5f15bSPeter Brune Y = snes->vec_sol_update; 128eed5f15bSPeter Brune if (!jac->Xorig) {ierr = VecDuplicate(X,&jac->Xorig);CHKERRQ(ierr);} 129eed5f15bSPeter Brune Xorig = jac->Xorig; 130eed5f15bSPeter Brune ierr = VecCopy(X,Xorig); 131eed5f15bSPeter Brune if (!next) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"No composite SNESes supplied via SNESCompositeAddSNES() or -snes_composite_sneses"); 13272edecb9SPeter Brune if (snes->normschedule == SNES_NORM_ALWAYS) { 133eed5f15bSPeter Brune ierr = SNESSetInitialFunction(next->snes,F);CHKERRQ(ierr); 134eed5f15bSPeter Brune if (fnorm) {ierr = SNESSetInitialFunctionNorm(next->snes,*fnorm);CHKERRQ(ierr);} 135eed5f15bSPeter Brune while (next->next) { 136eed5f15bSPeter Brune next = next->next; 137eed5f15bSPeter Brune ierr = SNESSetInitialFunction(next->snes,F);CHKERRQ(ierr); 138eed5f15bSPeter Brune if (fnorm) {ierr = SNESSetInitialFunctionNorm(next->snes,*fnorm);CHKERRQ(ierr);} 139eed5f15bSPeter Brune } 140eed5f15bSPeter Brune } 141eed5f15bSPeter Brune next = jac->head; 1428f626970SPeter Brune ierr = VecCopy(Xorig,Y);CHKERRQ(ierr); 143eed5f15bSPeter Brune ierr = SNESSolve(next->snes,B,Y);CHKERRQ(ierr); 144*d5a53f06SPeter Brune ierr = SNESGetConvergedReason(next->snes,&reason);CHKERRQ(ierr); 145*d5a53f06SPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 146*d5a53f06SPeter Brune snes->reason = SNES_DIVERGED_INNER; 147*d5a53f06SPeter Brune PetscFunctionReturn(0); 148*d5a53f06SPeter Brune } 149eed5f15bSPeter Brune ierr = VecAXPY(Y,-1.0,Xorig);CHKERRQ(ierr); 1508f626970SPeter Brune ierr = VecAXPY(X,next->dmp,Y);CHKERRQ(ierr); 1518f626970SPeter Brune while (next->next) { 1528f626970SPeter Brune next = next->next; 1538f626970SPeter Brune ierr = VecCopy(Xorig,Y);CHKERRQ(ierr); 1548f626970SPeter Brune ierr = SNESSolve(next->snes,B,Y);CHKERRQ(ierr); 155*d5a53f06SPeter Brune ierr = SNESGetConvergedReason(next->snes,&reason);CHKERRQ(ierr); 156*d5a53f06SPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 157*d5a53f06SPeter Brune snes->reason = SNES_DIVERGED_INNER; 158*d5a53f06SPeter Brune PetscFunctionReturn(0); 159*d5a53f06SPeter Brune } 1608f626970SPeter Brune ierr = VecAXPY(Y,-1.0,Xorig);CHKERRQ(ierr); 1618f626970SPeter Brune ierr = VecAXPY(X,next->dmp,Y);CHKERRQ(ierr); 162eed5f15bSPeter Brune } 16372edecb9SPeter Brune if (snes->normschedule == SNES_NORM_ALWAYS) { 164eed5f15bSPeter Brune ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 165eed5f15bSPeter Brune if (fnorm) {ierr = VecNorm(F,NORM_2,fnorm);CHKERRQ(ierr);} 166eed5f15bSPeter Brune } 167eed5f15bSPeter Brune PetscFunctionReturn(0); 168eed5f15bSPeter Brune } 16990a8ba9bSPeter Brune 17090a8ba9bSPeter Brune #undef __FUNCT__ 17190a8ba9bSPeter Brune #define __FUNCT__ "SNESCompositeApply_AdditiveOptimal" 17290a8ba9bSPeter Brune /* approximately solve the overdetermined system: 17390a8ba9bSPeter Brune 17490a8ba9bSPeter Brune 2*F(x_i)\cdot F(\x_j)\alpha_i = 0 17590a8ba9bSPeter Brune \alpha_i = 1 17690a8ba9bSPeter Brune 17790a8ba9bSPeter Brune Which minimizes the L2 norm of the linearization of: 17890a8ba9bSPeter Brune ||F(\sum_i \alpha_i*x_i)||^2 17990a8ba9bSPeter Brune 18090a8ba9bSPeter Brune With the constraint that \sum_i\alpha_i = 1 18190a8ba9bSPeter Brune Where x_i is the solution from the ith subsolver. 18290a8ba9bSPeter Brune */ 18390a8ba9bSPeter Brune static PetscErrorCode SNESCompositeApply_AdditiveOptimal(SNES snes,Vec X,Vec B,Vec F,PetscReal *fnorm) 18490a8ba9bSPeter Brune { 18590a8ba9bSPeter Brune PetscErrorCode ierr; 18690a8ba9bSPeter Brune SNES_Composite *jac = (SNES_Composite*)snes->data; 18790a8ba9bSPeter Brune SNES_CompositeLink next = jac->head; 18890a8ba9bSPeter Brune Vec *Xes = jac->Xes,*Fes = jac->Fes; 18990a8ba9bSPeter Brune PetscInt i,j; 1905e806d2eSPeter Brune PetscScalar tot,total,ftf; 1915e806d2eSPeter Brune PetscReal min_fnorm; 1925e806d2eSPeter Brune PetscInt min_i; 193*d5a53f06SPeter Brune SNESConvergedReason reason; 19490a8ba9bSPeter Brune 19590a8ba9bSPeter Brune PetscFunctionBegin; 19690a8ba9bSPeter Brune if (!next) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"No composite SNESes supplied via SNESCompositeAddSNES() or -snes_composite_sneses"); 19758bdfa14SPeter Brune 19858bdfa14SPeter Brune if (snes->normschedule == SNES_NORM_ALWAYS) { 19958bdfa14SPeter Brune next = jac->head; 20058bdfa14SPeter Brune ierr = SNESSetInitialFunction(next->snes,F);CHKERRQ(ierr); 20158bdfa14SPeter Brune if (fnorm) {ierr = SNESSetInitialFunctionNorm(next->snes,*fnorm);CHKERRQ(ierr);} 20258bdfa14SPeter Brune while (next->next) { 20358bdfa14SPeter Brune next = next->next; 20458bdfa14SPeter Brune ierr = SNESSetInitialFunction(next->snes,F);CHKERRQ(ierr); 20558bdfa14SPeter Brune if (fnorm) {ierr = SNESSetInitialFunctionNorm(next->snes,*fnorm);CHKERRQ(ierr);} 20658bdfa14SPeter Brune } 20758bdfa14SPeter Brune } 20858bdfa14SPeter Brune 20990a8ba9bSPeter Brune next = jac->head; 21090a8ba9bSPeter Brune i = 0; 21190a8ba9bSPeter Brune ierr = VecCopy(X,Xes[i]);CHKERRQ(ierr); 21290a8ba9bSPeter Brune ierr = SNESSolve(next->snes,B,Xes[i]);CHKERRQ(ierr); 213*d5a53f06SPeter Brune ierr = SNESGetConvergedReason(next->snes,&reason);CHKERRQ(ierr); 214*d5a53f06SPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 215*d5a53f06SPeter Brune snes->reason = SNES_DIVERGED_INNER; 216*d5a53f06SPeter Brune PetscFunctionReturn(0); 217*d5a53f06SPeter Brune } 21890a8ba9bSPeter Brune while (next->next) { 21990a8ba9bSPeter Brune i++; 22090a8ba9bSPeter Brune next = next->next; 22190a8ba9bSPeter Brune ierr = VecCopy(X,Xes[i]);CHKERRQ(ierr); 22290a8ba9bSPeter Brune ierr = SNESSolve(next->snes,B,Xes[i]);CHKERRQ(ierr); 223*d5a53f06SPeter Brune ierr = SNESGetConvergedReason(next->snes,&reason);CHKERRQ(ierr); 224*d5a53f06SPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 225*d5a53f06SPeter Brune snes->reason = SNES_DIVERGED_INNER; 226*d5a53f06SPeter Brune PetscFunctionReturn(0); 227*d5a53f06SPeter Brune } 22890a8ba9bSPeter Brune } 22990a8ba9bSPeter Brune 23090a8ba9bSPeter Brune /* all the solutions are collected; combine optimally */ 23190a8ba9bSPeter Brune for (i=0;i<jac->n;i++) { 23290a8ba9bSPeter Brune for (j=0;j<i+1;j++) { 2339c2f3473SPeter Brune ierr = VecDotBegin(Fes[i],Fes[j],&jac->h[i + j*jac->n]);CHKERRQ(ierr); 23490a8ba9bSPeter Brune } 2359c2f3473SPeter Brune ierr = VecDotBegin(Fes[i],F,&jac->g[i]);CHKERRQ(ierr); 23690a8ba9bSPeter Brune } 2375e806d2eSPeter Brune 23890a8ba9bSPeter Brune for (i=0;i<jac->n;i++) { 23990a8ba9bSPeter Brune for (j=0;j<i+1;j++) { 2409c2f3473SPeter Brune ierr = VecDotEnd(Fes[i],Fes[j],&jac->h[i + j*jac->n]);CHKERRQ(ierr); 2419c2f3473SPeter Brune if (i == j) jac->fnorms[i] = PetscSqrtReal(PetscRealPart(jac->h[i + j*jac->n])); 24290a8ba9bSPeter Brune } 2439c2f3473SPeter Brune ierr = VecDotEnd(Fes[i],F,&jac->g[i]);CHKERRQ(ierr); 24490a8ba9bSPeter Brune } 24590a8ba9bSPeter Brune 2469c2f3473SPeter Brune ftf = (*fnorm)*(*fnorm); 24790a8ba9bSPeter Brune 24890a8ba9bSPeter Brune for (i=0; i<jac->n; i++) { 24990a8ba9bSPeter Brune for (j=i+1;j<jac->n;j++) { 2509c2f3473SPeter Brune jac->h[i + j*jac->n] = jac->h[j + i*jac->n]; 25190a8ba9bSPeter Brune } 25290a8ba9bSPeter Brune } 25390a8ba9bSPeter Brune 25490a8ba9bSPeter Brune for (i=0; i<jac->n; i++) { 2559c2f3473SPeter Brune for (j=0; j<jac->n; j++) { 2569c2f3473SPeter Brune jac->h[i + j*jac->n] = jac->h[i + j*jac->n] - jac->g[j] - jac->g[i] + ftf; 25790a8ba9bSPeter Brune } 2589c2f3473SPeter Brune jac->beta[i] = ftf - jac->g[i]; 2599c2f3473SPeter Brune } 26090a8ba9bSPeter Brune 26190a8ba9bSPeter Brune #if defined(PETSC_MISSING_LAPACK_GELSS) 26290a8ba9bSPeter Brune SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_SUP,"SNESCOMPOSITE with ADDITIVEOPTIMAL requires the LAPACK GELSS routine."); 26390a8ba9bSPeter Brune #else 26490a8ba9bSPeter Brune jac->info = 0; 26590a8ba9bSPeter Brune jac->rcond = -1.; 26690a8ba9bSPeter Brune ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr); 26790a8ba9bSPeter Brune #if defined(PETSC_USE_COMPLEX) 2689c2f3473SPeter 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)); 26990a8ba9bSPeter Brune #else 2709c2f3473SPeter 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)); 27190a8ba9bSPeter Brune #endif 27290a8ba9bSPeter Brune ierr = PetscFPTrapPop();CHKERRQ(ierr); 27390a8ba9bSPeter Brune if (jac->info < 0) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_LIB,"Bad argument to GELSS"); 27490a8ba9bSPeter Brune if (jac->info > 0) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_LIB,"SVD failed to converge"); 27590a8ba9bSPeter Brune #endif 2769c2f3473SPeter Brune tot = 0.; 2775e806d2eSPeter Brune total = 0.; 27890a8ba9bSPeter Brune for (i=0; i<jac->n; i++) { 27990a8ba9bSPeter Brune if (PetscIsInfOrNanScalar(jac->beta[i])) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_LIB,"SVD generated inconsistent output"); 280b3000dc2SPeter Brune ierr = PetscInfo2(snes,"%d: %f\n",i,PetscRealPart(jac->beta[i]));CHKERRQ(ierr); 2819c2f3473SPeter Brune tot += jac->beta[i]; 2825e806d2eSPeter Brune total += PetscAbsScalar(jac->beta[i]); 28390a8ba9bSPeter Brune } 2849c2f3473SPeter Brune ierr = VecScale(X,(1. - tot));CHKERRQ(ierr); 28590a8ba9bSPeter Brune ierr = VecMAXPY(X,jac->n,jac->beta,Xes);CHKERRQ(ierr); 28690a8ba9bSPeter Brune ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 2879c2f3473SPeter Brune ierr = VecNorm(F,NORM_2,fnorm);CHKERRQ(ierr); 28890a8ba9bSPeter Brune 2895e806d2eSPeter Brune /* take the minimum-normed candidate if it beats the combination by a factor of rtol or the combination has stagnated */ 2905e806d2eSPeter Brune min_fnorm = jac->fnorms[0]; 2915e806d2eSPeter Brune min_i = 0; 2929c2f3473SPeter Brune for (i=0; i<jac->n; i++) { 2935e806d2eSPeter Brune if (jac->fnorms[i] < min_fnorm) { 2945e806d2eSPeter Brune min_fnorm = jac->fnorms[i]; 2955e806d2eSPeter Brune min_i = i; 2969c2f3473SPeter Brune } 2979c2f3473SPeter Brune } 2985e806d2eSPeter Brune 2995e806d2eSPeter Brune /* stagnation or divergence restart to the solution of the solver that failed the least */ 3005e806d2eSPeter Brune if (PetscRealPart(total) < jac->stol || min_fnorm*jac->rtol < *fnorm) { 30158bdfa14SPeter Brune ierr = VecCopy(jac->Xes[min_i],X);CHKERRQ(ierr); 30258bdfa14SPeter Brune ierr = VecCopy(jac->Fes[min_i],F);CHKERRQ(ierr); 3035e806d2eSPeter Brune *fnorm = min_fnorm; 3045e806d2eSPeter Brune } 30590a8ba9bSPeter Brune PetscFunctionReturn(0); 30690a8ba9bSPeter Brune } 30790a8ba9bSPeter Brune 308eed5f15bSPeter Brune #undef __FUNCT__ 309eed5f15bSPeter Brune #define __FUNCT__ "SNESSetUp_Composite" 310eed5f15bSPeter Brune static PetscErrorCode SNESSetUp_Composite(SNES snes) 311eed5f15bSPeter Brune { 312eed5f15bSPeter Brune PetscErrorCode ierr; 313dd515a93SPeter Brune DM dm; 314eed5f15bSPeter Brune SNES_Composite *jac = (SNES_Composite*)snes->data; 315eed5f15bSPeter Brune SNES_CompositeLink next = jac->head; 31690a8ba9bSPeter Brune PetscInt n=0,i; 31790a8ba9bSPeter Brune Vec F; 318eed5f15bSPeter Brune 319eed5f15bSPeter Brune PetscFunctionBegin; 320dd515a93SPeter Brune ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 321eed5f15bSPeter Brune while (next) { 32290a8ba9bSPeter Brune n++; 323dd515a93SPeter Brune ierr = SNESSetDM(next->snes,dm);CHKERRQ(ierr); 324eed5f15bSPeter Brune ierr = SNESSetFromOptions(next->snes);CHKERRQ(ierr); 325eed5f15bSPeter Brune next = next->next; 326eed5f15bSPeter Brune } 32790a8ba9bSPeter Brune jac->nsnes = n; 32890a8ba9bSPeter Brune ierr = SNESGetFunction(snes,&F,NULL,NULL);CHKERRQ(ierr); 32990a8ba9bSPeter Brune if (jac->type == SNES_COMPOSITE_ADDITIVEOPTIMAL) { 33090a8ba9bSPeter Brune ierr = VecDuplicateVecs(F,jac->nsnes,&jac->Xes);CHKERRQ(ierr); 33190a8ba9bSPeter Brune ierr = PetscMalloc(sizeof(Vec)*n,&jac->Fes);CHKERRQ(ierr); 3329c2f3473SPeter Brune ierr = PetscMalloc(sizeof(PetscReal)*n,&jac->fnorms);CHKERRQ(ierr); 33390a8ba9bSPeter Brune next = jac->head; 33490a8ba9bSPeter Brune i = 0; 33590a8ba9bSPeter Brune while (next) { 33690a8ba9bSPeter Brune ierr = SNESGetFunction(next->snes,&F,NULL,NULL);CHKERRQ(ierr); 33790a8ba9bSPeter Brune jac->Fes[i] = F; 33890a8ba9bSPeter Brune ierr = PetscObjectReference((PetscObject)F);CHKERRQ(ierr); 33990a8ba9bSPeter Brune next = next->next; 34090a8ba9bSPeter Brune i++; 34190a8ba9bSPeter Brune } 34290a8ba9bSPeter Brune /* allocate the subspace direct solve area */ 34390a8ba9bSPeter Brune jac->nrhs = 1; 3449c2f3473SPeter Brune jac->lda = jac->nsnes; 34590a8ba9bSPeter Brune jac->ldb = jac->nsnes; 34690a8ba9bSPeter Brune jac->n = jac->nsnes; 34790a8ba9bSPeter Brune 3489c2f3473SPeter Brune ierr = PetscMalloc(jac->n*jac->n*sizeof(PetscScalar),&jac->h);CHKERRQ(ierr); 3499c2f3473SPeter Brune ierr = PetscMalloc(jac->n*sizeof(PetscScalar),&jac->beta);CHKERRQ(ierr); 3509c2f3473SPeter Brune ierr = PetscMalloc(jac->n*sizeof(PetscScalar),&jac->s);CHKERRQ(ierr); 3519c2f3473SPeter Brune ierr = PetscMalloc(jac->n*sizeof(PetscScalar),&jac->g);CHKERRQ(ierr); 3529c2f3473SPeter Brune jac->lwork = 12*jac->n; 35390a8ba9bSPeter Brune #if PETSC_USE_COMPLEX 35490a8ba9bSPeter Brune ierr = PetscMalloc(sizeof(PetscReal)*jac->lwork,&jac->rwork);CHKERRQ(ierr); 35590a8ba9bSPeter Brune #endif 35690a8ba9bSPeter Brune ierr = PetscMalloc(sizeof(PetscScalar)*jac->lwork,&jac->work);CHKERRQ(ierr); 35790a8ba9bSPeter Brune } 35890a8ba9bSPeter Brune 359eed5f15bSPeter Brune PetscFunctionReturn(0); 360eed5f15bSPeter Brune } 361eed5f15bSPeter Brune 362eed5f15bSPeter Brune #undef __FUNCT__ 363eed5f15bSPeter Brune #define __FUNCT__ "SNESReset_Composite" 364eed5f15bSPeter Brune static PetscErrorCode SNESReset_Composite(SNES snes) 365eed5f15bSPeter Brune { 366eed5f15bSPeter Brune SNES_Composite *jac = (SNES_Composite*)snes->data; 367eed5f15bSPeter Brune PetscErrorCode ierr; 368eed5f15bSPeter Brune SNES_CompositeLink next = jac->head; 369eed5f15bSPeter Brune 370eed5f15bSPeter Brune PetscFunctionBegin; 371eed5f15bSPeter Brune while (next) { 372eed5f15bSPeter Brune ierr = SNESReset(next->snes);CHKERRQ(ierr); 373eed5f15bSPeter Brune next = next->next; 374eed5f15bSPeter Brune } 375eed5f15bSPeter Brune ierr = VecDestroy(&jac->Xorig);CHKERRQ(ierr); 37690a8ba9bSPeter Brune if (jac->Xes) {ierr = VecDestroyVecs(jac->nsnes,&jac->Xes);CHKERRQ(ierr);} 37790a8ba9bSPeter Brune if (jac->Fes) {ierr = VecDestroyVecs(jac->nsnes,&jac->Fes);CHKERRQ(ierr);} 3789c2f3473SPeter Brune ierr = PetscFree(jac->fnorms);CHKERRQ(ierr); 37990a8ba9bSPeter Brune ierr = PetscFree(jac->h);CHKERRQ(ierr); 38090a8ba9bSPeter Brune ierr = PetscFree(jac->s);CHKERRQ(ierr); 3819c2f3473SPeter Brune ierr = PetscFree(jac->g);CHKERRQ(ierr); 38290a8ba9bSPeter Brune ierr = PetscFree(jac->beta);CHKERRQ(ierr); 38390a8ba9bSPeter Brune ierr = PetscFree(jac->work);CHKERRQ(ierr); 38490a8ba9bSPeter Brune ierr = PetscFree(jac->rwork);CHKERRQ(ierr); 38590a8ba9bSPeter Brune 386eed5f15bSPeter Brune PetscFunctionReturn(0); 387eed5f15bSPeter Brune } 388eed5f15bSPeter Brune 389eed5f15bSPeter Brune #undef __FUNCT__ 390eed5f15bSPeter Brune #define __FUNCT__ "SNESDestroy_Composite" 391eed5f15bSPeter Brune static PetscErrorCode SNESDestroy_Composite(SNES snes) 392eed5f15bSPeter Brune { 393eed5f15bSPeter Brune SNES_Composite *jac = (SNES_Composite*)snes->data; 394eed5f15bSPeter Brune PetscErrorCode ierr; 395eed5f15bSPeter Brune SNES_CompositeLink next = jac->head,next_tmp; 396eed5f15bSPeter Brune 397eed5f15bSPeter Brune PetscFunctionBegin; 398eed5f15bSPeter Brune ierr = SNESReset_Composite(snes);CHKERRQ(ierr); 399eed5f15bSPeter Brune while (next) { 400eed5f15bSPeter Brune ierr = SNESDestroy(&next->snes);CHKERRQ(ierr); 401eed5f15bSPeter Brune next_tmp = next; 402eed5f15bSPeter Brune next = next->next; 403eed5f15bSPeter Brune ierr = PetscFree(next_tmp);CHKERRQ(ierr); 404eed5f15bSPeter Brune } 405eed5f15bSPeter Brune ierr = PetscFree(snes->data);CHKERRQ(ierr); 406eed5f15bSPeter Brune PetscFunctionReturn(0); 407eed5f15bSPeter Brune } 408eed5f15bSPeter Brune 409eed5f15bSPeter Brune #undef __FUNCT__ 410eed5f15bSPeter Brune #define __FUNCT__ "SNESSetFromOptions_Composite" 411eed5f15bSPeter Brune static PetscErrorCode SNESSetFromOptions_Composite(SNES snes) 412eed5f15bSPeter Brune { 413eed5f15bSPeter Brune SNES_Composite *jac = (SNES_Composite*)snes->data; 414eed5f15bSPeter Brune PetscErrorCode ierr; 415eed5f15bSPeter Brune PetscInt nmax = 8,i; 416eed5f15bSPeter Brune SNES_CompositeLink next; 417eed5f15bSPeter Brune char *sneses[8]; 4188f626970SPeter Brune PetscReal dmps[8]; 419eed5f15bSPeter Brune PetscBool flg; 420eed5f15bSPeter Brune 421eed5f15bSPeter Brune PetscFunctionBegin; 422eed5f15bSPeter Brune ierr = PetscOptionsHead("Composite preconditioner options");CHKERRQ(ierr); 423eed5f15bSPeter Brune ierr = PetscOptionsEnum("-snes_composite_type","Type of composition","SNESCompositeSetType",SNESCompositeTypes,(PetscEnum)jac->type,(PetscEnum*)&jac->type,&flg);CHKERRQ(ierr); 424eed5f15bSPeter Brune if (flg) { 425eed5f15bSPeter Brune ierr = SNESCompositeSetType(snes,jac->type);CHKERRQ(ierr); 426eed5f15bSPeter Brune } 427eed5f15bSPeter Brune ierr = PetscOptionsStringArray("-snes_composite_sneses","List of composite solvers","SNESCompositeAddSNES",sneses,&nmax,&flg);CHKERRQ(ierr); 428eed5f15bSPeter Brune if (flg) { 429eed5f15bSPeter Brune for (i=0; i<nmax; i++) { 430eed5f15bSPeter Brune ierr = SNESCompositeAddSNES(snes,sneses[i]);CHKERRQ(ierr); 431eed5f15bSPeter Brune ierr = PetscFree(sneses[i]);CHKERRQ(ierr); /* deallocate string sneses[i], which is allocated in PetscOptionsStringArray() */ 432eed5f15bSPeter Brune } 433eed5f15bSPeter Brune } 4348f626970SPeter Brune ierr = PetscOptionsRealArray("-snes_composite_damping","Damping of the additive composite solvers","SNESCompositeSetDamping",dmps,&nmax,&flg);CHKERRQ(ierr); 4358f626970SPeter Brune if (flg) { 4368f626970SPeter Brune for (i=0; i<nmax; i++) { 4378f626970SPeter Brune ierr = SNESCompositeSetDamping(snes,i,dmps[i]);CHKERRQ(ierr); 4388f626970SPeter Brune } 4398f626970SPeter Brune } 4405e806d2eSPeter Brune ierr = PetscOptionsReal("-snes_composite_stol","Step tolerance for restart on the additive composite solvers","",jac->stol,&jac->stol,NULL);CHKERRQ(ierr); 4415e806d2eSPeter Brune ierr = PetscOptionsReal("-snes_composite_rtol","Residual tolerance for the additive composite solvers","",jac->rtol,&jac->rtol,NULL);CHKERRQ(ierr); 442eed5f15bSPeter Brune ierr = PetscOptionsTail();CHKERRQ(ierr); 443eed5f15bSPeter Brune 444eed5f15bSPeter Brune next = jac->head; 445eed5f15bSPeter Brune while (next) { 446eed5f15bSPeter Brune ierr = SNESSetFromOptions(next->snes);CHKERRQ(ierr); 447eed5f15bSPeter Brune next = next->next; 448eed5f15bSPeter Brune } 449eed5f15bSPeter Brune PetscFunctionReturn(0); 450eed5f15bSPeter Brune } 451eed5f15bSPeter Brune 452eed5f15bSPeter Brune #undef __FUNCT__ 453eed5f15bSPeter Brune #define __FUNCT__ "SNESView_Composite" 454eed5f15bSPeter Brune static PetscErrorCode SNESView_Composite(SNES snes,PetscViewer viewer) 455eed5f15bSPeter Brune { 456eed5f15bSPeter Brune SNES_Composite *jac = (SNES_Composite*)snes->data; 457eed5f15bSPeter Brune PetscErrorCode ierr; 458eed5f15bSPeter Brune SNES_CompositeLink next = jac->head; 459eed5f15bSPeter Brune PetscBool iascii; 460eed5f15bSPeter Brune 461eed5f15bSPeter Brune PetscFunctionBegin; 462eed5f15bSPeter Brune ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 463eed5f15bSPeter Brune if (iascii) { 464eed5f15bSPeter Brune ierr = PetscViewerASCIIPrintf(viewer,"Composite SNES type - %s\n",SNESCompositeTypes[jac->type]);CHKERRQ(ierr); 465eed5f15bSPeter Brune ierr = PetscViewerASCIIPrintf(viewer,"SNESes on composite preconditioner follow\n");CHKERRQ(ierr); 466eed5f15bSPeter Brune ierr = PetscViewerASCIIPrintf(viewer,"---------------------------------\n");CHKERRQ(ierr); 467eed5f15bSPeter Brune } 468eed5f15bSPeter Brune if (iascii) { 469eed5f15bSPeter Brune ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 470eed5f15bSPeter Brune } 471eed5f15bSPeter Brune while (next) { 472eed5f15bSPeter Brune ierr = SNESView(next->snes,viewer);CHKERRQ(ierr); 473eed5f15bSPeter Brune next = next->next; 474eed5f15bSPeter Brune } 475eed5f15bSPeter Brune if (iascii) { 476eed5f15bSPeter Brune ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 477eed5f15bSPeter Brune ierr = PetscViewerASCIIPrintf(viewer,"---------------------------------\n");CHKERRQ(ierr); 478eed5f15bSPeter Brune } 479eed5f15bSPeter Brune PetscFunctionReturn(0); 480eed5f15bSPeter Brune } 481eed5f15bSPeter Brune 482eed5f15bSPeter Brune /* ------------------------------------------------------------------------------*/ 483eed5f15bSPeter Brune 484eed5f15bSPeter Brune #undef __FUNCT__ 485eed5f15bSPeter Brune #define __FUNCT__ "SNESCompositeSetType_Composite" 486eed5f15bSPeter Brune static PetscErrorCode SNESCompositeSetType_Composite(SNES snes,SNESCompositeType type) 487eed5f15bSPeter Brune { 488eed5f15bSPeter Brune SNES_Composite *jac = (SNES_Composite*)snes->data; 489eed5f15bSPeter Brune 490eed5f15bSPeter Brune PetscFunctionBegin; 491eed5f15bSPeter Brune jac->type = type; 492eed5f15bSPeter Brune PetscFunctionReturn(0); 493eed5f15bSPeter Brune } 494eed5f15bSPeter Brune 495eed5f15bSPeter Brune #undef __FUNCT__ 496eed5f15bSPeter Brune #define __FUNCT__ "SNESCompositeAddSNES_Composite" 497eed5f15bSPeter Brune static PetscErrorCode SNESCompositeAddSNES_Composite(SNES snes,SNESType type) 498eed5f15bSPeter Brune { 499eed5f15bSPeter Brune SNES_Composite *jac; 500eed5f15bSPeter Brune SNES_CompositeLink next,ilink; 501eed5f15bSPeter Brune PetscErrorCode ierr; 502eed5f15bSPeter Brune PetscInt cnt = 0; 503eed5f15bSPeter Brune const char *prefix; 504eed5f15bSPeter Brune char newprefix[8]; 505eed5f15bSPeter Brune DM dm; 506eed5f15bSPeter Brune 507eed5f15bSPeter Brune PetscFunctionBegin; 508eed5f15bSPeter Brune ierr = PetscNewLog(snes,struct _SNES_CompositeLink,&ilink);CHKERRQ(ierr); 509eed5f15bSPeter Brune ilink->next = 0; 510eed5f15bSPeter Brune ierr = SNESCreate(PetscObjectComm((PetscObject)snes),&ilink->snes);CHKERRQ(ierr); 511eed5f15bSPeter Brune ierr = PetscLogObjectParent((PetscObject)snes,(PetscObject)ilink->snes);CHKERRQ(ierr); 512eed5f15bSPeter Brune ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 513eed5f15bSPeter Brune ierr = SNESSetDM(ilink->snes,dm);CHKERRQ(ierr); 514cf5b3eb5SPeter Brune ierr = SNESSetTolerances(ilink->snes,ilink->snes->abstol,ilink->snes->rtol,ilink->snes->stol,1,ilink->snes->max_funcs);CHKERRQ(ierr); 515eed5f15bSPeter Brune jac = (SNES_Composite*)snes->data; 516eed5f15bSPeter Brune next = jac->head; 517eed5f15bSPeter Brune if (!next) { 518eed5f15bSPeter Brune jac->head = ilink; 519eed5f15bSPeter Brune ilink->previous = NULL; 520eed5f15bSPeter Brune } else { 521eed5f15bSPeter Brune cnt++; 522eed5f15bSPeter Brune while (next->next) { 523eed5f15bSPeter Brune next = next->next; 524eed5f15bSPeter Brune cnt++; 525eed5f15bSPeter Brune } 526eed5f15bSPeter Brune next->next = ilink; 527eed5f15bSPeter Brune ilink->previous = next; 528eed5f15bSPeter Brune } 529eed5f15bSPeter Brune ierr = SNESGetOptionsPrefix(snes,&prefix);CHKERRQ(ierr); 530eed5f15bSPeter Brune ierr = SNESSetOptionsPrefix(ilink->snes,prefix);CHKERRQ(ierr); 531eed5f15bSPeter Brune sprintf(newprefix,"sub_%d_",(int)cnt); 532eed5f15bSPeter Brune ierr = SNESAppendOptionsPrefix(ilink->snes,newprefix);CHKERRQ(ierr); 5338f626970SPeter Brune ierr = PetscObjectIncrementTabLevel((PetscObject)ilink->snes,(PetscObject)snes,1);CHKERRQ(ierr); 534eed5f15bSPeter Brune ierr = SNESSetType(ilink->snes,type);CHKERRQ(ierr); 53572edecb9SPeter Brune ierr = SNESSetNormSchedule(ilink->snes, SNES_NORM_FINAL_ONLY);CHKERRQ(ierr); 5368f626970SPeter Brune ilink->dmp = 1.0; 53790a8ba9bSPeter Brune jac->nsnes++; 538eed5f15bSPeter Brune PetscFunctionReturn(0); 539eed5f15bSPeter Brune } 540eed5f15bSPeter Brune 541eed5f15bSPeter Brune #undef __FUNCT__ 542eed5f15bSPeter Brune #define __FUNCT__ "SNESCompositeGetSNES_Composite" 543eed5f15bSPeter Brune static PetscErrorCode SNESCompositeGetSNES_Composite(SNES snes,PetscInt n,SNES *subsnes) 544eed5f15bSPeter Brune { 545eed5f15bSPeter Brune SNES_Composite *jac; 546eed5f15bSPeter Brune SNES_CompositeLink next; 547eed5f15bSPeter Brune PetscInt i; 548eed5f15bSPeter Brune 549eed5f15bSPeter Brune PetscFunctionBegin; 550eed5f15bSPeter Brune jac = (SNES_Composite*)snes->data; 551eed5f15bSPeter Brune next = jac->head; 552eed5f15bSPeter Brune for (i=0; i<n; i++) { 553eed5f15bSPeter Brune if (!next->next) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_INCOMP,"Not enough SNESes in composite preconditioner"); 554eed5f15bSPeter Brune next = next->next; 555eed5f15bSPeter Brune } 556eed5f15bSPeter Brune *subsnes = next->snes; 557eed5f15bSPeter Brune PetscFunctionReturn(0); 558eed5f15bSPeter Brune } 559eed5f15bSPeter Brune 560eed5f15bSPeter Brune /* -------------------------------------------------------------------------------- */ 561eed5f15bSPeter Brune #undef __FUNCT__ 562eed5f15bSPeter Brune #define __FUNCT__ "SNESCompositeSetType" 563e31ab4f9SPeter Brune /*@C 564eed5f15bSPeter Brune SNESCompositeSetType - Sets the type of composite preconditioner. 565eed5f15bSPeter Brune 566eed5f15bSPeter Brune Logically Collective on SNES 567eed5f15bSPeter Brune 568eed5f15bSPeter Brune Input Parameter: 569eed5f15bSPeter Brune + snes - the preconditioner context 570eed5f15bSPeter Brune - type - SNES_COMPOSITE_ADDITIVE (default), SNES_COMPOSITE_MULTIPLICATIVE 571eed5f15bSPeter Brune 572eed5f15bSPeter Brune Options Database Key: 573eed5f15bSPeter Brune . -snes_composite_type <type: one of multiplicative, additive, special> - Sets composite preconditioner type 574eed5f15bSPeter Brune 575eed5f15bSPeter Brune Level: Developer 576eed5f15bSPeter Brune 577eed5f15bSPeter Brune .keywords: SNES, set, type, composite preconditioner, additive, multiplicative 578eed5f15bSPeter Brune @*/ 579eed5f15bSPeter Brune PetscErrorCode SNESCompositeSetType(SNES snes,SNESCompositeType type) 580eed5f15bSPeter Brune { 581eed5f15bSPeter Brune PetscErrorCode ierr; 582eed5f15bSPeter Brune 583eed5f15bSPeter Brune PetscFunctionBegin; 584eed5f15bSPeter Brune PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 585eed5f15bSPeter Brune PetscValidLogicalCollectiveEnum(snes,type,2); 586eed5f15bSPeter Brune ierr = PetscTryMethod(snes,"SNESCompositeSetType_C",(SNES,SNESCompositeType),(snes,type));CHKERRQ(ierr); 587eed5f15bSPeter Brune PetscFunctionReturn(0); 588eed5f15bSPeter Brune } 589eed5f15bSPeter Brune 590eed5f15bSPeter Brune #undef __FUNCT__ 591eed5f15bSPeter Brune #define __FUNCT__ "SNESCompositeAddSNES" 592eed5f15bSPeter Brune /*@C 593eed5f15bSPeter Brune SNESCompositeAddSNES - Adds another SNES to the composite SNES. 594eed5f15bSPeter Brune 595eed5f15bSPeter Brune Collective on SNES 596eed5f15bSPeter Brune 597eed5f15bSPeter Brune Input Parameters: 598eed5f15bSPeter Brune + snes - the preconditioner context 599eed5f15bSPeter Brune - type - the type of the new preconditioner 600eed5f15bSPeter Brune 601eed5f15bSPeter Brune Level: Developer 602eed5f15bSPeter Brune 603eed5f15bSPeter Brune .keywords: SNES, composite preconditioner, add 604eed5f15bSPeter Brune @*/ 605eed5f15bSPeter Brune PetscErrorCode SNESCompositeAddSNES(SNES snes,SNESType type) 606eed5f15bSPeter Brune { 607eed5f15bSPeter Brune PetscErrorCode ierr; 608eed5f15bSPeter Brune 609eed5f15bSPeter Brune PetscFunctionBegin; 610eed5f15bSPeter Brune PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 611eed5f15bSPeter Brune ierr = PetscTryMethod(snes,"SNESCompositeAddSNES_C",(SNES,SNESType),(snes,type));CHKERRQ(ierr); 612eed5f15bSPeter Brune PetscFunctionReturn(0); 613eed5f15bSPeter Brune } 614eed5f15bSPeter Brune #undef __FUNCT__ 615eed5f15bSPeter Brune #define __FUNCT__ "SNESCompositeGetSNES" 616eed5f15bSPeter Brune /*@ 617eed5f15bSPeter Brune SNESCompositeGetSNES - Gets one of the SNES objects in the composite SNES. 618eed5f15bSPeter Brune 619eed5f15bSPeter Brune Not Collective 620eed5f15bSPeter Brune 621eed5f15bSPeter Brune Input Parameter: 622eed5f15bSPeter Brune + snes - the preconditioner context 623eed5f15bSPeter Brune - n - the number of the snes requested 624eed5f15bSPeter Brune 625eed5f15bSPeter Brune Output Parameters: 626eed5f15bSPeter Brune . subsnes - the SNES requested 627eed5f15bSPeter Brune 628eed5f15bSPeter Brune Level: Developer 629eed5f15bSPeter Brune 630eed5f15bSPeter Brune .keywords: SNES, get, composite preconditioner, sub preconditioner 631eed5f15bSPeter Brune 632eed5f15bSPeter Brune .seealso: SNESCompositeAddSNES() 633eed5f15bSPeter Brune @*/ 634eed5f15bSPeter Brune PetscErrorCode SNESCompositeGetSNES(SNES snes,PetscInt n,SNES *subsnes) 635eed5f15bSPeter Brune { 636eed5f15bSPeter Brune PetscErrorCode ierr; 637eed5f15bSPeter Brune 638eed5f15bSPeter Brune PetscFunctionBegin; 639eed5f15bSPeter Brune PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 640eed5f15bSPeter Brune PetscValidPointer(subsnes,3); 641eed5f15bSPeter Brune ierr = PetscUseMethod(snes,"SNESCompositeGetSNES_C",(SNES,PetscInt,SNES*),(snes,n,subsnes));CHKERRQ(ierr); 642eed5f15bSPeter Brune PetscFunctionReturn(0); 643eed5f15bSPeter Brune } 644eed5f15bSPeter Brune 645eed5f15bSPeter Brune #undef __FUNCT__ 6468f626970SPeter Brune #define __FUNCT__ "SNESCompositeSetDamping_Composite" 6478f626970SPeter Brune static PetscErrorCode SNESCompositeSetDamping_Composite(SNES snes,PetscInt n,PetscReal dmp) 6488f626970SPeter Brune { 6498f626970SPeter Brune SNES_Composite *jac; 6508f626970SPeter Brune SNES_CompositeLink next; 6518f626970SPeter Brune PetscInt i; 6528f626970SPeter Brune 6538f626970SPeter Brune PetscFunctionBegin; 6548f626970SPeter Brune jac = (SNES_Composite*)snes->data; 6558f626970SPeter Brune next = jac->head; 6568f626970SPeter Brune for (i=0; i<n; i++) { 6578f626970SPeter Brune if (!next->next) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_INCOMP,"Not enough SNESes in composite preconditioner"); 6588f626970SPeter Brune next = next->next; 6598f626970SPeter Brune } 6608f626970SPeter Brune next->dmp = dmp; 6618f626970SPeter Brune PetscFunctionReturn(0); 6628f626970SPeter Brune } 6638f626970SPeter Brune 6648f626970SPeter Brune #undef __FUNCT__ 6658f626970SPeter Brune #define __FUNCT__ "SNESCompositeSetDamping" 6668f626970SPeter Brune /*@ 6678f626970SPeter Brune SNESCompositeSetDamping - Sets the damping of a subsolver when using additive composite SNES. 6688f626970SPeter Brune 6698f626970SPeter Brune Not Collective 6708f626970SPeter Brune 6718f626970SPeter Brune Input Parameter: 6728f626970SPeter Brune + snes - the preconditioner context 6738f626970SPeter Brune . n - the number of the snes requested 6748f626970SPeter Brune - dmp - the damping 6758f626970SPeter Brune 6768f626970SPeter Brune Level: Developer 6778f626970SPeter Brune 6788f626970SPeter Brune .keywords: SNES, get, composite preconditioner, sub preconditioner 6798f626970SPeter Brune 6808f626970SPeter Brune .seealso: SNESCompositeAddSNES() 6818f626970SPeter Brune @*/ 6828f626970SPeter Brune PetscErrorCode SNESCompositeSetDamping(SNES snes,PetscInt n,PetscReal dmp) 6838f626970SPeter Brune { 6848f626970SPeter Brune PetscErrorCode ierr; 6858f626970SPeter Brune 6868f626970SPeter Brune PetscFunctionBegin; 6878f626970SPeter Brune PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 6888f626970SPeter Brune ierr = PetscUseMethod(snes,"SNESCompositeSetDamping_C",(SNES,PetscInt,PetscReal),(snes,n,dmp));CHKERRQ(ierr); 6898f626970SPeter Brune PetscFunctionReturn(0); 6908f626970SPeter Brune } 6918f626970SPeter Brune 6928f626970SPeter Brune #undef __FUNCT__ 693eed5f15bSPeter Brune #define __FUNCT__ "SNESSolve_Composite" 694eed5f15bSPeter Brune PetscErrorCode SNESSolve_Composite(SNES snes) 695eed5f15bSPeter Brune { 696eed5f15bSPeter Brune Vec F; 697eed5f15bSPeter Brune Vec X; 698eed5f15bSPeter Brune Vec B; 699eed5f15bSPeter Brune PetscInt i; 700eed5f15bSPeter Brune PetscReal fnorm = 0.0; 701eed5f15bSPeter Brune PetscErrorCode ierr; 70272edecb9SPeter Brune SNESNormSchedule normtype; 703eed5f15bSPeter Brune SNES_Composite *comp = (SNES_Composite*)snes->data; 704eed5f15bSPeter Brune 705eed5f15bSPeter Brune PetscFunctionBegin; 706eed5f15bSPeter Brune X = snes->vec_sol; 707eed5f15bSPeter Brune F = snes->vec_func; 708eed5f15bSPeter Brune B = snes->vec_rhs; 709eed5f15bSPeter Brune 710eed5f15bSPeter Brune ierr = PetscObjectAMSTakeAccess((PetscObject)snes);CHKERRQ(ierr); 711eed5f15bSPeter Brune snes->iter = 0; 712eed5f15bSPeter Brune snes->norm = 0.; 713eed5f15bSPeter Brune ierr = PetscObjectAMSGrantAccess((PetscObject)snes);CHKERRQ(ierr); 714eed5f15bSPeter Brune snes->reason = SNES_CONVERGED_ITERATING; 71572edecb9SPeter Brune ierr = SNESGetNormSchedule(snes, &normtype);CHKERRQ(ierr); 71672edecb9SPeter Brune if (normtype == SNES_NORM_ALWAYS || normtype == SNES_NORM_INITIAL_ONLY || normtype == SNES_NORM_INITIAL_FINAL_ONLY) { 717eed5f15bSPeter Brune if (!snes->vec_func_init_set) { 718eed5f15bSPeter Brune ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 719eed5f15bSPeter Brune if (snes->domainerror) { 720eed5f15bSPeter Brune snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 721eed5f15bSPeter Brune PetscFunctionReturn(0); 722eed5f15bSPeter Brune } 723eed5f15bSPeter Brune } else snes->vec_func_init_set = PETSC_FALSE; 724eed5f15bSPeter Brune 725eed5f15bSPeter Brune /* convergence test */ 726eed5f15bSPeter Brune if (!snes->norm_init_set) { 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 PetscFunctionReturn(0); 731eed5f15bSPeter Brune } 732eed5f15bSPeter Brune } else { 733eed5f15bSPeter Brune fnorm = snes->norm_init; 734eed5f15bSPeter Brune snes->norm_init_set = PETSC_FALSE; 735eed5f15bSPeter Brune } 736eed5f15bSPeter Brune ierr = PetscObjectAMSTakeAccess((PetscObject)snes);CHKERRQ(ierr); 737eed5f15bSPeter Brune snes->iter = 0; 738eed5f15bSPeter Brune snes->norm = fnorm; 739eed5f15bSPeter Brune ierr = PetscObjectAMSGrantAccess((PetscObject)snes);CHKERRQ(ierr); 740eed5f15bSPeter Brune ierr = SNESLogConvergenceHistory(snes,snes->norm,0);CHKERRQ(ierr); 741eed5f15bSPeter Brune ierr = SNESMonitor(snes,0,snes->norm);CHKERRQ(ierr); 742eed5f15bSPeter Brune 743eed5f15bSPeter Brune /* set parameter for default relative tolerance convergence test */ 744eed5f15bSPeter Brune snes->ttol = fnorm*snes->rtol; 745eed5f15bSPeter Brune 746eed5f15bSPeter Brune /* test convergence */ 747eed5f15bSPeter Brune ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 748eed5f15bSPeter Brune if (snes->reason) PetscFunctionReturn(0); 749eed5f15bSPeter Brune } else { 750eed5f15bSPeter Brune ierr = PetscObjectAMSGrantAccess((PetscObject)snes);CHKERRQ(ierr); 751eed5f15bSPeter Brune ierr = SNESLogConvergenceHistory(snes,snes->norm,0);CHKERRQ(ierr); 752eed5f15bSPeter Brune ierr = SNESMonitor(snes,0,snes->norm);CHKERRQ(ierr); 753eed5f15bSPeter Brune } 754eed5f15bSPeter Brune 755eed5f15bSPeter Brune /* Call general purpose update function */ 756eed5f15bSPeter Brune if (snes->ops->update) { 757eed5f15bSPeter Brune ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 758eed5f15bSPeter Brune } 759eed5f15bSPeter Brune 760eed5f15bSPeter Brune for (i = 0; i < snes->max_its; i++) { 761eed5f15bSPeter Brune if (comp->type == SNES_COMPOSITE_ADDITIVE) { 762eed5f15bSPeter Brune ierr = SNESCompositeApply_Additive(snes,X,B,F,&fnorm);CHKERRQ(ierr); 763eed5f15bSPeter Brune } else if (comp->type == SNES_COMPOSITE_MULTIPLICATIVE) { 764eed5f15bSPeter Brune ierr = SNESCompositeApply_Multiplicative(snes,X,B,F,&fnorm);CHKERRQ(ierr); 76590a8ba9bSPeter Brune } else if (comp->type == SNES_COMPOSITE_ADDITIVEOPTIMAL) { 76690a8ba9bSPeter Brune ierr = SNESCompositeApply_AdditiveOptimal(snes,X,B,F,&fnorm);CHKERRQ(ierr); 767eed5f15bSPeter Brune } else { 76890a8ba9bSPeter Brune SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"Unsupported SNESComposite type"); 769eed5f15bSPeter Brune } 770*d5a53f06SPeter Brune if (snes->reason < 0) break; 771*d5a53f06SPeter Brune 772eed5f15bSPeter Brune if ((i == snes->max_its - 1) && (normtype == SNES_NORM_INITIAL_FINAL_ONLY || normtype == SNES_NORM_FINAL_ONLY)) { 773eed5f15bSPeter Brune ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 774eed5f15bSPeter Brune if (snes->domainerror) { 775eed5f15bSPeter Brune snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 776eed5f15bSPeter Brune break; 777eed5f15bSPeter Brune } 778*d5a53f06SPeter Brune ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); 779eed5f15bSPeter Brune if (PetscIsInfOrNanReal(fnorm)) { 780eed5f15bSPeter Brune snes->reason = SNES_DIVERGED_FNORM_NAN; 781eed5f15bSPeter Brune break; 782eed5f15bSPeter Brune } 783eed5f15bSPeter Brune } 784eed5f15bSPeter Brune /* Monitor convergence */ 785eed5f15bSPeter Brune ierr = PetscObjectAMSTakeAccess((PetscObject)snes);CHKERRQ(ierr); 786eed5f15bSPeter Brune snes->iter = i+1; 787eed5f15bSPeter Brune snes->norm = fnorm; 788eed5f15bSPeter Brune ierr = PetscObjectAMSGrantAccess((PetscObject)snes);CHKERRQ(ierr); 789eed5f15bSPeter Brune ierr = SNESLogConvergenceHistory(snes,snes->norm,0);CHKERRQ(ierr); 790eed5f15bSPeter Brune ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr); 791eed5f15bSPeter Brune /* Test for convergence */ 79272edecb9SPeter Brune if (normtype == SNES_NORM_ALWAYS) {ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);} 793eed5f15bSPeter Brune if (snes->reason) break; 794eed5f15bSPeter Brune /* Call general purpose update function */ 795eed5f15bSPeter Brune if (snes->ops->update) {ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);} 796eed5f15bSPeter Brune } 79772edecb9SPeter Brune if (normtype == SNES_NORM_ALWAYS) { 798eed5f15bSPeter Brune if (i == snes->max_its) { 799eed5f15bSPeter Brune ierr = PetscInfo1(snes,"Maximum number of iterations has been reached: %D\n",snes->max_its);CHKERRQ(ierr); 800eed5f15bSPeter Brune if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT; 801eed5f15bSPeter Brune } 802eed5f15bSPeter Brune } else if (!snes->reason) snes->reason = SNES_CONVERGED_ITS; 803eed5f15bSPeter Brune PetscFunctionReturn(0); 804eed5f15bSPeter Brune } 805eed5f15bSPeter Brune 806eed5f15bSPeter Brune /* -------------------------------------------------------------------------------------------*/ 807eed5f15bSPeter Brune 808eed5f15bSPeter Brune /*MC 809eed5f15bSPeter Brune SNESCOMPOSITE - Build a preconditioner by composing together several nonlinear solvers 810eed5f15bSPeter Brune 811eed5f15bSPeter Brune Options Database Keys: 812eed5f15bSPeter Brune + -snes_composite_type <type: one of multiplicative, additive, symmetric_multiplicative, special> - Sets composite preconditioner type 813eed5f15bSPeter Brune - -snes_composite_sneses - <snes0,snes1,...> list of SNESes to compose 814eed5f15bSPeter Brune 815eed5f15bSPeter Brune Level: intermediate 816eed5f15bSPeter Brune 817eed5f15bSPeter Brune Concepts: composing solvers 818eed5f15bSPeter Brune 819eed5f15bSPeter Brune .seealso: SNESCreate(), SNESSetType(), SNESType (for list of available types), SNES, 820eed5f15bSPeter Brune SNESSHELL, SNESCompositeSetType(), SNESCompositeSpecialSetAlpha(), SNESCompositeAddSNES(), 821eed5f15bSPeter Brune SNESCompositeGetSNES() 822eed5f15bSPeter Brune 823eed5f15bSPeter Brune M*/ 824eed5f15bSPeter Brune 825eed5f15bSPeter Brune #undef __FUNCT__ 826eed5f15bSPeter Brune #define __FUNCT__ "SNESCreate_Composite" 827eed5f15bSPeter Brune PETSC_EXTERN PetscErrorCode SNESCreate_Composite(SNES snes) 828eed5f15bSPeter Brune { 829eed5f15bSPeter Brune PetscErrorCode ierr; 830eed5f15bSPeter Brune SNES_Composite *jac; 831eed5f15bSPeter Brune 832eed5f15bSPeter Brune PetscFunctionBegin; 833eed5f15bSPeter Brune ierr = PetscNewLog(snes,SNES_Composite,&jac);CHKERRQ(ierr); 834eed5f15bSPeter Brune 835eed5f15bSPeter Brune snes->ops->solve = SNESSolve_Composite; 836eed5f15bSPeter Brune snes->ops->setup = SNESSetUp_Composite; 837eed5f15bSPeter Brune snes->ops->reset = SNESReset_Composite; 838eed5f15bSPeter Brune snes->ops->destroy = SNESDestroy_Composite; 839eed5f15bSPeter Brune snes->ops->setfromoptions = SNESSetFromOptions_Composite; 840eed5f15bSPeter Brune snes->ops->view = SNESView_Composite; 841eed5f15bSPeter Brune 842eed5f15bSPeter Brune snes->data = (void*)jac; 84390a8ba9bSPeter Brune jac->type = SNES_COMPOSITE_ADDITIVEOPTIMAL; 84490a8ba9bSPeter Brune jac->Fes = NULL; 84590a8ba9bSPeter Brune jac->Xes = NULL; 8469c2f3473SPeter Brune jac->fnorms = NULL; 84790a8ba9bSPeter Brune jac->nsnes = 0; 848eed5f15bSPeter Brune jac->head = 0; 8495e806d2eSPeter Brune jac->stol = 0.1; 8505e806d2eSPeter Brune jac->rtol = 1.1; 851eed5f15bSPeter Brune 85290a8ba9bSPeter Brune jac->h = NULL; 85390a8ba9bSPeter Brune jac->s = NULL; 85490a8ba9bSPeter Brune jac->beta = NULL; 85590a8ba9bSPeter Brune jac->work = NULL; 85690a8ba9bSPeter Brune jac->rwork = NULL; 85790a8ba9bSPeter Brune 8588f626970SPeter Brune ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESCompositeSetType_C",SNESCompositeSetType_Composite);CHKERRQ(ierr); 8598f626970SPeter Brune ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESCompositeAddSNES_C",SNESCompositeAddSNES_Composite);CHKERRQ(ierr); 8608f626970SPeter Brune ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESCompositeGetSNES_C",SNESCompositeGetSNES_Composite);CHKERRQ(ierr); 8618f626970SPeter Brune ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESCompositeSetDamping_C",SNESCompositeSetDamping_Composite);CHKERRQ(ierr); 862eed5f15bSPeter Brune PetscFunctionReturn(0); 863eed5f15bSPeter Brune } 864eed5f15bSPeter Brune 865