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

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

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

Public Member Functions

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

Constructor & Destructor Documentation

◆ OpSpatialPhysical_du_domega()

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

Definition at line 755 of file EshelbianPlasticity.hpp.

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

Reimplemented from EshelbianPlasticity::OpAssembleBasic< VolUserDataOperator >.

Examples
EshelbianOperators.cpp.

Definition at line 2475 of file EshelbianOperators.cpp.

2476 {
2478
2480 auto t_L = symm_L_tensor();
2481
2482 int row_nb_dofs = row_data.getIndices().size();
2483 int col_nb_dofs = col_data.getIndices().size();
2484 auto get_ftensor3 = [](MatrixDouble &m, const int r, const int c) {
2486
2487 &m(r + 0, c + 0), &m(r + 0, c + 1), &m(r + 0, c + 2),
2488
2489 &m(r + 1, c + 0), &m(r + 1, c + 1), &m(r + 1, c + 2),
2490
2491 &m(r + 2, c + 0), &m(r + 2, c + 1), &m(r + 2, c + 2),
2492
2493 &m(r + 3, c + 0), &m(r + 3, c + 1), &m(r + 3, c + 2),
2494
2495 &m(r + 4, c + 0), &m(r + 4, c + 1), &m(r + 4, c + 2),
2496
2497 &m(r + 5, c + 0), &m(r + 5, c + 1), &m(r + 5, c + 2)
2498
2499 );
2500 };
2501 FTensor::Index<'i', 3> i;
2502 FTensor::Index<'j', 3> j;
2503 FTensor::Index<'k', 3> k;
2504 FTensor::Index<'m', 3> m;
2505 FTensor::Index<'n', 3> n;
2506
2507 auto v = getVolume();
2508 auto t_w = getFTensor0IntegrationWeight();
2509 auto t_approx_P_adjoint_log_du_domega =
2510 getFTensor2FromMat<3, size_symm>(dataAtPts->adjointPdUdOmegaAtPts);
2511
2512 int row_nb_base_functions = row_data.getN().size2();
2513 auto t_row_base_fun = row_data.getFTensor0N();
2514
2515 int nb_integration_pts = row_data.getN().size1();
2516
2517 for (int gg = 0; gg != nb_integration_pts; ++gg) {
2518 double a = v * t_w;
2519
2520 int rr = 0;
2521 for (; rr != row_nb_dofs / 6; ++rr) {
2522 auto t_col_base_fun = col_data.getFTensor0N(gg, 0);
2523 auto t_m = get_ftensor3(K, 6 * rr, 0);
2524 for (int cc = 0; cc != col_nb_dofs / 3; ++cc) {
2525 double v = a * t_row_base_fun * t_col_base_fun;
2526 t_m(L, k) -= v * t_approx_P_adjoint_log_du_domega(k, L);
2527 ++t_m;
2528 ++t_col_base_fun;
2529 }
2530 ++t_row_base_fun;
2531 }
2532
2533 for (; rr != row_nb_base_functions; ++rr)
2534 ++t_row_base_fun;
2535
2536 ++t_w;
2537 ++t_approx_P_adjoint_log_du_domega;
2538 }
2539
2541}
constexpr double a
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
FTensor::Index< 'i', SPACE_DIM > i
const double c
speed of light (cm/ns)
const double n
refractive index of diffusive medium
FTensor::Index< 'j', 3 > j
FTensor::Index< 'k', 3 > k
static constexpr auto size_symm
FTensor::Tensor2< FTensor::PackPtr< double *, 1 >, Tensor_Dim1, Tensor_Dim2 > getFTensor2FromMat(MatrixDouble &data)
Get tensor rank 2 (matrix) form data matrix.
int r
Definition sdf.py:8
FTensor::Index< 'm', 3 > m
FTensor::Tensor0< FTensor::PackPtr< double *, 1 > > getFTensor0N(const FieldApproximationBase base)
Get base function as Tensor0.
MatrixDouble & getN(const FieldApproximationBase base)
get base functions this return matrix (nb. of rows is equal to nb. of Gauss pts, nb....
const VectorInt & getIndices() const
Get global indices of dofs on entity.
auto getFTensor0IntegrationWeight()
Get integration weights.

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