xref: /libCEED/rust/libceed/src/operator.rs (revision 80a9ef0545a39c00cdcaab1ca26f8053604f3120)
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 // -----------------------------------------------------------------------------
269df49d7eSJed Brown pub(crate) struct OperatorCore<'a> {
279df49d7eSJed Brown     ceed: &'a crate::Ceed,
289df49d7eSJed Brown     ptr: bind_ceed::CeedOperator,
299df49d7eSJed Brown }
309df49d7eSJed Brown 
319df49d7eSJed Brown pub struct Operator<'a> {
329df49d7eSJed Brown     op_core: OperatorCore<'a>,
339df49d7eSJed Brown }
349df49d7eSJed Brown 
359df49d7eSJed Brown pub struct CompositeOperator<'a> {
369df49d7eSJed Brown     op_core: OperatorCore<'a>,
379df49d7eSJed Brown }
389df49d7eSJed Brown 
399df49d7eSJed Brown // -----------------------------------------------------------------------------
409df49d7eSJed Brown // Destructor
419df49d7eSJed Brown // -----------------------------------------------------------------------------
429df49d7eSJed Brown impl<'a> Drop for OperatorCore<'a> {
439df49d7eSJed Brown     fn drop(&mut self) {
449df49d7eSJed Brown         unsafe {
459df49d7eSJed Brown             bind_ceed::CeedOperatorDestroy(&mut self.ptr);
469df49d7eSJed Brown         }
479df49d7eSJed Brown     }
489df49d7eSJed Brown }
499df49d7eSJed Brown 
509df49d7eSJed Brown // -----------------------------------------------------------------------------
519df49d7eSJed Brown // Display
529df49d7eSJed Brown // -----------------------------------------------------------------------------
539df49d7eSJed Brown impl<'a> fmt::Display for OperatorCore<'a> {
549df49d7eSJed Brown     fn fmt(&self, f: &mut fmt::Formatter) -> fmt::Result {
559df49d7eSJed Brown         let mut ptr = std::ptr::null_mut();
569df49d7eSJed Brown         let mut sizeloc = crate::MAX_BUFFER_LENGTH;
579df49d7eSJed Brown         let cstring = unsafe {
589df49d7eSJed Brown             let file = bind_ceed::open_memstream(&mut ptr, &mut sizeloc);
599df49d7eSJed Brown             bind_ceed::CeedOperatorView(self.ptr, file);
609df49d7eSJed Brown             bind_ceed::fclose(file);
619df49d7eSJed Brown             CString::from_raw(ptr)
629df49d7eSJed Brown         };
639df49d7eSJed Brown         cstring.to_string_lossy().fmt(f)
649df49d7eSJed Brown     }
659df49d7eSJed Brown }
669df49d7eSJed Brown 
679df49d7eSJed Brown /// View an Operator
689df49d7eSJed Brown ///
699df49d7eSJed Brown /// ```
709df49d7eSJed Brown /// # use libceed::prelude::*;
719df49d7eSJed Brown /// # let ceed = libceed::Ceed::default_init();
729df49d7eSJed Brown /// let qf = ceed.q_function_interior_by_name("Mass1DBuild").unwrap();
739df49d7eSJed Brown ///
749df49d7eSJed Brown /// // Operator field arguments
759df49d7eSJed Brown /// let ne = 3;
769df49d7eSJed Brown /// let q = 4 as usize;
779df49d7eSJed Brown /// let mut ind: Vec<i32> = vec![0; 2 * ne];
789df49d7eSJed Brown /// for i in 0..ne {
799df49d7eSJed Brown ///     ind[2 * i + 0] = i as i32;
809df49d7eSJed Brown ///     ind[2 * i + 1] = (i + 1) as i32;
819df49d7eSJed Brown /// }
829df49d7eSJed Brown /// let r = ceed
839df49d7eSJed Brown ///     .elem_restriction(ne, 2, 1, 1, ne + 1, MemType::Host, &ind)
849df49d7eSJed Brown ///     .unwrap();
859df49d7eSJed Brown /// let strides: [i32; 3] = [1, q as i32, q as i32];
869df49d7eSJed Brown /// let rq = ceed
879df49d7eSJed Brown ///     .strided_elem_restriction(ne, 2, 1, q * ne, strides)
889df49d7eSJed Brown ///     .unwrap();
899df49d7eSJed Brown ///
909df49d7eSJed Brown /// let b = ceed
919df49d7eSJed Brown ///     .basis_tensor_H1_Lagrange(1, 1, 2, q, QuadMode::Gauss)
929df49d7eSJed Brown ///     .unwrap();
939df49d7eSJed Brown ///
949df49d7eSJed Brown /// // Operator fields
959df49d7eSJed Brown /// let op = ceed
969df49d7eSJed Brown ///     .operator(&qf, QFunctionOpt::None, QFunctionOpt::None)
979df49d7eSJed Brown ///     .unwrap()
989df49d7eSJed Brown ///     .field("dx", &r, &b, VectorOpt::Active)
999df49d7eSJed Brown ///     .unwrap()
1009df49d7eSJed Brown ///     .field("weights", ElemRestrictionOpt::None, &b, VectorOpt::None)
1019df49d7eSJed Brown ///     .unwrap()
1029df49d7eSJed Brown ///     .field("qdata", &rq, BasisOpt::Collocated, VectorOpt::Active)
1039df49d7eSJed Brown ///     .unwrap();
1049df49d7eSJed Brown ///
1059df49d7eSJed Brown /// println!("{}", op);
1069df49d7eSJed Brown /// ```
1079df49d7eSJed Brown impl<'a> fmt::Display for Operator<'a> {
1089df49d7eSJed Brown     fn fmt(&self, f: &mut fmt::Formatter) -> fmt::Result {
1099df49d7eSJed Brown         self.op_core.fmt(f)
1109df49d7eSJed Brown     }
1119df49d7eSJed Brown }
1129df49d7eSJed Brown 
1139df49d7eSJed Brown /// View a composite Operator
1149df49d7eSJed Brown ///
1159df49d7eSJed Brown /// ```
1169df49d7eSJed Brown /// # use libceed::prelude::*;
1179df49d7eSJed Brown /// # let ceed = libceed::Ceed::default_init();
1189df49d7eSJed Brown ///
1199df49d7eSJed Brown /// // Sub operator field arguments
1209df49d7eSJed Brown /// let ne = 3;
1219df49d7eSJed Brown /// let q = 4 as usize;
1229df49d7eSJed Brown /// let mut ind: Vec<i32> = vec![0; 2 * ne];
1239df49d7eSJed Brown /// for i in 0..ne {
1249df49d7eSJed Brown ///     ind[2 * i + 0] = i as i32;
1259df49d7eSJed Brown ///     ind[2 * i + 1] = (i + 1) as i32;
1269df49d7eSJed Brown /// }
1279df49d7eSJed Brown /// let r = ceed
1289df49d7eSJed Brown ///     .elem_restriction(ne, 2, 1, 1, ne + 1, MemType::Host, &ind)
1299df49d7eSJed Brown ///     .unwrap();
1309df49d7eSJed Brown /// let strides: [i32; 3] = [1, q as i32, q as i32];
1319df49d7eSJed Brown /// let rq = ceed
1329df49d7eSJed Brown ///     .strided_elem_restriction(ne, 2, 1, q * ne, strides)
1339df49d7eSJed Brown ///     .unwrap();
1349df49d7eSJed Brown ///
1359df49d7eSJed Brown /// let b = ceed
1369df49d7eSJed Brown ///     .basis_tensor_H1_Lagrange(1, 1, 2, q, QuadMode::Gauss)
1379df49d7eSJed Brown ///     .unwrap();
1389df49d7eSJed Brown ///
1399df49d7eSJed Brown /// let qdata_mass = ceed.vector(q * ne).unwrap();
1409df49d7eSJed Brown /// let qdata_diff = ceed.vector(q * ne).unwrap();
1419df49d7eSJed Brown ///
1429df49d7eSJed Brown /// let qf_mass = ceed.q_function_interior_by_name("MassApply").unwrap();
1439df49d7eSJed Brown /// let op_mass = ceed
1449df49d7eSJed Brown ///     .operator(&qf_mass, QFunctionOpt::None, QFunctionOpt::None)
1459df49d7eSJed Brown ///     .unwrap()
1469df49d7eSJed Brown ///     .field("u", &r, &b, VectorOpt::Active)
1479df49d7eSJed Brown ///     .unwrap()
1489df49d7eSJed Brown ///     .field("qdata", &rq, BasisOpt::Collocated, &qdata_mass)
1499df49d7eSJed Brown ///     .unwrap()
1509df49d7eSJed Brown ///     .field("v", &r, &b, VectorOpt::Active)
1519df49d7eSJed Brown ///     .unwrap();
1529df49d7eSJed Brown ///
1539df49d7eSJed Brown /// let qf_diff = ceed.q_function_interior_by_name("Poisson1DApply").unwrap();
1549df49d7eSJed Brown /// let op_diff = ceed
1559df49d7eSJed Brown ///     .operator(&qf_diff, QFunctionOpt::None, QFunctionOpt::None)
1569df49d7eSJed Brown ///     .unwrap()
1579df49d7eSJed Brown ///     .field("du", &r, &b, VectorOpt::Active)
1589df49d7eSJed Brown ///     .unwrap()
1599df49d7eSJed Brown ///     .field("qdata", &rq, BasisOpt::Collocated, &qdata_diff)
1609df49d7eSJed Brown ///     .unwrap()
1619df49d7eSJed Brown ///     .field("dv", &r, &b, VectorOpt::Active)
1629df49d7eSJed Brown ///     .unwrap();
1639df49d7eSJed Brown ///
1649df49d7eSJed Brown /// let op = ceed
1659df49d7eSJed Brown ///     .composite_operator()
1669df49d7eSJed Brown ///     .unwrap()
1679df49d7eSJed Brown ///     .sub_operator(&op_mass)
1689df49d7eSJed Brown ///     .unwrap()
1699df49d7eSJed Brown ///     .sub_operator(&op_diff)
1709df49d7eSJed Brown ///     .unwrap();
1719df49d7eSJed Brown ///
1729df49d7eSJed Brown /// println!("{}", op);
1739df49d7eSJed Brown /// ```
1749df49d7eSJed Brown impl<'a> fmt::Display for CompositeOperator<'a> {
1759df49d7eSJed Brown     fn fmt(&self, f: &mut fmt::Formatter) -> fmt::Result {
1769df49d7eSJed Brown         self.op_core.fmt(f)
1779df49d7eSJed Brown     }
1789df49d7eSJed Brown }
1799df49d7eSJed Brown 
1809df49d7eSJed Brown // -----------------------------------------------------------------------------
1819df49d7eSJed Brown // Core functionality
1829df49d7eSJed Brown // -----------------------------------------------------------------------------
1839df49d7eSJed Brown impl<'a> OperatorCore<'a> {
1849df49d7eSJed Brown     // Common implementations
1859df49d7eSJed Brown     pub fn apply(&self, input: &Vector, output: &mut Vector) -> crate::Result<i32> {
1869df49d7eSJed Brown         let ierr = unsafe {
1879df49d7eSJed Brown             bind_ceed::CeedOperatorApply(
1889df49d7eSJed Brown                 self.ptr,
1899df49d7eSJed Brown                 input.ptr,
1909df49d7eSJed Brown                 output.ptr,
1919df49d7eSJed Brown                 bind_ceed::CEED_REQUEST_IMMEDIATE,
1929df49d7eSJed Brown             )
1939df49d7eSJed Brown         };
1949df49d7eSJed Brown         self.ceed.check_error(ierr)
1959df49d7eSJed Brown     }
1969df49d7eSJed Brown 
1979df49d7eSJed Brown     pub fn apply_add(&self, input: &Vector, output: &mut Vector) -> crate::Result<i32> {
1989df49d7eSJed Brown         let ierr = unsafe {
1999df49d7eSJed Brown             bind_ceed::CeedOperatorApplyAdd(
2009df49d7eSJed Brown                 self.ptr,
2019df49d7eSJed Brown                 input.ptr,
2029df49d7eSJed Brown                 output.ptr,
2039df49d7eSJed Brown                 bind_ceed::CEED_REQUEST_IMMEDIATE,
2049df49d7eSJed Brown             )
2059df49d7eSJed Brown         };
2069df49d7eSJed Brown         self.ceed.check_error(ierr)
2079df49d7eSJed Brown     }
2089df49d7eSJed Brown 
2099df49d7eSJed Brown     pub fn linear_assemble_diagonal(&self, assembled: &mut Vector) -> crate::Result<i32> {
2109df49d7eSJed Brown         let ierr = unsafe {
2119df49d7eSJed Brown             bind_ceed::CeedOperatorLinearAssembleDiagonal(
2129df49d7eSJed Brown                 self.ptr,
2139df49d7eSJed Brown                 assembled.ptr,
2149df49d7eSJed Brown                 bind_ceed::CEED_REQUEST_IMMEDIATE,
2159df49d7eSJed Brown             )
2169df49d7eSJed Brown         };
2179df49d7eSJed Brown         self.ceed.check_error(ierr)
2189df49d7eSJed Brown     }
2199df49d7eSJed Brown 
2209df49d7eSJed Brown     pub fn linear_assemble_add_diagonal(&self, assembled: &mut Vector) -> crate::Result<i32> {
2219df49d7eSJed Brown         let ierr = unsafe {
2229df49d7eSJed Brown             bind_ceed::CeedOperatorLinearAssembleAddDiagonal(
2239df49d7eSJed Brown                 self.ptr,
2249df49d7eSJed Brown                 assembled.ptr,
2259df49d7eSJed Brown                 bind_ceed::CEED_REQUEST_IMMEDIATE,
2269df49d7eSJed Brown             )
2279df49d7eSJed Brown         };
2289df49d7eSJed Brown         self.ceed.check_error(ierr)
2299df49d7eSJed Brown     }
2309df49d7eSJed Brown 
2319df49d7eSJed Brown     pub fn linear_assemble_point_block_diagonal(
2329df49d7eSJed Brown         &self,
2339df49d7eSJed Brown         assembled: &mut Vector,
2349df49d7eSJed Brown     ) -> crate::Result<i32> {
2359df49d7eSJed Brown         let ierr = unsafe {
2369df49d7eSJed Brown             bind_ceed::CeedOperatorLinearAssemblePointBlockDiagonal(
2379df49d7eSJed Brown                 self.ptr,
2389df49d7eSJed Brown                 assembled.ptr,
2399df49d7eSJed Brown                 bind_ceed::CEED_REQUEST_IMMEDIATE,
2409df49d7eSJed Brown             )
2419df49d7eSJed Brown         };
2429df49d7eSJed Brown         self.ceed.check_error(ierr)
2439df49d7eSJed Brown     }
2449df49d7eSJed Brown 
2459df49d7eSJed Brown     pub fn linear_assemble_add_point_block_diagonal(
2469df49d7eSJed Brown         &self,
2479df49d7eSJed Brown         assembled: &mut Vector,
2489df49d7eSJed Brown     ) -> crate::Result<i32> {
2499df49d7eSJed Brown         let ierr = unsafe {
2509df49d7eSJed Brown             bind_ceed::CeedOperatorLinearAssembleAddPointBlockDiagonal(
2519df49d7eSJed Brown                 self.ptr,
2529df49d7eSJed Brown                 assembled.ptr,
2539df49d7eSJed Brown                 bind_ceed::CEED_REQUEST_IMMEDIATE,
2549df49d7eSJed Brown             )
2559df49d7eSJed Brown         };
2569df49d7eSJed Brown         self.ceed.check_error(ierr)
2579df49d7eSJed Brown     }
2589df49d7eSJed Brown }
2599df49d7eSJed Brown 
2609df49d7eSJed Brown // -----------------------------------------------------------------------------
2619df49d7eSJed Brown // Operator
2629df49d7eSJed Brown // -----------------------------------------------------------------------------
2639df49d7eSJed Brown impl<'a> Operator<'a> {
2649df49d7eSJed Brown     // Constructor
2659df49d7eSJed Brown     pub fn create<'b>(
2669df49d7eSJed Brown         ceed: &'a crate::Ceed,
2679df49d7eSJed Brown         qf: impl Into<QFunctionOpt<'b>>,
2689df49d7eSJed Brown         dqf: impl Into<QFunctionOpt<'b>>,
2699df49d7eSJed Brown         dqfT: impl Into<QFunctionOpt<'b>>,
2709df49d7eSJed Brown     ) -> crate::Result<Self> {
2719df49d7eSJed Brown         let mut ptr = std::ptr::null_mut();
2729df49d7eSJed Brown         let ierr = unsafe {
2739df49d7eSJed Brown             bind_ceed::CeedOperatorCreate(
2749df49d7eSJed Brown                 ceed.ptr,
2759df49d7eSJed Brown                 qf.into().to_raw(),
2769df49d7eSJed Brown                 dqf.into().to_raw(),
2779df49d7eSJed Brown                 dqfT.into().to_raw(),
2789df49d7eSJed Brown                 &mut ptr,
2799df49d7eSJed Brown             )
2809df49d7eSJed Brown         };
2819df49d7eSJed Brown         ceed.check_error(ierr)?;
2829df49d7eSJed Brown         Ok(Self {
2839df49d7eSJed Brown             op_core: OperatorCore { ceed, ptr },
2849df49d7eSJed Brown         })
2859df49d7eSJed Brown     }
2869df49d7eSJed Brown 
2879df49d7eSJed Brown     fn from_raw(ceed: &'a crate::Ceed, ptr: bind_ceed::CeedOperator) -> crate::Result<Self> {
2889df49d7eSJed Brown         Ok(Self {
2899df49d7eSJed Brown             op_core: OperatorCore { ceed, ptr },
2909df49d7eSJed Brown         })
2919df49d7eSJed Brown     }
2929df49d7eSJed Brown 
2939df49d7eSJed Brown     /// Apply Operator to a vector
2949df49d7eSJed Brown     ///
2959df49d7eSJed Brown     /// * `input`  - Input Vector
2969df49d7eSJed Brown     /// * `output` - Output Vector
2979df49d7eSJed Brown     ///
2989df49d7eSJed Brown     /// ```
2999df49d7eSJed Brown     /// # use libceed::prelude::*;
3009df49d7eSJed Brown     /// # let ceed = libceed::Ceed::default_init();
3019df49d7eSJed Brown     /// let ne = 4;
3029df49d7eSJed Brown     /// let p = 3;
3039df49d7eSJed Brown     /// let q = 4;
3049df49d7eSJed Brown     /// let ndofs = p * ne - ne + 1;
3059df49d7eSJed Brown     ///
3069df49d7eSJed Brown     /// // Vectors
3079df49d7eSJed Brown     /// let x = ceed.vector_from_slice(&[-1., -0.5, 0.0, 0.5, 1.0]).unwrap();
3089df49d7eSJed Brown     /// let mut qdata = ceed.vector(ne * q).unwrap();
3099df49d7eSJed Brown     /// qdata.set_value(0.0);
3109df49d7eSJed Brown     /// let u = ceed.vector_from_slice(&vec![1.0; ndofs]).unwrap();
3119df49d7eSJed Brown     /// let mut v = ceed.vector(ndofs).unwrap();
3129df49d7eSJed Brown     /// v.set_value(0.0);
3139df49d7eSJed Brown     ///
3149df49d7eSJed Brown     /// // Restrictions
3159df49d7eSJed Brown     /// let mut indx: Vec<i32> = vec![0; 2 * ne];
3169df49d7eSJed Brown     /// for i in 0..ne {
3179df49d7eSJed Brown     ///     indx[2 * i + 0] = i as i32;
3189df49d7eSJed Brown     ///     indx[2 * i + 1] = (i + 1) as i32;
3199df49d7eSJed Brown     /// }
3209df49d7eSJed Brown     /// let rx = ceed
3219df49d7eSJed Brown     ///     .elem_restriction(ne, 2, 1, 1, ne + 1, MemType::Host, &indx)
3229df49d7eSJed Brown     ///     .unwrap();
3239df49d7eSJed Brown     /// let mut indu: Vec<i32> = vec![0; p * ne];
3249df49d7eSJed Brown     /// for i in 0..ne {
3259df49d7eSJed Brown     ///     indu[p * i + 0] = i as i32;
3269df49d7eSJed Brown     ///     indu[p * i + 1] = (i + 1) as i32;
3279df49d7eSJed Brown     ///     indu[p * i + 2] = (i + 2) as i32;
3289df49d7eSJed Brown     /// }
3299df49d7eSJed Brown     /// let ru = ceed
3309df49d7eSJed Brown     ///     .elem_restriction(ne, 3, 1, 1, ndofs, MemType::Host, &indu)
3319df49d7eSJed Brown     ///     .unwrap();
3329df49d7eSJed Brown     /// let strides: [i32; 3] = [1, q as i32, q as i32];
3339df49d7eSJed Brown     /// let rq = ceed
3349df49d7eSJed Brown     ///     .strided_elem_restriction(ne, q, 1, q * ne, strides)
3359df49d7eSJed Brown     ///     .unwrap();
3369df49d7eSJed Brown     ///
3379df49d7eSJed Brown     /// // Bases
3389df49d7eSJed Brown     /// let bx = ceed
3399df49d7eSJed Brown     ///     .basis_tensor_H1_Lagrange(1, 1, 2, q, QuadMode::Gauss)
3409df49d7eSJed Brown     ///     .unwrap();
3419df49d7eSJed Brown     /// let bu = ceed
3429df49d7eSJed Brown     ///     .basis_tensor_H1_Lagrange(1, 1, p, q, QuadMode::Gauss)
3439df49d7eSJed Brown     ///     .unwrap();
3449df49d7eSJed Brown     ///
3459df49d7eSJed Brown     /// // Build quadrature data
3469df49d7eSJed Brown     /// let qf_build = ceed.q_function_interior_by_name("Mass1DBuild").unwrap();
3479df49d7eSJed Brown     /// ceed.operator(&qf_build, QFunctionOpt::None, QFunctionOpt::None)
3489df49d7eSJed Brown     ///     .unwrap()
3499df49d7eSJed Brown     ///     .field("dx", &rx, &bx, VectorOpt::Active)
3509df49d7eSJed Brown     ///     .unwrap()
3519df49d7eSJed Brown     ///     .field("weights", ElemRestrictionOpt::None, &bx, VectorOpt::None)
3529df49d7eSJed Brown     ///     .unwrap()
3539df49d7eSJed Brown     ///     .field("qdata", &rq, BasisOpt::Collocated, VectorOpt::Active)
3549df49d7eSJed Brown     ///     .unwrap()
3559df49d7eSJed Brown     ///     .apply(&x, &mut qdata)
3569df49d7eSJed Brown     ///     .unwrap();
3579df49d7eSJed Brown     ///
3589df49d7eSJed Brown     /// // Mass operator
3599df49d7eSJed Brown     /// let qf_mass = ceed.q_function_interior_by_name("MassApply").unwrap();
3609df49d7eSJed Brown     /// let op_mass = ceed
3619df49d7eSJed Brown     ///     .operator(&qf_mass, QFunctionOpt::None, QFunctionOpt::None)
3629df49d7eSJed Brown     ///     .unwrap()
3639df49d7eSJed Brown     ///     .field("u", &ru, &bu, VectorOpt::Active)
3649df49d7eSJed Brown     ///     .unwrap()
3659df49d7eSJed Brown     ///     .field("qdata", &rq, BasisOpt::Collocated, &qdata)
3669df49d7eSJed Brown     ///     .unwrap()
3679df49d7eSJed Brown     ///     .field("v", &ru, &bu, VectorOpt::Active)
3689df49d7eSJed Brown     ///     .unwrap();
3699df49d7eSJed Brown     ///
3709df49d7eSJed Brown     /// v.set_value(0.0);
3719df49d7eSJed Brown     /// op_mass.apply(&u, &mut v).unwrap();
3729df49d7eSJed Brown     ///
3739df49d7eSJed Brown     /// // Check
374*80a9ef05SNatalie Beams     /// let sum: Scalar = v.view().iter().sum();
3759df49d7eSJed Brown     /// assert!(
376*80a9ef05SNatalie Beams     ///     (sum - 2.0).abs() < 10.0 * libceed::EPSILON,
3779df49d7eSJed Brown     ///     "Incorrect interval length computed"
3789df49d7eSJed Brown     /// );
3799df49d7eSJed Brown     /// ```
3809df49d7eSJed Brown     pub fn apply(&self, input: &Vector, output: &mut Vector) -> crate::Result<i32> {
3819df49d7eSJed Brown         self.op_core.apply(input, output)
3829df49d7eSJed Brown     }
3839df49d7eSJed Brown 
3849df49d7eSJed Brown     /// Apply Operator to a vector and add result to output vector
3859df49d7eSJed Brown     ///
3869df49d7eSJed Brown     /// * `input`  - Input Vector
3879df49d7eSJed Brown     /// * `output` - Output Vector
3889df49d7eSJed Brown     ///
3899df49d7eSJed Brown     /// ```
3909df49d7eSJed Brown     /// # use libceed::prelude::*;
3919df49d7eSJed Brown     /// # let ceed = libceed::Ceed::default_init();
3929df49d7eSJed Brown     /// let ne = 4;
3939df49d7eSJed Brown     /// let p = 3;
3949df49d7eSJed Brown     /// let q = 4;
3959df49d7eSJed Brown     /// let ndofs = p * ne - ne + 1;
3969df49d7eSJed Brown     ///
3979df49d7eSJed Brown     /// // Vectors
3989df49d7eSJed Brown     /// let x = ceed.vector_from_slice(&[-1., -0.5, 0.0, 0.5, 1.0]).unwrap();
3999df49d7eSJed Brown     /// let mut qdata = ceed.vector(ne * q).unwrap();
4009df49d7eSJed Brown     /// qdata.set_value(0.0);
4019df49d7eSJed Brown     /// let u = ceed.vector_from_slice(&vec![1.0; ndofs]).unwrap();
4029df49d7eSJed Brown     /// let mut v = ceed.vector(ndofs).unwrap();
4039df49d7eSJed Brown     ///
4049df49d7eSJed Brown     /// // Restrictions
4059df49d7eSJed Brown     /// let mut indx: Vec<i32> = vec![0; 2 * ne];
4069df49d7eSJed Brown     /// for i in 0..ne {
4079df49d7eSJed Brown     ///     indx[2 * i + 0] = i as i32;
4089df49d7eSJed Brown     ///     indx[2 * i + 1] = (i + 1) as i32;
4099df49d7eSJed Brown     /// }
4109df49d7eSJed Brown     /// let rx = ceed
4119df49d7eSJed Brown     ///     .elem_restriction(ne, 2, 1, 1, ne + 1, MemType::Host, &indx)
4129df49d7eSJed Brown     ///     .unwrap();
4139df49d7eSJed Brown     /// let mut indu: Vec<i32> = vec![0; p * ne];
4149df49d7eSJed Brown     /// for i in 0..ne {
4159df49d7eSJed Brown     ///     indu[p * i + 0] = i as i32;
4169df49d7eSJed Brown     ///     indu[p * i + 1] = (i + 1) as i32;
4179df49d7eSJed Brown     ///     indu[p * i + 2] = (i + 2) as i32;
4189df49d7eSJed Brown     /// }
4199df49d7eSJed Brown     /// let ru = ceed
4209df49d7eSJed Brown     ///     .elem_restriction(ne, 3, 1, 1, ndofs, MemType::Host, &indu)
4219df49d7eSJed Brown     ///     .unwrap();
4229df49d7eSJed Brown     /// let strides: [i32; 3] = [1, q as i32, q as i32];
4239df49d7eSJed Brown     /// let rq = ceed
4249df49d7eSJed Brown     ///     .strided_elem_restriction(ne, q, 1, q * ne, strides)
4259df49d7eSJed Brown     ///     .unwrap();
4269df49d7eSJed Brown     ///
4279df49d7eSJed Brown     /// // Bases
4289df49d7eSJed Brown     /// let bx = ceed
4299df49d7eSJed Brown     ///     .basis_tensor_H1_Lagrange(1, 1, 2, q, QuadMode::Gauss)
4309df49d7eSJed Brown     ///     .unwrap();
4319df49d7eSJed Brown     /// let bu = ceed
4329df49d7eSJed Brown     ///     .basis_tensor_H1_Lagrange(1, 1, p, q, QuadMode::Gauss)
4339df49d7eSJed Brown     ///     .unwrap();
4349df49d7eSJed Brown     ///
4359df49d7eSJed Brown     /// // Build quadrature data
4369df49d7eSJed Brown     /// let qf_build = ceed.q_function_interior_by_name("Mass1DBuild").unwrap();
4379df49d7eSJed Brown     /// ceed.operator(&qf_build, QFunctionOpt::None, QFunctionOpt::None)
4389df49d7eSJed Brown     ///     .unwrap()
4399df49d7eSJed Brown     ///     .field("dx", &rx, &bx, VectorOpt::Active)
4409df49d7eSJed Brown     ///     .unwrap()
4419df49d7eSJed Brown     ///     .field("weights", ElemRestrictionOpt::None, &bx, VectorOpt::None)
4429df49d7eSJed Brown     ///     .unwrap()
4439df49d7eSJed Brown     ///     .field("qdata", &rq, BasisOpt::Collocated, VectorOpt::Active)
4449df49d7eSJed Brown     ///     .unwrap()
4459df49d7eSJed Brown     ///     .apply(&x, &mut qdata)
4469df49d7eSJed Brown     ///     .unwrap();
4479df49d7eSJed Brown     ///
4489df49d7eSJed Brown     /// // Mass operator
4499df49d7eSJed Brown     /// let qf_mass = ceed.q_function_interior_by_name("MassApply").unwrap();
4509df49d7eSJed Brown     /// let op_mass = ceed
4519df49d7eSJed Brown     ///     .operator(&qf_mass, QFunctionOpt::None, QFunctionOpt::None)
4529df49d7eSJed Brown     ///     .unwrap()
4539df49d7eSJed Brown     ///     .field("u", &ru, &bu, VectorOpt::Active)
4549df49d7eSJed Brown     ///     .unwrap()
4559df49d7eSJed Brown     ///     .field("qdata", &rq, BasisOpt::Collocated, &qdata)
4569df49d7eSJed Brown     ///     .unwrap()
4579df49d7eSJed Brown     ///     .field("v", &ru, &bu, VectorOpt::Active)
4589df49d7eSJed Brown     ///     .unwrap();
4599df49d7eSJed Brown     ///
4609df49d7eSJed Brown     /// v.set_value(1.0);
4619df49d7eSJed Brown     /// op_mass.apply_add(&u, &mut v).unwrap();
4629df49d7eSJed Brown     ///
4639df49d7eSJed Brown     /// // Check
464*80a9ef05SNatalie Beams     /// let sum: Scalar = v.view().iter().sum();
4659df49d7eSJed Brown     /// assert!(
466*80a9ef05SNatalie Beams     ///     (sum - (2.0 + ndofs as Scalar)).abs() < 10.0 * libceed::EPSILON,
4679df49d7eSJed Brown     ///     "Incorrect interval length computed and added"
4689df49d7eSJed Brown     /// );
4699df49d7eSJed Brown     /// ```
4709df49d7eSJed Brown     pub fn apply_add(&self, input: &Vector, output: &mut Vector) -> crate::Result<i32> {
4719df49d7eSJed Brown         self.op_core.apply_add(input, output)
4729df49d7eSJed Brown     }
4739df49d7eSJed Brown 
4749df49d7eSJed Brown     /// Provide a field to a Operator for use by its QFunction
4759df49d7eSJed Brown     ///
4769df49d7eSJed Brown     /// * `fieldname` - Name of the field (to be matched with the name used by
4779df49d7eSJed Brown     ///                   the QFunction)
4789df49d7eSJed Brown     /// * `r`         - ElemRestriction
4799df49d7eSJed Brown     /// * `b`         - Basis in which the field resides or `BasisCollocated` if
4809df49d7eSJed Brown     ///                   collocated with quadrature points
4819df49d7eSJed Brown     /// * `v`         - Vector to be used by Operator or `VectorActive` if field
4829df49d7eSJed Brown     ///                   is active or `VectorNone` if using `Weight` with the
4839df49d7eSJed Brown     ///                   QFunction
4849df49d7eSJed Brown     ///
4859df49d7eSJed Brown     ///
4869df49d7eSJed Brown     /// ```
4879df49d7eSJed Brown     /// # use libceed::prelude::*;
4889df49d7eSJed Brown     /// # let ceed = libceed::Ceed::default_init();
4899df49d7eSJed Brown     /// let qf = ceed.q_function_interior_by_name("Mass1DBuild").unwrap();
4909df49d7eSJed Brown     /// let mut op = ceed
4919df49d7eSJed Brown     ///     .operator(&qf, QFunctionOpt::None, QFunctionOpt::None)
4929df49d7eSJed Brown     ///     .unwrap();
4939df49d7eSJed Brown     ///
4949df49d7eSJed Brown     /// // Operator field arguments
4959df49d7eSJed Brown     /// let ne = 3;
4969df49d7eSJed Brown     /// let q = 4;
4979df49d7eSJed Brown     /// let mut ind: Vec<i32> = vec![0; 2 * ne];
4989df49d7eSJed Brown     /// for i in 0..ne {
4999df49d7eSJed Brown     ///     ind[2 * i + 0] = i as i32;
5009df49d7eSJed Brown     ///     ind[2 * i + 1] = (i + 1) as i32;
5019df49d7eSJed Brown     /// }
5029df49d7eSJed Brown     /// let r = ceed
5039df49d7eSJed Brown     ///     .elem_restriction(ne, 2, 1, 1, ne + 1, MemType::Host, &ind)
5049df49d7eSJed Brown     ///     .unwrap();
5059df49d7eSJed Brown     ///
5069df49d7eSJed Brown     /// let b = ceed
5079df49d7eSJed Brown     ///     .basis_tensor_H1_Lagrange(1, 1, 2, q, QuadMode::Gauss)
5089df49d7eSJed Brown     ///     .unwrap();
5099df49d7eSJed Brown     ///
5109df49d7eSJed Brown     /// // Operator field
5119df49d7eSJed Brown     /// op = op.field("dx", &r, &b, VectorOpt::Active).unwrap();
5129df49d7eSJed Brown     /// ```
5139df49d7eSJed Brown     #[allow(unused_mut)]
5149df49d7eSJed Brown     pub fn field<'b>(
5159df49d7eSJed Brown         mut self,
5169df49d7eSJed Brown         fieldname: &str,
5179df49d7eSJed Brown         r: impl Into<ElemRestrictionOpt<'b>>,
5189df49d7eSJed Brown         b: impl Into<BasisOpt<'b>>,
5199df49d7eSJed Brown         v: impl Into<VectorOpt<'b>>,
5209df49d7eSJed Brown     ) -> crate::Result<Self> {
5219df49d7eSJed Brown         let fieldname = CString::new(fieldname).expect("CString::new failed");
5229df49d7eSJed Brown         let fieldname = fieldname.as_ptr() as *const i8;
5239df49d7eSJed Brown         let ierr = unsafe {
5249df49d7eSJed Brown             bind_ceed::CeedOperatorSetField(
5259df49d7eSJed Brown                 self.op_core.ptr,
5269df49d7eSJed Brown                 fieldname,
5279df49d7eSJed Brown                 r.into().to_raw(),
5289df49d7eSJed Brown                 b.into().to_raw(),
5299df49d7eSJed Brown                 v.into().to_raw(),
5309df49d7eSJed Brown             )
5319df49d7eSJed Brown         };
5329df49d7eSJed Brown         self.op_core.ceed.check_error(ierr)?;
5339df49d7eSJed Brown         Ok(self)
5349df49d7eSJed Brown     }
5359df49d7eSJed Brown 
5369df49d7eSJed Brown     /// Assemble the diagonal of a square linear Operator
5379df49d7eSJed Brown     ///
5389df49d7eSJed Brown     /// This overwrites a Vector with the diagonal of a linear Operator.
5399df49d7eSJed Brown     ///
5409df49d7eSJed Brown     /// Note: Currently only non-composite Operators with a single field and
5419df49d7eSJed Brown     ///      composite Operators with single field sub-operators are supported.
5429df49d7eSJed Brown     ///
5439df49d7eSJed Brown     /// * `op`        - Operator to assemble QFunction
5449df49d7eSJed Brown     /// * `assembled` - Vector to store assembled Operator diagonal
5459df49d7eSJed Brown     ///
5469df49d7eSJed Brown     /// ```
5479df49d7eSJed Brown     /// # use libceed::prelude::*;
5489df49d7eSJed Brown     /// # let ceed = libceed::Ceed::default_init();
5499df49d7eSJed Brown     /// let ne = 4;
5509df49d7eSJed Brown     /// let p = 3;
5519df49d7eSJed Brown     /// let q = 4;
5529df49d7eSJed Brown     /// let ndofs = p * ne - ne + 1;
5539df49d7eSJed Brown     ///
5549df49d7eSJed Brown     /// // Vectors
5559df49d7eSJed Brown     /// let x = ceed.vector_from_slice(&[-1., -0.5, 0.0, 0.5, 1.0]).unwrap();
5569df49d7eSJed Brown     /// let mut qdata = ceed.vector(ne * q).unwrap();
5579df49d7eSJed Brown     /// qdata.set_value(0.0);
5589df49d7eSJed Brown     /// let mut u = ceed.vector(ndofs).unwrap();
5599df49d7eSJed Brown     /// u.set_value(1.0);
5609df49d7eSJed Brown     /// let mut v = ceed.vector(ndofs).unwrap();
5619df49d7eSJed Brown     /// v.set_value(0.0);
5629df49d7eSJed Brown     ///
5639df49d7eSJed Brown     /// // Restrictions
5649df49d7eSJed Brown     /// let mut indx: Vec<i32> = vec![0; 2 * ne];
5659df49d7eSJed Brown     /// for i in 0..ne {
5669df49d7eSJed Brown     ///     indx[2 * i + 0] = i as i32;
5679df49d7eSJed Brown     ///     indx[2 * i + 1] = (i + 1) as i32;
5689df49d7eSJed Brown     /// }
5699df49d7eSJed Brown     /// let rx = ceed
5709df49d7eSJed Brown     ///     .elem_restriction(ne, 2, 1, 1, ne + 1, MemType::Host, &indx)
5719df49d7eSJed Brown     ///     .unwrap();
5729df49d7eSJed Brown     /// let mut indu: Vec<i32> = vec![0; p * ne];
5739df49d7eSJed Brown     /// for i in 0..ne {
5749df49d7eSJed Brown     ///     indu[p * i + 0] = (2 * i) as i32;
5759df49d7eSJed Brown     ///     indu[p * i + 1] = (2 * i + 1) as i32;
5769df49d7eSJed Brown     ///     indu[p * i + 2] = (2 * i + 2) as i32;
5779df49d7eSJed Brown     /// }
5789df49d7eSJed Brown     /// let ru = ceed
5799df49d7eSJed Brown     ///     .elem_restriction(ne, p, 1, 1, ndofs, MemType::Host, &indu)
5809df49d7eSJed Brown     ///     .unwrap();
5819df49d7eSJed Brown     /// let strides: [i32; 3] = [1, q as i32, q as i32];
5829df49d7eSJed Brown     /// let rq = ceed
5839df49d7eSJed Brown     ///     .strided_elem_restriction(ne, q, 1, q * ne, strides)
5849df49d7eSJed Brown     ///     .unwrap();
5859df49d7eSJed Brown     ///
5869df49d7eSJed Brown     /// // Bases
5879df49d7eSJed Brown     /// let bx = ceed
5889df49d7eSJed Brown     ///     .basis_tensor_H1_Lagrange(1, 1, 2, q, QuadMode::Gauss)
5899df49d7eSJed Brown     ///     .unwrap();
5909df49d7eSJed Brown     /// let bu = ceed
5919df49d7eSJed Brown     ///     .basis_tensor_H1_Lagrange(1, 1, p, q, QuadMode::Gauss)
5929df49d7eSJed Brown     ///     .unwrap();
5939df49d7eSJed Brown     ///
5949df49d7eSJed Brown     /// // Build quadrature data
5959df49d7eSJed Brown     /// let qf_build = ceed.q_function_interior_by_name("Mass1DBuild").unwrap();
5969df49d7eSJed Brown     /// ceed.operator(&qf_build, QFunctionOpt::None, QFunctionOpt::None)
5979df49d7eSJed Brown     ///     .unwrap()
5989df49d7eSJed Brown     ///     .field("dx", &rx, &bx, VectorOpt::Active)
5999df49d7eSJed Brown     ///     .unwrap()
6009df49d7eSJed Brown     ///     .field("weights", ElemRestrictionOpt::None, &bx, VectorOpt::None)
6019df49d7eSJed Brown     ///     .unwrap()
6029df49d7eSJed Brown     ///     .field("qdata", &rq, BasisOpt::Collocated, VectorOpt::Active)
6039df49d7eSJed Brown     ///     .unwrap()
6049df49d7eSJed Brown     ///     .apply(&x, &mut qdata)
6059df49d7eSJed Brown     ///     .unwrap();
6069df49d7eSJed Brown     ///
6079df49d7eSJed Brown     /// // Mass operator
6089df49d7eSJed Brown     /// let qf_mass = ceed.q_function_interior_by_name("MassApply").unwrap();
6099df49d7eSJed Brown     /// let op_mass = ceed
6109df49d7eSJed Brown     ///     .operator(&qf_mass, QFunctionOpt::None, QFunctionOpt::None)
6119df49d7eSJed Brown     ///     .unwrap()
6129df49d7eSJed Brown     ///     .field("u", &ru, &bu, VectorOpt::Active)
6139df49d7eSJed Brown     ///     .unwrap()
6149df49d7eSJed Brown     ///     .field("qdata", &rq, BasisOpt::Collocated, &qdata)
6159df49d7eSJed Brown     ///     .unwrap()
6169df49d7eSJed Brown     ///     .field("v", &ru, &bu, VectorOpt::Active)
6179df49d7eSJed Brown     ///     .unwrap();
6189df49d7eSJed Brown     ///
6199df49d7eSJed Brown     /// // Diagonal
6209df49d7eSJed Brown     /// let mut diag = ceed.vector(ndofs).unwrap();
6219df49d7eSJed Brown     /// diag.set_value(0.0);
6229df49d7eSJed Brown     /// op_mass.linear_assemble_diagonal(&mut diag).unwrap();
6239df49d7eSJed Brown     ///
6249df49d7eSJed Brown     /// // Manual diagonal computation
6259df49d7eSJed Brown     /// let mut true_diag = ceed.vector(ndofs).unwrap();
6269df49d7eSJed Brown     /// for i in 0..ndofs {
6279df49d7eSJed Brown     ///     u.set_value(0.0);
6289df49d7eSJed Brown     ///     {
6299df49d7eSJed Brown     ///         let mut u_array = u.view_mut();
6309df49d7eSJed Brown     ///         u_array[i] = 1.;
6319df49d7eSJed Brown     ///     }
6329df49d7eSJed Brown     ///
6339df49d7eSJed Brown     ///     op_mass.apply(&u, &mut v).unwrap();
6349df49d7eSJed Brown     ///
6359df49d7eSJed Brown     ///     {
6369df49d7eSJed Brown     ///         let v_array = v.view_mut();
6379df49d7eSJed Brown     ///         let mut true_array = true_diag.view_mut();
6389df49d7eSJed Brown     ///         true_array[i] = v_array[i];
6399df49d7eSJed Brown     ///     }
6409df49d7eSJed Brown     /// }
6419df49d7eSJed Brown     ///
6429df49d7eSJed Brown     /// // Check
6439df49d7eSJed Brown     /// diag.view()
6449df49d7eSJed Brown     ///     .iter()
6459df49d7eSJed Brown     ///     .zip(true_diag.view().iter())
6469df49d7eSJed Brown     ///     .for_each(|(computed, actual)| {
6479df49d7eSJed Brown     ///         assert!(
648*80a9ef05SNatalie Beams     ///             (*computed - *actual).abs() < 10.0 * libceed::EPSILON,
6499df49d7eSJed Brown     ///             "Diagonal entry incorrect"
6509df49d7eSJed Brown     ///         );
6519df49d7eSJed Brown     ///     });
6529df49d7eSJed Brown     /// ```
6539df49d7eSJed Brown     pub fn linear_assemble_diagonal(&self, assembled: &mut Vector) -> crate::Result<i32> {
6549df49d7eSJed Brown         self.op_core.linear_assemble_diagonal(assembled)
6559df49d7eSJed Brown     }
6569df49d7eSJed Brown 
6579df49d7eSJed Brown     /// Assemble the diagonal of a square linear Operator
6589df49d7eSJed Brown     ///
6599df49d7eSJed Brown     /// This sums into a Vector with the diagonal of a linear Operator.
6609df49d7eSJed Brown     ///
6619df49d7eSJed Brown     /// Note: Currently only non-composite Operators with a single field and
6629df49d7eSJed Brown     ///      composite Operators with single field sub-operators are supported.
6639df49d7eSJed Brown     ///
6649df49d7eSJed Brown     /// * `op`        - Operator to assemble QFunction
6659df49d7eSJed Brown     /// * `assembled` - Vector to store assembled Operator diagonal
6669df49d7eSJed Brown     ///
6679df49d7eSJed Brown     ///
6689df49d7eSJed Brown     /// ```
6699df49d7eSJed Brown     /// # use libceed::prelude::*;
6709df49d7eSJed Brown     /// # let ceed = libceed::Ceed::default_init();
6719df49d7eSJed Brown     /// let ne = 4;
6729df49d7eSJed Brown     /// let p = 3;
6739df49d7eSJed Brown     /// let q = 4;
6749df49d7eSJed Brown     /// let ndofs = p * ne - ne + 1;
6759df49d7eSJed Brown     ///
6769df49d7eSJed Brown     /// // Vectors
6779df49d7eSJed Brown     /// let x = ceed.vector_from_slice(&[-1., -0.5, 0.0, 0.5, 1.0]).unwrap();
6789df49d7eSJed Brown     /// let mut qdata = ceed.vector(ne * q).unwrap();
6799df49d7eSJed Brown     /// qdata.set_value(0.0);
6809df49d7eSJed Brown     /// let mut u = ceed.vector(ndofs).unwrap();
6819df49d7eSJed Brown     /// u.set_value(1.0);
6829df49d7eSJed Brown     /// let mut v = ceed.vector(ndofs).unwrap();
6839df49d7eSJed Brown     /// v.set_value(0.0);
6849df49d7eSJed Brown     ///
6859df49d7eSJed Brown     /// // Restrictions
6869df49d7eSJed Brown     /// let mut indx: Vec<i32> = vec![0; 2 * ne];
6879df49d7eSJed Brown     /// for i in 0..ne {
6889df49d7eSJed Brown     ///     indx[2 * i + 0] = i as i32;
6899df49d7eSJed Brown     ///     indx[2 * i + 1] = (i + 1) as i32;
6909df49d7eSJed Brown     /// }
6919df49d7eSJed Brown     /// let rx = ceed
6929df49d7eSJed Brown     ///     .elem_restriction(ne, 2, 1, 1, ne + 1, MemType::Host, &indx)
6939df49d7eSJed Brown     ///     .unwrap();
6949df49d7eSJed Brown     /// let mut indu: Vec<i32> = vec![0; p * ne];
6959df49d7eSJed Brown     /// for i in 0..ne {
6969df49d7eSJed Brown     ///     indu[p * i + 0] = (2 * i) as i32;
6979df49d7eSJed Brown     ///     indu[p * i + 1] = (2 * i + 1) as i32;
6989df49d7eSJed Brown     ///     indu[p * i + 2] = (2 * i + 2) as i32;
6999df49d7eSJed Brown     /// }
7009df49d7eSJed Brown     /// let ru = ceed
7019df49d7eSJed Brown     ///     .elem_restriction(ne, p, 1, 1, ndofs, MemType::Host, &indu)
7029df49d7eSJed Brown     ///     .unwrap();
7039df49d7eSJed Brown     /// let strides: [i32; 3] = [1, q as i32, q as i32];
7049df49d7eSJed Brown     /// let rq = ceed
7059df49d7eSJed Brown     ///     .strided_elem_restriction(ne, q, 1, q * ne, strides)
7069df49d7eSJed Brown     ///     .unwrap();
7079df49d7eSJed Brown     ///
7089df49d7eSJed Brown     /// // Bases
7099df49d7eSJed Brown     /// let bx = ceed
7109df49d7eSJed Brown     ///     .basis_tensor_H1_Lagrange(1, 1, 2, q, QuadMode::Gauss)
7119df49d7eSJed Brown     ///     .unwrap();
7129df49d7eSJed Brown     /// let bu = ceed
7139df49d7eSJed Brown     ///     .basis_tensor_H1_Lagrange(1, 1, p, q, QuadMode::Gauss)
7149df49d7eSJed Brown     ///     .unwrap();
7159df49d7eSJed Brown     ///
7169df49d7eSJed Brown     /// // Build quadrature data
7179df49d7eSJed Brown     /// let qf_build = ceed.q_function_interior_by_name("Mass1DBuild").unwrap();
7189df49d7eSJed Brown     /// ceed.operator(&qf_build, QFunctionOpt::None, QFunctionOpt::None)
7199df49d7eSJed Brown     ///     .unwrap()
7209df49d7eSJed Brown     ///     .field("dx", &rx, &bx, VectorOpt::Active)
7219df49d7eSJed Brown     ///     .unwrap()
7229df49d7eSJed Brown     ///     .field("weights", ElemRestrictionOpt::None, &bx, VectorOpt::None)
7239df49d7eSJed Brown     ///     .unwrap()
7249df49d7eSJed Brown     ///     .field("qdata", &rq, BasisOpt::Collocated, VectorOpt::Active)
7259df49d7eSJed Brown     ///     .unwrap()
7269df49d7eSJed Brown     ///     .apply(&x, &mut qdata)
7279df49d7eSJed Brown     ///     .unwrap();
7289df49d7eSJed Brown     ///
7299df49d7eSJed Brown     /// // Mass operator
7309df49d7eSJed Brown     /// let qf_mass = ceed.q_function_interior_by_name("MassApply").unwrap();
7319df49d7eSJed Brown     /// let op_mass = ceed
7329df49d7eSJed Brown     ///     .operator(&qf_mass, QFunctionOpt::None, QFunctionOpt::None)
7339df49d7eSJed Brown     ///     .unwrap()
7349df49d7eSJed Brown     ///     .field("u", &ru, &bu, VectorOpt::Active)
7359df49d7eSJed Brown     ///     .unwrap()
7369df49d7eSJed Brown     ///     .field("qdata", &rq, BasisOpt::Collocated, &qdata)
7379df49d7eSJed Brown     ///     .unwrap()
7389df49d7eSJed Brown     ///     .field("v", &ru, &bu, VectorOpt::Active)
7399df49d7eSJed Brown     ///     .unwrap();
7409df49d7eSJed Brown     ///
7419df49d7eSJed Brown     /// // Diagonal
7429df49d7eSJed Brown     /// let mut diag = ceed.vector(ndofs).unwrap();
7439df49d7eSJed Brown     /// diag.set_value(1.0);
7449df49d7eSJed Brown     /// op_mass.linear_assemble_add_diagonal(&mut diag).unwrap();
7459df49d7eSJed Brown     ///
7469df49d7eSJed Brown     /// // Manual diagonal computation
7479df49d7eSJed Brown     /// let mut true_diag = ceed.vector(ndofs).unwrap();
7489df49d7eSJed Brown     /// for i in 0..ndofs {
7499df49d7eSJed Brown     ///     u.set_value(0.0);
7509df49d7eSJed Brown     ///     {
7519df49d7eSJed Brown     ///         let mut u_array = u.view_mut();
7529df49d7eSJed Brown     ///         u_array[i] = 1.;
7539df49d7eSJed Brown     ///     }
7549df49d7eSJed Brown     ///
7559df49d7eSJed Brown     ///     op_mass.apply(&u, &mut v).unwrap();
7569df49d7eSJed Brown     ///
7579df49d7eSJed Brown     ///     {
7589df49d7eSJed Brown     ///         let v_array = v.view_mut();
7599df49d7eSJed Brown     ///         let mut true_array = true_diag.view_mut();
7609df49d7eSJed Brown     ///         true_array[i] = v_array[i] + 1.0;
7619df49d7eSJed Brown     ///     }
7629df49d7eSJed Brown     /// }
7639df49d7eSJed Brown     ///
7649df49d7eSJed Brown     /// // Check
7659df49d7eSJed Brown     /// diag.view()
7669df49d7eSJed Brown     ///     .iter()
7679df49d7eSJed Brown     ///     .zip(true_diag.view().iter())
7689df49d7eSJed Brown     ///     .for_each(|(computed, actual)| {
7699df49d7eSJed Brown     ///         assert!(
770*80a9ef05SNatalie Beams     ///             (*computed - *actual).abs() < 10.0 * libceed::EPSILON,
7719df49d7eSJed Brown     ///             "Diagonal entry incorrect"
7729df49d7eSJed Brown     ///         );
7739df49d7eSJed Brown     ///     });
7749df49d7eSJed Brown     /// ```
7759df49d7eSJed Brown     pub fn linear_assemble_add_diagonal(&self, assembled: &mut Vector) -> crate::Result<i32> {
7769df49d7eSJed Brown         self.op_core.linear_assemble_add_diagonal(assembled)
7779df49d7eSJed Brown     }
7789df49d7eSJed Brown 
7799df49d7eSJed Brown     /// Assemble the point block diagonal of a square linear Operator
7809df49d7eSJed Brown     ///
7819df49d7eSJed Brown     /// This overwrites a Vector with the point block diagonal of a linear
7829df49d7eSJed Brown     /// Operator.
7839df49d7eSJed Brown     ///
7849df49d7eSJed Brown     /// Note: Currently only non-composite Operators with a single field and
7859df49d7eSJed Brown     ///      composite Operators with single field sub-operators are supported.
7869df49d7eSJed Brown     ///
7879df49d7eSJed Brown     /// * `op`        - Operator to assemble QFunction
7889df49d7eSJed Brown     /// * `assembled` - Vector to store assembled CeedOperator point block
7899df49d7eSJed Brown     ///                   diagonal, provided in row-major form with an
7909df49d7eSJed Brown     ///                   `ncomp * ncomp` block at each node. The dimensions of
7919df49d7eSJed Brown     ///                   this vector are derived from the active vector for
7929df49d7eSJed Brown     ///                   the CeedOperator. The array has shape
7939df49d7eSJed Brown     ///                   `[nodes, component out, component in]`.
7949df49d7eSJed Brown     ///
7959df49d7eSJed Brown     /// ```
7969df49d7eSJed Brown     /// # use libceed::prelude::*;
7979df49d7eSJed Brown     /// # let ceed = libceed::Ceed::default_init();
7989df49d7eSJed Brown     /// let ne = 4;
7999df49d7eSJed Brown     /// let p = 3;
8009df49d7eSJed Brown     /// let q = 4;
8019df49d7eSJed Brown     /// let ncomp = 2;
8029df49d7eSJed Brown     /// let ndofs = p * ne - ne + 1;
8039df49d7eSJed Brown     ///
8049df49d7eSJed Brown     /// // Vectors
8059df49d7eSJed Brown     /// let x = ceed.vector_from_slice(&[-1., -0.5, 0.0, 0.5, 1.0]).unwrap();
8069df49d7eSJed Brown     /// let mut qdata = ceed.vector(ne * q).unwrap();
8079df49d7eSJed Brown     /// qdata.set_value(0.0);
8089df49d7eSJed Brown     /// let mut u = ceed.vector(ncomp * ndofs).unwrap();
8099df49d7eSJed Brown     /// u.set_value(1.0);
8109df49d7eSJed Brown     /// let mut v = ceed.vector(ncomp * ndofs).unwrap();
8119df49d7eSJed Brown     /// v.set_value(0.0);
8129df49d7eSJed Brown     ///
8139df49d7eSJed Brown     /// // Restrictions
8149df49d7eSJed Brown     /// let mut indx: Vec<i32> = vec![0; 2 * ne];
8159df49d7eSJed Brown     /// for i in 0..ne {
8169df49d7eSJed Brown     ///     indx[2 * i + 0] = i as i32;
8179df49d7eSJed Brown     ///     indx[2 * i + 1] = (i + 1) as i32;
8189df49d7eSJed Brown     /// }
8199df49d7eSJed Brown     /// let rx = ceed
8209df49d7eSJed Brown     ///     .elem_restriction(ne, 2, 1, 1, ne + 1, MemType::Host, &indx)
8219df49d7eSJed Brown     ///     .unwrap();
8229df49d7eSJed Brown     /// let mut indu: Vec<i32> = vec![0; p * ne];
8239df49d7eSJed Brown     /// for i in 0..ne {
8249df49d7eSJed Brown     ///     indu[p * i + 0] = (2 * i) as i32;
8259df49d7eSJed Brown     ///     indu[p * i + 1] = (2 * i + 1) as i32;
8269df49d7eSJed Brown     ///     indu[p * i + 2] = (2 * i + 2) as i32;
8279df49d7eSJed Brown     /// }
8289df49d7eSJed Brown     /// let ru = ceed
8299df49d7eSJed Brown     ///     .elem_restriction(ne, p, ncomp, ndofs, ncomp * ndofs, MemType::Host, &indu)
8309df49d7eSJed Brown     ///     .unwrap();
8319df49d7eSJed Brown     /// let strides: [i32; 3] = [1, q as i32, q as i32];
8329df49d7eSJed Brown     /// let rq = ceed
8339df49d7eSJed Brown     ///     .strided_elem_restriction(ne, q, 1, q * ne, strides)
8349df49d7eSJed Brown     ///     .unwrap();
8359df49d7eSJed Brown     ///
8369df49d7eSJed Brown     /// // Bases
8379df49d7eSJed Brown     /// let bx = ceed
8389df49d7eSJed Brown     ///     .basis_tensor_H1_Lagrange(1, 1, 2, q, QuadMode::Gauss)
8399df49d7eSJed Brown     ///     .unwrap();
8409df49d7eSJed Brown     /// let bu = ceed
8419df49d7eSJed Brown     ///     .basis_tensor_H1_Lagrange(1, ncomp, p, q, QuadMode::Gauss)
8429df49d7eSJed Brown     ///     .unwrap();
8439df49d7eSJed Brown     ///
8449df49d7eSJed Brown     /// // Build quadrature data
8459df49d7eSJed Brown     /// let qf_build = ceed.q_function_interior_by_name("Mass1DBuild").unwrap();
8469df49d7eSJed Brown     /// ceed.operator(&qf_build, QFunctionOpt::None, QFunctionOpt::None)
8479df49d7eSJed Brown     ///     .unwrap()
8489df49d7eSJed Brown     ///     .field("dx", &rx, &bx, VectorOpt::Active)
8499df49d7eSJed Brown     ///     .unwrap()
8509df49d7eSJed Brown     ///     .field("weights", ElemRestrictionOpt::None, &bx, VectorOpt::None)
8519df49d7eSJed Brown     ///     .unwrap()
8529df49d7eSJed Brown     ///     .field("qdata", &rq, BasisOpt::Collocated, VectorOpt::Active)
8539df49d7eSJed Brown     ///     .unwrap()
8549df49d7eSJed Brown     ///     .apply(&x, &mut qdata)
8559df49d7eSJed Brown     ///     .unwrap();
8569df49d7eSJed Brown     ///
8579df49d7eSJed Brown     /// // Mass operator
8589df49d7eSJed Brown     /// let mut mass_2_comp = |[u, qdata, ..]: QFunctionInputs, [v, ..]: QFunctionOutputs| {
8599df49d7eSJed Brown     ///     // Number of quadrature points
8609df49d7eSJed Brown     ///     let q = qdata.len();
8619df49d7eSJed Brown     ///
8629df49d7eSJed Brown     ///     // Iterate over quadrature points
8639df49d7eSJed Brown     ///     for i in 0..q {
8649df49d7eSJed Brown     ///         v[i + 0 * q] = u[i + 1 * q] * qdata[i];
8659df49d7eSJed Brown     ///         v[i + 1 * q] = u[i + 0 * q] * qdata[i];
8669df49d7eSJed Brown     ///     }
8679df49d7eSJed Brown     ///
8689df49d7eSJed Brown     ///     // Return clean error code
8699df49d7eSJed Brown     ///     0
8709df49d7eSJed Brown     /// };
8719df49d7eSJed Brown     ///
8729df49d7eSJed Brown     /// let qf_mass = ceed
8739df49d7eSJed Brown     ///     .q_function_interior(1, Box::new(mass_2_comp))
8749df49d7eSJed Brown     ///     .unwrap()
8759df49d7eSJed Brown     ///     .input("u", 2, EvalMode::Interp)
8769df49d7eSJed Brown     ///     .unwrap()
8779df49d7eSJed Brown     ///     .input("qdata", 1, EvalMode::None)
8789df49d7eSJed Brown     ///     .unwrap()
8799df49d7eSJed Brown     ///     .output("v", 2, EvalMode::Interp)
8809df49d7eSJed Brown     ///     .unwrap();
8819df49d7eSJed Brown     ///
8829df49d7eSJed Brown     /// let op_mass = ceed
8839df49d7eSJed Brown     ///     .operator(&qf_mass, QFunctionOpt::None, QFunctionOpt::None)
8849df49d7eSJed Brown     ///     .unwrap()
8859df49d7eSJed Brown     ///     .field("u", &ru, &bu, VectorOpt::Active)
8869df49d7eSJed Brown     ///     .unwrap()
8879df49d7eSJed Brown     ///     .field("qdata", &rq, BasisOpt::Collocated, &qdata)
8889df49d7eSJed Brown     ///     .unwrap()
8899df49d7eSJed Brown     ///     .field("v", &ru, &bu, VectorOpt::Active)
8909df49d7eSJed Brown     ///     .unwrap();
8919df49d7eSJed Brown     ///
8929df49d7eSJed Brown     /// // Diagonal
8939df49d7eSJed Brown     /// let mut diag = ceed.vector(ncomp * ncomp * ndofs).unwrap();
8949df49d7eSJed Brown     /// diag.set_value(0.0);
8959df49d7eSJed Brown     /// op_mass
8969df49d7eSJed Brown     ///     .linear_assemble_point_block_diagonal(&mut diag)
8979df49d7eSJed Brown     ///     .unwrap();
8989df49d7eSJed Brown     ///
8999df49d7eSJed Brown     /// // Manual diagonal computation
9009df49d7eSJed Brown     /// let mut true_diag = ceed.vector(ncomp * ncomp * ndofs).unwrap();
9019df49d7eSJed Brown     /// for i in 0..ndofs {
9029df49d7eSJed Brown     ///     for j in 0..ncomp {
9039df49d7eSJed Brown     ///         u.set_value(0.0);
9049df49d7eSJed Brown     ///         {
9059df49d7eSJed Brown     ///             let mut u_array = u.view_mut();
9069df49d7eSJed Brown     ///             u_array[i + j * ndofs] = 1.;
9079df49d7eSJed Brown     ///         }
9089df49d7eSJed Brown     ///
9099df49d7eSJed Brown     ///         op_mass.apply(&u, &mut v).unwrap();
9109df49d7eSJed Brown     ///
9119df49d7eSJed Brown     ///         {
9129df49d7eSJed Brown     ///             let v_array = v.view_mut();
9139df49d7eSJed Brown     ///             let mut true_array = true_diag.view_mut();
9149df49d7eSJed Brown     ///             for k in 0..ncomp {
9159df49d7eSJed Brown     ///                 true_array[i * ncomp * ncomp + k * ncomp + j] = v_array[i + k * ndofs];
9169df49d7eSJed Brown     ///             }
9179df49d7eSJed Brown     ///         }
9189df49d7eSJed Brown     ///     }
9199df49d7eSJed Brown     /// }
9209df49d7eSJed Brown     ///
9219df49d7eSJed Brown     /// // Check
9229df49d7eSJed Brown     /// diag.view()
9239df49d7eSJed Brown     ///     .iter()
9249df49d7eSJed Brown     ///     .zip(true_diag.view().iter())
9259df49d7eSJed Brown     ///     .for_each(|(computed, actual)| {
9269df49d7eSJed Brown     ///         assert!(
927*80a9ef05SNatalie Beams     ///             (*computed - *actual).abs() < 10.0 * libceed::EPSILON,
9289df49d7eSJed Brown     ///             "Diagonal entry incorrect"
9299df49d7eSJed Brown     ///         );
9309df49d7eSJed Brown     ///     });
9319df49d7eSJed Brown     /// ```
9329df49d7eSJed Brown     pub fn linear_assemble_point_block_diagonal(
9339df49d7eSJed Brown         &self,
9349df49d7eSJed Brown         assembled: &mut Vector,
9359df49d7eSJed Brown     ) -> crate::Result<i32> {
9369df49d7eSJed Brown         self.op_core.linear_assemble_point_block_diagonal(assembled)
9379df49d7eSJed Brown     }
9389df49d7eSJed Brown 
9399df49d7eSJed Brown     /// Assemble the point block diagonal of a square linear Operator
9409df49d7eSJed Brown     ///
9419df49d7eSJed Brown     /// This sums into a Vector with the point block diagonal of a linear
9429df49d7eSJed Brown     /// Operator.
9439df49d7eSJed Brown     ///
9449df49d7eSJed Brown     /// Note: Currently only non-composite Operators with a single field and
9459df49d7eSJed Brown     ///      composite Operators with single field sub-operators are supported.
9469df49d7eSJed Brown     ///
9479df49d7eSJed Brown     /// * `op`        -     Operator to assemble QFunction
9489df49d7eSJed Brown     /// * `assembled` - Vector to store assembled CeedOperator point block
9499df49d7eSJed Brown     ///                   diagonal, provided in row-major form with an
9509df49d7eSJed Brown     ///                   `ncomp * ncomp` block at each node. The dimensions of
9519df49d7eSJed Brown     ///                   this vector are derived from the active vector for
9529df49d7eSJed Brown     ///                   the CeedOperator. The array has shape
9539df49d7eSJed Brown     ///                   `[nodes, component out, component in]`.
9549df49d7eSJed Brown     ///
9559df49d7eSJed Brown     /// ```
9569df49d7eSJed Brown     /// # use libceed::prelude::*;
9579df49d7eSJed Brown     /// # let ceed = libceed::Ceed::default_init();
9589df49d7eSJed Brown     /// let ne = 4;
9599df49d7eSJed Brown     /// let p = 3;
9609df49d7eSJed Brown     /// let q = 4;
9619df49d7eSJed Brown     /// let ncomp = 2;
9629df49d7eSJed Brown     /// let ndofs = p * ne - ne + 1;
9639df49d7eSJed Brown     ///
9649df49d7eSJed Brown     /// // Vectors
9659df49d7eSJed Brown     /// let x = ceed.vector_from_slice(&[-1., -0.5, 0.0, 0.5, 1.0]).unwrap();
9669df49d7eSJed Brown     /// let mut qdata = ceed.vector(ne * q).unwrap();
9679df49d7eSJed Brown     /// qdata.set_value(0.0);
9689df49d7eSJed Brown     /// let mut u = ceed.vector(ncomp * ndofs).unwrap();
9699df49d7eSJed Brown     /// u.set_value(1.0);
9709df49d7eSJed Brown     /// let mut v = ceed.vector(ncomp * ndofs).unwrap();
9719df49d7eSJed Brown     /// v.set_value(0.0);
9729df49d7eSJed Brown     ///
9739df49d7eSJed Brown     /// // Restrictions
9749df49d7eSJed Brown     /// let mut indx: Vec<i32> = vec![0; 2 * ne];
9759df49d7eSJed Brown     /// for i in 0..ne {
9769df49d7eSJed Brown     ///     indx[2 * i + 0] = i as i32;
9779df49d7eSJed Brown     ///     indx[2 * i + 1] = (i + 1) as i32;
9789df49d7eSJed Brown     /// }
9799df49d7eSJed Brown     /// let rx = ceed
9809df49d7eSJed Brown     ///     .elem_restriction(ne, 2, 1, 1, ne + 1, MemType::Host, &indx)
9819df49d7eSJed Brown     ///     .unwrap();
9829df49d7eSJed Brown     /// let mut indu: Vec<i32> = vec![0; p * ne];
9839df49d7eSJed Brown     /// for i in 0..ne {
9849df49d7eSJed Brown     ///     indu[p * i + 0] = (2 * i) as i32;
9859df49d7eSJed Brown     ///     indu[p * i + 1] = (2 * i + 1) as i32;
9869df49d7eSJed Brown     ///     indu[p * i + 2] = (2 * i + 2) as i32;
9879df49d7eSJed Brown     /// }
9889df49d7eSJed Brown     /// let ru = ceed
9899df49d7eSJed Brown     ///     .elem_restriction(ne, p, ncomp, ndofs, ncomp * ndofs, MemType::Host, &indu)
9909df49d7eSJed Brown     ///     .unwrap();
9919df49d7eSJed Brown     /// let strides: [i32; 3] = [1, q as i32, q as i32];
9929df49d7eSJed Brown     /// let rq = ceed
9939df49d7eSJed Brown     ///     .strided_elem_restriction(ne, q, 1, q * ne, strides)
9949df49d7eSJed Brown     ///     .unwrap();
9959df49d7eSJed Brown     ///
9969df49d7eSJed Brown     /// // Bases
9979df49d7eSJed Brown     /// let bx = ceed
9989df49d7eSJed Brown     ///     .basis_tensor_H1_Lagrange(1, 1, 2, q, QuadMode::Gauss)
9999df49d7eSJed Brown     ///     .unwrap();
10009df49d7eSJed Brown     /// let bu = ceed
10019df49d7eSJed Brown     ///     .basis_tensor_H1_Lagrange(1, ncomp, p, q, QuadMode::Gauss)
10029df49d7eSJed Brown     ///     .unwrap();
10039df49d7eSJed Brown     ///
10049df49d7eSJed Brown     /// // Build quadrature data
10059df49d7eSJed Brown     /// let qf_build = ceed.q_function_interior_by_name("Mass1DBuild").unwrap();
10069df49d7eSJed Brown     /// ceed.operator(&qf_build, QFunctionOpt::None, QFunctionOpt::None)
10079df49d7eSJed Brown     ///     .unwrap()
10089df49d7eSJed Brown     ///     .field("dx", &rx, &bx, VectorOpt::Active)
10099df49d7eSJed Brown     ///     .unwrap()
10109df49d7eSJed Brown     ///     .field("weights", ElemRestrictionOpt::None, &bx, VectorOpt::None)
10119df49d7eSJed Brown     ///     .unwrap()
10129df49d7eSJed Brown     ///     .field("qdata", &rq, BasisOpt::Collocated, VectorOpt::Active)
10139df49d7eSJed Brown     ///     .unwrap()
10149df49d7eSJed Brown     ///     .apply(&x, &mut qdata)
10159df49d7eSJed Brown     ///     .unwrap();
10169df49d7eSJed Brown     ///
10179df49d7eSJed Brown     /// // Mass operator
10189df49d7eSJed Brown     /// let mut mass_2_comp = |[u, qdata, ..]: QFunctionInputs, [v, ..]: QFunctionOutputs| {
10199df49d7eSJed Brown     ///     // Number of quadrature points
10209df49d7eSJed Brown     ///     let q = qdata.len();
10219df49d7eSJed Brown     ///
10229df49d7eSJed Brown     ///     // Iterate over quadrature points
10239df49d7eSJed Brown     ///     for i in 0..q {
10249df49d7eSJed Brown     ///         v[i + 0 * q] = u[i + 1 * q] * qdata[i];
10259df49d7eSJed Brown     ///         v[i + 1 * q] = u[i + 0 * q] * qdata[i];
10269df49d7eSJed Brown     ///     }
10279df49d7eSJed Brown     ///
10289df49d7eSJed Brown     ///     // Return clean error code
10299df49d7eSJed Brown     ///     0
10309df49d7eSJed Brown     /// };
10319df49d7eSJed Brown     ///
10329df49d7eSJed Brown     /// let qf_mass = ceed
10339df49d7eSJed Brown     ///     .q_function_interior(1, Box::new(mass_2_comp))
10349df49d7eSJed Brown     ///     .unwrap()
10359df49d7eSJed Brown     ///     .input("u", 2, EvalMode::Interp)
10369df49d7eSJed Brown     ///     .unwrap()
10379df49d7eSJed Brown     ///     .input("qdata", 1, EvalMode::None)
10389df49d7eSJed Brown     ///     .unwrap()
10399df49d7eSJed Brown     ///     .output("v", 2, EvalMode::Interp)
10409df49d7eSJed Brown     ///     .unwrap();
10419df49d7eSJed Brown     ///
10429df49d7eSJed Brown     /// let op_mass = ceed
10439df49d7eSJed Brown     ///     .operator(&qf_mass, QFunctionOpt::None, QFunctionOpt::None)
10449df49d7eSJed Brown     ///     .unwrap()
10459df49d7eSJed Brown     ///     .field("u", &ru, &bu, VectorOpt::Active)
10469df49d7eSJed Brown     ///     .unwrap()
10479df49d7eSJed Brown     ///     .field("qdata", &rq, BasisOpt::Collocated, &qdata)
10489df49d7eSJed Brown     ///     .unwrap()
10499df49d7eSJed Brown     ///     .field("v", &ru, &bu, VectorOpt::Active)
10509df49d7eSJed Brown     ///     .unwrap();
10519df49d7eSJed Brown     ///
10529df49d7eSJed Brown     /// // Diagonal
10539df49d7eSJed Brown     /// let mut diag = ceed.vector(ncomp * ncomp * ndofs).unwrap();
10549df49d7eSJed Brown     /// diag.set_value(1.0);
10559df49d7eSJed Brown     /// op_mass
10569df49d7eSJed Brown     ///     .linear_assemble_add_point_block_diagonal(&mut diag)
10579df49d7eSJed Brown     ///     .unwrap();
10589df49d7eSJed Brown     ///
10599df49d7eSJed Brown     /// // Manual diagonal computation
10609df49d7eSJed Brown     /// let mut true_diag = ceed.vector(ncomp * ncomp * ndofs).unwrap();
10619df49d7eSJed Brown     /// for i in 0..ndofs {
10629df49d7eSJed Brown     ///     for j in 0..ncomp {
10639df49d7eSJed Brown     ///         u.set_value(0.0);
10649df49d7eSJed Brown     ///         {
10659df49d7eSJed Brown     ///             let mut u_array = u.view_mut();
10669df49d7eSJed Brown     ///             u_array[i + j * ndofs] = 1.;
10679df49d7eSJed Brown     ///         }
10689df49d7eSJed Brown     ///
10699df49d7eSJed Brown     ///         op_mass.apply(&u, &mut v).unwrap();
10709df49d7eSJed Brown     ///
10719df49d7eSJed Brown     ///         {
10729df49d7eSJed Brown     ///             let v_array = v.view_mut();
10739df49d7eSJed Brown     ///             let mut true_array = true_diag.view_mut();
10749df49d7eSJed Brown     ///             for k in 0..ncomp {
10759df49d7eSJed Brown     ///                 true_array[i * ncomp * ncomp + k * ncomp + j] = v_array[i + k * ndofs];
10769df49d7eSJed Brown     ///             }
10779df49d7eSJed Brown     ///         }
10789df49d7eSJed Brown     ///     }
10799df49d7eSJed Brown     /// }
10809df49d7eSJed Brown     ///
10819df49d7eSJed Brown     /// // Check
10829df49d7eSJed Brown     /// diag.view()
10839df49d7eSJed Brown     ///     .iter()
10849df49d7eSJed Brown     ///     .zip(true_diag.view().iter())
10859df49d7eSJed Brown     ///     .for_each(|(computed, actual)| {
10869df49d7eSJed Brown     ///         assert!(
1087*80a9ef05SNatalie Beams     ///             (*computed - 1.0 - *actual).abs() < 10.0 * libceed::EPSILON,
10889df49d7eSJed Brown     ///             "Diagonal entry incorrect"
10899df49d7eSJed Brown     ///         );
10909df49d7eSJed Brown     ///     });
10919df49d7eSJed Brown     /// ```
10929df49d7eSJed Brown     pub fn linear_assemble_add_point_block_diagonal(
10939df49d7eSJed Brown         &self,
10949df49d7eSJed Brown         assembled: &mut Vector,
10959df49d7eSJed Brown     ) -> crate::Result<i32> {
10969df49d7eSJed Brown         self.op_core
10979df49d7eSJed Brown             .linear_assemble_add_point_block_diagonal(assembled)
10989df49d7eSJed Brown     }
10999df49d7eSJed Brown 
11009df49d7eSJed Brown     /// Create a multigrid coarse Operator and level transfer Operators for a
11019df49d7eSJed Brown     ///   given Operator, creating the prolongation basis from the fine and
11029df49d7eSJed Brown     ///   coarse grid interpolation
11039df49d7eSJed Brown     ///
11049df49d7eSJed Brown     /// * `p_mult_fine`  - Lvector multiplicity in parallel gather/scatter
11059df49d7eSJed Brown     /// * `rstr_coarse`  - Coarse grid restriction
11069df49d7eSJed Brown     /// * `basis_coarse` - Coarse grid active vector basis
11079df49d7eSJed Brown     ///
11089df49d7eSJed Brown     /// ```
11099df49d7eSJed Brown     /// # use libceed::prelude::*;
11109df49d7eSJed Brown     /// # let ceed = libceed::Ceed::default_init();
11119df49d7eSJed Brown     /// let ne = 15;
11129df49d7eSJed Brown     /// let p_coarse = 3;
11139df49d7eSJed Brown     /// let p_fine = 5;
11149df49d7eSJed Brown     /// let q = 6;
11159df49d7eSJed Brown     /// let ndofs_coarse = p_coarse * ne - ne + 1;
11169df49d7eSJed Brown     /// let ndofs_fine = p_fine * ne - ne + 1;
11179df49d7eSJed Brown     ///
11189df49d7eSJed Brown     /// // Vectors
11199df49d7eSJed Brown     /// let x_array = (0..ne + 1)
1120*80a9ef05SNatalie Beams     ///     .map(|i| 2.0 * i as Scalar / ne as Scalar - 1.0)
1121*80a9ef05SNatalie Beams     ///     .collect::<Vec<Scalar>>();
11229df49d7eSJed Brown     /// let x = ceed.vector_from_slice(&x_array).unwrap();
11239df49d7eSJed Brown     /// let mut qdata = ceed.vector(ne * q).unwrap();
11249df49d7eSJed Brown     /// qdata.set_value(0.0);
11259df49d7eSJed Brown     /// let mut u_coarse = ceed.vector(ndofs_coarse).unwrap();
11269df49d7eSJed Brown     /// u_coarse.set_value(1.0);
11279df49d7eSJed Brown     /// let mut u_fine = ceed.vector(ndofs_fine).unwrap();
11289df49d7eSJed Brown     /// u_fine.set_value(1.0);
11299df49d7eSJed Brown     /// let mut v_coarse = ceed.vector(ndofs_coarse).unwrap();
11309df49d7eSJed Brown     /// v_coarse.set_value(0.0);
11319df49d7eSJed Brown     /// let mut v_fine = ceed.vector(ndofs_fine).unwrap();
11329df49d7eSJed Brown     /// v_fine.set_value(0.0);
11339df49d7eSJed Brown     /// let mut multiplicity = ceed.vector(ndofs_fine).unwrap();
11349df49d7eSJed Brown     /// multiplicity.set_value(1.0);
11359df49d7eSJed Brown     ///
11369df49d7eSJed Brown     /// // Restrictions
11379df49d7eSJed Brown     /// let strides: [i32; 3] = [1, q as i32, q as i32];
11389df49d7eSJed Brown     /// let rq = ceed
11399df49d7eSJed Brown     ///     .strided_elem_restriction(ne, q, 1, q * ne, strides)
11409df49d7eSJed Brown     ///     .unwrap();
11419df49d7eSJed Brown     ///
11429df49d7eSJed Brown     /// let mut indx: Vec<i32> = vec![0; 2 * ne];
11439df49d7eSJed Brown     /// for i in 0..ne {
11449df49d7eSJed Brown     ///     indx[2 * i + 0] = i as i32;
11459df49d7eSJed Brown     ///     indx[2 * i + 1] = (i + 1) as i32;
11469df49d7eSJed Brown     /// }
11479df49d7eSJed Brown     /// let rx = ceed
11489df49d7eSJed Brown     ///     .elem_restriction(ne, 2, 1, 1, ne + 1, MemType::Host, &indx)
11499df49d7eSJed Brown     ///     .unwrap();
11509df49d7eSJed Brown     ///
11519df49d7eSJed Brown     /// let mut indu_coarse: Vec<i32> = vec![0; p_coarse * ne];
11529df49d7eSJed Brown     /// for i in 0..ne {
11539df49d7eSJed Brown     ///     for j in 0..p_coarse {
11549df49d7eSJed Brown     ///         indu_coarse[p_coarse * i + j] = (i + j) as i32;
11559df49d7eSJed Brown     ///     }
11569df49d7eSJed Brown     /// }
11579df49d7eSJed Brown     /// let ru_coarse = ceed
11589df49d7eSJed Brown     ///     .elem_restriction(
11599df49d7eSJed Brown     ///         ne,
11609df49d7eSJed Brown     ///         p_coarse,
11619df49d7eSJed Brown     ///         1,
11629df49d7eSJed Brown     ///         1,
11639df49d7eSJed Brown     ///         ndofs_coarse,
11649df49d7eSJed Brown     ///         MemType::Host,
11659df49d7eSJed Brown     ///         &indu_coarse,
11669df49d7eSJed Brown     ///     )
11679df49d7eSJed Brown     ///     .unwrap();
11689df49d7eSJed Brown     ///
11699df49d7eSJed Brown     /// let mut indu_fine: Vec<i32> = vec![0; p_fine * ne];
11709df49d7eSJed Brown     /// for i in 0..ne {
11719df49d7eSJed Brown     ///     for j in 0..p_fine {
11729df49d7eSJed Brown     ///         indu_fine[p_fine * i + j] = (i + j) as i32;
11739df49d7eSJed Brown     ///     }
11749df49d7eSJed Brown     /// }
11759df49d7eSJed Brown     /// let ru_fine = ceed
11769df49d7eSJed Brown     ///     .elem_restriction(ne, p_fine, 1, 1, ndofs_fine, MemType::Host, &indu_fine)
11779df49d7eSJed Brown     ///     .unwrap();
11789df49d7eSJed Brown     ///
11799df49d7eSJed Brown     /// // Bases
11809df49d7eSJed Brown     /// let bx = ceed
11819df49d7eSJed Brown     ///     .basis_tensor_H1_Lagrange(1, 1, 2, q, QuadMode::Gauss)
11829df49d7eSJed Brown     ///     .unwrap();
11839df49d7eSJed Brown     /// let bu_coarse = ceed
11849df49d7eSJed Brown     ///     .basis_tensor_H1_Lagrange(1, 1, p_coarse, q, QuadMode::Gauss)
11859df49d7eSJed Brown     ///     .unwrap();
11869df49d7eSJed Brown     /// let bu_fine = ceed
11879df49d7eSJed Brown     ///     .basis_tensor_H1_Lagrange(1, 1, p_fine, q, QuadMode::Gauss)
11889df49d7eSJed Brown     ///     .unwrap();
11899df49d7eSJed Brown     ///
11909df49d7eSJed Brown     /// // Build quadrature data
11919df49d7eSJed Brown     /// let qf_build = ceed.q_function_interior_by_name("Mass1DBuild").unwrap();
11929df49d7eSJed Brown     /// ceed.operator(&qf_build, QFunctionOpt::None, QFunctionOpt::None)
11939df49d7eSJed Brown     ///     .unwrap()
11949df49d7eSJed Brown     ///     .field("dx", &rx, &bx, VectorOpt::Active)
11959df49d7eSJed Brown     ///     .unwrap()
11969df49d7eSJed Brown     ///     .field("weights", ElemRestrictionOpt::None, &bx, VectorOpt::None)
11979df49d7eSJed Brown     ///     .unwrap()
11989df49d7eSJed Brown     ///     .field("qdata", &rq, BasisOpt::Collocated, VectorOpt::Active)
11999df49d7eSJed Brown     ///     .unwrap()
12009df49d7eSJed Brown     ///     .apply(&x, &mut qdata)
12019df49d7eSJed Brown     ///     .unwrap();
12029df49d7eSJed Brown     ///
12039df49d7eSJed Brown     /// // Mass operator
12049df49d7eSJed Brown     /// let qf_mass = ceed.q_function_interior_by_name("MassApply").unwrap();
12059df49d7eSJed Brown     /// let op_mass_fine = ceed
12069df49d7eSJed Brown     ///     .operator(&qf_mass, QFunctionOpt::None, QFunctionOpt::None)
12079df49d7eSJed Brown     ///     .unwrap()
12089df49d7eSJed Brown     ///     .field("u", &ru_fine, &bu_fine, VectorOpt::Active)
12099df49d7eSJed Brown     ///     .unwrap()
12109df49d7eSJed Brown     ///     .field("qdata", &rq, BasisOpt::Collocated, &qdata)
12119df49d7eSJed Brown     ///     .unwrap()
12129df49d7eSJed Brown     ///     .field("v", &ru_fine, &bu_fine, VectorOpt::Active)
12139df49d7eSJed Brown     ///     .unwrap();
12149df49d7eSJed Brown     ///
12159df49d7eSJed Brown     /// // Multigrid setup
12169df49d7eSJed Brown     /// let (op_mass_coarse, op_prolong, op_restrict) = op_mass_fine
12179df49d7eSJed Brown     ///     .create_multigrid_level(&multiplicity, &ru_coarse, &bu_coarse)
12189df49d7eSJed Brown     ///     .unwrap();
12199df49d7eSJed Brown     ///
12209df49d7eSJed Brown     /// // Coarse problem
12219df49d7eSJed Brown     /// u_coarse.set_value(1.0);
12229df49d7eSJed Brown     /// op_mass_coarse.apply(&u_coarse, &mut v_coarse).unwrap();
12239df49d7eSJed Brown     ///
12249df49d7eSJed Brown     /// // Check
1225*80a9ef05SNatalie Beams     /// let sum: Scalar = v_coarse.view().iter().sum();
12269df49d7eSJed Brown     /// assert!(
1227*80a9ef05SNatalie Beams     ///     (sum - 2.0).abs() < 50.0 * libceed::EPSILON,
12289df49d7eSJed Brown     ///     "Incorrect interval length computed"
12299df49d7eSJed Brown     /// );
12309df49d7eSJed Brown     ///
12319df49d7eSJed Brown     /// // Prolong
12329df49d7eSJed Brown     /// op_prolong.apply(&u_coarse, &mut u_fine).unwrap();
12339df49d7eSJed Brown     ///
12349df49d7eSJed Brown     /// // Fine problem
12359df49d7eSJed Brown     /// op_mass_fine.apply(&u_fine, &mut v_fine).unwrap();
12369df49d7eSJed Brown     ///
12379df49d7eSJed Brown     /// // Check
1238*80a9ef05SNatalie Beams     /// let sum: Scalar = v_fine.view().iter().sum();
12399df49d7eSJed Brown     /// assert!(
1240*80a9ef05SNatalie Beams     ///     (sum - 2.0).abs() < 50.0 * libceed::EPSILON,
12419df49d7eSJed Brown     ///     "Incorrect interval length computed"
12429df49d7eSJed Brown     /// );
12439df49d7eSJed Brown     ///
12449df49d7eSJed Brown     /// // Restrict
12459df49d7eSJed Brown     /// op_restrict.apply(&v_fine, &mut v_coarse).unwrap();
12469df49d7eSJed Brown     ///
12479df49d7eSJed Brown     /// // Check
1248*80a9ef05SNatalie Beams     /// let sum: Scalar = v_coarse.view().iter().sum();
12499df49d7eSJed Brown     /// assert!(
1250*80a9ef05SNatalie Beams     ///     (sum - 2.0).abs() < 50.0 * libceed::EPSILON,
12519df49d7eSJed Brown     ///     "Incorrect interval length computed"
12529df49d7eSJed Brown     /// );
12539df49d7eSJed Brown     /// ```
12549df49d7eSJed Brown     pub fn create_multigrid_level(
12559df49d7eSJed Brown         &self,
12569df49d7eSJed Brown         p_mult_fine: &Vector,
12579df49d7eSJed Brown         rstr_coarse: &ElemRestriction,
12589df49d7eSJed Brown         basis_coarse: &Basis,
12599df49d7eSJed Brown     ) -> crate::Result<(Operator, Operator, Operator)> {
12609df49d7eSJed Brown         let mut ptr_coarse = std::ptr::null_mut();
12619df49d7eSJed Brown         let mut ptr_prolong = std::ptr::null_mut();
12629df49d7eSJed Brown         let mut ptr_restrict = std::ptr::null_mut();
12639df49d7eSJed Brown         let ierr = unsafe {
12649df49d7eSJed Brown             bind_ceed::CeedOperatorMultigridLevelCreate(
12659df49d7eSJed Brown                 self.op_core.ptr,
12669df49d7eSJed Brown                 p_mult_fine.ptr,
12679df49d7eSJed Brown                 rstr_coarse.ptr,
12689df49d7eSJed Brown                 basis_coarse.ptr,
12699df49d7eSJed Brown                 &mut ptr_coarse,
12709df49d7eSJed Brown                 &mut ptr_prolong,
12719df49d7eSJed Brown                 &mut ptr_restrict,
12729df49d7eSJed Brown             )
12739df49d7eSJed Brown         };
12749df49d7eSJed Brown         self.op_core.ceed.check_error(ierr)?;
12759df49d7eSJed Brown         let op_coarse = Operator::from_raw(self.op_core.ceed, ptr_coarse)?;
12769df49d7eSJed Brown         let op_prolong = Operator::from_raw(self.op_core.ceed, ptr_prolong)?;
12779df49d7eSJed Brown         let op_restrict = Operator::from_raw(self.op_core.ceed, ptr_restrict)?;
12789df49d7eSJed Brown         Ok((op_coarse, op_prolong, op_restrict))
12799df49d7eSJed Brown     }
12809df49d7eSJed Brown 
12819df49d7eSJed Brown     /// Create a multigrid coarse Operator and level transfer Operators for a
12829df49d7eSJed Brown     ///   given Operator with a tensor basis for the active basis
12839df49d7eSJed Brown     ///
12849df49d7eSJed Brown     /// * `p_mult_fine`   - Lvector multiplicity in parallel gather/scatter
12859df49d7eSJed Brown     /// * `rstr_coarse`   - Coarse grid restriction
12869df49d7eSJed Brown     /// * `basis_coarse`  - Coarse grid active vector basis
12879df49d7eSJed Brown     /// * `interp_c_to_f` - Matrix for coarse to fine
12889df49d7eSJed Brown     ///
12899df49d7eSJed Brown     /// ```
12909df49d7eSJed Brown     /// # use libceed::prelude::*;
12919df49d7eSJed Brown     /// # let ceed = libceed::Ceed::default_init();
12929df49d7eSJed Brown     /// let ne = 15;
12939df49d7eSJed Brown     /// let p_coarse = 3;
12949df49d7eSJed Brown     /// let p_fine = 5;
12959df49d7eSJed Brown     /// let q = 6;
12969df49d7eSJed Brown     /// let ndofs_coarse = p_coarse * ne - ne + 1;
12979df49d7eSJed Brown     /// let ndofs_fine = p_fine * ne - ne + 1;
12989df49d7eSJed Brown     ///
12999df49d7eSJed Brown     /// // Vectors
13009df49d7eSJed Brown     /// let x_array = (0..ne + 1)
1301*80a9ef05SNatalie Beams     ///     .map(|i| 2.0 * i as Scalar / ne as Scalar - 1.0)
1302*80a9ef05SNatalie Beams     ///     .collect::<Vec<Scalar>>();
13039df49d7eSJed Brown     /// let x = ceed.vector_from_slice(&x_array).unwrap();
13049df49d7eSJed Brown     /// let mut qdata = ceed.vector(ne * q).unwrap();
13059df49d7eSJed Brown     /// qdata.set_value(0.0);
13069df49d7eSJed Brown     /// let mut u_coarse = ceed.vector(ndofs_coarse).unwrap();
13079df49d7eSJed Brown     /// u_coarse.set_value(1.0);
13089df49d7eSJed Brown     /// let mut u_fine = ceed.vector(ndofs_fine).unwrap();
13099df49d7eSJed Brown     /// u_fine.set_value(1.0);
13109df49d7eSJed Brown     /// let mut v_coarse = ceed.vector(ndofs_coarse).unwrap();
13119df49d7eSJed Brown     /// v_coarse.set_value(0.0);
13129df49d7eSJed Brown     /// let mut v_fine = ceed.vector(ndofs_fine).unwrap();
13139df49d7eSJed Brown     /// v_fine.set_value(0.0);
13149df49d7eSJed Brown     /// let mut multiplicity = ceed.vector(ndofs_fine).unwrap();
13159df49d7eSJed Brown     /// multiplicity.set_value(1.0);
13169df49d7eSJed Brown     ///
13179df49d7eSJed Brown     /// // Restrictions
13189df49d7eSJed Brown     /// let strides: [i32; 3] = [1, q as i32, q as i32];
13199df49d7eSJed Brown     /// let rq = ceed
13209df49d7eSJed Brown     ///     .strided_elem_restriction(ne, q, 1, q * ne, strides)
13219df49d7eSJed Brown     ///     .unwrap();
13229df49d7eSJed Brown     ///
13239df49d7eSJed Brown     /// let mut indx: Vec<i32> = vec![0; 2 * ne];
13249df49d7eSJed Brown     /// for i in 0..ne {
13259df49d7eSJed Brown     ///     indx[2 * i + 0] = i as i32;
13269df49d7eSJed Brown     ///     indx[2 * i + 1] = (i + 1) as i32;
13279df49d7eSJed Brown     /// }
13289df49d7eSJed Brown     /// let rx = ceed
13299df49d7eSJed Brown     ///     .elem_restriction(ne, 2, 1, 1, ne + 1, MemType::Host, &indx)
13309df49d7eSJed Brown     ///     .unwrap();
13319df49d7eSJed Brown     ///
13329df49d7eSJed Brown     /// let mut indu_coarse: Vec<i32> = vec![0; p_coarse * ne];
13339df49d7eSJed Brown     /// for i in 0..ne {
13349df49d7eSJed Brown     ///     for j in 0..p_coarse {
13359df49d7eSJed Brown     ///         indu_coarse[p_coarse * i + j] = (i + j) as i32;
13369df49d7eSJed Brown     ///     }
13379df49d7eSJed Brown     /// }
13389df49d7eSJed Brown     /// let ru_coarse = ceed
13399df49d7eSJed Brown     ///     .elem_restriction(
13409df49d7eSJed Brown     ///         ne,
13419df49d7eSJed Brown     ///         p_coarse,
13429df49d7eSJed Brown     ///         1,
13439df49d7eSJed Brown     ///         1,
13449df49d7eSJed Brown     ///         ndofs_coarse,
13459df49d7eSJed Brown     ///         MemType::Host,
13469df49d7eSJed Brown     ///         &indu_coarse,
13479df49d7eSJed Brown     ///     )
13489df49d7eSJed Brown     ///     .unwrap();
13499df49d7eSJed Brown     ///
13509df49d7eSJed Brown     /// let mut indu_fine: Vec<i32> = vec![0; p_fine * ne];
13519df49d7eSJed Brown     /// for i in 0..ne {
13529df49d7eSJed Brown     ///     for j in 0..p_fine {
13539df49d7eSJed Brown     ///         indu_fine[p_fine * i + j] = (i + j) as i32;
13549df49d7eSJed Brown     ///     }
13559df49d7eSJed Brown     /// }
13569df49d7eSJed Brown     /// let ru_fine = ceed
13579df49d7eSJed Brown     ///     .elem_restriction(ne, p_fine, 1, 1, ndofs_fine, MemType::Host, &indu_fine)
13589df49d7eSJed Brown     ///     .unwrap();
13599df49d7eSJed Brown     ///
13609df49d7eSJed Brown     /// // Bases
13619df49d7eSJed Brown     /// let bx = ceed
13629df49d7eSJed Brown     ///     .basis_tensor_H1_Lagrange(1, 1, 2, q, QuadMode::Gauss)
13639df49d7eSJed Brown     ///     .unwrap();
13649df49d7eSJed Brown     /// let bu_coarse = ceed
13659df49d7eSJed Brown     ///     .basis_tensor_H1_Lagrange(1, 1, p_coarse, q, QuadMode::Gauss)
13669df49d7eSJed Brown     ///     .unwrap();
13679df49d7eSJed Brown     /// let bu_fine = ceed
13689df49d7eSJed Brown     ///     .basis_tensor_H1_Lagrange(1, 1, p_fine, q, QuadMode::Gauss)
13699df49d7eSJed Brown     ///     .unwrap();
13709df49d7eSJed Brown     ///
13719df49d7eSJed Brown     /// // Build quadrature data
13729df49d7eSJed Brown     /// let qf_build = ceed.q_function_interior_by_name("Mass1DBuild").unwrap();
13739df49d7eSJed Brown     /// ceed.operator(&qf_build, QFunctionOpt::None, QFunctionOpt::None)
13749df49d7eSJed Brown     ///     .unwrap()
13759df49d7eSJed Brown     ///     .field("dx", &rx, &bx, VectorOpt::Active)
13769df49d7eSJed Brown     ///     .unwrap()
13779df49d7eSJed Brown     ///     .field("weights", ElemRestrictionOpt::None, &bx, VectorOpt::None)
13789df49d7eSJed Brown     ///     .unwrap()
13799df49d7eSJed Brown     ///     .field("qdata", &rq, BasisOpt::Collocated, VectorOpt::Active)
13809df49d7eSJed Brown     ///     .unwrap()
13819df49d7eSJed Brown     ///     .apply(&x, &mut qdata)
13829df49d7eSJed Brown     ///     .unwrap();
13839df49d7eSJed Brown     ///
13849df49d7eSJed Brown     /// // Mass operator
13859df49d7eSJed Brown     /// let qf_mass = ceed.q_function_interior_by_name("MassApply").unwrap();
13869df49d7eSJed Brown     /// let op_mass_fine = ceed
13879df49d7eSJed Brown     ///     .operator(&qf_mass, QFunctionOpt::None, QFunctionOpt::None)
13889df49d7eSJed Brown     ///     .unwrap()
13899df49d7eSJed Brown     ///     .field("u", &ru_fine, &bu_fine, VectorOpt::Active)
13909df49d7eSJed Brown     ///     .unwrap()
13919df49d7eSJed Brown     ///     .field("qdata", &rq, BasisOpt::Collocated, &qdata)
13929df49d7eSJed Brown     ///     .unwrap()
13939df49d7eSJed Brown     ///     .field("v", &ru_fine, &bu_fine, VectorOpt::Active)
13949df49d7eSJed Brown     ///     .unwrap();
13959df49d7eSJed Brown     ///
13969df49d7eSJed Brown     /// // Multigrid setup
1397*80a9ef05SNatalie Beams     /// let mut interp_c_to_f: Vec<Scalar> = vec![0.; p_coarse * p_fine];
13989df49d7eSJed Brown     /// {
13999df49d7eSJed Brown     ///     let mut coarse = ceed.vector(p_coarse).unwrap();
14009df49d7eSJed Brown     ///     let mut fine = ceed.vector(p_fine).unwrap();
14019df49d7eSJed Brown     ///     let basis_c_to_f = ceed
14029df49d7eSJed Brown     ///         .basis_tensor_H1_Lagrange(1, 1, p_coarse, p_fine, QuadMode::GaussLobatto)
14039df49d7eSJed Brown     ///         .unwrap();
14049df49d7eSJed Brown     ///     for i in 0..p_coarse {
14059df49d7eSJed Brown     ///         coarse.set_value(0.0);
14069df49d7eSJed Brown     ///         {
14079df49d7eSJed Brown     ///             let mut array = coarse.view_mut();
14089df49d7eSJed Brown     ///             array[i] = 1.;
14099df49d7eSJed Brown     ///         }
14109df49d7eSJed Brown     ///         basis_c_to_f
14119df49d7eSJed Brown     ///             .apply(
14129df49d7eSJed Brown     ///                 1,
14139df49d7eSJed Brown     ///                 TransposeMode::NoTranspose,
14149df49d7eSJed Brown     ///                 EvalMode::Interp,
14159df49d7eSJed Brown     ///                 &coarse,
14169df49d7eSJed Brown     ///                 &mut fine,
14179df49d7eSJed Brown     ///             )
14189df49d7eSJed Brown     ///             .unwrap();
14199df49d7eSJed Brown     ///         let array = fine.view();
14209df49d7eSJed Brown     ///         for j in 0..p_fine {
14219df49d7eSJed Brown     ///             interp_c_to_f[j * p_coarse + i] = array[j];
14229df49d7eSJed Brown     ///         }
14239df49d7eSJed Brown     ///     }
14249df49d7eSJed Brown     /// }
14259df49d7eSJed Brown     /// let (op_mass_coarse, op_prolong, op_restrict) = op_mass_fine
14269df49d7eSJed Brown     ///     .create_multigrid_level_tensor_H1(&multiplicity, &ru_coarse, &bu_coarse, &interp_c_to_f)
14279df49d7eSJed Brown     ///     .unwrap();
14289df49d7eSJed Brown     ///
14299df49d7eSJed Brown     /// // Coarse problem
14309df49d7eSJed Brown     /// u_coarse.set_value(1.0);
14319df49d7eSJed Brown     /// op_mass_coarse.apply(&u_coarse, &mut v_coarse).unwrap();
14329df49d7eSJed Brown     ///
14339df49d7eSJed Brown     /// // Check
1434*80a9ef05SNatalie Beams     /// let sum: Scalar = v_coarse.view().iter().sum();
14359df49d7eSJed Brown     /// assert!(
1436*80a9ef05SNatalie Beams     ///     (sum - 2.0).abs() < 10.0 * libceed::EPSILON,
14379df49d7eSJed Brown     ///     "Incorrect interval length computed"
14389df49d7eSJed Brown     /// );
14399df49d7eSJed Brown     ///
14409df49d7eSJed Brown     /// // Prolong
14419df49d7eSJed Brown     /// op_prolong.apply(&u_coarse, &mut u_fine).unwrap();
14429df49d7eSJed Brown     ///
14439df49d7eSJed Brown     /// // Fine problem
14449df49d7eSJed Brown     /// op_mass_fine.apply(&u_fine, &mut v_fine).unwrap();
14459df49d7eSJed Brown     ///
14469df49d7eSJed Brown     /// // Check
1447*80a9ef05SNatalie Beams     /// let sum: Scalar = v_fine.view().iter().sum();
14489df49d7eSJed Brown     /// assert!(
1449*80a9ef05SNatalie Beams     ///     (sum - 2.0).abs() < 10.0 * libceed::EPSILON,
14509df49d7eSJed Brown     ///     "Incorrect interval length computed"
14519df49d7eSJed Brown     /// );
14529df49d7eSJed Brown     ///
14539df49d7eSJed Brown     /// // Restrict
14549df49d7eSJed Brown     /// op_restrict.apply(&v_fine, &mut v_coarse).unwrap();
14559df49d7eSJed Brown     ///
14569df49d7eSJed Brown     /// // Check
1457*80a9ef05SNatalie Beams     /// let sum: Scalar = v_coarse.view().iter().sum();
14589df49d7eSJed Brown     /// assert!(
1459*80a9ef05SNatalie Beams     ///     (sum - 2.0).abs() < 10.0 * libceed::EPSILON,
14609df49d7eSJed Brown     ///     "Incorrect interval length computed"
14619df49d7eSJed Brown     /// );
14629df49d7eSJed Brown     /// ```
14639df49d7eSJed Brown     pub fn create_multigrid_level_tensor_H1(
14649df49d7eSJed Brown         &self,
14659df49d7eSJed Brown         p_mult_fine: &Vector,
14669df49d7eSJed Brown         rstr_coarse: &ElemRestriction,
14679df49d7eSJed Brown         basis_coarse: &Basis,
1468*80a9ef05SNatalie Beams         interpCtoF: &Vec<Scalar>,
14699df49d7eSJed Brown     ) -> crate::Result<(Operator, Operator, Operator)> {
14709df49d7eSJed Brown         let mut ptr_coarse = std::ptr::null_mut();
14719df49d7eSJed Brown         let mut ptr_prolong = std::ptr::null_mut();
14729df49d7eSJed Brown         let mut ptr_restrict = std::ptr::null_mut();
14739df49d7eSJed Brown         let ierr = unsafe {
14749df49d7eSJed Brown             bind_ceed::CeedOperatorMultigridLevelCreateTensorH1(
14759df49d7eSJed Brown                 self.op_core.ptr,
14769df49d7eSJed Brown                 p_mult_fine.ptr,
14779df49d7eSJed Brown                 rstr_coarse.ptr,
14789df49d7eSJed Brown                 basis_coarse.ptr,
14799df49d7eSJed Brown                 interpCtoF.as_ptr(),
14809df49d7eSJed Brown                 &mut ptr_coarse,
14819df49d7eSJed Brown                 &mut ptr_prolong,
14829df49d7eSJed Brown                 &mut ptr_restrict,
14839df49d7eSJed Brown             )
14849df49d7eSJed Brown         };
14859df49d7eSJed Brown         self.op_core.ceed.check_error(ierr)?;
14869df49d7eSJed Brown         let op_coarse = Operator::from_raw(self.op_core.ceed, ptr_coarse)?;
14879df49d7eSJed Brown         let op_prolong = Operator::from_raw(self.op_core.ceed, ptr_prolong)?;
14889df49d7eSJed Brown         let op_restrict = Operator::from_raw(self.op_core.ceed, ptr_restrict)?;
14899df49d7eSJed Brown         Ok((op_coarse, op_prolong, op_restrict))
14909df49d7eSJed Brown     }
14919df49d7eSJed Brown 
14929df49d7eSJed Brown     /// Create a multigrid coarse Operator and level transfer Operators for a
14939df49d7eSJed Brown     ///   given Operator with a non-tensor basis for the active basis
14949df49d7eSJed Brown     ///
14959df49d7eSJed Brown     /// * `p_mult_fine`   - Lvector multiplicity in parallel gather/scatter
14969df49d7eSJed Brown     /// * `rstr_coarse`   - Coarse grid restriction
14979df49d7eSJed Brown     /// * `basis_coarse`  - Coarse grid active vector basis
14989df49d7eSJed Brown     /// * `interp_c_to_f` - Matrix for coarse to fine
14999df49d7eSJed Brown     ///
15009df49d7eSJed Brown     /// ```
15019df49d7eSJed Brown     /// # use libceed::prelude::*;
15029df49d7eSJed Brown     /// # let ceed = libceed::Ceed::default_init();
15039df49d7eSJed Brown     /// let ne = 15;
15049df49d7eSJed Brown     /// let p_coarse = 3;
15059df49d7eSJed Brown     /// let p_fine = 5;
15069df49d7eSJed Brown     /// let q = 6;
15079df49d7eSJed Brown     /// let ndofs_coarse = p_coarse * ne - ne + 1;
15089df49d7eSJed Brown     /// let ndofs_fine = p_fine * ne - ne + 1;
15099df49d7eSJed Brown     ///
15109df49d7eSJed Brown     /// // Vectors
15119df49d7eSJed Brown     /// let x_array = (0..ne + 1)
1512*80a9ef05SNatalie Beams     ///     .map(|i| 2.0 * i as Scalar / ne as Scalar - 1.0)
1513*80a9ef05SNatalie Beams     ///     .collect::<Vec<Scalar>>();
15149df49d7eSJed Brown     /// let x = ceed.vector_from_slice(&x_array).unwrap();
15159df49d7eSJed Brown     /// let mut qdata = ceed.vector(ne * q).unwrap();
15169df49d7eSJed Brown     /// qdata.set_value(0.0);
15179df49d7eSJed Brown     /// let mut u_coarse = ceed.vector(ndofs_coarse).unwrap();
15189df49d7eSJed Brown     /// u_coarse.set_value(1.0);
15199df49d7eSJed Brown     /// let mut u_fine = ceed.vector(ndofs_fine).unwrap();
15209df49d7eSJed Brown     /// u_fine.set_value(1.0);
15219df49d7eSJed Brown     /// let mut v_coarse = ceed.vector(ndofs_coarse).unwrap();
15229df49d7eSJed Brown     /// v_coarse.set_value(0.0);
15239df49d7eSJed Brown     /// let mut v_fine = ceed.vector(ndofs_fine).unwrap();
15249df49d7eSJed Brown     /// v_fine.set_value(0.0);
15259df49d7eSJed Brown     /// let mut multiplicity = ceed.vector(ndofs_fine).unwrap();
15269df49d7eSJed Brown     /// multiplicity.set_value(1.0);
15279df49d7eSJed Brown     ///
15289df49d7eSJed Brown     /// // Restrictions
15299df49d7eSJed Brown     /// let strides: [i32; 3] = [1, q as i32, q as i32];
15309df49d7eSJed Brown     /// let rq = ceed
15319df49d7eSJed Brown     ///     .strided_elem_restriction(ne, q, 1, q * ne, strides)
15329df49d7eSJed Brown     ///     .unwrap();
15339df49d7eSJed Brown     ///
15349df49d7eSJed Brown     /// let mut indx: Vec<i32> = vec![0; 2 * ne];
15359df49d7eSJed Brown     /// for i in 0..ne {
15369df49d7eSJed Brown     ///     indx[2 * i + 0] = i as i32;
15379df49d7eSJed Brown     ///     indx[2 * i + 1] = (i + 1) as i32;
15389df49d7eSJed Brown     /// }
15399df49d7eSJed Brown     /// let rx = ceed
15409df49d7eSJed Brown     ///     .elem_restriction(ne, 2, 1, 1, ne + 1, MemType::Host, &indx)
15419df49d7eSJed Brown     ///     .unwrap();
15429df49d7eSJed Brown     ///
15439df49d7eSJed Brown     /// let mut indu_coarse: Vec<i32> = vec![0; p_coarse * ne];
15449df49d7eSJed Brown     /// for i in 0..ne {
15459df49d7eSJed Brown     ///     for j in 0..p_coarse {
15469df49d7eSJed Brown     ///         indu_coarse[p_coarse * i + j] = (i + j) as i32;
15479df49d7eSJed Brown     ///     }
15489df49d7eSJed Brown     /// }
15499df49d7eSJed Brown     /// let ru_coarse = ceed
15509df49d7eSJed Brown     ///     .elem_restriction(
15519df49d7eSJed Brown     ///         ne,
15529df49d7eSJed Brown     ///         p_coarse,
15539df49d7eSJed Brown     ///         1,
15549df49d7eSJed Brown     ///         1,
15559df49d7eSJed Brown     ///         ndofs_coarse,
15569df49d7eSJed Brown     ///         MemType::Host,
15579df49d7eSJed Brown     ///         &indu_coarse,
15589df49d7eSJed Brown     ///     )
15599df49d7eSJed Brown     ///     .unwrap();
15609df49d7eSJed Brown     ///
15619df49d7eSJed Brown     /// let mut indu_fine: Vec<i32> = vec![0; p_fine * ne];
15629df49d7eSJed Brown     /// for i in 0..ne {
15639df49d7eSJed Brown     ///     for j in 0..p_fine {
15649df49d7eSJed Brown     ///         indu_fine[p_fine * i + j] = (i + j) as i32;
15659df49d7eSJed Brown     ///     }
15669df49d7eSJed Brown     /// }
15679df49d7eSJed Brown     /// let ru_fine = ceed
15689df49d7eSJed Brown     ///     .elem_restriction(ne, p_fine, 1, 1, ndofs_fine, MemType::Host, &indu_fine)
15699df49d7eSJed Brown     ///     .unwrap();
15709df49d7eSJed Brown     ///
15719df49d7eSJed Brown     /// // Bases
15729df49d7eSJed Brown     /// let bx = ceed
15739df49d7eSJed Brown     ///     .basis_tensor_H1_Lagrange(1, 1, 2, q, QuadMode::Gauss)
15749df49d7eSJed Brown     ///     .unwrap();
15759df49d7eSJed Brown     /// let bu_coarse = ceed
15769df49d7eSJed Brown     ///     .basis_tensor_H1_Lagrange(1, 1, p_coarse, q, QuadMode::Gauss)
15779df49d7eSJed Brown     ///     .unwrap();
15789df49d7eSJed Brown     /// let bu_fine = ceed
15799df49d7eSJed Brown     ///     .basis_tensor_H1_Lagrange(1, 1, p_fine, q, QuadMode::Gauss)
15809df49d7eSJed Brown     ///     .unwrap();
15819df49d7eSJed Brown     ///
15829df49d7eSJed Brown     /// // Build quadrature data
15839df49d7eSJed Brown     /// let qf_build = ceed.q_function_interior_by_name("Mass1DBuild").unwrap();
15849df49d7eSJed Brown     /// ceed.operator(&qf_build, QFunctionOpt::None, QFunctionOpt::None)
15859df49d7eSJed Brown     ///     .unwrap()
15869df49d7eSJed Brown     ///     .field("dx", &rx, &bx, VectorOpt::Active)
15879df49d7eSJed Brown     ///     .unwrap()
15889df49d7eSJed Brown     ///     .field("weights", ElemRestrictionOpt::None, &bx, VectorOpt::None)
15899df49d7eSJed Brown     ///     .unwrap()
15909df49d7eSJed Brown     ///     .field("qdata", &rq, BasisOpt::Collocated, VectorOpt::Active)
15919df49d7eSJed Brown     ///     .unwrap()
15929df49d7eSJed Brown     ///     .apply(&x, &mut qdata)
15939df49d7eSJed Brown     ///     .unwrap();
15949df49d7eSJed Brown     ///
15959df49d7eSJed Brown     /// // Mass operator
15969df49d7eSJed Brown     /// let qf_mass = ceed.q_function_interior_by_name("MassApply").unwrap();
15979df49d7eSJed Brown     /// let op_mass_fine = ceed
15989df49d7eSJed Brown     ///     .operator(&qf_mass, QFunctionOpt::None, QFunctionOpt::None)
15999df49d7eSJed Brown     ///     .unwrap()
16009df49d7eSJed Brown     ///     .field("u", &ru_fine, &bu_fine, VectorOpt::Active)
16019df49d7eSJed Brown     ///     .unwrap()
16029df49d7eSJed Brown     ///     .field("qdata", &rq, BasisOpt::Collocated, &qdata)
16039df49d7eSJed Brown     ///     .unwrap()
16049df49d7eSJed Brown     ///     .field("v", &ru_fine, &bu_fine, VectorOpt::Active)
16059df49d7eSJed Brown     ///     .unwrap();
16069df49d7eSJed Brown     ///
16079df49d7eSJed Brown     /// // Multigrid setup
1608*80a9ef05SNatalie Beams     /// let mut interp_c_to_f: Vec<Scalar> = vec![0.; p_coarse * p_fine];
16099df49d7eSJed Brown     /// {
16109df49d7eSJed Brown     ///     let mut coarse = ceed.vector(p_coarse).unwrap();
16119df49d7eSJed Brown     ///     let mut fine = ceed.vector(p_fine).unwrap();
16129df49d7eSJed Brown     ///     let basis_c_to_f = ceed
16139df49d7eSJed Brown     ///         .basis_tensor_H1_Lagrange(1, 1, p_coarse, p_fine, QuadMode::GaussLobatto)
16149df49d7eSJed Brown     ///         .unwrap();
16159df49d7eSJed Brown     ///     for i in 0..p_coarse {
16169df49d7eSJed Brown     ///         coarse.set_value(0.0);
16179df49d7eSJed Brown     ///         {
16189df49d7eSJed Brown     ///             let mut array = coarse.view_mut();
16199df49d7eSJed Brown     ///             array[i] = 1.;
16209df49d7eSJed Brown     ///         }
16219df49d7eSJed Brown     ///         basis_c_to_f
16229df49d7eSJed Brown     ///             .apply(
16239df49d7eSJed Brown     ///                 1,
16249df49d7eSJed Brown     ///                 TransposeMode::NoTranspose,
16259df49d7eSJed Brown     ///                 EvalMode::Interp,
16269df49d7eSJed Brown     ///                 &coarse,
16279df49d7eSJed Brown     ///                 &mut fine,
16289df49d7eSJed Brown     ///             )
16299df49d7eSJed Brown     ///             .unwrap();
16309df49d7eSJed Brown     ///         let array = fine.view();
16319df49d7eSJed Brown     ///         for j in 0..p_fine {
16329df49d7eSJed Brown     ///             interp_c_to_f[j * p_coarse + i] = array[j];
16339df49d7eSJed Brown     ///         }
16349df49d7eSJed Brown     ///     }
16359df49d7eSJed Brown     /// }
16369df49d7eSJed Brown     /// let (op_mass_coarse, op_prolong, op_restrict) = op_mass_fine
16379df49d7eSJed Brown     ///     .create_multigrid_level_H1(&multiplicity, &ru_coarse, &bu_coarse, &interp_c_to_f)
16389df49d7eSJed Brown     ///     .unwrap();
16399df49d7eSJed Brown     ///
16409df49d7eSJed Brown     /// // Coarse problem
16419df49d7eSJed Brown     /// u_coarse.set_value(1.0);
16429df49d7eSJed Brown     /// op_mass_coarse.apply(&u_coarse, &mut v_coarse).unwrap();
16439df49d7eSJed Brown     ///
16449df49d7eSJed Brown     /// // Check
1645*80a9ef05SNatalie Beams     /// let sum: Scalar = v_coarse.view().iter().sum();
16469df49d7eSJed Brown     /// assert!(
1647*80a9ef05SNatalie Beams     ///     (sum - 2.0).abs() < 10.0 * libceed::EPSILON,
16489df49d7eSJed Brown     ///     "Incorrect interval length computed"
16499df49d7eSJed Brown     /// );
16509df49d7eSJed Brown     ///
16519df49d7eSJed Brown     /// // Prolong
16529df49d7eSJed Brown     /// op_prolong.apply(&u_coarse, &mut u_fine).unwrap();
16539df49d7eSJed Brown     ///
16549df49d7eSJed Brown     /// // Fine problem
16559df49d7eSJed Brown     /// op_mass_fine.apply(&u_fine, &mut v_fine).unwrap();
16569df49d7eSJed Brown     ///
16579df49d7eSJed Brown     /// // Check
1658*80a9ef05SNatalie Beams     /// let sum: Scalar = v_fine.view().iter().sum();
16599df49d7eSJed Brown     /// assert!(
1660*80a9ef05SNatalie Beams     ///     (sum - 2.0).abs() < 10.0 * libceed::EPSILON,
16619df49d7eSJed Brown     ///     "Incorrect interval length computed"
16629df49d7eSJed Brown     /// );
16639df49d7eSJed Brown     ///
16649df49d7eSJed Brown     /// // Restrict
16659df49d7eSJed Brown     /// op_restrict.apply(&v_fine, &mut v_coarse).unwrap();
16669df49d7eSJed Brown     ///
16679df49d7eSJed Brown     /// // Check
1668*80a9ef05SNatalie Beams     /// let sum: Scalar = v_coarse.view().iter().sum();
16699df49d7eSJed Brown     /// assert!(
1670*80a9ef05SNatalie Beams     ///     (sum - 2.0).abs() < 10.0 * libceed::EPSILON,
16719df49d7eSJed Brown     ///     "Incorrect interval length computed"
16729df49d7eSJed Brown     /// );
16739df49d7eSJed Brown     /// ```
16749df49d7eSJed Brown     pub fn create_multigrid_level_H1(
16759df49d7eSJed Brown         &self,
16769df49d7eSJed Brown         p_mult_fine: &Vector,
16779df49d7eSJed Brown         rstr_coarse: &ElemRestriction,
16789df49d7eSJed Brown         basis_coarse: &Basis,
1679*80a9ef05SNatalie Beams         interpCtoF: &[Scalar],
16809df49d7eSJed Brown     ) -> crate::Result<(Operator, Operator, Operator)> {
16819df49d7eSJed Brown         let mut ptr_coarse = std::ptr::null_mut();
16829df49d7eSJed Brown         let mut ptr_prolong = std::ptr::null_mut();
16839df49d7eSJed Brown         let mut ptr_restrict = std::ptr::null_mut();
16849df49d7eSJed Brown         let ierr = unsafe {
16859df49d7eSJed Brown             bind_ceed::CeedOperatorMultigridLevelCreateH1(
16869df49d7eSJed Brown                 self.op_core.ptr,
16879df49d7eSJed Brown                 p_mult_fine.ptr,
16889df49d7eSJed Brown                 rstr_coarse.ptr,
16899df49d7eSJed Brown                 basis_coarse.ptr,
16909df49d7eSJed Brown                 interpCtoF.as_ptr(),
16919df49d7eSJed Brown                 &mut ptr_coarse,
16929df49d7eSJed Brown                 &mut ptr_prolong,
16939df49d7eSJed Brown                 &mut ptr_restrict,
16949df49d7eSJed Brown             )
16959df49d7eSJed Brown         };
16969df49d7eSJed Brown         self.op_core.ceed.check_error(ierr)?;
16979df49d7eSJed Brown         let op_coarse = Operator::from_raw(self.op_core.ceed, ptr_coarse)?;
16989df49d7eSJed Brown         let op_prolong = Operator::from_raw(self.op_core.ceed, ptr_prolong)?;
16999df49d7eSJed Brown         let op_restrict = Operator::from_raw(self.op_core.ceed, ptr_restrict)?;
17009df49d7eSJed Brown         Ok((op_coarse, op_prolong, op_restrict))
17019df49d7eSJed Brown     }
17029df49d7eSJed Brown }
17039df49d7eSJed Brown 
17049df49d7eSJed Brown // -----------------------------------------------------------------------------
17059df49d7eSJed Brown // Composite Operator
17069df49d7eSJed Brown // -----------------------------------------------------------------------------
17079df49d7eSJed Brown impl<'a> CompositeOperator<'a> {
17089df49d7eSJed Brown     // Constructor
17099df49d7eSJed Brown     pub fn create(ceed: &'a crate::Ceed) -> crate::Result<Self> {
17109df49d7eSJed Brown         let mut ptr = std::ptr::null_mut();
17119df49d7eSJed Brown         let ierr = unsafe { bind_ceed::CeedCompositeOperatorCreate(ceed.ptr, &mut ptr) };
17129df49d7eSJed Brown         ceed.check_error(ierr)?;
17139df49d7eSJed Brown         Ok(Self {
17149df49d7eSJed Brown             op_core: OperatorCore { ceed, ptr },
17159df49d7eSJed Brown         })
17169df49d7eSJed Brown     }
17179df49d7eSJed Brown 
17189df49d7eSJed Brown     /// Apply Operator to a vector
17199df49d7eSJed Brown     ///
17209df49d7eSJed Brown     /// * `input`  - Input Vector
17219df49d7eSJed Brown     /// * `output` - Output Vector
17229df49d7eSJed Brown     ///
17239df49d7eSJed Brown     /// ```
17249df49d7eSJed Brown     /// # use libceed::prelude::*;
17259df49d7eSJed Brown     /// # let ceed = libceed::Ceed::default_init();
17269df49d7eSJed Brown     /// let ne = 4;
17279df49d7eSJed Brown     /// let p = 3;
17289df49d7eSJed Brown     /// let q = 4;
17299df49d7eSJed Brown     /// let ndofs = p * ne - ne + 1;
17309df49d7eSJed Brown     ///
17319df49d7eSJed Brown     /// // Vectors
17329df49d7eSJed Brown     /// let x = ceed.vector_from_slice(&[-1., -0.5, 0.0, 0.5, 1.0]).unwrap();
17339df49d7eSJed Brown     /// let mut qdata_mass = ceed.vector(ne * q).unwrap();
17349df49d7eSJed Brown     /// qdata_mass.set_value(0.0);
17359df49d7eSJed Brown     /// let mut qdata_diff = ceed.vector(ne * q).unwrap();
17369df49d7eSJed Brown     /// qdata_diff.set_value(0.0);
17379df49d7eSJed Brown     /// let mut u = ceed.vector(ndofs).unwrap();
17389df49d7eSJed Brown     /// u.set_value(1.0);
17399df49d7eSJed Brown     /// let mut v = ceed.vector(ndofs).unwrap();
17409df49d7eSJed Brown     /// v.set_value(0.0);
17419df49d7eSJed Brown     ///
17429df49d7eSJed Brown     /// // Restrictions
17439df49d7eSJed Brown     /// let mut indx: Vec<i32> = vec![0; 2 * ne];
17449df49d7eSJed Brown     /// for i in 0..ne {
17459df49d7eSJed Brown     ///     indx[2 * i + 0] = i as i32;
17469df49d7eSJed Brown     ///     indx[2 * i + 1] = (i + 1) as i32;
17479df49d7eSJed Brown     /// }
17489df49d7eSJed Brown     /// let rx = ceed
17499df49d7eSJed Brown     ///     .elem_restriction(ne, 2, 1, 1, ne + 1, MemType::Host, &indx)
17509df49d7eSJed Brown     ///     .unwrap();
17519df49d7eSJed Brown     /// let mut indu: Vec<i32> = vec![0; p * ne];
17529df49d7eSJed Brown     /// for i in 0..ne {
17539df49d7eSJed Brown     ///     indu[p * i + 0] = i as i32;
17549df49d7eSJed Brown     ///     indu[p * i + 1] = (i + 1) as i32;
17559df49d7eSJed Brown     ///     indu[p * i + 2] = (i + 2) as i32;
17569df49d7eSJed Brown     /// }
17579df49d7eSJed Brown     /// let ru = ceed
17589df49d7eSJed Brown     ///     .elem_restriction(ne, 3, 1, 1, ndofs, MemType::Host, &indu)
17599df49d7eSJed Brown     ///     .unwrap();
17609df49d7eSJed Brown     /// let strides: [i32; 3] = [1, q as i32, q as i32];
17619df49d7eSJed Brown     /// let rq = ceed
17629df49d7eSJed Brown     ///     .strided_elem_restriction(ne, q, 1, q * ne, strides)
17639df49d7eSJed Brown     ///     .unwrap();
17649df49d7eSJed Brown     ///
17659df49d7eSJed Brown     /// // Bases
17669df49d7eSJed Brown     /// let bx = ceed
17679df49d7eSJed Brown     ///     .basis_tensor_H1_Lagrange(1, 1, 2, q, QuadMode::Gauss)
17689df49d7eSJed Brown     ///     .unwrap();
17699df49d7eSJed Brown     /// let bu = ceed
17709df49d7eSJed Brown     ///     .basis_tensor_H1_Lagrange(1, 1, p, q, QuadMode::Gauss)
17719df49d7eSJed Brown     ///     .unwrap();
17729df49d7eSJed Brown     ///
17739df49d7eSJed Brown     /// // Build quadrature data
17749df49d7eSJed Brown     /// let qf_build_mass = ceed.q_function_interior_by_name("Mass1DBuild").unwrap();
17759df49d7eSJed Brown     /// ceed.operator(&qf_build_mass, QFunctionOpt::None, QFunctionOpt::None)
17769df49d7eSJed Brown     ///     .unwrap()
17779df49d7eSJed Brown     ///     .field("dx", &rx, &bx, VectorOpt::Active)
17789df49d7eSJed Brown     ///     .unwrap()
17799df49d7eSJed Brown     ///     .field("weights", ElemRestrictionOpt::None, &bx, VectorOpt::None)
17809df49d7eSJed Brown     ///     .unwrap()
17819df49d7eSJed Brown     ///     .field("qdata", &rq, BasisOpt::Collocated, VectorOpt::Active)
17829df49d7eSJed Brown     ///     .unwrap()
17839df49d7eSJed Brown     ///     .apply(&x, &mut qdata_mass)
17849df49d7eSJed Brown     ///     .unwrap();
17859df49d7eSJed Brown     ///
17869df49d7eSJed Brown     /// let qf_build_diff = ceed.q_function_interior_by_name("Poisson1DBuild").unwrap();
17879df49d7eSJed Brown     /// ceed.operator(&qf_build_diff, QFunctionOpt::None, QFunctionOpt::None)
17889df49d7eSJed Brown     ///     .unwrap()
17899df49d7eSJed Brown     ///     .field("dx", &rx, &bx, VectorOpt::Active)
17909df49d7eSJed Brown     ///     .unwrap()
17919df49d7eSJed Brown     ///     .field("weights", ElemRestrictionOpt::None, &bx, VectorOpt::None)
17929df49d7eSJed Brown     ///     .unwrap()
17939df49d7eSJed Brown     ///     .field("qdata", &rq, BasisOpt::Collocated, VectorOpt::Active)
17949df49d7eSJed Brown     ///     .unwrap()
17959df49d7eSJed Brown     ///     .apply(&x, &mut qdata_diff)
17969df49d7eSJed Brown     ///     .unwrap();
17979df49d7eSJed Brown     ///
17989df49d7eSJed Brown     /// // Application operator
17999df49d7eSJed Brown     /// let qf_mass = ceed.q_function_interior_by_name("MassApply").unwrap();
18009df49d7eSJed Brown     /// let op_mass = ceed
18019df49d7eSJed Brown     ///     .operator(&qf_mass, QFunctionOpt::None, QFunctionOpt::None)
18029df49d7eSJed Brown     ///     .unwrap()
18039df49d7eSJed Brown     ///     .field("u", &ru, &bu, VectorOpt::Active)
18049df49d7eSJed Brown     ///     .unwrap()
18059df49d7eSJed Brown     ///     .field("qdata", &rq, BasisOpt::Collocated, &qdata_mass)
18069df49d7eSJed Brown     ///     .unwrap()
18079df49d7eSJed Brown     ///     .field("v", &ru, &bu, VectorOpt::Active)
18089df49d7eSJed Brown     ///     .unwrap();
18099df49d7eSJed Brown     ///
18109df49d7eSJed Brown     /// let qf_diff = ceed.q_function_interior_by_name("Poisson1DApply").unwrap();
18119df49d7eSJed Brown     /// let op_diff = ceed
18129df49d7eSJed Brown     ///     .operator(&qf_diff, QFunctionOpt::None, QFunctionOpt::None)
18139df49d7eSJed Brown     ///     .unwrap()
18149df49d7eSJed Brown     ///     .field("du", &ru, &bu, VectorOpt::Active)
18159df49d7eSJed Brown     ///     .unwrap()
18169df49d7eSJed Brown     ///     .field("qdata", &rq, BasisOpt::Collocated, &qdata_diff)
18179df49d7eSJed Brown     ///     .unwrap()
18189df49d7eSJed Brown     ///     .field("dv", &ru, &bu, VectorOpt::Active)
18199df49d7eSJed Brown     ///     .unwrap();
18209df49d7eSJed Brown     ///
18219df49d7eSJed Brown     /// let op_composite = ceed
18229df49d7eSJed Brown     ///     .composite_operator()
18239df49d7eSJed Brown     ///     .unwrap()
18249df49d7eSJed Brown     ///     .sub_operator(&op_mass)
18259df49d7eSJed Brown     ///     .unwrap()
18269df49d7eSJed Brown     ///     .sub_operator(&op_diff)
18279df49d7eSJed Brown     ///     .unwrap();
18289df49d7eSJed Brown     ///
18299df49d7eSJed Brown     /// v.set_value(0.0);
18309df49d7eSJed Brown     /// op_composite.apply(&u, &mut v).unwrap();
18319df49d7eSJed Brown     ///
18329df49d7eSJed Brown     /// // Check
1833*80a9ef05SNatalie Beams     /// let sum: Scalar = v.view().iter().sum();
18349df49d7eSJed Brown     /// assert!(
1835*80a9ef05SNatalie Beams     ///     (sum - 2.0).abs() < 10.0 * libceed::EPSILON,
18369df49d7eSJed Brown     ///     "Incorrect interval length computed"
18379df49d7eSJed Brown     /// );
18389df49d7eSJed Brown     /// ```
18399df49d7eSJed Brown     pub fn apply(&self, input: &Vector, output: &mut Vector) -> crate::Result<i32> {
18409df49d7eSJed Brown         self.op_core.apply(input, output)
18419df49d7eSJed Brown     }
18429df49d7eSJed Brown 
18439df49d7eSJed Brown     /// Apply Operator to a vector and add result to output vector
18449df49d7eSJed Brown     ///
18459df49d7eSJed Brown     /// * `input`  - Input Vector
18469df49d7eSJed Brown     /// * `output` - Output Vector
18479df49d7eSJed Brown     ///
18489df49d7eSJed Brown     /// ```
18499df49d7eSJed Brown     /// # use libceed::prelude::*;
18509df49d7eSJed Brown     /// # let ceed = libceed::Ceed::default_init();
18519df49d7eSJed Brown     /// let ne = 4;
18529df49d7eSJed Brown     /// let p = 3;
18539df49d7eSJed Brown     /// let q = 4;
18549df49d7eSJed Brown     /// let ndofs = p * ne - ne + 1;
18559df49d7eSJed Brown     ///
18569df49d7eSJed Brown     /// // Vectors
18579df49d7eSJed Brown     /// let x = ceed.vector_from_slice(&[-1., -0.5, 0.0, 0.5, 1.0]).unwrap();
18589df49d7eSJed Brown     /// let mut qdata_mass = ceed.vector(ne * q).unwrap();
18599df49d7eSJed Brown     /// qdata_mass.set_value(0.0);
18609df49d7eSJed Brown     /// let mut qdata_diff = ceed.vector(ne * q).unwrap();
18619df49d7eSJed Brown     /// qdata_diff.set_value(0.0);
18629df49d7eSJed Brown     /// let mut u = ceed.vector(ndofs).unwrap();
18639df49d7eSJed Brown     /// u.set_value(1.0);
18649df49d7eSJed Brown     /// let mut v = ceed.vector(ndofs).unwrap();
18659df49d7eSJed Brown     /// v.set_value(0.0);
18669df49d7eSJed Brown     ///
18679df49d7eSJed Brown     /// // Restrictions
18689df49d7eSJed Brown     /// let mut indx: Vec<i32> = vec![0; 2 * ne];
18699df49d7eSJed Brown     /// for i in 0..ne {
18709df49d7eSJed Brown     ///     indx[2 * i + 0] = i as i32;
18719df49d7eSJed Brown     ///     indx[2 * i + 1] = (i + 1) as i32;
18729df49d7eSJed Brown     /// }
18739df49d7eSJed Brown     /// let rx = ceed
18749df49d7eSJed Brown     ///     .elem_restriction(ne, 2, 1, 1, ne + 1, MemType::Host, &indx)
18759df49d7eSJed Brown     ///     .unwrap();
18769df49d7eSJed Brown     /// let mut indu: Vec<i32> = vec![0; p * ne];
18779df49d7eSJed Brown     /// for i in 0..ne {
18789df49d7eSJed Brown     ///     indu[p * i + 0] = i as i32;
18799df49d7eSJed Brown     ///     indu[p * i + 1] = (i + 1) as i32;
18809df49d7eSJed Brown     ///     indu[p * i + 2] = (i + 2) as i32;
18819df49d7eSJed Brown     /// }
18829df49d7eSJed Brown     /// let ru = ceed
18839df49d7eSJed Brown     ///     .elem_restriction(ne, 3, 1, 1, ndofs, MemType::Host, &indu)
18849df49d7eSJed Brown     ///     .unwrap();
18859df49d7eSJed Brown     /// let strides: [i32; 3] = [1, q as i32, q as i32];
18869df49d7eSJed Brown     /// let rq = ceed
18879df49d7eSJed Brown     ///     .strided_elem_restriction(ne, q, 1, q * ne, strides)
18889df49d7eSJed Brown     ///     .unwrap();
18899df49d7eSJed Brown     ///
18909df49d7eSJed Brown     /// // Bases
18919df49d7eSJed Brown     /// let bx = ceed
18929df49d7eSJed Brown     ///     .basis_tensor_H1_Lagrange(1, 1, 2, q, QuadMode::Gauss)
18939df49d7eSJed Brown     ///     .unwrap();
18949df49d7eSJed Brown     /// let bu = ceed
18959df49d7eSJed Brown     ///     .basis_tensor_H1_Lagrange(1, 1, p, q, QuadMode::Gauss)
18969df49d7eSJed Brown     ///     .unwrap();
18979df49d7eSJed Brown     ///
18989df49d7eSJed Brown     /// // Build quadrature data
18999df49d7eSJed Brown     /// let qf_build_mass = ceed.q_function_interior_by_name("Mass1DBuild").unwrap();
19009df49d7eSJed Brown     /// ceed.operator(&qf_build_mass, QFunctionOpt::None, QFunctionOpt::None)
19019df49d7eSJed Brown     ///     .unwrap()
19029df49d7eSJed Brown     ///     .field("dx", &rx, &bx, VectorOpt::Active)
19039df49d7eSJed Brown     ///     .unwrap()
19049df49d7eSJed Brown     ///     .field("weights", ElemRestrictionOpt::None, &bx, VectorOpt::None)
19059df49d7eSJed Brown     ///     .unwrap()
19069df49d7eSJed Brown     ///     .field("qdata", &rq, BasisOpt::Collocated, VectorOpt::Active)
19079df49d7eSJed Brown     ///     .unwrap()
19089df49d7eSJed Brown     ///     .apply(&x, &mut qdata_mass)
19099df49d7eSJed Brown     ///     .unwrap();
19109df49d7eSJed Brown     ///
19119df49d7eSJed Brown     /// let qf_build_diff = ceed.q_function_interior_by_name("Poisson1DBuild").unwrap();
19129df49d7eSJed Brown     /// ceed.operator(&qf_build_diff, QFunctionOpt::None, QFunctionOpt::None)
19139df49d7eSJed Brown     ///     .unwrap()
19149df49d7eSJed Brown     ///     .field("dx", &rx, &bx, VectorOpt::Active)
19159df49d7eSJed Brown     ///     .unwrap()
19169df49d7eSJed Brown     ///     .field("weights", ElemRestrictionOpt::None, &bx, VectorOpt::None)
19179df49d7eSJed Brown     ///     .unwrap()
19189df49d7eSJed Brown     ///     .field("qdata", &rq, BasisOpt::Collocated, VectorOpt::Active)
19199df49d7eSJed Brown     ///     .unwrap()
19209df49d7eSJed Brown     ///     .apply(&x, &mut qdata_diff)
19219df49d7eSJed Brown     ///     .unwrap();
19229df49d7eSJed Brown     ///
19239df49d7eSJed Brown     /// // Application operator
19249df49d7eSJed Brown     /// let qf_mass = ceed.q_function_interior_by_name("MassApply").unwrap();
19259df49d7eSJed Brown     /// let op_mass = ceed
19269df49d7eSJed Brown     ///     .operator(&qf_mass, QFunctionOpt::None, QFunctionOpt::None)
19279df49d7eSJed Brown     ///     .unwrap()
19289df49d7eSJed Brown     ///     .field("u", &ru, &bu, VectorOpt::Active)
19299df49d7eSJed Brown     ///     .unwrap()
19309df49d7eSJed Brown     ///     .field("qdata", &rq, BasisOpt::Collocated, &qdata_mass)
19319df49d7eSJed Brown     ///     .unwrap()
19329df49d7eSJed Brown     ///     .field("v", &ru, &bu, VectorOpt::Active)
19339df49d7eSJed Brown     ///     .unwrap();
19349df49d7eSJed Brown     ///
19359df49d7eSJed Brown     /// let qf_diff = ceed.q_function_interior_by_name("Poisson1DApply").unwrap();
19369df49d7eSJed Brown     /// let op_diff = ceed
19379df49d7eSJed Brown     ///     .operator(&qf_diff, QFunctionOpt::None, QFunctionOpt::None)
19389df49d7eSJed Brown     ///     .unwrap()
19399df49d7eSJed Brown     ///     .field("du", &ru, &bu, VectorOpt::Active)
19409df49d7eSJed Brown     ///     .unwrap()
19419df49d7eSJed Brown     ///     .field("qdata", &rq, BasisOpt::Collocated, &qdata_diff)
19429df49d7eSJed Brown     ///     .unwrap()
19439df49d7eSJed Brown     ///     .field("dv", &ru, &bu, VectorOpt::Active)
19449df49d7eSJed Brown     ///     .unwrap();
19459df49d7eSJed Brown     ///
19469df49d7eSJed Brown     /// let op_composite = ceed
19479df49d7eSJed Brown     ///     .composite_operator()
19489df49d7eSJed Brown     ///     .unwrap()
19499df49d7eSJed Brown     ///     .sub_operator(&op_mass)
19509df49d7eSJed Brown     ///     .unwrap()
19519df49d7eSJed Brown     ///     .sub_operator(&op_diff)
19529df49d7eSJed Brown     ///     .unwrap();
19539df49d7eSJed Brown     ///
19549df49d7eSJed Brown     /// v.set_value(1.0);
19559df49d7eSJed Brown     /// op_composite.apply_add(&u, &mut v).unwrap();
19569df49d7eSJed Brown     ///
19579df49d7eSJed Brown     /// // Check
1958*80a9ef05SNatalie Beams     /// let sum: Scalar = v.view().iter().sum();
19599df49d7eSJed Brown     /// assert!(
1960*80a9ef05SNatalie Beams     ///     (sum - (2.0 + ndofs as Scalar)).abs() < 10.0 * libceed::EPSILON,
19619df49d7eSJed Brown     ///     "Incorrect interval length computed"
19629df49d7eSJed Brown     /// );
19639df49d7eSJed Brown     /// ```
19649df49d7eSJed Brown     pub fn apply_add(&self, input: &Vector, output: &mut Vector) -> crate::Result<i32> {
19659df49d7eSJed Brown         self.op_core.apply_add(input, output)
19669df49d7eSJed Brown     }
19679df49d7eSJed Brown 
19689df49d7eSJed Brown     /// Add a sub-Operator to a Composite Operator
19699df49d7eSJed Brown     ///
19709df49d7eSJed Brown     /// * `subop` - Sub-Operator
19719df49d7eSJed Brown     ///
19729df49d7eSJed Brown     /// ```
19739df49d7eSJed Brown     /// # use libceed::prelude::*;
19749df49d7eSJed Brown     /// # let ceed = libceed::Ceed::default_init();
19759df49d7eSJed Brown     /// let mut op = ceed.composite_operator().unwrap();
19769df49d7eSJed Brown     ///
19779df49d7eSJed Brown     /// let qf_mass = ceed.q_function_interior_by_name("MassApply").unwrap();
19789df49d7eSJed Brown     /// let op_mass = ceed
19799df49d7eSJed Brown     ///     .operator(&qf_mass, QFunctionOpt::None, QFunctionOpt::None)
19809df49d7eSJed Brown     ///     .unwrap();
19819df49d7eSJed Brown     /// op = op.sub_operator(&op_mass).unwrap();
19829df49d7eSJed Brown     ///
19839df49d7eSJed Brown     /// let qf_diff = ceed.q_function_interior_by_name("Poisson1DApply").unwrap();
19849df49d7eSJed Brown     /// let op_diff = ceed
19859df49d7eSJed Brown     ///     .operator(&qf_diff, QFunctionOpt::None, QFunctionOpt::None)
19869df49d7eSJed Brown     ///     .unwrap();
19879df49d7eSJed Brown     /// op = op.sub_operator(&op_diff).unwrap();
19889df49d7eSJed Brown     /// ```
19899df49d7eSJed Brown     #[allow(unused_mut)]
19909df49d7eSJed Brown     pub fn sub_operator(mut self, subop: &Operator) -> crate::Result<Self> {
19919df49d7eSJed Brown         let ierr =
19929df49d7eSJed Brown             unsafe { bind_ceed::CeedCompositeOperatorAddSub(self.op_core.ptr, subop.op_core.ptr) };
19939df49d7eSJed Brown         self.op_core.ceed.check_error(ierr)?;
19949df49d7eSJed Brown         Ok(self)
19959df49d7eSJed Brown     }
19969df49d7eSJed Brown 
19979df49d7eSJed Brown     /// Assemble the diagonal of a square linear Operator
19989df49d7eSJed Brown     ///
19999df49d7eSJed Brown     /// This overwrites a Vector with the diagonal of a linear Operator.
20009df49d7eSJed Brown     ///
20019df49d7eSJed Brown     /// Note: Currently only non-composite Operators with a single field and
20029df49d7eSJed Brown     ///      composite Operators with single field sub-operators are supported.
20039df49d7eSJed Brown     ///
20049df49d7eSJed Brown     /// * `op`        - Operator to assemble QFunction
20059df49d7eSJed Brown     /// * `assembled` - Vector to store assembled Operator diagonal
20069df49d7eSJed Brown     pub fn linear_asssemble_diagonal(&self, assembled: &mut Vector) -> crate::Result<i32> {
20079df49d7eSJed Brown         self.op_core.linear_assemble_add_diagonal(assembled)
20089df49d7eSJed Brown     }
20099df49d7eSJed Brown 
20109df49d7eSJed Brown     /// Assemble the point block diagonal of a square linear Operator
20119df49d7eSJed Brown     ///
20129df49d7eSJed Brown     /// This overwrites a Vector with the point block diagonal of a linear
20139df49d7eSJed Brown     ///   Operator.
20149df49d7eSJed Brown     ///
20159df49d7eSJed Brown     /// Note: Currently only non-composite Operators with a single field and
20169df49d7eSJed Brown     ///      composite Operators with single field sub-operators are supported.
20179df49d7eSJed Brown     ///
20189df49d7eSJed Brown     /// * `op`        - Operator to assemble QFunction
20199df49d7eSJed Brown     /// * `assembled` - Vector to store assembled Operator diagonal
20209df49d7eSJed Brown     pub fn linear_assemble_add_diagonal(&self, assembled: &mut Vector) -> crate::Result<i32> {
20219df49d7eSJed Brown         self.op_core.linear_assemble_add_diagonal(assembled)
20229df49d7eSJed Brown     }
20239df49d7eSJed Brown 
20249df49d7eSJed Brown     /// Assemble the diagonal of a square linear Operator
20259df49d7eSJed Brown     ///
20269df49d7eSJed Brown     /// This overwrites a Vector with the diagonal of a linear Operator.
20279df49d7eSJed Brown     ///
20289df49d7eSJed Brown     /// Note: Currently only non-composite Operators with a single field and
20299df49d7eSJed Brown     ///      composite Operators with single field sub-operators are supported.
20309df49d7eSJed Brown     ///
20319df49d7eSJed Brown     /// * `op`        - Operator to assemble QFunction
20329df49d7eSJed Brown     /// * `assembled` - Vector to store assembled CeedOperator point block
20339df49d7eSJed Brown     ///                   diagonal, provided in row-major form with an
20349df49d7eSJed Brown     ///                   `ncomp * ncomp` block at each node. The dimensions of
20359df49d7eSJed Brown     ///                   this vector are derived from the active vector for
20369df49d7eSJed Brown     ///                   the CeedOperator. The array has shape
20379df49d7eSJed Brown     ///                   `[nodes, component out, component in]`.
20389df49d7eSJed Brown     pub fn linear_assemble_point_block_diagonal(
20399df49d7eSJed Brown         &self,
20409df49d7eSJed Brown         assembled: &mut Vector,
20419df49d7eSJed Brown     ) -> crate::Result<i32> {
20429df49d7eSJed Brown         self.op_core.linear_assemble_add_diagonal(assembled)
20439df49d7eSJed Brown     }
20449df49d7eSJed Brown 
20459df49d7eSJed Brown     /// Assemble the diagonal of a square linear Operator
20469df49d7eSJed Brown     ///
20479df49d7eSJed Brown     /// This sums into a Vector with the diagonal of a linear Operator.
20489df49d7eSJed Brown     ///
20499df49d7eSJed Brown     /// Note: Currently only non-composite Operators with a single field and
20509df49d7eSJed Brown     ///      composite Operators with single field sub-operators are supported.
20519df49d7eSJed Brown     ///
20529df49d7eSJed Brown     /// * `op`        - Operator to assemble QFunction
20539df49d7eSJed Brown     /// * `assembled` - Vector to store assembled CeedOperator point block
20549df49d7eSJed Brown     ///                   diagonal, provided in row-major form with an
20559df49d7eSJed Brown     ///                   `ncomp * ncomp` block at each node. The dimensions of
20569df49d7eSJed Brown     ///                   this vector are derived from the active vector for
20579df49d7eSJed Brown     ///                   the CeedOperator. The array has shape
20589df49d7eSJed Brown     ///                   `[nodes, component out, component in]`.
20599df49d7eSJed Brown     pub fn linear_assemble_add_point_block_diagonal(
20609df49d7eSJed Brown         &self,
20619df49d7eSJed Brown         assembled: &mut Vector,
20629df49d7eSJed Brown     ) -> crate::Result<i32> {
20639df49d7eSJed Brown         self.op_core.linear_assemble_add_diagonal(assembled)
20649df49d7eSJed Brown     }
20659df49d7eSJed Brown }
20669df49d7eSJed Brown 
20679df49d7eSJed Brown // -----------------------------------------------------------------------------
2068