Page Menu
Home
c4science
Search
Configure Global Search
Log In
Files
F62348946
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
Sun, May 12, 14:21
Size
5 KB
Mime Type
text/x-c++
Expires
Tue, May 14, 14:21 (2 d)
Engine
blob
Format
Raw Data
Handle
17636205
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 "kelvin_helper.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
);
}
/* -------------------------------------------------------------------------- */
/* Operator implementation */
/* -------------------------------------------------------------------------- */
template
<
model_type
type
,
UInt
order
>
void
Kelvin
<
type
,
order
>::
applyIf
(
GridBase
<
Real
>&
source
,
GridBase
<
Real
>&
out
,
filter_t
pred
)
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
,
&
pred
](
UInt
node
,
decltype
(
this
->
source_buffers
)
&
source_buffers
,
decltype
(
this
->
disp_buffer
)
&
out_buffer
)
{
constexpr
UInt
derivative
=
order
-
2
;
const
Real
L
=
this
->
model
->
getSystemSize
().
front
();
const
UInt
N
=
this
->
model
->
getDiscretization
().
front
();
const
Real
dl
=
L
/
(
N
-
1
);
const
Real
xi
=
node
*
dl
;
detail
::
LinearElement
<
type
,
decltype
(
kelvin
)
>
element
{
thrust
::
make_pair
(
node
,
xi
),
kelvin
,
source_buffers
,
out_buffer
,
this
->
wavevectors
};
out_buffer
=
0
;
// Simple condition to ensure correct integration
auto
filter_out_elem
=
[
&
pred
](
UInt
e
)
{
return
not
pred
(
e
)
and
not
pred
(
e
+
1
);
};
// Loop over elements
for
(
UInt
e
:
Loop
::
range
(
N
-
1
))
{
if
(
filter_out_elem
(
e
))
continue
;
const
thrust
::
pair
<
Real
,
Real
>
node_positions
{
e
*
dl
,
(
e
+
1
)
*
dl
};
element
.
integrateElement
(
e
,
node_positions
);
}
// Setting fundamental frequency to zero
typename
detail
::
KelvinTrait
<
influence
::
Kelvin
<
trait
::
dimension
,
derivative
>>::
out_t
out_fundamental
(
out_buffer
(
0
));
out_fundamental
=
0
;
};
this
->
fourierApplyIf
(
apply
,
source
,
out
,
pred
);
}
/* -------------------------------------------------------------------------- */
template
<
model_type
type
,
UInt
order
>
void
Kelvin
<
type
,
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