Page Menu
Home
c4science
Search
Configure Global Search
Log In
Files
F68693336
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
Fri, Jun 28, 13:02
Size
10 KB
Mime Type
text/x-c++
Expires
Sun, Jun 30, 13:02 (2 d)
Engine
blob
Format
Raw Data
Handle
18621255
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 "influence.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
<>
void
Kelvin
<
model_type
::
volume_2d
,
2
>::
apply
(
GridBase
<
Real
>&
source
,
GridBase
<
Real
>&
out
)
const
{
Real
nu
=
model
->
getPoissonRatio
(),
mu
=
model
->
getShearModulus
();
VectorProxy
<
const
Real
,
trait
::
dimension
>
domain
(
model
->
getSystemSize
()[
0
]);
influence
::
Kelvin
<
trait
::
dimension
,
0
>
kelvin
(
mu
,
nu
);
auto
apply
=
[
this
,
&
kelvin
](
UInt
i
,
decltype
(
source_buffers
)
&
source_buffers
,
decltype
(
disp_buffer
)
&
displacement
)
{
constexpr
UInt
dim
=
trait
::
dimension
;
const
Real
L
=
this
->
model
->
getSystemSize
().
front
();
const
UInt
N
=
this
->
model
->
getDiscretization
().
front
();
const
Real
dl
=
L
/
(
N
-
1
);
// Compute displacements u_i
displacement
=
0
;
for
(
UInt
j
:
Loop
::
range
(
N
))
{
const
Real
dij
=
j
*
dl
-
i
*
dl
;
// don't factorize!
auto
&
source
=
source_buffers
[
j
];
#define POTENTIAL(yj_xi) \
Loop::stridedLoop( \
[&kelvin, dij, dl](VectorProxy<Complex, dim>&& u, \
VectorProxy<Complex, dim>&& f, \
VectorProxy<const Real, dim - 1>&& q) { \
/* Cutoff */
\
if (-q.l2norm() * std::abs(dij) < std::log(1e-2)) \
return; \
influence::KelvinIntegrator<1>::integrate<yj_xi>(u, f, kelvin, q, dij, \
dl); \
}, \
displacement, source, this->wavevectors)
if
(
j
>
i
)
{
POTENTIAL
(
1
);
}
else
if
(
j
==
i
)
{
POTENTIAL
(
0
);
}
else
{
POTENTIAL
(
-
1
);
}
#undef POTENTIAL
}
// Setting fundamental frequency to zero
VectorProxy
<
Complex
,
dim
>
u_fundamental
(
displacement
(
0
));
u_fundamental
=
0
;
};
this
->
fourierApply
(
apply
,
source
,
out
);
}
/* -------------------------------------------------------------------------- */
template
<>
void
Kelvin
<
model_type
::
volume_2d
,
3
>::
apply
(
GridBase
<
Real
>&
source
,
GridBase
<
Real
>&
out
)
const
{
Real
nu
=
model
->
getPoissonRatio
(),
mu
=
model
->
getShearModulus
();
VectorProxy
<
const
Real
,
trait
::
dimension
>
domain
(
model
->
getSystemSize
()[
0
]);
influence
::
Kelvin
<
trait
::
dimension
,
1
>
kelvin
(
mu
,
nu
);
auto
apply
=
[
this
,
&
kelvin
](
UInt
i
,
decltype
(
source_buffers
)
&
source_buffers
,
decltype
(
disp_buffer
)
&
displacement
)
{
constexpr
UInt
dim
=
trait
::
dimension
;
const
Real
L
=
this
->
model
->
getSystemSize
().
front
();
const
UInt
N
=
this
->
model
->
getDiscretization
().
front
();
const
Real
dl
=
L
/
(
N
-
1
);
// Compute displacements u_i
displacement
=
0
;
for
(
UInt
j
:
Loop
::
range
(
N
))
{
const
Real
dij
=
j
*
dl
-
i
*
dl
;
// don't factorize!
auto
&
source
=
source_buffers
[
j
];
#define POTENTIAL(yj_xi) \
Loop::stridedLoop( \
[&kelvin, dij, dl](VectorProxy<Complex, dim>&& u, \
MatrixProxy<Complex, dim, dim>&& f, \
VectorProxy<const Real, dim - 1>&& q) { \
/* Cutoff */
\
if (-q.l2norm() * std::abs(dij) < std::log(1e-2)) \
return; \
influence::KelvinIntegrator<1>::integrate<yj_xi>(u, f, kelvin, q, dij, \
dl); \
}, \
displacement, source, this->wavevectors)
if
(
j
>
i
)
{
POTENTIAL
(
1
);
}
else
if
(
j
==
i
)
{
POTENTIAL
(
0
);
}
else
{
POTENTIAL
(
-
1
);
}
#undef POTENTIAL
}
// Setting fundamental frequency to zero
VectorProxy
<
Complex
,
dim
>
u_fundamental
(
displacement
(
0
));
u_fundamental
=
0
;
};
this
->
fourierApply
(
apply
,
source
,
out
);
}
/* -------------------------------------------------------------------------- */
template
<>
void
Kelvin
<
model_type
::
volume_2d
,
4
>::
apply
(
GridBase
<
Real
>&
source
,
GridBase
<
Real
>&
out
)
const
{
Real
nu
=
model
->
getPoissonRatio
(),
mu
=
model
->
getShearModulus
();
VectorProxy
<
const
Real
,
trait
::
dimension
>
domain
(
model
->
getSystemSize
()[
0
]);
influence
::
Kelvin
<
trait
::
dimension
,
2
>
kelvin
(
mu
,
nu
);
auto
apply
=
[
this
,
&
kelvin
](
UInt
i
,
decltype
(
source_buffers
)
&
source_buffers
,
decltype
(
disp_buffer
)
&
gradu
)
{
constexpr
UInt
dim
=
trait
::
dimension
;
const
Real
L
=
this
->
model
->
getSystemSize
().
front
();
const
UInt
N
=
this
->
model
->
getDiscretization
().
front
();
const
Real
dl
=
L
/
(
N
-
1
);
// Compute displacement gradient
gradu
=
0
;
for
(
UInt
j
:
Loop
::
range
(
N
))
{
const
Real
dij
=
j
*
dl
-
i
*
dl
;
// don't factorize!
auto
&
source
=
source_buffers
[
j
];
#define POTENTIAL(yj_xi) \
Loop::stridedLoop( \
[&kelvin, dij, dl](MatrixProxy<Complex, dim, dim>&& u, \
MatrixProxy<Complex, dim, dim>&& f, \
VectorProxy<const Real, dim - 1>&& q) { \
/* Cutoff */
\
if (-q.l2norm() * std::abs(dij) < std::log(1e-2)) \
return; \
influence::KelvinIntegrator<0>::integrate<yj_xi>(u, f, kelvin, q, dij, \
dl); \
}, \
gradu, source, this->wavevectors)
if
(
j
>
i
)
{
POTENTIAL
(
1
);
}
else
if
(
j
==
i
)
{
POTENTIAL
(
0
);
// Additional free term from discontinous Kelvin gradient
Loop
::
stridedLoop
(
[
&
kelvin
](
MatrixProxy
<
Complex
,
dim
,
dim
>&&
u
,
MatrixProxy
<
Complex
,
dim
,
dim
>&&
f
,
VectorProxy
<
const
Real
,
dim
-
1
>&&
q
)
{
u
+=
kelvin
.
applyDiscontinuityTerm
(
q
,
f
);
},
gradu
,
source
,
this
->
wavevectors
);
}
else
{
POTENTIAL
(
-
1
);
}
#undef POTENTIAL
}
gradu
*=
-
1
;
// Setting fundamental frequency to zero
MatrixProxy
<
Complex
,
dim
,
dim
>
u_fundamental
(
gradu
(
0
));
u_fundamental
=
0
;
};
this
->
fourierApply
(
apply
,
source
,
out
);
}
/* -------------------------------------------------------------------------- */
template
<
model_type
type
,
UInt
tensor_order
>
void
Kelvin
<
type
,
tensor_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
<
model_type
type
,
UInt
tensor_order
>
void
Kelvin
<
type
,
tensor_order
>::
apply
(
GridBase
<
Real
>&
source
,
GridBase
<
Real
>&
out
)
const
{
TAMAAS_EXCEPTION
(
"The requested operator has not been implemented"
);
}
/* -------------------------------------------------------------------------- */
/* 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