Struct libceed::Ceed[][src]

pub struct Ceed { /* fields omitted */ }

A Ceed is a library context representing control of a logical hardware resource.

Implementations

impl Ceed[src]

pub fn init(resource: &str) -> Self[src]

Returns a Ceed context initialized with the specified resource

arguments

  • resource - Resource to use, e.g., “/cpu/self”
let ceed = libceed::Ceed::init("/cpu/self/ref/serial");

pub fn init_with_error_handler(
    resource: &str,
    handler: CeedErrorHandler
) -> Self
[src]

Returns a Ceed context initialized with the specified resource

arguments

  • resource - Resource to use, e.g., “/cpu/self”
let ceed = libceed::Ceed::init_with_error_handler(
    "/cpu/self/ref/serial",
    libceed::CeedErrorHandler::ErrorAbort,
);

pub fn resource(&self) -> String[src]

Returns full resource name for a Ceed context

let ceed = libceed::Ceed::init("/cpu/self/ref/serial");
let resource = ceed.resource();

assert_eq!(resource, "/cpu/self/ref/serial".to_string())

pub fn vector(&self, n: usize) -> Result<Vector<'_>, CeedError>[src]

Returns a CeedVector of the specified length (does not allocate memory)

arguments

  • n - Length of vector
let vec = ceed.vector(10).unwrap();

pub fn vector_from_slice(&self, slice: &[f64]) -> Result<Vector<'_>, CeedError>[src]

Create a Vector initialized with the data (copied) from a slice

arguments

  • slice - Slice containing data
let vec = ceed.vector_from_slice(&[1., 2., 3.]).unwrap();
assert_eq!(vec.length(), 3);

pub fn elem_restriction(
    &self,
    nelem: usize,
    elemsize: usize,
    ncomp: usize,
    compstride: usize,
    lsize: usize,
    mtype: MemType,
    offsets: &[i32]
) -> Result<ElemRestriction<'_>, CeedError>
[src]

Returns a ElemRestriction

arguments

  • nelem - Number of elements described in the offsets array
  • elemsize - Size (number of “nodes”) per element
  • ncomp - Number of field components per interpolation node (1 for scalar fields)
  • compstride - Stride between components for the same Lvector “node”. Data for node i, component j, element k can be found in the Lvector at index offsets[i + k*elemsize] + j*compstride.
  • lsize - The size of the Lvector. This vector may be larger than the elements and fields given by this restriction.
  • mtype - Memory type of the offsets array, see CeedMemType
  • offsets - Array of shape [nelem, elemsize]. Row i holds the ordered list of the offsets (into the input CeedVector) for the unknowns corresponding to element i, where 0 <= i < nelem. All offsets must be in the range [0, lsize - 1].
let nelem = 3;
let mut ind: Vec<i32> = vec![0; 2 * nelem];
for i in 0..nelem {
    ind[2 * i + 0] = i as i32;
    ind[2 * i + 1] = (i + 1) as i32;
}
let r = ceed
    .elem_restriction(nelem, 2, 1, 1, nelem + 1, MemType::Host, &ind)
    .unwrap();

pub fn strided_elem_restriction(
    &self,
    nelem: usize,
    elemsize: usize,
    ncomp: usize,
    lsize: usize,
    strides: [i32; 3]
) -> Result<ElemRestriction<'_>, CeedError>
[src]

Returns a ElemRestriction

arguments

  • nelem - Number of elements described in the offsets array
  • elemsize - Size (number of “nodes”) per element
  • ncomp - Number of field components per interpolation node (1 for scalar fields)
  • compstride - Stride between components for the same Lvector “node”. Data for node i, component j, element k can be found in the Lvector at index offsets[i + k*elemsize] + j*compstride.
  • lsize - The size of the Lvector. This vector may be larger than the elements and fields given by this restriction.
  • strides - Array for strides between [nodes, components, elements]. Data for node i, component j, element k can be found in the Lvector at index i*strides[0] + j*strides[1] + k*strides[2]. CEED_STRIDES_BACKEND may be used with vectors created by a Ceed backend.
let nelem = 3;
let strides: [i32; 3] = [1, 2, 2];
let r = ceed
    .strided_elem_restriction(nelem, 2, 1, nelem * 2, strides)
    .unwrap();

pub fn basis_tensor_H1(
    &self,
    dim: usize,
    ncomp: usize,
    P1d: usize,
    Q1d: usize,
    interp1d: &[f64],
    grad1d: &[f64],
    qref1d: &[f64],
    qweight1d: &[f64]
) -> Result<Basis<'_>, CeedError>
[src]

Returns a tensor-product basis

arguments

  • dim - Topological dimension of element
  • ncomp - Number of field components (1 for scalar fields)
  • P1d - Number of Gauss-Lobatto nodes in one dimension. The polynomial degree of the resulting Q_k element is k=P-1.
  • Q1d - Number of quadrature points in one dimension
  • interp1d - Row-major (Q1d * P1d) matrix expressing the values of nodal basis functions at quadrature points
  • grad1d - Row-major (Q1d * P1d) matrix expressing derivatives of nodal basis functions at quadrature points
  • qref1d - Array of length Q1d holding the locations of quadrature points on the 1D reference element [-1, 1]
  • qweight1d - Array of length Q1d holding the quadrature weights on the reference element
let interp1d  = [ 0.62994317,  0.47255875, -0.14950343,  0.04700152,
                 -0.07069480,  0.97297619,  0.13253993, -0.03482132,
                 -0.03482132,  0.13253993,  0.97297619, -0.07069480,
                  0.04700152, -0.14950343,  0.47255875,  0.62994317];
let grad1d    = [-2.34183742,  2.78794489, -0.63510411,  0.18899664,
                 -0.51670214, -0.48795249,  1.33790510, -0.33325047,
                 -0.18899664,  0.63510411, -2.78794489,  2.34183742];
let qref1d    = [-0.86113631, -0.33998104,  0.33998104,  0.86113631];
let qweight1d = [ 0.34785485,  0.65214515,  0.65214515,  0.34785485];
let b = ceed.
basis_tensor_H1(2, 1, 4, 4, &interp1d, &grad1d, &qref1d, &qweight1d).unwrap();

pub fn basis_tensor_H1_Lagrange(
    &self,
    dim: usize,
    ncomp: usize,
    P: usize,
    Q: usize,
    qmode: QuadMode
) -> Result<Basis<'_>, CeedError>
[src]

Returns a tensor-product Lagrange basis

arguments

  • dim - Topological dimension of element
  • ncomp - Number of field components (1 for scalar fields)
  • P - Number of Gauss-Lobatto nodes in one dimension. The polynomial degree of the resulting Q_k element is k=P-1.
  • Q - Number of quadrature points in one dimension
  • qmode - Distribution of the Q quadrature points (affects order of accuracy for the quadrature)
let b = ceed
    .basis_tensor_H1_Lagrange(2, 1, 3, 4, QuadMode::Gauss)
    .unwrap();

pub fn basis_H1(
    &self,
    topo: ElemTopology,
    ncomp: usize,
    nnodes: usize,
    nqpts: usize,
    interp: &[f64],
    grad: &[f64],
    qref: &[f64],
    qweight: &[f64]
) -> Result<Basis<'_>, CeedError>
[src]

Returns a tensor-product basis

arguments

  • topo - Topology of element, e.g. hypercube, simplex, ect
  • ncomp - Number of field components (1 for scalar fields)
  • nnodes - Total number of nodes
  • nqpts - Total number of quadrature points
  • interp - Row-major (nqpts * nnodes) matrix expressing the values of nodal basis functions at quadrature points
  • grad - Row-major (nqpts * dim * nnodes) matrix expressing derivatives of nodal basis functions at quadrature points
  • qref - Array of length nqpts holding the locations of quadrature points on the reference element [-1, 1]
  • qweight - Array of length nqpts holding the quadrature weights on the reference element
let interp = [
    0.12000000,
    0.48000000,
    -0.12000000,
    0.48000000,
    0.16000000,
    -0.12000000,
    -0.12000000,
    0.48000000,
    0.12000000,
    0.16000000,
    0.48000000,
    -0.12000000,
    -0.11111111,
    0.44444444,
    -0.11111111,
    0.44444444,
    0.44444444,
    -0.11111111,
    -0.12000000,
    0.16000000,
    -0.12000000,
    0.48000000,
    0.48000000,
    0.12000000,
];
let grad = [
    -1.40000000,
    1.60000000,
    -0.20000000,
    -0.80000000,
    0.80000000,
    0.00000000,
    0.20000000,
    -1.60000000,
    1.40000000,
    -0.80000000,
    0.80000000,
    0.00000000,
    -0.33333333,
    0.00000000,
    0.33333333,
    -1.33333333,
    1.33333333,
    0.00000000,
    0.20000000,
    0.00000000,
    -0.20000000,
    -2.40000000,
    2.40000000,
    0.00000000,
    -1.40000000,
    -0.80000000,
    0.00000000,
    1.60000000,
    0.80000000,
    -0.20000000,
    0.20000000,
    -2.40000000,
    0.00000000,
    0.00000000,
    2.40000000,
    -0.20000000,
    -0.33333333,
    -1.33333333,
    0.00000000,
    0.00000000,
    1.33333333,
    0.33333333,
    0.20000000,
    -0.80000000,
    0.00000000,
    -1.60000000,
    0.80000000,
    1.40000000,
];
let qref = [
    0.20000000, 0.60000000, 0.33333333, 0.20000000, 0.20000000, 0.20000000, 0.33333333,
    0.60000000,
];
let qweight = [0.26041667, 0.26041667, -0.28125000, 0.26041667];
let b = ceed
    .basis_H1(
        ElemTopology::Triangle,
        1,
        6,
        4,
        &interp,
        &grad,
        &qref,
        &qweight,
    )
    .unwrap();

pub fn q_function_interior(
    &self,
    vlength: usize,
    f: Box<QFunctionUserClosure>
) -> Result<QFunction<'_>, CeedError>
[src]

Returns a CeedQFunction for evaluating interior (volumetric) terms

arguments

  • vlength - Vector length. Caller must ensure that number of quadrature points is a multiple of vlength.
  • f - Boxed closure to evaluate action at quadrature points.
let mut user_f = |[u, weights, ..]: QFunctionInputs, [v, ..]: QFunctionOutputs| {
    // Iterate over quadrature points
    v.iter_mut()
        .zip(u.iter().zip(weights.iter()))
        .for_each(|(v, (u, w))| *v = u * w);

    // Return clean error code
    0
};

let qf = ceed.q_function_interior(1, Box::new(user_f)).unwrap();

pub fn q_function_interior_by_name(
    &self,
    name: &str
) -> Result<QFunctionByName<'_>, CeedError>
[src]

Returns a CeedQFunction for evaluating interior (volumetric) terms created by name

let qf = ceed.q_function_interior_by_name("Mass1DBuild").unwrap();

pub fn operator<'b>(
    &self,
    qf: impl Into<QFunctionOpt<'b>>,
    dqf: impl Into<QFunctionOpt<'b>>,
    dqfT: impl Into<QFunctionOpt<'b>>
) -> Result<Operator<'_>, CeedError>
[src]

Returns a Operator and associate a QFunction. A Basis and ElemRestriction can be associated with QFunction fields with set_field().

  • qf - QFunction defining the action of the operator at quadrature points
  • dqf - QFunction defining the action of the Jacobian of the qf (or qfunction_none)
  • dqfT - QFunction defining the action of the transpose of the Jacobian of the qf (or qfunction_none)
let qf = ceed.q_function_interior_by_name("Mass1DBuild").unwrap();
let op = ceed
    .operator(&qf, QFunctionOpt::None, QFunctionOpt::None)
    .unwrap();

pub fn composite_operator(&self) -> Result<CompositeOperator<'_>, CeedError>[src]

Returns an Operator that composes the action of several Operators

let op = ceed.composite_operator().unwrap();

Trait Implementations

impl Debug for Ceed[src]

impl Display for Ceed[src]

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

View a Ceed

let ceed = libceed::Ceed::init("/cpu/self/ref/serial");
println!("{}", ceed);

impl Drop for Ceed[src]

Auto Trait Implementations

impl RefUnwindSafe for Ceed

impl !Send for Ceed

impl !Sync for Ceed

impl Unpin for Ceed

impl UnwindSafe for Ceed

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.