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 // ----------------------------------------------------------------------------- 25c68be7a2SJeremy L Thompson #[derive(Debug)] 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 } 45*e03682afSJeremy L Thompson 46*e03682afSJeremy L Thompson /// Check if a BasisOpt is Some 47*e03682afSJeremy L Thompson /// 48*e03682afSJeremy L Thompson /// ``` 49*e03682afSJeremy L Thompson /// # use libceed::prelude::*; 50*e03682afSJeremy L Thompson /// # fn main() -> Result<(), libceed::CeedError> { 51*e03682afSJeremy L Thompson /// # let ceed = libceed::Ceed::default_init(); 52*e03682afSJeremy L Thompson /// let b = ceed.basis_tensor_H1_Lagrange(1, 2, 3, 4, QuadMode::Gauss)?; 53*e03682afSJeremy L Thompson /// let b_opt = BasisOpt::from(&b); 54*e03682afSJeremy L Thompson /// assert!(b_opt.is_some(), "Incorrect BasisOpt"); 55*e03682afSJeremy L Thompson /// 56*e03682afSJeremy L Thompson /// let b_opt = BasisOpt::Collocated; 57*e03682afSJeremy L Thompson /// assert!(!b_opt.is_some(), "Incorrect BasisOpt"); 58*e03682afSJeremy L Thompson /// # Ok(()) 59*e03682afSJeremy L Thompson /// # } 60*e03682afSJeremy L Thompson /// ``` 61*e03682afSJeremy L Thompson pub fn is_some(&self) -> bool { 62*e03682afSJeremy L Thompson match self { 63*e03682afSJeremy L Thompson Self::Some(_) => true, 64*e03682afSJeremy L Thompson Self::Collocated => false, 65*e03682afSJeremy L Thompson } 66*e03682afSJeremy L Thompson } 67*e03682afSJeremy L Thompson 68*e03682afSJeremy L Thompson /// Check if a BasisOpt is Collocated 69*e03682afSJeremy L Thompson /// 70*e03682afSJeremy L Thompson /// ``` 71*e03682afSJeremy L Thompson /// # use libceed::prelude::*; 72*e03682afSJeremy L Thompson /// # fn main() -> Result<(), libceed::CeedError> { 73*e03682afSJeremy L Thompson /// # let ceed = libceed::Ceed::default_init(); 74*e03682afSJeremy L Thompson /// let b = ceed.basis_tensor_H1_Lagrange(1, 2, 3, 4, QuadMode::Gauss)?; 75*e03682afSJeremy L Thompson /// let b_opt = BasisOpt::from(&b); 76*e03682afSJeremy L Thompson /// assert!(!b_opt.is_collocated(), "Incorrect BasisOpt"); 77*e03682afSJeremy L Thompson /// 78*e03682afSJeremy L Thompson /// let b_opt = BasisOpt::Collocated; 79*e03682afSJeremy L Thompson /// assert!(b_opt.is_collocated(), "Incorrect BasisOpt"); 80*e03682afSJeremy L Thompson /// # Ok(()) 81*e03682afSJeremy L Thompson /// # } 82*e03682afSJeremy L Thompson /// ``` 83*e03682afSJeremy L Thompson pub fn is_collocated(&self) -> bool { 84*e03682afSJeremy L Thompson match self { 85*e03682afSJeremy L Thompson Self::Some(_) => false, 86*e03682afSJeremy L Thompson Self::Collocated => true, 87*e03682afSJeremy L Thompson } 88*e03682afSJeremy L Thompson } 899df49d7eSJed Brown } 909df49d7eSJed Brown 919df49d7eSJed Brown // ----------------------------------------------------------------------------- 929df49d7eSJed Brown // CeedBasis context wrapper 939df49d7eSJed Brown // ----------------------------------------------------------------------------- 94c68be7a2SJeremy L Thompson #[derive(Debug)] 959df49d7eSJed Brown pub struct Basis<'a> { 969df49d7eSJed Brown ceed: &'a crate::Ceed, 979df49d7eSJed Brown pub(crate) ptr: bind_ceed::CeedBasis, 989df49d7eSJed Brown } 999df49d7eSJed Brown 1009df49d7eSJed Brown // ----------------------------------------------------------------------------- 1019df49d7eSJed Brown // Destructor 1029df49d7eSJed Brown // ----------------------------------------------------------------------------- 1039df49d7eSJed Brown impl<'a> Drop for Basis<'a> { 1049df49d7eSJed Brown fn drop(&mut self) { 1059df49d7eSJed Brown unsafe { 1069df49d7eSJed Brown if self.ptr != bind_ceed::CEED_BASIS_COLLOCATED { 1079df49d7eSJed Brown bind_ceed::CeedBasisDestroy(&mut self.ptr); 1089df49d7eSJed Brown } 1099df49d7eSJed Brown } 1109df49d7eSJed Brown } 1119df49d7eSJed Brown } 1129df49d7eSJed Brown 1139df49d7eSJed Brown // ----------------------------------------------------------------------------- 1149df49d7eSJed Brown // Display 1159df49d7eSJed Brown // ----------------------------------------------------------------------------- 1169df49d7eSJed Brown impl<'a> fmt::Display for Basis<'a> { 1179df49d7eSJed Brown /// View a Basis 1189df49d7eSJed Brown /// 1199df49d7eSJed Brown /// ``` 1209df49d7eSJed Brown /// # use libceed::prelude::*; 121c68be7a2SJeremy L Thompson /// # fn main() -> Result<(), libceed::CeedError> { 1229df49d7eSJed Brown /// # let ceed = libceed::Ceed::default_init(); 123c68be7a2SJeremy L Thompson /// let b = ceed.basis_tensor_H1_Lagrange(1, 2, 3, 4, QuadMode::Gauss)?; 1249df49d7eSJed Brown /// println!("{}", b); 125c68be7a2SJeremy L Thompson /// # Ok(()) 126c68be7a2SJeremy L Thompson /// # } 1279df49d7eSJed Brown /// ``` 1289df49d7eSJed Brown fn fmt(&self, f: &mut fmt::Formatter) -> fmt::Result { 1299df49d7eSJed Brown let mut ptr = std::ptr::null_mut(); 1309df49d7eSJed Brown let mut sizeloc = crate::MAX_BUFFER_LENGTH; 1319df49d7eSJed Brown let cstring = unsafe { 1329df49d7eSJed Brown let file = bind_ceed::open_memstream(&mut ptr, &mut sizeloc); 1339df49d7eSJed Brown bind_ceed::CeedBasisView(self.ptr, file); 1349df49d7eSJed Brown bind_ceed::fclose(file); 1359df49d7eSJed Brown CString::from_raw(ptr) 1369df49d7eSJed Brown }; 1379df49d7eSJed Brown cstring.to_string_lossy().fmt(f) 1389df49d7eSJed Brown } 1399df49d7eSJed Brown } 1409df49d7eSJed Brown 1419df49d7eSJed Brown // ----------------------------------------------------------------------------- 1429df49d7eSJed Brown // Implementations 1439df49d7eSJed Brown // ----------------------------------------------------------------------------- 1449df49d7eSJed Brown impl<'a> Basis<'a> { 1459df49d7eSJed Brown // Constructors 1469df49d7eSJed Brown pub fn create_tensor_H1( 1479df49d7eSJed Brown ceed: &'a crate::Ceed, 1489df49d7eSJed Brown dim: usize, 1499df49d7eSJed Brown ncomp: usize, 1509df49d7eSJed Brown P1d: usize, 1519df49d7eSJed Brown Q1d: usize, 15280a9ef05SNatalie Beams interp1d: &[crate::Scalar], 15380a9ef05SNatalie Beams grad1d: &[crate::Scalar], 15480a9ef05SNatalie Beams qref1d: &[crate::Scalar], 15580a9ef05SNatalie Beams qweight1d: &[crate::Scalar], 1569df49d7eSJed Brown ) -> crate::Result<Self> { 1579df49d7eSJed Brown let mut ptr = std::ptr::null_mut(); 1589df49d7eSJed Brown let (dim, ncomp, P1d, Q1d) = ( 1599df49d7eSJed Brown i32::try_from(dim).unwrap(), 1609df49d7eSJed Brown i32::try_from(ncomp).unwrap(), 1619df49d7eSJed Brown i32::try_from(P1d).unwrap(), 1629df49d7eSJed Brown i32::try_from(Q1d).unwrap(), 1639df49d7eSJed Brown ); 1649df49d7eSJed Brown let ierr = unsafe { 1659df49d7eSJed Brown bind_ceed::CeedBasisCreateTensorH1( 1669df49d7eSJed Brown ceed.ptr, 1679df49d7eSJed Brown dim, 1689df49d7eSJed Brown ncomp, 1699df49d7eSJed Brown P1d, 1709df49d7eSJed Brown Q1d, 1719df49d7eSJed Brown interp1d.as_ptr(), 1729df49d7eSJed Brown grad1d.as_ptr(), 1739df49d7eSJed Brown qref1d.as_ptr(), 1749df49d7eSJed Brown qweight1d.as_ptr(), 1759df49d7eSJed Brown &mut ptr, 1769df49d7eSJed Brown ) 1779df49d7eSJed Brown }; 1789df49d7eSJed Brown ceed.check_error(ierr)?; 1799df49d7eSJed Brown Ok(Self { ceed, ptr }) 1809df49d7eSJed Brown } 1819df49d7eSJed Brown 1829df49d7eSJed Brown pub fn create_tensor_H1_Lagrange( 1839df49d7eSJed Brown ceed: &'a crate::Ceed, 1849df49d7eSJed Brown dim: usize, 1859df49d7eSJed Brown ncomp: usize, 1869df49d7eSJed Brown P: usize, 1879df49d7eSJed Brown Q: usize, 1889df49d7eSJed Brown qmode: crate::QuadMode, 1899df49d7eSJed Brown ) -> crate::Result<Self> { 1909df49d7eSJed Brown let mut ptr = std::ptr::null_mut(); 1919df49d7eSJed Brown let (dim, ncomp, P, Q, qmode) = ( 1929df49d7eSJed Brown i32::try_from(dim).unwrap(), 1939df49d7eSJed Brown i32::try_from(ncomp).unwrap(), 1949df49d7eSJed Brown i32::try_from(P).unwrap(), 1959df49d7eSJed Brown i32::try_from(Q).unwrap(), 1969df49d7eSJed Brown qmode as bind_ceed::CeedQuadMode, 1979df49d7eSJed Brown ); 1989df49d7eSJed Brown let ierr = unsafe { 1999df49d7eSJed Brown bind_ceed::CeedBasisCreateTensorH1Lagrange(ceed.ptr, dim, ncomp, P, Q, qmode, &mut ptr) 2009df49d7eSJed Brown }; 2019df49d7eSJed Brown ceed.check_error(ierr)?; 2029df49d7eSJed Brown Ok(Self { ceed, ptr }) 2039df49d7eSJed Brown } 2049df49d7eSJed Brown 2059df49d7eSJed Brown pub fn create_H1( 2069df49d7eSJed Brown ceed: &'a crate::Ceed, 2079df49d7eSJed Brown topo: crate::ElemTopology, 2089df49d7eSJed Brown ncomp: usize, 2099df49d7eSJed Brown nnodes: usize, 2109df49d7eSJed Brown nqpts: usize, 21180a9ef05SNatalie Beams interp: &[crate::Scalar], 21280a9ef05SNatalie Beams grad: &[crate::Scalar], 21380a9ef05SNatalie Beams qref: &[crate::Scalar], 21480a9ef05SNatalie Beams qweight: &[crate::Scalar], 2159df49d7eSJed Brown ) -> crate::Result<Self> { 2169df49d7eSJed Brown let mut ptr = std::ptr::null_mut(); 2179df49d7eSJed Brown let (topo, ncomp, nnodes, nqpts) = ( 2189df49d7eSJed Brown topo as bind_ceed::CeedElemTopology, 2199df49d7eSJed Brown i32::try_from(ncomp).unwrap(), 2209df49d7eSJed Brown i32::try_from(nnodes).unwrap(), 2219df49d7eSJed Brown i32::try_from(nqpts).unwrap(), 2229df49d7eSJed Brown ); 2239df49d7eSJed Brown let ierr = unsafe { 2249df49d7eSJed Brown bind_ceed::CeedBasisCreateH1( 2259df49d7eSJed Brown ceed.ptr, 2269df49d7eSJed Brown topo, 2279df49d7eSJed Brown ncomp, 2289df49d7eSJed Brown nnodes, 2299df49d7eSJed Brown nqpts, 2309df49d7eSJed Brown interp.as_ptr(), 2319df49d7eSJed Brown grad.as_ptr(), 2329df49d7eSJed Brown qref.as_ptr(), 2339df49d7eSJed Brown qweight.as_ptr(), 2349df49d7eSJed Brown &mut ptr, 2359df49d7eSJed Brown ) 2369df49d7eSJed Brown }; 2379df49d7eSJed Brown ceed.check_error(ierr)?; 2389df49d7eSJed Brown Ok(Self { ceed, ptr }) 2399df49d7eSJed Brown } 2409df49d7eSJed Brown 2419df49d7eSJed Brown /// Apply basis evaluation from nodes to quadrature points or vice versa 2429df49d7eSJed Brown /// 2439df49d7eSJed Brown /// * `nelem` - The number of elements to apply the basis evaluation to 2449df49d7eSJed Brown /// * `tmode` - `TrasposeMode::NoTranspose` to evaluate from nodes to 2459df49d7eSJed Brown /// quadrature points, `TransposeMode::Transpose` to apply the 2469df49d7eSJed Brown /// transpose, mapping from quadrature points to nodes 2479df49d7eSJed Brown /// * `emode` - `EvalMode::None` to use values directly, `EvalMode::Interp` 2489df49d7eSJed Brown /// to use interpolated values, `EvalMode::Grad` to use 2499df49d7eSJed Brown /// gradients, `EvalMode::Weight` to use quadrature weights 2509df49d7eSJed Brown /// * `u` - Input Vector 2519df49d7eSJed Brown /// * `v` - Output Vector 2529df49d7eSJed Brown /// 2539df49d7eSJed Brown /// ``` 2549df49d7eSJed Brown /// # use libceed::prelude::*; 255c68be7a2SJeremy L Thompson /// # fn main() -> Result<(), libceed::CeedError> { 2569df49d7eSJed Brown /// # let ceed = libceed::Ceed::default_init(); 2579df49d7eSJed Brown /// const Q: usize = 6; 258c68be7a2SJeremy L Thompson /// let bu = ceed.basis_tensor_H1_Lagrange(1, 1, Q, Q, QuadMode::GaussLobatto)?; 259c68be7a2SJeremy L Thompson /// let bx = ceed.basis_tensor_H1_Lagrange(1, 1, 2, Q, QuadMode::Gauss)?; 2609df49d7eSJed Brown /// 261c68be7a2SJeremy L Thompson /// let x_corners = ceed.vector_from_slice(&[-1., 1.])?; 262c68be7a2SJeremy L Thompson /// let mut x_qpts = ceed.vector(Q)?; 263c68be7a2SJeremy L Thompson /// let mut x_nodes = ceed.vector(Q)?; 2649df49d7eSJed Brown /// bx.apply( 2659df49d7eSJed Brown /// 1, 2669df49d7eSJed Brown /// TransposeMode::NoTranspose, 2679df49d7eSJed Brown /// EvalMode::Interp, 2689df49d7eSJed Brown /// &x_corners, 2699df49d7eSJed Brown /// &mut x_nodes, 2709df49d7eSJed Brown /// ); 2719df49d7eSJed Brown /// bu.apply( 2729df49d7eSJed Brown /// 1, 2739df49d7eSJed Brown /// TransposeMode::NoTranspose, 2749df49d7eSJed Brown /// EvalMode::Interp, 2759df49d7eSJed Brown /// &x_nodes, 2769df49d7eSJed Brown /// &mut x_qpts, 2779df49d7eSJed Brown /// ); 2789df49d7eSJed Brown /// 2799df49d7eSJed Brown /// // Create function x^3 + 1 on Gauss Lobatto points 2809df49d7eSJed Brown /// let mut u_arr = [0.; Q]; 2819df49d7eSJed Brown /// u_arr 2829df49d7eSJed Brown /// .iter_mut() 2839df49d7eSJed Brown /// .zip(x_nodes.view().iter()) 2849df49d7eSJed Brown /// .for_each(|(u, x)| *u = x * x * x + 1.); 285c68be7a2SJeremy L Thompson /// let u = ceed.vector_from_slice(&u_arr)?; 2869df49d7eSJed Brown /// 2879df49d7eSJed Brown /// // Map function to Gauss points 288c68be7a2SJeremy L Thompson /// let mut v = ceed.vector(Q)?; 2899df49d7eSJed Brown /// v.set_value(0.); 290c68be7a2SJeremy L Thompson /// bu.apply(1, TransposeMode::NoTranspose, EvalMode::Interp, &u, &mut v)?; 2919df49d7eSJed Brown /// 2929df49d7eSJed Brown /// // Verify results 2939df49d7eSJed Brown /// v.view() 2949df49d7eSJed Brown /// .iter() 2959df49d7eSJed Brown /// .zip(x_qpts.view().iter()) 2969df49d7eSJed Brown /// .for_each(|(v, x)| { 2979df49d7eSJed Brown /// let true_value = x * x * x + 1.; 29880a9ef05SNatalie Beams /// assert!( 29980a9ef05SNatalie Beams /// (*v - true_value).abs() < 10.0 * libceed::EPSILON, 30080a9ef05SNatalie Beams /// "Incorrect basis application" 30180a9ef05SNatalie Beams /// ); 3029df49d7eSJed Brown /// }); 303c68be7a2SJeremy L Thompson /// # Ok(()) 304c68be7a2SJeremy L Thompson /// # } 3059df49d7eSJed Brown /// ``` 3069df49d7eSJed Brown pub fn apply( 3079df49d7eSJed Brown &self, 3089df49d7eSJed Brown nelem: usize, 3099df49d7eSJed Brown tmode: TransposeMode, 3109df49d7eSJed Brown emode: EvalMode, 3119df49d7eSJed Brown u: &Vector, 3129df49d7eSJed Brown v: &mut Vector, 3139df49d7eSJed Brown ) -> crate::Result<i32> { 3149df49d7eSJed Brown let (nelem, tmode, emode) = ( 3159df49d7eSJed Brown i32::try_from(nelem).unwrap(), 3169df49d7eSJed Brown tmode as bind_ceed::CeedTransposeMode, 3179df49d7eSJed Brown emode as bind_ceed::CeedEvalMode, 3189df49d7eSJed Brown ); 3199df49d7eSJed Brown let ierr = 3209df49d7eSJed Brown unsafe { bind_ceed::CeedBasisApply(self.ptr, nelem, tmode, emode, u.ptr, v.ptr) }; 3219df49d7eSJed Brown self.ceed.check_error(ierr) 3229df49d7eSJed Brown } 3239df49d7eSJed Brown 3249df49d7eSJed Brown /// Returns the dimension for given CeedBasis 3259df49d7eSJed Brown /// 3269df49d7eSJed Brown /// ``` 3279df49d7eSJed Brown /// # use libceed::prelude::*; 328c68be7a2SJeremy L Thompson /// # fn main() -> Result<(), libceed::CeedError> { 3299df49d7eSJed Brown /// # let ceed = libceed::Ceed::default_init(); 3309df49d7eSJed Brown /// let dim = 2; 331c68be7a2SJeremy L Thompson /// let b = ceed.basis_tensor_H1_Lagrange(dim, 1, 3, 4, QuadMode::Gauss)?; 3329df49d7eSJed Brown /// 3339df49d7eSJed Brown /// let d = b.dimension(); 3349df49d7eSJed Brown /// assert_eq!(d, dim, "Incorrect dimension"); 335c68be7a2SJeremy L Thompson /// # Ok(()) 336c68be7a2SJeremy L Thompson /// # } 3379df49d7eSJed Brown /// ``` 3389df49d7eSJed Brown pub fn dimension(&self) -> usize { 3399df49d7eSJed Brown let mut dim = 0; 3409df49d7eSJed Brown unsafe { bind_ceed::CeedBasisGetDimension(self.ptr, &mut dim) }; 3419df49d7eSJed Brown usize::try_from(dim).unwrap() 3429df49d7eSJed Brown } 3439df49d7eSJed Brown 3449df49d7eSJed Brown /// Returns number of components for given CeedBasis 3459df49d7eSJed Brown /// 3469df49d7eSJed Brown /// ``` 3479df49d7eSJed Brown /// # use libceed::prelude::*; 348c68be7a2SJeremy L Thompson /// # fn main() -> Result<(), libceed::CeedError> { 3499df49d7eSJed Brown /// # let ceed = libceed::Ceed::default_init(); 3509df49d7eSJed Brown /// let ncomp = 2; 351c68be7a2SJeremy L Thompson /// let b = ceed.basis_tensor_H1_Lagrange(1, ncomp, 3, 4, QuadMode::Gauss)?; 3529df49d7eSJed Brown /// 3539df49d7eSJed Brown /// let n = b.num_components(); 3549df49d7eSJed Brown /// assert_eq!(n, ncomp, "Incorrect number of components"); 355c68be7a2SJeremy L Thompson /// # Ok(()) 356c68be7a2SJeremy L Thompson /// # } 3579df49d7eSJed Brown /// ``` 3589df49d7eSJed Brown pub fn num_components(&self) -> usize { 3599df49d7eSJed Brown let mut ncomp = 0; 3609df49d7eSJed Brown unsafe { bind_ceed::CeedBasisGetNumComponents(self.ptr, &mut ncomp) }; 3619df49d7eSJed Brown usize::try_from(ncomp).unwrap() 3629df49d7eSJed Brown } 3639df49d7eSJed Brown 3649df49d7eSJed Brown /// Returns total number of nodes (in dim dimensions) of a CeedBasis 3659df49d7eSJed Brown /// 3669df49d7eSJed Brown /// ``` 3679df49d7eSJed Brown /// # use libceed::prelude::*; 368c68be7a2SJeremy L Thompson /// # fn main() -> Result<(), libceed::CeedError> { 3699df49d7eSJed Brown /// # let ceed = libceed::Ceed::default_init(); 3709df49d7eSJed Brown /// let p = 3; 371c68be7a2SJeremy L Thompson /// let b = ceed.basis_tensor_H1_Lagrange(2, 1, p, 4, QuadMode::Gauss)?; 3729df49d7eSJed Brown /// 3739df49d7eSJed Brown /// let nnodes = b.num_nodes(); 3749df49d7eSJed Brown /// assert_eq!(nnodes, p * p, "Incorrect number of nodes"); 375c68be7a2SJeremy L Thompson /// # Ok(()) 376c68be7a2SJeremy L Thompson /// # } 3779df49d7eSJed Brown /// ``` 3789df49d7eSJed Brown pub fn num_nodes(&self) -> usize { 3799df49d7eSJed Brown let mut nnodes = 0; 3809df49d7eSJed Brown unsafe { bind_ceed::CeedBasisGetNumNodes(self.ptr, &mut nnodes) }; 3819df49d7eSJed Brown usize::try_from(nnodes).unwrap() 3829df49d7eSJed Brown } 3839df49d7eSJed Brown 3849df49d7eSJed Brown /// Returns total number of quadrature points (in dim dimensions) of a 3859df49d7eSJed Brown /// CeedBasis 3869df49d7eSJed Brown /// 3879df49d7eSJed Brown /// ``` 3889df49d7eSJed Brown /// # use libceed::prelude::*; 389c68be7a2SJeremy L Thompson /// # fn main() -> Result<(), libceed::CeedError> { 3909df49d7eSJed Brown /// # let ceed = libceed::Ceed::default_init(); 3919df49d7eSJed Brown /// let q = 4; 392c68be7a2SJeremy L Thompson /// let b = ceed.basis_tensor_H1_Lagrange(2, 1, 3, q, QuadMode::Gauss)?; 3939df49d7eSJed Brown /// 3949df49d7eSJed Brown /// let nqpts = b.num_quadrature_points(); 3959df49d7eSJed Brown /// assert_eq!(nqpts, q * q, "Incorrect number of quadrature points"); 396c68be7a2SJeremy L Thompson /// # Ok(()) 397c68be7a2SJeremy L Thompson /// # } 3989df49d7eSJed Brown /// ``` 3999df49d7eSJed Brown pub fn num_quadrature_points(&self) -> usize { 4009df49d7eSJed Brown let mut Q = 0; 4019df49d7eSJed Brown unsafe { 4029df49d7eSJed Brown bind_ceed::CeedBasisGetNumQuadraturePoints(self.ptr, &mut Q); 4039df49d7eSJed Brown } 4049df49d7eSJed Brown usize::try_from(Q).unwrap() 4059df49d7eSJed Brown } 4069df49d7eSJed Brown } 4079df49d7eSJed Brown 4089df49d7eSJed Brown // ----------------------------------------------------------------------------- 409