v0.14.0
Loading...
Searching...
No Matches
Public Member Functions | Public Attributes | List of all members
MixTransport::UnsaturatedFlowElement::OpDivTauU_HdivL2 Struct Reference

#include <users_modules/tutorials/cor-0to1/src/UnsaturatedFlow.hpp>

Inheritance diagram for MixTransport::UnsaturatedFlowElement::OpDivTauU_HdivL2:
[legend]
Collaboration diagram for MixTransport::UnsaturatedFlowElement::OpDivTauU_HdivL2:
[legend]

Public Member Functions

 OpDivTauU_HdivL2 (UnsaturatedFlowElement &ctx, const std::string &flux_name_col, const std::string &val_name_row)
 Constructor. More...
 
MoFEMErrorCode doWork (int row_side, int col_side, EntityType row_type, EntityType col_type, EntitiesFieldData::EntData &row_data, EntitiesFieldData::EntData &col_data)
 Do calculations. More...
 
- 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 Attributes

UnsaturatedFlowElementcTx
 
MatrixDouble nN
 
VectorDouble divVec
 
FTensor::Index< 'i', 3 > i
 
- 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...
 

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)>
 
- 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 738 of file UnsaturatedFlow.hpp.

Constructor & Destructor Documentation

◆ OpDivTauU_HdivL2()

MixTransport::UnsaturatedFlowElement::OpDivTauU_HdivL2::OpDivTauU_HdivL2 ( UnsaturatedFlowElement ctx,
const std::string &  flux_name_col,
const std::string &  val_name_row 
)
inline

Constructor.

Definition at line 746 of file UnsaturatedFlow.hpp.

Member Function Documentation

◆ doWork()

MoFEMErrorCode MixTransport::UnsaturatedFlowElement::OpDivTauU_HdivL2::doWork ( int  row_side,
int  col_side,
EntityType  row_type,
EntityType  col_type,
EntitiesFieldData::EntData &  row_data,
EntitiesFieldData::EntData &  col_data 
)
inline

Do calculations.

Parameters
row_sidelocal index of entity on row
col_sidelocal index of entity on column
row_typetype of row entity
col_typetype of col entity
row_datarow data structure carrying information about base functions, DOFs indices, etc.
col_datacolumn data structure carrying information about base functions, DOFs indices, etc.
Returns
error code
Examples
UnsaturatedFlow.hpp.

Definition at line 769 of file UnsaturatedFlow.hpp.

772 {
774 int nb_row = row_data.getFieldData().size();
775 int nb_col = col_data.getFieldData().size();
776 if (nb_row == 0)
778 if (nb_col == 0)
780 nN.resize(nb_row, nb_col, false);
781 divVec.resize(nb_row, false);
782 nN.clear();
783 // Get EntityHandle of the finite element
784 EntityHandle fe_ent = getNumeredEntFiniteElementPtr()->getEnt();
785 // Get material block id
786 int block_id;
787 CHKERR cTx.getMaterial(fe_ent, block_id);
788 // Get material block
789 boost::shared_ptr<GenericMaterial> &block_data = cTx.dMatMap.at(block_id);
790 // Get pressure
791 auto t_h = getFTensor0FromVec(cTx.valuesAtGaussPts);
792 // Get flux
793 auto t_flux = getFTensor1FromMat<3>(cTx.fluxesAtGaussPts);
794 // Coords at integration points
795 auto t_coords = getFTensor1CoordsAtGaussPts();
796 // Get integration weight
797 auto t_w = getFTensor0IntegrationWeight();
798 // Get base function
799 auto t_n_hdiv_row = row_data.getFTensor1N<3>();
800 // Get derivatives of base functions
802 auto t_base_diff_hdiv = row_data.getFTensor2DiffN<3, 3>();
803 // Get volume
804 double vol = getVolume();
805 int nb_gauss_pts = row_data.getN().size1();
806 for (int gg = 0; gg != nb_gauss_pts; gg++) {
807 block_data->h = t_h;
808 block_data->x = t_coords(0);
809 block_data->y = t_coords(1);
810 block_data->z = t_coords(2);
811 CHKERR block_data->calK();
812 CHKERR block_data->calDiffK();
813 const double K = block_data->K;
814 // const double z = t_coords(2);
815 const double KK = K * K;
816 const double diffK = block_data->diffK;
817 double alpha = t_w * vol;
818 for (auto &v : divVec) {
819 v = t_base_diff_hdiv(i, i);
820 ++t_base_diff_hdiv;
821 }
822 noalias(nN) -= alpha * outer_prod(divVec, col_data.getN(gg));
823 FTensor::Tensor0<double *> t_a(&*nN.data().begin());
824 for (int rr = 0; rr != nb_row; rr++) {
825 double beta = alpha * (-diffK / KK) * (t_n_hdiv_row(i) * t_flux(i));
826 auto t_n_col = col_data.getFTensor0N(gg, 0);
827 for (int cc = 0; cc != nb_col; cc++) {
828 t_a += beta * t_n_col;
829 ++t_n_col;
830 ++t_a;
831 }
832 ++t_n_hdiv_row;
833 }
834 ++t_w;
835 ++t_coords;
836 ++t_h;
837 ++t_flux;
838 }
839 CHKERR MatSetValues(getFEMethod()->ts_B, nb_row,
840 &row_data.getIndices()[0], nb_col,
841 &col_data.getIndices()[0], &nN(0, 0), ADD_VALUES);
843 }
static Index< 'K', 3 > K
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:447
#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
const double v
phase velocity of light in medium (cm/ns)
MoFEMErrorCode MatSetValues(Mat M, const EntitiesFieldData::EntData &row_data, const EntitiesFieldData::EntData &col_data, const double *ptr, InsertMode iora)
Assemble PETSc matrix.
VectorDouble valuesAtGaussPts
values at integration points on element
MatrixDouble fluxesAtGaussPts
fluxes at integration points on element
MaterialsDoubleMap dMatMap
materials database
virtual MoFEMErrorCode getMaterial(const EntityHandle ent, int &block_id) const
For given element handle get material block Id.
auto getFTensor1CoordsAtGaussPts()
Get coordinates at integration points assuming linear geometry.
auto getFTensor0IntegrationWeight()
Get integration weights.
boost::shared_ptr< const NumeredEntFiniteElement > getNumeredEntFiniteElementPtr() const
Return raw pointer to NumeredEntFiniteElement.
const FEMethod * getFEMethod() const
Return raw pointer to Finite Element Method object.

Member Data Documentation

◆ cTx

UnsaturatedFlowElement& MixTransport::UnsaturatedFlowElement::OpDivTauU_HdivL2::cTx
Examples
UnsaturatedFlow.hpp.

Definition at line 741 of file UnsaturatedFlow.hpp.

◆ divVec

VectorDouble MixTransport::UnsaturatedFlowElement::OpDivTauU_HdivL2::divVec
Examples
UnsaturatedFlow.hpp.

Definition at line 754 of file UnsaturatedFlow.hpp.

◆ i

FTensor::Index<'i', 3> MixTransport::UnsaturatedFlowElement::OpDivTauU_HdivL2::i
Examples
UnsaturatedFlow.hpp.

Definition at line 755 of file UnsaturatedFlow.hpp.

◆ nN

MatrixDouble MixTransport::UnsaturatedFlowElement::OpDivTauU_HdivL2::nN
Examples
UnsaturatedFlow.hpp.

Definition at line 753 of file UnsaturatedFlow.hpp.


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