xref: /libCEED/rust/libceed/src/lib.rs (revision 4b61a5a04425d85ec42b770771dd05c355ddbfd8)
13d8e8822SJeremy L Thompson // Copyright (c) 2017-2022, Lawrence Livermore National Security, LLC and other CEED contributors.
23d8e8822SJeremy L Thompson // All Rights Reserved. See the top-level LICENSE and NOTICE files for details.
39df49d7eSJed Brown //
43d8e8822SJeremy L Thompson // SPDX-License-Identifier: BSD-2-Clause
59df49d7eSJed Brown //
63d8e8822SJeremy L Thompson // This file is part of CEED:  http://github.com/ceed
79df49d7eSJed Brown 
8a1cbad85SJed Brown // Fenced `rust` code blocks included from README.md are executed as part of doctests.
9a1cbad85SJed Brown #![doc = include_str!("../README.md")]
109df49d7eSJed Brown // -----------------------------------------------------------------------------
119df49d7eSJed Brown // Exceptions
129df49d7eSJed Brown // -----------------------------------------------------------------------------
139df49d7eSJed Brown #![allow(non_snake_case)]
149df49d7eSJed Brown 
159df49d7eSJed Brown // -----------------------------------------------------------------------------
169df49d7eSJed Brown // Crate prelude
179df49d7eSJed Brown // -----------------------------------------------------------------------------
189df49d7eSJed Brown use crate::prelude::*;
199df49d7eSJed Brown use std::sync::Once;
209df49d7eSJed Brown 
219df49d7eSJed Brown pub mod prelude {
229df49d7eSJed Brown     pub use crate::{
239df49d7eSJed Brown         basis::{self, Basis, BasisOpt},
249df49d7eSJed Brown         elem_restriction::{self, ElemRestriction, ElemRestrictionOpt},
2508778c6fSJeremy L Thompson         operator::{self, CompositeOperator, Operator, OperatorField},
269df49d7eSJed Brown         qfunction::{
2708778c6fSJeremy L Thompson             self, QFunction, QFunctionByName, QFunctionField, QFunctionInputs, QFunctionOpt,
2808778c6fSJeremy L Thompson             QFunctionOutputs,
299df49d7eSJed Brown         },
30486868d3SJeremy L Thompson         vector::{self, Vector, VectorOpt, VectorSliceWrapper},
312ba8e59cSJeremy L Thompson         ElemTopology, EvalMode, MemType, NormType, QuadMode, Scalar, TransposeMode,
3280a9ef05SNatalie Beams         CEED_STRIDES_BACKEND, EPSILON, MAX_QFUNCTION_FIELDS,
339df49d7eSJed Brown     };
34630ad4c9Sjeremylt     pub(crate) use libceed_sys::bind_ceed;
359df49d7eSJed Brown     pub(crate) use std::convert::TryFrom;
3608778c6fSJeremy L Thompson     pub(crate) use std::ffi::{CStr, CString};
379df49d7eSJed Brown     pub(crate) use std::fmt;
381142270cSJeremy L Thompson     pub(crate) use std::marker::PhantomData;
399df49d7eSJed Brown }
409df49d7eSJed Brown 
419df49d7eSJed Brown // -----------------------------------------------------------------------------
429df49d7eSJed Brown // Modules
439df49d7eSJed Brown // -----------------------------------------------------------------------------
449df49d7eSJed Brown pub mod basis;
459df49d7eSJed Brown pub mod elem_restriction;
469df49d7eSJed Brown pub mod operator;
479df49d7eSJed Brown pub mod qfunction;
489df49d7eSJed Brown pub mod vector;
499df49d7eSJed Brown 
509df49d7eSJed Brown // -----------------------------------------------------------------------------
5180a9ef05SNatalie Beams // Typedef for scalar
5280a9ef05SNatalie Beams // -----------------------------------------------------------------------------
5380a9ef05SNatalie Beams pub type Scalar = bind_ceed::CeedScalar;
5480a9ef05SNatalie Beams 
5580a9ef05SNatalie Beams // -----------------------------------------------------------------------------
569df49d7eSJed Brown // Constants for library
579df49d7eSJed Brown // -----------------------------------------------------------------------------
58257916f8SJed Brown const MAX_BUFFER_LENGTH: usize = 4096;
599df49d7eSJed Brown pub const MAX_QFUNCTION_FIELDS: usize = 16;
609df49d7eSJed Brown pub const CEED_STRIDES_BACKEND: [i32; 3] = [0; 3];
6180a9ef05SNatalie Beams pub const EPSILON: crate::Scalar = bind_ceed::CEED_EPSILON as crate::Scalar;
629df49d7eSJed Brown 
639df49d7eSJed Brown // -----------------------------------------------------------------------------
649df49d7eSJed Brown // Enums for libCEED
659df49d7eSJed Brown // -----------------------------------------------------------------------------
669df49d7eSJed Brown #[derive(Clone, Copy, PartialEq, Eq)]
679df49d7eSJed Brown /// Many Ceed interfaces take or return pointers to memory.  This enum is used to
689df49d7eSJed Brown /// specify where the memory being provided or requested must reside.
699df49d7eSJed Brown pub enum MemType {
709df49d7eSJed Brown     Host = bind_ceed::CeedMemType_CEED_MEM_HOST as isize,
719df49d7eSJed Brown     Device = bind_ceed::CeedMemType_CEED_MEM_DEVICE as isize,
729df49d7eSJed Brown }
739df49d7eSJed Brown 
749df49d7eSJed Brown #[derive(Clone, Copy, PartialEq, Eq)]
759df49d7eSJed Brown // OwnPointer will not be used by user but is included for internal use
769df49d7eSJed Brown #[allow(dead_code)]
779df49d7eSJed Brown /// Conveys ownership status of arrays passed to Ceed interfaces.
789df49d7eSJed Brown pub(crate) enum CopyMode {
799df49d7eSJed Brown     CopyValues = bind_ceed::CeedCopyMode_CEED_COPY_VALUES as isize,
809df49d7eSJed Brown     UsePointer = bind_ceed::CeedCopyMode_CEED_USE_POINTER as isize,
819df49d7eSJed Brown     OwnPointer = bind_ceed::CeedCopyMode_CEED_OWN_POINTER as isize,
829df49d7eSJed Brown }
839df49d7eSJed Brown 
849df49d7eSJed Brown #[derive(Clone, Copy, PartialEq, Eq)]
859df49d7eSJed Brown /// Denotes type of vector norm to be computed
869df49d7eSJed Brown pub enum NormType {
879df49d7eSJed Brown     One = bind_ceed::CeedNormType_CEED_NORM_1 as isize,
889df49d7eSJed Brown     Two = bind_ceed::CeedNormType_CEED_NORM_2 as isize,
899df49d7eSJed Brown     Max = bind_ceed::CeedNormType_CEED_NORM_MAX as isize,
909df49d7eSJed Brown }
919df49d7eSJed Brown 
929df49d7eSJed Brown #[derive(Clone, Copy, PartialEq, Eq)]
939df49d7eSJed Brown /// Denotes whether a linear transformation or its transpose should be applied
949df49d7eSJed Brown pub enum TransposeMode {
959df49d7eSJed Brown     NoTranspose = bind_ceed::CeedTransposeMode_CEED_NOTRANSPOSE as isize,
969df49d7eSJed Brown     Transpose = bind_ceed::CeedTransposeMode_CEED_TRANSPOSE as isize,
979df49d7eSJed Brown }
989df49d7eSJed Brown 
999df49d7eSJed Brown #[derive(Clone, Copy, PartialEq, Eq)]
1009df49d7eSJed Brown /// Type of quadrature; also used for location of nodes
1019df49d7eSJed Brown pub enum QuadMode {
1029df49d7eSJed Brown     Gauss = bind_ceed::CeedQuadMode_CEED_GAUSS as isize,
1039df49d7eSJed Brown     GaussLobatto = bind_ceed::CeedQuadMode_CEED_GAUSS_LOBATTO as isize,
1049df49d7eSJed Brown }
1059df49d7eSJed Brown 
1069df49d7eSJed Brown #[derive(Clone, Copy, PartialEq, Eq)]
1079df49d7eSJed Brown /// Type of basis shape to create non-tensor H1 element basis
1089df49d7eSJed Brown pub enum ElemTopology {
10950c301a5SRezgar Shakeri     Line = bind_ceed::CeedElemTopology_CEED_TOPOLOGY_LINE as isize,
11050c301a5SRezgar Shakeri     Triangle = bind_ceed::CeedElemTopology_CEED_TOPOLOGY_TRIANGLE as isize,
11150c301a5SRezgar Shakeri     Quad = bind_ceed::CeedElemTopology_CEED_TOPOLOGY_QUAD as isize,
11250c301a5SRezgar Shakeri     Tet = bind_ceed::CeedElemTopology_CEED_TOPOLOGY_TET as isize,
11350c301a5SRezgar Shakeri     Pyramid = bind_ceed::CeedElemTopology_CEED_TOPOLOGY_PYRAMID as isize,
11450c301a5SRezgar Shakeri     Prism = bind_ceed::CeedElemTopology_CEED_TOPOLOGY_PRISM as isize,
11550c301a5SRezgar Shakeri     Hex = bind_ceed::CeedElemTopology_CEED_TOPOLOGY_HEX as isize,
1169df49d7eSJed Brown }
1179df49d7eSJed Brown 
118c68be7a2SJeremy L Thompson #[derive(Clone, Copy, Debug, PartialEq, Eq)]
1199df49d7eSJed Brown /// Basis evaluation mode
1209df49d7eSJed Brown pub enum EvalMode {
1219df49d7eSJed Brown     None = bind_ceed::CeedEvalMode_CEED_EVAL_NONE as isize,
1229df49d7eSJed Brown     Interp = bind_ceed::CeedEvalMode_CEED_EVAL_INTERP as isize,
1239df49d7eSJed Brown     Grad = bind_ceed::CeedEvalMode_CEED_EVAL_GRAD as isize,
1249df49d7eSJed Brown     Div = bind_ceed::CeedEvalMode_CEED_EVAL_DIV as isize,
1259df49d7eSJed Brown     Curl = bind_ceed::CeedEvalMode_CEED_EVAL_CURL as isize,
1269df49d7eSJed Brown     Weight = bind_ceed::CeedEvalMode_CEED_EVAL_WEIGHT as isize,
1279df49d7eSJed Brown }
128c68be7a2SJeremy L Thompson impl EvalMode {
129c68be7a2SJeremy L Thompson     pub(crate) fn from_u32(value: u32) -> EvalMode {
130c68be7a2SJeremy L Thompson         match value {
131c68be7a2SJeremy L Thompson             bind_ceed::CeedEvalMode_CEED_EVAL_NONE => EvalMode::None,
132c68be7a2SJeremy L Thompson             bind_ceed::CeedEvalMode_CEED_EVAL_INTERP => EvalMode::Interp,
133c68be7a2SJeremy L Thompson             bind_ceed::CeedEvalMode_CEED_EVAL_GRAD => EvalMode::Grad,
134c68be7a2SJeremy L Thompson             bind_ceed::CeedEvalMode_CEED_EVAL_DIV => EvalMode::Div,
135c68be7a2SJeremy L Thompson             bind_ceed::CeedEvalMode_CEED_EVAL_CURL => EvalMode::Curl,
136c68be7a2SJeremy L Thompson             bind_ceed::CeedEvalMode_CEED_EVAL_WEIGHT => EvalMode::Weight,
137c68be7a2SJeremy L Thompson             _ => panic!("Unknown value: {}", value),
138c68be7a2SJeremy L Thompson         }
139c68be7a2SJeremy L Thompson     }
140c68be7a2SJeremy L Thompson }
1419df49d7eSJed Brown 
1429df49d7eSJed Brown // -----------------------------------------------------------------------------
1439df49d7eSJed Brown // Ceed error
1449df49d7eSJed Brown // -----------------------------------------------------------------------------
1452ba8e59cSJeremy L Thompson pub type Result<T> = std::result::Result<T, Error>;
1469df49d7eSJed Brown 
1477ed177dbSJed Brown /// libCEED error messages - returning an Error without painc!ing indicates
1487ed177dbSJed Brown ///   the function call failed but the data is not corrupted
1499df49d7eSJed Brown #[derive(Debug)]
1502ba8e59cSJeremy L Thompson pub struct Error {
1519df49d7eSJed Brown     pub message: String,
1529df49d7eSJed Brown }
1539df49d7eSJed Brown 
1542ba8e59cSJeremy L Thompson impl fmt::Display for Error {
1559df49d7eSJed Brown     fn fmt(&self, f: &mut fmt::Formatter) -> fmt::Result {
1569df49d7eSJed Brown         write!(f, "{}", self.message)
1579df49d7eSJed Brown     }
1589df49d7eSJed Brown }
1599df49d7eSJed Brown 
1609df49d7eSJed Brown // -----------------------------------------------------------------------------
1611142270cSJeremy L Thompson // Internal error checker
1621142270cSJeremy L Thompson // -----------------------------------------------------------------------------
1631142270cSJeremy L Thompson #[doc(hidden)]
1641142270cSJeremy L Thompson pub(crate) fn check_error(ceed_ptr: bind_ceed::Ceed, ierr: i32) -> Result<i32> {
1651142270cSJeremy L Thompson     // Return early if code is clean
1661142270cSJeremy L Thompson     if ierr == bind_ceed::CeedErrorType_CEED_ERROR_SUCCESS {
1671142270cSJeremy L Thompson         return Ok(ierr);
1681142270cSJeremy L Thompson     }
1691142270cSJeremy L Thompson     // Retrieve error message
1701142270cSJeremy L Thompson     let mut ptr: *const std::os::raw::c_char = std::ptr::null_mut();
1711142270cSJeremy L Thompson     let c_str = unsafe {
1721142270cSJeremy L Thompson         bind_ceed::CeedGetErrorMessage(ceed_ptr, &mut ptr);
1731142270cSJeremy L Thompson         std::ffi::CStr::from_ptr(ptr)
1741142270cSJeremy L Thompson     };
1751142270cSJeremy L Thompson     let message = c_str.to_string_lossy().to_string();
1761142270cSJeremy L Thompson     // Panic if negative code, otherwise return error
1771142270cSJeremy L Thompson     if ierr < bind_ceed::CeedErrorType_CEED_ERROR_SUCCESS {
1781142270cSJeremy L Thompson         panic!("{}", message);
1791142270cSJeremy L Thompson     }
1802ba8e59cSJeremy L Thompson     Err(Error { message })
1811142270cSJeremy L Thompson }
1821142270cSJeremy L Thompson 
1831142270cSJeremy L Thompson // -----------------------------------------------------------------------------
1849df49d7eSJed Brown // Ceed error handler
1859df49d7eSJed Brown // -----------------------------------------------------------------------------
1862ba8e59cSJeremy L Thompson pub enum ErrorHandler {
1879df49d7eSJed Brown     ErrorAbort,
1889df49d7eSJed Brown     ErrorExit,
1899df49d7eSJed Brown     ErrorReturn,
1909df49d7eSJed Brown     ErrorStore,
1919df49d7eSJed Brown }
1929df49d7eSJed Brown 
1939df49d7eSJed Brown // -----------------------------------------------------------------------------
1949df49d7eSJed Brown // Ceed context wrapper
1959df49d7eSJed Brown // -----------------------------------------------------------------------------
1969df49d7eSJed Brown /// A Ceed is a library context representing control of a logical hardware
1979df49d7eSJed Brown /// resource.
1989df49d7eSJed Brown #[derive(Debug)]
1999df49d7eSJed Brown pub struct Ceed {
2009df49d7eSJed Brown     ptr: bind_ceed::Ceed,
2019df49d7eSJed Brown }
2029df49d7eSJed Brown 
2039df49d7eSJed Brown // -----------------------------------------------------------------------------
2049df49d7eSJed Brown // Destructor
2059df49d7eSJed Brown // -----------------------------------------------------------------------------
2069df49d7eSJed Brown impl Drop for Ceed {
2079df49d7eSJed Brown     fn drop(&mut self) {
2089df49d7eSJed Brown         unsafe {
2099df49d7eSJed Brown             bind_ceed::CeedDestroy(&mut self.ptr);
2109df49d7eSJed Brown         }
2119df49d7eSJed Brown     }
2129df49d7eSJed Brown }
2139df49d7eSJed Brown 
2149df49d7eSJed Brown // -----------------------------------------------------------------------------
21559189cfaSJeremy L Thompson // Cloning
21659189cfaSJeremy L Thompson // -----------------------------------------------------------------------------
21759189cfaSJeremy L Thompson impl Clone for Ceed {
21859189cfaSJeremy L Thompson     /// Perform a shallow clone of a Ceed context
21959189cfaSJeremy L Thompson     ///
22059189cfaSJeremy L Thompson     /// ```
22159189cfaSJeremy L Thompson     /// let ceed = libceed::Ceed::init("/cpu/self/ref/serial");
22259189cfaSJeremy L Thompson     /// let ceed_clone = ceed.clone();
22359189cfaSJeremy L Thompson     ///
2247ed177dbSJed Brown     /// println!(" original:{} \n clone:{}", ceed, ceed_clone);
22559189cfaSJeremy L Thompson     /// ```
22659189cfaSJeremy L Thompson     fn clone(&self) -> Self {
22759189cfaSJeremy L Thompson         let mut ptr_clone = std::ptr::null_mut();
22859189cfaSJeremy L Thompson         let ierr = unsafe { bind_ceed::CeedReferenceCopy(self.ptr, &mut ptr_clone) };
22959189cfaSJeremy L Thompson         self.check_error(ierr).expect("failed to clone Ceed");
23059189cfaSJeremy L Thompson         Self { ptr: ptr_clone }
23159189cfaSJeremy L Thompson     }
23259189cfaSJeremy L Thompson }
23359189cfaSJeremy L Thompson 
23459189cfaSJeremy L Thompson // -----------------------------------------------------------------------------
2359df49d7eSJed Brown // Display
2369df49d7eSJed Brown // -----------------------------------------------------------------------------
2379df49d7eSJed Brown impl fmt::Display for Ceed {
2389df49d7eSJed Brown     /// View a Ceed
2399df49d7eSJed Brown     ///
2409df49d7eSJed Brown     /// ```
2419df49d7eSJed Brown     /// let ceed = libceed::Ceed::init("/cpu/self/ref/serial");
2429df49d7eSJed Brown     /// println!("{}", ceed);
2439df49d7eSJed Brown     /// ```
2449df49d7eSJed Brown     fn fmt(&self, f: &mut fmt::Formatter) -> fmt::Result {
2459df49d7eSJed Brown         let mut ptr = std::ptr::null_mut();
2469df49d7eSJed Brown         let mut sizeloc = crate::MAX_BUFFER_LENGTH;
2479df49d7eSJed Brown         let cstring = unsafe {
2489df49d7eSJed Brown             let file = bind_ceed::open_memstream(&mut ptr, &mut sizeloc);
2499df49d7eSJed Brown             bind_ceed::CeedView(self.ptr, file);
2509df49d7eSJed Brown             bind_ceed::fclose(file);
2519df49d7eSJed Brown             CString::from_raw(ptr)
2529df49d7eSJed Brown         };
2539df49d7eSJed Brown         cstring.to_string_lossy().fmt(f)
2549df49d7eSJed Brown     }
2559df49d7eSJed Brown }
2569df49d7eSJed Brown 
2579df49d7eSJed Brown static REGISTER: Once = Once::new();
2589df49d7eSJed Brown 
2599df49d7eSJed Brown // -----------------------------------------------------------------------------
2609df49d7eSJed Brown // Object constructors
2619df49d7eSJed Brown // -----------------------------------------------------------------------------
2629df49d7eSJed Brown impl Ceed {
2637ed177dbSJed Brown     #[cfg_attr(feature = "katexit", katexit::katexit)]
2649df49d7eSJed Brown     /// Returns a Ceed context initialized with the specified resource
2659df49d7eSJed Brown     ///
2669df49d7eSJed Brown     /// # arguments
2679df49d7eSJed Brown     ///
2689df49d7eSJed Brown     /// * `resource` - Resource to use, e.g., "/cpu/self"
2699df49d7eSJed Brown     ///
2709df49d7eSJed Brown     /// ```
2719df49d7eSJed Brown     /// let ceed = libceed::Ceed::init("/cpu/self/ref/serial");
2729df49d7eSJed Brown     /// ```
2739df49d7eSJed Brown     pub fn init(resource: &str) -> Self {
2742ba8e59cSJeremy L Thompson         Ceed::init_with_error_handler(resource, ErrorHandler::ErrorStore)
2759df49d7eSJed Brown     }
2769df49d7eSJed Brown 
2779df49d7eSJed Brown     /// Returns a Ceed context initialized with the specified resource
2789df49d7eSJed Brown     ///
2799df49d7eSJed Brown     /// # arguments
2809df49d7eSJed Brown     ///
2819df49d7eSJed Brown     /// * `resource` - Resource to use, e.g., "/cpu/self"
2829df49d7eSJed Brown     ///
2839df49d7eSJed Brown     /// ```
2849df49d7eSJed Brown     /// let ceed = libceed::Ceed::init_with_error_handler(
2859df49d7eSJed Brown     ///     "/cpu/self/ref/serial",
2862ba8e59cSJeremy L Thompson     ///     libceed::ErrorHandler::ErrorAbort,
2879df49d7eSJed Brown     /// );
2889df49d7eSJed Brown     /// ```
2892ba8e59cSJeremy L Thompson     pub fn init_with_error_handler(resource: &str, handler: ErrorHandler) -> Self {
2909df49d7eSJed Brown         REGISTER.call_once(|| unsafe {
2919df49d7eSJed Brown             bind_ceed::CeedRegisterAll();
2929df49d7eSJed Brown             bind_ceed::CeedQFunctionRegisterAll();
2939df49d7eSJed Brown         });
2949df49d7eSJed Brown 
2959df49d7eSJed Brown         // Convert to C string
2969df49d7eSJed Brown         let c_resource = CString::new(resource).expect("CString::new failed");
2979df49d7eSJed Brown 
2989df49d7eSJed Brown         // Get error handler pointer
2999df49d7eSJed Brown         let eh = match handler {
3002ba8e59cSJeremy L Thompson             ErrorHandler::ErrorAbort => bind_ceed::CeedErrorAbort,
3012ba8e59cSJeremy L Thompson             ErrorHandler::ErrorExit => bind_ceed::CeedErrorExit,
3022ba8e59cSJeremy L Thompson             ErrorHandler::ErrorReturn => bind_ceed::CeedErrorReturn,
3032ba8e59cSJeremy L Thompson             ErrorHandler::ErrorStore => bind_ceed::CeedErrorStore,
3049df49d7eSJed Brown         };
3059df49d7eSJed Brown 
3069df49d7eSJed Brown         // Call to libCEED
3079df49d7eSJed Brown         let mut ptr = std::ptr::null_mut();
3089df49d7eSJed Brown         let mut ierr = unsafe { bind_ceed::CeedInit(c_resource.as_ptr() as *const i8, &mut ptr) };
3099df49d7eSJed Brown         if ierr != 0 {
3109df49d7eSJed Brown             panic!("Error initializing backend resource: {}", resource)
3119df49d7eSJed Brown         }
3129df49d7eSJed Brown         ierr = unsafe { bind_ceed::CeedSetErrorHandler(ptr, Some(eh)) };
3139df49d7eSJed Brown         let ceed = Ceed { ptr };
3149df49d7eSJed Brown         ceed.check_error(ierr).unwrap();
3159df49d7eSJed Brown         ceed
3169df49d7eSJed Brown     }
3179df49d7eSJed Brown 
3189df49d7eSJed Brown     /// Default initializer for testing
3199df49d7eSJed Brown     #[doc(hidden)]
3209df49d7eSJed Brown     pub fn default_init() -> Self {
3219df49d7eSJed Brown         // Convert to C string
3229df49d7eSJed Brown         let resource = "/cpu/self/ref/serial";
3239df49d7eSJed Brown         crate::Ceed::init(resource)
3249df49d7eSJed Brown     }
3259df49d7eSJed Brown 
3269df49d7eSJed Brown     /// Internal error checker
3279df49d7eSJed Brown     #[doc(hidden)]
3289df49d7eSJed Brown     fn check_error(&self, ierr: i32) -> Result<i32> {
3299df49d7eSJed Brown         // Return early if code is clean
3309df49d7eSJed Brown         if ierr == bind_ceed::CeedErrorType_CEED_ERROR_SUCCESS {
3319df49d7eSJed Brown             return Ok(ierr);
3329df49d7eSJed Brown         }
3339df49d7eSJed Brown         // Retrieve error message
3349df49d7eSJed Brown         let mut ptr: *const std::os::raw::c_char = std::ptr::null_mut();
3359df49d7eSJed Brown         let c_str = unsafe {
3369df49d7eSJed Brown             bind_ceed::CeedGetErrorMessage(self.ptr, &mut ptr);
3379df49d7eSJed Brown             std::ffi::CStr::from_ptr(ptr)
3389df49d7eSJed Brown         };
3399df49d7eSJed Brown         let message = c_str.to_string_lossy().to_string();
3409df49d7eSJed Brown         // Panic if negative code, otherwise return error
3419df49d7eSJed Brown         if ierr < bind_ceed::CeedErrorType_CEED_ERROR_SUCCESS {
3429df49d7eSJed Brown             panic!("{}", message);
3439df49d7eSJed Brown         }
3442ba8e59cSJeremy L Thompson         Err(Error { message })
3459df49d7eSJed Brown     }
3469df49d7eSJed Brown 
3479df49d7eSJed Brown     /// Returns full resource name for a Ceed context
3489df49d7eSJed Brown     ///
3499df49d7eSJed Brown     /// ```
3509df49d7eSJed Brown     /// let ceed = libceed::Ceed::init("/cpu/self/ref/serial");
3519df49d7eSJed Brown     /// let resource = ceed.resource();
3529df49d7eSJed Brown     ///
3539df49d7eSJed Brown     /// assert_eq!(resource, "/cpu/self/ref/serial".to_string())
3549df49d7eSJed Brown     /// ```
3559df49d7eSJed Brown     pub fn resource(&self) -> String {
3569df49d7eSJed Brown         let mut ptr: *const std::os::raw::c_char = std::ptr::null_mut();
3579df49d7eSJed Brown         let c_str = unsafe {
3589df49d7eSJed Brown             bind_ceed::CeedGetResource(self.ptr, &mut ptr);
3599df49d7eSJed Brown             std::ffi::CStr::from_ptr(ptr)
3609df49d7eSJed Brown         };
3619df49d7eSJed Brown         c_str.to_string_lossy().to_string()
3629df49d7eSJed Brown     }
3639df49d7eSJed Brown 
3647ed177dbSJed Brown     /// Returns a Vector of the specified length (does not allocate memory)
3659df49d7eSJed Brown     ///
3669df49d7eSJed Brown     /// # arguments
3679df49d7eSJed Brown     ///
3689df49d7eSJed Brown     /// * `n` - Length of vector
3699df49d7eSJed Brown     ///
3709df49d7eSJed Brown     /// ```
3719df49d7eSJed Brown     /// # use libceed::prelude::*;
3724d27c890SJeremy L Thompson     /// # fn main() -> libceed::Result<()> {
3739df49d7eSJed Brown     /// # let ceed = libceed::Ceed::default_init();
374c68be7a2SJeremy L Thompson     /// let vec = ceed.vector(10)?;
375c68be7a2SJeremy L Thompson     /// # Ok(())
376c68be7a2SJeremy L Thompson     /// # }
3779df49d7eSJed Brown     /// ```
378594ef120SJeremy L Thompson     pub fn vector<'a>(&self, n: usize) -> Result<Vector<'a>> {
3799df49d7eSJed Brown         Vector::create(self, n)
3809df49d7eSJed Brown     }
3819df49d7eSJed Brown 
3829df49d7eSJed Brown     /// Create a Vector initialized with the data (copied) from a slice
3839df49d7eSJed Brown     ///
3849df49d7eSJed Brown     /// # arguments
3859df49d7eSJed Brown     ///
3869df49d7eSJed Brown     /// * `slice` - Slice containing data
3879df49d7eSJed Brown     ///
3889df49d7eSJed Brown     /// ```
3899df49d7eSJed Brown     /// # use libceed::prelude::*;
3904d27c890SJeremy L Thompson     /// # fn main() -> libceed::Result<()> {
3919df49d7eSJed Brown     /// # let ceed = libceed::Ceed::default_init();
392c68be7a2SJeremy L Thompson     /// let vec = ceed.vector_from_slice(&[1., 2., 3.])?;
3939df49d7eSJed Brown     /// assert_eq!(vec.length(), 3);
394c68be7a2SJeremy L Thompson     /// # Ok(())
395c68be7a2SJeremy L Thompson     /// # }
3969df49d7eSJed Brown     /// ```
397594ef120SJeremy L Thompson     pub fn vector_from_slice<'a>(&self, slice: &[crate::Scalar]) -> Result<Vector<'a>> {
3989df49d7eSJed Brown         Vector::from_slice(self, slice)
3999df49d7eSJed Brown     }
4009df49d7eSJed Brown 
4017ed177dbSJed Brown     /// Returns an ElemRestriction, $\mathcal{E}$, which extracts the degrees of
4027ed177dbSJed Brown     ///   freedom for each element from the local vector into the element vector
4037ed177dbSJed Brown     ///   or assembles contributions from each element in the element vector to
4047ed177dbSJed Brown     ///   the local vector
4059df49d7eSJed Brown     ///
4069df49d7eSJed Brown     /// # arguments
4079df49d7eSJed Brown     ///
4089df49d7eSJed Brown     /// * `nelem`      - Number of elements described in the offsets array
4099df49d7eSJed Brown     /// * `elemsize`   - Size (number of "nodes") per element
4109df49d7eSJed Brown     /// * `ncomp`      - Number of field components per interpolation node (1
4119df49d7eSJed Brown     ///                    for scalar fields)
4129df49d7eSJed Brown     /// * `compstride` - Stride between components for the same Lvector "node".
4139df49d7eSJed Brown     ///                    Data for node `i`, component `j`, element `k` can be
4149df49d7eSJed Brown     ///                    found in the Lvector at index
4159df49d7eSJed Brown     ///                    `offsets[i + k*elemsize] + j*compstride`.
4169df49d7eSJed Brown     /// * `lsize`      - The size of the Lvector. This vector may be larger
4179df49d7eSJed Brown     ///                    than the elements and fields given by this
4189df49d7eSJed Brown     ///                    restriction.
4199df49d7eSJed Brown     /// * `mtype`      - Memory type of the offsets array, see CeedMemType
4209df49d7eSJed Brown     /// * `offsets`    - Array of shape `[nelem, elemsize]`. Row `i` holds the
4219df49d7eSJed Brown     ///                    ordered list of the offsets (into the input CeedVector)
4229df49d7eSJed Brown     ///                    for the unknowns corresponding to element `i`, where
4239df49d7eSJed Brown     ///                    `0 <= i < nelem`. All offsets must be in the range
4249df49d7eSJed Brown     ///                    `[0, lsize - 1]`.
4259df49d7eSJed Brown     ///
4269df49d7eSJed Brown     /// ```
4279df49d7eSJed Brown     /// # use libceed::prelude::*;
4284d27c890SJeremy L Thompson     /// # fn main() -> libceed::Result<()> {
4299df49d7eSJed Brown     /// # let ceed = libceed::Ceed::default_init();
4309df49d7eSJed Brown     /// let nelem = 3;
4319df49d7eSJed Brown     /// let mut ind: Vec<i32> = vec![0; 2 * nelem];
4329df49d7eSJed Brown     /// for i in 0..nelem {
4339df49d7eSJed Brown     ///     ind[2 * i + 0] = i as i32;
4349df49d7eSJed Brown     ///     ind[2 * i + 1] = (i + 1) as i32;
4359df49d7eSJed Brown     /// }
436c68be7a2SJeremy L Thompson     /// let r = ceed.elem_restriction(nelem, 2, 1, 1, nelem + 1, MemType::Host, &ind)?;
437c68be7a2SJeremy L Thompson     /// # Ok(())
438c68be7a2SJeremy L Thompson     /// # }
4399df49d7eSJed Brown     /// ```
440594ef120SJeremy L Thompson     pub fn elem_restriction<'a>(
4419df49d7eSJed Brown         &self,
4429df49d7eSJed Brown         nelem: usize,
4439df49d7eSJed Brown         elemsize: usize,
4449df49d7eSJed Brown         ncomp: usize,
4459df49d7eSJed Brown         compstride: usize,
4469df49d7eSJed Brown         lsize: usize,
4479df49d7eSJed Brown         mtype: MemType,
4489df49d7eSJed Brown         offsets: &[i32],
449594ef120SJeremy L Thompson     ) -> Result<ElemRestriction<'a>> {
4509df49d7eSJed Brown         ElemRestriction::create(
4519df49d7eSJed Brown             self, nelem, elemsize, ncomp, compstride, lsize, mtype, offsets,
452709403c1SSebastian Grimberg         )
453709403c1SSebastian Grimberg     }
454709403c1SSebastian Grimberg 
455709403c1SSebastian Grimberg     /// Returns an oriented ElemRestriction, $\mathcal{E}$, which extracts the
456709403c1SSebastian Grimberg     ///   degrees of freedom for each element from the local vector into the
457709403c1SSebastian Grimberg     ///   element vector or assembles contributions from each element in the
458709403c1SSebastian Grimberg     ///   element vector to the local vector
459709403c1SSebastian Grimberg     ///
460709403c1SSebastian Grimberg     /// # arguments
461709403c1SSebastian Grimberg     ///
462709403c1SSebastian Grimberg     /// * `nelem`      - Number of elements described in the offsets array
463709403c1SSebastian Grimberg     /// * `elemsize`   - Size (number of "nodes") per element
464709403c1SSebastian Grimberg     /// * `ncomp`      - Number of field components per interpolation node (1
465709403c1SSebastian Grimberg     ///                    for scalar fields)
466709403c1SSebastian Grimberg     /// * `compstride` - Stride between components for the same Lvector "node".
467709403c1SSebastian Grimberg     ///                    Data for node `i`, component `j`, element `k` can be
468709403c1SSebastian Grimberg     ///                    found in the Lvector at index
469709403c1SSebastian Grimberg     ///                    `offsets[i + k*elemsize] + j*compstride`.
470709403c1SSebastian Grimberg     /// * `lsize`      - The size of the Lvector. This vector may be larger
471709403c1SSebastian Grimberg     ///                    than the elements and fields given by this
472709403c1SSebastian Grimberg     ///                    restriction.
473709403c1SSebastian Grimberg     /// * `mtype`      - Memory type of the offsets array, see CeedMemType
474709403c1SSebastian Grimberg     /// * `offsets`    - Array of shape `[nelem, elemsize]`. Row `i` holds the
475709403c1SSebastian Grimberg     ///                    ordered list of the offsets (into the input CeedVector)
476709403c1SSebastian Grimberg     ///                    for the unknowns corresponding to element `i`, where
477709403c1SSebastian Grimberg     ///                    `0 <= i < nelem`. All offsets must be in the range
478709403c1SSebastian Grimberg     ///                    `[0, lsize - 1]`.
479709403c1SSebastian Grimberg     /// * `orients`    - Array of shape `[nelem, elemsize]`. Row `i` holds the
480709403c1SSebastian Grimberg     ///                    ordered list of the orientations for the unknowns
481709403c1SSebastian Grimberg     ///                    corresponding to element `i`, with bool `false` used
482709403c1SSebastian Grimberg     ///                    for positively oriented and `true` to flip the
483709403c1SSebastian Grimberg     ///                    orientation.
484709403c1SSebastian Grimberg     ///
485709403c1SSebastian Grimberg     /// ```
486709403c1SSebastian Grimberg     /// # use libceed::prelude::*;
487709403c1SSebastian Grimberg     /// # fn main() -> libceed::Result<()> {
488709403c1SSebastian Grimberg     /// # let ceed = libceed::Ceed::default_init();
489709403c1SSebastian Grimberg     /// let nelem = 3;
490709403c1SSebastian Grimberg     /// let mut ind: Vec<i32> = vec![0; 2 * nelem];
491709403c1SSebastian Grimberg     /// let mut orients: Vec<bool> = vec![false; 2 * nelem];
492709403c1SSebastian Grimberg     /// for i in 0..nelem {
493709403c1SSebastian Grimberg     ///     ind[2 * i + 0] = i as i32;
494709403c1SSebastian Grimberg     ///     ind[2 * i + 1] = (i + 1) as i32;
495709403c1SSebastian Grimberg     ///     orients[2 * i + 0] = (i % 2) > 0; // flip the dofs on element 1, 3, ...
496709403c1SSebastian Grimberg     ///     orients[2 * i + 1] = (i % 2) > 0;
497709403c1SSebastian Grimberg     /// }
498709403c1SSebastian Grimberg     /// let r =
499709403c1SSebastian Grimberg     ///     ceed.oriented_elem_restriction(nelem, 2, 1, 1, nelem + 1, MemType::Host, &ind, &orients)?;
500709403c1SSebastian Grimberg     /// # Ok(())
501709403c1SSebastian Grimberg     /// # }
502709403c1SSebastian Grimberg     /// ```
503709403c1SSebastian Grimberg     pub fn oriented_elem_restriction<'a>(
504709403c1SSebastian Grimberg         &self,
505709403c1SSebastian Grimberg         nelem: usize,
506709403c1SSebastian Grimberg         elemsize: usize,
507709403c1SSebastian Grimberg         ncomp: usize,
508709403c1SSebastian Grimberg         compstride: usize,
509709403c1SSebastian Grimberg         lsize: usize,
510709403c1SSebastian Grimberg         mtype: MemType,
511709403c1SSebastian Grimberg         offsets: &[i32],
512709403c1SSebastian Grimberg         orients: &[bool],
513709403c1SSebastian Grimberg     ) -> Result<ElemRestriction<'a>> {
514709403c1SSebastian Grimberg         ElemRestriction::create_oriented(
515709403c1SSebastian Grimberg             self, nelem, elemsize, ncomp, compstride, lsize, mtype, offsets, orients,
516709403c1SSebastian Grimberg         )
517709403c1SSebastian Grimberg     }
518709403c1SSebastian Grimberg 
519709403c1SSebastian Grimberg     /// Returns a curl-oriented ElemRestriction, $\mathcal{E}$, which extracts
520709403c1SSebastian Grimberg     ///   the degrees of freedom for each element from the local vector into
521709403c1SSebastian Grimberg     ///   the element vector or assembles contributions from each element in
522709403c1SSebastian Grimberg     ///   the element vector to the local vector
523709403c1SSebastian Grimberg     ///
524709403c1SSebastian Grimberg     /// # arguments
525709403c1SSebastian Grimberg     ///
526709403c1SSebastian Grimberg     /// * `nelem`       - Number of elements described in the offsets array
527709403c1SSebastian Grimberg     /// * `elemsize`    - Size (number of "nodes") per element
528709403c1SSebastian Grimberg     /// * `ncomp`       - Number of field components per interpolation node (1
529709403c1SSebastian Grimberg     ///                     for scalar fields)
530709403c1SSebastian Grimberg     /// * `compstride`  - Stride between components for the same Lvector "node".
531709403c1SSebastian Grimberg     ///                     Data for node `i`, component `j`, element `k` can be
532709403c1SSebastian Grimberg     ///                     found in the Lvector at index
533709403c1SSebastian Grimberg     ///                     `offsets[i + k*elemsize] + j*compstride`.
534709403c1SSebastian Grimberg     /// * `lsize`       - The size of the Lvector. This vector may be larger
535709403c1SSebastian Grimberg     ///                     than the elements and fields given by this
536709403c1SSebastian Grimberg     ///                     restriction.
537709403c1SSebastian Grimberg     /// * `mtype`       - Memory type of the offsets array, see CeedMemType
538709403c1SSebastian Grimberg     /// * `offsets`     - Array of shape `[nelem, elemsize]`. Row `i` holds the
539709403c1SSebastian Grimberg     ///                     ordered list of the offsets (into the input CeedVector)
540709403c1SSebastian Grimberg     ///                     for the unknowns corresponding to element `i`, where
541709403c1SSebastian Grimberg     ///                     `0 <= i < nelem`. All offsets must be in the range
542709403c1SSebastian Grimberg     ///                     `[0, lsize - 1]`.
543709403c1SSebastian Grimberg     /// * `curlorients` - Array of shape `[nelem, 3 * elemsize]`. Row `i` holds
544709403c1SSebastian Grimberg     ///                     a row-major tridiagonal matrix (`curlorients[i, 0] =
545709403c1SSebastian Grimberg     ///                     curlorients[i, 3 * elemsize - 1] = 0`, where
546709403c1SSebastian Grimberg     ///                     `0 <= i < nelem`) which is applied to the element
547709403c1SSebastian Grimberg     ///                     unknowns upon restriction.
548709403c1SSebastian Grimberg     ///
549709403c1SSebastian Grimberg     /// ```
550709403c1SSebastian Grimberg     /// # use libceed::prelude::*;
551709403c1SSebastian Grimberg     /// # fn main() -> libceed::Result<()> {
552709403c1SSebastian Grimberg     /// # let ceed = libceed::Ceed::default_init();
553709403c1SSebastian Grimberg     /// let nelem = 3;
554709403c1SSebastian Grimberg     /// let mut ind: Vec<i32> = vec![0; 2 * nelem];
555709403c1SSebastian Grimberg     /// let mut curlorients: Vec<i8> = vec![0; 3 * 2 * nelem];
556709403c1SSebastian Grimberg     /// for i in 0..nelem {
557709403c1SSebastian Grimberg     ///     ind[2 * i + 0] = i as i32;
558709403c1SSebastian Grimberg     ///     ind[2 * i + 1] = (i + 1) as i32;
559709403c1SSebastian Grimberg     ///     curlorients[3 * 2 * i] = 0 as i8;
560709403c1SSebastian Grimberg     ///     curlorients[3 * 2 * (i + 1) - 1] = 0 as i8;
561709403c1SSebastian Grimberg     ///     if (i % 2 > 0) {
562709403c1SSebastian Grimberg     ///         // T = [0  -1]
563709403c1SSebastian Grimberg     ///         //     [-1  0]
564709403c1SSebastian Grimberg     ///         curlorients[3 * 2 * i + 1] = 0 as i8;
565709403c1SSebastian Grimberg     ///         curlorients[3 * 2 * i + 2] = -1 as i8;
566709403c1SSebastian Grimberg     ///         curlorients[3 * 2 * i + 3] = -1 as i8;
567709403c1SSebastian Grimberg     ///         curlorients[3 * 2 * i + 4] = 0 as i8;
568709403c1SSebastian Grimberg     ///     } else {
569709403c1SSebastian Grimberg     ///         // T = I
570709403c1SSebastian Grimberg     ///         curlorients[3 * 2 * i + 1] = 1 as i8;
571709403c1SSebastian Grimberg     ///         curlorients[3 * 2 * i + 2] = 0 as i8;
572709403c1SSebastian Grimberg     ///         curlorients[3 * 2 * i + 3] = 0 as i8;
573709403c1SSebastian Grimberg     ///         curlorients[3 * 2 * i + 4] = 1 as i8;
574709403c1SSebastian Grimberg     ///     }
575709403c1SSebastian Grimberg     /// }
576709403c1SSebastian Grimberg     /// let r = ceed.curl_oriented_elem_restriction(
577709403c1SSebastian Grimberg     ///     nelem,
578709403c1SSebastian Grimberg     ///     2,
579709403c1SSebastian Grimberg     ///     1,
580709403c1SSebastian Grimberg     ///     1,
581709403c1SSebastian Grimberg     ///     nelem + 1,
582709403c1SSebastian Grimberg     ///     MemType::Host,
583709403c1SSebastian Grimberg     ///     &ind,
584709403c1SSebastian Grimberg     ///     &curlorients,
585709403c1SSebastian Grimberg     /// )?;
586709403c1SSebastian Grimberg     /// # Ok(())
587709403c1SSebastian Grimberg     /// # }
588709403c1SSebastian Grimberg     /// ```
589709403c1SSebastian Grimberg     pub fn curl_oriented_elem_restriction<'a>(
590709403c1SSebastian Grimberg         &self,
591709403c1SSebastian Grimberg         nelem: usize,
592709403c1SSebastian Grimberg         elemsize: usize,
593709403c1SSebastian Grimberg         ncomp: usize,
594709403c1SSebastian Grimberg         compstride: usize,
595709403c1SSebastian Grimberg         lsize: usize,
596709403c1SSebastian Grimberg         mtype: MemType,
597709403c1SSebastian Grimberg         offsets: &[i32],
598709403c1SSebastian Grimberg         curlorients: &[i8],
599709403c1SSebastian Grimberg     ) -> Result<ElemRestriction<'a>> {
600709403c1SSebastian Grimberg         ElemRestriction::create_curl_oriented(
601709403c1SSebastian Grimberg             self,
602709403c1SSebastian Grimberg             nelem,
603709403c1SSebastian Grimberg             elemsize,
604709403c1SSebastian Grimberg             ncomp,
605709403c1SSebastian Grimberg             compstride,
606709403c1SSebastian Grimberg             lsize,
607709403c1SSebastian Grimberg             mtype,
608709403c1SSebastian Grimberg             offsets,
609709403c1SSebastian Grimberg             curlorients,
6109df49d7eSJed Brown         )
6119df49d7eSJed Brown     }
6129df49d7eSJed Brown 
6137ed177dbSJed Brown     /// Returns an ElemRestriction, $\mathcal{E}$, from an local vector to
6147ed177dbSJed Brown     ///   an element vector where data can be indexed from the `strides` array
6159df49d7eSJed Brown     ///
6169df49d7eSJed Brown     /// # arguments
6179df49d7eSJed Brown     ///
6189df49d7eSJed Brown     /// * `nelem`      - Number of elements described in the offsets array
6199df49d7eSJed Brown     /// * `elemsize`   - Size (number of "nodes") per element
6209df49d7eSJed Brown     /// * `ncomp`      - Number of field components per interpolation node (1
6219df49d7eSJed Brown     ///                    for scalar fields)
6229df49d7eSJed Brown     /// * `lsize`      - The size of the Lvector. This vector may be larger
6239df49d7eSJed Brown     ///   than the elements and fields given by this restriction.
6249df49d7eSJed Brown     /// * `strides`   - Array for strides between `[nodes, components, elements]`.
6259df49d7eSJed Brown     ///                   Data for node `i`, component `j`, element `k` can be
6269df49d7eSJed Brown     ///                   found in the Lvector at index
6279df49d7eSJed Brown     ///                   `i*strides[0] + j*strides[1] + k*strides[2]`.
6289df49d7eSJed Brown     ///                   CEED_STRIDES_BACKEND may be used with vectors created
6299df49d7eSJed Brown     ///                   by a Ceed backend.
6309df49d7eSJed Brown     ///
6319df49d7eSJed Brown     /// ```
6329df49d7eSJed Brown     /// # use libceed::prelude::*;
6334d27c890SJeremy L Thompson     /// # fn main() -> libceed::Result<()> {
6349df49d7eSJed Brown     /// # let ceed = libceed::Ceed::default_init();
6359df49d7eSJed Brown     /// let nelem = 3;
6369df49d7eSJed Brown     /// let strides: [i32; 3] = [1, 2, 2];
637c68be7a2SJeremy L Thompson     /// let r = ceed.strided_elem_restriction(nelem, 2, 1, nelem * 2, strides)?;
638c68be7a2SJeremy L Thompson     /// # Ok(())
639c68be7a2SJeremy L Thompson     /// # }
6409df49d7eSJed Brown     /// ```
641594ef120SJeremy L Thompson     pub fn strided_elem_restriction<'a>(
6429df49d7eSJed Brown         &self,
6439df49d7eSJed Brown         nelem: usize,
6449df49d7eSJed Brown         elemsize: usize,
6459df49d7eSJed Brown         ncomp: usize,
6469df49d7eSJed Brown         lsize: usize,
6479df49d7eSJed Brown         strides: [i32; 3],
648594ef120SJeremy L Thompson     ) -> Result<ElemRestriction<'a>> {
6499df49d7eSJed Brown         ElemRestriction::create_strided(self, nelem, elemsize, ncomp, lsize, strides)
6509df49d7eSJed Brown     }
6519df49d7eSJed Brown 
6527ed177dbSJed Brown     /// Returns an $H^1$ tensor-product Basis
6539df49d7eSJed Brown     ///
6549df49d7eSJed Brown     /// # arguments
6559df49d7eSJed Brown     ///
6569df49d7eSJed Brown     /// * `dim`       - Topological dimension of element
6579df49d7eSJed Brown     /// * `ncomp`     - Number of field components (1 for scalar fields)
6589df49d7eSJed Brown     /// * `P1d`       - Number of Gauss-Lobatto nodes in one dimension.  The
6599df49d7eSJed Brown     ///                   polynomial degree of the resulting `Q_k` element is
6609df49d7eSJed Brown     ///                   `k=P-1`.
6619df49d7eSJed Brown     /// * `Q1d`       - Number of quadrature points in one dimension
6629df49d7eSJed Brown     /// * `interp1d`  - Row-major `(Q1d * P1d)` matrix expressing the values of
6639df49d7eSJed Brown     ///                   nodal basis functions at quadrature points
6649df49d7eSJed Brown     /// * `grad1d`    - Row-major `(Q1d * P1d)` matrix expressing derivatives of
6659df49d7eSJed Brown     ///                   nodal basis functions at quadrature points
6669df49d7eSJed Brown     /// * `qref1d`    - Array of length `Q1d` holding the locations of quadrature
6679df49d7eSJed Brown     ///                   points on the 1D reference element `[-1, 1]`
6689df49d7eSJed Brown     /// * `qweight1d` - Array of length `Q1d` holding the quadrature weights on
6699df49d7eSJed Brown     ///                   the reference element
6709df49d7eSJed Brown     ///
6719df49d7eSJed Brown     /// ```
6729df49d7eSJed Brown     /// # use libceed::prelude::*;
6734d27c890SJeremy L Thompson     /// # fn main() -> libceed::Result<()> {
6749df49d7eSJed Brown     /// # let ceed = libceed::Ceed::default_init();
6759df49d7eSJed Brown     /// let interp1d  = [ 0.62994317,  0.47255875, -0.14950343,  0.04700152,
6769df49d7eSJed Brown     ///                  -0.07069480,  0.97297619,  0.13253993, -0.03482132,
6779df49d7eSJed Brown     ///                  -0.03482132,  0.13253993,  0.97297619, -0.07069480,
6789df49d7eSJed Brown     ///                   0.04700152, -0.14950343,  0.47255875,  0.62994317];
6799df49d7eSJed Brown     /// let grad1d    = [-2.34183742,  2.78794489, -0.63510411,  0.18899664,
6809df49d7eSJed Brown     ///                  -0.51670214, -0.48795249,  1.33790510, -0.33325047,
6819df49d7eSJed Brown     //                    0.33325047, -1.33790510,  0.48795249,  0.51670214,
6829df49d7eSJed Brown     ///                  -0.18899664,  0.63510411, -2.78794489,  2.34183742];
6839df49d7eSJed Brown     /// let qref1d    = [-0.86113631, -0.33998104,  0.33998104,  0.86113631];
6849df49d7eSJed Brown     /// let qweight1d = [ 0.34785485,  0.65214515,  0.65214515,  0.34785485];
6859df49d7eSJed Brown     /// let b = ceed.
686c68be7a2SJeremy L Thompson     /// basis_tensor_H1(2, 1, 4, 4, &interp1d, &grad1d, &qref1d, &qweight1d)?;
687c68be7a2SJeremy L Thompson     /// # Ok(())
688c68be7a2SJeremy L Thompson     /// # }
6899df49d7eSJed Brown     /// ```
690594ef120SJeremy L Thompson     pub fn basis_tensor_H1<'a>(
6919df49d7eSJed Brown         &self,
6929df49d7eSJed Brown         dim: usize,
6939df49d7eSJed Brown         ncomp: usize,
6949df49d7eSJed Brown         P1d: usize,
6959df49d7eSJed Brown         Q1d: usize,
69680a9ef05SNatalie Beams         interp1d: &[crate::Scalar],
69780a9ef05SNatalie Beams         grad1d: &[crate::Scalar],
69880a9ef05SNatalie Beams         qref1d: &[crate::Scalar],
69980a9ef05SNatalie Beams         qweight1d: &[crate::Scalar],
700594ef120SJeremy L Thompson     ) -> Result<Basis<'a>> {
7019df49d7eSJed Brown         Basis::create_tensor_H1(
7029df49d7eSJed Brown             self, dim, ncomp, P1d, Q1d, interp1d, grad1d, qref1d, qweight1d,
7039df49d7eSJed Brown         )
7049df49d7eSJed Brown     }
7059df49d7eSJed Brown 
7067ed177dbSJed Brown     /// Returns an $H^1$ Lagrange tensor-product Basis
7079df49d7eSJed Brown     ///
7089df49d7eSJed Brown     /// # arguments
7099df49d7eSJed Brown     ///
7109df49d7eSJed Brown     /// * `dim`   - Topological dimension of element
7119df49d7eSJed Brown     /// * `ncomp` - Number of field components (1 for scalar fields)
7129df49d7eSJed Brown     /// * `P`     - Number of Gauss-Lobatto nodes in one dimension.  The
7139df49d7eSJed Brown     ///               polynomial degree of the resulting `Q_k` element is `k=P-1`.
7149df49d7eSJed Brown     /// * `Q`     - Number of quadrature points in one dimension
7159df49d7eSJed Brown     /// * `qmode` - Distribution of the `Q` quadrature points (affects order of
7169df49d7eSJed Brown     ///               accuracy for the quadrature)
7179df49d7eSJed Brown     ///
7189df49d7eSJed Brown     /// ```
7199df49d7eSJed Brown     /// # use libceed::prelude::*;
7204d27c890SJeremy L Thompson     /// # fn main() -> libceed::Result<()> {
7219df49d7eSJed Brown     /// # let ceed = libceed::Ceed::default_init();
722c68be7a2SJeremy L Thompson     /// let b = ceed.basis_tensor_H1_Lagrange(2, 1, 3, 4, QuadMode::Gauss)?;
723c68be7a2SJeremy L Thompson     /// # Ok(())
724c68be7a2SJeremy L Thompson     /// # }
7259df49d7eSJed Brown     /// ```
726594ef120SJeremy L Thompson     pub fn basis_tensor_H1_Lagrange<'a>(
7279df49d7eSJed Brown         &self,
7289df49d7eSJed Brown         dim: usize,
7299df49d7eSJed Brown         ncomp: usize,
7309df49d7eSJed Brown         P: usize,
7319df49d7eSJed Brown         Q: usize,
7329df49d7eSJed Brown         qmode: QuadMode,
733594ef120SJeremy L Thompson     ) -> Result<Basis<'a>> {
7349df49d7eSJed Brown         Basis::create_tensor_H1_Lagrange(self, dim, ncomp, P, Q, qmode)
7359df49d7eSJed Brown     }
7369df49d7eSJed Brown 
7377ed177dbSJed Brown     /// Returns an $H-1$ Basis
7389df49d7eSJed Brown     ///
7399df49d7eSJed Brown     /// # arguments
7409df49d7eSJed Brown     ///
7419df49d7eSJed Brown     /// * `topo`    - Topology of element, e.g. hypercube, simplex, ect
7429df49d7eSJed Brown     /// * `ncomp`   - Number of field components (1 for scalar fields)
7439df49d7eSJed Brown     /// * `nnodes`  - Total number of nodes
7449df49d7eSJed Brown     /// * `nqpts`   - Total number of quadrature points
7459df49d7eSJed Brown     /// * `interp`  - Row-major `(nqpts * nnodes)` matrix expressing the values of
7469df49d7eSJed Brown     ///                 nodal basis functions at quadrature points
74797c1c57aSSebastian Grimberg     /// * `grad`    - Row-major `(dim * nqpts * nnodes)` matrix expressing
7489df49d7eSJed Brown     ///                 derivatives of nodal basis functions at quadrature points
7499df49d7eSJed Brown     /// * `qref`    - Array of length `nqpts` holding the locations of quadrature
7509df49d7eSJed Brown     ///                 points on the reference element `[-1, 1]`
7519df49d7eSJed Brown     /// * `qweight` - Array of length `nqpts` holding the quadrature weights on
7529df49d7eSJed Brown     ///                 the reference element
7539df49d7eSJed Brown     ///
7549df49d7eSJed Brown     /// ```
7559df49d7eSJed Brown     /// # use libceed::prelude::*;
7564d27c890SJeremy L Thompson     /// # fn main() -> libceed::Result<()> {
7579df49d7eSJed Brown     /// # let ceed = libceed::Ceed::default_init();
7589df49d7eSJed Brown     /// let interp = [
7599df49d7eSJed Brown     ///     0.12000000,
7609df49d7eSJed Brown     ///     0.48000000,
7619df49d7eSJed Brown     ///     -0.12000000,
7629df49d7eSJed Brown     ///     0.48000000,
7639df49d7eSJed Brown     ///     0.16000000,
7649df49d7eSJed Brown     ///     -0.12000000,
7659df49d7eSJed Brown     ///     -0.12000000,
7669df49d7eSJed Brown     ///     0.48000000,
7679df49d7eSJed Brown     ///     0.12000000,
7689df49d7eSJed Brown     ///     0.16000000,
7699df49d7eSJed Brown     ///     0.48000000,
7709df49d7eSJed Brown     ///     -0.12000000,
7719df49d7eSJed Brown     ///     -0.11111111,
7729df49d7eSJed Brown     ///     0.44444444,
7739df49d7eSJed Brown     ///     -0.11111111,
7749df49d7eSJed Brown     ///     0.44444444,
7759df49d7eSJed Brown     ///     0.44444444,
7769df49d7eSJed Brown     ///     -0.11111111,
7779df49d7eSJed Brown     ///     -0.12000000,
7789df49d7eSJed Brown     ///     0.16000000,
7799df49d7eSJed Brown     ///     -0.12000000,
7809df49d7eSJed Brown     ///     0.48000000,
7819df49d7eSJed Brown     ///     0.48000000,
7829df49d7eSJed Brown     ///     0.12000000,
7839df49d7eSJed Brown     /// ];
7849df49d7eSJed Brown     /// let grad = [
7859df49d7eSJed Brown     ///     -1.40000000,
7869df49d7eSJed Brown     ///     1.60000000,
7879df49d7eSJed Brown     ///     -0.20000000,
7889df49d7eSJed Brown     ///     -0.80000000,
7899df49d7eSJed Brown     ///     0.80000000,
7909df49d7eSJed Brown     ///     0.00000000,
7919df49d7eSJed Brown     ///     0.20000000,
7929df49d7eSJed Brown     ///     -1.60000000,
7939df49d7eSJed Brown     ///     1.40000000,
7949df49d7eSJed Brown     ///     -0.80000000,
7959df49d7eSJed Brown     ///     0.80000000,
7969df49d7eSJed Brown     ///     0.00000000,
7979df49d7eSJed Brown     ///     -0.33333333,
7989df49d7eSJed Brown     ///     0.00000000,
7999df49d7eSJed Brown     ///     0.33333333,
8009df49d7eSJed Brown     ///     -1.33333333,
8019df49d7eSJed Brown     ///     1.33333333,
8029df49d7eSJed Brown     ///     0.00000000,
8039df49d7eSJed Brown     ///     0.20000000,
8049df49d7eSJed Brown     ///     0.00000000,
8059df49d7eSJed Brown     ///     -0.20000000,
8069df49d7eSJed Brown     ///     -2.40000000,
8079df49d7eSJed Brown     ///     2.40000000,
8089df49d7eSJed Brown     ///     0.00000000,
8099df49d7eSJed Brown     ///     -1.40000000,
8109df49d7eSJed Brown     ///     -0.80000000,
8119df49d7eSJed Brown     ///     0.00000000,
8129df49d7eSJed Brown     ///     1.60000000,
8139df49d7eSJed Brown     ///     0.80000000,
8149df49d7eSJed Brown     ///     -0.20000000,
8159df49d7eSJed Brown     ///     0.20000000,
8169df49d7eSJed Brown     ///     -2.40000000,
8179df49d7eSJed Brown     ///     0.00000000,
8189df49d7eSJed Brown     ///     0.00000000,
8199df49d7eSJed Brown     ///     2.40000000,
8209df49d7eSJed Brown     ///     -0.20000000,
8219df49d7eSJed Brown     ///     -0.33333333,
8229df49d7eSJed Brown     ///     -1.33333333,
8239df49d7eSJed Brown     ///     0.00000000,
8249df49d7eSJed Brown     ///     0.00000000,
8259df49d7eSJed Brown     ///     1.33333333,
8269df49d7eSJed Brown     ///     0.33333333,
8279df49d7eSJed Brown     ///     0.20000000,
8289df49d7eSJed Brown     ///     -0.80000000,
8299df49d7eSJed Brown     ///     0.00000000,
8309df49d7eSJed Brown     ///     -1.60000000,
8319df49d7eSJed Brown     ///     0.80000000,
8329df49d7eSJed Brown     ///     1.40000000,
8339df49d7eSJed Brown     /// ];
8349df49d7eSJed Brown     /// let qref = [
8359df49d7eSJed Brown     ///     0.20000000, 0.60000000, 0.33333333, 0.20000000, 0.20000000, 0.20000000, 0.33333333,
8369df49d7eSJed Brown     ///     0.60000000,
8379df49d7eSJed Brown     /// ];
8389df49d7eSJed Brown     /// let qweight = [0.26041667, 0.26041667, -0.28125000, 0.26041667];
839c68be7a2SJeremy L Thompson     /// let b = ceed.basis_H1(
8409df49d7eSJed Brown     ///     ElemTopology::Triangle,
8419df49d7eSJed Brown     ///     1,
8429df49d7eSJed Brown     ///     6,
8439df49d7eSJed Brown     ///     4,
8449df49d7eSJed Brown     ///     &interp,
8459df49d7eSJed Brown     ///     &grad,
8469df49d7eSJed Brown     ///     &qref,
8479df49d7eSJed Brown     ///     &qweight,
848c68be7a2SJeremy L Thompson     /// )?;
849c68be7a2SJeremy L Thompson     /// # Ok(())
850c68be7a2SJeremy L Thompson     /// # }
8519df49d7eSJed Brown     /// ```
852594ef120SJeremy L Thompson     pub fn basis_H1<'a>(
8539df49d7eSJed Brown         &self,
8549df49d7eSJed Brown         topo: ElemTopology,
8559df49d7eSJed Brown         ncomp: usize,
8569df49d7eSJed Brown         nnodes: usize,
8579df49d7eSJed Brown         nqpts: usize,
85880a9ef05SNatalie Beams         interp: &[crate::Scalar],
85980a9ef05SNatalie Beams         grad: &[crate::Scalar],
86080a9ef05SNatalie Beams         qref: &[crate::Scalar],
86180a9ef05SNatalie Beams         qweight: &[crate::Scalar],
862594ef120SJeremy L Thompson     ) -> Result<Basis<'a>> {
8639df49d7eSJed Brown         Basis::create_H1(
8649df49d7eSJed Brown             self, topo, ncomp, nnodes, nqpts, interp, grad, qref, qweight,
8659df49d7eSJed Brown         )
8669df49d7eSJed Brown     }
8679df49d7eSJed Brown 
86897c1c57aSSebastian Grimberg     /// Returns an $H(div)$ Basis
86997c1c57aSSebastian Grimberg     ///
87097c1c57aSSebastian Grimberg     /// # arguments
87197c1c57aSSebastian Grimberg     ///
87297c1c57aSSebastian Grimberg     /// * `topo`    - Topology of element, e.g. hypercube, simplex, ect
87397c1c57aSSebastian Grimberg     /// * `ncomp`   - Number of field components (1 for scalar fields)
87497c1c57aSSebastian Grimberg     /// * `nnodes`  - Total number of nodes
87597c1c57aSSebastian Grimberg     /// * `nqpts`   - Total number of quadrature points
87697c1c57aSSebastian Grimberg     /// * `interp`  - Row-major `(dim * nqpts * nnodes)` matrix expressing the
87797c1c57aSSebastian Grimberg     ///                 values of basis functions at quadrature points
87897c1c57aSSebastian Grimberg     /// * `div`     - Row-major `(nqpts * nnodes)` matrix expressing the
87997c1c57aSSebastian Grimberg     ///                 divergence of basis functions at quadrature points
88097c1c57aSSebastian Grimberg     /// * `qref`    - Array of length `nqpts` holding the locations of quadrature
88197c1c57aSSebastian Grimberg     ///                 points on the reference element `[-1, 1]`
88297c1c57aSSebastian Grimberg     /// * `qweight` - Array of length `nqpts` holding the quadrature weights on
88397c1c57aSSebastian Grimberg     ///                 the reference element
88497c1c57aSSebastian Grimberg     ///
88597c1c57aSSebastian Grimberg     /// ```
88697c1c57aSSebastian Grimberg     /// # use libceed::prelude::*;
88797c1c57aSSebastian Grimberg     /// # fn main() -> libceed::Result<()> {
88897c1c57aSSebastian Grimberg     /// # let ceed = libceed::Ceed::default_init();
88997c1c57aSSebastian Grimberg     /// let interp = [
89097c1c57aSSebastian Grimberg     ///     0.00000000,
89197c1c57aSSebastian Grimberg     ///     -1.57735027,
89297c1c57aSSebastian Grimberg     ///     0.57735027,
89397c1c57aSSebastian Grimberg     ///     0.00000000,
89497c1c57aSSebastian Grimberg     ///     0.00000000,
89597c1c57aSSebastian Grimberg     ///     -1.57735027,
89697c1c57aSSebastian Grimberg     ///     0.57735027,
89797c1c57aSSebastian Grimberg     ///     0.00000000,
89897c1c57aSSebastian Grimberg     ///     0.00000000,
89997c1c57aSSebastian Grimberg     ///     -1.57735027,
90097c1c57aSSebastian Grimberg     ///     0.57735027,
90197c1c57aSSebastian Grimberg     ///     0.00000000,
90297c1c57aSSebastian Grimberg     ///     0.00000000,
90397c1c57aSSebastian Grimberg     ///     -0.42264973,
90497c1c57aSSebastian Grimberg     ///     -0.57735027,
90597c1c57aSSebastian Grimberg     ///     0.00000000,
90697c1c57aSSebastian Grimberg     ///     0.42264973,
90797c1c57aSSebastian Grimberg     ///     0.00000000,
90897c1c57aSSebastian Grimberg     ///     0.00000000,
90997c1c57aSSebastian Grimberg     ///     0.57735027,
91097c1c57aSSebastian Grimberg     ///     0.42264973,
91197c1c57aSSebastian Grimberg     ///     0.00000000,
91297c1c57aSSebastian Grimberg     ///     0.00000000,
91397c1c57aSSebastian Grimberg     ///     0.57735027,
91497c1c57aSSebastian Grimberg     ///     1.57735027,
91597c1c57aSSebastian Grimberg     ///     0.00000000,
91697c1c57aSSebastian Grimberg     ///     0.00000000,
91797c1c57aSSebastian Grimberg     ///     -0.57735027,
91897c1c57aSSebastian Grimberg     ///     1.57735027,
91997c1c57aSSebastian Grimberg     ///     0.00000000,
92097c1c57aSSebastian Grimberg     ///     0.00000000,
92197c1c57aSSebastian Grimberg     ///     -0.57735027,
92297c1c57aSSebastian Grimberg     /// ];
92397c1c57aSSebastian Grimberg     /// let div = [
92497c1c57aSSebastian Grimberg     ///     -1.00000000,
92597c1c57aSSebastian Grimberg     ///     1.00000000,
92697c1c57aSSebastian Grimberg     ///     -1.00000000,
92797c1c57aSSebastian Grimberg     ///     1.00000000,
92897c1c57aSSebastian Grimberg     ///     -1.00000000,
92997c1c57aSSebastian Grimberg     ///     1.00000000,
93097c1c57aSSebastian Grimberg     ///     -1.00000000,
93197c1c57aSSebastian Grimberg     ///     1.00000000,
93297c1c57aSSebastian Grimberg     ///     -1.00000000,
93397c1c57aSSebastian Grimberg     ///     1.00000000,
93497c1c57aSSebastian Grimberg     ///     -1.00000000,
93597c1c57aSSebastian Grimberg     ///     1.00000000,
93697c1c57aSSebastian Grimberg     ///     -1.00000000,
93797c1c57aSSebastian Grimberg     ///     1.00000000,
93897c1c57aSSebastian Grimberg     ///     -1.00000000,
93997c1c57aSSebastian Grimberg     ///     1.00000000,
94097c1c57aSSebastian Grimberg     /// ];
94197c1c57aSSebastian Grimberg     /// let qref = [
94297c1c57aSSebastian Grimberg     ///     0.57735026, 0.57735026, 0.57735026, 0.57735026, 0.57735026, 0.57735026, 0.57735026,
94397c1c57aSSebastian Grimberg     ///     0.57735026,
94497c1c57aSSebastian Grimberg     /// ];
94597c1c57aSSebastian Grimberg     /// let qweight = [1.00000000, 1.00000000, 1.00000000, 1.00000000];
94697c1c57aSSebastian Grimberg     /// let b = ceed.basis_Hdiv(ElemTopology::Quad, 1, 4, 4, &interp, &div, &qref, &qweight)?;
94797c1c57aSSebastian Grimberg     /// # Ok(())
94897c1c57aSSebastian Grimberg     /// # }
94997c1c57aSSebastian Grimberg     /// ```
95097c1c57aSSebastian Grimberg     pub fn basis_Hdiv<'a>(
95197c1c57aSSebastian Grimberg         &self,
95297c1c57aSSebastian Grimberg         topo: ElemTopology,
95397c1c57aSSebastian Grimberg         ncomp: usize,
95497c1c57aSSebastian Grimberg         nnodes: usize,
95597c1c57aSSebastian Grimberg         nqpts: usize,
95697c1c57aSSebastian Grimberg         interp: &[crate::Scalar],
95797c1c57aSSebastian Grimberg         div: &[crate::Scalar],
95897c1c57aSSebastian Grimberg         qref: &[crate::Scalar],
95997c1c57aSSebastian Grimberg         qweight: &[crate::Scalar],
96097c1c57aSSebastian Grimberg     ) -> Result<Basis<'a>> {
96197c1c57aSSebastian Grimberg         Basis::create_Hdiv(self, topo, ncomp, nnodes, nqpts, interp, div, qref, qweight)
96297c1c57aSSebastian Grimberg     }
96397c1c57aSSebastian Grimberg 
96497c1c57aSSebastian Grimberg     /// Returns an $H(curl)$ Basis
96597c1c57aSSebastian Grimberg     ///
96697c1c57aSSebastian Grimberg     /// # arguments
96797c1c57aSSebastian Grimberg     ///
96897c1c57aSSebastian Grimberg     /// * `topo`    - Topology of element, e.g. hypercube, simplex, ect
96997c1c57aSSebastian Grimberg     /// * `ncomp`   - Number of field components (1 for scalar fields)
97097c1c57aSSebastian Grimberg     /// * `nnodes`  - Total number of nodes
97197c1c57aSSebastian Grimberg     /// * `nqpts`   - Total number of quadrature points
97297c1c57aSSebastian Grimberg     /// * `interp`  - Row-major `(dim * nqpts * nnodes)` matrix expressing the
97397c1c57aSSebastian Grimberg     ///                 values of basis functions at quadrature points
97497c1c57aSSebastian Grimberg     /// * `curl`    - Row-major `(curl_comp * nqpts * nnodes)`, `curl_comp = 1 if
97597c1c57aSSebastian Grimberg     ///                 dim < 3 else dim` matrix expressing the curl of basis
97697c1c57aSSebastian Grimberg     ///                 functions at quadrature points
97797c1c57aSSebastian Grimberg     /// * `qref`    - Array of length `nqpts` holding the locations of quadrature
97897c1c57aSSebastian Grimberg     ///                 points on the reference element `[-1, 1]`
97997c1c57aSSebastian Grimberg     /// * `qweight` - Array of length `nqpts` holding the quadrature weights on
98097c1c57aSSebastian Grimberg     ///                 the reference element
98197c1c57aSSebastian Grimberg     ///
98297c1c57aSSebastian Grimberg     /// ```
98397c1c57aSSebastian Grimberg     /// # use libceed::prelude::*;
98497c1c57aSSebastian Grimberg     /// # fn main() -> libceed::Result<()> {
98597c1c57aSSebastian Grimberg     /// # let ceed = libceed::Ceed::default_init();
98697c1c57aSSebastian Grimberg     /// let interp = [
98797c1c57aSSebastian Grimberg     ///     -0.20000000,
98897c1c57aSSebastian Grimberg     ///     0.20000000,
98997c1c57aSSebastian Grimberg     ///     0.80000000,
99097c1c57aSSebastian Grimberg     ///     -0.20000000,
99197c1c57aSSebastian Grimberg     ///     0.20000000,
99297c1c57aSSebastian Grimberg     ///     0.80000000,
99397c1c57aSSebastian Grimberg     ///     -0.33333333,
99497c1c57aSSebastian Grimberg     ///     0.33333333,
99597c1c57aSSebastian Grimberg     ///     0.66666667,
99697c1c57aSSebastian Grimberg     ///     -0.60000000,
99797c1c57aSSebastian Grimberg     ///     0.60000000,
99897c1c57aSSebastian Grimberg     ///     0.40000000,
99997c1c57aSSebastian Grimberg     ///     0.20000000,
100097c1c57aSSebastian Grimberg     ///     0.80000000,
100197c1c57aSSebastian Grimberg     ///     0.20000000,
100297c1c57aSSebastian Grimberg     ///     0.60000000,
100397c1c57aSSebastian Grimberg     ///     0.40000000,
100497c1c57aSSebastian Grimberg     ///     0.60000000,
100597c1c57aSSebastian Grimberg     ///     0.33333333,
100697c1c57aSSebastian Grimberg     ///     0.66666667,
100797c1c57aSSebastian Grimberg     ///     0.33333333,
100897c1c57aSSebastian Grimberg     ///     0.20000000,
100997c1c57aSSebastian Grimberg     ///     0.80000000,
101097c1c57aSSebastian Grimberg     ///     0.20000000,
101197c1c57aSSebastian Grimberg     /// ];
101297c1c57aSSebastian Grimberg     /// let curl = [
101397c1c57aSSebastian Grimberg     ///     2.00000000,
101497c1c57aSSebastian Grimberg     ///     -2.00000000,
101597c1c57aSSebastian Grimberg     ///     -2.00000000,
101697c1c57aSSebastian Grimberg     ///     2.00000000,
101797c1c57aSSebastian Grimberg     ///     -2.00000000,
101897c1c57aSSebastian Grimberg     ///     -2.00000000,
101997c1c57aSSebastian Grimberg     ///     2.00000000,
102097c1c57aSSebastian Grimberg     ///     -2.00000000,
102197c1c57aSSebastian Grimberg     ///     -2.00000000,
102297c1c57aSSebastian Grimberg     ///     2.00000000,
102397c1c57aSSebastian Grimberg     ///     -2.00000000,
102497c1c57aSSebastian Grimberg     ///     -2.00000000,
102597c1c57aSSebastian Grimberg     /// ];
102697c1c57aSSebastian Grimberg     /// let qref = [
102797c1c57aSSebastian Grimberg     ///     0.20000000, 0.60000000, 0.33333333, 0.20000000, 0.20000000, 0.20000000, 0.33333333,
102897c1c57aSSebastian Grimberg     ///     0.60000000,
102997c1c57aSSebastian Grimberg     /// ];
103097c1c57aSSebastian Grimberg     /// let qweight = [0.26041667, 0.26041667, -0.28125000, 0.26041667];
103197c1c57aSSebastian Grimberg     /// let b = ceed.basis_Hcurl(
103297c1c57aSSebastian Grimberg     ///     ElemTopology::Triangle,
103397c1c57aSSebastian Grimberg     ///     1,
103497c1c57aSSebastian Grimberg     ///     3,
103597c1c57aSSebastian Grimberg     ///     4,
103697c1c57aSSebastian Grimberg     ///     &interp,
103797c1c57aSSebastian Grimberg     ///     &curl,
103897c1c57aSSebastian Grimberg     ///     &qref,
103997c1c57aSSebastian Grimberg     ///     &qweight,
104097c1c57aSSebastian Grimberg     /// )?;
104197c1c57aSSebastian Grimberg     /// # Ok(())
104297c1c57aSSebastian Grimberg     /// # }
104397c1c57aSSebastian Grimberg     /// ```
104497c1c57aSSebastian Grimberg     pub fn basis_Hcurl<'a>(
104597c1c57aSSebastian Grimberg         &self,
104697c1c57aSSebastian Grimberg         topo: ElemTopology,
104797c1c57aSSebastian Grimberg         ncomp: usize,
104897c1c57aSSebastian Grimberg         nnodes: usize,
104997c1c57aSSebastian Grimberg         nqpts: usize,
105097c1c57aSSebastian Grimberg         interp: &[crate::Scalar],
105197c1c57aSSebastian Grimberg         curl: &[crate::Scalar],
105297c1c57aSSebastian Grimberg         qref: &[crate::Scalar],
105397c1c57aSSebastian Grimberg         qweight: &[crate::Scalar],
105497c1c57aSSebastian Grimberg     ) -> Result<Basis<'a>> {
105597c1c57aSSebastian Grimberg         Basis::create_Hcurl(
105697c1c57aSSebastian Grimberg             self, topo, ncomp, nnodes, nqpts, interp, curl, qref, qweight,
105797c1c57aSSebastian Grimberg         )
105897c1c57aSSebastian Grimberg     }
105997c1c57aSSebastian Grimberg 
10607ed177dbSJed Brown     /// Returns a QFunction for evaluating interior (volumetric) terms
10617ed177dbSJed Brown     ///   of a weak form corresponding to the $L^2$ inner product
10627ed177dbSJed Brown     ///
10637ed177dbSJed Brown     /// $$
10647ed177dbSJed Brown     /// \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),
10657ed177dbSJed Brown     /// $$
10667ed177dbSJed Brown     ///
10677ed177dbSJed Brown     /// where $v \cdot f_0$ represents contraction over fields and $\nabla v : f_1$
10687ed177dbSJed Brown     ///   represents contraction over both fields and spatial dimensions.
10699df49d7eSJed Brown     ///
10709df49d7eSJed Brown     /// # arguments
10719df49d7eSJed Brown     ///
10729df49d7eSJed Brown     /// * `vlength` - Vector length. Caller must ensure that number of
10739df49d7eSJed Brown     ///                 quadrature points is a multiple of vlength.
10747ed177dbSJed Brown     /// * `f`       - Boxed closure to evaluate weak form at quadrature points.
10759df49d7eSJed Brown     ///
10769df49d7eSJed Brown     /// ```
10779df49d7eSJed Brown     /// # use libceed::prelude::*;
10784d27c890SJeremy L Thompson     /// # fn main() -> libceed::Result<()> {
10799df49d7eSJed Brown     /// # let ceed = libceed::Ceed::default_init();
10809df49d7eSJed Brown     /// let mut user_f = |[u, weights, ..]: QFunctionInputs, [v, ..]: QFunctionOutputs| {
10819df49d7eSJed Brown     ///     // Iterate over quadrature points
10829df49d7eSJed Brown     ///     v.iter_mut()
10839df49d7eSJed Brown     ///         .zip(u.iter().zip(weights.iter()))
10849df49d7eSJed Brown     ///         .for_each(|(v, (u, w))| *v = u * w);
10859df49d7eSJed Brown     ///
10869df49d7eSJed Brown     ///     // Return clean error code
10879df49d7eSJed Brown     ///     0
10889df49d7eSJed Brown     /// };
10899df49d7eSJed Brown     ///
1090c68be7a2SJeremy L Thompson     /// let qf = ceed.q_function_interior(1, Box::new(user_f))?;
1091c68be7a2SJeremy L Thompson     /// # Ok(())
1092c68be7a2SJeremy L Thompson     /// # }
10939df49d7eSJed Brown     /// ```
1094594ef120SJeremy L Thompson     pub fn q_function_interior<'a>(
10959df49d7eSJed Brown         &self,
10969df49d7eSJed Brown         vlength: usize,
10979df49d7eSJed Brown         f: Box<qfunction::QFunctionUserClosure>,
1098594ef120SJeremy L Thompson     ) -> Result<QFunction<'a>> {
10999df49d7eSJed Brown         QFunction::create(self, vlength, f)
11009df49d7eSJed Brown     }
11019df49d7eSJed Brown 
11027ed177dbSJed Brown     /// Returns a QFunction for evaluating interior (volumetric) terms
11039df49d7eSJed Brown     ///   created by name
11049df49d7eSJed Brown     ///
11057ed177dbSJed Brown     /// # arguments
11067ed177dbSJed Brown     ///
11077ed177dbSJed Brown     /// * `name` - name of QFunction from libCEED gallery
11087ed177dbSJed Brown     ///
11099df49d7eSJed Brown     /// ```
11109df49d7eSJed Brown     /// # use libceed::prelude::*;
11114d27c890SJeremy L Thompson     /// # fn main() -> libceed::Result<()> {
11129df49d7eSJed Brown     /// # let ceed = libceed::Ceed::default_init();
1113c68be7a2SJeremy L Thompson     /// let qf = ceed.q_function_interior_by_name("Mass1DBuild")?;
1114c68be7a2SJeremy L Thompson     /// # Ok(())
1115c68be7a2SJeremy L Thompson     /// # }
11169df49d7eSJed Brown     /// ```
1117594ef120SJeremy L Thompson     pub fn q_function_interior_by_name<'a>(&self, name: &str) -> Result<QFunctionByName<'a>> {
11189df49d7eSJed Brown         QFunctionByName::create(self, name)
11199df49d7eSJed Brown     }
11209df49d7eSJed Brown 
11217ed177dbSJed Brown     /// Returns an Operator and associate a QFunction. A Basis and
11227ed177dbSJed Brown     ///   ElemRestriction can be associated with QFunction fields via
11239df49d7eSJed Brown     ///   set_field().
11249df49d7eSJed Brown     ///
11257ed177dbSJed Brown     /// # arguments
11267ed177dbSJed Brown     ///
11279df49d7eSJed Brown     /// * `qf`   - QFunction defining the action of the operator at quadrature
11289df49d7eSJed Brown     ///              points
11299df49d7eSJed Brown     /// * `dqf`  - QFunction defining the action of the Jacobian of the qf (or
11309df49d7eSJed Brown     ///              qfunction_none)
11319df49d7eSJed Brown     /// * `dqfT` - QFunction defining the action of the transpose of the
11329df49d7eSJed Brown     ///              Jacobian of the qf (or qfunction_none)
11339df49d7eSJed Brown     ///
11349df49d7eSJed Brown     /// ```
11359df49d7eSJed Brown     /// # use libceed::prelude::*;
11364d27c890SJeremy L Thompson     /// # fn main() -> libceed::Result<()> {
11379df49d7eSJed Brown     /// # let ceed = libceed::Ceed::default_init();
1138c68be7a2SJeremy L Thompson     /// let qf = ceed.q_function_interior_by_name("Mass1DBuild")?;
1139c68be7a2SJeremy L Thompson     /// let op = ceed.operator(&qf, QFunctionOpt::None, QFunctionOpt::None)?;
1140c68be7a2SJeremy L Thompson     /// # Ok(())
1141c68be7a2SJeremy L Thompson     /// # }
11429df49d7eSJed Brown     /// ```
1143594ef120SJeremy L Thompson     pub fn operator<'a, 'b>(
11449df49d7eSJed Brown         &self,
11459df49d7eSJed Brown         qf: impl Into<QFunctionOpt<'b>>,
11469df49d7eSJed Brown         dqf: impl Into<QFunctionOpt<'b>>,
11479df49d7eSJed Brown         dqfT: impl Into<QFunctionOpt<'b>>,
1148594ef120SJeremy L Thompson     ) -> Result<Operator<'a>> {
11499df49d7eSJed Brown         Operator::create(self, qf, dqf, dqfT)
11509df49d7eSJed Brown     }
11519df49d7eSJed Brown 
11529df49d7eSJed Brown     /// Returns an Operator that composes the action of several Operators
11539df49d7eSJed Brown     ///
11549df49d7eSJed Brown     /// ```
11559df49d7eSJed Brown     /// # use libceed::prelude::*;
11564d27c890SJeremy L Thompson     /// # fn main() -> libceed::Result<()> {
11579df49d7eSJed Brown     /// # let ceed = libceed::Ceed::default_init();
1158c68be7a2SJeremy L Thompson     /// let op = ceed.composite_operator()?;
1159c68be7a2SJeremy L Thompson     /// # Ok(())
1160c68be7a2SJeremy L Thompson     /// # }
11619df49d7eSJed Brown     /// ```
1162594ef120SJeremy L Thompson     pub fn composite_operator<'a>(&self) -> Result<CompositeOperator<'a>> {
11639df49d7eSJed Brown         CompositeOperator::create(self)
11649df49d7eSJed Brown     }
11659df49d7eSJed Brown }
11669df49d7eSJed Brown 
11679df49d7eSJed Brown // -----------------------------------------------------------------------------
11689df49d7eSJed Brown // Tests
11699df49d7eSJed Brown // -----------------------------------------------------------------------------
11709df49d7eSJed Brown #[cfg(test)]
11719df49d7eSJed Brown mod tests {
11729df49d7eSJed Brown     use super::*;
11739df49d7eSJed Brown 
117489d15d5fSJeremy L Thompson     fn ceed_t501() -> Result<()> {
11759df49d7eSJed Brown         let resource = "/cpu/self/ref/blocked";
11769df49d7eSJed Brown         let ceed = Ceed::init(resource);
11779df49d7eSJed Brown         let nelem = 4;
11789df49d7eSJed Brown         let p = 3;
11799df49d7eSJed Brown         let q = 4;
11809df49d7eSJed Brown         let ndofs = p * nelem - nelem + 1;
11819df49d7eSJed Brown 
11829df49d7eSJed Brown         // Vectors
11839df49d7eSJed Brown         let x = ceed.vector_from_slice(&[-1., -0.5, 0.0, 0.5, 1.0])?;
11849df49d7eSJed Brown         let mut qdata = ceed.vector(nelem * q)?;
11859df49d7eSJed Brown         qdata.set_value(0.0)?;
11869df49d7eSJed Brown         let mut u = ceed.vector(ndofs)?;
11879df49d7eSJed Brown         u.set_value(1.0)?;
11889df49d7eSJed Brown         let mut v = ceed.vector(ndofs)?;
11899df49d7eSJed Brown         v.set_value(0.0)?;
11909df49d7eSJed Brown 
11919df49d7eSJed Brown         // Restrictions
11929df49d7eSJed Brown         let mut indx: Vec<i32> = vec![0; 2 * nelem];
11939df49d7eSJed Brown         for i in 0..nelem {
11941ce8139fSJeremy L Thompson             for j in 0..2 {
11951ce8139fSJeremy L Thompson                 indx[2 * i + j] = (i + j) as i32;
11961ce8139fSJeremy L Thompson             }
11979df49d7eSJed Brown         }
11989df49d7eSJed Brown         let rx = ceed.elem_restriction(nelem, 2, 1, 1, nelem + 1, MemType::Host, &indx)?;
11999df49d7eSJed Brown         let mut indu: Vec<i32> = vec![0; p * nelem];
12009df49d7eSJed Brown         for i in 0..nelem {
12011ce8139fSJeremy L Thompson             for j in 0..p {
12021ce8139fSJeremy L Thompson                 indu[p * i + j] = (i + j) as i32;
12039df49d7eSJed Brown             }
12041ce8139fSJeremy L Thompson         }
12051ce8139fSJeremy L Thompson         let ru = ceed.elem_restriction(nelem, p, 1, 1, ndofs, MemType::Host, &indu)?;
12069df49d7eSJed Brown         let strides: [i32; 3] = [1, q as i32, q as i32];
12079df49d7eSJed Brown         let rq = ceed.strided_elem_restriction(nelem, q, 1, q * nelem, strides)?;
12089df49d7eSJed Brown 
12099df49d7eSJed Brown         // Bases
12109df49d7eSJed Brown         let bx = ceed.basis_tensor_H1_Lagrange(1, 1, 2, q, QuadMode::Gauss)?;
12119df49d7eSJed Brown         let bu = ceed.basis_tensor_H1_Lagrange(1, 1, p, q, QuadMode::Gauss)?;
12129df49d7eSJed Brown 
12139df49d7eSJed Brown         // Build quadrature data
12149df49d7eSJed Brown         let qf_build = ceed.q_function_interior_by_name("Mass1DBuild")?;
12159df49d7eSJed Brown         ceed.operator(&qf_build, QFunctionOpt::None, QFunctionOpt::None)?
12169df49d7eSJed Brown             .field("dx", &rx, &bx, VectorOpt::Active)?
12179df49d7eSJed Brown             .field("weights", ElemRestrictionOpt::None, &bx, VectorOpt::None)?
1218356036faSJeremy L Thompson             .field("qdata", &rq, BasisOpt::None, VectorOpt::Active)?
12199df49d7eSJed Brown             .apply(&x, &mut qdata)?;
12209df49d7eSJed Brown 
12219df49d7eSJed Brown         // Mass operator
12229df49d7eSJed Brown         let qf_mass = ceed.q_function_interior_by_name("MassApply")?;
12239df49d7eSJed Brown         let op_mass = ceed
12249df49d7eSJed Brown             .operator(&qf_mass, QFunctionOpt::None, QFunctionOpt::None)?
12259df49d7eSJed Brown             .field("u", &ru, &bu, VectorOpt::Active)?
1226356036faSJeremy L Thompson             .field("qdata", &rq, BasisOpt::None, &qdata)?
12276f97ff0aSJeremy L Thompson             .field("v", &ru, &bu, VectorOpt::Active)?
12286f97ff0aSJeremy L Thompson             .check()?;
12299df49d7eSJed Brown 
12309df49d7eSJed Brown         v.set_value(0.0)?;
12319df49d7eSJed Brown         op_mass.apply(&u, &mut v)?;
12329df49d7eSJed Brown 
12339df49d7eSJed Brown         // Check
1234e78171edSJeremy L Thompson         let sum: Scalar = v.view()?.iter().sum();
12351ce8139fSJeremy L Thompson         let error: Scalar = (sum - 2.0).abs();
12369df49d7eSJed Brown         assert!(
1237*4b61a5a0SRezgar Shakeri             error < 50.0 * EPSILON,
12381ce8139fSJeremy L Thompson             "Incorrect interval length computed. Expected: 2.0, Found: {}, Error: {:.12e}",
12391ce8139fSJeremy L Thompson             sum,
12401ce8139fSJeremy L Thompson             error
12419df49d7eSJed Brown         );
124289d15d5fSJeremy L Thompson         Ok(())
12439df49d7eSJed Brown     }
12449df49d7eSJed Brown 
12459df49d7eSJed Brown     #[test]
12469df49d7eSJed Brown     fn test_ceed_t501() {
12479df49d7eSJed Brown         assert!(ceed_t501().is_ok());
12489df49d7eSJed Brown     }
12499df49d7eSJed Brown }
12509df49d7eSJed Brown 
12519df49d7eSJed Brown // -----------------------------------------------------------------------------
1252