Plugins
hiptmair_hybrid_smoother.h
Go to the documentation of this file.
1 /*
2  * Copyright (c) 2013: 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  * The hybrid smoother proposed by R. Hiptmair for the discretizations
35  * based on the Nedelec element. The smoother is implemented for the
36  * discretizations of the time-harmonic equations.
37  */
38 #ifndef __H__UG__PLUGINS__ELECTROMAGNETISM__HIPTMAIR_HYBRID_SMOOTHER__
39 #define __H__UG__PLUGINS__ELECTROMAGNETISM__HIPTMAIR_HYBRID_SMOOTHER__
40 
41 // basic ug4 headers
42 #include "common/common.h"
43 #include "lib_grid/lg_base.h"
44 
45 // library-specific ug4 headers
50 
51 // plugin headers
52 #include "../em_material.h"
53 
54 namespace ug{
55 namespace Electromagnetism{
56 
58 
78 template <typename TDomain, typename TAlgebra>
80  public ILinearIterator<typename TAlgebra::vector_type>,
81  public DebugWritingObject<TAlgebra>
82 {
84 
85 public:
92 
101 
106 
107 private:
108 
110  bool zero_mass_entry (double val) {return fabs (val) < 1e-64;} // this could be merely 0
111 
114 
117 
122 
125 
130 
133 
136 
139 
142 
144  struct tEdgeInfo
145  {
146  size_t vrt_index [2];
147 
148  char flags;
149 
150  bool is_Dirichlet () {return (flags & 1) != 0;}
151  bool conductive_vrt_0 () {return (flags & 2) != 0;}
152  bool conductive_vrt_1 () {return (flags & 4) != 0;}
153  bool is_init () {return (flags & 8) != 0;}
154 
155  void clear_flags () {flags = 0;}
156  void set_Dirichlet () {flags |= 1;}
157  void set_conductive_vrt_0 () {flags |= 2;}
158  void set_conductive_vrt_1 () {flags |= 4;}
159  void set_init () {flags |= 8;}
160  };
163 
165  bool m_bInit;
166 
169 
170 public:
173  (
174  SmartPtr<ApproximationSpace<TDomain> > vertApproxSpace,
175  SmartPtr<ILinearIterator<vector_type> > edgeSmoother,
177  )
179  m_pPotCorRe (NULL), m_pPotCorIm (NULL),
180  m_spVertApproxSpace (vertApproxSpace),
181  m_spEdgeSmoother (edgeSmoother), m_spVertSmoother (vertSmoother),
182  m_bInit (false),
183  m_bSkipEdge (false), m_bSkipVertex (false)
184  {
186  UG_THROW (name() << ": illegal vert.-centered approx. space");
187  if (m_spEdgeSmoother.invalid ())
188  UG_THROW (name() << ": illegal edge-centered smoother");
189  if (m_spVertSmoother.invalid ())
190  UG_THROW (name() << ": illegal vert.-centered smoother");
191  }
192 
195  {
196  delete m_pPotCorRe; delete m_pPotCorIm;
197  }
198 
200  const char* name() const {return "Hiptmair hybrid smoother for Whitney-1 elements";}
201 
203  bool supports_parallel() const {return true;}
204 
207 
210  {
211  UG_THROW(name () << ": Cannot initialize the hybrid smoother without the geometry. Specify the 2nd argument for init!");
212  return false;
213  }
214 
216  bool apply (vector_type & c, const vector_type & d);
217 
220  {
221  // compute the correction
222  if(! apply (c, d)) return false;
223 
224  // update the defect d := d - A*c
225  if(! m_spSysMat->matmul_minus (d, c))
226  UG_THROW (name() << "::apply_update_defect: failed to execute matmul_minus");
227 
228  return true;
229  }
230 
233  (
235  )
236  {
237  m_spDirichlet = spDirichlet;
238  }
239 
241  void set_skip_edge_smoother (bool skip_edge) {m_bSkipEdge = skip_edge;}
243  void set_skip_vertex_smoother (bool skip_vertex) {m_bSkipVertex = skip_vertex;}
244 
247  {
248  SmartPtr<this_type> newInst
250  m_spEdgeSmoother->clone (), m_spVertSmoother->clone ()));
251  newInst->set_Dirichlet (m_spDirichlet);
252  newInst->set_skip_edge_smoother (m_bSkipEdge);
253  newInst->set_skip_vertex_smoother (m_bSkipVertex);
254  newInst->set_debug (this->debug_writer ());
255  newInst->set_damp (this->damping ());
256  return newInst;
257  }
258 
259 private: // Auxiliary functions:
260 
263  (
264  const DoFDistribution * pEdgeDD,
265  const DoFDistribution * pVertDD
266  );
267 
270  (
271  const DoFDistribution * pEdgeDD,
272  const DoFDistribution * pVertDD
273  );
274 
276  void compute_GtMG ();
277 
280  (
281  const vector_type & d,
282  pot_vector_type & potDefRe,
283  pot_vector_type & potDefIm
284  );
285 
288  (
289  pot_vector_type & potCorRe,
290  pot_vector_type & potCorIm,
291  vector_type & c
292  );
293 
294 #ifdef UG_PARALLEL
296  struct t_red_op_or
297  {
298  static inline bool op (bool a, bool b) {return a || b;}
299  };
302  {
303  static inline bool op (bool a, bool b) {return a && b;}
304  };
305 #endif
306 
307 }; // class TimeHarmonicNedelecHybridSmoother
308 
309 } // end namespace Electromagnetism
310 } // end namespace ug
311 
312 // Implementation of the functions:
314 
315 #endif /* __H__UG__PLUGINS__ELECTROMAGNETISM__HIPTMAIR_HYBRID_SMOOTHER__ */
316 
317 /* End of File */
bool invalid() const
SmartPtr< IDebugWriter< algebra_type > > debug_writer()
Common interface to get the Dirichlet boundary conditions.
Definition: em_material.h:287
The hybrid smoother by R. Hiptmair.
Definition: hiptmair_hybrid_smoother.h:82
MatrixOperator< matrix_type, vector_type > matrix_operator_type
Matrix Operator type.
Definition: hiptmair_hybrid_smoother.h:91
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
~TimeHarmonicNedelecHybridSmoother()
Destructor.
Definition: hiptmair_hybrid_smoother.h:194
CPUAlgebra TPotAlgebra
The auxiliary algebra type for the space of the potential. (Note: It should be scalar....
Definition: hiptmair_hybrid_smoother.h:94
TPotGridFunc * m_pPotCorIm
Vertex-centered grid function for the Im-part of potential corrections (this is a GridFunction to all...
Definition: hiptmair_hybrid_smoother.h:121
GridFunction< TDomain, TAlgebra > TGridFunc
Grid function type for the solution.
Definition: hiptmair_hybrid_smoother.h:103
TPotAlgebra::matrix_type pot_matrix_type
Matrix type for the potential space.
Definition: hiptmair_hybrid_smoother.h:98
void set_Dirichlet(SmartPtr< EMDirichlet< TDomain, TAlgebra > > spDirichlet)
Sets the Dirichlet boundary.
Definition: hiptmair_hybrid_smoother.h:233
GridFunction< TDomain, TPotAlgebra > TPotGridFunc
Grid function type for the potential.
Definition: hiptmair_hybrid_smoother.h:105
bool apply_update_defect(vector_type &c, vector_type &d)
Computes the correction and updates the defect d := d - L*c.
Definition: hiptmair_hybrid_smoother.h:219
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
void set_skip_edge_smoother(bool skip_edge)
Skip flag the edge smoother.
Definition: hiptmair_hybrid_smoother.h:241
TPotAlgebra::vector_type pot_vector_type
Vector type for the potential space.
Definition: hiptmair_hybrid_smoother.h:96
SmartPtr< matrix_operator_type > m_spSysMat
Edge-centered matrix of the original equation.
Definition: hiptmair_hybrid_smoother.h:113
const char * name() const
Returns the name.
Definition: hiptmair_hybrid_smoother.h:200
MatrixOperator< pot_matrix_type, pot_vector_type > pot_matrix_operator_type
Matrix Operator type for the potential space.
Definition: hiptmair_hybrid_smoother.h:100
TimeHarmonicNedelecHybridSmoother< TDomain, TAlgebra > this_type
Definition: hiptmair_hybrid_smoother.h:83
SmartPtr< DoFDistribution > m_spEdgeDD
DoF distribution for the Nedelec elements.
Definition: hiptmair_hybrid_smoother.h:124
SmartPtr< ApproximationSpace< TDomain > > m_spVertApproxSpace
Approximation space for the potential (vertex-centered) space.
Definition: hiptmair_hybrid_smoother.h:127
bool zero_mass_entry(double val)
Measure for numerically zero entries in the mass matrix.
Definition: hiptmair_hybrid_smoother.h:110
bool apply(vector_type &c, const vector_type &d)
Computes the correction.
Definition: hiptmair_hybrid_smoother_impl.h:455
bool m_bInit
Whether initialized.
Definition: hiptmair_hybrid_smoother.h:165
TAlgebra::matrix_type matrix_type
Matrix type.
Definition: hiptmair_hybrid_smoother.h:89
SmartPtr< ILinearIterator< vector_type > > m_spEdgeSmoother
Smoother for the edge dofs.
Definition: hiptmair_hybrid_smoother.h:135
SmartPtr< ILinearIterator< vector_type, vector_type > > clone()
Clone the smoother by copying the data.
Definition: hiptmair_hybrid_smoother.h:246
void compute_GtMG()
Computes the product .
Definition: hiptmair_hybrid_smoother_impl.h:128
SmartPtr< ILinearIterator< pot_vector_type > > m_spVertSmoother
Smoother for the vertex dofs.
Definition: hiptmair_hybrid_smoother.h:138
bool m_bSkipVertex
Definition: hiptmair_hybrid_smoother.h:168
bool init(SmartPtr< ILinearOperator< vector_type > > L)
We cannot initialize without the geometry (this version of init cannot not work)
Definition: hiptmair_hybrid_smoother.h:209
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
SmartPtr< pot_matrix_operator_type > m_spPotMat
Vertex-centered matrix for the potential.
Definition: hiptmair_hybrid_smoother.h:116
SmartPtr< EMDirichlet< TDomain, TAlgebra > > m_spDirichlet
Dirichlet boundary.
Definition: hiptmair_hybrid_smoother.h:141
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
SmartPtr< DoFDistribution > m_spVertDD
Level DoF distribution for the vertex centered grid func.
Definition: hiptmair_hybrid_smoother.h:129
void set_skip_vertex_smoother(bool skip_vertex)
Skip flag the vertex smoother.
Definition: hiptmair_hybrid_smoother.h:243
bool supports_parallel() const
Currently returns false because the computation of the potentials is not purely parallel up to now.
Definition: hiptmair_hybrid_smoother.h:203
VariableArray1< tEdgeInfo > m_vEdgeInfo
Storage for the information about the edge-vertex interconnections.
Definition: hiptmair_hybrid_smoother.h:162
bool m_bSkipEdge
Needed mainly for debugging: Whether to skip one of the stages.
Definition: hiptmair_hybrid_smoother.h:168
GridLevel m_GridLevel
Grid level (or surface) where to smooth.
Definition: hiptmair_hybrid_smoother.h:132
TimeHarmonicNedelecHybridSmoother(SmartPtr< ApproximationSpace< TDomain > > vertApproxSpace, SmartPtr< ILinearIterator< vector_type > > edgeSmoother, SmartPtr< ILinearIterator< pot_vector_type > > vertSmoother)
Constructor setting the approx. spaces and the default subsmoothers.
Definition: hiptmair_hybrid_smoother.h:173
TPotGridFunc * m_pPotCorRe
Vertex-centered grid function for the Re-part of potential corrections (this is a GridFunction to all...
Definition: hiptmair_hybrid_smoother.h:119
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
SmartPtr< IDamping< X, Y > > damping()
ParallelMatrix< SparseMatrix< double > > matrix_type
ParallelVector< Vector< double > > vector_type
#define UG_THROW(msg)
"and" reduction operation class for the conductivity condition
Definition: hiptmair_hybrid_smoother.h:302
static bool op(bool a, bool b)
Definition: hiptmair_hybrid_smoother.h:303
"or" reduction operation class for the conductivity condition
Definition: hiptmair_hybrid_smoother.h:297
static bool op(bool a, bool b)
Definition: hiptmair_hybrid_smoother.h:298
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_Dirichlet()
Definition: hiptmair_hybrid_smoother.h:156
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
bool is_init()
Definition: hiptmair_hybrid_smoother.h:153
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
char flags
flags: (0) Dirichlet, (1) conductive vertex 0, (2) conductive vertex 1, (3) init-ed
Definition: hiptmair_hybrid_smoother.h:148
void clear_flags()
Definition: hiptmair_hybrid_smoother.h:155
function ProblemDisc new(problemDesc, dom)