37 return Selection(
"Solute_AdvectionDiffusion_BC_Type",
"Types of boundary conditions for advection-diffusion solute transport model.")
39 "Default transport boundary condition.\n" 40 "On water inflow (($(q_w \\le 0)$)), total flux is given by the reference concentration 'bc_conc'. " 41 "On water outflow we prescribe zero diffusive flux, " 42 "i.e. the mass flows out only due to advection.")
44 "Dirichlet boundary condition (($ c = c_D $)).\n" 45 "The prescribed concentration (($c_D$)) is specified by the field 'bc_conc'.")
47 "Total mass flux boundary condition.\n" 48 "The prescribed total incoming flux can have the general form (($\\delta(f_N+\\sigma_R(c_R-c) )$)), " 49 "where the absolute flux (($f_N$)) is specified by the field 'bc_flux', " 50 "the transition parameter (($\\sigma_R$)) by 'bc_robin_sigma', " 51 "and the reference concentration (($c_R$)) by 'bc_conc'.")
52 .
add_value(bc_diffusive_flux,
"diffusive_flux",
53 "Diffusive flux boundary condition.\n" 54 "The prescribed incoming mass flux due to diffusion can have the general form (($\\delta(f_N+\\sigma_R(c_R-c) )$)), " 55 "where the absolute flux (($f_N$)) is specified by the field 'bc_flux', " 56 "the transition parameter (($\\sigma_R$)) by 'bc_robin_sigma', " 57 "and the reference concentration (($c_R$)) by 'bc_conc'.")
72 "Type of boundary condition.")
80 .
description(
"Dirichlet boundary condition (for each substance).")
86 .description(
"Flux in Neumann boundary condition.")
87 .units(
UnitSI().kg().m().s(-1).md() )
92 .
name(
"bc_robin_sigma")
93 .description(
"Conductivity coefficient in Robin boundary condition.")
94 .units(
UnitSI().m(4).s(-1).md() )
100 .
description(
"Initial values for concentration of substances.")
104 .
description(
"Longitudinal dispersivity in the liquid (for each substance).")
110 .
description(
"Transverse dispersivity in the liquid (for each substance).")
116 .
description(
"Molecular diffusivity in the liquid (for each substance).")
121 .
name(
"rock_density")
127 .
name(
"sorption_coefficient")
146 description +
" for solute transport.")
149 "Density of the solvent [ (($kg.m^{-3}$)) ].");
158 "Selection of output fields for Diffusive Solute Transport DG model.");
180 vector<double> elem_csec(point_list.size()), wc(point_list.size());
186 for (
unsigned int i=0; i<point_list.size(); i++)
187 mm_coef[i] = elem_csec[i]*wc[i];
196 por_m(point_list.size()),
197 rho_l(point_list.size()),
198 rho_s(point_list.size()),
199 sorp_mult(point_list.size());
209 for (
unsigned int i=0; i<point_list.size(); i++)
211 ret_coef[sbi][i] = (1.-por_m[i])*rho_s[i]*sorp_mult[i]*elem_csec[i];
218 const arma::mat33 &Dm,
double alphaL,
double alphaT,
double water_content,
double porosity,
double cross_cut,
arma::mat33 &K)
220 double vnorm = arma::norm(velocity, 2);
225 double tortuosity = pow(water_content, 7.0 / 3.0)/ (porosity * porosity);
231 if (fabs(vnorm) > 0) {
237 K = ((alphaL - alphaT) / vnorm) * arma::kron(velocity.t(), velocity);
247 K += alphaT*vnorm*arma::eye(3,3) + Dm*(tortuosity*cross_cut*water_content);
260 const unsigned int qsize = point_list.size();
261 const unsigned int n_subst = dif_coef.size();
263 std::vector<double> alphaL(qsize), alphaT(qsize), por_m(qsize), csection(qsize), wc(qsize);
269 for (
unsigned int sbi=0; sbi<n_subst; sbi++)
271 data().
diff_m[sbi].value_list(point_list, ele_acc, Dm);
272 data().
disp_l[sbi].value_list(point_list, ele_acc, alphaL);
273 data().
disp_t[sbi].value_list(point_list, ele_acc, alphaT);
274 for (
unsigned int i=0; i<qsize; i++)
276 ad_coef[sbi][i] = velocity[i];
278 Dm[i], alphaL[i], alphaT[i], wc[i], por_m[i], csection[i],
290 data().
init_conc[sbi].value_list(point_list, ele_acc, init_values[sbi]);
294 arma::uvec &bc_types)
311 data().
bc_flux[index].value_list(point_list, ele_acc, bc_flux);
316 for (
auto f : bc_flux) f = -f;
334 const unsigned int qsize = point_list.size();
339 data().
sources_conc[sbi].value_list(point_list, ele_acc, sources_value[sbi]);
343 for (
unsigned int k=0; k<qsize; k++)
345 sources_density[sbi][k] *= csection[k];
346 sources_sigma[sbi][k] *= csection[k];
356 const unsigned int qsize = point_list.size();
363 for (
unsigned int k=0; k<qsize; k++)
364 sources_sigma[sbi][k] *= csection[k];
void get_bc_type(const ElementAccessor< 3 > &ele_acc, arma::uvec &bc_types) override
void compute_retardation_coefficient(const std::vector< arma::vec3 > &point_list, const ElementAccessor< 3 > &ele_acc, std::vector< std::vector< double > > &ret_coef) override
SubstanceList substances_
Transported substances.
static const Input::Type::Selection & get_bc_type_selection()
void set_balance_object(std::shared_ptr< Balance > balance) override
ConcentrationTransportModel(Mesh &mesh, const Input::Record &in_rec)
MultiField< 3, FieldValue< 3 >::Scalar > output_field
FieldCommon & input_selection(Input::Type::Selection element_selection)
auto disable_where(const MultiField< spacedim, typename FieldValue< spacedim >::Enum > &control_field, const vector< FieldEnum > &value_list) -> MultiField &
static constexpr Mask in_main_matrix
A field is part of main "stiffness matrix" of the equation.
MultiField< 3, FieldValue< 3 >::Scalar > init_conc
Initial concentrations.
const std::vector< std::string > & names()
void init_from_input(const Input::Record &in_rec) override
Read necessary data from input record.
FieldCommon & flags_add(FieldFlag::Flags::Mask mask)
MultiField< 3, FieldValue< 3 >::Scalar > sorption_coefficient
Coefficient of linear sorption.
void calculate_dispersivity_tensor(const arma::vec3 &velocity, const arma::mat33 &Dm, double alphaL, double alphaT, double water_content, double porosity, double cross_cut, arma::mat33 &K)
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) override
Return data for diffusive or total flux b.c.
MultiField< 3, FieldValue< 3 >::Scalar > sources_density
Concentration sources - density of substance source, only positive part is used.
void compute_init_cond(const std::vector< arma::vec3 > &point_list, const ElementAccessor< 3 > &ele_acc, std::vector< std::vector< double > > &init_values) override
MultiField< 3, FieldValue< 3 >::Scalar > sources_conc
bool flux_changed
Indicator of change in advection vector field.
static constexpr const char * name()
void compute_mass_matrix_coefficient(const std::vector< arma::vec3 > &point_list, const ElementAccessor< 3 > &ele_acc, std::vector< double > &mm_coef) override
static Input::Type::Abstract & get_input_type()
Common specification of the input record for secondary equations.
Discontinuous Galerkin method for equation of transport with dispersion.
FieldCommon & units(const UnitSI &units)
Set basic units of the field.
virtual void value_list(const std::vector< Point > &point_list, const ElementAccessor< spacedim > &elm, std::vector< typename Value::return_type > &value_list) const
unsigned int size() const
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).
arma::vec::fixed< spacedim > centre() const
Computes the barycenter.
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.
static constexpr Mask equation_result
Match result fields. These are never given by input or copy of input.
static constexpr Mask in_time_term
A field is part of time term of the equation.
FieldCommon & input_default(const string &input_default)
MultiField< 3, FieldValue< 3 >::Scalar > sources_sigma
Concentration sources - Robin type, in_flux = sources_sigma * (sources_conc - mobile_conc) ...
Field< 3, FieldValue< 3 >::Scalar > porosity
Mobile porosity - usually saturated water content in the case of unsaturated flow model...
std::shared_ptr< Balance > balance_
object for calculation and writing the mass balance to file.
static IT::Selection get_output_selection()
unsigned int n_substances() override
Returns number of transported substances.
static constexpr Mask in_rhs
A field is part of the right hand side of the equation.
BCMultiField< 3, FieldValue< 3 >::Scalar > bc_robin_sigma
Transition coefficient in total/diffusive flux b.c.
Field< 3, FieldValue< 3 >::Scalar > water_content
Water content - result of unsaturated water flow model or porosity.
const MH_DofHandler * mh_dh
FieldCommon & description(const string &description)
MultiField< 3, FieldValue< 3 >::Scalar > disp_l
Longitudal dispersivity (for each substance).
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) override
Return transition coefficient for flux b.c.
BCMultiField< 3, FieldValue< 3 >::Scalar > bc_flux
Flux value in total/diffusive flux b.c.
FieldCommon & name(const string &name)
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.
FieldCommon & flags(FieldFlag::Flags::Mask mask)
Field< 3, FieldValue< 3 >::Scalar > cross_section
Pointer to DarcyFlow field cross_section.
Class for representation SI units of Fields.
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) override
static UnitSI & dimensionless()
Returns dimensionless unit.
virtual ModelEqData & data()=0
Derived class should implement getter for ModelEqData instance.
void compute_sources_sigma(const std::vector< arma::vec3 > &point_list, const ElementAccessor< 3 > &ele_acc, std::vector< std::vector< double > > &sources_sigma) override
~ConcentrationTransportModel() override
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) override