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