Line data Source code
1 : /*
2 : * Copyright (c) 2010-2015: G-CSC, Goethe University Frankfurt
3 : * Author: Andreas Vogel
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__SPATIAL_DISC__CONSTRAINTS__CONSTRAINTS__INTERFACE__
34 : #define __H__UG__LIB_DISC__SPATIAL_DISC__CONSTRAINTS__CONSTRAINTS__INTERFACE__
35 :
36 : // extern headers
37 : #include <vector>
38 :
39 : // intern headers
40 : #include "lib_disc/assemble_interface.h"
41 : #include "lib_disc/common/local_algebra.h"
42 : #include "lib_disc/dof_manager/dof_distribution.h"
43 : #include "lib_disc/function_spaces/approximation_space.h"
44 :
45 : namespace ug{
46 :
47 : /// interface for adjustment of constraints
48 : /**
49 : * This class is the base class for the handling of constraints. Implementations
50 : * should adjust the jacobian/defect in order to guarantee the constraints. The
51 : * constraints can also be set in a solution vector by calling the adjust_solution
52 : * method.
53 : * Each implementation must specify the type of constraints since in some
54 : * situations it is important to apply the constraint modification in a certain
55 : * order.
56 : *
57 : * \tparam TDomain type of Domain
58 : * \tparam TDoFDistribution type of DoF Distribution
59 : * \tparam TAlgebra type of Algebra
60 : */
61 : template <typename TAlgebra>
62 : class IConstraint
63 : {
64 : public:
65 : /// Algebra type
66 : typedef TAlgebra algebra_type;
67 :
68 : /// Type of algebra matrix
69 : typedef typename algebra_type::matrix_type matrix_type;
70 :
71 : /// Type of algebra vector
72 : typedef typename algebra_type::vector_type vector_type;
73 :
74 : public:
75 : /// adapts jacobian to enforce constraints
76 : virtual void adjust_jacobian(matrix_type& J, const vector_type& u,
77 : ConstSmartPtr<DoFDistribution> dd, int type, number time = 0.0,
78 : ConstSmartPtr<VectorTimeSeries<vector_type> > vSol = SPNULL,
79 : const number s_a0 = 1.0) = 0;
80 :
81 : /// adapts defect to enforce constraints
82 : virtual void adjust_defect(vector_type& d, const vector_type& u,
83 : ConstSmartPtr<DoFDistribution> dd, int type, number time = 0.0,
84 : ConstSmartPtr<VectorTimeSeries<vector_type> > vSol = SPNULL,
85 : const std::vector<number>* vScaleMass = NULL,
86 : const std::vector<number>* vScaleStiff = NULL) = 0;
87 :
88 : /// adapts matrix and rhs (linear case) to enforce constraints
89 : virtual void adjust_linear(matrix_type& mat, vector_type& rhs,
90 : ConstSmartPtr<DoFDistribution> dd, int type, number time = 0.0) = 0;
91 :
92 : /// adapts a rhs to enforce constraints
93 : virtual void adjust_rhs(vector_type& rhs, const vector_type& u,
94 : ConstSmartPtr<DoFDistribution> dd, int type, number time = 0.0) = 0;
95 :
96 : /// sets the constraints in a solution vector
97 : virtual void adjust_solution(vector_type& u, ConstSmartPtr<DoFDistribution> dd, int type,
98 : number time = 0.0) = 0;
99 :
100 : /// adapts correction to enforce constraints
101 : /// for additive constraints (e.g. linear constraints, but NOT Dirichlet other than Dirichlet-0!),
102 : /// this is the same as adjust_solution, therefore default implementation
103 0 : virtual void adjust_correction(vector_type& c, ConstSmartPtr<DoFDistribution> dd, int type,
104 0 : number time = 0.0) {adjust_solution(c, dd, type, time);}
105 :
106 : /// adjust linear residual
107 0 : virtual void adjust_linear_residual(vector_type& d, const vector_type& u,
108 : ConstSmartPtr<DoFDistribution> dd, int type, number time = 0.0)
109 0 : {adjust_correction(d, dd, type, time);}
110 :
111 : /// adjusts the assembled error estimator values in the attachments according to the constraint
112 0 : virtual void adjust_error(const vector_type& u, ConstSmartPtr<DoFDistribution> dd, int type,
113 : number time = 0.0,
114 : ConstSmartPtr<VectorTimeSeries<vector_type> > vSol = SPNULL,
115 : const std::vector<number>* vScaleMass = NULL,
116 0 : const std::vector<number>* vScaleStiff = NULL) {};
117 :
118 : /// sets constraints in prolongation
119 0 : virtual void adjust_prolongation(matrix_type& P,
120 : ConstSmartPtr<DoFDistribution> ddFine,
121 : ConstSmartPtr<DoFDistribution> ddCoarse,
122 : int type,
123 0 : number time = 0.0) {};
124 :
125 : /// sets constraints in restriction
126 0 : virtual void adjust_restriction(matrix_type& R,
127 : ConstSmartPtr<DoFDistribution> ddCoarse,
128 : ConstSmartPtr<DoFDistribution> ddFine,
129 : int type,
130 0 : number time = 0.0) {};
131 :
132 : /// sets the constraints in a solution vector
133 0 : virtual void adjust_restriction(vector_type& uCoarse, GridLevel coarseLvl,
134 : const vector_type& uFine, GridLevel fineLvl,
135 0 : int type) {};
136 :
137 : /// sets the constraints in a solution vector
138 0 : virtual void adjust_prolongation(vector_type& uFine, GridLevel fineLvl,
139 : const vector_type& uCoarse, GridLevel coarseLvl,
140 0 : int type) {};
141 :
142 : /// modifies solution vector before calling the assembling routine
143 0 : virtual void modify_solution(vector_type& uMod, const vector_type& u,
144 0 : ConstSmartPtr<DoFDistribution> dd, int type) {};
145 :
146 : /// modify_solution for instationary case
147 0 : virtual void modify_solution(SmartPtr<VectorTimeSeries<vector_type> > vSolMod,
148 : ConstSmartPtr<VectorTimeSeries<vector_type> > vSol,
149 : ConstSmartPtr<DoFDistribution> dd,
150 0 : int type) {};
151 :
152 : /// returns the type of constraints
153 : virtual int type() const = 0;
154 :
155 : /// virtual destructor
156 : virtual ~IConstraint() {};
157 :
158 : };
159 :
160 : template <typename TDomain, typename TAlgebra>
161 : class IDomainConstraint : public IConstraint<TAlgebra>
162 : {
163 : public:
164 : /// Domain Type
165 : typedef TDomain domain_type;
166 :
167 : /// Algebra type
168 : typedef TAlgebra algebra_type;
169 :
170 : /// Type of algebra matrix
171 : typedef typename algebra_type::matrix_type matrix_type;
172 :
173 : /// Type of algebra vector
174 : typedef typename algebra_type::vector_type vector_type;
175 :
176 : public:
177 : /// constructor
178 0 : IDomainConstraint() : IConstraint<TAlgebra>(),
179 0 : m_bDoErrEst(false), m_spErrEstData(SPNULL) {};
180 :
181 : /// sets the approximation space
182 0 : virtual void set_approximation_space(SmartPtr<ApproximationSpace<TDomain> > approxSpace)
183 : {
184 0 : m_spApproxSpace = approxSpace;
185 0 : }
186 :
187 : /// returns approximation space
188 : SmartPtr<ApproximationSpace<TDomain> > approximation_space()
189 : {
190 : return m_spApproxSpace;
191 : }
192 :
193 : /// returns approximation space
194 : ConstSmartPtr<ApproximationSpace<TDomain> > approximation_space() const
195 : {
196 : return m_spApproxSpace;
197 : }
198 :
199 : /// returns the type of constraints
200 : virtual int type() const = 0;
201 :
202 : /// sets the assemble adapter for the constraints
203 : void set_ass_tuner(ConstSmartPtr<AssemblingTuner<TAlgebra> > spAssemblingTuner = NULL)
204 : {
205 0 : m_spAssTuner = spAssemblingTuner;
206 : }
207 :
208 : // //////////////////////////
209 : // Error estimator
210 : // //////////////////////////
211 : public:
212 : /// sets the pointer to an error estimator data object (or NULL)
213 : /**
214 : * This function sets the pointer to an error estimator data object
215 : * that should be used for this discretization. Note that the ElemDisc
216 : * object must use RTTI to try to convert this pointer to the type
217 : * of the objects accepted by it for this purpose. If the conversion
218 : * fails than an exception must be thrown since this situation is not
219 : * allowed.
220 : */
221 0 : void set_error_estimator(SmartPtr<IErrEstData<TDomain> > ee) {m_spErrEstData = ee; m_bDoErrEst = true;}
222 :
223 : /// find out whether or not a posteriori error estimation is to be performed for this disc
224 0 : bool err_est_enabled() const {return m_bDoErrEst;}
225 :
226 : /// returns the pointer to the error estimator data object (or NULL)
227 0 : virtual SmartPtr<IErrEstData<TDomain> > err_est_data() {return m_spErrEstData;}
228 :
229 : private:
230 : /// flag indicating whether or not a posteriori error estimation is to be performed for this disc
231 : bool m_bDoErrEst;
232 :
233 : protected:
234 : /// error estimation object associated to the element discretization
235 : SmartPtr<IErrEstData<TDomain> > m_spErrEstData;
236 :
237 : // ///////////////////////////
238 :
239 :
240 : protected:
241 : /// returns the level dof distribution
242 : ConstSmartPtr<DoFDistribution> dd(const GridLevel& gl) const
243 : {
244 : return m_spApproxSpace->dof_distribution(gl);
245 : }
246 :
247 : protected:
248 : /// Approximation Space
249 : SmartPtr<ApproximationSpace<TDomain> > m_spApproxSpace;
250 :
251 : /// Assemble adapter
252 : ConstSmartPtr<AssemblingTuner<TAlgebra> > m_spAssTuner;
253 :
254 : };
255 :
256 : } // end namespace ug
257 :
258 : #endif /* __H__UG__LIB_DISC__SPATIAL_DISC__CONSTRAINTS__CONSTRAINTS__INTERFACE__ */
|