use oxiphysics_core::math::Vec3;
use crate::assembly::{assemble_load_vector, assemble_stiffness};
use crate::boundary::{DirichletBc, NeumannBc, apply_dirichlet, apply_neumann};
use crate::constitutive::LinearElasticMaterial;
use crate::element::LinearTetrahedron;
use crate::mesh::TetrahedralMesh;
use crate::solvers::{jacobi_preconditioner, preconditioned_cg};
#[derive(Debug, Clone, PartialEq)]
pub enum FemError {
SingularSystem,
InvalidInput(String),
SolverNotConverged,
}
impl std::fmt::Display for FemError {
fn fmt(&self, f: &mut std::fmt::Formatter<'_>) -> std::fmt::Result {
match self {
FemError::SingularSystem => write!(f, "singular or ill-conditioned system"),
FemError::InvalidInput(msg) => write!(f, "invalid input: {msg}"),
FemError::SolverNotConverged => write!(f, "iterative solver did not converge"),
}
}
}
impl std::error::Error for FemError {}
#[derive(Debug, Clone)]
pub struct FemResult {
pub displacements: Vec<Vec3>,
pub stress: Vec<[f64; 6]>,
pub strain: Vec<[f64; 6]>,
}
pub struct LinearStaticAnalysis {
pub max_iter: usize,
pub tolerance: f64,
}
impl Default for LinearStaticAnalysis {
fn default() -> Self {
Self {
max_iter: 10000,
tolerance: 1e-10,
}
}
}
impl LinearStaticAnalysis {
pub fn new() -> Self {
Self::default()
}
pub fn solve(
&self,
mesh: &TetrahedralMesh,
material: &LinearElasticMaterial,
dirichlet_bcs: &[DirichletBc],
neumann_bcs: &[NeumannBc],
body_force: &Vec3,
) -> FemResult {
let ndof = mesh.num_nodes() * 3;
let mut k = assemble_stiffness(mesh, material);
let mut f = assemble_load_vector(mesh, body_force);
apply_neumann(&mut f, neumann_bcs);
apply_dirichlet(&mut k, &mut f, dirichlet_bcs);
let precond = jacobi_preconditioner(&k);
let u = preconditioned_cg(&k, &f, &precond, self.max_iter, self.tolerance);
let num_nodes = mesh.num_nodes();
let mut displacements = Vec::with_capacity(num_nodes);
for i in 0..num_nodes {
displacements.push(Vec3::new(u[i * 3], u[i * 3 + 1], u[i * 3 + 2]));
}
let d_matrix = material.constitutive_matrix();
let num_elements = mesh.num_elements();
let mut stress = Vec::with_capacity(num_elements);
let mut strain = Vec::with_capacity(num_elements);
for elem in &mesh.elements {
let nodes = [
mesh.nodes[elem[0]],
mesh.nodes[elem[1]],
mesh.nodes[elem[2]],
mesh.nodes[elem[3]],
];
let (b_mat, _vol) = LinearTetrahedron::b_matrix(&nodes);
let mut u_e = [0.0f64; 12];
for (local, &global) in elem.iter().enumerate() {
u_e[local * 3] = u[global * 3];
u_e[local * 3 + 1] = u[global * 3 + 1];
u_e[local * 3 + 2] = u[global * 3 + 2];
}
let mut eps = [0.0f64; 6];
for i in 0..6 {
for j in 0..12 {
eps[i] += b_mat[i * 12 + j] * u_e[j];
}
}
let mut sig = [0.0f64; 6];
for i in 0..6 {
for j in 0..6 {
sig[i] += d_matrix[i][j] * eps[j];
}
}
strain.push(eps);
stress.push(sig);
}
let _ = ndof;
FemResult {
displacements,
stress,
strain,
}
}
}
pub struct StaticAnalysis {
pub n_dofs: usize,
pub stiffness_matrix: Vec<f64>,
pub force_vector: Vec<f64>,
}
impl StaticAnalysis {
pub fn new(n_dofs: usize) -> Self {
Self {
n_dofs,
stiffness_matrix: vec![0.0; n_dofs * n_dofs],
force_vector: vec![0.0; n_dofs],
}
}
pub fn assemble_global_stiffness(&mut self, elements: &[(Vec<usize>, Vec<f64>)]) {
let n = self.n_dofs;
for (dofs, ke) in elements {
let ndof_e = dofs.len();
for i in 0..ndof_e {
for j in 0..ndof_e {
let gi = dofs[i];
let gj = dofs[j];
self.stiffness_matrix[gi * n + gj] += ke[i * ndof_e + j];
}
}
}
}
pub fn solve(&self) -> Result<Vec<f64>, FemError> {
let n = self.n_dofs;
if n == 0 {
return Ok(Vec::new());
}
let mut aug: Vec<f64> = Vec::with_capacity(n * (n + 1));
for i in 0..n {
for j in 0..n {
aug.push(self.stiffness_matrix[i * n + j]);
}
aug.push(self.force_vector[i]);
}
let cols = n + 1;
for col in 0..n {
let mut max_val = aug[col * cols + col].abs();
let mut max_row = col;
for row in (col + 1)..n {
let v = aug[row * cols + col].abs();
if v > max_val {
max_val = v;
max_row = row;
}
}
if max_val < 1e-14 {
return Err(FemError::SingularSystem);
}
if max_row != col {
for k in 0..cols {
aug.swap(col * cols + k, max_row * cols + k);
}
}
let pivot = aug[col * cols + col];
for row in (col + 1)..n {
let factor = aug[row * cols + col] / pivot;
for k in col..cols {
let sub = factor * aug[col * cols + k];
aug[row * cols + k] -= sub;
}
}
}
let mut u = vec![0.0; n];
for i in (0..n).rev() {
let mut sum = aug[i * cols + n];
for j in (i + 1)..n {
sum -= aug[i * cols + j] * u[j];
}
u[i] = sum / aug[i * cols + i];
}
Ok(u)
}
pub fn compute_strains(
&self,
displacements: &[f64],
elements: &[(Vec<usize>, Vec<f64>)],
) -> Vec<[f64; 6]> {
let mut strains = Vec::with_capacity(elements.len());
for (dofs, b_mat) in elements {
let ndof_e = dofs.len();
let mut u_e = vec![0.0; ndof_e];
for (k, &g) in dofs.iter().enumerate() {
u_e[k] = displacements[g];
}
let mut eps = [0.0f64; 6];
for i in 0..6 {
for j in 0..ndof_e {
eps[i] += b_mat[i * ndof_e + j] * u_e[j];
}
}
strains.push(eps);
}
strains
}
pub fn compute_stresses(strains: &[[f64; 6]], d_matrix: &[f64; 36]) -> Vec<[f64; 6]> {
strains
.iter()
.map(|eps| {
let mut sig = [0.0f64; 6];
for i in 0..6 {
for j in 0..6 {
sig[i] += d_matrix[i * 6 + j] * eps[j];
}
}
sig
})
.collect()
}
}
pub struct ModalAnalysis {
pub max_iter: usize,
pub tol: f64,
}
impl Default for ModalAnalysis {
fn default() -> Self {
Self {
max_iter: 500,
tol: 1e-10,
}
}
}
impl ModalAnalysis {
pub fn new() -> Self {
Self::default()
}
fn dot(a: &[f64], b: &[f64]) -> f64 {
a.iter().zip(b.iter()).map(|(x, y)| x * y).sum()
}
fn matvec_dense(a: &[f64], x: &[f64], n: usize) -> Vec<f64> {
let mut y = vec![0.0; n];
for i in 0..n {
for j in 0..n {
y[i] += a[i * n + j] * x[j];
}
}
y
}
fn solve_dense(a: &[f64], b: &[f64], n: usize) -> Option<Vec<f64>> {
let cols = n + 1;
let mut aug = vec![0.0; n * cols];
for i in 0..n {
for j in 0..n {
aug[i * cols + j] = a[i * n + j];
}
aug[i * cols + n] = b[i];
}
for col in 0..n {
let mut max_val = aug[col * cols + col].abs();
let mut max_row = col;
for row in (col + 1)..n {
let v = aug[row * cols + col].abs();
if v > max_val {
max_val = v;
max_row = row;
}
}
if max_val < 1e-15 {
return None;
}
if max_row != col {
for k in 0..cols {
aug.swap(col * cols + k, max_row * cols + k);
}
}
let pivot = aug[col * cols + col];
for row in (col + 1)..n {
let factor = aug[row * cols + col] / pivot;
for k in col..cols {
let sub = factor * aug[col * cols + k];
aug[row * cols + k] -= sub;
}
}
}
let mut x = vec![0.0; n];
for i in (0..n).rev() {
let mut s = aug[i * cols + n];
for j in (i + 1)..n {
s -= aug[i * cols + j] * x[j];
}
x[i] = s / aug[i * cols + i];
}
Some(x)
}
pub fn compute_natural_frequencies(
k: &[f64],
m: &[f64],
n: usize,
n_modes: usize,
max_iter: usize,
) -> Vec<f64> {
let mut freqs = Vec::with_capacity(n_modes);
let mut found_modes: Vec<Vec<f64>> = Vec::new();
for mode_idx in 0..n_modes {
let mut q: Vec<f64> = (0..n)
.map(|i| {
let v = ((i + mode_idx * 13 + 1) as f64 * 0.7).sin();
if i == mode_idx % n { v + 1.0 } else { v }
})
.collect();
let mq = Self::matvec_dense(m, &q, n);
let nm = Self::dot(&q, &mq).max(0.0).sqrt();
if nm > 1e-60 {
for qi in q.iter_mut() {
*qi /= nm;
}
}
let mut omega_sq = 0.0_f64;
for _iter in 0..max_iter {
for phi in &found_modes {
let mphi = Self::matvec_dense(m, phi, n);
let c = Self::dot(&q, &mphi);
for i in 0..n {
q[i] -= c * phi[i];
}
}
let kq = Self::matvec_dense(k, &q, n);
let mq2 = Self::matvec_dense(m, &q, n);
let num = Self::dot(&q, &kq);
let den = Self::dot(&q, &mq2);
let rq = if den.abs() > 1e-60 { num / den } else { 0.0 };
let mut km = vec![0.0; n * n];
for i in 0..n {
for j in 0..n {
km[i * n + j] = k[i * n + j] - rq * m[i * n + j];
}
}
let rhs = Self::matvec_dense(m, &q, n);
let y_opt = Self::solve_dense(&km, &rhs, n);
let y = match y_opt {
Some(y) => y,
None => {
let mut rhs2 = Self::matvec_dense(m, &q, n);
if let Some(sol) = Self::solve_dense(k, &rhs2, n) {
rhs2 = sol;
}
rhs2
}
};
let mut yd = y.clone();
for phi in &found_modes {
let mphi = Self::matvec_dense(m, phi, n);
let c = Self::dot(&yd, &mphi);
for i in 0..n {
yd[i] -= c * phi[i];
}
}
let my = Self::matvec_dense(m, &yd, n);
let nm2 = Self::dot(&yd, &my).max(0.0).sqrt();
if nm2 < 1e-60 {
break;
}
let q_new: Vec<f64> = yd.iter().map(|yi| yi / nm2).collect();
let kqn = Self::matvec_dense(k, &q_new, n);
let mqn = Self::matvec_dense(m, &q_new, n);
let num_new = Self::dot(&q_new, &kqn);
let den_new = Self::dot(&q_new, &mqn);
let rq_new = if den_new.abs() > 1e-60 {
num_new / den_new
} else {
0.0
};
let change = (rq_new - omega_sq).abs() / (rq_new.abs() + 1e-60);
q = q_new;
omega_sq = rq_new;
if change < 1e-10 && _iter > 0 {
break;
}
}
freqs.push(omega_sq.max(0.0).sqrt());
found_modes.push(q);
}
freqs
}
pub fn compute_mode_shapes(
k: &[f64],
m: &[f64],
n: usize,
omega: f64,
max_iter: usize,
) -> Vec<f64> {
let shift = omega * omega;
let mut km = vec![0.0; n * n];
for i in 0..n {
for j in 0..n {
km[i * n + j] = k[i * n + j] - shift * m[i * n + j];
}
}
let mut q: Vec<f64> = (0..n).map(|i| ((i + 1) as f64 * 0.5).sin()).collect();
let nrm: f64 = q.iter().map(|x| x * x).sum::<f64>().sqrt();
if nrm > 1e-60 {
for qi in q.iter_mut() {
*qi /= nrm;
}
}
for _ in 0..max_iter {
let rhs = Self::matvec_dense(m, &q, n);
let y = match Self::solve_dense(&km, &rhs, n) {
Some(y) => y,
None => break,
};
let ny: f64 = y.iter().map(|x| x * x).sum::<f64>().sqrt();
if ny < 1e-60 {
break;
}
let q_new: Vec<f64> = y.iter().map(|yi| yi / ny).collect();
let change: f64 = q
.iter()
.zip(q_new.iter())
.map(|(a, b)| (a - b).abs())
.sum::<f64>()
/ (n as f64);
q = q_new;
if change < 1e-12 {
break;
}
}
let mq = Self::matvec_dense(m, &q, n);
let mn: f64 = Self::dot(&q, &mq).max(0.0).sqrt();
if mn > 1e-60 {
for qi in q.iter_mut() {
*qi /= mn;
}
}
q
}
}
pub struct DynamicAnalysis {
pub m: Vec<f64>,
pub c: Vec<f64>,
pub k: Vec<f64>,
pub n: usize,
}
impl DynamicAnalysis {
pub fn new(m: Vec<f64>, c: Vec<f64>, k: Vec<f64>, n: usize) -> Self {
Self { m, c, k, n }
}
fn mv(a: &[f64], x: &[f64], n: usize) -> Vec<f64> {
let mut y = vec![0.0; n];
for i in 0..n {
for j in 0..n {
y[i] += a[i * n + j] * x[j];
}
}
y
}
fn solve(a: &[f64], b: &[f64], n: usize) -> Vec<f64> {
let cols = n + 1;
let mut aug = vec![0.0; n * cols];
for i in 0..n {
for j in 0..n {
aug[i * cols + j] = a[i * n + j];
}
aug[i * cols + n] = b[i];
}
for col in 0..n {
let mut max_val = aug[col * cols + col].abs();
let mut max_row = col;
for row in (col + 1)..n {
let v = aug[row * cols + col].abs();
if v > max_val {
max_val = v;
max_row = row;
}
}
if max_val < 1e-15 {
continue;
}
if max_row != col {
for k in 0..cols {
aug.swap(col * cols + k, max_row * cols + k);
}
}
let piv = aug[col * cols + col];
for row in (col + 1)..n {
let fac = aug[row * cols + col] / piv;
for k in col..cols {
let sub = fac * aug[col * cols + k];
aug[row * cols + k] -= sub;
}
}
}
let mut x = vec![0.0; n];
for i in (0..n).rev() {
let mut s = aug[i * cols + n];
for j in (i + 1)..n {
s -= aug[i * cols + j] * x[j];
}
x[i] = s / aug[i * cols + i];
}
x
}
pub fn newmark_beta_step(
&self,
u: &[f64],
v: &[f64],
a: &[f64],
f_ext: &[f64],
dt: f64,
beta: f64,
gamma: f64,
) -> (Vec<f64>, Vec<f64>, Vec<f64>) {
let n = self.n;
let u_pred: Vec<f64> = (0..n)
.map(|i| u[i] + dt * v[i] + dt * dt * (0.5 - beta) * a[i])
.collect();
let v_pred: Vec<f64> = (0..n).map(|i| v[i] + dt * (1.0 - gamma) * a[i]).collect();
let c1 = 1.0 / (beta * dt * dt);
let c2 = gamma / (beta * dt);
let mut k_eff = vec![0.0; n * n];
for i in 0..n {
for j in 0..n {
k_eff[i * n + j] =
c1 * self.m[i * n + j] + c2 * self.c[i * n + j] + self.k[i * n + j];
}
}
let ku_pred = Self::mv(&self.k, &u_pred, n);
let cv_pred = Self::mv(&self.c, &v_pred, n);
let rhs: Vec<f64> = (0..n).map(|i| f_ext[i] - ku_pred[i] - cv_pred[i]).collect();
let a_new = Self::solve(&k_eff, &rhs, n);
let u_new: Vec<f64> = (0..n)
.map(|i| u_pred[i] + beta * dt * dt * a_new[i])
.collect();
let v_new: Vec<f64> = (0..n).map(|i| v_pred[i] + gamma * dt * a_new[i]).collect();
(u_new, v_new, a_new)
}
pub fn wilson_theta_step(
&self,
u: &[f64],
v: &[f64],
a: &[f64],
f_ext: &[f64],
dt: f64,
theta: f64,
) -> (Vec<f64>, Vec<f64>, Vec<f64>) {
let n = self.n;
let tau = theta * dt;
let c1 = 6.0 / (tau * tau);
let c2 = 3.0 / tau;
let mut k_eff = vec![0.0; n * n];
for i in 0..n {
for j in 0..n {
k_eff[i * n + j] =
c1 * self.m[i * n + j] + c2 * self.c[i * n + j] + self.k[i * n + j];
}
}
let f_tau = f_ext.to_vec();
let mu = Self::mv(&self.m, u, n);
let mv_curr = Self::mv(&self.m, v, n);
let ma = Self::mv(&self.m, a, n);
let cu = Self::mv(&self.c, u, n);
let cv_curr = Self::mv(&self.c, v, n);
let rhs: Vec<f64> = (0..n)
.map(|i| {
f_tau[i] + c1 * mu[i] + c2 * mv_curr[i] + 2.0 * ma[i]
- self.k[i * n..i * n + n]
.iter()
.zip(u.iter())
.map(|(kij, uj)| kij * uj)
.sum::<f64>()
+ (c2 * cu[i] + cv_curr[i])
})
.collect();
let a_tau = Self::solve(&k_eff, &rhs, n);
let a_new: Vec<f64> = (0..n).map(|i| a[i] + (a_tau[i] - a[i]) / theta).collect();
let v_new: Vec<f64> = (0..n)
.map(|i| v[i] + dt / 2.0 * (a[i] + a_new[i]))
.collect();
let u_new: Vec<f64> = (0..n)
.map(|i| u[i] + dt * v[i] + dt * dt / 6.0 * (2.0 * a[i] + a_new[i]))
.collect();
(u_new, v_new, a_new)
}
}
#[derive(Debug, Clone, Copy)]
pub struct StressState {
pub sigma: [f64; 6],
}
impl StressState {
pub fn from_voigt(sigma: [f64; 6]) -> Self {
Self { sigma }
}
pub fn new(sxx: f64, syy: f64, szz: f64, txy: f64, tyz: f64, txz: f64) -> Self {
Self {
sigma: [sxx, syy, szz, txy, tyz, txz],
}
}
pub fn hydrostatic(&self) -> f64 {
(self.sigma[0] + self.sigma[1] + self.sigma[2]) / 3.0
}
pub fn i1(&self) -> f64 {
self.sigma[0] + self.sigma[1] + self.sigma[2]
}
pub fn j2(&self) -> f64 {
let p = self.hydrostatic();
let sx = self.sigma[0] - p;
let sy = self.sigma[1] - p;
let sz = self.sigma[2] - p;
let txy = self.sigma[3];
let tyz = self.sigma[4];
let txz = self.sigma[5];
0.5 * (sx * sx + sy * sy + sz * sz + 2.0 * (txy * txy + tyz * tyz + txz * txz))
}
pub fn j3(&self) -> f64 {
let p = self.hydrostatic();
let sx = self.sigma[0] - p;
let sy = self.sigma[1] - p;
let sz = self.sigma[2] - p;
let txy = self.sigma[3];
let tyz = self.sigma[4];
let txz = self.sigma[5];
sx * (sy * sz - tyz * tyz) - txy * (txy * sz - tyz * txz) + txz * (txy * tyz - sy * txz)
}
pub fn von_mises(&self) -> f64 {
(3.0 * self.j2()).max(0.0).sqrt()
}
pub fn tresca_stress(&self) -> f64 {
let [s1, _s2, s3] = self.principal_stresses();
s1 - s3
}
pub fn principal_stresses(&self) -> [f64; 3] {
let p = self.hydrostatic();
let j2 = self.j2();
let j3 = self.j3();
let r = (j2 / 3.0).max(0.0).sqrt();
if r < 1e-30 {
return [p, p, p];
}
let cos_arg = (j3 / (2.0 * r * r * r)).clamp(-1.0, 1.0);
let theta = cos_arg.acos() / 3.0;
use std::f64::consts::PI;
let s1 = p + 2.0 * r * theta.cos();
let s2 = p + 2.0 * r * (theta - 2.0 * PI / 3.0).cos();
let s3 = p + 2.0 * r * (theta - 4.0 * PI / 3.0).cos();
let mut arr = [s1, s2, s3];
arr.sort_by(|a, b| b.partial_cmp(a).unwrap_or(std::cmp::Ordering::Equal));
arr
}
pub fn lode_angle_deg(&self) -> f64 {
let j2 = self.j2();
let j3 = self.j3();
if j2 < 1e-30 {
return 0.0;
}
let cos_arg = (3.0_f64.sqrt() * 3.0 * j3 / (2.0 * j2.powf(1.5))).clamp(-1.0, 1.0);
let theta_rad = cos_arg.acos() / 3.0;
theta_rad.to_degrees()
}
pub fn safety_factor_von_mises(&self, yield_stress: f64) -> f64 {
let vm = self.von_mises();
if vm < 1e-30 {
return f64::INFINITY;
}
yield_stress / vm
}
pub fn safety_factor_tresca(&self, yield_stress: f64) -> f64 {
let tresca = self.tresca_stress();
if tresca < 1e-30 {
return f64::INFINITY;
}
yield_stress / tresca
}
pub fn drucker_prager(&self, alpha: f64) -> f64 {
alpha * self.i1() + self.j2().max(0.0).sqrt()
}
}
#[derive(Debug, Clone)]
pub struct FatigueLifeEstimator {
pub sigma_uts: f64,
pub n_ref: f64,
pub endurance_limit: f64,
}
impl FatigueLifeEstimator {
pub fn new(sigma_uts: f64, n_ref: f64, endurance_limit: f64) -> Self {
assert!(sigma_uts > 0.0);
assert!(n_ref > 0.0);
assert!(endurance_limit >= 0.0 && endurance_limit < sigma_uts);
Self {
sigma_uts,
n_ref,
endurance_limit,
}
}
pub fn cycles_to_failure(&self, sigma_a: f64) -> f64 {
if sigma_a <= self.endurance_limit {
return f64::INFINITY;
}
if sigma_a >= self.sigma_uts {
return 1.0;
}
let m = self.basquin_exponent();
self.n_ref * (self.sigma_uts / sigma_a).powf(m)
}
pub fn basquin_exponent(&self) -> f64 {
self.n_ref.log10() / 2.0_f64.log10()
}
pub fn goodman_corrected_amplitude(&self, sigma_a: f64, sigma_mean: f64) -> f64 {
let denom = 1.0 - sigma_mean / self.sigma_uts;
if denom <= 0.0 {
return f64::INFINITY;
}
sigma_a / denom
}
}
pub fn miner_damage(cycle_ratios: &[f64]) -> f64 {
cycle_ratios.iter().sum()
}
pub fn compute_element_strains(
displacements: &[f64],
element_dofs: &[Vec<usize>],
b_matrices: &[Vec<f64>],
) -> Vec<[f64; 6]> {
element_dofs
.iter()
.zip(b_matrices.iter())
.map(|(dofs, b)| {
let mut eps = [0.0f64; 6];
let n_dof = dofs.len();
for row in 0..6 {
let mut sum = 0.0;
for col in 0..n_dof {
if dofs[col] < displacements.len() {
sum += b[row * n_dof + col] * displacements[dofs[col]];
}
}
eps[row] = sum;
}
eps
})
.collect()
}
pub fn von_mises_field(stresses: &[[f64; 6]]) -> Vec<f64> {
stresses
.iter()
.map(|&s| StressState::from_voigt(s).von_mises())
.collect()
}
pub fn safety_factor_field(stresses: &[[f64; 6]], yield_stress: f64) -> Vec<f64> {
stresses
.iter()
.map(|&s| StressState::from_voigt(s).safety_factor_von_mises(yield_stress))
.collect()
}
pub fn smooth_field(values: &[f64], neighbors: &[Vec<usize>]) -> Vec<f64> {
values
.iter()
.enumerate()
.map(|(i, &v)| {
let nbrs = &neighbors[i];
if nbrs.is_empty() {
return v;
}
let sum: f64 = nbrs
.iter()
.map(|&j| if j < values.len() { values[j] } else { 0.0 })
.sum();
(v + sum) / (1 + nbrs.len()) as f64
})
.collect()
}
#[derive(Debug, Clone, Copy)]
pub struct RayleighDamping {
pub alpha: f64,
pub beta: f64,
}
impl RayleighDamping {
pub fn from_two_modes(omega_i: f64, omega_j: f64, zeta_i: f64, zeta_j: f64) -> Self {
let det = (omega_j * omega_j - omega_i * omega_i) / (2.0 * omega_i * omega_j);
if det.abs() < 1e-30 {
return Self {
alpha: 0.0,
beta: 0.0,
};
}
let alpha = 2.0 * omega_i * omega_j * (zeta_i * omega_j - zeta_j * omega_i)
/ (omega_j * omega_j - omega_i * omega_i);
let beta =
2.0 * (zeta_j * omega_j - zeta_i * omega_i) / (omega_j * omega_j - omega_i * omega_i);
Self { alpha, beta }
}
pub fn damping_ratio(&self, omega: f64) -> f64 {
if omega < 1e-30 {
return 0.0;
}
self.alpha / (2.0 * omega) + self.beta * omega / 2.0
}
pub fn assemble(&self, mass: &[f64], stiffness: &[f64], n: usize) -> Vec<f64> {
assert_eq!(mass.len(), n * n);
assert_eq!(stiffness.len(), n * n);
(0..n * n)
.map(|k| self.alpha * mass[k] + self.beta * stiffness[k])
.collect()
}
}
#[cfg(test)]
mod tests {
use super::*;
#[test]
fn test_static_solve_2dof() {
let mut sa = StaticAnalysis::new(2);
sa.stiffness_matrix = vec![4.0, -1.0, -1.0, 3.0];
sa.force_vector = vec![1.0, 2.0];
let u = sa.solve().unwrap();
let expected0 = 5.0 / 11.0;
let expected1 = 9.0 / 11.0;
assert!(
(u[0] - expected0).abs() < 1e-12,
"u[0] = {}, expected {expected0}",
u[0]
);
assert!(
(u[1] - expected1).abs() < 1e-12,
"u[1] = {}, expected {expected1}",
u[1]
);
}
#[test]
fn test_static_solve_singular() {
let mut sa = StaticAnalysis::new(2);
sa.stiffness_matrix = vec![1.0, 1.0, 1.0, 1.0];
sa.force_vector = vec![1.0, 1.0];
let result = sa.solve();
assert!(result.is_err(), "expected SingularSystem error");
}
#[test]
fn test_static_assemble_and_solve() {
let mut sa = StaticAnalysis::new(2);
let elements = vec![(vec![0usize, 1], vec![2.0, -2.0, -2.0, 2.0])];
sa.assemble_global_stiffness(&elements);
sa.stiffness_matrix[0] = 1.0;
sa.stiffness_matrix[1] = 0.0;
sa.stiffness_matrix[2] = 0.0;
sa.force_vector[0] = 0.0;
sa.force_vector[1] = 1.0; let u = sa.solve().unwrap();
assert!((u[0]).abs() < 1e-12, "u[0] should be 0, got {}", u[0]);
assert!(u[1] > 0.0, "u[1] should be positive, got {}", u[1]);
}
#[test]
fn test_compute_strains() {
let sa = StaticAnalysis::new(1);
let b_mat: Vec<f64> = {
let mut v = vec![0.0; 6];
v[0] = 1.0;
v
};
let elements = vec![(vec![0usize], b_mat)];
let u = vec![0.5_f64];
let strains = sa.compute_strains(&u, &elements);
assert_eq!(strains.len(), 1);
assert!((strains[0][0] - 0.5).abs() < 1e-15);
for (k, &val) in strains[0].iter().enumerate().skip(1) {
assert_eq!(val, 0.0, "strains[0][{k}] should be 0");
}
}
#[test]
fn test_compute_stresses_identity_d() {
let mut d = [0.0f64; 36];
for i in 0..6 {
d[i * 6 + i] = 1.0;
}
let eps = [1.0, 2.0, 3.0, 4.0, 5.0, 6.0];
let stresses = StaticAnalysis::compute_stresses(&[eps], &d);
assert_eq!(stresses[0], eps);
}
#[test]
fn test_modal_1dof() {
let k_val = 400.0f64;
let m_val = 1.0f64;
let k = vec![k_val];
let m = vec![m_val];
let _ma = ModalAnalysis::default();
let freqs = ModalAnalysis::compute_natural_frequencies(&k, &m, 1, 1, 200);
let expected = k_val.sqrt();
assert!(
(freqs[0] - expected).abs() / expected < 1e-5,
"ω = {}, expected {expected}",
freqs[0]
);
}
#[test]
fn test_modal_2dof_diagonal() {
let k = vec![1.0, 0.0, 0.0, 9.0];
let m = vec![1.0, 0.0, 0.0, 1.0];
let _ma = ModalAnalysis::default();
let mut freqs = ModalAnalysis::compute_natural_frequencies(&k, &m, 2, 2, 500);
freqs.sort_by(|a, b| a.partial_cmp(b).unwrap());
assert!((freqs[0] - 1.0).abs() < 0.05, "ω₁ = {}", freqs[0]);
assert!((freqs[1] - 3.0).abs() < 0.05, "ω₂ = {}", freqs[1]);
}
#[test]
fn test_mode_shape_1dof() {
let k = vec![25.0_f64];
let m = vec![1.0_f64];
let phi = ModalAnalysis::compute_mode_shapes(&k, &m, 1, 5.0, 100);
let norm: f64 = phi.iter().map(|x| x * x).sum::<f64>().sqrt();
assert!((norm - 1.0).abs() < 1e-6, "mass norm = {norm}");
}
#[test]
fn test_newmark_free_vibration() {
let m_val = 1.0_f64;
let k_val = 100.0_f64;
let dt = 0.001;
let beta = 0.25;
let gamma = 0.5;
let da = DynamicAnalysis::new(vec![m_val], vec![0.0], vec![k_val], 1);
let mut u = vec![1.0_f64];
let mut v = vec![0.0_f64];
let mut a = vec![-k_val / m_val];
for _ in 0..100 {
let (un, vn, an) = da.newmark_beta_step(&u, &v, &a, &[0.0], dt, beta, gamma);
u = un;
v = vn;
a = an;
}
assert!(u[0].abs() <= 1.01, "u should be bounded, got {}", u[0]);
}
#[test]
fn test_newmark_static_equilibrium() {
let da = DynamicAnalysis::new(vec![1.0], vec![0.1], vec![10.0], 1);
let mut u = vec![0.0_f64];
let mut v = vec![0.0_f64];
let mut a = vec![0.0_f64];
let f = vec![5.0_f64]; for _ in 0..5000 {
let (un, vn, an) = da.newmark_beta_step(&u, &v, &a, &f, 0.01, 0.25, 0.5);
u = un;
v = vn;
a = an;
}
assert!(u[0] > 0.0, "u should be positive, got {}", u[0]);
}
#[test]
fn test_wilson_theta_step_displacement_sign() {
let da = DynamicAnalysis::new(vec![1.0], vec![0.0], vec![1.0], 1);
let u = vec![0.0_f64];
let v = vec![0.0_f64];
let a = vec![0.0_f64];
let (u_new, _, _) = da.wilson_theta_step(&u, &v, &a, &[1.0], 0.01, 1.4);
assert!(
u_new[0] >= 0.0,
"displacement should be non-negative, got {}",
u_new[0]
);
}
#[test]
fn test_cantilever_deflection() {
let length: f64 = 0.5;
let width: f64 = 0.1;
let height: f64 = 0.1;
let e: f64 = 200.0e9;
let nu: f64 = 0.3;
let p: f64 = -500.0;
let mesh = TetrahedralMesh::generate_beam(length, width, height, 4, 1, 1);
let material = LinearElasticMaterial::new(e, nu);
let mut dirichlet_bcs = Vec::new();
for (i, node) in mesh.nodes.iter().enumerate() {
if node.x.abs() < 1e-10 {
dirichlet_bcs.push(DirichletBc::new(i, 0, 0.0));
dirichlet_bcs.push(DirichletBc::new(i, 1, 0.0));
dirichlet_bcs.push(DirichletBc::new(i, 2, 0.0));
}
}
let free_end_nodes: Vec<usize> = mesh
.nodes
.iter()
.enumerate()
.filter(|(_, n)| (n.x - length).abs() < 1e-10)
.map(|(i, _)| i)
.collect();
let num_free = free_end_nodes.len() as f64;
let neumann_bcs: Vec<NeumannBc> = free_end_nodes
.iter()
.map(|&i| NeumannBc::new(i, [0.0, p / num_free, 0.0]))
.collect();
let analysis = LinearStaticAnalysis {
max_iter: 20000,
tolerance: 1e-10,
};
let result = analysis.solve(
&mesh,
&material,
&dirichlet_bcs,
&neumann_bcs,
&Vec3::new(0.0, 0.0, 0.0),
);
let max_y_disp = free_end_nodes
.iter()
.map(|&i| result.displacements[i].y)
.fold(f64::INFINITY, f64::min);
assert!(
max_y_disp < 0.0,
"Tip deflection should be downward (negative y), got {max_y_disp:.6e}"
);
}
#[test]
fn test_hydrostatic_pure_hydrostatic() {
let sigma = [1.0, 1.0, 1.0, 0.0, 0.0, 0.0];
let s = StressState::from_voigt(sigma);
assert!(
(s.hydrostatic() - 1.0).abs() < 1e-12,
"hydrostatic={}",
s.hydrostatic()
);
}
#[test]
fn test_von_mises_zero_for_hydrostatic() {
let sigma = [5.0, 5.0, 5.0, 0.0, 0.0, 0.0];
let s = StressState::from_voigt(sigma);
assert!(
s.von_mises() < 1e-10,
"von Mises should be 0 for hydrostatic: {}",
s.von_mises()
);
}
#[test]
fn test_von_mises_uniaxial() {
let s0 = 300e6;
let sigma = [s0, 0.0, 0.0, 0.0, 0.0, 0.0];
let s = StressState::from_voigt(sigma);
assert!(
(s.von_mises() - s0).abs() < 1.0,
"von Mises uniaxial: {}",
s.von_mises()
);
}
#[test]
fn test_principal_stresses_sorted() {
let sigma = [200e6, 100e6, 50e6, 0.0, 0.0, 0.0];
let s = StressState::from_voigt(sigma);
let [s1, s2, s3] = s.principal_stresses();
assert!(
s1 >= s2 && s2 >= s3,
"principal stresses not sorted: {s1} {s2} {s3}"
);
}
#[test]
fn test_safety_factor_von_mises() {
let sigma = [200e6, 0.0, 0.0, 0.0, 0.0, 0.0];
let s = StressState::from_voigt(sigma);
let sf = s.safety_factor_von_mises(400e6);
assert!((sf - 2.0).abs() < 1e-6, "safety factor={sf}");
}
#[test]
fn test_lode_angle_uniaxial_tension() {
let sigma = [1.0, 0.0, 0.0, 0.0, 0.0, 0.0];
let s = StressState::from_voigt(sigma);
let theta = s.lode_angle_deg();
assert!(theta.is_finite(), "Lode angle should be finite: {theta}");
}
#[test]
fn test_tresca_criterion_uniaxial() {
let sigma = [200e6, 0.0, 0.0, 0.0, 0.0, 0.0];
let s = StressState::from_voigt(sigma);
let tresca = s.tresca_stress();
assert!((tresca - 200e6).abs() < 1.0, "Tresca={tresca}");
}
#[test]
fn test_fatigue_life_below_endurance() {
let est = FatigueLifeEstimator::new(500e6, 1e6, 200e6);
let life = est.cycles_to_failure(100e6);
assert!(
life == f64::INFINITY,
"below endurance should return infinity: {life}"
);
}
#[test]
fn test_fatigue_life_at_ultimate() {
let est = FatigueLifeEstimator::new(500e6, 1e6, 200e6);
let life = est.cycles_to_failure(500e6);
assert!(life <= 1e6 + 1.0, "at UTS life should be <= N_ref: {life}");
}
#[test]
fn test_damage_accumulation_miner() {
let d = miner_damage(&[0.3, 0.5, 0.2]);
assert!((d - 1.0).abs() < 1e-12, "Miner sum = {d}");
}
#[test]
fn test_cantilever_beam() {
let length: f64 = 1.0;
let width: f64 = 0.1;
let height: f64 = 0.1;
let e: f64 = 200.0e9;
let nu: f64 = 0.3;
let p: f64 = -1000.0;
let inertia = width * height.powi(3) / 12.0;
let analytical_deflection = p * length.powi(3) / (3.0 * e * inertia);
let nx = 10;
let ny = 2;
let nz = 2;
let mesh = TetrahedralMesh::generate_beam(length, width, height, nx, ny, nz);
let material = LinearElasticMaterial::new(e, nu);
let mut dirichlet_bcs = Vec::new();
for (i, node) in mesh.nodes.iter().enumerate() {
if node.x.abs() < 1e-10 {
dirichlet_bcs.push(DirichletBc::new(i, 0, 0.0));
dirichlet_bcs.push(DirichletBc::new(i, 1, 0.0));
dirichlet_bcs.push(DirichletBc::new(i, 2, 0.0));
}
}
let free_end_nodes: Vec<usize> = mesh
.nodes
.iter()
.enumerate()
.filter(|(_, n)| (n.x - length).abs() < 1e-10)
.map(|(i, _)| i)
.collect();
let num_free = free_end_nodes.len() as f64;
let neumann_bcs: Vec<NeumannBc> = free_end_nodes
.iter()
.map(|&i| NeumannBc::new(i, [0.0, p / num_free, 0.0]))
.collect();
let analysis = LinearStaticAnalysis {
max_iter: 50000,
tolerance: 1e-12,
};
let result = analysis.solve(
&mesh,
&material,
&dirichlet_bcs,
&neumann_bcs,
&Vec3::new(0.0, 0.0, 0.0),
);
let max_y_disp = free_end_nodes
.iter()
.map(|&i| result.displacements[i].y)
.fold(f64::INFINITY, f64::min);
let error = ((max_y_disp - analytical_deflection) / analytical_deflection).abs();
assert!(
max_y_disp < 0.0,
"Deflection should be negative (downward), got {max_y_disp:.6e}",
);
assert!(
error < 1.0,
"Cantilever deflection error too large: {:.1}%. \
FEM = {max_y_disp:.6e}, analytical = {analytical_deflection:.6e}",
error * 100.0,
);
}
}