Loading [MathJax]/extensions/tex2jax.js
Plugins
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Modules Pages
hiptmair_hybrid_smoother_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
33/*
34 * Implementation of the class functions templates of the class of the
35 * hybrid smoother by Hiptmair.
36 */
37
38// ug4 headers
40
41// plugin headers
42#ifdef UG_PARALLEL
43#include "../memreduce.h"
44#endif
45
46namespace ug{
47namespace Electromagnetism{
48
49//==== Initialization ====
50
59template <typename TDomain, typename TAlgebra>
61(
62 const DoFDistribution * pEdgeDD,
63 const DoFDistribution * pVertDD
64)
65{
67
68 m_vEdgeInfo.resize (pEdgeDD->num_indices (), false);
69 for (size_t i = 0; i < pEdgeDD->num_indices (); i++) m_vEdgeInfo[i].clear_flags ();
70
71 std::vector<size_t> vEdgeInd (1), vVertInd (1);
72
73// Get the edge-to-vertex correspondence
74 t_edge_iter edgeIterEnd = pEdgeDD->end<Edge> ();
75 for (t_edge_iter edgeIter = pEdgeDD->begin<Edge> (); edgeIter != edgeIterEnd; ++edgeIter)
76 {
77 Edge * pEdge = *edgeIter;
78 if (pEdgeDD->inner_algebra_indices (pEdge, vEdgeInd) != 1)
79 UG_THROW (name () << ": Edge DoF distribution mismatch. Not the Nedelec-Type-1 element?");
80 tEdgeInfo & EdgeInfo = m_vEdgeInfo [vEdgeInd [0]];
81
82 EdgeInfo.set_init ();
83
84 for (size_t i = 0; i < 2; i++)
85 {
86 Vertex * pVertex = pEdge->vertex (i);
87 if (pVertDD->inner_algebra_indices (pVertex, vVertInd) != 1)
88 UG_THROW (name () << ": Vertex DoF distribution mismatch. Not the Lagrange-Order-1 element?");
89 EdgeInfo.vrt_index [i] = vVertInd [0];
90 }
91 }
92
93 if (m_spDirichlet.invalid ()) return; // no Dirichlet boundaries specified
94
95// Mark Dirichlet edges
96 SubsetGroup dirichlet_ssgrp (pEdgeDD->subset_handler ());
97 m_spDirichlet->get_dirichlet_subsets (dirichlet_ssgrp);
98 for (size_t j = 0; j < dirichlet_ssgrp.size (); j++)
99 {
100 int si = dirichlet_ssgrp [j];
101 // Loop the edges in the subset
102 t_edge_iter iterEnd = pEdgeDD->end<Edge> (si);
103 for (t_edge_iter iter = pEdgeDD->begin<Edge> (si); iter != iterEnd; iter++)
104 {
105 Edge * pEdge = *iter;
106 if (pEdgeDD->inner_algebra_indices (pEdge, vEdgeInd) != 1)
107 UG_THROW (name () << ": Edge DoF distribution mismatch. Not the Nedelec-Type-1 element?");
108 m_vEdgeInfo[vEdgeInd[0]].set_Dirichlet ();
109 }
110 }
111}
112
127template <typename TDomain, typename TAlgebra>
129{
130 size_t N_edges = m_spSysMat->num_rows ();
131 size_t N_verts = m_spPotMat->num_rows ();
132
133 m_spPotMat->set (0.0);
134
135 // Flags of the 'conductive' vertices:
136 VariableArray1<bool> vConductiveVertex;
137 vConductiveVertex.resize (N_verts, false);
138 memset (& (vConductiveVertex.at (0)), 0, N_verts * sizeof (bool));
139
140 // 1. Compute G^T M^{(1)}_h G and mark the 'conductive' vertices.
141 // Loop over the rows:
142 for (size_t row = 0; row < N_edges; row++)
143 {
144 if (! m_vEdgeInfo[row].is_init ())
145 UG_THROW (name () << ": Uninitialized edge data in the computation of the potential.");
146
147 bool conductive_edge = false;
148 // Remark: For 'conductive edges', the mass matrix is non-zero. We check that.
149
150 size_t startVert [2];
151 memcpy (startVert, m_vEdgeInfo[row].vrt_index, 2 * sizeof (size_t));
152
153 // Loop over the connections:
154 for (typename matrix_type::row_iterator row_it = m_spSysMat->begin_row(row);
155 row_it != m_spSysMat->end_row(row); ++row_it)
156 {
157 double M_entry = BlockRef (row_it.value (), 1, 0); // the entry of M^{(1)}_h
158
159 // Skip zero entries
160 if (zero_mass_entry (M_entry))
161 continue; // to prevent extra connections in the potential stiffness matrix
162 else
163 conductive_edge = true;
164
165 size_t destVert [2];
166 memcpy (destVert, m_vEdgeInfo[row_it.index()].vrt_index, 2 * sizeof (size_t));
167
168 // Update the entries of the potential stiffness matrix:
169 (* m_spPotMat) (startVert[0], destVert[0]) += M_entry;
170 (* m_spPotMat) (startVert[0], destVert[1]) -= M_entry;
171 (* m_spPotMat) (startVert[1], destVert[0]) -= M_entry;
172 (* m_spPotMat) (startVert[1], destVert[1]) += M_entry;
173 }
174
175 // Mark the 'conductive' nodes:
176 if (conductive_edge)
177 vConductiveVertex[startVert[0]] = vConductiveVertex[startVert[1]] = true;
178 }
179# ifdef UG_PARALLEL
180 MemAllReduce<VariableArray1<bool>, IndexLayout, t_red_op_or>
181 (&vConductiveVertex, m_spPotMat->layouts()->master(),
182 m_spPotMat->layouts()->slave(), & (m_spPotMat->layouts()->comm()));
183# endif
184
185 // 2. Unmark vertices at the Dirichlet boundary:
186 // Loop over the edges:
187 for (size_t row = 0; row < N_edges; row++)
188 {
189 tEdgeInfo & EdgeInfo = m_vEdgeInfo [row];
190 if (EdgeInfo.is_Dirichlet ())
191 vConductiveVertex[EdgeInfo.vrt_index[0]]
192 = vConductiveVertex[EdgeInfo.vrt_index[1]] = false;
193 }
194# ifdef UG_PARALLEL
195 MemAllReduce<VariableArray1<bool>, IndexLayout, t_red_op_and>
196 (&vConductiveVertex, m_spPotMat->layouts()->master(),
197 m_spPotMat->layouts()->slave(), & (m_spPotMat->layouts()->comm()));
198# endif
199
200 // 3. Assemble the identity matrix for the 'non-conductive' vertices
201 // (otherwise we would have zero lines there).
202 // Loop over the vertices:
203 for (size_t row = 0; row < N_verts; row++)
204 if (! vConductiveVertex [row])
205 SetDirichletRow (* m_spPotMat.get (), row);
206
207 // 4. Move the conductive vertex marks into the EdgeInfo structores:
208 // Loop over the edges:
209 for (size_t row = 0; row < N_edges; row++)
210 {
211 tEdgeInfo & EdgeInfo = m_vEdgeInfo [row];
212 if (vConductiveVertex[EdgeInfo.vrt_index[0]])
213 EdgeInfo.set_conductive_vrt_0 ();
214 if (vConductiveVertex[EdgeInfo.vrt_index[1]])
215 EdgeInfo.set_conductive_vrt_1 ();
216 }
217}
218
241template <typename TDomain, typename TAlgebra>
243(
244 const DoFDistribution * pEdgeDD,
245 const DoFDistribution * pVertDD
246)
247{
248 size_t N_edges = pEdgeDD->num_indices ();
249 size_t N_verts = pVertDD->num_indices ();
250
251 // To be on the safe side:
252 if (N_edges != m_spSysMat->num_rows () || N_edges != m_spSysMat->num_cols ())
253 UG_THROW (name () << ": smoother not init. or illegal matrix type");
254
255 // Allocate the auxiliary matrices and vectors:
256 try
257 {
258 m_spPotMat->resize_and_clear (N_verts, N_verts);
259 }
260 UG_CATCH_THROW (name () << ": Failed to allocate the auxiliary matrix for the potential.");
261# ifdef UG_PARALLEL
262 m_spPotMat->set_layouts (pVertDD->layouts ());
263# endif
264
265 // Get the edge-vertex dof correspondence:
266 get_edge_vert_correspondence (pEdgeDD, pVertDD);
267
268 // Compute the potential stiffness matrix:
269 compute_GtMG ();
270}
271
281template <typename TDomain, typename TAlgebra>
283(
285 const vector_type & u
286)
287{
288 try
289 {
290 // Check the vertex approx. space:
291 if (m_spVertApproxSpace.invalid())
292 UG_THROW (name() << ", init: No approx. space for the potential specified.");
293
294 // Cast and check if a suitable operator
295 m_spSysMat = J.template cast_dynamic<MatrixOperator<matrix_type, vector_type> >();
296 if(m_spSysMat.invalid())
297 UG_THROW (name() << ", init: No or illegal lin. operator passed (not a matrix?).");
298
299 // Try to cast the vector to the grid function and get the edge DoF distribution:
300 const TGridFunc * pGridF;
301 if ((pGridF = dynamic_cast<const TGridFunc *> (&u)) == 0)
302 UG_THROW (name() << ", init: No DoFDistribution specified.");
303 m_spEdgeDD = ((TGridFunc *) pGridF)->dof_distribution ();
304 m_GridLevel = m_spEdgeDD->grid_level ();
305 m_spVertDD = m_spVertApproxSpace->dof_distribution (m_GridLevel);
306
307 // Prepare the smoother for the potential:
308 delete m_pPotCorRe; m_pPotCorRe = 0;
309 delete m_pPotCorIm; m_pPotCorIm = 0;
310
311 if (! m_bSkipVertex)
312 {
313 // Compute the matrix of the potential equation:
314 compute_potential_matrix (m_spEdgeDD.get (), m_spVertDD.get ());
315
316 // Allocate the grid functions for the defect in the potential space:
317 m_pPotCorRe = new TPotGridFunc (m_spVertApproxSpace, m_spVertDD);
318 m_pPotCorIm = new TPotGridFunc (m_spVertApproxSpace, m_spVertDD);
319# ifdef UG_PARALLEL
320 m_pPotCorRe->set_storage_type (PST_CONSISTENT);
321 m_pPotCorIm->set_storage_type (PST_CONSISTENT);
322# endif
323
324 // Initialize the subordinated smoother for the vertex dof:
325 m_pPotCorRe->set (0.0);
326 if (! m_spVertSmoother->init (m_spPotMat, *m_pPotCorRe))
327 UG_THROW (name() << ", init: Failed to initialize the subordinated vertex-based smoother.");
328 }
329
330 // Initialize the subordinated smoother for the edge dofs:
331 if (! m_spEdgeSmoother->init (J, u))
332 UG_THROW (name() << ", init: Failed to initialize the subordinated edge-based smoother.");
333 }
334 catch (...)
335 {
336 m_bInit = false;
337 delete m_pPotCorRe; m_pPotCorRe = 0;
338 delete m_pPotCorIm; m_pPotCorIm = 0;
339 m_vEdgeInfo.resize (0, false);
340 m_spPotMat->resize_and_clear (0, 0);
341 throw;
342 }
343
344 m_bInit = true;
345 return true;
346}
347
348//==== Computation of the correction ====
349
355template <typename TDomain, typename TAlgebra>
357(
358 const vector_type & d,
359 pot_vector_type & potDefRe,
360 pot_vector_type & potDefIm
361)
362{
363 size_t N_edges = m_vEdgeInfo.size ();
364
365 potDefRe.set (0.0); potDefIm.set (0.0);
366
367 // Loop over the edges:
368 for (size_t edge = 0; edge < N_edges; edge++)
369 {
370 tEdgeInfo & EdgeInfo = m_vEdgeInfo [edge];
371 number Re_d = BlockRef (d [edge], 0);
372 number Im_d = BlockRef (d [edge], 1);
373
374 if (EdgeInfo.conductive_vrt_0 ())
375 {
376 size_t i = EdgeInfo.vrt_index [0];
377 potDefRe [i] -= Re_d;
378 potDefIm [i] -= Im_d;
379 }
380
381 if (EdgeInfo.conductive_vrt_1 ())
382 {
383 size_t i = EdgeInfo.vrt_index [1];
384 potDefRe [i] += Re_d;
385 potDefIm [i] += Im_d;
386 }
387 }
388}
389
393template <typename TDomain, typename TAlgebra>
395(
396 pot_vector_type & potCorRe,
397 pot_vector_type & potCorIm,
398 vector_type & c
399)
400{
401 size_t N_edges = m_vEdgeInfo.size ();
402
403 // Loop over the edges:
404 for (size_t edge = 0; edge < N_edges; edge++)
405 {
406 tEdgeInfo & EdgeInfo = m_vEdgeInfo [edge];
407
408 BlockRef (c [edge], 0) += potCorRe [EdgeInfo.vrt_index[1]]
409 - potCorRe [EdgeInfo.vrt_index[0]];
410
411 BlockRef (c [edge], 1) += potCorIm [EdgeInfo.vrt_index[1]]
412 - potCorIm [EdgeInfo.vrt_index[0]];
413 }
414}
415
453template <typename TDomain, typename TAlgebra>
455(
456 vector_type & c,
457 const vector_type & d
458)
459{
460//---- Protection:
461
462// Check that object is initialized:
463 if(! m_bInit)
464 {
465 UG_LOG("ERROR in '"<<name()<<"::apply': Object not initialized.\n");
466 return false;
467 }
468
469// Check sizes:
470 if (d.size() != m_spSysMat->num_rows())
471 UG_THROW("ERROR in '"<<name()<<"::apply': Vector [size= "<<d.size()<<
472 "] and Row [size= "<<m_spSysMat->num_rows()<<"] sizes have to match!");
473 if (c.size() != m_spSysMat->num_cols())
474 UG_THROW("ERROR in '"<<name()<<"::apply': Vector [size= "<<c.size()<<
475 "] and Column [size= "<<m_spSysMat->num_cols()<<"] sizes have to match!");
476 if (d.size() != c.size())
477 UG_THROW("ERROR in '"<<name()<<"::apply': Vector [d size= "<<d.size()<<
478 ", c size = " << c.size() << "] sizes have to match!");
479
480//---- Numerics:
481
482// A temporary vectors for the defect:
483 vector_type auxEdgeDef (d.size ());
484
485 auxEdgeDef = d;
486
487// Apply the edge-based smoother:
488 if (! m_bSkipEdge)
489 {
490 if (! m_spEdgeSmoother->apply_update_defect (c, auxEdgeDef))
491 return false;
492# ifdef UG_PARALLEL
493 c.change_storage_type (PST_CONSISTENT);
494# endif
495 }
496 else c.set (0.0);
497
498 if (! m_bSkipVertex)
499 {
500 // Temporary vectors:
501 pot_vector_type potDefRe (m_pPotCorRe->size ());
502 pot_vector_type potDefIm (m_pPotCorIm->size ());
503# ifdef UG_PARALLEL
504 potDefRe.set_layouts (m_pPotCorRe->layouts ());
506 potDefIm.set_layouts (m_pPotCorIm->layouts ());
508# endif
509
510 // Compute the 'vertex defect' of the potential:
511 collect_edge_defect (auxEdgeDef, potDefRe, potDefIm);
512
513 // Apply the vertex-centered smoother:
514 if (! m_spVertSmoother->apply (*m_pPotCorRe, potDefIm))
515 return false;
516 if (! m_spVertSmoother->apply (*m_pPotCorIm, potDefRe))
517 return false;
518# ifdef UG_PARALLEL
519 m_pPotCorRe->change_storage_type (PST_CONSISTENT);
520 m_pPotCorIm->change_storage_type (PST_CONSISTENT);
521# endif
522
523 // Add the vertex-centered correction into the edge-centered one:
524 distribute_vertex_correction (*m_pPotCorRe, *m_pPotCorIm, c);
525 }
526
527// Damping (the standard way, like for IPreconditioner):
528 number kappa = this->damping()->damping (c, d, m_spSysMat);
529 if (kappa != 1.0)
530 c *= kappa;
531
532 return true;
533}
534
535} // end namespace Electromagnetism
536} // end namespace ug
537
538/* End of File */
traits< TElem >::iterator begin()
ConstSmartPtr< AlgebraLayouts > layouts() const
traits< TElem >::iterator end()
size_t inner_algebra_indices(GridObject *elem, std::vector< size_t > &ind, bool bClear=true) const
size_t num_indices() const
ConstSmartPtr< ISubsetHandler > subset_handler() const
virtual Vertex * vertex(size_t index) const
void collect_edge_defect(const vector_type &d, pot_vector_type &potDefRe, pot_vector_type &potDefIm)
Computes the vertex-centered defect .
Definition hiptmair_hybrid_smoother_impl.h:357
bool init(SmartPtr< ILinearOperator< vector_type > > J, const vector_type &u)
Initialization using a matrix and a GridFunction (not merely a vector!)
Definition hiptmair_hybrid_smoother_impl.h:283
bool apply(vector_type &c, const vector_type &d)
Computes the correction.
Definition hiptmair_hybrid_smoother_impl.h:455
void compute_GtMG()
Computes the product .
Definition hiptmair_hybrid_smoother_impl.h:128
void compute_potential_matrix(const DoFDistribution *pEdgeDD, const DoFDistribution *pVertDD)
Computes the matrix for the smoother in the potential space and marks the "conductive nodes".
Definition hiptmair_hybrid_smoother_impl.h:243
void distribute_vertex_correction(pot_vector_type &potCorRe, pot_vector_type &potCorIm, vector_type &c)
Adds the vertex-centered correction to the edge-centered one:
Definition hiptmair_hybrid_smoother_impl.h:395
TAlgebra::vector_type vector_type
Vector type.
Definition hiptmair_hybrid_smoother.h:87
void get_edge_vert_correspondence(const DoFDistribution *pEdgeDD, const DoFDistribution *pVertDD)
Gets the correspondence between the edges and the vertices.
Definition hiptmair_hybrid_smoother_impl.h:61
void set(number w)
void set_storage_type(uint type)
void set_layouts(ConstSmartPtr< AlgebraLayouts > layouts)
size_t size() const
bool resize(size_type n, bool bCopyValues=true)
T & at(size_type i)
pcl::SingleLevelLayout< pcl::OrderedInterface< size_t, std::vector > > IndexLayout
PST_CONSISTENT
PST_ADDITIVE
#define UG_CATCH_THROW(msg)
#define UG_THROW(msg)
#define UG_LOG(msg)
double number
void SetDirichletRow(TMatrix &mat, const DoFIndex &ind)
const number & BlockRef(const number &m, size_t i)
"and" reduction operation class for the conductivity condition
Definition hiptmair_hybrid_smoother.h:302
"or" reduction operation class for the conductivity condition
Definition hiptmair_hybrid_smoother.h:297
Structure of the storage for the information about the edge-vertex interconnections.
Definition hiptmair_hybrid_smoother.h:145
void set_init()
Definition hiptmair_hybrid_smoother.h:159
bool is_Dirichlet()
Definition hiptmair_hybrid_smoother.h:150
void set_conductive_vrt_1()
Definition hiptmair_hybrid_smoother.h:158
bool conductive_vrt_0()
Definition hiptmair_hybrid_smoother.h:151
void set_conductive_vrt_0()
Definition hiptmair_hybrid_smoother.h:157
size_t vrt_index[2]
vertex dof's of the beginning and the end of the edge
Definition hiptmair_hybrid_smoother.h:146
bool conductive_vrt_1()
Definition hiptmair_hybrid_smoother.h:152