Flow123d  JS_before_hm-1003-g4e68d2c
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 
26 namespace IT = Input::Type;
27 
28 /**
29  * AdvectionDiffusionModel is a base class for description of a physical process described
30  * by the advection-diffusion partial differential equation (PDE). The derived classes define input parameters
31  * and implement methods that calculate coefficients of the PDE. These methods are then used by a template
32  * class for numerical solution, whose specialization derives from the model class.
33  */
35 public:
36 
38 // abc_none,
43  };
44 
45  /// Read necessary data from input record.
46  virtual void init_from_input(const Input::Record &in_rec) = 0;
47 
48  /**
49  * Compute coefficients of mass matrix.
50  * @param point_list Points at which to evaluate.
51  * @param ele_acc Element accessor.
52  * @param mm_coef Coefficient vector (output).
53  */
54  virtual void compute_mass_matrix_coefficient(const Armor::array &point_list,
55  const ElementAccessor<3> &ele_acc,
56  std::vector<double> &mm_coef) = 0;
57 
58 
59  /**
60  * Compute retardation coefficients due to sorption.
61  * @param point_list Points at which to evaluate.
62  * @param ele_acc Element accessor.
63  * @param ret_coef Coefficient vector (output).
64  */
65  virtual void compute_retardation_coefficient(const Armor::array &point_list,
66  const ElementAccessor<3> &ele_acc,
67  std::vector<std::vector<double> > &ret_coef) = 0;
68 
69  /**
70  * Compute coefficients of stiffness matrix.
71  * @param point_list Points at which to evaluate.
72  * @param velocity Velocity field (input). Temporary solution before we can pass data from other
73  * equations.
74  * @param ele_acc Element accessor.
75  * @param ad_coef Coefficients of advection (output).
76  * @param dif_coef Coefficients of diffusion (output).
77  */
78  virtual void compute_advection_diffusion_coefficients(const Armor::array &point_list,
79  const std::vector<arma::vec3> &velocity,
80  const ElementAccessor<3> &ele_acc,
82  std::vector<std::vector<arma::mat33> > &dif_coef) = 0;
83 
84  /**
85  * Compute initial conditions.
86  * @param point_list Points at which to evaluate.
87  * @param ele_acc Element accessor.
88  * @param init_values Vector of intial values (output).
89  */
90  virtual void compute_init_cond(const Armor::array &point_list,
91  const ElementAccessor<3> &ele_acc,
92  std::vector<std::vector<double> > &init_values) = 0;
93 
94  /**
95  * Return types of boundary conditions for each solution component.
96  * @param ele_acc Element accessor.
97  * @param bc_types Vector of bc. types (output, see BC_Type)
98  */
99  virtual void get_bc_type(const ElementAccessor<3> &ele_acc,
100  arma::uvec &bc_types) = 0;
101 
102  /**
103  * \brief Return data for diffusive or total flux b.c.
104  *
105  * The flux can in general take the form
106  *
107  * cross_section*(flux + sigma*(solution - ref_value))
108  *
109  * @param index Component index.
110  * @param point_list Points at which to evaluate.
111  * @param ele_acc Element accessor.
112  * @param bc_flux Neumann flux (output).
113  * @param bc_sigma Transition parameter (output).
114  * @param bc_ref_value Reference value (output).
115  */
116  virtual void get_flux_bc_data(unsigned int index,
117  const Armor::array &point_list,
118  const ElementAccessor<3> &ele_acc,
119  std::vector< double > &bc_flux,
120  std::vector< double > &bc_sigma,
121  std::vector< double > &bc_ref_value) = 0;
122 
123  /**
124  * \brief Return transition coefficient for flux b.c.
125  *
126  * In assembly of system matrix one does not teed all data for total/diffusive flux b.c.
127  * This method therefore returns only the sigma coefficient.
128  *
129  * @param index Component index.
130  * @param point_list Points at which to evaluate.
131  * @param ele_acc Element accessor.
132  * @param bc_sigma Transition parameter (output).
133  */
134  virtual void get_flux_bc_sigma(unsigned int index,
135  const Armor::array &point_list,
136  const ElementAccessor<3> &ele_acc,
137  std::vector< double > &bc_sigma) = 0;
138 
139  /**
140  * Compute coefficients of volume sources.
141  * @param point_list Points at which to evaluate.
142  * @param ele_acc Element accessor.
143  * @param sources_conc Source concentrations (output).
144  * @param sources_density Source densities (output).
145  * @param sources_sigma Source sigmas (output).
146  */
147  virtual void compute_source_coefficients(const Armor::array &point_list,
148  const ElementAccessor<3> &ele_acc,
149  std::vector<std::vector<double> > &sources_conc,
150  std::vector<std::vector<double> > &sources_density,
151  std::vector<std::vector<double> > &sources_sigma) = 0;
152 
153  /**
154  * Compute coefficients of volume sources.
155  * @param point_list Points at which to evaluate.
156  * @param ele_acc Element accessor.
157  * @param sources_sigma Source sigmas (output).
158  */
159  virtual void compute_sources_sigma(const Armor::array &point_list,
160  const ElementAccessor<3> &ele_acc,
161  std::vector<std::vector<double> > &sources_sigma) = 0;
162 
163  /// Destructor.
165 
166 
167 
168 };
169 
170 
171 
172 
173 
174 
175 #endif /* AD_MODEL_HH_ */
virtual ~AdvectionDiffusionModel()
Destructor.
virtual void compute_retardation_coefficient(const Armor::array &point_list, const ElementAccessor< 3 > &ele_acc, std::vector< std::vector< double > > &ret_coef)=0
virtual void compute_source_coefficients(const Armor::array &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 compute_init_cond(const Armor::array &point_list, const ElementAccessor< 3 > &ele_acc, std::vector< std::vector< double > > &init_values)=0
virtual void get_flux_bc_sigma(unsigned int index, const Armor::array &point_list, const ElementAccessor< 3 > &ele_acc, std::vector< double > &bc_sigma)=0
Return transition coefficient for flux b.c.
virtual void compute_sources_sigma(const Armor::array &point_list, const ElementAccessor< 3 > &ele_acc, std::vector< std::vector< double > > &sources_sigma)=0
Accessor to the data with type Type::Record.
Definition: accessors.hh:291
virtual void init_from_input(const Input::Record &in_rec)=0
Read necessary data from input record.
virtual void get_flux_bc_data(unsigned int index, const Armor::array &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 Armor::array &point_list, const ElementAccessor< 3 > &ele_acc, std::vector< double > &mm_coef)=0
virtual void get_bc_type(const ElementAccessor< 3 > &ele_acc, arma::uvec &bc_types)=0
virtual void compute_advection_diffusion_coefficients(const Armor::array &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