gemlab 2.0.0

Geometry and meshes laboratory for finite element analyses
Documentation
//! Functions to perform numerical integration using Shapes
//!
//! # Definitions
//!
//! * `ngauss` -- number of integration (Gauss) points
//! * `|J|` -- the determinant of the Jacobian dx/dξ
//! * `||J||` -- the norm of the Jacobian vector for lines in multi-dimensions
//! * `ξ` -- ksi (or xi) -- coordinates in the reference space
//! * `ιᵖ := ξᵖ` -- (**iota**-p gets ksi-p (or xi-p)) --
//!    the coordinate of the integration point on the reference (natural) space
//! * `wᵖ` -- the weight of the p-th integration point
//! * Xᵀ -- Transposed matrix of coordinates
//! * N -- Shape functions
//! * L -- Derivatives of shape functions
//! * J -- Jacobian tensor
//! * J⁻¹ -- Inverse Jacobian matrix (only available if geo_ndim = space_ndim)
//! * B -- Gradients of shape functions (derivatives w.r.t real coordinates x)
//! * see also the subsection named **Expressions** below
//!
//! # Integration of scalar field over a geometric shape
//!
//! Function [scalar_field()]
//!
//! ```text
//!     ⌠   → →
//! I = │ s(x(ξ)) dΩ
//!//!     Ωₑ
//! ```
//!
//! # Vector results: Integration of some combinations involving N and B resulting in vectors
//!
//! ## VEC 01: Shape(N) times scalar(S)
//!
//! Function [vec_01_ns()]
//!
//! ```text
//!      ⌠    → →     →
//! aᵐ = │ Nᵐ(x(ξ)) s(x) dΩ
//!//!      Ωₑ
//! ```
//!
//! ## VEC 01(bry): Shape(N) times scalar(S) (boundary integral version)
//!
//! Function [vec_01_ns_bry()]
//!
//! ```text
//!      ⌠    → →     →
//! aᵐ = │ Nᵐ(x(ξ)) s(x) dΩ
//!//!      Ωₑ
//! ```
//!
//! ## VEC 02: Shape(N) times vector(V)
//!
//! Function [vec_02_nv()]
//!
//! ```text
//! →    ⌠    → →   → →
//! bᵐ = │ Nᵐ(x(ξ)) v(x) dΩ
//!//!      Ωₑ
//! ```
//!
//! ## VEC 02(bry): Shape(N) times vector(V) (boundary integral version)
//!
//! Function [vec_02_nv_bry()]
//!
//! ```text
//! →    ⌠    → →   → →
//! bᵐ = │ Nᵐ(x(ξ)) v(x) dΓ
//!//!      Γₑ
//! ```
//!
//! ## VEC 03: Vector(V) dot gradient(B)
//!
//! Function [vec_03_bv()]
//!
//! ```text
//!      ⌠ →  → →     → →
//! cᵐ = │ Bᵐ(x(ξ)) · w(x) dΩ
//!//!      Ωₑ
//! ```
//!
//! ## VEC 04: gradient(B) dot transpose tensor(T)
//!
//! Function [vec_04_bt()]
//!
//! ```text
//! →    ⌠ →  → →        →      ⌠   →    →  → →
//! dᵐ = │ Bᵐ(x(ξ)) · σᵀ(x) dΩ  │ σ(x) · Bᵐ(x(ξ)) dΩ
//!      ⌡            ▔         ⌡ ▔
//!      Ωₑ                     Ωₑ
//! ```
//!
//! # Matrix results: Integration of some combinations involving N, tensors, and B, resulting in matrices
//!
//! ![cases](https://github.com/cpmech/gemlab/raw/main/data/figures/integ-matrix-cases.png)
//!
//! ## MAT 01: Shape(N) times scalar(S) times shape(N) (e.g., diffusion matrix)
//!
//! Function [mat_01_nsn()]
//!
//! ```text
//!//! Kᵐⁿ = │ Nᵐ s Nⁿ dΩ
//!//!       Ωₑ
//! ```
//!
//! ## MAT 02: Gradient(B) dot vector(V) times shape(N) (e.g., compressibility matrix)
//!
//! Function [mat_02_bvn()]
//!
//! ```text
//!       ⌠ →    →
//! Kᵐⁿ = │ Bᵐ ⋅ v Nⁿ dΩ
//!//!       Ωₑ
//! ```
//!
//! ## MAT 03: Gradient(B) dot tensor(T) dot gradient(B) (e.g., conductivity matrix)
//!
//! Function [mat_03_btb()]
//!
//! ```text
//!       ⌠ →        →
//! Kᵐⁿ = │ Bᵐ ⋅ T ⋅ Bⁿ dΩ
//!       ⌡      ▔
//!       Ωₑ
//! ```
//!
//! ## MAT 04: shape(Nb) times scalar(S) times gradient(B) (e.g., coupling matrix)
//!
//! Function [mat_04_nsb()]
//!
//! ```text
//! →     ⌠       →
//! Kᵐⁿ = │ Nbᵐ s Bⁿ dΩ
//!//!       Ωₑ
//! ```
//!
//! ## MAT 05: Gradient(Bb) times tensor(T) times shape(N) (e.g., coupling matrix)
//!
//! Function [mat_05_btn()]
//!
//! ```text
//! →     ⌠ →
//! Kᵐⁿ = │ Bbᵐ ⋅ T Nⁿ dΩ
//!       ⌡       ▔
//!       Ωₑ
//! ```
//!
//! ## MAT 06: Shape(N) times vector(V) times shape(Nb) (e.g., coupling matrix)
//!
//! Function [mat_06_nvn()]
//!
//! ```text
//! →     ⌠    →
//! Kᵐⁿ = │ Nᵐ v Nbⁿ dΩ
//!//!       Ωₑ
//! ```
//!
//! ## MAT 07: Gradient(B) times scalar(S) times shape(Nb) (e.g., coupling matrix)
//!
//! Function [mat_07_bsn()]
//!
//! ```text
//! →     ⌠ →
//! Kᵐⁿ = │ Bᵐ s Nbⁿ dΩ
//!//!       Ωₑ
//! ```
//!
//! ## MAT 08: Shape(N) times tensor(T) times shape(N) (e.g., mass matrix)
//!
//! Function [mat_08_ntn()]
//!
//! ```text
//!//! Kᵐⁿ = │ Nᵐ T Nⁿ dΩ
//! ▔     ⌡    ▔
//!       Ωₑ
//! ```
//!
//! ## MAT 09: Shape(N) times vector(V) dot gradient(B) (e.g., variable density matrix)
//!
//! Function [mat_09_nvb()]
//!
//! ```text
//!       ⌠    →   →
//! Kᵐⁿ = │ Nᵐ v ⊗ Bⁿ dΩ
//! ▔     ⌡
//!       Ωₑ
//! ```
//!
//! ## MAT 10: Gradient(B) dot 4th-tensor(D) dot gradient(B) (e.g., stiffness matrix)
//!
//! Function [mat_10_bdb()]
//!
//! ```text
//!       ⌠                       →    →
//! Kᵐⁿ = │ Σ Σ Σ Σ Bᵐₖ Dᵢₖⱼₗ Bⁿₗ eᵢ ⊗ eⱼ dΩ
//! ▔     ⌡ i j k l
//!       Ωₑ
//! ```
//!
//! # Expressions
//!
//! See also [crate::shapes]
//!
//! ## Xᵀ: Transposed matrix of coordinates
//!
//! ```text
//!      ┌                              ┐  superscript = node
//!      | x⁰₀  x¹₀  x²₀  x³₀       xᴹ₀ |  subscript = space dimension
//! Xᵀ = | x⁰₁  x¹₁  x²₁  x³₁  ...  xᴹ₁ |
//!      | x⁰₂  x¹₂  x²₂  x³₂       xᴹ₂ |
//!      └                              ┘_(space_ndim,nnode)
//! where `M = nnode - 1`
//! ```
//!
//! ## N: Shape functions
//!
//! Shape function of node m at ξ (ksi; i.e., xi):
//!
//! ```text
//! Nᵐ(ξ)
//! matrix notation: N is an (nnode) vector
//! ```
//!
//! ## L: Derivatives of shape functions
//!
//! Derivatives of shape functions with respect to the natural coordinates:
//!
//! ```text
//!//! →  →    dNᵐ(ξ)
//! Lᵐ(ξ) = ——————
//!//!//! matrix notation: L is an (nnode,geo_ndim) matrix
//! ```
//!
//! ## Isoparametric formula
//!
//! ```text
//! → →         →  →
//! x(ξ) = Σ Nᵐ(ξ) xᵐ
//!        m
//! matrix notation: x = Xᵀ ⋅ N
//! ```
//!
//! ## J: Jacobian tensor (SOLID case with geo_ndim = space_ndim = 2 or 3)
//!
//! ```text
//!//!   →    dx     →    →
//! J(ξ) = —— = Σ xᵐ ⊗ Lᵐ
//!         →   m
//!//! matrix notation: J = Xᵀ · L
//! J is a (space_ndim,geo_ndim) matrix
//! ```
//!
//! ## J: Jacobian vector (CABLE case with geo_ndim = 1 and space_ndim = 2 or 3)
//!
//! ```text
//!//! →    →      →    →  →    dx
//! J := Jcable(ξ) = g₁(ξ) = ——
//!//! matrix notation: J = Jcable = Xᵀ · L
//! J is a (space_ndim,1) matrix; i.e., a vector
//! ```
//!
//! ## J: Jacobian matrix (SHELL case with geo_ndim = 2 and space_ndim = 3)
//!
//! ```text
//!                 dx
//! J(ξ) = Jshell = ——
//!//! matrix notation: J = Jshell = Xᵀ · L
//! J is a (3,2) matrix
//! ```
//!
//! ## B: Gradients of shape functions (derivatives w.r.t real coordinates x) (only for SOLID case)
//!
//! ```text
//!//! →  →    dNᵐ(ξ)
//! Bᵐ(ξ) = ——————
//!//!           dx
//! matrix notation: B = L · J⁻¹
//! B is an (nnode,space_ndim) matrix
//! ```

mod analytical_qua4;
mod analytical_qua8;
mod analytical_tet4;
mod analytical_tri3;
mod common_args;
mod gauss;
mod mat_01_nsn;
mod mat_01_nsn_bry;
mod mat_02_bvn;
mod mat_03_btb;
mod mat_04_nsb;
mod mat_05_btn;
mod mat_06_nvn;
mod mat_07_bsn;
mod mat_08_ntn;
mod mat_09_nvb;
mod mat_10_bdb;
mod scalar_field;
mod testing;
mod vec_01_ns;
mod vec_01_ns_bry;
mod vec_02_nv;
mod vec_02_nv_bry;
mod vec_03_bv;
mod vec_04_bt;

pub use analytical_qua4::*;
pub use analytical_qua8::*;
pub use analytical_tet4::*;
pub use analytical_tri3::*;
pub use common_args::*;
pub use gauss::*;
pub use mat_01_nsn::*;
pub use mat_01_nsn_bry::*;
pub use mat_02_bvn::*;
pub use mat_03_btb::*;
pub use mat_04_nsb::*;
pub use mat_05_btn::*;
pub use mat_06_nvn::*;
pub use mat_07_bsn::*;
pub use mat_08_ntn::*;
pub use mat_09_nvb::*;
pub use mat_10_bdb::*;
pub use scalar_field::*;
pub use vec_01_ns::*;
pub use vec_01_ns_bry::*;
pub use vec_02_nv::*;
pub use vec_02_nv_bry::*;
pub use vec_03_bv::*;
pub use vec_04_bt::*;