Page MenuHomec4science

diffusion.hpp
No OneTemporary

File Metadata

Created
Sun, Aug 25, 09:21

diffusion.hpp

/*-------------------------------------------------------
- Module : reactmicp/systems
- File : diffusion.hpp
- Author : Fabien Georget
Copyright (c) 2014, Fabien 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:
* 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_REACTMICP_SYSTEMS_DIFFUSION_HPP
#define SPECMICP_REACTMICP_SYSTEMS_DIFFUSION_HPP
#include <vector>
#include <memory>
#include "database/data_container.hpp"
#include "reactmicp/common.hpp"
#include "reactmicp/micpsolvers/parabolicprogram.hpp"
#include "reactmicp/meshes/mesh1d.hpp"
namespace specmicp {
namespace reactmicp {
namespace systems {
class DiffusionProgram: solvers::ParabolicProgram<DiffusionProgram>
{
public:
using Variables = solvers::ParabolicVariables;
const int no_equation = -1;
//! \brief Return the number of degree of freedoms
ind_t get_neq() const {return m_neq;}
//! \brief Return the residuals
void ebe_residuals(const Variables& variable,
Vector& residual);
//! \brief Residuals for an element
void elemental_residuals(int element,
const Variables& variable,
Vector& residual);
//! brief Return the jacobian
void ebe_jacobian(const Variables& variable,
Vector& residual,
SparseMatrix& jacobian);
//! \brief Return the index of the MiCP equations
const std::vector<bool>& micp_index() const { return m_is_micp;}
private:
// Element by element procedures
// =============================
// Degree of freedom numbering
// ===========================
//! \brief Return the number of the dof
ind_t get_dof(ind_t node, ind_t elemental_dof) const {
return elemental_dof + m_ndf*node;}
//! \brief Return the ndof of a diffusion equation
ind_t get_ndof_diffusion(ind_t ind_component) const {
return ind_diffusion;}
//! \brief Return the ndof of a mass balance equation
ind_t get_ndof_massbalance(ind_t ind_component) const {
return m_data->nb_component + ind_component;}
//! \brief Return the ndof of a mineral equation
ind_t get_ndof_mineral(ind_t ind_mineral) const {
return 2*m_data->nb_component + ind_mineral;}
//! \brief Return the dof of a diffusion equation
ind_t get_dof_diffusion(ind_t node, ind_t ind_component) const {
return get_dof(node, get_ndof_diffusion(ind_component));}
//! \brief Return the dof of a mass balance equation
ind_t get_dof_massbalance(ind_t node, ind_t ind_component) const {
return get_dof(node, get_ndof_massbalance(ind_component));}
//! \brief Return the dof of a mineral equation
ind_t get_dof_mineral(ind_t node, ind_t ind_mineral) const {
return get_dof(node, get_ndof_mineral(ind_mineral));}
//! \brief Return the index of the equation
//!
//! Return no_equation if there is no equation for this degree of freedom
ind_t ideq(ind_t index) const {return m_ideq[index];}
//! \brief Return the equation's index of a diffusion equation
ind_t get_ideq_diffusion(ind_t node, ind_t ind_component) const {
return get_ideq(get_dof_diffusion(node, ind_component));}
//! \brief Return the equation's index of a mass balance equation
ind_t get_ideq_massbalance(ind_t node, ind_t ind_component) const {
return get_ideq(get_dof_massbalance(node, ind_component));}
//! \brief Return the equation's index of a mineral equation
ind_t get_ideq_mineral(ind_t node, ind_t ind_mineral) const {
return get_ideq(get_dof_mineral(node, ind_mineral));}
// Secondary variables
// ===================
// Algorithms
// ----------
//! \brief Set secondary concentrations
void compute_secondary_concentrations(const Variables& variable);
//! \brief Set secondary concentrations for node 'node'
void nodal_secondary_concentrations(ind_t node, const Variables& variable);
//! \brief Solve for secondary variables
int solve_secondary_variables(const Variables& variable);
//! \brief Compute secondary variables
void compute_secondary_variables(const Variables& variable);
//! \brief compute the ionic strengths
void compute_ionic_strength(const Variables& variable);
//! \brief compute the ionic strength for node 'node'
void nodal_ionic_strength(ind_t node, const Variables& variable);
//! \brief compute the activity coefficients
void compute_loggamma();
//! \brief compute the logarithm of the activity coefficients for node 'node'
void nodal_loggamma(ind_t node);
//! \brief compute the activity coefficient for the components of node 'node'
void nodal_loggamma_component(ind_t node);
//! \brief compute the activity coefficient for the aqueous species of node 'node'
void nodal_loggamma_aqueous(ind_t node);
// Getters and setters
// -------------------
//! \brief Secondary species concentration
double& secondary_concentration(ind_t node, ind_t aqueous) {
return m_secondary_concentration.coeffRef(aqueous, node);
}
//! \brief Secondary species concentration
double secondary_concentration(ind_t node, ind_t aqueous) const {
return m_secondary_concentration.coeff(aqueous, node);
}
//! \brief Return the loggamma of a component
double& loggamma_component(ind_t node, ind_t component) {
return m_secondary_variables.coeffRef(component, node);
}
//! \brief Return the loggamma of a component
double loggamma_component(ind_t node, ind_t component) const {
return m_secondary_variables.coeff(component, node);
}
//! \brief Return the loggamma of a secondary species
double& loggamma_secondary(ind_t node, ind_t species) {
return m_secondary_variables.coeffRef(m_data->nb_component + species, node);
}
//! \brief Return the loggamma of a secondary species
double loggamma_secondary(ind_t node, ind_t species) const {
return m_secondary_variables.coeff(m_data->nb_component + species, node);
}
//! \brief Return ionic strength
double& ionic_strength(ind_t node) {
return m_secondary_variables.coeffRef(m_data->nb_component+m_data->nb_aqueous, node);
}
//! \brief Return ionic strength
double ionic_strength(ind_t node) const {
return m_secondary_variables.coeff(m_data->nb_component+m_data->nb_aqueous, node);
}
// Member
// =======
ind_t m_ndf;
ind_t m_neq;
std::vector<bool> m_is_micp;
std::vector<ind_t> m_ideq;
SparseVector m_boundary_condition;
Eigen::Matrix<double> m_secondary_variables;
Eigen::Matrix<double> m_secondary_concentration;
std::shared_ptr<mesh::Mesh1D> m_mesh;
std::shared_ptr<database::DataContainer> m_data;
};
} // end namespace systems
} // end namespace reactmicp
} // end namespace specmicp
#endif // SPECMICP_REACTMICP_SYSTEMS_DIFFUSION_HPP

Event Timeline