1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206 207 208 209 210 211 212 213 214 215 216 217 218 219 220 221 222 223 224 225 226 227 228 229 230 231 232 233 234 235 236 237 238 239 240 241 242 243 244 245 246 247 248 249 250 251 252 253 254 255 256 257 258 259 260 261 262 263 264 265 266 267 268 269 270 271 272 273 274 275 276 277 278 279 280 281 282 283 284 285 286 287 288 289 290 291 292 293 294 295 296 297 298 299 300 301 302 303 304 305 306 307 308 309 310 311 312 313 314 315 316 317 318 319 320 321 322 323 324 325 326 327 328 329 330 331 332 333 334 335 336 337 338 339 340 341 342 343 344 345 346 347 348 349 350 351 352 353 354 355 356 357 358 359 360 361 362 363 364 365 366 367 368 369 370 371 372 373 374 375 376 377 378 379 380 381 382 383 384 385 386 387 388 389 390 391 392 393 394 395 396 397 398 399 400 401 402 403 404 405 406 407 408 409 410 411 412 413 414 415 416 417 418 419 420 421 422 423 424 425 426 427 428 429 430 431 432 433 434 435 436 437 438 439 440 441 442 443 444 445 446 447 448 449 450 451 452 453 454 455 456 457 458 459 460 461 462 463 464 465 466 467 468 469 470 471 472 473 474 475 476 477 478 479 480 481 482 483 484 485 486 487 488 489 490 491 492 493 494 495 496 497 498 499 500 501 502 503 504 505 506 507 508 509 510 511 512 513 514 515 516 517 518 519 520 521 522 523 524 525 526 527 528 529 530 531 532 533 534 535
// Copyright (c) 2017-2022, Lawrence Livermore National Security, LLC and other CEED contributors.
// All Rights Reserved. See the top-level LICENSE and NOTICE files for details.
//
// SPDX-License-Identifier: BSD-2-Clause
//
// This file is part of CEED: http://github.com/ceed
//! A Ceed Basis defines the discrete finite element basis and associated
//! quadrature rule.
use crate::prelude::*;
// -----------------------------------------------------------------------------
// Basis option
// -----------------------------------------------------------------------------
#[derive(Debug)]
pub enum BasisOpt<'a> {
Some(&'a Basis<'a>),
None,
}
/// Construct a BasisOpt reference from a Basis reference
impl<'a> From<&'a Basis<'_>> for BasisOpt<'a> {
fn from(basis: &'a Basis) -> Self {
debug_assert!(basis.ptr != unsafe { bind_ceed::CEED_BASIS_NONE });
Self::Some(basis)
}
}
impl<'a> BasisOpt<'a> {
/// Transform a Rust libCEED BasisOpt into C libCEED CeedBasis
pub(crate) fn to_raw(self) -> bind_ceed::CeedBasis {
match self {
Self::Some(basis) => basis.ptr,
Self::None => unsafe { bind_ceed::CEED_BASIS_NONE },
}
}
/// Check if a BasisOpt is Some
///
/// ```
/// # use libceed::prelude::*;
/// # fn main() -> libceed::Result<()> {
/// # let ceed = libceed::Ceed::default_init();
/// let b = ceed.basis_tensor_H1_Lagrange(1, 2, 3, 4, QuadMode::Gauss)?;
/// let b_opt = BasisOpt::from(&b);
/// assert!(b_opt.is_some(), "Incorrect BasisOpt");
///
/// let b_opt = BasisOpt::None;
/// assert!(!b_opt.is_some(), "Incorrect BasisOpt");
/// # Ok(())
/// # }
/// ```
pub fn is_some(&self) -> bool {
match self {
Self::Some(_) => true,
Self::None => false,
}
}
/// Check if a BasisOpt is None
///
/// ```
/// # use libceed::prelude::*;
/// # fn main() -> libceed::Result<()> {
/// # let ceed = libceed::Ceed::default_init();
/// let b = ceed.basis_tensor_H1_Lagrange(1, 2, 3, 4, QuadMode::Gauss)?;
/// let b_opt = BasisOpt::from(&b);
/// assert!(!b_opt.is_none(), "Incorrect BasisOpt");
///
/// let b_opt = BasisOpt::None;
/// assert!(b_opt.is_none(), "Incorrect BasisOpt");
/// # Ok(())
/// # }
/// ```
pub fn is_none(&self) -> bool {
match self {
Self::Some(_) => false,
Self::None => true,
}
}
}
// -----------------------------------------------------------------------------
// Basis context wrapper
// -----------------------------------------------------------------------------
#[derive(Debug)]
pub struct Basis<'a> {
pub(crate) ptr: bind_ceed::CeedBasis,
_lifeline: PhantomData<&'a ()>,
}
// -----------------------------------------------------------------------------
// Destructor
// -----------------------------------------------------------------------------
impl<'a> Drop for Basis<'a> {
fn drop(&mut self) {
unsafe {
if self.ptr != bind_ceed::CEED_BASIS_NONE {
bind_ceed::CeedBasisDestroy(&mut self.ptr);
}
}
}
}
// -----------------------------------------------------------------------------
// Display
// -----------------------------------------------------------------------------
impl<'a> fmt::Display for Basis<'a> {
/// View a Basis
///
/// ```
/// # use libceed::prelude::*;
/// # fn main() -> libceed::Result<()> {
/// # let ceed = libceed::Ceed::default_init();
/// let b = ceed.basis_tensor_H1_Lagrange(1, 2, 3, 4, QuadMode::Gauss)?;
/// println!("{}", b);
/// # Ok(())
/// # }
/// ```
fn fmt(&self, f: &mut fmt::Formatter) -> fmt::Result {
let mut ptr = std::ptr::null_mut();
let mut sizeloc = crate::MAX_BUFFER_LENGTH;
let cstring = unsafe {
let file = bind_ceed::open_memstream(&mut ptr, &mut sizeloc);
bind_ceed::CeedBasisView(self.ptr, file);
bind_ceed::fclose(file);
CString::from_raw(ptr)
};
cstring.to_string_lossy().fmt(f)
}
}
// -----------------------------------------------------------------------------
// Implementations
// -----------------------------------------------------------------------------
impl<'a> Basis<'a> {
// Constructors
pub fn create_tensor_H1(
ceed: &crate::Ceed,
dim: usize,
ncomp: usize,
P1d: usize,
Q1d: usize,
interp1d: &[crate::Scalar],
grad1d: &[crate::Scalar],
qref1d: &[crate::Scalar],
qweight1d: &[crate::Scalar],
) -> crate::Result<Self> {
let mut ptr = std::ptr::null_mut();
let (dim, ncomp, P1d, Q1d) = (
i32::try_from(dim).unwrap(),
i32::try_from(ncomp).unwrap(),
i32::try_from(P1d).unwrap(),
i32::try_from(Q1d).unwrap(),
);
let ierr = unsafe {
bind_ceed::CeedBasisCreateTensorH1(
ceed.ptr,
dim,
ncomp,
P1d,
Q1d,
interp1d.as_ptr(),
grad1d.as_ptr(),
qref1d.as_ptr(),
qweight1d.as_ptr(),
&mut ptr,
)
};
ceed.check_error(ierr)?;
Ok(Self {
ptr,
_lifeline: PhantomData,
})
}
pub fn create_tensor_H1_Lagrange(
ceed: &crate::Ceed,
dim: usize,
ncomp: usize,
P: usize,
Q: usize,
qmode: crate::QuadMode,
) -> crate::Result<Self> {
let mut ptr = std::ptr::null_mut();
let (dim, ncomp, P, Q, qmode) = (
i32::try_from(dim).unwrap(),
i32::try_from(ncomp).unwrap(),
i32::try_from(P).unwrap(),
i32::try_from(Q).unwrap(),
qmode as bind_ceed::CeedQuadMode,
);
let ierr = unsafe {
bind_ceed::CeedBasisCreateTensorH1Lagrange(ceed.ptr, dim, ncomp, P, Q, qmode, &mut ptr)
};
ceed.check_error(ierr)?;
Ok(Self {
ptr,
_lifeline: PhantomData,
})
}
pub fn create_H1(
ceed: &crate::Ceed,
topo: crate::ElemTopology,
ncomp: usize,
nnodes: usize,
nqpts: usize,
interp: &[crate::Scalar],
grad: &[crate::Scalar],
qref: &[crate::Scalar],
qweight: &[crate::Scalar],
) -> crate::Result<Self> {
let mut ptr = std::ptr::null_mut();
let (topo, ncomp, nnodes, nqpts) = (
topo as bind_ceed::CeedElemTopology,
i32::try_from(ncomp).unwrap(),
i32::try_from(nnodes).unwrap(),
i32::try_from(nqpts).unwrap(),
);
let ierr = unsafe {
bind_ceed::CeedBasisCreateH1(
ceed.ptr,
topo,
ncomp,
nnodes,
nqpts,
interp.as_ptr(),
grad.as_ptr(),
qref.as_ptr(),
qweight.as_ptr(),
&mut ptr,
)
};
ceed.check_error(ierr)?;
Ok(Self {
ptr,
_lifeline: PhantomData,
})
}
pub fn create_Hdiv(
ceed: &crate::Ceed,
topo: crate::ElemTopology,
ncomp: usize,
nnodes: usize,
nqpts: usize,
interp: &[crate::Scalar],
div: &[crate::Scalar],
qref: &[crate::Scalar],
qweight: &[crate::Scalar],
) -> crate::Result<Self> {
let mut ptr = std::ptr::null_mut();
let (topo, ncomp, nnodes, nqpts) = (
topo as bind_ceed::CeedElemTopology,
i32::try_from(ncomp).unwrap(),
i32::try_from(nnodes).unwrap(),
i32::try_from(nqpts).unwrap(),
);
let ierr = unsafe {
bind_ceed::CeedBasisCreateHdiv(
ceed.ptr,
topo,
ncomp,
nnodes,
nqpts,
interp.as_ptr(),
div.as_ptr(),
qref.as_ptr(),
qweight.as_ptr(),
&mut ptr,
)
};
ceed.check_error(ierr)?;
Ok(Self {
ptr,
_lifeline: PhantomData,
})
}
pub fn create_Hcurl(
ceed: &crate::Ceed,
topo: crate::ElemTopology,
ncomp: usize,
nnodes: usize,
nqpts: usize,
interp: &[crate::Scalar],
curl: &[crate::Scalar],
qref: &[crate::Scalar],
qweight: &[crate::Scalar],
) -> crate::Result<Self> {
let mut ptr = std::ptr::null_mut();
let (topo, ncomp, nnodes, nqpts) = (
topo as bind_ceed::CeedElemTopology,
i32::try_from(ncomp).unwrap(),
i32::try_from(nnodes).unwrap(),
i32::try_from(nqpts).unwrap(),
);
let ierr = unsafe {
bind_ceed::CeedBasisCreateHcurl(
ceed.ptr,
topo,
ncomp,
nnodes,
nqpts,
interp.as_ptr(),
curl.as_ptr(),
qref.as_ptr(),
qweight.as_ptr(),
&mut ptr,
)
};
ceed.check_error(ierr)?;
Ok(Self {
ptr,
_lifeline: PhantomData,
})
}
// Error handling
#[doc(hidden)]
fn check_error(&self, ierr: i32) -> crate::Result<i32> {
let mut ptr = std::ptr::null_mut();
unsafe {
bind_ceed::CeedBasisGetCeed(self.ptr, &mut ptr);
}
crate::check_error(ptr, ierr)
}
/// Apply basis evaluation from nodes to quadrature points or vice versa
///
/// * `nelem` - The number of elements to apply the basis evaluation to
/// * `tmode` - `TrasposeMode::NoTranspose` to evaluate from nodes to
/// quadrature points, `TransposeMode::Transpose` to apply the
/// transpose, mapping from quadrature points to nodes
/// * `emode` - `EvalMode::None` to use values directly, `EvalMode::Interp`
/// to use interpolated values, `EvalMode::Grad` to use
/// gradients, `EvalMode::Weight` to use quadrature weights
/// * `u` - Input Vector
/// * `v` - Output Vector
///
/// ```
/// # use libceed::prelude::*;
/// # fn main() -> libceed::Result<()> {
/// # let ceed = libceed::Ceed::default_init();
/// const Q: usize = 6;
/// let bu = ceed.basis_tensor_H1_Lagrange(1, 1, Q, Q, QuadMode::GaussLobatto)?;
/// let bx = ceed.basis_tensor_H1_Lagrange(1, 1, 2, Q, QuadMode::Gauss)?;
///
/// let x_corners = ceed.vector_from_slice(&[-1., 1.])?;
/// let mut x_qpts = ceed.vector(Q)?;
/// let mut x_nodes = ceed.vector(Q)?;
/// bx.apply(
/// 1,
/// TransposeMode::NoTranspose,
/// EvalMode::Interp,
/// &x_corners,
/// &mut x_nodes,
/// )?;
/// bu.apply(
/// 1,
/// TransposeMode::NoTranspose,
/// EvalMode::Interp,
/// &x_nodes,
/// &mut x_qpts,
/// )?;
///
/// // Create function x^3 + 1 on Gauss Lobatto points
/// let mut u_arr = [0.; Q];
/// u_arr
/// .iter_mut()
/// .zip(x_nodes.view()?.iter())
/// .for_each(|(u, x)| *u = x * x * x + 1.);
/// let u = ceed.vector_from_slice(&u_arr)?;
///
/// // Map function to Gauss points
/// let mut v = ceed.vector(Q)?;
/// v.set_value(0.);
/// bu.apply(1, TransposeMode::NoTranspose, EvalMode::Interp, &u, &mut v)?;
///
/// // Verify results
/// v.view()?
/// .iter()
/// .zip(x_qpts.view()?.iter())
/// .for_each(|(v, x)| {
/// let true_value = x * x * x + 1.;
/// assert!(
/// (*v - true_value).abs() < 10.0 * libceed::EPSILON,
/// "Incorrect basis application"
/// );
/// });
/// # Ok(())
/// # }
/// ```
pub fn apply(
&self,
nelem: usize,
tmode: TransposeMode,
emode: EvalMode,
u: &Vector,
v: &mut Vector,
) -> crate::Result<i32> {
let (nelem, tmode, emode) = (
i32::try_from(nelem).unwrap(),
tmode as bind_ceed::CeedTransposeMode,
emode as bind_ceed::CeedEvalMode,
);
let ierr =
unsafe { bind_ceed::CeedBasisApply(self.ptr, nelem, tmode, emode, u.ptr, v.ptr) };
self.check_error(ierr)
}
/// Returns the dimension for given Basis
///
/// ```
/// # use libceed::prelude::*;
/// # fn main() -> libceed::Result<()> {
/// # let ceed = libceed::Ceed::default_init();
/// let dim = 2;
/// let b = ceed.basis_tensor_H1_Lagrange(dim, 1, 3, 4, QuadMode::Gauss)?;
///
/// let d = b.dimension();
/// assert_eq!(d, dim, "Incorrect dimension");
/// # Ok(())
/// # }
/// ```
pub fn dimension(&self) -> usize {
let mut dim = 0;
unsafe { bind_ceed::CeedBasisGetDimension(self.ptr, &mut dim) };
usize::try_from(dim).unwrap()
}
/// Returns number of components for given Basis
///
/// ```
/// # use libceed::prelude::*;
/// # fn main() -> libceed::Result<()> {
/// # let ceed = libceed::Ceed::default_init();
/// let ncomp = 2;
/// let b = ceed.basis_tensor_H1_Lagrange(1, ncomp, 3, 4, QuadMode::Gauss)?;
///
/// let n = b.num_components();
/// assert_eq!(n, ncomp, "Incorrect number of components");
/// # Ok(())
/// # }
/// ```
pub fn num_components(&self) -> usize {
let mut ncomp = 0;
unsafe { bind_ceed::CeedBasisGetNumComponents(self.ptr, &mut ncomp) };
usize::try_from(ncomp).unwrap()
}
/// Returns total number of nodes (in dim dimensions) of a Basis
///
/// ```
/// # use libceed::prelude::*;
/// # fn main() -> libceed::Result<()> {
/// # let ceed = libceed::Ceed::default_init();
/// let p = 3;
/// let b = ceed.basis_tensor_H1_Lagrange(2, 1, p, 4, QuadMode::Gauss)?;
///
/// let nnodes = b.num_nodes();
/// assert_eq!(nnodes, p * p, "Incorrect number of nodes");
/// # Ok(())
/// # }
/// ```
pub fn num_nodes(&self) -> usize {
let mut nnodes = 0;
unsafe { bind_ceed::CeedBasisGetNumNodes(self.ptr, &mut nnodes) };
usize::try_from(nnodes).unwrap()
}
/// Returns total number of quadrature points (in dim dimensions) of a
/// Basis
///
/// ```
/// # use libceed::prelude::*;
/// # fn main() -> libceed::Result<()> {
/// # let ceed = libceed::Ceed::default_init();
/// let q = 4;
/// let b = ceed.basis_tensor_H1_Lagrange(2, 1, 3, q, QuadMode::Gauss)?;
///
/// let nqpts = b.num_quadrature_points();
/// assert_eq!(nqpts, q * q, "Incorrect number of quadrature points");
/// # Ok(())
/// # }
/// ```
pub fn num_quadrature_points(&self) -> usize {
let mut Q = 0;
unsafe {
bind_ceed::CeedBasisGetNumQuadraturePoints(self.ptr, &mut Q);
}
usize::try_from(Q).unwrap()
}
/// Create projection from self to specified Basis.
///
/// Both bases must have the same quadrature space. The input bases need not
/// be nested as function spaces; this interface solves a least squares
/// problem to find a representation in the `to` basis that agrees at
/// quadrature points with the origin basis. Since the bases need not be
/// Lagrange, the resulting projection "basis" will have empty quadrature
/// points and weights.
///
/// ```
/// # use libceed::prelude::*;
/// # fn main() -> libceed::Result<()> {
/// # let ceed = libceed::Ceed::default_init();
/// let coarse = ceed.basis_tensor_H1_Lagrange(1, 1, 2, 3, QuadMode::Gauss)?;
/// let fine = ceed.basis_tensor_H1_Lagrange(1, 1, 3, 3, QuadMode::Gauss)?;
/// let proj = coarse.create_projection(&fine)?;
/// let u = ceed.vector_from_slice(&[1., 2.])?;
/// let mut v = ceed.vector(3)?;
/// proj.apply(1, TransposeMode::NoTranspose, EvalMode::Interp, &u, &mut v)?;
/// let expected = [1., 1.5, 2.];
/// for (a, b) in v.view()?.iter().zip(expected) {
/// assert!(
/// (a - b).abs() < 10.0 * libceed::EPSILON,
/// "Incorrect projection of linear Lagrange to quadratic Lagrange"
/// );
/// }
/// # Ok(())
/// # }
/// ```
pub fn create_projection(&self, to: &Self) -> crate::Result<Self> {
let mut ptr = std::ptr::null_mut();
let ierr = unsafe { bind_ceed::CeedBasisCreateProjection(self.ptr, to.ptr, &mut ptr) };
self.check_error(ierr)?;
Ok(Self {
ptr,
_lifeline: PhantomData,
})
}
}
// -----------------------------------------------------------------------------