19df49d7eSJed Brown // Copyright (c) 2017, Lawrence Livermore National Security, LLC. Produced at 29df49d7eSJed Brown // the Lawrence Livermore National Laboratory. LLNL-CODE-734707. All Rights 39df49d7eSJed Brown // reserved. See files LICENSE and NOTICE for details. 49df49d7eSJed Brown // 59df49d7eSJed Brown // This file is part of CEED, a collection of benchmarks, miniapps, software 69df49d7eSJed Brown // libraries and APIs for efficient high-order finite element and spectral 79df49d7eSJed Brown // element discretizations for exascale applications. For more information and 89df49d7eSJed Brown // source code availability see http://github.com/ceed. 99df49d7eSJed Brown // 109df49d7eSJed Brown // The CEED research is supported by the Exascale Computing Project 17-SC-20-SC, 119df49d7eSJed Brown // a collaborative effort of two U.S. Department of Energy organizations (Office 129df49d7eSJed Brown // of Science and the National Nuclear Security Administration) responsible for 139df49d7eSJed Brown // the planning and preparation of a capable exascale ecosystem, including 149df49d7eSJed Brown // software, applications, hardware, advanced system engineering and early 159df49d7eSJed Brown // testbed platforms, in support of the nation's exascale computing imperative 169df49d7eSJed Brown 179df49d7eSJed Brown //! A Ceed Operator defines the finite/spectral element operator associated to a 189df49d7eSJed Brown //! Ceed QFunction. A Ceed Operator connects Ceed ElemRestrictions, 199df49d7eSJed Brown //! Ceed Bases, and Ceed QFunctions. 209df49d7eSJed Brown 219df49d7eSJed Brown use crate::prelude::*; 229df49d7eSJed Brown 239df49d7eSJed Brown // ----------------------------------------------------------------------------- 249df49d7eSJed Brown // CeedOperator context wrapper 259df49d7eSJed Brown // ----------------------------------------------------------------------------- 26c68be7a2SJeremy L Thompson #[derive(Debug)] 279df49d7eSJed Brown pub(crate) struct OperatorCore<'a> { 289df49d7eSJed Brown ptr: bind_ceed::CeedOperator, 291142270cSJeremy L Thompson _lifeline: PhantomData<&'a ()>, 309df49d7eSJed Brown } 319df49d7eSJed Brown 32c68be7a2SJeremy L Thompson #[derive(Debug)] 339df49d7eSJed Brown pub struct Operator<'a> { 349df49d7eSJed Brown op_core: OperatorCore<'a>, 359df49d7eSJed Brown } 369df49d7eSJed Brown 37c68be7a2SJeremy L Thompson #[derive(Debug)] 389df49d7eSJed Brown pub struct CompositeOperator<'a> { 399df49d7eSJed Brown op_core: OperatorCore<'a>, 409df49d7eSJed Brown } 419df49d7eSJed Brown 429df49d7eSJed Brown // ----------------------------------------------------------------------------- 439df49d7eSJed Brown // Destructor 449df49d7eSJed Brown // ----------------------------------------------------------------------------- 459df49d7eSJed Brown impl<'a> Drop for OperatorCore<'a> { 469df49d7eSJed Brown fn drop(&mut self) { 479df49d7eSJed Brown unsafe { 489df49d7eSJed Brown bind_ceed::CeedOperatorDestroy(&mut self.ptr); 499df49d7eSJed Brown } 509df49d7eSJed Brown } 519df49d7eSJed Brown } 529df49d7eSJed Brown 539df49d7eSJed Brown // ----------------------------------------------------------------------------- 549df49d7eSJed Brown // Display 559df49d7eSJed Brown // ----------------------------------------------------------------------------- 569df49d7eSJed Brown impl<'a> fmt::Display for OperatorCore<'a> { 579df49d7eSJed Brown fn fmt(&self, f: &mut fmt::Formatter) -> fmt::Result { 589df49d7eSJed Brown let mut ptr = std::ptr::null_mut(); 599df49d7eSJed Brown let mut sizeloc = crate::MAX_BUFFER_LENGTH; 609df49d7eSJed Brown let cstring = unsafe { 619df49d7eSJed Brown let file = bind_ceed::open_memstream(&mut ptr, &mut sizeloc); 629df49d7eSJed Brown bind_ceed::CeedOperatorView(self.ptr, file); 639df49d7eSJed Brown bind_ceed::fclose(file); 649df49d7eSJed Brown CString::from_raw(ptr) 659df49d7eSJed Brown }; 669df49d7eSJed Brown cstring.to_string_lossy().fmt(f) 679df49d7eSJed Brown } 689df49d7eSJed Brown } 699df49d7eSJed Brown 709df49d7eSJed Brown /// View an Operator 719df49d7eSJed Brown /// 729df49d7eSJed Brown /// ``` 739df49d7eSJed Brown /// # use libceed::prelude::*; 74c68be7a2SJeremy L Thompson /// # fn main() -> Result<(), libceed::CeedError> { 759df49d7eSJed Brown /// # let ceed = libceed::Ceed::default_init(); 76c68be7a2SJeremy L Thompson /// let qf = ceed.q_function_interior_by_name("Mass1DBuild")?; 779df49d7eSJed Brown /// 789df49d7eSJed Brown /// // Operator field arguments 799df49d7eSJed Brown /// let ne = 3; 809df49d7eSJed Brown /// let q = 4 as usize; 819df49d7eSJed Brown /// let mut ind: Vec<i32> = vec![0; 2 * ne]; 829df49d7eSJed Brown /// for i in 0..ne { 839df49d7eSJed Brown /// ind[2 * i + 0] = i as i32; 849df49d7eSJed Brown /// ind[2 * i + 1] = (i + 1) as i32; 859df49d7eSJed Brown /// } 86c68be7a2SJeremy L Thompson /// let r = ceed.elem_restriction(ne, 2, 1, 1, ne + 1, MemType::Host, &ind)?; 879df49d7eSJed Brown /// let strides: [i32; 3] = [1, q as i32, q as i32]; 88c68be7a2SJeremy L Thompson /// let rq = ceed.strided_elem_restriction(ne, 2, 1, q * ne, strides)?; 899df49d7eSJed Brown /// 90c68be7a2SJeremy L Thompson /// let b = ceed.basis_tensor_H1_Lagrange(1, 1, 2, q, QuadMode::Gauss)?; 919df49d7eSJed Brown /// 929df49d7eSJed Brown /// // Operator fields 939df49d7eSJed Brown /// let op = ceed 94c68be7a2SJeremy L Thompson /// .operator(&qf, QFunctionOpt::None, QFunctionOpt::None)? 95c68be7a2SJeremy L Thompson /// .field("dx", &r, &b, VectorOpt::Active)? 96c68be7a2SJeremy L Thompson /// .field("weights", ElemRestrictionOpt::None, &b, VectorOpt::None)? 97c68be7a2SJeremy L Thompson /// .field("qdata", &rq, BasisOpt::Collocated, VectorOpt::Active)?; 989df49d7eSJed Brown /// 999df49d7eSJed Brown /// println!("{}", op); 100c68be7a2SJeremy L Thompson /// # Ok(()) 101c68be7a2SJeremy L Thompson /// # } 1029df49d7eSJed Brown /// ``` 1039df49d7eSJed Brown impl<'a> fmt::Display for Operator<'a> { 1049df49d7eSJed Brown fn fmt(&self, f: &mut fmt::Formatter) -> fmt::Result { 1059df49d7eSJed Brown self.op_core.fmt(f) 1069df49d7eSJed Brown } 1079df49d7eSJed Brown } 1089df49d7eSJed Brown 1099df49d7eSJed Brown /// View a composite Operator 1109df49d7eSJed Brown /// 1119df49d7eSJed Brown /// ``` 1129df49d7eSJed Brown /// # use libceed::prelude::*; 113c68be7a2SJeremy L Thompson /// # fn main() -> Result<(), libceed::CeedError> { 1149df49d7eSJed Brown /// # let ceed = libceed::Ceed::default_init(); 1159df49d7eSJed Brown /// 1169df49d7eSJed Brown /// // Sub operator field arguments 1179df49d7eSJed Brown /// let ne = 3; 1189df49d7eSJed Brown /// let q = 4 as usize; 1199df49d7eSJed Brown /// let mut ind: Vec<i32> = vec![0; 2 * ne]; 1209df49d7eSJed Brown /// for i in 0..ne { 1219df49d7eSJed Brown /// ind[2 * i + 0] = i as i32; 1229df49d7eSJed Brown /// ind[2 * i + 1] = (i + 1) as i32; 1239df49d7eSJed Brown /// } 124c68be7a2SJeremy L Thompson /// let r = ceed.elem_restriction(ne, 2, 1, 1, ne + 1, MemType::Host, &ind)?; 1259df49d7eSJed Brown /// let strides: [i32; 3] = [1, q as i32, q as i32]; 126c68be7a2SJeremy L Thompson /// let rq = ceed.strided_elem_restriction(ne, 2, 1, q * ne, strides)?; 1279df49d7eSJed Brown /// 128c68be7a2SJeremy L Thompson /// let b = ceed.basis_tensor_H1_Lagrange(1, 1, 2, q, QuadMode::Gauss)?; 1299df49d7eSJed Brown /// 130c68be7a2SJeremy L Thompson /// let qdata_mass = ceed.vector(q * ne)?; 131c68be7a2SJeremy L Thompson /// let qdata_diff = ceed.vector(q * ne)?; 1329df49d7eSJed Brown /// 133c68be7a2SJeremy L Thompson /// let qf_mass = ceed.q_function_interior_by_name("MassApply")?; 1349df49d7eSJed Brown /// let op_mass = ceed 135c68be7a2SJeremy L Thompson /// .operator(&qf_mass, QFunctionOpt::None, QFunctionOpt::None)? 136c68be7a2SJeremy L Thompson /// .field("u", &r, &b, VectorOpt::Active)? 137c68be7a2SJeremy L Thompson /// .field("qdata", &rq, BasisOpt::Collocated, &qdata_mass)? 138c68be7a2SJeremy L Thompson /// .field("v", &r, &b, VectorOpt::Active)?; 1399df49d7eSJed Brown /// 140c68be7a2SJeremy L Thompson /// let qf_diff = ceed.q_function_interior_by_name("Poisson1DApply")?; 1419df49d7eSJed Brown /// let op_diff = ceed 142c68be7a2SJeremy L Thompson /// .operator(&qf_diff, QFunctionOpt::None, QFunctionOpt::None)? 143c68be7a2SJeremy L Thompson /// .field("du", &r, &b, VectorOpt::Active)? 144c68be7a2SJeremy L Thompson /// .field("qdata", &rq, BasisOpt::Collocated, &qdata_diff)? 145c68be7a2SJeremy L Thompson /// .field("dv", &r, &b, VectorOpt::Active)?; 1469df49d7eSJed Brown /// 1479df49d7eSJed Brown /// let op = ceed 148c68be7a2SJeremy L Thompson /// .composite_operator()? 149c68be7a2SJeremy L Thompson /// .sub_operator(&op_mass)? 150c68be7a2SJeremy L Thompson /// .sub_operator(&op_diff)?; 1519df49d7eSJed Brown /// 1529df49d7eSJed Brown /// println!("{}", op); 153c68be7a2SJeremy L Thompson /// # Ok(()) 154c68be7a2SJeremy L Thompson /// # } 1559df49d7eSJed Brown /// ``` 1569df49d7eSJed Brown impl<'a> fmt::Display for CompositeOperator<'a> { 1579df49d7eSJed Brown fn fmt(&self, f: &mut fmt::Formatter) -> fmt::Result { 1589df49d7eSJed Brown self.op_core.fmt(f) 1599df49d7eSJed Brown } 1609df49d7eSJed Brown } 1619df49d7eSJed Brown 1629df49d7eSJed Brown // ----------------------------------------------------------------------------- 1639df49d7eSJed Brown // Core functionality 1649df49d7eSJed Brown // ----------------------------------------------------------------------------- 1659df49d7eSJed Brown impl<'a> OperatorCore<'a> { 1661142270cSJeremy L Thompson // Error handling 1671142270cSJeremy L Thompson #[doc(hidden)] 1681142270cSJeremy L Thompson fn check_error(&self, ierr: i32) -> crate::Result<i32> { 1691142270cSJeremy L Thompson let mut ptr = std::ptr::null_mut(); 1701142270cSJeremy L Thompson unsafe { 1711142270cSJeremy L Thompson bind_ceed::CeedOperatorGetCeed(self.ptr, &mut ptr); 1721142270cSJeremy L Thompson } 1731142270cSJeremy L Thompson crate::check_error(ptr, ierr) 1741142270cSJeremy L Thompson } 1751142270cSJeremy L Thompson 1769df49d7eSJed Brown // Common implementations 1779df49d7eSJed Brown pub fn apply(&self, input: &Vector, output: &mut Vector) -> crate::Result<i32> { 1789df49d7eSJed Brown let ierr = unsafe { 1799df49d7eSJed Brown bind_ceed::CeedOperatorApply( 1809df49d7eSJed Brown self.ptr, 1819df49d7eSJed Brown input.ptr, 1829df49d7eSJed Brown output.ptr, 1839df49d7eSJed Brown bind_ceed::CEED_REQUEST_IMMEDIATE, 1849df49d7eSJed Brown ) 1859df49d7eSJed Brown }; 1861142270cSJeremy L Thompson self.check_error(ierr) 1879df49d7eSJed Brown } 1889df49d7eSJed Brown 1899df49d7eSJed Brown pub fn apply_add(&self, input: &Vector, output: &mut Vector) -> crate::Result<i32> { 1909df49d7eSJed Brown let ierr = unsafe { 1919df49d7eSJed Brown bind_ceed::CeedOperatorApplyAdd( 1929df49d7eSJed Brown self.ptr, 1939df49d7eSJed Brown input.ptr, 1949df49d7eSJed Brown output.ptr, 1959df49d7eSJed Brown bind_ceed::CEED_REQUEST_IMMEDIATE, 1969df49d7eSJed Brown ) 1979df49d7eSJed Brown }; 1981142270cSJeremy L Thompson self.check_error(ierr) 1999df49d7eSJed Brown } 2009df49d7eSJed Brown 2019df49d7eSJed Brown pub fn linear_assemble_diagonal(&self, assembled: &mut Vector) -> crate::Result<i32> { 2029df49d7eSJed Brown let ierr = unsafe { 2039df49d7eSJed Brown bind_ceed::CeedOperatorLinearAssembleDiagonal( 2049df49d7eSJed Brown self.ptr, 2059df49d7eSJed Brown assembled.ptr, 2069df49d7eSJed Brown bind_ceed::CEED_REQUEST_IMMEDIATE, 2079df49d7eSJed Brown ) 2089df49d7eSJed Brown }; 2091142270cSJeremy L Thompson self.check_error(ierr) 2109df49d7eSJed Brown } 2119df49d7eSJed Brown 2129df49d7eSJed Brown pub fn linear_assemble_add_diagonal(&self, assembled: &mut Vector) -> crate::Result<i32> { 2139df49d7eSJed Brown let ierr = unsafe { 2149df49d7eSJed Brown bind_ceed::CeedOperatorLinearAssembleAddDiagonal( 2159df49d7eSJed Brown self.ptr, 2169df49d7eSJed Brown assembled.ptr, 2179df49d7eSJed Brown bind_ceed::CEED_REQUEST_IMMEDIATE, 2189df49d7eSJed Brown ) 2199df49d7eSJed Brown }; 2201142270cSJeremy L Thompson self.check_error(ierr) 2219df49d7eSJed Brown } 2229df49d7eSJed Brown 2239df49d7eSJed Brown pub fn linear_assemble_point_block_diagonal( 2249df49d7eSJed Brown &self, 2259df49d7eSJed Brown assembled: &mut Vector, 2269df49d7eSJed Brown ) -> crate::Result<i32> { 2279df49d7eSJed Brown let ierr = unsafe { 2289df49d7eSJed Brown bind_ceed::CeedOperatorLinearAssemblePointBlockDiagonal( 2299df49d7eSJed Brown self.ptr, 2309df49d7eSJed Brown assembled.ptr, 2319df49d7eSJed Brown bind_ceed::CEED_REQUEST_IMMEDIATE, 2329df49d7eSJed Brown ) 2339df49d7eSJed Brown }; 2341142270cSJeremy L Thompson self.check_error(ierr) 2359df49d7eSJed Brown } 2369df49d7eSJed Brown 2379df49d7eSJed Brown pub fn linear_assemble_add_point_block_diagonal( 2389df49d7eSJed Brown &self, 2399df49d7eSJed Brown assembled: &mut Vector, 2409df49d7eSJed Brown ) -> crate::Result<i32> { 2419df49d7eSJed Brown let ierr = unsafe { 2429df49d7eSJed Brown bind_ceed::CeedOperatorLinearAssembleAddPointBlockDiagonal( 2439df49d7eSJed Brown self.ptr, 2449df49d7eSJed Brown assembled.ptr, 2459df49d7eSJed Brown bind_ceed::CEED_REQUEST_IMMEDIATE, 2469df49d7eSJed Brown ) 2479df49d7eSJed Brown }; 2481142270cSJeremy L Thompson self.check_error(ierr) 2499df49d7eSJed Brown } 2509df49d7eSJed Brown } 2519df49d7eSJed Brown 2529df49d7eSJed Brown // ----------------------------------------------------------------------------- 2539df49d7eSJed Brown // Operator 2549df49d7eSJed Brown // ----------------------------------------------------------------------------- 2559df49d7eSJed Brown impl<'a> Operator<'a> { 2569df49d7eSJed Brown // Constructor 2579df49d7eSJed Brown pub fn create<'b>( 2589df49d7eSJed Brown ceed: &'a crate::Ceed, 2599df49d7eSJed Brown qf: impl Into<QFunctionOpt<'b>>, 2609df49d7eSJed Brown dqf: impl Into<QFunctionOpt<'b>>, 2619df49d7eSJed Brown dqfT: impl Into<QFunctionOpt<'b>>, 2629df49d7eSJed Brown ) -> crate::Result<Self> { 2639df49d7eSJed Brown let mut ptr = std::ptr::null_mut(); 2649df49d7eSJed Brown let ierr = unsafe { 2659df49d7eSJed Brown bind_ceed::CeedOperatorCreate( 2669df49d7eSJed Brown ceed.ptr, 2679df49d7eSJed Brown qf.into().to_raw(), 2689df49d7eSJed Brown dqf.into().to_raw(), 2699df49d7eSJed Brown dqfT.into().to_raw(), 2709df49d7eSJed Brown &mut ptr, 2719df49d7eSJed Brown ) 2729df49d7eSJed Brown }; 2739df49d7eSJed Brown ceed.check_error(ierr)?; 2749df49d7eSJed Brown Ok(Self { 2751142270cSJeremy L Thompson op_core: OperatorCore { 2761142270cSJeremy L Thompson ptr, 2771142270cSJeremy L Thompson _lifeline: PhantomData, 2781142270cSJeremy L Thompson }, 2799df49d7eSJed Brown }) 2809df49d7eSJed Brown } 2819df49d7eSJed Brown 2821142270cSJeremy L Thompson fn from_raw(ptr: bind_ceed::CeedOperator) -> crate::Result<Self> { 2839df49d7eSJed Brown Ok(Self { 2841142270cSJeremy L Thompson op_core: OperatorCore { 2851142270cSJeremy L Thompson ptr, 2861142270cSJeremy L Thompson _lifeline: PhantomData, 2871142270cSJeremy L Thompson }, 2889df49d7eSJed Brown }) 2899df49d7eSJed Brown } 2909df49d7eSJed Brown 2919df49d7eSJed Brown /// Apply Operator to a vector 2929df49d7eSJed Brown /// 2939df49d7eSJed Brown /// * `input` - Input Vector 2949df49d7eSJed Brown /// * `output` - Output Vector 2959df49d7eSJed Brown /// 2969df49d7eSJed Brown /// ``` 2979df49d7eSJed Brown /// # use libceed::prelude::*; 298c68be7a2SJeremy L Thompson /// # fn main() -> Result<(), libceed::CeedError> { 2999df49d7eSJed Brown /// # let ceed = libceed::Ceed::default_init(); 3009df49d7eSJed Brown /// let ne = 4; 3019df49d7eSJed Brown /// let p = 3; 3029df49d7eSJed Brown /// let q = 4; 3039df49d7eSJed Brown /// let ndofs = p * ne - ne + 1; 3049df49d7eSJed Brown /// 3059df49d7eSJed Brown /// // Vectors 306c68be7a2SJeremy L Thompson /// let x = ceed.vector_from_slice(&[-1., -0.5, 0.0, 0.5, 1.0])?; 307c68be7a2SJeremy L Thompson /// let mut qdata = ceed.vector(ne * q)?; 3089df49d7eSJed Brown /// qdata.set_value(0.0); 309c68be7a2SJeremy L Thompson /// let u = ceed.vector_from_slice(&vec![1.0; ndofs])?; 310c68be7a2SJeremy L Thompson /// let mut v = ceed.vector(ndofs)?; 3119df49d7eSJed Brown /// v.set_value(0.0); 3129df49d7eSJed Brown /// 3139df49d7eSJed Brown /// // Restrictions 3149df49d7eSJed Brown /// let mut indx: Vec<i32> = vec![0; 2 * ne]; 3159df49d7eSJed Brown /// for i in 0..ne { 3169df49d7eSJed Brown /// indx[2 * i + 0] = i as i32; 3179df49d7eSJed Brown /// indx[2 * i + 1] = (i + 1) as i32; 3189df49d7eSJed Brown /// } 319c68be7a2SJeremy L Thompson /// let rx = ceed.elem_restriction(ne, 2, 1, 1, ne + 1, MemType::Host, &indx)?; 3209df49d7eSJed Brown /// let mut indu: Vec<i32> = vec![0; p * ne]; 3219df49d7eSJed Brown /// for i in 0..ne { 3229df49d7eSJed Brown /// indu[p * i + 0] = i as i32; 3239df49d7eSJed Brown /// indu[p * i + 1] = (i + 1) as i32; 3249df49d7eSJed Brown /// indu[p * i + 2] = (i + 2) as i32; 3259df49d7eSJed Brown /// } 326c68be7a2SJeremy L Thompson /// let ru = ceed.elem_restriction(ne, 3, 1, 1, ndofs, MemType::Host, &indu)?; 3279df49d7eSJed Brown /// let strides: [i32; 3] = [1, q as i32, q as i32]; 328c68be7a2SJeremy L Thompson /// let rq = ceed.strided_elem_restriction(ne, q, 1, q * ne, strides)?; 3299df49d7eSJed Brown /// 3309df49d7eSJed Brown /// // Bases 331c68be7a2SJeremy L Thompson /// let bx = ceed.basis_tensor_H1_Lagrange(1, 1, 2, q, QuadMode::Gauss)?; 332c68be7a2SJeremy L Thompson /// let bu = ceed.basis_tensor_H1_Lagrange(1, 1, p, q, QuadMode::Gauss)?; 3339df49d7eSJed Brown /// 3349df49d7eSJed Brown /// // Build quadrature data 335c68be7a2SJeremy L Thompson /// let qf_build = ceed.q_function_interior_by_name("Mass1DBuild")?; 336c68be7a2SJeremy L Thompson /// ceed.operator(&qf_build, QFunctionOpt::None, QFunctionOpt::None)? 337c68be7a2SJeremy L Thompson /// .field("dx", &rx, &bx, VectorOpt::Active)? 338c68be7a2SJeremy L Thompson /// .field("weights", ElemRestrictionOpt::None, &bx, VectorOpt::None)? 339c68be7a2SJeremy L Thompson /// .field("qdata", &rq, BasisOpt::Collocated, VectorOpt::Active)? 340c68be7a2SJeremy L Thompson /// .apply(&x, &mut qdata)?; 3419df49d7eSJed Brown /// 3429df49d7eSJed Brown /// // Mass operator 343c68be7a2SJeremy L Thompson /// let qf_mass = ceed.q_function_interior_by_name("MassApply")?; 3449df49d7eSJed Brown /// let op_mass = ceed 345c68be7a2SJeremy L Thompson /// .operator(&qf_mass, QFunctionOpt::None, QFunctionOpt::None)? 346c68be7a2SJeremy L Thompson /// .field("u", &ru, &bu, VectorOpt::Active)? 347c68be7a2SJeremy L Thompson /// .field("qdata", &rq, BasisOpt::Collocated, &qdata)? 348c68be7a2SJeremy L Thompson /// .field("v", &ru, &bu, VectorOpt::Active)?; 3499df49d7eSJed Brown /// 3509df49d7eSJed Brown /// v.set_value(0.0); 351c68be7a2SJeremy L Thompson /// op_mass.apply(&u, &mut v)?; 3529df49d7eSJed Brown /// 3539df49d7eSJed Brown /// // Check 35480a9ef05SNatalie Beams /// let sum: Scalar = v.view().iter().sum(); 3559df49d7eSJed Brown /// assert!( 35680a9ef05SNatalie Beams /// (sum - 2.0).abs() < 10.0 * libceed::EPSILON, 3579df49d7eSJed Brown /// "Incorrect interval length computed" 3589df49d7eSJed Brown /// ); 359c68be7a2SJeremy L Thompson /// # Ok(()) 360c68be7a2SJeremy L Thompson /// # } 3619df49d7eSJed Brown /// ``` 3629df49d7eSJed Brown pub fn apply(&self, input: &Vector, output: &mut Vector) -> crate::Result<i32> { 3639df49d7eSJed Brown self.op_core.apply(input, output) 3649df49d7eSJed Brown } 3659df49d7eSJed Brown 3669df49d7eSJed Brown /// Apply Operator to a vector and add result to output vector 3679df49d7eSJed Brown /// 3689df49d7eSJed Brown /// * `input` - Input Vector 3699df49d7eSJed Brown /// * `output` - Output Vector 3709df49d7eSJed Brown /// 3719df49d7eSJed Brown /// ``` 3729df49d7eSJed Brown /// # use libceed::prelude::*; 373c68be7a2SJeremy L Thompson /// # fn main() -> Result<(), libceed::CeedError> { 3749df49d7eSJed Brown /// # let ceed = libceed::Ceed::default_init(); 3759df49d7eSJed Brown /// let ne = 4; 3769df49d7eSJed Brown /// let p = 3; 3779df49d7eSJed Brown /// let q = 4; 3789df49d7eSJed Brown /// let ndofs = p * ne - ne + 1; 3799df49d7eSJed Brown /// 3809df49d7eSJed Brown /// // Vectors 381c68be7a2SJeremy L Thompson /// let x = ceed.vector_from_slice(&[-1., -0.5, 0.0, 0.5, 1.0])?; 382c68be7a2SJeremy L Thompson /// let mut qdata = ceed.vector(ne * q)?; 3839df49d7eSJed Brown /// qdata.set_value(0.0); 384c68be7a2SJeremy L Thompson /// let u = ceed.vector_from_slice(&vec![1.0; ndofs])?; 385c68be7a2SJeremy L Thompson /// let mut v = ceed.vector(ndofs)?; 3869df49d7eSJed Brown /// 3879df49d7eSJed Brown /// // Restrictions 3889df49d7eSJed Brown /// let mut indx: Vec<i32> = vec![0; 2 * ne]; 3899df49d7eSJed Brown /// for i in 0..ne { 3909df49d7eSJed Brown /// indx[2 * i + 0] = i as i32; 3919df49d7eSJed Brown /// indx[2 * i + 1] = (i + 1) as i32; 3929df49d7eSJed Brown /// } 393c68be7a2SJeremy L Thompson /// let rx = ceed.elem_restriction(ne, 2, 1, 1, ne + 1, MemType::Host, &indx)?; 3949df49d7eSJed Brown /// let mut indu: Vec<i32> = vec![0; p * ne]; 3959df49d7eSJed Brown /// for i in 0..ne { 3969df49d7eSJed Brown /// indu[p * i + 0] = i as i32; 3979df49d7eSJed Brown /// indu[p * i + 1] = (i + 1) as i32; 3989df49d7eSJed Brown /// indu[p * i + 2] = (i + 2) as i32; 3999df49d7eSJed Brown /// } 400c68be7a2SJeremy L Thompson /// let ru = ceed.elem_restriction(ne, 3, 1, 1, ndofs, MemType::Host, &indu)?; 4019df49d7eSJed Brown /// let strides: [i32; 3] = [1, q as i32, q as i32]; 402c68be7a2SJeremy L Thompson /// let rq = ceed.strided_elem_restriction(ne, q, 1, q * ne, strides)?; 4039df49d7eSJed Brown /// 4049df49d7eSJed Brown /// // Bases 405c68be7a2SJeremy L Thompson /// let bx = ceed.basis_tensor_H1_Lagrange(1, 1, 2, q, QuadMode::Gauss)?; 406c68be7a2SJeremy L Thompson /// let bu = ceed.basis_tensor_H1_Lagrange(1, 1, p, q, QuadMode::Gauss)?; 4079df49d7eSJed Brown /// 4089df49d7eSJed Brown /// // Build quadrature data 409c68be7a2SJeremy L Thompson /// let qf_build = ceed.q_function_interior_by_name("Mass1DBuild")?; 410c68be7a2SJeremy L Thompson /// ceed.operator(&qf_build, QFunctionOpt::None, QFunctionOpt::None)? 411c68be7a2SJeremy L Thompson /// .field("dx", &rx, &bx, VectorOpt::Active)? 412c68be7a2SJeremy L Thompson /// .field("weights", ElemRestrictionOpt::None, &bx, VectorOpt::None)? 413c68be7a2SJeremy L Thompson /// .field("qdata", &rq, BasisOpt::Collocated, VectorOpt::Active)? 414c68be7a2SJeremy L Thompson /// .apply(&x, &mut qdata)?; 4159df49d7eSJed Brown /// 4169df49d7eSJed Brown /// // Mass operator 417c68be7a2SJeremy L Thompson /// let qf_mass = ceed.q_function_interior_by_name("MassApply")?; 4189df49d7eSJed Brown /// let op_mass = ceed 419c68be7a2SJeremy L Thompson /// .operator(&qf_mass, QFunctionOpt::None, QFunctionOpt::None)? 420c68be7a2SJeremy L Thompson /// .field("u", &ru, &bu, VectorOpt::Active)? 421c68be7a2SJeremy L Thompson /// .field("qdata", &rq, BasisOpt::Collocated, &qdata)? 422c68be7a2SJeremy L Thompson /// .field("v", &ru, &bu, VectorOpt::Active)?; 4239df49d7eSJed Brown /// 4249df49d7eSJed Brown /// v.set_value(1.0); 425c68be7a2SJeremy L Thompson /// op_mass.apply_add(&u, &mut v)?; 4269df49d7eSJed Brown /// 4279df49d7eSJed Brown /// // Check 42880a9ef05SNatalie Beams /// let sum: Scalar = v.view().iter().sum(); 4299df49d7eSJed Brown /// assert!( 43080a9ef05SNatalie Beams /// (sum - (2.0 + ndofs as Scalar)).abs() < 10.0 * libceed::EPSILON, 4319df49d7eSJed Brown /// "Incorrect interval length computed and added" 4329df49d7eSJed Brown /// ); 433c68be7a2SJeremy L Thompson /// # Ok(()) 434c68be7a2SJeremy L Thompson /// # } 4359df49d7eSJed Brown /// ``` 4369df49d7eSJed Brown pub fn apply_add(&self, input: &Vector, output: &mut Vector) -> crate::Result<i32> { 4379df49d7eSJed Brown self.op_core.apply_add(input, output) 4389df49d7eSJed Brown } 4399df49d7eSJed Brown 4409df49d7eSJed Brown /// Provide a field to a Operator for use by its QFunction 4419df49d7eSJed Brown /// 4429df49d7eSJed Brown /// * `fieldname` - Name of the field (to be matched with the name used by 4439df49d7eSJed Brown /// the QFunction) 4449df49d7eSJed Brown /// * `r` - ElemRestriction 4459df49d7eSJed Brown /// * `b` - Basis in which the field resides or `BasisCollocated` if 4469df49d7eSJed Brown /// collocated with quadrature points 4479df49d7eSJed Brown /// * `v` - Vector to be used by Operator or `VectorActive` if field 4489df49d7eSJed Brown /// is active or `VectorNone` if using `Weight` with the 4499df49d7eSJed Brown /// QFunction 4509df49d7eSJed Brown /// 4519df49d7eSJed Brown /// 4529df49d7eSJed Brown /// ``` 4539df49d7eSJed Brown /// # use libceed::prelude::*; 454c68be7a2SJeremy L Thompson /// # fn main() -> Result<(), libceed::CeedError> { 4559df49d7eSJed Brown /// # let ceed = libceed::Ceed::default_init(); 456c68be7a2SJeremy L Thompson /// let qf = ceed.q_function_interior_by_name("Mass1DBuild")?; 457c68be7a2SJeremy L Thompson /// let mut op = ceed.operator(&qf, QFunctionOpt::None, QFunctionOpt::None)?; 4589df49d7eSJed Brown /// 4599df49d7eSJed Brown /// // Operator field arguments 4609df49d7eSJed Brown /// let ne = 3; 4619df49d7eSJed Brown /// let q = 4; 4629df49d7eSJed Brown /// let mut ind: Vec<i32> = vec![0; 2 * ne]; 4639df49d7eSJed Brown /// for i in 0..ne { 4649df49d7eSJed Brown /// ind[2 * i + 0] = i as i32; 4659df49d7eSJed Brown /// ind[2 * i + 1] = (i + 1) as i32; 4669df49d7eSJed Brown /// } 467c68be7a2SJeremy L Thompson /// let r = ceed.elem_restriction(ne, 2, 1, 1, ne + 1, MemType::Host, &ind)?; 4689df49d7eSJed Brown /// 469c68be7a2SJeremy L Thompson /// let b = ceed.basis_tensor_H1_Lagrange(1, 1, 2, q, QuadMode::Gauss)?; 4709df49d7eSJed Brown /// 4719df49d7eSJed Brown /// // Operator field 472c68be7a2SJeremy L Thompson /// op = op.field("dx", &r, &b, VectorOpt::Active)?; 473c68be7a2SJeremy L Thompson /// # Ok(()) 474c68be7a2SJeremy L Thompson /// # } 4759df49d7eSJed Brown /// ``` 4769df49d7eSJed Brown #[allow(unused_mut)] 4779df49d7eSJed Brown pub fn field<'b>( 4789df49d7eSJed Brown mut self, 4799df49d7eSJed Brown fieldname: &str, 4809df49d7eSJed Brown r: impl Into<ElemRestrictionOpt<'b>>, 4819df49d7eSJed Brown b: impl Into<BasisOpt<'b>>, 4829df49d7eSJed Brown v: impl Into<VectorOpt<'b>>, 4839df49d7eSJed Brown ) -> crate::Result<Self> { 4849df49d7eSJed Brown let fieldname = CString::new(fieldname).expect("CString::new failed"); 4859df49d7eSJed Brown let fieldname = fieldname.as_ptr() as *const i8; 4869df49d7eSJed Brown let ierr = unsafe { 4879df49d7eSJed Brown bind_ceed::CeedOperatorSetField( 4889df49d7eSJed Brown self.op_core.ptr, 4899df49d7eSJed Brown fieldname, 4909df49d7eSJed Brown r.into().to_raw(), 4919df49d7eSJed Brown b.into().to_raw(), 4929df49d7eSJed Brown v.into().to_raw(), 4939df49d7eSJed Brown ) 4949df49d7eSJed Brown }; 4951142270cSJeremy L Thompson self.op_core.check_error(ierr)?; 4969df49d7eSJed Brown Ok(self) 4979df49d7eSJed Brown } 4989df49d7eSJed Brown 499*3f2759cfSJeremy L Thompson /// Get the number of elements associated with an Operator 500*3f2759cfSJeremy L Thompson /// 501*3f2759cfSJeremy L Thompson /// 502*3f2759cfSJeremy L Thompson /// ``` 503*3f2759cfSJeremy L Thompson /// # use libceed::prelude::*; 504*3f2759cfSJeremy L Thompson /// # fn main() -> Result<(), libceed::CeedError> { 505*3f2759cfSJeremy L Thompson /// # let ceed = libceed::Ceed::default_init(); 506*3f2759cfSJeremy L Thompson /// let qf = ceed.q_function_interior_by_name("Mass1DBuild")?; 507*3f2759cfSJeremy L Thompson /// let mut op = ceed.operator(&qf, QFunctionOpt::None, QFunctionOpt::None)?; 508*3f2759cfSJeremy L Thompson /// 509*3f2759cfSJeremy L Thompson /// // Operator field arguments 510*3f2759cfSJeremy L Thompson /// let ne = 3; 511*3f2759cfSJeremy L Thompson /// let q = 4; 512*3f2759cfSJeremy L Thompson /// let mut ind: Vec<i32> = vec![0; 2 * ne]; 513*3f2759cfSJeremy L Thompson /// for i in 0..ne { 514*3f2759cfSJeremy L Thompson /// ind[2 * i + 0] = i as i32; 515*3f2759cfSJeremy L Thompson /// ind[2 * i + 1] = (i + 1) as i32; 516*3f2759cfSJeremy L Thompson /// } 517*3f2759cfSJeremy L Thompson /// let r = ceed.elem_restriction(ne, 2, 1, 1, ne + 1, MemType::Host, &ind)?; 518*3f2759cfSJeremy L Thompson /// 519*3f2759cfSJeremy L Thompson /// let b = ceed.basis_tensor_H1_Lagrange(1, 1, 2, q, QuadMode::Gauss)?; 520*3f2759cfSJeremy L Thompson /// 521*3f2759cfSJeremy L Thompson /// // Operator field 522*3f2759cfSJeremy L Thompson /// op = op.field("dx", &r, &b, VectorOpt::Active)?; 523*3f2759cfSJeremy L Thompson /// 524*3f2759cfSJeremy L Thompson /// // Check number of elements 525*3f2759cfSJeremy L Thompson /// assert_eq!(op.num_elements(), ne, "Incorrect number of elements"); 526*3f2759cfSJeremy L Thompson /// # Ok(()) 527*3f2759cfSJeremy L Thompson /// # } 528*3f2759cfSJeremy L Thompson /// ``` 529*3f2759cfSJeremy L Thompson pub fn num_elements(&self) -> usize { 530*3f2759cfSJeremy L Thompson let mut nelem = 0; 531*3f2759cfSJeremy L Thompson unsafe { bind_ceed::CeedOperatorGetNumElements(self.op_core.ptr, &mut nelem) }; 532*3f2759cfSJeremy L Thompson usize::try_from(nelem).unwrap() 533*3f2759cfSJeremy L Thompson } 534*3f2759cfSJeremy L Thompson 535*3f2759cfSJeremy L Thompson /// Get the number of quadrature points associated with each element of 536*3f2759cfSJeremy L Thompson /// an Operator 537*3f2759cfSJeremy L Thompson /// 538*3f2759cfSJeremy L Thompson /// 539*3f2759cfSJeremy L Thompson /// ``` 540*3f2759cfSJeremy L Thompson /// # use libceed::prelude::*; 541*3f2759cfSJeremy L Thompson /// # fn main() -> Result<(), libceed::CeedError> { 542*3f2759cfSJeremy L Thompson /// # let ceed = libceed::Ceed::default_init(); 543*3f2759cfSJeremy L Thompson /// let qf = ceed.q_function_interior_by_name("Mass1DBuild")?; 544*3f2759cfSJeremy L Thompson /// let mut op = ceed.operator(&qf, QFunctionOpt::None, QFunctionOpt::None)?; 545*3f2759cfSJeremy L Thompson /// 546*3f2759cfSJeremy L Thompson /// // Operator field arguments 547*3f2759cfSJeremy L Thompson /// let ne = 3; 548*3f2759cfSJeremy L Thompson /// let q = 4; 549*3f2759cfSJeremy L Thompson /// let mut ind: Vec<i32> = vec![0; 2 * ne]; 550*3f2759cfSJeremy L Thompson /// for i in 0..ne { 551*3f2759cfSJeremy L Thompson /// ind[2 * i + 0] = i as i32; 552*3f2759cfSJeremy L Thompson /// ind[2 * i + 1] = (i + 1) as i32; 553*3f2759cfSJeremy L Thompson /// } 554*3f2759cfSJeremy L Thompson /// let r = ceed.elem_restriction(ne, 2, 1, 1, ne + 1, MemType::Host, &ind)?; 555*3f2759cfSJeremy L Thompson /// 556*3f2759cfSJeremy L Thompson /// let b = ceed.basis_tensor_H1_Lagrange(1, 1, 2, q, QuadMode::Gauss)?; 557*3f2759cfSJeremy L Thompson /// 558*3f2759cfSJeremy L Thompson /// // Operator field 559*3f2759cfSJeremy L Thompson /// op = op.field("dx", &r, &b, VectorOpt::Active)?; 560*3f2759cfSJeremy L Thompson /// 561*3f2759cfSJeremy L Thompson /// // Check number of quadrature points 562*3f2759cfSJeremy L Thompson /// assert_eq!( 563*3f2759cfSJeremy L Thompson /// op.num_quadrature_points(), 564*3f2759cfSJeremy L Thompson /// q, 565*3f2759cfSJeremy L Thompson /// "Incorrect number of quadrature points" 566*3f2759cfSJeremy L Thompson /// ); 567*3f2759cfSJeremy L Thompson /// # Ok(()) 568*3f2759cfSJeremy L Thompson /// # } 569*3f2759cfSJeremy L Thompson /// ``` 570*3f2759cfSJeremy L Thompson pub fn num_quadrature_points(&self) -> usize { 571*3f2759cfSJeremy L Thompson let mut Q = 0; 572*3f2759cfSJeremy L Thompson unsafe { bind_ceed::CeedOperatorGetNumQuadraturePoints(self.op_core.ptr, &mut Q) }; 573*3f2759cfSJeremy L Thompson usize::try_from(Q).unwrap() 574*3f2759cfSJeremy L Thompson } 575*3f2759cfSJeremy L Thompson 5769df49d7eSJed Brown /// Assemble the diagonal of a square linear Operator 5779df49d7eSJed Brown /// 5789df49d7eSJed Brown /// This overwrites a Vector with the diagonal of a linear Operator. 5799df49d7eSJed Brown /// 5809df49d7eSJed Brown /// Note: Currently only non-composite Operators with a single field and 5819df49d7eSJed Brown /// composite Operators with single field sub-operators are supported. 5829df49d7eSJed Brown /// 5839df49d7eSJed Brown /// * `op` - Operator to assemble QFunction 5849df49d7eSJed Brown /// * `assembled` - Vector to store assembled Operator diagonal 5859df49d7eSJed Brown /// 5869df49d7eSJed Brown /// ``` 5879df49d7eSJed Brown /// # use libceed::prelude::*; 588c68be7a2SJeremy L Thompson /// # fn main() -> Result<(), libceed::CeedError> { 5899df49d7eSJed Brown /// # let ceed = libceed::Ceed::default_init(); 5909df49d7eSJed Brown /// let ne = 4; 5919df49d7eSJed Brown /// let p = 3; 5929df49d7eSJed Brown /// let q = 4; 5939df49d7eSJed Brown /// let ndofs = p * ne - ne + 1; 5949df49d7eSJed Brown /// 5959df49d7eSJed Brown /// // Vectors 596c68be7a2SJeremy L Thompson /// let x = ceed.vector_from_slice(&[-1., -0.5, 0.0, 0.5, 1.0])?; 597c68be7a2SJeremy L Thompson /// let mut qdata = ceed.vector(ne * q)?; 5989df49d7eSJed Brown /// qdata.set_value(0.0); 599c68be7a2SJeremy L Thompson /// let mut u = ceed.vector(ndofs)?; 6009df49d7eSJed Brown /// u.set_value(1.0); 601c68be7a2SJeremy L Thompson /// let mut v = ceed.vector(ndofs)?; 6029df49d7eSJed Brown /// v.set_value(0.0); 6039df49d7eSJed Brown /// 6049df49d7eSJed Brown /// // Restrictions 6059df49d7eSJed Brown /// let mut indx: Vec<i32> = vec![0; 2 * ne]; 6069df49d7eSJed Brown /// for i in 0..ne { 6079df49d7eSJed Brown /// indx[2 * i + 0] = i as i32; 6089df49d7eSJed Brown /// indx[2 * i + 1] = (i + 1) as i32; 6099df49d7eSJed Brown /// } 610c68be7a2SJeremy L Thompson /// let rx = ceed.elem_restriction(ne, 2, 1, 1, ne + 1, MemType::Host, &indx)?; 6119df49d7eSJed Brown /// let mut indu: Vec<i32> = vec![0; p * ne]; 6129df49d7eSJed Brown /// for i in 0..ne { 6139df49d7eSJed Brown /// indu[p * i + 0] = (2 * i) as i32; 6149df49d7eSJed Brown /// indu[p * i + 1] = (2 * i + 1) as i32; 6159df49d7eSJed Brown /// indu[p * i + 2] = (2 * i + 2) as i32; 6169df49d7eSJed Brown /// } 617c68be7a2SJeremy L Thompson /// let ru = ceed.elem_restriction(ne, p, 1, 1, ndofs, MemType::Host, &indu)?; 6189df49d7eSJed Brown /// let strides: [i32; 3] = [1, q as i32, q as i32]; 619c68be7a2SJeremy L Thompson /// let rq = ceed.strided_elem_restriction(ne, q, 1, q * ne, strides)?; 6209df49d7eSJed Brown /// 6219df49d7eSJed Brown /// // Bases 622c68be7a2SJeremy L Thompson /// let bx = ceed.basis_tensor_H1_Lagrange(1, 1, 2, q, QuadMode::Gauss)?; 623c68be7a2SJeremy L Thompson /// let bu = ceed.basis_tensor_H1_Lagrange(1, 1, p, q, QuadMode::Gauss)?; 6249df49d7eSJed Brown /// 6259df49d7eSJed Brown /// // Build quadrature data 626c68be7a2SJeremy L Thompson /// let qf_build = ceed.q_function_interior_by_name("Mass1DBuild")?; 627c68be7a2SJeremy L Thompson /// ceed.operator(&qf_build, QFunctionOpt::None, QFunctionOpt::None)? 628c68be7a2SJeremy L Thompson /// .field("dx", &rx, &bx, VectorOpt::Active)? 629c68be7a2SJeremy L Thompson /// .field("weights", ElemRestrictionOpt::None, &bx, VectorOpt::None)? 630c68be7a2SJeremy L Thompson /// .field("qdata", &rq, BasisOpt::Collocated, VectorOpt::Active)? 631c68be7a2SJeremy L Thompson /// .apply(&x, &mut qdata)?; 6329df49d7eSJed Brown /// 6339df49d7eSJed Brown /// // Mass operator 634c68be7a2SJeremy L Thompson /// let qf_mass = ceed.q_function_interior_by_name("MassApply")?; 6359df49d7eSJed Brown /// let op_mass = ceed 636c68be7a2SJeremy L Thompson /// .operator(&qf_mass, QFunctionOpt::None, QFunctionOpt::None)? 637c68be7a2SJeremy L Thompson /// .field("u", &ru, &bu, VectorOpt::Active)? 638c68be7a2SJeremy L Thompson /// .field("qdata", &rq, BasisOpt::Collocated, &qdata)? 639c68be7a2SJeremy L Thompson /// .field("v", &ru, &bu, VectorOpt::Active)?; 6409df49d7eSJed Brown /// 6419df49d7eSJed Brown /// // Diagonal 642c68be7a2SJeremy L Thompson /// let mut diag = ceed.vector(ndofs)?; 6439df49d7eSJed Brown /// diag.set_value(0.0); 644c68be7a2SJeremy L Thompson /// op_mass.linear_assemble_diagonal(&mut diag)?; 6459df49d7eSJed Brown /// 6469df49d7eSJed Brown /// // Manual diagonal computation 647c68be7a2SJeremy L Thompson /// let mut true_diag = ceed.vector(ndofs)?; 6489df49d7eSJed Brown /// for i in 0..ndofs { 6499df49d7eSJed Brown /// u.set_value(0.0); 6509df49d7eSJed Brown /// { 6519df49d7eSJed Brown /// let mut u_array = u.view_mut(); 6529df49d7eSJed Brown /// u_array[i] = 1.; 6539df49d7eSJed Brown /// } 6549df49d7eSJed Brown /// 655c68be7a2SJeremy L Thompson /// op_mass.apply(&u, &mut v)?; 6569df49d7eSJed Brown /// 6579df49d7eSJed Brown /// { 6589df49d7eSJed Brown /// let v_array = v.view_mut(); 6599df49d7eSJed Brown /// let mut true_array = true_diag.view_mut(); 6609df49d7eSJed Brown /// true_array[i] = v_array[i]; 6619df49d7eSJed Brown /// } 6629df49d7eSJed Brown /// } 6639df49d7eSJed Brown /// 6649df49d7eSJed Brown /// // Check 6659df49d7eSJed Brown /// diag.view() 6669df49d7eSJed Brown /// .iter() 6679df49d7eSJed Brown /// .zip(true_diag.view().iter()) 6689df49d7eSJed Brown /// .for_each(|(computed, actual)| { 6699df49d7eSJed Brown /// assert!( 67080a9ef05SNatalie Beams /// (*computed - *actual).abs() < 10.0 * libceed::EPSILON, 6719df49d7eSJed Brown /// "Diagonal entry incorrect" 6729df49d7eSJed Brown /// ); 6739df49d7eSJed Brown /// }); 674c68be7a2SJeremy L Thompson /// # Ok(()) 675c68be7a2SJeremy L Thompson /// # } 6769df49d7eSJed Brown /// ``` 6779df49d7eSJed Brown pub fn linear_assemble_diagonal(&self, assembled: &mut Vector) -> crate::Result<i32> { 6789df49d7eSJed Brown self.op_core.linear_assemble_diagonal(assembled) 6799df49d7eSJed Brown } 6809df49d7eSJed Brown 6819df49d7eSJed Brown /// Assemble the diagonal of a square linear Operator 6829df49d7eSJed Brown /// 6839df49d7eSJed Brown /// This sums into a Vector with the diagonal of a linear Operator. 6849df49d7eSJed Brown /// 6859df49d7eSJed Brown /// Note: Currently only non-composite Operators with a single field and 6869df49d7eSJed Brown /// composite Operators with single field sub-operators are supported. 6879df49d7eSJed Brown /// 6889df49d7eSJed Brown /// * `op` - Operator to assemble QFunction 6899df49d7eSJed Brown /// * `assembled` - Vector to store assembled Operator diagonal 6909df49d7eSJed Brown /// 6919df49d7eSJed Brown /// 6929df49d7eSJed Brown /// ``` 6939df49d7eSJed Brown /// # use libceed::prelude::*; 694c68be7a2SJeremy L Thompson /// # fn main() -> Result<(), libceed::CeedError> { 6959df49d7eSJed Brown /// # let ceed = libceed::Ceed::default_init(); 6969df49d7eSJed Brown /// let ne = 4; 6979df49d7eSJed Brown /// let p = 3; 6989df49d7eSJed Brown /// let q = 4; 6999df49d7eSJed Brown /// let ndofs = p * ne - ne + 1; 7009df49d7eSJed Brown /// 7019df49d7eSJed Brown /// // Vectors 702c68be7a2SJeremy L Thompson /// let x = ceed.vector_from_slice(&[-1., -0.5, 0.0, 0.5, 1.0])?; 703c68be7a2SJeremy L Thompson /// let mut qdata = ceed.vector(ne * q)?; 7049df49d7eSJed Brown /// qdata.set_value(0.0); 705c68be7a2SJeremy L Thompson /// let mut u = ceed.vector(ndofs)?; 7069df49d7eSJed Brown /// u.set_value(1.0); 707c68be7a2SJeremy L Thompson /// let mut v = ceed.vector(ndofs)?; 7089df49d7eSJed Brown /// v.set_value(0.0); 7099df49d7eSJed Brown /// 7109df49d7eSJed Brown /// // Restrictions 7119df49d7eSJed Brown /// let mut indx: Vec<i32> = vec![0; 2 * ne]; 7129df49d7eSJed Brown /// for i in 0..ne { 7139df49d7eSJed Brown /// indx[2 * i + 0] = i as i32; 7149df49d7eSJed Brown /// indx[2 * i + 1] = (i + 1) as i32; 7159df49d7eSJed Brown /// } 716c68be7a2SJeremy L Thompson /// let rx = ceed.elem_restriction(ne, 2, 1, 1, ne + 1, MemType::Host, &indx)?; 7179df49d7eSJed Brown /// let mut indu: Vec<i32> = vec![0; p * ne]; 7189df49d7eSJed Brown /// for i in 0..ne { 7199df49d7eSJed Brown /// indu[p * i + 0] = (2 * i) as i32; 7209df49d7eSJed Brown /// indu[p * i + 1] = (2 * i + 1) as i32; 7219df49d7eSJed Brown /// indu[p * i + 2] = (2 * i + 2) as i32; 7229df49d7eSJed Brown /// } 723c68be7a2SJeremy L Thompson /// let ru = ceed.elem_restriction(ne, p, 1, 1, ndofs, MemType::Host, &indu)?; 7249df49d7eSJed Brown /// let strides: [i32; 3] = [1, q as i32, q as i32]; 725c68be7a2SJeremy L Thompson /// let rq = ceed.strided_elem_restriction(ne, q, 1, q * ne, strides)?; 7269df49d7eSJed Brown /// 7279df49d7eSJed Brown /// // Bases 728c68be7a2SJeremy L Thompson /// let bx = ceed.basis_tensor_H1_Lagrange(1, 1, 2, q, QuadMode::Gauss)?; 729c68be7a2SJeremy L Thompson /// let bu = ceed.basis_tensor_H1_Lagrange(1, 1, p, q, QuadMode::Gauss)?; 7309df49d7eSJed Brown /// 7319df49d7eSJed Brown /// // Build quadrature data 732c68be7a2SJeremy L Thompson /// let qf_build = ceed.q_function_interior_by_name("Mass1DBuild")?; 733c68be7a2SJeremy L Thompson /// ceed.operator(&qf_build, QFunctionOpt::None, QFunctionOpt::None)? 734c68be7a2SJeremy L Thompson /// .field("dx", &rx, &bx, VectorOpt::Active)? 735c68be7a2SJeremy L Thompson /// .field("weights", ElemRestrictionOpt::None, &bx, VectorOpt::None)? 736c68be7a2SJeremy L Thompson /// .field("qdata", &rq, BasisOpt::Collocated, VectorOpt::Active)? 737c68be7a2SJeremy L Thompson /// .apply(&x, &mut qdata)?; 7389df49d7eSJed Brown /// 7399df49d7eSJed Brown /// // Mass operator 740c68be7a2SJeremy L Thompson /// let qf_mass = ceed.q_function_interior_by_name("MassApply")?; 7419df49d7eSJed Brown /// let op_mass = ceed 742c68be7a2SJeremy L Thompson /// .operator(&qf_mass, QFunctionOpt::None, QFunctionOpt::None)? 743c68be7a2SJeremy L Thompson /// .field("u", &ru, &bu, VectorOpt::Active)? 744c68be7a2SJeremy L Thompson /// .field("qdata", &rq, BasisOpt::Collocated, &qdata)? 745c68be7a2SJeremy L Thompson /// .field("v", &ru, &bu, VectorOpt::Active)?; 7469df49d7eSJed Brown /// 7479df49d7eSJed Brown /// // Diagonal 748c68be7a2SJeremy L Thompson /// let mut diag = ceed.vector(ndofs)?; 7499df49d7eSJed Brown /// diag.set_value(1.0); 750c68be7a2SJeremy L Thompson /// op_mass.linear_assemble_add_diagonal(&mut diag)?; 7519df49d7eSJed Brown /// 7529df49d7eSJed Brown /// // Manual diagonal computation 753c68be7a2SJeremy L Thompson /// let mut true_diag = ceed.vector(ndofs)?; 7549df49d7eSJed Brown /// for i in 0..ndofs { 7559df49d7eSJed Brown /// u.set_value(0.0); 7569df49d7eSJed Brown /// { 7579df49d7eSJed Brown /// let mut u_array = u.view_mut(); 7589df49d7eSJed Brown /// u_array[i] = 1.; 7599df49d7eSJed Brown /// } 7609df49d7eSJed Brown /// 761c68be7a2SJeremy L Thompson /// op_mass.apply(&u, &mut v)?; 7629df49d7eSJed Brown /// 7639df49d7eSJed Brown /// { 7649df49d7eSJed Brown /// let v_array = v.view_mut(); 7659df49d7eSJed Brown /// let mut true_array = true_diag.view_mut(); 7669df49d7eSJed Brown /// true_array[i] = v_array[i] + 1.0; 7679df49d7eSJed Brown /// } 7689df49d7eSJed Brown /// } 7699df49d7eSJed Brown /// 7709df49d7eSJed Brown /// // Check 7719df49d7eSJed Brown /// diag.view() 7729df49d7eSJed Brown /// .iter() 7739df49d7eSJed Brown /// .zip(true_diag.view().iter()) 7749df49d7eSJed Brown /// .for_each(|(computed, actual)| { 7759df49d7eSJed Brown /// assert!( 77680a9ef05SNatalie Beams /// (*computed - *actual).abs() < 10.0 * libceed::EPSILON, 7779df49d7eSJed Brown /// "Diagonal entry incorrect" 7789df49d7eSJed Brown /// ); 7799df49d7eSJed Brown /// }); 780c68be7a2SJeremy L Thompson /// # Ok(()) 781c68be7a2SJeremy L Thompson /// # } 7829df49d7eSJed Brown /// ``` 7839df49d7eSJed Brown pub fn linear_assemble_add_diagonal(&self, assembled: &mut Vector) -> crate::Result<i32> { 7849df49d7eSJed Brown self.op_core.linear_assemble_add_diagonal(assembled) 7859df49d7eSJed Brown } 7869df49d7eSJed Brown 7879df49d7eSJed Brown /// Assemble the point block diagonal of a square linear Operator 7889df49d7eSJed Brown /// 7899df49d7eSJed Brown /// This overwrites a Vector with the point block diagonal of a linear 7909df49d7eSJed Brown /// Operator. 7919df49d7eSJed Brown /// 7929df49d7eSJed Brown /// Note: Currently only non-composite Operators with a single field and 7939df49d7eSJed Brown /// composite Operators with single field sub-operators are supported. 7949df49d7eSJed Brown /// 7959df49d7eSJed Brown /// * `op` - Operator to assemble QFunction 7969df49d7eSJed Brown /// * `assembled` - Vector to store assembled CeedOperator point block 7979df49d7eSJed Brown /// diagonal, provided in row-major form with an 7989df49d7eSJed Brown /// `ncomp * ncomp` block at each node. The dimensions of 7999df49d7eSJed Brown /// this vector are derived from the active vector for 8009df49d7eSJed Brown /// the CeedOperator. The array has shape 8019df49d7eSJed Brown /// `[nodes, component out, component in]`. 8029df49d7eSJed Brown /// 8039df49d7eSJed Brown /// ``` 8049df49d7eSJed Brown /// # use libceed::prelude::*; 805c68be7a2SJeremy L Thompson /// # fn main() -> Result<(), libceed::CeedError> { 8069df49d7eSJed Brown /// # let ceed = libceed::Ceed::default_init(); 8079df49d7eSJed Brown /// let ne = 4; 8089df49d7eSJed Brown /// let p = 3; 8099df49d7eSJed Brown /// let q = 4; 8109df49d7eSJed Brown /// let ncomp = 2; 8119df49d7eSJed Brown /// let ndofs = p * ne - ne + 1; 8129df49d7eSJed Brown /// 8139df49d7eSJed Brown /// // Vectors 814c68be7a2SJeremy L Thompson /// let x = ceed.vector_from_slice(&[-1., -0.5, 0.0, 0.5, 1.0])?; 815c68be7a2SJeremy L Thompson /// let mut qdata = ceed.vector(ne * q)?; 8169df49d7eSJed Brown /// qdata.set_value(0.0); 817c68be7a2SJeremy L Thompson /// let mut u = ceed.vector(ncomp * ndofs)?; 8189df49d7eSJed Brown /// u.set_value(1.0); 819c68be7a2SJeremy L Thompson /// let mut v = ceed.vector(ncomp * ndofs)?; 8209df49d7eSJed Brown /// v.set_value(0.0); 8219df49d7eSJed Brown /// 8229df49d7eSJed Brown /// // Restrictions 8239df49d7eSJed Brown /// let mut indx: Vec<i32> = vec![0; 2 * ne]; 8249df49d7eSJed Brown /// for i in 0..ne { 8259df49d7eSJed Brown /// indx[2 * i + 0] = i as i32; 8269df49d7eSJed Brown /// indx[2 * i + 1] = (i + 1) as i32; 8279df49d7eSJed Brown /// } 828c68be7a2SJeremy L Thompson /// let rx = ceed.elem_restriction(ne, 2, 1, 1, ne + 1, MemType::Host, &indx)?; 8299df49d7eSJed Brown /// let mut indu: Vec<i32> = vec![0; p * ne]; 8309df49d7eSJed Brown /// for i in 0..ne { 8319df49d7eSJed Brown /// indu[p * i + 0] = (2 * i) as i32; 8329df49d7eSJed Brown /// indu[p * i + 1] = (2 * i + 1) as i32; 8339df49d7eSJed Brown /// indu[p * i + 2] = (2 * i + 2) as i32; 8349df49d7eSJed Brown /// } 835c68be7a2SJeremy L Thompson /// let ru = ceed.elem_restriction(ne, p, ncomp, ndofs, ncomp * ndofs, MemType::Host, &indu)?; 8369df49d7eSJed Brown /// let strides: [i32; 3] = [1, q as i32, q as i32]; 837c68be7a2SJeremy L Thompson /// let rq = ceed.strided_elem_restriction(ne, q, 1, q * ne, strides)?; 8389df49d7eSJed Brown /// 8399df49d7eSJed Brown /// // Bases 840c68be7a2SJeremy L Thompson /// let bx = ceed.basis_tensor_H1_Lagrange(1, 1, 2, q, QuadMode::Gauss)?; 841c68be7a2SJeremy L Thompson /// let bu = ceed.basis_tensor_H1_Lagrange(1, ncomp, p, q, QuadMode::Gauss)?; 8429df49d7eSJed Brown /// 8439df49d7eSJed Brown /// // Build quadrature data 844c68be7a2SJeremy L Thompson /// let qf_build = ceed.q_function_interior_by_name("Mass1DBuild")?; 845c68be7a2SJeremy L Thompson /// ceed.operator(&qf_build, QFunctionOpt::None, QFunctionOpt::None)? 846c68be7a2SJeremy L Thompson /// .field("dx", &rx, &bx, VectorOpt::Active)? 847c68be7a2SJeremy L Thompson /// .field("weights", ElemRestrictionOpt::None, &bx, VectorOpt::None)? 848c68be7a2SJeremy L Thompson /// .field("qdata", &rq, BasisOpt::Collocated, VectorOpt::Active)? 849c68be7a2SJeremy L Thompson /// .apply(&x, &mut qdata)?; 8509df49d7eSJed Brown /// 8519df49d7eSJed Brown /// // Mass operator 8529df49d7eSJed Brown /// let mut mass_2_comp = |[u, qdata, ..]: QFunctionInputs, [v, ..]: QFunctionOutputs| { 8539df49d7eSJed Brown /// // Number of quadrature points 8549df49d7eSJed Brown /// let q = qdata.len(); 8559df49d7eSJed Brown /// 8569df49d7eSJed Brown /// // Iterate over quadrature points 8579df49d7eSJed Brown /// for i in 0..q { 8589df49d7eSJed Brown /// v[i + 0 * q] = u[i + 1 * q] * qdata[i]; 8599df49d7eSJed Brown /// v[i + 1 * q] = u[i + 0 * q] * qdata[i]; 8609df49d7eSJed Brown /// } 8619df49d7eSJed Brown /// 8629df49d7eSJed Brown /// // Return clean error code 8639df49d7eSJed Brown /// 0 8649df49d7eSJed Brown /// }; 8659df49d7eSJed Brown /// 8669df49d7eSJed Brown /// let qf_mass = ceed 867c68be7a2SJeremy L Thompson /// .q_function_interior(1, Box::new(mass_2_comp))? 868c68be7a2SJeremy L Thompson /// .input("u", 2, EvalMode::Interp)? 869c68be7a2SJeremy L Thompson /// .input("qdata", 1, EvalMode::None)? 870c68be7a2SJeremy L Thompson /// .output("v", 2, EvalMode::Interp)?; 8719df49d7eSJed Brown /// 8729df49d7eSJed Brown /// let op_mass = ceed 873c68be7a2SJeremy L Thompson /// .operator(&qf_mass, QFunctionOpt::None, QFunctionOpt::None)? 874c68be7a2SJeremy L Thompson /// .field("u", &ru, &bu, VectorOpt::Active)? 875c68be7a2SJeremy L Thompson /// .field("qdata", &rq, BasisOpt::Collocated, &qdata)? 876c68be7a2SJeremy L Thompson /// .field("v", &ru, &bu, VectorOpt::Active)?; 8779df49d7eSJed Brown /// 8789df49d7eSJed Brown /// // Diagonal 879c68be7a2SJeremy L Thompson /// let mut diag = ceed.vector(ncomp * ncomp * ndofs)?; 8809df49d7eSJed Brown /// diag.set_value(0.0); 881c68be7a2SJeremy L Thompson /// op_mass.linear_assemble_point_block_diagonal(&mut diag)?; 8829df49d7eSJed Brown /// 8839df49d7eSJed Brown /// // Manual diagonal computation 884c68be7a2SJeremy L Thompson /// let mut true_diag = ceed.vector(ncomp * ncomp * ndofs)?; 8859df49d7eSJed Brown /// for i in 0..ndofs { 8869df49d7eSJed Brown /// for j in 0..ncomp { 8879df49d7eSJed Brown /// u.set_value(0.0); 8889df49d7eSJed Brown /// { 8899df49d7eSJed Brown /// let mut u_array = u.view_mut(); 8909df49d7eSJed Brown /// u_array[i + j * ndofs] = 1.; 8919df49d7eSJed Brown /// } 8929df49d7eSJed Brown /// 893c68be7a2SJeremy L Thompson /// op_mass.apply(&u, &mut v)?; 8949df49d7eSJed Brown /// 8959df49d7eSJed Brown /// { 8969df49d7eSJed Brown /// let v_array = v.view_mut(); 8979df49d7eSJed Brown /// let mut true_array = true_diag.view_mut(); 8989df49d7eSJed Brown /// for k in 0..ncomp { 8999df49d7eSJed Brown /// true_array[i * ncomp * ncomp + k * ncomp + j] = v_array[i + k * ndofs]; 9009df49d7eSJed Brown /// } 9019df49d7eSJed Brown /// } 9029df49d7eSJed Brown /// } 9039df49d7eSJed Brown /// } 9049df49d7eSJed Brown /// 9059df49d7eSJed Brown /// // Check 9069df49d7eSJed Brown /// diag.view() 9079df49d7eSJed Brown /// .iter() 9089df49d7eSJed Brown /// .zip(true_diag.view().iter()) 9099df49d7eSJed Brown /// .for_each(|(computed, actual)| { 9109df49d7eSJed Brown /// assert!( 91180a9ef05SNatalie Beams /// (*computed - *actual).abs() < 10.0 * libceed::EPSILON, 9129df49d7eSJed Brown /// "Diagonal entry incorrect" 9139df49d7eSJed Brown /// ); 9149df49d7eSJed Brown /// }); 915c68be7a2SJeremy L Thompson /// # Ok(()) 916c68be7a2SJeremy L Thompson /// # } 9179df49d7eSJed Brown /// ``` 9189df49d7eSJed Brown pub fn linear_assemble_point_block_diagonal( 9199df49d7eSJed Brown &self, 9209df49d7eSJed Brown assembled: &mut Vector, 9219df49d7eSJed Brown ) -> crate::Result<i32> { 9229df49d7eSJed Brown self.op_core.linear_assemble_point_block_diagonal(assembled) 9239df49d7eSJed Brown } 9249df49d7eSJed Brown 9259df49d7eSJed Brown /// Assemble the point block diagonal of a square linear Operator 9269df49d7eSJed Brown /// 9279df49d7eSJed Brown /// This sums into a Vector with the point block diagonal of a linear 9289df49d7eSJed Brown /// Operator. 9299df49d7eSJed Brown /// 9309df49d7eSJed Brown /// Note: Currently only non-composite Operators with a single field and 9319df49d7eSJed Brown /// composite Operators with single field sub-operators are supported. 9329df49d7eSJed Brown /// 9339df49d7eSJed Brown /// * `op` - Operator to assemble QFunction 9349df49d7eSJed Brown /// * `assembled` - Vector to store assembled CeedOperator point block 9359df49d7eSJed Brown /// diagonal, provided in row-major form with an 9369df49d7eSJed Brown /// `ncomp * ncomp` block at each node. The dimensions of 9379df49d7eSJed Brown /// this vector are derived from the active vector for 9389df49d7eSJed Brown /// the CeedOperator. The array has shape 9399df49d7eSJed Brown /// `[nodes, component out, component in]`. 9409df49d7eSJed Brown /// 9419df49d7eSJed Brown /// ``` 9429df49d7eSJed Brown /// # use libceed::prelude::*; 943c68be7a2SJeremy L Thompson /// # fn main() -> Result<(), libceed::CeedError> { 9449df49d7eSJed Brown /// # let ceed = libceed::Ceed::default_init(); 9459df49d7eSJed Brown /// let ne = 4; 9469df49d7eSJed Brown /// let p = 3; 9479df49d7eSJed Brown /// let q = 4; 9489df49d7eSJed Brown /// let ncomp = 2; 9499df49d7eSJed Brown /// let ndofs = p * ne - ne + 1; 9509df49d7eSJed Brown /// 9519df49d7eSJed Brown /// // Vectors 952c68be7a2SJeremy L Thompson /// let x = ceed.vector_from_slice(&[-1., -0.5, 0.0, 0.5, 1.0])?; 953c68be7a2SJeremy L Thompson /// let mut qdata = ceed.vector(ne * q)?; 9549df49d7eSJed Brown /// qdata.set_value(0.0); 955c68be7a2SJeremy L Thompson /// let mut u = ceed.vector(ncomp * ndofs)?; 9569df49d7eSJed Brown /// u.set_value(1.0); 957c68be7a2SJeremy L Thompson /// let mut v = ceed.vector(ncomp * ndofs)?; 9589df49d7eSJed Brown /// v.set_value(0.0); 9599df49d7eSJed Brown /// 9609df49d7eSJed Brown /// // Restrictions 9619df49d7eSJed Brown /// let mut indx: Vec<i32> = vec![0; 2 * ne]; 9629df49d7eSJed Brown /// for i in 0..ne { 9639df49d7eSJed Brown /// indx[2 * i + 0] = i as i32; 9649df49d7eSJed Brown /// indx[2 * i + 1] = (i + 1) as i32; 9659df49d7eSJed Brown /// } 966c68be7a2SJeremy L Thompson /// let rx = ceed.elem_restriction(ne, 2, 1, 1, ne + 1, MemType::Host, &indx)?; 9679df49d7eSJed Brown /// let mut indu: Vec<i32> = vec![0; p * ne]; 9689df49d7eSJed Brown /// for i in 0..ne { 9699df49d7eSJed Brown /// indu[p * i + 0] = (2 * i) as i32; 9709df49d7eSJed Brown /// indu[p * i + 1] = (2 * i + 1) as i32; 9719df49d7eSJed Brown /// indu[p * i + 2] = (2 * i + 2) as i32; 9729df49d7eSJed Brown /// } 973c68be7a2SJeremy L Thompson /// let ru = ceed.elem_restriction(ne, p, ncomp, ndofs, ncomp * ndofs, MemType::Host, &indu)?; 9749df49d7eSJed Brown /// let strides: [i32; 3] = [1, q as i32, q as i32]; 975c68be7a2SJeremy L Thompson /// let rq = ceed.strided_elem_restriction(ne, q, 1, q * ne, strides)?; 9769df49d7eSJed Brown /// 9779df49d7eSJed Brown /// // Bases 978c68be7a2SJeremy L Thompson /// let bx = ceed.basis_tensor_H1_Lagrange(1, 1, 2, q, QuadMode::Gauss)?; 979c68be7a2SJeremy L Thompson /// let bu = ceed.basis_tensor_H1_Lagrange(1, ncomp, p, q, QuadMode::Gauss)?; 9809df49d7eSJed Brown /// 9819df49d7eSJed Brown /// // Build quadrature data 982c68be7a2SJeremy L Thompson /// let qf_build = ceed.q_function_interior_by_name("Mass1DBuild")?; 983c68be7a2SJeremy L Thompson /// ceed.operator(&qf_build, QFunctionOpt::None, QFunctionOpt::None)? 984c68be7a2SJeremy L Thompson /// .field("dx", &rx, &bx, VectorOpt::Active)? 985c68be7a2SJeremy L Thompson /// .field("weights", ElemRestrictionOpt::None, &bx, VectorOpt::None)? 986c68be7a2SJeremy L Thompson /// .field("qdata", &rq, BasisOpt::Collocated, VectorOpt::Active)? 987c68be7a2SJeremy L Thompson /// .apply(&x, &mut qdata)?; 9889df49d7eSJed Brown /// 9899df49d7eSJed Brown /// // Mass operator 9909df49d7eSJed Brown /// let mut mass_2_comp = |[u, qdata, ..]: QFunctionInputs, [v, ..]: QFunctionOutputs| { 9919df49d7eSJed Brown /// // Number of quadrature points 9929df49d7eSJed Brown /// let q = qdata.len(); 9939df49d7eSJed Brown /// 9949df49d7eSJed Brown /// // Iterate over quadrature points 9959df49d7eSJed Brown /// for i in 0..q { 9969df49d7eSJed Brown /// v[i + 0 * q] = u[i + 1 * q] * qdata[i]; 9979df49d7eSJed Brown /// v[i + 1 * q] = u[i + 0 * q] * qdata[i]; 9989df49d7eSJed Brown /// } 9999df49d7eSJed Brown /// 10009df49d7eSJed Brown /// // Return clean error code 10019df49d7eSJed Brown /// 0 10029df49d7eSJed Brown /// }; 10039df49d7eSJed Brown /// 10049df49d7eSJed Brown /// let qf_mass = ceed 1005c68be7a2SJeremy L Thompson /// .q_function_interior(1, Box::new(mass_2_comp))? 1006c68be7a2SJeremy L Thompson /// .input("u", 2, EvalMode::Interp)? 1007c68be7a2SJeremy L Thompson /// .input("qdata", 1, EvalMode::None)? 1008c68be7a2SJeremy L Thompson /// .output("v", 2, EvalMode::Interp)?; 10099df49d7eSJed Brown /// 10109df49d7eSJed Brown /// let op_mass = ceed 1011c68be7a2SJeremy L Thompson /// .operator(&qf_mass, QFunctionOpt::None, QFunctionOpt::None)? 1012c68be7a2SJeremy L Thompson /// .field("u", &ru, &bu, VectorOpt::Active)? 1013c68be7a2SJeremy L Thompson /// .field("qdata", &rq, BasisOpt::Collocated, &qdata)? 1014c68be7a2SJeremy L Thompson /// .field("v", &ru, &bu, VectorOpt::Active)?; 10159df49d7eSJed Brown /// 10169df49d7eSJed Brown /// // Diagonal 1017c68be7a2SJeremy L Thompson /// let mut diag = ceed.vector(ncomp * ncomp * ndofs)?; 10189df49d7eSJed Brown /// diag.set_value(1.0); 1019c68be7a2SJeremy L Thompson /// op_mass.linear_assemble_add_point_block_diagonal(&mut diag)?; 10209df49d7eSJed Brown /// 10219df49d7eSJed Brown /// // Manual diagonal computation 1022c68be7a2SJeremy L Thompson /// let mut true_diag = ceed.vector(ncomp * ncomp * ndofs)?; 10239df49d7eSJed Brown /// for i in 0..ndofs { 10249df49d7eSJed Brown /// for j in 0..ncomp { 10259df49d7eSJed Brown /// u.set_value(0.0); 10269df49d7eSJed Brown /// { 10279df49d7eSJed Brown /// let mut u_array = u.view_mut(); 10289df49d7eSJed Brown /// u_array[i + j * ndofs] = 1.; 10299df49d7eSJed Brown /// } 10309df49d7eSJed Brown /// 1031c68be7a2SJeremy L Thompson /// op_mass.apply(&u, &mut v)?; 10329df49d7eSJed Brown /// 10339df49d7eSJed Brown /// { 10349df49d7eSJed Brown /// let v_array = v.view_mut(); 10359df49d7eSJed Brown /// let mut true_array = true_diag.view_mut(); 10369df49d7eSJed Brown /// for k in 0..ncomp { 10379df49d7eSJed Brown /// true_array[i * ncomp * ncomp + k * ncomp + j] = v_array[i + k * ndofs]; 10389df49d7eSJed Brown /// } 10399df49d7eSJed Brown /// } 10409df49d7eSJed Brown /// } 10419df49d7eSJed Brown /// } 10429df49d7eSJed Brown /// 10439df49d7eSJed Brown /// // Check 10449df49d7eSJed Brown /// diag.view() 10459df49d7eSJed Brown /// .iter() 10469df49d7eSJed Brown /// .zip(true_diag.view().iter()) 10479df49d7eSJed Brown /// .for_each(|(computed, actual)| { 10489df49d7eSJed Brown /// assert!( 104980a9ef05SNatalie Beams /// (*computed - 1.0 - *actual).abs() < 10.0 * libceed::EPSILON, 10509df49d7eSJed Brown /// "Diagonal entry incorrect" 10519df49d7eSJed Brown /// ); 10529df49d7eSJed Brown /// }); 1053c68be7a2SJeremy L Thompson /// # Ok(()) 1054c68be7a2SJeremy L Thompson /// # } 10559df49d7eSJed Brown /// ``` 10569df49d7eSJed Brown pub fn linear_assemble_add_point_block_diagonal( 10579df49d7eSJed Brown &self, 10589df49d7eSJed Brown assembled: &mut Vector, 10599df49d7eSJed Brown ) -> crate::Result<i32> { 10609df49d7eSJed Brown self.op_core 10619df49d7eSJed Brown .linear_assemble_add_point_block_diagonal(assembled) 10629df49d7eSJed Brown } 10639df49d7eSJed Brown 10649df49d7eSJed Brown /// Create a multigrid coarse Operator and level transfer Operators for a 10659df49d7eSJed Brown /// given Operator, creating the prolongation basis from the fine and 10669df49d7eSJed Brown /// coarse grid interpolation 10679df49d7eSJed Brown /// 10689df49d7eSJed Brown /// * `p_mult_fine` - Lvector multiplicity in parallel gather/scatter 10699df49d7eSJed Brown /// * `rstr_coarse` - Coarse grid restriction 10709df49d7eSJed Brown /// * `basis_coarse` - Coarse grid active vector basis 10719df49d7eSJed Brown /// 10729df49d7eSJed Brown /// ``` 10739df49d7eSJed Brown /// # use libceed::prelude::*; 1074c68be7a2SJeremy L Thompson /// # fn main() -> Result<(), libceed::CeedError> { 10759df49d7eSJed Brown /// # let ceed = libceed::Ceed::default_init(); 10769df49d7eSJed Brown /// let ne = 15; 10779df49d7eSJed Brown /// let p_coarse = 3; 10789df49d7eSJed Brown /// let p_fine = 5; 10799df49d7eSJed Brown /// let q = 6; 10809df49d7eSJed Brown /// let ndofs_coarse = p_coarse * ne - ne + 1; 10819df49d7eSJed Brown /// let ndofs_fine = p_fine * ne - ne + 1; 10829df49d7eSJed Brown /// 10839df49d7eSJed Brown /// // Vectors 10849df49d7eSJed Brown /// let x_array = (0..ne + 1) 108580a9ef05SNatalie Beams /// .map(|i| 2.0 * i as Scalar / ne as Scalar - 1.0) 108680a9ef05SNatalie Beams /// .collect::<Vec<Scalar>>(); 1087c68be7a2SJeremy L Thompson /// let x = ceed.vector_from_slice(&x_array)?; 1088c68be7a2SJeremy L Thompson /// let mut qdata = ceed.vector(ne * q)?; 10899df49d7eSJed Brown /// qdata.set_value(0.0); 1090c68be7a2SJeremy L Thompson /// let mut u_coarse = ceed.vector(ndofs_coarse)?; 10919df49d7eSJed Brown /// u_coarse.set_value(1.0); 1092c68be7a2SJeremy L Thompson /// let mut u_fine = ceed.vector(ndofs_fine)?; 10939df49d7eSJed Brown /// u_fine.set_value(1.0); 1094c68be7a2SJeremy L Thompson /// let mut v_coarse = ceed.vector(ndofs_coarse)?; 10959df49d7eSJed Brown /// v_coarse.set_value(0.0); 1096c68be7a2SJeremy L Thompson /// let mut v_fine = ceed.vector(ndofs_fine)?; 10979df49d7eSJed Brown /// v_fine.set_value(0.0); 1098c68be7a2SJeremy L Thompson /// let mut multiplicity = ceed.vector(ndofs_fine)?; 10999df49d7eSJed Brown /// multiplicity.set_value(1.0); 11009df49d7eSJed Brown /// 11019df49d7eSJed Brown /// // Restrictions 11029df49d7eSJed Brown /// let strides: [i32; 3] = [1, q as i32, q as i32]; 1103c68be7a2SJeremy L Thompson /// let rq = ceed.strided_elem_restriction(ne, q, 1, q * ne, strides)?; 11049df49d7eSJed Brown /// 11059df49d7eSJed Brown /// let mut indx: Vec<i32> = vec![0; 2 * ne]; 11069df49d7eSJed Brown /// for i in 0..ne { 11079df49d7eSJed Brown /// indx[2 * i + 0] = i as i32; 11089df49d7eSJed Brown /// indx[2 * i + 1] = (i + 1) as i32; 11099df49d7eSJed Brown /// } 1110c68be7a2SJeremy L Thompson /// let rx = ceed.elem_restriction(ne, 2, 1, 1, ne + 1, MemType::Host, &indx)?; 11119df49d7eSJed Brown /// 11129df49d7eSJed Brown /// let mut indu_coarse: Vec<i32> = vec![0; p_coarse * ne]; 11139df49d7eSJed Brown /// for i in 0..ne { 11149df49d7eSJed Brown /// for j in 0..p_coarse { 11159df49d7eSJed Brown /// indu_coarse[p_coarse * i + j] = (i + j) as i32; 11169df49d7eSJed Brown /// } 11179df49d7eSJed Brown /// } 1118c68be7a2SJeremy L Thompson /// let ru_coarse = ceed.elem_restriction( 11199df49d7eSJed Brown /// ne, 11209df49d7eSJed Brown /// p_coarse, 11219df49d7eSJed Brown /// 1, 11229df49d7eSJed Brown /// 1, 11239df49d7eSJed Brown /// ndofs_coarse, 11249df49d7eSJed Brown /// MemType::Host, 11259df49d7eSJed Brown /// &indu_coarse, 1126c68be7a2SJeremy L Thompson /// )?; 11279df49d7eSJed Brown /// 11289df49d7eSJed Brown /// let mut indu_fine: Vec<i32> = vec![0; p_fine * ne]; 11299df49d7eSJed Brown /// for i in 0..ne { 11309df49d7eSJed Brown /// for j in 0..p_fine { 11319df49d7eSJed Brown /// indu_fine[p_fine * i + j] = (i + j) as i32; 11329df49d7eSJed Brown /// } 11339df49d7eSJed Brown /// } 1134c68be7a2SJeremy L Thompson /// let ru_fine = ceed.elem_restriction(ne, p_fine, 1, 1, ndofs_fine, MemType::Host, &indu_fine)?; 11359df49d7eSJed Brown /// 11369df49d7eSJed Brown /// // Bases 1137c68be7a2SJeremy L Thompson /// let bx = ceed.basis_tensor_H1_Lagrange(1, 1, 2, q, QuadMode::Gauss)?; 1138c68be7a2SJeremy L Thompson /// let bu_coarse = ceed.basis_tensor_H1_Lagrange(1, 1, p_coarse, q, QuadMode::Gauss)?; 1139c68be7a2SJeremy L Thompson /// let bu_fine = ceed.basis_tensor_H1_Lagrange(1, 1, p_fine, q, QuadMode::Gauss)?; 11409df49d7eSJed Brown /// 11419df49d7eSJed Brown /// // Build quadrature data 1142c68be7a2SJeremy L Thompson /// let qf_build = ceed.q_function_interior_by_name("Mass1DBuild")?; 1143c68be7a2SJeremy L Thompson /// ceed.operator(&qf_build, QFunctionOpt::None, QFunctionOpt::None)? 1144c68be7a2SJeremy L Thompson /// .field("dx", &rx, &bx, VectorOpt::Active)? 1145c68be7a2SJeremy L Thompson /// .field("weights", ElemRestrictionOpt::None, &bx, VectorOpt::None)? 1146c68be7a2SJeremy L Thompson /// .field("qdata", &rq, BasisOpt::Collocated, VectorOpt::Active)? 1147c68be7a2SJeremy L Thompson /// .apply(&x, &mut qdata)?; 11489df49d7eSJed Brown /// 11499df49d7eSJed Brown /// // Mass operator 1150c68be7a2SJeremy L Thompson /// let qf_mass = ceed.q_function_interior_by_name("MassApply")?; 11519df49d7eSJed Brown /// let op_mass_fine = ceed 1152c68be7a2SJeremy L Thompson /// .operator(&qf_mass, QFunctionOpt::None, QFunctionOpt::None)? 1153c68be7a2SJeremy L Thompson /// .field("u", &ru_fine, &bu_fine, VectorOpt::Active)? 1154c68be7a2SJeremy L Thompson /// .field("qdata", &rq, BasisOpt::Collocated, &qdata)? 1155c68be7a2SJeremy L Thompson /// .field("v", &ru_fine, &bu_fine, VectorOpt::Active)?; 11569df49d7eSJed Brown /// 11579df49d7eSJed Brown /// // Multigrid setup 1158c68be7a2SJeremy L Thompson /// let (op_mass_coarse, op_prolong, op_restrict) = 1159c68be7a2SJeremy L Thompson /// op_mass_fine.create_multigrid_level(&multiplicity, &ru_coarse, &bu_coarse)?; 11609df49d7eSJed Brown /// 11619df49d7eSJed Brown /// // Coarse problem 11629df49d7eSJed Brown /// u_coarse.set_value(1.0); 1163c68be7a2SJeremy L Thompson /// op_mass_coarse.apply(&u_coarse, &mut v_coarse)?; 11649df49d7eSJed Brown /// 11659df49d7eSJed Brown /// // Check 116680a9ef05SNatalie Beams /// let sum: Scalar = v_coarse.view().iter().sum(); 11679df49d7eSJed Brown /// assert!( 116880a9ef05SNatalie Beams /// (sum - 2.0).abs() < 50.0 * libceed::EPSILON, 11699df49d7eSJed Brown /// "Incorrect interval length computed" 11709df49d7eSJed Brown /// ); 11719df49d7eSJed Brown /// 11729df49d7eSJed Brown /// // Prolong 1173c68be7a2SJeremy L Thompson /// op_prolong.apply(&u_coarse, &mut u_fine)?; 11749df49d7eSJed Brown /// 11759df49d7eSJed Brown /// // Fine problem 1176c68be7a2SJeremy L Thompson /// op_mass_fine.apply(&u_fine, &mut v_fine)?; 11779df49d7eSJed Brown /// 11789df49d7eSJed Brown /// // Check 117980a9ef05SNatalie Beams /// let sum: Scalar = v_fine.view().iter().sum(); 11809df49d7eSJed Brown /// assert!( 118180a9ef05SNatalie Beams /// (sum - 2.0).abs() < 50.0 * libceed::EPSILON, 11829df49d7eSJed Brown /// "Incorrect interval length computed" 11839df49d7eSJed Brown /// ); 11849df49d7eSJed Brown /// 11859df49d7eSJed Brown /// // Restrict 1186c68be7a2SJeremy L Thompson /// op_restrict.apply(&v_fine, &mut v_coarse)?; 11879df49d7eSJed Brown /// 11889df49d7eSJed Brown /// // Check 118980a9ef05SNatalie Beams /// let sum: Scalar = v_coarse.view().iter().sum(); 11909df49d7eSJed Brown /// assert!( 119180a9ef05SNatalie Beams /// (sum - 2.0).abs() < 50.0 * libceed::EPSILON, 11929df49d7eSJed Brown /// "Incorrect interval length computed" 11939df49d7eSJed Brown /// ); 1194c68be7a2SJeremy L Thompson /// # Ok(()) 1195c68be7a2SJeremy L Thompson /// # } 11969df49d7eSJed Brown /// ``` 11979df49d7eSJed Brown pub fn create_multigrid_level( 11989df49d7eSJed Brown &self, 11999df49d7eSJed Brown p_mult_fine: &Vector, 12009df49d7eSJed Brown rstr_coarse: &ElemRestriction, 12019df49d7eSJed Brown basis_coarse: &Basis, 12029df49d7eSJed Brown ) -> crate::Result<(Operator, Operator, Operator)> { 12039df49d7eSJed Brown let mut ptr_coarse = std::ptr::null_mut(); 12049df49d7eSJed Brown let mut ptr_prolong = std::ptr::null_mut(); 12059df49d7eSJed Brown let mut ptr_restrict = std::ptr::null_mut(); 12069df49d7eSJed Brown let ierr = unsafe { 12079df49d7eSJed Brown bind_ceed::CeedOperatorMultigridLevelCreate( 12089df49d7eSJed Brown self.op_core.ptr, 12099df49d7eSJed Brown p_mult_fine.ptr, 12109df49d7eSJed Brown rstr_coarse.ptr, 12119df49d7eSJed Brown basis_coarse.ptr, 12129df49d7eSJed Brown &mut ptr_coarse, 12139df49d7eSJed Brown &mut ptr_prolong, 12149df49d7eSJed Brown &mut ptr_restrict, 12159df49d7eSJed Brown ) 12169df49d7eSJed Brown }; 12171142270cSJeremy L Thompson self.op_core.check_error(ierr)?; 12181142270cSJeremy L Thompson let op_coarse = Operator::from_raw(ptr_coarse)?; 12191142270cSJeremy L Thompson let op_prolong = Operator::from_raw(ptr_prolong)?; 12201142270cSJeremy L Thompson let op_restrict = Operator::from_raw(ptr_restrict)?; 12219df49d7eSJed Brown Ok((op_coarse, op_prolong, op_restrict)) 12229df49d7eSJed Brown } 12239df49d7eSJed Brown 12249df49d7eSJed Brown /// Create a multigrid coarse Operator and level transfer Operators for a 12259df49d7eSJed Brown /// given Operator with a tensor basis for the active basis 12269df49d7eSJed Brown /// 12279df49d7eSJed Brown /// * `p_mult_fine` - Lvector multiplicity in parallel gather/scatter 12289df49d7eSJed Brown /// * `rstr_coarse` - Coarse grid restriction 12299df49d7eSJed Brown /// * `basis_coarse` - Coarse grid active vector basis 12309df49d7eSJed Brown /// * `interp_c_to_f` - Matrix for coarse to fine 12319df49d7eSJed Brown /// 12329df49d7eSJed Brown /// ``` 12339df49d7eSJed Brown /// # use libceed::prelude::*; 1234c68be7a2SJeremy L Thompson /// # fn main() -> Result<(), libceed::CeedError> { 12359df49d7eSJed Brown /// # let ceed = libceed::Ceed::default_init(); 12369df49d7eSJed Brown /// let ne = 15; 12379df49d7eSJed Brown /// let p_coarse = 3; 12389df49d7eSJed Brown /// let p_fine = 5; 12399df49d7eSJed Brown /// let q = 6; 12409df49d7eSJed Brown /// let ndofs_coarse = p_coarse * ne - ne + 1; 12419df49d7eSJed Brown /// let ndofs_fine = p_fine * ne - ne + 1; 12429df49d7eSJed Brown /// 12439df49d7eSJed Brown /// // Vectors 12449df49d7eSJed Brown /// let x_array = (0..ne + 1) 124580a9ef05SNatalie Beams /// .map(|i| 2.0 * i as Scalar / ne as Scalar - 1.0) 124680a9ef05SNatalie Beams /// .collect::<Vec<Scalar>>(); 1247c68be7a2SJeremy L Thompson /// let x = ceed.vector_from_slice(&x_array)?; 1248c68be7a2SJeremy L Thompson /// let mut qdata = ceed.vector(ne * q)?; 12499df49d7eSJed Brown /// qdata.set_value(0.0); 1250c68be7a2SJeremy L Thompson /// let mut u_coarse = ceed.vector(ndofs_coarse)?; 12519df49d7eSJed Brown /// u_coarse.set_value(1.0); 1252c68be7a2SJeremy L Thompson /// let mut u_fine = ceed.vector(ndofs_fine)?; 12539df49d7eSJed Brown /// u_fine.set_value(1.0); 1254c68be7a2SJeremy L Thompson /// let mut v_coarse = ceed.vector(ndofs_coarse)?; 12559df49d7eSJed Brown /// v_coarse.set_value(0.0); 1256c68be7a2SJeremy L Thompson /// let mut v_fine = ceed.vector(ndofs_fine)?; 12579df49d7eSJed Brown /// v_fine.set_value(0.0); 1258c68be7a2SJeremy L Thompson /// let mut multiplicity = ceed.vector(ndofs_fine)?; 12599df49d7eSJed Brown /// multiplicity.set_value(1.0); 12609df49d7eSJed Brown /// 12619df49d7eSJed Brown /// // Restrictions 12629df49d7eSJed Brown /// let strides: [i32; 3] = [1, q as i32, q as i32]; 1263c68be7a2SJeremy L Thompson /// let rq = ceed.strided_elem_restriction(ne, q, 1, q * ne, strides)?; 12649df49d7eSJed Brown /// 12659df49d7eSJed Brown /// let mut indx: Vec<i32> = vec![0; 2 * ne]; 12669df49d7eSJed Brown /// for i in 0..ne { 12679df49d7eSJed Brown /// indx[2 * i + 0] = i as i32; 12689df49d7eSJed Brown /// indx[2 * i + 1] = (i + 1) as i32; 12699df49d7eSJed Brown /// } 1270c68be7a2SJeremy L Thompson /// let rx = ceed.elem_restriction(ne, 2, 1, 1, ne + 1, MemType::Host, &indx)?; 12719df49d7eSJed Brown /// 12729df49d7eSJed Brown /// let mut indu_coarse: Vec<i32> = vec![0; p_coarse * ne]; 12739df49d7eSJed Brown /// for i in 0..ne { 12749df49d7eSJed Brown /// for j in 0..p_coarse { 12759df49d7eSJed Brown /// indu_coarse[p_coarse * i + j] = (i + j) as i32; 12769df49d7eSJed Brown /// } 12779df49d7eSJed Brown /// } 1278c68be7a2SJeremy L Thompson /// let ru_coarse = ceed.elem_restriction( 12799df49d7eSJed Brown /// ne, 12809df49d7eSJed Brown /// p_coarse, 12819df49d7eSJed Brown /// 1, 12829df49d7eSJed Brown /// 1, 12839df49d7eSJed Brown /// ndofs_coarse, 12849df49d7eSJed Brown /// MemType::Host, 12859df49d7eSJed Brown /// &indu_coarse, 1286c68be7a2SJeremy L Thompson /// )?; 12879df49d7eSJed Brown /// 12889df49d7eSJed Brown /// let mut indu_fine: Vec<i32> = vec![0; p_fine * ne]; 12899df49d7eSJed Brown /// for i in 0..ne { 12909df49d7eSJed Brown /// for j in 0..p_fine { 12919df49d7eSJed Brown /// indu_fine[p_fine * i + j] = (i + j) as i32; 12929df49d7eSJed Brown /// } 12939df49d7eSJed Brown /// } 1294c68be7a2SJeremy L Thompson /// let ru_fine = ceed.elem_restriction(ne, p_fine, 1, 1, ndofs_fine, MemType::Host, &indu_fine)?; 12959df49d7eSJed Brown /// 12969df49d7eSJed Brown /// // Bases 1297c68be7a2SJeremy L Thompson /// let bx = ceed.basis_tensor_H1_Lagrange(1, 1, 2, q, QuadMode::Gauss)?; 1298c68be7a2SJeremy L Thompson /// let bu_coarse = ceed.basis_tensor_H1_Lagrange(1, 1, p_coarse, q, QuadMode::Gauss)?; 1299c68be7a2SJeremy L Thompson /// let bu_fine = ceed.basis_tensor_H1_Lagrange(1, 1, p_fine, q, QuadMode::Gauss)?; 13009df49d7eSJed Brown /// 13019df49d7eSJed Brown /// // Build quadrature data 1302c68be7a2SJeremy L Thompson /// let qf_build = ceed.q_function_interior_by_name("Mass1DBuild")?; 1303c68be7a2SJeremy L Thompson /// ceed.operator(&qf_build, QFunctionOpt::None, QFunctionOpt::None)? 1304c68be7a2SJeremy L Thompson /// .field("dx", &rx, &bx, VectorOpt::Active)? 1305c68be7a2SJeremy L Thompson /// .field("weights", ElemRestrictionOpt::None, &bx, VectorOpt::None)? 1306c68be7a2SJeremy L Thompson /// .field("qdata", &rq, BasisOpt::Collocated, VectorOpt::Active)? 1307c68be7a2SJeremy L Thompson /// .apply(&x, &mut qdata)?; 13089df49d7eSJed Brown /// 13099df49d7eSJed Brown /// // Mass operator 1310c68be7a2SJeremy L Thompson /// let qf_mass = ceed.q_function_interior_by_name("MassApply")?; 13119df49d7eSJed Brown /// let op_mass_fine = ceed 1312c68be7a2SJeremy L Thompson /// .operator(&qf_mass, QFunctionOpt::None, QFunctionOpt::None)? 1313c68be7a2SJeremy L Thompson /// .field("u", &ru_fine, &bu_fine, VectorOpt::Active)? 1314c68be7a2SJeremy L Thompson /// .field("qdata", &rq, BasisOpt::Collocated, &qdata)? 1315c68be7a2SJeremy L Thompson /// .field("v", &ru_fine, &bu_fine, VectorOpt::Active)?; 13169df49d7eSJed Brown /// 13179df49d7eSJed Brown /// // Multigrid setup 131880a9ef05SNatalie Beams /// let mut interp_c_to_f: Vec<Scalar> = vec![0.; p_coarse * p_fine]; 13199df49d7eSJed Brown /// { 1320c68be7a2SJeremy L Thompson /// let mut coarse = ceed.vector(p_coarse)?; 1321c68be7a2SJeremy L Thompson /// let mut fine = ceed.vector(p_fine)?; 1322c68be7a2SJeremy L Thompson /// let basis_c_to_f = 1323c68be7a2SJeremy L Thompson /// ceed.basis_tensor_H1_Lagrange(1, 1, p_coarse, p_fine, QuadMode::GaussLobatto)?; 13249df49d7eSJed Brown /// for i in 0..p_coarse { 13259df49d7eSJed Brown /// coarse.set_value(0.0); 13269df49d7eSJed Brown /// { 13279df49d7eSJed Brown /// let mut array = coarse.view_mut(); 13289df49d7eSJed Brown /// array[i] = 1.; 13299df49d7eSJed Brown /// } 1330c68be7a2SJeremy L Thompson /// basis_c_to_f.apply( 13319df49d7eSJed Brown /// 1, 13329df49d7eSJed Brown /// TransposeMode::NoTranspose, 13339df49d7eSJed Brown /// EvalMode::Interp, 13349df49d7eSJed Brown /// &coarse, 13359df49d7eSJed Brown /// &mut fine, 1336c68be7a2SJeremy L Thompson /// )?; 13379df49d7eSJed Brown /// let array = fine.view(); 13389df49d7eSJed Brown /// for j in 0..p_fine { 13399df49d7eSJed Brown /// interp_c_to_f[j * p_coarse + i] = array[j]; 13409df49d7eSJed Brown /// } 13419df49d7eSJed Brown /// } 13429df49d7eSJed Brown /// } 1343c68be7a2SJeremy L Thompson /// let (op_mass_coarse, op_prolong, op_restrict) = op_mass_fine.create_multigrid_level_tensor_H1( 1344c68be7a2SJeremy L Thompson /// &multiplicity, 1345c68be7a2SJeremy L Thompson /// &ru_coarse, 1346c68be7a2SJeremy L Thompson /// &bu_coarse, 1347c68be7a2SJeremy L Thompson /// &interp_c_to_f, 1348c68be7a2SJeremy L Thompson /// )?; 13499df49d7eSJed Brown /// 13509df49d7eSJed Brown /// // Coarse problem 13519df49d7eSJed Brown /// u_coarse.set_value(1.0); 1352c68be7a2SJeremy L Thompson /// op_mass_coarse.apply(&u_coarse, &mut v_coarse)?; 13539df49d7eSJed Brown /// 13549df49d7eSJed Brown /// // Check 135580a9ef05SNatalie Beams /// let sum: Scalar = v_coarse.view().iter().sum(); 13569df49d7eSJed Brown /// assert!( 135780a9ef05SNatalie Beams /// (sum - 2.0).abs() < 10.0 * libceed::EPSILON, 13589df49d7eSJed Brown /// "Incorrect interval length computed" 13599df49d7eSJed Brown /// ); 13609df49d7eSJed Brown /// 13619df49d7eSJed Brown /// // Prolong 1362c68be7a2SJeremy L Thompson /// op_prolong.apply(&u_coarse, &mut u_fine)?; 13639df49d7eSJed Brown /// 13649df49d7eSJed Brown /// // Fine problem 1365c68be7a2SJeremy L Thompson /// op_mass_fine.apply(&u_fine, &mut v_fine)?; 13669df49d7eSJed Brown /// 13679df49d7eSJed Brown /// // Check 136880a9ef05SNatalie Beams /// let sum: Scalar = v_fine.view().iter().sum(); 13699df49d7eSJed Brown /// assert!( 137080a9ef05SNatalie Beams /// (sum - 2.0).abs() < 10.0 * libceed::EPSILON, 13719df49d7eSJed Brown /// "Incorrect interval length computed" 13729df49d7eSJed Brown /// ); 13739df49d7eSJed Brown /// 13749df49d7eSJed Brown /// // Restrict 1375c68be7a2SJeremy L Thompson /// op_restrict.apply(&v_fine, &mut v_coarse)?; 13769df49d7eSJed Brown /// 13779df49d7eSJed Brown /// // Check 137880a9ef05SNatalie Beams /// let sum: Scalar = v_coarse.view().iter().sum(); 13799df49d7eSJed Brown /// assert!( 138080a9ef05SNatalie Beams /// (sum - 2.0).abs() < 10.0 * libceed::EPSILON, 13819df49d7eSJed Brown /// "Incorrect interval length computed" 13829df49d7eSJed Brown /// ); 1383c68be7a2SJeremy L Thompson /// # Ok(()) 1384c68be7a2SJeremy L Thompson /// # } 13859df49d7eSJed Brown /// ``` 13869df49d7eSJed Brown pub fn create_multigrid_level_tensor_H1( 13879df49d7eSJed Brown &self, 13889df49d7eSJed Brown p_mult_fine: &Vector, 13899df49d7eSJed Brown rstr_coarse: &ElemRestriction, 13909df49d7eSJed Brown basis_coarse: &Basis, 139180a9ef05SNatalie Beams interpCtoF: &Vec<Scalar>, 13929df49d7eSJed Brown ) -> crate::Result<(Operator, Operator, Operator)> { 13939df49d7eSJed Brown let mut ptr_coarse = std::ptr::null_mut(); 13949df49d7eSJed Brown let mut ptr_prolong = std::ptr::null_mut(); 13959df49d7eSJed Brown let mut ptr_restrict = std::ptr::null_mut(); 13969df49d7eSJed Brown let ierr = unsafe { 13979df49d7eSJed Brown bind_ceed::CeedOperatorMultigridLevelCreateTensorH1( 13989df49d7eSJed Brown self.op_core.ptr, 13999df49d7eSJed Brown p_mult_fine.ptr, 14009df49d7eSJed Brown rstr_coarse.ptr, 14019df49d7eSJed Brown basis_coarse.ptr, 14029df49d7eSJed Brown interpCtoF.as_ptr(), 14039df49d7eSJed Brown &mut ptr_coarse, 14049df49d7eSJed Brown &mut ptr_prolong, 14059df49d7eSJed Brown &mut ptr_restrict, 14069df49d7eSJed Brown ) 14079df49d7eSJed Brown }; 14081142270cSJeremy L Thompson self.op_core.check_error(ierr)?; 14091142270cSJeremy L Thompson let op_coarse = Operator::from_raw(ptr_coarse)?; 14101142270cSJeremy L Thompson let op_prolong = Operator::from_raw(ptr_prolong)?; 14111142270cSJeremy L Thompson let op_restrict = Operator::from_raw(ptr_restrict)?; 14129df49d7eSJed Brown Ok((op_coarse, op_prolong, op_restrict)) 14139df49d7eSJed Brown } 14149df49d7eSJed Brown 14159df49d7eSJed Brown /// Create a multigrid coarse Operator and level transfer Operators for a 14169df49d7eSJed Brown /// given Operator with a non-tensor basis for the active basis 14179df49d7eSJed Brown /// 14189df49d7eSJed Brown /// * `p_mult_fine` - Lvector multiplicity in parallel gather/scatter 14199df49d7eSJed Brown /// * `rstr_coarse` - Coarse grid restriction 14209df49d7eSJed Brown /// * `basis_coarse` - Coarse grid active vector basis 14219df49d7eSJed Brown /// * `interp_c_to_f` - Matrix for coarse to fine 14229df49d7eSJed Brown /// 14239df49d7eSJed Brown /// ``` 14249df49d7eSJed Brown /// # use libceed::prelude::*; 1425c68be7a2SJeremy L Thompson /// # fn main() -> Result<(), libceed::CeedError> { 14269df49d7eSJed Brown /// # let ceed = libceed::Ceed::default_init(); 14279df49d7eSJed Brown /// let ne = 15; 14289df49d7eSJed Brown /// let p_coarse = 3; 14299df49d7eSJed Brown /// let p_fine = 5; 14309df49d7eSJed Brown /// let q = 6; 14319df49d7eSJed Brown /// let ndofs_coarse = p_coarse * ne - ne + 1; 14329df49d7eSJed Brown /// let ndofs_fine = p_fine * ne - ne + 1; 14339df49d7eSJed Brown /// 14349df49d7eSJed Brown /// // Vectors 14359df49d7eSJed Brown /// let x_array = (0..ne + 1) 143680a9ef05SNatalie Beams /// .map(|i| 2.0 * i as Scalar / ne as Scalar - 1.0) 143780a9ef05SNatalie Beams /// .collect::<Vec<Scalar>>(); 1438c68be7a2SJeremy L Thompson /// let x = ceed.vector_from_slice(&x_array)?; 1439c68be7a2SJeremy L Thompson /// let mut qdata = ceed.vector(ne * q)?; 14409df49d7eSJed Brown /// qdata.set_value(0.0); 1441c68be7a2SJeremy L Thompson /// let mut u_coarse = ceed.vector(ndofs_coarse)?; 14429df49d7eSJed Brown /// u_coarse.set_value(1.0); 1443c68be7a2SJeremy L Thompson /// let mut u_fine = ceed.vector(ndofs_fine)?; 14449df49d7eSJed Brown /// u_fine.set_value(1.0); 1445c68be7a2SJeremy L Thompson /// let mut v_coarse = ceed.vector(ndofs_coarse)?; 14469df49d7eSJed Brown /// v_coarse.set_value(0.0); 1447c68be7a2SJeremy L Thompson /// let mut v_fine = ceed.vector(ndofs_fine)?; 14489df49d7eSJed Brown /// v_fine.set_value(0.0); 1449c68be7a2SJeremy L Thompson /// let mut multiplicity = ceed.vector(ndofs_fine)?; 14509df49d7eSJed Brown /// multiplicity.set_value(1.0); 14519df49d7eSJed Brown /// 14529df49d7eSJed Brown /// // Restrictions 14539df49d7eSJed Brown /// let strides: [i32; 3] = [1, q as i32, q as i32]; 1454c68be7a2SJeremy L Thompson /// let rq = ceed.strided_elem_restriction(ne, q, 1, q * ne, strides)?; 14559df49d7eSJed Brown /// 14569df49d7eSJed Brown /// let mut indx: Vec<i32> = vec![0; 2 * ne]; 14579df49d7eSJed Brown /// for i in 0..ne { 14589df49d7eSJed Brown /// indx[2 * i + 0] = i as i32; 14599df49d7eSJed Brown /// indx[2 * i + 1] = (i + 1) as i32; 14609df49d7eSJed Brown /// } 1461c68be7a2SJeremy L Thompson /// let rx = ceed.elem_restriction(ne, 2, 1, 1, ne + 1, MemType::Host, &indx)?; 14629df49d7eSJed Brown /// 14639df49d7eSJed Brown /// let mut indu_coarse: Vec<i32> = vec![0; p_coarse * ne]; 14649df49d7eSJed Brown /// for i in 0..ne { 14659df49d7eSJed Brown /// for j in 0..p_coarse { 14669df49d7eSJed Brown /// indu_coarse[p_coarse * i + j] = (i + j) as i32; 14679df49d7eSJed Brown /// } 14689df49d7eSJed Brown /// } 1469c68be7a2SJeremy L Thompson /// let ru_coarse = ceed.elem_restriction( 14709df49d7eSJed Brown /// ne, 14719df49d7eSJed Brown /// p_coarse, 14729df49d7eSJed Brown /// 1, 14739df49d7eSJed Brown /// 1, 14749df49d7eSJed Brown /// ndofs_coarse, 14759df49d7eSJed Brown /// MemType::Host, 14769df49d7eSJed Brown /// &indu_coarse, 1477c68be7a2SJeremy L Thompson /// )?; 14789df49d7eSJed Brown /// 14799df49d7eSJed Brown /// let mut indu_fine: Vec<i32> = vec![0; p_fine * ne]; 14809df49d7eSJed Brown /// for i in 0..ne { 14819df49d7eSJed Brown /// for j in 0..p_fine { 14829df49d7eSJed Brown /// indu_fine[p_fine * i + j] = (i + j) as i32; 14839df49d7eSJed Brown /// } 14849df49d7eSJed Brown /// } 1485c68be7a2SJeremy L Thompson /// let ru_fine = ceed.elem_restriction(ne, p_fine, 1, 1, ndofs_fine, MemType::Host, &indu_fine)?; 14869df49d7eSJed Brown /// 14879df49d7eSJed Brown /// // Bases 1488c68be7a2SJeremy L Thompson /// let bx = ceed.basis_tensor_H1_Lagrange(1, 1, 2, q, QuadMode::Gauss)?; 1489c68be7a2SJeremy L Thompson /// let bu_coarse = ceed.basis_tensor_H1_Lagrange(1, 1, p_coarse, q, QuadMode::Gauss)?; 1490c68be7a2SJeremy L Thompson /// let bu_fine = ceed.basis_tensor_H1_Lagrange(1, 1, p_fine, q, QuadMode::Gauss)?; 14919df49d7eSJed Brown /// 14929df49d7eSJed Brown /// // Build quadrature data 1493c68be7a2SJeremy L Thompson /// let qf_build = ceed.q_function_interior_by_name("Mass1DBuild")?; 1494c68be7a2SJeremy L Thompson /// ceed.operator(&qf_build, QFunctionOpt::None, QFunctionOpt::None)? 1495c68be7a2SJeremy L Thompson /// .field("dx", &rx, &bx, VectorOpt::Active)? 1496c68be7a2SJeremy L Thompson /// .field("weights", ElemRestrictionOpt::None, &bx, VectorOpt::None)? 1497c68be7a2SJeremy L Thompson /// .field("qdata", &rq, BasisOpt::Collocated, VectorOpt::Active)? 1498c68be7a2SJeremy L Thompson /// .apply(&x, &mut qdata)?; 14999df49d7eSJed Brown /// 15009df49d7eSJed Brown /// // Mass operator 1501c68be7a2SJeremy L Thompson /// let qf_mass = ceed.q_function_interior_by_name("MassApply")?; 15029df49d7eSJed Brown /// let op_mass_fine = ceed 1503c68be7a2SJeremy L Thompson /// .operator(&qf_mass, QFunctionOpt::None, QFunctionOpt::None)? 1504c68be7a2SJeremy L Thompson /// .field("u", &ru_fine, &bu_fine, VectorOpt::Active)? 1505c68be7a2SJeremy L Thompson /// .field("qdata", &rq, BasisOpt::Collocated, &qdata)? 1506c68be7a2SJeremy L Thompson /// .field("v", &ru_fine, &bu_fine, VectorOpt::Active)?; 15079df49d7eSJed Brown /// 15089df49d7eSJed Brown /// // Multigrid setup 150980a9ef05SNatalie Beams /// let mut interp_c_to_f: Vec<Scalar> = vec![0.; p_coarse * p_fine]; 15109df49d7eSJed Brown /// { 1511c68be7a2SJeremy L Thompson /// let mut coarse = ceed.vector(p_coarse)?; 1512c68be7a2SJeremy L Thompson /// let mut fine = ceed.vector(p_fine)?; 1513c68be7a2SJeremy L Thompson /// let basis_c_to_f = 1514c68be7a2SJeremy L Thompson /// ceed.basis_tensor_H1_Lagrange(1, 1, p_coarse, p_fine, QuadMode::GaussLobatto)?; 15159df49d7eSJed Brown /// for i in 0..p_coarse { 15169df49d7eSJed Brown /// coarse.set_value(0.0); 15179df49d7eSJed Brown /// { 15189df49d7eSJed Brown /// let mut array = coarse.view_mut(); 15199df49d7eSJed Brown /// array[i] = 1.; 15209df49d7eSJed Brown /// } 1521c68be7a2SJeremy L Thompson /// basis_c_to_f.apply( 15229df49d7eSJed Brown /// 1, 15239df49d7eSJed Brown /// TransposeMode::NoTranspose, 15249df49d7eSJed Brown /// EvalMode::Interp, 15259df49d7eSJed Brown /// &coarse, 15269df49d7eSJed Brown /// &mut fine, 1527c68be7a2SJeremy L Thompson /// )?; 15289df49d7eSJed Brown /// let array = fine.view(); 15299df49d7eSJed Brown /// for j in 0..p_fine { 15309df49d7eSJed Brown /// interp_c_to_f[j * p_coarse + i] = array[j]; 15319df49d7eSJed Brown /// } 15329df49d7eSJed Brown /// } 15339df49d7eSJed Brown /// } 1534c68be7a2SJeremy L Thompson /// let (op_mass_coarse, op_prolong, op_restrict) = op_mass_fine.create_multigrid_level_H1( 1535c68be7a2SJeremy L Thompson /// &multiplicity, 1536c68be7a2SJeremy L Thompson /// &ru_coarse, 1537c68be7a2SJeremy L Thompson /// &bu_coarse, 1538c68be7a2SJeremy L Thompson /// &interp_c_to_f, 1539c68be7a2SJeremy L Thompson /// )?; 15409df49d7eSJed Brown /// 15419df49d7eSJed Brown /// // Coarse problem 15429df49d7eSJed Brown /// u_coarse.set_value(1.0); 1543c68be7a2SJeremy L Thompson /// op_mass_coarse.apply(&u_coarse, &mut v_coarse)?; 15449df49d7eSJed Brown /// 15459df49d7eSJed Brown /// // Check 154680a9ef05SNatalie Beams /// let sum: Scalar = v_coarse.view().iter().sum(); 15479df49d7eSJed Brown /// assert!( 154880a9ef05SNatalie Beams /// (sum - 2.0).abs() < 10.0 * libceed::EPSILON, 15499df49d7eSJed Brown /// "Incorrect interval length computed" 15509df49d7eSJed Brown /// ); 15519df49d7eSJed Brown /// 15529df49d7eSJed Brown /// // Prolong 1553c68be7a2SJeremy L Thompson /// op_prolong.apply(&u_coarse, &mut u_fine)?; 15549df49d7eSJed Brown /// 15559df49d7eSJed Brown /// // Fine problem 1556c68be7a2SJeremy L Thompson /// op_mass_fine.apply(&u_fine, &mut v_fine)?; 15579df49d7eSJed Brown /// 15589df49d7eSJed Brown /// // Check 155980a9ef05SNatalie Beams /// let sum: Scalar = v_fine.view().iter().sum(); 15609df49d7eSJed Brown /// assert!( 156180a9ef05SNatalie Beams /// (sum - 2.0).abs() < 10.0 * libceed::EPSILON, 15629df49d7eSJed Brown /// "Incorrect interval length computed" 15639df49d7eSJed Brown /// ); 15649df49d7eSJed Brown /// 15659df49d7eSJed Brown /// // Restrict 1566c68be7a2SJeremy L Thompson /// op_restrict.apply(&v_fine, &mut v_coarse)?; 15679df49d7eSJed Brown /// 15689df49d7eSJed Brown /// // Check 156980a9ef05SNatalie Beams /// let sum: Scalar = v_coarse.view().iter().sum(); 15709df49d7eSJed Brown /// assert!( 157180a9ef05SNatalie Beams /// (sum - 2.0).abs() < 10.0 * libceed::EPSILON, 15729df49d7eSJed Brown /// "Incorrect interval length computed" 15739df49d7eSJed Brown /// ); 1574c68be7a2SJeremy L Thompson /// # Ok(()) 1575c68be7a2SJeremy L Thompson /// # } 15769df49d7eSJed Brown /// ``` 15779df49d7eSJed Brown pub fn create_multigrid_level_H1( 15789df49d7eSJed Brown &self, 15799df49d7eSJed Brown p_mult_fine: &Vector, 15809df49d7eSJed Brown rstr_coarse: &ElemRestriction, 15819df49d7eSJed Brown basis_coarse: &Basis, 158280a9ef05SNatalie Beams interpCtoF: &[Scalar], 15839df49d7eSJed Brown ) -> crate::Result<(Operator, Operator, Operator)> { 15849df49d7eSJed Brown let mut ptr_coarse = std::ptr::null_mut(); 15859df49d7eSJed Brown let mut ptr_prolong = std::ptr::null_mut(); 15869df49d7eSJed Brown let mut ptr_restrict = std::ptr::null_mut(); 15879df49d7eSJed Brown let ierr = unsafe { 15889df49d7eSJed Brown bind_ceed::CeedOperatorMultigridLevelCreateH1( 15899df49d7eSJed Brown self.op_core.ptr, 15909df49d7eSJed Brown p_mult_fine.ptr, 15919df49d7eSJed Brown rstr_coarse.ptr, 15929df49d7eSJed Brown basis_coarse.ptr, 15939df49d7eSJed Brown interpCtoF.as_ptr(), 15949df49d7eSJed Brown &mut ptr_coarse, 15959df49d7eSJed Brown &mut ptr_prolong, 15969df49d7eSJed Brown &mut ptr_restrict, 15979df49d7eSJed Brown ) 15989df49d7eSJed Brown }; 15991142270cSJeremy L Thompson self.op_core.check_error(ierr)?; 16001142270cSJeremy L Thompson let op_coarse = Operator::from_raw(ptr_coarse)?; 16011142270cSJeremy L Thompson let op_prolong = Operator::from_raw(ptr_prolong)?; 16021142270cSJeremy L Thompson let op_restrict = Operator::from_raw(ptr_restrict)?; 16039df49d7eSJed Brown Ok((op_coarse, op_prolong, op_restrict)) 16049df49d7eSJed Brown } 16059df49d7eSJed Brown } 16069df49d7eSJed Brown 16079df49d7eSJed Brown // ----------------------------------------------------------------------------- 16089df49d7eSJed Brown // Composite Operator 16099df49d7eSJed Brown // ----------------------------------------------------------------------------- 16109df49d7eSJed Brown impl<'a> CompositeOperator<'a> { 16119df49d7eSJed Brown // Constructor 16129df49d7eSJed Brown pub fn create(ceed: &'a crate::Ceed) -> crate::Result<Self> { 16139df49d7eSJed Brown let mut ptr = std::ptr::null_mut(); 16149df49d7eSJed Brown let ierr = unsafe { bind_ceed::CeedCompositeOperatorCreate(ceed.ptr, &mut ptr) }; 16159df49d7eSJed Brown ceed.check_error(ierr)?; 16169df49d7eSJed Brown Ok(Self { 16171142270cSJeremy L Thompson op_core: OperatorCore { 16181142270cSJeremy L Thompson ptr, 16191142270cSJeremy L Thompson _lifeline: PhantomData, 16201142270cSJeremy L Thompson }, 16219df49d7eSJed Brown }) 16229df49d7eSJed Brown } 16239df49d7eSJed Brown 16249df49d7eSJed Brown /// Apply Operator to a vector 16259df49d7eSJed Brown /// 16269df49d7eSJed Brown /// * `input` - Input Vector 16279df49d7eSJed Brown /// * `output` - Output Vector 16289df49d7eSJed Brown /// 16299df49d7eSJed Brown /// ``` 16309df49d7eSJed Brown /// # use libceed::prelude::*; 1631c68be7a2SJeremy L Thompson /// # fn main() -> Result<(), libceed::CeedError> { 16329df49d7eSJed Brown /// # let ceed = libceed::Ceed::default_init(); 16339df49d7eSJed Brown /// let ne = 4; 16349df49d7eSJed Brown /// let p = 3; 16359df49d7eSJed Brown /// let q = 4; 16369df49d7eSJed Brown /// let ndofs = p * ne - ne + 1; 16379df49d7eSJed Brown /// 16389df49d7eSJed Brown /// // Vectors 1639c68be7a2SJeremy L Thompson /// let x = ceed.vector_from_slice(&[-1., -0.5, 0.0, 0.5, 1.0])?; 1640c68be7a2SJeremy L Thompson /// let mut qdata_mass = ceed.vector(ne * q)?; 16419df49d7eSJed Brown /// qdata_mass.set_value(0.0); 1642c68be7a2SJeremy L Thompson /// let mut qdata_diff = ceed.vector(ne * q)?; 16439df49d7eSJed Brown /// qdata_diff.set_value(0.0); 1644c68be7a2SJeremy L Thompson /// let mut u = ceed.vector(ndofs)?; 16459df49d7eSJed Brown /// u.set_value(1.0); 1646c68be7a2SJeremy L Thompson /// let mut v = ceed.vector(ndofs)?; 16479df49d7eSJed Brown /// v.set_value(0.0); 16489df49d7eSJed Brown /// 16499df49d7eSJed Brown /// // Restrictions 16509df49d7eSJed Brown /// let mut indx: Vec<i32> = vec![0; 2 * ne]; 16519df49d7eSJed Brown /// for i in 0..ne { 16529df49d7eSJed Brown /// indx[2 * i + 0] = i as i32; 16539df49d7eSJed Brown /// indx[2 * i + 1] = (i + 1) as i32; 16549df49d7eSJed Brown /// } 1655c68be7a2SJeremy L Thompson /// let rx = ceed.elem_restriction(ne, 2, 1, 1, ne + 1, MemType::Host, &indx)?; 16569df49d7eSJed Brown /// let mut indu: Vec<i32> = vec![0; p * ne]; 16579df49d7eSJed Brown /// for i in 0..ne { 16589df49d7eSJed Brown /// indu[p * i + 0] = i as i32; 16599df49d7eSJed Brown /// indu[p * i + 1] = (i + 1) as i32; 16609df49d7eSJed Brown /// indu[p * i + 2] = (i + 2) as i32; 16619df49d7eSJed Brown /// } 1662c68be7a2SJeremy L Thompson /// let ru = ceed.elem_restriction(ne, 3, 1, 1, ndofs, MemType::Host, &indu)?; 16639df49d7eSJed Brown /// let strides: [i32; 3] = [1, q as i32, q as i32]; 1664c68be7a2SJeremy L Thompson /// let rq = ceed.strided_elem_restriction(ne, q, 1, q * ne, strides)?; 16659df49d7eSJed Brown /// 16669df49d7eSJed Brown /// // Bases 1667c68be7a2SJeremy L Thompson /// let bx = ceed.basis_tensor_H1_Lagrange(1, 1, 2, q, QuadMode::Gauss)?; 1668c68be7a2SJeremy L Thompson /// let bu = ceed.basis_tensor_H1_Lagrange(1, 1, p, q, QuadMode::Gauss)?; 16699df49d7eSJed Brown /// 16709df49d7eSJed Brown /// // Build quadrature data 1671c68be7a2SJeremy L Thompson /// let qf_build_mass = ceed.q_function_interior_by_name("Mass1DBuild")?; 1672c68be7a2SJeremy L Thompson /// ceed.operator(&qf_build_mass, QFunctionOpt::None, QFunctionOpt::None)? 1673c68be7a2SJeremy L Thompson /// .field("dx", &rx, &bx, VectorOpt::Active)? 1674c68be7a2SJeremy L Thompson /// .field("weights", ElemRestrictionOpt::None, &bx, VectorOpt::None)? 1675c68be7a2SJeremy L Thompson /// .field("qdata", &rq, BasisOpt::Collocated, VectorOpt::Active)? 1676c68be7a2SJeremy L Thompson /// .apply(&x, &mut qdata_mass)?; 16779df49d7eSJed Brown /// 1678c68be7a2SJeremy L Thompson /// let qf_build_diff = ceed.q_function_interior_by_name("Poisson1DBuild")?; 1679c68be7a2SJeremy L Thompson /// ceed.operator(&qf_build_diff, QFunctionOpt::None, QFunctionOpt::None)? 1680c68be7a2SJeremy L Thompson /// .field("dx", &rx, &bx, VectorOpt::Active)? 1681c68be7a2SJeremy L Thompson /// .field("weights", ElemRestrictionOpt::None, &bx, VectorOpt::None)? 1682c68be7a2SJeremy L Thompson /// .field("qdata", &rq, BasisOpt::Collocated, VectorOpt::Active)? 1683c68be7a2SJeremy L Thompson /// .apply(&x, &mut qdata_diff)?; 16849df49d7eSJed Brown /// 16859df49d7eSJed Brown /// // Application operator 1686c68be7a2SJeremy L Thompson /// let qf_mass = ceed.q_function_interior_by_name("MassApply")?; 16879df49d7eSJed Brown /// let op_mass = ceed 1688c68be7a2SJeremy L Thompson /// .operator(&qf_mass, QFunctionOpt::None, QFunctionOpt::None)? 1689c68be7a2SJeremy L Thompson /// .field("u", &ru, &bu, VectorOpt::Active)? 1690c68be7a2SJeremy L Thompson /// .field("qdata", &rq, BasisOpt::Collocated, &qdata_mass)? 1691c68be7a2SJeremy L Thompson /// .field("v", &ru, &bu, VectorOpt::Active)?; 16929df49d7eSJed Brown /// 1693c68be7a2SJeremy L Thompson /// let qf_diff = ceed.q_function_interior_by_name("Poisson1DApply")?; 16949df49d7eSJed Brown /// let op_diff = ceed 1695c68be7a2SJeremy L Thompson /// .operator(&qf_diff, QFunctionOpt::None, QFunctionOpt::None)? 1696c68be7a2SJeremy L Thompson /// .field("du", &ru, &bu, VectorOpt::Active)? 1697c68be7a2SJeremy L Thompson /// .field("qdata", &rq, BasisOpt::Collocated, &qdata_diff)? 1698c68be7a2SJeremy L Thompson /// .field("dv", &ru, &bu, VectorOpt::Active)?; 16999df49d7eSJed Brown /// 17009df49d7eSJed Brown /// let op_composite = ceed 1701c68be7a2SJeremy L Thompson /// .composite_operator()? 1702c68be7a2SJeremy L Thompson /// .sub_operator(&op_mass)? 1703c68be7a2SJeremy L Thompson /// .sub_operator(&op_diff)?; 17049df49d7eSJed Brown /// 17059df49d7eSJed Brown /// v.set_value(0.0); 1706c68be7a2SJeremy L Thompson /// op_composite.apply(&u, &mut v)?; 17079df49d7eSJed Brown /// 17089df49d7eSJed Brown /// // Check 170980a9ef05SNatalie Beams /// let sum: Scalar = v.view().iter().sum(); 17109df49d7eSJed Brown /// assert!( 171180a9ef05SNatalie Beams /// (sum - 2.0).abs() < 10.0 * libceed::EPSILON, 17129df49d7eSJed Brown /// "Incorrect interval length computed" 17139df49d7eSJed Brown /// ); 1714c68be7a2SJeremy L Thompson /// # Ok(()) 1715c68be7a2SJeremy L Thompson /// # } 17169df49d7eSJed Brown /// ``` 17179df49d7eSJed Brown pub fn apply(&self, input: &Vector, output: &mut Vector) -> crate::Result<i32> { 17189df49d7eSJed Brown self.op_core.apply(input, output) 17199df49d7eSJed Brown } 17209df49d7eSJed Brown 17219df49d7eSJed Brown /// Apply Operator to a vector and add result to output vector 17229df49d7eSJed Brown /// 17239df49d7eSJed Brown /// * `input` - Input Vector 17249df49d7eSJed Brown /// * `output` - Output Vector 17259df49d7eSJed Brown /// 17269df49d7eSJed Brown /// ``` 17279df49d7eSJed Brown /// # use libceed::prelude::*; 1728c68be7a2SJeremy L Thompson /// # fn main() -> Result<(), libceed::CeedError> { 17299df49d7eSJed Brown /// # let ceed = libceed::Ceed::default_init(); 17309df49d7eSJed Brown /// let ne = 4; 17319df49d7eSJed Brown /// let p = 3; 17329df49d7eSJed Brown /// let q = 4; 17339df49d7eSJed Brown /// let ndofs = p * ne - ne + 1; 17349df49d7eSJed Brown /// 17359df49d7eSJed Brown /// // Vectors 1736c68be7a2SJeremy L Thompson /// let x = ceed.vector_from_slice(&[-1., -0.5, 0.0, 0.5, 1.0])?; 1737c68be7a2SJeremy L Thompson /// let mut qdata_mass = ceed.vector(ne * q)?; 17389df49d7eSJed Brown /// qdata_mass.set_value(0.0); 1739c68be7a2SJeremy L Thompson /// let mut qdata_diff = ceed.vector(ne * q)?; 17409df49d7eSJed Brown /// qdata_diff.set_value(0.0); 1741c68be7a2SJeremy L Thompson /// let mut u = ceed.vector(ndofs)?; 17429df49d7eSJed Brown /// u.set_value(1.0); 1743c68be7a2SJeremy L Thompson /// let mut v = ceed.vector(ndofs)?; 17449df49d7eSJed Brown /// v.set_value(0.0); 17459df49d7eSJed Brown /// 17469df49d7eSJed Brown /// // Restrictions 17479df49d7eSJed Brown /// let mut indx: Vec<i32> = vec![0; 2 * ne]; 17489df49d7eSJed Brown /// for i in 0..ne { 17499df49d7eSJed Brown /// indx[2 * i + 0] = i as i32; 17509df49d7eSJed Brown /// indx[2 * i + 1] = (i + 1) as i32; 17519df49d7eSJed Brown /// } 1752c68be7a2SJeremy L Thompson /// let rx = ceed.elem_restriction(ne, 2, 1, 1, ne + 1, MemType::Host, &indx)?; 17539df49d7eSJed Brown /// let mut indu: Vec<i32> = vec![0; p * ne]; 17549df49d7eSJed Brown /// for i in 0..ne { 17559df49d7eSJed Brown /// indu[p * i + 0] = i as i32; 17569df49d7eSJed Brown /// indu[p * i + 1] = (i + 1) as i32; 17579df49d7eSJed Brown /// indu[p * i + 2] = (i + 2) as i32; 17589df49d7eSJed Brown /// } 1759c68be7a2SJeremy L Thompson /// let ru = ceed.elem_restriction(ne, 3, 1, 1, ndofs, MemType::Host, &indu)?; 17609df49d7eSJed Brown /// let strides: [i32; 3] = [1, q as i32, q as i32]; 1761c68be7a2SJeremy L Thompson /// let rq = ceed.strided_elem_restriction(ne, q, 1, q * ne, strides)?; 17629df49d7eSJed Brown /// 17639df49d7eSJed Brown /// // Bases 1764c68be7a2SJeremy L Thompson /// let bx = ceed.basis_tensor_H1_Lagrange(1, 1, 2, q, QuadMode::Gauss)?; 1765c68be7a2SJeremy L Thompson /// let bu = ceed.basis_tensor_H1_Lagrange(1, 1, p, q, QuadMode::Gauss)?; 17669df49d7eSJed Brown /// 17679df49d7eSJed Brown /// // Build quadrature data 1768c68be7a2SJeremy L Thompson /// let qf_build_mass = ceed.q_function_interior_by_name("Mass1DBuild")?; 1769c68be7a2SJeremy L Thompson /// ceed.operator(&qf_build_mass, QFunctionOpt::None, QFunctionOpt::None)? 1770c68be7a2SJeremy L Thompson /// .field("dx", &rx, &bx, VectorOpt::Active)? 1771c68be7a2SJeremy L Thompson /// .field("weights", ElemRestrictionOpt::None, &bx, VectorOpt::None)? 1772c68be7a2SJeremy L Thompson /// .field("qdata", &rq, BasisOpt::Collocated, VectorOpt::Active)? 1773c68be7a2SJeremy L Thompson /// .apply(&x, &mut qdata_mass)?; 17749df49d7eSJed Brown /// 1775c68be7a2SJeremy L Thompson /// let qf_build_diff = ceed.q_function_interior_by_name("Poisson1DBuild")?; 1776c68be7a2SJeremy L Thompson /// ceed.operator(&qf_build_diff, QFunctionOpt::None, QFunctionOpt::None)? 1777c68be7a2SJeremy L Thompson /// .field("dx", &rx, &bx, VectorOpt::Active)? 1778c68be7a2SJeremy L Thompson /// .field("weights", ElemRestrictionOpt::None, &bx, VectorOpt::None)? 1779c68be7a2SJeremy L Thompson /// .field("qdata", &rq, BasisOpt::Collocated, VectorOpt::Active)? 1780c68be7a2SJeremy L Thompson /// .apply(&x, &mut qdata_diff)?; 17819df49d7eSJed Brown /// 17829df49d7eSJed Brown /// // Application operator 1783c68be7a2SJeremy L Thompson /// let qf_mass = ceed.q_function_interior_by_name("MassApply")?; 17849df49d7eSJed Brown /// let op_mass = ceed 1785c68be7a2SJeremy L Thompson /// .operator(&qf_mass, QFunctionOpt::None, QFunctionOpt::None)? 1786c68be7a2SJeremy L Thompson /// .field("u", &ru, &bu, VectorOpt::Active)? 1787c68be7a2SJeremy L Thompson /// .field("qdata", &rq, BasisOpt::Collocated, &qdata_mass)? 1788c68be7a2SJeremy L Thompson /// .field("v", &ru, &bu, VectorOpt::Active)?; 17899df49d7eSJed Brown /// 1790c68be7a2SJeremy L Thompson /// let qf_diff = ceed.q_function_interior_by_name("Poisson1DApply")?; 17919df49d7eSJed Brown /// let op_diff = ceed 1792c68be7a2SJeremy L Thompson /// .operator(&qf_diff, QFunctionOpt::None, QFunctionOpt::None)? 1793c68be7a2SJeremy L Thompson /// .field("du", &ru, &bu, VectorOpt::Active)? 1794c68be7a2SJeremy L Thompson /// .field("qdata", &rq, BasisOpt::Collocated, &qdata_diff)? 1795c68be7a2SJeremy L Thompson /// .field("dv", &ru, &bu, VectorOpt::Active)?; 17969df49d7eSJed Brown /// 17979df49d7eSJed Brown /// let op_composite = ceed 1798c68be7a2SJeremy L Thompson /// .composite_operator()? 1799c68be7a2SJeremy L Thompson /// .sub_operator(&op_mass)? 1800c68be7a2SJeremy L Thompson /// .sub_operator(&op_diff)?; 18019df49d7eSJed Brown /// 18029df49d7eSJed Brown /// v.set_value(1.0); 1803c68be7a2SJeremy L Thompson /// op_composite.apply_add(&u, &mut v)?; 18049df49d7eSJed Brown /// 18059df49d7eSJed Brown /// // Check 180680a9ef05SNatalie Beams /// let sum: Scalar = v.view().iter().sum(); 18079df49d7eSJed Brown /// assert!( 180880a9ef05SNatalie Beams /// (sum - (2.0 + ndofs as Scalar)).abs() < 10.0 * libceed::EPSILON, 18099df49d7eSJed Brown /// "Incorrect interval length computed" 18109df49d7eSJed Brown /// ); 1811c68be7a2SJeremy L Thompson /// # Ok(()) 1812c68be7a2SJeremy L Thompson /// # } 18139df49d7eSJed Brown /// ``` 18149df49d7eSJed Brown pub fn apply_add(&self, input: &Vector, output: &mut Vector) -> crate::Result<i32> { 18159df49d7eSJed Brown self.op_core.apply_add(input, output) 18169df49d7eSJed Brown } 18179df49d7eSJed Brown 18189df49d7eSJed Brown /// Add a sub-Operator to a Composite Operator 18199df49d7eSJed Brown /// 18209df49d7eSJed Brown /// * `subop` - Sub-Operator 18219df49d7eSJed Brown /// 18229df49d7eSJed Brown /// ``` 18239df49d7eSJed Brown /// # use libceed::prelude::*; 1824c68be7a2SJeremy L Thompson /// # fn main() -> Result<(), libceed::CeedError> { 18259df49d7eSJed Brown /// # let ceed = libceed::Ceed::default_init(); 1826c68be7a2SJeremy L Thompson /// let mut op = ceed.composite_operator()?; 18279df49d7eSJed Brown /// 1828c68be7a2SJeremy L Thompson /// let qf_mass = ceed.q_function_interior_by_name("MassApply")?; 1829c68be7a2SJeremy L Thompson /// let op_mass = ceed.operator(&qf_mass, QFunctionOpt::None, QFunctionOpt::None)?; 1830c68be7a2SJeremy L Thompson /// op = op.sub_operator(&op_mass)?; 18319df49d7eSJed Brown /// 1832c68be7a2SJeremy L Thompson /// let qf_diff = ceed.q_function_interior_by_name("Poisson1DApply")?; 1833c68be7a2SJeremy L Thompson /// let op_diff = ceed.operator(&qf_diff, QFunctionOpt::None, QFunctionOpt::None)?; 1834c68be7a2SJeremy L Thompson /// op = op.sub_operator(&op_diff)?; 1835c68be7a2SJeremy L Thompson /// # Ok(()) 1836c68be7a2SJeremy L Thompson /// # } 18379df49d7eSJed Brown /// ``` 18389df49d7eSJed Brown #[allow(unused_mut)] 18399df49d7eSJed Brown pub fn sub_operator(mut self, subop: &Operator) -> crate::Result<Self> { 18409df49d7eSJed Brown let ierr = 18419df49d7eSJed Brown unsafe { bind_ceed::CeedCompositeOperatorAddSub(self.op_core.ptr, subop.op_core.ptr) }; 18421142270cSJeremy L Thompson self.op_core.check_error(ierr)?; 18439df49d7eSJed Brown Ok(self) 18449df49d7eSJed Brown } 18459df49d7eSJed Brown 18469df49d7eSJed Brown /// Assemble the diagonal of a square linear Operator 18479df49d7eSJed Brown /// 18489df49d7eSJed Brown /// This overwrites a Vector with the diagonal of a linear Operator. 18499df49d7eSJed Brown /// 18509df49d7eSJed Brown /// Note: Currently only non-composite Operators with a single field and 18519df49d7eSJed Brown /// composite Operators with single field sub-operators are supported. 18529df49d7eSJed Brown /// 18539df49d7eSJed Brown /// * `op` - Operator to assemble QFunction 18549df49d7eSJed Brown /// * `assembled` - Vector to store assembled Operator diagonal 18559df49d7eSJed Brown pub fn linear_asssemble_diagonal(&self, assembled: &mut Vector) -> crate::Result<i32> { 18569df49d7eSJed Brown self.op_core.linear_assemble_add_diagonal(assembled) 18579df49d7eSJed Brown } 18589df49d7eSJed Brown 18599df49d7eSJed Brown /// Assemble the point block diagonal of a square linear Operator 18609df49d7eSJed Brown /// 18619df49d7eSJed Brown /// This overwrites a Vector with the point block diagonal of a linear 18629df49d7eSJed Brown /// Operator. 18639df49d7eSJed Brown /// 18649df49d7eSJed Brown /// Note: Currently only non-composite Operators with a single field and 18659df49d7eSJed Brown /// composite Operators with single field sub-operators are supported. 18669df49d7eSJed Brown /// 18679df49d7eSJed Brown /// * `op` - Operator to assemble QFunction 18689df49d7eSJed Brown /// * `assembled` - Vector to store assembled Operator diagonal 18699df49d7eSJed Brown pub fn linear_assemble_add_diagonal(&self, assembled: &mut Vector) -> crate::Result<i32> { 18709df49d7eSJed Brown self.op_core.linear_assemble_add_diagonal(assembled) 18719df49d7eSJed Brown } 18729df49d7eSJed Brown 18739df49d7eSJed Brown /// Assemble the diagonal of a square linear Operator 18749df49d7eSJed Brown /// 18759df49d7eSJed Brown /// This overwrites a Vector with the diagonal of a linear Operator. 18769df49d7eSJed Brown /// 18779df49d7eSJed Brown /// Note: Currently only non-composite Operators with a single field and 18789df49d7eSJed Brown /// composite Operators with single field sub-operators are supported. 18799df49d7eSJed Brown /// 18809df49d7eSJed Brown /// * `op` - Operator to assemble QFunction 18819df49d7eSJed Brown /// * `assembled` - Vector to store assembled CeedOperator point block 18829df49d7eSJed Brown /// diagonal, provided in row-major form with an 18839df49d7eSJed Brown /// `ncomp * ncomp` block at each node. The dimensions of 18849df49d7eSJed Brown /// this vector are derived from the active vector for 18859df49d7eSJed Brown /// the CeedOperator. The array has shape 18869df49d7eSJed Brown /// `[nodes, component out, component in]`. 18879df49d7eSJed Brown pub fn linear_assemble_point_block_diagonal( 18889df49d7eSJed Brown &self, 18899df49d7eSJed Brown assembled: &mut Vector, 18909df49d7eSJed Brown ) -> crate::Result<i32> { 18919df49d7eSJed Brown self.op_core.linear_assemble_add_diagonal(assembled) 18929df49d7eSJed Brown } 18939df49d7eSJed Brown 18949df49d7eSJed Brown /// Assemble the diagonal of a square linear Operator 18959df49d7eSJed Brown /// 18969df49d7eSJed Brown /// This sums into a Vector with the diagonal of a linear Operator. 18979df49d7eSJed Brown /// 18989df49d7eSJed Brown /// Note: Currently only non-composite Operators with a single field and 18999df49d7eSJed Brown /// composite Operators with single field sub-operators are supported. 19009df49d7eSJed Brown /// 19019df49d7eSJed Brown /// * `op` - Operator to assemble QFunction 19029df49d7eSJed Brown /// * `assembled` - Vector to store assembled CeedOperator point block 19039df49d7eSJed Brown /// diagonal, provided in row-major form with an 19049df49d7eSJed Brown /// `ncomp * ncomp` block at each node. The dimensions of 19059df49d7eSJed Brown /// this vector are derived from the active vector for 19069df49d7eSJed Brown /// the CeedOperator. The array has shape 19079df49d7eSJed Brown /// `[nodes, component out, component in]`. 19089df49d7eSJed Brown pub fn linear_assemble_add_point_block_diagonal( 19099df49d7eSJed Brown &self, 19109df49d7eSJed Brown assembled: &mut Vector, 19119df49d7eSJed Brown ) -> crate::Result<i32> { 19129df49d7eSJed Brown self.op_core.linear_assemble_add_diagonal(assembled) 19139df49d7eSJed Brown } 19149df49d7eSJed Brown } 19159df49d7eSJed Brown 19169df49d7eSJed Brown // ----------------------------------------------------------------------------- 1917