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