Loading [MathJax]/extensions/AMSmath.js
v0.14.0
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Modules Pages
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 2808 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 2810 of file EshelbianPlasticity.cpp.

2811  {
2812 
2813 #ifdef PYTHON_SDF
2814  auto setup_sdf = [&]() {
2815  boost::shared_ptr<ContactOps::SDFPython> sdf_python_ptr;
2816 
2817  auto file_exists = [](std::string myfile) {
2818  std::ifstream file(myfile.c_str());
2819  if (file) {
2820  return true;
2821  }
2822  return false;
2823  };
2824 
2825  char sdf_file_name[255] = "sdf.py";
2826  CHKERR PetscOptionsGetString(PETSC_NULL, PETSC_NULL, "-sdf_file",
2827  sdf_file_name, 255, PETSC_NULL);
2828 
2829  if (file_exists(sdf_file_name)) {
2830  MOFEM_LOG("EP", Sev::inform) << sdf_file_name << " file found";
2831  sdf_python_ptr = boost::make_shared<ContactOps::SDFPython>();
2832  CHKERR sdf_python_ptr->sdfInit(sdf_file_name);
2833  ContactOps::sdfPythonWeakPtr = sdf_python_ptr;
2834  } else {
2835  MOFEM_LOG("EP", Sev::warning) << sdf_file_name << " file NOT found";
2836  }
2837  return sdf_python_ptr;
2838  };
2839  auto sdf_python_ptr = setup_sdf();
2840 #endif
2841 
2842  auto setup_ts_monitor = [&]() {
2843  boost::shared_ptr<TsCtx> ts_ctx;
2844  DMMoFEMGetTsCtx(ep_ptr->dmElastic, ts_ctx);
2845  if (set_ts_monitor) {
2846  CHKERR TSMonitorSet(ts, TsMonitorSet, ts_ctx.get(), PETSC_NULL);
2847  auto monitor_ptr = boost::make_shared<EshelbianMonitor>(*ep_ptr);
2848  ts_ctx->getLoopsMonitor().push_back(
2849  TsCtx::PairNameFEMethodPtr(ep_ptr->elementVolumeName, monitor_ptr));
2850  }
2851  return std::make_tuple(ts_ctx);
2852  };
2853 
2854  auto setup_snes_monitor = [&]() {
2856  SNES snes;
2857  CHKERR TSGetSNES(ts, &snes);
2858  PetscViewerAndFormat *vf;
2859  CHKERR PetscViewerAndFormatCreate(PETSC_VIEWER_STDOUT_WORLD,
2860  PETSC_VIEWER_DEFAULT, &vf);
2861  CHKERR SNESMonitorSet(
2862  snes,
2863  (MoFEMErrorCode(*)(SNES, PetscInt, PetscReal,
2864  void *))SNESMonitorFields,
2865  vf, (MoFEMErrorCode(*)(void **))PetscViewerAndFormatDestroy);
2867  };
2868 
2869  auto setup_section = [&]() {
2870  PetscSection section_raw;
2871  CHKERR DMGetSection(ep_ptr->dmElastic, &section_raw);
2872  int num_fields;
2873  CHKERR PetscSectionGetNumFields(section_raw, &num_fields);
2874  for (int ff = 0; ff != num_fields; ff++) {
2875  const char *field_name;
2876  CHKERR PetscSectionGetFieldName(section_raw, ff, &field_name);
2877  MOFEM_LOG_C("EP", Sev::inform, "Field %d name %s", ff, field_name);
2878  }
2879  return section_raw;
2880  };
2881 
2882  auto set_vector_on_mesh = [&]() {
2884  CHKERR DMoFEMMeshToLocalVector(ep_ptr->dmElastic, x, INSERT_VALUES,
2885  SCATTER_FORWARD);
2886  CHKERR VecGhostUpdateBegin(x, INSERT_VALUES, SCATTER_FORWARD);
2887  CHKERR VecGhostUpdateEnd(x, INSERT_VALUES, SCATTER_FORWARD);
2889  };
2890 
2891  auto setup_schur_block_solver = [&]() {
2892  CHKERR TSAppendOptionsPrefix(ts, "elastic_");
2893  CHKERR TSSetFromOptions(ts);
2894  CHKERR TSSetDM(ts, ep_ptr->dmElastic);
2895  // Adding field split solver
2896  boost::shared_ptr<EshelbianCore::SetUpSchur> schur_ptr;
2897  if constexpr (A == AssemblyType::BLOCK_MAT) {
2898  schur_ptr =
2900  CHK_THROW_MESSAGE(schur_ptr->setUp(ts), "setup schur");
2901  }
2902  return schur_ptr;
2903  };
2904 
2905  // Warning: sequence of construction is not guaranteed for tuple. You have
2906  // to enforce order by proper packaging.
2907 
2908 #ifdef PYTHON_SDF
2909  return std::make_tuple(setup_sdf(), setup_ts_monitor(),
2910  setup_snes_monitor(), setup_section(),
2911  set_vector_on_mesh(), setup_schur_block_solver());
2912 #else
2913  return std::make_tuple(setup_ts_monitor(), setup_snes_monitor(),
2914  setup_section(), set_vector_on_mesh(),
2915  setup_schur_block_solver());
2916 #endif
2917  }

The documentation for this struct was generated from the following file:
MoFEMFunctionReturnHot
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:460
EshelbianCore::SetUpSchur::createSetUpSchur
static boost::shared_ptr< SetUpSchur > createSetUpSchur(MoFEM::Interface &m_field, EshelbianCore *ep_core_ptr)
Definition: SetUpSchurImpl.cpp:578
CHK_THROW_MESSAGE
#define CHK_THROW_MESSAGE(err, msg)
Check and throw MoFEM exception.
Definition: definitions.h:609
MoFEM::TsCtx::getLoopsMonitor
FEMethodsSequence & getLoopsMonitor()
Get the loops to do Monitor object.
Definition: TsCtx.hpp:102
MoFEM::BLOCK_MAT
@ BLOCK_MAT
Definition: FormsIntegrators.hpp:107
MoFEM::Exceptions::MoFEMErrorCode
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
Definition: Exceptions.hpp:56
A
constexpr AssemblyType A
Definition: operators_tests.cpp:30
MoFEM::DMoFEMMeshToLocalVector
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:523
ts_ctx
MoFEM::TsCtx * ts_ctx
Definition: level_set.cpp:1932
MoFEM::TsMonitorSet
PetscErrorCode TsMonitorSet(TS ts, PetscInt step, PetscReal t, Vec u, void *ctx)
Set monitor for TS solver.
Definition: TsCtx.cpp:263
EshelbianCore::dmElastic
SmartPetscObj< DM > dmElastic
Elastic problem.
Definition: EshelbianCore.hpp:96
CHKERR
#define CHKERR
Inline error check.
Definition: definitions.h:548
MOFEM_LOG_C
#define MOFEM_LOG_C(channel, severity, format,...)
Definition: LogManager.hpp:311
MoFEM::PairNameFEMethodPtr
Definition: AuxPETSc.hpp:12
EshelbianCore::elementVolumeName
const std::string elementVolumeName
Definition: EshelbianCore.hpp:113
EshelbianCore::mField
MoFEM::Interface & mField
Definition: EshelbianCore.hpp:84
field_name
constexpr auto field_name
Definition: poisson_2d_homogeneous.cpp:13
MOFEM_LOG
#define MOFEM_LOG(channel, severity)
Log.
Definition: LogManager.hpp:308
MoFEM::PetscOptionsGetString
PetscErrorCode PetscOptionsGetString(PetscOptions *, const char pre[], const char name[], char str[], size_t size, PetscBool *set)
Definition: DeprecatedPetsc.hpp:172
MoFEM::DMMoFEMGetTsCtx
PetscErrorCode DMMoFEMGetTsCtx(DM dm, MoFEM::TsCtx **ts_ctx)
get MoFEM::TsCtx data structure
Definition: DMMoFEM.cpp:1142
MoFEMFunctionBeginHot
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:453