Page Menu
Home
c4science
Search
Configure Global Search
Log In
Files
F91421071
transport_program.hpp
No One
Temporary
Actions
Download File
Edit File
Delete File
View Transforms
Subscribe
Mute Notifications
Award Token
Subscribers
None
File Metadata
Details
File Info
Storage
Attached
Created
Sun, Nov 10, 23:00
Size
4 KB
Mime Type
text/x-c++
Expires
Tue, Nov 12, 23:00 (2 d)
Engine
blob
Format
Raw Data
Handle
22260835
Attached To
rSPECMICP SpecMiCP / ReactMiCP
transport_program.hpp
View Options
#ifndef SPECMICP_REACTMICP_SATURATED_DIFFUSION_TRANSPORTPROGRAM_HPP
#define SPECMICP_REACTMICP_SATURATED_DIFFUSION_TRANSPORTPROGRAM_HPP
#include <memory>
#include "database/data_container.hpp"
#include "dfpmsolver/parabolic_program.hpp"
#include "transport_parameters.hpp"
#include "reactmicp/meshes/mesh1d.hpp"
namespace
specmicp
{
namespace
reactmicp
{
namespace
systems
{
namespace
siasaturated
{
//! \brief Transport program for saturated diffusion
class
SaturatedDiffusionProgram
:
dfpmsolver
::
ParabolicProgram
<
SaturatedDiffusionProgram
>
{
public
:
using
triplet_t
=
Eigen
::
Triplet
<
scalar_t
>
;
//!< Triplet type for building the transport matrix
using
list_triplet_t
=
std
::
vector
<
triplet_t
>
;
//!< List of triplet for building the sparse matrix
SaturatedDiffusionProgram
(
std
::
shared_ptr
<
mesh
::
Mesh1D
>
the_mesh
,
std
::
shared_ptr
<
database
::
DataContainer
>
the_data
,
std
::
shared_ptr
<
SaturatedDiffusionTransportParameters
>
the_param
)
:
m_mesh
(
the_mesh
),
m_data
(
the_data
),
m_param
(
the_param
),
m_external_flux
(
Eigen
::
VectorXd
::
Zero
(
the_mesh
->
nnodes
()
*
(
the_data
->
nb_component
-
1
)))
{}
//! \brief Return the number of equations
index_t
get_neq
()
const
{
return
m_neq
;}
//! \brief Return the number of degrees of freedom per node
index_t
get_ndf
()
const
{
return
m_data
->
nb_component
-
1
;}
// number the equations
void
number_equations
(
std
::
vector
<
index_t
>
list_fixed_nodes
);
//! \brief Compute the residuals
void
compute_residuals
(
const
Vector
&
displacement
,
const
Vector
&
velocity
,
Vector
&
residual
);
//! \brief Compute the jacobian
void
compute_jacobian
(
Vector
&
displacement
,
Vector
&
velocity
,
Eigen
::
SparseMatrix
<
scalar_t
>
&
jacobian
,
scalar_t
alphadt
);
//! \brief Update the solutions
void
update_solution
(
const
Vector
&
update
,
scalar_t
lambda
,
scalar_t
alpha_dt
,
Vector
&
predictor
,
Vector
&
displacement
,
Vector
&
velocity
);
//! \brief Apply boundary conditions to the velocity vector
//!
//! by default do nothing.
void
apply_bc
(
scalar_t
dt
,
const
Vector
&
displacement
,
Vector
velocity
)
{}
//! Return the degree of freedom number of 'component' at 'node'
index_t
get_dof_component
(
index_t
node
,
index_t
component
)
const
{
specmicp_assert
(
component
>=
1
and
component
<
m_data
->
nb_component
);
specmicp_assert
(
node
>=
0
and
node
<
m_mesh
->
nnodes
());
return
(
component
-
1
)
+
get_ndf
()
*
node
;
}
//! Return the degree of freedom number of 'component' at 'node'
index_t
get_dof_mineral
(
index_t
node
,
index_t
mineral
)
const
{
specmicp_assert
(
mineral
>=
0
and
mineral
<
m_data
->
nb_mineral
);
specmicp_assert
(
node
>=
0
and
node
<
m_mesh
->
nnodes
());
return
mineral
+
m_data
->
nb_mineral
*
node
;
}
//! Return the total mobile concentration
scalar_t
&
get_total_mobile_concentration
(
index_t
node
,
index_t
component
,
Vector
&
concentrations
)
const
{
return
concentrations
.
coeffRef
(
get_dof_component
(
node
,
component
));
}
//! Return the total mobile concentration
scalar_t
get_total_mobile_concentration
(
index_t
node
,
index_t
component
,
const
Vector
&
concentrations
)
const
{
return
concentrations
.
coeff
(
get_dof_component
(
node
,
component
));
}
//! Return the external flux
Vector
&
external_flux
()
{
return
m_external_flux
;}
//! Return the external flux for 'component' at 'node'
scalar_t
&
external_flux
(
index_t
node
,
index_t
component
)
{
return
m_external_flux
(
get_dof_component
(
node
,
component
));}
private
:
void
element_residual_transport
(
index_t
element
,
const
Eigen
::
VectorXd
&
displacement
,
const
Eigen
::
VectorXd
&
velocity
,
Eigen
::
VectorXd
&
residual
);
void
element_residual_transport_component
(
index_t
element
,
index_t
component
,
const
Vector
&
displacement
,
const
Vector
&
velocity
,
Vector
&
element_residual
);
void
element_jacobian_transport
(
index_t
element
,
Vector
&
displacement
,
Vector
&
velocity
,
list_triplet_t
&
jacobian
,
scalar_t
alphadt
);
index_t
m_neq
;
Vector
m_ideq
;
std
::
shared_ptr
<
mesh
::
Mesh1D
>
m_mesh
;
std
::
shared_ptr
<
database
::
DataContainer
>
m_data
;
std
::
shared_ptr
<
SaturatedDiffusionTransportParameters
>
m_param
;
Vector
m_external_flux
;
};
}
// end namespace siasaturated
}
// end namespace systems
}
// end namespace reactmicp
}
// end namespace specmicp
#endif
// SPECMICP_REACTMICP_SATURATED_DIFFUSION_TRANSPORTPROGRAM_HPP
Event Timeline
Log In to Comment