v0.15.0
Loading...
Searching...
No Matches
EshelbianPlasticity::OpSpatialRotation_domega_domega Struct Reference

#include "users_modules/eshelbian_plasticity/src/EshelbianPlasticity.hpp"

Inheritance diagram for EshelbianPlasticity::OpSpatialRotation_domega_domega:
[legend]
Collaboration diagram for EshelbianPlasticity::OpSpatialRotation_domega_domega:
[legend]

Public Member Functions

 OpSpatialRotation_domega_domega (std::string row_field, std::string col_field, boost::shared_ptr< DataAtIntegrationPts > data_ptr, double alpha_omega)
 
MoFEMErrorCode integrate (EntData &row_data, EntData &col_data)
 
- Public Member Functions inherited from EshelbianPlasticity::OpAssembleVolume
MoFEMErrorCode assemble (int row_side, int col_side, EntityType row_type, EntityType col_type, EntData &row_data, EntData &col_data)
 
- Public Member Functions inherited from EshelbianPlasticity::OpAssembleBasic< VolUserDataOperator >
 OpAssembleBasic (const std::string &field_name, boost::shared_ptr< DataAtIntegrationPts > data_ptr, const char type)
 
 OpAssembleBasic (std::string row_field, std::string col_field, boost::shared_ptr< DataAtIntegrationPts > data_ptr, const char type, const bool assemble_symmetry, ScaleOff scale_off=[]() { return 1;})
 
 OpAssembleBasic (const FieldSpace space)
 
virtual MoFEMErrorCode integrate (EntData &data)
 
virtual MoFEMErrorCode integrate (int row_side, EntityType row_type, EntData &data)
 
virtual MoFEMErrorCode assemble (EntData &data)
 
virtual MoFEMErrorCode assemble (int row_side, EntityType row_type, EntData &data)
 
MoFEMErrorCode doWork (int side, EntityType type, EntData &data)
 
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
 
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
 

Private Attributes

double alphaOmega
 

Additional Inherited Members

- Public Types inherited from EshelbianPlasticity::OpAssembleVolume
using OP = OpAssembleBasic<VolUserDataOperator>
 
- Public Types inherited from EshelbianPlasticity::OpAssembleBasic< VolUserDataOperator >
using ScaleOff
 
- 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 EshelbianPlasticity::OpAssembleBasic< VolUserDataOperator >
const bool assembleSymmetry
 
boost::shared_ptr< DataAtIntegrationPtsdataAtPts
 data at integration pts
 
VectorDouble nF
 local right hand side vector
 
MatrixDouble K
 local tangent matrix
 
MatrixDouble transposeK
 
ScaleOff scaleOff
 
- 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 Member Functions inherited from MoFEM::VolumeElementForcesAndSourcesCore::UserDataOperator
MoFEMErrorCode setPtrFE (ForcesAndSourcesCore *ptr)
 
- Protected Member Functions inherited from MoFEM::ForcesAndSourcesCore::UserDataOperator
- Protected Attributes inherited from MoFEM::ForcesAndSourcesCore::UserDataOperator
ForcesAndSourcesCoreptrFE
 
- Static Protected Attributes inherited from EshelbianPlasticity::OpAssembleVolume
static std::map< std::pair< std::string, std::string >, MatrixDouble > mapMatrix
 

Detailed Description

Definition at line 799 of file EshelbianPlasticity.hpp.

Constructor & Destructor Documentation

◆ OpSpatialRotation_domega_domega()

EshelbianPlasticity::OpSpatialRotation_domega_domega::OpSpatialRotation_domega_domega ( std::string row_field,
std::string col_field,
boost::shared_ptr< DataAtIntegrationPts > data_ptr,
double alpha_omega )
inline

Definition at line 800 of file EshelbianPlasticity.hpp.

834 {
835
836struct CGGUserPolynomialBase : public TetPolynomialBase {
837
838 using CachePhi = boost::tuple<int, int, MatrixDouble>;
839
840 CGGUserPolynomialBase(boost::shared_ptr<CachePhi> cache_phi = nullptr);
841 ~CGGUserPolynomialBase() = default;
842
843 MoFEMErrorCode query_interface(boost::typeindex::type_index type_index,
844 BaseFunctionUnknownInterface **iface) const;
845 MoFEMErrorCode getValue(MatrixDouble &pts,
846 boost::shared_ptr<BaseFunctionCtx> ctx_ptr);
847
848private:
849 MatrixDouble shapeFun;
850 boost::shared_ptr<CachePhi> cachePhiPtr;
851
852 MoFEMErrorCode getValueHdivForCGGBubble(MatrixDouble &pts);
853};
854
855struct ContactTree;
856
859enum StretchSelector {
860 LINEAR,
861 LOG /*linar extenion*/,
862 LOG_QUADRATIC /*quadratic extension*/,
864};
866
867using MatrixPtr = boost::shared_ptr<MatrixDouble>;
868using VectorPtr = boost::shared_ptr<VectorDouble>;
869
870using EntData = EntitiesFieldData::EntData;
871using UserDataOperator = ForcesAndSourcesCore::UserDataOperator;
872using VolUserDataOperator = VolumeElementForcesAndSourcesCore::UserDataOperator;
873using FaceUserDataOperator = FaceElementForcesAndSourcesCore::UserDataOperator;
874using EleOnSide = PipelineManager::ElementsAndOpsByDim<SPACE_DIM>::FaceSideEle;
875using SideEleOp = EleOnSide::UserDataOperator;
876
877struct PhysicalEquations;
878
879struct AnalyticalExprPython;
880
882 : public boost::enable_shared_from_this<DataAtIntegrationPts> {
883
884 MatrixDouble approxPAtPts;
885 MatrixDouble approxSigmaAtPts;
886 MatrixDouble divPAtPts;
887 MatrixDouble divSigmaAtPts;
888 MatrixDouble wL2AtPts;
889 MatrixDouble wL2DotAtPts;
890 MatrixDouble wL2DotDotAtPts;
891 MatrixDouble logStretchTensorAtPts;
892 MatrixDouble stretchTensorAtPts;
893 VectorDouble jacobianAtPts;
894 MatrixDouble tractionAtPts;
895
896 MatrixDouble diffStretchTensorAtPts;
897 VectorDouble detStretchTensorAtPts;
898 MatrixDouble detStretchTensorAtPts_du;
899 MatrixDouble logStretchDotTensorAtPts;
900 MatrixDouble gradLogStretchDotTensorAtPts;
901 MatrixDouble rotAxisAtPts;
902 MatrixDouble rotAxisDotAtPts;
903 MatrixDouble rotAxisGradDotAtPts;
904 MatrixDouble rotAxis0AtPts;
905 MatrixDouble WAtPts;
906 MatrixDouble W0AtPts;
907 MatrixDouble GAtPts;
908 MatrixDouble G0AtPts;
909 MatrixDouble wH1AtPts;
910 MatrixDouble XH1AtPts;
911 MatrixDouble contactL2AtPts;
912 MatrixDouble wGradH1AtPts;
913 MatrixDouble logStretch2H1AtPts;
914 MatrixDouble logStretchTotalTensorAtPts;
915
916 MatrixDouble hAtPts;
917 MatrixDouble hdOmegaAtPts;
918 MatrixDouble hdLogStretchAtPts;
919 MatrixDouble leviKirchhoffAtPts;
920 MatrixDouble leviKirchhoffdOmegaAtPts;
921 MatrixDouble leviKirchhoffdLogStreatchAtPts;
922 MatrixDouble leviKirchhoffPAtPts;
923 MatrixDouble adjointPdstretchAtPts;
924 MatrixDouble adjointPdUAtPts;
925 MatrixDouble adjointPdUdOmegaAtPts;
926 MatrixDouble adjointPdUdPAtPts;
927 MatrixDouble rotMatAtPts;
928 MatrixDouble PAtPts;
929 MatrixDouble SigmaAtPts;
930 VectorDouble energyAtPts; //< this is density of energy at integration points
931 MatrixDouble flowL2AtPts;
932
933 MatrixDouble varRotAxis;
934 MatrixDouble varLogStreach;
935 MatrixDouble varPiola;
936 MatrixDouble varDivPiola;
937 MatrixDouble varWL2;
938
939 MatrixDouble P_du;
940
941 MatrixDouble eigenVals;
942 MatrixDouble eigenVecs;
943 VectorInt nbUniq;
944 MatrixDouble eigenValsC;
945 MatrixDouble eigenVecsC;
946 VectorInt nbUniqC;
947
948 MatrixDouble matD;
949 MatrixDouble matAxiatorD;
950 MatrixDouble matDeviatorD;
951 MatrixDouble matInvD;
952
953 MatrixDouble internalStressAtPts;
954
955 MatrixDouble facePiolaAtPts;
956 MatrixDouble gradHybridDispAtPts;
957 MatrixDouble faceMaterialForceAtPts;
958 VectorDouble normalPressureAtPts;
959
960 double mu;
961 double lambda;
962 double piolaScale = 1.;
963 double faceEnergy = 0.;
964
965 inline auto getPiolaScalePtr() {
966 return boost::shared_ptr<double>(shared_from_this(), &piolaScale);
967 }
968
969 inline MatrixPtr getApproxSigmaAtPts() {
970 return boost::shared_ptr<MatrixDouble>(shared_from_this(),
971 &approxSigmaAtPts);
972 }
973 inline MatrixPtr getApproxPAtPts() {
974 return boost::shared_ptr<MatrixDouble>(shared_from_this(), &approxPAtPts);
975 }
976
977 inline MatrixPtr getDivPAtPts() {
978 return boost::shared_ptr<MatrixDouble>(shared_from_this(), &divPAtPts);
979 }
980
981 inline MatrixPtr getDivSigmaAtPts() {
982 return boost::shared_ptr<MatrixDouble>(shared_from_this(), &divSigmaAtPts);
983 }
984
985 inline MatrixPtr getSmallWL2AtPts() {
986 return boost::shared_ptr<MatrixDouble>(shared_from_this(), &wL2AtPts);
987 }
988
989 inline MatrixPtr getSmallWL2DotAtPts() {
990 return boost::shared_ptr<MatrixDouble>(shared_from_this(), &wL2DotAtPts);
991 }
992
993 inline MatrixPtr getSmallWL2DotDotAtPts() {
994 return boost::shared_ptr<MatrixDouble>(shared_from_this(), &wL2DotDotAtPts);
995 }
996
997 inline MatrixPtr getLogStretchTensorAtPts() {
998 return boost::shared_ptr<MatrixDouble>(shared_from_this(),
999 &logStretchTensorAtPts);
1000 }
1001
1002 inline MatrixPtr getStretchTensorAtPts() {
1003 return boost::shared_ptr<MatrixDouble>(shared_from_this(),
1004 &stretchTensorAtPts);
1005 }
1006
1007 inline MatrixPtr getLogStretchDotTensorAtPts() {
1008 return boost::shared_ptr<MatrixDouble>(shared_from_this(),
1009 &logStretchDotTensorAtPts);
1010 }
1011
1012 inline MatrixPtr getGradLogStretchDotTensorAtPts() {
1013 return boost::shared_ptr<MatrixDouble>(shared_from_this(),
1014 &gradLogStretchDotTensorAtPts);
1015 }
1016
1017 inline MatrixPtr getRotAxisAtPts() {
1018 return boost::shared_ptr<MatrixDouble>(shared_from_this(), &rotAxisAtPts);
1019 }
1020
1021 inline MatrixPtr getRotAxis0AtPts() {
1022 return boost::shared_ptr<MatrixDouble>(shared_from_this(), &rotAxis0AtPts);
1023 }
1024
1025 inline MatrixPtr getRotAxisDotAtPts() {
1026 return boost::shared_ptr<MatrixDouble>(shared_from_this(),
1027 &rotAxisDotAtPts);
1028 }
1029
1030 inline MatrixPtr getRotAxisGradDotAtPts() {
1031 return boost::shared_ptr<MatrixDouble>(shared_from_this(),
1032 &rotAxisGradDotAtPts);
1033 }
1034
1035 inline MatrixPtr getBigGAtPts() {
1036 return boost::shared_ptr<MatrixDouble>(shared_from_this(), &GAtPts);
1037 }
1038
1039 inline MatrixPtr getBigG0AtPts() {
1040 return boost::shared_ptr<MatrixDouble>(shared_from_this(), &GAtPts);
1041 }
1042
1043 inline MatrixPtr getMatDPtr() {
1044 return boost::shared_ptr<MatrixDouble>(shared_from_this(), &matD);
1045 }
1046
1047 inline MatrixPtr getMatInvDPtr() {
1048 return boost::shared_ptr<MatrixDouble>(shared_from_this(), &matInvD);
1049 }
1050
1051 inline MatrixPtr getMatAxiatorDPtr() {
1052 return boost::shared_ptr<MatrixDouble>(shared_from_this(), &matAxiatorD);
1053 }
1054
1055 inline MatrixPtr getMatDeviatorDPtr() {
1056 return boost::shared_ptr<MatrixDouble>(shared_from_this(), &matDeviatorD);
1057 }
1058
1059 inline MatrixPtr getSmallWH1AtPts() {
1060 return boost::shared_ptr<MatrixDouble>(shared_from_this(), &wH1AtPts);
1061 }
1062
1063 inline MatrixPtr getLargeXH1AtPts() {
1064 return boost::shared_ptr<MatrixDouble>(shared_from_this(), &XH1AtPts);
1065 }
1066
1067 inline MatrixPtr getContactL2AtPts() {
1068 return boost::shared_ptr<MatrixDouble>(shared_from_this(), &contactL2AtPts);
1069 }
1070
1071 inline MatrixPtr getSmallWGradH1AtPts() {
1072 return boost::shared_ptr<MatrixDouble>(shared_from_this(), &wGradH1AtPts);
1073 }
1074
1075 inline VectorPtr getJacobianAtPts() {
1076 return boost::shared_ptr<VectorDouble>(shared_from_this(), &jacobianAtPts);
1077 }
1078
1079 inline MatrixPtr getLeviKirchhoffAtPts() {
1080 return boost::shared_ptr<MatrixDouble>(shared_from_this(),
1081 &leviKirchhoffAtPts);
1082 };
1083
1084 inline MatrixPtr getVarRotAxisPts() {
1085 return boost::shared_ptr<MatrixDouble>(shared_from_this(), &varRotAxis);
1086 }
1087
1088 inline MatrixPtr getVarLogStreachPts() {
1089 return boost::shared_ptr<MatrixDouble>(shared_from_this(), &varLogStreach);
1090 }
1091
1092 inline MatrixPtr getVarPiolaPts() {
1093 return boost::shared_ptr<MatrixDouble>(shared_from_this(), &varPiola);
1094 }
1095
1096 inline MatrixPtr getDivVarPiolaPts() {
1097 return boost::shared_ptr<MatrixDouble>(shared_from_this(), &varDivPiola);
1098 }
1099
1100 inline MatrixPtr getVarWL2Pts() {
1101 return boost::shared_ptr<MatrixDouble>(shared_from_this(), &varWL2);
1102 }
1103
1104 inline MatrixPtr getInternalStressAtPts() {
1105 return boost::shared_ptr<MatrixDouble>(shared_from_this(),
1106 &internalStressAtPts);
1107 }
1108
1109 inline MatrixPtr getFacePiolaAtPts() {
1110 return boost::shared_ptr<MatrixDouble>(shared_from_this(), &facePiolaAtPts);
1111 }
1112
1113 inline MatrixPtr getGradHybridDispAtPts() {
1114 return boost::shared_ptr<MatrixDouble>(shared_from_this(),
1115 &gradHybridDispAtPts);
1116 }
1117
1118 inline MatrixPtr getEigenValsAtPts() {
1119 return boost::shared_ptr<MatrixDouble>(shared_from_this(), &eigenVals);
1120 }
1121
1122 boost::shared_ptr<PhysicalEquations> physicsPtr;
1123};
1124
1125struct OpJacobian;
1126
1127// Forward declarations
1128struct ExternalStrain;
1129using ExternalStrainVec = std::vector<ExternalStrain>;
1130struct PhysicalEquations {
1131
1132 typedef FTensor::Tensor1<adouble, 3> ATensor1;
1133 typedef FTensor::Tensor2<adouble, 3, 3> ATensor2;
1134 typedef FTensor::Tensor3<adouble, 3, 3, 3> ATensor3;
1135 typedef FTensor::Tensor1<double, 3> DTensor1;
1136 typedef FTensor::Tensor2<double, 3, 3> DTensor2;
1137 typedef FTensor::Tensor3<double, 3, 3, 3> DTensor3;
1138
1140
1141 typedef FTensor::Tensor2<FTensor::PackPtr<double *, 1>, 3, 3> DTensor2Ptr;
1142 typedef FTensor::Tensor3<FTensor::PackPtr<double *, 1>, 3, 3, 3> DTensor3Ptr;
1143
1144 PhysicalEquations() = delete;
1145 PhysicalEquations(const int size_active, const int size_dependent)
1146 : activeVariables(size_active, 0),
1147 dependentVariablesPiola(size_dependent, 0),
1148 dependentVariablesPiolaDirevatives(size_dependent * size_active, 0) {}
1149 virtual ~PhysicalEquations() = default;
1150
1151 virtual MoFEMErrorCode recordTape(const int tag, DTensor2Ptr *t_h) = 0;
1152
1153 virtual UserDataOperator *
1154 returnOpJacobian(const int tag, const bool eval_rhs, const bool eval_lhs,
1155 boost::shared_ptr<DataAtIntegrationPts> data_ptr,
1156 boost::shared_ptr<PhysicalEquations> physics_ptr);
1157
1158 virtual VolUserDataOperator *
1159 returnOpSpatialPhysical(const std::string &field_name,
1160 boost::shared_ptr<DataAtIntegrationPts> data_ptr,
1161 const double alpha_u);
1162
1163 virtual VolUserDataOperator *
1164 returnOpSpatialPhysicalExternalStrain(
1165 const std::string &field_name,
1166 boost::shared_ptr<DataAtIntegrationPts> data_ptr,
1167 boost::shared_ptr<ExternalStrainVec> external_strain_vec_ptr,
1168 std::map<std::string, boost::shared_ptr<ScalingMethod>> smv);
1169
1170 virtual VolUserDataOperator *returnOpSpatialPhysical_du_du(
1171 std::string row_field, std::string col_field,
1172 boost::shared_ptr<DataAtIntegrationPts> data_ptr, const double alpha);
1173
1174 virtual VolUserDataOperator *
1175 returnOpCalculateEnergy(boost::shared_ptr<DataAtIntegrationPts> data_ptr,
1176 boost::shared_ptr<double> total_energy_ptr);
1177
1178 virtual VolUserDataOperator *returnOpCalculateStretchFromStress(
1179 boost::shared_ptr<DataAtIntegrationPts> data_ptr,
1180 boost::shared_ptr<PhysicalEquations> physics_ptr);
1181
1182 virtual VolUserDataOperator *returnOpCalculateVarStretchFromStress(
1183 boost::shared_ptr<DataAtIntegrationPts> data_ptr,
1184 boost::shared_ptr<PhysicalEquations> physics_ptr);
1185
1186 virtual VolUserDataOperator *
1187 returnOpSetScale(boost::shared_ptr<double> scale_ptr,
1188 boost::shared_ptr<PhysicalEquations> physics_ptr);
1189
1190 std::vector<double> activeVariables;
1191 std::vector<double> dependentVariablesPiola;
1192 std::vector<double> dependentVariablesPiolaDirevatives;
1193
1194 /** \name Active variables */
1195
1196 /**@{*/
1197
1198 template <int S>
1199 inline static DTensor2Ptr get_VecTensor2(std::vector<double> &v) {
1200 return DTensor2Ptr(&v[S + 0], &v[S + 1], &v[S + 2], &v[S + 3], &v[S + 4],
1201 &v[S + 5], &v[S + 6], &v[S + 7], &v[S + 8]);
1202 }
1203
1204 template <int S>
1205 inline static DTensor0Ptr get_VecTensor0(std::vector<double> &v) {
1206 return DTensor0Ptr(&v[S + 0]);
1207 }
1208
1209 template <int S0>
1210 inline static DTensor3Ptr get_vecTensor3(std::vector<double> &v,
1211 const int nba) {
1212
1213 const int A00 = nba * 0 + S0;
1214 const int A01 = nba * 1 + S0;
1215 const int A02 = nba * 2 + S0;
1216 const int A10 = nba * 3 + S0;
1217 const int A11 = nba * 4 + S0;
1218 const int A12 = nba * 5 + S0;
1219 const int A20 = nba * 6 + S0;
1220 const int A21 = nba * 7 + S0;
1221 const int A22 = nba * 8 + S0;
1222
1223 return DTensor3Ptr(
1224
1225 &v[A00 + 0], &v[A00 + 1], &v[A00 + 2], &v[A01 + 0], &v[A01 + 1],
1226 &v[A01 + 2], &v[A02 + 0], &v[A02 + 1], &v[A02 + 2],
1227
1228 &v[A10 + 0], &v[A10 + 1], &v[A10 + 2], &v[A11 + 0], &v[A11 + 1],
1229 &v[A11 + 2], &v[A12 + 0], &v[A12 + 1], &v[A12 + 2],
1230
1231 &v[A20 + 0], &v[A20 + 1], &v[A20 + 2], &v[A21 + 0], &v[A21 + 1],
1232 &v[A21 + 2], &v[A22 + 0], &v[A22 + 1], &v[A22 + 2]
1233
1234 );
1235 }
1236
1237 /**@}*/
1238
1239 /** \name Active variables */
1240
1241 /**@{*/
1242
1243 inline DTensor2Ptr get_h() { return get_VecTensor2<0>(activeVariables); }
1244
1245 /**@}*/
1246
1247 /** \name Dependent variables */
1248
1249 /**@{*/
1250
1251 inline DTensor2Ptr get_P() {
1252 return get_VecTensor2<0>(dependentVariablesPiola);
1253 }
1254
1255 /**@}*/
1256
1257 /** \name Derivatives of dependent variables */
1258
1259 /**@{*/
1260
1261 inline DTensor3Ptr get_P_dh0() {
1262 return get_vecTensor3<0>(dependentVariablesPiolaDirevatives,
1263 activeVariables.size());
1264 }
1265 inline DTensor3Ptr get_P_dh1() {
1266 return get_vecTensor3<3>(dependentVariablesPiolaDirevatives,
1267 activeVariables.size());
1268 }
1269 inline DTensor3Ptr get_P_dh2() {
1270 return get_vecTensor3<6>(dependentVariablesPiolaDirevatives,
1271 activeVariables.size());
1272 }
1273
1274 /**@}*/
1275};
1276
1277struct BcDisp {
1278 BcDisp(std::string name, std::vector<double> attr, Range faces);
1279 std::string blockName;
1280 Range faces;
1281 VectorDouble3 vals;
1282 VectorInt3 flags;
1283};
1284using BcDispVec = std::vector<BcDisp>;
1285
1286struct BcRot {
1287 BcRot(std::string name, std::vector<double> attr, Range faces);
1288 std::string blockName;
1289 Range faces;
1290 VectorDouble vals;
1291 double theta;
1292};
1293using BcRotVec = std::vector<BcRot>;
1294
1295typedef std::vector<Range> TractionFreeBc;
1296
1297struct TractionBc {
1298 TractionBc(std::string name, std::vector<double> attr, Range faces);
1299 std::string blockName;
1300 Range faces;
1301 VectorDouble3 vals;
1302 VectorInt3 flags;
1303};
1304using TractionBcVec = std::vector<TractionBc>;
1305
1306struct NormalDisplacementBc {
1307 NormalDisplacementBc(std::string name, std::vector<double> attr, Range faces);
1308 std::string blockName;
1309 Range faces;
1310 double val;
1311};
1312using NormalDisplacementBcVec = std::vector<NormalDisplacementBc>;
1313
1314struct AnalyticalDisplacementBc {
1315 AnalyticalDisplacementBc(std::string name, std::vector<double> attr,
1316 Range faces);
1317 std::string blockName;
1318 Range faces;
1319 VectorInt3 flags;
1320};
1321using AnalyticalDisplacementBcVec = std::vector<AnalyticalDisplacementBc>;
1322
1323struct AnalyticalTractionBc {
1324 AnalyticalTractionBc(std::string name, std::vector<double> attr, Range faces);
1325 std::string blockName;
1326 Range faces;
1327 VectorInt3 flags;
1328};
1329using AnalyticalTractionBcVec = std::vector<AnalyticalTractionBc>;
1330
1331struct PressureBc {
1332 PressureBc(std::string name, std::vector<double> attr, Range faces);
1333 std::string blockName;
1334 Range faces;
1335 double val;
1336};
1337using PressureBcVec = std::vector<PressureBc>;
1338
1339struct ExternalStrain {
1340 ExternalStrain(std::string name, std::vector<double> attr, Range ents);
1341 std::string blockName;
1342 Range ents;
1343 double val;
1344 double bulkModulusK;
1345};
1346
1347 #ifdef ENABLE_PYTHON_BINDING
1348struct AnalyticalExprPython {
1349 AnalyticalExprPython() = default;
1350 virtual ~AnalyticalExprPython() = default;
1351
1352 MoFEMErrorCode analyticalExprInit(const std::string py_file);
1353 MoFEMErrorCode evalAnalyticalDisp(double delta_t, double t, np::ndarray x,
1354 np::ndarray y, np::ndarray z,
1355 np::ndarray nx, np::ndarray ny,
1356 np::ndarray nz,
1357 const std::string &block_name,
1358 np::ndarray &analytical_expr);
1359 MoFEMErrorCode evalAnalyticalTraction(double delta_t, double t, np::ndarray x,
1360 np::ndarray y, np::ndarray z,
1361 np::ndarray nx, np::ndarray ny,
1362 np::ndarray nz,
1363 const std::string &block_name,
1364 np::ndarray &analytical_expr);
1365
1366 template <typename T>
1367 inline std::vector<T>
1368 py_list_to_std_vector(const boost::python::object &iterable) {
1369 return std::vector<T>(boost::python::stl_input_iterator<T>(iterable),
1370 boost::python::stl_input_iterator<T>());
1371 }
1372
1373private:
1374 bp::object mainNamespace;
1375 bp::object analyticalDispFun;
1376 bp::object analyticalTractionFun;
1377};
1378
1379extern boost::weak_ptr<AnalyticalExprPython> AnalyticalExprPythonWeakPtr;
1380 #endif
1381
1382 #include "EshelbianCore.hpp"
1383 #include "EshelbianOperators.hpp"
1384
1385} // namespace EshelbianPlasticity
1386
1387#endif //__ESHELBIAN_PLASTICITY_HPP__ge64_NodeáM
1388
ForcesAndSourcesCore::UserDataOperator UserDataOperator
static double lambda
const double v
phase velocity of light in medium (cm/ns)
boost::shared_ptr< VectorDouble > VectorPtr
std::vector< AnalyticalTractionBc > AnalyticalTractionBcVec
std::vector< BcDisp > BcDispVec
std::vector< NormalDisplacementBc > NormalDisplacementBcVec
std::vector< BcRot > BcRotVec
std::vector< ExternalStrain > ExternalStrainVec
std::vector< AnalyticalDisplacementBc > AnalyticalDisplacementBcVec
std::vector< Range > TractionFreeBc
std::vector< PressureBc > PressureBcVec
std::vector< TractionBc > TractionBcVec
boost::shared_ptr< MatrixDouble > MatrixPtr
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
constexpr double t
plate stiffness
Definition plate.cpp:58
constexpr auto field_name
PipelineManager::ElementsAndOpsByDim< SPACE_DIM >::FaceSideEle EleOnSide
Data on single entity (This is passed as argument to DataOperator::doWork)

Member Function Documentation

◆ integrate()

MoFEMErrorCode OpSpatialRotation_domega_domega::integrate ( EntData & row_data,
EntData & col_data )
virtual

Reimplemented from EshelbianPlasticity::OpAssembleBasic< VolUserDataOperator >.

Examples
EshelbianOperators.cpp.

Definition at line 2703 of file EshelbianOperators.cpp.

2704 {
2706 int nb_integration_pts = getGaussPts().size2();
2707 int row_nb_dofs = row_data.getIndices().size();
2708 int col_nb_dofs = col_data.getIndices().size();
2709 auto get_ftensor2 = [](MatrixDouble &m, const int r, const int c) {
2711
2712 &m(r + 0, c + 0), &m(r + 0, c + 1), &m(r + 0, c + 2),
2713
2714 &m(r + 1, c + 0), &m(r + 1, c + 1), &m(r + 1, c + 2),
2715
2716 &m(r + 2, c + 0), &m(r + 2, c + 1), &m(r + 2, c + 2)
2717
2718 );
2719 };
2720 FTensor::Index<'i', 3> i;
2721 FTensor::Index<'j', 3> j;
2722 FTensor::Index<'k', 3> k;
2723 FTensor::Index<'l', 3> l;
2724 FTensor::Index<'m', 3> m;
2725 FTensor::Index<'n', 3> n;
2726
2728
2729 auto v = getVolume();
2730 auto ts_a = getTSa();
2731 auto t_w = getFTensor0IntegrationWeight();
2732 auto t_levi_kirchhoff_domega =
2733 getFTensor2FromMat<3, 3>(dataAtPts->leviKirchhoffdOmegaAtPts);
2734 int row_nb_base_functions = row_data.getN().size2();
2735 auto t_row_base_fun = row_data.getFTensor0N();
2736 auto t_row_grad_fun = row_data.getFTensor1DiffN<3>();
2737
2738 auto time_step = getTStimeStep();
2739 for (int gg = 0; gg != nb_integration_pts; ++gg) {
2740 double a = v * t_w;
2741 double c = (a * alphaOmega) * (ts_a /*/ time_step*/);
2742
2743 int rr = 0;
2744 for (; rr != row_nb_dofs / 3; ++rr) {
2745 auto t_m = get_ftensor2(K, 3 * rr, 0);
2746 const double b = a * t_row_base_fun;
2747 auto t_col_base_fun = col_data.getFTensor0N(gg, 0);
2748 auto t_col_grad_fun = col_data.getFTensor1DiffN<3>(gg, 0);
2749 for (int cc = 0; cc != col_nb_dofs / 3; ++cc) {
2750 t_m(k, l) -= (b * t_col_base_fun) * t_levi_kirchhoff_domega(k, l);
2751 t_m(k, l) += t_kd(k, l) * (c * (t_row_grad_fun(i) * t_col_grad_fun(i)));
2752 ++t_m;
2753 ++t_col_base_fun;
2754 ++t_col_grad_fun;
2755 }
2756 ++t_row_base_fun;
2757 ++t_row_grad_fun;
2758 }
2759 for (; rr != row_nb_base_functions; ++rr) {
2760 ++t_row_base_fun;
2761 ++t_row_grad_fun;
2762 }
2763 ++t_w;
2764 ++t_levi_kirchhoff_domega;
2765 }
2767}
constexpr double a
Kronecker Delta class.
#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()
constexpr auto t_kd
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.
FTensor::Tensor0< FTensor::PackPtr< double *, 1 > > getFTensor0N(const FieldApproximationBase base)
Get base function as Tensor0.
MatrixDouble & getN(const FieldApproximationBase base)
get base functions this return matrix (nb. of rows is equal to nb. of Gauss pts, nb....
const VectorInt & getIndices() const
Get global indices of dofs on entity.
auto getFTensor0IntegrationWeight()
Get integration weights.
MatrixDouble & getGaussPts()
matrix of integration (Gauss) points for Volume Element

Member Data Documentation

◆ alphaOmega

double EshelbianPlasticity::OpSpatialRotation_domega_domega::alphaOmega
private
Examples
EshelbianOperators.cpp.

Definition at line 810 of file EshelbianPlasticity.hpp.


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