v0.14.0
Static Public Member Functions | List of all members
MetaNeumannForces Struct Reference

Set of high-level function declaring elements and setting operators to apply forces/fluxes. More...

#include <users_modules/basic_finite_elements/src/SurfacePressure.hpp>

Static Public Member Functions

static MoFEMErrorCode addNeumannBCElements (MoFEM::Interface &m_field, const std::string field_name, const std::string mesh_nodals_positions="MESH_NODE_POSITIONS", Range *intersect_ptr=NULL)
 Declare finite element. More...
 
static MoFEMErrorCode setMomentumFluxOperators (MoFEM::Interface &m_field, boost::ptr_map< std::string, NeumannForcesSurface > &neumann_forces, Vec F, const std::string field_name, const std::string mesh_nodals_positions="MESH_NODE_POSITIONS")
 Set operators to finite elements calculating right hand side vector. More...
 
static MoFEMErrorCode addNeumannFluxBCElements (MoFEM::Interface &m_field, const std::string field_name, const std::string mesh_nodals_positions="MESH_NODE_POSITIONS")
 
static MoFEMErrorCode setMassFluxOperators (MoFEM::Interface &m_field, boost::ptr_map< std::string, NeumannForcesSurface > &neumann_forces, Vec F, const std::string field_name, const std::string mesh_nodals_positions="MESH_NODE_POSITIONS")
 

Detailed Description

Set of high-level function declaring elements and setting operators to apply forces/fluxes.

Definition at line 1261 of file SurfacePressure.hpp.

Member Function Documentation

◆ addNeumannBCElements()

MoFEMErrorCode MetaNeumannForces::addNeumannBCElements ( MoFEM::Interface m_field,
const std::string  field_name,
const std::string  mesh_nodals_positions = "MESH_NODE_POSITIONS",
Range intersect_ptr = NULL 
)
static

Declare finite element.

Search cubit sidesets and blocksets with pressure bc and declare surface elemen

Block set has to have name “PRESSURE”. Can have name “PRESSURE_01” or any other name with prefix. The first attribute of block set is pressure value.

Parameters
m_fieldInterface insurance
field_nameField name (f.e. DISPLACEMENT)
mesh_nodals_positionsName of field on which ho-geometry is defined
intersect_ptrPointer to range to interect meshset entities
Returns
Error code
Examples
elasticity.cpp, elasticity_mixed_formulation.cpp, mortar_contact_thermal.cpp, navier_stokes.cpp, nonlinear_dynamics.cpp, simple_contact.cpp, and simple_contact_thermal.cpp.

Definition at line 1974 of file SurfacePressure.cpp.

1976  {
1978 
1979  // Define boundary element that operates on rows, columns and data of a
1980  // given field
1981  CHKERR m_field.add_finite_element("FORCE_FE", MF_ZERO);
1985  if (m_field.check_field(mesh_nodals_positions)) {
1986  CHKERR m_field.modify_finite_element_add_field_data("FORCE_FE",
1987  mesh_nodals_positions);
1988  }
1989  // Add entities to that element, here we add all triangles with FORCESET
1990  // from cubit
1992  it)) {
1993  Range tris;
1994  CHKERR m_field.get_moab().get_entities_by_dimension(it->meshset, 2, tris,
1995  true);
1996  if (intersect_ptr)
1997  tris = intersect(tris, *intersect_ptr);
1998  CHKERR m_field.add_ents_to_finite_element_by_dim(tris, 2, "FORCE_FE");
1999  }
2000 
2001  CHKERR m_field.add_finite_element("PRESSURE_FE", MF_ZERO);
2002  CHKERR m_field.modify_finite_element_add_field_row("PRESSURE_FE", field_name);
2003  CHKERR m_field.modify_finite_element_add_field_col("PRESSURE_FE", field_name);
2004  CHKERR m_field.modify_finite_element_add_field_data("PRESSURE_FE",
2005  field_name);
2006  if (m_field.check_field(mesh_nodals_positions)) {
2007  CHKERR m_field.modify_finite_element_add_field_data("PRESSURE_FE",
2008  mesh_nodals_positions);
2009  }
2010 
2012  SIDESET | PRESSURESET, it)) {
2013  Range tris;
2014  CHKERR m_field.get_moab().get_entities_by_dimension(it->meshset, 2, tris,
2015  true);
2016  if (intersect_ptr)
2017  tris = intersect(tris, *intersect_ptr);
2018  CHKERR m_field.add_ents_to_finite_element_by_dim(tris, 2, "PRESSURE_FE");
2019  }
2020 
2021  // Reading forces from BLOCKSET
2022 
2023  const string block_set_force_name("FORCE");
2024  // search for block named FORCE and add its attributes to FORCE_FE element
2026  if (it->getName().compare(0, block_set_force_name.length(),
2027  block_set_force_name) == 0) {
2028  Range tris;
2029  CHKERR m_field.get_moab().get_entities_by_dimension(it->meshset, 2, tris,
2030  true);
2031  if (intersect_ptr)
2032  tris = intersect(tris, *intersect_ptr);
2033  CHKERR m_field.add_ents_to_finite_element_by_dim(tris, 2, "FORCE_FE");
2034  }
2035  }
2036  // search for block named PRESSURE and add its attributes to PRESSURE_FE
2037  // element
2038  const string block_set_pressure_name("PRESSURE");
2040  if (it->getName().compare(0, block_set_pressure_name.length(),
2041  block_set_pressure_name) == 0) {
2042  Range tris;
2043  CHKERR m_field.get_moab().get_entities_by_dimension(it->meshset, 2, tris,
2044  true);
2045  if (intersect_ptr)
2046  tris = intersect(tris, *intersect_ptr);
2047  CHKERR m_field.add_ents_to_finite_element_by_dim(tris, 2, "PRESSURE_FE");
2048  }
2049  }
2050 
2051  // search for block named LINEAR_PRESSURE and add its attributes to
2052  // PRESSURE_FE element
2053  const string block_set_linear_pressure_name("LINEAR_PRESSURE");
2055  if (it->getName().compare(0, block_set_linear_pressure_name.length(),
2056  block_set_linear_pressure_name) == 0) {
2057  Range tris;
2058  CHKERR m_field.get_moab().get_entities_by_dimension(it->meshset, 2, tris,
2059  true);
2060  if (intersect_ptr)
2061  tris = intersect(tris, *intersect_ptr);
2062  CHKERR m_field.add_ents_to_finite_element_by_dim(tris, 2, "PRESSURE_FE");
2063  }
2064  }
2065 
2067 }

◆ addNeumannFluxBCElements()

MoFEMErrorCode MetaNeumannForces::addNeumannFluxBCElements ( MoFEM::Interface m_field,
const std::string  field_name,
const std::string  mesh_nodals_positions = "MESH_NODE_POSITIONS" 
)
static

Definition at line 2141 of file SurfacePressure.cpp.

2143  {
2145 
2146  CHKERR m_field.add_finite_element("FLUX_FE", MF_ZERO);
2150  if (m_field.check_field(mesh_nodals_positions)) {
2151  CHKERR m_field.modify_finite_element_add_field_data("FLUX_FE",
2152  mesh_nodals_positions);
2153  }
2154 
2156  SIDESET | PRESSURESET, it)) {
2157  Range tris;
2158  CHKERR m_field.get_moab().get_entities_by_dimension(it->meshset, 2, tris,
2159  true);
2160  CHKERR m_field.add_ents_to_finite_element_by_dim(tris, 2, "FLUX_FE");
2161  }
2162 
2164 }

◆ setMassFluxOperators()

MoFEMErrorCode MetaNeumannForces::setMassFluxOperators ( MoFEM::Interface m_field,
boost::ptr_map< std::string, NeumannForcesSurface > &  neumann_forces,
Vec  F,
const std::string  field_name,
const std::string  mesh_nodals_positions = "MESH_NODE_POSITIONS" 
)
static

Definition at line 2166 of file SurfacePressure.cpp.

2169  {
2171  bool ho_geometry = m_field.check_field(mesh_nodals_positions);
2172 
2173  string fe_name;
2174  fe_name = "FLUX_FE";
2175  neumann_forces.insert(fe_name, new NeumannForcesSurface(m_field));
2176 
2177  auto &nf = neumann_forces.at(fe_name);
2178  auto &fe = nf.getLoopFe();
2179 
2180  if (ho_geometry)
2181  CHKERR addHOOpsFace3D(mesh_nodals_positions, fe, false, false);
2182 
2184  SIDESET | PRESSURESET, it)) {
2185  CHKERR nf.addFlux(field_name, F, it->getMeshsetId(), ho_geometry);
2186  }
2188 }

◆ setMomentumFluxOperators()

MoFEMErrorCode MetaNeumannForces::setMomentumFluxOperators ( MoFEM::Interface m_field,
boost::ptr_map< std::string, NeumannForcesSurface > &  neumann_forces,
Vec  F,
const std::string  field_name,
const std::string  mesh_nodals_positions = "MESH_NODE_POSITIONS" 
)
static

Set operators to finite elements calculating right hand side vector.

Parameters
m_fieldInterface
neumann_forcesMap of pointers to force/pressure elements
FRight hand side vector
field_nameField name (f.e. DISPLACEMENT)
mesh_nodals_positionsName of field on which ho-geometry is defined
Returns
Error code
Examples
elasticity.cpp, elasticity_mixed_formulation.cpp, mortar_contact_thermal.cpp, navier_stokes.cpp, simple_contact.cpp, and simple_contact_thermal.cpp.

Definition at line 2069 of file SurfacePressure.cpp.

2072  {
2074  bool ho_geometry = m_field.check_field(mesh_nodals_positions);
2075  // Add forces
2076  {
2077  std::string fe_name = "FORCE_FE";
2078  neumann_forces.insert(fe_name, new NeumannForcesSurface(m_field));
2079  auto &nf = neumann_forces.at(fe_name);
2080  auto &fe = nf.getLoopFe();
2081 
2082  if (ho_geometry)
2083  CHKERR addHOOpsFace3D(mesh_nodals_positions, fe, false, false);
2084 
2086  it)) {
2087  CHKERR nf.addForce(field_name, F, it->getMeshsetId(), ho_geometry, false);
2088  }
2089  // Reading forces from BLOCKSET
2090  const string block_set_force_name("FORCE");
2092  if (it->getName().compare(0, block_set_force_name.length(),
2093  block_set_force_name) == 0) {
2094  CHKERR nf.addForce(field_name, F, it->getMeshsetId(), ho_geometry,
2095  true);
2096  }
2097  }
2098  }
2099 
2100  // Add pressures
2101  {
2102  std::string fe_name = "PRESSURE_FE";
2103  neumann_forces.insert(fe_name, new NeumannForcesSurface(m_field));
2104 
2105  auto &nf = neumann_forces.at(fe_name);
2106  auto &fe = nf.getLoopFe();
2107 
2108  if (ho_geometry)
2109  CHKERR addHOOpsFace3D(mesh_nodals_positions, fe, false, false);
2110 
2112  m_field, SIDESET | PRESSURESET, it)) {
2113  CHKERR nf.addPressure(field_name, F, it->getMeshsetId(), ho_geometry,
2114  false);
2115  }
2116 
2117  // Reading pressures from BLOCKSET
2118  const string block_set_pressure_name("PRESSURE");
2120  if (it->getName().compare(0, block_set_pressure_name.length(),
2121  block_set_pressure_name) == 0) {
2122  CHKERR nf.addPressure(field_name, F, it->getMeshsetId(), ho_geometry,
2123  true);
2124  }
2125  }
2126 
2127  // Reading pressures from BLOCKSET
2128  const string block_set_linear_pressure_name("LINEAR_PRESSURE");
2130  if (it->getName().compare(0, block_set_linear_pressure_name.length(),
2131  block_set_linear_pressure_name) == 0) {
2132  CHKERR nf.addLinearPressure(field_name, F, it->getMeshsetId(),
2133  ho_geometry);
2134  }
2135  }
2136  }
2137 
2139 }

The documentation for this struct was generated from the following files:
MoFEM::addHOOpsFace3D
MoFEMErrorCode addHOOpsFace3D(const std::string field, E &e, bool hcurl, bool hdiv)
Definition: HODataOperators.hpp:699
SIDESET
@ SIDESET
Definition: definitions.h:160
PRESSURESET
@ PRESSURESET
Definition: definitions.h:165
MoFEM::CoreInterface::modify_finite_element_add_field_row
virtual MoFEMErrorCode modify_finite_element_add_field_row(const std::string &fe_name, const std::string name_row)=0
set field row which finite element use
_IT_CUBITMESHSETS_BY_BCDATA_TYPE_FOR_LOOP_
#define _IT_CUBITMESHSETS_BY_BCDATA_TYPE_FOR_LOOP_(MESHSET_MANAGER, CUBITBCTYPE, IT)
Iterator that loops over a specific Cubit MeshSet in a moFEM field.
Definition: MeshsetsManager.hpp:49
NODESET
@ NODESET
Definition: definitions.h:159
CHKERR
#define CHKERR
Inline error check.
Definition: definitions.h:548
MoFEM::CoreInterface::add_finite_element
virtual MoFEMErrorCode add_finite_element(const std::string &fe_name, enum MoFEMTypes bh=MF_EXCL, int verb=DEFAULT_VERBOSITY)=0
add finite element
MoFEM::CoreInterface::get_moab
virtual moab::Interface & get_moab()=0
MoFEM::CoreInterface::modify_finite_element_add_field_col
virtual MoFEMErrorCode modify_finite_element_add_field_col(const std::string &fe_name, const std::string name_row)=0
set field col which finite element use
FORCESET
@ FORCESET
Definition: definitions.h:164
MoFEM::CoreInterface::check_field
virtual bool check_field(const std::string &name) const =0
check if field is in database
MoFEM::CoreInterface::add_ents_to_finite_element_by_dim
virtual MoFEMErrorCode add_ents_to_finite_element_by_dim(const EntityHandle entities, const int dim, const std::string &name, const bool recursive=true)=0
add entities to finite element
MoFEM::CoreInterface::modify_finite_element_add_field_data
virtual MoFEMErrorCode modify_finite_element_add_field_data(const std::string &fe_name, const std::string name_field)=0
set finite element field data
field_name
constexpr auto field_name
Definition: poisson_2d_homogeneous.cpp:13
Range
MF_ZERO
@ MF_ZERO
Definition: definitions.h:111
_IT_CUBITMESHSETS_BY_SET_TYPE_FOR_LOOP_
#define _IT_CUBITMESHSETS_BY_SET_TYPE_FOR_LOOP_(MESHSET_MANAGER, CUBITBCTYPE, IT)
Iterator that loops over a specific Cubit MeshSet having a particular BC meshset in a moFEM field.
Definition: MeshsetsManager.hpp:71
BLOCKSET
@ BLOCKSET
Definition: definitions.h:161
NeumannForcesSurface
Finite element and operators to apply force/pressures applied to surfaces.
Definition: SurfacePressure.hpp:14
MoFEMFunctionReturn
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:429
MoFEMFunctionBegin
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:359
F
@ F
Definition: free_surface.cpp:394