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

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

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

Public Member Functions

 OpSpatialConsistency_dBubble_dP (std::string row_field, std::string col_field, boost::shared_ptr< DataAtIntegrationPts > data_ptr)
 
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 Member Functions

template<int S>
MoFEMErrorCode integrateImpl (EntData &row_data, EntData &col_data)
 

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 840 of file EshelbianPlasticity.hpp.

Constructor & Destructor Documentation

◆ OpSpatialConsistency_dBubble_dP()

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

Definition at line 841 of file EshelbianPlasticity.hpp.

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

Reimplemented from EshelbianPlasticity::OpAssembleBasic< VolUserDataOperator >.

Examples
EshelbianOperators.cpp.

Definition at line 2899 of file EshelbianOperators.cpp.

2900 {
2902 if (dataAtPts->matInvD.size1() == size_symm &&
2903 dataAtPts->matInvD.size2() == size_symm) {
2904 MoFEMFunctionReturnHot(integrateImpl<0>(row_data, col_data));
2905 } else {
2907 }
2909};
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
static constexpr auto size_symm
MoFEMErrorCode integrateImpl(EntData &row_data, EntData &col_data)

◆ integrateImpl()

template<int S>
MoFEMErrorCode OpSpatialConsistency_dBubble_dP::integrateImpl ( EntData & row_data,
EntData & col_data )
private
Examples
EshelbianOperators.cpp.

Definition at line 2913 of file EshelbianOperators.cpp.

2914 {
2916
2917 auto get_ftensor1 = [](MatrixDouble &m, const int r, const int c) {
2919
2920 &m(r + 0, c + 0), &m(r + 0, c + 1), &m(r + 0, c + 2)
2921
2922 );
2923 };
2924
2925 int nb_integration_pts = getGaussPts().size2();
2926 int row_nb_dofs = row_data.getIndices().size();
2927 int col_nb_dofs = col_data.getIndices().size();
2928
2929 auto v = getVolume();
2930 auto t_w = getFTensor0IntegrationWeight();
2931 int row_nb_base_functions = row_data.getN().size2() / 9;
2932
2939
2941 &*dataAtPts->matInvD.data().begin());
2942
2943 auto t_row_base = row_data.getFTensor2N<3, 3>();
2944 for (int gg = 0; gg != nb_integration_pts; ++gg) {
2945 double a = v * t_w;
2946
2947 auto t_m = get_ftensor1(K, 0, 0);
2948
2949 int rr = 0;
2950 for (; rr != row_nb_dofs; ++rr) {
2951 auto t_col_base = col_data.getFTensor1N<3>(gg, 0);
2952 for (int cc = 0; cc != col_nb_dofs / 3; ++cc) {
2953 t_m(k) -= a * (t_row_base(i, j) * t_inv_D(i, j, k, l)) * t_col_base(l);
2954 ++t_col_base;
2955 ++t_m;
2956 }
2957
2958 ++t_row_base;
2959 }
2960
2961 for (; rr != row_nb_base_functions; ++rr)
2962 ++t_row_base;
2963 ++t_w;
2964 ++t_inv_D;
2965 }
2967}
#define FTENSOR_INDEX(DIM, I)
constexpr double a
constexpr int SPACE_DIM
#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
static auto getFTensor4DdgFromPtr(T *ptr)
int r
Definition sdf.py:8
FTensor::Index< 'm', 3 > m
FTensor::Tensor2< FTensor::PackPtr< double *, Tensor_Dim0 *Tensor_Dim1 >, Tensor_Dim0, Tensor_Dim1 > getFTensor2N(FieldApproximationBase base)
Get base functions for Hdiv/Hcurl spaces.
MatrixDouble & getN(const FieldApproximationBase base)
get base functions this return matrix (nb. of rows is equal to nb. of Gauss pts, nb....
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.
auto getFTensor0IntegrationWeight()
Get integration weights.
MatrixDouble & getGaussPts()
matrix of integration (Gauss) points for Volume Element

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