Skip to content
Snippets Groups Projects
Commit 60870867 authored by Oliver Sander's avatar Oliver Sander Committed by sander@PCPOOL.MI.FU-BERLIN.DE
Browse files

RodSolver now accepts a generic GeodesicFEAssembler instead of only a RodAssembler

[[Imported from SVN: r4037]]
parent 05f035da
No related branches found
No related tags found
No related merge requests found
...@@ -27,7 +27,7 @@ ...@@ -27,7 +27,7 @@
template <class GridType> template <class GridType>
void RodSolver<GridType>::setup(const GridType& grid, void RodSolver<GridType>::setup(const GridType& grid,
const RodAssembler<GridType>* rodAssembler, const GeodesicFEAssembler<typename GridType::LeafGridView, RigidBodyMotion<3> >* assembler,
const SolutionType& x, const SolutionType& x,
const Dune::BitSetVector<blocksize>& dirichletNodes, const Dune::BitSetVector<blocksize>& dirichletNodes,
double tolerance, double tolerance,
...@@ -45,7 +45,7 @@ void RodSolver<GridType>::setup(const GridType& grid, ...@@ -45,7 +45,7 @@ void RodSolver<GridType>::setup(const GridType& grid,
using namespace Dune; using namespace Dune;
grid_ = &grid; grid_ = &grid;
rodAssembler_ = rodAssembler; assembler_ = assembler;
x_ = x; x_ = x;
tolerance_ = tolerance; tolerance_ = tolerance;
maxTrustRegionSteps_ = maxTrustRegionSteps; maxTrustRegionSteps_ = maxTrustRegionSteps;
...@@ -125,7 +125,7 @@ void RodSolver<GridType>::setup(const GridType& grid, ...@@ -125,7 +125,7 @@ void RodSolver<GridType>::setup(const GridType& grid,
hessianMatrix_ = new MatrixType; hessianMatrix_ = new MatrixType;
MatrixIndexSet indices(grid_->size(1), grid_->size(1)); MatrixIndexSet indices(grid_->size(1), grid_->size(1));
rodAssembler_->getNeighborsPerVertex(indices); assembler_->getNeighborsPerVertex(indices);
indices.exportIdx(*hessianMatrix_); indices.exportIdx(*hessianMatrix_);
// ////////////////////////////////////////////////////////// // //////////////////////////////////////////////////////////
...@@ -180,7 +180,7 @@ void RodSolver<GridType>::solve() ...@@ -180,7 +180,7 @@ void RodSolver<GridType>::solve()
std::cout << "----------------------------------------------------" << std::endl; std::cout << "----------------------------------------------------" << std::endl;
std::cout << " Trust-Region Step Number: " << i std::cout << " Trust-Region Step Number: " << i
<< ", radius: " << trustRegion.radius() << ", radius: " << trustRegion.radius()
<< ", energy: " << rodAssembler_->computeEnergy(x_) << std::endl; << ", energy: " << assembler_->computeEnergy(x_) << std::endl;
std::cout << "----------------------------------------------------" << std::endl; std::cout << "----------------------------------------------------" << std::endl;
} }
...@@ -188,9 +188,8 @@ void RodSolver<GridType>::solve() ...@@ -188,9 +188,8 @@ void RodSolver<GridType>::solve()
CorrectionType corr(x_.size()); CorrectionType corr(x_.size());
corr = 0; corr = 0;
rodAssembler_->assembleGradient(x_, rhs); assembler_->assembleGradient(x_, rhs);
rodAssembler_->assembleMatrix(x_, *hessianMatrix_); assembler_->assembleMatrix(x_, *hessianMatrix_);
//rodAssembler_->assembleMatrixFD(x_, *hessianMatrix_);
//gradientFDCheck(x_, rhs, *rodAssembler_); //gradientFDCheck(x_, rhs, *rodAssembler_);
//hessianFDCheck(x_, *hessianMatrix_, *rodAssembler_); //hessianFDCheck(x_, *hessianMatrix_, *rodAssembler_);
...@@ -324,8 +323,8 @@ void RodSolver<GridType>::solve() ...@@ -324,8 +323,8 @@ void RodSolver<GridType>::solve()
} }
/** \todo Don't always recompute oldEnergy */ /** \todo Don't always recompute oldEnergy */
double oldEnergy = rodAssembler_->computeEnergy(x_); double oldEnergy = assembler_->computeEnergy(x_);
double energy = rodAssembler_->computeEnergy(newIterate); double energy = assembler_->computeEnergy(newIterate);
// compute the model decrease // compute the model decrease
// It is $ m(x) - m(x+s) = -<g,s> - 0.5 <s, Hs> // It is $ m(x) - m(x+s) = -<g,s> - 0.5 <s, Hs>
......
...@@ -12,7 +12,7 @@ ...@@ -12,7 +12,7 @@
#include <dune-solvers/norms/h1seminorm.hh> #include <dune-solvers/norms/h1seminorm.hh>
#include <dune-solvers/solvers/loopsolver.hh> #include <dune-solvers/solvers/loopsolver.hh>
#include "rodassembler.hh" #include "geodesicfeassembler.hh"
#include "rigidbodymotion.hh" #include "rigidbodymotion.hh"
...@@ -38,7 +38,7 @@ public: ...@@ -38,7 +38,7 @@ public:
{} {}
void setup(const GridType& grid, void setup(const GridType& grid,
const RodAssembler<GridType>* rodAssembler, const GeodesicFEAssembler<typename GridType::LeafGridView, RigidBodyMotion<3> >* rodAssembler,
const SolutionType& x, const SolutionType& x,
const Dune::BitSetVector<blocksize>& dirichletNodes, const Dune::BitSetVector<blocksize>& dirichletNodes,
double tolerance, double tolerance,
...@@ -99,7 +99,7 @@ protected: ...@@ -99,7 +99,7 @@ protected:
MatrixType* hessianMatrix_; MatrixType* hessianMatrix_;
/** \brief The assembler for the material law */ /** \brief The assembler for the material law */
const RodAssembler<GridType>* rodAssembler_; const GeodesicFEAssembler<typename GridType::LeafGridView, RigidBodyMotion<3> >* assembler_;
/** \brief The multigrid solver */ /** \brief The multigrid solver */
LoopSolver<CorrectionType>* mmgSolver_; LoopSolver<CorrectionType>* mmgSolver_;
......
0% Loading or .
You are about to add 0 people to the discussion. Proceed with caution.
Please register or to comment