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 8 Vec *x; /* solution vectors */ 9 Vec *xl; /* solution local vectors */ 10 Vec *y; /* step vectors */ 11 Vec *b; /* rhs vectors */ 12 13 VecScatter *oscatter; /* scatter from global space to the subdomain global space */ 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 17 PCASMType type; /* ASM type */ 18 19 PetscBool usesdm; /* use the DM for setting up the subproblems */ 20 } SNES_NASM; 21 22 #undef __FUNCT__ 23 #define __FUNCT__ "SNESReset_NASM" 24 PetscErrorCode SNESReset_NASM(SNES snes) 25 { 26 SNES_NASM *nasm = (SNES_NASM *)snes->data; 27 PetscErrorCode ierr; 28 PetscInt i; 29 PetscFunctionBegin; 30 for (i=0;i<nasm->n;i++) { 31 if (nasm->x) { ierr = VecDestroy(&nasm->x[i]);CHKERRQ(ierr); } 32 if (nasm->xl){ ierr = VecDestroy(&nasm->xl[i]);CHKERRQ(ierr); } 33 if (nasm->y) { ierr = VecDestroy(&nasm->y[i]);CHKERRQ(ierr); } 34 if (nasm->b) { ierr = VecDestroy(&nasm->b[i]);CHKERRQ(ierr); } 35 36 if (nasm->subsnes) { ierr = SNESDestroy(&nasm->subsnes[i]);CHKERRQ(ierr); } 37 if (nasm->oscatter) { ierr = VecScatterDestroy(&nasm->oscatter[i]);CHKERRQ(ierr); } 38 if (nasm->iscatter) { ierr = VecScatterDestroy(&nasm->iscatter[i]);CHKERRQ(ierr); } 39 if (nasm->gscatter) { ierr = VecScatterDestroy(&nasm->gscatter[i]);CHKERRQ(ierr); } 40 } 41 42 if (nasm->x) {ierr = PetscFree(nasm->x);CHKERRQ(ierr);} 43 if (nasm->xl) {ierr = PetscFree(nasm->xl);CHKERRQ(ierr);} 44 if (nasm->y) {ierr = PetscFree(nasm->y);CHKERRQ(ierr);} 45 if (nasm->b) {ierr = PetscFree(nasm->b);CHKERRQ(ierr);} 46 47 if (nasm->subsnes) {ierr = PetscFree(nasm->subsnes);CHKERRQ(ierr);} 48 if (nasm->oscatter) {ierr = PetscFree(nasm->oscatter);CHKERRQ(ierr);} 49 if (nasm->iscatter) {ierr = PetscFree(nasm->iscatter);CHKERRQ(ierr);} 50 if (nasm->gscatter) {ierr = PetscFree(nasm->gscatter);CHKERRQ(ierr);} 51 52 PetscFunctionReturn(0); 53 } 54 55 #undef __FUNCT__ 56 #define __FUNCT__ "SNESDestroy_NASM" 57 PetscErrorCode SNESDestroy_NASM(SNES snes) 58 { 59 PetscErrorCode ierr; 60 PetscFunctionBegin; 61 ierr = SNESReset_NASM(snes);CHKERRQ(ierr); 62 ierr = PetscFree(snes->data);CHKERRQ(ierr); 63 PetscFunctionReturn(0); 64 } 65 66 #undef __FUNCT__ 67 #define __FUNCT__ "DMGlobalToLocalSubDomainDirichletHook_Private" 68 PetscErrorCode DMGlobalToLocalSubDomainDirichletHook_Private(DM dm,Vec g,InsertMode mode,Vec l,void *ctx) { 69 PetscErrorCode ierr; 70 Vec bcs = (Vec)ctx; 71 PetscFunctionBegin; 72 ierr = VecCopy(bcs,l);CHKERRQ(ierr); 73 PetscFunctionReturn(0); 74 } 75 76 #undef __FUNCT__ 77 #define __FUNCT__ "SNESSetUp_NASM" 78 PetscErrorCode SNESSetUp_NASM(SNES snes) 79 { 80 SNES_NASM *nasm = (SNES_NASM *)snes->data; 81 PetscErrorCode ierr; 82 DM dm,ddm; 83 DM *subdms; 84 PetscInt i; 85 const char *optionsprefix; 86 Vec F; 87 88 PetscFunctionBegin; 89 90 if (!nasm->subsnes) { 91 ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 92 if (dm) { 93 nasm->usesdm = PETSC_TRUE; 94 ierr = DMCreateDomainDecomposition(dm,&nasm->n,PETSC_NULL,PETSC_NULL,PETSC_NULL,&subdms);CHKERRQ(ierr); 95 if (!subdms) { 96 ierr = DMCreateDomainDecompositionDM(dm,"default",&ddm);CHKERRQ(ierr); 97 if (!ddm) { 98 SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_ARG_WRONGSTATE,"DM has no default decomposition defined. Set subsolves manually with SNESNASMSetSubdomains()."); 99 } 100 ierr = SNESSetDM(snes,ddm);CHKERRQ(ierr); 101 ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 102 ierr = DMCreateDomainDecomposition(dm,&nasm->n,PETSC_NULL,PETSC_NULL,PETSC_NULL,&subdms);CHKERRQ(ierr); 103 } 104 if (!subdms) { 105 SETERRQ(((PetscObject)dm)->comm,PETSC_ERR_ARG_WRONGSTATE,"DM has no default decomposition defined. Set subsolves manually with SNESNASMSetSubdomains()."); 106 } 107 ierr = DMCreateDomainDecompositionScatters(dm,nasm->n,subdms,&nasm->iscatter,&nasm->oscatter,&nasm->gscatter);CHKERRQ(ierr); 108 109 ierr = SNESGetOptionsPrefix(snes, &optionsprefix);CHKERRQ(ierr); 110 ierr = PetscMalloc(nasm->n*sizeof(SNES),&nasm->subsnes);CHKERRQ(ierr); 111 112 for (i=0;i<nasm->n;i++) { 113 ierr = SNESCreate(PETSC_COMM_SELF,&nasm->subsnes[i]);CHKERRQ(ierr); 114 ierr = SNESAppendOptionsPrefix(nasm->subsnes[i],optionsprefix);CHKERRQ(ierr); 115 ierr = SNESAppendOptionsPrefix(nasm->subsnes[i],"sub_");CHKERRQ(ierr); 116 ierr = SNESSetDM(nasm->subsnes[i],subdms[i]);CHKERRQ(ierr); 117 ierr = SNESSetFromOptions(nasm->subsnes[i]);CHKERRQ(ierr); 118 ierr = DMDestroy(&subdms[i]);CHKERRQ(ierr); 119 } 120 ierr = PetscFree(subdms);CHKERRQ(ierr); 121 } else { 122 SETERRQ(((PetscObject)snes)->comm,PETSC_ERR_ARG_WRONGSTATE,"Cannot construct local problems automatically without a DM!"); 123 } 124 } else { 125 SETERRQ(((PetscObject)snes)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set subproblems manually if there is no DM!"); 126 } 127 /* allocate the global vectors */ 128 ierr = PetscMalloc(nasm->n*sizeof(Vec),&nasm->x);CHKERRQ(ierr); 129 ierr = PetscMalloc(nasm->n*sizeof(Vec),&nasm->xl);CHKERRQ(ierr); 130 ierr = PetscMalloc(nasm->n*sizeof(Vec),&nasm->y);CHKERRQ(ierr); 131 ierr = PetscMalloc(nasm->n*sizeof(Vec),&nasm->b);CHKERRQ(ierr); 132 133 for (i=0;i<nasm->n;i++) { 134 DM subdm; 135 ierr = SNESGetFunction(nasm->subsnes[i],&F,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); 136 ierr = VecDuplicate(F,&nasm->x[i]);CHKERRQ(ierr); 137 ierr = VecDuplicate(F,&nasm->y[i]);CHKERRQ(ierr); 138 ierr = VecDuplicate(F,&nasm->b[i]);CHKERRQ(ierr); 139 ierr = SNESGetDM(nasm->subsnes[i],&subdm);CHKERRQ(ierr); 140 ierr = DMCreateLocalVector(subdm,&nasm->xl[i]);CHKERRQ(ierr); 141 ierr = DMGlobalToLocalHookAdd(subdm,DMGlobalToLocalSubDomainDirichletHook_Private,PETSC_NULL,nasm->xl[i]);CHKERRQ(ierr); 142 } 143 144 PetscFunctionReturn(0); 145 } 146 147 #undef __FUNCT__ 148 #define __FUNCT__ "SNESSetFromOptions_NASM" 149 PetscErrorCode SNESSetFromOptions_NASM(SNES snes) 150 { 151 PetscErrorCode ierr; 152 DM dm,ddm; 153 char ddm_name[1024]; 154 PCASMType asmtype; 155 const char *const SNESNASMTypes[] = {"NONE","RESTRICT","INTERPOLATE","BASIC","PCASMType","PC_ASM_",0}; 156 PetscBool flg; 157 SNES_NASM *nasm = (SNES_NASM *)snes->data; 158 PetscFunctionBegin; 159 ierr = PetscOptionsHead("Nonlinear Additive Schwartz options");CHKERRQ(ierr); 160 ierr = PetscOptionsEnum("-snes_nasm_type","Type of restriction/extension","",SNESNASMTypes,(PetscEnum)nasm->type,(PetscEnum*)&asmtype,&flg);CHKERRQ(ierr); 161 if (flg) {nasm->type = asmtype;} 162 ierr = PetscOptionsString("-snes_nasm_decomposition", "Name of the DM defining the composition", "SNESSetDM", ddm_name, ddm_name,1024,&flg);CHKERRQ(ierr); 163 ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 164 if (flg) { 165 if (dm) { 166 ierr = DMCreateDomainDecompositionDM(dm, ddm_name, &ddm);CHKERRQ(ierr); 167 if (!ddm) { 168 SETERRQ1(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONGSTATE, "Unknown DM decomposition name %s", ddm_name); 169 } 170 ierr = PetscInfo(snes,"Using domain decomposition DM defined using options database\n");CHKERRQ(ierr); 171 ierr = SNESSetDM(snes,ddm);CHKERRQ(ierr); 172 } else { 173 SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONGSTATE, "No DM to decompose"); 174 } 175 } 176 ierr = PetscOptionsTail();CHKERRQ(ierr); 177 PetscFunctionReturn(0); 178 } 179 180 #undef __FUNCT__ 181 #define __FUNCT__ "SNESView_NASM" 182 PetscErrorCode SNESView_NASM(SNES snes, PetscViewer viewer) 183 { 184 PetscFunctionBegin; 185 186 187 PetscFunctionReturn(0); 188 } 189 190 #undef __FUNCT__ 191 #define __FUNCT__ "SNESNASMSetSubdomains" 192 PetscErrorCode SNESNASMSetSubdomains(SNES snes,PetscInt n,SNES subsnes[],VecScatter iscatter[],VecScatter oscatter[],VecScatter gscatter[]) { 193 PetscErrorCode ierr; 194 PetscErrorCode (*f)(SNES,PetscInt,SNES*,VecScatter*,VecScatter*,VecScatter*); 195 PetscFunctionBegin; 196 ierr = PetscObjectQueryFunction((PetscObject)snes,"SNESNASMSetSubdomains_C",(void (**)(void))&f);CHKERRQ(ierr); 197 ierr = (f)(snes,n,subsnes,iscatter,oscatter,gscatter);CHKERRQ(ierr); 198 PetscFunctionReturn(0); 199 } 200 201 EXTERN_C_BEGIN 202 #undef __FUNCT__ 203 #define __FUNCT__ "SNESNASMSetSubdomains_NASM" 204 PetscErrorCode SNESNASMSetSubdomains_NASM(SNES snes,PetscInt n,SNES subsnes[],VecScatter iscatter[],VecScatter oscatter[],VecScatter gscatter[]) { 205 PetscInt i; 206 PetscErrorCode ierr; 207 SNES_NASM *nasm = (SNES_NASM *)snes->data; 208 PetscFunctionBegin; 209 if (snes->setupcalled)SETERRQ(((PetscObject)snes)->comm,PETSC_ERR_ARG_WRONGSTATE,"SNESNASMSetSubdomains() should be called before calling SNESSetUp()."); 210 211 /* tear down the previously set things */ 212 ierr = SNESReset(snes);CHKERRQ(ierr); 213 214 nasm->n = n; 215 if (oscatter) { 216 for (i=0; i<n; i++) {ierr = PetscObjectReference((PetscObject)oscatter[i]);CHKERRQ(ierr);} 217 } 218 if (iscatter) { 219 for (i=0; i<n; i++) {ierr = PetscObjectReference((PetscObject)iscatter[i]);CHKERRQ(ierr);} 220 } 221 if (gscatter) { 222 for (i=0; i<n; i++) {ierr = PetscObjectReference((PetscObject)gscatter[i]);CHKERRQ(ierr);} 223 } 224 if (oscatter) { 225 ierr = PetscMalloc(n*sizeof(IS),&nasm->oscatter);CHKERRQ(ierr); 226 for (i=0; i<n; i++) { 227 nasm->oscatter[i] = oscatter[i]; 228 } 229 } 230 if (iscatter) { 231 ierr = PetscMalloc(n*sizeof(IS),&nasm->iscatter);CHKERRQ(ierr); 232 for (i=0; i<n; i++) { 233 nasm->iscatter[i] = iscatter[i]; 234 } 235 } 236 if (gscatter) { 237 ierr = PetscMalloc(n*sizeof(IS),&nasm->gscatter);CHKERRQ(ierr); 238 for (i=0; i<n; i++) { 239 nasm->gscatter[i] = gscatter[i]; 240 } 241 } 242 243 if (subsnes) { 244 ierr = PetscMalloc(n*sizeof(SNES),&nasm->subsnes);CHKERRQ(ierr); 245 for (i=0; i<n; i++) { 246 nasm->subsnes[i] = subsnes[i]; 247 } 248 } 249 PetscFunctionReturn(0); 250 } 251 EXTERN_C_END 252 253 #undef __FUNCT__ 254 #define __FUNCT__ "SNESNASMSolveLocal_Private" 255 PetscErrorCode SNESNASMSolveLocal_Private(SNES snes,Vec B,Vec Y,Vec X) { 256 SNES_NASM *nasm = (SNES_NASM *)snes->data; 257 SNES subsnes; 258 PetscInt i; 259 PetscErrorCode ierr; 260 Vec Xlloc,Xl,Bl,Yl; 261 VecScatter iscat,oscat,gscat; 262 DM dm,subdm; 263 PetscFunctionBegin; 264 ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 265 ierr = VecSet(Y,0);CHKERRQ(ierr); 266 for(i=0;i<nasm->n;i++) { 267 subsnes = nasm->subsnes[i]; 268 ierr = SNESGetDM(subsnes,&subdm);CHKERRQ(ierr); 269 iscat = nasm->iscatter[i]; 270 oscat = nasm->oscatter[i]; 271 gscat = nasm->gscatter[i]; 272 ierr = DMSubDomainRestrict(dm,oscat,gscat,subdm);CHKERRQ(ierr); 273 274 /* scatter the solution to the local solution */ 275 Xl = nasm->x[i]; 276 Xlloc = nasm->xl[i]; 277 Yl = nasm->y[i]; 278 279 ierr = VecScatterBegin(oscat,X,Xl,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 280 ierr = VecScatterEnd(oscat,X,Xl,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 281 282 ierr = VecScatterBegin(gscat,X,Xlloc,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 283 ierr = VecScatterEnd(gscat,X,Xlloc,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 284 285 ierr = VecCopy(Xl,Yl);CHKERRQ(ierr); 286 287 if (B) { 288 /* scatter the RHS to the local RHS */ 289 Bl = nasm->b[i]; 290 ierr = VecScatterBegin(oscat,B,Bl,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 291 ierr = VecScatterEnd(oscat,B,Bl,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 292 } else { 293 Bl = PETSC_NULL; 294 } 295 ierr = SNESSolve(subsnes,Bl,Yl);CHKERRQ(ierr); 296 297 ierr = VecAXPY(Yl,-1.0,Xl);CHKERRQ(ierr); 298 299 if (nasm->type == PC_ASM_BASIC) { 300 ierr = VecScatterBegin(oscat,Yl,Y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 301 ierr = VecScatterEnd(oscat,Yl,Y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 302 } else if (nasm->type == PC_ASM_RESTRICT) { 303 ierr = VecScatterBegin(iscat,Yl,Y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 304 ierr = VecScatterEnd(iscat,Yl,Y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 305 } else { 306 SETERRQ(((PetscObject)snes)->comm,PETSC_ERR_ARG_WRONGSTATE,"Only basic and interpolate types are supported for SNESNASM"); 307 } 308 } 309 310 ierr = VecAssemblyBegin(Y);CHKERRQ(ierr); 311 ierr = VecAssemblyEnd(Y);CHKERRQ(ierr); 312 313 ierr = VecAXPY(X,1.0,Y);CHKERRQ(ierr); 314 315 PetscFunctionReturn(0); 316 } 317 318 #undef __FUNCT__ 319 #define __FUNCT__ "SNESSolve_NASM" 320 PetscErrorCode SNESSolve_NASM(SNES snes) 321 { 322 Vec F; 323 Vec X; 324 Vec B; 325 Vec Y; 326 PetscInt i; 327 PetscReal fnorm; 328 PetscErrorCode ierr; 329 SNESNormType normtype; 330 331 PetscFunctionBegin; 332 333 X = snes->vec_sol; 334 Y = snes->vec_sol_update; 335 F = snes->vec_func; 336 B = snes->vec_rhs; 337 338 ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 339 snes->iter = 0; 340 snes->norm = 0.; 341 ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 342 snes->reason = SNES_CONVERGED_ITERATING; 343 ierr = SNESGetNormType(snes, &normtype);CHKERRQ(ierr); 344 if (normtype == SNES_NORM_FUNCTION || normtype == SNES_NORM_INITIAL_ONLY || normtype == SNES_NORM_INITIAL_FINAL_ONLY) { 345 /* compute the initial function and preconditioned update delX */ 346 if (!snes->vec_func_init_set) { 347 ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 348 if (snes->domainerror) { 349 snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 350 PetscFunctionReturn(0); 351 } 352 } else { 353 snes->vec_func_init_set = PETSC_FALSE; 354 } 355 356 /* convergence test */ 357 if (!snes->norm_init_set) { 358 ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F|| */ 359 if (PetscIsInfOrNanReal(fnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"Infinite or not-a-number generated in norm"); 360 } else { 361 fnorm = snes->norm_init; 362 snes->norm_init_set = PETSC_FALSE; 363 } 364 ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 365 snes->iter = 0; 366 snes->norm = fnorm; 367 ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 368 SNESLogConvHistory(snes,snes->norm,0); 369 ierr = SNESMonitor(snes,0,snes->norm);CHKERRQ(ierr); 370 371 /* set parameter for default relative tolerance convergence test */ 372 snes->ttol = fnorm*snes->rtol; 373 374 /* test convergence */ 375 ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 376 if (snes->reason) PetscFunctionReturn(0); 377 } else { 378 ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 379 SNESLogConvHistory(snes,snes->norm,0); 380 ierr = SNESMonitor(snes,0,snes->norm);CHKERRQ(ierr); 381 } 382 383 /* Call general purpose update function */ 384 if (snes->ops->update) { 385 ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 386 } 387 388 for (i = 0; i < snes->max_its; i++) { 389 ierr = SNESNASMSolveLocal_Private(snes,B,Y,X);CHKERRQ(ierr); 390 if (normtype == SNES_NORM_FUNCTION || ((i == snes->max_its - 1) && (normtype == SNES_NORM_INITIAL_FINAL_ONLY || normtype == SNES_NORM_FINAL_ONLY))) { 391 ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 392 if (snes->domainerror) { 393 snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 394 PetscFunctionReturn(0); 395 } 396 ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F|| */ 397 if (PetscIsInfOrNanReal(fnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"Infinite or not-a-number generated in norm"); 398 } 399 /* Monitor convergence */ 400 ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 401 snes->iter = i+1; 402 snes->norm = fnorm; 403 ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 404 SNESLogConvHistory(snes,snes->norm,0); 405 ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr); 406 /* Test for convergence */ 407 if (normtype == SNES_NORM_FUNCTION) ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 408 if (snes->reason) PetscFunctionReturn(0); 409 /* Call general purpose update function */ 410 if (snes->ops->update) { 411 ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 412 } 413 } 414 if (normtype == SNES_NORM_FUNCTION) { 415 if (i == snes->max_its) { 416 ierr = PetscInfo1(snes,"Maximum number of iterations has been reached: %D\n",snes->max_its);CHKERRQ(ierr); 417 if(!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT; 418 } 419 } else { 420 if (!snes->reason) snes->reason = SNES_CONVERGED_ITS; /* NASM is meant to be used as a preconditioner */ 421 } 422 PetscFunctionReturn(0); 423 } 424 425 /*MC 426 SNESNASM - Nonlinear Additive Schwartz 427 428 Level: advanced 429 430 .seealso: SNESCreate(), SNES, SNESSetType(), SNESType (for list of available types) 431 M*/ 432 433 EXTERN_C_BEGIN 434 #undef __FUNCT__ 435 #define __FUNCT__ "SNESCreate_NASM" 436 PetscErrorCode SNESCreate_NASM(SNES snes) 437 { 438 SNES_NASM *nasm; 439 PetscErrorCode ierr; 440 441 PetscFunctionBegin; 442 443 444 ierr = PetscNewLog(snes, SNES_NASM, &nasm);CHKERRQ(ierr); 445 snes->data = (void*)nasm; 446 447 nasm->n = PETSC_DECIDE; 448 nasm->subsnes = 0; 449 nasm->x = 0; 450 nasm->xl = 0; 451 nasm->y = 0; 452 nasm->b = 0; 453 nasm->oscatter = 0; 454 nasm->iscatter = 0; 455 nasm->gscatter = 0; 456 457 nasm->type = PC_ASM_BASIC; 458 459 snes->ops->destroy = SNESDestroy_NASM; 460 snes->ops->setup = SNESSetUp_NASM; 461 snes->ops->setfromoptions = SNESSetFromOptions_NASM; 462 snes->ops->view = SNESView_NASM; 463 snes->ops->solve = SNESSolve_NASM; 464 snes->ops->reset = SNESReset_NASM; 465 466 snes->usesksp = PETSC_FALSE; 467 snes->usespc = PETSC_FALSE; 468 469 if (!snes->tolerancesset) { 470 snes->max_its = 10000; 471 snes->max_funcs = 10000; 472 } 473 474 ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESNASMSetSubdomains_C","SNESNASMSetSubdomains_NASM", 475 SNESNASMSetSubdomains_NASM);CHKERRQ(ierr); 476 477 PetscFunctionReturn(0); 478 } 479 EXTERN_C_END 480