Page Menu
Home
c4science
Search
Configure Global Search
Log In
Files
F86425326
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, Oct 6, 10:23
Size
7 KB
Mime Type
text/x-c++
Expires
Tue, Oct 8, 10:23 (2 d)
Engine
blob
Format
Raw Data
Handle
21418649
Attached To
rSPECMICP SpecMiCP / ReactMiCP
transport_program.hpp
View Options
/* =============================================================================
Copyright (c) 2014-2017 F. Georget <fabieng@princeton.edu> Princeton University
Copyright (c) 2017 F. Georget <fabien.georget@epfl.ch> EPFL
All rights reserved.
Redistribution and use in source and binary forms, with or without modification,
are permitted provided that the following conditions are met:
1. Redistributions of source code must retain the above copyright notice,
this list of conditions and the following disclaimer.
2. 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.
3. Neither the name of the copyright holder 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 HOLDER 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 REACTMICP_CHLORIDE_TRANSPORTPROGRAM_HPP
#define REACTMICP_CHLORIDE_TRANSPORTPROGRAM_HPP
#include "dfpm/solver/parabolic_program.hpp"
#include "variables.hpp"
#include "dfpm/meshes/mesh1dfwd.hpp"
#include "dfpm/types.hpp"
#include <functional>
#include "specmicp_common/return_code.hpp"
#include "specmicp_database/database_holder.hpp"
#include "reactmicp/systems/chloride/transport_program_struct.hpp"
#include "chem_stagger_base.hpp"
#include <iostream>
namespace
specmicp
{
struct
AdimensionalSystemSolution
;
struct
AdimensionalSystemConstraints
;
namespace
reactmicp
{
namespace
solver
{
class
ChemistryStaggerBase
;
}
// namespace solver
namespace
systems
{
namespace
chloride
{
//! \brief Concept class for a program
class
ChlorideProgram
:
public
dfpmsolver
::
ParabolicProgram
<
ChlorideProgram
>
,
public
database
::
DatabaseHolder
{
public
:
constexpr
static
bool
USE_CHEMISTRY_RATE
=
true
;
constexpr
static
bool
NO_CHEMISTRY_RATE
=
false
;
ChlorideProgram
(
ChlorideVariablesPtr
vars
,
std
::
vector
<
index_t
>
list_fixed_dof
);
//! \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_vars
->
ndof
();}
//! \brief Return the total number of degrees of freedom
index_t
get_tot_ndf
()
const
{
return
m_vars
->
ndof
()
*
m_mesh
->
nb_nodes
();}
//! \brief Return the id of the equation corresponding to the degree of freedom 'id_dof'
//!
//! Return 'no_equation' if no equation exist
index_t
id_equation
(
index_t
id_dof
)
const
{
return
m_ideq
(
id_dof
);}
void
initialize_timestep
(
scalar_t
dt
)
{
m_dt
=
dt
;
}
//! \brief Compute the residuals
void
compute_residuals
(
const
Vector
&
displacement
,
const
Vector
&
velocity
,
Vector
&
residual
,
bool
use_chemistry_rate
,
bool
use_grad_psi
);
//! \brief Compute the residuals
void
compute_residuals
(
const
Vector
&
displacement
,
const
Vector
&
velocity
,
Vector
&
residual
)
{
compute_residuals
(
displacement
,
velocity
,
residual
,
USE_CHEMISTRY_RATE
,
true
);
//std::cout << "Res : " << residual.norm() << std::endl;
}
//! \brief Compute the initial residuals
//!
//! By default, it calls the same residuals as compute_residuals
void
compute_residuals_0
(
const
Vector
&
displacement
,
const
Vector
&
velocity
,
Vector
&
residual
)
{
compute_residuals
(
displacement
,
velocity
,
residual
,
NO_CHEMISTRY_RATE
,
false
);
//std::cout << "Res 0 : " << residual.norm() << std::endl;
}
//! \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
)
{}
//! \brief Register Chemistry stagger
void
register_chem_stagger
(
std
::
shared_ptr
<
ChlorideChemistryStaggerBase
>
chem_stagger
);
//! \brief Compute residuals element
void
compute_element_residuals
(
index_t
element
,
const
Vector
&
edisplacement
,
const
Vector
&
evelocity
,
const
AdimensionalSystemSolution
&
adim_sol_0
,
const
AdimensionalSystemSolution
&
adim_sol_1
,
Vector
&
eresidual
,
bool
use_chemistry_rate
,
bool
use_grad_psi
);
//! \brief Set the approximation method
void
set_approx_method
(
AdimSolutionPerturbationMethod
method
)
{
m_approx_method
=
method
;
}
//! \brief Compute the transport rates
void
compute_transport_rate
(
scalar_t
dt
,
const
Vector
&
displacement
);
void
set_scaling
(
Vector
scaling_vector
)
{
m_scaling
=
scaling_vector
;
}
Vector
get_current
()
{
return
m_current
;
}
private
:
void
number_equations
(
std
::
vector
<
index_t
>&
list_fixed_dof
);
void
compute_element_jacobian
(
index_t
element
,
Vector
&
displacement
,
Vector
&
velocity
,
dfpm
::
list_triplet_t
&
jacobian
,
scalar_t
alphadt
);
ChlorideVariablesPtr
m_vars
;
mesh
::
Mesh1DPtr
m_mesh
;
index_t
m_neq
;
Vector
m_ideq
;
// ndof*node
scalar_t
m_dt
;
//!< Current timestep
std
::
shared_ptr
<
ChlorideChemistryStaggerBase
>
m_chem_stagger
;
AdimSolutionPerturbationMethod
m_approx_method
{
AdimSolutionPerturbationMethod
::
full
};
AdimensionalSystemSolution
perturb_adim_solution
(
index_t
element
,
index_t
enode
,
const
Vector
&
edisplacement
);
AdimensionalSystemSolution
perturb_adim_solution
(
index_t
node
,
const
Vector
&
displacement
);
// Scaling vector (size : ndof)
Vector
m_scaling
;
Vector
m_current
;
};
}
// namespace chloride
}
// namespace systems
}
// namespace reactmicp
}
// namespace specmicp
#endif
// REACTMICP_CHLORIDE_TRANSPORTPROGRAM_HPP
Event Timeline
Log In to Comment