Flow123d  JS_before_hm-62-ge56f9d5
concentration_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 concentration_model.hh
15  * @brief Discontinuous Galerkin method for equation of transport with dispersion.
16  * @author Jan Stebel
17  */
18 
19 #ifndef CONC_TRANS_MODEL_HH_
20 #define CONC_TRANS_MODEL_HH_
21 
22 #include <boost/exception/info.hpp> // for operator<<, err...
23 #include <memory> // for shared_ptr
24 #include <string> // for string
25 #include <vector> // for vector
26 #include <armadillo>
28 #include "fields/field_values.hh" // for FieldValue<>::S...
29 #include "fields/bc_multi_field.hh"
30 #include "fields/field.hh"
31 #include "fields/multi_field.hh"
32 #include "input/type_base.hh" // for Array
33 #include "input/type_generic.hh" // for Instance
34 #include "input/type_record.hh" // for Record
35 #include "input/type_selection.hh" // for Selection
36 #include "transport/substance.hh" // for SubstanceList
37 #include "transport/transport_operator_splitting.hh" // for ConcentrationTr...
38 
39 class Balance;
40 class MH_DofHandler;
41 class Mesh;
42 class OutputTime;
43 namespace Input { class Record; }
44 template <int spacedim> class ElementAccessor;
45 
46 
47 
48 
50 public:
51 
52  class ModelEqData : public TransportEqData {
53  public:
54 
59  bc_diffusive_flux
60  };
61 
62  /// Type of boundary condition (see also BC_Type)
64  /// Prescribed concentration for Dirichlet/reference concentration for flux b.c.
66  /// Flux value in total/diffusive flux b.c.
68  /// Transition coefficient in total/diffusive flux b.c.
70  /// Initial concentrations.
72  /// Longitudal dispersivity (for each substance).
74  /// Transversal dispersivity (for each substance).
76  /// Molecular diffusivity (for each substance).
78 
79  Field<3, FieldValue<3>::Scalar > rock_density; ///< Rock matrix density.
80  MultiField<3, FieldValue<3>::Scalar > sorption_coefficient; ///< Coefficient of linear sorption.
81 
82 
84 
85 
86 
87  ModelEqData();
88 
89  static constexpr const char * name() { return "Solute_AdvectionDiffusion"; }
90 
91  static string default_output_field() { return "\"conc\""; }
92 
93  static const Input::Type::Selection & get_bc_type_selection();
94 
95  static IT::Selection get_output_selection();
96 
97  };
98 
99 
101 
102 
103  ConcentrationTransportModel(Mesh &mesh, const Input::Record &in_rec);
104 
105  void init_from_input(const Input::Record &in_rec) override;
106 
107  void compute_mass_matrix_coefficient(const std::vector<arma::vec3 > &point_list,
108  const ElementAccessor<3> &ele_acc,
109  std::vector<double> &mm_coef) override;
110 
111  void compute_retardation_coefficient(const std::vector<arma::vec3 > &point_list,
112  const ElementAccessor<3> &ele_acc,
113  std::vector<std::vector<double> > &ret_coef) override;
114 
115 
116 
117  void compute_advection_diffusion_coefficients(const std::vector<arma::vec3 > &point_list,
118  const std::vector<arma::vec3> &velocity,
119  const ElementAccessor<3> &ele_acc,
121  std::vector<std::vector<arma::mat33> > &dif_coef) override;
122 
123  void compute_init_cond(const std::vector<arma::vec3> &point_list,
124  const ElementAccessor<3> &ele_acc,
125  std::vector<std::vector<double> > &init_values) override;
126 
127  void get_bc_type(const ElementAccessor<3> &ele_acc,
128  arma::uvec &bc_types) override;
129 
130  void get_flux_bc_data(unsigned int index,
131  const std::vector<arma::vec3> &point_list,
132  const ElementAccessor<3> &ele_acc,
133  std::vector< double > &bc_flux,
134  std::vector< double > &bc_sigma,
135  std::vector< double > &bc_ref_value) override;
136 
137  void get_flux_bc_sigma(unsigned int index,
138  const std::vector<arma::vec3> &point_list,
139  const ElementAccessor<3> &ele_acc,
140  std::vector< double > &bc_sigma) override;
141 
142  void compute_source_coefficients(const std::vector<arma::vec3> &point_list,
143  const ElementAccessor<3> &ele_acc,
144  std::vector<std::vector<double> > &sources_conc,
145  std::vector<std::vector<double> > &sources_density,
146  std::vector<std::vector<double> > &sources_sigma) override;
147 
148  void compute_sources_sigma(const std::vector<arma::vec3> &point_list,
149  const ElementAccessor<3> &ele_acc,
150  std::vector<std::vector<double> > &sources_sigma) override;
151 
152  ~ConcentrationTransportModel() override;
153 
154 
155  /**
156  * @brief Updates the velocity field which determines some coefficients of the transport equation.
157  *
158  * @param dh mixed hybrid dof handler
159  *
160  * (So far it does not work since the flow module returns a vector of zeros.)
161  * @param velocity_vector Input array of velocity values.
162  */
163  inline void set_velocity_field(const MH_DofHandler &dh) override
164  {
165  mh_dh = &dh;
166  flux_changed = true;
167  }
168 
169  /// Returns number of transported substances.
170  inline unsigned int n_substances() override
171  { return substances_.size(); }
172 
173  /// Returns reference to the vector of substance names.
174  inline SubstanceList &substances() override
175  { return substances_; }
176 
177 
178  // Methods inherited from ConcentrationTransportBase:
179 
180  // Must be implemented in descendants.
181  void set_target_time(double) override {};
182 
183  void set_balance_object(std::shared_ptr<Balance> balance) override;
184 
186  { return subst_idx; }
187 
188  void set_output_stream(std::shared_ptr<OutputTime> stream)
189  { output_stream_ = stream; }
190 
191  std::shared_ptr<OutputTime> output_stream() override
192  { return output_stream_; }
193 
194 
195 
196 protected:
197 
198  /// Derived class should implement getter for ModelEqData instance.
199  virtual ModelEqData &data() = 0;
200 
201  /**
202  * Create input type that can be passed to the derived class.
203  * @param implementation String characterizing the numerical method, e.g. DG, FEM, FVM.
204  * @param description Comment used to describe the record key.
205  * @return
206  */
207  static IT::Record get_input_type(const string &implementation, const string &description);
208 
209  /**
210  * Formula to calculate the dispersivity tensor.
211  * @param velocity Fluid velocity.
212  * @param Dm Molecular diffusivity.
213  * @param alphaL Longitudal dispersivity.
214  * @param alphaT Transversal dispersivity.
215  * @param porosity Porosity.
216  * @param cross_cut Cross-section.
217  * @param K Dispersivity tensor (output).
218  */
219  void calculate_dispersivity_tensor(const arma::vec3 &velocity,
220  const arma::mat33 &Dm,
221  double alphaL,
222  double alphaT,
223  double water_content,
224  double porosity,
225  double cross_cut,
226  arma::mat33 &K);
227 
228  /// Indicator of change in advection vector field.
230 
231  /// Transported substances.
233 
234  /// List of indices used to call balance methods for a set of quantities.
236 
237  /// Density of liquid (a global constant).
239 
240  /**
241  * Temporary solution how to pass velocity field form the flow model.
242  * TODO: introduce FieldDiscrete -containing true DOFHandler and data vector and pass such object together with other
243  * data. Possibly make more general set_data method, allowing setting data given by name. needs support from EqDataBase.
244  */
246 
247  std::shared_ptr<OutputTime> output_stream_;
248 
249 
250 
251 };
252 
253 
254 
255 
256 
257 
258 
259 
260 #endif /* CONC_TRANS_MODEL_HH_ */
SubstanceList substances_
Transported substances.
MultiField< 3, FieldValue< 3 >::Scalar > output_field
MultiField< 3, FieldValue< 3 >::Scalar > init_conc
Initial concentrations.
void set_output_stream(std::shared_ptr< OutputTime > stream)
Setter for output stream.
MultiField< 3, FieldValue< 3 >::Scalar > sorption_coefficient
Coefficient of linear sorption.
Abstract linear system class.
Definition: balance.hh:37
Class template representing a field with values dependent on: point, element, and region...
Definition: field.hh:83
SubstanceList & substances() override
Returns reference to the vector of substance names.
Definition: mesh.h:76
ConcentrationTransportBase FactoryBaseType
bool flux_changed
Indicator of change in advection vector field.
static constexpr const char * name()
void set_velocity_field(const MH_DofHandler &dh) override
Updates the velocity field which determines some coefficients of the transport equation.
double solvent_density_
Density of liquid (a global constant).
Field< 3, FieldValue< 3 >::Scalar > rock_density
Rock matrix density.
MultiField< 3, FieldValue< 3 >::Scalar > disp_t
Transversal dispersivity (for each substance).
BCMultiField< 3, FieldValue< 3 >::Enum > bc_type
Type of boundary condition (see also BC_Type)
BCMultiField< 3, FieldValue< 3 >::Scalar > bc_dirichlet_value
Prescribed concentration for Dirichlet/reference concentration for flux b.c.
Accessor to the data with type Type::Record.
Definition: accessors.hh:292
unsigned int n_substances() override
Returns number of transported substances.
The class for outputting data during time.
Definition: output_time.hh:50
const vector< unsigned int > & get_subst_idx() override
Return substance indices used in balance.
BCMultiField< 3, FieldValue< 3 >::Scalar > bc_robin_sigma
Transition coefficient in total/diffusive flux b.c.
std::shared_ptr< OutputTime > output_stream() override
Getter for output stream.
MultiField< 3, FieldValue< 3 >::Scalar > disp_l
Longitudal dispersivity (for each substance).
Discontinuous Galerkin method for equation of transport with dispersion.
std::shared_ptr< OutputTime > output_stream_
BCMultiField< 3, FieldValue< 3 >::Scalar > bc_flux
Flux value in total/diffusive flux b.c.
Classes for storing substance data.
void set_target_time(double) override
MultiField< 3, FieldValue< 3 >::TensorFixed > diff_m
Molecular diffusivity (for each substance).
vector< unsigned int > subst_idx
List of indices used to call balance methods for a set of quantities.
Record type proxy class.
Definition: type_record.hh:182
Class for representation of a vector of fields of the same physical quantity.
Definition: multi_field.hh:87
Template for classes storing finite set of named values.