Flow123d  release_2.2.0-21-g2806cf6
advection_diffusion_model.hh
Go to the documentation of this file.
1 /*!
2  *
3  * Copyright (C) 2015 Technical University of Liberec. All rights reserved.
4  *
5  * This program is free software; you can redistribute it and/or modify it under
6  * the terms of the GNU General Public License version 3 as published by the
7  * Free Software Foundation. (http://www.gnu.org/licenses/gpl-3.0.en.html)
8  *
9  * This program is distributed in the hope that it will be useful, but WITHOUT
10  * ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
11  * FOR A PARTICULAR PURPOSE. See the GNU General Public License for more details.
12  *
13  *
14  * @file advection_diffusion_model.hh
15  * @brief Discontinuous Galerkin method for equation of transport with dispersion.
16  * @author Jan Stebel
17  */
18 
19 #ifndef AD_MODEL_HH_
20 #define AD_MODEL_HH_
21 
22 #include <armadillo>
23 #include <vector>
24 
25 class SubstanceList;
26 
27 namespace IT = Input::Type;
28 
29 /**
30  * AdvectionDiffusionModel is a base class for description of a physical process described
31  * by the advection-diffusion partial differential equation (PDE). The derived classes define input parameters
32  * and implement methods that calculate coefficients of the PDE. These methods are then used by a template
33  * class for numerical solution, whose specialization derives from the model class.
34  */
36 public:
37 
39 // abc_none,
44  };
45 
46  /// Read necessary data from input record.
47  virtual void init_from_input(const Input::Record &in_rec) = 0;
48 
49  /**
50  * Compute coefficients of mass matrix.
51  * @param point_list Points at which to evaluate.
52  * @param ele_acc Element accessor.
53  * @param mm_coef Coefficient vector (output).
54  */
55  virtual void compute_mass_matrix_coefficient(const std::vector<arma::vec3 > &point_list,
56  const ElementAccessor<3> &ele_acc,
57  std::vector<double> &mm_coef) = 0;
58 
59 
60  /**
61  * Compute retardation coefficients due to sorption.
62  * @param point_list Points at which to evaluate.
63  * @param ele_acc Element accessor.
64  * @param ret_coef Coefficient vector (output).
65  */
66  virtual void compute_retardation_coefficient(const std::vector<arma::vec3 > &point_list,
67  const ElementAccessor<3> &ele_acc,
68  std::vector<std::vector<double> > &ret_coef) = 0;
69 
70  /**
71  * Compute coefficients of stiffness matrix.
72  * @param point_list Points at which to evaluate.
73  * @param velocity Velocity field (input). Temporary solution before we can pass data from other
74  * equations.
75  * @param ele_acc Element accessor.
76  * @param ad_coef Coefficients of advection (output).
77  * @param dif_coef Coefficients of diffusion (output).
78  */
80  const std::vector<arma::vec3> &velocity,
81  const ElementAccessor<3> &ele_acc,
83  std::vector<std::vector<arma::mat33> > &dif_coef) = 0;
84 
85  /**
86  * Compute initial conditions.
87  * @param point_list Points at which to evaluate.
88  * @param ele_acc Element accessor.
89  * @param init_values Vector of intial values (output).
90  */
91  virtual void compute_init_cond(const std::vector<arma::vec3> &point_list,
92  const ElementAccessor<3> &ele_acc,
93  std::vector<std::vector<double> > &init_values) = 0;
94 
95  /**
96  * Return types of boundary conditions for each solution component.
97  * @param ele_acc Element accessor.
98  * @param bc_types Vector of bc. types (output, see BC_Type)
99  */
100  virtual void get_bc_type(const ElementAccessor<3> &ele_acc,
101  arma::uvec &bc_types) = 0;
102 
103  /**
104  * \brief Return data for diffusive or total flux b.c.
105  *
106  * The flux can in general take the form
107  *
108  * cross_section*(flux + sigma*(solution - ref_value))
109  *
110  * @param index Component index.
111  * @param point_list Points at which to evaluate.
112  * @param ele_acc Element accessor.
113  * @param bc_flux Neumann flux (output).
114  * @param bc_sigma Transition parameter (output).
115  * @param bc_ref_value Reference value (output).
116  */
117  virtual void get_flux_bc_data(unsigned int index,
118  const std::vector<arma::vec3> &point_list,
119  const ElementAccessor<3> &ele_acc,
120  std::vector< double > &bc_flux,
121  std::vector< double > &bc_sigma,
122  std::vector< double > &bc_ref_value) = 0;
123 
124  /**
125  * \brief Return transition coefficient for flux b.c.
126  *
127  * In assembly of system matrix one does not teed all data for total/diffusive flux b.c.
128  * This method therefore returns only the sigma coefficient.
129  *
130  * @param index Component index.
131  * @param point_list Points at which to evaluate.
132  * @param ele_acc Element accessor.
133  * @param bc_sigma Transition parameter (output).
134  */
135  virtual void get_flux_bc_sigma(unsigned int index,
136  const std::vector<arma::vec3> &point_list,
137  const ElementAccessor<3> &ele_acc,
138  std::vector< double > &bc_sigma) = 0;
139 
140  /**
141  * Compute coefficients of volume sources.
142  * @param point_list Points at which to evaluate.
143  * @param ele_acc Element accessor.
144  * @param sources_conc Source concentrations (output).
145  * @param sources_density Source densities (output).
146  * @param sources_sigma Source sigmas (output).
147  */
148  virtual void compute_source_coefficients(const std::vector<arma::vec3> &point_list,
149  const ElementAccessor<3> &ele_acc,
150  std::vector<std::vector<double> > &sources_conc,
151  std::vector<std::vector<double> > &sources_density,
152  std::vector<std::vector<double> > &sources_sigma) = 0;
153 
154  /**
155  * Compute coefficients of volume sources.
156  * @param point_list Points at which to evaluate.
157  * @param ele_acc Element accessor.
158  * @param sources_sigma Source sigmas (output).
159  */
160  virtual void compute_sources_sigma(const std::vector<arma::vec3> &point_list,
161  const ElementAccessor<3> &ele_acc,
162  std::vector<std::vector<double> > &sources_sigma) = 0;
163 
164  /// Destructor.
166 
167 
168 
169 };
170 
171 
172 
173 
174 
175 
176 #endif /* AD_MODEL_HH_ */
virtual ~AdvectionDiffusionModel()
Destructor.
virtual void compute_init_cond(const std::vector< arma::vec3 > &point_list, const ElementAccessor< 3 > &ele_acc, std::vector< std::vector< double > > &init_values)=0
virtual void compute_sources_sigma(const std::vector< arma::vec3 > &point_list, const ElementAccessor< 3 > &ele_acc, std::vector< std::vector< double > > &sources_sigma)=0
virtual void compute_source_coefficients(const std::vector< arma::vec3 > &point_list, const ElementAccessor< 3 > &ele_acc, std::vector< std::vector< double > > &sources_conc, std::vector< std::vector< double > > &sources_density, std::vector< std::vector< double > > &sources_sigma)=0
virtual void get_flux_bc_sigma(unsigned int index, const std::vector< arma::vec3 > &point_list, const ElementAccessor< 3 > &ele_acc, std::vector< double > &bc_sigma)=0
Return transition coefficient for flux b.c.
virtual void get_flux_bc_data(unsigned int index, const std::vector< arma::vec3 > &point_list, const ElementAccessor< 3 > &ele_acc, std::vector< double > &bc_flux, std::vector< double > &bc_sigma, std::vector< double > &bc_ref_value)=0
Return data for diffusive or total flux b.c.
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:292
virtual void init_from_input(const Input::Record &in_rec)=0
Read necessary data from input record.
virtual void get_bc_type(const ElementAccessor< 3 > &ele_acc, arma::uvec &bc_types)=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_retardation_coefficient(const std::vector< arma::vec3 > &point_list, const ElementAccessor< 3 > &ele_acc, std::vector< std::vector< double > > &ret_coef)=0