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