13d8e8822SJeremy L Thompson // Copyright (c) 2017-2022, 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)] 1641142270cSJeremy L Thompson pub(crate) fn check_error(ceed_ptr: bind_ceed::Ceed, ierr: i32) -> Result<i32> { 1651142270cSJeremy L Thompson // Return early if code is clean 1661142270cSJeremy L Thompson if ierr == bind_ceed::CeedErrorType_CEED_ERROR_SUCCESS { 1671142270cSJeremy L Thompson return Ok(ierr); 1681142270cSJeremy L Thompson } 1691142270cSJeremy L Thompson // Retrieve error message 1701142270cSJeremy L Thompson let mut ptr: *const std::os::raw::c_char = std::ptr::null_mut(); 1711142270cSJeremy L Thompson let c_str = unsafe { 1721142270cSJeremy L Thompson bind_ceed::CeedGetErrorMessage(ceed_ptr, &mut ptr); 1731142270cSJeremy L Thompson std::ffi::CStr::from_ptr(ptr) 1741142270cSJeremy L Thompson }; 1751142270cSJeremy L Thompson let message = c_str.to_string_lossy().to_string(); 1761142270cSJeremy L Thompson // Panic if negative code, otherwise return error 1771142270cSJeremy L Thompson if ierr < bind_ceed::CeedErrorType_CEED_ERROR_SUCCESS { 1781142270cSJeremy L Thompson panic!("{}", message); 1791142270cSJeremy L Thompson } 1802ba8e59cSJeremy L Thompson Err(Error { message }) 1811142270cSJeremy L Thompson } 1821142270cSJeremy L Thompson 1831142270cSJeremy L Thompson // ----------------------------------------------------------------------------- 1849df49d7eSJed Brown // Ceed error handler 1859df49d7eSJed Brown // ----------------------------------------------------------------------------- 1862ba8e59cSJeremy L Thompson pub enum ErrorHandler { 1879df49d7eSJed Brown ErrorAbort, 1889df49d7eSJed Brown ErrorExit, 1899df49d7eSJed Brown ErrorReturn, 1909df49d7eSJed Brown ErrorStore, 1919df49d7eSJed Brown } 1929df49d7eSJed Brown 1939df49d7eSJed Brown // ----------------------------------------------------------------------------- 1949df49d7eSJed Brown // Ceed context wrapper 1959df49d7eSJed Brown // ----------------------------------------------------------------------------- 1969df49d7eSJed Brown /// A Ceed is a library context representing control of a logical hardware 1979df49d7eSJed Brown /// resource. 1989df49d7eSJed Brown #[derive(Debug)] 1999df49d7eSJed Brown pub struct Ceed { 2009df49d7eSJed Brown ptr: bind_ceed::Ceed, 2019df49d7eSJed Brown } 2029df49d7eSJed Brown 2039df49d7eSJed Brown // ----------------------------------------------------------------------------- 2049df49d7eSJed Brown // Destructor 2059df49d7eSJed Brown // ----------------------------------------------------------------------------- 2069df49d7eSJed Brown impl Drop for Ceed { 2079df49d7eSJed Brown fn drop(&mut self) { 2089df49d7eSJed Brown unsafe { 2099df49d7eSJed Brown bind_ceed::CeedDestroy(&mut self.ptr); 2109df49d7eSJed Brown } 2119df49d7eSJed Brown } 2129df49d7eSJed Brown } 2139df49d7eSJed Brown 2149df49d7eSJed Brown // ----------------------------------------------------------------------------- 21559189cfaSJeremy L Thompson // Cloning 21659189cfaSJeremy L Thompson // ----------------------------------------------------------------------------- 21759189cfaSJeremy L Thompson impl Clone for Ceed { 21859189cfaSJeremy L Thompson /// Perform a shallow clone of a Ceed context 21959189cfaSJeremy L Thompson /// 22059189cfaSJeremy L Thompson /// ``` 22159189cfaSJeremy L Thompson /// let ceed = libceed::Ceed::init("/cpu/self/ref/serial"); 22259189cfaSJeremy L Thompson /// let ceed_clone = ceed.clone(); 22359189cfaSJeremy L Thompson /// 2247ed177dbSJed Brown /// println!(" original:{} \n clone:{}", ceed, ceed_clone); 22559189cfaSJeremy L Thompson /// ``` 22659189cfaSJeremy L Thompson fn clone(&self) -> Self { 22759189cfaSJeremy L Thompson let mut ptr_clone = std::ptr::null_mut(); 22859189cfaSJeremy L Thompson let ierr = unsafe { bind_ceed::CeedReferenceCopy(self.ptr, &mut ptr_clone) }; 22959189cfaSJeremy L Thompson self.check_error(ierr).expect("failed to clone Ceed"); 23059189cfaSJeremy L Thompson Self { ptr: ptr_clone } 23159189cfaSJeremy L Thompson } 23259189cfaSJeremy L Thompson } 23359189cfaSJeremy L Thompson 23459189cfaSJeremy L Thompson // ----------------------------------------------------------------------------- 2359df49d7eSJed Brown // Display 2369df49d7eSJed Brown // ----------------------------------------------------------------------------- 2379df49d7eSJed Brown impl fmt::Display for Ceed { 2389df49d7eSJed Brown /// View a Ceed 2399df49d7eSJed Brown /// 2409df49d7eSJed Brown /// ``` 2419df49d7eSJed Brown /// let ceed = libceed::Ceed::init("/cpu/self/ref/serial"); 2429df49d7eSJed Brown /// println!("{}", ceed); 2439df49d7eSJed Brown /// ``` 2449df49d7eSJed Brown fn fmt(&self, f: &mut fmt::Formatter) -> fmt::Result { 2459df49d7eSJed Brown let mut ptr = std::ptr::null_mut(); 2469df49d7eSJed Brown let mut sizeloc = crate::MAX_BUFFER_LENGTH; 2479df49d7eSJed Brown let cstring = unsafe { 2489df49d7eSJed Brown let file = bind_ceed::open_memstream(&mut ptr, &mut sizeloc); 2499df49d7eSJed Brown bind_ceed::CeedView(self.ptr, file); 2509df49d7eSJed Brown bind_ceed::fclose(file); 2519df49d7eSJed Brown CString::from_raw(ptr) 2529df49d7eSJed Brown }; 2539df49d7eSJed Brown cstring.to_string_lossy().fmt(f) 2549df49d7eSJed Brown } 2559df49d7eSJed Brown } 2569df49d7eSJed Brown 2579df49d7eSJed Brown static REGISTER: Once = Once::new(); 2589df49d7eSJed Brown 2599df49d7eSJed Brown // ----------------------------------------------------------------------------- 2609df49d7eSJed Brown // Object constructors 2619df49d7eSJed Brown // ----------------------------------------------------------------------------- 2629df49d7eSJed Brown impl Ceed { 2637ed177dbSJed Brown #[cfg_attr(feature = "katexit", katexit::katexit)] 2649df49d7eSJed Brown /// Returns a Ceed context initialized with the specified resource 2659df49d7eSJed Brown /// 2669df49d7eSJed Brown /// # arguments 2679df49d7eSJed Brown /// 2689df49d7eSJed Brown /// * `resource` - Resource to use, e.g., "/cpu/self" 2699df49d7eSJed Brown /// 2709df49d7eSJed Brown /// ``` 2719df49d7eSJed Brown /// let ceed = libceed::Ceed::init("/cpu/self/ref/serial"); 2729df49d7eSJed Brown /// ``` 2739df49d7eSJed Brown pub fn init(resource: &str) -> Self { 2742ba8e59cSJeremy L Thompson Ceed::init_with_error_handler(resource, ErrorHandler::ErrorStore) 2759df49d7eSJed Brown } 2769df49d7eSJed Brown 2779df49d7eSJed Brown /// Returns a Ceed context initialized with the specified resource 2789df49d7eSJed Brown /// 2799df49d7eSJed Brown /// # arguments 2809df49d7eSJed Brown /// 2819df49d7eSJed Brown /// * `resource` - Resource to use, e.g., "/cpu/self" 2829df49d7eSJed Brown /// 2839df49d7eSJed Brown /// ``` 2849df49d7eSJed Brown /// let ceed = libceed::Ceed::init_with_error_handler( 2859df49d7eSJed Brown /// "/cpu/self/ref/serial", 2862ba8e59cSJeremy L Thompson /// libceed::ErrorHandler::ErrorAbort, 2879df49d7eSJed Brown /// ); 2889df49d7eSJed Brown /// ``` 2892ba8e59cSJeremy L Thompson pub fn init_with_error_handler(resource: &str, handler: ErrorHandler) -> Self { 2909df49d7eSJed Brown REGISTER.call_once(|| unsafe { 2919df49d7eSJed Brown bind_ceed::CeedRegisterAll(); 2929df49d7eSJed Brown bind_ceed::CeedQFunctionRegisterAll(); 2939df49d7eSJed Brown }); 2949df49d7eSJed Brown 2959df49d7eSJed Brown // Convert to C string 2969df49d7eSJed Brown let c_resource = CString::new(resource).expect("CString::new failed"); 2979df49d7eSJed Brown 2989df49d7eSJed Brown // Get error handler pointer 2999df49d7eSJed Brown let eh = match handler { 3002ba8e59cSJeremy L Thompson ErrorHandler::ErrorAbort => bind_ceed::CeedErrorAbort, 3012ba8e59cSJeremy L Thompson ErrorHandler::ErrorExit => bind_ceed::CeedErrorExit, 3022ba8e59cSJeremy L Thompson ErrorHandler::ErrorReturn => bind_ceed::CeedErrorReturn, 3032ba8e59cSJeremy L Thompson ErrorHandler::ErrorStore => bind_ceed::CeedErrorStore, 3049df49d7eSJed Brown }; 3059df49d7eSJed Brown 3069df49d7eSJed Brown // Call to libCEED 3079df49d7eSJed Brown let mut ptr = std::ptr::null_mut(); 3089df49d7eSJed Brown let mut ierr = unsafe { bind_ceed::CeedInit(c_resource.as_ptr() as *const i8, &mut ptr) }; 3099df49d7eSJed Brown if ierr != 0 { 3109df49d7eSJed Brown panic!("Error initializing backend resource: {}", resource) 3119df49d7eSJed Brown } 3129df49d7eSJed Brown ierr = unsafe { bind_ceed::CeedSetErrorHandler(ptr, Some(eh)) }; 3139df49d7eSJed Brown let ceed = Ceed { ptr }; 3149df49d7eSJed Brown ceed.check_error(ierr).unwrap(); 3159df49d7eSJed Brown ceed 3169df49d7eSJed Brown } 3179df49d7eSJed Brown 3189df49d7eSJed Brown /// Default initializer for testing 3199df49d7eSJed Brown #[doc(hidden)] 3209df49d7eSJed Brown pub fn default_init() -> Self { 3219df49d7eSJed Brown // Convert to C string 3229df49d7eSJed Brown let resource = "/cpu/self/ref/serial"; 3239df49d7eSJed Brown crate::Ceed::init(resource) 3249df49d7eSJed Brown } 3259df49d7eSJed Brown 3269df49d7eSJed Brown /// Internal error checker 3279df49d7eSJed Brown #[doc(hidden)] 3289df49d7eSJed Brown fn check_error(&self, ierr: i32) -> Result<i32> { 3299df49d7eSJed Brown // Return early if code is clean 3309df49d7eSJed Brown if ierr == bind_ceed::CeedErrorType_CEED_ERROR_SUCCESS { 3319df49d7eSJed Brown return Ok(ierr); 3329df49d7eSJed Brown } 3339df49d7eSJed Brown // Retrieve error message 3349df49d7eSJed Brown let mut ptr: *const std::os::raw::c_char = std::ptr::null_mut(); 3359df49d7eSJed Brown let c_str = unsafe { 3369df49d7eSJed Brown bind_ceed::CeedGetErrorMessage(self.ptr, &mut ptr); 3379df49d7eSJed Brown std::ffi::CStr::from_ptr(ptr) 3389df49d7eSJed Brown }; 3399df49d7eSJed Brown let message = c_str.to_string_lossy().to_string(); 3409df49d7eSJed Brown // Panic if negative code, otherwise return error 3419df49d7eSJed Brown if ierr < bind_ceed::CeedErrorType_CEED_ERROR_SUCCESS { 3429df49d7eSJed Brown panic!("{}", message); 3439df49d7eSJed Brown } 3442ba8e59cSJeremy L Thompson Err(Error { message }) 3459df49d7eSJed Brown } 3469df49d7eSJed Brown 3479df49d7eSJed Brown /// Returns full resource name for a Ceed context 3489df49d7eSJed Brown /// 3499df49d7eSJed Brown /// ``` 3509df49d7eSJed Brown /// let ceed = libceed::Ceed::init("/cpu/self/ref/serial"); 3519df49d7eSJed Brown /// let resource = ceed.resource(); 3529df49d7eSJed Brown /// 3539df49d7eSJed Brown /// assert_eq!(resource, "/cpu/self/ref/serial".to_string()) 3549df49d7eSJed Brown /// ``` 3559df49d7eSJed Brown pub fn resource(&self) -> String { 3569df49d7eSJed Brown let mut ptr: *const std::os::raw::c_char = std::ptr::null_mut(); 3579df49d7eSJed Brown let c_str = unsafe { 3589df49d7eSJed Brown bind_ceed::CeedGetResource(self.ptr, &mut ptr); 3599df49d7eSJed Brown std::ffi::CStr::from_ptr(ptr) 3609df49d7eSJed Brown }; 3619df49d7eSJed Brown c_str.to_string_lossy().to_string() 3629df49d7eSJed Brown } 3639df49d7eSJed Brown 3647ed177dbSJed Brown /// Returns a Vector of the specified length (does not allocate memory) 3659df49d7eSJed Brown /// 3669df49d7eSJed Brown /// # arguments 3679df49d7eSJed Brown /// 3689df49d7eSJed Brown /// * `n` - Length of vector 3699df49d7eSJed Brown /// 3709df49d7eSJed Brown /// ``` 3719df49d7eSJed Brown /// # use libceed::prelude::*; 3724d27c890SJeremy L Thompson /// # fn main() -> libceed::Result<()> { 3739df49d7eSJed Brown /// # let ceed = libceed::Ceed::default_init(); 374c68be7a2SJeremy L Thompson /// let vec = ceed.vector(10)?; 375c68be7a2SJeremy L Thompson /// # Ok(()) 376c68be7a2SJeremy L Thompson /// # } 3779df49d7eSJed Brown /// ``` 378594ef120SJeremy L Thompson pub fn vector<'a>(&self, n: usize) -> Result<Vector<'a>> { 3799df49d7eSJed Brown Vector::create(self, n) 3809df49d7eSJed Brown } 3819df49d7eSJed Brown 3829df49d7eSJed Brown /// Create a Vector initialized with the data (copied) from a slice 3839df49d7eSJed Brown /// 3849df49d7eSJed Brown /// # arguments 3859df49d7eSJed Brown /// 3869df49d7eSJed Brown /// * `slice` - Slice containing data 3879df49d7eSJed Brown /// 3889df49d7eSJed Brown /// ``` 3899df49d7eSJed Brown /// # use libceed::prelude::*; 3904d27c890SJeremy L Thompson /// # fn main() -> libceed::Result<()> { 3919df49d7eSJed Brown /// # let ceed = libceed::Ceed::default_init(); 392c68be7a2SJeremy L Thompson /// let vec = ceed.vector_from_slice(&[1., 2., 3.])?; 3939df49d7eSJed Brown /// assert_eq!(vec.length(), 3); 394c68be7a2SJeremy L Thompson /// # Ok(()) 395c68be7a2SJeremy L Thompson /// # } 3969df49d7eSJed Brown /// ``` 397594ef120SJeremy L Thompson pub fn vector_from_slice<'a>(&self, slice: &[crate::Scalar]) -> Result<Vector<'a>> { 3989df49d7eSJed Brown Vector::from_slice(self, slice) 3999df49d7eSJed Brown } 4009df49d7eSJed Brown 4017ed177dbSJed Brown /// Returns an ElemRestriction, $\mathcal{E}$, which extracts the degrees of 4027ed177dbSJed Brown /// freedom for each element from the local vector into the element vector 4037ed177dbSJed Brown /// or assembles contributions from each element in the element vector to 4047ed177dbSJed Brown /// the local vector 4059df49d7eSJed Brown /// 4069df49d7eSJed Brown /// # arguments 4079df49d7eSJed Brown /// 4089df49d7eSJed Brown /// * `nelem` - Number of elements described in the offsets array 4099df49d7eSJed Brown /// * `elemsize` - Size (number of "nodes") per element 4109df49d7eSJed Brown /// * `ncomp` - Number of field components per interpolation node (1 4119df49d7eSJed Brown /// for scalar fields) 4129df49d7eSJed Brown /// * `compstride` - Stride between components for the same Lvector "node". 4139df49d7eSJed Brown /// Data for node `i`, component `j`, element `k` can be 4149df49d7eSJed Brown /// found in the Lvector at index 4159df49d7eSJed Brown /// `offsets[i + k*elemsize] + j*compstride`. 4169df49d7eSJed Brown /// * `lsize` - The size of the Lvector. This vector may be larger 4179df49d7eSJed Brown /// than the elements and fields given by this 4189df49d7eSJed Brown /// restriction. 4199df49d7eSJed Brown /// * `mtype` - Memory type of the offsets array, see CeedMemType 4209df49d7eSJed Brown /// * `offsets` - Array of shape `[nelem, elemsize]`. Row `i` holds the 4219df49d7eSJed Brown /// ordered list of the offsets (into the input CeedVector) 4229df49d7eSJed Brown /// for the unknowns corresponding to element `i`, where 4239df49d7eSJed Brown /// `0 <= i < nelem`. All offsets must be in the range 4249df49d7eSJed Brown /// `[0, lsize - 1]`. 4259df49d7eSJed Brown /// 4269df49d7eSJed Brown /// ``` 4279df49d7eSJed Brown /// # use libceed::prelude::*; 4284d27c890SJeremy L Thompson /// # fn main() -> libceed::Result<()> { 4299df49d7eSJed Brown /// # let ceed = libceed::Ceed::default_init(); 4309df49d7eSJed Brown /// let nelem = 3; 4319df49d7eSJed Brown /// let mut ind: Vec<i32> = vec![0; 2 * nelem]; 4329df49d7eSJed Brown /// for i in 0..nelem { 4339df49d7eSJed Brown /// ind[2 * i + 0] = i as i32; 4349df49d7eSJed Brown /// ind[2 * i + 1] = (i + 1) as i32; 4359df49d7eSJed Brown /// } 436c68be7a2SJeremy L Thompson /// let r = ceed.elem_restriction(nelem, 2, 1, 1, nelem + 1, MemType::Host, &ind)?; 437c68be7a2SJeremy L Thompson /// # Ok(()) 438c68be7a2SJeremy L Thompson /// # } 4399df49d7eSJed Brown /// ``` 440594ef120SJeremy L Thompson pub fn elem_restriction<'a>( 4419df49d7eSJed Brown &self, 4429df49d7eSJed Brown nelem: usize, 4439df49d7eSJed Brown elemsize: usize, 4449df49d7eSJed Brown ncomp: usize, 4459df49d7eSJed Brown compstride: usize, 4469df49d7eSJed Brown lsize: usize, 4479df49d7eSJed Brown mtype: MemType, 4489df49d7eSJed Brown offsets: &[i32], 449594ef120SJeremy L Thompson ) -> Result<ElemRestriction<'a>> { 4509df49d7eSJed Brown ElemRestriction::create( 4519df49d7eSJed Brown self, nelem, elemsize, ncomp, compstride, lsize, mtype, offsets, 4529df49d7eSJed Brown ) 4539df49d7eSJed Brown } 4549df49d7eSJed Brown 4557ed177dbSJed Brown /// Returns an ElemRestriction, $\mathcal{E}$, from an local vector to 4567ed177dbSJed Brown /// an element vector where data can be indexed from the `strides` array 4579df49d7eSJed Brown /// 4589df49d7eSJed Brown /// # arguments 4599df49d7eSJed Brown /// 4609df49d7eSJed Brown /// * `nelem` - Number of elements described in the offsets array 4619df49d7eSJed Brown /// * `elemsize` - Size (number of "nodes") per element 4629df49d7eSJed Brown /// * `ncomp` - Number of field components per interpolation node (1 4639df49d7eSJed Brown /// for scalar fields) 4649df49d7eSJed Brown /// * `lsize` - The size of the Lvector. This vector may be larger 4659df49d7eSJed Brown /// than the elements and fields given by this restriction. 4669df49d7eSJed Brown /// * `strides` - Array for strides between `[nodes, components, elements]`. 4679df49d7eSJed Brown /// Data for node `i`, component `j`, element `k` can be 4689df49d7eSJed Brown /// found in the Lvector at index 4699df49d7eSJed Brown /// `i*strides[0] + j*strides[1] + k*strides[2]`. 4709df49d7eSJed Brown /// CEED_STRIDES_BACKEND may be used with vectors created 4719df49d7eSJed Brown /// by a Ceed backend. 4729df49d7eSJed Brown /// 4739df49d7eSJed Brown /// ``` 4749df49d7eSJed Brown /// # use libceed::prelude::*; 4754d27c890SJeremy L Thompson /// # fn main() -> libceed::Result<()> { 4769df49d7eSJed Brown /// # let ceed = libceed::Ceed::default_init(); 4779df49d7eSJed Brown /// let nelem = 3; 4789df49d7eSJed Brown /// let strides: [i32; 3] = [1, 2, 2]; 479c68be7a2SJeremy L Thompson /// let r = ceed.strided_elem_restriction(nelem, 2, 1, nelem * 2, strides)?; 480c68be7a2SJeremy L Thompson /// # Ok(()) 481c68be7a2SJeremy L Thompson /// # } 4829df49d7eSJed Brown /// ``` 483594ef120SJeremy L Thompson pub fn strided_elem_restriction<'a>( 4849df49d7eSJed Brown &self, 4859df49d7eSJed Brown nelem: usize, 4869df49d7eSJed Brown elemsize: usize, 4879df49d7eSJed Brown ncomp: usize, 4889df49d7eSJed Brown lsize: usize, 4899df49d7eSJed Brown strides: [i32; 3], 490594ef120SJeremy L Thompson ) -> Result<ElemRestriction<'a>> { 4919df49d7eSJed Brown ElemRestriction::create_strided(self, nelem, elemsize, ncomp, lsize, strides) 4929df49d7eSJed Brown } 4939df49d7eSJed Brown 4947ed177dbSJed Brown /// Returns an $H^1$ tensor-product Basis 4959df49d7eSJed Brown /// 4969df49d7eSJed Brown /// # arguments 4979df49d7eSJed Brown /// 4989df49d7eSJed Brown /// * `dim` - Topological dimension of element 4999df49d7eSJed Brown /// * `ncomp` - Number of field components (1 for scalar fields) 5009df49d7eSJed Brown /// * `P1d` - Number of Gauss-Lobatto nodes in one dimension. The 5019df49d7eSJed Brown /// polynomial degree of the resulting `Q_k` element is 5029df49d7eSJed Brown /// `k=P-1`. 5039df49d7eSJed Brown /// * `Q1d` - Number of quadrature points in one dimension 5049df49d7eSJed Brown /// * `interp1d` - Row-major `(Q1d * P1d)` matrix expressing the values of 5059df49d7eSJed Brown /// nodal basis functions at quadrature points 5069df49d7eSJed Brown /// * `grad1d` - Row-major `(Q1d * P1d)` matrix expressing derivatives of 5079df49d7eSJed Brown /// nodal basis functions at quadrature points 5089df49d7eSJed Brown /// * `qref1d` - Array of length `Q1d` holding the locations of quadrature 5099df49d7eSJed Brown /// points on the 1D reference element `[-1, 1]` 5109df49d7eSJed Brown /// * `qweight1d` - Array of length `Q1d` holding the quadrature weights on 5119df49d7eSJed Brown /// the reference element 5129df49d7eSJed Brown /// 5139df49d7eSJed Brown /// ``` 5149df49d7eSJed Brown /// # use libceed::prelude::*; 5154d27c890SJeremy L Thompson /// # fn main() -> libceed::Result<()> { 5169df49d7eSJed Brown /// # let ceed = libceed::Ceed::default_init(); 5179df49d7eSJed Brown /// let interp1d = [ 0.62994317, 0.47255875, -0.14950343, 0.04700152, 5189df49d7eSJed Brown /// -0.07069480, 0.97297619, 0.13253993, -0.03482132, 5199df49d7eSJed Brown /// -0.03482132, 0.13253993, 0.97297619, -0.07069480, 5209df49d7eSJed Brown /// 0.04700152, -0.14950343, 0.47255875, 0.62994317]; 5219df49d7eSJed Brown /// let grad1d = [-2.34183742, 2.78794489, -0.63510411, 0.18899664, 5229df49d7eSJed Brown /// -0.51670214, -0.48795249, 1.33790510, -0.33325047, 5239df49d7eSJed Brown // 0.33325047, -1.33790510, 0.48795249, 0.51670214, 5249df49d7eSJed Brown /// -0.18899664, 0.63510411, -2.78794489, 2.34183742]; 5259df49d7eSJed Brown /// let qref1d = [-0.86113631, -0.33998104, 0.33998104, 0.86113631]; 5269df49d7eSJed Brown /// let qweight1d = [ 0.34785485, 0.65214515, 0.65214515, 0.34785485]; 5279df49d7eSJed Brown /// let b = ceed. 528c68be7a2SJeremy L Thompson /// basis_tensor_H1(2, 1, 4, 4, &interp1d, &grad1d, &qref1d, &qweight1d)?; 529c68be7a2SJeremy L Thompson /// # Ok(()) 530c68be7a2SJeremy L Thompson /// # } 5319df49d7eSJed Brown /// ``` 532594ef120SJeremy L Thompson pub fn basis_tensor_H1<'a>( 5339df49d7eSJed Brown &self, 5349df49d7eSJed Brown dim: usize, 5359df49d7eSJed Brown ncomp: usize, 5369df49d7eSJed Brown P1d: usize, 5379df49d7eSJed Brown Q1d: usize, 53880a9ef05SNatalie Beams interp1d: &[crate::Scalar], 53980a9ef05SNatalie Beams grad1d: &[crate::Scalar], 54080a9ef05SNatalie Beams qref1d: &[crate::Scalar], 54180a9ef05SNatalie Beams qweight1d: &[crate::Scalar], 542594ef120SJeremy L Thompson ) -> Result<Basis<'a>> { 5439df49d7eSJed Brown Basis::create_tensor_H1( 5449df49d7eSJed Brown self, dim, ncomp, P1d, Q1d, interp1d, grad1d, qref1d, qweight1d, 5459df49d7eSJed Brown ) 5469df49d7eSJed Brown } 5479df49d7eSJed Brown 5487ed177dbSJed Brown /// Returns an $H^1$ Lagrange tensor-product Basis 5499df49d7eSJed Brown /// 5509df49d7eSJed Brown /// # arguments 5519df49d7eSJed Brown /// 5529df49d7eSJed Brown /// * `dim` - Topological dimension of element 5539df49d7eSJed Brown /// * `ncomp` - Number of field components (1 for scalar fields) 5549df49d7eSJed Brown /// * `P` - Number of Gauss-Lobatto nodes in one dimension. The 5559df49d7eSJed Brown /// polynomial degree of the resulting `Q_k` element is `k=P-1`. 5569df49d7eSJed Brown /// * `Q` - Number of quadrature points in one dimension 5579df49d7eSJed Brown /// * `qmode` - Distribution of the `Q` quadrature points (affects order of 5589df49d7eSJed Brown /// accuracy for the quadrature) 5599df49d7eSJed Brown /// 5609df49d7eSJed Brown /// ``` 5619df49d7eSJed Brown /// # use libceed::prelude::*; 5624d27c890SJeremy L Thompson /// # fn main() -> libceed::Result<()> { 5639df49d7eSJed Brown /// # let ceed = libceed::Ceed::default_init(); 564c68be7a2SJeremy L Thompson /// let b = ceed.basis_tensor_H1_Lagrange(2, 1, 3, 4, QuadMode::Gauss)?; 565c68be7a2SJeremy L Thompson /// # Ok(()) 566c68be7a2SJeremy L Thompson /// # } 5679df49d7eSJed Brown /// ``` 568594ef120SJeremy L Thompson pub fn basis_tensor_H1_Lagrange<'a>( 5699df49d7eSJed Brown &self, 5709df49d7eSJed Brown dim: usize, 5719df49d7eSJed Brown ncomp: usize, 5729df49d7eSJed Brown P: usize, 5739df49d7eSJed Brown Q: usize, 5749df49d7eSJed Brown qmode: QuadMode, 575594ef120SJeremy L Thompson ) -> Result<Basis<'a>> { 5769df49d7eSJed Brown Basis::create_tensor_H1_Lagrange(self, dim, ncomp, P, Q, qmode) 5779df49d7eSJed Brown } 5789df49d7eSJed Brown 5797ed177dbSJed Brown /// Returns an $H-1$ Basis 5809df49d7eSJed Brown /// 5819df49d7eSJed Brown /// # arguments 5829df49d7eSJed Brown /// 5839df49d7eSJed Brown /// * `topo` - Topology of element, e.g. hypercube, simplex, ect 5849df49d7eSJed Brown /// * `ncomp` - Number of field components (1 for scalar fields) 5859df49d7eSJed Brown /// * `nnodes` - Total number of nodes 5869df49d7eSJed Brown /// * `nqpts` - Total number of quadrature points 5879df49d7eSJed Brown /// * `interp` - Row-major `(nqpts * nnodes)` matrix expressing the values of 5889df49d7eSJed Brown /// nodal basis functions at quadrature points 589*97c1c57aSSebastian Grimberg /// * `grad` - Row-major `(dim * nqpts * nnodes)` matrix expressing 5909df49d7eSJed Brown /// derivatives of nodal basis functions at quadrature points 5919df49d7eSJed Brown /// * `qref` - Array of length `nqpts` holding the locations of quadrature 5929df49d7eSJed Brown /// points on the reference element `[-1, 1]` 5939df49d7eSJed Brown /// * `qweight` - Array of length `nqpts` holding the quadrature weights on 5949df49d7eSJed Brown /// the reference element 5959df49d7eSJed Brown /// 5969df49d7eSJed Brown /// ``` 5979df49d7eSJed Brown /// # use libceed::prelude::*; 5984d27c890SJeremy L Thompson /// # fn main() -> libceed::Result<()> { 5999df49d7eSJed Brown /// # let ceed = libceed::Ceed::default_init(); 6009df49d7eSJed Brown /// let interp = [ 6019df49d7eSJed Brown /// 0.12000000, 6029df49d7eSJed Brown /// 0.48000000, 6039df49d7eSJed Brown /// -0.12000000, 6049df49d7eSJed Brown /// 0.48000000, 6059df49d7eSJed Brown /// 0.16000000, 6069df49d7eSJed Brown /// -0.12000000, 6079df49d7eSJed Brown /// -0.12000000, 6089df49d7eSJed Brown /// 0.48000000, 6099df49d7eSJed Brown /// 0.12000000, 6109df49d7eSJed Brown /// 0.16000000, 6119df49d7eSJed Brown /// 0.48000000, 6129df49d7eSJed Brown /// -0.12000000, 6139df49d7eSJed Brown /// -0.11111111, 6149df49d7eSJed Brown /// 0.44444444, 6159df49d7eSJed Brown /// -0.11111111, 6169df49d7eSJed Brown /// 0.44444444, 6179df49d7eSJed Brown /// 0.44444444, 6189df49d7eSJed Brown /// -0.11111111, 6199df49d7eSJed Brown /// -0.12000000, 6209df49d7eSJed Brown /// 0.16000000, 6219df49d7eSJed Brown /// -0.12000000, 6229df49d7eSJed Brown /// 0.48000000, 6239df49d7eSJed Brown /// 0.48000000, 6249df49d7eSJed Brown /// 0.12000000, 6259df49d7eSJed Brown /// ]; 6269df49d7eSJed Brown /// let grad = [ 6279df49d7eSJed Brown /// -1.40000000, 6289df49d7eSJed Brown /// 1.60000000, 6299df49d7eSJed Brown /// -0.20000000, 6309df49d7eSJed Brown /// -0.80000000, 6319df49d7eSJed Brown /// 0.80000000, 6329df49d7eSJed Brown /// 0.00000000, 6339df49d7eSJed Brown /// 0.20000000, 6349df49d7eSJed Brown /// -1.60000000, 6359df49d7eSJed Brown /// 1.40000000, 6369df49d7eSJed Brown /// -0.80000000, 6379df49d7eSJed Brown /// 0.80000000, 6389df49d7eSJed Brown /// 0.00000000, 6399df49d7eSJed Brown /// -0.33333333, 6409df49d7eSJed Brown /// 0.00000000, 6419df49d7eSJed Brown /// 0.33333333, 6429df49d7eSJed Brown /// -1.33333333, 6439df49d7eSJed Brown /// 1.33333333, 6449df49d7eSJed Brown /// 0.00000000, 6459df49d7eSJed Brown /// 0.20000000, 6469df49d7eSJed Brown /// 0.00000000, 6479df49d7eSJed Brown /// -0.20000000, 6489df49d7eSJed Brown /// -2.40000000, 6499df49d7eSJed Brown /// 2.40000000, 6509df49d7eSJed Brown /// 0.00000000, 6519df49d7eSJed Brown /// -1.40000000, 6529df49d7eSJed Brown /// -0.80000000, 6539df49d7eSJed Brown /// 0.00000000, 6549df49d7eSJed Brown /// 1.60000000, 6559df49d7eSJed Brown /// 0.80000000, 6569df49d7eSJed Brown /// -0.20000000, 6579df49d7eSJed Brown /// 0.20000000, 6589df49d7eSJed Brown /// -2.40000000, 6599df49d7eSJed Brown /// 0.00000000, 6609df49d7eSJed Brown /// 0.00000000, 6619df49d7eSJed Brown /// 2.40000000, 6629df49d7eSJed Brown /// -0.20000000, 6639df49d7eSJed Brown /// -0.33333333, 6649df49d7eSJed Brown /// -1.33333333, 6659df49d7eSJed Brown /// 0.00000000, 6669df49d7eSJed Brown /// 0.00000000, 6679df49d7eSJed Brown /// 1.33333333, 6689df49d7eSJed Brown /// 0.33333333, 6699df49d7eSJed Brown /// 0.20000000, 6709df49d7eSJed Brown /// -0.80000000, 6719df49d7eSJed Brown /// 0.00000000, 6729df49d7eSJed Brown /// -1.60000000, 6739df49d7eSJed Brown /// 0.80000000, 6749df49d7eSJed Brown /// 1.40000000, 6759df49d7eSJed Brown /// ]; 6769df49d7eSJed Brown /// let qref = [ 6779df49d7eSJed Brown /// 0.20000000, 0.60000000, 0.33333333, 0.20000000, 0.20000000, 0.20000000, 0.33333333, 6789df49d7eSJed Brown /// 0.60000000, 6799df49d7eSJed Brown /// ]; 6809df49d7eSJed Brown /// let qweight = [0.26041667, 0.26041667, -0.28125000, 0.26041667]; 681c68be7a2SJeremy L Thompson /// let b = ceed.basis_H1( 6829df49d7eSJed Brown /// ElemTopology::Triangle, 6839df49d7eSJed Brown /// 1, 6849df49d7eSJed Brown /// 6, 6859df49d7eSJed Brown /// 4, 6869df49d7eSJed Brown /// &interp, 6879df49d7eSJed Brown /// &grad, 6889df49d7eSJed Brown /// &qref, 6899df49d7eSJed Brown /// &qweight, 690c68be7a2SJeremy L Thompson /// )?; 691c68be7a2SJeremy L Thompson /// # Ok(()) 692c68be7a2SJeremy L Thompson /// # } 6939df49d7eSJed Brown /// ``` 694594ef120SJeremy L Thompson pub fn basis_H1<'a>( 6959df49d7eSJed Brown &self, 6969df49d7eSJed Brown topo: ElemTopology, 6979df49d7eSJed Brown ncomp: usize, 6989df49d7eSJed Brown nnodes: usize, 6999df49d7eSJed Brown nqpts: usize, 70080a9ef05SNatalie Beams interp: &[crate::Scalar], 70180a9ef05SNatalie Beams grad: &[crate::Scalar], 70280a9ef05SNatalie Beams qref: &[crate::Scalar], 70380a9ef05SNatalie Beams qweight: &[crate::Scalar], 704594ef120SJeremy L Thompson ) -> Result<Basis<'a>> { 7059df49d7eSJed Brown Basis::create_H1( 7069df49d7eSJed Brown self, topo, ncomp, nnodes, nqpts, interp, grad, qref, qweight, 7079df49d7eSJed Brown ) 7089df49d7eSJed Brown } 7099df49d7eSJed Brown 710*97c1c57aSSebastian Grimberg /// Returns an $H(div)$ Basis 711*97c1c57aSSebastian Grimberg /// 712*97c1c57aSSebastian Grimberg /// # arguments 713*97c1c57aSSebastian Grimberg /// 714*97c1c57aSSebastian Grimberg /// * `topo` - Topology of element, e.g. hypercube, simplex, ect 715*97c1c57aSSebastian Grimberg /// * `ncomp` - Number of field components (1 for scalar fields) 716*97c1c57aSSebastian Grimberg /// * `nnodes` - Total number of nodes 717*97c1c57aSSebastian Grimberg /// * `nqpts` - Total number of quadrature points 718*97c1c57aSSebastian Grimberg /// * `interp` - Row-major `(dim * nqpts * nnodes)` matrix expressing the 719*97c1c57aSSebastian Grimberg /// values of basis functions at quadrature points 720*97c1c57aSSebastian Grimberg /// * `div` - Row-major `(nqpts * nnodes)` matrix expressing the 721*97c1c57aSSebastian Grimberg /// divergence of basis functions at quadrature points 722*97c1c57aSSebastian Grimberg /// * `qref` - Array of length `nqpts` holding the locations of quadrature 723*97c1c57aSSebastian Grimberg /// points on the reference element `[-1, 1]` 724*97c1c57aSSebastian Grimberg /// * `qweight` - Array of length `nqpts` holding the quadrature weights on 725*97c1c57aSSebastian Grimberg /// the reference element 726*97c1c57aSSebastian Grimberg /// 727*97c1c57aSSebastian Grimberg /// ``` 728*97c1c57aSSebastian Grimberg /// # use libceed::prelude::*; 729*97c1c57aSSebastian Grimberg /// # fn main() -> libceed::Result<()> { 730*97c1c57aSSebastian Grimberg /// # let ceed = libceed::Ceed::default_init(); 731*97c1c57aSSebastian Grimberg /// let interp = [ 732*97c1c57aSSebastian Grimberg /// 0.00000000, 733*97c1c57aSSebastian Grimberg /// -1.57735027, 734*97c1c57aSSebastian Grimberg /// 0.57735027, 735*97c1c57aSSebastian Grimberg /// 0.00000000, 736*97c1c57aSSebastian Grimberg /// 0.00000000, 737*97c1c57aSSebastian Grimberg /// -1.57735027, 738*97c1c57aSSebastian Grimberg /// 0.57735027, 739*97c1c57aSSebastian Grimberg /// 0.00000000, 740*97c1c57aSSebastian Grimberg /// 0.00000000, 741*97c1c57aSSebastian Grimberg /// -1.57735027, 742*97c1c57aSSebastian Grimberg /// 0.57735027, 743*97c1c57aSSebastian Grimberg /// 0.00000000, 744*97c1c57aSSebastian Grimberg /// 0.00000000, 745*97c1c57aSSebastian Grimberg /// -0.42264973, 746*97c1c57aSSebastian Grimberg /// -0.57735027, 747*97c1c57aSSebastian Grimberg /// 0.00000000, 748*97c1c57aSSebastian Grimberg /// 0.42264973, 749*97c1c57aSSebastian Grimberg /// 0.00000000, 750*97c1c57aSSebastian Grimberg /// 0.00000000, 751*97c1c57aSSebastian Grimberg /// 0.57735027, 752*97c1c57aSSebastian Grimberg /// 0.42264973, 753*97c1c57aSSebastian Grimberg /// 0.00000000, 754*97c1c57aSSebastian Grimberg /// 0.00000000, 755*97c1c57aSSebastian Grimberg /// 0.57735027, 756*97c1c57aSSebastian Grimberg /// 1.57735027, 757*97c1c57aSSebastian Grimberg /// 0.00000000, 758*97c1c57aSSebastian Grimberg /// 0.00000000, 759*97c1c57aSSebastian Grimberg /// -0.57735027, 760*97c1c57aSSebastian Grimberg /// 1.57735027, 761*97c1c57aSSebastian Grimberg /// 0.00000000, 762*97c1c57aSSebastian Grimberg /// 0.00000000, 763*97c1c57aSSebastian Grimberg /// -0.57735027, 764*97c1c57aSSebastian Grimberg /// ]; 765*97c1c57aSSebastian Grimberg /// let div = [ 766*97c1c57aSSebastian Grimberg /// -1.00000000, 767*97c1c57aSSebastian Grimberg /// 1.00000000, 768*97c1c57aSSebastian Grimberg /// -1.00000000, 769*97c1c57aSSebastian Grimberg /// 1.00000000, 770*97c1c57aSSebastian Grimberg /// -1.00000000, 771*97c1c57aSSebastian Grimberg /// 1.00000000, 772*97c1c57aSSebastian Grimberg /// -1.00000000, 773*97c1c57aSSebastian Grimberg /// 1.00000000, 774*97c1c57aSSebastian Grimberg /// -1.00000000, 775*97c1c57aSSebastian Grimberg /// 1.00000000, 776*97c1c57aSSebastian Grimberg /// -1.00000000, 777*97c1c57aSSebastian Grimberg /// 1.00000000, 778*97c1c57aSSebastian Grimberg /// -1.00000000, 779*97c1c57aSSebastian Grimberg /// 1.00000000, 780*97c1c57aSSebastian Grimberg /// -1.00000000, 781*97c1c57aSSebastian Grimberg /// 1.00000000, 782*97c1c57aSSebastian Grimberg /// ]; 783*97c1c57aSSebastian Grimberg /// let qref = [ 784*97c1c57aSSebastian Grimberg /// 0.57735026, 0.57735026, 0.57735026, 0.57735026, 0.57735026, 0.57735026, 0.57735026, 785*97c1c57aSSebastian Grimberg /// 0.57735026, 786*97c1c57aSSebastian Grimberg /// ]; 787*97c1c57aSSebastian Grimberg /// let qweight = [1.00000000, 1.00000000, 1.00000000, 1.00000000]; 788*97c1c57aSSebastian Grimberg /// let b = ceed.basis_Hdiv(ElemTopology::Quad, 1, 4, 4, &interp, &div, &qref, &qweight)?; 789*97c1c57aSSebastian Grimberg /// # Ok(()) 790*97c1c57aSSebastian Grimberg /// # } 791*97c1c57aSSebastian Grimberg /// ``` 792*97c1c57aSSebastian Grimberg pub fn basis_Hdiv<'a>( 793*97c1c57aSSebastian Grimberg &self, 794*97c1c57aSSebastian Grimberg topo: ElemTopology, 795*97c1c57aSSebastian Grimberg ncomp: usize, 796*97c1c57aSSebastian Grimberg nnodes: usize, 797*97c1c57aSSebastian Grimberg nqpts: usize, 798*97c1c57aSSebastian Grimberg interp: &[crate::Scalar], 799*97c1c57aSSebastian Grimberg div: &[crate::Scalar], 800*97c1c57aSSebastian Grimberg qref: &[crate::Scalar], 801*97c1c57aSSebastian Grimberg qweight: &[crate::Scalar], 802*97c1c57aSSebastian Grimberg ) -> Result<Basis<'a>> { 803*97c1c57aSSebastian Grimberg Basis::create_Hdiv(self, topo, ncomp, nnodes, nqpts, interp, div, qref, qweight) 804*97c1c57aSSebastian Grimberg } 805*97c1c57aSSebastian Grimberg 806*97c1c57aSSebastian Grimberg /// Returns an $H(curl)$ Basis 807*97c1c57aSSebastian Grimberg /// 808*97c1c57aSSebastian Grimberg /// # arguments 809*97c1c57aSSebastian Grimberg /// 810*97c1c57aSSebastian Grimberg /// * `topo` - Topology of element, e.g. hypercube, simplex, ect 811*97c1c57aSSebastian Grimberg /// * `ncomp` - Number of field components (1 for scalar fields) 812*97c1c57aSSebastian Grimberg /// * `nnodes` - Total number of nodes 813*97c1c57aSSebastian Grimberg /// * `nqpts` - Total number of quadrature points 814*97c1c57aSSebastian Grimberg /// * `interp` - Row-major `(dim * nqpts * nnodes)` matrix expressing the 815*97c1c57aSSebastian Grimberg /// values of basis functions at quadrature points 816*97c1c57aSSebastian Grimberg /// * `curl` - Row-major `(curl_comp * nqpts * nnodes)`, `curl_comp = 1 if 817*97c1c57aSSebastian Grimberg /// dim < 3 else dim` matrix expressing the curl of basis 818*97c1c57aSSebastian Grimberg /// functions at quadrature points 819*97c1c57aSSebastian Grimberg /// * `qref` - Array of length `nqpts` holding the locations of quadrature 820*97c1c57aSSebastian Grimberg /// points on the reference element `[-1, 1]` 821*97c1c57aSSebastian Grimberg /// * `qweight` - Array of length `nqpts` holding the quadrature weights on 822*97c1c57aSSebastian Grimberg /// the reference element 823*97c1c57aSSebastian Grimberg /// 824*97c1c57aSSebastian Grimberg /// ``` 825*97c1c57aSSebastian Grimberg /// # use libceed::prelude::*; 826*97c1c57aSSebastian Grimberg /// # fn main() -> libceed::Result<()> { 827*97c1c57aSSebastian Grimberg /// # let ceed = libceed::Ceed::default_init(); 828*97c1c57aSSebastian Grimberg /// let interp = [ 829*97c1c57aSSebastian Grimberg /// -0.20000000, 830*97c1c57aSSebastian Grimberg /// 0.20000000, 831*97c1c57aSSebastian Grimberg /// 0.80000000, 832*97c1c57aSSebastian Grimberg /// -0.20000000, 833*97c1c57aSSebastian Grimberg /// 0.20000000, 834*97c1c57aSSebastian Grimberg /// 0.80000000, 835*97c1c57aSSebastian Grimberg /// -0.33333333, 836*97c1c57aSSebastian Grimberg /// 0.33333333, 837*97c1c57aSSebastian Grimberg /// 0.66666667, 838*97c1c57aSSebastian Grimberg /// -0.60000000, 839*97c1c57aSSebastian Grimberg /// 0.60000000, 840*97c1c57aSSebastian Grimberg /// 0.40000000, 841*97c1c57aSSebastian Grimberg /// 0.20000000, 842*97c1c57aSSebastian Grimberg /// 0.80000000, 843*97c1c57aSSebastian Grimberg /// 0.20000000, 844*97c1c57aSSebastian Grimberg /// 0.60000000, 845*97c1c57aSSebastian Grimberg /// 0.40000000, 846*97c1c57aSSebastian Grimberg /// 0.60000000, 847*97c1c57aSSebastian Grimberg /// 0.33333333, 848*97c1c57aSSebastian Grimberg /// 0.66666667, 849*97c1c57aSSebastian Grimberg /// 0.33333333, 850*97c1c57aSSebastian Grimberg /// 0.20000000, 851*97c1c57aSSebastian Grimberg /// 0.80000000, 852*97c1c57aSSebastian Grimberg /// 0.20000000, 853*97c1c57aSSebastian Grimberg /// ]; 854*97c1c57aSSebastian Grimberg /// let curl = [ 855*97c1c57aSSebastian Grimberg /// 2.00000000, 856*97c1c57aSSebastian Grimberg /// -2.00000000, 857*97c1c57aSSebastian Grimberg /// -2.00000000, 858*97c1c57aSSebastian Grimberg /// 2.00000000, 859*97c1c57aSSebastian Grimberg /// -2.00000000, 860*97c1c57aSSebastian Grimberg /// -2.00000000, 861*97c1c57aSSebastian Grimberg /// 2.00000000, 862*97c1c57aSSebastian Grimberg /// -2.00000000, 863*97c1c57aSSebastian Grimberg /// -2.00000000, 864*97c1c57aSSebastian Grimberg /// 2.00000000, 865*97c1c57aSSebastian Grimberg /// -2.00000000, 866*97c1c57aSSebastian Grimberg /// -2.00000000, 867*97c1c57aSSebastian Grimberg /// ]; 868*97c1c57aSSebastian Grimberg /// let qref = [ 869*97c1c57aSSebastian Grimberg /// 0.20000000, 0.60000000, 0.33333333, 0.20000000, 0.20000000, 0.20000000, 0.33333333, 870*97c1c57aSSebastian Grimberg /// 0.60000000, 871*97c1c57aSSebastian Grimberg /// ]; 872*97c1c57aSSebastian Grimberg /// let qweight = [0.26041667, 0.26041667, -0.28125000, 0.26041667]; 873*97c1c57aSSebastian Grimberg /// let b = ceed.basis_Hcurl( 874*97c1c57aSSebastian Grimberg /// ElemTopology::Triangle, 875*97c1c57aSSebastian Grimberg /// 1, 876*97c1c57aSSebastian Grimberg /// 3, 877*97c1c57aSSebastian Grimberg /// 4, 878*97c1c57aSSebastian Grimberg /// &interp, 879*97c1c57aSSebastian Grimberg /// &curl, 880*97c1c57aSSebastian Grimberg /// &qref, 881*97c1c57aSSebastian Grimberg /// &qweight, 882*97c1c57aSSebastian Grimberg /// )?; 883*97c1c57aSSebastian Grimberg /// # Ok(()) 884*97c1c57aSSebastian Grimberg /// # } 885*97c1c57aSSebastian Grimberg /// ``` 886*97c1c57aSSebastian Grimberg pub fn basis_Hcurl<'a>( 887*97c1c57aSSebastian Grimberg &self, 888*97c1c57aSSebastian Grimberg topo: ElemTopology, 889*97c1c57aSSebastian Grimberg ncomp: usize, 890*97c1c57aSSebastian Grimberg nnodes: usize, 891*97c1c57aSSebastian Grimberg nqpts: usize, 892*97c1c57aSSebastian Grimberg interp: &[crate::Scalar], 893*97c1c57aSSebastian Grimberg curl: &[crate::Scalar], 894*97c1c57aSSebastian Grimberg qref: &[crate::Scalar], 895*97c1c57aSSebastian Grimberg qweight: &[crate::Scalar], 896*97c1c57aSSebastian Grimberg ) -> Result<Basis<'a>> { 897*97c1c57aSSebastian Grimberg Basis::create_Hcurl( 898*97c1c57aSSebastian Grimberg self, topo, ncomp, nnodes, nqpts, interp, curl, qref, qweight, 899*97c1c57aSSebastian Grimberg ) 900*97c1c57aSSebastian Grimberg } 901*97c1c57aSSebastian Grimberg 9027ed177dbSJed Brown /// Returns a QFunction for evaluating interior (volumetric) terms 9037ed177dbSJed Brown /// of a weak form corresponding to the $L^2$ inner product 9047ed177dbSJed Brown /// 9057ed177dbSJed Brown /// $$ 9067ed177dbSJed 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), 9077ed177dbSJed Brown /// $$ 9087ed177dbSJed Brown /// 9097ed177dbSJed Brown /// where $v \cdot f_0$ represents contraction over fields and $\nabla v : f_1$ 9107ed177dbSJed Brown /// represents contraction over both fields and spatial dimensions. 9119df49d7eSJed Brown /// 9129df49d7eSJed Brown /// # arguments 9139df49d7eSJed Brown /// 9149df49d7eSJed Brown /// * `vlength` - Vector length. Caller must ensure that number of 9159df49d7eSJed Brown /// quadrature points is a multiple of vlength. 9167ed177dbSJed Brown /// * `f` - Boxed closure to evaluate weak form at quadrature points. 9179df49d7eSJed Brown /// 9189df49d7eSJed Brown /// ``` 9199df49d7eSJed Brown /// # use libceed::prelude::*; 9204d27c890SJeremy L Thompson /// # fn main() -> libceed::Result<()> { 9219df49d7eSJed Brown /// # let ceed = libceed::Ceed::default_init(); 9229df49d7eSJed Brown /// let mut user_f = |[u, weights, ..]: QFunctionInputs, [v, ..]: QFunctionOutputs| { 9239df49d7eSJed Brown /// // Iterate over quadrature points 9249df49d7eSJed Brown /// v.iter_mut() 9259df49d7eSJed Brown /// .zip(u.iter().zip(weights.iter())) 9269df49d7eSJed Brown /// .for_each(|(v, (u, w))| *v = u * w); 9279df49d7eSJed Brown /// 9289df49d7eSJed Brown /// // Return clean error code 9299df49d7eSJed Brown /// 0 9309df49d7eSJed Brown /// }; 9319df49d7eSJed Brown /// 932c68be7a2SJeremy L Thompson /// let qf = ceed.q_function_interior(1, Box::new(user_f))?; 933c68be7a2SJeremy L Thompson /// # Ok(()) 934c68be7a2SJeremy L Thompson /// # } 9359df49d7eSJed Brown /// ``` 936594ef120SJeremy L Thompson pub fn q_function_interior<'a>( 9379df49d7eSJed Brown &self, 9389df49d7eSJed Brown vlength: usize, 9399df49d7eSJed Brown f: Box<qfunction::QFunctionUserClosure>, 940594ef120SJeremy L Thompson ) -> Result<QFunction<'a>> { 9419df49d7eSJed Brown QFunction::create(self, vlength, f) 9429df49d7eSJed Brown } 9439df49d7eSJed Brown 9447ed177dbSJed Brown /// Returns a QFunction for evaluating interior (volumetric) terms 9459df49d7eSJed Brown /// created by name 9469df49d7eSJed Brown /// 9477ed177dbSJed Brown /// # arguments 9487ed177dbSJed Brown /// 9497ed177dbSJed Brown /// * `name` - name of QFunction from libCEED gallery 9507ed177dbSJed Brown /// 9519df49d7eSJed Brown /// ``` 9529df49d7eSJed Brown /// # use libceed::prelude::*; 9534d27c890SJeremy L Thompson /// # fn main() -> libceed::Result<()> { 9549df49d7eSJed Brown /// # let ceed = libceed::Ceed::default_init(); 955c68be7a2SJeremy L Thompson /// let qf = ceed.q_function_interior_by_name("Mass1DBuild")?; 956c68be7a2SJeremy L Thompson /// # Ok(()) 957c68be7a2SJeremy L Thompson /// # } 9589df49d7eSJed Brown /// ``` 959594ef120SJeremy L Thompson pub fn q_function_interior_by_name<'a>(&self, name: &str) -> Result<QFunctionByName<'a>> { 9609df49d7eSJed Brown QFunctionByName::create(self, name) 9619df49d7eSJed Brown } 9629df49d7eSJed Brown 9637ed177dbSJed Brown /// Returns an Operator and associate a QFunction. A Basis and 9647ed177dbSJed Brown /// ElemRestriction can be associated with QFunction fields via 9659df49d7eSJed Brown /// set_field(). 9669df49d7eSJed Brown /// 9677ed177dbSJed Brown /// # arguments 9687ed177dbSJed Brown /// 9699df49d7eSJed Brown /// * `qf` - QFunction defining the action of the operator at quadrature 9709df49d7eSJed Brown /// points 9719df49d7eSJed Brown /// * `dqf` - QFunction defining the action of the Jacobian of the qf (or 9729df49d7eSJed Brown /// qfunction_none) 9739df49d7eSJed Brown /// * `dqfT` - QFunction defining the action of the transpose of the 9749df49d7eSJed Brown /// Jacobian of the qf (or qfunction_none) 9759df49d7eSJed Brown /// 9769df49d7eSJed Brown /// ``` 9779df49d7eSJed Brown /// # use libceed::prelude::*; 9784d27c890SJeremy L Thompson /// # fn main() -> libceed::Result<()> { 9799df49d7eSJed Brown /// # let ceed = libceed::Ceed::default_init(); 980c68be7a2SJeremy L Thompson /// let qf = ceed.q_function_interior_by_name("Mass1DBuild")?; 981c68be7a2SJeremy L Thompson /// let op = ceed.operator(&qf, QFunctionOpt::None, QFunctionOpt::None)?; 982c68be7a2SJeremy L Thompson /// # Ok(()) 983c68be7a2SJeremy L Thompson /// # } 9849df49d7eSJed Brown /// ``` 985594ef120SJeremy L Thompson pub fn operator<'a, 'b>( 9869df49d7eSJed Brown &self, 9879df49d7eSJed Brown qf: impl Into<QFunctionOpt<'b>>, 9889df49d7eSJed Brown dqf: impl Into<QFunctionOpt<'b>>, 9899df49d7eSJed Brown dqfT: impl Into<QFunctionOpt<'b>>, 990594ef120SJeremy L Thompson ) -> Result<Operator<'a>> { 9919df49d7eSJed Brown Operator::create(self, qf, dqf, dqfT) 9929df49d7eSJed Brown } 9939df49d7eSJed Brown 9949df49d7eSJed Brown /// Returns an Operator that composes the action of several Operators 9959df49d7eSJed Brown /// 9969df49d7eSJed Brown /// ``` 9979df49d7eSJed Brown /// # use libceed::prelude::*; 9984d27c890SJeremy L Thompson /// # fn main() -> libceed::Result<()> { 9999df49d7eSJed Brown /// # let ceed = libceed::Ceed::default_init(); 1000c68be7a2SJeremy L Thompson /// let op = ceed.composite_operator()?; 1001c68be7a2SJeremy L Thompson /// # Ok(()) 1002c68be7a2SJeremy L Thompson /// # } 10039df49d7eSJed Brown /// ``` 1004594ef120SJeremy L Thompson pub fn composite_operator<'a>(&self) -> Result<CompositeOperator<'a>> { 10059df49d7eSJed Brown CompositeOperator::create(self) 10069df49d7eSJed Brown } 10079df49d7eSJed Brown } 10089df49d7eSJed Brown 10099df49d7eSJed Brown // ----------------------------------------------------------------------------- 10109df49d7eSJed Brown // Tests 10119df49d7eSJed Brown // ----------------------------------------------------------------------------- 10129df49d7eSJed Brown #[cfg(test)] 10139df49d7eSJed Brown mod tests { 10149df49d7eSJed Brown use super::*; 10159df49d7eSJed Brown 101689d15d5fSJeremy L Thompson fn ceed_t501() -> Result<()> { 10179df49d7eSJed Brown let resource = "/cpu/self/ref/blocked"; 10189df49d7eSJed Brown let ceed = Ceed::init(resource); 10199df49d7eSJed Brown let nelem = 4; 10209df49d7eSJed Brown let p = 3; 10219df49d7eSJed Brown let q = 4; 10229df49d7eSJed Brown let ndofs = p * nelem - nelem + 1; 10239df49d7eSJed Brown 10249df49d7eSJed Brown // Vectors 10259df49d7eSJed Brown let x = ceed.vector_from_slice(&[-1., -0.5, 0.0, 0.5, 1.0])?; 10269df49d7eSJed Brown let mut qdata = ceed.vector(nelem * q)?; 10279df49d7eSJed Brown qdata.set_value(0.0)?; 10289df49d7eSJed Brown let mut u = ceed.vector(ndofs)?; 10299df49d7eSJed Brown u.set_value(1.0)?; 10309df49d7eSJed Brown let mut v = ceed.vector(ndofs)?; 10319df49d7eSJed Brown v.set_value(0.0)?; 10329df49d7eSJed Brown 10339df49d7eSJed Brown // Restrictions 10349df49d7eSJed Brown let mut indx: Vec<i32> = vec![0; 2 * nelem]; 10359df49d7eSJed Brown for i in 0..nelem { 10369df49d7eSJed Brown indx[2 * i + 0] = i as i32; 10379df49d7eSJed Brown indx[2 * i + 1] = (i + 1) as i32; 10389df49d7eSJed Brown } 10399df49d7eSJed Brown let rx = ceed.elem_restriction(nelem, 2, 1, 1, nelem + 1, MemType::Host, &indx)?; 10409df49d7eSJed Brown let mut indu: Vec<i32> = vec![0; p * nelem]; 10419df49d7eSJed Brown for i in 0..nelem { 10429df49d7eSJed Brown indu[p * i + 0] = i as i32; 10439df49d7eSJed Brown indu[p * i + 1] = (i + 1) as i32; 10449df49d7eSJed Brown indu[p * i + 2] = (i + 2) as i32; 10459df49d7eSJed Brown } 10469df49d7eSJed Brown let ru = ceed.elem_restriction(nelem, 3, 1, 1, ndofs, MemType::Host, &indu)?; 10479df49d7eSJed Brown let strides: [i32; 3] = [1, q as i32, q as i32]; 10489df49d7eSJed Brown let rq = ceed.strided_elem_restriction(nelem, q, 1, q * nelem, strides)?; 10499df49d7eSJed Brown 10509df49d7eSJed Brown // Bases 10519df49d7eSJed Brown let bx = ceed.basis_tensor_H1_Lagrange(1, 1, 2, q, QuadMode::Gauss)?; 10529df49d7eSJed Brown let bu = ceed.basis_tensor_H1_Lagrange(1, 1, p, q, QuadMode::Gauss)?; 10539df49d7eSJed Brown 10549df49d7eSJed Brown // Build quadrature data 10559df49d7eSJed Brown let qf_build = ceed.q_function_interior_by_name("Mass1DBuild")?; 10569df49d7eSJed Brown ceed.operator(&qf_build, QFunctionOpt::None, QFunctionOpt::None)? 10579df49d7eSJed Brown .field("dx", &rx, &bx, VectorOpt::Active)? 10589df49d7eSJed Brown .field("weights", ElemRestrictionOpt::None, &bx, VectorOpt::None)? 10599df49d7eSJed Brown .field("qdata", &rq, BasisOpt::Collocated, VectorOpt::Active)? 10609df49d7eSJed Brown .apply(&x, &mut qdata)?; 10619df49d7eSJed Brown 10629df49d7eSJed Brown // Mass operator 10639df49d7eSJed Brown let qf_mass = ceed.q_function_interior_by_name("MassApply")?; 10649df49d7eSJed Brown let op_mass = ceed 10659df49d7eSJed Brown .operator(&qf_mass, QFunctionOpt::None, QFunctionOpt::None)? 10669df49d7eSJed Brown .field("u", &ru, &bu, VectorOpt::Active)? 10679df49d7eSJed Brown .field("qdata", &rq, BasisOpt::Collocated, &qdata)? 10686f97ff0aSJeremy L Thompson .field("v", &ru, &bu, VectorOpt::Active)? 10696f97ff0aSJeremy L Thompson .check()?; 10709df49d7eSJed Brown 10719df49d7eSJed Brown v.set_value(0.0)?; 10729df49d7eSJed Brown op_mass.apply(&u, &mut v)?; 10739df49d7eSJed Brown 10749df49d7eSJed Brown // Check 1075e78171edSJeremy L Thompson let sum: Scalar = v.view()?.iter().sum(); 10769df49d7eSJed Brown assert!( 10779df49d7eSJed Brown (sum - 2.0).abs() < 1e-15, 10789df49d7eSJed Brown "Incorrect interval length computed" 10799df49d7eSJed Brown ); 108089d15d5fSJeremy L Thompson Ok(()) 10819df49d7eSJed Brown } 10829df49d7eSJed Brown 10839df49d7eSJed Brown #[test] 10849df49d7eSJed Brown fn test_ceed_t501() { 10859df49d7eSJed Brown assert!(ceed_t501().is_ok()); 10869df49d7eSJed Brown } 10879df49d7eSJed Brown } 10889df49d7eSJed Brown 10899df49d7eSJed Brown // ----------------------------------------------------------------------------- 1090