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 = SNESGetFunctionNorm(next->snes,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 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); 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: %f\n",i,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 = PetscMalloc(sizeof(Vec)*n,&jac->Fes);CHKERRQ(ierr); 324 ierr = PetscMalloc(sizeof(PetscReal)*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 = PetscMalloc(jac->n*jac->n*sizeof(PetscScalar),&jac->h);CHKERRQ(ierr); 341 ierr = PetscMalloc(jac->n*sizeof(PetscScalar),&jac->beta);CHKERRQ(ierr); 342 ierr = PetscMalloc(jac->n*sizeof(PetscScalar),&jac->s);CHKERRQ(ierr); 343 ierr = PetscMalloc(jac->n*sizeof(PetscScalar),&jac->g);CHKERRQ(ierr); 344 jac->lwork = 12*jac->n; 345 #if PETSC_USE_COMPLEX 346 ierr = PetscMalloc(sizeof(PetscReal)*jac->lwork,&jac->rwork);CHKERRQ(ierr); 347 #endif 348 ierr = PetscMalloc(sizeof(PetscScalar)*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 378 PetscFunctionReturn(0); 379 } 380 381 #undef __FUNCT__ 382 #define __FUNCT__ "SNESDestroy_Composite" 383 static PetscErrorCode SNESDestroy_Composite(SNES snes) 384 { 385 SNES_Composite *jac = (SNES_Composite*)snes->data; 386 PetscErrorCode ierr; 387 SNES_CompositeLink next = jac->head,next_tmp; 388 389 PetscFunctionBegin; 390 ierr = SNESReset_Composite(snes);CHKERRQ(ierr); 391 while (next) { 392 ierr = SNESDestroy(&next->snes);CHKERRQ(ierr); 393 next_tmp = next; 394 next = next->next; 395 ierr = PetscFree(next_tmp);CHKERRQ(ierr); 396 } 397 ierr = PetscFree(snes->data);CHKERRQ(ierr); 398 PetscFunctionReturn(0); 399 } 400 401 #undef __FUNCT__ 402 #define __FUNCT__ "SNESSetFromOptions_Composite" 403 static PetscErrorCode SNESSetFromOptions_Composite(SNES snes) 404 { 405 SNES_Composite *jac = (SNES_Composite*)snes->data; 406 PetscErrorCode ierr; 407 PetscInt nmax = 8,i; 408 SNES_CompositeLink next; 409 char *sneses[8]; 410 PetscReal dmps[8]; 411 PetscBool flg; 412 413 PetscFunctionBegin; 414 ierr = PetscOptionsHead("Composite preconditioner options");CHKERRQ(ierr); 415 ierr = PetscOptionsEnum("-snes_composite_type","Type of composition","SNESCompositeSetType",SNESCompositeTypes,(PetscEnum)jac->type,(PetscEnum*)&jac->type,&flg);CHKERRQ(ierr); 416 if (flg) { 417 ierr = SNESCompositeSetType(snes,jac->type);CHKERRQ(ierr); 418 } 419 ierr = PetscOptionsStringArray("-snes_composite_sneses","List of composite solvers","SNESCompositeAddSNES",sneses,&nmax,&flg);CHKERRQ(ierr); 420 if (flg) { 421 for (i=0; i<nmax; i++) { 422 ierr = SNESCompositeAddSNES(snes,sneses[i]);CHKERRQ(ierr); 423 ierr = PetscFree(sneses[i]);CHKERRQ(ierr); /* deallocate string sneses[i], which is allocated in PetscOptionsStringArray() */ 424 } 425 } 426 ierr = PetscOptionsRealArray("-snes_composite_damping","Damping of the additive composite solvers","SNESCompositeSetDamping",dmps,&nmax,&flg);CHKERRQ(ierr); 427 if (flg) { 428 for (i=0; i<nmax; i++) { 429 ierr = SNESCompositeSetDamping(snes,i,dmps[i]);CHKERRQ(ierr); 430 } 431 } 432 ierr = PetscOptionsReal("-snes_composite_stol","Step tolerance for restart on the additive composite solvers","",jac->stol,&jac->stol,NULL);CHKERRQ(ierr); 433 ierr = PetscOptionsReal("-snes_composite_rtol","Residual tolerance for the additive composite solvers","",jac->rtol,&jac->rtol,NULL);CHKERRQ(ierr); 434 ierr = PetscOptionsTail();CHKERRQ(ierr); 435 436 next = jac->head; 437 while (next) { 438 ierr = SNESSetFromOptions(next->snes);CHKERRQ(ierr); 439 next = next->next; 440 } 441 PetscFunctionReturn(0); 442 } 443 444 #undef __FUNCT__ 445 #define __FUNCT__ "SNESView_Composite" 446 static PetscErrorCode SNESView_Composite(SNES snes,PetscViewer viewer) 447 { 448 SNES_Composite *jac = (SNES_Composite*)snes->data; 449 PetscErrorCode ierr; 450 SNES_CompositeLink next = jac->head; 451 PetscBool iascii; 452 453 PetscFunctionBegin; 454 ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 455 if (iascii) { 456 ierr = PetscViewerASCIIPrintf(viewer,"Composite SNES type - %s\n",SNESCompositeTypes[jac->type]);CHKERRQ(ierr); 457 ierr = PetscViewerASCIIPrintf(viewer,"SNESes on composite preconditioner follow\n");CHKERRQ(ierr); 458 ierr = PetscViewerASCIIPrintf(viewer,"---------------------------------\n");CHKERRQ(ierr); 459 } 460 if (iascii) { 461 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 462 } 463 while (next) { 464 ierr = SNESView(next->snes,viewer);CHKERRQ(ierr); 465 next = next->next; 466 } 467 if (iascii) { 468 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 469 ierr = PetscViewerASCIIPrintf(viewer,"---------------------------------\n");CHKERRQ(ierr); 470 } 471 PetscFunctionReturn(0); 472 } 473 474 /* ------------------------------------------------------------------------------*/ 475 476 #undef __FUNCT__ 477 #define __FUNCT__ "SNESCompositeSetType_Composite" 478 static PetscErrorCode SNESCompositeSetType_Composite(SNES snes,SNESCompositeType type) 479 { 480 SNES_Composite *jac = (SNES_Composite*)snes->data; 481 482 PetscFunctionBegin; 483 jac->type = type; 484 PetscFunctionReturn(0); 485 } 486 487 #undef __FUNCT__ 488 #define __FUNCT__ "SNESCompositeAddSNES_Composite" 489 static PetscErrorCode SNESCompositeAddSNES_Composite(SNES snes,SNESType type) 490 { 491 SNES_Composite *jac; 492 SNES_CompositeLink next,ilink; 493 PetscErrorCode ierr; 494 PetscInt cnt = 0; 495 const char *prefix; 496 char newprefix[8]; 497 DM dm; 498 499 PetscFunctionBegin; 500 ierr = PetscNewLog(snes,struct _SNES_CompositeLink,&ilink);CHKERRQ(ierr); 501 ilink->next = 0; 502 ierr = SNESCreate(PetscObjectComm((PetscObject)snes),&ilink->snes);CHKERRQ(ierr); 503 ierr = PetscLogObjectParent((PetscObject)snes,(PetscObject)ilink->snes);CHKERRQ(ierr); 504 ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 505 ierr = SNESSetDM(ilink->snes,dm);CHKERRQ(ierr); 506 ierr = SNESSetTolerances(ilink->snes,ilink->snes->abstol,ilink->snes->rtol,ilink->snes->stol,1,ilink->snes->max_funcs);CHKERRQ(ierr); 507 jac = (SNES_Composite*)snes->data; 508 next = jac->head; 509 if (!next) { 510 jac->head = ilink; 511 ilink->previous = NULL; 512 } else { 513 cnt++; 514 while (next->next) { 515 next = next->next; 516 cnt++; 517 } 518 next->next = ilink; 519 ilink->previous = next; 520 } 521 ierr = SNESGetOptionsPrefix(snes,&prefix);CHKERRQ(ierr); 522 ierr = SNESSetOptionsPrefix(ilink->snes,prefix);CHKERRQ(ierr); 523 sprintf(newprefix,"sub_%d_",(int)cnt); 524 ierr = SNESAppendOptionsPrefix(ilink->snes,newprefix);CHKERRQ(ierr); 525 ierr = PetscObjectIncrementTabLevel((PetscObject)ilink->snes,(PetscObject)snes,1);CHKERRQ(ierr); 526 ierr = SNESSetType(ilink->snes,type);CHKERRQ(ierr); 527 ierr = SNESSetNormSchedule(ilink->snes, SNES_NORM_FINAL_ONLY);CHKERRQ(ierr); 528 ilink->dmp = 1.0; 529 jac->nsnes++; 530 PetscFunctionReturn(0); 531 } 532 533 #undef __FUNCT__ 534 #define __FUNCT__ "SNESCompositeGetSNES_Composite" 535 static PetscErrorCode SNESCompositeGetSNES_Composite(SNES snes,PetscInt n,SNES *subsnes) 536 { 537 SNES_Composite *jac; 538 SNES_CompositeLink next; 539 PetscInt i; 540 541 PetscFunctionBegin; 542 jac = (SNES_Composite*)snes->data; 543 next = jac->head; 544 for (i=0; i<n; i++) { 545 if (!next->next) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_INCOMP,"Not enough SNESes in composite preconditioner"); 546 next = next->next; 547 } 548 *subsnes = next->snes; 549 PetscFunctionReturn(0); 550 } 551 552 /* -------------------------------------------------------------------------------- */ 553 #undef __FUNCT__ 554 #define __FUNCT__ "SNESCompositeSetType" 555 /*@C 556 SNESCompositeSetType - Sets the type of composite preconditioner. 557 558 Logically Collective on SNES 559 560 Input Parameter: 561 + snes - the preconditioner context 562 - type - SNES_COMPOSITE_ADDITIVE (default), SNES_COMPOSITE_MULTIPLICATIVE 563 564 Options Database Key: 565 . -snes_composite_type <type: one of multiplicative, additive, special> - Sets composite preconditioner type 566 567 Level: Developer 568 569 .keywords: SNES, set, type, composite preconditioner, additive, multiplicative 570 @*/ 571 PetscErrorCode SNESCompositeSetType(SNES snes,SNESCompositeType type) 572 { 573 PetscErrorCode ierr; 574 575 PetscFunctionBegin; 576 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 577 PetscValidLogicalCollectiveEnum(snes,type,2); 578 ierr = PetscTryMethod(snes,"SNESCompositeSetType_C",(SNES,SNESCompositeType),(snes,type));CHKERRQ(ierr); 579 PetscFunctionReturn(0); 580 } 581 582 #undef __FUNCT__ 583 #define __FUNCT__ "SNESCompositeAddSNES" 584 /*@C 585 SNESCompositeAddSNES - Adds another SNES to the composite SNES. 586 587 Collective on SNES 588 589 Input Parameters: 590 + snes - the preconditioner context 591 - type - the type of the new preconditioner 592 593 Level: Developer 594 595 .keywords: SNES, composite preconditioner, add 596 @*/ 597 PetscErrorCode SNESCompositeAddSNES(SNES snes,SNESType type) 598 { 599 PetscErrorCode ierr; 600 601 PetscFunctionBegin; 602 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 603 ierr = PetscTryMethod(snes,"SNESCompositeAddSNES_C",(SNES,SNESType),(snes,type));CHKERRQ(ierr); 604 PetscFunctionReturn(0); 605 } 606 #undef __FUNCT__ 607 #define __FUNCT__ "SNESCompositeGetSNES" 608 /*@ 609 SNESCompositeGetSNES - Gets one of the SNES objects in the composite SNES. 610 611 Not Collective 612 613 Input Parameter: 614 + snes - the preconditioner context 615 - n - the number of the snes requested 616 617 Output Parameters: 618 . subsnes - the SNES requested 619 620 Level: Developer 621 622 .keywords: SNES, get, composite preconditioner, sub preconditioner 623 624 .seealso: SNESCompositeAddSNES() 625 @*/ 626 PetscErrorCode SNESCompositeGetSNES(SNES snes,PetscInt n,SNES *subsnes) 627 { 628 PetscErrorCode ierr; 629 630 PetscFunctionBegin; 631 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 632 PetscValidPointer(subsnes,3); 633 ierr = PetscUseMethod(snes,"SNESCompositeGetSNES_C",(SNES,PetscInt,SNES*),(snes,n,subsnes));CHKERRQ(ierr); 634 PetscFunctionReturn(0); 635 } 636 637 #undef __FUNCT__ 638 #define __FUNCT__ "SNESCompositeSetDamping_Composite" 639 static PetscErrorCode SNESCompositeSetDamping_Composite(SNES snes,PetscInt n,PetscReal dmp) 640 { 641 SNES_Composite *jac; 642 SNES_CompositeLink next; 643 PetscInt i; 644 645 PetscFunctionBegin; 646 jac = (SNES_Composite*)snes->data; 647 next = jac->head; 648 for (i=0; i<n; i++) { 649 if (!next->next) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_INCOMP,"Not enough SNESes in composite preconditioner"); 650 next = next->next; 651 } 652 next->dmp = dmp; 653 PetscFunctionReturn(0); 654 } 655 656 #undef __FUNCT__ 657 #define __FUNCT__ "SNESCompositeSetDamping" 658 /*@ 659 SNESCompositeSetDamping - Sets the damping of a subsolver when using additive composite SNES. 660 661 Not Collective 662 663 Input Parameter: 664 + snes - the preconditioner context 665 . n - the number of the snes requested 666 - dmp - the damping 667 668 Level: Developer 669 670 .keywords: SNES, get, composite preconditioner, sub preconditioner 671 672 .seealso: SNESCompositeAddSNES() 673 @*/ 674 PetscErrorCode SNESCompositeSetDamping(SNES snes,PetscInt n,PetscReal dmp) 675 { 676 PetscErrorCode ierr; 677 678 PetscFunctionBegin; 679 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 680 ierr = PetscUseMethod(snes,"SNESCompositeSetDamping_C",(SNES,PetscInt,PetscReal),(snes,n,dmp));CHKERRQ(ierr); 681 PetscFunctionReturn(0); 682 } 683 684 #undef __FUNCT__ 685 #define __FUNCT__ "SNESSolve_Composite" 686 PetscErrorCode SNESSolve_Composite(SNES snes) 687 { 688 Vec F; 689 Vec X; 690 Vec B; 691 PetscInt i; 692 PetscReal fnorm = 0.0; 693 PetscErrorCode ierr; 694 SNESNormSchedule normtype; 695 SNES_Composite *comp = (SNES_Composite*)snes->data; 696 697 PetscFunctionBegin; 698 X = snes->vec_sol; 699 F = snes->vec_func; 700 B = snes->vec_rhs; 701 702 ierr = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr); 703 snes->iter = 0; 704 snes->norm = 0.; 705 ierr = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr); 706 snes->reason = SNES_CONVERGED_ITERATING; 707 ierr = SNESGetNormSchedule(snes, &normtype);CHKERRQ(ierr); 708 if (normtype == SNES_NORM_ALWAYS || normtype == SNES_NORM_INITIAL_ONLY || normtype == SNES_NORM_INITIAL_FINAL_ONLY) { 709 if (!snes->vec_func_init_set) { 710 ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 711 if (snes->domainerror) { 712 snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 713 PetscFunctionReturn(0); 714 } 715 } else snes->vec_func_init_set = PETSC_FALSE; 716 717 ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F|| */ 718 if (PetscIsInfOrNanReal(fnorm)) { 719 snes->reason = SNES_DIVERGED_FNORM_NAN; 720 PetscFunctionReturn(0); 721 } 722 ierr = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr); 723 snes->iter = 0; 724 snes->norm = fnorm; 725 ierr = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr); 726 ierr = SNESLogConvergenceHistory(snes,snes->norm,0);CHKERRQ(ierr); 727 ierr = SNESMonitor(snes,0,snes->norm);CHKERRQ(ierr); 728 729 /* set parameter for default relative tolerance convergence test */ 730 snes->ttol = fnorm*snes->rtol; 731 732 /* test convergence */ 733 ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 734 if (snes->reason) PetscFunctionReturn(0); 735 } else { 736 ierr = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr); 737 ierr = SNESLogConvergenceHistory(snes,snes->norm,0);CHKERRQ(ierr); 738 ierr = SNESMonitor(snes,0,snes->norm);CHKERRQ(ierr); 739 } 740 741 /* Call general purpose update function */ 742 if (snes->ops->update) { 743 ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 744 } 745 746 for (i = 0; i < snes->max_its; i++) { 747 if (comp->type == SNES_COMPOSITE_ADDITIVE) { 748 ierr = SNESCompositeApply_Additive(snes,X,B,F,&fnorm);CHKERRQ(ierr); 749 } else if (comp->type == SNES_COMPOSITE_MULTIPLICATIVE) { 750 ierr = SNESCompositeApply_Multiplicative(snes,X,B,F,&fnorm);CHKERRQ(ierr); 751 } else if (comp->type == SNES_COMPOSITE_ADDITIVEOPTIMAL) { 752 ierr = SNESCompositeApply_AdditiveOptimal(snes,X,B,F,&fnorm);CHKERRQ(ierr); 753 } else { 754 SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"Unsupported SNESComposite type"); 755 } 756 if (snes->reason < 0) break; 757 758 if ((i == snes->max_its - 1) && (normtype == SNES_NORM_INITIAL_FINAL_ONLY || normtype == SNES_NORM_FINAL_ONLY)) { 759 ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 760 if (snes->domainerror) { 761 snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 762 break; 763 } 764 ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); 765 if (PetscIsInfOrNanReal(fnorm)) { 766 snes->reason = SNES_DIVERGED_FNORM_NAN; 767 break; 768 } 769 } 770 /* Monitor convergence */ 771 ierr = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr); 772 snes->iter = i+1; 773 snes->norm = fnorm; 774 ierr = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr); 775 ierr = SNESLogConvergenceHistory(snes,snes->norm,0);CHKERRQ(ierr); 776 ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr); 777 /* Test for convergence */ 778 if (normtype == SNES_NORM_ALWAYS) {ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);} 779 if (snes->reason) break; 780 /* Call general purpose update function */ 781 if (snes->ops->update) {ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);} 782 } 783 if (normtype == SNES_NORM_ALWAYS) { 784 if (i == snes->max_its) { 785 ierr = PetscInfo1(snes,"Maximum number of iterations has been reached: %D\n",snes->max_its);CHKERRQ(ierr); 786 if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT; 787 } 788 } else if (!snes->reason) snes->reason = SNES_CONVERGED_ITS; 789 PetscFunctionReturn(0); 790 } 791 792 /* -------------------------------------------------------------------------------------------*/ 793 794 /*MC 795 SNESCOMPOSITE - Build a preconditioner by composing together several nonlinear solvers 796 797 Options Database Keys: 798 + -snes_composite_type <type: one of multiplicative, additive, symmetric_multiplicative, special> - Sets composite preconditioner type 799 - -snes_composite_sneses - <snes0,snes1,...> list of SNESes to compose 800 801 Level: intermediate 802 803 Concepts: composing solvers 804 805 .seealso: SNESCreate(), SNESSetType(), SNESType (for list of available types), SNES, 806 SNESSHELL, SNESCompositeSetType(), SNESCompositeSpecialSetAlpha(), SNESCompositeAddSNES(), 807 SNESCompositeGetSNES() 808 809 M*/ 810 811 #undef __FUNCT__ 812 #define __FUNCT__ "SNESCreate_Composite" 813 PETSC_EXTERN PetscErrorCode SNESCreate_Composite(SNES snes) 814 { 815 PetscErrorCode ierr; 816 SNES_Composite *jac; 817 818 PetscFunctionBegin; 819 ierr = PetscNewLog(snes,SNES_Composite,&jac);CHKERRQ(ierr); 820 821 snes->ops->solve = SNESSolve_Composite; 822 snes->ops->setup = SNESSetUp_Composite; 823 snes->ops->reset = SNESReset_Composite; 824 snes->ops->destroy = SNESDestroy_Composite; 825 snes->ops->setfromoptions = SNESSetFromOptions_Composite; 826 snes->ops->view = SNESView_Composite; 827 828 snes->data = (void*)jac; 829 jac->type = SNES_COMPOSITE_ADDITIVEOPTIMAL; 830 jac->Fes = NULL; 831 jac->Xes = NULL; 832 jac->fnorms = NULL; 833 jac->nsnes = 0; 834 jac->head = 0; 835 jac->stol = 0.1; 836 jac->rtol = 1.1; 837 838 jac->h = NULL; 839 jac->s = NULL; 840 jac->beta = NULL; 841 jac->work = NULL; 842 jac->rwork = NULL; 843 844 ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESCompositeSetType_C",SNESCompositeSetType_Composite);CHKERRQ(ierr); 845 ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESCompositeAddSNES_C",SNESCompositeAddSNES_Composite);CHKERRQ(ierr); 846 ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESCompositeGetSNES_C",SNESCompositeGetSNES_Composite);CHKERRQ(ierr); 847 ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESCompositeSetDamping_C",SNESCompositeSetDamping_Composite);CHKERRQ(ierr); 848 PetscFunctionReturn(0); 849 } 850 851