Struct libceed::basis::Basis[][src]

pub struct Basis<'a> { /* fields omitted */ }

Implementations

impl<'a> Basis<'a>[src]

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

pub fn create_tensor_H1_Lagrange(
    ceed: &'a Ceed,
    dim: usize,
    ncomp: usize,
    P: usize,
    Q: usize,
    qmode: QuadMode
) -> Result<Self, CeedError>
[src]

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

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

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)
    .unwrap();
let bx = ceed
    .basis_tensor_H1_Lagrange(1, 1, 2, Q, QuadMode::Gauss)
    .unwrap();

let x_corners = ceed.vector_from_slice(&[-1., 1.]).unwrap();
let mut x_qpts = ceed.vector(Q).unwrap();
let mut x_nodes = ceed.vector(Q).unwrap();
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).unwrap();

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

// Verify results
v.view()
    .iter()
    .zip(x_qpts.view().iter())
    .for_each(|(v, x)| {
        let true_value = x * x * x + 1.;
        assert_eq!(*v, true_value, "Incorrect basis application");
    });

pub fn dimension(&self) -> usize[src]

Returns the dimension for given CeedBasis

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

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

pub fn num_components(&self) -> usize[src]

Returns number of components for given CeedBasis

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

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

pub fn num_nodes(&self) -> usize[src]

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

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

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

pub fn num_quadrature_points(&self) -> usize[src]

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

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

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

Trait Implementations

impl<'a> Display for Basis<'a>[src]

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

View a Basis

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

impl<'a> Drop for Basis<'a>[src]

impl<'a> From<&'a Basis<'_>> for BasisOpt<'a>[src]

Construct a BasisOpt reference from a Basis reference

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

impl<T> Any for T where
    T: 'static + ?Sized
[src]

impl<T> Borrow<T> for T where
    T: ?Sized
[src]

impl<T> BorrowMut<T> for T where
    T: ?Sized
[src]

impl<T> From<T> for T[src]

impl<T, U> Into<U> for T where
    U: From<T>, 
[src]

impl<T> ToString for T where
    T: Display + ?Sized
[src]

impl<T, U> TryFrom<U> for T where
    U: Into<T>, 
[src]

type Error = Infallible

The type returned in the event of a conversion error.

impl<T, U> TryInto<U> for T where
    U: TryFrom<T>, 
[src]

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

The type returned in the event of a conversion error.