Page MenuHomec4science

reactive_transport_solver.hpp
No OneTemporary

File Metadata

Created
Wed, Aug 14, 21:22

reactive_transport_solver.hpp

#ifndef SPECMICP_REACTMICP_SATURATED_DIFFUSION_REACTIVETRANSPORTSOLVER_HPP
#define SPECMICP_REACTMICP_SATURATED_DIFFUSION_REACTIVETRANSPORTSOLVER_HPP
#include "common.hpp"
#include "options.hpp"
#include "variables.hpp"
#include "transport_program.hpp"
#include "dfpmsolver/parabolic_driver.hpp"
#include "boundary_conditions.hpp"
namespace specmicp {
namespace reactmicp {
using namespace systems::siasaturated;
namespace internals {
struct SIAResiduals; // forward declaration
} // end namespace internals
//! \brief SIA solver for saturated reactive transport problem
//!
class SIASaturatedReactiveTransportSolver
{
public:
using retcode_t = SIASaturatedReactiveTransportSolverReturnCode;
SIASaturatedReactiveTransportSolver(
std::shared_ptr<mesh::Mesh1D> the_mesh,
RawDatabasePtr the_database,
std::shared_ptr<SaturatedDiffusionTransportParameters> transport_parameters
):
m_mesh(the_mesh),
m_database(the_database),
m_variables(the_database, the_mesh, transport_parameters),
m_transport_program(the_mesh, the_database, transport_parameters),
m_transport_solver(m_transport_program),
m_bcs(the_mesh->nnodes())
{}
void apply_boundary_conditions(SIABoundaryConditions bcs);
SIASaturatedVariables& get_variables() {return m_variables;}
//! \brief Solve a timestep
SIASaturatedReactiveTransportSolverReturnCode solve_timestep(scalar_t dt);
//! \brief Return a read-write reference to the options
SIASaturatedReactiveTransportSolverOptions& get_options() {return m_options;}
//! \brief Return a read-only reference to the options
const SIASaturatedReactiveTransportSolverOptions& get_options() const {return m_options;}
//! \brief Use the sequential iterative algorithm
void use_sia(int nb_iterations) {get_options().use_sia(nb_iterations);}
//! \brief Use the sequential non-iterative algorithm
void use_snia() {get_options().use_snia();}
private:
//! \brief Initialize the problem
void initialize();
//! \brief Restart the transport problem
void restart_transport_iterations(internals::SIAResiduals& residuals_storage);
//! \brief Restart the speciation problems
void restart_speciation_iterations(scalar_t dt, internals::SIAResiduals& residuals_storage);
//! \brief Solve the speciation problem at one node
void solve_speciation_node(scalar_t dt, index_t node, internals::SIAResiduals& residuals_storage);
retcode_t check_convergence(internals::SIAResiduals& residuals_storage);
scalar_t residuals_transport();
SIASaturatedReactiveTransportSolverOptions m_options;
std::shared_ptr<mesh::Mesh1D> m_mesh;
std::shared_ptr<database::DataContainer> m_database;
SIASaturatedVariables m_variables;
Eigen::MatrixXd m_minerals_start_iteration;
// transport program
SaturatedDiffusionProgram m_transport_program;
dfpmsolver::ParabolicDriver<SaturatedDiffusionProgram> m_transport_solver;
//
Eigen::VectorXi m_bcs;
Eigen::VectorXi m_flag_compute_speciation;
};
} // end namespace reactmicp
} // end namespace specmicp
#endif // SPECMICP_REACTMICP_SATURATED_DIFFUSION_REACTIVETRANSPORTSOLVER_HPP

Event Timeline