Line data Source code
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 :
38 : #include "lib_algebra/operator/interface/operator_inverse.h"
39 : #include "lib_algebra/operator/interface/linear_operator_inverse.h"
40 : #include "lib_algebra/operator/debug_writer.h"
41 :
42 : // modul intern headers
43 : #include "lib_disc/assemble_interface.h"
44 : #include "lib_disc/operator/non_linear_operator/assembled_non_linear_operator.h"
45 : #include "lib_disc/operator/linear_operator/assembled_linear_operator.h"
46 : #include "lib_disc/assemble_interface.h"
47 :
48 : #include "lib_disc/spatial_disc/domain_disc.h"
49 : #include "lib_disc/function_spaces/metric_spaces.h"
50 : #include "lib_disc/function_spaces/error_elem_marking_strategy.h"
51 :
52 : #include "lib_algebra/operator/debug_writer.h"
53 :
54 : namespace ug {
55 :
56 : //! Nested iteration solver (e.g. for full multigrid)
57 : /*! If error estimators are available and the convergence rate of the FE/FV method is known,
58 : * these methods allow to construct an optimal order method.
59 : * */
60 : template <typename TDomain, typename TAlgebra>
61 : class NestedIterationSolver
62 : : public IOperatorInverse<typename TAlgebra::vector_type>,
63 : public DebugWritingObject<TAlgebra>
64 : {
65 : public:
66 : /// Algebra type
67 : typedef TAlgebra algebra_type;
68 :
69 : /// Vector type
70 : typedef typename TAlgebra::vector_type vector_type;
71 :
72 : /// Matrix type
73 : typedef typename TAlgebra::matrix_type matrix_type;
74 :
75 : /// GridFunction type
76 : typedef GridFunction<TDomain, TAlgebra> grid_function_type;
77 :
78 : /// Error function type
79 : typedef GridFunction<TDomain, CPUAlgebra> error_function_type;
80 :
81 : public:
82 : /// default constructor
83 : NestedIterationSolver();
84 :
85 : /// constructor setting operator
86 : NestedIterationSolver(SmartPtr<IOperator<vector_type> > N);
87 :
88 : /// constructor using assembling
89 : NestedIterationSolver(SmartPtr<IAssemble<TAlgebra> > spAss);
90 :
91 : /// constructor using assembling
92 : NestedIterationSolver(SmartPtr<IAssemble<TAlgebra> > spAss, SmartPtr<IAssemble<TAlgebra> > spDomErr);
93 :
94 : /// constructor
95 : NestedIterationSolver(SmartPtr<ILinearOperatorInverse<vector_type> > LinearSolver);
96 :
97 : /// sets the linear solver (this should be a fixed number of multi-grid cycles )
98 0 : void set_linear_solver(SmartPtr<ILinearOperatorInverse<vector_type> > LinearSolver) {m_spLinearSolver = LinearSolver;}
99 :
100 : /** @name IOperatorInverse interface*/
101 : ///@{
102 : /// This operator inverts the Operator N: Y -> X
103 : virtual bool init(SmartPtr<IOperator<vector_type> > N);
104 :
105 : /// prepare Operator
106 : virtual bool prepare(vector_type& u);
107 :
108 : /// apply Operator, i.e. N^{-1}(0) = u
109 : virtual bool apply(vector_type& u);
110 : ///@}
111 :
112 :
113 : /// returns information about configuration parameters
114 : /**
115 : * this should return necessary information about parameters and possibly
116 : * calling config_string of subcomponents.
117 : *
118 : * \returns std::string necessary information about configuration parameters
119 : */
120 :
121 : virtual std::string config_string() const;
122 :
123 :
124 :
125 :
126 : /** @name Adaptive refinement
127 : * for adaptive nested iteration
128 : */
129 : ///@{
130 :
131 : //! getter/setter for top level
132 0 : int top_level() const {return m_topLevel;}
133 0 : void set_top_level(int lev) { m_topLevel = lev;}
134 :
135 : //! getter/setter for base level
136 0 : int base_level() const {return m_baseLevel;}
137 0 : void set_base_level(int lev) { m_baseLevel = lev;}
138 :
139 : //! set grid refiner
140 0 : void set_refiner(SmartPtr<IRefiner> r) {m_spRefiner =r; }
141 :
142 : //! set marking strategy
143 0 : void set_refinement_marking(SmartPtr<IElementMarkingStrategy<TDomain> > m) {m_spRefinementMarking =m; }
144 : void set_coarsening_marking(SmartPtr<IElementMarkingStrategy<TDomain> > m) {m_spCoarseningMarking =m; }
145 0 : void set_tolerance(number tol)
146 : {
147 0 : m_TOL = tol;
148 0 : UG_LOG("NI::set_tolerance" << m_TOL <<std::endl);
149 0 : }
150 :
151 :
152 :
153 : //! indicates, if grids should be refined further (if desired accuracy has not been achieved)
154 0 : bool use_adaptive_refinement() const {return m_bUseAdaptiveRefinement;}
155 : //! disable grid refinement
156 0 : void disable_adaptive_refinement() {m_bUseAdaptiveRefinement= false;}
157 : //! enable grid refinement
158 0 : void enable_adaptive_refinement() {m_bUseAdaptiveRefinement= true;}
159 :
160 :
161 0 : void set_max_steps(int steps) {m_maxSteps = steps;}
162 0 : number last_error() const {return m_lastError;}
163 :
164 0 : void set_associated_space(SmartPtr<IGridFunctionSpace<grid_function_type> > spSpace)
165 0 : { m_spAssociatedSpace = spSpace;}
166 :
167 0 : void set_absolute_tolerance(number atol)
168 : {
169 :
170 0 : m_absTOL = atol;
171 0 : UG_LOG("NI::set_absolute_tolerance" << m_absTOL <<std::endl);
172 0 : }
173 :
174 : ///@}
175 :
176 : /// for debug output
177 : using DebugWritingObject<TAlgebra>::write_debug;
178 : using DebugWritingObject<TAlgebra>::vector_debug_writer;
179 0 : void set_debug_elem_error(SmartPtr<error_function_type> spErrEta)
180 0 : {m_spElemError = spErrEta;}
181 :
182 : protected:
183 : void estimate_and_mark_domain(const grid_function_type& u, SmartPtr<IElementMarkingStrategy<TDomain> > spMarking, bool bClearMarks = true);
184 : number refine_domain(const grid_function_type& u);
185 : number coarsen_domain(const grid_function_type& u);
186 :
187 :
188 : private:
189 : /// linear solver
190 : SmartPtr<ILinearOperatorInverse<vector_type> > m_spLinearSolver;
191 :
192 : /// assembling routine
193 : SmartPtr<AssembledOperator<algebra_type> > m_N;
194 : /// jacobi operator
195 : SmartPtr<AssembledLinearOperator<algebra_type> > m_J;
196 : /// assembling
197 : SmartPtr<IAssemble<TAlgebra> > m_spAss;
198 : SmartPtr<IAssemble<TAlgebra> > m_spDomErr;
199 :
200 : /// call counter
201 : int m_dgbCall;
202 : int m_lastNumSteps;
203 :
204 :
205 : int m_baseLevel;
206 : int m_topLevel;
207 :
208 : /// (adaptive) refinement
209 : SmartPtr<IRefiner> m_spRefiner;
210 :
211 : SmartPtr<IElementMarkingStrategy<TDomain> > m_spRefinementMarking;
212 : SmartPtr<IElementMarkingStrategy<TDomain> > m_spCoarseningMarking;
213 : bool m_bUseAdaptiveRefinement;
214 :
215 : int m_maxSteps;
216 : number m_lastError;
217 :
218 : /// tolerance
219 : number m_TOL;
220 :
221 : /// associated norm (for relative error)
222 : SmartPtr<IGridFunctionSpace<grid_function_type> > m_spAssociatedSpace;
223 :
224 : /// absolute tolerance (for relative error)
225 : number m_absTOL;
226 :
227 : /// (optional) debug info for adaptive refinement
228 : SmartPtr<error_function_type> m_spElemError;
229 : };
230 :
231 : }
232 :
233 : #include "nested_iteration_impl.h"
234 :
235 : #endif /* __H__UG__LIB_DISC__OPERATOR__NON_LINEAR_OPERATOR__NI_SOLVER__NI__ */
|