Page MenuHomec4science

equilibrium_stagger.cpp
No OneTemporary

File Metadata

Created
Tue, May 14, 13:57

equilibrium_stagger.cpp

/* =============================================================================
Copyright (c) 2014-2017 F. Georget <fabieng@princeton.edu> Princeton University
Copyright (c) 2017 F. Georget <fabien.georget@epfl.ch> EPFL
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 "equilibrium_stagger.hpp"
#include "specmicp_common/cached_vector.hpp"
#include "specmicp/adimensional/adimensional_system_solution.hpp"
#include "specmicp/adimensional/adimensional_system_structs.hpp"
#include "specmicp/adimensional/adimensional_system_solver_structs.hpp"
#include "specmicp/adimensional/adimensional_system_solver.hpp"
#include "reactmicp/solver/staggers_base/stagger_structs.hpp"
#include "variables.hpp"
#include "boundary_conditions.hpp"
#include "reactmicp/systems/common/equilibrium_options_vector.hpp"
#include "specmicp_common/log.hpp"
#include <iostream>
#include "specmicp/io/print.hpp"
namespace specmicp {
namespace reactmicp {
namespace systems {
namespace chloride {
struct ChlorideEquilibriumStagger::ChlorideEquilibriumStaggerImpl {
std::shared_ptr<EquilibriumOptionsVector> m_options;
std::shared_ptr<BoundaryConditions> m_bcs;
RawDatabasePtr m_data;
scalar_t m_dt;
ChlorideEquilibriumStaggerImpl(
std::shared_ptr<ChlorideVariables> var,
std::shared_ptr<BoundaryConditions> bcs,
std::shared_ptr<EquilibriumOptionsVector> opts
):
m_options(opts),
m_bcs(bcs),
m_data(var->get_database())
{}
void initialize(ChlorideVariables * const var);
void initialize_timestep(scalar_t dt, ChlorideVariables * const var);
AdimensionalSystemConstraints pre_filled_constraints(
index_t node,
ChlorideVariables * const vars
);
micpsolver::MiCPSolverReturnCode solve_one_node(
index_t node,
const AdimensionalSystemConstraints& constraints,
const AdimensionalSystemSolution& previous,
AdimensionalSystemSolution& output
);
AdimensionalSystemConstraints pre_filled_constraints_perturb(
index_t node,
ChlorideVariables * const vars,
const Eigen::Ref<const Vector>& nodal_displacement
);
int compute_one_node(
index_t node,
ChlorideVariables* const vars,
AdimensionalSystemSolution& out
);
void analyse_solution(
index_t node,
const AdimensionalSystemSolution& solution,
ChlorideVariables* const vars
);
ChlorideVariables * const cast_to_true_var(VariablesBase * const var) {
return static_cast<ChlorideVariables * const>(var);
}
};
// Main interface
// ==============
ChlorideEquilibriumStagger::ChlorideEquilibriumStagger(
std::shared_ptr<ChlorideVariables> var,
std::shared_ptr<BoundaryConditions> bcs,
std::shared_ptr<EquilibriumOptionsVector> opts
):
m_impl(utils::make_pimpl<ChlorideEquilibriumStaggerImpl>(var, bcs, opts))
{
}
ChlorideEquilibriumStagger::~ChlorideEquilibriumStagger() = default;
void ChlorideEquilibriumStagger::initialize(VariablesBase * const var)
{
return m_impl->initialize(m_impl->cast_to_true_var(var));
}
void ChlorideEquilibriumStagger::initialize_timestep(scalar_t dt, VariablesBase * const var)
{
return m_impl->initialize_timestep(dt, m_impl->cast_to_true_var(var));
}
solver::StaggerReturnCode ChlorideEquilibriumStagger::restart_timestep(VariablesBase * const var)
{
std::cout << "Restart eq timestep" << std::endl;
ChlorideVariables * const true_var = m_impl->cast_to_true_var(var);
int sum_retcode = 0;
#ifdef SPECMICP_HAVE_OPENMP
#pragma omp parallel for\
schedule(dynamic, 5)\
reduction(+: sum_retcode)
#endif // SPECMICP_HAVE_OPENMP
for (index_t n=0; n<true_var->get_mesh()->nb_nodes();++n)
{
AdimensionalSystemSolution sol;
if (static_cast<int>(solve_equilibrium_at_node(n, true_var, sol)) <= 0) {
sum_retcode += -1;
}
}
if (sum_retcode != 0)
{
return solver::StaggerReturnCode::UnknownError;
}
return solver::StaggerReturnCode::ResidualMinimized;
}
ReturnCode ChlorideEquilibriumStagger::solve_equilibrium_at_node(
index_t node,
solver::VariablesBase* const var,
AdimensionalSystemSolution& out)
{
ChlorideVariables * const vars = m_impl->cast_to_true_var(var);
AdimensionalSystemConstraints constraints = m_impl->pre_filled_constraints(node, vars);
micpsolver::MiCPSolverReturnCode retcode = m_impl->solve_one_node(node, constraints, vars->adim_solution(node), out);
if (retcode <= micpsolver::MiCPSolverReturnCode::Success) {
io::print_specmicp_constraints(m_impl->m_data, constraints);
scalar_t sum = 0.0;
for (auto idc: m_impl->m_data->range_aqueous_component()) {
sum += m_impl->m_data->charge_component(idc)*constraints.total_concentrations(idc);
}
std::cout << "Charge balance : " << sum << "\n";
ERROR << "Failed to solve equilibrium at node : " << node;
} else {
m_impl->analyse_solution(node, out, vars);
}
return micpsolver::to_generic_return_code(retcode);
}
// Implementation
// ==============
void ChlorideEquilibriumStagger::ChlorideEquilibriumStaggerImpl::initialize(ChlorideVariables * const var)
{
}
void ChlorideEquilibriumStagger::ChlorideEquilibriumStaggerImpl::initialize_timestep(scalar_t dt, ChlorideVariables * const var)
{
m_dt = dt;
for (auto node: var->get_mesh()->range_nodes()) {
auto sol = var->adim_solution(node);
if (not sol.is_valid) {
throw std::runtime_error("Invalid adim solution for node "
+std::to_string(node) + ". Whyyyyyyyy ?");
}
}
}
ReturnCode ChlorideEquilibriumStagger::perturb_adim_solution(
index_t node,
ChlorideVariables * const var,
const Eigen::Ref<const Vector>& nodal_displacement,
const AdimensionalSystemSolution& previous,
AdimensionalSystemSolution& out
)
{
AdimensionalSystemConstraints constraints = m_impl->pre_filled_constraints_perturb(
node, var, nodal_displacement);
//AdimensionalSystemSolutionExtractor extr(var->adim_solution(node), var->get_database(), var->get_units());
//std::cout << extr.volume_fraction_mineral(var->get_database()->get_id_mineral("Portlandite")) << std::endl;
//if (node == 1) {
//std::cout << node << "\n";
//io::print_specmicp_constraints(var->get_database(), constraints);
//}
AdimensionalSystemSolverOptions opts = (*m_impl->m_options)[node];
opts.system_options.non_ideality = false;
//opts.solver_options.set_maximum_step_length(0.1);
//opts.solver_options.set_tolerance(1e-2);
AdimensionalSystemSolver solver(m_impl->m_data, constraints, previous, opts);
auto x = previous.main_variables;
auto perf = solver.solve(x, false);
if (perf.return_code > micpsolver::MiCPSolverReturnCode::NotConvergedYet) {
out = solver.merge_raw_solution_with_immobile_species(x, previous);
}
else if (perf.return_code <= micpsolver::MiCPSolverReturnCode::Success) {
io::print_specmicp_constraints(m_impl->m_data, constraints);
io::print_specmicp_options(std::cout, opts);
ERROR << "Failed to solve equilibrium at node : " << node;
}
return micpsolver::to_generic_return_code(perf.return_code);
}
// ###
int ChlorideEquilibriumStagger::ChlorideEquilibriumStaggerImpl::compute_one_node(
index_t node,
ChlorideVariables* const vars,
AdimensionalSystemSolution& out
)
{
AdimensionalSystemConstraints constraints = pre_filled_constraints(node, vars);
micpsolver::MiCPSolverReturnCode retcode = solve_one_node(node, constraints, vars->adim_solution(node), out);
if (retcode < micpsolver::MiCPSolverReturnCode::Success) {
return -1;
}
analyse_solution(node, out, vars);
return 0;
}
void ChlorideEquilibriumStagger::ChlorideEquilibriumStaggerImpl::analyse_solution(
index_t node,
const AdimensionalSystemSolution& solution,
ChlorideVariables* const vars
)
{
vars->adim_solution(node) = solution;
AdimensionalSystemSolutionExtractor extr(solution, m_data, vars->get_units() );
if (node == 1) {
std::cout << extr.volume_fraction_mineral(m_data->get_id_mineral("Friedels_salt"))
<< "\n";
}
for (auto i: m_data->range_aqueous_component())
{
auto dof = vars->dof_component(node, i);
auto conc_aq = extr.total_aqueous_concentration(i)*extr.density_water();
auto conc_aq_vel = (conc_aq - vars->predictors()(dof))/m_dt;
vars->velocities()(dof) = conc_aq_vel;
vars->main_variables()(dof) = conc_aq;
auto conc_immobile = extr.total_immobile_concentration(i);
auto conc_immobile_vel = (conc_immobile - vars->predictor_immobile_concentrations()(dof))/m_dt;
vars->total_immobile_concentrations()(dof) = conc_immobile;
vars->velocity_immobile_concentrations()(dof) = conc_immobile_vel;
vars->chemistry_rates()(dof) = -conc_immobile_vel;
}
{
auto porosity = extr.porosity();
vars->velocities_porosity()(node) = (porosity - vars->predictors_porosity()(node))/m_dt;
vars->porosity()(node) = porosity;
}
}
micpsolver::MiCPSolverReturnCode ChlorideEquilibriumStagger::ChlorideEquilibriumStaggerImpl::solve_one_node(
index_t node,
const AdimensionalSystemConstraints& constraints,
const AdimensionalSystemSolution& previous,
AdimensionalSystemSolution& output
)
{
auto& opts = (*m_options)[node];
if (opts.system_options.solve_solid_solutions == false) {
std::cout << "Without solid solution \n";
AdimensionalSystemSolver solver(m_data, constraints, previous, opts);
output.main_variables = previous.main_variables;
auto perf = solver.solve(output.main_variables, false);
if (perf.return_code > micpsolver::MiCPSolverReturnCode::NotConvergedYet) {
output = solver.get_raw_solution(output.main_variables);
}
return perf.return_code;
} else {
std::cout << "With solid solution \n";
opts.system_options.solve_solid_solutions = false;
AdimensionalSystemSolver solver(m_data, constraints, opts);
output.main_variables = previous.main_variables;
auto perf = solver.solve(output.main_variables, false);
if (perf.return_code > micpsolver::MiCPSolverReturnCode::NotConvergedYet) {
std::cout << "Solve no solid solution \n";
output = solver.get_raw_solution(output.main_variables);
opts.system_options.solve_solid_solutions = true;
AdimensionalSystemSolver solver(m_data, constraints, output, opts);
perf = solver.solve(output.main_variables, false);
if (perf.return_code > micpsolver::MiCPSolverReturnCode::NotConvergedYet) {
std::cout << "Solve solid solution \n";
output = solver.get_raw_solution(output.main_variables);
}
}
return perf.return_code;
}
}
// ### implementation
AdimensionalSystemConstraints
ChlorideEquilibriumStagger::ChlorideEquilibriumStaggerImpl::pre_filled_constraints(
index_t node,
ChlorideVariables * const vars
)
{
AdimensionalSystemConstraints constraints = m_bcs->get_constraint(node);
constraints.inert_volume_fraction = vars->inert_volume_fraction(node);
Vector tot_concs = Vector::Zero(m_data->nb_component());
for (auto aqc: m_data->range_aqueous_component())
{
tot_concs(aqc) = vars->porosity(node)*vars->total_mobile_concentration(node, aqc)+vars->total_immobile_concentration(node, aqc);
}
constraints.set_total_concentrations(tot_concs);
return constraints;
}
AdimensionalSystemConstraints
ChlorideEquilibriumStagger::ChlorideEquilibriumStaggerImpl::pre_filled_constraints_perturb(
index_t node,
ChlorideVariables * const vars,
const Eigen::Ref<const Vector>& nodal_displacement
)
{
AdimensionalSystemConstraints constraints = m_bcs->get_constraint(node);
//constraints.inert_volume_fraction = vars->inert_volume_fraction(node);
Vector tot_concs = Vector::Zero(m_data->nb_component());
for (auto aqc: m_data->range_aqueous_component())
{
const index_t ndof = vars->ndof_component(aqc);
//tot_concs(aqc) = vars->porosity(node)*nodal_displacement(ndof);
tot_concs(aqc) = nodal_displacement(ndof);
}
constraints.set_total_concentrations(tot_concs);
constraints.disable_immobile_species();
return constraints;
}
} // namespace chloride
} // namespace systems
} // namespace reactmicp
} // namespace specmicp

Event Timeline