Struct libceed::basis::Basis

source ·
pub struct Basis<'a> { /* private fields */ }

Implementations§

source§

impl<'a> Basis<'a>

source

pub fn create_tensor_H1( ceed: &Ceed, dim: usize, ncomp: usize, P1d: usize, Q1d: usize, interp1d: &[Scalar], grad1d: &[Scalar], qref1d: &[Scalar], qweight1d: &[Scalar] ) -> Result<Self>

source

pub fn create_tensor_H1_Lagrange( ceed: &Ceed, dim: usize, ncomp: usize, P: usize, Q: usize, qmode: QuadMode ) -> Result<Self>

source

pub fn create_H1( ceed: &Ceed, topo: ElemTopology, ncomp: usize, nnodes: usize, nqpts: usize, interp: &[Scalar], grad: &[Scalar], qref: &[Scalar], qweight: &[Scalar] ) -> Result<Self>

source

pub fn create_Hdiv( ceed: &Ceed, topo: ElemTopology, ncomp: usize, nnodes: usize, nqpts: usize, interp: &[Scalar], div: &[Scalar], qref: &[Scalar], qweight: &[Scalar] ) -> Result<Self>

source

pub fn create_Hcurl( ceed: &Ceed, topo: ElemTopology, ncomp: usize, nnodes: usize, nqpts: usize, interp: &[Scalar], curl: &[Scalar], qref: &[Scalar], qweight: &[Scalar] ) -> Result<Self>

source

pub fn apply( &self, nelem: usize, tmode: TransposeMode, emode: EvalMode, u: &Vector<'_>, v: &mut Vector<'_> ) -> Result<i32>

Apply basis evaluation from nodes to quadrature points or vice versa

  • nelem - The number of elements to apply the basis evaluation to
  • tmode - TrasposeMode::NoTranspose to evaluate from nodes to quadrature points, TransposeMode::Transpose to apply the transpose, mapping from quadrature points to nodes
  • emode - EvalMode::None to use values directly, EvalMode::Interp to use interpolated values, EvalMode::Grad to use gradients, EvalMode::Weight to use quadrature weights
  • u - Input Vector
  • v - Output Vector
const Q: usize = 6;
let bu = ceed.basis_tensor_H1_Lagrange(1, 1, Q, Q, QuadMode::GaussLobatto)?;
let bx = ceed.basis_tensor_H1_Lagrange(1, 1, 2, Q, QuadMode::Gauss)?;

let x_corners = ceed.vector_from_slice(&[-1., 1.])?;
let mut x_qpts = ceed.vector(Q)?;
let mut x_nodes = ceed.vector(Q)?;
bx.apply(
    1,
    TransposeMode::NoTranspose,
    EvalMode::Interp,
    &x_corners,
    &mut x_nodes,
)?;
bu.apply(
    1,
    TransposeMode::NoTranspose,
    EvalMode::Interp,
    &x_nodes,
    &mut x_qpts,
)?;

// Create function x^3 + 1 on Gauss Lobatto points
let mut u_arr = [0.; Q];
u_arr
    .iter_mut()
    .zip(x_nodes.view()?.iter())
    .for_each(|(u, x)| *u = x * x * x + 1.);
let u = ceed.vector_from_slice(&u_arr)?;

// Map function to Gauss points
let mut v = ceed.vector(Q)?;
v.set_value(0.);
bu.apply(1, TransposeMode::NoTranspose, EvalMode::Interp, &u, &mut v)?;

// Verify results
v.view()?
    .iter()
    .zip(x_qpts.view()?.iter())
    .for_each(|(v, x)| {
        let true_value = x * x * x + 1.;
        assert!(
            (*v - true_value).abs() < 10.0 * libceed::EPSILON,
            "Incorrect basis application"
        );
    });
source

pub fn dimension(&self) -> usize

Returns the dimension for given Basis

let dim = 2;
let b = ceed.basis_tensor_H1_Lagrange(dim, 1, 3, 4, QuadMode::Gauss)?;

let d = b.dimension();
assert_eq!(d, dim, "Incorrect dimension");
source

pub fn num_components(&self) -> usize

Returns number of components for given Basis

let ncomp = 2;
let b = ceed.basis_tensor_H1_Lagrange(1, ncomp, 3, 4, QuadMode::Gauss)?;

let n = b.num_components();
assert_eq!(n, ncomp, "Incorrect number of components");
source

pub fn num_nodes(&self) -> usize

Returns total number of nodes (in dim dimensions) of a Basis

let p = 3;
let b = ceed.basis_tensor_H1_Lagrange(2, 1, p, 4, QuadMode::Gauss)?;

let nnodes = b.num_nodes();
assert_eq!(nnodes, p * p, "Incorrect number of nodes");
source

pub fn num_quadrature_points(&self) -> usize

Returns total number of quadrature points (in dim dimensions) of a Basis

let q = 4;
let b = ceed.basis_tensor_H1_Lagrange(2, 1, 3, q, QuadMode::Gauss)?;

let nqpts = b.num_quadrature_points();
assert_eq!(nqpts, q * q, "Incorrect number of quadrature points");
source

pub fn create_projection(&self, to: &Self) -> Result<Self>

Create projection from self to specified Basis.

Both bases must have the same quadrature space. The input bases need not be nested as function spaces; this interface solves a least squares problem to find a representation in the to basis that agrees at quadrature points with the origin basis. Since the bases need not be Lagrange, the resulting projection “basis” will have empty quadrature points and weights.

let coarse = ceed.basis_tensor_H1_Lagrange(1, 1, 2, 3, QuadMode::Gauss)?;
let fine = ceed.basis_tensor_H1_Lagrange(1, 1, 3, 3, QuadMode::Gauss)?;
let proj = coarse.create_projection(&fine)?;
let u = ceed.vector_from_slice(&[1., 2.])?;
let mut v = ceed.vector(3)?;
proj.apply(1, TransposeMode::NoTranspose, EvalMode::Interp, &u, &mut v)?;
let expected = [1., 1.5, 2.];
for (a, b) in v.view()?.iter().zip(expected) {
    assert!(
        (a - b).abs() < 10.0 * libceed::EPSILON,
        "Incorrect projection of linear Lagrange to quadratic Lagrange"
    );
}

Trait Implementations§

source§

impl<'a> Debug for Basis<'a>

source§

fn fmt(&self, f: &mut Formatter<'_>) -> Result

Formats the value using the given formatter. Read more
source§

impl<'a> Display for Basis<'a>

source§

fn fmt(&self, f: &mut Formatter<'_>) -> Result

View a Basis

let b = ceed.basis_tensor_H1_Lagrange(1, 2, 3, 4, QuadMode::Gauss)?;
println!("{}", b);
source§

impl<'a> Drop for Basis<'a>

source§

fn drop(&mut self)

Executes the destructor for this type. Read more
source§

impl<'a> From<&'a Basis<'_>> for BasisOpt<'a>

Construct a BasisOpt reference from a Basis reference

source§

fn from(basis: &'a Basis<'_>) -> Self

Converts to this type from the input type.

Auto Trait Implementations§

§

impl<'a> RefUnwindSafe for Basis<'a>

§

impl<'a> !Send for Basis<'a>

§

impl<'a> !Sync for Basis<'a>

§

impl<'a> Unpin for Basis<'a>

§

impl<'a> UnwindSafe for Basis<'a>

Blanket Implementations§

source§

impl<T> Any for T
where T: 'static + ?Sized,

source§

fn type_id(&self) -> TypeId

Gets the TypeId of self. Read more
source§

impl<T> Borrow<T> for T
where T: ?Sized,

source§

fn borrow(&self) -> &T

Immutably borrows from an owned value. Read more
source§

impl<T> BorrowMut<T> for T
where T: ?Sized,

source§

fn borrow_mut(&mut self) -> &mut T

Mutably borrows from an owned value. Read more
source§

impl<T> From<T> for T

source§

fn from(t: T) -> T

Returns the argument unchanged.

source§

impl<T, U> Into<U> for T
where U: From<T>,

source§

fn into(self) -> U

Calls U::from(self).

That is, this conversion is whatever the implementation of From<T> for U chooses to do.

source§

impl<T> ToString for T
where T: Display + ?Sized,

source§

default fn to_string(&self) -> String

Converts the given value to a String. Read more
source§

impl<T, U> TryFrom<U> for T
where U: Into<T>,

§

type Error = Infallible

The type returned in the event of a conversion error.
source§

fn try_from(value: U) -> Result<T, <T as TryFrom<U>>::Error>

Performs the conversion.
source§

impl<T, U> TryInto<U> for T
where U: TryFrom<T>,

§

type Error = <U as TryFrom<T>>::Error

The type returned in the event of a conversion error.
source§

fn try_into(self) -> Result<U, <U as TryFrom<T>>::Error>

Performs the conversion.