v0.15.0
Loading...
Searching...
No Matches
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.cpp, and Remodeling.hpp.

Definition at line 1386 of file Remodeling.cpp.

1388 : mField(m_field), commonData(common_data), postProc(m_field),
1389 postProcElastic(m_field),
1390 // densityMaps(m_field),
1391 iNit(false) {}
Remodeling::CommonData & commonData
PostProcVolumeOnRefinedMesh postProcElastic
PostProcVolumeOnRefinedMesh postProc

Member Function Documentation

◆ operator()()

MoFEMErrorCode MonitorPostProc::operator() ( )
Examples
Remodeling.cpp, and Remodeling.hpp.

Definition at line 1398 of file Remodeling.cpp.

1398 {
1401}
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...

◆ postProcess()

MoFEMErrorCode MonitorPostProc::postProcess ( )
Examples
Remodeling.cpp, and Remodeling.hpp.

Definition at line 1403 of file Remodeling.cpp.

1403 {
1405
1406 PetscBool mass_postproc = PETSC_FALSE;
1407 PetscBool equilibrium_flg = PETSC_FALSE;
1408 PetscBool analysis_mesh_flg = PETSC_FALSE;
1409 rate = 1;
1410
1411 PetscOptionsBegin(PETSC_COMM_WORLD, "", "Bone remodeling post-process",
1412 "none");
1413
1414 CHKERR PetscOptionsBool(
1415 "-mass_postproc",
1416 "is used for testing, calculates mass and energy at each time step", "",
1417 mass_postproc, &mass_postproc, PETSC_NULL);
1418
1419 CHKERR PetscOptionsBool("-analysis_mesh",
1420 "saves analysis mesh at each time step", "",
1421 analysis_mesh_flg, &analysis_mesh_flg, PETSC_NULL);
1422
1423 CHKERR PetscOptionsScalar(
1424 "-equilibrium_stop_rate",
1425 "is used to stop calculations when equilibium state is achieved", "",
1427 PetscOptionsEnd();
1428
1429 if (!iNit) {
1430
1431 PetscOptionsBegin(PETSC_COMM_WORLD, "",
1432 "Bone remodeling post-process", "none");
1433
1434 pRT = 1;
1435 CHKERR PetscOptionsInt("-my_output_prt",
1436 "frequncy how often results are dumped on hard disk",
1437 "", 1, &pRT, NULL);
1438 PetscOptionsEnd();
1439
1441 CHKERR postProc.addFieldValuesPostProc("DISPLACEMENTS");
1442 CHKERR postProc.addFieldValuesPostProc("MESH_NODE_POSITIONS");
1445 }
1446 // CHKERR postProc.addFieldValuesPostProc("RHO");
1447 // CHKERR postProc.addFieldValuesGradientPostProc("RHO");
1448
1449 {
1450 auto &list_of_operators = postProc.getOpPtrVector();
1451
1452 list_of_operators.push_back(
1454 list_of_operators.push_back(new OpCalculateScalarFieldGradient<3>(
1455 "RHO", commonData.data.gradRhoPtr));
1456 // Get displacement gradient at integration points
1457 list_of_operators.push_back(new OpCalculateVectorFieldGradient<3, 3>(
1458 "DISPLACEMENTS", commonData.data.gradDispPtr));
1459 list_of_operators.push_back(new OpCalculateStress(commonData));
1460 list_of_operators.push_back(new OpPostProcStress(
1462 }
1463
1466 CHKERR postProcElastic.addFieldValuesPostProc("MESH_NODE_POSITIONS");
1468 std::map<int, NonlinearElasticElement::BlockData>::iterator sit =
1469 commonData.elasticPtr->setOfBlocks.begin();
1470 for (; sit != commonData.elasticPtr->setOfBlocks.end(); sit++) {
1473 "DISPLACEMENTS", sit->second, postProcElastic.commonData));
1474 }
1475
1476 iNit = true;
1477 }
1478
1479 if (mass_postproc) {
1480 Vec mass_vec;
1481 Vec energ_vec;
1482 int mass_vec_ghost[] = {0};
1483 CHKERR VecCreateGhost(PETSC_COMM_WORLD, (!mField.get_comm_rank()) ? 1 : 0,
1484 1, 1, mass_vec_ghost, &mass_vec);
1485
1486 CHKERR VecZeroEntries(mass_vec);
1487 CHKERR VecDuplicate(mass_vec, &energ_vec);
1488
1489 Remodeling::Fe energyCalc(mField);
1490 CHKERR addHOOpsVol("MESH_NODE_POSITIONS", energyCalc, true, false, false,
1491 true);
1492 energyCalc.getOpPtrVector().push_back(
1494 energyCalc.getOpPtrVector().push_back(new OpCalculateScalarFieldGradient<3>(
1495 "RHO", commonData.data.gradRhoPtr));
1496 energyCalc.getOpPtrVector().push_back(
1497 new OpCalculateVectorFieldGradient<3, 3>("DISPLACEMENTS",
1499 energyCalc.getOpPtrVector().push_back(new OpCalculateStress(commonData));
1500 energyCalc.getOpPtrVector().push_back(
1501 new OpMassAndEnergyCalculation("RHO", commonData, energ_vec, mass_vec));
1502 CHKERR DMoFEMLoopFiniteElements(commonData.dm, "FE_REMODELLING",
1503 &energyCalc);
1504
1505 CHKERR VecAssemblyBegin(mass_vec);
1506 CHKERR VecAssemblyEnd(mass_vec);
1507 double mass_sum;
1508 CHKERR VecSum(mass_vec, &mass_sum);
1509
1510 CHKERR VecAssemblyBegin(energ_vec);
1511 CHKERR VecAssemblyEnd(energ_vec);
1512 double energ_sum;
1513 CHKERR VecSum(energ_vec, &energ_sum);
1514
1515 CHKERR PetscPrintf(PETSC_COMM_WORLD,
1516 "Time: %g Mass: %6.5g Elastic energy: %6.5g \n", ts_t,
1517 mass_sum, energ_sum);
1518 CHKERR VecDestroy(&mass_vec);
1519 CHKERR VecDestroy(&energ_vec);
1520
1521 if (equilibrium_flg) {
1522 double equilibrium_rate =
1523 fabs(energ_sum - commonData.cUrrent_psi) / energ_sum;
1524 double equilibrium_mass_rate =
1525 fabs(mass_sum - commonData.cUrrent_mass) / mass_sum;
1526 if (equilibrium_rate < rate) {
1527 CHKERR PetscPrintf(
1528 PETSC_COMM_WORLD,
1529 "Energy equilibrium state is achieved! Difference = %0.6g %%. \n",
1530 equilibrium_rate * 100);
1531 }
1532 if (equilibrium_mass_rate < rate) {
1533 CHKERR PetscPrintf(
1534 PETSC_COMM_WORLD,
1535 "Mass equilibrium state is achieved! Difference = %0.6g %%. \n",
1536 equilibrium_mass_rate * 100);
1537 }
1538 commonData.cUrrent_psi = energ_sum;
1539 commonData.cUrrent_mass = mass_sum;
1540 }
1541 }
1542
1543 int step;
1544#if PETSC_VERSION_LE(3, 8, 0)
1545 CHKERR TSGetTimeStepNumber(ts, &step);
1546#else
1547 CHKERR TSGetStepNumber(ts, &step);
1548#endif
1549
1550 if ((step) % pRT == 0) {
1551
1552 ostringstream sss;
1553 sss << "out_" << step << ".h5m";
1555 CHKERR postProc.postProcMesh.write_file(sss.str().c_str(), "MOAB",
1556 "PARALLEL=WRITE_PART");
1557 if (analysis_mesh_flg) {
1558 ostringstream ttt;
1559 ttt << "analysis_mesh_" << step << ".h5m";
1560 CHKERR mField.get_moab().write_file(ttt.str().c_str(), "MOAB",
1561 "PARALLEL=WRITE_PART");
1562 }
1563
1567 ostringstream ss;
1568 ss << "out_elastic_" << step << ".h5m";
1569 CHKERR postProcElastic.postProcMesh.write_file(ss.str().c_str(), "MOAB",
1570 "PARALLEL=WRITE_PART");
1571 }
1572 }
1574}
#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.
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:576
MoFEMErrorCode addFieldValuesPostProc(const std::string field_name, Vec v=PETSC_NULLPTR)
Add operator to post-process L2, H1, Hdiv, Hcurl field value.
MoFEMErrorCode addFieldValuesGradientPostProc(const std::string field_name, Vec v=PETSC_NULLPTR)
Add operator to post-process L2 or H1 field gradient.
MoFEMErrorCode addHOOpsVol(const std::string field, E &e, bool h1, bool hcurl, bool hdiv, bool l2)
boost::shared_ptr< MatrixDouble > gradDispPtr
Ptr to gradient of displacements matrix container.
boost::shared_ptr< VectorDouble > rhoPtr
Ptr to density matrix container.
boost::shared_ptr< MatrixDouble > gradRhoPtr
Gradient of density field.
PetscBool less_post_proc
reduce file size
double cUrrent_psi
current free energy for evaluating equilibrium state
boost::shared_ptr< NonlinearElasticElement > elasticPtr
double cUrrent_mass
current free energy for evaluating equilibrium state
virtual moab::Interface & get_moab()=0
virtual int get_comm_rank() const =0
boost::ptr_deque< UserDataOperator > & getOpPtrVector()
Use to push back operator for row operator.
Get field gradients at integration pts for scalar field rank 0, i.e. vector field.
Get value at integration points for scalar field.
Get field gradients at integration pts for scalar field rank 0, i.e. vector field.
std::vector< EntityHandle > mapGaussPts
MoFEMErrorCode generateReferenceElementMesh()
Generate reference mesh on single element.

◆ preProcess()

MoFEMErrorCode MonitorPostProc::preProcess ( )
Examples
Remodeling.cpp, and Remodeling.hpp.

Definition at line 1393 of file Remodeling.cpp.

1393 {
1396}

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: