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

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

Inheritance diagram for ReactionDiffusion::OpAssembleLhsTauTau< dim >:
[legend]
Collaboration diagram for ReactionDiffusion::OpAssembleLhsTauTau< dim >:
[legend]

Public Member Functions

 OpAssembleLhsTauTau (std::string flux_field, boost::shared_ptr< PreviousData > &commonData, 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::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 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...
 

Private Attributes

boost::shared_ptr< PreviousDatacommonData
 
MatrixDouble mat
 
MatrixDouble transMat
 
Range essential_bd_ents
 
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 []
 
virtual MoFEMErrorCode setPtrFE (ForcesAndSourcesCore *ptr)
 
- Protected Attributes inherited from MoFEM::ForcesAndSourcesCore::UserDataOperator
ForcesAndSourcesCoreptrFE
 

Detailed Description

template<int dim>
struct ReactionDiffusion::OpAssembleLhsTauTau< dim >

Definition at line 689 of file RDOperators.hpp.

Constructor & Destructor Documentation

◆ OpAssembleLhsTauTau()

template<int dim>
ReactionDiffusion::OpAssembleLhsTauTau< dim >::OpAssembleLhsTauTau ( std::string  flux_field,
boost::shared_ptr< PreviousData > &  commonData,
std::map< int, BlockData > &  block_map 
)
inline

Definition at line 691 of file RDOperators.hpp.

694 : OpFaceEle(flux_field, flux_field, OpFaceEle::OPROWCOL),
695 setOfBlock(block_map), commonData(commonData)
696 {
697 sYmm = true;
698 }
FaceEle::UserDataOperator OpFaceEle
Definition: RDOperators.hpp:13
bool sYmm
If true assume that matrix is symmetric structure.
@ OPROWCOL
operator doWork is executed on FE rows &columns
std::map< int, BlockData > setOfBlock
boost::shared_ptr< PreviousData > commonData

Member Function Documentation

◆ doWork()

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

Definition at line 700 of file RDOperators.hpp.

702 {
704 const int nb_row_dofs = row_data.getIndices().size();
705 const int nb_col_dofs = col_data.getIndices().size();
706
707 if (nb_row_dofs && nb_col_dofs) {
708 auto find_block_data = [&]() {
710 BlockData *block_raw_ptr = nullptr;
711 for (auto &m : setOfBlock) {
712 if (m.second.block_ents.find(fe_ent) != m.second.block_ents.end()) {
713 block_raw_ptr = &m.second;
714 break;
715 }
716 }
717 return block_raw_ptr;
718 };
719
720 auto block_data_ptr = find_block_data();
721 if (!block_data_ptr)
722 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY, "Block not found");
723 auto &block_data = *block_data_ptr;
724
725 mat.resize(nb_row_dofs, nb_col_dofs, false);
726 mat.clear();
727 const int nb_integration_pts = getGaussPts().size2();
728 auto t_mass_value = getFTensor0FromVec(commonData->mass_values);
729
730 auto t_row_tau_base = row_data.getFTensor1N<3>();
731
732 auto t_w = getFTensor0IntegrationWeight();
733 const double vol = getMeasure();
734
735 for (int gg = 0; gg != nb_integration_pts; ++gg) {
736 const double a = vol * t_w;
737 const double K = B_epsilon + (block_data.B0 + B * t_mass_value);
738 const double K_inv = 1. / K;
739 for (int rr = 0; rr != nb_row_dofs; ++rr) {
740 auto t_col_tau_base = col_data.getFTensor1N<3>(gg, 0);
741 for (int cc = 0; cc != nb_col_dofs; ++cc) {
742 mat(rr, cc) += (K_inv * t_row_tau_base(i) * t_col_tau_base(i)) * a;
743 ++t_col_tau_base;
744 }
745 ++t_row_tau_base;
746 }
747 ++t_mass_value;
748 ++t_w;
749 }
750 CHKERR MatSetValues(getFEMethod()->ts_B, row_data, col_data, &mat(0, 0),
751 ADD_VALUES);
752 if (row_side != col_side || row_type != col_type) {
753 transMat.resize(nb_col_dofs, nb_row_dofs, false);
754 noalias(transMat) = trans(mat);
755 CHKERR MatSetValues(getFEMethod()->ts_B, col_data, row_data,
756 &transMat(0, 0), ADD_VALUES);
757 }
758 }
760 }
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
@ 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
#define CHKERR
Inline error check.
Definition: definitions.h:535
FTensor::Index< 'm', SPACE_DIM > m
MoFEMErrorCode MatSetValues(Mat M, const EntitiesFieldData::EntData &row_data, const EntitiesFieldData::EntData &col_data, const double *ptr, InsertMode iora)
Assemble PETSc matrix.
static auto getFTensor0FromVec(ublas::vector< T, A > &data)
Get tensor rank 0 (scalar) form data vector.
Definition: Templates.hpp:135
FTensor::Index< 'i', 3 > i
Definition: RDOperators.hpp:29
const double B_epsilon
Definition: RDOperators.hpp:21
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.
EntityHandle getFEEntityHandle() const
Return finite element entity handle.
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

Member Data Documentation

◆ commonData

template<int dim>
boost::shared_ptr<PreviousData> ReactionDiffusion::OpAssembleLhsTauTau< dim >::commonData
private

Definition at line 763 of file RDOperators.hpp.

◆ essential_bd_ents

template<int dim>
Range ReactionDiffusion::OpAssembleLhsTauTau< dim >::essential_bd_ents
private

Definition at line 765 of file RDOperators.hpp.

◆ mat

template<int dim>
MatrixDouble ReactionDiffusion::OpAssembleLhsTauTau< dim >::mat
private

Definition at line 764 of file RDOperators.hpp.

◆ setOfBlock

template<int dim>
std::map<int, BlockData> ReactionDiffusion::OpAssembleLhsTauTau< dim >::setOfBlock
private

Definition at line 766 of file RDOperators.hpp.

◆ transMat

template<int dim>
MatrixDouble ReactionDiffusion::OpAssembleLhsTauTau< dim >::transMat
private

Definition at line 764 of file RDOperators.hpp.


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