Page MenuHomec4science

variables.hpp
No OneTemporary

File Metadata

Created
Thu, Aug 15, 22:01

variables.hpp

#ifndef SPECMICP_REACTMICP_SYSTEMS_SATURATED_VARIABLES_HPP
#define SPECMICP_REACTMICP_SYSTEMS_SATURATED_VARIABLES_HPP
#include "database.hpp"
#include "reactmicp/solver/staggers_base/variables_base.hpp"
#include "specmicp/adimensional/adimensional_system_solution.hpp"
#include <vector>
// forward declaration
// ===================
#include "dfpm/meshes/mesh1dfwd.hpp"
namespace specmicp {
namespace reactmicp {
namespace solver {
using VariablesBasePtr = std::shared_ptr<VariablesBase>;
}
namespace systems {
namespace satdiff {
class SaturatedVariablesFactory;
} // end namespace satdiff
} // end namespace systems
} // end namespace reactmicp
} // end namespace specmicp
// Class declaration
// =================
namespace specmicp {
namespace reactmicp {
namespace systems {
namespace satdiff {
//! \brief Variables for the saturated reactive transport system
//!
//! Contain all the variables that need to be shared between the staggers
class SaturatedVariables: public solver::VariablesBase
{
// SaturatedVariablesFactory should be the class to use to inialize
// the variables correctly
friend class SaturatedVariablesFactory;
public:
SaturatedVariables(mesh::Mesh1DPtr the_mesh,
RawDatabasePtr the_database);
//! \brief Return the mesh
mesh::Mesh1DPtr get_mesh() {return m_mesh;}
//! \brief Return the database
RawDatabasePtr get_database() {return m_database;}
//! \brief Return the number of components
index_t nb_component() {return m_database->nb_component();}
//! \brief Return the number of nodes
index_t nb_nodes() {return m_is_fixed_composition.size();}
//! \brief Return true if 'node' has a fixed composition
index_t is_fixed_composition(index_t node) {return m_is_fixed_composition[node];}
// Main variables
// ==============
//! \brief Return the main variable vector
Vector& displacement() {return m_displacement;}
//! \brief Return the main variable vector at the beginning of the timestep
Vector& predictor() {return m_predictor;}
//! \brief Return the velocity of the main variables
Vector& velocity() {return m_velocity;}
//! \brief Return the rate of change of the main variables due to the transport operator
Vector& transport_rate() {return m_transport_rate;}
//! \brief Return the rate of change of the main variables due to the chemistry operator
Vector& chemistry_rate() {return m_chemistry_rate;}
// Access to main variables
// ========================
//! \brief Return the number of degree of freedom (per node) in the main variables vector
index_t ndf() {return 2*m_database->nb_component();}
//! \brief Return the offset of 'node' in the main variables vector
index_t offset_node(index_t node) {return node*ndf();}
//! \brief Return the offset of the aqueous concentration variables in the main variables vector
index_t offset_aqueous_concentration() {return 0;}
//! \brief Return the offset of the aqueous concentrations variables in the main variables vector
index_t offset_aqueous_concentration(index_t node) {
return offset_aqueous_concentration()+offset_node(node);}
//! \brief Return the offset of the solid concentration variables in the main variables vector
index_t offset_solid_concentration() {return m_database->nb_component();}
//! \brief Return the offset of the solid concentrations variables in the main variables vector
index_t offset_solid_concentration(index_t node) {
return offset_solid_concentration()+offset_node(node);}
//! \brief Return the degree of freedom number for the aqueous concentration of 'component' at 'node'
index_t dof_aqueous_concentration(index_t node, index_t component) {
return (component + offset_aqueous_concentration(node));
}
//! \brief Return the degree of freedom number for the solid concentration of 'component' at 'node'
index_t dof_solid_concentration(index_t node, index_t component) {
return (component + offset_solid_concentration(node));
}
//! \brief Return the aqueous concentration of 'component' at 'node' in 'var'
//!
//! 'var' is any of the main variables vector, it may be a velocity vector
scalar_t& aqueous_concentration(index_t node, index_t component, Vector& var) {
return var(dof_aqueous_concentration(node, component));
}
//! \brief Return the aqueous concentration of 'component' at 'node' in main variables
//!
//! 'var' is any of the main variables vector, it may be a velocity vector
scalar_t& aqueous_concentration(index_t node, index_t component) {
return m_displacement(dof_aqueous_concentration(node, component));
}
//! \brief Return the solid concentration of 'component' at 'node' in 'var'
//!
//! 'var' is any of the main variables vector, it may be a velocity vector
scalar_t& solid_concentration(index_t node, index_t component, Vector& var){
return var(dof_solid_concentration(node, component));
}
//! \brief Return the solid concentration of 'component' at 'node' in main variables
//!
//! 'var' is any of the main variables vector, it may be a velocity vector
scalar_t& solid_concentration(index_t node, index_t component){
return m_displacement(dof_solid_concentration(node, component));
}
//! \brief Return a vector containing the total concentrations computed from the main variables
//!
//! This is to be used to restart the chemistry computation
Vector total_concentrations(index_t node);
// Equilibrium
// ===========
//! \brief Returh the solution of the speciation solver at 'node'
AdimensionalSystemSolution& equilibrium_solution(index_t node) {
return m_equilibrium_solutions[node];
}
// Upscaling
// =========
//! \brief Return the offset for 'node' in the upscaling variables vector
index_t offset_node_upscaling(index_t node) {return ndf_upscaling()*node;}
//! \brief Return the number fo degree of freedom (per node) for the upscaling vector
index_t ndf_upscaling() {return 5;}
//! \brief Return the degree of freedom for the porosity at 'node'
index_t dof_porosity(index_t node) {return 0 + offset_node_upscaling(node);}
//! \brief Return the degree of freedom of the porosity velocity at 'node'
index_t dof_vel_porosity(index_t node) {return 1 + offset_node_upscaling(node);}
//! \brief Return the degree of freedom of the diffusion coefficient at 'node'
index_t dof_diffusion_coefficient(index_t node) {return 2 + offset_node_upscaling(node);}
//! \brief Return the degree of freedom of the permeability at 'node'
index_t dof_permeability(index_t node) {return 3 + offset_node_upscaling(node);}
//! \brief Return the fluid velocity
index_t dof_fluid_velocity(index_t node) {return 4 + offset_node_upscaling(node);}
//! \brief Return the porosity at 'node'
scalar_t& porosity(index_t node) {return m_upscaling(dof_porosity(node));}
//! \brief Return the rate of change of the porosity at 'node'
scalar_t& vel_porosity(index_t node) {return m_upscaling(dof_vel_porosity(node));}
//! \brief Return the diffusion coefficient at 'node'
scalar_t& diffusion_coefficient(index_t node) {return m_upscaling(dof_diffusion_coefficient(node));}
//! \brief Return the permeability at 'node'
scalar_t& permeability(index_t node) {return m_upscaling(dof_permeability(node));}
//! \brief Return the fluid velocity at 'node'
scalar_t& fluid_velocity(index_t node) {return m_upscaling(dof_fluid_velocity(node));}
//! \brief Return the vector of upscaling variables
Vector& upscaling_variables() {return m_upscaling;}
//! \brief Reset the main variables
void reset_main_variables() override;
private:
// ############ //
// Attributes //
// ############ //
mesh::Mesh1DPtr m_mesh;
RawDatabasePtr m_database;
std::vector<bool> m_is_fixed_composition;
// Main variables
// ==============
Vector m_displacement;
Vector m_predictor;
Vector m_velocity;
Vector m_transport_rate;
Vector m_chemistry_rate;
// Equilibrium
// ===========
std::vector<AdimensionalSystemSolution> m_equilibrium_solutions;
// Upscaling
// =========
Vector m_upscaling;
};
//! \brief typedef of a shared pointer of a SaturatedVariables
using SaturatedVariablesPtr = std::shared_ptr<SaturatedVariables>;
// Casting function
// =================
//! \brief Static cast to a SaturatedVariablesPtr
inline SaturatedVariablesPtr cast_var_from_base(solver::VariablesBasePtr var)
{
return std::static_pointer_cast<SaturatedVariables>(var);
}
//! \brief Static cast from a SaturatedVariablesPtr
inline solver::VariablesBasePtr cast_var_to_base(SaturatedVariablesPtr var)
{
return std::static_pointer_cast<solver::VariablesBase>(var);
}
} // end namespace satdiff
} // end namespace systems
} // end namespace reactmicp
} // end namespace specmicp
#endif // SPECMICP_REACTMICP_SYSTEMS_SATURATED_VARIABLES_HPP

Event Timeline