1*5aed82e4SJeremy L Thompson // Copyright (c) 2017-2024, 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 119df49d7eSJed Brown use crate::prelude::*; 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 309df49d7eSJed Brown 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 /// ``` 40e03682afSJeremy L Thompson /// # use libceed::prelude::*; 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 /// ``` 62e03682afSJeremy L Thompson /// # use libceed::prelude::*; 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 /// ``` 1119df49d7eSJed Brown /// # use libceed::prelude::*; 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 1379df49d7eSJed Brown pub fn create_tensor_H1( 138594ef120SJeremy L Thompson ceed: &crate::Ceed, 1399df49d7eSJed Brown dim: usize, 1409df49d7eSJed Brown ncomp: usize, 1419df49d7eSJed Brown P1d: usize, 1429df49d7eSJed Brown Q1d: usize, 14380a9ef05SNatalie Beams interp1d: &[crate::Scalar], 14480a9ef05SNatalie Beams grad1d: &[crate::Scalar], 14580a9ef05SNatalie Beams qref1d: &[crate::Scalar], 14680a9ef05SNatalie Beams qweight1d: &[crate::Scalar], 1479df49d7eSJed Brown ) -> crate::Result<Self> { 1489df49d7eSJed Brown let mut ptr = std::ptr::null_mut(); 1499df49d7eSJed Brown let (dim, ncomp, P1d, Q1d) = ( 1509df49d7eSJed Brown i32::try_from(dim).unwrap(), 1519df49d7eSJed Brown i32::try_from(ncomp).unwrap(), 1529df49d7eSJed Brown i32::try_from(P1d).unwrap(), 1539df49d7eSJed Brown i32::try_from(Q1d).unwrap(), 1549df49d7eSJed Brown ); 1559df49d7eSJed Brown let ierr = unsafe { 1569df49d7eSJed Brown bind_ceed::CeedBasisCreateTensorH1( 1579df49d7eSJed Brown ceed.ptr, 1589df49d7eSJed Brown dim, 1599df49d7eSJed Brown ncomp, 1609df49d7eSJed Brown P1d, 1619df49d7eSJed Brown Q1d, 1629df49d7eSJed Brown interp1d.as_ptr(), 1639df49d7eSJed Brown grad1d.as_ptr(), 1649df49d7eSJed Brown qref1d.as_ptr(), 1659df49d7eSJed Brown qweight1d.as_ptr(), 1669df49d7eSJed Brown &mut ptr, 1679df49d7eSJed Brown ) 1689df49d7eSJed Brown }; 1699df49d7eSJed Brown ceed.check_error(ierr)?; 1701142270cSJeremy L Thompson Ok(Self { 1711142270cSJeremy L Thompson ptr, 1721142270cSJeremy L Thompson _lifeline: PhantomData, 1731142270cSJeremy L Thompson }) 1749df49d7eSJed Brown } 1759df49d7eSJed Brown 1769df49d7eSJed Brown pub fn create_tensor_H1_Lagrange( 177594ef120SJeremy L Thompson ceed: &crate::Ceed, 1789df49d7eSJed Brown dim: usize, 1799df49d7eSJed Brown ncomp: usize, 1809df49d7eSJed Brown P: usize, 1819df49d7eSJed Brown Q: usize, 1829df49d7eSJed Brown qmode: crate::QuadMode, 1839df49d7eSJed Brown ) -> crate::Result<Self> { 1849df49d7eSJed Brown let mut ptr = std::ptr::null_mut(); 1859df49d7eSJed Brown let (dim, ncomp, P, Q, qmode) = ( 1869df49d7eSJed Brown i32::try_from(dim).unwrap(), 1879df49d7eSJed Brown i32::try_from(ncomp).unwrap(), 1889df49d7eSJed Brown i32::try_from(P).unwrap(), 1899df49d7eSJed Brown i32::try_from(Q).unwrap(), 1909df49d7eSJed Brown qmode as bind_ceed::CeedQuadMode, 1919df49d7eSJed Brown ); 1929df49d7eSJed Brown let ierr = unsafe { 1939df49d7eSJed Brown bind_ceed::CeedBasisCreateTensorH1Lagrange(ceed.ptr, dim, ncomp, P, Q, qmode, &mut ptr) 1949df49d7eSJed Brown }; 1959df49d7eSJed Brown ceed.check_error(ierr)?; 1961142270cSJeremy L Thompson Ok(Self { 1971142270cSJeremy L Thompson ptr, 1981142270cSJeremy L Thompson _lifeline: PhantomData, 1991142270cSJeremy L Thompson }) 2009df49d7eSJed Brown } 2019df49d7eSJed Brown 2029df49d7eSJed Brown pub fn create_H1( 203594ef120SJeremy L Thompson ceed: &crate::Ceed, 2049df49d7eSJed Brown topo: crate::ElemTopology, 2059df49d7eSJed Brown ncomp: usize, 2069df49d7eSJed Brown nnodes: usize, 2079df49d7eSJed Brown nqpts: usize, 20880a9ef05SNatalie Beams interp: &[crate::Scalar], 20980a9ef05SNatalie Beams grad: &[crate::Scalar], 21080a9ef05SNatalie Beams qref: &[crate::Scalar], 21180a9ef05SNatalie Beams qweight: &[crate::Scalar], 2129df49d7eSJed Brown ) -> crate::Result<Self> { 2139df49d7eSJed Brown let mut ptr = std::ptr::null_mut(); 2149df49d7eSJed Brown let (topo, ncomp, nnodes, nqpts) = ( 2159df49d7eSJed Brown topo as bind_ceed::CeedElemTopology, 2169df49d7eSJed Brown i32::try_from(ncomp).unwrap(), 2179df49d7eSJed Brown i32::try_from(nnodes).unwrap(), 2189df49d7eSJed Brown i32::try_from(nqpts).unwrap(), 2199df49d7eSJed Brown ); 2209df49d7eSJed Brown let ierr = unsafe { 2219df49d7eSJed Brown bind_ceed::CeedBasisCreateH1( 2229df49d7eSJed Brown ceed.ptr, 2239df49d7eSJed Brown topo, 2249df49d7eSJed Brown ncomp, 2259df49d7eSJed Brown nnodes, 2269df49d7eSJed Brown nqpts, 2279df49d7eSJed Brown interp.as_ptr(), 2289df49d7eSJed Brown grad.as_ptr(), 2299df49d7eSJed Brown qref.as_ptr(), 2309df49d7eSJed Brown qweight.as_ptr(), 2319df49d7eSJed Brown &mut ptr, 2329df49d7eSJed Brown ) 2339df49d7eSJed Brown }; 2349df49d7eSJed Brown ceed.check_error(ierr)?; 2351142270cSJeremy L Thompson Ok(Self { 2361142270cSJeremy L Thompson ptr, 2371142270cSJeremy L Thompson _lifeline: PhantomData, 2381142270cSJeremy L Thompson }) 2391142270cSJeremy L Thompson } 2401142270cSJeremy L Thompson 24111b88ddaSSebastian Grimberg pub fn create_Hdiv( 24211b88ddaSSebastian Grimberg ceed: &crate::Ceed, 24311b88ddaSSebastian Grimberg topo: crate::ElemTopology, 24411b88ddaSSebastian Grimberg ncomp: usize, 24511b88ddaSSebastian Grimberg nnodes: usize, 24611b88ddaSSebastian Grimberg nqpts: usize, 24711b88ddaSSebastian Grimberg interp: &[crate::Scalar], 24811b88ddaSSebastian Grimberg div: &[crate::Scalar], 24911b88ddaSSebastian Grimberg qref: &[crate::Scalar], 25011b88ddaSSebastian Grimberg qweight: &[crate::Scalar], 25111b88ddaSSebastian Grimberg ) -> crate::Result<Self> { 25211b88ddaSSebastian Grimberg let mut ptr = std::ptr::null_mut(); 25311b88ddaSSebastian Grimberg let (topo, ncomp, nnodes, nqpts) = ( 25411b88ddaSSebastian Grimberg topo as bind_ceed::CeedElemTopology, 25511b88ddaSSebastian Grimberg i32::try_from(ncomp).unwrap(), 25611b88ddaSSebastian Grimberg i32::try_from(nnodes).unwrap(), 25711b88ddaSSebastian Grimberg i32::try_from(nqpts).unwrap(), 25811b88ddaSSebastian Grimberg ); 25911b88ddaSSebastian Grimberg let ierr = unsafe { 26011b88ddaSSebastian Grimberg bind_ceed::CeedBasisCreateHdiv( 26111b88ddaSSebastian Grimberg ceed.ptr, 26211b88ddaSSebastian Grimberg topo, 26311b88ddaSSebastian Grimberg ncomp, 26411b88ddaSSebastian Grimberg nnodes, 26511b88ddaSSebastian Grimberg nqpts, 26611b88ddaSSebastian Grimberg interp.as_ptr(), 26711b88ddaSSebastian Grimberg div.as_ptr(), 26811b88ddaSSebastian Grimberg qref.as_ptr(), 26911b88ddaSSebastian Grimberg qweight.as_ptr(), 27011b88ddaSSebastian Grimberg &mut ptr, 27111b88ddaSSebastian Grimberg ) 27211b88ddaSSebastian Grimberg }; 27311b88ddaSSebastian Grimberg ceed.check_error(ierr)?; 27411b88ddaSSebastian Grimberg Ok(Self { 27511b88ddaSSebastian Grimberg ptr, 27611b88ddaSSebastian Grimberg _lifeline: PhantomData, 27711b88ddaSSebastian Grimberg }) 27811b88ddaSSebastian Grimberg } 27911b88ddaSSebastian Grimberg 28011b88ddaSSebastian Grimberg pub fn create_Hcurl( 28111b88ddaSSebastian Grimberg ceed: &crate::Ceed, 28211b88ddaSSebastian Grimberg topo: crate::ElemTopology, 28311b88ddaSSebastian Grimberg ncomp: usize, 28411b88ddaSSebastian Grimberg nnodes: usize, 28511b88ddaSSebastian Grimberg nqpts: usize, 28611b88ddaSSebastian Grimberg interp: &[crate::Scalar], 28711b88ddaSSebastian Grimberg curl: &[crate::Scalar], 28811b88ddaSSebastian Grimberg qref: &[crate::Scalar], 28911b88ddaSSebastian Grimberg qweight: &[crate::Scalar], 29011b88ddaSSebastian Grimberg ) -> crate::Result<Self> { 29111b88ddaSSebastian Grimberg let mut ptr = std::ptr::null_mut(); 29211b88ddaSSebastian Grimberg let (topo, ncomp, nnodes, nqpts) = ( 29311b88ddaSSebastian Grimberg topo as bind_ceed::CeedElemTopology, 29411b88ddaSSebastian Grimberg i32::try_from(ncomp).unwrap(), 29511b88ddaSSebastian Grimberg i32::try_from(nnodes).unwrap(), 29611b88ddaSSebastian Grimberg i32::try_from(nqpts).unwrap(), 29711b88ddaSSebastian Grimberg ); 29811b88ddaSSebastian Grimberg let ierr = unsafe { 29911b88ddaSSebastian Grimberg bind_ceed::CeedBasisCreateHcurl( 30011b88ddaSSebastian Grimberg ceed.ptr, 30111b88ddaSSebastian Grimberg topo, 30211b88ddaSSebastian Grimberg ncomp, 30311b88ddaSSebastian Grimberg nnodes, 30411b88ddaSSebastian Grimberg nqpts, 30511b88ddaSSebastian Grimberg interp.as_ptr(), 30611b88ddaSSebastian Grimberg curl.as_ptr(), 30711b88ddaSSebastian Grimberg qref.as_ptr(), 30811b88ddaSSebastian Grimberg qweight.as_ptr(), 30911b88ddaSSebastian Grimberg &mut ptr, 31011b88ddaSSebastian Grimberg ) 31111b88ddaSSebastian Grimberg }; 31211b88ddaSSebastian Grimberg ceed.check_error(ierr)?; 31311b88ddaSSebastian Grimberg Ok(Self { 31411b88ddaSSebastian Grimberg ptr, 31511b88ddaSSebastian Grimberg _lifeline: PhantomData, 31611b88ddaSSebastian Grimberg }) 31711b88ddaSSebastian Grimberg } 31811b88ddaSSebastian Grimberg 3191142270cSJeremy L Thompson // Error handling 3201142270cSJeremy L Thompson #[doc(hidden)] 3211142270cSJeremy L Thompson fn check_error(&self, ierr: i32) -> crate::Result<i32> { 3221142270cSJeremy L Thompson let mut ptr = std::ptr::null_mut(); 3231142270cSJeremy L Thompson unsafe { 3241142270cSJeremy L Thompson bind_ceed::CeedBasisGetCeed(self.ptr, &mut ptr); 3251142270cSJeremy L Thompson } 3261142270cSJeremy L Thompson crate::check_error(ptr, ierr) 3279df49d7eSJed Brown } 3289df49d7eSJed Brown 3299df49d7eSJed Brown /// Apply basis evaluation from nodes to quadrature points or vice versa 3309df49d7eSJed Brown /// 3319df49d7eSJed Brown /// * `nelem` - The number of elements to apply the basis evaluation to 3329df49d7eSJed Brown /// * `tmode` - `TrasposeMode::NoTranspose` to evaluate from nodes to 3339df49d7eSJed Brown /// quadrature points, `TransposeMode::Transpose` to apply the 3349df49d7eSJed Brown /// transpose, mapping from quadrature points to nodes 3359df49d7eSJed Brown /// * `emode` - `EvalMode::None` to use values directly, `EvalMode::Interp` 3369df49d7eSJed Brown /// to use interpolated values, `EvalMode::Grad` to use 3379df49d7eSJed Brown /// gradients, `EvalMode::Weight` to use quadrature weights 3389df49d7eSJed Brown /// * `u` - Input Vector 3399df49d7eSJed Brown /// * `v` - Output Vector 3409df49d7eSJed Brown /// 3419df49d7eSJed Brown /// ``` 3429df49d7eSJed Brown /// # use libceed::prelude::*; 3434d27c890SJeremy L Thompson /// # fn main() -> libceed::Result<()> { 3449df49d7eSJed Brown /// # let ceed = libceed::Ceed::default_init(); 3459df49d7eSJed Brown /// const Q: usize = 6; 346c68be7a2SJeremy L Thompson /// let bu = ceed.basis_tensor_H1_Lagrange(1, 1, Q, Q, QuadMode::GaussLobatto)?; 347c68be7a2SJeremy L Thompson /// let bx = ceed.basis_tensor_H1_Lagrange(1, 1, 2, Q, QuadMode::Gauss)?; 3489df49d7eSJed Brown /// 349c68be7a2SJeremy L Thompson /// let x_corners = ceed.vector_from_slice(&[-1., 1.])?; 350c68be7a2SJeremy L Thompson /// let mut x_qpts = ceed.vector(Q)?; 351c68be7a2SJeremy L Thompson /// let mut x_nodes = ceed.vector(Q)?; 3529df49d7eSJed Brown /// bx.apply( 3539df49d7eSJed Brown /// 1, 3549df49d7eSJed Brown /// TransposeMode::NoTranspose, 3559df49d7eSJed Brown /// EvalMode::Interp, 3569df49d7eSJed Brown /// &x_corners, 3579df49d7eSJed Brown /// &mut x_nodes, 358d7f01795SJeremy L Thompson /// )?; 3599df49d7eSJed Brown /// bu.apply( 3609df49d7eSJed Brown /// 1, 3619df49d7eSJed Brown /// TransposeMode::NoTranspose, 3629df49d7eSJed Brown /// EvalMode::Interp, 3639df49d7eSJed Brown /// &x_nodes, 3649df49d7eSJed Brown /// &mut x_qpts, 365d7f01795SJeremy L Thompson /// )?; 3669df49d7eSJed Brown /// 3679df49d7eSJed Brown /// // Create function x^3 + 1 on Gauss Lobatto points 3689df49d7eSJed Brown /// let mut u_arr = [0.; Q]; 3699df49d7eSJed Brown /// u_arr 3709df49d7eSJed Brown /// .iter_mut() 371e78171edSJeremy L Thompson /// .zip(x_nodes.view()?.iter()) 3729df49d7eSJed Brown /// .for_each(|(u, x)| *u = x * x * x + 1.); 373c68be7a2SJeremy L Thompson /// let u = ceed.vector_from_slice(&u_arr)?; 3749df49d7eSJed Brown /// 3759df49d7eSJed Brown /// // Map function to Gauss points 376c68be7a2SJeremy L Thompson /// let mut v = ceed.vector(Q)?; 3779df49d7eSJed Brown /// v.set_value(0.); 378c68be7a2SJeremy L Thompson /// bu.apply(1, TransposeMode::NoTranspose, EvalMode::Interp, &u, &mut v)?; 3799df49d7eSJed Brown /// 3809df49d7eSJed Brown /// // Verify results 381e78171edSJeremy L Thompson /// v.view()? 3829df49d7eSJed Brown /// .iter() 383e78171edSJeremy L Thompson /// .zip(x_qpts.view()?.iter()) 3849df49d7eSJed Brown /// .for_each(|(v, x)| { 3859df49d7eSJed Brown /// let true_value = x * x * x + 1.; 38680a9ef05SNatalie Beams /// assert!( 38780a9ef05SNatalie Beams /// (*v - true_value).abs() < 10.0 * libceed::EPSILON, 38880a9ef05SNatalie Beams /// "Incorrect basis application" 38980a9ef05SNatalie Beams /// ); 3909df49d7eSJed Brown /// }); 391c68be7a2SJeremy L Thompson /// # Ok(()) 392c68be7a2SJeremy L Thompson /// # } 3939df49d7eSJed Brown /// ``` 3949df49d7eSJed Brown pub fn apply( 3959df49d7eSJed Brown &self, 3969df49d7eSJed Brown nelem: usize, 3979df49d7eSJed Brown tmode: TransposeMode, 3989df49d7eSJed Brown emode: EvalMode, 3999df49d7eSJed Brown u: &Vector, 4009df49d7eSJed Brown v: &mut Vector, 4019df49d7eSJed Brown ) -> crate::Result<i32> { 4029df49d7eSJed Brown let (nelem, tmode, emode) = ( 4039df49d7eSJed Brown i32::try_from(nelem).unwrap(), 4049df49d7eSJed Brown tmode as bind_ceed::CeedTransposeMode, 4059df49d7eSJed Brown emode as bind_ceed::CeedEvalMode, 4069df49d7eSJed Brown ); 4079df49d7eSJed Brown let ierr = 4089df49d7eSJed Brown unsafe { bind_ceed::CeedBasisApply(self.ptr, nelem, tmode, emode, u.ptr, v.ptr) }; 4091142270cSJeremy L Thompson self.check_error(ierr) 4109df49d7eSJed Brown } 4119df49d7eSJed Brown 4127ed177dbSJed Brown /// Returns the dimension for given Basis 4139df49d7eSJed Brown /// 4149df49d7eSJed Brown /// ``` 4159df49d7eSJed Brown /// # use libceed::prelude::*; 4164d27c890SJeremy L Thompson /// # fn main() -> libceed::Result<()> { 4179df49d7eSJed Brown /// # let ceed = libceed::Ceed::default_init(); 4189df49d7eSJed Brown /// let dim = 2; 419c68be7a2SJeremy L Thompson /// let b = ceed.basis_tensor_H1_Lagrange(dim, 1, 3, 4, QuadMode::Gauss)?; 4209df49d7eSJed Brown /// 4219df49d7eSJed Brown /// let d = b.dimension(); 4229df49d7eSJed Brown /// assert_eq!(d, dim, "Incorrect dimension"); 423c68be7a2SJeremy L Thompson /// # Ok(()) 424c68be7a2SJeremy L Thompson /// # } 4259df49d7eSJed Brown /// ``` 4269df49d7eSJed Brown pub fn dimension(&self) -> usize { 4279df49d7eSJed Brown let mut dim = 0; 4289df49d7eSJed Brown unsafe { bind_ceed::CeedBasisGetDimension(self.ptr, &mut dim) }; 4299df49d7eSJed Brown usize::try_from(dim).unwrap() 4309df49d7eSJed Brown } 4319df49d7eSJed Brown 4327ed177dbSJed Brown /// Returns number of components for given Basis 4339df49d7eSJed Brown /// 4349df49d7eSJed Brown /// ``` 4359df49d7eSJed Brown /// # use libceed::prelude::*; 4364d27c890SJeremy L Thompson /// # fn main() -> libceed::Result<()> { 4379df49d7eSJed Brown /// # let ceed = libceed::Ceed::default_init(); 4389df49d7eSJed Brown /// let ncomp = 2; 439c68be7a2SJeremy L Thompson /// let b = ceed.basis_tensor_H1_Lagrange(1, ncomp, 3, 4, QuadMode::Gauss)?; 4409df49d7eSJed Brown /// 4419df49d7eSJed Brown /// let n = b.num_components(); 4429df49d7eSJed Brown /// assert_eq!(n, ncomp, "Incorrect number of components"); 443c68be7a2SJeremy L Thompson /// # Ok(()) 444c68be7a2SJeremy L Thompson /// # } 4459df49d7eSJed Brown /// ``` 4469df49d7eSJed Brown pub fn num_components(&self) -> usize { 4479df49d7eSJed Brown let mut ncomp = 0; 4489df49d7eSJed Brown unsafe { bind_ceed::CeedBasisGetNumComponents(self.ptr, &mut ncomp) }; 4499df49d7eSJed Brown usize::try_from(ncomp).unwrap() 4509df49d7eSJed Brown } 4519df49d7eSJed Brown 4527ed177dbSJed Brown /// Returns total number of nodes (in dim dimensions) of a Basis 4539df49d7eSJed Brown /// 4549df49d7eSJed Brown /// ``` 4559df49d7eSJed Brown /// # use libceed::prelude::*; 4564d27c890SJeremy L Thompson /// # fn main() -> libceed::Result<()> { 4579df49d7eSJed Brown /// # let ceed = libceed::Ceed::default_init(); 4589df49d7eSJed Brown /// let p = 3; 459c68be7a2SJeremy L Thompson /// let b = ceed.basis_tensor_H1_Lagrange(2, 1, p, 4, QuadMode::Gauss)?; 4609df49d7eSJed Brown /// 4619df49d7eSJed Brown /// let nnodes = b.num_nodes(); 4629df49d7eSJed Brown /// assert_eq!(nnodes, p * p, "Incorrect number of nodes"); 463c68be7a2SJeremy L Thompson /// # Ok(()) 464c68be7a2SJeremy L Thompson /// # } 4659df49d7eSJed Brown /// ``` 4669df49d7eSJed Brown pub fn num_nodes(&self) -> usize { 4679df49d7eSJed Brown let mut nnodes = 0; 4689df49d7eSJed Brown unsafe { bind_ceed::CeedBasisGetNumNodes(self.ptr, &mut nnodes) }; 4699df49d7eSJed Brown usize::try_from(nnodes).unwrap() 4709df49d7eSJed Brown } 4719df49d7eSJed Brown 4729df49d7eSJed Brown /// Returns total number of quadrature points (in dim dimensions) of a 4737ed177dbSJed Brown /// Basis 4749df49d7eSJed Brown /// 4759df49d7eSJed Brown /// ``` 4769df49d7eSJed Brown /// # use libceed::prelude::*; 4774d27c890SJeremy L Thompson /// # fn main() -> libceed::Result<()> { 4789df49d7eSJed Brown /// # let ceed = libceed::Ceed::default_init(); 4799df49d7eSJed Brown /// let q = 4; 480c68be7a2SJeremy L Thompson /// let b = ceed.basis_tensor_H1_Lagrange(2, 1, 3, q, QuadMode::Gauss)?; 4819df49d7eSJed Brown /// 4829df49d7eSJed Brown /// let nqpts = b.num_quadrature_points(); 4839df49d7eSJed Brown /// assert_eq!(nqpts, q * q, "Incorrect number of quadrature points"); 484c68be7a2SJeremy L Thompson /// # Ok(()) 485c68be7a2SJeremy L Thompson /// # } 4869df49d7eSJed Brown /// ``` 4879df49d7eSJed Brown pub fn num_quadrature_points(&self) -> usize { 4889df49d7eSJed Brown let mut Q = 0; 4899df49d7eSJed Brown unsafe { 4909df49d7eSJed Brown bind_ceed::CeedBasisGetNumQuadraturePoints(self.ptr, &mut Q); 4919df49d7eSJed Brown } 4929df49d7eSJed Brown usize::try_from(Q).unwrap() 4939df49d7eSJed Brown } 49400d548f6SJed Brown 49500d548f6SJed Brown /// Create projection from self to specified Basis. 49600d548f6SJed Brown /// 49700d548f6SJed Brown /// Both bases must have the same quadrature space. The input bases need not 49800d548f6SJed Brown /// be nested as function spaces; this interface solves a least squares 499b748b478SJeremy L Thompson /// problem to find a representation in the `to` basis that agrees at 50000d548f6SJed Brown /// quadrature points with the origin basis. Since the bases need not be 50100d548f6SJed Brown /// Lagrange, the resulting projection "basis" will have empty quadrature 50200d548f6SJed Brown /// points and weights. 50300d548f6SJed Brown /// 50400d548f6SJed Brown /// ``` 50500d548f6SJed Brown /// # use libceed::prelude::*; 50600d548f6SJed Brown /// # fn main() -> libceed::Result<()> { 50700d548f6SJed Brown /// # let ceed = libceed::Ceed::default_init(); 50800d548f6SJed Brown /// let coarse = ceed.basis_tensor_H1_Lagrange(1, 1, 2, 3, QuadMode::Gauss)?; 50900d548f6SJed Brown /// let fine = ceed.basis_tensor_H1_Lagrange(1, 1, 3, 3, QuadMode::Gauss)?; 51000d548f6SJed Brown /// let proj = coarse.create_projection(&fine)?; 51100d548f6SJed Brown /// let u = ceed.vector_from_slice(&[1., 2.])?; 51200d548f6SJed Brown /// let mut v = ceed.vector(3)?; 51300d548f6SJed Brown /// proj.apply(1, TransposeMode::NoTranspose, EvalMode::Interp, &u, &mut v)?; 51400d548f6SJed Brown /// let expected = [1., 1.5, 2.]; 51500d548f6SJed Brown /// for (a, b) in v.view()?.iter().zip(expected) { 51600d548f6SJed Brown /// assert!( 51700d548f6SJed Brown /// (a - b).abs() < 10.0 * libceed::EPSILON, 51800d548f6SJed Brown /// "Incorrect projection of linear Lagrange to quadratic Lagrange" 51900d548f6SJed Brown /// ); 52000d548f6SJed Brown /// } 52100d548f6SJed Brown /// # Ok(()) 52200d548f6SJed Brown /// # } 52300d548f6SJed Brown /// ``` 52400d548f6SJed Brown pub fn create_projection(&self, to: &Self) -> crate::Result<Self> { 52500d548f6SJed Brown let mut ptr = std::ptr::null_mut(); 52600d548f6SJed Brown let ierr = unsafe { bind_ceed::CeedBasisCreateProjection(self.ptr, to.ptr, &mut ptr) }; 52700d548f6SJed Brown self.check_error(ierr)?; 52800d548f6SJed Brown Ok(Self { 52900d548f6SJed Brown ptr, 53000d548f6SJed Brown _lifeline: PhantomData, 53100d548f6SJed Brown }) 53200d548f6SJed Brown } 5339df49d7eSJed Brown } 5349df49d7eSJed Brown 5359df49d7eSJed Brown // ----------------------------------------------------------------------------- 536