Page MenuHomec4science

saturation_equation.cpp
No OneTemporary

File Metadata

Created
Sat, Jul 6, 14:21

saturation_equation.cpp

/*-------------------------------------------------------------------------------
Copyright (c) 2014,2015 F. Georget <fabieng@princeton.edu>, Princeton University
All rights reserved.
Redistribution and use in source and binary forms, with or without modification,
are permitted provided that the following conditions are met:
1. Redistributions of source code must retain the above copyright notice, this
list of conditions and the following disclaimer.
2. Redistributions in binary form must reproduce the above copyright notice,
this list of conditions and the following disclaimer in the documentation and/or
other materials provided with the distribution.
3. Neither the name of the copyright holder nor the names of its contributors
may be used to endorse or promote products derived from this software without
specific prior written permission.
THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS "AS IS" AND
ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE IMPLIED
WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE ARE
DISCLAIMED. IN NO EVENT SHALL THE COPYRIGHT HOLDER OR CONTRIBUTORS BE LIABLE FOR
ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL, EXEMPLARY, OR CONSEQUENTIAL DAMAGES
(INCLUDING, BUT NOT LIMITED TO, PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES;
LOSS OF USE, DATA, OR PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON
ANY THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT
(INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
-----------------------------------------------------------------------------*/
#include "saturation_equation.hpp"
#include "variables_box.hpp"
#include "transport_constraints.hpp"
#include "../../../dfpm/meshes/mesh1d.hpp"
#include "../../../dfpmsolver/parabolic_driver.hpp"
#include "../../../utils/compat.hpp"
#include <vector>
namespace specmicp {
namespace dfpmsolver {
// explicit template instanciation
template class dfpmsolver::ParabolicDriver<reactmicp::systems::unsaturated::SaturationEquation>;
} //end namespace dfpmsolver
} //end namespace specmicp
namespace specmicp {
namespace reactmicp {
namespace systems {
namespace unsaturated {
static constexpr index_t no_equation {-1};
static constexpr index_t no_eq_no_var {-2};
static constexpr index_t not_initialized {-5};
struct SPECMICP_DLL_LOCAL SaturationEquation::SaturationEquationImpl
{
mesh::Mesh1DPtr m_mesh;
SaturationVariableBox m_vars;
std::vector<index_t> m_ideq;
bool m_store_residual_info;
scalar_t m_scaling {1.0};
SaturationEquationImpl(
mesh::Mesh1DPtr the_mesh,
SaturationVariableBox vars):
m_mesh(the_mesh),
m_vars(vars),
m_ideq(the_mesh->nb_nodes(), not_initialized)
{}
index_t& id_equation(index_t node) {return m_ideq[node];}
bool node_can_flux(index_t node) {return m_ideq[node] != no_eq_no_var;}
bool node_has_eq(index_t node) {return m_ideq[node] > no_equation;}
void set_store_residual_info() {
m_store_residual_info = true;
}
void reset_store_residual_info() {
m_store_residual_info = false;
}
bool store_residual_info() {
return m_store_residual_info;
}
//! \brief Return a pointer to the mesh
mesh::Mesh1D* mesh() {return m_mesh.get();}
range_t range_nodes() {return m_mesh->range_nodes();}
void set_relative_variables(const Vector& displacement);
void set_relative_variables(index_t node, const Vector& displacement);
void compute_transport_rate(scalar_t dt, const Vector& displacement);
};
SaturationEquation::~SaturationEquation() = default;
SaturationEquation::SaturationEquation(mesh::Mesh1DPtr the_mesh,
SaturationVariableBox &variables,
const TransportConstraints& constraints
):
m_neq(-1),
m_tot_ndf(the_mesh->nb_nodes()),
m_impl(make_unique<SaturationEquationImpl>(the_mesh, variables))
{
number_equations(constraints);
}
void SaturationEquation::number_equations(const TransportConstraints& constraints)
{
for (index_t node: constraints.fixed_nodes())
{
m_impl->id_equation(node) = no_equation;
}
for (index_t node: constraints.gas_nodes())
{
m_impl->id_equation(node) = no_eq_no_var;
}
index_t neq = 0;
for (index_t node: m_impl->range_nodes()) {
if (m_impl->id_equation(node) == not_initialized)
{
m_impl->id_equation(node) = neq;
++neq;
}
}
m_neq = neq;
}
index_t SaturationEquation::id_equation(index_t id_dof)
{
return m_impl->id_equation(id_dof)<0?no_equation:m_impl->id_equation(id_dof);
}
void SaturationEquation::set_scaling(scalar_t value)
{
m_impl->m_scaling = value;
}
// Residuals
// =========
void SaturationEquation::residuals_element(
index_t element,
const Vector& displacement,
const Vector& velocity,
Eigen::Vector2d& element_residual,
bool use_chemistry_rate
)
{
element_residual.setZero();
mesh::Mesh1D* m_mesh = m_impl->mesh();
SaturationVariableBox& vars = m_impl->m_vars;
const scalar_t mass_coeff_0 = m_mesh->get_volume_cell_element(element, 0);
const scalar_t mass_coeff_1 = m_mesh->get_volume_cell_element(element, 1);
const index_t node_0 = m_mesh->get_node(element, 0);
const index_t node_1 = m_mesh->get_node(element, 1);
scalar_t flux_0 = 0.0;
scalar_t flux_1 = 0.0;
if (m_impl->node_can_flux(node_0) and m_impl->node_can_flux(node_1))
{
// Cap pressure gradient
const scalar_t perm_0 = vars.liquid_permeability(node_0)
* vars.relative_liquid_permeability(node_0);
const scalar_t perm_1 = vars.liquid_permeability(node_1)
* vars.relative_liquid_permeability(node_1);
const scalar_t perm = 2.0/(1.0/perm_0 + 1.0/perm_1);
const scalar_t aq_coefficient = (vars.aqueous_concentration(node_0)
+ vars.aqueous_concentration(node_1)
) / 2.0;
const scalar_t cap_pressure_gradient = ( vars.capillary_pressure(node_1)
- vars.capillary_pressure(node_0)
) / m_mesh->get_dx(element);
const scalar_t advection_flux = (perm/vars.constants.viscosity_liquid_water)*cap_pressure_gradient;
const scalar_t cappres_flux = -aq_coefficient*advection_flux;
// Diffusion Cw
const scalar_t coeff_diff_0 = vars.liquid_diffusivity(node_0)
* vars.relative_liquid_diffusivity(node_0);
const scalar_t coeff_diff_1 = vars.liquid_diffusivity(node_1)
* vars.relative_liquid_diffusivity(node_1);
const scalar_t coeff_diff = 2.0/(1.0/coeff_diff_0 + 1.0/coeff_diff_1);
const scalar_t aq_flux = coeff_diff*(vars.aqueous_concentration(node_1)
- vars.aqueous_concentration(node_0)
) / m_mesh->get_dx(element);
// Tot flux
const scalar_t tot_flux = m_mesh->get_face_area(element)*(cappres_flux + aq_flux);
flux_0 = tot_flux;
flux_1 = -tot_flux;
// Storage
if (m_impl->store_residual_info())
{
// advective flux stored by element
vars.advection_flux(element) = advection_flux;
// fluxes to compute exchange term
vars.liquid_saturation.transport_fluxes(node_0) += flux_0;
vars.liquid_saturation.transport_fluxes(node_1) += flux_1;
}
}
// transient
if (m_impl->node_has_eq(node_0))
{
const scalar_t porosity_0 = vars.porosity(node_0);
const scalar_t aq_tot_conc_0 = vars.aqueous_concentration(node_0);
const scalar_t saturation_0 = displacement(node_0);
scalar_t transient_0 = (
porosity_0 * aq_tot_conc_0 * velocity(node_0)
+ saturation_0 * aq_tot_conc_0 * vars.porosity.velocity(node_0)
+ porosity_0 * saturation_0 * vars.aqueous_concentration.velocity(node_0)
);
auto res = mass_coeff_0*transient_0 - flux_0;
if (use_chemistry_rate)
{
const scalar_t chemistry_0 = vars.liquid_saturation.chemistry_rate(node_0)
+ vars.solid_concentration.chemistry_rate(node_0)
+ vars.vapor_pressure.chemistry_rate(node_0)
;
res -= mass_coeff_0*chemistry_0;
}
element_residual(0) = res / m_impl->m_scaling;
}
if (m_impl->node_has_eq(node_1))
{
const scalar_t porosity_1 = vars.porosity(node_1);
const scalar_t aq_tot_conc_1 = vars.aqueous_concentration(node_1);
const scalar_t saturation_1 = displacement(node_1);
scalar_t transient_1 = (
porosity_1 * aq_tot_conc_1 * velocity(node_1)
+ saturation_1 * aq_tot_conc_1 * vars.porosity.velocity(node_1)
+ porosity_1 * saturation_1 * vars.aqueous_concentration.velocity(node_1)
);
auto res = mass_coeff_1*transient_1 - flux_1;
if (use_chemistry_rate)
{
const scalar_t chemistry_1 = vars.liquid_saturation.chemistry_rate(node_1)
+ vars.solid_concentration.chemistry_rate(node_1)
+ vars.vapor_pressure.chemistry_rate(node_1)
;
res -= mass_coeff_1*chemistry_1;
}
element_residual(1) = res / m_impl->m_scaling;
}
}
//! \brief Compute the residuals
void SaturationEquation::compute_residuals(
const Vector& displacement,
const Vector& velocity,
Vector& residuals,
bool use_chemistry_rate
)
{
mesh::Mesh1D* m_mesh = m_impl->mesh();
residuals.setZero(get_neq());
m_impl->set_store_residual_info();
set_relative_variables(displacement);
Eigen::Vector2d element_residual;
for (index_t element: m_mesh->range_elements())
{
residuals_element(element, displacement, velocity, element_residual, use_chemistry_rate);
const index_t node_0 = m_mesh->get_node(element, 0);
if (m_impl->node_has_eq(node_0))
{
residuals(m_impl->id_equation(node_0)) += element_residual(0);
}
const index_t node_1 = m_mesh->get_node(element, 1);
if (m_impl->node_has_eq(node_1))
{
residuals(m_impl->id_equation(node_1)) += element_residual(1);
}
}
m_impl->reset_store_residual_info();
}
void SaturationEquation::compute_jacobian(
Vector& displacement,
Vector& velocity,
Eigen::SparseMatrix<scalar_t>& jacobian,
scalar_t alphadt
)
{
mesh::Mesh1D* m_mesh = m_impl->mesh();
dfpm::list_triplet_t jacob;
const index_t estimation = 3*get_neq();
jacob.reserve(estimation);
// assume relative variables are set
for (index_t element: m_mesh->range_elements())
{
Eigen::Vector2d element_residual_orig;
residuals_element(element, displacement, velocity, element_residual_orig, false);
for (index_t enodec=0; enodec<2; ++enodec)
{
const index_t nodec = m_mesh->get_node(element, enodec);
if (not m_impl->node_has_eq(nodec)) continue;
const scalar_t tmp_d = displacement(nodec);
const scalar_t tmp_v = velocity(nodec);
scalar_t h = eps_jacobian*std::abs(tmp_v);
if (h<1e-4*eps_jacobian) h = eps_jacobian;
velocity(nodec) = tmp_v + h;
h = velocity(nodec) - tmp_v;
displacement(nodec) = tmp_d + alphadt*h;
m_impl->set_relative_variables(nodec, displacement);
Eigen::Vector2d element_residual;
residuals_element(element, displacement, velocity, element_residual, false);
displacement(nodec) = tmp_d;
velocity(nodec) = tmp_v;
m_impl->set_relative_variables(nodec, displacement);
for (index_t enoder=0; enoder<2; ++enoder)
{
const index_t noder = m_mesh->get_node(element, enoder);
if (not m_impl->node_has_eq(noder)) continue;
jacob.push_back(dfpm::triplet_t(
m_impl->id_equation(noder),
m_impl->id_equation(nodec),
(element_residual(enoder) - element_residual_orig(enoder))/h
));
}
}
}
jacobian = Eigen::SparseMatrix<scalar_t>(get_neq(), get_neq());
jacobian.setFromTriplets(jacob.begin(), jacob.end());
}
void SaturationEquation::update_solution(
const Vector& update,
scalar_t lambda,
scalar_t alpha_dt,
Vector& predictor,
Vector& displacement,
Vector& velocity
)
{
for (index_t node: m_impl->mesh()->range_nodes())
{
if (m_impl->node_has_eq(node))
{
velocity(node) += lambda*update(m_impl->id_equation(node));
}
}
displacement = predictor + alpha_dt*velocity;
m_impl->compute_transport_rate(alpha_dt, displacement);
}
void SaturationEquation::set_relative_variables(const Vector& displacement)
{
return m_impl->set_relative_variables(displacement);
}
void SaturationEquation::SaturationEquationImpl::set_relative_variables(
index_t node,
const Vector& displacement
)
{
if (not node_can_flux(node)) return;
const scalar_t saturation = displacement(node);
m_vars.relative_liquid_diffusivity(node) = m_vars.relative_liquid_diffusivity_f(node, saturation);
m_vars.relative_liquid_permeability(node) = m_vars.relative_liquid_permeability_f(node, saturation);
m_vars.capillary_pressure(node) = m_vars.capillary_pressure_f(node, saturation);
}
void SaturationEquation::SaturationEquationImpl::set_relative_variables(const Vector& displacement)
{
for (index_t node: m_mesh->range_nodes())
{
set_relative_variables(node, displacement);
}
}
void SaturationEquation::SaturationEquationImpl::compute_transport_rate(
scalar_t dt,
const Vector& displacement)
{
MainVariable& saturation = m_vars.liquid_saturation;
const MainVariable& solid_conc = m_vars.solid_concentration;
const MainVariable& pressure = m_vars.vapor_pressure;
const SecondaryTransientVariable& porosity = m_vars.porosity;
const SecondaryTransientVariable& aqueous_concentration = m_vars.aqueous_concentration;
for (index_t node: m_mesh->range_nodes())
{
if (! node_has_eq(node)) continue;
const scalar_t transient = (
(porosity(node)*aqueous_concentration(node)*displacement(node))
- (porosity.predictor(node)*aqueous_concentration.predictor(node)*saturation.predictor(node))
) / dt;
const scalar_t chem_rates = (
saturation.chemistry_rate(node)
+ solid_conc.chemistry_rate(node)
+ pressure.chemistry_rate(node)
);
saturation.transport_fluxes(node) = transient - chem_rates;
}
}
} //end namespace unsaturated
} //end namespace systems
} //end namespace reactmicp
} //end namespace specmicp

Event Timeline