1eed5f15bSPeter Brune 2eed5f15bSPeter Brune /* 3eed5f15bSPeter Brune Defines a SNES that can consist of a collection of SNESes 4eed5f15bSPeter Brune */ 5af0996ceSBarry Smith #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) { 8963cdc2bdSPatrick Farrell if (snes->xl && snes->xu) { 9063cdc2bdSPatrick Farrell ierr = SNESVIComputeInactiveSetFnorm(snes, F, X, fnorm);CHKERRQ(ierr); 9163cdc2bdSPatrick Farrell } else { 9271dbe336SPeter Brune ierr = VecNorm(F, NORM_2, fnorm);CHKERRQ(ierr); 9363cdc2bdSPatrick Farrell } 94*422a814eSBarry Smith SNESCheckFunctionNorm(snes,*fnorm); 95d5a53f06SPeter Brune } 9672edecb9SPeter Brune } else if (snes->normschedule == SNES_NORM_ALWAYS) { 97be95d8f1SBarry Smith ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 98d5a53f06SPeter Brune if (fnorm) { 99a6da83ecSPatrick Farrell if (snes->xl && snes->xu) { 100a6da83ecSPatrick Farrell ierr = SNESVIComputeInactiveSetFnorm(snes, F, X, fnorm);CHKERRQ(ierr); 101a6da83ecSPatrick Farrell } else { 102d5a53f06SPeter Brune ierr = VecNorm(F, NORM_2, fnorm);CHKERRQ(ierr); 103a6da83ecSPatrick Farrell } 104*422a814eSBarry Smith SNESCheckFunctionNorm(snes,*fnorm); 105d5a53f06SPeter Brune } 106eed5f15bSPeter Brune } 107eed5f15bSPeter Brune PetscFunctionReturn(0); 108eed5f15bSPeter Brune } 109eed5f15bSPeter Brune 110eed5f15bSPeter Brune #undef __FUNCT__ 111eed5f15bSPeter Brune #define __FUNCT__ "SNESCompositeApply_Additive" 112eed5f15bSPeter Brune static PetscErrorCode SNESCompositeApply_Additive(SNES snes,Vec X,Vec B,Vec F,PetscReal *fnorm) 113eed5f15bSPeter Brune { 114eed5f15bSPeter Brune PetscErrorCode ierr; 115eed5f15bSPeter Brune SNES_Composite *jac = (SNES_Composite*)snes->data; 116eed5f15bSPeter Brune SNES_CompositeLink next = jac->head; 117eed5f15bSPeter Brune Vec Y,Xorig; 118d5a53f06SPeter Brune SNESConvergedReason reason; 119eed5f15bSPeter Brune 120eed5f15bSPeter Brune PetscFunctionBegin; 121eed5f15bSPeter Brune Y = snes->vec_sol_update; 122eed5f15bSPeter Brune if (!jac->Xorig) {ierr = VecDuplicate(X,&jac->Xorig);CHKERRQ(ierr);} 123eed5f15bSPeter Brune Xorig = jac->Xorig; 124302440fdSBarry Smith ierr = VecCopy(X,Xorig);CHKERRQ(ierr); 125eed5f15bSPeter Brune if (!next) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"No composite SNESes supplied via SNESCompositeAddSNES() or -snes_composite_sneses"); 12672edecb9SPeter Brune if (snes->normschedule == SNES_NORM_ALWAYS) { 127eed5f15bSPeter Brune ierr = SNESSetInitialFunction(next->snes,F);CHKERRQ(ierr); 128eed5f15bSPeter Brune while (next->next) { 129eed5f15bSPeter Brune next = next->next; 130eed5f15bSPeter Brune ierr = SNESSetInitialFunction(next->snes,F);CHKERRQ(ierr); 131eed5f15bSPeter Brune } 132eed5f15bSPeter Brune } 133eed5f15bSPeter Brune next = jac->head; 1348f626970SPeter Brune ierr = VecCopy(Xorig,Y);CHKERRQ(ierr); 135eed5f15bSPeter Brune ierr = SNESSolve(next->snes,B,Y);CHKERRQ(ierr); 136d5a53f06SPeter Brune ierr = SNESGetConvergedReason(next->snes,&reason);CHKERRQ(ierr); 137d5a53f06SPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 138d5a53f06SPeter Brune snes->reason = SNES_DIVERGED_INNER; 139d5a53f06SPeter Brune PetscFunctionReturn(0); 140d5a53f06SPeter Brune } 141eed5f15bSPeter Brune ierr = VecAXPY(Y,-1.0,Xorig);CHKERRQ(ierr); 1428f626970SPeter Brune ierr = VecAXPY(X,next->dmp,Y);CHKERRQ(ierr); 1438f626970SPeter Brune while (next->next) { 1448f626970SPeter Brune next = next->next; 1458f626970SPeter Brune ierr = VecCopy(Xorig,Y);CHKERRQ(ierr); 1468f626970SPeter Brune ierr = SNESSolve(next->snes,B,Y);CHKERRQ(ierr); 147d5a53f06SPeter Brune ierr = SNESGetConvergedReason(next->snes,&reason);CHKERRQ(ierr); 148d5a53f06SPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 149d5a53f06SPeter Brune snes->reason = SNES_DIVERGED_INNER; 150d5a53f06SPeter Brune PetscFunctionReturn(0); 151d5a53f06SPeter Brune } 1528f626970SPeter Brune ierr = VecAXPY(Y,-1.0,Xorig);CHKERRQ(ierr); 1538f626970SPeter Brune ierr = VecAXPY(X,next->dmp,Y);CHKERRQ(ierr); 154eed5f15bSPeter Brune } 15572edecb9SPeter Brune if (snes->normschedule == SNES_NORM_ALWAYS) { 156eed5f15bSPeter Brune ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 157a6da83ecSPatrick Farrell if (fnorm) { 158a6da83ecSPatrick Farrell if (snes->xl && snes->xu) { 159a6da83ecSPatrick Farrell ierr = SNESVIComputeInactiveSetFnorm(snes, F, X, fnorm);CHKERRQ(ierr); 160a6da83ecSPatrick Farrell } else { 161a6da83ecSPatrick Farrell ierr = VecNorm(F, NORM_2, fnorm);CHKERRQ(ierr); 162a6da83ecSPatrick Farrell } 163*422a814eSBarry Smith SNESCheckFunctionNorm(snes,*fnorm); 164a6da83ecSPatrick Farrell } 165eed5f15bSPeter Brune } 166eed5f15bSPeter Brune PetscFunctionReturn(0); 167eed5f15bSPeter Brune } 16890a8ba9bSPeter Brune 16990a8ba9bSPeter Brune #undef __FUNCT__ 17090a8ba9bSPeter Brune #define __FUNCT__ "SNESCompositeApply_AdditiveOptimal" 17190a8ba9bSPeter Brune /* approximately solve the overdetermined system: 17290a8ba9bSPeter Brune 17390a8ba9bSPeter Brune 2*F(x_i)\cdot F(\x_j)\alpha_i = 0 17490a8ba9bSPeter Brune \alpha_i = 1 17590a8ba9bSPeter Brune 17690a8ba9bSPeter Brune Which minimizes the L2 norm of the linearization of: 17790a8ba9bSPeter Brune ||F(\sum_i \alpha_i*x_i)||^2 17890a8ba9bSPeter Brune 17990a8ba9bSPeter Brune With the constraint that \sum_i\alpha_i = 1 18090a8ba9bSPeter Brune Where x_i is the solution from the ith subsolver. 18190a8ba9bSPeter Brune */ 18290a8ba9bSPeter Brune static PetscErrorCode SNESCompositeApply_AdditiveOptimal(SNES snes,Vec X,Vec B,Vec F,PetscReal *fnorm) 18390a8ba9bSPeter Brune { 18490a8ba9bSPeter Brune PetscErrorCode ierr; 18590a8ba9bSPeter Brune SNES_Composite *jac = (SNES_Composite*)snes->data; 18690a8ba9bSPeter Brune SNES_CompositeLink next = jac->head; 18790a8ba9bSPeter Brune Vec *Xes = jac->Xes,*Fes = jac->Fes; 18890a8ba9bSPeter Brune PetscInt i,j; 1895e806d2eSPeter Brune PetscScalar tot,total,ftf; 1905e806d2eSPeter Brune PetscReal min_fnorm; 1915e806d2eSPeter Brune PetscInt min_i; 192d5a53f06SPeter Brune SNESConvergedReason reason; 19390a8ba9bSPeter Brune 19490a8ba9bSPeter Brune PetscFunctionBegin; 19590a8ba9bSPeter Brune if (!next) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"No composite SNESes supplied via SNESCompositeAddSNES() or -snes_composite_sneses"); 19658bdfa14SPeter Brune 19758bdfa14SPeter Brune if (snes->normschedule == SNES_NORM_ALWAYS) { 19858bdfa14SPeter Brune next = jac->head; 19958bdfa14SPeter Brune ierr = SNESSetInitialFunction(next->snes,F);CHKERRQ(ierr); 20058bdfa14SPeter Brune while (next->next) { 20158bdfa14SPeter Brune next = next->next; 20258bdfa14SPeter Brune ierr = SNESSetInitialFunction(next->snes,F);CHKERRQ(ierr); 20358bdfa14SPeter Brune } 20458bdfa14SPeter Brune } 20558bdfa14SPeter Brune 20690a8ba9bSPeter Brune next = jac->head; 20790a8ba9bSPeter Brune i = 0; 20890a8ba9bSPeter Brune ierr = VecCopy(X,Xes[i]);CHKERRQ(ierr); 20990a8ba9bSPeter Brune ierr = SNESSolve(next->snes,B,Xes[i]);CHKERRQ(ierr); 210d5a53f06SPeter Brune ierr = SNESGetConvergedReason(next->snes,&reason);CHKERRQ(ierr); 211d5a53f06SPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 212d5a53f06SPeter Brune snes->reason = SNES_DIVERGED_INNER; 213d5a53f06SPeter Brune PetscFunctionReturn(0); 214d5a53f06SPeter Brune } 21590a8ba9bSPeter Brune while (next->next) { 21690a8ba9bSPeter Brune i++; 21790a8ba9bSPeter Brune next = next->next; 21890a8ba9bSPeter Brune ierr = VecCopy(X,Xes[i]);CHKERRQ(ierr); 21990a8ba9bSPeter Brune ierr = SNESSolve(next->snes,B,Xes[i]);CHKERRQ(ierr); 220d5a53f06SPeter Brune ierr = SNESGetConvergedReason(next->snes,&reason);CHKERRQ(ierr); 221d5a53f06SPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 222d5a53f06SPeter Brune snes->reason = SNES_DIVERGED_INNER; 223d5a53f06SPeter Brune PetscFunctionReturn(0); 224d5a53f06SPeter Brune } 22590a8ba9bSPeter Brune } 22690a8ba9bSPeter Brune 22790a8ba9bSPeter Brune /* all the solutions are collected; combine optimally */ 22890a8ba9bSPeter Brune for (i=0;i<jac->n;i++) { 22990a8ba9bSPeter Brune for (j=0;j<i+1;j++) { 2309c2f3473SPeter Brune ierr = VecDotBegin(Fes[i],Fes[j],&jac->h[i + j*jac->n]);CHKERRQ(ierr); 23190a8ba9bSPeter Brune } 2329c2f3473SPeter Brune ierr = VecDotBegin(Fes[i],F,&jac->g[i]);CHKERRQ(ierr); 23390a8ba9bSPeter Brune } 2345e806d2eSPeter Brune 23590a8ba9bSPeter Brune for (i=0;i<jac->n;i++) { 23690a8ba9bSPeter Brune for (j=0;j<i+1;j++) { 2379c2f3473SPeter Brune ierr = VecDotEnd(Fes[i],Fes[j],&jac->h[i + j*jac->n]);CHKERRQ(ierr); 2389c2f3473SPeter Brune if (i == j) jac->fnorms[i] = PetscSqrtReal(PetscRealPart(jac->h[i + j*jac->n])); 23990a8ba9bSPeter Brune } 2409c2f3473SPeter Brune ierr = VecDotEnd(Fes[i],F,&jac->g[i]);CHKERRQ(ierr); 24190a8ba9bSPeter Brune } 24290a8ba9bSPeter Brune 2439c2f3473SPeter Brune ftf = (*fnorm)*(*fnorm); 24490a8ba9bSPeter Brune 24590a8ba9bSPeter Brune for (i=0; i<jac->n; i++) { 24690a8ba9bSPeter Brune for (j=i+1;j<jac->n;j++) { 2479c2f3473SPeter Brune jac->h[i + j*jac->n] = jac->h[j + i*jac->n]; 24890a8ba9bSPeter Brune } 24990a8ba9bSPeter Brune } 25090a8ba9bSPeter Brune 25190a8ba9bSPeter Brune for (i=0; i<jac->n; i++) { 2529c2f3473SPeter Brune for (j=0; j<jac->n; j++) { 2539c2f3473SPeter Brune jac->h[i + j*jac->n] = jac->h[i + j*jac->n] - jac->g[j] - jac->g[i] + ftf; 25490a8ba9bSPeter Brune } 2559c2f3473SPeter Brune jac->beta[i] = ftf - jac->g[i]; 2569c2f3473SPeter Brune } 25790a8ba9bSPeter Brune 25890a8ba9bSPeter Brune #if defined(PETSC_MISSING_LAPACK_GELSS) 25990a8ba9bSPeter Brune SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_SUP,"SNESCOMPOSITE with ADDITIVEOPTIMAL requires the LAPACK GELSS routine."); 26090a8ba9bSPeter Brune #else 26190a8ba9bSPeter Brune jac->info = 0; 26290a8ba9bSPeter Brune jac->rcond = -1.; 26390a8ba9bSPeter Brune ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr); 26490a8ba9bSPeter Brune #if defined(PETSC_USE_COMPLEX) 2659c2f3473SPeter 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)); 26690a8ba9bSPeter Brune #else 2679c2f3473SPeter 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)); 26890a8ba9bSPeter Brune #endif 26990a8ba9bSPeter Brune ierr = PetscFPTrapPop();CHKERRQ(ierr); 27090a8ba9bSPeter Brune if (jac->info < 0) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_LIB,"Bad argument to GELSS"); 27190a8ba9bSPeter Brune if (jac->info > 0) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_LIB,"SVD failed to converge"); 27290a8ba9bSPeter Brune #endif 2739c2f3473SPeter Brune tot = 0.; 2745e806d2eSPeter Brune total = 0.; 27590a8ba9bSPeter Brune for (i=0; i<jac->n; i++) { 276*422a814eSBarry Smith if (snes->errorifnotconverged && PetscIsInfOrNanScalar(jac->beta[i])) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_LIB,"SVD generated inconsistent output"); 277c783eb9eSBarry Smith ierr = PetscInfo2(snes,"%D: %g\n",i,(double)PetscRealPart(jac->beta[i]));CHKERRQ(ierr); 2789c2f3473SPeter Brune tot += jac->beta[i]; 2795e806d2eSPeter Brune total += PetscAbsScalar(jac->beta[i]); 28090a8ba9bSPeter Brune } 2819c2f3473SPeter Brune ierr = VecScale(X,(1. - tot));CHKERRQ(ierr); 28290a8ba9bSPeter Brune ierr = VecMAXPY(X,jac->n,jac->beta,Xes);CHKERRQ(ierr); 28390a8ba9bSPeter Brune ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 284a6da83ecSPatrick Farrell 285a6da83ecSPatrick Farrell if (snes->xl && snes->xu) { 286a6da83ecSPatrick Farrell ierr = SNESVIComputeInactiveSetFnorm(snes, F, X, fnorm);CHKERRQ(ierr); 287a6da83ecSPatrick Farrell } else { 2889c2f3473SPeter Brune ierr = VecNorm(F, NORM_2, fnorm);CHKERRQ(ierr); 289a6da83ecSPatrick Farrell } 29090a8ba9bSPeter Brune 2915e806d2eSPeter Brune /* take the minimum-normed candidate if it beats the combination by a factor of rtol or the combination has stagnated */ 2925e806d2eSPeter Brune min_fnorm = jac->fnorms[0]; 2935e806d2eSPeter Brune min_i = 0; 2949c2f3473SPeter Brune for (i=0; i<jac->n; i++) { 2955e806d2eSPeter Brune if (jac->fnorms[i] < min_fnorm) { 2965e806d2eSPeter Brune min_fnorm = jac->fnorms[i]; 2975e806d2eSPeter Brune min_i = i; 2989c2f3473SPeter Brune } 2999c2f3473SPeter Brune } 3005e806d2eSPeter Brune 3015e806d2eSPeter Brune /* stagnation or divergence restart to the solution of the solver that failed the least */ 3025e806d2eSPeter Brune if (PetscRealPart(total) < jac->stol || min_fnorm*jac->rtol < *fnorm) { 30358bdfa14SPeter Brune ierr = VecCopy(jac->Xes[min_i],X);CHKERRQ(ierr); 30458bdfa14SPeter Brune ierr = VecCopy(jac->Fes[min_i],F);CHKERRQ(ierr); 3055e806d2eSPeter Brune *fnorm = min_fnorm; 3065e806d2eSPeter Brune } 30790a8ba9bSPeter Brune PetscFunctionReturn(0); 30890a8ba9bSPeter Brune } 30990a8ba9bSPeter Brune 310eed5f15bSPeter Brune #undef __FUNCT__ 311eed5f15bSPeter Brune #define __FUNCT__ "SNESSetUp_Composite" 312eed5f15bSPeter Brune static PetscErrorCode SNESSetUp_Composite(SNES snes) 313eed5f15bSPeter Brune { 314eed5f15bSPeter Brune PetscErrorCode ierr; 315dd515a93SPeter Brune DM dm; 316eed5f15bSPeter Brune SNES_Composite *jac = (SNES_Composite*)snes->data; 317eed5f15bSPeter Brune SNES_CompositeLink next = jac->head; 31890a8ba9bSPeter Brune PetscInt n=0,i; 31990a8ba9bSPeter Brune Vec F; 320eed5f15bSPeter Brune 321eed5f15bSPeter Brune PetscFunctionBegin; 322dd515a93SPeter Brune ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 32363cdc2bdSPatrick Farrell 32463cdc2bdSPatrick Farrell if (snes->ops->computevariablebounds) { 32563cdc2bdSPatrick Farrell /* SNESVI only ever calls computevariablebounds once, so calling it once here is justified */ 32663cdc2bdSPatrick Farrell if (!snes->xl) {ierr = VecDuplicate(snes->vec_sol,&snes->xl);CHKERRQ(ierr);} 32763cdc2bdSPatrick Farrell if (!snes->xu) {ierr = VecDuplicate(snes->vec_sol,&snes->xu);CHKERRQ(ierr);} 32863cdc2bdSPatrick Farrell ierr = (*snes->ops->computevariablebounds)(snes,snes->xl,snes->xu);CHKERRQ(ierr); 32963cdc2bdSPatrick Farrell } 33063cdc2bdSPatrick Farrell 331eed5f15bSPeter Brune while (next) { 33290a8ba9bSPeter Brune n++; 333dd515a93SPeter Brune ierr = SNESSetDM(next->snes,dm);CHKERRQ(ierr); 334eed5f15bSPeter Brune ierr = SNESSetFromOptions(next->snes);CHKERRQ(ierr); 33563cdc2bdSPatrick Farrell ierr = SNESSetApplicationContext(next->snes, snes->user);CHKERRQ(ierr); 33663cdc2bdSPatrick Farrell if (snes->xl && snes->xu) { 33763cdc2bdSPatrick Farrell if (snes->ops->computevariablebounds) { 33863cdc2bdSPatrick Farrell ierr = SNESVISetComputeVariableBounds(next->snes, snes->ops->computevariablebounds);CHKERRQ(ierr); 33963cdc2bdSPatrick Farrell } else { 34063cdc2bdSPatrick Farrell ierr = SNESVISetVariableBounds(next->snes,snes->xl,snes->xu);CHKERRQ(ierr); 34163cdc2bdSPatrick Farrell } 34263cdc2bdSPatrick Farrell } 34363cdc2bdSPatrick Farrell 344eed5f15bSPeter Brune next = next->next; 345eed5f15bSPeter Brune } 34690a8ba9bSPeter Brune jac->nsnes = n; 34790a8ba9bSPeter Brune ierr = SNESGetFunction(snes,&F,NULL,NULL);CHKERRQ(ierr); 34890a8ba9bSPeter Brune if (jac->type == SNES_COMPOSITE_ADDITIVEOPTIMAL) { 34990a8ba9bSPeter Brune ierr = VecDuplicateVecs(F,jac->nsnes,&jac->Xes);CHKERRQ(ierr); 350854ce69bSBarry Smith ierr = PetscMalloc1(n,&jac->Fes);CHKERRQ(ierr); 351854ce69bSBarry Smith ierr = PetscMalloc1(n,&jac->fnorms);CHKERRQ(ierr); 35290a8ba9bSPeter Brune next = jac->head; 35390a8ba9bSPeter Brune i = 0; 35490a8ba9bSPeter Brune while (next) { 35590a8ba9bSPeter Brune ierr = SNESGetFunction(next->snes,&F,NULL,NULL);CHKERRQ(ierr); 35690a8ba9bSPeter Brune jac->Fes[i] = F; 35790a8ba9bSPeter Brune ierr = PetscObjectReference((PetscObject)F);CHKERRQ(ierr); 35890a8ba9bSPeter Brune next = next->next; 35990a8ba9bSPeter Brune i++; 36090a8ba9bSPeter Brune } 36190a8ba9bSPeter Brune /* allocate the subspace direct solve area */ 36290a8ba9bSPeter Brune jac->nrhs = 1; 3639c2f3473SPeter Brune jac->lda = jac->nsnes; 36490a8ba9bSPeter Brune jac->ldb = jac->nsnes; 36590a8ba9bSPeter Brune jac->n = jac->nsnes; 36690a8ba9bSPeter Brune 367785e854fSJed Brown ierr = PetscMalloc1(jac->n*jac->n,&jac->h);CHKERRQ(ierr); 368785e854fSJed Brown ierr = PetscMalloc1(jac->n,&jac->beta);CHKERRQ(ierr); 369785e854fSJed Brown ierr = PetscMalloc1(jac->n,&jac->s);CHKERRQ(ierr); 370785e854fSJed Brown ierr = PetscMalloc1(jac->n,&jac->g);CHKERRQ(ierr); 3719c2f3473SPeter Brune jac->lwork = 12*jac->n; 37290a8ba9bSPeter Brune #if PETSC_USE_COMPLEX 373854ce69bSBarry Smith ierr = PetscMalloc1(jac->lwork,&jac->rwork);CHKERRQ(ierr); 37490a8ba9bSPeter Brune #endif 375854ce69bSBarry Smith ierr = PetscMalloc1(jac->lwork,&jac->work);CHKERRQ(ierr); 37690a8ba9bSPeter Brune } 37790a8ba9bSPeter Brune 378eed5f15bSPeter Brune PetscFunctionReturn(0); 379eed5f15bSPeter Brune } 380eed5f15bSPeter Brune 381eed5f15bSPeter Brune #undef __FUNCT__ 382eed5f15bSPeter Brune #define __FUNCT__ "SNESReset_Composite" 383eed5f15bSPeter Brune static PetscErrorCode SNESReset_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; 388eed5f15bSPeter Brune 389eed5f15bSPeter Brune PetscFunctionBegin; 390eed5f15bSPeter Brune while (next) { 391eed5f15bSPeter Brune ierr = SNESReset(next->snes);CHKERRQ(ierr); 392eed5f15bSPeter Brune next = next->next; 393eed5f15bSPeter Brune } 394eed5f15bSPeter Brune ierr = VecDestroy(&jac->Xorig);CHKERRQ(ierr); 39590a8ba9bSPeter Brune if (jac->Xes) {ierr = VecDestroyVecs(jac->nsnes,&jac->Xes);CHKERRQ(ierr);} 39690a8ba9bSPeter Brune if (jac->Fes) {ierr = VecDestroyVecs(jac->nsnes,&jac->Fes);CHKERRQ(ierr);} 3979c2f3473SPeter Brune ierr = PetscFree(jac->fnorms);CHKERRQ(ierr); 39890a8ba9bSPeter Brune ierr = PetscFree(jac->h);CHKERRQ(ierr); 39990a8ba9bSPeter Brune ierr = PetscFree(jac->s);CHKERRQ(ierr); 4009c2f3473SPeter Brune ierr = PetscFree(jac->g);CHKERRQ(ierr); 40190a8ba9bSPeter Brune ierr = PetscFree(jac->beta);CHKERRQ(ierr); 40290a8ba9bSPeter Brune ierr = PetscFree(jac->work);CHKERRQ(ierr); 40390a8ba9bSPeter Brune ierr = PetscFree(jac->rwork);CHKERRQ(ierr); 404eed5f15bSPeter Brune PetscFunctionReturn(0); 405eed5f15bSPeter Brune } 406eed5f15bSPeter Brune 407eed5f15bSPeter Brune #undef __FUNCT__ 408eed5f15bSPeter Brune #define __FUNCT__ "SNESDestroy_Composite" 409eed5f15bSPeter Brune static PetscErrorCode SNESDestroy_Composite(SNES snes) 410eed5f15bSPeter Brune { 411eed5f15bSPeter Brune SNES_Composite *jac = (SNES_Composite*)snes->data; 412eed5f15bSPeter Brune PetscErrorCode ierr; 413eed5f15bSPeter Brune SNES_CompositeLink next = jac->head,next_tmp; 414eed5f15bSPeter Brune 415eed5f15bSPeter Brune PetscFunctionBegin; 416eed5f15bSPeter Brune ierr = SNESReset_Composite(snes);CHKERRQ(ierr); 417eed5f15bSPeter Brune while (next) { 418eed5f15bSPeter Brune ierr = SNESDestroy(&next->snes);CHKERRQ(ierr); 419eed5f15bSPeter Brune next_tmp = next; 420eed5f15bSPeter Brune next = next->next; 421eed5f15bSPeter Brune ierr = PetscFree(next_tmp);CHKERRQ(ierr); 422eed5f15bSPeter Brune } 423eed5f15bSPeter Brune ierr = PetscFree(snes->data);CHKERRQ(ierr); 424eed5f15bSPeter Brune PetscFunctionReturn(0); 425eed5f15bSPeter Brune } 426eed5f15bSPeter Brune 427eed5f15bSPeter Brune #undef __FUNCT__ 428eed5f15bSPeter Brune #define __FUNCT__ "SNESSetFromOptions_Composite" 4298c34d3f5SBarry Smith static PetscErrorCode SNESSetFromOptions_Composite(PetscOptions *PetscOptionsObject,SNES snes) 430eed5f15bSPeter Brune { 431eed5f15bSPeter Brune SNES_Composite *jac = (SNES_Composite*)snes->data; 432eed5f15bSPeter Brune PetscErrorCode ierr; 433eed5f15bSPeter Brune PetscInt nmax = 8,i; 434eed5f15bSPeter Brune SNES_CompositeLink next; 435eed5f15bSPeter Brune char *sneses[8]; 4368f626970SPeter Brune PetscReal dmps[8]; 437eed5f15bSPeter Brune PetscBool flg; 438eed5f15bSPeter Brune 439eed5f15bSPeter Brune PetscFunctionBegin; 440e55864a3SBarry Smith ierr = PetscOptionsHead(PetscOptionsObject,"Composite preconditioner options");CHKERRQ(ierr); 441eed5f15bSPeter Brune ierr = PetscOptionsEnum("-snes_composite_type","Type of composition","SNESCompositeSetType",SNESCompositeTypes,(PetscEnum)jac->type,(PetscEnum*)&jac->type,&flg);CHKERRQ(ierr); 442eed5f15bSPeter Brune if (flg) { 443eed5f15bSPeter Brune ierr = SNESCompositeSetType(snes,jac->type);CHKERRQ(ierr); 444eed5f15bSPeter Brune } 445eed5f15bSPeter Brune ierr = PetscOptionsStringArray("-snes_composite_sneses","List of composite solvers","SNESCompositeAddSNES",sneses,&nmax,&flg);CHKERRQ(ierr); 446eed5f15bSPeter Brune if (flg) { 447eed5f15bSPeter Brune for (i=0; i<nmax; i++) { 448eed5f15bSPeter Brune ierr = SNESCompositeAddSNES(snes,sneses[i]);CHKERRQ(ierr); 449eed5f15bSPeter Brune ierr = PetscFree(sneses[i]);CHKERRQ(ierr); /* deallocate string sneses[i], which is allocated in PetscOptionsStringArray() */ 450eed5f15bSPeter Brune } 451eed5f15bSPeter Brune } 4528f626970SPeter Brune ierr = PetscOptionsRealArray("-snes_composite_damping","Damping of the additive composite solvers","SNESCompositeSetDamping",dmps,&nmax,&flg);CHKERRQ(ierr); 4538f626970SPeter Brune if (flg) { 4548f626970SPeter Brune for (i=0; i<nmax; i++) { 4558f626970SPeter Brune ierr = SNESCompositeSetDamping(snes,i,dmps[i]);CHKERRQ(ierr); 4568f626970SPeter Brune } 4578f626970SPeter Brune } 4585e806d2eSPeter Brune ierr = PetscOptionsReal("-snes_composite_stol","Step tolerance for restart on the additive composite solvers","",jac->stol,&jac->stol,NULL);CHKERRQ(ierr); 4595e806d2eSPeter Brune ierr = PetscOptionsReal("-snes_composite_rtol","Residual tolerance for the additive composite solvers","",jac->rtol,&jac->rtol,NULL);CHKERRQ(ierr); 460eed5f15bSPeter Brune ierr = PetscOptionsTail();CHKERRQ(ierr); 461eed5f15bSPeter Brune 462eed5f15bSPeter Brune next = jac->head; 463eed5f15bSPeter Brune while (next) { 464eed5f15bSPeter Brune ierr = SNESSetFromOptions(next->snes);CHKERRQ(ierr); 465eed5f15bSPeter Brune next = next->next; 466eed5f15bSPeter Brune } 467eed5f15bSPeter Brune PetscFunctionReturn(0); 468eed5f15bSPeter Brune } 469eed5f15bSPeter Brune 470eed5f15bSPeter Brune #undef __FUNCT__ 471eed5f15bSPeter Brune #define __FUNCT__ "SNESView_Composite" 472eed5f15bSPeter Brune static PetscErrorCode SNESView_Composite(SNES snes,PetscViewer viewer) 473eed5f15bSPeter Brune { 474eed5f15bSPeter Brune SNES_Composite *jac = (SNES_Composite*)snes->data; 475eed5f15bSPeter Brune PetscErrorCode ierr; 476eed5f15bSPeter Brune SNES_CompositeLink next = jac->head; 477eed5f15bSPeter Brune PetscBool iascii; 478eed5f15bSPeter Brune 479eed5f15bSPeter Brune PetscFunctionBegin; 480eed5f15bSPeter Brune ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 481eed5f15bSPeter Brune if (iascii) { 482eed5f15bSPeter Brune ierr = PetscViewerASCIIPrintf(viewer,"Composite SNES type - %s\n",SNESCompositeTypes[jac->type]);CHKERRQ(ierr); 483eed5f15bSPeter Brune ierr = PetscViewerASCIIPrintf(viewer,"SNESes on composite preconditioner follow\n");CHKERRQ(ierr); 484eed5f15bSPeter Brune ierr = PetscViewerASCIIPrintf(viewer,"---------------------------------\n");CHKERRQ(ierr); 485eed5f15bSPeter Brune } 486eed5f15bSPeter Brune if (iascii) { 487eed5f15bSPeter Brune ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 488eed5f15bSPeter Brune } 489eed5f15bSPeter Brune while (next) { 490eed5f15bSPeter Brune ierr = SNESView(next->snes,viewer);CHKERRQ(ierr); 491eed5f15bSPeter Brune next = next->next; 492eed5f15bSPeter Brune } 493eed5f15bSPeter Brune if (iascii) { 494eed5f15bSPeter Brune ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 495eed5f15bSPeter Brune ierr = PetscViewerASCIIPrintf(viewer,"---------------------------------\n");CHKERRQ(ierr); 496eed5f15bSPeter Brune } 497eed5f15bSPeter Brune PetscFunctionReturn(0); 498eed5f15bSPeter Brune } 499eed5f15bSPeter Brune 500eed5f15bSPeter Brune /* ------------------------------------------------------------------------------*/ 501eed5f15bSPeter Brune 502eed5f15bSPeter Brune #undef __FUNCT__ 503eed5f15bSPeter Brune #define __FUNCT__ "SNESCompositeSetType_Composite" 504eed5f15bSPeter Brune static PetscErrorCode SNESCompositeSetType_Composite(SNES snes,SNESCompositeType type) 505eed5f15bSPeter Brune { 506eed5f15bSPeter Brune SNES_Composite *jac = (SNES_Composite*)snes->data; 507eed5f15bSPeter Brune 508eed5f15bSPeter Brune PetscFunctionBegin; 509eed5f15bSPeter Brune jac->type = type; 510eed5f15bSPeter Brune PetscFunctionReturn(0); 511eed5f15bSPeter Brune } 512eed5f15bSPeter Brune 513eed5f15bSPeter Brune #undef __FUNCT__ 514eed5f15bSPeter Brune #define __FUNCT__ "SNESCompositeAddSNES_Composite" 515eed5f15bSPeter Brune static PetscErrorCode SNESCompositeAddSNES_Composite(SNES snes,SNESType type) 516eed5f15bSPeter Brune { 517eed5f15bSPeter Brune SNES_Composite *jac; 518eed5f15bSPeter Brune SNES_CompositeLink next,ilink; 519eed5f15bSPeter Brune PetscErrorCode ierr; 520eed5f15bSPeter Brune PetscInt cnt = 0; 521eed5f15bSPeter Brune const char *prefix; 522eed5f15bSPeter Brune char newprefix[8]; 523eed5f15bSPeter Brune DM dm; 524eed5f15bSPeter Brune 525eed5f15bSPeter Brune PetscFunctionBegin; 526b00a9115SJed Brown ierr = PetscNewLog(snes,&ilink);CHKERRQ(ierr); 527eed5f15bSPeter Brune ilink->next = 0; 528eed5f15bSPeter Brune ierr = SNESCreate(PetscObjectComm((PetscObject)snes),&ilink->snes);CHKERRQ(ierr); 529eed5f15bSPeter Brune ierr = PetscLogObjectParent((PetscObject)snes,(PetscObject)ilink->snes);CHKERRQ(ierr); 530eed5f15bSPeter Brune ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 531eed5f15bSPeter Brune ierr = SNESSetDM(ilink->snes,dm);CHKERRQ(ierr); 532cf5b3eb5SPeter Brune ierr = SNESSetTolerances(ilink->snes,ilink->snes->abstol,ilink->snes->rtol,ilink->snes->stol,1,ilink->snes->max_funcs);CHKERRQ(ierr); 533eed5f15bSPeter Brune jac = (SNES_Composite*)snes->data; 534eed5f15bSPeter Brune next = jac->head; 535eed5f15bSPeter Brune if (!next) { 536eed5f15bSPeter Brune jac->head = ilink; 537eed5f15bSPeter Brune ilink->previous = NULL; 538eed5f15bSPeter Brune } else { 539eed5f15bSPeter Brune cnt++; 540eed5f15bSPeter Brune while (next->next) { 541eed5f15bSPeter Brune next = next->next; 542eed5f15bSPeter Brune cnt++; 543eed5f15bSPeter Brune } 544eed5f15bSPeter Brune next->next = ilink; 545eed5f15bSPeter Brune ilink->previous = next; 546eed5f15bSPeter Brune } 547eed5f15bSPeter Brune ierr = SNESGetOptionsPrefix(snes,&prefix);CHKERRQ(ierr); 548eed5f15bSPeter Brune ierr = SNESSetOptionsPrefix(ilink->snes,prefix);CHKERRQ(ierr); 549eed5f15bSPeter Brune sprintf(newprefix,"sub_%d_",(int)cnt); 550eed5f15bSPeter Brune ierr = SNESAppendOptionsPrefix(ilink->snes,newprefix);CHKERRQ(ierr); 5518f626970SPeter Brune ierr = PetscObjectIncrementTabLevel((PetscObject)ilink->snes,(PetscObject)snes,1);CHKERRQ(ierr); 552eed5f15bSPeter Brune ierr = SNESSetType(ilink->snes,type);CHKERRQ(ierr); 55372edecb9SPeter Brune ierr = SNESSetNormSchedule(ilink->snes, SNES_NORM_FINAL_ONLY);CHKERRQ(ierr); 55463cdc2bdSPatrick Farrell 5558f626970SPeter Brune ilink->dmp = 1.0; 55690a8ba9bSPeter Brune jac->nsnes++; 557eed5f15bSPeter Brune PetscFunctionReturn(0); 558eed5f15bSPeter Brune } 559eed5f15bSPeter Brune 560eed5f15bSPeter Brune #undef __FUNCT__ 561eed5f15bSPeter Brune #define __FUNCT__ "SNESCompositeGetSNES_Composite" 562eed5f15bSPeter Brune static PetscErrorCode SNESCompositeGetSNES_Composite(SNES snes,PetscInt n,SNES *subsnes) 563eed5f15bSPeter Brune { 564eed5f15bSPeter Brune SNES_Composite *jac; 565eed5f15bSPeter Brune SNES_CompositeLink next; 566eed5f15bSPeter Brune PetscInt i; 567eed5f15bSPeter Brune 568eed5f15bSPeter Brune PetscFunctionBegin; 569eed5f15bSPeter Brune jac = (SNES_Composite*)snes->data; 570eed5f15bSPeter Brune next = jac->head; 571eed5f15bSPeter Brune for (i=0; i<n; i++) { 572eed5f15bSPeter Brune if (!next->next) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_INCOMP,"Not enough SNESes in composite preconditioner"); 573eed5f15bSPeter Brune next = next->next; 574eed5f15bSPeter Brune } 575eed5f15bSPeter Brune *subsnes = next->snes; 576eed5f15bSPeter Brune PetscFunctionReturn(0); 577eed5f15bSPeter Brune } 578eed5f15bSPeter Brune 579eed5f15bSPeter Brune /* -------------------------------------------------------------------------------- */ 580eed5f15bSPeter Brune #undef __FUNCT__ 581eed5f15bSPeter Brune #define __FUNCT__ "SNESCompositeSetType" 582e31ab4f9SPeter Brune /*@C 583eed5f15bSPeter Brune SNESCompositeSetType - Sets the type of composite preconditioner. 584eed5f15bSPeter Brune 585eed5f15bSPeter Brune Logically Collective on SNES 586eed5f15bSPeter Brune 587eed5f15bSPeter Brune Input Parameter: 588eed5f15bSPeter Brune + snes - the preconditioner context 589eed5f15bSPeter Brune - type - SNES_COMPOSITE_ADDITIVE (default), SNES_COMPOSITE_MULTIPLICATIVE 590eed5f15bSPeter Brune 591eed5f15bSPeter Brune Options Database Key: 592eed5f15bSPeter Brune . -snes_composite_type <type: one of multiplicative, additive, special> - Sets composite preconditioner type 593eed5f15bSPeter Brune 594eed5f15bSPeter Brune Level: Developer 595eed5f15bSPeter Brune 596eed5f15bSPeter Brune .keywords: SNES, set, type, composite preconditioner, additive, multiplicative 597eed5f15bSPeter Brune @*/ 598eed5f15bSPeter Brune PetscErrorCode SNESCompositeSetType(SNES snes,SNESCompositeType type) 599eed5f15bSPeter Brune { 600eed5f15bSPeter Brune PetscErrorCode ierr; 601eed5f15bSPeter Brune 602eed5f15bSPeter Brune PetscFunctionBegin; 603eed5f15bSPeter Brune PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 604eed5f15bSPeter Brune PetscValidLogicalCollectiveEnum(snes,type,2); 605eed5f15bSPeter Brune ierr = PetscTryMethod(snes,"SNESCompositeSetType_C",(SNES,SNESCompositeType),(snes,type));CHKERRQ(ierr); 606eed5f15bSPeter Brune PetscFunctionReturn(0); 607eed5f15bSPeter Brune } 608eed5f15bSPeter Brune 609eed5f15bSPeter Brune #undef __FUNCT__ 610eed5f15bSPeter Brune #define __FUNCT__ "SNESCompositeAddSNES" 611eed5f15bSPeter Brune /*@C 612eed5f15bSPeter Brune SNESCompositeAddSNES - Adds another SNES to the composite SNES. 613eed5f15bSPeter Brune 614eed5f15bSPeter Brune Collective on SNES 615eed5f15bSPeter Brune 616eed5f15bSPeter Brune Input Parameters: 617eed5f15bSPeter Brune + snes - the preconditioner context 618eed5f15bSPeter Brune - type - the type of the new preconditioner 619eed5f15bSPeter Brune 620eed5f15bSPeter Brune Level: Developer 621eed5f15bSPeter Brune 622eed5f15bSPeter Brune .keywords: SNES, composite preconditioner, add 623eed5f15bSPeter Brune @*/ 624eed5f15bSPeter Brune PetscErrorCode SNESCompositeAddSNES(SNES snes,SNESType type) 625eed5f15bSPeter Brune { 626eed5f15bSPeter Brune PetscErrorCode ierr; 627eed5f15bSPeter Brune 628eed5f15bSPeter Brune PetscFunctionBegin; 629eed5f15bSPeter Brune PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 630eed5f15bSPeter Brune ierr = PetscTryMethod(snes,"SNESCompositeAddSNES_C",(SNES,SNESType),(snes,type));CHKERRQ(ierr); 631eed5f15bSPeter Brune PetscFunctionReturn(0); 632eed5f15bSPeter Brune } 633eed5f15bSPeter Brune #undef __FUNCT__ 634eed5f15bSPeter Brune #define __FUNCT__ "SNESCompositeGetSNES" 635eed5f15bSPeter Brune /*@ 636eed5f15bSPeter Brune SNESCompositeGetSNES - Gets one of the SNES objects in the composite SNES. 637eed5f15bSPeter Brune 638eed5f15bSPeter Brune Not Collective 639eed5f15bSPeter Brune 640eed5f15bSPeter Brune Input Parameter: 641eed5f15bSPeter Brune + snes - the preconditioner context 642eed5f15bSPeter Brune - n - the number of the snes requested 643eed5f15bSPeter Brune 644eed5f15bSPeter Brune Output Parameters: 645eed5f15bSPeter Brune . subsnes - the SNES requested 646eed5f15bSPeter Brune 647eed5f15bSPeter Brune Level: Developer 648eed5f15bSPeter Brune 649eed5f15bSPeter Brune .keywords: SNES, get, composite preconditioner, sub preconditioner 650eed5f15bSPeter Brune 651eed5f15bSPeter Brune .seealso: SNESCompositeAddSNES() 652eed5f15bSPeter Brune @*/ 653eed5f15bSPeter Brune PetscErrorCode SNESCompositeGetSNES(SNES snes,PetscInt n,SNES *subsnes) 654eed5f15bSPeter Brune { 655eed5f15bSPeter Brune PetscErrorCode ierr; 656eed5f15bSPeter Brune 657eed5f15bSPeter Brune PetscFunctionBegin; 658eed5f15bSPeter Brune PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 659eed5f15bSPeter Brune PetscValidPointer(subsnes,3); 660eed5f15bSPeter Brune ierr = PetscUseMethod(snes,"SNESCompositeGetSNES_C",(SNES,PetscInt,SNES*),(snes,n,subsnes));CHKERRQ(ierr); 661eed5f15bSPeter Brune PetscFunctionReturn(0); 662eed5f15bSPeter Brune } 663eed5f15bSPeter Brune 664eed5f15bSPeter Brune #undef __FUNCT__ 6656b2b7f7bSPatrick Farrell #define __FUNCT__ "SNESCompositeGetNumber" 6666b2b7f7bSPatrick Farrell /*@ 6676b2b7f7bSPatrick Farrell SNESCompositeGetNumber - Get the number of subsolvers in the composite SNES. 6686b2b7f7bSPatrick Farrell 6696b2b7f7bSPatrick Farrell Logically Collective on SNES 6706b2b7f7bSPatrick Farrell 6716b2b7f7bSPatrick Farrell Input Parameter: 6726b2b7f7bSPatrick Farrell snes - the preconditioner context 6736b2b7f7bSPatrick Farrell 6746b2b7f7bSPatrick Farrell Output Parameter: 6756b2b7f7bSPatrick Farrell n - the number of subsolvers 6766b2b7f7bSPatrick Farrell 6776b2b7f7bSPatrick Farrell Level: Developer 6786b2b7f7bSPatrick Farrell 6796b2b7f7bSPatrick Farrell .keywords: SNES, composite preconditioner 6806b2b7f7bSPatrick Farrell @*/ 6816b2b7f7bSPatrick Farrell PetscErrorCode SNESCompositeGetNumber(SNES snes,PetscInt *n) 6826b2b7f7bSPatrick Farrell { 6836b2b7f7bSPatrick Farrell SNES_Composite *jac; 6846b2b7f7bSPatrick Farrell SNES_CompositeLink next; 6856b2b7f7bSPatrick Farrell 6866b2b7f7bSPatrick Farrell PetscFunctionBegin; 6876b2b7f7bSPatrick Farrell jac = (SNES_Composite*)snes->data; 6886b2b7f7bSPatrick Farrell next = jac->head; 6896b2b7f7bSPatrick Farrell 6906b2b7f7bSPatrick Farrell *n = 0; 6916b2b7f7bSPatrick Farrell while (next) { 6926b2b7f7bSPatrick Farrell *n = *n + 1; 6936b2b7f7bSPatrick Farrell next = next->next; 6946b2b7f7bSPatrick Farrell } 6956b2b7f7bSPatrick Farrell PetscFunctionReturn(0); 6966b2b7f7bSPatrick Farrell } 6976b2b7f7bSPatrick Farrell 6986b2b7f7bSPatrick Farrell #undef __FUNCT__ 6998f626970SPeter Brune #define __FUNCT__ "SNESCompositeSetDamping_Composite" 7008f626970SPeter Brune static PetscErrorCode SNESCompositeSetDamping_Composite(SNES snes,PetscInt n,PetscReal dmp) 7018f626970SPeter Brune { 7028f626970SPeter Brune SNES_Composite *jac; 7038f626970SPeter Brune SNES_CompositeLink next; 7048f626970SPeter Brune PetscInt i; 7058f626970SPeter Brune 7068f626970SPeter Brune PetscFunctionBegin; 7078f626970SPeter Brune jac = (SNES_Composite*)snes->data; 7088f626970SPeter Brune next = jac->head; 7098f626970SPeter Brune for (i=0; i<n; i++) { 7108f626970SPeter Brune if (!next->next) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_INCOMP,"Not enough SNESes in composite preconditioner"); 7118f626970SPeter Brune next = next->next; 7128f626970SPeter Brune } 7138f626970SPeter Brune next->dmp = dmp; 7148f626970SPeter Brune PetscFunctionReturn(0); 7158f626970SPeter Brune } 7168f626970SPeter Brune 7178f626970SPeter Brune #undef __FUNCT__ 7188f626970SPeter Brune #define __FUNCT__ "SNESCompositeSetDamping" 7198f626970SPeter Brune /*@ 7208f626970SPeter Brune SNESCompositeSetDamping - Sets the damping of a subsolver when using additive composite SNES. 7218f626970SPeter Brune 7228f626970SPeter Brune Not Collective 7238f626970SPeter Brune 7248f626970SPeter Brune Input Parameter: 7258f626970SPeter Brune + snes - the preconditioner context 7268f626970SPeter Brune . n - the number of the snes requested 7278f626970SPeter Brune - dmp - the damping 7288f626970SPeter Brune 7298f626970SPeter Brune Level: Developer 7308f626970SPeter Brune 7318f626970SPeter Brune .keywords: SNES, get, composite preconditioner, sub preconditioner 7328f626970SPeter Brune 7338f626970SPeter Brune .seealso: SNESCompositeAddSNES() 7348f626970SPeter Brune @*/ 7358f626970SPeter Brune PetscErrorCode SNESCompositeSetDamping(SNES snes,PetscInt n,PetscReal dmp) 7368f626970SPeter Brune { 7378f626970SPeter Brune PetscErrorCode ierr; 7388f626970SPeter Brune 7398f626970SPeter Brune PetscFunctionBegin; 7408f626970SPeter Brune PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 7418f626970SPeter Brune ierr = PetscUseMethod(snes,"SNESCompositeSetDamping_C",(SNES,PetscInt,PetscReal),(snes,n,dmp));CHKERRQ(ierr); 7428f626970SPeter Brune PetscFunctionReturn(0); 7438f626970SPeter Brune } 7448f626970SPeter Brune 7458f626970SPeter Brune #undef __FUNCT__ 746eed5f15bSPeter Brune #define __FUNCT__ "SNESSolve_Composite" 747eed5f15bSPeter Brune PetscErrorCode SNESSolve_Composite(SNES snes) 748eed5f15bSPeter Brune { 749eed5f15bSPeter Brune Vec F; 750eed5f15bSPeter Brune Vec X; 751eed5f15bSPeter Brune Vec B; 752eed5f15bSPeter Brune PetscInt i; 753b22975d2SPatrick Farrell PetscReal fnorm = 0.0, xnorm = 0.0, snorm = 0.0; 754eed5f15bSPeter Brune PetscErrorCode ierr; 75572edecb9SPeter Brune SNESNormSchedule normtype; 756eed5f15bSPeter Brune SNES_Composite *comp = (SNES_Composite*)snes->data; 757eed5f15bSPeter Brune 758eed5f15bSPeter Brune PetscFunctionBegin; 759eed5f15bSPeter Brune X = snes->vec_sol; 760eed5f15bSPeter Brune F = snes->vec_func; 761eed5f15bSPeter Brune B = snes->vec_rhs; 762eed5f15bSPeter Brune 763e04113cfSBarry Smith ierr = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr); 764eed5f15bSPeter Brune snes->iter = 0; 765eed5f15bSPeter Brune snes->norm = 0.; 766e04113cfSBarry Smith ierr = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr); 767b22975d2SPatrick Farrell ierr = SNESSetWorkVecs(snes, 1);CHKERRQ(ierr); 768eed5f15bSPeter Brune snes->reason = SNES_CONVERGED_ITERATING; 76972edecb9SPeter Brune ierr = SNESGetNormSchedule(snes, &normtype);CHKERRQ(ierr); 77072edecb9SPeter Brune if (normtype == SNES_NORM_ALWAYS || normtype == SNES_NORM_INITIAL_ONLY || normtype == SNES_NORM_INITIAL_FINAL_ONLY) { 771eed5f15bSPeter Brune if (!snes->vec_func_init_set) { 772eed5f15bSPeter Brune ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 773eed5f15bSPeter Brune } else snes->vec_func_init_set = PETSC_FALSE; 774eed5f15bSPeter Brune 775a6da83ecSPatrick Farrell if (snes->xl && snes->xu) { 776a6da83ecSPatrick Farrell ierr = SNESVIComputeInactiveSetFnorm(snes, F, X, &fnorm);CHKERRQ(ierr); 777a6da83ecSPatrick Farrell } else { 778eed5f15bSPeter Brune ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F|| */ 779a6da83ecSPatrick Farrell } 780*422a814eSBarry Smith SNESCheckFunctionNorm(snes,fnorm); 781e04113cfSBarry Smith ierr = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr); 782eed5f15bSPeter Brune snes->iter = 0; 783eed5f15bSPeter Brune snes->norm = fnorm; 784e04113cfSBarry Smith ierr = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr); 785eed5f15bSPeter Brune ierr = SNESLogConvergenceHistory(snes,snes->norm,0);CHKERRQ(ierr); 786eed5f15bSPeter Brune ierr = SNESMonitor(snes,0,snes->norm);CHKERRQ(ierr); 787eed5f15bSPeter Brune 788eed5f15bSPeter Brune /* test convergence */ 789eed5f15bSPeter Brune ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 790eed5f15bSPeter Brune if (snes->reason) PetscFunctionReturn(0); 791eed5f15bSPeter Brune } else { 792e04113cfSBarry Smith ierr = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr); 793eed5f15bSPeter Brune ierr = SNESLogConvergenceHistory(snes,snes->norm,0);CHKERRQ(ierr); 794eed5f15bSPeter Brune ierr = SNESMonitor(snes,0,snes->norm);CHKERRQ(ierr); 795eed5f15bSPeter Brune } 796eed5f15bSPeter Brune 797eed5f15bSPeter Brune /* Call general purpose update function */ 798eed5f15bSPeter Brune if (snes->ops->update) { 799eed5f15bSPeter Brune ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 800eed5f15bSPeter Brune } 801eed5f15bSPeter Brune 802eed5f15bSPeter Brune for (i = 0; i < snes->max_its; i++) { 803b22975d2SPatrick Farrell /* Copy the state before modification by application of the composite solver; 804b22975d2SPatrick Farrell we will subtract the new state after application */ 805b22975d2SPatrick Farrell ierr = VecCopy(X, snes->work[0]);CHKERRQ(ierr); 806b22975d2SPatrick Farrell 807eed5f15bSPeter Brune if (comp->type == SNES_COMPOSITE_ADDITIVE) { 808eed5f15bSPeter Brune ierr = SNESCompositeApply_Additive(snes,X,B,F,&fnorm);CHKERRQ(ierr); 809eed5f15bSPeter Brune } else if (comp->type == SNES_COMPOSITE_MULTIPLICATIVE) { 810eed5f15bSPeter Brune ierr = SNESCompositeApply_Multiplicative(snes,X,B,F,&fnorm);CHKERRQ(ierr); 81190a8ba9bSPeter Brune } else if (comp->type == SNES_COMPOSITE_ADDITIVEOPTIMAL) { 81290a8ba9bSPeter Brune ierr = SNESCompositeApply_AdditiveOptimal(snes,X,B,F,&fnorm);CHKERRQ(ierr); 813eed5f15bSPeter Brune } else { 81490a8ba9bSPeter Brune SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"Unsupported SNESComposite type"); 815eed5f15bSPeter Brune } 816d5a53f06SPeter Brune if (snes->reason < 0) break; 817d5a53f06SPeter Brune 818b22975d2SPatrick Farrell /* Compute the solution update for convergence testing */ 819b22975d2SPatrick Farrell ierr = VecAXPY(snes->work[0], -1.0, X);CHKERRQ(ierr); 820b22975d2SPatrick Farrell ierr = VecScale(snes->work[0], -1.0);CHKERRQ(ierr); 821b22975d2SPatrick Farrell 822eed5f15bSPeter Brune if ((i == snes->max_its - 1) && (normtype == SNES_NORM_INITIAL_FINAL_ONLY || normtype == SNES_NORM_FINAL_ONLY)) { 823eed5f15bSPeter Brune ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 824b22975d2SPatrick Farrell 825a6da83ecSPatrick Farrell if (snes->xl && snes->xu) { 826a6da83ecSPatrick Farrell ierr = VecNormBegin(X, NORM_2, &xnorm);CHKERRQ(ierr); 827a6da83ecSPatrick Farrell ierr = VecNormBegin(snes->work[0], NORM_2, &snorm);CHKERRQ(ierr); 828a6da83ecSPatrick Farrell ierr = SNESVIComputeInactiveSetFnorm(snes, F, X, &fnorm);CHKERRQ(ierr); 829a6da83ecSPatrick Farrell ierr = VecNormEnd(X, NORM_2, &xnorm);CHKERRQ(ierr); 830a6da83ecSPatrick Farrell ierr = VecNormEnd(snes->work[0], NORM_2, &snorm);CHKERRQ(ierr); 831a6da83ecSPatrick Farrell } else { 832b22975d2SPatrick Farrell ierr = VecNormBegin(F, NORM_2, &fnorm);CHKERRQ(ierr); 833b22975d2SPatrick Farrell ierr = VecNormBegin(X, NORM_2, &xnorm);CHKERRQ(ierr); 834b22975d2SPatrick Farrell ierr = VecNormBegin(snes->work[0], NORM_2, &snorm);CHKERRQ(ierr); 835b22975d2SPatrick Farrell 836b22975d2SPatrick Farrell ierr = VecNormEnd(F, NORM_2, &fnorm);CHKERRQ(ierr); 837b22975d2SPatrick Farrell ierr = VecNormEnd(X, NORM_2, &xnorm);CHKERRQ(ierr); 838b22975d2SPatrick Farrell ierr = VecNormEnd(snes->work[0], NORM_2, &snorm);CHKERRQ(ierr); 839a6da83ecSPatrick Farrell } 840*422a814eSBarry Smith SNESCheckFunctionNorm(snes,fnorm); 841b22975d2SPatrick Farrell } else if (normtype == SNES_NORM_ALWAYS) { 842b22975d2SPatrick Farrell ierr = VecNormBegin(X, NORM_2, &xnorm);CHKERRQ(ierr); 843b22975d2SPatrick Farrell ierr = VecNormBegin(snes->work[0], NORM_2, &snorm);CHKERRQ(ierr); 844b22975d2SPatrick Farrell ierr = VecNormEnd(X, NORM_2, &xnorm);CHKERRQ(ierr); 845b22975d2SPatrick Farrell ierr = VecNormEnd(snes->work[0], NORM_2, &snorm);CHKERRQ(ierr); 846eed5f15bSPeter Brune } 847eed5f15bSPeter Brune /* Monitor convergence */ 848e04113cfSBarry Smith ierr = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr); 849eed5f15bSPeter Brune snes->iter = i+1; 850eed5f15bSPeter Brune snes->norm = fnorm; 851e04113cfSBarry Smith ierr = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr); 852eed5f15bSPeter Brune ierr = SNESLogConvergenceHistory(snes,snes->norm,0);CHKERRQ(ierr); 853eed5f15bSPeter Brune ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr); 854eed5f15bSPeter Brune /* Test for convergence */ 855b22975d2SPatrick Farrell if (normtype == SNES_NORM_ALWAYS) {ierr = (*snes->ops->converged)(snes,snes->iter,xnorm,snorm,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);} 856eed5f15bSPeter Brune if (snes->reason) break; 857eed5f15bSPeter Brune /* Call general purpose update function */ 858eed5f15bSPeter Brune if (snes->ops->update) {ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);} 859eed5f15bSPeter Brune } 86072edecb9SPeter Brune if (normtype == SNES_NORM_ALWAYS) { 861eed5f15bSPeter Brune if (i == snes->max_its) { 862eed5f15bSPeter Brune ierr = PetscInfo1(snes,"Maximum number of iterations has been reached: %D\n",snes->max_its);CHKERRQ(ierr); 863eed5f15bSPeter Brune if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT; 864eed5f15bSPeter Brune } 865eed5f15bSPeter Brune } else if (!snes->reason) snes->reason = SNES_CONVERGED_ITS; 866eed5f15bSPeter Brune PetscFunctionReturn(0); 867eed5f15bSPeter Brune } 868eed5f15bSPeter Brune 869eed5f15bSPeter Brune /* -------------------------------------------------------------------------------------------*/ 870eed5f15bSPeter Brune 871eed5f15bSPeter Brune /*MC 872eed5f15bSPeter Brune SNESCOMPOSITE - Build a preconditioner by composing together several nonlinear solvers 873eed5f15bSPeter Brune 874eed5f15bSPeter Brune Options Database Keys: 875eed5f15bSPeter Brune + -snes_composite_type <type: one of multiplicative, additive, symmetric_multiplicative, special> - Sets composite preconditioner type 876eed5f15bSPeter Brune - -snes_composite_sneses - <snes0,snes1,...> list of SNESes to compose 877eed5f15bSPeter Brune 878eed5f15bSPeter Brune Level: intermediate 879eed5f15bSPeter Brune 880eed5f15bSPeter Brune Concepts: composing solvers 881eed5f15bSPeter Brune 882eed5f15bSPeter Brune .seealso: SNESCreate(), SNESSetType(), SNESType (for list of available types), SNES, 883eed5f15bSPeter Brune SNESSHELL, SNESCompositeSetType(), SNESCompositeSpecialSetAlpha(), SNESCompositeAddSNES(), 884eed5f15bSPeter Brune SNESCompositeGetSNES() 885eed5f15bSPeter Brune 886eed5f15bSPeter Brune M*/ 887eed5f15bSPeter Brune 888eed5f15bSPeter Brune #undef __FUNCT__ 889eed5f15bSPeter Brune #define __FUNCT__ "SNESCreate_Composite" 890eed5f15bSPeter Brune PETSC_EXTERN PetscErrorCode SNESCreate_Composite(SNES snes) 891eed5f15bSPeter Brune { 892eed5f15bSPeter Brune PetscErrorCode ierr; 893eed5f15bSPeter Brune SNES_Composite *jac; 894eed5f15bSPeter Brune 895eed5f15bSPeter Brune PetscFunctionBegin; 896b00a9115SJed Brown ierr = PetscNewLog(snes,&jac);CHKERRQ(ierr); 897eed5f15bSPeter Brune 898eed5f15bSPeter Brune snes->ops->solve = SNESSolve_Composite; 899eed5f15bSPeter Brune snes->ops->setup = SNESSetUp_Composite; 900eed5f15bSPeter Brune snes->ops->reset = SNESReset_Composite; 901eed5f15bSPeter Brune snes->ops->destroy = SNESDestroy_Composite; 902eed5f15bSPeter Brune snes->ops->setfromoptions = SNESSetFromOptions_Composite; 903eed5f15bSPeter Brune snes->ops->view = SNESView_Composite; 904eed5f15bSPeter Brune 905eed5f15bSPeter Brune snes->data = (void*)jac; 90690a8ba9bSPeter Brune jac->type = SNES_COMPOSITE_ADDITIVEOPTIMAL; 90790a8ba9bSPeter Brune jac->Fes = NULL; 90890a8ba9bSPeter Brune jac->Xes = NULL; 9099c2f3473SPeter Brune jac->fnorms = NULL; 91090a8ba9bSPeter Brune jac->nsnes = 0; 911eed5f15bSPeter Brune jac->head = 0; 9125e806d2eSPeter Brune jac->stol = 0.1; 9135e806d2eSPeter Brune jac->rtol = 1.1; 914eed5f15bSPeter Brune 91590a8ba9bSPeter Brune jac->h = NULL; 91690a8ba9bSPeter Brune jac->s = NULL; 91790a8ba9bSPeter Brune jac->beta = NULL; 91890a8ba9bSPeter Brune jac->work = NULL; 91990a8ba9bSPeter Brune jac->rwork = NULL; 92090a8ba9bSPeter Brune 9218f626970SPeter Brune ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESCompositeSetType_C",SNESCompositeSetType_Composite);CHKERRQ(ierr); 9228f626970SPeter Brune ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESCompositeAddSNES_C",SNESCompositeAddSNES_Composite);CHKERRQ(ierr); 9238f626970SPeter Brune ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESCompositeGetSNES_C",SNESCompositeGetSNES_Composite);CHKERRQ(ierr); 9248f626970SPeter Brune ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESCompositeSetDamping_C",SNESCompositeSetDamping_Composite);CHKERRQ(ierr); 925eed5f15bSPeter Brune PetscFunctionReturn(0); 926eed5f15bSPeter Brune } 927eed5f15bSPeter Brune 928