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