Flow123d
Main Page
Related Pages
Modules
Namespaces
Classes
Files
File List
File Members
flow123d
src
transport
advection_diffusion_model.hh
Go to the documentation of this file.
1
/*!
2
*
3
* Copyright (C) 2007 Technical University of Liberec. All rights reserved.
4
*
5
* Please make a following refer to Flow123d on your project site if you use the program for any purpose,
6
* especially for academic research:
7
* Flow123d, Research Centre: Advanced Remedial Technologies, Technical University of Liberec, Czech Republic
8
*
9
* This program is free software; you can redistribute it and/or modify it under the terms
10
* of the GNU General Public License version 3 as published by the Free Software Foundation.
11
*
12
* This program is distributed in the hope that it will be useful, but WITHOUT ANY WARRANTY;
13
* without even the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.
14
* See the GNU General Public License for more details.
15
*
16
* You should have received a copy of the GNU General Public License along with this program; if not,
17
* write to the Free Software Foundation, Inc., 59 Temple Place - Suite 330, Boston, MA 021110-1307, USA.
18
*
19
*
20
* $Id$
21
* $Revision$
22
* $LastChangedBy$
23
* $LastChangedDate$
24
*
25
* @file
26
* @brief Discontinuous Galerkin method for equation of transport with dispersion.
27
* @author Jan Stebel
28
*/
29
30
#ifndef AD_MODEL_HH_
31
#define AD_MODEL_HH_
32
33
#include <armadillo>
34
#include <vector>
35
36
namespace
IT = Input::Type;
37
38
/**
39
* AdvectionDiffusionModel is a base class for description of a physical process described
40
* by the advection-diffusion partial differential equation (PDE). The derived classes define input parameters
41
* and implement methods that calculate coefficients of the PDE. These methods are then used by a template
42
* class for numerical solution, whose specialization derives from the model class.
43
*/
44
class
AdvectionDiffusionModel
{
45
public
:
46
47
/// Initialize model data. E.g. set vector field dimensions.
48
virtual
void
init_data
(
unsigned
int
n_subst_) = 0;
49
50
/// Temporary solution for sharing data from other equations.
51
//virtual void set_cross_section_field(const Field< 3, FieldValue<3>::Scalar > &cross_section) = 0;
52
53
/// Read or set names of solution components.
54
virtual
void
set_component_names
(
std::vector<string>
&names,
const
Input::Record
&in_rec) = 0;
55
56
/// Check if mass matrix coefficients have changed.
57
virtual
bool
mass_matrix_changed
() = 0;
58
59
/// Check if stiffness matrix coefficients have changed.
60
virtual
bool
stiffness_matrix_changed
() = 0;
61
62
/// Check if right hand side coefficients have changed.
63
virtual
bool
rhs_changed
() = 0;
64
65
/**
66
* Compute coefficients of mass matrix.
67
* @param point_list Points at which to evaluate.
68
* @param ele_acc Element accessor.
69
* @param mm_coef Coefficient vector (output).
70
*/
71
virtual
void
compute_mass_matrix_coefficient
(
const
std::vector<arma::vec3 >
&point_list,
72
const
ElementAccessor<3>
&ele_acc,
73
std::vector<double>
&mm_coef) = 0;
74
75
/**
76
* Compute coefficients of stiffness matrix.
77
* @param point_list Points at which to evaluate.
78
* @param velocity Velocity field (input). Temporary solution before we can pass data from other
79
* equations.
80
* @param ele_acc Element accessor.
81
* @param ad_coef Coefficients of advection (output).
82
* @param dif_coef Coefficients of diffusion (output).
83
*/
84
virtual
void
compute_advection_diffusion_coefficients
(
const
std::vector<arma::vec3>
&point_list,
85
const
std::vector<arma::vec3>
&velocity,
86
const
ElementAccessor<3>
&ele_acc,
87
std::vector
<
std::vector<arma::vec3>
> &ad_coef,
88
std::vector
<
std::vector<arma::mat33>
> &dif_coef) = 0;
89
90
/**
91
* Compute initial conditions.
92
* @param point_list Points at which to evaluate.
93
* @param ele_acc Element accessor.
94
* @param init_values Vector of intial values (output).
95
*/
96
virtual
void
compute_init_cond
(
const
std::vector<arma::vec3>
&point_list,
97
const
ElementAccessor<3>
&ele_acc,
98
std::vector< arma::vec >
&init_values) = 0;
99
100
/**
101
* Computes the Dirichlet boundary condition values.
102
* @param point_list Points at which to evaluate.
103
* @param ele_acc Element accessor.
104
* @param bc_values Vector of b.c. values (output).
105
*/
106
virtual
void
compute_dirichlet_bc
(
const
std::vector<arma::vec3>
&point_list,
107
const
ElementAccessor<3>
&ele_acc,
108
std::vector< arma::vec >
&bc_values) = 0;
109
110
/**
111
* Compute coefficients of volume sources.
112
* @param point_list Points at which to evaluate.
113
* @param ele_acc Element accessor.
114
* @param sources_conc Source concentrations (output).
115
* @param sources_density Source densities (output).
116
* @param sources_sigma Source sigmas (output).
117
*/
118
virtual
void
compute_source_coefficients
(
const
std::vector<arma::vec3>
&point_list,
119
const
ElementAccessor<3>
&ele_acc,
120
std::vector<arma::vec>
&sources_conc,
121
std::vector<arma::vec>
&sources_density,
122
std::vector<arma::vec>
&sources_sigma) = 0;
123
124
/**
125
* Compute coefficients of volume sources.
126
* @param point_list Points at which to evaluate.
127
* @param ele_acc Element accessor.
128
* @param sources_sigma Source sigmas (output).
129
*/
130
virtual
void
compute_sources_sigma
(
const
std::vector<arma::vec3>
&point_list,
131
const
ElementAccessor<3>
&ele_acc,
132
std::vector<arma::vec>
&sources_sigma) = 0;
133
134
/// Destructor.
135
virtual
~AdvectionDiffusionModel
() {};
136
137
138
139
};
140
141
142
143
144
145
146
#endif
/* AD_MODEL_HH_ */
Generated on Thu May 29 2014 23:14:49 for Flow123d by
1.8.4