1af0996ceSBarry Smith #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 */ 11f10b3e88SPatrick Farrell Vec weight; /* weighting for adding updates on overlaps, in global space */ 12111ade9eSPeter Brune VecScatter *oscatter; /* scatter from global space to the subdomain global space */ 13f10b3e88SPatrick Farrell VecScatter *oscatter_copy; /* copy of the above */ 14111ade9eSPeter Brune VecScatter *iscatter; /* scatter from global space to the nonoverlapping subdomain space */ 15111ade9eSPeter Brune VecScatter *gscatter; /* scatter from global space to the subdomain local space */ 16111ade9eSPeter Brune PCASMType type; /* ASM type */ 17111ade9eSPeter Brune PetscBool usesdm; /* use the DM for setting up the subproblems */ 18d728fb7dSPeter Brune PetscBool finaljacobian; /* compute the jacobian of the converged solution */ 19610116beSPeter Brune PetscReal damping; /* damping parameter for updates from the blocks */ 20f10b3e88SPatrick Farrell PetscBool weight_set; /* use a weight in the overlap updates */ 21b20c023fSPeter Brune 22b20c023fSPeter Brune /* logging events */ 23b20c023fSPeter Brune PetscLogEvent eventrestrictinterp; 24b20c023fSPeter Brune PetscLogEvent eventsubsolve; 25602bec5dSPeter Brune 26602bec5dSPeter Brune PetscInt fjtype; /* type of computed jacobian */ 27602bec5dSPeter Brune Vec xinit; /* initial solution in case the final jacobian type is computed as first */ 28eaedb033SPeter Brune } SNES_NASM; 29eaedb033SPeter Brune 309e5d0892SLisandro Dalcin const char *const SNESNASMTypes[] = {"NONE", "RESTRICT", "INTERPOLATE", "BASIC", "PCASMType", "PC_ASM_", NULL}; 31602bec5dSPeter Brune const char *const SNESNASMFJTypes[] = {"FINALOUTER", "FINALINNER", "INITIAL"}; 32b20c023fSPeter Brune 33d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESReset_NASM(SNES snes) 34d71ae5a4SJacob Faibussowitsch { 35eaedb033SPeter Brune SNES_NASM *nasm = (SNES_NASM *)snes->data; 36eaedb033SPeter Brune PetscInt i; 376e111a19SKarl Rupp 38eaedb033SPeter Brune PetscFunctionBegin; 39eaedb033SPeter Brune for (i = 0; i < nasm->n; i++) { 409566063dSJacob Faibussowitsch if (nasm->xl) PetscCall(VecDestroy(&nasm->xl[i])); 419566063dSJacob Faibussowitsch if (nasm->x) PetscCall(VecDestroy(&nasm->x[i])); 429566063dSJacob Faibussowitsch if (nasm->y) PetscCall(VecDestroy(&nasm->y[i])); 439566063dSJacob Faibussowitsch if (nasm->b) PetscCall(VecDestroy(&nasm->b[i])); 44eaedb033SPeter Brune 459566063dSJacob Faibussowitsch if (nasm->subsnes) PetscCall(SNESDestroy(&nasm->subsnes[i])); 469566063dSJacob Faibussowitsch if (nasm->oscatter) PetscCall(VecScatterDestroy(&nasm->oscatter[i])); 479566063dSJacob Faibussowitsch if (nasm->oscatter_copy) PetscCall(VecScatterDestroy(&nasm->oscatter_copy[i])); 489566063dSJacob Faibussowitsch if (nasm->iscatter) PetscCall(VecScatterDestroy(&nasm->iscatter[i])); 499566063dSJacob Faibussowitsch if (nasm->gscatter) PetscCall(VecScatterDestroy(&nasm->gscatter[i])); 50eaedb033SPeter Brune } 51111ade9eSPeter Brune 529566063dSJacob Faibussowitsch PetscCall(PetscFree(nasm->x)); 539566063dSJacob Faibussowitsch PetscCall(PetscFree(nasm->xl)); 549566063dSJacob Faibussowitsch PetscCall(PetscFree(nasm->y)); 559566063dSJacob Faibussowitsch PetscCall(PetscFree(nasm->b)); 56111ade9eSPeter Brune 579566063dSJacob Faibussowitsch if (nasm->xinit) PetscCall(VecDestroy(&nasm->xinit)); 58602bec5dSPeter Brune 599566063dSJacob Faibussowitsch PetscCall(PetscFree(nasm->subsnes)); 609566063dSJacob Faibussowitsch PetscCall(PetscFree(nasm->oscatter)); 619566063dSJacob Faibussowitsch PetscCall(PetscFree(nasm->oscatter_copy)); 629566063dSJacob Faibussowitsch PetscCall(PetscFree(nasm->iscatter)); 639566063dSJacob Faibussowitsch PetscCall(PetscFree(nasm->gscatter)); 64b20c023fSPeter Brune 6548a46eb9SPierre Jolivet if (nasm->weight_set) PetscCall(VecDestroy(&nasm->weight)); 66f10b3e88SPatrick Farrell 67b20c023fSPeter Brune nasm->eventrestrictinterp = 0; 68b20c023fSPeter Brune nasm->eventsubsolve = 0; 693ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 70eaedb033SPeter Brune } 71eaedb033SPeter Brune 72d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESDestroy_NASM(SNES snes) 73d71ae5a4SJacob Faibussowitsch { 74eaedb033SPeter Brune PetscFunctionBegin; 759566063dSJacob Faibussowitsch PetscCall(SNESReset_NASM(snes)); 762e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNASMSetType_C", NULL)); 772e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNASMGetType_C", NULL)); 782e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNASMSetSubdomains_C", NULL)); 792e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNASMGetSubdomains_C", NULL)); 802e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNASMSetDamping_C", NULL)); 812e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNASMGetDamping_C", NULL)); 822e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNASMGetSubdomainVecs_C", NULL)); 832e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNASMSetComputeFinalJacobian_C", NULL)); 849566063dSJacob Faibussowitsch PetscCall(PetscFree(snes->data)); 853ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 86eaedb033SPeter Brune } 87eaedb033SPeter Brune 88d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMGlobalToLocalSubDomainDirichletHook_Private(DM dm, Vec g, InsertMode mode, Vec l, void *ctx) 89d71ae5a4SJacob Faibussowitsch { 90111ade9eSPeter Brune Vec bcs = (Vec)ctx; 916e111a19SKarl Rupp 92111ade9eSPeter Brune PetscFunctionBegin; 939566063dSJacob Faibussowitsch PetscCall(VecCopy(bcs, l)); 943ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 95111ade9eSPeter Brune } 96111ade9eSPeter Brune 97d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESSetUp_NASM(SNES snes) 98d71ae5a4SJacob Faibussowitsch { 99eaedb033SPeter Brune SNES_NASM *nasm = (SNES_NASM *)snes->data; 10076857b2aSPeter Brune DM dm, subdm; 101111ade9eSPeter Brune DM *subdms; 102111ade9eSPeter Brune PetscInt i; 103eaedb033SPeter Brune const char *optionsprefix; 104111ade9eSPeter Brune Vec F; 105eaedb033SPeter Brune 106eaedb033SPeter Brune PetscFunctionBegin; 107eaedb033SPeter Brune if (!nasm->subsnes) { 1089566063dSJacob Faibussowitsch PetscCall(SNESGetDM(snes, &dm)); 1090a696f66SPeter Brune if (dm) { 110eaedb033SPeter Brune nasm->usesdm = PETSC_TRUE; 1119566063dSJacob Faibussowitsch PetscCall(DMCreateDomainDecomposition(dm, &nasm->n, NULL, NULL, NULL, &subdms)); 11228b400f6SJacob Faibussowitsch PetscCheck(subdms, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "DM has no default decomposition defined. Set subsolves manually with SNESNASMSetSubdomains()."); 1139566063dSJacob Faibussowitsch PetscCall(DMCreateDomainDecompositionScatters(dm, nasm->n, subdms, &nasm->iscatter, &nasm->oscatter, &nasm->gscatter)); 1149566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nasm->n, &nasm->oscatter_copy)); 11548a46eb9SPierre Jolivet for (i = 0; i < nasm->n; i++) PetscCall(VecScatterCopy(nasm->oscatter[i], &nasm->oscatter_copy[i])); 116eaedb033SPeter Brune 1179566063dSJacob Faibussowitsch PetscCall(SNESGetOptionsPrefix(snes, &optionsprefix)); 1189566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nasm->n, &nasm->subsnes)); 119111ade9eSPeter Brune for (i = 0; i < nasm->n; i++) { 120de54d9edSStefano Zampini PetscCall(SNESCreate(PetscObjectComm((PetscObject)subdms[i]), &nasm->subsnes[i])); 1219566063dSJacob Faibussowitsch PetscCall(PetscObjectIncrementTabLevel((PetscObject)nasm->subsnes[i], (PetscObject)snes, 1)); 1229566063dSJacob Faibussowitsch PetscCall(SNESAppendOptionsPrefix(nasm->subsnes[i], optionsprefix)); 1239566063dSJacob Faibussowitsch PetscCall(SNESAppendOptionsPrefix(nasm->subsnes[i], "sub_")); 1249566063dSJacob Faibussowitsch PetscCall(SNESSetDM(nasm->subsnes[i], subdms[i])); 1252891b2c6SStefano Zampini if (snes->ops->usercompute) { 1262891b2c6SStefano Zampini PetscCall(SNESSetComputeApplicationContext(nasm->subsnes[i], snes->ops->usercompute, snes->ops->userdestroy)); 1272891b2c6SStefano Zampini } else { 1282891b2c6SStefano Zampini void *ctx; 1292891b2c6SStefano Zampini 1302891b2c6SStefano Zampini PetscCall(SNESGetApplicationContext(snes, &ctx)); 1312891b2c6SStefano Zampini PetscCall(SNESSetApplicationContext(nasm->subsnes[i], ctx)); 1322891b2c6SStefano Zampini } 1339566063dSJacob Faibussowitsch PetscCall(SNESSetFromOptions(nasm->subsnes[i])); 1349566063dSJacob Faibussowitsch PetscCall(DMDestroy(&subdms[i])); 135111ade9eSPeter Brune } 1369566063dSJacob Faibussowitsch PetscCall(PetscFree(subdms)); 137ce94432eSBarry Smith } else SETERRQ(PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_WRONGSTATE, "Cannot construct local problems automatically without a DM!"); 138ce94432eSBarry Smith } else SETERRQ(PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_WRONGSTATE, "Must set subproblems manually if there is no DM!"); 139111ade9eSPeter Brune /* allocate the global vectors */ 14048a46eb9SPierre Jolivet if (!nasm->x) PetscCall(PetscCalloc1(nasm->n, &nasm->x)); 14148a46eb9SPierre Jolivet if (!nasm->xl) PetscCall(PetscCalloc1(nasm->n, &nasm->xl)); 14248a46eb9SPierre Jolivet if (!nasm->y) PetscCall(PetscCalloc1(nasm->n, &nasm->y)); 14348a46eb9SPierre Jolivet if (!nasm->b) PetscCall(PetscCalloc1(nasm->n, &nasm->b)); 144111ade9eSPeter Brune 145111ade9eSPeter Brune for (i = 0; i < nasm->n; i++) { 1469566063dSJacob Faibussowitsch PetscCall(SNESGetFunction(nasm->subsnes[i], &F, NULL, NULL)); 1479566063dSJacob Faibussowitsch if (!nasm->x[i]) PetscCall(VecDuplicate(F, &nasm->x[i])); 1489566063dSJacob Faibussowitsch if (!nasm->y[i]) PetscCall(VecDuplicate(F, &nasm->y[i])); 1499566063dSJacob Faibussowitsch if (!nasm->b[i]) PetscCall(VecDuplicate(F, &nasm->b[i])); 15076857b2aSPeter Brune if (!nasm->xl[i]) { 1519566063dSJacob Faibussowitsch PetscCall(SNESGetDM(nasm->subsnes[i], &subdm)); 1529566063dSJacob Faibussowitsch PetscCall(DMCreateLocalVector(subdm, &nasm->xl[i])); 1539566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalHookAdd(subdm, DMGlobalToLocalSubDomainDirichletHook_Private, NULL, nasm->xl[i])); 154111ade9eSPeter Brune } 15561ba4676SBarry Smith } 156602bec5dSPeter Brune if (nasm->finaljacobian) { 1579566063dSJacob Faibussowitsch PetscCall(SNESSetUpMatrices(snes)); 15848a46eb9SPierre Jolivet if (nasm->fjtype == 2) PetscCall(VecDuplicate(snes->vec_sol, &nasm->xinit)); 15948a46eb9SPierre Jolivet for (i = 0; i < nasm->n; i++) PetscCall(SNESSetUpMatrices(nasm->subsnes[i])); 160602bec5dSPeter Brune } 1613ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 162eaedb033SPeter Brune } 163eaedb033SPeter Brune 164d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESSetFromOptions_NASM(SNES snes, PetscOptionItems *PetscOptionsObject) 165d71ae5a4SJacob Faibussowitsch { 166111ade9eSPeter Brune PCASMType asmtype; 16783dc3634SPierre Jolivet PetscBool flg, monflg; 168111ade9eSPeter Brune SNES_NASM *nasm = (SNES_NASM *)snes->data; 1696e111a19SKarl Rupp 170eaedb033SPeter Brune PetscFunctionBegin; 171d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject, "Nonlinear Additive Schwarz options"); 1729566063dSJacob Faibussowitsch PetscCall(PetscOptionsEnum("-snes_nasm_type", "Type of restriction/extension", "", SNESNASMTypes, (PetscEnum)nasm->type, (PetscEnum *)&asmtype, &flg)); 1739566063dSJacob Faibussowitsch if (flg) PetscCall(SNESNASMSetType(snes, asmtype)); 174b20c023fSPeter Brune flg = PETSC_FALSE; 175b20c023fSPeter Brune monflg = PETSC_TRUE; 1769566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-snes_nasm_damping", "The new solution is obtained as old solution plus dmp times (sum of the solutions on the subdomains)", "SNESNASMSetDamping", nasm->damping, &nasm->damping, &flg)); 1779566063dSJacob Faibussowitsch if (flg) PetscCall(SNESNASMSetDamping(snes, nasm->damping)); 1789566063dSJacob Faibussowitsch PetscCall(PetscOptionsDeprecated("-snes_nasm_sub_view", NULL, "3.15", "Use -snes_view ::ascii_info_detail")); 1799566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-snes_nasm_finaljacobian", "Compute the global jacobian of the final iterate (for ASPIN)", "", nasm->finaljacobian, &nasm->finaljacobian, NULL)); 1809566063dSJacob Faibussowitsch PetscCall(PetscOptionsEList("-snes_nasm_finaljacobian_type", "The type of the final jacobian computed.", "", SNESNASMFJTypes, 3, SNESNASMFJTypes[0], &nasm->fjtype, NULL)); 1819566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-snes_nasm_log", "Log times for subSNES solves and restriction", "", monflg, &monflg, &flg)); 182b20c023fSPeter Brune if (flg) { 1839566063dSJacob Faibussowitsch PetscCall(PetscLogEventRegister("SNESNASMSubSolve", ((PetscObject)snes)->classid, &nasm->eventsubsolve)); 1849566063dSJacob Faibussowitsch PetscCall(PetscLogEventRegister("SNESNASMRestrict", ((PetscObject)snes)->classid, &nasm->eventrestrictinterp)); 185b20c023fSPeter Brune } 186d0609cedSBarry Smith PetscOptionsHeadEnd(); 1873ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 188eaedb033SPeter Brune } 189eaedb033SPeter Brune 190d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESView_NASM(SNES snes, PetscViewer viewer) 191d71ae5a4SJacob Faibussowitsch { 192b20c023fSPeter Brune SNES_NASM *nasm = (SNES_NASM *)snes->data; 193a4f17876SPeter Brune PetscMPIInt rank, size; 194dd2fa690SBarry Smith PetscInt i, N, bsz; 195b20c023fSPeter Brune PetscBool iascii, isstring; 196b20c023fSPeter Brune PetscViewer sviewer; 197ce94432eSBarry Smith MPI_Comm comm; 19883dc3634SPierre Jolivet PetscViewerFormat format; 19983dc3634SPierre Jolivet const char *prefix; 200b20c023fSPeter Brune 201b20c023fSPeter Brune PetscFunctionBegin; 2029566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)snes, &comm)); 2039566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii)); 2049566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERSTRING, &isstring)); 2059566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank)); 2069566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm, &size)); 2071c2dc1cbSBarry Smith PetscCall(MPIU_Allreduce(&nasm->n, &N, 1, MPIU_INT, MPI_SUM, comm)); 208b20c023fSPeter Brune if (iascii) { 20963a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " total subdomain blocks = %" PetscInt_FMT "\n", N)); 2109566063dSJacob Faibussowitsch PetscCall(PetscViewerGetFormat(viewer, &format)); 21183dc3634SPierre Jolivet if (format != PETSC_VIEWER_ASCII_INFO_DETAIL) { 212a4f17876SPeter Brune if (nasm->subsnes) { 2139566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " Local solver information for first block on rank 0:\n")); 2149566063dSJacob Faibussowitsch PetscCall(SNESGetOptionsPrefix(snes, &prefix)); 2159566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " Use -%ssnes_view ::ascii_info_detail to display information for all blocks\n", prefix ? prefix : "")); 2169566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushTab(viewer)); 2179566063dSJacob Faibussowitsch PetscCall(PetscViewerGetSubViewer(viewer, PETSC_COMM_SELF, &sviewer)); 218dd400576SPatrick Sanan if (rank == 0) { 219*b4025f61SBarry Smith PetscCall(PetscViewerASCIIPushTab(sviewer)); 2209566063dSJacob Faibussowitsch PetscCall(SNESView(nasm->subsnes[0], sviewer)); 221*b4025f61SBarry Smith PetscCall(PetscViewerASCIIPopTab(sviewer)); 222a4f17876SPeter Brune } 2239566063dSJacob Faibussowitsch PetscCall(PetscViewerRestoreSubViewer(viewer, PETSC_COMM_SELF, &sviewer)); 2249566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopTab(viewer)); 225a4f17876SPeter Brune } 226a4f17876SPeter Brune } else { 227a4f17876SPeter Brune /* print the solver on each block */ 2289566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushSynchronized(viewer)); 22963a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, " [%d] number of local blocks = %" PetscInt_FMT "\n", (int)rank, nasm->n)); 2309566063dSJacob Faibussowitsch PetscCall(PetscViewerFlush(viewer)); 2319566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopSynchronized(viewer)); 2329566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " Local solver information for each block is in the following SNES objects:\n")); 2339566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushTab(viewer)); 2349566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "- - - - - - - - - - - - - - - - - -\n")); 2359566063dSJacob Faibussowitsch PetscCall(PetscViewerGetSubViewer(viewer, PETSC_COMM_SELF, &sviewer)); 236a4f17876SPeter Brune for (i = 0; i < nasm->n; i++) { 2379566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(nasm->x[i], &bsz)); 23863a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(sviewer, "[%d] local block number %" PetscInt_FMT ", size = %" PetscInt_FMT "\n", (int)rank, i, bsz)); 2399566063dSJacob Faibussowitsch PetscCall(SNESView(nasm->subsnes[i], sviewer)); 2409566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(sviewer, "- - - - - - - - - - - - - - - - - -\n")); 241b20c023fSPeter Brune } 2429566063dSJacob Faibussowitsch PetscCall(PetscViewerRestoreSubViewer(viewer, PETSC_COMM_SELF, &sviewer)); 2439566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopTab(viewer)); 244a4f17876SPeter Brune } 245b20c023fSPeter Brune } else if (isstring) { 24663a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerStringSPrintf(viewer, " blocks=%" PetscInt_FMT ",type=%s", N, SNESNASMTypes[nasm->type])); 2479566063dSJacob Faibussowitsch PetscCall(PetscViewerGetSubViewer(viewer, PETSC_COMM_SELF, &sviewer)); 2489566063dSJacob Faibussowitsch if (nasm->subsnes && rank == 0) PetscCall(SNESView(nasm->subsnes[0], sviewer)); 2499566063dSJacob Faibussowitsch PetscCall(PetscViewerRestoreSubViewer(viewer, PETSC_COMM_SELF, &sviewer)); 250b20c023fSPeter Brune } 2513ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 252eaedb033SPeter Brune } 253eaedb033SPeter Brune 254e0331734SPeter Brune /*@ 255420bcc1bSBarry Smith SNESNASMSetType - Set the type of subdomain update used for the nonlinear additive Schwarz solver `SNESNASM` 256e0331734SPeter Brune 257c3339decSBarry Smith Logically Collective 258e0331734SPeter Brune 259e0331734SPeter Brune Input Parameters: 260f6dfbefdSBarry Smith + snes - the `SNES` context 261f6dfbefdSBarry Smith - type - the type of update, `PC_ASM_BASIC` or `PC_ASM_RESTRICT` 262e0331734SPeter Brune 263420bcc1bSBarry Smith Options Database Key: 264420bcc1bSBarry Smith . -snes_nasm_type <basic,restrict> - type of subdomain update used 265420bcc1bSBarry Smith 266e0331734SPeter Brune Level: intermediate 267e0331734SPeter Brune 268420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESNASM`, `SNESNASMGetType()`, `PCASMSetType()`, `PC_ASM_BASIC`, `PC_ASM_RESTRICT`, `PCASMType` 269e0331734SPeter Brune @*/ 270d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESNASMSetType(SNES snes, PCASMType type) 271d71ae5a4SJacob Faibussowitsch { 272e0331734SPeter Brune PetscErrorCode (*f)(SNES, PCASMType); 273e0331734SPeter Brune 274e0331734SPeter Brune PetscFunctionBegin; 2759566063dSJacob Faibussowitsch PetscCall(PetscObjectQueryFunction((PetscObject)snes, "SNESNASMSetType_C", &f)); 2769566063dSJacob Faibussowitsch if (f) PetscCall((f)(snes, type)); 2773ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 278e0331734SPeter Brune } 279e0331734SPeter Brune 280d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESNASMSetType_NASM(SNES snes, PCASMType type) 281d71ae5a4SJacob Faibussowitsch { 282e0331734SPeter Brune SNES_NASM *nasm = (SNES_NASM *)snes->data; 283e0331734SPeter Brune 284e0331734SPeter Brune PetscFunctionBegin; 2850b121fc5SBarry Smith PetscCheck(type == PC_ASM_BASIC || type == PC_ASM_RESTRICT, PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_OUTOFRANGE, "SNESNASM only supports basic and restrict types"); 286e0331734SPeter Brune nasm->type = type; 2873ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 288e0331734SPeter Brune } 289e0331734SPeter Brune 290e0331734SPeter Brune /*@ 291420bcc1bSBarry Smith SNESNASMGetType - Get the type of subdomain update used for the nonlinear additive Schwarz solver `SNESNASM` 292e0331734SPeter Brune 293c3339decSBarry Smith Logically Collective 294e0331734SPeter Brune 2952fe279fdSBarry Smith Input Parameter: 296f6dfbefdSBarry Smith . snes - the `SNES` context 297e0331734SPeter Brune 2982fe279fdSBarry Smith Output Parameter: 299e0331734SPeter Brune . type - the type of update 300e0331734SPeter Brune 301e0331734SPeter Brune Level: intermediate 302e0331734SPeter Brune 303420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESNASM`, `SNESNASMSetType()`, `PCASMGetType()`, `PC_ASM_BASIC`, `PC_ASM_RESTRICT`, `PCASMType` 304e0331734SPeter Brune @*/ 305d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESNASMGetType(SNES snes, PCASMType *type) 306d71ae5a4SJacob Faibussowitsch { 307e0331734SPeter Brune PetscFunctionBegin; 308cac4c232SBarry Smith PetscUseMethod(snes, "SNESNASMGetType_C", (SNES, PCASMType *), (snes, type)); 3093ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 310e0331734SPeter Brune } 311e0331734SPeter Brune 312d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESNASMGetType_NASM(SNES snes, PCASMType *type) 313d71ae5a4SJacob Faibussowitsch { 314e0331734SPeter Brune SNES_NASM *nasm = (SNES_NASM *)snes->data; 315e0331734SPeter Brune 316e0331734SPeter Brune PetscFunctionBegin; 317e0331734SPeter Brune *type = nasm->type; 3183ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 319e0331734SPeter Brune } 320e0331734SPeter Brune 32176857b2aSPeter Brune /*@ 322f6dfbefdSBarry Smith SNESNASMSetSubdomains - Manually Set the context required to restrict and solve subdomain problems in the nonlinear additive Schwarz solver 32376857b2aSPeter Brune 324f6dfbefdSBarry Smith Logically Collective 32576857b2aSPeter Brune 32676857b2aSPeter Brune Input Parameters: 327420bcc1bSBarry Smith + snes - the `SNES` context 32876857b2aSPeter Brune . n - the number of local subdomains 32976857b2aSPeter Brune . subsnes - solvers defined on the local subdomains 33076857b2aSPeter Brune . iscatter - scatters into the nonoverlapping portions of the local subdomains 33176857b2aSPeter Brune . oscatter - scatters into the overlapping portions of the local subdomains 33276857b2aSPeter Brune - gscatter - scatters into the (ghosted) local vector of the local subdomain 33376857b2aSPeter Brune 33476857b2aSPeter Brune Level: intermediate 33576857b2aSPeter Brune 336420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESNASM`, `SNESNASMGetSubdomains()` 33776857b2aSPeter Brune @*/ 338d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESNASMSetSubdomains(SNES snes, PetscInt n, SNES subsnes[], VecScatter iscatter[], VecScatter oscatter[], VecScatter gscatter[]) 339d71ae5a4SJacob Faibussowitsch { 340111ade9eSPeter Brune PetscErrorCode (*f)(SNES, PetscInt, SNES *, VecScatter *, VecScatter *, VecScatter *); 3416e111a19SKarl Rupp 342eaedb033SPeter Brune PetscFunctionBegin; 3439566063dSJacob Faibussowitsch PetscCall(PetscObjectQueryFunction((PetscObject)snes, "SNESNASMSetSubdomains_C", &f)); 3449566063dSJacob Faibussowitsch if (f) PetscCall((f)(snes, n, subsnes, iscatter, oscatter, gscatter)); 3453ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 346eaedb033SPeter Brune } 347eaedb033SPeter Brune 348d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESNASMSetSubdomains_NASM(SNES snes, PetscInt n, SNES subsnes[], VecScatter iscatter[], VecScatter oscatter[], VecScatter gscatter[]) 349d71ae5a4SJacob Faibussowitsch { 350eaedb033SPeter Brune PetscInt i; 351eaedb033SPeter Brune SNES_NASM *nasm = (SNES_NASM *)snes->data; 3526e111a19SKarl Rupp 353eaedb033SPeter Brune PetscFunctionBegin; 35428b400f6SJacob Faibussowitsch PetscCheck(!snes->setupcalled, PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_WRONGSTATE, "SNESNASMSetSubdomains() should be called before calling SNESSetUp()."); 355eaedb033SPeter Brune 356111ade9eSPeter Brune /* tear down the previously set things */ 3579566063dSJacob Faibussowitsch PetscCall(SNESReset(snes)); 358111ade9eSPeter Brune 359eaedb033SPeter Brune nasm->n = n; 360111ade9eSPeter Brune if (oscatter) { 3619566063dSJacob Faibussowitsch for (i = 0; i < n; i++) PetscCall(PetscObjectReference((PetscObject)oscatter[i])); 362eaedb033SPeter Brune } 363111ade9eSPeter Brune if (iscatter) { 3649566063dSJacob Faibussowitsch for (i = 0; i < n; i++) PetscCall(PetscObjectReference((PetscObject)iscatter[i])); 365eaedb033SPeter Brune } 366111ade9eSPeter Brune if (gscatter) { 3679566063dSJacob Faibussowitsch for (i = 0; i < n; i++) PetscCall(PetscObjectReference((PetscObject)gscatter[i])); 368111ade9eSPeter Brune } 369111ade9eSPeter Brune if (oscatter) { 3709566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &nasm->oscatter)); 3719566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &nasm->oscatter_copy)); 372eaedb033SPeter Brune for (i = 0; i < n; i++) { 373111ade9eSPeter Brune nasm->oscatter[i] = oscatter[i]; 3749566063dSJacob Faibussowitsch PetscCall(VecScatterCopy(oscatter[i], &nasm->oscatter_copy[i])); 375eaedb033SPeter Brune } 376111ade9eSPeter Brune } 377111ade9eSPeter Brune if (iscatter) { 3789566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &nasm->iscatter)); 379ad540459SPierre Jolivet for (i = 0; i < n; i++) nasm->iscatter[i] = iscatter[i]; 380eaedb033SPeter Brune } 381111ade9eSPeter Brune if (gscatter) { 3829566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &nasm->gscatter)); 383ad540459SPierre Jolivet for (i = 0; i < n; i++) nasm->gscatter[i] = gscatter[i]; 384eaedb033SPeter Brune } 385111ade9eSPeter Brune 386eaedb033SPeter Brune if (subsnes) { 3879566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &nasm->subsnes)); 388ad540459SPierre Jolivet for (i = 0; i < n; i++) nasm->subsnes[i] = subsnes[i]; 389eaedb033SPeter Brune } 3903ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 391eaedb033SPeter Brune } 392eaedb033SPeter Brune 39376857b2aSPeter Brune /*@ 394f6dfbefdSBarry Smith SNESNASMGetSubdomains - Get the local subdomain contexts for the nonlinear additive Schwarz solver 39576857b2aSPeter Brune 396f6dfbefdSBarry Smith Not Collective but some of the objects returned will be parallel 39776857b2aSPeter Brune 398f899ff85SJose E. Roman Input Parameter: 399f6dfbefdSBarry Smith . snes - the `SNES` context 40076857b2aSPeter Brune 40176857b2aSPeter Brune Output Parameters: 40276857b2aSPeter Brune + n - the number of local subdomains 40376857b2aSPeter Brune . subsnes - solvers defined on the local subdomains 40476857b2aSPeter Brune . iscatter - scatters into the nonoverlapping portions of the local subdomains 40576857b2aSPeter Brune . oscatter - scatters into the overlapping portions of the local subdomains 40676857b2aSPeter Brune - gscatter - scatters into the (ghosted) local vector of the local subdomain 40776857b2aSPeter Brune 40876857b2aSPeter Brune Level: intermediate 40976857b2aSPeter Brune 410420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESNASM`, `SNESNASMSetSubdomains()` 41176857b2aSPeter Brune @*/ 412d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESNASMGetSubdomains(SNES snes, PetscInt *n, SNES *subsnes[], VecScatter *iscatter[], VecScatter *oscatter[], VecScatter *gscatter[]) 413d71ae5a4SJacob Faibussowitsch { 41476857b2aSPeter Brune PetscErrorCode (*f)(SNES, PetscInt *, SNES **, VecScatter **, VecScatter **, VecScatter **); 41576857b2aSPeter Brune 41676857b2aSPeter Brune PetscFunctionBegin; 4179566063dSJacob Faibussowitsch PetscCall(PetscObjectQueryFunction((PetscObject)snes, "SNESNASMGetSubdomains_C", &f)); 4189566063dSJacob Faibussowitsch if (f) PetscCall((f)(snes, n, subsnes, iscatter, oscatter, gscatter)); 4193ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 42076857b2aSPeter Brune } 42176857b2aSPeter Brune 422d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESNASMGetSubdomains_NASM(SNES snes, PetscInt *n, SNES *subsnes[], VecScatter *iscatter[], VecScatter *oscatter[], VecScatter *gscatter[]) 423d71ae5a4SJacob Faibussowitsch { 42476857b2aSPeter Brune SNES_NASM *nasm = (SNES_NASM *)snes->data; 42576857b2aSPeter Brune 42676857b2aSPeter Brune PetscFunctionBegin; 42776857b2aSPeter Brune if (n) *n = nasm->n; 42876857b2aSPeter Brune if (oscatter) *oscatter = nasm->oscatter; 42976857b2aSPeter Brune if (iscatter) *iscatter = nasm->iscatter; 43076857b2aSPeter Brune if (gscatter) *gscatter = nasm->gscatter; 43183dc3634SPierre Jolivet if (subsnes) *subsnes = nasm->subsnes; 4323ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 43376857b2aSPeter Brune } 43476857b2aSPeter Brune 43576857b2aSPeter Brune /*@ 436f6dfbefdSBarry Smith SNESNASMGetSubdomainVecs - Get the processor-local subdomain vectors for the nonlinear additive Schwarz solver 43776857b2aSPeter Brune 43876857b2aSPeter Brune Not Collective 43976857b2aSPeter Brune 440f899ff85SJose E. Roman Input Parameter: 441f6dfbefdSBarry Smith . snes - the `SNES` context 44276857b2aSPeter Brune 44376857b2aSPeter Brune Output Parameters: 44476857b2aSPeter Brune + n - the number of local subdomains 44576857b2aSPeter Brune . x - The subdomain solution vector 44676857b2aSPeter Brune . y - The subdomain step vector 44776857b2aSPeter Brune . b - The subdomain RHS vector 44876857b2aSPeter Brune - xl - The subdomain local vectors (ghosted) 44976857b2aSPeter Brune 45076857b2aSPeter Brune Level: developer 45176857b2aSPeter Brune 452420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESNASM`, `SNESNASMGetSubdomains()` 45376857b2aSPeter Brune @*/ 454d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESNASMGetSubdomainVecs(SNES snes, PetscInt *n, Vec **x, Vec **y, Vec **b, Vec **xl) 455d71ae5a4SJacob Faibussowitsch { 45676857b2aSPeter Brune PetscErrorCode (*f)(SNES, PetscInt *, Vec **, Vec **, Vec **, Vec **); 45776857b2aSPeter Brune 45876857b2aSPeter Brune PetscFunctionBegin; 4599566063dSJacob Faibussowitsch PetscCall(PetscObjectQueryFunction((PetscObject)snes, "SNESNASMGetSubdomainVecs_C", &f)); 4609566063dSJacob Faibussowitsch if (f) PetscCall((f)(snes, n, x, y, b, xl)); 4613ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 46276857b2aSPeter Brune } 46376857b2aSPeter Brune 464d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESNASMGetSubdomainVecs_NASM(SNES snes, PetscInt *n, Vec **x, Vec **y, Vec **b, Vec **xl) 465d71ae5a4SJacob Faibussowitsch { 46676857b2aSPeter Brune SNES_NASM *nasm = (SNES_NASM *)snes->data; 46776857b2aSPeter Brune 46876857b2aSPeter Brune PetscFunctionBegin; 46976857b2aSPeter Brune if (n) *n = nasm->n; 47076857b2aSPeter Brune if (x) *x = nasm->x; 47176857b2aSPeter Brune if (y) *y = nasm->y; 47276857b2aSPeter Brune if (b) *b = nasm->b; 47376857b2aSPeter Brune if (xl) *xl = nasm->xl; 4743ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 47576857b2aSPeter Brune } 47676857b2aSPeter Brune 477d728fb7dSPeter Brune /*@ 478f6dfbefdSBarry Smith SNESNASMSetComputeFinalJacobian - Schedules the computation of the global and subdomain Jacobians upon convergence for the 479f6dfbefdSBarry Smith nonlinear additive Schwarz solver 480d728fb7dSPeter Brune 481c3339decSBarry Smith Collective 482d728fb7dSPeter Brune 483d728fb7dSPeter Brune Input Parameters: 484f6dfbefdSBarry Smith + snes - the SNES context 485f6dfbefdSBarry Smith - flg - `PETSC_TRUE` to compute the Jacobians 486d728fb7dSPeter Brune 487d728fb7dSPeter Brune Level: developer 488d728fb7dSPeter Brune 48995452b02SPatrick Sanan Notes: 490f6dfbefdSBarry Smith This is used almost exclusively in the implementation of `SNESASPIN`, where the converged subdomain and global Jacobian 491d728fb7dSPeter Brune is needed at each linear iteration. 492d728fb7dSPeter Brune 493420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESNASM`, `SNESNASMGetSubdomains()` 494d728fb7dSPeter Brune @*/ 495d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESNASMSetComputeFinalJacobian(SNES snes, PetscBool flg) 496d71ae5a4SJacob Faibussowitsch { 497d728fb7dSPeter Brune PetscErrorCode (*f)(SNES, PetscBool); 498d728fb7dSPeter Brune 499d728fb7dSPeter Brune PetscFunctionBegin; 5009566063dSJacob Faibussowitsch PetscCall(PetscObjectQueryFunction((PetscObject)snes, "SNESNASMSetComputeFinalJacobian_C", &f)); 5019566063dSJacob Faibussowitsch if (f) PetscCall((f)(snes, flg)); 5023ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 503d728fb7dSPeter Brune } 504d728fb7dSPeter Brune 505d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESNASMSetComputeFinalJacobian_NASM(SNES snes, PetscBool flg) 506d71ae5a4SJacob Faibussowitsch { 507d728fb7dSPeter Brune SNES_NASM *nasm = (SNES_NASM *)snes->data; 508d728fb7dSPeter Brune 509d728fb7dSPeter Brune PetscFunctionBegin; 510d728fb7dSPeter Brune nasm->finaljacobian = flg; 5113ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 512d728fb7dSPeter Brune } 51376857b2aSPeter Brune 514610116beSPeter Brune /*@ 515f6dfbefdSBarry Smith SNESNASMSetDamping - Sets the update damping for `SNESNASM` the nonlinear additive Schwarz solver 516610116beSPeter Brune 51720f4b53cSBarry Smith Logically Collective 518610116beSPeter Brune 519610116beSPeter Brune Input Parameters: 520f6dfbefdSBarry Smith + snes - the `SNES` context 521610116beSPeter Brune - dmp - damping 522610116beSPeter Brune 523420bcc1bSBarry Smith Options Database Key: 524420bcc1bSBarry Smith . -snes_nasm_damping <dmp> - the new solution is obtained as old solution plus `dmp` times (sum of the solutions on the subdomains) 525420bcc1bSBarry Smith 526610116beSPeter Brune Level: intermediate 527610116beSPeter Brune 528420bcc1bSBarry Smith Note: 52995452b02SPatrick Sanan The new solution is obtained as old solution plus dmp times (sum of the solutions on the subdomains) 5305dfa0f3bSBarry Smith 531420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESNASM`, `SNESNASMGetDamping()` 532610116beSPeter Brune @*/ 533d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESNASMSetDamping(SNES snes, PetscReal dmp) 534d71ae5a4SJacob Faibussowitsch { 535610116beSPeter Brune PetscErrorCode (*f)(SNES, PetscReal); 536610116beSPeter Brune 537610116beSPeter Brune PetscFunctionBegin; 5389566063dSJacob Faibussowitsch PetscCall(PetscObjectQueryFunction((PetscObject)snes, "SNESNASMSetDamping_C", (void (**)(void)) & f)); 5399566063dSJacob Faibussowitsch if (f) PetscCall((f)(snes, dmp)); 5403ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 541610116beSPeter Brune } 542610116beSPeter Brune 543d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESNASMSetDamping_NASM(SNES snes, PetscReal dmp) 544d71ae5a4SJacob Faibussowitsch { 545610116beSPeter Brune SNES_NASM *nasm = (SNES_NASM *)snes->data; 546610116beSPeter Brune 547610116beSPeter Brune PetscFunctionBegin; 548610116beSPeter Brune nasm->damping = dmp; 5493ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 550610116beSPeter Brune } 551610116beSPeter Brune 552610116beSPeter Brune /*@ 553f6dfbefdSBarry Smith SNESNASMGetDamping - Gets the update damping for `SNESNASM` the nonlinear additive Schwarz solver 554610116beSPeter Brune 555610116beSPeter Brune Not Collective 556610116beSPeter Brune 557f6dfbefdSBarry Smith Input Parameter: 558420bcc1bSBarry Smith . snes - the `SNES` context 559f6dfbefdSBarry Smith 560f6dfbefdSBarry Smith Output Parameter: 561f6dfbefdSBarry Smith . dmp - damping 562610116beSPeter Brune 563610116beSPeter Brune Level: intermediate 564610116beSPeter Brune 565420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESNASM`, `SNESNASMSetDamping()` 566610116beSPeter Brune @*/ 567d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESNASMGetDamping(SNES snes, PetscReal *dmp) 568d71ae5a4SJacob Faibussowitsch { 569610116beSPeter Brune PetscFunctionBegin; 570cac4c232SBarry Smith PetscUseMethod(snes, "SNESNASMGetDamping_C", (SNES, PetscReal *), (snes, dmp)); 5713ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 572610116beSPeter Brune } 573610116beSPeter Brune 574d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESNASMGetDamping_NASM(SNES snes, PetscReal *dmp) 575d71ae5a4SJacob Faibussowitsch { 576610116beSPeter Brune SNES_NASM *nasm = (SNES_NASM *)snes->data; 577610116beSPeter Brune 578610116beSPeter Brune PetscFunctionBegin; 579610116beSPeter Brune *dmp = nasm->damping; 5803ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 581610116beSPeter Brune } 582610116beSPeter Brune 58314eb1c5cSMatthew G. Knepley /* 58414eb1c5cSMatthew G. Knepley Input Parameters: 58514eb1c5cSMatthew G. Knepley + snes - The solver 58614eb1c5cSMatthew G. Knepley . B - The RHS vector 58714eb1c5cSMatthew G. Knepley - X - The initial guess 58814eb1c5cSMatthew G. Knepley 5892fe279fdSBarry Smith Output Parameter: 59014eb1c5cSMatthew G. Knepley . Y - The solution update 59114eb1c5cSMatthew G. Knepley 59214eb1c5cSMatthew G. Knepley TODO: All scatters should be packed into one 59314eb1c5cSMatthew G. Knepley */ 59466976f2fSJacob Faibussowitsch static PetscErrorCode SNESNASMSolveLocal_Private(SNES snes, Vec B, Vec Y, Vec X) 595d71ae5a4SJacob Faibussowitsch { 596eaedb033SPeter Brune SNES_NASM *nasm = (SNES_NASM *)snes->data; 597258e1594SPeter Brune SNES subsnes; 598eaedb033SPeter Brune PetscInt i; 599610116beSPeter Brune PetscReal dmp; 600f10b3e88SPatrick Farrell Vec Xl, Bl, Yl, Xlloc; 601f10b3e88SPatrick Farrell VecScatter iscat, oscat, gscat, oscat_copy; 602111ade9eSPeter Brune DM dm, subdm; 603e0331734SPeter Brune PCASMType type; 6040adebc6cSBarry Smith 605eaedb033SPeter Brune PetscFunctionBegin; 6069566063dSJacob Faibussowitsch PetscCall(SNESNASMGetType(snes, &type)); 6079566063dSJacob Faibussowitsch PetscCall(SNESGetDM(snes, &dm)); 6089566063dSJacob Faibussowitsch PetscCall(VecSet(Y, 0)); 6099566063dSJacob Faibussowitsch if (nasm->eventrestrictinterp) PetscCall(PetscLogEventBegin(nasm->eventrestrictinterp, snes, 0, 0, 0)); 610de54d9edSStefano Zampini for (i = 0; i < nasm->n; i++) { /* scatter the global solution to the overlap solution and the local solution */ 611f10b3e88SPatrick Farrell Xl = nasm->x[i]; 61270c78f05SPeter Brune Xlloc = nasm->xl[i]; 61370c78f05SPeter Brune oscat = nasm->oscatter[i]; 614f10b3e88SPatrick Farrell oscat_copy = nasm->oscatter_copy[i]; 615f10b3e88SPatrick Farrell gscat = nasm->gscatter[i]; 6169566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(oscat, X, Xl, INSERT_VALUES, SCATTER_FORWARD)); 6179566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(gscat, X, Xlloc, INSERT_VALUES, SCATTER_FORWARD)); 618de54d9edSStefano Zampini 61970c78f05SPeter Brune if (B) { 62070c78f05SPeter Brune /* scatter the RHS to the local RHS */ 62170c78f05SPeter Brune Bl = nasm->b[i]; 6229566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(oscat_copy, B, Bl, INSERT_VALUES, SCATTER_FORWARD)); 62370c78f05SPeter Brune } 62470c78f05SPeter Brune } 6259566063dSJacob Faibussowitsch if (nasm->eventrestrictinterp) PetscCall(PetscLogEventEnd(nasm->eventrestrictinterp, snes, 0, 0, 0)); 626b20c023fSPeter Brune 6279566063dSJacob Faibussowitsch if (nasm->eventsubsolve) PetscCall(PetscLogEventBegin(nasm->eventsubsolve, snes, 0, 0, 0)); 62870c78f05SPeter Brune for (i = 0; i < nasm->n; i++) { 629de54d9edSStefano Zampini PetscErrorCode (*bl)(DM, Vec, void *); 630de54d9edSStefano Zampini void *bctx; 631de54d9edSStefano Zampini 63270c78f05SPeter Brune Xl = nasm->x[i]; 63370c78f05SPeter Brune Xlloc = nasm->xl[i]; 63470c78f05SPeter Brune Yl = nasm->y[i]; 635258e1594SPeter Brune subsnes = nasm->subsnes[i]; 6369566063dSJacob Faibussowitsch PetscCall(SNESGetDM(subsnes, &subdm)); 637111ade9eSPeter Brune iscat = nasm->iscatter[i]; 638111ade9eSPeter Brune oscat = nasm->oscatter[i]; 639f10b3e88SPatrick Farrell oscat_copy = nasm->oscatter_copy[i]; 640111ade9eSPeter Brune gscat = nasm->gscatter[i]; 6419566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(oscat, X, Xl, INSERT_VALUES, SCATTER_FORWARD)); 6429566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(gscat, X, Xlloc, INSERT_VALUES, SCATTER_FORWARD)); 64324b7f281SPeter Brune if (B) { 64424b7f281SPeter Brune Bl = nasm->b[i]; 6459566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(oscat_copy, B, Bl, INSERT_VALUES, SCATTER_FORWARD)); 646ed3c11a9SPeter Brune } else Bl = NULL; 647f10b3e88SPatrick Farrell 648de54d9edSStefano Zampini PetscCall(SNESGetDM(subsnes, &subdm)); 649de54d9edSStefano Zampini PetscCall(DMSNESGetBoundaryLocal(subdm, &bl, &bctx)); 650de54d9edSStefano Zampini if (bl) PetscCall((*bl)(subdm, Xlloc, bctx)); 651de54d9edSStefano Zampini 6529566063dSJacob Faibussowitsch PetscCall(DMSubDomainRestrict(dm, oscat, gscat, subdm)); 6539566063dSJacob Faibussowitsch PetscCall(VecCopy(Xl, Yl)); 6549566063dSJacob Faibussowitsch PetscCall(SNESSolve(subsnes, Bl, Xl)); 6559566063dSJacob Faibussowitsch PetscCall(VecAYPX(Yl, -1.0, Xl)); 6569566063dSJacob Faibussowitsch PetscCall(VecScale(Yl, nasm->damping)); 657e0331734SPeter Brune if (type == PC_ASM_BASIC) { 6589566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(oscat, Yl, Y, ADD_VALUES, SCATTER_REVERSE)); 6599566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(oscat, Yl, Y, ADD_VALUES, SCATTER_REVERSE)); 660e0331734SPeter Brune } else if (type == PC_ASM_RESTRICT) { 6619566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(iscat, Yl, Y, ADD_VALUES, SCATTER_REVERSE)); 6629566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(iscat, Yl, Y, ADD_VALUES, SCATTER_REVERSE)); 663ce94432eSBarry Smith } else SETERRQ(PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_WRONGSTATE, "Only basic and restrict types are supported for SNESNASM"); 664eaedb033SPeter Brune } 6659566063dSJacob Faibussowitsch if (nasm->eventsubsolve) PetscCall(PetscLogEventEnd(nasm->eventsubsolve, snes, 0, 0, 0)); 6669566063dSJacob Faibussowitsch if (nasm->eventrestrictinterp) PetscCall(PetscLogEventBegin(nasm->eventrestrictinterp, snes, 0, 0, 0)); 6671baa6e33SBarry Smith if (nasm->weight_set) PetscCall(VecPointwiseMult(Y, Y, nasm->weight)); 6689566063dSJacob Faibussowitsch if (nasm->eventrestrictinterp) PetscCall(PetscLogEventEnd(nasm->eventrestrictinterp, snes, 0, 0, 0)); 6699566063dSJacob Faibussowitsch PetscCall(SNESNASMGetDamping(snes, &dmp)); 6709566063dSJacob Faibussowitsch PetscCall(VecAXPY(X, dmp, Y)); 6713ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 672eaedb033SPeter Brune } 673eaedb033SPeter Brune 674d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESNASMComputeFinalJacobian_Private(SNES snes, Vec Xfinal) 675d71ae5a4SJacob Faibussowitsch { 676602bec5dSPeter Brune Vec X = Xfinal; 677d728fb7dSPeter Brune SNES_NASM *nasm = (SNES_NASM *)snes->data; 678d728fb7dSPeter Brune SNES subsnes; 679ca44f815SPeter Brune PetscInt i, lag = 1; 680e59f0a30SPeter Brune Vec Xlloc, Xl, Fl, F; 681d728fb7dSPeter Brune VecScatter oscat, gscat; 682d728fb7dSPeter Brune DM dm, subdm; 683d1e9a80fSBarry Smith 684d728fb7dSPeter Brune PetscFunctionBegin; 685602bec5dSPeter Brune if (nasm->fjtype == 2) X = nasm->xinit; 686e59f0a30SPeter Brune F = snes->vec_func; 6879566063dSJacob Faibussowitsch if (snes->normschedule == SNES_NORM_NONE) PetscCall(SNESComputeFunction(snes, X, F)); 6889566063dSJacob Faibussowitsch PetscCall(SNESComputeJacobian(snes, X, snes->jacobian, snes->jacobian_pre)); 6899566063dSJacob Faibussowitsch PetscCall(SNESGetDM(snes, &dm)); 6909566063dSJacob Faibussowitsch if (nasm->eventrestrictinterp) PetscCall(PetscLogEventBegin(nasm->eventrestrictinterp, snes, 0, 0, 0)); 691602bec5dSPeter Brune if (nasm->fjtype != 1) { 692d728fb7dSPeter Brune for (i = 0; i < nasm->n; i++) { 693d728fb7dSPeter Brune Xlloc = nasm->xl[i]; 694d728fb7dSPeter Brune gscat = nasm->gscatter[i]; 6959566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(gscat, X, Xlloc, INSERT_VALUES, SCATTER_FORWARD)); 696602bec5dSPeter Brune } 697d728fb7dSPeter Brune } 6989566063dSJacob Faibussowitsch if (nasm->eventrestrictinterp) PetscCall(PetscLogEventEnd(nasm->eventrestrictinterp, snes, 0, 0, 0)); 699d728fb7dSPeter Brune for (i = 0; i < nasm->n; i++) { 700e59f0a30SPeter Brune Fl = nasm->subsnes[i]->vec_func; 701d728fb7dSPeter Brune Xl = nasm->x[i]; 702d728fb7dSPeter Brune Xlloc = nasm->xl[i]; 703d728fb7dSPeter Brune subsnes = nasm->subsnes[i]; 704d728fb7dSPeter Brune oscat = nasm->oscatter[i]; 705d728fb7dSPeter Brune gscat = nasm->gscatter[i]; 7069566063dSJacob Faibussowitsch if (nasm->fjtype != 1) PetscCall(VecScatterEnd(gscat, X, Xlloc, INSERT_VALUES, SCATTER_FORWARD)); 7079566063dSJacob Faibussowitsch PetscCall(SNESGetDM(subsnes, &subdm)); 7089566063dSJacob Faibussowitsch PetscCall(DMSubDomainRestrict(dm, oscat, gscat, subdm)); 709602bec5dSPeter Brune if (nasm->fjtype != 1) { 7109566063dSJacob Faibussowitsch PetscCall(DMLocalToGlobalBegin(subdm, Xlloc, INSERT_VALUES, Xl)); 7119566063dSJacob Faibussowitsch PetscCall(DMLocalToGlobalEnd(subdm, Xlloc, INSERT_VALUES, Xl)); 712602bec5dSPeter Brune } 713ca44f815SPeter Brune if (subsnes->lagjacobian == -1) subsnes->lagjacobian = -2; 714ca44f815SPeter Brune else if (subsnes->lagjacobian > 1) lag = subsnes->lagjacobian; 7159566063dSJacob Faibussowitsch PetscCall(SNESComputeFunction(subsnes, Xl, Fl)); 7169566063dSJacob Faibussowitsch PetscCall(SNESComputeJacobian(subsnes, Xl, subsnes->jacobian, subsnes->jacobian_pre)); 717ca44f815SPeter Brune if (lag > 1) subsnes->lagjacobian = lag; 718d728fb7dSPeter Brune } 7193ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 720d728fb7dSPeter Brune } 721d728fb7dSPeter Brune 722d71ae5a4SJacob Faibussowitsch static PetscErrorCode SNESSolve_NASM(SNES snes) 723d71ae5a4SJacob Faibussowitsch { 724eaedb033SPeter Brune Vec F; 725eaedb033SPeter Brune Vec X; 726eaedb033SPeter Brune Vec B; 727111ade9eSPeter Brune Vec Y; 728eaedb033SPeter Brune PetscInt i; 729ed3c11a9SPeter Brune PetscReal fnorm = 0.0; 730365a6726SPeter Brune SNESNormSchedule normschedule; 731d728fb7dSPeter Brune SNES_NASM *nasm = (SNES_NASM *)snes->data; 732eaedb033SPeter Brune 733eaedb033SPeter Brune PetscFunctionBegin; 734c579b300SPatrick Farrell 7350b121fc5SBarry Smith PetscCheck(!snes->xl & !snes->xu && !snes->ops->computevariablebounds, PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_WRONGSTATE, "SNES solver %s does not support bounds", ((PetscObject)snes)->type_name); 736c579b300SPatrick Farrell 7379566063dSJacob Faibussowitsch PetscCall(PetscCitationsRegister(SNESCitation, &SNEScite)); 738eaedb033SPeter Brune X = snes->vec_sol; 739111ade9eSPeter Brune Y = snes->vec_sol_update; 740eaedb033SPeter Brune F = snes->vec_func; 741eaedb033SPeter Brune B = snes->vec_rhs; 742eaedb033SPeter Brune 7439566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes)); 744eaedb033SPeter Brune snes->iter = 0; 745eaedb033SPeter Brune snes->norm = 0.; 7469566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 747eaedb033SPeter Brune snes->reason = SNES_CONVERGED_ITERATING; 7489566063dSJacob Faibussowitsch PetscCall(SNESGetNormSchedule(snes, &normschedule)); 749de54d9edSStefano Zampini if (normschedule == SNES_NORM_ALWAYS || normschedule == SNES_NORM_INITIAL_ONLY || normschedule == SNES_NORM_INITIAL_FINAL_ONLY || !snes->max_its) { 750eaedb033SPeter Brune /* compute the initial function and preconditioned update delX */ 751eaedb033SPeter Brune if (!snes->vec_func_init_set) { 7529566063dSJacob Faibussowitsch PetscCall(SNESComputeFunction(snes, X, F)); 7531aa26658SKarl Rupp } else snes->vec_func_init_set = PETSC_FALSE; 754eaedb033SPeter Brune 7559566063dSJacob Faibussowitsch PetscCall(VecNorm(F, NORM_2, &fnorm)); /* fnorm <- ||F|| */ 756422a814eSBarry Smith SNESCheckFunctionNorm(snes, fnorm); 7579566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes)); 758eaedb033SPeter Brune snes->iter = 0; 759eaedb033SPeter Brune snes->norm = fnorm; 7609566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 7619566063dSJacob Faibussowitsch PetscCall(SNESLogConvergenceHistory(snes, snes->norm, 0)); 762eaedb033SPeter Brune 763eaedb033SPeter Brune /* test convergence */ 7642d157150SStefano Zampini PetscCall(SNESConverged(snes, 0, 0.0, 0.0, fnorm)); 7652d157150SStefano Zampini PetscCall(SNESMonitor(snes, 0, snes->norm)); 7663ba16761SJacob Faibussowitsch if (snes->reason) PetscFunctionReturn(PETSC_SUCCESS); 767eaedb033SPeter Brune } else { 7689566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 7699566063dSJacob Faibussowitsch PetscCall(SNESLogConvergenceHistory(snes, snes->norm, 0)); 7702d157150SStefano Zampini PetscCall(SNESMonitor(snes, snes->iter, snes->norm)); 771eaedb033SPeter Brune } 772eaedb033SPeter Brune 773eaedb033SPeter Brune /* Call general purpose update function */ 774dbbe0bcdSBarry Smith PetscTryTypeMethod(snes, update, snes->iter); 775602bec5dSPeter Brune /* copy the initial solution over for later */ 7769566063dSJacob Faibussowitsch if (nasm->fjtype == 2) PetscCall(VecCopy(X, nasm->xinit)); 777eaedb033SPeter Brune 778eaedb033SPeter Brune for (i = 0; i < snes->max_its; i++) { 7799566063dSJacob Faibussowitsch PetscCall(SNESNASMSolveLocal_Private(snes, B, Y, X)); 780365a6726SPeter Brune if (normschedule == SNES_NORM_ALWAYS || ((i == snes->max_its - 1) && (normschedule == SNES_NORM_INITIAL_FINAL_ONLY || normschedule == SNES_NORM_FINAL_ONLY))) { 7819566063dSJacob Faibussowitsch PetscCall(SNESComputeFunction(snes, X, F)); 7829566063dSJacob Faibussowitsch PetscCall(VecNorm(F, NORM_2, &fnorm)); /* fnorm <- ||F|| */ 783422a814eSBarry Smith SNESCheckFunctionNorm(snes, fnorm); 784eaedb033SPeter Brune } 785eaedb033SPeter Brune /* Monitor convergence */ 7869566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes)); 787eaedb033SPeter Brune snes->iter = i + 1; 788eaedb033SPeter Brune snes->norm = fnorm; 7899566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 7909566063dSJacob Faibussowitsch PetscCall(SNESLogConvergenceHistory(snes, snes->norm, 0)); 791eaedb033SPeter Brune /* Test for convergence */ 7922d157150SStefano Zampini PetscCall(SNESConverged(snes, snes->iter, 0.0, 0.0, fnorm)); 7932d157150SStefano Zampini PetscCall(SNESMonitor(snes, snes->iter, snes->norm)); 794d728fb7dSPeter Brune if (snes->reason) break; 795eaedb033SPeter Brune /* Call general purpose update function */ 796dbbe0bcdSBarry Smith PetscTryTypeMethod(snes, update, snes->iter); 797eaedb033SPeter Brune } 79807b62357SFande Kong if (nasm->finaljacobian) { 7999566063dSJacob Faibussowitsch PetscCall(SNESNASMComputeFinalJacobian_Private(snes, X)); 80007b62357SFande Kong SNESCheckJacobianDomainerror(snes); 80107b62357SFande Kong } 8023ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 803eaedb033SPeter Brune } 804eaedb033SPeter Brune 805eaedb033SPeter Brune /*MC 8061d27aa22SBarry Smith SNESNASM - Nonlinear Additive Schwarz solver {cite}`ck02`, {cite}`bruneknepleysmithtu15` 807eaedb033SPeter Brune 808f6dfbefdSBarry Smith Options Database Keys: 80970c78603SPeter Brune + -snes_nasm_log - enable logging events for the communication and solve stages 81070c78603SPeter Brune . -snes_nasm_type <basic,restrict> - type of subdomain update used 811420bcc1bSBarry Smith . -snes_nasm_damping <dmp> - the new solution is obtained as old solution plus dmp times (sum of the solutions on the subdomains) 8121d27aa22SBarry Smith . -snes_nasm_finaljacobian - compute the local and global Jacobians of the final iterate 8131d27aa22SBarry Smith . -snes_nasm_finaljacobian_type <finalinner,finalouter,initial> - pick state the Jacobian is calculated at 81470c78603SPeter Brune . -sub_snes_ - options prefix of the subdomain nonlinear solves 81570c78603SPeter Brune . -sub_ksp_ - options prefix of the subdomain Krylov solver 81670c78603SPeter Brune - -sub_pc_ - options prefix of the subdomain preconditioner 81770c78603SPeter Brune 818eaedb033SPeter Brune Level: advanced 819eaedb033SPeter Brune 820f6dfbefdSBarry Smith Note: 821f6dfbefdSBarry Smith This is not often used directly as a solver, it converges too slowly. However it works well as a nonlinear preconditioner for 822f6dfbefdSBarry Smith the `SNESASPIN` solver 823f6dfbefdSBarry Smith 824f6dfbefdSBarry Smith Developer Note: 825f6dfbefdSBarry Smith This is a non-Newton based nonlinear solver that does not directly require a Jacobian; hence the flag snes->usesksp is set to 826f6dfbefdSBarry Smith false and `SNESView()` and -snes_view do not display a `KSP` object. However, if the flag nasm->finaljacobian is set (for example, if 827f6dfbefdSBarry Smith `SNESNASM` is used as a nonlinear preconditioner for `SNESASPIN`) then `SNESSetUpMatrices()` is called to generate the 828f6dfbefdSBarry Smith Jacobian (needed by `SNESASPIN`) 829f6dfbefdSBarry Smith and this utilizes the inner `KSP` object for storing the matrices, but the `KSP` is never used for solving a linear system. When `SNESNASM` is 830f6dfbefdSBarry Smith used by `SNESASPIN` they share the same Jacobian matrices because `SNESSetUp()` (called on the outer `SNESASPIN`) causes the inner `SNES` 831f6dfbefdSBarry Smith object (in this case `SNESNASM`) to inherit the outer Jacobian matrices. 832951fe5abSBarry Smith 833420bcc1bSBarry Smith .seealso: [](ch_snes), `SNESCreate()`, `SNES`, `SNESSetType()`, `SNESType`, `SNESNASMSetType()`, `SNESNASMGetType()`, `SNESNASMSetSubdomains()`, `SNESNASMGetSubdomains()`, 834420bcc1bSBarry Smith `SNESNASMGetSubdomainVecs()`, `SNESNASMSetComputeFinalJacobian()`, `SNESNASMSetDamping()`, `SNESNASMGetDamping()`, `SNESNASMSetWeight()`, 835420bcc1bSBarry Smith `SNESNASMGetSNES()`, `SNESNASMGetNumber()` 836eaedb033SPeter Brune M*/ 837eaedb033SPeter Brune 838d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode SNESCreate_NASM(SNES snes) 839d71ae5a4SJacob Faibussowitsch { 840eaedb033SPeter Brune SNES_NASM *nasm; 841eaedb033SPeter Brune 842eaedb033SPeter Brune PetscFunctionBegin; 8434dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&nasm)); 844eaedb033SPeter Brune snes->data = (void *)nasm; 845eaedb033SPeter Brune 846eaedb033SPeter Brune nasm->n = PETSC_DECIDE; 8479e5d0892SLisandro Dalcin nasm->subsnes = NULL; 8489e5d0892SLisandro Dalcin nasm->x = NULL; 8499e5d0892SLisandro Dalcin nasm->xl = NULL; 8509e5d0892SLisandro Dalcin nasm->y = NULL; 8519e5d0892SLisandro Dalcin nasm->b = NULL; 8529e5d0892SLisandro Dalcin nasm->oscatter = NULL; 8539e5d0892SLisandro Dalcin nasm->oscatter_copy = NULL; 8549e5d0892SLisandro Dalcin nasm->iscatter = NULL; 8559e5d0892SLisandro Dalcin nasm->gscatter = NULL; 856610116beSPeter Brune nasm->damping = 1.; 857111ade9eSPeter Brune 858111ade9eSPeter Brune nasm->type = PC_ASM_BASIC; 859d728fb7dSPeter Brune nasm->finaljacobian = PETSC_FALSE; 860f10b3e88SPatrick Farrell nasm->weight_set = PETSC_FALSE; 861eaedb033SPeter Brune 862eaedb033SPeter Brune snes->ops->destroy = SNESDestroy_NASM; 863eaedb033SPeter Brune snes->ops->setup = SNESSetUp_NASM; 864eaedb033SPeter Brune snes->ops->setfromoptions = SNESSetFromOptions_NASM; 865eaedb033SPeter Brune snes->ops->view = SNESView_NASM; 866eaedb033SPeter Brune snes->ops->solve = SNESSolve_NASM; 867eaedb033SPeter Brune snes->ops->reset = SNESReset_NASM; 868eaedb033SPeter Brune 869eaedb033SPeter Brune snes->usesksp = PETSC_FALSE; 870efd4aadfSBarry Smith snes->usesnpc = PETSC_FALSE; 871eaedb033SPeter Brune 8724fc747eaSLawrence Mitchell snes->alwayscomputesfinalresidual = PETSC_FALSE; 8734fc747eaSLawrence Mitchell 874602bec5dSPeter Brune nasm->fjtype = 0; 875602bec5dSPeter Brune nasm->xinit = NULL; 8760298fd71SBarry Smith nasm->eventrestrictinterp = 0; 8770298fd71SBarry Smith nasm->eventsubsolve = 0; 878b20c023fSPeter Brune 879eaedb033SPeter Brune if (!snes->tolerancesset) { 880eaedb033SPeter Brune snes->max_its = 10000; 881eaedb033SPeter Brune snes->max_funcs = 10000; 882eaedb033SPeter Brune } 883eaedb033SPeter Brune 8849566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNASMSetType_C", SNESNASMSetType_NASM)); 8859566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNASMGetType_C", SNESNASMGetType_NASM)); 8869566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNASMSetSubdomains_C", SNESNASMSetSubdomains_NASM)); 8879566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNASMGetSubdomains_C", SNESNASMGetSubdomains_NASM)); 8889566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNASMSetDamping_C", SNESNASMSetDamping_NASM)); 8899566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNASMGetDamping_C", SNESNASMGetDamping_NASM)); 8909566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNASMGetSubdomainVecs_C", SNESNASMGetSubdomainVecs_NASM)); 8919566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNASMSetComputeFinalJacobian_C", SNESNASMSetComputeFinalJacobian_NASM)); 8923ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 893eaedb033SPeter Brune } 89499e0435eSBarry Smith 895448b6425SPatrick Farrell /*@ 896448b6425SPatrick Farrell SNESNASMGetSNES - Gets a subsolver 897448b6425SPatrick Farrell 89820f4b53cSBarry Smith Not Collective 899448b6425SPatrick Farrell 900448b6425SPatrick Farrell Input Parameters: 90120f4b53cSBarry Smith + snes - the `SNES` context 902448b6425SPatrick Farrell - i - the number of the subsnes to get 903448b6425SPatrick Farrell 9042fe279fdSBarry Smith Output Parameter: 905448b6425SPatrick Farrell . subsnes - the subsolver context 906448b6425SPatrick Farrell 907448b6425SPatrick Farrell Level: intermediate 908448b6425SPatrick Farrell 909420bcc1bSBarry Smith .seealso: [](ch_snes), `SNESNASM`, `SNESNASMGetNumber()` 910448b6425SPatrick Farrell @*/ 911d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESNASMGetSNES(SNES snes, PetscInt i, SNES *subsnes) 912d71ae5a4SJacob Faibussowitsch { 913448b6425SPatrick Farrell SNES_NASM *nasm = (SNES_NASM *)snes->data; 914448b6425SPatrick Farrell 915448b6425SPatrick Farrell PetscFunctionBegin; 9160b121fc5SBarry Smith PetscCheck(i >= 0 && i < nasm->n, PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_OUTOFRANGE, "No such subsolver"); 917448b6425SPatrick Farrell *subsnes = nasm->subsnes[i]; 9183ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 919448b6425SPatrick Farrell } 920448b6425SPatrick Farrell 921448b6425SPatrick Farrell /*@ 922448b6425SPatrick Farrell SNESNASMGetNumber - Gets number of subsolvers 923448b6425SPatrick Farrell 92420f4b53cSBarry Smith Not Collective 925448b6425SPatrick Farrell 9262fe279fdSBarry Smith Input Parameter: 92720f4b53cSBarry Smith . snes - the `SNES` context 928448b6425SPatrick Farrell 9292fe279fdSBarry Smith Output Parameter: 930448b6425SPatrick Farrell . n - the number of subsolvers 931448b6425SPatrick Farrell 932448b6425SPatrick Farrell Level: intermediate 933448b6425SPatrick Farrell 934420bcc1bSBarry Smith .seealso: [](ch_snes), `SNESNASM`, `SNESNASMGetSNES()` 935448b6425SPatrick Farrell @*/ 936d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESNASMGetNumber(SNES snes, PetscInt *n) 937d71ae5a4SJacob Faibussowitsch { 938448b6425SPatrick Farrell SNES_NASM *nasm = (SNES_NASM *)snes->data; 939448b6425SPatrick Farrell 940448b6425SPatrick Farrell PetscFunctionBegin; 941448b6425SPatrick Farrell *n = nasm->n; 9423ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 943448b6425SPatrick Farrell } 944448b6425SPatrick Farrell 945f10b3e88SPatrick Farrell /*@ 946f10b3e88SPatrick Farrell SNESNASMSetWeight - Sets weight to use when adding overlapping updates 947f10b3e88SPatrick Farrell 948f10b3e88SPatrick Farrell Collective 949f10b3e88SPatrick Farrell 950f10b3e88SPatrick Farrell Input Parameters: 95120f4b53cSBarry Smith + snes - the `SNES` context 952f10b3e88SPatrick Farrell - weight - the weights to use (typically 1/N for each dof, where N is the number of patches it appears in) 953f10b3e88SPatrick Farrell 954f10b3e88SPatrick Farrell Level: intermediate 955f10b3e88SPatrick Farrell 956420bcc1bSBarry Smith .seealso: [](ch_snes), `SNESNASM` 957f10b3e88SPatrick Farrell @*/ 958d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESNASMSetWeight(SNES snes, Vec weight) 959d71ae5a4SJacob Faibussowitsch { 960f10b3e88SPatrick Farrell SNES_NASM *nasm = (SNES_NASM *)snes->data; 961f10b3e88SPatrick Farrell 962f10b3e88SPatrick Farrell PetscFunctionBegin; 963f10b3e88SPatrick Farrell 9649566063dSJacob Faibussowitsch PetscCall(VecDestroy(&nasm->weight)); 965f10b3e88SPatrick Farrell nasm->weight_set = PETSC_TRUE; 966f10b3e88SPatrick Farrell nasm->weight = weight; 9679566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)nasm->weight)); 968f10b3e88SPatrick Farrell 9693ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 970f10b3e88SPatrick Farrell } 971