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