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

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

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

Public Member Functions

 OpSpatialConsistency_dP_domega (std::string row_field, std::string col_field, boost::shared_ptr< DataAtIntegrationPts > data_ptr, const bool assemble_off)
 
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
 

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

Constructor & Destructor Documentation

◆ OpSpatialConsistency_dP_domega()

EshelbianPlasticity::OpSpatialConsistency_dP_domega::OpSpatialConsistency_dP_domega ( std::string row_field,
std::string col_field,
boost::shared_ptr< DataAtIntegrationPts > data_ptr,
const bool assemble_off )
inline

Definition at line 855 of file EshelbianPlasticity.hpp.

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

Reimplemented from EshelbianPlasticity::OpAssembleBasic< VolUserDataOperator >.

Examples
EshelbianOperators.cpp.

Definition at line 2969 of file EshelbianOperators.cpp.

2970 {
2972
2979
2980 int nb_integration_pts = row_data.getN().size1();
2981 int row_nb_dofs = row_data.getIndices().size();
2982 int col_nb_dofs = col_data.getIndices().size();
2983
2984 auto get_ftensor2 = [](MatrixDouble &m, const int r, const int c) {
2986
2987 &m(r + 0, c + 0), &m(r + 0, c + 1), &m(r + 0, c + 2),
2988
2989 &m(r + 1, c + 0), &m(r + 1, c + 1), &m(r + 1, c + 2),
2990
2991 &m(r + 2, c + 0), &m(r + 2, c + 1), &m(r + 2, c + 2)
2992
2993 );
2994 };
2995
2996 auto v = getVolume();
2997 auto t_w = getFTensor0IntegrationWeight();
2998 int row_nb_base_functions = row_data.getN().size2() / 3;
2999 auto t_row_base_fun = row_data.getFTensor1N<3>();
3000
3001 auto t_h_domega = getFTensor3FromMat<3, 3, 3>(dataAtPts->hdOmegaAtPts);
3002
3003 for (int gg = 0; gg != nb_integration_pts; ++gg) {
3004 double a = v * t_w;
3005
3006 int rr = 0;
3007 for (; rr != row_nb_dofs / 3; ++rr) {
3008
3010 t_PRT(i, k) = t_row_base_fun(j) * t_h_domega(i, j, k);
3011
3012 auto t_col_base_fun = col_data.getFTensor0N(gg, 0);
3013 auto t_m = get_ftensor2(K, 3 * rr, 0);
3014 for (int cc = 0; cc != col_nb_dofs / 3; ++cc) {
3015 t_m(i, j) -= (a * t_col_base_fun) * t_PRT(i, j);
3016 ++t_m;
3017 ++t_col_base_fun;
3018 }
3019
3020 ++t_row_base_fun;
3021 }
3022
3023 for (; rr != row_nb_base_functions; ++rr)
3024 ++t_row_base_fun;
3025 ++t_w;
3026 ++t_h_domega;
3027 }
3029}
#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 FTensor::Tensor3< FTensor::PackPtr< T *, 1 >, Tensor_Dim0, Tensor_Dim1, Tensor_Dim2 > getFTensor3FromMat(ublas::matrix< T, L, A > &data)
Get tensor rank 3 (non symmetries) form data matrix.
int r
Definition sdf.py:8
FTensor::Index< 'm', 3 > m
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....
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.

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