1*9df49d7eSJed Brown // Copyright (c) 2017, Lawrence Livermore National Security, LLC. Produced at 2*9df49d7eSJed Brown // the Lawrence Livermore National Laboratory. LLNL-CODE-734707. All Rights 3*9df49d7eSJed Brown // reserved. See files LICENSE and NOTICE for details. 4*9df49d7eSJed Brown // 5*9df49d7eSJed Brown // This file is part of CEED, a collection of benchmarks, miniapps, software 6*9df49d7eSJed Brown // libraries and APIs for efficient high-order finite element and spectral 7*9df49d7eSJed Brown // element discretizations for exascale applications. For more information and 8*9df49d7eSJed Brown // source code availability see http://github.com/ceed. 9*9df49d7eSJed Brown // 10*9df49d7eSJed Brown // The CEED research is supported by the Exascale Computing Project 17-SC-20-SC, 11*9df49d7eSJed Brown // a collaborative effort of two U.S. Department of Energy organizations (Office 12*9df49d7eSJed Brown // of Science and the National Nuclear Security Administration) responsible for 13*9df49d7eSJed Brown // the planning and preparation of a capable exascale ecosystem, including 14*9df49d7eSJed Brown // software, applications, hardware, advanced system engineering and early 15*9df49d7eSJed Brown // testbed platforms, in support of the nation's exascale computing imperative 16*9df49d7eSJed Brown 17*9df49d7eSJed Brown //! A Ceed Basis defines the discrete finite element basis and associated 18*9df49d7eSJed Brown //! quadrature rule. 19*9df49d7eSJed Brown 20*9df49d7eSJed Brown use crate::prelude::*; 21*9df49d7eSJed Brown 22*9df49d7eSJed Brown // ----------------------------------------------------------------------------- 23*9df49d7eSJed Brown // CeedBasis option 24*9df49d7eSJed Brown // ----------------------------------------------------------------------------- 25*9df49d7eSJed Brown #[derive(Clone, Copy)] 26*9df49d7eSJed Brown pub enum BasisOpt<'a> { 27*9df49d7eSJed Brown Some(&'a Basis<'a>), 28*9df49d7eSJed Brown Collocated, 29*9df49d7eSJed Brown } 30*9df49d7eSJed Brown /// Construct a BasisOpt reference from a Basis reference 31*9df49d7eSJed Brown impl<'a> From<&'a Basis<'_>> for BasisOpt<'a> { 32*9df49d7eSJed Brown fn from(basis: &'a Basis) -> Self { 33*9df49d7eSJed Brown debug_assert!(basis.ptr != unsafe { bind_ceed::CEED_BASIS_COLLOCATED }); 34*9df49d7eSJed Brown Self::Some(basis) 35*9df49d7eSJed Brown } 36*9df49d7eSJed Brown } 37*9df49d7eSJed Brown impl<'a> BasisOpt<'a> { 38*9df49d7eSJed Brown /// Transform a Rust libCEED BasisOpt into C libCEED CeedBasis 39*9df49d7eSJed Brown pub(crate) fn to_raw(self) -> bind_ceed::CeedBasis { 40*9df49d7eSJed Brown match self { 41*9df49d7eSJed Brown Self::Some(basis) => basis.ptr, 42*9df49d7eSJed Brown Self::Collocated => unsafe { bind_ceed::CEED_BASIS_COLLOCATED }, 43*9df49d7eSJed Brown } 44*9df49d7eSJed Brown } 45*9df49d7eSJed Brown } 46*9df49d7eSJed Brown 47*9df49d7eSJed Brown // ----------------------------------------------------------------------------- 48*9df49d7eSJed Brown // CeedBasis context wrapper 49*9df49d7eSJed Brown // ----------------------------------------------------------------------------- 50*9df49d7eSJed Brown pub struct Basis<'a> { 51*9df49d7eSJed Brown ceed: &'a crate::Ceed, 52*9df49d7eSJed Brown pub(crate) ptr: bind_ceed::CeedBasis, 53*9df49d7eSJed Brown } 54*9df49d7eSJed Brown 55*9df49d7eSJed Brown // ----------------------------------------------------------------------------- 56*9df49d7eSJed Brown // Destructor 57*9df49d7eSJed Brown // ----------------------------------------------------------------------------- 58*9df49d7eSJed Brown impl<'a> Drop for Basis<'a> { 59*9df49d7eSJed Brown fn drop(&mut self) { 60*9df49d7eSJed Brown unsafe { 61*9df49d7eSJed Brown if self.ptr != bind_ceed::CEED_BASIS_COLLOCATED { 62*9df49d7eSJed Brown bind_ceed::CeedBasisDestroy(&mut self.ptr); 63*9df49d7eSJed Brown } 64*9df49d7eSJed Brown } 65*9df49d7eSJed Brown } 66*9df49d7eSJed Brown } 67*9df49d7eSJed Brown 68*9df49d7eSJed Brown // ----------------------------------------------------------------------------- 69*9df49d7eSJed Brown // Display 70*9df49d7eSJed Brown // ----------------------------------------------------------------------------- 71*9df49d7eSJed Brown impl<'a> fmt::Display for Basis<'a> { 72*9df49d7eSJed Brown /// View a Basis 73*9df49d7eSJed Brown /// 74*9df49d7eSJed Brown /// ``` 75*9df49d7eSJed Brown /// # use libceed::prelude::*; 76*9df49d7eSJed Brown /// # let ceed = libceed::Ceed::default_init(); 77*9df49d7eSJed Brown /// let b = ceed 78*9df49d7eSJed Brown /// .basis_tensor_H1_Lagrange(1, 2, 3, 4, QuadMode::Gauss) 79*9df49d7eSJed Brown /// .unwrap(); 80*9df49d7eSJed Brown /// println!("{}", b); 81*9df49d7eSJed Brown /// ``` 82*9df49d7eSJed Brown fn fmt(&self, f: &mut fmt::Formatter) -> fmt::Result { 83*9df49d7eSJed Brown let mut ptr = std::ptr::null_mut(); 84*9df49d7eSJed Brown let mut sizeloc = crate::MAX_BUFFER_LENGTH; 85*9df49d7eSJed Brown let cstring = unsafe { 86*9df49d7eSJed Brown let file = bind_ceed::open_memstream(&mut ptr, &mut sizeloc); 87*9df49d7eSJed Brown bind_ceed::CeedBasisView(self.ptr, file); 88*9df49d7eSJed Brown bind_ceed::fclose(file); 89*9df49d7eSJed Brown CString::from_raw(ptr) 90*9df49d7eSJed Brown }; 91*9df49d7eSJed Brown cstring.to_string_lossy().fmt(f) 92*9df49d7eSJed Brown } 93*9df49d7eSJed Brown } 94*9df49d7eSJed Brown 95*9df49d7eSJed Brown // ----------------------------------------------------------------------------- 96*9df49d7eSJed Brown // Implementations 97*9df49d7eSJed Brown // ----------------------------------------------------------------------------- 98*9df49d7eSJed Brown impl<'a> Basis<'a> { 99*9df49d7eSJed Brown // Constructors 100*9df49d7eSJed Brown pub fn create_tensor_H1( 101*9df49d7eSJed Brown ceed: &'a crate::Ceed, 102*9df49d7eSJed Brown dim: usize, 103*9df49d7eSJed Brown ncomp: usize, 104*9df49d7eSJed Brown P1d: usize, 105*9df49d7eSJed Brown Q1d: usize, 106*9df49d7eSJed Brown interp1d: &[f64], 107*9df49d7eSJed Brown grad1d: &[f64], 108*9df49d7eSJed Brown qref1d: &[f64], 109*9df49d7eSJed Brown qweight1d: &[f64], 110*9df49d7eSJed Brown ) -> crate::Result<Self> { 111*9df49d7eSJed Brown let mut ptr = std::ptr::null_mut(); 112*9df49d7eSJed Brown let (dim, ncomp, P1d, Q1d) = ( 113*9df49d7eSJed Brown i32::try_from(dim).unwrap(), 114*9df49d7eSJed Brown i32::try_from(ncomp).unwrap(), 115*9df49d7eSJed Brown i32::try_from(P1d).unwrap(), 116*9df49d7eSJed Brown i32::try_from(Q1d).unwrap(), 117*9df49d7eSJed Brown ); 118*9df49d7eSJed Brown let ierr = unsafe { 119*9df49d7eSJed Brown bind_ceed::CeedBasisCreateTensorH1( 120*9df49d7eSJed Brown ceed.ptr, 121*9df49d7eSJed Brown dim, 122*9df49d7eSJed Brown ncomp, 123*9df49d7eSJed Brown P1d, 124*9df49d7eSJed Brown Q1d, 125*9df49d7eSJed Brown interp1d.as_ptr(), 126*9df49d7eSJed Brown grad1d.as_ptr(), 127*9df49d7eSJed Brown qref1d.as_ptr(), 128*9df49d7eSJed Brown qweight1d.as_ptr(), 129*9df49d7eSJed Brown &mut ptr, 130*9df49d7eSJed Brown ) 131*9df49d7eSJed Brown }; 132*9df49d7eSJed Brown ceed.check_error(ierr)?; 133*9df49d7eSJed Brown Ok(Self { ceed, ptr }) 134*9df49d7eSJed Brown } 135*9df49d7eSJed Brown 136*9df49d7eSJed Brown pub fn create_tensor_H1_Lagrange( 137*9df49d7eSJed Brown ceed: &'a crate::Ceed, 138*9df49d7eSJed Brown dim: usize, 139*9df49d7eSJed Brown ncomp: usize, 140*9df49d7eSJed Brown P: usize, 141*9df49d7eSJed Brown Q: usize, 142*9df49d7eSJed Brown qmode: crate::QuadMode, 143*9df49d7eSJed Brown ) -> crate::Result<Self> { 144*9df49d7eSJed Brown let mut ptr = std::ptr::null_mut(); 145*9df49d7eSJed Brown let (dim, ncomp, P, Q, qmode) = ( 146*9df49d7eSJed Brown i32::try_from(dim).unwrap(), 147*9df49d7eSJed Brown i32::try_from(ncomp).unwrap(), 148*9df49d7eSJed Brown i32::try_from(P).unwrap(), 149*9df49d7eSJed Brown i32::try_from(Q).unwrap(), 150*9df49d7eSJed Brown qmode as bind_ceed::CeedQuadMode, 151*9df49d7eSJed Brown ); 152*9df49d7eSJed Brown let ierr = unsafe { 153*9df49d7eSJed Brown bind_ceed::CeedBasisCreateTensorH1Lagrange(ceed.ptr, dim, ncomp, P, Q, qmode, &mut ptr) 154*9df49d7eSJed Brown }; 155*9df49d7eSJed Brown ceed.check_error(ierr)?; 156*9df49d7eSJed Brown Ok(Self { ceed, ptr }) 157*9df49d7eSJed Brown } 158*9df49d7eSJed Brown 159*9df49d7eSJed Brown pub fn create_H1( 160*9df49d7eSJed Brown ceed: &'a crate::Ceed, 161*9df49d7eSJed Brown topo: crate::ElemTopology, 162*9df49d7eSJed Brown ncomp: usize, 163*9df49d7eSJed Brown nnodes: usize, 164*9df49d7eSJed Brown nqpts: usize, 165*9df49d7eSJed Brown interp: &[f64], 166*9df49d7eSJed Brown grad: &[f64], 167*9df49d7eSJed Brown qref: &[f64], 168*9df49d7eSJed Brown qweight: &[f64], 169*9df49d7eSJed Brown ) -> crate::Result<Self> { 170*9df49d7eSJed Brown let mut ptr = std::ptr::null_mut(); 171*9df49d7eSJed Brown let (topo, ncomp, nnodes, nqpts) = ( 172*9df49d7eSJed Brown topo as bind_ceed::CeedElemTopology, 173*9df49d7eSJed Brown i32::try_from(ncomp).unwrap(), 174*9df49d7eSJed Brown i32::try_from(nnodes).unwrap(), 175*9df49d7eSJed Brown i32::try_from(nqpts).unwrap(), 176*9df49d7eSJed Brown ); 177*9df49d7eSJed Brown let ierr = unsafe { 178*9df49d7eSJed Brown bind_ceed::CeedBasisCreateH1( 179*9df49d7eSJed Brown ceed.ptr, 180*9df49d7eSJed Brown topo, 181*9df49d7eSJed Brown ncomp, 182*9df49d7eSJed Brown nnodes, 183*9df49d7eSJed Brown nqpts, 184*9df49d7eSJed Brown interp.as_ptr(), 185*9df49d7eSJed Brown grad.as_ptr(), 186*9df49d7eSJed Brown qref.as_ptr(), 187*9df49d7eSJed Brown qweight.as_ptr(), 188*9df49d7eSJed Brown &mut ptr, 189*9df49d7eSJed Brown ) 190*9df49d7eSJed Brown }; 191*9df49d7eSJed Brown ceed.check_error(ierr)?; 192*9df49d7eSJed Brown Ok(Self { ceed, ptr }) 193*9df49d7eSJed Brown } 194*9df49d7eSJed Brown 195*9df49d7eSJed Brown /// Apply basis evaluation from nodes to quadrature points or vice versa 196*9df49d7eSJed Brown /// 197*9df49d7eSJed Brown /// * `nelem` - The number of elements to apply the basis evaluation to 198*9df49d7eSJed Brown /// * `tmode` - `TrasposeMode::NoTranspose` to evaluate from nodes to 199*9df49d7eSJed Brown /// quadrature points, `TransposeMode::Transpose` to apply the 200*9df49d7eSJed Brown /// transpose, mapping from quadrature points to nodes 201*9df49d7eSJed Brown /// * `emode` - `EvalMode::None` to use values directly, `EvalMode::Interp` 202*9df49d7eSJed Brown /// to use interpolated values, `EvalMode::Grad` to use 203*9df49d7eSJed Brown /// gradients, `EvalMode::Weight` to use quadrature weights 204*9df49d7eSJed Brown /// * `u` - Input Vector 205*9df49d7eSJed Brown /// * `v` - Output Vector 206*9df49d7eSJed Brown /// 207*9df49d7eSJed Brown /// ``` 208*9df49d7eSJed Brown /// # use libceed::prelude::*; 209*9df49d7eSJed Brown /// # let ceed = libceed::Ceed::default_init(); 210*9df49d7eSJed Brown /// const Q: usize = 6; 211*9df49d7eSJed Brown /// let bu = ceed 212*9df49d7eSJed Brown /// .basis_tensor_H1_Lagrange(1, 1, Q, Q, QuadMode::GaussLobatto) 213*9df49d7eSJed Brown /// .unwrap(); 214*9df49d7eSJed Brown /// let bx = ceed 215*9df49d7eSJed Brown /// .basis_tensor_H1_Lagrange(1, 1, 2, Q, QuadMode::Gauss) 216*9df49d7eSJed Brown /// .unwrap(); 217*9df49d7eSJed Brown /// 218*9df49d7eSJed Brown /// let x_corners = ceed.vector_from_slice(&[-1., 1.]).unwrap(); 219*9df49d7eSJed Brown /// let mut x_qpts = ceed.vector(Q).unwrap(); 220*9df49d7eSJed Brown /// let mut x_nodes = ceed.vector(Q).unwrap(); 221*9df49d7eSJed Brown /// bx.apply( 222*9df49d7eSJed Brown /// 1, 223*9df49d7eSJed Brown /// TransposeMode::NoTranspose, 224*9df49d7eSJed Brown /// EvalMode::Interp, 225*9df49d7eSJed Brown /// &x_corners, 226*9df49d7eSJed Brown /// &mut x_nodes, 227*9df49d7eSJed Brown /// ); 228*9df49d7eSJed Brown /// bu.apply( 229*9df49d7eSJed Brown /// 1, 230*9df49d7eSJed Brown /// TransposeMode::NoTranspose, 231*9df49d7eSJed Brown /// EvalMode::Interp, 232*9df49d7eSJed Brown /// &x_nodes, 233*9df49d7eSJed Brown /// &mut x_qpts, 234*9df49d7eSJed Brown /// ); 235*9df49d7eSJed Brown /// 236*9df49d7eSJed Brown /// // Create function x^3 + 1 on Gauss Lobatto points 237*9df49d7eSJed Brown /// let mut u_arr = [0.; Q]; 238*9df49d7eSJed Brown /// u_arr 239*9df49d7eSJed Brown /// .iter_mut() 240*9df49d7eSJed Brown /// .zip(x_nodes.view().iter()) 241*9df49d7eSJed Brown /// .for_each(|(u, x)| *u = x * x * x + 1.); 242*9df49d7eSJed Brown /// let u = ceed.vector_from_slice(&u_arr).unwrap(); 243*9df49d7eSJed Brown /// 244*9df49d7eSJed Brown /// // Map function to Gauss points 245*9df49d7eSJed Brown /// let mut v = ceed.vector(Q).unwrap(); 246*9df49d7eSJed Brown /// v.set_value(0.); 247*9df49d7eSJed Brown /// bu.apply(1, TransposeMode::NoTranspose, EvalMode::Interp, &u, &mut v) 248*9df49d7eSJed Brown /// .unwrap(); 249*9df49d7eSJed Brown /// 250*9df49d7eSJed Brown /// // Verify results 251*9df49d7eSJed Brown /// v.view() 252*9df49d7eSJed Brown /// .iter() 253*9df49d7eSJed Brown /// .zip(x_qpts.view().iter()) 254*9df49d7eSJed Brown /// .for_each(|(v, x)| { 255*9df49d7eSJed Brown /// let true_value = x * x * x + 1.; 256*9df49d7eSJed Brown /// assert_eq!(*v, true_value, "Incorrect basis application"); 257*9df49d7eSJed Brown /// }); 258*9df49d7eSJed Brown /// ``` 259*9df49d7eSJed Brown pub fn apply( 260*9df49d7eSJed Brown &self, 261*9df49d7eSJed Brown nelem: usize, 262*9df49d7eSJed Brown tmode: TransposeMode, 263*9df49d7eSJed Brown emode: EvalMode, 264*9df49d7eSJed Brown u: &Vector, 265*9df49d7eSJed Brown v: &mut Vector, 266*9df49d7eSJed Brown ) -> crate::Result<i32> { 267*9df49d7eSJed Brown let (nelem, tmode, emode) = ( 268*9df49d7eSJed Brown i32::try_from(nelem).unwrap(), 269*9df49d7eSJed Brown tmode as bind_ceed::CeedTransposeMode, 270*9df49d7eSJed Brown emode as bind_ceed::CeedEvalMode, 271*9df49d7eSJed Brown ); 272*9df49d7eSJed Brown let ierr = 273*9df49d7eSJed Brown unsafe { bind_ceed::CeedBasisApply(self.ptr, nelem, tmode, emode, u.ptr, v.ptr) }; 274*9df49d7eSJed Brown self.ceed.check_error(ierr) 275*9df49d7eSJed Brown } 276*9df49d7eSJed Brown 277*9df49d7eSJed Brown /// Returns the dimension for given CeedBasis 278*9df49d7eSJed Brown /// 279*9df49d7eSJed Brown /// ``` 280*9df49d7eSJed Brown /// # use libceed::prelude::*; 281*9df49d7eSJed Brown /// # let ceed = libceed::Ceed::default_init(); 282*9df49d7eSJed Brown /// let dim = 2; 283*9df49d7eSJed Brown /// let b = ceed 284*9df49d7eSJed Brown /// .basis_tensor_H1_Lagrange(dim, 1, 3, 4, QuadMode::Gauss) 285*9df49d7eSJed Brown /// .unwrap(); 286*9df49d7eSJed Brown /// 287*9df49d7eSJed Brown /// let d = b.dimension(); 288*9df49d7eSJed Brown /// assert_eq!(d, dim, "Incorrect dimension"); 289*9df49d7eSJed Brown /// ``` 290*9df49d7eSJed Brown pub fn dimension(&self) -> usize { 291*9df49d7eSJed Brown let mut dim = 0; 292*9df49d7eSJed Brown unsafe { bind_ceed::CeedBasisGetDimension(self.ptr, &mut dim) }; 293*9df49d7eSJed Brown usize::try_from(dim).unwrap() 294*9df49d7eSJed Brown } 295*9df49d7eSJed Brown 296*9df49d7eSJed Brown /// Returns number of components for given CeedBasis 297*9df49d7eSJed Brown /// 298*9df49d7eSJed Brown /// ``` 299*9df49d7eSJed Brown /// # use libceed::prelude::*; 300*9df49d7eSJed Brown /// # let ceed = libceed::Ceed::default_init(); 301*9df49d7eSJed Brown /// let ncomp = 2; 302*9df49d7eSJed Brown /// let b = ceed 303*9df49d7eSJed Brown /// .basis_tensor_H1_Lagrange(1, ncomp, 3, 4, QuadMode::Gauss) 304*9df49d7eSJed Brown /// .unwrap(); 305*9df49d7eSJed Brown /// 306*9df49d7eSJed Brown /// let n = b.num_components(); 307*9df49d7eSJed Brown /// assert_eq!(n, ncomp, "Incorrect number of components"); 308*9df49d7eSJed Brown /// ``` 309*9df49d7eSJed Brown pub fn num_components(&self) -> usize { 310*9df49d7eSJed Brown let mut ncomp = 0; 311*9df49d7eSJed Brown unsafe { bind_ceed::CeedBasisGetNumComponents(self.ptr, &mut ncomp) }; 312*9df49d7eSJed Brown usize::try_from(ncomp).unwrap() 313*9df49d7eSJed Brown } 314*9df49d7eSJed Brown 315*9df49d7eSJed Brown /// Returns total number of nodes (in dim dimensions) of a CeedBasis 316*9df49d7eSJed Brown /// 317*9df49d7eSJed Brown /// ``` 318*9df49d7eSJed Brown /// # use libceed::prelude::*; 319*9df49d7eSJed Brown /// # let ceed = libceed::Ceed::default_init(); 320*9df49d7eSJed Brown /// let p = 3; 321*9df49d7eSJed Brown /// let b = ceed 322*9df49d7eSJed Brown /// .basis_tensor_H1_Lagrange(2, 1, p, 4, QuadMode::Gauss) 323*9df49d7eSJed Brown /// .unwrap(); 324*9df49d7eSJed Brown /// 325*9df49d7eSJed Brown /// let nnodes = b.num_nodes(); 326*9df49d7eSJed Brown /// assert_eq!(nnodes, p * p, "Incorrect number of nodes"); 327*9df49d7eSJed Brown /// ``` 328*9df49d7eSJed Brown pub fn num_nodes(&self) -> usize { 329*9df49d7eSJed Brown let mut nnodes = 0; 330*9df49d7eSJed Brown unsafe { bind_ceed::CeedBasisGetNumNodes(self.ptr, &mut nnodes) }; 331*9df49d7eSJed Brown usize::try_from(nnodes).unwrap() 332*9df49d7eSJed Brown } 333*9df49d7eSJed Brown 334*9df49d7eSJed Brown /// Returns total number of quadrature points (in dim dimensions) of a 335*9df49d7eSJed Brown /// CeedBasis 336*9df49d7eSJed Brown /// 337*9df49d7eSJed Brown /// ``` 338*9df49d7eSJed Brown /// # use libceed::prelude::*; 339*9df49d7eSJed Brown /// # let ceed = libceed::Ceed::default_init(); 340*9df49d7eSJed Brown /// let q = 4; 341*9df49d7eSJed Brown /// let b = ceed 342*9df49d7eSJed Brown /// .basis_tensor_H1_Lagrange(2, 1, 3, q, QuadMode::Gauss) 343*9df49d7eSJed Brown /// .unwrap(); 344*9df49d7eSJed Brown /// 345*9df49d7eSJed Brown /// let nqpts = b.num_quadrature_points(); 346*9df49d7eSJed Brown /// assert_eq!(nqpts, q * q, "Incorrect number of quadrature points"); 347*9df49d7eSJed Brown /// ``` 348*9df49d7eSJed Brown pub fn num_quadrature_points(&self) -> usize { 349*9df49d7eSJed Brown let mut Q = 0; 350*9df49d7eSJed Brown unsafe { 351*9df49d7eSJed Brown bind_ceed::CeedBasisGetNumQuadraturePoints(self.ptr, &mut Q); 352*9df49d7eSJed Brown } 353*9df49d7eSJed Brown usize::try_from(Q).unwrap() 354*9df49d7eSJed Brown } 355*9df49d7eSJed Brown } 356*9df49d7eSJed Brown 357*9df49d7eSJed Brown // ----------------------------------------------------------------------------- 358