ug4
domain_disc.h
Go to the documentation of this file.
1 /*
2  * Copyright (c) 2010-2015: G-CSC, Goethe University Frankfurt
3  * Author: Andreas Vogel
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__SPATIAL_DISC__DOMAIN_DISC__
34 #define __H__UG__LIB_DISC__SPATIAL_DISC__DOMAIN_DISC__
35 
36 // other ug4 modules
37 #include "common/common.h"
39 
40 // library intern headers
41 #include "subset_assemble_util.h"
42 #include "domain_disc_interface.h"
47 #include "disc_item.h"
51 
52 namespace ug {
53 
56 
58 
76 template <typename TDomain, typename TAlgebra, typename TGlobAssembler>
78 : public IDomainDiscretization<TAlgebra>,
79  public IDomainMarker<TDomain>,
80  protected TGlobAssembler
81 {
83  typedef TGlobAssembler gass_type;
84 
85  public:
87  typedef TDomain domain_type;
88 
90  typedef TAlgebra algebra_type;
91 
93  typedef typename algebra_type::matrix_type matrix_type;
94 
96  typedef typename algebra_type::vector_type vector_type;
97 
100 
103 
105  static const int dim = TDomain::dim;
106 
107  public:
110  m_bErrorCalculated(false),
111  m_spApproxSpace(pApproxSpace), m_spAssTuner(new AssemblingTuner<TAlgebra>)
112  {};
113 
116 
118  // Time independent part
120 
123  virtual void assemble_jacobian(matrix_type& J, const vector_type& u, const GridLevel& gl)
124  {assemble_jacobian(J, u, dd(gl));}
125 
128  virtual void assemble_defect(vector_type& d, const vector_type& u, const GridLevel& gl)
129  {assemble_defect(d, u, dd(gl));}
130 
133  virtual void assemble_linear(matrix_type& mat, vector_type& rhs, const GridLevel& gl)
134  {assemble_linear(mat, rhs, dd(gl));}
135 
138  virtual void assemble_rhs(vector_type& rhs, const vector_type& u, const GridLevel& gl)
139  {assemble_rhs(rhs, u, dd(gl));}
140 
143  {assemble_rhs(rhs, rhs, dd);}
144  virtual void assemble_rhs(vector_type& rhs, const GridLevel& gl)
145  {assemble_rhs(rhs, dd(gl));}
146 
149  virtual void adjust_solution(vector_type& u, const GridLevel& gl)
150  {adjust_solution(u, dd(gl));}
151 
155  {assemble_jacobian(J, u, u.dof_distribution());}
156 
158  {assemble_defect(d, u, u.dof_distribution());}
159 
161  {assemble_linear(A, rhs, rhs.dof_distribution());}
162 
164  {assemble_rhs(rhs, u, u.dof_distribution());}
165 
167  {assemble_rhs(b, b.dof_distribution());}
168 
172 
174  // Time dependent part
178  virtual void prepare_timestep(ConstSmartPtr<VectorTimeSeries<vector_type> > vSol, number future_time, const GridLevel& gl)
179  {prepare_timestep(vSol, future_time, dd(gl));}
180 
181 
186  {prepare_timestep_elem(vSol, dd(gl));}
187 
189  virtual void assemble_jacobian(matrix_type& J,
191  const number s_a0,
195  const number s_a, const GridLevel& gl)
196  {assemble_jacobian(J, vSol, s_a, dd(gl));}
197 
199  virtual void assemble_defect(vector_type& d,
201  const std::vector<number>& vScaleMass,
202  const std::vector<number>& vScaleStiff,
204  virtual void assemble_defect(vector_type& d,
206  const std::vector<number>& vScaleMass,
207  const std::vector<number>& vScaleStiff,
208  const GridLevel& gl)
209  {assemble_defect(d, vSol, vScaleMass, vScaleStiff, dd(gl));}
210 
212  virtual void assemble_linear(matrix_type& A, vector_type& b,
214  const std::vector<number>& vScaleMass,
215  const std::vector<number>& vScaleStiff,
219  const std::vector<number>& vScaleMass,
220  const std::vector<number>& vScaleStiff,
221  const GridLevel& gl)
222  {assemble_linear(A, b, vSol, vScaleMass, vScaleStiff, dd(gl));}
223 
225  virtual void assemble_rhs( vector_type& b,
227  const std::vector<number>& vScaleMass,
228  const std::vector<number>& vScaleStiff,
230  virtual void assemble_rhs( vector_type& b,
232  const std::vector<number>& vScaleMass,
233  const std::vector<number>& vScaleStiff,
234  const GridLevel& gl)
235  {assemble_rhs(b, vSol, vScaleMass, vScaleStiff, dd(gl));}
236 
239  virtual void adjust_solution(vector_type& u, number time, const GridLevel& gl)
240  {adjust_solution(u, time, dd(gl));}
241 
245  {finish_timestep(vSol, dd(gl));}
246 
250  {finish_timestep_elem(vSol, dd(gl));}
251 
253  // Mass and Stiffness Matrix
255 
257  virtual void assemble_mass_matrix(matrix_type& M, const vector_type& u,
259  virtual void assemble_mass_matrix(matrix_type& M, const vector_type& u,
260  const GridLevel& gl)
261  {assemble_mass_matrix(M, u, dd(gl));}
262 
264  virtual void assemble_stiffness_matrix(matrix_type& A, const vector_type& u,
267  const GridLevel& gl)
268  {assemble_stiffness_matrix(A, u, dd(gl));}
269 
271  // Error estimator ///
272 public:
274  // stationary
275  virtual void calc_error(const vector_type& u, ConstSmartPtr<DoFDistribution> dd, error_vector_type* u_vtk = NULL);
276  virtual void calc_error(const vector_type& u, const GridLevel& gl, error_vector_type* u_vtk = NULL)
277  {calc_error(u, dd(gl));}
278 
280  {calc_error(u, u.dd(), NULL);}
282  {calc_error(u, u.dd(), u_vtk);}
283 
284  // instationary
287  const std::vector<number>& vScaleMass,
288  const std::vector<number>& vScaleStiff,
289  error_vector_type* u_vtk);
291  const std::vector<number>& vScaleMass,
292  const std::vector<number>& vScaleStiff,
293  const GridLevel& gl,
294  error_vector_type* u_vtk)
295  {
297  vScaleMass, vScaleStiff, u_vtk);
298  }
299 
300  virtual void mark_with_strategy(IRefiner& refiner, SmartPtr <IElementMarkingStrategy<TDomain> > strategy);
301 
304  virtual void invalidate_error();
305 
307  virtual bool is_error_valid();
308 
309  protected:
312 
313  // Error estimator //
315 
316  public:
321 
322  public:
324 
334  {
335  // check that not already registered
336  for(size_t i = 0; i < m_vDomainElemDisc.size(); ++i)
337  if(m_vDomainElemDisc[i] == elem)
338  return;
339 
340  // set approximation space
341  elem->set_approximation_space(m_spApproxSpace);
342 
343  // add it
344  m_vDomainElemDisc.push_back(elem);
345  }
346 
348 
355  {
356  // check that already registered
357  for (size_t i = 0; i < m_vDomainElemDisc.size(); i++)
358  {
359  if (m_vDomainElemDisc[i] == elem)
360  {
361  // remove constraint
362  m_vDomainElemDisc.erase(m_vDomainElemDisc.begin()+i);
363  return;
364  }
365  }
366 
367  UG_LOG("Tried to remove ElemDisc from DomainDisc"
368  ", but could not find it there.");
369  }
370 
372 
382  {
383  // check that not already registered
384  for(size_t i = 0; i < m_vDomainElemError.size(); ++i)
385  if(m_vDomainElemError[i] == elem) return;
386 
387  // set approximation space
388  elem->set_approximation_space(m_spApproxSpace);
389 
390  // add it
391  m_vDomainElemError.push_back(elem);
392  }
393 
395 
402  {
403  // check that 'elem' already registered
404  for (size_t i = 0; i < m_vDomainElemError.size(); i++)
405  {
406  if (m_vDomainElemError[i] == elem)
407  {
408  // remove constraint
409  m_vDomainElemError.erase(m_vDomainElemError.begin()+i);
410  return;
411  }
412  }
413 
414  UG_LOG("Tried to remove ElemError from DomainDisc"
415  ", but could not find it there.");
416  }
417 
418 
420 
427  {
428  // check that not already registered
429  for(size_t i = 0; i < m_vConstraint.size(); ++i)
430  if(m_vConstraint[i] == pp)
431  return;
432 
433  // set approximation space
434  pp->set_approximation_space(m_spApproxSpace);
435 
436  // add constraint
437  m_vConstraint.push_back(pp);
438  }
439 
441 
448  {
449  // check that already registered
450  for (size_t i = 0; i < m_vConstraint.size(); i++)
451  {
452  if (m_vConstraint[i] == pp)
453  {
454  // remove constraint
455  m_vConstraint.erase(m_vConstraint.begin()+i);
456  return;
457  }
458  }
459 
460  UG_LOG("Tried to remove DomainConstraint from DomainDisc"
461  ", but could not find it there.");
462  }
463 
465 
472  {
473  // add elem discs
474  for(size_t i = 0; i < di->num_elem_disc(); ++i)
475  add(di->elem_disc(i));
476 
477  // add constraints
478  for(size_t i = 0; i < di->num_constraint(); ++i)
479  add(di->constraint(i));
480  }
481 
483  virtual size_t num_constraints() const {return m_vConstraint.size();}
484 
486  virtual SmartPtr<IConstraint<TAlgebra> > constraint(size_t i) {return m_vConstraint[i];}
487 
490 
491  protected:
493  void update_elem_discs();
494  void update_elem_errors();
495  void update_constraints();
496  void update_disc_items();
497  void update_error_items();
498 
499  protected:
501  ConstSmartPtr<DoFDistribution> dd(const GridLevel& gl) const{return m_spApproxSpace->dof_distribution(gl);}
502 
503  protected:
505  std::vector<SmartPtr<IElemDisc<TDomain> > > m_vDomainElemDisc;
506 
508  std::vector<IElemDisc<TDomain>*> m_vElemDisc;
509 
511  std::vector<SmartPtr<IElemError<TDomain> > > m_vDomainElemError;
512 
514  std::vector<IElemError<TDomain>*> m_vElemError;
515 
517  std::vector<SmartPtr<IDomainConstraint<TDomain, TAlgebra> > > m_vConstraint;
518 
521 
524 
525  private:
526  //---- Auxiliary function templates for the assembling ----//
527  // These functions call the corresponding functions from the global assembler for a composed list of elements:
528  //-- for stationary problems --//
529  template <typename TElem>
530  void AssembleMassMatrix( const std::vector<IElemDisc<domain_type>*>& vElemDisc,
532  int si, bool bNonRegularGrid,
533  matrix_type& M,
534  const vector_type& u);
535  template <typename TElem>
536  void AssembleStiffnessMatrix( const std::vector<IElemDisc<domain_type>*>& vElemDisc,
538  int si, bool bNonRegularGrid,
539  matrix_type& A,
540  const vector_type& u);
541  template <typename TElem>
542  void AssembleJacobian( const std::vector<IElemDisc<domain_type>*>& vElemDisc,
544  int si, bool bNonRegularGrid,
545  matrix_type& J,
546  const vector_type& u);
547  template <typename TElem>
548  void AssembleDefect( const std::vector<IElemDisc<domain_type>*>& vElemDisc,
550  int si, bool bNonRegularGrid,
551  vector_type& d,
552  const vector_type& u);
553  template <typename TElem>
554  void AssembleLinear( const std::vector<IElemDisc<domain_type>*>& vElemDisc,
556  int si, bool bNonRegularGrid,
557  matrix_type& A,
558  vector_type& rhs);
559  template <typename TElem>
560  void AssembleRhs( const std::vector<IElemDisc<domain_type>*>& vElemDisc,
562  int si, bool bNonRegularGrid,
563  vector_type& rhs,
564  const vector_type& u);
565  template <typename TElem>
566  void AssembleErrorEstimator( const std::vector<IElemError<domain_type>*>& vElemDisc,
568  int si, bool bNonRegularGrid,
569  const vector_type& u);
570  //-- for time-dependent problems --//
571  template <typename TElem>
572  void PrepareTimestepElem( const std::vector<IElemDisc<domain_type>*>& vElemDisc,
574  int si, bool bNonRegularGrid,
576  template <typename TElem>
577  void AssembleJacobian( const std::vector<IElemDisc<domain_type>*>& vElemDisc,
579  int si, bool bNonRegularGrid,
580  matrix_type& J,
582  number s_a0);
583  template <typename TElem>
584  void AssembleDefect( const std::vector<IElemDisc<domain_type>*>& vElemDisc,
586  int si, bool bNonRegularGrid,
587  vector_type& d,
589  const std::vector<number>& vScaleMass,
590  const std::vector<number>& vScaleStiff);
591  template <typename TElem>
592  void AssembleLinear( const std::vector<IElemDisc<domain_type>*>& vElemDisc,
594  int si, bool bNonRegularGrid,
595  matrix_type& A,
596  vector_type& rhs,
598  const std::vector<number>& vScaleMass,
599  const std::vector<number>& vScaleStiff);
600  template <typename TElem>
601  void AssembleRhs( const std::vector<IElemDisc<domain_type>*>& vElemDisc,
603  int si, bool bNonRegularGrid,
604  vector_type& rhs,
606  const std::vector<number>& vScaleMass,
607  const std::vector<number>& vScaleStiff);
608  template <typename TElem>
609  void FinishTimestepElem( const std::vector<IElemDisc<domain_type>*>& vElemDisc,
611  int si, bool bNonRegularGrid,
613  template <typename TElem>
614  void AssembleErrorEstimator( const std::vector<IElemError<domain_type>*>& vElemDisc,
616  int si, bool bNonRegularGrid,
617  const std::vector<number>& vScaleMass,
618  const std::vector<number>& vScaleStiff,
620 };
621 
623 
635 template <typename TDomain, typename TAlgebra>
637 : public DomainDiscretizationBase<TDomain, TAlgebra, StdGlobAssembler<TDomain, TAlgebra> >
638 {
641 
642  public:
644  typedef TDomain domain_type;
645 
647  typedef TAlgebra algebra_type;
648 
650  typedef typename algebra_type::matrix_type matrix_type;
651 
653  typedef typename algebra_type::vector_type vector_type;
654 
657 
659  static const int dim = TDomain::dim;
660 
661  public:
665  {};
666 
668  virtual ~DomainDiscretization() {};
669 };
670 
672 
673 } // end namespace ug
674 
675 // include documentation
676 #include "domain_disc_impl.h"
677 
678 #endif /* __H__UG__LIB_DISC__SPATIAL_DISC__DOMAIN_DISC__ */
Definition: smart_pointer.h:296
Definition: smart_pointer.h:108
base class for approximation spaces without type of algebra or dof distribution
Definition: approximation_space.h:279
The AssemblingTuner class combines tools to adapt the assembling routine.
Definition: ass_tuner.h:90
generic domain discretization implementing the interface
Definition: domain_disc.h:81
void update_constraints()
Definition: domain_disc_impl.h:112
void add_elem_error_indicator(SmartPtr< IElemError< TDomain > > elem)
adds an element error indicator to the assembling process
Definition: domain_disc.h:381
virtual void calc_error(const vector_type &u, ConstSmartPtr< DoFDistribution > dd, error_vector_type *u_vtk=NULL)
Definition: domain_disc_impl.h:2134
void remove(SmartPtr< IElemDisc< TDomain > > elem)
removes a element discretization from the assembling process
Definition: domain_disc.h:354
SmartPtr< AssemblingTuner< TAlgebra > > m_spAssTuner
this object provides tools to adapt the assemble routine
Definition: domain_disc.h:523
std::vector< SmartPtr< IElemDisc< TDomain > > > m_vDomainElemDisc
vector holding all registered elem discs
Definition: domain_disc.h:505
void add(SmartPtr< IDiscretizationItem< TDomain, TAlgebra > > di)
adds a disc item to the assembling process
Definition: domain_disc.h:471
void update_elem_discs()
set the approximation space in the elem discs and extract IElemDiscs
Definition: domain_disc_impl.h:68
virtual bool is_error_valid()
returns whether current error values are valid
Definition: domain_disc_impl.h:2705
virtual void calc_error(const vector_type &u, const GridLevel &gl, error_vector_type *u_vtk=NULL)
Definition: domain_disc.h:276
virtual size_t num_constraints() const
returns number of registered constraints
Definition: domain_disc.h:483
TDomain domain_type
Type of Domain.
Definition: domain_disc.h:87
void assemble_defect(vector_type &d, GridFunction< TDomain, TAlgebra > &u)
Definition: domain_disc.h:157
bool m_bErrorCalculated
Definition: domain_disc.h:311
virtual void invalidate_error()
Definition: domain_disc_impl.h:2690
void AssembleJacobian(const std::vector< IElemDisc< domain_type > * > &vElemDisc, ConstSmartPtr< DoFDistribution > dd, int si, bool bNonRegularGrid, matrix_type &J, const vector_type &u)
Definition: domain_disc_impl.h:615
void AssembleErrorEstimator(const std::vector< IElemError< domain_type > * > &vElemDisc, ConstSmartPtr< DoFDistribution > dd, int si, bool bNonRegularGrid, const vector_type &u)
Definition: domain_disc_impl.h:2397
std::vector< IElemDisc< TDomain > * > m_vElemDisc
vector holding all registered elem discs
Definition: domain_disc.h:508
ApproximationSpace< TDomain > approx_space_type
Type of approximation space.
Definition: domain_disc.h:102
void remove_elem_error_indicator(SmartPtr< IElemError< TDomain > > elem)
removes a element discretization from the assembling process
Definition: domain_disc.h:401
virtual void finish_timestep(ConstSmartPtr< VectorTimeSeries< vector_type > > vSol, ConstSmartPtr< DoFDistribution > dd)
finishes time step
Definition: domain_disc_impl.h:2717
TAlgebra algebra_type
Type of algebra.
Definition: domain_disc.h:90
algebra_type::matrix_type matrix_type
Type of algebra matrix.
Definition: domain_disc.h:93
virtual void assemble_rhs(vector_type &rhs, const GridLevel &gl)
assembles the rhs
Definition: domain_disc.h:144
virtual void prepare_timestep_elem(ConstSmartPtr< VectorTimeSeries< vector_type > > vSol, ConstSmartPtr< DoFDistribution > dd)
prepares time step element-wise
Definition: domain_disc_impl.h:1241
static const int dim
world dimension
Definition: domain_disc.h:105
void adjust_solution(GridFunction< TDomain, TAlgebra > &u)
Definition: domain_disc.h:169
virtual void adjust_solution(vector_type &u, number time, const GridLevel &gl)
sets dirichlet values in solution vector
Definition: domain_disc.h:239
virtual void calc_error(const GridFunction< TDomain, TAlgebra > &u)
Definition: domain_disc.h:279
virtual void prepare_timestep(ConstSmartPtr< VectorTimeSeries< vector_type > > vSol, number future_time, ConstSmartPtr< DoFDistribution > dd)
prepares time step
Definition: domain_disc_impl.h:1204
virtual void assemble_rhs(vector_type &rhs, const vector_type &u, const GridLevel &gl)
assembles the rhs
Definition: domain_disc.h:138
virtual void assemble_jacobian(matrix_type &J, ConstSmartPtr< VectorTimeSeries< vector_type > > vSol, const number s_a, const GridLevel &gl)
assembles Jacobian (or Approximation of Jacobian)
Definition: domain_disc.h:193
void update_disc_items()
Definition: domain_disc_impl.h:127
void AssembleMassMatrix(const std::vector< IElemDisc< domain_type > * > &vElemDisc, ConstSmartPtr< DoFDistribution > dd, int si, bool bNonRegularGrid, matrix_type &M, const vector_type &u)
Definition: domain_disc_impl.h:271
virtual void assemble_defect(vector_type &d, const vector_type &u, ConstSmartPtr< DoFDistribution > dd)
assembles Defect
Definition: domain_disc_impl.h:648
virtual void mark_with_strategy(IRefiner &refiner, SmartPtr< IElementMarkingStrategy< TDomain > > strategy)
Definition: domain_disc_impl.h:2671
void assemble_linear(matrix_type &A, GridFunction< TDomain, TAlgebra > &rhs)
Definition: domain_disc.h:160
void remove(SmartPtr< IDomainConstraint< TDomain, TAlgebra > > pp)
removes a constraint from the assembling process
Definition: domain_disc.h:447
virtual void prepare_timestep(ConstSmartPtr< VectorTimeSeries< vector_type > > vSol, number future_time, const GridLevel &gl)
Definition: domain_disc.h:178
void AssembleRhs(const std::vector< IElemDisc< domain_type > * > &vElemDisc, ConstSmartPtr< DoFDistribution > dd, int si, bool bNonRegularGrid, vector_type &rhs, const vector_type &u)
Definition: domain_disc_impl.h:1124
virtual void assemble_mass_matrix(matrix_type &M, const vector_type &u, const GridLevel &gl)
assembles the mass matrix
Definition: domain_disc.h:259
virtual void adjust_solution(vector_type &u, const GridLevel &gl)
sets dirichlet values in solution vector
Definition: domain_disc.h:149
virtual void assemble_mass_matrix(matrix_type &M, const vector_type &u, ConstSmartPtr< DoFDistribution > dd)
assembles the mass matrix
Definition: domain_disc_impl.h:145
DomainDiscretizationBase(SmartPtr< approx_space_type > pApproxSpace)
default Constructor
Definition: domain_disc.h:109
virtual void assemble_linear(matrix_type &mat, vector_type &rhs, const GridLevel &gl)
Assembles Matrix and Right-Hand-Side for a linear problem.
Definition: domain_disc.h:133
virtual void assemble_stiffness_matrix(matrix_type &A, const vector_type &u, const GridLevel &gl)
assembles the stiffness matrix
Definition: domain_disc.h:266
virtual void assemble_defect(vector_type &d, ConstSmartPtr< VectorTimeSeries< vector_type > > vSol, const std::vector< number > &vScaleMass, const std::vector< number > &vScaleStiff, const GridLevel &gl)
assembles Defect
Definition: domain_disc.h:204
IMultigridElementIndicators< TDomain > m_mgElemErrors
Definition: domain_disc.h:310
void assemble_rhs(vector_type &rhs, GridFunction< TDomain, TAlgebra > &u)
Definition: domain_disc.h:163
void AssembleDefect(const std::vector< IElemDisc< domain_type > * > &vElemDisc, ConstSmartPtr< DoFDistribution > dd, int si, bool bNonRegularGrid, vector_type &d, const vector_type &u)
Definition: domain_disc_impl.h:804
void add(SmartPtr< IElemDisc< TDomain > > elem)
adds an element discretization to the assembling process
Definition: domain_disc.h:333
std::vector< SmartPtr< IElemError< TDomain > > > m_vDomainElemError
vector holding all registered elem discs
Definition: domain_disc.h:511
void add(SmartPtr< IDomainConstraint< TDomain, TAlgebra > > pp)
adds a constraint to the assembling process
Definition: domain_disc.h:426
CPUAlgebra::vector_type error_vector_type
Type of error vector.
Definition: domain_disc.h:99
virtual void assemble_linear(matrix_type &A, vector_type &b, ConstSmartPtr< VectorTimeSeries< vector_type > > vSol, const std::vector< number > &vScaleMass, const std::vector< number > &vScaleStiff, const GridLevel &gl)
Assembles matrix_type and Right-Hand-Side for a linear problem.
Definition: domain_disc.h:217
ConstSmartPtr< DoFDistribution > dd(const GridLevel &gl) const
returns the level dof distribution
Definition: domain_disc.h:501
virtual ~DomainDiscretizationBase()
virtual destructor
Definition: domain_disc.h:115
void AssembleLinear(const std::vector< IElemDisc< domain_type > * > &vElemDisc, ConstSmartPtr< DoFDistribution > dd, int si, bool bNonRegularGrid, matrix_type &A, vector_type &rhs)
Definition: domain_disc_impl.h:965
std::vector< IElemError< TDomain > * > m_vElemError
vector holding all registered elem discs
Definition: domain_disc.h:514
virtual void assemble_stiffness_matrix(matrix_type &A, const vector_type &u, ConstSmartPtr< DoFDistribution > dd)
assembles the stiffness matrix
Definition: domain_disc_impl.h:305
virtual void assemble_rhs(vector_type &b, ConstSmartPtr< VectorTimeSeries< vector_type > > vSol, const std::vector< number > &vScaleMass, const std::vector< number > &vScaleStiff, const GridLevel &gl)
Assembles Right-Hand-Side for a linear problem.
Definition: domain_disc.h:230
virtual void assemble_rhs(vector_type &rhs, const vector_type &u, ConstSmartPtr< DoFDistribution > dd)
assembles rhs
Definition: domain_disc_impl.h:998
virtual void assemble_linear(matrix_type &A, vector_type &b, ConstSmartPtr< DoFDistribution > dd)
Assembles Matrix and Right-Hand-Side for a linear problem.
Definition: domain_disc_impl.h:837
TGlobAssembler gass_type
Type of the global assembler.
Definition: domain_disc.h:83
virtual ConstSmartPtr< AssemblingTuner< TAlgebra > > ass_tuner() const
Definition: domain_disc.h:319
std::vector< SmartPtr< IDomainConstraint< TDomain, TAlgebra > > > m_vConstraint
vector holding all registered constraints
Definition: domain_disc.h:517
void FinishTimestepElem(const std::vector< IElemDisc< domain_type > * > &vElemDisc, ConstSmartPtr< DoFDistribution > dd, int si, bool bNonRegularGrid, ConstSmartPtr< VectorTimeSeries< vector_type > > vSol)
Definition: domain_disc_impl.h:2855
void assemble_jacobian(matrix_type &J, GridFunction< TDomain, TAlgebra > &u)
Definition: domain_disc.h:154
virtual void finish_timestep_elem(ConstSmartPtr< VectorTimeSeries< vector_type > > vSol, const GridLevel &gl)
finishes timestep
Definition: domain_disc.h:249
algebra_type::vector_type vector_type
Type of algebra vector.
Definition: domain_disc.h:96
SmartPtr< approx_space_type > approximation_space()
Definition: domain_disc.h:488
virtual void finish_timestep_elem(ConstSmartPtr< VectorTimeSeries< vector_type > > vSol, ConstSmartPtr< DoFDistribution > dd)
finishes timestep
Definition: domain_disc_impl.h:2753
void AssembleStiffnessMatrix(const std::vector< IElemDisc< domain_type > * > &vElemDisc, ConstSmartPtr< DoFDistribution > dd, int si, bool bNonRegularGrid, matrix_type &A, const vector_type &u)
Definition: domain_disc_impl.h:431
ConstSmartPtr< approx_space_type > approximation_space() const
Definition: domain_disc.h:489
void assemble_rhs(GridFunction< TDomain, TAlgebra > &b)
Definition: domain_disc.h:166
virtual void assemble_jacobian(matrix_type &J, const vector_type &u, const GridLevel &gl)
assembles Jacobian (or Approximation of Jacobian)
Definition: domain_disc.h:123
SmartPtr< approx_space_type > m_spApproxSpace
current approximation space
Definition: domain_disc.h:520
virtual void calc_error(ConstSmartPtr< VectorTimeSeries< vector_type > > vSol, const std::vector< number > &vScaleMass, const std::vector< number > &vScaleStiff, const GridLevel &gl, error_vector_type *u_vtk)
Transient version.
Definition: domain_disc.h:290
virtual void calc_error(const GridFunction< TDomain, TAlgebra > &u, error_vector_type *u_vtk)
Definition: domain_disc.h:281
virtual void prepare_timestep_elem(ConstSmartPtr< VectorTimeSeries< vector_type > > vSol, const GridLevel &gl)
Definition: domain_disc.h:185
void PrepareTimestepElem(const std::vector< IElemDisc< domain_type > * > &vElemDisc, ConstSmartPtr< DoFDistribution > dd, int si, bool bNonRegularGrid, ConstSmartPtr< VectorTimeSeries< vector_type > > vSol)
Definition: domain_disc_impl.h:1344
virtual SmartPtr< AssemblingTuner< TAlgebra > > ass_tuner()
Definition: domain_disc.h:318
virtual SmartPtr< IConstraint< TAlgebra > > constraint(size_t i)
returns the i'th constraint
Definition: domain_disc.h:486
virtual void assemble_jacobian(matrix_type &J, const vector_type &u, ConstSmartPtr< DoFDistribution > dd)
assembles Jacobian (or Approximation of Jacobian)
Definition: domain_disc_impl.h:471
void update_error_items()
Definition: domain_disc_impl.h:134
void update_elem_errors()
Definition: domain_disc_impl.h:89
virtual void adjust_solution(vector_type &u, ConstSmartPtr< DoFDistribution > dd)
sets dirichlet values in solution vector
Definition: domain_disc_impl.h:1156
virtual void assemble_rhs(vector_type &rhs, ConstSmartPtr< DoFDistribution > dd)
assembles rhs
Definition: domain_disc.h:142
virtual void assemble_defect(vector_type &d, const vector_type &u, const GridLevel &gl)
assembles Defect
Definition: domain_disc.h:128
virtual void finish_timestep(ConstSmartPtr< VectorTimeSeries< vector_type > > vSol, const GridLevel &gl)
Definition: domain_disc.h:244
domain discretization implementing the interface
Definition: domain_disc.h:638
virtual ~DomainDiscretization()
virtual destructor
Definition: domain_disc.h:668
TDomain domain_type
Type of Domain.
Definition: domain_disc.h:644
TAlgebra algebra_type
Type of algebra.
Definition: domain_disc.h:647
static const int dim
world dimension
Definition: domain_disc.h:659
algebra_type::vector_type vector_type
Type of algebra vector.
Definition: domain_disc.h:653
StdGlobAssembler< TDomain, TAlgebra > gass_type
Type of the global assembler.
Definition: domain_disc.h:640
algebra_type::matrix_type matrix_type
Type of algebra matrix.
Definition: domain_disc.h:650
DomainDiscretization(SmartPtr< approx_space_type > pApproxSpace)
default Constructor
Definition: domain_disc.h:663
ApproximationSpace< TDomain > approx_space_type
Type of approximation space.
Definition: domain_disc.h:656
represents numerical solutions on a grid using an algebraic vector
Definition: grid_function.h:121
SmartPtr< DoFDistribution > dof_distribution()
Definition: grid_function.h:236
SmartPtr< DoFDistribution > dd()
Definition: grid_function.h:237
Definition: grid_level.h:42
TAlgebra algebra_type
Algebra type.
Definition: assemble_interface.h:113
TAlgebra::vector_type vector_type
Type of algebra vector.
Definition: assemble_interface.h:119
TAlgebra::matrix_type matrix_type
Type of algebra matrix.
Definition: assemble_interface.h:116
Definition: disc_item.h:43
Definition: constraint_interface.h:162
Interface for domain discretization.
Definition: domain_disc_interface.h:136
Interface for an object that can mark elements based on a strategy.
Definition: domain_disc_interface.h:107
Definition: elem_disc_interface.h:800
Definition: elem_disc_interface.h:760
Abstract base class for element marking (in adaptive refinement)
Definition: error_elem_marking_strategy.h:125
This class encapsulates the multi-grid attachments for error estimation.
Definition: error_elem_marking_strategy.h:54
The refiner interface allows to mark elements for refinement and to call refine.
Definition: refiner_interface.h:67
Wrapper for sequential matrices to handle them in parallel.
Definition: parallel_matrix.h:65
Global assembler based on the straightforward application of the local discretizations.
Definition: elem_disc_assemble_util.h:78
time series of solutions and corresponding time point
Definition: solution_time_series.h:59
static const int dim
#define UG_LOG(msg)
Definition: log.h:367
double number
Definition: types.h:124
CPUAlgebra::vector_type vector_type
the ug namespace