v0.14.0
Public Member Functions | Public Attributes | List of all members
BoneRemodeling::MonitorPostProc Struct Reference

#include <users_modules/bone_remodelling/src/Remodeling.hpp>

Inheritance diagram for BoneRemodeling::MonitorPostProc:
[legend]
Collaboration diagram for BoneRemodeling::MonitorPostProc:
[legend]

Public Member Functions

 MonitorPostProc (MoFEM::Interface &m_field, Remodeling::CommonData &common_data)
 
MoFEMErrorCode preProcess ()
 
MoFEMErrorCode operator() ()
 
MoFEMErrorCode postProcess ()
 

Public Attributes

MoFEM::InterfacemField
 
Remodeling::CommonDatacommonData
 
PostProcVolumeOnRefinedMesh postProc
 
PostProcVolumeOnRefinedMesh postProcElastic
 
PetscBool mass_postproc
 
PetscBool equilibrium_flg
 
double rate
 
bool iNit
 
int pRT
 

Detailed Description

Element used to post-process results at each time step

Definition at line 553 of file Remodeling.hpp.

Constructor & Destructor Documentation

◆ MonitorPostProc()

MonitorPostProc::MonitorPostProc ( MoFEM::Interface m_field,
Remodeling::CommonData common_data 
)
Examples
Remodeling.hpp.

Definition at line 1387 of file Remodeling.cpp.

1389  : mField(m_field), commonData(common_data), postProc(m_field),
1390  postProcElastic(m_field),
1391  // densityMaps(m_field),
1392  iNit(false) {}

Member Function Documentation

◆ operator()()

MoFEMErrorCode MonitorPostProc::operator() ( )
Examples
Remodeling.hpp.

Definition at line 1399 of file Remodeling.cpp.

1399  {
1402 }

◆ postProcess()

MoFEMErrorCode MonitorPostProc::postProcess ( )
Examples
Remodeling.hpp.

Definition at line 1404 of file Remodeling.cpp.

1404  {
1406 
1407  PetscBool mass_postproc = PETSC_FALSE;
1408  PetscBool equilibrium_flg = PETSC_FALSE;
1409  PetscBool analysis_mesh_flg = PETSC_FALSE;
1410  rate = 1;
1411 
1412  CHKERR PetscOptionsBegin(PETSC_COMM_WORLD, "", "Bone remodeling post-process",
1413  "none");
1414 
1415  CHKERR PetscOptionsBool(
1416  "-mass_postproc",
1417  "is used for testing, calculates mass and energy at each time step", "",
1418  mass_postproc, &mass_postproc, PETSC_NULL);
1419 
1420  CHKERR PetscOptionsBool("-analysis_mesh",
1421  "saves analysis mesh at each time step", "",
1422  analysis_mesh_flg, &analysis_mesh_flg, PETSC_NULL);
1423 
1424  CHKERR PetscOptionsScalar(
1425  "-equilibrium_stop_rate",
1426  "is used to stop calculations when equilibium state is achieved", "",
1427  rate, &rate, &equilibrium_flg);
1428  ierr = PetscOptionsEnd();
1429  CHKERRQ(ierr);
1430 
1431  if (!iNit) {
1432 
1433  CHKERR PetscOptionsBegin(PETSC_COMM_WORLD, "",
1434  "Bone remodeling post-process", "none");
1435 
1436  pRT = 1;
1437  CHKERR PetscOptionsInt("-my_output_prt",
1438  "frequncy how often results are dumped on hard disk",
1439  "", 1, &pRT, NULL);
1440  ierr = PetscOptionsEnd();
1441  CHKERRQ(ierr);
1442 
1444  CHKERR postProc.addFieldValuesPostProc("DISPLACEMENTS");
1445  CHKERR postProc.addFieldValuesPostProc("MESH_NODE_POSITIONS");
1446  if (!commonData.less_post_proc) {
1448  }
1449  // CHKERR postProc.addFieldValuesPostProc("RHO");
1450  // CHKERR postProc.addFieldValuesGradientPostProc("RHO");
1451 
1452  {
1453  auto &list_of_operators = postProc.getOpPtrVector();
1454 
1455  list_of_operators.push_back(
1457  list_of_operators.push_back(new OpCalculateScalarFieldGradient<3>(
1458  "RHO", commonData.data.gradRhoPtr));
1459  // Get displacement gradient at integration points
1460  list_of_operators.push_back(new OpCalculateVectorFieldGradient<3, 3>(
1461  "DISPLACEMENTS", commonData.data.gradDispPtr));
1462  list_of_operators.push_back(new OpCalculateStress(commonData));
1463  list_of_operators.push_back(new OpPostProcStress(
1465  }
1466 
1469  CHKERR postProcElastic.addFieldValuesPostProc("MESH_NODE_POSITIONS");
1471  std::map<int, NonlinearElasticElement::BlockData>::iterator sit =
1472  commonData.elasticPtr->setOfBlocks.begin();
1473  for (; sit != commonData.elasticPtr->setOfBlocks.end(); sit++) {
1476  "DISPLACEMENTS", sit->second, postProcElastic.commonData));
1477  }
1478 
1479  iNit = true;
1480  }
1481 
1482  if (mass_postproc) {
1483  Vec mass_vec;
1484  Vec energ_vec;
1485  int mass_vec_ghost[] = {0};
1486  CHKERR VecCreateGhost(PETSC_COMM_WORLD, (!mField.get_comm_rank()) ? 1 : 0,
1487  1, 1, mass_vec_ghost, &mass_vec);
1488 
1489  CHKERR VecZeroEntries(mass_vec);
1490  CHKERR VecDuplicate(mass_vec, &energ_vec);
1491 
1492  Remodeling::Fe energyCalc(mField);
1493  CHKERR addHOOpsVol("MESH_NODE_POSITIONS", energyCalc, true, false, false,
1494  true);
1495  energyCalc.getOpPtrVector().push_back(
1497  energyCalc.getOpPtrVector().push_back(new OpCalculateScalarFieldGradient<3>(
1498  "RHO", commonData.data.gradRhoPtr));
1499  energyCalc.getOpPtrVector().push_back(
1500  new OpCalculateVectorFieldGradient<3, 3>("DISPLACEMENTS",
1502  energyCalc.getOpPtrVector().push_back(new OpCalculateStress(commonData));
1503  energyCalc.getOpPtrVector().push_back(
1504  new OpMassAndEnergyCalculation("RHO", commonData, energ_vec, mass_vec));
1505  CHKERR DMoFEMLoopFiniteElements(commonData.dm, "FE_REMODELLING",
1506  &energyCalc);
1507 
1508  CHKERR VecAssemblyBegin(mass_vec);
1509  CHKERR VecAssemblyEnd(mass_vec);
1510  double mass_sum;
1511  CHKERR VecSum(mass_vec, &mass_sum);
1512 
1513  CHKERR VecAssemblyBegin(energ_vec);
1514  CHKERR VecAssemblyEnd(energ_vec);
1515  double energ_sum;
1516  CHKERR VecSum(energ_vec, &energ_sum);
1517 
1518  CHKERR PetscPrintf(PETSC_COMM_WORLD,
1519  "Time: %g Mass: %6.5g Elastic energy: %6.5g \n", ts_t,
1520  mass_sum, energ_sum);
1521  CHKERR VecDestroy(&mass_vec);
1522  CHKERR VecDestroy(&energ_vec);
1523 
1524  if (equilibrium_flg) {
1525  double equilibrium_rate =
1526  fabs(energ_sum - commonData.cUrrent_psi) / energ_sum;
1527  double equilibrium_mass_rate =
1528  fabs(mass_sum - commonData.cUrrent_mass) / mass_sum;
1529  if (equilibrium_rate < rate) {
1530  CHKERR PetscPrintf(
1531  PETSC_COMM_WORLD,
1532  "Energy equilibrium state is achieved! Difference = %0.6g %%. \n",
1533  equilibrium_rate * 100);
1534  }
1535  if (equilibrium_mass_rate < rate) {
1536  CHKERR PetscPrintf(
1537  PETSC_COMM_WORLD,
1538  "Mass equilibrium state is achieved! Difference = %0.6g %%. \n",
1539  equilibrium_mass_rate * 100);
1540  }
1541  commonData.cUrrent_psi = energ_sum;
1542  commonData.cUrrent_mass = mass_sum;
1543  }
1544  }
1545 
1546  int step;
1547 #if PETSC_VERSION_LE(3, 8, 0)
1548  CHKERR TSGetTimeStepNumber(ts, &step);
1549 #else
1550  CHKERR TSGetStepNumber(ts, &step);
1551 #endif
1552 
1553  if ((step) % pRT == 0) {
1554 
1555  ostringstream sss;
1556  sss << "out_" << step << ".h5m";
1557  CHKERR DMoFEMLoopFiniteElements(commonData.dm, "FE_REMODELLING", &postProc);
1558  CHKERR postProc.postProcMesh.write_file(sss.str().c_str(), "MOAB",
1559  "PARALLEL=WRITE_PART");
1560  if (analysis_mesh_flg) {
1561  ostringstream ttt;
1562  ttt << "analysis_mesh_" << step << ".h5m";
1563  CHKERR mField.get_moab().write_file(ttt.str().c_str(), "MOAB",
1564  "PARALLEL=WRITE_PART");
1565  }
1566 
1569  &postProcElastic);
1570  ostringstream ss;
1571  ss << "out_elastic_" << step << ".h5m";
1572  CHKERR postProcElastic.postProcMesh.write_file(ss.str().c_str(), "MOAB",
1573  "PARALLEL=WRITE_PART");
1574  }
1575  }
1577 }

◆ preProcess()

MoFEMErrorCode MonitorPostProc::preProcess ( )
Examples
Remodeling.hpp.

Definition at line 1394 of file Remodeling.cpp.

1394  {
1397 }

Member Data Documentation

◆ commonData

Remodeling::CommonData& BoneRemodeling::MonitorPostProc::commonData
Examples
Remodeling.hpp.

Definition at line 556 of file Remodeling.hpp.

◆ equilibrium_flg

PetscBool BoneRemodeling::MonitorPostProc::equilibrium_flg
Examples
Remodeling.hpp.

Definition at line 561 of file Remodeling.hpp.

◆ iNit

bool BoneRemodeling::MonitorPostProc::iNit
Examples
Remodeling.hpp.

Definition at line 563 of file Remodeling.hpp.

◆ mass_postproc

PetscBool BoneRemodeling::MonitorPostProc::mass_postproc
Examples
Remodeling.hpp.

Definition at line 560 of file Remodeling.hpp.

◆ mField

MoFEM::Interface& BoneRemodeling::MonitorPostProc::mField
Examples
Remodeling.hpp.

Definition at line 555 of file Remodeling.hpp.

◆ postProc

PostProcVolumeOnRefinedMesh BoneRemodeling::MonitorPostProc::postProc
Examples
Remodeling.hpp.

Definition at line 557 of file Remodeling.hpp.

◆ postProcElastic

PostProcVolumeOnRefinedMesh BoneRemodeling::MonitorPostProc::postProcElastic
Examples
Remodeling.hpp.

Definition at line 558 of file Remodeling.hpp.

◆ pRT

int BoneRemodeling::MonitorPostProc::pRT
Examples
Remodeling.hpp.

Definition at line 564 of file Remodeling.hpp.

◆ rate

double BoneRemodeling::MonitorPostProc::rate
Examples
Remodeling.hpp.

Definition at line 562 of file Remodeling.hpp.


The documentation for this struct was generated from the following files:
MoFEMFunctionReturnHot
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:447
BoneRemodeling::Remodeling::CommonData::nOremodellingBlock
bool nOremodellingBlock
Definition: Remodeling.hpp:164
BoneRemodeling::MonitorPostProc::mass_postproc
PetscBool mass_postproc
Definition: Remodeling.hpp:560
PostProcTemplateOnRefineMesh::postProcMesh
moab::Interface & postProcMesh
Definition: PostProcOnRefMesh.hpp:122
MoFEM::addHOOpsVol
MoFEMErrorCode addHOOpsVol(const std::string field, E &e, bool h1, bool hcurl, bool hdiv, bool l2)
Definition: HODataOperators.hpp:764
MoFEM::CoreInterface::get_comm_rank
virtual int get_comm_rank() const =0
BoneRemodeling::Remodeling::CommonData::cUrrent_psi
double cUrrent_psi
current free energy for evaluating equilibrium state
Definition: Remodeling.hpp:158
BoneRemodeling::Remodeling::CommonData::cUrrent_mass
double cUrrent_mass
current free energy for evaluating equilibrium state
Definition: Remodeling.hpp:160
BoneRemodeling::Remodeling::CommonData::DataContainers::rhoPtr
boost::shared_ptr< VectorDouble > rhoPtr
Ptr to density matrix container.
Definition: Remodeling.hpp:187
BoneRemodeling::Remodeling::CommonData::DataContainers::gradDispPtr
boost::shared_ptr< MatrixDouble > gradDispPtr
Ptr to gradient of displacements matrix container.
Definition: Remodeling.hpp:185
MoFEM::OpCalculateScalarFieldGradient
Get field gradients at integration pts for scalar filed rank 0, i.e. vector field.
Definition: UserDataOperators.hpp:1294
OpPostProcStress
Definition: ElasticityMixedFormulation.hpp:429
BoneRemodeling::Remodeling::CommonData::less_post_proc
PetscBool less_post_proc
reduce file size
Definition: Remodeling.hpp:163
CHKERR
#define CHKERR
Inline error check.
Definition: definitions.h:535
PostProcStress
Definition: PostProcStresses.hpp:17
MoFEM::CoreInterface::get_moab
virtual moab::Interface & get_moab()=0
BoneRemodeling::MonitorPostProc::postProcElastic
PostProcVolumeOnRefinedMesh postProcElastic
Definition: Remodeling.hpp:558
BoneRemodeling::Remodeling::CommonData::elasticPtr
boost::shared_ptr< NonlinearElasticElement > elasticPtr
Definition: Remodeling.hpp:140
BoneRemodeling::MonitorPostProc::equilibrium_flg
PetscBool equilibrium_flg
Definition: Remodeling.hpp:561
BoneRemodeling::Remodeling::CommonData::DataContainers::gradRhoPtr
boost::shared_ptr< MatrixDouble > gradRhoPtr
Gradient of density field.
Definition: Remodeling.hpp:188
PostProcTemplateVolumeOnRefinedMesh::commonData
CommonData commonData
Definition: PostProcOnRefMesh.hpp:287
MoFEM::OpCalculateScalarFieldValues
Get value at integration points for scalar field.
Definition: UserDataOperators.hpp:82
BoneRemodeling::MonitorPostProc::rate
double rate
Definition: Remodeling.hpp:562
BoneRemodeling::MonitorPostProc::pRT
int pRT
Definition: Remodeling.hpp:564
PostProcTemplateVolumeOnRefinedMesh::generateReferenceElementMesh
MoFEMErrorCode generateReferenceElementMesh()
Generate reference mesh on single element.
Definition: PostProcOnRefMesh.hpp:301
PostProcTemplateOnRefineMesh::addFieldValuesGradientPostProc
MoFEMErrorCode addFieldValuesGradientPostProc(const std::string field_name, Vec v=PETSC_NULL)
Add operator to post-process L2 or H1 field gradient.
Definition: PostProcOnRefMesh.hpp:195
MoFEM::OpCalculateVectorFieldGradient
Get field gradients at integration pts for scalar filed rank 0, i.e. vector field.
Definition: UserDataOperators.hpp:1536
BoneRemodeling::MonitorPostProc::postProc
PostProcVolumeOnRefinedMesh postProc
Definition: Remodeling.hpp:557
BoneRemodeling::MonitorPostProc::iNit
bool iNit
Definition: Remodeling.hpp:563
PostProcTemplateOnRefineMesh::mapGaussPts
std::vector< EntityHandle > mapGaussPts
Definition: PostProcOnRefMesh.hpp:125
BoneRemodeling::MonitorPostProc::mField
MoFEM::Interface & mField
Definition: Remodeling.hpp:555
BoneRemodeling::Remodeling::CommonData::dm
DM dm
Discretization manager.
Definition: Remodeling.hpp:127
MoFEM::Exceptions::ierr
static MoFEMErrorCodeGeneric< PetscErrorCode > ierr
Definition: Exceptions.hpp:76
MoFEM::ForcesAndSourcesCore::getOpPtrVector
boost::ptr_deque< UserDataOperator > & getOpPtrVector()
Use to push back operator for row operator.
Definition: ForcesAndSourcesCore.hpp:83
EigenMatrix::Vec
const FTensor::Tensor2< T, Dim, Dim > Vec
Definition: MatrixFunction.hpp:66
OpCalculateStress
Definition: HookeElement.hpp:153
MoFEMFunctionBeginHot
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:440
BoneRemodeling::Remodeling::CommonData::data
DataContainers data
Definition: Remodeling.hpp:203
MoFEM::DMoFEMLoopFiniteElements
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:590
MoFEMFunctionReturn
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:416
BoneRemodeling::MonitorPostProc::commonData
Remodeling::CommonData & commonData
Definition: Remodeling.hpp:556
MoFEMFunctionBegin
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:346
PostProcTemplateOnRefineMesh::addFieldValuesPostProc
MoFEMErrorCode addFieldValuesPostProc(const std::string field_name, Vec v=PETSC_NULL)
Add operator to post-process L2, H1, Hdiv, Hcurl field value.
Definition: PostProcOnRefMesh.hpp:153