Page Menu
Home
c4science
Search
Configure Global Search
Log In
Files
F102413754
material_anisotropic.cc
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
Thu, Feb 20, 11:36
Size
3 KB
Mime Type
text/x-c
Expires
Sat, Feb 22, 11:36 (2 d)
Engine
blob
Format
Raw Data
Handle
24351834
Attached To
rMUSPECTRE µSpectre
material_anisotropic.cc
View Options
/**
* @file material_anisotropic.cc
*
* @author Ali Falsafi<ali.falsafi@epfl.ch>
*
* @date 09 Jul 2018
*
* @brief Base class for materials (constitutive models)
*
* Copyright © 2017 Till Junge
*
* µSpectre is free software; you can redistribute it and/or
* modify it under the terms of the GNU General Public License as
* published by the Free Software Foundation, either version 3, or (at
* your option) any later version.
*
* µSpectre 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
* General Public License for more details.
18::00 *
* You should have received a copy of the GNU General Public License
* along with GNU Emacs; see the file COPYING. If not, write to the
* Free Software Foundation, Inc., 59 Temple Place - Suite 330,
* Boston, MA 02111-1307, USA.x
*/
#include "material_anisotropic.hh"
#include "libmugrid/T4_map_proxy.hh"
#include "common/muSpectre_common.hh"
#include "libmugrid/eigen_tools.hh"
#include "libmugrid/tensor_algebra.hh"
#include "common/voigt_conversion.hh"
#include "material_base.hh"
namespace muSpectre {
/* ---------------------------------------------------------------------- */
template <Dim_t DimS, Dim_t DimM>
MaterialAnisotropic<DimS, DimM>::MaterialAnisotropic(
std::string name, std::vector<Real> input_c)
: Parent(name), C{c_maker(input_c)} {};
/* ---------------------------------------------------------------------- */
template <Dim_t DimS, Dim_t DimM>
auto MaterialAnisotropic<DimS, DimM>::c_maker(std::vector<Real> input)
-> Stiffness_t {
std::array<Dim_t, 2> constexpr input_size{6, 21};
Stiffness_t C4{Stiffness_t::Zero()};
Stiffness_t C4_v{Stiffness_t::Zero()};
if (input.size() != input_size[DimM - 2]) {
std::stringstream err_str{};
err_str << "Number of the inputs should be" << input_size[DimM - 2]
<< std::endl;
throw std::runtime_error(err_str.str());
}
constexpr Dim_t v_diff{DimM * DimM - vsize(DimM)};
// memory order voigt -> col major
using t_t = VoigtConversion<DimM>;
auto v_order{t_t::get_vec_vec()}; // non-sym for C
int counter{0};
for (int i{0}; i < vsize(DimM); ++i) {
// diagonal terms
C4(i, i) = input[counter];
counter++;
for (int j{i + 1}; j < vsize(DimM); ++j) {
C4(i, j) = C4(j, i) = input[counter];
if (j >= DimM) {
C4(j + v_diff, i) = C4(j, i + v_diff) = C4(i + v_diff, j) =
C4(i, j + v_diff) = input[counter];
}
counter++;
}
}
C4.bottomRightCorner(v_diff, v_diff) =
C4.block(DimM + v_diff, DimM, v_diff, v_diff) =
C4.block(DimM, DimM + v_diff, v_diff, v_diff) =
C4.block(DimM, DimM, v_diff, v_diff);
for (int i = 0; i < DimM * DimM; i++) {
for (int j = 0; j < DimM * DimM; j++) {
C4_v(i, j) = C4(v_order[i], v_order[j]);
}
}
return C4_v;
}
/* ---------------------------------------------------------------------- */
template class MaterialAnisotropic<twoD, twoD>;
template class MaterialAnisotropic<twoD, threeD>;
template class MaterialAnisotropic<threeD, threeD>;
} // namespace muSpectre
Event Timeline
Log In to Comment