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 89e5d0892SLisandro Dalcin const char *const SNESCompositeTypes[] = {"ADDITIVE","MULTIPLICATIVE","ADDITIVEOPTIMAL","SNESCompositeType","SNES_COMPOSITE",NULL}; 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 static PetscErrorCode SNESCompositeApply_Multiplicative(SNES snes,Vec X,Vec B,Vec F,PetscReal *fnorm) 49eed5f15bSPeter Brune { 50eed5f15bSPeter Brune PetscErrorCode ierr; 51eed5f15bSPeter Brune SNES_Composite *jac = (SNES_Composite*)snes->data; 52eed5f15bSPeter Brune SNES_CompositeLink next = jac->head; 53eed5f15bSPeter Brune Vec FSub; 54d5a53f06SPeter Brune SNESConvergedReason reason; 55eed5f15bSPeter Brune 56eed5f15bSPeter Brune PetscFunctionBegin; 57eed5f15bSPeter Brune if (!next) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"No composite SNESes supplied via SNESCompositeAddSNES() or -snes_composite_sneses"); 5872edecb9SPeter Brune if (snes->normschedule == SNES_NORM_ALWAYS) { 59eed5f15bSPeter Brune ierr = SNESSetInitialFunction(next->snes,F);CHKERRQ(ierr); 60eed5f15bSPeter Brune } 61eed5f15bSPeter Brune ierr = SNESSolve(next->snes,B,X);CHKERRQ(ierr); 62d5a53f06SPeter Brune ierr = SNESGetConvergedReason(next->snes,&reason);CHKERRQ(ierr); 63d5a53f06SPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 640b762f1fSPatrick Farrell jac->innerFailures++; 650b762f1fSPatrick Farrell if (jac->innerFailures >= snes->maxFailures) { 66d5a53f06SPeter Brune snes->reason = SNES_DIVERGED_INNER; 67d5a53f06SPeter Brune PetscFunctionReturn(0); 68d5a53f06SPeter Brune } 690b762f1fSPatrick Farrell } 70eed5f15bSPeter Brune 71eed5f15bSPeter Brune while (next->next) { 72eed5f15bSPeter Brune /* only copy the function over in the case where the functions correspond */ 73efd4aadfSBarry Smith if (next->snes->npcside== PC_RIGHT && next->snes->normschedule != SNES_NORM_NONE) { 74eed5f15bSPeter Brune ierr = SNESGetFunction(next->snes,&FSub,NULL,NULL);CHKERRQ(ierr); 75eed5f15bSPeter Brune next = next->next; 76eed5f15bSPeter Brune ierr = SNESSetInitialFunction(next->snes,FSub);CHKERRQ(ierr); 77eed5f15bSPeter Brune } else { 78eed5f15bSPeter Brune next = next->next; 79eed5f15bSPeter Brune } 80eed5f15bSPeter Brune ierr = SNESSolve(next->snes,B,X);CHKERRQ(ierr); 81d5a53f06SPeter Brune ierr = SNESGetConvergedReason(next->snes,&reason);CHKERRQ(ierr); 82d5a53f06SPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 830b762f1fSPatrick Farrell jac->innerFailures++; 840b762f1fSPatrick Farrell if (jac->innerFailures >= snes->maxFailures) { 85d5a53f06SPeter Brune snes->reason = SNES_DIVERGED_INNER; 86d5a53f06SPeter Brune PetscFunctionReturn(0); 87d5a53f06SPeter Brune } 88eed5f15bSPeter Brune } 890b762f1fSPatrick Farrell } 90efd4aadfSBarry Smith if (next->snes->npcside== PC_RIGHT) { 91eed5f15bSPeter Brune ierr = SNESGetFunction(next->snes,&FSub,NULL,NULL);CHKERRQ(ierr); 92eed5f15bSPeter Brune ierr = VecCopy(FSub,F);CHKERRQ(ierr); 93d5a53f06SPeter Brune if (fnorm) { 9463cdc2bdSPatrick Farrell if (snes->xl && snes->xu) { 9563cdc2bdSPatrick Farrell ierr = SNESVIComputeInactiveSetFnorm(snes, F, X, fnorm);CHKERRQ(ierr); 9663cdc2bdSPatrick Farrell } else { 9771dbe336SPeter Brune ierr = VecNorm(F, NORM_2, fnorm);CHKERRQ(ierr); 9863cdc2bdSPatrick Farrell } 99422a814eSBarry Smith SNESCheckFunctionNorm(snes,*fnorm); 100d5a53f06SPeter Brune } 10172edecb9SPeter Brune } else if (snes->normschedule == SNES_NORM_ALWAYS) { 102be95d8f1SBarry Smith ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 103d5a53f06SPeter Brune if (fnorm) { 104a6da83ecSPatrick Farrell if (snes->xl && snes->xu) { 105a6da83ecSPatrick Farrell ierr = SNESVIComputeInactiveSetFnorm(snes, F, X, fnorm);CHKERRQ(ierr); 106a6da83ecSPatrick Farrell } else { 107d5a53f06SPeter Brune ierr = VecNorm(F, NORM_2, fnorm);CHKERRQ(ierr); 108a6da83ecSPatrick Farrell } 109422a814eSBarry Smith SNESCheckFunctionNorm(snes,*fnorm); 110d5a53f06SPeter Brune } 111eed5f15bSPeter Brune } 112eed5f15bSPeter Brune PetscFunctionReturn(0); 113eed5f15bSPeter Brune } 114eed5f15bSPeter Brune 115eed5f15bSPeter Brune static PetscErrorCode SNESCompositeApply_Additive(SNES snes,Vec X,Vec B,Vec F,PetscReal *fnorm) 116eed5f15bSPeter Brune { 117eed5f15bSPeter Brune PetscErrorCode ierr; 118eed5f15bSPeter Brune SNES_Composite *jac = (SNES_Composite*)snes->data; 119eed5f15bSPeter Brune SNES_CompositeLink next = jac->head; 120eed5f15bSPeter Brune Vec Y,Xorig; 121d5a53f06SPeter Brune SNESConvergedReason reason; 122eed5f15bSPeter Brune 123eed5f15bSPeter Brune PetscFunctionBegin; 124eed5f15bSPeter Brune Y = snes->vec_sol_update; 125eed5f15bSPeter Brune if (!jac->Xorig) {ierr = VecDuplicate(X,&jac->Xorig);CHKERRQ(ierr);} 126eed5f15bSPeter Brune Xorig = jac->Xorig; 127302440fdSBarry Smith ierr = VecCopy(X,Xorig);CHKERRQ(ierr); 128eed5f15bSPeter Brune if (!next) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"No composite SNESes supplied via SNESCompositeAddSNES() or -snes_composite_sneses"); 12972edecb9SPeter Brune if (snes->normschedule == SNES_NORM_ALWAYS) { 130eed5f15bSPeter Brune ierr = SNESSetInitialFunction(next->snes,F);CHKERRQ(ierr); 131eed5f15bSPeter Brune while (next->next) { 132eed5f15bSPeter Brune next = next->next; 133eed5f15bSPeter Brune ierr = SNESSetInitialFunction(next->snes,F);CHKERRQ(ierr); 134eed5f15bSPeter Brune } 135eed5f15bSPeter Brune } 136eed5f15bSPeter Brune next = jac->head; 1378f626970SPeter Brune ierr = VecCopy(Xorig,Y);CHKERRQ(ierr); 138eed5f15bSPeter Brune ierr = SNESSolve(next->snes,B,Y);CHKERRQ(ierr); 139d5a53f06SPeter Brune ierr = SNESGetConvergedReason(next->snes,&reason);CHKERRQ(ierr); 140d5a53f06SPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 1410b762f1fSPatrick Farrell jac->innerFailures++; 1420b762f1fSPatrick Farrell if (jac->innerFailures >= snes->maxFailures) { 143d5a53f06SPeter Brune snes->reason = SNES_DIVERGED_INNER; 144d5a53f06SPeter Brune PetscFunctionReturn(0); 145d5a53f06SPeter Brune } 1460b762f1fSPatrick Farrell } 147eed5f15bSPeter Brune ierr = VecAXPY(Y,-1.0,Xorig);CHKERRQ(ierr); 1488f626970SPeter Brune ierr = VecAXPY(X,next->dmp,Y);CHKERRQ(ierr); 1498f626970SPeter Brune while (next->next) { 1508f626970SPeter Brune next = next->next; 1518f626970SPeter Brune ierr = VecCopy(Xorig,Y);CHKERRQ(ierr); 1528f626970SPeter Brune ierr = SNESSolve(next->snes,B,Y);CHKERRQ(ierr); 153d5a53f06SPeter Brune ierr = SNESGetConvergedReason(next->snes,&reason);CHKERRQ(ierr); 154d5a53f06SPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 1550b762f1fSPatrick Farrell jac->innerFailures++; 1560b762f1fSPatrick Farrell if (jac->innerFailures >= snes->maxFailures) { 157d5a53f06SPeter Brune snes->reason = SNES_DIVERGED_INNER; 158d5a53f06SPeter Brune PetscFunctionReturn(0); 159d5a53f06SPeter Brune } 1600b762f1fSPatrick Farrell } 1618f626970SPeter Brune ierr = VecAXPY(Y,-1.0,Xorig);CHKERRQ(ierr); 1628f626970SPeter Brune ierr = VecAXPY(X,next->dmp,Y);CHKERRQ(ierr); 163eed5f15bSPeter Brune } 16472edecb9SPeter Brune if (snes->normschedule == SNES_NORM_ALWAYS) { 165eed5f15bSPeter Brune ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 166a6da83ecSPatrick Farrell if (fnorm) { 167a6da83ecSPatrick Farrell if (snes->xl && snes->xu) { 168a6da83ecSPatrick Farrell ierr = SNESVIComputeInactiveSetFnorm(snes, F, X, fnorm);CHKERRQ(ierr); 169a6da83ecSPatrick Farrell } else { 170a6da83ecSPatrick Farrell ierr = VecNorm(F, NORM_2, fnorm);CHKERRQ(ierr); 171a6da83ecSPatrick Farrell } 172422a814eSBarry Smith SNESCheckFunctionNorm(snes,*fnorm); 173a6da83ecSPatrick Farrell } 174eed5f15bSPeter Brune } 175eed5f15bSPeter Brune PetscFunctionReturn(0); 176eed5f15bSPeter Brune } 17790a8ba9bSPeter Brune 17890a8ba9bSPeter Brune /* approximately solve the overdetermined system: 17990a8ba9bSPeter Brune 18090a8ba9bSPeter Brune 2*F(x_i)\cdot F(\x_j)\alpha_i = 0 18190a8ba9bSPeter Brune \alpha_i = 1 18290a8ba9bSPeter Brune 18390a8ba9bSPeter Brune Which minimizes the L2 norm of the linearization of: 18490a8ba9bSPeter Brune ||F(\sum_i \alpha_i*x_i)||^2 18590a8ba9bSPeter Brune 18690a8ba9bSPeter Brune With the constraint that \sum_i\alpha_i = 1 18790a8ba9bSPeter Brune Where x_i is the solution from the ith subsolver. 18890a8ba9bSPeter Brune */ 18990a8ba9bSPeter Brune static PetscErrorCode SNESCompositeApply_AdditiveOptimal(SNES snes,Vec X,Vec B,Vec F,PetscReal *fnorm) 19090a8ba9bSPeter Brune { 19190a8ba9bSPeter Brune PetscErrorCode ierr; 19290a8ba9bSPeter Brune SNES_Composite *jac = (SNES_Composite*)snes->data; 19390a8ba9bSPeter Brune SNES_CompositeLink next = jac->head; 19490a8ba9bSPeter Brune Vec *Xes = jac->Xes,*Fes = jac->Fes; 19590a8ba9bSPeter Brune PetscInt i,j; 1965e806d2eSPeter Brune PetscScalar tot,total,ftf; 1975e806d2eSPeter Brune PetscReal min_fnorm; 1985e806d2eSPeter Brune PetscInt min_i; 199d5a53f06SPeter Brune SNESConvergedReason reason; 20090a8ba9bSPeter Brune 20190a8ba9bSPeter Brune PetscFunctionBegin; 20290a8ba9bSPeter Brune if (!next) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"No composite SNESes supplied via SNESCompositeAddSNES() or -snes_composite_sneses"); 20358bdfa14SPeter Brune 20458bdfa14SPeter Brune if (snes->normschedule == SNES_NORM_ALWAYS) { 20558bdfa14SPeter Brune next = jac->head; 20658bdfa14SPeter Brune ierr = SNESSetInitialFunction(next->snes,F);CHKERRQ(ierr); 20758bdfa14SPeter Brune while (next->next) { 20858bdfa14SPeter Brune next = next->next; 20958bdfa14SPeter Brune ierr = SNESSetInitialFunction(next->snes,F);CHKERRQ(ierr); 21058bdfa14SPeter Brune } 21158bdfa14SPeter Brune } 21258bdfa14SPeter Brune 21390a8ba9bSPeter Brune next = jac->head; 21490a8ba9bSPeter Brune i = 0; 21590a8ba9bSPeter Brune ierr = VecCopy(X,Xes[i]);CHKERRQ(ierr); 21690a8ba9bSPeter Brune ierr = SNESSolve(next->snes,B,Xes[i]);CHKERRQ(ierr); 217d5a53f06SPeter Brune ierr = SNESGetConvergedReason(next->snes,&reason);CHKERRQ(ierr); 218d5a53f06SPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 2190b762f1fSPatrick Farrell jac->innerFailures++; 2200b762f1fSPatrick Farrell if (jac->innerFailures >= snes->maxFailures) { 221d5a53f06SPeter Brune snes->reason = SNES_DIVERGED_INNER; 222d5a53f06SPeter Brune PetscFunctionReturn(0); 223d5a53f06SPeter Brune } 2240b762f1fSPatrick Farrell } 22590a8ba9bSPeter Brune while (next->next) { 22690a8ba9bSPeter Brune i++; 22790a8ba9bSPeter Brune next = next->next; 22890a8ba9bSPeter Brune ierr = VecCopy(X,Xes[i]);CHKERRQ(ierr); 22990a8ba9bSPeter Brune ierr = SNESSolve(next->snes,B,Xes[i]);CHKERRQ(ierr); 230d5a53f06SPeter Brune ierr = SNESGetConvergedReason(next->snes,&reason);CHKERRQ(ierr); 231d5a53f06SPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 2320b762f1fSPatrick Farrell jac->innerFailures++; 2330b762f1fSPatrick Farrell if (jac->innerFailures >= snes->maxFailures) { 234d5a53f06SPeter Brune snes->reason = SNES_DIVERGED_INNER; 235d5a53f06SPeter Brune PetscFunctionReturn(0); 236d5a53f06SPeter Brune } 23790a8ba9bSPeter Brune } 2380b762f1fSPatrick Farrell } 23990a8ba9bSPeter Brune 24090a8ba9bSPeter Brune /* all the solutions are collected; combine optimally */ 24190a8ba9bSPeter Brune for (i=0;i<jac->n;i++) { 24290a8ba9bSPeter Brune for (j=0;j<i+1;j++) { 2439c2f3473SPeter Brune ierr = VecDotBegin(Fes[i],Fes[j],&jac->h[i + j*jac->n]);CHKERRQ(ierr); 24490a8ba9bSPeter Brune } 2459c2f3473SPeter Brune ierr = VecDotBegin(Fes[i],F,&jac->g[i]);CHKERRQ(ierr); 24690a8ba9bSPeter Brune } 2475e806d2eSPeter Brune 24890a8ba9bSPeter Brune for (i=0;i<jac->n;i++) { 24990a8ba9bSPeter Brune for (j=0;j<i+1;j++) { 2509c2f3473SPeter Brune ierr = VecDotEnd(Fes[i],Fes[j],&jac->h[i + j*jac->n]);CHKERRQ(ierr); 2519c2f3473SPeter Brune if (i == j) jac->fnorms[i] = PetscSqrtReal(PetscRealPart(jac->h[i + j*jac->n])); 25290a8ba9bSPeter Brune } 2539c2f3473SPeter Brune ierr = VecDotEnd(Fes[i],F,&jac->g[i]);CHKERRQ(ierr); 25490a8ba9bSPeter Brune } 25590a8ba9bSPeter Brune 2569c2f3473SPeter Brune ftf = (*fnorm)*(*fnorm); 25790a8ba9bSPeter Brune 25890a8ba9bSPeter Brune for (i=0; i<jac->n; i++) { 25990a8ba9bSPeter Brune for (j=i+1;j<jac->n;j++) { 2609c2f3473SPeter Brune jac->h[i + j*jac->n] = jac->h[j + i*jac->n]; 26190a8ba9bSPeter Brune } 26290a8ba9bSPeter Brune } 26390a8ba9bSPeter Brune 26490a8ba9bSPeter Brune for (i=0; i<jac->n; i++) { 2659c2f3473SPeter Brune for (j=0; j<jac->n; j++) { 2669c2f3473SPeter Brune jac->h[i + j*jac->n] = jac->h[i + j*jac->n] - jac->g[j] - jac->g[i] + ftf; 26790a8ba9bSPeter Brune } 2689c2f3473SPeter Brune jac->beta[i] = ftf - jac->g[i]; 2699c2f3473SPeter Brune } 27090a8ba9bSPeter Brune 27190a8ba9bSPeter Brune jac->info = 0; 27290a8ba9bSPeter Brune jac->rcond = -1.; 27390a8ba9bSPeter Brune ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr); 27490a8ba9bSPeter Brune #if defined(PETSC_USE_COMPLEX) 27573cf7048SBarry Smith PetscStackCallBLAS("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)); 27690a8ba9bSPeter Brune #else 27773cf7048SBarry Smith PetscStackCallBLAS("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)); 27890a8ba9bSPeter Brune #endif 27990a8ba9bSPeter Brune ierr = PetscFPTrapPop();CHKERRQ(ierr); 28090a8ba9bSPeter Brune if (jac->info < 0) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_LIB,"Bad argument to GELSS"); 28190a8ba9bSPeter Brune if (jac->info > 0) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_LIB,"SVD failed to converge"); 2829c2f3473SPeter Brune tot = 0.; 2835e806d2eSPeter Brune total = 0.; 28490a8ba9bSPeter Brune for (i=0; i<jac->n; i++) { 285422a814eSBarry Smith if (snes->errorifnotconverged && PetscIsInfOrNanScalar(jac->beta[i])) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_LIB,"SVD generated inconsistent output"); 286c783eb9eSBarry Smith ierr = PetscInfo2(snes,"%D: %g\n",i,(double)PetscRealPart(jac->beta[i]));CHKERRQ(ierr); 2879c2f3473SPeter Brune tot += jac->beta[i]; 2885e806d2eSPeter Brune total += PetscAbsScalar(jac->beta[i]); 28990a8ba9bSPeter Brune } 2909c2f3473SPeter Brune ierr = VecScale(X,(1. - tot));CHKERRQ(ierr); 29190a8ba9bSPeter Brune ierr = VecMAXPY(X,jac->n,jac->beta,Xes);CHKERRQ(ierr); 29290a8ba9bSPeter Brune ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 293a6da83ecSPatrick Farrell 294a6da83ecSPatrick Farrell if (snes->xl && snes->xu) { 295a6da83ecSPatrick Farrell ierr = SNESVIComputeInactiveSetFnorm(snes, F, X, fnorm);CHKERRQ(ierr); 296a6da83ecSPatrick Farrell } else { 2979c2f3473SPeter Brune ierr = VecNorm(F, NORM_2, fnorm);CHKERRQ(ierr); 298a6da83ecSPatrick Farrell } 29990a8ba9bSPeter Brune 3005e806d2eSPeter Brune /* take the minimum-normed candidate if it beats the combination by a factor of rtol or the combination has stagnated */ 3015e806d2eSPeter Brune min_fnorm = jac->fnorms[0]; 3025e806d2eSPeter Brune min_i = 0; 3039c2f3473SPeter Brune for (i=0; i<jac->n; i++) { 3045e806d2eSPeter Brune if (jac->fnorms[i] < min_fnorm) { 3055e806d2eSPeter Brune min_fnorm = jac->fnorms[i]; 3065e806d2eSPeter Brune min_i = i; 3079c2f3473SPeter Brune } 3089c2f3473SPeter Brune } 3095e806d2eSPeter Brune 3105e806d2eSPeter Brune /* stagnation or divergence restart to the solution of the solver that failed the least */ 3115e806d2eSPeter Brune if (PetscRealPart(total) < jac->stol || min_fnorm*jac->rtol < *fnorm) { 31258bdfa14SPeter Brune ierr = VecCopy(jac->Xes[min_i],X);CHKERRQ(ierr); 31358bdfa14SPeter Brune ierr = VecCopy(jac->Fes[min_i],F);CHKERRQ(ierr); 3145e806d2eSPeter Brune *fnorm = min_fnorm; 3155e806d2eSPeter Brune } 31690a8ba9bSPeter Brune PetscFunctionReturn(0); 31790a8ba9bSPeter Brune } 31890a8ba9bSPeter Brune 319eed5f15bSPeter Brune static PetscErrorCode SNESSetUp_Composite(SNES snes) 320eed5f15bSPeter Brune { 321eed5f15bSPeter Brune PetscErrorCode ierr; 322dd515a93SPeter Brune DM dm; 323eed5f15bSPeter Brune SNES_Composite *jac = (SNES_Composite*)snes->data; 324eed5f15bSPeter Brune SNES_CompositeLink next = jac->head; 32590a8ba9bSPeter Brune PetscInt n=0,i; 32690a8ba9bSPeter Brune Vec F; 327eed5f15bSPeter Brune 328eed5f15bSPeter Brune PetscFunctionBegin; 329dd515a93SPeter Brune ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 33063cdc2bdSPatrick Farrell 33163cdc2bdSPatrick Farrell if (snes->ops->computevariablebounds) { 33263cdc2bdSPatrick Farrell /* SNESVI only ever calls computevariablebounds once, so calling it once here is justified */ 33363cdc2bdSPatrick Farrell if (!snes->xl) {ierr = VecDuplicate(snes->vec_sol,&snes->xl);CHKERRQ(ierr);} 33463cdc2bdSPatrick Farrell if (!snes->xu) {ierr = VecDuplicate(snes->vec_sol,&snes->xu);CHKERRQ(ierr);} 33563cdc2bdSPatrick Farrell ierr = (*snes->ops->computevariablebounds)(snes,snes->xl,snes->xu);CHKERRQ(ierr); 33663cdc2bdSPatrick Farrell } 33763cdc2bdSPatrick Farrell 338eed5f15bSPeter Brune while (next) { 33990a8ba9bSPeter Brune n++; 340dd515a93SPeter Brune ierr = SNESSetDM(next->snes,dm);CHKERRQ(ierr); 34163cdc2bdSPatrick Farrell ierr = SNESSetApplicationContext(next->snes, snes->user);CHKERRQ(ierr); 34263cdc2bdSPatrick Farrell if (snes->xl && snes->xu) { 34363cdc2bdSPatrick Farrell if (snes->ops->computevariablebounds) { 34463cdc2bdSPatrick Farrell ierr = SNESVISetComputeVariableBounds(next->snes, snes->ops->computevariablebounds);CHKERRQ(ierr); 34563cdc2bdSPatrick Farrell } else { 34663cdc2bdSPatrick Farrell ierr = SNESVISetVariableBounds(next->snes,snes->xl,snes->xu);CHKERRQ(ierr); 34763cdc2bdSPatrick Farrell } 34863cdc2bdSPatrick Farrell } 34963cdc2bdSPatrick Farrell 350eed5f15bSPeter Brune next = next->next; 351eed5f15bSPeter Brune } 35290a8ba9bSPeter Brune jac->nsnes = n; 35390a8ba9bSPeter Brune ierr = SNESGetFunction(snes,&F,NULL,NULL);CHKERRQ(ierr); 35490a8ba9bSPeter Brune if (jac->type == SNES_COMPOSITE_ADDITIVEOPTIMAL) { 35590a8ba9bSPeter Brune ierr = VecDuplicateVecs(F,jac->nsnes,&jac->Xes);CHKERRQ(ierr); 356854ce69bSBarry Smith ierr = PetscMalloc1(n,&jac->Fes);CHKERRQ(ierr); 357854ce69bSBarry Smith ierr = PetscMalloc1(n,&jac->fnorms);CHKERRQ(ierr); 35890a8ba9bSPeter Brune next = jac->head; 35990a8ba9bSPeter Brune i = 0; 36090a8ba9bSPeter Brune while (next) { 36190a8ba9bSPeter Brune ierr = SNESGetFunction(next->snes,&F,NULL,NULL);CHKERRQ(ierr); 36290a8ba9bSPeter Brune jac->Fes[i] = F; 36390a8ba9bSPeter Brune ierr = PetscObjectReference((PetscObject)F);CHKERRQ(ierr); 36490a8ba9bSPeter Brune next = next->next; 36590a8ba9bSPeter Brune i++; 36690a8ba9bSPeter Brune } 36790a8ba9bSPeter Brune /* allocate the subspace direct solve area */ 36890a8ba9bSPeter Brune jac->nrhs = 1; 3699c2f3473SPeter Brune jac->lda = jac->nsnes; 37090a8ba9bSPeter Brune jac->ldb = jac->nsnes; 37190a8ba9bSPeter Brune jac->n = jac->nsnes; 37290a8ba9bSPeter Brune 373785e854fSJed Brown ierr = PetscMalloc1(jac->n*jac->n,&jac->h);CHKERRQ(ierr); 374785e854fSJed Brown ierr = PetscMalloc1(jac->n,&jac->beta);CHKERRQ(ierr); 375785e854fSJed Brown ierr = PetscMalloc1(jac->n,&jac->s);CHKERRQ(ierr); 376785e854fSJed Brown ierr = PetscMalloc1(jac->n,&jac->g);CHKERRQ(ierr); 3779c2f3473SPeter Brune jac->lwork = 12*jac->n; 378dd63322aSSatish Balay #if defined(PETSC_USE_COMPLEX) 379854ce69bSBarry Smith ierr = PetscMalloc1(jac->lwork,&jac->rwork);CHKERRQ(ierr); 38090a8ba9bSPeter Brune #endif 381854ce69bSBarry Smith ierr = PetscMalloc1(jac->lwork,&jac->work);CHKERRQ(ierr); 38290a8ba9bSPeter Brune } 38390a8ba9bSPeter Brune 384eed5f15bSPeter Brune PetscFunctionReturn(0); 385eed5f15bSPeter Brune } 386eed5f15bSPeter Brune 387eed5f15bSPeter Brune static PetscErrorCode SNESReset_Composite(SNES snes) 388eed5f15bSPeter Brune { 389eed5f15bSPeter Brune SNES_Composite *jac = (SNES_Composite*)snes->data; 390eed5f15bSPeter Brune PetscErrorCode ierr; 391eed5f15bSPeter Brune SNES_CompositeLink next = jac->head; 392eed5f15bSPeter Brune 393eed5f15bSPeter Brune PetscFunctionBegin; 394eed5f15bSPeter Brune while (next) { 395eed5f15bSPeter Brune ierr = SNESReset(next->snes);CHKERRQ(ierr); 396eed5f15bSPeter Brune next = next->next; 397eed5f15bSPeter Brune } 398eed5f15bSPeter Brune ierr = VecDestroy(&jac->Xorig);CHKERRQ(ierr); 39990a8ba9bSPeter Brune if (jac->Xes) {ierr = VecDestroyVecs(jac->nsnes,&jac->Xes);CHKERRQ(ierr);} 40090a8ba9bSPeter Brune if (jac->Fes) {ierr = VecDestroyVecs(jac->nsnes,&jac->Fes);CHKERRQ(ierr);} 4019c2f3473SPeter Brune ierr = PetscFree(jac->fnorms);CHKERRQ(ierr); 40290a8ba9bSPeter Brune ierr = PetscFree(jac->h);CHKERRQ(ierr); 40390a8ba9bSPeter Brune ierr = PetscFree(jac->s);CHKERRQ(ierr); 4049c2f3473SPeter Brune ierr = PetscFree(jac->g);CHKERRQ(ierr); 40590a8ba9bSPeter Brune ierr = PetscFree(jac->beta);CHKERRQ(ierr); 40690a8ba9bSPeter Brune ierr = PetscFree(jac->work);CHKERRQ(ierr); 40790a8ba9bSPeter Brune ierr = PetscFree(jac->rwork);CHKERRQ(ierr); 408eed5f15bSPeter Brune PetscFunctionReturn(0); 409eed5f15bSPeter Brune } 410eed5f15bSPeter Brune 411eed5f15bSPeter Brune static PetscErrorCode SNESDestroy_Composite(SNES snes) 412eed5f15bSPeter Brune { 413eed5f15bSPeter Brune SNES_Composite *jac = (SNES_Composite*)snes->data; 414eed5f15bSPeter Brune PetscErrorCode ierr; 415eed5f15bSPeter Brune SNES_CompositeLink next = jac->head,next_tmp; 416eed5f15bSPeter Brune 417eed5f15bSPeter Brune PetscFunctionBegin; 418eed5f15bSPeter Brune ierr = SNESReset_Composite(snes);CHKERRQ(ierr); 419eed5f15bSPeter Brune while (next) { 420eed5f15bSPeter Brune ierr = SNESDestroy(&next->snes);CHKERRQ(ierr); 421eed5f15bSPeter Brune next_tmp = next; 422eed5f15bSPeter Brune next = next->next; 423eed5f15bSPeter Brune ierr = PetscFree(next_tmp);CHKERRQ(ierr); 424eed5f15bSPeter Brune } 425eed5f15bSPeter Brune ierr = PetscFree(snes->data);CHKERRQ(ierr); 426eed5f15bSPeter Brune PetscFunctionReturn(0); 427eed5f15bSPeter Brune } 428eed5f15bSPeter Brune 4294416b707SBarry Smith static PetscErrorCode SNESSetFromOptions_Composite(PetscOptionItems *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 static PetscErrorCode SNESView_Composite(SNES snes,PetscViewer viewer) 471eed5f15bSPeter Brune { 472eed5f15bSPeter Brune SNES_Composite *jac = (SNES_Composite*)snes->data; 473eed5f15bSPeter Brune PetscErrorCode ierr; 474eed5f15bSPeter Brune SNES_CompositeLink next = jac->head; 475eed5f15bSPeter Brune PetscBool iascii; 476eed5f15bSPeter Brune 477eed5f15bSPeter Brune PetscFunctionBegin; 478eed5f15bSPeter Brune ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 479eed5f15bSPeter Brune if (iascii) { 480efd4aadfSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," type - %s\n",SNESCompositeTypes[jac->type]);CHKERRQ(ierr); 481eed5f15bSPeter Brune ierr = PetscViewerASCIIPrintf(viewer," SNESes on composite preconditioner follow\n");CHKERRQ(ierr); 482eed5f15bSPeter Brune ierr = PetscViewerASCIIPrintf(viewer," ---------------------------------\n");CHKERRQ(ierr); 483eed5f15bSPeter Brune } 484eed5f15bSPeter Brune if (iascii) { 485eed5f15bSPeter Brune ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 486eed5f15bSPeter Brune } 487eed5f15bSPeter Brune while (next) { 488eed5f15bSPeter Brune ierr = SNESView(next->snes,viewer);CHKERRQ(ierr); 489eed5f15bSPeter Brune next = next->next; 490eed5f15bSPeter Brune } 491eed5f15bSPeter Brune if (iascii) { 492eed5f15bSPeter Brune ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 493eed5f15bSPeter Brune ierr = PetscViewerASCIIPrintf(viewer," ---------------------------------\n");CHKERRQ(ierr); 494eed5f15bSPeter Brune } 495eed5f15bSPeter Brune PetscFunctionReturn(0); 496eed5f15bSPeter Brune } 497eed5f15bSPeter Brune 498eed5f15bSPeter Brune /* ------------------------------------------------------------------------------*/ 499eed5f15bSPeter Brune 500eed5f15bSPeter Brune static PetscErrorCode SNESCompositeSetType_Composite(SNES snes,SNESCompositeType type) 501eed5f15bSPeter Brune { 502eed5f15bSPeter Brune SNES_Composite *jac = (SNES_Composite*)snes->data; 503eed5f15bSPeter Brune 504eed5f15bSPeter Brune PetscFunctionBegin; 505eed5f15bSPeter Brune jac->type = type; 506eed5f15bSPeter Brune PetscFunctionReturn(0); 507eed5f15bSPeter Brune } 508eed5f15bSPeter Brune 509eed5f15bSPeter Brune static PetscErrorCode SNESCompositeAddSNES_Composite(SNES snes,SNESType type) 510eed5f15bSPeter Brune { 511eed5f15bSPeter Brune SNES_Composite *jac; 512eed5f15bSPeter Brune SNES_CompositeLink next,ilink; 513eed5f15bSPeter Brune PetscErrorCode ierr; 514eed5f15bSPeter Brune PetscInt cnt = 0; 515eed5f15bSPeter Brune const char *prefix; 516d726e3a5SJed Brown char newprefix[20]; 517eed5f15bSPeter Brune DM dm; 518eed5f15bSPeter Brune 519eed5f15bSPeter Brune PetscFunctionBegin; 520b00a9115SJed Brown ierr = PetscNewLog(snes,&ilink);CHKERRQ(ierr); 5219e5d0892SLisandro Dalcin ilink->next = NULL; 522eed5f15bSPeter Brune ierr = SNESCreate(PetscObjectComm((PetscObject)snes),&ilink->snes);CHKERRQ(ierr); 523eed5f15bSPeter Brune ierr = PetscLogObjectParent((PetscObject)snes,(PetscObject)ilink->snes);CHKERRQ(ierr); 524eed5f15bSPeter Brune ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 525eed5f15bSPeter Brune ierr = SNESSetDM(ilink->snes,dm);CHKERRQ(ierr); 5264ae48641SStefano Zampini ierr = SNESSetTolerances(ilink->snes,snes->abstol,snes->rtol,snes->stol,1,snes->max_funcs);CHKERRQ(ierr); 5275a94abe6SSatish Balay ierr = PetscObjectCopyFortranFunctionPointers((PetscObject)snes,(PetscObject)ilink->snes);CHKERRQ(ierr); 528eed5f15bSPeter Brune jac = (SNES_Composite*)snes->data; 529eed5f15bSPeter Brune next = jac->head; 530eed5f15bSPeter Brune if (!next) { 531eed5f15bSPeter Brune jac->head = ilink; 532eed5f15bSPeter Brune ilink->previous = NULL; 533eed5f15bSPeter Brune } else { 534eed5f15bSPeter Brune cnt++; 535eed5f15bSPeter Brune while (next->next) { 536eed5f15bSPeter Brune next = next->next; 537eed5f15bSPeter Brune cnt++; 538eed5f15bSPeter Brune } 539eed5f15bSPeter Brune next->next = ilink; 540eed5f15bSPeter Brune ilink->previous = next; 541eed5f15bSPeter Brune } 542eed5f15bSPeter Brune ierr = SNESGetOptionsPrefix(snes,&prefix);CHKERRQ(ierr); 543eed5f15bSPeter Brune ierr = SNESSetOptionsPrefix(ilink->snes,prefix);CHKERRQ(ierr); 5444ae48641SStefano Zampini ierr = PetscSNPrintf(newprefix,sizeof(newprefix),"sub_%d_",(int)cnt);CHKERRQ(ierr); 545eed5f15bSPeter Brune ierr = SNESAppendOptionsPrefix(ilink->snes,newprefix);CHKERRQ(ierr); 5468f626970SPeter Brune ierr = PetscObjectIncrementTabLevel((PetscObject)ilink->snes,(PetscObject)snes,1);CHKERRQ(ierr); 547eed5f15bSPeter Brune ierr = SNESSetType(ilink->snes,type);CHKERRQ(ierr); 54872edecb9SPeter Brune ierr = SNESSetNormSchedule(ilink->snes, SNES_NORM_FINAL_ONLY);CHKERRQ(ierr); 54963cdc2bdSPatrick Farrell 5508f626970SPeter Brune ilink->dmp = 1.0; 55190a8ba9bSPeter Brune jac->nsnes++; 552eed5f15bSPeter Brune PetscFunctionReturn(0); 553eed5f15bSPeter Brune } 554eed5f15bSPeter Brune 555eed5f15bSPeter Brune static PetscErrorCode SNESCompositeGetSNES_Composite(SNES snes,PetscInt n,SNES *subsnes) 556eed5f15bSPeter Brune { 557eed5f15bSPeter Brune SNES_Composite *jac; 558eed5f15bSPeter Brune SNES_CompositeLink next; 559eed5f15bSPeter Brune PetscInt i; 560eed5f15bSPeter Brune 561eed5f15bSPeter Brune PetscFunctionBegin; 562eed5f15bSPeter Brune jac = (SNES_Composite*)snes->data; 563eed5f15bSPeter Brune next = jac->head; 564eed5f15bSPeter Brune for (i=0; i<n; i++) { 565eed5f15bSPeter Brune if (!next->next) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_INCOMP,"Not enough SNESes in composite preconditioner"); 566eed5f15bSPeter Brune next = next->next; 567eed5f15bSPeter Brune } 568eed5f15bSPeter Brune *subsnes = next->snes; 569eed5f15bSPeter Brune PetscFunctionReturn(0); 570eed5f15bSPeter Brune } 571eed5f15bSPeter Brune 572eed5f15bSPeter Brune /* -------------------------------------------------------------------------------- */ 573e31ab4f9SPeter Brune /*@C 574eed5f15bSPeter Brune SNESCompositeSetType - Sets the type of composite preconditioner. 575eed5f15bSPeter Brune 576eed5f15bSPeter Brune Logically Collective on SNES 577eed5f15bSPeter Brune 578*d8d19677SJose E. Roman Input Parameters: 579eed5f15bSPeter Brune + snes - the preconditioner context 580eed5f15bSPeter Brune - type - SNES_COMPOSITE_ADDITIVE (default), SNES_COMPOSITE_MULTIPLICATIVE 581eed5f15bSPeter Brune 582eed5f15bSPeter Brune Options Database Key: 583eed5f15bSPeter Brune . -snes_composite_type <type: one of multiplicative, additive, special> - Sets composite preconditioner type 584eed5f15bSPeter Brune 585eed5f15bSPeter Brune Level: Developer 586eed5f15bSPeter Brune 587eed5f15bSPeter Brune @*/ 588eed5f15bSPeter Brune PetscErrorCode SNESCompositeSetType(SNES snes,SNESCompositeType type) 589eed5f15bSPeter Brune { 590eed5f15bSPeter Brune PetscErrorCode ierr; 591eed5f15bSPeter Brune 592eed5f15bSPeter Brune PetscFunctionBegin; 593eed5f15bSPeter Brune PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 594eed5f15bSPeter Brune PetscValidLogicalCollectiveEnum(snes,type,2); 595eed5f15bSPeter Brune ierr = PetscTryMethod(snes,"SNESCompositeSetType_C",(SNES,SNESCompositeType),(snes,type));CHKERRQ(ierr); 596eed5f15bSPeter Brune PetscFunctionReturn(0); 597eed5f15bSPeter Brune } 598eed5f15bSPeter Brune 599eed5f15bSPeter Brune /*@C 600eed5f15bSPeter Brune SNESCompositeAddSNES - Adds another SNES to the composite SNES. 601eed5f15bSPeter Brune 602eed5f15bSPeter Brune Collective on SNES 603eed5f15bSPeter Brune 604eed5f15bSPeter Brune Input Parameters: 605eed5f15bSPeter Brune + snes - the preconditioner context 606eed5f15bSPeter Brune - type - the type of the new preconditioner 607eed5f15bSPeter Brune 608eed5f15bSPeter Brune Level: Developer 609eed5f15bSPeter Brune 610eed5f15bSPeter Brune @*/ 611eed5f15bSPeter Brune PetscErrorCode SNESCompositeAddSNES(SNES snes,SNESType type) 612eed5f15bSPeter Brune { 613eed5f15bSPeter Brune PetscErrorCode ierr; 614eed5f15bSPeter Brune 615eed5f15bSPeter Brune PetscFunctionBegin; 616eed5f15bSPeter Brune PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 617eed5f15bSPeter Brune ierr = PetscTryMethod(snes,"SNESCompositeAddSNES_C",(SNES,SNESType),(snes,type));CHKERRQ(ierr); 618eed5f15bSPeter Brune PetscFunctionReturn(0); 619eed5f15bSPeter Brune } 6204ae48641SStefano Zampini 621eed5f15bSPeter Brune /*@ 622eed5f15bSPeter Brune SNESCompositeGetSNES - Gets one of the SNES objects in the composite SNES. 623eed5f15bSPeter Brune 624eed5f15bSPeter Brune Not Collective 625eed5f15bSPeter Brune 626*d8d19677SJose E. Roman Input Parameters: 627eed5f15bSPeter Brune + snes - the preconditioner context 628eed5f15bSPeter Brune - n - the number of the snes requested 629eed5f15bSPeter Brune 630eed5f15bSPeter Brune Output Parameters: 631eed5f15bSPeter Brune . subsnes - the SNES requested 632eed5f15bSPeter Brune 633eed5f15bSPeter Brune Level: Developer 634eed5f15bSPeter Brune 635eed5f15bSPeter Brune .seealso: SNESCompositeAddSNES() 636eed5f15bSPeter Brune @*/ 637eed5f15bSPeter Brune PetscErrorCode SNESCompositeGetSNES(SNES snes,PetscInt n,SNES *subsnes) 638eed5f15bSPeter Brune { 639eed5f15bSPeter Brune PetscErrorCode ierr; 640eed5f15bSPeter Brune 641eed5f15bSPeter Brune PetscFunctionBegin; 642eed5f15bSPeter Brune PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 643eed5f15bSPeter Brune PetscValidPointer(subsnes,3); 644eed5f15bSPeter Brune ierr = PetscUseMethod(snes,"SNESCompositeGetSNES_C",(SNES,PetscInt,SNES*),(snes,n,subsnes));CHKERRQ(ierr); 645eed5f15bSPeter Brune PetscFunctionReturn(0); 646eed5f15bSPeter Brune } 647eed5f15bSPeter Brune 6486b2b7f7bSPatrick Farrell /*@ 6496b2b7f7bSPatrick Farrell SNESCompositeGetNumber - Get the number of subsolvers in the composite SNES. 6506b2b7f7bSPatrick Farrell 6516b2b7f7bSPatrick Farrell Logically Collective on SNES 6526b2b7f7bSPatrick Farrell 6536b2b7f7bSPatrick Farrell Input Parameter: 6546b2b7f7bSPatrick Farrell snes - the preconditioner context 6556b2b7f7bSPatrick Farrell 6566b2b7f7bSPatrick Farrell Output Parameter: 6576b2b7f7bSPatrick Farrell n - the number of subsolvers 6586b2b7f7bSPatrick Farrell 6596b2b7f7bSPatrick Farrell Level: Developer 6606b2b7f7bSPatrick Farrell 6616b2b7f7bSPatrick Farrell @*/ 6626b2b7f7bSPatrick Farrell PetscErrorCode SNESCompositeGetNumber(SNES snes,PetscInt *n) 6636b2b7f7bSPatrick Farrell { 6646b2b7f7bSPatrick Farrell SNES_Composite *jac; 6656b2b7f7bSPatrick Farrell SNES_CompositeLink next; 6666b2b7f7bSPatrick Farrell 6676b2b7f7bSPatrick Farrell PetscFunctionBegin; 6686b2b7f7bSPatrick Farrell jac = (SNES_Composite*)snes->data; 6696b2b7f7bSPatrick Farrell next = jac->head; 6706b2b7f7bSPatrick Farrell 6716b2b7f7bSPatrick Farrell *n = 0; 6726b2b7f7bSPatrick Farrell while (next) { 6736b2b7f7bSPatrick Farrell *n = *n + 1; 6746b2b7f7bSPatrick Farrell next = next->next; 6756b2b7f7bSPatrick Farrell } 6766b2b7f7bSPatrick Farrell PetscFunctionReturn(0); 6776b2b7f7bSPatrick Farrell } 6786b2b7f7bSPatrick Farrell 6798f626970SPeter Brune static PetscErrorCode SNESCompositeSetDamping_Composite(SNES snes,PetscInt n,PetscReal dmp) 6808f626970SPeter Brune { 6818f626970SPeter Brune SNES_Composite *jac; 6828f626970SPeter Brune SNES_CompositeLink next; 6838f626970SPeter Brune PetscInt i; 6848f626970SPeter Brune 6858f626970SPeter Brune PetscFunctionBegin; 6868f626970SPeter Brune jac = (SNES_Composite*)snes->data; 6878f626970SPeter Brune next = jac->head; 6888f626970SPeter Brune for (i=0; i<n; i++) { 6898f626970SPeter Brune if (!next->next) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_INCOMP,"Not enough SNESes in composite preconditioner"); 6908f626970SPeter Brune next = next->next; 6918f626970SPeter Brune } 6928f626970SPeter Brune next->dmp = dmp; 6938f626970SPeter Brune PetscFunctionReturn(0); 6948f626970SPeter Brune } 6958f626970SPeter Brune 6968f626970SPeter Brune /*@ 6978f626970SPeter Brune SNESCompositeSetDamping - Sets the damping of a subsolver when using additive composite SNES. 6988f626970SPeter Brune 6998f626970SPeter Brune Not Collective 7008f626970SPeter Brune 701*d8d19677SJose E. Roman Input Parameters: 7028f626970SPeter Brune + snes - the preconditioner context 7038f626970SPeter Brune . n - the number of the snes requested 7048f626970SPeter Brune - dmp - the damping 7058f626970SPeter Brune 7068f626970SPeter Brune Level: Developer 7078f626970SPeter Brune 7088f626970SPeter Brune .seealso: SNESCompositeAddSNES() 7098f626970SPeter Brune @*/ 7108f626970SPeter Brune PetscErrorCode SNESCompositeSetDamping(SNES snes,PetscInt n,PetscReal dmp) 7118f626970SPeter Brune { 7128f626970SPeter Brune PetscErrorCode ierr; 7138f626970SPeter Brune 7148f626970SPeter Brune PetscFunctionBegin; 7158f626970SPeter Brune PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 7168f626970SPeter Brune ierr = PetscUseMethod(snes,"SNESCompositeSetDamping_C",(SNES,PetscInt,PetscReal),(snes,n,dmp));CHKERRQ(ierr); 7178f626970SPeter Brune PetscFunctionReturn(0); 7188f626970SPeter Brune } 7198f626970SPeter Brune 72040244768SBarry Smith static PetscErrorCode SNESSolve_Composite(SNES snes) 721eed5f15bSPeter Brune { 7224ae48641SStefano Zampini Vec F,X,B,Y; 723eed5f15bSPeter Brune PetscInt i; 724b22975d2SPatrick Farrell PetscReal fnorm = 0.0, xnorm = 0.0, snorm = 0.0; 725eed5f15bSPeter Brune PetscErrorCode ierr; 72672edecb9SPeter Brune SNESNormSchedule normtype; 727eed5f15bSPeter Brune SNES_Composite *comp = (SNES_Composite*)snes->data; 728eed5f15bSPeter Brune 729eed5f15bSPeter Brune PetscFunctionBegin; 730eed5f15bSPeter Brune X = snes->vec_sol; 731eed5f15bSPeter Brune F = snes->vec_func; 732eed5f15bSPeter Brune B = snes->vec_rhs; 7334ae48641SStefano Zampini Y = snes->vec_sol_update; 734eed5f15bSPeter Brune 735e04113cfSBarry Smith ierr = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr); 736eed5f15bSPeter Brune snes->iter = 0; 737eed5f15bSPeter Brune snes->norm = 0.; 7380b762f1fSPatrick Farrell comp->innerFailures = 0; 739e04113cfSBarry Smith ierr = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr); 740eed5f15bSPeter Brune snes->reason = SNES_CONVERGED_ITERATING; 74172edecb9SPeter Brune ierr = SNESGetNormSchedule(snes, &normtype);CHKERRQ(ierr); 74272edecb9SPeter Brune if (normtype == SNES_NORM_ALWAYS || normtype == SNES_NORM_INITIAL_ONLY || normtype == SNES_NORM_INITIAL_FINAL_ONLY) { 743eed5f15bSPeter Brune if (!snes->vec_func_init_set) { 744eed5f15bSPeter Brune ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 745eed5f15bSPeter Brune } else snes->vec_func_init_set = PETSC_FALSE; 746eed5f15bSPeter Brune 747a6da83ecSPatrick Farrell if (snes->xl && snes->xu) { 748a6da83ecSPatrick Farrell ierr = SNESVIComputeInactiveSetFnorm(snes, F, X, &fnorm);CHKERRQ(ierr); 749a6da83ecSPatrick Farrell } else { 750eed5f15bSPeter Brune ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F|| */ 751a6da83ecSPatrick Farrell } 752422a814eSBarry Smith SNESCheckFunctionNorm(snes,fnorm); 753e04113cfSBarry Smith ierr = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr); 754eed5f15bSPeter Brune snes->iter = 0; 755eed5f15bSPeter Brune snes->norm = fnorm; 756e04113cfSBarry Smith ierr = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr); 757eed5f15bSPeter Brune ierr = SNESLogConvergenceHistory(snes,snes->norm,0);CHKERRQ(ierr); 758eed5f15bSPeter Brune ierr = SNESMonitor(snes,0,snes->norm);CHKERRQ(ierr); 759eed5f15bSPeter Brune 760eed5f15bSPeter Brune /* test convergence */ 761eed5f15bSPeter Brune ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 762eed5f15bSPeter Brune if (snes->reason) PetscFunctionReturn(0); 763eed5f15bSPeter Brune } else { 764e04113cfSBarry Smith ierr = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr); 765eed5f15bSPeter Brune ierr = SNESLogConvergenceHistory(snes,snes->norm,0);CHKERRQ(ierr); 766eed5f15bSPeter Brune ierr = SNESMonitor(snes,0,snes->norm);CHKERRQ(ierr); 767eed5f15bSPeter Brune } 768eed5f15bSPeter Brune 7694ae48641SStefano Zampini for (i = 0; i < snes->max_its; i++) { 770eed5f15bSPeter Brune /* Call general purpose update function */ 771eed5f15bSPeter Brune if (snes->ops->update) { 772eed5f15bSPeter Brune ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 773eed5f15bSPeter Brune } 774eed5f15bSPeter Brune 775b22975d2SPatrick Farrell /* Copy the state before modification by application of the composite solver; 776b22975d2SPatrick Farrell we will subtract the new state after application */ 7774ae48641SStefano Zampini ierr = VecCopy(X, Y);CHKERRQ(ierr); 778b22975d2SPatrick Farrell 779eed5f15bSPeter Brune if (comp->type == SNES_COMPOSITE_ADDITIVE) { 780eed5f15bSPeter Brune ierr = SNESCompositeApply_Additive(snes,X,B,F,&fnorm);CHKERRQ(ierr); 781eed5f15bSPeter Brune } else if (comp->type == SNES_COMPOSITE_MULTIPLICATIVE) { 782eed5f15bSPeter Brune ierr = SNESCompositeApply_Multiplicative(snes,X,B,F,&fnorm);CHKERRQ(ierr); 78390a8ba9bSPeter Brune } else if (comp->type == SNES_COMPOSITE_ADDITIVEOPTIMAL) { 78490a8ba9bSPeter Brune ierr = SNESCompositeApply_AdditiveOptimal(snes,X,B,F,&fnorm);CHKERRQ(ierr); 7856c4ed002SBarry Smith } else SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"Unsupported SNESComposite type"); 786d5a53f06SPeter Brune if (snes->reason < 0) break; 787d5a53f06SPeter Brune 788b22975d2SPatrick Farrell /* Compute the solution update for convergence testing */ 7894ae48641SStefano Zampini ierr = VecAYPX(Y, -1.0, X);CHKERRQ(ierr); 790b22975d2SPatrick Farrell 791eed5f15bSPeter Brune if ((i == snes->max_its - 1) && (normtype == SNES_NORM_INITIAL_FINAL_ONLY || normtype == SNES_NORM_FINAL_ONLY)) { 792eed5f15bSPeter Brune ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 793b22975d2SPatrick Farrell 794a6da83ecSPatrick Farrell if (snes->xl && snes->xu) { 795a6da83ecSPatrick Farrell ierr = VecNormBegin(X, NORM_2, &xnorm);CHKERRQ(ierr); 7964ae48641SStefano Zampini ierr = VecNormBegin(Y, NORM_2, &snorm);CHKERRQ(ierr); 797a6da83ecSPatrick Farrell ierr = SNESVIComputeInactiveSetFnorm(snes, F, X, &fnorm);CHKERRQ(ierr); 798a6da83ecSPatrick Farrell ierr = VecNormEnd(X, NORM_2, &xnorm);CHKERRQ(ierr); 7994ae48641SStefano Zampini ierr = VecNormEnd(Y, NORM_2, &snorm);CHKERRQ(ierr); 800a6da83ecSPatrick Farrell } else { 801b22975d2SPatrick Farrell ierr = VecNormBegin(F, NORM_2, &fnorm);CHKERRQ(ierr); 802b22975d2SPatrick Farrell ierr = VecNormBegin(X, NORM_2, &xnorm);CHKERRQ(ierr); 8034ae48641SStefano Zampini ierr = VecNormBegin(Y, NORM_2, &snorm);CHKERRQ(ierr); 804b22975d2SPatrick Farrell 805b22975d2SPatrick Farrell ierr = VecNormEnd(F, NORM_2, &fnorm);CHKERRQ(ierr); 806b22975d2SPatrick Farrell ierr = VecNormEnd(X, NORM_2, &xnorm);CHKERRQ(ierr); 8074ae48641SStefano Zampini ierr = VecNormEnd(Y, NORM_2, &snorm);CHKERRQ(ierr); 808a6da83ecSPatrick Farrell } 809422a814eSBarry Smith SNESCheckFunctionNorm(snes,fnorm); 810b22975d2SPatrick Farrell } else if (normtype == SNES_NORM_ALWAYS) { 811b22975d2SPatrick Farrell ierr = VecNormBegin(X, NORM_2, &xnorm);CHKERRQ(ierr); 8124ae48641SStefano Zampini ierr = VecNormBegin(Y, NORM_2, &snorm);CHKERRQ(ierr); 813b22975d2SPatrick Farrell ierr = VecNormEnd(X, NORM_2, &xnorm);CHKERRQ(ierr); 8144ae48641SStefano Zampini ierr = VecNormEnd(Y, NORM_2, &snorm);CHKERRQ(ierr); 815eed5f15bSPeter Brune } 816eed5f15bSPeter Brune /* Monitor convergence */ 817e04113cfSBarry Smith ierr = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr); 818eed5f15bSPeter Brune snes->iter = i+1; 819eed5f15bSPeter Brune snes->norm = fnorm; 820c1e67a49SFande Kong snes->xnorm = xnorm; 821c1e67a49SFande Kong snes->ynorm = snorm; 822e04113cfSBarry Smith ierr = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr); 823eed5f15bSPeter Brune ierr = SNESLogConvergenceHistory(snes,snes->norm,0);CHKERRQ(ierr); 824eed5f15bSPeter Brune ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr); 825eed5f15bSPeter Brune /* Test for convergence */ 826b22975d2SPatrick Farrell if (normtype == SNES_NORM_ALWAYS) {ierr = (*snes->ops->converged)(snes,snes->iter,xnorm,snorm,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);} 827eed5f15bSPeter Brune if (snes->reason) break; 828eed5f15bSPeter Brune } 82972edecb9SPeter Brune if (normtype == SNES_NORM_ALWAYS) { 830eed5f15bSPeter Brune if (i == snes->max_its) { 831eed5f15bSPeter Brune ierr = PetscInfo1(snes,"Maximum number of iterations has been reached: %D\n",snes->max_its);CHKERRQ(ierr); 832eed5f15bSPeter Brune if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT; 833eed5f15bSPeter Brune } 834eed5f15bSPeter Brune } else if (!snes->reason) snes->reason = SNES_CONVERGED_ITS; 835eed5f15bSPeter Brune PetscFunctionReturn(0); 836eed5f15bSPeter Brune } 837eed5f15bSPeter Brune 838eed5f15bSPeter Brune /* -------------------------------------------------------------------------------------------*/ 839eed5f15bSPeter Brune 840eed5f15bSPeter Brune /*MC 841eed5f15bSPeter Brune SNESCOMPOSITE - Build a preconditioner by composing together several nonlinear solvers 842eed5f15bSPeter Brune 843eed5f15bSPeter Brune Options Database Keys: 844eed5f15bSPeter Brune + -snes_composite_type <type: one of multiplicative, additive, symmetric_multiplicative, special> - Sets composite preconditioner type 845eed5f15bSPeter Brune - -snes_composite_sneses - <snes0,snes1,...> list of SNESes to compose 846eed5f15bSPeter Brune 847eed5f15bSPeter Brune Level: intermediate 848eed5f15bSPeter Brune 849eed5f15bSPeter Brune .seealso: SNESCreate(), SNESSetType(), SNESType (for list of available types), SNES, 850eed5f15bSPeter Brune SNESSHELL, SNESCompositeSetType(), SNESCompositeSpecialSetAlpha(), SNESCompositeAddSNES(), 851eed5f15bSPeter Brune SNESCompositeGetSNES() 852eed5f15bSPeter Brune 8534f02bc6aSBarry Smith References: 85496a0c994SBarry Smith . 1. - Peter R. Brune, Matthew G. Knepley, Barry F. Smith, and Xuemin Tu, "Composing Scalable Nonlinear Algebraic Solvers", 8554f02bc6aSBarry Smith SIAM Review, 57(4), 2015 8564f02bc6aSBarry Smith 857eed5f15bSPeter Brune M*/ 858eed5f15bSPeter Brune 859eed5f15bSPeter Brune PETSC_EXTERN PetscErrorCode SNESCreate_Composite(SNES snes) 860eed5f15bSPeter Brune { 861eed5f15bSPeter Brune PetscErrorCode ierr; 862eed5f15bSPeter Brune SNES_Composite *jac; 863eed5f15bSPeter Brune 864eed5f15bSPeter Brune PetscFunctionBegin; 865b00a9115SJed Brown ierr = PetscNewLog(snes,&jac);CHKERRQ(ierr); 866eed5f15bSPeter Brune 867eed5f15bSPeter Brune snes->ops->solve = SNESSolve_Composite; 868eed5f15bSPeter Brune snes->ops->setup = SNESSetUp_Composite; 869eed5f15bSPeter Brune snes->ops->reset = SNESReset_Composite; 870eed5f15bSPeter Brune snes->ops->destroy = SNESDestroy_Composite; 871eed5f15bSPeter Brune snes->ops->setfromoptions = SNESSetFromOptions_Composite; 872eed5f15bSPeter Brune snes->ops->view = SNESView_Composite; 873eed5f15bSPeter Brune 874d8d34be6SBarry Smith snes->usesksp = PETSC_FALSE; 875d8d34be6SBarry Smith 8764fc747eaSLawrence Mitchell snes->alwayscomputesfinalresidual = PETSC_FALSE; 8774fc747eaSLawrence Mitchell 878eed5f15bSPeter Brune snes->data = (void*)jac; 87990a8ba9bSPeter Brune jac->type = SNES_COMPOSITE_ADDITIVEOPTIMAL; 88090a8ba9bSPeter Brune jac->Fes = NULL; 88190a8ba9bSPeter Brune jac->Xes = NULL; 8829c2f3473SPeter Brune jac->fnorms = NULL; 88390a8ba9bSPeter Brune jac->nsnes = 0; 8849e5d0892SLisandro Dalcin jac->head = NULL; 8855e806d2eSPeter Brune jac->stol = 0.1; 8865e806d2eSPeter Brune jac->rtol = 1.1; 887eed5f15bSPeter Brune 88890a8ba9bSPeter Brune jac->h = NULL; 88990a8ba9bSPeter Brune jac->s = NULL; 89090a8ba9bSPeter Brune jac->beta = NULL; 89190a8ba9bSPeter Brune jac->work = NULL; 89290a8ba9bSPeter Brune jac->rwork = NULL; 89390a8ba9bSPeter Brune 8948f626970SPeter Brune ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESCompositeSetType_C",SNESCompositeSetType_Composite);CHKERRQ(ierr); 8958f626970SPeter Brune ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESCompositeAddSNES_C",SNESCompositeAddSNES_Composite);CHKERRQ(ierr); 8968f626970SPeter Brune ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESCompositeGetSNES_C",SNESCompositeGetSNES_Composite);CHKERRQ(ierr); 8978f626970SPeter Brune ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESCompositeSetDamping_C",SNESCompositeSetDamping_Composite);CHKERRQ(ierr); 898eed5f15bSPeter Brune PetscFunctionReturn(0); 899eed5f15bSPeter Brune } 900