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