xref: /libCEED/rust/libceed/src/operator.rs (revision c68be7a2e45631197b626561fad53d5b146fcd59)
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 // -----------------------------------------------------------------------------
26*c68be7a2SJeremy L Thompson #[derive(Debug)]
279df49d7eSJed Brown pub(crate) struct OperatorCore<'a> {
28*c68be7a2SJeremy L Thompson     pub(crate) ceed: &'a crate::Ceed,
299df49d7eSJed Brown     ptr: bind_ceed::CeedOperator,
309df49d7eSJed Brown }
319df49d7eSJed Brown 
32*c68be7a2SJeremy L Thompson #[derive(Debug)]
339df49d7eSJed Brown pub struct Operator<'a> {
349df49d7eSJed Brown     op_core: OperatorCore<'a>,
359df49d7eSJed Brown }
369df49d7eSJed Brown 
37*c68be7a2SJeremy L Thompson #[derive(Debug)]
389df49d7eSJed Brown pub struct CompositeOperator<'a> {
399df49d7eSJed Brown     op_core: OperatorCore<'a>,
409df49d7eSJed Brown }
419df49d7eSJed Brown 
429df49d7eSJed Brown // -----------------------------------------------------------------------------
439df49d7eSJed Brown // Destructor
449df49d7eSJed Brown // -----------------------------------------------------------------------------
459df49d7eSJed Brown impl<'a> Drop for OperatorCore<'a> {
469df49d7eSJed Brown     fn drop(&mut self) {
479df49d7eSJed Brown         unsafe {
489df49d7eSJed Brown             bind_ceed::CeedOperatorDestroy(&mut self.ptr);
499df49d7eSJed Brown         }
509df49d7eSJed Brown     }
519df49d7eSJed Brown }
529df49d7eSJed Brown 
539df49d7eSJed Brown // -----------------------------------------------------------------------------
549df49d7eSJed Brown // Display
559df49d7eSJed Brown // -----------------------------------------------------------------------------
569df49d7eSJed Brown impl<'a> fmt::Display for OperatorCore<'a> {
579df49d7eSJed Brown     fn fmt(&self, f: &mut fmt::Formatter) -> fmt::Result {
589df49d7eSJed Brown         let mut ptr = std::ptr::null_mut();
599df49d7eSJed Brown         let mut sizeloc = crate::MAX_BUFFER_LENGTH;
609df49d7eSJed Brown         let cstring = unsafe {
619df49d7eSJed Brown             let file = bind_ceed::open_memstream(&mut ptr, &mut sizeloc);
629df49d7eSJed Brown             bind_ceed::CeedOperatorView(self.ptr, file);
639df49d7eSJed Brown             bind_ceed::fclose(file);
649df49d7eSJed Brown             CString::from_raw(ptr)
659df49d7eSJed Brown         };
669df49d7eSJed Brown         cstring.to_string_lossy().fmt(f)
679df49d7eSJed Brown     }
689df49d7eSJed Brown }
699df49d7eSJed Brown 
709df49d7eSJed Brown /// View an Operator
719df49d7eSJed Brown ///
729df49d7eSJed Brown /// ```
739df49d7eSJed Brown /// # use libceed::prelude::*;
74*c68be7a2SJeremy L Thompson /// # fn main() -> Result<(), libceed::CeedError> {
759df49d7eSJed Brown /// # let ceed = libceed::Ceed::default_init();
76*c68be7a2SJeremy L Thompson /// let qf = ceed.q_function_interior_by_name("Mass1DBuild")?;
779df49d7eSJed Brown ///
789df49d7eSJed Brown /// // Operator field arguments
799df49d7eSJed Brown /// let ne = 3;
809df49d7eSJed Brown /// let q = 4 as usize;
819df49d7eSJed Brown /// let mut ind: Vec<i32> = vec![0; 2 * ne];
829df49d7eSJed Brown /// for i in 0..ne {
839df49d7eSJed Brown ///     ind[2 * i + 0] = i as i32;
849df49d7eSJed Brown ///     ind[2 * i + 1] = (i + 1) as i32;
859df49d7eSJed Brown /// }
86*c68be7a2SJeremy L Thompson /// let r = ceed.elem_restriction(ne, 2, 1, 1, ne + 1, MemType::Host, &ind)?;
879df49d7eSJed Brown /// let strides: [i32; 3] = [1, q as i32, q as i32];
88*c68be7a2SJeremy L Thompson /// let rq = ceed.strided_elem_restriction(ne, 2, 1, q * ne, strides)?;
899df49d7eSJed Brown ///
90*c68be7a2SJeremy L Thompson /// let b = ceed.basis_tensor_H1_Lagrange(1, 1, 2, q, QuadMode::Gauss)?;
919df49d7eSJed Brown ///
929df49d7eSJed Brown /// // Operator fields
939df49d7eSJed Brown /// let op = ceed
94*c68be7a2SJeremy L Thompson ///     .operator(&qf, QFunctionOpt::None, QFunctionOpt::None)?
95*c68be7a2SJeremy L Thompson ///     .field("dx", &r, &b, VectorOpt::Active)?
96*c68be7a2SJeremy L Thompson ///     .field("weights", ElemRestrictionOpt::None, &b, VectorOpt::None)?
97*c68be7a2SJeremy L Thompson ///     .field("qdata", &rq, BasisOpt::Collocated, VectorOpt::Active)?;
989df49d7eSJed Brown ///
999df49d7eSJed Brown /// println!("{}", op);
100*c68be7a2SJeremy L Thompson /// # Ok(())
101*c68be7a2SJeremy L Thompson /// # }
1029df49d7eSJed Brown /// ```
1039df49d7eSJed Brown impl<'a> fmt::Display for Operator<'a> {
1049df49d7eSJed Brown     fn fmt(&self, f: &mut fmt::Formatter) -> fmt::Result {
1059df49d7eSJed Brown         self.op_core.fmt(f)
1069df49d7eSJed Brown     }
1079df49d7eSJed Brown }
1089df49d7eSJed Brown 
1099df49d7eSJed Brown /// View a composite Operator
1109df49d7eSJed Brown ///
1119df49d7eSJed Brown /// ```
1129df49d7eSJed Brown /// # use libceed::prelude::*;
113*c68be7a2SJeremy L Thompson /// # fn main() -> Result<(), libceed::CeedError> {
1149df49d7eSJed Brown /// # let ceed = libceed::Ceed::default_init();
1159df49d7eSJed Brown ///
1169df49d7eSJed Brown /// // Sub operator field arguments
1179df49d7eSJed Brown /// let ne = 3;
1189df49d7eSJed Brown /// let q = 4 as usize;
1199df49d7eSJed Brown /// let mut ind: Vec<i32> = vec![0; 2 * ne];
1209df49d7eSJed Brown /// for i in 0..ne {
1219df49d7eSJed Brown ///     ind[2 * i + 0] = i as i32;
1229df49d7eSJed Brown ///     ind[2 * i + 1] = (i + 1) as i32;
1239df49d7eSJed Brown /// }
124*c68be7a2SJeremy L Thompson /// let r = ceed.elem_restriction(ne, 2, 1, 1, ne + 1, MemType::Host, &ind)?;
1259df49d7eSJed Brown /// let strides: [i32; 3] = [1, q as i32, q as i32];
126*c68be7a2SJeremy L Thompson /// let rq = ceed.strided_elem_restriction(ne, 2, 1, q * ne, strides)?;
1279df49d7eSJed Brown ///
128*c68be7a2SJeremy L Thompson /// let b = ceed.basis_tensor_H1_Lagrange(1, 1, 2, q, QuadMode::Gauss)?;
1299df49d7eSJed Brown ///
130*c68be7a2SJeremy L Thompson /// let qdata_mass = ceed.vector(q * ne)?;
131*c68be7a2SJeremy L Thompson /// let qdata_diff = ceed.vector(q * ne)?;
1329df49d7eSJed Brown ///
133*c68be7a2SJeremy L Thompson /// let qf_mass = ceed.q_function_interior_by_name("MassApply")?;
1349df49d7eSJed Brown /// let op_mass = ceed
135*c68be7a2SJeremy L Thompson ///     .operator(&qf_mass, QFunctionOpt::None, QFunctionOpt::None)?
136*c68be7a2SJeremy L Thompson ///     .field("u", &r, &b, VectorOpt::Active)?
137*c68be7a2SJeremy L Thompson ///     .field("qdata", &rq, BasisOpt::Collocated, &qdata_mass)?
138*c68be7a2SJeremy L Thompson ///     .field("v", &r, &b, VectorOpt::Active)?;
1399df49d7eSJed Brown ///
140*c68be7a2SJeremy L Thompson /// let qf_diff = ceed.q_function_interior_by_name("Poisson1DApply")?;
1419df49d7eSJed Brown /// let op_diff = ceed
142*c68be7a2SJeremy L Thompson ///     .operator(&qf_diff, QFunctionOpt::None, QFunctionOpt::None)?
143*c68be7a2SJeremy L Thompson ///     .field("du", &r, &b, VectorOpt::Active)?
144*c68be7a2SJeremy L Thompson ///     .field("qdata", &rq, BasisOpt::Collocated, &qdata_diff)?
145*c68be7a2SJeremy L Thompson ///     .field("dv", &r, &b, VectorOpt::Active)?;
1469df49d7eSJed Brown ///
1479df49d7eSJed Brown /// let op = ceed
148*c68be7a2SJeremy L Thompson ///     .composite_operator()?
149*c68be7a2SJeremy L Thompson ///     .sub_operator(&op_mass)?
150*c68be7a2SJeremy L Thompson ///     .sub_operator(&op_diff)?;
1519df49d7eSJed Brown ///
1529df49d7eSJed Brown /// println!("{}", op);
153*c68be7a2SJeremy L Thompson /// # Ok(())
154*c68be7a2SJeremy L Thompson /// # }
1559df49d7eSJed Brown /// ```
1569df49d7eSJed Brown impl<'a> fmt::Display for CompositeOperator<'a> {
1579df49d7eSJed Brown     fn fmt(&self, f: &mut fmt::Formatter) -> fmt::Result {
1589df49d7eSJed Brown         self.op_core.fmt(f)
1599df49d7eSJed Brown     }
1609df49d7eSJed Brown }
1619df49d7eSJed Brown 
1629df49d7eSJed Brown // -----------------------------------------------------------------------------
1639df49d7eSJed Brown // Core functionality
1649df49d7eSJed Brown // -----------------------------------------------------------------------------
1659df49d7eSJed Brown impl<'a> OperatorCore<'a> {
1669df49d7eSJed Brown     // Common implementations
1679df49d7eSJed Brown     pub fn apply(&self, input: &Vector, output: &mut Vector) -> crate::Result<i32> {
1689df49d7eSJed Brown         let ierr = unsafe {
1699df49d7eSJed Brown             bind_ceed::CeedOperatorApply(
1709df49d7eSJed Brown                 self.ptr,
1719df49d7eSJed Brown                 input.ptr,
1729df49d7eSJed Brown                 output.ptr,
1739df49d7eSJed Brown                 bind_ceed::CEED_REQUEST_IMMEDIATE,
1749df49d7eSJed Brown             )
1759df49d7eSJed Brown         };
1769df49d7eSJed Brown         self.ceed.check_error(ierr)
1779df49d7eSJed Brown     }
1789df49d7eSJed Brown 
1799df49d7eSJed Brown     pub fn apply_add(&self, input: &Vector, output: &mut Vector) -> crate::Result<i32> {
1809df49d7eSJed Brown         let ierr = unsafe {
1819df49d7eSJed Brown             bind_ceed::CeedOperatorApplyAdd(
1829df49d7eSJed Brown                 self.ptr,
1839df49d7eSJed Brown                 input.ptr,
1849df49d7eSJed Brown                 output.ptr,
1859df49d7eSJed Brown                 bind_ceed::CEED_REQUEST_IMMEDIATE,
1869df49d7eSJed Brown             )
1879df49d7eSJed Brown         };
1889df49d7eSJed Brown         self.ceed.check_error(ierr)
1899df49d7eSJed Brown     }
1909df49d7eSJed Brown 
1919df49d7eSJed Brown     pub fn linear_assemble_diagonal(&self, assembled: &mut Vector) -> crate::Result<i32> {
1929df49d7eSJed Brown         let ierr = unsafe {
1939df49d7eSJed Brown             bind_ceed::CeedOperatorLinearAssembleDiagonal(
1949df49d7eSJed Brown                 self.ptr,
1959df49d7eSJed Brown                 assembled.ptr,
1969df49d7eSJed Brown                 bind_ceed::CEED_REQUEST_IMMEDIATE,
1979df49d7eSJed Brown             )
1989df49d7eSJed Brown         };
1999df49d7eSJed Brown         self.ceed.check_error(ierr)
2009df49d7eSJed Brown     }
2019df49d7eSJed Brown 
2029df49d7eSJed Brown     pub fn linear_assemble_add_diagonal(&self, assembled: &mut Vector) -> crate::Result<i32> {
2039df49d7eSJed Brown         let ierr = unsafe {
2049df49d7eSJed Brown             bind_ceed::CeedOperatorLinearAssembleAddDiagonal(
2059df49d7eSJed Brown                 self.ptr,
2069df49d7eSJed Brown                 assembled.ptr,
2079df49d7eSJed Brown                 bind_ceed::CEED_REQUEST_IMMEDIATE,
2089df49d7eSJed Brown             )
2099df49d7eSJed Brown         };
2109df49d7eSJed Brown         self.ceed.check_error(ierr)
2119df49d7eSJed Brown     }
2129df49d7eSJed Brown 
2139df49d7eSJed Brown     pub fn linear_assemble_point_block_diagonal(
2149df49d7eSJed Brown         &self,
2159df49d7eSJed Brown         assembled: &mut Vector,
2169df49d7eSJed Brown     ) -> crate::Result<i32> {
2179df49d7eSJed Brown         let ierr = unsafe {
2189df49d7eSJed Brown             bind_ceed::CeedOperatorLinearAssemblePointBlockDiagonal(
2199df49d7eSJed Brown                 self.ptr,
2209df49d7eSJed Brown                 assembled.ptr,
2219df49d7eSJed Brown                 bind_ceed::CEED_REQUEST_IMMEDIATE,
2229df49d7eSJed Brown             )
2239df49d7eSJed Brown         };
2249df49d7eSJed Brown         self.ceed.check_error(ierr)
2259df49d7eSJed Brown     }
2269df49d7eSJed Brown 
2279df49d7eSJed Brown     pub fn linear_assemble_add_point_block_diagonal(
2289df49d7eSJed Brown         &self,
2299df49d7eSJed Brown         assembled: &mut Vector,
2309df49d7eSJed Brown     ) -> crate::Result<i32> {
2319df49d7eSJed Brown         let ierr = unsafe {
2329df49d7eSJed Brown             bind_ceed::CeedOperatorLinearAssembleAddPointBlockDiagonal(
2339df49d7eSJed Brown                 self.ptr,
2349df49d7eSJed Brown                 assembled.ptr,
2359df49d7eSJed Brown                 bind_ceed::CEED_REQUEST_IMMEDIATE,
2369df49d7eSJed Brown             )
2379df49d7eSJed Brown         };
2389df49d7eSJed Brown         self.ceed.check_error(ierr)
2399df49d7eSJed Brown     }
2409df49d7eSJed Brown }
2419df49d7eSJed Brown 
2429df49d7eSJed Brown // -----------------------------------------------------------------------------
2439df49d7eSJed Brown // Operator
2449df49d7eSJed Brown // -----------------------------------------------------------------------------
2459df49d7eSJed Brown impl<'a> Operator<'a> {
2469df49d7eSJed Brown     // Constructor
2479df49d7eSJed Brown     pub fn create<'b>(
2489df49d7eSJed Brown         ceed: &'a crate::Ceed,
2499df49d7eSJed Brown         qf: impl Into<QFunctionOpt<'b>>,
2509df49d7eSJed Brown         dqf: impl Into<QFunctionOpt<'b>>,
2519df49d7eSJed Brown         dqfT: impl Into<QFunctionOpt<'b>>,
2529df49d7eSJed Brown     ) -> crate::Result<Self> {
2539df49d7eSJed Brown         let mut ptr = std::ptr::null_mut();
2549df49d7eSJed Brown         let ierr = unsafe {
2559df49d7eSJed Brown             bind_ceed::CeedOperatorCreate(
2569df49d7eSJed Brown                 ceed.ptr,
2579df49d7eSJed Brown                 qf.into().to_raw(),
2589df49d7eSJed Brown                 dqf.into().to_raw(),
2599df49d7eSJed Brown                 dqfT.into().to_raw(),
2609df49d7eSJed Brown                 &mut ptr,
2619df49d7eSJed Brown             )
2629df49d7eSJed Brown         };
2639df49d7eSJed Brown         ceed.check_error(ierr)?;
2649df49d7eSJed Brown         Ok(Self {
2659df49d7eSJed Brown             op_core: OperatorCore { ceed, ptr },
2669df49d7eSJed Brown         })
2679df49d7eSJed Brown     }
2689df49d7eSJed Brown 
2699df49d7eSJed Brown     fn from_raw(ceed: &'a crate::Ceed, ptr: bind_ceed::CeedOperator) -> crate::Result<Self> {
2709df49d7eSJed Brown         Ok(Self {
2719df49d7eSJed Brown             op_core: OperatorCore { ceed, ptr },
2729df49d7eSJed Brown         })
2739df49d7eSJed Brown     }
2749df49d7eSJed Brown 
2759df49d7eSJed Brown     /// Apply Operator to a vector
2769df49d7eSJed Brown     ///
2779df49d7eSJed Brown     /// * `input`  - Input Vector
2789df49d7eSJed Brown     /// * `output` - Output Vector
2799df49d7eSJed Brown     ///
2809df49d7eSJed Brown     /// ```
2819df49d7eSJed Brown     /// # use libceed::prelude::*;
282*c68be7a2SJeremy L Thompson     /// # fn main() -> Result<(), libceed::CeedError> {
2839df49d7eSJed Brown     /// # let ceed = libceed::Ceed::default_init();
2849df49d7eSJed Brown     /// let ne = 4;
2859df49d7eSJed Brown     /// let p = 3;
2869df49d7eSJed Brown     /// let q = 4;
2879df49d7eSJed Brown     /// let ndofs = p * ne - ne + 1;
2889df49d7eSJed Brown     ///
2899df49d7eSJed Brown     /// // Vectors
290*c68be7a2SJeremy L Thompson     /// let x = ceed.vector_from_slice(&[-1., -0.5, 0.0, 0.5, 1.0])?;
291*c68be7a2SJeremy L Thompson     /// let mut qdata = ceed.vector(ne * q)?;
2929df49d7eSJed Brown     /// qdata.set_value(0.0);
293*c68be7a2SJeremy L Thompson     /// let u = ceed.vector_from_slice(&vec![1.0; ndofs])?;
294*c68be7a2SJeremy L Thompson     /// let mut v = ceed.vector(ndofs)?;
2959df49d7eSJed Brown     /// v.set_value(0.0);
2969df49d7eSJed Brown     ///
2979df49d7eSJed Brown     /// // Restrictions
2989df49d7eSJed Brown     /// let mut indx: Vec<i32> = vec![0; 2 * ne];
2999df49d7eSJed Brown     /// for i in 0..ne {
3009df49d7eSJed Brown     ///     indx[2 * i + 0] = i as i32;
3019df49d7eSJed Brown     ///     indx[2 * i + 1] = (i + 1) as i32;
3029df49d7eSJed Brown     /// }
303*c68be7a2SJeremy L Thompson     /// let rx = ceed.elem_restriction(ne, 2, 1, 1, ne + 1, MemType::Host, &indx)?;
3049df49d7eSJed Brown     /// let mut indu: Vec<i32> = vec![0; p * ne];
3059df49d7eSJed Brown     /// for i in 0..ne {
3069df49d7eSJed Brown     ///     indu[p * i + 0] = i as i32;
3079df49d7eSJed Brown     ///     indu[p * i + 1] = (i + 1) as i32;
3089df49d7eSJed Brown     ///     indu[p * i + 2] = (i + 2) as i32;
3099df49d7eSJed Brown     /// }
310*c68be7a2SJeremy L Thompson     /// let ru = ceed.elem_restriction(ne, 3, 1, 1, ndofs, MemType::Host, &indu)?;
3119df49d7eSJed Brown     /// let strides: [i32; 3] = [1, q as i32, q as i32];
312*c68be7a2SJeremy L Thompson     /// let rq = ceed.strided_elem_restriction(ne, q, 1, q * ne, strides)?;
3139df49d7eSJed Brown     ///
3149df49d7eSJed Brown     /// // Bases
315*c68be7a2SJeremy L Thompson     /// let bx = ceed.basis_tensor_H1_Lagrange(1, 1, 2, q, QuadMode::Gauss)?;
316*c68be7a2SJeremy L Thompson     /// let bu = ceed.basis_tensor_H1_Lagrange(1, 1, p, q, QuadMode::Gauss)?;
3179df49d7eSJed Brown     ///
3189df49d7eSJed Brown     /// // Build quadrature data
319*c68be7a2SJeremy L Thompson     /// let qf_build = ceed.q_function_interior_by_name("Mass1DBuild")?;
320*c68be7a2SJeremy L Thompson     /// ceed.operator(&qf_build, QFunctionOpt::None, QFunctionOpt::None)?
321*c68be7a2SJeremy L Thompson     ///     .field("dx", &rx, &bx, VectorOpt::Active)?
322*c68be7a2SJeremy L Thompson     ///     .field("weights", ElemRestrictionOpt::None, &bx, VectorOpt::None)?
323*c68be7a2SJeremy L Thompson     ///     .field("qdata", &rq, BasisOpt::Collocated, VectorOpt::Active)?
324*c68be7a2SJeremy L Thompson     ///     .apply(&x, &mut qdata)?;
3259df49d7eSJed Brown     ///
3269df49d7eSJed Brown     /// // Mass operator
327*c68be7a2SJeremy L Thompson     /// let qf_mass = ceed.q_function_interior_by_name("MassApply")?;
3289df49d7eSJed Brown     /// let op_mass = ceed
329*c68be7a2SJeremy L Thompson     ///     .operator(&qf_mass, QFunctionOpt::None, QFunctionOpt::None)?
330*c68be7a2SJeremy L Thompson     ///     .field("u", &ru, &bu, VectorOpt::Active)?
331*c68be7a2SJeremy L Thompson     ///     .field("qdata", &rq, BasisOpt::Collocated, &qdata)?
332*c68be7a2SJeremy L Thompson     ///     .field("v", &ru, &bu, VectorOpt::Active)?;
3339df49d7eSJed Brown     ///
3349df49d7eSJed Brown     /// v.set_value(0.0);
335*c68be7a2SJeremy L Thompson     /// op_mass.apply(&u, &mut v)?;
3369df49d7eSJed Brown     ///
3379df49d7eSJed Brown     /// // Check
33880a9ef05SNatalie Beams     /// let sum: Scalar = v.view().iter().sum();
3399df49d7eSJed Brown     /// assert!(
34080a9ef05SNatalie Beams     ///     (sum - 2.0).abs() < 10.0 * libceed::EPSILON,
3419df49d7eSJed Brown     ///     "Incorrect interval length computed"
3429df49d7eSJed Brown     /// );
343*c68be7a2SJeremy L Thompson     /// # Ok(())
344*c68be7a2SJeremy L Thompson     /// # }
3459df49d7eSJed Brown     /// ```
3469df49d7eSJed Brown     pub fn apply(&self, input: &Vector, output: &mut Vector) -> crate::Result<i32> {
3479df49d7eSJed Brown         self.op_core.apply(input, output)
3489df49d7eSJed Brown     }
3499df49d7eSJed Brown 
3509df49d7eSJed Brown     /// Apply Operator to a vector and add result to output vector
3519df49d7eSJed Brown     ///
3529df49d7eSJed Brown     /// * `input`  - Input Vector
3539df49d7eSJed Brown     /// * `output` - Output Vector
3549df49d7eSJed Brown     ///
3559df49d7eSJed Brown     /// ```
3569df49d7eSJed Brown     /// # use libceed::prelude::*;
357*c68be7a2SJeremy L Thompson     /// # fn main() -> Result<(), libceed::CeedError> {
3589df49d7eSJed Brown     /// # let ceed = libceed::Ceed::default_init();
3599df49d7eSJed Brown     /// let ne = 4;
3609df49d7eSJed Brown     /// let p = 3;
3619df49d7eSJed Brown     /// let q = 4;
3629df49d7eSJed Brown     /// let ndofs = p * ne - ne + 1;
3639df49d7eSJed Brown     ///
3649df49d7eSJed Brown     /// // Vectors
365*c68be7a2SJeremy L Thompson     /// let x = ceed.vector_from_slice(&[-1., -0.5, 0.0, 0.5, 1.0])?;
366*c68be7a2SJeremy L Thompson     /// let mut qdata = ceed.vector(ne * q)?;
3679df49d7eSJed Brown     /// qdata.set_value(0.0);
368*c68be7a2SJeremy L Thompson     /// let u = ceed.vector_from_slice(&vec![1.0; ndofs])?;
369*c68be7a2SJeremy L Thompson     /// let mut v = ceed.vector(ndofs)?;
3709df49d7eSJed Brown     ///
3719df49d7eSJed Brown     /// // Restrictions
3729df49d7eSJed Brown     /// let mut indx: Vec<i32> = vec![0; 2 * ne];
3739df49d7eSJed Brown     /// for i in 0..ne {
3749df49d7eSJed Brown     ///     indx[2 * i + 0] = i as i32;
3759df49d7eSJed Brown     ///     indx[2 * i + 1] = (i + 1) as i32;
3769df49d7eSJed Brown     /// }
377*c68be7a2SJeremy L Thompson     /// let rx = ceed.elem_restriction(ne, 2, 1, 1, ne + 1, MemType::Host, &indx)?;
3789df49d7eSJed Brown     /// let mut indu: Vec<i32> = vec![0; p * ne];
3799df49d7eSJed Brown     /// for i in 0..ne {
3809df49d7eSJed Brown     ///     indu[p * i + 0] = i as i32;
3819df49d7eSJed Brown     ///     indu[p * i + 1] = (i + 1) as i32;
3829df49d7eSJed Brown     ///     indu[p * i + 2] = (i + 2) as i32;
3839df49d7eSJed Brown     /// }
384*c68be7a2SJeremy L Thompson     /// let ru = ceed.elem_restriction(ne, 3, 1, 1, ndofs, MemType::Host, &indu)?;
3859df49d7eSJed Brown     /// let strides: [i32; 3] = [1, q as i32, q as i32];
386*c68be7a2SJeremy L Thompson     /// let rq = ceed.strided_elem_restriction(ne, q, 1, q * ne, strides)?;
3879df49d7eSJed Brown     ///
3889df49d7eSJed Brown     /// // Bases
389*c68be7a2SJeremy L Thompson     /// let bx = ceed.basis_tensor_H1_Lagrange(1, 1, 2, q, QuadMode::Gauss)?;
390*c68be7a2SJeremy L Thompson     /// let bu = ceed.basis_tensor_H1_Lagrange(1, 1, p, q, QuadMode::Gauss)?;
3919df49d7eSJed Brown     ///
3929df49d7eSJed Brown     /// // Build quadrature data
393*c68be7a2SJeremy L Thompson     /// let qf_build = ceed.q_function_interior_by_name("Mass1DBuild")?;
394*c68be7a2SJeremy L Thompson     /// ceed.operator(&qf_build, QFunctionOpt::None, QFunctionOpt::None)?
395*c68be7a2SJeremy L Thompson     ///     .field("dx", &rx, &bx, VectorOpt::Active)?
396*c68be7a2SJeremy L Thompson     ///     .field("weights", ElemRestrictionOpt::None, &bx, VectorOpt::None)?
397*c68be7a2SJeremy L Thompson     ///     .field("qdata", &rq, BasisOpt::Collocated, VectorOpt::Active)?
398*c68be7a2SJeremy L Thompson     ///     .apply(&x, &mut qdata)?;
3999df49d7eSJed Brown     ///
4009df49d7eSJed Brown     /// // Mass operator
401*c68be7a2SJeremy L Thompson     /// let qf_mass = ceed.q_function_interior_by_name("MassApply")?;
4029df49d7eSJed Brown     /// let op_mass = ceed
403*c68be7a2SJeremy L Thompson     ///     .operator(&qf_mass, QFunctionOpt::None, QFunctionOpt::None)?
404*c68be7a2SJeremy L Thompson     ///     .field("u", &ru, &bu, VectorOpt::Active)?
405*c68be7a2SJeremy L Thompson     ///     .field("qdata", &rq, BasisOpt::Collocated, &qdata)?
406*c68be7a2SJeremy L Thompson     ///     .field("v", &ru, &bu, VectorOpt::Active)?;
4079df49d7eSJed Brown     ///
4089df49d7eSJed Brown     /// v.set_value(1.0);
409*c68be7a2SJeremy L Thompson     /// op_mass.apply_add(&u, &mut v)?;
4109df49d7eSJed Brown     ///
4119df49d7eSJed Brown     /// // Check
41280a9ef05SNatalie Beams     /// let sum: Scalar = v.view().iter().sum();
4139df49d7eSJed Brown     /// assert!(
41480a9ef05SNatalie Beams     ///     (sum - (2.0 + ndofs as Scalar)).abs() < 10.0 * libceed::EPSILON,
4159df49d7eSJed Brown     ///     "Incorrect interval length computed and added"
4169df49d7eSJed Brown     /// );
417*c68be7a2SJeremy L Thompson     /// # Ok(())
418*c68be7a2SJeremy L Thompson     /// # }
4199df49d7eSJed Brown     /// ```
4209df49d7eSJed Brown     pub fn apply_add(&self, input: &Vector, output: &mut Vector) -> crate::Result<i32> {
4219df49d7eSJed Brown         self.op_core.apply_add(input, output)
4229df49d7eSJed Brown     }
4239df49d7eSJed Brown 
4249df49d7eSJed Brown     /// Provide a field to a Operator for use by its QFunction
4259df49d7eSJed Brown     ///
4269df49d7eSJed Brown     /// * `fieldname` - Name of the field (to be matched with the name used by
4279df49d7eSJed Brown     ///                   the QFunction)
4289df49d7eSJed Brown     /// * `r`         - ElemRestriction
4299df49d7eSJed Brown     /// * `b`         - Basis in which the field resides or `BasisCollocated` if
4309df49d7eSJed Brown     ///                   collocated with quadrature points
4319df49d7eSJed Brown     /// * `v`         - Vector to be used by Operator or `VectorActive` if field
4329df49d7eSJed Brown     ///                   is active or `VectorNone` if using `Weight` with the
4339df49d7eSJed Brown     ///                   QFunction
4349df49d7eSJed Brown     ///
4359df49d7eSJed Brown     ///
4369df49d7eSJed Brown     /// ```
4379df49d7eSJed Brown     /// # use libceed::prelude::*;
438*c68be7a2SJeremy L Thompson     /// # fn main() -> Result<(), libceed::CeedError> {
4399df49d7eSJed Brown     /// # let ceed = libceed::Ceed::default_init();
440*c68be7a2SJeremy L Thompson     /// let qf = ceed.q_function_interior_by_name("Mass1DBuild")?;
441*c68be7a2SJeremy L Thompson     /// let mut op = ceed.operator(&qf, QFunctionOpt::None, QFunctionOpt::None)?;
4429df49d7eSJed Brown     ///
4439df49d7eSJed Brown     /// // Operator field arguments
4449df49d7eSJed Brown     /// let ne = 3;
4459df49d7eSJed Brown     /// let q = 4;
4469df49d7eSJed Brown     /// let mut ind: Vec<i32> = vec![0; 2 * ne];
4479df49d7eSJed Brown     /// for i in 0..ne {
4489df49d7eSJed Brown     ///     ind[2 * i + 0] = i as i32;
4499df49d7eSJed Brown     ///     ind[2 * i + 1] = (i + 1) as i32;
4509df49d7eSJed Brown     /// }
451*c68be7a2SJeremy L Thompson     /// let r = ceed.elem_restriction(ne, 2, 1, 1, ne + 1, MemType::Host, &ind)?;
4529df49d7eSJed Brown     ///
453*c68be7a2SJeremy L Thompson     /// let b = ceed.basis_tensor_H1_Lagrange(1, 1, 2, q, QuadMode::Gauss)?;
4549df49d7eSJed Brown     ///
4559df49d7eSJed Brown     /// // Operator field
456*c68be7a2SJeremy L Thompson     /// op = op.field("dx", &r, &b, VectorOpt::Active)?;
457*c68be7a2SJeremy L Thompson     /// # Ok(())
458*c68be7a2SJeremy L Thompson     /// # }
4599df49d7eSJed Brown     /// ```
4609df49d7eSJed Brown     #[allow(unused_mut)]
4619df49d7eSJed Brown     pub fn field<'b>(
4629df49d7eSJed Brown         mut self,
4639df49d7eSJed Brown         fieldname: &str,
4649df49d7eSJed Brown         r: impl Into<ElemRestrictionOpt<'b>>,
4659df49d7eSJed Brown         b: impl Into<BasisOpt<'b>>,
4669df49d7eSJed Brown         v: impl Into<VectorOpt<'b>>,
4679df49d7eSJed Brown     ) -> crate::Result<Self> {
4689df49d7eSJed Brown         let fieldname = CString::new(fieldname).expect("CString::new failed");
4699df49d7eSJed Brown         let fieldname = fieldname.as_ptr() as *const i8;
4709df49d7eSJed Brown         let ierr = unsafe {
4719df49d7eSJed Brown             bind_ceed::CeedOperatorSetField(
4729df49d7eSJed Brown                 self.op_core.ptr,
4739df49d7eSJed Brown                 fieldname,
4749df49d7eSJed Brown                 r.into().to_raw(),
4759df49d7eSJed Brown                 b.into().to_raw(),
4769df49d7eSJed Brown                 v.into().to_raw(),
4779df49d7eSJed Brown             )
4789df49d7eSJed Brown         };
4799df49d7eSJed Brown         self.op_core.ceed.check_error(ierr)?;
4809df49d7eSJed Brown         Ok(self)
4819df49d7eSJed Brown     }
4829df49d7eSJed Brown 
4839df49d7eSJed Brown     /// Assemble the diagonal of a square linear Operator
4849df49d7eSJed Brown     ///
4859df49d7eSJed Brown     /// This overwrites a Vector with the diagonal of a linear Operator.
4869df49d7eSJed Brown     ///
4879df49d7eSJed Brown     /// Note: Currently only non-composite Operators with a single field and
4889df49d7eSJed Brown     ///      composite Operators with single field sub-operators are supported.
4899df49d7eSJed Brown     ///
4909df49d7eSJed Brown     /// * `op`        - Operator to assemble QFunction
4919df49d7eSJed Brown     /// * `assembled` - Vector to store assembled Operator diagonal
4929df49d7eSJed Brown     ///
4939df49d7eSJed Brown     /// ```
4949df49d7eSJed Brown     /// # use libceed::prelude::*;
495*c68be7a2SJeremy L Thompson     /// # fn main() -> Result<(), libceed::CeedError> {
4969df49d7eSJed Brown     /// # let ceed = libceed::Ceed::default_init();
4979df49d7eSJed Brown     /// let ne = 4;
4989df49d7eSJed Brown     /// let p = 3;
4999df49d7eSJed Brown     /// let q = 4;
5009df49d7eSJed Brown     /// let ndofs = p * ne - ne + 1;
5019df49d7eSJed Brown     ///
5029df49d7eSJed Brown     /// // Vectors
503*c68be7a2SJeremy L Thompson     /// let x = ceed.vector_from_slice(&[-1., -0.5, 0.0, 0.5, 1.0])?;
504*c68be7a2SJeremy L Thompson     /// let mut qdata = ceed.vector(ne * q)?;
5059df49d7eSJed Brown     /// qdata.set_value(0.0);
506*c68be7a2SJeremy L Thompson     /// let mut u = ceed.vector(ndofs)?;
5079df49d7eSJed Brown     /// u.set_value(1.0);
508*c68be7a2SJeremy L Thompson     /// let mut v = ceed.vector(ndofs)?;
5099df49d7eSJed Brown     /// v.set_value(0.0);
5109df49d7eSJed Brown     ///
5119df49d7eSJed Brown     /// // Restrictions
5129df49d7eSJed Brown     /// let mut indx: Vec<i32> = vec![0; 2 * ne];
5139df49d7eSJed Brown     /// for i in 0..ne {
5149df49d7eSJed Brown     ///     indx[2 * i + 0] = i as i32;
5159df49d7eSJed Brown     ///     indx[2 * i + 1] = (i + 1) as i32;
5169df49d7eSJed Brown     /// }
517*c68be7a2SJeremy L Thompson     /// let rx = ceed.elem_restriction(ne, 2, 1, 1, ne + 1, MemType::Host, &indx)?;
5189df49d7eSJed Brown     /// let mut indu: Vec<i32> = vec![0; p * ne];
5199df49d7eSJed Brown     /// for i in 0..ne {
5209df49d7eSJed Brown     ///     indu[p * i + 0] = (2 * i) as i32;
5219df49d7eSJed Brown     ///     indu[p * i + 1] = (2 * i + 1) as i32;
5229df49d7eSJed Brown     ///     indu[p * i + 2] = (2 * i + 2) as i32;
5239df49d7eSJed Brown     /// }
524*c68be7a2SJeremy L Thompson     /// let ru = ceed.elem_restriction(ne, p, 1, 1, ndofs, MemType::Host, &indu)?;
5259df49d7eSJed Brown     /// let strides: [i32; 3] = [1, q as i32, q as i32];
526*c68be7a2SJeremy L Thompson     /// let rq = ceed.strided_elem_restriction(ne, q, 1, q * ne, strides)?;
5279df49d7eSJed Brown     ///
5289df49d7eSJed Brown     /// // Bases
529*c68be7a2SJeremy L Thompson     /// let bx = ceed.basis_tensor_H1_Lagrange(1, 1, 2, q, QuadMode::Gauss)?;
530*c68be7a2SJeremy L Thompson     /// let bu = ceed.basis_tensor_H1_Lagrange(1, 1, p, q, QuadMode::Gauss)?;
5319df49d7eSJed Brown     ///
5329df49d7eSJed Brown     /// // Build quadrature data
533*c68be7a2SJeremy L Thompson     /// let qf_build = ceed.q_function_interior_by_name("Mass1DBuild")?;
534*c68be7a2SJeremy L Thompson     /// ceed.operator(&qf_build, QFunctionOpt::None, QFunctionOpt::None)?
535*c68be7a2SJeremy L Thompson     ///     .field("dx", &rx, &bx, VectorOpt::Active)?
536*c68be7a2SJeremy L Thompson     ///     .field("weights", ElemRestrictionOpt::None, &bx, VectorOpt::None)?
537*c68be7a2SJeremy L Thompson     ///     .field("qdata", &rq, BasisOpt::Collocated, VectorOpt::Active)?
538*c68be7a2SJeremy L Thompson     ///     .apply(&x, &mut qdata)?;
5399df49d7eSJed Brown     ///
5409df49d7eSJed Brown     /// // Mass operator
541*c68be7a2SJeremy L Thompson     /// let qf_mass = ceed.q_function_interior_by_name("MassApply")?;
5429df49d7eSJed Brown     /// let op_mass = ceed
543*c68be7a2SJeremy L Thompson     ///     .operator(&qf_mass, QFunctionOpt::None, QFunctionOpt::None)?
544*c68be7a2SJeremy L Thompson     ///     .field("u", &ru, &bu, VectorOpt::Active)?
545*c68be7a2SJeremy L Thompson     ///     .field("qdata", &rq, BasisOpt::Collocated, &qdata)?
546*c68be7a2SJeremy L Thompson     ///     .field("v", &ru, &bu, VectorOpt::Active)?;
5479df49d7eSJed Brown     ///
5489df49d7eSJed Brown     /// // Diagonal
549*c68be7a2SJeremy L Thompson     /// let mut diag = ceed.vector(ndofs)?;
5509df49d7eSJed Brown     /// diag.set_value(0.0);
551*c68be7a2SJeremy L Thompson     /// op_mass.linear_assemble_diagonal(&mut diag)?;
5529df49d7eSJed Brown     ///
5539df49d7eSJed Brown     /// // Manual diagonal computation
554*c68be7a2SJeremy L Thompson     /// let mut true_diag = ceed.vector(ndofs)?;
5559df49d7eSJed Brown     /// for i in 0..ndofs {
5569df49d7eSJed Brown     ///     u.set_value(0.0);
5579df49d7eSJed Brown     ///     {
5589df49d7eSJed Brown     ///         let mut u_array = u.view_mut();
5599df49d7eSJed Brown     ///         u_array[i] = 1.;
5609df49d7eSJed Brown     ///     }
5619df49d7eSJed Brown     ///
562*c68be7a2SJeremy L Thompson     ///     op_mass.apply(&u, &mut v)?;
5639df49d7eSJed Brown     ///
5649df49d7eSJed Brown     ///     {
5659df49d7eSJed Brown     ///         let v_array = v.view_mut();
5669df49d7eSJed Brown     ///         let mut true_array = true_diag.view_mut();
5679df49d7eSJed Brown     ///         true_array[i] = v_array[i];
5689df49d7eSJed Brown     ///     }
5699df49d7eSJed Brown     /// }
5709df49d7eSJed Brown     ///
5719df49d7eSJed Brown     /// // Check
5729df49d7eSJed Brown     /// diag.view()
5739df49d7eSJed Brown     ///     .iter()
5749df49d7eSJed Brown     ///     .zip(true_diag.view().iter())
5759df49d7eSJed Brown     ///     .for_each(|(computed, actual)| {
5769df49d7eSJed Brown     ///         assert!(
57780a9ef05SNatalie Beams     ///             (*computed - *actual).abs() < 10.0 * libceed::EPSILON,
5789df49d7eSJed Brown     ///             "Diagonal entry incorrect"
5799df49d7eSJed Brown     ///         );
5809df49d7eSJed Brown     ///     });
581*c68be7a2SJeremy L Thompson     /// # Ok(())
582*c68be7a2SJeremy L Thompson     /// # }
5839df49d7eSJed Brown     /// ```
5849df49d7eSJed Brown     pub fn linear_assemble_diagonal(&self, assembled: &mut Vector) -> crate::Result<i32> {
5859df49d7eSJed Brown         self.op_core.linear_assemble_diagonal(assembled)
5869df49d7eSJed Brown     }
5879df49d7eSJed Brown 
5889df49d7eSJed Brown     /// Assemble the diagonal of a square linear Operator
5899df49d7eSJed Brown     ///
5909df49d7eSJed Brown     /// This sums into a Vector with the diagonal of a linear Operator.
5919df49d7eSJed Brown     ///
5929df49d7eSJed Brown     /// Note: Currently only non-composite Operators with a single field and
5939df49d7eSJed Brown     ///      composite Operators with single field sub-operators are supported.
5949df49d7eSJed Brown     ///
5959df49d7eSJed Brown     /// * `op`        - Operator to assemble QFunction
5969df49d7eSJed Brown     /// * `assembled` - Vector to store assembled Operator diagonal
5979df49d7eSJed Brown     ///
5989df49d7eSJed Brown     ///
5999df49d7eSJed Brown     /// ```
6009df49d7eSJed Brown     /// # use libceed::prelude::*;
601*c68be7a2SJeremy L Thompson     /// # fn main() -> Result<(), libceed::CeedError> {
6029df49d7eSJed Brown     /// # let ceed = libceed::Ceed::default_init();
6039df49d7eSJed Brown     /// let ne = 4;
6049df49d7eSJed Brown     /// let p = 3;
6059df49d7eSJed Brown     /// let q = 4;
6069df49d7eSJed Brown     /// let ndofs = p * ne - ne + 1;
6079df49d7eSJed Brown     ///
6089df49d7eSJed Brown     /// // Vectors
609*c68be7a2SJeremy L Thompson     /// let x = ceed.vector_from_slice(&[-1., -0.5, 0.0, 0.5, 1.0])?;
610*c68be7a2SJeremy L Thompson     /// let mut qdata = ceed.vector(ne * q)?;
6119df49d7eSJed Brown     /// qdata.set_value(0.0);
612*c68be7a2SJeremy L Thompson     /// let mut u = ceed.vector(ndofs)?;
6139df49d7eSJed Brown     /// u.set_value(1.0);
614*c68be7a2SJeremy L Thompson     /// let mut v = ceed.vector(ndofs)?;
6159df49d7eSJed Brown     /// v.set_value(0.0);
6169df49d7eSJed Brown     ///
6179df49d7eSJed Brown     /// // Restrictions
6189df49d7eSJed Brown     /// let mut indx: Vec<i32> = vec![0; 2 * ne];
6199df49d7eSJed Brown     /// for i in 0..ne {
6209df49d7eSJed Brown     ///     indx[2 * i + 0] = i as i32;
6219df49d7eSJed Brown     ///     indx[2 * i + 1] = (i + 1) as i32;
6229df49d7eSJed Brown     /// }
623*c68be7a2SJeremy L Thompson     /// let rx = ceed.elem_restriction(ne, 2, 1, 1, ne + 1, MemType::Host, &indx)?;
6249df49d7eSJed Brown     /// let mut indu: Vec<i32> = vec![0; p * ne];
6259df49d7eSJed Brown     /// for i in 0..ne {
6269df49d7eSJed Brown     ///     indu[p * i + 0] = (2 * i) as i32;
6279df49d7eSJed Brown     ///     indu[p * i + 1] = (2 * i + 1) as i32;
6289df49d7eSJed Brown     ///     indu[p * i + 2] = (2 * i + 2) as i32;
6299df49d7eSJed Brown     /// }
630*c68be7a2SJeremy L Thompson     /// let ru = ceed.elem_restriction(ne, p, 1, 1, ndofs, MemType::Host, &indu)?;
6319df49d7eSJed Brown     /// let strides: [i32; 3] = [1, q as i32, q as i32];
632*c68be7a2SJeremy L Thompson     /// let rq = ceed.strided_elem_restriction(ne, q, 1, q * ne, strides)?;
6339df49d7eSJed Brown     ///
6349df49d7eSJed Brown     /// // Bases
635*c68be7a2SJeremy L Thompson     /// let bx = ceed.basis_tensor_H1_Lagrange(1, 1, 2, q, QuadMode::Gauss)?;
636*c68be7a2SJeremy L Thompson     /// let bu = ceed.basis_tensor_H1_Lagrange(1, 1, p, q, QuadMode::Gauss)?;
6379df49d7eSJed Brown     ///
6389df49d7eSJed Brown     /// // Build quadrature data
639*c68be7a2SJeremy L Thompson     /// let qf_build = ceed.q_function_interior_by_name("Mass1DBuild")?;
640*c68be7a2SJeremy L Thompson     /// ceed.operator(&qf_build, QFunctionOpt::None, QFunctionOpt::None)?
641*c68be7a2SJeremy L Thompson     ///     .field("dx", &rx, &bx, VectorOpt::Active)?
642*c68be7a2SJeremy L Thompson     ///     .field("weights", ElemRestrictionOpt::None, &bx, VectorOpt::None)?
643*c68be7a2SJeremy L Thompson     ///     .field("qdata", &rq, BasisOpt::Collocated, VectorOpt::Active)?
644*c68be7a2SJeremy L Thompson     ///     .apply(&x, &mut qdata)?;
6459df49d7eSJed Brown     ///
6469df49d7eSJed Brown     /// // Mass operator
647*c68be7a2SJeremy L Thompson     /// let qf_mass = ceed.q_function_interior_by_name("MassApply")?;
6489df49d7eSJed Brown     /// let op_mass = ceed
649*c68be7a2SJeremy L Thompson     ///     .operator(&qf_mass, QFunctionOpt::None, QFunctionOpt::None)?
650*c68be7a2SJeremy L Thompson     ///     .field("u", &ru, &bu, VectorOpt::Active)?
651*c68be7a2SJeremy L Thompson     ///     .field("qdata", &rq, BasisOpt::Collocated, &qdata)?
652*c68be7a2SJeremy L Thompson     ///     .field("v", &ru, &bu, VectorOpt::Active)?;
6539df49d7eSJed Brown     ///
6549df49d7eSJed Brown     /// // Diagonal
655*c68be7a2SJeremy L Thompson     /// let mut diag = ceed.vector(ndofs)?;
6569df49d7eSJed Brown     /// diag.set_value(1.0);
657*c68be7a2SJeremy L Thompson     /// op_mass.linear_assemble_add_diagonal(&mut diag)?;
6589df49d7eSJed Brown     ///
6599df49d7eSJed Brown     /// // Manual diagonal computation
660*c68be7a2SJeremy L Thompson     /// let mut true_diag = ceed.vector(ndofs)?;
6619df49d7eSJed Brown     /// for i in 0..ndofs {
6629df49d7eSJed Brown     ///     u.set_value(0.0);
6639df49d7eSJed Brown     ///     {
6649df49d7eSJed Brown     ///         let mut u_array = u.view_mut();
6659df49d7eSJed Brown     ///         u_array[i] = 1.;
6669df49d7eSJed Brown     ///     }
6679df49d7eSJed Brown     ///
668*c68be7a2SJeremy L Thompson     ///     op_mass.apply(&u, &mut v)?;
6699df49d7eSJed Brown     ///
6709df49d7eSJed Brown     ///     {
6719df49d7eSJed Brown     ///         let v_array = v.view_mut();
6729df49d7eSJed Brown     ///         let mut true_array = true_diag.view_mut();
6739df49d7eSJed Brown     ///         true_array[i] = v_array[i] + 1.0;
6749df49d7eSJed Brown     ///     }
6759df49d7eSJed Brown     /// }
6769df49d7eSJed Brown     ///
6779df49d7eSJed Brown     /// // Check
6789df49d7eSJed Brown     /// diag.view()
6799df49d7eSJed Brown     ///     .iter()
6809df49d7eSJed Brown     ///     .zip(true_diag.view().iter())
6819df49d7eSJed Brown     ///     .for_each(|(computed, actual)| {
6829df49d7eSJed Brown     ///         assert!(
68380a9ef05SNatalie Beams     ///             (*computed - *actual).abs() < 10.0 * libceed::EPSILON,
6849df49d7eSJed Brown     ///             "Diagonal entry incorrect"
6859df49d7eSJed Brown     ///         );
6869df49d7eSJed Brown     ///     });
687*c68be7a2SJeremy L Thompson     /// # Ok(())
688*c68be7a2SJeremy L Thompson     /// # }
6899df49d7eSJed Brown     /// ```
6909df49d7eSJed Brown     pub fn linear_assemble_add_diagonal(&self, assembled: &mut Vector) -> crate::Result<i32> {
6919df49d7eSJed Brown         self.op_core.linear_assemble_add_diagonal(assembled)
6929df49d7eSJed Brown     }
6939df49d7eSJed Brown 
6949df49d7eSJed Brown     /// Assemble the point block diagonal of a square linear Operator
6959df49d7eSJed Brown     ///
6969df49d7eSJed Brown     /// This overwrites a Vector with the point block diagonal of a linear
6979df49d7eSJed Brown     /// Operator.
6989df49d7eSJed Brown     ///
6999df49d7eSJed Brown     /// Note: Currently only non-composite Operators with a single field and
7009df49d7eSJed Brown     ///      composite Operators with single field sub-operators are supported.
7019df49d7eSJed Brown     ///
7029df49d7eSJed Brown     /// * `op`        - Operator to assemble QFunction
7039df49d7eSJed Brown     /// * `assembled` - Vector to store assembled CeedOperator point block
7049df49d7eSJed Brown     ///                   diagonal, provided in row-major form with an
7059df49d7eSJed Brown     ///                   `ncomp * ncomp` block at each node. The dimensions of
7069df49d7eSJed Brown     ///                   this vector are derived from the active vector for
7079df49d7eSJed Brown     ///                   the CeedOperator. The array has shape
7089df49d7eSJed Brown     ///                   `[nodes, component out, component in]`.
7099df49d7eSJed Brown     ///
7109df49d7eSJed Brown     /// ```
7119df49d7eSJed Brown     /// # use libceed::prelude::*;
712*c68be7a2SJeremy L Thompson     /// # fn main() -> Result<(), libceed::CeedError> {
7139df49d7eSJed Brown     /// # let ceed = libceed::Ceed::default_init();
7149df49d7eSJed Brown     /// let ne = 4;
7159df49d7eSJed Brown     /// let p = 3;
7169df49d7eSJed Brown     /// let q = 4;
7179df49d7eSJed Brown     /// let ncomp = 2;
7189df49d7eSJed Brown     /// let ndofs = p * ne - ne + 1;
7199df49d7eSJed Brown     ///
7209df49d7eSJed Brown     /// // Vectors
721*c68be7a2SJeremy L Thompson     /// let x = ceed.vector_from_slice(&[-1., -0.5, 0.0, 0.5, 1.0])?;
722*c68be7a2SJeremy L Thompson     /// let mut qdata = ceed.vector(ne * q)?;
7239df49d7eSJed Brown     /// qdata.set_value(0.0);
724*c68be7a2SJeremy L Thompson     /// let mut u = ceed.vector(ncomp * ndofs)?;
7259df49d7eSJed Brown     /// u.set_value(1.0);
726*c68be7a2SJeremy L Thompson     /// let mut v = ceed.vector(ncomp * ndofs)?;
7279df49d7eSJed Brown     /// v.set_value(0.0);
7289df49d7eSJed Brown     ///
7299df49d7eSJed Brown     /// // Restrictions
7309df49d7eSJed Brown     /// let mut indx: Vec<i32> = vec![0; 2 * ne];
7319df49d7eSJed Brown     /// for i in 0..ne {
7329df49d7eSJed Brown     ///     indx[2 * i + 0] = i as i32;
7339df49d7eSJed Brown     ///     indx[2 * i + 1] = (i + 1) as i32;
7349df49d7eSJed Brown     /// }
735*c68be7a2SJeremy L Thompson     /// let rx = ceed.elem_restriction(ne, 2, 1, 1, ne + 1, MemType::Host, &indx)?;
7369df49d7eSJed Brown     /// let mut indu: Vec<i32> = vec![0; p * ne];
7379df49d7eSJed Brown     /// for i in 0..ne {
7389df49d7eSJed Brown     ///     indu[p * i + 0] = (2 * i) as i32;
7399df49d7eSJed Brown     ///     indu[p * i + 1] = (2 * i + 1) as i32;
7409df49d7eSJed Brown     ///     indu[p * i + 2] = (2 * i + 2) as i32;
7419df49d7eSJed Brown     /// }
742*c68be7a2SJeremy L Thompson     /// let ru = ceed.elem_restriction(ne, p, ncomp, ndofs, ncomp * ndofs, MemType::Host, &indu)?;
7439df49d7eSJed Brown     /// let strides: [i32; 3] = [1, q as i32, q as i32];
744*c68be7a2SJeremy L Thompson     /// let rq = ceed.strided_elem_restriction(ne, q, 1, q * ne, strides)?;
7459df49d7eSJed Brown     ///
7469df49d7eSJed Brown     /// // Bases
747*c68be7a2SJeremy L Thompson     /// let bx = ceed.basis_tensor_H1_Lagrange(1, 1, 2, q, QuadMode::Gauss)?;
748*c68be7a2SJeremy L Thompson     /// let bu = ceed.basis_tensor_H1_Lagrange(1, ncomp, p, q, QuadMode::Gauss)?;
7499df49d7eSJed Brown     ///
7509df49d7eSJed Brown     /// // Build quadrature data
751*c68be7a2SJeremy L Thompson     /// let qf_build = ceed.q_function_interior_by_name("Mass1DBuild")?;
752*c68be7a2SJeremy L Thompson     /// ceed.operator(&qf_build, QFunctionOpt::None, QFunctionOpt::None)?
753*c68be7a2SJeremy L Thompson     ///     .field("dx", &rx, &bx, VectorOpt::Active)?
754*c68be7a2SJeremy L Thompson     ///     .field("weights", ElemRestrictionOpt::None, &bx, VectorOpt::None)?
755*c68be7a2SJeremy L Thompson     ///     .field("qdata", &rq, BasisOpt::Collocated, VectorOpt::Active)?
756*c68be7a2SJeremy L Thompson     ///     .apply(&x, &mut qdata)?;
7579df49d7eSJed Brown     ///
7589df49d7eSJed Brown     /// // Mass operator
7599df49d7eSJed Brown     /// let mut mass_2_comp = |[u, qdata, ..]: QFunctionInputs, [v, ..]: QFunctionOutputs| {
7609df49d7eSJed Brown     ///     // Number of quadrature points
7619df49d7eSJed Brown     ///     let q = qdata.len();
7629df49d7eSJed Brown     ///
7639df49d7eSJed Brown     ///     // Iterate over quadrature points
7649df49d7eSJed Brown     ///     for i in 0..q {
7659df49d7eSJed Brown     ///         v[i + 0 * q] = u[i + 1 * q] * qdata[i];
7669df49d7eSJed Brown     ///         v[i + 1 * q] = u[i + 0 * q] * qdata[i];
7679df49d7eSJed Brown     ///     }
7689df49d7eSJed Brown     ///
7699df49d7eSJed Brown     ///     // Return clean error code
7709df49d7eSJed Brown     ///     0
7719df49d7eSJed Brown     /// };
7729df49d7eSJed Brown     ///
7739df49d7eSJed Brown     /// let qf_mass = ceed
774*c68be7a2SJeremy L Thompson     ///     .q_function_interior(1, Box::new(mass_2_comp))?
775*c68be7a2SJeremy L Thompson     ///     .input("u", 2, EvalMode::Interp)?
776*c68be7a2SJeremy L Thompson     ///     .input("qdata", 1, EvalMode::None)?
777*c68be7a2SJeremy L Thompson     ///     .output("v", 2, EvalMode::Interp)?;
7789df49d7eSJed Brown     ///
7799df49d7eSJed Brown     /// let op_mass = ceed
780*c68be7a2SJeremy L Thompson     ///     .operator(&qf_mass, QFunctionOpt::None, QFunctionOpt::None)?
781*c68be7a2SJeremy L Thompson     ///     .field("u", &ru, &bu, VectorOpt::Active)?
782*c68be7a2SJeremy L Thompson     ///     .field("qdata", &rq, BasisOpt::Collocated, &qdata)?
783*c68be7a2SJeremy L Thompson     ///     .field("v", &ru, &bu, VectorOpt::Active)?;
7849df49d7eSJed Brown     ///
7859df49d7eSJed Brown     /// // Diagonal
786*c68be7a2SJeremy L Thompson     /// let mut diag = ceed.vector(ncomp * ncomp * ndofs)?;
7879df49d7eSJed Brown     /// diag.set_value(0.0);
788*c68be7a2SJeremy L Thompson     /// op_mass.linear_assemble_point_block_diagonal(&mut diag)?;
7899df49d7eSJed Brown     ///
7909df49d7eSJed Brown     /// // Manual diagonal computation
791*c68be7a2SJeremy L Thompson     /// let mut true_diag = ceed.vector(ncomp * ncomp * ndofs)?;
7929df49d7eSJed Brown     /// for i in 0..ndofs {
7939df49d7eSJed Brown     ///     for j in 0..ncomp {
7949df49d7eSJed Brown     ///         u.set_value(0.0);
7959df49d7eSJed Brown     ///         {
7969df49d7eSJed Brown     ///             let mut u_array = u.view_mut();
7979df49d7eSJed Brown     ///             u_array[i + j * ndofs] = 1.;
7989df49d7eSJed Brown     ///         }
7999df49d7eSJed Brown     ///
800*c68be7a2SJeremy L Thompson     ///         op_mass.apply(&u, &mut v)?;
8019df49d7eSJed Brown     ///
8029df49d7eSJed Brown     ///         {
8039df49d7eSJed Brown     ///             let v_array = v.view_mut();
8049df49d7eSJed Brown     ///             let mut true_array = true_diag.view_mut();
8059df49d7eSJed Brown     ///             for k in 0..ncomp {
8069df49d7eSJed Brown     ///                 true_array[i * ncomp * ncomp + k * ncomp + j] = v_array[i + k * ndofs];
8079df49d7eSJed Brown     ///             }
8089df49d7eSJed Brown     ///         }
8099df49d7eSJed Brown     ///     }
8109df49d7eSJed Brown     /// }
8119df49d7eSJed Brown     ///
8129df49d7eSJed Brown     /// // Check
8139df49d7eSJed Brown     /// diag.view()
8149df49d7eSJed Brown     ///     .iter()
8159df49d7eSJed Brown     ///     .zip(true_diag.view().iter())
8169df49d7eSJed Brown     ///     .for_each(|(computed, actual)| {
8179df49d7eSJed Brown     ///         assert!(
81880a9ef05SNatalie Beams     ///             (*computed - *actual).abs() < 10.0 * libceed::EPSILON,
8199df49d7eSJed Brown     ///             "Diagonal entry incorrect"
8209df49d7eSJed Brown     ///         );
8219df49d7eSJed Brown     ///     });
822*c68be7a2SJeremy L Thompson     /// # Ok(())
823*c68be7a2SJeremy L Thompson     /// # }
8249df49d7eSJed Brown     /// ```
8259df49d7eSJed Brown     pub fn linear_assemble_point_block_diagonal(
8269df49d7eSJed Brown         &self,
8279df49d7eSJed Brown         assembled: &mut Vector,
8289df49d7eSJed Brown     ) -> crate::Result<i32> {
8299df49d7eSJed Brown         self.op_core.linear_assemble_point_block_diagonal(assembled)
8309df49d7eSJed Brown     }
8319df49d7eSJed Brown 
8329df49d7eSJed Brown     /// Assemble the point block diagonal of a square linear Operator
8339df49d7eSJed Brown     ///
8349df49d7eSJed Brown     /// This sums into a Vector with the point block diagonal of a linear
8359df49d7eSJed Brown     /// Operator.
8369df49d7eSJed Brown     ///
8379df49d7eSJed Brown     /// Note: Currently only non-composite Operators with a single field and
8389df49d7eSJed Brown     ///      composite Operators with single field sub-operators are supported.
8399df49d7eSJed Brown     ///
8409df49d7eSJed Brown     /// * `op`        -     Operator to assemble QFunction
8419df49d7eSJed Brown     /// * `assembled` - Vector to store assembled CeedOperator point block
8429df49d7eSJed Brown     ///                   diagonal, provided in row-major form with an
8439df49d7eSJed Brown     ///                   `ncomp * ncomp` block at each node. The dimensions of
8449df49d7eSJed Brown     ///                   this vector are derived from the active vector for
8459df49d7eSJed Brown     ///                   the CeedOperator. The array has shape
8469df49d7eSJed Brown     ///                   `[nodes, component out, component in]`.
8479df49d7eSJed Brown     ///
8489df49d7eSJed Brown     /// ```
8499df49d7eSJed Brown     /// # use libceed::prelude::*;
850*c68be7a2SJeremy L Thompson     /// # fn main() -> Result<(), libceed::CeedError> {
8519df49d7eSJed Brown     /// # let ceed = libceed::Ceed::default_init();
8529df49d7eSJed Brown     /// let ne = 4;
8539df49d7eSJed Brown     /// let p = 3;
8549df49d7eSJed Brown     /// let q = 4;
8559df49d7eSJed Brown     /// let ncomp = 2;
8569df49d7eSJed Brown     /// let ndofs = p * ne - ne + 1;
8579df49d7eSJed Brown     ///
8589df49d7eSJed Brown     /// // Vectors
859*c68be7a2SJeremy L Thompson     /// let x = ceed.vector_from_slice(&[-1., -0.5, 0.0, 0.5, 1.0])?;
860*c68be7a2SJeremy L Thompson     /// let mut qdata = ceed.vector(ne * q)?;
8619df49d7eSJed Brown     /// qdata.set_value(0.0);
862*c68be7a2SJeremy L Thompson     /// let mut u = ceed.vector(ncomp * ndofs)?;
8639df49d7eSJed Brown     /// u.set_value(1.0);
864*c68be7a2SJeremy L Thompson     /// let mut v = ceed.vector(ncomp * ndofs)?;
8659df49d7eSJed Brown     /// v.set_value(0.0);
8669df49d7eSJed Brown     ///
8679df49d7eSJed Brown     /// // Restrictions
8689df49d7eSJed Brown     /// let mut indx: Vec<i32> = vec![0; 2 * ne];
8699df49d7eSJed Brown     /// for i in 0..ne {
8709df49d7eSJed Brown     ///     indx[2 * i + 0] = i as i32;
8719df49d7eSJed Brown     ///     indx[2 * i + 1] = (i + 1) as i32;
8729df49d7eSJed Brown     /// }
873*c68be7a2SJeremy L Thompson     /// let rx = ceed.elem_restriction(ne, 2, 1, 1, ne + 1, MemType::Host, &indx)?;
8749df49d7eSJed Brown     /// let mut indu: Vec<i32> = vec![0; p * ne];
8759df49d7eSJed Brown     /// for i in 0..ne {
8769df49d7eSJed Brown     ///     indu[p * i + 0] = (2 * i) as i32;
8779df49d7eSJed Brown     ///     indu[p * i + 1] = (2 * i + 1) as i32;
8789df49d7eSJed Brown     ///     indu[p * i + 2] = (2 * i + 2) as i32;
8799df49d7eSJed Brown     /// }
880*c68be7a2SJeremy L Thompson     /// let ru = ceed.elem_restriction(ne, p, ncomp, ndofs, ncomp * ndofs, MemType::Host, &indu)?;
8819df49d7eSJed Brown     /// let strides: [i32; 3] = [1, q as i32, q as i32];
882*c68be7a2SJeremy L Thompson     /// let rq = ceed.strided_elem_restriction(ne, q, 1, q * ne, strides)?;
8839df49d7eSJed Brown     ///
8849df49d7eSJed Brown     /// // Bases
885*c68be7a2SJeremy L Thompson     /// let bx = ceed.basis_tensor_H1_Lagrange(1, 1, 2, q, QuadMode::Gauss)?;
886*c68be7a2SJeremy L Thompson     /// let bu = ceed.basis_tensor_H1_Lagrange(1, ncomp, p, q, QuadMode::Gauss)?;
8879df49d7eSJed Brown     ///
8889df49d7eSJed Brown     /// // Build quadrature data
889*c68be7a2SJeremy L Thompson     /// let qf_build = ceed.q_function_interior_by_name("Mass1DBuild")?;
890*c68be7a2SJeremy L Thompson     /// ceed.operator(&qf_build, QFunctionOpt::None, QFunctionOpt::None)?
891*c68be7a2SJeremy L Thompson     ///     .field("dx", &rx, &bx, VectorOpt::Active)?
892*c68be7a2SJeremy L Thompson     ///     .field("weights", ElemRestrictionOpt::None, &bx, VectorOpt::None)?
893*c68be7a2SJeremy L Thompson     ///     .field("qdata", &rq, BasisOpt::Collocated, VectorOpt::Active)?
894*c68be7a2SJeremy L Thompson     ///     .apply(&x, &mut qdata)?;
8959df49d7eSJed Brown     ///
8969df49d7eSJed Brown     /// // Mass operator
8979df49d7eSJed Brown     /// let mut mass_2_comp = |[u, qdata, ..]: QFunctionInputs, [v, ..]: QFunctionOutputs| {
8989df49d7eSJed Brown     ///     // Number of quadrature points
8999df49d7eSJed Brown     ///     let q = qdata.len();
9009df49d7eSJed Brown     ///
9019df49d7eSJed Brown     ///     // Iterate over quadrature points
9029df49d7eSJed Brown     ///     for i in 0..q {
9039df49d7eSJed Brown     ///         v[i + 0 * q] = u[i + 1 * q] * qdata[i];
9049df49d7eSJed Brown     ///         v[i + 1 * q] = u[i + 0 * q] * qdata[i];
9059df49d7eSJed Brown     ///     }
9069df49d7eSJed Brown     ///
9079df49d7eSJed Brown     ///     // Return clean error code
9089df49d7eSJed Brown     ///     0
9099df49d7eSJed Brown     /// };
9109df49d7eSJed Brown     ///
9119df49d7eSJed Brown     /// let qf_mass = ceed
912*c68be7a2SJeremy L Thompson     ///     .q_function_interior(1, Box::new(mass_2_comp))?
913*c68be7a2SJeremy L Thompson     ///     .input("u", 2, EvalMode::Interp)?
914*c68be7a2SJeremy L Thompson     ///     .input("qdata", 1, EvalMode::None)?
915*c68be7a2SJeremy L Thompson     ///     .output("v", 2, EvalMode::Interp)?;
9169df49d7eSJed Brown     ///
9179df49d7eSJed Brown     /// let op_mass = ceed
918*c68be7a2SJeremy L Thompson     ///     .operator(&qf_mass, QFunctionOpt::None, QFunctionOpt::None)?
919*c68be7a2SJeremy L Thompson     ///     .field("u", &ru, &bu, VectorOpt::Active)?
920*c68be7a2SJeremy L Thompson     ///     .field("qdata", &rq, BasisOpt::Collocated, &qdata)?
921*c68be7a2SJeremy L Thompson     ///     .field("v", &ru, &bu, VectorOpt::Active)?;
9229df49d7eSJed Brown     ///
9239df49d7eSJed Brown     /// // Diagonal
924*c68be7a2SJeremy L Thompson     /// let mut diag = ceed.vector(ncomp * ncomp * ndofs)?;
9259df49d7eSJed Brown     /// diag.set_value(1.0);
926*c68be7a2SJeremy L Thompson     /// op_mass.linear_assemble_add_point_block_diagonal(&mut diag)?;
9279df49d7eSJed Brown     ///
9289df49d7eSJed Brown     /// // Manual diagonal computation
929*c68be7a2SJeremy L Thompson     /// let mut true_diag = ceed.vector(ncomp * ncomp * ndofs)?;
9309df49d7eSJed Brown     /// for i in 0..ndofs {
9319df49d7eSJed Brown     ///     for j in 0..ncomp {
9329df49d7eSJed Brown     ///         u.set_value(0.0);
9339df49d7eSJed Brown     ///         {
9349df49d7eSJed Brown     ///             let mut u_array = u.view_mut();
9359df49d7eSJed Brown     ///             u_array[i + j * ndofs] = 1.;
9369df49d7eSJed Brown     ///         }
9379df49d7eSJed Brown     ///
938*c68be7a2SJeremy L Thompson     ///         op_mass.apply(&u, &mut v)?;
9399df49d7eSJed Brown     ///
9409df49d7eSJed Brown     ///         {
9419df49d7eSJed Brown     ///             let v_array = v.view_mut();
9429df49d7eSJed Brown     ///             let mut true_array = true_diag.view_mut();
9439df49d7eSJed Brown     ///             for k in 0..ncomp {
9449df49d7eSJed Brown     ///                 true_array[i * ncomp * ncomp + k * ncomp + j] = v_array[i + k * ndofs];
9459df49d7eSJed Brown     ///             }
9469df49d7eSJed Brown     ///         }
9479df49d7eSJed Brown     ///     }
9489df49d7eSJed Brown     /// }
9499df49d7eSJed Brown     ///
9509df49d7eSJed Brown     /// // Check
9519df49d7eSJed Brown     /// diag.view()
9529df49d7eSJed Brown     ///     .iter()
9539df49d7eSJed Brown     ///     .zip(true_diag.view().iter())
9549df49d7eSJed Brown     ///     .for_each(|(computed, actual)| {
9559df49d7eSJed Brown     ///         assert!(
95680a9ef05SNatalie Beams     ///             (*computed - 1.0 - *actual).abs() < 10.0 * libceed::EPSILON,
9579df49d7eSJed Brown     ///             "Diagonal entry incorrect"
9589df49d7eSJed Brown     ///         );
9599df49d7eSJed Brown     ///     });
960*c68be7a2SJeremy L Thompson     /// # Ok(())
961*c68be7a2SJeremy L Thompson     /// # }
9629df49d7eSJed Brown     /// ```
9639df49d7eSJed Brown     pub fn linear_assemble_add_point_block_diagonal(
9649df49d7eSJed Brown         &self,
9659df49d7eSJed Brown         assembled: &mut Vector,
9669df49d7eSJed Brown     ) -> crate::Result<i32> {
9679df49d7eSJed Brown         self.op_core
9689df49d7eSJed Brown             .linear_assemble_add_point_block_diagonal(assembled)
9699df49d7eSJed Brown     }
9709df49d7eSJed Brown 
9719df49d7eSJed Brown     /// Create a multigrid coarse Operator and level transfer Operators for a
9729df49d7eSJed Brown     ///   given Operator, creating the prolongation basis from the fine and
9739df49d7eSJed Brown     ///   coarse grid interpolation
9749df49d7eSJed Brown     ///
9759df49d7eSJed Brown     /// * `p_mult_fine`  - Lvector multiplicity in parallel gather/scatter
9769df49d7eSJed Brown     /// * `rstr_coarse`  - Coarse grid restriction
9779df49d7eSJed Brown     /// * `basis_coarse` - Coarse grid active vector basis
9789df49d7eSJed Brown     ///
9799df49d7eSJed Brown     /// ```
9809df49d7eSJed Brown     /// # use libceed::prelude::*;
981*c68be7a2SJeremy L Thompson     /// # fn main() -> Result<(), libceed::CeedError> {
9829df49d7eSJed Brown     /// # let ceed = libceed::Ceed::default_init();
9839df49d7eSJed Brown     /// let ne = 15;
9849df49d7eSJed Brown     /// let p_coarse = 3;
9859df49d7eSJed Brown     /// let p_fine = 5;
9869df49d7eSJed Brown     /// let q = 6;
9879df49d7eSJed Brown     /// let ndofs_coarse = p_coarse * ne - ne + 1;
9889df49d7eSJed Brown     /// let ndofs_fine = p_fine * ne - ne + 1;
9899df49d7eSJed Brown     ///
9909df49d7eSJed Brown     /// // Vectors
9919df49d7eSJed Brown     /// let x_array = (0..ne + 1)
99280a9ef05SNatalie Beams     ///     .map(|i| 2.0 * i as Scalar / ne as Scalar - 1.0)
99380a9ef05SNatalie Beams     ///     .collect::<Vec<Scalar>>();
994*c68be7a2SJeremy L Thompson     /// let x = ceed.vector_from_slice(&x_array)?;
995*c68be7a2SJeremy L Thompson     /// let mut qdata = ceed.vector(ne * q)?;
9969df49d7eSJed Brown     /// qdata.set_value(0.0);
997*c68be7a2SJeremy L Thompson     /// let mut u_coarse = ceed.vector(ndofs_coarse)?;
9989df49d7eSJed Brown     /// u_coarse.set_value(1.0);
999*c68be7a2SJeremy L Thompson     /// let mut u_fine = ceed.vector(ndofs_fine)?;
10009df49d7eSJed Brown     /// u_fine.set_value(1.0);
1001*c68be7a2SJeremy L Thompson     /// let mut v_coarse = ceed.vector(ndofs_coarse)?;
10029df49d7eSJed Brown     /// v_coarse.set_value(0.0);
1003*c68be7a2SJeremy L Thompson     /// let mut v_fine = ceed.vector(ndofs_fine)?;
10049df49d7eSJed Brown     /// v_fine.set_value(0.0);
1005*c68be7a2SJeremy L Thompson     /// let mut multiplicity = ceed.vector(ndofs_fine)?;
10069df49d7eSJed Brown     /// multiplicity.set_value(1.0);
10079df49d7eSJed Brown     ///
10089df49d7eSJed Brown     /// // Restrictions
10099df49d7eSJed Brown     /// let strides: [i32; 3] = [1, q as i32, q as i32];
1010*c68be7a2SJeremy L Thompson     /// let rq = ceed.strided_elem_restriction(ne, q, 1, q * ne, strides)?;
10119df49d7eSJed Brown     ///
10129df49d7eSJed Brown     /// let mut indx: Vec<i32> = vec![0; 2 * ne];
10139df49d7eSJed Brown     /// for i in 0..ne {
10149df49d7eSJed Brown     ///     indx[2 * i + 0] = i as i32;
10159df49d7eSJed Brown     ///     indx[2 * i + 1] = (i + 1) as i32;
10169df49d7eSJed Brown     /// }
1017*c68be7a2SJeremy L Thompson     /// let rx = ceed.elem_restriction(ne, 2, 1, 1, ne + 1, MemType::Host, &indx)?;
10189df49d7eSJed Brown     ///
10199df49d7eSJed Brown     /// let mut indu_coarse: Vec<i32> = vec![0; p_coarse * ne];
10209df49d7eSJed Brown     /// for i in 0..ne {
10219df49d7eSJed Brown     ///     for j in 0..p_coarse {
10229df49d7eSJed Brown     ///         indu_coarse[p_coarse * i + j] = (i + j) as i32;
10239df49d7eSJed Brown     ///     }
10249df49d7eSJed Brown     /// }
1025*c68be7a2SJeremy L Thompson     /// let ru_coarse = ceed.elem_restriction(
10269df49d7eSJed Brown     ///     ne,
10279df49d7eSJed Brown     ///     p_coarse,
10289df49d7eSJed Brown     ///     1,
10299df49d7eSJed Brown     ///     1,
10309df49d7eSJed Brown     ///     ndofs_coarse,
10319df49d7eSJed Brown     ///     MemType::Host,
10329df49d7eSJed Brown     ///     &indu_coarse,
1033*c68be7a2SJeremy L Thompson     /// )?;
10349df49d7eSJed Brown     ///
10359df49d7eSJed Brown     /// let mut indu_fine: Vec<i32> = vec![0; p_fine * ne];
10369df49d7eSJed Brown     /// for i in 0..ne {
10379df49d7eSJed Brown     ///     for j in 0..p_fine {
10389df49d7eSJed Brown     ///         indu_fine[p_fine * i + j] = (i + j) as i32;
10399df49d7eSJed Brown     ///     }
10409df49d7eSJed Brown     /// }
1041*c68be7a2SJeremy L Thompson     /// let ru_fine = ceed.elem_restriction(ne, p_fine, 1, 1, ndofs_fine, MemType::Host, &indu_fine)?;
10429df49d7eSJed Brown     ///
10439df49d7eSJed Brown     /// // Bases
1044*c68be7a2SJeremy L Thompson     /// let bx = ceed.basis_tensor_H1_Lagrange(1, 1, 2, q, QuadMode::Gauss)?;
1045*c68be7a2SJeremy L Thompson     /// let bu_coarse = ceed.basis_tensor_H1_Lagrange(1, 1, p_coarse, q, QuadMode::Gauss)?;
1046*c68be7a2SJeremy L Thompson     /// let bu_fine = ceed.basis_tensor_H1_Lagrange(1, 1, p_fine, q, QuadMode::Gauss)?;
10479df49d7eSJed Brown     ///
10489df49d7eSJed Brown     /// // Build quadrature data
1049*c68be7a2SJeremy L Thompson     /// let qf_build = ceed.q_function_interior_by_name("Mass1DBuild")?;
1050*c68be7a2SJeremy L Thompson     /// ceed.operator(&qf_build, QFunctionOpt::None, QFunctionOpt::None)?
1051*c68be7a2SJeremy L Thompson     ///     .field("dx", &rx, &bx, VectorOpt::Active)?
1052*c68be7a2SJeremy L Thompson     ///     .field("weights", ElemRestrictionOpt::None, &bx, VectorOpt::None)?
1053*c68be7a2SJeremy L Thompson     ///     .field("qdata", &rq, BasisOpt::Collocated, VectorOpt::Active)?
1054*c68be7a2SJeremy L Thompson     ///     .apply(&x, &mut qdata)?;
10559df49d7eSJed Brown     ///
10569df49d7eSJed Brown     /// // Mass operator
1057*c68be7a2SJeremy L Thompson     /// let qf_mass = ceed.q_function_interior_by_name("MassApply")?;
10589df49d7eSJed Brown     /// let op_mass_fine = ceed
1059*c68be7a2SJeremy L Thompson     ///     .operator(&qf_mass, QFunctionOpt::None, QFunctionOpt::None)?
1060*c68be7a2SJeremy L Thompson     ///     .field("u", &ru_fine, &bu_fine, VectorOpt::Active)?
1061*c68be7a2SJeremy L Thompson     ///     .field("qdata", &rq, BasisOpt::Collocated, &qdata)?
1062*c68be7a2SJeremy L Thompson     ///     .field("v", &ru_fine, &bu_fine, VectorOpt::Active)?;
10639df49d7eSJed Brown     ///
10649df49d7eSJed Brown     /// // Multigrid setup
1065*c68be7a2SJeremy L Thompson     /// let (op_mass_coarse, op_prolong, op_restrict) =
1066*c68be7a2SJeremy L Thompson     ///     op_mass_fine.create_multigrid_level(&multiplicity, &ru_coarse, &bu_coarse)?;
10679df49d7eSJed Brown     ///
10689df49d7eSJed Brown     /// // Coarse problem
10699df49d7eSJed Brown     /// u_coarse.set_value(1.0);
1070*c68be7a2SJeremy L Thompson     /// op_mass_coarse.apply(&u_coarse, &mut v_coarse)?;
10719df49d7eSJed Brown     ///
10729df49d7eSJed Brown     /// // Check
107380a9ef05SNatalie Beams     /// let sum: Scalar = v_coarse.view().iter().sum();
10749df49d7eSJed Brown     /// assert!(
107580a9ef05SNatalie Beams     ///     (sum - 2.0).abs() < 50.0 * libceed::EPSILON,
10769df49d7eSJed Brown     ///     "Incorrect interval length computed"
10779df49d7eSJed Brown     /// );
10789df49d7eSJed Brown     ///
10799df49d7eSJed Brown     /// // Prolong
1080*c68be7a2SJeremy L Thompson     /// op_prolong.apply(&u_coarse, &mut u_fine)?;
10819df49d7eSJed Brown     ///
10829df49d7eSJed Brown     /// // Fine problem
1083*c68be7a2SJeremy L Thompson     /// op_mass_fine.apply(&u_fine, &mut v_fine)?;
10849df49d7eSJed Brown     ///
10859df49d7eSJed Brown     /// // Check
108680a9ef05SNatalie Beams     /// let sum: Scalar = v_fine.view().iter().sum();
10879df49d7eSJed Brown     /// assert!(
108880a9ef05SNatalie Beams     ///     (sum - 2.0).abs() < 50.0 * libceed::EPSILON,
10899df49d7eSJed Brown     ///     "Incorrect interval length computed"
10909df49d7eSJed Brown     /// );
10919df49d7eSJed Brown     ///
10929df49d7eSJed Brown     /// // Restrict
1093*c68be7a2SJeremy L Thompson     /// op_restrict.apply(&v_fine, &mut v_coarse)?;
10949df49d7eSJed Brown     ///
10959df49d7eSJed Brown     /// // Check
109680a9ef05SNatalie Beams     /// let sum: Scalar = v_coarse.view().iter().sum();
10979df49d7eSJed Brown     /// assert!(
109880a9ef05SNatalie Beams     ///     (sum - 2.0).abs() < 50.0 * libceed::EPSILON,
10999df49d7eSJed Brown     ///     "Incorrect interval length computed"
11009df49d7eSJed Brown     /// );
1101*c68be7a2SJeremy L Thompson     /// # Ok(())
1102*c68be7a2SJeremy L Thompson     /// # }
11039df49d7eSJed Brown     /// ```
11049df49d7eSJed Brown     pub fn create_multigrid_level(
11059df49d7eSJed Brown         &self,
11069df49d7eSJed Brown         p_mult_fine: &Vector,
11079df49d7eSJed Brown         rstr_coarse: &ElemRestriction,
11089df49d7eSJed Brown         basis_coarse: &Basis,
11099df49d7eSJed Brown     ) -> crate::Result<(Operator, Operator, Operator)> {
11109df49d7eSJed Brown         let mut ptr_coarse = std::ptr::null_mut();
11119df49d7eSJed Brown         let mut ptr_prolong = std::ptr::null_mut();
11129df49d7eSJed Brown         let mut ptr_restrict = std::ptr::null_mut();
11139df49d7eSJed Brown         let ierr = unsafe {
11149df49d7eSJed Brown             bind_ceed::CeedOperatorMultigridLevelCreate(
11159df49d7eSJed Brown                 self.op_core.ptr,
11169df49d7eSJed Brown                 p_mult_fine.ptr,
11179df49d7eSJed Brown                 rstr_coarse.ptr,
11189df49d7eSJed Brown                 basis_coarse.ptr,
11199df49d7eSJed Brown                 &mut ptr_coarse,
11209df49d7eSJed Brown                 &mut ptr_prolong,
11219df49d7eSJed Brown                 &mut ptr_restrict,
11229df49d7eSJed Brown             )
11239df49d7eSJed Brown         };
11249df49d7eSJed Brown         self.op_core.ceed.check_error(ierr)?;
11259df49d7eSJed Brown         let op_coarse = Operator::from_raw(self.op_core.ceed, ptr_coarse)?;
11269df49d7eSJed Brown         let op_prolong = Operator::from_raw(self.op_core.ceed, ptr_prolong)?;
11279df49d7eSJed Brown         let op_restrict = Operator::from_raw(self.op_core.ceed, ptr_restrict)?;
11289df49d7eSJed Brown         Ok((op_coarse, op_prolong, op_restrict))
11299df49d7eSJed Brown     }
11309df49d7eSJed Brown 
11319df49d7eSJed Brown     /// Create a multigrid coarse Operator and level transfer Operators for a
11329df49d7eSJed Brown     ///   given Operator with a tensor basis for the active basis
11339df49d7eSJed Brown     ///
11349df49d7eSJed Brown     /// * `p_mult_fine`   - Lvector multiplicity in parallel gather/scatter
11359df49d7eSJed Brown     /// * `rstr_coarse`   - Coarse grid restriction
11369df49d7eSJed Brown     /// * `basis_coarse`  - Coarse grid active vector basis
11379df49d7eSJed Brown     /// * `interp_c_to_f` - Matrix for coarse to fine
11389df49d7eSJed Brown     ///
11399df49d7eSJed Brown     /// ```
11409df49d7eSJed Brown     /// # use libceed::prelude::*;
1141*c68be7a2SJeremy L Thompson     /// # fn main() -> Result<(), libceed::CeedError> {
11429df49d7eSJed Brown     /// # let ceed = libceed::Ceed::default_init();
11439df49d7eSJed Brown     /// let ne = 15;
11449df49d7eSJed Brown     /// let p_coarse = 3;
11459df49d7eSJed Brown     /// let p_fine = 5;
11469df49d7eSJed Brown     /// let q = 6;
11479df49d7eSJed Brown     /// let ndofs_coarse = p_coarse * ne - ne + 1;
11489df49d7eSJed Brown     /// let ndofs_fine = p_fine * ne - ne + 1;
11499df49d7eSJed Brown     ///
11509df49d7eSJed Brown     /// // Vectors
11519df49d7eSJed Brown     /// let x_array = (0..ne + 1)
115280a9ef05SNatalie Beams     ///     .map(|i| 2.0 * i as Scalar / ne as Scalar - 1.0)
115380a9ef05SNatalie Beams     ///     .collect::<Vec<Scalar>>();
1154*c68be7a2SJeremy L Thompson     /// let x = ceed.vector_from_slice(&x_array)?;
1155*c68be7a2SJeremy L Thompson     /// let mut qdata = ceed.vector(ne * q)?;
11569df49d7eSJed Brown     /// qdata.set_value(0.0);
1157*c68be7a2SJeremy L Thompson     /// let mut u_coarse = ceed.vector(ndofs_coarse)?;
11589df49d7eSJed Brown     /// u_coarse.set_value(1.0);
1159*c68be7a2SJeremy L Thompson     /// let mut u_fine = ceed.vector(ndofs_fine)?;
11609df49d7eSJed Brown     /// u_fine.set_value(1.0);
1161*c68be7a2SJeremy L Thompson     /// let mut v_coarse = ceed.vector(ndofs_coarse)?;
11629df49d7eSJed Brown     /// v_coarse.set_value(0.0);
1163*c68be7a2SJeremy L Thompson     /// let mut v_fine = ceed.vector(ndofs_fine)?;
11649df49d7eSJed Brown     /// v_fine.set_value(0.0);
1165*c68be7a2SJeremy L Thompson     /// let mut multiplicity = ceed.vector(ndofs_fine)?;
11669df49d7eSJed Brown     /// multiplicity.set_value(1.0);
11679df49d7eSJed Brown     ///
11689df49d7eSJed Brown     /// // Restrictions
11699df49d7eSJed Brown     /// let strides: [i32; 3] = [1, q as i32, q as i32];
1170*c68be7a2SJeremy L Thompson     /// let rq = ceed.strided_elem_restriction(ne, q, 1, q * ne, strides)?;
11719df49d7eSJed Brown     ///
11729df49d7eSJed Brown     /// let mut indx: Vec<i32> = vec![0; 2 * ne];
11739df49d7eSJed Brown     /// for i in 0..ne {
11749df49d7eSJed Brown     ///     indx[2 * i + 0] = i as i32;
11759df49d7eSJed Brown     ///     indx[2 * i + 1] = (i + 1) as i32;
11769df49d7eSJed Brown     /// }
1177*c68be7a2SJeremy L Thompson     /// let rx = ceed.elem_restriction(ne, 2, 1, 1, ne + 1, MemType::Host, &indx)?;
11789df49d7eSJed Brown     ///
11799df49d7eSJed Brown     /// let mut indu_coarse: Vec<i32> = vec![0; p_coarse * ne];
11809df49d7eSJed Brown     /// for i in 0..ne {
11819df49d7eSJed Brown     ///     for j in 0..p_coarse {
11829df49d7eSJed Brown     ///         indu_coarse[p_coarse * i + j] = (i + j) as i32;
11839df49d7eSJed Brown     ///     }
11849df49d7eSJed Brown     /// }
1185*c68be7a2SJeremy L Thompson     /// let ru_coarse = ceed.elem_restriction(
11869df49d7eSJed Brown     ///     ne,
11879df49d7eSJed Brown     ///     p_coarse,
11889df49d7eSJed Brown     ///     1,
11899df49d7eSJed Brown     ///     1,
11909df49d7eSJed Brown     ///     ndofs_coarse,
11919df49d7eSJed Brown     ///     MemType::Host,
11929df49d7eSJed Brown     ///     &indu_coarse,
1193*c68be7a2SJeremy L Thompson     /// )?;
11949df49d7eSJed Brown     ///
11959df49d7eSJed Brown     /// let mut indu_fine: Vec<i32> = vec![0; p_fine * ne];
11969df49d7eSJed Brown     /// for i in 0..ne {
11979df49d7eSJed Brown     ///     for j in 0..p_fine {
11989df49d7eSJed Brown     ///         indu_fine[p_fine * i + j] = (i + j) as i32;
11999df49d7eSJed Brown     ///     }
12009df49d7eSJed Brown     /// }
1201*c68be7a2SJeremy L Thompson     /// let ru_fine = ceed.elem_restriction(ne, p_fine, 1, 1, ndofs_fine, MemType::Host, &indu_fine)?;
12029df49d7eSJed Brown     ///
12039df49d7eSJed Brown     /// // Bases
1204*c68be7a2SJeremy L Thompson     /// let bx = ceed.basis_tensor_H1_Lagrange(1, 1, 2, q, QuadMode::Gauss)?;
1205*c68be7a2SJeremy L Thompson     /// let bu_coarse = ceed.basis_tensor_H1_Lagrange(1, 1, p_coarse, q, QuadMode::Gauss)?;
1206*c68be7a2SJeremy L Thompson     /// let bu_fine = ceed.basis_tensor_H1_Lagrange(1, 1, p_fine, q, QuadMode::Gauss)?;
12079df49d7eSJed Brown     ///
12089df49d7eSJed Brown     /// // Build quadrature data
1209*c68be7a2SJeremy L Thompson     /// let qf_build = ceed.q_function_interior_by_name("Mass1DBuild")?;
1210*c68be7a2SJeremy L Thompson     /// ceed.operator(&qf_build, QFunctionOpt::None, QFunctionOpt::None)?
1211*c68be7a2SJeremy L Thompson     ///     .field("dx", &rx, &bx, VectorOpt::Active)?
1212*c68be7a2SJeremy L Thompson     ///     .field("weights", ElemRestrictionOpt::None, &bx, VectorOpt::None)?
1213*c68be7a2SJeremy L Thompson     ///     .field("qdata", &rq, BasisOpt::Collocated, VectorOpt::Active)?
1214*c68be7a2SJeremy L Thompson     ///     .apply(&x, &mut qdata)?;
12159df49d7eSJed Brown     ///
12169df49d7eSJed Brown     /// // Mass operator
1217*c68be7a2SJeremy L Thompson     /// let qf_mass = ceed.q_function_interior_by_name("MassApply")?;
12189df49d7eSJed Brown     /// let op_mass_fine = ceed
1219*c68be7a2SJeremy L Thompson     ///     .operator(&qf_mass, QFunctionOpt::None, QFunctionOpt::None)?
1220*c68be7a2SJeremy L Thompson     ///     .field("u", &ru_fine, &bu_fine, VectorOpt::Active)?
1221*c68be7a2SJeremy L Thompson     ///     .field("qdata", &rq, BasisOpt::Collocated, &qdata)?
1222*c68be7a2SJeremy L Thompson     ///     .field("v", &ru_fine, &bu_fine, VectorOpt::Active)?;
12239df49d7eSJed Brown     ///
12249df49d7eSJed Brown     /// // Multigrid setup
122580a9ef05SNatalie Beams     /// let mut interp_c_to_f: Vec<Scalar> = vec![0.; p_coarse * p_fine];
12269df49d7eSJed Brown     /// {
1227*c68be7a2SJeremy L Thompson     ///     let mut coarse = ceed.vector(p_coarse)?;
1228*c68be7a2SJeremy L Thompson     ///     let mut fine = ceed.vector(p_fine)?;
1229*c68be7a2SJeremy L Thompson     ///     let basis_c_to_f =
1230*c68be7a2SJeremy L Thompson     ///         ceed.basis_tensor_H1_Lagrange(1, 1, p_coarse, p_fine, QuadMode::GaussLobatto)?;
12319df49d7eSJed Brown     ///     for i in 0..p_coarse {
12329df49d7eSJed Brown     ///         coarse.set_value(0.0);
12339df49d7eSJed Brown     ///         {
12349df49d7eSJed Brown     ///             let mut array = coarse.view_mut();
12359df49d7eSJed Brown     ///             array[i] = 1.;
12369df49d7eSJed Brown     ///         }
1237*c68be7a2SJeremy L Thompson     ///         basis_c_to_f.apply(
12389df49d7eSJed Brown     ///             1,
12399df49d7eSJed Brown     ///             TransposeMode::NoTranspose,
12409df49d7eSJed Brown     ///             EvalMode::Interp,
12419df49d7eSJed Brown     ///             &coarse,
12429df49d7eSJed Brown     ///             &mut fine,
1243*c68be7a2SJeremy L Thompson     ///         )?;
12449df49d7eSJed Brown     ///         let array = fine.view();
12459df49d7eSJed Brown     ///         for j in 0..p_fine {
12469df49d7eSJed Brown     ///             interp_c_to_f[j * p_coarse + i] = array[j];
12479df49d7eSJed Brown     ///         }
12489df49d7eSJed Brown     ///     }
12499df49d7eSJed Brown     /// }
1250*c68be7a2SJeremy L Thompson     /// let (op_mass_coarse, op_prolong, op_restrict) = op_mass_fine.create_multigrid_level_tensor_H1(
1251*c68be7a2SJeremy L Thompson     ///     &multiplicity,
1252*c68be7a2SJeremy L Thompson     ///     &ru_coarse,
1253*c68be7a2SJeremy L Thompson     ///     &bu_coarse,
1254*c68be7a2SJeremy L Thompson     ///     &interp_c_to_f,
1255*c68be7a2SJeremy L Thompson     /// )?;
12569df49d7eSJed Brown     ///
12579df49d7eSJed Brown     /// // Coarse problem
12589df49d7eSJed Brown     /// u_coarse.set_value(1.0);
1259*c68be7a2SJeremy L Thompson     /// op_mass_coarse.apply(&u_coarse, &mut v_coarse)?;
12609df49d7eSJed Brown     ///
12619df49d7eSJed Brown     /// // Check
126280a9ef05SNatalie Beams     /// let sum: Scalar = v_coarse.view().iter().sum();
12639df49d7eSJed Brown     /// assert!(
126480a9ef05SNatalie Beams     ///     (sum - 2.0).abs() < 10.0 * libceed::EPSILON,
12659df49d7eSJed Brown     ///     "Incorrect interval length computed"
12669df49d7eSJed Brown     /// );
12679df49d7eSJed Brown     ///
12689df49d7eSJed Brown     /// // Prolong
1269*c68be7a2SJeremy L Thompson     /// op_prolong.apply(&u_coarse, &mut u_fine)?;
12709df49d7eSJed Brown     ///
12719df49d7eSJed Brown     /// // Fine problem
1272*c68be7a2SJeremy L Thompson     /// op_mass_fine.apply(&u_fine, &mut v_fine)?;
12739df49d7eSJed Brown     ///
12749df49d7eSJed Brown     /// // Check
127580a9ef05SNatalie Beams     /// let sum: Scalar = v_fine.view().iter().sum();
12769df49d7eSJed Brown     /// assert!(
127780a9ef05SNatalie Beams     ///     (sum - 2.0).abs() < 10.0 * libceed::EPSILON,
12789df49d7eSJed Brown     ///     "Incorrect interval length computed"
12799df49d7eSJed Brown     /// );
12809df49d7eSJed Brown     ///
12819df49d7eSJed Brown     /// // Restrict
1282*c68be7a2SJeremy L Thompson     /// op_restrict.apply(&v_fine, &mut v_coarse)?;
12839df49d7eSJed Brown     ///
12849df49d7eSJed Brown     /// // Check
128580a9ef05SNatalie Beams     /// let sum: Scalar = v_coarse.view().iter().sum();
12869df49d7eSJed Brown     /// assert!(
128780a9ef05SNatalie Beams     ///     (sum - 2.0).abs() < 10.0 * libceed::EPSILON,
12889df49d7eSJed Brown     ///     "Incorrect interval length computed"
12899df49d7eSJed Brown     /// );
1290*c68be7a2SJeremy L Thompson     /// # Ok(())
1291*c68be7a2SJeremy L Thompson     /// # }
12929df49d7eSJed Brown     /// ```
12939df49d7eSJed Brown     pub fn create_multigrid_level_tensor_H1(
12949df49d7eSJed Brown         &self,
12959df49d7eSJed Brown         p_mult_fine: &Vector,
12969df49d7eSJed Brown         rstr_coarse: &ElemRestriction,
12979df49d7eSJed Brown         basis_coarse: &Basis,
129880a9ef05SNatalie Beams         interpCtoF: &Vec<Scalar>,
12999df49d7eSJed Brown     ) -> crate::Result<(Operator, Operator, Operator)> {
13009df49d7eSJed Brown         let mut ptr_coarse = std::ptr::null_mut();
13019df49d7eSJed Brown         let mut ptr_prolong = std::ptr::null_mut();
13029df49d7eSJed Brown         let mut ptr_restrict = std::ptr::null_mut();
13039df49d7eSJed Brown         let ierr = unsafe {
13049df49d7eSJed Brown             bind_ceed::CeedOperatorMultigridLevelCreateTensorH1(
13059df49d7eSJed Brown                 self.op_core.ptr,
13069df49d7eSJed Brown                 p_mult_fine.ptr,
13079df49d7eSJed Brown                 rstr_coarse.ptr,
13089df49d7eSJed Brown                 basis_coarse.ptr,
13099df49d7eSJed Brown                 interpCtoF.as_ptr(),
13109df49d7eSJed Brown                 &mut ptr_coarse,
13119df49d7eSJed Brown                 &mut ptr_prolong,
13129df49d7eSJed Brown                 &mut ptr_restrict,
13139df49d7eSJed Brown             )
13149df49d7eSJed Brown         };
13159df49d7eSJed Brown         self.op_core.ceed.check_error(ierr)?;
13169df49d7eSJed Brown         let op_coarse = Operator::from_raw(self.op_core.ceed, ptr_coarse)?;
13179df49d7eSJed Brown         let op_prolong = Operator::from_raw(self.op_core.ceed, ptr_prolong)?;
13189df49d7eSJed Brown         let op_restrict = Operator::from_raw(self.op_core.ceed, ptr_restrict)?;
13199df49d7eSJed Brown         Ok((op_coarse, op_prolong, op_restrict))
13209df49d7eSJed Brown     }
13219df49d7eSJed Brown 
13229df49d7eSJed Brown     /// Create a multigrid coarse Operator and level transfer Operators for a
13239df49d7eSJed Brown     ///   given Operator with a non-tensor basis for the active basis
13249df49d7eSJed Brown     ///
13259df49d7eSJed Brown     /// * `p_mult_fine`   - Lvector multiplicity in parallel gather/scatter
13269df49d7eSJed Brown     /// * `rstr_coarse`   - Coarse grid restriction
13279df49d7eSJed Brown     /// * `basis_coarse`  - Coarse grid active vector basis
13289df49d7eSJed Brown     /// * `interp_c_to_f` - Matrix for coarse to fine
13299df49d7eSJed Brown     ///
13309df49d7eSJed Brown     /// ```
13319df49d7eSJed Brown     /// # use libceed::prelude::*;
1332*c68be7a2SJeremy L Thompson     /// # fn main() -> Result<(), libceed::CeedError> {
13339df49d7eSJed Brown     /// # let ceed = libceed::Ceed::default_init();
13349df49d7eSJed Brown     /// let ne = 15;
13359df49d7eSJed Brown     /// let p_coarse = 3;
13369df49d7eSJed Brown     /// let p_fine = 5;
13379df49d7eSJed Brown     /// let q = 6;
13389df49d7eSJed Brown     /// let ndofs_coarse = p_coarse * ne - ne + 1;
13399df49d7eSJed Brown     /// let ndofs_fine = p_fine * ne - ne + 1;
13409df49d7eSJed Brown     ///
13419df49d7eSJed Brown     /// // Vectors
13429df49d7eSJed Brown     /// let x_array = (0..ne + 1)
134380a9ef05SNatalie Beams     ///     .map(|i| 2.0 * i as Scalar / ne as Scalar - 1.0)
134480a9ef05SNatalie Beams     ///     .collect::<Vec<Scalar>>();
1345*c68be7a2SJeremy L Thompson     /// let x = ceed.vector_from_slice(&x_array)?;
1346*c68be7a2SJeremy L Thompson     /// let mut qdata = ceed.vector(ne * q)?;
13479df49d7eSJed Brown     /// qdata.set_value(0.0);
1348*c68be7a2SJeremy L Thompson     /// let mut u_coarse = ceed.vector(ndofs_coarse)?;
13499df49d7eSJed Brown     /// u_coarse.set_value(1.0);
1350*c68be7a2SJeremy L Thompson     /// let mut u_fine = ceed.vector(ndofs_fine)?;
13519df49d7eSJed Brown     /// u_fine.set_value(1.0);
1352*c68be7a2SJeremy L Thompson     /// let mut v_coarse = ceed.vector(ndofs_coarse)?;
13539df49d7eSJed Brown     /// v_coarse.set_value(0.0);
1354*c68be7a2SJeremy L Thompson     /// let mut v_fine = ceed.vector(ndofs_fine)?;
13559df49d7eSJed Brown     /// v_fine.set_value(0.0);
1356*c68be7a2SJeremy L Thompson     /// let mut multiplicity = ceed.vector(ndofs_fine)?;
13579df49d7eSJed Brown     /// multiplicity.set_value(1.0);
13589df49d7eSJed Brown     ///
13599df49d7eSJed Brown     /// // Restrictions
13609df49d7eSJed Brown     /// let strides: [i32; 3] = [1, q as i32, q as i32];
1361*c68be7a2SJeremy L Thompson     /// let rq = ceed.strided_elem_restriction(ne, q, 1, q * ne, strides)?;
13629df49d7eSJed Brown     ///
13639df49d7eSJed Brown     /// let mut indx: Vec<i32> = vec![0; 2 * ne];
13649df49d7eSJed Brown     /// for i in 0..ne {
13659df49d7eSJed Brown     ///     indx[2 * i + 0] = i as i32;
13669df49d7eSJed Brown     ///     indx[2 * i + 1] = (i + 1) as i32;
13679df49d7eSJed Brown     /// }
1368*c68be7a2SJeremy L Thompson     /// let rx = ceed.elem_restriction(ne, 2, 1, 1, ne + 1, MemType::Host, &indx)?;
13699df49d7eSJed Brown     ///
13709df49d7eSJed Brown     /// let mut indu_coarse: Vec<i32> = vec![0; p_coarse * ne];
13719df49d7eSJed Brown     /// for i in 0..ne {
13729df49d7eSJed Brown     ///     for j in 0..p_coarse {
13739df49d7eSJed Brown     ///         indu_coarse[p_coarse * i + j] = (i + j) as i32;
13749df49d7eSJed Brown     ///     }
13759df49d7eSJed Brown     /// }
1376*c68be7a2SJeremy L Thompson     /// let ru_coarse = ceed.elem_restriction(
13779df49d7eSJed Brown     ///     ne,
13789df49d7eSJed Brown     ///     p_coarse,
13799df49d7eSJed Brown     ///     1,
13809df49d7eSJed Brown     ///     1,
13819df49d7eSJed Brown     ///     ndofs_coarse,
13829df49d7eSJed Brown     ///     MemType::Host,
13839df49d7eSJed Brown     ///     &indu_coarse,
1384*c68be7a2SJeremy L Thompson     /// )?;
13859df49d7eSJed Brown     ///
13869df49d7eSJed Brown     /// let mut indu_fine: Vec<i32> = vec![0; p_fine * ne];
13879df49d7eSJed Brown     /// for i in 0..ne {
13889df49d7eSJed Brown     ///     for j in 0..p_fine {
13899df49d7eSJed Brown     ///         indu_fine[p_fine * i + j] = (i + j) as i32;
13909df49d7eSJed Brown     ///     }
13919df49d7eSJed Brown     /// }
1392*c68be7a2SJeremy L Thompson     /// let ru_fine = ceed.elem_restriction(ne, p_fine, 1, 1, ndofs_fine, MemType::Host, &indu_fine)?;
13939df49d7eSJed Brown     ///
13949df49d7eSJed Brown     /// // Bases
1395*c68be7a2SJeremy L Thompson     /// let bx = ceed.basis_tensor_H1_Lagrange(1, 1, 2, q, QuadMode::Gauss)?;
1396*c68be7a2SJeremy L Thompson     /// let bu_coarse = ceed.basis_tensor_H1_Lagrange(1, 1, p_coarse, q, QuadMode::Gauss)?;
1397*c68be7a2SJeremy L Thompson     /// let bu_fine = ceed.basis_tensor_H1_Lagrange(1, 1, p_fine, q, QuadMode::Gauss)?;
13989df49d7eSJed Brown     ///
13999df49d7eSJed Brown     /// // Build quadrature data
1400*c68be7a2SJeremy L Thompson     /// let qf_build = ceed.q_function_interior_by_name("Mass1DBuild")?;
1401*c68be7a2SJeremy L Thompson     /// ceed.operator(&qf_build, QFunctionOpt::None, QFunctionOpt::None)?
1402*c68be7a2SJeremy L Thompson     ///     .field("dx", &rx, &bx, VectorOpt::Active)?
1403*c68be7a2SJeremy L Thompson     ///     .field("weights", ElemRestrictionOpt::None, &bx, VectorOpt::None)?
1404*c68be7a2SJeremy L Thompson     ///     .field("qdata", &rq, BasisOpt::Collocated, VectorOpt::Active)?
1405*c68be7a2SJeremy L Thompson     ///     .apply(&x, &mut qdata)?;
14069df49d7eSJed Brown     ///
14079df49d7eSJed Brown     /// // Mass operator
1408*c68be7a2SJeremy L Thompson     /// let qf_mass = ceed.q_function_interior_by_name("MassApply")?;
14099df49d7eSJed Brown     /// let op_mass_fine = ceed
1410*c68be7a2SJeremy L Thompson     ///     .operator(&qf_mass, QFunctionOpt::None, QFunctionOpt::None)?
1411*c68be7a2SJeremy L Thompson     ///     .field("u", &ru_fine, &bu_fine, VectorOpt::Active)?
1412*c68be7a2SJeremy L Thompson     ///     .field("qdata", &rq, BasisOpt::Collocated, &qdata)?
1413*c68be7a2SJeremy L Thompson     ///     .field("v", &ru_fine, &bu_fine, VectorOpt::Active)?;
14149df49d7eSJed Brown     ///
14159df49d7eSJed Brown     /// // Multigrid setup
141680a9ef05SNatalie Beams     /// let mut interp_c_to_f: Vec<Scalar> = vec![0.; p_coarse * p_fine];
14179df49d7eSJed Brown     /// {
1418*c68be7a2SJeremy L Thompson     ///     let mut coarse = ceed.vector(p_coarse)?;
1419*c68be7a2SJeremy L Thompson     ///     let mut fine = ceed.vector(p_fine)?;
1420*c68be7a2SJeremy L Thompson     ///     let basis_c_to_f =
1421*c68be7a2SJeremy L Thompson     ///         ceed.basis_tensor_H1_Lagrange(1, 1, p_coarse, p_fine, QuadMode::GaussLobatto)?;
14229df49d7eSJed Brown     ///     for i in 0..p_coarse {
14239df49d7eSJed Brown     ///         coarse.set_value(0.0);
14249df49d7eSJed Brown     ///         {
14259df49d7eSJed Brown     ///             let mut array = coarse.view_mut();
14269df49d7eSJed Brown     ///             array[i] = 1.;
14279df49d7eSJed Brown     ///         }
1428*c68be7a2SJeremy L Thompson     ///         basis_c_to_f.apply(
14299df49d7eSJed Brown     ///             1,
14309df49d7eSJed Brown     ///             TransposeMode::NoTranspose,
14319df49d7eSJed Brown     ///             EvalMode::Interp,
14329df49d7eSJed Brown     ///             &coarse,
14339df49d7eSJed Brown     ///             &mut fine,
1434*c68be7a2SJeremy L Thompson     ///         )?;
14359df49d7eSJed Brown     ///         let array = fine.view();
14369df49d7eSJed Brown     ///         for j in 0..p_fine {
14379df49d7eSJed Brown     ///             interp_c_to_f[j * p_coarse + i] = array[j];
14389df49d7eSJed Brown     ///         }
14399df49d7eSJed Brown     ///     }
14409df49d7eSJed Brown     /// }
1441*c68be7a2SJeremy L Thompson     /// let (op_mass_coarse, op_prolong, op_restrict) = op_mass_fine.create_multigrid_level_H1(
1442*c68be7a2SJeremy L Thompson     ///     &multiplicity,
1443*c68be7a2SJeremy L Thompson     ///     &ru_coarse,
1444*c68be7a2SJeremy L Thompson     ///     &bu_coarse,
1445*c68be7a2SJeremy L Thompson     ///     &interp_c_to_f,
1446*c68be7a2SJeremy L Thompson     /// )?;
14479df49d7eSJed Brown     ///
14489df49d7eSJed Brown     /// // Coarse problem
14499df49d7eSJed Brown     /// u_coarse.set_value(1.0);
1450*c68be7a2SJeremy L Thompson     /// op_mass_coarse.apply(&u_coarse, &mut v_coarse)?;
14519df49d7eSJed Brown     ///
14529df49d7eSJed Brown     /// // Check
145380a9ef05SNatalie Beams     /// let sum: Scalar = v_coarse.view().iter().sum();
14549df49d7eSJed Brown     /// assert!(
145580a9ef05SNatalie Beams     ///     (sum - 2.0).abs() < 10.0 * libceed::EPSILON,
14569df49d7eSJed Brown     ///     "Incorrect interval length computed"
14579df49d7eSJed Brown     /// );
14589df49d7eSJed Brown     ///
14599df49d7eSJed Brown     /// // Prolong
1460*c68be7a2SJeremy L Thompson     /// op_prolong.apply(&u_coarse, &mut u_fine)?;
14619df49d7eSJed Brown     ///
14629df49d7eSJed Brown     /// // Fine problem
1463*c68be7a2SJeremy L Thompson     /// op_mass_fine.apply(&u_fine, &mut v_fine)?;
14649df49d7eSJed Brown     ///
14659df49d7eSJed Brown     /// // Check
146680a9ef05SNatalie Beams     /// let sum: Scalar = v_fine.view().iter().sum();
14679df49d7eSJed Brown     /// assert!(
146880a9ef05SNatalie Beams     ///     (sum - 2.0).abs() < 10.0 * libceed::EPSILON,
14699df49d7eSJed Brown     ///     "Incorrect interval length computed"
14709df49d7eSJed Brown     /// );
14719df49d7eSJed Brown     ///
14729df49d7eSJed Brown     /// // Restrict
1473*c68be7a2SJeremy L Thompson     /// op_restrict.apply(&v_fine, &mut v_coarse)?;
14749df49d7eSJed Brown     ///
14759df49d7eSJed Brown     /// // Check
147680a9ef05SNatalie Beams     /// let sum: Scalar = v_coarse.view().iter().sum();
14779df49d7eSJed Brown     /// assert!(
147880a9ef05SNatalie Beams     ///     (sum - 2.0).abs() < 10.0 * libceed::EPSILON,
14799df49d7eSJed Brown     ///     "Incorrect interval length computed"
14809df49d7eSJed Brown     /// );
1481*c68be7a2SJeremy L Thompson     /// # Ok(())
1482*c68be7a2SJeremy L Thompson     /// # }
14839df49d7eSJed Brown     /// ```
14849df49d7eSJed Brown     pub fn create_multigrid_level_H1(
14859df49d7eSJed Brown         &self,
14869df49d7eSJed Brown         p_mult_fine: &Vector,
14879df49d7eSJed Brown         rstr_coarse: &ElemRestriction,
14889df49d7eSJed Brown         basis_coarse: &Basis,
148980a9ef05SNatalie Beams         interpCtoF: &[Scalar],
14909df49d7eSJed Brown     ) -> crate::Result<(Operator, Operator, Operator)> {
14919df49d7eSJed Brown         let mut ptr_coarse = std::ptr::null_mut();
14929df49d7eSJed Brown         let mut ptr_prolong = std::ptr::null_mut();
14939df49d7eSJed Brown         let mut ptr_restrict = std::ptr::null_mut();
14949df49d7eSJed Brown         let ierr = unsafe {
14959df49d7eSJed Brown             bind_ceed::CeedOperatorMultigridLevelCreateH1(
14969df49d7eSJed Brown                 self.op_core.ptr,
14979df49d7eSJed Brown                 p_mult_fine.ptr,
14989df49d7eSJed Brown                 rstr_coarse.ptr,
14999df49d7eSJed Brown                 basis_coarse.ptr,
15009df49d7eSJed Brown                 interpCtoF.as_ptr(),
15019df49d7eSJed Brown                 &mut ptr_coarse,
15029df49d7eSJed Brown                 &mut ptr_prolong,
15039df49d7eSJed Brown                 &mut ptr_restrict,
15049df49d7eSJed Brown             )
15059df49d7eSJed Brown         };
15069df49d7eSJed Brown         self.op_core.ceed.check_error(ierr)?;
15079df49d7eSJed Brown         let op_coarse = Operator::from_raw(self.op_core.ceed, ptr_coarse)?;
15089df49d7eSJed Brown         let op_prolong = Operator::from_raw(self.op_core.ceed, ptr_prolong)?;
15099df49d7eSJed Brown         let op_restrict = Operator::from_raw(self.op_core.ceed, ptr_restrict)?;
15109df49d7eSJed Brown         Ok((op_coarse, op_prolong, op_restrict))
15119df49d7eSJed Brown     }
15129df49d7eSJed Brown }
15139df49d7eSJed Brown 
15149df49d7eSJed Brown // -----------------------------------------------------------------------------
15159df49d7eSJed Brown // Composite Operator
15169df49d7eSJed Brown // -----------------------------------------------------------------------------
15179df49d7eSJed Brown impl<'a> CompositeOperator<'a> {
15189df49d7eSJed Brown     // Constructor
15199df49d7eSJed Brown     pub fn create(ceed: &'a crate::Ceed) -> crate::Result<Self> {
15209df49d7eSJed Brown         let mut ptr = std::ptr::null_mut();
15219df49d7eSJed Brown         let ierr = unsafe { bind_ceed::CeedCompositeOperatorCreate(ceed.ptr, &mut ptr) };
15229df49d7eSJed Brown         ceed.check_error(ierr)?;
15239df49d7eSJed Brown         Ok(Self {
15249df49d7eSJed Brown             op_core: OperatorCore { ceed, ptr },
15259df49d7eSJed Brown         })
15269df49d7eSJed Brown     }
15279df49d7eSJed Brown 
15289df49d7eSJed Brown     /// Apply Operator to a vector
15299df49d7eSJed Brown     ///
15309df49d7eSJed Brown     /// * `input`  - Input Vector
15319df49d7eSJed Brown     /// * `output` - Output Vector
15329df49d7eSJed Brown     ///
15339df49d7eSJed Brown     /// ```
15349df49d7eSJed Brown     /// # use libceed::prelude::*;
1535*c68be7a2SJeremy L Thompson     /// # fn main() -> Result<(), libceed::CeedError> {
15369df49d7eSJed Brown     /// # let ceed = libceed::Ceed::default_init();
15379df49d7eSJed Brown     /// let ne = 4;
15389df49d7eSJed Brown     /// let p = 3;
15399df49d7eSJed Brown     /// let q = 4;
15409df49d7eSJed Brown     /// let ndofs = p * ne - ne + 1;
15419df49d7eSJed Brown     ///
15429df49d7eSJed Brown     /// // Vectors
1543*c68be7a2SJeremy L Thompson     /// let x = ceed.vector_from_slice(&[-1., -0.5, 0.0, 0.5, 1.0])?;
1544*c68be7a2SJeremy L Thompson     /// let mut qdata_mass = ceed.vector(ne * q)?;
15459df49d7eSJed Brown     /// qdata_mass.set_value(0.0);
1546*c68be7a2SJeremy L Thompson     /// let mut qdata_diff = ceed.vector(ne * q)?;
15479df49d7eSJed Brown     /// qdata_diff.set_value(0.0);
1548*c68be7a2SJeremy L Thompson     /// let mut u = ceed.vector(ndofs)?;
15499df49d7eSJed Brown     /// u.set_value(1.0);
1550*c68be7a2SJeremy L Thompson     /// let mut v = ceed.vector(ndofs)?;
15519df49d7eSJed Brown     /// v.set_value(0.0);
15529df49d7eSJed Brown     ///
15539df49d7eSJed Brown     /// // Restrictions
15549df49d7eSJed Brown     /// let mut indx: Vec<i32> = vec![0; 2 * ne];
15559df49d7eSJed Brown     /// for i in 0..ne {
15569df49d7eSJed Brown     ///     indx[2 * i + 0] = i as i32;
15579df49d7eSJed Brown     ///     indx[2 * i + 1] = (i + 1) as i32;
15589df49d7eSJed Brown     /// }
1559*c68be7a2SJeremy L Thompson     /// let rx = ceed.elem_restriction(ne, 2, 1, 1, ne + 1, MemType::Host, &indx)?;
15609df49d7eSJed Brown     /// let mut indu: Vec<i32> = vec![0; p * ne];
15619df49d7eSJed Brown     /// for i in 0..ne {
15629df49d7eSJed Brown     ///     indu[p * i + 0] = i as i32;
15639df49d7eSJed Brown     ///     indu[p * i + 1] = (i + 1) as i32;
15649df49d7eSJed Brown     ///     indu[p * i + 2] = (i + 2) as i32;
15659df49d7eSJed Brown     /// }
1566*c68be7a2SJeremy L Thompson     /// let ru = ceed.elem_restriction(ne, 3, 1, 1, ndofs, MemType::Host, &indu)?;
15679df49d7eSJed Brown     /// let strides: [i32; 3] = [1, q as i32, q as i32];
1568*c68be7a2SJeremy L Thompson     /// let rq = ceed.strided_elem_restriction(ne, q, 1, q * ne, strides)?;
15699df49d7eSJed Brown     ///
15709df49d7eSJed Brown     /// // Bases
1571*c68be7a2SJeremy L Thompson     /// let bx = ceed.basis_tensor_H1_Lagrange(1, 1, 2, q, QuadMode::Gauss)?;
1572*c68be7a2SJeremy L Thompson     /// let bu = ceed.basis_tensor_H1_Lagrange(1, 1, p, q, QuadMode::Gauss)?;
15739df49d7eSJed Brown     ///
15749df49d7eSJed Brown     /// // Build quadrature data
1575*c68be7a2SJeremy L Thompson     /// let qf_build_mass = ceed.q_function_interior_by_name("Mass1DBuild")?;
1576*c68be7a2SJeremy L Thompson     /// ceed.operator(&qf_build_mass, QFunctionOpt::None, QFunctionOpt::None)?
1577*c68be7a2SJeremy L Thompson     ///     .field("dx", &rx, &bx, VectorOpt::Active)?
1578*c68be7a2SJeremy L Thompson     ///     .field("weights", ElemRestrictionOpt::None, &bx, VectorOpt::None)?
1579*c68be7a2SJeremy L Thompson     ///     .field("qdata", &rq, BasisOpt::Collocated, VectorOpt::Active)?
1580*c68be7a2SJeremy L Thompson     ///     .apply(&x, &mut qdata_mass)?;
15819df49d7eSJed Brown     ///
1582*c68be7a2SJeremy L Thompson     /// let qf_build_diff = ceed.q_function_interior_by_name("Poisson1DBuild")?;
1583*c68be7a2SJeremy L Thompson     /// ceed.operator(&qf_build_diff, QFunctionOpt::None, QFunctionOpt::None)?
1584*c68be7a2SJeremy L Thompson     ///     .field("dx", &rx, &bx, VectorOpt::Active)?
1585*c68be7a2SJeremy L Thompson     ///     .field("weights", ElemRestrictionOpt::None, &bx, VectorOpt::None)?
1586*c68be7a2SJeremy L Thompson     ///     .field("qdata", &rq, BasisOpt::Collocated, VectorOpt::Active)?
1587*c68be7a2SJeremy L Thompson     ///     .apply(&x, &mut qdata_diff)?;
15889df49d7eSJed Brown     ///
15899df49d7eSJed Brown     /// // Application operator
1590*c68be7a2SJeremy L Thompson     /// let qf_mass = ceed.q_function_interior_by_name("MassApply")?;
15919df49d7eSJed Brown     /// let op_mass = ceed
1592*c68be7a2SJeremy L Thompson     ///     .operator(&qf_mass, QFunctionOpt::None, QFunctionOpt::None)?
1593*c68be7a2SJeremy L Thompson     ///     .field("u", &ru, &bu, VectorOpt::Active)?
1594*c68be7a2SJeremy L Thompson     ///     .field("qdata", &rq, BasisOpt::Collocated, &qdata_mass)?
1595*c68be7a2SJeremy L Thompson     ///     .field("v", &ru, &bu, VectorOpt::Active)?;
15969df49d7eSJed Brown     ///
1597*c68be7a2SJeremy L Thompson     /// let qf_diff = ceed.q_function_interior_by_name("Poisson1DApply")?;
15989df49d7eSJed Brown     /// let op_diff = ceed
1599*c68be7a2SJeremy L Thompson     ///     .operator(&qf_diff, QFunctionOpt::None, QFunctionOpt::None)?
1600*c68be7a2SJeremy L Thompson     ///     .field("du", &ru, &bu, VectorOpt::Active)?
1601*c68be7a2SJeremy L Thompson     ///     .field("qdata", &rq, BasisOpt::Collocated, &qdata_diff)?
1602*c68be7a2SJeremy L Thompson     ///     .field("dv", &ru, &bu, VectorOpt::Active)?;
16039df49d7eSJed Brown     ///
16049df49d7eSJed Brown     /// let op_composite = ceed
1605*c68be7a2SJeremy L Thompson     ///     .composite_operator()?
1606*c68be7a2SJeremy L Thompson     ///     .sub_operator(&op_mass)?
1607*c68be7a2SJeremy L Thompson     ///     .sub_operator(&op_diff)?;
16089df49d7eSJed Brown     ///
16099df49d7eSJed Brown     /// v.set_value(0.0);
1610*c68be7a2SJeremy L Thompson     /// op_composite.apply(&u, &mut v)?;
16119df49d7eSJed Brown     ///
16129df49d7eSJed Brown     /// // Check
161380a9ef05SNatalie Beams     /// let sum: Scalar = v.view().iter().sum();
16149df49d7eSJed Brown     /// assert!(
161580a9ef05SNatalie Beams     ///     (sum - 2.0).abs() < 10.0 * libceed::EPSILON,
16169df49d7eSJed Brown     ///     "Incorrect interval length computed"
16179df49d7eSJed Brown     /// );
1618*c68be7a2SJeremy L Thompson     /// # Ok(())
1619*c68be7a2SJeremy L Thompson     /// # }
16209df49d7eSJed Brown     /// ```
16219df49d7eSJed Brown     pub fn apply(&self, input: &Vector, output: &mut Vector) -> crate::Result<i32> {
16229df49d7eSJed Brown         self.op_core.apply(input, output)
16239df49d7eSJed Brown     }
16249df49d7eSJed Brown 
16259df49d7eSJed Brown     /// Apply Operator to a vector and add result to output vector
16269df49d7eSJed Brown     ///
16279df49d7eSJed Brown     /// * `input`  - Input Vector
16289df49d7eSJed Brown     /// * `output` - Output Vector
16299df49d7eSJed Brown     ///
16309df49d7eSJed Brown     /// ```
16319df49d7eSJed Brown     /// # use libceed::prelude::*;
1632*c68be7a2SJeremy L Thompson     /// # fn main() -> Result<(), libceed::CeedError> {
16339df49d7eSJed Brown     /// # let ceed = libceed::Ceed::default_init();
16349df49d7eSJed Brown     /// let ne = 4;
16359df49d7eSJed Brown     /// let p = 3;
16369df49d7eSJed Brown     /// let q = 4;
16379df49d7eSJed Brown     /// let ndofs = p * ne - ne + 1;
16389df49d7eSJed Brown     ///
16399df49d7eSJed Brown     /// // Vectors
1640*c68be7a2SJeremy L Thompson     /// let x = ceed.vector_from_slice(&[-1., -0.5, 0.0, 0.5, 1.0])?;
1641*c68be7a2SJeremy L Thompson     /// let mut qdata_mass = ceed.vector(ne * q)?;
16429df49d7eSJed Brown     /// qdata_mass.set_value(0.0);
1643*c68be7a2SJeremy L Thompson     /// let mut qdata_diff = ceed.vector(ne * q)?;
16449df49d7eSJed Brown     /// qdata_diff.set_value(0.0);
1645*c68be7a2SJeremy L Thompson     /// let mut u = ceed.vector(ndofs)?;
16469df49d7eSJed Brown     /// u.set_value(1.0);
1647*c68be7a2SJeremy L Thompson     /// let mut v = ceed.vector(ndofs)?;
16489df49d7eSJed Brown     /// v.set_value(0.0);
16499df49d7eSJed Brown     ///
16509df49d7eSJed Brown     /// // Restrictions
16519df49d7eSJed Brown     /// let mut indx: Vec<i32> = vec![0; 2 * ne];
16529df49d7eSJed Brown     /// for i in 0..ne {
16539df49d7eSJed Brown     ///     indx[2 * i + 0] = i as i32;
16549df49d7eSJed Brown     ///     indx[2 * i + 1] = (i + 1) as i32;
16559df49d7eSJed Brown     /// }
1656*c68be7a2SJeremy L Thompson     /// let rx = ceed.elem_restriction(ne, 2, 1, 1, ne + 1, MemType::Host, &indx)?;
16579df49d7eSJed Brown     /// let mut indu: Vec<i32> = vec![0; p * ne];
16589df49d7eSJed Brown     /// for i in 0..ne {
16599df49d7eSJed Brown     ///     indu[p * i + 0] = i as i32;
16609df49d7eSJed Brown     ///     indu[p * i + 1] = (i + 1) as i32;
16619df49d7eSJed Brown     ///     indu[p * i + 2] = (i + 2) as i32;
16629df49d7eSJed Brown     /// }
1663*c68be7a2SJeremy L Thompson     /// let ru = ceed.elem_restriction(ne, 3, 1, 1, ndofs, MemType::Host, &indu)?;
16649df49d7eSJed Brown     /// let strides: [i32; 3] = [1, q as i32, q as i32];
1665*c68be7a2SJeremy L Thompson     /// let rq = ceed.strided_elem_restriction(ne, q, 1, q * ne, strides)?;
16669df49d7eSJed Brown     ///
16679df49d7eSJed Brown     /// // Bases
1668*c68be7a2SJeremy L Thompson     /// let bx = ceed.basis_tensor_H1_Lagrange(1, 1, 2, q, QuadMode::Gauss)?;
1669*c68be7a2SJeremy L Thompson     /// let bu = ceed.basis_tensor_H1_Lagrange(1, 1, p, q, QuadMode::Gauss)?;
16709df49d7eSJed Brown     ///
16719df49d7eSJed Brown     /// // Build quadrature data
1672*c68be7a2SJeremy L Thompson     /// let qf_build_mass = ceed.q_function_interior_by_name("Mass1DBuild")?;
1673*c68be7a2SJeremy L Thompson     /// ceed.operator(&qf_build_mass, QFunctionOpt::None, QFunctionOpt::None)?
1674*c68be7a2SJeremy L Thompson     ///     .field("dx", &rx, &bx, VectorOpt::Active)?
1675*c68be7a2SJeremy L Thompson     ///     .field("weights", ElemRestrictionOpt::None, &bx, VectorOpt::None)?
1676*c68be7a2SJeremy L Thompson     ///     .field("qdata", &rq, BasisOpt::Collocated, VectorOpt::Active)?
1677*c68be7a2SJeremy L Thompson     ///     .apply(&x, &mut qdata_mass)?;
16789df49d7eSJed Brown     ///
1679*c68be7a2SJeremy L Thompson     /// let qf_build_diff = ceed.q_function_interior_by_name("Poisson1DBuild")?;
1680*c68be7a2SJeremy L Thompson     /// ceed.operator(&qf_build_diff, QFunctionOpt::None, QFunctionOpt::None)?
1681*c68be7a2SJeremy L Thompson     ///     .field("dx", &rx, &bx, VectorOpt::Active)?
1682*c68be7a2SJeremy L Thompson     ///     .field("weights", ElemRestrictionOpt::None, &bx, VectorOpt::None)?
1683*c68be7a2SJeremy L Thompson     ///     .field("qdata", &rq, BasisOpt::Collocated, VectorOpt::Active)?
1684*c68be7a2SJeremy L Thompson     ///     .apply(&x, &mut qdata_diff)?;
16859df49d7eSJed Brown     ///
16869df49d7eSJed Brown     /// // Application operator
1687*c68be7a2SJeremy L Thompson     /// let qf_mass = ceed.q_function_interior_by_name("MassApply")?;
16889df49d7eSJed Brown     /// let op_mass = ceed
1689*c68be7a2SJeremy L Thompson     ///     .operator(&qf_mass, QFunctionOpt::None, QFunctionOpt::None)?
1690*c68be7a2SJeremy L Thompson     ///     .field("u", &ru, &bu, VectorOpt::Active)?
1691*c68be7a2SJeremy L Thompson     ///     .field("qdata", &rq, BasisOpt::Collocated, &qdata_mass)?
1692*c68be7a2SJeremy L Thompson     ///     .field("v", &ru, &bu, VectorOpt::Active)?;
16939df49d7eSJed Brown     ///
1694*c68be7a2SJeremy L Thompson     /// let qf_diff = ceed.q_function_interior_by_name("Poisson1DApply")?;
16959df49d7eSJed Brown     /// let op_diff = ceed
1696*c68be7a2SJeremy L Thompson     ///     .operator(&qf_diff, QFunctionOpt::None, QFunctionOpt::None)?
1697*c68be7a2SJeremy L Thompson     ///     .field("du", &ru, &bu, VectorOpt::Active)?
1698*c68be7a2SJeremy L Thompson     ///     .field("qdata", &rq, BasisOpt::Collocated, &qdata_diff)?
1699*c68be7a2SJeremy L Thompson     ///     .field("dv", &ru, &bu, VectorOpt::Active)?;
17009df49d7eSJed Brown     ///
17019df49d7eSJed Brown     /// let op_composite = ceed
1702*c68be7a2SJeremy L Thompson     ///     .composite_operator()?
1703*c68be7a2SJeremy L Thompson     ///     .sub_operator(&op_mass)?
1704*c68be7a2SJeremy L Thompson     ///     .sub_operator(&op_diff)?;
17059df49d7eSJed Brown     ///
17069df49d7eSJed Brown     /// v.set_value(1.0);
1707*c68be7a2SJeremy L Thompson     /// op_composite.apply_add(&u, &mut v)?;
17089df49d7eSJed Brown     ///
17099df49d7eSJed Brown     /// // Check
171080a9ef05SNatalie Beams     /// let sum: Scalar = v.view().iter().sum();
17119df49d7eSJed Brown     /// assert!(
171280a9ef05SNatalie Beams     ///     (sum - (2.0 + ndofs as Scalar)).abs() < 10.0 * libceed::EPSILON,
17139df49d7eSJed Brown     ///     "Incorrect interval length computed"
17149df49d7eSJed Brown     /// );
1715*c68be7a2SJeremy L Thompson     /// # Ok(())
1716*c68be7a2SJeremy L Thompson     /// # }
17179df49d7eSJed Brown     /// ```
17189df49d7eSJed Brown     pub fn apply_add(&self, input: &Vector, output: &mut Vector) -> crate::Result<i32> {
17199df49d7eSJed Brown         self.op_core.apply_add(input, output)
17209df49d7eSJed Brown     }
17219df49d7eSJed Brown 
17229df49d7eSJed Brown     /// Add a sub-Operator to a Composite Operator
17239df49d7eSJed Brown     ///
17249df49d7eSJed Brown     /// * `subop` - Sub-Operator
17259df49d7eSJed Brown     ///
17269df49d7eSJed Brown     /// ```
17279df49d7eSJed Brown     /// # use libceed::prelude::*;
1728*c68be7a2SJeremy L Thompson     /// # fn main() -> Result<(), libceed::CeedError> {
17299df49d7eSJed Brown     /// # let ceed = libceed::Ceed::default_init();
1730*c68be7a2SJeremy L Thompson     /// let mut op = ceed.composite_operator()?;
17319df49d7eSJed Brown     ///
1732*c68be7a2SJeremy L Thompson     /// let qf_mass = ceed.q_function_interior_by_name("MassApply")?;
1733*c68be7a2SJeremy L Thompson     /// let op_mass = ceed.operator(&qf_mass, QFunctionOpt::None, QFunctionOpt::None)?;
1734*c68be7a2SJeremy L Thompson     /// op = op.sub_operator(&op_mass)?;
17359df49d7eSJed Brown     ///
1736*c68be7a2SJeremy L Thompson     /// let qf_diff = ceed.q_function_interior_by_name("Poisson1DApply")?;
1737*c68be7a2SJeremy L Thompson     /// let op_diff = ceed.operator(&qf_diff, QFunctionOpt::None, QFunctionOpt::None)?;
1738*c68be7a2SJeremy L Thompson     /// op = op.sub_operator(&op_diff)?;
1739*c68be7a2SJeremy L Thompson     /// # Ok(())
1740*c68be7a2SJeremy L Thompson     /// # }
17419df49d7eSJed Brown     /// ```
17429df49d7eSJed Brown     #[allow(unused_mut)]
17439df49d7eSJed Brown     pub fn sub_operator(mut self, subop: &Operator) -> crate::Result<Self> {
17449df49d7eSJed Brown         let ierr =
17459df49d7eSJed Brown             unsafe { bind_ceed::CeedCompositeOperatorAddSub(self.op_core.ptr, subop.op_core.ptr) };
17469df49d7eSJed Brown         self.op_core.ceed.check_error(ierr)?;
17479df49d7eSJed Brown         Ok(self)
17489df49d7eSJed Brown     }
17499df49d7eSJed Brown 
17509df49d7eSJed Brown     /// Assemble the diagonal of a square linear Operator
17519df49d7eSJed Brown     ///
17529df49d7eSJed Brown     /// This overwrites a Vector with the diagonal of a linear Operator.
17539df49d7eSJed Brown     ///
17549df49d7eSJed Brown     /// Note: Currently only non-composite Operators with a single field and
17559df49d7eSJed Brown     ///      composite Operators with single field sub-operators are supported.
17569df49d7eSJed Brown     ///
17579df49d7eSJed Brown     /// * `op`        - Operator to assemble QFunction
17589df49d7eSJed Brown     /// * `assembled` - Vector to store assembled Operator diagonal
17599df49d7eSJed Brown     pub fn linear_asssemble_diagonal(&self, assembled: &mut Vector) -> crate::Result<i32> {
17609df49d7eSJed Brown         self.op_core.linear_assemble_add_diagonal(assembled)
17619df49d7eSJed Brown     }
17629df49d7eSJed Brown 
17639df49d7eSJed Brown     /// Assemble the point block diagonal of a square linear Operator
17649df49d7eSJed Brown     ///
17659df49d7eSJed Brown     /// This overwrites a Vector with the point block diagonal of a linear
17669df49d7eSJed Brown     ///   Operator.
17679df49d7eSJed Brown     ///
17689df49d7eSJed Brown     /// Note: Currently only non-composite Operators with a single field and
17699df49d7eSJed Brown     ///      composite Operators with single field sub-operators are supported.
17709df49d7eSJed Brown     ///
17719df49d7eSJed Brown     /// * `op`        - Operator to assemble QFunction
17729df49d7eSJed Brown     /// * `assembled` - Vector to store assembled Operator diagonal
17739df49d7eSJed Brown     pub fn linear_assemble_add_diagonal(&self, assembled: &mut Vector) -> crate::Result<i32> {
17749df49d7eSJed Brown         self.op_core.linear_assemble_add_diagonal(assembled)
17759df49d7eSJed Brown     }
17769df49d7eSJed Brown 
17779df49d7eSJed Brown     /// Assemble the diagonal of a square linear Operator
17789df49d7eSJed Brown     ///
17799df49d7eSJed Brown     /// This overwrites a Vector with the diagonal of a linear Operator.
17809df49d7eSJed Brown     ///
17819df49d7eSJed Brown     /// Note: Currently only non-composite Operators with a single field and
17829df49d7eSJed Brown     ///      composite Operators with single field sub-operators are supported.
17839df49d7eSJed Brown     ///
17849df49d7eSJed Brown     /// * `op`        - Operator to assemble QFunction
17859df49d7eSJed Brown     /// * `assembled` - Vector to store assembled CeedOperator point block
17869df49d7eSJed Brown     ///                   diagonal, provided in row-major form with an
17879df49d7eSJed Brown     ///                   `ncomp * ncomp` block at each node. The dimensions of
17889df49d7eSJed Brown     ///                   this vector are derived from the active vector for
17899df49d7eSJed Brown     ///                   the CeedOperator. The array has shape
17909df49d7eSJed Brown     ///                   `[nodes, component out, component in]`.
17919df49d7eSJed Brown     pub fn linear_assemble_point_block_diagonal(
17929df49d7eSJed Brown         &self,
17939df49d7eSJed Brown         assembled: &mut Vector,
17949df49d7eSJed Brown     ) -> crate::Result<i32> {
17959df49d7eSJed Brown         self.op_core.linear_assemble_add_diagonal(assembled)
17969df49d7eSJed Brown     }
17979df49d7eSJed Brown 
17989df49d7eSJed Brown     /// Assemble the diagonal of a square linear Operator
17999df49d7eSJed Brown     ///
18009df49d7eSJed Brown     /// This sums into a Vector with the diagonal of a linear Operator.
18019df49d7eSJed Brown     ///
18029df49d7eSJed Brown     /// Note: Currently only non-composite Operators with a single field and
18039df49d7eSJed Brown     ///      composite Operators with single field sub-operators are supported.
18049df49d7eSJed Brown     ///
18059df49d7eSJed Brown     /// * `op`        - Operator to assemble QFunction
18069df49d7eSJed Brown     /// * `assembled` - Vector to store assembled CeedOperator point block
18079df49d7eSJed Brown     ///                   diagonal, provided in row-major form with an
18089df49d7eSJed Brown     ///                   `ncomp * ncomp` block at each node. The dimensions of
18099df49d7eSJed Brown     ///                   this vector are derived from the active vector for
18109df49d7eSJed Brown     ///                   the CeedOperator. The array has shape
18119df49d7eSJed Brown     ///                   `[nodes, component out, component in]`.
18129df49d7eSJed Brown     pub fn linear_assemble_add_point_block_diagonal(
18139df49d7eSJed Brown         &self,
18149df49d7eSJed Brown         assembled: &mut Vector,
18159df49d7eSJed Brown     ) -> crate::Result<i32> {
18169df49d7eSJed Brown         self.op_core.linear_assemble_add_diagonal(assembled)
18179df49d7eSJed Brown     }
18189df49d7eSJed Brown }
18199df49d7eSJed Brown 
18209df49d7eSJed Brown // -----------------------------------------------------------------------------
1821