ug4
domain_disc_interface.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_INTERFACE__
34 #define __H__UG__LIB_DISC__SPATIAL_DISC__DOMAIN_DISC_INTERFACE__
35 
42 
43 namespace ug {
44 
47 
49 template <typename TAlgebra>
51 {
52 public:
54  typedef typename TAlgebra::vector_type vector_type;
55 
58 
59  // (virtual) destructor
60  virtual ~IDomainErrorIndicator() {};
61 
68  virtual void calc_error
69  ( const vector_type& u,
70  const GridLevel& gl,
71  error_vector_type* u_vtk = NULL
72  ) = 0;
73  virtual void calc_error
74  ( const vector_type& u,
76  error_vector_type* u_vtk = NULL
77  ) = 0;
78 
80  virtual void calc_error
83  const std::vector<number>& vScaleMass,
84  const std::vector<number>& vScaleStiff,
85  error_vector_type* u_vtk
86  ) = 0;
87 
89  virtual void calc_error
91  const std::vector<number>& vScaleMass,
92  const std::vector<number>& vScaleStiff,
93  const GridLevel& gl,
94  error_vector_type* u_vtk
95  ) = 0;
96 
98  virtual void invalidate_error() = 0;
99 
101  virtual bool is_error_valid() = 0;
102 };
103 
105 template <typename TDomain>
107 {
108 public:
111 
112  // (virtual) destructor
113  virtual ~IDomainMarker() {};
114 
115  virtual void mark_with_strategy
116  ( IRefiner& refiner,
118  ) = 0 ;
119 };
120 
122 
134 template <typename TAlgebra>
135 class IDomainDiscretization : public IAssemble<TAlgebra>, public IDomainErrorIndicator<TAlgebra>
136 {
138 
139  public:
141  typedef TAlgebra algebra_type;
142 
144  typedef typename algebra_type::matrix_type matrix_type;
145 
147  typedef typename algebra_type::vector_type vector_type;
149 
150  public:
152 
159  virtual void assemble_jacobian(matrix_type& J, const vector_type& u, const GridLevel& gl) = 0;
161 
163 
170  virtual void assemble_defect(vector_type& d, const vector_type& u, const GridLevel& gl) = 0;
172 
174 
181  virtual void assemble_linear(matrix_type& A, vector_type& b, const GridLevel& gl) = 0;
183 
185  virtual void assemble_rhs(vector_type& rhs, const vector_type& u, const GridLevel& gl) = 0;
187 
189  virtual void assemble_rhs(vector_type& b, const GridLevel& gl) = 0;
191 
193 
200  virtual void adjust_solution(vector_type& u, const GridLevel& gl) = 0;
202 
204  virtual void assemble_mass_matrix(matrix_type& M, const vector_type& u, const GridLevel& gl) = 0;
206 
208  virtual void assemble_stiffness_matrix(matrix_type& A, const vector_type& u, const GridLevel& gl) = 0;
210 
211 
212  public:
214 
222  virtual void prepare_timestep(ConstSmartPtr<VectorTimeSeries<vector_type> > vSol, number future_time, const GridLevel& gl) = 0;
224  {prepare_timestep(vSol, future_time, GridLevel());}
225 
227 
235 
238  {prepare_timestep_elem(vSol, GridLevel());}
239 
241 
251  const number s_a, const GridLevel& gl) = 0;
254  const number s_a0,
256 
260  const number s_a)
261  {assemble_jacobian(J, vSol, s_a, GridLevel());}
262 
263 
265 
274  virtual void assemble_defect(vector_type& d,
276  const std::vector<number>& vScaleMass,
277  const std::vector<number>& vScaleStiff,
278  const GridLevel& gl) = 0;
279  virtual void assemble_defect(vector_type& d,
281  const std::vector<number>& vScaleMass,
282  const std::vector<number>& vScaleStiff,
284 
288  const std::vector<number>& vScaleMass,
289  const std::vector<number>& vScaleStiff)
290  {assemble_defect(d, vSol, vScaleMass, vScaleStiff, GridLevel());}
291 
293 
305  const std::vector<number>& vScaleMass,
306  const std::vector<number>& vScaleStiff,
307  const GridLevel& gl) = 0;
310  const std::vector<number>& vScaleMass,
311  const std::vector<number>& vScaleStiff,
313 
317  const std::vector<number>& vScaleMass,
318  const std::vector<number>& vScaleStiff)
319  {assemble_linear(A,b,vSol,vScaleMass,vScaleStiff, GridLevel());}
320 
322 
331  virtual void assemble_rhs( vector_type& b,
333  const std::vector<number>& vScaleMass,
334  const std::vector<number>& vScaleStiff,
335  const GridLevel& gl) = 0;
336  virtual void assemble_rhs( vector_type& b,
338  const std::vector<number>& vScaleMass,
339  const std::vector<number>& vScaleStiff,
341 
345  const std::vector<number>& vScaleMass,
346  const std::vector<number>& vScaleStiff)
347  {assemble_rhs(b,vSol,vScaleMass,vScaleStiff, GridLevel());}
348 
349 
351 
358  virtual void adjust_solution(vector_type& u, number time, const GridLevel& gl) = 0;
360 
363  {adjust_solution(u,time, GridLevel());}
364 
366 
376  {finish_timestep(vSol, GridLevel());}
377 
379 
387 
390  {finish_timestep_elem(vSol, GridLevel());}
391 
392 
393 
394 
396  virtual size_t num_constraints() const = 0;
397 
400 };
401 
403 
404 }; // namespace ug
405 
406 #endif /* __H__UG__LIB_DISC__SPATIAL_DISC__DOMAIN_DISC_INTERFACE__ */
Definition: smart_pointer.h:296
Definition: smart_pointer.h:108
Definition: grid_level.h:42
Interface providing Jacobian and Defect of a discretization.
Definition: assemble_interface.h:110
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
Interface for domain discretization.
Definition: domain_disc_interface.h:136
virtual SmartPtr< IConstraint< TAlgebra > > constraint(size_t i)=0
returns the i'th post process
virtual void assemble_linear(matrix_type &A, vector_type &b, const GridLevel &gl)=0
Assembles Matrix and Right-Hand-Side for a linear problem.
virtual void assemble_defect(vector_type &d, const vector_type &u, const GridLevel &gl)=0
assembles Defect
virtual void assemble_jacobian(matrix_type &J, ConstSmartPtr< VectorTimeSeries< vector_type > > vSol, const number s_a, const GridLevel &gl)=0
assembles Jacobian (or Approximation of Jacobian)
virtual void assemble_jacobian(matrix_type &J, const vector_type &u, ConstSmartPtr< DoFDistribution > dd)=0
virtual void prepare_timestep_elem(ConstSmartPtr< VectorTimeSeries< vector_type > > vSol, const GridLevel &gl)=0
virtual void prepare_timestep(ConstSmartPtr< VectorTimeSeries< vector_type > > vSol, number future_time)
Definition: domain_disc_interface.h:223
virtual void assemble_defect(vector_type &d, const vector_type &u, ConstSmartPtr< DoFDistribution > dd)=0
virtual void prepare_timestep(ConstSmartPtr< VectorTimeSeries< vector_type > > vSol, number future_time, ConstSmartPtr< DoFDistribution > dd)=0
prepares time step
virtual void assemble_mass_matrix(matrix_type &M, const vector_type &u, const GridLevel &gl)=0
assembles the mass matrix
algebra_type::vector_type vector_type
Type of algebra vector.
Definition: domain_disc_interface.h:147
virtual void assemble_rhs(vector_type &b, ConstSmartPtr< VectorTimeSeries< vector_type > > vSol, const std::vector< number > &vScaleMass, const std::vector< number > &vScaleStiff, ConstSmartPtr< DoFDistribution > dd)=0
virtual void finish_timestep_elem(ConstSmartPtr< VectorTimeSeries< vector_type > > vSol, ConstSmartPtr< DoFDistribution > dd)=0
virtual void assemble_mass_matrix(matrix_type &M, const vector_type &u, ConstSmartPtr< DoFDistribution > dd)=0
virtual void finish_timestep_elem(ConstSmartPtr< VectorTimeSeries< vector_type > > vSol, const GridLevel &gl)=0
finishes timestep
virtual void assemble_rhs(vector_type &b, ConstSmartPtr< DoFDistribution > dd)=0
void prepare_timestep_elem(ConstSmartPtr< VectorTimeSeries< vector_type > > vSol)
prepares timestep on surface level
Definition: domain_disc_interface.h:237
void assemble_jacobian(matrix_type &J, ConstSmartPtr< VectorTimeSeries< vector_type > > vSol, const number s_a)
assembles jacobian on surface level
Definition: domain_disc_interface.h:258
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)=0
Assembles Right-Hand-Side for a linear problem.
void assemble_linear(matrix_type &A, vector_type &b, ConstSmartPtr< VectorTimeSeries< vector_type > > vSol, const std::vector< number > &vScaleMass, const std::vector< number > &vScaleStiff)
assembles linear on surface level
Definition: domain_disc_interface.h:315
virtual void prepare_timestep_elem(ConstSmartPtr< VectorTimeSeries< vector_type > > vSol, ConstSmartPtr< DoFDistribution > dd)=0
prepares time step element-wise
virtual void assemble_linear(matrix_type &A, vector_type &b, ConstSmartPtr< DoFDistribution > dd)=0
virtual void assemble_defect(vector_type &d, ConstSmartPtr< VectorTimeSeries< vector_type > > vSol, const std::vector< number > &vScaleMass, const std::vector< number > &vScaleStiff, ConstSmartPtr< DoFDistribution > dd)=0
virtual size_t num_constraints() const =0
returns the number of post processes
virtual ~IDomainDiscretization()
Definition: domain_disc_interface.h:148
TAlgebra algebra_type
Algebra type.
Definition: domain_disc_interface.h:141
void assemble_rhs(vector_type &b, ConstSmartPtr< VectorTimeSeries< vector_type > > vSol, const std::vector< number > &vScaleMass, const std::vector< number > &vScaleStiff)
assembles rhs on surface level
Definition: domain_disc_interface.h:343
void finish_timestep_elem(ConstSmartPtr< VectorTimeSeries< vector_type > > vSol)
finishes timestep on surface level
Definition: domain_disc_interface.h:389
virtual void finish_timestep(ConstSmartPtr< VectorTimeSeries< vector_type > > vSol)
Definition: domain_disc_interface.h:375
virtual void adjust_solution(vector_type &u, const GridLevel &gl)=0
sets dirichlet values in solution vector
void assemble_defect(vector_type &d, ConstSmartPtr< VectorTimeSeries< vector_type > > vSol, const std::vector< number > &vScaleMass, const std::vector< number > &vScaleStiff)
assembles defect on surface level
Definition: domain_disc_interface.h:286
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, ConstSmartPtr< DoFDistribution > dd)=0
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)=0
assembles Defect
virtual void finish_timestep(ConstSmartPtr< VectorTimeSeries< vector_type > > vSol, ConstSmartPtr< DoFDistribution > dd)=0
finishes time step
virtual void assemble_rhs(vector_type &b, const GridLevel &gl)=0
assembles the rhs
virtual void assemble_jacobian(matrix_type &J, ConstSmartPtr< VectorTimeSeries< vector_type > > vSol, const number s_a0, ConstSmartPtr< DoFDistribution > dd)=0
virtual void assemble_rhs(vector_type &rhs, const vector_type &u, const GridLevel &gl)=0
assembles the rhs
virtual void adjust_solution(vector_type &u, number time, ConstSmartPtr< DoFDistribution > dd)=0
algebra_type::matrix_type matrix_type
Type of algebra matrix.
Definition: domain_disc_interface.h:144
virtual void assemble_jacobian(matrix_type &J, const vector_type &u, const GridLevel &gl)=0
assembles Jacobian (or Approximation of Jacobian)
virtual void assemble_rhs(vector_type &rhs, const vector_type &u, ConstSmartPtr< DoFDistribution > dd)=0
virtual void assemble_stiffness_matrix(matrix_type &A, const vector_type &u, const GridLevel &gl)=0
assembles the stiffness matrix
virtual void prepare_timestep(ConstSmartPtr< VectorTimeSeries< vector_type > > vSol, number future_time, const GridLevel &gl)=0
virtual void adjust_solution(vector_type &u, ConstSmartPtr< DoFDistribution > dd)=0
virtual void adjust_solution(vector_type &u, number time, const GridLevel &gl)=0
sets dirichlet values in solution vector
virtual void finish_timestep(ConstSmartPtr< VectorTimeSeries< vector_type > > vSol, const GridLevel &gl)=0
virtual void assemble_stiffness_matrix(matrix_type &A, const vector_type &u, ConstSmartPtr< DoFDistribution > dd)=0
void adjust_solution(vector_type &u, number time)
adjust solution on surface level
Definition: domain_disc_interface.h:362
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)=0
Assembles matrix_type and Right-Hand-Side for a linear problem.
Interface for an object that can estimate the (global) error.
Definition: domain_disc_interface.h:51
virtual void calc_error(const vector_type &u, ConstSmartPtr< DoFDistribution > dd, error_vector_type *u_vtk=NULL)=0
virtual void calc_error(const vector_type &u, const GridLevel &gl, error_vector_type *u_vtk=NULL)=0
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)=0
Transient version.
virtual ~IDomainErrorIndicator()
Definition: domain_disc_interface.h:60
virtual bool is_error_valid()=0
returns whether current error values are valid
virtual void invalidate_error()=0
marks error indicators as invalid
TAlgebra::vector_type vector_type
Type of algebra vector.
Definition: domain_disc_interface.h:54
CPUAlgebra::vector_type error_vector_type
Type of error vector.
Definition: domain_disc_interface.h:57
virtual void calc_error(ConstSmartPtr< VectorTimeSeries< vector_type > > vSol, ConstSmartPtr< DoFDistribution > dd, const std::vector< number > &vScaleMass, const std::vector< number > &vScaleStiff, error_vector_type *u_vtk)=0
Transient version.
Interface for an object that can mark elements based on a strategy.
Definition: domain_disc_interface.h:107
virtual void mark_with_strategy(IRefiner &refiner, SmartPtr< element_marking_strategy_type > spMarkingStrategy)=0
virtual ~IDomainMarker()
Definition: domain_disc_interface.h:113
IElementMarkingStrategy< TDomain > element_marking_strategy_type
Type of algebra vector.
Definition: domain_disc_interface.h:110
Abstract base class for element marking (in adaptive refinement)
Definition: error_elem_marking_strategy.h:125
The refiner interface allows to mark elements for refinement and to call refine.
Definition: refiner_interface.h:67
time series of solutions and corresponding time point
Definition: solution_time_series.h:59
double number
Definition: types.h:124
the ug namespace