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 //! # libCEED Rust Interface 18*9df49d7eSJed Brown //! 19*9df49d7eSJed Brown //! This is the documentation for the high level libCEED Rust interface. 20*9df49d7eSJed Brown //! See the full libCEED user manual [here](https://libceed.readthedocs.io). 21*9df49d7eSJed Brown //! 22*9df49d7eSJed Brown //! libCEED is a low-level API for for the efficient high-order discretization methods 23*9df49d7eSJed Brown //! developed by the ECP co-design Center for Efficient Exascale Discretizations (CEED). 24*9df49d7eSJed Brown //! While our focus is on high-order finite elements, the approach is mostly algebraic 25*9df49d7eSJed Brown //! and thus applicable to other discretizations in factored form. 26*9df49d7eSJed Brown //! 27*9df49d7eSJed Brown //! ## Usage 28*9df49d7eSJed Brown //! 29*9df49d7eSJed Brown //! To call libCEED from a Rust package, the following `Cargo.toml` can be used. 30*9df49d7eSJed Brown //! ```toml 31*9df49d7eSJed Brown //! [dependencies] 32*9df49d7eSJed Brown //! libceed = { git = "https://github.com/CEED/libCEED", branch = "main" } 33*9df49d7eSJed Brown //! ``` 34*9df49d7eSJed Brown //! 35*9df49d7eSJed Brown //! Supported features: 36*9df49d7eSJed Brown //! * `static` (default): link to static libceed.a 37*9df49d7eSJed Brown //! * `system`: use libceed from a system directory (otherwise, install from source) 38*9df49d7eSJed Brown //! 39*9df49d7eSJed Brown //! ``` 40*9df49d7eSJed Brown //! extern crate libceed; 41*9df49d7eSJed Brown //! 42*9df49d7eSJed Brown //! fn main() { 43*9df49d7eSJed Brown //! let ceed = libceed::Ceed::init("/cpu/self/ref"); 44*9df49d7eSJed Brown //! let xc = ceed.vector_from_slice(&[0., 0.5, 1.0]).unwrap(); 45*9df49d7eSJed Brown //! let xs = xc.view(); 46*9df49d7eSJed Brown //! assert_eq!(xs[..], [0., 0.5, 1.0]); 47*9df49d7eSJed Brown //! } 48*9df49d7eSJed Brown //! ``` 49*9df49d7eSJed Brown //! 50*9df49d7eSJed Brown //! ## Examples 51*9df49d7eSJed Brown //! 52*9df49d7eSJed Brown //! Examples of libCEED can be found in the libCEED GitHub repository under the 53*9df49d7eSJed Brown //! `examples/rust` folder. 54*9df49d7eSJed Brown //! 55*9df49d7eSJed Brown //! ## Development 56*9df49d7eSJed Brown //! 57*9df49d7eSJed Brown //! To develop libCEED, use `cargo build` in the `rust/` directory to install a 58*9df49d7eSJed Brown //! local copy and build the bindings. If you need custom flags for the C 59*9df49d7eSJed Brown //! project, we recommend using `make configure` to cache arguments. If you 60*9df49d7eSJed Brown //! disable the `static` feature, then you'll need to set `LD_LIBRARY_PATH` for 61*9df49d7eSJed Brown //! doctests to be able to find it. You can do this in `$CEED_DIR/lib` and set 62*9df49d7eSJed Brown //! `PKG_CONFIG_PATH`. 63*9df49d7eSJed Brown //! 64*9df49d7eSJed Brown //! Note: the `LD_LIBRARY_PATH` workarounds will become unnecessary if [this 65*9df49d7eSJed Brown //! issue](https://github.com/rust-lang/cargo/issues/1592) is resolved -- it's 66*9df49d7eSJed Brown //! currently closed, but the problem still exists. 67*9df49d7eSJed Brown 68*9df49d7eSJed Brown // ----------------------------------------------------------------------------- 69*9df49d7eSJed Brown // Exceptions 70*9df49d7eSJed Brown // ----------------------------------------------------------------------------- 71*9df49d7eSJed Brown #![allow(non_snake_case)] 72*9df49d7eSJed Brown 73*9df49d7eSJed Brown // ----------------------------------------------------------------------------- 74*9df49d7eSJed Brown // Crate prelude 75*9df49d7eSJed Brown // ----------------------------------------------------------------------------- 76*9df49d7eSJed Brown use crate::prelude::*; 77*9df49d7eSJed Brown use std::sync::Once; 78*9df49d7eSJed Brown 79*9df49d7eSJed Brown pub mod prelude { 80*9df49d7eSJed Brown pub(crate) use libceed_sys::bind_ceed; 81*9df49d7eSJed Brown pub use crate::{ 82*9df49d7eSJed Brown basis::{self, Basis, BasisOpt}, 83*9df49d7eSJed Brown elem_restriction::{self, ElemRestriction, ElemRestrictionOpt}, 84*9df49d7eSJed Brown operator::{self, CompositeOperator, Operator}, 85*9df49d7eSJed Brown qfunction::{ 86*9df49d7eSJed Brown self, QFunction, QFunctionByName, QFunctionInputs, QFunctionOpt, QFunctionOutputs, 87*9df49d7eSJed Brown }, 88*9df49d7eSJed Brown vector::{self, Vector, VectorOpt}, 89*9df49d7eSJed Brown ElemTopology, EvalMode, MemType, NormType, QuadMode, TransposeMode, CEED_STRIDES_BACKEND, 90*9df49d7eSJed Brown MAX_QFUNCTION_FIELDS, 91*9df49d7eSJed Brown }; 92*9df49d7eSJed Brown pub(crate) use std::convert::TryFrom; 93*9df49d7eSJed Brown pub(crate) use std::ffi::CString; 94*9df49d7eSJed Brown pub(crate) use std::fmt; 95*9df49d7eSJed Brown } 96*9df49d7eSJed Brown 97*9df49d7eSJed Brown // ----------------------------------------------------------------------------- 98*9df49d7eSJed Brown // Modules 99*9df49d7eSJed Brown // ----------------------------------------------------------------------------- 100*9df49d7eSJed Brown pub mod basis; 101*9df49d7eSJed Brown pub mod elem_restriction; 102*9df49d7eSJed Brown pub mod operator; 103*9df49d7eSJed Brown pub mod qfunction; 104*9df49d7eSJed Brown pub mod vector; 105*9df49d7eSJed Brown 106*9df49d7eSJed Brown // ----------------------------------------------------------------------------- 107*9df49d7eSJed Brown // Constants for library 108*9df49d7eSJed Brown // ----------------------------------------------------------------------------- 109*9df49d7eSJed Brown const MAX_BUFFER_LENGTH: u64 = 4096; 110*9df49d7eSJed Brown pub const MAX_QFUNCTION_FIELDS: usize = 16; 111*9df49d7eSJed Brown pub const CEED_STRIDES_BACKEND: [i32; 3] = [0; 3]; 112*9df49d7eSJed Brown 113*9df49d7eSJed Brown // ----------------------------------------------------------------------------- 114*9df49d7eSJed Brown // Enums for libCEED 115*9df49d7eSJed Brown // ----------------------------------------------------------------------------- 116*9df49d7eSJed Brown #[derive(Clone, Copy, PartialEq, Eq)] 117*9df49d7eSJed Brown /// Many Ceed interfaces take or return pointers to memory. This enum is used to 118*9df49d7eSJed Brown /// specify where the memory being provided or requested must reside. 119*9df49d7eSJed Brown pub enum MemType { 120*9df49d7eSJed Brown Host = bind_ceed::CeedMemType_CEED_MEM_HOST as isize, 121*9df49d7eSJed Brown Device = bind_ceed::CeedMemType_CEED_MEM_DEVICE as isize, 122*9df49d7eSJed Brown } 123*9df49d7eSJed Brown 124*9df49d7eSJed Brown #[derive(Clone, Copy, PartialEq, Eq)] 125*9df49d7eSJed Brown // OwnPointer will not be used by user but is included for internal use 126*9df49d7eSJed Brown #[allow(dead_code)] 127*9df49d7eSJed Brown /// Conveys ownership status of arrays passed to Ceed interfaces. 128*9df49d7eSJed Brown pub(crate) enum CopyMode { 129*9df49d7eSJed Brown CopyValues = bind_ceed::CeedCopyMode_CEED_COPY_VALUES as isize, 130*9df49d7eSJed Brown UsePointer = bind_ceed::CeedCopyMode_CEED_USE_POINTER as isize, 131*9df49d7eSJed Brown OwnPointer = bind_ceed::CeedCopyMode_CEED_OWN_POINTER as isize, 132*9df49d7eSJed Brown } 133*9df49d7eSJed Brown 134*9df49d7eSJed Brown #[derive(Clone, Copy, PartialEq, Eq)] 135*9df49d7eSJed Brown /// Denotes type of vector norm to be computed 136*9df49d7eSJed Brown pub enum NormType { 137*9df49d7eSJed Brown One = bind_ceed::CeedNormType_CEED_NORM_1 as isize, 138*9df49d7eSJed Brown Two = bind_ceed::CeedNormType_CEED_NORM_2 as isize, 139*9df49d7eSJed Brown Max = bind_ceed::CeedNormType_CEED_NORM_MAX as isize, 140*9df49d7eSJed Brown } 141*9df49d7eSJed Brown 142*9df49d7eSJed Brown #[derive(Clone, Copy, PartialEq, Eq)] 143*9df49d7eSJed Brown /// Denotes whether a linear transformation or its transpose should be applied 144*9df49d7eSJed Brown pub enum TransposeMode { 145*9df49d7eSJed Brown NoTranspose = bind_ceed::CeedTransposeMode_CEED_NOTRANSPOSE as isize, 146*9df49d7eSJed Brown Transpose = bind_ceed::CeedTransposeMode_CEED_TRANSPOSE as isize, 147*9df49d7eSJed Brown } 148*9df49d7eSJed Brown 149*9df49d7eSJed Brown #[derive(Clone, Copy, PartialEq, Eq)] 150*9df49d7eSJed Brown /// Type of quadrature; also used for location of nodes 151*9df49d7eSJed Brown pub enum QuadMode { 152*9df49d7eSJed Brown Gauss = bind_ceed::CeedQuadMode_CEED_GAUSS as isize, 153*9df49d7eSJed Brown GaussLobatto = bind_ceed::CeedQuadMode_CEED_GAUSS_LOBATTO as isize, 154*9df49d7eSJed Brown } 155*9df49d7eSJed Brown 156*9df49d7eSJed Brown #[derive(Clone, Copy, PartialEq, Eq)] 157*9df49d7eSJed Brown /// Type of basis shape to create non-tensor H1 element basis 158*9df49d7eSJed Brown pub enum ElemTopology { 159*9df49d7eSJed Brown Line = bind_ceed::CeedElemTopology_CEED_LINE as isize, 160*9df49d7eSJed Brown Triangle = bind_ceed::CeedElemTopology_CEED_TRIANGLE as isize, 161*9df49d7eSJed Brown Quad = bind_ceed::CeedElemTopology_CEED_QUAD as isize, 162*9df49d7eSJed Brown Tet = bind_ceed::CeedElemTopology_CEED_TET as isize, 163*9df49d7eSJed Brown Pyramid = bind_ceed::CeedElemTopology_CEED_PYRAMID as isize, 164*9df49d7eSJed Brown Prism = bind_ceed::CeedElemTopology_CEED_PRISM as isize, 165*9df49d7eSJed Brown Hex = bind_ceed::CeedElemTopology_CEED_HEX as isize, 166*9df49d7eSJed Brown } 167*9df49d7eSJed Brown 168*9df49d7eSJed Brown #[derive(Clone, Copy, PartialEq, Eq)] 169*9df49d7eSJed Brown /// Basis evaluation mode 170*9df49d7eSJed Brown pub enum EvalMode { 171*9df49d7eSJed Brown None = bind_ceed::CeedEvalMode_CEED_EVAL_NONE as isize, 172*9df49d7eSJed Brown Interp = bind_ceed::CeedEvalMode_CEED_EVAL_INTERP as isize, 173*9df49d7eSJed Brown Grad = bind_ceed::CeedEvalMode_CEED_EVAL_GRAD as isize, 174*9df49d7eSJed Brown Div = bind_ceed::CeedEvalMode_CEED_EVAL_DIV as isize, 175*9df49d7eSJed Brown Curl = bind_ceed::CeedEvalMode_CEED_EVAL_CURL as isize, 176*9df49d7eSJed Brown Weight = bind_ceed::CeedEvalMode_CEED_EVAL_WEIGHT as isize, 177*9df49d7eSJed Brown } 178*9df49d7eSJed Brown 179*9df49d7eSJed Brown // ----------------------------------------------------------------------------- 180*9df49d7eSJed Brown // Ceed error 181*9df49d7eSJed Brown // ----------------------------------------------------------------------------- 182*9df49d7eSJed Brown type Result<T> = std::result::Result<T, CeedError>; 183*9df49d7eSJed Brown 184*9df49d7eSJed Brown #[derive(Debug)] 185*9df49d7eSJed Brown pub struct CeedError { 186*9df49d7eSJed Brown pub message: String, 187*9df49d7eSJed Brown } 188*9df49d7eSJed Brown 189*9df49d7eSJed Brown impl fmt::Display for CeedError { 190*9df49d7eSJed Brown fn fmt(&self, f: &mut fmt::Formatter) -> fmt::Result { 191*9df49d7eSJed Brown write!(f, "{}", self.message) 192*9df49d7eSJed Brown } 193*9df49d7eSJed Brown } 194*9df49d7eSJed Brown 195*9df49d7eSJed Brown // ----------------------------------------------------------------------------- 196*9df49d7eSJed Brown // Ceed error handler 197*9df49d7eSJed Brown // ----------------------------------------------------------------------------- 198*9df49d7eSJed Brown pub enum CeedErrorHandler { 199*9df49d7eSJed Brown ErrorAbort, 200*9df49d7eSJed Brown ErrorExit, 201*9df49d7eSJed Brown ErrorReturn, 202*9df49d7eSJed Brown ErrorStore, 203*9df49d7eSJed Brown } 204*9df49d7eSJed Brown 205*9df49d7eSJed Brown // ----------------------------------------------------------------------------- 206*9df49d7eSJed Brown // Ceed context wrapper 207*9df49d7eSJed Brown // ----------------------------------------------------------------------------- 208*9df49d7eSJed Brown /// A Ceed is a library context representing control of a logical hardware 209*9df49d7eSJed Brown /// resource. 210*9df49d7eSJed Brown #[derive(Debug)] 211*9df49d7eSJed Brown pub struct Ceed { 212*9df49d7eSJed Brown ptr: bind_ceed::Ceed, 213*9df49d7eSJed Brown } 214*9df49d7eSJed Brown 215*9df49d7eSJed Brown // ----------------------------------------------------------------------------- 216*9df49d7eSJed Brown // Destructor 217*9df49d7eSJed Brown // ----------------------------------------------------------------------------- 218*9df49d7eSJed Brown impl Drop for Ceed { 219*9df49d7eSJed Brown fn drop(&mut self) { 220*9df49d7eSJed Brown unsafe { 221*9df49d7eSJed Brown bind_ceed::CeedDestroy(&mut self.ptr); 222*9df49d7eSJed Brown } 223*9df49d7eSJed Brown } 224*9df49d7eSJed Brown } 225*9df49d7eSJed Brown 226*9df49d7eSJed Brown // ----------------------------------------------------------------------------- 227*9df49d7eSJed Brown // Display 228*9df49d7eSJed Brown // ----------------------------------------------------------------------------- 229*9df49d7eSJed Brown impl fmt::Display for Ceed { 230*9df49d7eSJed Brown /// View a Ceed 231*9df49d7eSJed Brown /// 232*9df49d7eSJed Brown /// ``` 233*9df49d7eSJed Brown /// let ceed = libceed::Ceed::init("/cpu/self/ref/serial"); 234*9df49d7eSJed Brown /// println!("{}", ceed); 235*9df49d7eSJed Brown /// ``` 236*9df49d7eSJed Brown fn fmt(&self, f: &mut fmt::Formatter) -> fmt::Result { 237*9df49d7eSJed Brown let mut ptr = std::ptr::null_mut(); 238*9df49d7eSJed Brown let mut sizeloc = crate::MAX_BUFFER_LENGTH; 239*9df49d7eSJed Brown let cstring = unsafe { 240*9df49d7eSJed Brown let file = bind_ceed::open_memstream(&mut ptr, &mut sizeloc); 241*9df49d7eSJed Brown bind_ceed::CeedView(self.ptr, file); 242*9df49d7eSJed Brown bind_ceed::fclose(file); 243*9df49d7eSJed Brown CString::from_raw(ptr) 244*9df49d7eSJed Brown }; 245*9df49d7eSJed Brown cstring.to_string_lossy().fmt(f) 246*9df49d7eSJed Brown } 247*9df49d7eSJed Brown } 248*9df49d7eSJed Brown 249*9df49d7eSJed Brown static REGISTER: Once = Once::new(); 250*9df49d7eSJed Brown 251*9df49d7eSJed Brown // ----------------------------------------------------------------------------- 252*9df49d7eSJed Brown // Object constructors 253*9df49d7eSJed Brown // ----------------------------------------------------------------------------- 254*9df49d7eSJed Brown impl Ceed { 255*9df49d7eSJed Brown /// Returns a Ceed context initialized with the specified resource 256*9df49d7eSJed Brown /// 257*9df49d7eSJed Brown /// # arguments 258*9df49d7eSJed Brown /// 259*9df49d7eSJed Brown /// * `resource` - Resource to use, e.g., "/cpu/self" 260*9df49d7eSJed Brown /// 261*9df49d7eSJed Brown /// ``` 262*9df49d7eSJed Brown /// let ceed = libceed::Ceed::init("/cpu/self/ref/serial"); 263*9df49d7eSJed Brown /// ``` 264*9df49d7eSJed Brown pub fn init(resource: &str) -> Self { 265*9df49d7eSJed Brown Ceed::init_with_error_handler(resource, CeedErrorHandler::ErrorStore) 266*9df49d7eSJed Brown } 267*9df49d7eSJed Brown 268*9df49d7eSJed Brown /// Returns a Ceed context initialized with the specified resource 269*9df49d7eSJed Brown /// 270*9df49d7eSJed Brown /// # arguments 271*9df49d7eSJed Brown /// 272*9df49d7eSJed Brown /// * `resource` - Resource to use, e.g., "/cpu/self" 273*9df49d7eSJed Brown /// 274*9df49d7eSJed Brown /// ``` 275*9df49d7eSJed Brown /// let ceed = libceed::Ceed::init_with_error_handler( 276*9df49d7eSJed Brown /// "/cpu/self/ref/serial", 277*9df49d7eSJed Brown /// libceed::CeedErrorHandler::ErrorAbort, 278*9df49d7eSJed Brown /// ); 279*9df49d7eSJed Brown /// ``` 280*9df49d7eSJed Brown pub fn init_with_error_handler(resource: &str, handler: CeedErrorHandler) -> Self { 281*9df49d7eSJed Brown REGISTER.call_once(|| unsafe { 282*9df49d7eSJed Brown bind_ceed::CeedRegisterAll(); 283*9df49d7eSJed Brown bind_ceed::CeedQFunctionRegisterAll(); 284*9df49d7eSJed Brown }); 285*9df49d7eSJed Brown 286*9df49d7eSJed Brown // Convert to C string 287*9df49d7eSJed Brown let c_resource = CString::new(resource).expect("CString::new failed"); 288*9df49d7eSJed Brown 289*9df49d7eSJed Brown // Get error handler pointer 290*9df49d7eSJed Brown let eh = match handler { 291*9df49d7eSJed Brown CeedErrorHandler::ErrorAbort => bind_ceed::CeedErrorAbort, 292*9df49d7eSJed Brown CeedErrorHandler::ErrorExit => bind_ceed::CeedErrorExit, 293*9df49d7eSJed Brown CeedErrorHandler::ErrorReturn => bind_ceed::CeedErrorReturn, 294*9df49d7eSJed Brown CeedErrorHandler::ErrorStore => bind_ceed::CeedErrorStore, 295*9df49d7eSJed Brown }; 296*9df49d7eSJed Brown 297*9df49d7eSJed Brown // Call to libCEED 298*9df49d7eSJed Brown let mut ptr = std::ptr::null_mut(); 299*9df49d7eSJed Brown let mut ierr = unsafe { bind_ceed::CeedInit(c_resource.as_ptr() as *const i8, &mut ptr) }; 300*9df49d7eSJed Brown if ierr != 0 { 301*9df49d7eSJed Brown panic!("Error initializing backend resource: {}", resource) 302*9df49d7eSJed Brown } 303*9df49d7eSJed Brown ierr = unsafe { bind_ceed::CeedSetErrorHandler(ptr, Some(eh)) }; 304*9df49d7eSJed Brown let ceed = Ceed { ptr }; 305*9df49d7eSJed Brown ceed.check_error(ierr).unwrap(); 306*9df49d7eSJed Brown ceed 307*9df49d7eSJed Brown } 308*9df49d7eSJed Brown 309*9df49d7eSJed Brown /// Default initializer for testing 310*9df49d7eSJed Brown #[doc(hidden)] 311*9df49d7eSJed Brown pub fn default_init() -> Self { 312*9df49d7eSJed Brown // Convert to C string 313*9df49d7eSJed Brown let resource = "/cpu/self/ref/serial"; 314*9df49d7eSJed Brown crate::Ceed::init(resource) 315*9df49d7eSJed Brown } 316*9df49d7eSJed Brown 317*9df49d7eSJed Brown /// Internal error checker 318*9df49d7eSJed Brown #[doc(hidden)] 319*9df49d7eSJed Brown fn check_error(&self, ierr: i32) -> Result<i32> { 320*9df49d7eSJed Brown // Return early if code is clean 321*9df49d7eSJed Brown if ierr == bind_ceed::CeedErrorType_CEED_ERROR_SUCCESS { 322*9df49d7eSJed Brown return Ok(ierr); 323*9df49d7eSJed Brown } 324*9df49d7eSJed Brown // Retrieve error message 325*9df49d7eSJed Brown let mut ptr: *const std::os::raw::c_char = std::ptr::null_mut(); 326*9df49d7eSJed Brown let c_str = unsafe { 327*9df49d7eSJed Brown bind_ceed::CeedGetErrorMessage(self.ptr, &mut ptr); 328*9df49d7eSJed Brown std::ffi::CStr::from_ptr(ptr) 329*9df49d7eSJed Brown }; 330*9df49d7eSJed Brown let message = c_str.to_string_lossy().to_string(); 331*9df49d7eSJed Brown // Panic if negative code, otherwise return error 332*9df49d7eSJed Brown if ierr < bind_ceed::CeedErrorType_CEED_ERROR_SUCCESS { 333*9df49d7eSJed Brown panic!("{}", message); 334*9df49d7eSJed Brown } 335*9df49d7eSJed Brown Err(CeedError { message }) 336*9df49d7eSJed Brown } 337*9df49d7eSJed Brown 338*9df49d7eSJed Brown /// Returns full resource name for a Ceed context 339*9df49d7eSJed Brown /// 340*9df49d7eSJed Brown /// ``` 341*9df49d7eSJed Brown /// let ceed = libceed::Ceed::init("/cpu/self/ref/serial"); 342*9df49d7eSJed Brown /// let resource = ceed.resource(); 343*9df49d7eSJed Brown /// 344*9df49d7eSJed Brown /// assert_eq!(resource, "/cpu/self/ref/serial".to_string()) 345*9df49d7eSJed Brown /// ``` 346*9df49d7eSJed Brown pub fn resource(&self) -> String { 347*9df49d7eSJed Brown let mut ptr: *const std::os::raw::c_char = std::ptr::null_mut(); 348*9df49d7eSJed Brown let c_str = unsafe { 349*9df49d7eSJed Brown bind_ceed::CeedGetResource(self.ptr, &mut ptr); 350*9df49d7eSJed Brown std::ffi::CStr::from_ptr(ptr) 351*9df49d7eSJed Brown }; 352*9df49d7eSJed Brown c_str.to_string_lossy().to_string() 353*9df49d7eSJed Brown } 354*9df49d7eSJed Brown 355*9df49d7eSJed Brown /// Returns a CeedVector of the specified length (does not allocate memory) 356*9df49d7eSJed Brown /// 357*9df49d7eSJed Brown /// # arguments 358*9df49d7eSJed Brown /// 359*9df49d7eSJed Brown /// * `n` - Length of vector 360*9df49d7eSJed Brown /// 361*9df49d7eSJed Brown /// ``` 362*9df49d7eSJed Brown /// # use libceed::prelude::*; 363*9df49d7eSJed Brown /// # let ceed = libceed::Ceed::default_init(); 364*9df49d7eSJed Brown /// let vec = ceed.vector(10).unwrap(); 365*9df49d7eSJed Brown /// ``` 366*9df49d7eSJed Brown pub fn vector(&self, n: usize) -> Result<Vector> { 367*9df49d7eSJed Brown Vector::create(self, n) 368*9df49d7eSJed Brown } 369*9df49d7eSJed Brown 370*9df49d7eSJed Brown /// Create a Vector initialized with the data (copied) from a slice 371*9df49d7eSJed Brown /// 372*9df49d7eSJed Brown /// # arguments 373*9df49d7eSJed Brown /// 374*9df49d7eSJed Brown /// * `slice` - Slice containing data 375*9df49d7eSJed Brown /// 376*9df49d7eSJed Brown /// ``` 377*9df49d7eSJed Brown /// # use libceed::prelude::*; 378*9df49d7eSJed Brown /// # let ceed = libceed::Ceed::default_init(); 379*9df49d7eSJed Brown /// let vec = ceed.vector_from_slice(&[1., 2., 3.]).unwrap(); 380*9df49d7eSJed Brown /// assert_eq!(vec.length(), 3); 381*9df49d7eSJed Brown /// ``` 382*9df49d7eSJed Brown pub fn vector_from_slice(&self, slice: &[f64]) -> Result<Vector> { 383*9df49d7eSJed Brown Vector::from_slice(self, slice) 384*9df49d7eSJed Brown } 385*9df49d7eSJed Brown 386*9df49d7eSJed Brown /// Returns a ElemRestriction 387*9df49d7eSJed Brown /// 388*9df49d7eSJed Brown /// # arguments 389*9df49d7eSJed Brown /// 390*9df49d7eSJed Brown /// * `nelem` - Number of elements described in the offsets array 391*9df49d7eSJed Brown /// * `elemsize` - Size (number of "nodes") per element 392*9df49d7eSJed Brown /// * `ncomp` - Number of field components per interpolation node (1 393*9df49d7eSJed Brown /// for scalar fields) 394*9df49d7eSJed Brown /// * `compstride` - Stride between components for the same Lvector "node". 395*9df49d7eSJed Brown /// Data for node `i`, component `j`, element `k` can be 396*9df49d7eSJed Brown /// found in the Lvector at index 397*9df49d7eSJed Brown /// `offsets[i + k*elemsize] + j*compstride`. 398*9df49d7eSJed Brown /// * `lsize` - The size of the Lvector. This vector may be larger 399*9df49d7eSJed Brown /// than the elements and fields given by this 400*9df49d7eSJed Brown /// restriction. 401*9df49d7eSJed Brown /// * `mtype` - Memory type of the offsets array, see CeedMemType 402*9df49d7eSJed Brown /// * `offsets` - Array of shape `[nelem, elemsize]`. Row `i` holds the 403*9df49d7eSJed Brown /// ordered list of the offsets (into the input CeedVector) 404*9df49d7eSJed Brown /// for the unknowns corresponding to element `i`, where 405*9df49d7eSJed Brown /// `0 <= i < nelem`. All offsets must be in the range 406*9df49d7eSJed Brown /// `[0, lsize - 1]`. 407*9df49d7eSJed Brown /// 408*9df49d7eSJed Brown /// ``` 409*9df49d7eSJed Brown /// # use libceed::prelude::*; 410*9df49d7eSJed Brown /// # let ceed = libceed::Ceed::default_init(); 411*9df49d7eSJed Brown /// let nelem = 3; 412*9df49d7eSJed Brown /// let mut ind: Vec<i32> = vec![0; 2 * nelem]; 413*9df49d7eSJed Brown /// for i in 0..nelem { 414*9df49d7eSJed Brown /// ind[2 * i + 0] = i as i32; 415*9df49d7eSJed Brown /// ind[2 * i + 1] = (i + 1) as i32; 416*9df49d7eSJed Brown /// } 417*9df49d7eSJed Brown /// let r = ceed 418*9df49d7eSJed Brown /// .elem_restriction(nelem, 2, 1, 1, nelem + 1, MemType::Host, &ind) 419*9df49d7eSJed Brown /// .unwrap(); 420*9df49d7eSJed Brown /// ``` 421*9df49d7eSJed Brown pub fn elem_restriction( 422*9df49d7eSJed Brown &self, 423*9df49d7eSJed Brown nelem: usize, 424*9df49d7eSJed Brown elemsize: usize, 425*9df49d7eSJed Brown ncomp: usize, 426*9df49d7eSJed Brown compstride: usize, 427*9df49d7eSJed Brown lsize: usize, 428*9df49d7eSJed Brown mtype: MemType, 429*9df49d7eSJed Brown offsets: &[i32], 430*9df49d7eSJed Brown ) -> Result<ElemRestriction> { 431*9df49d7eSJed Brown ElemRestriction::create( 432*9df49d7eSJed Brown self, nelem, elemsize, ncomp, compstride, lsize, mtype, offsets, 433*9df49d7eSJed Brown ) 434*9df49d7eSJed Brown } 435*9df49d7eSJed Brown 436*9df49d7eSJed Brown /// Returns a ElemRestriction 437*9df49d7eSJed Brown /// 438*9df49d7eSJed Brown /// # arguments 439*9df49d7eSJed Brown /// 440*9df49d7eSJed Brown /// * `nelem` - Number of elements described in the offsets array 441*9df49d7eSJed Brown /// * `elemsize` - Size (number of "nodes") per element 442*9df49d7eSJed Brown /// * `ncomp` - Number of field components per interpolation node (1 443*9df49d7eSJed Brown /// for scalar fields) 444*9df49d7eSJed Brown /// * `compstride` - Stride between components for the same Lvector "node". 445*9df49d7eSJed Brown /// Data for node `i`, component `j`, element `k` can be 446*9df49d7eSJed Brown /// found in the Lvector at index 447*9df49d7eSJed Brown /// `offsets[i + k*elemsize] + j*compstride`. 448*9df49d7eSJed Brown /// * `lsize` - The size of the Lvector. This vector may be larger 449*9df49d7eSJed Brown /// than the elements and fields given by this restriction. 450*9df49d7eSJed Brown /// * `strides` - Array for strides between `[nodes, components, elements]`. 451*9df49d7eSJed Brown /// Data for node `i`, component `j`, element `k` can be 452*9df49d7eSJed Brown /// found in the Lvector at index 453*9df49d7eSJed Brown /// `i*strides[0] + j*strides[1] + k*strides[2]`. 454*9df49d7eSJed Brown /// CEED_STRIDES_BACKEND may be used with vectors created 455*9df49d7eSJed Brown /// by a Ceed backend. 456*9df49d7eSJed Brown /// 457*9df49d7eSJed Brown /// ``` 458*9df49d7eSJed Brown /// # use libceed::prelude::*; 459*9df49d7eSJed Brown /// # let ceed = libceed::Ceed::default_init(); 460*9df49d7eSJed Brown /// let nelem = 3; 461*9df49d7eSJed Brown /// let strides: [i32; 3] = [1, 2, 2]; 462*9df49d7eSJed Brown /// let r = ceed 463*9df49d7eSJed Brown /// .strided_elem_restriction(nelem, 2, 1, nelem * 2, strides) 464*9df49d7eSJed Brown /// .unwrap(); 465*9df49d7eSJed Brown /// ``` 466*9df49d7eSJed Brown pub fn strided_elem_restriction( 467*9df49d7eSJed Brown &self, 468*9df49d7eSJed Brown nelem: usize, 469*9df49d7eSJed Brown elemsize: usize, 470*9df49d7eSJed Brown ncomp: usize, 471*9df49d7eSJed Brown lsize: usize, 472*9df49d7eSJed Brown strides: [i32; 3], 473*9df49d7eSJed Brown ) -> Result<ElemRestriction> { 474*9df49d7eSJed Brown ElemRestriction::create_strided(self, nelem, elemsize, ncomp, lsize, strides) 475*9df49d7eSJed Brown } 476*9df49d7eSJed Brown 477*9df49d7eSJed Brown /// Returns a tensor-product basis 478*9df49d7eSJed Brown /// 479*9df49d7eSJed Brown /// # arguments 480*9df49d7eSJed Brown /// 481*9df49d7eSJed Brown /// * `dim` - Topological dimension of element 482*9df49d7eSJed Brown /// * `ncomp` - Number of field components (1 for scalar fields) 483*9df49d7eSJed Brown /// * `P1d` - Number of Gauss-Lobatto nodes in one dimension. The 484*9df49d7eSJed Brown /// polynomial degree of the resulting `Q_k` element is 485*9df49d7eSJed Brown /// `k=P-1`. 486*9df49d7eSJed Brown /// * `Q1d` - Number of quadrature points in one dimension 487*9df49d7eSJed Brown /// * `interp1d` - Row-major `(Q1d * P1d)` matrix expressing the values of 488*9df49d7eSJed Brown /// nodal basis functions at quadrature points 489*9df49d7eSJed Brown /// * `grad1d` - Row-major `(Q1d * P1d)` matrix expressing derivatives of 490*9df49d7eSJed Brown /// nodal basis functions at quadrature points 491*9df49d7eSJed Brown /// * `qref1d` - Array of length `Q1d` holding the locations of quadrature 492*9df49d7eSJed Brown /// points on the 1D reference element `[-1, 1]` 493*9df49d7eSJed Brown /// * `qweight1d` - Array of length `Q1d` holding the quadrature weights on 494*9df49d7eSJed Brown /// the reference element 495*9df49d7eSJed Brown /// 496*9df49d7eSJed Brown /// ``` 497*9df49d7eSJed Brown /// # use libceed::prelude::*; 498*9df49d7eSJed Brown /// # let ceed = libceed::Ceed::default_init(); 499*9df49d7eSJed Brown /// let interp1d = [ 0.62994317, 0.47255875, -0.14950343, 0.04700152, 500*9df49d7eSJed Brown /// -0.07069480, 0.97297619, 0.13253993, -0.03482132, 501*9df49d7eSJed Brown /// -0.03482132, 0.13253993, 0.97297619, -0.07069480, 502*9df49d7eSJed Brown /// 0.04700152, -0.14950343, 0.47255875, 0.62994317]; 503*9df49d7eSJed Brown /// let grad1d = [-2.34183742, 2.78794489, -0.63510411, 0.18899664, 504*9df49d7eSJed Brown /// -0.51670214, -0.48795249, 1.33790510, -0.33325047, 505*9df49d7eSJed Brown // 0.33325047, -1.33790510, 0.48795249, 0.51670214, 506*9df49d7eSJed Brown /// -0.18899664, 0.63510411, -2.78794489, 2.34183742]; 507*9df49d7eSJed Brown /// let qref1d = [-0.86113631, -0.33998104, 0.33998104, 0.86113631]; 508*9df49d7eSJed Brown /// let qweight1d = [ 0.34785485, 0.65214515, 0.65214515, 0.34785485]; 509*9df49d7eSJed Brown /// let b = ceed. 510*9df49d7eSJed Brown /// basis_tensor_H1(2, 1, 4, 4, &interp1d, &grad1d, &qref1d, &qweight1d).unwrap(); 511*9df49d7eSJed Brown /// ``` 512*9df49d7eSJed Brown pub fn basis_tensor_H1( 513*9df49d7eSJed Brown &self, 514*9df49d7eSJed Brown dim: usize, 515*9df49d7eSJed Brown ncomp: usize, 516*9df49d7eSJed Brown P1d: usize, 517*9df49d7eSJed Brown Q1d: usize, 518*9df49d7eSJed Brown interp1d: &[f64], 519*9df49d7eSJed Brown grad1d: &[f64], 520*9df49d7eSJed Brown qref1d: &[f64], 521*9df49d7eSJed Brown qweight1d: &[f64], 522*9df49d7eSJed Brown ) -> Result<Basis> { 523*9df49d7eSJed Brown Basis::create_tensor_H1( 524*9df49d7eSJed Brown self, dim, ncomp, P1d, Q1d, interp1d, grad1d, qref1d, qweight1d, 525*9df49d7eSJed Brown ) 526*9df49d7eSJed Brown } 527*9df49d7eSJed Brown 528*9df49d7eSJed Brown /// Returns a tensor-product Lagrange basis 529*9df49d7eSJed Brown /// 530*9df49d7eSJed Brown /// # arguments 531*9df49d7eSJed Brown /// 532*9df49d7eSJed Brown /// * `dim` - Topological dimension of element 533*9df49d7eSJed Brown /// * `ncomp` - Number of field components (1 for scalar fields) 534*9df49d7eSJed Brown /// * `P` - Number of Gauss-Lobatto nodes in one dimension. The 535*9df49d7eSJed Brown /// polynomial degree of the resulting `Q_k` element is `k=P-1`. 536*9df49d7eSJed Brown /// * `Q` - Number of quadrature points in one dimension 537*9df49d7eSJed Brown /// * `qmode` - Distribution of the `Q` quadrature points (affects order of 538*9df49d7eSJed Brown /// accuracy for the quadrature) 539*9df49d7eSJed Brown /// 540*9df49d7eSJed Brown /// ``` 541*9df49d7eSJed Brown /// # use libceed::prelude::*; 542*9df49d7eSJed Brown /// # let ceed = libceed::Ceed::default_init(); 543*9df49d7eSJed Brown /// let b = ceed 544*9df49d7eSJed Brown /// .basis_tensor_H1_Lagrange(2, 1, 3, 4, QuadMode::Gauss) 545*9df49d7eSJed Brown /// .unwrap(); 546*9df49d7eSJed Brown /// ``` 547*9df49d7eSJed Brown pub fn basis_tensor_H1_Lagrange( 548*9df49d7eSJed Brown &self, 549*9df49d7eSJed Brown dim: usize, 550*9df49d7eSJed Brown ncomp: usize, 551*9df49d7eSJed Brown P: usize, 552*9df49d7eSJed Brown Q: usize, 553*9df49d7eSJed Brown qmode: QuadMode, 554*9df49d7eSJed Brown ) -> Result<Basis> { 555*9df49d7eSJed Brown Basis::create_tensor_H1_Lagrange(self, dim, ncomp, P, Q, qmode) 556*9df49d7eSJed Brown } 557*9df49d7eSJed Brown 558*9df49d7eSJed Brown /// Returns a tensor-product basis 559*9df49d7eSJed Brown /// 560*9df49d7eSJed Brown /// # arguments 561*9df49d7eSJed Brown /// 562*9df49d7eSJed Brown /// * `topo` - Topology of element, e.g. hypercube, simplex, ect 563*9df49d7eSJed Brown /// * `ncomp` - Number of field components (1 for scalar fields) 564*9df49d7eSJed Brown /// * `nnodes` - Total number of nodes 565*9df49d7eSJed Brown /// * `nqpts` - Total number of quadrature points 566*9df49d7eSJed Brown /// * `interp` - Row-major `(nqpts * nnodes)` matrix expressing the values of 567*9df49d7eSJed Brown /// nodal basis functions at quadrature points 568*9df49d7eSJed Brown /// * `grad` - Row-major `(nqpts * dim * nnodes)` matrix expressing 569*9df49d7eSJed Brown /// derivatives of nodal basis functions at quadrature points 570*9df49d7eSJed Brown /// * `qref` - Array of length `nqpts` holding the locations of quadrature 571*9df49d7eSJed Brown /// points on the reference element `[-1, 1]` 572*9df49d7eSJed Brown /// * `qweight` - Array of length `nqpts` holding the quadrature weights on 573*9df49d7eSJed Brown /// the reference element 574*9df49d7eSJed Brown /// 575*9df49d7eSJed Brown /// ``` 576*9df49d7eSJed Brown /// # use libceed::prelude::*; 577*9df49d7eSJed Brown /// # let ceed = libceed::Ceed::default_init(); 578*9df49d7eSJed Brown /// let interp = [ 579*9df49d7eSJed Brown /// 0.12000000, 580*9df49d7eSJed Brown /// 0.48000000, 581*9df49d7eSJed Brown /// -0.12000000, 582*9df49d7eSJed Brown /// 0.48000000, 583*9df49d7eSJed Brown /// 0.16000000, 584*9df49d7eSJed Brown /// -0.12000000, 585*9df49d7eSJed Brown /// -0.12000000, 586*9df49d7eSJed Brown /// 0.48000000, 587*9df49d7eSJed Brown /// 0.12000000, 588*9df49d7eSJed Brown /// 0.16000000, 589*9df49d7eSJed Brown /// 0.48000000, 590*9df49d7eSJed Brown /// -0.12000000, 591*9df49d7eSJed Brown /// -0.11111111, 592*9df49d7eSJed Brown /// 0.44444444, 593*9df49d7eSJed Brown /// -0.11111111, 594*9df49d7eSJed Brown /// 0.44444444, 595*9df49d7eSJed Brown /// 0.44444444, 596*9df49d7eSJed Brown /// -0.11111111, 597*9df49d7eSJed Brown /// -0.12000000, 598*9df49d7eSJed Brown /// 0.16000000, 599*9df49d7eSJed Brown /// -0.12000000, 600*9df49d7eSJed Brown /// 0.48000000, 601*9df49d7eSJed Brown /// 0.48000000, 602*9df49d7eSJed Brown /// 0.12000000, 603*9df49d7eSJed Brown /// ]; 604*9df49d7eSJed Brown /// let grad = [ 605*9df49d7eSJed Brown /// -1.40000000, 606*9df49d7eSJed Brown /// 1.60000000, 607*9df49d7eSJed Brown /// -0.20000000, 608*9df49d7eSJed Brown /// -0.80000000, 609*9df49d7eSJed Brown /// 0.80000000, 610*9df49d7eSJed Brown /// 0.00000000, 611*9df49d7eSJed Brown /// 0.20000000, 612*9df49d7eSJed Brown /// -1.60000000, 613*9df49d7eSJed Brown /// 1.40000000, 614*9df49d7eSJed Brown /// -0.80000000, 615*9df49d7eSJed Brown /// 0.80000000, 616*9df49d7eSJed Brown /// 0.00000000, 617*9df49d7eSJed Brown /// -0.33333333, 618*9df49d7eSJed Brown /// 0.00000000, 619*9df49d7eSJed Brown /// 0.33333333, 620*9df49d7eSJed Brown /// -1.33333333, 621*9df49d7eSJed Brown /// 1.33333333, 622*9df49d7eSJed Brown /// 0.00000000, 623*9df49d7eSJed Brown /// 0.20000000, 624*9df49d7eSJed Brown /// 0.00000000, 625*9df49d7eSJed Brown /// -0.20000000, 626*9df49d7eSJed Brown /// -2.40000000, 627*9df49d7eSJed Brown /// 2.40000000, 628*9df49d7eSJed Brown /// 0.00000000, 629*9df49d7eSJed Brown /// -1.40000000, 630*9df49d7eSJed Brown /// -0.80000000, 631*9df49d7eSJed Brown /// 0.00000000, 632*9df49d7eSJed Brown /// 1.60000000, 633*9df49d7eSJed Brown /// 0.80000000, 634*9df49d7eSJed Brown /// -0.20000000, 635*9df49d7eSJed Brown /// 0.20000000, 636*9df49d7eSJed Brown /// -2.40000000, 637*9df49d7eSJed Brown /// 0.00000000, 638*9df49d7eSJed Brown /// 0.00000000, 639*9df49d7eSJed Brown /// 2.40000000, 640*9df49d7eSJed Brown /// -0.20000000, 641*9df49d7eSJed Brown /// -0.33333333, 642*9df49d7eSJed Brown /// -1.33333333, 643*9df49d7eSJed Brown /// 0.00000000, 644*9df49d7eSJed Brown /// 0.00000000, 645*9df49d7eSJed Brown /// 1.33333333, 646*9df49d7eSJed Brown /// 0.33333333, 647*9df49d7eSJed Brown /// 0.20000000, 648*9df49d7eSJed Brown /// -0.80000000, 649*9df49d7eSJed Brown /// 0.00000000, 650*9df49d7eSJed Brown /// -1.60000000, 651*9df49d7eSJed Brown /// 0.80000000, 652*9df49d7eSJed Brown /// 1.40000000, 653*9df49d7eSJed Brown /// ]; 654*9df49d7eSJed Brown /// let qref = [ 655*9df49d7eSJed Brown /// 0.20000000, 0.60000000, 0.33333333, 0.20000000, 0.20000000, 0.20000000, 0.33333333, 656*9df49d7eSJed Brown /// 0.60000000, 657*9df49d7eSJed Brown /// ]; 658*9df49d7eSJed Brown /// let qweight = [0.26041667, 0.26041667, -0.28125000, 0.26041667]; 659*9df49d7eSJed Brown /// let b = ceed 660*9df49d7eSJed Brown /// .basis_H1( 661*9df49d7eSJed Brown /// ElemTopology::Triangle, 662*9df49d7eSJed Brown /// 1, 663*9df49d7eSJed Brown /// 6, 664*9df49d7eSJed Brown /// 4, 665*9df49d7eSJed Brown /// &interp, 666*9df49d7eSJed Brown /// &grad, 667*9df49d7eSJed Brown /// &qref, 668*9df49d7eSJed Brown /// &qweight, 669*9df49d7eSJed Brown /// ) 670*9df49d7eSJed Brown /// .unwrap(); 671*9df49d7eSJed Brown /// ``` 672*9df49d7eSJed Brown pub fn basis_H1( 673*9df49d7eSJed Brown &self, 674*9df49d7eSJed Brown topo: ElemTopology, 675*9df49d7eSJed Brown ncomp: usize, 676*9df49d7eSJed Brown nnodes: usize, 677*9df49d7eSJed Brown nqpts: usize, 678*9df49d7eSJed Brown interp: &[f64], 679*9df49d7eSJed Brown grad: &[f64], 680*9df49d7eSJed Brown qref: &[f64], 681*9df49d7eSJed Brown qweight: &[f64], 682*9df49d7eSJed Brown ) -> Result<Basis> { 683*9df49d7eSJed Brown Basis::create_H1( 684*9df49d7eSJed Brown self, topo, ncomp, nnodes, nqpts, interp, grad, qref, qweight, 685*9df49d7eSJed Brown ) 686*9df49d7eSJed Brown } 687*9df49d7eSJed Brown 688*9df49d7eSJed Brown /// Returns a CeedQFunction for evaluating interior (volumetric) terms 689*9df49d7eSJed Brown /// 690*9df49d7eSJed Brown /// # arguments 691*9df49d7eSJed Brown /// 692*9df49d7eSJed Brown /// * `vlength` - Vector length. Caller must ensure that number of 693*9df49d7eSJed Brown /// quadrature points is a multiple of vlength. 694*9df49d7eSJed Brown /// * `f` - Boxed closure to evaluate action at quadrature points. 695*9df49d7eSJed Brown /// 696*9df49d7eSJed Brown /// ``` 697*9df49d7eSJed Brown /// # use libceed::prelude::*; 698*9df49d7eSJed Brown /// # let ceed = libceed::Ceed::default_init(); 699*9df49d7eSJed Brown /// let mut user_f = |[u, weights, ..]: QFunctionInputs, [v, ..]: QFunctionOutputs| { 700*9df49d7eSJed Brown /// // Iterate over quadrature points 701*9df49d7eSJed Brown /// v.iter_mut() 702*9df49d7eSJed Brown /// .zip(u.iter().zip(weights.iter())) 703*9df49d7eSJed Brown /// .for_each(|(v, (u, w))| *v = u * w); 704*9df49d7eSJed Brown /// 705*9df49d7eSJed Brown /// // Return clean error code 706*9df49d7eSJed Brown /// 0 707*9df49d7eSJed Brown /// }; 708*9df49d7eSJed Brown /// 709*9df49d7eSJed Brown /// let qf = ceed.q_function_interior(1, Box::new(user_f)).unwrap(); 710*9df49d7eSJed Brown /// ``` 711*9df49d7eSJed Brown pub fn q_function_interior( 712*9df49d7eSJed Brown &self, 713*9df49d7eSJed Brown vlength: usize, 714*9df49d7eSJed Brown f: Box<qfunction::QFunctionUserClosure>, 715*9df49d7eSJed Brown ) -> Result<QFunction> { 716*9df49d7eSJed Brown QFunction::create(self, vlength, f) 717*9df49d7eSJed Brown } 718*9df49d7eSJed Brown 719*9df49d7eSJed Brown /// Returns a CeedQFunction for evaluating interior (volumetric) terms 720*9df49d7eSJed Brown /// created by name 721*9df49d7eSJed Brown /// 722*9df49d7eSJed Brown /// ``` 723*9df49d7eSJed Brown /// # use libceed::prelude::*; 724*9df49d7eSJed Brown /// # let ceed = libceed::Ceed::default_init(); 725*9df49d7eSJed Brown /// let qf = ceed.q_function_interior_by_name("Mass1DBuild").unwrap(); 726*9df49d7eSJed Brown /// ``` 727*9df49d7eSJed Brown pub fn q_function_interior_by_name(&self, name: &str) -> Result<QFunctionByName> { 728*9df49d7eSJed Brown QFunctionByName::create(self, name) 729*9df49d7eSJed Brown } 730*9df49d7eSJed Brown 731*9df49d7eSJed Brown /// Returns a Operator and associate a QFunction. A Basis and 732*9df49d7eSJed Brown /// ElemRestriction can be associated with QFunction fields with 733*9df49d7eSJed Brown /// set_field(). 734*9df49d7eSJed Brown /// 735*9df49d7eSJed Brown /// * `qf` - QFunction defining the action of the operator at quadrature 736*9df49d7eSJed Brown /// points 737*9df49d7eSJed Brown /// * `dqf` - QFunction defining the action of the Jacobian of the qf (or 738*9df49d7eSJed Brown /// qfunction_none) 739*9df49d7eSJed Brown /// * `dqfT` - QFunction defining the action of the transpose of the 740*9df49d7eSJed Brown /// Jacobian of the qf (or qfunction_none) 741*9df49d7eSJed Brown /// 742*9df49d7eSJed Brown /// ``` 743*9df49d7eSJed Brown /// # use libceed::prelude::*; 744*9df49d7eSJed Brown /// # let ceed = libceed::Ceed::default_init(); 745*9df49d7eSJed Brown /// let qf = ceed.q_function_interior_by_name("Mass1DBuild").unwrap(); 746*9df49d7eSJed Brown /// let op = ceed 747*9df49d7eSJed Brown /// .operator(&qf, QFunctionOpt::None, QFunctionOpt::None) 748*9df49d7eSJed Brown /// .unwrap(); 749*9df49d7eSJed Brown /// ``` 750*9df49d7eSJed Brown pub fn operator<'b>( 751*9df49d7eSJed Brown &self, 752*9df49d7eSJed Brown qf: impl Into<QFunctionOpt<'b>>, 753*9df49d7eSJed Brown dqf: impl Into<QFunctionOpt<'b>>, 754*9df49d7eSJed Brown dqfT: impl Into<QFunctionOpt<'b>>, 755*9df49d7eSJed Brown ) -> Result<Operator> { 756*9df49d7eSJed Brown Operator::create(self, qf, dqf, dqfT) 757*9df49d7eSJed Brown } 758*9df49d7eSJed Brown 759*9df49d7eSJed Brown /// Returns an Operator that composes the action of several Operators 760*9df49d7eSJed Brown /// 761*9df49d7eSJed Brown /// ``` 762*9df49d7eSJed Brown /// # use libceed::prelude::*; 763*9df49d7eSJed Brown /// # let ceed = libceed::Ceed::default_init(); 764*9df49d7eSJed Brown /// let op = ceed.composite_operator().unwrap(); 765*9df49d7eSJed Brown /// ``` 766*9df49d7eSJed Brown pub fn composite_operator(&self) -> Result<CompositeOperator> { 767*9df49d7eSJed Brown CompositeOperator::create(self) 768*9df49d7eSJed Brown } 769*9df49d7eSJed Brown } 770*9df49d7eSJed Brown 771*9df49d7eSJed Brown // ----------------------------------------------------------------------------- 772*9df49d7eSJed Brown // Tests 773*9df49d7eSJed Brown // ----------------------------------------------------------------------------- 774*9df49d7eSJed Brown #[cfg(test)] 775*9df49d7eSJed Brown mod tests { 776*9df49d7eSJed Brown use super::*; 777*9df49d7eSJed Brown 778*9df49d7eSJed Brown fn ceed_t501() -> Result<i32> { 779*9df49d7eSJed Brown let resource = "/cpu/self/ref/blocked"; 780*9df49d7eSJed Brown let ceed = Ceed::init(resource); 781*9df49d7eSJed Brown let nelem = 4; 782*9df49d7eSJed Brown let p = 3; 783*9df49d7eSJed Brown let q = 4; 784*9df49d7eSJed Brown let ndofs = p * nelem - nelem + 1; 785*9df49d7eSJed Brown 786*9df49d7eSJed Brown // Vectors 787*9df49d7eSJed Brown let x = ceed.vector_from_slice(&[-1., -0.5, 0.0, 0.5, 1.0])?; 788*9df49d7eSJed Brown let mut qdata = ceed.vector(nelem * q)?; 789*9df49d7eSJed Brown qdata.set_value(0.0)?; 790*9df49d7eSJed Brown let mut u = ceed.vector(ndofs)?; 791*9df49d7eSJed Brown u.set_value(1.0)?; 792*9df49d7eSJed Brown let mut v = ceed.vector(ndofs)?; 793*9df49d7eSJed Brown v.set_value(0.0)?; 794*9df49d7eSJed Brown 795*9df49d7eSJed Brown // Restrictions 796*9df49d7eSJed Brown let mut indx: Vec<i32> = vec![0; 2 * nelem]; 797*9df49d7eSJed Brown for i in 0..nelem { 798*9df49d7eSJed Brown indx[2 * i + 0] = i as i32; 799*9df49d7eSJed Brown indx[2 * i + 1] = (i + 1) as i32; 800*9df49d7eSJed Brown } 801*9df49d7eSJed Brown let rx = ceed.elem_restriction(nelem, 2, 1, 1, nelem + 1, MemType::Host, &indx)?; 802*9df49d7eSJed Brown let mut indu: Vec<i32> = vec![0; p * nelem]; 803*9df49d7eSJed Brown for i in 0..nelem { 804*9df49d7eSJed Brown indu[p * i + 0] = i as i32; 805*9df49d7eSJed Brown indu[p * i + 1] = (i + 1) as i32; 806*9df49d7eSJed Brown indu[p * i + 2] = (i + 2) as i32; 807*9df49d7eSJed Brown } 808*9df49d7eSJed Brown let ru = ceed.elem_restriction(nelem, 3, 1, 1, ndofs, MemType::Host, &indu)?; 809*9df49d7eSJed Brown let strides: [i32; 3] = [1, q as i32, q as i32]; 810*9df49d7eSJed Brown let rq = ceed.strided_elem_restriction(nelem, q, 1, q * nelem, strides)?; 811*9df49d7eSJed Brown 812*9df49d7eSJed Brown // Bases 813*9df49d7eSJed Brown let bx = ceed.basis_tensor_H1_Lagrange(1, 1, 2, q, QuadMode::Gauss)?; 814*9df49d7eSJed Brown let bu = ceed.basis_tensor_H1_Lagrange(1, 1, p, q, QuadMode::Gauss)?; 815*9df49d7eSJed Brown 816*9df49d7eSJed Brown // Build quadrature data 817*9df49d7eSJed Brown let qf_build = ceed.q_function_interior_by_name("Mass1DBuild")?; 818*9df49d7eSJed Brown ceed.operator(&qf_build, QFunctionOpt::None, QFunctionOpt::None)? 819*9df49d7eSJed Brown .field("dx", &rx, &bx, VectorOpt::Active)? 820*9df49d7eSJed Brown .field("weights", ElemRestrictionOpt::None, &bx, VectorOpt::None)? 821*9df49d7eSJed Brown .field("qdata", &rq, BasisOpt::Collocated, VectorOpt::Active)? 822*9df49d7eSJed Brown .apply(&x, &mut qdata)?; 823*9df49d7eSJed Brown 824*9df49d7eSJed Brown // Mass operator 825*9df49d7eSJed Brown let qf_mass = ceed.q_function_interior_by_name("MassApply")?; 826*9df49d7eSJed Brown let op_mass = ceed 827*9df49d7eSJed Brown .operator(&qf_mass, QFunctionOpt::None, QFunctionOpt::None)? 828*9df49d7eSJed Brown .field("u", &ru, &bu, VectorOpt::Active)? 829*9df49d7eSJed Brown .field("qdata", &rq, BasisOpt::Collocated, &qdata)? 830*9df49d7eSJed Brown .field("v", &ru, &bu, VectorOpt::Active)?; 831*9df49d7eSJed Brown 832*9df49d7eSJed Brown v.set_value(0.0)?; 833*9df49d7eSJed Brown op_mass.apply(&u, &mut v)?; 834*9df49d7eSJed Brown 835*9df49d7eSJed Brown // Check 836*9df49d7eSJed Brown let sum: f64 = v.view().iter().sum(); 837*9df49d7eSJed Brown assert!( 838*9df49d7eSJed Brown (sum - 2.0).abs() < 1e-15, 839*9df49d7eSJed Brown "Incorrect interval length computed" 840*9df49d7eSJed Brown ); 841*9df49d7eSJed Brown Ok(0) 842*9df49d7eSJed Brown } 843*9df49d7eSJed Brown 844*9df49d7eSJed Brown #[test] 845*9df49d7eSJed Brown fn test_ceed_t501() { 846*9df49d7eSJed Brown assert!(ceed_t501().is_ok()); 847*9df49d7eSJed Brown } 848*9df49d7eSJed Brown } 849*9df49d7eSJed Brown 850*9df49d7eSJed Brown // ----------------------------------------------------------------------------- 851