Page Menu
Home
c4science
Search
Configure Global Search
Log In
Files
F69659906
kelvin.cpp
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
Tue, Jul 2, 21:34
Size
7 KB
Mime Type
text/x-c++
Expires
Thu, Jul 4, 21:34 (2 d)
Engine
blob
Format
Raw Data
Handle
18738730
Attached To
rTAMAAS tamaas
kelvin.cpp
View Options
/**
* @file
*
* @author Lucas Frérot <lucas.frerot@epfl.ch>
*
* @section LICENSE
*
* Copyright (©) 2017 EPFL (Ecole Polytechnique Fédérale de
* Lausanne) Laboratory (LSMS - Laboratoire de Simulation en Mécanique des
* Solides)
*
* Tamaas is free software: you can redistribute it and/or modify it under the
* terms of the GNU Lesser General Public License as published by the Free
* Software Foundation, either version 3 of the License, or (at your option) any
* later version.
*
* Tamaas is distributed in the hope that it will be useful, but WITHOUT ANY
* WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR
* A PARTICULAR PURPOSE. See the GNU Lesser General Public License for more
* details.
*
* You should have received a copy of the GNU Lesser General Public License
* along with Tamaas. If not, see <http://www.gnu.org/licenses/>.
*
*/
/* -------------------------------------------------------------------------- */
#include "kelvin.hh"
#include "elasto_plastic_model.hh"
#include "influence.hh"
/* -------------------------------------------------------------------------- */
__BEGIN_TAMAAS__
/* -------------------------------------------------------------------------- */
/* -------------------------------------------------------------------------- */
/* Constructors */
/* -------------------------------------------------------------------------- */
template <>
Kelvin<model_type::volume_2d, 2>::Kelvin(Model* model)
: VolumePotential(model) {
initialize(trait::components, trait::components);
}
template <>
Kelvin<model_type::volume_2d, 3>::Kelvin(Model* model)
: VolumePotential(model) {
initialize(trait::components * trait::components, trait::components);
}
template <>
Kelvin<model_type::volume_2d, 4>::Kelvin(Model* model)
: VolumePotential(model) {
initialize(trait::components * trait::components,
trait::components * trait::components);
}
/* -------------------------------------------------------------------------- */
namespace detail {
template <typename T>
struct KelvinTrait;
template <UInt dim>
struct KelvinTrait<influence::Kelvin<dim, 0>> {
using source_t = VectorProxy<Complex, dim>;
using out_t = VectorProxy<Complex, dim>;
};
template <UInt dim>
struct KelvinTrait<influence::Kelvin<dim, 1>> {
using source_t = MatrixProxy<Complex, dim, dim>;
using out_t = VectorProxy<Complex, dim>;
};
template <UInt dim>
struct KelvinTrait<influence::Kelvin<dim, 2>> {
using source_t = MatrixProxy<Complex, dim, dim>;
using out_t = MatrixProxy<Complex, dim, dim>;
};
template <model_type type, typename T>
struct KelvinHelper {
using kelvin_t = T;
using trait = model_type_traits<type>;
static constexpr UInt dim = trait::dimension;
static constexpr UInt bdim = trait::boundary_dimension;
KelvinHelper(const T& kelvin, Real dij, Real dl, Real cutoff)
: kelvin(kelvin), dij(dij), dl(dl), cutoff(cutoff) {}
template <Int yj_xi>
inline void applyIntegral(GridHermitian<Real, bdim>& out,
GridHermitian<Real, bdim>& source,
const Grid<Real, bdim>& wavevectors) const {
Loop::stridedLoop(
[&](typename KelvinTrait<T>::out_t&& out_local,
typename KelvinTrait<T>::source_t&& source_local,
VectorProxy<const Real, bdim>&& q) {
// Cutoff
if (-q.l2norm() * std::abs(dij) < std::log(cutoff))
return;
influence::KelvinIntegrator<1>::integrate<yj_xi>(
out_local, source_local, kelvin, q, dij, dl);
},
out, source, wavevectors);
}
inline void applyFreeTerm(GridHermitian<Real, bdim>& out,
GridHermitian<Real, bdim>& source,
const Grid<Real, bdim>& wavevectors) const {}
const T& kelvin;
const Real dij, dl, cutoff;
};
template <>
inline void
KelvinHelper<model_type::volume_2d, influence::Kelvin<3, 2>>::applyFreeTerm(
GridHermitian<Real, bdim>& out, GridHermitian<Real, bdim>& source,
const Grid<Real, bdim>& wavevectors) const {
Loop::stridedLoop(
[&](KelvinTrait<kelvin_t>::out_t&& out_local,
KelvinTrait<kelvin_t>::source_t&& source_local,
VectorProxy<const Real, bdim>&& q) {
out_local += kelvin.applyDiscontinuityTerm(q, source_local);
},
out, source, wavevectors);
}
} // namespace detail
/* -------------------------------------------------------------------------- */
/* Operator implementation */
/* -------------------------------------------------------------------------- */
template <model_type type, UInt order>
void Kelvin<type, order>::apply(GridBase<Real>& source,
GridBase<Real>& out) const {
constexpr UInt derivative = order - 2;
Real nu = this->model->getPoissonRatio(), mu = this->model->getShearModulus();
influence::Kelvin<trait::dimension, derivative> kelvin(mu, nu);
auto apply = [this, &kelvin](UInt i, decltype(this->source_buffers)& source_buffers,
decltype(this->disp_buffer)& out_buffer) {
constexpr UInt derivative = order - 2;
constexpr UInt dim = trait::dimension;
const Real L = this->model->getSystemSize().front();
const UInt N = this->model->getDiscretization().front();
const Real dl = L / (N - 1);
out_buffer = 0;
for (UInt j : Loop::range(N)) {
const Real dij = j * dl - i * dl; // don't factorize!
auto& source = source_buffers[j];
detail::KelvinHelper<type, influence::Kelvin<dim, derivative>> helper(
kelvin, dij, dl, 1e-2);
if (j > i) {
helper.template applyIntegral<1>(out_buffer, source, this->wavevectors);
} else if (j == i) {
helper.template applyIntegral<0>(out_buffer, source, this->wavevectors);
helper.applyFreeTerm(out_buffer, source, this->wavevectors);
} else {
helper.template applyIntegral<-1>(out_buffer, source, this->wavevectors);
}
}
// Setting fundamental frequency to zero
typename detail::KelvinTrait<
influence::Kelvin<trait::dimension, derivative>>::out_t
out_fundamental(out_buffer(0));
out_fundamental = 0;
};
this->fourierApply(apply, source, out);
}
/* -------------------------------------------------------------------------- */
template <model_type type, UInt tensor_order>
void Kelvin<type, tensor_order>::initialize(UInt source_components,
UInt out_components) {
// Copy horizontal sizes
std::array<UInt, trait::boundary_dimension> sizes;
std::copy(this->model->getDiscretization().begin() + 1,
this->model->getDiscretization().end(), sizes.begin());
auto hermitian_sizes =
GridHermitian<Real, trait::boundary_dimension>::hermitianDimensions(
sizes);
/// Initializing buffers
this->source_buffers.resize(this->model->getDiscretization()[0]);
std::for_each(this->source_buffers.begin(), this->source_buffers.end(),
[&](typename parent::BufferType& buffer) {
buffer.setNbComponents(source_components);
buffer.resize(hermitian_sizes);
});
this->disp_buffer.setNbComponents(out_components);
this->disp_buffer.resize(hermitian_sizes);
}
/* -------------------------------------------------------------------------- */
/* Template instanciation */
/* -------------------------------------------------------------------------- */
template class Kelvin<model_type::volume_2d, 2>;
template class Kelvin<model_type::volume_2d, 3>;
template class Kelvin<model_type::volume_2d, 4>;
/* -------------------------------------------------------------------------- */
__END_TAMAAS__
Event Timeline
Log In to Comment