Loading [MathJax]/extensions/tex2jax.js
ug4
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Modules Pages
composite_time_disc.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__UG__LIB_DISC__TIME_DISC__COMBINED_TIME_DISC__
34#define __H__UG__LIB_DISC__TIME_DISC__COMBINED_TIME_DISC__
35
36#include <cstddef> // for size_t
37#include <vector> // for vector
38
39#include "lib_disc/time_disc/time_disc_interface.h" // for ITimeDiscretization
40
41
42namespace ug {
43
46
48template <typename TAlgebra>
50 : public IAssemble<TAlgebra>
51{
52 public:
53 typedef TAlgebra algebra_type;
54 typedef typename algebra_type::vector_type vector_type;
55 typedef typename algebra_type::matrix_type matrix_type;
56
57 public:
60
61 public:
69 : public AssemblingTuner<TAlgebra>
70 {
71 public:
73 {m_vAssTuner.push_back(assTuner);}
74
75 void set_force_regular_grid(bool bForce)
76 {
77 for (size_t i = 0; i < m_vAssTuner.size(); ++i)
79 }
80
81 protected:
82 std::vector<SmartPtr<AssemblingTuner<TAlgebra>, FreeDelete> > m_vAssTuner;
83 };
84
85 // inherited from ITimeDiscretization
86 public:
89
91 virtual void prepare_step_elem
92 (
94 number dt,
95 const GridLevel& gl
96 );
97
100
102 virtual void finish_step_elem
103 (
105 const GridLevel& gl
106 );
107
109 virtual number future_time() const;
110
112 virtual size_t num_prev_steps() const;
113
115 virtual size_t num_stages() const;
116
118 virtual void set_stage(size_t stage);
119
120 // inherited from IAssemble
121 public:
123 void assemble_jacobian(matrix_type& J, const vector_type& u, const GridLevel& gl);
124
126 void assemble_defect(vector_type& d, const vector_type& u, const GridLevel& gl);
127
129 void assemble_linear(matrix_type& A, vector_type& b, const GridLevel& gl);
130
132 void assemble_rhs(vector_type& b, const vector_type& u, const GridLevel& gl);
133
135 void assemble_rhs(vector_type& b, const GridLevel& gl);
136
138 void adjust_solution(vector_type& u, const GridLevel& gl);
139
145
147 virtual size_t num_constraints() const;
148
150 virtual SmartPtr<IConstraint<TAlgebra> > constraint(size_t i);
151
152 protected:
153 std::vector<SmartPtr<ITimeDiscretization<TAlgebra> > > m_vTimeDisc;
154};
155
156} // end namespace ug
157
158
160
161// include implementation
163
164#endif // __H__UG__LIB_DISC__TIME_DISC__COMBINED_TIME_DISC__
Definition smart_pointer.h:296
Definition smart_pointer.h:51
Definition smart_pointer.h:108
The AssemblingTuner class combines tools to adapt the assembling routine.
Definition ass_tuner.h:90
Definition composite_time_disc.h:70
std::vector< SmartPtr< AssemblingTuner< TAlgebra >, FreeDelete > > m_vAssTuner
Definition composite_time_disc.h:82
void add_ass_tuner(SmartPtr< AssemblingTuner< TAlgebra > > assTuner)
Definition composite_time_disc.h:72
void set_force_regular_grid(bool bForce)
forces the assembling to consider the grid as regular
Definition composite_time_disc.h:75
combine several time discretizations into one
Definition composite_time_disc.h:51
virtual void prepare_step(SmartPtr< VectorTimeSeries< vector_type > > prevSol, number dt)
prepares the assembling of Defect/Jacobian for a time step
Definition composite_time_disc_impl.h:40
virtual number future_time() const
returns the future time point (i.e. the one that will be computed)
Definition composite_time_disc_impl.h:78
virtual size_t num_constraints() const
returns the number of constraints
Definition composite_time_disc_impl.h:237
virtual SmartPtr< IConstraint< TAlgebra > > constraint(size_t i)
returns the i'th constraint
Definition composite_time_disc_impl.h:247
void assemble_rhs(vector_type &b, const vector_type &u, const GridLevel &gl)
Definition composite_time_disc_impl.h:178
std::vector< SmartPtr< ITimeDiscretization< TAlgebra > > > m_vTimeDisc
Definition composite_time_disc.h:153
virtual size_t num_stages() const
returns the number of stages
Definition composite_time_disc_impl.h:96
virtual void finish_step(SmartPtr< VectorTimeSeries< vector_type > > currSol)
Definition composite_time_disc_impl.h:60
void assemble_jacobian(matrix_type &J, const vector_type &u, const GridLevel &gl)
assembles Jacobian (or Approximation of Jacobian)
Definition composite_time_disc_impl.h:115
virtual void finish_step_elem(SmartPtr< VectorTimeSeries< vector_type > > currSol, const GridLevel &gl)
Definition composite_time_disc_impl.h:68
void assemble_defect(vector_type &d, const vector_type &u, const GridLevel &gl)
assembles Defect
Definition composite_time_disc_impl.h:136
void assemble_linear(matrix_type &A, vector_type &b, const GridLevel &gl)
Assembles Matrix and Right-Hand-Side for a linear problem.
Definition composite_time_disc_impl.h:157
virtual size_t num_prev_steps() const
returns number of previous time steps needed
Definition composite_time_disc_impl.h:86
virtual SmartPtr< AssemblingTuner< TAlgebra > > ass_tuner()
Definition composite_time_disc_impl.h:221
algebra_type::matrix_type matrix_type
Definition composite_time_disc.h:55
algebra_type::vector_type vector_type
Definition composite_time_disc.h:54
void adjust_solution(vector_type &u, const GridLevel &gl)
sets dirichlet values in solution vector
Definition composite_time_disc_impl.h:214
virtual void prepare_step_elem(SmartPtr< VectorTimeSeries< vector_type > > prevSol, number dt, const GridLevel &gl)
prepares the assembling of Defect/Jacobian for a time step
Definition composite_time_disc_impl.h:48
TAlgebra algebra_type
Definition composite_time_disc.h:53
virtual void set_stage(size_t stage)
sets the stage
Definition composite_time_disc_impl.h:106
void add_time_disc(SmartPtr< ITimeDiscretization< TAlgebra > > tDisc)
Definition composite_time_disc.h:58
Definition grid_level.h:42
Interface providing Jacobian and Defect of a discretization.
Definition assemble_interface.h:110
Time Discretization Interface.
Definition time_disc_interface.h:65
time series of solutions and corresponding time point
Definition solution_time_series.h:59
double number
Definition types.h:124
the ug namespace