Page Menu
Home
c4science
Search
Configure Global Search
Log In
Files
F90846980
material_orthotropic.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, Nov 5, 07:33
Size
4 KB
Mime Type
text/x-c
Expires
Thu, Nov 7, 07:33 (2 d)
Engine
blob
Format
Raw Data
Handle
22145604
Attached To
rMUSPECTRE µSpectre
material_orthotropic.hh
View Options
/**
* @file material_anisotropic.hh
*
* @author Ali Falsafi<ali.falsafi@epfl.ch>
*
* @date 11 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.
*
* 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 MATERIAL_ORTHOTROPIC_H
#define MATERIAL_ORTHOTROPIC_H
#include "materials/material_base.hh"
#include "materials/material_muSpectre_base.hh"
#include "materials/material_anisotropic.hh"
#include "common/common.hh"
#include "common/T4_map_proxy.hh"
#include "cell/cell_base.hh"
namespace
muSpectre
{
// Forward declaration for factory function
//template <Dim_t DimS, Dim_t DimM>
//class CellBase;
template
<
Dim_t
DimS
,
Dim_t
DimM
>
class
MaterialOrthotropic
;
//traits for orthotropic material
template
<
Dim_t
DimS
,
Dim_t
DimM
>
struct
MaterialMuSpectre_traits
<
MaterialOrthotropic
<
DimS
,
DimM
>>
{
using
Parent
=
MaterialMuSpectre_traits
<
void
>
;
//!< base for elasticity
//! global field collection
using
GFieldCollection_t
=
typename
MaterialBase
<
DimS
,
DimM
>::
GFieldCollection_t
;
//! expected map type for strain fields
using
StrainMap_t
=
MatrixFieldMap
<
GFieldCollection_t
,
Real
,
DimM
,
DimM
,
true
>
;
//! expected map type for stress fields
using
StressMap_t
=
MatrixFieldMap
<
GFieldCollection_t
,
Real
,
DimM
,
DimM
>
;
//! expected map type for tangent stiffness fields
using
TangentMap_t
=
T4MatrixFieldMap
<
GFieldCollection_t
,
Real
,
DimM
>
;
//! declare what type of strain measure your law takes as input
constexpr
static
auto
strain_measure
{
StrainMeasure
::
GreenLagrange
};
//! declare what type of stress measure your law yields as output
constexpr
static
auto
stress_measure
{
StressMeasure
::
PK2
};
//! anisotropicity without internal variables
using
InternalVariables
=
std
::
tuple
<>
;
};
/**
* Material implementation for orthotropic constitutive law
*/
template
<
Dim_t
DimS
,
Dim_t
DimM
=
DimS
>
class
MaterialOrthotropic
:
public
MaterialAnisotropic
<
DimS
,
DimM
>
{
//! base class
using
Parent
=
MaterialAnisotropic
<
DimS
,
DimM
>
;
using
Stiffness_t
=
T4Mat
<
Real
,
DimM
>
;
//! traits of this material
using
traits
=
MaterialMuSpectre_traits
<
MaterialOrthotropic
>
;
//! this law does not have any internal variables
using
InternalVariables
=
typename
traits
::
InternalVariables
;
public
:
//! global field collection
using
GFieldCollection_t
=
typename
MaterialBase
<
DimS
,
DimM
>::
GFieldCollection_t
;
//! expected map type for tangent stiffness fields
using
Tangent_t
=
T4MatrixFieldMap
<
GFieldCollection_t
,
Real
,
DimM
>
;
//! Default constructor
MaterialOrthotropic
()
=
delete
;
MaterialOrthotropic
(
std
::
string
name
,
std
::
vector
<
Real
>
input
);
//! Copy constructor
MaterialOrthotropic
(
const
MaterialOrthotropic
&
other
)
=
delete
;
//! Move constructor
MaterialOrthotropic
(
MaterialOrthotropic
&&
other
)
=
delete
;
//! Destructor
virtual
~
MaterialOrthotropic
()
=
default
;
/* overloaded make function in order to make pyrhon binding
able to make an object of materila orthotropic*/
static
MaterialOrthotropic
<
DimS
,
DimM
>
&
make
(
CellBase
<
DimS
,
DimM
>
&
cell
,
std
::
string
name
,
std
::
vector
<
Real
>
input
);
std
::
vector
<
Real
>
input_c_maker
(
std
::
vector
<
Real
>
input
);
protected
:
/**
* these variable are used to determine which elements of the
* stiffness matrix should be replaced with the inpts for the
* orthotropic material
*/
constexpr
static
std
::
array
<
std
::
size_t
,
2
>
output_size
{
6
,
21
};
static
std
::
array
<
bool
,
output_size
[
DimM
-
2
]
>
ret_flag
;
};
/* ---------------------------------------------------------------------- */
template
<
Dim_t
DimS
,
Dim_t
DimM
>
MaterialOrthotropic
<
DimS
,
DimM
>
&
MaterialOrthotropic
<
DimS
,
DimM
>::
make
(
CellBase
<
DimS
,
DimM
>
&
cell
,
std
::
string
name
,
std
::
vector
<
Real
>
input
){
auto
mat
=
std
::
make_unique
<
MaterialOrthotropic
<
DimS
,
DimM
>>
(
name
,
input
);
auto
&
mat_ref
=
*
mat
;
cell
.
add_material
(
std
::
move
(
mat
));
return
mat_ref
;
}
}
//muSpectre
#endif
/* MATERIAL_ORTHOTROPIC_H */
Event Timeline
Log In to Comment