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 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 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 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 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 } 160 ierr = DMGlobalToLocalHookAdd(subdm,DMGlobalToLocalSubDomainDirichletHook_Private,NULL,nasm->xl[i]);CHKERRQ(ierr); 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 PetscErrorCode SNESSetFromOptions_NASM(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("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 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 = MPI_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 = PetscViewerGetSingleton(viewer,&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 = PetscViewerRestoreSingleton(viewer,&sviewer);CHKERRQ(ierr); 242 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 243 } 244 } else { 245 /* print the solver on each block */ 246 ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_TRUE);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 = PetscViewerASCIISynchronizedAllow(viewer,PETSC_FALSE);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 = PetscViewerGetSingleton(viewer,&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 = PetscViewerRestoreSingleton(viewer,&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 = PetscViewerGetSingleton(viewer,&sviewer);CHKERRQ(ierr); 267 if (nasm->subsnes && !rank) {ierr = SNESView(nasm->subsnes[0],sviewer);CHKERRQ(ierr);} 268 ierr = PetscViewerRestoreSingleton(viewer,&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 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 PetscErrorCode (*f)(SNES,PCASMType*); 336 337 PetscFunctionBegin; 338 ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESNASMGetType_C",&f);CHKERRQ(ierr); 339 if (f) {ierr = (f)(snes,type);CHKERRQ(ierr);} 340 PetscFunctionReturn(0); 341 } 342 343 #undef __FUNCT__ 344 #define __FUNCT__ "SNESNASMGetType_NASM" 345 PetscErrorCode SNESNASMGetType_NASM(SNES snes,PCASMType *type) 346 { 347 SNES_NASM *nasm = (SNES_NASM*)snes->data; 348 349 PetscFunctionBegin; 350 *type = nasm->type; 351 PetscFunctionReturn(0); 352 } 353 354 #undef __FUNCT__ 355 #define __FUNCT__ "SNESNASMSetSubdomains" 356 /*@ 357 SNESNASMSetSubdomains - Manually Set the context required to restrict and solve subdomain problems. 358 359 Not Collective 360 361 Input Parameters: 362 + SNES - the SNES context 363 . n - the number of local subdomains 364 . subsnes - solvers defined on the local subdomains 365 . iscatter - scatters into the nonoverlapping portions of the local subdomains 366 . oscatter - scatters into the overlapping portions of the local subdomains 367 - gscatter - scatters into the (ghosted) local vector of the local subdomain 368 369 Level: intermediate 370 371 .keywords: SNES, NASM 372 373 .seealso: SNESNASM, SNESNASMGetSubdomains() 374 @*/ 375 PetscErrorCode SNESNASMSetSubdomains(SNES snes,PetscInt n,SNES subsnes[],VecScatter iscatter[],VecScatter oscatter[],VecScatter gscatter[]) 376 { 377 PetscErrorCode ierr; 378 PetscErrorCode (*f)(SNES,PetscInt,SNES*,VecScatter*,VecScatter*,VecScatter*); 379 380 PetscFunctionBegin; 381 ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESNASMSetSubdomains_C",&f);CHKERRQ(ierr); 382 if (f) {ierr = (f)(snes,n,subsnes,iscatter,oscatter,gscatter);CHKERRQ(ierr);} 383 PetscFunctionReturn(0); 384 } 385 386 #undef __FUNCT__ 387 #define __FUNCT__ "SNESNASMSetSubdomains_NASM" 388 PetscErrorCode SNESNASMSetSubdomains_NASM(SNES snes,PetscInt n,SNES subsnes[],VecScatter iscatter[],VecScatter oscatter[],VecScatter gscatter[]) 389 { 390 PetscInt i; 391 PetscErrorCode ierr; 392 SNES_NASM *nasm = (SNES_NASM*)snes->data; 393 394 PetscFunctionBegin; 395 if (snes->setupcalled) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"SNESNASMSetSubdomains() should be called before calling SNESSetUp()."); 396 397 /* tear down the previously set things */ 398 ierr = SNESReset(snes);CHKERRQ(ierr); 399 400 nasm->n = n; 401 if (oscatter) { 402 for (i=0; i<n; i++) {ierr = PetscObjectReference((PetscObject)oscatter[i]);CHKERRQ(ierr);} 403 } 404 if (iscatter) { 405 for (i=0; i<n; i++) {ierr = PetscObjectReference((PetscObject)iscatter[i]);CHKERRQ(ierr);} 406 } 407 if (gscatter) { 408 for (i=0; i<n; i++) {ierr = PetscObjectReference((PetscObject)gscatter[i]);CHKERRQ(ierr);} 409 } 410 if (oscatter) { 411 ierr = PetscMalloc1(n,&nasm->oscatter);CHKERRQ(ierr); 412 for (i=0; i<n; i++) { 413 nasm->oscatter[i] = oscatter[i]; 414 } 415 } 416 if (iscatter) { 417 ierr = PetscMalloc1(n,&nasm->iscatter);CHKERRQ(ierr); 418 for (i=0; i<n; i++) { 419 nasm->iscatter[i] = iscatter[i]; 420 } 421 } 422 if (gscatter) { 423 ierr = PetscMalloc1(n,&nasm->gscatter);CHKERRQ(ierr); 424 for (i=0; i<n; i++) { 425 nasm->gscatter[i] = gscatter[i]; 426 } 427 } 428 429 if (subsnes) { 430 ierr = PetscMalloc1(n,&nasm->subsnes);CHKERRQ(ierr); 431 for (i=0; i<n; i++) { 432 nasm->subsnes[i] = subsnes[i]; 433 } 434 nasm->same_local_solves = PETSC_FALSE; 435 } 436 PetscFunctionReturn(0); 437 } 438 439 #undef __FUNCT__ 440 #define __FUNCT__ "SNESNASMGetSubdomains" 441 /*@ 442 SNESNASMGetSubdomains - Get the local subdomain context. 443 444 Not Collective 445 446 Input Parameters: 447 . SNES - the SNES context 448 449 Output Parameters: 450 + n - the number of local subdomains 451 . subsnes - solvers defined on the local subdomains 452 . iscatter - scatters into the nonoverlapping portions of the local subdomains 453 . oscatter - scatters into the overlapping portions of the local subdomains 454 - gscatter - scatters into the (ghosted) local vector of the local subdomain 455 456 Level: intermediate 457 458 .keywords: SNES, NASM 459 460 .seealso: SNESNASM, SNESNASMSetSubdomains() 461 @*/ 462 PetscErrorCode SNESNASMGetSubdomains(SNES snes,PetscInt *n,SNES *subsnes[],VecScatter *iscatter[],VecScatter *oscatter[],VecScatter *gscatter[]) 463 { 464 PetscErrorCode ierr; 465 PetscErrorCode (*f)(SNES,PetscInt*,SNES**,VecScatter**,VecScatter**,VecScatter**); 466 467 PetscFunctionBegin; 468 ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESNASMGetSubdomains_C",&f);CHKERRQ(ierr); 469 if (f) {ierr = (f)(snes,n,subsnes,iscatter,oscatter,gscatter);CHKERRQ(ierr);} 470 PetscFunctionReturn(0); 471 } 472 473 #undef __FUNCT__ 474 #define __FUNCT__ "SNESNASMGetSubdomains_NASM" 475 PetscErrorCode SNESNASMGetSubdomains_NASM(SNES snes,PetscInt *n,SNES *subsnes[],VecScatter *iscatter[],VecScatter *oscatter[],VecScatter *gscatter[]) 476 { 477 SNES_NASM *nasm = (SNES_NASM*)snes->data; 478 479 PetscFunctionBegin; 480 if (n) *n = nasm->n; 481 if (oscatter) *oscatter = nasm->oscatter; 482 if (iscatter) *iscatter = nasm->iscatter; 483 if (gscatter) *gscatter = nasm->gscatter; 484 if (subsnes) { 485 *subsnes = nasm->subsnes; 486 nasm->same_local_solves = PETSC_FALSE; 487 } 488 PetscFunctionReturn(0); 489 } 490 491 #undef __FUNCT__ 492 #define __FUNCT__ "SNESNASMGetSubdomainVecs" 493 /*@ 494 SNESNASMGetSubdomainVecs - Get the processor-local subdomain vectors 495 496 Not Collective 497 498 Input Parameters: 499 . SNES - the SNES context 500 501 Output Parameters: 502 + n - the number of local subdomains 503 . x - The subdomain solution vector 504 . y - The subdomain step vector 505 . b - The subdomain RHS vector 506 - xl - The subdomain local vectors (ghosted) 507 508 Level: developer 509 510 .keywords: SNES, NASM 511 512 .seealso: SNESNASM, SNESNASMGetSubdomains() 513 @*/ 514 PetscErrorCode SNESNASMGetSubdomainVecs(SNES snes,PetscInt *n,Vec **x,Vec **y,Vec **b, Vec **xl) 515 { 516 PetscErrorCode ierr; 517 PetscErrorCode (*f)(SNES,PetscInt*,Vec**,Vec**,Vec**,Vec**); 518 519 PetscFunctionBegin; 520 ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESNASMGetSubdomainVecs_C",&f);CHKERRQ(ierr); 521 if (f) {ierr = (f)(snes,n,x,y,b,xl);CHKERRQ(ierr);} 522 PetscFunctionReturn(0); 523 } 524 525 #undef __FUNCT__ 526 #define __FUNCT__ "SNESNASMGetSubdomainVecs_NASM" 527 PetscErrorCode SNESNASMGetSubdomainVecs_NASM(SNES snes,PetscInt *n,Vec **x,Vec **y,Vec **b,Vec **xl) 528 { 529 SNES_NASM *nasm = (SNES_NASM*)snes->data; 530 531 PetscFunctionBegin; 532 if (n) *n = nasm->n; 533 if (x) *x = nasm->x; 534 if (y) *y = nasm->y; 535 if (b) *b = nasm->b; 536 if (xl) *xl = nasm->xl; 537 PetscFunctionReturn(0); 538 } 539 540 #undef __FUNCT__ 541 #define __FUNCT__ "SNESNASMSetComputeFinalJacobian" 542 /*@ 543 SNESNASMSetComputeFinalJacobian - Schedules the computation of the global and subdomain jacobians upon convergence 544 545 Collective on SNES 546 547 Input Parameters: 548 + SNES - the SNES context 549 - flg - indication of whether to compute the jacobians or not 550 551 Level: developer 552 553 Notes: This is used almost exclusively in the implementation of ASPIN, where the converged subdomain and global jacobian 554 is needed at each linear iteration. 555 556 .keywords: SNES, NASM, ASPIN 557 558 .seealso: SNESNASM, SNESNASMGetSubdomains() 559 @*/ 560 PetscErrorCode SNESNASMSetComputeFinalJacobian(SNES snes,PetscBool flg) 561 { 562 PetscErrorCode (*f)(SNES,PetscBool); 563 PetscErrorCode ierr; 564 565 PetscFunctionBegin; 566 ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESNASMSetComputeFinalJacobian_C",&f);CHKERRQ(ierr); 567 if (f) {ierr = (f)(snes,flg);CHKERRQ(ierr);} 568 PetscFunctionReturn(0); 569 } 570 571 #undef __FUNCT__ 572 #define __FUNCT__ "SNESNASMSetComputeFinalJacobian_NASM" 573 PetscErrorCode SNESNASMSetComputeFinalJacobian_NASM(SNES snes,PetscBool flg) 574 { 575 SNES_NASM *nasm = (SNES_NASM*)snes->data; 576 577 PetscFunctionBegin; 578 nasm->finaljacobian = flg; 579 if (flg) snes->usesksp = PETSC_TRUE; 580 PetscFunctionReturn(0); 581 } 582 583 #undef __FUNCT__ 584 #define __FUNCT__ "SNESNASMSetDamping" 585 /*@ 586 SNESNASMSetDamping - Sets the update damping for NASM 587 588 Logically collective on SNES 589 590 Input Parameters: 591 + SNES - the SNES context 592 - dmp - damping 593 594 Level: intermediate 595 596 .keywords: SNES, NASM, damping 597 598 .seealso: SNESNASM, SNESNASMGetDamping() 599 @*/ 600 PetscErrorCode SNESNASMSetDamping(SNES snes,PetscReal dmp) 601 { 602 PetscErrorCode (*f)(SNES,PetscReal); 603 PetscErrorCode ierr; 604 605 PetscFunctionBegin; 606 ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESNASMSetDamping_C",(void (**)(void))&f);CHKERRQ(ierr); 607 if (f) {ierr = (f)(snes,dmp);CHKERRQ(ierr);} 608 PetscFunctionReturn(0); 609 } 610 611 #undef __FUNCT__ 612 #define __FUNCT__ "SNESNASMSetDamping_NASM" 613 PetscErrorCode SNESNASMSetDamping_NASM(SNES snes,PetscReal dmp) 614 { 615 SNES_NASM *nasm = (SNES_NASM*)snes->data; 616 617 PetscFunctionBegin; 618 nasm->damping = dmp; 619 PetscFunctionReturn(0); 620 } 621 622 #undef __FUNCT__ 623 #define __FUNCT__ "SNESNASMGetDamping" 624 /*@ 625 SNESNASMGetDamping - Gets the update damping for NASM 626 627 Not Collective 628 629 Input Parameters: 630 + SNES - the SNES context 631 - dmp - damping 632 633 Level: intermediate 634 635 .keywords: SNES, NASM, damping 636 637 .seealso: SNESNASM, SNESNASMSetDamping() 638 @*/ 639 PetscErrorCode SNESNASMGetDamping(SNES snes,PetscReal *dmp) 640 { 641 PetscErrorCode (*f)(SNES,PetscReal*); 642 PetscErrorCode ierr; 643 644 PetscFunctionBegin; 645 ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESNASMGetDamping_C",(void (**)(void))&f);CHKERRQ(ierr); 646 if (f) {ierr = (f)(snes,dmp);CHKERRQ(ierr);} 647 PetscFunctionReturn(0); 648 } 649 650 #undef __FUNCT__ 651 #define __FUNCT__ "SNESNASMGetDamping_NASM" 652 PetscErrorCode SNESNASMGetDamping_NASM(SNES snes,PetscReal *dmp) 653 { 654 SNES_NASM *nasm = (SNES_NASM*)snes->data; 655 656 PetscFunctionBegin; 657 *dmp = nasm->damping; 658 PetscFunctionReturn(0); 659 } 660 661 662 #undef __FUNCT__ 663 #define __FUNCT__ "SNESNASMSolveLocal_Private" 664 PetscErrorCode SNESNASMSolveLocal_Private(SNES snes,Vec B,Vec Y,Vec X) 665 { 666 SNES_NASM *nasm = (SNES_NASM*)snes->data; 667 SNES subsnes; 668 PetscInt i; 669 PetscReal dmp; 670 PetscErrorCode ierr; 671 Vec Xlloc,Xl,Bl,Yl; 672 VecScatter iscat,oscat,gscat; 673 DM dm,subdm; 674 PCASMType type; 675 676 PetscFunctionBegin; 677 ierr = SNESNASMGetType(snes,&type);CHKERRQ(ierr); 678 ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 679 ierr = SNESNASMGetDamping(snes,&dmp);CHKERRQ(ierr); 680 ierr = VecSet(Y,0);CHKERRQ(ierr); 681 if (nasm->eventrestrictinterp) {ierr = PetscLogEventBegin(nasm->eventrestrictinterp,snes,0,0,0);CHKERRQ(ierr);} 682 for (i=0; i<nasm->n; i++) { 683 /* scatter the solution to the local solution */ 684 Xlloc = nasm->xl[i]; 685 gscat = nasm->gscatter[i]; 686 oscat = nasm->oscatter[i]; 687 ierr = VecScatterBegin(gscat,X,Xlloc,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 688 if (B) { 689 /* scatter the RHS to the local RHS */ 690 Bl = nasm->b[i]; 691 ierr = VecScatterBegin(oscat,B,Bl,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 692 } 693 } 694 if (nasm->eventrestrictinterp) {ierr = PetscLogEventEnd(nasm->eventrestrictinterp,snes,0,0,0);CHKERRQ(ierr);} 695 696 697 if (nasm->eventsubsolve) {ierr = PetscLogEventBegin(nasm->eventsubsolve,snes,0,0,0);CHKERRQ(ierr);} 698 for (i=0; i<nasm->n; i++) { 699 Xl = nasm->x[i]; 700 Xlloc = nasm->xl[i]; 701 Yl = nasm->y[i]; 702 subsnes = nasm->subsnes[i]; 703 ierr = SNESGetDM(subsnes,&subdm);CHKERRQ(ierr); 704 iscat = nasm->iscatter[i]; 705 oscat = nasm->oscatter[i]; 706 gscat = nasm->gscatter[i]; 707 ierr = VecScatterEnd(gscat,X,Xlloc,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 708 if (B) { 709 Bl = nasm->b[i]; 710 ierr = VecScatterEnd(oscat,B,Bl,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 711 } else Bl = NULL; 712 ierr = DMSubDomainRestrict(dm,oscat,gscat,subdm);CHKERRQ(ierr); 713 ierr = DMLocalToGlobalBegin(subdm,Xlloc,INSERT_VALUES,Xl);CHKERRQ(ierr); 714 ierr = DMLocalToGlobalEnd(subdm,Xlloc,INSERT_VALUES,Xl);CHKERRQ(ierr); 715 ierr = VecCopy(Xl,Yl);CHKERRQ(ierr); 716 ierr = SNESSolve(subsnes,Bl,Xl);CHKERRQ(ierr); 717 ierr = VecAYPX(Yl,-1.0,Xl);CHKERRQ(ierr); 718 if (type == PC_ASM_BASIC) { 719 ierr = VecScatterBegin(oscat,Yl,Y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 720 } else if (type == PC_ASM_RESTRICT) { 721 ierr = VecScatterBegin(iscat,Yl,Y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 722 } else SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"Only basic and restrict types are supported for SNESNASM"); 723 } 724 if (nasm->eventsubsolve) {ierr = PetscLogEventEnd(nasm->eventsubsolve,snes,0,0,0);CHKERRQ(ierr);} 725 if (nasm->eventrestrictinterp) {ierr = PetscLogEventBegin(nasm->eventrestrictinterp,snes,0,0,0);CHKERRQ(ierr);} 726 for (i=0; i<nasm->n; i++) { 727 Yl = nasm->y[i]; 728 iscat = nasm->iscatter[i]; 729 oscat = nasm->oscatter[i]; 730 if (type == PC_ASM_BASIC) { 731 ierr = VecScatterEnd(oscat,Yl,Y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 732 } else if (type == PC_ASM_RESTRICT) { 733 ierr = VecScatterEnd(iscat,Yl,Y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 734 } else SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"Only basic and restrict types are supported for SNESNASM"); 735 } 736 if (nasm->eventrestrictinterp) {ierr = PetscLogEventEnd(nasm->eventrestrictinterp,snes,0,0,0);CHKERRQ(ierr);} 737 ierr = VecAXPY(X,dmp,Y);CHKERRQ(ierr); 738 PetscFunctionReturn(0); 739 } 740 741 #undef __FUNCT__ 742 #define __FUNCT__ "SNESNASMComputeFinalJacobian_Private" 743 PetscErrorCode SNESNASMComputeFinalJacobian_Private(SNES snes, Vec Xfinal) 744 { 745 Vec X = Xfinal; 746 SNES_NASM *nasm = (SNES_NASM*)snes->data; 747 SNES subsnes; 748 PetscInt i,lag = 1; 749 PetscErrorCode ierr; 750 Vec Xlloc,Xl,Fl,F; 751 VecScatter oscat,gscat; 752 DM dm,subdm; 753 754 PetscFunctionBegin; 755 if (nasm->fjtype == 2) X = nasm->xinit; 756 F = snes->vec_func; 757 if (snes->normschedule == SNES_NORM_NONE) {ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);} 758 ierr = SNESComputeJacobian(snes,X,snes->jacobian,snes->jacobian_pre);CHKERRQ(ierr); 759 ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 760 if (nasm->eventrestrictinterp) {ierr = PetscLogEventBegin(nasm->eventrestrictinterp,snes,0,0,0);CHKERRQ(ierr);} 761 if (nasm->fjtype != 1) { 762 for (i=0; i<nasm->n; i++) { 763 Xlloc = nasm->xl[i]; 764 gscat = nasm->gscatter[i]; 765 oscat = nasm->oscatter[i]; 766 ierr = VecScatterBegin(gscat,X,Xlloc,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 767 } 768 } 769 if (nasm->eventrestrictinterp) {ierr = PetscLogEventEnd(nasm->eventrestrictinterp,snes,0,0,0);CHKERRQ(ierr);} 770 for (i=0; i<nasm->n; i++) { 771 Fl = nasm->subsnes[i]->vec_func; 772 Xl = nasm->x[i]; 773 Xlloc = nasm->xl[i]; 774 subsnes = nasm->subsnes[i]; 775 oscat = nasm->oscatter[i]; 776 gscat = nasm->gscatter[i]; 777 if (nasm->fjtype != 1) {ierr = VecScatterEnd(gscat,X,Xlloc,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);} 778 ierr = SNESGetDM(subsnes,&subdm);CHKERRQ(ierr); 779 ierr = DMSubDomainRestrict(dm,oscat,gscat,subdm);CHKERRQ(ierr); 780 if (nasm->fjtype != 1) { 781 ierr = DMLocalToGlobalBegin(subdm,Xlloc,INSERT_VALUES,Xl);CHKERRQ(ierr); 782 ierr = DMLocalToGlobalEnd(subdm,Xlloc,INSERT_VALUES,Xl);CHKERRQ(ierr); 783 } 784 if (subsnes->lagjacobian == -1) subsnes->lagjacobian = -2; 785 else if (subsnes->lagjacobian > 1) lag = subsnes->lagjacobian; 786 ierr = SNESComputeFunction(subsnes,Xl,Fl);CHKERRQ(ierr); 787 ierr = SNESComputeJacobian(subsnes,Xl,subsnes->jacobian,subsnes->jacobian_pre);CHKERRQ(ierr); 788 if (lag > 1) subsnes->lagjacobian = lag; 789 } 790 PetscFunctionReturn(0); 791 } 792 793 #undef __FUNCT__ 794 #define __FUNCT__ "SNESSolve_NASM" 795 PetscErrorCode SNESSolve_NASM(SNES snes) 796 { 797 Vec F; 798 Vec X; 799 Vec B; 800 Vec Y; 801 PetscInt i; 802 PetscReal fnorm = 0.0; 803 PetscErrorCode ierr; 804 SNESNormSchedule normschedule; 805 SNES_NASM *nasm = (SNES_NASM*)snes->data; 806 807 PetscFunctionBegin; 808 ierr = PetscCitationsRegister(SNESCitation,&SNEScite);CHKERRQ(ierr); 809 X = snes->vec_sol; 810 Y = snes->vec_sol_update; 811 F = snes->vec_func; 812 B = snes->vec_rhs; 813 814 ierr = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr); 815 snes->iter = 0; 816 snes->norm = 0.; 817 ierr = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr); 818 snes->reason = SNES_CONVERGED_ITERATING; 819 ierr = SNESGetNormSchedule(snes, &normschedule);CHKERRQ(ierr); 820 if (normschedule == SNES_NORM_ALWAYS || normschedule == SNES_NORM_INITIAL_ONLY || normschedule == SNES_NORM_INITIAL_FINAL_ONLY) { 821 /* compute the initial function and preconditioned update delX */ 822 if (!snes->vec_func_init_set) { 823 ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 824 if (snes->domainerror) { 825 snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 826 PetscFunctionReturn(0); 827 } 828 } else snes->vec_func_init_set = PETSC_FALSE; 829 830 ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F|| */ 831 if (PetscIsInfOrNanReal(fnorm)) { 832 snes->reason = SNES_DIVERGED_FNORM_NAN; 833 PetscFunctionReturn(0); 834 } 835 ierr = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr); 836 snes->iter = 0; 837 snes->norm = fnorm; 838 ierr = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr); 839 ierr = SNESLogConvergenceHistory(snes,snes->norm,0);CHKERRQ(ierr); 840 ierr = SNESMonitor(snes,0,snes->norm);CHKERRQ(ierr); 841 842 /* test convergence */ 843 ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 844 if (snes->reason) PetscFunctionReturn(0); 845 } else { 846 ierr = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr); 847 ierr = SNESLogConvergenceHistory(snes,snes->norm,0);CHKERRQ(ierr); 848 ierr = SNESMonitor(snes,0,snes->norm);CHKERRQ(ierr); 849 } 850 851 /* Call general purpose update function */ 852 if (snes->ops->update) { 853 ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 854 } 855 /* copy the initial solution over for later */ 856 if (nasm->fjtype == 2) {ierr = VecCopy(X,nasm->xinit);CHKERRQ(ierr);} 857 858 for (i = 0; i < snes->max_its; i++) { 859 ierr = SNESNASMSolveLocal_Private(snes,B,Y,X);CHKERRQ(ierr); 860 if (normschedule == SNES_NORM_ALWAYS || ((i == snes->max_its - 1) && (normschedule == SNES_NORM_INITIAL_FINAL_ONLY || normschedule == SNES_NORM_FINAL_ONLY))) { 861 ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 862 if (snes->domainerror) { 863 snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 864 break; 865 } 866 ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F|| */ 867 if (PetscIsInfOrNanReal(fnorm)) { 868 snes->reason = SNES_DIVERGED_FNORM_NAN; 869 break; 870 } 871 } 872 /* Monitor convergence */ 873 ierr = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr); 874 snes->iter = i+1; 875 snes->norm = fnorm; 876 ierr = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr); 877 ierr = SNESLogConvergenceHistory(snes,snes->norm,0);CHKERRQ(ierr); 878 ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr); 879 /* Test for convergence */ 880 if (normschedule == SNES_NORM_ALWAYS) {ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);} 881 if (snes->reason) break; 882 /* Call general purpose update function */ 883 if (snes->ops->update) {ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);} 884 } 885 if (nasm->finaljacobian) {ierr = SNESNASMComputeFinalJacobian_Private(snes,X);CHKERRQ(ierr);} 886 if (normschedule == SNES_NORM_ALWAYS) { 887 if (i == snes->max_its) { 888 ierr = PetscInfo1(snes,"Maximum number of iterations has been reached: %D\n",snes->max_its);CHKERRQ(ierr); 889 if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT; 890 } 891 } else if (!snes->reason) snes->reason = SNES_CONVERGED_ITS; /* NASM is meant to be used as a preconditioner */ 892 PetscFunctionReturn(0); 893 } 894 895 /*MC 896 SNESNASM - Nonlinear Additive Schwartz 897 898 Options Database: 899 + -snes_nasm_log - enable logging events for the communication and solve stages 900 . -snes_nasm_type <basic,restrict> - type of subdomain update used 901 . -snes_nasm_finaljacobian - compute the local and global jacobians of the final iterate 902 . -snes_nasm_finaljacobian_type <finalinner,finalouter,initial> pick state the jacobian is calculated at 903 . -sub_snes_ - options prefix of the subdomain nonlinear solves 904 . -sub_ksp_ - options prefix of the subdomain Krylov solver 905 - -sub_pc_ - options prefix of the subdomain preconditioner 906 907 Level: advanced 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