Page MenuHomec4science

diffusion_numbering.hpp
No OneTemporary

File Metadata

Created
Tue, Jul 23, 05:28

diffusion_numbering.hpp

/*-------------------------------------------------------
- Module : reactmicp/systems/diffusion
- File : diffusion_numbering.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_NUMBERING_HPP
#define SPECMICP_REACTMICP_SYSTEMS_DIFFUSION_NUMBERING_HPP
#include <memory>
#include "reactmicp/common.hpp"
#include "specmicp/equilibrium_data.hpp"
namespace specmicp {
namespace reactmicp {
namespace mesh {
// forward declaration
class Mesh1D;
} // end namespace mesh
namespace systems {
namespace diffusion {
struct BoundaryCondition
{
int node;
EquilibriumState composition;
std::vector<int> tot_mobile_concentration;
std::vector<int> component;
std::vector<int> mineral;
void fix_all();
};
using ListBoundaryCondition = std::vector<BoundaryCondition>;
//! \brief Equation numbering scheme
class DiffusionNumbering
{
public:
const int no_equation = -1;
using Variables = ParabolicVariables;
DiffusionNumbering(ind_t number_nodes,
std::shared_ptr<database::DataContainer> data,
const ListBoundaryCondition& list_bc,
Variables& variable):
DiffusionNumbering(number_nodes, data->nb_component-1, data->nb_mineral, list_bc, variable)
{}
DiffusionNumbering(ind_t number_nodes,
ind_t number_aqueous_component,
ind_t number_mineral,
const ListBoundaryCondition& list_bc,
Variables& variable);
//! \brief Return the number of equations
ind_t get_neq() const {return nb_neq;}
//! \brief Return the number of degree of freedom (per node)
ind_t get_ndf() const {return nb_ndf;}
// 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 + nb_ndf*node;}
//! \brief Return the ndof of a diffusion equation
ind_t get_ndof_diffusion(ind_t ind_component) const {
return ind_component-1;}
//! \brief Return the ndof of a mass balance equation
ind_t get_ndof_massbalance(ind_t ind_component) const {
return nb_aq_component + ind_component - 1;}
//! \brief Return the ndof of a mineral equation
ind_t get_ndof_mineral(ind_t ind_mineral) const {
return 2*nb_aq_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 id_equation(ind_t node, ind_t elemental_dof) const {return m_ideq(elemental_dof, node);}
//! \brief Return the equation's index of a diffusion equation
ind_t id_equation_diffusion(ind_t node, ind_t ind_component) const {
return id_equation(node, get_ndof_diffusion(ind_component));}
//! \brief Return the equation's index of a mass balance equation
ind_t id_equation_massbalance(ind_t node, ind_t ind_component) const {
return id_equation(node, get_ndof_massbalance(ind_component));}
//! \brief Return the equation's index of a mineral equation
ind_t id_equation_mineral(ind_t node, ind_t ind_mineral) const {
return id_equation(node, get_ndof_mineral(ind_mineral));}
// Access primary variable
// =======================
//! \brief Return the total mobile concentration for 'component' at 'node'
double mobile_total_concentration(ind_t node, ind_t component, const Variables& var) const {
return var.displacement(get_dof_diffusion(node, component));
}
//! \brief Return a reference to the total mobile concentration for 'component' at 'node'
double& mobile_total_concentration(ind_t node, ind_t component, Variables& var) {
return var.displacement.coeffRef(get_dof_diffusion(node, component));
}
//! \brief Return the component concentration for 'component' at 'node'
double component_concentration(ind_t node, ind_t component, const Variables& var) const {
return var.displacement(get_dof_massbalance(node, component));
}
//! \brief Return a reference to the component concentration for 'component' at 'node'
double& component_concentration(ind_t node, ind_t component, Variables& var) {
return var.displacement.coeffRef(get_dof_massbalance(node, component));
}
//! \brief Return the mineral mole number for 'mineral' at 'node'
double mole_mineral(ind_t node, ind_t mineral, const Variables& var) const {
return var.displacement(get_dof_mineral(node, mineral));
}
//! \brief Return a reference to the mineral mole number for 'mineral' at 'node'
double& mole_mineral(ind_t node, ind_t mineral, Variables& var) {
return var.displacement.coeffRef(get_dof_mineral(node, mineral));
}
//! \brief Return true if the node has a boundary condition
double node_has_bc(ind_t node);
private:
//! \brief Initialize the variable
void initialize_variable(Variables& variable);
//! \brief Parse the BCs and number the equations
void number_equations(const ListBoundaryCondition& list_bc,
Variables& variable);
//! \brief Parse the BCS, set the corresponding values in variable
void parse_bcs(const ListBoundaryCondition& list_bc,
Variables& variable);
//! \brief number the equations
void do_numbering();
ind_t nb_aq_component;
ind_t nb_mineral;
ind_t nb_ndf;
ind_t nb_neq;
Eigen::MatrixXi m_ideq;
std::vector<ind_t> m_nodes_with_bc;
};
} // end namespace diffusion
} // end namespace systems
} // end namespace reactmicp
} // end namespace specmicp
#endif // SPECMICP_REACTMICP_SYSTEMS_DIFFUSION_NUMBERING_HPP

Event Timeline