Page MenuHomec4science

bem_uzawa.cpp
No OneTemporary

File Metadata

Created
Sun, May 19, 13:57

bem_uzawa.cpp

/**
*
* @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_uzawa.hh"
#include <iostream>
#include <sstream>
#include <fstream>
#include <iomanip>
#include <sstream>
/* -------------------------------------------------------------------------- */
#define TIMER
#include "surface_timer.hh"
/* -------------------------------------------------------------------------- */
__BEGIN_TAMAAS__
void BemUzawa::computeTractionsFromDisplacements() {
this->applyInverseInfluenceFunctions(this->true_displacements,this->surface_tractions);
}
Real BemUzawa::computeEquilibrium(Real epsilon,
Real mean_displacement, Real penalization_parameter) {
this->computeSpectralInfluenceOverDisplacement();
this->search_direction = 0.;
this->true_displacements = 0;
this->true_displacements = mean_displacement;
this->computeGaps();
this->multipliers = 0;
Real f_2 = 1e300;
Real new_pen = 1e300;
Real old_pen = 1e300;
convergence_iterations.clear();
UInt nb_iterations_1 = 0;
UInt nb_iterations_2 = 0;
std::ofstream file ("output.txt");
Real max_iterations_1 = 500;
Real max_iterations_2 = 1000;
while (new_pen> epsilon && nb_iterations_1++ < max_iterations_1) {
std::cout << "iter ext "<< nb_iterations_1 << std::endl;
while (f_2> 1e-12 && nb_iterations_2++ < max_iterations_2) {
std::cout << "iter int "<< nb_iterations_2 << std::endl;
this->functional->computeGradFU();
this->computeSearchDirection(mean_displacement,penalization_parameter);
Real alpha=1/(10*penalization_parameter);
this->old_displacements=this->true_displacements;
this->updateUnknown( alpha, mean_displacement);
this->computeGaps();
f_2 = computeStoppingCriterion();
if (nb_iterations_2 % dump_freq == 0) {
std::cout << std::scientific << std::setprecision(10)
<< nb_iterations << " "
<< f_2 << std::fixed << std::endl;
this->computePressures();
Real orth=0.;
UInt n = surface.size();
UInt size = n*n;
#pragma omp parallel for reduction(+:orth)
for (UInt i = 0; i < size; ++i) {
orth += std::abs(surface_tractions(i) * ( this->true_displacements(i)-surface(i)));
}
file<< std::scientific << std::setprecision(10)<< nb_iterations_2 << " "
<<f_2<< " "<<orth<< std::endl;
}
}
this->updateMultipliers(penalization_parameter);
old_pen=computeInterpenetration(this->old_displacements);
new_pen=computeInterpenetration(this->true_displacements);
std::cout << "old penetration is "<< old_pen << std::endl;
std::cout << "new penetration is "<< new_pen << std::endl;
penalization_parameter=this->updatePenalization(penalization_parameter, old_pen, new_pen);
//to avoid ill conditioned system
if (penalization_parameter>1.0e8) new_pen=0.;
nb_iterations_2 = 0;
f_2 = 1e300;
std::cout << "penalization is "<< penalization_parameter << std::endl;
}
this->computeGaps();
this->computePressures();
return new_pen;
}
/* -------------------------------------------------------------------------- */
Real BemUzawa::computeStoppingCriterion() {
Real crit=0.;
Real disp_norm = 0.;
UInt n = surface.size();
UInt size = n*n;
#pragma omp parallel for reduction(+:crit, disp_norm)
for (UInt i = 0; i < size; ++i) {
crit += (this->search_direction(i))*(this->search_direction(i));
disp_norm += (true_displacements(i)*true_displacements(i));
}
return crit / disp_norm;
}
/* -------------------------------------------------------------------------- */
void BemUzawa::computeSearchDirection(Real mean_displacement,Real penalization_parameter) {
STARTTIMER("computeOptimalStep");
UInt n = surface.size();
UInt size = n*n;
const Surface<Real> & gradF = this->functional->getGradF();
#pragma omp parallel for
for (UInt i = 1; i < size; ++i) {
this->search_direction(i)=gradF(i);
if (this->multipliers(i)+penalization_parameter*gap(i)<=0){
this->search_direction(i)=this->search_direction(i)+ this->multipliers(i)+penalization_parameter*gap(i);
}
}
}
/* -------------------------------------------------------------------------- */
Real BemUzawa::computeOptimalStep() {
STARTTIMER("computeOptimalStep");
this->applyInverseInfluenceFunctions(search_direction, surface_r);
UInt n = surface.size();
UInt size = n*n;
Real numerator = 0., denominator = 0.;
#pragma omp parallel for reduction(+: numerator, denominator)
for (UInt i = 0; i < size; ++i) {
numerator += search_direction(i) * search_direction(i);
denominator += surface_r(i) * search_direction(i);
}
Real alpha = numerator / denominator;
STOPTIMER("computeOptimalStep");
return alpha;
}
/* -------------------------------------------------------------------------- */
void BemUzawa::updateMultipliers(Real penalization_parameter) {
STARTTIMER("updateMultipliers");
UInt n = surface.size();
UInt size = n*n;
#pragma omp parallel for
for (UInt i = 0; i < size; ++i) {
if (this->multipliers(i)+penalization_parameter*gap(i)>0.){
this->multipliers(i) = 0;
}
else{
this->multipliers(i)+=penalization_parameter*gap(i);
}
}
STOPTIMER("updateMultipliers");
}
/* -------------------------------------------------------------------------- */
Real BemUzawa::updatePenalization(Real penalization_parameter, Real old_pen, Real new_pen) {
STARTTIMER("updatePenalization");
if (new_pen>0.25*old_pen){
penalization_parameter*=5;
}
return penalization_parameter;
STOPTIMER("updatePenalization");
}
/* -------------------------------------------------------------------------- */
Real BemUzawa::computeInterpenetration( Surface<Real> & displacements ) {
STARTTIMER("updateMultipliers");
UInt n = surface.size();
UInt size = n*n;
Real res=0.;
#pragma omp parallel for
for (UInt i = 0; i < size; ++i) {
if (0>this->true_displacements(i)-surface(i)){
res+=(this->true_displacements(i)-surface(i))*(this->true_displacements(i)-surface(i));
}
}
return sqrt(res);
STOPTIMER("updateMultipliers");
}
/* -------------------------------------------------------------------------- */
void BemUzawa::updateUnknown(Real alpha, Real mean_displacement) {
STARTTIMER("updateDisplacements");
UInt n = surface.size();
UInt size = n*n;
#pragma omp parallel for
for (UInt i = 0; i < size; ++i) {
this->true_displacements(i) -= alpha*this->search_direction(i);
}
Real moyenne=SurfaceStatistics::computeAverage(this->true_displacements, 0);
for (UInt i = 0; i < size; ++i) {
this->true_displacements(i) =this->true_displacements(i)-moyenne+mean_displacement;
}
STOPTIMER("updateDisplacements");
}
/* -------------------------------------------------------------------------- */
void BemUzawa::computePressures() {
this->computeTractionsFromDisplacements();
this->functional->computeGradFU();
//const Surface<Real> & gradF = this->functional->getGradF();
//Real min = SurfaceStatistics::computeMinimum(gradF);
this->surface_tractions -= this->surface_tractions(0);
}
/* -------------------------------------------------------------------------- */
__END_TAMAAS__

Event Timeline