Loading [MathJax]/extensions/tex2jax.js
ug4
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Modules Pages
std_injection_impl.h
Go to the documentation of this file.
1/*
2 * Copyright (c) 2012-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__OPERATOR__LINEAR_OPERATOR__STD_INJECTION_IMPL__
34#define __H__UG__LIB_DISC__OPERATOR__LINEAR_OPERATOR__STD_INJECTION_IMPL__
35
36#include "std_injection.h"
37
38namespace ug{
39
49template <typename TAlgebra>
50void AssembleInjectionForP1Lagrange(typename TAlgebra::matrix_type& mat,
51 const DoFDistribution& coarseDD,
52 const DoFDistribution& fineDD)
53{
54 PROFILE_FUNC_GROUP("gmg");
55// Allow only lagrange P1 functions
56 for(size_t fct = 0; fct < fineDD.num_fct(); ++fct)
57 if(fineDD.local_finite_element_id(fct).type() != LFEID::LAGRANGE ||
58 fineDD.local_finite_element_id(fct).order() != 1)
59 UG_THROW("AssembleInjectionForP1Lagrange: "
60 "Interpolation only implemented for Lagrange P1 functions.");
61
62// get MultiGrid
63 const MultiGrid& grid = *coarseDD.multi_grid();
64
65// get number of dofs on different levels
66 const size_t numFineDoFs = fineDD.num_indices();
67 const size_t numCoarseDoFs = coarseDD.num_indices();
68
69// resize matrix
70 mat.resize_and_clear(numCoarseDoFs, numFineDoFs);
71
72 std::vector<size_t> coarseInd, fineInd;
73
74// RegularVertex iterators
76 const_iterator iter, iterBegin, iterEnd;
77
78 iterBegin = fineDD.template begin<RegularVertex>();
79 iterEnd = fineDD.template end<RegularVertex>();
80
81// loop nodes of fine subset
82 for(iter = iterBegin; iter != iterEnd; ++iter)
83 {
84 // get father
85 GridObject* geomObj = grid.get_parent(*iter);
86 Vertex* vert = dynamic_cast<Vertex*>(geomObj);
87
88 // Check if father is RegularVertex
89 if(vert != NULL)
90 {
91 // get global indices
92 coarseDD.inner_algebra_indices(vert, coarseInd);
93 }
94 else continue;
95
96 // get global indices
97 fineDD.inner_algebra_indices(*iter, fineInd);
98
99 for(size_t i = 0; i < coarseInd.size(); ++i)
100 mat(coarseInd[i], fineInd[i]) = 1.0;
101 }
102}
103
104template <int dim, typename TAlgebra>
105void AssembleInjectionByAverageOfChildren(typename TAlgebra::matrix_type& mat,
106 const DoFDistribution& coarseDD,
107 const DoFDistribution& fineDD)
108{
109 PROFILE_FUNC_GROUP("gmg");
110// get MultiGrid
111 const MultiGrid& grid = *coarseDD.multi_grid();
112
113 std::vector<size_t> coarseInd, fineInd;
114
115// RegularVertex iterators
116 typedef typename DoFDistribution::dim_traits<dim>::const_iterator const_iterator;
118 const_iterator iter, iterBegin, iterEnd;
119
120 iterBegin = coarseDD.template begin<Element>();
121 iterEnd = coarseDD.template end<Element>();
122
123// loop elements of coarse subset
124 for(iter = iterBegin; iter != iterEnd; ++iter)
125 {
126 // get element
127 Element* coarseElem = *iter;
128
129 // get children
130 const size_t numChild = grid.num_children<Element, Element>(coarseElem);
131 if(numChild == 0) continue;
132
133 // get global indices
134 coarseDD.inner_algebra_indices(coarseElem, coarseInd);
135
136 for(size_t c = 0; c < numChild; ++c)
137 {
138 Element* child = grid.get_child<Element, Element>(coarseElem, c);
139
140 fineDD.inner_algebra_indices(child, fineInd);
141
142 for(size_t i = 0; i < coarseInd.size(); ++i)
143 mat(coarseInd[i], fineInd[i]) = 1.0/(numChild);
144 }
145 }
146}
147
148template <typename TAlgebra>
149void AssembleInjectionByAverageOfChildren(typename TAlgebra::matrix_type& mat,
150 const DoFDistribution& coarseDD,
151 const DoFDistribution& fineDD)
152{
153// get number of dofs on different levels
154 const size_t numFineDoFs = fineDD.num_indices();
155 const size_t numCoarseDoFs = coarseDD.num_indices();
156
157// resize matrix
158 mat.resize_and_clear(numCoarseDoFs, numFineDoFs);
159
160 if(coarseDD.max_dofs(VERTEX)) AssembleInjectionByAverageOfChildren<0, TAlgebra>(mat, coarseDD, fineDD);
161 if(coarseDD.max_dofs(EDGE)) AssembleInjectionByAverageOfChildren<1, TAlgebra>(mat, coarseDD, fineDD);
162 if(coarseDD.max_dofs(FACE)) AssembleInjectionByAverageOfChildren<2, TAlgebra>(mat, coarseDD, fineDD);
163 if(coarseDD.max_dofs(VOLUME)) AssembleInjectionByAverageOfChildren<3, TAlgebra>(mat, coarseDD, fineDD);
164}
165
166template <typename TDomain, typename TAlgebra>
167template <typename TElem>
170 const DoFDistribution& coarseDD, const DoFDistribution& fineDD)
171{
172 PROFILE_FUNC_GROUP("gmg");
173
174 std::vector<size_t> vCoarseIndex, vFineIndex;
175 const MultiGrid& mg = *coarseDD.multi_grid();
176
177// iterators
178 typedef typename DoFDistribution::traits<TElem>::const_iterator const_iterator;
179 const_iterator iter, iterBegin, iterEnd;
180
181// loop subsets on fine level
182 for(int si = 0; si < coarseDD.num_subsets(); ++si)
183 {
184 iterBegin = coarseDD.template begin<TElem>(si);
185 iterEnd = coarseDD.template end<TElem>(si);
186
187 // loop vertices for fine level subset
188 for(iter = iterBegin; iter != iterEnd; ++iter)
189 {
190 // get element
191 TElem* coarseElem = *iter;
192
193 const size_t numChild = mg.num_children<TElem, TElem>(coarseElem);
194 if(numChild != 0) continue;
195
196 // get indices
197 coarseDD.inner_algebra_indices(coarseElem, vCoarseIndex);
198 fineDD.inner_algebra_indices(coarseElem, vFineIndex);
199 UG_ASSERT(vCoarseIndex.size() == vFineIndex.size(), "Size mismatch");
200
201 // set identity
202 for(size_t i = 0; i < vCoarseIndex.size(); ++i)
203 mat(vCoarseIndex[i], vFineIndex[i]) = 1.0;
204 }
205 }
206}
207
208template <typename TDomain, typename TAlgebra>
211 const DoFDistribution& coarseDD, const DoFDistribution& fineDD)
212{
213 if(coarseDD.max_dofs(VERTEX)) set_identity_on_pure_surface<Vertex>(mat, coarseDD, fineDD);
214 if(coarseDD.max_dofs(EDGE)) set_identity_on_pure_surface<Edge>(mat, coarseDD, fineDD);
215 if(coarseDD.max_dofs(FACE)) set_identity_on_pure_surface<Face>(mat, coarseDD, fineDD);
216 if(coarseDD.max_dofs(VOLUME)) set_identity_on_pure_surface<Volume>(mat, coarseDD, fineDD);
217}
218
221// StdInjection
224
225
226template <typename TDomain, typename TAlgebra>
229{
230 m_spApproxSpace = approxSpace;
231}
232
233template <typename TDomain, typename TAlgebra>
235set_levels(GridLevel coarseLevel, GridLevel fineLevel)
236{
237 m_fineLevel = fineLevel;
238 m_coarseLevel = coarseLevel;
239
240 if(m_fineLevel.level() - m_coarseLevel.level() != 1)
241 UG_THROW("StdInjection::set_levels:"
242 " Can only project between successive level.");
243}
244
245template <typename TDomain, typename TAlgebra>
247{
248 PROFILE_FUNC_GROUP("gmg");
249 if(!m_spApproxSpace.valid())
250 UG_THROW("StdInjection::init: "
251 "Approximation Space not set. Cannot init Projection.");
252
253// check only lagrange P1 functions
254 bool P1LagrangeOnly = true;
255 for(size_t fct = 0; fct < m_spApproxSpace->num_fct(); ++fct)
256 if(m_spApproxSpace->local_finite_element_id(fct).type() != LFEID::LAGRANGE ||
257 m_spApproxSpace->local_finite_element_id(fct).order() != 1)
258 P1LagrangeOnly = false;
259
260 try{
261 if(P1LagrangeOnly)
262 {
263 AssembleInjectionForP1Lagrange<TAlgebra>
264 (m_matrix,
265 *m_spApproxSpace->dof_distribution(m_coarseLevel),
266 *m_spApproxSpace->dof_distribution(m_fineLevel));
267 }
268 else
269 {
270 AssembleInjectionByAverageOfChildren<TAlgebra>
271 (m_matrix,
272 *m_spApproxSpace->dof_distribution(m_coarseLevel),
273 *m_spApproxSpace->dof_distribution(m_fineLevel));
274 }
275
276 } UG_CATCH_THROW("StdInjection::init():"
277 " Cannot assemble interpolation matrix.");
278
279 if(m_coarseLevel.is_surface()){
280 set_identity_on_pure_surface(m_matrix, *m_spApproxSpace->dof_distribution(m_coarseLevel), *m_spApproxSpace->dof_distribution(m_fineLevel));
281 }
282
283 #ifdef UG_PARALLEL
284 m_matrix.set_storage_type(PST_CONSISTENT);
285 #endif
286
287 m_bInit = true;
288}
289
290template <typename TDomain, typename TAlgebra>
292prolongate(vector_type& uFine, const vector_type& uCoarse)
293{
294 PROFILE_FUNC_GROUP("gmg");
295// Check, that operator is initiallized
296 if(!m_bInit)
297 UG_THROW("StdInjection::apply:"
298 " Operator not initialized.");
299
300// Some Assertions
301 UG_ASSERT(uFine.size() >= m_matrix.num_rows(),
302 "Vector [size= " << uFine.size() << "] and Rows [size= "
303 << m_matrix.num_rows() <<"] sizes have to match!");
304 UG_ASSERT(uCoarse.size() >= m_matrix.num_cols(), "Vector [size= "
305 << uCoarse.size() << "] and Cols [size= " <<
306 m_matrix.num_cols() <<"] sizes have to match!");
307
308// Apply matrix
309 if(!m_matrix.apply_transposed(uFine, uCoarse))
310 UG_THROW("StdInjection::apply: Cannot apply matrix.");
311}
312
313template <typename TDomain, typename TAlgebra>
315do_restrict(vector_type& uCoarse, const vector_type& uFine)
316{
317 PROFILE_FUNC_GROUP("gmg");
318// Check, that operator is initialized
319 if(!m_bInit)
320 UG_THROW("StdInjection::apply_transposed:"
321 "Operator not initialized.");
322
323// Some Assertions
324 UG_ASSERT(uFine.size() >= m_matrix.num_cols(),
325 "Vector [size= " << uFine.size() << "] and Cols [size= "
326 << m_matrix.num_cols() <<"] sizes have to match!");
327 UG_ASSERT(uCoarse.size() >= m_matrix.num_rows(), "Vector [size= "
328 << uCoarse.size() << "] and Rows [size= " <<
329 m_matrix.num_rows() <<"] sizes have to match!");
330
331// Apply matrix
332 try{
333 m_matrix.apply_ignore_zero_rows(uCoarse, 1.0, uFine);
334 }
335 UG_CATCH_THROW("StdInjection::apply_transposed:"
336 " Cannot apply transposed matrix.");
337}
338
339
340template <typename TDomain, typename TAlgebra>
343{
345 op->set_approximation_space(m_spApproxSpace);
346 return op;
347}
348
349} // end namespace ug
350
351#endif /* __H__UG__LIB_DISC__OPERATOR__LINEAR_OPERATOR__STD_INJECTION_IMPL__ */
Definition smart_pointer.h:108
base class for approximation spaces without type of algebra or dof distribution
Definition approximation_space.h:279
Definition dof_distribution.h:51
SmartPtr< MultiGrid > multi_grid()
Definition dof_distribution.h:69
size_t inner_algebra_indices(GridObject *elem, std::vector< size_t > &ind, bool bClear=true) const
extracts all algebra indices in the inner of the element (not sorted)
Definition dof_distribution.cpp:1276
size_t num_indices() const
return the number of dofs distributed
Definition dof_distribution.h:218
const LFEID & local_finite_element_id(size_t fct) const
Definition dof_distribution_info.h:225
int num_subsets() const
returns number of subsets
Definition dof_distribution_info.h:192
size_t max_dofs(const int dim) const
returns the maximum number of dofs on grid objects in a dimension
Definition dof_distribution_info.h:231
size_t num_fct() const
number of discrete functions on subset si
Definition dof_distribution_info.h:176
Definition grid_level.h:42
The base class for all geometric objects, such as vertices, edges, faces, volumes,...
Definition grid_base_objects.h:157
SpaceType type() const
returns the type of the local finite element
Definition local_finite_element_id.h:132
int order() const
returns the order of the local finite element
Definition local_finite_element_id.h:126
@ LAGRANGE
Definition local_finite_element_id.h:104
Definition multi_grid.h:72
TChild * get_child(TElem *elem, size_t ind) const
returns the i-th child of the given child-type
Definition multi_grid.h:268
GridObject * get_parent(GridObject *parent) const
Definition multi_grid.cpp:180
size_t num_children(TElem *elem) const
returns the number of children of the given child-type
Definition multi_grid.h:225
Definition std_injection.h:61
virtual SmartPtr< ITransferOperator< TDomain, TAlgebra > > clone()
clones the operator
Definition std_injection_impl.h:342
void set_identity_on_pure_surface(matrix_type &mat, const DoFDistribution &coarseDD, const DoFDistribution &fineDD)
Definition std_injection_impl.h:169
virtual void prolongate(vector_type &uFine, const vector_type &uCoarse)
Project uFine to uCoarse; uCoarse = P(uFine);.
Definition std_injection_impl.h:292
virtual void init()
Init operator.
Definition std_injection_impl.h:246
virtual void do_restrict(vector_type &uCoarse, const vector_type &uFine)
Apply Transposed Operator u = L^T*f.
Definition std_injection_impl.h:315
TAlgebra::vector_type vector_type
Type of Vector.
Definition std_injection.h:70
void set_approximation_space(SmartPtr< ApproximationSpace< TDomain > > approxSpace)
Set Approximation Space.
Definition std_injection_impl.h:228
TAlgebra::matrix_type matrix_type
Type of Vector.
Definition std_injection.h:73
void set_levels(GridLevel coarseLevel, GridLevel fineLevel)
Set approximation level.
Definition std_injection_impl.h:235
Base-class for all vertex-types.
Definition grid_base_objects.h:231
@ PST_CONSISTENT
Definition parallel_storage_type.h:68
#define UG_ASSERT(expr, msg)
Definition assert.h:70
#define UG_CATCH_THROW(msg)
Definition error.h:64
#define UG_THROW(msg)
Definition error.h:57
the ug namespace
void AssembleInjectionForP1Lagrange(typename TAlgebra::matrix_type &mat, const DoFDistribution &coarseDD, const DoFDistribution &fineDD)
Definition std_injection_impl.h:50
@ VOLUME
Definition grid_base_objects.h:63
@ VERTEX
Definition grid_base_objects.h:60
@ EDGE
Definition grid_base_objects.h:61
@ FACE
Definition grid_base_objects.h:62
void AssembleInjectionByAverageOfChildren(typename TAlgebra::matrix_type &mat, const DoFDistribution &coarseDD, const DoFDistribution &fineDD)
Definition std_injection_impl.h:105
#define PROFILE_FUNC_GROUP(groups)
Definition profiler.h:258
domain_traits< dim >::grid_base_object grid_base_object
Definition dof_distribution.h:88
SurfaceView::traits< grid_base_object >::const_iterator const_iterator
Definition dof_distribution.h:90
SurfaceView::traits< TElem >::const_iterator const_iterator
Definition dof_distribution.h:82