ug4
constrained_linear_iterator.h
Go to the documentation of this file.
1 /*
2  * Copyright (c) 2010-2015: G-CSC, Goethe University Frankfurt
3  * Author: Markus Breit
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__LIB_ALGEBRA__OPERATOR__INTERFACE__CONSTRAINED_LINEAR_ITERATOR__
34 #define __H__LIB_ALGEBRA__OPERATOR__INTERFACE__CONSTRAINED_LINEAR_ITERATOR__
35 
40 #include "preconditioner.h"
41 
42 #include <boost/core/enable_if.hpp>
43 #include <boost/type_traits/is_base_of.hpp>
44 
45 namespace ug {
46 
47 
63 template <typename TDomain, typename TAlgebra, typename TLinIt, typename = void>
64 class ConstrainedLinearIterator : public TLinIt
65 {
66  public:
68  typedef typename TAlgebra::vector_type vector_type;
69  typedef TLinIt base_type;
70 
71  using TLinIt::init;
72 
73  public:
76  : base_type(), m_spDomDisc(domDisc), m_time(0.0)
77  {}
78 
81  : base_type(parent), m_spDomDisc(parent.m_spDomDisc), m_time(parent.m_time)
82  {}
83 
86  {
88  }
89 
92 
93  protected:
95  virtual const char* name() const
96  {
97  return base_type::name();
98  }
99 
101  {
102  return base_type::init(J, u);
103  }
104 
106  {
107  return base_type::init(L);
108  }
109 
110  virtual bool apply(vector_type& c, const vector_type& d)
111  {
112  // perform normal step
113  if (!base_type::apply(c, d)) return false;
114 
115  // apply constraints on correction
116  if (!m_spDomDisc.valid())
117  UG_THROW("Domain discretization passed to ConstrainedLinearIterator is not valid.\n"
118  "Cannot apply any constraints.")
119 
120  gf_type* gf = dynamic_cast<gf_type*>(&c);
121  if (gf && gf->dof_distribution().valid())
122  {
123  size_t nConstr = m_spDomDisc->num_constraints();
124 
125  // first Dirichlet (hanging nodes might be constrained by Dirichlet nodes)
126  for (size_t i = 0; i < nConstr; i++)
127  if (m_spDomDisc->constraint(i)->type() & CT_DIRICHLET)
128  m_spDomDisc->constraint(i)->adjust_correction(c, gf->dof_distribution(), CT_DIRICHLET, m_time);
129 
130  // then other interpolation
131  for (size_t i = 0; i < nConstr; i++)
132  if (m_spDomDisc->constraint(i)->type() & CT_CONSTRAINTS)
133  m_spDomDisc->constraint(i)->adjust_correction(c, gf->dof_distribution(), CT_CONSTRAINTS, m_time);
134 
135  // then hanging nodes
136  for (size_t i = 0; i < nConstr; i++)
137  if (m_spDomDisc->constraint(i)->type() & CT_HANGING)
138  m_spDomDisc->constraint(i)->adjust_correction(c, gf->dof_distribution(), CT_HANGING, m_time);
139 
140  // and Dirichlet again (Dirichlet nodes might also be hanging)
141  for (size_t i = 0; i < nConstr; i++)
142  if (m_spDomDisc->constraint(i)->type() & CT_DIRICHLET)
143  m_spDomDisc->constraint(i)->adjust_correction(c, gf->dof_distribution(), CT_DIRICHLET, m_time);
144  }
145  else UG_THROW("Vector is not a grid function. This is not supported.")
146 
147  return true;
148  }
149 
151  {
153  return impl(c, d);
154  }
155 
156  private:
157  template <typename S, typename = void>//, typename otherDummy = void>
159  {
161 
163  : cli(_cli) {}
164 
166  {
167  // FIXME: This implementation is not correct.
168  // As the defect is calculated from possibly erroneous corrections, there is no telling
169  // how to adjust the defect. Except for setting it to zero everywhere, which, by the way,
170  // is NOT what adjust_defect does, in general. We would need something like
171  // adjust_defect_linearSolver() here to distinguish between the non-linear defect and
172  // the residual in the linear solver.
173  // We CANNOT use our own apply() here and then update the defect ourselves as the
174  // ILinearIterator does not have an update method nor any way of applying the underlying
175  // linear operator.
176 
177  // perform normal step
178  if (!cli.base_type::apply_update_defect(c, d)) return false;
179 
180  // apply constraints on correction and defect
181  if (!cli.m_spDomDisc.valid())
182  UG_THROW("Domain discretization passed to ConstrainedLinearIterator is not valid.\n"
183  "Cannot apply any constraints.")
184 
185  gf_type* gf = dynamic_cast<gf_type*>(&c);
186  if (gf && gf->dof_distribution().valid())
187  {
188  size_t nConstr = cli.m_spDomDisc->num_constraints();
189 
190  // corrections
191  // first Dirichlet (hanging nodes might be constrained by Dirichlet nodes)
192  for (size_t i = 0; i < nConstr; i++)
193  if (cli.m_spDomDisc->constraint(i)->type() & CT_DIRICHLET)
194  cli.m_spDomDisc->constraint(i)->adjust_correction(c, gf->dof_distribution(), CT_DIRICHLET, cli.m_time);
195 
196  // then other interpolation
197  for (size_t i = 0; i < nConstr; i++)
198  if (cli.m_spDomDisc->constraint(i)->type() & CT_CONSTRAINTS)
199  cli.m_spDomDisc->constraint(i)->adjust_correction(c, gf->dof_distribution(), CT_CONSTRAINTS, cli.m_time);
200 
201  // then hanging nodes
202  for (size_t i = 0; i < nConstr; i++)
203  if (cli.m_spDomDisc->constraint(i)->type() & CT_HANGING)
204  cli.m_spDomDisc->constraint(i)->adjust_correction(c, gf->dof_distribution(), CT_HANGING, cli.m_time);
205 
206  // and Dirichlet again (Dirichlet nodes might also be hanging)
207  for (size_t i = 0; i < nConstr; i++)
208  if (cli.m_spDomDisc->constraint(i)->type() & CT_DIRICHLET)
209  cli.m_spDomDisc->constraint(i)->adjust_correction(c, gf->dof_distribution(), CT_DIRICHLET, cli.m_time);
210 
211  // residuals
212  // hanging nodes first
213  for (size_t i = 0; i < nConstr; i++)
214  if (cli.m_spDomDisc->constraint(i)->type() & CT_HANGING)
215  cli.m_spDomDisc->constraint(i)->adjust_linear_residual(d, c, gf->dof_distribution(), CT_HANGING, cli.m_time);
216 
217  // then other constraints
218  for (size_t i = 0; i < nConstr; i++)
219  if (cli.m_spDomDisc->constraint(i)->type() & CT_CONSTRAINTS)
220  cli.m_spDomDisc->constraint(i)->adjust_linear_residual(d, c, gf->dof_distribution(), CT_CONSTRAINTS, cli.m_time);
221 
222  // and Dirichlet last
223  for (size_t i = 0; i < nConstr; i++)
224  if (cli.m_spDomDisc->constraint(i)->type() & CT_DIRICHLET)
225  cli.m_spDomDisc->constraint(i)->adjust_linear_residual(d, c, gf->dof_distribution(), CT_DIRICHLET, cli.m_time);
226 
227  }
228  else UG_THROW("Vector is not a grid function. This is not supported.")
229 
230  return true;
231  }
232  };
233 
234 
235  // special implementation for IPreconditioner
236  template <typename S>
237  struct apply_update_defect_impl<S, typename boost::enable_if<boost::is_base_of<IPreconditioner<TAlgebra>, S> >::type>
238  {
240 
242  : cli(_cli) {}
243 
245  {
246  // apply precond and constraints
247  if (!cli.apply(c, d)) return false;
248 
249  // calculate new defect from preconditioner defect operator
250  cli.m_spDefectOperator->apply_sub(d, c);
251 
252  return true;
253  }
254  };
255 
256  friend struct apply_update_defect_impl<TLinIt>;
257 
258 
259  public:
261  virtual bool supports_parallel() const
262  {
263  return base_type::supports_parallel();
264  }
265 
267  void set_time(number time)
268  {
269  m_time = time;
270  };
271 
272  protected:
276 };
277 
278 
279 
280 
281 } // end namespace ug
282 
283 
284 #endif // __H__LIB_ALGEBRA__OPERATOR__INTERFACE__CONSTRAINED_LINEAR_ITERATOR__
location name
Definition: checkpoint_util.lua:128
Definition: smart_pointer.h:108
Definition: constrained_linear_iterator.h:65
virtual bool init(SmartPtr< ILinearOperator< vector_type, vector_type > > J, const vector_type &u)
Definition: constrained_linear_iterator.h:100
ConstrainedLinearIterator(const ConstrainedLinearIterator< TDomain, TAlgebra, TLinIt > &parent)
clone constructor
Definition: constrained_linear_iterator.h:80
virtual const char * name() const
returns the name of iterator
Definition: constrained_linear_iterator.h:95
GridFunction< TDomain, TAlgebra > gf_type
Definition: constrained_linear_iterator.h:67
virtual bool supports_parallel() const
returns if parallel solving is supported
Definition: constrained_linear_iterator.h:261
virtual bool init(SmartPtr< ILinearOperator< vector_type, vector_type > > L)
Definition: constrained_linear_iterator.h:105
SmartPtr< ILinearIterator< vector_type > > m_spLinIt
Definition: constrained_linear_iterator.h:270
TLinIt base_type
Definition: constrained_linear_iterator.h:69
virtual bool apply(vector_type &c, const vector_type &d)
Definition: constrained_linear_iterator.h:110
virtual ~ConstrainedLinearIterator()
destructor
Definition: constrained_linear_iterator.h:91
TAlgebra::vector_type vector_type
Definition: constrained_linear_iterator.h:68
number m_time
Definition: constrained_linear_iterator.h:275
SmartPtr< IDomainDiscretization< TAlgebra > > m_spDomDisc
Definition: constrained_linear_iterator.h:274
virtual bool apply_update_defect(vector_type &c, vector_type &d)
Definition: constrained_linear_iterator.h:150
void set_time(number time)
setter for time
Definition: constrained_linear_iterator.h:267
ConstrainedLinearIterator(SmartPtr< IDomainDiscretization< TAlgebra > > domDisc)
constructor
Definition: constrained_linear_iterator.h:75
virtual SmartPtr< ILinearIterator< vector_type > > clone()
clone
Definition: constrained_linear_iterator.h:85
represents numerical solutions on a grid using an algebraic vector
Definition: grid_function.h:121
SmartPtr< DoFDistribution > dof_distribution()
Definition: grid_function.h:236
Interface for domain discretization.
Definition: domain_disc_interface.h:136
describes a linear mapping X->Y
Definition: linear_operator.h:80
virtual void init()
#define UG_THROW(msg)
Definition: error.h:57
double number
Definition: types.h:124
Definition: boost_serialization_routines.h:49
the ug namespace
@ CT_DIRICHLET
Definition: ass_tuner.h:59
@ CT_HANGING
Definition: ass_tuner.h:57
@ CT_CONSTRAINTS
Definition: ass_tuner.h:58
SmartPtr< T, FreePolicy > make_sp(T *inst)
returns a SmartPtr for the passed raw pointer
Definition: smart_pointer.h:836
Definition: constrained_linear_iterator.h:159
ConstrainedLinearIterator & cli
Definition: constrained_linear_iterator.h:160
apply_update_defect_impl(ConstrainedLinearIterator &_cli)
Definition: constrained_linear_iterator.h:162
bool operator()(vector_type &c, vector_type &d)
Definition: constrained_linear_iterator.h:165