v0.14.0
Loading...
Searching...
No Matches
Public Member Functions | Private Attributes | List of all members
ElectroPhysiology::OpAssembleLhsTauV< dim > Struct Template Reference

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

Inheritance diagram for ElectroPhysiology::OpAssembleLhsTauV< dim >:
[legend]
Collaboration diagram for ElectroPhysiology::OpAssembleLhsTauV< dim >:
[legend]

Public Member Functions

 OpAssembleLhsTauV (std::string flux_field, std::string mass_field)
 
MoFEMErrorCode doWork (int row_side, int col_side, EntityType row_type, EntityType col_type, EntData &row_data, EntData &col_data)
 
 OpAssembleLhsTauV (std::string flux_field, std::string mass_field, boost::shared_ptr< PreviousData > &data, std::map< int, BlockData > &block_map)
 
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 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...
 
- 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
 
boost::shared_ptr< PreviousDatacommonData
 
std::map< int, BlockDatasetOfBlock
 

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

template<int dim>
struct ElectroPhysiology::OpAssembleLhsTauV< dim >

Definition at line 604 of file ElecPhysOperators.hpp.

Constructor & Destructor Documentation

◆ OpAssembleLhsTauV() [1/2]

template<int dim>
ElectroPhysiology::OpAssembleLhsTauV< dim >::OpAssembleLhsTauV ( std::string  flux_field,
std::string  mass_field 
)
inline

Definition at line 606 of file ElecPhysOperators.hpp.

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

◆ OpAssembleLhsTauV() [2/2]

template<int dim>
ElectroPhysiology::OpAssembleLhsTauV< dim >::OpAssembleLhsTauV ( std::string  flux_field,
std::string  mass_field,
boost::shared_ptr< PreviousData > &  data,
std::map< int, BlockData > &  block_map 
)
inline

Definition at line 675 of file ElecPhysOperators2D.hpp.

678 : OpFaceEle(flux_field, mass_field, OpFaceEle::OPROWCOL),
679 commonData(data), setOfBlock(block_map) {
680 sYmm = false;
681 }
FaceEle::UserDataOperator OpFaceEle
boost::shared_ptr< PreviousData > commonData

Member Function Documentation

◆ doWork() [1/2]

template<int dim>
MoFEMErrorCode ElectroPhysiology::OpAssembleLhsTauV< dim >::doWork ( int  row_side,
int  col_side,
EntityType  row_type,
EntityType  col_type,
EntData row_data,
EntData col_data 
)
inline

Definition at line 611 of file ElecPhysOperators.hpp.

613 {
615 const int nb_row_dofs = row_data.getIndices().size();
616 const int nb_col_dofs = col_data.getIndices().size();
617
618 if (nb_row_dofs && nb_col_dofs) {
619
620 mat.resize(nb_row_dofs, nb_col_dofs, false);
621 mat.clear();
622 const int nb_integration_pts = getGaussPts().size2();
623 auto t_w = getFTensor0IntegrationWeight();
624 auto t_row_tau_base = row_data.getFTensor1N<3>();
625
626 auto t_row_tau_grad = row_data.getFTensor2DiffN<3, 3>();
627 const double vol = getMeasure();
628
629 for (int gg = 0; gg != nb_integration_pts; ++gg) {
630 const double a = vol * t_w;
631
632 for (int rr = 0; rr != nb_row_dofs; ++rr) {
633 auto t_col_v_base = col_data.getFTensor0N(gg, 0);
634 for (int cc = 0; cc != nb_col_dofs; ++cc) {
635 double div_row_base = t_row_tau_grad(0, 0) + t_row_tau_grad(1, 1) +
636 t_row_tau_grad(2, 2);
637 mat(rr, cc) += -(div_row_base * t_col_v_base) * a;
638 ++t_col_v_base;
639 }
640 ++t_row_tau_base;
641 ++t_row_tau_grad;
642 }
643 ++t_w;
644 }
645 CHKERR MatSetValues(getFEMethod()->ts_B, row_data, col_data, &mat(0, 0),
646 ADD_VALUES);
647 }
649 }
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.
FTensor::Tensor1< FTensor::PackPtr< double *, Tensor_Dim >, Tensor_Dim > getFTensor1N(FieldApproximationBase base)
Get base functions for Hdiv/Hcurl spaces.
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]

template<int dim>
MoFEMErrorCode ElectroPhysiology::OpAssembleLhsTauV< dim >::doWork ( int  row_side,
int  col_side,
EntityType  row_type,
EntityType  col_type,
EntData row_data,
EntData col_data 
)
inline

Definition at line 683 of file ElecPhysOperators2D.hpp.

685 {
687 const int nb_row_dofs = row_data.getIndices().size();
688 const int nb_col_dofs = col_data.getIndices().size();
689
690 if (nb_row_dofs && nb_col_dofs) {
691 // auto find_block_data = [&]() {
692 // EntityHandle fe_ent = getFEEntityHandle();
693 // BlockData *block_raw_ptr = nullptr;
694 // for (auto &m : setOfBlock) {
695 // if (m.second.block_ents.find(fe_ent) != m.second.block_ents.end()) {
696 // block_raw_ptr = &m.second;
697 // break;
698 // }
699 // }
700 // return block_raw_ptr;
701 // };
702
703 // auto block_data_ptr = find_block_data();
704 // if (!block_data_ptr)
705 // SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY, "Block not found");
706 // auto &block_data = *block_data_ptr;
707 mat.resize(nb_row_dofs, nb_col_dofs, false);
708 mat.clear();
709 const int nb_integration_pts = getGaussPts().size2();
710 auto t_w = getFTensor0IntegrationWeight();
711 auto t_row_tau_base = row_data.getFTensor1N<3>();
712
713 auto t_row_tau_grad = row_data.getFTensor2DiffN<3, 2>();
714 auto t_mass_value = getFTensor0FromVec(commonData->mass_values);
715 auto t_flux_value = getFTensor1FromMat<3>(commonData->flux_values);
716 const double vol = getMeasure();
717
718 for (int gg = 0; gg != nb_integration_pts; ++gg) {
719 const double a = vol * t_w;
720 const double K = B_epsilon + (block.B0 + B * t_mass_value);
721 const double K_inv = 1. / K;
722 const double K_diff = B;
723
724 for (int rr = 0; rr != nb_row_dofs; ++rr) {
725 auto t_col_v_base = col_data.getFTensor0N(gg, 0);
726 for (int cc = 0; cc != nb_col_dofs; ++cc) {
727 double div_row_base = t_row_tau_grad(0, 0) + t_row_tau_grad(1, 1);
728 mat(rr, cc) += (-(t_row_tau_base(i) * t_flux_value(i) * K_inv *
729 K_inv * K_diff * t_col_v_base) -
730 (div_row_base * t_col_v_base)) *
731 a;
732 ++t_col_v_base;
733 }
734 ++t_row_tau_base;
735 ++t_row_tau_grad;
736 }
737 ++t_w;
738 ++t_mass_value;
739 ++t_flux_value;
740 }
741 CHKERR MatSetValues(getFEMethod()->ts_B, row_data, col_data, &mat(0, 0),
742 ADD_VALUES);
743 }
745 }
FTensor::Index< 'i', 3 > i
static auto getFTensor0FromVec(ublas::vector< T, A > &data)
Get tensor rank 0 (scalar) form data vector.
Definition: Templates.hpp:135

Member Data Documentation

◆ commonData

template<int dim>
boost::shared_ptr<PreviousData> ElectroPhysiology::OpAssembleLhsTauV< dim >::commonData
private

Definition at line 748 of file ElecPhysOperators2D.hpp.

◆ mat

template<int dim>
MatrixDouble ElectroPhysiology::OpAssembleLhsTauV< dim >::mat
private

Definition at line 652 of file ElecPhysOperators.hpp.

◆ setOfBlock

template<int dim>
std::map<int, BlockData> ElectroPhysiology::OpAssembleLhsTauV< dim >::setOfBlock
private

Definition at line 750 of file ElecPhysOperators2D.hpp.


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