1 #ifndef DUNE_PDELAB_STATIONARY_LINEARPROBLEM_HH
2 #define DUNE_PDELAB_STATIONARY_LINEARPROBLEM_HH
6 #include <dune/common/timer.hh>
7 #include <dune/common/parametertree.hh>
24 template<
class RFType>
43 template<
typename GO,
typename LS,
typename V>
46 typedef typename Dune::template FieldTraits<typename V::ElementType >::real_type Real;
47 typedef typename GO::Traits::Jacobian M;
48 typedef typename GO::Traits::TrialGridFunctionSpace TrialGridFunctionSpace;
50 typedef GO GridOperator;
60 , _min_defect(min_defect)
61 , _hanging_node_modifications(false)
71 , _min_defect(min_defect)
72 , _hanging_node_modifications(false)
99 , _reduction(params.get<Real>(
"reduction"))
100 , _min_defect(params.get<Real>(
"min_defect",1
e-99))
101 , _hanging_node_modifications(params.get<bool>(
"hanging_node_modifications",false))
102 , _keep_matrix(params.get<bool>(
"keep_matrix",true))
103 , _verbose(params.get<int>(
"verbosity",1))
128 , _reduction(params.get<Real>(
"reduction"))
129 , _min_defect(params.get<Real>(
"min_defect",1
e-99))
130 , _hanging_node_modifications(params.get<bool>(
"hanging_node_modifications",false))
131 , _keep_matrix(params.get<bool>(
"keep_matrix",true))
132 , _verbose(params.get<int>(
"verbosity",1))
138 _hanging_node_modifications=b;
144 return _hanging_node_modifications;
164 void apply(V& x,
bool reuse_matrix =
false) {
169 void apply (
bool reuse_matrix =
false)
172 double timing,assembler_time=0;
179 _jacobian = std::make_shared<M>(_go);
180 timing = watch.elapsed();
181 if (_go.trialGridFunctionSpace().gridView().comm().rank()==0 && _verbose>=1)
182 std::cout <<
"=== matrix setup (max) " << timing <<
" s" << std::endl;
184 assembler_time += timing;
186 else if (_go.trialGridFunctionSpace().gridView().comm().rank()==0 && _verbose>=1)
187 std::cout <<
"=== matrix setup skipped (matrix already allocated)" << std::endl;
189 if (_hanging_node_modifications)
192 _go.localAssembler().backtransform(*_x);
197 (*_jacobian) = Real(0.0);
198 _go.jacobian(*_x,*_jacobian);
201 timing = watch.elapsed();
203 if (_go.trialGridFunctionSpace().gridView().comm().rank()==0 && _verbose>=1)
206 std::cout <<
"=== matrix assembly SKIPPED" << std::endl;
208 std::cout <<
"=== matrix assembly (max) " << timing <<
" s" << std::endl;
211 assembler_time += timing;
216 W r(_go.testGridFunctionSpace(),0.0);
219 timing = watch.elapsed();
221 if (_go.trialGridFunctionSpace().gridView().comm().rank()==0 && _verbose>=1)
222 std::cout <<
"=== residual assembly (max) " << timing <<
" s" << std::endl;
223 assembler_time += timing;
226 auto defect = _ls.norm(r);
230 V z(_go.trialGridFunctionSpace(),0.0);
231 auto red = std::max(_reduction,_min_defect/defect);
232 if (_go.trialGridFunctionSpace().gridView().comm().rank()==0 && _verbose>=1)
234 std::cout <<
"=== solving (reduction: " << red <<
") ";
236 std::cout << std::flush;
238 std::cout << std::endl;
240 _ls.apply(*_jacobian,z,r,red);
241 _linear_solver_result = _ls.result();
242 timing = watch.elapsed();
244 if (_go.trialGridFunctionSpace().gridView().comm().rank()==0 && _verbose>=1)
245 std::cout << timing <<
" s" << std::endl;
254 _res.
defect =
static_cast<double>(defect)*_linear_solver_result.
reduction;
258 if (_hanging_node_modifications)
261 if (_hanging_node_modifications)
262 _go.localAssembler().backtransform(*_x);
276 return _linear_solver_result;
294 std::shared_ptr<M> _jacobian;
299 bool _hanging_node_modifications;
const Entity & e
Definition: localfunctionspace.hh:121
void set_shifted_dofs(const CG &cg, typename XG::ElementType x, XG &xg)
Definition: constraints.hh:1014
For backward compatibility – Do not use this!
Definition: adaptivity.hh:28
typename impl::BackendVectorSelector< GridFunctionSpace, FieldType >::Type Vector
alias of the return type of BackendVectorSelector
Definition: backend/interface.hh:106
RFType conv_rate
Definition: solver.hh:36
bool converged
Definition: solver.hh:32
RFType reduction
Definition: solver.hh:35
unsigned int iterations
Definition: solver.hh:33
double elapsed
Definition: solver.hh:34
Definition: linearproblem.hh:26
int linear_solver_iterations
Definition: linearproblem.hh:31
double linear_solver_time
Definition: linearproblem.hh:30
StationaryLinearProblemSolverResult()
Definition: linearproblem.hh:33
RFType defect
Definition: linearproblem.hh:28
double assembler_time
Definition: linearproblem.hh:29
RFType first_defect
Definition: linearproblem.hh:27
Definition: linearproblem.hh:45
void setKeepMatrix(bool b)
Set whether the jacobian matrix should be kept across calls to apply().
Definition: linearproblem.hh:148
Real reduction() const
Definition: linearproblem.hh:279
void apply(bool reuse_matrix=false)
Definition: linearproblem.hh:169
StationaryLinearProblemSolver(const GO &go, LS &ls, Real reduction, Real min_defect=1e-99, int verbose=1)
Definition: linearproblem.hh:66
void setHangingNodeModifications(bool b)
Set whether the solver should apply the necessary transformations for calculations on hanging nodes.
Definition: linearproblem.hh:136
bool hangingNodeModifications() const
Return whether the solver performs the necessary transformations for calculations on hanging nodes.
Definition: linearproblem.hh:142
StationaryLinearProblemSolver(const GO &go, LS &ls, V &x, const ParameterTree ¶ms)
Construct a StationaryLinearProblemSolver for the given objects and read parameters from a ParameterT...
Definition: linearproblem.hh:95
StationaryLinearProblemSolver(const GO &go, LS &ls, V &x, Real reduction, Real min_defect=1e-99, int verbose=1)
Definition: linearproblem.hh:55
StationaryLinearProblemSolverResult< double > Result
Definition: linearproblem.hh:53
const Dune::PDELab::LinearSolverResult< double > & ls_result() const
Definition: linearproblem.hh:275
void apply(V &x, bool reuse_matrix=false)
Definition: linearproblem.hh:164
void discardMatrix()
Discard the stored Jacobian matrix.
Definition: linearproblem.hh:269
void setReduction(Real reduction)
Definition: linearproblem.hh:284
const Result & result() const
Definition: linearproblem.hh:159
bool keepMatrix() const
Return whether the jacobian matrix is kept across calls to apply().
Definition: linearproblem.hh:154
StationaryLinearProblemSolver(const GO &go, LS &ls, const ParameterTree ¶ms)
Construct a StationaryLinearProblemSolver for the given objects and read parameters from a ParameterT...
Definition: linearproblem.hh:124