13d8e8822SJeremy L Thompson // Copyright (c) 2017-2022, 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, 2008778c6fSJeremy L Thompson _lifeline: PhantomData<&'a ()>, 2108778c6fSJeremy L Thompson } 2208778c6fSJeremy L Thompson 2308778c6fSJeremy L Thompson // ----------------------------------------------------------------------------- 2408778c6fSJeremy L Thompson // Implementations 2508778c6fSJeremy L Thompson // ----------------------------------------------------------------------------- 2608778c6fSJeremy L Thompson impl<'a> OperatorField<'a> { 2708778c6fSJeremy L Thompson /// Get the name of an OperatorField 2808778c6fSJeremy L Thompson /// 2908778c6fSJeremy L Thompson /// ``` 3008778c6fSJeremy L Thompson /// # use libceed::prelude::*; 314d27c890SJeremy L Thompson /// # fn main() -> libceed::Result<()> { 3208778c6fSJeremy L Thompson /// # let ceed = libceed::Ceed::default_init(); 3308778c6fSJeremy L Thompson /// let qf = ceed.q_function_interior_by_name("Mass1DBuild")?; 3408778c6fSJeremy L Thompson /// 3508778c6fSJeremy L Thompson /// // Operator field arguments 3608778c6fSJeremy L Thompson /// let ne = 3; 3708778c6fSJeremy L Thompson /// let q = 4 as usize; 3808778c6fSJeremy L Thompson /// let mut ind: Vec<i32> = vec![0; 2 * ne]; 3908778c6fSJeremy L Thompson /// for i in 0..ne { 4008778c6fSJeremy L Thompson /// ind[2 * i + 0] = i as i32; 4108778c6fSJeremy L Thompson /// ind[2 * i + 1] = (i + 1) as i32; 4208778c6fSJeremy L Thompson /// } 4308778c6fSJeremy L Thompson /// let r = ceed.elem_restriction(ne, 2, 1, 1, ne + 1, MemType::Host, &ind)?; 4408778c6fSJeremy L Thompson /// let strides: [i32; 3] = [1, q as i32, q as i32]; 45d9267b76SJeremy L Thompson /// let rq = ceed.strided_elem_restriction(ne, q, 1, q * ne, strides)?; 4608778c6fSJeremy L Thompson /// 4708778c6fSJeremy L Thompson /// let b = ceed.basis_tensor_H1_Lagrange(1, 1, 2, q, QuadMode::Gauss)?; 4808778c6fSJeremy L Thompson /// 4908778c6fSJeremy L Thompson /// // Operator fields 5008778c6fSJeremy L Thompson /// let op = ceed 5108778c6fSJeremy L Thompson /// .operator(&qf, QFunctionOpt::None, QFunctionOpt::None)? 5208778c6fSJeremy L Thompson /// .field("dx", &r, &b, VectorOpt::Active)? 5308778c6fSJeremy L Thompson /// .field("weights", ElemRestrictionOpt::None, &b, VectorOpt::None)? 54356036faSJeremy L Thompson /// .field("qdata", &rq, BasisOpt::None, VectorOpt::Active)?; 5508778c6fSJeremy L Thompson /// 5608778c6fSJeremy L Thompson /// let inputs = op.inputs()?; 5708778c6fSJeremy L Thompson /// 5808778c6fSJeremy L Thompson /// assert_eq!(inputs[0].name(), "dx", "Incorrect input name"); 5908778c6fSJeremy L Thompson /// assert_eq!(inputs[1].name(), "weights", "Incorrect input name"); 6008778c6fSJeremy L Thompson /// # Ok(()) 6108778c6fSJeremy L Thompson /// # } 6208778c6fSJeremy L Thompson /// ``` 6308778c6fSJeremy L Thompson pub fn name(&self) -> &str { 6408778c6fSJeremy L Thompson let mut name_ptr: *mut std::os::raw::c_char = std::ptr::null_mut(); 6508778c6fSJeremy L Thompson unsafe { 66*6f8994e9SJeremy L Thompson bind_ceed::CeedOperatorFieldGetName( 67*6f8994e9SJeremy L Thompson self.ptr, 68*6f8994e9SJeremy L Thompson &mut name_ptr as *const _ as *mut *const _, 69*6f8994e9SJeremy L Thompson ); 7008778c6fSJeremy L Thompson } 7108778c6fSJeremy L Thompson unsafe { CStr::from_ptr(name_ptr) }.to_str().unwrap() 7208778c6fSJeremy L Thompson } 7308778c6fSJeremy L Thompson 7408778c6fSJeremy L Thompson /// Get the ElemRestriction of an OperatorField 7508778c6fSJeremy L Thompson /// 7608778c6fSJeremy L Thompson /// ``` 7708778c6fSJeremy L Thompson /// # use libceed::prelude::*; 784d27c890SJeremy L Thompson /// # fn main() -> libceed::Result<()> { 7908778c6fSJeremy L Thompson /// # let ceed = libceed::Ceed::default_init(); 8008778c6fSJeremy L Thompson /// let qf = ceed.q_function_interior_by_name("Mass1DBuild")?; 8108778c6fSJeremy L Thompson /// 8208778c6fSJeremy L Thompson /// // Operator field arguments 8308778c6fSJeremy L Thompson /// let ne = 3; 8408778c6fSJeremy L Thompson /// let q = 4 as usize; 8508778c6fSJeremy L Thompson /// let mut ind: Vec<i32> = vec![0; 2 * ne]; 8608778c6fSJeremy L Thompson /// for i in 0..ne { 8708778c6fSJeremy L Thompson /// ind[2 * i + 0] = i as i32; 8808778c6fSJeremy L Thompson /// ind[2 * i + 1] = (i + 1) as i32; 8908778c6fSJeremy L Thompson /// } 9008778c6fSJeremy L Thompson /// let r = ceed.elem_restriction(ne, 2, 1, 1, ne + 1, MemType::Host, &ind)?; 9108778c6fSJeremy L Thompson /// let strides: [i32; 3] = [1, q as i32, q as i32]; 92d9267b76SJeremy L Thompson /// let rq = ceed.strided_elem_restriction(ne, q, 1, q * ne, strides)?; 9308778c6fSJeremy L Thompson /// 9408778c6fSJeremy L Thompson /// let b = ceed.basis_tensor_H1_Lagrange(1, 1, 2, q, QuadMode::Gauss)?; 9508778c6fSJeremy L Thompson /// 9608778c6fSJeremy L Thompson /// // Operator fields 9708778c6fSJeremy L Thompson /// let op = ceed 9808778c6fSJeremy L Thompson /// .operator(&qf, QFunctionOpt::None, QFunctionOpt::None)? 9908778c6fSJeremy L Thompson /// .field("dx", &r, &b, VectorOpt::Active)? 10008778c6fSJeremy L Thompson /// .field("weights", ElemRestrictionOpt::None, &b, VectorOpt::None)? 101356036faSJeremy L Thompson /// .field("qdata", &rq, BasisOpt::None, VectorOpt::Active)?; 10208778c6fSJeremy L Thompson /// 10308778c6fSJeremy L Thompson /// let inputs = op.inputs()?; 10408778c6fSJeremy L Thompson /// 10508778c6fSJeremy L Thompson /// assert!( 10608778c6fSJeremy L Thompson /// inputs[0].elem_restriction().is_some(), 10708778c6fSJeremy L Thompson /// "Incorrect field ElemRestriction" 10808778c6fSJeremy L Thompson /// ); 10908778c6fSJeremy L Thompson /// assert!( 11008778c6fSJeremy L Thompson /// inputs[1].elem_restriction().is_none(), 11108778c6fSJeremy L Thompson /// "Incorrect field ElemRestriction" 11208778c6fSJeremy L Thompson /// ); 11308778c6fSJeremy L Thompson /// # Ok(()) 11408778c6fSJeremy L Thompson /// # } 11508778c6fSJeremy L Thompson /// ``` 11608778c6fSJeremy L Thompson pub fn elem_restriction(&self) -> ElemRestrictionOpt { 11708778c6fSJeremy L Thompson let mut ptr = std::ptr::null_mut(); 11808778c6fSJeremy L Thompson unsafe { 11908778c6fSJeremy L Thompson bind_ceed::CeedOperatorFieldGetElemRestriction(self.ptr, &mut ptr); 12008778c6fSJeremy L Thompson } 12108778c6fSJeremy L Thompson if ptr == unsafe { bind_ceed::CEED_ELEMRESTRICTION_NONE } { 12208778c6fSJeremy L Thompson ElemRestrictionOpt::None 12308778c6fSJeremy L Thompson } else { 12408778c6fSJeremy L Thompson let slice = unsafe { 12508778c6fSJeremy L Thompson std::slice::from_raw_parts( 12608778c6fSJeremy L Thompson &ptr as *const bind_ceed::CeedElemRestriction as *const crate::ElemRestriction, 12708778c6fSJeremy L Thompson 1 as usize, 12808778c6fSJeremy L Thompson ) 12908778c6fSJeremy L Thompson }; 13008778c6fSJeremy L Thompson ElemRestrictionOpt::Some(&slice[0]) 13108778c6fSJeremy L Thompson } 13208778c6fSJeremy L Thompson } 13308778c6fSJeremy L Thompson 13408778c6fSJeremy L Thompson /// Get the Basis of an OperatorField 13508778c6fSJeremy L Thompson /// 13608778c6fSJeremy L Thompson /// ``` 13708778c6fSJeremy L Thompson /// # use libceed::prelude::*; 1384d27c890SJeremy L Thompson /// # fn main() -> libceed::Result<()> { 13908778c6fSJeremy L Thompson /// # let ceed = libceed::Ceed::default_init(); 14008778c6fSJeremy L Thompson /// let qf = ceed.q_function_interior_by_name("Mass1DBuild")?; 14108778c6fSJeremy L Thompson /// 14208778c6fSJeremy L Thompson /// // Operator field arguments 14308778c6fSJeremy L Thompson /// let ne = 3; 14408778c6fSJeremy L Thompson /// let q = 4 as usize; 14508778c6fSJeremy L Thompson /// let mut ind: Vec<i32> = vec![0; 2 * ne]; 14608778c6fSJeremy L Thompson /// for i in 0..ne { 14708778c6fSJeremy L Thompson /// ind[2 * i + 0] = i as i32; 14808778c6fSJeremy L Thompson /// ind[2 * i + 1] = (i + 1) as i32; 14908778c6fSJeremy L Thompson /// } 15008778c6fSJeremy L Thompson /// let r = ceed.elem_restriction(ne, 2, 1, 1, ne + 1, MemType::Host, &ind)?; 15108778c6fSJeremy L Thompson /// let strides: [i32; 3] = [1, q as i32, q as i32]; 152d9267b76SJeremy L Thompson /// let rq = ceed.strided_elem_restriction(ne, q, 1, q * ne, strides)?; 15308778c6fSJeremy L Thompson /// 15408778c6fSJeremy L Thompson /// let b = ceed.basis_tensor_H1_Lagrange(1, 1, 2, q, QuadMode::Gauss)?; 15508778c6fSJeremy L Thompson /// 15608778c6fSJeremy L Thompson /// // Operator fields 15708778c6fSJeremy L Thompson /// let op = ceed 15808778c6fSJeremy L Thompson /// .operator(&qf, QFunctionOpt::None, QFunctionOpt::None)? 15908778c6fSJeremy L Thompson /// .field("dx", &r, &b, VectorOpt::Active)? 16008778c6fSJeremy L Thompson /// .field("weights", ElemRestrictionOpt::None, &b, VectorOpt::None)? 161356036faSJeremy L Thompson /// .field("qdata", &rq, BasisOpt::None, VectorOpt::Active)?; 16208778c6fSJeremy L Thompson /// 16308778c6fSJeremy L Thompson /// let inputs = op.inputs()?; 16408778c6fSJeremy L Thompson /// 16508778c6fSJeremy L Thompson /// assert!(inputs[0].basis().is_some(), "Incorrect field Basis"); 16608778c6fSJeremy L Thompson /// assert!(inputs[1].basis().is_some(), "Incorrect field Basis"); 16708778c6fSJeremy L Thompson /// 16808778c6fSJeremy L Thompson /// let outputs = op.outputs()?; 16908778c6fSJeremy L Thompson /// 170356036faSJeremy L Thompson /// assert!(outputs[0].basis().is_none(), "Incorrect field Basis"); 17108778c6fSJeremy L Thompson /// # Ok(()) 17208778c6fSJeremy L Thompson /// # } 17308778c6fSJeremy L Thompson /// ``` 17408778c6fSJeremy L Thompson pub fn basis(&self) -> BasisOpt { 17508778c6fSJeremy L Thompson let mut ptr = std::ptr::null_mut(); 17608778c6fSJeremy L Thompson unsafe { 17708778c6fSJeremy L Thompson bind_ceed::CeedOperatorFieldGetBasis(self.ptr, &mut ptr); 17808778c6fSJeremy L Thompson } 179356036faSJeremy L Thompson if ptr == unsafe { bind_ceed::CEED_BASIS_NONE } { 180356036faSJeremy L Thompson BasisOpt::None 18108778c6fSJeremy L Thompson } else { 18208778c6fSJeremy L Thompson let slice = unsafe { 18308778c6fSJeremy L Thompson std::slice::from_raw_parts( 18408778c6fSJeremy L Thompson &ptr as *const bind_ceed::CeedBasis as *const crate::Basis, 18508778c6fSJeremy L Thompson 1 as usize, 18608778c6fSJeremy L Thompson ) 18708778c6fSJeremy L Thompson }; 18808778c6fSJeremy L Thompson BasisOpt::Some(&slice[0]) 18908778c6fSJeremy L Thompson } 19008778c6fSJeremy L Thompson } 19108778c6fSJeremy L Thompson 19208778c6fSJeremy L Thompson /// Get the Vector of an OperatorField 19308778c6fSJeremy L Thompson /// 19408778c6fSJeremy L Thompson /// ``` 19508778c6fSJeremy L Thompson /// # use libceed::prelude::*; 1964d27c890SJeremy L Thompson /// # fn main() -> libceed::Result<()> { 19708778c6fSJeremy L Thompson /// # let ceed = libceed::Ceed::default_init(); 19808778c6fSJeremy L Thompson /// let qf = ceed.q_function_interior_by_name("Mass1DBuild")?; 19908778c6fSJeremy L Thompson /// 20008778c6fSJeremy L Thompson /// // Operator field arguments 20108778c6fSJeremy L Thompson /// let ne = 3; 20208778c6fSJeremy L Thompson /// let q = 4 as usize; 20308778c6fSJeremy L Thompson /// let mut ind: Vec<i32> = vec![0; 2 * ne]; 20408778c6fSJeremy L Thompson /// for i in 0..ne { 20508778c6fSJeremy L Thompson /// ind[2 * i + 0] = i as i32; 20608778c6fSJeremy L Thompson /// ind[2 * i + 1] = (i + 1) as i32; 20708778c6fSJeremy L Thompson /// } 20808778c6fSJeremy L Thompson /// let r = ceed.elem_restriction(ne, 2, 1, 1, ne + 1, MemType::Host, &ind)?; 20908778c6fSJeremy L Thompson /// let strides: [i32; 3] = [1, q as i32, q as i32]; 210d9267b76SJeremy L Thompson /// let rq = ceed.strided_elem_restriction(ne, q, 1, q * ne, strides)?; 21108778c6fSJeremy L Thompson /// 21208778c6fSJeremy L Thompson /// let b = ceed.basis_tensor_H1_Lagrange(1, 1, 2, q, QuadMode::Gauss)?; 21308778c6fSJeremy L Thompson /// 21408778c6fSJeremy L Thompson /// // Operator fields 21508778c6fSJeremy L Thompson /// let op = ceed 21608778c6fSJeremy L Thompson /// .operator(&qf, QFunctionOpt::None, QFunctionOpt::None)? 21708778c6fSJeremy L Thompson /// .field("dx", &r, &b, VectorOpt::Active)? 21808778c6fSJeremy L Thompson /// .field("weights", ElemRestrictionOpt::None, &b, VectorOpt::None)? 219356036faSJeremy L Thompson /// .field("qdata", &rq, BasisOpt::None, VectorOpt::Active)?; 22008778c6fSJeremy L Thompson /// 22108778c6fSJeremy L Thompson /// let inputs = op.inputs()?; 22208778c6fSJeremy L Thompson /// 22308778c6fSJeremy L Thompson /// assert!(inputs[0].vector().is_active(), "Incorrect field Vector"); 22408778c6fSJeremy L Thompson /// assert!(inputs[1].vector().is_none(), "Incorrect field Vector"); 22508778c6fSJeremy L Thompson /// # Ok(()) 22608778c6fSJeremy L Thompson /// # } 22708778c6fSJeremy L Thompson /// ``` 22808778c6fSJeremy L Thompson pub fn vector(&self) -> VectorOpt { 22908778c6fSJeremy L Thompson let mut ptr = std::ptr::null_mut(); 23008778c6fSJeremy L Thompson unsafe { 23108778c6fSJeremy L Thompson bind_ceed::CeedOperatorFieldGetVector(self.ptr, &mut ptr); 23208778c6fSJeremy L Thompson } 23308778c6fSJeremy L Thompson if ptr == unsafe { bind_ceed::CEED_VECTOR_ACTIVE } { 23408778c6fSJeremy L Thompson VectorOpt::Active 23508778c6fSJeremy L Thompson } else if ptr == unsafe { bind_ceed::CEED_VECTOR_NONE } { 23608778c6fSJeremy L Thompson VectorOpt::None 23708778c6fSJeremy L Thompson } else { 23808778c6fSJeremy L Thompson let slice = unsafe { 23908778c6fSJeremy L Thompson std::slice::from_raw_parts( 24008778c6fSJeremy L Thompson &ptr as *const bind_ceed::CeedVector as *const crate::Vector, 24108778c6fSJeremy L Thompson 1 as usize, 24208778c6fSJeremy L Thompson ) 24308778c6fSJeremy L Thompson }; 24408778c6fSJeremy L Thompson VectorOpt::Some(&slice[0]) 24508778c6fSJeremy L Thompson } 24608778c6fSJeremy L Thompson } 24708778c6fSJeremy L Thompson } 24808778c6fSJeremy L Thompson 24908778c6fSJeremy L Thompson // ----------------------------------------------------------------------------- 2507ed177dbSJed Brown // Operator context wrapper 2519df49d7eSJed Brown // ----------------------------------------------------------------------------- 252c68be7a2SJeremy L Thompson #[derive(Debug)] 2539df49d7eSJed Brown pub(crate) struct OperatorCore<'a> { 2549df49d7eSJed Brown ptr: bind_ceed::CeedOperator, 2551142270cSJeremy L Thompson _lifeline: PhantomData<&'a ()>, 2569df49d7eSJed Brown } 2579df49d7eSJed Brown 258c68be7a2SJeremy L Thompson #[derive(Debug)] 2599df49d7eSJed Brown pub struct Operator<'a> { 2609df49d7eSJed Brown op_core: OperatorCore<'a>, 2619df49d7eSJed Brown } 2629df49d7eSJed Brown 263c68be7a2SJeremy L Thompson #[derive(Debug)] 2649df49d7eSJed Brown pub struct CompositeOperator<'a> { 2659df49d7eSJed Brown op_core: OperatorCore<'a>, 2669df49d7eSJed Brown } 2679df49d7eSJed Brown 2689df49d7eSJed Brown // ----------------------------------------------------------------------------- 2699df49d7eSJed Brown // Destructor 2709df49d7eSJed Brown // ----------------------------------------------------------------------------- 2719df49d7eSJed Brown impl<'a> Drop for OperatorCore<'a> { 2729df49d7eSJed Brown fn drop(&mut self) { 2739df49d7eSJed Brown unsafe { 2749df49d7eSJed Brown bind_ceed::CeedOperatorDestroy(&mut self.ptr); 2759df49d7eSJed Brown } 2769df49d7eSJed Brown } 2779df49d7eSJed Brown } 2789df49d7eSJed Brown 2799df49d7eSJed Brown // ----------------------------------------------------------------------------- 2809df49d7eSJed Brown // Display 2819df49d7eSJed Brown // ----------------------------------------------------------------------------- 2829df49d7eSJed Brown impl<'a> fmt::Display for OperatorCore<'a> { 2839df49d7eSJed Brown fn fmt(&self, f: &mut fmt::Formatter) -> fmt::Result { 2849df49d7eSJed Brown let mut ptr = std::ptr::null_mut(); 2859df49d7eSJed Brown let mut sizeloc = crate::MAX_BUFFER_LENGTH; 2869df49d7eSJed Brown let cstring = unsafe { 2879df49d7eSJed Brown let file = bind_ceed::open_memstream(&mut ptr, &mut sizeloc); 2889df49d7eSJed Brown bind_ceed::CeedOperatorView(self.ptr, file); 2899df49d7eSJed Brown bind_ceed::fclose(file); 2909df49d7eSJed Brown CString::from_raw(ptr) 2919df49d7eSJed Brown }; 2929df49d7eSJed Brown cstring.to_string_lossy().fmt(f) 2939df49d7eSJed Brown } 2949df49d7eSJed Brown } 2959df49d7eSJed Brown 2969df49d7eSJed Brown /// View an Operator 2979df49d7eSJed Brown /// 2989df49d7eSJed Brown /// ``` 2999df49d7eSJed Brown /// # use libceed::prelude::*; 3004d27c890SJeremy L Thompson /// # fn main() -> libceed::Result<()> { 3019df49d7eSJed Brown /// # let ceed = libceed::Ceed::default_init(); 302c68be7a2SJeremy L Thompson /// let qf = ceed.q_function_interior_by_name("Mass1DBuild")?; 3039df49d7eSJed Brown /// 3049df49d7eSJed Brown /// // Operator field arguments 3059df49d7eSJed Brown /// let ne = 3; 3069df49d7eSJed Brown /// let q = 4 as usize; 3079df49d7eSJed Brown /// let mut ind: Vec<i32> = vec![0; 2 * ne]; 3089df49d7eSJed Brown /// for i in 0..ne { 3099df49d7eSJed Brown /// ind[2 * i + 0] = i as i32; 3109df49d7eSJed Brown /// ind[2 * i + 1] = (i + 1) as i32; 3119df49d7eSJed Brown /// } 312c68be7a2SJeremy L Thompson /// let r = ceed.elem_restriction(ne, 2, 1, 1, ne + 1, MemType::Host, &ind)?; 3139df49d7eSJed Brown /// let strides: [i32; 3] = [1, q as i32, q as i32]; 314d9267b76SJeremy L Thompson /// let rq = ceed.strided_elem_restriction(ne, q, 1, q * ne, strides)?; 3159df49d7eSJed Brown /// 316c68be7a2SJeremy L Thompson /// let b = ceed.basis_tensor_H1_Lagrange(1, 1, 2, q, QuadMode::Gauss)?; 3179df49d7eSJed Brown /// 3189df49d7eSJed Brown /// // Operator fields 3199df49d7eSJed Brown /// let op = ceed 320c68be7a2SJeremy L Thompson /// .operator(&qf, QFunctionOpt::None, QFunctionOpt::None)? 321ea6b5821SJeremy L Thompson /// .name("mass")? 322c68be7a2SJeremy L Thompson /// .field("dx", &r, &b, VectorOpt::Active)? 323c68be7a2SJeremy L Thompson /// .field("weights", ElemRestrictionOpt::None, &b, VectorOpt::None)? 324356036faSJeremy L Thompson /// .field("qdata", &rq, BasisOpt::None, VectorOpt::Active)?; 3259df49d7eSJed Brown /// 3269df49d7eSJed Brown /// println!("{}", op); 327c68be7a2SJeremy L Thompson /// # Ok(()) 328c68be7a2SJeremy L Thompson /// # } 3299df49d7eSJed Brown /// ``` 3309df49d7eSJed Brown impl<'a> fmt::Display for Operator<'a> { 3319df49d7eSJed Brown fn fmt(&self, f: &mut fmt::Formatter) -> fmt::Result { 3329df49d7eSJed Brown self.op_core.fmt(f) 3339df49d7eSJed Brown } 3349df49d7eSJed Brown } 3359df49d7eSJed Brown 3369df49d7eSJed Brown /// View a composite Operator 3379df49d7eSJed Brown /// 3389df49d7eSJed Brown /// ``` 3399df49d7eSJed Brown /// # use libceed::prelude::*; 3404d27c890SJeremy L Thompson /// # fn main() -> libceed::Result<()> { 3419df49d7eSJed Brown /// # let ceed = libceed::Ceed::default_init(); 3429df49d7eSJed Brown /// 3439df49d7eSJed Brown /// // Sub operator field arguments 3449df49d7eSJed Brown /// let ne = 3; 3459df49d7eSJed Brown /// let q = 4 as usize; 3469df49d7eSJed Brown /// let mut ind: Vec<i32> = vec![0; 2 * ne]; 3479df49d7eSJed Brown /// for i in 0..ne { 3489df49d7eSJed Brown /// ind[2 * i + 0] = i as i32; 3499df49d7eSJed Brown /// ind[2 * i + 1] = (i + 1) as i32; 3509df49d7eSJed Brown /// } 351c68be7a2SJeremy L Thompson /// let r = ceed.elem_restriction(ne, 2, 1, 1, ne + 1, MemType::Host, &ind)?; 3529df49d7eSJed Brown /// let strides: [i32; 3] = [1, q as i32, q as i32]; 353d9267b76SJeremy L Thompson /// let rq = ceed.strided_elem_restriction(ne, q, 1, q * ne, strides)?; 3549df49d7eSJed Brown /// 355c68be7a2SJeremy L Thompson /// let b = ceed.basis_tensor_H1_Lagrange(1, 1, 2, q, QuadMode::Gauss)?; 3569df49d7eSJed Brown /// 357c68be7a2SJeremy L Thompson /// let qdata_mass = ceed.vector(q * ne)?; 358c68be7a2SJeremy L Thompson /// let qdata_diff = ceed.vector(q * ne)?; 3599df49d7eSJed Brown /// 360c68be7a2SJeremy L Thompson /// let qf_mass = ceed.q_function_interior_by_name("MassApply")?; 3619df49d7eSJed Brown /// let op_mass = ceed 362c68be7a2SJeremy L Thompson /// .operator(&qf_mass, QFunctionOpt::None, QFunctionOpt::None)? 363ea6b5821SJeremy L Thompson /// .name("Mass term")? 364c68be7a2SJeremy L Thompson /// .field("u", &r, &b, VectorOpt::Active)? 365356036faSJeremy L Thompson /// .field("qdata", &rq, BasisOpt::None, &qdata_mass)? 366c68be7a2SJeremy L Thompson /// .field("v", &r, &b, VectorOpt::Active)?; 3679df49d7eSJed Brown /// 368c68be7a2SJeremy L Thompson /// let qf_diff = ceed.q_function_interior_by_name("Poisson1DApply")?; 3699df49d7eSJed Brown /// let op_diff = ceed 370c68be7a2SJeremy L Thompson /// .operator(&qf_diff, QFunctionOpt::None, QFunctionOpt::None)? 371ea6b5821SJeremy L Thompson /// .name("Poisson term")? 372c68be7a2SJeremy L Thompson /// .field("du", &r, &b, VectorOpt::Active)? 373356036faSJeremy L Thompson /// .field("qdata", &rq, BasisOpt::None, &qdata_diff)? 374c68be7a2SJeremy L Thompson /// .field("dv", &r, &b, VectorOpt::Active)?; 3759df49d7eSJed Brown /// 3769df49d7eSJed Brown /// let op = ceed 377c68be7a2SJeremy L Thompson /// .composite_operator()? 378ea6b5821SJeremy L Thompson /// .name("Screened Poisson")? 379c68be7a2SJeremy L Thompson /// .sub_operator(&op_mass)? 380c68be7a2SJeremy L Thompson /// .sub_operator(&op_diff)?; 3819df49d7eSJed Brown /// 3829df49d7eSJed Brown /// println!("{}", op); 383c68be7a2SJeremy L Thompson /// # Ok(()) 384c68be7a2SJeremy L Thompson /// # } 3859df49d7eSJed Brown /// ``` 3869df49d7eSJed Brown impl<'a> fmt::Display for CompositeOperator<'a> { 3879df49d7eSJed Brown fn fmt(&self, f: &mut fmt::Formatter) -> fmt::Result { 3889df49d7eSJed Brown self.op_core.fmt(f) 3899df49d7eSJed Brown } 3909df49d7eSJed Brown } 3919df49d7eSJed Brown 3929df49d7eSJed Brown // ----------------------------------------------------------------------------- 3939df49d7eSJed Brown // Core functionality 3949df49d7eSJed Brown // ----------------------------------------------------------------------------- 3959df49d7eSJed Brown impl<'a> OperatorCore<'a> { 3961142270cSJeremy L Thompson // Error handling 3971142270cSJeremy L Thompson #[doc(hidden)] 3981142270cSJeremy L Thompson fn check_error(&self, ierr: i32) -> crate::Result<i32> { 3991142270cSJeremy L Thompson let mut ptr = std::ptr::null_mut(); 4001142270cSJeremy L Thompson unsafe { 4011142270cSJeremy L Thompson bind_ceed::CeedOperatorGetCeed(self.ptr, &mut ptr); 4021142270cSJeremy L Thompson } 4031142270cSJeremy L Thompson crate::check_error(ptr, ierr) 4041142270cSJeremy L Thompson } 4051142270cSJeremy L Thompson 4069df49d7eSJed Brown // Common implementations 4076f97ff0aSJeremy L Thompson pub fn check(&self) -> crate::Result<i32> { 4086f97ff0aSJeremy L Thompson let ierr = unsafe { bind_ceed::CeedOperatorCheckReady(self.ptr) }; 4096f97ff0aSJeremy L Thompson self.check_error(ierr) 4106f97ff0aSJeremy L Thompson } 4116f97ff0aSJeremy L Thompson 412ea6b5821SJeremy L Thompson pub fn name(&self, name: &str) -> crate::Result<i32> { 413ea6b5821SJeremy L Thompson let name_c = CString::new(name).expect("CString::new failed"); 414ea6b5821SJeremy L Thompson let ierr = unsafe { bind_ceed::CeedOperatorSetName(self.ptr, name_c.as_ptr()) }; 415ea6b5821SJeremy L Thompson self.check_error(ierr) 416ea6b5821SJeremy L Thompson } 417ea6b5821SJeremy L Thompson 4189df49d7eSJed Brown pub fn apply(&self, input: &Vector, output: &mut Vector) -> crate::Result<i32> { 4199df49d7eSJed Brown let ierr = unsafe { 4209df49d7eSJed Brown bind_ceed::CeedOperatorApply( 4219df49d7eSJed Brown self.ptr, 4229df49d7eSJed Brown input.ptr, 4239df49d7eSJed Brown output.ptr, 4249df49d7eSJed Brown bind_ceed::CEED_REQUEST_IMMEDIATE, 4259df49d7eSJed Brown ) 4269df49d7eSJed Brown }; 4271142270cSJeremy L Thompson self.check_error(ierr) 4289df49d7eSJed Brown } 4299df49d7eSJed Brown 4309df49d7eSJed Brown pub fn apply_add(&self, input: &Vector, output: &mut Vector) -> crate::Result<i32> { 4319df49d7eSJed Brown let ierr = unsafe { 4329df49d7eSJed Brown bind_ceed::CeedOperatorApplyAdd( 4339df49d7eSJed Brown self.ptr, 4349df49d7eSJed Brown input.ptr, 4359df49d7eSJed Brown output.ptr, 4369df49d7eSJed Brown bind_ceed::CEED_REQUEST_IMMEDIATE, 4379df49d7eSJed Brown ) 4389df49d7eSJed Brown }; 4391142270cSJeremy L Thompson self.check_error(ierr) 4409df49d7eSJed Brown } 4419df49d7eSJed Brown 4429df49d7eSJed Brown pub fn linear_assemble_diagonal(&self, assembled: &mut Vector) -> crate::Result<i32> { 4439df49d7eSJed Brown let ierr = unsafe { 4449df49d7eSJed Brown bind_ceed::CeedOperatorLinearAssembleDiagonal( 4459df49d7eSJed Brown self.ptr, 4469df49d7eSJed Brown assembled.ptr, 4479df49d7eSJed Brown bind_ceed::CEED_REQUEST_IMMEDIATE, 4489df49d7eSJed Brown ) 4499df49d7eSJed Brown }; 4501142270cSJeremy L Thompson self.check_error(ierr) 4519df49d7eSJed Brown } 4529df49d7eSJed Brown 4539df49d7eSJed Brown pub fn linear_assemble_add_diagonal(&self, assembled: &mut Vector) -> crate::Result<i32> { 4549df49d7eSJed Brown let ierr = unsafe { 4559df49d7eSJed Brown bind_ceed::CeedOperatorLinearAssembleAddDiagonal( 4569df49d7eSJed Brown self.ptr, 4579df49d7eSJed Brown assembled.ptr, 4589df49d7eSJed Brown bind_ceed::CEED_REQUEST_IMMEDIATE, 4599df49d7eSJed Brown ) 4609df49d7eSJed Brown }; 4611142270cSJeremy L Thompson self.check_error(ierr) 4629df49d7eSJed Brown } 4639df49d7eSJed Brown 4649df49d7eSJed Brown pub fn linear_assemble_point_block_diagonal( 4659df49d7eSJed Brown &self, 4669df49d7eSJed Brown assembled: &mut Vector, 4679df49d7eSJed Brown ) -> crate::Result<i32> { 4689df49d7eSJed Brown let ierr = unsafe { 4699df49d7eSJed Brown bind_ceed::CeedOperatorLinearAssemblePointBlockDiagonal( 4709df49d7eSJed Brown self.ptr, 4719df49d7eSJed Brown assembled.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 linear_assemble_add_point_block_diagonal( 4799df49d7eSJed Brown &self, 4809df49d7eSJed Brown assembled: &mut Vector, 4819df49d7eSJed Brown ) -> crate::Result<i32> { 4829df49d7eSJed Brown let ierr = unsafe { 4839df49d7eSJed Brown bind_ceed::CeedOperatorLinearAssembleAddPointBlockDiagonal( 4849df49d7eSJed Brown self.ptr, 4859df49d7eSJed Brown assembled.ptr, 4869df49d7eSJed Brown bind_ceed::CEED_REQUEST_IMMEDIATE, 4879df49d7eSJed Brown ) 4889df49d7eSJed Brown }; 4891142270cSJeremy L Thompson self.check_error(ierr) 4909df49d7eSJed Brown } 4919df49d7eSJed Brown } 4929df49d7eSJed Brown 4939df49d7eSJed Brown // ----------------------------------------------------------------------------- 4949df49d7eSJed Brown // Operator 4959df49d7eSJed Brown // ----------------------------------------------------------------------------- 4969df49d7eSJed Brown impl<'a> Operator<'a> { 4979df49d7eSJed Brown // Constructor 4989df49d7eSJed Brown pub fn create<'b>( 499594ef120SJeremy L Thompson ceed: &crate::Ceed, 5009df49d7eSJed Brown qf: impl Into<QFunctionOpt<'b>>, 5019df49d7eSJed Brown dqf: impl Into<QFunctionOpt<'b>>, 5029df49d7eSJed Brown dqfT: impl Into<QFunctionOpt<'b>>, 5039df49d7eSJed Brown ) -> crate::Result<Self> { 5049df49d7eSJed Brown let mut ptr = std::ptr::null_mut(); 5059df49d7eSJed Brown let ierr = unsafe { 5069df49d7eSJed Brown bind_ceed::CeedOperatorCreate( 5079df49d7eSJed Brown ceed.ptr, 5089df49d7eSJed Brown qf.into().to_raw(), 5099df49d7eSJed Brown dqf.into().to_raw(), 5109df49d7eSJed Brown dqfT.into().to_raw(), 5119df49d7eSJed Brown &mut ptr, 5129df49d7eSJed Brown ) 5139df49d7eSJed Brown }; 5149df49d7eSJed Brown ceed.check_error(ierr)?; 5159df49d7eSJed Brown Ok(Self { 5161142270cSJeremy L Thompson op_core: OperatorCore { 5171142270cSJeremy L Thompson ptr, 5181142270cSJeremy L Thompson _lifeline: PhantomData, 5191142270cSJeremy L Thompson }, 5209df49d7eSJed Brown }) 5219df49d7eSJed Brown } 5229df49d7eSJed Brown 5231142270cSJeremy L Thompson fn from_raw(ptr: bind_ceed::CeedOperator) -> crate::Result<Self> { 5249df49d7eSJed Brown Ok(Self { 5251142270cSJeremy L Thompson op_core: OperatorCore { 5261142270cSJeremy L Thompson ptr, 5271142270cSJeremy L Thompson _lifeline: PhantomData, 5281142270cSJeremy L Thompson }, 5299df49d7eSJed Brown }) 5309df49d7eSJed Brown } 5319df49d7eSJed Brown 532ea6b5821SJeremy L Thompson /// Set name for Operator printing 533ea6b5821SJeremy L Thompson /// 534ea6b5821SJeremy L Thompson /// * 'name' - Name to set 535ea6b5821SJeremy L Thompson /// 536ea6b5821SJeremy L Thompson /// ``` 537ea6b5821SJeremy L Thompson /// # use libceed::prelude::*; 538ea6b5821SJeremy L Thompson /// # fn main() -> libceed::Result<()> { 539ea6b5821SJeremy L Thompson /// # let ceed = libceed::Ceed::default_init(); 540ea6b5821SJeremy L Thompson /// let qf = ceed.q_function_interior_by_name("Mass1DBuild")?; 541ea6b5821SJeremy L Thompson /// 542ea6b5821SJeremy L Thompson /// // Operator field arguments 543ea6b5821SJeremy L Thompson /// let ne = 3; 544ea6b5821SJeremy L Thompson /// let q = 4 as usize; 545ea6b5821SJeremy L Thompson /// let mut ind: Vec<i32> = vec![0; 2 * ne]; 546ea6b5821SJeremy L Thompson /// for i in 0..ne { 547ea6b5821SJeremy L Thompson /// ind[2 * i + 0] = i as i32; 548ea6b5821SJeremy L Thompson /// ind[2 * i + 1] = (i + 1) as i32; 549ea6b5821SJeremy L Thompson /// } 550ea6b5821SJeremy L Thompson /// let r = ceed.elem_restriction(ne, 2, 1, 1, ne + 1, MemType::Host, &ind)?; 551ea6b5821SJeremy L Thompson /// let strides: [i32; 3] = [1, q as i32, q as i32]; 552d9267b76SJeremy L Thompson /// let rq = ceed.strided_elem_restriction(ne, q, 1, q * ne, strides)?; 553ea6b5821SJeremy L Thompson /// 554ea6b5821SJeremy L Thompson /// let b = ceed.basis_tensor_H1_Lagrange(1, 1, 2, q, QuadMode::Gauss)?; 555ea6b5821SJeremy L Thompson /// 556ea6b5821SJeremy L Thompson /// // Operator fields 557ea6b5821SJeremy L Thompson /// let op = ceed 558ea6b5821SJeremy L Thompson /// .operator(&qf, QFunctionOpt::None, QFunctionOpt::None)? 559ea6b5821SJeremy L Thompson /// .name("mass")? 560ea6b5821SJeremy L Thompson /// .field("dx", &r, &b, VectorOpt::Active)? 561ea6b5821SJeremy L Thompson /// .field("weights", ElemRestrictionOpt::None, &b, VectorOpt::None)? 562356036faSJeremy L Thompson /// .field("qdata", &rq, BasisOpt::None, VectorOpt::Active)?; 563ea6b5821SJeremy L Thompson /// # Ok(()) 564ea6b5821SJeremy L Thompson /// # } 565ea6b5821SJeremy L Thompson /// ``` 566ea6b5821SJeremy L Thompson #[allow(unused_mut)] 567ea6b5821SJeremy L Thompson pub fn name(mut self, name: &str) -> crate::Result<Self> { 568ea6b5821SJeremy L Thompson self.op_core.name(name)?; 569ea6b5821SJeremy L Thompson Ok(self) 570ea6b5821SJeremy L Thompson } 571ea6b5821SJeremy L Thompson 5729df49d7eSJed Brown /// Apply Operator to a vector 5739df49d7eSJed Brown /// 5749df49d7eSJed Brown /// * `input` - Input Vector 5759df49d7eSJed Brown /// * `output` - Output Vector 5769df49d7eSJed Brown /// 5779df49d7eSJed Brown /// ``` 5789df49d7eSJed Brown /// # use libceed::prelude::*; 5794d27c890SJeremy L Thompson /// # fn main() -> libceed::Result<()> { 5809df49d7eSJed Brown /// # let ceed = libceed::Ceed::default_init(); 5819df49d7eSJed Brown /// let ne = 4; 5829df49d7eSJed Brown /// let p = 3; 5839df49d7eSJed Brown /// let q = 4; 5849df49d7eSJed Brown /// let ndofs = p * ne - ne + 1; 5859df49d7eSJed Brown /// 5869df49d7eSJed Brown /// // Vectors 587c68be7a2SJeremy L Thompson /// let x = ceed.vector_from_slice(&[-1., -0.5, 0.0, 0.5, 1.0])?; 588c68be7a2SJeremy L Thompson /// let mut qdata = ceed.vector(ne * q)?; 5899df49d7eSJed Brown /// qdata.set_value(0.0); 590c68be7a2SJeremy L Thompson /// let u = ceed.vector_from_slice(&vec![1.0; ndofs])?; 591c68be7a2SJeremy L Thompson /// let mut v = ceed.vector(ndofs)?; 5929df49d7eSJed Brown /// v.set_value(0.0); 5939df49d7eSJed Brown /// 5949df49d7eSJed Brown /// // Restrictions 5959df49d7eSJed Brown /// let mut indx: Vec<i32> = vec![0; 2 * ne]; 5969df49d7eSJed Brown /// for i in 0..ne { 5979df49d7eSJed Brown /// indx[2 * i + 0] = i as i32; 5989df49d7eSJed Brown /// indx[2 * i + 1] = (i + 1) as i32; 5999df49d7eSJed Brown /// } 600c68be7a2SJeremy L Thompson /// let rx = ceed.elem_restriction(ne, 2, 1, 1, ne + 1, MemType::Host, &indx)?; 6019df49d7eSJed Brown /// let mut indu: Vec<i32> = vec![0; p * ne]; 6029df49d7eSJed Brown /// for i in 0..ne { 6039df49d7eSJed Brown /// indu[p * i + 0] = i as i32; 6049df49d7eSJed Brown /// indu[p * i + 1] = (i + 1) as i32; 6059df49d7eSJed Brown /// indu[p * i + 2] = (i + 2) as i32; 6069df49d7eSJed Brown /// } 607c68be7a2SJeremy L Thompson /// let ru = ceed.elem_restriction(ne, 3, 1, 1, ndofs, MemType::Host, &indu)?; 6089df49d7eSJed Brown /// let strides: [i32; 3] = [1, q as i32, q as i32]; 609c68be7a2SJeremy L Thompson /// let rq = ceed.strided_elem_restriction(ne, q, 1, q * ne, strides)?; 6109df49d7eSJed Brown /// 6119df49d7eSJed Brown /// // Bases 612c68be7a2SJeremy L Thompson /// let bx = ceed.basis_tensor_H1_Lagrange(1, 1, 2, q, QuadMode::Gauss)?; 613c68be7a2SJeremy L Thompson /// let bu = ceed.basis_tensor_H1_Lagrange(1, 1, p, q, QuadMode::Gauss)?; 6149df49d7eSJed Brown /// 6159df49d7eSJed Brown /// // Build quadrature data 616c68be7a2SJeremy L Thompson /// let qf_build = ceed.q_function_interior_by_name("Mass1DBuild")?; 617c68be7a2SJeremy L Thompson /// ceed.operator(&qf_build, QFunctionOpt::None, QFunctionOpt::None)? 618c68be7a2SJeremy L Thompson /// .field("dx", &rx, &bx, VectorOpt::Active)? 619c68be7a2SJeremy L Thompson /// .field("weights", ElemRestrictionOpt::None, &bx, VectorOpt::None)? 620356036faSJeremy L Thompson /// .field("qdata", &rq, BasisOpt::None, VectorOpt::Active)? 621c68be7a2SJeremy L Thompson /// .apply(&x, &mut qdata)?; 6229df49d7eSJed Brown /// 6239df49d7eSJed Brown /// // Mass operator 624c68be7a2SJeremy L Thompson /// let qf_mass = ceed.q_function_interior_by_name("MassApply")?; 6259df49d7eSJed Brown /// let op_mass = ceed 626c68be7a2SJeremy L Thompson /// .operator(&qf_mass, QFunctionOpt::None, QFunctionOpt::None)? 627c68be7a2SJeremy L Thompson /// .field("u", &ru, &bu, VectorOpt::Active)? 628356036faSJeremy L Thompson /// .field("qdata", &rq, BasisOpt::None, &qdata)? 629c68be7a2SJeremy L Thompson /// .field("v", &ru, &bu, VectorOpt::Active)?; 6309df49d7eSJed Brown /// 6319df49d7eSJed Brown /// v.set_value(0.0); 632c68be7a2SJeremy L Thompson /// op_mass.apply(&u, &mut v)?; 6339df49d7eSJed Brown /// 6349df49d7eSJed Brown /// // Check 635e78171edSJeremy L Thompson /// let sum: Scalar = v.view()?.iter().sum(); 6364b61a5a0SRezgar Shakeri /// let error: Scalar = (sum - 2.0).abs(); 6379df49d7eSJed Brown /// assert!( 6384b61a5a0SRezgar Shakeri /// error < 50.0 * libceed::EPSILON, 6394b61a5a0SRezgar Shakeri /// "Incorrect interval length computed. Expected: 2.0, Found: {}, Error: {:.12e}", 6404b61a5a0SRezgar Shakeri /// sum, 6414b61a5a0SRezgar Shakeri /// error 6429df49d7eSJed Brown /// ); 643c68be7a2SJeremy L Thompson /// # Ok(()) 644c68be7a2SJeremy L Thompson /// # } 6459df49d7eSJed Brown /// ``` 6469df49d7eSJed Brown pub fn apply(&self, input: &Vector, output: &mut Vector) -> crate::Result<i32> { 6479df49d7eSJed Brown self.op_core.apply(input, output) 6489df49d7eSJed Brown } 6499df49d7eSJed Brown 6509df49d7eSJed Brown /// Apply Operator to a vector and add result to output vector 6519df49d7eSJed Brown /// 6529df49d7eSJed Brown /// * `input` - Input Vector 6539df49d7eSJed Brown /// * `output` - Output Vector 6549df49d7eSJed Brown /// 6559df49d7eSJed Brown /// ``` 6569df49d7eSJed Brown /// # use libceed::prelude::*; 6574d27c890SJeremy L Thompson /// # fn main() -> libceed::Result<()> { 6589df49d7eSJed Brown /// # let ceed = libceed::Ceed::default_init(); 6599df49d7eSJed Brown /// let ne = 4; 6609df49d7eSJed Brown /// let p = 3; 6619df49d7eSJed Brown /// let q = 4; 6629df49d7eSJed Brown /// let ndofs = p * ne - ne + 1; 6639df49d7eSJed Brown /// 6649df49d7eSJed Brown /// // Vectors 665c68be7a2SJeremy L Thompson /// let x = ceed.vector_from_slice(&[-1., -0.5, 0.0, 0.5, 1.0])?; 666c68be7a2SJeremy L Thompson /// let mut qdata = ceed.vector(ne * q)?; 6679df49d7eSJed Brown /// qdata.set_value(0.0); 668c68be7a2SJeremy L Thompson /// let u = ceed.vector_from_slice(&vec![1.0; ndofs])?; 669c68be7a2SJeremy L Thompson /// let mut v = ceed.vector(ndofs)?; 6709df49d7eSJed Brown /// 6719df49d7eSJed Brown /// // Restrictions 6729df49d7eSJed Brown /// let mut indx: Vec<i32> = vec![0; 2 * ne]; 6739df49d7eSJed Brown /// for i in 0..ne { 6749df49d7eSJed Brown /// indx[2 * i + 0] = i as i32; 6759df49d7eSJed Brown /// indx[2 * i + 1] = (i + 1) as i32; 6769df49d7eSJed Brown /// } 677c68be7a2SJeremy L Thompson /// let rx = ceed.elem_restriction(ne, 2, 1, 1, ne + 1, MemType::Host, &indx)?; 6789df49d7eSJed Brown /// let mut indu: Vec<i32> = vec![0; p * ne]; 6799df49d7eSJed Brown /// for i in 0..ne { 6809df49d7eSJed Brown /// indu[p * i + 0] = i as i32; 6819df49d7eSJed Brown /// indu[p * i + 1] = (i + 1) as i32; 6829df49d7eSJed Brown /// indu[p * i + 2] = (i + 2) as i32; 6839df49d7eSJed Brown /// } 684c68be7a2SJeremy L Thompson /// let ru = ceed.elem_restriction(ne, 3, 1, 1, ndofs, MemType::Host, &indu)?; 6859df49d7eSJed Brown /// let strides: [i32; 3] = [1, q as i32, q as i32]; 686c68be7a2SJeremy L Thompson /// let rq = ceed.strided_elem_restriction(ne, q, 1, q * ne, strides)?; 6879df49d7eSJed Brown /// 6889df49d7eSJed Brown /// // Bases 689c68be7a2SJeremy L Thompson /// let bx = ceed.basis_tensor_H1_Lagrange(1, 1, 2, q, QuadMode::Gauss)?; 690c68be7a2SJeremy L Thompson /// let bu = ceed.basis_tensor_H1_Lagrange(1, 1, p, q, QuadMode::Gauss)?; 6919df49d7eSJed Brown /// 6929df49d7eSJed Brown /// // Build quadrature data 693c68be7a2SJeremy L Thompson /// let qf_build = ceed.q_function_interior_by_name("Mass1DBuild")?; 694c68be7a2SJeremy L Thompson /// ceed.operator(&qf_build, QFunctionOpt::None, QFunctionOpt::None)? 695c68be7a2SJeremy L Thompson /// .field("dx", &rx, &bx, VectorOpt::Active)? 696c68be7a2SJeremy L Thompson /// .field("weights", ElemRestrictionOpt::None, &bx, VectorOpt::None)? 697356036faSJeremy L Thompson /// .field("qdata", &rq, BasisOpt::None, VectorOpt::Active)? 698c68be7a2SJeremy L Thompson /// .apply(&x, &mut qdata)?; 6999df49d7eSJed Brown /// 7009df49d7eSJed Brown /// // Mass operator 701c68be7a2SJeremy L Thompson /// let qf_mass = ceed.q_function_interior_by_name("MassApply")?; 7029df49d7eSJed Brown /// let op_mass = ceed 703c68be7a2SJeremy L Thompson /// .operator(&qf_mass, QFunctionOpt::None, QFunctionOpt::None)? 704c68be7a2SJeremy L Thompson /// .field("u", &ru, &bu, VectorOpt::Active)? 705356036faSJeremy L Thompson /// .field("qdata", &rq, BasisOpt::None, &qdata)? 706c68be7a2SJeremy L Thompson /// .field("v", &ru, &bu, VectorOpt::Active)?; 7079df49d7eSJed Brown /// 7089df49d7eSJed Brown /// v.set_value(1.0); 709c68be7a2SJeremy L Thompson /// op_mass.apply_add(&u, &mut v)?; 7109df49d7eSJed Brown /// 7119df49d7eSJed Brown /// // Check 712e78171edSJeremy L Thompson /// let sum: Scalar = v.view()?.iter().sum(); 7139df49d7eSJed Brown /// assert!( 71480a9ef05SNatalie Beams /// (sum - (2.0 + ndofs as Scalar)).abs() < 10.0 * libceed::EPSILON, 7159df49d7eSJed Brown /// "Incorrect interval length computed and added" 7169df49d7eSJed Brown /// ); 717c68be7a2SJeremy L Thompson /// # Ok(()) 718c68be7a2SJeremy L Thompson /// # } 7199df49d7eSJed Brown /// ``` 7209df49d7eSJed Brown pub fn apply_add(&self, input: &Vector, output: &mut Vector) -> crate::Result<i32> { 7219df49d7eSJed Brown self.op_core.apply_add(input, output) 7229df49d7eSJed Brown } 7239df49d7eSJed Brown 7249df49d7eSJed Brown /// Provide a field to a Operator for use by its QFunction 7259df49d7eSJed Brown /// 7269df49d7eSJed Brown /// * `fieldname` - Name of the field (to be matched with the name used by 7279df49d7eSJed Brown /// the QFunction) 7289df49d7eSJed Brown /// * `r` - ElemRestriction 729356036faSJeremy L Thompson /// * `b` - Basis in which the field resides or `BasisOpt::None` 730356036faSJeremy L Thompson /// * `v` - Vector to be used by Operator or `VectorOpt::Active` if field 731356036faSJeremy L Thompson /// is active or `VectorOpt::None` if using `Weight` with the 7329df49d7eSJed Brown /// QFunction 7339df49d7eSJed Brown /// 7349df49d7eSJed Brown /// 7359df49d7eSJed Brown /// ``` 7369df49d7eSJed Brown /// # use libceed::prelude::*; 7374d27c890SJeremy L Thompson /// # fn main() -> libceed::Result<()> { 7389df49d7eSJed Brown /// # let ceed = libceed::Ceed::default_init(); 739c68be7a2SJeremy L Thompson /// let qf = ceed.q_function_interior_by_name("Mass1DBuild")?; 740c68be7a2SJeremy L Thompson /// let mut op = ceed.operator(&qf, QFunctionOpt::None, QFunctionOpt::None)?; 7419df49d7eSJed Brown /// 7429df49d7eSJed Brown /// // Operator field arguments 7439df49d7eSJed Brown /// let ne = 3; 7449df49d7eSJed Brown /// let q = 4; 7459df49d7eSJed Brown /// let mut ind: Vec<i32> = vec![0; 2 * ne]; 7469df49d7eSJed Brown /// for i in 0..ne { 7479df49d7eSJed Brown /// ind[2 * i + 0] = i as i32; 7489df49d7eSJed Brown /// ind[2 * i + 1] = (i + 1) as i32; 7499df49d7eSJed Brown /// } 750c68be7a2SJeremy L Thompson /// let r = ceed.elem_restriction(ne, 2, 1, 1, ne + 1, MemType::Host, &ind)?; 7519df49d7eSJed Brown /// 752c68be7a2SJeremy L Thompson /// let b = ceed.basis_tensor_H1_Lagrange(1, 1, 2, q, QuadMode::Gauss)?; 7539df49d7eSJed Brown /// 7549df49d7eSJed Brown /// // Operator field 755c68be7a2SJeremy L Thompson /// op = op.field("dx", &r, &b, VectorOpt::Active)?; 756c68be7a2SJeremy L Thompson /// # Ok(()) 757c68be7a2SJeremy L Thompson /// # } 7589df49d7eSJed Brown /// ``` 7599df49d7eSJed Brown #[allow(unused_mut)] 7609df49d7eSJed Brown pub fn field<'b>( 7619df49d7eSJed Brown mut self, 7629df49d7eSJed Brown fieldname: &str, 7639df49d7eSJed Brown r: impl Into<ElemRestrictionOpt<'b>>, 7649df49d7eSJed Brown b: impl Into<BasisOpt<'b>>, 7659df49d7eSJed Brown v: impl Into<VectorOpt<'b>>, 7669df49d7eSJed Brown ) -> crate::Result<Self> { 7679df49d7eSJed Brown let fieldname = CString::new(fieldname).expect("CString::new failed"); 7689df49d7eSJed Brown let fieldname = fieldname.as_ptr() as *const i8; 7699df49d7eSJed Brown let ierr = unsafe { 7709df49d7eSJed Brown bind_ceed::CeedOperatorSetField( 7719df49d7eSJed Brown self.op_core.ptr, 7729df49d7eSJed Brown fieldname, 7739df49d7eSJed Brown r.into().to_raw(), 7749df49d7eSJed Brown b.into().to_raw(), 7759df49d7eSJed Brown v.into().to_raw(), 7769df49d7eSJed Brown ) 7779df49d7eSJed Brown }; 7781142270cSJeremy L Thompson self.op_core.check_error(ierr)?; 7799df49d7eSJed Brown Ok(self) 7809df49d7eSJed Brown } 7819df49d7eSJed Brown 78208778c6fSJeremy L Thompson /// Get a slice of Operator inputs 78308778c6fSJeremy L Thompson /// 78408778c6fSJeremy L Thompson /// ``` 78508778c6fSJeremy L Thompson /// # use libceed::prelude::*; 7864d27c890SJeremy L Thompson /// # fn main() -> libceed::Result<()> { 78708778c6fSJeremy L Thompson /// # let ceed = libceed::Ceed::default_init(); 78808778c6fSJeremy L Thompson /// let qf = ceed.q_function_interior_by_name("Mass1DBuild")?; 78908778c6fSJeremy L Thompson /// 79008778c6fSJeremy L Thompson /// // Operator field arguments 79108778c6fSJeremy L Thompson /// let ne = 3; 79208778c6fSJeremy L Thompson /// let q = 4 as usize; 79308778c6fSJeremy L Thompson /// let mut ind: Vec<i32> = vec![0; 2 * ne]; 79408778c6fSJeremy L Thompson /// for i in 0..ne { 79508778c6fSJeremy L Thompson /// ind[2 * i + 0] = i as i32; 79608778c6fSJeremy L Thompson /// ind[2 * i + 1] = (i + 1) as i32; 79708778c6fSJeremy L Thompson /// } 79808778c6fSJeremy L Thompson /// let r = ceed.elem_restriction(ne, 2, 1, 1, ne + 1, MemType::Host, &ind)?; 79908778c6fSJeremy L Thompson /// let strides: [i32; 3] = [1, q as i32, q as i32]; 800d9267b76SJeremy L Thompson /// let rq = ceed.strided_elem_restriction(ne, q, 1, q * ne, strides)?; 80108778c6fSJeremy L Thompson /// 80208778c6fSJeremy L Thompson /// let b = ceed.basis_tensor_H1_Lagrange(1, 1, 2, q, QuadMode::Gauss)?; 80308778c6fSJeremy L Thompson /// 80408778c6fSJeremy L Thompson /// // Operator fields 80508778c6fSJeremy L Thompson /// let op = ceed 80608778c6fSJeremy L Thompson /// .operator(&qf, QFunctionOpt::None, QFunctionOpt::None)? 80708778c6fSJeremy L Thompson /// .field("dx", &r, &b, VectorOpt::Active)? 80808778c6fSJeremy L Thompson /// .field("weights", ElemRestrictionOpt::None, &b, VectorOpt::None)? 809356036faSJeremy L Thompson /// .field("qdata", &rq, BasisOpt::None, VectorOpt::Active)?; 81008778c6fSJeremy L Thompson /// 81108778c6fSJeremy L Thompson /// let inputs = op.inputs()?; 81208778c6fSJeremy L Thompson /// 81308778c6fSJeremy L Thompson /// assert_eq!(inputs.len(), 2, "Incorrect inputs array"); 81408778c6fSJeremy L Thompson /// # Ok(()) 81508778c6fSJeremy L Thompson /// # } 81608778c6fSJeremy L Thompson /// ``` 81708778c6fSJeremy L Thompson pub fn inputs(&self) -> crate::Result<&[crate::OperatorField]> { 81808778c6fSJeremy L Thompson // Get array of raw C pointers for inputs 81908778c6fSJeremy L Thompson let mut num_inputs = 0; 82008778c6fSJeremy L Thompson let mut inputs_ptr = std::ptr::null_mut(); 82108778c6fSJeremy L Thompson let ierr = unsafe { 82208778c6fSJeremy L Thompson bind_ceed::CeedOperatorGetFields( 82308778c6fSJeremy L Thompson self.op_core.ptr, 82408778c6fSJeremy L Thompson &mut num_inputs, 82508778c6fSJeremy L Thompson &mut inputs_ptr, 82608778c6fSJeremy L Thompson std::ptr::null_mut() as *mut bind_ceed::CeedInt, 82708778c6fSJeremy L Thompson std::ptr::null_mut() as *mut *mut bind_ceed::CeedOperatorField, 82808778c6fSJeremy L Thompson ) 82908778c6fSJeremy L Thompson }; 83008778c6fSJeremy L Thompson self.op_core.check_error(ierr)?; 83108778c6fSJeremy L Thompson // Convert raw C pointers to fixed length slice 83208778c6fSJeremy L Thompson let inputs_slice = unsafe { 83308778c6fSJeremy L Thompson std::slice::from_raw_parts( 83408778c6fSJeremy L Thompson inputs_ptr as *const crate::OperatorField, 83508778c6fSJeremy L Thompson num_inputs as usize, 83608778c6fSJeremy L Thompson ) 83708778c6fSJeremy L Thompson }; 83808778c6fSJeremy L Thompson Ok(inputs_slice) 83908778c6fSJeremy L Thompson } 84008778c6fSJeremy L Thompson 84108778c6fSJeremy L Thompson /// Get a slice of Operator outputs 84208778c6fSJeremy L Thompson /// 84308778c6fSJeremy L Thompson /// ``` 84408778c6fSJeremy L Thompson /// # use libceed::prelude::*; 8454d27c890SJeremy L Thompson /// # fn main() -> libceed::Result<()> { 84608778c6fSJeremy L Thompson /// # let ceed = libceed::Ceed::default_init(); 84708778c6fSJeremy L Thompson /// let qf = ceed.q_function_interior_by_name("Mass1DBuild")?; 84808778c6fSJeremy L Thompson /// 84908778c6fSJeremy L Thompson /// // Operator field arguments 85008778c6fSJeremy L Thompson /// let ne = 3; 85108778c6fSJeremy L Thompson /// let q = 4 as usize; 85208778c6fSJeremy L Thompson /// let mut ind: Vec<i32> = vec![0; 2 * ne]; 85308778c6fSJeremy L Thompson /// for i in 0..ne { 85408778c6fSJeremy L Thompson /// ind[2 * i + 0] = i as i32; 85508778c6fSJeremy L Thompson /// ind[2 * i + 1] = (i + 1) as i32; 85608778c6fSJeremy L Thompson /// } 85708778c6fSJeremy L Thompson /// let r = ceed.elem_restriction(ne, 2, 1, 1, ne + 1, MemType::Host, &ind)?; 85808778c6fSJeremy L Thompson /// let strides: [i32; 3] = [1, q as i32, q as i32]; 859d9267b76SJeremy L Thompson /// let rq = ceed.strided_elem_restriction(ne, q, 1, q * ne, strides)?; 86008778c6fSJeremy L Thompson /// 86108778c6fSJeremy L Thompson /// let b = ceed.basis_tensor_H1_Lagrange(1, 1, 2, q, QuadMode::Gauss)?; 86208778c6fSJeremy L Thompson /// 86308778c6fSJeremy L Thompson /// // Operator fields 86408778c6fSJeremy L Thompson /// let op = ceed 86508778c6fSJeremy L Thompson /// .operator(&qf, QFunctionOpt::None, QFunctionOpt::None)? 86608778c6fSJeremy L Thompson /// .field("dx", &r, &b, VectorOpt::Active)? 86708778c6fSJeremy L Thompson /// .field("weights", ElemRestrictionOpt::None, &b, VectorOpt::None)? 868356036faSJeremy L Thompson /// .field("qdata", &rq, BasisOpt::None, VectorOpt::Active)?; 86908778c6fSJeremy L Thompson /// 87008778c6fSJeremy L Thompson /// let outputs = op.outputs()?; 87108778c6fSJeremy L Thompson /// 87208778c6fSJeremy L Thompson /// assert_eq!(outputs.len(), 1, "Incorrect outputs array"); 87308778c6fSJeremy L Thompson /// # Ok(()) 87408778c6fSJeremy L Thompson /// # } 87508778c6fSJeremy L Thompson /// ``` 87608778c6fSJeremy L Thompson pub fn outputs(&self) -> crate::Result<&[crate::OperatorField]> { 87708778c6fSJeremy L Thompson // Get array of raw C pointers for outputs 87808778c6fSJeremy L Thompson let mut num_outputs = 0; 87908778c6fSJeremy L Thompson let mut outputs_ptr = std::ptr::null_mut(); 88008778c6fSJeremy L Thompson let ierr = unsafe { 88108778c6fSJeremy L Thompson bind_ceed::CeedOperatorGetFields( 88208778c6fSJeremy L Thompson self.op_core.ptr, 88308778c6fSJeremy L Thompson std::ptr::null_mut() as *mut bind_ceed::CeedInt, 88408778c6fSJeremy L Thompson std::ptr::null_mut() as *mut *mut bind_ceed::CeedOperatorField, 88508778c6fSJeremy L Thompson &mut num_outputs, 88608778c6fSJeremy L Thompson &mut outputs_ptr, 88708778c6fSJeremy L Thompson ) 88808778c6fSJeremy L Thompson }; 88908778c6fSJeremy L Thompson self.op_core.check_error(ierr)?; 89008778c6fSJeremy L Thompson // Convert raw C pointers to fixed length slice 89108778c6fSJeremy L Thompson let outputs_slice = unsafe { 89208778c6fSJeremy L Thompson std::slice::from_raw_parts( 89308778c6fSJeremy L Thompson outputs_ptr as *const crate::OperatorField, 89408778c6fSJeremy L Thompson num_outputs as usize, 89508778c6fSJeremy L Thompson ) 89608778c6fSJeremy L Thompson }; 89708778c6fSJeremy L Thompson Ok(outputs_slice) 89808778c6fSJeremy L Thompson } 89908778c6fSJeremy L Thompson 9006f97ff0aSJeremy L Thompson /// Check if Operator is setup correctly 9016f97ff0aSJeremy L Thompson /// 9026f97ff0aSJeremy L Thompson /// ``` 9036f97ff0aSJeremy L Thompson /// # use libceed::prelude::*; 9046f97ff0aSJeremy L Thompson /// # fn main() -> libceed::Result<()> { 9056f97ff0aSJeremy L Thompson /// # let ceed = libceed::Ceed::default_init(); 9066f97ff0aSJeremy L Thompson /// let ne = 4; 9076f97ff0aSJeremy L Thompson /// let p = 3; 9086f97ff0aSJeremy L Thompson /// let q = 4; 9096f97ff0aSJeremy L Thompson /// let ndofs = p * ne - ne + 1; 9106f97ff0aSJeremy L Thompson /// 9116f97ff0aSJeremy L Thompson /// // Restrictions 9126f97ff0aSJeremy L Thompson /// let mut indx: Vec<i32> = vec![0; 2 * ne]; 9136f97ff0aSJeremy L Thompson /// for i in 0..ne { 9146f97ff0aSJeremy L Thompson /// indx[2 * i + 0] = i as i32; 9156f97ff0aSJeremy L Thompson /// indx[2 * i + 1] = (i + 1) as i32; 9166f97ff0aSJeremy L Thompson /// } 9176f97ff0aSJeremy L Thompson /// let rx = ceed.elem_restriction(ne, 2, 1, 1, ne + 1, MemType::Host, &indx)?; 9186f97ff0aSJeremy L Thompson /// let strides: [i32; 3] = [1, q as i32, q as i32]; 9196f97ff0aSJeremy L Thompson /// let rq = ceed.strided_elem_restriction(ne, q, 1, q * ne, strides)?; 9206f97ff0aSJeremy L Thompson /// 9216f97ff0aSJeremy L Thompson /// // Bases 9226f97ff0aSJeremy L Thompson /// let bx = ceed.basis_tensor_H1_Lagrange(1, 1, 2, q, QuadMode::Gauss)?; 9236f97ff0aSJeremy L Thompson /// 9246f97ff0aSJeremy L Thompson /// // Build quadrature data 9256f97ff0aSJeremy L Thompson /// let qf_build = ceed.q_function_interior_by_name("Mass1DBuild")?; 9266f97ff0aSJeremy L Thompson /// let op_build = ceed 9276f97ff0aSJeremy L Thompson /// .operator(&qf_build, QFunctionOpt::None, QFunctionOpt::None)? 9286f97ff0aSJeremy L Thompson /// .field("dx", &rx, &bx, VectorOpt::Active)? 9296f97ff0aSJeremy L Thompson /// .field("weights", ElemRestrictionOpt::None, &bx, VectorOpt::None)? 930356036faSJeremy L Thompson /// .field("qdata", &rq, BasisOpt::None, VectorOpt::Active)?; 9316f97ff0aSJeremy L Thompson /// 9326f97ff0aSJeremy L Thompson /// // Check 9336f97ff0aSJeremy L Thompson /// assert!(op_build.check().is_ok(), "Operator setup incorrectly"); 9346f97ff0aSJeremy L Thompson /// # Ok(()) 9356f97ff0aSJeremy L Thompson /// # } 9366f97ff0aSJeremy L Thompson /// ``` 9376f97ff0aSJeremy L Thompson pub fn check(self) -> crate::Result<Self> { 9386f97ff0aSJeremy L Thompson self.op_core.check()?; 9396f97ff0aSJeremy L Thompson Ok(self) 9406f97ff0aSJeremy L Thompson } 9416f97ff0aSJeremy L Thompson 9423f2759cfSJeremy L Thompson /// Get the number of elements associated with an Operator 9433f2759cfSJeremy L Thompson /// 9443f2759cfSJeremy L Thompson /// 9453f2759cfSJeremy L Thompson /// ``` 9463f2759cfSJeremy L Thompson /// # use libceed::prelude::*; 9474d27c890SJeremy L Thompson /// # fn main() -> libceed::Result<()> { 9483f2759cfSJeremy L Thompson /// # let ceed = libceed::Ceed::default_init(); 9493f2759cfSJeremy L Thompson /// let qf = ceed.q_function_interior_by_name("Mass1DBuild")?; 9503f2759cfSJeremy L Thompson /// let mut op = ceed.operator(&qf, QFunctionOpt::None, QFunctionOpt::None)?; 9513f2759cfSJeremy L Thompson /// 9523f2759cfSJeremy L Thompson /// // Operator field arguments 9533f2759cfSJeremy L Thompson /// let ne = 3; 9543f2759cfSJeremy L Thompson /// let q = 4; 9553f2759cfSJeremy L Thompson /// let mut ind: Vec<i32> = vec![0; 2 * ne]; 9563f2759cfSJeremy L Thompson /// for i in 0..ne { 9573f2759cfSJeremy L Thompson /// ind[2 * i + 0] = i as i32; 9583f2759cfSJeremy L Thompson /// ind[2 * i + 1] = (i + 1) as i32; 9593f2759cfSJeremy L Thompson /// } 9603f2759cfSJeremy L Thompson /// let r = ceed.elem_restriction(ne, 2, 1, 1, ne + 1, MemType::Host, &ind)?; 9613f2759cfSJeremy L Thompson /// 9623f2759cfSJeremy L Thompson /// let b = ceed.basis_tensor_H1_Lagrange(1, 1, 2, q, QuadMode::Gauss)?; 9633f2759cfSJeremy L Thompson /// 9643f2759cfSJeremy L Thompson /// // Operator field 9653f2759cfSJeremy L Thompson /// op = op.field("dx", &r, &b, VectorOpt::Active)?; 9663f2759cfSJeremy L Thompson /// 9673f2759cfSJeremy L Thompson /// // Check number of elements 9683f2759cfSJeremy L Thompson /// assert_eq!(op.num_elements(), ne, "Incorrect number of elements"); 9693f2759cfSJeremy L Thompson /// # Ok(()) 9703f2759cfSJeremy L Thompson /// # } 9713f2759cfSJeremy L Thompson /// ``` 9723f2759cfSJeremy L Thompson pub fn num_elements(&self) -> usize { 9733f2759cfSJeremy L Thompson let mut nelem = 0; 9743f2759cfSJeremy L Thompson unsafe { bind_ceed::CeedOperatorGetNumElements(self.op_core.ptr, &mut nelem) }; 9753f2759cfSJeremy L Thompson usize::try_from(nelem).unwrap() 9763f2759cfSJeremy L Thompson } 9773f2759cfSJeremy L Thompson 9783f2759cfSJeremy L Thompson /// Get the number of quadrature points associated with each element of 9793f2759cfSJeremy L Thompson /// an Operator 9803f2759cfSJeremy L Thompson /// 9813f2759cfSJeremy L Thompson /// 9823f2759cfSJeremy L Thompson /// ``` 9833f2759cfSJeremy L Thompson /// # use libceed::prelude::*; 9844d27c890SJeremy L Thompson /// # fn main() -> libceed::Result<()> { 9853f2759cfSJeremy L Thompson /// # let ceed = libceed::Ceed::default_init(); 9863f2759cfSJeremy L Thompson /// let qf = ceed.q_function_interior_by_name("Mass1DBuild")?; 9873f2759cfSJeremy L Thompson /// let mut op = ceed.operator(&qf, QFunctionOpt::None, QFunctionOpt::None)?; 9883f2759cfSJeremy L Thompson /// 9893f2759cfSJeremy L Thompson /// // Operator field arguments 9903f2759cfSJeremy L Thompson /// let ne = 3; 9913f2759cfSJeremy L Thompson /// let q = 4; 9923f2759cfSJeremy L Thompson /// let mut ind: Vec<i32> = vec![0; 2 * ne]; 9933f2759cfSJeremy L Thompson /// for i in 0..ne { 9943f2759cfSJeremy L Thompson /// ind[2 * i + 0] = i as i32; 9953f2759cfSJeremy L Thompson /// ind[2 * i + 1] = (i + 1) as i32; 9963f2759cfSJeremy L Thompson /// } 9973f2759cfSJeremy L Thompson /// let r = ceed.elem_restriction(ne, 2, 1, 1, ne + 1, MemType::Host, &ind)?; 9983f2759cfSJeremy L Thompson /// 9993f2759cfSJeremy L Thompson /// let b = ceed.basis_tensor_H1_Lagrange(1, 1, 2, q, QuadMode::Gauss)?; 10003f2759cfSJeremy L Thompson /// 10013f2759cfSJeremy L Thompson /// // Operator field 10023f2759cfSJeremy L Thompson /// op = op.field("dx", &r, &b, VectorOpt::Active)?; 10033f2759cfSJeremy L Thompson /// 10043f2759cfSJeremy L Thompson /// // Check number of quadrature points 10053f2759cfSJeremy L Thompson /// assert_eq!( 10063f2759cfSJeremy L Thompson /// op.num_quadrature_points(), 10073f2759cfSJeremy L Thompson /// q, 10083f2759cfSJeremy L Thompson /// "Incorrect number of quadrature points" 10093f2759cfSJeremy L Thompson /// ); 10103f2759cfSJeremy L Thompson /// # Ok(()) 10113f2759cfSJeremy L Thompson /// # } 10123f2759cfSJeremy L Thompson /// ``` 10133f2759cfSJeremy L Thompson pub fn num_quadrature_points(&self) -> usize { 10143f2759cfSJeremy L Thompson let mut Q = 0; 10153f2759cfSJeremy L Thompson unsafe { bind_ceed::CeedOperatorGetNumQuadraturePoints(self.op_core.ptr, &mut Q) }; 10163f2759cfSJeremy L Thompson usize::try_from(Q).unwrap() 10173f2759cfSJeremy L Thompson } 10183f2759cfSJeremy L Thompson 10199df49d7eSJed Brown /// Assemble the diagonal of a square linear Operator 10209df49d7eSJed Brown /// 10219df49d7eSJed Brown /// This overwrites a Vector with the diagonal of a linear Operator. 10229df49d7eSJed Brown /// 10239df49d7eSJed Brown /// Note: Currently only non-composite Operators with a single field and 10249df49d7eSJed Brown /// composite Operators with single field sub-operators are supported. 10259df49d7eSJed Brown /// 10269df49d7eSJed Brown /// * `op` - Operator to assemble QFunction 10279df49d7eSJed Brown /// * `assembled` - Vector to store assembled Operator diagonal 10289df49d7eSJed Brown /// 10299df49d7eSJed Brown /// ``` 10309df49d7eSJed Brown /// # use libceed::prelude::*; 10314d27c890SJeremy L Thompson /// # fn main() -> libceed::Result<()> { 10329df49d7eSJed Brown /// # let ceed = libceed::Ceed::default_init(); 10339df49d7eSJed Brown /// let ne = 4; 10349df49d7eSJed Brown /// let p = 3; 10359df49d7eSJed Brown /// let q = 4; 10369df49d7eSJed Brown /// let ndofs = p * ne - ne + 1; 10379df49d7eSJed Brown /// 10389df49d7eSJed Brown /// // Vectors 1039c68be7a2SJeremy L Thompson /// let x = ceed.vector_from_slice(&[-1., -0.5, 0.0, 0.5, 1.0])?; 1040c68be7a2SJeremy L Thompson /// let mut qdata = ceed.vector(ne * q)?; 10419df49d7eSJed Brown /// qdata.set_value(0.0); 1042c68be7a2SJeremy L Thompson /// let mut u = ceed.vector(ndofs)?; 10439df49d7eSJed Brown /// u.set_value(1.0); 1044c68be7a2SJeremy L Thompson /// let mut v = ceed.vector(ndofs)?; 10459df49d7eSJed Brown /// v.set_value(0.0); 10469df49d7eSJed Brown /// 10479df49d7eSJed Brown /// // Restrictions 10489df49d7eSJed Brown /// let mut indx: Vec<i32> = vec![0; 2 * ne]; 10499df49d7eSJed Brown /// for i in 0..ne { 10509df49d7eSJed Brown /// indx[2 * i + 0] = i as i32; 10519df49d7eSJed Brown /// indx[2 * i + 1] = (i + 1) as i32; 10529df49d7eSJed Brown /// } 1053c68be7a2SJeremy L Thompson /// let rx = ceed.elem_restriction(ne, 2, 1, 1, ne + 1, MemType::Host, &indx)?; 10549df49d7eSJed Brown /// let mut indu: Vec<i32> = vec![0; p * ne]; 10559df49d7eSJed Brown /// for i in 0..ne { 10569df49d7eSJed Brown /// indu[p * i + 0] = (2 * i) as i32; 10579df49d7eSJed Brown /// indu[p * i + 1] = (2 * i + 1) as i32; 10589df49d7eSJed Brown /// indu[p * i + 2] = (2 * i + 2) as i32; 10599df49d7eSJed Brown /// } 1060c68be7a2SJeremy L Thompson /// let ru = ceed.elem_restriction(ne, p, 1, 1, ndofs, MemType::Host, &indu)?; 10619df49d7eSJed Brown /// let strides: [i32; 3] = [1, q as i32, q as i32]; 1062c68be7a2SJeremy L Thompson /// let rq = ceed.strided_elem_restriction(ne, q, 1, q * ne, strides)?; 10639df49d7eSJed Brown /// 10649df49d7eSJed Brown /// // Bases 1065c68be7a2SJeremy L Thompson /// let bx = ceed.basis_tensor_H1_Lagrange(1, 1, 2, q, QuadMode::Gauss)?; 1066c68be7a2SJeremy L Thompson /// let bu = ceed.basis_tensor_H1_Lagrange(1, 1, p, q, QuadMode::Gauss)?; 10679df49d7eSJed Brown /// 10689df49d7eSJed Brown /// // Build quadrature data 1069c68be7a2SJeremy L Thompson /// let qf_build = ceed.q_function_interior_by_name("Mass1DBuild")?; 1070c68be7a2SJeremy L Thompson /// ceed.operator(&qf_build, QFunctionOpt::None, QFunctionOpt::None)? 1071c68be7a2SJeremy L Thompson /// .field("dx", &rx, &bx, VectorOpt::Active)? 1072c68be7a2SJeremy L Thompson /// .field("weights", ElemRestrictionOpt::None, &bx, VectorOpt::None)? 1073356036faSJeremy L Thompson /// .field("qdata", &rq, BasisOpt::None, VectorOpt::Active)? 1074c68be7a2SJeremy L Thompson /// .apply(&x, &mut qdata)?; 10759df49d7eSJed Brown /// 10769df49d7eSJed Brown /// // Mass operator 1077c68be7a2SJeremy L Thompson /// let qf_mass = ceed.q_function_interior_by_name("MassApply")?; 10789df49d7eSJed Brown /// let op_mass = ceed 1079c68be7a2SJeremy L Thompson /// .operator(&qf_mass, QFunctionOpt::None, QFunctionOpt::None)? 1080c68be7a2SJeremy L Thompson /// .field("u", &ru, &bu, VectorOpt::Active)? 1081356036faSJeremy L Thompson /// .field("qdata", &rq, BasisOpt::None, &qdata)? 1082c68be7a2SJeremy L Thompson /// .field("v", &ru, &bu, VectorOpt::Active)?; 10839df49d7eSJed Brown /// 10849df49d7eSJed Brown /// // Diagonal 1085c68be7a2SJeremy L Thompson /// let mut diag = ceed.vector(ndofs)?; 10869df49d7eSJed Brown /// diag.set_value(0.0); 1087c68be7a2SJeremy L Thompson /// op_mass.linear_assemble_diagonal(&mut diag)?; 10889df49d7eSJed Brown /// 10899df49d7eSJed Brown /// // Manual diagonal computation 1090c68be7a2SJeremy L Thompson /// let mut true_diag = ceed.vector(ndofs)?; 10919c774eddSJeremy L Thompson /// true_diag.set_value(0.0)?; 10929df49d7eSJed Brown /// for i in 0..ndofs { 10939df49d7eSJed Brown /// u.set_value(0.0); 10949df49d7eSJed Brown /// { 1095e78171edSJeremy L Thompson /// let mut u_array = u.view_mut()?; 10969df49d7eSJed Brown /// u_array[i] = 1.; 10979df49d7eSJed Brown /// } 10989df49d7eSJed Brown /// 1099c68be7a2SJeremy L Thompson /// op_mass.apply(&u, &mut v)?; 11009df49d7eSJed Brown /// 11019df49d7eSJed Brown /// { 11029c774eddSJeremy L Thompson /// let v_array = v.view()?; 1103e78171edSJeremy L Thompson /// let mut true_array = true_diag.view_mut()?; 11049df49d7eSJed Brown /// true_array[i] = v_array[i]; 11059df49d7eSJed Brown /// } 11069df49d7eSJed Brown /// } 11079df49d7eSJed Brown /// 11089df49d7eSJed Brown /// // Check 1109e78171edSJeremy L Thompson /// diag.view()? 11109df49d7eSJed Brown /// .iter() 1111e78171edSJeremy L Thompson /// .zip(true_diag.view()?.iter()) 11129df49d7eSJed Brown /// .for_each(|(computed, actual)| { 11139df49d7eSJed Brown /// assert!( 111480a9ef05SNatalie Beams /// (*computed - *actual).abs() < 10.0 * libceed::EPSILON, 11159df49d7eSJed Brown /// "Diagonal entry incorrect" 11169df49d7eSJed Brown /// ); 11179df49d7eSJed Brown /// }); 1118c68be7a2SJeremy L Thompson /// # Ok(()) 1119c68be7a2SJeremy L Thompson /// # } 11209df49d7eSJed Brown /// ``` 11219df49d7eSJed Brown pub fn linear_assemble_diagonal(&self, assembled: &mut Vector) -> crate::Result<i32> { 11229df49d7eSJed Brown self.op_core.linear_assemble_diagonal(assembled) 11239df49d7eSJed Brown } 11249df49d7eSJed Brown 11259df49d7eSJed Brown /// Assemble the diagonal of a square linear Operator 11269df49d7eSJed Brown /// 11279df49d7eSJed Brown /// This sums into a Vector with the diagonal of a linear Operator. 11289df49d7eSJed Brown /// 11299df49d7eSJed Brown /// Note: Currently only non-composite Operators with a single field and 11309df49d7eSJed Brown /// composite Operators with single field sub-operators are supported. 11319df49d7eSJed Brown /// 11329df49d7eSJed Brown /// * `op` - Operator to assemble QFunction 11339df49d7eSJed Brown /// * `assembled` - Vector to store assembled Operator diagonal 11349df49d7eSJed Brown /// 11359df49d7eSJed Brown /// 11369df49d7eSJed Brown /// ``` 11379df49d7eSJed Brown /// # use libceed::prelude::*; 11384d27c890SJeremy L Thompson /// # fn main() -> libceed::Result<()> { 11399df49d7eSJed Brown /// # let ceed = libceed::Ceed::default_init(); 11409df49d7eSJed Brown /// let ne = 4; 11419df49d7eSJed Brown /// let p = 3; 11429df49d7eSJed Brown /// let q = 4; 11439df49d7eSJed Brown /// let ndofs = p * ne - ne + 1; 11449df49d7eSJed Brown /// 11459df49d7eSJed Brown /// // Vectors 1146c68be7a2SJeremy L Thompson /// let x = ceed.vector_from_slice(&[-1., -0.5, 0.0, 0.5, 1.0])?; 1147c68be7a2SJeremy L Thompson /// let mut qdata = ceed.vector(ne * q)?; 11489df49d7eSJed Brown /// qdata.set_value(0.0); 1149c68be7a2SJeremy L Thompson /// let mut u = ceed.vector(ndofs)?; 11509df49d7eSJed Brown /// u.set_value(1.0); 1151c68be7a2SJeremy L Thompson /// let mut v = ceed.vector(ndofs)?; 11529df49d7eSJed Brown /// v.set_value(0.0); 11539df49d7eSJed Brown /// 11549df49d7eSJed Brown /// // Restrictions 11559df49d7eSJed Brown /// let mut indx: Vec<i32> = vec![0; 2 * ne]; 11569df49d7eSJed Brown /// for i in 0..ne { 11579df49d7eSJed Brown /// indx[2 * i + 0] = i as i32; 11589df49d7eSJed Brown /// indx[2 * i + 1] = (i + 1) as i32; 11599df49d7eSJed Brown /// } 1160c68be7a2SJeremy L Thompson /// let rx = ceed.elem_restriction(ne, 2, 1, 1, ne + 1, MemType::Host, &indx)?; 11619df49d7eSJed Brown /// let mut indu: Vec<i32> = vec![0; p * ne]; 11629df49d7eSJed Brown /// for i in 0..ne { 11639df49d7eSJed Brown /// indu[p * i + 0] = (2 * i) as i32; 11649df49d7eSJed Brown /// indu[p * i + 1] = (2 * i + 1) as i32; 11659df49d7eSJed Brown /// indu[p * i + 2] = (2 * i + 2) as i32; 11669df49d7eSJed Brown /// } 1167c68be7a2SJeremy L Thompson /// let ru = ceed.elem_restriction(ne, p, 1, 1, ndofs, MemType::Host, &indu)?; 11689df49d7eSJed Brown /// let strides: [i32; 3] = [1, q as i32, q as i32]; 1169c68be7a2SJeremy L Thompson /// let rq = ceed.strided_elem_restriction(ne, q, 1, q * ne, strides)?; 11709df49d7eSJed Brown /// 11719df49d7eSJed Brown /// // Bases 1172c68be7a2SJeremy L Thompson /// let bx = ceed.basis_tensor_H1_Lagrange(1, 1, 2, q, QuadMode::Gauss)?; 1173c68be7a2SJeremy L Thompson /// let bu = ceed.basis_tensor_H1_Lagrange(1, 1, p, q, QuadMode::Gauss)?; 11749df49d7eSJed Brown /// 11759df49d7eSJed Brown /// // Build quadrature data 1176c68be7a2SJeremy L Thompson /// let qf_build = ceed.q_function_interior_by_name("Mass1DBuild")?; 1177c68be7a2SJeremy L Thompson /// ceed.operator(&qf_build, QFunctionOpt::None, QFunctionOpt::None)? 1178c68be7a2SJeremy L Thompson /// .field("dx", &rx, &bx, VectorOpt::Active)? 1179c68be7a2SJeremy L Thompson /// .field("weights", ElemRestrictionOpt::None, &bx, VectorOpt::None)? 1180356036faSJeremy L Thompson /// .field("qdata", &rq, BasisOpt::None, VectorOpt::Active)? 1181c68be7a2SJeremy L Thompson /// .apply(&x, &mut qdata)?; 11829df49d7eSJed Brown /// 11839df49d7eSJed Brown /// // Mass operator 1184c68be7a2SJeremy L Thompson /// let qf_mass = ceed.q_function_interior_by_name("MassApply")?; 11859df49d7eSJed Brown /// let op_mass = ceed 1186c68be7a2SJeremy L Thompson /// .operator(&qf_mass, QFunctionOpt::None, QFunctionOpt::None)? 1187c68be7a2SJeremy L Thompson /// .field("u", &ru, &bu, VectorOpt::Active)? 1188356036faSJeremy L Thompson /// .field("qdata", &rq, BasisOpt::None, &qdata)? 1189c68be7a2SJeremy L Thompson /// .field("v", &ru, &bu, VectorOpt::Active)?; 11909df49d7eSJed Brown /// 11919df49d7eSJed Brown /// // Diagonal 1192c68be7a2SJeremy L Thompson /// let mut diag = ceed.vector(ndofs)?; 11939df49d7eSJed Brown /// diag.set_value(1.0); 1194c68be7a2SJeremy L Thompson /// op_mass.linear_assemble_add_diagonal(&mut diag)?; 11959df49d7eSJed Brown /// 11969df49d7eSJed Brown /// // Manual diagonal computation 1197c68be7a2SJeremy L Thompson /// let mut true_diag = ceed.vector(ndofs)?; 11989c774eddSJeremy L Thompson /// true_diag.set_value(0.0)?; 11999df49d7eSJed Brown /// for i in 0..ndofs { 12009df49d7eSJed Brown /// u.set_value(0.0); 12019df49d7eSJed Brown /// { 1202e78171edSJeremy L Thompson /// let mut u_array = u.view_mut()?; 12039df49d7eSJed Brown /// u_array[i] = 1.; 12049df49d7eSJed Brown /// } 12059df49d7eSJed Brown /// 1206c68be7a2SJeremy L Thompson /// op_mass.apply(&u, &mut v)?; 12079df49d7eSJed Brown /// 12089df49d7eSJed Brown /// { 12099c774eddSJeremy L Thompson /// let v_array = v.view()?; 1210e78171edSJeremy L Thompson /// let mut true_array = true_diag.view_mut()?; 12119df49d7eSJed Brown /// true_array[i] = v_array[i] + 1.0; 12129df49d7eSJed Brown /// } 12139df49d7eSJed Brown /// } 12149df49d7eSJed Brown /// 12159df49d7eSJed Brown /// // Check 1216e78171edSJeremy L Thompson /// diag.view()? 12179df49d7eSJed Brown /// .iter() 1218e78171edSJeremy L Thompson /// .zip(true_diag.view()?.iter()) 12199df49d7eSJed Brown /// .for_each(|(computed, actual)| { 12209df49d7eSJed Brown /// assert!( 122180a9ef05SNatalie Beams /// (*computed - *actual).abs() < 10.0 * libceed::EPSILON, 12229df49d7eSJed Brown /// "Diagonal entry incorrect" 12239df49d7eSJed Brown /// ); 12249df49d7eSJed Brown /// }); 1225c68be7a2SJeremy L Thompson /// # Ok(()) 1226c68be7a2SJeremy L Thompson /// # } 12279df49d7eSJed Brown /// ``` 12289df49d7eSJed Brown pub fn linear_assemble_add_diagonal(&self, assembled: &mut Vector) -> crate::Result<i32> { 12299df49d7eSJed Brown self.op_core.linear_assemble_add_diagonal(assembled) 12309df49d7eSJed Brown } 12319df49d7eSJed Brown 12329df49d7eSJed Brown /// Assemble the point block diagonal of a square linear Operator 12339df49d7eSJed Brown /// 12349df49d7eSJed Brown /// This overwrites a Vector with the point block diagonal of a linear 12359df49d7eSJed Brown /// Operator. 12369df49d7eSJed Brown /// 12379df49d7eSJed Brown /// Note: Currently only non-composite Operators with a single field and 12389df49d7eSJed Brown /// composite Operators with single field sub-operators are supported. 12399df49d7eSJed Brown /// 12409df49d7eSJed Brown /// * `op` - Operator to assemble QFunction 12419df49d7eSJed Brown /// * `assembled` - Vector to store assembled CeedOperator point block 12429df49d7eSJed Brown /// diagonal, provided in row-major form with an 12439df49d7eSJed Brown /// `ncomp * ncomp` block at each node. The dimensions of 12449df49d7eSJed Brown /// this vector are derived from the active vector for 12459df49d7eSJed Brown /// the CeedOperator. The array has shape 12469df49d7eSJed Brown /// `[nodes, component out, component in]`. 12479df49d7eSJed Brown /// 12489df49d7eSJed Brown /// ``` 12499df49d7eSJed Brown /// # use libceed::prelude::*; 12504d27c890SJeremy L Thompson /// # fn main() -> libceed::Result<()> { 12519df49d7eSJed Brown /// # let ceed = libceed::Ceed::default_init(); 12529df49d7eSJed Brown /// let ne = 4; 12539df49d7eSJed Brown /// let p = 3; 12549df49d7eSJed Brown /// let q = 4; 12559df49d7eSJed Brown /// let ncomp = 2; 12569df49d7eSJed Brown /// let ndofs = p * ne - ne + 1; 12579df49d7eSJed Brown /// 12589df49d7eSJed Brown /// // Vectors 1259c68be7a2SJeremy L Thompson /// let x = ceed.vector_from_slice(&[-1., -0.5, 0.0, 0.5, 1.0])?; 1260c68be7a2SJeremy L Thompson /// let mut qdata = ceed.vector(ne * q)?; 12619df49d7eSJed Brown /// qdata.set_value(0.0); 1262c68be7a2SJeremy L Thompson /// let mut u = ceed.vector(ncomp * ndofs)?; 12639df49d7eSJed Brown /// u.set_value(1.0); 1264c68be7a2SJeremy L Thompson /// let mut v = ceed.vector(ncomp * ndofs)?; 12659df49d7eSJed Brown /// v.set_value(0.0); 12669df49d7eSJed Brown /// 12679df49d7eSJed Brown /// // Restrictions 12689df49d7eSJed Brown /// let mut indx: Vec<i32> = vec![0; 2 * ne]; 12699df49d7eSJed Brown /// for i in 0..ne { 12709df49d7eSJed Brown /// indx[2 * i + 0] = i as i32; 12719df49d7eSJed Brown /// indx[2 * i + 1] = (i + 1) as i32; 12729df49d7eSJed Brown /// } 1273c68be7a2SJeremy L Thompson /// let rx = ceed.elem_restriction(ne, 2, 1, 1, ne + 1, MemType::Host, &indx)?; 12749df49d7eSJed Brown /// let mut indu: Vec<i32> = vec![0; p * ne]; 12759df49d7eSJed Brown /// for i in 0..ne { 12769df49d7eSJed Brown /// indu[p * i + 0] = (2 * i) as i32; 12779df49d7eSJed Brown /// indu[p * i + 1] = (2 * i + 1) as i32; 12789df49d7eSJed Brown /// indu[p * i + 2] = (2 * i + 2) as i32; 12799df49d7eSJed Brown /// } 1280c68be7a2SJeremy L Thompson /// let ru = ceed.elem_restriction(ne, p, ncomp, ndofs, ncomp * ndofs, MemType::Host, &indu)?; 12819df49d7eSJed Brown /// let strides: [i32; 3] = [1, q as i32, q as i32]; 1282c68be7a2SJeremy L Thompson /// let rq = ceed.strided_elem_restriction(ne, q, 1, q * ne, strides)?; 12839df49d7eSJed Brown /// 12849df49d7eSJed Brown /// // Bases 1285c68be7a2SJeremy L Thompson /// let bx = ceed.basis_tensor_H1_Lagrange(1, 1, 2, q, QuadMode::Gauss)?; 1286c68be7a2SJeremy L Thompson /// let bu = ceed.basis_tensor_H1_Lagrange(1, ncomp, p, q, QuadMode::Gauss)?; 12879df49d7eSJed Brown /// 12889df49d7eSJed Brown /// // Build quadrature data 1289c68be7a2SJeremy L Thompson /// let qf_build = ceed.q_function_interior_by_name("Mass1DBuild")?; 1290c68be7a2SJeremy L Thompson /// ceed.operator(&qf_build, QFunctionOpt::None, QFunctionOpt::None)? 1291c68be7a2SJeremy L Thompson /// .field("dx", &rx, &bx, VectorOpt::Active)? 1292c68be7a2SJeremy L Thompson /// .field("weights", ElemRestrictionOpt::None, &bx, VectorOpt::None)? 1293356036faSJeremy L Thompson /// .field("qdata", &rq, BasisOpt::None, VectorOpt::Active)? 1294c68be7a2SJeremy L Thompson /// .apply(&x, &mut qdata)?; 12959df49d7eSJed Brown /// 12969df49d7eSJed Brown /// // Mass operator 12979df49d7eSJed Brown /// let mut mass_2_comp = |[u, qdata, ..]: QFunctionInputs, [v, ..]: QFunctionOutputs| { 12989df49d7eSJed Brown /// // Number of quadrature points 12999df49d7eSJed Brown /// let q = qdata.len(); 13009df49d7eSJed Brown /// 13019df49d7eSJed Brown /// // Iterate over quadrature points 13029df49d7eSJed Brown /// for i in 0..q { 13039df49d7eSJed Brown /// v[i + 0 * q] = u[i + 1 * q] * qdata[i]; 13049df49d7eSJed Brown /// v[i + 1 * q] = u[i + 0 * q] * qdata[i]; 13059df49d7eSJed Brown /// } 13069df49d7eSJed Brown /// 13079df49d7eSJed Brown /// // Return clean error code 13089df49d7eSJed Brown /// 0 13099df49d7eSJed Brown /// }; 13109df49d7eSJed Brown /// 13119df49d7eSJed Brown /// let qf_mass = ceed 1312c68be7a2SJeremy L Thompson /// .q_function_interior(1, Box::new(mass_2_comp))? 1313c68be7a2SJeremy L Thompson /// .input("u", 2, EvalMode::Interp)? 1314c68be7a2SJeremy L Thompson /// .input("qdata", 1, EvalMode::None)? 1315c68be7a2SJeremy L Thompson /// .output("v", 2, EvalMode::Interp)?; 13169df49d7eSJed Brown /// 13179df49d7eSJed Brown /// let op_mass = ceed 1318c68be7a2SJeremy L Thompson /// .operator(&qf_mass, QFunctionOpt::None, QFunctionOpt::None)? 1319c68be7a2SJeremy L Thompson /// .field("u", &ru, &bu, VectorOpt::Active)? 1320356036faSJeremy L Thompson /// .field("qdata", &rq, BasisOpt::None, &qdata)? 1321c68be7a2SJeremy L Thompson /// .field("v", &ru, &bu, VectorOpt::Active)?; 13229df49d7eSJed Brown /// 13239df49d7eSJed Brown /// // Diagonal 1324c68be7a2SJeremy L Thompson /// let mut diag = ceed.vector(ncomp * ncomp * ndofs)?; 13259df49d7eSJed Brown /// diag.set_value(0.0); 1326c68be7a2SJeremy L Thompson /// op_mass.linear_assemble_point_block_diagonal(&mut diag)?; 13279df49d7eSJed Brown /// 13289df49d7eSJed Brown /// // Manual diagonal computation 1329c68be7a2SJeremy L Thompson /// let mut true_diag = ceed.vector(ncomp * ncomp * ndofs)?; 13309c774eddSJeremy L Thompson /// true_diag.set_value(0.0)?; 13319df49d7eSJed Brown /// for i in 0..ndofs { 13329df49d7eSJed Brown /// for j in 0..ncomp { 13339df49d7eSJed Brown /// u.set_value(0.0); 13349df49d7eSJed Brown /// { 1335e78171edSJeremy L Thompson /// let mut u_array = u.view_mut()?; 13369df49d7eSJed Brown /// u_array[i + j * ndofs] = 1.; 13379df49d7eSJed Brown /// } 13389df49d7eSJed Brown /// 1339c68be7a2SJeremy L Thompson /// op_mass.apply(&u, &mut v)?; 13409df49d7eSJed Brown /// 13419df49d7eSJed Brown /// { 13429c774eddSJeremy L Thompson /// let v_array = v.view()?; 1343e78171edSJeremy L Thompson /// let mut true_array = true_diag.view_mut()?; 13449df49d7eSJed Brown /// for k in 0..ncomp { 13459df49d7eSJed Brown /// true_array[i * ncomp * ncomp + k * ncomp + j] = v_array[i + k * ndofs]; 13469df49d7eSJed Brown /// } 13479df49d7eSJed Brown /// } 13489df49d7eSJed Brown /// } 13499df49d7eSJed Brown /// } 13509df49d7eSJed Brown /// 13519df49d7eSJed Brown /// // Check 1352e78171edSJeremy L Thompson /// diag.view()? 13539df49d7eSJed Brown /// .iter() 1354e78171edSJeremy L Thompson /// .zip(true_diag.view()?.iter()) 13559df49d7eSJed Brown /// .for_each(|(computed, actual)| { 13569df49d7eSJed Brown /// assert!( 135780a9ef05SNatalie Beams /// (*computed - *actual).abs() < 10.0 * libceed::EPSILON, 13589df49d7eSJed Brown /// "Diagonal entry incorrect" 13599df49d7eSJed Brown /// ); 13609df49d7eSJed Brown /// }); 1361c68be7a2SJeremy L Thompson /// # Ok(()) 1362c68be7a2SJeremy L Thompson /// # } 13639df49d7eSJed Brown /// ``` 13649df49d7eSJed Brown pub fn linear_assemble_point_block_diagonal( 13659df49d7eSJed Brown &self, 13669df49d7eSJed Brown assembled: &mut Vector, 13679df49d7eSJed Brown ) -> crate::Result<i32> { 13689df49d7eSJed Brown self.op_core.linear_assemble_point_block_diagonal(assembled) 13699df49d7eSJed Brown } 13709df49d7eSJed Brown 13719df49d7eSJed Brown /// Assemble the point block diagonal of a square linear Operator 13729df49d7eSJed Brown /// 13739df49d7eSJed Brown /// This sums into a Vector with the point block diagonal of a linear 13749df49d7eSJed Brown /// Operator. 13759df49d7eSJed Brown /// 13769df49d7eSJed Brown /// Note: Currently only non-composite Operators with a single field and 13779df49d7eSJed Brown /// composite Operators with single field sub-operators are supported. 13789df49d7eSJed Brown /// 13799df49d7eSJed Brown /// * `op` - Operator to assemble QFunction 13809df49d7eSJed Brown /// * `assembled` - Vector to store assembled CeedOperator point block 13819df49d7eSJed Brown /// diagonal, provided in row-major form with an 13829df49d7eSJed Brown /// `ncomp * ncomp` block at each node. The dimensions of 13839df49d7eSJed Brown /// this vector are derived from the active vector for 13849df49d7eSJed Brown /// the CeedOperator. The array has shape 13859df49d7eSJed Brown /// `[nodes, component out, component in]`. 13869df49d7eSJed Brown /// 13879df49d7eSJed Brown /// ``` 13889df49d7eSJed Brown /// # use libceed::prelude::*; 13894d27c890SJeremy L Thompson /// # fn main() -> libceed::Result<()> { 13909df49d7eSJed Brown /// # let ceed = libceed::Ceed::default_init(); 13919df49d7eSJed Brown /// let ne = 4; 13929df49d7eSJed Brown /// let p = 3; 13939df49d7eSJed Brown /// let q = 4; 13949df49d7eSJed Brown /// let ncomp = 2; 13959df49d7eSJed Brown /// let ndofs = p * ne - ne + 1; 13969df49d7eSJed Brown /// 13979df49d7eSJed Brown /// // Vectors 1398c68be7a2SJeremy L Thompson /// let x = ceed.vector_from_slice(&[-1., -0.5, 0.0, 0.5, 1.0])?; 1399c68be7a2SJeremy L Thompson /// let mut qdata = ceed.vector(ne * q)?; 14009df49d7eSJed Brown /// qdata.set_value(0.0); 1401c68be7a2SJeremy L Thompson /// let mut u = ceed.vector(ncomp * ndofs)?; 14029df49d7eSJed Brown /// u.set_value(1.0); 1403c68be7a2SJeremy L Thompson /// let mut v = ceed.vector(ncomp * ndofs)?; 14049df49d7eSJed Brown /// v.set_value(0.0); 14059df49d7eSJed Brown /// 14069df49d7eSJed Brown /// // Restrictions 14079df49d7eSJed Brown /// let mut indx: Vec<i32> = vec![0; 2 * ne]; 14089df49d7eSJed Brown /// for i in 0..ne { 14099df49d7eSJed Brown /// indx[2 * i + 0] = i as i32; 14109df49d7eSJed Brown /// indx[2 * i + 1] = (i + 1) as i32; 14119df49d7eSJed Brown /// } 1412c68be7a2SJeremy L Thompson /// let rx = ceed.elem_restriction(ne, 2, 1, 1, ne + 1, MemType::Host, &indx)?; 14139df49d7eSJed Brown /// let mut indu: Vec<i32> = vec![0; p * ne]; 14149df49d7eSJed Brown /// for i in 0..ne { 14159df49d7eSJed Brown /// indu[p * i + 0] = (2 * i) as i32; 14169df49d7eSJed Brown /// indu[p * i + 1] = (2 * i + 1) as i32; 14179df49d7eSJed Brown /// indu[p * i + 2] = (2 * i + 2) as i32; 14189df49d7eSJed Brown /// } 1419c68be7a2SJeremy L Thompson /// let ru = ceed.elem_restriction(ne, p, ncomp, ndofs, ncomp * ndofs, MemType::Host, &indu)?; 14209df49d7eSJed Brown /// let strides: [i32; 3] = [1, q as i32, q as i32]; 1421c68be7a2SJeremy L Thompson /// let rq = ceed.strided_elem_restriction(ne, q, 1, q * ne, strides)?; 14229df49d7eSJed Brown /// 14239df49d7eSJed Brown /// // Bases 1424c68be7a2SJeremy L Thompson /// let bx = ceed.basis_tensor_H1_Lagrange(1, 1, 2, q, QuadMode::Gauss)?; 1425c68be7a2SJeremy L Thompson /// let bu = ceed.basis_tensor_H1_Lagrange(1, ncomp, p, q, QuadMode::Gauss)?; 14269df49d7eSJed Brown /// 14279df49d7eSJed Brown /// // Build quadrature data 1428c68be7a2SJeremy L Thompson /// let qf_build = ceed.q_function_interior_by_name("Mass1DBuild")?; 1429c68be7a2SJeremy L Thompson /// ceed.operator(&qf_build, QFunctionOpt::None, QFunctionOpt::None)? 1430c68be7a2SJeremy L Thompson /// .field("dx", &rx, &bx, VectorOpt::Active)? 1431c68be7a2SJeremy L Thompson /// .field("weights", ElemRestrictionOpt::None, &bx, VectorOpt::None)? 1432356036faSJeremy L Thompson /// .field("qdata", &rq, BasisOpt::None, VectorOpt::Active)? 1433c68be7a2SJeremy L Thompson /// .apply(&x, &mut qdata)?; 14349df49d7eSJed Brown /// 14359df49d7eSJed Brown /// // Mass operator 14369df49d7eSJed Brown /// let mut mass_2_comp = |[u, qdata, ..]: QFunctionInputs, [v, ..]: QFunctionOutputs| { 14379df49d7eSJed Brown /// // Number of quadrature points 14389df49d7eSJed Brown /// let q = qdata.len(); 14399df49d7eSJed Brown /// 14409df49d7eSJed Brown /// // Iterate over quadrature points 14419df49d7eSJed Brown /// for i in 0..q { 14429df49d7eSJed Brown /// v[i + 0 * q] = u[i + 1 * q] * qdata[i]; 14439df49d7eSJed Brown /// v[i + 1 * q] = u[i + 0 * q] * qdata[i]; 14449df49d7eSJed Brown /// } 14459df49d7eSJed Brown /// 14469df49d7eSJed Brown /// // Return clean error code 14479df49d7eSJed Brown /// 0 14489df49d7eSJed Brown /// }; 14499df49d7eSJed Brown /// 14509df49d7eSJed Brown /// let qf_mass = ceed 1451c68be7a2SJeremy L Thompson /// .q_function_interior(1, Box::new(mass_2_comp))? 1452c68be7a2SJeremy L Thompson /// .input("u", 2, EvalMode::Interp)? 1453c68be7a2SJeremy L Thompson /// .input("qdata", 1, EvalMode::None)? 1454c68be7a2SJeremy L Thompson /// .output("v", 2, EvalMode::Interp)?; 14559df49d7eSJed Brown /// 14569df49d7eSJed Brown /// let op_mass = ceed 1457c68be7a2SJeremy L Thompson /// .operator(&qf_mass, QFunctionOpt::None, QFunctionOpt::None)? 1458c68be7a2SJeremy L Thompson /// .field("u", &ru, &bu, VectorOpt::Active)? 1459356036faSJeremy L Thompson /// .field("qdata", &rq, BasisOpt::None, &qdata)? 1460c68be7a2SJeremy L Thompson /// .field("v", &ru, &bu, VectorOpt::Active)?; 14619df49d7eSJed Brown /// 14629df49d7eSJed Brown /// // Diagonal 1463c68be7a2SJeremy L Thompson /// let mut diag = ceed.vector(ncomp * ncomp * ndofs)?; 14649df49d7eSJed Brown /// diag.set_value(1.0); 1465c68be7a2SJeremy L Thompson /// op_mass.linear_assemble_add_point_block_diagonal(&mut diag)?; 14669df49d7eSJed Brown /// 14679df49d7eSJed Brown /// // Manual diagonal computation 1468c68be7a2SJeremy L Thompson /// let mut true_diag = ceed.vector(ncomp * ncomp * ndofs)?; 14699c774eddSJeremy L Thompson /// true_diag.set_value(0.0)?; 14709df49d7eSJed Brown /// for i in 0..ndofs { 14719df49d7eSJed Brown /// for j in 0..ncomp { 14729df49d7eSJed Brown /// u.set_value(0.0); 14739df49d7eSJed Brown /// { 1474e78171edSJeremy L Thompson /// let mut u_array = u.view_mut()?; 14759df49d7eSJed Brown /// u_array[i + j * ndofs] = 1.; 14769df49d7eSJed Brown /// } 14779df49d7eSJed Brown /// 1478c68be7a2SJeremy L Thompson /// op_mass.apply(&u, &mut v)?; 14799df49d7eSJed Brown /// 14809df49d7eSJed Brown /// { 14819c774eddSJeremy L Thompson /// let v_array = v.view()?; 1482e78171edSJeremy L Thompson /// let mut true_array = true_diag.view_mut()?; 14839df49d7eSJed Brown /// for k in 0..ncomp { 14849df49d7eSJed Brown /// true_array[i * ncomp * ncomp + k * ncomp + j] = v_array[i + k * ndofs]; 14859df49d7eSJed Brown /// } 14869df49d7eSJed Brown /// } 14879df49d7eSJed Brown /// } 14889df49d7eSJed Brown /// } 14899df49d7eSJed Brown /// 14909df49d7eSJed Brown /// // Check 1491e78171edSJeremy L Thompson /// diag.view()? 14929df49d7eSJed Brown /// .iter() 1493e78171edSJeremy L Thompson /// .zip(true_diag.view()?.iter()) 14949df49d7eSJed Brown /// .for_each(|(computed, actual)| { 14959df49d7eSJed Brown /// assert!( 149680a9ef05SNatalie Beams /// (*computed - 1.0 - *actual).abs() < 10.0 * libceed::EPSILON, 14979df49d7eSJed Brown /// "Diagonal entry incorrect" 14989df49d7eSJed Brown /// ); 14999df49d7eSJed Brown /// }); 1500c68be7a2SJeremy L Thompson /// # Ok(()) 1501c68be7a2SJeremy L Thompson /// # } 15029df49d7eSJed Brown /// ``` 15039df49d7eSJed Brown pub fn linear_assemble_add_point_block_diagonal( 15049df49d7eSJed Brown &self, 15059df49d7eSJed Brown assembled: &mut Vector, 15069df49d7eSJed Brown ) -> crate::Result<i32> { 15079df49d7eSJed Brown self.op_core 15089df49d7eSJed Brown .linear_assemble_add_point_block_diagonal(assembled) 15099df49d7eSJed Brown } 15109df49d7eSJed Brown 15119df49d7eSJed Brown /// Create a multigrid coarse Operator and level transfer Operators for a 15129df49d7eSJed Brown /// given Operator, creating the prolongation basis from the fine and 15139df49d7eSJed Brown /// coarse grid interpolation 15149df49d7eSJed Brown /// 15159df49d7eSJed Brown /// * `p_mult_fine` - Lvector multiplicity in parallel gather/scatter 15169df49d7eSJed Brown /// * `rstr_coarse` - Coarse grid restriction 15179df49d7eSJed Brown /// * `basis_coarse` - Coarse grid active vector basis 15189df49d7eSJed Brown /// 15199df49d7eSJed Brown /// ``` 15209df49d7eSJed Brown /// # use libceed::prelude::*; 15214d27c890SJeremy L Thompson /// # fn main() -> libceed::Result<()> { 15229df49d7eSJed Brown /// # let ceed = libceed::Ceed::default_init(); 15239df49d7eSJed Brown /// let ne = 15; 15249df49d7eSJed Brown /// let p_coarse = 3; 15259df49d7eSJed Brown /// let p_fine = 5; 15269df49d7eSJed Brown /// let q = 6; 15279df49d7eSJed Brown /// let ndofs_coarse = p_coarse * ne - ne + 1; 15289df49d7eSJed Brown /// let ndofs_fine = p_fine * ne - ne + 1; 15299df49d7eSJed Brown /// 15309df49d7eSJed Brown /// // Vectors 15319df49d7eSJed Brown /// let x_array = (0..ne + 1) 153280a9ef05SNatalie Beams /// .map(|i| 2.0 * i as Scalar / ne as Scalar - 1.0) 153380a9ef05SNatalie Beams /// .collect::<Vec<Scalar>>(); 1534c68be7a2SJeremy L Thompson /// let x = ceed.vector_from_slice(&x_array)?; 1535c68be7a2SJeremy L Thompson /// let mut qdata = ceed.vector(ne * q)?; 15369df49d7eSJed Brown /// qdata.set_value(0.0); 1537c68be7a2SJeremy L Thompson /// let mut u_coarse = ceed.vector(ndofs_coarse)?; 15389df49d7eSJed Brown /// u_coarse.set_value(1.0); 1539c68be7a2SJeremy L Thompson /// let mut u_fine = ceed.vector(ndofs_fine)?; 15409df49d7eSJed Brown /// u_fine.set_value(1.0); 1541c68be7a2SJeremy L Thompson /// let mut v_coarse = ceed.vector(ndofs_coarse)?; 15429df49d7eSJed Brown /// v_coarse.set_value(0.0); 1543c68be7a2SJeremy L Thompson /// let mut v_fine = ceed.vector(ndofs_fine)?; 15449df49d7eSJed Brown /// v_fine.set_value(0.0); 1545c68be7a2SJeremy L Thompson /// let mut multiplicity = ceed.vector(ndofs_fine)?; 15469df49d7eSJed Brown /// multiplicity.set_value(1.0); 15479df49d7eSJed Brown /// 15489df49d7eSJed Brown /// // Restrictions 15499df49d7eSJed Brown /// let strides: [i32; 3] = [1, q as i32, q as i32]; 1550c68be7a2SJeremy L Thompson /// let rq = ceed.strided_elem_restriction(ne, q, 1, q * ne, strides)?; 15519df49d7eSJed Brown /// 15529df49d7eSJed Brown /// let mut indx: Vec<i32> = vec![0; 2 * ne]; 15539df49d7eSJed Brown /// for i in 0..ne { 15549df49d7eSJed Brown /// indx[2 * i + 0] = i as i32; 15559df49d7eSJed Brown /// indx[2 * i + 1] = (i + 1) as i32; 15569df49d7eSJed Brown /// } 1557c68be7a2SJeremy L Thompson /// let rx = ceed.elem_restriction(ne, 2, 1, 1, ne + 1, MemType::Host, &indx)?; 15589df49d7eSJed Brown /// 15599df49d7eSJed Brown /// let mut indu_coarse: Vec<i32> = vec![0; p_coarse * ne]; 15609df49d7eSJed Brown /// for i in 0..ne { 15619df49d7eSJed Brown /// for j in 0..p_coarse { 15629df49d7eSJed Brown /// indu_coarse[p_coarse * i + j] = (i + j) as i32; 15639df49d7eSJed Brown /// } 15649df49d7eSJed Brown /// } 1565c68be7a2SJeremy L Thompson /// let ru_coarse = ceed.elem_restriction( 15669df49d7eSJed Brown /// ne, 15679df49d7eSJed Brown /// p_coarse, 15689df49d7eSJed Brown /// 1, 15699df49d7eSJed Brown /// 1, 15709df49d7eSJed Brown /// ndofs_coarse, 15719df49d7eSJed Brown /// MemType::Host, 15729df49d7eSJed Brown /// &indu_coarse, 1573c68be7a2SJeremy L Thompson /// )?; 15749df49d7eSJed Brown /// 15759df49d7eSJed Brown /// let mut indu_fine: Vec<i32> = vec![0; p_fine * ne]; 15769df49d7eSJed Brown /// for i in 0..ne { 15779df49d7eSJed Brown /// for j in 0..p_fine { 15789df49d7eSJed Brown /// indu_fine[p_fine * i + j] = (i + j) as i32; 15799df49d7eSJed Brown /// } 15809df49d7eSJed Brown /// } 1581c68be7a2SJeremy L Thompson /// let ru_fine = ceed.elem_restriction(ne, p_fine, 1, 1, ndofs_fine, MemType::Host, &indu_fine)?; 15829df49d7eSJed Brown /// 15839df49d7eSJed Brown /// // Bases 1584c68be7a2SJeremy L Thompson /// let bx = ceed.basis_tensor_H1_Lagrange(1, 1, 2, q, QuadMode::Gauss)?; 1585c68be7a2SJeremy L Thompson /// let bu_coarse = ceed.basis_tensor_H1_Lagrange(1, 1, p_coarse, q, QuadMode::Gauss)?; 1586c68be7a2SJeremy L Thompson /// let bu_fine = ceed.basis_tensor_H1_Lagrange(1, 1, p_fine, q, QuadMode::Gauss)?; 15879df49d7eSJed Brown /// 15889df49d7eSJed Brown /// // Build quadrature data 1589c68be7a2SJeremy L Thompson /// let qf_build = ceed.q_function_interior_by_name("Mass1DBuild")?; 1590c68be7a2SJeremy L Thompson /// ceed.operator(&qf_build, QFunctionOpt::None, QFunctionOpt::None)? 1591c68be7a2SJeremy L Thompson /// .field("dx", &rx, &bx, VectorOpt::Active)? 1592c68be7a2SJeremy L Thompson /// .field("weights", ElemRestrictionOpt::None, &bx, VectorOpt::None)? 1593356036faSJeremy L Thompson /// .field("qdata", &rq, BasisOpt::None, VectorOpt::Active)? 1594c68be7a2SJeremy L Thompson /// .apply(&x, &mut qdata)?; 15959df49d7eSJed Brown /// 15969df49d7eSJed Brown /// // Mass operator 1597c68be7a2SJeremy L Thompson /// let qf_mass = ceed.q_function_interior_by_name("MassApply")?; 15989df49d7eSJed Brown /// let op_mass_fine = ceed 1599c68be7a2SJeremy L Thompson /// .operator(&qf_mass, QFunctionOpt::None, QFunctionOpt::None)? 1600c68be7a2SJeremy L Thompson /// .field("u", &ru_fine, &bu_fine, VectorOpt::Active)? 1601356036faSJeremy L Thompson /// .field("qdata", &rq, BasisOpt::None, &qdata)? 1602c68be7a2SJeremy L Thompson /// .field("v", &ru_fine, &bu_fine, VectorOpt::Active)?; 16039df49d7eSJed Brown /// 16049df49d7eSJed Brown /// // Multigrid setup 1605c68be7a2SJeremy L Thompson /// let (op_mass_coarse, op_prolong, op_restrict) = 1606c68be7a2SJeremy L Thompson /// op_mass_fine.create_multigrid_level(&multiplicity, &ru_coarse, &bu_coarse)?; 16079df49d7eSJed Brown /// 16089df49d7eSJed Brown /// // Coarse problem 16099df49d7eSJed Brown /// u_coarse.set_value(1.0); 1610c68be7a2SJeremy L Thompson /// op_mass_coarse.apply(&u_coarse, &mut v_coarse)?; 16119df49d7eSJed Brown /// 16129df49d7eSJed Brown /// // Check 1613e78171edSJeremy L Thompson /// let sum: Scalar = v_coarse.view()?.iter().sum(); 16149df49d7eSJed Brown /// assert!( 161580a9ef05SNatalie Beams /// (sum - 2.0).abs() < 50.0 * libceed::EPSILON, 16169df49d7eSJed Brown /// "Incorrect interval length computed" 16179df49d7eSJed Brown /// ); 16189df49d7eSJed Brown /// 16199df49d7eSJed Brown /// // Prolong 1620c68be7a2SJeremy L Thompson /// op_prolong.apply(&u_coarse, &mut u_fine)?; 16219df49d7eSJed Brown /// 16229df49d7eSJed Brown /// // Fine problem 1623c68be7a2SJeremy L Thompson /// op_mass_fine.apply(&u_fine, &mut v_fine)?; 16249df49d7eSJed Brown /// 16259df49d7eSJed Brown /// // Check 1626e78171edSJeremy L Thompson /// let sum: Scalar = v_fine.view()?.iter().sum(); 16279df49d7eSJed Brown /// assert!( 162880a9ef05SNatalie Beams /// (sum - 2.0).abs() < 50.0 * libceed::EPSILON, 16299df49d7eSJed Brown /// "Incorrect interval length computed" 16309df49d7eSJed Brown /// ); 16319df49d7eSJed Brown /// 16329df49d7eSJed Brown /// // Restrict 1633c68be7a2SJeremy L Thompson /// op_restrict.apply(&v_fine, &mut v_coarse)?; 16349df49d7eSJed Brown /// 16359df49d7eSJed Brown /// // Check 1636e78171edSJeremy L Thompson /// let sum: Scalar = v_coarse.view()?.iter().sum(); 16379df49d7eSJed Brown /// assert!( 163880a9ef05SNatalie Beams /// (sum - 2.0).abs() < 50.0 * libceed::EPSILON, 16399df49d7eSJed Brown /// "Incorrect interval length computed" 16409df49d7eSJed Brown /// ); 1641c68be7a2SJeremy L Thompson /// # Ok(()) 1642c68be7a2SJeremy L Thompson /// # } 16439df49d7eSJed Brown /// ``` 1644594ef120SJeremy L Thompson pub fn create_multigrid_level<'b>( 16459df49d7eSJed Brown &self, 16469df49d7eSJed Brown p_mult_fine: &Vector, 16479df49d7eSJed Brown rstr_coarse: &ElemRestriction, 16489df49d7eSJed Brown basis_coarse: &Basis, 1649594ef120SJeremy L Thompson ) -> crate::Result<(Operator<'b>, Operator<'b>, Operator<'b>)> { 16509df49d7eSJed Brown let mut ptr_coarse = std::ptr::null_mut(); 16519df49d7eSJed Brown let mut ptr_prolong = std::ptr::null_mut(); 16529df49d7eSJed Brown let mut ptr_restrict = std::ptr::null_mut(); 16539df49d7eSJed Brown let ierr = unsafe { 16549df49d7eSJed Brown bind_ceed::CeedOperatorMultigridLevelCreate( 16559df49d7eSJed Brown self.op_core.ptr, 16569df49d7eSJed Brown p_mult_fine.ptr, 16579df49d7eSJed Brown rstr_coarse.ptr, 16589df49d7eSJed Brown basis_coarse.ptr, 16599df49d7eSJed Brown &mut ptr_coarse, 16609df49d7eSJed Brown &mut ptr_prolong, 16619df49d7eSJed Brown &mut ptr_restrict, 16629df49d7eSJed Brown ) 16639df49d7eSJed Brown }; 16641142270cSJeremy L Thompson self.op_core.check_error(ierr)?; 16651142270cSJeremy L Thompson let op_coarse = Operator::from_raw(ptr_coarse)?; 16661142270cSJeremy L Thompson let op_prolong = Operator::from_raw(ptr_prolong)?; 16671142270cSJeremy L Thompson let op_restrict = Operator::from_raw(ptr_restrict)?; 16689df49d7eSJed Brown Ok((op_coarse, op_prolong, op_restrict)) 16699df49d7eSJed Brown } 16709df49d7eSJed Brown 16719df49d7eSJed Brown /// Create a multigrid coarse Operator and level transfer Operators for a 16729df49d7eSJed Brown /// given Operator with a tensor basis for the active basis 16739df49d7eSJed Brown /// 16749df49d7eSJed Brown /// * `p_mult_fine` - Lvector multiplicity in parallel gather/scatter 16759df49d7eSJed Brown /// * `rstr_coarse` - Coarse grid restriction 16769df49d7eSJed Brown /// * `basis_coarse` - Coarse grid active vector basis 16779df49d7eSJed Brown /// * `interp_c_to_f` - Matrix for coarse to fine 16789df49d7eSJed Brown /// 16799df49d7eSJed Brown /// ``` 16809df49d7eSJed Brown /// # use libceed::prelude::*; 16814d27c890SJeremy L Thompson /// # fn main() -> libceed::Result<()> { 16829df49d7eSJed Brown /// # let ceed = libceed::Ceed::default_init(); 16839df49d7eSJed Brown /// let ne = 15; 16849df49d7eSJed Brown /// let p_coarse = 3; 16859df49d7eSJed Brown /// let p_fine = 5; 16869df49d7eSJed Brown /// let q = 6; 16879df49d7eSJed Brown /// let ndofs_coarse = p_coarse * ne - ne + 1; 16889df49d7eSJed Brown /// let ndofs_fine = p_fine * ne - ne + 1; 16899df49d7eSJed Brown /// 16909df49d7eSJed Brown /// // Vectors 16919df49d7eSJed Brown /// let x_array = (0..ne + 1) 169280a9ef05SNatalie Beams /// .map(|i| 2.0 * i as Scalar / ne as Scalar - 1.0) 169380a9ef05SNatalie Beams /// .collect::<Vec<Scalar>>(); 1694c68be7a2SJeremy L Thompson /// let x = ceed.vector_from_slice(&x_array)?; 1695c68be7a2SJeremy L Thompson /// let mut qdata = ceed.vector(ne * q)?; 16969df49d7eSJed Brown /// qdata.set_value(0.0); 1697c68be7a2SJeremy L Thompson /// let mut u_coarse = ceed.vector(ndofs_coarse)?; 16989df49d7eSJed Brown /// u_coarse.set_value(1.0); 1699c68be7a2SJeremy L Thompson /// let mut u_fine = ceed.vector(ndofs_fine)?; 17009df49d7eSJed Brown /// u_fine.set_value(1.0); 1701c68be7a2SJeremy L Thompson /// let mut v_coarse = ceed.vector(ndofs_coarse)?; 17029df49d7eSJed Brown /// v_coarse.set_value(0.0); 1703c68be7a2SJeremy L Thompson /// let mut v_fine = ceed.vector(ndofs_fine)?; 17049df49d7eSJed Brown /// v_fine.set_value(0.0); 1705c68be7a2SJeremy L Thompson /// let mut multiplicity = ceed.vector(ndofs_fine)?; 17069df49d7eSJed Brown /// multiplicity.set_value(1.0); 17079df49d7eSJed Brown /// 17089df49d7eSJed Brown /// // Restrictions 17099df49d7eSJed Brown /// let strides: [i32; 3] = [1, q as i32, q as i32]; 1710c68be7a2SJeremy L Thompson /// let rq = ceed.strided_elem_restriction(ne, q, 1, q * ne, strides)?; 17119df49d7eSJed Brown /// 17129df49d7eSJed Brown /// let mut indx: Vec<i32> = vec![0; 2 * ne]; 17139df49d7eSJed Brown /// for i in 0..ne { 17149df49d7eSJed Brown /// indx[2 * i + 0] = i as i32; 17159df49d7eSJed Brown /// indx[2 * i + 1] = (i + 1) as i32; 17169df49d7eSJed Brown /// } 1717c68be7a2SJeremy L Thompson /// let rx = ceed.elem_restriction(ne, 2, 1, 1, ne + 1, MemType::Host, &indx)?; 17189df49d7eSJed Brown /// 17199df49d7eSJed Brown /// let mut indu_coarse: Vec<i32> = vec![0; p_coarse * ne]; 17209df49d7eSJed Brown /// for i in 0..ne { 17219df49d7eSJed Brown /// for j in 0..p_coarse { 17229df49d7eSJed Brown /// indu_coarse[p_coarse * i + j] = (i + j) as i32; 17239df49d7eSJed Brown /// } 17249df49d7eSJed Brown /// } 1725c68be7a2SJeremy L Thompson /// let ru_coarse = ceed.elem_restriction( 17269df49d7eSJed Brown /// ne, 17279df49d7eSJed Brown /// p_coarse, 17289df49d7eSJed Brown /// 1, 17299df49d7eSJed Brown /// 1, 17309df49d7eSJed Brown /// ndofs_coarse, 17319df49d7eSJed Brown /// MemType::Host, 17329df49d7eSJed Brown /// &indu_coarse, 1733c68be7a2SJeremy L Thompson /// )?; 17349df49d7eSJed Brown /// 17359df49d7eSJed Brown /// let mut indu_fine: Vec<i32> = vec![0; p_fine * ne]; 17369df49d7eSJed Brown /// for i in 0..ne { 17379df49d7eSJed Brown /// for j in 0..p_fine { 17389df49d7eSJed Brown /// indu_fine[p_fine * i + j] = (i + j) as i32; 17399df49d7eSJed Brown /// } 17409df49d7eSJed Brown /// } 1741c68be7a2SJeremy L Thompson /// let ru_fine = ceed.elem_restriction(ne, p_fine, 1, 1, ndofs_fine, MemType::Host, &indu_fine)?; 17429df49d7eSJed Brown /// 17439df49d7eSJed Brown /// // Bases 1744c68be7a2SJeremy L Thompson /// let bx = ceed.basis_tensor_H1_Lagrange(1, 1, 2, q, QuadMode::Gauss)?; 1745c68be7a2SJeremy L Thompson /// let bu_coarse = ceed.basis_tensor_H1_Lagrange(1, 1, p_coarse, q, QuadMode::Gauss)?; 1746c68be7a2SJeremy L Thompson /// let bu_fine = ceed.basis_tensor_H1_Lagrange(1, 1, p_fine, q, QuadMode::Gauss)?; 17479df49d7eSJed Brown /// 17489df49d7eSJed Brown /// // Build quadrature data 1749c68be7a2SJeremy L Thompson /// let qf_build = ceed.q_function_interior_by_name("Mass1DBuild")?; 1750c68be7a2SJeremy L Thompson /// ceed.operator(&qf_build, QFunctionOpt::None, QFunctionOpt::None)? 1751c68be7a2SJeremy L Thompson /// .field("dx", &rx, &bx, VectorOpt::Active)? 1752c68be7a2SJeremy L Thompson /// .field("weights", ElemRestrictionOpt::None, &bx, VectorOpt::None)? 1753356036faSJeremy L Thompson /// .field("qdata", &rq, BasisOpt::None, VectorOpt::Active)? 1754c68be7a2SJeremy L Thompson /// .apply(&x, &mut qdata)?; 17559df49d7eSJed Brown /// 17569df49d7eSJed Brown /// // Mass operator 1757c68be7a2SJeremy L Thompson /// let qf_mass = ceed.q_function_interior_by_name("MassApply")?; 17589df49d7eSJed Brown /// let op_mass_fine = ceed 1759c68be7a2SJeremy L Thompson /// .operator(&qf_mass, QFunctionOpt::None, QFunctionOpt::None)? 1760c68be7a2SJeremy L Thompson /// .field("u", &ru_fine, &bu_fine, VectorOpt::Active)? 1761356036faSJeremy L Thompson /// .field("qdata", &rq, BasisOpt::None, &qdata)? 1762c68be7a2SJeremy L Thompson /// .field("v", &ru_fine, &bu_fine, VectorOpt::Active)?; 17639df49d7eSJed Brown /// 17649df49d7eSJed Brown /// // Multigrid setup 176580a9ef05SNatalie Beams /// let mut interp_c_to_f: Vec<Scalar> = vec![0.; p_coarse * p_fine]; 17669df49d7eSJed Brown /// { 1767c68be7a2SJeremy L Thompson /// let mut coarse = ceed.vector(p_coarse)?; 1768c68be7a2SJeremy L Thompson /// let mut fine = ceed.vector(p_fine)?; 1769c68be7a2SJeremy L Thompson /// let basis_c_to_f = 1770c68be7a2SJeremy L Thompson /// ceed.basis_tensor_H1_Lagrange(1, 1, p_coarse, p_fine, QuadMode::GaussLobatto)?; 17719df49d7eSJed Brown /// for i in 0..p_coarse { 17729df49d7eSJed Brown /// coarse.set_value(0.0); 17739df49d7eSJed Brown /// { 1774e78171edSJeremy L Thompson /// let mut array = coarse.view_mut()?; 17759df49d7eSJed Brown /// array[i] = 1.; 17769df49d7eSJed Brown /// } 1777c68be7a2SJeremy L Thompson /// basis_c_to_f.apply( 17789df49d7eSJed Brown /// 1, 17799df49d7eSJed Brown /// TransposeMode::NoTranspose, 17809df49d7eSJed Brown /// EvalMode::Interp, 17819df49d7eSJed Brown /// &coarse, 17829df49d7eSJed Brown /// &mut fine, 1783c68be7a2SJeremy L Thompson /// )?; 1784e78171edSJeremy L Thompson /// let array = fine.view()?; 17859df49d7eSJed Brown /// for j in 0..p_fine { 17869df49d7eSJed Brown /// interp_c_to_f[j * p_coarse + i] = array[j]; 17879df49d7eSJed Brown /// } 17889df49d7eSJed Brown /// } 17899df49d7eSJed Brown /// } 1790c68be7a2SJeremy L Thompson /// let (op_mass_coarse, op_prolong, op_restrict) = op_mass_fine.create_multigrid_level_tensor_H1( 1791c68be7a2SJeremy L Thompson /// &multiplicity, 1792c68be7a2SJeremy L Thompson /// &ru_coarse, 1793c68be7a2SJeremy L Thompson /// &bu_coarse, 1794c68be7a2SJeremy L Thompson /// &interp_c_to_f, 1795c68be7a2SJeremy L Thompson /// )?; 17969df49d7eSJed Brown /// 17979df49d7eSJed Brown /// // Coarse problem 17989df49d7eSJed Brown /// u_coarse.set_value(1.0); 1799c68be7a2SJeremy L Thompson /// op_mass_coarse.apply(&u_coarse, &mut v_coarse)?; 18009df49d7eSJed Brown /// 18019df49d7eSJed Brown /// // Check 1802e78171edSJeremy L Thompson /// let sum: Scalar = v_coarse.view()?.iter().sum(); 18039df49d7eSJed Brown /// assert!( 180480a9ef05SNatalie Beams /// (sum - 2.0).abs() < 10.0 * libceed::EPSILON, 18059df49d7eSJed Brown /// "Incorrect interval length computed" 18069df49d7eSJed Brown /// ); 18079df49d7eSJed Brown /// 18089df49d7eSJed Brown /// // Prolong 1809c68be7a2SJeremy L Thompson /// op_prolong.apply(&u_coarse, &mut u_fine)?; 18109df49d7eSJed Brown /// 18119df49d7eSJed Brown /// // Fine problem 1812c68be7a2SJeremy L Thompson /// op_mass_fine.apply(&u_fine, &mut v_fine)?; 18139df49d7eSJed Brown /// 18149df49d7eSJed Brown /// // Check 1815e78171edSJeremy L Thompson /// let sum: Scalar = v_fine.view()?.iter().sum(); 18169df49d7eSJed Brown /// assert!( 181780a9ef05SNatalie Beams /// (sum - 2.0).abs() < 10.0 * libceed::EPSILON, 18189df49d7eSJed Brown /// "Incorrect interval length computed" 18199df49d7eSJed Brown /// ); 18209df49d7eSJed Brown /// 18219df49d7eSJed Brown /// // Restrict 1822c68be7a2SJeremy L Thompson /// op_restrict.apply(&v_fine, &mut v_coarse)?; 18239df49d7eSJed Brown /// 18249df49d7eSJed Brown /// // Check 1825e78171edSJeremy L Thompson /// let sum: Scalar = v_coarse.view()?.iter().sum(); 18269df49d7eSJed Brown /// assert!( 182780a9ef05SNatalie Beams /// (sum - 2.0).abs() < 10.0 * libceed::EPSILON, 18289df49d7eSJed Brown /// "Incorrect interval length computed" 18299df49d7eSJed Brown /// ); 1830c68be7a2SJeremy L Thompson /// # Ok(()) 1831c68be7a2SJeremy L Thompson /// # } 18329df49d7eSJed Brown /// ``` 1833594ef120SJeremy L Thompson pub fn create_multigrid_level_tensor_H1<'b>( 18349df49d7eSJed Brown &self, 18359df49d7eSJed Brown p_mult_fine: &Vector, 18369df49d7eSJed Brown rstr_coarse: &ElemRestriction, 18379df49d7eSJed Brown basis_coarse: &Basis, 183880a9ef05SNatalie Beams interpCtoF: &Vec<Scalar>, 1839594ef120SJeremy L Thompson ) -> crate::Result<(Operator<'b>, Operator<'b>, Operator<'b>)> { 18409df49d7eSJed Brown let mut ptr_coarse = std::ptr::null_mut(); 18419df49d7eSJed Brown let mut ptr_prolong = std::ptr::null_mut(); 18429df49d7eSJed Brown let mut ptr_restrict = std::ptr::null_mut(); 18439df49d7eSJed Brown let ierr = unsafe { 18449df49d7eSJed Brown bind_ceed::CeedOperatorMultigridLevelCreateTensorH1( 18459df49d7eSJed Brown self.op_core.ptr, 18469df49d7eSJed Brown p_mult_fine.ptr, 18479df49d7eSJed Brown rstr_coarse.ptr, 18489df49d7eSJed Brown basis_coarse.ptr, 18499df49d7eSJed Brown interpCtoF.as_ptr(), 18509df49d7eSJed Brown &mut ptr_coarse, 18519df49d7eSJed Brown &mut ptr_prolong, 18529df49d7eSJed Brown &mut ptr_restrict, 18539df49d7eSJed Brown ) 18549df49d7eSJed Brown }; 18551142270cSJeremy L Thompson self.op_core.check_error(ierr)?; 18561142270cSJeremy L Thompson let op_coarse = Operator::from_raw(ptr_coarse)?; 18571142270cSJeremy L Thompson let op_prolong = Operator::from_raw(ptr_prolong)?; 18581142270cSJeremy L Thompson let op_restrict = Operator::from_raw(ptr_restrict)?; 18599df49d7eSJed Brown Ok((op_coarse, op_prolong, op_restrict)) 18609df49d7eSJed Brown } 18619df49d7eSJed Brown 18629df49d7eSJed Brown /// Create a multigrid coarse Operator and level transfer Operators for a 18639df49d7eSJed Brown /// given Operator with a non-tensor basis for the active basis 18649df49d7eSJed Brown /// 18659df49d7eSJed Brown /// * `p_mult_fine` - Lvector multiplicity in parallel gather/scatter 18669df49d7eSJed Brown /// * `rstr_coarse` - Coarse grid restriction 18679df49d7eSJed Brown /// * `basis_coarse` - Coarse grid active vector basis 18689df49d7eSJed Brown /// * `interp_c_to_f` - Matrix for coarse to fine 18699df49d7eSJed Brown /// 18709df49d7eSJed Brown /// ``` 18719df49d7eSJed Brown /// # use libceed::prelude::*; 18724d27c890SJeremy L Thompson /// # fn main() -> libceed::Result<()> { 18739df49d7eSJed Brown /// # let ceed = libceed::Ceed::default_init(); 18749df49d7eSJed Brown /// let ne = 15; 18759df49d7eSJed Brown /// let p_coarse = 3; 18769df49d7eSJed Brown /// let p_fine = 5; 18779df49d7eSJed Brown /// let q = 6; 18789df49d7eSJed Brown /// let ndofs_coarse = p_coarse * ne - ne + 1; 18799df49d7eSJed Brown /// let ndofs_fine = p_fine * ne - ne + 1; 18809df49d7eSJed Brown /// 18819df49d7eSJed Brown /// // Vectors 18829df49d7eSJed Brown /// let x_array = (0..ne + 1) 188380a9ef05SNatalie Beams /// .map(|i| 2.0 * i as Scalar / ne as Scalar - 1.0) 188480a9ef05SNatalie Beams /// .collect::<Vec<Scalar>>(); 1885c68be7a2SJeremy L Thompson /// let x = ceed.vector_from_slice(&x_array)?; 1886c68be7a2SJeremy L Thompson /// let mut qdata = ceed.vector(ne * q)?; 18879df49d7eSJed Brown /// qdata.set_value(0.0); 1888c68be7a2SJeremy L Thompson /// let mut u_coarse = ceed.vector(ndofs_coarse)?; 18899df49d7eSJed Brown /// u_coarse.set_value(1.0); 1890c68be7a2SJeremy L Thompson /// let mut u_fine = ceed.vector(ndofs_fine)?; 18919df49d7eSJed Brown /// u_fine.set_value(1.0); 1892c68be7a2SJeremy L Thompson /// let mut v_coarse = ceed.vector(ndofs_coarse)?; 18939df49d7eSJed Brown /// v_coarse.set_value(0.0); 1894c68be7a2SJeremy L Thompson /// let mut v_fine = ceed.vector(ndofs_fine)?; 18959df49d7eSJed Brown /// v_fine.set_value(0.0); 1896c68be7a2SJeremy L Thompson /// let mut multiplicity = ceed.vector(ndofs_fine)?; 18979df49d7eSJed Brown /// multiplicity.set_value(1.0); 18989df49d7eSJed Brown /// 18999df49d7eSJed Brown /// // Restrictions 19009df49d7eSJed Brown /// let strides: [i32; 3] = [1, q as i32, q as i32]; 1901c68be7a2SJeremy L Thompson /// let rq = ceed.strided_elem_restriction(ne, q, 1, q * ne, strides)?; 19029df49d7eSJed Brown /// 19039df49d7eSJed Brown /// let mut indx: Vec<i32> = vec![0; 2 * ne]; 19049df49d7eSJed Brown /// for i in 0..ne { 19059df49d7eSJed Brown /// indx[2 * i + 0] = i as i32; 19069df49d7eSJed Brown /// indx[2 * i + 1] = (i + 1) as i32; 19079df49d7eSJed Brown /// } 1908c68be7a2SJeremy L Thompson /// let rx = ceed.elem_restriction(ne, 2, 1, 1, ne + 1, MemType::Host, &indx)?; 19099df49d7eSJed Brown /// 19109df49d7eSJed Brown /// let mut indu_coarse: Vec<i32> = vec![0; p_coarse * ne]; 19119df49d7eSJed Brown /// for i in 0..ne { 19129df49d7eSJed Brown /// for j in 0..p_coarse { 19139df49d7eSJed Brown /// indu_coarse[p_coarse * i + j] = (i + j) as i32; 19149df49d7eSJed Brown /// } 19159df49d7eSJed Brown /// } 1916c68be7a2SJeremy L Thompson /// let ru_coarse = ceed.elem_restriction( 19179df49d7eSJed Brown /// ne, 19189df49d7eSJed Brown /// p_coarse, 19199df49d7eSJed Brown /// 1, 19209df49d7eSJed Brown /// 1, 19219df49d7eSJed Brown /// ndofs_coarse, 19229df49d7eSJed Brown /// MemType::Host, 19239df49d7eSJed Brown /// &indu_coarse, 1924c68be7a2SJeremy L Thompson /// )?; 19259df49d7eSJed Brown /// 19269df49d7eSJed Brown /// let mut indu_fine: Vec<i32> = vec![0; p_fine * ne]; 19279df49d7eSJed Brown /// for i in 0..ne { 19289df49d7eSJed Brown /// for j in 0..p_fine { 19299df49d7eSJed Brown /// indu_fine[p_fine * i + j] = (i + j) as i32; 19309df49d7eSJed Brown /// } 19319df49d7eSJed Brown /// } 1932c68be7a2SJeremy L Thompson /// let ru_fine = ceed.elem_restriction(ne, p_fine, 1, 1, ndofs_fine, MemType::Host, &indu_fine)?; 19339df49d7eSJed Brown /// 19349df49d7eSJed Brown /// // Bases 1935c68be7a2SJeremy L Thompson /// let bx = ceed.basis_tensor_H1_Lagrange(1, 1, 2, q, QuadMode::Gauss)?; 1936c68be7a2SJeremy L Thompson /// let bu_coarse = ceed.basis_tensor_H1_Lagrange(1, 1, p_coarse, q, QuadMode::Gauss)?; 1937c68be7a2SJeremy L Thompson /// let bu_fine = ceed.basis_tensor_H1_Lagrange(1, 1, p_fine, q, QuadMode::Gauss)?; 19389df49d7eSJed Brown /// 19399df49d7eSJed Brown /// // Build quadrature data 1940c68be7a2SJeremy L Thompson /// let qf_build = ceed.q_function_interior_by_name("Mass1DBuild")?; 1941c68be7a2SJeremy L Thompson /// ceed.operator(&qf_build, QFunctionOpt::None, QFunctionOpt::None)? 1942c68be7a2SJeremy L Thompson /// .field("dx", &rx, &bx, VectorOpt::Active)? 1943c68be7a2SJeremy L Thompson /// .field("weights", ElemRestrictionOpt::None, &bx, VectorOpt::None)? 1944356036faSJeremy L Thompson /// .field("qdata", &rq, BasisOpt::None, VectorOpt::Active)? 1945c68be7a2SJeremy L Thompson /// .apply(&x, &mut qdata)?; 19469df49d7eSJed Brown /// 19479df49d7eSJed Brown /// // Mass operator 1948c68be7a2SJeremy L Thompson /// let qf_mass = ceed.q_function_interior_by_name("MassApply")?; 19499df49d7eSJed Brown /// let op_mass_fine = ceed 1950c68be7a2SJeremy L Thompson /// .operator(&qf_mass, QFunctionOpt::None, QFunctionOpt::None)? 1951c68be7a2SJeremy L Thompson /// .field("u", &ru_fine, &bu_fine, VectorOpt::Active)? 1952356036faSJeremy L Thompson /// .field("qdata", &rq, BasisOpt::None, &qdata)? 1953c68be7a2SJeremy L Thompson /// .field("v", &ru_fine, &bu_fine, VectorOpt::Active)?; 19549df49d7eSJed Brown /// 19559df49d7eSJed Brown /// // Multigrid setup 195680a9ef05SNatalie Beams /// let mut interp_c_to_f: Vec<Scalar> = vec![0.; p_coarse * p_fine]; 19579df49d7eSJed Brown /// { 1958c68be7a2SJeremy L Thompson /// let mut coarse = ceed.vector(p_coarse)?; 1959c68be7a2SJeremy L Thompson /// let mut fine = ceed.vector(p_fine)?; 1960c68be7a2SJeremy L Thompson /// let basis_c_to_f = 1961c68be7a2SJeremy L Thompson /// ceed.basis_tensor_H1_Lagrange(1, 1, p_coarse, p_fine, QuadMode::GaussLobatto)?; 19629df49d7eSJed Brown /// for i in 0..p_coarse { 19639df49d7eSJed Brown /// coarse.set_value(0.0); 19649df49d7eSJed Brown /// { 1965e78171edSJeremy L Thompson /// let mut array = coarse.view_mut()?; 19669df49d7eSJed Brown /// array[i] = 1.; 19679df49d7eSJed Brown /// } 1968c68be7a2SJeremy L Thompson /// basis_c_to_f.apply( 19699df49d7eSJed Brown /// 1, 19709df49d7eSJed Brown /// TransposeMode::NoTranspose, 19719df49d7eSJed Brown /// EvalMode::Interp, 19729df49d7eSJed Brown /// &coarse, 19739df49d7eSJed Brown /// &mut fine, 1974c68be7a2SJeremy L Thompson /// )?; 1975e78171edSJeremy L Thompson /// let array = fine.view()?; 19769df49d7eSJed Brown /// for j in 0..p_fine { 19779df49d7eSJed Brown /// interp_c_to_f[j * p_coarse + i] = array[j]; 19789df49d7eSJed Brown /// } 19799df49d7eSJed Brown /// } 19809df49d7eSJed Brown /// } 1981c68be7a2SJeremy L Thompson /// let (op_mass_coarse, op_prolong, op_restrict) = op_mass_fine.create_multigrid_level_H1( 1982c68be7a2SJeremy L Thompson /// &multiplicity, 1983c68be7a2SJeremy L Thompson /// &ru_coarse, 1984c68be7a2SJeremy L Thompson /// &bu_coarse, 1985c68be7a2SJeremy L Thompson /// &interp_c_to_f, 1986c68be7a2SJeremy L Thompson /// )?; 19879df49d7eSJed Brown /// 19889df49d7eSJed Brown /// // Coarse problem 19899df49d7eSJed Brown /// u_coarse.set_value(1.0); 1990c68be7a2SJeremy L Thompson /// op_mass_coarse.apply(&u_coarse, &mut v_coarse)?; 19919df49d7eSJed Brown /// 19929df49d7eSJed Brown /// // Check 1993e78171edSJeremy L Thompson /// let sum: Scalar = v_coarse.view()?.iter().sum(); 19949df49d7eSJed Brown /// assert!( 199580a9ef05SNatalie Beams /// (sum - 2.0).abs() < 10.0 * libceed::EPSILON, 19969df49d7eSJed Brown /// "Incorrect interval length computed" 19979df49d7eSJed Brown /// ); 19989df49d7eSJed Brown /// 19999df49d7eSJed Brown /// // Prolong 2000c68be7a2SJeremy L Thompson /// op_prolong.apply(&u_coarse, &mut u_fine)?; 20019df49d7eSJed Brown /// 20029df49d7eSJed Brown /// // Fine problem 2003c68be7a2SJeremy L Thompson /// op_mass_fine.apply(&u_fine, &mut v_fine)?; 20049df49d7eSJed Brown /// 20059df49d7eSJed Brown /// // Check 2006e78171edSJeremy L Thompson /// let sum: Scalar = v_fine.view()?.iter().sum(); 20079df49d7eSJed Brown /// assert!( 200880a9ef05SNatalie Beams /// (sum - 2.0).abs() < 10.0 * libceed::EPSILON, 20099df49d7eSJed Brown /// "Incorrect interval length computed" 20109df49d7eSJed Brown /// ); 20119df49d7eSJed Brown /// 20129df49d7eSJed Brown /// // Restrict 2013c68be7a2SJeremy L Thompson /// op_restrict.apply(&v_fine, &mut v_coarse)?; 20149df49d7eSJed Brown /// 20159df49d7eSJed Brown /// // Check 2016e78171edSJeremy L Thompson /// let sum: Scalar = v_coarse.view()?.iter().sum(); 20179df49d7eSJed Brown /// assert!( 201880a9ef05SNatalie Beams /// (sum - 2.0).abs() < 10.0 * libceed::EPSILON, 20199df49d7eSJed Brown /// "Incorrect interval length computed" 20209df49d7eSJed Brown /// ); 2021c68be7a2SJeremy L Thompson /// # Ok(()) 2022c68be7a2SJeremy L Thompson /// # } 20239df49d7eSJed Brown /// ``` 2024594ef120SJeremy L Thompson pub fn create_multigrid_level_H1<'b>( 20259df49d7eSJed Brown &self, 20269df49d7eSJed Brown p_mult_fine: &Vector, 20279df49d7eSJed Brown rstr_coarse: &ElemRestriction, 20289df49d7eSJed Brown basis_coarse: &Basis, 202980a9ef05SNatalie Beams interpCtoF: &[Scalar], 2030594ef120SJeremy L Thompson ) -> crate::Result<(Operator<'b>, Operator<'b>, Operator<'b>)> { 20319df49d7eSJed Brown let mut ptr_coarse = std::ptr::null_mut(); 20329df49d7eSJed Brown let mut ptr_prolong = std::ptr::null_mut(); 20339df49d7eSJed Brown let mut ptr_restrict = std::ptr::null_mut(); 20349df49d7eSJed Brown let ierr = unsafe { 20359df49d7eSJed Brown bind_ceed::CeedOperatorMultigridLevelCreateH1( 20369df49d7eSJed Brown self.op_core.ptr, 20379df49d7eSJed Brown p_mult_fine.ptr, 20389df49d7eSJed Brown rstr_coarse.ptr, 20399df49d7eSJed Brown basis_coarse.ptr, 20409df49d7eSJed Brown interpCtoF.as_ptr(), 20419df49d7eSJed Brown &mut ptr_coarse, 20429df49d7eSJed Brown &mut ptr_prolong, 20439df49d7eSJed Brown &mut ptr_restrict, 20449df49d7eSJed Brown ) 20459df49d7eSJed Brown }; 20461142270cSJeremy L Thompson self.op_core.check_error(ierr)?; 20471142270cSJeremy L Thompson let op_coarse = Operator::from_raw(ptr_coarse)?; 20481142270cSJeremy L Thompson let op_prolong = Operator::from_raw(ptr_prolong)?; 20491142270cSJeremy L Thompson let op_restrict = Operator::from_raw(ptr_restrict)?; 20509df49d7eSJed Brown Ok((op_coarse, op_prolong, op_restrict)) 20519df49d7eSJed Brown } 20529df49d7eSJed Brown } 20539df49d7eSJed Brown 20549df49d7eSJed Brown // ----------------------------------------------------------------------------- 20559df49d7eSJed Brown // Composite Operator 20569df49d7eSJed Brown // ----------------------------------------------------------------------------- 20579df49d7eSJed Brown impl<'a> CompositeOperator<'a> { 20589df49d7eSJed Brown // Constructor 2059594ef120SJeremy L Thompson pub fn create(ceed: &crate::Ceed) -> crate::Result<Self> { 20609df49d7eSJed Brown let mut ptr = std::ptr::null_mut(); 20619df49d7eSJed Brown let ierr = unsafe { bind_ceed::CeedCompositeOperatorCreate(ceed.ptr, &mut ptr) }; 20629df49d7eSJed Brown ceed.check_error(ierr)?; 20639df49d7eSJed Brown Ok(Self { 20641142270cSJeremy L Thompson op_core: OperatorCore { 20651142270cSJeremy L Thompson ptr, 20661142270cSJeremy L Thompson _lifeline: PhantomData, 20671142270cSJeremy L Thompson }, 20689df49d7eSJed Brown }) 20699df49d7eSJed Brown } 20709df49d7eSJed Brown 2071ea6b5821SJeremy L Thompson /// Set name for CompositeOperator printing 2072ea6b5821SJeremy L Thompson /// 2073ea6b5821SJeremy L Thompson /// * 'name' - Name to set 2074ea6b5821SJeremy L Thompson /// 2075ea6b5821SJeremy L Thompson /// ``` 2076ea6b5821SJeremy L Thompson /// # use libceed::prelude::*; 2077ea6b5821SJeremy L Thompson /// # fn main() -> libceed::Result<()> { 2078ea6b5821SJeremy L Thompson /// # let ceed = libceed::Ceed::default_init(); 2079ea6b5821SJeremy L Thompson /// 2080ea6b5821SJeremy L Thompson /// // Sub operator field arguments 2081ea6b5821SJeremy L Thompson /// let ne = 3; 2082ea6b5821SJeremy L Thompson /// let q = 4 as usize; 2083ea6b5821SJeremy L Thompson /// let mut ind: Vec<i32> = vec![0; 2 * ne]; 2084ea6b5821SJeremy L Thompson /// for i in 0..ne { 2085ea6b5821SJeremy L Thompson /// ind[2 * i + 0] = i as i32; 2086ea6b5821SJeremy L Thompson /// ind[2 * i + 1] = (i + 1) as i32; 2087ea6b5821SJeremy L Thompson /// } 2088ea6b5821SJeremy L Thompson /// let r = ceed.elem_restriction(ne, 2, 1, 1, ne + 1, MemType::Host, &ind)?; 2089ea6b5821SJeremy L Thompson /// let strides: [i32; 3] = [1, q as i32, q as i32]; 2090d9267b76SJeremy L Thompson /// let rq = ceed.strided_elem_restriction(ne, q, 1, q * ne, strides)?; 2091ea6b5821SJeremy L Thompson /// 2092ea6b5821SJeremy L Thompson /// let b = ceed.basis_tensor_H1_Lagrange(1, 1, 2, q, QuadMode::Gauss)?; 2093ea6b5821SJeremy L Thompson /// 2094ea6b5821SJeremy L Thompson /// let qdata_mass = ceed.vector(q * ne)?; 2095ea6b5821SJeremy L Thompson /// let qdata_diff = ceed.vector(q * ne)?; 2096ea6b5821SJeremy L Thompson /// 2097ea6b5821SJeremy L Thompson /// let qf_mass = ceed.q_function_interior_by_name("MassApply")?; 2098ea6b5821SJeremy L Thompson /// let op_mass = ceed 2099ea6b5821SJeremy L Thompson /// .operator(&qf_mass, QFunctionOpt::None, QFunctionOpt::None)? 2100ea6b5821SJeremy L Thompson /// .name("Mass term")? 2101ea6b5821SJeremy L Thompson /// .field("u", &r, &b, VectorOpt::Active)? 2102356036faSJeremy L Thompson /// .field("qdata", &rq, BasisOpt::None, &qdata_mass)? 2103ea6b5821SJeremy L Thompson /// .field("v", &r, &b, VectorOpt::Active)?; 2104ea6b5821SJeremy L Thompson /// 2105ea6b5821SJeremy L Thompson /// let qf_diff = ceed.q_function_interior_by_name("Poisson1DApply")?; 2106ea6b5821SJeremy L Thompson /// let op_diff = ceed 2107ea6b5821SJeremy L Thompson /// .operator(&qf_diff, QFunctionOpt::None, QFunctionOpt::None)? 2108ea6b5821SJeremy L Thompson /// .name("Poisson term")? 2109ea6b5821SJeremy L Thompson /// .field("du", &r, &b, VectorOpt::Active)? 2110356036faSJeremy L Thompson /// .field("qdata", &rq, BasisOpt::None, &qdata_diff)? 2111ea6b5821SJeremy L Thompson /// .field("dv", &r, &b, VectorOpt::Active)?; 2112ea6b5821SJeremy L Thompson /// 2113ea6b5821SJeremy L Thompson /// let op = ceed 2114ea6b5821SJeremy L Thompson /// .composite_operator()? 2115ea6b5821SJeremy L Thompson /// .name("Screened Poisson")? 2116ea6b5821SJeremy L Thompson /// .sub_operator(&op_mass)? 2117ea6b5821SJeremy L Thompson /// .sub_operator(&op_diff)?; 2118ea6b5821SJeremy L Thompson /// # Ok(()) 2119ea6b5821SJeremy L Thompson /// # } 2120ea6b5821SJeremy L Thompson /// ``` 2121ea6b5821SJeremy L Thompson #[allow(unused_mut)] 2122ea6b5821SJeremy L Thompson pub fn name(mut self, name: &str) -> crate::Result<Self> { 2123ea6b5821SJeremy L Thompson self.op_core.name(name)?; 2124ea6b5821SJeremy L Thompson Ok(self) 2125ea6b5821SJeremy L Thompson } 2126ea6b5821SJeremy L Thompson 21279df49d7eSJed Brown /// Apply Operator to a vector 21289df49d7eSJed Brown /// 2129ea6b5821SJeremy L Thompson /// * `input` - Inpuht Vector 21309df49d7eSJed Brown /// * `output` - Output Vector 21319df49d7eSJed Brown /// 21329df49d7eSJed Brown /// ``` 21339df49d7eSJed Brown /// # use libceed::prelude::*; 21344d27c890SJeremy L Thompson /// # fn main() -> libceed::Result<()> { 21359df49d7eSJed Brown /// # let ceed = libceed::Ceed::default_init(); 21369df49d7eSJed Brown /// let ne = 4; 21379df49d7eSJed Brown /// let p = 3; 21389df49d7eSJed Brown /// let q = 4; 21399df49d7eSJed Brown /// let ndofs = p * ne - ne + 1; 21409df49d7eSJed Brown /// 21419df49d7eSJed Brown /// // Vectors 2142c68be7a2SJeremy L Thompson /// let x = ceed.vector_from_slice(&[-1., -0.5, 0.0, 0.5, 1.0])?; 2143c68be7a2SJeremy L Thompson /// let mut qdata_mass = ceed.vector(ne * q)?; 21449df49d7eSJed Brown /// qdata_mass.set_value(0.0); 2145c68be7a2SJeremy L Thompson /// let mut qdata_diff = ceed.vector(ne * q)?; 21469df49d7eSJed Brown /// qdata_diff.set_value(0.0); 2147c68be7a2SJeremy L Thompson /// let mut u = ceed.vector(ndofs)?; 21489df49d7eSJed Brown /// u.set_value(1.0); 2149c68be7a2SJeremy L Thompson /// let mut v = ceed.vector(ndofs)?; 21509df49d7eSJed Brown /// v.set_value(0.0); 21519df49d7eSJed Brown /// 21529df49d7eSJed Brown /// // Restrictions 21539df49d7eSJed Brown /// let mut indx: Vec<i32> = vec![0; 2 * ne]; 21549df49d7eSJed Brown /// for i in 0..ne { 21559df49d7eSJed Brown /// indx[2 * i + 0] = i as i32; 21569df49d7eSJed Brown /// indx[2 * i + 1] = (i + 1) as i32; 21579df49d7eSJed Brown /// } 2158c68be7a2SJeremy L Thompson /// let rx = ceed.elem_restriction(ne, 2, 1, 1, ne + 1, MemType::Host, &indx)?; 21599df49d7eSJed Brown /// let mut indu: Vec<i32> = vec![0; p * ne]; 21609df49d7eSJed Brown /// for i in 0..ne { 21619df49d7eSJed Brown /// indu[p * i + 0] = i as i32; 21629df49d7eSJed Brown /// indu[p * i + 1] = (i + 1) as i32; 21639df49d7eSJed Brown /// indu[p * i + 2] = (i + 2) as i32; 21649df49d7eSJed Brown /// } 2165c68be7a2SJeremy L Thompson /// let ru = ceed.elem_restriction(ne, 3, 1, 1, ndofs, MemType::Host, &indu)?; 21669df49d7eSJed Brown /// let strides: [i32; 3] = [1, q as i32, q as i32]; 2167c68be7a2SJeremy L Thompson /// let rq = ceed.strided_elem_restriction(ne, q, 1, q * ne, strides)?; 21689df49d7eSJed Brown /// 21699df49d7eSJed Brown /// // Bases 2170c68be7a2SJeremy L Thompson /// let bx = ceed.basis_tensor_H1_Lagrange(1, 1, 2, q, QuadMode::Gauss)?; 2171c68be7a2SJeremy L Thompson /// let bu = ceed.basis_tensor_H1_Lagrange(1, 1, p, q, QuadMode::Gauss)?; 21729df49d7eSJed Brown /// 21739df49d7eSJed Brown /// // Build quadrature data 2174c68be7a2SJeremy L Thompson /// let qf_build_mass = ceed.q_function_interior_by_name("Mass1DBuild")?; 2175c68be7a2SJeremy L Thompson /// ceed.operator(&qf_build_mass, QFunctionOpt::None, QFunctionOpt::None)? 2176c68be7a2SJeremy L Thompson /// .field("dx", &rx, &bx, VectorOpt::Active)? 2177c68be7a2SJeremy L Thompson /// .field("weights", ElemRestrictionOpt::None, &bx, VectorOpt::None)? 2178356036faSJeremy L Thompson /// .field("qdata", &rq, BasisOpt::None, VectorOpt::Active)? 2179c68be7a2SJeremy L Thompson /// .apply(&x, &mut qdata_mass)?; 21809df49d7eSJed Brown /// 2181c68be7a2SJeremy L Thompson /// let qf_build_diff = ceed.q_function_interior_by_name("Poisson1DBuild")?; 2182c68be7a2SJeremy L Thompson /// ceed.operator(&qf_build_diff, QFunctionOpt::None, QFunctionOpt::None)? 2183c68be7a2SJeremy L Thompson /// .field("dx", &rx, &bx, VectorOpt::Active)? 2184c68be7a2SJeremy L Thompson /// .field("weights", ElemRestrictionOpt::None, &bx, VectorOpt::None)? 2185356036faSJeremy L Thompson /// .field("qdata", &rq, BasisOpt::None, VectorOpt::Active)? 2186c68be7a2SJeremy L Thompson /// .apply(&x, &mut qdata_diff)?; 21879df49d7eSJed Brown /// 21889df49d7eSJed Brown /// // Application operator 2189c68be7a2SJeremy L Thompson /// let qf_mass = ceed.q_function_interior_by_name("MassApply")?; 21909df49d7eSJed Brown /// let op_mass = ceed 2191c68be7a2SJeremy L Thompson /// .operator(&qf_mass, QFunctionOpt::None, QFunctionOpt::None)? 2192c68be7a2SJeremy L Thompson /// .field("u", &ru, &bu, VectorOpt::Active)? 2193356036faSJeremy L Thompson /// .field("qdata", &rq, BasisOpt::None, &qdata_mass)? 2194c68be7a2SJeremy L Thompson /// .field("v", &ru, &bu, VectorOpt::Active)?; 21959df49d7eSJed Brown /// 2196c68be7a2SJeremy L Thompson /// let qf_diff = ceed.q_function_interior_by_name("Poisson1DApply")?; 21979df49d7eSJed Brown /// let op_diff = ceed 2198c68be7a2SJeremy L Thompson /// .operator(&qf_diff, QFunctionOpt::None, QFunctionOpt::None)? 2199c68be7a2SJeremy L Thompson /// .field("du", &ru, &bu, VectorOpt::Active)? 2200356036faSJeremy L Thompson /// .field("qdata", &rq, BasisOpt::None, &qdata_diff)? 2201c68be7a2SJeremy L Thompson /// .field("dv", &ru, &bu, VectorOpt::Active)?; 22029df49d7eSJed Brown /// 22039df49d7eSJed Brown /// let op_composite = ceed 2204c68be7a2SJeremy L Thompson /// .composite_operator()? 2205c68be7a2SJeremy L Thompson /// .sub_operator(&op_mass)? 2206c68be7a2SJeremy L Thompson /// .sub_operator(&op_diff)?; 22079df49d7eSJed Brown /// 22089df49d7eSJed Brown /// v.set_value(0.0); 2209c68be7a2SJeremy L Thompson /// op_composite.apply(&u, &mut v)?; 22109df49d7eSJed Brown /// 22119df49d7eSJed Brown /// // Check 2212e78171edSJeremy L Thompson /// let sum: Scalar = v.view()?.iter().sum(); 22139df49d7eSJed Brown /// assert!( 221480a9ef05SNatalie Beams /// (sum - 2.0).abs() < 10.0 * libceed::EPSILON, 22159df49d7eSJed Brown /// "Incorrect interval length computed" 22169df49d7eSJed Brown /// ); 2217c68be7a2SJeremy L Thompson /// # Ok(()) 2218c68be7a2SJeremy L Thompson /// # } 22199df49d7eSJed Brown /// ``` 22209df49d7eSJed Brown pub fn apply(&self, input: &Vector, output: &mut Vector) -> crate::Result<i32> { 22219df49d7eSJed Brown self.op_core.apply(input, output) 22229df49d7eSJed Brown } 22239df49d7eSJed Brown 22249df49d7eSJed Brown /// Apply Operator to a vector and add result to output vector 22259df49d7eSJed Brown /// 22269df49d7eSJed Brown /// * `input` - Input Vector 22279df49d7eSJed Brown /// * `output` - Output Vector 22289df49d7eSJed Brown /// 22299df49d7eSJed Brown /// ``` 22309df49d7eSJed Brown /// # use libceed::prelude::*; 22314d27c890SJeremy L Thompson /// # fn main() -> libceed::Result<()> { 22329df49d7eSJed Brown /// # let ceed = libceed::Ceed::default_init(); 22339df49d7eSJed Brown /// let ne = 4; 22349df49d7eSJed Brown /// let p = 3; 22359df49d7eSJed Brown /// let q = 4; 22369df49d7eSJed Brown /// let ndofs = p * ne - ne + 1; 22379df49d7eSJed Brown /// 22389df49d7eSJed Brown /// // Vectors 2239c68be7a2SJeremy L Thompson /// let x = ceed.vector_from_slice(&[-1., -0.5, 0.0, 0.5, 1.0])?; 2240c68be7a2SJeremy L Thompson /// let mut qdata_mass = ceed.vector(ne * q)?; 22419df49d7eSJed Brown /// qdata_mass.set_value(0.0); 2242c68be7a2SJeremy L Thompson /// let mut qdata_diff = ceed.vector(ne * q)?; 22439df49d7eSJed Brown /// qdata_diff.set_value(0.0); 2244c68be7a2SJeremy L Thompson /// let mut u = ceed.vector(ndofs)?; 22459df49d7eSJed Brown /// u.set_value(1.0); 2246c68be7a2SJeremy L Thompson /// let mut v = ceed.vector(ndofs)?; 22479df49d7eSJed Brown /// v.set_value(0.0); 22489df49d7eSJed Brown /// 22499df49d7eSJed Brown /// // Restrictions 22509df49d7eSJed Brown /// let mut indx: Vec<i32> = vec![0; 2 * ne]; 22519df49d7eSJed Brown /// for i in 0..ne { 22529df49d7eSJed Brown /// indx[2 * i + 0] = i as i32; 22539df49d7eSJed Brown /// indx[2 * i + 1] = (i + 1) as i32; 22549df49d7eSJed Brown /// } 2255c68be7a2SJeremy L Thompson /// let rx = ceed.elem_restriction(ne, 2, 1, 1, ne + 1, MemType::Host, &indx)?; 22569df49d7eSJed Brown /// let mut indu: Vec<i32> = vec![0; p * ne]; 22579df49d7eSJed Brown /// for i in 0..ne { 22589df49d7eSJed Brown /// indu[p * i + 0] = i as i32; 22599df49d7eSJed Brown /// indu[p * i + 1] = (i + 1) as i32; 22609df49d7eSJed Brown /// indu[p * i + 2] = (i + 2) as i32; 22619df49d7eSJed Brown /// } 2262c68be7a2SJeremy L Thompson /// let ru = ceed.elem_restriction(ne, 3, 1, 1, ndofs, MemType::Host, &indu)?; 22639df49d7eSJed Brown /// let strides: [i32; 3] = [1, q as i32, q as i32]; 2264c68be7a2SJeremy L Thompson /// let rq = ceed.strided_elem_restriction(ne, q, 1, q * ne, strides)?; 22659df49d7eSJed Brown /// 22669df49d7eSJed Brown /// // Bases 2267c68be7a2SJeremy L Thompson /// let bx = ceed.basis_tensor_H1_Lagrange(1, 1, 2, q, QuadMode::Gauss)?; 2268c68be7a2SJeremy L Thompson /// let bu = ceed.basis_tensor_H1_Lagrange(1, 1, p, q, QuadMode::Gauss)?; 22699df49d7eSJed Brown /// 22709df49d7eSJed Brown /// // Build quadrature data 2271c68be7a2SJeremy L Thompson /// let qf_build_mass = ceed.q_function_interior_by_name("Mass1DBuild")?; 2272c68be7a2SJeremy L Thompson /// ceed.operator(&qf_build_mass, QFunctionOpt::None, QFunctionOpt::None)? 2273c68be7a2SJeremy L Thompson /// .field("dx", &rx, &bx, VectorOpt::Active)? 2274c68be7a2SJeremy L Thompson /// .field("weights", ElemRestrictionOpt::None, &bx, VectorOpt::None)? 2275356036faSJeremy L Thompson /// .field("qdata", &rq, BasisOpt::None, VectorOpt::Active)? 2276c68be7a2SJeremy L Thompson /// .apply(&x, &mut qdata_mass)?; 22779df49d7eSJed Brown /// 2278c68be7a2SJeremy L Thompson /// let qf_build_diff = ceed.q_function_interior_by_name("Poisson1DBuild")?; 2279c68be7a2SJeremy L Thompson /// ceed.operator(&qf_build_diff, QFunctionOpt::None, QFunctionOpt::None)? 2280c68be7a2SJeremy L Thompson /// .field("dx", &rx, &bx, VectorOpt::Active)? 2281c68be7a2SJeremy L Thompson /// .field("weights", ElemRestrictionOpt::None, &bx, VectorOpt::None)? 2282356036faSJeremy L Thompson /// .field("qdata", &rq, BasisOpt::None, VectorOpt::Active)? 2283c68be7a2SJeremy L Thompson /// .apply(&x, &mut qdata_diff)?; 22849df49d7eSJed Brown /// 22859df49d7eSJed Brown /// // Application operator 2286c68be7a2SJeremy L Thompson /// let qf_mass = ceed.q_function_interior_by_name("MassApply")?; 22879df49d7eSJed Brown /// let op_mass = ceed 2288c68be7a2SJeremy L Thompson /// .operator(&qf_mass, QFunctionOpt::None, QFunctionOpt::None)? 2289c68be7a2SJeremy L Thompson /// .field("u", &ru, &bu, VectorOpt::Active)? 2290356036faSJeremy L Thompson /// .field("qdata", &rq, BasisOpt::None, &qdata_mass)? 2291c68be7a2SJeremy L Thompson /// .field("v", &ru, &bu, VectorOpt::Active)?; 22929df49d7eSJed Brown /// 2293c68be7a2SJeremy L Thompson /// let qf_diff = ceed.q_function_interior_by_name("Poisson1DApply")?; 22949df49d7eSJed Brown /// let op_diff = ceed 2295c68be7a2SJeremy L Thompson /// .operator(&qf_diff, QFunctionOpt::None, QFunctionOpt::None)? 2296c68be7a2SJeremy L Thompson /// .field("du", &ru, &bu, VectorOpt::Active)? 2297356036faSJeremy L Thompson /// .field("qdata", &rq, BasisOpt::None, &qdata_diff)? 2298c68be7a2SJeremy L Thompson /// .field("dv", &ru, &bu, VectorOpt::Active)?; 22999df49d7eSJed Brown /// 23009df49d7eSJed Brown /// let op_composite = ceed 2301c68be7a2SJeremy L Thompson /// .composite_operator()? 2302c68be7a2SJeremy L Thompson /// .sub_operator(&op_mass)? 2303c68be7a2SJeremy L Thompson /// .sub_operator(&op_diff)?; 23049df49d7eSJed Brown /// 23059df49d7eSJed Brown /// v.set_value(1.0); 2306c68be7a2SJeremy L Thompson /// op_composite.apply_add(&u, &mut v)?; 23079df49d7eSJed Brown /// 23089df49d7eSJed Brown /// // Check 2309e78171edSJeremy L Thompson /// let sum: Scalar = v.view()?.iter().sum(); 23109df49d7eSJed Brown /// assert!( 231180a9ef05SNatalie Beams /// (sum - (2.0 + ndofs as Scalar)).abs() < 10.0 * libceed::EPSILON, 23129df49d7eSJed Brown /// "Incorrect interval length computed" 23139df49d7eSJed Brown /// ); 2314c68be7a2SJeremy L Thompson /// # Ok(()) 2315c68be7a2SJeremy L Thompson /// # } 23169df49d7eSJed Brown /// ``` 23179df49d7eSJed Brown pub fn apply_add(&self, input: &Vector, output: &mut Vector) -> crate::Result<i32> { 23189df49d7eSJed Brown self.op_core.apply_add(input, output) 23199df49d7eSJed Brown } 23209df49d7eSJed Brown 23219df49d7eSJed Brown /// Add a sub-Operator to a Composite Operator 23229df49d7eSJed Brown /// 23239df49d7eSJed Brown /// * `subop` - Sub-Operator 23249df49d7eSJed Brown /// 23259df49d7eSJed Brown /// ``` 23269df49d7eSJed Brown /// # use libceed::prelude::*; 23274d27c890SJeremy L Thompson /// # fn main() -> libceed::Result<()> { 23289df49d7eSJed Brown /// # let ceed = libceed::Ceed::default_init(); 2329c68be7a2SJeremy L Thompson /// let mut op = ceed.composite_operator()?; 23309df49d7eSJed Brown /// 2331c68be7a2SJeremy L Thompson /// let qf_mass = ceed.q_function_interior_by_name("MassApply")?; 2332c68be7a2SJeremy L Thompson /// let op_mass = ceed.operator(&qf_mass, QFunctionOpt::None, QFunctionOpt::None)?; 2333c68be7a2SJeremy L Thompson /// op = op.sub_operator(&op_mass)?; 23349df49d7eSJed Brown /// 2335c68be7a2SJeremy L Thompson /// let qf_diff = ceed.q_function_interior_by_name("Poisson1DApply")?; 2336c68be7a2SJeremy L Thompson /// let op_diff = ceed.operator(&qf_diff, QFunctionOpt::None, QFunctionOpt::None)?; 2337c68be7a2SJeremy L Thompson /// op = op.sub_operator(&op_diff)?; 2338c68be7a2SJeremy L Thompson /// # Ok(()) 2339c68be7a2SJeremy L Thompson /// # } 23409df49d7eSJed Brown /// ``` 23419df49d7eSJed Brown #[allow(unused_mut)] 23429df49d7eSJed Brown pub fn sub_operator(mut self, subop: &Operator) -> crate::Result<Self> { 23439df49d7eSJed Brown let ierr = 23449df49d7eSJed Brown unsafe { bind_ceed::CeedCompositeOperatorAddSub(self.op_core.ptr, subop.op_core.ptr) }; 23451142270cSJeremy L Thompson self.op_core.check_error(ierr)?; 23469df49d7eSJed Brown Ok(self) 23479df49d7eSJed Brown } 23489df49d7eSJed Brown 23496f97ff0aSJeremy L Thompson /// Check if CompositeOperator is setup correctly 23506f97ff0aSJeremy L Thompson /// 23516f97ff0aSJeremy L Thompson /// ``` 23526f97ff0aSJeremy L Thompson /// # use libceed::prelude::*; 23536f97ff0aSJeremy L Thompson /// # fn main() -> libceed::Result<()> { 23546f97ff0aSJeremy L Thompson /// # let ceed = libceed::Ceed::default_init(); 23556f97ff0aSJeremy L Thompson /// let ne = 4; 23566f97ff0aSJeremy L Thompson /// let p = 3; 23576f97ff0aSJeremy L Thompson /// let q = 4; 23586f97ff0aSJeremy L Thompson /// let ndofs = p * ne - ne + 1; 23596f97ff0aSJeremy L Thompson /// 23606f97ff0aSJeremy L Thompson /// // Restrictions 23616f97ff0aSJeremy L Thompson /// let mut indx: Vec<i32> = vec![0; 2 * ne]; 23626f97ff0aSJeremy L Thompson /// for i in 0..ne { 23636f97ff0aSJeremy L Thompson /// indx[2 * i + 0] = i as i32; 23646f97ff0aSJeremy L Thompson /// indx[2 * i + 1] = (i + 1) as i32; 23656f97ff0aSJeremy L Thompson /// } 23666f97ff0aSJeremy L Thompson /// let rx = ceed.elem_restriction(ne, 2, 1, 1, ne + 1, MemType::Host, &indx)?; 23676f97ff0aSJeremy L Thompson /// let strides: [i32; 3] = [1, q as i32, q as i32]; 23686f97ff0aSJeremy L Thompson /// let rq = ceed.strided_elem_restriction(ne, q, 1, q * ne, strides)?; 23696f97ff0aSJeremy L Thompson /// 23706f97ff0aSJeremy L Thompson /// // Bases 23716f97ff0aSJeremy L Thompson /// let bx = ceed.basis_tensor_H1_Lagrange(1, 1, 2, q, QuadMode::Gauss)?; 23726f97ff0aSJeremy L Thompson /// 23736f97ff0aSJeremy L Thompson /// // Build quadrature data 23746f97ff0aSJeremy L Thompson /// let qf_build_mass = ceed.q_function_interior_by_name("Mass1DBuild")?; 23756f97ff0aSJeremy L Thompson /// let op_build_mass = ceed 23766f97ff0aSJeremy L Thompson /// .operator(&qf_build_mass, QFunctionOpt::None, QFunctionOpt::None)? 23776f97ff0aSJeremy L Thompson /// .field("dx", &rx, &bx, VectorOpt::Active)? 23786f97ff0aSJeremy L Thompson /// .field("weights", ElemRestrictionOpt::None, &bx, VectorOpt::None)? 2379356036faSJeremy L Thompson /// .field("qdata", &rq, BasisOpt::None, VectorOpt::Active)?; 23806f97ff0aSJeremy L Thompson /// 23816f97ff0aSJeremy L Thompson /// let qf_build_diff = ceed.q_function_interior_by_name("Poisson1DBuild")?; 23826f97ff0aSJeremy L Thompson /// let op_build_diff = ceed 23836f97ff0aSJeremy L Thompson /// .operator(&qf_build_diff, QFunctionOpt::None, QFunctionOpt::None)? 23846f97ff0aSJeremy L Thompson /// .field("dx", &rx, &bx, VectorOpt::Active)? 23856f97ff0aSJeremy L Thompson /// .field("weights", ElemRestrictionOpt::None, &bx, VectorOpt::None)? 2386356036faSJeremy L Thompson /// .field("qdata", &rq, BasisOpt::None, VectorOpt::Active)?; 23876f97ff0aSJeremy L Thompson /// 23886f97ff0aSJeremy L Thompson /// let op_build = ceed 23896f97ff0aSJeremy L Thompson /// .composite_operator()? 23906f97ff0aSJeremy L Thompson /// .sub_operator(&op_build_mass)? 23916f97ff0aSJeremy L Thompson /// .sub_operator(&op_build_diff)?; 23926f97ff0aSJeremy L Thompson /// 23936f97ff0aSJeremy L Thompson /// // Check 23946f97ff0aSJeremy L Thompson /// assert!(op_build.check().is_ok(), "Operator setup incorrectly"); 23956f97ff0aSJeremy L Thompson /// # Ok(()) 23966f97ff0aSJeremy L Thompson /// # } 23976f97ff0aSJeremy L Thompson /// ``` 23986f97ff0aSJeremy L Thompson pub fn check(self) -> crate::Result<Self> { 23996f97ff0aSJeremy L Thompson self.op_core.check()?; 24006f97ff0aSJeremy L Thompson Ok(self) 24016f97ff0aSJeremy L Thompson } 24026f97ff0aSJeremy L Thompson 24039df49d7eSJed Brown /// Assemble the diagonal of a square linear Operator 24049df49d7eSJed Brown /// 24059df49d7eSJed Brown /// This overwrites a Vector with the diagonal of a linear Operator. 24069df49d7eSJed Brown /// 24079df49d7eSJed Brown /// Note: Currently only non-composite Operators with a single field and 24089df49d7eSJed Brown /// composite Operators with single field sub-operators are supported. 24099df49d7eSJed Brown /// 24109df49d7eSJed Brown /// * `op` - Operator to assemble QFunction 24119df49d7eSJed Brown /// * `assembled` - Vector to store assembled Operator diagonal 2412b748b478SJeremy L Thompson pub fn linear_assemble_diagonal(&self, assembled: &mut Vector) -> crate::Result<i32> { 24139df49d7eSJed Brown self.op_core.linear_assemble_add_diagonal(assembled) 24149df49d7eSJed Brown } 24159df49d7eSJed Brown 24169df49d7eSJed Brown /// Assemble the point block diagonal of a square linear Operator 24179df49d7eSJed Brown /// 24189df49d7eSJed Brown /// This overwrites a Vector with the point block diagonal of a linear 24199df49d7eSJed Brown /// Operator. 24209df49d7eSJed Brown /// 24219df49d7eSJed Brown /// Note: Currently only non-composite Operators with a single field and 24229df49d7eSJed Brown /// composite Operators with single field sub-operators are supported. 24239df49d7eSJed Brown /// 24249df49d7eSJed Brown /// * `op` - Operator to assemble QFunction 24259df49d7eSJed Brown /// * `assembled` - Vector to store assembled Operator diagonal 24269df49d7eSJed Brown pub fn linear_assemble_add_diagonal(&self, assembled: &mut Vector) -> crate::Result<i32> { 24279df49d7eSJed Brown self.op_core.linear_assemble_add_diagonal(assembled) 24289df49d7eSJed Brown } 24299df49d7eSJed Brown 24309df49d7eSJed Brown /// Assemble the diagonal of a square linear Operator 24319df49d7eSJed Brown /// 24329df49d7eSJed Brown /// This overwrites a Vector with the diagonal of a linear Operator. 24339df49d7eSJed Brown /// 24349df49d7eSJed Brown /// Note: Currently only non-composite Operators with a single field and 24359df49d7eSJed Brown /// composite Operators with single field sub-operators are supported. 24369df49d7eSJed Brown /// 24379df49d7eSJed Brown /// * `op` - Operator to assemble QFunction 24389df49d7eSJed Brown /// * `assembled` - Vector to store assembled CeedOperator point block 24399df49d7eSJed Brown /// diagonal, provided in row-major form with an 24409df49d7eSJed Brown /// `ncomp * ncomp` block at each node. The dimensions of 24419df49d7eSJed Brown /// this vector are derived from the active vector for 24429df49d7eSJed Brown /// the CeedOperator. The array has shape 24439df49d7eSJed Brown /// `[nodes, component out, component in]`. 24449df49d7eSJed Brown pub fn linear_assemble_point_block_diagonal( 24459df49d7eSJed Brown &self, 24469df49d7eSJed Brown assembled: &mut Vector, 24479df49d7eSJed Brown ) -> crate::Result<i32> { 24489df49d7eSJed Brown self.op_core.linear_assemble_add_diagonal(assembled) 24499df49d7eSJed Brown } 24509df49d7eSJed Brown 24519df49d7eSJed Brown /// Assemble the diagonal of a square linear Operator 24529df49d7eSJed Brown /// 24539df49d7eSJed Brown /// This sums into a Vector with the diagonal of a linear Operator. 24549df49d7eSJed Brown /// 24559df49d7eSJed Brown /// Note: Currently only non-composite Operators with a single field and 24569df49d7eSJed Brown /// composite Operators with single field sub-operators are supported. 24579df49d7eSJed Brown /// 24589df49d7eSJed Brown /// * `op` - Operator to assemble QFunction 24599df49d7eSJed Brown /// * `assembled` - Vector to store assembled CeedOperator point block 24609df49d7eSJed Brown /// diagonal, provided in row-major form with an 24619df49d7eSJed Brown /// `ncomp * ncomp` block at each node. The dimensions of 24629df49d7eSJed Brown /// this vector are derived from the active vector for 24639df49d7eSJed Brown /// the CeedOperator. The array has shape 24649df49d7eSJed Brown /// `[nodes, component out, component in]`. 24659df49d7eSJed Brown pub fn linear_assemble_add_point_block_diagonal( 24669df49d7eSJed Brown &self, 24679df49d7eSJed Brown assembled: &mut Vector, 24689df49d7eSJed Brown ) -> crate::Result<i32> { 24699df49d7eSJed Brown self.op_core.linear_assemble_add_diagonal(assembled) 24709df49d7eSJed Brown } 24719df49d7eSJed Brown } 24729df49d7eSJed Brown 24739df49d7eSJed Brown // ----------------------------------------------------------------------------- 2474