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

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

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

Public Member Functions

 OpSpatialEquilibrium_dw_dw (std::string row_field, std::string col_field, boost::shared_ptr< DataAtIntegrationPts > data_ptr, const double alpha, const double rho)
 
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
 

Public Attributes

const double alphaW
 
const double alphaRho
 
- 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
 

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

Constructor & Destructor Documentation

◆ OpSpatialEquilibrium_dw_dw()

EshelbianPlasticity::OpSpatialEquilibrium_dw_dw::OpSpatialEquilibrium_dw_dw ( std::string row_field,
std::string col_field,
boost::shared_ptr< DataAtIntegrationPts > data_ptr,
const double alpha,
const double rho )
inline

Definition at line 719 of file EshelbianPlasticity.hpp.

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

Reimplemented from EshelbianPlasticity::OpAssembleBasic< VolUserDataOperator >.

Examples
EshelbianOperators.cpp.

Definition at line 2300 of file EshelbianOperators.cpp.

2301 {
2303
2304 if (alphaW < std::numeric_limits<double>::epsilon() &&
2305 alphaRho < std::numeric_limits<double>::epsilon())
2307
2308 const int nb_integration_pts = row_data.getN().size1();
2309 const int row_nb_dofs = row_data.getIndices().size();
2310 auto get_ftensor1 = [](MatrixDouble &m, const int r, const int c) {
2312 &m(r + 0, c + 0), &m(r + 1, c + 1), &m(r + 2, c + 2)
2313
2314 );
2315 };
2316 FTensor::Index<'i', 3> i;
2317
2318 auto v = getVolume();
2319 auto t_w = getFTensor0IntegrationWeight();
2320
2321 auto piola_scale = dataAtPts->piolaScale;
2322 auto alpha_w = alphaW / piola_scale;
2323 auto alpha_rho = alphaRho / piola_scale;
2324
2325 int row_nb_base_functions = row_data.getN().size2();
2326 auto t_row_base_fun = row_data.getFTensor0N();
2327
2328 double ts_scale = alpha_w * getTSa();
2329 if (std::abs(alphaRho) > std::numeric_limits<double>::epsilon())
2330 ts_scale += alpha_rho * getTSaa();
2331
2332 for (int gg = 0; gg != nb_integration_pts; ++gg) {
2333 double a = v * t_w * ts_scale;
2334
2335 int rr = 0;
2336 for (; rr != row_nb_dofs / 3; ++rr) {
2337
2338 auto t_col_base_fun = row_data.getFTensor0N(gg, 0);
2339 auto t_m = get_ftensor1(K, 3 * rr, 0);
2340 for (int cc = 0; cc != row_nb_dofs / 3; ++cc) {
2341 const double b = a * t_row_base_fun * t_col_base_fun;
2342 t_m(i) += b;
2343 ++t_m;
2344 ++t_col_base_fun;
2345 }
2346
2347 ++t_row_base_fun;
2348 }
2349
2350 for (; rr != row_nb_base_functions; ++rr)
2351 ++t_row_base_fun;
2352
2353 ++t_w;
2354 }
2355
2357}
constexpr double a
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
#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)
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.

Member Data Documentation

◆ alphaRho

const double EshelbianPlasticity::OpSpatialEquilibrium_dw_dw::alphaRho
Examples
EshelbianOperators.cpp.

Definition at line 718 of file EshelbianPlasticity.hpp.

◆ alphaW

const double EshelbianPlasticity::OpSpatialEquilibrium_dw_dw::alphaW
Examples
EshelbianOperators.cpp.

Definition at line 717 of file EshelbianPlasticity.hpp.


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