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