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 PetscReal damping; /* damping parameter for updates from the blocks */ 18 PetscBool same_local_solves; /* flag to determine if the solvers have been individually modified */ 19 20 /* logging events */ 21 PetscLogEvent eventrestrictinterp; 22 PetscLogEvent eventsubsolve; 23 24 PetscInt fjtype; /* type of computed jacobian */ 25 Vec xinit; /* initial solution in case the final jacobian type is computed as first */ 26 } SNES_NASM; 27 28 const char *const SNESNASMTypes[] = {"NONE","RESTRICT","INTERPOLATE","BASIC","PCASMType","PC_ASM_",0}; 29 const char *const SNESNASMFJTypes[] = {"FINALOUTER","FINALINNER","INITIAL"}; 30 31 #undef __FUNCT__ 32 #define __FUNCT__ "SNESReset_NASM" 33 static PetscErrorCode SNESReset_NASM(SNES snes) 34 { 35 SNES_NASM *nasm = (SNES_NASM*)snes->data; 36 PetscErrorCode ierr; 37 PetscInt i; 38 39 PetscFunctionBegin; 40 for (i=0; i<nasm->n; i++) { 41 if (nasm->xl) { ierr = VecDestroy(&nasm->xl[i]);CHKERRQ(ierr); } 42 if (nasm->x) { ierr = VecDestroy(&nasm->x[i]);CHKERRQ(ierr); } 43 if (nasm->y) { ierr = VecDestroy(&nasm->y[i]);CHKERRQ(ierr); } 44 if (nasm->b) { ierr = VecDestroy(&nasm->b[i]);CHKERRQ(ierr); } 45 46 if (nasm->subsnes) { ierr = SNESDestroy(&nasm->subsnes[i]);CHKERRQ(ierr); } 47 if (nasm->oscatter) { ierr = VecScatterDestroy(&nasm->oscatter[i]);CHKERRQ(ierr); } 48 if (nasm->iscatter) { ierr = VecScatterDestroy(&nasm->iscatter[i]);CHKERRQ(ierr); } 49 if (nasm->gscatter) { ierr = VecScatterDestroy(&nasm->gscatter[i]);CHKERRQ(ierr); } 50 } 51 52 if (nasm->x) {ierr = PetscFree(nasm->x);CHKERRQ(ierr);} 53 if (nasm->xl) {ierr = PetscFree(nasm->xl);CHKERRQ(ierr);} 54 if (nasm->y) {ierr = PetscFree(nasm->y);CHKERRQ(ierr);} 55 if (nasm->b) {ierr = PetscFree(nasm->b);CHKERRQ(ierr);} 56 57 if (nasm->xinit) {ierr = VecDestroy(&nasm->xinit);CHKERRQ(ierr);} 58 59 if (nasm->subsnes) {ierr = PetscFree(nasm->subsnes);CHKERRQ(ierr);} 60 if (nasm->oscatter) {ierr = PetscFree(nasm->oscatter);CHKERRQ(ierr);} 61 if (nasm->iscatter) {ierr = PetscFree(nasm->iscatter);CHKERRQ(ierr);} 62 if (nasm->gscatter) {ierr = PetscFree(nasm->gscatter);CHKERRQ(ierr);} 63 64 nasm->eventrestrictinterp = 0; 65 nasm->eventsubsolve = 0; 66 PetscFunctionReturn(0); 67 } 68 69 #undef __FUNCT__ 70 #define __FUNCT__ "SNESDestroy_NASM" 71 static PetscErrorCode SNESDestroy_NASM(SNES snes) 72 { 73 PetscErrorCode ierr; 74 75 PetscFunctionBegin; 76 ierr = SNESReset_NASM(snes);CHKERRQ(ierr); 77 ierr = PetscFree(snes->data);CHKERRQ(ierr); 78 PetscFunctionReturn(0); 79 } 80 81 #undef __FUNCT__ 82 #define __FUNCT__ "DMGlobalToLocalSubDomainDirichletHook_Private" 83 static PetscErrorCode DMGlobalToLocalSubDomainDirichletHook_Private(DM dm,Vec g,InsertMode mode,Vec l,void *ctx) 84 { 85 PetscErrorCode ierr; 86 Vec bcs = (Vec)ctx; 87 88 PetscFunctionBegin; 89 ierr = VecCopy(bcs,l);CHKERRQ(ierr); 90 PetscFunctionReturn(0); 91 } 92 93 #undef __FUNCT__ 94 #define __FUNCT__ "SNESSetUp_NASM" 95 static PetscErrorCode SNESSetUp_NASM(SNES snes) 96 { 97 SNES_NASM *nasm = (SNES_NASM*)snes->data; 98 PetscErrorCode ierr; 99 DM dm,subdm; 100 DM *subdms; 101 PetscInt i; 102 const char *optionsprefix; 103 Vec F; 104 PetscMPIInt size; 105 KSP ksp; 106 PC pc; 107 108 PetscFunctionBegin; 109 if (!nasm->subsnes) { 110 ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 111 if (dm) { 112 nasm->usesdm = PETSC_TRUE; 113 ierr = DMCreateDomainDecomposition(dm,&nasm->n,NULL,NULL,NULL,&subdms);CHKERRQ(ierr); 114 if (!subdms) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"DM has no default decomposition defined. Set subsolves manually with SNESNASMSetSubdomains()."); 115 ierr = DMCreateDomainDecompositionScatters(dm,nasm->n,subdms,&nasm->iscatter,&nasm->oscatter,&nasm->gscatter);CHKERRQ(ierr); 116 117 ierr = SNESGetOptionsPrefix(snes, &optionsprefix);CHKERRQ(ierr); 118 ierr = PetscMalloc1(nasm->n,&nasm->subsnes);CHKERRQ(ierr); 119 for (i=0; i<nasm->n; i++) { 120 ierr = SNESCreate(PETSC_COMM_SELF,&nasm->subsnes[i]);CHKERRQ(ierr); 121 ierr = SNESAppendOptionsPrefix(nasm->subsnes[i],optionsprefix);CHKERRQ(ierr); 122 ierr = SNESAppendOptionsPrefix(nasm->subsnes[i],"sub_");CHKERRQ(ierr); 123 ierr = SNESSetDM(nasm->subsnes[i],subdms[i]);CHKERRQ(ierr); 124 ierr = MPI_Comm_size(PetscObjectComm((PetscObject)nasm->subsnes[i]),&size);CHKERRQ(ierr); 125 if (size == 1) { 126 ierr = SNESGetKSP(nasm->subsnes[i],&ksp);CHKERRQ(ierr); 127 ierr = KSPGetPC(ksp,&pc);CHKERRQ(ierr); 128 ierr = KSPSetType(ksp,KSPPREONLY);CHKERRQ(ierr); 129 ierr = PCSetType(pc,PCLU);CHKERRQ(ierr); 130 } 131 ierr = SNESSetFromOptions(nasm->subsnes[i]);CHKERRQ(ierr); 132 ierr = DMDestroy(&subdms[i]);CHKERRQ(ierr); 133 } 134 ierr = PetscFree(subdms);CHKERRQ(ierr); 135 } else SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"Cannot construct local problems automatically without a DM!"); 136 } else SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"Must set subproblems manually if there is no DM!"); 137 /* allocate the global vectors */ 138 if (!nasm->x) { 139 ierr = PetscCalloc1(nasm->n,&nasm->x);CHKERRQ(ierr); 140 } 141 if (!nasm->xl) { 142 ierr = PetscCalloc1(nasm->n,&nasm->xl);CHKERRQ(ierr); 143 } 144 if (!nasm->y) { 145 ierr = PetscCalloc1(nasm->n,&nasm->y);CHKERRQ(ierr); 146 } 147 if (!nasm->b) { 148 ierr = PetscCalloc1(nasm->n,&nasm->b);CHKERRQ(ierr); 149 } 150 151 for (i=0; i<nasm->n; i++) { 152 ierr = SNESGetFunction(nasm->subsnes[i],&F,NULL,NULL);CHKERRQ(ierr); 153 if (!nasm->x[i]) {ierr = VecDuplicate(F,&nasm->x[i]);CHKERRQ(ierr);} 154 if (!nasm->y[i]) {ierr = VecDuplicate(F,&nasm->y[i]);CHKERRQ(ierr);} 155 if (!nasm->b[i]) {ierr = VecDuplicate(F,&nasm->b[i]);CHKERRQ(ierr);} 156 if (!nasm->xl[i]) { 157 ierr = SNESGetDM(nasm->subsnes[i],&subdm);CHKERRQ(ierr); 158 ierr = DMCreateLocalVector(subdm,&nasm->xl[i]);CHKERRQ(ierr); 159 ierr = DMGlobalToLocalHookAdd(subdm,DMGlobalToLocalSubDomainDirichletHook_Private,NULL,nasm->xl[i]);CHKERRQ(ierr); 160 } 161 } 162 if (nasm->finaljacobian) { 163 ierr = SNESSetUpMatrices(snes);CHKERRQ(ierr); 164 if (nasm->fjtype == 2) { 165 ierr = VecDuplicate(snes->vec_sol,&nasm->xinit);CHKERRQ(ierr); 166 } 167 for (i=0; i<nasm->n;i++) { 168 ierr = SNESSetUpMatrices(nasm->subsnes[i]);CHKERRQ(ierr); 169 } 170 } 171 PetscFunctionReturn(0); 172 } 173 174 #undef __FUNCT__ 175 #define __FUNCT__ "SNESSetFromOptions_NASM" 176 static PetscErrorCode SNESSetFromOptions_NASM(PetscOptionItems *PetscOptionsObject,SNES snes) 177 { 178 PetscErrorCode ierr; 179 PCASMType asmtype; 180 PetscBool flg,monflg,subviewflg; 181 SNES_NASM *nasm = (SNES_NASM*)snes->data; 182 183 PetscFunctionBegin; 184 ierr = PetscOptionsHead(PetscOptionsObject,"Nonlinear Additive Schwartz options");CHKERRQ(ierr); 185 ierr = PetscOptionsEnum("-snes_nasm_type","Type of restriction/extension","",SNESNASMTypes,(PetscEnum)nasm->type,(PetscEnum*)&asmtype,&flg);CHKERRQ(ierr); 186 if (flg) {ierr = SNESNASMSetType(snes,asmtype);CHKERRQ(ierr);} 187 flg = PETSC_FALSE; 188 monflg = PETSC_TRUE; 189 ierr = PetscOptionsReal("-snes_nasm_damping","Log times for subSNES solves and restriction","SNESNASMSetDamping",nasm->damping,&nasm->damping,&flg);CHKERRQ(ierr); 190 if (flg) {ierr = SNESNASMSetDamping(snes,nasm->damping);CHKERRQ(ierr);} 191 subviewflg = PETSC_FALSE; 192 ierr = PetscOptionsBool("-snes_nasm_sub_view","Print detailed information for every processor when using -snes_view","",subviewflg,&subviewflg,&flg);CHKERRQ(ierr); 193 if (flg) { 194 nasm->same_local_solves = PETSC_FALSE; 195 if (!subviewflg) { 196 nasm->same_local_solves = PETSC_TRUE; 197 } 198 } 199 ierr = PetscOptionsBool("-snes_nasm_finaljacobian","Compute the global jacobian of the final iterate (for ASPIN)","",nasm->finaljacobian,&nasm->finaljacobian,NULL);CHKERRQ(ierr); 200 ierr = PetscOptionsEList("-snes_nasm_finaljacobian_type","The type of the final jacobian computed.","",SNESNASMFJTypes,3,SNESNASMFJTypes[0],&nasm->fjtype,NULL);CHKERRQ(ierr); 201 ierr = PetscOptionsBool("-snes_nasm_log","Log times for subSNES solves and restriction","",monflg,&monflg,&flg);CHKERRQ(ierr); 202 if (flg) { 203 ierr = PetscLogEventRegister("SNESNASMSubSolve",((PetscObject)snes)->classid,&nasm->eventsubsolve);CHKERRQ(ierr); 204 ierr = PetscLogEventRegister("SNESNASMRestrict",((PetscObject)snes)->classid,&nasm->eventrestrictinterp);CHKERRQ(ierr); 205 } 206 ierr = PetscOptionsTail();CHKERRQ(ierr); 207 PetscFunctionReturn(0); 208 } 209 210 #undef __FUNCT__ 211 #define __FUNCT__ "SNESView_NASM" 212 static PetscErrorCode SNESView_NASM(SNES snes, PetscViewer viewer) 213 { 214 SNES_NASM *nasm = (SNES_NASM*)snes->data; 215 PetscErrorCode ierr; 216 PetscMPIInt rank,size; 217 PetscInt i,N,bsz; 218 PetscBool iascii,isstring; 219 PetscViewer sviewer; 220 MPI_Comm comm; 221 222 PetscFunctionBegin; 223 ierr = PetscObjectGetComm((PetscObject)snes,&comm);CHKERRQ(ierr); 224 ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 225 ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring);CHKERRQ(ierr); 226 ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 227 ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 228 ierr = MPIU_Allreduce(&nasm->n,&N,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(ierr); 229 if (iascii) { 230 ierr = PetscViewerASCIIPrintf(viewer, " Nonlinear Additive Schwarz: total subdomain blocks = %D\n",N);CHKERRQ(ierr); 231 if (nasm->same_local_solves) { 232 if (nasm->subsnes) { 233 ierr = PetscViewerASCIIPrintf(viewer," Local solve is the same for all blocks:\n");CHKERRQ(ierr); 234 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 235 ierr = PetscViewerGetSubViewer(viewer,PETSC_COMM_SELF,&sviewer);CHKERRQ(ierr); 236 if (!rank) { 237 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 238 ierr = SNESView(nasm->subsnes[0],sviewer);CHKERRQ(ierr); 239 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 240 } 241 ierr = PetscViewerRestoreSubViewer(viewer,PETSC_COMM_SELF,&sviewer);CHKERRQ(ierr); 242 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 243 } 244 } else { 245 /* print the solver on each block */ 246 ierr = PetscViewerASCIIPushSynchronized(viewer);CHKERRQ(ierr); 247 ierr = PetscViewerASCIISynchronizedPrintf(viewer," [%d] number of local blocks = %D\n",(int)rank,nasm->n);CHKERRQ(ierr); 248 ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 249 ierr = PetscViewerASCIIPopSynchronized(viewer);CHKERRQ(ierr); 250 ierr = PetscViewerASCIIPrintf(viewer," Local solve info for each block is in the following SNES objects:\n");CHKERRQ(ierr); 251 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 252 ierr = PetscViewerASCIIPrintf(viewer,"- - - - - - - - - - - - - - - - - -\n");CHKERRQ(ierr); 253 ierr = PetscViewerGetSubViewer(viewer,PETSC_COMM_SELF,&sviewer);CHKERRQ(ierr); 254 for (i=0; i<nasm->n; i++) { 255 ierr = VecGetLocalSize(nasm->x[i],&bsz);CHKERRQ(ierr); 256 ierr = PetscViewerASCIIPrintf(sviewer,"[%d] local block number %D, size = %D\n",(int)rank,i,bsz);CHKERRQ(ierr); 257 ierr = SNESView(nasm->subsnes[i],sviewer);CHKERRQ(ierr); 258 ierr = PetscViewerASCIIPrintf(sviewer,"- - - - - - - - - - - - - - - - - -\n");CHKERRQ(ierr); 259 } 260 ierr = PetscViewerRestoreSubViewer(viewer,PETSC_COMM_SELF,&sviewer);CHKERRQ(ierr); 261 ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 262 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 263 } 264 } else if (isstring) { 265 ierr = PetscViewerStringSPrintf(viewer," blocks=%D,type=%s",N,SNESNASMTypes[nasm->type]);CHKERRQ(ierr); 266 ierr = PetscViewerGetSubViewer(viewer,PETSC_COMM_SELF,&sviewer);CHKERRQ(ierr); 267 if (nasm->subsnes && !rank) {ierr = SNESView(nasm->subsnes[0],sviewer);CHKERRQ(ierr);} 268 ierr = PetscViewerRestoreSubViewer(viewer,PETSC_COMM_SELF,&sviewer);CHKERRQ(ierr); 269 } 270 PetscFunctionReturn(0); 271 } 272 273 #undef __FUNCT__ 274 #define __FUNCT__ "SNESNASMSetType" 275 /*@ 276 SNESNASMSetType - Set the type of subdomain update used 277 278 Logically Collective on SNES 279 280 Input Parameters: 281 + SNES - the SNES context 282 - type - the type of update, PC_ASM_BASIC or PC_ASM_RESTRICT 283 284 Level: intermediate 285 286 .keywords: SNES, NASM 287 288 .seealso: SNESNASM, SNESNASMGetType(), PCASMSetType() 289 @*/ 290 PetscErrorCode SNESNASMSetType(SNES snes,PCASMType type) 291 { 292 PetscErrorCode ierr; 293 PetscErrorCode (*f)(SNES,PCASMType); 294 295 PetscFunctionBegin; 296 ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESNASMSetType_C",&f);CHKERRQ(ierr); 297 if (f) {ierr = (f)(snes,type);CHKERRQ(ierr);} 298 PetscFunctionReturn(0); 299 } 300 301 #undef __FUNCT__ 302 #define __FUNCT__ "SNESNASMSetType_NASM" 303 static PetscErrorCode SNESNASMSetType_NASM(SNES snes,PCASMType type) 304 { 305 SNES_NASM *nasm = (SNES_NASM*)snes->data; 306 307 PetscFunctionBegin; 308 if (type != PC_ASM_BASIC && type != PC_ASM_RESTRICT) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_OUTOFRANGE,"SNESNASM only supports basic and restrict types"); 309 nasm->type = type; 310 PetscFunctionReturn(0); 311 } 312 313 #undef __FUNCT__ 314 #define __FUNCT__ "SNESNASMGetType" 315 /*@ 316 SNESNASMGetType - Get the type of subdomain update used 317 318 Logically Collective on SNES 319 320 Input Parameters: 321 . SNES - the SNES context 322 323 Output Parameters: 324 . type - the type of update 325 326 Level: intermediate 327 328 .keywords: SNES, NASM 329 330 .seealso: SNESNASM, SNESNASMSetType(), PCASMGetType() 331 @*/ 332 PetscErrorCode SNESNASMGetType(SNES snes,PCASMType *type) 333 { 334 PetscErrorCode ierr; 335 336 PetscFunctionBegin; 337 ierr = PetscUseMethod(snes,"SNESNASMGetType_C",(SNES,PCASMType*),(snes,type));CHKERRQ(ierr); 338 PetscFunctionReturn(0); 339 } 340 341 #undef __FUNCT__ 342 #define __FUNCT__ "SNESNASMGetType_NASM" 343 static PetscErrorCode SNESNASMGetType_NASM(SNES snes,PCASMType *type) 344 { 345 SNES_NASM *nasm = (SNES_NASM*)snes->data; 346 347 PetscFunctionBegin; 348 *type = nasm->type; 349 PetscFunctionReturn(0); 350 } 351 352 #undef __FUNCT__ 353 #define __FUNCT__ "SNESNASMSetSubdomains" 354 /*@ 355 SNESNASMSetSubdomains - Manually Set the context required to restrict and solve subdomain problems. 356 357 Not Collective 358 359 Input Parameters: 360 + SNES - the SNES context 361 . n - the number of local subdomains 362 . subsnes - solvers defined on the local subdomains 363 . iscatter - scatters into the nonoverlapping portions of the local subdomains 364 . oscatter - scatters into the overlapping portions of the local subdomains 365 - gscatter - scatters into the (ghosted) local vector of the local subdomain 366 367 Level: intermediate 368 369 .keywords: SNES, NASM 370 371 .seealso: SNESNASM, SNESNASMGetSubdomains() 372 @*/ 373 PetscErrorCode SNESNASMSetSubdomains(SNES snes,PetscInt n,SNES subsnes[],VecScatter iscatter[],VecScatter oscatter[],VecScatter gscatter[]) 374 { 375 PetscErrorCode ierr; 376 PetscErrorCode (*f)(SNES,PetscInt,SNES*,VecScatter*,VecScatter*,VecScatter*); 377 378 PetscFunctionBegin; 379 ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESNASMSetSubdomains_C",&f);CHKERRQ(ierr); 380 if (f) {ierr = (f)(snes,n,subsnes,iscatter,oscatter,gscatter);CHKERRQ(ierr);} 381 PetscFunctionReturn(0); 382 } 383 384 #undef __FUNCT__ 385 #define __FUNCT__ "SNESNASMSetSubdomains_NASM" 386 static PetscErrorCode SNESNASMSetSubdomains_NASM(SNES snes,PetscInt n,SNES subsnes[],VecScatter iscatter[],VecScatter oscatter[],VecScatter gscatter[]) 387 { 388 PetscInt i; 389 PetscErrorCode ierr; 390 SNES_NASM *nasm = (SNES_NASM*)snes->data; 391 392 PetscFunctionBegin; 393 if (snes->setupcalled) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"SNESNASMSetSubdomains() should be called before calling SNESSetUp()."); 394 395 /* tear down the previously set things */ 396 ierr = SNESReset(snes);CHKERRQ(ierr); 397 398 nasm->n = n; 399 if (oscatter) { 400 for (i=0; i<n; i++) {ierr = PetscObjectReference((PetscObject)oscatter[i]);CHKERRQ(ierr);} 401 } 402 if (iscatter) { 403 for (i=0; i<n; i++) {ierr = PetscObjectReference((PetscObject)iscatter[i]);CHKERRQ(ierr);} 404 } 405 if (gscatter) { 406 for (i=0; i<n; i++) {ierr = PetscObjectReference((PetscObject)gscatter[i]);CHKERRQ(ierr);} 407 } 408 if (oscatter) { 409 ierr = PetscMalloc1(n,&nasm->oscatter);CHKERRQ(ierr); 410 for (i=0; i<n; i++) { 411 nasm->oscatter[i] = oscatter[i]; 412 } 413 } 414 if (iscatter) { 415 ierr = PetscMalloc1(n,&nasm->iscatter);CHKERRQ(ierr); 416 for (i=0; i<n; i++) { 417 nasm->iscatter[i] = iscatter[i]; 418 } 419 } 420 if (gscatter) { 421 ierr = PetscMalloc1(n,&nasm->gscatter);CHKERRQ(ierr); 422 for (i=0; i<n; i++) { 423 nasm->gscatter[i] = gscatter[i]; 424 } 425 } 426 427 if (subsnes) { 428 ierr = PetscMalloc1(n,&nasm->subsnes);CHKERRQ(ierr); 429 for (i=0; i<n; i++) { 430 nasm->subsnes[i] = subsnes[i]; 431 } 432 nasm->same_local_solves = PETSC_FALSE; 433 } 434 PetscFunctionReturn(0); 435 } 436 437 #undef __FUNCT__ 438 #define __FUNCT__ "SNESNASMGetSubdomains" 439 /*@ 440 SNESNASMGetSubdomains - Get the local subdomain context. 441 442 Not Collective 443 444 Input Parameters: 445 . SNES - the SNES context 446 447 Output Parameters: 448 + n - the number of local subdomains 449 . subsnes - solvers defined on the local subdomains 450 . iscatter - scatters into the nonoverlapping portions of the local subdomains 451 . oscatter - scatters into the overlapping portions of the local subdomains 452 - gscatter - scatters into the (ghosted) local vector of the local subdomain 453 454 Level: intermediate 455 456 .keywords: SNES, NASM 457 458 .seealso: SNESNASM, SNESNASMSetSubdomains() 459 @*/ 460 PetscErrorCode SNESNASMGetSubdomains(SNES snes,PetscInt *n,SNES *subsnes[],VecScatter *iscatter[],VecScatter *oscatter[],VecScatter *gscatter[]) 461 { 462 PetscErrorCode ierr; 463 PetscErrorCode (*f)(SNES,PetscInt*,SNES**,VecScatter**,VecScatter**,VecScatter**); 464 465 PetscFunctionBegin; 466 ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESNASMGetSubdomains_C",&f);CHKERRQ(ierr); 467 if (f) {ierr = (f)(snes,n,subsnes,iscatter,oscatter,gscatter);CHKERRQ(ierr);} 468 PetscFunctionReturn(0); 469 } 470 471 #undef __FUNCT__ 472 #define __FUNCT__ "SNESNASMGetSubdomains_NASM" 473 static PetscErrorCode SNESNASMGetSubdomains_NASM(SNES snes,PetscInt *n,SNES *subsnes[],VecScatter *iscatter[],VecScatter *oscatter[],VecScatter *gscatter[]) 474 { 475 SNES_NASM *nasm = (SNES_NASM*)snes->data; 476 477 PetscFunctionBegin; 478 if (n) *n = nasm->n; 479 if (oscatter) *oscatter = nasm->oscatter; 480 if (iscatter) *iscatter = nasm->iscatter; 481 if (gscatter) *gscatter = nasm->gscatter; 482 if (subsnes) { 483 *subsnes = nasm->subsnes; 484 nasm->same_local_solves = PETSC_FALSE; 485 } 486 PetscFunctionReturn(0); 487 } 488 489 #undef __FUNCT__ 490 #define __FUNCT__ "SNESNASMGetSubdomainVecs" 491 /*@ 492 SNESNASMGetSubdomainVecs - Get the processor-local subdomain vectors 493 494 Not Collective 495 496 Input Parameters: 497 . SNES - the SNES context 498 499 Output Parameters: 500 + n - the number of local subdomains 501 . x - The subdomain solution vector 502 . y - The subdomain step vector 503 . b - The subdomain RHS vector 504 - xl - The subdomain local vectors (ghosted) 505 506 Level: developer 507 508 .keywords: SNES, NASM 509 510 .seealso: SNESNASM, SNESNASMGetSubdomains() 511 @*/ 512 PetscErrorCode SNESNASMGetSubdomainVecs(SNES snes,PetscInt *n,Vec **x,Vec **y,Vec **b, Vec **xl) 513 { 514 PetscErrorCode ierr; 515 PetscErrorCode (*f)(SNES,PetscInt*,Vec**,Vec**,Vec**,Vec**); 516 517 PetscFunctionBegin; 518 ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESNASMGetSubdomainVecs_C",&f);CHKERRQ(ierr); 519 if (f) {ierr = (f)(snes,n,x,y,b,xl);CHKERRQ(ierr);} 520 PetscFunctionReturn(0); 521 } 522 523 #undef __FUNCT__ 524 #define __FUNCT__ "SNESNASMGetSubdomainVecs_NASM" 525 static PetscErrorCode SNESNASMGetSubdomainVecs_NASM(SNES snes,PetscInt *n,Vec **x,Vec **y,Vec **b,Vec **xl) 526 { 527 SNES_NASM *nasm = (SNES_NASM*)snes->data; 528 529 PetscFunctionBegin; 530 if (n) *n = nasm->n; 531 if (x) *x = nasm->x; 532 if (y) *y = nasm->y; 533 if (b) *b = nasm->b; 534 if (xl) *xl = nasm->xl; 535 PetscFunctionReturn(0); 536 } 537 538 #undef __FUNCT__ 539 #define __FUNCT__ "SNESNASMSetComputeFinalJacobian" 540 /*@ 541 SNESNASMSetComputeFinalJacobian - Schedules the computation of the global and subdomain jacobians upon convergence 542 543 Collective on SNES 544 545 Input Parameters: 546 + SNES - the SNES context 547 - flg - indication of whether to compute the jacobians or not 548 549 Level: developer 550 551 Notes: This is used almost exclusively in the implementation of ASPIN, where the converged subdomain and global jacobian 552 is needed at each linear iteration. 553 554 .keywords: SNES, NASM, ASPIN 555 556 .seealso: SNESNASM, SNESNASMGetSubdomains() 557 @*/ 558 PetscErrorCode SNESNASMSetComputeFinalJacobian(SNES snes,PetscBool flg) 559 { 560 PetscErrorCode (*f)(SNES,PetscBool); 561 PetscErrorCode ierr; 562 563 PetscFunctionBegin; 564 ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESNASMSetComputeFinalJacobian_C",&f);CHKERRQ(ierr); 565 if (f) {ierr = (f)(snes,flg);CHKERRQ(ierr);} 566 PetscFunctionReturn(0); 567 } 568 569 #undef __FUNCT__ 570 #define __FUNCT__ "SNESNASMSetComputeFinalJacobian_NASM" 571 static PetscErrorCode SNESNASMSetComputeFinalJacobian_NASM(SNES snes,PetscBool flg) 572 { 573 SNES_NASM *nasm = (SNES_NASM*)snes->data; 574 575 PetscFunctionBegin; 576 nasm->finaljacobian = flg; 577 if (flg) snes->usesksp = PETSC_TRUE; 578 PetscFunctionReturn(0); 579 } 580 581 #undef __FUNCT__ 582 #define __FUNCT__ "SNESNASMSetDamping" 583 /*@ 584 SNESNASMSetDamping - Sets the update damping for NASM 585 586 Logically collective on SNES 587 588 Input Parameters: 589 + SNES - the SNES context 590 - dmp - damping 591 592 Level: intermediate 593 594 .keywords: SNES, NASM, damping 595 596 .seealso: SNESNASM, SNESNASMGetDamping() 597 @*/ 598 PetscErrorCode SNESNASMSetDamping(SNES snes,PetscReal dmp) 599 { 600 PetscErrorCode (*f)(SNES,PetscReal); 601 PetscErrorCode ierr; 602 603 PetscFunctionBegin; 604 ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESNASMSetDamping_C",(void (**)(void))&f);CHKERRQ(ierr); 605 if (f) {ierr = (f)(snes,dmp);CHKERRQ(ierr);} 606 PetscFunctionReturn(0); 607 } 608 609 #undef __FUNCT__ 610 #define __FUNCT__ "SNESNASMSetDamping_NASM" 611 static PetscErrorCode SNESNASMSetDamping_NASM(SNES snes,PetscReal dmp) 612 { 613 SNES_NASM *nasm = (SNES_NASM*)snes->data; 614 615 PetscFunctionBegin; 616 nasm->damping = dmp; 617 PetscFunctionReturn(0); 618 } 619 620 #undef __FUNCT__ 621 #define __FUNCT__ "SNESNASMGetDamping" 622 /*@ 623 SNESNASMGetDamping - Gets the update damping for NASM 624 625 Not Collective 626 627 Input Parameters: 628 + SNES - the SNES context 629 - dmp - damping 630 631 Level: intermediate 632 633 .keywords: SNES, NASM, damping 634 635 .seealso: SNESNASM, SNESNASMSetDamping() 636 @*/ 637 PetscErrorCode SNESNASMGetDamping(SNES snes,PetscReal *dmp) 638 { 639 PetscErrorCode ierr; 640 641 PetscFunctionBegin; 642 ierr = PetscUseMethod(snes,"SNESNASMGetDamping_C",(SNES,PetscReal*),(snes,dmp));CHKERRQ(ierr); 643 PetscFunctionReturn(0); 644 } 645 646 #undef __FUNCT__ 647 #define __FUNCT__ "SNESNASMGetDamping_NASM" 648 static PetscErrorCode SNESNASMGetDamping_NASM(SNES snes,PetscReal *dmp) 649 { 650 SNES_NASM *nasm = (SNES_NASM*)snes->data; 651 652 PetscFunctionBegin; 653 *dmp = nasm->damping; 654 PetscFunctionReturn(0); 655 } 656 657 658 #undef __FUNCT__ 659 #define __FUNCT__ "SNESNASMSolveLocal_Private" 660 /* 661 Input Parameters: 662 + snes - The solver 663 . B - The RHS vector 664 - X - The initial guess 665 666 Output Parameters: 667 . Y - The solution update 668 669 TODO: All scatters should be packed into one 670 */ 671 PetscErrorCode SNESNASMSolveLocal_Private(SNES snes,Vec B,Vec Y,Vec X) 672 { 673 SNES_NASM *nasm = (SNES_NASM*)snes->data; 674 SNES subsnes; 675 PetscInt i; 676 PetscReal dmp; 677 PetscErrorCode ierr; 678 Vec Xlloc,Xl,Bl,Yl; 679 VecScatter iscat,oscat,gscat; 680 DM dm,subdm; 681 PCASMType type; 682 683 PetscFunctionBegin; 684 ierr = SNESNASMGetType(snes,&type);CHKERRQ(ierr); 685 ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 686 ierr = SNESNASMGetDamping(snes,&dmp);CHKERRQ(ierr); 687 ierr = VecSet(Y,0);CHKERRQ(ierr); 688 if (nasm->eventrestrictinterp) {ierr = PetscLogEventBegin(nasm->eventrestrictinterp,snes,0,0,0);CHKERRQ(ierr);} 689 for (i=0; i<nasm->n; i++) { 690 /* scatter the solution to the local solution */ 691 Xlloc = nasm->xl[i]; 692 gscat = nasm->gscatter[i]; 693 oscat = nasm->oscatter[i]; 694 ierr = VecScatterBegin(gscat,X,Xlloc,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 695 if (B) { 696 /* scatter the RHS to the local RHS */ 697 Bl = nasm->b[i]; 698 ierr = VecScatterBegin(oscat,B,Bl,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 699 } 700 } 701 if (nasm->eventrestrictinterp) {ierr = PetscLogEventEnd(nasm->eventrestrictinterp,snes,0,0,0);CHKERRQ(ierr);} 702 703 704 if (nasm->eventsubsolve) {ierr = PetscLogEventBegin(nasm->eventsubsolve,snes,0,0,0);CHKERRQ(ierr);} 705 for (i=0; i<nasm->n; i++) { 706 Xl = nasm->x[i]; 707 Xlloc = nasm->xl[i]; 708 Yl = nasm->y[i]; 709 subsnes = nasm->subsnes[i]; 710 ierr = SNESGetDM(subsnes,&subdm);CHKERRQ(ierr); 711 iscat = nasm->iscatter[i]; 712 oscat = nasm->oscatter[i]; 713 gscat = nasm->gscatter[i]; 714 ierr = VecScatterEnd(gscat,X,Xlloc,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 715 if (B) { 716 Bl = nasm->b[i]; 717 ierr = VecScatterEnd(oscat,B,Bl,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 718 } else Bl = NULL; 719 ierr = DMSubDomainRestrict(dm,oscat,gscat,subdm);CHKERRQ(ierr); 720 /* Could scatter directly from X */ 721 ierr = DMLocalToGlobalBegin(subdm,Xlloc,INSERT_VALUES,Xl);CHKERRQ(ierr); 722 ierr = DMLocalToGlobalEnd(subdm,Xlloc,INSERT_VALUES,Xl);CHKERRQ(ierr); 723 ierr = VecCopy(Xl,Yl);CHKERRQ(ierr); 724 ierr = SNESSolve(subsnes,Bl,Xl);CHKERRQ(ierr); 725 ierr = VecAYPX(Yl,-1.0,Xl);CHKERRQ(ierr); 726 if (type == PC_ASM_BASIC) { 727 ierr = VecScatterBegin(oscat,Yl,Y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 728 } else if (type == PC_ASM_RESTRICT) { 729 ierr = VecScatterBegin(iscat,Yl,Y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 730 } else SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"Only basic and restrict types are supported for SNESNASM"); 731 } 732 if (nasm->eventsubsolve) {ierr = PetscLogEventEnd(nasm->eventsubsolve,snes,0,0,0);CHKERRQ(ierr);} 733 if (nasm->eventrestrictinterp) {ierr = PetscLogEventBegin(nasm->eventrestrictinterp,snes,0,0,0);CHKERRQ(ierr);} 734 for (i=0; i<nasm->n; i++) { 735 Yl = nasm->y[i]; 736 iscat = nasm->iscatter[i]; 737 oscat = nasm->oscatter[i]; 738 if (type == PC_ASM_BASIC) { 739 ierr = VecScatterEnd(oscat,Yl,Y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 740 } else if (type == PC_ASM_RESTRICT) { 741 ierr = VecScatterEnd(iscat,Yl,Y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 742 } else SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"Only basic and restrict types are supported for SNESNASM"); 743 } 744 if (nasm->eventrestrictinterp) {ierr = PetscLogEventEnd(nasm->eventrestrictinterp,snes,0,0,0);CHKERRQ(ierr);} 745 ierr = VecAXPY(X,dmp,Y);CHKERRQ(ierr); 746 PetscFunctionReturn(0); 747 } 748 749 #undef __FUNCT__ 750 #define __FUNCT__ "SNESNASMComputeFinalJacobian_Private" 751 static PetscErrorCode SNESNASMComputeFinalJacobian_Private(SNES snes, Vec Xfinal) 752 { 753 Vec X = Xfinal; 754 SNES_NASM *nasm = (SNES_NASM*)snes->data; 755 SNES subsnes; 756 PetscInt i,lag = 1; 757 PetscErrorCode ierr; 758 Vec Xlloc,Xl,Fl,F; 759 VecScatter oscat,gscat; 760 DM dm,subdm; 761 762 PetscFunctionBegin; 763 if (nasm->fjtype == 2) X = nasm->xinit; 764 F = snes->vec_func; 765 if (snes->normschedule == SNES_NORM_NONE) {ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);} 766 ierr = SNESComputeJacobian(snes,X,snes->jacobian,snes->jacobian_pre);CHKERRQ(ierr); 767 ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 768 if (nasm->eventrestrictinterp) {ierr = PetscLogEventBegin(nasm->eventrestrictinterp,snes,0,0,0);CHKERRQ(ierr);} 769 if (nasm->fjtype != 1) { 770 for (i=0; i<nasm->n; i++) { 771 Xlloc = nasm->xl[i]; 772 gscat = nasm->gscatter[i]; 773 ierr = VecScatterBegin(gscat,X,Xlloc,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 774 } 775 } 776 if (nasm->eventrestrictinterp) {ierr = PetscLogEventEnd(nasm->eventrestrictinterp,snes,0,0,0);CHKERRQ(ierr);} 777 for (i=0; i<nasm->n; i++) { 778 Fl = nasm->subsnes[i]->vec_func; 779 Xl = nasm->x[i]; 780 Xlloc = nasm->xl[i]; 781 subsnes = nasm->subsnes[i]; 782 oscat = nasm->oscatter[i]; 783 gscat = nasm->gscatter[i]; 784 if (nasm->fjtype != 1) {ierr = VecScatterEnd(gscat,X,Xlloc,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);} 785 ierr = SNESGetDM(subsnes,&subdm);CHKERRQ(ierr); 786 ierr = DMSubDomainRestrict(dm,oscat,gscat,subdm);CHKERRQ(ierr); 787 if (nasm->fjtype != 1) { 788 ierr = DMLocalToGlobalBegin(subdm,Xlloc,INSERT_VALUES,Xl);CHKERRQ(ierr); 789 ierr = DMLocalToGlobalEnd(subdm,Xlloc,INSERT_VALUES,Xl);CHKERRQ(ierr); 790 } 791 if (subsnes->lagjacobian == -1) subsnes->lagjacobian = -2; 792 else if (subsnes->lagjacobian > 1) lag = subsnes->lagjacobian; 793 ierr = SNESComputeFunction(subsnes,Xl,Fl);CHKERRQ(ierr); 794 ierr = SNESComputeJacobian(subsnes,Xl,subsnes->jacobian,subsnes->jacobian_pre);CHKERRQ(ierr); 795 if (lag > 1) subsnes->lagjacobian = lag; 796 } 797 PetscFunctionReturn(0); 798 } 799 800 #undef __FUNCT__ 801 #define __FUNCT__ "SNESSolve_NASM" 802 static PetscErrorCode SNESSolve_NASM(SNES snes) 803 { 804 Vec F; 805 Vec X; 806 Vec B; 807 Vec Y; 808 PetscInt i; 809 PetscReal fnorm = 0.0; 810 PetscErrorCode ierr; 811 SNESNormSchedule normschedule; 812 SNES_NASM *nasm = (SNES_NASM*)snes->data; 813 814 PetscFunctionBegin; 815 816 if (snes->xl || snes->xu || snes->ops->computevariablebounds) SETERRQ1(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE, "SNES solver %s does not support bounds", ((PetscObject)snes)->type_name); 817 818 ierr = PetscCitationsRegister(SNESCitation,&SNEScite);CHKERRQ(ierr); 819 X = snes->vec_sol; 820 Y = snes->vec_sol_update; 821 F = snes->vec_func; 822 B = snes->vec_rhs; 823 824 ierr = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr); 825 snes->iter = 0; 826 snes->norm = 0.; 827 ierr = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr); 828 snes->reason = SNES_CONVERGED_ITERATING; 829 ierr = SNESGetNormSchedule(snes, &normschedule);CHKERRQ(ierr); 830 if (normschedule == SNES_NORM_ALWAYS || normschedule == SNES_NORM_INITIAL_ONLY || normschedule == SNES_NORM_INITIAL_FINAL_ONLY) { 831 /* compute the initial function and preconditioned update delX */ 832 if (!snes->vec_func_init_set) { 833 ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 834 } else snes->vec_func_init_set = PETSC_FALSE; 835 836 ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F|| */ 837 SNESCheckFunctionNorm(snes,fnorm); 838 ierr = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr); 839 snes->iter = 0; 840 snes->norm = fnorm; 841 ierr = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr); 842 ierr = SNESLogConvergenceHistory(snes,snes->norm,0);CHKERRQ(ierr); 843 ierr = SNESMonitor(snes,0,snes->norm);CHKERRQ(ierr); 844 845 /* test convergence */ 846 ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 847 if (snes->reason) PetscFunctionReturn(0); 848 } else { 849 ierr = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr); 850 ierr = SNESLogConvergenceHistory(snes,snes->norm,0);CHKERRQ(ierr); 851 ierr = SNESMonitor(snes,0,snes->norm);CHKERRQ(ierr); 852 } 853 854 /* Call general purpose update function */ 855 if (snes->ops->update) { 856 ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 857 } 858 /* copy the initial solution over for later */ 859 if (nasm->fjtype == 2) {ierr = VecCopy(X,nasm->xinit);CHKERRQ(ierr);} 860 861 for (i = 0; i < snes->max_its; i++) { 862 ierr = SNESNASMSolveLocal_Private(snes,B,Y,X);CHKERRQ(ierr); 863 if (normschedule == SNES_NORM_ALWAYS || ((i == snes->max_its - 1) && (normschedule == SNES_NORM_INITIAL_FINAL_ONLY || normschedule == SNES_NORM_FINAL_ONLY))) { 864 ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 865 ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F|| */ 866 SNESCheckFunctionNorm(snes,fnorm); 867 } 868 /* Monitor convergence */ 869 ierr = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr); 870 snes->iter = i+1; 871 snes->norm = fnorm; 872 ierr = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr); 873 ierr = SNESLogConvergenceHistory(snes,snes->norm,0);CHKERRQ(ierr); 874 ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr); 875 /* Test for convergence */ 876 if (normschedule == SNES_NORM_ALWAYS) {ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);} 877 if (snes->reason) break; 878 /* Call general purpose update function */ 879 if (snes->ops->update) {ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);} 880 } 881 if (nasm->finaljacobian) {ierr = SNESNASMComputeFinalJacobian_Private(snes,X);CHKERRQ(ierr);} 882 if (normschedule == SNES_NORM_ALWAYS) { 883 if (i == snes->max_its) { 884 ierr = PetscInfo1(snes,"Maximum number of iterations has been reached: %D\n",snes->max_its);CHKERRQ(ierr); 885 if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT; 886 } 887 } else if (!snes->reason) snes->reason = SNES_CONVERGED_ITS; /* NASM is meant to be used as a preconditioner */ 888 PetscFunctionReturn(0); 889 } 890 891 /*MC 892 SNESNASM - Nonlinear Additive Schwartz 893 894 Options Database: 895 + -snes_nasm_log - enable logging events for the communication and solve stages 896 . -snes_nasm_type <basic,restrict> - type of subdomain update used 897 . -snes_nasm_finaljacobian - compute the local and global jacobians of the final iterate 898 . -snes_nasm_finaljacobian_type <finalinner,finalouter,initial> - pick state the jacobian is calculated at 899 . -sub_snes_ - options prefix of the subdomain nonlinear solves 900 . -sub_ksp_ - options prefix of the subdomain Krylov solver 901 - -sub_pc_ - options prefix of the subdomain preconditioner 902 903 Level: advanced 904 905 References: 906 . 1. - Peter R. Brune, Matthew G. Knepley, Barry F. Smith, and Xuemin Tu, "Composing Scalable Nonlinear Algebraic Solvers", 907 SIAM Review, 57(4), 2015 908 909 .seealso: SNESCreate(), SNES, SNESSetType(), SNESType (for list of available types) 910 M*/ 911 912 #undef __FUNCT__ 913 #define __FUNCT__ "SNESCreate_NASM" 914 PETSC_EXTERN PetscErrorCode SNESCreate_NASM(SNES snes) 915 { 916 SNES_NASM *nasm; 917 PetscErrorCode ierr; 918 919 PetscFunctionBegin; 920 ierr = PetscNewLog(snes,&nasm);CHKERRQ(ierr); 921 snes->data = (void*)nasm; 922 923 nasm->n = PETSC_DECIDE; 924 nasm->subsnes = 0; 925 nasm->x = 0; 926 nasm->xl = 0; 927 nasm->y = 0; 928 nasm->b = 0; 929 nasm->oscatter = 0; 930 nasm->iscatter = 0; 931 nasm->gscatter = 0; 932 nasm->damping = 1.; 933 934 nasm->type = PC_ASM_BASIC; 935 nasm->finaljacobian = PETSC_FALSE; 936 nasm->same_local_solves = PETSC_TRUE; 937 938 snes->ops->destroy = SNESDestroy_NASM; 939 snes->ops->setup = SNESSetUp_NASM; 940 snes->ops->setfromoptions = SNESSetFromOptions_NASM; 941 snes->ops->view = SNESView_NASM; 942 snes->ops->solve = SNESSolve_NASM; 943 snes->ops->reset = SNESReset_NASM; 944 945 snes->usesksp = PETSC_FALSE; 946 snes->usespc = PETSC_FALSE; 947 948 nasm->fjtype = 0; 949 nasm->xinit = NULL; 950 nasm->eventrestrictinterp = 0; 951 nasm->eventsubsolve = 0; 952 953 if (!snes->tolerancesset) { 954 snes->max_its = 10000; 955 snes->max_funcs = 10000; 956 } 957 958 ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESNASMSetType_C",SNESNASMSetType_NASM);CHKERRQ(ierr); 959 ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESNASMGetType_C",SNESNASMGetType_NASM);CHKERRQ(ierr); 960 ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESNASMSetSubdomains_C",SNESNASMSetSubdomains_NASM);CHKERRQ(ierr); 961 ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESNASMGetSubdomains_C",SNESNASMGetSubdomains_NASM);CHKERRQ(ierr); 962 ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESNASMSetDamping_C",SNESNASMSetDamping_NASM);CHKERRQ(ierr); 963 ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESNASMGetDamping_C",SNESNASMGetDamping_NASM);CHKERRQ(ierr); 964 ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESNASMGetSubdomainVecs_C",SNESNASMGetSubdomainVecs_NASM);CHKERRQ(ierr); 965 ierr = PetscObjectComposeFunction((PetscObject)snes,"SNESNASMSetComputeFinalJacobian_C",SNESNASMSetComputeFinalJacobian_NASM);CHKERRQ(ierr); 966 PetscFunctionReturn(0); 967 } 968 969