diff --git a/src/specmicp/reduced_system.cpp b/src/specmicp/reduced_system.cpp index e8fbad8..fed4668 100644 --- a/src/specmicp/reduced_system.cpp +++ b/src/specmicp/reduced_system.cpp @@ -1,495 +1,525 @@ /*------------------------------------------------------- - Module : specmicp - File : specmicp.cpp - Author : Fabien Georget Copyright (c) 2014, Fabien Georget , 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: * Redistributions of source code must retain the above copyright notice, this list of conditions and the following disclaimer. * 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. * Neither the name of the Princeton University 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 OWNER 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 #include "reduced_system.hpp" #include "utils/log.hpp" #include "equilibrium_data.hpp" #include "physics/constants.hpp" #include "physics/laws.hpp" namespace specmicp { inline double pow10(double x) { return std::pow(10.0, x); } ReducedSystem::ReducedSystem( RawDatabasePtr ptrdata, const Vector &totconc, bool conservation_water ): m_data(ptrdata), m_tot_conc(totconc), m_secondary_conc(ptrdata->nb_aqueous), m_loggamma(ptrdata->nb_component+ptrdata->nb_aqueous) { if (not m_simulbox.is_fixed_temperature() or m_simulbox.get_temperature() != units::celsius(25)) { throw std::runtime_error("Only a fixed temperature of 25°C is supported"); } number_eq(conservation_water); m_secondary_conc.setZero(); m_loggamma.setZero(); } ReducedSystem::ReducedSystem( RawDatabasePtr ptrdata, const Vector &totconc, const SimulationBox &simul_box, bool conservation_water ): m_simulbox(simul_box), m_data(ptrdata), m_tot_conc(totconc), m_secondary_conc(ptrdata->nb_aqueous), m_loggamma(ptrdata->nb_component+ptrdata->nb_aqueous) { if (not m_simulbox.is_fixed_temperature() or m_simulbox.get_temperature() != units::celsius(25)) { throw std::runtime_error("Only a fixed temperature of 25°C is supported"); } number_eq(conservation_water); m_secondary_conc.setZero(); m_loggamma.setZero(); } ReducedSystem::ReducedSystem( RawDatabasePtr ptrdata, const Vector& totconc, const SimulationBox& simul_box, const EquilibriumState& previous_solution, bool conservation_water ): m_simulbox(simul_box), m_data(ptrdata), m_tot_conc(totconc), m_secondary_conc(previous_solution.molality_secondary()), m_loggamma(previous_solution.activity_coefficient()) { if (not m_simulbox.is_fixed_temperature() or m_simulbox.get_temperature() != units::celsius(25)) { throw std::runtime_error("Only a fixed temperature of 25°C is supported"); } number_eq(conservation_water); } void ReducedSystem::number_eq(bool water_equation) { index_t neq = 0; m_ideq.reserve(m_data->nb_component+m_data->nb_mineral); if (water_equation) { m_ideq.push_back(neq); ++neq; } else { m_ideq.push_back(no_equation); } for (index_t i: m_data->range_aqueous_component()) { // Remove components that does not exist if (m_tot_conc(i) == 0 and i!= 1) { INFO << "Component " << m_data->labels_basis[i] << "has total concentration equal to zero, we desactivate it" << std::endl; m_nonactive_component.push_back(i); m_ideq.push_back(no_equation); continue; } else { m_ideq.push_back(neq); ++neq; } } m_nb_free_vars = neq; for (index_t m: m_data->range_mineral()) { bool can_precipitate = true; // Remove minerals that cannot precipitate for (auto it=m_nonactive_component.begin(); it!=m_nonactive_component.end(); ++it) { if (m_data->nu_mineral(m, *it) != 0) { can_precipitate = false; break; // this is not a mineral that can precipitate } } if (can_precipitate) { m_ideq.push_back(neq); ++neq; } else { m_ideq.push_back(no_equation); } } m_nb_tot_vars = neq; m_nb_compl_vars = m_nb_tot_vars - m_nb_free_vars; // aqueous species m_active_aqueous.reserve(m_data->nb_aqueous); for (index_t j: m_data->range_aqueous()) { bool can_exist = true; for (auto it=m_nonactive_component.begin(); it!=m_nonactive_component.end(); ++it) { if (m_data->nu_aqueous(j,*it) != 0) { can_exist = false; } } m_active_aqueous.push_back(can_exist); } } scalar_t ReducedSystem::residual_water(const Vector& x) const { scalar_t res = m_tot_conc(0) - mass_water(x)/m_data->molar_mass_basis_si(0); for (index_t j: m_data->range_aqueous()) { if (m_data->nu_aqueous(j, 0) == 0) continue; res -= mass_water(x)*m_data->nu_aqueous(j, 0)*m_secondary_conc(j); } for (index_t m: m_data->range_mineral()) { if (ideq_min(m) == no_equation or m_data->nu_mineral(m, 0) == 0) continue; res -= m_data->nu_mineral(m, 0)*x(ideq_min(m)); } return res; } double ReducedSystem::residual_component(const Vector &x, index_t i) const { const index_t idp = ideq_paq(i); scalar_t res = m_tot_conc(i) - mass_water(x)*pow10(x(idp)); for (index_t j: m_data->range_aqueous()) { if (not m_active_aqueous[j]) continue; res -= mass_water(x)*m_data->nu_aqueous(j, i)*m_secondary_conc(j); } for (index_t m: m_data->range_mineral()) { if (ideq_min(m) == no_equation) continue; res -= m_data->nu_mineral(m, i)*x(ideq_min(m)); } return res; } double ReducedSystem::residual_mineral(const Vector& x, index_t m) const { scalar_t res = m_data->logk_mineral(m); for (index_t i: m_data->range_aqueous_component()) { if (m_data->nu_mineral(m, i) != 0) res -= m_data->nu_mineral(m, i)*(x(ideq_paq(i)) + m_loggamma(i)); } return res; } void ReducedSystem::get_residuals(const Vector& x, Vector& residual) { set_secondary_concentration(x); if (ideq_w() != no_equation) residual(ideq_w()) = residual_water(x); for (index_t i: m_data->range_aqueous_component()) { if (ideq_paq(i) != no_equation) residual(ideq_paq(i)) = residual_component(x, i); } for (index_t m: m_data->range_mineral()) { if (ideq_min(m) != no_equation) residual(ideq_min(m)) = residual_mineral(x, m); } } //non-optimized Finite difference, for test only //void ReducedSystem::get_jacobian(Eigen::VectorXd& x, // Eigen::MatrixXd& jacobian) //{ // const int neq = total_variables(); // Eigen::VectorXd res(total_variables()); // Eigen::VectorXd perturbed_res(total_variables()); // get_residuals(x, res); // for (int j=0; jmolar_mass_basis_si(0); for (index_t j: m_data->range_aqueous()) { if ( m_data->nu_aqueous(j, 0) == 0 ) continue; tmp -= m_data->nu_aqueous(j, 0)*m_secondary_conc(j); } jacobian(0, 0) = tmp; for (index_t k: m_data->range_aqueous_component()) { if (ideq_paq(k) == no_equation) continue; scalar_t tmp = 0; for (index_t j: m_data->range_aqueous()) { tmp -= m_data->nu_aqueous(j,0)*m_data->nu_aqueous(j,k)*m_secondary_conc(j); } jacobian(0, ideq_paq(k)) = x(ideq_w())*log10 * tmp; } for (index_t m: m_data->range_mineral()) { if (ideq_min(m) == no_equation) continue; jacobian(0, ideq_min(m)) = -m_data->nu_mineral(m, 0); } } // aqueous component for (index_t i: m_data->range_aqueous_component()) { const index_t idp = ideq_paq(i); if (idp == no_equation) continue; // aqueous components for (index_t k: m_data->range_aqueous_component()) { if (ideq_paq(k) == no_equation) continue; scalar_t tmp_iip = 0; if (k == i) tmp_iip -= pow10(x(ideq_paq(i)))*log10; for (index_t j: m_data->range_aqueous()) { tmp_iip -= log10*m_data->nu_aqueous(j, i)*m_data->nu_aqueous(j, k)*m_secondary_conc(j); } jacobian(idp, ideq_paq(k)) = mass_water(x)*tmp_iip; } // minerals for (index_t m: m_data->range_mineral()) { if (ideq_min(m) == no_equation) continue; jacobian(idp, ideq_min(m)) = - m_data->nu_mineral(m, i); } if (ideq_w() != no_equation) // Water { scalar_t tmp_iw = -std::pow(10.0, x(idp)); for (index_t j: m_data->range_aqueous()) { if ( m_data->nu_aqueous(j, i) == 0 ) continue; tmp_iw -= m_data->nu_aqueous(j, i)*m_secondary_conc(j); } jacobian(idp, ideq_w()) = tmp_iw; } } // mineral equilibrium for (index_t m: m_data->range_mineral()) { const index_t idm = ideq_min(m); if (idm == no_equation) continue; for (index_t i: m_data->range_aqueous_component()) { if (ideq_paq(i) == no_equation) continue; jacobian(idm, ideq_paq(i)) = -m_data->nu_mineral(m, i); } } } void ReducedSystem::set_secondary_concentration(const Vector& x) { for (index_t j: m_data->range_aqueous()) { if (m_active_aqueous[j] == false) { m_secondary_conc(j) = 0; continue; } scalar_t conc = -m_data->logk_aqueous(j) - m_loggamma(j+m_data->nb_component); for (index_t k: m_data->range_aqueous_component()) { if (m_data->nu_aqueous(j, k) == 0) continue; conc += m_data->nu_aqueous(j, k)*(x(ideq_paq(k)) + m_loggamma(k)); } m_secondary_conc(j) = pow10(conc); } } void ReducedSystem::set_ionic_strength(const Vector& x) { scalar_t ionic = 0; for (index_t i: m_data->range_aqueous_component()) { if (ideq_paq(i) == no_equation or m_data->charge_component(i) == 0) continue; ionic += pow10(x(ideq_paq(i)))*std::pow(m_data->charge_component(i),2); } for (index_t j: m_data->range_aqueous()) { if (not m_active_aqueous[j] or m_data->charge_aqueous(j) == 0) continue; ionic += m_secondary_conc(j)*std::pow(m_data->charge_aqueous(j),2); } m_ionic_strength = ionic/2; } void ReducedSystem::compute_log_gamma(const Vector& x) { set_ionic_strength(x); const scalar_t sqrti = std::sqrt(m_ionic_strength); for (index_t i: m_data->range_aqueous_component()) { if (ideq_paq(i) == no_equation) { m_loggamma(i) = 0; continue; } m_loggamma(i) = laws::extended_debye_huckel(m_ionic_strength, sqrti, m_data->charge_component(i), m_data->a_debye_component(i), m_data->b_debye_component(i)); } for (index_t j: m_data->range_aqueous()) { if (not m_active_aqueous[j]) { m_loggamma(m_data->nb_component + j) = 0; continue; } m_loggamma(m_data->nb_component + j) = laws::extended_debye_huckel( m_ionic_strength, sqrti, m_data->charge_aqueous(j), m_data->a_debye_aqueous(j), m_data->b_debye_aqueous(j)); } } bool ReducedSystem::hook_start_iteration(const Vector& x, scalar_t norm_residual) { not_in_linesearch = true; scalar_t previous_norm = m_loggamma.norm(); bool may_have_converged = false; if (norm_residual < nb_free_variables()) { // Use fix point iterations for non-ideality for (int i=0; i<6; ++i) { set_secondary_concentration(x); compute_log_gamma(x); if (std::abs(previous_norm - m_loggamma.norm()) < 1e-6) {may_have_converged = true; break;} previous_norm = m_loggamma.norm(); } } return may_have_converged; } double ReducedSystem::max_lambda(const Vector& x, const Vector& update) { if (ideq_w() != no_equation) { return 1.0/std::max(1.0, -update(0)/(0.9*x(0))); } else { return 1.0; } } void ReducedSystem::reasonable_starting_guess(Vector &x, bool for_min) { x.resize(m_data->nb_component+ m_data->nb_mineral); x(0) = 2*m_tot_conc(0)*m_data->molar_mass_basis_si(0); for (index_t i: m_data->range_aqueous_component()) { x(i) = -4.0; } if (for_min) x.block(m_data->nb_component, 0, m_data->nb_mineral, 1).setZero(); else x.block(m_data->nb_component, 0, m_data->nb_mineral, 1).setZero(); m_loggamma.setZero(); m_secondary_conc.setZero(); } void ReducedSystem::reasonable_restarting_guess(Vector& x) { x(0) = m_tot_conc(0)*m_data->molar_mass_basis_si(0); for (index_t i=m_data->nb_component; inb_component+m_data->nb_mineral; ++i) { x(i) = 0.0; } m_loggamma.setZero(); m_secondary_conc.setZero(); } EquilibriumState ReducedSystem::get_solution(const Vector& xtot, const Vector& x) { double previous_norm = m_loggamma.norm(); set_secondary_concentration(x); compute_log_gamma(x); if (std::abs(previous_norm - m_loggamma.norm()) > 1e-6) { WARNING << "Activity coefficient have not converged !" << std::endl << "output can not be trusted\n Difference "+std::to_string(std::abs(previous_norm - m_loggamma.norm())); } return EquilibriumState(xtot, m_secondary_conc, m_loggamma, m_ionic_strength, m_data); } +void ReducedSystem::reduce_mineral_problem(Vector& true_var) +{ + for (index_t mineral: m_data->range_mineral()) + { + if (ideq_min(mineral) == no_equation) continue; + if (m_data->stability_mineral(mineral) == database::MineralStabilityClass::cannot_dissolve) + { + if (true_var(ideq_min(mineral)) == 0.0) continue; + for (index_t component: m_data->range_component()) + { + if (m_data->nu_mineral(mineral, component) == 0.0) continue; + m_tot_conc(component) -= m_data->nu_mineral(mineral, component)*true_var(ideq_min(mineral)); + } + true_var(ideq_min(mineral)) = 0.0; + } + } +} + +void ReducedSystem::reset_mineral_system(Vector& true_var, Vector& output_var) +{ + for (index_t mineral: m_data->range_mineral()) + { + if (ideq_min(mineral) == no_equation) continue; + if (m_data->stability_mineral(mineral) == database::MineralStabilityClass::cannot_dissolve) + { + output_var(mineral) += true_var(ideq_min(mineral)); + } + } +} + } // end namespace specmicp diff --git a/src/specmicp/reduced_system.hpp b/src/specmicp/reduced_system.hpp index 4d2fc47..e5a1d43 100644 --- a/src/specmicp/reduced_system.hpp +++ b/src/specmicp/reduced_system.hpp @@ -1,206 +1,213 @@ /*------------------------------------------------------- - Module : specmicp - File : specmicp.hpp - Author : Fabien Georget Copyright (c) 2014, Fabien Georget , 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: * Redistributions of source code must retain the above copyright notice, this list of conditions and the following disclaimer. * 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. * Neither the name of the Princeton University 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 OWNER 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. ---------------------------------------------------------*/ #ifndef SPECMICP_REDUCEDSYSTEM_HPP #define SPECMICP_REDUCEDSYSTEM_HPP //! \file reduced_system.hpp The MiCP program to solve speciation #include "common.hpp" #include "database.hpp" #include "kineticsdata.hpp" #include "micpsolver/micpprog.hpp" #include "simulation_box.hpp" //! \namespace specmicp main namespace for the SpecMiCP solver namespace specmicp { class EquilibriumState; //! \brief The equilibrium speciation problem //! //! Represent the equilibrium speciation problem as a Mixed Complementarity program //! It is called reduced becaused it does not solve explicitely for the secondary species. //! //! Equations //! ========== //! //! Main equations //! -------------- //! //! - Water conservation equation : \f[ T_w = m_w \left( \frac{1}{M_w} + \sum_{j=1}^{N_r} \nu_{wj} m_j \right) + \sum_{l=1}^{N_{min}} \nu_{lw} n_l \f] //! - Component conservation equations : \f[ T_i = m_w \left( m_i + \sum_{j=1}^{N_r} \nu_{ij} m_j \right) + \sum_{l=1}^{N_{min}} \nu_{li} n_l \f] //! - Mineral equilibrium conditions : \f[ n_l \geq 0, \; 1 - \mathrm{SI}_l \geq 0 , \; \mathrm{and} \; n_l ( 1 - \mathrm{SI}_l ) = 0 \f] //! //! Secondary equations //! ------------------- //! //! - Secondary aqueous species equilibrium conditions : \f[ m_j = \frac{\prod_{i=2}^{N_c} (\gamma_i m_i)^{\nu_{ji}}}{K_j \gamma_j} \f] //! - Ionic Strength : \f[ I = \sum_{k=2}^{N_c+N_r} z_k^2 m_k \f] //! - Activity coefficients : \f[ \log \gamma_k = \frac{-A z_k^2 \sqrt{I}}{1 + B \dot{a}_k \sqrt{I}} + \dot{b}_k I \f] //! class ReducedSystem : public micpsolver::MiCPProg { public: //! \brief Create a reduced system //! //! \param ptrthermo Shared pointer to the thermodynamic data //! \param totconc Vector containing the total concentrations (in mol) for each components ReducedSystem(RawDatabasePtr ptrdata, const Vector& totconc, bool conservation_water=true ); //! \brief Create a reduced system //! //! \param ptrthermo Shared pointer to the thermodynamic data //! \param totconc Vector containing the total concentrations (in mol) for each components ReducedSystem(RawDatabasePtr ptrdata, const Vector& totconc, const SimulationBox& simul_box, bool conservation_water=true ); ReducedSystem(RawDatabasePtr ptrdata, const Vector& totconc, const SimulationBox& simul_box, const EquilibriumState& previous_solution, bool conservation_water=true ); // Variables // ========== //! \brief Return the total number of variables index_t total_variables() const {return m_nb_tot_vars;} //! \brief Return the number of 'free' variables (not subject to complementarity conditions) index_t nb_free_variables() const {return m_nb_free_vars;} //! \brief Return the number of variables subjected to the complementarity conditions index_t nb_complementarity_variables() const {return m_nb_compl_vars;} // The linear system // ================== //! \brief Return the residual for the water conservation equation scalar_t residual_water(const Vector& x) const; //! \brief Return the residual for the conservation equation of component (i) //! //! For (i==0), water use the method : 'residual_water' double residual_component(const Vector& x, index_t i) const; //! \brief Equilibrium condition for the minerals \f$ 1 - \mathrm{SI}_l \f$ double residual_mineral(const Vector& x, index_t m) const; //! \brief Return the residuals void get_residuals(const Vector& x, Vector& residual); //! \brief Return the jacobian void get_jacobian(Vector& x, Matrix& jacobian); //! \brief Return the equation id index_t ideq(index_t i) const {return m_ideq[i];} //! \brief Return the equation number for conservation of water index_t ideq_w() const {return m_ideq[0];} //! \brief Return the equation number of component 'i' index_t ideq_paq(index_t i) const { specmicp_assert(i < m_data->nb_component); return m_ideq[i]; } //! \brief Return the equation number of mineral 'm' index_t ideq_min(index_t m) const { specmicp_assert(m < m_data->nb_mineral); return m_ideq[m+m_data->nb_component]; } //! \brief Return the mass of water scalar_t mass_water(const Vector& x) const { if (ideq_w() != no_equation) return x(ideq_w()); else return 1.0; //FIXME } //! \brief Compute the ionic strength void set_ionic_strength(const Vector& x); //! \brief Compute the activity coefficients void compute_log_gamma(const Vector& x); //! \brief Compute the activity model, called by the solver at the beginning of an iteration bool hook_start_iteration(const Vector &x, scalar_t norm_residual); //! \brief Return a scale factor to avoid negative mass during Newton's iteration double max_lambda(const Vector &x, const Vector &update); //! \brief Return the component that does not exist in the solution (inactive dof) const std::vector& get_non_active_component() {return m_nonactive_component;} //! \brief A reasonable (well... maybe...) starting guess void reasonable_starting_guess(Vector& x, bool for_min=false); //! \brief A reasonable (maybe...) restarting guess void reasonable_restarting_guess(Vector& x); //! Return the equilibrium state of the system, the 'Solution' of the speciation problem EquilibriumState get_solution(const Vector& xtot, const Vector& x); + //! \brief Return the total concentration of 'component' + scalar_t total_concentration(index_t component) const {return m_tot_conc(component);} + + //! \brief Reduce the mineral system if one cannot dissolve + void reduce_mineral_problem(Vector& true_var); + //! \brief Reset the mineral system to the true problem + void reset_mineral_system(Vector& true_var, Vector& output_var); + private: // For the normal algorithm void set_secondary_concentration(const Vector& x); // For the tot aq_conc pb //void unsafe_set_secondary_concentration(const Eigen::VectorXd& x); void number_eq(bool water_equation); SimulationBox m_simulbox; RawDatabasePtr m_data; Vector m_tot_conc; - index_t m_nb_tot_vars; index_t m_nb_free_vars; index_t m_nb_compl_vars; std::vector m_ideq; Vector m_secondary_conc; Vector m_gas_pressure; scalar_t m_ionic_strength; Vector m_loggamma; bool not_in_linesearch; std::vector m_nonactive_component; std::vector m_active_aqueous; std::vector m_active_gas; }; } // end namespace specmicp #endif // SPECMICP_REDUCEDSYSTEM_HPP diff --git a/src/specmicp/reduced_system_solver.cpp b/src/specmicp/reduced_system_solver.cpp index c2adacb..3a814ad 100644 --- a/src/specmicp/reduced_system_solver.cpp +++ b/src/specmicp/reduced_system_solver.cpp @@ -1,256 +1,258 @@ /*------------------------------------------------------- - Module : specmicp - File : reduced_system_solver - Author : Fabien Georget Copyright (c) 2014, Fabien Georget , 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: * Redistributions of source code must retain the above copyright notice, this list of conditions and the following disclaimer. * 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. * Neither the name of the Princeton University 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 OWNER 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 "reduced_system_solver.hpp" #include "micpsolver/micpsolver.hpp" #include "micpsolver/micpsolverold.hpp" //#include "micpsolver/micpsolver_min.hpp" #include "equilibrium_data.hpp" #include "simulation_box.hpp" #include namespace specmicp { ReducedSystemSolver::ReducedSystemSolver( RawDatabasePtr data, const Vector& totconc, const SimulationBox& simulbox, ReducedSystemSolverOptions options ): m_options(options), m_data(data), m_system(std::make_shared( data, totconc, simulbox, options.conservation_water)), m_var(Vector::Zero(data->nb_component+data->nb_mineral)) {} ReducedSystemSolver::ReducedSystemSolver( RawDatabasePtr data, const Vector& totconc, const SimulationBox& simulbox, const EquilibriumState& previous_solution, ReducedSystemSolverOptions options ): m_options(options), m_data(data), m_system(std::make_shared( data, totconc, simulbox, previous_solution, options.conservation_water)), m_var(Vector::Zero(data->nb_component+data->nb_mineral)) {} ReducedSystemSolver::ReducedSystemSolver( RawDatabasePtr data, const Vector& totconc, ReducedSystemSolverOptions options ): m_options(options), m_data(data), m_system(std::make_shared( data, totconc, options.conservation_water)), m_var(Vector::Zero(data->nb_component+data->nb_mineral)) {} ReducedSystemSolver::ReducedSystemSolver( RawDatabasePtr data, const Vector& totconc, const EquilibriumState& previous_solution, ReducedSystemSolverOptions options ): m_options(options), m_data(data), m_system(std::make_shared( data, totconc, SimulationBox(), previous_solution, options.conservation_water)), m_var(Vector::Zero(data->nb_component+data->nb_mineral)) {} EquilibriumState ReducedSystemSolver::get_solution(const Eigen::VectorXd& x) { set_true_variable_vector(x); return m_system->get_solution(x, m_var); } micpsolver::MiCPPerformance ReducedSystemSolver::solve(Vector& x, bool init) { if (init) { if (get_options().ncp_function == micpsolver::NCPfunction::min) m_system->reasonable_starting_guess(x, true); else m_system->reasonable_starting_guess(x, false); } set_true_variable_vector(x); micpsolver::MiCPPerformance perf = solve(); if (perf.return_code != micpsolver::MiCPSolverReturnCode::ResidualMinimized and get_options().allow_restart) { WARNING << "Failed to solve the system ! - We shake it up and start again"; const micpsolver::NCPfunction save_ncp = get_options().ncp_function; const int save_penalization_factor = get_options().solver_options.penalization_factor; if (save_ncp == micpsolver::NCPfunction::min) { get_options().ncp_function = micpsolver::NCPfunction::penalizedFB; m_system->reasonable_restarting_guess(x); } else { if (save_penalization_factor == 1) get_options().solver_options.penalization_factor = 0.8; set_return_vector(x); m_system->reasonable_restarting_guess(x); } set_true_variable_vector(x); micpsolver::MiCPPerformance perf2 = solve(); get_options().ncp_function = save_ncp; get_options().solver_options.penalization_factor = save_penalization_factor; perf += perf2; } if (perf.return_code > micpsolver::MiCPSolverReturnCode::NotConvergedYet) set_return_vector(x); return perf; } void ReducedSystemSolver::set_true_variable_vector(const Vector& x) { const std::vector& non_active_component = m_system->get_non_active_component(); if (non_active_component.size() == 0) { // we still copy the data, if we failed to solve the problem, we can restart if (m_options.conservation_water) m_var = x; // direct copy else m_var = x.block(1, 0, x.rows()-1, 1); for (int i=0; irange_aqueous_component()) { auto it = std::find(non_active_component.begin(), non_active_component.end(),i); if (it != non_active_component.end()) continue; m_var(new_i) = x(i); ++new_i; } for (index_t m: m_data->range_mineral()) { for (auto it=non_active_component.begin(); it!=non_active_component.end(); ++it) { if (m_data->nu_mineral(m, *it) != 0) continue; m_var(new_i) = x(m_data->nb_component+m); ++new_i; } } m_var.conservativeResize(new_i); } + m_system->reduce_mineral_problem(m_var); } micpsolver::MiCPPerformance ReducedSystemSolver::solve() { if (get_options().ncp_function == micpsolver::NCPfunction::penalizedFB) { // the default micpsolver::MiCPSolver solver(m_system); solver.set_options(get_options().solver_options); solver.solve(m_var); return solver.get_performance(); } else { micpsolver::MiCPSolverOLD solver(m_system); solver.set_options(get_options().solver_options); solver.solve(m_var); return solver.get_performance(); } } void ReducedSystemSolver::set_return_vector(Vector& x) { const std::vector& non_active_component = m_system->get_non_active_component(); if (non_active_component.size() == 0) { if (m_options.conservation_water) x = m_var; //direct copy else x.block(1, 0, x.rows()-1, 1) = m_var; // at that point we should have the correct solution } else { uindex_t new_i = 0; if (m_options.conservation_water) { x(0) = m_var(new_i); ++new_i; } for (index_t i: m_data->range_aqueous_component()) { auto it = std::find(non_active_component.begin(), non_active_component.end(),i); if (it != non_active_component.end()) { x(i) = -HUGE_VAL; continue; } x(i) = m_var(new_i) ; ++new_i; } for (index_t m: m_data->range_mineral()) { for (auto it=non_active_component.begin(); it!=non_active_component.end(); ++it) { if (m_data->nu_mineral(m, *it) != 0) { x(m_data->nb_component+m) = 0; continue; } x(m_data->nb_component+m) =m_var(new_i); ++new_i; } } } + m_system->reset_mineral_system(m_var, x); } } // end namespace specmicp