Page Menu
Home
c4science
Search
Configure Global Search
Log In
Files
F92897606
laminate_homogenisation.hh
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, Nov 24, 15:03
Size
8 KB
Mime Type
text/x-c
Expires
Tue, Nov 26, 15:03 (1 d, 21 h)
Engine
blob
Format
Raw Data
Handle
22534385
Attached To
rMUSPECTRE µSpectre
laminate_homogenisation.hh
View Options
/**
* @file lamminate_hemogenization.hh
*
* @author Ali Falsafi <ali.falsafi@epfl.ch>
*
* @date 28 Sep 2018
*
* @brief Laminatehomogenisation enables one to obtain the resulting stress
* and stiffness tensors of a laminate pixel that is consisted of two
* materialswith a certain normal vector of their interface plane.
* note that it is supposed to be used in static way. so it does note
* any data member. It is merely a collection of functions used to
* calculate effective stress and stiffness.
*
*
* Copyright © 2017 Till Junge, Ali Falsafi
*
* µ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.
*
* 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.
*/
#ifndef SRC_MATERIALS_LAMINATE_HOMOGENISATION_HH_
#define SRC_MATERIALS_LAMINATE_HOMOGENISATION_HH_
#include "common/geometry.hh"
#include "common/muSpectre_common.hh"
#include "libmugrid/field_map.hh"
namespace
muSpectre
{
// TODO(faslafi): to template with the formulation as well
template
<
Dim_t
Dim
>
struct
LamHomogen
{
//! typedefs for data handled by this interface
using
Vec_t
=
Eigen
::
Matrix
<
Real
,
Dim
,
1
>
;
using
Stiffness_t
=
muGrid
::
T4Mat
<
Real
,
Dim
>
;
using
Serial_stiffness_t
=
Eigen
::
Matrix
<
Real
,
2
*
Dim
-
1
,
2
*
Dim
-
1
>
;
using
Parallel_stiffness_t
=
muGrid
::
T4Mat
<
Real
,
(
Dim
-
1
)
>
;
using
Equation_stiffness_t
=
Eigen
::
Matrix
<
Real
,
Dim
,
Dim
>
;
using
Strain_t
=
Eigen
::
Matrix
<
Real
,
Dim
*
Dim
,
1
>
;
using
Serial_strain_t
=
Eigen
::
Matrix
<
Real
,
2
*
Dim
-
1
,
1
>
;
using
Parallel_strain_t
=
Eigen
::
Matrix
<
Real
,
(
Dim
-
1
)
*
(
Dim
-
1
),
1
>
;
using
Equation_strain_t
=
Eigen
::
Matrix
<
Real
,
Dim
,
1
>
;
using
Stress_t
=
Strain_t
;
using
Serial_stress_t
=
Serial_strain_t
;
using
Parallel_stress_t
=
Parallel_strain_t
;
using
Equation_stress_t
=
Equation_strain_t
;
using
StrainMat_t
=
Eigen
::
Matrix
<
Real
,
Dim
,
Dim
>
;
using
StressMat_t
=
StrainMat_t
;
using
Function_t
=
std
::
function
<
std
::
tuple
<
StressMat_t
,
Stiffness_t
>
(
const
Eigen
::
Ref
<
StrainMat_t
>
&
)
>
;
// these two functions return the indexes in column major strain and stress
// tensors that the behavior is either serial or parallel concerning
// combining the laminate layers
static
constexpr
auto
get_serial_indexes
()
->
std
::
array
<
Dim_t
,
2
*
Dim
-
1
>
;
static
constexpr
auto
get_parallel_indexes
()
->
std
::
array
<
Dim_t
,
(
Dim
-
1
)
*
(
Dim
-
1
)
>
;
static
constexpr
auto
get_if_serial
()
->
std
::
array
<
bool
,
Dim
*
Dim
>
;
// these functions return the parts of a stress or strain tensor that
// behave either serial or parallel in a laminate structure
// they are used to obtain a certian part of tensor which is needed by the
// solver used in this struct to calculate the effective stress and
// stiffness or compose the complete tensor from its sub sets.
// obtain the serial part of stress or strain tensor
static
auto
get_serial_strain
(
const
Eigen
::
Ref
<
Strain_t
>
&
F
)
->
Serial_strain_t
;
// obtain the serial part of stress or strain tensor
static
auto
get_equation_strain_from_serial
(
const
Eigen
::
Ref
<
Serial_strain_t
>
&
F
)
->
Equation_strain_t
;
static
auto
get_serial_strain_from_equation
(
const
Eigen
::
Ref
<
Equation_strain_t
>
&
Fs
)
->
Serial_strain_t
;
// obtain the parallel part of stress or strain tensor
static
auto
get_parallel_strain
(
const
Eigen
::
Ref
<
Strain_t
>
&
F
)
->
Parallel_strain_t
;
// compose the complete strain or stress tensor from
// its serial and parallel parts
static
auto
get_total_strain
(
const
Eigen
::
Ref
<
Serial_strain_t
>
&
F_seri
,
const
Eigen
::
Ref
<
Parallel_strain_t
>
&
F_para
)
->
Strain_t
;
// obtain the serial part of stifness tensor
static
auto
get_serial_stiffness
(
const
Eigen
::
Ref
<
Stiffness_t
>
&
K
)
->
Serial_stiffness_t
;
// obtain the parallel part of stifness tensor
static
auto
get_parallel_stiffness
(
const
Eigen
::
Ref
<
Stiffness_t
>
&
K
)
->
Parallel_stiffness_t
;
static
auto
get_equation_stiffness_from_serial
(
const
Eigen
::
Ref
<
Serial_stiffness_t
>
&
Ks
)
->
Equation_stiffness_t
;
// compose the complete stiffness tensor from
// its serial and parallel parts
static
auto
get_total_stiffness
(
const
Eigen
::
Ref
<
Serial_stiffness_t
>
&
K_seri
,
const
Eigen
::
Ref
<
Parallel_stiffness_t
>
&
K_para
)
->
Stiffness_t
;
/**
* this functions calculate the strain in the second materials considering
* the fact that:
* 1. In the directions that the laminate layers are in series the weighted
* average of the strains should be equal to the total strain of the
* laminate
* 2. In the directions that the laminate layers are in parallel the strain
* is constant and equal in all of the layers.
*/
static
auto
linear_eqs
(
Real
ratio
,
const
Eigen
::
Ref
<
Strain_t
>
&
F_0
,
const
Eigen
::
Ref
<
Strain_t
>
&
F_1
)
->
Strain_t
;
static
auto
linear_eqs_seri
(
Real
ratio
,
const
Eigen
::
Ref
<
Strain_t
>
&
F_0
,
const
Eigen
::
Ref
<
Serial_strain_t
>
&
F_1
)
->
Serial_strain_t
;
/**
* the objective in hemgenisation of a single laminate pixel is equating the
* stress in the serial directions so the difference of stress between their
* layers should tend to zero. this function return the stress difference
* and the difference of Stiffness matrices which is used as the Jacobian in
* the solution process
*/
static
auto
delta_stress_eval
(
const
Function_t
&
mat_1_stress_eval
,
const
Function_t
&
mat_2_stress_eval
,
const
Eigen
::
Ref
<
StrainMat_t
>
&
F_1
,
const
Eigen
::
Ref
<
StrainMat_t
>
&
F_2
,
RotatorNormal
<
Dim
>
rotator
,
Real
ratio
)
->
std
::
tuple
<
Serial_stress_t
,
Serial_stiffness_t
>
;
static
auto
delta_stress_eval_F_1
(
const
Function_t
&
mat_1_stress_eval
,
const
Function_t
&
mat_2_stress_eval
,
Eigen
::
Ref
<
Strain_t
>
F_0
,
Eigen
::
Ref
<
Strain_t
>
F_1_rot
,
RotatorNormal
<
Dim
>
rotator
,
Real
ratio
)
->
std
::
tuple
<
Serial_stress_t
,
Serial_stiffness_t
>
;
/**
* the following function s claculates the energy dimension error of the
* solution. it will beused in each step of the solution to detremine the
* relevant difference that implementation of that step has had on
* convergence to the solution.
*/
inline
static
auto
del_energy_eval
(
const
Real
del_F_norm
,
const
Real
delta_P_norm
)
->
Real
;
/**
* This functions calculate the resultant stress and tangent matrices
* according to the computed F_1 and F_2 from the solver.
*
*/
static
auto
lam_P_combine
(
const
Eigen
::
Ref
<
StressMat_t
>
&
P_1
,
const
Eigen
::
Ref
<
StressMat_t
>
&
P_2
,
const
Real
ratio
)
->
StressMat_t
;
static
auto
lam_K_combine
(
const
Eigen
::
Ref
<
Stiffness_t
>
&
K_1
,
const
Eigen
::
Ref
<
Stiffness_t
>
&
K_2
,
const
Real
ratio
)
->
Stiffness_t
;
/**
* This is the main solver function that might be called staically from
* an external file. this will return the resultant stress and stiffness
* tensor according to interanl "equilibrium" of the lamiante.
* The inputs are :
* 1- global Strain
* 2- stress calculation function of the layer 1
* 3- stress calculation function of the layer 2
* 4- the ratio of the first material in the laminate sturucture of the
* pixel 5- the normal vector of the interface of two layers 6- the
* tolerance error for the internal solution of the laminate pixel 7- the
* maximum iterations for the internal solution of the laminate pixel
*/
static
auto
laminate_solver
(
Eigen
::
Ref
<
Strain_t
>
F_coord
,
Function_t
mat_1_stress_eval
,
Function_t
mat_2_stress_eval
,
Real
ratio
,
const
Eigen
::
Ref
<
Vec_t
>
&
normal_vec
,
Real
tol
=
1e-6
,
Dim_t
max_iter
=
1000
)
->
std
::
tuple
<
Dim_t
,
Real
,
StressMat_t
,
Stiffness_t
>
;
};
// LamHomogen
}
// namespace muSpectre
#endif
// SRC_MATERIALS_LAMINATE_HOMOGENISATION_HH_
Event Timeline
Log In to Comment