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