v0.15.0
Loading...
Searching...
No Matches
Static Public Member Functions | List of all members
EshelbianPlasticity::solve_elastic_setup Struct Reference

Static Public Member Functions

static auto setup (EshelbianCore *ep_ptr, TS ts, Vec x, bool set_ts_monitor)
 

Detailed Description

Definition at line 3113 of file EshelbianPlasticity.cpp.

Member Function Documentation

◆ setup()

static auto EshelbianPlasticity::solve_elastic_setup::setup ( EshelbianCore ep_ptr,
TS  ts,
Vec  x,
bool  set_ts_monitor 
)
inlinestatic

Definition at line 3115 of file EshelbianPlasticity.cpp.

3116 {
3117
3118#ifdef ENABLE_PYTHON_BINDING
3119 auto setup_sdf = [&]() {
3120 boost::shared_ptr<ContactOps::SDFPython> sdf_python_ptr;
3121
3122 auto file_exists = [](std::string myfile) {
3123 std::ifstream file(myfile.c_str());
3124 if (file) {
3125 return true;
3126 }
3127 return false;
3128 };
3129
3130 char sdf_file_name[255] = "sdf.py";
3131 CHKERR PetscOptionsGetString(PETSC_NULLPTR, PETSC_NULLPTR, "-sdf_file",
3132 sdf_file_name, 255, PETSC_NULLPTR);
3133
3134 if (file_exists(sdf_file_name)) {
3135 MOFEM_LOG("EP", Sev::inform) << sdf_file_name << " file found";
3136 sdf_python_ptr = boost::make_shared<ContactOps::SDFPython>();
3137 CHKERR sdf_python_ptr->sdfInit(sdf_file_name);
3138 ContactOps::sdfPythonWeakPtr = sdf_python_ptr;
3139 MOFEM_LOG("EP", Sev::inform) << "SdfPython initialized";
3140 } else {
3141 MOFEM_LOG("EP", Sev::warning) << sdf_file_name << " file NOT found";
3142 }
3143 return sdf_python_ptr;
3144 };
3145 auto sdf_python_ptr = setup_sdf();
3146#endif
3147
3148 auto setup_ts_monitor = [&]() {
3149 boost::shared_ptr<TsCtx> ts_ctx;
3151 if (set_ts_monitor) {
3152 CHKERR TSMonitorSet(ts, TsMonitorSet, ts_ctx.get(), PETSC_NULLPTR);
3153 auto monitor_ptr = boost::make_shared<EshelbianMonitor>(*ep_ptr);
3154 ts_ctx->getLoopsMonitor().push_back(
3155 TsCtx::PairNameFEMethodPtr(ep_ptr->elementVolumeName, monitor_ptr));
3156 }
3157 MOFEM_LOG("EP", Sev::inform) << "TS monitor setup";
3158 return std::make_tuple(ts_ctx);
3159 };
3160
3161 auto setup_snes_monitor = [&]() {
3163 SNES snes;
3164 CHKERR TSGetSNES(ts, &snes);
3165 auto snes_ctx = getDMSnesCtx(ep_ptr->dmElastic);
3166 CHKERR SNESMonitorSet(snes,
3167 (MoFEMErrorCode (*)(SNES, PetscInt, PetscReal,
3168 void *))MoFEMSNESMonitorEnergy,
3169 (void *)(snes_ctx.get()), PETSC_NULLPTR);
3170 MOFEM_LOG("EP", Sev::inform) << "SNES monitor setup";
3172 };
3173
3174 auto setup_snes_conergence_test = [&]() {
3176
3177 auto snes_convergence_test = [](SNES snes, PetscInt it, PetscReal xnorm,
3178 PetscReal snorm, PetscReal fnorm,
3179 SNESConvergedReason *reason, void *cctx) {
3181 // EshelbianCore *ep_ptr = (EshelbianCore *)cctx;
3182 CHKERR SNESConvergedDefault(snes, it, xnorm, snorm, fnorm, reason,
3183 PETSC_NULLPTR);
3184
3185 Vec x_update, r;
3186 CHKERR SNESGetSolutionUpdate(snes, &x_update);
3187 CHKERR SNESGetFunction(snes, &r, PETSC_NULLPTR, PETSC_NULLPTR);
3188
3189 // *reason = SNES_CONVERGED_ITERATING;
3190 // if (!it) {
3191 // /* set parameter for default relative tolerance convergence test */
3192 // snes->ttol = fnorm * snes->rtol;
3193 // snes->rnorm0 = fnorm;
3194 // }
3195 // if (PetscIsInfOrNanReal(fnorm)) {
3196 // MOFEM_LOG("EP", Sev::error)
3197 // << "SNES convergence test: function norm is NaN";
3198 // *reason = SNES_DIVERGED_FNORM_NAN;
3199 // } else if (fnorm < snes->abstol && (it || !snes->forceiteration)) {
3200 // MOFEM_LOG("EP", Sev::inform)
3201 // << "SNES convergence test: Converged due to function norm "
3202 // << std::setw(14) << std::setprecision(12) << std::scientific
3203 // << fnorm << " < " << std::setw(14) << std::setprecision(12)
3204 // << std::scientific << snes->abstol;
3205
3206 // PetscCall(PetscInfo(
3207 // snes, "Converged due to function norm %14.12e < %14.12e\n",
3208 // (double)fnorm, (double)snes->abstol));
3209 // *reason = SNES_CONVERGED_FNORM_ABS;
3210 // } else if (snes->nfuncs >= snes->max_funcs && snes->max_funcs >= 0) {
3211 // PetscCall(PetscInfo(
3212 // snes,
3213 // "Exceeded maximum number of function evaluations: %" PetscInt_FMT
3214 // " > %" PetscInt_FMT "\n",
3215 // snes->nfuncs, snes->max_funcs));
3216 // *reason = SNES_DIVERGED_FUNCTION_COUNT;
3217 // }
3218
3219 // if (it && !*reason) {
3220 // if (fnorm <= snes->ttol) {
3221 // PetscCall(PetscInfo(snes,
3222 // "Converged due to function norm %14.12e < "
3223 // "%14.12e (relative tolerance)\n",
3224 // (double)fnorm, (double)snes->ttol));
3225 // *reason = SNES_CONVERGED_FNORM_RELATIVE;
3226 // } else if (snorm < snes->stol * xnorm) {
3227 // PetscCall(PetscInfo(snes,
3228 // "Converged due to small update length: %14.12e "
3229 // "< %14.12e * %14.12e\n",
3230 // (double)snorm, (double)snes->stol,
3231 // (double)xnorm));
3232 // *reason = SNES_CONVERGED_SNORM_RELATIVE;
3233 // } else if (snes->divtol != PETSC_UNLIMITED &&
3234 // (fnorm > snes->divtol * snes->rnorm0)) {
3235 // PetscCall(PetscInfo(snes,
3236 // "Diverged due to increase in function norm: "
3237 // "%14.12e > %14.12e * %14.12e\n",
3238 // (double)fnorm, (double)snes->divtol,
3239 // (double)snes->rnorm0));
3240 // *reason = SNES_DIVERGED_DTOL;
3241 // }
3242 // }
3243
3245 };
3246
3247 // SNES snes;
3248 // CHKERR TSGetSNES(ts, &snes);
3249 // CHKERR SNESSetConvergenceTest(snes, snes_convergence_test, ep_ptr,
3250 // PETSC_NULLPTR);
3251 // MOFEM_LOG("EP", Sev::inform) << "SNES convergence test setup";
3253 };
3254
3255 auto setup_section = [&]() {
3256 PetscSection section_raw;
3257 CHKERR DMGetSection(ep_ptr->dmElastic, &section_raw);
3258 int num_fields;
3259 CHKERR PetscSectionGetNumFields(section_raw, &num_fields);
3260 for (int ff = 0; ff != num_fields; ff++) {
3261 const char *field_name;
3262 CHKERR PetscSectionGetFieldName(section_raw, ff, &field_name);
3263 MOFEM_LOG_C("EP", Sev::inform, "Field %d name %s", ff, field_name);
3264 }
3265 return section_raw;
3266 };
3267
3268 auto set_vector_on_mesh = [&]() {
3270 CHKERR DMoFEMMeshToLocalVector(ep_ptr->dmElastic, x, INSERT_VALUES,
3271 SCATTER_FORWARD);
3272 CHKERR VecGhostUpdateBegin(x, INSERT_VALUES, SCATTER_FORWARD);
3273 CHKERR VecGhostUpdateEnd(x, INSERT_VALUES, SCATTER_FORWARD);
3274 MOFEM_LOG("EP", Sev::inform) << "Vector set on mesh";
3276 };
3277
3278 auto setup_schur_block_solver = [&]() {
3279 MOFEM_LOG("EP", Sev::inform) << "Setting up Schur block solver";
3280 CHKERR TSAppendOptionsPrefix(ts, "elastic_");
3281 CHKERR TSSetFromOptions(ts);
3282 CHKERR TSSetDM(ts, ep_ptr->dmElastic);
3283 // Adding field split solver
3284 boost::shared_ptr<EshelbianCore::SetUpSchur> schur_ptr;
3285 if constexpr (A == AssemblyType::BLOCK_MAT) {
3286 schur_ptr =
3288 CHK_THROW_MESSAGE(schur_ptr->setUp(ts), "setup schur");
3289 }
3290 MOFEM_LOG("EP", Sev::inform) << "Setting up Schur block solver done";
3291 return schur_ptr;
3292 };
3293
3294 // Warning: sequence of construction is not guaranteed for tuple. You have
3295 // to enforce order by proper packaging.
3296
3297#ifdef ENABLE_PYTHON_BINDING
3298 return std::make_tuple(setup_sdf(), setup_ts_monitor(),
3299 setup_snes_monitor(), setup_snes_conergence_test(),
3300 setup_section(), set_vector_on_mesh(),
3301 setup_schur_block_solver());
3302#else
3303 return std::make_tuple(setup_ts_monitor(), setup_snes_monitor(),
3304 setup_snes_conergence_test(), setup_section(),
3305 set_vector_on_mesh(), setup_schur_block_solver());
3306#endif
3307 }
#define MOFEM_LOG_C(channel, severity, format,...)
#define CHK_THROW_MESSAGE(err, msg)
Check and throw MoFEM exception.
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
#define CHKERR
Inline error check.
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
PetscErrorCode DMoFEMMeshToLocalVector(DM dm, Vec l, InsertMode mode, ScatterMode scatter_mode)
set local (or ghosted) vector values on mesh for partition only
Definition DMMoFEM.cpp:514
PetscErrorCode DMMoFEMGetTsCtx(DM dm, MoFEM::TsCtx **ts_ctx)
get MoFEM::TsCtx data structure
Definition DMMoFEM.cpp:1132
#define MOFEM_LOG(channel, severity)
Log.
MoFEM::TsCtx * ts_ctx
const FTensor::Tensor2< T, Dim, Dim > Vec
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
PetscErrorCode TsMonitorSet(TS ts, PetscInt step, PetscReal t, Vec u, void *ctx)
Set monitor for TS solver.
Definition TsCtx.cpp:263
MoFEMErrorCode MoFEMSNESMonitorEnergy(SNES snes, PetscInt its, PetscReal fgnorm, SnesCtx *ctx)
Sens monitor printing residual field by field.
Definition SnesCtx.cpp:648
PetscErrorCode PetscOptionsGetString(PetscOptions *, const char pre[], const char name[], char str[], size_t size, PetscBool *set)
auto getDMSnesCtx(DM dm)
Get SNES context data structure used by DM.
Definition DMMoFEM.hpp:1262
int r
Definition sdf.py:8
constexpr AssemblyType A
constexpr auto field_name
static boost::shared_ptr< SetUpSchur > createSetUpSchur(MoFEM::Interface &m_field, EshelbianCore *ep_core_ptr)
MoFEM::Interface & mField
const std::string elementVolumeName
SmartPetscObj< DM > dmElastic
Elastic problem.
FEMethodsSequence & getLoopsMonitor()
Get the loops to do Monitor object.
Definition TsCtx.hpp:102

The documentation for this struct was generated from the following file: