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 // ----------------------------------------------------------------------------- 25*c68be7a2SJeremy 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 } 459df49d7eSJed Brown } 469df49d7eSJed Brown 479df49d7eSJed Brown // ----------------------------------------------------------------------------- 489df49d7eSJed Brown // CeedBasis context wrapper 499df49d7eSJed Brown // ----------------------------------------------------------------------------- 50*c68be7a2SJeremy L Thompson #[derive(Debug)] 519df49d7eSJed Brown pub struct Basis<'a> { 529df49d7eSJed Brown ceed: &'a crate::Ceed, 539df49d7eSJed Brown pub(crate) ptr: bind_ceed::CeedBasis, 549df49d7eSJed Brown } 559df49d7eSJed Brown 569df49d7eSJed Brown // ----------------------------------------------------------------------------- 579df49d7eSJed Brown // Destructor 589df49d7eSJed Brown // ----------------------------------------------------------------------------- 599df49d7eSJed Brown impl<'a> Drop for Basis<'a> { 609df49d7eSJed Brown fn drop(&mut self) { 619df49d7eSJed Brown unsafe { 629df49d7eSJed Brown if self.ptr != bind_ceed::CEED_BASIS_COLLOCATED { 639df49d7eSJed Brown bind_ceed::CeedBasisDestroy(&mut self.ptr); 649df49d7eSJed Brown } 659df49d7eSJed Brown } 669df49d7eSJed Brown } 679df49d7eSJed Brown } 689df49d7eSJed Brown 699df49d7eSJed Brown // ----------------------------------------------------------------------------- 709df49d7eSJed Brown // Display 719df49d7eSJed Brown // ----------------------------------------------------------------------------- 729df49d7eSJed Brown impl<'a> fmt::Display for Basis<'a> { 739df49d7eSJed Brown /// View a Basis 749df49d7eSJed Brown /// 759df49d7eSJed Brown /// ``` 769df49d7eSJed Brown /// # use libceed::prelude::*; 77*c68be7a2SJeremy L Thompson /// # fn main() -> Result<(), libceed::CeedError> { 789df49d7eSJed Brown /// # let ceed = libceed::Ceed::default_init(); 79*c68be7a2SJeremy L Thompson /// let b = ceed.basis_tensor_H1_Lagrange(1, 2, 3, 4, QuadMode::Gauss)?; 809df49d7eSJed Brown /// println!("{}", b); 81*c68be7a2SJeremy L Thompson /// # Ok(()) 82*c68be7a2SJeremy L Thompson /// # } 839df49d7eSJed Brown /// ``` 849df49d7eSJed Brown fn fmt(&self, f: &mut fmt::Formatter) -> fmt::Result { 859df49d7eSJed Brown let mut ptr = std::ptr::null_mut(); 869df49d7eSJed Brown let mut sizeloc = crate::MAX_BUFFER_LENGTH; 879df49d7eSJed Brown let cstring = unsafe { 889df49d7eSJed Brown let file = bind_ceed::open_memstream(&mut ptr, &mut sizeloc); 899df49d7eSJed Brown bind_ceed::CeedBasisView(self.ptr, file); 909df49d7eSJed Brown bind_ceed::fclose(file); 919df49d7eSJed Brown CString::from_raw(ptr) 929df49d7eSJed Brown }; 939df49d7eSJed Brown cstring.to_string_lossy().fmt(f) 949df49d7eSJed Brown } 959df49d7eSJed Brown } 969df49d7eSJed Brown 979df49d7eSJed Brown // ----------------------------------------------------------------------------- 989df49d7eSJed Brown // Implementations 999df49d7eSJed Brown // ----------------------------------------------------------------------------- 1009df49d7eSJed Brown impl<'a> Basis<'a> { 1019df49d7eSJed Brown // Constructors 1029df49d7eSJed Brown pub fn create_tensor_H1( 1039df49d7eSJed Brown ceed: &'a crate::Ceed, 1049df49d7eSJed Brown dim: usize, 1059df49d7eSJed Brown ncomp: usize, 1069df49d7eSJed Brown P1d: usize, 1079df49d7eSJed Brown Q1d: usize, 10880a9ef05SNatalie Beams interp1d: &[crate::Scalar], 10980a9ef05SNatalie Beams grad1d: &[crate::Scalar], 11080a9ef05SNatalie Beams qref1d: &[crate::Scalar], 11180a9ef05SNatalie Beams qweight1d: &[crate::Scalar], 1129df49d7eSJed Brown ) -> crate::Result<Self> { 1139df49d7eSJed Brown let mut ptr = std::ptr::null_mut(); 1149df49d7eSJed Brown let (dim, ncomp, P1d, Q1d) = ( 1159df49d7eSJed Brown i32::try_from(dim).unwrap(), 1169df49d7eSJed Brown i32::try_from(ncomp).unwrap(), 1179df49d7eSJed Brown i32::try_from(P1d).unwrap(), 1189df49d7eSJed Brown i32::try_from(Q1d).unwrap(), 1199df49d7eSJed Brown ); 1209df49d7eSJed Brown let ierr = unsafe { 1219df49d7eSJed Brown bind_ceed::CeedBasisCreateTensorH1( 1229df49d7eSJed Brown ceed.ptr, 1239df49d7eSJed Brown dim, 1249df49d7eSJed Brown ncomp, 1259df49d7eSJed Brown P1d, 1269df49d7eSJed Brown Q1d, 1279df49d7eSJed Brown interp1d.as_ptr(), 1289df49d7eSJed Brown grad1d.as_ptr(), 1299df49d7eSJed Brown qref1d.as_ptr(), 1309df49d7eSJed Brown qweight1d.as_ptr(), 1319df49d7eSJed Brown &mut ptr, 1329df49d7eSJed Brown ) 1339df49d7eSJed Brown }; 1349df49d7eSJed Brown ceed.check_error(ierr)?; 1359df49d7eSJed Brown Ok(Self { ceed, ptr }) 1369df49d7eSJed Brown } 1379df49d7eSJed Brown 1389df49d7eSJed Brown pub fn create_tensor_H1_Lagrange( 1399df49d7eSJed Brown ceed: &'a crate::Ceed, 1409df49d7eSJed Brown dim: usize, 1419df49d7eSJed Brown ncomp: usize, 1429df49d7eSJed Brown P: usize, 1439df49d7eSJed Brown Q: usize, 1449df49d7eSJed Brown qmode: crate::QuadMode, 1459df49d7eSJed Brown ) -> crate::Result<Self> { 1469df49d7eSJed Brown let mut ptr = std::ptr::null_mut(); 1479df49d7eSJed Brown let (dim, ncomp, P, Q, qmode) = ( 1489df49d7eSJed Brown i32::try_from(dim).unwrap(), 1499df49d7eSJed Brown i32::try_from(ncomp).unwrap(), 1509df49d7eSJed Brown i32::try_from(P).unwrap(), 1519df49d7eSJed Brown i32::try_from(Q).unwrap(), 1529df49d7eSJed Brown qmode as bind_ceed::CeedQuadMode, 1539df49d7eSJed Brown ); 1549df49d7eSJed Brown let ierr = unsafe { 1559df49d7eSJed Brown bind_ceed::CeedBasisCreateTensorH1Lagrange(ceed.ptr, dim, ncomp, P, Q, qmode, &mut ptr) 1569df49d7eSJed Brown }; 1579df49d7eSJed Brown ceed.check_error(ierr)?; 1589df49d7eSJed Brown Ok(Self { ceed, ptr }) 1599df49d7eSJed Brown } 1609df49d7eSJed Brown 1619df49d7eSJed Brown pub fn create_H1( 1629df49d7eSJed Brown ceed: &'a crate::Ceed, 1639df49d7eSJed Brown topo: crate::ElemTopology, 1649df49d7eSJed Brown ncomp: usize, 1659df49d7eSJed Brown nnodes: usize, 1669df49d7eSJed Brown nqpts: usize, 16780a9ef05SNatalie Beams interp: &[crate::Scalar], 16880a9ef05SNatalie Beams grad: &[crate::Scalar], 16980a9ef05SNatalie Beams qref: &[crate::Scalar], 17080a9ef05SNatalie Beams qweight: &[crate::Scalar], 1719df49d7eSJed Brown ) -> crate::Result<Self> { 1729df49d7eSJed Brown let mut ptr = std::ptr::null_mut(); 1739df49d7eSJed Brown let (topo, ncomp, nnodes, nqpts) = ( 1749df49d7eSJed Brown topo as bind_ceed::CeedElemTopology, 1759df49d7eSJed Brown i32::try_from(ncomp).unwrap(), 1769df49d7eSJed Brown i32::try_from(nnodes).unwrap(), 1779df49d7eSJed Brown i32::try_from(nqpts).unwrap(), 1789df49d7eSJed Brown ); 1799df49d7eSJed Brown let ierr = unsafe { 1809df49d7eSJed Brown bind_ceed::CeedBasisCreateH1( 1819df49d7eSJed Brown ceed.ptr, 1829df49d7eSJed Brown topo, 1839df49d7eSJed Brown ncomp, 1849df49d7eSJed Brown nnodes, 1859df49d7eSJed Brown nqpts, 1869df49d7eSJed Brown interp.as_ptr(), 1879df49d7eSJed Brown grad.as_ptr(), 1889df49d7eSJed Brown qref.as_ptr(), 1899df49d7eSJed Brown qweight.as_ptr(), 1909df49d7eSJed Brown &mut ptr, 1919df49d7eSJed Brown ) 1929df49d7eSJed Brown }; 1939df49d7eSJed Brown ceed.check_error(ierr)?; 1949df49d7eSJed Brown Ok(Self { ceed, ptr }) 1959df49d7eSJed Brown } 1969df49d7eSJed Brown 1979df49d7eSJed Brown /// Apply basis evaluation from nodes to quadrature points or vice versa 1989df49d7eSJed Brown /// 1999df49d7eSJed Brown /// * `nelem` - The number of elements to apply the basis evaluation to 2009df49d7eSJed Brown /// * `tmode` - `TrasposeMode::NoTranspose` to evaluate from nodes to 2019df49d7eSJed Brown /// quadrature points, `TransposeMode::Transpose` to apply the 2029df49d7eSJed Brown /// transpose, mapping from quadrature points to nodes 2039df49d7eSJed Brown /// * `emode` - `EvalMode::None` to use values directly, `EvalMode::Interp` 2049df49d7eSJed Brown /// to use interpolated values, `EvalMode::Grad` to use 2059df49d7eSJed Brown /// gradients, `EvalMode::Weight` to use quadrature weights 2069df49d7eSJed Brown /// * `u` - Input Vector 2079df49d7eSJed Brown /// * `v` - Output Vector 2089df49d7eSJed Brown /// 2099df49d7eSJed Brown /// ``` 2109df49d7eSJed Brown /// # use libceed::prelude::*; 211*c68be7a2SJeremy L Thompson /// # fn main() -> Result<(), libceed::CeedError> { 2129df49d7eSJed Brown /// # let ceed = libceed::Ceed::default_init(); 2139df49d7eSJed Brown /// const Q: usize = 6; 214*c68be7a2SJeremy L Thompson /// let bu = ceed.basis_tensor_H1_Lagrange(1, 1, Q, Q, QuadMode::GaussLobatto)?; 215*c68be7a2SJeremy L Thompson /// let bx = ceed.basis_tensor_H1_Lagrange(1, 1, 2, Q, QuadMode::Gauss)?; 2169df49d7eSJed Brown /// 217*c68be7a2SJeremy L Thompson /// let x_corners = ceed.vector_from_slice(&[-1., 1.])?; 218*c68be7a2SJeremy L Thompson /// let mut x_qpts = ceed.vector(Q)?; 219*c68be7a2SJeremy L Thompson /// let mut x_nodes = ceed.vector(Q)?; 2209df49d7eSJed Brown /// bx.apply( 2219df49d7eSJed Brown /// 1, 2229df49d7eSJed Brown /// TransposeMode::NoTranspose, 2239df49d7eSJed Brown /// EvalMode::Interp, 2249df49d7eSJed Brown /// &x_corners, 2259df49d7eSJed Brown /// &mut x_nodes, 2269df49d7eSJed Brown /// ); 2279df49d7eSJed Brown /// bu.apply( 2289df49d7eSJed Brown /// 1, 2299df49d7eSJed Brown /// TransposeMode::NoTranspose, 2309df49d7eSJed Brown /// EvalMode::Interp, 2319df49d7eSJed Brown /// &x_nodes, 2329df49d7eSJed Brown /// &mut x_qpts, 2339df49d7eSJed Brown /// ); 2349df49d7eSJed Brown /// 2359df49d7eSJed Brown /// // Create function x^3 + 1 on Gauss Lobatto points 2369df49d7eSJed Brown /// let mut u_arr = [0.; Q]; 2379df49d7eSJed Brown /// u_arr 2389df49d7eSJed Brown /// .iter_mut() 2399df49d7eSJed Brown /// .zip(x_nodes.view().iter()) 2409df49d7eSJed Brown /// .for_each(|(u, x)| *u = x * x * x + 1.); 241*c68be7a2SJeremy L Thompson /// let u = ceed.vector_from_slice(&u_arr)?; 2429df49d7eSJed Brown /// 2439df49d7eSJed Brown /// // Map function to Gauss points 244*c68be7a2SJeremy L Thompson /// let mut v = ceed.vector(Q)?; 2459df49d7eSJed Brown /// v.set_value(0.); 246*c68be7a2SJeremy L Thompson /// bu.apply(1, TransposeMode::NoTranspose, EvalMode::Interp, &u, &mut v)?; 2479df49d7eSJed Brown /// 2489df49d7eSJed Brown /// // Verify results 2499df49d7eSJed Brown /// v.view() 2509df49d7eSJed Brown /// .iter() 2519df49d7eSJed Brown /// .zip(x_qpts.view().iter()) 2529df49d7eSJed Brown /// .for_each(|(v, x)| { 2539df49d7eSJed Brown /// let true_value = x * x * x + 1.; 25480a9ef05SNatalie Beams /// assert!( 25580a9ef05SNatalie Beams /// (*v - true_value).abs() < 10.0 * libceed::EPSILON, 25680a9ef05SNatalie Beams /// "Incorrect basis application" 25780a9ef05SNatalie Beams /// ); 2589df49d7eSJed Brown /// }); 259*c68be7a2SJeremy L Thompson /// # Ok(()) 260*c68be7a2SJeremy L Thompson /// # } 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::*; 284*c68be7a2SJeremy L Thompson /// # fn main() -> Result<(), libceed::CeedError> { 2859df49d7eSJed Brown /// # let ceed = libceed::Ceed::default_init(); 2869df49d7eSJed Brown /// let dim = 2; 287*c68be7a2SJeremy L Thompson /// let b = ceed.basis_tensor_H1_Lagrange(dim, 1, 3, 4, QuadMode::Gauss)?; 2889df49d7eSJed Brown /// 2899df49d7eSJed Brown /// let d = b.dimension(); 2909df49d7eSJed Brown /// assert_eq!(d, dim, "Incorrect dimension"); 291*c68be7a2SJeremy L Thompson /// # Ok(()) 292*c68be7a2SJeremy L Thompson /// # } 2939df49d7eSJed Brown /// ``` 2949df49d7eSJed Brown pub fn dimension(&self) -> usize { 2959df49d7eSJed Brown let mut dim = 0; 2969df49d7eSJed Brown unsafe { bind_ceed::CeedBasisGetDimension(self.ptr, &mut dim) }; 2979df49d7eSJed Brown usize::try_from(dim).unwrap() 2989df49d7eSJed Brown } 2999df49d7eSJed Brown 3009df49d7eSJed Brown /// Returns number of components for given CeedBasis 3019df49d7eSJed Brown /// 3029df49d7eSJed Brown /// ``` 3039df49d7eSJed Brown /// # use libceed::prelude::*; 304*c68be7a2SJeremy L Thompson /// # fn main() -> Result<(), libceed::CeedError> { 3059df49d7eSJed Brown /// # let ceed = libceed::Ceed::default_init(); 3069df49d7eSJed Brown /// let ncomp = 2; 307*c68be7a2SJeremy L Thompson /// let b = ceed.basis_tensor_H1_Lagrange(1, ncomp, 3, 4, QuadMode::Gauss)?; 3089df49d7eSJed Brown /// 3099df49d7eSJed Brown /// let n = b.num_components(); 3109df49d7eSJed Brown /// assert_eq!(n, ncomp, "Incorrect number of components"); 311*c68be7a2SJeremy L Thompson /// # Ok(()) 312*c68be7a2SJeremy L Thompson /// # } 3139df49d7eSJed Brown /// ``` 3149df49d7eSJed Brown pub fn num_components(&self) -> usize { 3159df49d7eSJed Brown let mut ncomp = 0; 3169df49d7eSJed Brown unsafe { bind_ceed::CeedBasisGetNumComponents(self.ptr, &mut ncomp) }; 3179df49d7eSJed Brown usize::try_from(ncomp).unwrap() 3189df49d7eSJed Brown } 3199df49d7eSJed Brown 3209df49d7eSJed Brown /// Returns total number of nodes (in dim dimensions) of a CeedBasis 3219df49d7eSJed Brown /// 3229df49d7eSJed Brown /// ``` 3239df49d7eSJed Brown /// # use libceed::prelude::*; 324*c68be7a2SJeremy L Thompson /// # fn main() -> Result<(), libceed::CeedError> { 3259df49d7eSJed Brown /// # let ceed = libceed::Ceed::default_init(); 3269df49d7eSJed Brown /// let p = 3; 327*c68be7a2SJeremy L Thompson /// let b = ceed.basis_tensor_H1_Lagrange(2, 1, p, 4, QuadMode::Gauss)?; 3289df49d7eSJed Brown /// 3299df49d7eSJed Brown /// let nnodes = b.num_nodes(); 3309df49d7eSJed Brown /// assert_eq!(nnodes, p * p, "Incorrect number of nodes"); 331*c68be7a2SJeremy L Thompson /// # Ok(()) 332*c68be7a2SJeremy L Thompson /// # } 3339df49d7eSJed Brown /// ``` 3349df49d7eSJed Brown pub fn num_nodes(&self) -> usize { 3359df49d7eSJed Brown let mut nnodes = 0; 3369df49d7eSJed Brown unsafe { bind_ceed::CeedBasisGetNumNodes(self.ptr, &mut nnodes) }; 3379df49d7eSJed Brown usize::try_from(nnodes).unwrap() 3389df49d7eSJed Brown } 3399df49d7eSJed Brown 3409df49d7eSJed Brown /// Returns total number of quadrature points (in dim dimensions) of a 3419df49d7eSJed Brown /// CeedBasis 3429df49d7eSJed Brown /// 3439df49d7eSJed Brown /// ``` 3449df49d7eSJed Brown /// # use libceed::prelude::*; 345*c68be7a2SJeremy L Thompson /// # fn main() -> Result<(), libceed::CeedError> { 3469df49d7eSJed Brown /// # let ceed = libceed::Ceed::default_init(); 3479df49d7eSJed Brown /// let q = 4; 348*c68be7a2SJeremy L Thompson /// let b = ceed.basis_tensor_H1_Lagrange(2, 1, 3, q, QuadMode::Gauss)?; 3499df49d7eSJed Brown /// 3509df49d7eSJed Brown /// let nqpts = b.num_quadrature_points(); 3519df49d7eSJed Brown /// assert_eq!(nqpts, q * q, "Incorrect number of quadrature points"); 352*c68be7a2SJeremy L Thompson /// # Ok(()) 353*c68be7a2SJeremy L Thompson /// # } 3549df49d7eSJed Brown /// ``` 3559df49d7eSJed Brown pub fn num_quadrature_points(&self) -> usize { 3569df49d7eSJed Brown let mut Q = 0; 3579df49d7eSJed Brown unsafe { 3589df49d7eSJed Brown bind_ceed::CeedBasisGetNumQuadraturePoints(self.ptr, &mut Q); 3599df49d7eSJed Brown } 3609df49d7eSJed Brown usize::try_from(Q).unwrap() 3619df49d7eSJed Brown } 3629df49d7eSJed Brown } 3639df49d7eSJed Brown 3649df49d7eSJed Brown // ----------------------------------------------------------------------------- 365