37 return Selection(
"ConvectionDiffusion_BC_Type",
"Types of boundary conditions for 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'.")
70 "Type of boundary condition.")
78 .
description(
"Dirichlet boundary condition (for each substance).")
84 .description(
"Flux in Neumann boundary condition.")
85 .units(
UnitSI().kg().m().s(-1).md() )
90 .
name(
"bc_robin_sigma")
91 .description(
"Conductivity coefficient in Robin boundary condition.")
92 .units(
UnitSI().m(4).s(-1).md() )
102 .
description(
"Longitudal dispersivity (for each substance).")
108 .
description(
"Transversal dispersivity (for each substance).")
114 .
description(
"Molecular diffusivity (for each substance).")
131 description +
" for solute transport.")
141 "Selection of output fields for Diffusive Solute Transport DG model.");
162 vector<double> elem_csec(point_list.size()), por_m(point_list.size());
167 for (
unsigned int i=0; i<point_list.size(); i++)
168 mm_coef[i] = elem_csec[i]*por_m[i];
173 double Dm,
double alphaL,
double alphaT,
double porosity,
double cross_cut,
arma::mat33 &K)
175 double vnorm = arma::norm(velocity, 2);
178 for (
int i=0; i<3; i++)
179 for (
int j=0; j<3; j++)
180 K(i,j) = (velocity[i]/vnorm)*(velocity[j]/vnorm)*(alphaL-alphaT) + alphaT*(i==j?1:0);
186 K = (vnorm*K + (Dm*pow(porosity, 1./3)*porosity*cross_cut)*arma::eye(3,3));
198 const unsigned int qsize = point_list.size();
199 const unsigned int n_subst = dif_coef.size();
200 std::vector<arma::vec> Dm(qsize, arma::vec(n_subst) ), alphaL(qsize, arma::vec(n_subst) ), alphaT(qsize, arma::vec(n_subst) );
209 for (
unsigned int i=0; i<qsize; i++) {
210 for (
unsigned int sbi=0; sbi<n_subst; sbi++) {
211 ad_coef[sbi][i] = velocity[i];
213 Dm[i][sbi], alphaL[i][sbi], alphaT[i][sbi], por_m[i], csection[i],
228 arma::uvec &bc_types)
248 for (
auto f : bc_flux) f = -f;
266 const unsigned int qsize = point_list.size();
273 for (
unsigned int k=0; k<qsize; k++)
275 sources_density[k] *= csection[k];
276 sources_sigma[k] *= csection[k];
285 const unsigned int qsize = point_list.size();
290 for (
unsigned int k=0; k<qsize; k++)
291 sources_sigma[k] *= csection[k];
void get_bc_type(const ElementAccessor< 3 > &ele_acc, arma::uvec &bc_types) override
SubstanceList substances_
Transported substances.
static const Input::Type::Selection & get_bc_type_selection()
ConcentrationTransportModel(Mesh &mesh, const Input::Record &in_rec)
MultiField< 3, FieldValue< 3 >::Scalar > output_field
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()
FieldCommon & flags_add(FieldFlag::Flags::Mask mask)
SubstanceList & substances() override
Returns reference to the vector of substance names.
MultiField< 3, FieldValue< 3 >::Scalar > diff_m
Molecular diffusivity (for each substance).
void initialize(const Input::Array &in_array)
Read from input array.
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.
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
void set_components(SubstanceList &substances, const Input::Record &in_rec) override
Read or set names of solution components.
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.
void compute_init_cond(const std::vector< arma::vec3 > &point_list, const ElementAccessor< 3 > &ele_acc, std::vector< arma::vec > &init_values) override
void compute_sources_sigma(const std::vector< arma::vec3 > &point_list, const ElementAccessor< 3 > &ele_acc, std::vector< arma::vec > &sources_sigma) override
static constexpr Mask equation_result
Match result fields. These are never given by input or copy of input.
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) ...
FieldCommon & input_selection(const Input::Type::Selection *element_selection)
void calculate_dispersivity_tensor(const arma::vec3 &velocity, double Dm, double alphaL, double alphaT, double porosity, double cross_cut, arma::mat33 &K)
Field< 3, FieldValue< 3 >::Scalar > porosity
Mobile porosity.
static IT::Selection get_output_selection()
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.
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) override
const MH_DofHandler * mh_dh
void set_balance_object(boost::shared_ptr< Balance > balance) override
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.
virtual MultiFieldValue::return_type value(const Point &p, const ElementAccessor< spacedim > &elm) const
BCMultiField< 3, FieldValue< 3 >::Scalar > bc_flux
Flux value in total/diffusive flux b.c.
FieldCommon & name(const string &name)
vector< unsigned int > subst_idx
List of indices used to call balance methods for a set of quantities.
arma::vec::fixed< spacedim > centre() const
boost::shared_ptr< Balance > balance_
object for calculation and writing the mass balance to file.
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.
virtual void value_list(const std::vector< Point > &point_list, const ElementAccessor< spacedim > &elm, std::vector< typename MultiFieldValue::return_type > &value_list) const
~ConcentrationTransportModel() override