1eaedb033SPeter Brune #include <petsc-private/snesimpl.h> /*I "petscsnes.h" I*/ 2111ade9eSPeter Brune #include <petscdm.h> 3eaedb033SPeter Brune 4eaedb033SPeter Brune typedef struct { 5eaedb033SPeter Brune PetscInt n; /* local subdomains */ 6eaedb033SPeter Brune SNES *subsnes; /* nonlinear solvers for each subdomain */ 7eaedb033SPeter Brune Vec *x; /* solution vectors */ 8111ade9eSPeter Brune Vec *xl; /* solution local vectors */ 9111ade9eSPeter Brune Vec *y; /* step vectors */ 10eaedb033SPeter Brune Vec *b; /* rhs vectors */ 11111ade9eSPeter Brune VecScatter *oscatter; /* scatter from global space to the subdomain global space */ 12111ade9eSPeter Brune VecScatter *iscatter; /* scatter from global space to the nonoverlapping subdomain space */ 13111ade9eSPeter Brune VecScatter *gscatter; /* scatter from global space to the subdomain local space */ 14111ade9eSPeter Brune PCASMType type; /* ASM type */ 15111ade9eSPeter Brune PetscBool usesdm; /* use the DM for setting up the subproblems */ 16d728fb7dSPeter Brune PetscBool finaljacobian; /* compute the jacobian of the converged solution */ 17610116beSPeter Brune PetscReal damping; /* damping parameter for updates from the blocks */ 18a4f17876SPeter Brune PetscBool same_local_solves; /* flag to determine if the solvers have been individually modified */ 19b20c023fSPeter Brune 20b20c023fSPeter Brune /* logging events */ 21b20c023fSPeter Brune PetscLogEvent eventrestrictinterp; 22b20c023fSPeter Brune PetscLogEvent eventsubsolve; 23602bec5dSPeter Brune 24602bec5dSPeter Brune PetscInt fjtype; /* type of computed jacobian */ 25602bec5dSPeter Brune Vec xinit; /* initial solution in case the final jacobian type is computed as first */ 26eaedb033SPeter Brune } SNES_NASM; 27eaedb033SPeter Brune 28b20c023fSPeter Brune const char *const SNESNASMTypes[] = {"NONE","RESTRICT","INTERPOLATE","BASIC","PCASMType","PC_ASM_",0}; 29602bec5dSPeter Brune const char *const SNESNASMFJTypes[] = {"FINALOUTER","FINALINNER","INITIAL"}; 30b20c023fSPeter Brune 31eaedb033SPeter Brune #undef __FUNCT__ 32eaedb033SPeter Brune #define __FUNCT__ "SNESReset_NASM" 33eaedb033SPeter Brune PetscErrorCode SNESReset_NASM(SNES snes) 34eaedb033SPeter Brune { 35eaedb033SPeter Brune SNES_NASM *nasm = (SNES_NASM*)snes->data; 36eaedb033SPeter Brune PetscErrorCode ierr; 37eaedb033SPeter Brune PetscInt i; 386e111a19SKarl Rupp 39eaedb033SPeter Brune PetscFunctionBegin; 40eaedb033SPeter Brune for (i=0; i<nasm->n; i++) { 41111ade9eSPeter Brune if (nasm->xl) { ierr = VecDestroy(&nasm->xl[i]);CHKERRQ(ierr); } 42f5f7c1b9SKarl Rupp if (nasm->x) { ierr = VecDestroy(&nasm->x[i]);CHKERRQ(ierr); } 43111ade9eSPeter Brune if (nasm->y) { ierr = VecDestroy(&nasm->y[i]);CHKERRQ(ierr); } 44bc8c1f72SJose Roman if (nasm->b) { ierr = VecDestroy(&nasm->b[i]);CHKERRQ(ierr); } 45eaedb033SPeter Brune 46bc8c1f72SJose Roman if (nasm->subsnes) { ierr = SNESDestroy(&nasm->subsnes[i]);CHKERRQ(ierr); } 47111ade9eSPeter Brune if (nasm->oscatter) { ierr = VecScatterDestroy(&nasm->oscatter[i]);CHKERRQ(ierr); } 48111ade9eSPeter Brune if (nasm->iscatter) { ierr = VecScatterDestroy(&nasm->iscatter[i]);CHKERRQ(ierr); } 49111ade9eSPeter Brune if (nasm->gscatter) { ierr = VecScatterDestroy(&nasm->gscatter[i]);CHKERRQ(ierr); } 50eaedb033SPeter Brune } 51111ade9eSPeter Brune 52111ade9eSPeter Brune if (nasm->x) {ierr = PetscFree(nasm->x);CHKERRQ(ierr);} 53111ade9eSPeter Brune if (nasm->xl) {ierr = PetscFree(nasm->xl);CHKERRQ(ierr);} 54111ade9eSPeter Brune if (nasm->y) {ierr = PetscFree(nasm->y);CHKERRQ(ierr);} 55111ade9eSPeter Brune if (nasm->b) {ierr = PetscFree(nasm->b);CHKERRQ(ierr);} 56111ade9eSPeter Brune 57602bec5dSPeter Brune if (nasm->xinit) {ierr = VecDestroy(&nasm->xinit);CHKERRQ(ierr);} 58602bec5dSPeter Brune 59111ade9eSPeter Brune if (nasm->subsnes) {ierr = PetscFree(nasm->subsnes);CHKERRQ(ierr);} 60111ade9eSPeter Brune if (nasm->oscatter) {ierr = PetscFree(nasm->oscatter);CHKERRQ(ierr);} 61111ade9eSPeter Brune if (nasm->iscatter) {ierr = PetscFree(nasm->iscatter);CHKERRQ(ierr);} 62111ade9eSPeter Brune if (nasm->gscatter) {ierr = PetscFree(nasm->gscatter);CHKERRQ(ierr);} 63b20c023fSPeter Brune 64b20c023fSPeter Brune nasm->eventrestrictinterp = 0; 65b20c023fSPeter Brune nasm->eventsubsolve = 0; 66eaedb033SPeter Brune PetscFunctionReturn(0); 67eaedb033SPeter Brune } 68eaedb033SPeter Brune 69eaedb033SPeter Brune #undef __FUNCT__ 70eaedb033SPeter Brune #define __FUNCT__ "SNESDestroy_NASM" 71eaedb033SPeter Brune PetscErrorCode SNESDestroy_NASM(SNES snes) 72eaedb033SPeter Brune { 73eaedb033SPeter Brune PetscErrorCode ierr; 746e111a19SKarl Rupp 75eaedb033SPeter Brune PetscFunctionBegin; 76eaedb033SPeter Brune ierr = SNESReset_NASM(snes);CHKERRQ(ierr); 7722d28d08SBarry Smith ierr = PetscFree(snes->data);CHKERRQ(ierr); 78eaedb033SPeter Brune PetscFunctionReturn(0); 79eaedb033SPeter Brune } 80eaedb033SPeter Brune 81eaedb033SPeter Brune #undef __FUNCT__ 82111ade9eSPeter Brune #define __FUNCT__ "DMGlobalToLocalSubDomainDirichletHook_Private" 830adebc6cSBarry Smith PetscErrorCode DMGlobalToLocalSubDomainDirichletHook_Private(DM dm,Vec g,InsertMode mode,Vec l,void *ctx) 840adebc6cSBarry Smith { 85111ade9eSPeter Brune PetscErrorCode ierr; 86111ade9eSPeter Brune Vec bcs = (Vec)ctx; 876e111a19SKarl Rupp 88111ade9eSPeter Brune PetscFunctionBegin; 89111ade9eSPeter Brune ierr = VecCopy(bcs,l);CHKERRQ(ierr); 90111ade9eSPeter Brune PetscFunctionReturn(0); 91111ade9eSPeter Brune } 92111ade9eSPeter Brune 93111ade9eSPeter Brune #undef __FUNCT__ 94eaedb033SPeter Brune #define __FUNCT__ "SNESSetUp_NASM" 95eaedb033SPeter Brune PetscErrorCode SNESSetUp_NASM(SNES snes) 96eaedb033SPeter Brune { 97eaedb033SPeter Brune SNES_NASM *nasm = (SNES_NASM*)snes->data; 98eaedb033SPeter Brune PetscErrorCode ierr; 9976857b2aSPeter Brune DM dm,subdm; 100111ade9eSPeter Brune DM *subdms; 101111ade9eSPeter Brune PetscInt i; 102eaedb033SPeter Brune const char *optionsprefix; 103111ade9eSPeter Brune Vec F; 104ed3c11a9SPeter Brune PetscMPIInt size; 105ed3c11a9SPeter Brune KSP ksp; 106ed3c11a9SPeter Brune PC pc; 107eaedb033SPeter Brune 108eaedb033SPeter Brune PetscFunctionBegin; 109eaedb033SPeter Brune if (!nasm->subsnes) { 110eaedb033SPeter Brune ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 1110a696f66SPeter Brune if (dm) { 112eaedb033SPeter Brune nasm->usesdm = PETSC_TRUE; 1130298fd71SBarry Smith ierr = DMCreateDomainDecomposition(dm,&nasm->n,NULL,NULL,NULL,&subdms);CHKERRQ(ierr); 114ce94432eSBarry Smith if (!subdms) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"DM has no default decomposition defined. Set subsolves manually with SNESNASMSetSubdomains()."); 115111ade9eSPeter Brune ierr = DMCreateDomainDecompositionScatters(dm,nasm->n,subdms,&nasm->iscatter,&nasm->oscatter,&nasm->gscatter);CHKERRQ(ierr); 116eaedb033SPeter Brune 117eaedb033SPeter Brune ierr = SNESGetOptionsPrefix(snes, &optionsprefix);CHKERRQ(ierr); 118111ade9eSPeter Brune ierr = PetscMalloc(nasm->n*sizeof(SNES),&nasm->subsnes);CHKERRQ(ierr); 119111ade9eSPeter Brune for (i=0; i<nasm->n; i++) { 120cdb298fcSPeter Brune ierr = SNESCreate(PETSC_COMM_SELF,&nasm->subsnes[i]);CHKERRQ(ierr); 121cdb298fcSPeter Brune ierr = SNESAppendOptionsPrefix(nasm->subsnes[i],optionsprefix);CHKERRQ(ierr); 122cdb298fcSPeter Brune ierr = SNESAppendOptionsPrefix(nasm->subsnes[i],"sub_");CHKERRQ(ierr); 123cdb298fcSPeter Brune ierr = SNESSetDM(nasm->subsnes[i],subdms[i]);CHKERRQ(ierr); 124ed3c11a9SPeter Brune ierr = MPI_Comm_size(PetscObjectComm((PetscObject)nasm->subsnes[i]),&size);CHKERRQ(ierr); 125ed3c11a9SPeter Brune if (size == 1) { 126ed3c11a9SPeter Brune ierr = SNESGetKSP(nasm->subsnes[i],&ksp);CHKERRQ(ierr); 127ed3c11a9SPeter Brune ierr = KSPGetPC(ksp,&pc);CHKERRQ(ierr); 128ed3c11a9SPeter Brune ierr = KSPSetType(ksp,KSPPREONLY);CHKERRQ(ierr); 129ed3c11a9SPeter Brune ierr = PCSetType(pc,PCLU);CHKERRQ(ierr); 130ed3c11a9SPeter Brune } 131cdb298fcSPeter Brune ierr = SNESSetFromOptions(nasm->subsnes[i]);CHKERRQ(ierr); 132111ade9eSPeter Brune ierr = DMDestroy(&subdms[i]);CHKERRQ(ierr); 133111ade9eSPeter Brune } 134111ade9eSPeter Brune ierr = PetscFree(subdms);CHKERRQ(ierr); 135ce94432eSBarry Smith } else SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"Cannot construct local problems automatically without a DM!"); 136ce94432eSBarry Smith } else SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"Must set subproblems manually if there is no DM!"); 137111ade9eSPeter Brune /* allocate the global vectors */ 138e245e0aaSPeter Brune if (!nasm->x) { 139111ade9eSPeter Brune ierr = PetscMalloc(nasm->n*sizeof(Vec),&nasm->x);CHKERRQ(ierr); 140e245e0aaSPeter Brune ierr = PetscMemzero(nasm->x,nasm->n*sizeof(Vec));CHKERRQ(ierr); 141e245e0aaSPeter Brune } 142e245e0aaSPeter Brune if (!nasm->xl) { 143111ade9eSPeter Brune ierr = PetscMalloc(nasm->n*sizeof(Vec),&nasm->xl);CHKERRQ(ierr); 144e245e0aaSPeter Brune ierr = PetscMemzero(nasm->xl,nasm->n*sizeof(Vec));CHKERRQ(ierr); 145e245e0aaSPeter Brune } 146e245e0aaSPeter Brune if (!nasm->y) { 147111ade9eSPeter Brune ierr = PetscMalloc(nasm->n*sizeof(Vec),&nasm->y);CHKERRQ(ierr); 148e245e0aaSPeter Brune ierr = PetscMemzero(nasm->y,nasm->n*sizeof(Vec));CHKERRQ(ierr); 149e245e0aaSPeter Brune } 150e245e0aaSPeter Brune if (!nasm->b) { 151111ade9eSPeter Brune ierr = PetscMalloc(nasm->n*sizeof(Vec),&nasm->b);CHKERRQ(ierr); 152e245e0aaSPeter Brune ierr = PetscMemzero(nasm->b,nasm->n*sizeof(Vec));CHKERRQ(ierr); 153e245e0aaSPeter Brune } 154111ade9eSPeter Brune 155111ade9eSPeter Brune for (i=0; i<nasm->n; i++) { 1560298fd71SBarry Smith ierr = SNESGetFunction(nasm->subsnes[i],&F,NULL,NULL);CHKERRQ(ierr); 15776857b2aSPeter Brune if (!nasm->x[i]) {ierr = VecDuplicate(F,&nasm->x[i]);CHKERRQ(ierr);} 15876857b2aSPeter Brune if (!nasm->y[i]) {ierr = VecDuplicate(F,&nasm->y[i]);CHKERRQ(ierr);} 15976857b2aSPeter Brune if (!nasm->b[i]) {ierr = VecDuplicate(F,&nasm->b[i]);CHKERRQ(ierr);} 16076857b2aSPeter Brune if (!nasm->xl[i]) { 161111ade9eSPeter Brune ierr = SNESGetDM(nasm->subsnes[i],&subdm);CHKERRQ(ierr); 162111ade9eSPeter Brune ierr = DMCreateLocalVector(subdm,&nasm->xl[i]);CHKERRQ(ierr); 16376857b2aSPeter Brune } 1640298fd71SBarry Smith ierr = DMGlobalToLocalHookAdd(subdm,DMGlobalToLocalSubDomainDirichletHook_Private,NULL,nasm->xl[i]);CHKERRQ(ierr); 165111ade9eSPeter Brune } 166602bec5dSPeter Brune if (nasm->finaljacobian) { 167602bec5dSPeter Brune ierr = SNESSetUpMatrices(snes);CHKERRQ(ierr); 168602bec5dSPeter Brune if (nasm->fjtype == 2) { 169602bec5dSPeter Brune ierr = VecDuplicate(snes->vec_sol,&nasm->xinit);CHKERRQ(ierr); 170602bec5dSPeter Brune } 171602bec5dSPeter Brune for (i=0; i<nasm->n;i++) { 172602bec5dSPeter Brune ierr = SNESSetUpMatrices(nasm->subsnes[i]);CHKERRQ(ierr); 173602bec5dSPeter Brune } 174602bec5dSPeter Brune } 175eaedb033SPeter Brune PetscFunctionReturn(0); 176eaedb033SPeter Brune } 177eaedb033SPeter Brune 178eaedb033SPeter Brune #undef __FUNCT__ 179eaedb033SPeter Brune #define __FUNCT__ "SNESSetFromOptions_NASM" 180eaedb033SPeter Brune PetscErrorCode SNESSetFromOptions_NASM(SNES snes) 181eaedb033SPeter Brune { 182eaedb033SPeter Brune PetscErrorCode ierr; 183111ade9eSPeter Brune PCASMType asmtype; 184a4f17876SPeter Brune PetscBool flg,monflg,subviewflg; 185111ade9eSPeter Brune SNES_NASM *nasm = (SNES_NASM*)snes->data; 1866e111a19SKarl Rupp 187eaedb033SPeter Brune PetscFunctionBegin; 188111ade9eSPeter Brune ierr = PetscOptionsHead("Nonlinear Additive Schwartz options");CHKERRQ(ierr); 189111ade9eSPeter Brune ierr = PetscOptionsEnum("-snes_nasm_type","Type of restriction/extension","",SNESNASMTypes,(PetscEnum)nasm->type,(PetscEnum*)&asmtype,&flg);CHKERRQ(ierr); 1901aa26658SKarl Rupp if (flg) nasm->type = asmtype; 191b20c023fSPeter Brune flg = PETSC_FALSE; 192b20c023fSPeter Brune monflg = PETSC_TRUE; 193610116beSPeter Brune ierr = PetscOptionsReal("-snes_nasm_damping","Log times for subSNES solves and restriction","SNESNASMSetDamping",nasm->damping,&nasm->damping,&flg);CHKERRQ(ierr); 194610116beSPeter Brune if (flg) {ierr = SNESNASMSetDamping(snes,nasm->damping);CHKERRQ(ierr);} 195a4f17876SPeter Brune subviewflg = PETSC_FALSE; 196a4f17876SPeter Brune ierr = PetscOptionsBool("-snes_nasm_sub_view","Print detailed information for every processor when using -snes_view","",subviewflg,&subviewflg,&flg);CHKERRQ(ierr); 197a4f17876SPeter Brune if (flg) { 198a4f17876SPeter Brune nasm->same_local_solves = PETSC_FALSE; 199a4f17876SPeter Brune if (!subviewflg) { 200a4f17876SPeter Brune nasm->same_local_solves = PETSC_TRUE; 201a4f17876SPeter Brune } 202a4f17876SPeter Brune } 203d728fb7dSPeter Brune ierr = PetscOptionsBool("-snes_nasm_finaljacobian","Compute the global jacobian of the final iterate (for ASPIN)","",nasm->finaljacobian,&nasm->finaljacobian,NULL);CHKERRQ(ierr); 204602bec5dSPeter Brune ierr = PetscOptionsEList("-snes_nasm_finaljacobian_type","The type of the final jacobian computed.","",SNESNASMFJTypes,3,SNESNASMFJTypes[0],&nasm->fjtype,NULL);CHKERRQ(ierr); 205a4f17876SPeter Brune ierr = PetscOptionsBool("-snes_nasm_log","Log times for subSNES solves and restriction","",monflg,&monflg,&flg);CHKERRQ(ierr); 206b20c023fSPeter Brune if (flg) { 207b20c023fSPeter Brune ierr = PetscLogEventRegister("SNESNASMSubSolve",((PetscObject)snes)->classid,&nasm->eventsubsolve);CHKERRQ(ierr); 208b20c023fSPeter Brune ierr = PetscLogEventRegister("SNESNASMRestrict",((PetscObject)snes)->classid,&nasm->eventrestrictinterp);CHKERRQ(ierr); 209b20c023fSPeter Brune } 210eaedb033SPeter Brune ierr = PetscOptionsTail();CHKERRQ(ierr); 211eaedb033SPeter Brune PetscFunctionReturn(0); 212eaedb033SPeter Brune } 213eaedb033SPeter Brune 214eaedb033SPeter Brune #undef __FUNCT__ 215eaedb033SPeter Brune #define __FUNCT__ "SNESView_NASM" 216eaedb033SPeter Brune PetscErrorCode SNESView_NASM(SNES snes, PetscViewer viewer) 217eaedb033SPeter Brune { 218b20c023fSPeter Brune SNES_NASM *nasm = (SNES_NASM*)snes->data; 219b20c023fSPeter Brune PetscErrorCode ierr; 220a4f17876SPeter Brune PetscMPIInt rank,size; 22179f8d099SPeter Brune PetscInt i,j,N,bsz; 222b20c023fSPeter Brune PetscBool iascii,isstring; 223b20c023fSPeter Brune PetscViewer sviewer; 224ce94432eSBarry Smith MPI_Comm comm; 225b20c023fSPeter Brune 226b20c023fSPeter Brune PetscFunctionBegin; 227ce94432eSBarry Smith ierr = PetscObjectGetComm((PetscObject)snes,&comm);CHKERRQ(ierr); 228b20c023fSPeter Brune ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 229b20c023fSPeter Brune ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring);CHKERRQ(ierr); 230b20c023fSPeter Brune ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 231a4f17876SPeter Brune ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 23279f8d099SPeter Brune ierr = MPI_Allreduce(&nasm->n,&N,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(ierr); 233b20c023fSPeter Brune if (iascii) { 234b20c023fSPeter Brune ierr = PetscViewerASCIIPrintf(viewer, " Nonlinear Additive Schwarz: total subdomain blocks = %D\n",N);CHKERRQ(ierr); 235a4f17876SPeter Brune if (nasm->same_local_solves) { 236a4f17876SPeter Brune if (nasm->subsnes) { 237a4f17876SPeter Brune ierr = PetscViewerASCIIPrintf(viewer," Local solve is the same for all blocks:\n");CHKERRQ(ierr); 238a4f17876SPeter Brune ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 239a4f17876SPeter Brune ierr = PetscViewerGetSingleton(viewer,&sviewer);CHKERRQ(ierr); 240a4f17876SPeter Brune if (!rank) { 241a4f17876SPeter Brune ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 242a4f17876SPeter Brune ierr = SNESView(nasm->subsnes[0],sviewer);CHKERRQ(ierr); 243a4f17876SPeter Brune ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 244a4f17876SPeter Brune } 245a4f17876SPeter Brune ierr = PetscViewerRestoreSingleton(viewer,&sviewer);CHKERRQ(ierr); 246a4f17876SPeter Brune ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 247a4f17876SPeter Brune } 248a4f17876SPeter Brune } else { 249a4f17876SPeter Brune /* print the solver on each block */ 250b20c023fSPeter Brune ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_TRUE);CHKERRQ(ierr); 251b20c023fSPeter Brune ierr = PetscViewerASCIISynchronizedPrintf(viewer," [%d] number of local blocks = %D\n",(int)rank,nasm->n);CHKERRQ(ierr); 252b20c023fSPeter Brune ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 253a4f17876SPeter Brune ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_FALSE);CHKERRQ(ierr); 254a4f17876SPeter Brune ierr = PetscViewerASCIIPrintf(viewer," Local solve info for each block is in the following SNES objects:\n");CHKERRQ(ierr); 255b20c023fSPeter Brune ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 256a4f17876SPeter Brune ierr = PetscViewerASCIIPrintf(viewer,"- - - - - - - - - - - - - - - - - -\n");CHKERRQ(ierr); 257a4f17876SPeter Brune for (j=0; j<size; j++) { 258b20c023fSPeter Brune ierr = PetscViewerGetSingleton(viewer,&sviewer);CHKERRQ(ierr); 259a4f17876SPeter Brune if (rank == j) { 260a4f17876SPeter Brune for (i=0; i<nasm->n; i++) { 261a4f17876SPeter Brune ierr = VecGetLocalSize(nasm->x[i],&bsz);CHKERRQ(ierr); 262a4f17876SPeter Brune ierr = PetscViewerASCIIPrintf(sviewer,"[%d] local block number %D, size = %D\n",(int)rank,i,bsz);CHKERRQ(ierr); 263b20c023fSPeter Brune ierr = SNESView(nasm->subsnes[i],sviewer);CHKERRQ(ierr); 264a4f17876SPeter Brune ierr = PetscViewerASCIIPrintf(sviewer,"- - - - - - - - - - - - - - - - - -\n");CHKERRQ(ierr); 265b20c023fSPeter Brune } 266b20c023fSPeter Brune } 267a4f17876SPeter Brune ierr = PetscViewerRestoreSingleton(viewer,&sviewer);CHKERRQ(ierr); 268a4f17876SPeter Brune ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 269a4f17876SPeter Brune } 270b20c023fSPeter Brune ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 271a4f17876SPeter Brune } 272b20c023fSPeter Brune } else if (isstring) { 273b20c023fSPeter Brune ierr = PetscViewerStringSPrintf(viewer," blocks=%D,type=%s",N,SNESNASMTypes[nasm->type]);CHKERRQ(ierr); 274b20c023fSPeter Brune ierr = PetscViewerGetSingleton(viewer,&sviewer);CHKERRQ(ierr); 275b20c023fSPeter Brune if (nasm->subsnes && !rank) {ierr = SNESView(nasm->subsnes[0],sviewer);CHKERRQ(ierr);} 276b20c023fSPeter Brune ierr = PetscViewerRestoreSingleton(viewer,&sviewer);CHKERRQ(ierr); 277b20c023fSPeter Brune } 278eaedb033SPeter Brune PetscFunctionReturn(0); 279eaedb033SPeter Brune } 280eaedb033SPeter Brune 281eaedb033SPeter Brune #undef __FUNCT__ 282eaedb033SPeter Brune #define __FUNCT__ "SNESNASMSetSubdomains" 28376857b2aSPeter Brune /*@ 28476857b2aSPeter Brune SNESNASMSetSubdomains - Manually Set the context required to restrict and solve subdomain problems. 28576857b2aSPeter Brune 28676857b2aSPeter Brune Not Collective 28776857b2aSPeter Brune 28876857b2aSPeter Brune Input Parameters: 28976857b2aSPeter Brune + SNES - the SNES context 29076857b2aSPeter Brune . n - the number of local subdomains 29176857b2aSPeter Brune . subsnes - solvers defined on the local subdomains 29276857b2aSPeter Brune . iscatter - scatters into the nonoverlapping portions of the local subdomains 29376857b2aSPeter Brune . oscatter - scatters into the overlapping portions of the local subdomains 29476857b2aSPeter Brune - gscatter - scatters into the (ghosted) local vector of the local subdomain 29576857b2aSPeter Brune 29676857b2aSPeter Brune Level: intermediate 29776857b2aSPeter Brune 29876857b2aSPeter Brune .keywords: SNES, NASM 29976857b2aSPeter Brune 30076857b2aSPeter Brune .seealso: SNESNASM, SNESNASMGetSubdomains() 30176857b2aSPeter Brune @*/ 302a6dfd86eSKarl Rupp PetscErrorCode SNESNASMSetSubdomains(SNES snes,PetscInt n,SNES subsnes[],VecScatter iscatter[],VecScatter oscatter[],VecScatter gscatter[]) 303a6dfd86eSKarl Rupp { 304eaedb033SPeter Brune PetscErrorCode ierr; 305111ade9eSPeter Brune PetscErrorCode (*f)(SNES,PetscInt,SNES*,VecScatter*,VecScatter*,VecScatter*); 3066e111a19SKarl Rupp 307eaedb033SPeter Brune PetscFunctionBegin; 3080005d66cSJed Brown ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESNASMSetSubdomains_C",&f);CHKERRQ(ierr); 3094b838e8fSPeter Brune if (f) {ierr = (f)(snes,n,subsnes,iscatter,oscatter,gscatter);CHKERRQ(ierr);} 310eaedb033SPeter Brune PetscFunctionReturn(0); 311eaedb033SPeter Brune } 312eaedb033SPeter Brune 313eaedb033SPeter Brune #undef __FUNCT__ 314eaedb033SPeter Brune #define __FUNCT__ "SNESNASMSetSubdomains_NASM" 315a6dfd86eSKarl Rupp PetscErrorCode SNESNASMSetSubdomains_NASM(SNES snes,PetscInt n,SNES subsnes[],VecScatter iscatter[],VecScatter oscatter[],VecScatter gscatter[]) 316a6dfd86eSKarl Rupp { 317eaedb033SPeter Brune PetscInt i; 318eaedb033SPeter Brune PetscErrorCode ierr; 319eaedb033SPeter Brune SNES_NASM *nasm = (SNES_NASM*)snes->data; 3206e111a19SKarl Rupp 321eaedb033SPeter Brune PetscFunctionBegin; 322ce94432eSBarry Smith if (snes->setupcalled) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"SNESNASMSetSubdomains() should be called before calling SNESSetUp()."); 323eaedb033SPeter Brune 324111ade9eSPeter Brune /* tear down the previously set things */ 325111ade9eSPeter Brune ierr = SNESReset(snes);CHKERRQ(ierr); 326111ade9eSPeter Brune 327eaedb033SPeter Brune nasm->n = n; 328111ade9eSPeter Brune if (oscatter) { 329111ade9eSPeter Brune for (i=0; i<n; i++) {ierr = PetscObjectReference((PetscObject)oscatter[i]);CHKERRQ(ierr);} 330eaedb033SPeter Brune } 331111ade9eSPeter Brune if (iscatter) { 332111ade9eSPeter Brune for (i=0; i<n; i++) {ierr = PetscObjectReference((PetscObject)iscatter[i]);CHKERRQ(ierr);} 333eaedb033SPeter Brune } 334111ade9eSPeter Brune if (gscatter) { 335111ade9eSPeter Brune for (i=0; i<n; i++) {ierr = PetscObjectReference((PetscObject)gscatter[i]);CHKERRQ(ierr);} 336111ade9eSPeter Brune } 337111ade9eSPeter Brune if (oscatter) { 338111ade9eSPeter Brune ierr = PetscMalloc(n*sizeof(IS),&nasm->oscatter);CHKERRQ(ierr); 339eaedb033SPeter Brune for (i=0; i<n; i++) { 340111ade9eSPeter Brune nasm->oscatter[i] = oscatter[i]; 341eaedb033SPeter Brune } 342111ade9eSPeter Brune } 343111ade9eSPeter Brune if (iscatter) { 344111ade9eSPeter Brune ierr = PetscMalloc(n*sizeof(IS),&nasm->iscatter);CHKERRQ(ierr); 345eaedb033SPeter Brune for (i=0; i<n; i++) { 346111ade9eSPeter Brune nasm->iscatter[i] = iscatter[i]; 347eaedb033SPeter Brune } 348eaedb033SPeter Brune } 349111ade9eSPeter Brune if (gscatter) { 350111ade9eSPeter Brune ierr = PetscMalloc(n*sizeof(IS),&nasm->gscatter);CHKERRQ(ierr); 351eaedb033SPeter Brune for (i=0; i<n; i++) { 352111ade9eSPeter Brune nasm->gscatter[i] = gscatter[i]; 353eaedb033SPeter Brune } 354eaedb033SPeter Brune } 355111ade9eSPeter Brune 356eaedb033SPeter Brune if (subsnes) { 357eaedb033SPeter Brune ierr = PetscMalloc(n*sizeof(SNES),&nasm->subsnes);CHKERRQ(ierr); 358eaedb033SPeter Brune for (i=0; i<n; i++) { 359eaedb033SPeter Brune nasm->subsnes[i] = subsnes[i]; 360eaedb033SPeter Brune } 361a4f17876SPeter Brune nasm->same_local_solves = PETSC_FALSE; 362eaedb033SPeter Brune } 363eaedb033SPeter Brune PetscFunctionReturn(0); 364eaedb033SPeter Brune } 365eaedb033SPeter Brune 366eaedb033SPeter Brune #undef __FUNCT__ 36776857b2aSPeter Brune #define __FUNCT__ "SNESNASMGetSubdomains" 36876857b2aSPeter Brune /*@ 36976857b2aSPeter Brune SNESNASMGetSubdomains - Get the local subdomain context. 37076857b2aSPeter Brune 37176857b2aSPeter Brune Not Collective 37276857b2aSPeter Brune 37376857b2aSPeter Brune Input Parameters: 37476857b2aSPeter Brune . SNES - the SNES context 37576857b2aSPeter Brune 37676857b2aSPeter Brune Output Parameters: 37776857b2aSPeter Brune + n - the number of local subdomains 37876857b2aSPeter Brune . subsnes - solvers defined on the local subdomains 37976857b2aSPeter Brune . iscatter - scatters into the nonoverlapping portions of the local subdomains 38076857b2aSPeter Brune . oscatter - scatters into the overlapping portions of the local subdomains 38176857b2aSPeter Brune - gscatter - scatters into the (ghosted) local vector of the local subdomain 38276857b2aSPeter Brune 38376857b2aSPeter Brune Level: intermediate 38476857b2aSPeter Brune 38576857b2aSPeter Brune .keywords: SNES, NASM 38676857b2aSPeter Brune 38776857b2aSPeter Brune .seealso: SNESNASM, SNESNASMSetSubdomains() 38876857b2aSPeter Brune @*/ 38976857b2aSPeter Brune PetscErrorCode SNESNASMGetSubdomains(SNES snes,PetscInt *n,SNES *subsnes[],VecScatter *iscatter[],VecScatter *oscatter[],VecScatter *gscatter[]) 39076857b2aSPeter Brune { 39176857b2aSPeter Brune PetscErrorCode ierr; 39276857b2aSPeter Brune PetscErrorCode (*f)(SNES,PetscInt*,SNES**,VecScatter**,VecScatter**,VecScatter**); 39376857b2aSPeter Brune 39476857b2aSPeter Brune PetscFunctionBegin; 3950005d66cSJed Brown ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESNASMGetSubdomains_C",&f);CHKERRQ(ierr); 3964b838e8fSPeter Brune if (f) {ierr = (f)(snes,n,subsnes,iscatter,oscatter,gscatter);CHKERRQ(ierr);} 39776857b2aSPeter Brune PetscFunctionReturn(0); 39876857b2aSPeter Brune } 39976857b2aSPeter Brune 40076857b2aSPeter Brune #undef __FUNCT__ 40176857b2aSPeter Brune #define __FUNCT__ "SNESNASMGetSubdomains_NASM" 40276857b2aSPeter Brune PetscErrorCode SNESNASMGetSubdomains_NASM(SNES snes,PetscInt *n,SNES *subsnes[],VecScatter *iscatter[],VecScatter *oscatter[],VecScatter *gscatter[]) 40376857b2aSPeter Brune { 40476857b2aSPeter Brune SNES_NASM *nasm = (SNES_NASM*)snes->data; 40576857b2aSPeter Brune 40676857b2aSPeter Brune PetscFunctionBegin; 40776857b2aSPeter Brune if (n) *n = nasm->n; 40876857b2aSPeter Brune if (oscatter) *oscatter = nasm->oscatter; 40976857b2aSPeter Brune if (iscatter) *iscatter = nasm->iscatter; 41076857b2aSPeter Brune if (gscatter) *gscatter = nasm->gscatter; 411a4f17876SPeter Brune if (subsnes) { 412a4f17876SPeter Brune *subsnes = nasm->subsnes; 413a4f17876SPeter Brune nasm->same_local_solves = PETSC_FALSE; 414a4f17876SPeter Brune } 41576857b2aSPeter Brune PetscFunctionReturn(0); 41676857b2aSPeter Brune } 41776857b2aSPeter Brune 41876857b2aSPeter Brune #undef __FUNCT__ 41976857b2aSPeter Brune #define __FUNCT__ "SNESNASMGetSubdomainVecs" 42076857b2aSPeter Brune /*@ 42176857b2aSPeter Brune SNESNASMGetSubdomainVecs - Get the processor-local subdomain vectors 42276857b2aSPeter Brune 42376857b2aSPeter Brune Not Collective 42476857b2aSPeter Brune 42576857b2aSPeter Brune Input Parameters: 42676857b2aSPeter Brune . SNES - the SNES context 42776857b2aSPeter Brune 42876857b2aSPeter Brune Output Parameters: 42976857b2aSPeter Brune + n - the number of local subdomains 43076857b2aSPeter Brune . x - The subdomain solution vector 43176857b2aSPeter Brune . y - The subdomain step vector 43276857b2aSPeter Brune . b - The subdomain RHS vector 43376857b2aSPeter Brune - xl - The subdomain local vectors (ghosted) 43476857b2aSPeter Brune 43576857b2aSPeter Brune Level: developer 43676857b2aSPeter Brune 43776857b2aSPeter Brune .keywords: SNES, NASM 43876857b2aSPeter Brune 43976857b2aSPeter Brune .seealso: SNESNASM, SNESNASMGetSubdomains() 44076857b2aSPeter Brune @*/ 44176857b2aSPeter Brune PetscErrorCode SNESNASMGetSubdomainVecs(SNES snes,PetscInt *n,Vec **x,Vec **y,Vec **b, Vec **xl) 44276857b2aSPeter Brune { 44376857b2aSPeter Brune PetscErrorCode ierr; 44476857b2aSPeter Brune PetscErrorCode (*f)(SNES,PetscInt*,Vec**,Vec**,Vec**,Vec**); 44576857b2aSPeter Brune 44676857b2aSPeter Brune PetscFunctionBegin; 4470005d66cSJed Brown ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESNASMGetSubdomainVecs_C",&f);CHKERRQ(ierr); 4484b838e8fSPeter Brune if (f) {ierr = (f)(snes,n,x,y,b,xl);CHKERRQ(ierr);} 44976857b2aSPeter Brune PetscFunctionReturn(0); 45076857b2aSPeter Brune } 45176857b2aSPeter Brune 45276857b2aSPeter Brune #undef __FUNCT__ 45376857b2aSPeter Brune #define __FUNCT__ "SNESNASMGetSubdomainVecs_NASM" 45476857b2aSPeter Brune PetscErrorCode SNESNASMGetSubdomainVecs_NASM(SNES snes,PetscInt *n,Vec **x,Vec **y,Vec **b,Vec **xl) 45576857b2aSPeter Brune { 45676857b2aSPeter Brune SNES_NASM *nasm = (SNES_NASM*)snes->data; 45776857b2aSPeter Brune 45876857b2aSPeter Brune PetscFunctionBegin; 45976857b2aSPeter Brune if (n) *n = nasm->n; 46076857b2aSPeter Brune if (x) *x = nasm->x; 46176857b2aSPeter Brune if (y) *y = nasm->y; 46276857b2aSPeter Brune if (b) *b = nasm->b; 46376857b2aSPeter Brune if (xl) *xl = nasm->xl; 46476857b2aSPeter Brune PetscFunctionReturn(0); 46576857b2aSPeter Brune } 46676857b2aSPeter Brune 467d728fb7dSPeter Brune #undef __FUNCT__ 468d728fb7dSPeter Brune #define __FUNCT__ "SNESNASMSetComputeFinalJacobian" 469d728fb7dSPeter Brune /*@ 470d728fb7dSPeter Brune SNESNASMSetComputeFinalJacobian - Schedules the computation of the global and subdomain jacobians upon convergence 471d728fb7dSPeter Brune 472d728fb7dSPeter Brune Collective on SNES 473d728fb7dSPeter Brune 474d728fb7dSPeter Brune Input Parameters: 475d728fb7dSPeter Brune + SNES - the SNES context 476d728fb7dSPeter Brune - flg - indication of whether to compute the jacobians or not 477d728fb7dSPeter Brune 478d728fb7dSPeter Brune Level: developer 479d728fb7dSPeter Brune 480d728fb7dSPeter Brune Notes: This is used almost exclusively in the implementation of ASPIN, where the converged subdomain and global jacobian 481d728fb7dSPeter Brune is needed at each linear iteration. 482d728fb7dSPeter Brune 483d728fb7dSPeter Brune .keywords: SNES, NASM, ASPIN 484d728fb7dSPeter Brune 485d728fb7dSPeter Brune .seealso: SNESNASM, SNESNASMGetSubdomains() 486d728fb7dSPeter Brune @*/ 487d728fb7dSPeter Brune PetscErrorCode SNESNASMSetComputeFinalJacobian(SNES snes,PetscBool flg) 488d728fb7dSPeter Brune { 489d728fb7dSPeter Brune PetscErrorCode (*f)(SNES,PetscBool); 490d728fb7dSPeter Brune PetscErrorCode ierr; 491d728fb7dSPeter Brune 492d728fb7dSPeter Brune PetscFunctionBegin; 4930005d66cSJed Brown ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESNASMSetComputeFinalJacobian_C",&f);CHKERRQ(ierr); 4944b838e8fSPeter Brune if (f) {ierr = (f)(snes,flg);CHKERRQ(ierr);} 495d728fb7dSPeter Brune PetscFunctionReturn(0); 496d728fb7dSPeter Brune } 497d728fb7dSPeter Brune 498d728fb7dSPeter Brune #undef __FUNCT__ 499d728fb7dSPeter Brune #define __FUNCT__ "SNESNASMSetComputeFinalJacobian_NASM" 500d728fb7dSPeter Brune PetscErrorCode SNESNASMSetComputeFinalJacobian_NASM(SNES snes,PetscBool flg) 501d728fb7dSPeter Brune { 502d728fb7dSPeter Brune SNES_NASM *nasm = (SNES_NASM*)snes->data; 503d728fb7dSPeter Brune 504d728fb7dSPeter Brune PetscFunctionBegin; 505d728fb7dSPeter Brune nasm->finaljacobian = flg; 506602bec5dSPeter Brune if (flg) snes->usesksp = PETSC_TRUE; 507d728fb7dSPeter Brune PetscFunctionReturn(0); 508d728fb7dSPeter Brune } 50976857b2aSPeter Brune 51076857b2aSPeter Brune #undef __FUNCT__ 511610116beSPeter Brune #define __FUNCT__ "SNESNASMSetDamping" 512610116beSPeter Brune /*@ 513610116beSPeter Brune SNESNASMSetDamping - Sets the update damping for NASM 514610116beSPeter Brune 515610116beSPeter Brune Logically collective on SNES 516610116beSPeter Brune 517610116beSPeter Brune Input Parameters: 518610116beSPeter Brune + SNES - the SNES context 519610116beSPeter Brune - dmp - damping 520610116beSPeter Brune 521610116beSPeter Brune Level: intermediate 522610116beSPeter Brune 523610116beSPeter Brune .keywords: SNES, NASM, damping 524610116beSPeter Brune 525610116beSPeter Brune .seealso: SNESNASM, SNESNASMGetDamping() 526610116beSPeter Brune @*/ 527610116beSPeter Brune PetscErrorCode SNESNASMSetDamping(SNES snes,PetscReal dmp) 528610116beSPeter Brune { 529610116beSPeter Brune PetscErrorCode (*f)(SNES,PetscReal); 530610116beSPeter Brune PetscErrorCode ierr; 531610116beSPeter Brune 532610116beSPeter Brune PetscFunctionBegin; 533610116beSPeter Brune ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESNASMSetDamping_C",(void (**)(void))&f);CHKERRQ(ierr); 534e24b0d94SPeter Brune if (f) {ierr = (f)(snes,dmp);CHKERRQ(ierr);} 535610116beSPeter Brune PetscFunctionReturn(0); 536610116beSPeter Brune } 537610116beSPeter Brune 538610116beSPeter Brune #undef __FUNCT__ 539610116beSPeter Brune #define __FUNCT__ "SNESNASMSetDamping_NASM" 540610116beSPeter Brune PetscErrorCode SNESNASMSetDamping_NASM(SNES snes,PetscReal dmp) 541610116beSPeter Brune { 542610116beSPeter Brune SNES_NASM *nasm = (SNES_NASM*)snes->data; 543610116beSPeter Brune 544610116beSPeter Brune PetscFunctionBegin; 545610116beSPeter Brune nasm->damping = dmp; 546610116beSPeter Brune PetscFunctionReturn(0); 547610116beSPeter Brune } 548610116beSPeter Brune 549610116beSPeter Brune #undef __FUNCT__ 550610116beSPeter Brune #define __FUNCT__ "SNESNASMGetDamping" 551610116beSPeter Brune /*@ 552610116beSPeter Brune SNESNASMGetDamping - Gets the update damping for NASM 553610116beSPeter Brune 554610116beSPeter Brune Not Collective 555610116beSPeter Brune 556610116beSPeter Brune Input Parameters: 557610116beSPeter Brune + SNES - the SNES context 558610116beSPeter Brune - dmp - damping 559610116beSPeter Brune 560610116beSPeter Brune Level: intermediate 561610116beSPeter Brune 562610116beSPeter Brune .keywords: SNES, NASM, damping 563610116beSPeter Brune 564610116beSPeter Brune .seealso: SNESNASM, SNESNASMSetDamping() 565610116beSPeter Brune @*/ 566610116beSPeter Brune PetscErrorCode SNESNASMGetDamping(SNES snes,PetscReal *dmp) 567610116beSPeter Brune { 568610116beSPeter Brune PetscErrorCode (*f)(SNES,PetscReal*); 569610116beSPeter Brune PetscErrorCode ierr; 570610116beSPeter Brune 571610116beSPeter Brune PetscFunctionBegin; 572610116beSPeter Brune ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESNASMGetDamping_C",(void (**)(void))&f);CHKERRQ(ierr); 573e24b0d94SPeter Brune if (f) {ierr = (f)(snes,dmp);CHKERRQ(ierr);} 574610116beSPeter Brune PetscFunctionReturn(0); 575610116beSPeter Brune } 576610116beSPeter Brune 577610116beSPeter Brune #undef __FUNCT__ 578610116beSPeter Brune #define __FUNCT__ "SNESNASMGetDamping_NASM" 579610116beSPeter Brune PetscErrorCode SNESNASMGetDamping_NASM(SNES snes,PetscReal *dmp) 580610116beSPeter Brune { 581610116beSPeter Brune SNES_NASM *nasm = (SNES_NASM*)snes->data; 582610116beSPeter Brune 583610116beSPeter Brune PetscFunctionBegin; 584610116beSPeter Brune *dmp = nasm->damping; 585610116beSPeter Brune PetscFunctionReturn(0); 586610116beSPeter Brune } 587610116beSPeter Brune 588610116beSPeter Brune 589610116beSPeter Brune #undef __FUNCT__ 590eaedb033SPeter Brune #define __FUNCT__ "SNESNASMSolveLocal_Private" 5910adebc6cSBarry Smith PetscErrorCode SNESNASMSolveLocal_Private(SNES snes,Vec B,Vec Y,Vec X) 5920adebc6cSBarry Smith { 593eaedb033SPeter Brune SNES_NASM *nasm = (SNES_NASM*)snes->data; 594258e1594SPeter Brune SNES subsnes; 595eaedb033SPeter Brune PetscInt i; 596610116beSPeter Brune PetscReal dmp; 597eaedb033SPeter Brune PetscErrorCode ierr; 598111ade9eSPeter Brune Vec Xlloc,Xl,Bl,Yl; 599111ade9eSPeter Brune VecScatter iscat,oscat,gscat; 600111ade9eSPeter Brune DM dm,subdm; 6010adebc6cSBarry Smith 602eaedb033SPeter Brune PetscFunctionBegin; 603eaedb033SPeter Brune ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 604610116beSPeter Brune ierr = SNESNASMGetDamping(snes,&dmp);CHKERRQ(ierr); 605111ade9eSPeter Brune ierr = VecSet(Y,0);CHKERRQ(ierr); 606b20c023fSPeter Brune if (nasm->eventrestrictinterp) {ierr = PetscLogEventBegin(nasm->eventrestrictinterp,snes,0,0,0);CHKERRQ(ierr);} 607eaedb033SPeter Brune for (i=0; i<nasm->n; i++) { 60870c78f05SPeter Brune /* scatter the solution to the local solution */ 60970c78f05SPeter Brune Xlloc = nasm->xl[i]; 61070c78f05SPeter Brune gscat = nasm->gscatter[i]; 61170c78f05SPeter Brune oscat = nasm->oscatter[i]; 61270c78f05SPeter Brune ierr = VecScatterBegin(gscat,X,Xlloc,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 61370c78f05SPeter Brune if (B) { 61470c78f05SPeter Brune /* scatter the RHS to the local RHS */ 61570c78f05SPeter Brune Bl = nasm->b[i]; 61670c78f05SPeter Brune ierr = VecScatterBegin(oscat,B,Bl,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 61770c78f05SPeter Brune } 61870c78f05SPeter Brune } 619b20c023fSPeter Brune if (nasm->eventrestrictinterp) {ierr = PetscLogEventEnd(nasm->eventrestrictinterp,snes,0,0,0);CHKERRQ(ierr);} 620b20c023fSPeter Brune 621b20c023fSPeter Brune 622b20c023fSPeter Brune if (nasm->eventsubsolve) {ierr = PetscLogEventBegin(nasm->eventsubsolve,snes,0,0,0);CHKERRQ(ierr);} 62370c78f05SPeter Brune for (i=0; i<nasm->n; i++) { 62470c78f05SPeter Brune Xl = nasm->x[i]; 62570c78f05SPeter Brune Xlloc = nasm->xl[i]; 62670c78f05SPeter Brune Yl = nasm->y[i]; 627258e1594SPeter Brune subsnes = nasm->subsnes[i]; 628258e1594SPeter Brune ierr = SNESGetDM(subsnes,&subdm);CHKERRQ(ierr); 629111ade9eSPeter Brune iscat = nasm->iscatter[i]; 630111ade9eSPeter Brune oscat = nasm->oscatter[i]; 631111ade9eSPeter Brune gscat = nasm->gscatter[i]; 632ed3c11a9SPeter Brune ierr = VecScatterEnd(gscat,X,Xlloc,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 63324b7f281SPeter Brune if (B) { 63424b7f281SPeter Brune Bl = nasm->b[i]; 635ed3c11a9SPeter Brune ierr = VecScatterEnd(oscat,B,Bl,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 636ed3c11a9SPeter Brune } else Bl = NULL; 637ed3c11a9SPeter Brune ierr = DMSubDomainRestrict(dm,oscat,gscat,subdm);CHKERRQ(ierr); 63870c78f05SPeter Brune ierr = DMLocalToGlobalBegin(subdm,Xlloc,INSERT_VALUES,Xl);CHKERRQ(ierr); 63970c78f05SPeter Brune ierr = DMLocalToGlobalEnd(subdm,Xlloc,INSERT_VALUES,Xl);CHKERRQ(ierr); 640111ade9eSPeter Brune ierr = VecCopy(Xl,Yl);CHKERRQ(ierr); 641ed3c11a9SPeter Brune ierr = SNESSolve(subsnes,Bl,Xl);CHKERRQ(ierr); 642ed3c11a9SPeter Brune ierr = VecAYPX(Yl,-1.0,Xl);CHKERRQ(ierr); 643111ade9eSPeter Brune if (nasm->type == PC_ASM_BASIC) { 644111ade9eSPeter Brune ierr = VecScatterBegin(oscat,Yl,Y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 645111ade9eSPeter Brune } else if (nasm->type == PC_ASM_RESTRICT) { 646111ade9eSPeter Brune ierr = VecScatterBegin(iscat,Yl,Y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 647ce94432eSBarry Smith } else SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"Only basic and restrict types are supported for SNESNASM"); 648eaedb033SPeter Brune } 649ed3c11a9SPeter Brune if (nasm->eventsubsolve) {ierr = PetscLogEventEnd(nasm->eventsubsolve,snes,0,0,0);CHKERRQ(ierr);} 650ed3c11a9SPeter Brune if (nasm->eventrestrictinterp) {ierr = PetscLogEventBegin(nasm->eventrestrictinterp,snes,0,0,0);CHKERRQ(ierr);} 65170c78f05SPeter Brune for (i=0; i<nasm->n; i++) { 65270c78f05SPeter Brune Yl = nasm->y[i]; 65370c78f05SPeter Brune iscat = nasm->iscatter[i]; 65470c78f05SPeter Brune oscat = nasm->oscatter[i]; 65570c78f05SPeter Brune if (nasm->type == PC_ASM_BASIC) { 65670c78f05SPeter Brune ierr = VecScatterEnd(oscat,Yl,Y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 65770c78f05SPeter Brune } else if (nasm->type == PC_ASM_RESTRICT) { 65870c78f05SPeter Brune ierr = VecScatterEnd(iscat,Yl,Y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 659ce94432eSBarry Smith } else SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"Only basic and restrict types are supported for SNESNASM"); 66070c78f05SPeter Brune } 661b20c023fSPeter Brune if (nasm->eventrestrictinterp) {ierr = PetscLogEventEnd(nasm->eventrestrictinterp,snes,0,0,0);CHKERRQ(ierr);} 662610116beSPeter Brune ierr = VecAXPY(X,dmp,Y);CHKERRQ(ierr); 663eaedb033SPeter Brune PetscFunctionReturn(0); 664eaedb033SPeter Brune } 665eaedb033SPeter Brune 666eaedb033SPeter Brune #undef __FUNCT__ 667d728fb7dSPeter Brune #define __FUNCT__ "SNESNASMComputeFinalJacobian_Private" 668602bec5dSPeter Brune PetscErrorCode SNESNASMComputeFinalJacobian_Private(SNES snes, Vec Xfinal) 669d728fb7dSPeter Brune { 670602bec5dSPeter Brune Vec X = Xfinal; 671d728fb7dSPeter Brune SNES_NASM *nasm = (SNES_NASM*)snes->data; 672d728fb7dSPeter Brune SNES subsnes; 673ca44f815SPeter Brune PetscInt i,lag = 1; 674d728fb7dSPeter Brune PetscErrorCode ierr; 675e59f0a30SPeter Brune Vec Xlloc,Xl,Fl,F; 676d728fb7dSPeter Brune VecScatter oscat,gscat; 677d728fb7dSPeter Brune DM dm,subdm; 678d728fb7dSPeter Brune MatStructure flg = DIFFERENT_NONZERO_PATTERN; 679d728fb7dSPeter Brune PetscFunctionBegin; 680602bec5dSPeter Brune if (nasm->fjtype == 2) X = nasm->xinit; 681e59f0a30SPeter Brune F = snes->vec_func; 682365a6726SPeter Brune if (snes->normschedule == SNES_NORM_NONE) {ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);} 683d728fb7dSPeter Brune ierr = SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg);CHKERRQ(ierr); 684d728fb7dSPeter Brune ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 685d728fb7dSPeter Brune if (nasm->eventrestrictinterp) {ierr = PetscLogEventBegin(nasm->eventrestrictinterp,snes,0,0,0);CHKERRQ(ierr);} 686602bec5dSPeter Brune if (nasm->fjtype != 1) { 687d728fb7dSPeter Brune for (i=0; i<nasm->n; i++) { 688d728fb7dSPeter Brune Xlloc = nasm->xl[i]; 689d728fb7dSPeter Brune gscat = nasm->gscatter[i]; 690d728fb7dSPeter Brune oscat = nasm->oscatter[i]; 691d728fb7dSPeter Brune ierr = VecScatterBegin(gscat,X,Xlloc,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 692602bec5dSPeter Brune } 693d728fb7dSPeter Brune } 694d728fb7dSPeter Brune if (nasm->eventrestrictinterp) {ierr = PetscLogEventEnd(nasm->eventrestrictinterp,snes,0,0,0);CHKERRQ(ierr);} 695d728fb7dSPeter Brune for (i=0; i<nasm->n; i++) { 696e59f0a30SPeter Brune Fl = nasm->subsnes[i]->vec_func; 697d728fb7dSPeter Brune Xl = nasm->x[i]; 698d728fb7dSPeter Brune Xlloc = nasm->xl[i]; 699d728fb7dSPeter Brune subsnes = nasm->subsnes[i]; 700d728fb7dSPeter Brune oscat = nasm->oscatter[i]; 701d728fb7dSPeter Brune gscat = nasm->gscatter[i]; 702602bec5dSPeter Brune if (nasm->fjtype != 1) {ierr = VecScatterEnd(gscat,X,Xlloc,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);} 703ed3c11a9SPeter Brune ierr = SNESGetDM(subsnes,&subdm);CHKERRQ(ierr); 704d728fb7dSPeter Brune ierr = DMSubDomainRestrict(dm,oscat,gscat,subdm);CHKERRQ(ierr); 705602bec5dSPeter Brune if (nasm->fjtype != 1) { 706d728fb7dSPeter Brune ierr = DMLocalToGlobalBegin(subdm,Xlloc,INSERT_VALUES,Xl);CHKERRQ(ierr); 707d728fb7dSPeter Brune ierr = DMLocalToGlobalEnd(subdm,Xlloc,INSERT_VALUES,Xl);CHKERRQ(ierr); 708602bec5dSPeter Brune } 709ca44f815SPeter Brune if (subsnes->lagjacobian == -1) subsnes->lagjacobian = -2; 710ca44f815SPeter Brune else if (subsnes->lagjacobian > 1) lag = subsnes->lagjacobian; 711602bec5dSPeter Brune ierr = SNESComputeFunction(subsnes,Xl,Fl);CHKERRQ(ierr); 7122f543b25SPeter Brune ierr = SNESComputeJacobian(subsnes,Xl,&subsnes->jacobian,&subsnes->jacobian_pre,&flg);CHKERRQ(ierr); 713ca44f815SPeter Brune if (lag > 1) subsnes->lagjacobian = lag; 714d728fb7dSPeter Brune ierr = KSPSetOperators(subsnes->ksp,subsnes->jacobian,subsnes->jacobian_pre,flg);CHKERRQ(ierr); 715d728fb7dSPeter Brune } 716d728fb7dSPeter Brune PetscFunctionReturn(0); 717d728fb7dSPeter Brune } 718d728fb7dSPeter Brune 719d728fb7dSPeter Brune #undef __FUNCT__ 720eaedb033SPeter Brune #define __FUNCT__ "SNESSolve_NASM" 721eaedb033SPeter Brune PetscErrorCode SNESSolve_NASM(SNES snes) 722eaedb033SPeter Brune { 723eaedb033SPeter Brune Vec F; 724eaedb033SPeter Brune Vec X; 725eaedb033SPeter Brune Vec B; 726111ade9eSPeter Brune Vec Y; 727eaedb033SPeter Brune PetscInt i; 728ed3c11a9SPeter Brune PetscReal fnorm = 0.0; 729eaedb033SPeter Brune PetscErrorCode ierr; 730365a6726SPeter Brune SNESNormSchedule normschedule; 731d728fb7dSPeter Brune SNES_NASM *nasm = (SNES_NASM*)snes->data; 732eaedb033SPeter Brune 733eaedb033SPeter Brune PetscFunctionBegin; 734eaedb033SPeter Brune X = snes->vec_sol; 735111ade9eSPeter Brune Y = snes->vec_sol_update; 736eaedb033SPeter Brune F = snes->vec_func; 737eaedb033SPeter Brune B = snes->vec_rhs; 738eaedb033SPeter Brune 739*e04113cfSBarry Smith ierr = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr); 740eaedb033SPeter Brune snes->iter = 0; 741eaedb033SPeter Brune snes->norm = 0.; 742*e04113cfSBarry Smith ierr = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr); 743eaedb033SPeter Brune snes->reason = SNES_CONVERGED_ITERATING; 744365a6726SPeter Brune ierr = SNESGetNormSchedule(snes, &normschedule);CHKERRQ(ierr); 745365a6726SPeter Brune if (normschedule == SNES_NORM_ALWAYS || normschedule == SNES_NORM_INITIAL_ONLY || normschedule == SNES_NORM_INITIAL_FINAL_ONLY) { 746eaedb033SPeter Brune /* compute the initial function and preconditioned update delX */ 747eaedb033SPeter Brune if (!snes->vec_func_init_set) { 748eaedb033SPeter Brune ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 749eaedb033SPeter Brune if (snes->domainerror) { 750eaedb033SPeter Brune snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 751eaedb033SPeter Brune PetscFunctionReturn(0); 752eaedb033SPeter Brune } 7531aa26658SKarl Rupp } else snes->vec_func_init_set = PETSC_FALSE; 754eaedb033SPeter Brune 755eaedb033SPeter Brune /* convergence test */ 756eaedb033SPeter Brune if (!snes->norm_init_set) { 757eaedb033SPeter Brune ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F|| */ 758189a9710SBarry Smith if (PetscIsInfOrNanReal(fnorm)) { 759189a9710SBarry Smith snes->reason = SNES_DIVERGED_FNORM_NAN; 760189a9710SBarry Smith PetscFunctionReturn(0); 761189a9710SBarry Smith } 762eaedb033SPeter Brune } else { 763eaedb033SPeter Brune fnorm = snes->norm_init; 764eaedb033SPeter Brune snes->norm_init_set = PETSC_FALSE; 765eaedb033SPeter Brune } 766*e04113cfSBarry Smith ierr = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr); 767eaedb033SPeter Brune snes->iter = 0; 768eaedb033SPeter Brune snes->norm = fnorm; 769*e04113cfSBarry Smith ierr = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr); 770a71f0d7dSBarry Smith ierr = SNESLogConvergenceHistory(snes,snes->norm,0);CHKERRQ(ierr); 771eaedb033SPeter Brune ierr = SNESMonitor(snes,0,snes->norm);CHKERRQ(ierr); 772eaedb033SPeter Brune 773eaedb033SPeter Brune /* set parameter for default relative tolerance convergence test */ 774eaedb033SPeter Brune snes->ttol = fnorm*snes->rtol; 775eaedb033SPeter Brune 776eaedb033SPeter Brune /* test convergence */ 777eaedb033SPeter Brune ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 778eaedb033SPeter Brune if (snes->reason) PetscFunctionReturn(0); 779eaedb033SPeter Brune } else { 780*e04113cfSBarry Smith ierr = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr); 781a71f0d7dSBarry Smith ierr = SNESLogConvergenceHistory(snes,snes->norm,0);CHKERRQ(ierr); 782eaedb033SPeter Brune ierr = SNESMonitor(snes,0,snes->norm);CHKERRQ(ierr); 783eaedb033SPeter Brune } 784eaedb033SPeter Brune 785eaedb033SPeter Brune /* Call general purpose update function */ 786eaedb033SPeter Brune if (snes->ops->update) { 787eaedb033SPeter Brune ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 788eaedb033SPeter Brune } 789602bec5dSPeter Brune /* copy the initial solution over for later */ 790602bec5dSPeter Brune if (nasm->fjtype == 2) {ierr = VecCopy(X,nasm->xinit);CHKERRQ(ierr);} 791eaedb033SPeter Brune 792eaedb033SPeter Brune for (i = 0; i < snes->max_its; i++) { 793111ade9eSPeter Brune ierr = SNESNASMSolveLocal_Private(snes,B,Y,X);CHKERRQ(ierr); 794365a6726SPeter Brune if (normschedule == SNES_NORM_ALWAYS || ((i == snes->max_its - 1) && (normschedule == SNES_NORM_INITIAL_FINAL_ONLY || normschedule == SNES_NORM_FINAL_ONLY))) { 795eaedb033SPeter Brune ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 796eaedb033SPeter Brune if (snes->domainerror) { 797eaedb033SPeter Brune snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 798d728fb7dSPeter Brune break; 799eaedb033SPeter Brune } 800eaedb033SPeter Brune ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F|| */ 801189a9710SBarry Smith if (PetscIsInfOrNanReal(fnorm)) { 802189a9710SBarry Smith snes->reason = SNES_DIVERGED_FNORM_NAN; 803189a9710SBarry Smith break; 804189a9710SBarry Smith } 805eaedb033SPeter Brune } 806eaedb033SPeter Brune /* Monitor convergence */ 807*e04113cfSBarry Smith ierr = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr); 808eaedb033SPeter Brune snes->iter = i+1; 809eaedb033SPeter Brune snes->norm = fnorm; 810*e04113cfSBarry Smith ierr = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr); 811a71f0d7dSBarry Smith ierr = SNESLogConvergenceHistory(snes,snes->norm,0);CHKERRQ(ierr); 812eaedb033SPeter Brune ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr); 813eaedb033SPeter Brune /* Test for convergence */ 814365a6726SPeter Brune if (normschedule == SNES_NORM_ALWAYS) {ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);} 815d728fb7dSPeter Brune if (snes->reason) break; 816eaedb033SPeter Brune /* Call general purpose update function */ 817e59f0a30SPeter Brune if (snes->ops->update) {ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);} 818eaedb033SPeter Brune } 819d728fb7dSPeter Brune if (nasm->finaljacobian) {ierr = SNESNASMComputeFinalJacobian_Private(snes,X);CHKERRQ(ierr);} 820365a6726SPeter Brune if (normschedule == SNES_NORM_ALWAYS) { 821eaedb033SPeter Brune if (i == snes->max_its) { 822eaedb033SPeter Brune ierr = PetscInfo1(snes,"Maximum number of iterations has been reached: %D\n",snes->max_its);CHKERRQ(ierr); 823eaedb033SPeter Brune if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT; 824eaedb033SPeter Brune } 8251aa26658SKarl Rupp } else if (!snes->reason) snes->reason = SNES_CONVERGED_ITS; /* NASM is meant to be used as a preconditioner */ 826eaedb033SPeter Brune PetscFunctionReturn(0); 827eaedb033SPeter Brune } 828eaedb033SPeter Brune 829eaedb033SPeter Brune /*MC 830eaedb033SPeter Brune SNESNASM - Nonlinear Additive Schwartz 831eaedb033SPeter Brune 83270c78603SPeter Brune Options Database: 83370c78603SPeter Brune + -snes_nasm_log - enable logging events for the communication and solve stages 83470c78603SPeter Brune . -snes_nasm_type <basic,restrict> - type of subdomain update used 83570c78603SPeter Brune . -snes_nasm_finaljacobian - compute the local and global jacobians of the final iterate 836602bec5dSPeter Brune . -snes_nasm_finaljacobian_type <finalinner,finalouter,initial> pick state the jacobian is calculated at 83770c78603SPeter Brune . -sub_snes_ - options prefix of the subdomain nonlinear solves 83870c78603SPeter Brune . -sub_ksp_ - options prefix of the subdomain Krylov solver 83970c78603SPeter Brune - -sub_pc_ - options prefix of the subdomain preconditioner 84070c78603SPeter Brune 841eaedb033SPeter Brune Level: advanced 842eaedb033SPeter Brune 843eaedb033SPeter Brune .seealso: SNESCreate(), SNES, SNESSetType(), SNESType (for list of available types) 844eaedb033SPeter Brune M*/ 845eaedb033SPeter Brune 846eaedb033SPeter Brune #undef __FUNCT__ 847eaedb033SPeter Brune #define __FUNCT__ "SNESCreate_NASM" 8488cc058d9SJed Brown PETSC_EXTERN PetscErrorCode SNESCreate_NASM(SNES snes) 849eaedb033SPeter Brune { 850eaedb033SPeter Brune SNES_NASM *nasm; 851eaedb033SPeter Brune PetscErrorCode ierr; 852eaedb033SPeter Brune 853eaedb033SPeter Brune PetscFunctionBegin; 854eaedb033SPeter Brune ierr = PetscNewLog(snes, SNES_NASM, &nasm);CHKERRQ(ierr); 855eaedb033SPeter Brune snes->data = (void*)nasm; 856eaedb033SPeter Brune 857eaedb033SPeter Brune nasm->n = PETSC_DECIDE; 858eaedb033SPeter Brune nasm->subsnes = 0; 859eaedb033SPeter Brune nasm->x = 0; 860111ade9eSPeter Brune nasm->xl = 0; 861111ade9eSPeter Brune nasm->y = 0; 862eaedb033SPeter Brune nasm->b = 0; 863111ade9eSPeter Brune nasm->oscatter = 0; 864111ade9eSPeter Brune nasm->iscatter = 0; 865111ade9eSPeter Brune nasm->gscatter = 0; 866610116beSPeter Brune nasm->damping = 1.; 867111ade9eSPeter Brune 868111ade9eSPeter Brune nasm->type = PC_ASM_BASIC; 869d728fb7dSPeter Brune nasm->finaljacobian = PETSC_FALSE; 870a4f17876SPeter Brune nasm->same_local_solves = PETSC_TRUE; 871eaedb033SPeter Brune 872eaedb033SPeter Brune snes->ops->destroy = SNESDestroy_NASM; 873eaedb033SPeter Brune snes->ops->setup = SNESSetUp_NASM; 874eaedb033SPeter Brune snes->ops->setfromoptions = SNESSetFromOptions_NASM; 875eaedb033SPeter Brune snes->ops->view = SNESView_NASM; 876eaedb033SPeter Brune snes->ops->solve = SNESSolve_NASM; 877eaedb033SPeter Brune snes->ops->reset = SNESReset_NASM; 878eaedb033SPeter Brune 879eaedb033SPeter Brune snes->usesksp = PETSC_FALSE; 880eaedb033SPeter Brune snes->usespc = PETSC_FALSE; 881eaedb033SPeter Brune 882602bec5dSPeter Brune nasm->fjtype = 0; 883602bec5dSPeter Brune nasm->xinit = NULL; 8840298fd71SBarry Smith nasm->eventrestrictinterp = 0; 8850298fd71SBarry Smith nasm->eventsubsolve = 0; 886b20c023fSPeter Brune 887eaedb033SPeter Brune if (!snes->tolerancesset) { 888eaedb033SPeter Brune snes->max_its = 10000; 889eaedb033SPeter Brune snes->max_funcs = 10000; 890eaedb033SPeter Brune } 891eaedb033SPeter Brune 892bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESNASMSetSubdomains_C",SNESNASMSetSubdomains_NASM);CHKERRQ(ierr); 893bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESNASMGetSubdomains_C",SNESNASMGetSubdomains_NASM);CHKERRQ(ierr); 89490bcee39SPeter Brune ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESNASMSetDamping_C",SNESNASMSetDamping_NASM);CHKERRQ(ierr); 89590bcee39SPeter Brune ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESNASMGetDamping_C",SNESNASMGetDamping_NASM);CHKERRQ(ierr); 896bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESNASMGetSubdomainVecs_C",SNESNASMGetSubdomainVecs_NASM);CHKERRQ(ierr); 897bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESNASMSetComputeFinalJacobian_C",SNESNASMSetComputeFinalJacobian_NASM);CHKERRQ(ierr); 898eaedb033SPeter Brune PetscFunctionReturn(0); 899eaedb033SPeter Brune } 90099e0435eSBarry Smith 901