1 #include <petsc-private/snesimpl.h> /*I "petscsnes.h" I*/ 2 #include <petscdm.h> 3 4 typedef struct { 5 PetscInt n; /* local subdomains */ 6 SNES *subsnes; /* nonlinear solvers for each subdomain */ 7 Vec *x; /* solution vectors */ 8 Vec *xl; /* solution local vectors */ 9 Vec *y; /* step vectors */ 10 Vec *b; /* rhs vectors */ 11 VecScatter *oscatter; /* scatter from global space to the subdomain global space */ 12 VecScatter *iscatter; /* scatter from global space to the nonoverlapping subdomain space */ 13 VecScatter *gscatter; /* scatter from global space to the subdomain local space */ 14 PCASMType type; /* ASM type */ 15 PetscBool usesdm; /* use the DM for setting up the subproblems */ 16 PetscBool finaljacobian; /* compute the jacobian of the converged solution */ 17 18 /* logging events */ 19 PetscLogEvent eventrestrictinterp; 20 PetscLogEvent eventsubsolve; 21 22 PetscInt fjtype; /* type of computed jacobian */ 23 Vec xinit; /* initial solution in case the final jacobian type is computed as first */ 24 } SNES_NASM; 25 26 const char *const SNESNASMTypes[] = {"NONE","RESTRICT","INTERPOLATE","BASIC","PCASMType","PC_ASM_",0}; 27 const char *const SNESNASMFJTypes[] = {"FINALOUTER","FINALINNER","INITIAL"}; 28 29 #undef __FUNCT__ 30 #define __FUNCT__ "SNESReset_NASM" 31 PetscErrorCode SNESReset_NASM(SNES snes) 32 { 33 SNES_NASM *nasm = (SNES_NASM*)snes->data; 34 PetscErrorCode ierr; 35 PetscInt i; 36 37 PetscFunctionBegin; 38 for (i=0; i<nasm->n; i++) { 39 if (nasm->xl) { ierr = VecDestroy(&nasm->xl[i]);CHKERRQ(ierr); } 40 if (nasm->x) { ierr = VecDestroy(&nasm->x[i]);CHKERRQ(ierr); } 41 if (nasm->y) { ierr = VecDestroy(&nasm->y[i]);CHKERRQ(ierr); } 42 if (nasm->b) { ierr = VecDestroy(&nasm->b[i]);CHKERRQ(ierr); } 43 44 if (nasm->subsnes) { ierr = SNESDestroy(&nasm->subsnes[i]);CHKERRQ(ierr); } 45 if (nasm->oscatter) { ierr = VecScatterDestroy(&nasm->oscatter[i]);CHKERRQ(ierr); } 46 if (nasm->iscatter) { ierr = VecScatterDestroy(&nasm->iscatter[i]);CHKERRQ(ierr); } 47 if (nasm->gscatter) { ierr = VecScatterDestroy(&nasm->gscatter[i]);CHKERRQ(ierr); } 48 } 49 50 if (nasm->x) {ierr = PetscFree(nasm->x);CHKERRQ(ierr);} 51 if (nasm->xl) {ierr = PetscFree(nasm->xl);CHKERRQ(ierr);} 52 if (nasm->y) {ierr = PetscFree(nasm->y);CHKERRQ(ierr);} 53 if (nasm->b) {ierr = PetscFree(nasm->b);CHKERRQ(ierr);} 54 55 if (nasm->xinit) {ierr = VecDestroy(&nasm->xinit);CHKERRQ(ierr);} 56 57 if (nasm->subsnes) {ierr = PetscFree(nasm->subsnes);CHKERRQ(ierr);} 58 if (nasm->oscatter) {ierr = PetscFree(nasm->oscatter);CHKERRQ(ierr);} 59 if (nasm->iscatter) {ierr = PetscFree(nasm->iscatter);CHKERRQ(ierr);} 60 if (nasm->gscatter) {ierr = PetscFree(nasm->gscatter);CHKERRQ(ierr);} 61 62 nasm->eventrestrictinterp = 0; 63 nasm->eventsubsolve = 0; 64 PetscFunctionReturn(0); 65 } 66 67 #undef __FUNCT__ 68 #define __FUNCT__ "SNESDestroy_NASM" 69 PetscErrorCode SNESDestroy_NASM(SNES snes) 70 { 71 PetscErrorCode ierr; 72 73 PetscFunctionBegin; 74 ierr = SNESReset_NASM(snes);CHKERRQ(ierr); 75 ierr = PetscFree(snes->data);CHKERRQ(ierr); 76 PetscFunctionReturn(0); 77 } 78 79 #undef __FUNCT__ 80 #define __FUNCT__ "DMGlobalToLocalSubDomainDirichletHook_Private" 81 PetscErrorCode DMGlobalToLocalSubDomainDirichletHook_Private(DM dm,Vec g,InsertMode mode,Vec l,void *ctx) 82 { 83 PetscErrorCode ierr; 84 Vec bcs = (Vec)ctx; 85 86 PetscFunctionBegin; 87 ierr = VecCopy(bcs,l);CHKERRQ(ierr); 88 PetscFunctionReturn(0); 89 } 90 91 #undef __FUNCT__ 92 #define __FUNCT__ "SNESSetUp_NASM" 93 PetscErrorCode SNESSetUp_NASM(SNES snes) 94 { 95 SNES_NASM *nasm = (SNES_NASM*)snes->data; 96 PetscErrorCode ierr; 97 DM dm,subdm; 98 DM *subdms; 99 PetscInt i; 100 const char *optionsprefix; 101 Vec F; 102 PetscMPIInt size; 103 KSP ksp; 104 PC pc; 105 106 PetscFunctionBegin; 107 if (!nasm->subsnes) { 108 ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 109 if (dm) { 110 nasm->usesdm = PETSC_TRUE; 111 ierr = DMCreateDomainDecomposition(dm,&nasm->n,NULL,NULL,NULL,&subdms);CHKERRQ(ierr); 112 if (!subdms) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"DM has no default decomposition defined. Set subsolves manually with SNESNASMSetSubdomains()."); 113 ierr = DMCreateDomainDecompositionScatters(dm,nasm->n,subdms,&nasm->iscatter,&nasm->oscatter,&nasm->gscatter);CHKERRQ(ierr); 114 115 ierr = SNESGetOptionsPrefix(snes, &optionsprefix);CHKERRQ(ierr); 116 ierr = PetscMalloc(nasm->n*sizeof(SNES),&nasm->subsnes);CHKERRQ(ierr); 117 for (i=0; i<nasm->n; i++) { 118 ierr = SNESCreate(PETSC_COMM_SELF,&nasm->subsnes[i]);CHKERRQ(ierr); 119 ierr = SNESAppendOptionsPrefix(nasm->subsnes[i],optionsprefix);CHKERRQ(ierr); 120 ierr = SNESAppendOptionsPrefix(nasm->subsnes[i],"sub_");CHKERRQ(ierr); 121 ierr = SNESSetDM(nasm->subsnes[i],subdms[i]);CHKERRQ(ierr); 122 ierr = MPI_Comm_size(PetscObjectComm((PetscObject)nasm->subsnes[i]),&size);CHKERRQ(ierr); 123 if (size == 1) { 124 ierr = SNESGetKSP(nasm->subsnes[i],&ksp);CHKERRQ(ierr); 125 ierr = KSPGetPC(ksp,&pc);CHKERRQ(ierr); 126 ierr = KSPSetType(ksp,KSPPREONLY);CHKERRQ(ierr); 127 ierr = PCSetType(pc,PCLU);CHKERRQ(ierr); 128 } 129 ierr = SNESSetFromOptions(nasm->subsnes[i]);CHKERRQ(ierr); 130 ierr = DMDestroy(&subdms[i]);CHKERRQ(ierr); 131 } 132 ierr = PetscFree(subdms);CHKERRQ(ierr); 133 } else SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"Cannot construct local problems automatically without a DM!"); 134 } else SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"Must set subproblems manually if there is no DM!"); 135 /* allocate the global vectors */ 136 if (!nasm->x) { 137 ierr = PetscMalloc(nasm->n*sizeof(Vec),&nasm->x);CHKERRQ(ierr); 138 ierr = PetscMemzero(nasm->x,nasm->n*sizeof(Vec));CHKERRQ(ierr); 139 } 140 if (!nasm->xl) { 141 ierr = PetscMalloc(nasm->n*sizeof(Vec),&nasm->xl);CHKERRQ(ierr); 142 ierr = PetscMemzero(nasm->xl,nasm->n*sizeof(Vec));CHKERRQ(ierr); 143 } 144 if (!nasm->y) { 145 ierr = PetscMalloc(nasm->n*sizeof(Vec),&nasm->y);CHKERRQ(ierr); 146 ierr = PetscMemzero(nasm->y,nasm->n*sizeof(Vec));CHKERRQ(ierr); 147 } 148 if (!nasm->b) { 149 ierr = PetscMalloc(nasm->n*sizeof(Vec),&nasm->b);CHKERRQ(ierr); 150 ierr = PetscMemzero(nasm->b,nasm->n*sizeof(Vec));CHKERRQ(ierr); 151 } 152 153 for (i=0; i<nasm->n; i++) { 154 ierr = SNESGetFunction(nasm->subsnes[i],&F,NULL,NULL);CHKERRQ(ierr); 155 if (!nasm->x[i]) {ierr = VecDuplicate(F,&nasm->x[i]);CHKERRQ(ierr);} 156 if (!nasm->y[i]) {ierr = VecDuplicate(F,&nasm->y[i]);CHKERRQ(ierr);} 157 if (!nasm->b[i]) {ierr = VecDuplicate(F,&nasm->b[i]);CHKERRQ(ierr);} 158 if (!nasm->xl[i]) { 159 ierr = SNESGetDM(nasm->subsnes[i],&subdm);CHKERRQ(ierr); 160 ierr = DMCreateLocalVector(subdm,&nasm->xl[i]);CHKERRQ(ierr); 161 } 162 ierr = DMGlobalToLocalHookAdd(subdm,DMGlobalToLocalSubDomainDirichletHook_Private,NULL,nasm->xl[i]);CHKERRQ(ierr); 163 } 164 if (nasm->finaljacobian) { 165 ierr = SNESSetUpMatrices(snes);CHKERRQ(ierr); 166 if (nasm->fjtype == 2) { 167 ierr = VecDuplicate(snes->vec_sol,&nasm->xinit);CHKERRQ(ierr); 168 } 169 for (i=0; i<nasm->n;i++) { 170 ierr = SNESSetUpMatrices(nasm->subsnes[i]);CHKERRQ(ierr); 171 } 172 } 173 PetscFunctionReturn(0); 174 } 175 176 #undef __FUNCT__ 177 #define __FUNCT__ "SNESSetFromOptions_NASM" 178 PetscErrorCode SNESSetFromOptions_NASM(SNES snes) 179 { 180 PetscErrorCode ierr; 181 PCASMType asmtype; 182 PetscBool flg,monflg; 183 SNES_NASM *nasm = (SNES_NASM*)snes->data; 184 185 PetscFunctionBegin; 186 ierr = PetscOptionsHead("Nonlinear Additive Schwartz options");CHKERRQ(ierr); 187 ierr = PetscOptionsEnum("-snes_nasm_type","Type of restriction/extension","",SNESNASMTypes,(PetscEnum)nasm->type,(PetscEnum*)&asmtype,&flg);CHKERRQ(ierr); 188 if (flg) nasm->type = asmtype; 189 flg = PETSC_FALSE; 190 monflg = PETSC_TRUE; 191 ierr = PetscOptionsBool("-snes_nasm_log","Log times for subSNES solves and restriction","",monflg,&monflg,&flg);CHKERRQ(ierr); 192 ierr = PetscOptionsBool("-snes_nasm_finaljacobian","Compute the global jacobian of the final iterate (for ASPIN)","",nasm->finaljacobian,&nasm->finaljacobian,NULL);CHKERRQ(ierr); 193 ierr = PetscOptionsEList("-snes_nasm_finaljacobian_type","The type of the final jacobian computed.","",SNESNASMFJTypes,3,SNESNASMFJTypes[0],&nasm->fjtype,NULL);CHKERRQ(ierr); 194 if (flg) { 195 ierr = PetscLogEventRegister("SNESNASMSubSolve",((PetscObject)snes)->classid,&nasm->eventsubsolve);CHKERRQ(ierr); 196 ierr = PetscLogEventRegister("SNESNASMRestrict",((PetscObject)snes)->classid,&nasm->eventrestrictinterp);CHKERRQ(ierr); 197 } 198 ierr = PetscOptionsTail();CHKERRQ(ierr); 199 PetscFunctionReturn(0); 200 } 201 202 #undef __FUNCT__ 203 #define __FUNCT__ "SNESView_NASM" 204 PetscErrorCode SNESView_NASM(SNES snes, PetscViewer viewer) 205 { 206 SNES_NASM *nasm = (SNES_NASM*)snes->data; 207 PetscErrorCode ierr; 208 PetscMPIInt rank; 209 PetscInt i,N; 210 PetscBool iascii,isstring; 211 PetscViewer sviewer; 212 MPI_Comm comm; 213 214 PetscFunctionBegin; 215 ierr = PetscObjectGetComm((PetscObject)snes,&comm);CHKERRQ(ierr); 216 ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 217 ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring);CHKERRQ(ierr); 218 ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 219 ierr = MPI_Reduce(&nasm->n,&N,1,MPIU_INT,MPIU_SUM,0,comm);CHKERRQ(ierr); 220 if (iascii) { 221 ierr = PetscViewerASCIIPrintf(viewer, " Nonlinear Additive Schwarz: total subdomain blocks = %D\n",N);CHKERRQ(ierr); 222 ierr = PetscViewerASCIIPrintf(viewer," Nonlinear Additive Schwarz: restriction/interpolation type - %s\n",SNESNASMTypes[nasm->type]);CHKERRQ(ierr); 223 ierr = PetscViewerASCIIPrintf(viewer," Nonlinear Additive Schwarz: subSNES iterations: %D subKSP iterations: %D\n",SNESNASMTypes[nasm->type]);CHKERRQ(ierr); 224 ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_TRUE);CHKERRQ(ierr); 225 ierr = PetscViewerASCIISynchronizedPrintf(viewer," [%d] number of local blocks = %D\n",(int)rank,nasm->n);CHKERRQ(ierr); 226 ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 227 ierr = PetscViewerASCIIPrintf(viewer," Local SNES objects:\n");CHKERRQ(ierr); 228 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 229 if (!rank) { 230 for (i=0; i<nasm->n; i++) { 231 ierr = PetscViewerGetSingleton(viewer,&sviewer);CHKERRQ(ierr); 232 ierr = SNESView(nasm->subsnes[i],sviewer);CHKERRQ(ierr); 233 ierr = PetscViewerRestoreSingleton(viewer,&sviewer);CHKERRQ(ierr); 234 } 235 } 236 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 237 ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 238 ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_FALSE);CHKERRQ(ierr); 239 } else if (isstring) { 240 ierr = PetscViewerStringSPrintf(viewer," blocks=%D,type=%s",N,SNESNASMTypes[nasm->type]);CHKERRQ(ierr); 241 ierr = PetscViewerGetSingleton(viewer,&sviewer);CHKERRQ(ierr); 242 if (nasm->subsnes && !rank) {ierr = SNESView(nasm->subsnes[0],sviewer);CHKERRQ(ierr);} 243 ierr = PetscViewerRestoreSingleton(viewer,&sviewer);CHKERRQ(ierr); 244 } 245 PetscFunctionReturn(0); 246 } 247 248 #undef __FUNCT__ 249 #define __FUNCT__ "SNESNASMSetSubdomains" 250 /*@ 251 SNESNASMSetSubdomains - Manually Set the context required to restrict and solve subdomain problems. 252 253 Not Collective 254 255 Input Parameters: 256 + SNES - the SNES context 257 . n - the number of local subdomains 258 . subsnes - solvers defined on the local subdomains 259 . iscatter - scatters into the nonoverlapping portions of the local subdomains 260 . oscatter - scatters into the overlapping portions of the local subdomains 261 - gscatter - scatters into the (ghosted) local vector of the local subdomain 262 263 Level: intermediate 264 265 .keywords: SNES, NASM 266 267 .seealso: SNESNASM, SNESNASMGetSubdomains() 268 @*/ 269 PetscErrorCode SNESNASMSetSubdomains(SNES snes,PetscInt n,SNES subsnes[],VecScatter iscatter[],VecScatter oscatter[],VecScatter gscatter[]) 270 { 271 PetscErrorCode ierr; 272 PetscErrorCode (*f)(SNES,PetscInt,SNES*,VecScatter*,VecScatter*,VecScatter*); 273 274 PetscFunctionBegin; 275 ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESNASMSetSubdomains_C",(void (**)(void))&f);CHKERRQ(ierr); 276 ierr = (f)(snes,n,subsnes,iscatter,oscatter,gscatter);CHKERRQ(ierr); 277 PetscFunctionReturn(0); 278 } 279 280 #undef __FUNCT__ 281 #define __FUNCT__ "SNESNASMSetSubdomains_NASM" 282 PetscErrorCode SNESNASMSetSubdomains_NASM(SNES snes,PetscInt n,SNES subsnes[],VecScatter iscatter[],VecScatter oscatter[],VecScatter gscatter[]) 283 { 284 PetscInt i; 285 PetscErrorCode ierr; 286 SNES_NASM *nasm = (SNES_NASM*)snes->data; 287 288 PetscFunctionBegin; 289 if (snes->setupcalled) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"SNESNASMSetSubdomains() should be called before calling SNESSetUp()."); 290 291 /* tear down the previously set things */ 292 ierr = SNESReset(snes);CHKERRQ(ierr); 293 294 nasm->n = n; 295 if (oscatter) { 296 for (i=0; i<n; i++) {ierr = PetscObjectReference((PetscObject)oscatter[i]);CHKERRQ(ierr);} 297 } 298 if (iscatter) { 299 for (i=0; i<n; i++) {ierr = PetscObjectReference((PetscObject)iscatter[i]);CHKERRQ(ierr);} 300 } 301 if (gscatter) { 302 for (i=0; i<n; i++) {ierr = PetscObjectReference((PetscObject)gscatter[i]);CHKERRQ(ierr);} 303 } 304 if (oscatter) { 305 ierr = PetscMalloc(n*sizeof(IS),&nasm->oscatter);CHKERRQ(ierr); 306 for (i=0; i<n; i++) { 307 nasm->oscatter[i] = oscatter[i]; 308 } 309 } 310 if (iscatter) { 311 ierr = PetscMalloc(n*sizeof(IS),&nasm->iscatter);CHKERRQ(ierr); 312 for (i=0; i<n; i++) { 313 nasm->iscatter[i] = iscatter[i]; 314 } 315 } 316 if (gscatter) { 317 ierr = PetscMalloc(n*sizeof(IS),&nasm->gscatter);CHKERRQ(ierr); 318 for (i=0; i<n; i++) { 319 nasm->gscatter[i] = gscatter[i]; 320 } 321 } 322 323 if (subsnes) { 324 ierr = PetscMalloc(n*sizeof(SNES),&nasm->subsnes);CHKERRQ(ierr); 325 for (i=0; i<n; i++) { 326 nasm->subsnes[i] = subsnes[i]; 327 } 328 } 329 PetscFunctionReturn(0); 330 } 331 332 #undef __FUNCT__ 333 #define __FUNCT__ "SNESNASMGetSubdomains" 334 /*@ 335 SNESNASMGetSubdomains - Get the local subdomain context. 336 337 Not Collective 338 339 Input Parameters: 340 . SNES - the SNES context 341 342 Output Parameters: 343 + n - the number of local subdomains 344 . subsnes - solvers defined on the local subdomains 345 . iscatter - scatters into the nonoverlapping portions of the local subdomains 346 . oscatter - scatters into the overlapping portions of the local subdomains 347 - gscatter - scatters into the (ghosted) local vector of the local subdomain 348 349 Level: intermediate 350 351 .keywords: SNES, NASM 352 353 .seealso: SNESNASM, SNESNASMSetSubdomains() 354 @*/ 355 PetscErrorCode SNESNASMGetSubdomains(SNES snes,PetscInt *n,SNES *subsnes[],VecScatter *iscatter[],VecScatter *oscatter[],VecScatter *gscatter[]) 356 { 357 PetscErrorCode ierr; 358 PetscErrorCode (*f)(SNES,PetscInt*,SNES**,VecScatter**,VecScatter**,VecScatter**); 359 360 PetscFunctionBegin; 361 ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESNASMGetSubdomains_C",(void (**)(void))&f);CHKERRQ(ierr); 362 ierr = (f)(snes,n,subsnes,iscatter,oscatter,gscatter);CHKERRQ(ierr); 363 PetscFunctionReturn(0); 364 } 365 366 #undef __FUNCT__ 367 #define __FUNCT__ "SNESNASMGetSubdomains_NASM" 368 PetscErrorCode SNESNASMGetSubdomains_NASM(SNES snes,PetscInt *n,SNES *subsnes[],VecScatter *iscatter[],VecScatter *oscatter[],VecScatter *gscatter[]) 369 { 370 SNES_NASM *nasm = (SNES_NASM*)snes->data; 371 372 PetscFunctionBegin; 373 if (n) *n = nasm->n; 374 if (oscatter) *oscatter = nasm->oscatter; 375 if (iscatter) *iscatter = nasm->iscatter; 376 if (gscatter) *gscatter = nasm->gscatter; 377 if (subsnes) *subsnes = nasm->subsnes; 378 PetscFunctionReturn(0); 379 } 380 381 #undef __FUNCT__ 382 #define __FUNCT__ "SNESNASMGetSubdomainVecs" 383 /*@ 384 SNESNASMGetSubdomainVecs - Get the processor-local subdomain vectors 385 386 Not Collective 387 388 Input Parameters: 389 . SNES - the SNES context 390 391 Output Parameters: 392 + n - the number of local subdomains 393 . x - The subdomain solution vector 394 . y - The subdomain step vector 395 . b - The subdomain RHS vector 396 - xl - The subdomain local vectors (ghosted) 397 398 Level: developer 399 400 .keywords: SNES, NASM 401 402 .seealso: SNESNASM, SNESNASMGetSubdomains() 403 @*/ 404 PetscErrorCode SNESNASMGetSubdomainVecs(SNES snes,PetscInt *n,Vec **x,Vec **y,Vec **b, Vec **xl) 405 { 406 PetscErrorCode ierr; 407 PetscErrorCode (*f)(SNES,PetscInt*,Vec**,Vec**,Vec**,Vec**); 408 409 PetscFunctionBegin; 410 ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESNASMGetSubdomainVecs_C",(void (**)(void))&f);CHKERRQ(ierr); 411 ierr = (f)(snes,n,x,y,b,xl);CHKERRQ(ierr); 412 PetscFunctionReturn(0); 413 } 414 415 #undef __FUNCT__ 416 #define __FUNCT__ "SNESNASMGetSubdomainVecs_NASM" 417 PetscErrorCode SNESNASMGetSubdomainVecs_NASM(SNES snes,PetscInt *n,Vec **x,Vec **y,Vec **b,Vec **xl) 418 { 419 SNES_NASM *nasm = (SNES_NASM*)snes->data; 420 421 PetscFunctionBegin; 422 if (n) *n = nasm->n; 423 if (x) *x = nasm->x; 424 if (y) *y = nasm->y; 425 if (b) *b = nasm->b; 426 if (xl) *xl = nasm->xl; 427 PetscFunctionReturn(0); 428 } 429 430 #undef __FUNCT__ 431 #define __FUNCT__ "SNESNASMSetComputeFinalJacobian" 432 /*@ 433 SNESNASMSetComputeFinalJacobian - Schedules the computation of the global and subdomain jacobians upon convergence 434 435 Collective on SNES 436 437 Input Parameters: 438 + SNES - the SNES context 439 - flg - indication of whether to compute the jacobians or not 440 441 Level: developer 442 443 Notes: This is used almost exclusively in the implementation of ASPIN, where the converged subdomain and global jacobian 444 is needed at each linear iteration. 445 446 .keywords: SNES, NASM, ASPIN 447 448 .seealso: SNESNASM, SNESNASMGetSubdomains() 449 @*/ 450 PetscErrorCode SNESNASMSetComputeFinalJacobian(SNES snes,PetscBool flg) 451 { 452 PetscErrorCode (*f)(SNES,PetscBool); 453 PetscErrorCode ierr; 454 455 PetscFunctionBegin; 456 ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESNASMSetComputeFinalJacobian_C",(void (**)(void))&f);CHKERRQ(ierr); 457 ierr = (f)(snes,flg);CHKERRQ(ierr); 458 PetscFunctionReturn(0); 459 } 460 461 #undef __FUNCT__ 462 #define __FUNCT__ "SNESNASMSetComputeFinalJacobian_NASM" 463 PetscErrorCode SNESNASMSetComputeFinalJacobian_NASM(SNES snes,PetscBool flg) 464 { 465 SNES_NASM *nasm = (SNES_NASM*)snes->data; 466 467 PetscFunctionBegin; 468 nasm->finaljacobian = flg; 469 if (flg) snes->usesksp = PETSC_TRUE; 470 PetscFunctionReturn(0); 471 } 472 473 #undef __FUNCT__ 474 #define __FUNCT__ "SNESNASMSolveLocal_Private" 475 PetscErrorCode SNESNASMSolveLocal_Private(SNES snes,Vec B,Vec Y,Vec X) 476 { 477 SNES_NASM *nasm = (SNES_NASM*)snes->data; 478 SNES subsnes; 479 PetscInt i; 480 PetscErrorCode ierr; 481 Vec Xlloc,Xl,Bl,Yl; 482 VecScatter iscat,oscat,gscat; 483 DM dm,subdm; 484 485 PetscFunctionBegin; 486 ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 487 ierr = VecSet(Y,0);CHKERRQ(ierr); 488 if (nasm->eventrestrictinterp) {ierr = PetscLogEventBegin(nasm->eventrestrictinterp,snes,0,0,0);CHKERRQ(ierr);} 489 for (i=0; i<nasm->n; i++) { 490 /* scatter the solution to the local solution */ 491 Xlloc = nasm->xl[i]; 492 gscat = nasm->gscatter[i]; 493 oscat = nasm->oscatter[i]; 494 ierr = VecScatterBegin(gscat,X,Xlloc,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 495 if (B) { 496 /* scatter the RHS to the local RHS */ 497 Bl = nasm->b[i]; 498 ierr = VecScatterBegin(oscat,B,Bl,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 499 } 500 } 501 if (nasm->eventrestrictinterp) {ierr = PetscLogEventEnd(nasm->eventrestrictinterp,snes,0,0,0);CHKERRQ(ierr);} 502 503 504 if (nasm->eventsubsolve) {ierr = PetscLogEventBegin(nasm->eventsubsolve,snes,0,0,0);CHKERRQ(ierr);} 505 for (i=0; i<nasm->n; i++) { 506 Xl = nasm->x[i]; 507 Xlloc = nasm->xl[i]; 508 Yl = nasm->y[i]; 509 subsnes = nasm->subsnes[i]; 510 ierr = SNESGetDM(subsnes,&subdm);CHKERRQ(ierr); 511 iscat = nasm->iscatter[i]; 512 oscat = nasm->oscatter[i]; 513 gscat = nasm->gscatter[i]; 514 ierr = VecScatterEnd(gscat,X,Xlloc,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 515 if (B) { 516 Bl = nasm->b[i]; 517 ierr = VecScatterEnd(oscat,B,Bl,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 518 } else Bl = NULL; 519 ierr = DMSubDomainRestrict(dm,oscat,gscat,subdm);CHKERRQ(ierr); 520 ierr = DMLocalToGlobalBegin(subdm,Xlloc,INSERT_VALUES,Xl);CHKERRQ(ierr); 521 ierr = DMLocalToGlobalEnd(subdm,Xlloc,INSERT_VALUES,Xl);CHKERRQ(ierr); 522 ierr = VecCopy(Xl,Yl);CHKERRQ(ierr); 523 ierr = SNESSolve(subsnes,Bl,Xl);CHKERRQ(ierr); 524 ierr = VecAYPX(Yl,-1.0,Xl);CHKERRQ(ierr); 525 if (nasm->type == PC_ASM_BASIC) { 526 ierr = VecScatterBegin(oscat,Yl,Y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 527 } else if (nasm->type == PC_ASM_RESTRICT) { 528 ierr = VecScatterBegin(iscat,Yl,Y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 529 } else SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"Only basic and restrict types are supported for SNESNASM"); 530 } 531 if (nasm->eventsubsolve) {ierr = PetscLogEventEnd(nasm->eventsubsolve,snes,0,0,0);CHKERRQ(ierr);} 532 if (nasm->eventrestrictinterp) {ierr = PetscLogEventBegin(nasm->eventrestrictinterp,snes,0,0,0);CHKERRQ(ierr);} 533 for (i=0; i<nasm->n; i++) { 534 Yl = nasm->y[i]; 535 iscat = nasm->iscatter[i]; 536 oscat = nasm->oscatter[i]; 537 if (nasm->type == PC_ASM_BASIC) { 538 ierr = VecScatterEnd(oscat,Yl,Y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 539 } else if (nasm->type == PC_ASM_RESTRICT) { 540 ierr = VecScatterEnd(iscat,Yl,Y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 541 } else SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"Only basic and restrict types are supported for SNESNASM"); 542 } 543 if (nasm->eventrestrictinterp) {ierr = PetscLogEventEnd(nasm->eventrestrictinterp,snes,0,0,0);CHKERRQ(ierr);} 544 ierr = VecAXPY(X,1.0,Y);CHKERRQ(ierr); 545 PetscFunctionReturn(0); 546 } 547 548 #undef __FUNCT__ 549 #define __FUNCT__ "SNESNASMComputeFinalJacobian_Private" 550 PetscErrorCode SNESNASMComputeFinalJacobian_Private(SNES snes, Vec Xfinal) 551 { 552 Vec X = Xfinal; 553 SNES_NASM *nasm = (SNES_NASM*)snes->data; 554 SNES subsnes; 555 PetscInt i,lag = 1; 556 PetscErrorCode ierr; 557 Vec Xlloc,Xl,Fl,F; 558 VecScatter oscat,gscat; 559 DM dm,subdm; 560 MatStructure flg = DIFFERENT_NONZERO_PATTERN; 561 PetscFunctionBegin; 562 if (nasm->fjtype == 2) X = nasm->xinit; 563 F = snes->vec_func; 564 if (snes->normtype == SNES_NORM_NONE) {ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);} 565 ierr = SNESComputeJacobian(snes,X,&snes->jacobian,&snes->jacobian_pre,&flg);CHKERRQ(ierr); 566 ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 567 if (nasm->eventrestrictinterp) {ierr = PetscLogEventBegin(nasm->eventrestrictinterp,snes,0,0,0);CHKERRQ(ierr);} 568 if (nasm->fjtype != 1) { 569 for (i=0; i<nasm->n; i++) { 570 Xlloc = nasm->xl[i]; 571 gscat = nasm->gscatter[i]; 572 oscat = nasm->oscatter[i]; 573 ierr = VecScatterBegin(gscat,X,Xlloc,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 574 } 575 } 576 if (nasm->eventrestrictinterp) {ierr = PetscLogEventEnd(nasm->eventrestrictinterp,snes,0,0,0);CHKERRQ(ierr);} 577 for (i=0; i<nasm->n; i++) { 578 Fl = nasm->subsnes[i]->vec_func; 579 Xl = nasm->x[i]; 580 Xlloc = nasm->xl[i]; 581 subsnes = nasm->subsnes[i]; 582 oscat = nasm->oscatter[i]; 583 gscat = nasm->gscatter[i]; 584 if (nasm->fjtype != 1) {ierr = VecScatterEnd(gscat,X,Xlloc,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);} 585 ierr = SNESGetDM(subsnes,&subdm);CHKERRQ(ierr); 586 ierr = DMSubDomainRestrict(dm,oscat,gscat,subdm);CHKERRQ(ierr); 587 if (nasm->fjtype != 1) { 588 ierr = DMLocalToGlobalBegin(subdm,Xlloc,INSERT_VALUES,Xl);CHKERRQ(ierr); 589 ierr = DMLocalToGlobalEnd(subdm,Xlloc,INSERT_VALUES,Xl);CHKERRQ(ierr); 590 } 591 if (subsnes->lagjacobian == -1) subsnes->lagjacobian = -2; 592 else if (subsnes->lagjacobian > 1) lag = subsnes->lagjacobian; 593 ierr = SNESComputeFunction(subsnes,Xl,Fl);CHKERRQ(ierr); 594 ierr = SNESComputeJacobian(subsnes,Xl,&subsnes->jacobian,&subsnes->jacobian_pre,&flg);CHKERRQ(ierr); 595 if (lag > 1) subsnes->lagjacobian = lag; 596 ierr = KSPSetOperators(subsnes->ksp,subsnes->jacobian,subsnes->jacobian_pre,flg);CHKERRQ(ierr); 597 } 598 PetscFunctionReturn(0); 599 } 600 601 #undef __FUNCT__ 602 #define __FUNCT__ "SNESSolve_NASM" 603 PetscErrorCode SNESSolve_NASM(SNES snes) 604 { 605 Vec F; 606 Vec X; 607 Vec B; 608 Vec Y; 609 PetscInt i; 610 PetscReal fnorm = 0.0; 611 PetscErrorCode ierr; 612 SNESNormType normtype; 613 SNES_NASM *nasm = (SNES_NASM*)snes->data; 614 615 PetscFunctionBegin; 616 X = snes->vec_sol; 617 Y = snes->vec_sol_update; 618 F = snes->vec_func; 619 B = snes->vec_rhs; 620 621 ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 622 snes->iter = 0; 623 snes->norm = 0.; 624 ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 625 snes->reason = SNES_CONVERGED_ITERATING; 626 ierr = SNESGetNormType(snes, &normtype);CHKERRQ(ierr); 627 if (normtype == SNES_NORM_FUNCTION || normtype == SNES_NORM_INITIAL_ONLY || normtype == SNES_NORM_INITIAL_FINAL_ONLY) { 628 /* compute the initial function and preconditioned update delX */ 629 if (!snes->vec_func_init_set) { 630 ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 631 if (snes->domainerror) { 632 snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 633 PetscFunctionReturn(0); 634 } 635 } else snes->vec_func_init_set = PETSC_FALSE; 636 637 /* convergence test */ 638 if (!snes->norm_init_set) { 639 ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F|| */ 640 if (PetscIsInfOrNanReal(fnorm)) { 641 snes->reason = SNES_DIVERGED_FNORM_NAN; 642 PetscFunctionReturn(0); 643 } 644 } else { 645 fnorm = snes->norm_init; 646 snes->norm_init_set = PETSC_FALSE; 647 } 648 ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 649 snes->iter = 0; 650 snes->norm = fnorm; 651 ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 652 ierr = SNESLogConvergenceHistory(snes,snes->norm,0);CHKERRQ(ierr); 653 ierr = SNESMonitor(snes,0,snes->norm);CHKERRQ(ierr); 654 655 /* set parameter for default relative tolerance convergence test */ 656 snes->ttol = fnorm*snes->rtol; 657 658 /* test convergence */ 659 ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 660 if (snes->reason) PetscFunctionReturn(0); 661 } else { 662 ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 663 ierr = SNESLogConvergenceHistory(snes,snes->norm,0);CHKERRQ(ierr); 664 ierr = SNESMonitor(snes,0,snes->norm);CHKERRQ(ierr); 665 } 666 667 /* Call general purpose update function */ 668 if (snes->ops->update) { 669 ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 670 } 671 /* copy the initial solution over for later */ 672 if (nasm->fjtype == 2) {ierr = VecCopy(X,nasm->xinit);CHKERRQ(ierr);} 673 674 for (i = 0; i < snes->max_its; i++) { 675 ierr = SNESNASMSolveLocal_Private(snes,B,Y,X);CHKERRQ(ierr); 676 if (normtype == SNES_NORM_FUNCTION || ((i == snes->max_its - 1) && (normtype == SNES_NORM_INITIAL_FINAL_ONLY || normtype == SNES_NORM_FINAL_ONLY))) { 677 ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 678 if (snes->domainerror) { 679 snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 680 break; 681 } 682 ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F|| */ 683 if (PetscIsInfOrNanReal(fnorm)) { 684 snes->reason = SNES_DIVERGED_FNORM_NAN; 685 break; 686 } 687 } 688 /* Monitor convergence */ 689 ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 690 snes->iter = i+1; 691 snes->norm = fnorm; 692 ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 693 ierr = SNESLogConvergenceHistory(snes,snes->norm,0);CHKERRQ(ierr); 694 ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr); 695 /* Test for convergence */ 696 if (normtype == SNES_NORM_FUNCTION) {ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);} 697 if (snes->reason) break; 698 /* Call general purpose update function */ 699 if (snes->ops->update) {ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);} 700 } 701 if (nasm->finaljacobian) {ierr = SNESNASMComputeFinalJacobian_Private(snes,X);CHKERRQ(ierr);} 702 if (normtype == SNES_NORM_FUNCTION) { 703 if (i == snes->max_its) { 704 ierr = PetscInfo1(snes,"Maximum number of iterations has been reached: %D\n",snes->max_its);CHKERRQ(ierr); 705 if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT; 706 } 707 } else if (!snes->reason) snes->reason = SNES_CONVERGED_ITS; /* NASM is meant to be used as a preconditioner */ 708 PetscFunctionReturn(0); 709 } 710 711 /*MC 712 SNESNASM - Nonlinear Additive Schwartz 713 714 Options Database: 715 + -snes_nasm_log - enable logging events for the communication and solve stages 716 . -snes_nasm_type <basic,restrict> - type of subdomain update used 717 . -snes_nasm_finaljacobian - compute the local and global jacobians of the final iterate 718 . -snes_nasm_finaljacobian_type <finalinner,finalouter,initial> pick state the jacobian is calculated at 719 . -sub_snes_ - options prefix of the subdomain nonlinear solves 720 . -sub_ksp_ - options prefix of the subdomain Krylov solver 721 - -sub_pc_ - options prefix of the subdomain preconditioner 722 723 Level: advanced 724 725 .seealso: SNESCreate(), SNES, SNESSetType(), SNESType (for list of available types) 726 M*/ 727 728 #undef __FUNCT__ 729 #define __FUNCT__ "SNESCreate_NASM" 730 PETSC_EXTERN PetscErrorCode SNESCreate_NASM(SNES snes) 731 { 732 SNES_NASM *nasm; 733 PetscErrorCode ierr; 734 735 PetscFunctionBegin; 736 ierr = PetscNewLog(snes, SNES_NASM, &nasm);CHKERRQ(ierr); 737 snes->data = (void*)nasm; 738 739 nasm->n = PETSC_DECIDE; 740 nasm->subsnes = 0; 741 nasm->x = 0; 742 nasm->xl = 0; 743 nasm->y = 0; 744 nasm->b = 0; 745 nasm->oscatter = 0; 746 nasm->iscatter = 0; 747 nasm->gscatter = 0; 748 749 nasm->type = PC_ASM_BASIC; 750 nasm->finaljacobian = PETSC_FALSE; 751 752 snes->ops->destroy = SNESDestroy_NASM; 753 snes->ops->setup = SNESSetUp_NASM; 754 snes->ops->setfromoptions = SNESSetFromOptions_NASM; 755 snes->ops->view = SNESView_NASM; 756 snes->ops->solve = SNESSolve_NASM; 757 snes->ops->reset = SNESReset_NASM; 758 759 snes->usesksp = PETSC_FALSE; 760 snes->usespc = PETSC_FALSE; 761 762 nasm->fjtype = 0; 763 nasm->xinit = NULL; 764 nasm->eventrestrictinterp = 0; 765 nasm->eventsubsolve = 0; 766 767 if (!snes->tolerancesset) { 768 snes->max_its = 10000; 769 snes->max_funcs = 10000; 770 } 771 772 ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESNASMSetSubdomains_C","SNESNASMSetSubdomains_NASM",SNESNASMSetSubdomains_NASM);CHKERRQ(ierr); 773 ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESNASMGetSubdomains_C","SNESNASMGetSubdomains_NASM",SNESNASMGetSubdomains_NASM);CHKERRQ(ierr); 774 ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESNASMGetSubdomainVecs_C","SNESNASMGetSubdomainVecs_NASM",SNESNASMGetSubdomainVecs_NASM);CHKERRQ(ierr); 775 ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESNASMSetComputeFinalJacobian_C","SNESNASMSetComputeFinalJacobian_NASM",SNESNASMSetComputeFinalJacobian_NASM);CHKERRQ(ierr); 776 PetscFunctionReturn(0); 777 } 778 779