xref: /libCEED/rust/libceed/src/lib.rs (revision a1cbad85bd2cf51a2e81228682bd7ae54f5d0f99)
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 
17*a1cbad85SJed Brown // Fenced `rust` code blocks included from README.md are executed as part of doctests.
18*a1cbad85SJed 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},
399df49d7eSJed Brown         ElemTopology, EvalMode, MemType, NormType, QuadMode, TransposeMode, CEED_STRIDES_BACKEND,
409df49d7eSJed Brown         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 // -----------------------------------------------------------------------------
589df49d7eSJed Brown // Constants for library
599df49d7eSJed Brown // -----------------------------------------------------------------------------
609df49d7eSJed Brown const MAX_BUFFER_LENGTH: u64 = 4096;
619df49d7eSJed Brown pub const MAX_QFUNCTION_FIELDS: usize = 16;
629df49d7eSJed Brown pub const CEED_STRIDES_BACKEND: [i32; 3] = [0; 3];
639df49d7eSJed Brown 
649df49d7eSJed Brown // -----------------------------------------------------------------------------
659df49d7eSJed Brown // Enums for libCEED
669df49d7eSJed Brown // -----------------------------------------------------------------------------
679df49d7eSJed Brown #[derive(Clone, Copy, PartialEq, Eq)]
689df49d7eSJed Brown /// Many Ceed interfaces take or return pointers to memory.  This enum is used to
699df49d7eSJed Brown /// specify where the memory being provided or requested must reside.
709df49d7eSJed Brown pub enum MemType {
719df49d7eSJed Brown     Host = bind_ceed::CeedMemType_CEED_MEM_HOST as isize,
729df49d7eSJed Brown     Device = bind_ceed::CeedMemType_CEED_MEM_DEVICE as isize,
739df49d7eSJed Brown }
749df49d7eSJed Brown 
759df49d7eSJed Brown #[derive(Clone, Copy, PartialEq, Eq)]
769df49d7eSJed Brown // OwnPointer will not be used by user but is included for internal use
779df49d7eSJed Brown #[allow(dead_code)]
789df49d7eSJed Brown /// Conveys ownership status of arrays passed to Ceed interfaces.
799df49d7eSJed Brown pub(crate) enum CopyMode {
809df49d7eSJed Brown     CopyValues = bind_ceed::CeedCopyMode_CEED_COPY_VALUES as isize,
819df49d7eSJed Brown     UsePointer = bind_ceed::CeedCopyMode_CEED_USE_POINTER as isize,
829df49d7eSJed Brown     OwnPointer = bind_ceed::CeedCopyMode_CEED_OWN_POINTER as isize,
839df49d7eSJed Brown }
849df49d7eSJed Brown 
859df49d7eSJed Brown #[derive(Clone, Copy, PartialEq, Eq)]
869df49d7eSJed Brown /// Denotes type of vector norm to be computed
879df49d7eSJed Brown pub enum NormType {
889df49d7eSJed Brown     One = bind_ceed::CeedNormType_CEED_NORM_1 as isize,
899df49d7eSJed Brown     Two = bind_ceed::CeedNormType_CEED_NORM_2 as isize,
909df49d7eSJed Brown     Max = bind_ceed::CeedNormType_CEED_NORM_MAX as isize,
919df49d7eSJed Brown }
929df49d7eSJed Brown 
939df49d7eSJed Brown #[derive(Clone, Copy, PartialEq, Eq)]
949df49d7eSJed Brown /// Denotes whether a linear transformation or its transpose should be applied
959df49d7eSJed Brown pub enum TransposeMode {
969df49d7eSJed Brown     NoTranspose = bind_ceed::CeedTransposeMode_CEED_NOTRANSPOSE as isize,
979df49d7eSJed Brown     Transpose = bind_ceed::CeedTransposeMode_CEED_TRANSPOSE as isize,
989df49d7eSJed Brown }
999df49d7eSJed Brown 
1009df49d7eSJed Brown #[derive(Clone, Copy, PartialEq, Eq)]
1019df49d7eSJed Brown /// Type of quadrature; also used for location of nodes
1029df49d7eSJed Brown pub enum QuadMode {
1039df49d7eSJed Brown     Gauss = bind_ceed::CeedQuadMode_CEED_GAUSS as isize,
1049df49d7eSJed Brown     GaussLobatto = bind_ceed::CeedQuadMode_CEED_GAUSS_LOBATTO as isize,
1059df49d7eSJed Brown }
1069df49d7eSJed Brown 
1079df49d7eSJed Brown #[derive(Clone, Copy, PartialEq, Eq)]
1089df49d7eSJed Brown /// Type of basis shape to create non-tensor H1 element basis
1099df49d7eSJed Brown pub enum ElemTopology {
1109df49d7eSJed Brown     Line = bind_ceed::CeedElemTopology_CEED_LINE as isize,
1119df49d7eSJed Brown     Triangle = bind_ceed::CeedElemTopology_CEED_TRIANGLE as isize,
1129df49d7eSJed Brown     Quad = bind_ceed::CeedElemTopology_CEED_QUAD as isize,
1139df49d7eSJed Brown     Tet = bind_ceed::CeedElemTopology_CEED_TET as isize,
1149df49d7eSJed Brown     Pyramid = bind_ceed::CeedElemTopology_CEED_PYRAMID as isize,
1159df49d7eSJed Brown     Prism = bind_ceed::CeedElemTopology_CEED_PRISM as isize,
1169df49d7eSJed Brown     Hex = bind_ceed::CeedElemTopology_CEED_HEX as isize,
1179df49d7eSJed Brown }
1189df49d7eSJed Brown 
1199df49d7eSJed Brown #[derive(Clone, Copy, PartialEq, Eq)]
1209df49d7eSJed Brown /// Basis evaluation mode
1219df49d7eSJed Brown pub enum EvalMode {
1229df49d7eSJed Brown     None = bind_ceed::CeedEvalMode_CEED_EVAL_NONE as isize,
1239df49d7eSJed Brown     Interp = bind_ceed::CeedEvalMode_CEED_EVAL_INTERP as isize,
1249df49d7eSJed Brown     Grad = bind_ceed::CeedEvalMode_CEED_EVAL_GRAD as isize,
1259df49d7eSJed Brown     Div = bind_ceed::CeedEvalMode_CEED_EVAL_DIV as isize,
1269df49d7eSJed Brown     Curl = bind_ceed::CeedEvalMode_CEED_EVAL_CURL as isize,
1279df49d7eSJed Brown     Weight = bind_ceed::CeedEvalMode_CEED_EVAL_WEIGHT as isize,
1289df49d7eSJed Brown }
1299df49d7eSJed Brown 
1309df49d7eSJed Brown // -----------------------------------------------------------------------------
1319df49d7eSJed Brown // Ceed error
1329df49d7eSJed Brown // -----------------------------------------------------------------------------
1339df49d7eSJed Brown type Result<T> = std::result::Result<T, CeedError>;
1349df49d7eSJed Brown 
1359df49d7eSJed Brown #[derive(Debug)]
1369df49d7eSJed Brown pub struct CeedError {
1379df49d7eSJed Brown     pub message: String,
1389df49d7eSJed Brown }
1399df49d7eSJed Brown 
1409df49d7eSJed Brown impl fmt::Display for CeedError {
1419df49d7eSJed Brown     fn fmt(&self, f: &mut fmt::Formatter) -> fmt::Result {
1429df49d7eSJed Brown         write!(f, "{}", self.message)
1439df49d7eSJed Brown     }
1449df49d7eSJed Brown }
1459df49d7eSJed Brown 
1469df49d7eSJed Brown // -----------------------------------------------------------------------------
1479df49d7eSJed Brown // Ceed error handler
1489df49d7eSJed Brown // -----------------------------------------------------------------------------
1499df49d7eSJed Brown pub enum CeedErrorHandler {
1509df49d7eSJed Brown     ErrorAbort,
1519df49d7eSJed Brown     ErrorExit,
1529df49d7eSJed Brown     ErrorReturn,
1539df49d7eSJed Brown     ErrorStore,
1549df49d7eSJed Brown }
1559df49d7eSJed Brown 
1569df49d7eSJed Brown // -----------------------------------------------------------------------------
1579df49d7eSJed Brown // Ceed context wrapper
1589df49d7eSJed Brown // -----------------------------------------------------------------------------
1599df49d7eSJed Brown /// A Ceed is a library context representing control of a logical hardware
1609df49d7eSJed Brown /// resource.
1619df49d7eSJed Brown #[derive(Debug)]
1629df49d7eSJed Brown pub struct Ceed {
1639df49d7eSJed Brown     ptr: bind_ceed::Ceed,
1649df49d7eSJed Brown }
1659df49d7eSJed Brown 
1669df49d7eSJed Brown // -----------------------------------------------------------------------------
1679df49d7eSJed Brown // Destructor
1689df49d7eSJed Brown // -----------------------------------------------------------------------------
1699df49d7eSJed Brown impl Drop for Ceed {
1709df49d7eSJed Brown     fn drop(&mut self) {
1719df49d7eSJed Brown         unsafe {
1729df49d7eSJed Brown             bind_ceed::CeedDestroy(&mut self.ptr);
1739df49d7eSJed Brown         }
1749df49d7eSJed Brown     }
1759df49d7eSJed Brown }
1769df49d7eSJed Brown 
1779df49d7eSJed Brown // -----------------------------------------------------------------------------
1789df49d7eSJed Brown // Display
1799df49d7eSJed Brown // -----------------------------------------------------------------------------
1809df49d7eSJed Brown impl fmt::Display for Ceed {
1819df49d7eSJed Brown     /// View a Ceed
1829df49d7eSJed Brown     ///
1839df49d7eSJed Brown     /// ```
1849df49d7eSJed Brown     /// let ceed = libceed::Ceed::init("/cpu/self/ref/serial");
1859df49d7eSJed Brown     /// println!("{}", ceed);
1869df49d7eSJed Brown     /// ```
1879df49d7eSJed Brown     fn fmt(&self, f: &mut fmt::Formatter) -> fmt::Result {
1889df49d7eSJed Brown         let mut ptr = std::ptr::null_mut();
1899df49d7eSJed Brown         let mut sizeloc = crate::MAX_BUFFER_LENGTH;
1909df49d7eSJed Brown         let cstring = unsafe {
1919df49d7eSJed Brown             let file = bind_ceed::open_memstream(&mut ptr, &mut sizeloc);
1929df49d7eSJed Brown             bind_ceed::CeedView(self.ptr, file);
1939df49d7eSJed Brown             bind_ceed::fclose(file);
1949df49d7eSJed Brown             CString::from_raw(ptr)
1959df49d7eSJed Brown         };
1969df49d7eSJed Brown         cstring.to_string_lossy().fmt(f)
1979df49d7eSJed Brown     }
1989df49d7eSJed Brown }
1999df49d7eSJed Brown 
2009df49d7eSJed Brown static REGISTER: Once = Once::new();
2019df49d7eSJed Brown 
2029df49d7eSJed Brown // -----------------------------------------------------------------------------
2039df49d7eSJed Brown // Object constructors
2049df49d7eSJed Brown // -----------------------------------------------------------------------------
2059df49d7eSJed Brown impl Ceed {
2069df49d7eSJed Brown     /// Returns a Ceed context initialized with the specified resource
2079df49d7eSJed Brown     ///
2089df49d7eSJed Brown     /// # arguments
2099df49d7eSJed Brown     ///
2109df49d7eSJed Brown     /// * `resource` - Resource to use, e.g., "/cpu/self"
2119df49d7eSJed Brown     ///
2129df49d7eSJed Brown     /// ```
2139df49d7eSJed Brown     /// let ceed = libceed::Ceed::init("/cpu/self/ref/serial");
2149df49d7eSJed Brown     /// ```
2159df49d7eSJed Brown     pub fn init(resource: &str) -> Self {
2169df49d7eSJed Brown         Ceed::init_with_error_handler(resource, CeedErrorHandler::ErrorStore)
2179df49d7eSJed Brown     }
2189df49d7eSJed Brown 
2199df49d7eSJed Brown     /// Returns a Ceed context initialized with the specified resource
2209df49d7eSJed Brown     ///
2219df49d7eSJed Brown     /// # arguments
2229df49d7eSJed Brown     ///
2239df49d7eSJed Brown     /// * `resource` - Resource to use, e.g., "/cpu/self"
2249df49d7eSJed Brown     ///
2259df49d7eSJed Brown     /// ```
2269df49d7eSJed Brown     /// let ceed = libceed::Ceed::init_with_error_handler(
2279df49d7eSJed Brown     ///     "/cpu/self/ref/serial",
2289df49d7eSJed Brown     ///     libceed::CeedErrorHandler::ErrorAbort,
2299df49d7eSJed Brown     /// );
2309df49d7eSJed Brown     /// ```
2319df49d7eSJed Brown     pub fn init_with_error_handler(resource: &str, handler: CeedErrorHandler) -> Self {
2329df49d7eSJed Brown         REGISTER.call_once(|| unsafe {
2339df49d7eSJed Brown             bind_ceed::CeedRegisterAll();
2349df49d7eSJed Brown             bind_ceed::CeedQFunctionRegisterAll();
2359df49d7eSJed Brown         });
2369df49d7eSJed Brown 
2379df49d7eSJed Brown         // Convert to C string
2389df49d7eSJed Brown         let c_resource = CString::new(resource).expect("CString::new failed");
2399df49d7eSJed Brown 
2409df49d7eSJed Brown         // Get error handler pointer
2419df49d7eSJed Brown         let eh = match handler {
2429df49d7eSJed Brown             CeedErrorHandler::ErrorAbort => bind_ceed::CeedErrorAbort,
2439df49d7eSJed Brown             CeedErrorHandler::ErrorExit => bind_ceed::CeedErrorExit,
2449df49d7eSJed Brown             CeedErrorHandler::ErrorReturn => bind_ceed::CeedErrorReturn,
2459df49d7eSJed Brown             CeedErrorHandler::ErrorStore => bind_ceed::CeedErrorStore,
2469df49d7eSJed Brown         };
2479df49d7eSJed Brown 
2489df49d7eSJed Brown         // Call to libCEED
2499df49d7eSJed Brown         let mut ptr = std::ptr::null_mut();
2509df49d7eSJed Brown         let mut ierr = unsafe { bind_ceed::CeedInit(c_resource.as_ptr() as *const i8, &mut ptr) };
2519df49d7eSJed Brown         if ierr != 0 {
2529df49d7eSJed Brown             panic!("Error initializing backend resource: {}", resource)
2539df49d7eSJed Brown         }
2549df49d7eSJed Brown         ierr = unsafe { bind_ceed::CeedSetErrorHandler(ptr, Some(eh)) };
2559df49d7eSJed Brown         let ceed = Ceed { ptr };
2569df49d7eSJed Brown         ceed.check_error(ierr).unwrap();
2579df49d7eSJed Brown         ceed
2589df49d7eSJed Brown     }
2599df49d7eSJed Brown 
2609df49d7eSJed Brown     /// Default initializer for testing
2619df49d7eSJed Brown     #[doc(hidden)]
2629df49d7eSJed Brown     pub fn default_init() -> Self {
2639df49d7eSJed Brown         // Convert to C string
2649df49d7eSJed Brown         let resource = "/cpu/self/ref/serial";
2659df49d7eSJed Brown         crate::Ceed::init(resource)
2669df49d7eSJed Brown     }
2679df49d7eSJed Brown 
2689df49d7eSJed Brown     /// Internal error checker
2699df49d7eSJed Brown     #[doc(hidden)]
2709df49d7eSJed Brown     fn check_error(&self, ierr: i32) -> Result<i32> {
2719df49d7eSJed Brown         // Return early if code is clean
2729df49d7eSJed Brown         if ierr == bind_ceed::CeedErrorType_CEED_ERROR_SUCCESS {
2739df49d7eSJed Brown             return Ok(ierr);
2749df49d7eSJed Brown         }
2759df49d7eSJed Brown         // Retrieve error message
2769df49d7eSJed Brown         let mut ptr: *const std::os::raw::c_char = std::ptr::null_mut();
2779df49d7eSJed Brown         let c_str = unsafe {
2789df49d7eSJed Brown             bind_ceed::CeedGetErrorMessage(self.ptr, &mut ptr);
2799df49d7eSJed Brown             std::ffi::CStr::from_ptr(ptr)
2809df49d7eSJed Brown         };
2819df49d7eSJed Brown         let message = c_str.to_string_lossy().to_string();
2829df49d7eSJed Brown         // Panic if negative code, otherwise return error
2839df49d7eSJed Brown         if ierr < bind_ceed::CeedErrorType_CEED_ERROR_SUCCESS {
2849df49d7eSJed Brown             panic!("{}", message);
2859df49d7eSJed Brown         }
2869df49d7eSJed Brown         Err(CeedError { message })
2879df49d7eSJed Brown     }
2889df49d7eSJed Brown 
2899df49d7eSJed Brown     /// Returns full resource name for a Ceed context
2909df49d7eSJed Brown     ///
2919df49d7eSJed Brown     /// ```
2929df49d7eSJed Brown     /// let ceed = libceed::Ceed::init("/cpu/self/ref/serial");
2939df49d7eSJed Brown     /// let resource = ceed.resource();
2949df49d7eSJed Brown     ///
2959df49d7eSJed Brown     /// assert_eq!(resource, "/cpu/self/ref/serial".to_string())
2969df49d7eSJed Brown     /// ```
2979df49d7eSJed Brown     pub fn resource(&self) -> String {
2989df49d7eSJed Brown         let mut ptr: *const std::os::raw::c_char = std::ptr::null_mut();
2999df49d7eSJed Brown         let c_str = unsafe {
3009df49d7eSJed Brown             bind_ceed::CeedGetResource(self.ptr, &mut ptr);
3019df49d7eSJed Brown             std::ffi::CStr::from_ptr(ptr)
3029df49d7eSJed Brown         };
3039df49d7eSJed Brown         c_str.to_string_lossy().to_string()
3049df49d7eSJed Brown     }
3059df49d7eSJed Brown 
3069df49d7eSJed Brown     /// Returns a CeedVector of the specified length (does not allocate memory)
3079df49d7eSJed Brown     ///
3089df49d7eSJed Brown     /// # arguments
3099df49d7eSJed Brown     ///
3109df49d7eSJed Brown     /// * `n` - Length of vector
3119df49d7eSJed Brown     ///
3129df49d7eSJed Brown     /// ```
3139df49d7eSJed Brown     /// # use libceed::prelude::*;
3149df49d7eSJed Brown     /// # let ceed = libceed::Ceed::default_init();
3159df49d7eSJed Brown     /// let vec = ceed.vector(10).unwrap();
3169df49d7eSJed Brown     /// ```
3179df49d7eSJed Brown     pub fn vector(&self, n: usize) -> Result<Vector> {
3189df49d7eSJed Brown         Vector::create(self, n)
3199df49d7eSJed Brown     }
3209df49d7eSJed Brown 
3219df49d7eSJed Brown     /// Create a Vector initialized with the data (copied) from a slice
3229df49d7eSJed Brown     ///
3239df49d7eSJed Brown     /// # arguments
3249df49d7eSJed Brown     ///
3259df49d7eSJed Brown     /// * `slice` - Slice containing data
3269df49d7eSJed Brown     ///
3279df49d7eSJed Brown     /// ```
3289df49d7eSJed Brown     /// # use libceed::prelude::*;
3299df49d7eSJed Brown     /// # let ceed = libceed::Ceed::default_init();
3309df49d7eSJed Brown     /// let vec = ceed.vector_from_slice(&[1., 2., 3.]).unwrap();
3319df49d7eSJed Brown     /// assert_eq!(vec.length(), 3);
3329df49d7eSJed Brown     /// ```
3339df49d7eSJed Brown     pub fn vector_from_slice(&self, slice: &[f64]) -> Result<Vector> {
3349df49d7eSJed Brown         Vector::from_slice(self, slice)
3359df49d7eSJed Brown     }
3369df49d7eSJed Brown 
3379df49d7eSJed Brown     /// Returns a ElemRestriction
3389df49d7eSJed Brown     ///
3399df49d7eSJed Brown     /// # arguments
3409df49d7eSJed Brown     ///
3419df49d7eSJed Brown     /// * `nelem`      - Number of elements described in the offsets array
3429df49d7eSJed Brown     /// * `elemsize`   - Size (number of "nodes") per element
3439df49d7eSJed Brown     /// * `ncomp`      - Number of field components per interpolation node (1
3449df49d7eSJed Brown     ///                    for scalar fields)
3459df49d7eSJed Brown     /// * `compstride` - Stride between components for the same Lvector "node".
3469df49d7eSJed Brown     ///                    Data for node `i`, component `j`, element `k` can be
3479df49d7eSJed Brown     ///                    found in the Lvector at index
3489df49d7eSJed Brown     ///                    `offsets[i + k*elemsize] + j*compstride`.
3499df49d7eSJed Brown     /// * `lsize`      - The size of the Lvector. This vector may be larger
3509df49d7eSJed Brown     ///                    than the elements and fields given by this
3519df49d7eSJed Brown     ///                    restriction.
3529df49d7eSJed Brown     /// * `mtype`     - Memory type of the offsets array, see CeedMemType
3539df49d7eSJed Brown     /// * `offsets`    - Array of shape `[nelem, elemsize]`. Row `i` holds the
3549df49d7eSJed Brown     ///                    ordered list of the offsets (into the input CeedVector)
3559df49d7eSJed Brown     ///                    for the unknowns corresponding to element `i`, where
3569df49d7eSJed Brown     ///                    `0 <= i < nelem`. All offsets must be in the range
3579df49d7eSJed Brown     ///                    `[0, lsize - 1]`.
3589df49d7eSJed Brown     ///
3599df49d7eSJed Brown     /// ```
3609df49d7eSJed Brown     /// # use libceed::prelude::*;
3619df49d7eSJed Brown     /// # let ceed = libceed::Ceed::default_init();
3629df49d7eSJed Brown     /// let nelem = 3;
3639df49d7eSJed Brown     /// let mut ind: Vec<i32> = vec![0; 2 * nelem];
3649df49d7eSJed Brown     /// for i in 0..nelem {
3659df49d7eSJed Brown     ///     ind[2 * i + 0] = i as i32;
3669df49d7eSJed Brown     ///     ind[2 * i + 1] = (i + 1) as i32;
3679df49d7eSJed Brown     /// }
3689df49d7eSJed Brown     /// let r = ceed
3699df49d7eSJed Brown     ///     .elem_restriction(nelem, 2, 1, 1, nelem + 1, MemType::Host, &ind)
3709df49d7eSJed Brown     ///     .unwrap();
3719df49d7eSJed Brown     /// ```
3729df49d7eSJed Brown     pub fn elem_restriction(
3739df49d7eSJed Brown         &self,
3749df49d7eSJed Brown         nelem: usize,
3759df49d7eSJed Brown         elemsize: usize,
3769df49d7eSJed Brown         ncomp: usize,
3779df49d7eSJed Brown         compstride: usize,
3789df49d7eSJed Brown         lsize: usize,
3799df49d7eSJed Brown         mtype: MemType,
3809df49d7eSJed Brown         offsets: &[i32],
3819df49d7eSJed Brown     ) -> Result<ElemRestriction> {
3829df49d7eSJed Brown         ElemRestriction::create(
3839df49d7eSJed Brown             self, nelem, elemsize, ncomp, compstride, lsize, mtype, offsets,
3849df49d7eSJed Brown         )
3859df49d7eSJed Brown     }
3869df49d7eSJed Brown 
3879df49d7eSJed Brown     /// Returns a ElemRestriction
3889df49d7eSJed Brown     ///
3899df49d7eSJed Brown     /// # arguments
3909df49d7eSJed Brown     ///
3919df49d7eSJed Brown     /// * `nelem`      - Number of elements described in the offsets array
3929df49d7eSJed Brown     /// * `elemsize`   - Size (number of "nodes") per element
3939df49d7eSJed Brown     /// * `ncomp`      - Number of field components per interpolation node (1
3949df49d7eSJed Brown     ///                    for scalar fields)
3959df49d7eSJed Brown     /// * `compstride` - Stride between components for the same Lvector "node".
3969df49d7eSJed Brown     ///                    Data for node `i`, component `j`, element `k` can be
3979df49d7eSJed Brown     ///                    found in the Lvector at index
3989df49d7eSJed Brown     ///                    `offsets[i + k*elemsize] + j*compstride`.
3999df49d7eSJed Brown     /// * `lsize`      - The size of the Lvector. This vector may be larger
4009df49d7eSJed Brown     ///   than the elements and fields given by this restriction.
4019df49d7eSJed Brown     /// * `strides`   - Array for strides between `[nodes, components, elements]`.
4029df49d7eSJed Brown     ///                   Data for node `i`, component `j`, element `k` can be
4039df49d7eSJed Brown     ///                   found in the Lvector at index
4049df49d7eSJed Brown     ///                   `i*strides[0] + j*strides[1] + k*strides[2]`.
4059df49d7eSJed Brown     ///                   CEED_STRIDES_BACKEND may be used with vectors created
4069df49d7eSJed Brown     ///                   by a Ceed backend.
4079df49d7eSJed Brown     ///
4089df49d7eSJed Brown     /// ```
4099df49d7eSJed Brown     /// # use libceed::prelude::*;
4109df49d7eSJed Brown     /// # let ceed = libceed::Ceed::default_init();
4119df49d7eSJed Brown     /// let nelem = 3;
4129df49d7eSJed Brown     /// let strides: [i32; 3] = [1, 2, 2];
4139df49d7eSJed Brown     /// let r = ceed
4149df49d7eSJed Brown     ///     .strided_elem_restriction(nelem, 2, 1, nelem * 2, strides)
4159df49d7eSJed Brown     ///     .unwrap();
4169df49d7eSJed Brown     /// ```
4179df49d7eSJed Brown     pub fn strided_elem_restriction(
4189df49d7eSJed Brown         &self,
4199df49d7eSJed Brown         nelem: usize,
4209df49d7eSJed Brown         elemsize: usize,
4219df49d7eSJed Brown         ncomp: usize,
4229df49d7eSJed Brown         lsize: usize,
4239df49d7eSJed Brown         strides: [i32; 3],
4249df49d7eSJed Brown     ) -> Result<ElemRestriction> {
4259df49d7eSJed Brown         ElemRestriction::create_strided(self, nelem, elemsize, ncomp, lsize, strides)
4269df49d7eSJed Brown     }
4279df49d7eSJed Brown 
4289df49d7eSJed Brown     /// Returns a tensor-product basis
4299df49d7eSJed Brown     ///
4309df49d7eSJed Brown     /// # arguments
4319df49d7eSJed Brown     ///
4329df49d7eSJed Brown     /// * `dim`       - Topological dimension of element
4339df49d7eSJed Brown     /// * `ncomp`     - Number of field components (1 for scalar fields)
4349df49d7eSJed Brown     /// * `P1d`       - Number of Gauss-Lobatto nodes in one dimension.  The
4359df49d7eSJed Brown     ///                   polynomial degree of the resulting `Q_k` element is
4369df49d7eSJed Brown     ///                   `k=P-1`.
4379df49d7eSJed Brown     /// * `Q1d`       - Number of quadrature points in one dimension
4389df49d7eSJed Brown     /// * `interp1d`  - Row-major `(Q1d * P1d)` matrix expressing the values of
4399df49d7eSJed Brown     ///                   nodal basis functions at quadrature points
4409df49d7eSJed Brown     /// * `grad1d`    - Row-major `(Q1d * P1d)` matrix expressing derivatives of
4419df49d7eSJed Brown     ///                   nodal basis functions at quadrature points
4429df49d7eSJed Brown     /// * `qref1d`    - Array of length `Q1d` holding the locations of quadrature
4439df49d7eSJed Brown     ///                   points on the 1D reference element `[-1, 1]`
4449df49d7eSJed Brown     /// * `qweight1d` - Array of length `Q1d` holding the quadrature weights on
4459df49d7eSJed Brown     ///                   the reference element
4469df49d7eSJed Brown     ///
4479df49d7eSJed Brown     /// ```
4489df49d7eSJed Brown     /// # use libceed::prelude::*;
4499df49d7eSJed Brown     /// # let ceed = libceed::Ceed::default_init();
4509df49d7eSJed Brown     /// let interp1d  = [ 0.62994317,  0.47255875, -0.14950343,  0.04700152,
4519df49d7eSJed Brown     ///                  -0.07069480,  0.97297619,  0.13253993, -0.03482132,
4529df49d7eSJed Brown     ///                  -0.03482132,  0.13253993,  0.97297619, -0.07069480,
4539df49d7eSJed Brown     ///                   0.04700152, -0.14950343,  0.47255875,  0.62994317];
4549df49d7eSJed Brown     /// let grad1d    = [-2.34183742,  2.78794489, -0.63510411,  0.18899664,
4559df49d7eSJed Brown     ///                  -0.51670214, -0.48795249,  1.33790510, -0.33325047,
4569df49d7eSJed Brown     //                    0.33325047, -1.33790510,  0.48795249,  0.51670214,
4579df49d7eSJed Brown     ///                  -0.18899664,  0.63510411, -2.78794489,  2.34183742];
4589df49d7eSJed Brown     /// let qref1d    = [-0.86113631, -0.33998104,  0.33998104,  0.86113631];
4599df49d7eSJed Brown     /// let qweight1d = [ 0.34785485,  0.65214515,  0.65214515,  0.34785485];
4609df49d7eSJed Brown     /// let b = ceed.
4619df49d7eSJed Brown     /// basis_tensor_H1(2, 1, 4, 4, &interp1d, &grad1d, &qref1d, &qweight1d).unwrap();
4629df49d7eSJed Brown     /// ```
4639df49d7eSJed Brown     pub fn basis_tensor_H1(
4649df49d7eSJed Brown         &self,
4659df49d7eSJed Brown         dim: usize,
4669df49d7eSJed Brown         ncomp: usize,
4679df49d7eSJed Brown         P1d: usize,
4689df49d7eSJed Brown         Q1d: usize,
4699df49d7eSJed Brown         interp1d: &[f64],
4709df49d7eSJed Brown         grad1d: &[f64],
4719df49d7eSJed Brown         qref1d: &[f64],
4729df49d7eSJed Brown         qweight1d: &[f64],
4739df49d7eSJed Brown     ) -> Result<Basis> {
4749df49d7eSJed Brown         Basis::create_tensor_H1(
4759df49d7eSJed Brown             self, dim, ncomp, P1d, Q1d, interp1d, grad1d, qref1d, qweight1d,
4769df49d7eSJed Brown         )
4779df49d7eSJed Brown     }
4789df49d7eSJed Brown 
4799df49d7eSJed Brown     /// Returns a tensor-product Lagrange basis
4809df49d7eSJed Brown     ///
4819df49d7eSJed Brown     /// # arguments
4829df49d7eSJed Brown     ///
4839df49d7eSJed Brown     /// * `dim`   - Topological dimension of element
4849df49d7eSJed Brown     /// * `ncomp` - Number of field components (1 for scalar fields)
4859df49d7eSJed Brown     /// * `P`     - Number of Gauss-Lobatto nodes in one dimension.  The
4869df49d7eSJed Brown     ///               polynomial degree of the resulting `Q_k` element is `k=P-1`.
4879df49d7eSJed Brown     /// * `Q`     - Number of quadrature points in one dimension
4889df49d7eSJed Brown     /// * `qmode` - Distribution of the `Q` quadrature points (affects order of
4899df49d7eSJed Brown     ///               accuracy for the quadrature)
4909df49d7eSJed Brown     ///
4919df49d7eSJed Brown     /// ```
4929df49d7eSJed Brown     /// # use libceed::prelude::*;
4939df49d7eSJed Brown     /// # let ceed = libceed::Ceed::default_init();
4949df49d7eSJed Brown     /// let b = ceed
4959df49d7eSJed Brown     ///     .basis_tensor_H1_Lagrange(2, 1, 3, 4, QuadMode::Gauss)
4969df49d7eSJed Brown     ///     .unwrap();
4979df49d7eSJed Brown     /// ```
4989df49d7eSJed Brown     pub fn basis_tensor_H1_Lagrange(
4999df49d7eSJed Brown         &self,
5009df49d7eSJed Brown         dim: usize,
5019df49d7eSJed Brown         ncomp: usize,
5029df49d7eSJed Brown         P: usize,
5039df49d7eSJed Brown         Q: usize,
5049df49d7eSJed Brown         qmode: QuadMode,
5059df49d7eSJed Brown     ) -> Result<Basis> {
5069df49d7eSJed Brown         Basis::create_tensor_H1_Lagrange(self, dim, ncomp, P, Q, qmode)
5079df49d7eSJed Brown     }
5089df49d7eSJed Brown 
5099df49d7eSJed Brown     /// Returns a tensor-product basis
5109df49d7eSJed Brown     ///
5119df49d7eSJed Brown     /// # arguments
5129df49d7eSJed Brown     ///
5139df49d7eSJed Brown     /// * `topo`    - Topology of element, e.g. hypercube, simplex, ect
5149df49d7eSJed Brown     /// * `ncomp`   - Number of field components (1 for scalar fields)
5159df49d7eSJed Brown     /// * `nnodes`  - Total number of nodes
5169df49d7eSJed Brown     /// * `nqpts`   - Total number of quadrature points
5179df49d7eSJed Brown     /// * `interp`  - Row-major `(nqpts * nnodes)` matrix expressing the values of
5189df49d7eSJed Brown     ///                 nodal basis functions at quadrature points
5199df49d7eSJed Brown     /// * `grad`    - Row-major `(nqpts * dim * nnodes)` matrix expressing
5209df49d7eSJed Brown     ///                 derivatives of nodal basis functions at quadrature points
5219df49d7eSJed Brown     /// * `qref`    - Array of length `nqpts` holding the locations of quadrature
5229df49d7eSJed Brown     ///                 points on the reference element `[-1, 1]`
5239df49d7eSJed Brown     /// * `qweight` - Array of length `nqpts` holding the quadrature weights on
5249df49d7eSJed Brown     ///                 the reference element
5259df49d7eSJed Brown     ///
5269df49d7eSJed Brown     /// ```
5279df49d7eSJed Brown     /// # use libceed::prelude::*;
5289df49d7eSJed Brown     /// # let ceed = libceed::Ceed::default_init();
5299df49d7eSJed Brown     /// let interp = [
5309df49d7eSJed Brown     ///     0.12000000,
5319df49d7eSJed Brown     ///     0.48000000,
5329df49d7eSJed Brown     ///     -0.12000000,
5339df49d7eSJed Brown     ///     0.48000000,
5349df49d7eSJed Brown     ///     0.16000000,
5359df49d7eSJed Brown     ///     -0.12000000,
5369df49d7eSJed Brown     ///     -0.12000000,
5379df49d7eSJed Brown     ///     0.48000000,
5389df49d7eSJed Brown     ///     0.12000000,
5399df49d7eSJed Brown     ///     0.16000000,
5409df49d7eSJed Brown     ///     0.48000000,
5419df49d7eSJed Brown     ///     -0.12000000,
5429df49d7eSJed Brown     ///     -0.11111111,
5439df49d7eSJed Brown     ///     0.44444444,
5449df49d7eSJed Brown     ///     -0.11111111,
5459df49d7eSJed Brown     ///     0.44444444,
5469df49d7eSJed Brown     ///     0.44444444,
5479df49d7eSJed Brown     ///     -0.11111111,
5489df49d7eSJed Brown     ///     -0.12000000,
5499df49d7eSJed Brown     ///     0.16000000,
5509df49d7eSJed Brown     ///     -0.12000000,
5519df49d7eSJed Brown     ///     0.48000000,
5529df49d7eSJed Brown     ///     0.48000000,
5539df49d7eSJed Brown     ///     0.12000000,
5549df49d7eSJed Brown     /// ];
5559df49d7eSJed Brown     /// let grad = [
5569df49d7eSJed Brown     ///     -1.40000000,
5579df49d7eSJed Brown     ///     1.60000000,
5589df49d7eSJed Brown     ///     -0.20000000,
5599df49d7eSJed Brown     ///     -0.80000000,
5609df49d7eSJed Brown     ///     0.80000000,
5619df49d7eSJed Brown     ///     0.00000000,
5629df49d7eSJed Brown     ///     0.20000000,
5639df49d7eSJed Brown     ///     -1.60000000,
5649df49d7eSJed Brown     ///     1.40000000,
5659df49d7eSJed Brown     ///     -0.80000000,
5669df49d7eSJed Brown     ///     0.80000000,
5679df49d7eSJed Brown     ///     0.00000000,
5689df49d7eSJed Brown     ///     -0.33333333,
5699df49d7eSJed Brown     ///     0.00000000,
5709df49d7eSJed Brown     ///     0.33333333,
5719df49d7eSJed Brown     ///     -1.33333333,
5729df49d7eSJed Brown     ///     1.33333333,
5739df49d7eSJed Brown     ///     0.00000000,
5749df49d7eSJed Brown     ///     0.20000000,
5759df49d7eSJed Brown     ///     0.00000000,
5769df49d7eSJed Brown     ///     -0.20000000,
5779df49d7eSJed Brown     ///     -2.40000000,
5789df49d7eSJed Brown     ///     2.40000000,
5799df49d7eSJed Brown     ///     0.00000000,
5809df49d7eSJed Brown     ///     -1.40000000,
5819df49d7eSJed Brown     ///     -0.80000000,
5829df49d7eSJed Brown     ///     0.00000000,
5839df49d7eSJed Brown     ///     1.60000000,
5849df49d7eSJed Brown     ///     0.80000000,
5859df49d7eSJed Brown     ///     -0.20000000,
5869df49d7eSJed Brown     ///     0.20000000,
5879df49d7eSJed Brown     ///     -2.40000000,
5889df49d7eSJed Brown     ///     0.00000000,
5899df49d7eSJed Brown     ///     0.00000000,
5909df49d7eSJed Brown     ///     2.40000000,
5919df49d7eSJed Brown     ///     -0.20000000,
5929df49d7eSJed Brown     ///     -0.33333333,
5939df49d7eSJed Brown     ///     -1.33333333,
5949df49d7eSJed Brown     ///     0.00000000,
5959df49d7eSJed Brown     ///     0.00000000,
5969df49d7eSJed Brown     ///     1.33333333,
5979df49d7eSJed Brown     ///     0.33333333,
5989df49d7eSJed Brown     ///     0.20000000,
5999df49d7eSJed Brown     ///     -0.80000000,
6009df49d7eSJed Brown     ///     0.00000000,
6019df49d7eSJed Brown     ///     -1.60000000,
6029df49d7eSJed Brown     ///     0.80000000,
6039df49d7eSJed Brown     ///     1.40000000,
6049df49d7eSJed Brown     /// ];
6059df49d7eSJed Brown     /// let qref = [
6069df49d7eSJed Brown     ///     0.20000000, 0.60000000, 0.33333333, 0.20000000, 0.20000000, 0.20000000, 0.33333333,
6079df49d7eSJed Brown     ///     0.60000000,
6089df49d7eSJed Brown     /// ];
6099df49d7eSJed Brown     /// let qweight = [0.26041667, 0.26041667, -0.28125000, 0.26041667];
6109df49d7eSJed Brown     /// let b = ceed
6119df49d7eSJed Brown     ///     .basis_H1(
6129df49d7eSJed Brown     ///         ElemTopology::Triangle,
6139df49d7eSJed Brown     ///         1,
6149df49d7eSJed Brown     ///         6,
6159df49d7eSJed Brown     ///         4,
6169df49d7eSJed Brown     ///         &interp,
6179df49d7eSJed Brown     ///         &grad,
6189df49d7eSJed Brown     ///         &qref,
6199df49d7eSJed Brown     ///         &qweight,
6209df49d7eSJed Brown     ///     )
6219df49d7eSJed Brown     ///     .unwrap();
6229df49d7eSJed Brown     /// ```
6239df49d7eSJed Brown     pub fn basis_H1(
6249df49d7eSJed Brown         &self,
6259df49d7eSJed Brown         topo: ElemTopology,
6269df49d7eSJed Brown         ncomp: usize,
6279df49d7eSJed Brown         nnodes: usize,
6289df49d7eSJed Brown         nqpts: usize,
6299df49d7eSJed Brown         interp: &[f64],
6309df49d7eSJed Brown         grad: &[f64],
6319df49d7eSJed Brown         qref: &[f64],
6329df49d7eSJed Brown         qweight: &[f64],
6339df49d7eSJed Brown     ) -> Result<Basis> {
6349df49d7eSJed Brown         Basis::create_H1(
6359df49d7eSJed Brown             self, topo, ncomp, nnodes, nqpts, interp, grad, qref, qweight,
6369df49d7eSJed Brown         )
6379df49d7eSJed Brown     }
6389df49d7eSJed Brown 
6399df49d7eSJed Brown     /// Returns a CeedQFunction for evaluating interior (volumetric) terms
6409df49d7eSJed Brown     ///
6419df49d7eSJed Brown     /// # arguments
6429df49d7eSJed Brown     ///
6439df49d7eSJed Brown     /// * `vlength` - Vector length. Caller must ensure that number of
6449df49d7eSJed Brown     ///                 quadrature points is a multiple of vlength.
6459df49d7eSJed Brown     /// * `f`       - Boxed closure to evaluate action at quadrature points.
6469df49d7eSJed Brown     ///
6479df49d7eSJed Brown     /// ```
6489df49d7eSJed Brown     /// # use libceed::prelude::*;
6499df49d7eSJed Brown     /// # let ceed = libceed::Ceed::default_init();
6509df49d7eSJed Brown     /// let mut user_f = |[u, weights, ..]: QFunctionInputs, [v, ..]: QFunctionOutputs| {
6519df49d7eSJed Brown     ///     // Iterate over quadrature points
6529df49d7eSJed Brown     ///     v.iter_mut()
6539df49d7eSJed Brown     ///         .zip(u.iter().zip(weights.iter()))
6549df49d7eSJed Brown     ///         .for_each(|(v, (u, w))| *v = u * w);
6559df49d7eSJed Brown     ///
6569df49d7eSJed Brown     ///     // Return clean error code
6579df49d7eSJed Brown     ///     0
6589df49d7eSJed Brown     /// };
6599df49d7eSJed Brown     ///
6609df49d7eSJed Brown     /// let qf = ceed.q_function_interior(1, Box::new(user_f)).unwrap();
6619df49d7eSJed Brown     /// ```
6629df49d7eSJed Brown     pub fn q_function_interior(
6639df49d7eSJed Brown         &self,
6649df49d7eSJed Brown         vlength: usize,
6659df49d7eSJed Brown         f: Box<qfunction::QFunctionUserClosure>,
6669df49d7eSJed Brown     ) -> Result<QFunction> {
6679df49d7eSJed Brown         QFunction::create(self, vlength, f)
6689df49d7eSJed Brown     }
6699df49d7eSJed Brown 
6709df49d7eSJed Brown     /// Returns a CeedQFunction for evaluating interior (volumetric) terms
6719df49d7eSJed Brown     /// created by name
6729df49d7eSJed Brown     ///
6739df49d7eSJed Brown     /// ```
6749df49d7eSJed Brown     /// # use libceed::prelude::*;
6759df49d7eSJed Brown     /// # let ceed = libceed::Ceed::default_init();
6769df49d7eSJed Brown     /// let qf = ceed.q_function_interior_by_name("Mass1DBuild").unwrap();
6779df49d7eSJed Brown     /// ```
6789df49d7eSJed Brown     pub fn q_function_interior_by_name(&self, name: &str) -> Result<QFunctionByName> {
6799df49d7eSJed Brown         QFunctionByName::create(self, name)
6809df49d7eSJed Brown     }
6819df49d7eSJed Brown 
6829df49d7eSJed Brown     /// Returns a Operator and associate a QFunction. A Basis and
6839df49d7eSJed Brown     /// ElemRestriction can be   associated with QFunction fields with
6849df49d7eSJed Brown     /// set_field().
6859df49d7eSJed Brown     ///
6869df49d7eSJed Brown     /// * `qf`   - QFunction defining the action of the operator at quadrature
6879df49d7eSJed Brown     ///              points
6889df49d7eSJed Brown     /// * `dqf`  - QFunction defining the action of the Jacobian of the qf (or
6899df49d7eSJed Brown     ///              qfunction_none)
6909df49d7eSJed Brown     /// * `dqfT` - QFunction defining the action of the transpose of the
6919df49d7eSJed Brown     ///              Jacobian of the qf (or qfunction_none)
6929df49d7eSJed Brown     ///
6939df49d7eSJed Brown     /// ```
6949df49d7eSJed Brown     /// # use libceed::prelude::*;
6959df49d7eSJed Brown     /// # let ceed = libceed::Ceed::default_init();
6969df49d7eSJed Brown     /// let qf = ceed.q_function_interior_by_name("Mass1DBuild").unwrap();
6979df49d7eSJed Brown     /// let op = ceed
6989df49d7eSJed Brown     ///     .operator(&qf, QFunctionOpt::None, QFunctionOpt::None)
6999df49d7eSJed Brown     ///     .unwrap();
7009df49d7eSJed Brown     /// ```
7019df49d7eSJed Brown     pub fn operator<'b>(
7029df49d7eSJed Brown         &self,
7039df49d7eSJed Brown         qf: impl Into<QFunctionOpt<'b>>,
7049df49d7eSJed Brown         dqf: impl Into<QFunctionOpt<'b>>,
7059df49d7eSJed Brown         dqfT: impl Into<QFunctionOpt<'b>>,
7069df49d7eSJed Brown     ) -> Result<Operator> {
7079df49d7eSJed Brown         Operator::create(self, qf, dqf, dqfT)
7089df49d7eSJed Brown     }
7099df49d7eSJed Brown 
7109df49d7eSJed Brown     /// Returns an Operator that composes the action of several Operators
7119df49d7eSJed Brown     ///
7129df49d7eSJed Brown     /// ```
7139df49d7eSJed Brown     /// # use libceed::prelude::*;
7149df49d7eSJed Brown     /// # let ceed = libceed::Ceed::default_init();
7159df49d7eSJed Brown     /// let op = ceed.composite_operator().unwrap();
7169df49d7eSJed Brown     /// ```
7179df49d7eSJed Brown     pub fn composite_operator(&self) -> Result<CompositeOperator> {
7189df49d7eSJed Brown         CompositeOperator::create(self)
7199df49d7eSJed Brown     }
7209df49d7eSJed Brown }
7219df49d7eSJed Brown 
7229df49d7eSJed Brown // -----------------------------------------------------------------------------
7239df49d7eSJed Brown // Tests
7249df49d7eSJed Brown // -----------------------------------------------------------------------------
7259df49d7eSJed Brown #[cfg(test)]
7269df49d7eSJed Brown mod tests {
7279df49d7eSJed Brown     use super::*;
7289df49d7eSJed Brown 
7299df49d7eSJed Brown     fn ceed_t501() -> Result<i32> {
7309df49d7eSJed Brown         let resource = "/cpu/self/ref/blocked";
7319df49d7eSJed Brown         let ceed = Ceed::init(resource);
7329df49d7eSJed Brown         let nelem = 4;
7339df49d7eSJed Brown         let p = 3;
7349df49d7eSJed Brown         let q = 4;
7359df49d7eSJed Brown         let ndofs = p * nelem - nelem + 1;
7369df49d7eSJed Brown 
7379df49d7eSJed Brown         // Vectors
7389df49d7eSJed Brown         let x = ceed.vector_from_slice(&[-1., -0.5, 0.0, 0.5, 1.0])?;
7399df49d7eSJed Brown         let mut qdata = ceed.vector(nelem * q)?;
7409df49d7eSJed Brown         qdata.set_value(0.0)?;
7419df49d7eSJed Brown         let mut u = ceed.vector(ndofs)?;
7429df49d7eSJed Brown         u.set_value(1.0)?;
7439df49d7eSJed Brown         let mut v = ceed.vector(ndofs)?;
7449df49d7eSJed Brown         v.set_value(0.0)?;
7459df49d7eSJed Brown 
7469df49d7eSJed Brown         // Restrictions
7479df49d7eSJed Brown         let mut indx: Vec<i32> = vec![0; 2 * nelem];
7489df49d7eSJed Brown         for i in 0..nelem {
7499df49d7eSJed Brown             indx[2 * i + 0] = i as i32;
7509df49d7eSJed Brown             indx[2 * i + 1] = (i + 1) as i32;
7519df49d7eSJed Brown         }
7529df49d7eSJed Brown         let rx = ceed.elem_restriction(nelem, 2, 1, 1, nelem + 1, MemType::Host, &indx)?;
7539df49d7eSJed Brown         let mut indu: Vec<i32> = vec![0; p * nelem];
7549df49d7eSJed Brown         for i in 0..nelem {
7559df49d7eSJed Brown             indu[p * i + 0] = i as i32;
7569df49d7eSJed Brown             indu[p * i + 1] = (i + 1) as i32;
7579df49d7eSJed Brown             indu[p * i + 2] = (i + 2) as i32;
7589df49d7eSJed Brown         }
7599df49d7eSJed Brown         let ru = ceed.elem_restriction(nelem, 3, 1, 1, ndofs, MemType::Host, &indu)?;
7609df49d7eSJed Brown         let strides: [i32; 3] = [1, q as i32, q as i32];
7619df49d7eSJed Brown         let rq = ceed.strided_elem_restriction(nelem, q, 1, q * nelem, strides)?;
7629df49d7eSJed Brown 
7639df49d7eSJed Brown         // Bases
7649df49d7eSJed Brown         let bx = ceed.basis_tensor_H1_Lagrange(1, 1, 2, q, QuadMode::Gauss)?;
7659df49d7eSJed Brown         let bu = ceed.basis_tensor_H1_Lagrange(1, 1, p, q, QuadMode::Gauss)?;
7669df49d7eSJed Brown 
7679df49d7eSJed Brown         // Build quadrature data
7689df49d7eSJed Brown         let qf_build = ceed.q_function_interior_by_name("Mass1DBuild")?;
7699df49d7eSJed Brown         ceed.operator(&qf_build, QFunctionOpt::None, QFunctionOpt::None)?
7709df49d7eSJed Brown             .field("dx", &rx, &bx, VectorOpt::Active)?
7719df49d7eSJed Brown             .field("weights", ElemRestrictionOpt::None, &bx, VectorOpt::None)?
7729df49d7eSJed Brown             .field("qdata", &rq, BasisOpt::Collocated, VectorOpt::Active)?
7739df49d7eSJed Brown             .apply(&x, &mut qdata)?;
7749df49d7eSJed Brown 
7759df49d7eSJed Brown         // Mass operator
7769df49d7eSJed Brown         let qf_mass = ceed.q_function_interior_by_name("MassApply")?;
7779df49d7eSJed Brown         let op_mass = ceed
7789df49d7eSJed Brown             .operator(&qf_mass, QFunctionOpt::None, QFunctionOpt::None)?
7799df49d7eSJed Brown             .field("u", &ru, &bu, VectorOpt::Active)?
7809df49d7eSJed Brown             .field("qdata", &rq, BasisOpt::Collocated, &qdata)?
7819df49d7eSJed Brown             .field("v", &ru, &bu, VectorOpt::Active)?;
7829df49d7eSJed Brown 
7839df49d7eSJed Brown         v.set_value(0.0)?;
7849df49d7eSJed Brown         op_mass.apply(&u, &mut v)?;
7859df49d7eSJed Brown 
7869df49d7eSJed Brown         // Check
7879df49d7eSJed Brown         let sum: f64 = v.view().iter().sum();
7889df49d7eSJed Brown         assert!(
7899df49d7eSJed Brown             (sum - 2.0).abs() < 1e-15,
7909df49d7eSJed Brown             "Incorrect interval length computed"
7919df49d7eSJed Brown         );
7929df49d7eSJed Brown         Ok(0)
7939df49d7eSJed Brown     }
7949df49d7eSJed Brown 
7959df49d7eSJed Brown     #[test]
7969df49d7eSJed Brown     fn test_ceed_t501() {
7979df49d7eSJed Brown         assert!(ceed_t501().is_ok());
7989df49d7eSJed Brown     }
7999df49d7eSJed Brown }
8009df49d7eSJed Brown 
8019df49d7eSJed Brown // -----------------------------------------------------------------------------
802