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 17a1cbad85SJed Brown // Fenced `rust` code blocks included from README.md are executed as part of doctests. 18a1cbad85SJed Brown #![doc = include_str!("../README.md")] 199df49d7eSJed Brown // ----------------------------------------------------------------------------- 209df49d7eSJed Brown // Exceptions 219df49d7eSJed Brown // ----------------------------------------------------------------------------- 229df49d7eSJed Brown #![allow(non_snake_case)] 239df49d7eSJed Brown 249df49d7eSJed Brown // ----------------------------------------------------------------------------- 259df49d7eSJed Brown // Crate prelude 269df49d7eSJed Brown // ----------------------------------------------------------------------------- 279df49d7eSJed Brown use crate::prelude::*; 289df49d7eSJed Brown use std::sync::Once; 299df49d7eSJed Brown 309df49d7eSJed Brown pub mod prelude { 319df49d7eSJed Brown pub use crate::{ 329df49d7eSJed Brown basis::{self, Basis, BasisOpt}, 339df49d7eSJed Brown elem_restriction::{self, ElemRestriction, ElemRestrictionOpt}, 3408778c6fSJeremy L Thompson operator::{self, CompositeOperator, Operator, OperatorField}, 359df49d7eSJed Brown qfunction::{ 3608778c6fSJeremy L Thompson self, QFunction, QFunctionByName, QFunctionField, QFunctionInputs, QFunctionOpt, 3708778c6fSJeremy L Thompson QFunctionOutputs, 389df49d7eSJed Brown }, 39486868d3SJeremy L Thompson vector::{self, Vector, VectorOpt, VectorSliceWrapper}, 402ba8e59cSJeremy L Thompson ElemTopology, EvalMode, MemType, NormType, QuadMode, Scalar, TransposeMode, 4180a9ef05SNatalie Beams CEED_STRIDES_BACKEND, EPSILON, MAX_QFUNCTION_FIELDS, 429df49d7eSJed Brown }; 43630ad4c9Sjeremylt pub(crate) use libceed_sys::bind_ceed; 449df49d7eSJed Brown pub(crate) use std::convert::TryFrom; 4508778c6fSJeremy L Thompson pub(crate) use std::ffi::{CStr, CString}; 469df49d7eSJed Brown pub(crate) use std::fmt; 471142270cSJeremy L Thompson pub(crate) use std::marker::PhantomData; 489df49d7eSJed Brown } 499df49d7eSJed Brown 509df49d7eSJed Brown // ----------------------------------------------------------------------------- 519df49d7eSJed Brown // Modules 529df49d7eSJed Brown // ----------------------------------------------------------------------------- 539df49d7eSJed Brown pub mod basis; 549df49d7eSJed Brown pub mod elem_restriction; 559df49d7eSJed Brown pub mod operator; 569df49d7eSJed Brown pub mod qfunction; 579df49d7eSJed Brown pub mod vector; 589df49d7eSJed Brown 599df49d7eSJed Brown // ----------------------------------------------------------------------------- 6080a9ef05SNatalie Beams // Typedef for scalar 6180a9ef05SNatalie Beams // ----------------------------------------------------------------------------- 6280a9ef05SNatalie Beams pub type Scalar = bind_ceed::CeedScalar; 6380a9ef05SNatalie Beams 6480a9ef05SNatalie Beams // ----------------------------------------------------------------------------- 659df49d7eSJed Brown // Constants for library 669df49d7eSJed Brown // ----------------------------------------------------------------------------- 679df49d7eSJed Brown const MAX_BUFFER_LENGTH: u64 = 4096; 689df49d7eSJed Brown pub const MAX_QFUNCTION_FIELDS: usize = 16; 699df49d7eSJed Brown pub const CEED_STRIDES_BACKEND: [i32; 3] = [0; 3]; 7080a9ef05SNatalie Beams pub const EPSILON: crate::Scalar = bind_ceed::CEED_EPSILON as crate::Scalar; 719df49d7eSJed Brown 729df49d7eSJed Brown // ----------------------------------------------------------------------------- 739df49d7eSJed Brown // Enums for libCEED 749df49d7eSJed Brown // ----------------------------------------------------------------------------- 759df49d7eSJed Brown #[derive(Clone, Copy, PartialEq, Eq)] 769df49d7eSJed Brown /// Many Ceed interfaces take or return pointers to memory. This enum is used to 779df49d7eSJed Brown /// specify where the memory being provided or requested must reside. 789df49d7eSJed Brown pub enum MemType { 799df49d7eSJed Brown Host = bind_ceed::CeedMemType_CEED_MEM_HOST as isize, 809df49d7eSJed Brown Device = bind_ceed::CeedMemType_CEED_MEM_DEVICE as isize, 819df49d7eSJed Brown } 829df49d7eSJed Brown 839df49d7eSJed Brown #[derive(Clone, Copy, PartialEq, Eq)] 849df49d7eSJed Brown // OwnPointer will not be used by user but is included for internal use 859df49d7eSJed Brown #[allow(dead_code)] 869df49d7eSJed Brown /// Conveys ownership status of arrays passed to Ceed interfaces. 879df49d7eSJed Brown pub(crate) enum CopyMode { 889df49d7eSJed Brown CopyValues = bind_ceed::CeedCopyMode_CEED_COPY_VALUES as isize, 899df49d7eSJed Brown UsePointer = bind_ceed::CeedCopyMode_CEED_USE_POINTER as isize, 909df49d7eSJed Brown OwnPointer = bind_ceed::CeedCopyMode_CEED_OWN_POINTER as isize, 919df49d7eSJed Brown } 929df49d7eSJed Brown 939df49d7eSJed Brown #[derive(Clone, Copy, PartialEq, Eq)] 949df49d7eSJed Brown /// Denotes type of vector norm to be computed 959df49d7eSJed Brown pub enum NormType { 969df49d7eSJed Brown One = bind_ceed::CeedNormType_CEED_NORM_1 as isize, 979df49d7eSJed Brown Two = bind_ceed::CeedNormType_CEED_NORM_2 as isize, 989df49d7eSJed Brown Max = bind_ceed::CeedNormType_CEED_NORM_MAX as isize, 999df49d7eSJed Brown } 1009df49d7eSJed Brown 1019df49d7eSJed Brown #[derive(Clone, Copy, PartialEq, Eq)] 1029df49d7eSJed Brown /// Denotes whether a linear transformation or its transpose should be applied 1039df49d7eSJed Brown pub enum TransposeMode { 1049df49d7eSJed Brown NoTranspose = bind_ceed::CeedTransposeMode_CEED_NOTRANSPOSE as isize, 1059df49d7eSJed Brown Transpose = bind_ceed::CeedTransposeMode_CEED_TRANSPOSE as isize, 1069df49d7eSJed Brown } 1079df49d7eSJed Brown 1089df49d7eSJed Brown #[derive(Clone, Copy, PartialEq, Eq)] 1099df49d7eSJed Brown /// Type of quadrature; also used for location of nodes 1109df49d7eSJed Brown pub enum QuadMode { 1119df49d7eSJed Brown Gauss = bind_ceed::CeedQuadMode_CEED_GAUSS as isize, 1129df49d7eSJed Brown GaussLobatto = bind_ceed::CeedQuadMode_CEED_GAUSS_LOBATTO as isize, 1139df49d7eSJed Brown } 1149df49d7eSJed Brown 1159df49d7eSJed Brown #[derive(Clone, Copy, PartialEq, Eq)] 1169df49d7eSJed Brown /// Type of basis shape to create non-tensor H1 element basis 1179df49d7eSJed Brown pub enum ElemTopology { 1189df49d7eSJed Brown Line = bind_ceed::CeedElemTopology_CEED_LINE as isize, 1199df49d7eSJed Brown Triangle = bind_ceed::CeedElemTopology_CEED_TRIANGLE as isize, 1209df49d7eSJed Brown Quad = bind_ceed::CeedElemTopology_CEED_QUAD as isize, 1219df49d7eSJed Brown Tet = bind_ceed::CeedElemTopology_CEED_TET as isize, 1229df49d7eSJed Brown Pyramid = bind_ceed::CeedElemTopology_CEED_PYRAMID as isize, 1239df49d7eSJed Brown Prism = bind_ceed::CeedElemTopology_CEED_PRISM as isize, 1249df49d7eSJed Brown Hex = bind_ceed::CeedElemTopology_CEED_HEX as isize, 1259df49d7eSJed Brown } 1269df49d7eSJed Brown 127c68be7a2SJeremy L Thompson #[derive(Clone, Copy, Debug, PartialEq, Eq)] 1289df49d7eSJed Brown /// Basis evaluation mode 1299df49d7eSJed Brown pub enum EvalMode { 1309df49d7eSJed Brown None = bind_ceed::CeedEvalMode_CEED_EVAL_NONE as isize, 1319df49d7eSJed Brown Interp = bind_ceed::CeedEvalMode_CEED_EVAL_INTERP as isize, 1329df49d7eSJed Brown Grad = bind_ceed::CeedEvalMode_CEED_EVAL_GRAD as isize, 1339df49d7eSJed Brown Div = bind_ceed::CeedEvalMode_CEED_EVAL_DIV as isize, 1349df49d7eSJed Brown Curl = bind_ceed::CeedEvalMode_CEED_EVAL_CURL as isize, 1359df49d7eSJed Brown Weight = bind_ceed::CeedEvalMode_CEED_EVAL_WEIGHT as isize, 1369df49d7eSJed Brown } 137c68be7a2SJeremy L Thompson impl EvalMode { 138c68be7a2SJeremy L Thompson pub(crate) fn from_u32(value: u32) -> EvalMode { 139c68be7a2SJeremy L Thompson match value { 140c68be7a2SJeremy L Thompson bind_ceed::CeedEvalMode_CEED_EVAL_NONE => EvalMode::None, 141c68be7a2SJeremy L Thompson bind_ceed::CeedEvalMode_CEED_EVAL_INTERP => EvalMode::Interp, 142c68be7a2SJeremy L Thompson bind_ceed::CeedEvalMode_CEED_EVAL_GRAD => EvalMode::Grad, 143c68be7a2SJeremy L Thompson bind_ceed::CeedEvalMode_CEED_EVAL_DIV => EvalMode::Div, 144c68be7a2SJeremy L Thompson bind_ceed::CeedEvalMode_CEED_EVAL_CURL => EvalMode::Curl, 145c68be7a2SJeremy L Thompson bind_ceed::CeedEvalMode_CEED_EVAL_WEIGHT => EvalMode::Weight, 146c68be7a2SJeremy L Thompson _ => panic!("Unknown value: {}", value), 147c68be7a2SJeremy L Thompson } 148c68be7a2SJeremy L Thompson } 149c68be7a2SJeremy L Thompson } 1509df49d7eSJed Brown 1519df49d7eSJed Brown // ----------------------------------------------------------------------------- 1529df49d7eSJed Brown // Ceed error 1539df49d7eSJed Brown // ----------------------------------------------------------------------------- 1542ba8e59cSJeremy L Thompson pub type Result<T> = std::result::Result<T, Error>; 1559df49d7eSJed Brown 1569df49d7eSJed Brown #[derive(Debug)] 1572ba8e59cSJeremy L Thompson pub struct Error { 1589df49d7eSJed Brown pub message: String, 1599df49d7eSJed Brown } 1609df49d7eSJed Brown 1612ba8e59cSJeremy L Thompson impl fmt::Display for Error { 1629df49d7eSJed Brown fn fmt(&self, f: &mut fmt::Formatter) -> fmt::Result { 1639df49d7eSJed Brown write!(f, "{}", self.message) 1649df49d7eSJed Brown } 1659df49d7eSJed Brown } 1669df49d7eSJed Brown 1679df49d7eSJed Brown // ----------------------------------------------------------------------------- 1681142270cSJeremy L Thompson // Internal error checker 1691142270cSJeremy L Thompson // ----------------------------------------------------------------------------- 1701142270cSJeremy L Thompson #[doc(hidden)] 1711142270cSJeremy L Thompson pub(crate) fn check_error(ceed_ptr: bind_ceed::Ceed, ierr: i32) -> Result<i32> { 1721142270cSJeremy L Thompson // Return early if code is clean 1731142270cSJeremy L Thompson if ierr == bind_ceed::CeedErrorType_CEED_ERROR_SUCCESS { 1741142270cSJeremy L Thompson return Ok(ierr); 1751142270cSJeremy L Thompson } 1761142270cSJeremy L Thompson // Retrieve error message 1771142270cSJeremy L Thompson let mut ptr: *const std::os::raw::c_char = std::ptr::null_mut(); 1781142270cSJeremy L Thompson let c_str = unsafe { 1791142270cSJeremy L Thompson bind_ceed::CeedGetErrorMessage(ceed_ptr, &mut ptr); 1801142270cSJeremy L Thompson std::ffi::CStr::from_ptr(ptr) 1811142270cSJeremy L Thompson }; 1821142270cSJeremy L Thompson let message = c_str.to_string_lossy().to_string(); 1831142270cSJeremy L Thompson // Panic if negative code, otherwise return error 1841142270cSJeremy L Thompson if ierr < bind_ceed::CeedErrorType_CEED_ERROR_SUCCESS { 1851142270cSJeremy L Thompson panic!("{}", message); 1861142270cSJeremy L Thompson } 1872ba8e59cSJeremy L Thompson Err(Error { message }) 1881142270cSJeremy L Thompson } 1891142270cSJeremy L Thompson 1901142270cSJeremy L Thompson // ----------------------------------------------------------------------------- 1919df49d7eSJed Brown // Ceed error handler 1929df49d7eSJed Brown // ----------------------------------------------------------------------------- 1932ba8e59cSJeremy L Thompson pub enum ErrorHandler { 1949df49d7eSJed Brown ErrorAbort, 1959df49d7eSJed Brown ErrorExit, 1969df49d7eSJed Brown ErrorReturn, 1979df49d7eSJed Brown ErrorStore, 1989df49d7eSJed Brown } 1999df49d7eSJed Brown 2009df49d7eSJed Brown // ----------------------------------------------------------------------------- 2019df49d7eSJed Brown // Ceed context wrapper 2029df49d7eSJed Brown // ----------------------------------------------------------------------------- 2039df49d7eSJed Brown /// A Ceed is a library context representing control of a logical hardware 2049df49d7eSJed Brown /// resource. 2059df49d7eSJed Brown #[derive(Debug)] 2069df49d7eSJed Brown pub struct Ceed { 2079df49d7eSJed Brown ptr: bind_ceed::Ceed, 2089df49d7eSJed Brown } 2099df49d7eSJed Brown 2109df49d7eSJed Brown // ----------------------------------------------------------------------------- 2119df49d7eSJed Brown // Destructor 2129df49d7eSJed Brown // ----------------------------------------------------------------------------- 2139df49d7eSJed Brown impl Drop for Ceed { 2149df49d7eSJed Brown fn drop(&mut self) { 2159df49d7eSJed Brown unsafe { 2169df49d7eSJed Brown bind_ceed::CeedDestroy(&mut self.ptr); 2179df49d7eSJed Brown } 2189df49d7eSJed Brown } 2199df49d7eSJed Brown } 2209df49d7eSJed Brown 2219df49d7eSJed Brown // ----------------------------------------------------------------------------- 222*59189cfaSJeremy L Thompson // Cloning 223*59189cfaSJeremy L Thompson // ----------------------------------------------------------------------------- 224*59189cfaSJeremy L Thompson impl Clone for Ceed { 225*59189cfaSJeremy L Thompson /// Perform a shallow clone of a Ceed context 226*59189cfaSJeremy L Thompson /// 227*59189cfaSJeremy L Thompson /// ``` 228*59189cfaSJeremy L Thompson /// let ceed = libceed::Ceed::init("/cpu/self/ref/serial"); 229*59189cfaSJeremy L Thompson /// let ceed_clone = ceed.clone(); 230*59189cfaSJeremy L Thompson /// 231*59189cfaSJeremy L Thompson /// println!("{}", ceed); 232*59189cfaSJeremy L Thompson /// println!("{}", ceed_clone); 233*59189cfaSJeremy L Thompson /// ``` 234*59189cfaSJeremy L Thompson fn clone(&self) -> Self { 235*59189cfaSJeremy L Thompson let mut ptr_clone = std::ptr::null_mut(); 236*59189cfaSJeremy L Thompson let ierr = unsafe { bind_ceed::CeedReferenceCopy(self.ptr, &mut ptr_clone) }; 237*59189cfaSJeremy L Thompson self.check_error(ierr).expect("failed to clone Ceed"); 238*59189cfaSJeremy L Thompson Self { ptr: ptr_clone } 239*59189cfaSJeremy L Thompson } 240*59189cfaSJeremy L Thompson } 241*59189cfaSJeremy L Thompson 242*59189cfaSJeremy L Thompson // ----------------------------------------------------------------------------- 2439df49d7eSJed Brown // Display 2449df49d7eSJed Brown // ----------------------------------------------------------------------------- 2459df49d7eSJed Brown impl fmt::Display for Ceed { 2469df49d7eSJed Brown /// View a Ceed 2479df49d7eSJed Brown /// 2489df49d7eSJed Brown /// ``` 2499df49d7eSJed Brown /// let ceed = libceed::Ceed::init("/cpu/self/ref/serial"); 2509df49d7eSJed Brown /// println!("{}", ceed); 2519df49d7eSJed Brown /// ``` 2529df49d7eSJed Brown fn fmt(&self, f: &mut fmt::Formatter) -> fmt::Result { 2539df49d7eSJed Brown let mut ptr = std::ptr::null_mut(); 2549df49d7eSJed Brown let mut sizeloc = crate::MAX_BUFFER_LENGTH; 2559df49d7eSJed Brown let cstring = unsafe { 2569df49d7eSJed Brown let file = bind_ceed::open_memstream(&mut ptr, &mut sizeloc); 2579df49d7eSJed Brown bind_ceed::CeedView(self.ptr, file); 2589df49d7eSJed Brown bind_ceed::fclose(file); 2599df49d7eSJed Brown CString::from_raw(ptr) 2609df49d7eSJed Brown }; 2619df49d7eSJed Brown cstring.to_string_lossy().fmt(f) 2629df49d7eSJed Brown } 2639df49d7eSJed Brown } 2649df49d7eSJed Brown 2659df49d7eSJed Brown static REGISTER: Once = Once::new(); 2669df49d7eSJed Brown 2679df49d7eSJed Brown // ----------------------------------------------------------------------------- 2689df49d7eSJed Brown // Object constructors 2699df49d7eSJed Brown // ----------------------------------------------------------------------------- 2709df49d7eSJed Brown impl Ceed { 2719df49d7eSJed Brown /// Returns a Ceed context initialized with the specified resource 2729df49d7eSJed Brown /// 2739df49d7eSJed Brown /// # arguments 2749df49d7eSJed Brown /// 2759df49d7eSJed Brown /// * `resource` - Resource to use, e.g., "/cpu/self" 2769df49d7eSJed Brown /// 2779df49d7eSJed Brown /// ``` 2789df49d7eSJed Brown /// let ceed = libceed::Ceed::init("/cpu/self/ref/serial"); 2799df49d7eSJed Brown /// ``` 2809df49d7eSJed Brown pub fn init(resource: &str) -> Self { 2812ba8e59cSJeremy L Thompson Ceed::init_with_error_handler(resource, ErrorHandler::ErrorStore) 2829df49d7eSJed Brown } 2839df49d7eSJed Brown 2849df49d7eSJed Brown /// Returns a Ceed context initialized with the specified resource 2859df49d7eSJed Brown /// 2869df49d7eSJed Brown /// # arguments 2879df49d7eSJed Brown /// 2889df49d7eSJed Brown /// * `resource` - Resource to use, e.g., "/cpu/self" 2899df49d7eSJed Brown /// 2909df49d7eSJed Brown /// ``` 2919df49d7eSJed Brown /// let ceed = libceed::Ceed::init_with_error_handler( 2929df49d7eSJed Brown /// "/cpu/self/ref/serial", 2932ba8e59cSJeremy L Thompson /// libceed::ErrorHandler::ErrorAbort, 2949df49d7eSJed Brown /// ); 2959df49d7eSJed Brown /// ``` 2962ba8e59cSJeremy L Thompson pub fn init_with_error_handler(resource: &str, handler: ErrorHandler) -> Self { 2979df49d7eSJed Brown REGISTER.call_once(|| unsafe { 2989df49d7eSJed Brown bind_ceed::CeedRegisterAll(); 2999df49d7eSJed Brown bind_ceed::CeedQFunctionRegisterAll(); 3009df49d7eSJed Brown }); 3019df49d7eSJed Brown 3029df49d7eSJed Brown // Convert to C string 3039df49d7eSJed Brown let c_resource = CString::new(resource).expect("CString::new failed"); 3049df49d7eSJed Brown 3059df49d7eSJed Brown // Get error handler pointer 3069df49d7eSJed Brown let eh = match handler { 3072ba8e59cSJeremy L Thompson ErrorHandler::ErrorAbort => bind_ceed::CeedErrorAbort, 3082ba8e59cSJeremy L Thompson ErrorHandler::ErrorExit => bind_ceed::CeedErrorExit, 3092ba8e59cSJeremy L Thompson ErrorHandler::ErrorReturn => bind_ceed::CeedErrorReturn, 3102ba8e59cSJeremy L Thompson ErrorHandler::ErrorStore => bind_ceed::CeedErrorStore, 3119df49d7eSJed Brown }; 3129df49d7eSJed Brown 3139df49d7eSJed Brown // Call to libCEED 3149df49d7eSJed Brown let mut ptr = std::ptr::null_mut(); 3159df49d7eSJed Brown let mut ierr = unsafe { bind_ceed::CeedInit(c_resource.as_ptr() as *const i8, &mut ptr) }; 3169df49d7eSJed Brown if ierr != 0 { 3179df49d7eSJed Brown panic!("Error initializing backend resource: {}", resource) 3189df49d7eSJed Brown } 3199df49d7eSJed Brown ierr = unsafe { bind_ceed::CeedSetErrorHandler(ptr, Some(eh)) }; 3209df49d7eSJed Brown let ceed = Ceed { ptr }; 3219df49d7eSJed Brown ceed.check_error(ierr).unwrap(); 3229df49d7eSJed Brown ceed 3239df49d7eSJed Brown } 3249df49d7eSJed Brown 3259df49d7eSJed Brown /// Default initializer for testing 3269df49d7eSJed Brown #[doc(hidden)] 3279df49d7eSJed Brown pub fn default_init() -> Self { 3289df49d7eSJed Brown // Convert to C string 3299df49d7eSJed Brown let resource = "/cpu/self/ref/serial"; 3309df49d7eSJed Brown crate::Ceed::init(resource) 3319df49d7eSJed Brown } 3329df49d7eSJed Brown 3339df49d7eSJed Brown /// Internal error checker 3349df49d7eSJed Brown #[doc(hidden)] 3359df49d7eSJed Brown fn check_error(&self, ierr: i32) -> Result<i32> { 3369df49d7eSJed Brown // Return early if code is clean 3379df49d7eSJed Brown if ierr == bind_ceed::CeedErrorType_CEED_ERROR_SUCCESS { 3389df49d7eSJed Brown return Ok(ierr); 3399df49d7eSJed Brown } 3409df49d7eSJed Brown // Retrieve error message 3419df49d7eSJed Brown let mut ptr: *const std::os::raw::c_char = std::ptr::null_mut(); 3429df49d7eSJed Brown let c_str = unsafe { 3439df49d7eSJed Brown bind_ceed::CeedGetErrorMessage(self.ptr, &mut ptr); 3449df49d7eSJed Brown std::ffi::CStr::from_ptr(ptr) 3459df49d7eSJed Brown }; 3469df49d7eSJed Brown let message = c_str.to_string_lossy().to_string(); 3479df49d7eSJed Brown // Panic if negative code, otherwise return error 3489df49d7eSJed Brown if ierr < bind_ceed::CeedErrorType_CEED_ERROR_SUCCESS { 3499df49d7eSJed Brown panic!("{}", message); 3509df49d7eSJed Brown } 3512ba8e59cSJeremy L Thompson Err(Error { message }) 3529df49d7eSJed Brown } 3539df49d7eSJed Brown 3549df49d7eSJed Brown /// Returns full resource name for a Ceed context 3559df49d7eSJed Brown /// 3569df49d7eSJed Brown /// ``` 3579df49d7eSJed Brown /// let ceed = libceed::Ceed::init("/cpu/self/ref/serial"); 3589df49d7eSJed Brown /// let resource = ceed.resource(); 3599df49d7eSJed Brown /// 3609df49d7eSJed Brown /// assert_eq!(resource, "/cpu/self/ref/serial".to_string()) 3619df49d7eSJed Brown /// ``` 3629df49d7eSJed Brown pub fn resource(&self) -> String { 3639df49d7eSJed Brown let mut ptr: *const std::os::raw::c_char = std::ptr::null_mut(); 3649df49d7eSJed Brown let c_str = unsafe { 3659df49d7eSJed Brown bind_ceed::CeedGetResource(self.ptr, &mut ptr); 3669df49d7eSJed Brown std::ffi::CStr::from_ptr(ptr) 3679df49d7eSJed Brown }; 3689df49d7eSJed Brown c_str.to_string_lossy().to_string() 3699df49d7eSJed Brown } 3709df49d7eSJed Brown 3719df49d7eSJed Brown /// Returns a CeedVector of the specified length (does not allocate memory) 3729df49d7eSJed Brown /// 3739df49d7eSJed Brown /// # arguments 3749df49d7eSJed Brown /// 3759df49d7eSJed Brown /// * `n` - Length of vector 3769df49d7eSJed Brown /// 3779df49d7eSJed Brown /// ``` 3789df49d7eSJed Brown /// # use libceed::prelude::*; 3794d27c890SJeremy L Thompson /// # fn main() -> libceed::Result<()> { 3809df49d7eSJed Brown /// # let ceed = libceed::Ceed::default_init(); 381c68be7a2SJeremy L Thompson /// let vec = ceed.vector(10)?; 382c68be7a2SJeremy L Thompson /// # Ok(()) 383c68be7a2SJeremy L Thompson /// # } 3849df49d7eSJed Brown /// ``` 385594ef120SJeremy L Thompson pub fn vector<'a>(&self, n: usize) -> Result<Vector<'a>> { 3869df49d7eSJed Brown Vector::create(self, n) 3879df49d7eSJed Brown } 3889df49d7eSJed Brown 3899df49d7eSJed Brown /// Create a Vector initialized with the data (copied) from a slice 3909df49d7eSJed Brown /// 3919df49d7eSJed Brown /// # arguments 3929df49d7eSJed Brown /// 3939df49d7eSJed Brown /// * `slice` - Slice containing data 3949df49d7eSJed Brown /// 3959df49d7eSJed Brown /// ``` 3969df49d7eSJed Brown /// # use libceed::prelude::*; 3974d27c890SJeremy L Thompson /// # fn main() -> libceed::Result<()> { 3989df49d7eSJed Brown /// # let ceed = libceed::Ceed::default_init(); 399c68be7a2SJeremy L Thompson /// let vec = ceed.vector_from_slice(&[1., 2., 3.])?; 4009df49d7eSJed Brown /// assert_eq!(vec.length(), 3); 401c68be7a2SJeremy L Thompson /// # Ok(()) 402c68be7a2SJeremy L Thompson /// # } 4039df49d7eSJed Brown /// ``` 404594ef120SJeremy L Thompson pub fn vector_from_slice<'a>(&self, slice: &[crate::Scalar]) -> Result<Vector<'a>> { 4059df49d7eSJed Brown Vector::from_slice(self, slice) 4069df49d7eSJed Brown } 4079df49d7eSJed Brown 4089df49d7eSJed Brown /// Returns a ElemRestriction 4099df49d7eSJed Brown /// 4109df49d7eSJed Brown /// # arguments 4119df49d7eSJed Brown /// 4129df49d7eSJed Brown /// * `nelem` - Number of elements described in the offsets array 4139df49d7eSJed Brown /// * `elemsize` - Size (number of "nodes") per element 4149df49d7eSJed Brown /// * `ncomp` - Number of field components per interpolation node (1 4159df49d7eSJed Brown /// for scalar fields) 4169df49d7eSJed Brown /// * `compstride` - Stride between components for the same Lvector "node". 4179df49d7eSJed Brown /// Data for node `i`, component `j`, element `k` can be 4189df49d7eSJed Brown /// found in the Lvector at index 4199df49d7eSJed Brown /// `offsets[i + k*elemsize] + j*compstride`. 4209df49d7eSJed Brown /// * `lsize` - The size of the Lvector. This vector may be larger 4219df49d7eSJed Brown /// than the elements and fields given by this 4229df49d7eSJed Brown /// restriction. 4239df49d7eSJed Brown /// * `mtype` - Memory type of the offsets array, see CeedMemType 4249df49d7eSJed Brown /// * `offsets` - Array of shape `[nelem, elemsize]`. Row `i` holds the 4259df49d7eSJed Brown /// ordered list of the offsets (into the input CeedVector) 4269df49d7eSJed Brown /// for the unknowns corresponding to element `i`, where 4279df49d7eSJed Brown /// `0 <= i < nelem`. All offsets must be in the range 4289df49d7eSJed Brown /// `[0, lsize - 1]`. 4299df49d7eSJed Brown /// 4309df49d7eSJed Brown /// ``` 4319df49d7eSJed Brown /// # use libceed::prelude::*; 4324d27c890SJeremy L Thompson /// # fn main() -> libceed::Result<()> { 4339df49d7eSJed Brown /// # let ceed = libceed::Ceed::default_init(); 4349df49d7eSJed Brown /// let nelem = 3; 4359df49d7eSJed Brown /// let mut ind: Vec<i32> = vec![0; 2 * nelem]; 4369df49d7eSJed Brown /// for i in 0..nelem { 4379df49d7eSJed Brown /// ind[2 * i + 0] = i as i32; 4389df49d7eSJed Brown /// ind[2 * i + 1] = (i + 1) as i32; 4399df49d7eSJed Brown /// } 440c68be7a2SJeremy L Thompson /// let r = ceed.elem_restriction(nelem, 2, 1, 1, nelem + 1, MemType::Host, &ind)?; 441c68be7a2SJeremy L Thompson /// # Ok(()) 442c68be7a2SJeremy L Thompson /// # } 4439df49d7eSJed Brown /// ``` 444594ef120SJeremy L Thompson pub fn elem_restriction<'a>( 4459df49d7eSJed Brown &self, 4469df49d7eSJed Brown nelem: usize, 4479df49d7eSJed Brown elemsize: usize, 4489df49d7eSJed Brown ncomp: usize, 4499df49d7eSJed Brown compstride: usize, 4509df49d7eSJed Brown lsize: usize, 4519df49d7eSJed Brown mtype: MemType, 4529df49d7eSJed Brown offsets: &[i32], 453594ef120SJeremy L Thompson ) -> Result<ElemRestriction<'a>> { 4549df49d7eSJed Brown ElemRestriction::create( 4559df49d7eSJed Brown self, nelem, elemsize, ncomp, compstride, lsize, mtype, offsets, 4569df49d7eSJed Brown ) 4579df49d7eSJed Brown } 4589df49d7eSJed Brown 4599df49d7eSJed Brown /// Returns a ElemRestriction 4609df49d7eSJed Brown /// 4619df49d7eSJed Brown /// # arguments 4629df49d7eSJed Brown /// 4639df49d7eSJed Brown /// * `nelem` - Number of elements described in the offsets array 4649df49d7eSJed Brown /// * `elemsize` - Size (number of "nodes") per element 4659df49d7eSJed Brown /// * `ncomp` - Number of field components per interpolation node (1 4669df49d7eSJed Brown /// for scalar fields) 4679df49d7eSJed Brown /// * `compstride` - Stride between components for the same Lvector "node". 4689df49d7eSJed Brown /// Data for node `i`, component `j`, element `k` can be 4699df49d7eSJed Brown /// found in the Lvector at index 4709df49d7eSJed Brown /// `offsets[i + k*elemsize] + j*compstride`. 4719df49d7eSJed Brown /// * `lsize` - The size of the Lvector. This vector may be larger 4729df49d7eSJed Brown /// than the elements and fields given by this restriction. 4739df49d7eSJed Brown /// * `strides` - Array for strides between `[nodes, components, elements]`. 4749df49d7eSJed Brown /// Data for node `i`, component `j`, element `k` can be 4759df49d7eSJed Brown /// found in the Lvector at index 4769df49d7eSJed Brown /// `i*strides[0] + j*strides[1] + k*strides[2]`. 4779df49d7eSJed Brown /// CEED_STRIDES_BACKEND may be used with vectors created 4789df49d7eSJed Brown /// by a Ceed backend. 4799df49d7eSJed Brown /// 4809df49d7eSJed Brown /// ``` 4819df49d7eSJed Brown /// # use libceed::prelude::*; 4824d27c890SJeremy L Thompson /// # fn main() -> libceed::Result<()> { 4839df49d7eSJed Brown /// # let ceed = libceed::Ceed::default_init(); 4849df49d7eSJed Brown /// let nelem = 3; 4859df49d7eSJed Brown /// let strides: [i32; 3] = [1, 2, 2]; 486c68be7a2SJeremy L Thompson /// let r = ceed.strided_elem_restriction(nelem, 2, 1, nelem * 2, strides)?; 487c68be7a2SJeremy L Thompson /// # Ok(()) 488c68be7a2SJeremy L Thompson /// # } 4899df49d7eSJed Brown /// ``` 490594ef120SJeremy L Thompson pub fn strided_elem_restriction<'a>( 4919df49d7eSJed Brown &self, 4929df49d7eSJed Brown nelem: usize, 4939df49d7eSJed Brown elemsize: usize, 4949df49d7eSJed Brown ncomp: usize, 4959df49d7eSJed Brown lsize: usize, 4969df49d7eSJed Brown strides: [i32; 3], 497594ef120SJeremy L Thompson ) -> Result<ElemRestriction<'a>> { 4989df49d7eSJed Brown ElemRestriction::create_strided(self, nelem, elemsize, ncomp, lsize, strides) 4999df49d7eSJed Brown } 5009df49d7eSJed Brown 5019df49d7eSJed Brown /// Returns a tensor-product basis 5029df49d7eSJed Brown /// 5039df49d7eSJed Brown /// # arguments 5049df49d7eSJed Brown /// 5059df49d7eSJed Brown /// * `dim` - Topological dimension of element 5069df49d7eSJed Brown /// * `ncomp` - Number of field components (1 for scalar fields) 5079df49d7eSJed Brown /// * `P1d` - Number of Gauss-Lobatto nodes in one dimension. The 5089df49d7eSJed Brown /// polynomial degree of the resulting `Q_k` element is 5099df49d7eSJed Brown /// `k=P-1`. 5109df49d7eSJed Brown /// * `Q1d` - Number of quadrature points in one dimension 5119df49d7eSJed Brown /// * `interp1d` - Row-major `(Q1d * P1d)` matrix expressing the values of 5129df49d7eSJed Brown /// nodal basis functions at quadrature points 5139df49d7eSJed Brown /// * `grad1d` - Row-major `(Q1d * P1d)` matrix expressing derivatives of 5149df49d7eSJed Brown /// nodal basis functions at quadrature points 5159df49d7eSJed Brown /// * `qref1d` - Array of length `Q1d` holding the locations of quadrature 5169df49d7eSJed Brown /// points on the 1D reference element `[-1, 1]` 5179df49d7eSJed Brown /// * `qweight1d` - Array of length `Q1d` holding the quadrature weights on 5189df49d7eSJed Brown /// the reference element 5199df49d7eSJed Brown /// 5209df49d7eSJed Brown /// ``` 5219df49d7eSJed Brown /// # use libceed::prelude::*; 5224d27c890SJeremy L Thompson /// # fn main() -> libceed::Result<()> { 5239df49d7eSJed Brown /// # let ceed = libceed::Ceed::default_init(); 5249df49d7eSJed Brown /// let interp1d = [ 0.62994317, 0.47255875, -0.14950343, 0.04700152, 5259df49d7eSJed Brown /// -0.07069480, 0.97297619, 0.13253993, -0.03482132, 5269df49d7eSJed Brown /// -0.03482132, 0.13253993, 0.97297619, -0.07069480, 5279df49d7eSJed Brown /// 0.04700152, -0.14950343, 0.47255875, 0.62994317]; 5289df49d7eSJed Brown /// let grad1d = [-2.34183742, 2.78794489, -0.63510411, 0.18899664, 5299df49d7eSJed Brown /// -0.51670214, -0.48795249, 1.33790510, -0.33325047, 5309df49d7eSJed Brown // 0.33325047, -1.33790510, 0.48795249, 0.51670214, 5319df49d7eSJed Brown /// -0.18899664, 0.63510411, -2.78794489, 2.34183742]; 5329df49d7eSJed Brown /// let qref1d = [-0.86113631, -0.33998104, 0.33998104, 0.86113631]; 5339df49d7eSJed Brown /// let qweight1d = [ 0.34785485, 0.65214515, 0.65214515, 0.34785485]; 5349df49d7eSJed Brown /// let b = ceed. 535c68be7a2SJeremy L Thompson /// basis_tensor_H1(2, 1, 4, 4, &interp1d, &grad1d, &qref1d, &qweight1d)?; 536c68be7a2SJeremy L Thompson /// # Ok(()) 537c68be7a2SJeremy L Thompson /// # } 5389df49d7eSJed Brown /// ``` 539594ef120SJeremy L Thompson pub fn basis_tensor_H1<'a>( 5409df49d7eSJed Brown &self, 5419df49d7eSJed Brown dim: usize, 5429df49d7eSJed Brown ncomp: usize, 5439df49d7eSJed Brown P1d: usize, 5449df49d7eSJed Brown Q1d: usize, 54580a9ef05SNatalie Beams interp1d: &[crate::Scalar], 54680a9ef05SNatalie Beams grad1d: &[crate::Scalar], 54780a9ef05SNatalie Beams qref1d: &[crate::Scalar], 54880a9ef05SNatalie Beams qweight1d: &[crate::Scalar], 549594ef120SJeremy L Thompson ) -> Result<Basis<'a>> { 5509df49d7eSJed Brown Basis::create_tensor_H1( 5519df49d7eSJed Brown self, dim, ncomp, P1d, Q1d, interp1d, grad1d, qref1d, qweight1d, 5529df49d7eSJed Brown ) 5539df49d7eSJed Brown } 5549df49d7eSJed Brown 5559df49d7eSJed Brown /// Returns a tensor-product Lagrange basis 5569df49d7eSJed Brown /// 5579df49d7eSJed Brown /// # arguments 5589df49d7eSJed Brown /// 5599df49d7eSJed Brown /// * `dim` - Topological dimension of element 5609df49d7eSJed Brown /// * `ncomp` - Number of field components (1 for scalar fields) 5619df49d7eSJed Brown /// * `P` - Number of Gauss-Lobatto nodes in one dimension. The 5629df49d7eSJed Brown /// polynomial degree of the resulting `Q_k` element is `k=P-1`. 5639df49d7eSJed Brown /// * `Q` - Number of quadrature points in one dimension 5649df49d7eSJed Brown /// * `qmode` - Distribution of the `Q` quadrature points (affects order of 5659df49d7eSJed Brown /// accuracy for the quadrature) 5669df49d7eSJed Brown /// 5679df49d7eSJed Brown /// ``` 5689df49d7eSJed Brown /// # use libceed::prelude::*; 5694d27c890SJeremy L Thompson /// # fn main() -> libceed::Result<()> { 5709df49d7eSJed Brown /// # let ceed = libceed::Ceed::default_init(); 571c68be7a2SJeremy L Thompson /// let b = ceed.basis_tensor_H1_Lagrange(2, 1, 3, 4, QuadMode::Gauss)?; 572c68be7a2SJeremy L Thompson /// # Ok(()) 573c68be7a2SJeremy L Thompson /// # } 5749df49d7eSJed Brown /// ``` 575594ef120SJeremy L Thompson pub fn basis_tensor_H1_Lagrange<'a>( 5769df49d7eSJed Brown &self, 5779df49d7eSJed Brown dim: usize, 5789df49d7eSJed Brown ncomp: usize, 5799df49d7eSJed Brown P: usize, 5809df49d7eSJed Brown Q: usize, 5819df49d7eSJed Brown qmode: QuadMode, 582594ef120SJeremy L Thompson ) -> Result<Basis<'a>> { 5839df49d7eSJed Brown Basis::create_tensor_H1_Lagrange(self, dim, ncomp, P, Q, qmode) 5849df49d7eSJed Brown } 5859df49d7eSJed Brown 5869df49d7eSJed Brown /// Returns a tensor-product basis 5879df49d7eSJed Brown /// 5889df49d7eSJed Brown /// # arguments 5899df49d7eSJed Brown /// 5909df49d7eSJed Brown /// * `topo` - Topology of element, e.g. hypercube, simplex, ect 5919df49d7eSJed Brown /// * `ncomp` - Number of field components (1 for scalar fields) 5929df49d7eSJed Brown /// * `nnodes` - Total number of nodes 5939df49d7eSJed Brown /// * `nqpts` - Total number of quadrature points 5949df49d7eSJed Brown /// * `interp` - Row-major `(nqpts * nnodes)` matrix expressing the values of 5959df49d7eSJed Brown /// nodal basis functions at quadrature points 5969df49d7eSJed Brown /// * `grad` - Row-major `(nqpts * dim * nnodes)` matrix expressing 5979df49d7eSJed Brown /// derivatives of nodal basis functions at quadrature points 5989df49d7eSJed Brown /// * `qref` - Array of length `nqpts` holding the locations of quadrature 5999df49d7eSJed Brown /// points on the reference element `[-1, 1]` 6009df49d7eSJed Brown /// * `qweight` - Array of length `nqpts` holding the quadrature weights on 6019df49d7eSJed Brown /// the reference element 6029df49d7eSJed Brown /// 6039df49d7eSJed Brown /// ``` 6049df49d7eSJed Brown /// # use libceed::prelude::*; 6054d27c890SJeremy L Thompson /// # fn main() -> libceed::Result<()> { 6069df49d7eSJed Brown /// # let ceed = libceed::Ceed::default_init(); 6079df49d7eSJed Brown /// let interp = [ 6089df49d7eSJed Brown /// 0.12000000, 6099df49d7eSJed Brown /// 0.48000000, 6109df49d7eSJed Brown /// -0.12000000, 6119df49d7eSJed Brown /// 0.48000000, 6129df49d7eSJed Brown /// 0.16000000, 6139df49d7eSJed Brown /// -0.12000000, 6149df49d7eSJed Brown /// -0.12000000, 6159df49d7eSJed Brown /// 0.48000000, 6169df49d7eSJed Brown /// 0.12000000, 6179df49d7eSJed Brown /// 0.16000000, 6189df49d7eSJed Brown /// 0.48000000, 6199df49d7eSJed Brown /// -0.12000000, 6209df49d7eSJed Brown /// -0.11111111, 6219df49d7eSJed Brown /// 0.44444444, 6229df49d7eSJed Brown /// -0.11111111, 6239df49d7eSJed Brown /// 0.44444444, 6249df49d7eSJed Brown /// 0.44444444, 6259df49d7eSJed Brown /// -0.11111111, 6269df49d7eSJed Brown /// -0.12000000, 6279df49d7eSJed Brown /// 0.16000000, 6289df49d7eSJed Brown /// -0.12000000, 6299df49d7eSJed Brown /// 0.48000000, 6309df49d7eSJed Brown /// 0.48000000, 6319df49d7eSJed Brown /// 0.12000000, 6329df49d7eSJed Brown /// ]; 6339df49d7eSJed Brown /// let grad = [ 6349df49d7eSJed Brown /// -1.40000000, 6359df49d7eSJed Brown /// 1.60000000, 6369df49d7eSJed Brown /// -0.20000000, 6379df49d7eSJed Brown /// -0.80000000, 6389df49d7eSJed Brown /// 0.80000000, 6399df49d7eSJed Brown /// 0.00000000, 6409df49d7eSJed Brown /// 0.20000000, 6419df49d7eSJed Brown /// -1.60000000, 6429df49d7eSJed Brown /// 1.40000000, 6439df49d7eSJed Brown /// -0.80000000, 6449df49d7eSJed Brown /// 0.80000000, 6459df49d7eSJed Brown /// 0.00000000, 6469df49d7eSJed Brown /// -0.33333333, 6479df49d7eSJed Brown /// 0.00000000, 6489df49d7eSJed Brown /// 0.33333333, 6499df49d7eSJed Brown /// -1.33333333, 6509df49d7eSJed Brown /// 1.33333333, 6519df49d7eSJed Brown /// 0.00000000, 6529df49d7eSJed Brown /// 0.20000000, 6539df49d7eSJed Brown /// 0.00000000, 6549df49d7eSJed Brown /// -0.20000000, 6559df49d7eSJed Brown /// -2.40000000, 6569df49d7eSJed Brown /// 2.40000000, 6579df49d7eSJed Brown /// 0.00000000, 6589df49d7eSJed Brown /// -1.40000000, 6599df49d7eSJed Brown /// -0.80000000, 6609df49d7eSJed Brown /// 0.00000000, 6619df49d7eSJed Brown /// 1.60000000, 6629df49d7eSJed Brown /// 0.80000000, 6639df49d7eSJed Brown /// -0.20000000, 6649df49d7eSJed Brown /// 0.20000000, 6659df49d7eSJed Brown /// -2.40000000, 6669df49d7eSJed Brown /// 0.00000000, 6679df49d7eSJed Brown /// 0.00000000, 6689df49d7eSJed Brown /// 2.40000000, 6699df49d7eSJed Brown /// -0.20000000, 6709df49d7eSJed Brown /// -0.33333333, 6719df49d7eSJed Brown /// -1.33333333, 6729df49d7eSJed Brown /// 0.00000000, 6739df49d7eSJed Brown /// 0.00000000, 6749df49d7eSJed Brown /// 1.33333333, 6759df49d7eSJed Brown /// 0.33333333, 6769df49d7eSJed Brown /// 0.20000000, 6779df49d7eSJed Brown /// -0.80000000, 6789df49d7eSJed Brown /// 0.00000000, 6799df49d7eSJed Brown /// -1.60000000, 6809df49d7eSJed Brown /// 0.80000000, 6819df49d7eSJed Brown /// 1.40000000, 6829df49d7eSJed Brown /// ]; 6839df49d7eSJed Brown /// let qref = [ 6849df49d7eSJed Brown /// 0.20000000, 0.60000000, 0.33333333, 0.20000000, 0.20000000, 0.20000000, 0.33333333, 6859df49d7eSJed Brown /// 0.60000000, 6869df49d7eSJed Brown /// ]; 6879df49d7eSJed Brown /// let qweight = [0.26041667, 0.26041667, -0.28125000, 0.26041667]; 688c68be7a2SJeremy L Thompson /// let b = ceed.basis_H1( 6899df49d7eSJed Brown /// ElemTopology::Triangle, 6909df49d7eSJed Brown /// 1, 6919df49d7eSJed Brown /// 6, 6929df49d7eSJed Brown /// 4, 6939df49d7eSJed Brown /// &interp, 6949df49d7eSJed Brown /// &grad, 6959df49d7eSJed Brown /// &qref, 6969df49d7eSJed Brown /// &qweight, 697c68be7a2SJeremy L Thompson /// )?; 698c68be7a2SJeremy L Thompson /// # Ok(()) 699c68be7a2SJeremy L Thompson /// # } 7009df49d7eSJed Brown /// ``` 701594ef120SJeremy L Thompson pub fn basis_H1<'a>( 7029df49d7eSJed Brown &self, 7039df49d7eSJed Brown topo: ElemTopology, 7049df49d7eSJed Brown ncomp: usize, 7059df49d7eSJed Brown nnodes: usize, 7069df49d7eSJed Brown nqpts: usize, 70780a9ef05SNatalie Beams interp: &[crate::Scalar], 70880a9ef05SNatalie Beams grad: &[crate::Scalar], 70980a9ef05SNatalie Beams qref: &[crate::Scalar], 71080a9ef05SNatalie Beams qweight: &[crate::Scalar], 711594ef120SJeremy L Thompson ) -> Result<Basis<'a>> { 7129df49d7eSJed Brown Basis::create_H1( 7139df49d7eSJed Brown self, topo, ncomp, nnodes, nqpts, interp, grad, qref, qweight, 7149df49d7eSJed Brown ) 7159df49d7eSJed Brown } 7169df49d7eSJed Brown 7179df49d7eSJed Brown /// Returns a CeedQFunction for evaluating interior (volumetric) terms 7189df49d7eSJed Brown /// 7199df49d7eSJed Brown /// # arguments 7209df49d7eSJed Brown /// 7219df49d7eSJed Brown /// * `vlength` - Vector length. Caller must ensure that number of 7229df49d7eSJed Brown /// quadrature points is a multiple of vlength. 7239df49d7eSJed Brown /// * `f` - Boxed closure to evaluate action at quadrature points. 7249df49d7eSJed Brown /// 7259df49d7eSJed Brown /// ``` 7269df49d7eSJed Brown /// # use libceed::prelude::*; 7274d27c890SJeremy L Thompson /// # fn main() -> libceed::Result<()> { 7289df49d7eSJed Brown /// # let ceed = libceed::Ceed::default_init(); 7299df49d7eSJed Brown /// let mut user_f = |[u, weights, ..]: QFunctionInputs, [v, ..]: QFunctionOutputs| { 7309df49d7eSJed Brown /// // Iterate over quadrature points 7319df49d7eSJed Brown /// v.iter_mut() 7329df49d7eSJed Brown /// .zip(u.iter().zip(weights.iter())) 7339df49d7eSJed Brown /// .for_each(|(v, (u, w))| *v = u * w); 7349df49d7eSJed Brown /// 7359df49d7eSJed Brown /// // Return clean error code 7369df49d7eSJed Brown /// 0 7379df49d7eSJed Brown /// }; 7389df49d7eSJed Brown /// 739c68be7a2SJeremy L Thompson /// let qf = ceed.q_function_interior(1, Box::new(user_f))?; 740c68be7a2SJeremy L Thompson /// # Ok(()) 741c68be7a2SJeremy L Thompson /// # } 7429df49d7eSJed Brown /// ``` 743594ef120SJeremy L Thompson pub fn q_function_interior<'a>( 7449df49d7eSJed Brown &self, 7459df49d7eSJed Brown vlength: usize, 7469df49d7eSJed Brown f: Box<qfunction::QFunctionUserClosure>, 747594ef120SJeremy L Thompson ) -> Result<QFunction<'a>> { 7489df49d7eSJed Brown QFunction::create(self, vlength, f) 7499df49d7eSJed Brown } 7509df49d7eSJed Brown 7519df49d7eSJed Brown /// Returns a CeedQFunction for evaluating interior (volumetric) terms 7529df49d7eSJed Brown /// created by name 7539df49d7eSJed Brown /// 7549df49d7eSJed Brown /// ``` 7559df49d7eSJed Brown /// # use libceed::prelude::*; 7564d27c890SJeremy L Thompson /// # fn main() -> libceed::Result<()> { 7579df49d7eSJed Brown /// # let ceed = libceed::Ceed::default_init(); 758c68be7a2SJeremy L Thompson /// let qf = ceed.q_function_interior_by_name("Mass1DBuild")?; 759c68be7a2SJeremy L Thompson /// # Ok(()) 760c68be7a2SJeremy L Thompson /// # } 7619df49d7eSJed Brown /// ``` 762594ef120SJeremy L Thompson pub fn q_function_interior_by_name<'a>(&self, name: &str) -> Result<QFunctionByName<'a>> { 7639df49d7eSJed Brown QFunctionByName::create(self, name) 7649df49d7eSJed Brown } 7659df49d7eSJed Brown 7669df49d7eSJed Brown /// Returns a Operator and associate a QFunction. A Basis and 7679df49d7eSJed Brown /// ElemRestriction can be associated with QFunction fields with 7689df49d7eSJed Brown /// set_field(). 7699df49d7eSJed Brown /// 7709df49d7eSJed Brown /// * `qf` - QFunction defining the action of the operator at quadrature 7719df49d7eSJed Brown /// points 7729df49d7eSJed Brown /// * `dqf` - QFunction defining the action of the Jacobian of the qf (or 7739df49d7eSJed Brown /// qfunction_none) 7749df49d7eSJed Brown /// * `dqfT` - QFunction defining the action of the transpose of the 7759df49d7eSJed Brown /// Jacobian of the qf (or qfunction_none) 7769df49d7eSJed Brown /// 7779df49d7eSJed Brown /// ``` 7789df49d7eSJed Brown /// # use libceed::prelude::*; 7794d27c890SJeremy L Thompson /// # fn main() -> libceed::Result<()> { 7809df49d7eSJed Brown /// # let ceed = libceed::Ceed::default_init(); 781c68be7a2SJeremy L Thompson /// let qf = ceed.q_function_interior_by_name("Mass1DBuild")?; 782c68be7a2SJeremy L Thompson /// let op = ceed.operator(&qf, QFunctionOpt::None, QFunctionOpt::None)?; 783c68be7a2SJeremy L Thompson /// # Ok(()) 784c68be7a2SJeremy L Thompson /// # } 7859df49d7eSJed Brown /// ``` 786594ef120SJeremy L Thompson pub fn operator<'a, 'b>( 7879df49d7eSJed Brown &self, 7889df49d7eSJed Brown qf: impl Into<QFunctionOpt<'b>>, 7899df49d7eSJed Brown dqf: impl Into<QFunctionOpt<'b>>, 7909df49d7eSJed Brown dqfT: impl Into<QFunctionOpt<'b>>, 791594ef120SJeremy L Thompson ) -> Result<Operator<'a>> { 7929df49d7eSJed Brown Operator::create(self, qf, dqf, dqfT) 7939df49d7eSJed Brown } 7949df49d7eSJed Brown 7959df49d7eSJed Brown /// Returns an Operator that composes the action of several Operators 7969df49d7eSJed Brown /// 7979df49d7eSJed Brown /// ``` 7989df49d7eSJed Brown /// # use libceed::prelude::*; 7994d27c890SJeremy L Thompson /// # fn main() -> libceed::Result<()> { 8009df49d7eSJed Brown /// # let ceed = libceed::Ceed::default_init(); 801c68be7a2SJeremy L Thompson /// let op = ceed.composite_operator()?; 802c68be7a2SJeremy L Thompson /// # Ok(()) 803c68be7a2SJeremy L Thompson /// # } 8049df49d7eSJed Brown /// ``` 805594ef120SJeremy L Thompson pub fn composite_operator<'a>(&self) -> Result<CompositeOperator<'a>> { 8069df49d7eSJed Brown CompositeOperator::create(self) 8079df49d7eSJed Brown } 8089df49d7eSJed Brown } 8099df49d7eSJed Brown 8109df49d7eSJed Brown // ----------------------------------------------------------------------------- 8119df49d7eSJed Brown // Tests 8129df49d7eSJed Brown // ----------------------------------------------------------------------------- 8139df49d7eSJed Brown #[cfg(test)] 8149df49d7eSJed Brown mod tests { 8159df49d7eSJed Brown use super::*; 8169df49d7eSJed Brown 81789d15d5fSJeremy L Thompson fn ceed_t501() -> Result<()> { 8189df49d7eSJed Brown let resource = "/cpu/self/ref/blocked"; 8199df49d7eSJed Brown let ceed = Ceed::init(resource); 8209df49d7eSJed Brown let nelem = 4; 8219df49d7eSJed Brown let p = 3; 8229df49d7eSJed Brown let q = 4; 8239df49d7eSJed Brown let ndofs = p * nelem - nelem + 1; 8249df49d7eSJed Brown 8259df49d7eSJed Brown // Vectors 8269df49d7eSJed Brown let x = ceed.vector_from_slice(&[-1., -0.5, 0.0, 0.5, 1.0])?; 8279df49d7eSJed Brown let mut qdata = ceed.vector(nelem * q)?; 8289df49d7eSJed Brown qdata.set_value(0.0)?; 8299df49d7eSJed Brown let mut u = ceed.vector(ndofs)?; 8309df49d7eSJed Brown u.set_value(1.0)?; 8319df49d7eSJed Brown let mut v = ceed.vector(ndofs)?; 8329df49d7eSJed Brown v.set_value(0.0)?; 8339df49d7eSJed Brown 8349df49d7eSJed Brown // Restrictions 8359df49d7eSJed Brown let mut indx: Vec<i32> = vec![0; 2 * nelem]; 8369df49d7eSJed Brown for i in 0..nelem { 8379df49d7eSJed Brown indx[2 * i + 0] = i as i32; 8389df49d7eSJed Brown indx[2 * i + 1] = (i + 1) as i32; 8399df49d7eSJed Brown } 8409df49d7eSJed Brown let rx = ceed.elem_restriction(nelem, 2, 1, 1, nelem + 1, MemType::Host, &indx)?; 8419df49d7eSJed Brown let mut indu: Vec<i32> = vec![0; p * nelem]; 8429df49d7eSJed Brown for i in 0..nelem { 8439df49d7eSJed Brown indu[p * i + 0] = i as i32; 8449df49d7eSJed Brown indu[p * i + 1] = (i + 1) as i32; 8459df49d7eSJed Brown indu[p * i + 2] = (i + 2) as i32; 8469df49d7eSJed Brown } 8479df49d7eSJed Brown let ru = ceed.elem_restriction(nelem, 3, 1, 1, ndofs, MemType::Host, &indu)?; 8489df49d7eSJed Brown let strides: [i32; 3] = [1, q as i32, q as i32]; 8499df49d7eSJed Brown let rq = ceed.strided_elem_restriction(nelem, q, 1, q * nelem, strides)?; 8509df49d7eSJed Brown 8519df49d7eSJed Brown // Bases 8529df49d7eSJed Brown let bx = ceed.basis_tensor_H1_Lagrange(1, 1, 2, q, QuadMode::Gauss)?; 8539df49d7eSJed Brown let bu = ceed.basis_tensor_H1_Lagrange(1, 1, p, q, QuadMode::Gauss)?; 8549df49d7eSJed Brown 8559df49d7eSJed Brown // Build quadrature data 8569df49d7eSJed Brown let qf_build = ceed.q_function_interior_by_name("Mass1DBuild")?; 8579df49d7eSJed Brown ceed.operator(&qf_build, QFunctionOpt::None, QFunctionOpt::None)? 8589df49d7eSJed Brown .field("dx", &rx, &bx, VectorOpt::Active)? 8599df49d7eSJed Brown .field("weights", ElemRestrictionOpt::None, &bx, VectorOpt::None)? 8609df49d7eSJed Brown .field("qdata", &rq, BasisOpt::Collocated, VectorOpt::Active)? 8619df49d7eSJed Brown .apply(&x, &mut qdata)?; 8629df49d7eSJed Brown 8639df49d7eSJed Brown // Mass operator 8649df49d7eSJed Brown let qf_mass = ceed.q_function_interior_by_name("MassApply")?; 8659df49d7eSJed Brown let op_mass = ceed 8669df49d7eSJed Brown .operator(&qf_mass, QFunctionOpt::None, QFunctionOpt::None)? 8679df49d7eSJed Brown .field("u", &ru, &bu, VectorOpt::Active)? 8689df49d7eSJed Brown .field("qdata", &rq, BasisOpt::Collocated, &qdata)? 8696f97ff0aSJeremy L Thompson .field("v", &ru, &bu, VectorOpt::Active)? 8706f97ff0aSJeremy L Thompson .check()?; 8719df49d7eSJed Brown 8729df49d7eSJed Brown v.set_value(0.0)?; 8739df49d7eSJed Brown op_mass.apply(&u, &mut v)?; 8749df49d7eSJed Brown 8759df49d7eSJed Brown // Check 876e78171edSJeremy L Thompson let sum: Scalar = v.view()?.iter().sum(); 8779df49d7eSJed Brown assert!( 8789df49d7eSJed Brown (sum - 2.0).abs() < 1e-15, 8799df49d7eSJed Brown "Incorrect interval length computed" 8809df49d7eSJed Brown ); 88189d15d5fSJeremy L Thompson Ok(()) 8829df49d7eSJed Brown } 8839df49d7eSJed Brown 8849df49d7eSJed Brown #[test] 8859df49d7eSJed Brown fn test_ceed_t501() { 8869df49d7eSJed Brown assert!(ceed_t501().is_ok()); 8879df49d7eSJed Brown } 8889df49d7eSJed Brown } 8899df49d7eSJed Brown 8909df49d7eSJed Brown // ----------------------------------------------------------------------------- 891