Flow123d  jenkins-Flow123d-windows-release-multijob-285
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 class SubstanceList;
37 
38 namespace IT = Input::Type;
39 
40 /**
41  * AdvectionDiffusionModel is a base class for description of a physical process described
42  * by the advection-diffusion partial differential equation (PDE). The derived classes define input parameters
43  * and implement methods that calculate coefficients of the PDE. These methods are then used by a template
44  * class for numerical solution, whose specialization derives from the model class.
45  */
47 public:
48 
49  /// Read or set names of solution components.
50  virtual void set_components(SubstanceList &substances, const Input::Record &in_rec) = 0;
51 
52  /**
53  * Compute coefficients of mass matrix.
54  * @param point_list Points at which to evaluate.
55  * @param ele_acc Element accessor.
56  * @param mm_coef Coefficient vector (output).
57  */
58  virtual void compute_mass_matrix_coefficient(const std::vector<arma::vec3 > &point_list,
59  const ElementAccessor<3> &ele_acc,
60  std::vector<double> &mm_coef) = 0;
61 
62  /**
63  * Compute coefficients of stiffness matrix.
64  * @param point_list Points at which to evaluate.
65  * @param velocity Velocity field (input). Temporary solution before we can pass data from other
66  * equations.
67  * @param ele_acc Element accessor.
68  * @param ad_coef Coefficients of advection (output).
69  * @param dif_coef Coefficients of diffusion (output).
70  */
72  const std::vector<arma::vec3> &velocity,
73  const ElementAccessor<3> &ele_acc,
75  std::vector<std::vector<arma::mat33> > &dif_coef) = 0;
76 
77  /**
78  * Compute initial conditions.
79  * @param point_list Points at which to evaluate.
80  * @param ele_acc Element accessor.
81  * @param init_values Vector of intial values (output).
82  */
83  virtual void compute_init_cond(const std::vector<arma::vec3> &point_list,
84  const ElementAccessor<3> &ele_acc,
85  std::vector< arma::vec > &init_values) = 0;
86 
87  /**
88  * Computes the Dirichlet boundary condition values.
89  * @param point_list Points at which to evaluate.
90  * @param ele_acc Element accessor.
91  * @param bc_values Vector of b.c. values (output).
92  */
93  virtual void compute_dirichlet_bc(const std::vector<arma::vec3> &point_list,
94  const ElementAccessor<3> &ele_acc,
95  std::vector< arma::vec > &bc_values) = 0;
96 
97  /**
98  * Compute coefficients of volume sources.
99  * @param point_list Points at which to evaluate.
100  * @param ele_acc Element accessor.
101  * @param sources_conc Source concentrations (output).
102  * @param sources_density Source densities (output).
103  * @param sources_sigma Source sigmas (output).
104  */
105  virtual void compute_source_coefficients(const std::vector<arma::vec3> &point_list,
106  const ElementAccessor<3> &ele_acc,
107  std::vector<arma::vec> &sources_conc,
108  std::vector<arma::vec> &sources_density,
109  std::vector<arma::vec> &sources_sigma) = 0;
110 
111  /**
112  * Compute coefficients of volume sources.
113  * @param point_list Points at which to evaluate.
114  * @param ele_acc Element accessor.
115  * @param sources_sigma Source sigmas (output).
116  */
117  virtual void compute_sources_sigma(const std::vector<arma::vec3> &point_list,
118  const ElementAccessor<3> &ele_acc,
119  std::vector<arma::vec> &sources_sigma) = 0;
120 
121  /// Destructor.
123 
124 
125 
126 };
127 
128 
129 
130 
131 
132 
133 #endif /* AD_MODEL_HH_ */
virtual ~AdvectionDiffusionModel()
Destructor.
virtual void set_components(SubstanceList &substances, const Input::Record &in_rec)=0
Read or set names of solution components.
virtual void compute_sources_sigma(const std::vector< arma::vec3 > &point_list, const ElementAccessor< 3 > &ele_acc, std::vector< arma::vec > &sources_sigma)=0
virtual void compute_mass_matrix_coefficient(const std::vector< arma::vec3 > &point_list, const ElementAccessor< 3 > &ele_acc, std::vector< double > &mm_coef)=0
Accessor to the data with type Type::Record.
Definition: accessors.hh:327
virtual void compute_init_cond(const std::vector< arma::vec3 > &point_list, const ElementAccessor< 3 > &ele_acc, std::vector< arma::vec > &init_values)=0
virtual void compute_advection_diffusion_coefficients(const std::vector< arma::vec3 > &point_list, const std::vector< arma::vec3 > &velocity, const ElementAccessor< 3 > &ele_acc, std::vector< std::vector< arma::vec3 > > &ad_coef, std::vector< std::vector< arma::mat33 > > &dif_coef)=0
virtual void compute_source_coefficients(const std::vector< arma::vec3 > &point_list, const ElementAccessor< 3 > &ele_acc, std::vector< arma::vec > &sources_conc, std::vector< arma::vec > &sources_density, std::vector< arma::vec > &sources_sigma)=0
virtual void compute_dirichlet_bc(const std::vector< arma::vec3 > &point_list, const ElementAccessor< 3 > &ele_acc, std::vector< arma::vec > &bc_values)=0