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