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