15aed82e4SJeremy L Thompson // Copyright (c) 2017-2024, Lawrence Livermore National Security, LLC and other CEED contributors. 23d8e8822SJeremy L Thompson // All Rights Reserved. See the top-level LICENSE and NOTICE files for details. 39df49d7eSJed Brown // 43d8e8822SJeremy L Thompson // SPDX-License-Identifier: BSD-2-Clause 59df49d7eSJed Brown // 63d8e8822SJeremy L Thompson // This file is part of CEED: http://github.com/ceed 79df49d7eSJed Brown 89df49d7eSJed Brown //! A Ceed Operator defines the finite/spectral element operator associated to a 99df49d7eSJed Brown //! Ceed QFunction. A Ceed Operator connects Ceed ElemRestrictions, 109df49d7eSJed Brown //! Ceed Bases, and Ceed QFunctions. 119df49d7eSJed Brown 129df49d7eSJed Brown use crate::prelude::*; 139df49d7eSJed Brown 149df49d7eSJed Brown // ----------------------------------------------------------------------------- 157ed177dbSJed Brown // Operator Field context wrapper 1608778c6fSJeremy L Thompson // ----------------------------------------------------------------------------- 1708778c6fSJeremy L Thompson #[derive(Debug)] 1808778c6fSJeremy L Thompson pub struct OperatorField<'a> { 1908778c6fSJeremy L Thompson pub(crate) ptr: bind_ceed::CeedOperatorField, 20e03fef56SJeremy L Thompson pub(crate) vector: crate::Vector<'a>, 21e03fef56SJeremy L Thompson pub(crate) elem_restriction: crate::ElemRestriction<'a>, 22e03fef56SJeremy L Thompson pub(crate) basis: crate::Basis<'a>, 2308778c6fSJeremy L Thompson _lifeline: PhantomData<&'a ()>, 2408778c6fSJeremy L Thompson } 2508778c6fSJeremy L Thompson 2608778c6fSJeremy L Thompson // ----------------------------------------------------------------------------- 2708778c6fSJeremy L Thompson // Implementations 2808778c6fSJeremy L Thompson // ----------------------------------------------------------------------------- 2908778c6fSJeremy L Thompson impl<'a> OperatorField<'a> { 30e03fef56SJeremy L Thompson pub(crate) fn from_raw( 31e03fef56SJeremy L Thompson ptr: bind_ceed::CeedOperatorField, 32e03fef56SJeremy L Thompson ceed: crate::Ceed, 33e03fef56SJeremy L Thompson ) -> crate::Result<Self> { 34e03fef56SJeremy L Thompson let vector = { 35e03fef56SJeremy L Thompson let mut vector_ptr = std::ptr::null_mut(); 36e03fef56SJeremy L Thompson let ierr = unsafe { bind_ceed::CeedOperatorFieldGetVector(ptr, &mut vector_ptr) }; 37e03fef56SJeremy L Thompson ceed.check_error(ierr)?; 38e03fef56SJeremy L Thompson crate::Vector::from_raw(vector_ptr)? 39e03fef56SJeremy L Thompson }; 40e03fef56SJeremy L Thompson let elem_restriction = { 41e03fef56SJeremy L Thompson let mut elem_restriction_ptr = std::ptr::null_mut(); 42e03fef56SJeremy L Thompson let ierr = unsafe { 43e03fef56SJeremy L Thompson bind_ceed::CeedOperatorFieldGetElemRestriction(ptr, &mut elem_restriction_ptr) 44e03fef56SJeremy L Thompson }; 45e03fef56SJeremy L Thompson ceed.check_error(ierr)?; 46e03fef56SJeremy L Thompson crate::ElemRestriction::from_raw(elem_restriction_ptr)? 47e03fef56SJeremy L Thompson }; 48e03fef56SJeremy L Thompson let basis = { 49e03fef56SJeremy L Thompson let mut basis_ptr = std::ptr::null_mut(); 50e03fef56SJeremy L Thompson let ierr = unsafe { bind_ceed::CeedOperatorFieldGetBasis(ptr, &mut basis_ptr) }; 51e03fef56SJeremy L Thompson ceed.check_error(ierr)?; 52e03fef56SJeremy L Thompson crate::Basis::from_raw(basis_ptr)? 53e03fef56SJeremy L Thompson }; 54e03fef56SJeremy L Thompson Ok(Self { 55e03fef56SJeremy L Thompson ptr, 56e03fef56SJeremy L Thompson vector, 57e03fef56SJeremy L Thompson elem_restriction, 58e03fef56SJeremy L Thompson basis, 59e03fef56SJeremy L Thompson _lifeline: PhantomData, 60e03fef56SJeremy L Thompson }) 61e03fef56SJeremy L Thompson } 62e03fef56SJeremy L Thompson 6308778c6fSJeremy L Thompson /// Get the name of an OperatorField 6408778c6fSJeremy L Thompson /// 6508778c6fSJeremy L Thompson /// ``` 6608778c6fSJeremy L Thompson /// # use libceed::prelude::*; 674d27c890SJeremy L Thompson /// # fn main() -> libceed::Result<()> { 6808778c6fSJeremy L Thompson /// # let ceed = libceed::Ceed::default_init(); 6908778c6fSJeremy L Thompson /// let qf = ceed.q_function_interior_by_name("Mass1DBuild")?; 7008778c6fSJeremy L Thompson /// 7108778c6fSJeremy L Thompson /// // Operator field arguments 7208778c6fSJeremy L Thompson /// let ne = 3; 7308778c6fSJeremy L Thompson /// let q = 4 as usize; 7408778c6fSJeremy L Thompson /// let mut ind: Vec<i32> = vec![0; 2 * ne]; 7508778c6fSJeremy L Thompson /// for i in 0..ne { 7608778c6fSJeremy L Thompson /// ind[2 * i + 0] = i as i32; 7708778c6fSJeremy L Thompson /// ind[2 * i + 1] = (i + 1) as i32; 7808778c6fSJeremy L Thompson /// } 7908778c6fSJeremy L Thompson /// let r = ceed.elem_restriction(ne, 2, 1, 1, ne + 1, MemType::Host, &ind)?; 8008778c6fSJeremy L Thompson /// let strides: [i32; 3] = [1, q as i32, q as i32]; 81d9267b76SJeremy L Thompson /// let rq = ceed.strided_elem_restriction(ne, q, 1, q * ne, strides)?; 8208778c6fSJeremy L Thompson /// 8308778c6fSJeremy L Thompson /// let b = ceed.basis_tensor_H1_Lagrange(1, 1, 2, q, QuadMode::Gauss)?; 8408778c6fSJeremy L Thompson /// 8508778c6fSJeremy L Thompson /// // Operator fields 8608778c6fSJeremy L Thompson /// let op = ceed 8708778c6fSJeremy L Thompson /// .operator(&qf, QFunctionOpt::None, QFunctionOpt::None)? 8808778c6fSJeremy L Thompson /// .field("dx", &r, &b, VectorOpt::Active)? 8908778c6fSJeremy L Thompson /// .field("weights", ElemRestrictionOpt::None, &b, VectorOpt::None)? 90356036faSJeremy L Thompson /// .field("qdata", &rq, BasisOpt::None, VectorOpt::Active)?; 9108778c6fSJeremy L Thompson /// 9208778c6fSJeremy L Thompson /// let inputs = op.inputs()?; 9308778c6fSJeremy L Thompson /// 9408778c6fSJeremy L Thompson /// assert_eq!(inputs[0].name(), "dx", "Incorrect input name"); 9508778c6fSJeremy L Thompson /// assert_eq!(inputs[1].name(), "weights", "Incorrect input name"); 9608778c6fSJeremy L Thompson /// # Ok(()) 9708778c6fSJeremy L Thompson /// # } 9808778c6fSJeremy L Thompson /// ``` 9908778c6fSJeremy L Thompson pub fn name(&self) -> &str { 10008778c6fSJeremy L Thompson let mut name_ptr: *mut std::os::raw::c_char = std::ptr::null_mut(); 10108778c6fSJeremy L Thompson unsafe { 1026f8994e9SJeremy L Thompson bind_ceed::CeedOperatorFieldGetName( 1036f8994e9SJeremy L Thompson self.ptr, 1046f8994e9SJeremy L Thompson &mut name_ptr as *const _ as *mut *const _, 1056f8994e9SJeremy L Thompson ); 10608778c6fSJeremy L Thompson } 10708778c6fSJeremy L Thompson unsafe { CStr::from_ptr(name_ptr) }.to_str().unwrap() 10808778c6fSJeremy L Thompson } 10908778c6fSJeremy L Thompson 11008778c6fSJeremy L Thompson /// Get the ElemRestriction of an OperatorField 11108778c6fSJeremy L Thompson /// 11208778c6fSJeremy L Thompson /// ``` 11308778c6fSJeremy L Thompson /// # use libceed::prelude::*; 1144d27c890SJeremy L Thompson /// # fn main() -> libceed::Result<()> { 11508778c6fSJeremy L Thompson /// # let ceed = libceed::Ceed::default_init(); 11608778c6fSJeremy L Thompson /// let qf = ceed.q_function_interior_by_name("Mass1DBuild")?; 11708778c6fSJeremy L Thompson /// 11808778c6fSJeremy L Thompson /// // Operator field arguments 11908778c6fSJeremy L Thompson /// let ne = 3; 12008778c6fSJeremy L Thompson /// let q = 4 as usize; 12108778c6fSJeremy L Thompson /// let mut ind: Vec<i32> = vec![0; 2 * ne]; 12208778c6fSJeremy L Thompson /// for i in 0..ne { 12308778c6fSJeremy L Thompson /// ind[2 * i + 0] = i as i32; 12408778c6fSJeremy L Thompson /// ind[2 * i + 1] = (i + 1) as i32; 12508778c6fSJeremy L Thompson /// } 12608778c6fSJeremy L Thompson /// let r = ceed.elem_restriction(ne, 2, 1, 1, ne + 1, MemType::Host, &ind)?; 12708778c6fSJeremy L Thompson /// let strides: [i32; 3] = [1, q as i32, q as i32]; 128d9267b76SJeremy L Thompson /// let rq = ceed.strided_elem_restriction(ne, q, 1, q * ne, strides)?; 12908778c6fSJeremy L Thompson /// 13008778c6fSJeremy L Thompson /// let b = ceed.basis_tensor_H1_Lagrange(1, 1, 2, q, QuadMode::Gauss)?; 13108778c6fSJeremy L Thompson /// 13208778c6fSJeremy L Thompson /// // Operator fields 13308778c6fSJeremy L Thompson /// let op = ceed 13408778c6fSJeremy L Thompson /// .operator(&qf, QFunctionOpt::None, QFunctionOpt::None)? 13508778c6fSJeremy L Thompson /// .field("dx", &r, &b, VectorOpt::Active)? 13608778c6fSJeremy L Thompson /// .field("weights", ElemRestrictionOpt::None, &b, VectorOpt::None)? 137356036faSJeremy L Thompson /// .field("qdata", &rq, BasisOpt::None, VectorOpt::Active)?; 13808778c6fSJeremy L Thompson /// 13908778c6fSJeremy L Thompson /// let inputs = op.inputs()?; 14008778c6fSJeremy L Thompson /// 14108778c6fSJeremy L Thompson /// assert!( 14208778c6fSJeremy L Thompson /// inputs[0].elem_restriction().is_some(), 14308778c6fSJeremy L Thompson /// "Incorrect field ElemRestriction" 14408778c6fSJeremy L Thompson /// ); 145567c3c29SJeremy L Thompson /// if let ElemRestrictionOpt::Some(r) = inputs[0].elem_restriction() { 146567c3c29SJeremy L Thompson /// assert_eq!( 147567c3c29SJeremy L Thompson /// r.num_elements(), 148567c3c29SJeremy L Thompson /// ne, 149567c3c29SJeremy L Thompson /// "Incorrect field ElemRestriction number of elements" 150567c3c29SJeremy L Thompson /// ); 151567c3c29SJeremy L Thompson /// } 152567c3c29SJeremy L Thompson /// 15308778c6fSJeremy L Thompson /// assert!( 15408778c6fSJeremy L Thompson /// inputs[1].elem_restriction().is_none(), 15508778c6fSJeremy L Thompson /// "Incorrect field ElemRestriction" 15608778c6fSJeremy L Thompson /// ); 157e03fef56SJeremy L Thompson /// 158e03fef56SJeremy L Thompson /// let outputs = op.outputs()?; 159e03fef56SJeremy L Thompson /// 160e03fef56SJeremy L Thompson /// assert!( 161e03fef56SJeremy L Thompson /// outputs[0].elem_restriction().is_some(), 162e03fef56SJeremy L Thompson /// "Incorrect field ElemRestriction" 163e03fef56SJeremy L Thompson /// ); 164567c3c29SJeremy L Thompson /// if let ElemRestrictionOpt::Some(r) = outputs[0].elem_restriction() { 165567c3c29SJeremy L Thompson /// assert_eq!( 166567c3c29SJeremy L Thompson /// r.num_elements(), 167567c3c29SJeremy L Thompson /// ne, 168567c3c29SJeremy L Thompson /// "Incorrect field ElemRestriction number of elements" 169567c3c29SJeremy L Thompson /// ); 170567c3c29SJeremy L Thompson /// } 17108778c6fSJeremy L Thompson /// # Ok(()) 17208778c6fSJeremy L Thompson /// # } 17308778c6fSJeremy L Thompson /// ``` 17408778c6fSJeremy L Thompson pub fn elem_restriction(&self) -> ElemRestrictionOpt { 175e03fef56SJeremy L Thompson if self.elem_restriction.ptr == unsafe { bind_ceed::CEED_ELEMRESTRICTION_NONE } { 17608778c6fSJeremy L Thompson ElemRestrictionOpt::None 17708778c6fSJeremy L Thompson } else { 178e03fef56SJeremy L Thompson ElemRestrictionOpt::Some(&self.elem_restriction) 17908778c6fSJeremy L Thompson } 18008778c6fSJeremy L Thompson } 18108778c6fSJeremy L Thompson 18208778c6fSJeremy L Thompson /// Get the Basis of an OperatorField 18308778c6fSJeremy L Thompson /// 18408778c6fSJeremy L Thompson /// ``` 18508778c6fSJeremy L Thompson /// # use libceed::prelude::*; 1864d27c890SJeremy L Thompson /// # fn main() -> libceed::Result<()> { 18708778c6fSJeremy L Thompson /// # let ceed = libceed::Ceed::default_init(); 18808778c6fSJeremy L Thompson /// let qf = ceed.q_function_interior_by_name("Mass1DBuild")?; 18908778c6fSJeremy L Thompson /// 19008778c6fSJeremy L Thompson /// // Operator field arguments 19108778c6fSJeremy L Thompson /// let ne = 3; 19208778c6fSJeremy L Thompson /// let q = 4 as usize; 19308778c6fSJeremy L Thompson /// let mut ind: Vec<i32> = vec![0; 2 * ne]; 19408778c6fSJeremy L Thompson /// for i in 0..ne { 19508778c6fSJeremy L Thompson /// ind[2 * i + 0] = i as i32; 19608778c6fSJeremy L Thompson /// ind[2 * i + 1] = (i + 1) as i32; 19708778c6fSJeremy L Thompson /// } 19808778c6fSJeremy L Thompson /// let r = ceed.elem_restriction(ne, 2, 1, 1, ne + 1, MemType::Host, &ind)?; 19908778c6fSJeremy L Thompson /// let strides: [i32; 3] = [1, q as i32, q as i32]; 200d9267b76SJeremy L Thompson /// let rq = ceed.strided_elem_restriction(ne, q, 1, q * ne, strides)?; 20108778c6fSJeremy L Thompson /// 20208778c6fSJeremy L Thompson /// let b = ceed.basis_tensor_H1_Lagrange(1, 1, 2, q, QuadMode::Gauss)?; 20308778c6fSJeremy L Thompson /// 20408778c6fSJeremy L Thompson /// // Operator fields 20508778c6fSJeremy L Thompson /// let op = ceed 20608778c6fSJeremy L Thompson /// .operator(&qf, QFunctionOpt::None, QFunctionOpt::None)? 20708778c6fSJeremy L Thompson /// .field("dx", &r, &b, VectorOpt::Active)? 20808778c6fSJeremy L Thompson /// .field("weights", ElemRestrictionOpt::None, &b, VectorOpt::None)? 209356036faSJeremy L Thompson /// .field("qdata", &rq, BasisOpt::None, VectorOpt::Active)?; 21008778c6fSJeremy L Thompson /// 21108778c6fSJeremy L Thompson /// let inputs = op.inputs()?; 21208778c6fSJeremy L Thompson /// 21308778c6fSJeremy L Thompson /// assert!(inputs[0].basis().is_some(), "Incorrect field Basis"); 214567c3c29SJeremy L Thompson /// if let BasisOpt::Some(b) = inputs[0].basis() { 215567c3c29SJeremy L Thompson /// assert_eq!( 216567c3c29SJeremy L Thompson /// b.num_quadrature_points(), 217567c3c29SJeremy L Thompson /// q, 218567c3c29SJeremy L Thompson /// "Incorrect field Basis number of quadrature points" 219567c3c29SJeremy L Thompson /// ); 220567c3c29SJeremy L Thompson /// } 22108778c6fSJeremy L Thompson /// assert!(inputs[1].basis().is_some(), "Incorrect field Basis"); 222567c3c29SJeremy L Thompson /// if let BasisOpt::Some(b) = inputs[1].basis() { 223567c3c29SJeremy L Thompson /// assert_eq!( 224567c3c29SJeremy L Thompson /// b.num_quadrature_points(), 225567c3c29SJeremy L Thompson /// q, 226567c3c29SJeremy L Thompson /// "Incorrect field Basis number of quadrature points" 227567c3c29SJeremy L Thompson /// ); 228567c3c29SJeremy L Thompson /// } 22908778c6fSJeremy L Thompson /// 23008778c6fSJeremy L Thompson /// let outputs = op.outputs()?; 23108778c6fSJeremy L Thompson /// 232356036faSJeremy L Thompson /// assert!(outputs[0].basis().is_none(), "Incorrect field Basis"); 23308778c6fSJeremy L Thompson /// # Ok(()) 23408778c6fSJeremy L Thompson /// # } 23508778c6fSJeremy L Thompson /// ``` 23608778c6fSJeremy L Thompson pub fn basis(&self) -> BasisOpt { 237e03fef56SJeremy L Thompson if self.basis.ptr == unsafe { bind_ceed::CEED_BASIS_NONE } { 238356036faSJeremy L Thompson BasisOpt::None 23908778c6fSJeremy L Thompson } else { 240e03fef56SJeremy L Thompson BasisOpt::Some(&self.basis) 24108778c6fSJeremy L Thompson } 24208778c6fSJeremy L Thompson } 24308778c6fSJeremy L Thompson 24408778c6fSJeremy L Thompson /// Get the Vector of an OperatorField 24508778c6fSJeremy L Thompson /// 24608778c6fSJeremy L Thompson /// ``` 24708778c6fSJeremy L Thompson /// # use libceed::prelude::*; 2484d27c890SJeremy L Thompson /// # fn main() -> libceed::Result<()> { 24908778c6fSJeremy L Thompson /// # let ceed = libceed::Ceed::default_init(); 25008778c6fSJeremy L Thompson /// let qf = ceed.q_function_interior_by_name("Mass1DBuild")?; 25108778c6fSJeremy L Thompson /// 25208778c6fSJeremy L Thompson /// // Operator field arguments 25308778c6fSJeremy L Thompson /// let ne = 3; 25408778c6fSJeremy L Thompson /// let q = 4 as usize; 25508778c6fSJeremy L Thompson /// let mut ind: Vec<i32> = vec![0; 2 * ne]; 25608778c6fSJeremy L Thompson /// for i in 0..ne { 25708778c6fSJeremy L Thompson /// ind[2 * i + 0] = i as i32; 25808778c6fSJeremy L Thompson /// ind[2 * i + 1] = (i + 1) as i32; 25908778c6fSJeremy L Thompson /// } 26008778c6fSJeremy L Thompson /// let r = ceed.elem_restriction(ne, 2, 1, 1, ne + 1, MemType::Host, &ind)?; 26108778c6fSJeremy L Thompson /// let strides: [i32; 3] = [1, q as i32, q as i32]; 262d9267b76SJeremy L Thompson /// let rq = ceed.strided_elem_restriction(ne, q, 1, q * ne, strides)?; 26308778c6fSJeremy L Thompson /// 26408778c6fSJeremy L Thompson /// let b = ceed.basis_tensor_H1_Lagrange(1, 1, 2, q, QuadMode::Gauss)?; 26508778c6fSJeremy L Thompson /// 26608778c6fSJeremy L Thompson /// // Operator fields 26708778c6fSJeremy L Thompson /// let op = ceed 26808778c6fSJeremy L Thompson /// .operator(&qf, QFunctionOpt::None, QFunctionOpt::None)? 26908778c6fSJeremy L Thompson /// .field("dx", &r, &b, VectorOpt::Active)? 27008778c6fSJeremy L Thompson /// .field("weights", ElemRestrictionOpt::None, &b, VectorOpt::None)? 271356036faSJeremy L Thompson /// .field("qdata", &rq, BasisOpt::None, VectorOpt::Active)?; 27208778c6fSJeremy L Thompson /// 27308778c6fSJeremy L Thompson /// let inputs = op.inputs()?; 27408778c6fSJeremy L Thompson /// 27508778c6fSJeremy L Thompson /// assert!(inputs[0].vector().is_active(), "Incorrect field Vector"); 27608778c6fSJeremy L Thompson /// assert!(inputs[1].vector().is_none(), "Incorrect field Vector"); 277e03fef56SJeremy L Thompson /// 278e03fef56SJeremy L Thompson /// let outputs = op.outputs()?; 279e03fef56SJeremy L Thompson /// 280e03fef56SJeremy L Thompson /// assert!(outputs[0].vector().is_active(), "Incorrect field Vector"); 28108778c6fSJeremy L Thompson /// # Ok(()) 28208778c6fSJeremy L Thompson /// # } 28308778c6fSJeremy L Thompson /// ``` 28408778c6fSJeremy L Thompson pub fn vector(&self) -> VectorOpt { 285e03fef56SJeremy L Thompson if self.vector.ptr == unsafe { bind_ceed::CEED_VECTOR_ACTIVE } { 28608778c6fSJeremy L Thompson VectorOpt::Active 287e03fef56SJeremy L Thompson } else if self.vector.ptr == unsafe { bind_ceed::CEED_VECTOR_NONE } { 28808778c6fSJeremy L Thompson VectorOpt::None 28908778c6fSJeremy L Thompson } else { 290e03fef56SJeremy L Thompson VectorOpt::Some(&self.vector) 29108778c6fSJeremy L Thompson } 29208778c6fSJeremy L Thompson } 29308778c6fSJeremy L Thompson } 29408778c6fSJeremy L Thompson 29508778c6fSJeremy L Thompson // ----------------------------------------------------------------------------- 2967ed177dbSJed Brown // Operator context wrapper 2979df49d7eSJed Brown // ----------------------------------------------------------------------------- 298c68be7a2SJeremy L Thompson #[derive(Debug)] 2999df49d7eSJed Brown pub(crate) struct OperatorCore<'a> { 3009df49d7eSJed Brown ptr: bind_ceed::CeedOperator, 3011142270cSJeremy L Thompson _lifeline: PhantomData<&'a ()>, 3029df49d7eSJed Brown } 3039df49d7eSJed Brown 304c68be7a2SJeremy L Thompson #[derive(Debug)] 3059df49d7eSJed Brown pub struct Operator<'a> { 3069df49d7eSJed Brown op_core: OperatorCore<'a>, 3079df49d7eSJed Brown } 3089df49d7eSJed Brown 309c68be7a2SJeremy L Thompson #[derive(Debug)] 3109df49d7eSJed Brown pub struct CompositeOperator<'a> { 3119df49d7eSJed Brown op_core: OperatorCore<'a>, 3129df49d7eSJed Brown } 3139df49d7eSJed Brown 3149df49d7eSJed Brown // ----------------------------------------------------------------------------- 3159df49d7eSJed Brown // Destructor 3169df49d7eSJed Brown // ----------------------------------------------------------------------------- 3179df49d7eSJed Brown impl<'a> Drop for OperatorCore<'a> { 3189df49d7eSJed Brown fn drop(&mut self) { 3199df49d7eSJed Brown unsafe { 3209df49d7eSJed Brown bind_ceed::CeedOperatorDestroy(&mut self.ptr); 3219df49d7eSJed Brown } 3229df49d7eSJed Brown } 3239df49d7eSJed Brown } 3249df49d7eSJed Brown 3259df49d7eSJed Brown // ----------------------------------------------------------------------------- 3269df49d7eSJed Brown // Display 3279df49d7eSJed Brown // ----------------------------------------------------------------------------- 3289df49d7eSJed Brown impl<'a> fmt::Display for OperatorCore<'a> { 3299df49d7eSJed Brown fn fmt(&self, f: &mut fmt::Formatter) -> fmt::Result { 3309df49d7eSJed Brown let mut ptr = std::ptr::null_mut(); 3319df49d7eSJed Brown let mut sizeloc = crate::MAX_BUFFER_LENGTH; 3329df49d7eSJed Brown let cstring = unsafe { 3339df49d7eSJed Brown let file = bind_ceed::open_memstream(&mut ptr, &mut sizeloc); 3349df49d7eSJed Brown bind_ceed::CeedOperatorView(self.ptr, file); 3359df49d7eSJed Brown bind_ceed::fclose(file); 3369df49d7eSJed Brown CString::from_raw(ptr) 3379df49d7eSJed Brown }; 3389df49d7eSJed Brown cstring.to_string_lossy().fmt(f) 3399df49d7eSJed Brown } 3409df49d7eSJed Brown } 3419df49d7eSJed Brown 3429df49d7eSJed Brown /// View an Operator 3439df49d7eSJed Brown /// 3449df49d7eSJed Brown /// ``` 3459df49d7eSJed Brown /// # use libceed::prelude::*; 3464d27c890SJeremy L Thompson /// # fn main() -> libceed::Result<()> { 3479df49d7eSJed Brown /// # let ceed = libceed::Ceed::default_init(); 348c68be7a2SJeremy L Thompson /// let qf = ceed.q_function_interior_by_name("Mass1DBuild")?; 3499df49d7eSJed Brown /// 3509df49d7eSJed Brown /// // Operator field arguments 3519df49d7eSJed Brown /// let ne = 3; 3529df49d7eSJed Brown /// let q = 4 as usize; 3539df49d7eSJed Brown /// let mut ind: Vec<i32> = vec![0; 2 * ne]; 3549df49d7eSJed Brown /// for i in 0..ne { 3559df49d7eSJed Brown /// ind[2 * i + 0] = i as i32; 3569df49d7eSJed Brown /// ind[2 * i + 1] = (i + 1) as i32; 3579df49d7eSJed Brown /// } 358c68be7a2SJeremy L Thompson /// let r = ceed.elem_restriction(ne, 2, 1, 1, ne + 1, MemType::Host, &ind)?; 3599df49d7eSJed Brown /// let strides: [i32; 3] = [1, q as i32, q as i32]; 360d9267b76SJeremy L Thompson /// let rq = ceed.strided_elem_restriction(ne, q, 1, q * ne, strides)?; 3619df49d7eSJed Brown /// 362c68be7a2SJeremy L Thompson /// let b = ceed.basis_tensor_H1_Lagrange(1, 1, 2, q, QuadMode::Gauss)?; 3639df49d7eSJed Brown /// 3649df49d7eSJed Brown /// // Operator fields 3659df49d7eSJed Brown /// let op = ceed 366c68be7a2SJeremy L Thompson /// .operator(&qf, QFunctionOpt::None, QFunctionOpt::None)? 367ea6b5821SJeremy L Thompson /// .name("mass")? 368c68be7a2SJeremy L Thompson /// .field("dx", &r, &b, VectorOpt::Active)? 369c68be7a2SJeremy L Thompson /// .field("weights", ElemRestrictionOpt::None, &b, VectorOpt::None)? 370356036faSJeremy L Thompson /// .field("qdata", &rq, BasisOpt::None, VectorOpt::Active)?; 3719df49d7eSJed Brown /// 3729df49d7eSJed Brown /// println!("{}", op); 373c68be7a2SJeremy L Thompson /// # Ok(()) 374c68be7a2SJeremy L Thompson /// # } 3759df49d7eSJed Brown /// ``` 3769df49d7eSJed Brown impl<'a> fmt::Display for Operator<'a> { 3779df49d7eSJed Brown fn fmt(&self, f: &mut fmt::Formatter) -> fmt::Result { 3789df49d7eSJed Brown self.op_core.fmt(f) 3799df49d7eSJed Brown } 3809df49d7eSJed Brown } 3819df49d7eSJed Brown 3829df49d7eSJed Brown /// View a composite Operator 3839df49d7eSJed Brown /// 3849df49d7eSJed Brown /// ``` 3859df49d7eSJed Brown /// # use libceed::prelude::*; 3864d27c890SJeremy L Thompson /// # fn main() -> libceed::Result<()> { 3879df49d7eSJed Brown /// # let ceed = libceed::Ceed::default_init(); 3889df49d7eSJed Brown /// 3899df49d7eSJed Brown /// // Sub operator field arguments 3909df49d7eSJed Brown /// let ne = 3; 3919df49d7eSJed Brown /// let q = 4 as usize; 3929df49d7eSJed Brown /// let mut ind: Vec<i32> = vec![0; 2 * ne]; 3939df49d7eSJed Brown /// for i in 0..ne { 3949df49d7eSJed Brown /// ind[2 * i + 0] = i as i32; 3959df49d7eSJed Brown /// ind[2 * i + 1] = (i + 1) as i32; 3969df49d7eSJed Brown /// } 397c68be7a2SJeremy L Thompson /// let r = ceed.elem_restriction(ne, 2, 1, 1, ne + 1, MemType::Host, &ind)?; 3989df49d7eSJed Brown /// let strides: [i32; 3] = [1, q as i32, q as i32]; 399d9267b76SJeremy L Thompson /// let rq = ceed.strided_elem_restriction(ne, q, 1, q * ne, strides)?; 4009df49d7eSJed Brown /// 401c68be7a2SJeremy L Thompson /// let b = ceed.basis_tensor_H1_Lagrange(1, 1, 2, q, QuadMode::Gauss)?; 4029df49d7eSJed Brown /// 403c68be7a2SJeremy L Thompson /// let qdata_mass = ceed.vector(q * ne)?; 404c68be7a2SJeremy L Thompson /// let qdata_diff = ceed.vector(q * ne)?; 4059df49d7eSJed Brown /// 406c68be7a2SJeremy L Thompson /// let qf_mass = ceed.q_function_interior_by_name("MassApply")?; 4079df49d7eSJed Brown /// let op_mass = ceed 408c68be7a2SJeremy L Thompson /// .operator(&qf_mass, QFunctionOpt::None, QFunctionOpt::None)? 409ea6b5821SJeremy L Thompson /// .name("Mass term")? 410c68be7a2SJeremy L Thompson /// .field("u", &r, &b, VectorOpt::Active)? 411356036faSJeremy L Thompson /// .field("qdata", &rq, BasisOpt::None, &qdata_mass)? 412c68be7a2SJeremy L Thompson /// .field("v", &r, &b, VectorOpt::Active)?; 4139df49d7eSJed Brown /// 414c68be7a2SJeremy L Thompson /// let qf_diff = ceed.q_function_interior_by_name("Poisson1DApply")?; 4159df49d7eSJed Brown /// let op_diff = ceed 416c68be7a2SJeremy L Thompson /// .operator(&qf_diff, QFunctionOpt::None, QFunctionOpt::None)? 417ea6b5821SJeremy L Thompson /// .name("Poisson term")? 418c68be7a2SJeremy L Thompson /// .field("du", &r, &b, VectorOpt::Active)? 419356036faSJeremy L Thompson /// .field("qdata", &rq, BasisOpt::None, &qdata_diff)? 420c68be7a2SJeremy L Thompson /// .field("dv", &r, &b, VectorOpt::Active)?; 4219df49d7eSJed Brown /// 4229df49d7eSJed Brown /// let op = ceed 423c68be7a2SJeremy L Thompson /// .composite_operator()? 424ea6b5821SJeremy L Thompson /// .name("Screened Poisson")? 425c68be7a2SJeremy L Thompson /// .sub_operator(&op_mass)? 426c68be7a2SJeremy L Thompson /// .sub_operator(&op_diff)?; 4279df49d7eSJed Brown /// 4289df49d7eSJed Brown /// println!("{}", op); 429c68be7a2SJeremy L Thompson /// # Ok(()) 430c68be7a2SJeremy L Thompson /// # } 4319df49d7eSJed Brown /// ``` 4329df49d7eSJed Brown impl<'a> fmt::Display for CompositeOperator<'a> { 4339df49d7eSJed Brown fn fmt(&self, f: &mut fmt::Formatter) -> fmt::Result { 4349df49d7eSJed Brown self.op_core.fmt(f) 4359df49d7eSJed Brown } 4369df49d7eSJed Brown } 4379df49d7eSJed Brown 4389df49d7eSJed Brown // ----------------------------------------------------------------------------- 4399df49d7eSJed Brown // Core functionality 4409df49d7eSJed Brown // ----------------------------------------------------------------------------- 4419df49d7eSJed Brown impl<'a> OperatorCore<'a> { 442*11544396SJeremy L Thompson // Raw Ceed for error handling 443*11544396SJeremy L Thompson #[doc(hidden)] 444*11544396SJeremy L Thompson fn ceed(&self) -> bind_ceed::Ceed { 445*11544396SJeremy L Thompson unsafe { bind_ceed::CeedOperatorReturnCeed(self.ptr) } 446*11544396SJeremy L Thompson } 447*11544396SJeremy L Thompson 4481142270cSJeremy L Thompson // Error handling 4491142270cSJeremy L Thompson #[doc(hidden)] 4501142270cSJeremy L Thompson fn check_error(&self, ierr: i32) -> crate::Result<i32> { 451*11544396SJeremy L Thompson crate::check_error(|| self.ceed(), ierr) 4521142270cSJeremy L Thompson } 4531142270cSJeremy L Thompson 4549df49d7eSJed Brown // Common implementations 4556f97ff0aSJeremy L Thompson pub fn check(&self) -> crate::Result<i32> { 4566f97ff0aSJeremy L Thompson let ierr = unsafe { bind_ceed::CeedOperatorCheckReady(self.ptr) }; 4576f97ff0aSJeremy L Thompson self.check_error(ierr) 4586f97ff0aSJeremy L Thompson } 4596f97ff0aSJeremy L Thompson 460ea6b5821SJeremy L Thompson pub fn name(&self, name: &str) -> crate::Result<i32> { 461ea6b5821SJeremy L Thompson let name_c = CString::new(name).expect("CString::new failed"); 462ea6b5821SJeremy L Thompson let ierr = unsafe { bind_ceed::CeedOperatorSetName(self.ptr, name_c.as_ptr()) }; 463ea6b5821SJeremy L Thompson self.check_error(ierr) 464ea6b5821SJeremy L Thompson } 465ea6b5821SJeremy L Thompson 4669df49d7eSJed Brown pub fn apply(&self, input: &Vector, output: &mut Vector) -> crate::Result<i32> { 4679df49d7eSJed Brown let ierr = unsafe { 4689df49d7eSJed Brown bind_ceed::CeedOperatorApply( 4699df49d7eSJed Brown self.ptr, 4709df49d7eSJed Brown input.ptr, 4719df49d7eSJed Brown output.ptr, 4729df49d7eSJed Brown bind_ceed::CEED_REQUEST_IMMEDIATE, 4739df49d7eSJed Brown ) 4749df49d7eSJed Brown }; 4751142270cSJeremy L Thompson self.check_error(ierr) 4769df49d7eSJed Brown } 4779df49d7eSJed Brown 4789df49d7eSJed Brown pub fn apply_add(&self, input: &Vector, output: &mut Vector) -> crate::Result<i32> { 4799df49d7eSJed Brown let ierr = unsafe { 4809df49d7eSJed Brown bind_ceed::CeedOperatorApplyAdd( 4819df49d7eSJed Brown self.ptr, 4829df49d7eSJed Brown input.ptr, 4839df49d7eSJed Brown output.ptr, 4849df49d7eSJed Brown bind_ceed::CEED_REQUEST_IMMEDIATE, 4859df49d7eSJed Brown ) 4869df49d7eSJed Brown }; 4871142270cSJeremy L Thompson self.check_error(ierr) 4889df49d7eSJed Brown } 4899df49d7eSJed Brown 4909df49d7eSJed Brown pub fn linear_assemble_diagonal(&self, assembled: &mut Vector) -> crate::Result<i32> { 4919df49d7eSJed Brown let ierr = unsafe { 4929df49d7eSJed Brown bind_ceed::CeedOperatorLinearAssembleDiagonal( 4939df49d7eSJed Brown self.ptr, 4949df49d7eSJed Brown assembled.ptr, 4959df49d7eSJed Brown bind_ceed::CEED_REQUEST_IMMEDIATE, 4969df49d7eSJed Brown ) 4979df49d7eSJed Brown }; 4981142270cSJeremy L Thompson self.check_error(ierr) 4999df49d7eSJed Brown } 5009df49d7eSJed Brown 5019df49d7eSJed Brown pub fn linear_assemble_add_diagonal(&self, assembled: &mut Vector) -> crate::Result<i32> { 5029df49d7eSJed Brown let ierr = unsafe { 5039df49d7eSJed Brown bind_ceed::CeedOperatorLinearAssembleAddDiagonal( 5049df49d7eSJed Brown self.ptr, 5059df49d7eSJed Brown assembled.ptr, 5069df49d7eSJed Brown bind_ceed::CEED_REQUEST_IMMEDIATE, 5079df49d7eSJed Brown ) 5089df49d7eSJed Brown }; 5091142270cSJeremy L Thompson self.check_error(ierr) 5109df49d7eSJed Brown } 5119df49d7eSJed Brown 5129df49d7eSJed Brown pub fn linear_assemble_point_block_diagonal( 5139df49d7eSJed Brown &self, 5149df49d7eSJed Brown assembled: &mut Vector, 5159df49d7eSJed Brown ) -> crate::Result<i32> { 5169df49d7eSJed Brown let ierr = unsafe { 5179df49d7eSJed Brown bind_ceed::CeedOperatorLinearAssemblePointBlockDiagonal( 5189df49d7eSJed Brown self.ptr, 5199df49d7eSJed Brown assembled.ptr, 5209df49d7eSJed Brown bind_ceed::CEED_REQUEST_IMMEDIATE, 5219df49d7eSJed Brown ) 5229df49d7eSJed Brown }; 5231142270cSJeremy L Thompson self.check_error(ierr) 5249df49d7eSJed Brown } 5259df49d7eSJed Brown 5269df49d7eSJed Brown pub fn linear_assemble_add_point_block_diagonal( 5279df49d7eSJed Brown &self, 5289df49d7eSJed Brown assembled: &mut Vector, 5299df49d7eSJed Brown ) -> crate::Result<i32> { 5309df49d7eSJed Brown let ierr = unsafe { 5319df49d7eSJed Brown bind_ceed::CeedOperatorLinearAssembleAddPointBlockDiagonal( 5329df49d7eSJed Brown self.ptr, 5339df49d7eSJed Brown assembled.ptr, 5349df49d7eSJed Brown bind_ceed::CEED_REQUEST_IMMEDIATE, 5359df49d7eSJed Brown ) 5369df49d7eSJed Brown }; 5371142270cSJeremy L Thompson self.check_error(ierr) 5389df49d7eSJed Brown } 5399df49d7eSJed Brown } 5409df49d7eSJed Brown 5419df49d7eSJed Brown // ----------------------------------------------------------------------------- 5429df49d7eSJed Brown // Operator 5439df49d7eSJed Brown // ----------------------------------------------------------------------------- 5449df49d7eSJed Brown impl<'a> Operator<'a> { 5459df49d7eSJed Brown // Constructor 5469df49d7eSJed Brown pub fn create<'b>( 547594ef120SJeremy L Thompson ceed: &crate::Ceed, 5489df49d7eSJed Brown qf: impl Into<QFunctionOpt<'b>>, 5499df49d7eSJed Brown dqf: impl Into<QFunctionOpt<'b>>, 5509df49d7eSJed Brown dqfT: impl Into<QFunctionOpt<'b>>, 5519df49d7eSJed Brown ) -> crate::Result<Self> { 5529df49d7eSJed Brown let mut ptr = std::ptr::null_mut(); 5539df49d7eSJed Brown let ierr = unsafe { 5549df49d7eSJed Brown bind_ceed::CeedOperatorCreate( 5559df49d7eSJed Brown ceed.ptr, 5569df49d7eSJed Brown qf.into().to_raw(), 5579df49d7eSJed Brown dqf.into().to_raw(), 5589df49d7eSJed Brown dqfT.into().to_raw(), 5599df49d7eSJed Brown &mut ptr, 5609df49d7eSJed Brown ) 5619df49d7eSJed Brown }; 5629df49d7eSJed Brown ceed.check_error(ierr)?; 5639df49d7eSJed Brown Ok(Self { 5641142270cSJeremy L Thompson op_core: OperatorCore { 5651142270cSJeremy L Thompson ptr, 5661142270cSJeremy L Thompson _lifeline: PhantomData, 5671142270cSJeremy L Thompson }, 5689df49d7eSJed Brown }) 5699df49d7eSJed Brown } 5709df49d7eSJed Brown 5711142270cSJeremy L Thompson fn from_raw(ptr: bind_ceed::CeedOperator) -> crate::Result<Self> { 5729df49d7eSJed Brown Ok(Self { 5731142270cSJeremy L Thompson op_core: OperatorCore { 5741142270cSJeremy L Thompson ptr, 5751142270cSJeremy L Thompson _lifeline: PhantomData, 5761142270cSJeremy L Thompson }, 5779df49d7eSJed Brown }) 5789df49d7eSJed Brown } 5799df49d7eSJed Brown 580ea6b5821SJeremy L Thompson /// Set name for Operator printing 581ea6b5821SJeremy L Thompson /// 582ea6b5821SJeremy L Thompson /// * 'name' - Name to set 583ea6b5821SJeremy L Thompson /// 584ea6b5821SJeremy L Thompson /// ``` 585ea6b5821SJeremy L Thompson /// # use libceed::prelude::*; 586ea6b5821SJeremy L Thompson /// # fn main() -> libceed::Result<()> { 587ea6b5821SJeremy L Thompson /// # let ceed = libceed::Ceed::default_init(); 588ea6b5821SJeremy L Thompson /// let qf = ceed.q_function_interior_by_name("Mass1DBuild")?; 589ea6b5821SJeremy L Thompson /// 590ea6b5821SJeremy L Thompson /// // Operator field arguments 591ea6b5821SJeremy L Thompson /// let ne = 3; 592ea6b5821SJeremy L Thompson /// let q = 4 as usize; 593ea6b5821SJeremy L Thompson /// let mut ind: Vec<i32> = vec![0; 2 * ne]; 594ea6b5821SJeremy L Thompson /// for i in 0..ne { 595ea6b5821SJeremy L Thompson /// ind[2 * i + 0] = i as i32; 596ea6b5821SJeremy L Thompson /// ind[2 * i + 1] = (i + 1) as i32; 597ea6b5821SJeremy L Thompson /// } 598ea6b5821SJeremy L Thompson /// let r = ceed.elem_restriction(ne, 2, 1, 1, ne + 1, MemType::Host, &ind)?; 599ea6b5821SJeremy L Thompson /// let strides: [i32; 3] = [1, q as i32, q as i32]; 600d9267b76SJeremy L Thompson /// let rq = ceed.strided_elem_restriction(ne, q, 1, q * ne, strides)?; 601ea6b5821SJeremy L Thompson /// 602ea6b5821SJeremy L Thompson /// let b = ceed.basis_tensor_H1_Lagrange(1, 1, 2, q, QuadMode::Gauss)?; 603ea6b5821SJeremy L Thompson /// 604ea6b5821SJeremy L Thompson /// // Operator fields 605ea6b5821SJeremy L Thompson /// let op = ceed 606ea6b5821SJeremy L Thompson /// .operator(&qf, QFunctionOpt::None, QFunctionOpt::None)? 607ea6b5821SJeremy L Thompson /// .name("mass")? 608ea6b5821SJeremy L Thompson /// .field("dx", &r, &b, VectorOpt::Active)? 609ea6b5821SJeremy L Thompson /// .field("weights", ElemRestrictionOpt::None, &b, VectorOpt::None)? 610356036faSJeremy L Thompson /// .field("qdata", &rq, BasisOpt::None, VectorOpt::Active)?; 611ea6b5821SJeremy L Thompson /// # Ok(()) 612ea6b5821SJeremy L Thompson /// # } 613ea6b5821SJeremy L Thompson /// ``` 614ea6b5821SJeremy L Thompson #[allow(unused_mut)] 615ea6b5821SJeremy L Thompson pub fn name(mut self, name: &str) -> crate::Result<Self> { 616ea6b5821SJeremy L Thompson self.op_core.name(name)?; 617ea6b5821SJeremy L Thompson Ok(self) 618ea6b5821SJeremy L Thompson } 619ea6b5821SJeremy L Thompson 6209df49d7eSJed Brown /// Apply Operator to a vector 6219df49d7eSJed Brown /// 6229df49d7eSJed Brown /// * `input` - Input Vector 6239df49d7eSJed Brown /// * `output` - Output Vector 6249df49d7eSJed Brown /// 6259df49d7eSJed Brown /// ``` 6269df49d7eSJed Brown /// # use libceed::prelude::*; 6274d27c890SJeremy L Thompson /// # fn main() -> libceed::Result<()> { 6289df49d7eSJed Brown /// # let ceed = libceed::Ceed::default_init(); 6299df49d7eSJed Brown /// let ne = 4; 6309df49d7eSJed Brown /// let p = 3; 6319df49d7eSJed Brown /// let q = 4; 6329df49d7eSJed Brown /// let ndofs = p * ne - ne + 1; 6339df49d7eSJed Brown /// 6349df49d7eSJed Brown /// // Vectors 635c68be7a2SJeremy L Thompson /// let x = ceed.vector_from_slice(&[-1., -0.5, 0.0, 0.5, 1.0])?; 636c68be7a2SJeremy L Thompson /// let mut qdata = ceed.vector(ne * q)?; 6379df49d7eSJed Brown /// qdata.set_value(0.0); 638c68be7a2SJeremy L Thompson /// let u = ceed.vector_from_slice(&vec![1.0; ndofs])?; 639c68be7a2SJeremy L Thompson /// let mut v = ceed.vector(ndofs)?; 6409df49d7eSJed Brown /// v.set_value(0.0); 6419df49d7eSJed Brown /// 6429df49d7eSJed Brown /// // Restrictions 6439df49d7eSJed Brown /// let mut indx: Vec<i32> = vec![0; 2 * ne]; 6449df49d7eSJed Brown /// for i in 0..ne { 6459df49d7eSJed Brown /// indx[2 * i + 0] = i as i32; 6469df49d7eSJed Brown /// indx[2 * i + 1] = (i + 1) as i32; 6479df49d7eSJed Brown /// } 648c68be7a2SJeremy L Thompson /// let rx = ceed.elem_restriction(ne, 2, 1, 1, ne + 1, MemType::Host, &indx)?; 6499df49d7eSJed Brown /// let mut indu: Vec<i32> = vec![0; p * ne]; 6509df49d7eSJed Brown /// for i in 0..ne { 6519df49d7eSJed Brown /// indu[p * i + 0] = i as i32; 6529df49d7eSJed Brown /// indu[p * i + 1] = (i + 1) as i32; 6539df49d7eSJed Brown /// indu[p * i + 2] = (i + 2) as i32; 6549df49d7eSJed Brown /// } 655c68be7a2SJeremy L Thompson /// let ru = ceed.elem_restriction(ne, 3, 1, 1, ndofs, MemType::Host, &indu)?; 6569df49d7eSJed Brown /// let strides: [i32; 3] = [1, q as i32, q as i32]; 657c68be7a2SJeremy L Thompson /// let rq = ceed.strided_elem_restriction(ne, q, 1, q * ne, strides)?; 6589df49d7eSJed Brown /// 6599df49d7eSJed Brown /// // Bases 660c68be7a2SJeremy L Thompson /// let bx = ceed.basis_tensor_H1_Lagrange(1, 1, 2, q, QuadMode::Gauss)?; 661c68be7a2SJeremy L Thompson /// let bu = ceed.basis_tensor_H1_Lagrange(1, 1, p, q, QuadMode::Gauss)?; 6629df49d7eSJed Brown /// 6639df49d7eSJed Brown /// // Build quadrature data 664c68be7a2SJeremy L Thompson /// let qf_build = ceed.q_function_interior_by_name("Mass1DBuild")?; 665c68be7a2SJeremy L Thompson /// ceed.operator(&qf_build, QFunctionOpt::None, QFunctionOpt::None)? 666c68be7a2SJeremy L Thompson /// .field("dx", &rx, &bx, VectorOpt::Active)? 667c68be7a2SJeremy L Thompson /// .field("weights", ElemRestrictionOpt::None, &bx, VectorOpt::None)? 668356036faSJeremy L Thompson /// .field("qdata", &rq, BasisOpt::None, VectorOpt::Active)? 669c68be7a2SJeremy L Thompson /// .apply(&x, &mut qdata)?; 6709df49d7eSJed Brown /// 6719df49d7eSJed Brown /// // Mass operator 672c68be7a2SJeremy L Thompson /// let qf_mass = ceed.q_function_interior_by_name("MassApply")?; 6739df49d7eSJed Brown /// let op_mass = ceed 674c68be7a2SJeremy L Thompson /// .operator(&qf_mass, QFunctionOpt::None, QFunctionOpt::None)? 675c68be7a2SJeremy L Thompson /// .field("u", &ru, &bu, VectorOpt::Active)? 676356036faSJeremy L Thompson /// .field("qdata", &rq, BasisOpt::None, &qdata)? 677c68be7a2SJeremy L Thompson /// .field("v", &ru, &bu, VectorOpt::Active)?; 6789df49d7eSJed Brown /// 6799df49d7eSJed Brown /// v.set_value(0.0); 680c68be7a2SJeremy L Thompson /// op_mass.apply(&u, &mut v)?; 6819df49d7eSJed Brown /// 6829df49d7eSJed Brown /// // Check 683e78171edSJeremy L Thompson /// let sum: Scalar = v.view()?.iter().sum(); 6844b61a5a0SRezgar Shakeri /// let error: Scalar = (sum - 2.0).abs(); 6859df49d7eSJed Brown /// assert!( 6864b61a5a0SRezgar Shakeri /// error < 50.0 * libceed::EPSILON, 6874b61a5a0SRezgar Shakeri /// "Incorrect interval length computed. Expected: 2.0, Found: {}, Error: {:.12e}", 6884b61a5a0SRezgar Shakeri /// sum, 6894b61a5a0SRezgar Shakeri /// error 6909df49d7eSJed Brown /// ); 691c68be7a2SJeremy L Thompson /// # Ok(()) 692c68be7a2SJeremy L Thompson /// # } 6939df49d7eSJed Brown /// ``` 6949df49d7eSJed Brown pub fn apply(&self, input: &Vector, output: &mut Vector) -> crate::Result<i32> { 6959df49d7eSJed Brown self.op_core.apply(input, output) 6969df49d7eSJed Brown } 6979df49d7eSJed Brown 6989df49d7eSJed Brown /// Apply Operator to a vector and add result to output vector 6999df49d7eSJed Brown /// 7009df49d7eSJed Brown /// * `input` - Input Vector 7019df49d7eSJed Brown /// * `output` - Output Vector 7029df49d7eSJed Brown /// 7039df49d7eSJed Brown /// ``` 7049df49d7eSJed Brown /// # use libceed::prelude::*; 7054d27c890SJeremy L Thompson /// # fn main() -> libceed::Result<()> { 7069df49d7eSJed Brown /// # let ceed = libceed::Ceed::default_init(); 7079df49d7eSJed Brown /// let ne = 4; 7089df49d7eSJed Brown /// let p = 3; 7099df49d7eSJed Brown /// let q = 4; 7109df49d7eSJed Brown /// let ndofs = p * ne - ne + 1; 7119df49d7eSJed Brown /// 7129df49d7eSJed Brown /// // Vectors 713c68be7a2SJeremy L Thompson /// let x = ceed.vector_from_slice(&[-1., -0.5, 0.0, 0.5, 1.0])?; 714c68be7a2SJeremy L Thompson /// let mut qdata = ceed.vector(ne * q)?; 7159df49d7eSJed Brown /// qdata.set_value(0.0); 716c68be7a2SJeremy L Thompson /// let u = ceed.vector_from_slice(&vec![1.0; ndofs])?; 717c68be7a2SJeremy L Thompson /// let mut v = ceed.vector(ndofs)?; 7189df49d7eSJed Brown /// 7199df49d7eSJed Brown /// // Restrictions 7209df49d7eSJed Brown /// let mut indx: Vec<i32> = vec![0; 2 * ne]; 7219df49d7eSJed Brown /// for i in 0..ne { 7229df49d7eSJed Brown /// indx[2 * i + 0] = i as i32; 7239df49d7eSJed Brown /// indx[2 * i + 1] = (i + 1) as i32; 7249df49d7eSJed Brown /// } 725c68be7a2SJeremy L Thompson /// let rx = ceed.elem_restriction(ne, 2, 1, 1, ne + 1, MemType::Host, &indx)?; 7269df49d7eSJed Brown /// let mut indu: Vec<i32> = vec![0; p * ne]; 7279df49d7eSJed Brown /// for i in 0..ne { 7289df49d7eSJed Brown /// indu[p * i + 0] = i as i32; 7299df49d7eSJed Brown /// indu[p * i + 1] = (i + 1) as i32; 7309df49d7eSJed Brown /// indu[p * i + 2] = (i + 2) as i32; 7319df49d7eSJed Brown /// } 732c68be7a2SJeremy L Thompson /// let ru = ceed.elem_restriction(ne, 3, 1, 1, ndofs, MemType::Host, &indu)?; 7339df49d7eSJed Brown /// let strides: [i32; 3] = [1, q as i32, q as i32]; 734c68be7a2SJeremy L Thompson /// let rq = ceed.strided_elem_restriction(ne, q, 1, q * ne, strides)?; 7359df49d7eSJed Brown /// 7369df49d7eSJed Brown /// // Bases 737c68be7a2SJeremy L Thompson /// let bx = ceed.basis_tensor_H1_Lagrange(1, 1, 2, q, QuadMode::Gauss)?; 738c68be7a2SJeremy L Thompson /// let bu = ceed.basis_tensor_H1_Lagrange(1, 1, p, q, QuadMode::Gauss)?; 7399df49d7eSJed Brown /// 7409df49d7eSJed Brown /// // Build quadrature data 741c68be7a2SJeremy L Thompson /// let qf_build = ceed.q_function_interior_by_name("Mass1DBuild")?; 742c68be7a2SJeremy L Thompson /// ceed.operator(&qf_build, QFunctionOpt::None, QFunctionOpt::None)? 743c68be7a2SJeremy L Thompson /// .field("dx", &rx, &bx, VectorOpt::Active)? 744c68be7a2SJeremy L Thompson /// .field("weights", ElemRestrictionOpt::None, &bx, VectorOpt::None)? 745356036faSJeremy L Thompson /// .field("qdata", &rq, BasisOpt::None, VectorOpt::Active)? 746c68be7a2SJeremy L Thompson /// .apply(&x, &mut qdata)?; 7479df49d7eSJed Brown /// 7489df49d7eSJed Brown /// // Mass operator 749c68be7a2SJeremy L Thompson /// let qf_mass = ceed.q_function_interior_by_name("MassApply")?; 7509df49d7eSJed Brown /// let op_mass = ceed 751c68be7a2SJeremy L Thompson /// .operator(&qf_mass, QFunctionOpt::None, QFunctionOpt::None)? 752c68be7a2SJeremy L Thompson /// .field("u", &ru, &bu, VectorOpt::Active)? 753356036faSJeremy L Thompson /// .field("qdata", &rq, BasisOpt::None, &qdata)? 754c68be7a2SJeremy L Thompson /// .field("v", &ru, &bu, VectorOpt::Active)?; 7559df49d7eSJed Brown /// 7569df49d7eSJed Brown /// v.set_value(1.0); 757c68be7a2SJeremy L Thompson /// op_mass.apply_add(&u, &mut v)?; 7589df49d7eSJed Brown /// 7599df49d7eSJed Brown /// // Check 760e78171edSJeremy L Thompson /// let sum: Scalar = v.view()?.iter().sum(); 7619df49d7eSJed Brown /// assert!( 76280a9ef05SNatalie Beams /// (sum - (2.0 + ndofs as Scalar)).abs() < 10.0 * libceed::EPSILON, 7639df49d7eSJed Brown /// "Incorrect interval length computed and added" 7649df49d7eSJed Brown /// ); 765c68be7a2SJeremy L Thompson /// # Ok(()) 766c68be7a2SJeremy L Thompson /// # } 7679df49d7eSJed Brown /// ``` 7689df49d7eSJed Brown pub fn apply_add(&self, input: &Vector, output: &mut Vector) -> crate::Result<i32> { 7699df49d7eSJed Brown self.op_core.apply_add(input, output) 7709df49d7eSJed Brown } 7719df49d7eSJed Brown 7729df49d7eSJed Brown /// Provide a field to a Operator for use by its QFunction 7739df49d7eSJed Brown /// 7749df49d7eSJed Brown /// * `fieldname` - Name of the field (to be matched with the name used by 7759df49d7eSJed Brown /// the QFunction) 7769df49d7eSJed Brown /// * `r` - ElemRestriction 777356036faSJeremy L Thompson /// * `b` - Basis in which the field resides or `BasisOpt::None` 778356036faSJeremy L Thompson /// * `v` - Vector to be used by Operator or `VectorOpt::Active` if field 779356036faSJeremy L Thompson /// is active or `VectorOpt::None` if using `Weight` with the 7809df49d7eSJed Brown /// QFunction 7819df49d7eSJed Brown /// 7829df49d7eSJed Brown /// 7839df49d7eSJed Brown /// ``` 7849df49d7eSJed Brown /// # use libceed::prelude::*; 7854d27c890SJeremy L Thompson /// # fn main() -> libceed::Result<()> { 7869df49d7eSJed Brown /// # let ceed = libceed::Ceed::default_init(); 787c68be7a2SJeremy L Thompson /// let qf = ceed.q_function_interior_by_name("Mass1DBuild")?; 788c68be7a2SJeremy L Thompson /// let mut op = ceed.operator(&qf, QFunctionOpt::None, QFunctionOpt::None)?; 7899df49d7eSJed Brown /// 7909df49d7eSJed Brown /// // Operator field arguments 7919df49d7eSJed Brown /// let ne = 3; 7929df49d7eSJed Brown /// let q = 4; 7939df49d7eSJed Brown /// let mut ind: Vec<i32> = vec![0; 2 * ne]; 7949df49d7eSJed Brown /// for i in 0..ne { 7959df49d7eSJed Brown /// ind[2 * i + 0] = i as i32; 7969df49d7eSJed Brown /// ind[2 * i + 1] = (i + 1) as i32; 7979df49d7eSJed Brown /// } 798c68be7a2SJeremy L Thompson /// let r = ceed.elem_restriction(ne, 2, 1, 1, ne + 1, MemType::Host, &ind)?; 7999df49d7eSJed Brown /// 800c68be7a2SJeremy L Thompson /// let b = ceed.basis_tensor_H1_Lagrange(1, 1, 2, q, QuadMode::Gauss)?; 8019df49d7eSJed Brown /// 8029df49d7eSJed Brown /// // Operator field 803c68be7a2SJeremy L Thompson /// op = op.field("dx", &r, &b, VectorOpt::Active)?; 804c68be7a2SJeremy L Thompson /// # Ok(()) 805c68be7a2SJeremy L Thompson /// # } 8069df49d7eSJed Brown /// ``` 8079df49d7eSJed Brown #[allow(unused_mut)] 8089df49d7eSJed Brown pub fn field<'b>( 8099df49d7eSJed Brown mut self, 8109df49d7eSJed Brown fieldname: &str, 8119df49d7eSJed Brown r: impl Into<ElemRestrictionOpt<'b>>, 8129df49d7eSJed Brown b: impl Into<BasisOpt<'b>>, 8139df49d7eSJed Brown v: impl Into<VectorOpt<'b>>, 8149df49d7eSJed Brown ) -> crate::Result<Self> { 8159df49d7eSJed Brown let fieldname = CString::new(fieldname).expect("CString::new failed"); 8169df49d7eSJed Brown let fieldname = fieldname.as_ptr() as *const i8; 8179df49d7eSJed Brown let ierr = unsafe { 8189df49d7eSJed Brown bind_ceed::CeedOperatorSetField( 8199df49d7eSJed Brown self.op_core.ptr, 8209df49d7eSJed Brown fieldname, 8219df49d7eSJed Brown r.into().to_raw(), 8229df49d7eSJed Brown b.into().to_raw(), 8239df49d7eSJed Brown v.into().to_raw(), 8249df49d7eSJed Brown ) 8259df49d7eSJed Brown }; 8261142270cSJeremy L Thompson self.op_core.check_error(ierr)?; 8279df49d7eSJed Brown Ok(self) 8289df49d7eSJed Brown } 8299df49d7eSJed Brown 83008778c6fSJeremy L Thompson /// Get a slice of Operator inputs 83108778c6fSJeremy L Thompson /// 83208778c6fSJeremy L Thompson /// ``` 83308778c6fSJeremy L Thompson /// # use libceed::prelude::*; 8344d27c890SJeremy L Thompson /// # fn main() -> libceed::Result<()> { 83508778c6fSJeremy L Thompson /// # let ceed = libceed::Ceed::default_init(); 83608778c6fSJeremy L Thompson /// let qf = ceed.q_function_interior_by_name("Mass1DBuild")?; 83708778c6fSJeremy L Thompson /// 83808778c6fSJeremy L Thompson /// // Operator field arguments 83908778c6fSJeremy L Thompson /// let ne = 3; 84008778c6fSJeremy L Thompson /// let q = 4 as usize; 84108778c6fSJeremy L Thompson /// let mut ind: Vec<i32> = vec![0; 2 * ne]; 84208778c6fSJeremy L Thompson /// for i in 0..ne { 84308778c6fSJeremy L Thompson /// ind[2 * i + 0] = i as i32; 84408778c6fSJeremy L Thompson /// ind[2 * i + 1] = (i + 1) as i32; 84508778c6fSJeremy L Thompson /// } 84608778c6fSJeremy L Thompson /// let r = ceed.elem_restriction(ne, 2, 1, 1, ne + 1, MemType::Host, &ind)?; 84708778c6fSJeremy L Thompson /// let strides: [i32; 3] = [1, q as i32, q as i32]; 848d9267b76SJeremy L Thompson /// let rq = ceed.strided_elem_restriction(ne, q, 1, q * ne, strides)?; 84908778c6fSJeremy L Thompson /// 85008778c6fSJeremy L Thompson /// let b = ceed.basis_tensor_H1_Lagrange(1, 1, 2, q, QuadMode::Gauss)?; 85108778c6fSJeremy L Thompson /// 85208778c6fSJeremy L Thompson /// // Operator fields 85308778c6fSJeremy L Thompson /// let op = ceed 85408778c6fSJeremy L Thompson /// .operator(&qf, QFunctionOpt::None, QFunctionOpt::None)? 85508778c6fSJeremy L Thompson /// .field("dx", &r, &b, VectorOpt::Active)? 85608778c6fSJeremy L Thompson /// .field("weights", ElemRestrictionOpt::None, &b, VectorOpt::None)? 857356036faSJeremy L Thompson /// .field("qdata", &rq, BasisOpt::None, VectorOpt::Active)?; 85808778c6fSJeremy L Thompson /// 85908778c6fSJeremy L Thompson /// let inputs = op.inputs()?; 86008778c6fSJeremy L Thompson /// 86108778c6fSJeremy L Thompson /// assert_eq!(inputs.len(), 2, "Incorrect inputs array"); 86208778c6fSJeremy L Thompson /// # Ok(()) 86308778c6fSJeremy L Thompson /// # } 86408778c6fSJeremy L Thompson /// ``` 865e03fef56SJeremy L Thompson pub fn inputs(&self) -> crate::Result<Vec<crate::OperatorField>> { 86608778c6fSJeremy L Thompson // Get array of raw C pointers for inputs 86708778c6fSJeremy L Thompson let mut num_inputs = 0; 86808778c6fSJeremy L Thompson let mut inputs_ptr = std::ptr::null_mut(); 86908778c6fSJeremy L Thompson let ierr = unsafe { 87008778c6fSJeremy L Thompson bind_ceed::CeedOperatorGetFields( 87108778c6fSJeremy L Thompson self.op_core.ptr, 87208778c6fSJeremy L Thompson &mut num_inputs, 87308778c6fSJeremy L Thompson &mut inputs_ptr, 87408778c6fSJeremy L Thompson std::ptr::null_mut() as *mut bind_ceed::CeedInt, 87508778c6fSJeremy L Thompson std::ptr::null_mut() as *mut *mut bind_ceed::CeedOperatorField, 87608778c6fSJeremy L Thompson ) 87708778c6fSJeremy L Thompson }; 87808778c6fSJeremy L Thompson self.op_core.check_error(ierr)?; 87908778c6fSJeremy L Thompson // Convert raw C pointers to fixed length slice 88008778c6fSJeremy L Thompson let inputs_slice = unsafe { 88108778c6fSJeremy L Thompson std::slice::from_raw_parts( 882e03fef56SJeremy L Thompson inputs_ptr as *mut bind_ceed::CeedOperatorField, 88308778c6fSJeremy L Thompson num_inputs as usize, 88408778c6fSJeremy L Thompson ) 88508778c6fSJeremy L Thompson }; 886e03fef56SJeremy L Thompson // And finally build vec 887e03fef56SJeremy L Thompson let ceed = { 888e03fef56SJeremy L Thompson let mut ptr = std::ptr::null_mut(); 889e03fef56SJeremy L Thompson let mut ptr_copy = std::ptr::null_mut(); 890e03fef56SJeremy L Thompson unsafe { 891e03fef56SJeremy L Thompson bind_ceed::CeedOperatorGetCeed(self.op_core.ptr, &mut ptr); 892e03fef56SJeremy L Thompson bind_ceed::CeedReferenceCopy(ptr, &mut ptr_copy); // refcount 893e03fef56SJeremy L Thompson } 894e03fef56SJeremy L Thompson crate::Ceed { ptr } 895e03fef56SJeremy L Thompson }; 896e03fef56SJeremy L Thompson let inputs = (0..num_inputs as usize) 897e03fef56SJeremy L Thompson .map(|i| crate::OperatorField::from_raw(inputs_slice[i], ceed.clone())) 898e03fef56SJeremy L Thompson .collect::<crate::Result<Vec<_>>>()?; 899e03fef56SJeremy L Thompson Ok(inputs) 90008778c6fSJeremy L Thompson } 90108778c6fSJeremy L Thompson 90208778c6fSJeremy L Thompson /// Get a slice of Operator outputs 90308778c6fSJeremy L Thompson /// 90408778c6fSJeremy L Thompson /// ``` 90508778c6fSJeremy L Thompson /// # use libceed::prelude::*; 9064d27c890SJeremy L Thompson /// # fn main() -> libceed::Result<()> { 90708778c6fSJeremy L Thompson /// # let ceed = libceed::Ceed::default_init(); 90808778c6fSJeremy L Thompson /// let qf = ceed.q_function_interior_by_name("Mass1DBuild")?; 90908778c6fSJeremy L Thompson /// 91008778c6fSJeremy L Thompson /// // Operator field arguments 91108778c6fSJeremy L Thompson /// let ne = 3; 91208778c6fSJeremy L Thompson /// let q = 4 as usize; 91308778c6fSJeremy L Thompson /// let mut ind: Vec<i32> = vec![0; 2 * ne]; 91408778c6fSJeremy L Thompson /// for i in 0..ne { 91508778c6fSJeremy L Thompson /// ind[2 * i + 0] = i as i32; 91608778c6fSJeremy L Thompson /// ind[2 * i + 1] = (i + 1) as i32; 91708778c6fSJeremy L Thompson /// } 91808778c6fSJeremy L Thompson /// let r = ceed.elem_restriction(ne, 2, 1, 1, ne + 1, MemType::Host, &ind)?; 91908778c6fSJeremy L Thompson /// let strides: [i32; 3] = [1, q as i32, q as i32]; 920d9267b76SJeremy L Thompson /// let rq = ceed.strided_elem_restriction(ne, q, 1, q * ne, strides)?; 92108778c6fSJeremy L Thompson /// 92208778c6fSJeremy L Thompson /// let b = ceed.basis_tensor_H1_Lagrange(1, 1, 2, q, QuadMode::Gauss)?; 92308778c6fSJeremy L Thompson /// 92408778c6fSJeremy L Thompson /// // Operator fields 92508778c6fSJeremy L Thompson /// let op = ceed 92608778c6fSJeremy L Thompson /// .operator(&qf, QFunctionOpt::None, QFunctionOpt::None)? 92708778c6fSJeremy L Thompson /// .field("dx", &r, &b, VectorOpt::Active)? 92808778c6fSJeremy L Thompson /// .field("weights", ElemRestrictionOpt::None, &b, VectorOpt::None)? 929356036faSJeremy L Thompson /// .field("qdata", &rq, BasisOpt::None, VectorOpt::Active)?; 93008778c6fSJeremy L Thompson /// 93108778c6fSJeremy L Thompson /// let outputs = op.outputs()?; 93208778c6fSJeremy L Thompson /// 93308778c6fSJeremy L Thompson /// assert_eq!(outputs.len(), 1, "Incorrect outputs array"); 93408778c6fSJeremy L Thompson /// # Ok(()) 93508778c6fSJeremy L Thompson /// # } 93608778c6fSJeremy L Thompson /// ``` 937e03fef56SJeremy L Thompson pub fn outputs(&self) -> crate::Result<Vec<crate::OperatorField>> { 93808778c6fSJeremy L Thompson // Get array of raw C pointers for outputs 93908778c6fSJeremy L Thompson let mut num_outputs = 0; 94008778c6fSJeremy L Thompson let mut outputs_ptr = std::ptr::null_mut(); 94108778c6fSJeremy L Thompson let ierr = unsafe { 94208778c6fSJeremy L Thompson bind_ceed::CeedOperatorGetFields( 94308778c6fSJeremy L Thompson self.op_core.ptr, 94408778c6fSJeremy L Thompson std::ptr::null_mut() as *mut bind_ceed::CeedInt, 94508778c6fSJeremy L Thompson std::ptr::null_mut() as *mut *mut bind_ceed::CeedOperatorField, 94608778c6fSJeremy L Thompson &mut num_outputs, 94708778c6fSJeremy L Thompson &mut outputs_ptr, 94808778c6fSJeremy L Thompson ) 94908778c6fSJeremy L Thompson }; 95008778c6fSJeremy L Thompson self.op_core.check_error(ierr)?; 95108778c6fSJeremy L Thompson // Convert raw C pointers to fixed length slice 95208778c6fSJeremy L Thompson let outputs_slice = unsafe { 95308778c6fSJeremy L Thompson std::slice::from_raw_parts( 954e03fef56SJeremy L Thompson outputs_ptr as *mut bind_ceed::CeedOperatorField, 95508778c6fSJeremy L Thompson num_outputs as usize, 95608778c6fSJeremy L Thompson ) 95708778c6fSJeremy L Thompson }; 958e03fef56SJeremy L Thompson // And finally build vec 959e03fef56SJeremy L Thompson let ceed = { 960e03fef56SJeremy L Thompson let mut ptr = std::ptr::null_mut(); 961e03fef56SJeremy L Thompson let mut ptr_copy = std::ptr::null_mut(); 962e03fef56SJeremy L Thompson unsafe { 963e03fef56SJeremy L Thompson bind_ceed::CeedOperatorGetCeed(self.op_core.ptr, &mut ptr); 964e03fef56SJeremy L Thompson bind_ceed::CeedReferenceCopy(ptr, &mut ptr_copy); // refcount 965e03fef56SJeremy L Thompson } 966e03fef56SJeremy L Thompson crate::Ceed { ptr } 967e03fef56SJeremy L Thompson }; 968e03fef56SJeremy L Thompson let outputs = (0..num_outputs as usize) 969e03fef56SJeremy L Thompson .map(|i| crate::OperatorField::from_raw(outputs_slice[i], ceed.clone())) 970e03fef56SJeremy L Thompson .collect::<crate::Result<Vec<_>>>()?; 971e03fef56SJeremy L Thompson Ok(outputs) 97208778c6fSJeremy L Thompson } 97308778c6fSJeremy L Thompson 9746f97ff0aSJeremy L Thompson /// Check if Operator is setup correctly 9756f97ff0aSJeremy L Thompson /// 9766f97ff0aSJeremy L Thompson /// ``` 9776f97ff0aSJeremy L Thompson /// # use libceed::prelude::*; 9786f97ff0aSJeremy L Thompson /// # fn main() -> libceed::Result<()> { 9796f97ff0aSJeremy L Thompson /// # let ceed = libceed::Ceed::default_init(); 9806f97ff0aSJeremy L Thompson /// let ne = 4; 9816f97ff0aSJeremy L Thompson /// let p = 3; 9826f97ff0aSJeremy L Thompson /// let q = 4; 9836f97ff0aSJeremy L Thompson /// let ndofs = p * ne - ne + 1; 9846f97ff0aSJeremy L Thompson /// 9856f97ff0aSJeremy L Thompson /// // Restrictions 9866f97ff0aSJeremy L Thompson /// let mut indx: Vec<i32> = vec![0; 2 * ne]; 9876f97ff0aSJeremy L Thompson /// for i in 0..ne { 9886f97ff0aSJeremy L Thompson /// indx[2 * i + 0] = i as i32; 9896f97ff0aSJeremy L Thompson /// indx[2 * i + 1] = (i + 1) as i32; 9906f97ff0aSJeremy L Thompson /// } 9916f97ff0aSJeremy L Thompson /// let rx = ceed.elem_restriction(ne, 2, 1, 1, ne + 1, MemType::Host, &indx)?; 9926f97ff0aSJeremy L Thompson /// let strides: [i32; 3] = [1, q as i32, q as i32]; 9936f97ff0aSJeremy L Thompson /// let rq = ceed.strided_elem_restriction(ne, q, 1, q * ne, strides)?; 9946f97ff0aSJeremy L Thompson /// 9956f97ff0aSJeremy L Thompson /// // Bases 9966f97ff0aSJeremy L Thompson /// let bx = ceed.basis_tensor_H1_Lagrange(1, 1, 2, q, QuadMode::Gauss)?; 9976f97ff0aSJeremy L Thompson /// 9986f97ff0aSJeremy L Thompson /// // Build quadrature data 9996f97ff0aSJeremy L Thompson /// let qf_build = ceed.q_function_interior_by_name("Mass1DBuild")?; 10006f97ff0aSJeremy L Thompson /// let op_build = ceed 10016f97ff0aSJeremy L Thompson /// .operator(&qf_build, QFunctionOpt::None, QFunctionOpt::None)? 10026f97ff0aSJeremy L Thompson /// .field("dx", &rx, &bx, VectorOpt::Active)? 10036f97ff0aSJeremy L Thompson /// .field("weights", ElemRestrictionOpt::None, &bx, VectorOpt::None)? 1004356036faSJeremy L Thompson /// .field("qdata", &rq, BasisOpt::None, VectorOpt::Active)?; 10056f97ff0aSJeremy L Thompson /// 10066f97ff0aSJeremy L Thompson /// // Check 10076f97ff0aSJeremy L Thompson /// assert!(op_build.check().is_ok(), "Operator setup incorrectly"); 10086f97ff0aSJeremy L Thompson /// # Ok(()) 10096f97ff0aSJeremy L Thompson /// # } 10106f97ff0aSJeremy L Thompson /// ``` 10116f97ff0aSJeremy L Thompson pub fn check(self) -> crate::Result<Self> { 10126f97ff0aSJeremy L Thompson self.op_core.check()?; 10136f97ff0aSJeremy L Thompson Ok(self) 10146f97ff0aSJeremy L Thompson } 10156f97ff0aSJeremy L Thompson 10163f2759cfSJeremy L Thompson /// Get the number of elements associated with an Operator 10173f2759cfSJeremy L Thompson /// 10183f2759cfSJeremy L Thompson /// 10193f2759cfSJeremy L Thompson /// ``` 10203f2759cfSJeremy L Thompson /// # use libceed::prelude::*; 10214d27c890SJeremy L Thompson /// # fn main() -> libceed::Result<()> { 10223f2759cfSJeremy L Thompson /// # let ceed = libceed::Ceed::default_init(); 10233f2759cfSJeremy L Thompson /// let qf = ceed.q_function_interior_by_name("Mass1DBuild")?; 10243f2759cfSJeremy L Thompson /// let mut op = ceed.operator(&qf, QFunctionOpt::None, QFunctionOpt::None)?; 10253f2759cfSJeremy L Thompson /// 10263f2759cfSJeremy L Thompson /// // Operator field arguments 10273f2759cfSJeremy L Thompson /// let ne = 3; 10283f2759cfSJeremy L Thompson /// let q = 4; 10293f2759cfSJeremy L Thompson /// let mut ind: Vec<i32> = vec![0; 2 * ne]; 10303f2759cfSJeremy L Thompson /// for i in 0..ne { 10313f2759cfSJeremy L Thompson /// ind[2 * i + 0] = i as i32; 10323f2759cfSJeremy L Thompson /// ind[2 * i + 1] = (i + 1) as i32; 10333f2759cfSJeremy L Thompson /// } 10343f2759cfSJeremy L Thompson /// let r = ceed.elem_restriction(ne, 2, 1, 1, ne + 1, MemType::Host, &ind)?; 10353f2759cfSJeremy L Thompson /// 10363f2759cfSJeremy L Thompson /// let b = ceed.basis_tensor_H1_Lagrange(1, 1, 2, q, QuadMode::Gauss)?; 10373f2759cfSJeremy L Thompson /// 10383f2759cfSJeremy L Thompson /// // Operator field 10393f2759cfSJeremy L Thompson /// op = op.field("dx", &r, &b, VectorOpt::Active)?; 10403f2759cfSJeremy L Thompson /// 10413f2759cfSJeremy L Thompson /// // Check number of elements 10423f2759cfSJeremy L Thompson /// assert_eq!(op.num_elements(), ne, "Incorrect number of elements"); 10433f2759cfSJeremy L Thompson /// # Ok(()) 10443f2759cfSJeremy L Thompson /// # } 10453f2759cfSJeremy L Thompson /// ``` 10463f2759cfSJeremy L Thompson pub fn num_elements(&self) -> usize { 10473f2759cfSJeremy L Thompson let mut nelem = 0; 10483f2759cfSJeremy L Thompson unsafe { bind_ceed::CeedOperatorGetNumElements(self.op_core.ptr, &mut nelem) }; 10493f2759cfSJeremy L Thompson usize::try_from(nelem).unwrap() 10503f2759cfSJeremy L Thompson } 10513f2759cfSJeremy L Thompson 10523f2759cfSJeremy L Thompson /// Get the number of quadrature points associated with each element of 10533f2759cfSJeremy L Thompson /// an Operator 10543f2759cfSJeremy L Thompson /// 10553f2759cfSJeremy L Thompson /// 10563f2759cfSJeremy L Thompson /// ``` 10573f2759cfSJeremy L Thompson /// # use libceed::prelude::*; 10584d27c890SJeremy L Thompson /// # fn main() -> libceed::Result<()> { 10593f2759cfSJeremy L Thompson /// # let ceed = libceed::Ceed::default_init(); 10603f2759cfSJeremy L Thompson /// let qf = ceed.q_function_interior_by_name("Mass1DBuild")?; 10613f2759cfSJeremy L Thompson /// let mut op = ceed.operator(&qf, QFunctionOpt::None, QFunctionOpt::None)?; 10623f2759cfSJeremy L Thompson /// 10633f2759cfSJeremy L Thompson /// // Operator field arguments 10643f2759cfSJeremy L Thompson /// let ne = 3; 10653f2759cfSJeremy L Thompson /// let q = 4; 10663f2759cfSJeremy L Thompson /// let mut ind: Vec<i32> = vec![0; 2 * ne]; 10673f2759cfSJeremy L Thompson /// for i in 0..ne { 10683f2759cfSJeremy L Thompson /// ind[2 * i + 0] = i as i32; 10693f2759cfSJeremy L Thompson /// ind[2 * i + 1] = (i + 1) as i32; 10703f2759cfSJeremy L Thompson /// } 10713f2759cfSJeremy L Thompson /// let r = ceed.elem_restriction(ne, 2, 1, 1, ne + 1, MemType::Host, &ind)?; 10723f2759cfSJeremy L Thompson /// 10733f2759cfSJeremy L Thompson /// let b = ceed.basis_tensor_H1_Lagrange(1, 1, 2, q, QuadMode::Gauss)?; 10743f2759cfSJeremy L Thompson /// 10753f2759cfSJeremy L Thompson /// // Operator field 10763f2759cfSJeremy L Thompson /// op = op.field("dx", &r, &b, VectorOpt::Active)?; 10773f2759cfSJeremy L Thompson /// 10783f2759cfSJeremy L Thompson /// // Check number of quadrature points 10793f2759cfSJeremy L Thompson /// assert_eq!( 10803f2759cfSJeremy L Thompson /// op.num_quadrature_points(), 10813f2759cfSJeremy L Thompson /// q, 10823f2759cfSJeremy L Thompson /// "Incorrect number of quadrature points" 10833f2759cfSJeremy L Thompson /// ); 10843f2759cfSJeremy L Thompson /// # Ok(()) 10853f2759cfSJeremy L Thompson /// # } 10863f2759cfSJeremy L Thompson /// ``` 10873f2759cfSJeremy L Thompson pub fn num_quadrature_points(&self) -> usize { 10883f2759cfSJeremy L Thompson let mut Q = 0; 10893f2759cfSJeremy L Thompson unsafe { bind_ceed::CeedOperatorGetNumQuadraturePoints(self.op_core.ptr, &mut Q) }; 10903f2759cfSJeremy L Thompson usize::try_from(Q).unwrap() 10913f2759cfSJeremy L Thompson } 10923f2759cfSJeremy L Thompson 10939df49d7eSJed Brown /// Assemble the diagonal of a square linear Operator 10949df49d7eSJed Brown /// 10959df49d7eSJed Brown /// This overwrites a Vector with the diagonal of a linear Operator. 10969df49d7eSJed Brown /// 10979df49d7eSJed Brown /// Note: Currently only non-composite Operators with a single field and 10989df49d7eSJed Brown /// composite Operators with single field sub-operators are supported. 10999df49d7eSJed Brown /// 11009df49d7eSJed Brown /// * `op` - Operator to assemble QFunction 11019df49d7eSJed Brown /// * `assembled` - Vector to store assembled Operator diagonal 11029df49d7eSJed Brown /// 11039df49d7eSJed Brown /// ``` 11049df49d7eSJed Brown /// # use libceed::prelude::*; 11054d27c890SJeremy L Thompson /// # fn main() -> libceed::Result<()> { 11069df49d7eSJed Brown /// # let ceed = libceed::Ceed::default_init(); 11079df49d7eSJed Brown /// let ne = 4; 11089df49d7eSJed Brown /// let p = 3; 11099df49d7eSJed Brown /// let q = 4; 11109df49d7eSJed Brown /// let ndofs = p * ne - ne + 1; 11119df49d7eSJed Brown /// 11129df49d7eSJed Brown /// // Vectors 1113c68be7a2SJeremy L Thompson /// let x = ceed.vector_from_slice(&[-1., -0.5, 0.0, 0.5, 1.0])?; 1114c68be7a2SJeremy L Thompson /// let mut qdata = ceed.vector(ne * q)?; 11159df49d7eSJed Brown /// qdata.set_value(0.0); 1116c68be7a2SJeremy L Thompson /// let mut u = ceed.vector(ndofs)?; 11179df49d7eSJed Brown /// u.set_value(1.0); 1118c68be7a2SJeremy L Thompson /// let mut v = ceed.vector(ndofs)?; 11199df49d7eSJed Brown /// v.set_value(0.0); 11209df49d7eSJed Brown /// 11219df49d7eSJed Brown /// // Restrictions 11229df49d7eSJed Brown /// let mut indx: Vec<i32> = vec![0; 2 * ne]; 11239df49d7eSJed Brown /// for i in 0..ne { 11249df49d7eSJed Brown /// indx[2 * i + 0] = i as i32; 11259df49d7eSJed Brown /// indx[2 * i + 1] = (i + 1) as i32; 11269df49d7eSJed Brown /// } 1127c68be7a2SJeremy L Thompson /// let rx = ceed.elem_restriction(ne, 2, 1, 1, ne + 1, MemType::Host, &indx)?; 11289df49d7eSJed Brown /// let mut indu: Vec<i32> = vec![0; p * ne]; 11299df49d7eSJed Brown /// for i in 0..ne { 11309df49d7eSJed Brown /// indu[p * i + 0] = (2 * i) as i32; 11319df49d7eSJed Brown /// indu[p * i + 1] = (2 * i + 1) as i32; 11329df49d7eSJed Brown /// indu[p * i + 2] = (2 * i + 2) as i32; 11339df49d7eSJed Brown /// } 1134c68be7a2SJeremy L Thompson /// let ru = ceed.elem_restriction(ne, p, 1, 1, ndofs, MemType::Host, &indu)?; 11359df49d7eSJed Brown /// let strides: [i32; 3] = [1, q as i32, q as i32]; 1136c68be7a2SJeremy L Thompson /// let rq = ceed.strided_elem_restriction(ne, q, 1, q * ne, strides)?; 11379df49d7eSJed Brown /// 11389df49d7eSJed Brown /// // Bases 1139c68be7a2SJeremy L Thompson /// let bx = ceed.basis_tensor_H1_Lagrange(1, 1, 2, q, QuadMode::Gauss)?; 1140c68be7a2SJeremy L Thompson /// let bu = ceed.basis_tensor_H1_Lagrange(1, 1, p, q, QuadMode::Gauss)?; 11419df49d7eSJed Brown /// 11429df49d7eSJed Brown /// // Build quadrature data 1143c68be7a2SJeremy L Thompson /// let qf_build = ceed.q_function_interior_by_name("Mass1DBuild")?; 1144c68be7a2SJeremy L Thompson /// ceed.operator(&qf_build, QFunctionOpt::None, QFunctionOpt::None)? 1145c68be7a2SJeremy L Thompson /// .field("dx", &rx, &bx, VectorOpt::Active)? 1146c68be7a2SJeremy L Thompson /// .field("weights", ElemRestrictionOpt::None, &bx, VectorOpt::None)? 1147356036faSJeremy L Thompson /// .field("qdata", &rq, BasisOpt::None, VectorOpt::Active)? 1148c68be7a2SJeremy L Thompson /// .apply(&x, &mut qdata)?; 11499df49d7eSJed Brown /// 11509df49d7eSJed Brown /// // Mass operator 1151c68be7a2SJeremy L Thompson /// let qf_mass = ceed.q_function_interior_by_name("MassApply")?; 11529df49d7eSJed Brown /// let op_mass = ceed 1153c68be7a2SJeremy L Thompson /// .operator(&qf_mass, QFunctionOpt::None, QFunctionOpt::None)? 1154c68be7a2SJeremy L Thompson /// .field("u", &ru, &bu, VectorOpt::Active)? 1155356036faSJeremy L Thompson /// .field("qdata", &rq, BasisOpt::None, &qdata)? 1156c68be7a2SJeremy L Thompson /// .field("v", &ru, &bu, VectorOpt::Active)?; 11579df49d7eSJed Brown /// 11589df49d7eSJed Brown /// // Diagonal 1159c68be7a2SJeremy L Thompson /// let mut diag = ceed.vector(ndofs)?; 11609df49d7eSJed Brown /// diag.set_value(0.0); 1161c68be7a2SJeremy L Thompson /// op_mass.linear_assemble_diagonal(&mut diag)?; 11629df49d7eSJed Brown /// 11639df49d7eSJed Brown /// // Manual diagonal computation 1164c68be7a2SJeremy L Thompson /// let mut true_diag = ceed.vector(ndofs)?; 11659c774eddSJeremy L Thompson /// true_diag.set_value(0.0)?; 11669df49d7eSJed Brown /// for i in 0..ndofs { 11679df49d7eSJed Brown /// u.set_value(0.0); 11689df49d7eSJed Brown /// { 1169e78171edSJeremy L Thompson /// let mut u_array = u.view_mut()?; 11709df49d7eSJed Brown /// u_array[i] = 1.; 11719df49d7eSJed Brown /// } 11729df49d7eSJed Brown /// 1173c68be7a2SJeremy L Thompson /// op_mass.apply(&u, &mut v)?; 11749df49d7eSJed Brown /// 11759df49d7eSJed Brown /// { 11769c774eddSJeremy L Thompson /// let v_array = v.view()?; 1177e78171edSJeremy L Thompson /// let mut true_array = true_diag.view_mut()?; 11789df49d7eSJed Brown /// true_array[i] = v_array[i]; 11799df49d7eSJed Brown /// } 11809df49d7eSJed Brown /// } 11819df49d7eSJed Brown /// 11829df49d7eSJed Brown /// // Check 1183e78171edSJeremy L Thompson /// diag.view()? 11849df49d7eSJed Brown /// .iter() 1185e78171edSJeremy L Thompson /// .zip(true_diag.view()?.iter()) 11869df49d7eSJed Brown /// .for_each(|(computed, actual)| { 11879df49d7eSJed Brown /// assert!( 118880a9ef05SNatalie Beams /// (*computed - *actual).abs() < 10.0 * libceed::EPSILON, 11899df49d7eSJed Brown /// "Diagonal entry incorrect" 11909df49d7eSJed Brown /// ); 11919df49d7eSJed Brown /// }); 1192c68be7a2SJeremy L Thompson /// # Ok(()) 1193c68be7a2SJeremy L Thompson /// # } 11949df49d7eSJed Brown /// ``` 11959df49d7eSJed Brown pub fn linear_assemble_diagonal(&self, assembled: &mut Vector) -> crate::Result<i32> { 11969df49d7eSJed Brown self.op_core.linear_assemble_diagonal(assembled) 11979df49d7eSJed Brown } 11989df49d7eSJed Brown 11999df49d7eSJed Brown /// Assemble the diagonal of a square linear Operator 12009df49d7eSJed Brown /// 12019df49d7eSJed Brown /// This sums into a Vector with the diagonal of a linear Operator. 12029df49d7eSJed Brown /// 12039df49d7eSJed Brown /// Note: Currently only non-composite Operators with a single field and 12049df49d7eSJed Brown /// composite Operators with single field sub-operators are supported. 12059df49d7eSJed Brown /// 12069df49d7eSJed Brown /// * `op` - Operator to assemble QFunction 12079df49d7eSJed Brown /// * `assembled` - Vector to store assembled Operator diagonal 12089df49d7eSJed Brown /// 12099df49d7eSJed Brown /// 12109df49d7eSJed Brown /// ``` 12119df49d7eSJed Brown /// # use libceed::prelude::*; 12124d27c890SJeremy L Thompson /// # fn main() -> libceed::Result<()> { 12139df49d7eSJed Brown /// # let ceed = libceed::Ceed::default_init(); 12149df49d7eSJed Brown /// let ne = 4; 12159df49d7eSJed Brown /// let p = 3; 12169df49d7eSJed Brown /// let q = 4; 12179df49d7eSJed Brown /// let ndofs = p * ne - ne + 1; 12189df49d7eSJed Brown /// 12199df49d7eSJed Brown /// // Vectors 1220c68be7a2SJeremy L Thompson /// let x = ceed.vector_from_slice(&[-1., -0.5, 0.0, 0.5, 1.0])?; 1221c68be7a2SJeremy L Thompson /// let mut qdata = ceed.vector(ne * q)?; 12229df49d7eSJed Brown /// qdata.set_value(0.0); 1223c68be7a2SJeremy L Thompson /// let mut u = ceed.vector(ndofs)?; 12249df49d7eSJed Brown /// u.set_value(1.0); 1225c68be7a2SJeremy L Thompson /// let mut v = ceed.vector(ndofs)?; 12269df49d7eSJed Brown /// v.set_value(0.0); 12279df49d7eSJed Brown /// 12289df49d7eSJed Brown /// // Restrictions 12299df49d7eSJed Brown /// let mut indx: Vec<i32> = vec![0; 2 * ne]; 12309df49d7eSJed Brown /// for i in 0..ne { 12319df49d7eSJed Brown /// indx[2 * i + 0] = i as i32; 12329df49d7eSJed Brown /// indx[2 * i + 1] = (i + 1) as i32; 12339df49d7eSJed Brown /// } 1234c68be7a2SJeremy L Thompson /// let rx = ceed.elem_restriction(ne, 2, 1, 1, ne + 1, MemType::Host, &indx)?; 12359df49d7eSJed Brown /// let mut indu: Vec<i32> = vec![0; p * ne]; 12369df49d7eSJed Brown /// for i in 0..ne { 12379df49d7eSJed Brown /// indu[p * i + 0] = (2 * i) as i32; 12389df49d7eSJed Brown /// indu[p * i + 1] = (2 * i + 1) as i32; 12399df49d7eSJed Brown /// indu[p * i + 2] = (2 * i + 2) as i32; 12409df49d7eSJed Brown /// } 1241c68be7a2SJeremy L Thompson /// let ru = ceed.elem_restriction(ne, p, 1, 1, ndofs, MemType::Host, &indu)?; 12429df49d7eSJed Brown /// let strides: [i32; 3] = [1, q as i32, q as i32]; 1243c68be7a2SJeremy L Thompson /// let rq = ceed.strided_elem_restriction(ne, q, 1, q * ne, strides)?; 12449df49d7eSJed Brown /// 12459df49d7eSJed Brown /// // Bases 1246c68be7a2SJeremy L Thompson /// let bx = ceed.basis_tensor_H1_Lagrange(1, 1, 2, q, QuadMode::Gauss)?; 1247c68be7a2SJeremy L Thompson /// let bu = ceed.basis_tensor_H1_Lagrange(1, 1, p, q, QuadMode::Gauss)?; 12489df49d7eSJed Brown /// 12499df49d7eSJed Brown /// // Build quadrature data 1250c68be7a2SJeremy L Thompson /// let qf_build = ceed.q_function_interior_by_name("Mass1DBuild")?; 1251c68be7a2SJeremy L Thompson /// ceed.operator(&qf_build, QFunctionOpt::None, QFunctionOpt::None)? 1252c68be7a2SJeremy L Thompson /// .field("dx", &rx, &bx, VectorOpt::Active)? 1253c68be7a2SJeremy L Thompson /// .field("weights", ElemRestrictionOpt::None, &bx, VectorOpt::None)? 1254356036faSJeremy L Thompson /// .field("qdata", &rq, BasisOpt::None, VectorOpt::Active)? 1255c68be7a2SJeremy L Thompson /// .apply(&x, &mut qdata)?; 12569df49d7eSJed Brown /// 12579df49d7eSJed Brown /// // Mass operator 1258c68be7a2SJeremy L Thompson /// let qf_mass = ceed.q_function_interior_by_name("MassApply")?; 12599df49d7eSJed Brown /// let op_mass = ceed 1260c68be7a2SJeremy L Thompson /// .operator(&qf_mass, QFunctionOpt::None, QFunctionOpt::None)? 1261c68be7a2SJeremy L Thompson /// .field("u", &ru, &bu, VectorOpt::Active)? 1262356036faSJeremy L Thompson /// .field("qdata", &rq, BasisOpt::None, &qdata)? 1263c68be7a2SJeremy L Thompson /// .field("v", &ru, &bu, VectorOpt::Active)?; 12649df49d7eSJed Brown /// 12659df49d7eSJed Brown /// // Diagonal 1266c68be7a2SJeremy L Thompson /// let mut diag = ceed.vector(ndofs)?; 12679df49d7eSJed Brown /// diag.set_value(1.0); 1268c68be7a2SJeremy L Thompson /// op_mass.linear_assemble_add_diagonal(&mut diag)?; 12699df49d7eSJed Brown /// 12709df49d7eSJed Brown /// // Manual diagonal computation 1271c68be7a2SJeremy L Thompson /// let mut true_diag = ceed.vector(ndofs)?; 12729c774eddSJeremy L Thompson /// true_diag.set_value(0.0)?; 12739df49d7eSJed Brown /// for i in 0..ndofs { 12749df49d7eSJed Brown /// u.set_value(0.0); 12759df49d7eSJed Brown /// { 1276e78171edSJeremy L Thompson /// let mut u_array = u.view_mut()?; 12779df49d7eSJed Brown /// u_array[i] = 1.; 12789df49d7eSJed Brown /// } 12799df49d7eSJed Brown /// 1280c68be7a2SJeremy L Thompson /// op_mass.apply(&u, &mut v)?; 12819df49d7eSJed Brown /// 12829df49d7eSJed Brown /// { 12839c774eddSJeremy L Thompson /// let v_array = v.view()?; 1284e78171edSJeremy L Thompson /// let mut true_array = true_diag.view_mut()?; 12859df49d7eSJed Brown /// true_array[i] = v_array[i] + 1.0; 12869df49d7eSJed Brown /// } 12879df49d7eSJed Brown /// } 12889df49d7eSJed Brown /// 12899df49d7eSJed Brown /// // Check 1290e78171edSJeremy L Thompson /// diag.view()? 12919df49d7eSJed Brown /// .iter() 1292e78171edSJeremy L Thompson /// .zip(true_diag.view()?.iter()) 12939df49d7eSJed Brown /// .for_each(|(computed, actual)| { 12949df49d7eSJed Brown /// assert!( 129580a9ef05SNatalie Beams /// (*computed - *actual).abs() < 10.0 * libceed::EPSILON, 12969df49d7eSJed Brown /// "Diagonal entry incorrect" 12979df49d7eSJed Brown /// ); 12989df49d7eSJed Brown /// }); 1299c68be7a2SJeremy L Thompson /// # Ok(()) 1300c68be7a2SJeremy L Thompson /// # } 13019df49d7eSJed Brown /// ``` 13029df49d7eSJed Brown pub fn linear_assemble_add_diagonal(&self, assembled: &mut Vector) -> crate::Result<i32> { 13039df49d7eSJed Brown self.op_core.linear_assemble_add_diagonal(assembled) 13049df49d7eSJed Brown } 13059df49d7eSJed Brown 13069df49d7eSJed Brown /// Assemble the point block diagonal of a square linear Operator 13079df49d7eSJed Brown /// 13089df49d7eSJed Brown /// This overwrites a Vector with the point block diagonal of a linear 13099df49d7eSJed Brown /// Operator. 13109df49d7eSJed Brown /// 13119df49d7eSJed Brown /// Note: Currently only non-composite Operators with a single field and 13129df49d7eSJed Brown /// composite Operators with single field sub-operators are supported. 13139df49d7eSJed Brown /// 13149df49d7eSJed Brown /// * `op` - Operator to assemble QFunction 13159df49d7eSJed Brown /// * `assembled` - Vector to store assembled CeedOperator point block 13169df49d7eSJed Brown /// diagonal, provided in row-major form with an 13179df49d7eSJed Brown /// `ncomp * ncomp` block at each node. The dimensions of 13189df49d7eSJed Brown /// this vector are derived from the active vector for 13199df49d7eSJed Brown /// the CeedOperator. The array has shape 13209df49d7eSJed Brown /// `[nodes, component out, component in]`. 13219df49d7eSJed Brown /// 13229df49d7eSJed Brown /// ``` 13239df49d7eSJed Brown /// # use libceed::prelude::*; 13244d27c890SJeremy L Thompson /// # fn main() -> libceed::Result<()> { 13259df49d7eSJed Brown /// # let ceed = libceed::Ceed::default_init(); 13269df49d7eSJed Brown /// let ne = 4; 13279df49d7eSJed Brown /// let p = 3; 13289df49d7eSJed Brown /// let q = 4; 13299df49d7eSJed Brown /// let ncomp = 2; 13309df49d7eSJed Brown /// let ndofs = p * ne - ne + 1; 13319df49d7eSJed Brown /// 13329df49d7eSJed Brown /// // Vectors 1333c68be7a2SJeremy L Thompson /// let x = ceed.vector_from_slice(&[-1., -0.5, 0.0, 0.5, 1.0])?; 1334c68be7a2SJeremy L Thompson /// let mut qdata = ceed.vector(ne * q)?; 13359df49d7eSJed Brown /// qdata.set_value(0.0); 1336c68be7a2SJeremy L Thompson /// let mut u = ceed.vector(ncomp * ndofs)?; 13379df49d7eSJed Brown /// u.set_value(1.0); 1338c68be7a2SJeremy L Thompson /// let mut v = ceed.vector(ncomp * ndofs)?; 13399df49d7eSJed Brown /// v.set_value(0.0); 13409df49d7eSJed Brown /// 13419df49d7eSJed Brown /// // Restrictions 13429df49d7eSJed Brown /// let mut indx: Vec<i32> = vec![0; 2 * ne]; 13439df49d7eSJed Brown /// for i in 0..ne { 13449df49d7eSJed Brown /// indx[2 * i + 0] = i as i32; 13459df49d7eSJed Brown /// indx[2 * i + 1] = (i + 1) as i32; 13469df49d7eSJed Brown /// } 1347c68be7a2SJeremy L Thompson /// let rx = ceed.elem_restriction(ne, 2, 1, 1, ne + 1, MemType::Host, &indx)?; 13489df49d7eSJed Brown /// let mut indu: Vec<i32> = vec![0; p * ne]; 13499df49d7eSJed Brown /// for i in 0..ne { 13509df49d7eSJed Brown /// indu[p * i + 0] = (2 * i) as i32; 13519df49d7eSJed Brown /// indu[p * i + 1] = (2 * i + 1) as i32; 13529df49d7eSJed Brown /// indu[p * i + 2] = (2 * i + 2) as i32; 13539df49d7eSJed Brown /// } 1354c68be7a2SJeremy L Thompson /// let ru = ceed.elem_restriction(ne, p, ncomp, ndofs, ncomp * ndofs, MemType::Host, &indu)?; 13559df49d7eSJed Brown /// let strides: [i32; 3] = [1, q as i32, q as i32]; 1356c68be7a2SJeremy L Thompson /// let rq = ceed.strided_elem_restriction(ne, q, 1, q * ne, strides)?; 13579df49d7eSJed Brown /// 13589df49d7eSJed Brown /// // Bases 1359c68be7a2SJeremy L Thompson /// let bx = ceed.basis_tensor_H1_Lagrange(1, 1, 2, q, QuadMode::Gauss)?; 1360c68be7a2SJeremy L Thompson /// let bu = ceed.basis_tensor_H1_Lagrange(1, ncomp, p, q, QuadMode::Gauss)?; 13619df49d7eSJed Brown /// 13629df49d7eSJed Brown /// // Build quadrature data 1363c68be7a2SJeremy L Thompson /// let qf_build = ceed.q_function_interior_by_name("Mass1DBuild")?; 1364c68be7a2SJeremy L Thompson /// ceed.operator(&qf_build, QFunctionOpt::None, QFunctionOpt::None)? 1365c68be7a2SJeremy L Thompson /// .field("dx", &rx, &bx, VectorOpt::Active)? 1366c68be7a2SJeremy L Thompson /// .field("weights", ElemRestrictionOpt::None, &bx, VectorOpt::None)? 1367356036faSJeremy L Thompson /// .field("qdata", &rq, BasisOpt::None, VectorOpt::Active)? 1368c68be7a2SJeremy L Thompson /// .apply(&x, &mut qdata)?; 13699df49d7eSJed Brown /// 13709df49d7eSJed Brown /// // Mass operator 13719df49d7eSJed Brown /// let mut mass_2_comp = |[u, qdata, ..]: QFunctionInputs, [v, ..]: QFunctionOutputs| { 13729df49d7eSJed Brown /// // Number of quadrature points 13739df49d7eSJed Brown /// let q = qdata.len(); 13749df49d7eSJed Brown /// 13759df49d7eSJed Brown /// // Iterate over quadrature points 13769df49d7eSJed Brown /// for i in 0..q { 13779df49d7eSJed Brown /// v[i + 0 * q] = u[i + 1 * q] * qdata[i]; 13789df49d7eSJed Brown /// v[i + 1 * q] = u[i + 0 * q] * qdata[i]; 13799df49d7eSJed Brown /// } 13809df49d7eSJed Brown /// 13819df49d7eSJed Brown /// // Return clean error code 13829df49d7eSJed Brown /// 0 13839df49d7eSJed Brown /// }; 13849df49d7eSJed Brown /// 13859df49d7eSJed Brown /// let qf_mass = ceed 1386c68be7a2SJeremy L Thompson /// .q_function_interior(1, Box::new(mass_2_comp))? 1387c68be7a2SJeremy L Thompson /// .input("u", 2, EvalMode::Interp)? 1388c68be7a2SJeremy L Thompson /// .input("qdata", 1, EvalMode::None)? 1389c68be7a2SJeremy L Thompson /// .output("v", 2, EvalMode::Interp)?; 13909df49d7eSJed Brown /// 13919df49d7eSJed Brown /// let op_mass = ceed 1392c68be7a2SJeremy L Thompson /// .operator(&qf_mass, QFunctionOpt::None, QFunctionOpt::None)? 1393c68be7a2SJeremy L Thompson /// .field("u", &ru, &bu, VectorOpt::Active)? 1394356036faSJeremy L Thompson /// .field("qdata", &rq, BasisOpt::None, &qdata)? 1395c68be7a2SJeremy L Thompson /// .field("v", &ru, &bu, VectorOpt::Active)?; 13969df49d7eSJed Brown /// 13979df49d7eSJed Brown /// // Diagonal 1398c68be7a2SJeremy L Thompson /// let mut diag = ceed.vector(ncomp * ncomp * ndofs)?; 13999df49d7eSJed Brown /// diag.set_value(0.0); 1400c68be7a2SJeremy L Thompson /// op_mass.linear_assemble_point_block_diagonal(&mut diag)?; 14019df49d7eSJed Brown /// 14029df49d7eSJed Brown /// // Manual diagonal computation 1403c68be7a2SJeremy L Thompson /// let mut true_diag = ceed.vector(ncomp * ncomp * ndofs)?; 14049c774eddSJeremy L Thompson /// true_diag.set_value(0.0)?; 14059df49d7eSJed Brown /// for i in 0..ndofs { 14069df49d7eSJed Brown /// for j in 0..ncomp { 14079df49d7eSJed Brown /// u.set_value(0.0); 14089df49d7eSJed Brown /// { 1409e78171edSJeremy L Thompson /// let mut u_array = u.view_mut()?; 14109df49d7eSJed Brown /// u_array[i + j * ndofs] = 1.; 14119df49d7eSJed Brown /// } 14129df49d7eSJed Brown /// 1413c68be7a2SJeremy L Thompson /// op_mass.apply(&u, &mut v)?; 14149df49d7eSJed Brown /// 14159df49d7eSJed Brown /// { 14169c774eddSJeremy L Thompson /// let v_array = v.view()?; 1417e78171edSJeremy L Thompson /// let mut true_array = true_diag.view_mut()?; 14189df49d7eSJed Brown /// for k in 0..ncomp { 14199df49d7eSJed Brown /// true_array[i * ncomp * ncomp + k * ncomp + j] = v_array[i + k * ndofs]; 14209df49d7eSJed Brown /// } 14219df49d7eSJed Brown /// } 14229df49d7eSJed Brown /// } 14239df49d7eSJed Brown /// } 14249df49d7eSJed Brown /// 14259df49d7eSJed Brown /// // Check 1426e78171edSJeremy L Thompson /// diag.view()? 14279df49d7eSJed Brown /// .iter() 1428e78171edSJeremy L Thompson /// .zip(true_diag.view()?.iter()) 14299df49d7eSJed Brown /// .for_each(|(computed, actual)| { 14309df49d7eSJed Brown /// assert!( 143180a9ef05SNatalie Beams /// (*computed - *actual).abs() < 10.0 * libceed::EPSILON, 14329df49d7eSJed Brown /// "Diagonal entry incorrect" 14339df49d7eSJed Brown /// ); 14349df49d7eSJed Brown /// }); 1435c68be7a2SJeremy L Thompson /// # Ok(()) 1436c68be7a2SJeremy L Thompson /// # } 14379df49d7eSJed Brown /// ``` 14389df49d7eSJed Brown pub fn linear_assemble_point_block_diagonal( 14399df49d7eSJed Brown &self, 14409df49d7eSJed Brown assembled: &mut Vector, 14419df49d7eSJed Brown ) -> crate::Result<i32> { 14429df49d7eSJed Brown self.op_core.linear_assemble_point_block_diagonal(assembled) 14439df49d7eSJed Brown } 14449df49d7eSJed Brown 14459df49d7eSJed Brown /// Assemble the point block diagonal of a square linear Operator 14469df49d7eSJed Brown /// 14479df49d7eSJed Brown /// This sums into a Vector with the point block diagonal of a linear 14489df49d7eSJed Brown /// Operator. 14499df49d7eSJed Brown /// 14509df49d7eSJed Brown /// Note: Currently only non-composite Operators with a single field and 14519df49d7eSJed Brown /// composite Operators with single field sub-operators are supported. 14529df49d7eSJed Brown /// 14539df49d7eSJed Brown /// * `op` - Operator to assemble QFunction 14549df49d7eSJed Brown /// * `assembled` - Vector to store assembled CeedOperator point block 14559df49d7eSJed Brown /// diagonal, provided in row-major form with an 14569df49d7eSJed Brown /// `ncomp * ncomp` block at each node. The dimensions of 14579df49d7eSJed Brown /// this vector are derived from the active vector for 14589df49d7eSJed Brown /// the CeedOperator. The array has shape 14599df49d7eSJed Brown /// `[nodes, component out, component in]`. 14609df49d7eSJed Brown /// 14619df49d7eSJed Brown /// ``` 14629df49d7eSJed Brown /// # use libceed::prelude::*; 14634d27c890SJeremy L Thompson /// # fn main() -> libceed::Result<()> { 14649df49d7eSJed Brown /// # let ceed = libceed::Ceed::default_init(); 14659df49d7eSJed Brown /// let ne = 4; 14669df49d7eSJed Brown /// let p = 3; 14679df49d7eSJed Brown /// let q = 4; 14689df49d7eSJed Brown /// let ncomp = 2; 14699df49d7eSJed Brown /// let ndofs = p * ne - ne + 1; 14709df49d7eSJed Brown /// 14719df49d7eSJed Brown /// // Vectors 1472c68be7a2SJeremy L Thompson /// let x = ceed.vector_from_slice(&[-1., -0.5, 0.0, 0.5, 1.0])?; 1473c68be7a2SJeremy L Thompson /// let mut qdata = ceed.vector(ne * q)?; 14749df49d7eSJed Brown /// qdata.set_value(0.0); 1475c68be7a2SJeremy L Thompson /// let mut u = ceed.vector(ncomp * ndofs)?; 14769df49d7eSJed Brown /// u.set_value(1.0); 1477c68be7a2SJeremy L Thompson /// let mut v = ceed.vector(ncomp * ndofs)?; 14789df49d7eSJed Brown /// v.set_value(0.0); 14799df49d7eSJed Brown /// 14809df49d7eSJed Brown /// // Restrictions 14819df49d7eSJed Brown /// let mut indx: Vec<i32> = vec![0; 2 * ne]; 14829df49d7eSJed Brown /// for i in 0..ne { 14839df49d7eSJed Brown /// indx[2 * i + 0] = i as i32; 14849df49d7eSJed Brown /// indx[2 * i + 1] = (i + 1) as i32; 14859df49d7eSJed Brown /// } 1486c68be7a2SJeremy L Thompson /// let rx = ceed.elem_restriction(ne, 2, 1, 1, ne + 1, MemType::Host, &indx)?; 14879df49d7eSJed Brown /// let mut indu: Vec<i32> = vec![0; p * ne]; 14889df49d7eSJed Brown /// for i in 0..ne { 14899df49d7eSJed Brown /// indu[p * i + 0] = (2 * i) as i32; 14909df49d7eSJed Brown /// indu[p * i + 1] = (2 * i + 1) as i32; 14919df49d7eSJed Brown /// indu[p * i + 2] = (2 * i + 2) as i32; 14929df49d7eSJed Brown /// } 1493c68be7a2SJeremy L Thompson /// let ru = ceed.elem_restriction(ne, p, ncomp, ndofs, ncomp * ndofs, MemType::Host, &indu)?; 14949df49d7eSJed Brown /// let strides: [i32; 3] = [1, q as i32, q as i32]; 1495c68be7a2SJeremy L Thompson /// let rq = ceed.strided_elem_restriction(ne, q, 1, q * ne, strides)?; 14969df49d7eSJed Brown /// 14979df49d7eSJed Brown /// // Bases 1498c68be7a2SJeremy L Thompson /// let bx = ceed.basis_tensor_H1_Lagrange(1, 1, 2, q, QuadMode::Gauss)?; 1499c68be7a2SJeremy L Thompson /// let bu = ceed.basis_tensor_H1_Lagrange(1, ncomp, p, q, QuadMode::Gauss)?; 15009df49d7eSJed Brown /// 15019df49d7eSJed Brown /// // Build quadrature data 1502c68be7a2SJeremy L Thompson /// let qf_build = ceed.q_function_interior_by_name("Mass1DBuild")?; 1503c68be7a2SJeremy L Thompson /// ceed.operator(&qf_build, QFunctionOpt::None, QFunctionOpt::None)? 1504c68be7a2SJeremy L Thompson /// .field("dx", &rx, &bx, VectorOpt::Active)? 1505c68be7a2SJeremy L Thompson /// .field("weights", ElemRestrictionOpt::None, &bx, VectorOpt::None)? 1506356036faSJeremy L Thompson /// .field("qdata", &rq, BasisOpt::None, VectorOpt::Active)? 1507c68be7a2SJeremy L Thompson /// .apply(&x, &mut qdata)?; 15089df49d7eSJed Brown /// 15099df49d7eSJed Brown /// // Mass operator 15109df49d7eSJed Brown /// let mut mass_2_comp = |[u, qdata, ..]: QFunctionInputs, [v, ..]: QFunctionOutputs| { 15119df49d7eSJed Brown /// // Number of quadrature points 15129df49d7eSJed Brown /// let q = qdata.len(); 15139df49d7eSJed Brown /// 15149df49d7eSJed Brown /// // Iterate over quadrature points 15159df49d7eSJed Brown /// for i in 0..q { 15169df49d7eSJed Brown /// v[i + 0 * q] = u[i + 1 * q] * qdata[i]; 15179df49d7eSJed Brown /// v[i + 1 * q] = u[i + 0 * q] * qdata[i]; 15189df49d7eSJed Brown /// } 15199df49d7eSJed Brown /// 15209df49d7eSJed Brown /// // Return clean error code 15219df49d7eSJed Brown /// 0 15229df49d7eSJed Brown /// }; 15239df49d7eSJed Brown /// 15249df49d7eSJed Brown /// let qf_mass = ceed 1525c68be7a2SJeremy L Thompson /// .q_function_interior(1, Box::new(mass_2_comp))? 1526c68be7a2SJeremy L Thompson /// .input("u", 2, EvalMode::Interp)? 1527c68be7a2SJeremy L Thompson /// .input("qdata", 1, EvalMode::None)? 1528c68be7a2SJeremy L Thompson /// .output("v", 2, EvalMode::Interp)?; 15299df49d7eSJed Brown /// 15309df49d7eSJed Brown /// let op_mass = ceed 1531c68be7a2SJeremy L Thompson /// .operator(&qf_mass, QFunctionOpt::None, QFunctionOpt::None)? 1532c68be7a2SJeremy L Thompson /// .field("u", &ru, &bu, VectorOpt::Active)? 1533356036faSJeremy L Thompson /// .field("qdata", &rq, BasisOpt::None, &qdata)? 1534c68be7a2SJeremy L Thompson /// .field("v", &ru, &bu, VectorOpt::Active)?; 15359df49d7eSJed Brown /// 15369df49d7eSJed Brown /// // Diagonal 1537c68be7a2SJeremy L Thompson /// let mut diag = ceed.vector(ncomp * ncomp * ndofs)?; 15389df49d7eSJed Brown /// diag.set_value(1.0); 1539c68be7a2SJeremy L Thompson /// op_mass.linear_assemble_add_point_block_diagonal(&mut diag)?; 15409df49d7eSJed Brown /// 15419df49d7eSJed Brown /// // Manual diagonal computation 1542c68be7a2SJeremy L Thompson /// let mut true_diag = ceed.vector(ncomp * ncomp * ndofs)?; 15439c774eddSJeremy L Thompson /// true_diag.set_value(0.0)?; 15449df49d7eSJed Brown /// for i in 0..ndofs { 15459df49d7eSJed Brown /// for j in 0..ncomp { 15469df49d7eSJed Brown /// u.set_value(0.0); 15479df49d7eSJed Brown /// { 1548e78171edSJeremy L Thompson /// let mut u_array = u.view_mut()?; 15499df49d7eSJed Brown /// u_array[i + j * ndofs] = 1.; 15509df49d7eSJed Brown /// } 15519df49d7eSJed Brown /// 1552c68be7a2SJeremy L Thompson /// op_mass.apply(&u, &mut v)?; 15539df49d7eSJed Brown /// 15549df49d7eSJed Brown /// { 15559c774eddSJeremy L Thompson /// let v_array = v.view()?; 1556e78171edSJeremy L Thompson /// let mut true_array = true_diag.view_mut()?; 15579df49d7eSJed Brown /// for k in 0..ncomp { 15589df49d7eSJed Brown /// true_array[i * ncomp * ncomp + k * ncomp + j] = v_array[i + k * ndofs]; 15599df49d7eSJed Brown /// } 15609df49d7eSJed Brown /// } 15619df49d7eSJed Brown /// } 15629df49d7eSJed Brown /// } 15639df49d7eSJed Brown /// 15649df49d7eSJed Brown /// // Check 1565e78171edSJeremy L Thompson /// diag.view()? 15669df49d7eSJed Brown /// .iter() 1567e78171edSJeremy L Thompson /// .zip(true_diag.view()?.iter()) 15689df49d7eSJed Brown /// .for_each(|(computed, actual)| { 15699df49d7eSJed Brown /// assert!( 157080a9ef05SNatalie Beams /// (*computed - 1.0 - *actual).abs() < 10.0 * libceed::EPSILON, 15719df49d7eSJed Brown /// "Diagonal entry incorrect" 15729df49d7eSJed Brown /// ); 15739df49d7eSJed Brown /// }); 1574c68be7a2SJeremy L Thompson /// # Ok(()) 1575c68be7a2SJeremy L Thompson /// # } 15769df49d7eSJed Brown /// ``` 15779df49d7eSJed Brown pub fn linear_assemble_add_point_block_diagonal( 15789df49d7eSJed Brown &self, 15799df49d7eSJed Brown assembled: &mut Vector, 15809df49d7eSJed Brown ) -> crate::Result<i32> { 15819df49d7eSJed Brown self.op_core 15829df49d7eSJed Brown .linear_assemble_add_point_block_diagonal(assembled) 15839df49d7eSJed Brown } 15849df49d7eSJed Brown 15859df49d7eSJed Brown /// Create a multigrid coarse Operator and level transfer Operators for a 15869df49d7eSJed Brown /// given Operator, creating the prolongation basis from the fine and 15879df49d7eSJed Brown /// coarse grid interpolation 15889df49d7eSJed Brown /// 15899df49d7eSJed Brown /// * `p_mult_fine` - Lvector multiplicity in parallel gather/scatter 15909df49d7eSJed Brown /// * `rstr_coarse` - Coarse grid restriction 15919df49d7eSJed Brown /// * `basis_coarse` - Coarse grid active vector basis 15929df49d7eSJed Brown /// 15939df49d7eSJed Brown /// ``` 15949df49d7eSJed Brown /// # use libceed::prelude::*; 15954d27c890SJeremy L Thompson /// # fn main() -> libceed::Result<()> { 15969df49d7eSJed Brown /// # let ceed = libceed::Ceed::default_init(); 15979df49d7eSJed Brown /// let ne = 15; 15989df49d7eSJed Brown /// let p_coarse = 3; 15999df49d7eSJed Brown /// let p_fine = 5; 16009df49d7eSJed Brown /// let q = 6; 16019df49d7eSJed Brown /// let ndofs_coarse = p_coarse * ne - ne + 1; 16029df49d7eSJed Brown /// let ndofs_fine = p_fine * ne - ne + 1; 16039df49d7eSJed Brown /// 16049df49d7eSJed Brown /// // Vectors 16059df49d7eSJed Brown /// let x_array = (0..ne + 1) 160680a9ef05SNatalie Beams /// .map(|i| 2.0 * i as Scalar / ne as Scalar - 1.0) 160780a9ef05SNatalie Beams /// .collect::<Vec<Scalar>>(); 1608c68be7a2SJeremy L Thompson /// let x = ceed.vector_from_slice(&x_array)?; 1609c68be7a2SJeremy L Thompson /// let mut qdata = ceed.vector(ne * q)?; 16109df49d7eSJed Brown /// qdata.set_value(0.0); 1611c68be7a2SJeremy L Thompson /// let mut u_coarse = ceed.vector(ndofs_coarse)?; 16129df49d7eSJed Brown /// u_coarse.set_value(1.0); 1613c68be7a2SJeremy L Thompson /// let mut u_fine = ceed.vector(ndofs_fine)?; 16149df49d7eSJed Brown /// u_fine.set_value(1.0); 1615c68be7a2SJeremy L Thompson /// let mut v_coarse = ceed.vector(ndofs_coarse)?; 16169df49d7eSJed Brown /// v_coarse.set_value(0.0); 1617c68be7a2SJeremy L Thompson /// let mut v_fine = ceed.vector(ndofs_fine)?; 16189df49d7eSJed Brown /// v_fine.set_value(0.0); 1619c68be7a2SJeremy L Thompson /// let mut multiplicity = ceed.vector(ndofs_fine)?; 16209df49d7eSJed Brown /// multiplicity.set_value(1.0); 16219df49d7eSJed Brown /// 16229df49d7eSJed Brown /// // Restrictions 16239df49d7eSJed Brown /// let strides: [i32; 3] = [1, q as i32, q as i32]; 1624c68be7a2SJeremy L Thompson /// let rq = ceed.strided_elem_restriction(ne, q, 1, q * ne, strides)?; 16259df49d7eSJed Brown /// 16269df49d7eSJed Brown /// let mut indx: Vec<i32> = vec![0; 2 * ne]; 16279df49d7eSJed Brown /// for i in 0..ne { 16289df49d7eSJed Brown /// indx[2 * i + 0] = i as i32; 16299df49d7eSJed Brown /// indx[2 * i + 1] = (i + 1) as i32; 16309df49d7eSJed Brown /// } 1631c68be7a2SJeremy L Thompson /// let rx = ceed.elem_restriction(ne, 2, 1, 1, ne + 1, MemType::Host, &indx)?; 16329df49d7eSJed Brown /// 16339df49d7eSJed Brown /// let mut indu_coarse: Vec<i32> = vec![0; p_coarse * ne]; 16349df49d7eSJed Brown /// for i in 0..ne { 16359df49d7eSJed Brown /// for j in 0..p_coarse { 16369df49d7eSJed Brown /// indu_coarse[p_coarse * i + j] = (i + j) as i32; 16379df49d7eSJed Brown /// } 16389df49d7eSJed Brown /// } 1639c68be7a2SJeremy L Thompson /// let ru_coarse = ceed.elem_restriction( 16409df49d7eSJed Brown /// ne, 16419df49d7eSJed Brown /// p_coarse, 16429df49d7eSJed Brown /// 1, 16439df49d7eSJed Brown /// 1, 16449df49d7eSJed Brown /// ndofs_coarse, 16459df49d7eSJed Brown /// MemType::Host, 16469df49d7eSJed Brown /// &indu_coarse, 1647c68be7a2SJeremy L Thompson /// )?; 16489df49d7eSJed Brown /// 16499df49d7eSJed Brown /// let mut indu_fine: Vec<i32> = vec![0; p_fine * ne]; 16509df49d7eSJed Brown /// for i in 0..ne { 16519df49d7eSJed Brown /// for j in 0..p_fine { 16529df49d7eSJed Brown /// indu_fine[p_fine * i + j] = (i + j) as i32; 16539df49d7eSJed Brown /// } 16549df49d7eSJed Brown /// } 1655c68be7a2SJeremy L Thompson /// let ru_fine = ceed.elem_restriction(ne, p_fine, 1, 1, ndofs_fine, MemType::Host, &indu_fine)?; 16569df49d7eSJed Brown /// 16579df49d7eSJed Brown /// // Bases 1658c68be7a2SJeremy L Thompson /// let bx = ceed.basis_tensor_H1_Lagrange(1, 1, 2, q, QuadMode::Gauss)?; 1659c68be7a2SJeremy L Thompson /// let bu_coarse = ceed.basis_tensor_H1_Lagrange(1, 1, p_coarse, q, QuadMode::Gauss)?; 1660c68be7a2SJeremy L Thompson /// let bu_fine = ceed.basis_tensor_H1_Lagrange(1, 1, p_fine, q, QuadMode::Gauss)?; 16619df49d7eSJed Brown /// 16629df49d7eSJed Brown /// // Build quadrature data 1663c68be7a2SJeremy L Thompson /// let qf_build = ceed.q_function_interior_by_name("Mass1DBuild")?; 1664c68be7a2SJeremy L Thompson /// ceed.operator(&qf_build, QFunctionOpt::None, QFunctionOpt::None)? 1665c68be7a2SJeremy L Thompson /// .field("dx", &rx, &bx, VectorOpt::Active)? 1666c68be7a2SJeremy L Thompson /// .field("weights", ElemRestrictionOpt::None, &bx, VectorOpt::None)? 1667356036faSJeremy L Thompson /// .field("qdata", &rq, BasisOpt::None, VectorOpt::Active)? 1668c68be7a2SJeremy L Thompson /// .apply(&x, &mut qdata)?; 16699df49d7eSJed Brown /// 16709df49d7eSJed Brown /// // Mass operator 1671c68be7a2SJeremy L Thompson /// let qf_mass = ceed.q_function_interior_by_name("MassApply")?; 16729df49d7eSJed Brown /// let op_mass_fine = ceed 1673c68be7a2SJeremy L Thompson /// .operator(&qf_mass, QFunctionOpt::None, QFunctionOpt::None)? 1674c68be7a2SJeremy L Thompson /// .field("u", &ru_fine, &bu_fine, VectorOpt::Active)? 1675356036faSJeremy L Thompson /// .field("qdata", &rq, BasisOpt::None, &qdata)? 1676c68be7a2SJeremy L Thompson /// .field("v", &ru_fine, &bu_fine, VectorOpt::Active)?; 16779df49d7eSJed Brown /// 16789df49d7eSJed Brown /// // Multigrid setup 1679c68be7a2SJeremy L Thompson /// let (op_mass_coarse, op_prolong, op_restrict) = 1680c68be7a2SJeremy L Thompson /// op_mass_fine.create_multigrid_level(&multiplicity, &ru_coarse, &bu_coarse)?; 16819df49d7eSJed Brown /// 16829df49d7eSJed Brown /// // Coarse problem 16839df49d7eSJed Brown /// u_coarse.set_value(1.0); 1684c68be7a2SJeremy L Thompson /// op_mass_coarse.apply(&u_coarse, &mut v_coarse)?; 16859df49d7eSJed Brown /// 16869df49d7eSJed Brown /// // Check 1687e78171edSJeremy L Thompson /// let sum: Scalar = v_coarse.view()?.iter().sum(); 16889df49d7eSJed Brown /// assert!( 168980a9ef05SNatalie Beams /// (sum - 2.0).abs() < 50.0 * libceed::EPSILON, 16909df49d7eSJed Brown /// "Incorrect interval length computed" 16919df49d7eSJed Brown /// ); 16929df49d7eSJed Brown /// 16939df49d7eSJed Brown /// // Prolong 1694c68be7a2SJeremy L Thompson /// op_prolong.apply(&u_coarse, &mut u_fine)?; 16959df49d7eSJed Brown /// 16969df49d7eSJed Brown /// // Fine problem 1697c68be7a2SJeremy L Thompson /// op_mass_fine.apply(&u_fine, &mut v_fine)?; 16989df49d7eSJed Brown /// 16999df49d7eSJed Brown /// // Check 1700e78171edSJeremy L Thompson /// let sum: Scalar = v_fine.view()?.iter().sum(); 17019df49d7eSJed Brown /// assert!( 1702392d940cSJeremy L Thompson /// (sum - 2.0).abs() < 200.0 * libceed::EPSILON, 17039df49d7eSJed Brown /// "Incorrect interval length computed" 17049df49d7eSJed Brown /// ); 17059df49d7eSJed Brown /// 17069df49d7eSJed Brown /// // Restrict 1707c68be7a2SJeremy L Thompson /// op_restrict.apply(&v_fine, &mut v_coarse)?; 17089df49d7eSJed Brown /// 17099df49d7eSJed Brown /// // Check 1710e78171edSJeremy L Thompson /// let sum: Scalar = v_coarse.view()?.iter().sum(); 17119df49d7eSJed Brown /// assert!( 1712392d940cSJeremy L Thompson /// (sum - 2.0).abs() < 200.0 * libceed::EPSILON, 17139df49d7eSJed Brown /// "Incorrect interval length computed" 17149df49d7eSJed Brown /// ); 1715c68be7a2SJeremy L Thompson /// # Ok(()) 1716c68be7a2SJeremy L Thompson /// # } 17179df49d7eSJed Brown /// ``` 1718594ef120SJeremy L Thompson pub fn create_multigrid_level<'b>( 17199df49d7eSJed Brown &self, 17209df49d7eSJed Brown p_mult_fine: &Vector, 17219df49d7eSJed Brown rstr_coarse: &ElemRestriction, 17229df49d7eSJed Brown basis_coarse: &Basis, 1723594ef120SJeremy L Thompson ) -> crate::Result<(Operator<'b>, Operator<'b>, Operator<'b>)> { 17249df49d7eSJed Brown let mut ptr_coarse = std::ptr::null_mut(); 17259df49d7eSJed Brown let mut ptr_prolong = std::ptr::null_mut(); 17269df49d7eSJed Brown let mut ptr_restrict = std::ptr::null_mut(); 17279df49d7eSJed Brown let ierr = unsafe { 17289df49d7eSJed Brown bind_ceed::CeedOperatorMultigridLevelCreate( 17299df49d7eSJed Brown self.op_core.ptr, 17309df49d7eSJed Brown p_mult_fine.ptr, 17319df49d7eSJed Brown rstr_coarse.ptr, 17329df49d7eSJed Brown basis_coarse.ptr, 17339df49d7eSJed Brown &mut ptr_coarse, 17349df49d7eSJed Brown &mut ptr_prolong, 17359df49d7eSJed Brown &mut ptr_restrict, 17369df49d7eSJed Brown ) 17379df49d7eSJed Brown }; 17381142270cSJeremy L Thompson self.op_core.check_error(ierr)?; 17391142270cSJeremy L Thompson let op_coarse = Operator::from_raw(ptr_coarse)?; 17401142270cSJeremy L Thompson let op_prolong = Operator::from_raw(ptr_prolong)?; 17411142270cSJeremy L Thompson let op_restrict = Operator::from_raw(ptr_restrict)?; 17429df49d7eSJed Brown Ok((op_coarse, op_prolong, op_restrict)) 17439df49d7eSJed Brown } 17449df49d7eSJed Brown 17459df49d7eSJed Brown /// Create a multigrid coarse Operator and level transfer Operators for a 17469df49d7eSJed Brown /// given Operator with a tensor basis for the active basis 17479df49d7eSJed Brown /// 17489df49d7eSJed Brown /// * `p_mult_fine` - Lvector multiplicity in parallel gather/scatter 17499df49d7eSJed Brown /// * `rstr_coarse` - Coarse grid restriction 17509df49d7eSJed Brown /// * `basis_coarse` - Coarse grid active vector basis 17519df49d7eSJed Brown /// * `interp_c_to_f` - Matrix for coarse to fine 17529df49d7eSJed Brown /// 17539df49d7eSJed Brown /// ``` 17549df49d7eSJed Brown /// # use libceed::prelude::*; 17554d27c890SJeremy L Thompson /// # fn main() -> libceed::Result<()> { 17569df49d7eSJed Brown /// # let ceed = libceed::Ceed::default_init(); 17579df49d7eSJed Brown /// let ne = 15; 17589df49d7eSJed Brown /// let p_coarse = 3; 17599df49d7eSJed Brown /// let p_fine = 5; 17609df49d7eSJed Brown /// let q = 6; 17619df49d7eSJed Brown /// let ndofs_coarse = p_coarse * ne - ne + 1; 17629df49d7eSJed Brown /// let ndofs_fine = p_fine * ne - ne + 1; 17639df49d7eSJed Brown /// 17649df49d7eSJed Brown /// // Vectors 17659df49d7eSJed Brown /// let x_array = (0..ne + 1) 176680a9ef05SNatalie Beams /// .map(|i| 2.0 * i as Scalar / ne as Scalar - 1.0) 176780a9ef05SNatalie Beams /// .collect::<Vec<Scalar>>(); 1768c68be7a2SJeremy L Thompson /// let x = ceed.vector_from_slice(&x_array)?; 1769c68be7a2SJeremy L Thompson /// let mut qdata = ceed.vector(ne * q)?; 17709df49d7eSJed Brown /// qdata.set_value(0.0); 1771c68be7a2SJeremy L Thompson /// let mut u_coarse = ceed.vector(ndofs_coarse)?; 17729df49d7eSJed Brown /// u_coarse.set_value(1.0); 1773c68be7a2SJeremy L Thompson /// let mut u_fine = ceed.vector(ndofs_fine)?; 17749df49d7eSJed Brown /// u_fine.set_value(1.0); 1775c68be7a2SJeremy L Thompson /// let mut v_coarse = ceed.vector(ndofs_coarse)?; 17769df49d7eSJed Brown /// v_coarse.set_value(0.0); 1777c68be7a2SJeremy L Thompson /// let mut v_fine = ceed.vector(ndofs_fine)?; 17789df49d7eSJed Brown /// v_fine.set_value(0.0); 1779c68be7a2SJeremy L Thompson /// let mut multiplicity = ceed.vector(ndofs_fine)?; 17809df49d7eSJed Brown /// multiplicity.set_value(1.0); 17819df49d7eSJed Brown /// 17829df49d7eSJed Brown /// // Restrictions 17839df49d7eSJed Brown /// let strides: [i32; 3] = [1, q as i32, q as i32]; 1784c68be7a2SJeremy L Thompson /// let rq = ceed.strided_elem_restriction(ne, q, 1, q * ne, strides)?; 17859df49d7eSJed Brown /// 17869df49d7eSJed Brown /// let mut indx: Vec<i32> = vec![0; 2 * ne]; 17879df49d7eSJed Brown /// for i in 0..ne { 17889df49d7eSJed Brown /// indx[2 * i + 0] = i as i32; 17899df49d7eSJed Brown /// indx[2 * i + 1] = (i + 1) as i32; 17909df49d7eSJed Brown /// } 1791c68be7a2SJeremy L Thompson /// let rx = ceed.elem_restriction(ne, 2, 1, 1, ne + 1, MemType::Host, &indx)?; 17929df49d7eSJed Brown /// 17939df49d7eSJed Brown /// let mut indu_coarse: Vec<i32> = vec![0; p_coarse * ne]; 17949df49d7eSJed Brown /// for i in 0..ne { 17959df49d7eSJed Brown /// for j in 0..p_coarse { 17969df49d7eSJed Brown /// indu_coarse[p_coarse * i + j] = (i + j) as i32; 17979df49d7eSJed Brown /// } 17989df49d7eSJed Brown /// } 1799c68be7a2SJeremy L Thompson /// let ru_coarse = ceed.elem_restriction( 18009df49d7eSJed Brown /// ne, 18019df49d7eSJed Brown /// p_coarse, 18029df49d7eSJed Brown /// 1, 18039df49d7eSJed Brown /// 1, 18049df49d7eSJed Brown /// ndofs_coarse, 18059df49d7eSJed Brown /// MemType::Host, 18069df49d7eSJed Brown /// &indu_coarse, 1807c68be7a2SJeremy L Thompson /// )?; 18089df49d7eSJed Brown /// 18099df49d7eSJed Brown /// let mut indu_fine: Vec<i32> = vec![0; p_fine * ne]; 18109df49d7eSJed Brown /// for i in 0..ne { 18119df49d7eSJed Brown /// for j in 0..p_fine { 18129df49d7eSJed Brown /// indu_fine[p_fine * i + j] = (i + j) as i32; 18139df49d7eSJed Brown /// } 18149df49d7eSJed Brown /// } 1815c68be7a2SJeremy L Thompson /// let ru_fine = ceed.elem_restriction(ne, p_fine, 1, 1, ndofs_fine, MemType::Host, &indu_fine)?; 18169df49d7eSJed Brown /// 18179df49d7eSJed Brown /// // Bases 1818c68be7a2SJeremy L Thompson /// let bx = ceed.basis_tensor_H1_Lagrange(1, 1, 2, q, QuadMode::Gauss)?; 1819c68be7a2SJeremy L Thompson /// let bu_coarse = ceed.basis_tensor_H1_Lagrange(1, 1, p_coarse, q, QuadMode::Gauss)?; 1820c68be7a2SJeremy L Thompson /// let bu_fine = ceed.basis_tensor_H1_Lagrange(1, 1, p_fine, q, QuadMode::Gauss)?; 18219df49d7eSJed Brown /// 18229df49d7eSJed Brown /// // Build quadrature data 1823c68be7a2SJeremy L Thompson /// let qf_build = ceed.q_function_interior_by_name("Mass1DBuild")?; 1824c68be7a2SJeremy L Thompson /// ceed.operator(&qf_build, QFunctionOpt::None, QFunctionOpt::None)? 1825c68be7a2SJeremy L Thompson /// .field("dx", &rx, &bx, VectorOpt::Active)? 1826c68be7a2SJeremy L Thompson /// .field("weights", ElemRestrictionOpt::None, &bx, VectorOpt::None)? 1827356036faSJeremy L Thompson /// .field("qdata", &rq, BasisOpt::None, VectorOpt::Active)? 1828c68be7a2SJeremy L Thompson /// .apply(&x, &mut qdata)?; 18299df49d7eSJed Brown /// 18309df49d7eSJed Brown /// // Mass operator 1831c68be7a2SJeremy L Thompson /// let qf_mass = ceed.q_function_interior_by_name("MassApply")?; 18329df49d7eSJed Brown /// let op_mass_fine = ceed 1833c68be7a2SJeremy L Thompson /// .operator(&qf_mass, QFunctionOpt::None, QFunctionOpt::None)? 1834c68be7a2SJeremy L Thompson /// .field("u", &ru_fine, &bu_fine, VectorOpt::Active)? 1835356036faSJeremy L Thompson /// .field("qdata", &rq, BasisOpt::None, &qdata)? 1836c68be7a2SJeremy L Thompson /// .field("v", &ru_fine, &bu_fine, VectorOpt::Active)?; 18379df49d7eSJed Brown /// 18389df49d7eSJed Brown /// // Multigrid setup 183980a9ef05SNatalie Beams /// let mut interp_c_to_f: Vec<Scalar> = vec![0.; p_coarse * p_fine]; 18409df49d7eSJed Brown /// { 1841c68be7a2SJeremy L Thompson /// let mut coarse = ceed.vector(p_coarse)?; 1842c68be7a2SJeremy L Thompson /// let mut fine = ceed.vector(p_fine)?; 1843c68be7a2SJeremy L Thompson /// let basis_c_to_f = 1844c68be7a2SJeremy L Thompson /// ceed.basis_tensor_H1_Lagrange(1, 1, p_coarse, p_fine, QuadMode::GaussLobatto)?; 18459df49d7eSJed Brown /// for i in 0..p_coarse { 18469df49d7eSJed Brown /// coarse.set_value(0.0); 18479df49d7eSJed Brown /// { 1848e78171edSJeremy L Thompson /// let mut array = coarse.view_mut()?; 18499df49d7eSJed Brown /// array[i] = 1.; 18509df49d7eSJed Brown /// } 1851c68be7a2SJeremy L Thompson /// basis_c_to_f.apply( 18529df49d7eSJed Brown /// 1, 18539df49d7eSJed Brown /// TransposeMode::NoTranspose, 18549df49d7eSJed Brown /// EvalMode::Interp, 18559df49d7eSJed Brown /// &coarse, 18569df49d7eSJed Brown /// &mut fine, 1857c68be7a2SJeremy L Thompson /// )?; 1858e78171edSJeremy L Thompson /// let array = fine.view()?; 18599df49d7eSJed Brown /// for j in 0..p_fine { 18609df49d7eSJed Brown /// interp_c_to_f[j * p_coarse + i] = array[j]; 18619df49d7eSJed Brown /// } 18629df49d7eSJed Brown /// } 18639df49d7eSJed Brown /// } 1864c68be7a2SJeremy L Thompson /// let (op_mass_coarse, op_prolong, op_restrict) = op_mass_fine.create_multigrid_level_tensor_H1( 1865c68be7a2SJeremy L Thompson /// &multiplicity, 1866c68be7a2SJeremy L Thompson /// &ru_coarse, 1867c68be7a2SJeremy L Thompson /// &bu_coarse, 1868c68be7a2SJeremy L Thompson /// &interp_c_to_f, 1869c68be7a2SJeremy L Thompson /// )?; 18709df49d7eSJed Brown /// 18719df49d7eSJed Brown /// // Coarse problem 18729df49d7eSJed Brown /// u_coarse.set_value(1.0); 1873c68be7a2SJeremy L Thompson /// op_mass_coarse.apply(&u_coarse, &mut v_coarse)?; 18749df49d7eSJed Brown /// 18759df49d7eSJed Brown /// // Check 1876e78171edSJeremy L Thompson /// let sum: Scalar = v_coarse.view()?.iter().sum(); 18779df49d7eSJed Brown /// assert!( 187880a9ef05SNatalie Beams /// (sum - 2.0).abs() < 10.0 * libceed::EPSILON, 18799df49d7eSJed Brown /// "Incorrect interval length computed" 18809df49d7eSJed Brown /// ); 18819df49d7eSJed Brown /// 18829df49d7eSJed Brown /// // Prolong 1883c68be7a2SJeremy L Thompson /// op_prolong.apply(&u_coarse, &mut u_fine)?; 18849df49d7eSJed Brown /// 18859df49d7eSJed Brown /// // Fine problem 1886c68be7a2SJeremy L Thompson /// op_mass_fine.apply(&u_fine, &mut v_fine)?; 18879df49d7eSJed Brown /// 18889df49d7eSJed Brown /// // Check 1889e78171edSJeremy L Thompson /// let sum: Scalar = v_fine.view()?.iter().sum(); 18909df49d7eSJed Brown /// assert!( 1891392d940cSJeremy L Thompson /// (sum - 2.0).abs() < 200.0 * libceed::EPSILON, 18929df49d7eSJed Brown /// "Incorrect interval length computed" 18939df49d7eSJed Brown /// ); 18949df49d7eSJed Brown /// 18959df49d7eSJed Brown /// // Restrict 1896c68be7a2SJeremy L Thompson /// op_restrict.apply(&v_fine, &mut v_coarse)?; 18979df49d7eSJed Brown /// 18989df49d7eSJed Brown /// // Check 1899e78171edSJeremy L Thompson /// let sum: Scalar = v_coarse.view()?.iter().sum(); 19009df49d7eSJed Brown /// assert!( 1901392d940cSJeremy L Thompson /// (sum - 2.0).abs() < 200.0 * libceed::EPSILON, 19029df49d7eSJed Brown /// "Incorrect interval length computed" 19039df49d7eSJed Brown /// ); 1904c68be7a2SJeremy L Thompson /// # Ok(()) 1905c68be7a2SJeremy L Thompson /// # } 19069df49d7eSJed Brown /// ``` 1907594ef120SJeremy L Thompson pub fn create_multigrid_level_tensor_H1<'b>( 19089df49d7eSJed Brown &self, 19099df49d7eSJed Brown p_mult_fine: &Vector, 19109df49d7eSJed Brown rstr_coarse: &ElemRestriction, 19119df49d7eSJed Brown basis_coarse: &Basis, 191280a9ef05SNatalie Beams interpCtoF: &Vec<Scalar>, 1913594ef120SJeremy L Thompson ) -> crate::Result<(Operator<'b>, Operator<'b>, Operator<'b>)> { 19149df49d7eSJed Brown let mut ptr_coarse = std::ptr::null_mut(); 19159df49d7eSJed Brown let mut ptr_prolong = std::ptr::null_mut(); 19169df49d7eSJed Brown let mut ptr_restrict = std::ptr::null_mut(); 19179df49d7eSJed Brown let ierr = unsafe { 19189df49d7eSJed Brown bind_ceed::CeedOperatorMultigridLevelCreateTensorH1( 19199df49d7eSJed Brown self.op_core.ptr, 19209df49d7eSJed Brown p_mult_fine.ptr, 19219df49d7eSJed Brown rstr_coarse.ptr, 19229df49d7eSJed Brown basis_coarse.ptr, 19239df49d7eSJed Brown interpCtoF.as_ptr(), 19249df49d7eSJed Brown &mut ptr_coarse, 19259df49d7eSJed Brown &mut ptr_prolong, 19269df49d7eSJed Brown &mut ptr_restrict, 19279df49d7eSJed Brown ) 19289df49d7eSJed Brown }; 19291142270cSJeremy L Thompson self.op_core.check_error(ierr)?; 19301142270cSJeremy L Thompson let op_coarse = Operator::from_raw(ptr_coarse)?; 19311142270cSJeremy L Thompson let op_prolong = Operator::from_raw(ptr_prolong)?; 19321142270cSJeremy L Thompson let op_restrict = Operator::from_raw(ptr_restrict)?; 19339df49d7eSJed Brown Ok((op_coarse, op_prolong, op_restrict)) 19349df49d7eSJed Brown } 19359df49d7eSJed Brown 19369df49d7eSJed Brown /// Create a multigrid coarse Operator and level transfer Operators for a 19379df49d7eSJed Brown /// given Operator with a non-tensor basis for the active basis 19389df49d7eSJed Brown /// 19399df49d7eSJed Brown /// * `p_mult_fine` - Lvector multiplicity in parallel gather/scatter 19409df49d7eSJed Brown /// * `rstr_coarse` - Coarse grid restriction 19419df49d7eSJed Brown /// * `basis_coarse` - Coarse grid active vector basis 19429df49d7eSJed Brown /// * `interp_c_to_f` - Matrix for coarse to fine 19439df49d7eSJed Brown /// 19449df49d7eSJed Brown /// ``` 19459df49d7eSJed Brown /// # use libceed::prelude::*; 19464d27c890SJeremy L Thompson /// # fn main() -> libceed::Result<()> { 19479df49d7eSJed Brown /// # let ceed = libceed::Ceed::default_init(); 19489df49d7eSJed Brown /// let ne = 15; 19499df49d7eSJed Brown /// let p_coarse = 3; 19509df49d7eSJed Brown /// let p_fine = 5; 19519df49d7eSJed Brown /// let q = 6; 19529df49d7eSJed Brown /// let ndofs_coarse = p_coarse * ne - ne + 1; 19539df49d7eSJed Brown /// let ndofs_fine = p_fine * ne - ne + 1; 19549df49d7eSJed Brown /// 19559df49d7eSJed Brown /// // Vectors 19569df49d7eSJed Brown /// let x_array = (0..ne + 1) 195780a9ef05SNatalie Beams /// .map(|i| 2.0 * i as Scalar / ne as Scalar - 1.0) 195880a9ef05SNatalie Beams /// .collect::<Vec<Scalar>>(); 1959c68be7a2SJeremy L Thompson /// let x = ceed.vector_from_slice(&x_array)?; 1960c68be7a2SJeremy L Thompson /// let mut qdata = ceed.vector(ne * q)?; 19619df49d7eSJed Brown /// qdata.set_value(0.0); 1962c68be7a2SJeremy L Thompson /// let mut u_coarse = ceed.vector(ndofs_coarse)?; 19639df49d7eSJed Brown /// u_coarse.set_value(1.0); 1964c68be7a2SJeremy L Thompson /// let mut u_fine = ceed.vector(ndofs_fine)?; 19659df49d7eSJed Brown /// u_fine.set_value(1.0); 1966c68be7a2SJeremy L Thompson /// let mut v_coarse = ceed.vector(ndofs_coarse)?; 19679df49d7eSJed Brown /// v_coarse.set_value(0.0); 1968c68be7a2SJeremy L Thompson /// let mut v_fine = ceed.vector(ndofs_fine)?; 19699df49d7eSJed Brown /// v_fine.set_value(0.0); 1970c68be7a2SJeremy L Thompson /// let mut multiplicity = ceed.vector(ndofs_fine)?; 19719df49d7eSJed Brown /// multiplicity.set_value(1.0); 19729df49d7eSJed Brown /// 19739df49d7eSJed Brown /// // Restrictions 19749df49d7eSJed Brown /// let strides: [i32; 3] = [1, q as i32, q as i32]; 1975c68be7a2SJeremy L Thompson /// let rq = ceed.strided_elem_restriction(ne, q, 1, q * ne, strides)?; 19769df49d7eSJed Brown /// 19779df49d7eSJed Brown /// let mut indx: Vec<i32> = vec![0; 2 * ne]; 19789df49d7eSJed Brown /// for i in 0..ne { 19799df49d7eSJed Brown /// indx[2 * i + 0] = i as i32; 19809df49d7eSJed Brown /// indx[2 * i + 1] = (i + 1) as i32; 19819df49d7eSJed Brown /// } 1982c68be7a2SJeremy L Thompson /// let rx = ceed.elem_restriction(ne, 2, 1, 1, ne + 1, MemType::Host, &indx)?; 19839df49d7eSJed Brown /// 19849df49d7eSJed Brown /// let mut indu_coarse: Vec<i32> = vec![0; p_coarse * ne]; 19859df49d7eSJed Brown /// for i in 0..ne { 19869df49d7eSJed Brown /// for j in 0..p_coarse { 19879df49d7eSJed Brown /// indu_coarse[p_coarse * i + j] = (i + j) as i32; 19889df49d7eSJed Brown /// } 19899df49d7eSJed Brown /// } 1990c68be7a2SJeremy L Thompson /// let ru_coarse = ceed.elem_restriction( 19919df49d7eSJed Brown /// ne, 19929df49d7eSJed Brown /// p_coarse, 19939df49d7eSJed Brown /// 1, 19949df49d7eSJed Brown /// 1, 19959df49d7eSJed Brown /// ndofs_coarse, 19969df49d7eSJed Brown /// MemType::Host, 19979df49d7eSJed Brown /// &indu_coarse, 1998c68be7a2SJeremy L Thompson /// )?; 19999df49d7eSJed Brown /// 20009df49d7eSJed Brown /// let mut indu_fine: Vec<i32> = vec![0; p_fine * ne]; 20019df49d7eSJed Brown /// for i in 0..ne { 20029df49d7eSJed Brown /// for j in 0..p_fine { 20039df49d7eSJed Brown /// indu_fine[p_fine * i + j] = (i + j) as i32; 20049df49d7eSJed Brown /// } 20059df49d7eSJed Brown /// } 2006c68be7a2SJeremy L Thompson /// let ru_fine = ceed.elem_restriction(ne, p_fine, 1, 1, ndofs_fine, MemType::Host, &indu_fine)?; 20079df49d7eSJed Brown /// 20089df49d7eSJed Brown /// // Bases 2009c68be7a2SJeremy L Thompson /// let bx = ceed.basis_tensor_H1_Lagrange(1, 1, 2, q, QuadMode::Gauss)?; 2010c68be7a2SJeremy L Thompson /// let bu_coarse = ceed.basis_tensor_H1_Lagrange(1, 1, p_coarse, q, QuadMode::Gauss)?; 2011c68be7a2SJeremy L Thompson /// let bu_fine = ceed.basis_tensor_H1_Lagrange(1, 1, p_fine, q, QuadMode::Gauss)?; 20129df49d7eSJed Brown /// 20139df49d7eSJed Brown /// // Build quadrature data 2014c68be7a2SJeremy L Thompson /// let qf_build = ceed.q_function_interior_by_name("Mass1DBuild")?; 2015c68be7a2SJeremy L Thompson /// ceed.operator(&qf_build, QFunctionOpt::None, QFunctionOpt::None)? 2016c68be7a2SJeremy L Thompson /// .field("dx", &rx, &bx, VectorOpt::Active)? 2017c68be7a2SJeremy L Thompson /// .field("weights", ElemRestrictionOpt::None, &bx, VectorOpt::None)? 2018356036faSJeremy L Thompson /// .field("qdata", &rq, BasisOpt::None, VectorOpt::Active)? 2019c68be7a2SJeremy L Thompson /// .apply(&x, &mut qdata)?; 20209df49d7eSJed Brown /// 20219df49d7eSJed Brown /// // Mass operator 2022c68be7a2SJeremy L Thompson /// let qf_mass = ceed.q_function_interior_by_name("MassApply")?; 20239df49d7eSJed Brown /// let op_mass_fine = ceed 2024c68be7a2SJeremy L Thompson /// .operator(&qf_mass, QFunctionOpt::None, QFunctionOpt::None)? 2025c68be7a2SJeremy L Thompson /// .field("u", &ru_fine, &bu_fine, VectorOpt::Active)? 2026356036faSJeremy L Thompson /// .field("qdata", &rq, BasisOpt::None, &qdata)? 2027c68be7a2SJeremy L Thompson /// .field("v", &ru_fine, &bu_fine, VectorOpt::Active)?; 20289df49d7eSJed Brown /// 20299df49d7eSJed Brown /// // Multigrid setup 203080a9ef05SNatalie Beams /// let mut interp_c_to_f: Vec<Scalar> = vec![0.; p_coarse * p_fine]; 20319df49d7eSJed Brown /// { 2032c68be7a2SJeremy L Thompson /// let mut coarse = ceed.vector(p_coarse)?; 2033c68be7a2SJeremy L Thompson /// let mut fine = ceed.vector(p_fine)?; 2034c68be7a2SJeremy L Thompson /// let basis_c_to_f = 2035c68be7a2SJeremy L Thompson /// ceed.basis_tensor_H1_Lagrange(1, 1, p_coarse, p_fine, QuadMode::GaussLobatto)?; 20369df49d7eSJed Brown /// for i in 0..p_coarse { 20379df49d7eSJed Brown /// coarse.set_value(0.0); 20389df49d7eSJed Brown /// { 2039e78171edSJeremy L Thompson /// let mut array = coarse.view_mut()?; 20409df49d7eSJed Brown /// array[i] = 1.; 20419df49d7eSJed Brown /// } 2042c68be7a2SJeremy L Thompson /// basis_c_to_f.apply( 20439df49d7eSJed Brown /// 1, 20449df49d7eSJed Brown /// TransposeMode::NoTranspose, 20459df49d7eSJed Brown /// EvalMode::Interp, 20469df49d7eSJed Brown /// &coarse, 20479df49d7eSJed Brown /// &mut fine, 2048c68be7a2SJeremy L Thompson /// )?; 2049e78171edSJeremy L Thompson /// let array = fine.view()?; 20509df49d7eSJed Brown /// for j in 0..p_fine { 20519df49d7eSJed Brown /// interp_c_to_f[j * p_coarse + i] = array[j]; 20529df49d7eSJed Brown /// } 20539df49d7eSJed Brown /// } 20549df49d7eSJed Brown /// } 2055c68be7a2SJeremy L Thompson /// let (op_mass_coarse, op_prolong, op_restrict) = op_mass_fine.create_multigrid_level_H1( 2056c68be7a2SJeremy L Thompson /// &multiplicity, 2057c68be7a2SJeremy L Thompson /// &ru_coarse, 2058c68be7a2SJeremy L Thompson /// &bu_coarse, 2059c68be7a2SJeremy L Thompson /// &interp_c_to_f, 2060c68be7a2SJeremy L Thompson /// )?; 20619df49d7eSJed Brown /// 20629df49d7eSJed Brown /// // Coarse problem 20639df49d7eSJed Brown /// u_coarse.set_value(1.0); 2064c68be7a2SJeremy L Thompson /// op_mass_coarse.apply(&u_coarse, &mut v_coarse)?; 20659df49d7eSJed Brown /// 20669df49d7eSJed Brown /// // Check 2067e78171edSJeremy L Thompson /// let sum: Scalar = v_coarse.view()?.iter().sum(); 20689df49d7eSJed Brown /// assert!( 206980a9ef05SNatalie Beams /// (sum - 2.0).abs() < 10.0 * libceed::EPSILON, 20709df49d7eSJed Brown /// "Incorrect interval length computed" 20719df49d7eSJed Brown /// ); 20729df49d7eSJed Brown /// 20739df49d7eSJed Brown /// // Prolong 2074c68be7a2SJeremy L Thompson /// op_prolong.apply(&u_coarse, &mut u_fine)?; 20759df49d7eSJed Brown /// 20769df49d7eSJed Brown /// // Fine problem 2077c68be7a2SJeremy L Thompson /// op_mass_fine.apply(&u_fine, &mut v_fine)?; 20789df49d7eSJed Brown /// 20799df49d7eSJed Brown /// // Check 2080e78171edSJeremy L Thompson /// let sum: Scalar = v_fine.view()?.iter().sum(); 20819df49d7eSJed Brown /// assert!( 2082392d940cSJeremy L Thompson /// (sum - 2.0).abs() < 200.0 * libceed::EPSILON, 20839df49d7eSJed Brown /// "Incorrect interval length computed" 20849df49d7eSJed Brown /// ); 20859df49d7eSJed Brown /// 20869df49d7eSJed Brown /// // Restrict 2087c68be7a2SJeremy L Thompson /// op_restrict.apply(&v_fine, &mut v_coarse)?; 20889df49d7eSJed Brown /// 20899df49d7eSJed Brown /// // Check 2090e78171edSJeremy L Thompson /// let sum: Scalar = v_coarse.view()?.iter().sum(); 20919df49d7eSJed Brown /// assert!( 2092392d940cSJeremy L Thompson /// (sum - 2.0).abs() < 200.0 * libceed::EPSILON, 20939df49d7eSJed Brown /// "Incorrect interval length computed" 20949df49d7eSJed Brown /// ); 2095c68be7a2SJeremy L Thompson /// # Ok(()) 2096c68be7a2SJeremy L Thompson /// # } 20979df49d7eSJed Brown /// ``` 2098594ef120SJeremy L Thompson pub fn create_multigrid_level_H1<'b>( 20999df49d7eSJed Brown &self, 21009df49d7eSJed Brown p_mult_fine: &Vector, 21019df49d7eSJed Brown rstr_coarse: &ElemRestriction, 21029df49d7eSJed Brown basis_coarse: &Basis, 210380a9ef05SNatalie Beams interpCtoF: &[Scalar], 2104594ef120SJeremy L Thompson ) -> crate::Result<(Operator<'b>, Operator<'b>, Operator<'b>)> { 21059df49d7eSJed Brown let mut ptr_coarse = std::ptr::null_mut(); 21069df49d7eSJed Brown let mut ptr_prolong = std::ptr::null_mut(); 21079df49d7eSJed Brown let mut ptr_restrict = std::ptr::null_mut(); 21089df49d7eSJed Brown let ierr = unsafe { 21099df49d7eSJed Brown bind_ceed::CeedOperatorMultigridLevelCreateH1( 21109df49d7eSJed Brown self.op_core.ptr, 21119df49d7eSJed Brown p_mult_fine.ptr, 21129df49d7eSJed Brown rstr_coarse.ptr, 21139df49d7eSJed Brown basis_coarse.ptr, 21149df49d7eSJed Brown interpCtoF.as_ptr(), 21159df49d7eSJed Brown &mut ptr_coarse, 21169df49d7eSJed Brown &mut ptr_prolong, 21179df49d7eSJed Brown &mut ptr_restrict, 21189df49d7eSJed Brown ) 21199df49d7eSJed Brown }; 21201142270cSJeremy L Thompson self.op_core.check_error(ierr)?; 21211142270cSJeremy L Thompson let op_coarse = Operator::from_raw(ptr_coarse)?; 21221142270cSJeremy L Thompson let op_prolong = Operator::from_raw(ptr_prolong)?; 21231142270cSJeremy L Thompson let op_restrict = Operator::from_raw(ptr_restrict)?; 21249df49d7eSJed Brown Ok((op_coarse, op_prolong, op_restrict)) 21259df49d7eSJed Brown } 21269df49d7eSJed Brown } 21279df49d7eSJed Brown 21289df49d7eSJed Brown // ----------------------------------------------------------------------------- 21299df49d7eSJed Brown // Composite Operator 21309df49d7eSJed Brown // ----------------------------------------------------------------------------- 21319df49d7eSJed Brown impl<'a> CompositeOperator<'a> { 21329df49d7eSJed Brown // Constructor 2133594ef120SJeremy L Thompson pub fn create(ceed: &crate::Ceed) -> crate::Result<Self> { 21349df49d7eSJed Brown let mut ptr = std::ptr::null_mut(); 21359df49d7eSJed Brown let ierr = unsafe { bind_ceed::CeedCompositeOperatorCreate(ceed.ptr, &mut ptr) }; 21369df49d7eSJed Brown ceed.check_error(ierr)?; 21379df49d7eSJed Brown Ok(Self { 21381142270cSJeremy L Thompson op_core: OperatorCore { 21391142270cSJeremy L Thompson ptr, 21401142270cSJeremy L Thompson _lifeline: PhantomData, 21411142270cSJeremy L Thompson }, 21429df49d7eSJed Brown }) 21439df49d7eSJed Brown } 21449df49d7eSJed Brown 2145ea6b5821SJeremy L Thompson /// Set name for CompositeOperator printing 2146ea6b5821SJeremy L Thompson /// 2147ea6b5821SJeremy L Thompson /// * 'name' - Name to set 2148ea6b5821SJeremy L Thompson /// 2149ea6b5821SJeremy L Thompson /// ``` 2150ea6b5821SJeremy L Thompson /// # use libceed::prelude::*; 2151ea6b5821SJeremy L Thompson /// # fn main() -> libceed::Result<()> { 2152ea6b5821SJeremy L Thompson /// # let ceed = libceed::Ceed::default_init(); 2153ea6b5821SJeremy L Thompson /// 2154ea6b5821SJeremy L Thompson /// // Sub operator field arguments 2155ea6b5821SJeremy L Thompson /// let ne = 3; 2156ea6b5821SJeremy L Thompson /// let q = 4 as usize; 2157ea6b5821SJeremy L Thompson /// let mut ind: Vec<i32> = vec![0; 2 * ne]; 2158ea6b5821SJeremy L Thompson /// for i in 0..ne { 2159ea6b5821SJeremy L Thompson /// ind[2 * i + 0] = i as i32; 2160ea6b5821SJeremy L Thompson /// ind[2 * i + 1] = (i + 1) as i32; 2161ea6b5821SJeremy L Thompson /// } 2162ea6b5821SJeremy L Thompson /// let r = ceed.elem_restriction(ne, 2, 1, 1, ne + 1, MemType::Host, &ind)?; 2163ea6b5821SJeremy L Thompson /// let strides: [i32; 3] = [1, q as i32, q as i32]; 2164d9267b76SJeremy L Thompson /// let rq = ceed.strided_elem_restriction(ne, q, 1, q * ne, strides)?; 2165ea6b5821SJeremy L Thompson /// 2166ea6b5821SJeremy L Thompson /// let b = ceed.basis_tensor_H1_Lagrange(1, 1, 2, q, QuadMode::Gauss)?; 2167ea6b5821SJeremy L Thompson /// 2168ea6b5821SJeremy L Thompson /// let qdata_mass = ceed.vector(q * ne)?; 2169ea6b5821SJeremy L Thompson /// let qdata_diff = ceed.vector(q * ne)?; 2170ea6b5821SJeremy L Thompson /// 2171ea6b5821SJeremy L Thompson /// let qf_mass = ceed.q_function_interior_by_name("MassApply")?; 2172ea6b5821SJeremy L Thompson /// let op_mass = ceed 2173ea6b5821SJeremy L Thompson /// .operator(&qf_mass, QFunctionOpt::None, QFunctionOpt::None)? 2174ea6b5821SJeremy L Thompson /// .name("Mass term")? 2175ea6b5821SJeremy L Thompson /// .field("u", &r, &b, VectorOpt::Active)? 2176356036faSJeremy L Thompson /// .field("qdata", &rq, BasisOpt::None, &qdata_mass)? 2177ea6b5821SJeremy L Thompson /// .field("v", &r, &b, VectorOpt::Active)?; 2178ea6b5821SJeremy L Thompson /// 2179ea6b5821SJeremy L Thompson /// let qf_diff = ceed.q_function_interior_by_name("Poisson1DApply")?; 2180ea6b5821SJeremy L Thompson /// let op_diff = ceed 2181ea6b5821SJeremy L Thompson /// .operator(&qf_diff, QFunctionOpt::None, QFunctionOpt::None)? 2182ea6b5821SJeremy L Thompson /// .name("Poisson term")? 2183ea6b5821SJeremy L Thompson /// .field("du", &r, &b, VectorOpt::Active)? 2184356036faSJeremy L Thompson /// .field("qdata", &rq, BasisOpt::None, &qdata_diff)? 2185ea6b5821SJeremy L Thompson /// .field("dv", &r, &b, VectorOpt::Active)?; 2186ea6b5821SJeremy L Thompson /// 2187ea6b5821SJeremy L Thompson /// let op = ceed 2188ea6b5821SJeremy L Thompson /// .composite_operator()? 2189ea6b5821SJeremy L Thompson /// .name("Screened Poisson")? 2190ea6b5821SJeremy L Thompson /// .sub_operator(&op_mass)? 2191ea6b5821SJeremy L Thompson /// .sub_operator(&op_diff)?; 2192ea6b5821SJeremy L Thompson /// # Ok(()) 2193ea6b5821SJeremy L Thompson /// # } 2194ea6b5821SJeremy L Thompson /// ``` 2195ea6b5821SJeremy L Thompson #[allow(unused_mut)] 2196ea6b5821SJeremy L Thompson pub fn name(mut self, name: &str) -> crate::Result<Self> { 2197ea6b5821SJeremy L Thompson self.op_core.name(name)?; 2198ea6b5821SJeremy L Thompson Ok(self) 2199ea6b5821SJeremy L Thompson } 2200ea6b5821SJeremy L Thompson 22019df49d7eSJed Brown /// Apply Operator to a vector 22029df49d7eSJed Brown /// 2203ea6b5821SJeremy L Thompson /// * `input` - Inpuht Vector 22049df49d7eSJed Brown /// * `output` - Output Vector 22059df49d7eSJed Brown /// 22069df49d7eSJed Brown /// ``` 22079df49d7eSJed Brown /// # use libceed::prelude::*; 22084d27c890SJeremy L Thompson /// # fn main() -> libceed::Result<()> { 22099df49d7eSJed Brown /// # let ceed = libceed::Ceed::default_init(); 22109df49d7eSJed Brown /// let ne = 4; 22119df49d7eSJed Brown /// let p = 3; 22129df49d7eSJed Brown /// let q = 4; 22139df49d7eSJed Brown /// let ndofs = p * ne - ne + 1; 22149df49d7eSJed Brown /// 22159df49d7eSJed Brown /// // Vectors 2216c68be7a2SJeremy L Thompson /// let x = ceed.vector_from_slice(&[-1., -0.5, 0.0, 0.5, 1.0])?; 2217c68be7a2SJeremy L Thompson /// let mut qdata_mass = ceed.vector(ne * q)?; 22189df49d7eSJed Brown /// qdata_mass.set_value(0.0); 2219c68be7a2SJeremy L Thompson /// let mut qdata_diff = ceed.vector(ne * q)?; 22209df49d7eSJed Brown /// qdata_diff.set_value(0.0); 2221c68be7a2SJeremy L Thompson /// let mut u = ceed.vector(ndofs)?; 22229df49d7eSJed Brown /// u.set_value(1.0); 2223c68be7a2SJeremy L Thompson /// let mut v = ceed.vector(ndofs)?; 22249df49d7eSJed Brown /// v.set_value(0.0); 22259df49d7eSJed Brown /// 22269df49d7eSJed Brown /// // Restrictions 22279df49d7eSJed Brown /// let mut indx: Vec<i32> = vec![0; 2 * ne]; 22289df49d7eSJed Brown /// for i in 0..ne { 22299df49d7eSJed Brown /// indx[2 * i + 0] = i as i32; 22309df49d7eSJed Brown /// indx[2 * i + 1] = (i + 1) as i32; 22319df49d7eSJed Brown /// } 2232c68be7a2SJeremy L Thompson /// let rx = ceed.elem_restriction(ne, 2, 1, 1, ne + 1, MemType::Host, &indx)?; 22339df49d7eSJed Brown /// let mut indu: Vec<i32> = vec![0; p * ne]; 22349df49d7eSJed Brown /// for i in 0..ne { 22359df49d7eSJed Brown /// indu[p * i + 0] = i as i32; 22369df49d7eSJed Brown /// indu[p * i + 1] = (i + 1) as i32; 22379df49d7eSJed Brown /// indu[p * i + 2] = (i + 2) as i32; 22389df49d7eSJed Brown /// } 2239c68be7a2SJeremy L Thompson /// let ru = ceed.elem_restriction(ne, 3, 1, 1, ndofs, MemType::Host, &indu)?; 22409df49d7eSJed Brown /// let strides: [i32; 3] = [1, q as i32, q as i32]; 2241c68be7a2SJeremy L Thompson /// let rq = ceed.strided_elem_restriction(ne, q, 1, q * ne, strides)?; 22429df49d7eSJed Brown /// 22439df49d7eSJed Brown /// // Bases 2244c68be7a2SJeremy L Thompson /// let bx = ceed.basis_tensor_H1_Lagrange(1, 1, 2, q, QuadMode::Gauss)?; 2245c68be7a2SJeremy L Thompson /// let bu = ceed.basis_tensor_H1_Lagrange(1, 1, p, q, QuadMode::Gauss)?; 22469df49d7eSJed Brown /// 22479df49d7eSJed Brown /// // Build quadrature data 2248c68be7a2SJeremy L Thompson /// let qf_build_mass = ceed.q_function_interior_by_name("Mass1DBuild")?; 2249c68be7a2SJeremy L Thompson /// ceed.operator(&qf_build_mass, QFunctionOpt::None, QFunctionOpt::None)? 2250c68be7a2SJeremy L Thompson /// .field("dx", &rx, &bx, VectorOpt::Active)? 2251c68be7a2SJeremy L Thompson /// .field("weights", ElemRestrictionOpt::None, &bx, VectorOpt::None)? 2252356036faSJeremy L Thompson /// .field("qdata", &rq, BasisOpt::None, VectorOpt::Active)? 2253c68be7a2SJeremy L Thompson /// .apply(&x, &mut qdata_mass)?; 22549df49d7eSJed Brown /// 2255c68be7a2SJeremy L Thompson /// let qf_build_diff = ceed.q_function_interior_by_name("Poisson1DBuild")?; 2256c68be7a2SJeremy L Thompson /// ceed.operator(&qf_build_diff, QFunctionOpt::None, QFunctionOpt::None)? 2257c68be7a2SJeremy L Thompson /// .field("dx", &rx, &bx, VectorOpt::Active)? 2258c68be7a2SJeremy L Thompson /// .field("weights", ElemRestrictionOpt::None, &bx, VectorOpt::None)? 2259356036faSJeremy L Thompson /// .field("qdata", &rq, BasisOpt::None, VectorOpt::Active)? 2260c68be7a2SJeremy L Thompson /// .apply(&x, &mut qdata_diff)?; 22619df49d7eSJed Brown /// 22629df49d7eSJed Brown /// // Application operator 2263c68be7a2SJeremy L Thompson /// let qf_mass = ceed.q_function_interior_by_name("MassApply")?; 22649df49d7eSJed Brown /// let op_mass = ceed 2265c68be7a2SJeremy L Thompson /// .operator(&qf_mass, QFunctionOpt::None, QFunctionOpt::None)? 2266c68be7a2SJeremy L Thompson /// .field("u", &ru, &bu, VectorOpt::Active)? 2267356036faSJeremy L Thompson /// .field("qdata", &rq, BasisOpt::None, &qdata_mass)? 2268c68be7a2SJeremy L Thompson /// .field("v", &ru, &bu, VectorOpt::Active)?; 22699df49d7eSJed Brown /// 2270c68be7a2SJeremy L Thompson /// let qf_diff = ceed.q_function_interior_by_name("Poisson1DApply")?; 22719df49d7eSJed Brown /// let op_diff = ceed 2272c68be7a2SJeremy L Thompson /// .operator(&qf_diff, QFunctionOpt::None, QFunctionOpt::None)? 2273c68be7a2SJeremy L Thompson /// .field("du", &ru, &bu, VectorOpt::Active)? 2274356036faSJeremy L Thompson /// .field("qdata", &rq, BasisOpt::None, &qdata_diff)? 2275c68be7a2SJeremy L Thompson /// .field("dv", &ru, &bu, VectorOpt::Active)?; 22769df49d7eSJed Brown /// 22779df49d7eSJed Brown /// let op_composite = ceed 2278c68be7a2SJeremy L Thompson /// .composite_operator()? 2279c68be7a2SJeremy L Thompson /// .sub_operator(&op_mass)? 2280c68be7a2SJeremy L Thompson /// .sub_operator(&op_diff)?; 22819df49d7eSJed Brown /// 22829df49d7eSJed Brown /// v.set_value(0.0); 2283c68be7a2SJeremy L Thompson /// op_composite.apply(&u, &mut v)?; 22849df49d7eSJed Brown /// 22859df49d7eSJed Brown /// // Check 2286e78171edSJeremy L Thompson /// let sum: Scalar = v.view()?.iter().sum(); 22879df49d7eSJed Brown /// assert!( 228880a9ef05SNatalie Beams /// (sum - 2.0).abs() < 10.0 * libceed::EPSILON, 22899df49d7eSJed Brown /// "Incorrect interval length computed" 22909df49d7eSJed Brown /// ); 2291c68be7a2SJeremy L Thompson /// # Ok(()) 2292c68be7a2SJeremy L Thompson /// # } 22939df49d7eSJed Brown /// ``` 22949df49d7eSJed Brown pub fn apply(&self, input: &Vector, output: &mut Vector) -> crate::Result<i32> { 22959df49d7eSJed Brown self.op_core.apply(input, output) 22969df49d7eSJed Brown } 22979df49d7eSJed Brown 22989df49d7eSJed Brown /// Apply Operator to a vector and add result to output vector 22999df49d7eSJed Brown /// 23009df49d7eSJed Brown /// * `input` - Input Vector 23019df49d7eSJed Brown /// * `output` - Output Vector 23029df49d7eSJed Brown /// 23039df49d7eSJed Brown /// ``` 23049df49d7eSJed Brown /// # use libceed::prelude::*; 23054d27c890SJeremy L Thompson /// # fn main() -> libceed::Result<()> { 23069df49d7eSJed Brown /// # let ceed = libceed::Ceed::default_init(); 23079df49d7eSJed Brown /// let ne = 4; 23089df49d7eSJed Brown /// let p = 3; 23099df49d7eSJed Brown /// let q = 4; 23109df49d7eSJed Brown /// let ndofs = p * ne - ne + 1; 23119df49d7eSJed Brown /// 23129df49d7eSJed Brown /// // Vectors 2313c68be7a2SJeremy L Thompson /// let x = ceed.vector_from_slice(&[-1., -0.5, 0.0, 0.5, 1.0])?; 2314c68be7a2SJeremy L Thompson /// let mut qdata_mass = ceed.vector(ne * q)?; 23159df49d7eSJed Brown /// qdata_mass.set_value(0.0); 2316c68be7a2SJeremy L Thompson /// let mut qdata_diff = ceed.vector(ne * q)?; 23179df49d7eSJed Brown /// qdata_diff.set_value(0.0); 2318c68be7a2SJeremy L Thompson /// let mut u = ceed.vector(ndofs)?; 23199df49d7eSJed Brown /// u.set_value(1.0); 2320c68be7a2SJeremy L Thompson /// let mut v = ceed.vector(ndofs)?; 23219df49d7eSJed Brown /// v.set_value(0.0); 23229df49d7eSJed Brown /// 23239df49d7eSJed Brown /// // Restrictions 23249df49d7eSJed Brown /// let mut indx: Vec<i32> = vec![0; 2 * ne]; 23259df49d7eSJed Brown /// for i in 0..ne { 23269df49d7eSJed Brown /// indx[2 * i + 0] = i as i32; 23279df49d7eSJed Brown /// indx[2 * i + 1] = (i + 1) as i32; 23289df49d7eSJed Brown /// } 2329c68be7a2SJeremy L Thompson /// let rx = ceed.elem_restriction(ne, 2, 1, 1, ne + 1, MemType::Host, &indx)?; 23309df49d7eSJed Brown /// let mut indu: Vec<i32> = vec![0; p * ne]; 23319df49d7eSJed Brown /// for i in 0..ne { 23329df49d7eSJed Brown /// indu[p * i + 0] = i as i32; 23339df49d7eSJed Brown /// indu[p * i + 1] = (i + 1) as i32; 23349df49d7eSJed Brown /// indu[p * i + 2] = (i + 2) as i32; 23359df49d7eSJed Brown /// } 2336c68be7a2SJeremy L Thompson /// let ru = ceed.elem_restriction(ne, 3, 1, 1, ndofs, MemType::Host, &indu)?; 23379df49d7eSJed Brown /// let strides: [i32; 3] = [1, q as i32, q as i32]; 2338c68be7a2SJeremy L Thompson /// let rq = ceed.strided_elem_restriction(ne, q, 1, q * ne, strides)?; 23399df49d7eSJed Brown /// 23409df49d7eSJed Brown /// // Bases 2341c68be7a2SJeremy L Thompson /// let bx = ceed.basis_tensor_H1_Lagrange(1, 1, 2, q, QuadMode::Gauss)?; 2342c68be7a2SJeremy L Thompson /// let bu = ceed.basis_tensor_H1_Lagrange(1, 1, p, q, QuadMode::Gauss)?; 23439df49d7eSJed Brown /// 23449df49d7eSJed Brown /// // Build quadrature data 2345c68be7a2SJeremy L Thompson /// let qf_build_mass = ceed.q_function_interior_by_name("Mass1DBuild")?; 2346c68be7a2SJeremy L Thompson /// ceed.operator(&qf_build_mass, QFunctionOpt::None, QFunctionOpt::None)? 2347c68be7a2SJeremy L Thompson /// .field("dx", &rx, &bx, VectorOpt::Active)? 2348c68be7a2SJeremy L Thompson /// .field("weights", ElemRestrictionOpt::None, &bx, VectorOpt::None)? 2349356036faSJeremy L Thompson /// .field("qdata", &rq, BasisOpt::None, VectorOpt::Active)? 2350c68be7a2SJeremy L Thompson /// .apply(&x, &mut qdata_mass)?; 23519df49d7eSJed Brown /// 2352c68be7a2SJeremy L Thompson /// let qf_build_diff = ceed.q_function_interior_by_name("Poisson1DBuild")?; 2353c68be7a2SJeremy L Thompson /// ceed.operator(&qf_build_diff, QFunctionOpt::None, QFunctionOpt::None)? 2354c68be7a2SJeremy L Thompson /// .field("dx", &rx, &bx, VectorOpt::Active)? 2355c68be7a2SJeremy L Thompson /// .field("weights", ElemRestrictionOpt::None, &bx, VectorOpt::None)? 2356356036faSJeremy L Thompson /// .field("qdata", &rq, BasisOpt::None, VectorOpt::Active)? 2357c68be7a2SJeremy L Thompson /// .apply(&x, &mut qdata_diff)?; 23589df49d7eSJed Brown /// 23599df49d7eSJed Brown /// // Application operator 2360c68be7a2SJeremy L Thompson /// let qf_mass = ceed.q_function_interior_by_name("MassApply")?; 23619df49d7eSJed Brown /// let op_mass = ceed 2362c68be7a2SJeremy L Thompson /// .operator(&qf_mass, QFunctionOpt::None, QFunctionOpt::None)? 2363c68be7a2SJeremy L Thompson /// .field("u", &ru, &bu, VectorOpt::Active)? 2364356036faSJeremy L Thompson /// .field("qdata", &rq, BasisOpt::None, &qdata_mass)? 2365c68be7a2SJeremy L Thompson /// .field("v", &ru, &bu, VectorOpt::Active)?; 23669df49d7eSJed Brown /// 2367c68be7a2SJeremy L Thompson /// let qf_diff = ceed.q_function_interior_by_name("Poisson1DApply")?; 23689df49d7eSJed Brown /// let op_diff = ceed 2369c68be7a2SJeremy L Thompson /// .operator(&qf_diff, QFunctionOpt::None, QFunctionOpt::None)? 2370c68be7a2SJeremy L Thompson /// .field("du", &ru, &bu, VectorOpt::Active)? 2371356036faSJeremy L Thompson /// .field("qdata", &rq, BasisOpt::None, &qdata_diff)? 2372c68be7a2SJeremy L Thompson /// .field("dv", &ru, &bu, VectorOpt::Active)?; 23739df49d7eSJed Brown /// 23749df49d7eSJed Brown /// let op_composite = ceed 2375c68be7a2SJeremy L Thompson /// .composite_operator()? 2376c68be7a2SJeremy L Thompson /// .sub_operator(&op_mass)? 2377c68be7a2SJeremy L Thompson /// .sub_operator(&op_diff)?; 23789df49d7eSJed Brown /// 23799df49d7eSJed Brown /// v.set_value(1.0); 2380c68be7a2SJeremy L Thompson /// op_composite.apply_add(&u, &mut v)?; 23819df49d7eSJed Brown /// 23829df49d7eSJed Brown /// // Check 2383e78171edSJeremy L Thompson /// let sum: Scalar = v.view()?.iter().sum(); 23849df49d7eSJed Brown /// assert!( 238580a9ef05SNatalie Beams /// (sum - (2.0 + ndofs as Scalar)).abs() < 10.0 * libceed::EPSILON, 23869df49d7eSJed Brown /// "Incorrect interval length computed" 23879df49d7eSJed Brown /// ); 2388c68be7a2SJeremy L Thompson /// # Ok(()) 2389c68be7a2SJeremy L Thompson /// # } 23909df49d7eSJed Brown /// ``` 23919df49d7eSJed Brown pub fn apply_add(&self, input: &Vector, output: &mut Vector) -> crate::Result<i32> { 23929df49d7eSJed Brown self.op_core.apply_add(input, output) 23939df49d7eSJed Brown } 23949df49d7eSJed Brown 23959df49d7eSJed Brown /// Add a sub-Operator to a Composite Operator 23969df49d7eSJed Brown /// 23979df49d7eSJed Brown /// * `subop` - Sub-Operator 23989df49d7eSJed Brown /// 23999df49d7eSJed Brown /// ``` 24009df49d7eSJed Brown /// # use libceed::prelude::*; 24014d27c890SJeremy L Thompson /// # fn main() -> libceed::Result<()> { 24029df49d7eSJed Brown /// # let ceed = libceed::Ceed::default_init(); 2403c68be7a2SJeremy L Thompson /// let mut op = ceed.composite_operator()?; 24049df49d7eSJed Brown /// 2405c68be7a2SJeremy L Thompson /// let qf_mass = ceed.q_function_interior_by_name("MassApply")?; 2406c68be7a2SJeremy L Thompson /// let op_mass = ceed.operator(&qf_mass, QFunctionOpt::None, QFunctionOpt::None)?; 2407c68be7a2SJeremy L Thompson /// op = op.sub_operator(&op_mass)?; 24089df49d7eSJed Brown /// 2409c68be7a2SJeremy L Thompson /// let qf_diff = ceed.q_function_interior_by_name("Poisson1DApply")?; 2410c68be7a2SJeremy L Thompson /// let op_diff = ceed.operator(&qf_diff, QFunctionOpt::None, QFunctionOpt::None)?; 2411c68be7a2SJeremy L Thompson /// op = op.sub_operator(&op_diff)?; 2412c68be7a2SJeremy L Thompson /// # Ok(()) 2413c68be7a2SJeremy L Thompson /// # } 24149df49d7eSJed Brown /// ``` 24159df49d7eSJed Brown #[allow(unused_mut)] 24169df49d7eSJed Brown pub fn sub_operator(mut self, subop: &Operator) -> crate::Result<Self> { 24179df49d7eSJed Brown let ierr = 24189df49d7eSJed Brown unsafe { bind_ceed::CeedCompositeOperatorAddSub(self.op_core.ptr, subop.op_core.ptr) }; 24191142270cSJeremy L Thompson self.op_core.check_error(ierr)?; 24209df49d7eSJed Brown Ok(self) 24219df49d7eSJed Brown } 24229df49d7eSJed Brown 24236f97ff0aSJeremy L Thompson /// Check if CompositeOperator is setup correctly 24246f97ff0aSJeremy L Thompson /// 24256f97ff0aSJeremy L Thompson /// ``` 24266f97ff0aSJeremy L Thompson /// # use libceed::prelude::*; 24276f97ff0aSJeremy L Thompson /// # fn main() -> libceed::Result<()> { 24286f97ff0aSJeremy L Thompson /// # let ceed = libceed::Ceed::default_init(); 24296f97ff0aSJeremy L Thompson /// let ne = 4; 24306f97ff0aSJeremy L Thompson /// let p = 3; 24316f97ff0aSJeremy L Thompson /// let q = 4; 24326f97ff0aSJeremy L Thompson /// let ndofs = p * ne - ne + 1; 24336f97ff0aSJeremy L Thompson /// 24346f97ff0aSJeremy L Thompson /// // Restrictions 24356f97ff0aSJeremy L Thompson /// let mut indx: Vec<i32> = vec![0; 2 * ne]; 24366f97ff0aSJeremy L Thompson /// for i in 0..ne { 24376f97ff0aSJeremy L Thompson /// indx[2 * i + 0] = i as i32; 24386f97ff0aSJeremy L Thompson /// indx[2 * i + 1] = (i + 1) as i32; 24396f97ff0aSJeremy L Thompson /// } 24406f97ff0aSJeremy L Thompson /// let rx = ceed.elem_restriction(ne, 2, 1, 1, ne + 1, MemType::Host, &indx)?; 24416f97ff0aSJeremy L Thompson /// let strides: [i32; 3] = [1, q as i32, q as i32]; 24426f97ff0aSJeremy L Thompson /// let rq = ceed.strided_elem_restriction(ne, q, 1, q * ne, strides)?; 24436f97ff0aSJeremy L Thompson /// 24446f97ff0aSJeremy L Thompson /// // Bases 24456f97ff0aSJeremy L Thompson /// let bx = ceed.basis_tensor_H1_Lagrange(1, 1, 2, q, QuadMode::Gauss)?; 24466f97ff0aSJeremy L Thompson /// 24476f97ff0aSJeremy L Thompson /// // Build quadrature data 24486f97ff0aSJeremy L Thompson /// let qf_build_mass = ceed.q_function_interior_by_name("Mass1DBuild")?; 24496f97ff0aSJeremy L Thompson /// let op_build_mass = ceed 24506f97ff0aSJeremy L Thompson /// .operator(&qf_build_mass, QFunctionOpt::None, QFunctionOpt::None)? 24516f97ff0aSJeremy L Thompson /// .field("dx", &rx, &bx, VectorOpt::Active)? 24526f97ff0aSJeremy L Thompson /// .field("weights", ElemRestrictionOpt::None, &bx, VectorOpt::None)? 2453356036faSJeremy L Thompson /// .field("qdata", &rq, BasisOpt::None, VectorOpt::Active)?; 24546f97ff0aSJeremy L Thompson /// 24556f97ff0aSJeremy L Thompson /// let qf_build_diff = ceed.q_function_interior_by_name("Poisson1DBuild")?; 24566f97ff0aSJeremy L Thompson /// let op_build_diff = ceed 24576f97ff0aSJeremy L Thompson /// .operator(&qf_build_diff, QFunctionOpt::None, QFunctionOpt::None)? 24586f97ff0aSJeremy L Thompson /// .field("dx", &rx, &bx, VectorOpt::Active)? 24596f97ff0aSJeremy L Thompson /// .field("weights", ElemRestrictionOpt::None, &bx, VectorOpt::None)? 2460356036faSJeremy L Thompson /// .field("qdata", &rq, BasisOpt::None, VectorOpt::Active)?; 24616f97ff0aSJeremy L Thompson /// 24626f97ff0aSJeremy L Thompson /// let op_build = ceed 24636f97ff0aSJeremy L Thompson /// .composite_operator()? 24646f97ff0aSJeremy L Thompson /// .sub_operator(&op_build_mass)? 24656f97ff0aSJeremy L Thompson /// .sub_operator(&op_build_diff)?; 24666f97ff0aSJeremy L Thompson /// 24676f97ff0aSJeremy L Thompson /// // Check 24686f97ff0aSJeremy L Thompson /// assert!(op_build.check().is_ok(), "Operator setup incorrectly"); 24696f97ff0aSJeremy L Thompson /// # Ok(()) 24706f97ff0aSJeremy L Thompson /// # } 24716f97ff0aSJeremy L Thompson /// ``` 24726f97ff0aSJeremy L Thompson pub fn check(self) -> crate::Result<Self> { 24736f97ff0aSJeremy L Thompson self.op_core.check()?; 24746f97ff0aSJeremy L Thompson Ok(self) 24756f97ff0aSJeremy L Thompson } 24766f97ff0aSJeremy L Thompson 24779df49d7eSJed Brown /// Assemble the diagonal of a square linear Operator 24789df49d7eSJed Brown /// 24799df49d7eSJed Brown /// This overwrites a Vector with the diagonal of a linear Operator. 24809df49d7eSJed Brown /// 24819df49d7eSJed Brown /// Note: Currently only non-composite Operators with a single field and 24829df49d7eSJed Brown /// composite Operators with single field sub-operators are supported. 24839df49d7eSJed Brown /// 24849df49d7eSJed Brown /// * `op` - Operator to assemble QFunction 24859df49d7eSJed Brown /// * `assembled` - Vector to store assembled Operator diagonal 2486b748b478SJeremy L Thompson pub fn linear_assemble_diagonal(&self, assembled: &mut Vector) -> crate::Result<i32> { 24879df49d7eSJed Brown self.op_core.linear_assemble_add_diagonal(assembled) 24889df49d7eSJed Brown } 24899df49d7eSJed Brown 24909df49d7eSJed Brown /// Assemble the point block diagonal of a square linear Operator 24919df49d7eSJed Brown /// 24929df49d7eSJed Brown /// This overwrites a Vector with the point block diagonal of a linear 24939df49d7eSJed Brown /// Operator. 24949df49d7eSJed Brown /// 24959df49d7eSJed Brown /// Note: Currently only non-composite Operators with a single field and 24969df49d7eSJed Brown /// composite Operators with single field sub-operators are supported. 24979df49d7eSJed Brown /// 24989df49d7eSJed Brown /// * `op` - Operator to assemble QFunction 24999df49d7eSJed Brown /// * `assembled` - Vector to store assembled Operator diagonal 25009df49d7eSJed Brown pub fn linear_assemble_add_diagonal(&self, assembled: &mut Vector) -> crate::Result<i32> { 25019df49d7eSJed Brown self.op_core.linear_assemble_add_diagonal(assembled) 25029df49d7eSJed Brown } 25039df49d7eSJed Brown 25049df49d7eSJed Brown /// Assemble the diagonal of a square linear Operator 25059df49d7eSJed Brown /// 25069df49d7eSJed Brown /// This overwrites a Vector with the diagonal of a linear Operator. 25079df49d7eSJed Brown /// 25089df49d7eSJed Brown /// Note: Currently only non-composite Operators with a single field and 25099df49d7eSJed Brown /// composite Operators with single field sub-operators are supported. 25109df49d7eSJed Brown /// 25119df49d7eSJed Brown /// * `op` - Operator to assemble QFunction 25129df49d7eSJed Brown /// * `assembled` - Vector to store assembled CeedOperator point block 25139df49d7eSJed Brown /// diagonal, provided in row-major form with an 25149df49d7eSJed Brown /// `ncomp * ncomp` block at each node. The dimensions of 25159df49d7eSJed Brown /// this vector are derived from the active vector for 25169df49d7eSJed Brown /// the CeedOperator. The array has shape 25179df49d7eSJed Brown /// `[nodes, component out, component in]`. 25189df49d7eSJed Brown pub fn linear_assemble_point_block_diagonal( 25199df49d7eSJed Brown &self, 25209df49d7eSJed Brown assembled: &mut Vector, 25219df49d7eSJed Brown ) -> crate::Result<i32> { 25229df49d7eSJed Brown self.op_core.linear_assemble_add_diagonal(assembled) 25239df49d7eSJed Brown } 25249df49d7eSJed Brown 25259df49d7eSJed Brown /// Assemble the diagonal of a square linear Operator 25269df49d7eSJed Brown /// 25279df49d7eSJed Brown /// This sums into a Vector with the diagonal of a linear Operator. 25289df49d7eSJed Brown /// 25299df49d7eSJed Brown /// Note: Currently only non-composite Operators with a single field and 25309df49d7eSJed Brown /// composite Operators with single field sub-operators are supported. 25319df49d7eSJed Brown /// 25329df49d7eSJed Brown /// * `op` - Operator to assemble QFunction 25339df49d7eSJed Brown /// * `assembled` - Vector to store assembled CeedOperator point block 25349df49d7eSJed Brown /// diagonal, provided in row-major form with an 25359df49d7eSJed Brown /// `ncomp * ncomp` block at each node. The dimensions of 25369df49d7eSJed Brown /// this vector are derived from the active vector for 25379df49d7eSJed Brown /// the CeedOperator. The array has shape 25389df49d7eSJed Brown /// `[nodes, component out, component in]`. 25399df49d7eSJed Brown pub fn linear_assemble_add_point_block_diagonal( 25409df49d7eSJed Brown &self, 25419df49d7eSJed Brown assembled: &mut Vector, 25429df49d7eSJed Brown ) -> crate::Result<i32> { 25439df49d7eSJed Brown self.op_core.linear_assemble_add_diagonal(assembled) 25449df49d7eSJed Brown } 25459df49d7eSJed Brown } 25469df49d7eSJed Brown 25479df49d7eSJed Brown // ----------------------------------------------------------------------------- 2548