v0.15.0
Loading...
Searching...
No Matches
OpAleLhs_dX_dx Struct Reference

#include "users_modules/basic_finite_elements/src/HookeElement.hpp"

Inheritance diagram for OpAleLhs_dX_dx:
[legend]
Collaboration diagram for OpAleLhs_dX_dx:
[legend]

Public Member Functions

 OpAleLhs_dX_dx (const std::string row_field, const std::string col_field, boost::shared_ptr< DataAtIntegrationPts > &data_at_pts)
 
- Public Member Functions inherited from OpAssemble
 OpAssemble (const std::string row_field, const std::string col_field, boost::shared_ptr< DataAtIntegrationPts > data_at_pts, const char type, bool symm=false)
 
MoFEMErrorCode doWork (int row_side, int col_side, EntityType row_type, EntityType col_type, EntData &row_data, EntData &col_data)
 Do calculations for give operator.
 
MoFEMErrorCode doWork (int row_side, EntityType row_type, EntData &row_data)
 Operator for linear form, usually to calculate values on right hand side.
 
- Public Member Functions inherited from MoFEM::VolumeElementForcesAndSourcesCore::UserDataOperator
int getNumNodes ()
 get element number of nodes
 
const EntityHandlegetConn ()
 get element connectivity
 
double getVolume () const
 element volume (linear geometry)
 
doublegetVolume ()
 element volume (linear geometry)
 
FTensor::Tensor2< double *, 3, 3 > & getJac ()
 get element Jacobian
 
FTensor::Tensor2< double *, 3, 3 > & getInvJac ()
 get element inverse Jacobian
 
VectorDoublegetCoords ()
 nodal coordinates
 
VolumeElementForcesAndSourcesCoregetVolumeFE () const
 return pointer to Generic Volume Finite Element object
 
- 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.
 
EntityHandle getFEEntityHandle () const
 Return finite element entity handle.
 
int getFEDim () const
 Get dimension of finite element.
 
EntityType getFEType () const
 Get dimension of finite element.
 
boost::weak_ptr< SideNumbergetSideNumberPtr (const int side_number, const EntityType type)
 Get the side number pointer.
 
EntityHandle getSideEntity (const int side_number, const EntityType type)
 Get the side entity.
 
int getNumberOfNodesOnElement () const
 Get the number of nodes on finite element.
 
MoFEMErrorCode getProblemRowIndices (const std::string filed_name, const EntityType type, const int side, VectorInt &indices) const
 Get row indices.
 
MoFEMErrorCode getProblemColIndices (const std::string filed_name, const EntityType type, const int side, VectorInt &indices) const
 Get col indices.
 
const FEMethodgetFEMethod () const
 Return raw pointer to Finite Element Method object.
 
int getOpType () const
 Get operator types.
 
void setOpType (const OpType type)
 Set operator type.
 
void addOpType (const OpType type)
 Add operator type.
 
int getNinTheLoop () const
 get number of finite element in the loop
 
int getLoopSize () const
 get size of elements in the loop
 
std::string getFEName () const
 Get name of the element.
 
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
 
auto getFTensor0IntegrationWeight ()
 Get integration weights.
 
MatrixDoublegetCoordsAtGaussPts ()
 Gauss points and weight, matrix (nb. of points x 3)
 
auto getFTensor1CoordsAtGaussPts ()
 Get coordinates at integration points assuming linear geometry.
 
double getMeasure () const
 get measure of element
 
doublegetMeasure ()
 get measure of element
 
MoFEMErrorCode loopSide (const string &fe_name, ForcesAndSourcesCore *side_fe, const size_t dim, const EntityHandle ent_for_side=0, boost::shared_ptr< Range > fe_range=nullptr, 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.
 
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.
 
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.
 
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.
 
MoFEMErrorCode loopRange (const string &fe_name, ForcesAndSourcesCore *range_fe, boost::shared_ptr< Range > fe_range, const int verb=QUIET, const LogManager::SeverityLevel sev=Sev::noisy)
 Iterate over range of elements.
 
- Public Member Functions inherited from MoFEM::DataOperator
 DataOperator (const bool symm=true)
 
virtual ~DataOperator ()=default
 
virtual MoFEMErrorCode opLhs (EntitiesFieldData &row_data, EntitiesFieldData &col_data)
 
virtual MoFEMErrorCode opRhs (EntitiesFieldData &data, const bool error_if_no_base=false)
 
bool getSymm () const
 Get if operator uses symmetry of DOFs or not.
 
void setSymm ()
 set if operator is executed taking in account symmetry
 
void unSetSymm ()
 unset if operator is executed for non symmetric problem
 

Protected Member Functions

MoFEMErrorCode iNtegrate (EntData &row_data, EntData &col_data)
 Integrate tangent stiffness for material momentum.
 
- Protected Member Functions inherited from OpAssemble
virtual MoFEMErrorCode iNtegrate (EntData &row_data)
 
MoFEMErrorCode aSsemble (EntData &row_data, EntData &col_data)
 Assemble local entity block matrix.
 
MoFEMErrorCode aSsemble (EntData &row_data)
 Assemble local entity right-hand vector.
 
- Protected Member Functions inherited from MoFEM::VolumeElementForcesAndSourcesCore::UserDataOperator
MoFEMErrorCode setPtrFE (ForcesAndSourcesCore *ptr)
 
- Protected Member Functions inherited from MoFEM::ForcesAndSourcesCore::UserDataOperator

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
 
- Public Types inherited from MoFEM::DataOperator
using DoWorkLhsHookFunType
 
using DoWorkRhsHookFunType
 
- 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.
 
std::array< bool, MBMAXTYPE > doEntities
 If true operator is executed for entity.
 
booldoVertices
 \deprectaed If false skip vertices
 
booldoEdges
 \deprectaed If false skip edges
 
booldoQuads
 \deprectaed
 
booldoTris
 \deprectaed
 
booldoTets
 \deprectaed
 
booldoPrisms
 \deprectaed
 
- Static Public Attributes inherited from MoFEM::ForcesAndSourcesCore::UserDataOperator
static const char *const OpTypeNames []
 
- Protected Attributes inherited from OpAssemble
MatrixDouble K
 
MatrixDouble transK
 
VectorDouble nF
 
boost::shared_ptr< DataAtIntegrationPtsdataAtPts
 
VectorInt rowIndices
 
VectorInt colIndices
 
int nbRows
 number of dofs on rows
 
int nbCols
 number if dof on column
 
int nbIntegrationPts
 number of integration points
 
bool isDiag
 true if this block is on diagonal
 
- Protected Attributes inherited from MoFEM::ForcesAndSourcesCore::UserDataOperator
ForcesAndSourcesCoreptrFE
 

Detailed Description

Examples
HookeElement.cpp, and HookeElement.hpp.

Definition at line 547 of file HookeElement.hpp.

Constructor & Destructor Documentation

◆ OpAleLhs_dX_dx()

OpAleLhs_dX_dx::OpAleLhs_dX_dx ( const std::string row_field,
const std::string col_field,
boost::shared_ptr< DataAtIntegrationPts > & data_at_pts )
inline

Definition at line 549 of file HookeElement.hpp.

551 : OpAssemble(row_field, col_field, data_at_pts, OPROWCOL, false) {}
@ OPROWCOL
operator doWork is executed on FE rows &columns
OpAssemble(const std::string row_field, const std::string col_field, boost::shared_ptr< DataAtIntegrationPts > data_at_pts, const char type, bool symm=false)

Member Function Documentation

◆ iNtegrate()

MoFEMErrorCode OpAleLhs_dX_dx::iNtegrate ( EntData & row_data,
EntData & col_data )
protectedvirtual

Integrate tangent stiffness for material momentum.

Parameters
row_datarow data (consist base functions on row entity)
col_datacolumn data (consist base functions on column entity)
Returns
error code

Reimplemented from OpAssemble.

Definition at line 787 of file HookeElement.cpp.

788 {
790
791 // get sub-block (3x3) of local stiffens matrix, here represented by
792 // second order tensor
793 auto get_tensor2 = [](MatrixDouble &m, const int r, const int c) {
795 &m(r + 0, c + 0), &m(r + 0, c + 1), &m(r + 0, c + 2), &m(r + 1, c + 0),
796 &m(r + 1, c + 1), &m(r + 1, c + 2), &m(r + 2, c + 0), &m(r + 2, c + 1),
797 &m(r + 2, c + 2));
798 };
799
800 FTensor::Index<'i', 3> i;
801 FTensor::Index<'j', 3> j;
802 FTensor::Index<'k', 3> k;
803 FTensor::Index<'l', 3> l;
804 FTensor::Index<'m', 3> m;
805 FTensor::Index<'n', 3> n;
806
807 // get element volume
808 double vol = getVolume();
809
810 // get intergrayion weights
811 auto t_w = getFTensor0IntegrationWeight();
812
813 // get derivatives of base functions on rows
814 auto t_row_diff_base = row_data.getFTensor1DiffN<3>();
815 const int row_nb_base_fun = row_data.getN().size2();
816
817 auto t_invH = getFTensor2FromMat<3, 3>(*dataAtPts->invHMat);
818 auto &det_H = *dataAtPts->detHVec;
819
820 auto get_eshelby_stress_dx = [this]() {
822 t_eshelby_stress_dx;
823 int mm = 0;
824 for (int ii = 0; ii != 3; ++ii)
825 for (int jj = 0; jj != 3; ++jj)
826 for (int kk = 0; kk != 3; ++kk)
827 for (int ll = 0; ll != 3; ++ll)
828 t_eshelby_stress_dx.ptr(ii, jj, kk, ll) =
829 &(*dataAtPts->eshelbyStress_dx)(mm++, 0);
830 return t_eshelby_stress_dx;
831 };
832
833 auto t_eshelby_stress_dx = get_eshelby_stress_dx();
834
835 // iterate over integration points
836 for (int gg = 0; gg != nbIntegrationPts; ++gg) {
837
838 // calculate scalar weight times element volume
839 double a = t_w * vol * det_H[gg];
840
841 // iterate over row base functions
842 int rr = 0;
843 for (; rr != nbRows / 3; ++rr) {
844
845 // get sub matrix for the row
846 auto t_m = get_tensor2(K, 3 * rr, 0);
847
848 FTensor::Tensor1<double, 3> t_row_diff_base_pulled;
849 t_row_diff_base_pulled(i) = t_row_diff_base(j) * t_invH(j, i);
850
851 FTensor::Tensor3<double, 3, 3, 3> t_row_stress_dx;
852 t_row_stress_dx(i, k, l) =
853 a * t_row_diff_base_pulled(j) * t_eshelby_stress_dx(i, j, k, l);
854
855 // get derivatives of base functions for columns
856 auto t_col_diff_base = col_data.getFTensor1DiffN<3>(gg, 0);
857
858 // iterate column base functions
859 for (int cc = 0; cc != nbCols / 3; ++cc) {
860
861 t_m(i, k) += t_row_stress_dx(i, k, l) * t_col_diff_base(l);
862
863 // move to next column base function
864 ++t_col_diff_base;
865
866 // move to next block of local stiffens matrix
867 ++t_m;
868 }
869
870 // move to next row base function
871 ++t_row_diff_base;
872 }
873
874 for (; rr != row_nb_base_fun; ++rr)
875 ++t_row_diff_base;
876
877 ++t_w;
878 ++t_invH;
879 ++t_eshelby_stress_dx;
880 }
881
883}
constexpr double a
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
FTensor::Index< 'i', SPACE_DIM > i
const double c
speed of light (cm/ns)
const double n
refractive index of diffusive medium
FTensor::Index< 'l', 3 > l
FTensor::Index< 'j', 3 > j
FTensor::Index< 'k', 3 > k
FTensor::Tensor2< FTensor::PackPtr< double *, 1 >, Tensor_Dim1, Tensor_Dim2 > getFTensor2FromMat(MatrixDouble &data)
Get tensor rank 2 (matrix) form data matrix.
int r
Definition sdf.py:8
FTensor::Index< 'm', 3 > m
FTensor::Tensor1< FTensor::PackPtr< double *, Tensor_Dim >, Tensor_Dim > getFTensor1DiffN(const FieldApproximationBase base)
Get derivatives of base functions.
MatrixDouble & getN(const FieldApproximationBase base)
get base functions this return matrix (nb. of rows is equal to nb. of Gauss pts, nb....
auto getFTensor0IntegrationWeight()
Get integration weights.
int nbCols
number if dof on column
boost::shared_ptr< DataAtIntegrationPts > dataAtPts
int nbIntegrationPts
number of integration points
int nbRows
number of dofs on rows
MatrixDouble K

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