v0.14.0
Loading...
Searching...
No Matches
Public Member Functions | Private Attributes | List of all members
ContactPlasticity::MMonitor Struct Reference
Inheritance diagram for ContactPlasticity::MMonitor:
[legend]
Collaboration diagram for ContactPlasticity::MMonitor:
[legend]

Public Member Functions

 MMonitor (SmartPetscObj< DM > &dm, MoFEM::Interface &m_field, boost::shared_ptr< OpContactTools::CommonData > common_data_ptr, boost::shared_ptr< PostProcVolumeOnRefinedMesh > post_proc_fe, std::tuple< SmartPetscObj< Vec >, SmartPetscObj< VecScatter > > ux_scatter, std::tuple< SmartPetscObj< Vec >, SmartPetscObj< VecScatter > > uy_scatter, std::tuple< SmartPetscObj< Vec >, SmartPetscObj< VecScatter > > uz_scatter)
 
MoFEMErrorCode preProcess ()
 function is run at the beginning of loop More...
 
MoFEMErrorCode operator() ()
 function is run for every finite element More...
 
MoFEMErrorCode postProcess ()
 function is run at the end of loop More...
 
- Public Member Functions inherited from MoFEM::FEMethod
MoFEMErrorCode query_interface (boost::typeindex::type_index type_index, UnknownInterface **iface) const
 
 FEMethod ()=default
 
auto getFEName () const
 get finite element name More...
 
auto getDataDofsPtr () const
 
auto getDataVectorDofsPtr () const
 
const FieldEntity_vector_viewgetDataFieldEnts () const
 
boost::shared_ptr< FieldEntity_vector_view > & getDataFieldEntsPtr () const
 
auto & getRowFieldEnts () const
 
auto & getRowFieldEntsPtr () const
 
auto & getColFieldEnts () const
 
auto & getColFieldEntsPtr () const
 
auto getRowDofsPtr () const
 
auto getColDofsPtr () const
 
auto getNumberOfNodes () const
 
EntityHandle getFEEntityHandle () const
 
MoFEMErrorCode getNodeData (const std::string field_name, VectorDouble &data, const bool reset_dofs=true)
 
- Public Member Functions inherited from MoFEM::BasicMethod
MoFEMErrorCode query_interface (boost::typeindex::type_index type_index, UnknownInterface **iface) const
 
 BasicMethod ()
 
virtual ~BasicMethod ()=default
 
int getNinTheLoop () const
 get number of evaluated element in the loop More...
 
int getLoopSize () const
 get loop size More...
 
auto getLoHiFERank () const
 Get lo and hi processor rank of iterated entities. More...
 
auto getLoFERank () const
 Get upper rank in loop for iterating elements. More...
 
auto getHiFERank () const
 Get upper rank in loop for iterating elements. More...
 
unsigned int getFieldBitNumber (std::string field_name) const
 
MoFEMErrorCode copyBasicMethod (const BasicMethod &basic)
 Copy data from other base method to this base method. More...
 
virtual MoFEMErrorCode preProcess ()
 function is run at the beginning of loop More...
 
virtual MoFEMErrorCode operator() ()
 function is run for every finite element More...
 
virtual MoFEMErrorCode postProcess ()
 function is run at the end of loop More...
 
boost::weak_ptr< CacheTuplegetCacheWeakPtr () const
 Get the cache weak ptr object. More...
 
- Public Member Functions inherited from MoFEM::KspMethod
MoFEMErrorCode query_interface (boost::typeindex::type_index type_index, UnknownInterface **iface) const
 
 KspMethod ()
 
virtual ~KspMethod ()=default
 
MoFEMErrorCode copyKsp (const KspMethod &ksp)
 copy data form another method More...
 
- Public Member Functions inherited from MoFEM::PetscData
MoFEMErrorCode query_interface (boost::typeindex::type_index type_index, UnknownInterface **iface) const
 
 PetscData ()
 
virtual ~PetscData ()=default
 
MoFEMErrorCode copyPetscData (const PetscData &petsc_data)
 
- Public Member Functions inherited from MoFEM::UnknownInterface
virtual MoFEMErrorCode query_interface (boost::typeindex::type_index type_index, UnknownInterface **iface) const =0
 
template<class IFACE >
MoFEMErrorCode registerInterface (bool error_if_registration_failed=true)
 Register interface. More...
 
template<class IFACE >
MoFEMErrorCode getInterface (IFACE *&iface) const
 Get interface refernce to pointer of interface. More...
 
template<class IFACE >
MoFEMErrorCode getInterface (IFACE **const iface) const
 Get interface pointer to pointer of interface. More...
 
template<class IFACE , typename boost::enable_if< boost::is_pointer< IFACE >, int >::type = 0>
IFACE getInterface () const
 Get interface pointer to pointer of interface. More...
 
template<class IFACE , typename boost::enable_if< boost::is_reference< IFACE >, int >::type = 0>
IFACE getInterface () const
 Get reference to interface. More...
 
template<class IFACE >
IFACE * getInterface () const
 Function returning pointer to interface. More...
 
virtual ~UnknownInterface ()=default
 
- Public Member Functions inherited from MoFEM::SnesMethod
MoFEMErrorCode query_interface (boost::typeindex::type_index type_index, UnknownInterface **iface) const
 
 SnesMethod ()
 
virtual ~SnesMethod ()=default
 
MoFEMErrorCode copySnes (const SnesMethod &snes)
 Copy snes data. More...
 
- Public Member Functions inherited from MoFEM::TSMethod
MoFEMErrorCode query_interface (boost::typeindex::type_index type_index, UnknownInterface **iface) const
 
 TSMethod ()
 
virtual ~TSMethod ()=default
 
MoFEMErrorCode copyTs (const TSMethod &ts)
 Copy TS solver data. More...
 

Private Attributes

SmartPetscObj< DM > dM
 
MoFEM::InterfacemField
 
boost::shared_ptr< OpContactTools::CommonDatacommonDataPtr
 
boost::shared_ptr< PostProcVolumeOnRefinedMeshpostProcFe
 
boost::shared_ptr< DomainSideElevolSideFe
 
boost::shared_ptr< PostProcFaceOnRefinedMeshpostProcFeSkeleton
 
boost::shared_ptr< PostProcFaceOnRefinedMeshpostProcFeSkin
 
std::tuple< SmartPetscObj< Vec >, SmartPetscObj< VecScatter > > uXScatter
 
std::tuple< SmartPetscObj< Vec >, SmartPetscObj< VecScatter > > uYScatter
 
std::tuple< SmartPetscObj< Vec >, SmartPetscObj< VecScatter > > uZScatter
 

Additional Inherited Members

- Public Types inherited from MoFEM::KspMethod
enum  KSPContext { CTX_SETFUNCTION , CTX_OPERATORS , CTX_KSPNONE }
 pass information about context of KSP/DM for with finite element is computed More...
 
- Public Types inherited from MoFEM::PetscData
enum  DataContext {
  CTX_SET_NONE = 0 , CTX_SET_F = 1 << 0 , CTX_SET_A = 1 << 1 , CTX_SET_B = 1 << 2 ,
  CTX_SET_X = 1 << 3 , CTX_SET_X_T = 1 << 4 , CTX_SET_X_TT = 1 << 6 , CTX_SET_TIME = 1 << 7
}
 
using Switches = std::bitset< 8 >
 
- Public Types inherited from MoFEM::SnesMethod
enum  SNESContext { CTX_SNESSETFUNCTION , CTX_SNESSETJACOBIAN , CTX_SNESNONE }
 
- Public Types inherited from MoFEM::TSMethod
enum  TSContext {
  CTX_TSSETRHSFUNCTION , CTX_TSSETRHSJACOBIAN , CTX_TSSETIFUNCTION , CTX_TSSETIJACOBIAN ,
  CTX_TSTSMONITORSET , CTX_TSNONE
}
 
- Static Public Member Functions inherited from MoFEM::UnknownInterface
static MoFEMErrorCode getLibVersion (Version &version)
 Get library version. More...
 
static MoFEMErrorCode getFileVersion (moab::Interface &moab, Version &version)
 Get database major version. More...
 
static MoFEMErrorCode setFileVersion (moab::Interface &moab, Version version=Version(MoFEM_VERSION_MAJOR, MoFEM_VERSION_MINOR, MoFEM_VERSION_BUILD))
 Get database major version. More...
 
static MoFEMErrorCode getInterfaceVersion (Version &version)
 Get database major version. More...
 
- Public Attributes inherited from MoFEM::FEMethod
std::string feName
 Name of finite element. More...
 
boost::shared_ptr< const NumeredEntFiniteElementnumeredEntFiniteElementPtr
 
boost::function< bool(FEMethod *fe_method_ptr)> exeTestHook
 Tet if element to skip element. More...
 
- Public Attributes inherited from MoFEM::BasicMethod
int nInTheLoop
 number currently of processed method More...
 
int loopSize
 local number oe methods to process More...
 
std::pair< int, int > loHiFERank
 Llo and hi processor rank of iterated entities. More...
 
int rAnk
 processor rank More...
 
int sIze
 number of processors in communicator More...
 
const RefEntity_multiIndexrefinedEntitiesPtr
 container of mofem dof entities More...
 
const RefElement_multiIndexrefinedFiniteElementsPtr
 container of mofem finite element entities More...
 
const ProblemproblemPtr
 raw pointer to problem More...
 
const Field_multiIndexfieldsPtr
 raw pointer to fields container More...
 
const FieldEntity_multiIndexentitiesPtr
 raw pointer to container of field entities More...
 
const DofEntity_multiIndexdofsPtr
 raw pointer container of dofs More...
 
const FiniteElement_multiIndexfiniteElementsPtr
 raw pointer to container finite elements More...
 
const EntFiniteElement_multiIndexfiniteElementsEntitiesPtr
 
const FieldEntityEntFiniteElementAdjacencyMap_multiIndexadjacenciesPtr
 
boost::function< MoFEMErrorCode()> preProcessHook
 Hook function for pre-processing. More...
 
boost::function< MoFEMErrorCode()> postProcessHook
 Hook function for post-processing. More...
 
boost::function< MoFEMErrorCode()> operatorHook
 Hook function for operator. More...
 
boost::movelib::unique_ptr< boolvecAssembleSwitch
 
boost::movelib::unique_ptr< boolmatAssembleSwitch
 
boost::weak_ptr< CacheTuplecacheWeakPtr
 
- Public Attributes inherited from MoFEM::KspMethod
KSPContext ksp_ctx
 Context. More...
 
KSP ksp
 KSP solver. More...
 
Vec & ksp_f
 
Mat & ksp_A
 
Mat & ksp_B
 
- Public Attributes inherited from MoFEM::PetscData
Switches data_ctx
 
Vec f
 
Mat A
 
Mat B
 
Vec x
 
Vec x_t
 
Vec x_tt
 
- Public Attributes inherited from MoFEM::SnesMethod
SNESContext snes_ctx
 
SNES snes
 snes solver More...
 
Vec & snes_x
 state vector More...
 
Vec & snes_f
 residual More...
 
Mat & snes_A
 jacobian matrix More...
 
Mat & snes_B
 preconditioner of jacobian matrix More...
 
- Public Attributes inherited from MoFEM::TSMethod
TS ts
 time solver More...
 
TSContext ts_ctx
 
PetscInt ts_step
 time step number More...
 
PetscReal ts_a
 shift for U_t (see PETSc Time Solver) More...
 
PetscReal ts_aa
 shift for U_tt shift for U_tt More...
 
PetscReal ts_t
 time More...
 
PetscReal ts_dt
 time step size More...
 
Vec & ts_u
 state vector More...
 
Vec & ts_u_t
 time derivative of state vector More...
 
Vec & ts_u_tt
 second time derivative of state vector More...
 
Vec & ts_F
 residual vector More...
 
Mat & ts_A
 
Mat & ts_B
 Preconditioner for ts_A. More...
 
- Static Public Attributes inherited from MoFEM::PetscData
static constexpr Switches CtxSetNone = PetscData::Switches(CTX_SET_NONE)
 
static constexpr Switches CtxSetF = PetscData::Switches(CTX_SET_F)
 
static constexpr Switches CtxSetA = PetscData::Switches(CTX_SET_A)
 
static constexpr Switches CtxSetB = PetscData::Switches(CTX_SET_B)
 
static constexpr Switches CtxSetX = PetscData::Switches(CTX_SET_X)
 
static constexpr Switches CtxSetX_T = PetscData::Switches(CTX_SET_X_T)
 
static constexpr Switches CtxSetX_TT = PetscData::Switches(CTX_SET_X_TT)
 
static constexpr Switches CtxSetTime = PetscData::Switches(CTX_SET_TIME)
 

Detailed Description

Definition at line 246 of file multifield_plasticity.cpp.

Constructor & Destructor Documentation

◆ MMonitor()

ContactPlasticity::MMonitor::MMonitor ( SmartPetscObj< DM > &  dm,
MoFEM::Interface m_field,
boost::shared_ptr< OpContactTools::CommonData common_data_ptr,
boost::shared_ptr< PostProcVolumeOnRefinedMesh post_proc_fe,
std::tuple< SmartPetscObj< Vec >, SmartPetscObj< VecScatter > >  ux_scatter,
std::tuple< SmartPetscObj< Vec >, SmartPetscObj< VecScatter > >  uy_scatter,
std::tuple< SmartPetscObj< Vec >, SmartPetscObj< VecScatter > >  uz_scatter 
)
inline

Definition at line 248 of file multifield_plasticity.cpp.

255 : dM(dm), mField(m_field), commonDataPtr(common_data_ptr),
256 postProcFe(post_proc_fe), uXScatter(ux_scatter),
257 uYScatter(uy_scatter), uZScatter(uz_scatter) {
258 volSideFe = boost::make_shared<MoFEM::DomainSideEle>(mField);
259 volSideFe->getOpPtrVector().push_back(
261 volSideFe->getOpPtrVector().push_back(
264 boost::make_shared<PostProcFaceOnRefinedMesh>(mField);
265 postProcFeSkeleton->generateReferenceElementMesh();
266 postProcFeSkeleton->getOpPtrVector().push_back(
267 new OpCalculateJumpOnSkeleton("SKELETON_LAMBDA", commonDataPtr,
268 volSideFe));
269 postProcFeSkeleton->getOpPtrVector().push_back(
271 "SKELETON_LAMBDA", postProcFeSkeleton->postProcMesh,
272 postProcFeSkeleton->mapGaussPts, commonDataPtr));
273
274 postProcFeSkin = boost::make_shared<PostProcFaceOnRefinedMesh>(mField);
275 CHKERR postProcFeSkin->generateReferenceElementMesh();
276 CHKERR postProcFeSkin->addFieldValuesPostProc("U", "DISPLACEMENT");
278 postProcFeSkin->getUserPolynomialBase() =
279 boost::shared_ptr<BaseFunction>(new BubbleTauPolynomialBase());
280
282 postProcFeSkin->addFieldValuesPostProc("MESH_NODE_POSITIONS");
283
284 if (data.with_plasticity) {
285 CHKERR postProcFeSkin->addFieldValuesPostProcOnSkin("TAU");
286 CHKERR postProcFeSkin->addFieldValuesPostProcOnSkin("EP");
287 }
288
289 if (data.with_contact) {
290 CHKERR postProcFeSkin->addFieldValuesPostProcOnSkin("SIGMA");
291 }
292 };
#define CHKERR
Inline error check.
Definition: definitions.h:535
boost::shared_ptr< DomainSideEle > volSideFe
std::tuple< SmartPetscObj< Vec >, SmartPetscObj< VecScatter > > uZScatter
std::tuple< SmartPetscObj< Vec >, SmartPetscObj< VecScatter > > uYScatter
boost::shared_ptr< OpContactTools::CommonData > commonDataPtr
boost::shared_ptr< PostProcVolumeOnRefinedMesh > postProcFe
boost::shared_ptr< PostProcFaceOnRefinedMesh > postProcFeSkin
boost::shared_ptr< PostProcFaceOnRefinedMesh > postProcFeSkeleton
std::tuple< SmartPetscObj< Vec >, SmartPetscObj< VecScatter > > uXScatter
static ProblemData data

Member Function Documentation

◆ operator()()

MoFEMErrorCode ContactPlasticity::MMonitor::operator() ( )
inlinevirtual

function is run for every finite element

It is used to calculate element local matrices and assembly. It can be used for post-processing.

Reimplemented from MoFEM::BasicMethod.

Definition at line 315 of file multifield_plasticity.cpp.

315{ return 0; }

◆ postProcess()

MoFEMErrorCode ContactPlasticity::MMonitor::postProcess ( )
inlinevirtual

function is run at the end of loop

It is used to assembly matrices and vectors, calculating global variables, f.e. total internal energy, ect.

Iterating over dofs: Example1 iterating over dofs in row by name of the field for(IT_GET_FEROW_BY_NAME_DOFS_FOR_LOOP(this,"DISPLACEMENT",it)) { ... }

Reimplemented from MoFEM::BasicMethod.

Definition at line 317 of file multifield_plasticity.cpp.

317 {
319
320 auto make_vtk = [&]() {
322 if (data.no_output || ts_step % data.output_freq != 0)
325 CHKERR postProcFe->writeFile("out_plastic_" +
326 boost::lexical_cast<std::string>(ts_step) +
327 ".h5m");
329 };
330
331 auto save_skeleton = [&]() {
333
335 CHKERR postProcFeSkeleton->writeFile(
336 "out_skeleton_" + boost::lexical_cast<std::string>(ts_step) +
337 ".h5m");
339 };
340
341 auto save_skin = [&]() {
343 if (data.no_output || ts_step % data.output_freq != 0)
345
347 CHKERR postProcFeSkin->writeFile(
348 "out_skin_" + boost::lexical_cast<std::string>(ts_step) + ".h5m");
350 };
351
352 double max_disp, min_disp;
353
354 auto print_max_min = [&](auto &tuple, const std::string msg) {
356 CHKERR VecScatterBegin(std::get<1>(tuple), ts_u, std::get<0>(tuple),
357 INSERT_VALUES, SCATTER_FORWARD);
358 CHKERR VecScatterEnd(std::get<1>(tuple), ts_u, std::get<0>(tuple),
359 INSERT_VALUES, SCATTER_FORWARD);
360 CHKERR VecMax(std::get<0>(tuple), PETSC_NULL, &max_disp);
361 CHKERR VecMin(std::get<0>(tuple), PETSC_NULL, &min_disp);
362 PetscPrintf(PETSC_COMM_WORLD, "%s time %3.4e min %3.4e max %3.4e\n",
363 msg.c_str(), ts_t, min_disp, max_disp);
365 };
366
367 double gauss_area;
368
369 if (data.with_contact) {
370 auto &l_reactions = commonDataPtr->bodyReactions;
371 auto &lgauss_pts = commonDataPtr->stateArrayArea;
372
373 int roller_nb = (*cache).rollerDataVec.size();
374 l_reactions.resize(roller_nb * 3);
375 std::fill(l_reactions.begin(), l_reactions.end(), 0);
376
378
379 if (true) {
380 std::vector<double> gauss_pts(2, 0);
381 auto dm = mField.getInterface<Simple>()->getDM();
382 CHKERR MPIU_Allreduce(lgauss_pts.data(), gauss_pts.data(), 2,
383 MPIU_REAL, MPIU_SUM,
384 PetscObjectComm((PetscObject)dm));
385 MOFEM_LOG_C("WORLD", Sev::inform,
386 "\t TS %d Total contact area (ref): %g / %g", ts_step,
387 gauss_pts[0], gauss_pts[1]);
388 gauss_area = gauss_pts[0];
389 lgauss_pts[0] = lgauss_pts[1] = 0;
390
391 vector<double> g_reactions(l_reactions.size() * 3, 0);
392 CHKERR MPIU_Allreduce(l_reactions.data(), g_reactions.data(),
393 roller_nb * 3, MPIU_REAL, MPIU_SUM,
394 PetscObjectComm((PetscObject)dm));
395 for (int dd = 0; dd != roller_nb; ++dd)
396 MOFEM_LOG_C("WORLD", Sev::inform,
397 "\t Body %d reactions Time %g Force X: %3.4e Y: "
398 "%3.4e Z: %3.4e",
399 (*cache).rollerDataVec[dd].iD, ts_t,
400 g_reactions[3 * dd + 0], g_reactions[3 * dd + 1],
401 g_reactions[3 * dd + 2]);
402 }
403 }
404
405 if (ts_step % data.output_freq == 0 && !data.no_output) {
406 std::ostringstream ostrm, ostrm2;
407 ostrm << "out_debug_" << ts_step << ".vtk";
408 ostrm2 << "out_debug_" << ts_step << ".h5m";
410 if (mField.get_comm_size() == 1)
411 CHKERR data.moab_debug->write_file(ostrm.str().c_str(), "VTK", "");
412 else
413 CHKERR data.moab_debug->write_file(ostrm2.str().c_str(), "MOAB",
414 "PARALLEL=WRITE_PART");
415
416 data.moab_debug->delete_mesh();
417 }
418 }
419
421 auto &vec = commonDataPtr->reactionsVec;
422 CHKERR VecZeroEntries(vec);
423 CHKERR VecGhostUpdateBegin(vec, INSERT_VALUES, SCATTER_FORWARD);
424 CHKERR VecGhostUpdateEnd(vec, INSERT_VALUES, SCATTER_FORWARD);
427 Range ents;
428 CHKERR mField.get_moab().get_entities_by_dimension(0, 0, ents, true);
429
430 ParallelComm *pcomm =
431 ParallelComm::get_pcomm(&mField.get_moab(), MYPCOMM_INDEX);
432 CHKERR pcomm->reduce_tags(commonDataPtr->reactionTag, MPI_SUM, ents);
433
434 CHKERR make_vtk();
435 CHKERR save_skin();
436
437 // clear reactions tag
438 double def_value = 0.;
439 CHKERR mField.get_moab().tag_clear_data(commonDataPtr->reactionTag,
440 ents, &def_value);
441
442 // CHKERR VecGhostUpdateBegin(vec, ADD_VALUES, SCATTER_REVERSE);
443 // CHKERR VecGhostUpdateEnd(vec, ADD_VALUES, SCATTER_REVERSE);
444 CHKERR VecAssemblyBegin(vec);
445 CHKERR VecAssemblyEnd(vec);
446 CHKERR VecGhostUpdateBegin(vec, ADD_VALUES, SCATTER_REVERSE);
447 CHKERR VecGhostUpdateEnd(vec, ADD_VALUES, SCATTER_REVERSE);
448 CHKERR VecGhostUpdateBegin(vec, INSERT_VALUES, SCATTER_FORWARD);
449 CHKERR VecGhostUpdateEnd(vec, INSERT_VALUES, SCATTER_FORWARD);
450
451 const double *react_ptr;
452 CHKERR VecGetArrayRead(vec, &react_ptr);
453
454 PetscPrintf(PETSC_COMM_WORLD,
455 "Reactions time %3.4e X: %3.4e Y: %3.4e Z: %3.4e \n", ts_t,
456 react_ptr[0], react_ptr[1], react_ptr[2]);
457
458 CHKERR VecRestoreArrayRead(vec, &react_ptr);
459 } else {
460 CHKERR make_vtk();
461 CHKERR save_skin();
462 }
463 // FIXME: this does not work properly
465 CHKERR save_skeleton();
466
469 auto &m_field = postProcFe->mField;
470 data.moab_debug->delete_mesh();
471 if (m_field.get_comm_rank() == 0) {
472 std::ostringstream ost;
473 for (auto &rol : (*cache).rollerDataVec) {
474 EntityHandle vertex;
475 VectorDouble roller_coords = rol.BodyDispScaled;
476 roller_coords += rol.originCoords;
477
478 CHKERR data.moab_debug->create_vertex(&roller_coords(0), vertex);
479 CHKERR rol.saveBasicDataOnTag(*data.moab_debug, vertex);
480 }
481
482 ost << "out_roller_" << ts_step << ".vtk";
483 CHKERR data.moab_debug->write_file(ost.str().c_str(), "VTK", "");
484 data.moab_debug->delete_mesh();
485 }
486 }
487
488 if (data.save_restarts && ts_step % data.output_freq == 0) {
489
490 CHKERR VecGhostUpdateBegin(ts_u, INSERT_VALUES, SCATTER_FORWARD);
491 CHKERR VecGhostUpdateEnd(ts_u, INSERT_VALUES, SCATTER_FORWARD);
492
493 string rest_name =
494 "restart_" + to_string(ts_step) + "_" + to_string(ts_t) + ".dat";
495 PetscViewer viewer;
496 CHKERR PetscViewerBinaryOpen(PETSC_COMM_WORLD, rest_name.c_str(),
497 FILE_MODE_WRITE, &viewer);
498 CHKERR PetscViewerBinarySetSkipInfo(viewer, PETSC_TRUE);
499 VecView(ts_u, viewer);
500 PetscViewerDestroy(&viewer);
501 }
502
503 CHKERR print_max_min(uXScatter, "Ux");
504 CHKERR print_max_min(uYScatter, "Uy");
505 CHKERR print_max_min(uZScatter, "Uz");
506
507 switch (data.atom_test_nb) {
508 case 1: {
509 if (ts_step == 10) {
510 const double *react_ptr;
511 CHKERR VecGetArrayRead(commonDataPtr->reactionsVec, &react_ptr);
512 if (fabs(react_ptr[0] + 0.12519621572) > 1e-6)
513 SETERRQ(PETSC_COMM_SELF, MOFEM_ATOM_TEST_INVALID,
514 "atom test diverged!");
515 CHKERR VecRestoreArrayRead(commonDataPtr->reactionsVec, &react_ptr);
516 }
517 break;
518 }
519 case 2: {
520 if (ts_step == 5) {
521 const double *react_ptr;
522 CHKERR VecGetArrayRead(commonDataPtr->reactionsVec, &react_ptr);
523 // hertz solution 0.25*(4/3)*(1/(1-0.3*0.3))*(10^0.5)*0.1^(3/2)
524 if ((0.03663003663 - react_ptr[2]) / 0.03663003663 > 1e-2)
525 SETERRQ(PETSC_COMM_SELF, MOFEM_ATOM_TEST_INVALID,
526 "atom test diverged!");
527 CHKERR VecRestoreArrayRead(commonDataPtr->reactionsVec, &react_ptr);
528 }
529 break;
530 }
531 case 3: {
532 if (ts_t == 1) {
533 double min_disp;
534 CHKERR VecMin(std::get<0>(uXScatter), PETSC_NULL, &min_disp);
535 if (fabs((min_disp + 2.70) / 2.70) > 5e-2)
536 SETERRQ(PETSC_COMM_SELF, MOFEM_ATOM_TEST_INVALID,
537 "atom test diverged!");
538 }
539 break;
540 }
541 case 4: {
542 if (ts_step == 10) {
543 const double *react_ptr;
544 CHKERR VecGetArrayRead(commonDataPtr->reactionsVec, &react_ptr);
545 // Reaction Force at first yield Johnson = Contact Mechanics = Eq6.10
546 // Fy = (yield_stress/Em)^2 * yield_stress *(Pi^3 * 1.6^3 * R^2)/6
547 // Em = E /(1-poisson_ratio^2)
548 if (abs(0.056090728 - (react_ptr[2] * 4.0)) / 0.056090728 > 5.0e-2)
549 SETERRQ(PETSC_COMM_SELF, MOFEM_ATOM_TEST_INVALID,
550 "atom test diverged!");
551 }
552 case 5: {
553 if (ts_step == 14) {
554 const double *react_ptr;
555 //comparison with data from miehe.dat file
556 CHKERR VecGetArrayRead(commonDataPtr->reactionsVec, &react_ptr);
557 if (abs(19030.7 + react_ptr[2]) / 19030.7 > 1e-3)
558 SETERRQ(PETSC_COMM_SELF, MOFEM_ATOM_TEST_INVALID,
559 "atom test diverged!");
560 }
561 break;
562 }
563 break;
564 }
565
566 default:
567 break;
568 }
569
571 }
#define MOFEM_LOG_C(channel, severity, format,...)
Definition: LogManager.hpp:311
void simple(double P1[], double P2[], double P3[], double c[], const int N)
Definition: acoustic.cpp:69
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:447
#define MYPCOMM_INDEX
default communicator number PCOMM
Definition: definitions.h:215
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:346
@ MOFEM_ATOM_TEST_INVALID
Definition: definitions.h:40
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:416
PetscErrorCode DMoFEMLoopFiniteElements(DM dm, const char fe_name[], MoFEM::FEMethod *method, CacheTupleWeakPtr cache_ptr=CacheTupleSharedPtr())
Executes FEMethod for finite elements in DM.
Definition: DMMoFEM.cpp:572
const Tensor2_symmetric_Expr< const ddTensor0< T, Dim, i, j >, typename promote< T, double >::V, Dim, i, j > dd(const Tensor0< T * > &a, const Index< i, Dim > index1, const Index< j, Dim > index2, const Tensor1< int, Dim > &d_ijk, const Tensor1< double, Dim > &d_xyz)
Definition: ddTensor0.hpp:33
boost::shared_ptr< DomainEle > calc_reactions
boost::shared_ptr< BoundaryEle > update_roller
virtual int get_comm_size() const =0
virtual moab::Interface & get_moab()=0
Simple interface for fast problem set-up.
Definition: Simple.hpp:27
PetscReal ts_t
time
PetscInt ts_step
time step number
Vec & ts_u
state vector
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface refernce to pointer of interface.

◆ preProcess()

MoFEMErrorCode ContactPlasticity::MMonitor::preProcess ( )
inlinevirtual

function is run at the beginning of loop

It is used to zeroing matrices and vectors, calculation of shape functions on reference element, preprocessing boundary conditions, etc.

Reimplemented from MoFEM::BasicMethod.

Definition at line 294 of file multifield_plasticity.cpp.

294 {
296 CHKERR TSGetTimeStep(ts, &(commonDataPtr->timeStepSize));
297 if (!ts_step || ts_step == data.restart_step) {
298 for (auto &roller : (*cache).rollerDataVec) {
299 if (roller.methodOpForRollerPosition) {
300 roller.BodyDispScaled.clear();
302 this, roller.methodOpForRollerPosition, roller.BodyDispScaled);
303 }
304 if (roller.methodOpForRollerDirection) {
305 roller.BodyDirectionScaled.clear();
307 this, roller.methodOpForRollerDirection,
308 roller.BodyDirectionScaled);
309 }
310 }
311 }
312
314 }
static MoFEMErrorCode applyScale(const FEMethod *fe, boost::ptr_vector< MethodForForceScaling > &methods_op, VectorDouble &nf)
TS ts
time solver

Member Data Documentation

◆ commonDataPtr

boost::shared_ptr<OpContactTools::CommonData> ContactPlasticity::MMonitor::commonDataPtr
private

Definition at line 576 of file multifield_plasticity.cpp.

◆ dM

SmartPetscObj<DM> ContactPlasticity::MMonitor::dM
private

Definition at line 574 of file multifield_plasticity.cpp.

◆ mField

MoFEM::Interface& ContactPlasticity::MMonitor::mField
private

Definition at line 575 of file multifield_plasticity.cpp.

◆ postProcFe

boost::shared_ptr<PostProcVolumeOnRefinedMesh> ContactPlasticity::MMonitor::postProcFe
private

Definition at line 577 of file multifield_plasticity.cpp.

◆ postProcFeSkeleton

boost::shared_ptr<PostProcFaceOnRefinedMesh> ContactPlasticity::MMonitor::postProcFeSkeleton
private

Definition at line 580 of file multifield_plasticity.cpp.

◆ postProcFeSkin

boost::shared_ptr<PostProcFaceOnRefinedMesh> ContactPlasticity::MMonitor::postProcFeSkin
private

Definition at line 581 of file multifield_plasticity.cpp.

◆ uXScatter

std::tuple<SmartPetscObj<Vec>, SmartPetscObj<VecScatter> > ContactPlasticity::MMonitor::uXScatter
private

Definition at line 583 of file multifield_plasticity.cpp.

◆ uYScatter

std::tuple<SmartPetscObj<Vec>, SmartPetscObj<VecScatter> > ContactPlasticity::MMonitor::uYScatter
private

Definition at line 584 of file multifield_plasticity.cpp.

◆ uZScatter

std::tuple<SmartPetscObj<Vec>, SmartPetscObj<VecScatter> > ContactPlasticity::MMonitor::uZScatter
private

Definition at line 585 of file multifield_plasticity.cpp.

◆ volSideFe

boost::shared_ptr<DomainSideEle> ContactPlasticity::MMonitor::volSideFe
private

Definition at line 579 of file multifield_plasticity.cpp.


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