1eed5f15bSPeter Brune 2eed5f15bSPeter Brune /* 3eed5f15bSPeter Brune Defines a SNES that can consist of a collection of SNESes 4eed5f15bSPeter Brune */ 5eed5f15bSPeter Brune #include <petsc-private/snesimpl.h> /*I "petscsnes.h" I*/ 690a8ba9bSPeter Brune #include <petscblaslapack.h> 7eed5f15bSPeter Brune 890a8ba9bSPeter Brune const char *const SNESCompositeTypes[] = {"ADDITIVE","MULTIPLICATIVE","ADDITIVEOPTIMAL","SNESCompositeType","SNES_COMPOSITE",0}; 9eed5f15bSPeter Brune 10eed5f15bSPeter Brune typedef struct _SNES_CompositeLink *SNES_CompositeLink; 11eed5f15bSPeter Brune struct _SNES_CompositeLink { 12eed5f15bSPeter Brune SNES snes; 138f626970SPeter Brune PetscReal dmp; 1490a8ba9bSPeter Brune Vec X; 15eed5f15bSPeter Brune SNES_CompositeLink next; 16eed5f15bSPeter Brune SNES_CompositeLink previous; 17eed5f15bSPeter Brune }; 18eed5f15bSPeter Brune 19eed5f15bSPeter Brune typedef struct { 20eed5f15bSPeter Brune SNES_CompositeLink head; 2190a8ba9bSPeter Brune PetscInt nsnes; 22eed5f15bSPeter Brune SNESCompositeType type; 23eed5f15bSPeter Brune Vec Xorig; 2490a8ba9bSPeter Brune 2590a8ba9bSPeter Brune /* context for ADDITIVEOPTIMAL */ 2690a8ba9bSPeter Brune Vec *Xes,*Fes; /* solution and residual vectors for the subsolvers */ 279c2f3473SPeter Brune PetscReal *fnorms; /* norms of the solutions */ 2890a8ba9bSPeter Brune PetscScalar *h; /* the matrix formed as q_ij = (rdot_i, rdot_j) */ 299c2f3473SPeter Brune PetscScalar *g; /* the dotproducts of the previous function with the candidate functions */ 309c2f3473SPeter Brune PetscBLASInt n; /* matrix dimension -- nsnes */ 3190a8ba9bSPeter Brune PetscBLASInt nrhs; /* the number of right hand sides */ 3290a8ba9bSPeter Brune PetscBLASInt lda; /* the padded matrix dimension */ 3390a8ba9bSPeter Brune PetscBLASInt ldb; /* the padded vector dimension */ 3490a8ba9bSPeter Brune PetscReal *s; /* the singular values */ 359c2f3473SPeter Brune PetscScalar *beta; /* the RHS and combination */ 3690a8ba9bSPeter Brune PetscReal rcond; /* the exit condition */ 3790a8ba9bSPeter Brune PetscBLASInt rank; /* the effective rank */ 3890a8ba9bSPeter Brune PetscScalar *work; /* the work vector */ 3990a8ba9bSPeter Brune PetscReal *rwork; /* the real work vector used for complex */ 4090a8ba9bSPeter Brune PetscBLASInt lwork; /* the size of the work vector */ 4190a8ba9bSPeter Brune PetscBLASInt info; /* the output condition */ 4290a8ba9bSPeter Brune 435e806d2eSPeter Brune PetscReal rtol; /* restart tolerance for accepting the combination */ 445e806d2eSPeter Brune PetscReal stol; /* restart tolerance for the combination */ 45eed5f15bSPeter Brune } SNES_Composite; 46eed5f15bSPeter Brune 47eed5f15bSPeter Brune #undef __FUNCT__ 48eed5f15bSPeter Brune #define __FUNCT__ "SNESCompositeApply_Multiplicative" 49eed5f15bSPeter Brune static PetscErrorCode SNESCompositeApply_Multiplicative(SNES snes,Vec X,Vec B,Vec F,PetscReal *fnorm) 50eed5f15bSPeter Brune { 51eed5f15bSPeter Brune PetscErrorCode ierr; 52eed5f15bSPeter Brune SNES_Composite *jac = (SNES_Composite*)snes->data; 53eed5f15bSPeter Brune SNES_CompositeLink next = jac->head; 54eed5f15bSPeter Brune Vec FSub; 55eed5f15bSPeter Brune PetscReal fsubnorm; 56eed5f15bSPeter Brune 57eed5f15bSPeter Brune PetscFunctionBegin; 58eed5f15bSPeter Brune if (!next) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"No composite SNESes supplied via SNESCompositeAddSNES() or -snes_composite_sneses"); 59eed5f15bSPeter Brune if (snes->normtype == SNES_NORM_FUNCTION) { 60eed5f15bSPeter Brune ierr = SNESSetInitialFunction(next->snes,F);CHKERRQ(ierr); 61eed5f15bSPeter Brune if (fnorm) {ierr = SNESSetInitialFunctionNorm(next->snes,*fnorm);CHKERRQ(ierr);} 62eed5f15bSPeter Brune } 63eed5f15bSPeter Brune ierr = SNESSolve(next->snes,B,X);CHKERRQ(ierr); 64eed5f15bSPeter Brune 65eed5f15bSPeter Brune while (next->next) { 66eed5f15bSPeter Brune /* only copy the function over in the case where the functions correspond */ 67eed5f15bSPeter Brune if (next->snes->pcside == PC_RIGHT && next->snes->normtype != SNES_NORM_NONE) { 68eed5f15bSPeter Brune ierr = SNESGetFunction(next->snes,&FSub,NULL,NULL);CHKERRQ(ierr); 69eed5f15bSPeter Brune ierr = SNESGetFunctionNorm(next->snes,&fsubnorm);CHKERRQ(ierr); 70eed5f15bSPeter Brune next = next->next; 71eed5f15bSPeter Brune ierr = SNESSetInitialFunction(next->snes,FSub);CHKERRQ(ierr); 72eed5f15bSPeter Brune ierr = SNESSetInitialFunctionNorm(next->snes,fsubnorm);CHKERRQ(ierr); 73eed5f15bSPeter Brune } else { 74eed5f15bSPeter Brune next = next->next; 75eed5f15bSPeter Brune } 76eed5f15bSPeter Brune ierr = SNESSolve(next->snes,B,X);CHKERRQ(ierr); 77eed5f15bSPeter Brune } 78eed5f15bSPeter Brune if (next->snes->pcside == PC_RIGHT) { 79eed5f15bSPeter Brune ierr = SNESGetFunction(next->snes,&FSub,NULL,NULL);CHKERRQ(ierr); 80eed5f15bSPeter Brune ierr = VecCopy(FSub,F);CHKERRQ(ierr); 81eed5f15bSPeter Brune if (fnorm) {ierr = SNESGetFunctionNorm(next->snes,fnorm);CHKERRQ(ierr);} 82eed5f15bSPeter Brune } else if (snes->normtype == SNES_NORM_FUNCTION) { 83eed5f15bSPeter Brune SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 84eed5f15bSPeter Brune if (fnorm) {ierr = VecNorm(F,NORM_2,fnorm);CHKERRQ(ierr);} 85eed5f15bSPeter Brune } 86eed5f15bSPeter Brune PetscFunctionReturn(0); 87eed5f15bSPeter Brune } 88eed5f15bSPeter Brune 89eed5f15bSPeter Brune #undef __FUNCT__ 90eed5f15bSPeter Brune #define __FUNCT__ "SNESCompositeApply_Additive" 91eed5f15bSPeter Brune static PetscErrorCode SNESCompositeApply_Additive(SNES snes,Vec X,Vec B,Vec F,PetscReal *fnorm) 92eed5f15bSPeter Brune { 93eed5f15bSPeter Brune PetscErrorCode ierr; 94eed5f15bSPeter Brune SNES_Composite *jac = (SNES_Composite*)snes->data; 95eed5f15bSPeter Brune SNES_CompositeLink next = jac->head; 96eed5f15bSPeter Brune Vec Y,Xorig; 97eed5f15bSPeter Brune 98eed5f15bSPeter Brune PetscFunctionBegin; 99eed5f15bSPeter Brune Y = snes->vec_sol_update; 100eed5f15bSPeter Brune if (!jac->Xorig) {ierr = VecDuplicate(X,&jac->Xorig);CHKERRQ(ierr);} 101eed5f15bSPeter Brune Xorig = jac->Xorig; 102eed5f15bSPeter Brune ierr = VecCopy(X,Xorig); 103eed5f15bSPeter Brune if (!next) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"No composite SNESes supplied via SNESCompositeAddSNES() or -snes_composite_sneses"); 104eed5f15bSPeter Brune if (snes->normtype == SNES_NORM_FUNCTION) { 105eed5f15bSPeter Brune ierr = SNESSetInitialFunction(next->snes,F);CHKERRQ(ierr); 106eed5f15bSPeter Brune if (fnorm) {ierr = SNESSetInitialFunctionNorm(next->snes,*fnorm);CHKERRQ(ierr);} 107eed5f15bSPeter Brune while (next->next) { 108eed5f15bSPeter Brune next = next->next; 109eed5f15bSPeter Brune ierr = SNESSetInitialFunction(next->snes,F);CHKERRQ(ierr); 110eed5f15bSPeter Brune if (fnorm) {ierr = SNESSetInitialFunctionNorm(next->snes,*fnorm);CHKERRQ(ierr);} 111eed5f15bSPeter Brune } 112eed5f15bSPeter Brune } 113eed5f15bSPeter Brune next = jac->head; 1148f626970SPeter Brune ierr = VecCopy(Xorig,Y);CHKERRQ(ierr); 115eed5f15bSPeter Brune ierr = SNESSolve(next->snes,B,Y);CHKERRQ(ierr); 116eed5f15bSPeter Brune ierr = VecAXPY(Y,-1.0,Xorig);CHKERRQ(ierr); 1178f626970SPeter Brune ierr = VecAXPY(X,next->dmp,Y);CHKERRQ(ierr); 1188f626970SPeter Brune while (next->next) { 1198f626970SPeter Brune next = next->next; 1208f626970SPeter Brune ierr = VecCopy(Xorig,Y);CHKERRQ(ierr); 1218f626970SPeter Brune ierr = SNESSolve(next->snes,B,Y);CHKERRQ(ierr); 1228f626970SPeter Brune ierr = VecAXPY(Y,-1.0,Xorig);CHKERRQ(ierr); 1238f626970SPeter Brune ierr = VecAXPY(X,next->dmp,Y);CHKERRQ(ierr); 124eed5f15bSPeter Brune } 125eed5f15bSPeter Brune if (snes->normtype == SNES_NORM_FUNCTION) { 126eed5f15bSPeter Brune ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 127eed5f15bSPeter Brune if (fnorm) {ierr = VecNorm(F,NORM_2,fnorm);CHKERRQ(ierr);} 128eed5f15bSPeter Brune } 129eed5f15bSPeter Brune PetscFunctionReturn(0); 130eed5f15bSPeter Brune } 13190a8ba9bSPeter Brune 13290a8ba9bSPeter Brune #undef __FUNCT__ 13390a8ba9bSPeter Brune #define __FUNCT__ "SNESCompositeApply_AdditiveOptimal" 13490a8ba9bSPeter Brune /* approximately solve the overdetermined system: 13590a8ba9bSPeter Brune 13690a8ba9bSPeter Brune 2*F(x_i)\cdot F(\x_j)\alpha_i = 0 13790a8ba9bSPeter Brune \alpha_i = 1 13890a8ba9bSPeter Brune 13990a8ba9bSPeter Brune Which minimizes the L2 norm of the linearization of: 14090a8ba9bSPeter Brune ||F(\sum_i \alpha_i*x_i)||^2 14190a8ba9bSPeter Brune 14290a8ba9bSPeter Brune With the constraint that \sum_i\alpha_i = 1 14390a8ba9bSPeter Brune Where x_i is the solution from the ith subsolver. 14490a8ba9bSPeter Brune */ 14590a8ba9bSPeter Brune static PetscErrorCode SNESCompositeApply_AdditiveOptimal(SNES snes,Vec X,Vec B,Vec F,PetscReal *fnorm) 14690a8ba9bSPeter Brune { 14790a8ba9bSPeter Brune PetscErrorCode ierr; 14890a8ba9bSPeter Brune SNES_Composite *jac = (SNES_Composite*)snes->data; 14990a8ba9bSPeter Brune SNES_CompositeLink next = jac->head; 15090a8ba9bSPeter Brune Vec *Xes = jac->Xes,*Fes = jac->Fes; 15190a8ba9bSPeter Brune PetscInt i,j; 1525e806d2eSPeter Brune PetscScalar tot,total,ftf; 1535e806d2eSPeter Brune PetscReal min_fnorm; 1545e806d2eSPeter Brune PetscInt min_i; 15590a8ba9bSPeter Brune 15690a8ba9bSPeter Brune PetscFunctionBegin; 15790a8ba9bSPeter Brune if (!next) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"No composite SNESes supplied via SNESCompositeAddSNES() or -snes_composite_sneses"); 15890a8ba9bSPeter Brune next = jac->head; 15990a8ba9bSPeter Brune i = 0; 16090a8ba9bSPeter Brune ierr = VecCopy(X,Xes[i]);CHKERRQ(ierr); 16190a8ba9bSPeter Brune ierr = SNESSolve(next->snes,B,Xes[i]);CHKERRQ(ierr); 16290a8ba9bSPeter Brune while (next->next) { 16390a8ba9bSPeter Brune i++; 16490a8ba9bSPeter Brune next = next->next; 16590a8ba9bSPeter Brune ierr = VecCopy(X,Xes[i]);CHKERRQ(ierr); 16690a8ba9bSPeter Brune ierr = SNESSolve(next->snes,B,Xes[i]);CHKERRQ(ierr); 16790a8ba9bSPeter Brune } 16890a8ba9bSPeter Brune 16990a8ba9bSPeter Brune /* all the solutions are collected; combine optimally */ 17090a8ba9bSPeter Brune for (i=0;i<jac->n;i++) { 17190a8ba9bSPeter Brune for (j=0;j<i+1;j++) { 1729c2f3473SPeter Brune ierr = VecDotBegin(Fes[i],Fes[j],&jac->h[i + j*jac->n]);CHKERRQ(ierr); 17390a8ba9bSPeter Brune } 1749c2f3473SPeter Brune ierr = VecDotBegin(Fes[i],F,&jac->g[i]);CHKERRQ(ierr); 17590a8ba9bSPeter Brune } 1765e806d2eSPeter Brune 17790a8ba9bSPeter Brune for (i=0;i<jac->n;i++) { 17890a8ba9bSPeter Brune for (j=0;j<i+1;j++) { 1799c2f3473SPeter Brune ierr = VecDotEnd(Fes[i],Fes[j],&jac->h[i + j*jac->n]);CHKERRQ(ierr); 1809c2f3473SPeter Brune if (i == j) jac->fnorms[i] = PetscSqrtReal(PetscRealPart(jac->h[i + j*jac->n])); 18190a8ba9bSPeter Brune } 1829c2f3473SPeter Brune ierr = VecDotEnd(Fes[i],F,&jac->g[i]);CHKERRQ(ierr); 18390a8ba9bSPeter Brune } 18490a8ba9bSPeter Brune 1859c2f3473SPeter Brune ftf = (*fnorm)*(*fnorm); 18690a8ba9bSPeter Brune 18790a8ba9bSPeter Brune for (i=0; i<jac->n; i++) { 18890a8ba9bSPeter Brune for (j=i+1;j<jac->n;j++) { 1899c2f3473SPeter Brune jac->h[i + j*jac->n] = jac->h[j + i*jac->n]; 19090a8ba9bSPeter Brune } 19190a8ba9bSPeter Brune } 19290a8ba9bSPeter Brune 19390a8ba9bSPeter Brune for (i=0; i<jac->n; i++) { 1949c2f3473SPeter Brune for (j=0; j<jac->n; j++) { 1959c2f3473SPeter Brune jac->h[i + j*jac->n] = jac->h[i + j*jac->n] - jac->g[j] - jac->g[i] + ftf; 19690a8ba9bSPeter Brune } 1979c2f3473SPeter Brune jac->beta[i] = ftf - jac->g[i]; 1989c2f3473SPeter Brune } 19990a8ba9bSPeter Brune 20090a8ba9bSPeter Brune #if defined(PETSC_MISSING_LAPACK_GELSS) 20190a8ba9bSPeter Brune SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_SUP,"SNESCOMPOSITE with ADDITIVEOPTIMAL requires the LAPACK GELSS routine."); 20290a8ba9bSPeter Brune #else 20390a8ba9bSPeter Brune jac->info = 0; 20490a8ba9bSPeter Brune jac->rcond = -1.; 20590a8ba9bSPeter Brune ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr); 20690a8ba9bSPeter Brune #if defined(PETSC_USE_COMPLEX) 2079c2f3473SPeter 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)); 20890a8ba9bSPeter Brune #else 2099c2f3473SPeter 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)); 21090a8ba9bSPeter Brune #endif 21190a8ba9bSPeter Brune ierr = PetscFPTrapPop();CHKERRQ(ierr); 21290a8ba9bSPeter Brune if (jac->info < 0) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_LIB,"Bad argument to GELSS"); 21390a8ba9bSPeter Brune if (jac->info > 0) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_LIB,"SVD failed to converge"); 21490a8ba9bSPeter Brune #endif 2159c2f3473SPeter Brune tot = 0.; 2165e806d2eSPeter Brune total = 0.; 21790a8ba9bSPeter Brune for (i=0; i<jac->n; i++) { 21890a8ba9bSPeter Brune if (PetscIsInfOrNanScalar(jac->beta[i])) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_LIB,"SVD generated inconsistent output"); 219b3000dc2SPeter Brune ierr = PetscInfo2(snes,"%d: %f\n",i,PetscRealPart(jac->beta[i]));CHKERRQ(ierr); 2209c2f3473SPeter Brune tot += jac->beta[i]; 2215e806d2eSPeter Brune total += PetscAbsScalar(jac->beta[i]); 22290a8ba9bSPeter Brune } 2239c2f3473SPeter Brune ierr = VecScale(X,(1. - tot));CHKERRQ(ierr); 22490a8ba9bSPeter Brune ierr = VecMAXPY(X,jac->n,jac->beta,Xes);CHKERRQ(ierr); 22590a8ba9bSPeter Brune ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 2269c2f3473SPeter Brune ierr = VecNorm(F,NORM_2,fnorm);CHKERRQ(ierr); 22790a8ba9bSPeter Brune 2285e806d2eSPeter Brune /* take the minimum-normed candidate if it beats the combination by a factor of rtol or the combination has stagnated */ 2295e806d2eSPeter Brune min_fnorm = jac->fnorms[0]; 2305e806d2eSPeter Brune min_i = 0; 2319c2f3473SPeter Brune for (i=0; i<jac->n; i++) { 2325e806d2eSPeter Brune if (jac->fnorms[i] < min_fnorm) { 2335e806d2eSPeter Brune min_fnorm = jac->fnorms[i]; 2345e806d2eSPeter Brune min_i = i; 2359c2f3473SPeter Brune } 2369c2f3473SPeter Brune } 2375e806d2eSPeter Brune 2385e806d2eSPeter Brune /* stagnation or divergence restart to the solution of the solver that failed the least */ 2395e806d2eSPeter Brune if (PetscRealPart(total) < jac->stol || min_fnorm*jac->rtol < *fnorm) { 2405e806d2eSPeter Brune ierr = VecCopy(X,jac->Xes[min_i]);CHKERRQ(ierr); 2415e806d2eSPeter Brune ierr = VecCopy(F,jac->Fes[min_i]);CHKERRQ(ierr); 2425e806d2eSPeter Brune *fnorm = min_fnorm; 2435e806d2eSPeter Brune } 24490a8ba9bSPeter Brune PetscFunctionReturn(0); 24590a8ba9bSPeter Brune } 24690a8ba9bSPeter Brune 247eed5f15bSPeter Brune #undef __FUNCT__ 248eed5f15bSPeter Brune #define __FUNCT__ "SNESSetUp_Composite" 249eed5f15bSPeter Brune static PetscErrorCode SNESSetUp_Composite(SNES snes) 250eed5f15bSPeter Brune { 251eed5f15bSPeter Brune PetscErrorCode ierr; 252eed5f15bSPeter Brune SNES_Composite *jac = (SNES_Composite*)snes->data; 253eed5f15bSPeter Brune SNES_CompositeLink next = jac->head; 25490a8ba9bSPeter Brune PetscInt n=0,i; 25590a8ba9bSPeter Brune Vec F; 256eed5f15bSPeter Brune 257eed5f15bSPeter Brune PetscFunctionBegin; 258eed5f15bSPeter Brune while (next) { 25990a8ba9bSPeter Brune n++; 260eed5f15bSPeter Brune ierr = SNESSetFromOptions(next->snes);CHKERRQ(ierr); 261eed5f15bSPeter Brune next = next->next; 262eed5f15bSPeter Brune } 26390a8ba9bSPeter Brune jac->nsnes = n; 26490a8ba9bSPeter Brune ierr = SNESGetFunction(snes,&F,NULL,NULL);CHKERRQ(ierr); 26590a8ba9bSPeter Brune if (jac->type == SNES_COMPOSITE_ADDITIVEOPTIMAL) { 26690a8ba9bSPeter Brune ierr = VecDuplicateVecs(F,jac->nsnes,&jac->Xes);CHKERRQ(ierr); 26790a8ba9bSPeter Brune ierr = PetscMalloc(sizeof(Vec)*n,&jac->Fes);CHKERRQ(ierr); 2689c2f3473SPeter Brune ierr = PetscMalloc(sizeof(PetscReal)*n,&jac->fnorms);CHKERRQ(ierr); 26990a8ba9bSPeter Brune next = jac->head; 27090a8ba9bSPeter Brune i = 0; 27190a8ba9bSPeter Brune while (next) { 27290a8ba9bSPeter Brune ierr = SNESGetFunction(next->snes,&F,NULL,NULL);CHKERRQ(ierr); 27390a8ba9bSPeter Brune jac->Fes[i] = F; 27490a8ba9bSPeter Brune ierr = PetscObjectReference((PetscObject)F);CHKERRQ(ierr); 27590a8ba9bSPeter Brune next = next->next; 27690a8ba9bSPeter Brune i++; 27790a8ba9bSPeter Brune } 27890a8ba9bSPeter Brune /* allocate the subspace direct solve area */ 27990a8ba9bSPeter Brune jac->nrhs = 1; 2809c2f3473SPeter Brune jac->lda = jac->nsnes; 28190a8ba9bSPeter Brune jac->ldb = jac->nsnes; 28290a8ba9bSPeter Brune jac->n = jac->nsnes; 28390a8ba9bSPeter Brune 2849c2f3473SPeter Brune ierr = PetscMalloc(jac->n*jac->n*sizeof(PetscScalar),&jac->h);CHKERRQ(ierr); 2859c2f3473SPeter Brune ierr = PetscMalloc(jac->n*sizeof(PetscScalar),&jac->beta);CHKERRQ(ierr); 2869c2f3473SPeter Brune ierr = PetscMalloc(jac->n*sizeof(PetscScalar),&jac->s);CHKERRQ(ierr); 2879c2f3473SPeter Brune ierr = PetscMalloc(jac->n*sizeof(PetscScalar),&jac->g);CHKERRQ(ierr); 2889c2f3473SPeter Brune jac->lwork = 12*jac->n; 28990a8ba9bSPeter Brune #if PETSC_USE_COMPLEX 29090a8ba9bSPeter Brune ierr = PetscMalloc(sizeof(PetscReal)*jac->lwork,&jac->rwork);CHKERRQ(ierr); 29190a8ba9bSPeter Brune #endif 29290a8ba9bSPeter Brune ierr = PetscMalloc(sizeof(PetscScalar)*jac->lwork,&jac->work);CHKERRQ(ierr); 29390a8ba9bSPeter Brune } 29490a8ba9bSPeter Brune 295eed5f15bSPeter Brune PetscFunctionReturn(0); 296eed5f15bSPeter Brune } 297eed5f15bSPeter Brune 298eed5f15bSPeter Brune #undef __FUNCT__ 299eed5f15bSPeter Brune #define __FUNCT__ "SNESReset_Composite" 300eed5f15bSPeter Brune static PetscErrorCode SNESReset_Composite(SNES snes) 301eed5f15bSPeter Brune { 302eed5f15bSPeter Brune SNES_Composite *jac = (SNES_Composite*)snes->data; 303eed5f15bSPeter Brune PetscErrorCode ierr; 304eed5f15bSPeter Brune SNES_CompositeLink next = jac->head; 305eed5f15bSPeter Brune 306eed5f15bSPeter Brune PetscFunctionBegin; 307eed5f15bSPeter Brune while (next) { 308eed5f15bSPeter Brune ierr = SNESReset(next->snes);CHKERRQ(ierr); 309eed5f15bSPeter Brune next = next->next; 310eed5f15bSPeter Brune } 311eed5f15bSPeter Brune ierr = VecDestroy(&jac->Xorig);CHKERRQ(ierr); 31290a8ba9bSPeter Brune if (jac->Xes) {ierr = VecDestroyVecs(jac->nsnes,&jac->Xes);CHKERRQ(ierr);} 31390a8ba9bSPeter Brune if (jac->Fes) {ierr = VecDestroyVecs(jac->nsnes,&jac->Fes);CHKERRQ(ierr);} 3149c2f3473SPeter Brune ierr = PetscFree(jac->fnorms);CHKERRQ(ierr); 31590a8ba9bSPeter Brune ierr = PetscFree(jac->h);CHKERRQ(ierr); 31690a8ba9bSPeter Brune ierr = PetscFree(jac->s);CHKERRQ(ierr); 3179c2f3473SPeter Brune ierr = PetscFree(jac->g);CHKERRQ(ierr); 31890a8ba9bSPeter Brune ierr = PetscFree(jac->beta);CHKERRQ(ierr); 31990a8ba9bSPeter Brune ierr = PetscFree(jac->work);CHKERRQ(ierr); 32090a8ba9bSPeter Brune ierr = PetscFree(jac->rwork);CHKERRQ(ierr); 32190a8ba9bSPeter Brune 322eed5f15bSPeter Brune PetscFunctionReturn(0); 323eed5f15bSPeter Brune } 324eed5f15bSPeter Brune 325eed5f15bSPeter Brune #undef __FUNCT__ 326eed5f15bSPeter Brune #define __FUNCT__ "SNESDestroy_Composite" 327eed5f15bSPeter Brune static PetscErrorCode SNESDestroy_Composite(SNES snes) 328eed5f15bSPeter Brune { 329eed5f15bSPeter Brune SNES_Composite *jac = (SNES_Composite*)snes->data; 330eed5f15bSPeter Brune PetscErrorCode ierr; 331eed5f15bSPeter Brune SNES_CompositeLink next = jac->head,next_tmp; 332eed5f15bSPeter Brune 333eed5f15bSPeter Brune PetscFunctionBegin; 334eed5f15bSPeter Brune ierr = SNESReset_Composite(snes);CHKERRQ(ierr); 335eed5f15bSPeter Brune while (next) { 336eed5f15bSPeter Brune ierr = SNESDestroy(&next->snes);CHKERRQ(ierr); 337eed5f15bSPeter Brune next_tmp = next; 338eed5f15bSPeter Brune next = next->next; 339eed5f15bSPeter Brune ierr = PetscFree(next_tmp);CHKERRQ(ierr); 340eed5f15bSPeter Brune } 341eed5f15bSPeter Brune ierr = PetscFree(snes->data);CHKERRQ(ierr); 342eed5f15bSPeter Brune PetscFunctionReturn(0); 343eed5f15bSPeter Brune } 344eed5f15bSPeter Brune 345eed5f15bSPeter Brune #undef __FUNCT__ 346eed5f15bSPeter Brune #define __FUNCT__ "SNESSetFromOptions_Composite" 347eed5f15bSPeter Brune static PetscErrorCode SNESSetFromOptions_Composite(SNES snes) 348eed5f15bSPeter Brune { 349eed5f15bSPeter Brune SNES_Composite *jac = (SNES_Composite*)snes->data; 350eed5f15bSPeter Brune PetscErrorCode ierr; 351eed5f15bSPeter Brune PetscInt nmax = 8,i; 352eed5f15bSPeter Brune SNES_CompositeLink next; 353eed5f15bSPeter Brune char *sneses[8]; 3548f626970SPeter Brune PetscReal dmps[8]; 355eed5f15bSPeter Brune PetscBool flg; 356eed5f15bSPeter Brune 357eed5f15bSPeter Brune PetscFunctionBegin; 358eed5f15bSPeter Brune ierr = PetscOptionsHead("Composite preconditioner options");CHKERRQ(ierr); 359eed5f15bSPeter Brune ierr = PetscOptionsEnum("-snes_composite_type","Type of composition","SNESCompositeSetType",SNESCompositeTypes,(PetscEnum)jac->type,(PetscEnum*)&jac->type,&flg);CHKERRQ(ierr); 360eed5f15bSPeter Brune if (flg) { 361eed5f15bSPeter Brune ierr = SNESCompositeSetType(snes,jac->type);CHKERRQ(ierr); 362eed5f15bSPeter Brune } 363eed5f15bSPeter Brune ierr = PetscOptionsStringArray("-snes_composite_sneses","List of composite solvers","SNESCompositeAddSNES",sneses,&nmax,&flg);CHKERRQ(ierr); 364eed5f15bSPeter Brune if (flg) { 365eed5f15bSPeter Brune for (i=0; i<nmax; i++) { 366eed5f15bSPeter Brune ierr = SNESCompositeAddSNES(snes,sneses[i]);CHKERRQ(ierr); 367eed5f15bSPeter Brune ierr = PetscFree(sneses[i]);CHKERRQ(ierr); /* deallocate string sneses[i], which is allocated in PetscOptionsStringArray() */ 368eed5f15bSPeter Brune } 369eed5f15bSPeter Brune } 3708f626970SPeter Brune ierr = PetscOptionsRealArray("-snes_composite_damping","Damping of the additive composite solvers","SNESCompositeSetDamping",dmps,&nmax,&flg);CHKERRQ(ierr); 3718f626970SPeter Brune if (flg) { 3728f626970SPeter Brune for (i=0; i<nmax; i++) { 3738f626970SPeter Brune ierr = SNESCompositeSetDamping(snes,i,dmps[i]);CHKERRQ(ierr); 3748f626970SPeter Brune } 3758f626970SPeter Brune } 3765e806d2eSPeter Brune ierr = PetscOptionsReal("-snes_composite_stol","Step tolerance for restart on the additive composite solvers","",jac->stol,&jac->stol,NULL);CHKERRQ(ierr); 3775e806d2eSPeter Brune ierr = PetscOptionsReal("-snes_composite_rtol","Residual tolerance for the additive composite solvers","",jac->rtol,&jac->rtol,NULL);CHKERRQ(ierr); 378eed5f15bSPeter Brune ierr = PetscOptionsTail();CHKERRQ(ierr); 379eed5f15bSPeter Brune 380eed5f15bSPeter Brune next = jac->head; 381eed5f15bSPeter Brune while (next) { 382eed5f15bSPeter Brune ierr = SNESSetFromOptions(next->snes);CHKERRQ(ierr); 383eed5f15bSPeter Brune next = next->next; 384eed5f15bSPeter Brune } 385eed5f15bSPeter Brune PetscFunctionReturn(0); 386eed5f15bSPeter Brune } 387eed5f15bSPeter Brune 388eed5f15bSPeter Brune #undef __FUNCT__ 389eed5f15bSPeter Brune #define __FUNCT__ "SNESView_Composite" 390eed5f15bSPeter Brune static PetscErrorCode SNESView_Composite(SNES snes,PetscViewer viewer) 391eed5f15bSPeter Brune { 392eed5f15bSPeter Brune SNES_Composite *jac = (SNES_Composite*)snes->data; 393eed5f15bSPeter Brune PetscErrorCode ierr; 394eed5f15bSPeter Brune SNES_CompositeLink next = jac->head; 395eed5f15bSPeter Brune PetscBool iascii; 396eed5f15bSPeter Brune 397eed5f15bSPeter Brune PetscFunctionBegin; 398eed5f15bSPeter Brune ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 399eed5f15bSPeter Brune if (iascii) { 400eed5f15bSPeter Brune ierr = PetscViewerASCIIPrintf(viewer,"Composite SNES type - %s\n",SNESCompositeTypes[jac->type]);CHKERRQ(ierr); 401eed5f15bSPeter Brune ierr = PetscViewerASCIIPrintf(viewer,"SNESes on composite preconditioner follow\n");CHKERRQ(ierr); 402eed5f15bSPeter Brune ierr = PetscViewerASCIIPrintf(viewer,"---------------------------------\n");CHKERRQ(ierr); 403eed5f15bSPeter Brune } 404eed5f15bSPeter Brune if (iascii) { 405eed5f15bSPeter Brune ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 406eed5f15bSPeter Brune } 407eed5f15bSPeter Brune while (next) { 408eed5f15bSPeter Brune ierr = SNESView(next->snes,viewer);CHKERRQ(ierr); 409eed5f15bSPeter Brune next = next->next; 410eed5f15bSPeter Brune } 411eed5f15bSPeter Brune if (iascii) { 412eed5f15bSPeter Brune ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 413eed5f15bSPeter Brune ierr = PetscViewerASCIIPrintf(viewer,"---------------------------------\n");CHKERRQ(ierr); 414eed5f15bSPeter Brune } 415eed5f15bSPeter Brune PetscFunctionReturn(0); 416eed5f15bSPeter Brune } 417eed5f15bSPeter Brune 418eed5f15bSPeter Brune /* ------------------------------------------------------------------------------*/ 419eed5f15bSPeter Brune 420eed5f15bSPeter Brune #undef __FUNCT__ 421eed5f15bSPeter Brune #define __FUNCT__ "SNESCompositeSetType_Composite" 422eed5f15bSPeter Brune static PetscErrorCode SNESCompositeSetType_Composite(SNES snes,SNESCompositeType type) 423eed5f15bSPeter Brune { 424eed5f15bSPeter Brune SNES_Composite *jac = (SNES_Composite*)snes->data; 425eed5f15bSPeter Brune 426eed5f15bSPeter Brune PetscFunctionBegin; 427eed5f15bSPeter Brune jac->type = type; 428eed5f15bSPeter Brune PetscFunctionReturn(0); 429eed5f15bSPeter Brune } 430eed5f15bSPeter Brune 431eed5f15bSPeter Brune #undef __FUNCT__ 432eed5f15bSPeter Brune #define __FUNCT__ "SNESCompositeAddSNES_Composite" 433eed5f15bSPeter Brune static PetscErrorCode SNESCompositeAddSNES_Composite(SNES snes,SNESType type) 434eed5f15bSPeter Brune { 435eed5f15bSPeter Brune SNES_Composite *jac; 436eed5f15bSPeter Brune SNES_CompositeLink next,ilink; 437eed5f15bSPeter Brune PetscErrorCode ierr; 438eed5f15bSPeter Brune PetscInt cnt = 0; 439eed5f15bSPeter Brune const char *prefix; 440eed5f15bSPeter Brune char newprefix[8]; 441eed5f15bSPeter Brune DM dm; 442eed5f15bSPeter Brune 443eed5f15bSPeter Brune PetscFunctionBegin; 444eed5f15bSPeter Brune ierr = PetscNewLog(snes,struct _SNES_CompositeLink,&ilink);CHKERRQ(ierr); 445eed5f15bSPeter Brune ilink->next = 0; 446eed5f15bSPeter Brune ierr = SNESCreate(PetscObjectComm((PetscObject)snes),&ilink->snes);CHKERRQ(ierr); 447eed5f15bSPeter Brune ierr = PetscLogObjectParent((PetscObject)snes,(PetscObject)ilink->snes);CHKERRQ(ierr); 448eed5f15bSPeter Brune ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 449eed5f15bSPeter Brune ierr = SNESSetDM(ilink->snes,dm);CHKERRQ(ierr); 450*cf5b3eb5SPeter Brune ierr = SNESSetTolerances(ilink->snes,ilink->snes->abstol,ilink->snes->rtol,ilink->snes->stol,1,ilink->snes->max_funcs);CHKERRQ(ierr); 451eed5f15bSPeter Brune jac = (SNES_Composite*)snes->data; 452eed5f15bSPeter Brune next = jac->head; 453eed5f15bSPeter Brune if (!next) { 454eed5f15bSPeter Brune jac->head = ilink; 455eed5f15bSPeter Brune ilink->previous = NULL; 456eed5f15bSPeter Brune } else { 457eed5f15bSPeter Brune cnt++; 458eed5f15bSPeter Brune while (next->next) { 459eed5f15bSPeter Brune next = next->next; 460eed5f15bSPeter Brune cnt++; 461eed5f15bSPeter Brune } 462eed5f15bSPeter Brune next->next = ilink; 463eed5f15bSPeter Brune ilink->previous = next; 464eed5f15bSPeter Brune } 465eed5f15bSPeter Brune ierr = SNESGetOptionsPrefix(snes,&prefix);CHKERRQ(ierr); 466eed5f15bSPeter Brune ierr = SNESSetOptionsPrefix(ilink->snes,prefix);CHKERRQ(ierr); 467eed5f15bSPeter Brune sprintf(newprefix,"sub_%d_",(int)cnt); 468eed5f15bSPeter Brune ierr = SNESAppendOptionsPrefix(ilink->snes,newprefix);CHKERRQ(ierr); 4698f626970SPeter Brune ierr = PetscObjectIncrementTabLevel((PetscObject)ilink->snes,(PetscObject)snes,1);CHKERRQ(ierr); 470eed5f15bSPeter Brune ierr = SNESSetType(ilink->snes,type);CHKERRQ(ierr); 4719c2f3473SPeter Brune ierr = SNESSetNormType(ilink->snes, SNES_NORM_FINAL_ONLY);CHKERRQ(ierr); 4728f626970SPeter Brune ilink->dmp = 1.0; 47390a8ba9bSPeter Brune jac->nsnes++; 474eed5f15bSPeter Brune PetscFunctionReturn(0); 475eed5f15bSPeter Brune } 476eed5f15bSPeter Brune 477eed5f15bSPeter Brune #undef __FUNCT__ 478eed5f15bSPeter Brune #define __FUNCT__ "SNESCompositeGetSNES_Composite" 479eed5f15bSPeter Brune static PetscErrorCode SNESCompositeGetSNES_Composite(SNES snes,PetscInt n,SNES *subsnes) 480eed5f15bSPeter Brune { 481eed5f15bSPeter Brune SNES_Composite *jac; 482eed5f15bSPeter Brune SNES_CompositeLink next; 483eed5f15bSPeter Brune PetscInt i; 484eed5f15bSPeter Brune 485eed5f15bSPeter Brune PetscFunctionBegin; 486eed5f15bSPeter Brune jac = (SNES_Composite*)snes->data; 487eed5f15bSPeter Brune next = jac->head; 488eed5f15bSPeter Brune for (i=0; i<n; i++) { 489eed5f15bSPeter Brune if (!next->next) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_INCOMP,"Not enough SNESes in composite preconditioner"); 490eed5f15bSPeter Brune next = next->next; 491eed5f15bSPeter Brune } 492eed5f15bSPeter Brune *subsnes = next->snes; 493eed5f15bSPeter Brune PetscFunctionReturn(0); 494eed5f15bSPeter Brune } 495eed5f15bSPeter Brune 496eed5f15bSPeter Brune /* -------------------------------------------------------------------------------- */ 497eed5f15bSPeter Brune #undef __FUNCT__ 498eed5f15bSPeter Brune #define __FUNCT__ "SNESCompositeSetType" 499e31ab4f9SPeter Brune /*@C 500eed5f15bSPeter Brune SNESCompositeSetType - Sets the type of composite preconditioner. 501eed5f15bSPeter Brune 502eed5f15bSPeter Brune Logically Collective on SNES 503eed5f15bSPeter Brune 504eed5f15bSPeter Brune Input Parameter: 505eed5f15bSPeter Brune + snes - the preconditioner context 506eed5f15bSPeter Brune - type - SNES_COMPOSITE_ADDITIVE (default), SNES_COMPOSITE_MULTIPLICATIVE 507eed5f15bSPeter Brune 508eed5f15bSPeter Brune Options Database Key: 509eed5f15bSPeter Brune . -snes_composite_type <type: one of multiplicative, additive, special> - Sets composite preconditioner type 510eed5f15bSPeter Brune 511eed5f15bSPeter Brune Level: Developer 512eed5f15bSPeter Brune 513eed5f15bSPeter Brune .keywords: SNES, set, type, composite preconditioner, additive, multiplicative 514eed5f15bSPeter Brune @*/ 515eed5f15bSPeter Brune PetscErrorCode SNESCompositeSetType(SNES snes,SNESCompositeType type) 516eed5f15bSPeter Brune { 517eed5f15bSPeter Brune PetscErrorCode ierr; 518eed5f15bSPeter Brune 519eed5f15bSPeter Brune PetscFunctionBegin; 520eed5f15bSPeter Brune PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 521eed5f15bSPeter Brune PetscValidLogicalCollectiveEnum(snes,type,2); 522eed5f15bSPeter Brune ierr = PetscTryMethod(snes,"SNESCompositeSetType_C",(SNES,SNESCompositeType),(snes,type));CHKERRQ(ierr); 523eed5f15bSPeter Brune PetscFunctionReturn(0); 524eed5f15bSPeter Brune } 525eed5f15bSPeter Brune 526eed5f15bSPeter Brune #undef __FUNCT__ 527eed5f15bSPeter Brune #define __FUNCT__ "SNESCompositeAddSNES" 528eed5f15bSPeter Brune /*@C 529eed5f15bSPeter Brune SNESCompositeAddSNES - Adds another SNES to the composite SNES. 530eed5f15bSPeter Brune 531eed5f15bSPeter Brune Collective on SNES 532eed5f15bSPeter Brune 533eed5f15bSPeter Brune Input Parameters: 534eed5f15bSPeter Brune + snes - the preconditioner context 535eed5f15bSPeter Brune - type - the type of the new preconditioner 536eed5f15bSPeter Brune 537eed5f15bSPeter Brune Level: Developer 538eed5f15bSPeter Brune 539eed5f15bSPeter Brune .keywords: SNES, composite preconditioner, add 540eed5f15bSPeter Brune @*/ 541eed5f15bSPeter Brune PetscErrorCode SNESCompositeAddSNES(SNES snes,SNESType type) 542eed5f15bSPeter Brune { 543eed5f15bSPeter Brune PetscErrorCode ierr; 544eed5f15bSPeter Brune 545eed5f15bSPeter Brune PetscFunctionBegin; 546eed5f15bSPeter Brune PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 547eed5f15bSPeter Brune ierr = PetscTryMethod(snes,"SNESCompositeAddSNES_C",(SNES,SNESType),(snes,type));CHKERRQ(ierr); 548eed5f15bSPeter Brune PetscFunctionReturn(0); 549eed5f15bSPeter Brune } 550eed5f15bSPeter Brune #undef __FUNCT__ 551eed5f15bSPeter Brune #define __FUNCT__ "SNESCompositeGetSNES" 552eed5f15bSPeter Brune /*@ 553eed5f15bSPeter Brune SNESCompositeGetSNES - Gets one of the SNES objects in the composite SNES. 554eed5f15bSPeter Brune 555eed5f15bSPeter Brune Not Collective 556eed5f15bSPeter Brune 557eed5f15bSPeter Brune Input Parameter: 558eed5f15bSPeter Brune + snes - the preconditioner context 559eed5f15bSPeter Brune - n - the number of the snes requested 560eed5f15bSPeter Brune 561eed5f15bSPeter Brune Output Parameters: 562eed5f15bSPeter Brune . subsnes - the SNES requested 563eed5f15bSPeter Brune 564eed5f15bSPeter Brune Level: Developer 565eed5f15bSPeter Brune 566eed5f15bSPeter Brune .keywords: SNES, get, composite preconditioner, sub preconditioner 567eed5f15bSPeter Brune 568eed5f15bSPeter Brune .seealso: SNESCompositeAddSNES() 569eed5f15bSPeter Brune @*/ 570eed5f15bSPeter Brune PetscErrorCode SNESCompositeGetSNES(SNES snes,PetscInt n,SNES *subsnes) 571eed5f15bSPeter Brune { 572eed5f15bSPeter Brune PetscErrorCode ierr; 573eed5f15bSPeter Brune 574eed5f15bSPeter Brune PetscFunctionBegin; 575eed5f15bSPeter Brune PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 576eed5f15bSPeter Brune PetscValidPointer(subsnes,3); 577eed5f15bSPeter Brune ierr = PetscUseMethod(snes,"SNESCompositeGetSNES_C",(SNES,PetscInt,SNES*),(snes,n,subsnes));CHKERRQ(ierr); 578eed5f15bSPeter Brune PetscFunctionReturn(0); 579eed5f15bSPeter Brune } 580eed5f15bSPeter Brune 581eed5f15bSPeter Brune #undef __FUNCT__ 5828f626970SPeter Brune #define __FUNCT__ "SNESCompositeSetDamping_Composite" 5838f626970SPeter Brune static PetscErrorCode SNESCompositeSetDamping_Composite(SNES snes,PetscInt n,PetscReal dmp) 5848f626970SPeter Brune { 5858f626970SPeter Brune SNES_Composite *jac; 5868f626970SPeter Brune SNES_CompositeLink next; 5878f626970SPeter Brune PetscInt i; 5888f626970SPeter Brune 5898f626970SPeter Brune PetscFunctionBegin; 5908f626970SPeter Brune jac = (SNES_Composite*)snes->data; 5918f626970SPeter Brune next = jac->head; 5928f626970SPeter Brune for (i=0; i<n; i++) { 5938f626970SPeter Brune if (!next->next) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_INCOMP,"Not enough SNESes in composite preconditioner"); 5948f626970SPeter Brune next = next->next; 5958f626970SPeter Brune } 5968f626970SPeter Brune next->dmp = dmp; 5978f626970SPeter Brune PetscFunctionReturn(0); 5988f626970SPeter Brune } 5998f626970SPeter Brune 6008f626970SPeter Brune #undef __FUNCT__ 6018f626970SPeter Brune #define __FUNCT__ "SNESCompositeSetDamping" 6028f626970SPeter Brune /*@ 6038f626970SPeter Brune SNESCompositeSetDamping - Sets the damping of a subsolver when using additive composite SNES. 6048f626970SPeter Brune 6058f626970SPeter Brune Not Collective 6068f626970SPeter Brune 6078f626970SPeter Brune Input Parameter: 6088f626970SPeter Brune + snes - the preconditioner context 6098f626970SPeter Brune . n - the number of the snes requested 6108f626970SPeter Brune - dmp - the damping 6118f626970SPeter Brune 6128f626970SPeter Brune Level: Developer 6138f626970SPeter Brune 6148f626970SPeter Brune .keywords: SNES, get, composite preconditioner, sub preconditioner 6158f626970SPeter Brune 6168f626970SPeter Brune .seealso: SNESCompositeAddSNES() 6178f626970SPeter Brune @*/ 6188f626970SPeter Brune PetscErrorCode SNESCompositeSetDamping(SNES snes,PetscInt n,PetscReal dmp) 6198f626970SPeter Brune { 6208f626970SPeter Brune PetscErrorCode ierr; 6218f626970SPeter Brune 6228f626970SPeter Brune PetscFunctionBegin; 6238f626970SPeter Brune PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 6248f626970SPeter Brune ierr = PetscUseMethod(snes,"SNESCompositeSetDamping_C",(SNES,PetscInt,PetscReal),(snes,n,dmp));CHKERRQ(ierr); 6258f626970SPeter Brune PetscFunctionReturn(0); 6268f626970SPeter Brune } 6278f626970SPeter Brune 6288f626970SPeter Brune #undef __FUNCT__ 629eed5f15bSPeter Brune #define __FUNCT__ "SNESSolve_Composite" 630eed5f15bSPeter Brune PetscErrorCode SNESSolve_Composite(SNES snes) 631eed5f15bSPeter Brune { 632eed5f15bSPeter Brune Vec F; 633eed5f15bSPeter Brune Vec X; 634eed5f15bSPeter Brune Vec B; 635eed5f15bSPeter Brune PetscInt i; 636eed5f15bSPeter Brune PetscReal fnorm = 0.0; 637eed5f15bSPeter Brune PetscErrorCode ierr; 638eed5f15bSPeter Brune SNESNormType normtype; 639eed5f15bSPeter Brune SNES_Composite *comp = (SNES_Composite*)snes->data; 640eed5f15bSPeter Brune 641eed5f15bSPeter Brune PetscFunctionBegin; 642eed5f15bSPeter Brune X = snes->vec_sol; 643eed5f15bSPeter Brune F = snes->vec_func; 644eed5f15bSPeter Brune B = snes->vec_rhs; 645eed5f15bSPeter Brune 646eed5f15bSPeter Brune ierr = PetscObjectAMSTakeAccess((PetscObject)snes);CHKERRQ(ierr); 647eed5f15bSPeter Brune snes->iter = 0; 648eed5f15bSPeter Brune snes->norm = 0.; 649eed5f15bSPeter Brune ierr = PetscObjectAMSGrantAccess((PetscObject)snes);CHKERRQ(ierr); 650eed5f15bSPeter Brune snes->reason = SNES_CONVERGED_ITERATING; 651eed5f15bSPeter Brune ierr = SNESGetNormType(snes, &normtype);CHKERRQ(ierr); 652eed5f15bSPeter Brune if (normtype == SNES_NORM_FUNCTION || normtype == SNES_NORM_INITIAL_ONLY || normtype == SNES_NORM_INITIAL_FINAL_ONLY) { 653eed5f15bSPeter Brune if (!snes->vec_func_init_set) { 654eed5f15bSPeter Brune ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 655eed5f15bSPeter Brune if (snes->domainerror) { 656eed5f15bSPeter Brune snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 657eed5f15bSPeter Brune PetscFunctionReturn(0); 658eed5f15bSPeter Brune } 659eed5f15bSPeter Brune } else snes->vec_func_init_set = PETSC_FALSE; 660eed5f15bSPeter Brune 661eed5f15bSPeter Brune /* convergence test */ 662eed5f15bSPeter Brune if (!snes->norm_init_set) { 663eed5f15bSPeter Brune ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F|| */ 664eed5f15bSPeter Brune if (PetscIsInfOrNanReal(fnorm)) { 665eed5f15bSPeter Brune snes->reason = SNES_DIVERGED_FNORM_NAN; 666eed5f15bSPeter Brune PetscFunctionReturn(0); 667eed5f15bSPeter Brune } 668eed5f15bSPeter Brune } else { 669eed5f15bSPeter Brune fnorm = snes->norm_init; 670eed5f15bSPeter Brune snes->norm_init_set = PETSC_FALSE; 671eed5f15bSPeter Brune } 672eed5f15bSPeter Brune ierr = PetscObjectAMSTakeAccess((PetscObject)snes);CHKERRQ(ierr); 673eed5f15bSPeter Brune snes->iter = 0; 674eed5f15bSPeter Brune snes->norm = fnorm; 675eed5f15bSPeter Brune ierr = PetscObjectAMSGrantAccess((PetscObject)snes);CHKERRQ(ierr); 676eed5f15bSPeter Brune ierr = SNESLogConvergenceHistory(snes,snes->norm,0);CHKERRQ(ierr); 677eed5f15bSPeter Brune ierr = SNESMonitor(snes,0,snes->norm);CHKERRQ(ierr); 678eed5f15bSPeter Brune 679eed5f15bSPeter Brune /* set parameter for default relative tolerance convergence test */ 680eed5f15bSPeter Brune snes->ttol = fnorm*snes->rtol; 681eed5f15bSPeter Brune 682eed5f15bSPeter Brune /* test convergence */ 683eed5f15bSPeter Brune ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 684eed5f15bSPeter Brune if (snes->reason) PetscFunctionReturn(0); 685eed5f15bSPeter Brune } else { 686eed5f15bSPeter Brune ierr = PetscObjectAMSGrantAccess((PetscObject)snes);CHKERRQ(ierr); 687eed5f15bSPeter Brune ierr = SNESLogConvergenceHistory(snes,snes->norm,0);CHKERRQ(ierr); 688eed5f15bSPeter Brune ierr = SNESMonitor(snes,0,snes->norm);CHKERRQ(ierr); 689eed5f15bSPeter Brune } 690eed5f15bSPeter Brune 691eed5f15bSPeter Brune /* Call general purpose update function */ 692eed5f15bSPeter Brune if (snes->ops->update) { 693eed5f15bSPeter Brune ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 694eed5f15bSPeter Brune } 695eed5f15bSPeter Brune 696eed5f15bSPeter Brune for (i = 0; i < snes->max_its; i++) { 697eed5f15bSPeter Brune if (comp->type == SNES_COMPOSITE_ADDITIVE) { 698eed5f15bSPeter Brune ierr = SNESCompositeApply_Additive(snes,X,B,F,&fnorm);CHKERRQ(ierr); 699eed5f15bSPeter Brune } else if (comp->type == SNES_COMPOSITE_MULTIPLICATIVE) { 700eed5f15bSPeter Brune ierr = SNESCompositeApply_Multiplicative(snes,X,B,F,&fnorm);CHKERRQ(ierr); 70190a8ba9bSPeter Brune } else if (comp->type == SNES_COMPOSITE_ADDITIVEOPTIMAL) { 70290a8ba9bSPeter Brune ierr = SNESCompositeApply_AdditiveOptimal(snes,X,B,F,&fnorm);CHKERRQ(ierr); 703eed5f15bSPeter Brune } else { 70490a8ba9bSPeter Brune SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"Unsupported SNESComposite type"); 705eed5f15bSPeter Brune } 706eed5f15bSPeter Brune if ((i == snes->max_its - 1) && (normtype == SNES_NORM_INITIAL_FINAL_ONLY || normtype == SNES_NORM_FINAL_ONLY)) { 707eed5f15bSPeter Brune ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 708eed5f15bSPeter Brune if (snes->domainerror) { 709eed5f15bSPeter Brune snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 710eed5f15bSPeter Brune break; 711eed5f15bSPeter Brune } 712eed5f15bSPeter Brune ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F|| */ 713eed5f15bSPeter Brune if (PetscIsInfOrNanReal(fnorm)) { 714eed5f15bSPeter Brune snes->reason = SNES_DIVERGED_FNORM_NAN; 715eed5f15bSPeter Brune break; 716eed5f15bSPeter Brune } 717eed5f15bSPeter Brune } 718eed5f15bSPeter Brune /* Monitor convergence */ 719eed5f15bSPeter Brune ierr = PetscObjectAMSTakeAccess((PetscObject)snes);CHKERRQ(ierr); 720eed5f15bSPeter Brune snes->iter = i+1; 721eed5f15bSPeter Brune snes->norm = fnorm; 722eed5f15bSPeter Brune ierr = PetscObjectAMSGrantAccess((PetscObject)snes);CHKERRQ(ierr); 723eed5f15bSPeter Brune ierr = SNESLogConvergenceHistory(snes,snes->norm,0);CHKERRQ(ierr); 724eed5f15bSPeter Brune ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr); 725eed5f15bSPeter Brune /* Test for convergence */ 726eed5f15bSPeter Brune if (normtype == SNES_NORM_FUNCTION) {ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);} 727eed5f15bSPeter Brune if (snes->reason) break; 728eed5f15bSPeter Brune /* Call general purpose update function */ 729eed5f15bSPeter Brune if (snes->ops->update) {ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);} 730eed5f15bSPeter Brune } 731eed5f15bSPeter Brune if (normtype == SNES_NORM_FUNCTION) { 732eed5f15bSPeter Brune if (i == snes->max_its) { 733eed5f15bSPeter Brune ierr = PetscInfo1(snes,"Maximum number of iterations has been reached: %D\n",snes->max_its);CHKERRQ(ierr); 734eed5f15bSPeter Brune if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT; 735eed5f15bSPeter Brune } 736eed5f15bSPeter Brune } else if (!snes->reason) snes->reason = SNES_CONVERGED_ITS; 737eed5f15bSPeter Brune PetscFunctionReturn(0); 738eed5f15bSPeter Brune } 739eed5f15bSPeter Brune 740eed5f15bSPeter Brune /* -------------------------------------------------------------------------------------------*/ 741eed5f15bSPeter Brune 742eed5f15bSPeter Brune /*MC 743eed5f15bSPeter Brune SNESCOMPOSITE - Build a preconditioner by composing together several nonlinear solvers 744eed5f15bSPeter Brune 745eed5f15bSPeter Brune Options Database Keys: 746eed5f15bSPeter Brune + -snes_composite_type <type: one of multiplicative, additive, symmetric_multiplicative, special> - Sets composite preconditioner type 747eed5f15bSPeter Brune - -snes_composite_sneses - <snes0,snes1,...> list of SNESes to compose 748eed5f15bSPeter Brune 749eed5f15bSPeter Brune Level: intermediate 750eed5f15bSPeter Brune 751eed5f15bSPeter Brune Concepts: composing solvers 752eed5f15bSPeter Brune 753eed5f15bSPeter Brune .seealso: SNESCreate(), SNESSetType(), SNESType (for list of available types), SNES, 754eed5f15bSPeter Brune SNESSHELL, SNESCompositeSetType(), SNESCompositeSpecialSetAlpha(), SNESCompositeAddSNES(), 755eed5f15bSPeter Brune SNESCompositeGetSNES() 756eed5f15bSPeter Brune 757eed5f15bSPeter Brune M*/ 758eed5f15bSPeter Brune 759eed5f15bSPeter Brune #undef __FUNCT__ 760eed5f15bSPeter Brune #define __FUNCT__ "SNESCreate_Composite" 761eed5f15bSPeter Brune PETSC_EXTERN PetscErrorCode SNESCreate_Composite(SNES snes) 762eed5f15bSPeter Brune { 763eed5f15bSPeter Brune PetscErrorCode ierr; 764eed5f15bSPeter Brune SNES_Composite *jac; 765eed5f15bSPeter Brune 766eed5f15bSPeter Brune PetscFunctionBegin; 767eed5f15bSPeter Brune ierr = PetscNewLog(snes,SNES_Composite,&jac);CHKERRQ(ierr); 768eed5f15bSPeter Brune 769eed5f15bSPeter Brune snes->ops->solve = SNESSolve_Composite; 770eed5f15bSPeter Brune snes->ops->setup = SNESSetUp_Composite; 771eed5f15bSPeter Brune snes->ops->reset = SNESReset_Composite; 772eed5f15bSPeter Brune snes->ops->destroy = SNESDestroy_Composite; 773eed5f15bSPeter Brune snes->ops->setfromoptions = SNESSetFromOptions_Composite; 774eed5f15bSPeter Brune snes->ops->view = SNESView_Composite; 775eed5f15bSPeter Brune 776eed5f15bSPeter Brune snes->data = (void*)jac; 77790a8ba9bSPeter Brune jac->type = SNES_COMPOSITE_ADDITIVEOPTIMAL; 77890a8ba9bSPeter Brune jac->Fes = NULL; 77990a8ba9bSPeter Brune jac->Xes = NULL; 7809c2f3473SPeter Brune jac->fnorms = NULL; 78190a8ba9bSPeter Brune jac->nsnes = 0; 782eed5f15bSPeter Brune jac->head = 0; 7835e806d2eSPeter Brune jac->stol = 0.1; 7845e806d2eSPeter Brune jac->rtol = 1.1; 785eed5f15bSPeter Brune 78690a8ba9bSPeter Brune jac->h = NULL; 78790a8ba9bSPeter Brune jac->s = NULL; 78890a8ba9bSPeter Brune jac->beta = NULL; 78990a8ba9bSPeter Brune jac->work = NULL; 79090a8ba9bSPeter Brune jac->rwork = NULL; 79190a8ba9bSPeter Brune 7928f626970SPeter Brune ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESCompositeSetType_C",SNESCompositeSetType_Composite);CHKERRQ(ierr); 7938f626970SPeter Brune ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESCompositeAddSNES_C",SNESCompositeAddSNES_Composite);CHKERRQ(ierr); 7948f626970SPeter Brune ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESCompositeGetSNES_C",SNESCompositeGetSNES_Composite);CHKERRQ(ierr); 7958f626970SPeter Brune ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESCompositeSetDamping_C",SNESCompositeSetDamping_Composite);CHKERRQ(ierr); 796eed5f15bSPeter Brune PetscFunctionReturn(0); 797eed5f15bSPeter Brune } 798eed5f15bSPeter Brune 799