1 // Copyright (c) 2017-2025, Lawrence Livermore National Security, LLC and other CEED contributors. 2 // All Rights Reserved. See the top-level LICENSE and NOTICE files for details. 3 // 4 // SPDX-License-Identifier: BSD-2-Clause 5 // 6 // This file is part of CEED: http://github.com/ceed 7 8 // Fenced `rust` code blocks included from README.md are executed as part of doctests. 9 #![doc = include_str!("../README.md")] 10 // ----------------------------------------------------------------------------- 11 // Exceptions 12 // ----------------------------------------------------------------------------- 13 #![allow(non_snake_case)] 14 15 // ----------------------------------------------------------------------------- 16 // Crate prelude 17 // ----------------------------------------------------------------------------- 18 use crate::prelude::*; 19 use std::sync::Once; 20 21 pub mod prelude { 22 pub(crate) use libceed_sys::bind_ceed; 23 pub(crate) use std::convert::TryFrom; 24 pub(crate) use std::ffi::{CStr, CString}; 25 pub(crate) use std::fmt; 26 pub(crate) use std::marker::PhantomData; 27 } 28 29 // ----------------------------------------------------------------------------- 30 // Modules 31 // ----------------------------------------------------------------------------- 32 pub mod basis; 33 pub mod elem_restriction; 34 pub mod operator; 35 pub mod qfunction; 36 pub mod vector; 37 38 // ----------------------------------------------------------------------------- 39 // Typedef for scalar 40 // ----------------------------------------------------------------------------- 41 pub type Scalar = bind_ceed::CeedScalar; 42 43 // ----------------------------------------------------------------------------- 44 // Constants for library 45 // ----------------------------------------------------------------------------- 46 const MAX_BUFFER_LENGTH: usize = 4096; 47 pub const MAX_QFUNCTION_FIELDS: usize = 16; 48 pub const CEED_STRIDES_BACKEND: [i32; 3] = [0; 3]; 49 pub const EPSILON: crate::Scalar = bind_ceed::CEED_EPSILON as crate::Scalar; 50 51 // ----------------------------------------------------------------------------- 52 // Enums for libCEED 53 // ----------------------------------------------------------------------------- 54 #[derive(Clone, Copy, PartialEq, Eq)] 55 /// Many Ceed interfaces take or return pointers to memory. This enum is used to 56 /// specify where the memory being provided or requested must reside. 57 pub enum MemType { 58 Host = bind_ceed::CeedMemType_CEED_MEM_HOST as isize, 59 Device = bind_ceed::CeedMemType_CEED_MEM_DEVICE as isize, 60 } 61 62 #[derive(Clone, Copy, PartialEq, Eq)] 63 // OwnPointer will not be used by user but is included for internal use 64 #[allow(dead_code)] 65 /// Conveys ownership status of arrays passed to Ceed interfaces. 66 pub(crate) enum CopyMode { 67 CopyValues = bind_ceed::CeedCopyMode_CEED_COPY_VALUES as isize, 68 UsePointer = bind_ceed::CeedCopyMode_CEED_USE_POINTER as isize, 69 OwnPointer = bind_ceed::CeedCopyMode_CEED_OWN_POINTER as isize, 70 } 71 72 #[derive(Clone, Copy, PartialEq, Eq)] 73 /// Denotes type of vector norm to be computed 74 pub enum NormType { 75 One = bind_ceed::CeedNormType_CEED_NORM_1 as isize, 76 Two = bind_ceed::CeedNormType_CEED_NORM_2 as isize, 77 Max = bind_ceed::CeedNormType_CEED_NORM_MAX as isize, 78 } 79 80 #[derive(Clone, Copy, PartialEq, Eq)] 81 /// Denotes whether a linear transformation or its transpose should be applied 82 pub enum TransposeMode { 83 NoTranspose = bind_ceed::CeedTransposeMode_CEED_NOTRANSPOSE as isize, 84 Transpose = bind_ceed::CeedTransposeMode_CEED_TRANSPOSE as isize, 85 } 86 87 #[derive(Clone, Copy, PartialEq, Eq)] 88 /// Type of quadrature; also used for location of nodes 89 pub enum QuadMode { 90 Gauss = bind_ceed::CeedQuadMode_CEED_GAUSS as isize, 91 GaussLobatto = bind_ceed::CeedQuadMode_CEED_GAUSS_LOBATTO as isize, 92 } 93 94 #[derive(Clone, Copy, PartialEq, Eq)] 95 /// Type of basis shape to create non-tensor H1 element basis 96 pub enum ElemTopology { 97 Line = bind_ceed::CeedElemTopology_CEED_TOPOLOGY_LINE as isize, 98 Triangle = bind_ceed::CeedElemTopology_CEED_TOPOLOGY_TRIANGLE as isize, 99 Quad = bind_ceed::CeedElemTopology_CEED_TOPOLOGY_QUAD as isize, 100 Tet = bind_ceed::CeedElemTopology_CEED_TOPOLOGY_TET as isize, 101 Pyramid = bind_ceed::CeedElemTopology_CEED_TOPOLOGY_PYRAMID as isize, 102 Prism = bind_ceed::CeedElemTopology_CEED_TOPOLOGY_PRISM as isize, 103 Hex = bind_ceed::CeedElemTopology_CEED_TOPOLOGY_HEX as isize, 104 } 105 106 #[derive(Clone, Copy, Debug, PartialEq, Eq)] 107 /// Basis evaluation mode 108 pub enum EvalMode { 109 None = bind_ceed::CeedEvalMode_CEED_EVAL_NONE as isize, 110 Interp = bind_ceed::CeedEvalMode_CEED_EVAL_INTERP as isize, 111 Grad = bind_ceed::CeedEvalMode_CEED_EVAL_GRAD as isize, 112 Div = bind_ceed::CeedEvalMode_CEED_EVAL_DIV as isize, 113 Curl = bind_ceed::CeedEvalMode_CEED_EVAL_CURL as isize, 114 Weight = bind_ceed::CeedEvalMode_CEED_EVAL_WEIGHT as isize, 115 } 116 impl EvalMode { 117 pub(crate) fn from_u32(value: u32) -> EvalMode { 118 match value { 119 bind_ceed::CeedEvalMode_CEED_EVAL_NONE => EvalMode::None, 120 bind_ceed::CeedEvalMode_CEED_EVAL_INTERP => EvalMode::Interp, 121 bind_ceed::CeedEvalMode_CEED_EVAL_GRAD => EvalMode::Grad, 122 bind_ceed::CeedEvalMode_CEED_EVAL_DIV => EvalMode::Div, 123 bind_ceed::CeedEvalMode_CEED_EVAL_CURL => EvalMode::Curl, 124 bind_ceed::CeedEvalMode_CEED_EVAL_WEIGHT => EvalMode::Weight, 125 _ => panic!("Unknown value: {}", value), 126 } 127 } 128 } 129 130 // ----------------------------------------------------------------------------- 131 // Ceed error 132 // ----------------------------------------------------------------------------- 133 pub type Result<T> = std::result::Result<T, Error>; 134 135 /// libCEED error messages - returning an Error without painc!ing indicates 136 /// the function call failed but the data is not corrupted 137 #[derive(Debug)] 138 pub struct Error { 139 pub message: String, 140 } 141 142 impl fmt::Display for Error { 143 fn fmt(&self, f: &mut fmt::Formatter) -> fmt::Result { 144 write!(f, "{}", self.message) 145 } 146 } 147 148 // ----------------------------------------------------------------------------- 149 // Internal crate contents 150 // ----------------------------------------------------------------------------- 151 pub use crate::{ 152 basis::{Basis, BasisOpt}, 153 elem_restriction::{ElemRestriction, ElemRestrictionOpt}, 154 operator::{CompositeOperator, Operator, OperatorField}, 155 qfunction::{ 156 QFunction, QFunctionByName, QFunctionField, QFunctionInputs, QFunctionOpt, QFunctionOutputs, 157 }, 158 vector::{Vector, VectorOpt, VectorSliceWrapper}, 159 }; 160 161 // ----------------------------------------------------------------------------- 162 // Internal error checker 163 // ----------------------------------------------------------------------------- 164 #[doc(hidden)] 165 pub(crate) fn check_error<F>(ceed_ptr: F, ierr: i32) -> Result<i32> 166 where 167 F: FnOnce() -> bind_ceed::Ceed, 168 { 169 // Return early if code is clean 170 if ierr == bind_ceed::CeedErrorType_CEED_ERROR_SUCCESS { 171 return Ok(ierr); 172 } 173 // Retrieve error message 174 let mut ptr: *const std::os::raw::c_char = std::ptr::null_mut(); 175 let c_str = unsafe { 176 bind_ceed::CeedGetErrorMessage(ceed_ptr(), &mut ptr); 177 std::ffi::CStr::from_ptr(ptr) 178 }; 179 let message = c_str.to_string_lossy().to_string(); 180 // Panic if negative code, otherwise return error 181 if ierr < bind_ceed::CeedErrorType_CEED_ERROR_SUCCESS { 182 panic!("{}", message); 183 } 184 Err(Error { message }) 185 } 186 187 // ----------------------------------------------------------------------------- 188 // Ceed error handler 189 // ----------------------------------------------------------------------------- 190 pub enum ErrorHandler { 191 ErrorAbort, 192 ErrorExit, 193 ErrorReturn, 194 ErrorStore, 195 } 196 197 // ----------------------------------------------------------------------------- 198 // Ceed context wrapper 199 // ----------------------------------------------------------------------------- 200 /// A Ceed is a library context representing control of a logical hardware 201 /// resource. 202 #[derive(Debug)] 203 pub struct Ceed { 204 ptr: bind_ceed::Ceed, 205 } 206 207 // ----------------------------------------------------------------------------- 208 // Destructor 209 // ----------------------------------------------------------------------------- 210 impl Drop for Ceed { 211 fn drop(&mut self) { 212 unsafe { 213 bind_ceed::CeedDestroy(&mut self.ptr); 214 } 215 } 216 } 217 218 // ----------------------------------------------------------------------------- 219 // Cloning 220 // ----------------------------------------------------------------------------- 221 impl Clone for Ceed { 222 /// Perform a shallow clone of a Ceed context 223 /// 224 /// ``` 225 /// let ceed = libceed::Ceed::init("/cpu/self/ref/serial"); 226 /// let ceed_clone = ceed.clone(); 227 /// 228 /// println!(" original:{} \n clone:{}", ceed, ceed_clone); 229 /// ``` 230 fn clone(&self) -> Self { 231 let mut ptr_clone = std::ptr::null_mut(); 232 self.check_error(unsafe { bind_ceed::CeedReferenceCopy(self.ptr, &mut ptr_clone) }) 233 .expect("failed to clone Ceed"); 234 Self { ptr: ptr_clone } 235 } 236 } 237 238 // ----------------------------------------------------------------------------- 239 // Display 240 // ----------------------------------------------------------------------------- 241 impl fmt::Display for Ceed { 242 /// View a Ceed 243 /// 244 /// ``` 245 /// let ceed = libceed::Ceed::init("/cpu/self/ref/serial"); 246 /// println!("{}", ceed); 247 /// ``` 248 fn fmt(&self, f: &mut fmt::Formatter) -> fmt::Result { 249 let mut ptr = std::ptr::null_mut(); 250 let mut sizeloc = crate::MAX_BUFFER_LENGTH; 251 let cstring = unsafe { 252 let file = bind_ceed::open_memstream(&mut ptr, &mut sizeloc); 253 bind_ceed::CeedView(self.ptr, file); 254 bind_ceed::fclose(file); 255 CString::from_raw(ptr) 256 }; 257 cstring.to_string_lossy().fmt(f) 258 } 259 } 260 261 static REGISTER: Once = Once::new(); 262 263 // ----------------------------------------------------------------------------- 264 // Object constructors 265 // ----------------------------------------------------------------------------- 266 impl Ceed { 267 #[cfg_attr(feature = "katexit", katexit::katexit)] 268 /// Returns a Ceed context initialized with the specified resource 269 /// 270 /// # arguments 271 /// 272 /// * `resource` - Resource to use, e.g., "/cpu/self" 273 /// 274 /// ``` 275 /// let ceed = libceed::Ceed::init("/cpu/self/ref/serial"); 276 /// ``` 277 pub fn init(resource: &str) -> Self { 278 Ceed::init_with_error_handler(resource, ErrorHandler::ErrorStore) 279 } 280 281 /// Returns a Ceed context initialized with the specified resource 282 /// 283 /// # arguments 284 /// 285 /// * `resource` - Resource to use, e.g., "/cpu/self" 286 /// 287 /// ``` 288 /// let ceed = libceed::Ceed::init_with_error_handler( 289 /// "/cpu/self/ref/serial", 290 /// libceed::ErrorHandler::ErrorAbort, 291 /// ); 292 /// ``` 293 pub fn init_with_error_handler(resource: &str, handler: ErrorHandler) -> Self { 294 REGISTER.call_once(|| unsafe { 295 bind_ceed::CeedRegisterAll(); 296 bind_ceed::CeedQFunctionRegisterAll(); 297 }); 298 299 // Convert to C string 300 let c_resource = CString::new(resource).expect("CString::new failed"); 301 302 // Get error handler pointer 303 let eh = match handler { 304 ErrorHandler::ErrorAbort => bind_ceed::CeedErrorAbort, 305 ErrorHandler::ErrorExit => bind_ceed::CeedErrorExit, 306 ErrorHandler::ErrorReturn => bind_ceed::CeedErrorReturn, 307 ErrorHandler::ErrorStore => bind_ceed::CeedErrorStore, 308 }; 309 310 // Call to libCEED 311 let mut ptr = std::ptr::null_mut(); 312 let mut ierr = unsafe { bind_ceed::CeedInit(c_resource.as_ptr() as *const i8, &mut ptr) }; 313 if ierr != 0 { 314 panic!("Error initializing backend resource: {}", resource) 315 } 316 ierr = unsafe { bind_ceed::CeedSetErrorHandler(ptr, Some(eh)) }; 317 let ceed = Ceed { ptr }; 318 ceed.check_error(ierr).unwrap(); 319 ceed 320 } 321 322 /// Default initializer for testing 323 #[doc(hidden)] 324 pub fn default_init() -> Self { 325 // Convert to C string 326 let resource = "/cpu/self/ref/serial"; 327 crate::Ceed::init(resource) 328 } 329 330 /// Internal error checker 331 #[doc(hidden)] 332 fn check_error(&self, ierr: i32) -> Result<i32> { 333 // Return early if code is clean 334 if ierr == bind_ceed::CeedErrorType_CEED_ERROR_SUCCESS { 335 return Ok(ierr); 336 } 337 // Retrieve error message 338 let mut ptr: *const std::os::raw::c_char = std::ptr::null_mut(); 339 let c_str = unsafe { 340 bind_ceed::CeedGetErrorMessage(self.ptr, &mut ptr); 341 std::ffi::CStr::from_ptr(ptr) 342 }; 343 let message = c_str.to_string_lossy().to_string(); 344 // Panic if negative code, otherwise return error 345 if ierr < bind_ceed::CeedErrorType_CEED_ERROR_SUCCESS { 346 panic!("{}", message); 347 } 348 Err(Error { message }) 349 } 350 351 /// Returns full resource name for a Ceed context 352 /// 353 /// ``` 354 /// let ceed = libceed::Ceed::init("/cpu/self/ref/serial"); 355 /// let resource = ceed.resource(); 356 /// 357 /// assert_eq!(resource, "/cpu/self/ref/serial".to_string()) 358 /// ``` 359 pub fn resource(&self) -> String { 360 let mut ptr: *const std::os::raw::c_char = std::ptr::null_mut(); 361 let c_str = unsafe { 362 bind_ceed::CeedGetResource(self.ptr, &mut ptr); 363 std::ffi::CStr::from_ptr(ptr) 364 }; 365 c_str.to_string_lossy().to_string() 366 } 367 368 /// Returns a Vector of the specified length (does not allocate memory) 369 /// 370 /// # arguments 371 /// 372 /// * `n` - Length of vector 373 /// 374 /// ``` 375 /// # use libceed::prelude::*; 376 /// # fn main() -> libceed::Result<()> { 377 /// # let ceed = libceed::Ceed::default_init(); 378 /// let vec = ceed.vector(10)?; 379 /// # Ok(()) 380 /// # } 381 /// ``` 382 pub fn vector<'a>(&self, n: usize) -> Result<Vector<'a>> { 383 Vector::create(self, n) 384 } 385 386 /// Create a Vector initialized with the data (copied) from a slice 387 /// 388 /// # arguments 389 /// 390 /// * `slice` - Slice containing data 391 /// 392 /// ``` 393 /// # use libceed::prelude::*; 394 /// # fn main() -> libceed::Result<()> { 395 /// # let ceed = libceed::Ceed::default_init(); 396 /// let vec = ceed.vector_from_slice(&[1., 2., 3.])?; 397 /// assert_eq!(vec.length(), 3); 398 /// # Ok(()) 399 /// # } 400 /// ``` 401 pub fn vector_from_slice<'a>(&self, slice: &[crate::Scalar]) -> Result<Vector<'a>> { 402 Vector::from_slice(self, slice) 403 } 404 405 /// Returns an ElemRestriction, $\mathcal{E}$, which extracts the degrees of 406 /// freedom for each element from the local vector into the element vector 407 /// or assembles contributions from each element in the element vector to 408 /// the local vector 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::*, MemType}; 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 an oriented ElemRestriction, $\mathcal{E}$, which extracts the 460 /// degrees of freedom for each element from the local vector into the 461 /// element vector or assembles contributions from each element in the 462 /// element vector to the local vector 463 /// 464 /// # arguments 465 /// 466 /// * `nelem` - Number of elements described in the offsets array 467 /// * `elemsize` - Size (number of "nodes") per element 468 /// * `ncomp` - Number of field components per interpolation node (1 469 /// for scalar fields) 470 /// * `compstride` - Stride between components for the same Lvector "node". 471 /// Data for node `i`, component `j`, element `k` can be 472 /// found in the Lvector at index 473 /// `offsets[i + k*elemsize] + j*compstride`. 474 /// * `lsize` - The size of the Lvector. This vector may be larger 475 /// than the elements and fields given by this 476 /// restriction. 477 /// * `mtype` - Memory type of the offsets array, see CeedMemType 478 /// * `offsets` - Array of shape `[nelem, elemsize]`. Row `i` holds the 479 /// ordered list of the offsets (into the input CeedVector) 480 /// for the unknowns corresponding to element `i`, where 481 /// `0 <= i < nelem`. All offsets must be in the range 482 /// `[0, lsize - 1]`. 483 /// * `orients` - Array of shape `[nelem, elemsize]`. Row `i` holds the 484 /// ordered list of the orientations for the unknowns 485 /// corresponding to element `i`, with bool `false` used 486 /// for positively oriented and `true` to flip the 487 /// orientation. 488 /// 489 /// ``` 490 /// # use libceed::{prelude::*, MemType}; 491 /// # fn main() -> libceed::Result<()> { 492 /// # let ceed = libceed::Ceed::default_init(); 493 /// let nelem = 3; 494 /// let mut ind: Vec<i32> = vec![0; 2 * nelem]; 495 /// let mut orients: Vec<bool> = vec![false; 2 * nelem]; 496 /// for i in 0..nelem { 497 /// ind[2 * i + 0] = i as i32; 498 /// ind[2 * i + 1] = (i + 1) as i32; 499 /// orients[2 * i + 0] = (i % 2) > 0; // flip the dofs on element 1, 3, ... 500 /// orients[2 * i + 1] = (i % 2) > 0; 501 /// } 502 /// let r = 503 /// ceed.oriented_elem_restriction(nelem, 2, 1, 1, nelem + 1, MemType::Host, &ind, &orients)?; 504 /// # Ok(()) 505 /// # } 506 /// ``` 507 pub fn oriented_elem_restriction<'a>( 508 &self, 509 nelem: usize, 510 elemsize: usize, 511 ncomp: usize, 512 compstride: usize, 513 lsize: usize, 514 mtype: MemType, 515 offsets: &[i32], 516 orients: &[bool], 517 ) -> Result<ElemRestriction<'a>> { 518 ElemRestriction::create_oriented( 519 self, nelem, elemsize, ncomp, compstride, lsize, mtype, offsets, orients, 520 ) 521 } 522 523 /// Returns a curl-oriented ElemRestriction, $\mathcal{E}$, which extracts 524 /// the degrees of freedom for each element from the local vector into 525 /// the element vector or assembles contributions from each element in 526 /// the element vector to the local vector 527 /// 528 /// # arguments 529 /// 530 /// * `nelem` - Number of elements described in the offsets array 531 /// * `elemsize` - Size (number of "nodes") per element 532 /// * `ncomp` - Number of field components per interpolation node (1 533 /// for scalar fields) 534 /// * `compstride` - Stride between components for the same Lvector "node". 535 /// Data for node `i`, component `j`, element `k` can be 536 /// found in the Lvector at index 537 /// `offsets[i + k*elemsize] + j*compstride`. 538 /// * `lsize` - The size of the Lvector. This vector may be larger 539 /// than the elements and fields given by this 540 /// restriction. 541 /// * `mtype` - Memory type of the offsets array, see CeedMemType 542 /// * `offsets` - Array of shape `[nelem, elemsize]`. Row `i` holds the 543 /// ordered list of the offsets (into the input CeedVector) 544 /// for the unknowns corresponding to element `i`, where 545 /// `0 <= i < nelem`. All offsets must be in the range 546 /// `[0, lsize - 1]`. 547 /// * `curlorients` - Array of shape `[nelem, 3 * elemsize]`. Row `i` holds 548 /// a row-major tridiagonal matrix (`curlorients[i, 0] = 549 /// curlorients[i, 3 * elemsize - 1] = 0`, where 550 /// `0 <= i < nelem`) which is applied to the element 551 /// unknowns upon restriction. 552 /// 553 /// ``` 554 /// # use libceed::{prelude::*, MemType}; 555 /// # fn main() -> libceed::Result<()> { 556 /// # let ceed = libceed::Ceed::default_init(); 557 /// let nelem = 3; 558 /// let mut ind: Vec<i32> = vec![0; 2 * nelem]; 559 /// let mut curlorients: Vec<i8> = vec![0; 3 * 2 * nelem]; 560 /// for i in 0..nelem { 561 /// ind[2 * i + 0] = i as i32; 562 /// ind[2 * i + 1] = (i + 1) as i32; 563 /// curlorients[3 * 2 * i] = 0 as i8; 564 /// curlorients[3 * 2 * (i + 1) - 1] = 0 as i8; 565 /// if (i % 2 > 0) { 566 /// // T = [0 -1] 567 /// // [-1 0] 568 /// curlorients[3 * 2 * i + 1] = 0 as i8; 569 /// curlorients[3 * 2 * i + 2] = -1 as i8; 570 /// curlorients[3 * 2 * i + 3] = -1 as i8; 571 /// curlorients[3 * 2 * i + 4] = 0 as i8; 572 /// } else { 573 /// // T = I 574 /// curlorients[3 * 2 * i + 1] = 1 as i8; 575 /// curlorients[3 * 2 * i + 2] = 0 as i8; 576 /// curlorients[3 * 2 * i + 3] = 0 as i8; 577 /// curlorients[3 * 2 * i + 4] = 1 as i8; 578 /// } 579 /// } 580 /// let r = ceed.curl_oriented_elem_restriction( 581 /// nelem, 582 /// 2, 583 /// 1, 584 /// 1, 585 /// nelem + 1, 586 /// MemType::Host, 587 /// &ind, 588 /// &curlorients, 589 /// )?; 590 /// # Ok(()) 591 /// # } 592 /// ``` 593 pub fn curl_oriented_elem_restriction<'a>( 594 &self, 595 nelem: usize, 596 elemsize: usize, 597 ncomp: usize, 598 compstride: usize, 599 lsize: usize, 600 mtype: MemType, 601 offsets: &[i32], 602 curlorients: &[i8], 603 ) -> Result<ElemRestriction<'a>> { 604 ElemRestriction::create_curl_oriented( 605 self, 606 nelem, 607 elemsize, 608 ncomp, 609 compstride, 610 lsize, 611 mtype, 612 offsets, 613 curlorients, 614 ) 615 } 616 617 /// Returns an ElemRestriction, $\mathcal{E}$, from an local vector to 618 /// an element vector where data can be indexed from the `strides` array 619 /// 620 /// # arguments 621 /// 622 /// * `nelem` - Number of elements described in the offsets array 623 /// * `elemsize` - Size (number of "nodes") per element 624 /// * `ncomp` - Number of field components per interpolation node (1 625 /// for scalar fields) 626 /// * `lsize` - The size of the Lvector. This vector may be larger 627 /// than the elements and fields given by this restriction. 628 /// * `strides` - Array for strides between `[nodes, components, elements]`. 629 /// Data for node `i`, component `j`, element `k` can be 630 /// found in the Lvector at index 631 /// `i*strides[0] + j*strides[1] + k*strides[2]`. 632 /// CEED_STRIDES_BACKEND may be used with vectors created 633 /// by a Ceed backend. 634 /// 635 /// ``` 636 /// # use libceed::prelude::*; 637 /// # fn main() -> libceed::Result<()> { 638 /// # let ceed = libceed::Ceed::default_init(); 639 /// let nelem = 3; 640 /// let strides: [i32; 3] = [1, 2, 2]; 641 /// let r = ceed.strided_elem_restriction(nelem, 2, 1, nelem * 2, strides)?; 642 /// # Ok(()) 643 /// # } 644 /// ``` 645 pub fn strided_elem_restriction<'a>( 646 &self, 647 nelem: usize, 648 elemsize: usize, 649 ncomp: usize, 650 lsize: usize, 651 strides: [i32; 3], 652 ) -> Result<ElemRestriction<'a>> { 653 ElemRestriction::create_strided(self, nelem, elemsize, ncomp, lsize, strides) 654 } 655 656 /// Returns an $H^1$ tensor-product Basis 657 /// 658 /// # arguments 659 /// 660 /// * `dim` - Topological dimension of element 661 /// * `ncomp` - Number of field components (1 for scalar fields) 662 /// * `P1d` - Number of Gauss-Lobatto nodes in one dimension. The 663 /// polynomial degree of the resulting `Q_k` element is 664 /// `k=P-1`. 665 /// * `Q1d` - Number of quadrature points in one dimension 666 /// * `interp1d` - Row-major `(Q1d * P1d)` matrix expressing the values of 667 /// nodal basis functions at quadrature points 668 /// * `grad1d` - Row-major `(Q1d * P1d)` matrix expressing derivatives of 669 /// nodal basis functions at quadrature points 670 /// * `qref1d` - Array of length `Q1d` holding the locations of quadrature 671 /// points on the 1D reference element `[-1, 1]` 672 /// * `qweight1d` - Array of length `Q1d` holding the quadrature weights on 673 /// the reference element 674 /// 675 /// ``` 676 /// # use libceed::prelude::*; 677 /// # fn main() -> libceed::Result<()> { 678 /// # let ceed = libceed::Ceed::default_init(); 679 /// let interp1d = [ 0.62994317, 0.47255875, -0.14950343, 0.04700152, 680 /// -0.07069480, 0.97297619, 0.13253993, -0.03482132, 681 /// -0.03482132, 0.13253993, 0.97297619, -0.07069480, 682 /// 0.04700152, -0.14950343, 0.47255875, 0.62994317]; 683 /// let grad1d = [-2.34183742, 2.78794489, -0.63510411, 0.18899664, 684 /// -0.51670214, -0.48795249, 1.33790510, -0.33325047, 685 // 0.33325047, -1.33790510, 0.48795249, 0.51670214, 686 /// -0.18899664, 0.63510411, -2.78794489, 2.34183742]; 687 /// let qref1d = [-0.86113631, -0.33998104, 0.33998104, 0.86113631]; 688 /// let qweight1d = [ 0.34785485, 0.65214515, 0.65214515, 0.34785485]; 689 /// let b = ceed. 690 /// basis_tensor_H1(2, 1, 4, 4, &interp1d, &grad1d, &qref1d, &qweight1d)?; 691 /// # Ok(()) 692 /// # } 693 /// ``` 694 pub fn basis_tensor_H1<'a>( 695 &self, 696 dim: usize, 697 ncomp: usize, 698 P1d: usize, 699 Q1d: usize, 700 interp1d: &[crate::Scalar], 701 grad1d: &[crate::Scalar], 702 qref1d: &[crate::Scalar], 703 qweight1d: &[crate::Scalar], 704 ) -> Result<Basis<'a>> { 705 Basis::create_tensor_H1( 706 self, dim, ncomp, P1d, Q1d, interp1d, grad1d, qref1d, qweight1d, 707 ) 708 } 709 710 /// Returns an $H^1$ Lagrange tensor-product Basis 711 /// 712 /// # arguments 713 /// 714 /// * `dim` - Topological dimension of element 715 /// * `ncomp` - Number of field components (1 for scalar fields) 716 /// * `P` - Number of Gauss-Lobatto nodes in one dimension. The 717 /// polynomial degree of the resulting `Q_k` element is `k=P-1`. 718 /// * `Q` - Number of quadrature points in one dimension 719 /// * `qmode` - Distribution of the `Q` quadrature points (affects order of 720 /// accuracy for the quadrature) 721 /// 722 /// ``` 723 /// # use libceed::{prelude::*, QuadMode}; 724 /// # fn main() -> libceed::Result<()> { 725 /// # let ceed = libceed::Ceed::default_init(); 726 /// let b = ceed.basis_tensor_H1_Lagrange(2, 1, 3, 4, QuadMode::Gauss)?; 727 /// # Ok(()) 728 /// # } 729 /// ``` 730 pub fn basis_tensor_H1_Lagrange<'a>( 731 &self, 732 dim: usize, 733 ncomp: usize, 734 P: usize, 735 Q: usize, 736 qmode: QuadMode, 737 ) -> Result<Basis<'a>> { 738 Basis::create_tensor_H1_Lagrange(self, dim, ncomp, P, Q, qmode) 739 } 740 741 /// Returns an $H-1$ Basis 742 /// 743 /// # arguments 744 /// 745 /// * `topo` - Topology of element, e.g. hypercube, simplex, ect 746 /// * `ncomp` - Number of field components (1 for scalar fields) 747 /// * `nnodes` - Total number of nodes 748 /// * `nqpts` - Total number of quadrature points 749 /// * `interp` - Row-major `(nqpts * nnodes)` matrix expressing the values of 750 /// nodal basis functions at quadrature points 751 /// * `grad` - Row-major `(dim * nqpts * nnodes)` matrix expressing 752 /// derivatives of nodal basis functions at quadrature points 753 /// * `qref` - Array of length `nqpts` holding the locations of quadrature 754 /// points on the reference element `[-1, 1]` 755 /// * `qweight` - Array of length `nqpts` holding the quadrature weights on 756 /// the reference element 757 /// 758 /// ``` 759 /// # use libceed::{prelude::*, ElemTopology}; 760 /// # fn main() -> libceed::Result<()> { 761 /// # let ceed = libceed::Ceed::default_init(); 762 /// let interp = [ 763 /// 0.12000000, 764 /// 0.48000000, 765 /// -0.12000000, 766 /// 0.48000000, 767 /// 0.16000000, 768 /// -0.12000000, 769 /// -0.12000000, 770 /// 0.48000000, 771 /// 0.12000000, 772 /// 0.16000000, 773 /// 0.48000000, 774 /// -0.12000000, 775 /// -0.11111111, 776 /// 0.44444444, 777 /// -0.11111111, 778 /// 0.44444444, 779 /// 0.44444444, 780 /// -0.11111111, 781 /// -0.12000000, 782 /// 0.16000000, 783 /// -0.12000000, 784 /// 0.48000000, 785 /// 0.48000000, 786 /// 0.12000000, 787 /// ]; 788 /// let grad = [ 789 /// -1.40000000, 790 /// 1.60000000, 791 /// -0.20000000, 792 /// -0.80000000, 793 /// 0.80000000, 794 /// 0.00000000, 795 /// 0.20000000, 796 /// -1.60000000, 797 /// 1.40000000, 798 /// -0.80000000, 799 /// 0.80000000, 800 /// 0.00000000, 801 /// -0.33333333, 802 /// 0.00000000, 803 /// 0.33333333, 804 /// -1.33333333, 805 /// 1.33333333, 806 /// 0.00000000, 807 /// 0.20000000, 808 /// 0.00000000, 809 /// -0.20000000, 810 /// -2.40000000, 811 /// 2.40000000, 812 /// 0.00000000, 813 /// -1.40000000, 814 /// -0.80000000, 815 /// 0.00000000, 816 /// 1.60000000, 817 /// 0.80000000, 818 /// -0.20000000, 819 /// 0.20000000, 820 /// -2.40000000, 821 /// 0.00000000, 822 /// 0.00000000, 823 /// 2.40000000, 824 /// -0.20000000, 825 /// -0.33333333, 826 /// -1.33333333, 827 /// 0.00000000, 828 /// 0.00000000, 829 /// 1.33333333, 830 /// 0.33333333, 831 /// 0.20000000, 832 /// -0.80000000, 833 /// 0.00000000, 834 /// -1.60000000, 835 /// 0.80000000, 836 /// 1.40000000, 837 /// ]; 838 /// let qref = [ 839 /// 0.20000000, 0.60000000, 0.33333333, 0.20000000, 0.20000000, 0.20000000, 0.33333333, 840 /// 0.60000000, 841 /// ]; 842 /// let qweight = [0.26041667, 0.26041667, -0.28125000, 0.26041667]; 843 /// let b = ceed.basis_H1( 844 /// ElemTopology::Triangle, 845 /// 1, 846 /// 6, 847 /// 4, 848 /// &interp, 849 /// &grad, 850 /// &qref, 851 /// &qweight, 852 /// )?; 853 /// # Ok(()) 854 /// # } 855 /// ``` 856 pub fn basis_H1<'a>( 857 &self, 858 topo: ElemTopology, 859 ncomp: usize, 860 nnodes: usize, 861 nqpts: usize, 862 interp: &[crate::Scalar], 863 grad: &[crate::Scalar], 864 qref: &[crate::Scalar], 865 qweight: &[crate::Scalar], 866 ) -> Result<Basis<'a>> { 867 Basis::create_H1( 868 self, topo, ncomp, nnodes, nqpts, interp, grad, qref, qweight, 869 ) 870 } 871 872 /// Returns an $H(div)$ Basis 873 /// 874 /// # arguments 875 /// 876 /// * `topo` - Topology of element, e.g. hypercube, simplex, ect 877 /// * `ncomp` - Number of field components (1 for scalar fields) 878 /// * `nnodes` - Total number of nodes 879 /// * `nqpts` - Total number of quadrature points 880 /// * `interp` - Row-major `(dim * nqpts * nnodes)` matrix expressing the 881 /// values of basis functions at quadrature points 882 /// * `div` - Row-major `(nqpts * nnodes)` matrix expressing the 883 /// divergence of basis functions at quadrature points 884 /// * `qref` - Array of length `nqpts` holding the locations of quadrature 885 /// points on the reference element `[-1, 1]` 886 /// * `qweight` - Array of length `nqpts` holding the quadrature weights on 887 /// the reference element 888 /// 889 /// ``` 890 /// # use libceed::{prelude::*, ElemTopology}; 891 /// # fn main() -> libceed::Result<()> { 892 /// # let ceed = libceed::Ceed::default_init(); 893 /// let interp = [ 894 /// 0.00000000, 895 /// -1.57735027, 896 /// 0.57735027, 897 /// 0.00000000, 898 /// 0.00000000, 899 /// -1.57735027, 900 /// 0.57735027, 901 /// 0.00000000, 902 /// 0.00000000, 903 /// -1.57735027, 904 /// 0.57735027, 905 /// 0.00000000, 906 /// 0.00000000, 907 /// -0.42264973, 908 /// -0.57735027, 909 /// 0.00000000, 910 /// 0.42264973, 911 /// 0.00000000, 912 /// 0.00000000, 913 /// 0.57735027, 914 /// 0.42264973, 915 /// 0.00000000, 916 /// 0.00000000, 917 /// 0.57735027, 918 /// 1.57735027, 919 /// 0.00000000, 920 /// 0.00000000, 921 /// -0.57735027, 922 /// 1.57735027, 923 /// 0.00000000, 924 /// 0.00000000, 925 /// -0.57735027, 926 /// ]; 927 /// let div = [ 928 /// -1.00000000, 929 /// 1.00000000, 930 /// -1.00000000, 931 /// 1.00000000, 932 /// -1.00000000, 933 /// 1.00000000, 934 /// -1.00000000, 935 /// 1.00000000, 936 /// -1.00000000, 937 /// 1.00000000, 938 /// -1.00000000, 939 /// 1.00000000, 940 /// -1.00000000, 941 /// 1.00000000, 942 /// -1.00000000, 943 /// 1.00000000, 944 /// ]; 945 /// let qref = [ 946 /// 0.57735026, 0.57735026, 0.57735026, 0.57735026, 0.57735026, 0.57735026, 0.57735026, 947 /// 0.57735026, 948 /// ]; 949 /// let qweight = [1.00000000, 1.00000000, 1.00000000, 1.00000000]; 950 /// let b = ceed.basis_Hdiv(ElemTopology::Quad, 1, 4, 4, &interp, &div, &qref, &qweight)?; 951 /// # Ok(()) 952 /// # } 953 /// ``` 954 pub fn basis_Hdiv<'a>( 955 &self, 956 topo: ElemTopology, 957 ncomp: usize, 958 nnodes: usize, 959 nqpts: usize, 960 interp: &[crate::Scalar], 961 div: &[crate::Scalar], 962 qref: &[crate::Scalar], 963 qweight: &[crate::Scalar], 964 ) -> Result<Basis<'a>> { 965 Basis::create_Hdiv(self, topo, ncomp, nnodes, nqpts, interp, div, qref, qweight) 966 } 967 968 /// Returns an $H(curl)$ Basis 969 /// 970 /// # arguments 971 /// 972 /// * `topo` - Topology of element, e.g. hypercube, simplex, ect 973 /// * `ncomp` - Number of field components (1 for scalar fields) 974 /// * `nnodes` - Total number of nodes 975 /// * `nqpts` - Total number of quadrature points 976 /// * `interp` - Row-major `(dim * nqpts * nnodes)` matrix expressing the 977 /// values of basis functions at quadrature points 978 /// * `curl` - Row-major `(curl_comp * nqpts * nnodes)`, `curl_comp = 1 if 979 /// dim < 3 else dim` matrix expressing the curl of basis 980 /// functions at quadrature points 981 /// * `qref` - Array of length `nqpts` holding the locations of quadrature 982 /// points on the reference element `[-1, 1]` 983 /// * `qweight` - Array of length `nqpts` holding the quadrature weights on 984 /// the reference element 985 /// 986 /// ``` 987 /// # use libceed::{prelude::*, ElemTopology}; 988 /// # fn main() -> libceed::Result<()> { 989 /// # let ceed = libceed::Ceed::default_init(); 990 /// let interp = [ 991 /// -0.20000000, 992 /// 0.20000000, 993 /// 0.80000000, 994 /// -0.20000000, 995 /// 0.20000000, 996 /// 0.80000000, 997 /// -0.33333333, 998 /// 0.33333333, 999 /// 0.66666667, 1000 /// -0.60000000, 1001 /// 0.60000000, 1002 /// 0.40000000, 1003 /// 0.20000000, 1004 /// 0.80000000, 1005 /// 0.20000000, 1006 /// 0.60000000, 1007 /// 0.40000000, 1008 /// 0.60000000, 1009 /// 0.33333333, 1010 /// 0.66666667, 1011 /// 0.33333333, 1012 /// 0.20000000, 1013 /// 0.80000000, 1014 /// 0.20000000, 1015 /// ]; 1016 /// let curl = [ 1017 /// 2.00000000, 1018 /// -2.00000000, 1019 /// -2.00000000, 1020 /// 2.00000000, 1021 /// -2.00000000, 1022 /// -2.00000000, 1023 /// 2.00000000, 1024 /// -2.00000000, 1025 /// -2.00000000, 1026 /// 2.00000000, 1027 /// -2.00000000, 1028 /// -2.00000000, 1029 /// ]; 1030 /// let qref = [ 1031 /// 0.20000000, 0.60000000, 0.33333333, 0.20000000, 0.20000000, 0.20000000, 0.33333333, 1032 /// 0.60000000, 1033 /// ]; 1034 /// let qweight = [0.26041667, 0.26041667, -0.28125000, 0.26041667]; 1035 /// let b = ceed.basis_Hcurl( 1036 /// ElemTopology::Triangle, 1037 /// 1, 1038 /// 3, 1039 /// 4, 1040 /// &interp, 1041 /// &curl, 1042 /// &qref, 1043 /// &qweight, 1044 /// )?; 1045 /// # Ok(()) 1046 /// # } 1047 /// ``` 1048 pub fn basis_Hcurl<'a>( 1049 &self, 1050 topo: ElemTopology, 1051 ncomp: usize, 1052 nnodes: usize, 1053 nqpts: usize, 1054 interp: &[crate::Scalar], 1055 curl: &[crate::Scalar], 1056 qref: &[crate::Scalar], 1057 qweight: &[crate::Scalar], 1058 ) -> Result<Basis<'a>> { 1059 Basis::create_Hcurl( 1060 self, topo, ncomp, nnodes, nqpts, interp, curl, qref, qweight, 1061 ) 1062 } 1063 1064 /// Returns a QFunction for evaluating interior (volumetric) terms 1065 /// of a weak form corresponding to the $L^2$ inner product 1066 /// 1067 /// $$ 1068 /// \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), 1069 /// $$ 1070 /// 1071 /// where $v \cdot f_0$ represents contraction over fields and $\nabla v : f_1$ 1072 /// represents contraction over both fields and spatial dimensions. 1073 /// 1074 /// # arguments 1075 /// 1076 /// * `vlength` - Vector length. Caller must ensure that number of 1077 /// quadrature points is a multiple of vlength. 1078 /// * `f` - Boxed closure to evaluate weak form at quadrature points. 1079 /// 1080 /// ``` 1081 /// # use libceed::{prelude::*, QFunctionInputs, QFunctionOutputs}; 1082 /// # fn main() -> libceed::Result<()> { 1083 /// # let ceed = libceed::Ceed::default_init(); 1084 /// let mut user_f = |[u, weights, ..]: QFunctionInputs, [v, ..]: QFunctionOutputs| { 1085 /// // Iterate over quadrature points 1086 /// v.iter_mut() 1087 /// .zip(u.iter().zip(weights.iter())) 1088 /// .for_each(|(v, (u, w))| *v = u * w); 1089 /// 1090 /// // Return clean error code 1091 /// 0 1092 /// }; 1093 /// 1094 /// let qf = ceed.q_function_interior(1, Box::new(user_f))?; 1095 /// # Ok(()) 1096 /// # } 1097 /// ``` 1098 pub fn q_function_interior<'a>( 1099 &self, 1100 vlength: usize, 1101 f: Box<qfunction::QFunctionUserClosure>, 1102 ) -> Result<QFunction<'a>> { 1103 QFunction::create(self, vlength, f) 1104 } 1105 1106 /// Returns a QFunction for evaluating interior (volumetric) terms 1107 /// created by name 1108 /// 1109 /// # arguments 1110 /// 1111 /// * `name` - name of QFunction from libCEED gallery 1112 /// 1113 /// ``` 1114 /// # use libceed::prelude::*; 1115 /// # fn main() -> libceed::Result<()> { 1116 /// # let ceed = libceed::Ceed::default_init(); 1117 /// let qf = ceed.q_function_interior_by_name("Mass1DBuild")?; 1118 /// # Ok(()) 1119 /// # } 1120 /// ``` 1121 pub fn q_function_interior_by_name<'a>(&self, name: &str) -> Result<QFunctionByName<'a>> { 1122 QFunctionByName::create(self, name) 1123 } 1124 1125 /// Returns an Operator and associate a QFunction. A Basis and 1126 /// ElemRestriction can be associated with QFunction fields via 1127 /// set_field(). 1128 /// 1129 /// # arguments 1130 /// 1131 /// * `qf` - QFunction defining the action of the operator at quadrature 1132 /// points 1133 /// * `dqf` - QFunction defining the action of the Jacobian of the qf (or 1134 /// qfunction_none) 1135 /// * `dqfT` - QFunction defining the action of the transpose of the 1136 /// Jacobian of the qf (or qfunction_none) 1137 /// 1138 /// ``` 1139 /// # use libceed::{prelude::*, QFunctionOpt}; 1140 /// # fn main() -> libceed::Result<()> { 1141 /// # let ceed = libceed::Ceed::default_init(); 1142 /// let qf = ceed.q_function_interior_by_name("Mass1DBuild")?; 1143 /// let op = ceed.operator(&qf, QFunctionOpt::None, QFunctionOpt::None)?; 1144 /// # Ok(()) 1145 /// # } 1146 /// ``` 1147 pub fn operator<'a, 'b>( 1148 &self, 1149 qf: impl Into<QFunctionOpt<'b>>, 1150 dqf: impl Into<QFunctionOpt<'b>>, 1151 dqfT: impl Into<QFunctionOpt<'b>>, 1152 ) -> Result<Operator<'a>> { 1153 Operator::create(self, qf, dqf, dqfT) 1154 } 1155 1156 /// Returns an Operator that composes the action of several Operators 1157 /// 1158 /// ``` 1159 /// # use libceed::prelude::*; 1160 /// # fn main() -> libceed::Result<()> { 1161 /// # let ceed = libceed::Ceed::default_init(); 1162 /// let op = ceed.composite_operator()?; 1163 /// # Ok(()) 1164 /// # } 1165 /// ``` 1166 pub fn composite_operator<'a>(&self) -> Result<CompositeOperator<'a>> { 1167 CompositeOperator::create(self) 1168 } 1169 } 1170 1171 // ----------------------------------------------------------------------------- 1172 // Tests 1173 // ----------------------------------------------------------------------------- 1174 #[cfg(test)] 1175 mod tests { 1176 use super::*; 1177 1178 fn ceed_t501() -> Result<()> { 1179 let resource = "/cpu/self/ref/blocked"; 1180 let ceed = Ceed::init(resource); 1181 let nelem = 4; 1182 let p = 3; 1183 let q = 4; 1184 let ndofs = p * nelem - nelem + 1; 1185 1186 // Vectors 1187 let x = ceed.vector_from_slice(&[-1., -0.5, 0.0, 0.5, 1.0])?; 1188 let mut qdata = ceed.vector(nelem * q)?; 1189 qdata.set_value(0.0)?; 1190 let mut u = ceed.vector(ndofs)?; 1191 u.set_value(1.0)?; 1192 let mut v = ceed.vector(ndofs)?; 1193 v.set_value(0.0)?; 1194 1195 // Restrictions 1196 let mut indx: Vec<i32> = vec![0; 2 * nelem]; 1197 for i in 0..nelem { 1198 for j in 0..2 { 1199 indx[2 * i + j] = (i + j) as i32; 1200 } 1201 } 1202 let rx = ceed.elem_restriction(nelem, 2, 1, 1, nelem + 1, MemType::Host, &indx)?; 1203 let mut indu: Vec<i32> = vec![0; p * nelem]; 1204 for i in 0..nelem { 1205 for j in 0..p { 1206 indu[p * i + j] = (i + j) as i32; 1207 } 1208 } 1209 let ru = ceed.elem_restriction(nelem, p, 1, 1, ndofs, MemType::Host, &indu)?; 1210 let strides: [i32; 3] = [1, q as i32, q as i32]; 1211 let rq = ceed.strided_elem_restriction(nelem, q, 1, q * nelem, strides)?; 1212 1213 // Bases 1214 let bx = ceed.basis_tensor_H1_Lagrange(1, 1, 2, q, QuadMode::Gauss)?; 1215 let bu = ceed.basis_tensor_H1_Lagrange(1, 1, p, q, QuadMode::Gauss)?; 1216 1217 // Build quadrature data 1218 let qf_build = ceed.q_function_interior_by_name("Mass1DBuild")?; 1219 ceed.operator(&qf_build, QFunctionOpt::None, QFunctionOpt::None)? 1220 .field("dx", &rx, &bx, VectorOpt::Active)? 1221 .field("weights", ElemRestrictionOpt::None, &bx, VectorOpt::None)? 1222 .field("qdata", &rq, BasisOpt::None, VectorOpt::Active)? 1223 .apply(&x, &mut qdata)?; 1224 1225 // Mass operator 1226 let qf_mass = ceed.q_function_interior_by_name("MassApply")?; 1227 let op_mass = ceed 1228 .operator(&qf_mass, QFunctionOpt::None, QFunctionOpt::None)? 1229 .field("u", &ru, &bu, VectorOpt::Active)? 1230 .field("qdata", &rq, BasisOpt::None, &qdata)? 1231 .field("v", &ru, &bu, VectorOpt::Active)? 1232 .check()?; 1233 1234 v.set_value(0.0)?; 1235 op_mass.apply(&u, &mut v)?; 1236 1237 // Check 1238 let sum: Scalar = v.view()?.iter().sum(); 1239 let error: Scalar = (sum - 2.0).abs(); 1240 assert!( 1241 error < 50.0 * EPSILON, 1242 "Incorrect interval length computed. Expected: 2.0, Found: {}, Error: {:.12e}", 1243 sum, 1244 error 1245 ); 1246 Ok(()) 1247 } 1248 1249 #[test] 1250 fn test_ceed_t501() { 1251 assert!(ceed_t501().is_ok()); 1252 } 1253 } 1254 1255 // ----------------------------------------------------------------------------- 1256