Struct libceed::Ceed

source ·
pub struct Ceed { /* private fields */ }
Expand description

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

Implementations§

source§

impl Ceed

source

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

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");
source

pub fn init_with_error_handler(resource: &str, handler: ErrorHandler) -> Self

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::ErrorHandler::ErrorAbort,
);
source

pub fn resource(&self) -> String

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())
source

pub fn vector<'a>(&self, n: usize) -> Result<Vector<'a>>

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

arguments
  • n - Length of vector
let vec = ceed.vector(10)?;
source

pub fn vector_from_slice<'a>(&self, slice: &[Scalar]) -> Result<Vector<'a>>

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.])?;
assert_eq!(vec.length(), 3);
source

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

Returns an ElemRestriction, $\mathcal{E}$, which extracts the degrees of freedom for each element from the local vector into the element vector or assembles contributions from each element in the element vector to the local vector

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)?;
source

pub fn oriented_elem_restriction<'a>( &self, nelem: usize, elemsize: usize, ncomp: usize, compstride: usize, lsize: usize, mtype: MemType, offsets: &[i32], orients: &[bool] ) -> Result<ElemRestriction<'a>>

Returns an oriented ElemRestriction, $\mathcal{E}$, which extracts the degrees of freedom for each element from the local vector into the element vector or assembles contributions from each element in the element vector to the local vector

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].
  • orients - Array of shape [nelem, elemsize]. Row i holds the ordered list of the orientations for the unknowns corresponding to element i, with bool false used for positively oriented and true to flip the orientation.
let nelem = 3;
let mut ind: Vec<i32> = vec![0; 2 * nelem];
let mut orients: Vec<bool> = vec![false; 2 * nelem];
for i in 0..nelem {
    ind[2 * i + 0] = i as i32;
    ind[2 * i + 1] = (i + 1) as i32;
    orients[2 * i + 0] = (i % 2) > 0; // flip the dofs on element 1, 3, ...
    orients[2 * i + 1] = (i % 2) > 0;
}
let r =
    ceed.oriented_elem_restriction(nelem, 2, 1, 1, nelem + 1, MemType::Host, &ind, &orients)?;
source

pub fn curl_oriented_elem_restriction<'a>( &self, nelem: usize, elemsize: usize, ncomp: usize, compstride: usize, lsize: usize, mtype: MemType, offsets: &[i32], curlorients: &[i8] ) -> Result<ElemRestriction<'a>>

Returns a curl-oriented ElemRestriction, $\mathcal{E}$, which extracts the degrees of freedom for each element from the local vector into the element vector or assembles contributions from each element in the element vector to the local vector

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].
  • curlorients - Array of shape [nelem, 3 * elemsize]. Row i holds a row-major tridiagonal matrix (curlorients[i, 0] = curlorients[i, 3 * elemsize - 1] = 0, where 0 <= i < nelem) which is applied to the element unknowns upon restriction.
let nelem = 3;
let mut ind: Vec<i32> = vec![0; 2 * nelem];
let mut curlorients: Vec<i8> = vec![0; 3 * 2 * nelem];
for i in 0..nelem {
    ind[2 * i + 0] = i as i32;
    ind[2 * i + 1] = (i + 1) as i32;
    curlorients[3 * 2 * i] = 0 as i8;
    curlorients[3 * 2 * (i + 1) - 1] = 0 as i8;
    if (i % 2 > 0) {
        // T = [0  -1]
        //     [-1  0]
        curlorients[3 * 2 * i + 1] = 0 as i8;
        curlorients[3 * 2 * i + 2] = -1 as i8;
        curlorients[3 * 2 * i + 3] = -1 as i8;
        curlorients[3 * 2 * i + 4] = 0 as i8;
    } else {
        // T = I
        curlorients[3 * 2 * i + 1] = 1 as i8;
        curlorients[3 * 2 * i + 2] = 0 as i8;
        curlorients[3 * 2 * i + 3] = 0 as i8;
        curlorients[3 * 2 * i + 4] = 1 as i8;
    }
}
let r = ceed.curl_oriented_elem_restriction(
    nelem,
    2,
    1,
    1,
    nelem + 1,
    MemType::Host,
    &ind,
    &curlorients,
)?;
source

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

Returns an ElemRestriction, $\mathcal{E}$, from an local vector to an element vector where data can be indexed from the strides array

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)
  • 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)?;
source

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

Returns an $H^1$ 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)?;
source

pub fn basis_tensor_H1_Lagrange<'a>( &self, dim: usize, ncomp: usize, P: usize, Q: usize, qmode: QuadMode ) -> Result<Basis<'a>>

Returns an $H^1$ Lagrange tensor-product 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)?;
source

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

Returns an $H-1$ 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 (dim * nqpts * 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,
)?;
source

pub fn basis_Hdiv<'a>( &self, topo: ElemTopology, ncomp: usize, nnodes: usize, nqpts: usize, interp: &[Scalar], div: &[Scalar], qref: &[Scalar], qweight: &[Scalar] ) -> Result<Basis<'a>>

Returns an $H(div)$ 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 (dim * nqpts * nnodes) matrix expressing the values of basis functions at quadrature points
  • div - Row-major (nqpts * nnodes) matrix expressing the divergence of 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.00000000,
    -1.57735027,
    0.57735027,
    0.00000000,
    0.00000000,
    -1.57735027,
    0.57735027,
    0.00000000,
    0.00000000,
    -1.57735027,
    0.57735027,
    0.00000000,
    0.00000000,
    -0.42264973,
    -0.57735027,
    0.00000000,
    0.42264973,
    0.00000000,
    0.00000000,
    0.57735027,
    0.42264973,
    0.00000000,
    0.00000000,
    0.57735027,
    1.57735027,
    0.00000000,
    0.00000000,
    -0.57735027,
    1.57735027,
    0.00000000,
    0.00000000,
    -0.57735027,
];
let div = [
    -1.00000000,
    1.00000000,
    -1.00000000,
    1.00000000,
    -1.00000000,
    1.00000000,
    -1.00000000,
    1.00000000,
    -1.00000000,
    1.00000000,
    -1.00000000,
    1.00000000,
    -1.00000000,
    1.00000000,
    -1.00000000,
    1.00000000,
];
let qref = [
    0.57735026, 0.57735026, 0.57735026, 0.57735026, 0.57735026, 0.57735026, 0.57735026,
    0.57735026,
];
let qweight = [1.00000000, 1.00000000, 1.00000000, 1.00000000];
let b = ceed.basis_Hdiv(ElemTopology::Quad, 1, 4, 4, &interp, &div, &qref, &qweight)?;
source

pub fn basis_Hcurl<'a>( &self, topo: ElemTopology, ncomp: usize, nnodes: usize, nqpts: usize, interp: &[Scalar], curl: &[Scalar], qref: &[Scalar], qweight: &[Scalar] ) -> Result<Basis<'a>>

Returns an $H(curl)$ 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 (dim * nqpts * nnodes) matrix expressing the values of basis functions at quadrature points
  • curl - Row-major (curl_comp * nqpts * nnodes), curl_comp = 1 if dim < 3 else dim matrix expressing the curl of 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.20000000,
    0.20000000,
    0.80000000,
    -0.20000000,
    0.20000000,
    0.80000000,
    -0.33333333,
    0.33333333,
    0.66666667,
    -0.60000000,
    0.60000000,
    0.40000000,
    0.20000000,
    0.80000000,
    0.20000000,
    0.60000000,
    0.40000000,
    0.60000000,
    0.33333333,
    0.66666667,
    0.33333333,
    0.20000000,
    0.80000000,
    0.20000000,
];
let curl = [
    2.00000000,
    -2.00000000,
    -2.00000000,
    2.00000000,
    -2.00000000,
    -2.00000000,
    2.00000000,
    -2.00000000,
    -2.00000000,
    2.00000000,
    -2.00000000,
    -2.00000000,
];
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_Hcurl(
    ElemTopology::Triangle,
    1,
    3,
    4,
    &interp,
    &curl,
    &qref,
    &qweight,
)?;
source

pub fn q_function_interior<'a>( &self, vlength: usize, f: Box<QFunctionUserClosure> ) -> Result<QFunction<'a>>

Returns a QFunction for evaluating interior (volumetric) terms of a weak form corresponding to the $L^2$ inner product

$$ \langle v, F(u) \rangle = \int_\Omega v \cdot f_0 \left( u, \nabla u \right) + \left( \nabla v \right) : f_1 \left( u, \nabla u \right), $$

where $v \cdot f_0$ represents contraction over fields and $\nabla v : f_1$ represents contraction over both fields and spatial dimensions.

arguments
  • vlength - Vector length. Caller must ensure that number of quadrature points is a multiple of vlength.
  • f - Boxed closure to evaluate weak form 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))?;
source

pub fn q_function_interior_by_name<'a>( &self, name: &str ) -> Result<QFunctionByName<'a>>

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

arguments
  • name - name of QFunction from libCEED gallery
let qf = ceed.q_function_interior_by_name("Mass1DBuild")?;
source

pub fn operator<'a, 'b>( &self, qf: impl Into<QFunctionOpt<'b>>, dqf: impl Into<QFunctionOpt<'b>>, dqfT: impl Into<QFunctionOpt<'b>> ) -> Result<Operator<'a>>

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

arguments
  • 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")?;
let op = ceed.operator(&qf, QFunctionOpt::None, QFunctionOpt::None)?;
source

pub fn composite_operator<'a>(&self) -> Result<CompositeOperator<'a>>

Returns an Operator that composes the action of several Operators

let op = ceed.composite_operator()?;

Trait Implementations§

source§

impl Clone for Ceed

source§

fn clone(&self) -> Self

Perform a shallow clone of a Ceed context

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

println!(" original:{} \n clone:{}", ceed, ceed_clone);
1.0.0 · source§

fn clone_from(&mut self, source: &Self)

Performs copy-assignment from source. Read more
source§

impl Debug for Ceed

source§

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

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

impl Display for Ceed

source§

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

View a Ceed

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

impl Drop for Ceed

source§

fn drop(&mut self)

Executes the destructor for this type. Read more

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§

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> ToOwned for T
where T: Clone,

§

type Owned = T

The resulting type after obtaining ownership.
source§

fn to_owned(&self) -> T

Creates owned data from borrowed data, usually by cloning. Read more
source§

fn clone_into(&self, target: &mut T)

Uses borrowed data to replace owned data, usually by cloning. Read more
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.