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