Skip to content
Snippets Groups Projects
Commit f48b1b19 authored by Klaus Böhnlein's avatar Klaus Böhnlein Committed by Sander, Oliver
Browse files

Add setup method to Riemannian TR and PN methods that gets a ParameterSet object

parent 4559aa1d
No related branches found
No related tags found
1 merge request!130Add setup method to Riemannian TR and PN methods
......@@ -27,6 +27,24 @@
#include <dune/gfe/parallel/vectorcommunicator.hh>
#endif
template <class Basis, class TargetSpace, class Assembler>
void RiemannianProximalNewtonSolver<Basis, TargetSpace, Assembler>::
setup(const GridType& grid,
const Assembler* assembler,
const SolutionType& x,
const Dune::BitSetVector<blocksize>& dirichletNodes,
const Dune::ParameterTree& parameterSet)
{
setup(grid,
assembler,
x,
dirichletNodes,
parameterSet.get<double>("tolerance"),
parameterSet.get<int>("maxProximalNewtonSteps "),
parameterSet.get<double>("initialRegularization"),
parameterSet.get<bool>("instrumented", 0));
}
template <class Basis, class TargetSpace, class Assembler>
void RiemannianProximalNewtonSolver<Basis,TargetSpace,Assembler>::
setup(const GridType& grid,
......@@ -52,7 +70,7 @@ setup(const GridType& grid,
//////////////////////////////////////////////////////////////////
// Create global numbering for matrix and vector transfer
//////////////////////////////////////////////////////////////////
globalMapper_ = std::make_unique<GlobalMapper>(grid_->leafGridView());
// Transfer all Dirichlet data to the master processor
VectorCommunicator<GlobalMapper, typename GridType::LeafGridView::CollectiveCommunication, Dune::BitSetVector<blocksize> > vectorComm(*globalMapper_,
......
......@@ -4,6 +4,7 @@
#include <vector>
#include <dune/common/bitsetvector.hh>
#include <dune/common/parametertree.hh>
#include <dune/istl/bcrsmatrix.hh>
#include <dune/istl/bvector.hh>
......@@ -81,6 +82,13 @@ public:
double initialRegularization,
bool instrumented);
/** \brief Set up the solver using a monotone multigrid method as the inner solver */
void setup(const GridType& grid,
const Assembler* assembler,
const SolutionType& x,
const Dune::BitSetVector<blocksize>& dirichletNodes,
const Dune::ParameterTree& parameterSet);
void setScaling(const Dune::FieldVector<double,blocksize>& scaling)
{
scaling_ = scaling;
......
......@@ -28,6 +28,29 @@
#include <dune/gfe/parallel/vectorcommunicator.hh>
#endif
template <class Basis, class TargetSpace, class Assembler>
void RiemannianTrustRegionSolver<Basis, TargetSpace, Assembler>::
setup(const GridType& grid,
const Assembler* assembler,
const SolutionType& x,
const Dune::BitSetVector<blocksize>& dirichletNodes,
const Dune::ParameterTree& parameterSet)
{
setup(grid,
assembler,
x,
dirichletNodes,
parameterSet.get<double>("tolerance"),
parameterSet.get<int>("maxTrustRegionSteps"),
parameterSet.get<double>("initialTrustRegionRadius"),
parameterSet.get<int>("numIt"),
parameterSet.get<double>("mgTolerance"),
parameterSet.get<int>("mu"), parameterSet.get<int>("nu1"), parameterSet.get<int>("nu2"),
parameterSet.get<int>("baseIt"),
parameterSet.get<double>("baseTolerance"),
parameterSet.get<bool>("instrumented", 0));
}
template <class Basis, class TargetSpace, class Assembler>
void RiemannianTrustRegionSolver<Basis,TargetSpace,Assembler>::
setup(const GridType& grid,
......@@ -692,7 +715,7 @@ void RiemannianTrustRegionSolver<Basis,TargetSpace,Assembler>::solve()
fwrite(&x_[j], sizeof(TargetSpace), 1, fpIterate);
fclose(fpIterate);
}
if (rank==0)
......
......@@ -4,6 +4,7 @@
#include <vector>
#include <dune/common/bitsetvector.hh>
#include <dune/common/parametertree.hh>
#include <dune/istl/bcrsmatrix.hh>
#include <dune/istl/bvector.hh>
......@@ -119,6 +120,13 @@ public:
double baseTolerance,
bool instrumented);
/** \brief Set up the solver using a monotone multigrid method as the inner solver */
void setup(const GridType& grid,
const Assembler* assembler,
const SolutionType& x,
const Dune::BitSetVector<blocksize>& dirichletNodes,
const Dune::ParameterTree& parameterSet);
void setScaling(const Dune::FieldVector<double,blocksize>& scaling)
{
scaling_ = scaling;
......
0% Loading or .
You are about to add 0 people to the discussion. Proceed with caution.
Please register or to comment