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