xref: /libCEED/rust/libceed/src/basis.rs (revision 9ba83ac0e4b1fca39d6fa6737a318a9f0cbc172d)
1*9ba83ac0SJeremy L Thompson // Copyright (c) 2017-2026, 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 
89df49d7eSJed Brown //! A Ceed Basis defines the discrete finite element basis and associated
99df49d7eSJed Brown //! quadrature rule.
109df49d7eSJed Brown 
11eb07d68fSJeremy L Thompson use crate::{prelude::*, vector::Vector, EvalMode, TransposeMode};
129df49d7eSJed Brown 
139df49d7eSJed Brown // -----------------------------------------------------------------------------
147ed177dbSJed Brown // Basis option
159df49d7eSJed Brown // -----------------------------------------------------------------------------
16c68be7a2SJeremy L Thompson #[derive(Debug)]
179df49d7eSJed Brown pub enum BasisOpt<'a> {
189df49d7eSJed Brown     Some(&'a Basis<'a>),
19356036faSJeremy L Thompson     None,
209df49d7eSJed Brown }
219df49d7eSJed Brown /// Construct a BasisOpt reference from a Basis reference
229df49d7eSJed Brown impl<'a> From<&'a Basis<'_>> for BasisOpt<'a> {
239df49d7eSJed Brown     fn from(basis: &'a Basis) -> Self {
24356036faSJeremy L Thompson         debug_assert!(basis.ptr != unsafe { bind_ceed::CEED_BASIS_NONE });
259df49d7eSJed Brown         Self::Some(basis)
269df49d7eSJed Brown     }
279df49d7eSJed Brown }
289df49d7eSJed Brown impl<'a> BasisOpt<'a> {
299df49d7eSJed Brown     /// Transform a Rust libCEED BasisOpt into C libCEED CeedBasis
3078c2cefaSJeremy L Thompson     pub(crate) fn to_raw(&self) -> bind_ceed::CeedBasis {
319df49d7eSJed Brown         match self {
329df49d7eSJed Brown             Self::Some(basis) => basis.ptr,
33356036faSJeremy L Thompson             Self::None => unsafe { bind_ceed::CEED_BASIS_NONE },
349df49d7eSJed Brown         }
359df49d7eSJed Brown     }
36e03682afSJeremy L Thompson 
37e03682afSJeremy L Thompson     /// Check if a BasisOpt is Some
38e03682afSJeremy L Thompson     ///
39e03682afSJeremy L Thompson     /// ```
40eb07d68fSJeremy L Thompson     /// # use libceed::{prelude::*, BasisOpt, QuadMode};
414d27c890SJeremy L Thompson     /// # fn main() -> libceed::Result<()> {
42e03682afSJeremy L Thompson     /// # let ceed = libceed::Ceed::default_init();
43e03682afSJeremy L Thompson     /// let b = ceed.basis_tensor_H1_Lagrange(1, 2, 3, 4, QuadMode::Gauss)?;
44e03682afSJeremy L Thompson     /// let b_opt = BasisOpt::from(&b);
45e03682afSJeremy L Thompson     /// assert!(b_opt.is_some(), "Incorrect BasisOpt");
46e03682afSJeremy L Thompson     ///
47356036faSJeremy L Thompson     /// let b_opt = BasisOpt::None;
48e03682afSJeremy L Thompson     /// assert!(!b_opt.is_some(), "Incorrect BasisOpt");
49e03682afSJeremy L Thompson     /// # Ok(())
50e03682afSJeremy L Thompson     /// # }
51e03682afSJeremy L Thompson     /// ```
52e03682afSJeremy L Thompson     pub fn is_some(&self) -> bool {
53e03682afSJeremy L Thompson         match self {
54e03682afSJeremy L Thompson             Self::Some(_) => true,
55356036faSJeremy L Thompson             Self::None => false,
56e03682afSJeremy L Thompson         }
57e03682afSJeremy L Thompson     }
58e03682afSJeremy L Thompson 
59356036faSJeremy L Thompson     /// Check if a BasisOpt is None
60e03682afSJeremy L Thompson     ///
61e03682afSJeremy L Thompson     /// ```
62eb07d68fSJeremy L Thompson     /// # use libceed::{prelude::*, BasisOpt, QuadMode};
634d27c890SJeremy L Thompson     /// # fn main() -> libceed::Result<()> {
64e03682afSJeremy L Thompson     /// # let ceed = libceed::Ceed::default_init();
65e03682afSJeremy L Thompson     /// let b = ceed.basis_tensor_H1_Lagrange(1, 2, 3, 4, QuadMode::Gauss)?;
66e03682afSJeremy L Thompson     /// let b_opt = BasisOpt::from(&b);
67356036faSJeremy L Thompson     /// assert!(!b_opt.is_none(), "Incorrect BasisOpt");
68e03682afSJeremy L Thompson     ///
69356036faSJeremy L Thompson     /// let b_opt = BasisOpt::None;
70356036faSJeremy L Thompson     /// assert!(b_opt.is_none(), "Incorrect BasisOpt");
71e03682afSJeremy L Thompson     /// # Ok(())
72e03682afSJeremy L Thompson     /// # }
73e03682afSJeremy L Thompson     /// ```
74356036faSJeremy L Thompson     pub fn is_none(&self) -> bool {
75e03682afSJeremy L Thompson         match self {
76e03682afSJeremy L Thompson             Self::Some(_) => false,
77356036faSJeremy L Thompson             Self::None => true,
78e03682afSJeremy L Thompson         }
79e03682afSJeremy L Thompson     }
809df49d7eSJed Brown }
819df49d7eSJed Brown 
829df49d7eSJed Brown // -----------------------------------------------------------------------------
837ed177dbSJed Brown // Basis context wrapper
849df49d7eSJed Brown // -----------------------------------------------------------------------------
85c68be7a2SJeremy L Thompson #[derive(Debug)]
869df49d7eSJed Brown pub struct Basis<'a> {
879df49d7eSJed Brown     pub(crate) ptr: bind_ceed::CeedBasis,
881142270cSJeremy L Thompson     _lifeline: PhantomData<&'a ()>,
899df49d7eSJed Brown }
909df49d7eSJed Brown 
919df49d7eSJed Brown // -----------------------------------------------------------------------------
929df49d7eSJed Brown // Destructor
939df49d7eSJed Brown // -----------------------------------------------------------------------------
949df49d7eSJed Brown impl<'a> Drop for Basis<'a> {
959df49d7eSJed Brown     fn drop(&mut self) {
969df49d7eSJed Brown         unsafe {
97356036faSJeremy L Thompson             if self.ptr != bind_ceed::CEED_BASIS_NONE {
989df49d7eSJed Brown                 bind_ceed::CeedBasisDestroy(&mut self.ptr);
999df49d7eSJed Brown             }
1009df49d7eSJed Brown         }
1019df49d7eSJed Brown     }
1029df49d7eSJed Brown }
1039df49d7eSJed Brown 
1049df49d7eSJed Brown // -----------------------------------------------------------------------------
1059df49d7eSJed Brown // Display
1069df49d7eSJed Brown // -----------------------------------------------------------------------------
1079df49d7eSJed Brown impl<'a> fmt::Display for Basis<'a> {
1089df49d7eSJed Brown     /// View a Basis
1099df49d7eSJed Brown     ///
1109df49d7eSJed Brown     /// ```
111eb07d68fSJeremy L Thompson     /// # use libceed::{prelude::*, QuadMode};
1124d27c890SJeremy L Thompson     /// # fn main() -> libceed::Result<()> {
1139df49d7eSJed Brown     /// # let ceed = libceed::Ceed::default_init();
114c68be7a2SJeremy L Thompson     /// let b = ceed.basis_tensor_H1_Lagrange(1, 2, 3, 4, QuadMode::Gauss)?;
1159df49d7eSJed Brown     /// println!("{}", b);
116c68be7a2SJeremy L Thompson     /// # Ok(())
117c68be7a2SJeremy L Thompson     /// # }
1189df49d7eSJed Brown     /// ```
1199df49d7eSJed Brown     fn fmt(&self, f: &mut fmt::Formatter) -> fmt::Result {
1209df49d7eSJed Brown         let mut ptr = std::ptr::null_mut();
1219df49d7eSJed Brown         let mut sizeloc = crate::MAX_BUFFER_LENGTH;
1229df49d7eSJed Brown         let cstring = unsafe {
1239df49d7eSJed Brown             let file = bind_ceed::open_memstream(&mut ptr, &mut sizeloc);
1249df49d7eSJed Brown             bind_ceed::CeedBasisView(self.ptr, file);
1259df49d7eSJed Brown             bind_ceed::fclose(file);
1269df49d7eSJed Brown             CString::from_raw(ptr)
1279df49d7eSJed Brown         };
1289df49d7eSJed Brown         cstring.to_string_lossy().fmt(f)
1299df49d7eSJed Brown     }
1309df49d7eSJed Brown }
1319df49d7eSJed Brown 
1329df49d7eSJed Brown // -----------------------------------------------------------------------------
1339df49d7eSJed Brown // Implementations
1349df49d7eSJed Brown // -----------------------------------------------------------------------------
1359df49d7eSJed Brown impl<'a> Basis<'a> {
1369df49d7eSJed Brown     // Constructors
13778c2cefaSJeremy L Thompson     #[allow(clippy::too_many_arguments)]
1389df49d7eSJed Brown     pub fn create_tensor_H1(
139594ef120SJeremy L Thompson         ceed: &crate::Ceed,
1409df49d7eSJed Brown         dim: usize,
1419df49d7eSJed Brown         ncomp: usize,
1429df49d7eSJed Brown         P1d: usize,
1439df49d7eSJed Brown         Q1d: usize,
14480a9ef05SNatalie Beams         interp1d: &[crate::Scalar],
14580a9ef05SNatalie Beams         grad1d: &[crate::Scalar],
14680a9ef05SNatalie Beams         qref1d: &[crate::Scalar],
14780a9ef05SNatalie Beams         qweight1d: &[crate::Scalar],
1489df49d7eSJed Brown     ) -> crate::Result<Self> {
1499df49d7eSJed Brown         let mut ptr = std::ptr::null_mut();
1509df49d7eSJed Brown         let (dim, ncomp, P1d, Q1d) = (
1519df49d7eSJed Brown             i32::try_from(dim).unwrap(),
1529df49d7eSJed Brown             i32::try_from(ncomp).unwrap(),
1539df49d7eSJed Brown             i32::try_from(P1d).unwrap(),
1549df49d7eSJed Brown             i32::try_from(Q1d).unwrap(),
1559df49d7eSJed Brown         );
156656ef1e5SJeremy L Thompson         ceed.check_error(unsafe {
1579df49d7eSJed Brown             bind_ceed::CeedBasisCreateTensorH1(
1589df49d7eSJed Brown                 ceed.ptr,
1599df49d7eSJed Brown                 dim,
1609df49d7eSJed Brown                 ncomp,
1619df49d7eSJed Brown                 P1d,
1629df49d7eSJed Brown                 Q1d,
1639df49d7eSJed Brown                 interp1d.as_ptr(),
1649df49d7eSJed Brown                 grad1d.as_ptr(),
1659df49d7eSJed Brown                 qref1d.as_ptr(),
1669df49d7eSJed Brown                 qweight1d.as_ptr(),
1679df49d7eSJed Brown                 &mut ptr,
1689df49d7eSJed Brown             )
169656ef1e5SJeremy L Thompson         })?;
1701142270cSJeremy L Thompson         Ok(Self {
1711142270cSJeremy L Thompson             ptr,
1721142270cSJeremy L Thompson             _lifeline: PhantomData,
1731142270cSJeremy L Thompson         })
1749df49d7eSJed Brown     }
1759df49d7eSJed Brown 
1762b671a0aSJeremy L Thompson     pub(crate) unsafe fn from_raw(ptr: bind_ceed::CeedBasis) -> crate::Result<Self> {
177e03fef56SJeremy L Thompson         Ok(Self {
178e03fef56SJeremy L Thompson             ptr,
179e03fef56SJeremy L Thompson             _lifeline: PhantomData,
180e03fef56SJeremy L Thompson         })
181e03fef56SJeremy L Thompson     }
182e03fef56SJeremy L Thompson 
1839df49d7eSJed Brown     pub fn create_tensor_H1_Lagrange(
184594ef120SJeremy L Thompson         ceed: &crate::Ceed,
1859df49d7eSJed Brown         dim: usize,
1869df49d7eSJed Brown         ncomp: usize,
1879df49d7eSJed Brown         P: usize,
1889df49d7eSJed Brown         Q: usize,
1899df49d7eSJed Brown         qmode: crate::QuadMode,
1909df49d7eSJed Brown     ) -> crate::Result<Self> {
1919df49d7eSJed Brown         let mut ptr = std::ptr::null_mut();
1929df49d7eSJed Brown         let (dim, ncomp, P, Q, qmode) = (
1939df49d7eSJed Brown             i32::try_from(dim).unwrap(),
1949df49d7eSJed Brown             i32::try_from(ncomp).unwrap(),
1959df49d7eSJed Brown             i32::try_from(P).unwrap(),
1969df49d7eSJed Brown             i32::try_from(Q).unwrap(),
1979df49d7eSJed Brown             qmode as bind_ceed::CeedQuadMode,
1989df49d7eSJed Brown         );
199656ef1e5SJeremy L Thompson         ceed.check_error(unsafe {
2009df49d7eSJed Brown             bind_ceed::CeedBasisCreateTensorH1Lagrange(ceed.ptr, dim, ncomp, P, Q, qmode, &mut ptr)
201656ef1e5SJeremy L Thompson         })?;
2021142270cSJeremy L Thompson         Ok(Self {
2031142270cSJeremy L Thompson             ptr,
2041142270cSJeremy L Thompson             _lifeline: PhantomData,
2051142270cSJeremy L Thompson         })
2069df49d7eSJed Brown     }
2079df49d7eSJed Brown 
20878c2cefaSJeremy L Thompson     #[allow(clippy::too_many_arguments)]
2099df49d7eSJed Brown     pub fn create_H1(
210594ef120SJeremy L Thompson         ceed: &crate::Ceed,
2119df49d7eSJed Brown         topo: crate::ElemTopology,
2129df49d7eSJed Brown         ncomp: usize,
2139df49d7eSJed Brown         nnodes: usize,
2149df49d7eSJed Brown         nqpts: usize,
21580a9ef05SNatalie Beams         interp: &[crate::Scalar],
21680a9ef05SNatalie Beams         grad: &[crate::Scalar],
21780a9ef05SNatalie Beams         qref: &[crate::Scalar],
21880a9ef05SNatalie Beams         qweight: &[crate::Scalar],
2199df49d7eSJed Brown     ) -> crate::Result<Self> {
2209df49d7eSJed Brown         let mut ptr = std::ptr::null_mut();
2219df49d7eSJed Brown         let (topo, ncomp, nnodes, nqpts) = (
2229df49d7eSJed Brown             topo as bind_ceed::CeedElemTopology,
2239df49d7eSJed Brown             i32::try_from(ncomp).unwrap(),
2249df49d7eSJed Brown             i32::try_from(nnodes).unwrap(),
2259df49d7eSJed Brown             i32::try_from(nqpts).unwrap(),
2269df49d7eSJed Brown         );
227656ef1e5SJeremy L Thompson         ceed.check_error(unsafe {
2289df49d7eSJed Brown             bind_ceed::CeedBasisCreateH1(
2299df49d7eSJed Brown                 ceed.ptr,
2309df49d7eSJed Brown                 topo,
2319df49d7eSJed Brown                 ncomp,
2329df49d7eSJed Brown                 nnodes,
2339df49d7eSJed Brown                 nqpts,
2349df49d7eSJed Brown                 interp.as_ptr(),
2359df49d7eSJed Brown                 grad.as_ptr(),
2369df49d7eSJed Brown                 qref.as_ptr(),
2379df49d7eSJed Brown                 qweight.as_ptr(),
2389df49d7eSJed Brown                 &mut ptr,
2399df49d7eSJed Brown             )
240656ef1e5SJeremy L Thompson         })?;
2411142270cSJeremy L Thompson         Ok(Self {
2421142270cSJeremy L Thompson             ptr,
2431142270cSJeremy L Thompson             _lifeline: PhantomData,
2441142270cSJeremy L Thompson         })
2451142270cSJeremy L Thompson     }
2461142270cSJeremy L Thompson 
24778c2cefaSJeremy L Thompson     #[allow(clippy::too_many_arguments)]
24811b88ddaSSebastian Grimberg     pub fn create_Hdiv(
24911b88ddaSSebastian Grimberg         ceed: &crate::Ceed,
25011b88ddaSSebastian Grimberg         topo: crate::ElemTopology,
25111b88ddaSSebastian Grimberg         ncomp: usize,
25211b88ddaSSebastian Grimberg         nnodes: usize,
25311b88ddaSSebastian Grimberg         nqpts: usize,
25411b88ddaSSebastian Grimberg         interp: &[crate::Scalar],
25511b88ddaSSebastian Grimberg         div: &[crate::Scalar],
25611b88ddaSSebastian Grimberg         qref: &[crate::Scalar],
25711b88ddaSSebastian Grimberg         qweight: &[crate::Scalar],
25811b88ddaSSebastian Grimberg     ) -> crate::Result<Self> {
25911b88ddaSSebastian Grimberg         let mut ptr = std::ptr::null_mut();
26011b88ddaSSebastian Grimberg         let (topo, ncomp, nnodes, nqpts) = (
26111b88ddaSSebastian Grimberg             topo as bind_ceed::CeedElemTopology,
26211b88ddaSSebastian Grimberg             i32::try_from(ncomp).unwrap(),
26311b88ddaSSebastian Grimberg             i32::try_from(nnodes).unwrap(),
26411b88ddaSSebastian Grimberg             i32::try_from(nqpts).unwrap(),
26511b88ddaSSebastian Grimberg         );
266656ef1e5SJeremy L Thompson         ceed.check_error(unsafe {
26711b88ddaSSebastian Grimberg             bind_ceed::CeedBasisCreateHdiv(
26811b88ddaSSebastian Grimberg                 ceed.ptr,
26911b88ddaSSebastian Grimberg                 topo,
27011b88ddaSSebastian Grimberg                 ncomp,
27111b88ddaSSebastian Grimberg                 nnodes,
27211b88ddaSSebastian Grimberg                 nqpts,
27311b88ddaSSebastian Grimberg                 interp.as_ptr(),
27411b88ddaSSebastian Grimberg                 div.as_ptr(),
27511b88ddaSSebastian Grimberg                 qref.as_ptr(),
27611b88ddaSSebastian Grimberg                 qweight.as_ptr(),
27711b88ddaSSebastian Grimberg                 &mut ptr,
27811b88ddaSSebastian Grimberg             )
279656ef1e5SJeremy L Thompson         })?;
28011b88ddaSSebastian Grimberg         Ok(Self {
28111b88ddaSSebastian Grimberg             ptr,
28211b88ddaSSebastian Grimberg             _lifeline: PhantomData,
28311b88ddaSSebastian Grimberg         })
28411b88ddaSSebastian Grimberg     }
28511b88ddaSSebastian Grimberg 
28678c2cefaSJeremy L Thompson     #[allow(clippy::too_many_arguments)]
28711b88ddaSSebastian Grimberg     pub fn create_Hcurl(
28811b88ddaSSebastian Grimberg         ceed: &crate::Ceed,
28911b88ddaSSebastian Grimberg         topo: crate::ElemTopology,
29011b88ddaSSebastian Grimberg         ncomp: usize,
29111b88ddaSSebastian Grimberg         nnodes: usize,
29211b88ddaSSebastian Grimberg         nqpts: usize,
29311b88ddaSSebastian Grimberg         interp: &[crate::Scalar],
29411b88ddaSSebastian Grimberg         curl: &[crate::Scalar],
29511b88ddaSSebastian Grimberg         qref: &[crate::Scalar],
29611b88ddaSSebastian Grimberg         qweight: &[crate::Scalar],
29711b88ddaSSebastian Grimberg     ) -> crate::Result<Self> {
29811b88ddaSSebastian Grimberg         let mut ptr = std::ptr::null_mut();
29911b88ddaSSebastian Grimberg         let (topo, ncomp, nnodes, nqpts) = (
30011b88ddaSSebastian Grimberg             topo as bind_ceed::CeedElemTopology,
30111b88ddaSSebastian Grimberg             i32::try_from(ncomp).unwrap(),
30211b88ddaSSebastian Grimberg             i32::try_from(nnodes).unwrap(),
30311b88ddaSSebastian Grimberg             i32::try_from(nqpts).unwrap(),
30411b88ddaSSebastian Grimberg         );
305656ef1e5SJeremy L Thompson         ceed.check_error(unsafe {
30611b88ddaSSebastian Grimberg             bind_ceed::CeedBasisCreateHcurl(
30711b88ddaSSebastian Grimberg                 ceed.ptr,
30811b88ddaSSebastian Grimberg                 topo,
30911b88ddaSSebastian Grimberg                 ncomp,
31011b88ddaSSebastian Grimberg                 nnodes,
31111b88ddaSSebastian Grimberg                 nqpts,
31211b88ddaSSebastian Grimberg                 interp.as_ptr(),
31311b88ddaSSebastian Grimberg                 curl.as_ptr(),
31411b88ddaSSebastian Grimberg                 qref.as_ptr(),
31511b88ddaSSebastian Grimberg                 qweight.as_ptr(),
31611b88ddaSSebastian Grimberg                 &mut ptr,
31711b88ddaSSebastian Grimberg             )
318656ef1e5SJeremy L Thompson         })?;
31911b88ddaSSebastian Grimberg         Ok(Self {
32011b88ddaSSebastian Grimberg             ptr,
32111b88ddaSSebastian Grimberg             _lifeline: PhantomData,
32211b88ddaSSebastian Grimberg         })
32311b88ddaSSebastian Grimberg     }
32411b88ddaSSebastian Grimberg 
32511544396SJeremy L Thompson     // Raw Ceed for error handling
32611544396SJeremy L Thompson     #[doc(hidden)]
32711544396SJeremy L Thompson     fn ceed(&self) -> bind_ceed::Ceed {
32811544396SJeremy L Thompson         unsafe { bind_ceed::CeedBasisReturnCeed(self.ptr) }
32911544396SJeremy L Thompson     }
33011544396SJeremy L Thompson 
3311142270cSJeremy L Thompson     // Error handling
3321142270cSJeremy L Thompson     #[doc(hidden)]
3331142270cSJeremy L Thompson     fn check_error(&self, ierr: i32) -> crate::Result<i32> {
33411544396SJeremy L Thompson         crate::check_error(|| self.ceed(), ierr)
3359df49d7eSJed Brown     }
3369df49d7eSJed Brown 
3379df49d7eSJed Brown     /// Apply basis evaluation from nodes to quadrature points or vice versa
3389df49d7eSJed Brown     ///
3399df49d7eSJed Brown     /// * `nelem` - The number of elements to apply the basis evaluation to
3409df49d7eSJed Brown     /// * `tmode` - `TrasposeMode::NoTranspose` to evaluate from nodes to
3419df49d7eSJed Brown     ///               quadrature points, `TransposeMode::Transpose` to apply the
3429df49d7eSJed Brown     ///               transpose, mapping from quadrature points to nodes
3439df49d7eSJed Brown     /// * `emode` - `EvalMode::None` to use values directly, `EvalMode::Interp`
3449df49d7eSJed Brown     ///               to use interpolated values, `EvalMode::Grad` to use
3459df49d7eSJed Brown     ///               gradients, `EvalMode::Weight` to use quadrature weights
3469df49d7eSJed Brown     /// * `u`     - Input Vector
3479df49d7eSJed Brown     /// * `v`     - Output Vector
3489df49d7eSJed Brown     ///
3499df49d7eSJed Brown     /// ```
350eb07d68fSJeremy L Thompson     /// # use libceed::{prelude::*, EvalMode, TransposeMode, QuadMode};
3514d27c890SJeremy L Thompson     /// # fn main() -> libceed::Result<()> {
3529df49d7eSJed Brown     /// # let ceed = libceed::Ceed::default_init();
3539df49d7eSJed Brown     /// const Q: usize = 6;
354c68be7a2SJeremy L Thompson     /// let bu = ceed.basis_tensor_H1_Lagrange(1, 1, Q, Q, QuadMode::GaussLobatto)?;
355c68be7a2SJeremy L Thompson     /// let bx = ceed.basis_tensor_H1_Lagrange(1, 1, 2, Q, QuadMode::Gauss)?;
3569df49d7eSJed Brown     ///
357c68be7a2SJeremy L Thompson     /// let x_corners = ceed.vector_from_slice(&[-1., 1.])?;
358c68be7a2SJeremy L Thompson     /// let mut x_qpts = ceed.vector(Q)?;
359c68be7a2SJeremy L Thompson     /// let mut x_nodes = ceed.vector(Q)?;
3609df49d7eSJed Brown     /// bx.apply(
3619df49d7eSJed Brown     ///     1,
3629df49d7eSJed Brown     ///     TransposeMode::NoTranspose,
3639df49d7eSJed Brown     ///     EvalMode::Interp,
3649df49d7eSJed Brown     ///     &x_corners,
3659df49d7eSJed Brown     ///     &mut x_nodes,
366d7f01795SJeremy L Thompson     /// )?;
3679df49d7eSJed Brown     /// bu.apply(
3689df49d7eSJed Brown     ///     1,
3699df49d7eSJed Brown     ///     TransposeMode::NoTranspose,
3709df49d7eSJed Brown     ///     EvalMode::Interp,
3719df49d7eSJed Brown     ///     &x_nodes,
3729df49d7eSJed Brown     ///     &mut x_qpts,
373d7f01795SJeremy L Thompson     /// )?;
3749df49d7eSJed Brown     ///
3759df49d7eSJed Brown     /// // Create function x^3 + 1 on Gauss Lobatto points
3769df49d7eSJed Brown     /// let mut u_arr = [0.; Q];
3779df49d7eSJed Brown     /// u_arr
3789df49d7eSJed Brown     ///     .iter_mut()
379e78171edSJeremy L Thompson     ///     .zip(x_nodes.view()?.iter())
3809df49d7eSJed Brown     ///     .for_each(|(u, x)| *u = x * x * x + 1.);
381c68be7a2SJeremy L Thompson     /// let u = ceed.vector_from_slice(&u_arr)?;
3829df49d7eSJed Brown     ///
3839df49d7eSJed Brown     /// // Map function to Gauss points
384c68be7a2SJeremy L Thompson     /// let mut v = ceed.vector(Q)?;
3859df49d7eSJed Brown     /// v.set_value(0.);
386c68be7a2SJeremy L Thompson     /// bu.apply(1, TransposeMode::NoTranspose, EvalMode::Interp, &u, &mut v)?;
3879df49d7eSJed Brown     ///
3889df49d7eSJed Brown     /// // Verify results
389e78171edSJeremy L Thompson     /// v.view()?
3909df49d7eSJed Brown     ///     .iter()
391e78171edSJeremy L Thompson     ///     .zip(x_qpts.view()?.iter())
3929df49d7eSJed Brown     ///     .for_each(|(v, x)| {
3939df49d7eSJed Brown     ///         let true_value = x * x * x + 1.;
39480a9ef05SNatalie Beams     ///         assert!(
39580a9ef05SNatalie Beams     ///             (*v - true_value).abs() < 10.0 * libceed::EPSILON,
39680a9ef05SNatalie Beams     ///             "Incorrect basis application"
39780a9ef05SNatalie Beams     ///         );
3989df49d7eSJed Brown     ///     });
399c68be7a2SJeremy L Thompson     /// # Ok(())
400c68be7a2SJeremy L Thompson     /// # }
4019df49d7eSJed Brown     /// ```
4029df49d7eSJed Brown     pub fn apply(
4039df49d7eSJed Brown         &self,
4049df49d7eSJed Brown         nelem: usize,
4059df49d7eSJed Brown         tmode: TransposeMode,
4069df49d7eSJed Brown         emode: EvalMode,
4079df49d7eSJed Brown         u: &Vector,
4089df49d7eSJed Brown         v: &mut Vector,
4099df49d7eSJed Brown     ) -> crate::Result<i32> {
4109df49d7eSJed Brown         let (nelem, tmode, emode) = (
4119df49d7eSJed Brown             i32::try_from(nelem).unwrap(),
4129df49d7eSJed Brown             tmode as bind_ceed::CeedTransposeMode,
4139df49d7eSJed Brown             emode as bind_ceed::CeedEvalMode,
4149df49d7eSJed Brown         );
415656ef1e5SJeremy L Thompson         self.check_error(unsafe {
416656ef1e5SJeremy L Thompson             bind_ceed::CeedBasisApply(self.ptr, nelem, tmode, emode, u.ptr, v.ptr)
417656ef1e5SJeremy L Thompson         })
4189df49d7eSJed Brown     }
4199df49d7eSJed Brown 
4207ed177dbSJed Brown     /// Returns the dimension for given Basis
4219df49d7eSJed Brown     ///
4229df49d7eSJed Brown     /// ```
423eb07d68fSJeremy L Thompson     /// # use libceed::{prelude::*, QuadMode};
4244d27c890SJeremy L Thompson     /// # fn main() -> libceed::Result<()> {
4259df49d7eSJed Brown     /// # let ceed = libceed::Ceed::default_init();
4269df49d7eSJed Brown     /// let dim = 2;
427c68be7a2SJeremy L Thompson     /// let b = ceed.basis_tensor_H1_Lagrange(dim, 1, 3, 4, QuadMode::Gauss)?;
4289df49d7eSJed Brown     ///
4299df49d7eSJed Brown     /// let d = b.dimension();
4309df49d7eSJed Brown     /// assert_eq!(d, dim, "Incorrect dimension");
431c68be7a2SJeremy L Thompson     /// # Ok(())
432c68be7a2SJeremy L Thompson     /// # }
4339df49d7eSJed Brown     /// ```
4349df49d7eSJed Brown     pub fn dimension(&self) -> usize {
4359df49d7eSJed Brown         let mut dim = 0;
4369df49d7eSJed Brown         unsafe { bind_ceed::CeedBasisGetDimension(self.ptr, &mut dim) };
4379df49d7eSJed Brown         usize::try_from(dim).unwrap()
4389df49d7eSJed Brown     }
4399df49d7eSJed Brown 
4407ed177dbSJed Brown     /// Returns number of components for given Basis
4419df49d7eSJed Brown     ///
4429df49d7eSJed Brown     /// ```
443eb07d68fSJeremy L Thompson     /// # use libceed::{prelude::*, QuadMode};
4444d27c890SJeremy L Thompson     /// # fn main() -> libceed::Result<()> {
4459df49d7eSJed Brown     /// # let ceed = libceed::Ceed::default_init();
4469df49d7eSJed Brown     /// let ncomp = 2;
447c68be7a2SJeremy L Thompson     /// let b = ceed.basis_tensor_H1_Lagrange(1, ncomp, 3, 4, QuadMode::Gauss)?;
4489df49d7eSJed Brown     ///
4499df49d7eSJed Brown     /// let n = b.num_components();
4509df49d7eSJed Brown     /// assert_eq!(n, ncomp, "Incorrect number of components");
451c68be7a2SJeremy L Thompson     /// # Ok(())
452c68be7a2SJeremy L Thompson     /// # }
4539df49d7eSJed Brown     /// ```
4549df49d7eSJed Brown     pub fn num_components(&self) -> usize {
4559df49d7eSJed Brown         let mut ncomp = 0;
4569df49d7eSJed Brown         unsafe { bind_ceed::CeedBasisGetNumComponents(self.ptr, &mut ncomp) };
4579df49d7eSJed Brown         usize::try_from(ncomp).unwrap()
4589df49d7eSJed Brown     }
4599df49d7eSJed Brown 
4607ed177dbSJed Brown     /// Returns total number of nodes (in dim dimensions) of a Basis
4619df49d7eSJed Brown     ///
4629df49d7eSJed Brown     /// ```
463eb07d68fSJeremy L Thompson     /// # use libceed::{prelude::*, QuadMode};
4644d27c890SJeremy L Thompson     /// # fn main() -> libceed::Result<()> {
4659df49d7eSJed Brown     /// # let ceed = libceed::Ceed::default_init();
4669df49d7eSJed Brown     /// let p = 3;
467c68be7a2SJeremy L Thompson     /// let b = ceed.basis_tensor_H1_Lagrange(2, 1, p, 4, QuadMode::Gauss)?;
4689df49d7eSJed Brown     ///
4699df49d7eSJed Brown     /// let nnodes = b.num_nodes();
4709df49d7eSJed Brown     /// assert_eq!(nnodes, p * p, "Incorrect number of nodes");
471c68be7a2SJeremy L Thompson     /// # Ok(())
472c68be7a2SJeremy L Thompson     /// # }
4739df49d7eSJed Brown     /// ```
4749df49d7eSJed Brown     pub fn num_nodes(&self) -> usize {
4759df49d7eSJed Brown         let mut nnodes = 0;
4769df49d7eSJed Brown         unsafe { bind_ceed::CeedBasisGetNumNodes(self.ptr, &mut nnodes) };
4779df49d7eSJed Brown         usize::try_from(nnodes).unwrap()
4789df49d7eSJed Brown     }
4799df49d7eSJed Brown 
4809df49d7eSJed Brown     /// Returns total number of quadrature points (in dim dimensions) of a
4817ed177dbSJed Brown     /// Basis
4829df49d7eSJed Brown     ///
4839df49d7eSJed Brown     /// ```
484eb07d68fSJeremy L Thompson     /// # use libceed::{prelude::*, QuadMode};
4854d27c890SJeremy L Thompson     /// # fn main() -> libceed::Result<()> {
4869df49d7eSJed Brown     /// # let ceed = libceed::Ceed::default_init();
4879df49d7eSJed Brown     /// let q = 4;
488c68be7a2SJeremy L Thompson     /// let b = ceed.basis_tensor_H1_Lagrange(2, 1, 3, q, QuadMode::Gauss)?;
4899df49d7eSJed Brown     ///
4909df49d7eSJed Brown     /// let nqpts = b.num_quadrature_points();
4919df49d7eSJed Brown     /// assert_eq!(nqpts, q * q, "Incorrect number of quadrature points");
492c68be7a2SJeremy L Thompson     /// # Ok(())
493c68be7a2SJeremy L Thompson     /// # }
4949df49d7eSJed Brown     /// ```
4959df49d7eSJed Brown     pub fn num_quadrature_points(&self) -> usize {
4969df49d7eSJed Brown         let mut Q = 0;
4979df49d7eSJed Brown         unsafe {
4989df49d7eSJed Brown             bind_ceed::CeedBasisGetNumQuadraturePoints(self.ptr, &mut Q);
4999df49d7eSJed Brown         }
5009df49d7eSJed Brown         usize::try_from(Q).unwrap()
5019df49d7eSJed Brown     }
50200d548f6SJed Brown 
50300d548f6SJed Brown     /// Create projection from self to specified Basis.
50400d548f6SJed Brown     ///
50500d548f6SJed Brown     /// Both bases must have the same quadrature space. The input bases need not
50600d548f6SJed Brown     /// be nested as function spaces; this interface solves a least squares
507b748b478SJeremy L Thompson     /// problem to find a representation in the `to` basis that agrees at
50800d548f6SJed Brown     /// quadrature points with the origin basis. Since the bases need not be
50900d548f6SJed Brown     /// Lagrange, the resulting projection "basis" will have empty quadrature
51000d548f6SJed Brown     /// points and weights.
51100d548f6SJed Brown     ///
51200d548f6SJed Brown     /// ```
513eb07d68fSJeremy L Thompson     /// # use libceed::{prelude::*, EvalMode, TransposeMode, QuadMode};
51400d548f6SJed Brown     /// # fn main() -> libceed::Result<()> {
51500d548f6SJed Brown     /// # let ceed = libceed::Ceed::default_init();
51600d548f6SJed Brown     /// let coarse = ceed.basis_tensor_H1_Lagrange(1, 1, 2, 3, QuadMode::Gauss)?;
51700d548f6SJed Brown     /// let fine = ceed.basis_tensor_H1_Lagrange(1, 1, 3, 3, QuadMode::Gauss)?;
51800d548f6SJed Brown     /// let proj = coarse.create_projection(&fine)?;
51900d548f6SJed Brown     /// let u = ceed.vector_from_slice(&[1., 2.])?;
52000d548f6SJed Brown     /// let mut v = ceed.vector(3)?;
52100d548f6SJed Brown     /// proj.apply(1, TransposeMode::NoTranspose, EvalMode::Interp, &u, &mut v)?;
52200d548f6SJed Brown     /// let expected = [1., 1.5, 2.];
52300d548f6SJed Brown     /// for (a, b) in v.view()?.iter().zip(expected) {
52400d548f6SJed Brown     ///     assert!(
52500d548f6SJed Brown     ///         (a - b).abs() < 10.0 * libceed::EPSILON,
52600d548f6SJed Brown     ///         "Incorrect projection of linear Lagrange to quadratic Lagrange"
52700d548f6SJed Brown     ///     );
52800d548f6SJed Brown     /// }
52900d548f6SJed Brown     /// # Ok(())
53000d548f6SJed Brown     /// # }
53100d548f6SJed Brown     /// ```
53200d548f6SJed Brown     pub fn create_projection(&self, to: &Self) -> crate::Result<Self> {
53300d548f6SJed Brown         let mut ptr = std::ptr::null_mut();
534656ef1e5SJeremy L Thompson         self.check_error(unsafe {
535656ef1e5SJeremy L Thompson             bind_ceed::CeedBasisCreateProjection(self.ptr, to.ptr, &mut ptr)
536656ef1e5SJeremy L Thompson         })?;
53700d548f6SJed Brown         Ok(Self {
53800d548f6SJed Brown             ptr,
53900d548f6SJed Brown             _lifeline: PhantomData,
54000d548f6SJed Brown         })
54100d548f6SJed Brown     }
5429df49d7eSJed Brown }
5439df49d7eSJed Brown 
5449df49d7eSJed Brown // -----------------------------------------------------------------------------
545