Page Menu
Home
c4science
Search
Configure Global Search
Log In
Files
F71184321
bem_gigi.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
Wed, Jul 10, 04:10
Size
11 KB
Mime Type
text/x-c
Expires
Fri, Jul 12, 04:10 (1 d, 23 h)
Engine
blob
Format
Raw Data
Handle
18918084
Attached To
rTAMAAS tamaas
bem_gigi.cpp
View Options
/**
*
* @author Guillaume Anciaux <guillaume.anciaux@epfl.ch>
*
* @section LICENSE
*
* Copyright (©) 2016 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 <vector>
#include "surface.hh"
#include "bem_gigi.hh"
#include <iostream>
#include <sstream>
#include <fstream>
#include <iomanip>
#include <sstream>
#include <cmath>
/* -------------------------------------------------------------------------- */
#define TIMER
#include "surface_timer.hh"
/* -------------------------------------------------------------------------- */
Real
BemGigi
::
computeEquilibrium
(
Real
epsilon
,
Real
pressure
,
int
nthreads
){
this
->
setNumberOfThreads
(
nthreads
);
this
->
computeSpectralInfluenceOverDisplacement
();
this
->
surface_t
=
0.
;
this
->
surface_r
=
0.
;
this
->
surface_spectral_r
=
0.
;
this
->
surface_w
=
0.
;
this
->
surface_q
=
0.
;
this
->
surface_spectral_q
=
0.
;
this
->
pold
=
0.
;
this
->
surface_displacements
=
0.
;
Real
Gold
=
1.
;
Real
Rold
=
1.
;
Real
tau
=
0.
;
Real
f
=
1e300
;
Real
fPrevious
=
1e300
;
// Real current_pressure = SurfaceStatistics::computeAverage(surface_tractions);
surface_tractions
=
pressure
;
this
->
enforcePressureBalance
(
pressure
);
convergence_iterations
.
clear
();
nb_iterations
=
0
;
while
(
f
>
epsilon
&&
nb_iterations
<
max_iterations
)
{
++
nb_iterations
;
fPrevious
=
f
;
this
->
computeDisplacements
();
this
->
functional
->
computeGradF
();
this
->
functional
->
centerGradFOnContactArea
();
Real
R
=
this
->
computeR
();
this
->
updateW
(
R
,
Rold
);
Real
alpha
=
this
->
computeAlpha
();
pold
=
this
->
surface_tractions
;
Rold
=
R
;
Gold
=
R
;
this
->
updatePressurea
(
alpha
);
tau
=
alpha
;
this
->
computeGaps
();
UInt
Iol
=
this
->
computeIol
();
Real
G
=
this
->
computeG
();
//Inner loop
while
(
G
>
epsilon
)
{
++
nb_iterations
;
this
->
emptyoverlap
(
tau
);
this
->
enforcePressureBalance
(
pressure
);
this
->
functional
->
computeF
();
f
=
this
->
functional
->
getF
();
fPrevious
=
f
;
this
->
computeDisplacements
();
this
->
computeGaps
();
Real
gbar
=
this
->
computeMeanGapsInContact
();
this
->
gap
-=
gbar
;
G
=
this
->
computeG
();
this
->
updateT
(
G
,
Gold
);
Real
tau
=
this
->
computeTau
();
pold
=
this
->
surface_tractions
;
Gold
=
G
;
this
->
updatePressureb
(
tau
);
Iol
=
this
->
computeIol
();
}
this
->
enforcePressureBalance
(
pressure
);
this
->
functional
->
computeF
();
f
=
this
->
functional
->
getF
();
if
(
nb_iterations
%
dump_freq
==
0
){
Real
A
=
SurfaceStatistics
::
computeContactAreaRatio
(
surface_tractions
);
std
::
cout
<<
std
::
scientific
<<
std
::
setprecision
(
10
)
<<
nb_iterations
<<
" "
<<
f
<<
" "
<<
f
-
fPrevious
<<
" "
<<
A
<<
std
::
endl
;
}
convergence_iterations
.
push_back
(
f
);
}
return
f
;
}
/* -------------------------------------------------------------------------- */
void
BemGigi
::
truncatePressure
(){
STARTTIMER
(
"truncatePressure"
);
unsigned
int
n
=
surface
.
size
();
unsigned
int
size
=
n
*
n
;
#pragma omp parallel for
for
(
unsigned
int
i
=
0
;
i
<
size
;
++
i
)
{
if
(
this
->
surface_tractions
(
i
)
<
0.
)
this
->
surface_tractions
(
i
)
=
0
;
}
STOPTIMER
(
"truncatePressure"
);
}
/* -------------------------------------------------------------------------- */
Real
BemGigi
::
computeR
(){
STARTTIMER
(
"computeR"
);
unsigned
int
n
=
surface
.
size
();
unsigned
int
size
=
n
*
n
;
Real
res
=
0.
;
Surface
<
Real
>
&
gradF
=
this
->
functional
->
getGradF
();
#pragma omp parallel for reduction(+: res)
for
(
unsigned
int
i
=
0
;
i
<
size
;
++
i
)
{
Real
val
=
gradF
(
i
);
res
+=
val
*
val
;
}
STOPTIMER
(
"computeR"
);
return
res
;
}
/* -------------------------------------------------------------------------- */
Real
BemGigi
::
computeIol
(){
STARTTIMER
(
"computeIol"
);
unsigned
int
n
=
surface
.
size
();
unsigned
int
size
=
n
*
n
;
UInt
Iol
=
0
;
#pragma omp parallel for
for
(
unsigned
int
i
=
0
;
i
<
size
;
++
i
)
{
if
(
this
->
surface_tractions
(
i
)
<=
0.
&&
this
->
gap
(
i
)
<
0.
){
++
Iol
;}
}
STOPTIMER
(
"Iol"
);
return
Iol
;
}
/* -------------------------------------------------------------------------- */
void
BemGigi
::
updateW
(
Real
R
,
Real
Rold
){
STARTTIMER
(
"updateW"
);
unsigned
int
n
=
surface
.
size
();
unsigned
int
size
=
n
*
n
;
Real
factor
=
R
/
Rold
;
Surface
<
Real
>
&
gradF
=
this
->
functional
->
getGradF
();
#pragma omp parallel for
for
(
unsigned
int
i
=
0
;
i
<
size
;
++
i
)
{
this
->
surface_w
(
i
)
*=
factor
;
this
->
surface_w
(
i
)
+=
gradF
(
i
);
}
STOPTIMER
(
"updateW"
);
}
/* -------------------------------------------------------------------------- */
void
BemGigi
::
updateT
(
Real
G
,
Real
Gold
){
STARTTIMER
(
"updateT"
);
unsigned
int
n
=
surface
.
size
();
unsigned
int
size
=
n
*
n
;
Real
factor
=
G
/
Gold
;
#pragma omp parallel for
for
(
unsigned
int
i
=
0
;
i
<
size
;
++
i
)
{
if
(
this
->
surface_tractions
(
i
)
<=
0
)
this
->
surface_t
(
i
)
=
0.
;
else
{
this
->
surface_t
(
i
)
*=
factor
;
this
->
surface_t
(
i
)
+=
this
->
gap
(
i
);
}
}
STOPTIMER
(
"updateT"
);
}
/* -------------------------------------------------------------------------- */
Real
BemGigi
::
computeAlpha
(){
STARTTIMER
(
"computeAlpha"
);
surface_spectral_q
=
surface_w
;
surface_spectral_q
.
FFTTransform
(
nthreads
);
unsigned
int
n
=
surface
.
size
();
unsigned
int
size
=
n
*
n
;
#pragma omp parallel for
for
(
unsigned
int
i
=
0
;
i
<
size
;
++
i
)
{
surface_spectral_q
(
i
)
*=
this
->
surface_spectral_influence_disp
(
i
);
}
surface_spectral_q
.
FFTITransform
(
surface_q
,
nthreads
);
Real
qbar
=
0
;
UInt
nb_contact
=
0
;
#pragma omp parallel for reduction(+: nb_contact,qbar)
for
(
unsigned
int
i
=
0
;
i
<
size
;
++
i
)
{
++
nb_contact
;
qbar
+=
surface_q
(
i
);
}
qbar
/=
nb_contact
;
surface_q
-=
qbar
;
Real
alpha_sum1
=
0.
,
alpha_sum2
=
0.
;
Surface
<
Real
>
&
gradF
=
this
->
functional
->
getGradF
();
#pragma omp parallel for reduction(+: alpha_sum1, alpha_sum2)
for
(
unsigned
int
i
=
0
;
i
<
size
;
++
i
)
{
alpha_sum1
+=
gradF
(
i
)
*
surface_w
(
i
);
alpha_sum2
+=
surface_q
(
i
)
*
surface_w
(
i
);
}
Real
alpha
=
alpha_sum1
/
alpha_sum2
;
STOPTIMER
(
"computeAlpha"
);
return
alpha
;
}
/* -------------------------------------------------------------------------- */
void
BemGigi
::
updatePressurea
(
Real
alpha
){
STARTTIMER
(
"updatePressurea"
);
unsigned
int
n
=
surface
.
size
();
unsigned
int
size
=
n
*
n
;
#pragma omp parallel for
for
(
unsigned
int
i
=
0
;
i
<
size
;
++
i
)
{
this
->
surface_tractions
(
i
)
-=
alpha
*
this
->
surface_w
(
i
);
}
STOPTIMER
(
"updatePressurea"
);
}
/* -------------------------------------------------------------------------- */
void
BemGigi
::
updatePressureb
(
Real
tau
){
STARTTIMER
(
"updatePressureb"
);
unsigned
int
n
=
surface
.
size
();
unsigned
int
size
=
n
*
n
;
#pragma omp parallel for
for
(
unsigned
int
i
=
0
;
i
<
size
;
++
i
)
{
if
(
this
->
surface_tractions
(
i
)
<=
0
)
continue
;
this
->
surface_tractions
(
i
)
-=
tau
*
this
->
surface_t
(
i
);
}
STOPTIMER
(
"updatePressureb"
);
}
/* -------------------------------------------------------------------------- */
void
BemGigi
::
emptyoverlap
(
Real
tau
){
STARTTIMER
(
"emptyoverlap"
);
unsigned
int
n
=
surface
.
size
();
unsigned
int
size
=
n
*
n
;
#pragma omp parallel for
for
(
unsigned
int
i
=
0
;
i
<
size
;
++
i
)
{
if
(
this
->
surface_tractions
(
i
)
<=
0.
&&
this
->
gap
(
i
)
<
0.
){
this
->
surface_tractions
(
i
)
=
-
this
->
surface_tractions
(
i
)
-
tau
*
this
->
gap
(
i
);
}
}
STOPTIMER
(
"emptyoverlap"
);
}
/* -------------------------------------------------------------------------- */
void
BemGigi
::
enforcePressureBalance
(
Real
applied_pressure
){
STARTTIMER
(
"enforcePressureBalance"
);
unsigned
int
n
=
surface
.
size
();
unsigned
int
size
=
n
*
n
;
Real
pressure
=
0.
;
#pragma omp parallel for reduction(+: pressure)
for
(
unsigned
int
i
=
0
;
i
<
size
;
++
i
)
{
pressure
+=
this
->
surface_tractions
(
i
);
}
pressure
*=
1.
/
size
;
// std::cerr << std::scientific << std::setprecision(10)
// << "pressure " << pressure << std::endl;
#pragma omp parallel for
for
(
unsigned
int
i
=
0
;
i
<
size
;
++
i
)
{
this
->
surface_tractions
(
i
)
*=
applied_pressure
/
pressure
;
}
STOPTIMER
(
"enforcePressureBalance"
);
}
/* -------------------------------------------------------------------------- */
/* -------------------------------------------------------------------------- */
void
BemGigi
::
computeGaps
()
{
UInt
size
=
surface
.
size
();
#pragma omp parallel for
for
(
UInt
i
=
0
;
i
<
size
*
size
;
i
++
)
{
gap
(
i
)
=
surface_displacements
(
i
)
-
surface
(
i
);
}
}
/* -------------------------------------------------------------------------- */
Real
BemGigi
::
computeMeanGapsInContact
(){
STARTTIMER
(
"computeMeanGapsInContact"
);
unsigned
int
n
=
surface
.
size
();
unsigned
int
size
=
n
*
n
;
Real
res
=
0.
;
UInt
nb_contact
=
0
;
#pragma omp parallel for reduction(+: nb_contact,res)
for
(
unsigned
int
i
=
0
;
i
<
size
;
++
i
)
{
if
(
this
->
surface_tractions
(
i
)
<=
0
)
continue
;
++
nb_contact
;
res
+=
this
->
gap
(
i
);
}
res
/=
nb_contact
;
STOPTIMER
(
"computeMeanGapsInContact"
);
return
res
;
}
/* -------------------------------------------------------------------------- */
Real
BemGigi
::
computeG
(){
STARTTIMER
(
"computeG"
);
unsigned
int
n
=
surface
.
size
();
unsigned
int
size
=
n
*
n
;
Real
res
=
0.
;
#pragma omp parallel for reduction(+: res)
for
(
unsigned
int
i
=
0
;
i
<
size
;
++
i
)
{
if
(
this
->
surface_tractions
(
i
)
<=
0
)
continue
;
Real
val
=
this
->
gap
(
i
);
res
+=
val
*
val
;
}
// std::cout << res << std::endl;
STOPTIMER
(
"computeG"
);
return
res
;
}
/* -------------------------------------------------------------------------- */
Real
BemGigi
::
computeTau
(){
STARTTIMER
(
"computeTau"
);
surface_spectral_r
=
surface_t
;
surface_spectral_r
.
FFTTransform
(
nthreads
);
unsigned
int
n
=
surface
.
size
();
unsigned
int
size
=
n
*
n
;
#pragma omp parallel for
for
(
unsigned
int
i
=
0
;
i
<
size
;
++
i
)
{
surface_spectral_r
(
i
)
*=
this
->
surface_spectral_influence_disp
(
i
);
}
surface_spectral_r
.
FFTITransform
(
surface_r
,
nthreads
);
Real
rbar
=
0
;
UInt
nb_contact
=
0
;
#pragma omp parallel for reduction(+: nb_contact,rbar)
for
(
unsigned
int
i
=
0
;
i
<
size
;
++
i
)
{
if
(
this
->
surface_tractions
(
i
)
<=
0
)
continue
;
++
nb_contact
;
rbar
+=
surface_r
(
i
);
}
rbar
/=
nb_contact
;
surface_r
-=
rbar
;
Real
tau_sum1
=
0.
,
tau_sum2
=
0.
;
#pragma omp parallel for reduction(+: tau_sum1, tau_sum2)
for
(
unsigned
int
i
=
0
;
i
<
size
;
++
i
)
{
if
(
this
->
surface_tractions
(
i
)
<=
0
)
continue
;
tau_sum1
+=
this
->
gap
(
i
)
*
surface_t
(
i
);
tau_sum2
+=
surface_r
(
i
)
*
surface_t
(
i
);
}
Real
tau
=
tau_sum1
/
tau_sum2
;
STOPTIMER
(
"computeTau"
);
return
tau
;
}
Event Timeline
Log In to Comment