1 /* Defines the basic SNES object */ 2 #include <../src/snes/impls/fas/fasimpls.h> 3 4 /*MC 5 Full Approximation Scheme nonlinear multigrid solver. 6 7 The nonlinear problem is solved via the repeated application of nonlinear preconditioners and coarse-grid corrections 8 9 .seealso: SNESCreate(), SNES, SNESSetType(), SNESType (for list of available types) 10 M*/ 11 12 extern PetscErrorCode SNESDestroy_FAS(SNES snes); 13 extern PetscErrorCode SNESSetUp_FAS(SNES snes); 14 extern PetscErrorCode SNESSetFromOptions_FAS(SNES snes); 15 extern PetscErrorCode SNESView_FAS(SNES snes, PetscViewer viewer); 16 extern PetscErrorCode SNESSolve_FAS(SNES snes); 17 extern PetscErrorCode SNESReset_FAS(SNES snes); 18 19 EXTERN_C_BEGIN 20 21 #undef __FUNCT__ 22 #define __FUNCT__ "SNESCreate_FAS" 23 PetscErrorCode SNESCreate_FAS(SNES snes) 24 { 25 SNES_FAS * fas; 26 PetscErrorCode ierr; 27 28 PetscFunctionBegin; 29 snes->ops->destroy = SNESDestroy_FAS; 30 snes->ops->setup = SNESSetUp_FAS; 31 snes->ops->setfromoptions = SNESSetFromOptions_FAS; 32 snes->ops->view = SNESView_FAS; 33 snes->ops->solve = SNESSolve_FAS; 34 snes->ops->reset = SNESReset_FAS; 35 36 ierr = PetscNewLog(snes, SNES_FAS, &fas);CHKERRQ(ierr); 37 snes->data = (void*) fas; 38 fas->level = 0; 39 fas->presmooth = PETSC_NULL; 40 fas->postsmooth = PETSC_NULL; 41 fas->next = PETSC_NULL; 42 fas->interpolate = PETSC_NULL; 43 fas->restrct = PETSC_NULL; 44 45 PetscFunctionReturn(0); 46 } 47 EXTERN_C_END 48 49 #undef __FUNCT__ 50 #define __FUNCT__ "SNESFASGetLevels" 51 PetscErrorCode SNESFASGetLevels(SNES snes, PetscInt * levels) { 52 SNES_FAS * fas = (SNES_FAS *)snes->data; 53 PetscFunctionBegin; 54 *levels = fas->level; 55 PetscFunctionReturn(0); 56 } 57 58 #undef __FUNCT__ 59 #define __FUNCT__ "SNESFASGetSNES" 60 PetscErrorCode SNESFASGetSNES(SNES snes, PetscInt level, SNES * lsnes) { 61 SNES_FAS * fas = (SNES_FAS *)snes->data; 62 PetscInt levels = fas->level; 63 PetscInt i; 64 PetscFunctionBegin; 65 *lsnes = snes; 66 if (fas->level < level) { 67 SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "SNESFASGetSNES should only be called on a finer SNESFAS instance than the level."); 68 } 69 if (level > levels - 1) { 70 SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Level %d doesn't exist in the SNESFAS."); 71 } 72 for (i = fas->level; i > level; i--) { 73 *lsnes = fas->next; 74 fas = (SNES_FAS *)(*lsnes)->data; 75 } 76 if (fas->level != level) SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "SNESFASGetSNES didn't return the right level!"); 77 PetscFunctionReturn(0); 78 } 79 80 #undef __FUNCT__ 81 #define __FUNCT__ "SNESFASSetLevels" 82 PetscErrorCode SNESFASSetLevels(SNES snes, PetscInt levels, MPI_Comm * comms) { 83 PetscErrorCode ierr; 84 PetscInt i; 85 SNES_FAS * fas = (SNES_FAS *)snes->data; 86 MPI_Comm comm; 87 88 PetscFunctionBegin; 89 comm = PETSC_COMM_WORLD; 90 /* user has changed the number of levels; reset */ 91 ierr = SNESReset(snes);CHKERRQ(ierr); 92 /* destroy any coarser levels if necessary */ 93 if (fas->next) SNESDestroy(&fas->next);CHKERRQ(ierr); 94 /* setup the finest level */ 95 for (i = levels - 1; i >= 0; i--) { 96 if (comms) comm = comms[i]; 97 fas->level = i; 98 fas->levels = levels; 99 if (i > 0) { 100 ierr = SNESCreate(comm, &fas->next);CHKERRQ(ierr); 101 ierr = SNESSetType(fas->next, SNESFAS);CHKERRQ(ierr); 102 fas = (SNES_FAS *)fas->next->data; 103 } 104 } 105 PetscFunctionReturn(0); 106 } 107 108 #undef __FUNCT__ 109 #define __FUNCT__ "SNESFASSetInterpolation" 110 PetscErrorCode SNESFASSetInterpolation(SNES snes, PetscInt level, Mat mat) { 111 SNES_FAS * fas = (SNES_FAS *)snes->data; 112 PetscInt top_level = fas->level,i; 113 114 PetscFunctionBegin; 115 if (level > top_level) 116 SETERRQ1(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Bad level number %d in SNESFASSetInterpolation", level); 117 /* get to the correct level */ 118 for (i = fas->level; i > level; i--) { 119 fas = (SNES_FAS *)fas->next->data; 120 } 121 if (fas->level != level) 122 SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONG, "Inconsistent level labelling in SNESFASSetInterpolation"); 123 fas->interpolate = mat; 124 PetscFunctionReturn(0); 125 } 126 127 #undef __FUNCT__ 128 #define __FUNCT__ "SNESFASSetRestriction" 129 PetscErrorCode SNESFASSetRestriction(SNES snes, PetscInt level, Mat mat) { 130 SNES_FAS * fas = (SNES_FAS *)snes->data; 131 PetscInt top_level = fas->level,i; 132 133 PetscFunctionBegin; 134 if (level > top_level) 135 SETERRQ1(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Bad level number %d in SNESFASSetRestriction", level); 136 /* get to the correct level */ 137 for (i = fas->level; i > level; i--) { 138 fas = (SNES_FAS *)fas->next->data; 139 } 140 if (fas->level != level) 141 SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONG, "Inconsistent level labelling in SNESFASSetRestriction"); 142 fas->restrct = mat; 143 PetscFunctionReturn(0); 144 } 145 146 #undef __FUNCT__ 147 #define __FUNCT__ "SNESFASSetRScale" 148 PetscErrorCode SNESFASSetRScale(SNES snes, PetscInt level, Vec rscale) { 149 SNES_FAS * fas = (SNES_FAS *)snes->data; 150 PetscInt top_level = fas->level,i; 151 152 PetscFunctionBegin; 153 if (level > top_level) 154 SETERRQ1(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Bad level number %d in SNESFASSetRestriction", level); 155 /* get to the correct level */ 156 for (i = fas->level; i > level; i--) { 157 fas = (SNES_FAS *)fas->next->data; 158 } 159 if (fas->level != level) 160 SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONG, "Inconsistent level labelling in SNESFASSetRestriction"); 161 fas->rscale = rscale; 162 PetscFunctionReturn(0); 163 } 164 165 #undef __FUNCT__ 166 #define __FUNCT__ "SNESReset_FAS" 167 PetscErrorCode SNESReset_FAS(SNES snes) 168 { 169 PetscErrorCode ierr = 0; 170 SNES_FAS * fas = (SNES_FAS *)snes->data; 171 172 PetscFunctionBegin; 173 /* destroy local data created in SNESSetup_FAS */ 174 if (fas->presmooth) ierr = SNESDestroy(&fas->presmooth);CHKERRQ(ierr); 175 if (fas->postsmooth) ierr = SNESDestroy(&fas->postsmooth);CHKERRQ(ierr); 176 if (fas->interpolate) ierr = MatDestroy(&fas->interpolate);CHKERRQ(ierr); 177 if (fas->restrct) ierr = MatDestroy(&fas->restrct);CHKERRQ(ierr); 178 if (fas->rscale) ierr = VecDestroy(&fas->rscale);CHKERRQ(ierr); 179 180 181 /* recurse -- reset should NOT destroy the structures -- destroy should destroy the structures recursively */ 182 if (fas->next) ierr = SNESReset(fas->next);CHKERRQ(ierr); 183 if (snes->work) {ierr = VecDestroyVecs(snes->nwork,&snes->work);CHKERRQ(ierr);} 184 PetscFunctionReturn(0); 185 } 186 187 #undef __FUNCT__ 188 #define __FUNCT__ "SNESDestroy_FAS" 189 PetscErrorCode SNESDestroy_FAS(SNES snes) 190 { 191 SNES_FAS * fas = (SNES_FAS *)snes->data; 192 PetscErrorCode ierr; 193 194 PetscFunctionBegin; 195 /* recursively resets and then destroys */ 196 ierr = SNESReset_FAS(snes);CHKERRQ(ierr); 197 if (fas->next) ierr = SNESDestroy(&fas->next);CHKERRQ(ierr); 198 ierr = PetscFree(fas);CHKERRQ(ierr); 199 PetscFunctionReturn(0); 200 } 201 202 #undef __FUNCT__ 203 #define __FUNCT__ "SNESSetUp_FAS" 204 PetscErrorCode SNESSetUp_FAS(SNES snes) 205 { 206 SNES_FAS *fas = (SNES_FAS *) snes->data; 207 PetscErrorCode ierr; 208 209 PetscFunctionBegin; 210 ierr = SNESDefaultGetWork(snes, 4);CHKERRQ(ierr); /* the four work vectors are used to transfer stuff BACK */ 211 /* gets the solver ready for solution */ 212 if (snes->dm) { 213 /* construct EVERYTHING from the DM -- including the progressive set of smoothers */ 214 if (fas->next) { 215 /* for now -- assume the DM and the evaluation functions have been set externally */ 216 if (!fas->interpolate) { 217 if (!fas->next->dm) SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONGSTATE, "All levels must be assigned a DM"); 218 /* set the interpolation and restriction from the DM */ 219 ierr = DMGetInterpolation(fas->next->dm, snes->dm, &fas->interpolate, &fas->rscale);CHKERRQ(ierr); 220 fas->restrct = fas->interpolate; 221 } 222 /* TODO LATER: Preconditioner setup goes here */ 223 } 224 } 225 if (fas->next) { 226 /* gotta set up the solution vector for this to work */ 227 ierr = VecDuplicate(fas->rscale, &fas->next->vec_sol);CHKERRQ(ierr); 228 ierr = SNESSetUp(fas->next);CHKERRQ(ierr); 229 } 230 /* got to set them all up at once */ 231 PetscFunctionReturn(0); 232 } 233 234 #undef __FUNCT__ 235 #define __FUNCT__ "SNESSetFromOptions_FAS" 236 PetscErrorCode SNESSetFromOptions_FAS(SNES snes) 237 { 238 239 PetscErrorCode ierr; 240 241 PetscFunctionBegin; 242 /* types for the pre and postsmoothers */ 243 ierr = PetscOptionsHead("SNESFAS Options");CHKERRQ(ierr); 244 ierr = PetscOptionsTail();CHKERRQ(ierr); 245 PetscFunctionReturn(0); 246 } 247 248 #undef __FUNCT__ 249 #define __FUNCT__ "SNESView_FAS" 250 PetscErrorCode SNESView_FAS(SNES snes, PetscViewer viewer) 251 { 252 SNES_FAS *fas = (SNES_FAS *) snes->data; 253 PetscBool iascii; 254 PetscErrorCode ierr; 255 PetscInt levels = fas->levels; 256 PetscInt i; 257 258 PetscFunctionBegin; 259 ierr = PetscTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr); 260 if (iascii) { 261 ierr = PetscViewerASCIIPrintf(viewer, "FAS, levels = %d\n", fas->levels);CHKERRQ(ierr); 262 ierr = PetscViewerASCIIPushTab(viewer); 263 for (i = levels - 1; i >= 0; i--) { 264 ierr = PetscViewerASCIIPrintf(viewer, "level: %d\n", fas->level);CHKERRQ(ierr); 265 if (fas->presmooth) { 266 ierr = PetscViewerASCIIPrintf(viewer, "pre-smoother on level %D\n", fas->level);CHKERRQ(ierr); 267 ierr = PetscViewerASCIIPushTab(viewer); 268 ierr = SNESView(fas->presmooth, viewer);CHKERRQ(ierr); 269 ierr = PetscViewerASCIIPopTab(viewer); 270 } else { 271 ierr = PetscViewerASCIIPrintf(viewer, "no pre-smoother on level %D\n", fas->level);CHKERRQ(ierr); 272 } 273 if (fas->postsmooth) { 274 ierr = PetscViewerASCIIPrintf(viewer, "post-smoother on level %D\n", fas->level);CHKERRQ(ierr); 275 ierr = PetscViewerASCIIPushTab(viewer); 276 ierr = SNESView(fas->postsmooth, viewer);CHKERRQ(ierr); 277 ierr = PetscViewerASCIIPopTab(viewer); 278 } else { 279 ierr = PetscViewerASCIIPrintf(viewer, "no post-smoother on level %D\n", fas->level);CHKERRQ(ierr); 280 } 281 if (fas->next) fas = (SNES_FAS *)fas->next->data; 282 } 283 ierr = PetscViewerASCIIPopTab(viewer); 284 } else { 285 SETERRQ1(((PetscObject)snes)->comm,PETSC_ERR_SUP,"Viewer type %s not supported for SNESFAS",((PetscObject)viewer)->type_name); 286 } 287 PetscFunctionReturn(0); 288 } 289 290 /* 291 292 Defines the FAS cycle as: 293 294 fine problem: F(x) = 0 295 coarse problem: F^c(x) = b^c 296 297 b^c = F^c(I^c_fx^f - I^c_fF(x)) 298 299 correction: 300 301 x = x + I(x^c - Rx) 302 303 */ 304 305 #undef __FUNCT__ 306 #define __FUNCT__ "FASCycle_Private" 307 PetscErrorCode FASCycle_Private(SNES snes, Vec F, Vec X) { 308 309 PetscErrorCode ierr; 310 Vec X_c, Xo_c, F_c, B_c; 311 SNES_FAS * fas = (SNES_FAS *)snes->data; 312 313 PetscFunctionBegin; 314 /* pre-smooth -- just update using the pre-smoother */ 315 316 if (fas->next) { 317 X_c = fas->next->work[0]; 318 Xo_c = fas->next->work[1]; 319 F_c = fas->next->work[2]; 320 B_c = fas->next->work[3]; 321 /* project to coarse */ 322 ierr = MatRestrict(fas->restrct, X, X_c);CHKERRQ(ierr); 323 ierr = MatRestrict(fas->restrct, F, F_c);CHKERRQ(ierr); 324 ierr = VecPointwiseMult(Xo_c, fas->rscale, Xo_c);CHKERRQ(ierr); 325 ierr = VecPointwiseMult(F_c, fas->rscale, F_c);CHKERRQ(ierr); 326 /* solve the coarse problem corresponding to F^c(x^c) = b^c = F^c(Rx) - RF(x) */ 327 ierr = SNESComputeFunction(fas->next, Xo_c, B_c);CHKERRQ(ierr); /* B_c = F(X_c) */ 328 ierr = VecAXPY(B_c, -1.0, F_c);CHKERRQ(ierr); /* B_c = F(X_c) - F_C */ 329 ierr = VecCopy(Xo_c, X_c);CHKERRQ(ierr); 330 ierr = SNESSolve(fas->next, B_c, X_c);CHKERRQ(ierr); 331 /* correct as x <- x + I(x^c - Rx)*/ 332 ierr = VecAXPY(X_c, -1.0, Xo_c);CHKERRQ(ierr); 333 ierr = MatInterpolate(fas->interpolate, X_c, F);CHKERRQ(ierr); 334 ierr = VecAXPY(X, -1.0, F);CHKERRQ(ierr); 335 } 336 ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr); 337 /* post-smooth -- just update using the post-smoother */ 338 339 PetscFunctionReturn(0); 340 } 341 342 #undef __FUNCT__ 343 #define __FUNCT__ "SNESSolve_FAS" 344 345 PetscErrorCode SNESSolve_FAS(SNES snes) 346 { 347 PetscErrorCode ierr; 348 Vec X, F; 349 PetscFunctionBegin; 350 X = snes->vec_sol; 351 F = snes->vec_func; 352 ierr = FASCycle_Private(snes, snes->vec_func, snes->vec_sol);CHKERRQ(ierr); 353 PetscFunctionReturn(0); 354 } 355