1*d275d636SJeremy L Thompson // Copyright (c) 2017-2025, Lawrence Livermore National Security, LLC and other CEED contributors. 23d8e8822SJeremy L Thompson // All Rights Reserved. See the top-level LICENSE and NOTICE files for details. 39df49d7eSJed Brown // 43d8e8822SJeremy L Thompson // SPDX-License-Identifier: BSD-2-Clause 59df49d7eSJed Brown // 63d8e8822SJeremy L Thompson // This file is part of CEED: http://github.com/ceed 79df49d7eSJed Brown 8a1cbad85SJed Brown // Fenced `rust` code blocks included from README.md are executed as part of doctests. 9a1cbad85SJed Brown #![doc = include_str!("../README.md")] 109df49d7eSJed Brown // ----------------------------------------------------------------------------- 119df49d7eSJed Brown // Exceptions 129df49d7eSJed Brown // ----------------------------------------------------------------------------- 139df49d7eSJed Brown #![allow(non_snake_case)] 149df49d7eSJed Brown 159df49d7eSJed Brown // ----------------------------------------------------------------------------- 169df49d7eSJed Brown // Crate prelude 179df49d7eSJed Brown // ----------------------------------------------------------------------------- 189df49d7eSJed Brown use crate::prelude::*; 199df49d7eSJed Brown use std::sync::Once; 209df49d7eSJed Brown 219df49d7eSJed Brown pub mod prelude { 229df49d7eSJed Brown pub use crate::{ 239df49d7eSJed Brown basis::{self, Basis, BasisOpt}, 249df49d7eSJed Brown elem_restriction::{self, ElemRestriction, ElemRestrictionOpt}, 2508778c6fSJeremy L Thompson operator::{self, CompositeOperator, Operator, OperatorField}, 269df49d7eSJed Brown qfunction::{ 2708778c6fSJeremy L Thompson self, QFunction, QFunctionByName, QFunctionField, QFunctionInputs, QFunctionOpt, 2808778c6fSJeremy L Thompson QFunctionOutputs, 299df49d7eSJed Brown }, 30486868d3SJeremy L Thompson vector::{self, Vector, VectorOpt, VectorSliceWrapper}, 312ba8e59cSJeremy L Thompson ElemTopology, EvalMode, MemType, NormType, QuadMode, Scalar, TransposeMode, 3280a9ef05SNatalie Beams CEED_STRIDES_BACKEND, EPSILON, MAX_QFUNCTION_FIELDS, 339df49d7eSJed Brown }; 34630ad4c9Sjeremylt pub(crate) use libceed_sys::bind_ceed; 359df49d7eSJed Brown pub(crate) use std::convert::TryFrom; 3608778c6fSJeremy L Thompson pub(crate) use std::ffi::{CStr, CString}; 379df49d7eSJed Brown pub(crate) use std::fmt; 381142270cSJeremy L Thompson pub(crate) use std::marker::PhantomData; 399df49d7eSJed Brown } 409df49d7eSJed Brown 419df49d7eSJed Brown // ----------------------------------------------------------------------------- 429df49d7eSJed Brown // Modules 439df49d7eSJed Brown // ----------------------------------------------------------------------------- 449df49d7eSJed Brown pub mod basis; 459df49d7eSJed Brown pub mod elem_restriction; 469df49d7eSJed Brown pub mod operator; 479df49d7eSJed Brown pub mod qfunction; 489df49d7eSJed Brown pub mod vector; 499df49d7eSJed Brown 509df49d7eSJed Brown // ----------------------------------------------------------------------------- 5180a9ef05SNatalie Beams // Typedef for scalar 5280a9ef05SNatalie Beams // ----------------------------------------------------------------------------- 5380a9ef05SNatalie Beams pub type Scalar = bind_ceed::CeedScalar; 5480a9ef05SNatalie Beams 5580a9ef05SNatalie Beams // ----------------------------------------------------------------------------- 569df49d7eSJed Brown // Constants for library 579df49d7eSJed Brown // ----------------------------------------------------------------------------- 58257916f8SJed Brown const MAX_BUFFER_LENGTH: usize = 4096; 599df49d7eSJed Brown pub const MAX_QFUNCTION_FIELDS: usize = 16; 609df49d7eSJed Brown pub const CEED_STRIDES_BACKEND: [i32; 3] = [0; 3]; 6180a9ef05SNatalie Beams pub const EPSILON: crate::Scalar = bind_ceed::CEED_EPSILON as crate::Scalar; 629df49d7eSJed Brown 639df49d7eSJed Brown // ----------------------------------------------------------------------------- 649df49d7eSJed Brown // Enums for libCEED 659df49d7eSJed Brown // ----------------------------------------------------------------------------- 669df49d7eSJed Brown #[derive(Clone, Copy, PartialEq, Eq)] 679df49d7eSJed Brown /// Many Ceed interfaces take or return pointers to memory. This enum is used to 689df49d7eSJed Brown /// specify where the memory being provided or requested must reside. 699df49d7eSJed Brown pub enum MemType { 709df49d7eSJed Brown Host = bind_ceed::CeedMemType_CEED_MEM_HOST as isize, 719df49d7eSJed Brown Device = bind_ceed::CeedMemType_CEED_MEM_DEVICE as isize, 729df49d7eSJed Brown } 739df49d7eSJed Brown 749df49d7eSJed Brown #[derive(Clone, Copy, PartialEq, Eq)] 759df49d7eSJed Brown // OwnPointer will not be used by user but is included for internal use 769df49d7eSJed Brown #[allow(dead_code)] 779df49d7eSJed Brown /// Conveys ownership status of arrays passed to Ceed interfaces. 789df49d7eSJed Brown pub(crate) enum CopyMode { 799df49d7eSJed Brown CopyValues = bind_ceed::CeedCopyMode_CEED_COPY_VALUES as isize, 809df49d7eSJed Brown UsePointer = bind_ceed::CeedCopyMode_CEED_USE_POINTER as isize, 819df49d7eSJed Brown OwnPointer = bind_ceed::CeedCopyMode_CEED_OWN_POINTER as isize, 829df49d7eSJed Brown } 839df49d7eSJed Brown 849df49d7eSJed Brown #[derive(Clone, Copy, PartialEq, Eq)] 859df49d7eSJed Brown /// Denotes type of vector norm to be computed 869df49d7eSJed Brown pub enum NormType { 879df49d7eSJed Brown One = bind_ceed::CeedNormType_CEED_NORM_1 as isize, 889df49d7eSJed Brown Two = bind_ceed::CeedNormType_CEED_NORM_2 as isize, 899df49d7eSJed Brown Max = bind_ceed::CeedNormType_CEED_NORM_MAX as isize, 909df49d7eSJed Brown } 919df49d7eSJed Brown 929df49d7eSJed Brown #[derive(Clone, Copy, PartialEq, Eq)] 939df49d7eSJed Brown /// Denotes whether a linear transformation or its transpose should be applied 949df49d7eSJed Brown pub enum TransposeMode { 959df49d7eSJed Brown NoTranspose = bind_ceed::CeedTransposeMode_CEED_NOTRANSPOSE as isize, 969df49d7eSJed Brown Transpose = bind_ceed::CeedTransposeMode_CEED_TRANSPOSE as isize, 979df49d7eSJed Brown } 989df49d7eSJed Brown 999df49d7eSJed Brown #[derive(Clone, Copy, PartialEq, Eq)] 1009df49d7eSJed Brown /// Type of quadrature; also used for location of nodes 1019df49d7eSJed Brown pub enum QuadMode { 1029df49d7eSJed Brown Gauss = bind_ceed::CeedQuadMode_CEED_GAUSS as isize, 1039df49d7eSJed Brown GaussLobatto = bind_ceed::CeedQuadMode_CEED_GAUSS_LOBATTO as isize, 1049df49d7eSJed Brown } 1059df49d7eSJed Brown 1069df49d7eSJed Brown #[derive(Clone, Copy, PartialEq, Eq)] 1079df49d7eSJed Brown /// Type of basis shape to create non-tensor H1 element basis 1089df49d7eSJed Brown pub enum ElemTopology { 10950c301a5SRezgar Shakeri Line = bind_ceed::CeedElemTopology_CEED_TOPOLOGY_LINE as isize, 11050c301a5SRezgar Shakeri Triangle = bind_ceed::CeedElemTopology_CEED_TOPOLOGY_TRIANGLE as isize, 11150c301a5SRezgar Shakeri Quad = bind_ceed::CeedElemTopology_CEED_TOPOLOGY_QUAD as isize, 11250c301a5SRezgar Shakeri Tet = bind_ceed::CeedElemTopology_CEED_TOPOLOGY_TET as isize, 11350c301a5SRezgar Shakeri Pyramid = bind_ceed::CeedElemTopology_CEED_TOPOLOGY_PYRAMID as isize, 11450c301a5SRezgar Shakeri Prism = bind_ceed::CeedElemTopology_CEED_TOPOLOGY_PRISM as isize, 11550c301a5SRezgar Shakeri Hex = bind_ceed::CeedElemTopology_CEED_TOPOLOGY_HEX as isize, 1169df49d7eSJed Brown } 1179df49d7eSJed Brown 118c68be7a2SJeremy L Thompson #[derive(Clone, Copy, Debug, PartialEq, Eq)] 1199df49d7eSJed Brown /// Basis evaluation mode 1209df49d7eSJed Brown pub enum EvalMode { 1219df49d7eSJed Brown None = bind_ceed::CeedEvalMode_CEED_EVAL_NONE as isize, 1229df49d7eSJed Brown Interp = bind_ceed::CeedEvalMode_CEED_EVAL_INTERP as isize, 1239df49d7eSJed Brown Grad = bind_ceed::CeedEvalMode_CEED_EVAL_GRAD as isize, 1249df49d7eSJed Brown Div = bind_ceed::CeedEvalMode_CEED_EVAL_DIV as isize, 1259df49d7eSJed Brown Curl = bind_ceed::CeedEvalMode_CEED_EVAL_CURL as isize, 1269df49d7eSJed Brown Weight = bind_ceed::CeedEvalMode_CEED_EVAL_WEIGHT as isize, 1279df49d7eSJed Brown } 128c68be7a2SJeremy L Thompson impl EvalMode { 129c68be7a2SJeremy L Thompson pub(crate) fn from_u32(value: u32) -> EvalMode { 130c68be7a2SJeremy L Thompson match value { 131c68be7a2SJeremy L Thompson bind_ceed::CeedEvalMode_CEED_EVAL_NONE => EvalMode::None, 132c68be7a2SJeremy L Thompson bind_ceed::CeedEvalMode_CEED_EVAL_INTERP => EvalMode::Interp, 133c68be7a2SJeremy L Thompson bind_ceed::CeedEvalMode_CEED_EVAL_GRAD => EvalMode::Grad, 134c68be7a2SJeremy L Thompson bind_ceed::CeedEvalMode_CEED_EVAL_DIV => EvalMode::Div, 135c68be7a2SJeremy L Thompson bind_ceed::CeedEvalMode_CEED_EVAL_CURL => EvalMode::Curl, 136c68be7a2SJeremy L Thompson bind_ceed::CeedEvalMode_CEED_EVAL_WEIGHT => EvalMode::Weight, 137c68be7a2SJeremy L Thompson _ => panic!("Unknown value: {}", value), 138c68be7a2SJeremy L Thompson } 139c68be7a2SJeremy L Thompson } 140c68be7a2SJeremy L Thompson } 1419df49d7eSJed Brown 1429df49d7eSJed Brown // ----------------------------------------------------------------------------- 1439df49d7eSJed Brown // Ceed error 1449df49d7eSJed Brown // ----------------------------------------------------------------------------- 1452ba8e59cSJeremy L Thompson pub type Result<T> = std::result::Result<T, Error>; 1469df49d7eSJed Brown 1477ed177dbSJed Brown /// libCEED error messages - returning an Error without painc!ing indicates 1487ed177dbSJed Brown /// the function call failed but the data is not corrupted 1499df49d7eSJed Brown #[derive(Debug)] 1502ba8e59cSJeremy L Thompson pub struct Error { 1519df49d7eSJed Brown pub message: String, 1529df49d7eSJed Brown } 1539df49d7eSJed Brown 1542ba8e59cSJeremy L Thompson impl fmt::Display for Error { 1559df49d7eSJed Brown fn fmt(&self, f: &mut fmt::Formatter) -> fmt::Result { 1569df49d7eSJed Brown write!(f, "{}", self.message) 1579df49d7eSJed Brown } 1589df49d7eSJed Brown } 1599df49d7eSJed Brown 1609df49d7eSJed Brown // ----------------------------------------------------------------------------- 1611142270cSJeremy L Thompson // Internal error checker 1621142270cSJeremy L Thompson // ----------------------------------------------------------------------------- 1631142270cSJeremy L Thompson #[doc(hidden)] 16411544396SJeremy L Thompson pub(crate) fn check_error<F>(ceed_ptr: F, ierr: i32) -> Result<i32> 16511544396SJeremy L Thompson where 16611544396SJeremy L Thompson F: FnOnce() -> bind_ceed::Ceed, 16711544396SJeremy L Thompson { 1681142270cSJeremy L Thompson // Return early if code is clean 1691142270cSJeremy L Thompson if ierr == bind_ceed::CeedErrorType_CEED_ERROR_SUCCESS { 1701142270cSJeremy L Thompson return Ok(ierr); 1711142270cSJeremy L Thompson } 1721142270cSJeremy L Thompson // Retrieve error message 1731142270cSJeremy L Thompson let mut ptr: *const std::os::raw::c_char = std::ptr::null_mut(); 1741142270cSJeremy L Thompson let c_str = unsafe { 17511544396SJeremy L Thompson bind_ceed::CeedGetErrorMessage(ceed_ptr(), &mut ptr); 1761142270cSJeremy L Thompson std::ffi::CStr::from_ptr(ptr) 1771142270cSJeremy L Thompson }; 1781142270cSJeremy L Thompson let message = c_str.to_string_lossy().to_string(); 1791142270cSJeremy L Thompson // Panic if negative code, otherwise return error 1801142270cSJeremy L Thompson if ierr < bind_ceed::CeedErrorType_CEED_ERROR_SUCCESS { 1811142270cSJeremy L Thompson panic!("{}", message); 1821142270cSJeremy L Thompson } 1832ba8e59cSJeremy L Thompson Err(Error { message }) 1841142270cSJeremy L Thompson } 1851142270cSJeremy L Thompson 1861142270cSJeremy L Thompson // ----------------------------------------------------------------------------- 1879df49d7eSJed Brown // Ceed error handler 1889df49d7eSJed Brown // ----------------------------------------------------------------------------- 1892ba8e59cSJeremy L Thompson pub enum ErrorHandler { 1909df49d7eSJed Brown ErrorAbort, 1919df49d7eSJed Brown ErrorExit, 1929df49d7eSJed Brown ErrorReturn, 1939df49d7eSJed Brown ErrorStore, 1949df49d7eSJed Brown } 1959df49d7eSJed Brown 1969df49d7eSJed Brown // ----------------------------------------------------------------------------- 1979df49d7eSJed Brown // Ceed context wrapper 1989df49d7eSJed Brown // ----------------------------------------------------------------------------- 1999df49d7eSJed Brown /// A Ceed is a library context representing control of a logical hardware 2009df49d7eSJed Brown /// resource. 2019df49d7eSJed Brown #[derive(Debug)] 2029df49d7eSJed Brown pub struct Ceed { 2039df49d7eSJed Brown ptr: bind_ceed::Ceed, 2049df49d7eSJed Brown } 2059df49d7eSJed Brown 2069df49d7eSJed Brown // ----------------------------------------------------------------------------- 2079df49d7eSJed Brown // Destructor 2089df49d7eSJed Brown // ----------------------------------------------------------------------------- 2099df49d7eSJed Brown impl Drop for Ceed { 2109df49d7eSJed Brown fn drop(&mut self) { 2119df49d7eSJed Brown unsafe { 2129df49d7eSJed Brown bind_ceed::CeedDestroy(&mut self.ptr); 2139df49d7eSJed Brown } 2149df49d7eSJed Brown } 2159df49d7eSJed Brown } 2169df49d7eSJed Brown 2179df49d7eSJed Brown // ----------------------------------------------------------------------------- 21859189cfaSJeremy L Thompson // Cloning 21959189cfaSJeremy L Thompson // ----------------------------------------------------------------------------- 22059189cfaSJeremy L Thompson impl Clone for Ceed { 22159189cfaSJeremy L Thompson /// Perform a shallow clone of a Ceed context 22259189cfaSJeremy L Thompson /// 22359189cfaSJeremy L Thompson /// ``` 22459189cfaSJeremy L Thompson /// let ceed = libceed::Ceed::init("/cpu/self/ref/serial"); 22559189cfaSJeremy L Thompson /// let ceed_clone = ceed.clone(); 22659189cfaSJeremy L Thompson /// 2277ed177dbSJed Brown /// println!(" original:{} \n clone:{}", ceed, ceed_clone); 22859189cfaSJeremy L Thompson /// ``` 22959189cfaSJeremy L Thompson fn clone(&self) -> Self { 23059189cfaSJeremy L Thompson let mut ptr_clone = std::ptr::null_mut(); 231656ef1e5SJeremy L Thompson self.check_error(unsafe { bind_ceed::CeedReferenceCopy(self.ptr, &mut ptr_clone) }) 232656ef1e5SJeremy L Thompson .expect("failed to clone Ceed"); 23359189cfaSJeremy L Thompson Self { ptr: ptr_clone } 23459189cfaSJeremy L Thompson } 23559189cfaSJeremy L Thompson } 23659189cfaSJeremy L Thompson 23759189cfaSJeremy L Thompson // ----------------------------------------------------------------------------- 2389df49d7eSJed Brown // Display 2399df49d7eSJed Brown // ----------------------------------------------------------------------------- 2409df49d7eSJed Brown impl fmt::Display for Ceed { 2419df49d7eSJed Brown /// View a Ceed 2429df49d7eSJed Brown /// 2439df49d7eSJed Brown /// ``` 2449df49d7eSJed Brown /// let ceed = libceed::Ceed::init("/cpu/self/ref/serial"); 2459df49d7eSJed Brown /// println!("{}", ceed); 2469df49d7eSJed Brown /// ``` 2479df49d7eSJed Brown fn fmt(&self, f: &mut fmt::Formatter) -> fmt::Result { 2489df49d7eSJed Brown let mut ptr = std::ptr::null_mut(); 2499df49d7eSJed Brown let mut sizeloc = crate::MAX_BUFFER_LENGTH; 2509df49d7eSJed Brown let cstring = unsafe { 2519df49d7eSJed Brown let file = bind_ceed::open_memstream(&mut ptr, &mut sizeloc); 2529df49d7eSJed Brown bind_ceed::CeedView(self.ptr, file); 2539df49d7eSJed Brown bind_ceed::fclose(file); 2549df49d7eSJed Brown CString::from_raw(ptr) 2559df49d7eSJed Brown }; 2569df49d7eSJed Brown cstring.to_string_lossy().fmt(f) 2579df49d7eSJed Brown } 2589df49d7eSJed Brown } 2599df49d7eSJed Brown 2609df49d7eSJed Brown static REGISTER: Once = Once::new(); 2619df49d7eSJed Brown 2629df49d7eSJed Brown // ----------------------------------------------------------------------------- 2639df49d7eSJed Brown // Object constructors 2649df49d7eSJed Brown // ----------------------------------------------------------------------------- 2659df49d7eSJed Brown impl Ceed { 2667ed177dbSJed Brown #[cfg_attr(feature = "katexit", katexit::katexit)] 2679df49d7eSJed Brown /// Returns a Ceed context initialized with the specified resource 2689df49d7eSJed Brown /// 2699df49d7eSJed Brown /// # arguments 2709df49d7eSJed Brown /// 2719df49d7eSJed Brown /// * `resource` - Resource to use, e.g., "/cpu/self" 2729df49d7eSJed Brown /// 2739df49d7eSJed Brown /// ``` 2749df49d7eSJed Brown /// let ceed = libceed::Ceed::init("/cpu/self/ref/serial"); 2759df49d7eSJed Brown /// ``` 2769df49d7eSJed Brown pub fn init(resource: &str) -> Self { 2772ba8e59cSJeremy L Thompson Ceed::init_with_error_handler(resource, ErrorHandler::ErrorStore) 2789df49d7eSJed Brown } 2799df49d7eSJed Brown 2809df49d7eSJed Brown /// Returns a Ceed context initialized with the specified resource 2819df49d7eSJed Brown /// 2829df49d7eSJed Brown /// # arguments 2839df49d7eSJed Brown /// 2849df49d7eSJed Brown /// * `resource` - Resource to use, e.g., "/cpu/self" 2859df49d7eSJed Brown /// 2869df49d7eSJed Brown /// ``` 2879df49d7eSJed Brown /// let ceed = libceed::Ceed::init_with_error_handler( 2889df49d7eSJed Brown /// "/cpu/self/ref/serial", 2892ba8e59cSJeremy L Thompson /// libceed::ErrorHandler::ErrorAbort, 2909df49d7eSJed Brown /// ); 2919df49d7eSJed Brown /// ``` 2922ba8e59cSJeremy L Thompson pub fn init_with_error_handler(resource: &str, handler: ErrorHandler) -> Self { 2939df49d7eSJed Brown REGISTER.call_once(|| unsafe { 2949df49d7eSJed Brown bind_ceed::CeedRegisterAll(); 2959df49d7eSJed Brown bind_ceed::CeedQFunctionRegisterAll(); 2969df49d7eSJed Brown }); 2979df49d7eSJed Brown 2989df49d7eSJed Brown // Convert to C string 2999df49d7eSJed Brown let c_resource = CString::new(resource).expect("CString::new failed"); 3009df49d7eSJed Brown 3019df49d7eSJed Brown // Get error handler pointer 3029df49d7eSJed Brown let eh = match handler { 3032ba8e59cSJeremy L Thompson ErrorHandler::ErrorAbort => bind_ceed::CeedErrorAbort, 3042ba8e59cSJeremy L Thompson ErrorHandler::ErrorExit => bind_ceed::CeedErrorExit, 3052ba8e59cSJeremy L Thompson ErrorHandler::ErrorReturn => bind_ceed::CeedErrorReturn, 3062ba8e59cSJeremy L Thompson ErrorHandler::ErrorStore => bind_ceed::CeedErrorStore, 3079df49d7eSJed Brown }; 3089df49d7eSJed Brown 3099df49d7eSJed Brown // Call to libCEED 3109df49d7eSJed Brown let mut ptr = std::ptr::null_mut(); 3119df49d7eSJed Brown let mut ierr = unsafe { bind_ceed::CeedInit(c_resource.as_ptr() as *const i8, &mut ptr) }; 3129df49d7eSJed Brown if ierr != 0 { 3139df49d7eSJed Brown panic!("Error initializing backend resource: {}", resource) 3149df49d7eSJed Brown } 3159df49d7eSJed Brown ierr = unsafe { bind_ceed::CeedSetErrorHandler(ptr, Some(eh)) }; 3169df49d7eSJed Brown let ceed = Ceed { ptr }; 3179df49d7eSJed Brown ceed.check_error(ierr).unwrap(); 3189df49d7eSJed Brown ceed 3199df49d7eSJed Brown } 3209df49d7eSJed Brown 3219df49d7eSJed Brown /// Default initializer for testing 3229df49d7eSJed Brown #[doc(hidden)] 3239df49d7eSJed Brown pub fn default_init() -> Self { 3249df49d7eSJed Brown // Convert to C string 3259df49d7eSJed Brown let resource = "/cpu/self/ref/serial"; 3269df49d7eSJed Brown crate::Ceed::init(resource) 3279df49d7eSJed Brown } 3289df49d7eSJed Brown 3299df49d7eSJed Brown /// Internal error checker 3309df49d7eSJed Brown #[doc(hidden)] 3319df49d7eSJed Brown fn check_error(&self, ierr: i32) -> Result<i32> { 3329df49d7eSJed Brown // Return early if code is clean 3339df49d7eSJed Brown if ierr == bind_ceed::CeedErrorType_CEED_ERROR_SUCCESS { 3349df49d7eSJed Brown return Ok(ierr); 3359df49d7eSJed Brown } 3369df49d7eSJed Brown // Retrieve error message 3379df49d7eSJed Brown let mut ptr: *const std::os::raw::c_char = std::ptr::null_mut(); 3389df49d7eSJed Brown let c_str = unsafe { 3399df49d7eSJed Brown bind_ceed::CeedGetErrorMessage(self.ptr, &mut ptr); 3409df49d7eSJed Brown std::ffi::CStr::from_ptr(ptr) 3419df49d7eSJed Brown }; 3429df49d7eSJed Brown let message = c_str.to_string_lossy().to_string(); 3439df49d7eSJed Brown // Panic if negative code, otherwise return error 3449df49d7eSJed Brown if ierr < bind_ceed::CeedErrorType_CEED_ERROR_SUCCESS { 3459df49d7eSJed Brown panic!("{}", message); 3469df49d7eSJed Brown } 3472ba8e59cSJeremy L Thompson Err(Error { message }) 3489df49d7eSJed Brown } 3499df49d7eSJed Brown 3509df49d7eSJed Brown /// Returns full resource name for a Ceed context 3519df49d7eSJed Brown /// 3529df49d7eSJed Brown /// ``` 3539df49d7eSJed Brown /// let ceed = libceed::Ceed::init("/cpu/self/ref/serial"); 3549df49d7eSJed Brown /// let resource = ceed.resource(); 3559df49d7eSJed Brown /// 3569df49d7eSJed Brown /// assert_eq!(resource, "/cpu/self/ref/serial".to_string()) 3579df49d7eSJed Brown /// ``` 3589df49d7eSJed Brown pub fn resource(&self) -> String { 3599df49d7eSJed Brown let mut ptr: *const std::os::raw::c_char = std::ptr::null_mut(); 3609df49d7eSJed Brown let c_str = unsafe { 3619df49d7eSJed Brown bind_ceed::CeedGetResource(self.ptr, &mut ptr); 3629df49d7eSJed Brown std::ffi::CStr::from_ptr(ptr) 3639df49d7eSJed Brown }; 3649df49d7eSJed Brown c_str.to_string_lossy().to_string() 3659df49d7eSJed Brown } 3669df49d7eSJed Brown 3677ed177dbSJed Brown /// Returns a Vector of the specified length (does not allocate memory) 3689df49d7eSJed Brown /// 3699df49d7eSJed Brown /// # arguments 3709df49d7eSJed Brown /// 3719df49d7eSJed Brown /// * `n` - Length of vector 3729df49d7eSJed Brown /// 3739df49d7eSJed Brown /// ``` 3749df49d7eSJed Brown /// # use libceed::prelude::*; 3754d27c890SJeremy L Thompson /// # fn main() -> libceed::Result<()> { 3769df49d7eSJed Brown /// # let ceed = libceed::Ceed::default_init(); 377c68be7a2SJeremy L Thompson /// let vec = ceed.vector(10)?; 378c68be7a2SJeremy L Thompson /// # Ok(()) 379c68be7a2SJeremy L Thompson /// # } 3809df49d7eSJed Brown /// ``` 381594ef120SJeremy L Thompson pub fn vector<'a>(&self, n: usize) -> Result<Vector<'a>> { 3829df49d7eSJed Brown Vector::create(self, n) 3839df49d7eSJed Brown } 3849df49d7eSJed Brown 3859df49d7eSJed Brown /// Create a Vector initialized with the data (copied) from a slice 3869df49d7eSJed Brown /// 3879df49d7eSJed Brown /// # arguments 3889df49d7eSJed Brown /// 3899df49d7eSJed Brown /// * `slice` - Slice containing data 3909df49d7eSJed Brown /// 3919df49d7eSJed Brown /// ``` 3929df49d7eSJed Brown /// # use libceed::prelude::*; 3934d27c890SJeremy L Thompson /// # fn main() -> libceed::Result<()> { 3949df49d7eSJed Brown /// # let ceed = libceed::Ceed::default_init(); 395c68be7a2SJeremy L Thompson /// let vec = ceed.vector_from_slice(&[1., 2., 3.])?; 3969df49d7eSJed Brown /// assert_eq!(vec.length(), 3); 397c68be7a2SJeremy L Thompson /// # Ok(()) 398c68be7a2SJeremy L Thompson /// # } 3999df49d7eSJed Brown /// ``` 400594ef120SJeremy L Thompson pub fn vector_from_slice<'a>(&self, slice: &[crate::Scalar]) -> Result<Vector<'a>> { 4019df49d7eSJed Brown Vector::from_slice(self, slice) 4029df49d7eSJed Brown } 4039df49d7eSJed Brown 4047ed177dbSJed Brown /// Returns an ElemRestriction, $\mathcal{E}$, which extracts the degrees of 4057ed177dbSJed Brown /// freedom for each element from the local vector into the element vector 4067ed177dbSJed Brown /// or assembles contributions from each element in the element vector to 4077ed177dbSJed Brown /// the local vector 4089df49d7eSJed Brown /// 4099df49d7eSJed Brown /// # arguments 4109df49d7eSJed Brown /// 4119df49d7eSJed Brown /// * `nelem` - Number of elements described in the offsets array 4129df49d7eSJed Brown /// * `elemsize` - Size (number of "nodes") per element 4139df49d7eSJed Brown /// * `ncomp` - Number of field components per interpolation node (1 4149df49d7eSJed Brown /// for scalar fields) 4159df49d7eSJed Brown /// * `compstride` - Stride between components for the same Lvector "node". 4169df49d7eSJed Brown /// Data for node `i`, component `j`, element `k` can be 4179df49d7eSJed Brown /// found in the Lvector at index 4189df49d7eSJed Brown /// `offsets[i + k*elemsize] + j*compstride`. 4199df49d7eSJed Brown /// * `lsize` - The size of the Lvector. This vector may be larger 4209df49d7eSJed Brown /// than the elements and fields given by this 4219df49d7eSJed Brown /// restriction. 4229df49d7eSJed Brown /// * `mtype` - Memory type of the offsets array, see CeedMemType 4239df49d7eSJed Brown /// * `offsets` - Array of shape `[nelem, elemsize]`. Row `i` holds the 4249df49d7eSJed Brown /// ordered list of the offsets (into the input CeedVector) 4259df49d7eSJed Brown /// for the unknowns corresponding to element `i`, where 4269df49d7eSJed Brown /// `0 <= i < nelem`. All offsets must be in the range 4279df49d7eSJed Brown /// `[0, lsize - 1]`. 4289df49d7eSJed Brown /// 4299df49d7eSJed Brown /// ``` 4309df49d7eSJed Brown /// # use libceed::prelude::*; 4314d27c890SJeremy L Thompson /// # fn main() -> libceed::Result<()> { 4329df49d7eSJed Brown /// # let ceed = libceed::Ceed::default_init(); 4339df49d7eSJed Brown /// let nelem = 3; 4349df49d7eSJed Brown /// let mut ind: Vec<i32> = vec![0; 2 * nelem]; 4359df49d7eSJed Brown /// for i in 0..nelem { 4369df49d7eSJed Brown /// ind[2 * i + 0] = i as i32; 4379df49d7eSJed Brown /// ind[2 * i + 1] = (i + 1) as i32; 4389df49d7eSJed Brown /// } 439c68be7a2SJeremy L Thompson /// let r = ceed.elem_restriction(nelem, 2, 1, 1, nelem + 1, MemType::Host, &ind)?; 440c68be7a2SJeremy L Thompson /// # Ok(()) 441c68be7a2SJeremy L Thompson /// # } 4429df49d7eSJed Brown /// ``` 443594ef120SJeremy L Thompson pub fn elem_restriction<'a>( 4449df49d7eSJed Brown &self, 4459df49d7eSJed Brown nelem: usize, 4469df49d7eSJed Brown elemsize: usize, 4479df49d7eSJed Brown ncomp: usize, 4489df49d7eSJed Brown compstride: usize, 4499df49d7eSJed Brown lsize: usize, 4509df49d7eSJed Brown mtype: MemType, 4519df49d7eSJed Brown offsets: &[i32], 452594ef120SJeremy L Thompson ) -> Result<ElemRestriction<'a>> { 4539df49d7eSJed Brown ElemRestriction::create( 4549df49d7eSJed Brown self, nelem, elemsize, ncomp, compstride, lsize, mtype, offsets, 455709403c1SSebastian Grimberg ) 456709403c1SSebastian Grimberg } 457709403c1SSebastian Grimberg 458709403c1SSebastian Grimberg /// Returns an oriented ElemRestriction, $\mathcal{E}$, which extracts the 459709403c1SSebastian Grimberg /// degrees of freedom for each element from the local vector into the 460709403c1SSebastian Grimberg /// element vector or assembles contributions from each element in the 461709403c1SSebastian Grimberg /// element vector to the local vector 462709403c1SSebastian Grimberg /// 463709403c1SSebastian Grimberg /// # arguments 464709403c1SSebastian Grimberg /// 465709403c1SSebastian Grimberg /// * `nelem` - Number of elements described in the offsets array 466709403c1SSebastian Grimberg /// * `elemsize` - Size (number of "nodes") per element 467709403c1SSebastian Grimberg /// * `ncomp` - Number of field components per interpolation node (1 468709403c1SSebastian Grimberg /// for scalar fields) 469709403c1SSebastian Grimberg /// * `compstride` - Stride between components for the same Lvector "node". 470709403c1SSebastian Grimberg /// Data for node `i`, component `j`, element `k` can be 471709403c1SSebastian Grimberg /// found in the Lvector at index 472709403c1SSebastian Grimberg /// `offsets[i + k*elemsize] + j*compstride`. 473709403c1SSebastian Grimberg /// * `lsize` - The size of the Lvector. This vector may be larger 474709403c1SSebastian Grimberg /// than the elements and fields given by this 475709403c1SSebastian Grimberg /// restriction. 476709403c1SSebastian Grimberg /// * `mtype` - Memory type of the offsets array, see CeedMemType 477709403c1SSebastian Grimberg /// * `offsets` - Array of shape `[nelem, elemsize]`. Row `i` holds the 478709403c1SSebastian Grimberg /// ordered list of the offsets (into the input CeedVector) 479709403c1SSebastian Grimberg /// for the unknowns corresponding to element `i`, where 480709403c1SSebastian Grimberg /// `0 <= i < nelem`. All offsets must be in the range 481709403c1SSebastian Grimberg /// `[0, lsize - 1]`. 482709403c1SSebastian Grimberg /// * `orients` - Array of shape `[nelem, elemsize]`. Row `i` holds the 483709403c1SSebastian Grimberg /// ordered list of the orientations for the unknowns 484709403c1SSebastian Grimberg /// corresponding to element `i`, with bool `false` used 485709403c1SSebastian Grimberg /// for positively oriented and `true` to flip the 486709403c1SSebastian Grimberg /// orientation. 487709403c1SSebastian Grimberg /// 488709403c1SSebastian Grimberg /// ``` 489709403c1SSebastian Grimberg /// # use libceed::prelude::*; 490709403c1SSebastian Grimberg /// # fn main() -> libceed::Result<()> { 491709403c1SSebastian Grimberg /// # let ceed = libceed::Ceed::default_init(); 492709403c1SSebastian Grimberg /// let nelem = 3; 493709403c1SSebastian Grimberg /// let mut ind: Vec<i32> = vec![0; 2 * nelem]; 494709403c1SSebastian Grimberg /// let mut orients: Vec<bool> = vec![false; 2 * nelem]; 495709403c1SSebastian Grimberg /// for i in 0..nelem { 496709403c1SSebastian Grimberg /// ind[2 * i + 0] = i as i32; 497709403c1SSebastian Grimberg /// ind[2 * i + 1] = (i + 1) as i32; 498709403c1SSebastian Grimberg /// orients[2 * i + 0] = (i % 2) > 0; // flip the dofs on element 1, 3, ... 499709403c1SSebastian Grimberg /// orients[2 * i + 1] = (i % 2) > 0; 500709403c1SSebastian Grimberg /// } 501709403c1SSebastian Grimberg /// let r = 502709403c1SSebastian Grimberg /// ceed.oriented_elem_restriction(nelem, 2, 1, 1, nelem + 1, MemType::Host, &ind, &orients)?; 503709403c1SSebastian Grimberg /// # Ok(()) 504709403c1SSebastian Grimberg /// # } 505709403c1SSebastian Grimberg /// ``` 506709403c1SSebastian Grimberg pub fn oriented_elem_restriction<'a>( 507709403c1SSebastian Grimberg &self, 508709403c1SSebastian Grimberg nelem: usize, 509709403c1SSebastian Grimberg elemsize: usize, 510709403c1SSebastian Grimberg ncomp: usize, 511709403c1SSebastian Grimberg compstride: usize, 512709403c1SSebastian Grimberg lsize: usize, 513709403c1SSebastian Grimberg mtype: MemType, 514709403c1SSebastian Grimberg offsets: &[i32], 515709403c1SSebastian Grimberg orients: &[bool], 516709403c1SSebastian Grimberg ) -> Result<ElemRestriction<'a>> { 517709403c1SSebastian Grimberg ElemRestriction::create_oriented( 518709403c1SSebastian Grimberg self, nelem, elemsize, ncomp, compstride, lsize, mtype, offsets, orients, 519709403c1SSebastian Grimberg ) 520709403c1SSebastian Grimberg } 521709403c1SSebastian Grimberg 522709403c1SSebastian Grimberg /// Returns a curl-oriented ElemRestriction, $\mathcal{E}$, which extracts 523709403c1SSebastian Grimberg /// the degrees of freedom for each element from the local vector into 524709403c1SSebastian Grimberg /// the element vector or assembles contributions from each element in 525709403c1SSebastian Grimberg /// the element vector to the local vector 526709403c1SSebastian Grimberg /// 527709403c1SSebastian Grimberg /// # arguments 528709403c1SSebastian Grimberg /// 529709403c1SSebastian Grimberg /// * `nelem` - Number of elements described in the offsets array 530709403c1SSebastian Grimberg /// * `elemsize` - Size (number of "nodes") per element 531709403c1SSebastian Grimberg /// * `ncomp` - Number of field components per interpolation node (1 532709403c1SSebastian Grimberg /// for scalar fields) 533709403c1SSebastian Grimberg /// * `compstride` - Stride between components for the same Lvector "node". 534709403c1SSebastian Grimberg /// Data for node `i`, component `j`, element `k` can be 535709403c1SSebastian Grimberg /// found in the Lvector at index 536709403c1SSebastian Grimberg /// `offsets[i + k*elemsize] + j*compstride`. 537709403c1SSebastian Grimberg /// * `lsize` - The size of the Lvector. This vector may be larger 538709403c1SSebastian Grimberg /// than the elements and fields given by this 539709403c1SSebastian Grimberg /// restriction. 540709403c1SSebastian Grimberg /// * `mtype` - Memory type of the offsets array, see CeedMemType 541709403c1SSebastian Grimberg /// * `offsets` - Array of shape `[nelem, elemsize]`. Row `i` holds the 542709403c1SSebastian Grimberg /// ordered list of the offsets (into the input CeedVector) 543709403c1SSebastian Grimberg /// for the unknowns corresponding to element `i`, where 544709403c1SSebastian Grimberg /// `0 <= i < nelem`. All offsets must be in the range 545709403c1SSebastian Grimberg /// `[0, lsize - 1]`. 546709403c1SSebastian Grimberg /// * `curlorients` - Array of shape `[nelem, 3 * elemsize]`. Row `i` holds 547709403c1SSebastian Grimberg /// a row-major tridiagonal matrix (`curlorients[i, 0] = 548709403c1SSebastian Grimberg /// curlorients[i, 3 * elemsize - 1] = 0`, where 549709403c1SSebastian Grimberg /// `0 <= i < nelem`) which is applied to the element 550709403c1SSebastian Grimberg /// unknowns upon restriction. 551709403c1SSebastian Grimberg /// 552709403c1SSebastian Grimberg /// ``` 553709403c1SSebastian Grimberg /// # use libceed::prelude::*; 554709403c1SSebastian Grimberg /// # fn main() -> libceed::Result<()> { 555709403c1SSebastian Grimberg /// # let ceed = libceed::Ceed::default_init(); 556709403c1SSebastian Grimberg /// let nelem = 3; 557709403c1SSebastian Grimberg /// let mut ind: Vec<i32> = vec![0; 2 * nelem]; 558709403c1SSebastian Grimberg /// let mut curlorients: Vec<i8> = vec![0; 3 * 2 * nelem]; 559709403c1SSebastian Grimberg /// for i in 0..nelem { 560709403c1SSebastian Grimberg /// ind[2 * i + 0] = i as i32; 561709403c1SSebastian Grimberg /// ind[2 * i + 1] = (i + 1) as i32; 562709403c1SSebastian Grimberg /// curlorients[3 * 2 * i] = 0 as i8; 563709403c1SSebastian Grimberg /// curlorients[3 * 2 * (i + 1) - 1] = 0 as i8; 564709403c1SSebastian Grimberg /// if (i % 2 > 0) { 565709403c1SSebastian Grimberg /// // T = [0 -1] 566709403c1SSebastian Grimberg /// // [-1 0] 567709403c1SSebastian Grimberg /// curlorients[3 * 2 * i + 1] = 0 as i8; 568709403c1SSebastian Grimberg /// curlorients[3 * 2 * i + 2] = -1 as i8; 569709403c1SSebastian Grimberg /// curlorients[3 * 2 * i + 3] = -1 as i8; 570709403c1SSebastian Grimberg /// curlorients[3 * 2 * i + 4] = 0 as i8; 571709403c1SSebastian Grimberg /// } else { 572709403c1SSebastian Grimberg /// // T = I 573709403c1SSebastian Grimberg /// curlorients[3 * 2 * i + 1] = 1 as i8; 574709403c1SSebastian Grimberg /// curlorients[3 * 2 * i + 2] = 0 as i8; 575709403c1SSebastian Grimberg /// curlorients[3 * 2 * i + 3] = 0 as i8; 576709403c1SSebastian Grimberg /// curlorients[3 * 2 * i + 4] = 1 as i8; 577709403c1SSebastian Grimberg /// } 578709403c1SSebastian Grimberg /// } 579709403c1SSebastian Grimberg /// let r = ceed.curl_oriented_elem_restriction( 580709403c1SSebastian Grimberg /// nelem, 581709403c1SSebastian Grimberg /// 2, 582709403c1SSebastian Grimberg /// 1, 583709403c1SSebastian Grimberg /// 1, 584709403c1SSebastian Grimberg /// nelem + 1, 585709403c1SSebastian Grimberg /// MemType::Host, 586709403c1SSebastian Grimberg /// &ind, 587709403c1SSebastian Grimberg /// &curlorients, 588709403c1SSebastian Grimberg /// )?; 589709403c1SSebastian Grimberg /// # Ok(()) 590709403c1SSebastian Grimberg /// # } 591709403c1SSebastian Grimberg /// ``` 592709403c1SSebastian Grimberg pub fn curl_oriented_elem_restriction<'a>( 593709403c1SSebastian Grimberg &self, 594709403c1SSebastian Grimberg nelem: usize, 595709403c1SSebastian Grimberg elemsize: usize, 596709403c1SSebastian Grimberg ncomp: usize, 597709403c1SSebastian Grimberg compstride: usize, 598709403c1SSebastian Grimberg lsize: usize, 599709403c1SSebastian Grimberg mtype: MemType, 600709403c1SSebastian Grimberg offsets: &[i32], 601709403c1SSebastian Grimberg curlorients: &[i8], 602709403c1SSebastian Grimberg ) -> Result<ElemRestriction<'a>> { 603709403c1SSebastian Grimberg ElemRestriction::create_curl_oriented( 604709403c1SSebastian Grimberg self, 605709403c1SSebastian Grimberg nelem, 606709403c1SSebastian Grimberg elemsize, 607709403c1SSebastian Grimberg ncomp, 608709403c1SSebastian Grimberg compstride, 609709403c1SSebastian Grimberg lsize, 610709403c1SSebastian Grimberg mtype, 611709403c1SSebastian Grimberg offsets, 612709403c1SSebastian Grimberg curlorients, 6139df49d7eSJed Brown ) 6149df49d7eSJed Brown } 6159df49d7eSJed Brown 6167ed177dbSJed Brown /// Returns an ElemRestriction, $\mathcal{E}$, from an local vector to 6177ed177dbSJed Brown /// an element vector where data can be indexed from the `strides` array 6189df49d7eSJed Brown /// 6199df49d7eSJed Brown /// # arguments 6209df49d7eSJed Brown /// 6219df49d7eSJed Brown /// * `nelem` - Number of elements described in the offsets array 6229df49d7eSJed Brown /// * `elemsize` - Size (number of "nodes") per element 6239df49d7eSJed Brown /// * `ncomp` - Number of field components per interpolation node (1 6249df49d7eSJed Brown /// for scalar fields) 6259df49d7eSJed Brown /// * `lsize` - The size of the Lvector. This vector may be larger 6269df49d7eSJed Brown /// than the elements and fields given by this restriction. 6279df49d7eSJed Brown /// * `strides` - Array for strides between `[nodes, components, elements]`. 6289df49d7eSJed Brown /// Data for node `i`, component `j`, element `k` can be 6299df49d7eSJed Brown /// found in the Lvector at index 6309df49d7eSJed Brown /// `i*strides[0] + j*strides[1] + k*strides[2]`. 6319df49d7eSJed Brown /// CEED_STRIDES_BACKEND may be used with vectors created 6329df49d7eSJed Brown /// by a Ceed backend. 6339df49d7eSJed Brown /// 6349df49d7eSJed Brown /// ``` 6359df49d7eSJed Brown /// # use libceed::prelude::*; 6364d27c890SJeremy L Thompson /// # fn main() -> libceed::Result<()> { 6379df49d7eSJed Brown /// # let ceed = libceed::Ceed::default_init(); 6389df49d7eSJed Brown /// let nelem = 3; 6399df49d7eSJed Brown /// let strides: [i32; 3] = [1, 2, 2]; 640c68be7a2SJeremy L Thompson /// let r = ceed.strided_elem_restriction(nelem, 2, 1, nelem * 2, strides)?; 641c68be7a2SJeremy L Thompson /// # Ok(()) 642c68be7a2SJeremy L Thompson /// # } 6439df49d7eSJed Brown /// ``` 644594ef120SJeremy L Thompson pub fn strided_elem_restriction<'a>( 6459df49d7eSJed Brown &self, 6469df49d7eSJed Brown nelem: usize, 6479df49d7eSJed Brown elemsize: usize, 6489df49d7eSJed Brown ncomp: usize, 6499df49d7eSJed Brown lsize: usize, 6509df49d7eSJed Brown strides: [i32; 3], 651594ef120SJeremy L Thompson ) -> Result<ElemRestriction<'a>> { 6529df49d7eSJed Brown ElemRestriction::create_strided(self, nelem, elemsize, ncomp, lsize, strides) 6539df49d7eSJed Brown } 6549df49d7eSJed Brown 6557ed177dbSJed Brown /// Returns an $H^1$ tensor-product Basis 6569df49d7eSJed Brown /// 6579df49d7eSJed Brown /// # arguments 6589df49d7eSJed Brown /// 6599df49d7eSJed Brown /// * `dim` - Topological dimension of element 6609df49d7eSJed Brown /// * `ncomp` - Number of field components (1 for scalar fields) 6619df49d7eSJed Brown /// * `P1d` - Number of Gauss-Lobatto nodes in one dimension. The 6629df49d7eSJed Brown /// polynomial degree of the resulting `Q_k` element is 6639df49d7eSJed Brown /// `k=P-1`. 6649df49d7eSJed Brown /// * `Q1d` - Number of quadrature points in one dimension 6659df49d7eSJed Brown /// * `interp1d` - Row-major `(Q1d * P1d)` matrix expressing the values of 6669df49d7eSJed Brown /// nodal basis functions at quadrature points 6679df49d7eSJed Brown /// * `grad1d` - Row-major `(Q1d * P1d)` matrix expressing derivatives of 6689df49d7eSJed Brown /// nodal basis functions at quadrature points 6699df49d7eSJed Brown /// * `qref1d` - Array of length `Q1d` holding the locations of quadrature 6709df49d7eSJed Brown /// points on the 1D reference element `[-1, 1]` 6719df49d7eSJed Brown /// * `qweight1d` - Array of length `Q1d` holding the quadrature weights on 6729df49d7eSJed Brown /// the reference element 6739df49d7eSJed Brown /// 6749df49d7eSJed Brown /// ``` 6759df49d7eSJed Brown /// # use libceed::prelude::*; 6764d27c890SJeremy L Thompson /// # fn main() -> libceed::Result<()> { 6779df49d7eSJed Brown /// # let ceed = libceed::Ceed::default_init(); 6789df49d7eSJed Brown /// let interp1d = [ 0.62994317, 0.47255875, -0.14950343, 0.04700152, 6799df49d7eSJed Brown /// -0.07069480, 0.97297619, 0.13253993, -0.03482132, 6809df49d7eSJed Brown /// -0.03482132, 0.13253993, 0.97297619, -0.07069480, 6819df49d7eSJed Brown /// 0.04700152, -0.14950343, 0.47255875, 0.62994317]; 6829df49d7eSJed Brown /// let grad1d = [-2.34183742, 2.78794489, -0.63510411, 0.18899664, 6839df49d7eSJed Brown /// -0.51670214, -0.48795249, 1.33790510, -0.33325047, 6849df49d7eSJed Brown // 0.33325047, -1.33790510, 0.48795249, 0.51670214, 6859df49d7eSJed Brown /// -0.18899664, 0.63510411, -2.78794489, 2.34183742]; 6869df49d7eSJed Brown /// let qref1d = [-0.86113631, -0.33998104, 0.33998104, 0.86113631]; 6879df49d7eSJed Brown /// let qweight1d = [ 0.34785485, 0.65214515, 0.65214515, 0.34785485]; 6889df49d7eSJed Brown /// let b = ceed. 689c68be7a2SJeremy L Thompson /// basis_tensor_H1(2, 1, 4, 4, &interp1d, &grad1d, &qref1d, &qweight1d)?; 690c68be7a2SJeremy L Thompson /// # Ok(()) 691c68be7a2SJeremy L Thompson /// # } 6929df49d7eSJed Brown /// ``` 693594ef120SJeremy L Thompson pub fn basis_tensor_H1<'a>( 6949df49d7eSJed Brown &self, 6959df49d7eSJed Brown dim: usize, 6969df49d7eSJed Brown ncomp: usize, 6979df49d7eSJed Brown P1d: usize, 6989df49d7eSJed Brown Q1d: usize, 69980a9ef05SNatalie Beams interp1d: &[crate::Scalar], 70080a9ef05SNatalie Beams grad1d: &[crate::Scalar], 70180a9ef05SNatalie Beams qref1d: &[crate::Scalar], 70280a9ef05SNatalie Beams qweight1d: &[crate::Scalar], 703594ef120SJeremy L Thompson ) -> Result<Basis<'a>> { 7049df49d7eSJed Brown Basis::create_tensor_H1( 7059df49d7eSJed Brown self, dim, ncomp, P1d, Q1d, interp1d, grad1d, qref1d, qweight1d, 7069df49d7eSJed Brown ) 7079df49d7eSJed Brown } 7089df49d7eSJed Brown 7097ed177dbSJed Brown /// Returns an $H^1$ Lagrange tensor-product Basis 7109df49d7eSJed Brown /// 7119df49d7eSJed Brown /// # arguments 7129df49d7eSJed Brown /// 7139df49d7eSJed Brown /// * `dim` - Topological dimension of element 7149df49d7eSJed Brown /// * `ncomp` - Number of field components (1 for scalar fields) 7159df49d7eSJed Brown /// * `P` - Number of Gauss-Lobatto nodes in one dimension. The 7169df49d7eSJed Brown /// polynomial degree of the resulting `Q_k` element is `k=P-1`. 7179df49d7eSJed Brown /// * `Q` - Number of quadrature points in one dimension 7189df49d7eSJed Brown /// * `qmode` - Distribution of the `Q` quadrature points (affects order of 7199df49d7eSJed Brown /// accuracy for the quadrature) 7209df49d7eSJed Brown /// 7219df49d7eSJed Brown /// ``` 7229df49d7eSJed Brown /// # use libceed::prelude::*; 7234d27c890SJeremy L Thompson /// # fn main() -> libceed::Result<()> { 7249df49d7eSJed Brown /// # let ceed = libceed::Ceed::default_init(); 725c68be7a2SJeremy L Thompson /// let b = ceed.basis_tensor_H1_Lagrange(2, 1, 3, 4, QuadMode::Gauss)?; 726c68be7a2SJeremy L Thompson /// # Ok(()) 727c68be7a2SJeremy L Thompson /// # } 7289df49d7eSJed Brown /// ``` 729594ef120SJeremy L Thompson pub fn basis_tensor_H1_Lagrange<'a>( 7309df49d7eSJed Brown &self, 7319df49d7eSJed Brown dim: usize, 7329df49d7eSJed Brown ncomp: usize, 7339df49d7eSJed Brown P: usize, 7349df49d7eSJed Brown Q: usize, 7359df49d7eSJed Brown qmode: QuadMode, 736594ef120SJeremy L Thompson ) -> Result<Basis<'a>> { 7379df49d7eSJed Brown Basis::create_tensor_H1_Lagrange(self, dim, ncomp, P, Q, qmode) 7389df49d7eSJed Brown } 7399df49d7eSJed Brown 7407ed177dbSJed Brown /// Returns an $H-1$ Basis 7419df49d7eSJed Brown /// 7429df49d7eSJed Brown /// # arguments 7439df49d7eSJed Brown /// 7449df49d7eSJed Brown /// * `topo` - Topology of element, e.g. hypercube, simplex, ect 7459df49d7eSJed Brown /// * `ncomp` - Number of field components (1 for scalar fields) 7469df49d7eSJed Brown /// * `nnodes` - Total number of nodes 7479df49d7eSJed Brown /// * `nqpts` - Total number of quadrature points 7489df49d7eSJed Brown /// * `interp` - Row-major `(nqpts * nnodes)` matrix expressing the values of 7499df49d7eSJed Brown /// nodal basis functions at quadrature points 75097c1c57aSSebastian Grimberg /// * `grad` - Row-major `(dim * nqpts * nnodes)` matrix expressing 7519df49d7eSJed Brown /// derivatives of nodal basis functions at quadrature points 7529df49d7eSJed Brown /// * `qref` - Array of length `nqpts` holding the locations of quadrature 7539df49d7eSJed Brown /// points on the reference element `[-1, 1]` 7549df49d7eSJed Brown /// * `qweight` - Array of length `nqpts` holding the quadrature weights on 7559df49d7eSJed Brown /// the reference element 7569df49d7eSJed Brown /// 7579df49d7eSJed Brown /// ``` 7589df49d7eSJed Brown /// # use libceed::prelude::*; 7594d27c890SJeremy L Thompson /// # fn main() -> libceed::Result<()> { 7609df49d7eSJed Brown /// # let ceed = libceed::Ceed::default_init(); 7619df49d7eSJed Brown /// let interp = [ 7629df49d7eSJed Brown /// 0.12000000, 7639df49d7eSJed Brown /// 0.48000000, 7649df49d7eSJed Brown /// -0.12000000, 7659df49d7eSJed Brown /// 0.48000000, 7669df49d7eSJed Brown /// 0.16000000, 7679df49d7eSJed Brown /// -0.12000000, 7689df49d7eSJed Brown /// -0.12000000, 7699df49d7eSJed Brown /// 0.48000000, 7709df49d7eSJed Brown /// 0.12000000, 7719df49d7eSJed Brown /// 0.16000000, 7729df49d7eSJed Brown /// 0.48000000, 7739df49d7eSJed Brown /// -0.12000000, 7749df49d7eSJed Brown /// -0.11111111, 7759df49d7eSJed Brown /// 0.44444444, 7769df49d7eSJed Brown /// -0.11111111, 7779df49d7eSJed Brown /// 0.44444444, 7789df49d7eSJed Brown /// 0.44444444, 7799df49d7eSJed Brown /// -0.11111111, 7809df49d7eSJed Brown /// -0.12000000, 7819df49d7eSJed Brown /// 0.16000000, 7829df49d7eSJed Brown /// -0.12000000, 7839df49d7eSJed Brown /// 0.48000000, 7849df49d7eSJed Brown /// 0.48000000, 7859df49d7eSJed Brown /// 0.12000000, 7869df49d7eSJed Brown /// ]; 7879df49d7eSJed Brown /// let grad = [ 7889df49d7eSJed Brown /// -1.40000000, 7899df49d7eSJed Brown /// 1.60000000, 7909df49d7eSJed Brown /// -0.20000000, 7919df49d7eSJed Brown /// -0.80000000, 7929df49d7eSJed Brown /// 0.80000000, 7939df49d7eSJed Brown /// 0.00000000, 7949df49d7eSJed Brown /// 0.20000000, 7959df49d7eSJed Brown /// -1.60000000, 7969df49d7eSJed Brown /// 1.40000000, 7979df49d7eSJed Brown /// -0.80000000, 7989df49d7eSJed Brown /// 0.80000000, 7999df49d7eSJed Brown /// 0.00000000, 8009df49d7eSJed Brown /// -0.33333333, 8019df49d7eSJed Brown /// 0.00000000, 8029df49d7eSJed Brown /// 0.33333333, 8039df49d7eSJed Brown /// -1.33333333, 8049df49d7eSJed Brown /// 1.33333333, 8059df49d7eSJed Brown /// 0.00000000, 8069df49d7eSJed Brown /// 0.20000000, 8079df49d7eSJed Brown /// 0.00000000, 8089df49d7eSJed Brown /// -0.20000000, 8099df49d7eSJed Brown /// -2.40000000, 8109df49d7eSJed Brown /// 2.40000000, 8119df49d7eSJed Brown /// 0.00000000, 8129df49d7eSJed Brown /// -1.40000000, 8139df49d7eSJed Brown /// -0.80000000, 8149df49d7eSJed Brown /// 0.00000000, 8159df49d7eSJed Brown /// 1.60000000, 8169df49d7eSJed Brown /// 0.80000000, 8179df49d7eSJed Brown /// -0.20000000, 8189df49d7eSJed Brown /// 0.20000000, 8199df49d7eSJed Brown /// -2.40000000, 8209df49d7eSJed Brown /// 0.00000000, 8219df49d7eSJed Brown /// 0.00000000, 8229df49d7eSJed Brown /// 2.40000000, 8239df49d7eSJed Brown /// -0.20000000, 8249df49d7eSJed Brown /// -0.33333333, 8259df49d7eSJed Brown /// -1.33333333, 8269df49d7eSJed Brown /// 0.00000000, 8279df49d7eSJed Brown /// 0.00000000, 8289df49d7eSJed Brown /// 1.33333333, 8299df49d7eSJed Brown /// 0.33333333, 8309df49d7eSJed Brown /// 0.20000000, 8319df49d7eSJed Brown /// -0.80000000, 8329df49d7eSJed Brown /// 0.00000000, 8339df49d7eSJed Brown /// -1.60000000, 8349df49d7eSJed Brown /// 0.80000000, 8359df49d7eSJed Brown /// 1.40000000, 8369df49d7eSJed Brown /// ]; 8379df49d7eSJed Brown /// let qref = [ 8389df49d7eSJed Brown /// 0.20000000, 0.60000000, 0.33333333, 0.20000000, 0.20000000, 0.20000000, 0.33333333, 8399df49d7eSJed Brown /// 0.60000000, 8409df49d7eSJed Brown /// ]; 8419df49d7eSJed Brown /// let qweight = [0.26041667, 0.26041667, -0.28125000, 0.26041667]; 842c68be7a2SJeremy L Thompson /// let b = ceed.basis_H1( 8439df49d7eSJed Brown /// ElemTopology::Triangle, 8449df49d7eSJed Brown /// 1, 8459df49d7eSJed Brown /// 6, 8469df49d7eSJed Brown /// 4, 8479df49d7eSJed Brown /// &interp, 8489df49d7eSJed Brown /// &grad, 8499df49d7eSJed Brown /// &qref, 8509df49d7eSJed Brown /// &qweight, 851c68be7a2SJeremy L Thompson /// )?; 852c68be7a2SJeremy L Thompson /// # Ok(()) 853c68be7a2SJeremy L Thompson /// # } 8549df49d7eSJed Brown /// ``` 855594ef120SJeremy L Thompson pub fn basis_H1<'a>( 8569df49d7eSJed Brown &self, 8579df49d7eSJed Brown topo: ElemTopology, 8589df49d7eSJed Brown ncomp: usize, 8599df49d7eSJed Brown nnodes: usize, 8609df49d7eSJed Brown nqpts: usize, 86180a9ef05SNatalie Beams interp: &[crate::Scalar], 86280a9ef05SNatalie Beams grad: &[crate::Scalar], 86380a9ef05SNatalie Beams qref: &[crate::Scalar], 86480a9ef05SNatalie Beams qweight: &[crate::Scalar], 865594ef120SJeremy L Thompson ) -> Result<Basis<'a>> { 8669df49d7eSJed Brown Basis::create_H1( 8679df49d7eSJed Brown self, topo, ncomp, nnodes, nqpts, interp, grad, qref, qweight, 8689df49d7eSJed Brown ) 8699df49d7eSJed Brown } 8709df49d7eSJed Brown 87197c1c57aSSebastian Grimberg /// Returns an $H(div)$ Basis 87297c1c57aSSebastian Grimberg /// 87397c1c57aSSebastian Grimberg /// # arguments 87497c1c57aSSebastian Grimberg /// 87597c1c57aSSebastian Grimberg /// * `topo` - Topology of element, e.g. hypercube, simplex, ect 87697c1c57aSSebastian Grimberg /// * `ncomp` - Number of field components (1 for scalar fields) 87797c1c57aSSebastian Grimberg /// * `nnodes` - Total number of nodes 87897c1c57aSSebastian Grimberg /// * `nqpts` - Total number of quadrature points 87997c1c57aSSebastian Grimberg /// * `interp` - Row-major `(dim * nqpts * nnodes)` matrix expressing the 88097c1c57aSSebastian Grimberg /// values of basis functions at quadrature points 88197c1c57aSSebastian Grimberg /// * `div` - Row-major `(nqpts * nnodes)` matrix expressing the 88297c1c57aSSebastian Grimberg /// divergence of basis functions at quadrature points 88397c1c57aSSebastian Grimberg /// * `qref` - Array of length `nqpts` holding the locations of quadrature 88497c1c57aSSebastian Grimberg /// points on the reference element `[-1, 1]` 88597c1c57aSSebastian Grimberg /// * `qweight` - Array of length `nqpts` holding the quadrature weights on 88697c1c57aSSebastian Grimberg /// the reference element 88797c1c57aSSebastian Grimberg /// 88897c1c57aSSebastian Grimberg /// ``` 88997c1c57aSSebastian Grimberg /// # use libceed::prelude::*; 89097c1c57aSSebastian Grimberg /// # fn main() -> libceed::Result<()> { 89197c1c57aSSebastian Grimberg /// # let ceed = libceed::Ceed::default_init(); 89297c1c57aSSebastian Grimberg /// let interp = [ 89397c1c57aSSebastian Grimberg /// 0.00000000, 89497c1c57aSSebastian Grimberg /// -1.57735027, 89597c1c57aSSebastian Grimberg /// 0.57735027, 89697c1c57aSSebastian Grimberg /// 0.00000000, 89797c1c57aSSebastian Grimberg /// 0.00000000, 89897c1c57aSSebastian Grimberg /// -1.57735027, 89997c1c57aSSebastian Grimberg /// 0.57735027, 90097c1c57aSSebastian Grimberg /// 0.00000000, 90197c1c57aSSebastian Grimberg /// 0.00000000, 90297c1c57aSSebastian Grimberg /// -1.57735027, 90397c1c57aSSebastian Grimberg /// 0.57735027, 90497c1c57aSSebastian Grimberg /// 0.00000000, 90597c1c57aSSebastian Grimberg /// 0.00000000, 90697c1c57aSSebastian Grimberg /// -0.42264973, 90797c1c57aSSebastian Grimberg /// -0.57735027, 90897c1c57aSSebastian Grimberg /// 0.00000000, 90997c1c57aSSebastian Grimberg /// 0.42264973, 91097c1c57aSSebastian Grimberg /// 0.00000000, 91197c1c57aSSebastian Grimberg /// 0.00000000, 91297c1c57aSSebastian Grimberg /// 0.57735027, 91397c1c57aSSebastian Grimberg /// 0.42264973, 91497c1c57aSSebastian Grimberg /// 0.00000000, 91597c1c57aSSebastian Grimberg /// 0.00000000, 91697c1c57aSSebastian Grimberg /// 0.57735027, 91797c1c57aSSebastian Grimberg /// 1.57735027, 91897c1c57aSSebastian Grimberg /// 0.00000000, 91997c1c57aSSebastian Grimberg /// 0.00000000, 92097c1c57aSSebastian Grimberg /// -0.57735027, 92197c1c57aSSebastian Grimberg /// 1.57735027, 92297c1c57aSSebastian Grimberg /// 0.00000000, 92397c1c57aSSebastian Grimberg /// 0.00000000, 92497c1c57aSSebastian Grimberg /// -0.57735027, 92597c1c57aSSebastian Grimberg /// ]; 92697c1c57aSSebastian Grimberg /// let div = [ 92797c1c57aSSebastian Grimberg /// -1.00000000, 92897c1c57aSSebastian Grimberg /// 1.00000000, 92997c1c57aSSebastian Grimberg /// -1.00000000, 93097c1c57aSSebastian Grimberg /// 1.00000000, 93197c1c57aSSebastian Grimberg /// -1.00000000, 93297c1c57aSSebastian Grimberg /// 1.00000000, 93397c1c57aSSebastian Grimberg /// -1.00000000, 93497c1c57aSSebastian Grimberg /// 1.00000000, 93597c1c57aSSebastian Grimberg /// -1.00000000, 93697c1c57aSSebastian Grimberg /// 1.00000000, 93797c1c57aSSebastian Grimberg /// -1.00000000, 93897c1c57aSSebastian Grimberg /// 1.00000000, 93997c1c57aSSebastian Grimberg /// -1.00000000, 94097c1c57aSSebastian Grimberg /// 1.00000000, 94197c1c57aSSebastian Grimberg /// -1.00000000, 94297c1c57aSSebastian Grimberg /// 1.00000000, 94397c1c57aSSebastian Grimberg /// ]; 94497c1c57aSSebastian Grimberg /// let qref = [ 94597c1c57aSSebastian Grimberg /// 0.57735026, 0.57735026, 0.57735026, 0.57735026, 0.57735026, 0.57735026, 0.57735026, 94697c1c57aSSebastian Grimberg /// 0.57735026, 94797c1c57aSSebastian Grimberg /// ]; 94897c1c57aSSebastian Grimberg /// let qweight = [1.00000000, 1.00000000, 1.00000000, 1.00000000]; 94997c1c57aSSebastian Grimberg /// let b = ceed.basis_Hdiv(ElemTopology::Quad, 1, 4, 4, &interp, &div, &qref, &qweight)?; 95097c1c57aSSebastian Grimberg /// # Ok(()) 95197c1c57aSSebastian Grimberg /// # } 95297c1c57aSSebastian Grimberg /// ``` 95397c1c57aSSebastian Grimberg pub fn basis_Hdiv<'a>( 95497c1c57aSSebastian Grimberg &self, 95597c1c57aSSebastian Grimberg topo: ElemTopology, 95697c1c57aSSebastian Grimberg ncomp: usize, 95797c1c57aSSebastian Grimberg nnodes: usize, 95897c1c57aSSebastian Grimberg nqpts: usize, 95997c1c57aSSebastian Grimberg interp: &[crate::Scalar], 96097c1c57aSSebastian Grimberg div: &[crate::Scalar], 96197c1c57aSSebastian Grimberg qref: &[crate::Scalar], 96297c1c57aSSebastian Grimberg qweight: &[crate::Scalar], 96397c1c57aSSebastian Grimberg ) -> Result<Basis<'a>> { 96497c1c57aSSebastian Grimberg Basis::create_Hdiv(self, topo, ncomp, nnodes, nqpts, interp, div, qref, qweight) 96597c1c57aSSebastian Grimberg } 96697c1c57aSSebastian Grimberg 96797c1c57aSSebastian Grimberg /// Returns an $H(curl)$ Basis 96897c1c57aSSebastian Grimberg /// 96997c1c57aSSebastian Grimberg /// # arguments 97097c1c57aSSebastian Grimberg /// 97197c1c57aSSebastian Grimberg /// * `topo` - Topology of element, e.g. hypercube, simplex, ect 97297c1c57aSSebastian Grimberg /// * `ncomp` - Number of field components (1 for scalar fields) 97397c1c57aSSebastian Grimberg /// * `nnodes` - Total number of nodes 97497c1c57aSSebastian Grimberg /// * `nqpts` - Total number of quadrature points 97597c1c57aSSebastian Grimberg /// * `interp` - Row-major `(dim * nqpts * nnodes)` matrix expressing the 97697c1c57aSSebastian Grimberg /// values of basis functions at quadrature points 97797c1c57aSSebastian Grimberg /// * `curl` - Row-major `(curl_comp * nqpts * nnodes)`, `curl_comp = 1 if 97897c1c57aSSebastian Grimberg /// dim < 3 else dim` matrix expressing the curl of basis 97997c1c57aSSebastian Grimberg /// functions at quadrature points 98097c1c57aSSebastian Grimberg /// * `qref` - Array of length `nqpts` holding the locations of quadrature 98197c1c57aSSebastian Grimberg /// points on the reference element `[-1, 1]` 98297c1c57aSSebastian Grimberg /// * `qweight` - Array of length `nqpts` holding the quadrature weights on 98397c1c57aSSebastian Grimberg /// the reference element 98497c1c57aSSebastian Grimberg /// 98597c1c57aSSebastian Grimberg /// ``` 98697c1c57aSSebastian Grimberg /// # use libceed::prelude::*; 98797c1c57aSSebastian Grimberg /// # fn main() -> libceed::Result<()> { 98897c1c57aSSebastian Grimberg /// # let ceed = libceed::Ceed::default_init(); 98997c1c57aSSebastian Grimberg /// let interp = [ 99097c1c57aSSebastian Grimberg /// -0.20000000, 99197c1c57aSSebastian Grimberg /// 0.20000000, 99297c1c57aSSebastian Grimberg /// 0.80000000, 99397c1c57aSSebastian Grimberg /// -0.20000000, 99497c1c57aSSebastian Grimberg /// 0.20000000, 99597c1c57aSSebastian Grimberg /// 0.80000000, 99697c1c57aSSebastian Grimberg /// -0.33333333, 99797c1c57aSSebastian Grimberg /// 0.33333333, 99897c1c57aSSebastian Grimberg /// 0.66666667, 99997c1c57aSSebastian Grimberg /// -0.60000000, 100097c1c57aSSebastian Grimberg /// 0.60000000, 100197c1c57aSSebastian Grimberg /// 0.40000000, 100297c1c57aSSebastian Grimberg /// 0.20000000, 100397c1c57aSSebastian Grimberg /// 0.80000000, 100497c1c57aSSebastian Grimberg /// 0.20000000, 100597c1c57aSSebastian Grimberg /// 0.60000000, 100697c1c57aSSebastian Grimberg /// 0.40000000, 100797c1c57aSSebastian Grimberg /// 0.60000000, 100897c1c57aSSebastian Grimberg /// 0.33333333, 100997c1c57aSSebastian Grimberg /// 0.66666667, 101097c1c57aSSebastian Grimberg /// 0.33333333, 101197c1c57aSSebastian Grimberg /// 0.20000000, 101297c1c57aSSebastian Grimberg /// 0.80000000, 101397c1c57aSSebastian Grimberg /// 0.20000000, 101497c1c57aSSebastian Grimberg /// ]; 101597c1c57aSSebastian Grimberg /// let curl = [ 101697c1c57aSSebastian Grimberg /// 2.00000000, 101797c1c57aSSebastian Grimberg /// -2.00000000, 101897c1c57aSSebastian Grimberg /// -2.00000000, 101997c1c57aSSebastian Grimberg /// 2.00000000, 102097c1c57aSSebastian Grimberg /// -2.00000000, 102197c1c57aSSebastian Grimberg /// -2.00000000, 102297c1c57aSSebastian Grimberg /// 2.00000000, 102397c1c57aSSebastian Grimberg /// -2.00000000, 102497c1c57aSSebastian Grimberg /// -2.00000000, 102597c1c57aSSebastian Grimberg /// 2.00000000, 102697c1c57aSSebastian Grimberg /// -2.00000000, 102797c1c57aSSebastian Grimberg /// -2.00000000, 102897c1c57aSSebastian Grimberg /// ]; 102997c1c57aSSebastian Grimberg /// let qref = [ 103097c1c57aSSebastian Grimberg /// 0.20000000, 0.60000000, 0.33333333, 0.20000000, 0.20000000, 0.20000000, 0.33333333, 103197c1c57aSSebastian Grimberg /// 0.60000000, 103297c1c57aSSebastian Grimberg /// ]; 103397c1c57aSSebastian Grimberg /// let qweight = [0.26041667, 0.26041667, -0.28125000, 0.26041667]; 103497c1c57aSSebastian Grimberg /// let b = ceed.basis_Hcurl( 103597c1c57aSSebastian Grimberg /// ElemTopology::Triangle, 103697c1c57aSSebastian Grimberg /// 1, 103797c1c57aSSebastian Grimberg /// 3, 103897c1c57aSSebastian Grimberg /// 4, 103997c1c57aSSebastian Grimberg /// &interp, 104097c1c57aSSebastian Grimberg /// &curl, 104197c1c57aSSebastian Grimberg /// &qref, 104297c1c57aSSebastian Grimberg /// &qweight, 104397c1c57aSSebastian Grimberg /// )?; 104497c1c57aSSebastian Grimberg /// # Ok(()) 104597c1c57aSSebastian Grimberg /// # } 104697c1c57aSSebastian Grimberg /// ``` 104797c1c57aSSebastian Grimberg pub fn basis_Hcurl<'a>( 104897c1c57aSSebastian Grimberg &self, 104997c1c57aSSebastian Grimberg topo: ElemTopology, 105097c1c57aSSebastian Grimberg ncomp: usize, 105197c1c57aSSebastian Grimberg nnodes: usize, 105297c1c57aSSebastian Grimberg nqpts: usize, 105397c1c57aSSebastian Grimberg interp: &[crate::Scalar], 105497c1c57aSSebastian Grimberg curl: &[crate::Scalar], 105597c1c57aSSebastian Grimberg qref: &[crate::Scalar], 105697c1c57aSSebastian Grimberg qweight: &[crate::Scalar], 105797c1c57aSSebastian Grimberg ) -> Result<Basis<'a>> { 105897c1c57aSSebastian Grimberg Basis::create_Hcurl( 105997c1c57aSSebastian Grimberg self, topo, ncomp, nnodes, nqpts, interp, curl, qref, qweight, 106097c1c57aSSebastian Grimberg ) 106197c1c57aSSebastian Grimberg } 106297c1c57aSSebastian Grimberg 10637ed177dbSJed Brown /// Returns a QFunction for evaluating interior (volumetric) terms 10647ed177dbSJed Brown /// of a weak form corresponding to the $L^2$ inner product 10657ed177dbSJed Brown /// 10667ed177dbSJed Brown /// $$ 10677ed177dbSJed Brown /// \langle v, F(u) \rangle = \int_\Omega v \cdot f_0 \left( u, \nabla u \right) + \left( \nabla v \right) : f_1 \left( u, \nabla u \right), 10687ed177dbSJed Brown /// $$ 10697ed177dbSJed Brown /// 10707ed177dbSJed Brown /// where $v \cdot f_0$ represents contraction over fields and $\nabla v : f_1$ 10717ed177dbSJed Brown /// represents contraction over both fields and spatial dimensions. 10729df49d7eSJed Brown /// 10739df49d7eSJed Brown /// # arguments 10749df49d7eSJed Brown /// 10759df49d7eSJed Brown /// * `vlength` - Vector length. Caller must ensure that number of 10769df49d7eSJed Brown /// quadrature points is a multiple of vlength. 10777ed177dbSJed Brown /// * `f` - Boxed closure to evaluate weak form at quadrature points. 10789df49d7eSJed Brown /// 10799df49d7eSJed Brown /// ``` 10809df49d7eSJed Brown /// # use libceed::prelude::*; 10814d27c890SJeremy L Thompson /// # fn main() -> libceed::Result<()> { 10829df49d7eSJed Brown /// # let ceed = libceed::Ceed::default_init(); 10839df49d7eSJed Brown /// let mut user_f = |[u, weights, ..]: QFunctionInputs, [v, ..]: QFunctionOutputs| { 10849df49d7eSJed Brown /// // Iterate over quadrature points 10859df49d7eSJed Brown /// v.iter_mut() 10869df49d7eSJed Brown /// .zip(u.iter().zip(weights.iter())) 10879df49d7eSJed Brown /// .for_each(|(v, (u, w))| *v = u * w); 10889df49d7eSJed Brown /// 10899df49d7eSJed Brown /// // Return clean error code 10909df49d7eSJed Brown /// 0 10919df49d7eSJed Brown /// }; 10929df49d7eSJed Brown /// 1093c68be7a2SJeremy L Thompson /// let qf = ceed.q_function_interior(1, Box::new(user_f))?; 1094c68be7a2SJeremy L Thompson /// # Ok(()) 1095c68be7a2SJeremy L Thompson /// # } 10969df49d7eSJed Brown /// ``` 1097594ef120SJeremy L Thompson pub fn q_function_interior<'a>( 10989df49d7eSJed Brown &self, 10999df49d7eSJed Brown vlength: usize, 11009df49d7eSJed Brown f: Box<qfunction::QFunctionUserClosure>, 1101594ef120SJeremy L Thompson ) -> Result<QFunction<'a>> { 11029df49d7eSJed Brown QFunction::create(self, vlength, f) 11039df49d7eSJed Brown } 11049df49d7eSJed Brown 11057ed177dbSJed Brown /// Returns a QFunction for evaluating interior (volumetric) terms 11069df49d7eSJed Brown /// created by name 11079df49d7eSJed Brown /// 11087ed177dbSJed Brown /// # arguments 11097ed177dbSJed Brown /// 11107ed177dbSJed Brown /// * `name` - name of QFunction from libCEED gallery 11117ed177dbSJed Brown /// 11129df49d7eSJed Brown /// ``` 11139df49d7eSJed Brown /// # use libceed::prelude::*; 11144d27c890SJeremy L Thompson /// # fn main() -> libceed::Result<()> { 11159df49d7eSJed Brown /// # let ceed = libceed::Ceed::default_init(); 1116c68be7a2SJeremy L Thompson /// let qf = ceed.q_function_interior_by_name("Mass1DBuild")?; 1117c68be7a2SJeremy L Thompson /// # Ok(()) 1118c68be7a2SJeremy L Thompson /// # } 11199df49d7eSJed Brown /// ``` 1120594ef120SJeremy L Thompson pub fn q_function_interior_by_name<'a>(&self, name: &str) -> Result<QFunctionByName<'a>> { 11219df49d7eSJed Brown QFunctionByName::create(self, name) 11229df49d7eSJed Brown } 11239df49d7eSJed Brown 11247ed177dbSJed Brown /// Returns an Operator and associate a QFunction. A Basis and 11257ed177dbSJed Brown /// ElemRestriction can be associated with QFunction fields via 11269df49d7eSJed Brown /// set_field(). 11279df49d7eSJed Brown /// 11287ed177dbSJed Brown /// # arguments 11297ed177dbSJed Brown /// 11309df49d7eSJed Brown /// * `qf` - QFunction defining the action of the operator at quadrature 11319df49d7eSJed Brown /// points 11329df49d7eSJed Brown /// * `dqf` - QFunction defining the action of the Jacobian of the qf (or 11339df49d7eSJed Brown /// qfunction_none) 11349df49d7eSJed Brown /// * `dqfT` - QFunction defining the action of the transpose of the 11359df49d7eSJed Brown /// Jacobian of the qf (or qfunction_none) 11369df49d7eSJed Brown /// 11379df49d7eSJed Brown /// ``` 11389df49d7eSJed Brown /// # use libceed::prelude::*; 11394d27c890SJeremy L Thompson /// # fn main() -> libceed::Result<()> { 11409df49d7eSJed Brown /// # let ceed = libceed::Ceed::default_init(); 1141c68be7a2SJeremy L Thompson /// let qf = ceed.q_function_interior_by_name("Mass1DBuild")?; 1142c68be7a2SJeremy L Thompson /// let op = ceed.operator(&qf, QFunctionOpt::None, QFunctionOpt::None)?; 1143c68be7a2SJeremy L Thompson /// # Ok(()) 1144c68be7a2SJeremy L Thompson /// # } 11459df49d7eSJed Brown /// ``` 1146594ef120SJeremy L Thompson pub fn operator<'a, 'b>( 11479df49d7eSJed Brown &self, 11489df49d7eSJed Brown qf: impl Into<QFunctionOpt<'b>>, 11499df49d7eSJed Brown dqf: impl Into<QFunctionOpt<'b>>, 11509df49d7eSJed Brown dqfT: impl Into<QFunctionOpt<'b>>, 1151594ef120SJeremy L Thompson ) -> Result<Operator<'a>> { 11529df49d7eSJed Brown Operator::create(self, qf, dqf, dqfT) 11539df49d7eSJed Brown } 11549df49d7eSJed Brown 11559df49d7eSJed Brown /// Returns an Operator that composes the action of several Operators 11569df49d7eSJed Brown /// 11579df49d7eSJed Brown /// ``` 11589df49d7eSJed Brown /// # use libceed::prelude::*; 11594d27c890SJeremy L Thompson /// # fn main() -> libceed::Result<()> { 11609df49d7eSJed Brown /// # let ceed = libceed::Ceed::default_init(); 1161c68be7a2SJeremy L Thompson /// let op = ceed.composite_operator()?; 1162c68be7a2SJeremy L Thompson /// # Ok(()) 1163c68be7a2SJeremy L Thompson /// # } 11649df49d7eSJed Brown /// ``` 1165594ef120SJeremy L Thompson pub fn composite_operator<'a>(&self) -> Result<CompositeOperator<'a>> { 11669df49d7eSJed Brown CompositeOperator::create(self) 11679df49d7eSJed Brown } 11689df49d7eSJed Brown } 11699df49d7eSJed Brown 11709df49d7eSJed Brown // ----------------------------------------------------------------------------- 11719df49d7eSJed Brown // Tests 11729df49d7eSJed Brown // ----------------------------------------------------------------------------- 11739df49d7eSJed Brown #[cfg(test)] 11749df49d7eSJed Brown mod tests { 11759df49d7eSJed Brown use super::*; 11769df49d7eSJed Brown 117789d15d5fSJeremy L Thompson fn ceed_t501() -> Result<()> { 11789df49d7eSJed Brown let resource = "/cpu/self/ref/blocked"; 11799df49d7eSJed Brown let ceed = Ceed::init(resource); 11809df49d7eSJed Brown let nelem = 4; 11819df49d7eSJed Brown let p = 3; 11829df49d7eSJed Brown let q = 4; 11839df49d7eSJed Brown let ndofs = p * nelem - nelem + 1; 11849df49d7eSJed Brown 11859df49d7eSJed Brown // Vectors 11869df49d7eSJed Brown let x = ceed.vector_from_slice(&[-1., -0.5, 0.0, 0.5, 1.0])?; 11879df49d7eSJed Brown let mut qdata = ceed.vector(nelem * q)?; 11889df49d7eSJed Brown qdata.set_value(0.0)?; 11899df49d7eSJed Brown let mut u = ceed.vector(ndofs)?; 11909df49d7eSJed Brown u.set_value(1.0)?; 11919df49d7eSJed Brown let mut v = ceed.vector(ndofs)?; 11929df49d7eSJed Brown v.set_value(0.0)?; 11939df49d7eSJed Brown 11949df49d7eSJed Brown // Restrictions 11959df49d7eSJed Brown let mut indx: Vec<i32> = vec![0; 2 * nelem]; 11969df49d7eSJed Brown for i in 0..nelem { 11971ce8139fSJeremy L Thompson for j in 0..2 { 11981ce8139fSJeremy L Thompson indx[2 * i + j] = (i + j) as i32; 11991ce8139fSJeremy L Thompson } 12009df49d7eSJed Brown } 12019df49d7eSJed Brown let rx = ceed.elem_restriction(nelem, 2, 1, 1, nelem + 1, MemType::Host, &indx)?; 12029df49d7eSJed Brown let mut indu: Vec<i32> = vec![0; p * nelem]; 12039df49d7eSJed Brown for i in 0..nelem { 12041ce8139fSJeremy L Thompson for j in 0..p { 12051ce8139fSJeremy L Thompson indu[p * i + j] = (i + j) as i32; 12069df49d7eSJed Brown } 12071ce8139fSJeremy L Thompson } 12081ce8139fSJeremy L Thompson let ru = ceed.elem_restriction(nelem, p, 1, 1, ndofs, MemType::Host, &indu)?; 12099df49d7eSJed Brown let strides: [i32; 3] = [1, q as i32, q as i32]; 12109df49d7eSJed Brown let rq = ceed.strided_elem_restriction(nelem, q, 1, q * nelem, strides)?; 12119df49d7eSJed Brown 12129df49d7eSJed Brown // Bases 12139df49d7eSJed Brown let bx = ceed.basis_tensor_H1_Lagrange(1, 1, 2, q, QuadMode::Gauss)?; 12149df49d7eSJed Brown let bu = ceed.basis_tensor_H1_Lagrange(1, 1, p, q, QuadMode::Gauss)?; 12159df49d7eSJed Brown 12169df49d7eSJed Brown // Build quadrature data 12179df49d7eSJed Brown let qf_build = ceed.q_function_interior_by_name("Mass1DBuild")?; 12189df49d7eSJed Brown ceed.operator(&qf_build, QFunctionOpt::None, QFunctionOpt::None)? 12199df49d7eSJed Brown .field("dx", &rx, &bx, VectorOpt::Active)? 12209df49d7eSJed Brown .field("weights", ElemRestrictionOpt::None, &bx, VectorOpt::None)? 1221356036faSJeremy L Thompson .field("qdata", &rq, BasisOpt::None, VectorOpt::Active)? 12229df49d7eSJed Brown .apply(&x, &mut qdata)?; 12239df49d7eSJed Brown 12249df49d7eSJed Brown // Mass operator 12259df49d7eSJed Brown let qf_mass = ceed.q_function_interior_by_name("MassApply")?; 12269df49d7eSJed Brown let op_mass = ceed 12279df49d7eSJed Brown .operator(&qf_mass, QFunctionOpt::None, QFunctionOpt::None)? 12289df49d7eSJed Brown .field("u", &ru, &bu, VectorOpt::Active)? 1229356036faSJeremy L Thompson .field("qdata", &rq, BasisOpt::None, &qdata)? 12306f97ff0aSJeremy L Thompson .field("v", &ru, &bu, VectorOpt::Active)? 12316f97ff0aSJeremy L Thompson .check()?; 12329df49d7eSJed Brown 12339df49d7eSJed Brown v.set_value(0.0)?; 12349df49d7eSJed Brown op_mass.apply(&u, &mut v)?; 12359df49d7eSJed Brown 12369df49d7eSJed Brown // Check 1237e78171edSJeremy L Thompson let sum: Scalar = v.view()?.iter().sum(); 12381ce8139fSJeremy L Thompson let error: Scalar = (sum - 2.0).abs(); 12399df49d7eSJed Brown assert!( 12404b61a5a0SRezgar Shakeri error < 50.0 * EPSILON, 12411ce8139fSJeremy L Thompson "Incorrect interval length computed. Expected: 2.0, Found: {}, Error: {:.12e}", 12421ce8139fSJeremy L Thompson sum, 12431ce8139fSJeremy L Thompson error 12449df49d7eSJed Brown ); 124589d15d5fSJeremy L Thompson Ok(()) 12469df49d7eSJed Brown } 12479df49d7eSJed Brown 12489df49d7eSJed Brown #[test] 12499df49d7eSJed Brown fn test_ceed_t501() { 12509df49d7eSJed Brown assert!(ceed_t501().is_ok()); 12519df49d7eSJed Brown } 12529df49d7eSJed Brown } 12539df49d7eSJed Brown 12549df49d7eSJed Brown // ----------------------------------------------------------------------------- 1255