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

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

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

Public Member Functions

 OpSpatialPhysical_du_dP (std::string row_field, std::string col_field, boost::shared_ptr< DataAtIntegrationPts > data_ptr, const bool assemble_off=false)
 
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 729 of file EshelbianPlasticity.hpp.

Constructor & Destructor Documentation

◆ OpSpatialPhysical_du_dP()

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

Definition at line 730 of file EshelbianPlasticity.hpp.

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

Reimplemented from EshelbianPlasticity::OpAssembleBasic< VolUserDataOperator >.

Examples
EshelbianOperators.cpp.

Definition at line 2359 of file EshelbianOperators.cpp.

2360 {
2362
2368
2369 int nb_integration_pts = row_data.getN().size1();
2370 int row_nb_dofs = row_data.getIndices().size();
2371 int col_nb_dofs = col_data.getIndices().size();
2372 auto get_ftensor3 = [](MatrixDouble &m, const int r, const int c) {
2374
2375 &m(r + 0, c + 0), &m(r + 0, c + 1), &m(r + 0, c + 2),
2376
2377 &m(r + 1, c + 0), &m(r + 1, c + 1), &m(r + 1, c + 2),
2378
2379 &m(r + 2, c + 0), &m(r + 2, c + 1), &m(r + 2, c + 2),
2380
2381 &m(r + 3, c + 0), &m(r + 3, c + 1), &m(r + 3, c + 2),
2382
2383 &m(r + 4, c + 0), &m(r + 4, c + 1), &m(r + 4, c + 2),
2384
2385 &m(r + 5, c + 0), &m(r + 5, c + 1), &m(r + 5, c + 2));
2386 };
2387
2388 auto v = getVolume();
2389 auto t_w = getFTensor0IntegrationWeight();
2390
2391 int row_nb_base_functions = row_data.getN().size2();
2392 auto t_row_base_fun = row_data.getFTensor0N();
2393
2394 auto t_approx_P_adjoint_log_du_dP =
2396
2397 for (int gg = 0; gg != nb_integration_pts; ++gg) {
2398 double a = v * t_w;
2399 int rr = 0;
2400 for (; rr != row_nb_dofs / 6; ++rr) {
2401
2402 auto t_col_base_fun = col_data.getFTensor1N<3>(gg, 0);
2403 auto t_m = get_ftensor3(K, 6 * rr, 0);
2404
2405 for (int cc = 0; cc != col_nb_dofs / 3; ++cc) {
2406 t_m(L, i) -=
2407 a * (t_approx_P_adjoint_log_du_dP(i, j, L) * t_col_base_fun(j)) *
2408 t_row_base_fun;
2409 ++t_col_base_fun;
2410 ++t_m;
2411 }
2412
2413 ++t_row_base_fun;
2414 }
2415 for (; rr != row_nb_base_functions; ++rr)
2416 ++t_row_base_fun;
2417 ++t_w;
2418 ++t_approx_P_adjoint_log_du_dP;
2419 }
2420
2422}
#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)
FTensor::Index< 'l', 3 > l
FTensor::Index< 'j', 3 > j
FTensor::Index< 'k', 3 > k
static constexpr auto size_symm
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: