v0.14.0
Loading...
Searching...
No Matches
Public Member Functions | Private Attributes | List of all members
ElectroPhysiology::OpEssentialBC Struct Reference

#include <users_modules/softmech/chemo_mech/src/ElecPhysOperators.hpp>

Inheritance diagram for ElectroPhysiology::OpEssentialBC:
[legend]
Collaboration diagram for ElectroPhysiology::OpEssentialBC:
[legend]

Public Member Functions

 OpEssentialBC (std::string flux_field, Range &essential_bd_ents)
 
MoFEMErrorCode doWork (int side, EntityType type, EntData &data)
 
 OpEssentialBC (const std::string &flux_field, Range &essential_bd_ents)
 
MoFEMErrorCode doWork (int side, EntityType type, EntData &data)
 
- Public Member Functions inherited from MoFEM::FaceElementForcesAndSourcesCore::UserDataOperator
double getArea ()
 get area of face More...
 
VectorDoublegetNormal ()
 get triangle normal More...
 
VectorDoublegetTangent1 ()
 get triangle tangent 1 More...
 
VectorDoublegetTangent2 ()
 get triangle tangent 2 More...
 
auto getFTensor1Normal ()
 get normal as tensor More...
 
auto getFTensor1Tangent1 ()
 get tangentOne as tensor More...
 
auto getFTensor1Tangent2 ()
 get tangentTwo as tensor More...
 
int getNumNodes ()
 get element number of nodes More...
 
const EntityHandlegetConn ()
 get element connectivity More...
 
VectorDoublegetCoords ()
 get triangle coordinates More...
 
auto getFTensor1Coords ()
 get get coords at gauss points More...
 
MatrixDoublegetNormalsAtGaussPts ()
 if higher order geometry return normals at Gauss pts. More...
 
ublas::matrix_row< MatrixDoublegetNormalsAtGaussPts (const int gg)
 if higher order geometry return normals at Gauss pts. More...
 
MatrixDoublegetTangent1AtGaussPts ()
 if higher order geometry return tangent vector to triangle at Gauss pts. More...
 
MatrixDoublegetTangent2AtGaussPts ()
 if higher order geometry return tangent vector to triangle at Gauss pts. More...
 
auto getFTensor1NormalsAtGaussPts ()
 get normal at integration points More...
 
auto getFTensor1Tangent1AtGaussPts ()
 get tangent 1 at integration points More...
 
auto getFTensor1Tangent2AtGaussPts ()
 get tangent 2 at integration points More...
 
FaceElementForcesAndSourcesCoregetFaceFE ()
 return pointer to Generic Triangle Finite Element object More...
 
MoFEMErrorCode loopSideVolumes (const string fe_name, VolumeElementForcesAndSourcesCoreOnSide &fe_method)
 
- Public Member Functions inherited from MoFEM::ForcesAndSourcesCore::UserDataOperator
 UserDataOperator (const FieldSpace space, const char type=OPSPACE, const bool symm=true)
 
 UserDataOperator (const std::string field_name, const char type, const bool symm=true)
 
 UserDataOperator (const std::string row_field_name, const std::string col_field_name, const char type, const bool symm=true)
 
boost::shared_ptr< const NumeredEntFiniteElementgetNumeredEntFiniteElementPtr () const
 Return raw pointer to NumeredEntFiniteElement. More...
 
EntityHandle getFEEntityHandle () const
 Return finite element entity handle. More...
 
int getFEDim () const
 Get dimension of finite element. More...
 
EntityType getFEType () const
 Get dimension of finite element. More...
 
boost::weak_ptr< SideNumbergetSideNumberPtr (const int side_number, const EntityType type)
 Get the side number pointer. More...
 
EntityHandle getSideEntity (const int side_number, const EntityType type)
 Get the side entity. More...
 
int getNumberOfNodesOnElement () const
 Get the number of nodes on finite element. More...
 
MoFEMErrorCode getProblemRowIndices (const std::string filed_name, const EntityType type, const int side, VectorInt &indices) const
 Get row indices. More...
 
MoFEMErrorCode getProblemColIndices (const std::string filed_name, const EntityType type, const int side, VectorInt &indices) const
 Get col indices. More...
 
const FEMethodgetFEMethod () const
 Return raw pointer to Finite Element Method object. More...
 
int getOpType () const
 Get operator types. More...
 
void setOpType (const OpType type)
 Set operator type. More...
 
void addOpType (const OpType type)
 Add operator type. More...
 
int getNinTheLoop () const
 get number of finite element in the loop More...
 
int getLoopSize () const
 get size of elements in the loop More...
 
std::string getFEName () const
 Get name of the element. More...
 
ForcesAndSourcesCoregetPtrFE () const
 
ForcesAndSourcesCoregetSidePtrFE () const
 
ForcesAndSourcesCoregetRefinePtrFE () const
 
const PetscData::SwitchesgetDataCtx () const
 
const KspMethod::KSPContext getKSPCtx () const
 
const SnesMethod::SNESContext getSNESCtx () const
 
const TSMethod::TSContext getTSCtx () const
 
Vec getKSPf () const
 
Mat getKSPA () const
 
Mat getKSPB () const
 
Vec getSNESf () const
 
Vec getSNESx () const
 
Mat getSNESA () const
 
Mat getSNESB () const
 
Vec getTSu () const
 
Vec getTSu_t () const
 
Vec getTSu_tt () const
 
Vec getTSf () const
 
Mat getTSA () const
 
Mat getTSB () const
 
int getTSstep () const
 
double getTStime () const
 
double getTStimeStep () const
 
double getTSa () const
 
double getTSaa () const
 
MatrixDoublegetGaussPts ()
 matrix of integration (Gauss) points for Volume Element More...
 
auto getFTensor0IntegrationWeight ()
 Get integration weights. More...
 
MatrixDoublegetCoordsAtGaussPts ()
 Gauss points and weight, matrix (nb. of points x 3) More...
 
auto getFTensor1CoordsAtGaussPts ()
 Get coordinates at integration points assuming linear geometry. More...
 
double getMeasure () const
 get measure of element More...
 
doublegetMeasure ()
 get measure of element More...
 
MoFEMErrorCode loopSide (const string &fe_name, ForcesAndSourcesCore *side_fe, const size_t dim, const EntityHandle ent_for_side=0, const int verb=QUIET, const LogManager::SeverityLevel sev=Sev::noisy, AdjCache *adj_cache=nullptr)
 User call this function to loop over elements on the side of face. This function calls finite element with is operator to do calculations. More...
 
MoFEMErrorCode loopThis (const string &fe_name, ForcesAndSourcesCore *parent_fe, const int verb=QUIET, const LogManager::SeverityLevel sev=Sev::noisy)
 User call this function to loop over parent elements. This function calls finite element with is operator to do calculations. More...
 
MoFEMErrorCode loopParent (const string &fe_name, ForcesAndSourcesCore *parent_fe, const int verb=QUIET, const LogManager::SeverityLevel sev=Sev::noisy)
 User call this function to loop over parent elements. This function calls finite element with is operator to do calculations. More...
 
MoFEMErrorCode loopChildren (const string &fe_name, ForcesAndSourcesCore *child_fe, const int verb=QUIET, const LogManager::SeverityLevel sev=Sev::noisy)
 User call this function to loop over parent elements. This function calls finite element with is operator to do calculations. More...
 
- Public Member Functions inherited from MoFEM::DataOperator
 DataOperator (const bool symm=true)
 
virtual ~DataOperator ()=default
 
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. More...
 
virtual MoFEMErrorCode opLhs (EntitiesFieldData &row_data, EntitiesFieldData &col_data)
 
virtual MoFEMErrorCode doWork (int side, EntityType type, EntitiesFieldData::EntData &data)
 Operator for linear form, usually to calculate values on right hand side. More...
 
virtual MoFEMErrorCode opRhs (EntitiesFieldData &data, const bool error_if_no_base=false)
 
bool getSymm () const
 Get if operator uses symmetry of DOFs or not. More...
 
void setSymm ()
 set if operator is executed taking in account symmetry More...
 
void unSetSymm ()
 unset if operator is executed for non symmetric problem More...
 

Private Attributes

MatrixDouble nN
 
VectorDouble nF
 
Rangeessential_bd_ents
 

Additional Inherited Members

- Public Types inherited from MoFEM::ForcesAndSourcesCore::UserDataOperator
enum  OpType {
  OPROW = 1 << 0 , OPCOL = 1 << 1 , OPROWCOL = 1 << 2 , OPSPACE = 1 << 3 ,
  OPLAST = 1 << 3
}
 Controls loop over entities on element. More...
 
using AdjCache = std::map< EntityHandle, std::vector< boost::weak_ptr< NumeredEntFiniteElement > > >
 
- Public Types inherited from MoFEM::DataOperator
using DoWorkLhsHookFunType = boost::function< MoFEMErrorCode(DataOperator *op_ptr, int row_side, int col_side, EntityType row_type, EntityType col_type, EntitiesFieldData::EntData &row_data, EntitiesFieldData::EntData &col_data)>
 
using DoWorkRhsHookFunType = boost::function< MoFEMErrorCode(DataOperator *op_ptr, int side, EntityType type, EntitiesFieldData::EntData &data)>
 
- Public Attributes inherited from MoFEM::ForcesAndSourcesCore::UserDataOperator
char opType
 
std::string rowFieldName
 
std::string colFieldName
 
FieldSpace sPace
 
- Public Attributes inherited from MoFEM::DataOperator
DoWorkLhsHookFunType doWorkLhsHook
 
DoWorkRhsHookFunType doWorkRhsHook
 
bool sYmm
 If true assume that matrix is symmetric structure. More...
 
std::array< bool, MBMAXTYPE > doEntities
 If true operator is executed for entity. More...
 
booldoVertices
 \deprectaed If false skip vertices More...
 
booldoEdges
 \deprectaed If false skip edges More...
 
booldoQuads
 \deprectaed More...
 
booldoTris
 \deprectaed More...
 
booldoTets
 \deprectaed More...
 
booldoPrisms
 \deprectaed More...
 
- Static Public Attributes inherited from MoFEM::ForcesAndSourcesCore::UserDataOperator
static const char *const OpTypeNames []
 
virtual MoFEMErrorCode setPtrFE (ForcesAndSourcesCore *ptr)
 
- Protected Attributes inherited from MoFEM::ForcesAndSourcesCore::UserDataOperator
ForcesAndSourcesCoreptrFE
 

Detailed Description

Definition at line 58 of file ElecPhysOperators.hpp.

Constructor & Destructor Documentation

◆ OpEssentialBC() [1/2]

ElectroPhysiology::OpEssentialBC::OpEssentialBC ( std::string  flux_field,
Range essential_bd_ents 
)
inline

Definition at line 59 of file ElecPhysOperators.hpp.

61 }
FaceEle::UserDataOperator OpFaceEle
@ OPROW
operator doWork function is executed on FE rows

◆ OpEssentialBC() [2/2]

ElectroPhysiology::OpEssentialBC::OpEssentialBC ( const std::string &  flux_field,
Range essential_bd_ents 
)
inline

Definition at line 72 of file ElecPhysOperators2D.hpp.

BoundaryEle::UserDataOperator OpBoundaryEle

Member Function Documentation

◆ doWork() [1/2]

MoFEMErrorCode ElectroPhysiology::OpEssentialBC::doWork ( int  side,
EntityType  type,
EntData data 
)
inline

Definition at line 62 of file ElecPhysOperators.hpp.

62 {
64 int nb_dofs = data.getIndices().size();
65
66 if (nb_dofs) {
68 bool is_essential =
69 (essential_bd_ents.find(fe_ent) != essential_bd_ents.end());
70 if (is_essential) {
71 int nb_gauss_pts = getGaussPts().size2();
72 int size2 = data.getN().size2();
73 if (3 * nb_dofs != static_cast<int>(data.getN().size2()))
74 SETERRQ(PETSC_COMM_WORLD, MOFEM_DATA_INCONSISTENCY,
75 "wrong number of dofs");
76 nN.resize(nb_dofs, nb_dofs, false);
77 nF.resize(nb_dofs, false);
78 nN.clear();
79 nF.clear();
80
81 auto t_row_tau = data.getFTensor1N<3>();
82
83 double *normal_ptr;
84 if (getNormalsAtGaussPts().size1() == (unsigned int)nb_gauss_pts) {
85 // HO geometry
86 normal_ptr = &getNormalsAtGaussPts(0)[0];
87 } else {
88 // Linear geometry, i.e. constant normal on face
89 normal_ptr = &getNormal()[0];
90 }
91 // set tensor from pointer
92 FTensor::Tensor1<const double *, 3> t_normal(normal_ptr, &normal_ptr[1],
93 &normal_ptr[2], 3);
94
96 const double vol = getMeasure();
97 double nrm2 = 0;
98 for (int gg = 0; gg < nb_gauss_pts; gg++) {
99 if (gg == 0) {
100 nrm2 = sqrt(t_normal(i) * t_normal(i));
101 }
102 const double a = t_w * vol;
103 for (int rr = 0; rr != nb_dofs; rr++) {
105 &data.getVectorN<3>(gg)(0, HVEC0),
106 &data.getVectorN<3>(gg)(0, HVEC1),
107 &data.getVectorN<3>(gg)(0, HVEC2), 3);
108 nF[rr] += a * essen_value * t_row_tau(i) * t_normal(i) / nrm2;
109 for (int cc = 0; cc != nb_dofs; cc++) {
110 nN(rr, cc) += a * (t_row_tau(i) * t_normal(i)) *
111 (t_col_tau(i) * t_normal(i)) / (nrm2 * nrm2);
112 ++t_col_tau;
113 }
114 ++t_row_tau;
115 }
116 // If HO geometry increment t_normal to next integration point
117 if (getNormalsAtGaussPts().size1() == (unsigned int)nb_gauss_pts) {
118 ++t_normal;
119 nrm2 = sqrt(t_normal(i) * t_normal(i));
120 }
121 ++t_w;
122 }
123
125 cholesky_solve(nN, nF, ublas::lower());
126
127 for (auto &dof : data.getFieldDofs()) {
128 dof->getFieldData() = nF[dof->getEntDofIdx()];
129 }
130 }
131 }
133 }
constexpr double a
void cholesky_solve(const TRIA &L, VEC &x, ublas::lower)
solve system L L^T x = b inplace
Definition: cholesky.hpp:221
size_t cholesky_decompose(const MATRIX &A, TRIA &L)
decompose the symmetric positive definit matrix A into product L L^T.
Definition: cholesky.hpp:52
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:346
@ HVEC0
Definition: definitions.h:186
@ HVEC1
Definition: definitions.h:186
@ HVEC2
Definition: definitions.h:186
@ MOFEM_DATA_INCONSISTENCY
Definition: definitions.h:31
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:416
const double essen_value
FTensor::Index< 'i', 3 > i
MatrixDouble & getN(const FieldApproximationBase base)
get base functions this return matrix (nb. of rows is equal to nb. of Gauss pts, nb....
const MatrixAdaptor getVectorN(const FieldApproximationBase base, const int gg)
get Hdiv of base functions at Gauss pts
FTensor::Tensor1< FTensor::PackPtr< double *, Tensor_Dim >, Tensor_Dim > getFTensor1N(FieldApproximationBase base)
Get base functions for Hdiv/Hcurl spaces.
const VectorDofs & getFieldDofs() const
get dofs data stature FEDofEntity
const VectorInt & getIndices() const
Get global indices of dofs on entity.
MatrixDouble & getNormalsAtGaussPts()
if higher order geometry return normals at Gauss pts.
EntityHandle getFEEntityHandle() const
Return finite element entity handle.
auto getFTensor0IntegrationWeight()
Get integration weights.
double getMeasure() const
get measure of element
MatrixDouble & getGaussPts()
matrix of integration (Gauss) points for Volume Element

◆ doWork() [2/2]

MoFEMErrorCode ElectroPhysiology::OpEssentialBC::doWork ( int  side,
EntityType  type,
EntData data 
)
inline

Definition at line 76 of file ElecPhysOperators2D.hpp.

76 {
78 int nb_dofs = data.getIndices().size();
79 if (nb_dofs) {
81 bool is_essential =
82 (essential_bd_ents.find(fe_ent) != essential_bd_ents.end());
83 if (is_essential) {
84 int nb_gauss_pts = getGaussPts().size2();
85 int size2 = data.getN().size2();
86 if (3 * nb_dofs != static_cast<int>(data.getN().size2()))
87 SETERRQ(PETSC_COMM_WORLD, MOFEM_DATA_INCONSISTENCY,
88 "wrong number of dofs");
89 nN.resize(nb_dofs, nb_dofs, false);
90 nF.resize(nb_dofs, false);
91 nN.clear();
92 nF.clear();
93
94 auto t_row_tau = data.getFTensor1N<3>();
95
96 auto dir = getDirection();
97 double len = sqrt(dir[0] * dir[0] + dir[1] * dir[1] + dir[2] * dir[2]);
98
99 FTensor::Tensor1<double, 3> t_normal(-dir[1] / len, dir[0] / len,
100 dir[2] / len);
101
102 auto t_w = getFTensor0IntegrationWeight();
103 const double vol = getMeasure();
104
105 for (int gg = 0; gg < nb_gauss_pts; gg++) {
106 const double a = t_w * vol;
107 for (int rr = 0; rr != nb_dofs; rr++) {
108 auto t_col_tau = data.getFTensor1N<3>(gg, 0);
109 nF[rr] += a * essen_value * t_row_tau(i) * t_normal(i);
110 for (int cc = 0; cc != nb_dofs; cc++) {
111 nN(rr, cc) += a * (t_row_tau(i) * t_normal(i)) *
112 (t_col_tau(i) * t_normal(i));
113 ++t_col_tau;
114 }
115 ++t_row_tau;
116 }
117 ++t_w;
118 }
119
121 cholesky_solve(nN, nF, ublas::lower());
122
123 for (auto &dof : data.getFieldDofs()) {
124 dof->getFieldData() = nF[dof->getEntDofIdx()];
125 }
126 }
127 }
129 }

Member Data Documentation

◆ essential_bd_ents

Range & ElectroPhysiology::OpEssentialBC::essential_bd_ents
private

Definition at line 138 of file ElecPhysOperators.hpp.

◆ nF

VectorDouble ElectroPhysiology::OpEssentialBC::nF
private

Definition at line 137 of file ElecPhysOperators.hpp.

◆ nN

MatrixDouble ElectroPhysiology::OpEssentialBC::nN
private

Definition at line 136 of file ElecPhysOperators.hpp.


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