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

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

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

Public Types

typedef boost::function< double(const double, const double)> FUval
 
- 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 Member Functions

 OpAssembleLhsVV (std::string mass_field)
 
MoFEMErrorCode doWork (int row_side, int col_side, EntityType row_type, EntityType col_type, EntData &row_data, EntData &col_data)
 
 OpAssembleLhsVV (std::string mass_field, boost::shared_ptr< PreviousData > &datau, boost::shared_ptr< PreviousData > &datav, FUval Drhs_u)
 
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)
 

Public Attributes

boost::shared_ptr< PreviousData > & dataU
 
boost::shared_ptr< PreviousData > & dataV
 
FUval DRhs_u
 
- 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...
 

Private Attributes

MatrixDouble mat
 
MatrixDouble transMat
 

Additional Inherited Members

- 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 702 of file ElecPhysOperators.hpp.

Member Typedef Documentation

◆ FUval

typedef boost::function<double(const double, const double)> ElectroPhysiology::OpAssembleLhsVV::FUval

Definition at line 802 of file ElecPhysOperators2D.hpp.

Constructor & Destructor Documentation

◆ OpAssembleLhsVV() [1/2]

ElectroPhysiology::OpAssembleLhsVV::OpAssembleLhsVV ( std::string  mass_field)
inline

Definition at line 704 of file ElecPhysOperators.hpp.

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

◆ OpAssembleLhsVV() [2/2]

ElectroPhysiology::OpAssembleLhsVV::OpAssembleLhsVV ( std::string  mass_field,
boost::shared_ptr< PreviousData > &  datau,
boost::shared_ptr< PreviousData > &  datav,
FUval  Drhs_u 
)
inline

Definition at line 803 of file ElecPhysOperators2D.hpp.

807 : OpFaceEle(mass_field, mass_field, OpFaceEle::OPROWCOL)
808 , DRhs_u(Drhs_u)
809 , dataU(datau)
810 , dataV(datav){
811 sYmm = true;
812 }
FaceEle::UserDataOperator OpFaceEle
boost::shared_ptr< PreviousData > & dataV
boost::shared_ptr< PreviousData > & dataU

Member Function Documentation

◆ doWork() [1/2]

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

Definition at line 709 of file ElecPhysOperators.hpp.

711 {
713
714 const int nb_row_dofs = row_data.getIndices().size();
715 const int nb_col_dofs = col_data.getIndices().size();
716 if (nb_row_dofs && nb_col_dofs) {
717
718 mat.resize(nb_row_dofs, nb_col_dofs, false);
719 mat.clear();
720 const int nb_integration_pts = getGaussPts().size2();
721
722 auto t_row_v_base = row_data.getFTensor0N();
723
724 auto t_w = getFTensor0IntegrationWeight();
725 const double ts_a = getFEMethod()->ts_a;
726 const double vol = getMeasure();
727
728 for (int gg = 0; gg != nb_integration_pts; ++gg) {
729 const double a = vol * t_w;
730
731 for (int rr = 0; rr != nb_row_dofs; ++rr) {
732 auto t_col_v_base = col_data.getFTensor0N(gg, 0);
733 for (int cc = 0; cc != nb_col_dofs; ++cc) {
734 mat(rr, cc) += (ts_a * t_row_v_base * t_col_v_base) * a;
735
736 ++t_col_v_base;
737 }
738 ++t_row_v_base;
739 }
740 ++t_w;
741 }
742 CHKERR MatSetValues(getFEMethod()->ts_B, row_data, col_data, &mat(0, 0),
743 ADD_VALUES);
744 if (row_side != col_side || row_type != col_type) {
745 transMat.resize(nb_col_dofs, nb_row_dofs, false);
746 noalias(transMat) = trans(mat);
747 CHKERR MatSetValues(getFEMethod()->ts_B, col_data, row_data,
748 &transMat(0, 0), ADD_VALUES);
749 }
750 }
752 }
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::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
PetscReal ts_a
shift for U_t (see PETSc Time Solver)

◆ doWork() [2/2]

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

Definition at line 819 of file ElecPhysOperators2D.hpp.

821 {
823
824 const int nb_row_dofs = row_data.getIndices().size();
825 const int nb_col_dofs = col_data.getIndices().size();
826 if (nb_row_dofs && nb_col_dofs) {
827
828 mat.resize(nb_row_dofs, nb_col_dofs, false);
829 mat.clear();
830 const int nb_integration_pts = getGaussPts().size2();
831
832 auto t_row_v_base = row_data.getFTensor0N();
833 auto t_u_value = getFTensor0FromVec(dataU->mass_values);
834 auto t_v_value = getFTensor0FromVec(dataV->mass_values);
835
836 auto t_w = getFTensor0IntegrationWeight();
837 const double ts_a = getFEMethod()->ts_a;
838 const double vol = getMeasure();
839
840 for (int gg = 0; gg != nb_integration_pts; ++gg) {
841 const double a = vol * t_w;
842 double dfu = DRhs_u(t_u_value, t_v_value);
843 for (int rr = 0; rr != nb_row_dofs; ++rr) {
844 auto t_col_v_base = col_data.getFTensor0N(gg, 0);
845 for (int cc = 0; cc != nb_col_dofs; ++cc) {
846 mat(rr, cc) += ((ts_a - 0*dfu) * t_row_v_base * t_col_v_base) * a;
847
848 ++t_col_v_base;
849 }
850 ++t_row_v_base;
851
852 }
853 ++t_w;
854 ++t_u_value;
855 ++t_v_value;
856 }
857 CHKERR MatSetValues(getFEMethod()->ts_B, row_data, col_data, &mat(0, 0),
858 ADD_VALUES);
859 if (row_side != col_side || row_type != col_type) {
860 transMat.resize(nb_col_dofs, nb_row_dofs, false);
861 noalias(transMat) = trans(mat);
862 CHKERR MatSetValues(getFEMethod()->ts_B, col_data, row_data,
863 &transMat(0, 0), ADD_VALUES);
864 }
865 }
867 }
static auto getFTensor0FromVec(ublas::vector< T, A > &data)
Get tensor rank 0 (scalar) form data vector.
Definition: Templates.hpp:135

Member Data Documentation

◆ dataU

boost::shared_ptr<PreviousData>& ElectroPhysiology::OpAssembleLhsVV::dataU

Definition at line 814 of file ElecPhysOperators2D.hpp.

◆ dataV

boost::shared_ptr<PreviousData>& ElectroPhysiology::OpAssembleLhsVV::dataV

Definition at line 815 of file ElecPhysOperators2D.hpp.

◆ DRhs_u

FUval ElectroPhysiology::OpAssembleLhsVV::DRhs_u

Definition at line 817 of file ElecPhysOperators2D.hpp.

◆ mat

MatrixDouble ElectroPhysiology::OpAssembleLhsVV::mat
private

Definition at line 755 of file ElecPhysOperators.hpp.

◆ transMat

MatrixDouble ElectroPhysiology::OpAssembleLhsVV::transMat
private

Definition at line 755 of file ElecPhysOperators.hpp.


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