Flow123d  jenkins-Flow123d-windows32-release-multijob-51
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  */
45 public:
46 
47  /// Read or set names of solution components.
48  virtual void set_component_names(std::vector<string> &names, const Input::Record &in_rec) = 0;
49 
50  /**
51  * Compute coefficients of mass matrix.
52  * @param point_list Points at which to evaluate.
53  * @param ele_acc Element accessor.
54  * @param mm_coef Coefficient vector (output).
55  */
56  virtual void compute_mass_matrix_coefficient(const std::vector<arma::vec3 > &point_list,
57  const ElementAccessor<3> &ele_acc,
58  std::vector<double> &mm_coef) = 0;
59 
60  /**
61  * Compute coefficients of stiffness matrix.
62  * @param point_list Points at which to evaluate.
63  * @param velocity Velocity field (input). Temporary solution before we can pass data from other
64  * equations.
65  * @param ele_acc Element accessor.
66  * @param ad_coef Coefficients of advection (output).
67  * @param dif_coef Coefficients of diffusion (output).
68  */
70  const std::vector<arma::vec3> &velocity,
71  const ElementAccessor<3> &ele_acc,
73  std::vector<std::vector<arma::mat33> > &dif_coef) = 0;
74 
75  /**
76  * Compute initial conditions.
77  * @param point_list Points at which to evaluate.
78  * @param ele_acc Element accessor.
79  * @param init_values Vector of intial values (output).
80  */
81  virtual void compute_init_cond(const std::vector<arma::vec3> &point_list,
82  const ElementAccessor<3> &ele_acc,
83  std::vector< arma::vec > &init_values) = 0;
84 
85  /**
86  * Computes the Dirichlet boundary condition values.
87  * @param point_list Points at which to evaluate.
88  * @param ele_acc Element accessor.
89  * @param bc_values Vector of b.c. values (output).
90  */
91  virtual void compute_dirichlet_bc(const std::vector<arma::vec3> &point_list,
92  const ElementAccessor<3> &ele_acc,
93  std::vector< arma::vec > &bc_values) = 0;
94 
95  /**
96  * Compute coefficients of volume sources.
97  * @param point_list Points at which to evaluate.
98  * @param ele_acc Element accessor.
99  * @param sources_conc Source concentrations (output).
100  * @param sources_density Source densities (output).
101  * @param sources_sigma Source sigmas (output).
102  */
103  virtual void compute_source_coefficients(const std::vector<arma::vec3> &point_list,
104  const ElementAccessor<3> &ele_acc,
105  std::vector<arma::vec> &sources_conc,
106  std::vector<arma::vec> &sources_density,
107  std::vector<arma::vec> &sources_sigma) = 0;
108 
109  /**
110  * Compute coefficients of volume sources.
111  * @param point_list Points at which to evaluate.
112  * @param ele_acc Element accessor.
113  * @param sources_sigma Source sigmas (output).
114  */
115  virtual void compute_sources_sigma(const std::vector<arma::vec3> &point_list,
116  const ElementAccessor<3> &ele_acc,
117  std::vector<arma::vec> &sources_sigma) = 0;
118 
119  /// Destructor.
121 
122 
123 
124 };
125 
126 
127 
128 
129 
130 
131 #endif /* AD_MODEL_HH_ */
virtual ~AdvectionDiffusionModel()
Destructor.
virtual void set_component_names(std::vector< string > &names, 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:308
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