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

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

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

Public Member Functions

 OpAssembleLhsVTau (std::string mass_field, std::string flux_field)
 
MoFEMErrorCode doWork (int row_side, int col_side, EntityType row_type, EntityType col_type, EntData &row_data, EntData &col_data)
 
 OpAssembleLhsVTau (std::string mass_field, std::string flux_field)
 
MoFEMErrorCode doWork (int row_side, int col_side, EntityType row_type, EntityType col_type, EntData &row_data, EntData &col_data)
 
- Public Member Functions inherited from MoFEM::VolumeElementForcesAndSourcesCore::UserDataOperator
int getNumNodes ()
 get element number of nodes More...
 
const EntityHandlegetConn ()
 get element connectivity More...
 
double getVolume () const
 element volume (linear geometry) More...
 
doublegetVolume ()
 element volume (linear geometry) More...
 
FTensor::Tensor2< double *, 3, 3 > & getJac ()
 get element Jacobian More...
 
FTensor::Tensor2< double *, 3, 3 > & getInvJac ()
 get element inverse Jacobian More...
 
VectorDoublegetCoords ()
 nodal coordinates More...
 
VolumeElementForcesAndSourcesCoregetVolumeFE () const
 return pointer to Generic Volume Finite Element object More...
 
- 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 calls this function to loop over elements on the side of face. This function calls finite element with its operator to do calculations. More...
 
MoFEMErrorCode loopThis (const string &fe_name, ForcesAndSourcesCore *this_fe, const int verb=QUIET, const LogManager::SeverityLevel sev=Sev::noisy)
 User calls this function to loop over the same element using a different set of integration points. This function calls finite element with its 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 calls this function to loop over parent elements. This function calls finite element with its 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 calls this function to loop over parent elements. This function calls finite element with its 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...
 
- 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)
 

Private Attributes

MatrixDouble mat
 

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 []
 
- Protected Member Functions inherited from MoFEM::VolumeElementForcesAndSourcesCore::UserDataOperator
MoFEMErrorCode setPtrFE (ForcesAndSourcesCore *ptr)
 
virtual MoFEMErrorCode setPtrFE (ForcesAndSourcesCore *ptr)
 
- Protected Attributes inherited from MoFEM::ForcesAndSourcesCore::UserDataOperator
ForcesAndSourcesCoreptrFE
 

Detailed Description

Definition at line 656 of file ElecPhysOperators.hpp.

Constructor & Destructor Documentation

◆ OpAssembleLhsVTau() [1/2]

ElectroPhysiology::OpAssembleLhsVTau::OpAssembleLhsVTau ( std::string  mass_field,
std::string  flux_field 
)
inline

Definition at line 658 of file ElecPhysOperators.hpp.

659 : OpVolEle(mass_field, flux_field, OpVolEle::OPROWCOL) {
660 sYmm = false;
661 }
VolEle::UserDataOperator OpVolEle
bool sYmm
If true assume that matrix is symmetric structure.
@ OPROWCOL
operator doWork is executed on FE rows &columns

◆ OpAssembleLhsVTau() [2/2]

ElectroPhysiology::OpAssembleLhsVTau::OpAssembleLhsVTau ( std::string  mass_field,
std::string  flux_field 
)
inline

Definition at line 756 of file ElecPhysOperators2D.hpp.

757 : OpFaceEle(mass_field, flux_field, OpFaceEle::OPROWCOL) {
758 sYmm = false;
759 }
FaceEle::UserDataOperator OpFaceEle

Member Function Documentation

◆ doWork() [1/2]

MoFEMErrorCode ElectroPhysiology::OpAssembleLhsVTau::doWork ( int  row_side,
int  col_side,
EntityType  row_type,
EntityType  col_type,
EntData row_data,
EntData col_data 
)
inline

Definition at line 663 of file ElecPhysOperators.hpp.

665 {
667 const int nb_row_dofs = row_data.getIndices().size();
668 const int nb_col_dofs = col_data.getIndices().size();
669
670 if (nb_row_dofs && nb_col_dofs) {
671 mat.resize(nb_row_dofs, nb_col_dofs, false);
672 mat.clear();
673 const int nb_integration_pts = getGaussPts().size2();
674 auto t_w = getFTensor0IntegrationWeight();
675 auto t_row_v_base = row_data.getFTensor0N();
676 const double vol = getMeasure();
677
678 for (int gg = 0; gg != nb_integration_pts; ++gg) {
679 const double a = vol * t_w;
680 for (int rr = 0; rr != nb_row_dofs; ++rr) {
681 auto t_col_tau_grad = col_data.getFTensor2DiffN<3, 3>(gg, 0);
682 for (int cc = 0; cc != nb_col_dofs; ++cc) {
683 double div_col_base = t_col_tau_grad(0, 0) + t_col_tau_grad(1, 1) + t_col_tau_grad(2, 2);
684 mat(rr, cc) += (t_row_v_base * div_col_base) * a;
685 ++t_col_tau_grad;
686 }
687 ++t_row_v_base;
688 }
689 ++t_w;
690 }
691 CHKERR MatSetValues(getFEMethod()->ts_B, row_data, col_data, &mat(0, 0),
692 ADD_VALUES);
693 }
695 }
constexpr double a
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:346
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:416
#define CHKERR
Inline error check.
Definition: definitions.h:535
MoFEMErrorCode MatSetValues(Mat M, const EntitiesFieldData::EntData &row_data, const EntitiesFieldData::EntData &col_data, const double *ptr, InsertMode iora)
Assemble PETSc matrix.
FTensor::Tensor2< FTensor::PackPtr< double *, Tensor_Dim0 *Tensor_Dim1 >, Tensor_Dim0, Tensor_Dim1 > getFTensor2DiffN(FieldApproximationBase base)
Get derivatives of base functions for Hdiv space.
FTensor::Tensor0< FTensor::PackPtr< double *, 1 > > getFTensor0N(const FieldApproximationBase base)
Get base function as Tensor0.
const VectorInt & getIndices() const
Get global indices of dofs on entity.
auto getFTensor0IntegrationWeight()
Get integration weights.
double getMeasure() const
get measure of element
const FEMethod * getFEMethod() const
Return raw pointer to Finite Element Method object.
MatrixDouble & getGaussPts()
matrix of integration (Gauss) points for Volume Element

◆ doWork() [2/2]

MoFEMErrorCode ElectroPhysiology::OpAssembleLhsVTau::doWork ( int  row_side,
int  col_side,
EntityType  row_type,
EntityType  col_type,
EntData row_data,
EntData col_data 
)
inline

Definition at line 761 of file ElecPhysOperators2D.hpp.

763 {
765 const int nb_row_dofs = row_data.getIndices().size();
766 const int nb_col_dofs = col_data.getIndices().size();
767
768 if (nb_row_dofs && nb_col_dofs) {
769 mat.resize(nb_row_dofs, nb_col_dofs, false);
770 mat.clear();
771 const int nb_integration_pts = getGaussPts().size2();
772 auto t_w = getFTensor0IntegrationWeight();
773 auto t_row_v_base = row_data.getFTensor0N();
774 const double vol = getMeasure();
775
776 for (int gg = 0; gg != nb_integration_pts; ++gg) {
777 const double a = vol * t_w;
778 for (int rr = 0; rr != nb_row_dofs; ++rr) {
779 auto t_col_tau_grad = col_data.getFTensor2DiffN<3, 2>(gg, 0);
780 for (int cc = 0; cc != nb_col_dofs; ++cc) {
781 double div_col_base = t_col_tau_grad(0, 0) + t_col_tau_grad(1, 1);
782 mat(rr, cc) += (t_row_v_base * div_col_base) * a;
783 ++t_col_tau_grad;
784 }
785 ++t_row_v_base;
786 }
787 ++t_w;
788 }
789 CHKERR MatSetValues(getFEMethod()->ts_B, row_data, col_data, &mat(0, 0),
790 ADD_VALUES);
791 }
793 }

Member Data Documentation

◆ mat

MatrixDouble ElectroPhysiology::OpAssembleLhsVTau::mat
private

Definition at line 698 of file ElecPhysOperators.hpp.


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