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