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