1 /* 2 Defines a SNES that can consist of a collection of SNESes 3 */ 4 #include <petsc/private/snesimpl.h> /*I "petscsnes.h" I*/ 5 #include <petscblaslapack.h> 6 7 const char *const SNESCompositeTypes[] = {"ADDITIVE", "MULTIPLICATIVE", "ADDITIVEOPTIMAL", "SNESCompositeType", "SNES_COMPOSITE", NULL}; 8 9 typedef struct _SNES_CompositeLink *SNES_CompositeLink; 10 struct _SNES_CompositeLink { 11 SNES snes; 12 PetscReal dmp; 13 Vec X; 14 SNES_CompositeLink next; 15 SNES_CompositeLink previous; 16 }; 17 18 typedef struct { 19 SNES_CompositeLink head; 20 PetscInt nsnes; 21 SNESCompositeType type; 22 Vec Xorig; 23 PetscInt innerFailures; /* the number of inner failures we've seen */ 24 25 /* context for ADDITIVEOPTIMAL */ 26 Vec *Xes, *Fes; /* solution and residual vectors for the subsolvers */ 27 PetscReal *fnorms; /* norms of the solutions */ 28 PetscScalar *h; /* the matrix formed as q_ij = (rdot_i, rdot_j) */ 29 PetscScalar *g; /* the dotproducts of the previous function with the candidate functions */ 30 PetscBLASInt n; /* matrix dimension -- nsnes */ 31 PetscBLASInt nrhs; /* the number of right-hand sides */ 32 PetscBLASInt lda; /* the padded matrix dimension */ 33 PetscBLASInt ldb; /* the padded vector dimension */ 34 PetscReal *s; /* the singular values */ 35 PetscScalar *beta; /* the RHS and combination */ 36 PetscReal rcond; /* the exit condition */ 37 PetscBLASInt rank; /* the effective rank */ 38 PetscScalar *work; /* the work vector */ 39 PetscReal *rwork; /* the real work vector used for complex */ 40 PetscBLASInt lwork; /* the size of the work vector */ 41 PetscBLASInt info; /* the output condition */ 42 43 PetscReal rtol; /* restart tolerance for accepting the combination */ 44 PetscReal stol; /* restart tolerance for the combination */ 45 } SNES_Composite; 46 47 static PetscErrorCode SNESCompositeApply_Multiplicative(SNES snes, Vec X, Vec B, Vec F, PetscReal *fnorm) 48 { 49 SNES_Composite *jac = (SNES_Composite *)snes->data; 50 SNES_CompositeLink next = jac->head; 51 Vec FSub; 52 SNESConvergedReason reason; 53 54 PetscFunctionBegin; 55 PetscCheck(next, PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_WRONGSTATE, "No composite SNESes supplied via SNESCompositeAddSNES() or -snes_composite_sneses"); 56 if (snes->normschedule == SNES_NORM_ALWAYS) PetscCall(SNESSetInitialFunction(next->snes, F)); 57 PetscCall(SNESSolve(next->snes, B, X)); 58 PetscCall(SNESGetConvergedReason(next->snes, &reason)); 59 if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 60 jac->innerFailures++; 61 if (jac->innerFailures >= snes->maxFailures) { 62 snes->reason = SNES_DIVERGED_INNER; 63 PetscFunctionReturn(PETSC_SUCCESS); 64 } 65 } 66 67 while (next->next) { 68 /* only copy the function over in the case where the functions correspond */ 69 if (next->snes->npcside == PC_RIGHT && next->snes->normschedule != SNES_NORM_NONE) { 70 PetscCall(SNESGetFunction(next->snes, &FSub, NULL, NULL)); 71 next = next->next; 72 PetscCall(SNESSetInitialFunction(next->snes, FSub)); 73 } else { 74 next = next->next; 75 } 76 PetscCall(SNESSolve(next->snes, B, X)); 77 PetscCall(SNESGetConvergedReason(next->snes, &reason)); 78 if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 79 jac->innerFailures++; 80 if (jac->innerFailures >= snes->maxFailures) { 81 snes->reason = SNES_DIVERGED_INNER; 82 PetscFunctionReturn(PETSC_SUCCESS); 83 } 84 } 85 } 86 if (next->snes->npcside == PC_RIGHT) { 87 PetscCall(SNESGetFunction(next->snes, &FSub, NULL, NULL)); 88 PetscCall(VecCopy(FSub, F)); 89 if (fnorm) { 90 if (snes->xl && snes->xu) { 91 PetscCall(SNESVIComputeInactiveSetFnorm(snes, F, X, fnorm)); 92 } else { 93 PetscCall(VecNorm(F, NORM_2, fnorm)); 94 } 95 SNESCheckFunctionNorm(snes, *fnorm); 96 } 97 } else if (snes->normschedule == SNES_NORM_ALWAYS) { 98 PetscCall(SNESComputeFunction(snes, X, F)); 99 if (fnorm) { 100 if (snes->xl && snes->xu) { 101 PetscCall(SNESVIComputeInactiveSetFnorm(snes, F, X, fnorm)); 102 } else { 103 PetscCall(VecNorm(F, NORM_2, fnorm)); 104 } 105 SNESCheckFunctionNorm(snes, *fnorm); 106 } 107 } 108 PetscFunctionReturn(PETSC_SUCCESS); 109 } 110 111 static PetscErrorCode SNESCompositeApply_Additive(SNES snes, Vec X, Vec B, Vec F, PetscReal *fnorm) 112 { 113 SNES_Composite *jac = (SNES_Composite *)snes->data; 114 SNES_CompositeLink next = jac->head; 115 Vec Y, Xorig; 116 SNESConvergedReason reason; 117 118 PetscFunctionBegin; 119 Y = snes->vec_sol_update; 120 if (!jac->Xorig) PetscCall(VecDuplicate(X, &jac->Xorig)); 121 Xorig = jac->Xorig; 122 PetscCall(VecCopy(X, Xorig)); 123 PetscCheck(next, PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_WRONGSTATE, "No composite SNESes supplied via SNESCompositeAddSNES() or -snes_composite_sneses"); 124 if (snes->normschedule == SNES_NORM_ALWAYS) { 125 PetscCall(SNESSetInitialFunction(next->snes, F)); 126 while (next->next) { 127 next = next->next; 128 PetscCall(SNESSetInitialFunction(next->snes, F)); 129 } 130 } 131 next = jac->head; 132 PetscCall(VecCopy(Xorig, Y)); 133 PetscCall(SNESSolve(next->snes, B, Y)); 134 PetscCall(SNESGetConvergedReason(next->snes, &reason)); 135 if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 136 jac->innerFailures++; 137 if (jac->innerFailures >= snes->maxFailures) { 138 snes->reason = SNES_DIVERGED_INNER; 139 PetscFunctionReturn(PETSC_SUCCESS); 140 } 141 } 142 PetscCall(VecAXPY(Y, -1.0, Xorig)); 143 PetscCall(VecAXPY(X, next->dmp, Y)); 144 while (next->next) { 145 next = next->next; 146 PetscCall(VecCopy(Xorig, Y)); 147 PetscCall(SNESSolve(next->snes, B, Y)); 148 PetscCall(SNESGetConvergedReason(next->snes, &reason)); 149 if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 150 jac->innerFailures++; 151 if (jac->innerFailures >= snes->maxFailures) { 152 snes->reason = SNES_DIVERGED_INNER; 153 PetscFunctionReturn(PETSC_SUCCESS); 154 } 155 } 156 PetscCall(VecAXPY(Y, -1.0, Xorig)); 157 PetscCall(VecAXPY(X, next->dmp, Y)); 158 } 159 if (snes->normschedule == SNES_NORM_ALWAYS) { 160 PetscCall(SNESComputeFunction(snes, X, F)); 161 if (fnorm) { 162 if (snes->xl && snes->xu) { 163 PetscCall(SNESVIComputeInactiveSetFnorm(snes, F, X, fnorm)); 164 } else { 165 PetscCall(VecNorm(F, NORM_2, fnorm)); 166 } 167 SNESCheckFunctionNorm(snes, *fnorm); 168 } 169 } 170 PetscFunctionReturn(PETSC_SUCCESS); 171 } 172 173 /* approximately solve the overdetermined system: 174 175 2*F(x_i)\cdot F(\x_j)\alpha_i = 0 176 \alpha_i = 1 177 178 Which minimizes the L2 norm of the linearization of: 179 ||F(\sum_i \alpha_i*x_i)||^2 180 181 With the constraint that \sum_i\alpha_i = 1 182 Where x_i is the solution from the ith subsolver. 183 */ 184 static PetscErrorCode SNESCompositeApply_AdditiveOptimal(SNES snes, Vec X, Vec B, Vec F, PetscReal *fnorm) 185 { 186 SNES_Composite *jac = (SNES_Composite *)snes->data; 187 SNES_CompositeLink next = jac->head; 188 Vec *Xes = jac->Xes, *Fes = jac->Fes; 189 PetscInt i, j; 190 PetscScalar tot, total, ftf; 191 PetscReal min_fnorm; 192 PetscInt min_i; 193 SNESConvergedReason reason; 194 195 PetscFunctionBegin; 196 PetscCheck(next, PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_WRONGSTATE, "No composite SNESes supplied via SNESCompositeAddSNES() or -snes_composite_sneses"); 197 198 if (snes->normschedule == SNES_NORM_ALWAYS) { 199 next = jac->head; 200 PetscCall(SNESSetInitialFunction(next->snes, F)); 201 while (next->next) { 202 next = next->next; 203 PetscCall(SNESSetInitialFunction(next->snes, F)); 204 } 205 } 206 207 next = jac->head; 208 i = 0; 209 PetscCall(VecCopy(X, Xes[i])); 210 PetscCall(SNESSolve(next->snes, B, Xes[i])); 211 PetscCall(SNESGetConvergedReason(next->snes, &reason)); 212 if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 213 jac->innerFailures++; 214 if (jac->innerFailures >= snes->maxFailures) { 215 snes->reason = SNES_DIVERGED_INNER; 216 PetscFunctionReturn(PETSC_SUCCESS); 217 } 218 } 219 while (next->next) { 220 i++; 221 next = next->next; 222 PetscCall(VecCopy(X, Xes[i])); 223 PetscCall(SNESSolve(next->snes, B, Xes[i])); 224 PetscCall(SNESGetConvergedReason(next->snes, &reason)); 225 if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 226 jac->innerFailures++; 227 if (jac->innerFailures >= snes->maxFailures) { 228 snes->reason = SNES_DIVERGED_INNER; 229 PetscFunctionReturn(PETSC_SUCCESS); 230 } 231 } 232 } 233 234 /* all the solutions are collected; combine optimally */ 235 for (i = 0; i < jac->n; i++) { 236 for (j = 0; j < i + 1; j++) PetscCall(VecDotBegin(Fes[i], Fes[j], &jac->h[i + j * jac->n])); 237 PetscCall(VecDotBegin(Fes[i], F, &jac->g[i])); 238 } 239 240 for (i = 0; i < jac->n; i++) { 241 for (j = 0; j < i + 1; j++) { 242 PetscCall(VecDotEnd(Fes[i], Fes[j], &jac->h[i + j * jac->n])); 243 if (i == j) jac->fnorms[i] = PetscSqrtReal(PetscRealPart(jac->h[i + j * jac->n])); 244 } 245 PetscCall(VecDotEnd(Fes[i], F, &jac->g[i])); 246 } 247 248 ftf = (*fnorm) * (*fnorm); 249 250 for (i = 0; i < jac->n; i++) { 251 for (j = i + 1; j < jac->n; j++) jac->h[i + j * jac->n] = jac->h[j + i * jac->n]; 252 } 253 254 for (i = 0; i < jac->n; i++) { 255 for (j = 0; j < jac->n; j++) jac->h[i + j * jac->n] = jac->h[i + j * jac->n] - jac->g[j] - jac->g[i] + ftf; 256 jac->beta[i] = ftf - jac->g[i]; 257 } 258 259 jac->info = 0; 260 jac->rcond = -1.; 261 PetscCall(PetscFPTrapPush(PETSC_FP_TRAP_OFF)); 262 #if defined(PETSC_USE_COMPLEX) 263 PetscCallBLAS("LAPACKgelss", LAPACKgelss_(&jac->n, &jac->n, &jac->nrhs, jac->h, &jac->lda, jac->beta, &jac->lda, jac->s, &jac->rcond, &jac->rank, jac->work, &jac->lwork, jac->rwork, &jac->info)); 264 #else 265 PetscCallBLAS("LAPACKgelss", LAPACKgelss_(&jac->n, &jac->n, &jac->nrhs, jac->h, &jac->lda, jac->beta, &jac->lda, jac->s, &jac->rcond, &jac->rank, jac->work, &jac->lwork, &jac->info)); 266 #endif 267 PetscCall(PetscFPTrapPop()); 268 PetscCheck(jac->info >= 0, PetscObjectComm((PetscObject)snes), PETSC_ERR_LIB, "Bad argument to GELSS"); 269 PetscCheck(jac->info <= 0, PetscObjectComm((PetscObject)snes), PETSC_ERR_LIB, "SVD failed to converge"); 270 tot = 0.; 271 total = 0.; 272 for (i = 0; i < jac->n; i++) { 273 PetscCheck(!snes->errorifnotconverged || !PetscIsInfOrNanScalar(jac->beta[i]), PetscObjectComm((PetscObject)snes), PETSC_ERR_LIB, "SVD generated inconsistent output"); 274 PetscCall(PetscInfo(snes, "%" PetscInt_FMT ": %g\n", i, (double)PetscRealPart(jac->beta[i]))); 275 tot += jac->beta[i]; 276 total += PetscAbsScalar(jac->beta[i]); 277 } 278 PetscCall(VecScale(X, 1. - tot)); 279 PetscCall(VecMAXPY(X, jac->n, jac->beta, Xes)); 280 PetscCall(SNESComputeFunction(snes, X, F)); 281 282 if (snes->xl && snes->xu) { 283 PetscCall(SNESVIComputeInactiveSetFnorm(snes, F, X, fnorm)); 284 } else { 285 PetscCall(VecNorm(F, NORM_2, fnorm)); 286 } 287 288 /* take the minimum-normed candidate if it beats the combination by a factor of rtol or the combination has stagnated */ 289 min_fnorm = jac->fnorms[0]; 290 min_i = 0; 291 for (i = 0; i < jac->n; i++) { 292 if (jac->fnorms[i] < min_fnorm) { 293 min_fnorm = jac->fnorms[i]; 294 min_i = i; 295 } 296 } 297 298 /* stagnation or divergence restart to the solution of the solver that failed the least */ 299 if (PetscRealPart(total) < jac->stol || min_fnorm * jac->rtol < *fnorm) { 300 PetscCall(VecCopy(jac->Xes[min_i], X)); 301 PetscCall(VecCopy(jac->Fes[min_i], F)); 302 *fnorm = min_fnorm; 303 } 304 PetscFunctionReturn(PETSC_SUCCESS); 305 } 306 307 static PetscErrorCode SNESSetUp_Composite(SNES snes) 308 { 309 DM dm; 310 SNES_Composite *jac = (SNES_Composite *)snes->data; 311 SNES_CompositeLink next = jac->head; 312 PetscInt n = 0, i; 313 Vec F; 314 315 PetscFunctionBegin; 316 PetscCall(SNESGetDM(snes, &dm)); 317 318 if (snes->ops->computevariablebounds) { 319 /* SNESVI only ever calls computevariablebounds once, so calling it once here is justified */ 320 if (!snes->xl) PetscCall(VecDuplicate(snes->vec_sol, &snes->xl)); 321 if (!snes->xu) PetscCall(VecDuplicate(snes->vec_sol, &snes->xu)); 322 PetscUseTypeMethod(snes, computevariablebounds, snes->xl, snes->xu); 323 } 324 325 while (next) { 326 n++; 327 PetscCall(SNESSetDM(next->snes, dm)); 328 PetscCall(SNESSetJacobian(next->snes, snes->jacobian, snes->jacobian_pre, NULL, NULL)); 329 PetscCall(SNESSetApplicationContext(next->snes, snes->ctx)); 330 if (snes->xl && snes->xu) { 331 if (snes->ops->computevariablebounds) { 332 PetscCall(SNESVISetComputeVariableBounds(next->snes, snes->ops->computevariablebounds)); 333 } else { 334 PetscCall(SNESVISetVariableBounds(next->snes, snes->xl, snes->xu)); 335 } 336 } 337 338 next = next->next; 339 } 340 jac->nsnes = n; 341 PetscCall(SNESGetFunction(snes, &F, NULL, NULL)); 342 if (jac->type == SNES_COMPOSITE_ADDITIVEOPTIMAL) { 343 PetscCall(VecDuplicateVecs(F, jac->nsnes, &jac->Xes)); 344 PetscCall(PetscMalloc1(n, &jac->Fes)); 345 PetscCall(PetscMalloc1(n, &jac->fnorms)); 346 next = jac->head; 347 i = 0; 348 while (next) { 349 PetscCall(SNESGetFunction(next->snes, &F, NULL, NULL)); 350 jac->Fes[i] = F; 351 PetscCall(PetscObjectReference((PetscObject)F)); 352 next = next->next; 353 i++; 354 } 355 /* allocate the subspace direct solve area */ 356 jac->nrhs = 1; 357 PetscCall(PetscBLASIntCast(jac->nsnes, &jac->lda)); 358 PetscCall(PetscBLASIntCast(jac->nsnes, &jac->ldb)); 359 PetscCall(PetscBLASIntCast(jac->nsnes, &jac->n)); 360 361 PetscCall(PetscMalloc4(jac->n * jac->n, &jac->h, jac->n, &jac->beta, jac->n, &jac->s, jac->n, &jac->g)); 362 jac->lwork = 12 * jac->n; 363 #if defined(PETSC_USE_COMPLEX) 364 PetscCall(PetscMalloc1(jac->lwork, &jac->rwork)); 365 #endif 366 PetscCall(PetscMalloc1(jac->lwork, &jac->work)); 367 } 368 PetscFunctionReturn(PETSC_SUCCESS); 369 } 370 371 static PetscErrorCode SNESReset_Composite(SNES snes) 372 { 373 SNES_Composite *jac = (SNES_Composite *)snes->data; 374 SNES_CompositeLink next = jac->head; 375 376 PetscFunctionBegin; 377 while (next) { 378 PetscCall(SNESReset(next->snes)); 379 next = next->next; 380 } 381 PetscCall(VecDestroy(&jac->Xorig)); 382 if (jac->Xes) PetscCall(VecDestroyVecs(jac->nsnes, &jac->Xes)); 383 if (jac->Fes) PetscCall(VecDestroyVecs(jac->nsnes, &jac->Fes)); 384 PetscCall(PetscFree(jac->fnorms)); 385 PetscCall(PetscFree4(jac->h, jac->beta, jac->s, jac->g)); 386 PetscCall(PetscFree(jac->work)); 387 PetscCall(PetscFree(jac->rwork)); 388 PetscFunctionReturn(PETSC_SUCCESS); 389 } 390 391 static PetscErrorCode SNESDestroy_Composite(SNES snes) 392 { 393 SNES_Composite *jac = (SNES_Composite *)snes->data; 394 SNES_CompositeLink next = jac->head, next_tmp; 395 396 PetscFunctionBegin; 397 PetscCall(SNESReset_Composite(snes)); 398 while (next) { 399 PetscCall(SNESDestroy(&next->snes)); 400 next_tmp = next; 401 next = next->next; 402 PetscCall(PetscFree(next_tmp)); 403 } 404 PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESCompositeSetType_C", NULL)); 405 PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESCompositeAddSNES_C", NULL)); 406 PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESCompositeGetSNES_C", NULL)); 407 PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESCompositeSetDamping_C", NULL)); 408 PetscCall(PetscFree(snes->data)); 409 PetscFunctionReturn(PETSC_SUCCESS); 410 } 411 412 static PetscErrorCode SNESSetFromOptions_Composite(SNES snes, PetscOptionItems PetscOptionsObject) 413 { 414 SNES_Composite *jac = (SNES_Composite *)snes->data; 415 PetscInt nmax = 8, i; 416 SNES_CompositeLink next; 417 char *sneses[8]; 418 PetscReal dmps[8]; 419 PetscBool flg; 420 421 PetscFunctionBegin; 422 PetscOptionsHeadBegin(PetscOptionsObject, "Composite preconditioner options"); 423 PetscCall(PetscOptionsEnum("-snes_composite_type", "Type of composition", "SNESCompositeSetType", SNESCompositeTypes, (PetscEnum)jac->type, (PetscEnum *)&jac->type, &flg)); 424 if (flg) PetscCall(SNESCompositeSetType(snes, jac->type)); 425 PetscCall(PetscOptionsStringArray("-snes_composite_sneses", "List of composite solvers", "SNESCompositeAddSNES", sneses, &nmax, &flg)); 426 if (flg) { 427 for (i = 0; i < nmax; i++) { 428 PetscCall(SNESCompositeAddSNES(snes, sneses[i])); 429 PetscCall(PetscFree(sneses[i])); /* deallocate string sneses[i], which is allocated in PetscOptionsStringArray() */ 430 } 431 } 432 PetscCall(PetscOptionsRealArray("-snes_composite_damping", "Damping of the additive composite solvers", "SNESCompositeSetDamping", dmps, &nmax, &flg)); 433 if (flg) { 434 for (i = 0; i < nmax; i++) PetscCall(SNESCompositeSetDamping(snes, i, dmps[i])); 435 } 436 PetscCall(PetscOptionsReal("-snes_composite_stol", "Step tolerance for restart on the additive composite solvers", "", jac->stol, &jac->stol, NULL)); 437 PetscCall(PetscOptionsReal("-snes_composite_rtol", "Residual tolerance for the additive composite solvers", "", jac->rtol, &jac->rtol, NULL)); 438 PetscOptionsHeadEnd(); 439 440 next = jac->head; 441 while (next) { 442 PetscCall(SNESSetFromOptions(next->snes)); 443 next = next->next; 444 } 445 PetscFunctionReturn(PETSC_SUCCESS); 446 } 447 448 static PetscErrorCode SNESView_Composite(SNES snes, PetscViewer viewer) 449 { 450 SNES_Composite *jac = (SNES_Composite *)snes->data; 451 SNES_CompositeLink next = jac->head; 452 PetscBool isascii; 453 454 PetscFunctionBegin; 455 PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii)); 456 if (isascii) { 457 PetscCall(PetscViewerASCIIPrintf(viewer, " type - %s\n", SNESCompositeTypes[jac->type])); 458 PetscCall(PetscViewerASCIIPrintf(viewer, " SNESes on composite preconditioner follow\n")); 459 PetscCall(PetscViewerASCIIPrintf(viewer, " ---------------------------------\n")); 460 } 461 if (isascii) PetscCall(PetscViewerASCIIPushTab(viewer)); 462 while (next) { 463 PetscCall(SNESView(next->snes, viewer)); 464 next = next->next; 465 } 466 if (isascii) { 467 PetscCall(PetscViewerASCIIPopTab(viewer)); 468 PetscCall(PetscViewerASCIIPrintf(viewer, " ---------------------------------\n")); 469 } 470 PetscFunctionReturn(PETSC_SUCCESS); 471 } 472 473 static PetscErrorCode SNESCompositeSetType_Composite(SNES snes, SNESCompositeType type) 474 { 475 SNES_Composite *jac = (SNES_Composite *)snes->data; 476 477 PetscFunctionBegin; 478 jac->type = type; 479 PetscFunctionReturn(PETSC_SUCCESS); 480 } 481 482 static PetscErrorCode SNESCompositeAddSNES_Composite(SNES snes, SNESType type) 483 { 484 SNES_Composite *jac; 485 SNES_CompositeLink next, ilink; 486 PetscInt cnt = 0; 487 const char *prefix; 488 char newprefix[20]; 489 DM dm; 490 491 PetscFunctionBegin; 492 PetscCall(PetscNew(&ilink)); 493 ilink->next = NULL; 494 PetscCall(SNESCreate(PetscObjectComm((PetscObject)snes), &ilink->snes)); 495 PetscCall(PetscObjectIncrementTabLevel((PetscObject)ilink->snes, (PetscObject)snes, 1)); 496 PetscCall(SNESGetDM(snes, &dm)); 497 PetscCall(SNESSetDM(ilink->snes, dm)); 498 PetscCall(SNESSetTolerances(ilink->snes, snes->abstol, snes->rtol, snes->stol, 1, snes->max_funcs)); 499 PetscCall(PetscObjectCopyFortranFunctionPointers((PetscObject)snes, (PetscObject)ilink->snes)); 500 jac = (SNES_Composite *)snes->data; 501 next = jac->head; 502 if (!next) { 503 jac->head = ilink; 504 ilink->previous = NULL; 505 } else { 506 cnt++; 507 while (next->next) { 508 next = next->next; 509 cnt++; 510 } 511 next->next = ilink; 512 ilink->previous = next; 513 } 514 PetscCall(SNESGetOptionsPrefix(snes, &prefix)); 515 PetscCall(SNESSetOptionsPrefix(ilink->snes, prefix)); 516 PetscCall(PetscSNPrintf(newprefix, sizeof(newprefix), "sub_%" PetscInt_FMT "_", cnt)); 517 PetscCall(SNESAppendOptionsPrefix(ilink->snes, newprefix)); 518 PetscCall(SNESSetType(ilink->snes, type)); 519 PetscCall(SNESSetNormSchedule(ilink->snes, SNES_NORM_FINAL_ONLY)); 520 521 ilink->dmp = 1.0; 522 jac->nsnes++; 523 PetscFunctionReturn(PETSC_SUCCESS); 524 } 525 526 static PetscErrorCode SNESCompositeGetSNES_Composite(SNES snes, PetscInt n, SNES *subsnes) 527 { 528 SNES_Composite *jac; 529 SNES_CompositeLink next; 530 PetscInt i; 531 532 PetscFunctionBegin; 533 jac = (SNES_Composite *)snes->data; 534 next = jac->head; 535 for (i = 0; i < n; i++) { 536 PetscCheck(next->next, PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_INCOMP, "Not enough SNESes in composite preconditioner"); 537 next = next->next; 538 } 539 *subsnes = next->snes; 540 PetscFunctionReturn(PETSC_SUCCESS); 541 } 542 543 /*@ 544 SNESCompositeSetType - Sets the type of composite preconditioner. 545 546 Logically Collective 547 548 Input Parameters: 549 + snes - the preconditioner context 550 - type - `SNES_COMPOSITE_ADDITIVE` (default), `SNES_COMPOSITE_MULTIPLICATIVE`, or `SNES_COMPOSITE_ADDITIVEOPTIMAL` 551 552 Options Database Key: 553 . -snes_composite_type <type: one of multiplicative, additive, additiveoptimal> - Sets composite preconditioner type 554 555 Level: developer 556 557 .seealso: [](ch_snes), `SNES_COMPOSITE_ADDITIVE`, `SNES_COMPOSITE_MULTIPLICATIVE`, `SNESCompositeType`, `SNESCOMPOSITE`, `SNES_COMPOSITE_ADDITIVEOPTIMAL`, 558 `PCCompositeType` 559 @*/ 560 PetscErrorCode SNESCompositeSetType(SNES snes, SNESCompositeType type) 561 { 562 PetscFunctionBegin; 563 PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 564 PetscValidLogicalCollectiveEnum(snes, type, 2); 565 PetscTryMethod(snes, "SNESCompositeSetType_C", (SNES, SNESCompositeType), (snes, type)); 566 PetscFunctionReturn(PETSC_SUCCESS); 567 } 568 569 /*@ 570 SNESCompositeAddSNES - Adds another `SNES` to the `SNESCOMPOSITE` 571 572 Collective 573 574 Input Parameters: 575 + snes - the `SNES` context of type `SNESCOMPOSITE` 576 - type - the `SNESType` of the new solver 577 578 Level: developer 579 580 .seealso: [](ch_snes), `SNES`, `SNESCOMPOSITE`, `SNESCompositeGetSNES()` 581 @*/ 582 PetscErrorCode SNESCompositeAddSNES(SNES snes, SNESType type) 583 { 584 PetscFunctionBegin; 585 PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 586 PetscTryMethod(snes, "SNESCompositeAddSNES_C", (SNES, SNESType), (snes, type)); 587 PetscFunctionReturn(PETSC_SUCCESS); 588 } 589 590 /*@ 591 SNESCompositeGetSNES - Gets one of the `SNES` objects in the `SNES` of `SNESType` `SNESCOMPOSITE` 592 593 Not Collective 594 595 Input Parameters: 596 + snes - the `SNES` context 597 - n - the number of the composed `SNES` requested 598 599 Output Parameter: 600 . subsnes - the `SNES` requested 601 602 Level: developer 603 604 .seealso: [](ch_snes), `SNES`, `SNESCOMPOSITE`, `SNESCompositeAddSNES()`, `SNESCompositeGetNumber()` 605 @*/ 606 PetscErrorCode SNESCompositeGetSNES(SNES snes, PetscInt n, SNES *subsnes) 607 { 608 PetscFunctionBegin; 609 PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 610 PetscAssertPointer(subsnes, 3); 611 PetscUseMethod(snes, "SNESCompositeGetSNES_C", (SNES, PetscInt, SNES *), (snes, n, subsnes)); 612 PetscFunctionReturn(PETSC_SUCCESS); 613 } 614 615 /*@ 616 SNESCompositeGetNumber - Get the number of subsolvers in the `SNESCOMPOSITE` 617 618 Logically Collective 619 620 Input Parameter: 621 . snes - the `SNES` context 622 623 Output Parameter: 624 . n - the number of subsolvers 625 626 Level: developer 627 628 .seealso: [](ch_snes), `SNES`, `SNESCOMPOSITE`, `SNESCompositeAddSNES()`, `SNESCompositeGetSNES()` 629 @*/ 630 PetscErrorCode SNESCompositeGetNumber(SNES snes, PetscInt *n) 631 { 632 SNES_Composite *jac; 633 SNES_CompositeLink next; 634 635 PetscFunctionBegin; 636 jac = (SNES_Composite *)snes->data; 637 next = jac->head; 638 639 *n = 0; 640 while (next) { 641 *n = *n + 1; 642 next = next->next; 643 } 644 PetscFunctionReturn(PETSC_SUCCESS); 645 } 646 647 static PetscErrorCode SNESCompositeSetDamping_Composite(SNES snes, PetscInt n, PetscReal dmp) 648 { 649 SNES_Composite *jac; 650 SNES_CompositeLink next; 651 PetscInt i; 652 653 PetscFunctionBegin; 654 jac = (SNES_Composite *)snes->data; 655 next = jac->head; 656 for (i = 0; i < n; i++) { 657 PetscCheck(next->next, PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_INCOMP, "Not enough SNESes in composite preconditioner"); 658 next = next->next; 659 } 660 next->dmp = dmp; 661 PetscFunctionReturn(PETSC_SUCCESS); 662 } 663 664 /*@ 665 SNESCompositeSetDamping - Sets the damping of a subsolver when using `SNES_COMPOSITE_ADDITIVE` with a `SNES` of `SNESType` `SNESCOMPOSITE` 666 667 Not Collective 668 669 Input Parameters: 670 + snes - the `SNES` context 671 . n - the number of the sub-`SNES` object requested 672 - dmp - the damping 673 674 Level: intermediate 675 676 .seealso: [](ch_snes), `SNES`, `SNESCOMPOSITE`, `SNESCompositeAddSNES()`, `SNESCompositeGetSNES()`, 677 `SNES_COMPOSITE_ADDITIVE`, `SNES_COMPOSITE_MULTIPLICATIVE`, `SNESCompositeType`, `SNESCompositeSetType()` 678 @*/ 679 PetscErrorCode SNESCompositeSetDamping(SNES snes, PetscInt n, PetscReal dmp) 680 { 681 PetscFunctionBegin; 682 PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 683 PetscUseMethod(snes, "SNESCompositeSetDamping_C", (SNES, PetscInt, PetscReal), (snes, n, dmp)); 684 PetscFunctionReturn(PETSC_SUCCESS); 685 } 686 687 static PetscErrorCode SNESSolve_Composite(SNES snes) 688 { 689 Vec F, X, B, Y; 690 PetscInt i; 691 PetscReal fnorm = 0.0, xnorm = 0.0, snorm = 0.0; 692 SNESNormSchedule normtype; 693 SNES_Composite *comp = (SNES_Composite *)snes->data; 694 695 PetscFunctionBegin; 696 X = snes->vec_sol; 697 F = snes->vec_func; 698 B = snes->vec_rhs; 699 Y = snes->vec_sol_update; 700 701 PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes)); 702 snes->iter = 0; 703 snes->norm = 0.; 704 comp->innerFailures = 0; 705 PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 706 snes->reason = SNES_CONVERGED_ITERATING; 707 PetscCall(SNESGetNormSchedule(snes, &normtype)); 708 if (normtype == SNES_NORM_ALWAYS || normtype == SNES_NORM_INITIAL_ONLY || normtype == SNES_NORM_INITIAL_FINAL_ONLY) { 709 if (!snes->vec_func_init_set) PetscCall(SNESComputeFunction(snes, X, F)); 710 else snes->vec_func_init_set = PETSC_FALSE; 711 712 if (snes->xl && snes->xu) { 713 PetscCall(SNESVIComputeInactiveSetFnorm(snes, F, X, &fnorm)); 714 } else { 715 PetscCall(VecNorm(F, NORM_2, &fnorm)); /* fnorm <- ||F|| */ 716 } 717 SNESCheckFunctionNorm(snes, fnorm); 718 PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes)); 719 snes->iter = 0; 720 snes->norm = fnorm; 721 PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 722 PetscCall(SNESLogConvergenceHistory(snes, snes->norm, 0)); 723 724 /* test convergence */ 725 PetscCall(SNESConverged(snes, 0, 0.0, 0.0, fnorm)); 726 PetscCall(SNESMonitor(snes, 0, snes->norm)); 727 if (snes->reason) PetscFunctionReturn(PETSC_SUCCESS); 728 } else { 729 PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 730 PetscCall(SNESLogConvergenceHistory(snes, snes->norm, 0)); 731 PetscCall(SNESMonitor(snes, 0, snes->norm)); 732 } 733 734 for (i = 0; i < snes->max_its; i++) { 735 /* Call general purpose update function */ 736 PetscTryTypeMethod(snes, update, snes->iter); 737 738 /* Copy the state before modification by application of the composite solver; 739 we will subtract the new state after application */ 740 PetscCall(VecCopy(X, Y)); 741 742 if (comp->type == SNES_COMPOSITE_ADDITIVE) { 743 PetscCall(SNESCompositeApply_Additive(snes, X, B, F, &fnorm)); 744 } else if (comp->type == SNES_COMPOSITE_MULTIPLICATIVE) { 745 PetscCall(SNESCompositeApply_Multiplicative(snes, X, B, F, &fnorm)); 746 } else if (comp->type == SNES_COMPOSITE_ADDITIVEOPTIMAL) { 747 PetscCall(SNESCompositeApply_AdditiveOptimal(snes, X, B, F, &fnorm)); 748 } else SETERRQ(PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_WRONGSTATE, "Unsupported SNESComposite type"); 749 if (snes->reason < 0) break; 750 751 /* Compute the solution update for convergence testing */ 752 PetscCall(VecAYPX(Y, -1.0, X)); 753 754 if ((i == snes->max_its - 1) && (normtype == SNES_NORM_INITIAL_FINAL_ONLY || normtype == SNES_NORM_FINAL_ONLY)) { 755 PetscCall(SNESComputeFunction(snes, X, F)); 756 757 if (snes->xl && snes->xu) { 758 PetscCall(VecNormBegin(X, NORM_2, &xnorm)); 759 PetscCall(VecNormBegin(Y, NORM_2, &snorm)); 760 PetscCall(SNESVIComputeInactiveSetFnorm(snes, F, X, &fnorm)); 761 PetscCall(VecNormEnd(X, NORM_2, &xnorm)); 762 PetscCall(VecNormEnd(Y, NORM_2, &snorm)); 763 } else { 764 PetscCall(VecNormBegin(F, NORM_2, &fnorm)); 765 PetscCall(VecNormBegin(X, NORM_2, &xnorm)); 766 PetscCall(VecNormBegin(Y, NORM_2, &snorm)); 767 768 PetscCall(VecNormEnd(F, NORM_2, &fnorm)); 769 PetscCall(VecNormEnd(X, NORM_2, &xnorm)); 770 PetscCall(VecNormEnd(Y, NORM_2, &snorm)); 771 } 772 SNESCheckFunctionNorm(snes, fnorm); 773 } else if (normtype == SNES_NORM_ALWAYS) { 774 PetscCall(VecNormBegin(X, NORM_2, &xnorm)); 775 PetscCall(VecNormBegin(Y, NORM_2, &snorm)); 776 PetscCall(VecNormEnd(X, NORM_2, &xnorm)); 777 PetscCall(VecNormEnd(Y, NORM_2, &snorm)); 778 } 779 /* Monitor convergence */ 780 PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes)); 781 snes->iter = i + 1; 782 snes->norm = fnorm; 783 snes->xnorm = xnorm; 784 snes->ynorm = snorm; 785 PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 786 PetscCall(SNESLogConvergenceHistory(snes, snes->norm, 0)); 787 /* Test for convergence */ 788 PetscCall(SNESConverged(snes, snes->iter, xnorm, snorm, fnorm)); 789 PetscCall(SNESMonitor(snes, snes->iter, snes->norm)); 790 if (snes->reason) break; 791 } 792 PetscFunctionReturn(PETSC_SUCCESS); 793 } 794 795 /*MC 796 SNESCOMPOSITE - Builds a nonlinear solver/preconditioner by composing together several `SNES` nonlinear solvers {cite}`bruneknepleysmithtu15` 797 798 Options Database Keys: 799 + -snes_composite_type <type: one of multiplicative, additive, additiveoptimal> - Sets composite preconditioner type 800 - -snes_composite_sneses - <snes0,snes1,...> list of `SNES` to compose 801 802 Level: intermediate 803 804 .seealso: [](ch_snes), `SNES`, `SNESCOMPOSITE`, `SNESCompositeAddSNES()`, `SNESCompositeGetSNES()`, 805 `SNES_COMPOSITE_ADDITIVE`, `SNES_COMPOSITE_ADDITIVEOPTIMAL`, `SNES_COMPOSITE_MULTIPLICATIVE`, `SNESCompositeType`, 806 `SNESCompositeSetType()`, `SNESCompositeSetDamping()`, `SNESCompositeGetNumber()`, `PCCompositeType` 807 M*/ 808 809 PETSC_EXTERN PetscErrorCode SNESCreate_Composite(SNES snes) 810 { 811 SNES_Composite *jac; 812 813 PetscFunctionBegin; 814 PetscCall(PetscNew(&jac)); 815 816 snes->ops->solve = SNESSolve_Composite; 817 snes->ops->setup = SNESSetUp_Composite; 818 snes->ops->reset = SNESReset_Composite; 819 snes->ops->destroy = SNESDestroy_Composite; 820 snes->ops->setfromoptions = SNESSetFromOptions_Composite; 821 snes->ops->view = SNESView_Composite; 822 823 snes->usesksp = PETSC_FALSE; 824 825 snes->alwayscomputesfinalresidual = PETSC_FALSE; 826 827 PetscCall(SNESParametersInitialize(snes)); 828 829 snes->data = (void *)jac; 830 jac->type = SNES_COMPOSITE_ADDITIVEOPTIMAL; 831 jac->Fes = NULL; 832 jac->Xes = NULL; 833 jac->fnorms = NULL; 834 jac->nsnes = 0; 835 jac->head = NULL; 836 jac->stol = 0.1; 837 jac->rtol = 1.1; 838 839 jac->h = NULL; 840 jac->s = NULL; 841 jac->beta = NULL; 842 jac->work = NULL; 843 jac->rwork = NULL; 844 845 PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESCompositeSetType_C", SNESCompositeSetType_Composite)); 846 PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESCompositeAddSNES_C", SNESCompositeAddSNES_Composite)); 847 PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESCompositeGetSNES_C", SNESCompositeGetSNES_Composite)); 848 PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESCompositeSetDamping_C", SNESCompositeSetDamping_Composite)); 849 PetscFunctionReturn(PETSC_SUCCESS); 850 } 851