1 #include <petsc/private/snesimpl.h> /*I "petscsnes.h" I*/ 2 #include <petscdm.h> 3 4 typedef struct { 5 PetscInt n; /* local subdomains */ 6 SNES *subsnes; /* nonlinear solvers for each subdomain */ 7 Vec *x; /* solution vectors */ 8 Vec *xl; /* solution local vectors */ 9 Vec *y; /* step vectors */ 10 Vec *b; /* rhs vectors */ 11 Vec weight; /* weighting for adding updates on overlaps, in global space */ 12 VecScatter *oscatter; /* scatter from global space to the subdomain global space */ 13 VecScatter *oscatter_copy; /* copy of the above */ 14 VecScatter *iscatter; /* scatter from global space to the nonoverlapping subdomain space */ 15 VecScatter *gscatter; /* scatter from global space to the subdomain local space */ 16 PCASMType type; /* ASM type */ 17 PetscBool usesdm; /* use the DM for setting up the subproblems */ 18 PetscBool finaljacobian; /* compute the jacobian of the converged solution */ 19 PetscReal damping; /* damping parameter for updates from the blocks */ 20 PetscBool weight_set; /* use a weight in the overlap updates */ 21 22 /* logging events */ 23 PetscLogEvent eventrestrictinterp; 24 PetscLogEvent eventsubsolve; 25 26 PetscInt fjtype; /* type of computed jacobian */ 27 Vec xinit; /* initial solution in case the final jacobian type is computed as first */ 28 } SNES_NASM; 29 30 const char *const SNESNASMTypes[] = {"NONE", "RESTRICT", "INTERPOLATE", "BASIC", "PCASMType", "PC_ASM_", NULL}; 31 const char *const SNESNASMFJTypes[] = {"FINALOUTER", "FINALINNER", "INITIAL"}; 32 33 static PetscErrorCode SNESReset_NASM(SNES snes) 34 { 35 SNES_NASM *nasm = (SNES_NASM *)snes->data; 36 PetscInt i; 37 38 PetscFunctionBegin; 39 for (i = 0; i < nasm->n; i++) { 40 if (nasm->xl) PetscCall(VecDestroy(&nasm->xl[i])); 41 if (nasm->x) PetscCall(VecDestroy(&nasm->x[i])); 42 if (nasm->y) PetscCall(VecDestroy(&nasm->y[i])); 43 if (nasm->b) PetscCall(VecDestroy(&nasm->b[i])); 44 45 if (nasm->subsnes) PetscCall(SNESDestroy(&nasm->subsnes[i])); 46 if (nasm->oscatter) PetscCall(VecScatterDestroy(&nasm->oscatter[i])); 47 if (nasm->oscatter_copy) PetscCall(VecScatterDestroy(&nasm->oscatter_copy[i])); 48 if (nasm->iscatter) PetscCall(VecScatterDestroy(&nasm->iscatter[i])); 49 if (nasm->gscatter) PetscCall(VecScatterDestroy(&nasm->gscatter[i])); 50 } 51 52 PetscCall(PetscFree(nasm->x)); 53 PetscCall(PetscFree(nasm->xl)); 54 PetscCall(PetscFree(nasm->y)); 55 PetscCall(PetscFree(nasm->b)); 56 57 if (nasm->xinit) PetscCall(VecDestroy(&nasm->xinit)); 58 59 PetscCall(PetscFree(nasm->subsnes)); 60 PetscCall(PetscFree(nasm->oscatter)); 61 PetscCall(PetscFree(nasm->oscatter_copy)); 62 PetscCall(PetscFree(nasm->iscatter)); 63 PetscCall(PetscFree(nasm->gscatter)); 64 65 if (nasm->weight_set) PetscCall(VecDestroy(&nasm->weight)); 66 67 nasm->eventrestrictinterp = 0; 68 nasm->eventsubsolve = 0; 69 PetscFunctionReturn(PETSC_SUCCESS); 70 } 71 72 static PetscErrorCode SNESDestroy_NASM(SNES snes) 73 { 74 PetscFunctionBegin; 75 PetscCall(SNESReset_NASM(snes)); 76 PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNASMSetType_C", NULL)); 77 PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNASMGetType_C", NULL)); 78 PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNASMSetSubdomains_C", NULL)); 79 PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNASMGetSubdomains_C", NULL)); 80 PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNASMSetDamping_C", NULL)); 81 PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNASMGetDamping_C", NULL)); 82 PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNASMGetSubdomainVecs_C", NULL)); 83 PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNASMSetComputeFinalJacobian_C", NULL)); 84 PetscCall(PetscFree(snes->data)); 85 PetscFunctionReturn(PETSC_SUCCESS); 86 } 87 88 static PetscErrorCode DMGlobalToLocalSubDomainDirichletHook_Private(DM dm, Vec g, InsertMode mode, Vec l, void *ctx) 89 { 90 Vec bcs = (Vec)ctx; 91 92 PetscFunctionBegin; 93 PetscCall(VecCopy(bcs, l)); 94 PetscFunctionReturn(PETSC_SUCCESS); 95 } 96 97 static PetscErrorCode SNESSetUp_NASM(SNES snes) 98 { 99 SNES_NASM *nasm = (SNES_NASM *)snes->data; 100 DM dm, subdm; 101 DM *subdms; 102 PetscInt i; 103 const char *optionsprefix; 104 Vec F; 105 106 PetscFunctionBegin; 107 if (!nasm->subsnes) { 108 PetscCall(SNESGetDM(snes, &dm)); 109 if (dm) { 110 nasm->usesdm = PETSC_TRUE; 111 PetscCall(DMCreateDomainDecomposition(dm, &nasm->n, NULL, NULL, NULL, &subdms)); 112 PetscCheck(subdms, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "DM has no default decomposition defined. Set subsolves manually with SNESNASMSetSubdomains()."); 113 PetscCall(DMCreateDomainDecompositionScatters(dm, nasm->n, subdms, &nasm->iscatter, &nasm->oscatter, &nasm->gscatter)); 114 PetscCall(PetscMalloc1(nasm->n, &nasm->oscatter_copy)); 115 for (i = 0; i < nasm->n; i++) PetscCall(VecScatterCopy(nasm->oscatter[i], &nasm->oscatter_copy[i])); 116 117 PetscCall(SNESGetOptionsPrefix(snes, &optionsprefix)); 118 PetscCall(PetscMalloc1(nasm->n, &nasm->subsnes)); 119 for (i = 0; i < nasm->n; i++) { 120 PetscCall(SNESCreate(PetscObjectComm((PetscObject)subdms[i]), &nasm->subsnes[i])); 121 PetscCall(PetscObjectIncrementTabLevel((PetscObject)nasm->subsnes[i], (PetscObject)snes, 1)); 122 PetscCall(SNESAppendOptionsPrefix(nasm->subsnes[i], optionsprefix)); 123 PetscCall(SNESAppendOptionsPrefix(nasm->subsnes[i], "sub_")); 124 PetscCall(SNESSetDM(nasm->subsnes[i], subdms[i])); 125 if (snes->ops->usercompute) { 126 PetscCall(SNESSetComputeApplicationContext(nasm->subsnes[i], snes->ops->usercompute, snes->ops->userdestroy)); 127 } else { 128 void *ctx; 129 130 PetscCall(SNESGetApplicationContext(snes, &ctx)); 131 PetscCall(SNESSetApplicationContext(nasm->subsnes[i], ctx)); 132 } 133 PetscCall(SNESSetFromOptions(nasm->subsnes[i])); 134 PetscCall(DMDestroy(&subdms[i])); 135 } 136 PetscCall(PetscFree(subdms)); 137 } else SETERRQ(PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_WRONGSTATE, "Cannot construct local problems automatically without a DM!"); 138 } else SETERRQ(PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_WRONGSTATE, "Must set subproblems manually if there is no DM!"); 139 /* allocate the global vectors */ 140 if (!nasm->x) PetscCall(PetscCalloc1(nasm->n, &nasm->x)); 141 if (!nasm->xl) PetscCall(PetscCalloc1(nasm->n, &nasm->xl)); 142 if (!nasm->y) PetscCall(PetscCalloc1(nasm->n, &nasm->y)); 143 if (!nasm->b) PetscCall(PetscCalloc1(nasm->n, &nasm->b)); 144 145 for (i = 0; i < nasm->n; i++) { 146 PetscCall(SNESGetFunction(nasm->subsnes[i], &F, NULL, NULL)); 147 if (!nasm->x[i]) PetscCall(VecDuplicate(F, &nasm->x[i])); 148 if (!nasm->y[i]) PetscCall(VecDuplicate(F, &nasm->y[i])); 149 if (!nasm->b[i]) PetscCall(VecDuplicate(F, &nasm->b[i])); 150 if (!nasm->xl[i]) { 151 PetscCall(SNESGetDM(nasm->subsnes[i], &subdm)); 152 PetscCall(DMCreateLocalVector(subdm, &nasm->xl[i])); 153 PetscCall(DMGlobalToLocalHookAdd(subdm, DMGlobalToLocalSubDomainDirichletHook_Private, NULL, nasm->xl[i])); 154 } 155 } 156 if (nasm->finaljacobian) { 157 PetscCall(SNESSetUpMatrices(snes)); 158 if (nasm->fjtype == 2) PetscCall(VecDuplicate(snes->vec_sol, &nasm->xinit)); 159 for (i = 0; i < nasm->n; i++) PetscCall(SNESSetUpMatrices(nasm->subsnes[i])); 160 } 161 PetscFunctionReturn(PETSC_SUCCESS); 162 } 163 164 static PetscErrorCode SNESSetFromOptions_NASM(SNES snes, PetscOptionItems *PetscOptionsObject) 165 { 166 PCASMType asmtype; 167 PetscBool flg, monflg; 168 SNES_NASM *nasm = (SNES_NASM *)snes->data; 169 170 PetscFunctionBegin; 171 PetscOptionsHeadBegin(PetscOptionsObject, "Nonlinear Additive Schwarz options"); 172 PetscCall(PetscOptionsEnum("-snes_nasm_type", "Type of restriction/extension", "", SNESNASMTypes, (PetscEnum)nasm->type, (PetscEnum *)&asmtype, &flg)); 173 if (flg) PetscCall(SNESNASMSetType(snes, asmtype)); 174 flg = PETSC_FALSE; 175 monflg = PETSC_TRUE; 176 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)); 177 if (flg) PetscCall(SNESNASMSetDamping(snes, nasm->damping)); 178 PetscCall(PetscOptionsDeprecated("-snes_nasm_sub_view", NULL, "3.15", "Use -snes_view ::ascii_info_detail")); 179 PetscCall(PetscOptionsBool("-snes_nasm_finaljacobian", "Compute the global jacobian of the final iterate (for ASPIN)", "", nasm->finaljacobian, &nasm->finaljacobian, NULL)); 180 PetscCall(PetscOptionsEList("-snes_nasm_finaljacobian_type", "The type of the final jacobian computed.", "", SNESNASMFJTypes, 3, SNESNASMFJTypes[0], &nasm->fjtype, NULL)); 181 PetscCall(PetscOptionsBool("-snes_nasm_log", "Log times for subSNES solves and restriction", "", monflg, &monflg, &flg)); 182 if (flg) { 183 PetscCall(PetscLogEventRegister("SNESNASMSubSolve", ((PetscObject)snes)->classid, &nasm->eventsubsolve)); 184 PetscCall(PetscLogEventRegister("SNESNASMRestrict", ((PetscObject)snes)->classid, &nasm->eventrestrictinterp)); 185 } 186 PetscOptionsHeadEnd(); 187 PetscFunctionReturn(PETSC_SUCCESS); 188 } 189 190 static PetscErrorCode SNESView_NASM(SNES snes, PetscViewer viewer) 191 { 192 SNES_NASM *nasm = (SNES_NASM *)snes->data; 193 PetscMPIInt rank, size; 194 PetscInt i, N, bsz; 195 PetscBool iascii, isstring; 196 PetscViewer sviewer; 197 MPI_Comm comm; 198 PetscViewerFormat format; 199 const char *prefix; 200 201 PetscFunctionBegin; 202 PetscCall(PetscObjectGetComm((PetscObject)snes, &comm)); 203 PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii)); 204 PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERSTRING, &isstring)); 205 PetscCallMPI(MPI_Comm_rank(comm, &rank)); 206 PetscCallMPI(MPI_Comm_size(comm, &size)); 207 PetscCall(MPIU_Allreduce(&nasm->n, &N, 1, MPIU_INT, MPI_SUM, comm)); 208 if (iascii) { 209 PetscCall(PetscViewerASCIIPrintf(viewer, " total subdomain blocks = %" PetscInt_FMT "\n", N)); 210 PetscCall(PetscViewerGetFormat(viewer, &format)); 211 if (format != PETSC_VIEWER_ASCII_INFO_DETAIL) { 212 if (nasm->subsnes) { 213 PetscCall(PetscViewerASCIIPrintf(viewer, " Local solver information for first block on rank 0:\n")); 214 PetscCall(SNESGetOptionsPrefix(snes, &prefix)); 215 PetscCall(PetscViewerASCIIPrintf(viewer, " Use -%ssnes_view ::ascii_info_detail to display information for all blocks\n", prefix ? prefix : "")); 216 PetscCall(PetscViewerASCIIPushTab(viewer)); 217 PetscCall(PetscViewerGetSubViewer(viewer, PETSC_COMM_SELF, &sviewer)); 218 if (rank == 0) { 219 PetscCall(PetscViewerASCIIPushTab(sviewer)); 220 PetscCall(SNESView(nasm->subsnes[0], sviewer)); 221 PetscCall(PetscViewerASCIIPopTab(sviewer)); 222 } 223 PetscCall(PetscViewerRestoreSubViewer(viewer, PETSC_COMM_SELF, &sviewer)); 224 PetscCall(PetscViewerASCIIPopTab(viewer)); 225 } 226 } else { 227 /* print the solver on each block */ 228 PetscCall(PetscViewerASCIIPushSynchronized(viewer)); 229 PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, " [%d] number of local blocks = %" PetscInt_FMT "\n", (int)rank, nasm->n)); 230 PetscCall(PetscViewerFlush(viewer)); 231 PetscCall(PetscViewerASCIIPopSynchronized(viewer)); 232 PetscCall(PetscViewerASCIIPrintf(viewer, " Local solver information for each block is in the following SNES objects:\n")); 233 PetscCall(PetscViewerASCIIPushTab(viewer)); 234 PetscCall(PetscViewerASCIIPrintf(viewer, "- - - - - - - - - - - - - - - - - -\n")); 235 PetscCall(PetscViewerGetSubViewer(viewer, PETSC_COMM_SELF, &sviewer)); 236 for (i = 0; i < nasm->n; i++) { 237 PetscCall(VecGetLocalSize(nasm->x[i], &bsz)); 238 PetscCall(PetscViewerASCIIPrintf(sviewer, "[%d] local block number %" PetscInt_FMT ", size = %" PetscInt_FMT "\n", (int)rank, i, bsz)); 239 PetscCall(SNESView(nasm->subsnes[i], sviewer)); 240 PetscCall(PetscViewerASCIIPrintf(sviewer, "- - - - - - - - - - - - - - - - - -\n")); 241 } 242 PetscCall(PetscViewerRestoreSubViewer(viewer, PETSC_COMM_SELF, &sviewer)); 243 PetscCall(PetscViewerASCIIPopTab(viewer)); 244 } 245 } else if (isstring) { 246 PetscCall(PetscViewerStringSPrintf(viewer, " blocks=%" PetscInt_FMT ",type=%s", N, SNESNASMTypes[nasm->type])); 247 PetscCall(PetscViewerGetSubViewer(viewer, PETSC_COMM_SELF, &sviewer)); 248 if (nasm->subsnes && rank == 0) PetscCall(SNESView(nasm->subsnes[0], sviewer)); 249 PetscCall(PetscViewerRestoreSubViewer(viewer, PETSC_COMM_SELF, &sviewer)); 250 } 251 PetscFunctionReturn(PETSC_SUCCESS); 252 } 253 254 /*@ 255 SNESNASMSetType - Set the type of subdomain update used for the nonlinear additive Schwarz solver `SNESNASM` 256 257 Logically Collective 258 259 Input Parameters: 260 + snes - the `SNES` context 261 - type - the type of update, `PC_ASM_BASIC` or `PC_ASM_RESTRICT` 262 263 Options Database Key: 264 . -snes_nasm_type <basic,restrict> - type of subdomain update used 265 266 Level: intermediate 267 268 .seealso: [](ch_snes), `SNES`, `SNESNASM`, `SNESNASMGetType()`, `PCASMSetType()`, `PC_ASM_BASIC`, `PC_ASM_RESTRICT`, `PCASMType` 269 @*/ 270 PetscErrorCode SNESNASMSetType(SNES snes, PCASMType type) 271 { 272 PetscErrorCode (*f)(SNES, PCASMType); 273 274 PetscFunctionBegin; 275 PetscCall(PetscObjectQueryFunction((PetscObject)snes, "SNESNASMSetType_C", &f)); 276 if (f) PetscCall((f)(snes, type)); 277 PetscFunctionReturn(PETSC_SUCCESS); 278 } 279 280 static PetscErrorCode SNESNASMSetType_NASM(SNES snes, PCASMType type) 281 { 282 SNES_NASM *nasm = (SNES_NASM *)snes->data; 283 284 PetscFunctionBegin; 285 PetscCheck(type == PC_ASM_BASIC || type == PC_ASM_RESTRICT, PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_OUTOFRANGE, "SNESNASM only supports basic and restrict types"); 286 nasm->type = type; 287 PetscFunctionReturn(PETSC_SUCCESS); 288 } 289 290 /*@ 291 SNESNASMGetType - Get the type of subdomain update used for the nonlinear additive Schwarz solver `SNESNASM` 292 293 Logically Collective 294 295 Input Parameter: 296 . snes - the `SNES` context 297 298 Output Parameter: 299 . type - the type of update 300 301 Level: intermediate 302 303 .seealso: [](ch_snes), `SNES`, `SNESNASM`, `SNESNASMSetType()`, `PCASMGetType()`, `PC_ASM_BASIC`, `PC_ASM_RESTRICT`, `PCASMType` 304 @*/ 305 PetscErrorCode SNESNASMGetType(SNES snes, PCASMType *type) 306 { 307 PetscFunctionBegin; 308 PetscUseMethod(snes, "SNESNASMGetType_C", (SNES, PCASMType *), (snes, type)); 309 PetscFunctionReturn(PETSC_SUCCESS); 310 } 311 312 static PetscErrorCode SNESNASMGetType_NASM(SNES snes, PCASMType *type) 313 { 314 SNES_NASM *nasm = (SNES_NASM *)snes->data; 315 316 PetscFunctionBegin; 317 *type = nasm->type; 318 PetscFunctionReturn(PETSC_SUCCESS); 319 } 320 321 /*@ 322 SNESNASMSetSubdomains - Manually Set the context required to restrict and solve subdomain problems in the nonlinear additive Schwarz solver 323 324 Logically Collective 325 326 Input Parameters: 327 + snes - the `SNES` context 328 . n - the number of local subdomains 329 . subsnes - solvers defined on the local subdomains 330 . iscatter - scatters into the nonoverlapping portions of the local subdomains 331 . oscatter - scatters into the overlapping portions of the local subdomains 332 - gscatter - scatters into the (ghosted) local vector of the local subdomain 333 334 Level: intermediate 335 336 .seealso: [](ch_snes), `SNES`, `SNESNASM`, `SNESNASMGetSubdomains()` 337 @*/ 338 PetscErrorCode SNESNASMSetSubdomains(SNES snes, PetscInt n, SNES subsnes[], VecScatter iscatter[], VecScatter oscatter[], VecScatter gscatter[]) 339 { 340 PetscErrorCode (*f)(SNES, PetscInt, SNES *, VecScatter *, VecScatter *, VecScatter *); 341 342 PetscFunctionBegin; 343 PetscCall(PetscObjectQueryFunction((PetscObject)snes, "SNESNASMSetSubdomains_C", &f)); 344 if (f) PetscCall((f)(snes, n, subsnes, iscatter, oscatter, gscatter)); 345 PetscFunctionReturn(PETSC_SUCCESS); 346 } 347 348 static PetscErrorCode SNESNASMSetSubdomains_NASM(SNES snes, PetscInt n, SNES subsnes[], VecScatter iscatter[], VecScatter oscatter[], VecScatter gscatter[]) 349 { 350 PetscInt i; 351 SNES_NASM *nasm = (SNES_NASM *)snes->data; 352 353 PetscFunctionBegin; 354 PetscCheck(!snes->setupcalled, PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_WRONGSTATE, "SNESNASMSetSubdomains() should be called before calling SNESSetUp()."); 355 356 /* tear down the previously set things */ 357 PetscCall(SNESReset(snes)); 358 359 nasm->n = n; 360 if (oscatter) { 361 for (i = 0; i < n; i++) PetscCall(PetscObjectReference((PetscObject)oscatter[i])); 362 } 363 if (iscatter) { 364 for (i = 0; i < n; i++) PetscCall(PetscObjectReference((PetscObject)iscatter[i])); 365 } 366 if (gscatter) { 367 for (i = 0; i < n; i++) PetscCall(PetscObjectReference((PetscObject)gscatter[i])); 368 } 369 if (oscatter) { 370 PetscCall(PetscMalloc1(n, &nasm->oscatter)); 371 PetscCall(PetscMalloc1(n, &nasm->oscatter_copy)); 372 for (i = 0; i < n; i++) { 373 nasm->oscatter[i] = oscatter[i]; 374 PetscCall(VecScatterCopy(oscatter[i], &nasm->oscatter_copy[i])); 375 } 376 } 377 if (iscatter) { 378 PetscCall(PetscMalloc1(n, &nasm->iscatter)); 379 for (i = 0; i < n; i++) nasm->iscatter[i] = iscatter[i]; 380 } 381 if (gscatter) { 382 PetscCall(PetscMalloc1(n, &nasm->gscatter)); 383 for (i = 0; i < n; i++) nasm->gscatter[i] = gscatter[i]; 384 } 385 386 if (subsnes) { 387 PetscCall(PetscMalloc1(n, &nasm->subsnes)); 388 for (i = 0; i < n; i++) nasm->subsnes[i] = subsnes[i]; 389 } 390 PetscFunctionReturn(PETSC_SUCCESS); 391 } 392 393 /*@ 394 SNESNASMGetSubdomains - Get the local subdomain contexts for the nonlinear additive Schwarz solver 395 396 Not Collective but some of the objects returned will be parallel 397 398 Input Parameter: 399 . snes - the `SNES` context 400 401 Output Parameters: 402 + n - the number of local subdomains 403 . subsnes - solvers defined on the local subdomains 404 . iscatter - scatters into the nonoverlapping portions of the local subdomains 405 . oscatter - scatters into the overlapping portions of the local subdomains 406 - gscatter - scatters into the (ghosted) local vector of the local subdomain 407 408 Level: intermediate 409 410 .seealso: [](ch_snes), `SNES`, `SNESNASM`, `SNESNASMSetSubdomains()` 411 @*/ 412 PetscErrorCode SNESNASMGetSubdomains(SNES snes, PetscInt *n, SNES *subsnes[], VecScatter *iscatter[], VecScatter *oscatter[], VecScatter *gscatter[]) 413 { 414 PetscErrorCode (*f)(SNES, PetscInt *, SNES **, VecScatter **, VecScatter **, VecScatter **); 415 416 PetscFunctionBegin; 417 PetscCall(PetscObjectQueryFunction((PetscObject)snes, "SNESNASMGetSubdomains_C", &f)); 418 if (f) PetscCall((f)(snes, n, subsnes, iscatter, oscatter, gscatter)); 419 PetscFunctionReturn(PETSC_SUCCESS); 420 } 421 422 static PetscErrorCode SNESNASMGetSubdomains_NASM(SNES snes, PetscInt *n, SNES *subsnes[], VecScatter *iscatter[], VecScatter *oscatter[], VecScatter *gscatter[]) 423 { 424 SNES_NASM *nasm = (SNES_NASM *)snes->data; 425 426 PetscFunctionBegin; 427 if (n) *n = nasm->n; 428 if (oscatter) *oscatter = nasm->oscatter; 429 if (iscatter) *iscatter = nasm->iscatter; 430 if (gscatter) *gscatter = nasm->gscatter; 431 if (subsnes) *subsnes = nasm->subsnes; 432 PetscFunctionReturn(PETSC_SUCCESS); 433 } 434 435 /*@ 436 SNESNASMGetSubdomainVecs - Get the processor-local subdomain vectors for the nonlinear additive Schwarz solver 437 438 Not Collective 439 440 Input Parameter: 441 . snes - the `SNES` context 442 443 Output Parameters: 444 + n - the number of local subdomains 445 . x - The subdomain solution vector 446 . y - The subdomain step vector 447 . b - The subdomain RHS vector 448 - xl - The subdomain local vectors (ghosted) 449 450 Level: developer 451 452 .seealso: [](ch_snes), `SNES`, `SNESNASM`, `SNESNASMGetSubdomains()` 453 @*/ 454 PetscErrorCode SNESNASMGetSubdomainVecs(SNES snes, PetscInt *n, Vec **x, Vec **y, Vec **b, Vec **xl) 455 { 456 PetscErrorCode (*f)(SNES, PetscInt *, Vec **, Vec **, Vec **, Vec **); 457 458 PetscFunctionBegin; 459 PetscCall(PetscObjectQueryFunction((PetscObject)snes, "SNESNASMGetSubdomainVecs_C", &f)); 460 if (f) PetscCall((f)(snes, n, x, y, b, xl)); 461 PetscFunctionReturn(PETSC_SUCCESS); 462 } 463 464 static PetscErrorCode SNESNASMGetSubdomainVecs_NASM(SNES snes, PetscInt *n, Vec **x, Vec **y, Vec **b, Vec **xl) 465 { 466 SNES_NASM *nasm = (SNES_NASM *)snes->data; 467 468 PetscFunctionBegin; 469 if (n) *n = nasm->n; 470 if (x) *x = nasm->x; 471 if (y) *y = nasm->y; 472 if (b) *b = nasm->b; 473 if (xl) *xl = nasm->xl; 474 PetscFunctionReturn(PETSC_SUCCESS); 475 } 476 477 /*@ 478 SNESNASMSetComputeFinalJacobian - Schedules the computation of the global and subdomain Jacobians upon convergence for the 479 nonlinear additive Schwarz solver 480 481 Collective 482 483 Input Parameters: 484 + snes - the SNES context 485 - flg - `PETSC_TRUE` to compute the Jacobians 486 487 Level: developer 488 489 Notes: 490 This is used almost exclusively in the implementation of `SNESASPIN`, where the converged subdomain and global Jacobian 491 is needed at each linear iteration. 492 493 .seealso: [](ch_snes), `SNES`, `SNESNASM`, `SNESNASMGetSubdomains()` 494 @*/ 495 PetscErrorCode SNESNASMSetComputeFinalJacobian(SNES snes, PetscBool flg) 496 { 497 PetscErrorCode (*f)(SNES, PetscBool); 498 499 PetscFunctionBegin; 500 PetscCall(PetscObjectQueryFunction((PetscObject)snes, "SNESNASMSetComputeFinalJacobian_C", &f)); 501 if (f) PetscCall((f)(snes, flg)); 502 PetscFunctionReturn(PETSC_SUCCESS); 503 } 504 505 static PetscErrorCode SNESNASMSetComputeFinalJacobian_NASM(SNES snes, PetscBool flg) 506 { 507 SNES_NASM *nasm = (SNES_NASM *)snes->data; 508 509 PetscFunctionBegin; 510 nasm->finaljacobian = flg; 511 PetscFunctionReturn(PETSC_SUCCESS); 512 } 513 514 /*@ 515 SNESNASMSetDamping - Sets the update damping for `SNESNASM` the nonlinear additive Schwarz solver 516 517 Logically Collective 518 519 Input Parameters: 520 + snes - the `SNES` context 521 - dmp - damping 522 523 Options Database Key: 524 . -snes_nasm_damping <dmp> - the new solution is obtained as old solution plus `dmp` times (sum of the solutions on the subdomains) 525 526 Level: intermediate 527 528 Note: 529 The new solution is obtained as old solution plus dmp times (sum of the solutions on the subdomains) 530 531 .seealso: [](ch_snes), `SNES`, `SNESNASM`, `SNESNASMGetDamping()` 532 @*/ 533 PetscErrorCode SNESNASMSetDamping(SNES snes, PetscReal dmp) 534 { 535 PetscErrorCode (*f)(SNES, PetscReal); 536 537 PetscFunctionBegin; 538 PetscCall(PetscObjectQueryFunction((PetscObject)snes, "SNESNASMSetDamping_C", (void (**)(void)) & f)); 539 if (f) PetscCall((f)(snes, dmp)); 540 PetscFunctionReturn(PETSC_SUCCESS); 541 } 542 543 static PetscErrorCode SNESNASMSetDamping_NASM(SNES snes, PetscReal dmp) 544 { 545 SNES_NASM *nasm = (SNES_NASM *)snes->data; 546 547 PetscFunctionBegin; 548 nasm->damping = dmp; 549 PetscFunctionReturn(PETSC_SUCCESS); 550 } 551 552 /*@ 553 SNESNASMGetDamping - Gets the update damping for `SNESNASM` the nonlinear additive Schwarz solver 554 555 Not Collective 556 557 Input Parameter: 558 . snes - the `SNES` context 559 560 Output Parameter: 561 . dmp - damping 562 563 Level: intermediate 564 565 .seealso: [](ch_snes), `SNES`, `SNESNASM`, `SNESNASMSetDamping()` 566 @*/ 567 PetscErrorCode SNESNASMGetDamping(SNES snes, PetscReal *dmp) 568 { 569 PetscFunctionBegin; 570 PetscUseMethod(snes, "SNESNASMGetDamping_C", (SNES, PetscReal *), (snes, dmp)); 571 PetscFunctionReturn(PETSC_SUCCESS); 572 } 573 574 static PetscErrorCode SNESNASMGetDamping_NASM(SNES snes, PetscReal *dmp) 575 { 576 SNES_NASM *nasm = (SNES_NASM *)snes->data; 577 578 PetscFunctionBegin; 579 *dmp = nasm->damping; 580 PetscFunctionReturn(PETSC_SUCCESS); 581 } 582 583 /* 584 Input Parameters: 585 + snes - The solver 586 . B - The RHS vector 587 - X - The initial guess 588 589 Output Parameter: 590 . Y - The solution update 591 592 TODO: All scatters should be packed into one 593 */ 594 static PetscErrorCode SNESNASMSolveLocal_Private(SNES snes, Vec B, Vec Y, Vec X) 595 { 596 SNES_NASM *nasm = (SNES_NASM *)snes->data; 597 SNES subsnes; 598 PetscInt i; 599 PetscReal dmp; 600 Vec Xl, Bl, Yl, Xlloc; 601 VecScatter iscat, oscat, gscat, oscat_copy; 602 DM dm, subdm; 603 PCASMType type; 604 605 PetscFunctionBegin; 606 PetscCall(SNESNASMGetType(snes, &type)); 607 PetscCall(SNESGetDM(snes, &dm)); 608 PetscCall(VecSet(Y, 0)); 609 if (nasm->eventrestrictinterp) PetscCall(PetscLogEventBegin(nasm->eventrestrictinterp, snes, 0, 0, 0)); 610 for (i = 0; i < nasm->n; i++) { /* scatter the global solution to the overlap solution and the local solution */ 611 Xl = nasm->x[i]; 612 Xlloc = nasm->xl[i]; 613 oscat = nasm->oscatter[i]; 614 oscat_copy = nasm->oscatter_copy[i]; 615 gscat = nasm->gscatter[i]; 616 PetscCall(VecScatterBegin(oscat, X, Xl, INSERT_VALUES, SCATTER_FORWARD)); 617 PetscCall(VecScatterBegin(gscat, X, Xlloc, INSERT_VALUES, SCATTER_FORWARD)); 618 619 if (B) { 620 /* scatter the RHS to the local RHS */ 621 Bl = nasm->b[i]; 622 PetscCall(VecScatterBegin(oscat_copy, B, Bl, INSERT_VALUES, SCATTER_FORWARD)); 623 } 624 } 625 if (nasm->eventrestrictinterp) PetscCall(PetscLogEventEnd(nasm->eventrestrictinterp, snes, 0, 0, 0)); 626 627 if (nasm->eventsubsolve) PetscCall(PetscLogEventBegin(nasm->eventsubsolve, snes, 0, 0, 0)); 628 for (i = 0; i < nasm->n; i++) { 629 PetscErrorCode (*bl)(DM, Vec, void *); 630 void *bctx; 631 632 Xl = nasm->x[i]; 633 Xlloc = nasm->xl[i]; 634 Yl = nasm->y[i]; 635 subsnes = nasm->subsnes[i]; 636 PetscCall(SNESGetDM(subsnes, &subdm)); 637 iscat = nasm->iscatter[i]; 638 oscat = nasm->oscatter[i]; 639 oscat_copy = nasm->oscatter_copy[i]; 640 gscat = nasm->gscatter[i]; 641 PetscCall(VecScatterEnd(oscat, X, Xl, INSERT_VALUES, SCATTER_FORWARD)); 642 PetscCall(VecScatterEnd(gscat, X, Xlloc, INSERT_VALUES, SCATTER_FORWARD)); 643 if (B) { 644 Bl = nasm->b[i]; 645 PetscCall(VecScatterEnd(oscat_copy, B, Bl, INSERT_VALUES, SCATTER_FORWARD)); 646 } else Bl = NULL; 647 648 PetscCall(SNESGetDM(subsnes, &subdm)); 649 PetscCall(DMSNESGetBoundaryLocal(subdm, &bl, &bctx)); 650 if (bl) PetscCall((*bl)(subdm, Xlloc, bctx)); 651 652 PetscCall(DMSubDomainRestrict(dm, oscat, gscat, subdm)); 653 PetscCall(VecCopy(Xl, Yl)); 654 PetscCall(SNESSolve(subsnes, Bl, Xl)); 655 PetscCall(VecAYPX(Yl, -1.0, Xl)); 656 PetscCall(VecScale(Yl, nasm->damping)); 657 if (type == PC_ASM_BASIC) { 658 PetscCall(VecScatterBegin(oscat, Yl, Y, ADD_VALUES, SCATTER_REVERSE)); 659 PetscCall(VecScatterEnd(oscat, Yl, Y, ADD_VALUES, SCATTER_REVERSE)); 660 } else if (type == PC_ASM_RESTRICT) { 661 PetscCall(VecScatterBegin(iscat, Yl, Y, ADD_VALUES, SCATTER_REVERSE)); 662 PetscCall(VecScatterEnd(iscat, Yl, Y, ADD_VALUES, SCATTER_REVERSE)); 663 } else SETERRQ(PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_WRONGSTATE, "Only basic and restrict types are supported for SNESNASM"); 664 } 665 if (nasm->eventsubsolve) PetscCall(PetscLogEventEnd(nasm->eventsubsolve, snes, 0, 0, 0)); 666 if (nasm->eventrestrictinterp) PetscCall(PetscLogEventBegin(nasm->eventrestrictinterp, snes, 0, 0, 0)); 667 if (nasm->weight_set) PetscCall(VecPointwiseMult(Y, Y, nasm->weight)); 668 if (nasm->eventrestrictinterp) PetscCall(PetscLogEventEnd(nasm->eventrestrictinterp, snes, 0, 0, 0)); 669 PetscCall(SNESNASMGetDamping(snes, &dmp)); 670 PetscCall(VecAXPY(X, dmp, Y)); 671 PetscFunctionReturn(PETSC_SUCCESS); 672 } 673 674 static PetscErrorCode SNESNASMComputeFinalJacobian_Private(SNES snes, Vec Xfinal) 675 { 676 Vec X = Xfinal; 677 SNES_NASM *nasm = (SNES_NASM *)snes->data; 678 SNES subsnes; 679 PetscInt i, lag = 1; 680 Vec Xlloc, Xl, Fl, F; 681 VecScatter oscat, gscat; 682 DM dm, subdm; 683 684 PetscFunctionBegin; 685 if (nasm->fjtype == 2) X = nasm->xinit; 686 F = snes->vec_func; 687 if (snes->normschedule == SNES_NORM_NONE) PetscCall(SNESComputeFunction(snes, X, F)); 688 PetscCall(SNESComputeJacobian(snes, X, snes->jacobian, snes->jacobian_pre)); 689 PetscCall(SNESGetDM(snes, &dm)); 690 if (nasm->eventrestrictinterp) PetscCall(PetscLogEventBegin(nasm->eventrestrictinterp, snes, 0, 0, 0)); 691 if (nasm->fjtype != 1) { 692 for (i = 0; i < nasm->n; i++) { 693 Xlloc = nasm->xl[i]; 694 gscat = nasm->gscatter[i]; 695 PetscCall(VecScatterBegin(gscat, X, Xlloc, INSERT_VALUES, SCATTER_FORWARD)); 696 } 697 } 698 if (nasm->eventrestrictinterp) PetscCall(PetscLogEventEnd(nasm->eventrestrictinterp, snes, 0, 0, 0)); 699 for (i = 0; i < nasm->n; i++) { 700 Fl = nasm->subsnes[i]->vec_func; 701 Xl = nasm->x[i]; 702 Xlloc = nasm->xl[i]; 703 subsnes = nasm->subsnes[i]; 704 oscat = nasm->oscatter[i]; 705 gscat = nasm->gscatter[i]; 706 if (nasm->fjtype != 1) PetscCall(VecScatterEnd(gscat, X, Xlloc, INSERT_VALUES, SCATTER_FORWARD)); 707 PetscCall(SNESGetDM(subsnes, &subdm)); 708 PetscCall(DMSubDomainRestrict(dm, oscat, gscat, subdm)); 709 if (nasm->fjtype != 1) { 710 PetscCall(DMLocalToGlobalBegin(subdm, Xlloc, INSERT_VALUES, Xl)); 711 PetscCall(DMLocalToGlobalEnd(subdm, Xlloc, INSERT_VALUES, Xl)); 712 } 713 if (subsnes->lagjacobian == -1) subsnes->lagjacobian = -2; 714 else if (subsnes->lagjacobian > 1) lag = subsnes->lagjacobian; 715 PetscCall(SNESComputeFunction(subsnes, Xl, Fl)); 716 PetscCall(SNESComputeJacobian(subsnes, Xl, subsnes->jacobian, subsnes->jacobian_pre)); 717 if (lag > 1) subsnes->lagjacobian = lag; 718 } 719 PetscFunctionReturn(PETSC_SUCCESS); 720 } 721 722 static PetscErrorCode SNESSolve_NASM(SNES snes) 723 { 724 Vec F; 725 Vec X; 726 Vec B; 727 Vec Y; 728 PetscInt i; 729 PetscReal fnorm = 0.0; 730 SNESNormSchedule normschedule; 731 SNES_NASM *nasm = (SNES_NASM *)snes->data; 732 733 PetscFunctionBegin; 734 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); 735 736 PetscCall(PetscCitationsRegister(SNESCitation, &SNEScite)); 737 X = snes->vec_sol; 738 Y = snes->vec_sol_update; 739 F = snes->vec_func; 740 B = snes->vec_rhs; 741 742 PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes)); 743 snes->iter = 0; 744 snes->norm = 0.; 745 PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 746 snes->reason = SNES_CONVERGED_ITERATING; 747 PetscCall(SNESGetNormSchedule(snes, &normschedule)); 748 if (normschedule == SNES_NORM_ALWAYS || normschedule == SNES_NORM_INITIAL_ONLY || normschedule == SNES_NORM_INITIAL_FINAL_ONLY || !snes->max_its) { 749 /* compute the initial function and preconditioned update delX */ 750 if (!snes->vec_func_init_set) { 751 PetscCall(SNESComputeFunction(snes, X, F)); 752 } else snes->vec_func_init_set = PETSC_FALSE; 753 754 PetscCall(VecNorm(F, NORM_2, &fnorm)); /* fnorm <- ||F|| */ 755 SNESCheckFunctionNorm(snes, fnorm); 756 PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes)); 757 snes->iter = 0; 758 snes->norm = fnorm; 759 PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 760 PetscCall(SNESLogConvergenceHistory(snes, snes->norm, 0)); 761 762 /* test convergence */ 763 PetscCall(SNESConverged(snes, 0, 0.0, 0.0, fnorm)); 764 PetscCall(SNESMonitor(snes, 0, snes->norm)); 765 if (snes->reason) PetscFunctionReturn(PETSC_SUCCESS); 766 } else { 767 PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 768 PetscCall(SNESLogConvergenceHistory(snes, snes->norm, 0)); 769 PetscCall(SNESMonitor(snes, snes->iter, snes->norm)); 770 } 771 772 /* Call general purpose update function */ 773 PetscTryTypeMethod(snes, update, snes->iter); 774 /* copy the initial solution over for later */ 775 if (nasm->fjtype == 2) PetscCall(VecCopy(X, nasm->xinit)); 776 777 for (i = 0; i < snes->max_its; i++) { 778 PetscCall(SNESNASMSolveLocal_Private(snes, B, Y, X)); 779 if (normschedule == SNES_NORM_ALWAYS || ((i == snes->max_its - 1) && (normschedule == SNES_NORM_INITIAL_FINAL_ONLY || normschedule == SNES_NORM_FINAL_ONLY))) { 780 PetscCall(SNESComputeFunction(snes, X, F)); 781 PetscCall(VecNorm(F, NORM_2, &fnorm)); /* fnorm <- ||F|| */ 782 SNESCheckFunctionNorm(snes, fnorm); 783 } 784 /* Monitor convergence */ 785 PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes)); 786 snes->iter = i + 1; 787 snes->norm = fnorm; 788 PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 789 PetscCall(SNESLogConvergenceHistory(snes, snes->norm, 0)); 790 /* Test for convergence */ 791 PetscCall(SNESConverged(snes, snes->iter, 0.0, 0.0, fnorm)); 792 PetscCall(SNESMonitor(snes, snes->iter, snes->norm)); 793 if (snes->reason) break; 794 /* Call general purpose update function */ 795 PetscTryTypeMethod(snes, update, snes->iter); 796 } 797 if (nasm->finaljacobian) { 798 PetscCall(SNESNASMComputeFinalJacobian_Private(snes, X)); 799 SNESCheckJacobianDomainerror(snes); 800 } 801 PetscFunctionReturn(PETSC_SUCCESS); 802 } 803 804 /*MC 805 SNESNASM - Nonlinear Additive Schwarz solver {cite}`ck02`, {cite}`bruneknepleysmithtu15` 806 807 Options Database Keys: 808 + -snes_nasm_log - enable logging events for the communication and solve stages 809 . -snes_nasm_type <basic,restrict> - type of subdomain update used 810 . -snes_nasm_damping <dmp> - the new solution is obtained as old solution plus dmp times (sum of the solutions on the subdomains) 811 . -snes_nasm_finaljacobian - compute the local and global Jacobians of the final iterate 812 . -snes_nasm_finaljacobian_type <finalinner,finalouter,initial> - pick state the Jacobian is calculated at 813 . -sub_snes_ - options prefix of the subdomain nonlinear solves 814 . -sub_ksp_ - options prefix of the subdomain Krylov solver 815 - -sub_pc_ - options prefix of the subdomain preconditioner 816 817 Level: advanced 818 819 Note: 820 This is not often used directly as a solver, it converges too slowly. However it works well as a nonlinear preconditioner for 821 the `SNESASPIN` solver 822 823 Developer Note: 824 This is a non-Newton based nonlinear solver that does not directly require a Jacobian; hence the flag snes->usesksp is set to 825 false and `SNESView()` and -snes_view do not display a `KSP` object. However, if the flag nasm->finaljacobian is set (for example, if 826 `SNESNASM` is used as a nonlinear preconditioner for `SNESASPIN`) then `SNESSetUpMatrices()` is called to generate the 827 Jacobian (needed by `SNESASPIN`) 828 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 829 used by `SNESASPIN` they share the same Jacobian matrices because `SNESSetUp()` (called on the outer `SNESASPIN`) causes the inner `SNES` 830 object (in this case `SNESNASM`) to inherit the outer Jacobian matrices. 831 832 .seealso: [](ch_snes), `SNESCreate()`, `SNES`, `SNESSetType()`, `SNESType`, `SNESNASMSetType()`, `SNESNASMGetType()`, `SNESNASMSetSubdomains()`, `SNESNASMGetSubdomains()`, 833 `SNESNASMGetSubdomainVecs()`, `SNESNASMSetComputeFinalJacobian()`, `SNESNASMSetDamping()`, `SNESNASMGetDamping()`, `SNESNASMSetWeight()`, 834 `SNESNASMGetSNES()`, `SNESNASMGetNumber()` 835 M*/ 836 837 PETSC_EXTERN PetscErrorCode SNESCreate_NASM(SNES snes) 838 { 839 SNES_NASM *nasm; 840 841 PetscFunctionBegin; 842 PetscCall(PetscNew(&nasm)); 843 snes->data = (void *)nasm; 844 845 nasm->n = PETSC_DECIDE; 846 nasm->subsnes = NULL; 847 nasm->x = NULL; 848 nasm->xl = NULL; 849 nasm->y = NULL; 850 nasm->b = NULL; 851 nasm->oscatter = NULL; 852 nasm->oscatter_copy = NULL; 853 nasm->iscatter = NULL; 854 nasm->gscatter = NULL; 855 nasm->damping = 1.; 856 857 nasm->type = PC_ASM_BASIC; 858 nasm->finaljacobian = PETSC_FALSE; 859 nasm->weight_set = PETSC_FALSE; 860 861 snes->ops->destroy = SNESDestroy_NASM; 862 snes->ops->setup = SNESSetUp_NASM; 863 snes->ops->setfromoptions = SNESSetFromOptions_NASM; 864 snes->ops->view = SNESView_NASM; 865 snes->ops->solve = SNESSolve_NASM; 866 snes->ops->reset = SNESReset_NASM; 867 868 snes->usesksp = PETSC_FALSE; 869 snes->usesnpc = PETSC_FALSE; 870 871 snes->alwayscomputesfinalresidual = PETSC_FALSE; 872 873 nasm->fjtype = 0; 874 nasm->xinit = NULL; 875 nasm->eventrestrictinterp = 0; 876 nasm->eventsubsolve = 0; 877 878 if (!snes->tolerancesset) { 879 snes->max_its = 10000; 880 snes->max_funcs = 10000; 881 } 882 883 PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNASMSetType_C", SNESNASMSetType_NASM)); 884 PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNASMGetType_C", SNESNASMGetType_NASM)); 885 PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNASMSetSubdomains_C", SNESNASMSetSubdomains_NASM)); 886 PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNASMGetSubdomains_C", SNESNASMGetSubdomains_NASM)); 887 PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNASMSetDamping_C", SNESNASMSetDamping_NASM)); 888 PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNASMGetDamping_C", SNESNASMGetDamping_NASM)); 889 PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNASMGetSubdomainVecs_C", SNESNASMGetSubdomainVecs_NASM)); 890 PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESNASMSetComputeFinalJacobian_C", SNESNASMSetComputeFinalJacobian_NASM)); 891 PetscFunctionReturn(PETSC_SUCCESS); 892 } 893 894 /*@ 895 SNESNASMGetSNES - Gets a subsolver 896 897 Not Collective 898 899 Input Parameters: 900 + snes - the `SNES` context 901 - i - the number of the subsnes to get 902 903 Output Parameter: 904 . subsnes - the subsolver context 905 906 Level: intermediate 907 908 .seealso: [](ch_snes), `SNESNASM`, `SNESNASMGetNumber()` 909 @*/ 910 PetscErrorCode SNESNASMGetSNES(SNES snes, PetscInt i, SNES *subsnes) 911 { 912 SNES_NASM *nasm = (SNES_NASM *)snes->data; 913 914 PetscFunctionBegin; 915 PetscCheck(i >= 0 && i < nasm->n, PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_OUTOFRANGE, "No such subsolver"); 916 *subsnes = nasm->subsnes[i]; 917 PetscFunctionReturn(PETSC_SUCCESS); 918 } 919 920 /*@ 921 SNESNASMGetNumber - Gets number of subsolvers 922 923 Not Collective 924 925 Input Parameter: 926 . snes - the `SNES` context 927 928 Output Parameter: 929 . n - the number of subsolvers 930 931 Level: intermediate 932 933 .seealso: [](ch_snes), `SNESNASM`, `SNESNASMGetSNES()` 934 @*/ 935 PetscErrorCode SNESNASMGetNumber(SNES snes, PetscInt *n) 936 { 937 SNES_NASM *nasm = (SNES_NASM *)snes->data; 938 939 PetscFunctionBegin; 940 *n = nasm->n; 941 PetscFunctionReturn(PETSC_SUCCESS); 942 } 943 944 /*@ 945 SNESNASMSetWeight - Sets weight to use when adding overlapping updates 946 947 Collective 948 949 Input Parameters: 950 + snes - the `SNES` context 951 - weight - the weights to use (typically 1/N for each dof, where N is the number of patches it appears in) 952 953 Level: intermediate 954 955 .seealso: [](ch_snes), `SNESNASM` 956 @*/ 957 PetscErrorCode SNESNASMSetWeight(SNES snes, Vec weight) 958 { 959 SNES_NASM *nasm = (SNES_NASM *)snes->data; 960 961 PetscFunctionBegin; 962 PetscCall(VecDestroy(&nasm->weight)); 963 nasm->weight_set = PETSC_TRUE; 964 nasm->weight = weight; 965 PetscCall(PetscObjectReference((PetscObject)nasm->weight)); 966 PetscFunctionReturn(PETSC_SUCCESS); 967 } 968