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*/ 6eed5f15bSPeter Brune 7eed5f15bSPeter Brune const char *const SNESCompositeTypes[] = {"ADDITIVE","MULTIPLICATIVE","SNESCompositeType","SNES_COMPOSITE",0}; 8eed5f15bSPeter Brune 9eed5f15bSPeter Brune typedef struct _SNES_CompositeLink *SNES_CompositeLink; 10eed5f15bSPeter Brune struct _SNES_CompositeLink { 11eed5f15bSPeter Brune SNES snes; 128f626970SPeter Brune PetscReal dmp; 13eed5f15bSPeter Brune SNES_CompositeLink next; 14eed5f15bSPeter Brune SNES_CompositeLink previous; 15eed5f15bSPeter Brune }; 16eed5f15bSPeter Brune 17eed5f15bSPeter Brune typedef struct { 18eed5f15bSPeter Brune SNES_CompositeLink head; 19eed5f15bSPeter Brune SNESCompositeType type; 20eed5f15bSPeter Brune Vec Xorig; 21eed5f15bSPeter Brune } SNES_Composite; 22eed5f15bSPeter Brune 23eed5f15bSPeter Brune #undef __FUNCT__ 24eed5f15bSPeter Brune #define __FUNCT__ "SNESCompositeApply_Multiplicative" 25eed5f15bSPeter Brune static PetscErrorCode SNESCompositeApply_Multiplicative(SNES snes,Vec X,Vec B,Vec F,PetscReal *fnorm) 26eed5f15bSPeter Brune { 27eed5f15bSPeter Brune PetscErrorCode ierr; 28eed5f15bSPeter Brune SNES_Composite *jac = (SNES_Composite*)snes->data; 29eed5f15bSPeter Brune SNES_CompositeLink next = jac->head; 30eed5f15bSPeter Brune Vec FSub; 31eed5f15bSPeter Brune PetscReal fsubnorm; 32eed5f15bSPeter Brune 33eed5f15bSPeter Brune PetscFunctionBegin; 34eed5f15bSPeter Brune if (!next) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"No composite SNESes supplied via SNESCompositeAddSNES() or -snes_composite_sneses"); 35eed5f15bSPeter Brune if (snes->normtype == SNES_NORM_FUNCTION) { 36eed5f15bSPeter Brune ierr = SNESSetInitialFunction(next->snes,F);CHKERRQ(ierr); 37eed5f15bSPeter Brune if (fnorm) {ierr = SNESSetInitialFunctionNorm(next->snes,*fnorm);CHKERRQ(ierr);} 38eed5f15bSPeter Brune } 39eed5f15bSPeter Brune ierr = SNESSolve(next->snes,B,X);CHKERRQ(ierr); 40eed5f15bSPeter Brune 41eed5f15bSPeter Brune while (next->next) { 42eed5f15bSPeter Brune /* only copy the function over in the case where the functions correspond */ 43eed5f15bSPeter Brune if (next->snes->pcside == PC_RIGHT && next->snes->normtype != SNES_NORM_NONE) { 44eed5f15bSPeter Brune ierr = SNESGetFunction(next->snes,&FSub,NULL,NULL);CHKERRQ(ierr); 45eed5f15bSPeter Brune ierr = SNESGetFunctionNorm(next->snes,&fsubnorm);CHKERRQ(ierr); 46eed5f15bSPeter Brune next = next->next; 47eed5f15bSPeter Brune ierr = SNESSetInitialFunction(next->snes,FSub);CHKERRQ(ierr); 48eed5f15bSPeter Brune ierr = SNESSetInitialFunctionNorm(next->snes,fsubnorm);CHKERRQ(ierr); 49eed5f15bSPeter Brune } else { 50eed5f15bSPeter Brune next = next->next; 51eed5f15bSPeter Brune } 52eed5f15bSPeter Brune ierr = SNESSolve(next->snes,B,X);CHKERRQ(ierr); 53eed5f15bSPeter Brune } 54eed5f15bSPeter Brune if (next->snes->pcside == PC_RIGHT) { 55eed5f15bSPeter Brune ierr = SNESGetFunction(next->snes,&FSub,NULL,NULL);CHKERRQ(ierr); 56eed5f15bSPeter Brune ierr = VecCopy(FSub,F);CHKERRQ(ierr); 57eed5f15bSPeter Brune if (fnorm) {ierr = SNESGetFunctionNorm(next->snes,fnorm);CHKERRQ(ierr);} 58eed5f15bSPeter Brune } else if (snes->normtype == SNES_NORM_FUNCTION) { 59eed5f15bSPeter Brune SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 60eed5f15bSPeter Brune if (fnorm) {ierr = VecNorm(F,NORM_2,fnorm);CHKERRQ(ierr);} 61eed5f15bSPeter Brune } 62eed5f15bSPeter Brune PetscFunctionReturn(0); 63eed5f15bSPeter Brune } 64eed5f15bSPeter Brune 65eed5f15bSPeter Brune #undef __FUNCT__ 66eed5f15bSPeter Brune #define __FUNCT__ "SNESCompositeApply_Additive" 67eed5f15bSPeter Brune static PetscErrorCode SNESCompositeApply_Additive(SNES snes,Vec X,Vec B,Vec F,PetscReal *fnorm) 68eed5f15bSPeter Brune { 69eed5f15bSPeter Brune PetscErrorCode ierr; 70eed5f15bSPeter Brune SNES_Composite *jac = (SNES_Composite*)snes->data; 71eed5f15bSPeter Brune SNES_CompositeLink next = jac->head; 72eed5f15bSPeter Brune Vec Y,Xorig; 73eed5f15bSPeter Brune 74eed5f15bSPeter Brune PetscFunctionBegin; 75eed5f15bSPeter Brune Y = snes->vec_sol_update; 76eed5f15bSPeter Brune if (!jac->Xorig) {ierr = VecDuplicate(X,&jac->Xorig);CHKERRQ(ierr);} 77eed5f15bSPeter Brune Xorig = jac->Xorig; 78eed5f15bSPeter Brune ierr = VecCopy(X,Xorig); 79eed5f15bSPeter Brune if (!next) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"No composite SNESes supplied via SNESCompositeAddSNES() or -snes_composite_sneses"); 80eed5f15bSPeter Brune if (snes->normtype == SNES_NORM_FUNCTION) { 81eed5f15bSPeter Brune ierr = SNESSetInitialFunction(next->snes,F);CHKERRQ(ierr); 82eed5f15bSPeter Brune if (fnorm) {ierr = SNESSetInitialFunctionNorm(next->snes,*fnorm);CHKERRQ(ierr);} 83eed5f15bSPeter Brune while (next->next) { 84eed5f15bSPeter Brune next = next->next; 85eed5f15bSPeter Brune ierr = SNESSetInitialFunction(next->snes,F);CHKERRQ(ierr); 86eed5f15bSPeter Brune if (fnorm) {ierr = SNESSetInitialFunctionNorm(next->snes,*fnorm);CHKERRQ(ierr);} 87eed5f15bSPeter Brune } 88eed5f15bSPeter Brune } 89eed5f15bSPeter Brune next = jac->head; 908f626970SPeter Brune ierr = VecCopy(Xorig,Y);CHKERRQ(ierr); 91eed5f15bSPeter Brune ierr = SNESSolve(next->snes,B,Y);CHKERRQ(ierr); 92eed5f15bSPeter Brune ierr = VecAXPY(Y,-1.0,Xorig);CHKERRQ(ierr); 938f626970SPeter Brune ierr = VecAXPY(X,next->dmp,Y);CHKERRQ(ierr); 948f626970SPeter Brune while (next->next) { 958f626970SPeter Brune next = next->next; 968f626970SPeter Brune ierr = VecCopy(Xorig,Y);CHKERRQ(ierr); 978f626970SPeter Brune ierr = SNESSolve(next->snes,B,Y);CHKERRQ(ierr); 988f626970SPeter Brune ierr = VecAXPY(Y,-1.0,Xorig);CHKERRQ(ierr); 998f626970SPeter Brune ierr = VecAXPY(X,next->dmp,Y);CHKERRQ(ierr); 100eed5f15bSPeter Brune } 101eed5f15bSPeter Brune if (snes->normtype == SNES_NORM_FUNCTION) { 102eed5f15bSPeter Brune ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 103eed5f15bSPeter Brune if (fnorm) {ierr = VecNorm(F,NORM_2,fnorm);CHKERRQ(ierr);} 104eed5f15bSPeter Brune } 105eed5f15bSPeter Brune PetscFunctionReturn(0); 106eed5f15bSPeter Brune } 107eed5f15bSPeter Brune #undef __FUNCT__ 108eed5f15bSPeter Brune #define __FUNCT__ "SNESSetUp_Composite" 109eed5f15bSPeter Brune static PetscErrorCode SNESSetUp_Composite(SNES snes) 110eed5f15bSPeter Brune { 111eed5f15bSPeter Brune PetscErrorCode ierr; 112eed5f15bSPeter Brune SNES_Composite *jac = (SNES_Composite*)snes->data; 113eed5f15bSPeter Brune SNES_CompositeLink next = jac->head; 114eed5f15bSPeter Brune 115eed5f15bSPeter Brune PetscFunctionBegin; 116eed5f15bSPeter Brune while (next) { 117eed5f15bSPeter Brune ierr = SNESSetFromOptions(next->snes);CHKERRQ(ierr); 118eed5f15bSPeter Brune next = next->next; 119eed5f15bSPeter Brune } 120eed5f15bSPeter Brune PetscFunctionReturn(0); 121eed5f15bSPeter Brune } 122eed5f15bSPeter Brune 123eed5f15bSPeter Brune #undef __FUNCT__ 124eed5f15bSPeter Brune #define __FUNCT__ "SNESReset_Composite" 125eed5f15bSPeter Brune static PetscErrorCode SNESReset_Composite(SNES snes) 126eed5f15bSPeter Brune { 127eed5f15bSPeter Brune SNES_Composite *jac = (SNES_Composite*)snes->data; 128eed5f15bSPeter Brune PetscErrorCode ierr; 129eed5f15bSPeter Brune SNES_CompositeLink next = jac->head; 130eed5f15bSPeter Brune 131eed5f15bSPeter Brune PetscFunctionBegin; 132eed5f15bSPeter Brune while (next) { 133eed5f15bSPeter Brune ierr = SNESReset(next->snes);CHKERRQ(ierr); 134eed5f15bSPeter Brune next = next->next; 135eed5f15bSPeter Brune } 136eed5f15bSPeter Brune ierr = VecDestroy(&jac->Xorig);CHKERRQ(ierr); 137eed5f15bSPeter Brune PetscFunctionReturn(0); 138eed5f15bSPeter Brune } 139eed5f15bSPeter Brune 140eed5f15bSPeter Brune #undef __FUNCT__ 141eed5f15bSPeter Brune #define __FUNCT__ "SNESDestroy_Composite" 142eed5f15bSPeter Brune static PetscErrorCode SNESDestroy_Composite(SNES snes) 143eed5f15bSPeter Brune { 144eed5f15bSPeter Brune SNES_Composite *jac = (SNES_Composite*)snes->data; 145eed5f15bSPeter Brune PetscErrorCode ierr; 146eed5f15bSPeter Brune SNES_CompositeLink next = jac->head,next_tmp; 147eed5f15bSPeter Brune 148eed5f15bSPeter Brune PetscFunctionBegin; 149eed5f15bSPeter Brune ierr = SNESReset_Composite(snes);CHKERRQ(ierr); 150eed5f15bSPeter Brune while (next) { 151eed5f15bSPeter Brune ierr = SNESDestroy(&next->snes);CHKERRQ(ierr); 152eed5f15bSPeter Brune next_tmp = next; 153eed5f15bSPeter Brune next = next->next; 154eed5f15bSPeter Brune ierr = PetscFree(next_tmp);CHKERRQ(ierr); 155eed5f15bSPeter Brune } 156eed5f15bSPeter Brune ierr = PetscFree(snes->data);CHKERRQ(ierr); 157eed5f15bSPeter Brune PetscFunctionReturn(0); 158eed5f15bSPeter Brune } 159eed5f15bSPeter Brune 160eed5f15bSPeter Brune #undef __FUNCT__ 161eed5f15bSPeter Brune #define __FUNCT__ "SNESSetFromOptions_Composite" 162eed5f15bSPeter Brune static PetscErrorCode SNESSetFromOptions_Composite(SNES snes) 163eed5f15bSPeter Brune { 164eed5f15bSPeter Brune SNES_Composite *jac = (SNES_Composite*)snes->data; 165eed5f15bSPeter Brune PetscErrorCode ierr; 166eed5f15bSPeter Brune PetscInt nmax = 8,i; 167eed5f15bSPeter Brune SNES_CompositeLink next; 168eed5f15bSPeter Brune char *sneses[8]; 1698f626970SPeter Brune PetscReal dmps[8]; 170eed5f15bSPeter Brune PetscBool flg; 171eed5f15bSPeter Brune 172eed5f15bSPeter Brune PetscFunctionBegin; 173eed5f15bSPeter Brune ierr = PetscOptionsHead("Composite preconditioner options");CHKERRQ(ierr); 174eed5f15bSPeter Brune ierr = PetscOptionsEnum("-snes_composite_type","Type of composition","SNESCompositeSetType",SNESCompositeTypes,(PetscEnum)jac->type,(PetscEnum*)&jac->type,&flg);CHKERRQ(ierr); 175eed5f15bSPeter Brune if (flg) { 176eed5f15bSPeter Brune ierr = SNESCompositeSetType(snes,jac->type);CHKERRQ(ierr); 177eed5f15bSPeter Brune } 178eed5f15bSPeter Brune ierr = PetscOptionsStringArray("-snes_composite_sneses","List of composite solvers","SNESCompositeAddSNES",sneses,&nmax,&flg);CHKERRQ(ierr); 179eed5f15bSPeter Brune if (flg) { 180eed5f15bSPeter Brune for (i=0; i<nmax; i++) { 181eed5f15bSPeter Brune ierr = SNESCompositeAddSNES(snes,sneses[i]);CHKERRQ(ierr); 182eed5f15bSPeter Brune ierr = PetscFree(sneses[i]);CHKERRQ(ierr); /* deallocate string sneses[i], which is allocated in PetscOptionsStringArray() */ 183eed5f15bSPeter Brune } 184eed5f15bSPeter Brune } 1858f626970SPeter Brune ierr = PetscOptionsRealArray("-snes_composite_damping","Damping of the additive composite solvers","SNESCompositeSetDamping",dmps,&nmax,&flg);CHKERRQ(ierr); 1868f626970SPeter Brune if (flg) { 1878f626970SPeter Brune for (i=0; i<nmax; i++) { 1888f626970SPeter Brune ierr = SNESCompositeSetDamping(snes,i,dmps[i]);CHKERRQ(ierr); 1898f626970SPeter Brune } 1908f626970SPeter Brune } 191eed5f15bSPeter Brune ierr = PetscOptionsTail();CHKERRQ(ierr); 192eed5f15bSPeter Brune 193eed5f15bSPeter Brune next = jac->head; 194eed5f15bSPeter Brune while (next) { 195eed5f15bSPeter Brune ierr = SNESSetFromOptions(next->snes);CHKERRQ(ierr); 196eed5f15bSPeter Brune next = next->next; 197eed5f15bSPeter Brune } 198eed5f15bSPeter Brune PetscFunctionReturn(0); 199eed5f15bSPeter Brune } 200eed5f15bSPeter Brune 201eed5f15bSPeter Brune #undef __FUNCT__ 202eed5f15bSPeter Brune #define __FUNCT__ "SNESView_Composite" 203eed5f15bSPeter Brune static PetscErrorCode SNESView_Composite(SNES snes,PetscViewer viewer) 204eed5f15bSPeter Brune { 205eed5f15bSPeter Brune SNES_Composite *jac = (SNES_Composite*)snes->data; 206eed5f15bSPeter Brune PetscErrorCode ierr; 207eed5f15bSPeter Brune SNES_CompositeLink next = jac->head; 208eed5f15bSPeter Brune PetscBool iascii; 209eed5f15bSPeter Brune 210eed5f15bSPeter Brune PetscFunctionBegin; 211eed5f15bSPeter Brune ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 212eed5f15bSPeter Brune if (iascii) { 213eed5f15bSPeter Brune ierr = PetscViewerASCIIPrintf(viewer,"Composite SNES type - %s\n",SNESCompositeTypes[jac->type]);CHKERRQ(ierr); 214eed5f15bSPeter Brune ierr = PetscViewerASCIIPrintf(viewer,"SNESes on composite preconditioner follow\n");CHKERRQ(ierr); 215eed5f15bSPeter Brune ierr = PetscViewerASCIIPrintf(viewer,"---------------------------------\n");CHKERRQ(ierr); 216eed5f15bSPeter Brune } 217eed5f15bSPeter Brune if (iascii) { 218eed5f15bSPeter Brune ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 219eed5f15bSPeter Brune } 220eed5f15bSPeter Brune while (next) { 221eed5f15bSPeter Brune ierr = SNESView(next->snes,viewer);CHKERRQ(ierr); 222eed5f15bSPeter Brune next = next->next; 223eed5f15bSPeter Brune } 224eed5f15bSPeter Brune if (iascii) { 225eed5f15bSPeter Brune ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 226eed5f15bSPeter Brune ierr = PetscViewerASCIIPrintf(viewer,"---------------------------------\n");CHKERRQ(ierr); 227eed5f15bSPeter Brune } 228eed5f15bSPeter Brune PetscFunctionReturn(0); 229eed5f15bSPeter Brune } 230eed5f15bSPeter Brune 231eed5f15bSPeter Brune /* ------------------------------------------------------------------------------*/ 232eed5f15bSPeter Brune 233eed5f15bSPeter Brune #undef __FUNCT__ 234eed5f15bSPeter Brune #define __FUNCT__ "SNESCompositeSetType_Composite" 235eed5f15bSPeter Brune static PetscErrorCode SNESCompositeSetType_Composite(SNES snes,SNESCompositeType type) 236eed5f15bSPeter Brune { 237eed5f15bSPeter Brune SNES_Composite *jac = (SNES_Composite*)snes->data; 238eed5f15bSPeter Brune 239eed5f15bSPeter Brune PetscFunctionBegin; 240eed5f15bSPeter Brune jac->type = type; 241eed5f15bSPeter Brune PetscFunctionReturn(0); 242eed5f15bSPeter Brune } 243eed5f15bSPeter Brune 244eed5f15bSPeter Brune #undef __FUNCT__ 245eed5f15bSPeter Brune #define __FUNCT__ "SNESCompositeAddSNES_Composite" 246eed5f15bSPeter Brune static PetscErrorCode SNESCompositeAddSNES_Composite(SNES snes,SNESType type) 247eed5f15bSPeter Brune { 248eed5f15bSPeter Brune SNES_Composite *jac; 249eed5f15bSPeter Brune SNES_CompositeLink next,ilink; 250eed5f15bSPeter Brune PetscErrorCode ierr; 251eed5f15bSPeter Brune PetscInt cnt = 0; 252eed5f15bSPeter Brune const char *prefix; 253eed5f15bSPeter Brune char newprefix[8]; 254eed5f15bSPeter Brune DM dm; 255eed5f15bSPeter Brune 256eed5f15bSPeter Brune PetscFunctionBegin; 257eed5f15bSPeter Brune ierr = PetscNewLog(snes,struct _SNES_CompositeLink,&ilink);CHKERRQ(ierr); 258eed5f15bSPeter Brune ilink->next = 0; 259eed5f15bSPeter Brune ierr = SNESCreate(PetscObjectComm((PetscObject)snes),&ilink->snes);CHKERRQ(ierr); 260eed5f15bSPeter Brune ierr = PetscLogObjectParent((PetscObject)snes,(PetscObject)ilink->snes);CHKERRQ(ierr); 261eed5f15bSPeter Brune ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 262eed5f15bSPeter Brune ierr = SNESSetDM(ilink->snes,dm);CHKERRQ(ierr); 263eed5f15bSPeter Brune 264eed5f15bSPeter Brune jac = (SNES_Composite*)snes->data; 265eed5f15bSPeter Brune next = jac->head; 266eed5f15bSPeter Brune if (!next) { 267eed5f15bSPeter Brune jac->head = ilink; 268eed5f15bSPeter Brune ilink->previous = NULL; 269eed5f15bSPeter Brune } else { 270eed5f15bSPeter Brune cnt++; 271eed5f15bSPeter Brune while (next->next) { 272eed5f15bSPeter Brune next = next->next; 273eed5f15bSPeter Brune cnt++; 274eed5f15bSPeter Brune } 275eed5f15bSPeter Brune next->next = ilink; 276eed5f15bSPeter Brune ilink->previous = next; 277eed5f15bSPeter Brune } 278eed5f15bSPeter Brune ierr = SNESGetOptionsPrefix(snes,&prefix);CHKERRQ(ierr); 279eed5f15bSPeter Brune ierr = SNESSetOptionsPrefix(ilink->snes,prefix);CHKERRQ(ierr); 280eed5f15bSPeter Brune sprintf(newprefix,"sub_%d_",(int)cnt); 281eed5f15bSPeter Brune ierr = SNESAppendOptionsPrefix(ilink->snes,newprefix);CHKERRQ(ierr); 2828f626970SPeter Brune ierr = PetscObjectIncrementTabLevel((PetscObject)ilink->snes,(PetscObject)snes,1);CHKERRQ(ierr); 283eed5f15bSPeter Brune ierr = SNESSetType(ilink->snes,type);CHKERRQ(ierr); 2848f626970SPeter Brune ilink->dmp = 1.0; 285eed5f15bSPeter Brune PetscFunctionReturn(0); 286eed5f15bSPeter Brune } 287eed5f15bSPeter Brune 288eed5f15bSPeter Brune #undef __FUNCT__ 289eed5f15bSPeter Brune #define __FUNCT__ "SNESCompositeGetSNES_Composite" 290eed5f15bSPeter Brune static PetscErrorCode SNESCompositeGetSNES_Composite(SNES snes,PetscInt n,SNES *subsnes) 291eed5f15bSPeter Brune { 292eed5f15bSPeter Brune SNES_Composite *jac; 293eed5f15bSPeter Brune SNES_CompositeLink next; 294eed5f15bSPeter Brune PetscInt i; 295eed5f15bSPeter Brune 296eed5f15bSPeter Brune PetscFunctionBegin; 297eed5f15bSPeter Brune jac = (SNES_Composite*)snes->data; 298eed5f15bSPeter Brune next = jac->head; 299eed5f15bSPeter Brune for (i=0; i<n; i++) { 300eed5f15bSPeter Brune if (!next->next) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_INCOMP,"Not enough SNESes in composite preconditioner"); 301eed5f15bSPeter Brune next = next->next; 302eed5f15bSPeter Brune } 303eed5f15bSPeter Brune *subsnes = next->snes; 304eed5f15bSPeter Brune PetscFunctionReturn(0); 305eed5f15bSPeter Brune } 306eed5f15bSPeter Brune 307eed5f15bSPeter Brune /* -------------------------------------------------------------------------------- */ 308eed5f15bSPeter Brune #undef __FUNCT__ 309eed5f15bSPeter Brune #define __FUNCT__ "SNESCompositeSetType" 310*e31ab4f9SPeter Brune /*@C 311eed5f15bSPeter Brune SNESCompositeSetType - Sets the type of composite preconditioner. 312eed5f15bSPeter Brune 313eed5f15bSPeter Brune Logically Collective on SNES 314eed5f15bSPeter Brune 315eed5f15bSPeter Brune Input Parameter: 316eed5f15bSPeter Brune + snes - the preconditioner context 317eed5f15bSPeter Brune - type - SNES_COMPOSITE_ADDITIVE (default), SNES_COMPOSITE_MULTIPLICATIVE 318eed5f15bSPeter Brune 319eed5f15bSPeter Brune Options Database Key: 320eed5f15bSPeter Brune . -snes_composite_type <type: one of multiplicative, additive, special> - Sets composite preconditioner type 321eed5f15bSPeter Brune 322eed5f15bSPeter Brune Level: Developer 323eed5f15bSPeter Brune 324eed5f15bSPeter Brune .keywords: SNES, set, type, composite preconditioner, additive, multiplicative 325eed5f15bSPeter Brune @*/ 326eed5f15bSPeter Brune PetscErrorCode SNESCompositeSetType(SNES snes,SNESCompositeType type) 327eed5f15bSPeter Brune { 328eed5f15bSPeter Brune PetscErrorCode ierr; 329eed5f15bSPeter Brune 330eed5f15bSPeter Brune PetscFunctionBegin; 331eed5f15bSPeter Brune PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 332eed5f15bSPeter Brune PetscValidLogicalCollectiveEnum(snes,type,2); 333eed5f15bSPeter Brune ierr = PetscTryMethod(snes,"SNESCompositeSetType_C",(SNES,SNESCompositeType),(snes,type));CHKERRQ(ierr); 334eed5f15bSPeter Brune PetscFunctionReturn(0); 335eed5f15bSPeter Brune } 336eed5f15bSPeter Brune 337eed5f15bSPeter Brune #undef __FUNCT__ 338eed5f15bSPeter Brune #define __FUNCT__ "SNESCompositeAddSNES" 339eed5f15bSPeter Brune /*@C 340eed5f15bSPeter Brune SNESCompositeAddSNES - Adds another SNES to the composite SNES. 341eed5f15bSPeter Brune 342eed5f15bSPeter Brune Collective on SNES 343eed5f15bSPeter Brune 344eed5f15bSPeter Brune Input Parameters: 345eed5f15bSPeter Brune + snes - the preconditioner context 346eed5f15bSPeter Brune - type - the type of the new preconditioner 347eed5f15bSPeter Brune 348eed5f15bSPeter Brune Level: Developer 349eed5f15bSPeter Brune 350eed5f15bSPeter Brune .keywords: SNES, composite preconditioner, add 351eed5f15bSPeter Brune @*/ 352eed5f15bSPeter Brune PetscErrorCode SNESCompositeAddSNES(SNES snes,SNESType type) 353eed5f15bSPeter Brune { 354eed5f15bSPeter Brune PetscErrorCode ierr; 355eed5f15bSPeter Brune 356eed5f15bSPeter Brune PetscFunctionBegin; 357eed5f15bSPeter Brune PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 358eed5f15bSPeter Brune ierr = PetscTryMethod(snes,"SNESCompositeAddSNES_C",(SNES,SNESType),(snes,type));CHKERRQ(ierr); 359eed5f15bSPeter Brune PetscFunctionReturn(0); 360eed5f15bSPeter Brune } 361eed5f15bSPeter Brune #undef __FUNCT__ 362eed5f15bSPeter Brune #define __FUNCT__ "SNESCompositeGetSNES" 363eed5f15bSPeter Brune /*@ 364eed5f15bSPeter Brune SNESCompositeGetSNES - Gets one of the SNES objects in the composite SNES. 365eed5f15bSPeter Brune 366eed5f15bSPeter Brune Not Collective 367eed5f15bSPeter Brune 368eed5f15bSPeter Brune Input Parameter: 369eed5f15bSPeter Brune + snes - the preconditioner context 370eed5f15bSPeter Brune - n - the number of the snes requested 371eed5f15bSPeter Brune 372eed5f15bSPeter Brune Output Parameters: 373eed5f15bSPeter Brune . subsnes - the SNES requested 374eed5f15bSPeter Brune 375eed5f15bSPeter Brune Level: Developer 376eed5f15bSPeter Brune 377eed5f15bSPeter Brune .keywords: SNES, get, composite preconditioner, sub preconditioner 378eed5f15bSPeter Brune 379eed5f15bSPeter Brune .seealso: SNESCompositeAddSNES() 380eed5f15bSPeter Brune @*/ 381eed5f15bSPeter Brune PetscErrorCode SNESCompositeGetSNES(SNES snes,PetscInt n,SNES *subsnes) 382eed5f15bSPeter Brune { 383eed5f15bSPeter Brune PetscErrorCode ierr; 384eed5f15bSPeter Brune 385eed5f15bSPeter Brune PetscFunctionBegin; 386eed5f15bSPeter Brune PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 387eed5f15bSPeter Brune PetscValidPointer(subsnes,3); 388eed5f15bSPeter Brune ierr = PetscUseMethod(snes,"SNESCompositeGetSNES_C",(SNES,PetscInt,SNES*),(snes,n,subsnes));CHKERRQ(ierr); 389eed5f15bSPeter Brune PetscFunctionReturn(0); 390eed5f15bSPeter Brune } 391eed5f15bSPeter Brune 392eed5f15bSPeter Brune #undef __FUNCT__ 3938f626970SPeter Brune #define __FUNCT__ "SNESCompositeSetDamping_Composite" 3948f626970SPeter Brune static PetscErrorCode SNESCompositeSetDamping_Composite(SNES snes,PetscInt n,PetscReal dmp) 3958f626970SPeter Brune { 3968f626970SPeter Brune SNES_Composite *jac; 3978f626970SPeter Brune SNES_CompositeLink next; 3988f626970SPeter Brune PetscInt i; 3998f626970SPeter Brune 4008f626970SPeter Brune PetscFunctionBegin; 4018f626970SPeter Brune jac = (SNES_Composite*)snes->data; 4028f626970SPeter Brune next = jac->head; 4038f626970SPeter Brune for (i=0; i<n; i++) { 4048f626970SPeter Brune if (!next->next) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_INCOMP,"Not enough SNESes in composite preconditioner"); 4058f626970SPeter Brune next = next->next; 4068f626970SPeter Brune } 4078f626970SPeter Brune next->dmp = dmp; 4088f626970SPeter Brune PetscFunctionReturn(0); 4098f626970SPeter Brune } 4108f626970SPeter Brune 4118f626970SPeter Brune #undef __FUNCT__ 4128f626970SPeter Brune #define __FUNCT__ "SNESCompositeSetDamping" 4138f626970SPeter Brune /*@ 4148f626970SPeter Brune SNESCompositeSetDamping - Sets the damping of a subsolver when using additive composite SNES. 4158f626970SPeter Brune 4168f626970SPeter Brune Not Collective 4178f626970SPeter Brune 4188f626970SPeter Brune Input Parameter: 4198f626970SPeter Brune + snes - the preconditioner context 4208f626970SPeter Brune . n - the number of the snes requested 4218f626970SPeter Brune - dmp - the damping 4228f626970SPeter Brune 4238f626970SPeter Brune Level: Developer 4248f626970SPeter Brune 4258f626970SPeter Brune .keywords: SNES, get, composite preconditioner, sub preconditioner 4268f626970SPeter Brune 4278f626970SPeter Brune .seealso: SNESCompositeAddSNES() 4288f626970SPeter Brune @*/ 4298f626970SPeter Brune PetscErrorCode SNESCompositeSetDamping(SNES snes,PetscInt n,PetscReal dmp) 4308f626970SPeter Brune { 4318f626970SPeter Brune PetscErrorCode ierr; 4328f626970SPeter Brune 4338f626970SPeter Brune PetscFunctionBegin; 4348f626970SPeter Brune PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 4358f626970SPeter Brune ierr = PetscUseMethod(snes,"SNESCompositeSetDamping_C",(SNES,PetscInt,PetscReal),(snes,n,dmp));CHKERRQ(ierr); 4368f626970SPeter Brune PetscFunctionReturn(0); 4378f626970SPeter Brune } 4388f626970SPeter Brune 4398f626970SPeter Brune #undef __FUNCT__ 440eed5f15bSPeter Brune #define __FUNCT__ "SNESSolve_Composite" 441eed5f15bSPeter Brune PetscErrorCode SNESSolve_Composite(SNES snes) 442eed5f15bSPeter Brune { 443eed5f15bSPeter Brune Vec F; 444eed5f15bSPeter Brune Vec X; 445eed5f15bSPeter Brune Vec B; 446eed5f15bSPeter Brune PetscInt i; 447eed5f15bSPeter Brune PetscReal fnorm = 0.0; 448eed5f15bSPeter Brune PetscErrorCode ierr; 449eed5f15bSPeter Brune SNESNormType normtype; 450eed5f15bSPeter Brune SNES_Composite *comp = (SNES_Composite*)snes->data; 451eed5f15bSPeter Brune 452eed5f15bSPeter Brune PetscFunctionBegin; 453eed5f15bSPeter Brune X = snes->vec_sol; 454eed5f15bSPeter Brune F = snes->vec_func; 455eed5f15bSPeter Brune B = snes->vec_rhs; 456eed5f15bSPeter Brune 457eed5f15bSPeter Brune ierr = PetscObjectAMSTakeAccess((PetscObject)snes);CHKERRQ(ierr); 458eed5f15bSPeter Brune snes->iter = 0; 459eed5f15bSPeter Brune snes->norm = 0.; 460eed5f15bSPeter Brune ierr = PetscObjectAMSGrantAccess((PetscObject)snes);CHKERRQ(ierr); 461eed5f15bSPeter Brune snes->reason = SNES_CONVERGED_ITERATING; 462eed5f15bSPeter Brune ierr = SNESGetNormType(snes, &normtype);CHKERRQ(ierr); 463eed5f15bSPeter Brune if (normtype == SNES_NORM_FUNCTION || normtype == SNES_NORM_INITIAL_ONLY || normtype == SNES_NORM_INITIAL_FINAL_ONLY) { 464eed5f15bSPeter Brune if (!snes->vec_func_init_set) { 465eed5f15bSPeter Brune ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 466eed5f15bSPeter Brune if (snes->domainerror) { 467eed5f15bSPeter Brune snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 468eed5f15bSPeter Brune PetscFunctionReturn(0); 469eed5f15bSPeter Brune } 470eed5f15bSPeter Brune } else snes->vec_func_init_set = PETSC_FALSE; 471eed5f15bSPeter Brune 472eed5f15bSPeter Brune /* convergence test */ 473eed5f15bSPeter Brune if (!snes->norm_init_set) { 474eed5f15bSPeter Brune ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F|| */ 475eed5f15bSPeter Brune if (PetscIsInfOrNanReal(fnorm)) { 476eed5f15bSPeter Brune snes->reason = SNES_DIVERGED_FNORM_NAN; 477eed5f15bSPeter Brune PetscFunctionReturn(0); 478eed5f15bSPeter Brune } 479eed5f15bSPeter Brune } else { 480eed5f15bSPeter Brune fnorm = snes->norm_init; 481eed5f15bSPeter Brune snes->norm_init_set = PETSC_FALSE; 482eed5f15bSPeter Brune } 483eed5f15bSPeter Brune ierr = PetscObjectAMSTakeAccess((PetscObject)snes);CHKERRQ(ierr); 484eed5f15bSPeter Brune snes->iter = 0; 485eed5f15bSPeter Brune snes->norm = fnorm; 486eed5f15bSPeter Brune ierr = PetscObjectAMSGrantAccess((PetscObject)snes);CHKERRQ(ierr); 487eed5f15bSPeter Brune ierr = SNESLogConvergenceHistory(snes,snes->norm,0);CHKERRQ(ierr); 488eed5f15bSPeter Brune ierr = SNESMonitor(snes,0,snes->norm);CHKERRQ(ierr); 489eed5f15bSPeter Brune 490eed5f15bSPeter Brune /* set parameter for default relative tolerance convergence test */ 491eed5f15bSPeter Brune snes->ttol = fnorm*snes->rtol; 492eed5f15bSPeter Brune 493eed5f15bSPeter Brune /* test convergence */ 494eed5f15bSPeter Brune ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 495eed5f15bSPeter Brune if (snes->reason) PetscFunctionReturn(0); 496eed5f15bSPeter Brune } else { 497eed5f15bSPeter Brune ierr = PetscObjectAMSGrantAccess((PetscObject)snes);CHKERRQ(ierr); 498eed5f15bSPeter Brune ierr = SNESLogConvergenceHistory(snes,snes->norm,0);CHKERRQ(ierr); 499eed5f15bSPeter Brune ierr = SNESMonitor(snes,0,snes->norm);CHKERRQ(ierr); 500eed5f15bSPeter Brune } 501eed5f15bSPeter Brune 502eed5f15bSPeter Brune /* Call general purpose update function */ 503eed5f15bSPeter Brune if (snes->ops->update) { 504eed5f15bSPeter Brune ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 505eed5f15bSPeter Brune } 506eed5f15bSPeter Brune 507eed5f15bSPeter Brune for (i = 0; i < snes->max_its; i++) { 508eed5f15bSPeter Brune if (comp->type == SNES_COMPOSITE_ADDITIVE) { 509eed5f15bSPeter Brune ierr = SNESCompositeApply_Additive(snes,X,B,F,&fnorm);CHKERRQ(ierr); 510eed5f15bSPeter Brune } else if (comp->type == SNES_COMPOSITE_MULTIPLICATIVE) { 511eed5f15bSPeter Brune ierr = SNESCompositeApply_Multiplicative(snes,X,B,F,&fnorm);CHKERRQ(ierr); 512eed5f15bSPeter Brune } else { 513eed5f15bSPeter Brune } 514eed5f15bSPeter Brune if ((i == snes->max_its - 1) && (normtype == SNES_NORM_INITIAL_FINAL_ONLY || normtype == SNES_NORM_FINAL_ONLY)) { 515eed5f15bSPeter Brune ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 516eed5f15bSPeter Brune if (snes->domainerror) { 517eed5f15bSPeter Brune snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 518eed5f15bSPeter Brune break; 519eed5f15bSPeter Brune } 520eed5f15bSPeter Brune ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F|| */ 521eed5f15bSPeter Brune if (PetscIsInfOrNanReal(fnorm)) { 522eed5f15bSPeter Brune snes->reason = SNES_DIVERGED_FNORM_NAN; 523eed5f15bSPeter Brune break; 524eed5f15bSPeter Brune } 525eed5f15bSPeter Brune } 526eed5f15bSPeter Brune /* Monitor convergence */ 527eed5f15bSPeter Brune ierr = PetscObjectAMSTakeAccess((PetscObject)snes);CHKERRQ(ierr); 528eed5f15bSPeter Brune snes->iter = i+1; 529eed5f15bSPeter Brune snes->norm = fnorm; 530eed5f15bSPeter Brune ierr = PetscObjectAMSGrantAccess((PetscObject)snes);CHKERRQ(ierr); 531eed5f15bSPeter Brune ierr = SNESLogConvergenceHistory(snes,snes->norm,0);CHKERRQ(ierr); 532eed5f15bSPeter Brune ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr); 533eed5f15bSPeter Brune /* Test for convergence */ 534eed5f15bSPeter Brune if (normtype == SNES_NORM_FUNCTION) {ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);} 535eed5f15bSPeter Brune if (snes->reason) break; 536eed5f15bSPeter Brune /* Call general purpose update function */ 537eed5f15bSPeter Brune if (snes->ops->update) {ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);} 538eed5f15bSPeter Brune } 539eed5f15bSPeter Brune if (normtype == SNES_NORM_FUNCTION) { 540eed5f15bSPeter Brune if (i == snes->max_its) { 541eed5f15bSPeter Brune ierr = PetscInfo1(snes,"Maximum number of iterations has been reached: %D\n",snes->max_its);CHKERRQ(ierr); 542eed5f15bSPeter Brune if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT; 543eed5f15bSPeter Brune } 544eed5f15bSPeter Brune } else if (!snes->reason) snes->reason = SNES_CONVERGED_ITS; 545eed5f15bSPeter Brune PetscFunctionReturn(0); 546eed5f15bSPeter Brune } 547eed5f15bSPeter Brune 548eed5f15bSPeter Brune /* -------------------------------------------------------------------------------------------*/ 549eed5f15bSPeter Brune 550eed5f15bSPeter Brune /*MC 551eed5f15bSPeter Brune SNESCOMPOSITE - Build a preconditioner by composing together several nonlinear solvers 552eed5f15bSPeter Brune 553eed5f15bSPeter Brune Options Database Keys: 554eed5f15bSPeter Brune + -snes_composite_type <type: one of multiplicative, additive, symmetric_multiplicative, special> - Sets composite preconditioner type 555eed5f15bSPeter Brune - -snes_composite_sneses - <snes0,snes1,...> list of SNESes to compose 556eed5f15bSPeter Brune 557eed5f15bSPeter Brune Level: intermediate 558eed5f15bSPeter Brune 559eed5f15bSPeter Brune Concepts: composing solvers 560eed5f15bSPeter Brune 561eed5f15bSPeter Brune .seealso: SNESCreate(), SNESSetType(), SNESType (for list of available types), SNES, 562eed5f15bSPeter Brune SNESSHELL, SNESCompositeSetType(), SNESCompositeSpecialSetAlpha(), SNESCompositeAddSNES(), 563eed5f15bSPeter Brune SNESCompositeGetSNES() 564eed5f15bSPeter Brune 565eed5f15bSPeter Brune M*/ 566eed5f15bSPeter Brune 567eed5f15bSPeter Brune #undef __FUNCT__ 568eed5f15bSPeter Brune #define __FUNCT__ "SNESCreate_Composite" 569eed5f15bSPeter Brune PETSC_EXTERN PetscErrorCode SNESCreate_Composite(SNES snes) 570eed5f15bSPeter Brune { 571eed5f15bSPeter Brune PetscErrorCode ierr; 572eed5f15bSPeter Brune SNES_Composite *jac; 573eed5f15bSPeter Brune 574eed5f15bSPeter Brune PetscFunctionBegin; 575eed5f15bSPeter Brune ierr = PetscNewLog(snes,SNES_Composite,&jac);CHKERRQ(ierr); 576eed5f15bSPeter Brune 577eed5f15bSPeter Brune snes->ops->solve = SNESSolve_Composite; 578eed5f15bSPeter Brune snes->ops->setup = SNESSetUp_Composite; 579eed5f15bSPeter Brune snes->ops->reset = SNESReset_Composite; 580eed5f15bSPeter Brune snes->ops->destroy = SNESDestroy_Composite; 581eed5f15bSPeter Brune snes->ops->setfromoptions = SNESSetFromOptions_Composite; 582eed5f15bSPeter Brune snes->ops->view = SNESView_Composite; 583eed5f15bSPeter Brune 584eed5f15bSPeter Brune snes->data = (void*)jac; 585eed5f15bSPeter Brune jac->type = SNES_COMPOSITE_ADDITIVE; 586eed5f15bSPeter Brune jac->head = 0; 587eed5f15bSPeter Brune 5888f626970SPeter Brune ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESCompositeSetType_C",SNESCompositeSetType_Composite);CHKERRQ(ierr); 5898f626970SPeter Brune ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESCompositeAddSNES_C",SNESCompositeAddSNES_Composite);CHKERRQ(ierr); 5908f626970SPeter Brune ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESCompositeGetSNES_C",SNESCompositeGetSNES_Composite);CHKERRQ(ierr); 5918f626970SPeter Brune ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESCompositeSetDamping_C",SNESCompositeSetDamping_Composite);CHKERRQ(ierr); 592eed5f15bSPeter Brune PetscFunctionReturn(0); 593eed5f15bSPeter Brune } 594eed5f15bSPeter Brune 595