Page Menu
Home
c4science
Search
Configure Global Search
Log In
Files
F74939234
material_thermal.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
Tue, Jul 30, 12:19
Size
4 KB
Mime Type
text/x-c++
Expires
Thu, Aug 1, 12:19 (2 d)
Engine
blob
Format
Raw Data
Handle
19459119
Attached To
rAKA akantu
material_thermal.hh
View Options
/**
* @file material_thermal.hh
*
* @author Lucas Frerot <lucas.frerot@epfl.ch>
*
* @date creation: Fri Jun 18 2010
* @date last modification: Fri Apr 09 2021
*
* @brief Material isotropic thermo-elastic
*
*
* @section LICENSE
*
* Copyright (©) 2010-2021 EPFL (Ecole Polytechnique Fédérale de Lausanne)
* Laboratory (LSMS - Laboratoire de Simulation en Mécanique des Solides)
*
* Akantu 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.
*
* Akantu 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 Akantu. If not, see <http://www.gnu.org/licenses/>.
*
*/
/* -------------------------------------------------------------------------- */
#include "aka_common.hh"
#include "material.hh"
/* -------------------------------------------------------------------------- */
#ifndef AKANTU_MATERIAL_THERMAL_HH_
#define AKANTU_MATERIAL_THERMAL_HH_
namespace
akantu
{
template
<
Int
dim
>
class
MaterialThermal
:
public
Material
{
/* ------------------------------------------------------------------------ */
/* Constructors/Destructors */
/* ------------------------------------------------------------------------ */
public
:
MaterialThermal
(
SolidMechanicsModel
&
model
,
const
ID
&
id
=
""
);
MaterialThermal
(
SolidMechanicsModel
&
model
,
Int
spatial_dimension
,
const
Mesh
&
mesh
,
FEEngine
&
fe_engine
,
const
ID
&
id
=
""
);
~
MaterialThermal
()
override
=
default
;
protected
:
void
initialize
();
/* ------------------------------------------------------------------------ */
/* Methods */
/* ------------------------------------------------------------------------ */
public
:
void
initMaterial
()
override
;
/// constitutive law for all element of a type
void
computeStress
(
ElementType
el_type
,
GhostType
ghost_type
)
override
;
/// local computation of thermal stress
template
<
class
Args
>
inline
void
computeStressOnQuad
(
Args
&&
args
);
/* ------------------------------------------------------------------------ */
decltype
(
auto
)
getArguments
(
ElementType
el_type
,
GhostType
ghost_type
)
{
return
zip_append
(
Material
::
getArguments
<
dim
>
(
el_type
,
ghost_type
),
tuple
::
get
<
"delta_T"
_h
>
()
=
make_view
(
this
->
delta_T
(
el_type
,
ghost_type
)),
tuple
::
get
<
"sigma_th"
_h
>
()
=
make_view
(
this
->
sigma_th
(
el_type
,
ghost_type
)));
}
/* ------------------------------------------------------------------------ */
/* Class Members */
/* ------------------------------------------------------------------------ */
protected
:
/// Young modulus
Real
E
;
/// Poisson ratio
Real
nu
;
/// Thermal expansion coefficient
/// TODO : implement alpha as a matrix
Real
alpha
;
/// Temperature field
InternalField
<
Real
>
delta_T
;
/// Current thermal stress
InternalField
<
Real
>
sigma_th
;
/// Tell if we need to use the previous thermal stress
bool
use_previous_stress_thermal
;
};
/* ------------------------------------------------------------------------ */
/* Inline impl */
/* ------------------------------------------------------------------------ */
template
<
Int
dim
>
template
<
class
Args
>
inline
void
MaterialThermal
<
dim
>::
computeStressOnQuad
(
Args
&&
args
)
{
auto
&&
sigma
=
tuple
::
get
<
"sigma_th"
_h
>
(
args
);
auto
&&
deltaT
=
tuple
::
get
<
"delta_T"
_h
>
(
args
);
sigma
=
-
this
->
E
/
(
1.
-
2.
*
this
->
nu
)
*
this
->
alpha
*
deltaT
;
}
template
<>
template
<
class
Args
>
inline
void
MaterialThermal
<
1
>::
computeStressOnQuad
(
Args
&&
args
)
{
auto
&&
sigma
=
tuple
::
get
<
"sigma_th"
_h
>
(
args
);
auto
&&
deltaT
=
tuple
::
get
<
"delta_T"
_h
>
(
args
);
sigma
=
-
this
->
E
*
this
->
alpha
*
deltaT
;
}
}
// namespace akantu
#endif
/* AKANTU_MATERIAL_THERMAL_HH_ */
Event Timeline
Log In to Comment