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

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

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

Public Member Functions

 OpSpatialConsistency_dBubble_dBubble (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 826 of file EshelbianPlasticity.hpp.

Constructor & Destructor Documentation

◆ OpSpatialConsistency_dBubble_dBubble()

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

Definition at line 827 of file EshelbianPlasticity.hpp.

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

Reimplemented from EshelbianPlasticity::OpAssembleBasic< VolUserDataOperator >.

Examples
EshelbianOperators.cpp.

Definition at line 2841 of file EshelbianOperators.cpp.

2842 {
2844 if (dataAtPts->matInvD.size1() == size_symm &&
2845 dataAtPts->matInvD.size2() == size_symm) {
2846 MoFEMFunctionReturnHot(integrateImpl<0>(row_data, col_data));
2847 } else {
2849 }
2851};
#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_dBubble::integrateImpl ( EntData & row_data,
EntData & col_data )
private
Examples
EshelbianOperators.cpp.

Definition at line 2855 of file EshelbianOperators.cpp.

2856 {
2858
2859 int nb_integration_pts = getGaussPts().size2();
2860 int row_nb_dofs = row_data.getIndices().size();
2861 int col_nb_dofs = col_data.getIndices().size();
2862
2863 auto v = getVolume();
2864 auto t_w = getFTensor0IntegrationWeight();
2865 int row_nb_base_functions = row_data.getN().size2() / 9;
2866
2871
2873 &*dataAtPts->matInvD.data().begin());
2874
2875 auto t_row_base = row_data.getFTensor2N<3, 3>();
2876 for (int gg = 0; gg != nb_integration_pts; ++gg) {
2877 double a = v * t_w;
2878
2879 int rr = 0;
2880 for (; rr != row_nb_dofs; ++rr) {
2881 auto t_col_base = col_data.getFTensor2N<3, 3>(gg, 0);
2882 for (int cc = 0; cc != col_nb_dofs; ++cc) {
2883 K(rr, cc) -=
2884 a * (t_row_base(i, j) * (t_inv_D(i, j, k, l) * t_col_base(k, l)));
2885 ++t_col_base;
2886 }
2887
2888 ++t_row_base;
2889 }
2890
2891 for (; rr != row_nb_base_functions; ++rr)
2892 ++t_row_base;
2893 ++t_w;
2894 ++t_inv_D;
2895 }
2897}
#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
FTensor::Index< 'l', 3 > l
FTensor::Index< 'j', 3 > j
FTensor::Index< 'k', 3 > k
static auto getFTensor4DdgFromPtr(T *ptr)
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....
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: