Struct dzahui::solvers::fem::diffusion_solver::time_independent::DiffussionSolverTimeIndependent
source · pub struct DiffussionSolverTimeIndependent {
pub boundary_conditions: [f64; 2],
pub(crate) stiffness_matrix: Array2<f64>,
pub(crate) b_vector: Array1<f64>,
pub gauss_step: usize,
pub mu: f64,
pub b: f64,
}
Expand description
General Information
A diffusion solver with time-independence abstracts the equation: “- μu_xx + bu_x = 0” and contains boundary conditions along with mesh, “b” and “μ”
Fields
boundary_conditions
- Original boundary conditions (Only Dirichlet is supported for now).stiffness_matrix
- Left-side matrix of the resulting discrete equation.b_vector
- Right-side vector of the resulting discrete equation.gauss_step
- Precision of quadrature.mu
- First ot two needed constants.b
- Second of two needed constants.
Fields§
§boundary_conditions: [f64; 2]
§stiffness_matrix: Array2<f64>
§b_vector: Array1<f64>
§gauss_step: usize
§mu: f64
§b: f64
Implementations§
source§impl DiffussionSolverTimeIndependent
impl DiffussionSolverTimeIndependent
sourcepub fn new(
params: &DiffussionParamsTimeIndependent,
mesh: Vec<f64>,
gauss_step: usize
) -> Result<Self, Error>
pub fn new(
params: &DiffussionParamsTimeIndependent,
mesh: Vec<f64>,
gauss_step: usize
) -> Result<Self, Error>
Creates new instance
sourcepub fn gauss_legendre_integration(
boundary_conditions: [f64; 2],
mu: f64,
b: f64,
mesh: &Vec<f64>,
gauss_step: usize
) -> Result<(Array2<f64>, Array1<f64>), Error>
pub fn gauss_legendre_integration(
boundary_conditions: [f64; 2],
mu: f64,
b: f64,
mesh: &Vec<f64>,
gauss_step: usize
) -> Result<(Array2<f64>, Array1<f64>), Error>
General Information
First, it generates the basis for a solver from the linear basis constructor. Then the stiffnes matrix and vector b are generated based on linear basis integration via Gauss-Legendre and returned. Note that vector and matrix will have one on their diagonals’ boundaries and zero on other boundary elements to make boundary conditions permanent.
Parameters
boundary_conditions
- Conditions to guarantee system solution.mu
- Movement term.b
- Velocity term.mesh
- Vector of f64 representing a line.gauss_step
- How many nodes will be calculated for a given integration.
Returns
A tuple with both the stiffness matrix and the vector b.