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

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

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

Public Member Functions

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

Constructor & Destructor Documentation

◆ OpSpatialEquilibrium_dw_dP()

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

Definition at line 706 of file EshelbianPlasticity.hpp.

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

Reimplemented from EshelbianPlasticity::OpAssembleBasic< VolUserDataOperator >.

Examples
EshelbianOperators.cpp.

Definition at line 2264 of file EshelbianOperators.cpp.

2265 {
2267 int nb_integration_pts = row_data.getN().size1();
2268 int row_nb_dofs = row_data.getIndices().size();
2269 int col_nb_dofs = col_data.getIndices().size();
2270 auto get_ftensor1 = [](MatrixDouble &m, const int r, const int c) {
2272 &m(r + 0, c + 0), &m(r + 1, c + 1), &m(r + 2, c + 2));
2273 };
2274 FTensor::Index<'i', 3> i;
2275 auto v = getVolume();
2276 auto t_w = getFTensor0IntegrationWeight();
2277 int row_nb_base_functions = row_data.getN().size2();
2278 auto t_row_base_fun = row_data.getFTensor0N();
2279 for (int gg = 0; gg != nb_integration_pts; ++gg) {
2280 double a = v * t_w;
2281 int rr = 0;
2282 for (; rr != row_nb_dofs / 3; ++rr) {
2283 auto t_col_diff_base_fun = col_data.getFTensor2DiffN<3, 3>(gg, 0);
2284 auto t_m = get_ftensor1(K, 3 * rr, 0);
2285 for (int cc = 0; cc != col_nb_dofs / 3; ++cc) {
2286 double div_col_base = t_col_diff_base_fun(i, i);
2287 t_m(i) -= a * t_row_base_fun * div_col_base;
2288 ++t_m;
2289 ++t_col_diff_base_fun;
2290 }
2291 ++t_row_base_fun;
2292 }
2293 for (; rr != row_nb_base_functions; ++rr)
2294 ++t_row_base_fun;
2295 ++t_w;
2296 }
2298}
constexpr double a
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
FTensor::Index< 'i', SPACE_DIM > i
const double c
speed of light (cm/ns)
int r
Definition sdf.py:8
FTensor::Index< 'm', 3 > m
FTensor::Tensor2< FTensor::PackPtr< double *, Tensor_Dim0 *Tensor_Dim1 >, Tensor_Dim0, Tensor_Dim1 > getFTensor2DiffN(FieldApproximationBase base)
Get derivatives of base functions for Hdiv space.
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: