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