Struct dzahui::solvers::fem::diffusion_solver::time_dependent::DiffussionSolverTimeDependent
source · pub struct DiffussionSolverTimeDependent {
pub boundary_conditions: [f64; 2],
pub(crate) stiffness_matrix: Array2<f64>,
pub initial_conditions: Vec<f64>,
pub(crate) mass_matrix: Array2<f64>,
pub integration_step: usize,
pub(crate) state: Array1<f64>,
pub mu: f64,
pub b: f64,
}
Expand description
General Information
A diffusion solver with time-dependence abstracts the equation: “u_t - μu_xx + bu_x = 0” and contains boundary conditions, initial conditions, mesh, “stiffness_matrix” and “μ”.
Fields
boundary_conditions
- Boundary conditions (Only dirichlet is supported for now, Neumann is being worked on)stiffness_matrix
- Matrix of elements that is multiplied by timeinitial_conditions
- Every internal point needs an initial condition to advance the solution in timemass_matrix
- A matrix that pertains only to elements that are not multiplied by timeintegration_step
- Amount of terms to sum over to make the integralstate
- The state of every point at time tmu
- First ot two needed constantsb
- Second of two needed constants
Fields§
§boundary_conditions: [f64; 2]
§stiffness_matrix: Array2<f64>
§initial_conditions: Vec<f64>
§mass_matrix: Array2<f64>
§integration_step: usize
§state: Array1<f64>
§mu: f64
§b: f64
Implementations§
source§impl DiffussionSolverTimeDependent
impl DiffussionSolverTimeDependent
sourcepub fn new(
params: &DiffussionParamsTimeDependent,
mesh: Vec<f64>,
integration_step: usize
) -> Result<Self, Error>
pub fn new(
params: &DiffussionParamsTimeDependent,
mesh: Vec<f64>,
integration_step: usize
) -> Result<Self, Error>
Creates new instance checking initial conditions are the size they should be.
sourcefn gauss_legendre_integration(
mu: f64,
b: f64,
mesh: &Vec<f64>,
gauss_step: usize
) -> Result<(Array2<f64>, Array2<f64>), Error>
fn gauss_legendre_integration(
mu: f64,
b: f64,
mesh: &Vec<f64>,
gauss_step: usize
) -> Result<(Array2<f64>, Array2<f64>), Error>
General Information
Compĺete integration of linear basis to obtain mass matrix and stiffness matrix. Corners of every element have special values to attone for boundary conditions being constant. Matrices serve to solve the resulting problem: M(u_ti+1) = M(u_ti) + S(delta_t * u_ti) where M is mass matrix and S is stiffness matrix.
Parameters
mu
- First of two terms to solve equationb
- Second of two terms to solve equationmesh
- Vector of f64 representing a meshgauss_step
- Amount of nodes to compute for integration.