xref: /libCEED/rust/libceed/src/lib.rs (revision c68be7a2e45631197b626561fad53d5b146fcd59)
19df49d7eSJed Brown // Copyright (c) 2017, Lawrence Livermore National Security, LLC. Produced at
29df49d7eSJed Brown // the Lawrence Livermore National Laboratory. LLNL-CODE-734707. All Rights
39df49d7eSJed Brown // reserved. See files LICENSE and NOTICE for details.
49df49d7eSJed Brown //
59df49d7eSJed Brown // This file is part of CEED, a collection of benchmarks, miniapps, software
69df49d7eSJed Brown // libraries and APIs for efficient high-order finite element and spectral
79df49d7eSJed Brown // element discretizations for exascale applications. For more information and
89df49d7eSJed Brown // source code availability see http://github.com/ceed.
99df49d7eSJed Brown //
109df49d7eSJed Brown // The CEED research is supported by the Exascale Computing Project 17-SC-20-SC,
119df49d7eSJed Brown // a collaborative effort of two U.S. Department of Energy organizations (Office
129df49d7eSJed Brown // of Science and the National Nuclear Security Administration) responsible for
139df49d7eSJed Brown // the planning and preparation of a capable exascale ecosystem, including
149df49d7eSJed Brown // software, applications, hardware, advanced system engineering and early
159df49d7eSJed Brown // testbed platforms, in support of the nation's exascale computing imperative
169df49d7eSJed Brown 
17a1cbad85SJed Brown // Fenced `rust` code blocks included from README.md are executed as part of doctests.
18a1cbad85SJed Brown #![doc = include_str!("../README.md")]
199df49d7eSJed Brown // -----------------------------------------------------------------------------
209df49d7eSJed Brown // Exceptions
219df49d7eSJed Brown // -----------------------------------------------------------------------------
229df49d7eSJed Brown #![allow(non_snake_case)]
239df49d7eSJed Brown 
249df49d7eSJed Brown // -----------------------------------------------------------------------------
259df49d7eSJed Brown // Crate prelude
269df49d7eSJed Brown // -----------------------------------------------------------------------------
279df49d7eSJed Brown use crate::prelude::*;
289df49d7eSJed Brown use std::sync::Once;
299df49d7eSJed Brown 
309df49d7eSJed Brown pub mod prelude {
319df49d7eSJed Brown     pub use crate::{
329df49d7eSJed Brown         basis::{self, Basis, BasisOpt},
339df49d7eSJed Brown         elem_restriction::{self, ElemRestriction, ElemRestrictionOpt},
349df49d7eSJed Brown         operator::{self, CompositeOperator, Operator},
359df49d7eSJed Brown         qfunction::{
369df49d7eSJed Brown             self, QFunction, QFunctionByName, QFunctionInputs, QFunctionOpt, QFunctionOutputs,
379df49d7eSJed Brown         },
389df49d7eSJed Brown         vector::{self, Vector, VectorOpt},
3980a9ef05SNatalie Beams         ElemTopology, EvalMode, MemType, NormType, QuadMode, Scalar, TransposeMode,
4080a9ef05SNatalie Beams         CEED_STRIDES_BACKEND, EPSILON, MAX_QFUNCTION_FIELDS,
419df49d7eSJed Brown     };
42630ad4c9Sjeremylt     pub(crate) use libceed_sys::bind_ceed;
439df49d7eSJed Brown     pub(crate) use std::convert::TryFrom;
449df49d7eSJed Brown     pub(crate) use std::ffi::CString;
459df49d7eSJed Brown     pub(crate) use std::fmt;
469df49d7eSJed Brown }
479df49d7eSJed Brown 
489df49d7eSJed Brown // -----------------------------------------------------------------------------
499df49d7eSJed Brown // Modules
509df49d7eSJed Brown // -----------------------------------------------------------------------------
519df49d7eSJed Brown pub mod basis;
529df49d7eSJed Brown pub mod elem_restriction;
539df49d7eSJed Brown pub mod operator;
549df49d7eSJed Brown pub mod qfunction;
559df49d7eSJed Brown pub mod vector;
569df49d7eSJed Brown 
579df49d7eSJed Brown // -----------------------------------------------------------------------------
5880a9ef05SNatalie Beams // Typedef for scalar
5980a9ef05SNatalie Beams // -----------------------------------------------------------------------------
6080a9ef05SNatalie Beams pub type Scalar = bind_ceed::CeedScalar;
6180a9ef05SNatalie Beams 
6280a9ef05SNatalie Beams // -----------------------------------------------------------------------------
639df49d7eSJed Brown // Constants for library
649df49d7eSJed Brown // -----------------------------------------------------------------------------
659df49d7eSJed Brown const MAX_BUFFER_LENGTH: u64 = 4096;
669df49d7eSJed Brown pub const MAX_QFUNCTION_FIELDS: usize = 16;
679df49d7eSJed Brown pub const CEED_STRIDES_BACKEND: [i32; 3] = [0; 3];
6880a9ef05SNatalie Beams pub const EPSILON: crate::Scalar = bind_ceed::CEED_EPSILON as crate::Scalar;
699df49d7eSJed Brown 
709df49d7eSJed Brown // -----------------------------------------------------------------------------
719df49d7eSJed Brown // Enums for libCEED
729df49d7eSJed Brown // -----------------------------------------------------------------------------
739df49d7eSJed Brown #[derive(Clone, Copy, PartialEq, Eq)]
749df49d7eSJed Brown /// Many Ceed interfaces take or return pointers to memory.  This enum is used to
759df49d7eSJed Brown /// specify where the memory being provided or requested must reside.
769df49d7eSJed Brown pub enum MemType {
779df49d7eSJed Brown     Host = bind_ceed::CeedMemType_CEED_MEM_HOST as isize,
789df49d7eSJed Brown     Device = bind_ceed::CeedMemType_CEED_MEM_DEVICE as isize,
799df49d7eSJed Brown }
809df49d7eSJed Brown 
819df49d7eSJed Brown #[derive(Clone, Copy, PartialEq, Eq)]
829df49d7eSJed Brown // OwnPointer will not be used by user but is included for internal use
839df49d7eSJed Brown #[allow(dead_code)]
849df49d7eSJed Brown /// Conveys ownership status of arrays passed to Ceed interfaces.
859df49d7eSJed Brown pub(crate) enum CopyMode {
869df49d7eSJed Brown     CopyValues = bind_ceed::CeedCopyMode_CEED_COPY_VALUES as isize,
879df49d7eSJed Brown     UsePointer = bind_ceed::CeedCopyMode_CEED_USE_POINTER as isize,
889df49d7eSJed Brown     OwnPointer = bind_ceed::CeedCopyMode_CEED_OWN_POINTER as isize,
899df49d7eSJed Brown }
909df49d7eSJed Brown 
919df49d7eSJed Brown #[derive(Clone, Copy, PartialEq, Eq)]
929df49d7eSJed Brown /// Denotes type of vector norm to be computed
939df49d7eSJed Brown pub enum NormType {
949df49d7eSJed Brown     One = bind_ceed::CeedNormType_CEED_NORM_1 as isize,
959df49d7eSJed Brown     Two = bind_ceed::CeedNormType_CEED_NORM_2 as isize,
969df49d7eSJed Brown     Max = bind_ceed::CeedNormType_CEED_NORM_MAX as isize,
979df49d7eSJed Brown }
989df49d7eSJed Brown 
999df49d7eSJed Brown #[derive(Clone, Copy, PartialEq, Eq)]
1009df49d7eSJed Brown /// Denotes whether a linear transformation or its transpose should be applied
1019df49d7eSJed Brown pub enum TransposeMode {
1029df49d7eSJed Brown     NoTranspose = bind_ceed::CeedTransposeMode_CEED_NOTRANSPOSE as isize,
1039df49d7eSJed Brown     Transpose = bind_ceed::CeedTransposeMode_CEED_TRANSPOSE as isize,
1049df49d7eSJed Brown }
1059df49d7eSJed Brown 
1069df49d7eSJed Brown #[derive(Clone, Copy, PartialEq, Eq)]
1079df49d7eSJed Brown /// Type of quadrature; also used for location of nodes
1089df49d7eSJed Brown pub enum QuadMode {
1099df49d7eSJed Brown     Gauss = bind_ceed::CeedQuadMode_CEED_GAUSS as isize,
1109df49d7eSJed Brown     GaussLobatto = bind_ceed::CeedQuadMode_CEED_GAUSS_LOBATTO as isize,
1119df49d7eSJed Brown }
1129df49d7eSJed Brown 
1139df49d7eSJed Brown #[derive(Clone, Copy, PartialEq, Eq)]
1149df49d7eSJed Brown /// Type of basis shape to create non-tensor H1 element basis
1159df49d7eSJed Brown pub enum ElemTopology {
1169df49d7eSJed Brown     Line = bind_ceed::CeedElemTopology_CEED_LINE as isize,
1179df49d7eSJed Brown     Triangle = bind_ceed::CeedElemTopology_CEED_TRIANGLE as isize,
1189df49d7eSJed Brown     Quad = bind_ceed::CeedElemTopology_CEED_QUAD as isize,
1199df49d7eSJed Brown     Tet = bind_ceed::CeedElemTopology_CEED_TET as isize,
1209df49d7eSJed Brown     Pyramid = bind_ceed::CeedElemTopology_CEED_PYRAMID as isize,
1219df49d7eSJed Brown     Prism = bind_ceed::CeedElemTopology_CEED_PRISM as isize,
1229df49d7eSJed Brown     Hex = bind_ceed::CeedElemTopology_CEED_HEX as isize,
1239df49d7eSJed Brown }
1249df49d7eSJed Brown 
125*c68be7a2SJeremy L Thompson #[derive(Clone, Copy, Debug, PartialEq, Eq)]
1269df49d7eSJed Brown /// Basis evaluation mode
1279df49d7eSJed Brown pub enum EvalMode {
1289df49d7eSJed Brown     None = bind_ceed::CeedEvalMode_CEED_EVAL_NONE as isize,
1299df49d7eSJed Brown     Interp = bind_ceed::CeedEvalMode_CEED_EVAL_INTERP as isize,
1309df49d7eSJed Brown     Grad = bind_ceed::CeedEvalMode_CEED_EVAL_GRAD as isize,
1319df49d7eSJed Brown     Div = bind_ceed::CeedEvalMode_CEED_EVAL_DIV as isize,
1329df49d7eSJed Brown     Curl = bind_ceed::CeedEvalMode_CEED_EVAL_CURL as isize,
1339df49d7eSJed Brown     Weight = bind_ceed::CeedEvalMode_CEED_EVAL_WEIGHT as isize,
1349df49d7eSJed Brown }
135*c68be7a2SJeremy L Thompson impl EvalMode {
136*c68be7a2SJeremy L Thompson     pub(crate) fn from_u32(value: u32) -> EvalMode {
137*c68be7a2SJeremy L Thompson         match value {
138*c68be7a2SJeremy L Thompson             bind_ceed::CeedEvalMode_CEED_EVAL_NONE => EvalMode::None,
139*c68be7a2SJeremy L Thompson             bind_ceed::CeedEvalMode_CEED_EVAL_INTERP => EvalMode::Interp,
140*c68be7a2SJeremy L Thompson             bind_ceed::CeedEvalMode_CEED_EVAL_GRAD => EvalMode::Grad,
141*c68be7a2SJeremy L Thompson             bind_ceed::CeedEvalMode_CEED_EVAL_DIV => EvalMode::Div,
142*c68be7a2SJeremy L Thompson             bind_ceed::CeedEvalMode_CEED_EVAL_CURL => EvalMode::Curl,
143*c68be7a2SJeremy L Thompson             bind_ceed::CeedEvalMode_CEED_EVAL_WEIGHT => EvalMode::Weight,
144*c68be7a2SJeremy L Thompson             _ => panic!("Unknown value: {}", value),
145*c68be7a2SJeremy L Thompson         }
146*c68be7a2SJeremy L Thompson     }
147*c68be7a2SJeremy L Thompson }
1489df49d7eSJed Brown 
1499df49d7eSJed Brown // -----------------------------------------------------------------------------
1509df49d7eSJed Brown // Ceed error
1519df49d7eSJed Brown // -----------------------------------------------------------------------------
1529df49d7eSJed Brown type Result<T> = std::result::Result<T, CeedError>;
1539df49d7eSJed Brown 
1549df49d7eSJed Brown #[derive(Debug)]
1559df49d7eSJed Brown pub struct CeedError {
1569df49d7eSJed Brown     pub message: String,
1579df49d7eSJed Brown }
1589df49d7eSJed Brown 
1599df49d7eSJed Brown impl fmt::Display for CeedError {
1609df49d7eSJed Brown     fn fmt(&self, f: &mut fmt::Formatter) -> fmt::Result {
1619df49d7eSJed Brown         write!(f, "{}", self.message)
1629df49d7eSJed Brown     }
1639df49d7eSJed Brown }
1649df49d7eSJed Brown 
1659df49d7eSJed Brown // -----------------------------------------------------------------------------
1669df49d7eSJed Brown // Ceed error handler
1679df49d7eSJed Brown // -----------------------------------------------------------------------------
1689df49d7eSJed Brown pub enum CeedErrorHandler {
1699df49d7eSJed Brown     ErrorAbort,
1709df49d7eSJed Brown     ErrorExit,
1719df49d7eSJed Brown     ErrorReturn,
1729df49d7eSJed Brown     ErrorStore,
1739df49d7eSJed Brown }
1749df49d7eSJed Brown 
1759df49d7eSJed Brown // -----------------------------------------------------------------------------
1769df49d7eSJed Brown // Ceed context wrapper
1779df49d7eSJed Brown // -----------------------------------------------------------------------------
1789df49d7eSJed Brown /// A Ceed is a library context representing control of a logical hardware
1799df49d7eSJed Brown /// resource.
1809df49d7eSJed Brown #[derive(Debug)]
1819df49d7eSJed Brown pub struct Ceed {
1829df49d7eSJed Brown     ptr: bind_ceed::Ceed,
1839df49d7eSJed Brown }
1849df49d7eSJed Brown 
1859df49d7eSJed Brown // -----------------------------------------------------------------------------
1869df49d7eSJed Brown // Destructor
1879df49d7eSJed Brown // -----------------------------------------------------------------------------
1889df49d7eSJed Brown impl Drop for Ceed {
1899df49d7eSJed Brown     fn drop(&mut self) {
1909df49d7eSJed Brown         unsafe {
1919df49d7eSJed Brown             bind_ceed::CeedDestroy(&mut self.ptr);
1929df49d7eSJed Brown         }
1939df49d7eSJed Brown     }
1949df49d7eSJed Brown }
1959df49d7eSJed Brown 
1969df49d7eSJed Brown // -----------------------------------------------------------------------------
1979df49d7eSJed Brown // Display
1989df49d7eSJed Brown // -----------------------------------------------------------------------------
1999df49d7eSJed Brown impl fmt::Display for Ceed {
2009df49d7eSJed Brown     /// View a Ceed
2019df49d7eSJed Brown     ///
2029df49d7eSJed Brown     /// ```
2039df49d7eSJed Brown     /// let ceed = libceed::Ceed::init("/cpu/self/ref/serial");
2049df49d7eSJed Brown     /// println!("{}", ceed);
2059df49d7eSJed Brown     /// ```
2069df49d7eSJed Brown     fn fmt(&self, f: &mut fmt::Formatter) -> fmt::Result {
2079df49d7eSJed Brown         let mut ptr = std::ptr::null_mut();
2089df49d7eSJed Brown         let mut sizeloc = crate::MAX_BUFFER_LENGTH;
2099df49d7eSJed Brown         let cstring = unsafe {
2109df49d7eSJed Brown             let file = bind_ceed::open_memstream(&mut ptr, &mut sizeloc);
2119df49d7eSJed Brown             bind_ceed::CeedView(self.ptr, file);
2129df49d7eSJed Brown             bind_ceed::fclose(file);
2139df49d7eSJed Brown             CString::from_raw(ptr)
2149df49d7eSJed Brown         };
2159df49d7eSJed Brown         cstring.to_string_lossy().fmt(f)
2169df49d7eSJed Brown     }
2179df49d7eSJed Brown }
2189df49d7eSJed Brown 
2199df49d7eSJed Brown static REGISTER: Once = Once::new();
2209df49d7eSJed Brown 
2219df49d7eSJed Brown // -----------------------------------------------------------------------------
2229df49d7eSJed Brown // Object constructors
2239df49d7eSJed Brown // -----------------------------------------------------------------------------
2249df49d7eSJed Brown impl Ceed {
2259df49d7eSJed Brown     /// Returns a Ceed context initialized with the specified resource
2269df49d7eSJed Brown     ///
2279df49d7eSJed Brown     /// # arguments
2289df49d7eSJed Brown     ///
2299df49d7eSJed Brown     /// * `resource` - Resource to use, e.g., "/cpu/self"
2309df49d7eSJed Brown     ///
2319df49d7eSJed Brown     /// ```
2329df49d7eSJed Brown     /// let ceed = libceed::Ceed::init("/cpu/self/ref/serial");
2339df49d7eSJed Brown     /// ```
2349df49d7eSJed Brown     pub fn init(resource: &str) -> Self {
2359df49d7eSJed Brown         Ceed::init_with_error_handler(resource, CeedErrorHandler::ErrorStore)
2369df49d7eSJed Brown     }
2379df49d7eSJed Brown 
2389df49d7eSJed Brown     /// Returns a Ceed context initialized with the specified resource
2399df49d7eSJed Brown     ///
2409df49d7eSJed Brown     /// # arguments
2419df49d7eSJed Brown     ///
2429df49d7eSJed Brown     /// * `resource` - Resource to use, e.g., "/cpu/self"
2439df49d7eSJed Brown     ///
2449df49d7eSJed Brown     /// ```
2459df49d7eSJed Brown     /// let ceed = libceed::Ceed::init_with_error_handler(
2469df49d7eSJed Brown     ///     "/cpu/self/ref/serial",
2479df49d7eSJed Brown     ///     libceed::CeedErrorHandler::ErrorAbort,
2489df49d7eSJed Brown     /// );
2499df49d7eSJed Brown     /// ```
2509df49d7eSJed Brown     pub fn init_with_error_handler(resource: &str, handler: CeedErrorHandler) -> Self {
2519df49d7eSJed Brown         REGISTER.call_once(|| unsafe {
2529df49d7eSJed Brown             bind_ceed::CeedRegisterAll();
2539df49d7eSJed Brown             bind_ceed::CeedQFunctionRegisterAll();
2549df49d7eSJed Brown         });
2559df49d7eSJed Brown 
2569df49d7eSJed Brown         // Convert to C string
2579df49d7eSJed Brown         let c_resource = CString::new(resource).expect("CString::new failed");
2589df49d7eSJed Brown 
2599df49d7eSJed Brown         // Get error handler pointer
2609df49d7eSJed Brown         let eh = match handler {
2619df49d7eSJed Brown             CeedErrorHandler::ErrorAbort => bind_ceed::CeedErrorAbort,
2629df49d7eSJed Brown             CeedErrorHandler::ErrorExit => bind_ceed::CeedErrorExit,
2639df49d7eSJed Brown             CeedErrorHandler::ErrorReturn => bind_ceed::CeedErrorReturn,
2649df49d7eSJed Brown             CeedErrorHandler::ErrorStore => bind_ceed::CeedErrorStore,
2659df49d7eSJed Brown         };
2669df49d7eSJed Brown 
2679df49d7eSJed Brown         // Call to libCEED
2689df49d7eSJed Brown         let mut ptr = std::ptr::null_mut();
2699df49d7eSJed Brown         let mut ierr = unsafe { bind_ceed::CeedInit(c_resource.as_ptr() as *const i8, &mut ptr) };
2709df49d7eSJed Brown         if ierr != 0 {
2719df49d7eSJed Brown             panic!("Error initializing backend resource: {}", resource)
2729df49d7eSJed Brown         }
2739df49d7eSJed Brown         ierr = unsafe { bind_ceed::CeedSetErrorHandler(ptr, Some(eh)) };
2749df49d7eSJed Brown         let ceed = Ceed { ptr };
2759df49d7eSJed Brown         ceed.check_error(ierr).unwrap();
2769df49d7eSJed Brown         ceed
2779df49d7eSJed Brown     }
2789df49d7eSJed Brown 
2799df49d7eSJed Brown     /// Default initializer for testing
2809df49d7eSJed Brown     #[doc(hidden)]
2819df49d7eSJed Brown     pub fn default_init() -> Self {
2829df49d7eSJed Brown         // Convert to C string
2839df49d7eSJed Brown         let resource = "/cpu/self/ref/serial";
2849df49d7eSJed Brown         crate::Ceed::init(resource)
2859df49d7eSJed Brown     }
2869df49d7eSJed Brown 
2879df49d7eSJed Brown     /// Internal error checker
2889df49d7eSJed Brown     #[doc(hidden)]
2899df49d7eSJed Brown     fn check_error(&self, ierr: i32) -> Result<i32> {
2909df49d7eSJed Brown         // Return early if code is clean
2919df49d7eSJed Brown         if ierr == bind_ceed::CeedErrorType_CEED_ERROR_SUCCESS {
2929df49d7eSJed Brown             return Ok(ierr);
2939df49d7eSJed Brown         }
2949df49d7eSJed Brown         // Retrieve error message
2959df49d7eSJed Brown         let mut ptr: *const std::os::raw::c_char = std::ptr::null_mut();
2969df49d7eSJed Brown         let c_str = unsafe {
2979df49d7eSJed Brown             bind_ceed::CeedGetErrorMessage(self.ptr, &mut ptr);
2989df49d7eSJed Brown             std::ffi::CStr::from_ptr(ptr)
2999df49d7eSJed Brown         };
3009df49d7eSJed Brown         let message = c_str.to_string_lossy().to_string();
3019df49d7eSJed Brown         // Panic if negative code, otherwise return error
3029df49d7eSJed Brown         if ierr < bind_ceed::CeedErrorType_CEED_ERROR_SUCCESS {
3039df49d7eSJed Brown             panic!("{}", message);
3049df49d7eSJed Brown         }
3059df49d7eSJed Brown         Err(CeedError { message })
3069df49d7eSJed Brown     }
3079df49d7eSJed Brown 
3089df49d7eSJed Brown     /// Returns full resource name for a Ceed context
3099df49d7eSJed Brown     ///
3109df49d7eSJed Brown     /// ```
3119df49d7eSJed Brown     /// let ceed = libceed::Ceed::init("/cpu/self/ref/serial");
3129df49d7eSJed Brown     /// let resource = ceed.resource();
3139df49d7eSJed Brown     ///
3149df49d7eSJed Brown     /// assert_eq!(resource, "/cpu/self/ref/serial".to_string())
3159df49d7eSJed Brown     /// ```
3169df49d7eSJed Brown     pub fn resource(&self) -> String {
3179df49d7eSJed Brown         let mut ptr: *const std::os::raw::c_char = std::ptr::null_mut();
3189df49d7eSJed Brown         let c_str = unsafe {
3199df49d7eSJed Brown             bind_ceed::CeedGetResource(self.ptr, &mut ptr);
3209df49d7eSJed Brown             std::ffi::CStr::from_ptr(ptr)
3219df49d7eSJed Brown         };
3229df49d7eSJed Brown         c_str.to_string_lossy().to_string()
3239df49d7eSJed Brown     }
3249df49d7eSJed Brown 
3259df49d7eSJed Brown     /// Returns a CeedVector of the specified length (does not allocate memory)
3269df49d7eSJed Brown     ///
3279df49d7eSJed Brown     /// # arguments
3289df49d7eSJed Brown     ///
3299df49d7eSJed Brown     /// * `n` - Length of vector
3309df49d7eSJed Brown     ///
3319df49d7eSJed Brown     /// ```
3329df49d7eSJed Brown     /// # use libceed::prelude::*;
333*c68be7a2SJeremy L Thompson     /// # fn main() -> Result<(), libceed::CeedError> {
3349df49d7eSJed Brown     /// # let ceed = libceed::Ceed::default_init();
335*c68be7a2SJeremy L Thompson     /// let vec = ceed.vector(10)?;
336*c68be7a2SJeremy L Thompson     /// # Ok(())
337*c68be7a2SJeremy L Thompson     /// # }
3389df49d7eSJed Brown     /// ```
3399df49d7eSJed Brown     pub fn vector(&self, n: usize) -> Result<Vector> {
3409df49d7eSJed Brown         Vector::create(self, n)
3419df49d7eSJed Brown     }
3429df49d7eSJed Brown 
3439df49d7eSJed Brown     /// Create a Vector initialized with the data (copied) from a slice
3449df49d7eSJed Brown     ///
3459df49d7eSJed Brown     /// # arguments
3469df49d7eSJed Brown     ///
3479df49d7eSJed Brown     /// * `slice` - Slice containing data
3489df49d7eSJed Brown     ///
3499df49d7eSJed Brown     /// ```
3509df49d7eSJed Brown     /// # use libceed::prelude::*;
351*c68be7a2SJeremy L Thompson     /// # fn main() -> Result<(), libceed::CeedError> {
3529df49d7eSJed Brown     /// # let ceed = libceed::Ceed::default_init();
353*c68be7a2SJeremy L Thompson     /// let vec = ceed.vector_from_slice(&[1., 2., 3.])?;
3549df49d7eSJed Brown     /// assert_eq!(vec.length(), 3);
355*c68be7a2SJeremy L Thompson     /// # Ok(())
356*c68be7a2SJeremy L Thompson     /// # }
3579df49d7eSJed Brown     /// ```
35880a9ef05SNatalie Beams     pub fn vector_from_slice(&self, slice: &[crate::Scalar]) -> Result<Vector> {
3599df49d7eSJed Brown         Vector::from_slice(self, slice)
3609df49d7eSJed Brown     }
3619df49d7eSJed Brown 
3629df49d7eSJed Brown     /// Returns a ElemRestriction
3639df49d7eSJed Brown     ///
3649df49d7eSJed Brown     /// # arguments
3659df49d7eSJed Brown     ///
3669df49d7eSJed Brown     /// * `nelem`      - Number of elements described in the offsets array
3679df49d7eSJed Brown     /// * `elemsize`   - Size (number of "nodes") per element
3689df49d7eSJed Brown     /// * `ncomp`      - Number of field components per interpolation node (1
3699df49d7eSJed Brown     ///                    for scalar fields)
3709df49d7eSJed Brown     /// * `compstride` - Stride between components for the same Lvector "node".
3719df49d7eSJed Brown     ///                    Data for node `i`, component `j`, element `k` can be
3729df49d7eSJed Brown     ///                    found in the Lvector at index
3739df49d7eSJed Brown     ///                    `offsets[i + k*elemsize] + j*compstride`.
3749df49d7eSJed Brown     /// * `lsize`      - The size of the Lvector. This vector may be larger
3759df49d7eSJed Brown     ///                    than the elements and fields given by this
3769df49d7eSJed Brown     ///                    restriction.
3779df49d7eSJed Brown     /// * `mtype`     - Memory type of the offsets array, see CeedMemType
3789df49d7eSJed Brown     /// * `offsets`    - Array of shape `[nelem, elemsize]`. Row `i` holds the
3799df49d7eSJed Brown     ///                    ordered list of the offsets (into the input CeedVector)
3809df49d7eSJed Brown     ///                    for the unknowns corresponding to element `i`, where
3819df49d7eSJed Brown     ///                    `0 <= i < nelem`. All offsets must be in the range
3829df49d7eSJed Brown     ///                    `[0, lsize - 1]`.
3839df49d7eSJed Brown     ///
3849df49d7eSJed Brown     /// ```
3859df49d7eSJed Brown     /// # use libceed::prelude::*;
386*c68be7a2SJeremy L Thompson     /// # fn main() -> Result<(), libceed::CeedError> {
3879df49d7eSJed Brown     /// # let ceed = libceed::Ceed::default_init();
3889df49d7eSJed Brown     /// let nelem = 3;
3899df49d7eSJed Brown     /// let mut ind: Vec<i32> = vec![0; 2 * nelem];
3909df49d7eSJed Brown     /// for i in 0..nelem {
3919df49d7eSJed Brown     ///     ind[2 * i + 0] = i as i32;
3929df49d7eSJed Brown     ///     ind[2 * i + 1] = (i + 1) as i32;
3939df49d7eSJed Brown     /// }
394*c68be7a2SJeremy L Thompson     /// let r = ceed.elem_restriction(nelem, 2, 1, 1, nelem + 1, MemType::Host, &ind)?;
395*c68be7a2SJeremy L Thompson     /// # Ok(())
396*c68be7a2SJeremy L Thompson     /// # }
3979df49d7eSJed Brown     /// ```
3989df49d7eSJed Brown     pub fn elem_restriction(
3999df49d7eSJed Brown         &self,
4009df49d7eSJed Brown         nelem: usize,
4019df49d7eSJed Brown         elemsize: usize,
4029df49d7eSJed Brown         ncomp: usize,
4039df49d7eSJed Brown         compstride: usize,
4049df49d7eSJed Brown         lsize: usize,
4059df49d7eSJed Brown         mtype: MemType,
4069df49d7eSJed Brown         offsets: &[i32],
4079df49d7eSJed Brown     ) -> Result<ElemRestriction> {
4089df49d7eSJed Brown         ElemRestriction::create(
4099df49d7eSJed Brown             self, nelem, elemsize, ncomp, compstride, lsize, mtype, offsets,
4109df49d7eSJed Brown         )
4119df49d7eSJed Brown     }
4129df49d7eSJed Brown 
4139df49d7eSJed Brown     /// Returns a ElemRestriction
4149df49d7eSJed Brown     ///
4159df49d7eSJed Brown     /// # arguments
4169df49d7eSJed Brown     ///
4179df49d7eSJed Brown     /// * `nelem`      - Number of elements described in the offsets array
4189df49d7eSJed Brown     /// * `elemsize`   - Size (number of "nodes") per element
4199df49d7eSJed Brown     /// * `ncomp`      - Number of field components per interpolation node (1
4209df49d7eSJed Brown     ///                    for scalar fields)
4219df49d7eSJed Brown     /// * `compstride` - Stride between components for the same Lvector "node".
4229df49d7eSJed Brown     ///                    Data for node `i`, component `j`, element `k` can be
4239df49d7eSJed Brown     ///                    found in the Lvector at index
4249df49d7eSJed Brown     ///                    `offsets[i + k*elemsize] + j*compstride`.
4259df49d7eSJed Brown     /// * `lsize`      - The size of the Lvector. This vector may be larger
4269df49d7eSJed Brown     ///   than the elements and fields given by this restriction.
4279df49d7eSJed Brown     /// * `strides`   - Array for strides between `[nodes, components, elements]`.
4289df49d7eSJed Brown     ///                   Data for node `i`, component `j`, element `k` can be
4299df49d7eSJed Brown     ///                   found in the Lvector at index
4309df49d7eSJed Brown     ///                   `i*strides[0] + j*strides[1] + k*strides[2]`.
4319df49d7eSJed Brown     ///                   CEED_STRIDES_BACKEND may be used with vectors created
4329df49d7eSJed Brown     ///                   by a Ceed backend.
4339df49d7eSJed Brown     ///
4349df49d7eSJed Brown     /// ```
4359df49d7eSJed Brown     /// # use libceed::prelude::*;
436*c68be7a2SJeremy L Thompson     /// # fn main() -> Result<(), libceed::CeedError> {
4379df49d7eSJed Brown     /// # let ceed = libceed::Ceed::default_init();
4389df49d7eSJed Brown     /// let nelem = 3;
4399df49d7eSJed Brown     /// let strides: [i32; 3] = [1, 2, 2];
440*c68be7a2SJeremy L Thompson     /// let r = ceed.strided_elem_restriction(nelem, 2, 1, nelem * 2, strides)?;
441*c68be7a2SJeremy L Thompson     /// # Ok(())
442*c68be7a2SJeremy L Thompson     /// # }
4439df49d7eSJed Brown     /// ```
4449df49d7eSJed Brown     pub fn strided_elem_restriction(
4459df49d7eSJed Brown         &self,
4469df49d7eSJed Brown         nelem: usize,
4479df49d7eSJed Brown         elemsize: usize,
4489df49d7eSJed Brown         ncomp: usize,
4499df49d7eSJed Brown         lsize: usize,
4509df49d7eSJed Brown         strides: [i32; 3],
4519df49d7eSJed Brown     ) -> Result<ElemRestriction> {
4529df49d7eSJed Brown         ElemRestriction::create_strided(self, nelem, elemsize, ncomp, lsize, strides)
4539df49d7eSJed Brown     }
4549df49d7eSJed Brown 
4559df49d7eSJed Brown     /// Returns a tensor-product basis
4569df49d7eSJed Brown     ///
4579df49d7eSJed Brown     /// # arguments
4589df49d7eSJed Brown     ///
4599df49d7eSJed Brown     /// * `dim`       - Topological dimension of element
4609df49d7eSJed Brown     /// * `ncomp`     - Number of field components (1 for scalar fields)
4619df49d7eSJed Brown     /// * `P1d`       - Number of Gauss-Lobatto nodes in one dimension.  The
4629df49d7eSJed Brown     ///                   polynomial degree of the resulting `Q_k` element is
4639df49d7eSJed Brown     ///                   `k=P-1`.
4649df49d7eSJed Brown     /// * `Q1d`       - Number of quadrature points in one dimension
4659df49d7eSJed Brown     /// * `interp1d`  - Row-major `(Q1d * P1d)` matrix expressing the values of
4669df49d7eSJed Brown     ///                   nodal basis functions at quadrature points
4679df49d7eSJed Brown     /// * `grad1d`    - Row-major `(Q1d * P1d)` matrix expressing derivatives of
4689df49d7eSJed Brown     ///                   nodal basis functions at quadrature points
4699df49d7eSJed Brown     /// * `qref1d`    - Array of length `Q1d` holding the locations of quadrature
4709df49d7eSJed Brown     ///                   points on the 1D reference element `[-1, 1]`
4719df49d7eSJed Brown     /// * `qweight1d` - Array of length `Q1d` holding the quadrature weights on
4729df49d7eSJed Brown     ///                   the reference element
4739df49d7eSJed Brown     ///
4749df49d7eSJed Brown     /// ```
4759df49d7eSJed Brown     /// # use libceed::prelude::*;
476*c68be7a2SJeremy L Thompson     /// # fn main() -> Result<(), libceed::CeedError> {
4779df49d7eSJed Brown     /// # let ceed = libceed::Ceed::default_init();
4789df49d7eSJed Brown     /// let interp1d  = [ 0.62994317,  0.47255875, -0.14950343,  0.04700152,
4799df49d7eSJed Brown     ///                  -0.07069480,  0.97297619,  0.13253993, -0.03482132,
4809df49d7eSJed Brown     ///                  -0.03482132,  0.13253993,  0.97297619, -0.07069480,
4819df49d7eSJed Brown     ///                   0.04700152, -0.14950343,  0.47255875,  0.62994317];
4829df49d7eSJed Brown     /// let grad1d    = [-2.34183742,  2.78794489, -0.63510411,  0.18899664,
4839df49d7eSJed Brown     ///                  -0.51670214, -0.48795249,  1.33790510, -0.33325047,
4849df49d7eSJed Brown     //                    0.33325047, -1.33790510,  0.48795249,  0.51670214,
4859df49d7eSJed Brown     ///                  -0.18899664,  0.63510411, -2.78794489,  2.34183742];
4869df49d7eSJed Brown     /// let qref1d    = [-0.86113631, -0.33998104,  0.33998104,  0.86113631];
4879df49d7eSJed Brown     /// let qweight1d = [ 0.34785485,  0.65214515,  0.65214515,  0.34785485];
4889df49d7eSJed Brown     /// let b = ceed.
489*c68be7a2SJeremy L Thompson     /// basis_tensor_H1(2, 1, 4, 4, &interp1d, &grad1d, &qref1d, &qweight1d)?;
490*c68be7a2SJeremy L Thompson     /// # Ok(())
491*c68be7a2SJeremy L Thompson     /// # }
4929df49d7eSJed Brown     /// ```
4939df49d7eSJed Brown     pub fn basis_tensor_H1(
4949df49d7eSJed Brown         &self,
4959df49d7eSJed Brown         dim: usize,
4969df49d7eSJed Brown         ncomp: usize,
4979df49d7eSJed Brown         P1d: usize,
4989df49d7eSJed Brown         Q1d: usize,
49980a9ef05SNatalie Beams         interp1d: &[crate::Scalar],
50080a9ef05SNatalie Beams         grad1d: &[crate::Scalar],
50180a9ef05SNatalie Beams         qref1d: &[crate::Scalar],
50280a9ef05SNatalie Beams         qweight1d: &[crate::Scalar],
5039df49d7eSJed Brown     ) -> Result<Basis> {
5049df49d7eSJed Brown         Basis::create_tensor_H1(
5059df49d7eSJed Brown             self, dim, ncomp, P1d, Q1d, interp1d, grad1d, qref1d, qweight1d,
5069df49d7eSJed Brown         )
5079df49d7eSJed Brown     }
5089df49d7eSJed Brown 
5099df49d7eSJed Brown     /// Returns a tensor-product Lagrange basis
5109df49d7eSJed Brown     ///
5119df49d7eSJed Brown     /// # arguments
5129df49d7eSJed Brown     ///
5139df49d7eSJed Brown     /// * `dim`   - Topological dimension of element
5149df49d7eSJed Brown     /// * `ncomp` - Number of field components (1 for scalar fields)
5159df49d7eSJed Brown     /// * `P`     - Number of Gauss-Lobatto nodes in one dimension.  The
5169df49d7eSJed Brown     ///               polynomial degree of the resulting `Q_k` element is `k=P-1`.
5179df49d7eSJed Brown     /// * `Q`     - Number of quadrature points in one dimension
5189df49d7eSJed Brown     /// * `qmode` - Distribution of the `Q` quadrature points (affects order of
5199df49d7eSJed Brown     ///               accuracy for the quadrature)
5209df49d7eSJed Brown     ///
5219df49d7eSJed Brown     /// ```
5229df49d7eSJed Brown     /// # use libceed::prelude::*;
523*c68be7a2SJeremy L Thompson     /// # fn main() -> Result<(), libceed::CeedError> {
5249df49d7eSJed Brown     /// # let ceed = libceed::Ceed::default_init();
525*c68be7a2SJeremy L Thompson     /// let b = ceed.basis_tensor_H1_Lagrange(2, 1, 3, 4, QuadMode::Gauss)?;
526*c68be7a2SJeremy L Thompson     /// # Ok(())
527*c68be7a2SJeremy L Thompson     /// # }
5289df49d7eSJed Brown     /// ```
5299df49d7eSJed Brown     pub fn basis_tensor_H1_Lagrange(
5309df49d7eSJed Brown         &self,
5319df49d7eSJed Brown         dim: usize,
5329df49d7eSJed Brown         ncomp: usize,
5339df49d7eSJed Brown         P: usize,
5349df49d7eSJed Brown         Q: usize,
5359df49d7eSJed Brown         qmode: QuadMode,
5369df49d7eSJed Brown     ) -> Result<Basis> {
5379df49d7eSJed Brown         Basis::create_tensor_H1_Lagrange(self, dim, ncomp, P, Q, qmode)
5389df49d7eSJed Brown     }
5399df49d7eSJed Brown 
5409df49d7eSJed Brown     /// Returns a tensor-product basis
5419df49d7eSJed Brown     ///
5429df49d7eSJed Brown     /// # arguments
5439df49d7eSJed Brown     ///
5449df49d7eSJed Brown     /// * `topo`    - Topology of element, e.g. hypercube, simplex, ect
5459df49d7eSJed Brown     /// * `ncomp`   - Number of field components (1 for scalar fields)
5469df49d7eSJed Brown     /// * `nnodes`  - Total number of nodes
5479df49d7eSJed Brown     /// * `nqpts`   - Total number of quadrature points
5489df49d7eSJed Brown     /// * `interp`  - Row-major `(nqpts * nnodes)` matrix expressing the values of
5499df49d7eSJed Brown     ///                 nodal basis functions at quadrature points
5509df49d7eSJed Brown     /// * `grad`    - Row-major `(nqpts * dim * nnodes)` matrix expressing
5519df49d7eSJed Brown     ///                 derivatives of nodal basis functions at quadrature points
5529df49d7eSJed Brown     /// * `qref`    - Array of length `nqpts` holding the locations of quadrature
5539df49d7eSJed Brown     ///                 points on the reference element `[-1, 1]`
5549df49d7eSJed Brown     /// * `qweight` - Array of length `nqpts` holding the quadrature weights on
5559df49d7eSJed Brown     ///                 the reference element
5569df49d7eSJed Brown     ///
5579df49d7eSJed Brown     /// ```
5589df49d7eSJed Brown     /// # use libceed::prelude::*;
559*c68be7a2SJeremy L Thompson     /// # fn main() -> Result<(), libceed::CeedError> {
5609df49d7eSJed Brown     /// # let ceed = libceed::Ceed::default_init();
5619df49d7eSJed Brown     /// let interp = [
5629df49d7eSJed Brown     ///     0.12000000,
5639df49d7eSJed Brown     ///     0.48000000,
5649df49d7eSJed Brown     ///     -0.12000000,
5659df49d7eSJed Brown     ///     0.48000000,
5669df49d7eSJed Brown     ///     0.16000000,
5679df49d7eSJed Brown     ///     -0.12000000,
5689df49d7eSJed Brown     ///     -0.12000000,
5699df49d7eSJed Brown     ///     0.48000000,
5709df49d7eSJed Brown     ///     0.12000000,
5719df49d7eSJed Brown     ///     0.16000000,
5729df49d7eSJed Brown     ///     0.48000000,
5739df49d7eSJed Brown     ///     -0.12000000,
5749df49d7eSJed Brown     ///     -0.11111111,
5759df49d7eSJed Brown     ///     0.44444444,
5769df49d7eSJed Brown     ///     -0.11111111,
5779df49d7eSJed Brown     ///     0.44444444,
5789df49d7eSJed Brown     ///     0.44444444,
5799df49d7eSJed Brown     ///     -0.11111111,
5809df49d7eSJed Brown     ///     -0.12000000,
5819df49d7eSJed Brown     ///     0.16000000,
5829df49d7eSJed Brown     ///     -0.12000000,
5839df49d7eSJed Brown     ///     0.48000000,
5849df49d7eSJed Brown     ///     0.48000000,
5859df49d7eSJed Brown     ///     0.12000000,
5869df49d7eSJed Brown     /// ];
5879df49d7eSJed Brown     /// let grad = [
5889df49d7eSJed Brown     ///     -1.40000000,
5899df49d7eSJed Brown     ///     1.60000000,
5909df49d7eSJed Brown     ///     -0.20000000,
5919df49d7eSJed Brown     ///     -0.80000000,
5929df49d7eSJed Brown     ///     0.80000000,
5939df49d7eSJed Brown     ///     0.00000000,
5949df49d7eSJed Brown     ///     0.20000000,
5959df49d7eSJed Brown     ///     -1.60000000,
5969df49d7eSJed Brown     ///     1.40000000,
5979df49d7eSJed Brown     ///     -0.80000000,
5989df49d7eSJed Brown     ///     0.80000000,
5999df49d7eSJed Brown     ///     0.00000000,
6009df49d7eSJed Brown     ///     -0.33333333,
6019df49d7eSJed Brown     ///     0.00000000,
6029df49d7eSJed Brown     ///     0.33333333,
6039df49d7eSJed Brown     ///     -1.33333333,
6049df49d7eSJed Brown     ///     1.33333333,
6059df49d7eSJed Brown     ///     0.00000000,
6069df49d7eSJed Brown     ///     0.20000000,
6079df49d7eSJed Brown     ///     0.00000000,
6089df49d7eSJed Brown     ///     -0.20000000,
6099df49d7eSJed Brown     ///     -2.40000000,
6109df49d7eSJed Brown     ///     2.40000000,
6119df49d7eSJed Brown     ///     0.00000000,
6129df49d7eSJed Brown     ///     -1.40000000,
6139df49d7eSJed Brown     ///     -0.80000000,
6149df49d7eSJed Brown     ///     0.00000000,
6159df49d7eSJed Brown     ///     1.60000000,
6169df49d7eSJed Brown     ///     0.80000000,
6179df49d7eSJed Brown     ///     -0.20000000,
6189df49d7eSJed Brown     ///     0.20000000,
6199df49d7eSJed Brown     ///     -2.40000000,
6209df49d7eSJed Brown     ///     0.00000000,
6219df49d7eSJed Brown     ///     0.00000000,
6229df49d7eSJed Brown     ///     2.40000000,
6239df49d7eSJed Brown     ///     -0.20000000,
6249df49d7eSJed Brown     ///     -0.33333333,
6259df49d7eSJed Brown     ///     -1.33333333,
6269df49d7eSJed Brown     ///     0.00000000,
6279df49d7eSJed Brown     ///     0.00000000,
6289df49d7eSJed Brown     ///     1.33333333,
6299df49d7eSJed Brown     ///     0.33333333,
6309df49d7eSJed Brown     ///     0.20000000,
6319df49d7eSJed Brown     ///     -0.80000000,
6329df49d7eSJed Brown     ///     0.00000000,
6339df49d7eSJed Brown     ///     -1.60000000,
6349df49d7eSJed Brown     ///     0.80000000,
6359df49d7eSJed Brown     ///     1.40000000,
6369df49d7eSJed Brown     /// ];
6379df49d7eSJed Brown     /// let qref = [
6389df49d7eSJed Brown     ///     0.20000000, 0.60000000, 0.33333333, 0.20000000, 0.20000000, 0.20000000, 0.33333333,
6399df49d7eSJed Brown     ///     0.60000000,
6409df49d7eSJed Brown     /// ];
6419df49d7eSJed Brown     /// let qweight = [0.26041667, 0.26041667, -0.28125000, 0.26041667];
642*c68be7a2SJeremy L Thompson     /// let b = ceed.basis_H1(
6439df49d7eSJed Brown     ///     ElemTopology::Triangle,
6449df49d7eSJed Brown     ///     1,
6459df49d7eSJed Brown     ///     6,
6469df49d7eSJed Brown     ///     4,
6479df49d7eSJed Brown     ///     &interp,
6489df49d7eSJed Brown     ///     &grad,
6499df49d7eSJed Brown     ///     &qref,
6509df49d7eSJed Brown     ///     &qweight,
651*c68be7a2SJeremy L Thompson     /// )?;
652*c68be7a2SJeremy L Thompson     /// # Ok(())
653*c68be7a2SJeremy L Thompson     /// # }
6549df49d7eSJed Brown     /// ```
6559df49d7eSJed Brown     pub fn basis_H1(
6569df49d7eSJed Brown         &self,
6579df49d7eSJed Brown         topo: ElemTopology,
6589df49d7eSJed Brown         ncomp: usize,
6599df49d7eSJed Brown         nnodes: usize,
6609df49d7eSJed Brown         nqpts: usize,
66180a9ef05SNatalie Beams         interp: &[crate::Scalar],
66280a9ef05SNatalie Beams         grad: &[crate::Scalar],
66380a9ef05SNatalie Beams         qref: &[crate::Scalar],
66480a9ef05SNatalie Beams         qweight: &[crate::Scalar],
6659df49d7eSJed Brown     ) -> Result<Basis> {
6669df49d7eSJed Brown         Basis::create_H1(
6679df49d7eSJed Brown             self, topo, ncomp, nnodes, nqpts, interp, grad, qref, qweight,
6689df49d7eSJed Brown         )
6699df49d7eSJed Brown     }
6709df49d7eSJed Brown 
6719df49d7eSJed Brown     /// Returns a CeedQFunction for evaluating interior (volumetric) terms
6729df49d7eSJed Brown     ///
6739df49d7eSJed Brown     /// # arguments
6749df49d7eSJed Brown     ///
6759df49d7eSJed Brown     /// * `vlength` - Vector length. Caller must ensure that number of
6769df49d7eSJed Brown     ///                 quadrature points is a multiple of vlength.
6779df49d7eSJed Brown     /// * `f`       - Boxed closure to evaluate action at quadrature points.
6789df49d7eSJed Brown     ///
6799df49d7eSJed Brown     /// ```
6809df49d7eSJed Brown     /// # use libceed::prelude::*;
681*c68be7a2SJeremy L Thompson     /// # fn main() -> Result<(), libceed::CeedError> {
6829df49d7eSJed Brown     /// # let ceed = libceed::Ceed::default_init();
6839df49d7eSJed Brown     /// let mut user_f = |[u, weights, ..]: QFunctionInputs, [v, ..]: QFunctionOutputs| {
6849df49d7eSJed Brown     ///     // Iterate over quadrature points
6859df49d7eSJed Brown     ///     v.iter_mut()
6869df49d7eSJed Brown     ///         .zip(u.iter().zip(weights.iter()))
6879df49d7eSJed Brown     ///         .for_each(|(v, (u, w))| *v = u * w);
6889df49d7eSJed Brown     ///
6899df49d7eSJed Brown     ///     // Return clean error code
6909df49d7eSJed Brown     ///     0
6919df49d7eSJed Brown     /// };
6929df49d7eSJed Brown     ///
693*c68be7a2SJeremy L Thompson     /// let qf = ceed.q_function_interior(1, Box::new(user_f))?;
694*c68be7a2SJeremy L Thompson     /// # Ok(())
695*c68be7a2SJeremy L Thompson     /// # }
6969df49d7eSJed Brown     /// ```
6979df49d7eSJed Brown     pub fn q_function_interior(
6989df49d7eSJed Brown         &self,
6999df49d7eSJed Brown         vlength: usize,
7009df49d7eSJed Brown         f: Box<qfunction::QFunctionUserClosure>,
7019df49d7eSJed Brown     ) -> Result<QFunction> {
7029df49d7eSJed Brown         QFunction::create(self, vlength, f)
7039df49d7eSJed Brown     }
7049df49d7eSJed Brown 
7059df49d7eSJed Brown     /// Returns a CeedQFunction for evaluating interior (volumetric) terms
7069df49d7eSJed Brown     /// created by name
7079df49d7eSJed Brown     ///
7089df49d7eSJed Brown     /// ```
7099df49d7eSJed Brown     /// # use libceed::prelude::*;
710*c68be7a2SJeremy L Thompson     /// # fn main() -> Result<(), libceed::CeedError> {
7119df49d7eSJed Brown     /// # let ceed = libceed::Ceed::default_init();
712*c68be7a2SJeremy L Thompson     /// let qf = ceed.q_function_interior_by_name("Mass1DBuild")?;
713*c68be7a2SJeremy L Thompson     /// # Ok(())
714*c68be7a2SJeremy L Thompson     /// # }
7159df49d7eSJed Brown     /// ```
7169df49d7eSJed Brown     pub fn q_function_interior_by_name(&self, name: &str) -> Result<QFunctionByName> {
7179df49d7eSJed Brown         QFunctionByName::create(self, name)
7189df49d7eSJed Brown     }
7199df49d7eSJed Brown 
7209df49d7eSJed Brown     /// Returns a Operator and associate a QFunction. A Basis and
7219df49d7eSJed Brown     /// ElemRestriction can be   associated with QFunction fields with
7229df49d7eSJed Brown     /// set_field().
7239df49d7eSJed Brown     ///
7249df49d7eSJed Brown     /// * `qf`   - QFunction defining the action of the operator at quadrature
7259df49d7eSJed Brown     ///              points
7269df49d7eSJed Brown     /// * `dqf`  - QFunction defining the action of the Jacobian of the qf (or
7279df49d7eSJed Brown     ///              qfunction_none)
7289df49d7eSJed Brown     /// * `dqfT` - QFunction defining the action of the transpose of the
7299df49d7eSJed Brown     ///              Jacobian of the qf (or qfunction_none)
7309df49d7eSJed Brown     ///
7319df49d7eSJed Brown     /// ```
7329df49d7eSJed Brown     /// # use libceed::prelude::*;
733*c68be7a2SJeremy L Thompson     /// # fn main() -> Result<(), libceed::CeedError> {
7349df49d7eSJed Brown     /// # let ceed = libceed::Ceed::default_init();
735*c68be7a2SJeremy L Thompson     /// let qf = ceed.q_function_interior_by_name("Mass1DBuild")?;
736*c68be7a2SJeremy L Thompson     /// let op = ceed.operator(&qf, QFunctionOpt::None, QFunctionOpt::None)?;
737*c68be7a2SJeremy L Thompson     /// # Ok(())
738*c68be7a2SJeremy L Thompson     /// # }
7399df49d7eSJed Brown     /// ```
7409df49d7eSJed Brown     pub fn operator<'b>(
7419df49d7eSJed Brown         &self,
7429df49d7eSJed Brown         qf: impl Into<QFunctionOpt<'b>>,
7439df49d7eSJed Brown         dqf: impl Into<QFunctionOpt<'b>>,
7449df49d7eSJed Brown         dqfT: impl Into<QFunctionOpt<'b>>,
7459df49d7eSJed Brown     ) -> Result<Operator> {
7469df49d7eSJed Brown         Operator::create(self, qf, dqf, dqfT)
7479df49d7eSJed Brown     }
7489df49d7eSJed Brown 
7499df49d7eSJed Brown     /// Returns an Operator that composes the action of several Operators
7509df49d7eSJed Brown     ///
7519df49d7eSJed Brown     /// ```
7529df49d7eSJed Brown     /// # use libceed::prelude::*;
753*c68be7a2SJeremy L Thompson     /// # fn main() -> Result<(), libceed::CeedError> {
7549df49d7eSJed Brown     /// # let ceed = libceed::Ceed::default_init();
755*c68be7a2SJeremy L Thompson     /// let op = ceed.composite_operator()?;
756*c68be7a2SJeremy L Thompson     /// # Ok(())
757*c68be7a2SJeremy L Thompson     /// # }
7589df49d7eSJed Brown     /// ```
7599df49d7eSJed Brown     pub fn composite_operator(&self) -> Result<CompositeOperator> {
7609df49d7eSJed Brown         CompositeOperator::create(self)
7619df49d7eSJed Brown     }
7629df49d7eSJed Brown }
7639df49d7eSJed Brown 
7649df49d7eSJed Brown // -----------------------------------------------------------------------------
7659df49d7eSJed Brown // Tests
7669df49d7eSJed Brown // -----------------------------------------------------------------------------
7679df49d7eSJed Brown #[cfg(test)]
7689df49d7eSJed Brown mod tests {
7699df49d7eSJed Brown     use super::*;
7709df49d7eSJed Brown 
7719df49d7eSJed Brown     fn ceed_t501() -> Result<i32> {
7729df49d7eSJed Brown         let resource = "/cpu/self/ref/blocked";
7739df49d7eSJed Brown         let ceed = Ceed::init(resource);
7749df49d7eSJed Brown         let nelem = 4;
7759df49d7eSJed Brown         let p = 3;
7769df49d7eSJed Brown         let q = 4;
7779df49d7eSJed Brown         let ndofs = p * nelem - nelem + 1;
7789df49d7eSJed Brown 
7799df49d7eSJed Brown         // Vectors
7809df49d7eSJed Brown         let x = ceed.vector_from_slice(&[-1., -0.5, 0.0, 0.5, 1.0])?;
7819df49d7eSJed Brown         let mut qdata = ceed.vector(nelem * q)?;
7829df49d7eSJed Brown         qdata.set_value(0.0)?;
7839df49d7eSJed Brown         let mut u = ceed.vector(ndofs)?;
7849df49d7eSJed Brown         u.set_value(1.0)?;
7859df49d7eSJed Brown         let mut v = ceed.vector(ndofs)?;
7869df49d7eSJed Brown         v.set_value(0.0)?;
7879df49d7eSJed Brown 
7889df49d7eSJed Brown         // Restrictions
7899df49d7eSJed Brown         let mut indx: Vec<i32> = vec![0; 2 * nelem];
7909df49d7eSJed Brown         for i in 0..nelem {
7919df49d7eSJed Brown             indx[2 * i + 0] = i as i32;
7929df49d7eSJed Brown             indx[2 * i + 1] = (i + 1) as i32;
7939df49d7eSJed Brown         }
7949df49d7eSJed Brown         let rx = ceed.elem_restriction(nelem, 2, 1, 1, nelem + 1, MemType::Host, &indx)?;
7959df49d7eSJed Brown         let mut indu: Vec<i32> = vec![0; p * nelem];
7969df49d7eSJed Brown         for i in 0..nelem {
7979df49d7eSJed Brown             indu[p * i + 0] = i as i32;
7989df49d7eSJed Brown             indu[p * i + 1] = (i + 1) as i32;
7999df49d7eSJed Brown             indu[p * i + 2] = (i + 2) as i32;
8009df49d7eSJed Brown         }
8019df49d7eSJed Brown         let ru = ceed.elem_restriction(nelem, 3, 1, 1, ndofs, MemType::Host, &indu)?;
8029df49d7eSJed Brown         let strides: [i32; 3] = [1, q as i32, q as i32];
8039df49d7eSJed Brown         let rq = ceed.strided_elem_restriction(nelem, q, 1, q * nelem, strides)?;
8049df49d7eSJed Brown 
8059df49d7eSJed Brown         // Bases
8069df49d7eSJed Brown         let bx = ceed.basis_tensor_H1_Lagrange(1, 1, 2, q, QuadMode::Gauss)?;
8079df49d7eSJed Brown         let bu = ceed.basis_tensor_H1_Lagrange(1, 1, p, q, QuadMode::Gauss)?;
8089df49d7eSJed Brown 
8099df49d7eSJed Brown         // Build quadrature data
8109df49d7eSJed Brown         let qf_build = ceed.q_function_interior_by_name("Mass1DBuild")?;
8119df49d7eSJed Brown         ceed.operator(&qf_build, QFunctionOpt::None, QFunctionOpt::None)?
8129df49d7eSJed Brown             .field("dx", &rx, &bx, VectorOpt::Active)?
8139df49d7eSJed Brown             .field("weights", ElemRestrictionOpt::None, &bx, VectorOpt::None)?
8149df49d7eSJed Brown             .field("qdata", &rq, BasisOpt::Collocated, VectorOpt::Active)?
8159df49d7eSJed Brown             .apply(&x, &mut qdata)?;
8169df49d7eSJed Brown 
8179df49d7eSJed Brown         // Mass operator
8189df49d7eSJed Brown         let qf_mass = ceed.q_function_interior_by_name("MassApply")?;
8199df49d7eSJed Brown         let op_mass = ceed
8209df49d7eSJed Brown             .operator(&qf_mass, QFunctionOpt::None, QFunctionOpt::None)?
8219df49d7eSJed Brown             .field("u", &ru, &bu, VectorOpt::Active)?
8229df49d7eSJed Brown             .field("qdata", &rq, BasisOpt::Collocated, &qdata)?
8239df49d7eSJed Brown             .field("v", &ru, &bu, VectorOpt::Active)?;
8249df49d7eSJed Brown 
8259df49d7eSJed Brown         v.set_value(0.0)?;
8269df49d7eSJed Brown         op_mass.apply(&u, &mut v)?;
8279df49d7eSJed Brown 
8289df49d7eSJed Brown         // Check
82980a9ef05SNatalie Beams         let sum: Scalar = v.view().iter().sum();
8309df49d7eSJed Brown         assert!(
8319df49d7eSJed Brown             (sum - 2.0).abs() < 1e-15,
8329df49d7eSJed Brown             "Incorrect interval length computed"
8339df49d7eSJed Brown         );
8349df49d7eSJed Brown         Ok(0)
8359df49d7eSJed Brown     }
8369df49d7eSJed Brown 
8379df49d7eSJed Brown     #[test]
8389df49d7eSJed Brown     fn test_ceed_t501() {
8399df49d7eSJed Brown         assert!(ceed_t501().is_ok());
8409df49d7eSJed Brown     }
8419df49d7eSJed Brown }
8429df49d7eSJed Brown 
8439df49d7eSJed Brown // -----------------------------------------------------------------------------
844