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 3352 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 3354 of file EshelbianPlasticity.cpp.

3355 {
3356
3357#ifdef ENABLE_PYTHON_BINDING
3358 auto setup_sdf = [&]() {
3359 boost::shared_ptr<ContactOps::SDFPython> sdf_python_ptr;
3360
3361 auto file_exists = [](std::string myfile) {
3362 std::ifstream file(myfile.c_str());
3363 if (file) {
3364 return true;
3365 }
3366 return false;
3367 };
3368
3369 char sdf_file_name[255] = "sdf.py";
3370 CHKERR PetscOptionsGetString(PETSC_NULLPTR, PETSC_NULLPTR, "-sdf_file",
3371 sdf_file_name, 255, PETSC_NULLPTR);
3372
3373 if (file_exists(sdf_file_name)) {
3374 MOFEM_LOG("EP", Sev::inform) << sdf_file_name << " file found";
3375 sdf_python_ptr = boost::make_shared<ContactOps::SDFPython>();
3376 CHKERR sdf_python_ptr->sdfInit(sdf_file_name);
3377 ContactOps::sdfPythonWeakPtr = sdf_python_ptr;
3378 MOFEM_LOG("EP", Sev::inform) << "SdfPython initialized";
3379 } else {
3380 MOFEM_LOG("EP", Sev::warning) << sdf_file_name << " file NOT found";
3381 }
3382 return sdf_python_ptr;
3383 };
3384 auto sdf_python_ptr = setup_sdf();
3385#endif
3386
3387 auto setup_ts_monitor = [&]() {
3388 boost::shared_ptr<TsCtx> ts_ctx;
3390 if (set_ts_monitor) {
3391 CHKERR TSMonitorSet(ts, TsMonitorSet, ts_ctx.get(), PETSC_NULLPTR);
3392 auto monitor_ptr = boost::make_shared<EshelbianMonitor>(*ep_ptr);
3393 ts_ctx->getLoopsMonitor().push_back(
3394 TsCtx::PairNameFEMethodPtr(ep_ptr->elementVolumeName, monitor_ptr));
3395 }
3396 MOFEM_LOG("EP", Sev::inform) << "TS monitor setup";
3397 return std::make_tuple(ts_ctx);
3398 };
3399
3400 auto setup_snes_monitor = [&]() {
3402 SNES snes;
3403 CHKERR TSGetSNES(ts, &snes);
3404 auto snes_ctx = getDMSnesCtx(ep_ptr->dmElastic);
3405 CHKERR SNESMonitorSet(snes,
3406 (MoFEMErrorCode (*)(SNES, PetscInt, PetscReal,
3407 void *))MoFEMSNESMonitorEnergy,
3408 (void *)(snes_ctx.get()), PETSC_NULLPTR);
3409 MOFEM_LOG("EP", Sev::inform) << "SNES monitor setup";
3411 };
3412
3413 auto setup_snes_conergence_test = [&]() {
3415
3416 auto snes_convergence_test = [](SNES snes, PetscInt it, PetscReal xnorm,
3417 PetscReal snorm, PetscReal fnorm,
3418 SNESConvergedReason *reason, void *cctx) {
3420 // EshelbianCore *ep_ptr = (EshelbianCore *)cctx;
3421 CHKERR SNESConvergedDefault(snes, it, xnorm, snorm, fnorm, reason,
3422 PETSC_NULLPTR);
3423
3424 Vec x_update, r;
3425 CHKERR SNESGetSolutionUpdate(snes, &x_update);
3426 CHKERR SNESGetFunction(snes, &r, PETSC_NULLPTR, PETSC_NULLPTR);
3427
3428 // *reason = SNES_CONVERGED_ITERATING;
3429 // if (!it) {
3430 // /* set parameter for default relative tolerance convergence test */
3431 // snes->ttol = fnorm * snes->rtol;
3432 // snes->rnorm0 = fnorm;
3433 // }
3434 // if (PetscIsInfOrNanReal(fnorm)) {
3435 // MOFEM_LOG("EP", Sev::error)
3436 // << "SNES convergence test: function norm is NaN";
3437 // *reason = SNES_DIVERGED_FNORM_NAN;
3438 // } else if (fnorm < snes->abstol && (it || !snes->forceiteration)) {
3439 // MOFEM_LOG("EP", Sev::inform)
3440 // << "SNES convergence test: Converged due to function norm "
3441 // << std::setw(14) << std::setprecision(12) << std::scientific
3442 // << fnorm << " < " << std::setw(14) << std::setprecision(12)
3443 // << std::scientific << snes->abstol;
3444
3445 // PetscCall(PetscInfo(
3446 // snes, "Converged due to function norm %14.12e < %14.12e\n",
3447 // (double)fnorm, (double)snes->abstol));
3448 // *reason = SNES_CONVERGED_FNORM_ABS;
3449 // } else if (snes->nfuncs >= snes->max_funcs && snes->max_funcs >= 0) {
3450 // PetscCall(PetscInfo(
3451 // snes,
3452 // "Exceeded maximum number of function evaluations: %" PetscInt_FMT
3453 // " > %" PetscInt_FMT "\n",
3454 // snes->nfuncs, snes->max_funcs));
3455 // *reason = SNES_DIVERGED_FUNCTION_COUNT;
3456 // }
3457
3458 // if (it && !*reason) {
3459 // if (fnorm <= snes->ttol) {
3460 // PetscCall(PetscInfo(snes,
3461 // "Converged due to function norm %14.12e < "
3462 // "%14.12e (relative tolerance)\n",
3463 // (double)fnorm, (double)snes->ttol));
3464 // *reason = SNES_CONVERGED_FNORM_RELATIVE;
3465 // } else if (snorm < snes->stol * xnorm) {
3466 // PetscCall(PetscInfo(snes,
3467 // "Converged due to small update length: %14.12e "
3468 // "< %14.12e * %14.12e\n",
3469 // (double)snorm, (double)snes->stol,
3470 // (double)xnorm));
3471 // *reason = SNES_CONVERGED_SNORM_RELATIVE;
3472 // } else if (snes->divtol != PETSC_UNLIMITED &&
3473 // (fnorm > snes->divtol * snes->rnorm0)) {
3474 // PetscCall(PetscInfo(snes,
3475 // "Diverged due to increase in function norm: "
3476 // "%14.12e > %14.12e * %14.12e\n",
3477 // (double)fnorm, (double)snes->divtol,
3478 // (double)snes->rnorm0));
3479 // *reason = SNES_DIVERGED_DTOL;
3480 // }
3481 // }
3482
3484 };
3485
3486 // SNES snes;
3487 // CHKERR TSGetSNES(ts, &snes);
3488 // CHKERR SNESSetConvergenceTest(snes, snes_convergence_test, ep_ptr,
3489 // PETSC_NULLPTR);
3490 // MOFEM_LOG("EP", Sev::inform) << "SNES convergence test setup";
3492 };
3493
3494 auto setup_section = [&]() {
3495 PetscSection section_raw;
3496 CHKERR DMGetSection(ep_ptr->dmElastic, &section_raw);
3497 int num_fields;
3498 CHKERR PetscSectionGetNumFields(section_raw, &num_fields);
3499 for (int ff = 0; ff != num_fields; ff++) {
3500 const char *field_name;
3501 CHKERR PetscSectionGetFieldName(section_raw, ff, &field_name);
3502 MOFEM_LOG_C("EP", Sev::inform, "Field %d name %s", ff, field_name);
3503 }
3504 return section_raw;
3505 };
3506
3507 auto set_vector_on_mesh = [&]() {
3509 CHKERR DMoFEMMeshToLocalVector(ep_ptr->dmElastic, x, INSERT_VALUES,
3510 SCATTER_FORWARD);
3511 CHKERR VecGhostUpdateBegin(x, INSERT_VALUES, SCATTER_FORWARD);
3512 CHKERR VecGhostUpdateEnd(x, INSERT_VALUES, SCATTER_FORWARD);
3513 MOFEM_LOG("EP", Sev::inform) << "Vector set on mesh";
3515 };
3516
3517 auto setup_schur_block_solver = [&]() {
3518 MOFEM_LOG("EP", Sev::inform) << "Setting up Schur block solver";
3519 CHKERR TSAppendOptionsPrefix(ts, "elastic_");
3520 CHKERR TSSetFromOptions(ts);
3521 CHKERR TSSetDM(ts, ep_ptr->dmElastic);
3522 // Adding field split solver
3523 boost::shared_ptr<EshelbianCore::SetUpSchur> schur_ptr;
3524 if constexpr (A == AssemblyType::BLOCK_MAT) {
3525 schur_ptr =
3527 CHK_THROW_MESSAGE(schur_ptr->setUp(ts), "setup schur");
3528 }
3529 MOFEM_LOG("EP", Sev::inform) << "Setting up Schur block solver done";
3530 return schur_ptr;
3531 };
3532
3533 // Warning: sequence of construction is not guaranteed for tuple. You have
3534 // to enforce order by proper packaging.
3535
3536#ifdef ENABLE_PYTHON_BINDING
3537 return std::make_tuple(setup_sdf(), setup_ts_monitor(),
3538 setup_snes_monitor(), setup_snes_conergence_test(),
3539 setup_section(), set_vector_on_mesh(),
3540 setup_schur_block_solver());
3541#else
3542 return std::make_tuple(setup_ts_monitor(), setup_snes_monitor(),
3543 setup_snes_conergence_test(), setup_section(),
3544 set_vector_on_mesh(), setup_schur_block_solver());
3545#endif
3546 }
#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:652
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
[Define dimension]
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: