ug4
metric_spaces.h
Go to the documentation of this file.
1 /*
2  * Copyright (c) 2014-2017: G-CSC, Goethe University Frankfurt
3  * Author: Arne Naegel
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__FUNCTION_SPACE__METRIC_SPACES_H_
34 #define __H__UG__LIB_DISC__FUNCTION_SPACE__METRIC_SPACES_H_
35 
36 // C++ includes
37 #include <vector>
38 #include <cmath>
39 
40 // UG includes
41 #include "common/common.h"
43 
46 
51 
52 #include "lib_disc/spatial_disc/user_data/linker/darcy_velocity_linker.h" // required for Darcy velocity linker
53 
54 namespace ug {
56 template<typename TVector>
58 {
59 
60 public:
61  virtual ~IBanachSpace() {}
62 
64  virtual double norm(TVector &x)
65  { return x.norm(); }
66 
67  virtual double distance(TVector& x, TVector& y)
68  {
69  SmartPtr<TVector> delta = x.clone();
70  *delta -= y;
71  return norm(*delta);
72  }
73 
74 };
75 
76 
77 
79 template<typename TGridFunction>
81  public IBanachSpace<typename TGridFunction::vector_type>
82 {
83 
84 public:
85  typedef typename TGridFunction::vector_type vector_type;
86  typedef TGridFunction grid_function_type;
87 
89  virtual ~IGridFunctionSpace() {}
90 
92  virtual double norm(TGridFunction& x) = 0;
93  virtual double norm2(TGridFunction& x) = 0;
94 
95  // { return x.norm(); }
96 
98  virtual double distance(TGridFunction& x, TGridFunction& y) = 0;
99  virtual double distance2(TGridFunction& x, TGridFunction& y) = 0;
100  /*{
101  SmartPtr<TGridFunction> delta = x.clone();
102  *delta -= y;
103  return norm(*delta);
104  }*/
105 
107  virtual double norm(vector_type &x)
108  {
109  TGridFunction* gfX=dynamic_cast< TGridFunction*>(&x);
110  UG_ASSERT(gfX!=NULL, "Huhh: GridFunction required!");
111  return norm(*gfX);
112  }
113 
115  virtual double distance(vector_type &x, vector_type &y)
116  { return distance(static_cast<TGridFunction &>(x), static_cast<TGridFunction &>(y)); }
117 
118  virtual double scaling() const
119  { return 1.0; }
120 
121  virtual std::string config_string() const
122  { return std::string("IGridFunctionSpace"); }
123 };
124 
125 template<typename TGridFunction>
126 class AlgebraicSpace : public IGridFunctionSpace<TGridFunction>
127 {
128  typedef typename TGridFunction::vector_type vector_type;
129  virtual ~AlgebraicSpace() {}
130 
131  virtual double norm(TGridFunction& x)
132  {return x.norm();}
133 
134  virtual double norm2(TGridFunction& x)
135  { double n = this->norm(x); return n*n; }
136 
137  virtual double distance(TGridFunction& x, TGridFunction& y)
138  { SmartPtr<TGridFunction> delta = x.clone(); *delta -= y; return delta->norm(); }
139 
140  virtual double distance2(TGridFunction& x, TGridFunction& y)
141  { double d = this->distance(x,y);; return d*d;}
142 };
143 
144 
145 
147 template <typename W>
149 {
150 public:
151  typedef W weight_type;
152 
154  : m_spWeight(SPNULL) {}
155 
157  : m_spWeight(spW) {}
158 
161  { m_spWeight = spWeight; }
162 
164  { return m_spWeight; }
165 
166 protected:
168 };
169 
170 
172 /*
173 class ITimeData {
174 
175 public:
176  virtual ~ITimeData() {};
177 
179  virtual void update_time_data(number dt) = 0;
180 
181  virtual bool is_time_dependent() const
182  { return false; }
183 
184 };
185 
186 */
187 
189 template <typename TGridFunction>
191  public IGridFunctionSpace<TGridFunction>
192 {
193 protected:
194  std::string m_fctNames;
195  const char* m_ssNames;
197 
198 public:
200  static const int dim=TGridFunction::dim;
201 
202  IComponentSpace(const char *fctNames)
203  : m_fctNames(fctNames), m_ssNames(NULL), m_quadorder(3){}
204 
205  IComponentSpace(const char *fctNames, int order)
206  : m_fctNames(fctNames), m_ssNames(NULL), m_quadorder(order) {}
207 
208  IComponentSpace(const char *fctNames, const char* ssNames, int order)
209  : m_fctNames(fctNames), m_ssNames(ssNames), m_quadorder(order) {}
210 
211  virtual ~IComponentSpace() {};
212 
213  // per convention, norm must return sqrt of norm2
214  virtual double norm(TGridFunction& uFine)
215  { return sqrt(norm2(uFine)); }
216 
217  // per convention, distance must return sqrt of distance2
218  virtual double distance(TGridFunction& uFine, TGridFunction& uCoarse)
219  { return sqrt(distance2(uFine, uCoarse)); }
220 
221  virtual double norm2(TGridFunction& uFine) = 0;
222  virtual double distance2(TGridFunction& uFine, TGridFunction& uCoarse) = 0;
223 
224 
225 public:
227  virtual std::string config_string() const
228  {
229  std::stringstream ss;
230 
231  if (this->m_ssNames)
232  ss << this->m_fctNames << ", " << this->m_ssNames << ", " << this->m_quadorder
233  << ", type=" << std::endl;
234 
235  else
236  ss << this->m_fctNames << ", (no name), " << this->m_quadorder
237  << ", type=" << std::endl;
238  return ss.str();
239  }
240 
241 };
242 
243 
246 template <typename TGridFunction>
248 : public IComponentSpace<TGridFunction>
249 {
250  public:
251  GridFunctionComponentSpace(const char* fctNames)
252  : IComponentSpace<TGridFunction>(fctNames) {}
253 
254  GridFunctionComponentSpace(const char* fctNames, const char* ssNames)
255  : IComponentSpace<TGridFunction>(fctNames, ssNames, 1) {}
256 
258 
259  virtual double norm2(TGridFunction& uFine)
260  {
261  ConstSmartPtr<DoFDistribution> dd = uFine.dof_distribution();
262 
263  // find function indices
264  FunctionGroup fg(dd->function_pattern());
265  try {fg.add(TokenizeTrimString(this->m_fctNames));}
266  UG_CATCH_THROW("Could not convert function names to function indices.");
267 
268  // find subset indices
269  SubsetGroup sg(dd->subset_handler());
270  if (this->m_ssNames)
271  {
272  try {sg.add(TokenizeTrimString(this->m_ssNames));}
273  UG_CATCH_THROW("Could not convert subset names to subset indices.");
274  }
275  else
276  sg.add_all();
277 
278 
279  // loop subsets
280  number sum = 0.0;
281  const size_t sgSz = sg.size();
282  for (size_t i = 0; i < sgSz; ++i)
283  {
284  int si = sg[i];
285  if (dd->max_dofs(VERTEX, si)) add_norm_values<Vertex>(sum, uFine, dd, sg[i], fg);
286  if (dd->max_dofs(EDGE, si)) add_norm_values<Edge>(sum, uFine, dd, sg[i], fg);
287  if (dd->max_dofs(FACE, si)) add_norm_values<Face>(sum, uFine, dd, sg[i], fg);
288  if (dd->max_dofs(VOLUME, si)) add_norm_values<Volume>(sum, uFine, dd, sg[i], fg);
289  }
290 
291 #ifdef UG_PARALLEL
292  if (pcl::NumProcs() > 1)
293  {
295  sum = pc.allreduce(sum, PCL_RO_SUM);
296  }
297 #endif
298 
299  return sum;
300  }
301 
302  virtual double distance2(TGridFunction& uFine, TGridFunction& uCoarse)
303  {
304  ConstSmartPtr<DoFDistribution> dd = uFine.dof_distribution();
305  UG_COND_THROW(dd != uCoarse.dof_distribution(),
306  "GridFunctionComponentSpace::distance2: GF1 DoF distro is not the same as for GF2.\n"
307  "This case is not implemented.");
308 
309  // find function indices
310  FunctionGroup fg(dd->function_pattern());
311  try {fg.add(TokenizeTrimString(this->m_fctNames));}
312  UG_CATCH_THROW("Could not convert function names to function indices.");
313 
314  // find subset indices
315  SubsetGroup sg(dd->subset_handler());
316  if (this->m_ssNames)
317  {
318  try {sg.add(TokenizeTrimString(this->m_ssNames));}
319  UG_CATCH_THROW("Could not convert subset names to subset indices.");
320  }
321  else
322  sg.add_all();
323 
324 
325  // loop subsets
326  number sum = 0.0;
327  const size_t sgSz = sg.size();
328  for (size_t i = 0; i < sgSz; ++i)
329  {
330  int si = sg[i];
331  if (dd->max_dofs(VERTEX, si)) add_distance_values<Vertex>(sum, uFine, uCoarse, dd, sg[i], fg);
332  if (dd->max_dofs(EDGE, si)) add_distance_values<Edge>(sum, uFine, uCoarse, dd, sg[i], fg);
333  if (dd->max_dofs(FACE, si)) add_distance_values<Face>(sum, uFine, uCoarse, dd, sg[i], fg);
334  if (dd->max_dofs(VOLUME, si)) add_distance_values<Volume>(sum, uFine, uCoarse, dd, sg[i], fg);
335  }
336 
337 #ifdef UG_PARALLEL
338  if (pcl::NumProcs() > 1)
339  {
341  sum = pc.allreduce(sum, PCL_RO_SUM);
342  }
343 #endif
344 
345  return sum;
346  }
347 
348  private:
349  template <typename TBaseElem>
351  (
352  number& sum,
353  const TGridFunction& uFine,
354  const TGridFunction& uCoarse,
356  int si,
357  const FunctionGroup& fg
358  ) const
359  {
360  const SurfaceView& sv = *dd->surface_view();
361  const MultiGrid& mg = *dd->multi_grid();
362 
363  const number nFct = fg.size();
364 
365  // iterate all elements (including SHADOW_RIM_COPY!)
367  iter = dd->template begin<TBaseElem>(si, SurfaceView::ALL);
368  iterEnd = dd->template end<TBaseElem>(si, SurfaceView::ALL);
369  std::vector<DoFIndex> vInd;
370  for (; iter != iterEnd; ++iter)
371  {
372  TBaseElem* elem = *iter;
373  if (sv.is_contained(elem, dd->grid_level(), SurfaceView::SHADOW_RIM_COPY))
374  {
375  if (mg.num_children<TBaseElem>(elem) > 0)
376  {
377  TBaseElem* child = mg.get_child<TBaseElem>(elem, 0);
378  if (sv.is_contained(child, dd->grid_level(), SurfaceView::SURFACE_RIM))
379  continue;
380  }
381  }
382 
383  for (size_t fi = 0; fi < nFct; ++fi)
384  {
385  dd->inner_dof_indices(elem, fg[fi], vInd);
386 
387  const size_t nDof = vInd.size();
388  for (size_t dof = 0; dof < nDof; ++dof)
389  {
390  const number tmp = DoFRef(uFine, vInd[dof]) - DoFRef(uCoarse, vInd[dof]);
391  sum += tmp*tmp;
392  }
393  }
394  }
395  // note: no duplicate additions possible
396  }
397 
398  template <typename TBaseElem>
400  (
401  number& sum,
402  const TGridFunction& uFine,
404  int si,
405  const FunctionGroup& fg
406  ) const
407  {
408  const SurfaceView& sv = *dd->surface_view();
409  const MultiGrid& mg = *dd->multi_grid();
410 
411  const number nFct = fg.size();
412 
413  // iterate all elements (including SHADOW_RIM_COPY!)
415  iter = dd->template begin<TBaseElem>(si, SurfaceView::ALL);
416  iterEnd = dd->template end<TBaseElem>(si, SurfaceView::ALL);
417  std::vector<DoFIndex> vInd;
418  for (; iter != iterEnd; ++iter)
419  {
420  TBaseElem* elem = *iter;
421  if (sv.is_contained(elem, dd->grid_level(), SurfaceView::SHADOW_RIM_COPY))
422  {
423  if (mg.num_children<TBaseElem>(elem) > 0)
424  {
425  TBaseElem* child = mg.get_child<TBaseElem>(elem, 0);
426  if (sv.is_contained(child, dd->grid_level(), SurfaceView::SURFACE_RIM))
427  continue;
428  }
429  }
430 
431  for (size_t fi = 0; fi < nFct; ++fi)
432  {
433  dd->inner_dof_indices(elem, fg[fi], vInd);
434 
435  const size_t nDof = vInd.size();
436  for (size_t dof = 0; dof < nDof; ++dof)
437  {
438  const number tmp = DoFRef(uFine, vInd[dof]);
439  sum += tmp*tmp;
440  }
441  }
442  }
443  // note: no duplicate additions possible
444  }
445 };
446 
447 
448 
450 template <typename TGridFunction>
452  public IGridFunctionSpace<TGridFunction>
453 {
454 public:
457 
460  : m_spSpatialSpace(spSpace), m_tScale(tScale) {};
461 
463  virtual ~TimeDependentSpace() {};
464 
465  using base_type::norm;
466  using base_type::distance;
467  using base_type::scaling;
468 
470  double scaling() const
471  { return (m_spSpatialSpace->scaling()*m_tScale); }
472 
475  { m_tScale = tScale; }
476 
478  std::string config_string() const
479  {
480  std::stringstream ss;
481  ss << "TimeDependentSpace for " << std::endl;
482  ss << m_spSpatialSpace->config_string() << std::endl;
483  return ss.str();
484  }
485 
486 protected:
489 };
490 
491 
493 template <typename TGridFunction>
495  public IComponentSpace<TGridFunction>,
496  public IObjectWithWeights<typename L2DistIntegrand<TGridFunction>::weight_type >
497 {
498 public:
502 
504  L2ComponentSpace(const char *fctNames)
505  : base_type(fctNames), weighted_obj_type(make_sp(new ConstUserNumber<TGridFunction::dim>(1.0))) {};
506 
507  L2ComponentSpace(const char *fctNames, int order)
508  : base_type(fctNames, order), weighted_obj_type(make_sp(new ConstUserNumber<TGridFunction::dim>(1.0))) {};
509 
510  L2ComponentSpace(const char *fctNames, int order, double weight, const char* ssNames=0)
511  : base_type(fctNames, ssNames, order), weighted_obj_type(make_sp(new ConstUserNumber<TGridFunction::dim>(weight))) {};
512 
513  L2ComponentSpace(const char *fctNames, int order, ConstSmartPtr<weight_type> spWeight, const char* ssNames=0)
514  : base_type(fctNames, ssNames, order), weighted_obj_type(spWeight) {};
515 
518 
521 
526 
528  double norm2(TGridFunction& uFine)
530 
532  double distance2(TGridFunction& uFine, TGridFunction& uCoarse)
533  { return L2Distance2(uFine, base_type::m_fctNames.c_str(), uCoarse, base_type::m_fctNames.c_str(),
535 
536 
537 
538 };
539 
541 template <typename TGridFunction>
543  public IComponentSpace<TGridFunction>,
544  public IObjectWithWeights<typename L2DistIntegrand<TGridFunction>::weight_type >
545 {
546 public:
550 
552  L2QuotientSpace(const char *fctNames)
553  : base_type(fctNames), weighted_obj_type(make_sp(new ConstUserNumber<TGridFunction::dim>(1.0))) {};
554 
555  L2QuotientSpace(const char *fctNames, int order)
556  : base_type(fctNames, order), weighted_obj_type(make_sp(new ConstUserNumber<TGridFunction::dim>(1.0))) {};
557 
558  L2QuotientSpace(const char *fctNames, int order, double weight, const char* ssNames=0)
559  : base_type(fctNames, ssNames, order), weighted_obj_type(make_sp(new ConstUserNumber<TGridFunction::dim>(weight))) {};
560 
561  L2QuotientSpace(const char *fctNames, int order, ConstSmartPtr<weight_type> spWeight, const char* ssNames=0)
562  : base_type(fctNames, ssNames, order), weighted_obj_type(spWeight) {};
563 
566 
569 
574 
576  double norm2(TGridFunction& u)
577  {
578  typedef ConstUserNumber<TGridFunction::dim> MyConstUserData;
580 
581  SPUserData spConst= make_sp(new MyConstUserData(1.0));
582  number Meas = Integral(spConst, u, base_type::m_ssNames, 0.0, 1, "best");
584 
585  std::cerr << "Average:=" << uAvg <<"/" << Meas << " = " << uAvg/Meas << std::endl;
586  SPUserData spAvg = make_sp(new MyConstUserData(uAvg/Meas));
587  double qnorm = L2Error(spAvg, u, base_type::m_fctNames.c_str(), 0.0, base_type::m_quadorder, base_type::m_ssNames);
588 
589  return qnorm*qnorm;
590  }
591 
593  double distance2(TGridFunction& uFine, TGridFunction& uCoarse)
594  {
595  typedef ConstUserNumber<TGridFunction::dim> MyConstUserData;
597 
598  SPUserData spConst= make_sp(new MyConstUserData(1.0));
599  number Meas = Integral(spConst, uFine, base_type::m_ssNames, 0.0, 1, "best");
602 
603  std::cerr << "Average:=(" << avgFine << "-" << avgCoarse<<")/" << Meas << " = " << (avgFine-avgCoarse)/Meas << std::endl;
604  return L2Distance2(uFine, base_type::m_fctNames.c_str(),
605  uCoarse, base_type::m_fctNames.c_str(),
607  (avgFine-avgCoarse)/Meas);
608  }
609 
610 
611 
612 };
613 
614 
616 template <typename TGridFunction>
618 : public IComponentSpace<TGridFunction>,
619  public IObjectWithWeights<typename H1SemiDistIntegrand<TGridFunction>::weight_type >
620 {
621 public:
625 
626 
627  H1SemiComponentSpace(const char *fctNames)
628  : base_type(fctNames), weighted_obj_type(make_sp(new ConstUserMatrix<TGridFunction::dim>(1.0))) {};
629 
630  H1SemiComponentSpace(const char *fctNames, int order)
631  : base_type(fctNames, order), weighted_obj_type(make_sp(new ConstUserMatrix<TGridFunction::dim>(1.0))) {};
632 
633  H1SemiComponentSpace(const char *fctNames, int order, number weight, const char* ssNames=0)
634  : base_type(fctNames, ssNames, order), weighted_obj_type(make_sp(new ConstUserMatrix<TGridFunction::dim>(weight))) {};
635 
636  H1SemiComponentSpace(const char *fctNames, int order, ConstSmartPtr<weight_type> spWeight, const char* ssNames=0)
637  : base_type(fctNames, ssNames, order), weighted_obj_type(spWeight) {};
638 
639  H1SemiComponentSpace(const char *fctNames, int order, const char* ssNames, ConstSmartPtr<weight_type> spWeight)
640  : base_type(fctNames, ssNames, order), weighted_obj_type(spWeight) {};
641 
644 
645 
650 
652  double norm2(TGridFunction& uFine)
653  { return H1SemiNorm2<TGridFunction>(uFine, base_type::m_fctNames.c_str(), base_type::m_quadorder, NULL, weighted_obj_type::m_spWeight); }
654 
656  double distance2(TGridFunction& uFine, TGridFunction& uCoarse)
657  { return H1SemiDistance2<TGridFunction>(uFine, base_type::m_fctNames.c_str(), uCoarse, base_type::m_fctNames.c_str(), base_type::m_quadorder, m_spWeight); }
658 
663 
664 };
665 
666 
667 
668 /*
670 template <typename TGridFunction>
671 class KineticEnergyComponentSpace
672 : public IComponentSpace<TGridFunction>,
673  public IObjectWithWeights<typename L2Integrand<TGridFunction>::weight_type >
674 {
675 public:
676  typedef IComponentSpace<TGridFunction> base_type;
677  typedef typename L2Integrand<TGridFunction>::weight_type weight_type;
678  typedef IObjectWithWeights<weight_type> weighted_obj_type;
679  typedef DarcyVelocityLinker<TGridFunction::dim> velocity_type;
680 
681 
682  KineticEnergyComponentSpace(SmartPtr<velocity_type> spVelocity, const char *fctNames)
683  : base_type(fctNames), weighted_obj_type(make_sp(new ConstUserNumber<TGridFunction::dim>(1.0))), m_spVelocity(spVelocity) {};
684 
685  KineticEnergyComponentSpace(SmartPtr<velocity_type> spVelocity, const char *fctNames, int order)
686  : base_type(fctNames, order), weighted_obj_type(make_sp(new ConstUserNumber<TGridFunction::dim>(1.0))), m_spVelocity(spVelocity) {};
687 
688  KineticEnergyComponentSpace(SmartPtr<velocity_type> spVelocity, const char *fctNames, int order, number weight, const char* ssNames=0)
689  : base_type(fctNames, ssNames, order), weighted_obj_type(make_sp(new ConstUserNumber<TGridFunction::dim>(weight))), m_spVelocity(spVelocity) {};
690 
691  KineticEnergyComponentSpace(SmartPtr<velocity_type> spVelocity, const char *fctNames, int order, ConstSmartPtr<weight_type> spWeight, const char* ssNames=0)
692  : base_type(fctNames, ssNames, order), weighted_obj_type(spWeight), m_spVelocity(spVelocity) {};
693 
694  KineticEnergyComponentSpace(SmartPtr<velocity_type> spVelocity, const char *fctNames, int order, const char* ssNames, ConstSmartPtr<weight_type> spWeight)
695  : base_type(fctNames, ssNames, order), weighted_obj_type(spWeight), m_spVelocity(spVelocity) {};
696 
698  ~KineticEnergyComponentSpace() {};
699 
700 
701  using IComponentSpace<TGridFunction>::m_quadorder;
702  using IComponentSpace<TGridFunction>::norm;
703  using IComponentSpace<TGridFunction>::distance;
704 
706  double norm2(TGridFunction& u)
707  {
708  const char* subsets = NULL;
709  UserDataIntegrandSq<MathVector<TGridFunction::dim>, TGridFunction> integrand2(m_spVelocity, &u, 0.0);
710  return IntegrateSubsets(integrand2, u, subsets, m_quadorder);
711  }
712 
714  double distance2(TGridFunction& uFine, TGridFunction& uCoarse)
715  {
716  // UG_THROW("Not implemented!");
717  const char* subsets = NULL;
718  UserDataDistIntegrandSq<MathVector<TGridFunction::dim>, TGridFunction> integrand2(m_spVelocity, uFine, 0, uCoarse, 0);
719  return IntegrateSubsets(integrand2, uFine, subsets, m_quadorder);
720  }
721 
723  using weighted_obj_type::set_weight;
724  using weighted_obj_type::get_weight;
725  using weighted_obj_type::m_spWeight;
726 
727  void set_velocity(SmartPtr<velocity_type> spVelocity)
728  { m_spVelocity = spVelocity;}
729 
730 protected:
731  SmartPtr<velocity_type> m_spVelocity;
732 
733 };
734 */
735 
737 
738 template <typename TGridFunction>
740 : public IComponentSpace<TGridFunction>,
741  public IObjectWithWeights<typename H1EnergyDistIntegrand<TGridFunction>::weight_type >
742 {
743 public:
747  // typedef DarcyVelocityLinker<TGridFunction::dim> velocity_type;
748  static const int dim = TGridFunction::dim;
750 
751  H1EnergyComponentSpace(const char *fctNames)
752  : base_type(fctNames), weighted_obj_type(make_sp(new ConstUserMatrix<TGridFunction::dim>(1.0))), m_spVelocity(SPNULL) {};
753 
754  H1EnergyComponentSpace(const char *fctNames, int order)
755  : base_type(fctNames, order), weighted_obj_type(make_sp(new ConstUserMatrix<TGridFunction::dim>(1.0))), m_spVelocity(SPNULL) {};
756 
757  H1EnergyComponentSpace(const char *fctNames, int order, number weight, const char* ssNames=0)
758  : base_type(fctNames, ssNames, order), weighted_obj_type(make_sp(new ConstUserMatrix<TGridFunction::dim>(weight))), m_spVelocity(SPNULL) {};
759 
760  H1EnergyComponentSpace(const char *fctNames, int order, ConstSmartPtr<weight_type> spWeight, const char* ssNames=0)
761  : base_type(fctNames, ssNames, order), weighted_obj_type(spWeight), m_spVelocity(SPNULL) {};
762 
763  /*H1EnergyComponentSpace(const char *fctNames, int order, const char* ssNames, ConstSmartPtr<weight_type> spWeight)
764  : base_type(fctNames, ssNames, order), weighted_obj_type(spWeight), m_spVelocity(SPNULL) {};*/
765 
768 
769 
774 
776  double norm2(TGridFunction& uFine)
777  {
778  if (m_spVelocity.valid()) {
779  //const char* subsets = NULL; // [q^2]
780  UserDataIntegrandSq<MathVector<TGridFunction::dim>, TGridFunction> integrand2(m_spVelocity, &uFine, 0.0);
782  } else {
783  return H1EnergyNorm2<TGridFunction>(uFine, base_type::m_fctNames.c_str(), base_type::m_quadorder, base_type::m_ssNames, weighted_obj_type::m_spWeight);
784  }
785  }
786 
788  double distance2(TGridFunction& uFine, TGridFunction& uCoarse)
789  { return H1EnergyDistance2<TGridFunction>(uFine, base_type::m_fctNames.c_str(), uCoarse, base_type::m_fctNames.c_str(), base_type::m_quadorder,base_type::m_ssNames, weighted_obj_type::m_spWeight); }
790 
794 
796  { m_spVelocity = spVelocity;}
797 
798 protected:
800 
801 };
802 
804 template <typename TGridFunction>
806  public IComponentSpace<TGridFunction>
807 {
808 public:
810 
811  H1ComponentSpace(const char *fctNames) : base_type(fctNames) {};
812  H1ComponentSpace(const char *fctNames, int order) : base_type(fctNames, order) {};
813  H1ComponentSpace(const char *fctNames, const char* ssNames, int order) : base_type(fctNames, ssNames, order) {};
815 
818 
820  double norm2(TGridFunction& uFine)
821  { return H1Norm2<TGridFunction>(uFine, base_type::m_fctNames.c_str(), base_type::m_quadorder, base_type::m_ssNames); }
822 
824  double distance2(TGridFunction& uFine, TGridFunction& uCoarse)
825  { return H1Distance2<TGridFunction>(uFine, base_type::m_fctNames.c_str(), uCoarse, base_type::m_fctNames.c_str(), base_type::m_quadorder, base_type::m_ssNames); }
826 
827 };
828 
829 
830 
832 
835 template <typename TGridFunction>
836 class CompositeSpace : public IGridFunctionSpace<TGridFunction>
837 {
838 public:
842  typedef std::pair<SmartPtr<obj_type>, number> weighted_obj_type;
843 
845  // virtual ~CompositeSpace() {};
846 
847  using base_type::norm;
848  using base_type::distance;
849 
851  double norm(TGridFunction& uFine)
852  { return(sqrt(norm2(uFine))); }
853 
855  double norm2(TGridFunction& uFine)
856  {
857  number unorm2 = 0.0;
858  for (typename std::vector<weighted_obj_type>::iterator it = m_spWeightedSubspaces.begin();
859  it!= m_spWeightedSubspaces.end(); ++it)
860  {
861  double snorm2 = it->first->norm2(uFine);
862  unorm2 += it->second * snorm2; // scaling
863  UG_LOG("composite-norm2:\t" << snorm2 << "\t*\t" << it->second
864  << "\t=\t" << it->second * snorm2 << std::endl);
865  }
866  UG_LOG("composite-norm2-final:\t" << unorm2 << std::endl);
867  return unorm2;
868  }
869 
871  double distance2(TGridFunction& uFine, TGridFunction& uCoarse)
872  {
873  number unorm2 = 0.0;
874  for (typename std::vector<weighted_obj_type>::iterator it = m_spWeightedSubspaces.begin();
875  it!= m_spWeightedSubspaces.end(); ++it)
876  {
877  double sdist2 = it->first->distance2(uFine, uCoarse);
878  unorm2 += it->second * sdist2; // scaling
879  UG_LOG("composite-dist2:\t" << sdist2 << "\t*\t" << it->second
880  << "\t=\t" << it->second * sdist2 << std::endl);
881  }
882  UG_LOG("composite-dist2-final:\t" << unorm2 << std::endl);
883  return unorm2;
884  }
885 
887  double distance(TGridFunction& uFine, TGridFunction& uCoarse)
888  { return sqrt(distance2(uFine, uCoarse)); }
889 
891  void add(SmartPtr<obj_type> spSubSpace)
892  { m_spWeightedSubspaces.push_back(std::make_pair(spSubSpace, 1.0)); }
893 
895  void add(SmartPtr<obj_type> spSubSpace, number sigma)
896  { m_spWeightedSubspaces.push_back(std::make_pair(spSubSpace, sigma)); }
897 
899  std::string config_string() const
900  {
901  std::stringstream ss;
902  ss << "CompositeSpace:" << std::endl;
903 
904  for (typename std::vector<weighted_obj_type>::const_iterator it = m_spWeightedSubspaces.begin();
905  it!= m_spWeightedSubspaces.end(); ++it)
906  { ss << it->first->config_string(); }
907 
908  return ss.str();
909  }
910 
911 
914  {
915  for (typename std::vector<weighted_obj_type>::iterator it = m_spWeightedSubspaces.begin();
916  it!= m_spWeightedSubspaces.end(); ++it)
917  {
918  SmartPtr<time_dependent_obj_type> spSpaceT = it->first.template cast_dynamic<time_dependent_obj_type> ();
919  if (spSpaceT.valid()) spSpaceT->update_time_data(t);
920  }
921  }
922 
924  bool is_time_dependent() const
925  {
926  for (typename std::vector<weighted_obj_type>::const_iterator it = m_spWeightedSubspaces.begin();
927  it!= m_spWeightedSubspaces.end(); ++it)
928  {
929  SmartPtr<time_dependent_obj_type> spSpaceT = it->first.template cast_dynamic<time_dependent_obj_type>();
930  if (spSpaceT.valid()) return true;
931  }
932  return false;
933  }
934 
935  const std::vector<weighted_obj_type> &get_subspaces() const
936  { return m_spWeightedSubspaces; }
937 
938 protected:
939  std::vector<weighted_obj_type> m_spWeightedSubspaces;
940 
941 };
942 
943 /*
944 template <typename TGridFunction>
945 class TimeDependentCompositeSpace
946  : public CompositeSpace<TGridFunction>, public ITimeDependentSpace<TGridFunction>
947 {
948 protected:
949  using CompositeSpace<TGridFunction>::m_spSubspaces;
950 
951 public:
952  typedef typename IComponentSpace<TGridFunction> obj_type;
953 
954 
956  void add(SmartPtr<obj_type> spSubSpace)
957  { m_spSubspaces.push_back(make_sp(new ITimeDependentSpace(spSubSpace))); }
958 
960  void update_time_data(number t)
961  {
962  for (typename std::vector<SmartPtr<obj_type> >::iterator it = m_spSubspaces.begin();
963  it!= m_spSubspaces.end(); ++it)
964  { if ((*it)->is_time_dependent()) (*it)->update_time_data(t); }
965  }
966 
968  bool is_time_dependent() const
969  {
970  for (typename std::vector<SmartPtr<obj_type> >::const_iterator it = m_spSubspaces.begin();
971  it!= m_spSubspaces.end(); ++it)
972  { if ((*it)->is_time_dependent()) return true; }
973  return false;
974  }
975 };
976 
977  */
978 
979 } // namespace ug
980 #endif
bool valid() const
returns true if the pointer is valid, false if not.
Definition: smart_pointer.h:206
Definition: pcl_process_communicator.h:70
void allreduce(const void *sendBuf, void *recBuf, int count, DataType type, ReduceOperation op) const
performs MPI_Allreduce on the processes of the communicator.
Definition: pcl_process_communicator.cpp:318
Definition: metric_spaces.h:127
virtual double distance(TGridFunction &x, TGridFunction &y)
distance (for grid functions)
Definition: metric_spaces.h:137
TGridFunction::vector_type vector_type
Definition: metric_spaces.h:128
virtual ~AlgebraicSpace()
Definition: metric_spaces.h:129
virtual double distance2(TGridFunction &x, TGridFunction &y)
Definition: metric_spaces.h:140
virtual double norm2(TGridFunction &x)
Definition: metric_spaces.h:134
virtual double norm(TGridFunction &x)
norm (for grid functions)
Definition: metric_spaces.h:131
Defines a composite space, (i.e., additive composition from other spaces)
Definition: metric_spaces.h:837
IGridFunctionSpace< TGridFunction > base_type
Definition: metric_spaces.h:839
TimeDependentSpace< TGridFunction > time_dependent_obj_type
Definition: metric_spaces.h:841
std::string config_string() const
print config string
Definition: metric_spaces.h:899
double norm(TGridFunction &uFine)
norm (for grid functions)
Definition: metric_spaces.h:851
double norm2(TGridFunction &uFine)
Definition: metric_spaces.h:855
bool is_time_dependent() const
Check, if any object is time-dependent.
Definition: metric_spaces.h:924
double distance2(TGridFunction &uFine, TGridFunction &uCoarse)
Definition: metric_spaces.h:871
std::vector< weighted_obj_type > m_spWeightedSubspaces
Definition: metric_spaces.h:939
void add(SmartPtr< obj_type > spSubSpace)
add space to composite (with weight 1.0)
Definition: metric_spaces.h:891
std::pair< SmartPtr< obj_type >, number > weighted_obj_type
Definition: metric_spaces.h:842
void update_time_data(number t)
Forward update to all members.
Definition: metric_spaces.h:913
const std::vector< weighted_obj_type > & get_subspaces() const
Definition: metric_spaces.h:935
void add(SmartPtr< obj_type > spSubSpace, number sigma)
add space to composite (with variable weight)
Definition: metric_spaces.h:895
double distance(TGridFunction &uFine, TGridFunction &uCoarse)
Definition: metric_spaces.h:887
CompositeSpace()
Definition: metric_spaces.h:844
IComponentSpace< TGridFunction > obj_type
Definition: metric_spaces.h:840
constant matrix user data
Definition: const_user_data.h:232
constant scalar user data
Definition: const_user_data.h:153
Definition: function_group.h:52
void add(size_t fct)
adds a function by id to this group
Definition: function_group.cpp:79
size_t size() const
number of functions in this group
Definition: function_group.h:116
Definition: metric_spaces.h:249
virtual double distance2(TGridFunction &uFine, TGridFunction &uCoarse)
Definition: metric_spaces.h:302
GridFunctionComponentSpace(const char *fctNames, const char *ssNames)
Definition: metric_spaces.h:254
void add_distance_values(number &sum, const TGridFunction &uFine, const TGridFunction &uCoarse, ConstSmartPtr< DoFDistribution > dd, int si, const FunctionGroup &fg) const
Definition: metric_spaces.h:351
virtual double norm2(TGridFunction &uFine)
Definition: metric_spaces.h:259
void add_norm_values(number &sum, const TGridFunction &uFine, ConstSmartPtr< DoFDistribution > dd, int si, const FunctionGroup &fg) const
Definition: metric_spaces.h:400
GridFunctionComponentSpace(const char *fctNames)
Definition: metric_spaces.h:251
virtual ~GridFunctionComponentSpace()
Definition: metric_spaces.h:257
Definition: metric_spaces.h:807
~H1ComponentSpace()
Definition: metric_spaces.h:814
H1ComponentSpace(const char *fctNames, const char *ssNames, int order)
Definition: metric_spaces.h:813
H1ComponentSpace(const char *fctNames)
Definition: metric_spaces.h:811
double norm2(TGridFunction &uFine)
norm (for grid functions)
Definition: metric_spaces.h:820
H1ComponentSpace(const char *fctNames, int order)
Definition: metric_spaces.h:812
IComponentSpace< TGridFunction > base_type
Definition: metric_spaces.h:809
double distance2(TGridFunction &uFine, TGridFunction &uCoarse)
norm (for grid functions)
Definition: metric_spaces.h:824
Definition: metric_spaces.h:742
IObjectWithWeights< weight_type > weighted_obj_type
Definition: metric_spaces.h:746
H1SemiDistIntegrand< TGridFunction >::weight_type weight_type
Definition: metric_spaces.h:745
H1EnergyComponentSpace(const char *fctNames, int order, ConstSmartPtr< weight_type > spWeight, const char *ssNames=0)
Definition: metric_spaces.h:760
double norm2(TGridFunction &uFine)
Definition: metric_spaces.h:776
SmartPtr< velocity_type > m_spVelocity
Definition: metric_spaces.h:799
double distance2(TGridFunction &uFine, TGridFunction &uCoarse)
Definition: metric_spaces.h:788
IComponentSpace< TGridFunction > base_type
Definition: metric_spaces.h:744
UserData< MathVector< dim >, dim > velocity_type
Definition: metric_spaces.h:749
static const int dim
Definition: metric_spaces.h:748
~H1EnergyComponentSpace()
DTOR.
Definition: metric_spaces.h:767
H1EnergyComponentSpace(const char *fctNames)
Definition: metric_spaces.h:751
H1EnergyComponentSpace(const char *fctNames, int order)
Definition: metric_spaces.h:754
H1EnergyComponentSpace(const char *fctNames, int order, number weight, const char *ssNames=0)
Definition: metric_spaces.h:757
void set_velocity(SmartPtr< velocity_type > spVelocity)
Definition: metric_spaces.h:795
Definition: metric_spaces.h:620
H1SemiDistIntegrand< TGridFunction >::weight_type weight_type
Definition: metric_spaces.h:623
IComponentSpace< TGridFunction > base_type
Definition: metric_spaces.h:622
ConstSmartPtr< weight_type > m_spWeight
Definition: metric_spaces.h:167
H1SemiComponentSpace(const char *fctNames, int order, ConstSmartPtr< weight_type > spWeight, const char *ssNames=0)
Definition: metric_spaces.h:636
IObjectWithWeights< weight_type > weighted_obj_type
Definition: metric_spaces.h:624
H1SemiComponentSpace(const char *fctNames, int order, const char *ssNames, ConstSmartPtr< weight_type > spWeight)
Definition: metric_spaces.h:639
H1SemiComponentSpace(const char *fctNames, int order)
Definition: metric_spaces.h:630
~H1SemiComponentSpace()
DTOR.
Definition: metric_spaces.h:643
H1SemiComponentSpace(const char *fctNames, int order, number weight, const char *ssNames=0)
Definition: metric_spaces.h:633
H1SemiComponentSpace(const char *fctNames)
Definition: metric_spaces.h:627
double norm2(TGridFunction &uFine)
Definition: metric_spaces.h:652
double distance2(TGridFunction &uFine, TGridFunction &uCoarse)
Definition: metric_spaces.h:656
Abstract base class for (algebraic) vectors.
Definition: metric_spaces.h:58
virtual double norm(TVector &x)
euclidean norm (default)
Definition: metric_spaces.h:64
virtual double distance(TVector &x, TVector &y)
Definition: metric_spaces.h:67
virtual ~IBanachSpace()
Definition: metric_spaces.h:61
Estimate the error (based on the difference between two grid functions)
Definition: metric_spaces.h:192
virtual std::string config_string() const
print config string
Definition: metric_spaces.h:227
int m_quadorder
Definition: metric_spaces.h:196
IGridFunctionSpace< TGridFunction > base_type
Definition: metric_spaces.h:199
virtual double distance2(TGridFunction &uFine, TGridFunction &uCoarse)=0
IComponentSpace(const char *fctNames, int order)
Definition: metric_spaces.h:205
IComponentSpace(const char *fctNames, const char *ssNames, int order)
Definition: metric_spaces.h:208
IComponentSpace(const char *fctNames)
Definition: metric_spaces.h:202
virtual double norm(TGridFunction &uFine)
norm (for grid functions)
Definition: metric_spaces.h:214
virtual double norm2(TGridFunction &uFine)=0
virtual ~IComponentSpace()
Definition: metric_spaces.h:211
std::string m_fctNames
Definition: metric_spaces.h:194
virtual double distance(TGridFunction &uFine, TGridFunction &uCoarse)
distance (for grid functions)
Definition: metric_spaces.h:218
const char * m_ssNames
Definition: metric_spaces.h:195
static const int dim
Definition: metric_spaces.h:200
Abstract base class for grid functions.
Definition: metric_spaces.h:82
virtual double distance(vector_type &x, vector_type &y)
OVERRIDE distance (for vectors)
Definition: metric_spaces.h:115
virtual ~IGridFunctionSpace()
DTOR.
Definition: metric_spaces.h:89
virtual std::string config_string() const
Definition: metric_spaces.h:121
virtual double distance2(TGridFunction &x, TGridFunction &y)=0
virtual double norm(vector_type &x)
OVERRIDE norm (for vectors)
Definition: metric_spaces.h:107
TGridFunction grid_function_type
Definition: metric_spaces.h:86
virtual double norm(TGridFunction &x)=0
norm (for grid functions)
TGridFunction::vector_type vector_type
Definition: metric_spaces.h:85
virtual double distance(TGridFunction &x, TGridFunction &y)=0
distance (for grid functions)
virtual double scaling() const
Definition: metric_spaces.h:118
virtual double norm2(TGridFunction &x)=0
Definition: metric_spaces.h:149
W weight_type
Definition: metric_spaces.h:151
ConstSmartPtr< weight_type > m_spWeight
Definition: metric_spaces.h:167
void set_weight(ConstSmartPtr< weight_type > spWeight)
for weighted norms
Definition: metric_spaces.h:160
IObjectWithWeights()
Definition: metric_spaces.h:153
IObjectWithWeights(ConstSmartPtr< weight_type > spW)
Definition: metric_spaces.h:156
ConstSmartPtr< weight_type > get_weight()
Definition: metric_spaces.h:163
Definition: metric_spaces.h:497
L2DistIntegrand< TGridFunction >::weight_type weight_type
Definition: metric_spaces.h:500
IComponentSpace< TGridFunction > base_type
Definition: metric_spaces.h:499
L2ComponentSpace(const char *fctNames, int order, double weight, const char *ssNames=0)
Definition: metric_spaces.h:510
L2ComponentSpace(const char *fctNames, int order, ConstSmartPtr< weight_type > spWeight, const char *ssNames=0)
Definition: metric_spaces.h:513
L2ComponentSpace(const char *fctNames, int order)
Definition: metric_spaces.h:507
~L2ComponentSpace()
DTOR.
Definition: metric_spaces.h:517
double distance2(TGridFunction &uFine, TGridFunction &uCoarse)
distance (for grid functions)
Definition: metric_spaces.h:532
IObjectWithWeights< weight_type > weighted_obj_type
Definition: metric_spaces.h:501
L2ComponentSpace(const char *fctNames)
CTOR.
Definition: metric_spaces.h:504
double norm2(TGridFunction &uFine)
norm (for grid functions)
Definition: metric_spaces.h:528
Definition: metric_spaces.h:545
L2QuotientSpace(const char *fctNames, int order)
Definition: metric_spaces.h:555
~L2QuotientSpace()
DTOR.
Definition: metric_spaces.h:565
IComponentSpace< TGridFunction > base_type
Definition: metric_spaces.h:547
L2QuotientSpace(const char *fctNames)
CTOR.
Definition: metric_spaces.h:552
L2QuotientSpace(const char *fctNames, int order, ConstSmartPtr< weight_type > spWeight, const char *ssNames=0)
Definition: metric_spaces.h:561
double distance2(TGridFunction &uFine, TGridFunction &uCoarse)
distance (for grid functions)
Definition: metric_spaces.h:593
double norm2(TGridFunction &u)
norm (for grid functions)
Definition: metric_spaces.h:576
IObjectWithWeights< weight_type > weighted_obj_type
Definition: metric_spaces.h:549
L2DistIntegrand< TGridFunction >::weight_type weight_type
Definition: metric_spaces.h:548
L2QuotientSpace(const char *fctNames, int order, double weight, const char *ssNames=0)
Definition: metric_spaces.h:558
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
size_t num_children(TElem *elem) const
returns the number of children of the given child-type
Definition: multi_grid.h:225
Group of subsets.
Definition: subset_group.h:51
size_t size() const
number of subsets in this group
Definition: subset_group.h:122
void add_all()
select all subsets of underlying subset handler
Definition: subset_group.cpp:133
void add(int si)
adds a subset by number to this group
Definition: subset_group.cpp:64
Const iterator to traverse the surface of a multi-grid hierarchy.
Definition: surface_view.h:237
Represents the surface view of a multi-grid hierarchy.
Definition: surface_view.h:61
bool is_contained(TGeomObj *obj, const GridLevel &gl, SurfaceState validStates=ALL) const
returns if the element is contained in the surface view
Definition: surface_view_impl.hpp:442
@ SURFACE_RIM
Definition: surface_view.h:94
@ ALL
Definition: surface_view.h:102
@ SHADOW_RIM_COPY
Definition: surface_view.h:95
Wrapper class for time dependence.
Definition: metric_spaces.h:453
void update_time_data(number tScale)
characteristic time
Definition: metric_spaces.h:474
virtual ~TimeDependentSpace()
DTOR.
Definition: metric_spaces.h:463
double scaling() const
scaling (OVERRIDE)
Definition: metric_spaces.h:470
IGridFunctionSpace< TGridFunction > base_type
Definition: metric_spaces.h:455
SmartPtr< comp_space_type > m_spSpatialSpace
Definition: metric_spaces.h:487
TimeDependentSpace(SmartPtr< comp_space_type > spSpace, number tScale)
time dependent CTOR
Definition: metric_spaces.h:459
IComponentSpace< TGridFunction > comp_space_type
Definition: metric_spaces.h:456
number m_tScale
Definition: metric_spaces.h:488
std::string config_string() const
print config string
Definition: metric_spaces.h:478
Type based UserData.
Definition: user_data.h:143
For arbitrary UserData (of type TData), this class defines the integrand .
Definition: integrate.h:557
#define PCL_RO_SUM
Definition: pcl_methods.h:63
int NumProcs()
returns the number of processes
Definition: pcl_base.cpp:91
static const int dim
base_type::TBaseElem TBaseElem
const NullSmartPtr SPNULL
The equivalent to NULL for smart pointers.
Definition: smart_pointer.h:90
#define UG_ASSERT(expr, msg)
Definition: assert.h:70
#define UG_CATCH_THROW(msg)
Definition: error.h:64
#define UG_LOG(msg)
Definition: log.h:367
#define UG_COND_THROW(cond, msg)
UG_COND_THROW(cond, msg) : performs a UG_THROW(msg) if cond == true.
Definition: error.h:61
double number
Definition: types.h:124
the ug namespace
number L2Distance2(TGridFunction &spGridFct1, const char *cmp1, TGridFunction &spGridFct2, const char *cmp2, int quadOrder, const char *subsets, ConstSmartPtr< typename L2Integrand< TGridFunction >::weight_type > spWeight, number avgDist12=0.0)
computes the squared l2 distance between two functions
Definition: integrate.h:1911
number & DoFRef(TMatrix &mat, const DoFIndex &iInd, const DoFIndex &jInd)
Definition: multi_index.h:276
number L2Norm2(TGridFunction &u, const char *cmp, int quadOrder, const char *subsets, ConstSmartPtr< typename L2Integrand< TGridFunction >::weight_type > spWeight)
Definition: integrate.h:1680
number IntegrateSubsets(IIntegrand< number, TGridFunction::dim > &spIntegrand, TGridFunction &spGridFct, const char *subsets, int quadOrder, std::string quadType=std::string())
Definition: integrate.h:393
void TokenizeTrimString(const string &str, vector< string > &vToken, const char delimiter)
Definition: string_util.cpp:83
@ 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
number L2Error(SmartPtr< UserData< number, TGridFunction::dim > > spExactSol, TGridFunction &gridFct, const char *cmp, number time, int quadOrder, const char *subsets)
computes the l2 error function on the whole domain or on some subsets
Definition: integrate.h:1185
number Integral(SmartPtr< UserData< number, TGridFunction::dim > > spData, TGridFunction &spGridFct, const char *subsets, number time, int quadOrder, std::string quadType)
Definition: integrate.h:856
SmartPtr< T, FreePolicy > make_sp(T *inst)
returns a SmartPtr for the passed raw pointer
Definition: smart_pointer.h:836
function ProblemDisc new(problemDesc, dom)