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 } 45e03682afSJeremy L Thompson 46e03682afSJeremy L Thompson /// Check if a BasisOpt is Some 47e03682afSJeremy L Thompson /// 48e03682afSJeremy L Thompson /// ``` 49e03682afSJeremy L Thompson /// # use libceed::prelude::*; 5089d15d5fSJeremy L Thompson /// # fn main() -> Result<()> { 51e03682afSJeremy L Thompson /// # let ceed = libceed::Ceed::default_init(); 52e03682afSJeremy L Thompson /// let b = ceed.basis_tensor_H1_Lagrange(1, 2, 3, 4, QuadMode::Gauss)?; 53e03682afSJeremy L Thompson /// let b_opt = BasisOpt::from(&b); 54e03682afSJeremy L Thompson /// assert!(b_opt.is_some(), "Incorrect BasisOpt"); 55e03682afSJeremy L Thompson /// 56e03682afSJeremy L Thompson /// let b_opt = BasisOpt::Collocated; 57e03682afSJeremy L Thompson /// assert!(!b_opt.is_some(), "Incorrect BasisOpt"); 58e03682afSJeremy L Thompson /// # Ok(()) 59e03682afSJeremy L Thompson /// # } 60e03682afSJeremy L Thompson /// ``` 61e03682afSJeremy L Thompson pub fn is_some(&self) -> bool { 62e03682afSJeremy L Thompson match self { 63e03682afSJeremy L Thompson Self::Some(_) => true, 64e03682afSJeremy L Thompson Self::Collocated => false, 65e03682afSJeremy L Thompson } 66e03682afSJeremy L Thompson } 67e03682afSJeremy L Thompson 68e03682afSJeremy L Thompson /// Check if a BasisOpt is Collocated 69e03682afSJeremy L Thompson /// 70e03682afSJeremy L Thompson /// ``` 71e03682afSJeremy L Thompson /// # use libceed::prelude::*; 7289d15d5fSJeremy L Thompson /// # fn main() -> Result<()> { 73e03682afSJeremy L Thompson /// # let ceed = libceed::Ceed::default_init(); 74e03682afSJeremy L Thompson /// let b = ceed.basis_tensor_H1_Lagrange(1, 2, 3, 4, QuadMode::Gauss)?; 75e03682afSJeremy L Thompson /// let b_opt = BasisOpt::from(&b); 76e03682afSJeremy L Thompson /// assert!(!b_opt.is_collocated(), "Incorrect BasisOpt"); 77e03682afSJeremy L Thompson /// 78e03682afSJeremy L Thompson /// let b_opt = BasisOpt::Collocated; 79e03682afSJeremy L Thompson /// assert!(b_opt.is_collocated(), "Incorrect BasisOpt"); 80e03682afSJeremy L Thompson /// # Ok(()) 81e03682afSJeremy L Thompson /// # } 82e03682afSJeremy L Thompson /// ``` 83e03682afSJeremy L Thompson pub fn is_collocated(&self) -> bool { 84e03682afSJeremy L Thompson match self { 85e03682afSJeremy L Thompson Self::Some(_) => false, 86e03682afSJeremy L Thompson Self::Collocated => true, 87e03682afSJeremy L Thompson } 88e03682afSJeremy 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 pub(crate) ptr: bind_ceed::CeedBasis, 971142270cSJeremy L Thompson _lifeline: PhantomData<&'a ()>, 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::*; 12189d15d5fSJeremy L Thompson /// # fn main() -> Result<()> { 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)?; 1791142270cSJeremy L Thompson Ok(Self { 1801142270cSJeremy L Thompson ptr, 1811142270cSJeremy L Thompson _lifeline: PhantomData, 1821142270cSJeremy L Thompson }) 1839df49d7eSJed Brown } 1849df49d7eSJed Brown 1859df49d7eSJed Brown pub fn create_tensor_H1_Lagrange( 1869df49d7eSJed Brown ceed: &'a crate::Ceed, 1879df49d7eSJed Brown dim: usize, 1889df49d7eSJed Brown ncomp: usize, 1899df49d7eSJed Brown P: usize, 1909df49d7eSJed Brown Q: usize, 1919df49d7eSJed Brown qmode: crate::QuadMode, 1929df49d7eSJed Brown ) -> crate::Result<Self> { 1939df49d7eSJed Brown let mut ptr = std::ptr::null_mut(); 1949df49d7eSJed Brown let (dim, ncomp, P, Q, qmode) = ( 1959df49d7eSJed Brown i32::try_from(dim).unwrap(), 1969df49d7eSJed Brown i32::try_from(ncomp).unwrap(), 1979df49d7eSJed Brown i32::try_from(P).unwrap(), 1989df49d7eSJed Brown i32::try_from(Q).unwrap(), 1999df49d7eSJed Brown qmode as bind_ceed::CeedQuadMode, 2009df49d7eSJed Brown ); 2019df49d7eSJed Brown let ierr = unsafe { 2029df49d7eSJed Brown bind_ceed::CeedBasisCreateTensorH1Lagrange(ceed.ptr, dim, ncomp, P, Q, qmode, &mut ptr) 2039df49d7eSJed Brown }; 2049df49d7eSJed Brown ceed.check_error(ierr)?; 2051142270cSJeremy L Thompson Ok(Self { 2061142270cSJeremy L Thompson ptr, 2071142270cSJeremy L Thompson _lifeline: PhantomData, 2081142270cSJeremy L Thompson }) 2099df49d7eSJed Brown } 2109df49d7eSJed Brown 2119df49d7eSJed Brown pub fn create_H1( 2129df49d7eSJed Brown ceed: &'a crate::Ceed, 2139df49d7eSJed Brown topo: crate::ElemTopology, 2149df49d7eSJed Brown ncomp: usize, 2159df49d7eSJed Brown nnodes: usize, 2169df49d7eSJed Brown nqpts: usize, 21780a9ef05SNatalie Beams interp: &[crate::Scalar], 21880a9ef05SNatalie Beams grad: &[crate::Scalar], 21980a9ef05SNatalie Beams qref: &[crate::Scalar], 22080a9ef05SNatalie Beams qweight: &[crate::Scalar], 2219df49d7eSJed Brown ) -> crate::Result<Self> { 2229df49d7eSJed Brown let mut ptr = std::ptr::null_mut(); 2239df49d7eSJed Brown let (topo, ncomp, nnodes, nqpts) = ( 2249df49d7eSJed Brown topo as bind_ceed::CeedElemTopology, 2259df49d7eSJed Brown i32::try_from(ncomp).unwrap(), 2269df49d7eSJed Brown i32::try_from(nnodes).unwrap(), 2279df49d7eSJed Brown i32::try_from(nqpts).unwrap(), 2289df49d7eSJed Brown ); 2299df49d7eSJed Brown let ierr = unsafe { 2309df49d7eSJed Brown bind_ceed::CeedBasisCreateH1( 2319df49d7eSJed Brown ceed.ptr, 2329df49d7eSJed Brown topo, 2339df49d7eSJed Brown ncomp, 2349df49d7eSJed Brown nnodes, 2359df49d7eSJed Brown nqpts, 2369df49d7eSJed Brown interp.as_ptr(), 2379df49d7eSJed Brown grad.as_ptr(), 2389df49d7eSJed Brown qref.as_ptr(), 2399df49d7eSJed Brown qweight.as_ptr(), 2409df49d7eSJed Brown &mut ptr, 2419df49d7eSJed Brown ) 2429df49d7eSJed Brown }; 2439df49d7eSJed Brown ceed.check_error(ierr)?; 2441142270cSJeremy L Thompson Ok(Self { 2451142270cSJeremy L Thompson ptr, 2461142270cSJeremy L Thompson _lifeline: PhantomData, 2471142270cSJeremy L Thompson }) 2481142270cSJeremy L Thompson } 2491142270cSJeremy L Thompson 2501142270cSJeremy L Thompson // Error handling 2511142270cSJeremy L Thompson #[doc(hidden)] 2521142270cSJeremy L Thompson fn check_error(&self, ierr: i32) -> crate::Result<i32> { 2531142270cSJeremy L Thompson let mut ptr = std::ptr::null_mut(); 2541142270cSJeremy L Thompson unsafe { 2551142270cSJeremy L Thompson bind_ceed::CeedBasisGetCeed(self.ptr, &mut ptr); 2561142270cSJeremy L Thompson } 2571142270cSJeremy L Thompson crate::check_error(ptr, ierr) 2589df49d7eSJed Brown } 2599df49d7eSJed Brown 2609df49d7eSJed Brown /// Apply basis evaluation from nodes to quadrature points or vice versa 2619df49d7eSJed Brown /// 2629df49d7eSJed Brown /// * `nelem` - The number of elements to apply the basis evaluation to 2639df49d7eSJed Brown /// * `tmode` - `TrasposeMode::NoTranspose` to evaluate from nodes to 2649df49d7eSJed Brown /// quadrature points, `TransposeMode::Transpose` to apply the 2659df49d7eSJed Brown /// transpose, mapping from quadrature points to nodes 2669df49d7eSJed Brown /// * `emode` - `EvalMode::None` to use values directly, `EvalMode::Interp` 2679df49d7eSJed Brown /// to use interpolated values, `EvalMode::Grad` to use 2689df49d7eSJed Brown /// gradients, `EvalMode::Weight` to use quadrature weights 2699df49d7eSJed Brown /// * `u` - Input Vector 2709df49d7eSJed Brown /// * `v` - Output Vector 2719df49d7eSJed Brown /// 2729df49d7eSJed Brown /// ``` 2739df49d7eSJed Brown /// # use libceed::prelude::*; 27489d15d5fSJeremy L Thompson /// # fn main() -> Result<()> { 2759df49d7eSJed Brown /// # let ceed = libceed::Ceed::default_init(); 2769df49d7eSJed Brown /// const Q: usize = 6; 277c68be7a2SJeremy L Thompson /// let bu = ceed.basis_tensor_H1_Lagrange(1, 1, Q, Q, QuadMode::GaussLobatto)?; 278c68be7a2SJeremy L Thompson /// let bx = ceed.basis_tensor_H1_Lagrange(1, 1, 2, Q, QuadMode::Gauss)?; 2799df49d7eSJed Brown /// 280c68be7a2SJeremy L Thompson /// let x_corners = ceed.vector_from_slice(&[-1., 1.])?; 281c68be7a2SJeremy L Thompson /// let mut x_qpts = ceed.vector(Q)?; 282c68be7a2SJeremy L Thompson /// let mut x_nodes = ceed.vector(Q)?; 2839df49d7eSJed Brown /// bx.apply( 2849df49d7eSJed Brown /// 1, 2859df49d7eSJed Brown /// TransposeMode::NoTranspose, 2869df49d7eSJed Brown /// EvalMode::Interp, 2879df49d7eSJed Brown /// &x_corners, 2889df49d7eSJed Brown /// &mut x_nodes, 289*d7f01795SJeremy L Thompson /// )?; 2909df49d7eSJed Brown /// bu.apply( 2919df49d7eSJed Brown /// 1, 2929df49d7eSJed Brown /// TransposeMode::NoTranspose, 2939df49d7eSJed Brown /// EvalMode::Interp, 2949df49d7eSJed Brown /// &x_nodes, 2959df49d7eSJed Brown /// &mut x_qpts, 296*d7f01795SJeremy L Thompson /// )?; 2979df49d7eSJed Brown /// 2989df49d7eSJed Brown /// // Create function x^3 + 1 on Gauss Lobatto points 2999df49d7eSJed Brown /// let mut u_arr = [0.; Q]; 3009df49d7eSJed Brown /// u_arr 3019df49d7eSJed Brown /// .iter_mut() 3029df49d7eSJed Brown /// .zip(x_nodes.view().iter()) 3039df49d7eSJed Brown /// .for_each(|(u, x)| *u = x * x * x + 1.); 304c68be7a2SJeremy L Thompson /// let u = ceed.vector_from_slice(&u_arr)?; 3059df49d7eSJed Brown /// 3069df49d7eSJed Brown /// // Map function to Gauss points 307c68be7a2SJeremy L Thompson /// let mut v = ceed.vector(Q)?; 3089df49d7eSJed Brown /// v.set_value(0.); 309c68be7a2SJeremy L Thompson /// bu.apply(1, TransposeMode::NoTranspose, EvalMode::Interp, &u, &mut v)?; 3109df49d7eSJed Brown /// 3119df49d7eSJed Brown /// // Verify results 3129df49d7eSJed Brown /// v.view() 3139df49d7eSJed Brown /// .iter() 3149df49d7eSJed Brown /// .zip(x_qpts.view().iter()) 3159df49d7eSJed Brown /// .for_each(|(v, x)| { 3169df49d7eSJed Brown /// let true_value = x * x * x + 1.; 31780a9ef05SNatalie Beams /// assert!( 31880a9ef05SNatalie Beams /// (*v - true_value).abs() < 10.0 * libceed::EPSILON, 31980a9ef05SNatalie Beams /// "Incorrect basis application" 32080a9ef05SNatalie Beams /// ); 3219df49d7eSJed Brown /// }); 322c68be7a2SJeremy L Thompson /// # Ok(()) 323c68be7a2SJeremy L Thompson /// # } 3249df49d7eSJed Brown /// ``` 3259df49d7eSJed Brown pub fn apply( 3269df49d7eSJed Brown &self, 3279df49d7eSJed Brown nelem: usize, 3289df49d7eSJed Brown tmode: TransposeMode, 3299df49d7eSJed Brown emode: EvalMode, 3309df49d7eSJed Brown u: &Vector, 3319df49d7eSJed Brown v: &mut Vector, 3329df49d7eSJed Brown ) -> crate::Result<i32> { 3339df49d7eSJed Brown let (nelem, tmode, emode) = ( 3349df49d7eSJed Brown i32::try_from(nelem).unwrap(), 3359df49d7eSJed Brown tmode as bind_ceed::CeedTransposeMode, 3369df49d7eSJed Brown emode as bind_ceed::CeedEvalMode, 3379df49d7eSJed Brown ); 3389df49d7eSJed Brown let ierr = 3399df49d7eSJed Brown unsafe { bind_ceed::CeedBasisApply(self.ptr, nelem, tmode, emode, u.ptr, v.ptr) }; 3401142270cSJeremy L Thompson self.check_error(ierr) 3419df49d7eSJed Brown } 3429df49d7eSJed Brown 3439df49d7eSJed Brown /// Returns the dimension for given CeedBasis 3449df49d7eSJed Brown /// 3459df49d7eSJed Brown /// ``` 3469df49d7eSJed Brown /// # use libceed::prelude::*; 34789d15d5fSJeremy L Thompson /// # fn main() -> Result<()> { 3489df49d7eSJed Brown /// # let ceed = libceed::Ceed::default_init(); 3499df49d7eSJed Brown /// let dim = 2; 350c68be7a2SJeremy L Thompson /// let b = ceed.basis_tensor_H1_Lagrange(dim, 1, 3, 4, QuadMode::Gauss)?; 3519df49d7eSJed Brown /// 3529df49d7eSJed Brown /// let d = b.dimension(); 3539df49d7eSJed Brown /// assert_eq!(d, dim, "Incorrect dimension"); 354c68be7a2SJeremy L Thompson /// # Ok(()) 355c68be7a2SJeremy L Thompson /// # } 3569df49d7eSJed Brown /// ``` 3579df49d7eSJed Brown pub fn dimension(&self) -> usize { 3589df49d7eSJed Brown let mut dim = 0; 3599df49d7eSJed Brown unsafe { bind_ceed::CeedBasisGetDimension(self.ptr, &mut dim) }; 3609df49d7eSJed Brown usize::try_from(dim).unwrap() 3619df49d7eSJed Brown } 3629df49d7eSJed Brown 3639df49d7eSJed Brown /// Returns number of components for given CeedBasis 3649df49d7eSJed Brown /// 3659df49d7eSJed Brown /// ``` 3669df49d7eSJed Brown /// # use libceed::prelude::*; 36789d15d5fSJeremy L Thompson /// # fn main() -> Result<()> { 3689df49d7eSJed Brown /// # let ceed = libceed::Ceed::default_init(); 3699df49d7eSJed Brown /// let ncomp = 2; 370c68be7a2SJeremy L Thompson /// let b = ceed.basis_tensor_H1_Lagrange(1, ncomp, 3, 4, QuadMode::Gauss)?; 3719df49d7eSJed Brown /// 3729df49d7eSJed Brown /// let n = b.num_components(); 3739df49d7eSJed Brown /// assert_eq!(n, ncomp, "Incorrect number of components"); 374c68be7a2SJeremy L Thompson /// # Ok(()) 375c68be7a2SJeremy L Thompson /// # } 3769df49d7eSJed Brown /// ``` 3779df49d7eSJed Brown pub fn num_components(&self) -> usize { 3789df49d7eSJed Brown let mut ncomp = 0; 3799df49d7eSJed Brown unsafe { bind_ceed::CeedBasisGetNumComponents(self.ptr, &mut ncomp) }; 3809df49d7eSJed Brown usize::try_from(ncomp).unwrap() 3819df49d7eSJed Brown } 3829df49d7eSJed Brown 3839df49d7eSJed Brown /// Returns total number of nodes (in dim dimensions) of a CeedBasis 3849df49d7eSJed Brown /// 3859df49d7eSJed Brown /// ``` 3869df49d7eSJed Brown /// # use libceed::prelude::*; 38789d15d5fSJeremy L Thompson /// # fn main() -> Result<()> { 3889df49d7eSJed Brown /// # let ceed = libceed::Ceed::default_init(); 3899df49d7eSJed Brown /// let p = 3; 390c68be7a2SJeremy L Thompson /// let b = ceed.basis_tensor_H1_Lagrange(2, 1, p, 4, QuadMode::Gauss)?; 3919df49d7eSJed Brown /// 3929df49d7eSJed Brown /// let nnodes = b.num_nodes(); 3939df49d7eSJed Brown /// assert_eq!(nnodes, p * p, "Incorrect number of nodes"); 394c68be7a2SJeremy L Thompson /// # Ok(()) 395c68be7a2SJeremy L Thompson /// # } 3969df49d7eSJed Brown /// ``` 3979df49d7eSJed Brown pub fn num_nodes(&self) -> usize { 3989df49d7eSJed Brown let mut nnodes = 0; 3999df49d7eSJed Brown unsafe { bind_ceed::CeedBasisGetNumNodes(self.ptr, &mut nnodes) }; 4009df49d7eSJed Brown usize::try_from(nnodes).unwrap() 4019df49d7eSJed Brown } 4029df49d7eSJed Brown 4039df49d7eSJed Brown /// Returns total number of quadrature points (in dim dimensions) of a 4049df49d7eSJed Brown /// CeedBasis 4059df49d7eSJed Brown /// 4069df49d7eSJed Brown /// ``` 4079df49d7eSJed Brown /// # use libceed::prelude::*; 40889d15d5fSJeremy L Thompson /// # fn main() -> Result<()> { 4099df49d7eSJed Brown /// # let ceed = libceed::Ceed::default_init(); 4109df49d7eSJed Brown /// let q = 4; 411c68be7a2SJeremy L Thompson /// let b = ceed.basis_tensor_H1_Lagrange(2, 1, 3, q, QuadMode::Gauss)?; 4129df49d7eSJed Brown /// 4139df49d7eSJed Brown /// let nqpts = b.num_quadrature_points(); 4149df49d7eSJed Brown /// assert_eq!(nqpts, q * q, "Incorrect number of quadrature points"); 415c68be7a2SJeremy L Thompson /// # Ok(()) 416c68be7a2SJeremy L Thompson /// # } 4179df49d7eSJed Brown /// ``` 4189df49d7eSJed Brown pub fn num_quadrature_points(&self) -> usize { 4199df49d7eSJed Brown let mut Q = 0; 4209df49d7eSJed Brown unsafe { 4219df49d7eSJed Brown bind_ceed::CeedBasisGetNumQuadraturePoints(self.ptr, &mut Q); 4229df49d7eSJed Brown } 4239df49d7eSJed Brown usize::try_from(Q).unwrap() 4249df49d7eSJed Brown } 4259df49d7eSJed Brown } 4269df49d7eSJed Brown 4279df49d7eSJed Brown // ----------------------------------------------------------------------------- 428