Page Menu
Home
c4science
Search
Configure Global Search
Log In
Files
F86511054
ElementQuad4Planar.hpp
No One
Temporary
Actions
Download File
Edit File
Delete File
View Transforms
Subscribe
Mute Notifications
Award Token
Subscribers
None
File Metadata
Details
File Info
Storage
Attached
Created
Sun, Oct 6, 22:37
Size
9 KB
Mime Type
text/x-c
Expires
Tue, Oct 8, 22:37 (1 d, 23 h)
Engine
blob
Format
Raw Data
Handle
21399206
Attached To
rGOOSEFEM GooseFEM
ElementQuad4Planar.hpp
View Options
/**
Implementation of ElementQuad4Planar.h
\file ElementQuad4Planar.hpp
\copyright Copyright 2017. Tom de Geus. All rights reserved.
\license This project is released under the GNU Public License (GPLv3).
*/
#ifndef GOOSEFEM_ELEMENTQUAD4PLANAR_HPP
#define GOOSEFEM_ELEMENTQUAD4PLANAR_HPP
#include "ElementQuad4Planar.h"
#include "detail.hpp"
namespace GooseFEM {
namespace Element {
namespace Quad4 {
inline QuadraturePlanar::QuadraturePlanar(const xt::xtensor<double, 3>& x, double thick)
: QuadraturePlanar(x, Gauss::xi(), Gauss::w(), thick)
{
}
inline QuadraturePlanar::QuadraturePlanar(
const xt::xtensor<double, 3>& x,
const xt::xtensor<double, 2>& xi,
const xt::xtensor<double, 1>& w,
double thick)
{
m_thick = thick;
size_t nip = w.size();
xt::xtensor<double, 2> N = xt::empty<double>({nip, m_nne});
xt::xtensor<double, 3> dNxi = xt::empty<double>({nip, m_nne, m_ndim});
for (size_t q = 0; q < nip; ++q) {
N(q, 0) = 0.25 * (1.0 - xi(q, 0)) * (1.0 - xi(q, 1));
N(q, 1) = 0.25 * (1.0 + xi(q, 0)) * (1.0 - xi(q, 1));
N(q, 2) = 0.25 * (1.0 + xi(q, 0)) * (1.0 + xi(q, 1));
N(q, 3) = 0.25 * (1.0 - xi(q, 0)) * (1.0 + xi(q, 1));
}
for (size_t q = 0; q < nip; ++q) {
// - dN / dxi_0
dNxi(q, 0, 0) = -0.25 * (1.0 - xi(q, 1));
dNxi(q, 1, 0) = +0.25 * (1.0 - xi(q, 1));
dNxi(q, 2, 0) = +0.25 * (1.0 + xi(q, 1));
dNxi(q, 3, 0) = -0.25 * (1.0 + xi(q, 1));
// - dN / dxi_1
dNxi(q, 0, 1) = -0.25 * (1.0 - xi(q, 0));
dNxi(q, 1, 1) = -0.25 * (1.0 + xi(q, 0));
dNxi(q, 2, 1) = +0.25 * (1.0 + xi(q, 0));
dNxi(q, 3, 1) = +0.25 * (1.0 - xi(q, 0));
}
this->initQuadratureBaseCartesian(x, xi, w, N, dNxi);
}
inline void QuadraturePlanar::compute_dN()
{
#pragma omp parallel
{
xt::xtensor<double, 2> J = xt::empty<double>({2, 2});
xt::xtensor<double, 2> Jinv = xt::empty<double>({2, 2});
#pragma omp for
for (size_t e = 0; e < m_nelem; ++e) {
auto x = xt::adapt(&m_x(e, 0, 0), xt::xshape<m_nne, m_ndim>());
for (size_t q = 0; q < m_nip; ++q) {
auto dNxi = xt::adapt(&m_dNxi(q, 0, 0), xt::xshape<m_nne, m_ndim>());
auto dNx = xt::adapt(&m_dNx(e, q, 0, 0), xt::xshape<m_nne, m_ndim>());
// J(i,j) += dNxi(m,i) * x(m,j);
J(0, 0) = dNxi(0, 0) * x(0, 0) + dNxi(1, 0) * x(1, 0) + dNxi(2, 0) * x(2, 0) +
dNxi(3, 0) * x(3, 0);
J(0, 1) = dNxi(0, 0) * x(0, 1) + dNxi(1, 0) * x(1, 1) + dNxi(2, 0) * x(2, 1) +
dNxi(3, 0) * x(3, 1);
J(1, 0) = dNxi(0, 1) * x(0, 0) + dNxi(1, 1) * x(1, 0) + dNxi(2, 1) * x(2, 0) +
dNxi(3, 1) * x(3, 0);
J(1, 1) = dNxi(0, 1) * x(0, 1) + dNxi(1, 1) * x(1, 1) + dNxi(2, 1) * x(2, 1) +
dNxi(3, 1) * x(3, 1);
double Jdet = detail::tensor<2>::inv(J, Jinv);
// dNx(m,i) += Jinv(i,j) * dNxi(m,j);
for (size_t m = 0; m < m_nne; ++m) {
dNx(m, 0) = Jinv(0, 0) * dNxi(m, 0) + Jinv(0, 1) * dNxi(m, 1);
dNx(m, 1) = Jinv(1, 0) * dNxi(m, 0) + Jinv(1, 1) * dNxi(m, 1);
}
m_vol(e, q) = m_w(q) * Jdet * m_thick;
}
}
}
}
inline void QuadraturePlanar::gradN_vector(
const xt::xtensor<double, 3>& elemvec, xt::xtensor<double, 4>& qtensor) const
{
GOOSEFEM_ASSERT(xt::has_shape(elemvec, {m_nelem, m_nne, m_ndim}));
GOOSEFEM_ASSERT(xt::has_shape(qtensor, {m_nelem, m_nip, m_tdim, m_tdim}));
qtensor.fill(0.0);
#pragma omp parallel for
for (size_t e = 0; e < m_nelem; ++e) {
auto u = xt::adapt(&elemvec(e, 0, 0), xt::xshape<m_nne, m_ndim>());
for (size_t q = 0; q < m_nip; ++q) {
auto dNx = xt::adapt(&m_dNx(e, q, 0, 0), xt::xshape<m_nne, m_ndim>());
auto gradu = xt::adapt(&qtensor(e, q, 0, 0), xt::xshape<m_tdim, m_tdim>());
// gradu(i,j) += dNx(m,i) * u(m,j)
gradu(0, 0) = dNx(0, 0) * u(0, 0) + dNx(1, 0) * u(1, 0) + dNx(2, 0) * u(2, 0) +
dNx(3, 0) * u(3, 0);
gradu(0, 1) = dNx(0, 0) * u(0, 1) + dNx(1, 0) * u(1, 1) + dNx(2, 0) * u(2, 1) +
dNx(3, 0) * u(3, 1);
gradu(1, 0) = dNx(0, 1) * u(0, 0) + dNx(1, 1) * u(1, 0) + dNx(2, 1) * u(2, 0) +
dNx(3, 1) * u(3, 0);
gradu(1, 1) = dNx(0, 1) * u(0, 1) + dNx(1, 1) * u(1, 1) + dNx(2, 1) * u(2, 1) +
dNx(3, 1) * u(3, 1);
}
}
}
inline void QuadraturePlanar::gradN_vector_T(
const xt::xtensor<double, 3>& elemvec, xt::xtensor<double, 4>& qtensor) const
{
GOOSEFEM_ASSERT(xt::has_shape(elemvec, {m_nelem, m_nne, m_ndim}));
GOOSEFEM_ASSERT(xt::has_shape(qtensor, {m_nelem, m_nip, m_tdim, m_tdim}));
qtensor.fill(0.0);
#pragma omp parallel for
for (size_t e = 0; e < m_nelem; ++e) {
auto u = xt::adapt(&elemvec(e, 0, 0), xt::xshape<m_nne, m_ndim>());
for (size_t q = 0; q < m_nip; ++q) {
auto dNx = xt::adapt(&m_dNx(e, q, 0, 0), xt::xshape<m_nne, m_ndim>());
auto gradu = xt::adapt(&qtensor(e, q, 0, 0), xt::xshape<m_tdim, m_tdim>());
// gradu(j,i) += dNx(m,i) * u(m,j)
gradu(0, 0) = dNx(0, 0) * u(0, 0) + dNx(1, 0) * u(1, 0) + dNx(2, 0) * u(2, 0) +
dNx(3, 0) * u(3, 0);
gradu(1, 0) = dNx(0, 0) * u(0, 1) + dNx(1, 0) * u(1, 1) + dNx(2, 0) * u(2, 1) +
dNx(3, 0) * u(3, 1);
gradu(0, 1) = dNx(0, 1) * u(0, 0) + dNx(1, 1) * u(1, 0) + dNx(2, 1) * u(2, 0) +
dNx(3, 1) * u(3, 0);
gradu(1, 1) = dNx(0, 1) * u(0, 1) + dNx(1, 1) * u(1, 1) + dNx(2, 1) * u(2, 1) +
dNx(3, 1) * u(3, 1);
}
}
}
inline void QuadraturePlanar::symGradN_vector(
const xt::xtensor<double, 3>& elemvec, xt::xtensor<double, 4>& qtensor) const
{
GOOSEFEM_ASSERT(xt::has_shape(elemvec, {m_nelem, m_nne, m_ndim}));
GOOSEFEM_ASSERT(xt::has_shape(qtensor, {m_nelem, m_nip, m_tdim, m_tdim}));
qtensor.fill(0.0);
#pragma omp parallel for
for (size_t e = 0; e < m_nelem; ++e) {
auto u = xt::adapt(&elemvec(e, 0, 0), xt::xshape<m_nne, m_ndim>());
for (size_t q = 0; q < m_nip; ++q) {
auto dNx = xt::adapt(&m_dNx(e, q, 0, 0), xt::xshape<m_nne, m_ndim>());
auto eps = xt::adapt(&qtensor(e, q, 0, 0), xt::xshape<m_tdim, m_tdim>());
// gradu(i,j) += dNx(m,i) * u(m,j)
// eps(j,i) = 0.5 * (gradu(i,j) + gradu(j,i))
eps(0, 0) = dNx(0, 0) * u(0, 0) + dNx(1, 0) * u(1, 0) + dNx(2, 0) * u(2, 0) +
dNx(3, 0) * u(3, 0);
eps(1, 1) = dNx(0, 1) * u(0, 1) + dNx(1, 1) * u(1, 1) + dNx(2, 1) * u(2, 1) +
dNx(3, 1) * u(3, 1);
eps(0, 1) = 0.5 * (dNx(0, 0) * u(0, 1) + dNx(1, 0) * u(1, 1) + dNx(2, 0) * u(2, 1) +
dNx(3, 0) * u(3, 1) + dNx(0, 1) * u(0, 0) + dNx(1, 1) * u(1, 0) +
dNx(2, 1) * u(2, 0) + dNx(3, 1) * u(3, 0));
eps(1, 0) = eps(0, 1);
}
}
}
inline void QuadraturePlanar::int_N_scalar_NT_dV(
const xt::xtensor<double, 2>& qscalar, xt::xtensor<double, 3>& elemmat) const
{
GOOSEFEM_ASSERT(xt::has_shape(qscalar, {m_nelem, m_nip}));
GOOSEFEM_ASSERT(xt::has_shape(elemmat, {m_nelem, m_nne * m_ndim, m_nne * m_ndim}));
elemmat.fill(0.0);
#pragma omp parallel for
for (size_t e = 0; e < m_nelem; ++e) {
auto M = xt::adapt(&elemmat(e, 0, 0), xt::xshape<m_nne * m_ndim, m_nne * m_ndim>());
for (size_t q = 0; q < m_nip; ++q) {
auto N = xt::adapt(&m_N(q, 0), xt::xshape<m_nne>());
auto& vol = m_vol(e, q);
auto& rho = qscalar(e, q);
// M(m*ndim+i,n*ndim+i) += N(m) * scalar * N(n) * dV
for (size_t m = 0; m < m_nne; ++m) {
for (size_t n = 0; n < m_nne; ++n) {
M(m * m_ndim + 0, n * m_ndim + 0) += N(m) * rho * N(n) * vol;
M(m * m_ndim + 1, n * m_ndim + 1) += N(m) * rho * N(n) * vol;
}
}
}
}
}
inline void QuadraturePlanar::int_gradN_dot_tensor2_dV(
const xt::xtensor<double, 4>& qtensor, xt::xtensor<double, 3>& elemvec) const
{
GOOSEFEM_ASSERT(xt::has_shape(qtensor, {m_nelem, m_nip, m_tdim, m_tdim}));
GOOSEFEM_ASSERT(xt::has_shape(elemvec, {m_nelem, m_nne, m_ndim}));
elemvec.fill(0.0);
#pragma omp parallel for
for (size_t e = 0; e < m_nelem; ++e) {
auto f = xt::adapt(&elemvec(e, 0, 0), xt::xshape<m_nne, m_ndim>());
for (size_t q = 0; q < m_nip; ++q) {
auto dNx = xt::adapt(&m_dNx(e, q, 0, 0), xt::xshape<m_nne, m_ndim>());
auto sig = xt::adapt(&qtensor(e, q, 0, 0), xt::xshape<m_tdim, m_tdim>());
auto& vol = m_vol(e, q);
for (size_t m = 0; m < m_nne; ++m) {
f(m, 0) += (dNx(m, 0) * sig(0, 0) + dNx(m, 1) * sig(1, 0)) * vol;
f(m, 1) += (dNx(m, 0) * sig(0, 1) + dNx(m, 1) * sig(1, 1)) * vol;
}
}
}
}
} // namespace Quad4
} // namespace Element
} // namespace GooseFEM
#endif
Event Timeline
Log In to Comment