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