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; 240b762f1fSPatrick Farrell PetscInt innerFailures; /* the number of inner failures we've seen */ 2590a8ba9bSPeter Brune 2690a8ba9bSPeter Brune /* context for ADDITIVEOPTIMAL */ 2790a8ba9bSPeter Brune Vec *Xes,*Fes; /* solution and residual vectors for the subsolvers */ 289c2f3473SPeter Brune PetscReal *fnorms; /* norms of the solutions */ 2990a8ba9bSPeter Brune PetscScalar *h; /* the matrix formed as q_ij = (rdot_i, rdot_j) */ 309c2f3473SPeter Brune PetscScalar *g; /* the dotproducts of the previous function with the candidate functions */ 319c2f3473SPeter Brune PetscBLASInt n; /* matrix dimension -- nsnes */ 3290a8ba9bSPeter Brune PetscBLASInt nrhs; /* the number of right hand sides */ 3390a8ba9bSPeter Brune PetscBLASInt lda; /* the padded matrix dimension */ 3490a8ba9bSPeter Brune PetscBLASInt ldb; /* the padded vector dimension */ 3590a8ba9bSPeter Brune PetscReal *s; /* the singular values */ 369c2f3473SPeter Brune PetscScalar *beta; /* the RHS and combination */ 3790a8ba9bSPeter Brune PetscReal rcond; /* the exit condition */ 3890a8ba9bSPeter Brune PetscBLASInt rank; /* the effective rank */ 3990a8ba9bSPeter Brune PetscScalar *work; /* the work vector */ 4090a8ba9bSPeter Brune PetscReal *rwork; /* the real work vector used for complex */ 4190a8ba9bSPeter Brune PetscBLASInt lwork; /* the size of the work vector */ 4290a8ba9bSPeter Brune PetscBLASInt info; /* the output condition */ 4390a8ba9bSPeter Brune 445e806d2eSPeter Brune PetscReal rtol; /* restart tolerance for accepting the combination */ 455e806d2eSPeter Brune PetscReal stol; /* restart tolerance for the combination */ 46eed5f15bSPeter Brune } SNES_Composite; 47eed5f15bSPeter Brune 48eed5f15bSPeter Brune #undef __FUNCT__ 49eed5f15bSPeter Brune #define __FUNCT__ "SNESCompositeApply_Multiplicative" 50eed5f15bSPeter Brune static PetscErrorCode SNESCompositeApply_Multiplicative(SNES snes,Vec X,Vec B,Vec F,PetscReal *fnorm) 51eed5f15bSPeter Brune { 52eed5f15bSPeter Brune PetscErrorCode ierr; 53eed5f15bSPeter Brune SNES_Composite *jac = (SNES_Composite*)snes->data; 54eed5f15bSPeter Brune SNES_CompositeLink next = jac->head; 55eed5f15bSPeter Brune Vec FSub; 56d5a53f06SPeter 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 } 63eed5f15bSPeter Brune ierr = SNESSolve(next->snes,B,X);CHKERRQ(ierr); 64d5a53f06SPeter Brune ierr = SNESGetConvergedReason(next->snes,&reason);CHKERRQ(ierr); 65d5a53f06SPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 660b762f1fSPatrick Farrell jac->innerFailures++; 670b762f1fSPatrick Farrell if (jac->innerFailures >= snes->maxFailures) { 68d5a53f06SPeter Brune snes->reason = SNES_DIVERGED_INNER; 69d5a53f06SPeter Brune PetscFunctionReturn(0); 70d5a53f06SPeter Brune } 710b762f1fSPatrick Farrell } 72eed5f15bSPeter Brune 73eed5f15bSPeter Brune while (next->next) { 74eed5f15bSPeter Brune /* only copy the function over in the case where the functions correspond */ 7572edecb9SPeter Brune if (next->snes->pcside == PC_RIGHT && next->snes->normschedule != SNES_NORM_NONE) { 76eed5f15bSPeter Brune ierr = SNESGetFunction(next->snes,&FSub,NULL,NULL);CHKERRQ(ierr); 77eed5f15bSPeter Brune next = next->next; 78eed5f15bSPeter Brune ierr = SNESSetInitialFunction(next->snes,FSub);CHKERRQ(ierr); 79eed5f15bSPeter Brune } else { 80eed5f15bSPeter Brune next = next->next; 81eed5f15bSPeter Brune } 82eed5f15bSPeter Brune ierr = SNESSolve(next->snes,B,X);CHKERRQ(ierr); 83d5a53f06SPeter Brune ierr = SNESGetConvergedReason(next->snes,&reason);CHKERRQ(ierr); 84d5a53f06SPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 850b762f1fSPatrick Farrell jac->innerFailures++; 860b762f1fSPatrick Farrell if (jac->innerFailures >= snes->maxFailures) { 87d5a53f06SPeter Brune snes->reason = SNES_DIVERGED_INNER; 88d5a53f06SPeter Brune PetscFunctionReturn(0); 89d5a53f06SPeter Brune } 90eed5f15bSPeter Brune } 910b762f1fSPatrick Farrell } 92eed5f15bSPeter Brune if (next->snes->pcside == PC_RIGHT) { 93eed5f15bSPeter Brune ierr = SNESGetFunction(next->snes,&FSub,NULL,NULL);CHKERRQ(ierr); 94eed5f15bSPeter Brune ierr = VecCopy(FSub,F);CHKERRQ(ierr); 95d5a53f06SPeter Brune if (fnorm) { 9663cdc2bdSPatrick Farrell if (snes->xl && snes->xu) { 9763cdc2bdSPatrick Farrell ierr = SNESVIComputeInactiveSetFnorm(snes, F, X, fnorm);CHKERRQ(ierr); 9863cdc2bdSPatrick Farrell } else { 9971dbe336SPeter Brune ierr = VecNorm(F, NORM_2, fnorm);CHKERRQ(ierr); 10063cdc2bdSPatrick Farrell } 101422a814eSBarry Smith SNESCheckFunctionNorm(snes,*fnorm); 102d5a53f06SPeter Brune } 10372edecb9SPeter Brune } else if (snes->normschedule == SNES_NORM_ALWAYS) { 104be95d8f1SBarry Smith ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 105d5a53f06SPeter Brune if (fnorm) { 106a6da83ecSPatrick Farrell if (snes->xl && snes->xu) { 107a6da83ecSPatrick Farrell ierr = SNESVIComputeInactiveSetFnorm(snes, F, X, fnorm);CHKERRQ(ierr); 108a6da83ecSPatrick Farrell } else { 109d5a53f06SPeter Brune ierr = VecNorm(F, NORM_2, fnorm);CHKERRQ(ierr); 110a6da83ecSPatrick Farrell } 111422a814eSBarry Smith SNESCheckFunctionNorm(snes,*fnorm); 112d5a53f06SPeter Brune } 113eed5f15bSPeter Brune } 114eed5f15bSPeter Brune PetscFunctionReturn(0); 115eed5f15bSPeter Brune } 116eed5f15bSPeter Brune 117eed5f15bSPeter Brune #undef __FUNCT__ 118eed5f15bSPeter Brune #define __FUNCT__ "SNESCompositeApply_Additive" 119eed5f15bSPeter Brune static PetscErrorCode SNESCompositeApply_Additive(SNES snes,Vec X,Vec B,Vec F,PetscReal *fnorm) 120eed5f15bSPeter Brune { 121eed5f15bSPeter Brune PetscErrorCode ierr; 122eed5f15bSPeter Brune SNES_Composite *jac = (SNES_Composite*)snes->data; 123eed5f15bSPeter Brune SNES_CompositeLink next = jac->head; 124eed5f15bSPeter Brune Vec Y,Xorig; 125d5a53f06SPeter Brune SNESConvergedReason reason; 126eed5f15bSPeter Brune 127eed5f15bSPeter Brune PetscFunctionBegin; 128eed5f15bSPeter Brune Y = snes->vec_sol_update; 129eed5f15bSPeter Brune if (!jac->Xorig) {ierr = VecDuplicate(X,&jac->Xorig);CHKERRQ(ierr);} 130eed5f15bSPeter Brune Xorig = jac->Xorig; 131302440fdSBarry Smith ierr = VecCopy(X,Xorig);CHKERRQ(ierr); 132eed5f15bSPeter Brune if (!next) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"No composite SNESes supplied via SNESCompositeAddSNES() or -snes_composite_sneses"); 13372edecb9SPeter Brune if (snes->normschedule == SNES_NORM_ALWAYS) { 134eed5f15bSPeter Brune ierr = SNESSetInitialFunction(next->snes,F);CHKERRQ(ierr); 135eed5f15bSPeter Brune while (next->next) { 136eed5f15bSPeter Brune next = next->next; 137eed5f15bSPeter Brune ierr = SNESSetInitialFunction(next->snes,F);CHKERRQ(ierr); 138eed5f15bSPeter Brune } 139eed5f15bSPeter Brune } 140eed5f15bSPeter Brune next = jac->head; 1418f626970SPeter Brune ierr = VecCopy(Xorig,Y);CHKERRQ(ierr); 142eed5f15bSPeter Brune ierr = SNESSolve(next->snes,B,Y);CHKERRQ(ierr); 143d5a53f06SPeter Brune ierr = SNESGetConvergedReason(next->snes,&reason);CHKERRQ(ierr); 144d5a53f06SPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 1450b762f1fSPatrick Farrell jac->innerFailures++; 1460b762f1fSPatrick Farrell if (jac->innerFailures >= snes->maxFailures) { 147d5a53f06SPeter Brune snes->reason = SNES_DIVERGED_INNER; 148d5a53f06SPeter Brune PetscFunctionReturn(0); 149d5a53f06SPeter Brune } 1500b762f1fSPatrick Farrell } 151eed5f15bSPeter Brune ierr = VecAXPY(Y,-1.0,Xorig);CHKERRQ(ierr); 1528f626970SPeter Brune ierr = VecAXPY(X,next->dmp,Y);CHKERRQ(ierr); 1538f626970SPeter Brune while (next->next) { 1548f626970SPeter Brune next = next->next; 1558f626970SPeter Brune ierr = VecCopy(Xorig,Y);CHKERRQ(ierr); 1568f626970SPeter Brune ierr = SNESSolve(next->snes,B,Y);CHKERRQ(ierr); 157d5a53f06SPeter Brune ierr = SNESGetConvergedReason(next->snes,&reason);CHKERRQ(ierr); 158d5a53f06SPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 1590b762f1fSPatrick Farrell jac->innerFailures++; 1600b762f1fSPatrick Farrell if (jac->innerFailures >= snes->maxFailures) { 161d5a53f06SPeter Brune snes->reason = SNES_DIVERGED_INNER; 162d5a53f06SPeter Brune PetscFunctionReturn(0); 163d5a53f06SPeter Brune } 1640b762f1fSPatrick Farrell } 1658f626970SPeter Brune ierr = VecAXPY(Y,-1.0,Xorig);CHKERRQ(ierr); 1668f626970SPeter Brune ierr = VecAXPY(X,next->dmp,Y);CHKERRQ(ierr); 167eed5f15bSPeter Brune } 16872edecb9SPeter Brune if (snes->normschedule == SNES_NORM_ALWAYS) { 169eed5f15bSPeter Brune ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 170a6da83ecSPatrick Farrell if (fnorm) { 171a6da83ecSPatrick Farrell if (snes->xl && snes->xu) { 172a6da83ecSPatrick Farrell ierr = SNESVIComputeInactiveSetFnorm(snes, F, X, fnorm);CHKERRQ(ierr); 173a6da83ecSPatrick Farrell } else { 174a6da83ecSPatrick Farrell ierr = VecNorm(F, NORM_2, fnorm);CHKERRQ(ierr); 175a6da83ecSPatrick Farrell } 176422a814eSBarry Smith SNESCheckFunctionNorm(snes,*fnorm); 177a6da83ecSPatrick Farrell } 178eed5f15bSPeter Brune } 179eed5f15bSPeter Brune PetscFunctionReturn(0); 180eed5f15bSPeter Brune } 18190a8ba9bSPeter Brune 18290a8ba9bSPeter Brune #undef __FUNCT__ 18390a8ba9bSPeter Brune #define __FUNCT__ "SNESCompositeApply_AdditiveOptimal" 18490a8ba9bSPeter Brune /* approximately solve the overdetermined system: 18590a8ba9bSPeter Brune 18690a8ba9bSPeter Brune 2*F(x_i)\cdot F(\x_j)\alpha_i = 0 18790a8ba9bSPeter Brune \alpha_i = 1 18890a8ba9bSPeter Brune 18990a8ba9bSPeter Brune Which minimizes the L2 norm of the linearization of: 19090a8ba9bSPeter Brune ||F(\sum_i \alpha_i*x_i)||^2 19190a8ba9bSPeter Brune 19290a8ba9bSPeter Brune With the constraint that \sum_i\alpha_i = 1 19390a8ba9bSPeter Brune Where x_i is the solution from the ith subsolver. 19490a8ba9bSPeter Brune */ 19590a8ba9bSPeter Brune static PetscErrorCode SNESCompositeApply_AdditiveOptimal(SNES snes,Vec X,Vec B,Vec F,PetscReal *fnorm) 19690a8ba9bSPeter Brune { 19790a8ba9bSPeter Brune PetscErrorCode ierr; 19890a8ba9bSPeter Brune SNES_Composite *jac = (SNES_Composite*)snes->data; 19990a8ba9bSPeter Brune SNES_CompositeLink next = jac->head; 20090a8ba9bSPeter Brune Vec *Xes = jac->Xes,*Fes = jac->Fes; 20190a8ba9bSPeter Brune PetscInt i,j; 2025e806d2eSPeter Brune PetscScalar tot,total,ftf; 2035e806d2eSPeter Brune PetscReal min_fnorm; 2045e806d2eSPeter Brune PetscInt min_i; 205d5a53f06SPeter Brune SNESConvergedReason reason; 20690a8ba9bSPeter Brune 20790a8ba9bSPeter Brune PetscFunctionBegin; 20890a8ba9bSPeter Brune if (!next) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"No composite SNESes supplied via SNESCompositeAddSNES() or -snes_composite_sneses"); 20958bdfa14SPeter Brune 21058bdfa14SPeter Brune if (snes->normschedule == SNES_NORM_ALWAYS) { 21158bdfa14SPeter Brune next = jac->head; 21258bdfa14SPeter Brune ierr = SNESSetInitialFunction(next->snes,F);CHKERRQ(ierr); 21358bdfa14SPeter Brune while (next->next) { 21458bdfa14SPeter Brune next = next->next; 21558bdfa14SPeter Brune ierr = SNESSetInitialFunction(next->snes,F);CHKERRQ(ierr); 21658bdfa14SPeter Brune } 21758bdfa14SPeter Brune } 21858bdfa14SPeter Brune 21990a8ba9bSPeter Brune next = jac->head; 22090a8ba9bSPeter Brune i = 0; 22190a8ba9bSPeter Brune ierr = VecCopy(X,Xes[i]);CHKERRQ(ierr); 22290a8ba9bSPeter Brune ierr = SNESSolve(next->snes,B,Xes[i]);CHKERRQ(ierr); 223d5a53f06SPeter Brune ierr = SNESGetConvergedReason(next->snes,&reason);CHKERRQ(ierr); 224d5a53f06SPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 2250b762f1fSPatrick Farrell jac->innerFailures++; 2260b762f1fSPatrick Farrell if (jac->innerFailures >= snes->maxFailures) { 227d5a53f06SPeter Brune snes->reason = SNES_DIVERGED_INNER; 228d5a53f06SPeter Brune PetscFunctionReturn(0); 229d5a53f06SPeter Brune } 2300b762f1fSPatrick Farrell } 23190a8ba9bSPeter Brune while (next->next) { 23290a8ba9bSPeter Brune i++; 23390a8ba9bSPeter Brune next = next->next; 23490a8ba9bSPeter Brune ierr = VecCopy(X,Xes[i]);CHKERRQ(ierr); 23590a8ba9bSPeter Brune ierr = SNESSolve(next->snes,B,Xes[i]);CHKERRQ(ierr); 236d5a53f06SPeter Brune ierr = SNESGetConvergedReason(next->snes,&reason);CHKERRQ(ierr); 237d5a53f06SPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 2380b762f1fSPatrick Farrell jac->innerFailures++; 2390b762f1fSPatrick Farrell if (jac->innerFailures >= snes->maxFailures) { 240d5a53f06SPeter Brune snes->reason = SNES_DIVERGED_INNER; 241d5a53f06SPeter Brune PetscFunctionReturn(0); 242d5a53f06SPeter Brune } 24390a8ba9bSPeter Brune } 2440b762f1fSPatrick Farrell } 24590a8ba9bSPeter Brune 24690a8ba9bSPeter Brune /* all the solutions are collected; combine optimally */ 24790a8ba9bSPeter Brune for (i=0;i<jac->n;i++) { 24890a8ba9bSPeter Brune for (j=0;j<i+1;j++) { 2499c2f3473SPeter Brune ierr = VecDotBegin(Fes[i],Fes[j],&jac->h[i + j*jac->n]);CHKERRQ(ierr); 25090a8ba9bSPeter Brune } 2519c2f3473SPeter Brune ierr = VecDotBegin(Fes[i],F,&jac->g[i]);CHKERRQ(ierr); 25290a8ba9bSPeter Brune } 2535e806d2eSPeter Brune 25490a8ba9bSPeter Brune for (i=0;i<jac->n;i++) { 25590a8ba9bSPeter Brune for (j=0;j<i+1;j++) { 2569c2f3473SPeter Brune ierr = VecDotEnd(Fes[i],Fes[j],&jac->h[i + j*jac->n]);CHKERRQ(ierr); 2579c2f3473SPeter Brune if (i == j) jac->fnorms[i] = PetscSqrtReal(PetscRealPart(jac->h[i + j*jac->n])); 25890a8ba9bSPeter Brune } 2599c2f3473SPeter Brune ierr = VecDotEnd(Fes[i],F,&jac->g[i]);CHKERRQ(ierr); 26090a8ba9bSPeter Brune } 26190a8ba9bSPeter Brune 2629c2f3473SPeter Brune ftf = (*fnorm)*(*fnorm); 26390a8ba9bSPeter Brune 26490a8ba9bSPeter Brune for (i=0; i<jac->n; i++) { 26590a8ba9bSPeter Brune for (j=i+1;j<jac->n;j++) { 2669c2f3473SPeter Brune jac->h[i + j*jac->n] = jac->h[j + i*jac->n]; 26790a8ba9bSPeter Brune } 26890a8ba9bSPeter Brune } 26990a8ba9bSPeter Brune 27090a8ba9bSPeter Brune for (i=0; i<jac->n; i++) { 2719c2f3473SPeter Brune for (j=0; j<jac->n; j++) { 2729c2f3473SPeter Brune jac->h[i + j*jac->n] = jac->h[i + j*jac->n] - jac->g[j] - jac->g[i] + ftf; 27390a8ba9bSPeter Brune } 2749c2f3473SPeter Brune jac->beta[i] = ftf - jac->g[i]; 2759c2f3473SPeter Brune } 27690a8ba9bSPeter Brune 27790a8ba9bSPeter Brune #if defined(PETSC_MISSING_LAPACK_GELSS) 27890a8ba9bSPeter Brune SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_SUP,"SNESCOMPOSITE with ADDITIVEOPTIMAL requires the LAPACK GELSS routine."); 27990a8ba9bSPeter Brune #else 28090a8ba9bSPeter Brune jac->info = 0; 28190a8ba9bSPeter Brune jac->rcond = -1.; 28290a8ba9bSPeter Brune ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr); 28390a8ba9bSPeter Brune #if defined(PETSC_USE_COMPLEX) 2849c2f3473SPeter 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)); 28590a8ba9bSPeter Brune #else 2869c2f3473SPeter 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)); 28790a8ba9bSPeter Brune #endif 28890a8ba9bSPeter Brune ierr = PetscFPTrapPop();CHKERRQ(ierr); 28990a8ba9bSPeter Brune if (jac->info < 0) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_LIB,"Bad argument to GELSS"); 29090a8ba9bSPeter Brune if (jac->info > 0) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_LIB,"SVD failed to converge"); 29190a8ba9bSPeter Brune #endif 2929c2f3473SPeter Brune tot = 0.; 2935e806d2eSPeter Brune total = 0.; 29490a8ba9bSPeter Brune for (i=0; i<jac->n; i++) { 295422a814eSBarry Smith if (snes->errorifnotconverged && PetscIsInfOrNanScalar(jac->beta[i])) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_LIB,"SVD generated inconsistent output"); 296c783eb9eSBarry Smith ierr = PetscInfo2(snes,"%D: %g\n",i,(double)PetscRealPart(jac->beta[i]));CHKERRQ(ierr); 2979c2f3473SPeter Brune tot += jac->beta[i]; 2985e806d2eSPeter Brune total += PetscAbsScalar(jac->beta[i]); 29990a8ba9bSPeter Brune } 3009c2f3473SPeter Brune ierr = VecScale(X,(1. - tot));CHKERRQ(ierr); 30190a8ba9bSPeter Brune ierr = VecMAXPY(X,jac->n,jac->beta,Xes);CHKERRQ(ierr); 30290a8ba9bSPeter Brune ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 303a6da83ecSPatrick Farrell 304a6da83ecSPatrick Farrell if (snes->xl && snes->xu) { 305a6da83ecSPatrick Farrell ierr = SNESVIComputeInactiveSetFnorm(snes, F, X, fnorm);CHKERRQ(ierr); 306a6da83ecSPatrick Farrell } else { 3079c2f3473SPeter Brune ierr = VecNorm(F, NORM_2, fnorm);CHKERRQ(ierr); 308a6da83ecSPatrick Farrell } 30990a8ba9bSPeter Brune 3105e806d2eSPeter Brune /* take the minimum-normed candidate if it beats the combination by a factor of rtol or the combination has stagnated */ 3115e806d2eSPeter Brune min_fnorm = jac->fnorms[0]; 3125e806d2eSPeter Brune min_i = 0; 3139c2f3473SPeter Brune for (i=0; i<jac->n; i++) { 3145e806d2eSPeter Brune if (jac->fnorms[i] < min_fnorm) { 3155e806d2eSPeter Brune min_fnorm = jac->fnorms[i]; 3165e806d2eSPeter Brune min_i = i; 3179c2f3473SPeter Brune } 3189c2f3473SPeter Brune } 3195e806d2eSPeter Brune 3205e806d2eSPeter Brune /* stagnation or divergence restart to the solution of the solver that failed the least */ 3215e806d2eSPeter Brune if (PetscRealPart(total) < jac->stol || min_fnorm*jac->rtol < *fnorm) { 32258bdfa14SPeter Brune ierr = VecCopy(jac->Xes[min_i],X);CHKERRQ(ierr); 32358bdfa14SPeter Brune ierr = VecCopy(jac->Fes[min_i],F);CHKERRQ(ierr); 3245e806d2eSPeter Brune *fnorm = min_fnorm; 3255e806d2eSPeter Brune } 32690a8ba9bSPeter Brune PetscFunctionReturn(0); 32790a8ba9bSPeter Brune } 32890a8ba9bSPeter Brune 329eed5f15bSPeter Brune #undef __FUNCT__ 330eed5f15bSPeter Brune #define __FUNCT__ "SNESSetUp_Composite" 331eed5f15bSPeter Brune static PetscErrorCode SNESSetUp_Composite(SNES snes) 332eed5f15bSPeter Brune { 333eed5f15bSPeter Brune PetscErrorCode ierr; 334dd515a93SPeter Brune DM dm; 335eed5f15bSPeter Brune SNES_Composite *jac = (SNES_Composite*)snes->data; 336eed5f15bSPeter Brune SNES_CompositeLink next = jac->head; 33790a8ba9bSPeter Brune PetscInt n=0,i; 33890a8ba9bSPeter Brune Vec F; 339eed5f15bSPeter Brune 340eed5f15bSPeter Brune PetscFunctionBegin; 341dd515a93SPeter Brune ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 34263cdc2bdSPatrick Farrell 34363cdc2bdSPatrick Farrell if (snes->ops->computevariablebounds) { 34463cdc2bdSPatrick Farrell /* SNESVI only ever calls computevariablebounds once, so calling it once here is justified */ 34563cdc2bdSPatrick Farrell if (!snes->xl) {ierr = VecDuplicate(snes->vec_sol,&snes->xl);CHKERRQ(ierr);} 34663cdc2bdSPatrick Farrell if (!snes->xu) {ierr = VecDuplicate(snes->vec_sol,&snes->xu);CHKERRQ(ierr);} 34763cdc2bdSPatrick Farrell ierr = (*snes->ops->computevariablebounds)(snes,snes->xl,snes->xu);CHKERRQ(ierr); 34863cdc2bdSPatrick Farrell } 34963cdc2bdSPatrick Farrell 350eed5f15bSPeter Brune while (next) { 35190a8ba9bSPeter Brune n++; 352dd515a93SPeter Brune ierr = SNESSetDM(next->snes,dm);CHKERRQ(ierr); 35363cdc2bdSPatrick Farrell ierr = SNESSetApplicationContext(next->snes, snes->user);CHKERRQ(ierr); 35463cdc2bdSPatrick Farrell if (snes->xl && snes->xu) { 35563cdc2bdSPatrick Farrell if (snes->ops->computevariablebounds) { 35663cdc2bdSPatrick Farrell ierr = SNESVISetComputeVariableBounds(next->snes, snes->ops->computevariablebounds);CHKERRQ(ierr); 35763cdc2bdSPatrick Farrell } else { 35863cdc2bdSPatrick Farrell ierr = SNESVISetVariableBounds(next->snes,snes->xl,snes->xu);CHKERRQ(ierr); 35963cdc2bdSPatrick Farrell } 36063cdc2bdSPatrick Farrell } 36163cdc2bdSPatrick Farrell 362eed5f15bSPeter Brune next = next->next; 363eed5f15bSPeter Brune } 36490a8ba9bSPeter Brune jac->nsnes = n; 36590a8ba9bSPeter Brune ierr = SNESGetFunction(snes,&F,NULL,NULL);CHKERRQ(ierr); 36690a8ba9bSPeter Brune if (jac->type == SNES_COMPOSITE_ADDITIVEOPTIMAL) { 36790a8ba9bSPeter Brune ierr = VecDuplicateVecs(F,jac->nsnes,&jac->Xes);CHKERRQ(ierr); 368854ce69bSBarry Smith ierr = PetscMalloc1(n,&jac->Fes);CHKERRQ(ierr); 369854ce69bSBarry Smith ierr = PetscMalloc1(n,&jac->fnorms);CHKERRQ(ierr); 37090a8ba9bSPeter Brune next = jac->head; 37190a8ba9bSPeter Brune i = 0; 37290a8ba9bSPeter Brune while (next) { 37390a8ba9bSPeter Brune ierr = SNESGetFunction(next->snes,&F,NULL,NULL);CHKERRQ(ierr); 37490a8ba9bSPeter Brune jac->Fes[i] = F; 37590a8ba9bSPeter Brune ierr = PetscObjectReference((PetscObject)F);CHKERRQ(ierr); 37690a8ba9bSPeter Brune next = next->next; 37790a8ba9bSPeter Brune i++; 37890a8ba9bSPeter Brune } 37990a8ba9bSPeter Brune /* allocate the subspace direct solve area */ 38090a8ba9bSPeter Brune jac->nrhs = 1; 3819c2f3473SPeter Brune jac->lda = jac->nsnes; 38290a8ba9bSPeter Brune jac->ldb = jac->nsnes; 38390a8ba9bSPeter Brune jac->n = jac->nsnes; 38490a8ba9bSPeter Brune 385785e854fSJed Brown ierr = PetscMalloc1(jac->n*jac->n,&jac->h);CHKERRQ(ierr); 386785e854fSJed Brown ierr = PetscMalloc1(jac->n,&jac->beta);CHKERRQ(ierr); 387785e854fSJed Brown ierr = PetscMalloc1(jac->n,&jac->s);CHKERRQ(ierr); 388785e854fSJed Brown ierr = PetscMalloc1(jac->n,&jac->g);CHKERRQ(ierr); 3899c2f3473SPeter Brune jac->lwork = 12*jac->n; 39090a8ba9bSPeter Brune #if PETSC_USE_COMPLEX 391854ce69bSBarry Smith ierr = PetscMalloc1(jac->lwork,&jac->rwork);CHKERRQ(ierr); 39290a8ba9bSPeter Brune #endif 393854ce69bSBarry Smith ierr = PetscMalloc1(jac->lwork,&jac->work);CHKERRQ(ierr); 39490a8ba9bSPeter Brune } 39590a8ba9bSPeter Brune 396eed5f15bSPeter Brune PetscFunctionReturn(0); 397eed5f15bSPeter Brune } 398eed5f15bSPeter Brune 399eed5f15bSPeter Brune #undef __FUNCT__ 400eed5f15bSPeter Brune #define __FUNCT__ "SNESReset_Composite" 401eed5f15bSPeter Brune static PetscErrorCode SNESReset_Composite(SNES snes) 402eed5f15bSPeter Brune { 403eed5f15bSPeter Brune SNES_Composite *jac = (SNES_Composite*)snes->data; 404eed5f15bSPeter Brune PetscErrorCode ierr; 405eed5f15bSPeter Brune SNES_CompositeLink next = jac->head; 406eed5f15bSPeter Brune 407eed5f15bSPeter Brune PetscFunctionBegin; 408eed5f15bSPeter Brune while (next) { 409eed5f15bSPeter Brune ierr = SNESReset(next->snes);CHKERRQ(ierr); 410eed5f15bSPeter Brune next = next->next; 411eed5f15bSPeter Brune } 412eed5f15bSPeter Brune ierr = VecDestroy(&jac->Xorig);CHKERRQ(ierr); 41390a8ba9bSPeter Brune if (jac->Xes) {ierr = VecDestroyVecs(jac->nsnes,&jac->Xes);CHKERRQ(ierr);} 41490a8ba9bSPeter Brune if (jac->Fes) {ierr = VecDestroyVecs(jac->nsnes,&jac->Fes);CHKERRQ(ierr);} 4159c2f3473SPeter Brune ierr = PetscFree(jac->fnorms);CHKERRQ(ierr); 41690a8ba9bSPeter Brune ierr = PetscFree(jac->h);CHKERRQ(ierr); 41790a8ba9bSPeter Brune ierr = PetscFree(jac->s);CHKERRQ(ierr); 4189c2f3473SPeter Brune ierr = PetscFree(jac->g);CHKERRQ(ierr); 41990a8ba9bSPeter Brune ierr = PetscFree(jac->beta);CHKERRQ(ierr); 42090a8ba9bSPeter Brune ierr = PetscFree(jac->work);CHKERRQ(ierr); 42190a8ba9bSPeter Brune ierr = PetscFree(jac->rwork);CHKERRQ(ierr); 422eed5f15bSPeter Brune PetscFunctionReturn(0); 423eed5f15bSPeter Brune } 424eed5f15bSPeter Brune 425eed5f15bSPeter Brune #undef __FUNCT__ 426eed5f15bSPeter Brune #define __FUNCT__ "SNESDestroy_Composite" 427eed5f15bSPeter Brune static PetscErrorCode SNESDestroy_Composite(SNES snes) 428eed5f15bSPeter Brune { 429eed5f15bSPeter Brune SNES_Composite *jac = (SNES_Composite*)snes->data; 430eed5f15bSPeter Brune PetscErrorCode ierr; 431eed5f15bSPeter Brune SNES_CompositeLink next = jac->head,next_tmp; 432eed5f15bSPeter Brune 433eed5f15bSPeter Brune PetscFunctionBegin; 434eed5f15bSPeter Brune ierr = SNESReset_Composite(snes);CHKERRQ(ierr); 435eed5f15bSPeter Brune while (next) { 436eed5f15bSPeter Brune ierr = SNESDestroy(&next->snes);CHKERRQ(ierr); 437eed5f15bSPeter Brune next_tmp = next; 438eed5f15bSPeter Brune next = next->next; 439eed5f15bSPeter Brune ierr = PetscFree(next_tmp);CHKERRQ(ierr); 440eed5f15bSPeter Brune } 441eed5f15bSPeter Brune ierr = PetscFree(snes->data);CHKERRQ(ierr); 442eed5f15bSPeter Brune PetscFunctionReturn(0); 443eed5f15bSPeter Brune } 444eed5f15bSPeter Brune 445eed5f15bSPeter Brune #undef __FUNCT__ 446eed5f15bSPeter Brune #define __FUNCT__ "SNESSetFromOptions_Composite" 4474416b707SBarry Smith static PetscErrorCode SNESSetFromOptions_Composite(PetscOptionItems *PetscOptionsObject,SNES snes) 448eed5f15bSPeter Brune { 449eed5f15bSPeter Brune SNES_Composite *jac = (SNES_Composite*)snes->data; 450eed5f15bSPeter Brune PetscErrorCode ierr; 451eed5f15bSPeter Brune PetscInt nmax = 8,i; 452eed5f15bSPeter Brune SNES_CompositeLink next; 453eed5f15bSPeter Brune char *sneses[8]; 4548f626970SPeter Brune PetscReal dmps[8]; 455eed5f15bSPeter Brune PetscBool flg; 456eed5f15bSPeter Brune 457eed5f15bSPeter Brune PetscFunctionBegin; 458e55864a3SBarry Smith ierr = PetscOptionsHead(PetscOptionsObject,"Composite preconditioner options");CHKERRQ(ierr); 459eed5f15bSPeter Brune ierr = PetscOptionsEnum("-snes_composite_type","Type of composition","SNESCompositeSetType",SNESCompositeTypes,(PetscEnum)jac->type,(PetscEnum*)&jac->type,&flg);CHKERRQ(ierr); 460eed5f15bSPeter Brune if (flg) { 461eed5f15bSPeter Brune ierr = SNESCompositeSetType(snes,jac->type);CHKERRQ(ierr); 462eed5f15bSPeter Brune } 463eed5f15bSPeter Brune ierr = PetscOptionsStringArray("-snes_composite_sneses","List of composite solvers","SNESCompositeAddSNES",sneses,&nmax,&flg);CHKERRQ(ierr); 464eed5f15bSPeter Brune if (flg) { 465eed5f15bSPeter Brune for (i=0; i<nmax; i++) { 466eed5f15bSPeter Brune ierr = SNESCompositeAddSNES(snes,sneses[i]);CHKERRQ(ierr); 467eed5f15bSPeter Brune ierr = PetscFree(sneses[i]);CHKERRQ(ierr); /* deallocate string sneses[i], which is allocated in PetscOptionsStringArray() */ 468eed5f15bSPeter Brune } 469eed5f15bSPeter Brune } 4708f626970SPeter Brune ierr = PetscOptionsRealArray("-snes_composite_damping","Damping of the additive composite solvers","SNESCompositeSetDamping",dmps,&nmax,&flg);CHKERRQ(ierr); 4718f626970SPeter Brune if (flg) { 4728f626970SPeter Brune for (i=0; i<nmax; i++) { 4738f626970SPeter Brune ierr = SNESCompositeSetDamping(snes,i,dmps[i]);CHKERRQ(ierr); 4748f626970SPeter Brune } 4758f626970SPeter Brune } 4765e806d2eSPeter Brune ierr = PetscOptionsReal("-snes_composite_stol","Step tolerance for restart on the additive composite solvers","",jac->stol,&jac->stol,NULL);CHKERRQ(ierr); 4775e806d2eSPeter Brune ierr = PetscOptionsReal("-snes_composite_rtol","Residual tolerance for the additive composite solvers","",jac->rtol,&jac->rtol,NULL);CHKERRQ(ierr); 478eed5f15bSPeter Brune ierr = PetscOptionsTail();CHKERRQ(ierr); 479eed5f15bSPeter Brune 480eed5f15bSPeter Brune next = jac->head; 481eed5f15bSPeter Brune while (next) { 482eed5f15bSPeter Brune ierr = SNESSetFromOptions(next->snes);CHKERRQ(ierr); 483eed5f15bSPeter Brune next = next->next; 484eed5f15bSPeter Brune } 485eed5f15bSPeter Brune PetscFunctionReturn(0); 486eed5f15bSPeter Brune } 487eed5f15bSPeter Brune 488eed5f15bSPeter Brune #undef __FUNCT__ 489eed5f15bSPeter Brune #define __FUNCT__ "SNESView_Composite" 490eed5f15bSPeter Brune static PetscErrorCode SNESView_Composite(SNES snes,PetscViewer viewer) 491eed5f15bSPeter Brune { 492eed5f15bSPeter Brune SNES_Composite *jac = (SNES_Composite*)snes->data; 493eed5f15bSPeter Brune PetscErrorCode ierr; 494eed5f15bSPeter Brune SNES_CompositeLink next = jac->head; 495eed5f15bSPeter Brune PetscBool iascii; 496eed5f15bSPeter Brune 497eed5f15bSPeter Brune PetscFunctionBegin; 498eed5f15bSPeter Brune ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 499eed5f15bSPeter Brune if (iascii) { 500eed5f15bSPeter Brune ierr = PetscViewerASCIIPrintf(viewer,"Composite SNES type - %s\n",SNESCompositeTypes[jac->type]);CHKERRQ(ierr); 501eed5f15bSPeter Brune ierr = PetscViewerASCIIPrintf(viewer,"SNESes on composite preconditioner follow\n");CHKERRQ(ierr); 502eed5f15bSPeter Brune ierr = PetscViewerASCIIPrintf(viewer,"---------------------------------\n");CHKERRQ(ierr); 503eed5f15bSPeter Brune } 504eed5f15bSPeter Brune if (iascii) { 505eed5f15bSPeter Brune ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 506eed5f15bSPeter Brune } 507eed5f15bSPeter Brune while (next) { 508eed5f15bSPeter Brune ierr = SNESView(next->snes,viewer);CHKERRQ(ierr); 509eed5f15bSPeter Brune next = next->next; 510eed5f15bSPeter Brune } 511eed5f15bSPeter Brune if (iascii) { 512eed5f15bSPeter Brune ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 513eed5f15bSPeter Brune ierr = PetscViewerASCIIPrintf(viewer,"---------------------------------\n");CHKERRQ(ierr); 514eed5f15bSPeter Brune } 515eed5f15bSPeter Brune PetscFunctionReturn(0); 516eed5f15bSPeter Brune } 517eed5f15bSPeter Brune 518eed5f15bSPeter Brune /* ------------------------------------------------------------------------------*/ 519eed5f15bSPeter Brune 520eed5f15bSPeter Brune #undef __FUNCT__ 521eed5f15bSPeter Brune #define __FUNCT__ "SNESCompositeSetType_Composite" 522eed5f15bSPeter Brune static PetscErrorCode SNESCompositeSetType_Composite(SNES snes,SNESCompositeType type) 523eed5f15bSPeter Brune { 524eed5f15bSPeter Brune SNES_Composite *jac = (SNES_Composite*)snes->data; 525eed5f15bSPeter Brune 526eed5f15bSPeter Brune PetscFunctionBegin; 527eed5f15bSPeter Brune jac->type = type; 528eed5f15bSPeter Brune PetscFunctionReturn(0); 529eed5f15bSPeter Brune } 530eed5f15bSPeter Brune 531eed5f15bSPeter Brune #undef __FUNCT__ 532eed5f15bSPeter Brune #define __FUNCT__ "SNESCompositeAddSNES_Composite" 533eed5f15bSPeter Brune static PetscErrorCode SNESCompositeAddSNES_Composite(SNES snes,SNESType type) 534eed5f15bSPeter Brune { 535eed5f15bSPeter Brune SNES_Composite *jac; 536eed5f15bSPeter Brune SNES_CompositeLink next,ilink; 537eed5f15bSPeter Brune PetscErrorCode ierr; 538eed5f15bSPeter Brune PetscInt cnt = 0; 539eed5f15bSPeter Brune const char *prefix; 540eed5f15bSPeter Brune char newprefix[8]; 541eed5f15bSPeter Brune DM dm; 542eed5f15bSPeter Brune 543eed5f15bSPeter Brune PetscFunctionBegin; 544b00a9115SJed Brown ierr = PetscNewLog(snes,&ilink);CHKERRQ(ierr); 545eed5f15bSPeter Brune ilink->next = 0; 546eed5f15bSPeter Brune ierr = SNESCreate(PetscObjectComm((PetscObject)snes),&ilink->snes);CHKERRQ(ierr); 547eed5f15bSPeter Brune ierr = PetscLogObjectParent((PetscObject)snes,(PetscObject)ilink->snes);CHKERRQ(ierr); 548eed5f15bSPeter Brune ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 549eed5f15bSPeter Brune ierr = SNESSetDM(ilink->snes,dm);CHKERRQ(ierr); 550cf5b3eb5SPeter Brune ierr = SNESSetTolerances(ilink->snes,ilink->snes->abstol,ilink->snes->rtol,ilink->snes->stol,1,ilink->snes->max_funcs);CHKERRQ(ierr); 551eed5f15bSPeter Brune jac = (SNES_Composite*)snes->data; 552eed5f15bSPeter Brune next = jac->head; 553eed5f15bSPeter Brune if (!next) { 554eed5f15bSPeter Brune jac->head = ilink; 555eed5f15bSPeter Brune ilink->previous = NULL; 556eed5f15bSPeter Brune } else { 557eed5f15bSPeter Brune cnt++; 558eed5f15bSPeter Brune while (next->next) { 559eed5f15bSPeter Brune next = next->next; 560eed5f15bSPeter Brune cnt++; 561eed5f15bSPeter Brune } 562eed5f15bSPeter Brune next->next = ilink; 563eed5f15bSPeter Brune ilink->previous = next; 564eed5f15bSPeter Brune } 565eed5f15bSPeter Brune ierr = SNESGetOptionsPrefix(snes,&prefix);CHKERRQ(ierr); 566eed5f15bSPeter Brune ierr = SNESSetOptionsPrefix(ilink->snes,prefix);CHKERRQ(ierr); 567eed5f15bSPeter Brune sprintf(newprefix,"sub_%d_",(int)cnt); 568eed5f15bSPeter Brune ierr = SNESAppendOptionsPrefix(ilink->snes,newprefix);CHKERRQ(ierr); 5698f626970SPeter Brune ierr = PetscObjectIncrementTabLevel((PetscObject)ilink->snes,(PetscObject)snes,1);CHKERRQ(ierr); 570eed5f15bSPeter Brune ierr = SNESSetType(ilink->snes,type);CHKERRQ(ierr); 57172edecb9SPeter Brune ierr = SNESSetNormSchedule(ilink->snes, SNES_NORM_FINAL_ONLY);CHKERRQ(ierr); 57263cdc2bdSPatrick Farrell 5738f626970SPeter Brune ilink->dmp = 1.0; 57490a8ba9bSPeter Brune jac->nsnes++; 575eed5f15bSPeter Brune PetscFunctionReturn(0); 576eed5f15bSPeter Brune } 577eed5f15bSPeter Brune 578eed5f15bSPeter Brune #undef __FUNCT__ 579eed5f15bSPeter Brune #define __FUNCT__ "SNESCompositeGetSNES_Composite" 580eed5f15bSPeter Brune static PetscErrorCode SNESCompositeGetSNES_Composite(SNES snes,PetscInt n,SNES *subsnes) 581eed5f15bSPeter Brune { 582eed5f15bSPeter Brune SNES_Composite *jac; 583eed5f15bSPeter Brune SNES_CompositeLink next; 584eed5f15bSPeter Brune PetscInt i; 585eed5f15bSPeter Brune 586eed5f15bSPeter Brune PetscFunctionBegin; 587eed5f15bSPeter Brune jac = (SNES_Composite*)snes->data; 588eed5f15bSPeter Brune next = jac->head; 589eed5f15bSPeter Brune for (i=0; i<n; i++) { 590eed5f15bSPeter Brune if (!next->next) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_INCOMP,"Not enough SNESes in composite preconditioner"); 591eed5f15bSPeter Brune next = next->next; 592eed5f15bSPeter Brune } 593eed5f15bSPeter Brune *subsnes = next->snes; 594eed5f15bSPeter Brune PetscFunctionReturn(0); 595eed5f15bSPeter Brune } 596eed5f15bSPeter Brune 597eed5f15bSPeter Brune /* -------------------------------------------------------------------------------- */ 598eed5f15bSPeter Brune #undef __FUNCT__ 599eed5f15bSPeter Brune #define __FUNCT__ "SNESCompositeSetType" 600e31ab4f9SPeter Brune /*@C 601eed5f15bSPeter Brune SNESCompositeSetType - Sets the type of composite preconditioner. 602eed5f15bSPeter Brune 603eed5f15bSPeter Brune Logically Collective on SNES 604eed5f15bSPeter Brune 605eed5f15bSPeter Brune Input Parameter: 606eed5f15bSPeter Brune + snes - the preconditioner context 607eed5f15bSPeter Brune - type - SNES_COMPOSITE_ADDITIVE (default), SNES_COMPOSITE_MULTIPLICATIVE 608eed5f15bSPeter Brune 609eed5f15bSPeter Brune Options Database Key: 610eed5f15bSPeter Brune . -snes_composite_type <type: one of multiplicative, additive, special> - Sets composite preconditioner type 611eed5f15bSPeter Brune 612eed5f15bSPeter Brune Level: Developer 613eed5f15bSPeter Brune 614eed5f15bSPeter Brune .keywords: SNES, set, type, composite preconditioner, additive, multiplicative 615eed5f15bSPeter Brune @*/ 616eed5f15bSPeter Brune PetscErrorCode SNESCompositeSetType(SNES snes,SNESCompositeType type) 617eed5f15bSPeter Brune { 618eed5f15bSPeter Brune PetscErrorCode ierr; 619eed5f15bSPeter Brune 620eed5f15bSPeter Brune PetscFunctionBegin; 621eed5f15bSPeter Brune PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 622eed5f15bSPeter Brune PetscValidLogicalCollectiveEnum(snes,type,2); 623eed5f15bSPeter Brune ierr = PetscTryMethod(snes,"SNESCompositeSetType_C",(SNES,SNESCompositeType),(snes,type));CHKERRQ(ierr); 624eed5f15bSPeter Brune PetscFunctionReturn(0); 625eed5f15bSPeter Brune } 626eed5f15bSPeter Brune 627eed5f15bSPeter Brune #undef __FUNCT__ 628eed5f15bSPeter Brune #define __FUNCT__ "SNESCompositeAddSNES" 629eed5f15bSPeter Brune /*@C 630eed5f15bSPeter Brune SNESCompositeAddSNES - Adds another SNES to the composite SNES. 631eed5f15bSPeter Brune 632eed5f15bSPeter Brune Collective on SNES 633eed5f15bSPeter Brune 634eed5f15bSPeter Brune Input Parameters: 635eed5f15bSPeter Brune + snes - the preconditioner context 636eed5f15bSPeter Brune - type - the type of the new preconditioner 637eed5f15bSPeter Brune 638eed5f15bSPeter Brune Level: Developer 639eed5f15bSPeter Brune 640eed5f15bSPeter Brune .keywords: SNES, composite preconditioner, add 641eed5f15bSPeter Brune @*/ 642eed5f15bSPeter Brune PetscErrorCode SNESCompositeAddSNES(SNES snes,SNESType type) 643eed5f15bSPeter Brune { 644eed5f15bSPeter Brune PetscErrorCode ierr; 645eed5f15bSPeter Brune 646eed5f15bSPeter Brune PetscFunctionBegin; 647eed5f15bSPeter Brune PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 648eed5f15bSPeter Brune ierr = PetscTryMethod(snes,"SNESCompositeAddSNES_C",(SNES,SNESType),(snes,type));CHKERRQ(ierr); 649eed5f15bSPeter Brune PetscFunctionReturn(0); 650eed5f15bSPeter Brune } 651eed5f15bSPeter Brune #undef __FUNCT__ 652eed5f15bSPeter Brune #define __FUNCT__ "SNESCompositeGetSNES" 653eed5f15bSPeter Brune /*@ 654eed5f15bSPeter Brune SNESCompositeGetSNES - Gets one of the SNES objects in the composite SNES. 655eed5f15bSPeter Brune 656eed5f15bSPeter Brune Not Collective 657eed5f15bSPeter Brune 658eed5f15bSPeter Brune Input Parameter: 659eed5f15bSPeter Brune + snes - the preconditioner context 660eed5f15bSPeter Brune - n - the number of the snes requested 661eed5f15bSPeter Brune 662eed5f15bSPeter Brune Output Parameters: 663eed5f15bSPeter Brune . subsnes - the SNES requested 664eed5f15bSPeter Brune 665eed5f15bSPeter Brune Level: Developer 666eed5f15bSPeter Brune 667eed5f15bSPeter Brune .keywords: SNES, get, composite preconditioner, sub preconditioner 668eed5f15bSPeter Brune 669eed5f15bSPeter Brune .seealso: SNESCompositeAddSNES() 670eed5f15bSPeter Brune @*/ 671eed5f15bSPeter Brune PetscErrorCode SNESCompositeGetSNES(SNES snes,PetscInt n,SNES *subsnes) 672eed5f15bSPeter Brune { 673eed5f15bSPeter Brune PetscErrorCode ierr; 674eed5f15bSPeter Brune 675eed5f15bSPeter Brune PetscFunctionBegin; 676eed5f15bSPeter Brune PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 677eed5f15bSPeter Brune PetscValidPointer(subsnes,3); 678eed5f15bSPeter Brune ierr = PetscUseMethod(snes,"SNESCompositeGetSNES_C",(SNES,PetscInt,SNES*),(snes,n,subsnes));CHKERRQ(ierr); 679eed5f15bSPeter Brune PetscFunctionReturn(0); 680eed5f15bSPeter Brune } 681eed5f15bSPeter Brune 682eed5f15bSPeter Brune #undef __FUNCT__ 6836b2b7f7bSPatrick Farrell #define __FUNCT__ "SNESCompositeGetNumber" 6846b2b7f7bSPatrick Farrell /*@ 6856b2b7f7bSPatrick Farrell SNESCompositeGetNumber - Get the number of subsolvers in the composite SNES. 6866b2b7f7bSPatrick Farrell 6876b2b7f7bSPatrick Farrell Logically Collective on SNES 6886b2b7f7bSPatrick Farrell 6896b2b7f7bSPatrick Farrell Input Parameter: 6906b2b7f7bSPatrick Farrell snes - the preconditioner context 6916b2b7f7bSPatrick Farrell 6926b2b7f7bSPatrick Farrell Output Parameter: 6936b2b7f7bSPatrick Farrell n - the number of subsolvers 6946b2b7f7bSPatrick Farrell 6956b2b7f7bSPatrick Farrell Level: Developer 6966b2b7f7bSPatrick Farrell 6976b2b7f7bSPatrick Farrell .keywords: SNES, composite preconditioner 6986b2b7f7bSPatrick Farrell @*/ 6996b2b7f7bSPatrick Farrell PetscErrorCode SNESCompositeGetNumber(SNES snes,PetscInt *n) 7006b2b7f7bSPatrick Farrell { 7016b2b7f7bSPatrick Farrell SNES_Composite *jac; 7026b2b7f7bSPatrick Farrell SNES_CompositeLink next; 7036b2b7f7bSPatrick Farrell 7046b2b7f7bSPatrick Farrell PetscFunctionBegin; 7056b2b7f7bSPatrick Farrell jac = (SNES_Composite*)snes->data; 7066b2b7f7bSPatrick Farrell next = jac->head; 7076b2b7f7bSPatrick Farrell 7086b2b7f7bSPatrick Farrell *n = 0; 7096b2b7f7bSPatrick Farrell while (next) { 7106b2b7f7bSPatrick Farrell *n = *n + 1; 7116b2b7f7bSPatrick Farrell next = next->next; 7126b2b7f7bSPatrick Farrell } 7136b2b7f7bSPatrick Farrell PetscFunctionReturn(0); 7146b2b7f7bSPatrick Farrell } 7156b2b7f7bSPatrick Farrell 7166b2b7f7bSPatrick Farrell #undef __FUNCT__ 7178f626970SPeter Brune #define __FUNCT__ "SNESCompositeSetDamping_Composite" 7188f626970SPeter Brune static PetscErrorCode SNESCompositeSetDamping_Composite(SNES snes,PetscInt n,PetscReal dmp) 7198f626970SPeter Brune { 7208f626970SPeter Brune SNES_Composite *jac; 7218f626970SPeter Brune SNES_CompositeLink next; 7228f626970SPeter Brune PetscInt i; 7238f626970SPeter Brune 7248f626970SPeter Brune PetscFunctionBegin; 7258f626970SPeter Brune jac = (SNES_Composite*)snes->data; 7268f626970SPeter Brune next = jac->head; 7278f626970SPeter Brune for (i=0; i<n; i++) { 7288f626970SPeter Brune if (!next->next) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_INCOMP,"Not enough SNESes in composite preconditioner"); 7298f626970SPeter Brune next = next->next; 7308f626970SPeter Brune } 7318f626970SPeter Brune next->dmp = dmp; 7328f626970SPeter Brune PetscFunctionReturn(0); 7338f626970SPeter Brune } 7348f626970SPeter Brune 7358f626970SPeter Brune #undef __FUNCT__ 7368f626970SPeter Brune #define __FUNCT__ "SNESCompositeSetDamping" 7378f626970SPeter Brune /*@ 7388f626970SPeter Brune SNESCompositeSetDamping - Sets the damping of a subsolver when using additive composite SNES. 7398f626970SPeter Brune 7408f626970SPeter Brune Not Collective 7418f626970SPeter Brune 7428f626970SPeter Brune Input Parameter: 7438f626970SPeter Brune + snes - the preconditioner context 7448f626970SPeter Brune . n - the number of the snes requested 7458f626970SPeter Brune - dmp - the damping 7468f626970SPeter Brune 7478f626970SPeter Brune Level: Developer 7488f626970SPeter Brune 7498f626970SPeter Brune .keywords: SNES, get, composite preconditioner, sub preconditioner 7508f626970SPeter Brune 7518f626970SPeter Brune .seealso: SNESCompositeAddSNES() 7528f626970SPeter Brune @*/ 7538f626970SPeter Brune PetscErrorCode SNESCompositeSetDamping(SNES snes,PetscInt n,PetscReal dmp) 7548f626970SPeter Brune { 7558f626970SPeter Brune PetscErrorCode ierr; 7568f626970SPeter Brune 7578f626970SPeter Brune PetscFunctionBegin; 7588f626970SPeter Brune PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 7598f626970SPeter Brune ierr = PetscUseMethod(snes,"SNESCompositeSetDamping_C",(SNES,PetscInt,PetscReal),(snes,n,dmp));CHKERRQ(ierr); 7608f626970SPeter Brune PetscFunctionReturn(0); 7618f626970SPeter Brune } 7628f626970SPeter Brune 7638f626970SPeter Brune #undef __FUNCT__ 764eed5f15bSPeter Brune #define __FUNCT__ "SNESSolve_Composite" 765eed5f15bSPeter Brune PetscErrorCode SNESSolve_Composite(SNES snes) 766eed5f15bSPeter Brune { 767eed5f15bSPeter Brune Vec F; 768eed5f15bSPeter Brune Vec X; 769eed5f15bSPeter Brune Vec B; 770eed5f15bSPeter Brune PetscInt i; 771b22975d2SPatrick Farrell PetscReal fnorm = 0.0, xnorm = 0.0, snorm = 0.0; 772eed5f15bSPeter Brune PetscErrorCode ierr; 77372edecb9SPeter Brune SNESNormSchedule normtype; 774eed5f15bSPeter Brune SNES_Composite *comp = (SNES_Composite*)snes->data; 775eed5f15bSPeter Brune 776eed5f15bSPeter Brune PetscFunctionBegin; 777eed5f15bSPeter Brune X = snes->vec_sol; 778eed5f15bSPeter Brune F = snes->vec_func; 779eed5f15bSPeter Brune B = snes->vec_rhs; 780eed5f15bSPeter Brune 781e04113cfSBarry Smith ierr = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr); 782eed5f15bSPeter Brune snes->iter = 0; 783eed5f15bSPeter Brune snes->norm = 0.; 7840b762f1fSPatrick Farrell comp->innerFailures = 0; 785e04113cfSBarry Smith ierr = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr); 786b22975d2SPatrick Farrell ierr = SNESSetWorkVecs(snes, 1);CHKERRQ(ierr); 787eed5f15bSPeter Brune snes->reason = SNES_CONVERGED_ITERATING; 78872edecb9SPeter Brune ierr = SNESGetNormSchedule(snes, &normtype);CHKERRQ(ierr); 78972edecb9SPeter Brune if (normtype == SNES_NORM_ALWAYS || normtype == SNES_NORM_INITIAL_ONLY || normtype == SNES_NORM_INITIAL_FINAL_ONLY) { 790eed5f15bSPeter Brune if (!snes->vec_func_init_set) { 791eed5f15bSPeter Brune ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 792eed5f15bSPeter Brune } else snes->vec_func_init_set = PETSC_FALSE; 793eed5f15bSPeter Brune 794a6da83ecSPatrick Farrell if (snes->xl && snes->xu) { 795a6da83ecSPatrick Farrell ierr = SNESVIComputeInactiveSetFnorm(snes, F, X, &fnorm);CHKERRQ(ierr); 796a6da83ecSPatrick Farrell } else { 797eed5f15bSPeter Brune ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F|| */ 798a6da83ecSPatrick Farrell } 799422a814eSBarry Smith SNESCheckFunctionNorm(snes,fnorm); 800e04113cfSBarry Smith ierr = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr); 801eed5f15bSPeter Brune snes->iter = 0; 802eed5f15bSPeter Brune snes->norm = fnorm; 803e04113cfSBarry Smith ierr = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr); 804eed5f15bSPeter Brune ierr = SNESLogConvergenceHistory(snes,snes->norm,0);CHKERRQ(ierr); 805eed5f15bSPeter Brune ierr = SNESMonitor(snes,0,snes->norm);CHKERRQ(ierr); 806eed5f15bSPeter Brune 807eed5f15bSPeter Brune /* test convergence */ 808eed5f15bSPeter Brune ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 809eed5f15bSPeter Brune if (snes->reason) PetscFunctionReturn(0); 810eed5f15bSPeter Brune } else { 811e04113cfSBarry Smith ierr = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr); 812eed5f15bSPeter Brune ierr = SNESLogConvergenceHistory(snes,snes->norm,0);CHKERRQ(ierr); 813eed5f15bSPeter Brune ierr = SNESMonitor(snes,0,snes->norm);CHKERRQ(ierr); 814eed5f15bSPeter Brune } 815eed5f15bSPeter Brune 816eed5f15bSPeter Brune /* Call general purpose update function */ 817eed5f15bSPeter Brune if (snes->ops->update) { 818eed5f15bSPeter Brune ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 819eed5f15bSPeter Brune } 820eed5f15bSPeter Brune 821eed5f15bSPeter Brune for (i = 0; i < snes->max_its; i++) { 822b22975d2SPatrick Farrell /* Copy the state before modification by application of the composite solver; 823b22975d2SPatrick Farrell we will subtract the new state after application */ 824b22975d2SPatrick Farrell ierr = VecCopy(X, snes->work[0]);CHKERRQ(ierr); 825b22975d2SPatrick Farrell 826eed5f15bSPeter Brune if (comp->type == SNES_COMPOSITE_ADDITIVE) { 827eed5f15bSPeter Brune ierr = SNESCompositeApply_Additive(snes,X,B,F,&fnorm);CHKERRQ(ierr); 828eed5f15bSPeter Brune } else if (comp->type == SNES_COMPOSITE_MULTIPLICATIVE) { 829eed5f15bSPeter Brune ierr = SNESCompositeApply_Multiplicative(snes,X,B,F,&fnorm);CHKERRQ(ierr); 83090a8ba9bSPeter Brune } else if (comp->type == SNES_COMPOSITE_ADDITIVEOPTIMAL) { 83190a8ba9bSPeter Brune ierr = SNESCompositeApply_AdditiveOptimal(snes,X,B,F,&fnorm);CHKERRQ(ierr); 8326c4ed002SBarry Smith } else SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"Unsupported SNESComposite type"); 833d5a53f06SPeter Brune if (snes->reason < 0) break; 834d5a53f06SPeter Brune 835b22975d2SPatrick Farrell /* Compute the solution update for convergence testing */ 836b22975d2SPatrick Farrell ierr = VecAXPY(snes->work[0], -1.0, X);CHKERRQ(ierr); 837b22975d2SPatrick Farrell ierr = VecScale(snes->work[0], -1.0);CHKERRQ(ierr); 838b22975d2SPatrick Farrell 839eed5f15bSPeter Brune if ((i == snes->max_its - 1) && (normtype == SNES_NORM_INITIAL_FINAL_ONLY || normtype == SNES_NORM_FINAL_ONLY)) { 840eed5f15bSPeter Brune ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 841b22975d2SPatrick Farrell 842a6da83ecSPatrick Farrell if (snes->xl && snes->xu) { 843a6da83ecSPatrick Farrell ierr = VecNormBegin(X, NORM_2, &xnorm);CHKERRQ(ierr); 844a6da83ecSPatrick Farrell ierr = VecNormBegin(snes->work[0], NORM_2, &snorm);CHKERRQ(ierr); 845a6da83ecSPatrick Farrell ierr = SNESVIComputeInactiveSetFnorm(snes, F, X, &fnorm);CHKERRQ(ierr); 846a6da83ecSPatrick Farrell ierr = VecNormEnd(X, NORM_2, &xnorm);CHKERRQ(ierr); 847a6da83ecSPatrick Farrell ierr = VecNormEnd(snes->work[0], NORM_2, &snorm);CHKERRQ(ierr); 848a6da83ecSPatrick Farrell } else { 849b22975d2SPatrick Farrell ierr = VecNormBegin(F, NORM_2, &fnorm);CHKERRQ(ierr); 850b22975d2SPatrick Farrell ierr = VecNormBegin(X, NORM_2, &xnorm);CHKERRQ(ierr); 851b22975d2SPatrick Farrell ierr = VecNormBegin(snes->work[0], NORM_2, &snorm);CHKERRQ(ierr); 852b22975d2SPatrick Farrell 853b22975d2SPatrick Farrell ierr = VecNormEnd(F, NORM_2, &fnorm);CHKERRQ(ierr); 854b22975d2SPatrick Farrell ierr = VecNormEnd(X, NORM_2, &xnorm);CHKERRQ(ierr); 855b22975d2SPatrick Farrell ierr = VecNormEnd(snes->work[0], NORM_2, &snorm);CHKERRQ(ierr); 856a6da83ecSPatrick Farrell } 857422a814eSBarry Smith SNESCheckFunctionNorm(snes,fnorm); 858b22975d2SPatrick Farrell } else if (normtype == SNES_NORM_ALWAYS) { 859b22975d2SPatrick Farrell ierr = VecNormBegin(X, NORM_2, &xnorm);CHKERRQ(ierr); 860b22975d2SPatrick Farrell ierr = VecNormBegin(snes->work[0], NORM_2, &snorm);CHKERRQ(ierr); 861b22975d2SPatrick Farrell ierr = VecNormEnd(X, NORM_2, &xnorm);CHKERRQ(ierr); 862b22975d2SPatrick Farrell ierr = VecNormEnd(snes->work[0], NORM_2, &snorm);CHKERRQ(ierr); 863eed5f15bSPeter Brune } 864eed5f15bSPeter Brune /* Monitor convergence */ 865e04113cfSBarry Smith ierr = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr); 866eed5f15bSPeter Brune snes->iter = i+1; 867eed5f15bSPeter Brune snes->norm = fnorm; 868e04113cfSBarry Smith ierr = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr); 869eed5f15bSPeter Brune ierr = SNESLogConvergenceHistory(snes,snes->norm,0);CHKERRQ(ierr); 870eed5f15bSPeter Brune ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr); 871eed5f15bSPeter Brune /* Test for convergence */ 872b22975d2SPatrick Farrell if (normtype == SNES_NORM_ALWAYS) {ierr = (*snes->ops->converged)(snes,snes->iter,xnorm,snorm,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);} 873eed5f15bSPeter Brune if (snes->reason) break; 874eed5f15bSPeter Brune /* Call general purpose update function */ 875eed5f15bSPeter Brune if (snes->ops->update) {ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);} 876eed5f15bSPeter Brune } 87772edecb9SPeter Brune if (normtype == SNES_NORM_ALWAYS) { 878eed5f15bSPeter Brune if (i == snes->max_its) { 879eed5f15bSPeter Brune ierr = PetscInfo1(snes,"Maximum number of iterations has been reached: %D\n",snes->max_its);CHKERRQ(ierr); 880eed5f15bSPeter Brune if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT; 881eed5f15bSPeter Brune } 882eed5f15bSPeter Brune } else if (!snes->reason) snes->reason = SNES_CONVERGED_ITS; 883eed5f15bSPeter Brune PetscFunctionReturn(0); 884eed5f15bSPeter Brune } 885eed5f15bSPeter Brune 886eed5f15bSPeter Brune /* -------------------------------------------------------------------------------------------*/ 887eed5f15bSPeter Brune 888eed5f15bSPeter Brune /*MC 889eed5f15bSPeter Brune SNESCOMPOSITE - Build a preconditioner by composing together several nonlinear solvers 890eed5f15bSPeter Brune 891eed5f15bSPeter Brune Options Database Keys: 892eed5f15bSPeter Brune + -snes_composite_type <type: one of multiplicative, additive, symmetric_multiplicative, special> - Sets composite preconditioner type 893eed5f15bSPeter Brune - -snes_composite_sneses - <snes0,snes1,...> list of SNESes to compose 894eed5f15bSPeter Brune 895eed5f15bSPeter Brune Level: intermediate 896eed5f15bSPeter Brune 897eed5f15bSPeter Brune Concepts: composing solvers 898eed5f15bSPeter Brune 899eed5f15bSPeter Brune .seealso: SNESCreate(), SNESSetType(), SNESType (for list of available types), SNES, 900eed5f15bSPeter Brune SNESSHELL, SNESCompositeSetType(), SNESCompositeSpecialSetAlpha(), SNESCompositeAddSNES(), 901eed5f15bSPeter Brune SNESCompositeGetSNES() 902eed5f15bSPeter Brune 903*4f02bc6aSBarry Smith References: 904*4f02bc6aSBarry Smith 905*4f02bc6aSBarry Smith Peter R. Brune, Matthew G. Knepley, Barry F. Smith, and Xuemin Tu, "Composing Scalable Nonlinear Algebraic Solvers", 906*4f02bc6aSBarry Smith SIAM Review, 57(4), 2015 907*4f02bc6aSBarry Smith 908eed5f15bSPeter Brune M*/ 909eed5f15bSPeter Brune 910eed5f15bSPeter Brune #undef __FUNCT__ 911eed5f15bSPeter Brune #define __FUNCT__ "SNESCreate_Composite" 912eed5f15bSPeter Brune PETSC_EXTERN PetscErrorCode SNESCreate_Composite(SNES snes) 913eed5f15bSPeter Brune { 914eed5f15bSPeter Brune PetscErrorCode ierr; 915eed5f15bSPeter Brune SNES_Composite *jac; 916eed5f15bSPeter Brune 917eed5f15bSPeter Brune PetscFunctionBegin; 918b00a9115SJed Brown ierr = PetscNewLog(snes,&jac);CHKERRQ(ierr); 919eed5f15bSPeter Brune 920eed5f15bSPeter Brune snes->ops->solve = SNESSolve_Composite; 921eed5f15bSPeter Brune snes->ops->setup = SNESSetUp_Composite; 922eed5f15bSPeter Brune snes->ops->reset = SNESReset_Composite; 923eed5f15bSPeter Brune snes->ops->destroy = SNESDestroy_Composite; 924eed5f15bSPeter Brune snes->ops->setfromoptions = SNESSetFromOptions_Composite; 925eed5f15bSPeter Brune snes->ops->view = SNESView_Composite; 926eed5f15bSPeter Brune 927eed5f15bSPeter Brune snes->data = (void*)jac; 92890a8ba9bSPeter Brune jac->type = SNES_COMPOSITE_ADDITIVEOPTIMAL; 92990a8ba9bSPeter Brune jac->Fes = NULL; 93090a8ba9bSPeter Brune jac->Xes = NULL; 9319c2f3473SPeter Brune jac->fnorms = NULL; 93290a8ba9bSPeter Brune jac->nsnes = 0; 933eed5f15bSPeter Brune jac->head = 0; 9345e806d2eSPeter Brune jac->stol = 0.1; 9355e806d2eSPeter Brune jac->rtol = 1.1; 936eed5f15bSPeter Brune 93790a8ba9bSPeter Brune jac->h = NULL; 93890a8ba9bSPeter Brune jac->s = NULL; 93990a8ba9bSPeter Brune jac->beta = NULL; 94090a8ba9bSPeter Brune jac->work = NULL; 94190a8ba9bSPeter Brune jac->rwork = NULL; 94290a8ba9bSPeter Brune 9438f626970SPeter Brune ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESCompositeSetType_C",SNESCompositeSetType_Composite);CHKERRQ(ierr); 9448f626970SPeter Brune ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESCompositeAddSNES_C",SNESCompositeAddSNES_Composite);CHKERRQ(ierr); 9458f626970SPeter Brune ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESCompositeGetSNES_C",SNESCompositeGetSNES_Composite);CHKERRQ(ierr); 9468f626970SPeter Brune ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESCompositeSetDamping_C",SNESCompositeSetDamping_Composite);CHKERRQ(ierr); 947eed5f15bSPeter Brune PetscFunctionReturn(0); 948eed5f15bSPeter Brune } 949