Loading [MathJax]/extensions/tex2jax.js
Plugins
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Modules Pages
nedelec_transfer_impl.h
Go to the documentation of this file.
1/*
2 * Copyright (c) 2013-2014: G-CSC, Goethe University Frankfurt
3 * Author: Dmitry Logashenko
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
39#include "nedelec_local_ass.h"
40
41namespace ug{
42namespace Electromagnetism{
43
55template <typename TDomain, typename TAlgebra, typename TElem>
57(
58 const MultiGrid * mg,
59 Vertex * v,
60 TElem * base,
62)
63{
65
66// Get the parent of the vertex (note that this is typically not the 'base')
67 GridObject * parent = mg->get_parent (v);
68
69// Get the vertices of the parent
70 Grid::vertex_traits::secure_container vrts;
71 ((MultiGrid *) mg)->associated_elements (vrts, parent);
72 size_t n_co = vrts.size ();
73 UG_ASSERT (n_co != 0, "GetRegularLocalCoordinate: No associated vertices.")
74
75// Get the reference element for the base
77
78// Average the local coordinates of the vertices
79 local = 0.0;
80 for (size_t co = 0; co < n_co; co++)
81 {
82 int i = GetVertexIndex (base, vrts[co]);
83 if (i < 0)
84 UG_THROW ("GetRegularLocalCoordinate: No (implicit) parent-child relation between the vertex and the base element.");
85 local += rRefElem.corner (i);
86 }
87 local /= n_co;
88}
89
100template <typename TDomain, typename TAlgebra>
102(
103 const domain_type & domain,
104 const DoFDistribution & coarseDD,
105 const DoFDistribution & fineDD,
106 matrix_type & mat,
107 std::vector<bool> & vIsRestricted
108)
109{
111
112// Multiindices to access the components
113 std::vector<DoFIndex> c_ind (1), f_ind (1);
114
115// Get the multigrid:
116 const MultiGrid * grid = coarseDD.multi_grid().get ();
117
118// Get the number of the functions:
119 size_t num_fct = coarseDD.num_fct ();
120 UG_ASSERT (num_fct == fineDD.num_fct (), "NedelecTransfer: Coarse and find DD mismatch.");
121
122// Loop over all subsets:
123 for (int si = 0; si < coarseDD.num_subsets (); si++)
124 {
125 // Loop over all edges in the subset
126 iterator e_end = coarseDD.template end<RegularEdge> (si);
127 for (iterator edge_iter = coarseDD.template begin<RegularEdge> (si);
128 edge_iter != e_end; ++edge_iter)
129 {
130 number coef;
131 Edge * c_edge = * edge_iter;
132 const size_t n_children = grid->num_children<Edge, Edge> (c_edge);
133
134 if (n_children == 0)
135 continue;
136 else if (n_children == 1)
137 {
138 Edge * f_edge = grid->get_child<Edge, Edge> (c_edge, 0);
139
140 // Check the edge orientation:
141 GridObject * corner_0 = grid->get_parent (f_edge->vertex(0));
142 if (corner_0 == (GridObject *) (c_edge->vertex (0)))
143 coef = 1; // the edges should have the same orientation
144 else if (corner_0 == (GridObject *) (c_edge->vertex (1)))
145 coef = -1; // the edges have the inverse orientation
146 else
147 UG_THROW ("NedelecTransfer: Cannot find out the edge orientation.");
148
149 // Add the values:
150 for (size_t fct = 0; fct < num_fct; fct++)
151 {
152 // Get the DoFs:
153 if (coarseDD.inner_dof_indices (c_edge, fct, c_ind) != 1
154 || fineDD.inner_dof_indices (f_edge, fct, f_ind) != 1)
155 UG_THROW ("NedelecTransfer:"
156 "More than one DoF per edge. Not the Nedelec-type-1 element?");
157 // Assemble the + or - identity matrix here:
158 DoFRef (mat, f_ind[0], c_ind[0]) += coef;
159 vIsRestricted [c_ind[0][0]] = true;
160 }
161 }
162 else if (n_children == 2)
163 {
164 for (size_t child = 0; child < 2; child++)
165 {
166 Edge * f_edge = grid->get_child<Edge, Edge> (c_edge, child);
167
168 // Check the edge orientation:
169 GridObject * corner_0 = grid->get_parent (f_edge->vertex(0));
170 GridObject * corner_1 = grid->get_parent (f_edge->vertex(1));
171 if (corner_0 == (GridObject *) (c_edge->vertex (0))
172 || corner_1 == (GridObject *) (c_edge->vertex (1)))
173 coef = 0.5; // the edges should have the same orientation
174 else if (corner_0 == (GridObject *) (c_edge->vertex (1))
175 || corner_1 == (GridObject *) (c_edge->vertex (0)))
176 coef = -0.5; // the edges have the inverse orientation
177 else
178 UG_THROW ("NedelecTransfer: Cannot find out the edge orientation.");
179
180 // Add the values:
181 for (size_t fct = 0; fct < num_fct; fct++)
182 {
183 // Get the DoFs:
184 if (coarseDD.inner_dof_indices (c_edge, fct, c_ind) != 1
185 || fineDD.inner_dof_indices (f_edge, fct, f_ind) != 1)
186 UG_THROW ("NedelecTransfer:"
187 "More than one DoF per edge. Not the Nedelec-type-1 element?");
188 // Assemble the matrix entry:
189 DoFRef (mat, f_ind[0], c_ind[0]) += coef;
190 vIsRestricted [c_ind[0][0]] = true;
191 }
192 }
193 }
194 else
195 UG_THROW ("NedelecTransfer: Only regular refinement is supported,"
196 " but some edge is subdivided into " << n_children << " ones.");
197 }
198 }
199}
200
207template <typename TDomain, typename TAlgebra, typename TElem>
209(
210 const domain_type & domain,
211 const DoFDistribution & coarseDD,
212 const DoFDistribution & fineDD,
213 matrix_type & mat,
214 std::vector<bool> & vIsRestricted
215)
216{
217 typedef typename DoFDistribution::traits<TElem>::const_iterator iterator;
219
220// Multiindices to access the components
221 std::vector<DoFIndex> c_ind (ref_elem_type::numEdges), f_ind (1);
222
223// Get the multigrid:
224 const MultiGrid * grid = coarseDD.multi_grid().get ();
225
226// Get position accessor for the integration:
227 const typename TDomain::position_accessor_type & aaPos = domain.position_accessor ();
228
229// Get the number of the functions:
230 size_t num_fct = coarseDD.num_fct ();
231 UG_ASSERT (num_fct == fineDD.num_fct (), "NedelecTransfer: Coarse and find DD mismatch.");
232
233// Loop over all subsets:
234 for (int si = 0; si < coarseDD.num_subsets (); si++)
235 {
236 // Loop over all the elements of the given type in the subset
237 iterator e_end = coarseDD.template end<TElem> (si);
238 for (iterator elem_iter = coarseDD.template begin<TElem> (si);
239 elem_iter != e_end; ++elem_iter)
240 {
241 TElem * c_elem = * elem_iter;
242 const size_t n_children = grid->num_children<Edge, TElem> (c_elem);
243
244 if (n_children == 0)
245 continue;
246
247 // Get the corner positions:
248 position_type aCorners [ref_elem_type::numCorners];
249 for (size_t co = 0; co < (size_t) ref_elem_type::numCorners; co++)
250 aCorners [co] = aaPos [c_elem->vertex (co)];
251
252 // Loop over the child elements:
253 for (size_t child = 0; child < n_children; child++)
254 {
255 // Get the fine grid edge:
256 Edge * edge = grid->get_child<Edge, TElem> (c_elem, child);
257
258 // Get the local coordinates of the edge center w.r.t. the parent element:
259 MathVector<TElem::dim> loc_0, loc_1, loc_center;
260 GetRegularLocalCoordinate (grid, edge->vertex (0), c_elem, loc_0);
261 GetRegularLocalCoordinate (grid, edge->vertex (1), c_elem, loc_1);
262 VecAdd (loc_center, loc_0, loc_1);
263 loc_center /= 2;
264
265 // Get the length of the edge according to the local coordinates
266 MathVector<WDim> corner_0, corner_1, edge_vec;
267 DimReferenceMapping<dim, WDim> & map = ReferenceMappingProvider::get<dim, WDim> (roid);
268 map.update (aCorners);
269 map.local_to_global (corner_0, loc_0);
270 map.local_to_global (corner_1, loc_1);
271 VecSubtract (edge_vec, corner_1, corner_0);
272
273 // Get the shapes and assemble the contribution to the interpolation matrix:
274 MathVector<WDim> shape [ref_elem_type::numEdges];
275 NedelecT1_LDisc<TDomain, TElem>::get_shapes (&domain, c_elem, aCorners, loc_center, shape);
276 for (size_t fct = 0; fct < num_fct; fct++)
277 {
278 // Get the DoFs:
279 if (coarseDD.dof_indices (c_elem, fct, c_ind) != (size_t) ref_elem_type::numEdges
280 || fineDD.inner_dof_indices (edge, fct, f_ind) != 1)
281 UG_THROW ("NedelecTransfer:"
282 "More than one DoF per edge. Not the Nedelec-type-1 element?");
283 // Add the contributions to the matrix entry:
284 for (size_t i = 0; i < (size_t) ref_elem_type::numEdges; i++)
285 {
286 number coef = VecDot (shape [i], edge_vec);
287 DoFRef (mat, f_ind [0], c_ind [i]) += coef;
288 vIsRestricted [c_ind [i] [0]] = true;
289 }
290 }
291 }
292 }
293 }
294}
295
300template <typename TDomain, typename TAlgebra>
302{
303// The grid element types for this dimension:
304 typedef typename domain_traits<WDim>::AllElemList ElemList;
305
306// Verify the approximation space:
307 check_approximation_space ();
308
309// Get the DoF distributions:
310 const DoFDistribution & coarseDD = * m_spApproxSpace->dof_distribution (m_coarseLevel);
311 const DoFDistribution & fineDD = * m_spApproxSpace->dof_distribution (m_fineLevel);
312
313// Get number of dofs on the grid levels:
314 const size_t numFineDoFs = fineDD.num_indices ();
315 const size_t numCoarseDoFs = coarseDD.num_indices ();
316
317// Assemble the prolongation matrix:
318 m_vIsRestricted.clear (); m_vIsRestricted.resize (numCoarseDoFs, false);
319 try
320 {
321 m_prolongation_matrix.resize_and_clear (numFineDoFs, numCoarseDoFs);
322 }
323 UG_CATCH_THROW ("NedelecTransfer: Failed to allocate the prolongation matrix.");
324 boost::mpl::for_each<ElemList> (AssembleProlongationMatrix (this, coarseDD, fineDD));
325
326# ifdef UG_PARALLEL
327 m_prolongation_matrix.set_storage_type (PST_CONSISTENT);
328# endif
329
330// Done:
331 m_bInit = true;
332}
333
338template <typename TDomain, typename TAlgebra>
340(
341 vector_type & uFine,
342 const vector_type & uCoarse
343)
344{
345 if (! m_bInit)
346 UG_THROW("NedelecTransfer: Operator not initialized.");
347
348// Some assertions
349 UG_ASSERT(uFine.size () >= m_prolongation_matrix.num_rows (), "Vector ["
350 << uFine.size() << "] must be >= Row size " << m_prolongation_matrix.num_rows ());
351 UG_ASSERT(uCoarse.size () >= m_prolongation_matrix.num_cols (), "Vector ["
352 << uCoarse.size() << "] must be >= Col size " << m_prolongation_matrix.num_cols ());
353
354// Apply the prolongation matrix
355 if (! m_prolongation_matrix.apply (uFine, uCoarse))
356 {
357 std::stringstream ss;
358 ss << "NedelecTransfer::prolongate: Cannot apply matrix. ";
359#ifdef UG_PARALLEL
360 ss << "(Type uCoarse = " << uCoarse.get_storage_mask ();
361#endif
362 UG_THROW (ss.str().c_str ());
363 }
364
365// Set Dirichlet nodes to zero again
366 try
367 {
368 for (size_t i = 0; i < m_vConstraint.size (); ++i)
369 if (m_vConstraint[i]->type () & CT_DIRICHLET)
370 m_vConstraint[i]->adjust_defect (uFine, uFine,
371 m_spApproxSpace->dof_distribution (m_fineLevel), CT_DIRICHLET);
372 }
373 UG_CATCH_THROW("NedelecTransfer::prolongate: Error while setting dirichlet defect to zero.");
374
375// Call further constraints constraints (= adjust_restrict, member of class constraint)
376 try
377 {
378 for (int type = 1; type < CT_ALL; type = type << 1)
379 for (size_t i = 0; i < m_vConstraint.size (); ++i)
380 if (m_vConstraint[i]->type() & type)
381 m_vConstraint[i]->adjust_prolongation (uFine, m_fineLevel, uCoarse, m_coarseLevel, type);
382 }
383 UG_CATCH_THROW("NedelecTransfer::prolongate: Error while setting dirichlet defect to zero.");
384}
385
390template <typename TDomain, typename TAlgebra>
392(
393 vector_type & uCoarse,
394 const vector_type & uFine
395)
396{
397 if (! m_bInit)
398 UG_THROW("NedelecTransfer: Operator not initialized.");
399
400 vector_type uTmp; uTmp.resize(uCoarse.size());
401
402// Some assertions
403 UG_ASSERT (uFine.size () >= m_prolongation_matrix.num_rows (), "Vector ["
404 << uFine.size () << "] must be >= Row size " << m_prolongation_matrix.num_rows ());
405 UG_ASSERT (uCoarse.size () >= m_prolongation_matrix.num_cols (), "Vector ["
406 << uCoarse.size() << "] must be >= Col size " << m_prolongation_matrix.num_cols ());
407
408// Apply the transposed prolongation matrix
409 if (! m_prolongation_matrix.apply_transposed (uTmp, uFine))
410 UG_THROW ("NedelecTransfer::restrict: Cannot apply transposed matrix.");
411
412// Copy only restricted values:
413// This is needed for adaptive grids, where the defect that has been
414// projected from the surface should remain and only hidden (i.e.
415// indices with children) should be changed.
416 for (size_t i = 0; i < uTmp.size (); ++i)
417 if (m_vIsRestricted [i])
418 uCoarse[i] = uTmp[i];
419
420// Set dirichlet nodes to zero again
421 try
422 {
423 for(size_t i = 0; i < m_vConstraint.size (); ++i)
424 if (m_vConstraint[i]->type () & CT_DIRICHLET)
425 m_vConstraint[i]->adjust_defect (uCoarse, uCoarse,
426 m_spApproxSpace->dof_distribution (m_coarseLevel), CT_DIRICHLET);
427 }
428 UG_CATCH_THROW("NedelecTransfer::restrict: Error while setting dirichlet defect to zero.");
429
430// call restrictions due to added constraints (= adjust_restrict, member of class constraint)
431 try
432 {
433 for (int type = 1; type < CT_ALL; type = type << 1)
434 for (size_t i = 0; i < m_vConstraint.size (); ++i)
435 if (m_vConstraint[i]->type() & type)
436 m_vConstraint[i]->adjust_restriction (uCoarse, m_coarseLevel, uFine, m_fineLevel, type);
437 }
438 UG_CATCH_THROW("NedelecTransfer::restrict: Error while setting dirichlet defect to zero.");
439}
440
444template <typename TDomain, typename TAlgebra>
446{
447 if (m_spApproxSpace.invalid ())
448 UG_THROW ("NedelecTransfer: Approximation space not set.");
449
450 for (size_t fct = 0; fct < m_spApproxSpace->num_fct (); fct++)
451 if (m_spApproxSpace->local_finite_element_id (fct). type() != LFEID::NEDELEC)
452 UG_THROW ("NedelecTransfer: "
453 "All the functions should be based on the Nedelec element "
454 "but function #" << fct << " is not.");
455}
456
460template <typename TDomain, typename TAlgebra>
462(
463 GridLevel coarseLevel,
464 GridLevel fineLevel
465)
466{
467 m_bInit = false;
468 m_fineLevel = fineLevel;
469 m_coarseLevel = coarseLevel;
470
471 if(m_fineLevel.level () - m_coarseLevel.level () != 1)
472 UG_THROW("NedelecTransfer: Can only project between successive level.");
473
474 if(m_fineLevel.type () != GridLevel::LEVEL || m_coarseLevel.type () != GridLevel::LEVEL)
475 UG_THROW("NedelecTransfer: Can only project between level dof distributions, but fine="
476 << m_fineLevel << ", coarse=" << m_coarseLevel);
477}
478
482template <typename TDomain, typename TAlgebra>
484{
485 SmartPtr<NedelecTransfer> op (new NedelecTransfer (m_spApproxSpace));
486
487 for(size_t i = 0; i < m_vConstraint.size (); ++i)
488 op->add_constraint (m_vConstraint[i]);
489
490 return op;
491}
492
493
494} // end namespace Electromagnetism
495} // end namespace ug
496
497/* End of File */
virtual void local_to_global(MathVector< worldDim > &globPos, const MathVector< dim > &locPos) const=0
virtual void update(const MathVector< worldDim > *vCornerCoord)=0
size_t dof_indices(GridObject *elem, size_t fct, std::vector< DoFIndex > &ind, bool bHang=false, bool bClear=true) const
SmartPtr< MultiGrid > multi_grid()
size_t inner_dof_indices(GridObject *elem, size_t fct, std::vector< DoFIndex > &ind, bool bClear=true) const
size_t num_indices() const
virtual Vertex * vertex(size_t index) const
TAlgebra::matrix_type matrix_type
Type of Vector.
Definition nedelec_transfer.h:112
TDomain domain_type
Type of Domain.
Definition nedelec_transfer.h:106
TDomain domain_type
Type of Domain.
Definition nedelec_transfer.h:57
static void GetRegularLocalCoordinate(const MultiGrid *mg, Vertex *vrt, TElem *base, MathVector< TElem::dim > &local)
computes the local coordinates of a vertex according to the assumption of the regular refinement
Definition nedelec_transfer_impl.h:57
reference_element_traits< TElem >::reference_element_type ref_elem_type
reference element type
Definition nedelec_transfer.h:69
static void assemble_prolongation_matrix(const domain_type &domain, const DoFDistribution &coarseDD, const DoFDistribution &fineDD, matrix_type &mat, std::vector< bool > &vIsRestricted)
assembles the prolongation matrix for one type of the grid elements
Definition nedelec_transfer_impl.h:209
TDomain::position_type position_type
position type in the domain
Definition nedelec_transfer.h:75
TAlgebra::matrix_type matrix_type
Type of Vector.
Definition nedelec_transfer.h:63
static void get_shapes(const TDomain *domain, TElem *elem, const position_type *corners, const MathVector< dim > local, MathVector< WDim > shapes[])
computes the Nedelec shapes at a given point
Definition nedelec_local_ass.h:156
Class of the prolongation and the restriction of the Nedelec DoFs.
Definition nedelec_transfer.h:168
void prolongate(vector_type &uFineOut, const vector_type &uCoarse)
applies the prolongation
Definition nedelec_transfer_impl.h:340
TAlgebra::vector_type vector_type
Type of Vector.
Definition nedelec_transfer.h:183
void do_restrict(vector_type &uCoarse, const vector_type &uFine)
apples the restriction = transposed prolongation
Definition nedelec_transfer_impl.h:392
void init()
initializes the operator (computes the prolongation matrix etc)
Definition nedelec_transfer_impl.h:301
void check_approximation_space()
checks the approximation space
Definition nedelec_transfer_impl.h:445
virtual void set_levels(GridLevel coarseLevel, GridLevel fineLevel)
Set levels.
Definition nedelec_transfer_impl.h:462
SmartPtr< ITransferOperator< TDomain, TAlgebra > > clone()
returns new instance with same setting
Definition nedelec_transfer_impl.h:483
GridObject * get_parent(Edge *o) const
static TClass & get()
PST_CONSISTENT
SmartPtr< TGrid > grid()
int GetVertexIndex(EdgeVertices *e, Vertex *v)
Variant::Type type()
#define UG_ASSERT(expr, msg)
#define UG_CATCH_THROW(msg)
#define UG_THROW(msg)
double number
void VecAdd(vector_t &vOut, const vector_t &v, typename vector_t::value_type s)
void VecSubtract(vector_t &vOut, const vector_t &v, typename vector_t::value_type s)
vector_t::value_type VecDot(const vector_t &v1, const vector_t &v2)
CT_DIRICHLET
ReferenceObjectID
const number & DoFRef(const TMatrix &mat, const DoFIndex &iInd, const DoFIndex &jInd)
a helper class to call all the type-dependent assembling functions
Definition nedelec_transfer.h:230