v0.15.0
Loading...
Searching...
No Matches
Public Member Functions | Static Public Member Functions | List of all members
MoFEM::AddFluxToRhsPipelineImpl< MoFEM::OpFluxRhsImpl< ThermoElasticOps::SetTargetTemperature, 1, 1, A, I, OpBase >, A, I, OpBase > Struct Template Reference

#include "users_modules/multifield-thermoplasticity-private/tutorials/adv-2/src/ThermoElasticOps.hpp"

Public Member Functions

 AddFluxToRhsPipelineImpl ()=delete
 
 AddFluxToRhsPipelineImpl ()=delete
 

Static Public Member Functions

static MoFEMErrorCode add (boost::ptr_deque< ForcesAndSourcesCore::UserDataOperator > &pipeline, MoFEM::Interface &m_field, const std::string field_name, boost::shared_ptr< VectorDouble > temp_ptr, std::string block_name, Sev sev)
 
static MoFEMErrorCode add (boost::ptr_deque< ForcesAndSourcesCore::UserDataOperator > &pipeline, MoFEM::Interface &m_field, const std::string field_name, boost::shared_ptr< VectorDouble > temp_ptr, std::string block_name, Sev sev)
 

Detailed Description

template<AssemblyType A, IntegrationType I, typename OpBase>
struct MoFEM::AddFluxToRhsPipelineImpl< MoFEM::OpFluxRhsImpl< ThermoElasticOps::SetTargetTemperature, 1, 1, A, I, OpBase >, A, I, OpBase >

Definition at line 191 of file ThermoElasticOps.hpp.

Constructor & Destructor Documentation

◆ AddFluxToRhsPipelineImpl() [1/2]

template<AssemblyType A, IntegrationType I, typename OpBase >
MoFEM::AddFluxToRhsPipelineImpl< MoFEM::OpFluxRhsImpl< ThermoElasticOps::SetTargetTemperature, 1, 1, A, I, OpBase >, A, I, OpBase >::AddFluxToRhsPipelineImpl ( )
delete

◆ AddFluxToRhsPipelineImpl() [2/2]

template<AssemblyType A, IntegrationType I, typename OpBase >
MoFEM::AddFluxToRhsPipelineImpl< MoFEM::OpFluxRhsImpl< ThermoElasticOps::SetTargetTemperature, 1, 1, A, I, OpBase >, A, I, OpBase >::AddFluxToRhsPipelineImpl ( )
delete

Member Function Documentation

◆ add() [1/2]

template<AssemblyType A, IntegrationType I, typename OpBase >
static MoFEMErrorCode MoFEM::AddFluxToRhsPipelineImpl< MoFEM::OpFluxRhsImpl< ThermoElasticOps::SetTargetTemperature, 1, 1, A, I, OpBase >, A, I, OpBase >::add ( boost::ptr_deque< ForcesAndSourcesCore::UserDataOperator > &  pipeline,
MoFEM::Interface m_field,
const std::string  field_name,
boost::shared_ptr< VectorDouble temp_ptr,
std::string  block_name,
Sev  sev 
)
inlinestatic

Definition at line 201 of file ThermoElasticOps.hpp.

201 {
202 return boost::shared_ptr<double>(shared_from_this(), &refTemp);
203 }
204};
205
207 MoFEM::Interface &m_field,
208 boost::ptr_deque<ForcesAndSourcesCore::UserDataOperator> &pipeline,
209 std::string block_name,
210 boost::shared_ptr<BlockedThermoElasticParameters> blockedParamsPtr,
211 double default_coeff_expansion, double default_ref_temp, Sev sev) {
213
214 double local_coeff_expansion = default_coeff_expansion;
215 double local_ref_temp = default_ref_temp;
216
217 struct OpMatThermoElasticBlocks : public DomainEleOp {
218 OpMatThermoElasticBlocks(boost::shared_ptr<VectorDouble> expansion_ptr,
219 boost::shared_ptr<double> ref_temp_ptr,
220 double local_coeff_expansion,
221 double local_ref_temp, MoFEM::Interface &m_field,
222 Sev sev,
223 std::vector<const CubitMeshSets *> meshset_vec_ptr)
224 : DomainEleOp(NOSPACE, DomainEleOp::OPSPACE),
225 expansionPtr(expansion_ptr), refTempPtr(ref_temp_ptr),
226 defaultCoeffExpansion(local_coeff_expansion),
227 defaultRefTemp(local_ref_temp) {
229 extractThermoElasticBlockData(m_field, meshset_vec_ptr, sev),
230 "Cannot get data from thermoelastic block");
231 }
232
233 MoFEMErrorCode doWork(int side, EntityType type,
234 EntitiesFieldData::EntData &data) {
236
237 for (auto &b : blockData) {
238
239 if (b.blockEnts.find(getFEEntityHandle()) != b.blockEnts.end()) {
240 *expansionPtr = b.expansion;
241 *refTempPtr = b.ref_temp;
243 }
244 }
245
246 *expansionPtr = VectorDouble(SPACE_DIM, defaultCoeffExpansion);
247 *refTempPtr = defaultRefTemp;
248
250 }
251
252 private:
253 double defaultCoeffExpansion;
254 double defaultRefTemp;
255 struct BlockData {
256 double ref_temp;
257 VectorDouble expansion;
258 Range blockEnts;
259 };
260
261 std::vector<BlockData> blockData;
constexpr int SPACE_DIM
[Define dimension]
#define CHK_THROW_MESSAGE(err, msg)
Check and throw MoFEM exception.
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
@ NOSPACE
Definition definitions.h:83
#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()
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
UBlasVector< double > VectorDouble
Definition Types.hpp:68
MoFEM::LogManager::SeverityLevel Sev
MoFEMErrorCode addMatThermoElasticBlockOps(MoFEM::Interface &m_field, boost::ptr_deque< ForcesAndSourcesCore::UserDataOperator > &pipeline, std::string block_name, boost::shared_ptr< BlockedThermoElasticParameters > blockedParamsPtr, double default_coeff_expansion, double default_ref_temp, Sev sev)
Deprecated interface functions.
double default_ref_temp
double default_coeff_expansion

◆ add() [2/2]

template<AssemblyType A, IntegrationType I, typename OpBase >
static MoFEMErrorCode MoFEM::AddFluxToRhsPipelineImpl< MoFEM::OpFluxRhsImpl< ThermoElasticOps::SetTargetTemperature, 1, 1, A, I, OpBase >, A, I, OpBase >::add ( boost::ptr_deque< ForcesAndSourcesCore::UserDataOperator > &  pipeline,
MoFEM::Interface m_field,
const std::string  field_name,
boost::shared_ptr< VectorDouble temp_ptr,
std::string  block_name,
Sev  sev 
)
inlinestatic

Definition at line 519 of file ThermoElasticOps.hpp.

525 {
527
528 using OP_SOURCE = typename FormsIntegrators<OpBase>::template Assembly<
529 A>::template LinearForm<I>::template OpSource<1, 1>;
530 using OP_TEMP = typename FormsIntegrators<OpBase>::template Assembly<
531 A>::template LinearForm<I>::template OpBaseTimesScalar<1>;
532
533 auto add_op = [&](auto &&meshset_vec_ptr) {
535 for (auto m : meshset_vec_ptr) {
536 std::vector<double> block_data;
537 m->getAttributes(block_data);
538 if (block_data.size() < 2) {
539 SETERRQ(PETSC_COMM_WORLD, MOFEM_DATA_INCONSISTENCY,
540 "Expected two parameters");
541 }
542 double target_temperature = block_data[0];
543 double beta =
544 block_data[1]; // Set temperature parameter [ W/K * (1/m^3)]
545 auto ents_ptr = boost::make_shared<Range>();
546 CHKERR m_field.get_moab().get_entities_by_handle(m->meshset,
547 *(ents_ptr), true);
548
549 MOFEM_TAG_AND_LOG("WORLD", sev, "SetTargetTemperature")
550 << "Add " << *m << " target temperature " << target_temperature
551 << " penalty " << beta;
552
553 pipeline.push_back(new OP_SOURCE(
555 [target_temperature, beta](double, double, double) {
556 return target_temperature * beta;
557 },
558 ents_ptr));
559 pipeline.push_back(new OP_TEMP(
560 field_name, temp_ptr,
561 [beta](double, double, double) { return -beta; }, ents_ptr));
562 }
564 };
565
566 CHKERR add_op(
567
568 m_field.getInterface<MeshsetsManager>()->getCubitMeshsetPtr(std::regex(
569
570 (boost::format("%s(.*)") % block_name).str()
571
572 ))
573
574 );
575
576 MOFEM_LOG_CHANNEL("WORLD");
577
579 }
#define MOFEM_TAG_AND_LOG(channel, severity, tag)
Tag and log in channel.
@ MOFEM_DATA_INCONSISTENCY
Definition definitions.h:31
#define CHKERR
Inline error check.
#define MOFEM_LOG_CHANNEL(channel)
Set and reset channel.
constexpr AssemblyType A
[Define dimension]
constexpr auto field_name
FTensor::Index< 'm', 3 > m
virtual moab::Interface & get_moab()=0
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface reference to pointer of interface.

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