v0.15.0
Loading...
Searching...
No Matches
ConvectiveMassElement Struct Reference

structure grouping operators and data used for calculation of mass (convective) element \ nonlinear_elastic_elem More...

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

Collaboration diagram for ConvectiveMassElement:
[legend]

Classes

struct  BlockData
 data for calculation inertia forces More...
 
struct  CommonData
 common data used by volume elements More...
 
struct  CommonFunctions
 
struct  MatShellCtx
 
struct  MyVolumeFE
 definition of volume element More...
 
struct  OpEnergy
 
struct  OpEshelbyDynamicMaterialMomentumJacobian
 
struct  OpEshelbyDynamicMaterialMomentumLhs_dv
 
struct  OpEshelbyDynamicMaterialMomentumLhs_dX
 
struct  OpEshelbyDynamicMaterialMomentumLhs_dx
 
struct  OpEshelbyDynamicMaterialMomentumRhs
 
struct  OpGetCommonDataAtGaussPts
 
struct  OpGetDataAtGaussPts
 
struct  OpMassJacobian
 
struct  OpMassLhs_dM_dv
 
struct  OpMassLhs_dM_dX
 
struct  OpMassLhs_dM_dx
 
struct  OpMassRhs
 
struct  OpVelocityJacobian
 
struct  OpVelocityLhs_dV_dv
 
struct  OpVelocityLhs_dV_dX
 
struct  OpVelocityLhs_dV_dx
 
struct  OpVelocityRhs
 
struct  PCShellCtx
 
struct  ShellResidualElement
 
struct  UpdateAndControl
 Set fields DOT_. More...
 

Public Member Functions

MyVolumeFEgetLoopFeMassRhs ()
 get rhs volume element
 
MyVolumeFEgetLoopFeMassLhs ()
 get lhs volume element
 
MyVolumeFEgetLoopFeMassAuxLhs ()
 get lhs volume element for Kuu shell matrix
 
MyVolumeFEgetLoopFeVelRhs ()
 get rhs volume element
 
MyVolumeFEgetLoopFeVelLhs ()
 get lhs volume element
 
MyVolumeFEgetLoopFeTRhs ()
 get rhs volume element
 
MyVolumeFEgetLoopFeTLhs ()
 get lhs volume element
 
MyVolumeFEgetLoopFeEnergy ()
 get kinetic energy element
 
MoFEMErrorCode addHOOpsVol ()
 
 ConvectiveMassElement (MoFEM::Interface &m_field, short int tag)
 
MoFEMErrorCode setBlocks ()
 
MoFEMErrorCode addConvectiveMassElement (string element_name, string velocity_field_name, string spatial_position_field_name, string material_position_field_name="MESH_NODE_POSITIONS", bool ale=false, BitRefLevel bit=BitRefLevel())
 
MoFEMErrorCode addVelocityElement (string element_name, string velocity_field_name, string spatial_position_field_name, string material_position_field_name="MESH_NODE_POSITIONS", bool ale=false, BitRefLevel bit=BitRefLevel())
 
MoFEMErrorCode addEshelbyDynamicMaterialMomentum (string element_name, string velocity_field_name, string spatial_position_field_name, string material_position_field_name="MESH_NODE_POSITIONS", bool ale=false, BitRefLevel bit=BitRefLevel(), Range *intersected=NULL)
 
MoFEMErrorCode setConvectiveMassOperators (string velocity_field_name, string spatial_position_field_name, string material_position_field_name="MESH_NODE_POSITIONS", bool ale=false, bool linear=false)
 
MoFEMErrorCode setVelocityOperators (string velocity_field_name, string spatial_position_field_name, string material_position_field_name="MESH_NODE_POSITIONS", bool ale=false)
 
MoFEMErrorCode setKinematicEshelbyOperators (string velocity_field_name, string spatial_position_field_name, string material_position_field_name="MESH_NODE_POSITIONS", Range *forces_on_entities_ptr=NULL)
 
MoFEMErrorCode setShellMatrixMassOperators (string velocity_field_name, string spatial_position_field_name, string material_position_field_name="MESH_NODE_POSITIONS", bool linear=false)
 

Static Public Member Functions

static MoFEMErrorCode setBlocks (MoFEM::Interface &m_field, boost::shared_ptr< map< int, BlockData > > &block_sets_ptr)
 
static MoFEMErrorCode MultOpA (Mat A, Vec x, Vec f)
 Mult operator for shell matrix.
 
static MoFEMErrorCode ZeroEntriesOp (Mat A)
 
static MoFEMErrorCode PCShellSetUpOp (PC pc)
 
static MoFEMErrorCode PCShellDestroy (PC pc)
 
static MoFEMErrorCode PCShellApplyOp (PC pc, Vec f, Vec x)
 apply pre-conditioner for shell matrix
 

Public Attributes

MyVolumeFE feMassRhs
 calculate right hand side for tetrahedral elements
 
MyVolumeFE feMassLhs
 
MyVolumeFE feMassAuxLhs
 
MyVolumeFE feVelRhs
 calculate right hand side for tetrahedral elements
 
MyVolumeFE feVelLhs
 calculate left hand side for tetrahedral elements
 
MyVolumeFE feTRhs
 calculate right hand side for tetrahedral elements
 
MyVolumeFE feTLhs
 calculate left hand side for tetrahedral elements
 
MyVolumeFE feEnergy
 calculate kinetic energy
 
MoFEM::InterfacemField
 
short int tAg
 
std::map< int, BlockDatasetOfBlocks
 maps block set id with appropriate BlockData
 
CommonData commonData
 
boost::ptr_vector< MethodForForceScalingmethodsOp
 

Detailed Description

structure grouping operators and data used for calculation of mass (convective) element \ nonlinear_elastic_elem

structure grouping operators and data used for calculation of nonlinear elastic element

In order to assemble matrices and right hand vectors, the loops over elements, entities over that elements and finally loop over integration points are executed.

Following implementation separate those three celeries of loops and to each loop attach operator.

In order to assemble matrices and right hand vectors, the loops over elements, entities over that elements and finally loop over integration points are executed.

Following implementation separate those three categories of loops and to each loop attach operator.

Examples
HookeElement.hpp, and nonlinear_dynamics.cpp.

Definition at line 28 of file ConvectiveMassElement.hpp.

Constructor & Destructor Documentation

◆ ConvectiveMassElement()

ConvectiveMassElement::ConvectiveMassElement ( MoFEM::Interface & m_field,
short int tag )

Definition at line 75 of file ConvectiveMassElement.cpp.

77 : feMassRhs(m_field), feMassLhs(m_field), feMassAuxLhs(m_field),
78 feVelRhs(m_field), feVelLhs(m_field), feTRhs(m_field), feTLhs(m_field),
79 feEnergy(m_field), mField(m_field), tAg(tag) {}
MyVolumeFE feEnergy
calculate kinetic energy
MyVolumeFE feVelRhs
calculate right hand side for tetrahedral elements
MyVolumeFE feTRhs
calculate right hand side for tetrahedral elements
MyVolumeFE feMassRhs
calculate right hand side for tetrahedral elements
MyVolumeFE feTLhs
calculate left hand side for tetrahedral elements
MyVolumeFE feVelLhs
calculate left hand side for tetrahedral elements

Member Function Documentation

◆ addConvectiveMassElement()

MoFEMErrorCode ConvectiveMassElement::addConvectiveMassElement ( string element_name,
string velocity_field_name,
string spatial_position_field_name,
string material_position_field_name = "MESH_NODE_POSITIONS",
bool ale = false,
BitRefLevel bit = BitRefLevel() )

Definition at line 1836 of file ConvectiveMassElement.cpp.

1839 {
1841
1842 //
1843
1844 CHKERR mField.add_finite_element(element_name, MF_ZERO);
1846 velocity_field_name);
1848 velocity_field_name);
1850 velocity_field_name);
1852 element_name, spatial_position_field_name);
1854 element_name, spatial_position_field_name);
1856 element_name, spatial_position_field_name);
1857 if (mField.check_field(material_position_field_name)) {
1858 if (ale) {
1860 element_name, material_position_field_name);
1862 element_name, material_position_field_name);
1864 element_name, "DOT_" + material_position_field_name);
1865 }
1867 element_name, material_position_field_name);
1868 }
1870 element_name, "DOT_" + velocity_field_name);
1872 element_name, "DOT_" + spatial_position_field_name);
1873
1874 Range tets;
1875 if (bit.any()) {
1876 CHKERR mField.getInterface<BitRefManager>()->getEntitiesByTypeAndRefLevel(
1877 bit, BitRefLevel().set(), MBTET, tets);
1878 }
1879
1880 std::map<int, BlockData>::iterator sit = setOfBlocks.begin();
1881 for (; sit != setOfBlocks.end(); sit++) {
1882 Range add_tets = sit->second.tEts;
1883 if (!tets.empty()) {
1884 add_tets = intersect(add_tets, tets);
1885 }
1887 element_name);
1888 }
1889
1891}
@ MF_ZERO
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
#define CHKERR
Inline error check.
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
virtual MoFEMErrorCode add_finite_element(const std::string &fe_name, enum MoFEMTypes bh=MF_EXCL, int verb=DEFAULT_VERBOSITY)=0
add finite element
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
virtual MoFEMErrorCode add_ents_to_finite_element_by_type(const EntityHandle entities, const EntityType type, const std::string name, const bool recursive=true)=0
add entities to finite element
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
virtual MoFEMErrorCode modify_finite_element_add_field_data(const std::string &fe_name, const std::string name_field)=0
set finite element field data
virtual bool check_field(const std::string &name) const =0
check if field is in database
auto bit
set bit
std::bitset< BITREFLEVEL_SIZE > BitRefLevel
Bit structure attached to each entity identifying to what mesh entity is attached.
Definition Types.hpp:40
std::map< int, BlockData > setOfBlocks
maps block set id with appropriate BlockData
Managing BitRefLevels.
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface reference to pointer of interface.

◆ addEshelbyDynamicMaterialMomentum()

MoFEMErrorCode ConvectiveMassElement::addEshelbyDynamicMaterialMomentum ( string element_name,
string velocity_field_name,
string spatial_position_field_name,
string material_position_field_name = "MESH_NODE_POSITIONS",
bool ale = false,
BitRefLevel bit = BitRefLevel(),
Range * intersected = NULL )

Definition at line 1946 of file ConvectiveMassElement.cpp.

1949 {
1951
1952 //
1953
1954 CHKERR mField.add_finite_element(element_name, MF_ZERO);
1956 velocity_field_name);
1958 velocity_field_name);
1960 element_name, spatial_position_field_name);
1962 element_name, spatial_position_field_name);
1963 if (mField.check_field(material_position_field_name)) {
1964 if (ale) {
1966 element_name, material_position_field_name);
1968 element_name, material_position_field_name);
1970 element_name, "DOT_" + material_position_field_name);
1971 }
1973 element_name, material_position_field_name);
1974 }
1976 element_name, "DOT_" + velocity_field_name);
1978 element_name, "DOT_" + spatial_position_field_name);
1979
1980 Range tets;
1981 if (bit.any()) {
1982 CHKERR mField.getInterface<BitRefManager>()->getEntitiesByTypeAndRefLevel(
1983 bit, BitRefLevel().set(), MBTET, tets);
1984 }
1985 if (intersected != NULL) {
1986 if (tets.empty()) {
1987 tets = *intersected;
1988 } else {
1989 tets = intersect(*intersected, tets);
1990 }
1991 }
1992
1993 std::map<int, BlockData>::iterator sit = setOfBlocks.begin();
1994 for (; sit != setOfBlocks.end(); sit++) {
1995 Range add_tets = sit->second.tEts;
1996 if (!tets.empty()) {
1997 add_tets = intersect(add_tets, tets);
1998 }
2000 element_name);
2001 }
2002
2004}

◆ addHOOpsVol()

MoFEMErrorCode ConvectiveMassElement::addHOOpsVol ( )
inline

Definition at line 92 of file ConvectiveMassElement.hpp.

92 {
94 auto add_ops = [&](auto &fe) {
95 return MoFEM::addHOOpsVol("MESH_NODE_POSITIONS", fe, true, false, false,
96 false);
97 };
98 CHKERR add_ops(feMassRhs);
99 CHKERR add_ops(feMassLhs);
100 CHKERR add_ops(feMassAuxLhs);
101 CHKERR add_ops(feVelRhs);
102 CHKERR add_ops(feTRhs);
103 CHKERR add_ops(feTLhs);
104 CHKERR add_ops(feEnergy);
106 }
#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()
MoFEMErrorCode addHOOpsVol(const std::string field, E &e, bool h1, bool hcurl, bool hdiv, bool l2)

◆ addVelocityElement()

MoFEMErrorCode ConvectiveMassElement::addVelocityElement ( string element_name,
string velocity_field_name,
string spatial_position_field_name,
string material_position_field_name = "MESH_NODE_POSITIONS",
bool ale = false,
BitRefLevel bit = BitRefLevel() )

Definition at line 1893 of file ConvectiveMassElement.cpp.

1896 {
1898
1899 //
1900
1901 CHKERR mField.add_finite_element(element_name, MF_ZERO);
1903 velocity_field_name);
1905 velocity_field_name);
1907 velocity_field_name);
1909 element_name, spatial_position_field_name);
1911 element_name, spatial_position_field_name);
1912 if (mField.check_field(material_position_field_name)) {
1913 if (ale) {
1915 element_name, material_position_field_name);
1917 element_name, "DOT_" + material_position_field_name);
1918 }
1920 element_name, material_position_field_name);
1921 }
1923 element_name, "DOT_" + velocity_field_name);
1925 element_name, "DOT_" + spatial_position_field_name);
1926
1927 Range tets;
1928 if (bit.any()) {
1929 CHKERR mField.getInterface<BitRefManager>()->getEntitiesByTypeAndRefLevel(
1930 bit, BitRefLevel().set(), MBTET, tets);
1931 }
1932
1933 std::map<int, BlockData>::iterator sit = setOfBlocks.begin();
1934 for (; sit != setOfBlocks.end(); sit++) {
1935 Range add_tets = sit->second.tEts;
1936 if (!tets.empty()) {
1937 add_tets = intersect(add_tets, tets);
1938 }
1940 element_name);
1941 }
1942
1944}

◆ getLoopFeEnergy()

MyVolumeFE & ConvectiveMassElement::getLoopFeEnergy ( )
inline

get kinetic energy element

Definition at line 88 of file ConvectiveMassElement.hpp.

◆ getLoopFeMassAuxLhs()

MyVolumeFE & ConvectiveMassElement::getLoopFeMassAuxLhs ( )
inline

get lhs volume element for Kuu shell matrix

Definition at line 73 of file ConvectiveMassElement.hpp.

◆ getLoopFeMassLhs()

MyVolumeFE & ConvectiveMassElement::getLoopFeMassLhs ( )
inline

get lhs volume element

Definition at line 68 of file ConvectiveMassElement.hpp.

◆ getLoopFeMassRhs()

MyVolumeFE & ConvectiveMassElement::getLoopFeMassRhs ( )
inline

get rhs volume element

Definition at line 63 of file ConvectiveMassElement.hpp.

◆ getLoopFeTLhs()

MyVolumeFE & ConvectiveMassElement::getLoopFeTLhs ( )
inline

get lhs volume element

Definition at line 85 of file ConvectiveMassElement.hpp.

◆ getLoopFeTRhs()

MyVolumeFE & ConvectiveMassElement::getLoopFeTRhs ( )
inline

get rhs volume element

Definition at line 83 of file ConvectiveMassElement.hpp.

◆ getLoopFeVelLhs()

MyVolumeFE & ConvectiveMassElement::getLoopFeVelLhs ( )
inline

get lhs volume element

Definition at line 80 of file ConvectiveMassElement.hpp.

◆ getLoopFeVelRhs()

MyVolumeFE & ConvectiveMassElement::getLoopFeVelRhs ( )
inline

get rhs volume element

Definition at line 78 of file ConvectiveMassElement.hpp.

◆ MultOpA()

static MoFEMErrorCode ConvectiveMassElement::MultOpA ( Mat A,
Vec x,
Vec f )
inlinestatic

Mult operator for shell matrix.

\[ \left[ \begin{array}{cc} \mathbf{M} & \mathbf{K} \\ \mathbf{I} & -\mathbf{I}a \end{array} \right] \left[ \begin{array}{c} \mathbf{v} \\ \mathbf{u} \end{array} \right] = \left[ \begin{array}{c} \mathbf{r}_u \\ \mathbf{r}_v \end{array} \right] \]

Examples
nonlinear_dynamics.cpp.

Definition at line 546 of file ConvectiveMassElement.hpp.

546 {
548
549 void *void_ctx;
550 CHKERR MatShellGetContext(A, &void_ctx);
551 MatShellCtx *ctx = (MatShellCtx *)void_ctx;
552 if (!ctx->iNitialized) {
553 CHKERR ctx->iNit();
554 }
555 CHKERR VecZeroEntries(f);
556 // Mult Ku
557 CHKERR VecScatterBegin(ctx->scatterU, x, ctx->u, INSERT_VALUES,
558 SCATTER_FORWARD);
559 CHKERR VecScatterEnd(ctx->scatterU, x, ctx->u, INSERT_VALUES,
560 SCATTER_FORWARD);
561 CHKERR MatMult(ctx->K, ctx->u, ctx->Ku);
562 CHKERR VecScatterBegin(ctx->scatterU, ctx->Ku, f, INSERT_VALUES,
563 SCATTER_REVERSE);
564 CHKERR VecScatterEnd(ctx->scatterU, ctx->Ku, f, INSERT_VALUES,
565 SCATTER_REVERSE);
566 // Mult Mv
567 CHKERR VecScatterBegin(ctx->scatterV, x, ctx->v, INSERT_VALUES,
568 SCATTER_FORWARD);
569 CHKERR VecScatterEnd(ctx->scatterV, x, ctx->v, INSERT_VALUES,
570 SCATTER_FORWARD);
571 CHKERR MatMult(ctx->M, ctx->v, ctx->Mv);
572 CHKERR VecScatterBegin(ctx->scatterU, ctx->Mv, f, ADD_VALUES,
573 SCATTER_REVERSE);
574 CHKERR VecScatterEnd(ctx->scatterU, ctx->Mv, f, ADD_VALUES,
575 SCATTER_REVERSE);
576 // Velocities
577 CHKERR VecAXPY(ctx->v, -ctx->ts_a, ctx->u);
578 // CHKERR VecScale(ctx->v,ctx->scale);
579 CHKERR VecScatterBegin(ctx->scatterV, ctx->v, f, INSERT_VALUES,
580 SCATTER_REVERSE);
581 CHKERR VecScatterEnd(ctx->scatterV, ctx->v, f, INSERT_VALUES,
582 SCATTER_REVERSE);
583 // Assemble
584 CHKERR VecAssemblyBegin(f);
585 CHKERR VecAssemblyEnd(f);
587 }

◆ PCShellApplyOp()

static MoFEMErrorCode ConvectiveMassElement::PCShellApplyOp ( PC pc,
Vec f,
Vec x )
inlinestatic

apply pre-conditioner for shell matrix

\[ \left[ \begin{array}{cc} \mathbf{M} & \mathbf{K} \\ \mathbf{I} & -\mathbf{I}a \end{array} \right] \left[ \begin{array}{c} \mathbf{v} \\ \mathbf{u} \end{array} \right] = \left[ \begin{array}{c} \mathbf{r}_u \\ \mathbf{r}_v \end{array} \right] \]

where \(\mathbf{v} = \mathbf{r}_v + a\mathbf{u}\) and \(\mathbf{u}=(a\mathbf{M}+\mathbf{K})^{-1}(\mathbf{r}_u - \mathbf{M}\mathbf{r}_v\).

Examples
nonlinear_dynamics.cpp.

Definition at line 671 of file ConvectiveMassElement.hpp.

671 {
673
674 void *void_ctx;
675 CHKERR PCShellGetContext(pc, &void_ctx);
676 PCShellCtx *ctx = (PCShellCtx *)void_ctx;
677 MatShellCtx *shell_mat_ctx;
678 CHKERR MatShellGetContext(ctx->shellMat, &shell_mat_ctx);
679 // forward
680 CHKERR VecScatterBegin(shell_mat_ctx->scatterU, f, shell_mat_ctx->Ku,
681 INSERT_VALUES, SCATTER_FORWARD);
682 CHKERR VecScatterEnd(shell_mat_ctx->scatterU, f, shell_mat_ctx->Ku,
683 INSERT_VALUES, SCATTER_FORWARD);
684 CHKERR VecScatterBegin(shell_mat_ctx->scatterV, f, shell_mat_ctx->v,
685 INSERT_VALUES, SCATTER_FORWARD);
686 CHKERR VecScatterEnd(shell_mat_ctx->scatterV, f, shell_mat_ctx->v,
687 INSERT_VALUES, SCATTER_FORWARD);
688 // CHKERR VecScale(shell_mat_ctx->v,1/shell_mat_ctx->scale);
689 // apply pre-conditioner and calculate u
690 CHKERR MatMult(shell_mat_ctx->M, shell_mat_ctx->v,
691 shell_mat_ctx->Mv); // Mrv
692 CHKERR VecAXPY(shell_mat_ctx->Ku, -1, shell_mat_ctx->Mv); // f-Mrv
693 CHKERR PCApply(ctx->pC, shell_mat_ctx->Ku,
694 shell_mat_ctx->u); // u = (aM+K)^(-1)(ru-Mrv)
695 // VecView(shell_mat_ctx->u,PETSC_VIEWER_STDOUT_WORLD);
696 // calculate velocities
697 CHKERR VecAXPY(shell_mat_ctx->v, shell_mat_ctx->ts_a,
698 shell_mat_ctx->u); // v = v + a*u
699 // VecView(shell_mat_ctx->v,PETSC_VIEWER_STDOUT_WORLD);
700 // reverse
701 CHKERR VecZeroEntries(x);
702 CHKERR VecScatterBegin(shell_mat_ctx->scatterU, shell_mat_ctx->u, x,
703 INSERT_VALUES, SCATTER_REVERSE);
704 CHKERR VecScatterEnd(shell_mat_ctx->scatterU, shell_mat_ctx->u, x,
705 INSERT_VALUES, SCATTER_REVERSE);
706 CHKERR VecScatterBegin(shell_mat_ctx->scatterV, shell_mat_ctx->v, x,
707 INSERT_VALUES, SCATTER_REVERSE);
708 CHKERR VecScatterEnd(shell_mat_ctx->scatterV, shell_mat_ctx->v, x,
709 INSERT_VALUES, SCATTER_REVERSE);
710 CHKERR VecAssemblyBegin(x);
711 CHKERR VecAssemblyEnd(x);
713 }

◆ PCShellDestroy()

static MoFEMErrorCode ConvectiveMassElement::PCShellDestroy ( PC pc)
inlinestatic
Examples
nonlinear_dynamics.cpp.

Definition at line 633 of file ConvectiveMassElement.hpp.

633 {
635
636 void *void_ctx;
637 CHKERR PCShellGetContext(pc, &void_ctx);
638 PCShellCtx *ctx = (PCShellCtx *)void_ctx;
639 CHKERR ctx->dEstroy();
641 }

◆ PCShellSetUpOp()

static MoFEMErrorCode ConvectiveMassElement::PCShellSetUpOp ( PC pc)
inlinestatic
Examples
nonlinear_dynamics.cpp.

Definition at line 618 of file ConvectiveMassElement.hpp.

618 {
620
621 void *void_ctx;
622 CHKERR PCShellGetContext(pc, &void_ctx);
623 PCShellCtx *ctx = (PCShellCtx *)void_ctx;
624 CHKERR ctx->iNit();
625 MatShellCtx *shell_mat_ctx;
626 CHKERR MatShellGetContext(ctx->shellMat, &shell_mat_ctx);
627 CHKERR PCSetFromOptions(ctx->pC);
628 CHKERR PCSetOperators(ctx->pC, shell_mat_ctx->barK, shell_mat_ctx->barK);
629 CHKERR PCSetUp(ctx->pC);
631 }

◆ setBlocks() [1/2]

MoFEMErrorCode ConvectiveMassElement::setBlocks ( )
Examples
elasticity.cpp.

Definition at line 1762 of file ConvectiveMassElement.cpp.

1762 {
1764
1765 Range added_tets;
1767 mField, BLOCKSET | BODYFORCESSET, it)) {
1768 int id = it->getMeshsetId();
1769 EntityHandle meshset = it->getMeshset();
1770 CHKERR mField.get_moab().get_entities_by_type(meshset, MBTET,
1771 setOfBlocks[id].tEts, true);
1772 added_tets.merge(setOfBlocks[id].tEts);
1773 Block_BodyForces mydata;
1774 CHKERR it->getAttributeDataStructure(mydata);
1775 setOfBlocks[id].rho0 = mydata.data.density;
1776 setOfBlocks[id].a0.resize(3);
1777 setOfBlocks[id].a0[0] = mydata.data.acceleration_x;
1778 setOfBlocks[id].a0[1] = mydata.data.acceleration_y;
1779 setOfBlocks[id].a0[2] = mydata.data.acceleration_z;
1780 // std::cerr << setOfBlocks[id].tEts << std::endl;
1781 }
1782
1784 mField, BLOCKSET | MAT_ELASTICSET, it)) {
1785 Mat_Elastic mydata;
1786 CHKERR it->getAttributeDataStructure(mydata);
1787 if (mydata.data.User1 == 0)
1788 continue;
1789 Range tets;
1790 EntityHandle meshset = it->getMeshset();
1791 CHKERR mField.get_moab().get_entities_by_type(meshset, MBTET, tets, true);
1792 tets = subtract(tets, added_tets);
1793 if (tets.empty())
1794 continue;
1795 int id = it->getMeshsetId();
1796 setOfBlocks[-id].tEts = tets;
1797 setOfBlocks[-id].rho0 = mydata.data.User1;
1798 setOfBlocks[-id].a0.resize(3);
1799 setOfBlocks[-id].a0[0] = mydata.data.User2;
1800 setOfBlocks[-id].a0[1] = mydata.data.User3;
1801 setOfBlocks[-id].a0[2] = mydata.data.User4;
1802 // std::cerr << setOfBlocks[id].tEts << std::endl;
1803 }
1804
1806}
@ BODYFORCESSET
block name is "BODY_FORCES"
@ MAT_ELASTICSET
block name is "MAT_ELASTIC"
@ BLOCKSET
#define _IT_CUBITMESHSETS_BY_BCDATA_TYPE_FOR_LOOP_(MESHSET_MANAGER, CUBITBCTYPE, IT)
Iterator that loops over a specific Cubit MeshSet in a moFEM field.
Body force data structure.
virtual moab::Interface & get_moab()=0
Elastic material data structure.

◆ setBlocks() [2/2]

MoFEMErrorCode ConvectiveMassElement::setBlocks ( MoFEM::Interface & m_field,
boost::shared_ptr< map< int, BlockData > > & block_sets_ptr )
static

Definition at line 1808 of file ConvectiveMassElement.cpp.

1810 {
1812
1813 if (!block_sets_ptr)
1814 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
1815 "Pointer to block of sets is null");
1816
1818 m_field, BLOCKSET | BODYFORCESSET, it)) {
1819 Block_BodyForces mydata;
1820 CHKERR it->getAttributeDataStructure(mydata);
1821 int id = it->getMeshsetId();
1822 auto &block_data = (*block_sets_ptr)[id];
1823 EntityHandle meshset = it->getMeshset();
1824 CHKERR m_field.get_moab().get_entities_by_dimension(meshset, 3,
1825 block_data.tEts, true);
1826 block_data.rho0 = mydata.data.density;
1827 block_data.a0.resize(3);
1828 block_data.a0[0] = mydata.data.acceleration_x;
1829 block_data.a0[1] = mydata.data.acceleration_y;
1830 block_data.a0[2] = mydata.data.acceleration_z;
1831 }
1832
1834}
@ MOFEM_DATA_INCONSISTENCY
Definition definitions.h:31

◆ setConvectiveMassOperators()

MoFEMErrorCode ConvectiveMassElement::setConvectiveMassOperators ( string velocity_field_name,
string spatial_position_field_name,
string material_position_field_name = "MESH_NODE_POSITIONS",
bool ale = false,
bool linear = false )

Definition at line 2006 of file ConvectiveMassElement.cpp.

2008 {
2010
2011 commonData.spatialPositions = spatial_position_field_name;
2012 commonData.meshPositions = material_position_field_name;
2013 commonData.spatialVelocities = velocity_field_name;
2014 commonData.lInear = linear;
2015
2016 // Rhs
2017 feMassRhs.getOpPtrVector().push_back(
2018 new OpGetCommonDataAtGaussPts(velocity_field_name, commonData));
2019 feMassRhs.getOpPtrVector().push_back(
2020 new OpGetCommonDataAtGaussPts(spatial_position_field_name, commonData));
2021 feMassRhs.getOpPtrVector().push_back(
2022 new OpGetCommonDataAtGaussPts("DOT_" + velocity_field_name, commonData));
2023 feMassRhs.getOpPtrVector().push_back(new OpGetCommonDataAtGaussPts(
2024 "DOT_" + spatial_position_field_name, commonData));
2025 if (mField.check_field(material_position_field_name)) {
2026 feMassRhs.getOpPtrVector().push_back(new OpGetCommonDataAtGaussPts(
2027 material_position_field_name, commonData));
2028 if (ale) {
2029 feMassRhs.getOpPtrVector().push_back(new OpGetCommonDataAtGaussPts(
2030 "DOT_" + material_position_field_name, commonData));
2031 } else {
2032 feMassRhs.meshPositionsFieldName = material_position_field_name;
2033 }
2034 }
2035 std::map<int, BlockData>::iterator sit = setOfBlocks.begin();
2036 for (; sit != setOfBlocks.end(); sit++) {
2037 feMassRhs.getOpPtrVector().push_back(
2038 new OpMassJacobian(spatial_position_field_name, sit->second, commonData,
2039 methodsOp, tAg, false));
2040 feMassRhs.getOpPtrVector().push_back(
2041 new OpMassRhs(spatial_position_field_name, sit->second, commonData));
2042 }
2043
2044 // Lhs
2045 feMassLhs.getOpPtrVector().push_back(
2046 new OpGetCommonDataAtGaussPts(velocity_field_name, commonData));
2047 feMassLhs.getOpPtrVector().push_back(
2048 new OpGetCommonDataAtGaussPts(spatial_position_field_name, commonData));
2049 feMassLhs.getOpPtrVector().push_back(
2050 new OpGetCommonDataAtGaussPts("DOT_" + velocity_field_name, commonData));
2051 feMassLhs.getOpPtrVector().push_back(new OpGetCommonDataAtGaussPts(
2052 "DOT_" + spatial_position_field_name, commonData));
2053 if (mField.check_field(material_position_field_name)) {
2054 feMassLhs.getOpPtrVector().push_back(new OpGetCommonDataAtGaussPts(
2055 material_position_field_name, commonData));
2056 if (ale) {
2057 feMassLhs.getOpPtrVector().push_back(new OpGetCommonDataAtGaussPts(
2058 "DOT_" + material_position_field_name, commonData));
2059 } else {
2060 feMassLhs.meshPositionsFieldName = material_position_field_name;
2061 }
2062 }
2063 sit = setOfBlocks.begin();
2064 for (; sit != setOfBlocks.end(); sit++) {
2065 feMassLhs.getOpPtrVector().push_back(
2066 new OpMassJacobian(spatial_position_field_name, sit->second, commonData,
2067 methodsOp, tAg, true));
2068 feMassLhs.getOpPtrVector().push_back(
2069 new OpMassLhs_dM_dv(spatial_position_field_name, velocity_field_name,
2070 sit->second, commonData));
2071 feMassLhs.getOpPtrVector().push_back(new OpMassLhs_dM_dx(
2072 spatial_position_field_name, spatial_position_field_name, sit->second,
2073 commonData));
2074 if (mField.check_field(material_position_field_name)) {
2075 if (ale) {
2076 feMassLhs.getOpPtrVector().push_back(new OpMassLhs_dM_dX(
2077 spatial_position_field_name, material_position_field_name,
2078 sit->second, commonData));
2079 } else {
2080 feMassLhs.meshPositionsFieldName = material_position_field_name;
2081 }
2082 }
2083 }
2084
2085 // Energy
2086 feEnergy.getOpPtrVector().push_back(
2087 new OpGetCommonDataAtGaussPts(velocity_field_name, commonData));
2088 feEnergy.getOpPtrVector().push_back(
2089 new OpGetCommonDataAtGaussPts(spatial_position_field_name, commonData));
2090 if (mField.check_field(material_position_field_name)) {
2091 feEnergy.getOpPtrVector().push_back(new OpGetCommonDataAtGaussPts(
2092 material_position_field_name, commonData));
2093 feEnergy.meshPositionsFieldName = material_position_field_name;
2094 }
2095 sit = setOfBlocks.begin();
2096 for (; sit != setOfBlocks.end(); sit++) {
2097 feEnergy.getOpPtrVector().push_back(new OpEnergy(
2098 spatial_position_field_name, sit->second, commonData, feEnergy.V));
2099 }
2100
2102}
boost::ptr_vector< MethodForForceScaling > methodsOp
boost::ptr_deque< UserDataOperator > & getOpPtrVector()
Use to push back operator for row operator.

◆ setKinematicEshelbyOperators()

MoFEMErrorCode ConvectiveMassElement::setKinematicEshelbyOperators ( string velocity_field_name,
string spatial_position_field_name,
string material_position_field_name = "MESH_NODE_POSITIONS",
Range * forces_on_entities_ptr = NULL )

Definition at line 2182 of file ConvectiveMassElement.cpp.

2184 {
2186
2187 commonData.spatialPositions = spatial_position_field_name;
2188 commonData.meshPositions = material_position_field_name;
2189 commonData.spatialVelocities = velocity_field_name;
2190
2191 // Rhs
2192 feTRhs.getOpPtrVector().push_back(
2193 new OpGetCommonDataAtGaussPts(velocity_field_name, commonData));
2194 feTRhs.getOpPtrVector().push_back(
2195 new OpGetCommonDataAtGaussPts(spatial_position_field_name, commonData));
2196 feTRhs.getOpPtrVector().push_back(
2197 new OpGetCommonDataAtGaussPts(material_position_field_name, commonData));
2198 feTRhs.getOpPtrVector().push_back(
2199 new OpGetCommonDataAtGaussPts("DOT_" + velocity_field_name, commonData));
2200
2201 std::map<int, BlockData>::iterator sit = setOfBlocks.begin();
2202 for (; sit != setOfBlocks.end(); sit++) {
2203 feTRhs.getOpPtrVector().push_back(
2204 new OpEshelbyDynamicMaterialMomentumJacobian(
2205 material_position_field_name, sit->second, commonData, tAg, false));
2206 feTRhs.getOpPtrVector().push_back(new OpEshelbyDynamicMaterialMomentumRhs(
2207 material_position_field_name, sit->second, commonData,
2208 forces_on_entities_ptr));
2209 }
2210
2211 // Lhs
2212 feTLhs.getOpPtrVector().push_back(
2213 new OpGetCommonDataAtGaussPts(velocity_field_name, commonData));
2214 feTLhs.getOpPtrVector().push_back(
2215 new OpGetCommonDataAtGaussPts(spatial_position_field_name, commonData));
2216 feTLhs.getOpPtrVector().push_back(
2217 new OpGetCommonDataAtGaussPts("DOT_" + velocity_field_name, commonData));
2218 if (mField.check_field(material_position_field_name)) {
2219 feTLhs.getOpPtrVector().push_back(new OpGetCommonDataAtGaussPts(
2220 material_position_field_name, commonData));
2221 }
2222 sit = setOfBlocks.begin();
2223 for (; sit != setOfBlocks.end(); sit++) {
2224 feTLhs.getOpPtrVector().push_back(
2225 new OpEshelbyDynamicMaterialMomentumJacobian(
2226 material_position_field_name, sit->second, commonData, tAg));
2227 feTLhs.getOpPtrVector().push_back(
2228 new OpEshelbyDynamicMaterialMomentumLhs_dv(
2229 material_position_field_name, velocity_field_name, sit->second,
2230 commonData, forces_on_entities_ptr));
2231 feTLhs.getOpPtrVector().push_back(
2232 new OpEshelbyDynamicMaterialMomentumLhs_dx(
2233 material_position_field_name, spatial_position_field_name,
2234 sit->second, commonData, forces_on_entities_ptr));
2235 feTLhs.getOpPtrVector().push_back(
2236 new OpEshelbyDynamicMaterialMomentumLhs_dX(
2237 material_position_field_name, material_position_field_name,
2238 sit->second, commonData, forces_on_entities_ptr));
2239 }
2240
2242}

◆ setShellMatrixMassOperators()

MoFEMErrorCode ConvectiveMassElement::setShellMatrixMassOperators ( string velocity_field_name,
string spatial_position_field_name,
string material_position_field_name = "MESH_NODE_POSITIONS",
bool linear = false )

Definition at line 2244 of file ConvectiveMassElement.cpp.

2246 {
2248
2249 commonData.spatialPositions = spatial_position_field_name;
2250 commonData.meshPositions = material_position_field_name;
2251 commonData.spatialVelocities = velocity_field_name;
2252 commonData.lInear = linear;
2253
2254 // Rhs
2255 feMassRhs.getOpPtrVector().push_back(
2256 new OpGetCommonDataAtGaussPts(velocity_field_name, commonData));
2257 feMassRhs.getOpPtrVector().push_back(
2258 new OpGetCommonDataAtGaussPts(spatial_position_field_name, commonData));
2259 feMassRhs.getOpPtrVector().push_back(
2260 new OpGetCommonDataAtGaussPts("DOT_" + velocity_field_name, commonData));
2261 if (mField.check_field(material_position_field_name)) {
2262 feMassRhs.getOpPtrVector().push_back(new OpGetCommonDataAtGaussPts(
2263 material_position_field_name, commonData));
2264 feMassRhs.meshPositionsFieldName = material_position_field_name;
2265 }
2266 std::map<int, BlockData>::iterator sit = setOfBlocks.begin();
2267 for (; sit != setOfBlocks.end(); sit++) {
2268 feMassRhs.getOpPtrVector().push_back(
2269 new OpMassJacobian(spatial_position_field_name, sit->second, commonData,
2270 methodsOp, tAg, false));
2271 feMassRhs.getOpPtrVector().push_back(
2272 new OpMassRhs(spatial_position_field_name, sit->second, commonData));
2273 }
2274
2275 // Lhs
2276 feMassLhs.getOpPtrVector().push_back(
2277 new OpGetCommonDataAtGaussPts(velocity_field_name, commonData));
2278 feMassLhs.getOpPtrVector().push_back(
2279 new OpGetCommonDataAtGaussPts(spatial_position_field_name, commonData));
2280 feMassLhs.getOpPtrVector().push_back(
2281 new OpGetCommonDataAtGaussPts("DOT_" + velocity_field_name, commonData));
2282 if (mField.check_field(material_position_field_name)) {
2283 feMassLhs.getOpPtrVector().push_back(new OpGetCommonDataAtGaussPts(
2284 material_position_field_name, commonData));
2285 feMassLhs.meshPositionsFieldName = material_position_field_name;
2286 }
2287 sit = setOfBlocks.begin();
2288 for (; sit != setOfBlocks.end(); sit++) {
2289 feMassLhs.getOpPtrVector().push_back(
2290 new OpMassJacobian(spatial_position_field_name, sit->second, commonData,
2291 methodsOp, tAg, true));
2292 feMassLhs.getOpPtrVector().push_back(new OpMassLhs_dM_dv(
2293 spatial_position_field_name, spatial_position_field_name, sit->second,
2294 commonData));
2295 if (mField.check_field(material_position_field_name)) {
2296 feMassLhs.meshPositionsFieldName = material_position_field_name;
2297 }
2298 }
2299
2300 // Aux Lhs
2301 feMassAuxLhs.getOpPtrVector().push_back(
2302 new OpGetCommonDataAtGaussPts(velocity_field_name, commonData));
2303 feMassAuxLhs.getOpPtrVector().push_back(
2304 new OpGetCommonDataAtGaussPts(spatial_position_field_name, commonData));
2305 feMassAuxLhs.getOpPtrVector().push_back(
2306 new OpGetCommonDataAtGaussPts("DOT_" + velocity_field_name, commonData));
2307 if (mField.check_field(material_position_field_name)) {
2308 feMassAuxLhs.getOpPtrVector().push_back(new OpGetCommonDataAtGaussPts(
2309 material_position_field_name, commonData));
2310 feMassAuxLhs.meshPositionsFieldName = material_position_field_name;
2311 }
2312 sit = setOfBlocks.begin();
2313 for (; sit != setOfBlocks.end(); sit++) {
2314 feMassAuxLhs.getOpPtrVector().push_back(
2315 new OpMassJacobian(spatial_position_field_name, sit->second, commonData,
2316 methodsOp, tAg, true));
2317 feMassAuxLhs.getOpPtrVector().push_back(new OpMassLhs_dM_dx(
2318 spatial_position_field_name, spatial_position_field_name, sit->second,
2319 commonData));
2320 if (mField.check_field(material_position_field_name)) {
2321 feMassAuxLhs.meshPositionsFieldName = material_position_field_name;
2322 }
2323 }
2324
2325 // Energy E=0.5*rho*v*v
2326 feEnergy.getOpPtrVector().push_back(
2327 new OpGetCommonDataAtGaussPts(velocity_field_name, commonData));
2328 feEnergy.getOpPtrVector().push_back(
2329 new OpGetCommonDataAtGaussPts(spatial_position_field_name, commonData));
2330 if (mField.check_field(material_position_field_name)) {
2331 feEnergy.getOpPtrVector().push_back(new OpGetCommonDataAtGaussPts(
2332 material_position_field_name, commonData));
2333 feEnergy.meshPositionsFieldName = material_position_field_name;
2334 }
2335 sit = setOfBlocks.begin();
2336 for (; sit != setOfBlocks.end(); sit++) {
2337 feEnergy.getOpPtrVector().push_back(new OpEnergy(
2338 spatial_position_field_name, sit->second, commonData, feEnergy.V));
2339 }
2340
2342}

◆ setVelocityOperators()

MoFEMErrorCode ConvectiveMassElement::setVelocityOperators ( string velocity_field_name,
string spatial_position_field_name,
string material_position_field_name = "MESH_NODE_POSITIONS",
bool ale = false )

Definition at line 2104 of file ConvectiveMassElement.cpp.

2106 {
2108
2109 commonData.spatialPositions = spatial_position_field_name;
2110 commonData.meshPositions = material_position_field_name;
2111 commonData.spatialVelocities = velocity_field_name;
2112
2113 // Rhs
2114 feVelRhs.getOpPtrVector().push_back(
2115 new OpGetCommonDataAtGaussPts(velocity_field_name, commonData));
2116 feVelRhs.getOpPtrVector().push_back(
2117 new OpGetCommonDataAtGaussPts(spatial_position_field_name, commonData));
2118 feVelRhs.getOpPtrVector().push_back(
2119 new OpGetCommonDataAtGaussPts("DOT_" + velocity_field_name, commonData));
2120 if (mField.check_field(material_position_field_name)) {
2121 feVelRhs.getOpPtrVector().push_back(new OpGetCommonDataAtGaussPts(
2122 "DOT_" + spatial_position_field_name, commonData));
2123 feVelRhs.getOpPtrVector().push_back(new OpGetCommonDataAtGaussPts(
2124 material_position_field_name, commonData));
2125 if (ale) {
2126 feVelRhs.getOpPtrVector().push_back(new OpGetCommonDataAtGaussPts(
2127 "DOT_" + material_position_field_name, commonData));
2128 } else {
2129 feVelRhs.meshPositionsFieldName = material_position_field_name;
2130 }
2131 }
2132 std::map<int, BlockData>::iterator sit = setOfBlocks.begin();
2133 for (; sit != setOfBlocks.end(); sit++) {
2134 feVelRhs.getOpPtrVector().push_back(new OpVelocityJacobian(
2135 velocity_field_name, sit->second, commonData, tAg, false));
2136 feVelRhs.getOpPtrVector().push_back(
2137 new OpVelocityRhs(velocity_field_name, sit->second, commonData));
2138 }
2139
2140 // Lhs
2141 feVelLhs.getOpPtrVector().push_back(
2142 new OpGetCommonDataAtGaussPts(velocity_field_name, commonData));
2143 feVelLhs.getOpPtrVector().push_back(
2144 new OpGetCommonDataAtGaussPts(spatial_position_field_name, commonData));
2145 feVelLhs.getOpPtrVector().push_back(
2146 new OpGetCommonDataAtGaussPts("DOT_" + velocity_field_name, commonData));
2147 if (mField.check_field(material_position_field_name)) {
2148 feVelLhs.getOpPtrVector().push_back(new OpGetCommonDataAtGaussPts(
2149 "DOT_" + spatial_position_field_name, commonData));
2150 feVelLhs.getOpPtrVector().push_back(new OpGetCommonDataAtGaussPts(
2151 material_position_field_name, commonData));
2152 if (ale) {
2153 feVelLhs.getOpPtrVector().push_back(new OpGetCommonDataAtGaussPts(
2154 "DOT_" + material_position_field_name, commonData));
2155 } else {
2156 feVelLhs.meshPositionsFieldName = material_position_field_name;
2157 }
2158 }
2159 sit = setOfBlocks.begin();
2160 for (; sit != setOfBlocks.end(); sit++) {
2161 feVelLhs.getOpPtrVector().push_back(new OpVelocityJacobian(
2162 velocity_field_name, sit->second, commonData, tAg));
2163 feVelLhs.getOpPtrVector().push_back(new OpVelocityLhs_dV_dv(
2164 velocity_field_name, velocity_field_name, sit->second, commonData));
2165 feVelLhs.getOpPtrVector().push_back(new OpVelocityLhs_dV_dx(
2166 velocity_field_name, spatial_position_field_name, sit->second,
2167 commonData));
2168 if (mField.check_field(material_position_field_name)) {
2169 if (ale) {
2170 feVelLhs.getOpPtrVector().push_back(new OpVelocityLhs_dV_dX(
2171 velocity_field_name, material_position_field_name, sit->second,
2172 commonData));
2173 } else {
2174 feVelLhs.meshPositionsFieldName = material_position_field_name;
2175 }
2176 }
2177 }
2178
2180}

◆ ZeroEntriesOp()

static MoFEMErrorCode ConvectiveMassElement::ZeroEntriesOp ( Mat A)
inlinestatic
Examples
nonlinear_dynamics.cpp.

Definition at line 589 of file ConvectiveMassElement.hpp.

589 {
591
592 void *void_ctx;
593 CHKERR MatShellGetContext(A, &void_ctx);
594 MatShellCtx *ctx = (MatShellCtx *)void_ctx;
595 CHKERR MatZeroEntries(ctx->K);
596 CHKERR MatZeroEntries(ctx->M);
598 }

Member Data Documentation

◆ commonData

CommonData ConvectiveMassElement::commonData

Definition at line 148 of file ConvectiveMassElement.hpp.

◆ feEnergy

MyVolumeFE ConvectiveMassElement::feEnergy

calculate kinetic energy

Definition at line 87 of file ConvectiveMassElement.hpp.

◆ feMassAuxLhs

MyVolumeFE ConvectiveMassElement::feMassAuxLhs

calculate left hand side for tetrahedral elements for Kuu shell matrix

Definition at line 71 of file ConvectiveMassElement.hpp.

◆ feMassLhs

MyVolumeFE ConvectiveMassElement::feMassLhs

calculate left hand side for tetrahedral elements,i.e. mass element

Definition at line 66 of file ConvectiveMassElement.hpp.

◆ feMassRhs

MyVolumeFE ConvectiveMassElement::feMassRhs

calculate right hand side for tetrahedral elements

Definition at line 62 of file ConvectiveMassElement.hpp.

◆ feTLhs

MyVolumeFE ConvectiveMassElement::feTLhs

calculate left hand side for tetrahedral elements

Definition at line 84 of file ConvectiveMassElement.hpp.

◆ feTRhs

MyVolumeFE ConvectiveMassElement::feTRhs

calculate right hand side for tetrahedral elements

Definition at line 82 of file ConvectiveMassElement.hpp.

◆ feVelLhs

MyVolumeFE ConvectiveMassElement::feVelLhs

calculate left hand side for tetrahedral elements

Definition at line 79 of file ConvectiveMassElement.hpp.

◆ feVelRhs

MyVolumeFE ConvectiveMassElement::feVelRhs

calculate right hand side for tetrahedral elements

Definition at line 77 of file ConvectiveMassElement.hpp.

◆ methodsOp

boost::ptr_vector<MethodForForceScaling> ConvectiveMassElement::methodsOp

Definition at line 150 of file ConvectiveMassElement.hpp.

◆ mField

MoFEM::Interface& ConvectiveMassElement::mField

Definition at line 108 of file ConvectiveMassElement.hpp.

◆ setOfBlocks

std::map<int, BlockData> ConvectiveMassElement::setOfBlocks

maps block set id with appropriate BlockData

Definition at line 122 of file ConvectiveMassElement.hpp.

◆ tAg

short int ConvectiveMassElement::tAg

Definition at line 109 of file ConvectiveMassElement.hpp.


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