xref: /libCEED/rust/libceed/src/basis.rs (revision 80a9ef0545a39c00cdcaab1ca26f8053604f3120)
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 
179df49d7eSJed Brown //! A Ceed Basis defines the discrete finite element basis and associated
189df49d7eSJed Brown //! quadrature rule.
199df49d7eSJed Brown 
209df49d7eSJed Brown use crate::prelude::*;
219df49d7eSJed Brown 
229df49d7eSJed Brown // -----------------------------------------------------------------------------
239df49d7eSJed Brown // CeedBasis option
249df49d7eSJed Brown // -----------------------------------------------------------------------------
259df49d7eSJed Brown #[derive(Clone, Copy)]
269df49d7eSJed Brown pub enum BasisOpt<'a> {
279df49d7eSJed Brown     Some(&'a Basis<'a>),
289df49d7eSJed Brown     Collocated,
299df49d7eSJed Brown }
309df49d7eSJed Brown /// Construct a BasisOpt reference from a Basis reference
319df49d7eSJed Brown impl<'a> From<&'a Basis<'_>> for BasisOpt<'a> {
329df49d7eSJed Brown     fn from(basis: &'a Basis) -> Self {
339df49d7eSJed Brown         debug_assert!(basis.ptr != unsafe { bind_ceed::CEED_BASIS_COLLOCATED });
349df49d7eSJed Brown         Self::Some(basis)
359df49d7eSJed Brown     }
369df49d7eSJed Brown }
379df49d7eSJed Brown impl<'a> BasisOpt<'a> {
389df49d7eSJed Brown     /// Transform a Rust libCEED BasisOpt into C libCEED CeedBasis
399df49d7eSJed Brown     pub(crate) fn to_raw(self) -> bind_ceed::CeedBasis {
409df49d7eSJed Brown         match self {
419df49d7eSJed Brown             Self::Some(basis) => basis.ptr,
429df49d7eSJed Brown             Self::Collocated => unsafe { bind_ceed::CEED_BASIS_COLLOCATED },
439df49d7eSJed Brown         }
449df49d7eSJed Brown     }
459df49d7eSJed Brown }
469df49d7eSJed Brown 
479df49d7eSJed Brown // -----------------------------------------------------------------------------
489df49d7eSJed Brown // CeedBasis context wrapper
499df49d7eSJed Brown // -----------------------------------------------------------------------------
509df49d7eSJed Brown pub struct Basis<'a> {
519df49d7eSJed Brown     ceed: &'a crate::Ceed,
529df49d7eSJed Brown     pub(crate) ptr: bind_ceed::CeedBasis,
539df49d7eSJed Brown }
549df49d7eSJed Brown 
559df49d7eSJed Brown // -----------------------------------------------------------------------------
569df49d7eSJed Brown // Destructor
579df49d7eSJed Brown // -----------------------------------------------------------------------------
589df49d7eSJed Brown impl<'a> Drop for Basis<'a> {
599df49d7eSJed Brown     fn drop(&mut self) {
609df49d7eSJed Brown         unsafe {
619df49d7eSJed Brown             if self.ptr != bind_ceed::CEED_BASIS_COLLOCATED {
629df49d7eSJed Brown                 bind_ceed::CeedBasisDestroy(&mut self.ptr);
639df49d7eSJed Brown             }
649df49d7eSJed Brown         }
659df49d7eSJed Brown     }
669df49d7eSJed Brown }
679df49d7eSJed Brown 
689df49d7eSJed Brown // -----------------------------------------------------------------------------
699df49d7eSJed Brown // Display
709df49d7eSJed Brown // -----------------------------------------------------------------------------
719df49d7eSJed Brown impl<'a> fmt::Display for Basis<'a> {
729df49d7eSJed Brown     /// View a Basis
739df49d7eSJed Brown     ///
749df49d7eSJed Brown     /// ```
759df49d7eSJed Brown     /// # use libceed::prelude::*;
769df49d7eSJed Brown     /// # let ceed = libceed::Ceed::default_init();
779df49d7eSJed Brown     /// let b = ceed
789df49d7eSJed Brown     ///     .basis_tensor_H1_Lagrange(1, 2, 3, 4, QuadMode::Gauss)
799df49d7eSJed Brown     ///     .unwrap();
809df49d7eSJed Brown     /// println!("{}", b);
819df49d7eSJed Brown     /// ```
829df49d7eSJed Brown     fn fmt(&self, f: &mut fmt::Formatter) -> fmt::Result {
839df49d7eSJed Brown         let mut ptr = std::ptr::null_mut();
849df49d7eSJed Brown         let mut sizeloc = crate::MAX_BUFFER_LENGTH;
859df49d7eSJed Brown         let cstring = unsafe {
869df49d7eSJed Brown             let file = bind_ceed::open_memstream(&mut ptr, &mut sizeloc);
879df49d7eSJed Brown             bind_ceed::CeedBasisView(self.ptr, file);
889df49d7eSJed Brown             bind_ceed::fclose(file);
899df49d7eSJed Brown             CString::from_raw(ptr)
909df49d7eSJed Brown         };
919df49d7eSJed Brown         cstring.to_string_lossy().fmt(f)
929df49d7eSJed Brown     }
939df49d7eSJed Brown }
949df49d7eSJed Brown 
959df49d7eSJed Brown // -----------------------------------------------------------------------------
969df49d7eSJed Brown // Implementations
979df49d7eSJed Brown // -----------------------------------------------------------------------------
989df49d7eSJed Brown impl<'a> Basis<'a> {
999df49d7eSJed Brown     // Constructors
1009df49d7eSJed Brown     pub fn create_tensor_H1(
1019df49d7eSJed Brown         ceed: &'a crate::Ceed,
1029df49d7eSJed Brown         dim: usize,
1039df49d7eSJed Brown         ncomp: usize,
1049df49d7eSJed Brown         P1d: usize,
1059df49d7eSJed Brown         Q1d: usize,
106*80a9ef05SNatalie Beams         interp1d: &[crate::Scalar],
107*80a9ef05SNatalie Beams         grad1d: &[crate::Scalar],
108*80a9ef05SNatalie Beams         qref1d: &[crate::Scalar],
109*80a9ef05SNatalie Beams         qweight1d: &[crate::Scalar],
1109df49d7eSJed Brown     ) -> crate::Result<Self> {
1119df49d7eSJed Brown         let mut ptr = std::ptr::null_mut();
1129df49d7eSJed Brown         let (dim, ncomp, P1d, Q1d) = (
1139df49d7eSJed Brown             i32::try_from(dim).unwrap(),
1149df49d7eSJed Brown             i32::try_from(ncomp).unwrap(),
1159df49d7eSJed Brown             i32::try_from(P1d).unwrap(),
1169df49d7eSJed Brown             i32::try_from(Q1d).unwrap(),
1179df49d7eSJed Brown         );
1189df49d7eSJed Brown         let ierr = unsafe {
1199df49d7eSJed Brown             bind_ceed::CeedBasisCreateTensorH1(
1209df49d7eSJed Brown                 ceed.ptr,
1219df49d7eSJed Brown                 dim,
1229df49d7eSJed Brown                 ncomp,
1239df49d7eSJed Brown                 P1d,
1249df49d7eSJed Brown                 Q1d,
1259df49d7eSJed Brown                 interp1d.as_ptr(),
1269df49d7eSJed Brown                 grad1d.as_ptr(),
1279df49d7eSJed Brown                 qref1d.as_ptr(),
1289df49d7eSJed Brown                 qweight1d.as_ptr(),
1299df49d7eSJed Brown                 &mut ptr,
1309df49d7eSJed Brown             )
1319df49d7eSJed Brown         };
1329df49d7eSJed Brown         ceed.check_error(ierr)?;
1339df49d7eSJed Brown         Ok(Self { ceed, ptr })
1349df49d7eSJed Brown     }
1359df49d7eSJed Brown 
1369df49d7eSJed Brown     pub fn create_tensor_H1_Lagrange(
1379df49d7eSJed Brown         ceed: &'a crate::Ceed,
1389df49d7eSJed Brown         dim: usize,
1399df49d7eSJed Brown         ncomp: usize,
1409df49d7eSJed Brown         P: usize,
1419df49d7eSJed Brown         Q: usize,
1429df49d7eSJed Brown         qmode: crate::QuadMode,
1439df49d7eSJed Brown     ) -> crate::Result<Self> {
1449df49d7eSJed Brown         let mut ptr = std::ptr::null_mut();
1459df49d7eSJed Brown         let (dim, ncomp, P, Q, qmode) = (
1469df49d7eSJed Brown             i32::try_from(dim).unwrap(),
1479df49d7eSJed Brown             i32::try_from(ncomp).unwrap(),
1489df49d7eSJed Brown             i32::try_from(P).unwrap(),
1499df49d7eSJed Brown             i32::try_from(Q).unwrap(),
1509df49d7eSJed Brown             qmode as bind_ceed::CeedQuadMode,
1519df49d7eSJed Brown         );
1529df49d7eSJed Brown         let ierr = unsafe {
1539df49d7eSJed Brown             bind_ceed::CeedBasisCreateTensorH1Lagrange(ceed.ptr, dim, ncomp, P, Q, qmode, &mut ptr)
1549df49d7eSJed Brown         };
1559df49d7eSJed Brown         ceed.check_error(ierr)?;
1569df49d7eSJed Brown         Ok(Self { ceed, ptr })
1579df49d7eSJed Brown     }
1589df49d7eSJed Brown 
1599df49d7eSJed Brown     pub fn create_H1(
1609df49d7eSJed Brown         ceed: &'a crate::Ceed,
1619df49d7eSJed Brown         topo: crate::ElemTopology,
1629df49d7eSJed Brown         ncomp: usize,
1639df49d7eSJed Brown         nnodes: usize,
1649df49d7eSJed Brown         nqpts: usize,
165*80a9ef05SNatalie Beams         interp: &[crate::Scalar],
166*80a9ef05SNatalie Beams         grad: &[crate::Scalar],
167*80a9ef05SNatalie Beams         qref: &[crate::Scalar],
168*80a9ef05SNatalie Beams         qweight: &[crate::Scalar],
1699df49d7eSJed Brown     ) -> crate::Result<Self> {
1709df49d7eSJed Brown         let mut ptr = std::ptr::null_mut();
1719df49d7eSJed Brown         let (topo, ncomp, nnodes, nqpts) = (
1729df49d7eSJed Brown             topo as bind_ceed::CeedElemTopology,
1739df49d7eSJed Brown             i32::try_from(ncomp).unwrap(),
1749df49d7eSJed Brown             i32::try_from(nnodes).unwrap(),
1759df49d7eSJed Brown             i32::try_from(nqpts).unwrap(),
1769df49d7eSJed Brown         );
1779df49d7eSJed Brown         let ierr = unsafe {
1789df49d7eSJed Brown             bind_ceed::CeedBasisCreateH1(
1799df49d7eSJed Brown                 ceed.ptr,
1809df49d7eSJed Brown                 topo,
1819df49d7eSJed Brown                 ncomp,
1829df49d7eSJed Brown                 nnodes,
1839df49d7eSJed Brown                 nqpts,
1849df49d7eSJed Brown                 interp.as_ptr(),
1859df49d7eSJed Brown                 grad.as_ptr(),
1869df49d7eSJed Brown                 qref.as_ptr(),
1879df49d7eSJed Brown                 qweight.as_ptr(),
1889df49d7eSJed Brown                 &mut ptr,
1899df49d7eSJed Brown             )
1909df49d7eSJed Brown         };
1919df49d7eSJed Brown         ceed.check_error(ierr)?;
1929df49d7eSJed Brown         Ok(Self { ceed, ptr })
1939df49d7eSJed Brown     }
1949df49d7eSJed Brown 
1959df49d7eSJed Brown     /// Apply basis evaluation from nodes to quadrature points or vice versa
1969df49d7eSJed Brown     ///
1979df49d7eSJed Brown     /// * `nelem` - The number of elements to apply the basis evaluation to
1989df49d7eSJed Brown     /// * `tmode` - `TrasposeMode::NoTranspose` to evaluate from nodes to
1999df49d7eSJed Brown     ///               quadrature points, `TransposeMode::Transpose` to apply the
2009df49d7eSJed Brown     ///               transpose, mapping from quadrature points to nodes
2019df49d7eSJed Brown     /// * `emode` - `EvalMode::None` to use values directly, `EvalMode::Interp`
2029df49d7eSJed Brown     ///               to use interpolated values, `EvalMode::Grad` to use
2039df49d7eSJed Brown     ///               gradients, `EvalMode::Weight` to use quadrature weights
2049df49d7eSJed Brown     /// * `u`     - Input Vector
2059df49d7eSJed Brown     /// * `v`     - Output Vector
2069df49d7eSJed Brown     ///
2079df49d7eSJed Brown     /// ```
2089df49d7eSJed Brown     /// # use libceed::prelude::*;
2099df49d7eSJed Brown     /// # let ceed = libceed::Ceed::default_init();
2109df49d7eSJed Brown     /// const Q: usize = 6;
2119df49d7eSJed Brown     /// let bu = ceed
2129df49d7eSJed Brown     ///     .basis_tensor_H1_Lagrange(1, 1, Q, Q, QuadMode::GaussLobatto)
2139df49d7eSJed Brown     ///     .unwrap();
2149df49d7eSJed Brown     /// let bx = ceed
2159df49d7eSJed Brown     ///     .basis_tensor_H1_Lagrange(1, 1, 2, Q, QuadMode::Gauss)
2169df49d7eSJed Brown     ///     .unwrap();
2179df49d7eSJed Brown     ///
2189df49d7eSJed Brown     /// let x_corners = ceed.vector_from_slice(&[-1., 1.]).unwrap();
2199df49d7eSJed Brown     /// let mut x_qpts = ceed.vector(Q).unwrap();
2209df49d7eSJed Brown     /// let mut x_nodes = ceed.vector(Q).unwrap();
2219df49d7eSJed Brown     /// bx.apply(
2229df49d7eSJed Brown     ///     1,
2239df49d7eSJed Brown     ///     TransposeMode::NoTranspose,
2249df49d7eSJed Brown     ///     EvalMode::Interp,
2259df49d7eSJed Brown     ///     &x_corners,
2269df49d7eSJed Brown     ///     &mut x_nodes,
2279df49d7eSJed Brown     /// );
2289df49d7eSJed Brown     /// bu.apply(
2299df49d7eSJed Brown     ///     1,
2309df49d7eSJed Brown     ///     TransposeMode::NoTranspose,
2319df49d7eSJed Brown     ///     EvalMode::Interp,
2329df49d7eSJed Brown     ///     &x_nodes,
2339df49d7eSJed Brown     ///     &mut x_qpts,
2349df49d7eSJed Brown     /// );
2359df49d7eSJed Brown     ///
2369df49d7eSJed Brown     /// // Create function x^3 + 1 on Gauss Lobatto points
2379df49d7eSJed Brown     /// let mut u_arr = [0.; Q];
2389df49d7eSJed Brown     /// u_arr
2399df49d7eSJed Brown     ///     .iter_mut()
2409df49d7eSJed Brown     ///     .zip(x_nodes.view().iter())
2419df49d7eSJed Brown     ///     .for_each(|(u, x)| *u = x * x * x + 1.);
2429df49d7eSJed Brown     /// let u = ceed.vector_from_slice(&u_arr).unwrap();
2439df49d7eSJed Brown     ///
2449df49d7eSJed Brown     /// // Map function to Gauss points
2459df49d7eSJed Brown     /// let mut v = ceed.vector(Q).unwrap();
2469df49d7eSJed Brown     /// v.set_value(0.);
2479df49d7eSJed Brown     /// bu.apply(1, TransposeMode::NoTranspose, EvalMode::Interp, &u, &mut v)
2489df49d7eSJed Brown     ///     .unwrap();
2499df49d7eSJed Brown     ///
2509df49d7eSJed Brown     /// // Verify results
2519df49d7eSJed Brown     /// v.view()
2529df49d7eSJed Brown     ///     .iter()
2539df49d7eSJed Brown     ///     .zip(x_qpts.view().iter())
2549df49d7eSJed Brown     ///     .for_each(|(v, x)| {
2559df49d7eSJed Brown     ///         let true_value = x * x * x + 1.;
256*80a9ef05SNatalie Beams     ///         assert!(
257*80a9ef05SNatalie Beams     ///             (*v - true_value).abs() < 10.0 * libceed::EPSILON,
258*80a9ef05SNatalie Beams     ///             "Incorrect basis application"
259*80a9ef05SNatalie Beams     ///         );
2609df49d7eSJed Brown     ///     });
2619df49d7eSJed Brown     /// ```
2629df49d7eSJed Brown     pub fn apply(
2639df49d7eSJed Brown         &self,
2649df49d7eSJed Brown         nelem: usize,
2659df49d7eSJed Brown         tmode: TransposeMode,
2669df49d7eSJed Brown         emode: EvalMode,
2679df49d7eSJed Brown         u: &Vector,
2689df49d7eSJed Brown         v: &mut Vector,
2699df49d7eSJed Brown     ) -> crate::Result<i32> {
2709df49d7eSJed Brown         let (nelem, tmode, emode) = (
2719df49d7eSJed Brown             i32::try_from(nelem).unwrap(),
2729df49d7eSJed Brown             tmode as bind_ceed::CeedTransposeMode,
2739df49d7eSJed Brown             emode as bind_ceed::CeedEvalMode,
2749df49d7eSJed Brown         );
2759df49d7eSJed Brown         let ierr =
2769df49d7eSJed Brown             unsafe { bind_ceed::CeedBasisApply(self.ptr, nelem, tmode, emode, u.ptr, v.ptr) };
2779df49d7eSJed Brown         self.ceed.check_error(ierr)
2789df49d7eSJed Brown     }
2799df49d7eSJed Brown 
2809df49d7eSJed Brown     /// Returns the dimension for given CeedBasis
2819df49d7eSJed Brown     ///
2829df49d7eSJed Brown     /// ```
2839df49d7eSJed Brown     /// # use libceed::prelude::*;
2849df49d7eSJed Brown     /// # let ceed = libceed::Ceed::default_init();
2859df49d7eSJed Brown     /// let dim = 2;
2869df49d7eSJed Brown     /// let b = ceed
2879df49d7eSJed Brown     ///     .basis_tensor_H1_Lagrange(dim, 1, 3, 4, QuadMode::Gauss)
2889df49d7eSJed Brown     ///     .unwrap();
2899df49d7eSJed Brown     ///
2909df49d7eSJed Brown     /// let d = b.dimension();
2919df49d7eSJed Brown     /// assert_eq!(d, dim, "Incorrect dimension");
2929df49d7eSJed Brown     /// ```
2939df49d7eSJed Brown     pub fn dimension(&self) -> usize {
2949df49d7eSJed Brown         let mut dim = 0;
2959df49d7eSJed Brown         unsafe { bind_ceed::CeedBasisGetDimension(self.ptr, &mut dim) };
2969df49d7eSJed Brown         usize::try_from(dim).unwrap()
2979df49d7eSJed Brown     }
2989df49d7eSJed Brown 
2999df49d7eSJed Brown     /// Returns number of components for given CeedBasis
3009df49d7eSJed Brown     ///
3019df49d7eSJed Brown     /// ```
3029df49d7eSJed Brown     /// # use libceed::prelude::*;
3039df49d7eSJed Brown     /// # let ceed = libceed::Ceed::default_init();
3049df49d7eSJed Brown     /// let ncomp = 2;
3059df49d7eSJed Brown     /// let b = ceed
3069df49d7eSJed Brown     ///     .basis_tensor_H1_Lagrange(1, ncomp, 3, 4, QuadMode::Gauss)
3079df49d7eSJed Brown     ///     .unwrap();
3089df49d7eSJed Brown     ///
3099df49d7eSJed Brown     /// let n = b.num_components();
3109df49d7eSJed Brown     /// assert_eq!(n, ncomp, "Incorrect number of components");
3119df49d7eSJed Brown     /// ```
3129df49d7eSJed Brown     pub fn num_components(&self) -> usize {
3139df49d7eSJed Brown         let mut ncomp = 0;
3149df49d7eSJed Brown         unsafe { bind_ceed::CeedBasisGetNumComponents(self.ptr, &mut ncomp) };
3159df49d7eSJed Brown         usize::try_from(ncomp).unwrap()
3169df49d7eSJed Brown     }
3179df49d7eSJed Brown 
3189df49d7eSJed Brown     /// Returns total number of nodes (in dim dimensions) of a CeedBasis
3199df49d7eSJed Brown     ///
3209df49d7eSJed Brown     /// ```
3219df49d7eSJed Brown     /// # use libceed::prelude::*;
3229df49d7eSJed Brown     /// # let ceed = libceed::Ceed::default_init();
3239df49d7eSJed Brown     /// let p = 3;
3249df49d7eSJed Brown     /// let b = ceed
3259df49d7eSJed Brown     ///     .basis_tensor_H1_Lagrange(2, 1, p, 4, QuadMode::Gauss)
3269df49d7eSJed Brown     ///     .unwrap();
3279df49d7eSJed Brown     ///
3289df49d7eSJed Brown     /// let nnodes = b.num_nodes();
3299df49d7eSJed Brown     /// assert_eq!(nnodes, p * p, "Incorrect number of nodes");
3309df49d7eSJed Brown     /// ```
3319df49d7eSJed Brown     pub fn num_nodes(&self) -> usize {
3329df49d7eSJed Brown         let mut nnodes = 0;
3339df49d7eSJed Brown         unsafe { bind_ceed::CeedBasisGetNumNodes(self.ptr, &mut nnodes) };
3349df49d7eSJed Brown         usize::try_from(nnodes).unwrap()
3359df49d7eSJed Brown     }
3369df49d7eSJed Brown 
3379df49d7eSJed Brown     /// Returns total number of quadrature points (in dim dimensions) of a
3389df49d7eSJed Brown     /// CeedBasis
3399df49d7eSJed Brown     ///
3409df49d7eSJed Brown     /// ```
3419df49d7eSJed Brown     /// # use libceed::prelude::*;
3429df49d7eSJed Brown     /// # let ceed = libceed::Ceed::default_init();
3439df49d7eSJed Brown     /// let q = 4;
3449df49d7eSJed Brown     /// let b = ceed
3459df49d7eSJed Brown     ///     .basis_tensor_H1_Lagrange(2, 1, 3, q, QuadMode::Gauss)
3469df49d7eSJed Brown     ///     .unwrap();
3479df49d7eSJed Brown     ///
3489df49d7eSJed Brown     /// let nqpts = b.num_quadrature_points();
3499df49d7eSJed Brown     /// assert_eq!(nqpts, q * q, "Incorrect number of quadrature points");
3509df49d7eSJed Brown     /// ```
3519df49d7eSJed Brown     pub fn num_quadrature_points(&self) -> usize {
3529df49d7eSJed Brown         let mut Q = 0;
3539df49d7eSJed Brown         unsafe {
3549df49d7eSJed Brown             bind_ceed::CeedBasisGetNumQuadraturePoints(self.ptr, &mut Q);
3559df49d7eSJed Brown         }
3569df49d7eSJed Brown         usize::try_from(Q).unwrap()
3579df49d7eSJed Brown     }
3589df49d7eSJed Brown }
3599df49d7eSJed Brown 
3609df49d7eSJed Brown // -----------------------------------------------------------------------------
361