-
Notifications
You must be signed in to change notification settings - Fork 1
/
Cartesian2d_Cusp.hpp
163 lines (126 loc) · 3.56 KB
/
Cartesian2d_Cusp.hpp
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
/*
(c - MIT) T.W.J. de Geus (Tom) | www.geus.me | github.com/tdegeus/GMatElastoPlasticQPot
*/
#ifndef GMATELASTOPLASTICQPOT_CARTESIAN2D_CUSP_HPP
#define GMATELASTOPLASTICQPOT_CARTESIAN2D_CUSP_HPP
#include "Cartesian2d.h"
namespace GMatElastoPlasticQPot {
namespace Cartesian2d {
inline Cusp::Cusp(double K, double G, const xt::xtensor<double, 1>& epsy, bool init_elastic)
: m_K(K), m_G(G)
{
xt::xtensor<double, 1> y = xt::sort(epsy);
if (init_elastic) {
if (y(0) != -y(1)) {
y = xt::concatenate(xt::xtuple(xt::xtensor<double, 1>({-y(0)}), y));
}
}
GMATELASTOPLASTICQPOT_ASSERT(y.size() > 1);
m_yield = QPot::Static(0.0, y);
}
inline double Cusp::K() const
{
return m_K;
}
inline double Cusp::G() const
{
return m_G;
}
inline xt::xtensor<double, 1> Cusp::epsy() const
{
return m_yield.yield();
}
inline size_t Cusp::currentIndex() const
{
return m_yield.currentIndex();
}
inline double Cusp::currentYieldLeft() const
{
return m_yield.currentYieldLeft();
}
inline double Cusp::currentYieldRight() const
{
return m_yield.currentYieldRight();
}
inline double Cusp::epsp() const
{
return 0.5 * (m_yield.currentYieldLeft() + m_yield.currentYieldRight());
}
template <class T>
inline void Cusp::setStrain(const T& a)
{
GMATELASTOPLASTICQPOT_ASSERT(xt::has_shape(a, {2, 2}));
return this->setStrainIterator(a.cbegin());
}
template <class T>
inline void Cusp::setStrainIterator(const T& begin)
{
std::copy(begin, begin + 4, m_Eps.begin());
std::array<double, 4> Epsd;
double epsm = detail::hydrostatic_deviator(m_Eps, Epsd);
double epsd = std::sqrt(0.5 * detail::A2_ddot_B2(Epsd, Epsd));
m_yield.setPosition(epsd);
m_Sig[0] = m_Sig[3] = m_K * epsm;
if (epsd <= 0.0) {
m_Sig[1] = m_Sig[2] = 0.0;
return;
}
double eps_min = 0.5 * (m_yield.currentYieldRight() + m_yield.currentYieldLeft());
double g = m_G * (1.0 - eps_min / epsd);
m_Sig[0] += g * Epsd[0];
m_Sig[1] = g * Epsd[1];
m_Sig[2] = g * Epsd[2];
m_Sig[3] += g * Epsd[3];
}
template <class T>
inline void Cusp::stress(T& a) const
{
GMATELASTOPLASTICQPOT_ASSERT(xt::has_shape(a, {2, 2}));
return this->stressIterator(a.begin());
}
template <class T>
inline void Cusp::stressIterator(const T& begin) const
{
std::copy(m_Sig.begin(), m_Sig.end(), begin);
}
inline Tensor2 Cusp::Stress() const
{
auto ret = Tensor2::from_shape({2, 2});
this->stressIterator(ret.begin());
return ret;
}
template <class T>
inline void Cusp::tangent(T& C) const
{
auto II = Cartesian2d::II();
auto I4d = Cartesian2d::I4d();
xt::noalias(C) = 0.5 * m_K * II + m_G * I4d;
}
inline Tensor4 Cusp::Tangent() const
{
auto ret = Tensor4::from_shape({2, 2, 2, 2});
this->tangent(ret);
return ret;
}
inline double Cusp::energy() const
{
std::array<double, 4> Epsd;
double epsm = detail::hydrostatic_deviator(m_Eps, Epsd);
double epsd = std::sqrt(0.5 * detail::A2_ddot_B2(Epsd, Epsd));
double U = m_K * std::pow(epsm, 2.0);
double eps_min = 0.5 * (m_yield.currentYieldRight() + m_yield.currentYieldLeft());
double deps_y = 0.5 * (m_yield.currentYieldRight() - m_yield.currentYieldLeft());
double V = m_G * (std::pow(epsd - eps_min, 2.0) - std::pow(deps_y, 2.0));
return U + V;
}
inline bool Cusp::checkYieldBoundLeft(size_t n) const
{
return m_yield.checkYieldBoundLeft(n);
}
inline bool Cusp::checkYieldBoundRight(size_t n) const
{
return m_yield.checkYieldBoundRight(n);
}
} // namespace Cartesian2d
} // namespace GMatElastoPlasticQPot
#endif