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

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

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

Public Member Functions

 OpSpatialPhysical_du_dBubble (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 741 of file EshelbianPlasticity.hpp.

Constructor & Destructor Documentation

◆ OpSpatialPhysical_du_dBubble()

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

Definition at line 742 of file EshelbianPlasticity.hpp.

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

Reimplemented from EshelbianPlasticity::OpAssembleBasic< VolUserDataOperator >.

Examples
EshelbianOperators.cpp.

Definition at line 2424 of file EshelbianOperators.cpp.

2425 {
2427
2433
2434 int nb_integration_pts = row_data.getN().size1();
2435 int row_nb_dofs = row_data.getIndices().size();
2436 int col_nb_dofs = col_data.getIndices().size();
2437 auto get_ftensor2 = [](MatrixDouble &m, const int r, const int c) {
2439 &m(r + 0, c), &m(r + 1, c), &m(r + 2, c), &m(r + 3, c), &m(r + 4, c),
2440 &m(r + 5, c));
2441 };
2442
2443 auto v = getVolume();
2444 auto t_w = getFTensor0IntegrationWeight();
2445 auto t_row_base_fun = row_data.getFTensor0N();
2446
2447 int row_nb_base_functions = row_data.getN().size2();
2448
2449 auto t_approx_P_adjoint_log_du_dP =
2451
2452 for (int gg = 0; gg != nb_integration_pts; ++gg) {
2453 double a = v * t_w;
2454 int rr = 0;
2455 for (; rr != row_nb_dofs / 6; ++rr) {
2456 auto t_m = get_ftensor2(K, 6 * rr, 0);
2457 auto t_col_base_fun = col_data.getFTensor2N<3, 3>(gg, 0);
2458 for (int cc = 0; cc != col_nb_dofs; ++cc) {
2459 t_m(L) -=
2460 a * (t_approx_P_adjoint_log_du_dP(i, j, L) * t_col_base_fun(i, j)) *
2461 t_row_base_fun;
2462 ++t_m;
2463 ++t_col_base_fun;
2464 }
2465 ++t_row_base_fun;
2466 }
2467 for (; rr != row_nb_base_functions; ++rr)
2468 ++t_row_base_fun;
2469 ++t_w;
2470 ++t_approx_P_adjoint_log_du_dP;
2471 }
2473}
#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::Tensor2< FTensor::PackPtr< double *, Tensor_Dim0 *Tensor_Dim1 >, Tensor_Dim0, Tensor_Dim1 > getFTensor2N(FieldApproximationBase base)
Get base functions for Hdiv/Hcurl spaces.
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: