v0.15.0
Loading...
Searching...
No Matches
ThermalElement::TimeSeriesMonitor Struct Reference

TS monitore it records temperature at time steps. More...

#include "users_modules/basic_finite_elements/src/ThermalElement.hpp"

Inheritance diagram for ThermalElement::TimeSeriesMonitor:
[legend]
Collaboration diagram for ThermalElement::TimeSeriesMonitor:
[legend]

Public Types

using Ele = ForcesAndSourcesCore
 
using VolEle = VolumeElementForcesAndSourcesCore
 
using VolOp = VolumeElementForcesAndSourcesCore::UserDataOperator
 
using SetPtsData = FieldEvaluatorInterface::SetPtsData
 

Public Member Functions

 TimeSeriesMonitor (MoFEM::Interface &m_field, const std::string series_name, const std::string temp_name, std::vector< std::array< double, 3 > > eval_points={})
 
MoFEMErrorCode postProcess ()
 

Public Attributes

MoFEM::InterfacemField
 
const std::string seriesName
 
const std::string tempName
 
BitRefLevel mask
 
boost::shared_ptr< SetPtsDatadataFieldEval
 
boost::shared_ptr< VectorDouble > tempPtr
 
std::vector< std::array< double, 3 > > evalPoints
 

Detailed Description

TS monitore it records temperature at time steps.

Definition at line 580 of file ThermalElement.hpp.

Member Typedef Documentation

◆ Ele

using ThermalElement::TimeSeriesMonitor::Ele = ForcesAndSourcesCore

Definition at line 587 of file ThermalElement.hpp.

◆ SetPtsData

using ThermalElement::TimeSeriesMonitor::SetPtsData = FieldEvaluatorInterface::SetPtsData

Definition at line 590 of file ThermalElement.hpp.

◆ VolEle

using ThermalElement::TimeSeriesMonitor::VolEle = VolumeElementForcesAndSourcesCore

Definition at line 588 of file ThermalElement.hpp.

◆ VolOp

using ThermalElement::TimeSeriesMonitor::VolOp = VolumeElementForcesAndSourcesCore::UserDataOperator

Definition at line 589 of file ThermalElement.hpp.

Constructor & Destructor Documentation

◆ TimeSeriesMonitor()

ThermalElement::TimeSeriesMonitor::TimeSeriesMonitor ( MoFEM::Interface & m_field,
const std::string series_name,
const std::string temp_name,
std::vector< std::array< double, 3 > > eval_points = {} )
inline

Definition at line 597 of file ThermalElement.hpp.

599 {})
600 : mField(m_field), seriesName(series_name), tempName(temp_name),
601 evalPoints(eval_points),
602 dataFieldEval(m_field.getInterface<FieldEvaluatorInterface>()
603 ->getData<VolEle>()) {
604 mask.set();
605
606 if (!evalPoints.empty()) {
607 ierr = m_field.getInterface<FieldEvaluatorInterface>()->buildTree<3>(
608 dataFieldEval, "THERMAL_FE");
609 CHKERRABORT(PETSC_COMM_WORLD, ierr);
610
611 auto no_rule = [](int, int, int) { return -1; };
612
613 tempPtr = boost::make_shared<VectorDouble>();
614
615 boost::shared_ptr<Ele> vol_ele(dataFieldEval->feMethodPtr.lock());
616 vol_ele->getRuleHook = no_rule;
617
618 vol_ele->getOpPtrVector().push_back(
619 new OpCalculateScalarFieldValues("TEMP", tempPtr));
620 }
621 }
boost::shared_ptr< SetPtsData > dataFieldEval
std::vector< std::array< double, 3 > > evalPoints
boost::shared_ptr< VectorDouble > tempPtr

Member Function Documentation

◆ postProcess()

MoFEMErrorCode ThermalElement::TimeSeriesMonitor::postProcess ( )

Definition at line 454 of file ThermalElement.cpp.

454 {
456
457 CHKERR mField.getInterface<VecManager>()->setGlobalGhostVector(
458 problemPtr, ROW, ts_u, INSERT_VALUES, SCATTER_REVERSE);
459
460 BitRefLevel proble_bit_level = problemPtr->getBitRefLevel();
461
462 SeriesRecorder *recorder_ptr = NULL;
463 CHKERR mField.getInterface(recorder_ptr);
464 CHKERR recorder_ptr->record_begin(seriesName);
465 CHKERR recorder_ptr->record_field(seriesName, tempName, proble_bit_level,
466 mask);
467 CHKERR recorder_ptr->record_end(seriesName, ts_t);
468
469 auto post_proc_at_points = [&](std::array<double, 3> point, int num) {
471
472 dataFieldEval->setEvalPoints(point.data(), point.size() / 3);
473
474 struct OpPrint : public VolOp {
475
476 std::array<double, 3> pointCoords;
477 int pointNum;
478 boost::shared_ptr<VectorDouble> tempPtr;
479
480 OpPrint(boost::shared_ptr<VectorDouble> temp_ptr,
481 std::array<double, 3> &point_coords, int point_num)
482 : VolOp("TEMP", VolOp::OPROW), tempPtr(temp_ptr),
483 pointCoords(point_coords), pointNum(point_num) {}
484
485 MoFEMErrorCode doWork(int side, EntityType type,
488 if (type == MBVERTEX) {
489 if (getGaussPts().size2()) {
490
491 auto t_p = getFTensor0FromVec(*tempPtr);
492
493 MOFEM_LOG("THERMALSYNC", Sev::inform)
494 << "Pnt: " << std::to_string(pointNum)
495 << " Crd: " << getVectorAdaptor(pointCoords.data(), 3)
496 << " Tmp: " << t_p;
497 }
498 }
500 }
501 };
502
503 if (auto fe_ptr = dataFieldEval->feMethodPtr.lock()) {
504 fe_ptr->getOpPtrVector().push_back(new OpPrint(tempPtr, point, num));
506 ->evalFEAtThePoint<3>(point.data(), 1e-12, "DMTHERMAL", "THERMAL_FE",
508 mField.get_comm_rank(), nullptr, MF_EXIST,
509 QUIET);
510 fe_ptr->getOpPtrVector().pop_back();
511 }
512
514 };
515
516 if (!evalPoints.empty()) {
517 int num = 0;
518 for (auto p : evalPoints)
519 CHKERR post_proc_at_points(p, num++);
520 }
521
523
525}
#define MOFEM_LOG_SYNCHRONISE(comm)
Synchronise "SYNC" channel.
@ QUIET
@ ROW
@ MF_EXIST
#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 MOFEM_LOG(channel, severity)
Log.
virtual MoFEMErrorCode record_field(const std::string &serie_name, const std::string &field_name, const BitRefLevel &bit, const BitRefLevel &mask)
virtual MoFEMErrorCode record_begin(const std::string &serie_name)
virtual MoFEMErrorCode record_end(const std::string &serie_name, double time=0)
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
std::bitset< BITREFLEVEL_SIZE > BitRefLevel
Bit structure attached to each entity identifying to what mesh entity is attached.
Definition Types.hpp:40
auto getVectorAdaptor(T1 ptr, const size_t n)
Get Vector adaptor.
Definition Templates.hpp:31
static auto getFTensor0FromVec(ublas::vector< T, A > &data)
Get tensor rank 0 (scalar) form data vector.
virtual MPI_Comm & get_comm() const =0
virtual int get_comm_rank() const =0
virtual MoFEMErrorCode doWork(int row_side, int col_side, EntityType row_type, EntityType col_type, EntitiesFieldData::EntData &row_data, EntitiesFieldData::EntData &col_data)
Operator for bi-linear form, usually to calculate values on left hand side.
Data on single entity (This is passed as argument to DataOperator::doWork)
Field evaluator interface.
@ OPROW
operator doWork function is executed on FE rows
MatrixDouble & getGaussPts()
matrix of integration (Gauss) points for Volume Element
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface reference to pointer of interface.
Vector manager is used to create vectors \mofem_vectors.

Member Data Documentation

◆ dataFieldEval

boost::shared_ptr<SetPtsData> ThermalElement::TimeSeriesMonitor::dataFieldEval

Definition at line 592 of file ThermalElement.hpp.

◆ evalPoints

std::vector<std::array<double, 3> > ThermalElement::TimeSeriesMonitor::evalPoints

Definition at line 595 of file ThermalElement.hpp.

◆ mask

BitRefLevel ThermalElement::TimeSeriesMonitor::mask

Definition at line 585 of file ThermalElement.hpp.

◆ mField

MoFEM::Interface& ThermalElement::TimeSeriesMonitor::mField

Definition at line 582 of file ThermalElement.hpp.

◆ seriesName

const std::string ThermalElement::TimeSeriesMonitor::seriesName

Definition at line 583 of file ThermalElement.hpp.

◆ tempName

const std::string ThermalElement::TimeSeriesMonitor::tempName

Definition at line 584 of file ThermalElement.hpp.

◆ tempPtr

boost::shared_ptr<VectorDouble> ThermalElement::TimeSeriesMonitor::tempPtr

Definition at line 593 of file ThermalElement.hpp.


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