19df49d7eSJed Brown // Copyright (c) 2017, Lawrence Livermore National Security, LLC. Produced at 29df49d7eSJed Brown // the Lawrence Livermore National Laboratory. LLNL-CODE-734707. All Rights 39df49d7eSJed Brown // reserved. See files LICENSE and NOTICE for details. 49df49d7eSJed Brown // 59df49d7eSJed Brown // This file is part of CEED, a collection of benchmarks, miniapps, software 69df49d7eSJed Brown // libraries and APIs for efficient high-order finite element and spectral 79df49d7eSJed Brown // element discretizations for exascale applications. For more information and 89df49d7eSJed Brown // source code availability see http://github.com/ceed. 99df49d7eSJed Brown // 109df49d7eSJed Brown // The CEED research is supported by the Exascale Computing Project 17-SC-20-SC, 119df49d7eSJed Brown // a collaborative effort of two U.S. Department of Energy organizations (Office 129df49d7eSJed Brown // of Science and the National Nuclear Security Administration) responsible for 139df49d7eSJed Brown // the planning and preparation of a capable exascale ecosystem, including 149df49d7eSJed Brown // software, applications, hardware, advanced system engineering and early 159df49d7eSJed Brown // testbed platforms, in support of the nation's exascale computing imperative 169df49d7eSJed Brown 17a1cbad85SJed Brown // Fenced `rust` code blocks included from README.md are executed as part of doctests. 18a1cbad85SJed Brown #![doc = include_str!("../README.md")] 199df49d7eSJed Brown // ----------------------------------------------------------------------------- 209df49d7eSJed Brown // Exceptions 219df49d7eSJed Brown // ----------------------------------------------------------------------------- 229df49d7eSJed Brown #![allow(non_snake_case)] 239df49d7eSJed Brown 249df49d7eSJed Brown // ----------------------------------------------------------------------------- 259df49d7eSJed Brown // Crate prelude 269df49d7eSJed Brown // ----------------------------------------------------------------------------- 279df49d7eSJed Brown use crate::prelude::*; 289df49d7eSJed Brown use std::sync::Once; 299df49d7eSJed Brown 309df49d7eSJed Brown pub mod prelude { 319df49d7eSJed Brown pub use crate::{ 329df49d7eSJed Brown basis::{self, Basis, BasisOpt}, 339df49d7eSJed Brown elem_restriction::{self, ElemRestriction, ElemRestrictionOpt}, 349df49d7eSJed Brown operator::{self, CompositeOperator, Operator}, 359df49d7eSJed Brown qfunction::{ 369df49d7eSJed Brown self, QFunction, QFunctionByName, QFunctionInputs, QFunctionOpt, QFunctionOutputs, 379df49d7eSJed Brown }, 389df49d7eSJed Brown vector::{self, Vector, VectorOpt}, 3980a9ef05SNatalie Beams ElemTopology, EvalMode, MemType, NormType, QuadMode, Scalar, TransposeMode, 4080a9ef05SNatalie Beams CEED_STRIDES_BACKEND, EPSILON, MAX_QFUNCTION_FIELDS, 419df49d7eSJed Brown }; 42630ad4c9Sjeremylt pub(crate) use libceed_sys::bind_ceed; 439df49d7eSJed Brown pub(crate) use std::convert::TryFrom; 449df49d7eSJed Brown pub(crate) use std::ffi::CString; 459df49d7eSJed Brown pub(crate) use std::fmt; 46*1142270cSJeremy L Thompson pub(crate) use std::marker::PhantomData; 479df49d7eSJed Brown } 489df49d7eSJed Brown 499df49d7eSJed Brown // ----------------------------------------------------------------------------- 509df49d7eSJed Brown // Modules 519df49d7eSJed Brown // ----------------------------------------------------------------------------- 529df49d7eSJed Brown pub mod basis; 539df49d7eSJed Brown pub mod elem_restriction; 549df49d7eSJed Brown pub mod operator; 559df49d7eSJed Brown pub mod qfunction; 569df49d7eSJed Brown pub mod vector; 579df49d7eSJed Brown 589df49d7eSJed Brown // ----------------------------------------------------------------------------- 5980a9ef05SNatalie Beams // Typedef for scalar 6080a9ef05SNatalie Beams // ----------------------------------------------------------------------------- 6180a9ef05SNatalie Beams pub type Scalar = bind_ceed::CeedScalar; 6280a9ef05SNatalie Beams 6380a9ef05SNatalie Beams // ----------------------------------------------------------------------------- 649df49d7eSJed Brown // Constants for library 659df49d7eSJed Brown // ----------------------------------------------------------------------------- 669df49d7eSJed Brown const MAX_BUFFER_LENGTH: u64 = 4096; 679df49d7eSJed Brown pub const MAX_QFUNCTION_FIELDS: usize = 16; 689df49d7eSJed Brown pub const CEED_STRIDES_BACKEND: [i32; 3] = [0; 3]; 6980a9ef05SNatalie Beams pub const EPSILON: crate::Scalar = bind_ceed::CEED_EPSILON as crate::Scalar; 709df49d7eSJed Brown 719df49d7eSJed Brown // ----------------------------------------------------------------------------- 729df49d7eSJed Brown // Enums for libCEED 739df49d7eSJed Brown // ----------------------------------------------------------------------------- 749df49d7eSJed Brown #[derive(Clone, Copy, PartialEq, Eq)] 759df49d7eSJed Brown /// Many Ceed interfaces take or return pointers to memory. This enum is used to 769df49d7eSJed Brown /// specify where the memory being provided or requested must reside. 779df49d7eSJed Brown pub enum MemType { 789df49d7eSJed Brown Host = bind_ceed::CeedMemType_CEED_MEM_HOST as isize, 799df49d7eSJed Brown Device = bind_ceed::CeedMemType_CEED_MEM_DEVICE as isize, 809df49d7eSJed Brown } 819df49d7eSJed Brown 829df49d7eSJed Brown #[derive(Clone, Copy, PartialEq, Eq)] 839df49d7eSJed Brown // OwnPointer will not be used by user but is included for internal use 849df49d7eSJed Brown #[allow(dead_code)] 859df49d7eSJed Brown /// Conveys ownership status of arrays passed to Ceed interfaces. 869df49d7eSJed Brown pub(crate) enum CopyMode { 879df49d7eSJed Brown CopyValues = bind_ceed::CeedCopyMode_CEED_COPY_VALUES as isize, 889df49d7eSJed Brown UsePointer = bind_ceed::CeedCopyMode_CEED_USE_POINTER as isize, 899df49d7eSJed Brown OwnPointer = bind_ceed::CeedCopyMode_CEED_OWN_POINTER as isize, 909df49d7eSJed Brown } 919df49d7eSJed Brown 929df49d7eSJed Brown #[derive(Clone, Copy, PartialEq, Eq)] 939df49d7eSJed Brown /// Denotes type of vector norm to be computed 949df49d7eSJed Brown pub enum NormType { 959df49d7eSJed Brown One = bind_ceed::CeedNormType_CEED_NORM_1 as isize, 969df49d7eSJed Brown Two = bind_ceed::CeedNormType_CEED_NORM_2 as isize, 979df49d7eSJed Brown Max = bind_ceed::CeedNormType_CEED_NORM_MAX as isize, 989df49d7eSJed Brown } 999df49d7eSJed Brown 1009df49d7eSJed Brown #[derive(Clone, Copy, PartialEq, Eq)] 1019df49d7eSJed Brown /// Denotes whether a linear transformation or its transpose should be applied 1029df49d7eSJed Brown pub enum TransposeMode { 1039df49d7eSJed Brown NoTranspose = bind_ceed::CeedTransposeMode_CEED_NOTRANSPOSE as isize, 1049df49d7eSJed Brown Transpose = bind_ceed::CeedTransposeMode_CEED_TRANSPOSE as isize, 1059df49d7eSJed Brown } 1069df49d7eSJed Brown 1079df49d7eSJed Brown #[derive(Clone, Copy, PartialEq, Eq)] 1089df49d7eSJed Brown /// Type of quadrature; also used for location of nodes 1099df49d7eSJed Brown pub enum QuadMode { 1109df49d7eSJed Brown Gauss = bind_ceed::CeedQuadMode_CEED_GAUSS as isize, 1119df49d7eSJed Brown GaussLobatto = bind_ceed::CeedQuadMode_CEED_GAUSS_LOBATTO as isize, 1129df49d7eSJed Brown } 1139df49d7eSJed Brown 1149df49d7eSJed Brown #[derive(Clone, Copy, PartialEq, Eq)] 1159df49d7eSJed Brown /// Type of basis shape to create non-tensor H1 element basis 1169df49d7eSJed Brown pub enum ElemTopology { 1179df49d7eSJed Brown Line = bind_ceed::CeedElemTopology_CEED_LINE as isize, 1189df49d7eSJed Brown Triangle = bind_ceed::CeedElemTopology_CEED_TRIANGLE as isize, 1199df49d7eSJed Brown Quad = bind_ceed::CeedElemTopology_CEED_QUAD as isize, 1209df49d7eSJed Brown Tet = bind_ceed::CeedElemTopology_CEED_TET as isize, 1219df49d7eSJed Brown Pyramid = bind_ceed::CeedElemTopology_CEED_PYRAMID as isize, 1229df49d7eSJed Brown Prism = bind_ceed::CeedElemTopology_CEED_PRISM as isize, 1239df49d7eSJed Brown Hex = bind_ceed::CeedElemTopology_CEED_HEX as isize, 1249df49d7eSJed Brown } 1259df49d7eSJed Brown 126c68be7a2SJeremy L Thompson #[derive(Clone, Copy, Debug, PartialEq, Eq)] 1279df49d7eSJed Brown /// Basis evaluation mode 1289df49d7eSJed Brown pub enum EvalMode { 1299df49d7eSJed Brown None = bind_ceed::CeedEvalMode_CEED_EVAL_NONE as isize, 1309df49d7eSJed Brown Interp = bind_ceed::CeedEvalMode_CEED_EVAL_INTERP as isize, 1319df49d7eSJed Brown Grad = bind_ceed::CeedEvalMode_CEED_EVAL_GRAD as isize, 1329df49d7eSJed Brown Div = bind_ceed::CeedEvalMode_CEED_EVAL_DIV as isize, 1339df49d7eSJed Brown Curl = bind_ceed::CeedEvalMode_CEED_EVAL_CURL as isize, 1349df49d7eSJed Brown Weight = bind_ceed::CeedEvalMode_CEED_EVAL_WEIGHT as isize, 1359df49d7eSJed Brown } 136c68be7a2SJeremy L Thompson impl EvalMode { 137c68be7a2SJeremy L Thompson pub(crate) fn from_u32(value: u32) -> EvalMode { 138c68be7a2SJeremy L Thompson match value { 139c68be7a2SJeremy L Thompson bind_ceed::CeedEvalMode_CEED_EVAL_NONE => EvalMode::None, 140c68be7a2SJeremy L Thompson bind_ceed::CeedEvalMode_CEED_EVAL_INTERP => EvalMode::Interp, 141c68be7a2SJeremy L Thompson bind_ceed::CeedEvalMode_CEED_EVAL_GRAD => EvalMode::Grad, 142c68be7a2SJeremy L Thompson bind_ceed::CeedEvalMode_CEED_EVAL_DIV => EvalMode::Div, 143c68be7a2SJeremy L Thompson bind_ceed::CeedEvalMode_CEED_EVAL_CURL => EvalMode::Curl, 144c68be7a2SJeremy L Thompson bind_ceed::CeedEvalMode_CEED_EVAL_WEIGHT => EvalMode::Weight, 145c68be7a2SJeremy L Thompson _ => panic!("Unknown value: {}", value), 146c68be7a2SJeremy L Thompson } 147c68be7a2SJeremy L Thompson } 148c68be7a2SJeremy L Thompson } 1499df49d7eSJed Brown 1509df49d7eSJed Brown // ----------------------------------------------------------------------------- 1519df49d7eSJed Brown // Ceed error 1529df49d7eSJed Brown // ----------------------------------------------------------------------------- 1539df49d7eSJed Brown type Result<T> = std::result::Result<T, CeedError>; 1549df49d7eSJed Brown 1559df49d7eSJed Brown #[derive(Debug)] 1569df49d7eSJed Brown pub struct CeedError { 1579df49d7eSJed Brown pub message: String, 1589df49d7eSJed Brown } 1599df49d7eSJed Brown 1609df49d7eSJed Brown impl fmt::Display for CeedError { 1619df49d7eSJed Brown fn fmt(&self, f: &mut fmt::Formatter) -> fmt::Result { 1629df49d7eSJed Brown write!(f, "{}", self.message) 1639df49d7eSJed Brown } 1649df49d7eSJed Brown } 1659df49d7eSJed Brown 1669df49d7eSJed Brown // ----------------------------------------------------------------------------- 167*1142270cSJeremy L Thompson // Internal error checker 168*1142270cSJeremy L Thompson // ----------------------------------------------------------------------------- 169*1142270cSJeremy L Thompson #[doc(hidden)] 170*1142270cSJeremy L Thompson pub(crate) fn check_error(ceed_ptr: bind_ceed::Ceed, ierr: i32) -> Result<i32> { 171*1142270cSJeremy L Thompson // Return early if code is clean 172*1142270cSJeremy L Thompson if ierr == bind_ceed::CeedErrorType_CEED_ERROR_SUCCESS { 173*1142270cSJeremy L Thompson return Ok(ierr); 174*1142270cSJeremy L Thompson } 175*1142270cSJeremy L Thompson // Retrieve error message 176*1142270cSJeremy L Thompson let mut ptr: *const std::os::raw::c_char = std::ptr::null_mut(); 177*1142270cSJeremy L Thompson let c_str = unsafe { 178*1142270cSJeremy L Thompson bind_ceed::CeedGetErrorMessage(ceed_ptr, &mut ptr); 179*1142270cSJeremy L Thompson std::ffi::CStr::from_ptr(ptr) 180*1142270cSJeremy L Thompson }; 181*1142270cSJeremy L Thompson let message = c_str.to_string_lossy().to_string(); 182*1142270cSJeremy L Thompson // Panic if negative code, otherwise return error 183*1142270cSJeremy L Thompson if ierr < bind_ceed::CeedErrorType_CEED_ERROR_SUCCESS { 184*1142270cSJeremy L Thompson panic!("{}", message); 185*1142270cSJeremy L Thompson } 186*1142270cSJeremy L Thompson Err(CeedError { message }) 187*1142270cSJeremy L Thompson } 188*1142270cSJeremy L Thompson 189*1142270cSJeremy L Thompson // ----------------------------------------------------------------------------- 1909df49d7eSJed Brown // Ceed error handler 1919df49d7eSJed Brown // ----------------------------------------------------------------------------- 1929df49d7eSJed Brown pub enum CeedErrorHandler { 1939df49d7eSJed Brown ErrorAbort, 1949df49d7eSJed Brown ErrorExit, 1959df49d7eSJed Brown ErrorReturn, 1969df49d7eSJed Brown ErrorStore, 1979df49d7eSJed Brown } 1989df49d7eSJed Brown 1999df49d7eSJed Brown // ----------------------------------------------------------------------------- 2009df49d7eSJed Brown // Ceed context wrapper 2019df49d7eSJed Brown // ----------------------------------------------------------------------------- 2029df49d7eSJed Brown /// A Ceed is a library context representing control of a logical hardware 2039df49d7eSJed Brown /// resource. 2049df49d7eSJed Brown #[derive(Debug)] 2059df49d7eSJed Brown pub struct Ceed { 2069df49d7eSJed Brown ptr: bind_ceed::Ceed, 2079df49d7eSJed Brown } 2089df49d7eSJed Brown 2099df49d7eSJed Brown // ----------------------------------------------------------------------------- 2109df49d7eSJed Brown // Destructor 2119df49d7eSJed Brown // ----------------------------------------------------------------------------- 2129df49d7eSJed Brown impl Drop for Ceed { 2139df49d7eSJed Brown fn drop(&mut self) { 2149df49d7eSJed Brown unsafe { 2159df49d7eSJed Brown bind_ceed::CeedDestroy(&mut self.ptr); 2169df49d7eSJed Brown } 2179df49d7eSJed Brown } 2189df49d7eSJed Brown } 2199df49d7eSJed Brown 2209df49d7eSJed Brown // ----------------------------------------------------------------------------- 2219df49d7eSJed Brown // Display 2229df49d7eSJed Brown // ----------------------------------------------------------------------------- 2239df49d7eSJed Brown impl fmt::Display for Ceed { 2249df49d7eSJed Brown /// View a Ceed 2259df49d7eSJed Brown /// 2269df49d7eSJed Brown /// ``` 2279df49d7eSJed Brown /// let ceed = libceed::Ceed::init("/cpu/self/ref/serial"); 2289df49d7eSJed Brown /// println!("{}", ceed); 2299df49d7eSJed Brown /// ``` 2309df49d7eSJed Brown fn fmt(&self, f: &mut fmt::Formatter) -> fmt::Result { 2319df49d7eSJed Brown let mut ptr = std::ptr::null_mut(); 2329df49d7eSJed Brown let mut sizeloc = crate::MAX_BUFFER_LENGTH; 2339df49d7eSJed Brown let cstring = unsafe { 2349df49d7eSJed Brown let file = bind_ceed::open_memstream(&mut ptr, &mut sizeloc); 2359df49d7eSJed Brown bind_ceed::CeedView(self.ptr, file); 2369df49d7eSJed Brown bind_ceed::fclose(file); 2379df49d7eSJed Brown CString::from_raw(ptr) 2389df49d7eSJed Brown }; 2399df49d7eSJed Brown cstring.to_string_lossy().fmt(f) 2409df49d7eSJed Brown } 2419df49d7eSJed Brown } 2429df49d7eSJed Brown 2439df49d7eSJed Brown static REGISTER: Once = Once::new(); 2449df49d7eSJed Brown 2459df49d7eSJed Brown // ----------------------------------------------------------------------------- 2469df49d7eSJed Brown // Object constructors 2479df49d7eSJed Brown // ----------------------------------------------------------------------------- 2489df49d7eSJed Brown impl Ceed { 2499df49d7eSJed Brown /// Returns a Ceed context initialized with the specified resource 2509df49d7eSJed Brown /// 2519df49d7eSJed Brown /// # arguments 2529df49d7eSJed Brown /// 2539df49d7eSJed Brown /// * `resource` - Resource to use, e.g., "/cpu/self" 2549df49d7eSJed Brown /// 2559df49d7eSJed Brown /// ``` 2569df49d7eSJed Brown /// let ceed = libceed::Ceed::init("/cpu/self/ref/serial"); 2579df49d7eSJed Brown /// ``` 2589df49d7eSJed Brown pub fn init(resource: &str) -> Self { 2599df49d7eSJed Brown Ceed::init_with_error_handler(resource, CeedErrorHandler::ErrorStore) 2609df49d7eSJed Brown } 2619df49d7eSJed Brown 2629df49d7eSJed Brown /// Returns a Ceed context initialized with the specified resource 2639df49d7eSJed Brown /// 2649df49d7eSJed Brown /// # arguments 2659df49d7eSJed Brown /// 2669df49d7eSJed Brown /// * `resource` - Resource to use, e.g., "/cpu/self" 2679df49d7eSJed Brown /// 2689df49d7eSJed Brown /// ``` 2699df49d7eSJed Brown /// let ceed = libceed::Ceed::init_with_error_handler( 2709df49d7eSJed Brown /// "/cpu/self/ref/serial", 2719df49d7eSJed Brown /// libceed::CeedErrorHandler::ErrorAbort, 2729df49d7eSJed Brown /// ); 2739df49d7eSJed Brown /// ``` 2749df49d7eSJed Brown pub fn init_with_error_handler(resource: &str, handler: CeedErrorHandler) -> Self { 2759df49d7eSJed Brown REGISTER.call_once(|| unsafe { 2769df49d7eSJed Brown bind_ceed::CeedRegisterAll(); 2779df49d7eSJed Brown bind_ceed::CeedQFunctionRegisterAll(); 2789df49d7eSJed Brown }); 2799df49d7eSJed Brown 2809df49d7eSJed Brown // Convert to C string 2819df49d7eSJed Brown let c_resource = CString::new(resource).expect("CString::new failed"); 2829df49d7eSJed Brown 2839df49d7eSJed Brown // Get error handler pointer 2849df49d7eSJed Brown let eh = match handler { 2859df49d7eSJed Brown CeedErrorHandler::ErrorAbort => bind_ceed::CeedErrorAbort, 2869df49d7eSJed Brown CeedErrorHandler::ErrorExit => bind_ceed::CeedErrorExit, 2879df49d7eSJed Brown CeedErrorHandler::ErrorReturn => bind_ceed::CeedErrorReturn, 2889df49d7eSJed Brown CeedErrorHandler::ErrorStore => bind_ceed::CeedErrorStore, 2899df49d7eSJed Brown }; 2909df49d7eSJed Brown 2919df49d7eSJed Brown // Call to libCEED 2929df49d7eSJed Brown let mut ptr = std::ptr::null_mut(); 2939df49d7eSJed Brown let mut ierr = unsafe { bind_ceed::CeedInit(c_resource.as_ptr() as *const i8, &mut ptr) }; 2949df49d7eSJed Brown if ierr != 0 { 2959df49d7eSJed Brown panic!("Error initializing backend resource: {}", resource) 2969df49d7eSJed Brown } 2979df49d7eSJed Brown ierr = unsafe { bind_ceed::CeedSetErrorHandler(ptr, Some(eh)) }; 2989df49d7eSJed Brown let ceed = Ceed { ptr }; 2999df49d7eSJed Brown ceed.check_error(ierr).unwrap(); 3009df49d7eSJed Brown ceed 3019df49d7eSJed Brown } 3029df49d7eSJed Brown 3039df49d7eSJed Brown /// Default initializer for testing 3049df49d7eSJed Brown #[doc(hidden)] 3059df49d7eSJed Brown pub fn default_init() -> Self { 3069df49d7eSJed Brown // Convert to C string 3079df49d7eSJed Brown let resource = "/cpu/self/ref/serial"; 3089df49d7eSJed Brown crate::Ceed::init(resource) 3099df49d7eSJed Brown } 3109df49d7eSJed Brown 3119df49d7eSJed Brown /// Internal error checker 3129df49d7eSJed Brown #[doc(hidden)] 3139df49d7eSJed Brown fn check_error(&self, ierr: i32) -> Result<i32> { 3149df49d7eSJed Brown // Return early if code is clean 3159df49d7eSJed Brown if ierr == bind_ceed::CeedErrorType_CEED_ERROR_SUCCESS { 3169df49d7eSJed Brown return Ok(ierr); 3179df49d7eSJed Brown } 3189df49d7eSJed Brown // Retrieve error message 3199df49d7eSJed Brown let mut ptr: *const std::os::raw::c_char = std::ptr::null_mut(); 3209df49d7eSJed Brown let c_str = unsafe { 3219df49d7eSJed Brown bind_ceed::CeedGetErrorMessage(self.ptr, &mut ptr); 3229df49d7eSJed Brown std::ffi::CStr::from_ptr(ptr) 3239df49d7eSJed Brown }; 3249df49d7eSJed Brown let message = c_str.to_string_lossy().to_string(); 3259df49d7eSJed Brown // Panic if negative code, otherwise return error 3269df49d7eSJed Brown if ierr < bind_ceed::CeedErrorType_CEED_ERROR_SUCCESS { 3279df49d7eSJed Brown panic!("{}", message); 3289df49d7eSJed Brown } 3299df49d7eSJed Brown Err(CeedError { message }) 3309df49d7eSJed Brown } 3319df49d7eSJed Brown 3329df49d7eSJed Brown /// Returns full resource name for a Ceed context 3339df49d7eSJed Brown /// 3349df49d7eSJed Brown /// ``` 3359df49d7eSJed Brown /// let ceed = libceed::Ceed::init("/cpu/self/ref/serial"); 3369df49d7eSJed Brown /// let resource = ceed.resource(); 3379df49d7eSJed Brown /// 3389df49d7eSJed Brown /// assert_eq!(resource, "/cpu/self/ref/serial".to_string()) 3399df49d7eSJed Brown /// ``` 3409df49d7eSJed Brown pub fn resource(&self) -> String { 3419df49d7eSJed Brown let mut ptr: *const std::os::raw::c_char = std::ptr::null_mut(); 3429df49d7eSJed Brown let c_str = unsafe { 3439df49d7eSJed Brown bind_ceed::CeedGetResource(self.ptr, &mut ptr); 3449df49d7eSJed Brown std::ffi::CStr::from_ptr(ptr) 3459df49d7eSJed Brown }; 3469df49d7eSJed Brown c_str.to_string_lossy().to_string() 3479df49d7eSJed Brown } 3489df49d7eSJed Brown 3499df49d7eSJed Brown /// Returns a CeedVector of the specified length (does not allocate memory) 3509df49d7eSJed Brown /// 3519df49d7eSJed Brown /// # arguments 3529df49d7eSJed Brown /// 3539df49d7eSJed Brown /// * `n` - Length of vector 3549df49d7eSJed Brown /// 3559df49d7eSJed Brown /// ``` 3569df49d7eSJed Brown /// # use libceed::prelude::*; 357c68be7a2SJeremy L Thompson /// # fn main() -> Result<(), libceed::CeedError> { 3589df49d7eSJed Brown /// # let ceed = libceed::Ceed::default_init(); 359c68be7a2SJeremy L Thompson /// let vec = ceed.vector(10)?; 360c68be7a2SJeremy L Thompson /// # Ok(()) 361c68be7a2SJeremy L Thompson /// # } 3629df49d7eSJed Brown /// ``` 3639df49d7eSJed Brown pub fn vector(&self, n: usize) -> Result<Vector> { 3649df49d7eSJed Brown Vector::create(self, n) 3659df49d7eSJed Brown } 3669df49d7eSJed Brown 3679df49d7eSJed Brown /// Create a Vector initialized with the data (copied) from a slice 3689df49d7eSJed Brown /// 3699df49d7eSJed Brown /// # arguments 3709df49d7eSJed Brown /// 3719df49d7eSJed Brown /// * `slice` - Slice containing data 3729df49d7eSJed Brown /// 3739df49d7eSJed Brown /// ``` 3749df49d7eSJed Brown /// # use libceed::prelude::*; 375c68be7a2SJeremy L Thompson /// # fn main() -> Result<(), libceed::CeedError> { 3769df49d7eSJed Brown /// # let ceed = libceed::Ceed::default_init(); 377c68be7a2SJeremy L Thompson /// let vec = ceed.vector_from_slice(&[1., 2., 3.])?; 3789df49d7eSJed Brown /// assert_eq!(vec.length(), 3); 379c68be7a2SJeremy L Thompson /// # Ok(()) 380c68be7a2SJeremy L Thompson /// # } 3819df49d7eSJed Brown /// ``` 38280a9ef05SNatalie Beams pub fn vector_from_slice(&self, slice: &[crate::Scalar]) -> Result<Vector> { 3839df49d7eSJed Brown Vector::from_slice(self, slice) 3849df49d7eSJed Brown } 3859df49d7eSJed Brown 3869df49d7eSJed Brown /// Returns a ElemRestriction 3879df49d7eSJed Brown /// 3889df49d7eSJed Brown /// # arguments 3899df49d7eSJed Brown /// 3909df49d7eSJed Brown /// * `nelem` - Number of elements described in the offsets array 3919df49d7eSJed Brown /// * `elemsize` - Size (number of "nodes") per element 3929df49d7eSJed Brown /// * `ncomp` - Number of field components per interpolation node (1 3939df49d7eSJed Brown /// for scalar fields) 3949df49d7eSJed Brown /// * `compstride` - Stride between components for the same Lvector "node". 3959df49d7eSJed Brown /// Data for node `i`, component `j`, element `k` can be 3969df49d7eSJed Brown /// found in the Lvector at index 3979df49d7eSJed Brown /// `offsets[i + k*elemsize] + j*compstride`. 3989df49d7eSJed Brown /// * `lsize` - The size of the Lvector. This vector may be larger 3999df49d7eSJed Brown /// than the elements and fields given by this 4009df49d7eSJed Brown /// restriction. 4019df49d7eSJed Brown /// * `mtype` - Memory type of the offsets array, see CeedMemType 4029df49d7eSJed Brown /// * `offsets` - Array of shape `[nelem, elemsize]`. Row `i` holds the 4039df49d7eSJed Brown /// ordered list of the offsets (into the input CeedVector) 4049df49d7eSJed Brown /// for the unknowns corresponding to element `i`, where 4059df49d7eSJed Brown /// `0 <= i < nelem`. All offsets must be in the range 4069df49d7eSJed Brown /// `[0, lsize - 1]`. 4079df49d7eSJed Brown /// 4089df49d7eSJed Brown /// ``` 4099df49d7eSJed Brown /// # use libceed::prelude::*; 410c68be7a2SJeremy L Thompson /// # fn main() -> Result<(), libceed::CeedError> { 4119df49d7eSJed Brown /// # let ceed = libceed::Ceed::default_init(); 4129df49d7eSJed Brown /// let nelem = 3; 4139df49d7eSJed Brown /// let mut ind: Vec<i32> = vec![0; 2 * nelem]; 4149df49d7eSJed Brown /// for i in 0..nelem { 4159df49d7eSJed Brown /// ind[2 * i + 0] = i as i32; 4169df49d7eSJed Brown /// ind[2 * i + 1] = (i + 1) as i32; 4179df49d7eSJed Brown /// } 418c68be7a2SJeremy L Thompson /// let r = ceed.elem_restriction(nelem, 2, 1, 1, nelem + 1, MemType::Host, &ind)?; 419c68be7a2SJeremy L Thompson /// # Ok(()) 420c68be7a2SJeremy L Thompson /// # } 4219df49d7eSJed Brown /// ``` 4229df49d7eSJed Brown pub fn elem_restriction( 4239df49d7eSJed Brown &self, 4249df49d7eSJed Brown nelem: usize, 4259df49d7eSJed Brown elemsize: usize, 4269df49d7eSJed Brown ncomp: usize, 4279df49d7eSJed Brown compstride: usize, 4289df49d7eSJed Brown lsize: usize, 4299df49d7eSJed Brown mtype: MemType, 4309df49d7eSJed Brown offsets: &[i32], 4319df49d7eSJed Brown ) -> Result<ElemRestriction> { 4329df49d7eSJed Brown ElemRestriction::create( 4339df49d7eSJed Brown self, nelem, elemsize, ncomp, compstride, lsize, mtype, offsets, 4349df49d7eSJed Brown ) 4359df49d7eSJed Brown } 4369df49d7eSJed Brown 4379df49d7eSJed Brown /// Returns a ElemRestriction 4389df49d7eSJed Brown /// 4399df49d7eSJed Brown /// # arguments 4409df49d7eSJed Brown /// 4419df49d7eSJed Brown /// * `nelem` - Number of elements described in the offsets array 4429df49d7eSJed Brown /// * `elemsize` - Size (number of "nodes") per element 4439df49d7eSJed Brown /// * `ncomp` - Number of field components per interpolation node (1 4449df49d7eSJed Brown /// for scalar fields) 4459df49d7eSJed Brown /// * `compstride` - Stride between components for the same Lvector "node". 4469df49d7eSJed Brown /// Data for node `i`, component `j`, element `k` can be 4479df49d7eSJed Brown /// found in the Lvector at index 4489df49d7eSJed Brown /// `offsets[i + k*elemsize] + j*compstride`. 4499df49d7eSJed Brown /// * `lsize` - The size of the Lvector. This vector may be larger 4509df49d7eSJed Brown /// than the elements and fields given by this restriction. 4519df49d7eSJed Brown /// * `strides` - Array for strides between `[nodes, components, elements]`. 4529df49d7eSJed Brown /// Data for node `i`, component `j`, element `k` can be 4539df49d7eSJed Brown /// found in the Lvector at index 4549df49d7eSJed Brown /// `i*strides[0] + j*strides[1] + k*strides[2]`. 4559df49d7eSJed Brown /// CEED_STRIDES_BACKEND may be used with vectors created 4569df49d7eSJed Brown /// by a Ceed backend. 4579df49d7eSJed Brown /// 4589df49d7eSJed Brown /// ``` 4599df49d7eSJed Brown /// # use libceed::prelude::*; 460c68be7a2SJeremy L Thompson /// # fn main() -> Result<(), libceed::CeedError> { 4619df49d7eSJed Brown /// # let ceed = libceed::Ceed::default_init(); 4629df49d7eSJed Brown /// let nelem = 3; 4639df49d7eSJed Brown /// let strides: [i32; 3] = [1, 2, 2]; 464c68be7a2SJeremy L Thompson /// let r = ceed.strided_elem_restriction(nelem, 2, 1, nelem * 2, strides)?; 465c68be7a2SJeremy L Thompson /// # Ok(()) 466c68be7a2SJeremy L Thompson /// # } 4679df49d7eSJed Brown /// ``` 4689df49d7eSJed Brown pub fn strided_elem_restriction( 4699df49d7eSJed Brown &self, 4709df49d7eSJed Brown nelem: usize, 4719df49d7eSJed Brown elemsize: usize, 4729df49d7eSJed Brown ncomp: usize, 4739df49d7eSJed Brown lsize: usize, 4749df49d7eSJed Brown strides: [i32; 3], 4759df49d7eSJed Brown ) -> Result<ElemRestriction> { 4769df49d7eSJed Brown ElemRestriction::create_strided(self, nelem, elemsize, ncomp, lsize, strides) 4779df49d7eSJed Brown } 4789df49d7eSJed Brown 4799df49d7eSJed Brown /// Returns a tensor-product basis 4809df49d7eSJed Brown /// 4819df49d7eSJed Brown /// # arguments 4829df49d7eSJed Brown /// 4839df49d7eSJed Brown /// * `dim` - Topological dimension of element 4849df49d7eSJed Brown /// * `ncomp` - Number of field components (1 for scalar fields) 4859df49d7eSJed Brown /// * `P1d` - Number of Gauss-Lobatto nodes in one dimension. The 4869df49d7eSJed Brown /// polynomial degree of the resulting `Q_k` element is 4879df49d7eSJed Brown /// `k=P-1`. 4889df49d7eSJed Brown /// * `Q1d` - Number of quadrature points in one dimension 4899df49d7eSJed Brown /// * `interp1d` - Row-major `(Q1d * P1d)` matrix expressing the values of 4909df49d7eSJed Brown /// nodal basis functions at quadrature points 4919df49d7eSJed Brown /// * `grad1d` - Row-major `(Q1d * P1d)` matrix expressing derivatives of 4929df49d7eSJed Brown /// nodal basis functions at quadrature points 4939df49d7eSJed Brown /// * `qref1d` - Array of length `Q1d` holding the locations of quadrature 4949df49d7eSJed Brown /// points on the 1D reference element `[-1, 1]` 4959df49d7eSJed Brown /// * `qweight1d` - Array of length `Q1d` holding the quadrature weights on 4969df49d7eSJed Brown /// the reference element 4979df49d7eSJed Brown /// 4989df49d7eSJed Brown /// ``` 4999df49d7eSJed Brown /// # use libceed::prelude::*; 500c68be7a2SJeremy L Thompson /// # fn main() -> Result<(), libceed::CeedError> { 5019df49d7eSJed Brown /// # let ceed = libceed::Ceed::default_init(); 5029df49d7eSJed Brown /// let interp1d = [ 0.62994317, 0.47255875, -0.14950343, 0.04700152, 5039df49d7eSJed Brown /// -0.07069480, 0.97297619, 0.13253993, -0.03482132, 5049df49d7eSJed Brown /// -0.03482132, 0.13253993, 0.97297619, -0.07069480, 5059df49d7eSJed Brown /// 0.04700152, -0.14950343, 0.47255875, 0.62994317]; 5069df49d7eSJed Brown /// let grad1d = [-2.34183742, 2.78794489, -0.63510411, 0.18899664, 5079df49d7eSJed Brown /// -0.51670214, -0.48795249, 1.33790510, -0.33325047, 5089df49d7eSJed Brown // 0.33325047, -1.33790510, 0.48795249, 0.51670214, 5099df49d7eSJed Brown /// -0.18899664, 0.63510411, -2.78794489, 2.34183742]; 5109df49d7eSJed Brown /// let qref1d = [-0.86113631, -0.33998104, 0.33998104, 0.86113631]; 5119df49d7eSJed Brown /// let qweight1d = [ 0.34785485, 0.65214515, 0.65214515, 0.34785485]; 5129df49d7eSJed Brown /// let b = ceed. 513c68be7a2SJeremy L Thompson /// basis_tensor_H1(2, 1, 4, 4, &interp1d, &grad1d, &qref1d, &qweight1d)?; 514c68be7a2SJeremy L Thompson /// # Ok(()) 515c68be7a2SJeremy L Thompson /// # } 5169df49d7eSJed Brown /// ``` 5179df49d7eSJed Brown pub fn basis_tensor_H1( 5189df49d7eSJed Brown &self, 5199df49d7eSJed Brown dim: usize, 5209df49d7eSJed Brown ncomp: usize, 5219df49d7eSJed Brown P1d: usize, 5229df49d7eSJed Brown Q1d: usize, 52380a9ef05SNatalie Beams interp1d: &[crate::Scalar], 52480a9ef05SNatalie Beams grad1d: &[crate::Scalar], 52580a9ef05SNatalie Beams qref1d: &[crate::Scalar], 52680a9ef05SNatalie Beams qweight1d: &[crate::Scalar], 5279df49d7eSJed Brown ) -> Result<Basis> { 5289df49d7eSJed Brown Basis::create_tensor_H1( 5299df49d7eSJed Brown self, dim, ncomp, P1d, Q1d, interp1d, grad1d, qref1d, qweight1d, 5309df49d7eSJed Brown ) 5319df49d7eSJed Brown } 5329df49d7eSJed Brown 5339df49d7eSJed Brown /// Returns a tensor-product Lagrange basis 5349df49d7eSJed Brown /// 5359df49d7eSJed Brown /// # arguments 5369df49d7eSJed Brown /// 5379df49d7eSJed Brown /// * `dim` - Topological dimension of element 5389df49d7eSJed Brown /// * `ncomp` - Number of field components (1 for scalar fields) 5399df49d7eSJed Brown /// * `P` - Number of Gauss-Lobatto nodes in one dimension. The 5409df49d7eSJed Brown /// polynomial degree of the resulting `Q_k` element is `k=P-1`. 5419df49d7eSJed Brown /// * `Q` - Number of quadrature points in one dimension 5429df49d7eSJed Brown /// * `qmode` - Distribution of the `Q` quadrature points (affects order of 5439df49d7eSJed Brown /// accuracy for the quadrature) 5449df49d7eSJed Brown /// 5459df49d7eSJed Brown /// ``` 5469df49d7eSJed Brown /// # use libceed::prelude::*; 547c68be7a2SJeremy L Thompson /// # fn main() -> Result<(), libceed::CeedError> { 5489df49d7eSJed Brown /// # let ceed = libceed::Ceed::default_init(); 549c68be7a2SJeremy L Thompson /// let b = ceed.basis_tensor_H1_Lagrange(2, 1, 3, 4, QuadMode::Gauss)?; 550c68be7a2SJeremy L Thompson /// # Ok(()) 551c68be7a2SJeremy L Thompson /// # } 5529df49d7eSJed Brown /// ``` 5539df49d7eSJed Brown pub fn basis_tensor_H1_Lagrange( 5549df49d7eSJed Brown &self, 5559df49d7eSJed Brown dim: usize, 5569df49d7eSJed Brown ncomp: usize, 5579df49d7eSJed Brown P: usize, 5589df49d7eSJed Brown Q: usize, 5599df49d7eSJed Brown qmode: QuadMode, 5609df49d7eSJed Brown ) -> Result<Basis> { 5619df49d7eSJed Brown Basis::create_tensor_H1_Lagrange(self, dim, ncomp, P, Q, qmode) 5629df49d7eSJed Brown } 5639df49d7eSJed Brown 5649df49d7eSJed Brown /// Returns a tensor-product basis 5659df49d7eSJed Brown /// 5669df49d7eSJed Brown /// # arguments 5679df49d7eSJed Brown /// 5689df49d7eSJed Brown /// * `topo` - Topology of element, e.g. hypercube, simplex, ect 5699df49d7eSJed Brown /// * `ncomp` - Number of field components (1 for scalar fields) 5709df49d7eSJed Brown /// * `nnodes` - Total number of nodes 5719df49d7eSJed Brown /// * `nqpts` - Total number of quadrature points 5729df49d7eSJed Brown /// * `interp` - Row-major `(nqpts * nnodes)` matrix expressing the values of 5739df49d7eSJed Brown /// nodal basis functions at quadrature points 5749df49d7eSJed Brown /// * `grad` - Row-major `(nqpts * dim * nnodes)` matrix expressing 5759df49d7eSJed Brown /// derivatives of nodal basis functions at quadrature points 5769df49d7eSJed Brown /// * `qref` - Array of length `nqpts` holding the locations of quadrature 5779df49d7eSJed Brown /// points on the reference element `[-1, 1]` 5789df49d7eSJed Brown /// * `qweight` - Array of length `nqpts` holding the quadrature weights on 5799df49d7eSJed Brown /// the reference element 5809df49d7eSJed Brown /// 5819df49d7eSJed Brown /// ``` 5829df49d7eSJed Brown /// # use libceed::prelude::*; 583c68be7a2SJeremy L Thompson /// # fn main() -> Result<(), libceed::CeedError> { 5849df49d7eSJed Brown /// # let ceed = libceed::Ceed::default_init(); 5859df49d7eSJed Brown /// let interp = [ 5869df49d7eSJed Brown /// 0.12000000, 5879df49d7eSJed Brown /// 0.48000000, 5889df49d7eSJed Brown /// -0.12000000, 5899df49d7eSJed Brown /// 0.48000000, 5909df49d7eSJed Brown /// 0.16000000, 5919df49d7eSJed Brown /// -0.12000000, 5929df49d7eSJed Brown /// -0.12000000, 5939df49d7eSJed Brown /// 0.48000000, 5949df49d7eSJed Brown /// 0.12000000, 5959df49d7eSJed Brown /// 0.16000000, 5969df49d7eSJed Brown /// 0.48000000, 5979df49d7eSJed Brown /// -0.12000000, 5989df49d7eSJed Brown /// -0.11111111, 5999df49d7eSJed Brown /// 0.44444444, 6009df49d7eSJed Brown /// -0.11111111, 6019df49d7eSJed Brown /// 0.44444444, 6029df49d7eSJed Brown /// 0.44444444, 6039df49d7eSJed Brown /// -0.11111111, 6049df49d7eSJed Brown /// -0.12000000, 6059df49d7eSJed Brown /// 0.16000000, 6069df49d7eSJed Brown /// -0.12000000, 6079df49d7eSJed Brown /// 0.48000000, 6089df49d7eSJed Brown /// 0.48000000, 6099df49d7eSJed Brown /// 0.12000000, 6109df49d7eSJed Brown /// ]; 6119df49d7eSJed Brown /// let grad = [ 6129df49d7eSJed Brown /// -1.40000000, 6139df49d7eSJed Brown /// 1.60000000, 6149df49d7eSJed Brown /// -0.20000000, 6159df49d7eSJed Brown /// -0.80000000, 6169df49d7eSJed Brown /// 0.80000000, 6179df49d7eSJed Brown /// 0.00000000, 6189df49d7eSJed Brown /// 0.20000000, 6199df49d7eSJed Brown /// -1.60000000, 6209df49d7eSJed Brown /// 1.40000000, 6219df49d7eSJed Brown /// -0.80000000, 6229df49d7eSJed Brown /// 0.80000000, 6239df49d7eSJed Brown /// 0.00000000, 6249df49d7eSJed Brown /// -0.33333333, 6259df49d7eSJed Brown /// 0.00000000, 6269df49d7eSJed Brown /// 0.33333333, 6279df49d7eSJed Brown /// -1.33333333, 6289df49d7eSJed Brown /// 1.33333333, 6299df49d7eSJed Brown /// 0.00000000, 6309df49d7eSJed Brown /// 0.20000000, 6319df49d7eSJed Brown /// 0.00000000, 6329df49d7eSJed Brown /// -0.20000000, 6339df49d7eSJed Brown /// -2.40000000, 6349df49d7eSJed Brown /// 2.40000000, 6359df49d7eSJed Brown /// 0.00000000, 6369df49d7eSJed Brown /// -1.40000000, 6379df49d7eSJed Brown /// -0.80000000, 6389df49d7eSJed Brown /// 0.00000000, 6399df49d7eSJed Brown /// 1.60000000, 6409df49d7eSJed Brown /// 0.80000000, 6419df49d7eSJed Brown /// -0.20000000, 6429df49d7eSJed Brown /// 0.20000000, 6439df49d7eSJed Brown /// -2.40000000, 6449df49d7eSJed Brown /// 0.00000000, 6459df49d7eSJed Brown /// 0.00000000, 6469df49d7eSJed Brown /// 2.40000000, 6479df49d7eSJed Brown /// -0.20000000, 6489df49d7eSJed Brown /// -0.33333333, 6499df49d7eSJed Brown /// -1.33333333, 6509df49d7eSJed Brown /// 0.00000000, 6519df49d7eSJed Brown /// 0.00000000, 6529df49d7eSJed Brown /// 1.33333333, 6539df49d7eSJed Brown /// 0.33333333, 6549df49d7eSJed Brown /// 0.20000000, 6559df49d7eSJed Brown /// -0.80000000, 6569df49d7eSJed Brown /// 0.00000000, 6579df49d7eSJed Brown /// -1.60000000, 6589df49d7eSJed Brown /// 0.80000000, 6599df49d7eSJed Brown /// 1.40000000, 6609df49d7eSJed Brown /// ]; 6619df49d7eSJed Brown /// let qref = [ 6629df49d7eSJed Brown /// 0.20000000, 0.60000000, 0.33333333, 0.20000000, 0.20000000, 0.20000000, 0.33333333, 6639df49d7eSJed Brown /// 0.60000000, 6649df49d7eSJed Brown /// ]; 6659df49d7eSJed Brown /// let qweight = [0.26041667, 0.26041667, -0.28125000, 0.26041667]; 666c68be7a2SJeremy L Thompson /// let b = ceed.basis_H1( 6679df49d7eSJed Brown /// ElemTopology::Triangle, 6689df49d7eSJed Brown /// 1, 6699df49d7eSJed Brown /// 6, 6709df49d7eSJed Brown /// 4, 6719df49d7eSJed Brown /// &interp, 6729df49d7eSJed Brown /// &grad, 6739df49d7eSJed Brown /// &qref, 6749df49d7eSJed Brown /// &qweight, 675c68be7a2SJeremy L Thompson /// )?; 676c68be7a2SJeremy L Thompson /// # Ok(()) 677c68be7a2SJeremy L Thompson /// # } 6789df49d7eSJed Brown /// ``` 6799df49d7eSJed Brown pub fn basis_H1( 6809df49d7eSJed Brown &self, 6819df49d7eSJed Brown topo: ElemTopology, 6829df49d7eSJed Brown ncomp: usize, 6839df49d7eSJed Brown nnodes: usize, 6849df49d7eSJed Brown nqpts: usize, 68580a9ef05SNatalie Beams interp: &[crate::Scalar], 68680a9ef05SNatalie Beams grad: &[crate::Scalar], 68780a9ef05SNatalie Beams qref: &[crate::Scalar], 68880a9ef05SNatalie Beams qweight: &[crate::Scalar], 6899df49d7eSJed Brown ) -> Result<Basis> { 6909df49d7eSJed Brown Basis::create_H1( 6919df49d7eSJed Brown self, topo, ncomp, nnodes, nqpts, interp, grad, qref, qweight, 6929df49d7eSJed Brown ) 6939df49d7eSJed Brown } 6949df49d7eSJed Brown 6959df49d7eSJed Brown /// Returns a CeedQFunction for evaluating interior (volumetric) terms 6969df49d7eSJed Brown /// 6979df49d7eSJed Brown /// # arguments 6989df49d7eSJed Brown /// 6999df49d7eSJed Brown /// * `vlength` - Vector length. Caller must ensure that number of 7009df49d7eSJed Brown /// quadrature points is a multiple of vlength. 7019df49d7eSJed Brown /// * `f` - Boxed closure to evaluate action at quadrature points. 7029df49d7eSJed Brown /// 7039df49d7eSJed Brown /// ``` 7049df49d7eSJed Brown /// # use libceed::prelude::*; 705c68be7a2SJeremy L Thompson /// # fn main() -> Result<(), libceed::CeedError> { 7069df49d7eSJed Brown /// # let ceed = libceed::Ceed::default_init(); 7079df49d7eSJed Brown /// let mut user_f = |[u, weights, ..]: QFunctionInputs, [v, ..]: QFunctionOutputs| { 7089df49d7eSJed Brown /// // Iterate over quadrature points 7099df49d7eSJed Brown /// v.iter_mut() 7109df49d7eSJed Brown /// .zip(u.iter().zip(weights.iter())) 7119df49d7eSJed Brown /// .for_each(|(v, (u, w))| *v = u * w); 7129df49d7eSJed Brown /// 7139df49d7eSJed Brown /// // Return clean error code 7149df49d7eSJed Brown /// 0 7159df49d7eSJed Brown /// }; 7169df49d7eSJed Brown /// 717c68be7a2SJeremy L Thompson /// let qf = ceed.q_function_interior(1, Box::new(user_f))?; 718c68be7a2SJeremy L Thompson /// # Ok(()) 719c68be7a2SJeremy L Thompson /// # } 7209df49d7eSJed Brown /// ``` 7219df49d7eSJed Brown pub fn q_function_interior( 7229df49d7eSJed Brown &self, 7239df49d7eSJed Brown vlength: usize, 7249df49d7eSJed Brown f: Box<qfunction::QFunctionUserClosure>, 7259df49d7eSJed Brown ) -> Result<QFunction> { 7269df49d7eSJed Brown QFunction::create(self, vlength, f) 7279df49d7eSJed Brown } 7289df49d7eSJed Brown 7299df49d7eSJed Brown /// Returns a CeedQFunction for evaluating interior (volumetric) terms 7309df49d7eSJed Brown /// created by name 7319df49d7eSJed Brown /// 7329df49d7eSJed Brown /// ``` 7339df49d7eSJed Brown /// # use libceed::prelude::*; 734c68be7a2SJeremy L Thompson /// # fn main() -> Result<(), libceed::CeedError> { 7359df49d7eSJed Brown /// # let ceed = libceed::Ceed::default_init(); 736c68be7a2SJeremy L Thompson /// let qf = ceed.q_function_interior_by_name("Mass1DBuild")?; 737c68be7a2SJeremy L Thompson /// # Ok(()) 738c68be7a2SJeremy L Thompson /// # } 7399df49d7eSJed Brown /// ``` 7409df49d7eSJed Brown pub fn q_function_interior_by_name(&self, name: &str) -> Result<QFunctionByName> { 7419df49d7eSJed Brown QFunctionByName::create(self, name) 7429df49d7eSJed Brown } 7439df49d7eSJed Brown 7449df49d7eSJed Brown /// Returns a Operator and associate a QFunction. A Basis and 7459df49d7eSJed Brown /// ElemRestriction can be associated with QFunction fields with 7469df49d7eSJed Brown /// set_field(). 7479df49d7eSJed Brown /// 7489df49d7eSJed Brown /// * `qf` - QFunction defining the action of the operator at quadrature 7499df49d7eSJed Brown /// points 7509df49d7eSJed Brown /// * `dqf` - QFunction defining the action of the Jacobian of the qf (or 7519df49d7eSJed Brown /// qfunction_none) 7529df49d7eSJed Brown /// * `dqfT` - QFunction defining the action of the transpose of the 7539df49d7eSJed Brown /// Jacobian of the qf (or qfunction_none) 7549df49d7eSJed Brown /// 7559df49d7eSJed Brown /// ``` 7569df49d7eSJed Brown /// # use libceed::prelude::*; 757c68be7a2SJeremy L Thompson /// # fn main() -> Result<(), libceed::CeedError> { 7589df49d7eSJed Brown /// # let ceed = libceed::Ceed::default_init(); 759c68be7a2SJeremy L Thompson /// let qf = ceed.q_function_interior_by_name("Mass1DBuild")?; 760c68be7a2SJeremy L Thompson /// let op = ceed.operator(&qf, QFunctionOpt::None, QFunctionOpt::None)?; 761c68be7a2SJeremy L Thompson /// # Ok(()) 762c68be7a2SJeremy L Thompson /// # } 7639df49d7eSJed Brown /// ``` 7649df49d7eSJed Brown pub fn operator<'b>( 7659df49d7eSJed Brown &self, 7669df49d7eSJed Brown qf: impl Into<QFunctionOpt<'b>>, 7679df49d7eSJed Brown dqf: impl Into<QFunctionOpt<'b>>, 7689df49d7eSJed Brown dqfT: impl Into<QFunctionOpt<'b>>, 7699df49d7eSJed Brown ) -> Result<Operator> { 7709df49d7eSJed Brown Operator::create(self, qf, dqf, dqfT) 7719df49d7eSJed Brown } 7729df49d7eSJed Brown 7739df49d7eSJed Brown /// Returns an Operator that composes the action of several Operators 7749df49d7eSJed Brown /// 7759df49d7eSJed Brown /// ``` 7769df49d7eSJed Brown /// # use libceed::prelude::*; 777c68be7a2SJeremy L Thompson /// # fn main() -> Result<(), libceed::CeedError> { 7789df49d7eSJed Brown /// # let ceed = libceed::Ceed::default_init(); 779c68be7a2SJeremy L Thompson /// let op = ceed.composite_operator()?; 780c68be7a2SJeremy L Thompson /// # Ok(()) 781c68be7a2SJeremy L Thompson /// # } 7829df49d7eSJed Brown /// ``` 7839df49d7eSJed Brown pub fn composite_operator(&self) -> Result<CompositeOperator> { 7849df49d7eSJed Brown CompositeOperator::create(self) 7859df49d7eSJed Brown } 7869df49d7eSJed Brown } 7879df49d7eSJed Brown 7889df49d7eSJed Brown // ----------------------------------------------------------------------------- 7899df49d7eSJed Brown // Tests 7909df49d7eSJed Brown // ----------------------------------------------------------------------------- 7919df49d7eSJed Brown #[cfg(test)] 7929df49d7eSJed Brown mod tests { 7939df49d7eSJed Brown use super::*; 7949df49d7eSJed Brown 7959df49d7eSJed Brown fn ceed_t501() -> Result<i32> { 7969df49d7eSJed Brown let resource = "/cpu/self/ref/blocked"; 7979df49d7eSJed Brown let ceed = Ceed::init(resource); 7989df49d7eSJed Brown let nelem = 4; 7999df49d7eSJed Brown let p = 3; 8009df49d7eSJed Brown let q = 4; 8019df49d7eSJed Brown let ndofs = p * nelem - nelem + 1; 8029df49d7eSJed Brown 8039df49d7eSJed Brown // Vectors 8049df49d7eSJed Brown let x = ceed.vector_from_slice(&[-1., -0.5, 0.0, 0.5, 1.0])?; 8059df49d7eSJed Brown let mut qdata = ceed.vector(nelem * q)?; 8069df49d7eSJed Brown qdata.set_value(0.0)?; 8079df49d7eSJed Brown let mut u = ceed.vector(ndofs)?; 8089df49d7eSJed Brown u.set_value(1.0)?; 8099df49d7eSJed Brown let mut v = ceed.vector(ndofs)?; 8109df49d7eSJed Brown v.set_value(0.0)?; 8119df49d7eSJed Brown 8129df49d7eSJed Brown // Restrictions 8139df49d7eSJed Brown let mut indx: Vec<i32> = vec![0; 2 * nelem]; 8149df49d7eSJed Brown for i in 0..nelem { 8159df49d7eSJed Brown indx[2 * i + 0] = i as i32; 8169df49d7eSJed Brown indx[2 * i + 1] = (i + 1) as i32; 8179df49d7eSJed Brown } 8189df49d7eSJed Brown let rx = ceed.elem_restriction(nelem, 2, 1, 1, nelem + 1, MemType::Host, &indx)?; 8199df49d7eSJed Brown let mut indu: Vec<i32> = vec![0; p * nelem]; 8209df49d7eSJed Brown for i in 0..nelem { 8219df49d7eSJed Brown indu[p * i + 0] = i as i32; 8229df49d7eSJed Brown indu[p * i + 1] = (i + 1) as i32; 8239df49d7eSJed Brown indu[p * i + 2] = (i + 2) as i32; 8249df49d7eSJed Brown } 8259df49d7eSJed Brown let ru = ceed.elem_restriction(nelem, 3, 1, 1, ndofs, MemType::Host, &indu)?; 8269df49d7eSJed Brown let strides: [i32; 3] = [1, q as i32, q as i32]; 8279df49d7eSJed Brown let rq = ceed.strided_elem_restriction(nelem, q, 1, q * nelem, strides)?; 8289df49d7eSJed Brown 8299df49d7eSJed Brown // Bases 8309df49d7eSJed Brown let bx = ceed.basis_tensor_H1_Lagrange(1, 1, 2, q, QuadMode::Gauss)?; 8319df49d7eSJed Brown let bu = ceed.basis_tensor_H1_Lagrange(1, 1, p, q, QuadMode::Gauss)?; 8329df49d7eSJed Brown 8339df49d7eSJed Brown // Build quadrature data 8349df49d7eSJed Brown let qf_build = ceed.q_function_interior_by_name("Mass1DBuild")?; 8359df49d7eSJed Brown ceed.operator(&qf_build, QFunctionOpt::None, QFunctionOpt::None)? 8369df49d7eSJed Brown .field("dx", &rx, &bx, VectorOpt::Active)? 8379df49d7eSJed Brown .field("weights", ElemRestrictionOpt::None, &bx, VectorOpt::None)? 8389df49d7eSJed Brown .field("qdata", &rq, BasisOpt::Collocated, VectorOpt::Active)? 8399df49d7eSJed Brown .apply(&x, &mut qdata)?; 8409df49d7eSJed Brown 8419df49d7eSJed Brown // Mass operator 8429df49d7eSJed Brown let qf_mass = ceed.q_function_interior_by_name("MassApply")?; 8439df49d7eSJed Brown let op_mass = ceed 8449df49d7eSJed Brown .operator(&qf_mass, QFunctionOpt::None, QFunctionOpt::None)? 8459df49d7eSJed Brown .field("u", &ru, &bu, VectorOpt::Active)? 8469df49d7eSJed Brown .field("qdata", &rq, BasisOpt::Collocated, &qdata)? 8479df49d7eSJed Brown .field("v", &ru, &bu, VectorOpt::Active)?; 8489df49d7eSJed Brown 8499df49d7eSJed Brown v.set_value(0.0)?; 8509df49d7eSJed Brown op_mass.apply(&u, &mut v)?; 8519df49d7eSJed Brown 8529df49d7eSJed Brown // Check 85380a9ef05SNatalie Beams let sum: Scalar = v.view().iter().sum(); 8549df49d7eSJed Brown assert!( 8559df49d7eSJed Brown (sum - 2.0).abs() < 1e-15, 8569df49d7eSJed Brown "Incorrect interval length computed" 8579df49d7eSJed Brown ); 8589df49d7eSJed Brown Ok(0) 8599df49d7eSJed Brown } 8609df49d7eSJed Brown 8619df49d7eSJed Brown #[test] 8629df49d7eSJed Brown fn test_ceed_t501() { 8639df49d7eSJed Brown assert!(ceed_t501().is_ok()); 8649df49d7eSJed Brown } 8659df49d7eSJed Brown } 8669df49d7eSJed Brown 8679df49d7eSJed Brown // ----------------------------------------------------------------------------- 868