ug4
nested_iteration.h
Go to the documentation of this file.
1 /*
2  * Copyright (c) 2010-2015: G-CSC, Goethe University Frankfurt
3  * Author: Arne Naegel
4  *
5  * This file is part of UG4.
6  *
7  * UG4 is free software: you can redistribute it and/or modify it under the
8  * terms of the GNU Lesser General Public License version 3 (as published by the
9  * Free Software Foundation) with the following additional attribution
10  * requirements (according to LGPL/GPL v3 §7):
11  *
12  * (1) The following notice must be displayed in the Appropriate Legal Notices
13  * of covered and combined works: "Based on UG4 (www.ug4.org/license)".
14  *
15  * (2) The following notice must be displayed at a prominent place in the
16  * terminal output of covered works: "Based on UG4 (www.ug4.org/license)".
17  *
18  * (3) The following bibliography is recommended for citation and must be
19  * preserved in all covered files:
20  * "Reiter, S., Vogel, A., Heppner, I., Rupp, M., and Wittum, G. A massively
21  * parallel geometric multigrid solver on hierarchically distributed grids.
22  * Computing and visualization in science 16, 4 (2013), 151-164"
23  * "Vogel, A., Reiter, S., Rupp, M., Nägel, A., and Wittum, G. UG4 -- a novel
24  * flexible software system for simulating pde based models on high performance
25  * computers. Computing and visualization in science 16, 4 (2013), 165-179"
26  *
27  * This program is distributed in the hope that it will be useful,
28  * but WITHOUT ANY WARRANTY; without even the implied warranty of
29  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
30  * GNU Lesser General Public License for more details.
31  */
32 
33 #ifndef __H__UG__LIB_DISC__OPERATOR__NON_LINEAR_OPERATOR__NI_SOLVER__NI__
34 #define __H__UG__LIB_DISC__OPERATOR__NON_LINEAR_OPERATOR__NI_SOLVER__NI__
35 
36 #include <cmath>
37 
41 
42 // modul intern headers
47 
51 
53 
54 namespace ug {
55 
57 
60 template <typename TDomain, typename TAlgebra>
62  : public IOperatorInverse<typename TAlgebra::vector_type>,
63  public DebugWritingObject<TAlgebra>
64 {
65  public:
67  typedef TAlgebra algebra_type;
68 
70  typedef typename TAlgebra::vector_type vector_type;
71 
73  typedef typename TAlgebra::matrix_type matrix_type;
74 
77 
80 
81  public:
84 
87 
90 
93 
96 
99 
103  virtual bool init(SmartPtr<IOperator<vector_type> > N);
104 
106  virtual bool prepare(vector_type& u);
107 
109  virtual bool apply(vector_type& u);
111 
112 
114 
121  virtual std::string config_string() const;
122 
123 
124 
125 
130 
132  int top_level() const {return m_topLevel;}
133  void set_top_level(int lev) { m_topLevel = lev;}
134 
136  int base_level() const {return m_baseLevel;}
137  void set_base_level(int lev) { m_baseLevel = lev;}
138 
141 
146  {
147  m_TOL = tol;
148  UG_LOG("NI::set_tolerance" << m_TOL <<std::endl);
149  }
150 
151 
152 
159 
160 
161  void set_max_steps(int steps) {m_maxSteps = steps;}
162  number last_error() const {return m_lastError;}
163 
165  { m_spAssociatedSpace = spSpace;}
166 
168  {
169 
170  m_absTOL = atol;
171  UG_LOG("NI::set_absolute_tolerance" << m_absTOL <<std::endl);
172  }
173 
175 
180  {m_spElemError = spErrEta;}
181 
182  protected:
183  void estimate_and_mark_domain(const grid_function_type& u, SmartPtr<IElementMarkingStrategy<TDomain> > spMarking, bool bClearMarks = true);
186 
187 
188  private:
191 
199 
203 
204 
207 
210 
214 
217 
220 
223 
226 
229 };
230 
231 }
232 
233 #include "nested_iteration_impl.h"
234 
235 #endif /* __H__UG__LIB_DISC__OPERATOR__NON_LINEAR_OPERATOR__NI_SOLVER__NI__ */
Definition: smart_pointer.h:108
Definition: debug_writer.h:354
represents numerical solutions on a grid using an algebraic vector
Definition: grid_function.h:121
Interface providing Jacobian and Defect of a discretization.
Definition: assemble_interface.h:110
Abstract base class for element marking (in adaptive refinement)
Definition: error_elem_marking_strategy.h:125
Abstract base class for grid functions.
Definition: metric_spaces.h:82
describes an inverse linear mapping X->Y
Definition: linear_operator_inverse.h:80
describes a mapping X->Y
Definition: operator.h:86
describes an inverse mapping X->Y
Definition: operator_inverse.h:79
linear solver using abstract preconditioner interface
Definition: linear_solver.h:57
Nested iteration solver (e.g. for full multigrid)
Definition: nested_iteration.h:64
SmartPtr< AssembledOperator< algebra_type > > m_N
assembling routine
Definition: nested_iteration.h:193
virtual std::string config_string() const
returns information about configuration parameters
Definition: nested_iteration_impl.h:365
SmartPtr< IAssemble< TAlgebra > > m_spDomErr
Definition: nested_iteration.h:198
bool use_adaptive_refinement() const
indicates, if grids should be refined further (if desired accuracy has not been achieved)
Definition: nested_iteration.h:154
number m_TOL
tolerance
Definition: nested_iteration.h:219
SmartPtr< IElementMarkingStrategy< TDomain > > m_spCoarseningMarking
Definition: nested_iteration.h:212
void set_coarsening_marking(SmartPtr< IElementMarkingStrategy< TDomain > > m)
getter/setter for top level
Definition: nested_iteration.h:144
void set_max_steps(int steps)
getter/setter for top level
Definition: nested_iteration.h:161
SmartPtr< IRefiner > m_spRefiner
(adaptive) refinement
Definition: nested_iteration.h:209
TAlgebra algebra_type
Algebra type.
Definition: nested_iteration.h:67
void set_linear_solver(SmartPtr< ILinearOperatorInverse< vector_type > > LinearSolver)
sets the linear solver (this should be a fixed number of multi-grid cycles )
Definition: nested_iteration.h:98
void set_base_level(int lev)
getter/setter for top level
Definition: nested_iteration.h:137
int base_level() const
getter/setter for base level
Definition: nested_iteration.h:136
SmartPtr< IGridFunctionSpace< grid_function_type > > m_spAssociatedSpace
associated norm (for relative error)
Definition: nested_iteration.h:222
void set_debug_elem_error(SmartPtr< error_function_type > spErrEta)
Definition: nested_iteration.h:179
number last_error() const
getter/setter for top level
Definition: nested_iteration.h:162
GridFunction< TDomain, TAlgebra > grid_function_type
GridFunction type.
Definition: nested_iteration.h:76
virtual bool apply(vector_type &u)
apply Operator, i.e. N^{-1}(0) = u
Definition: nested_iteration_impl.h:208
number m_absTOL
absolute tolerance (for relative error)
Definition: nested_iteration.h:225
number refine_domain(const grid_function_type &u)
void set_refiner(SmartPtr< IRefiner > r)
set grid refiner
Definition: nested_iteration.h:140
SmartPtr< IElementMarkingStrategy< TDomain > > m_spRefinementMarking
Definition: nested_iteration.h:211
void set_absolute_tolerance(number atol)
getter/setter for top level
Definition: nested_iteration.h:167
void disable_adaptive_refinement()
disable grid refinement
Definition: nested_iteration.h:156
void set_refinement_marking(SmartPtr< IElementMarkingStrategy< TDomain > > m)
set marking strategy
Definition: nested_iteration.h:143
SmartPtr< error_function_type > m_spElemError
(optional) debug info for adaptive refinement
Definition: nested_iteration.h:228
int m_baseLevel
Definition: nested_iteration.h:205
SmartPtr< IAssemble< TAlgebra > > m_spAss
assembling
Definition: nested_iteration.h:197
SmartPtr< AssembledLinearOperator< algebra_type > > m_J
jacobi operator
Definition: nested_iteration.h:195
GridFunction< TDomain, CPUAlgebra > error_function_type
Error function type.
Definition: nested_iteration.h:79
int m_lastNumSteps
Definition: nested_iteration.h:202
number coarsen_domain(const grid_function_type &u)
Coarsen all elements.
Definition: nested_iteration_impl.h:192
void set_tolerance(number tol)
getter/setter for top level
Definition: nested_iteration.h:145
SmartPtr< ILinearOperatorInverse< vector_type > > m_spLinearSolver
linear solver
Definition: nested_iteration.h:190
TAlgebra::matrix_type matrix_type
Matrix type.
Definition: nested_iteration.h:73
virtual bool init(SmartPtr< IOperator< vector_type > > N)
Definition: nested_iteration_impl.h:134
int m_dgbCall
call counter
Definition: nested_iteration.h:201
void set_top_level(int lev)
getter/setter for top level
Definition: nested_iteration.h:133
void estimate_and_mark_domain(const grid_function_type &u, SmartPtr< IElementMarkingStrategy< TDomain > > spMarking, bool bClearMarks=true)
Refines domain and provides error estimate.
Definition: nested_iteration_impl.h:157
void enable_adaptive_refinement()
enable grid refinement
Definition: nested_iteration.h:158
void set_associated_space(SmartPtr< IGridFunctionSpace< grid_function_type > > spSpace)
getter/setter for top level
Definition: nested_iteration.h:164
TAlgebra::vector_type vector_type
Vector type.
Definition: nested_iteration.h:70
NestedIterationSolver()
default constructor
Definition: nested_iteration_impl.h:71
bool m_bUseAdaptiveRefinement
Definition: nested_iteration.h:213
int top_level() const
getter/setter for top level
Definition: nested_iteration.h:132
int m_topLevel
Definition: nested_iteration.h:206
virtual bool prepare(vector_type &u)
prepare Operator
Definition: nested_iteration_impl.h:147
int m_maxSteps
Definition: nested_iteration.h:215
number m_lastError
Definition: nested_iteration.h:216
#define UG_LOG(msg)
Definition: log.h:367
double number
Definition: types.h:124
CPUAlgebra::vector_type vector_type
the ug namespace