v0.14.0
Loading...
Searching...
No Matches
FluidPressure.hpp
Go to the documentation of this file.
1/* \file FluidPressure.hpp
2 *
3 * \brief Implementation of fluid pressure element
4 *
5 * \todo Implement nonlinear case (consrvative force, i.e. normal follows
6 * surface normal)
7 *
8 */
9
10
11
12#ifndef __FLUID_PRESSURE_HPP
13#define __FLUID_PRESSURE_HPP
14
15/** \brief Fluid pressure forces
16
17\todo Implementation for large displacements
18
19*/
21
24
27 int getRule(int order) { return order + 1; };
28
29 MoFEMErrorCode preProcess() {
32 }
33 };
35 MyTriangleFE &getLoopFe() { return fe; }
36
37 FluidPressure(MoFEM::Interface &m_field) : mField(m_field), fe(mField) {}
38
39 typedef int MeshSetId;
40 struct FluidData {
41 double dEnsity; ///< fluid density [kg/m^2] or any consistent unit
42 VectorDouble aCCeleration; ///< acceleration [m/s^2]
43 VectorDouble zEroPressure; ///< fluid level of reference zero pressure.
44 Range
45 tRis; ///< range of surface element to which fluid pressure is applied
46 friend std::ostream &operator<<(std::ostream &os,
48 };
49 std::map<MeshSetId, FluidData> setOfFluids;
50
51 boost::ptr_vector<MethodForForceScaling> methodsOp;
52
55
56 Vec F;
58 boost::ptr_vector<MethodForForceScaling> &methodsOp;
59 bool allowNegativePressure; ///< allows for negative pressures
61
62 OpCalculatePressure(const std::string field_name, Vec _F, FluidData &data,
63 boost::ptr_vector<MethodForForceScaling> &methods_op,
64 bool allow_negative_pressure, bool ho_geometry)
67 F(_F), dAta(data), methodsOp(methods_op),
68 allowNegativePressure(allow_negative_pressure),
69 hoGeometry(ho_geometry) {}
70
71 VectorDouble Nf;
72
73 MoFEMErrorCode doWork(int side, EntityType type,
74 EntitiesFieldData::EntData &data);
75
76 };
77
79 const std::string field_name,
80 const std::string mesh_nodals_positions = "MESH_NODE_POSITIONS");
81
83 string field_name, Vec F, bool allow_negative_pressure = true,
84 bool ho_geometry = false);
85};
86
87std::ostream &operator<<(std::ostream &os, const FluidPressure::FluidData &e);
88
89#endif //__FLUID_PRESSSURE_HPP
ForcesAndSourcesCore::UserDataOperator UserDataOperator
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:447
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:440
@ F
implementation of Data Operators for Forces and Sources
Definition: Common.hpp:10
constexpr auto field_name
#define _F(n)
Definition: quad.c:25
friend std::ostream & operator<<(std::ostream &os, const FluidPressure::FluidData &e)
VectorDouble zEroPressure
fluid level of reference zero pressure.
Range tRis
range of surface element to which fluid pressure is applied
double dEnsity
fluid density [kg/m^2] or any consistent unit
VectorDouble aCCeleration
acceleration [m/s^2]
MoFEMErrorCode preProcess()
function is run at the beginning of loop
MyTriangleFE(MoFEM::Interface &m_field)
boost::ptr_vector< MethodForForceScaling > & methodsOp
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
bool allowNegativePressure
allows for negative pressures
OpCalculatePressure(const std::string field_name, Vec _F, FluidData &data, boost::ptr_vector< MethodForForceScaling > &methods_op, bool allow_negative_pressure, bool ho_geometry)
Fluid pressure forces.
boost::ptr_vector< MethodForForceScaling > methodsOp
std::map< MeshSetId, FluidData > setOfFluids
MyTriangleFE & getLoopFe()
MoFEMErrorCode setNeumannFluidPressureFiniteElementOperators(string field_name, Vec F, bool allow_negative_pressure=true, bool ho_geometry=false)
MyTriangleFE fe
FluidPressure(MoFEM::Interface &m_field)
MoFEM::Interface & mField
MoFEMErrorCode addNeumannFluidPressureBCElements(const std::string field_name, const std::string mesh_nodals_positions="MESH_NODE_POSITIONS")
Deprecated interface functions.
@ OPROW
operator doWork function is executed on FE rows