use super::lagrange::*;
use crate::mesh::ElementType;
#[derive(Debug, Clone)]
pub struct ShapeValues {
pub values: Vec<f64>,
pub gradients: Vec<Vec<f64>>,
}
pub fn evaluate_shape(
element_type: ElementType,
degree: PolynomialDegree,
xi: f64,
eta: f64,
zeta: f64,
) -> ShapeValues {
match (element_type, degree) {
(ElementType::Triangle, PolynomialDegree::P1) => {
let vals = p1_triangle(xi, eta);
let grads = p1_triangle_grad();
ShapeValues {
values: vals.to_vec(),
gradients: grads.iter().map(|g| g.to_vec()).collect(),
}
}
(ElementType::Triangle, PolynomialDegree::P2) => {
let vals = p2_triangle(xi, eta);
let grads = p2_triangle_grad(xi, eta);
ShapeValues {
values: vals.to_vec(),
gradients: grads.iter().map(|g| g.to_vec()).collect(),
}
}
(ElementType::Quadrilateral, PolynomialDegree::P1) => {
let vals = q1_quadrilateral(xi, eta);
let grads = q1_quadrilateral_grad(xi, eta);
ShapeValues {
values: vals.to_vec(),
gradients: grads.iter().map(|g| g.to_vec()).collect(),
}
}
(ElementType::Tetrahedron, PolynomialDegree::P1) => {
let vals = p1_tetrahedron(xi, eta, zeta);
let grads = p1_tetrahedron_grad();
ShapeValues {
values: vals.to_vec(),
gradients: grads.iter().map(|g| g.to_vec()).collect(),
}
}
(ElementType::Hexahedron, PolynomialDegree::P1) => {
let vals = q1_hexahedron(xi, eta, zeta);
let grads = q1_hexahedron_grad(xi, eta, zeta);
ShapeValues {
values: vals.to_vec(),
gradients: grads.iter().map(|g| g.to_vec()).collect(),
}
}
_ => {
log::warn!(
"Shape functions for {:?} degree {:?} not implemented, using P1",
element_type,
degree
);
evaluate_shape(element_type, PolynomialDegree::P1, xi, eta, zeta)
}
}
}
#[derive(Debug, Clone)]
pub struct Jacobian {
pub matrix: Vec<Vec<f64>>,
pub det: f64,
pub inverse: Vec<Vec<f64>>,
}
impl Jacobian {
pub fn from_2d(grad_ref: &[Vec<f64>], coords: &[[f64; 2]]) -> Self {
let mut j = [[0.0; 2]; 2];
for (i, g) in grad_ref.iter().enumerate() {
j[0][0] += g[0] * coords[i][0]; j[0][1] += g[1] * coords[i][0]; j[1][0] += g[0] * coords[i][1]; j[1][1] += g[1] * coords[i][1]; }
let det = j[0][0] * j[1][1] - j[0][1] * j[1][0];
let inv_det = 1.0 / det;
let inverse = vec![
vec![j[1][1] * inv_det, -j[0][1] * inv_det],
vec![-j[1][0] * inv_det, j[0][0] * inv_det],
];
Self {
matrix: vec![vec![j[0][0], j[0][1]], vec![j[1][0], j[1][1]]],
det,
inverse,
}
}
pub fn from_3d(grad_ref: &[Vec<f64>], coords: &[[f64; 3]]) -> Self {
let mut j = [[0.0; 3]; 3];
for (i, g) in grad_ref.iter().enumerate() {
for k in 0..3 {
j[0][k] += g[k] * coords[i][0]; j[1][k] += g[k] * coords[i][1]; j[2][k] += g[k] * coords[i][2]; }
}
let det = j[0][0] * (j[1][1] * j[2][2] - j[1][2] * j[2][1])
- j[0][1] * (j[1][0] * j[2][2] - j[1][2] * j[2][0])
+ j[0][2] * (j[1][0] * j[2][1] - j[1][1] * j[2][0]);
let inv_det = 1.0 / det;
let inverse = vec![
vec![
(j[1][1] * j[2][2] - j[1][2] * j[2][1]) * inv_det,
(j[0][2] * j[2][1] - j[0][1] * j[2][2]) * inv_det,
(j[0][1] * j[1][2] - j[0][2] * j[1][1]) * inv_det,
],
vec![
(j[1][2] * j[2][0] - j[1][0] * j[2][2]) * inv_det,
(j[0][0] * j[2][2] - j[0][2] * j[2][0]) * inv_det,
(j[0][2] * j[1][0] - j[0][0] * j[1][2]) * inv_det,
],
vec![
(j[1][0] * j[2][1] - j[1][1] * j[2][0]) * inv_det,
(j[0][1] * j[2][0] - j[0][0] * j[2][1]) * inv_det,
(j[0][0] * j[1][1] - j[0][1] * j[1][0]) * inv_det,
],
];
Self {
matrix: vec![
vec![j[0][0], j[0][1], j[0][2]],
vec![j[1][0], j[1][1], j[1][2]],
vec![j[2][0], j[2][1], j[2][2]],
],
det,
inverse,
}
}
pub fn transform_gradient(&self, grad_ref: &[f64]) -> Vec<f64> {
let dim = self.inverse.len();
let mut result = vec![0.0; dim];
for i in 0..dim {
for j in 0..dim {
result[i] += self.inverse[j][i] * grad_ref[j];
}
}
result
}
}
#[cfg(test)]
mod tests {
use super::*;
#[test]
fn test_evaluate_shape_p1_triangle() {
let shape = evaluate_shape(ElementType::Triangle, PolynomialDegree::P1, 0.25, 0.25, 0.0);
assert_eq!(shape.values.len(), 3);
assert_eq!(shape.gradients.len(), 3);
let sum: f64 = shape.values.iter().sum();
assert!((sum - 1.0).abs() < 1e-14);
}
#[test]
fn test_jacobian_2d_unit_triangle() {
let coords = [[0.0, 0.0], [1.0, 0.0], [0.0, 1.0]];
let grads = vec![vec![-1.0, -1.0], vec![1.0, 0.0], vec![0.0, 1.0]];
let jac = Jacobian::from_2d(&grads, &coords);
assert!((jac.matrix[0][0] - 1.0).abs() < 1e-14);
assert!(jac.matrix[0][1].abs() < 1e-14);
assert!(jac.matrix[1][0].abs() < 1e-14);
assert!((jac.matrix[1][1] - 1.0).abs() < 1e-14);
assert!((jac.det - 1.0).abs() < 1e-14);
}
#[test]
fn test_jacobian_2d_scaled_triangle() {
let coords = [[0.0, 0.0], [2.0, 0.0], [0.0, 2.0]];
let grads = vec![vec![-1.0, -1.0], vec![1.0, 0.0], vec![0.0, 1.0]];
let jac = Jacobian::from_2d(&grads, &coords);
assert!((jac.matrix[0][0] - 2.0).abs() < 1e-14);
assert!((jac.matrix[1][1] - 2.0).abs() < 1e-14);
assert!((jac.det - 4.0).abs() < 1e-14);
}
#[test]
fn test_gradient_transformation() {
let coords = [[0.0, 0.0], [1.0, 0.0], [0.0, 1.0]];
let grads = vec![vec![-1.0, -1.0], vec![1.0, 0.0], vec![0.0, 1.0]];
let jac = Jacobian::from_2d(&grads, &coords);
let ref_grad = vec![1.0, 2.0];
let phys_grad = jac.transform_gradient(&ref_grad);
assert!((phys_grad[0] - 1.0).abs() < 1e-14);
assert!((phys_grad[1] - 2.0).abs() < 1e-14);
}
}