xref: /libCEED/rust/libceed/src/operator.rs (revision 17b5e52fb5f6c474a1e9fc79dd3f8a824f1b6293)
1 // Copyright (c) 2017-2022, Lawrence Livermore National Security, LLC and other CEED contributors.
2 // All Rights Reserved. See the top-level LICENSE and NOTICE files for details.
3 //
4 // SPDX-License-Identifier: BSD-2-Clause
5 //
6 // This file is part of CEED:  http://github.com/ceed
7 
8 //! A Ceed Operator defines the finite/spectral element operator associated to a
9 //! Ceed QFunction. A Ceed Operator connects Ceed ElemRestrictions,
10 //! Ceed Bases, and Ceed QFunctions.
11 
12 use crate::prelude::*;
13 
14 // -----------------------------------------------------------------------------
15 // Operator Field context wrapper
16 // -----------------------------------------------------------------------------
17 #[derive(Debug)]
18 pub struct OperatorField<'a> {
19     pub(crate) ptr: bind_ceed::CeedOperatorField,
20     _lifeline: PhantomData<&'a ()>,
21 }
22 
23 // -----------------------------------------------------------------------------
24 // Implementations
25 // -----------------------------------------------------------------------------
26 impl<'a> OperatorField<'a> {
27     /// Get the name of an OperatorField
28     ///
29     /// ```
30     /// # use libceed::prelude::*;
31     /// # fn main() -> libceed::Result<()> {
32     /// # let ceed = libceed::Ceed::default_init();
33     /// let qf = ceed.q_function_interior_by_name("Mass1DBuild")?;
34     ///
35     /// // Operator field arguments
36     /// let ne = 3;
37     /// let q = 4 as usize;
38     /// let mut ind: Vec<i32> = vec![0; 2 * ne];
39     /// for i in 0..ne {
40     ///     ind[2 * i + 0] = i as i32;
41     ///     ind[2 * i + 1] = (i + 1) as i32;
42     /// }
43     /// let r = ceed.elem_restriction(ne, 2, 1, 1, ne + 1, MemType::Host, &ind)?;
44     /// let strides: [i32; 3] = [1, q as i32, q as i32];
45     /// let rq = ceed.strided_elem_restriction(ne, q, 1, q * ne, strides)?;
46     ///
47     /// let b = ceed.basis_tensor_H1_Lagrange(1, 1, 2, q, QuadMode::Gauss)?;
48     ///
49     /// // Operator fields
50     /// let op = ceed
51     ///     .operator(&qf, QFunctionOpt::None, QFunctionOpt::None)?
52     ///     .field("dx", &r, &b, VectorOpt::Active)?
53     ///     .field("weights", ElemRestrictionOpt::None, &b, VectorOpt::None)?
54     ///     .field("qdata", &rq, BasisOpt::None, VectorOpt::Active)?;
55     ///
56     /// let inputs = op.inputs()?;
57     ///
58     /// assert_eq!(inputs[0].name(), "dx", "Incorrect input name");
59     /// assert_eq!(inputs[1].name(), "weights", "Incorrect input name");
60     /// # Ok(())
61     /// # }
62     /// ```
63     pub fn name(&self) -> &str {
64         let mut name_ptr: *mut std::os::raw::c_char = std::ptr::null_mut();
65         unsafe {
66             bind_ceed::CeedOperatorFieldGetName(self.ptr, &mut name_ptr);
67         }
68         unsafe { CStr::from_ptr(name_ptr) }.to_str().unwrap()
69     }
70 
71     /// Get the ElemRestriction of an OperatorField
72     ///
73     /// ```
74     /// # use libceed::prelude::*;
75     /// # fn main() -> libceed::Result<()> {
76     /// # let ceed = libceed::Ceed::default_init();
77     /// let qf = ceed.q_function_interior_by_name("Mass1DBuild")?;
78     ///
79     /// // Operator field arguments
80     /// let ne = 3;
81     /// let q = 4 as usize;
82     /// let mut ind: Vec<i32> = vec![0; 2 * ne];
83     /// for i in 0..ne {
84     ///     ind[2 * i + 0] = i as i32;
85     ///     ind[2 * i + 1] = (i + 1) as i32;
86     /// }
87     /// let r = ceed.elem_restriction(ne, 2, 1, 1, ne + 1, MemType::Host, &ind)?;
88     /// let strides: [i32; 3] = [1, q as i32, q as i32];
89     /// let rq = ceed.strided_elem_restriction(ne, q, 1, q * ne, strides)?;
90     ///
91     /// let b = ceed.basis_tensor_H1_Lagrange(1, 1, 2, q, QuadMode::Gauss)?;
92     ///
93     /// // Operator fields
94     /// let op = ceed
95     ///     .operator(&qf, QFunctionOpt::None, QFunctionOpt::None)?
96     ///     .field("dx", &r, &b, VectorOpt::Active)?
97     ///     .field("weights", ElemRestrictionOpt::None, &b, VectorOpt::None)?
98     ///     .field("qdata", &rq, BasisOpt::None, VectorOpt::Active)?;
99     ///
100     /// let inputs = op.inputs()?;
101     ///
102     /// assert!(
103     ///     inputs[0].elem_restriction().is_some(),
104     ///     "Incorrect field ElemRestriction"
105     /// );
106     /// if let ElemRestrictionOpt::Some(r) = inputs[0].elem_restriction() {
107     ///     assert_eq!(
108     ///         r.num_elements(),
109     ///         ne,
110     ///         "Incorrect field ElemRestriction number of elements"
111     ///     );
112     /// }
113     ///
114     /// assert!(
115     ///     inputs[1].elem_restriction().is_none(),
116     ///     "Incorrect field ElemRestriction"
117     /// );
118     /// # Ok(())
119     /// # }
120     /// ```
121     pub fn elem_restriction(&self) -> ElemRestrictionOpt {
122         let mut ptr = std::ptr::null_mut();
123         unsafe {
124             bind_ceed::CeedOperatorFieldGetElemRestriction(self.ptr, &mut ptr);
125         }
126         if ptr == unsafe { bind_ceed::CEED_ELEMRESTRICTION_NONE } {
127             ElemRestrictionOpt::None
128         } else {
129             let slice = unsafe {
130                 std::slice::from_raw_parts(
131                     &ptr as *const bind_ceed::CeedElemRestriction as *const crate::ElemRestriction,
132                     1 as usize,
133                 )
134             };
135             ElemRestrictionOpt::Some(&slice[0])
136         }
137     }
138 
139     /// Get the Basis of an OperatorField
140     ///
141     /// ```
142     /// # use libceed::prelude::*;
143     /// # fn main() -> libceed::Result<()> {
144     /// # let ceed = libceed::Ceed::default_init();
145     /// let qf = ceed.q_function_interior_by_name("Mass1DBuild")?;
146     ///
147     /// // Operator field arguments
148     /// let ne = 3;
149     /// let q = 4 as usize;
150     /// let mut ind: Vec<i32> = vec![0; 2 * ne];
151     /// for i in 0..ne {
152     ///     ind[2 * i + 0] = i as i32;
153     ///     ind[2 * i + 1] = (i + 1) as i32;
154     /// }
155     /// let r = ceed.elem_restriction(ne, 2, 1, 1, ne + 1, MemType::Host, &ind)?;
156     /// let strides: [i32; 3] = [1, q as i32, q as i32];
157     /// let rq = ceed.strided_elem_restriction(ne, q, 1, q * ne, strides)?;
158     ///
159     /// let b = ceed.basis_tensor_H1_Lagrange(1, 1, 2, q, QuadMode::Gauss)?;
160     ///
161     /// // Operator fields
162     /// let op = ceed
163     ///     .operator(&qf, QFunctionOpt::None, QFunctionOpt::None)?
164     ///     .field("dx", &r, &b, VectorOpt::Active)?
165     ///     .field("weights", ElemRestrictionOpt::None, &b, VectorOpt::None)?
166     ///     .field("qdata", &rq, BasisOpt::None, VectorOpt::Active)?;
167     ///
168     /// let inputs = op.inputs()?;
169     ///
170     /// assert!(inputs[0].basis().is_some(), "Incorrect field Basis");
171     /// if let BasisOpt::Some(b) = inputs[0].basis() {
172     ///     assert_eq!(
173     ///         b.num_quadrature_points(),
174     ///         q,
175     ///         "Incorrect field Basis number of quadrature points"
176     ///     );
177     /// }
178     /// assert!(inputs[1].basis().is_some(), "Incorrect field Basis");
179     /// if let BasisOpt::Some(b) = inputs[0].basis() {
180     ///     assert_eq!(
181     ///         b.num_quadrature_points(),
182     ///         q,
183     ///         "Incorrect field Basis number of quadrature points"
184     ///     );
185     /// }
186     ///
187     /// let outputs = op.outputs()?;
188     ///
189     /// assert!(outputs[0].basis().is_none(), "Incorrect field Basis");
190     /// # Ok(())
191     /// # }
192     /// ```
193     pub fn basis(&self) -> BasisOpt {
194         let mut ptr = std::ptr::null_mut();
195         unsafe {
196             bind_ceed::CeedOperatorFieldGetBasis(self.ptr, &mut ptr);
197         }
198         if ptr == unsafe { bind_ceed::CEED_BASIS_NONE } {
199             BasisOpt::None
200         } else {
201             let slice = unsafe {
202                 std::slice::from_raw_parts(
203                     &ptr as *const bind_ceed::CeedBasis as *const crate::Basis,
204                     1 as usize,
205                 )
206             };
207             BasisOpt::Some(&slice[0])
208         }
209     }
210 
211     /// Get the Vector of an OperatorField
212     ///
213     /// ```
214     /// # use libceed::prelude::*;
215     /// # fn main() -> libceed::Result<()> {
216     /// # let ceed = libceed::Ceed::default_init();
217     /// let qf = ceed.q_function_interior_by_name("Mass1DBuild")?;
218     ///
219     /// // Operator field arguments
220     /// let ne = 3;
221     /// let q = 4 as usize;
222     /// let mut ind: Vec<i32> = vec![0; 2 * ne];
223     /// for i in 0..ne {
224     ///     ind[2 * i + 0] = i as i32;
225     ///     ind[2 * i + 1] = (i + 1) as i32;
226     /// }
227     /// let r = ceed.elem_restriction(ne, 2, 1, 1, ne + 1, MemType::Host, &ind)?;
228     /// let strides: [i32; 3] = [1, q as i32, q as i32];
229     /// let rq = ceed.strided_elem_restriction(ne, q, 1, q * ne, strides)?;
230     ///
231     /// let b = ceed.basis_tensor_H1_Lagrange(1, 1, 2, q, QuadMode::Gauss)?;
232     ///
233     /// // Operator fields
234     /// let op = ceed
235     ///     .operator(&qf, QFunctionOpt::None, QFunctionOpt::None)?
236     ///     .field("dx", &r, &b, VectorOpt::Active)?
237     ///     .field("weights", ElemRestrictionOpt::None, &b, VectorOpt::None)?
238     ///     .field("qdata", &rq, BasisOpt::None, VectorOpt::Active)?;
239     ///
240     /// let inputs = op.inputs()?;
241     ///
242     /// assert!(inputs[0].vector().is_active(), "Incorrect field Vector");
243     /// assert!(inputs[1].vector().is_none(), "Incorrect field Vector");
244     /// # Ok(())
245     /// # }
246     /// ```
247     pub fn vector(&self) -> VectorOpt {
248         let mut ptr = std::ptr::null_mut();
249         unsafe {
250             bind_ceed::CeedOperatorFieldGetVector(self.ptr, &mut ptr);
251         }
252         if ptr == unsafe { bind_ceed::CEED_VECTOR_ACTIVE } {
253             VectorOpt::Active
254         } else if ptr == unsafe { bind_ceed::CEED_VECTOR_NONE } {
255             VectorOpt::None
256         } else {
257             let slice = unsafe {
258                 std::slice::from_raw_parts(
259                     &ptr as *const bind_ceed::CeedVector as *const crate::Vector,
260                     1 as usize,
261                 )
262             };
263             VectorOpt::Some(&slice[0])
264         }
265     }
266 }
267 
268 // -----------------------------------------------------------------------------
269 // Operator context wrapper
270 // -----------------------------------------------------------------------------
271 #[derive(Debug)]
272 pub(crate) struct OperatorCore<'a> {
273     ptr: bind_ceed::CeedOperator,
274     _lifeline: PhantomData<&'a ()>,
275 }
276 
277 #[derive(Debug)]
278 pub struct Operator<'a> {
279     op_core: OperatorCore<'a>,
280 }
281 
282 #[derive(Debug)]
283 pub struct CompositeOperator<'a> {
284     op_core: OperatorCore<'a>,
285 }
286 
287 // -----------------------------------------------------------------------------
288 // Destructor
289 // -----------------------------------------------------------------------------
290 impl<'a> Drop for OperatorCore<'a> {
291     fn drop(&mut self) {
292         unsafe {
293             bind_ceed::CeedOperatorDestroy(&mut self.ptr);
294         }
295     }
296 }
297 
298 // -----------------------------------------------------------------------------
299 // Display
300 // -----------------------------------------------------------------------------
301 impl<'a> fmt::Display for OperatorCore<'a> {
302     fn fmt(&self, f: &mut fmt::Formatter) -> fmt::Result {
303         let mut ptr = std::ptr::null_mut();
304         let mut sizeloc = crate::MAX_BUFFER_LENGTH;
305         let cstring = unsafe {
306             let file = bind_ceed::open_memstream(&mut ptr, &mut sizeloc);
307             bind_ceed::CeedOperatorView(self.ptr, file);
308             bind_ceed::fclose(file);
309             CString::from_raw(ptr)
310         };
311         cstring.to_string_lossy().fmt(f)
312     }
313 }
314 
315 /// View an Operator
316 ///
317 /// ```
318 /// # use libceed::prelude::*;
319 /// # fn main() -> libceed::Result<()> {
320 /// # let ceed = libceed::Ceed::default_init();
321 /// let qf = ceed.q_function_interior_by_name("Mass1DBuild")?;
322 ///
323 /// // Operator field arguments
324 /// let ne = 3;
325 /// let q = 4 as usize;
326 /// let mut ind: Vec<i32> = vec![0; 2 * ne];
327 /// for i in 0..ne {
328 ///     ind[2 * i + 0] = i as i32;
329 ///     ind[2 * i + 1] = (i + 1) as i32;
330 /// }
331 /// let r = ceed.elem_restriction(ne, 2, 1, 1, ne + 1, MemType::Host, &ind)?;
332 /// let strides: [i32; 3] = [1, q as i32, q as i32];
333 /// let rq = ceed.strided_elem_restriction(ne, q, 1, q * ne, strides)?;
334 ///
335 /// let b = ceed.basis_tensor_H1_Lagrange(1, 1, 2, q, QuadMode::Gauss)?;
336 ///
337 /// // Operator fields
338 /// let op = ceed
339 ///     .operator(&qf, QFunctionOpt::None, QFunctionOpt::None)?
340 ///     .name("mass")?
341 ///     .field("dx", &r, &b, VectorOpt::Active)?
342 ///     .field("weights", ElemRestrictionOpt::None, &b, VectorOpt::None)?
343 ///     .field("qdata", &rq, BasisOpt::None, VectorOpt::Active)?;
344 ///
345 /// println!("{}", op);
346 /// # Ok(())
347 /// # }
348 /// ```
349 impl<'a> fmt::Display for Operator<'a> {
350     fn fmt(&self, f: &mut fmt::Formatter) -> fmt::Result {
351         self.op_core.fmt(f)
352     }
353 }
354 
355 /// View a composite Operator
356 ///
357 /// ```
358 /// # use libceed::prelude::*;
359 /// # fn main() -> libceed::Result<()> {
360 /// # let ceed = libceed::Ceed::default_init();
361 ///
362 /// // Sub operator field arguments
363 /// let ne = 3;
364 /// let q = 4 as usize;
365 /// let mut ind: Vec<i32> = vec![0; 2 * ne];
366 /// for i in 0..ne {
367 ///     ind[2 * i + 0] = i as i32;
368 ///     ind[2 * i + 1] = (i + 1) as i32;
369 /// }
370 /// let r = ceed.elem_restriction(ne, 2, 1, 1, ne + 1, MemType::Host, &ind)?;
371 /// let strides: [i32; 3] = [1, q as i32, q as i32];
372 /// let rq = ceed.strided_elem_restriction(ne, q, 1, q * ne, strides)?;
373 ///
374 /// let b = ceed.basis_tensor_H1_Lagrange(1, 1, 2, q, QuadMode::Gauss)?;
375 ///
376 /// let qdata_mass = ceed.vector(q * ne)?;
377 /// let qdata_diff = ceed.vector(q * ne)?;
378 ///
379 /// let qf_mass = ceed.q_function_interior_by_name("MassApply")?;
380 /// let op_mass = ceed
381 ///     .operator(&qf_mass, QFunctionOpt::None, QFunctionOpt::None)?
382 ///     .name("Mass term")?
383 ///     .field("u", &r, &b, VectorOpt::Active)?
384 ///     .field("qdata", &rq, BasisOpt::None, &qdata_mass)?
385 ///     .field("v", &r, &b, VectorOpt::Active)?;
386 ///
387 /// let qf_diff = ceed.q_function_interior_by_name("Poisson1DApply")?;
388 /// let op_diff = ceed
389 ///     .operator(&qf_diff, QFunctionOpt::None, QFunctionOpt::None)?
390 ///     .name("Poisson term")?
391 ///     .field("du", &r, &b, VectorOpt::Active)?
392 ///     .field("qdata", &rq, BasisOpt::None, &qdata_diff)?
393 ///     .field("dv", &r, &b, VectorOpt::Active)?;
394 ///
395 /// let op = ceed
396 ///     .composite_operator()?
397 ///     .name("Screened Poisson")?
398 ///     .sub_operator(&op_mass)?
399 ///     .sub_operator(&op_diff)?;
400 ///
401 /// println!("{}", op);
402 /// # Ok(())
403 /// # }
404 /// ```
405 impl<'a> fmt::Display for CompositeOperator<'a> {
406     fn fmt(&self, f: &mut fmt::Formatter) -> fmt::Result {
407         self.op_core.fmt(f)
408     }
409 }
410 
411 // -----------------------------------------------------------------------------
412 // Core functionality
413 // -----------------------------------------------------------------------------
414 impl<'a> OperatorCore<'a> {
415     // Error handling
416     #[doc(hidden)]
417     fn check_error(&self, ierr: i32) -> crate::Result<i32> {
418         let mut ptr = std::ptr::null_mut();
419         unsafe {
420             bind_ceed::CeedOperatorGetCeed(self.ptr, &mut ptr);
421         }
422         crate::check_error(ptr, ierr)
423     }
424 
425     // Common implementations
426     pub fn check(&self) -> crate::Result<i32> {
427         let ierr = unsafe { bind_ceed::CeedOperatorCheckReady(self.ptr) };
428         self.check_error(ierr)
429     }
430 
431     pub fn name(&self, name: &str) -> crate::Result<i32> {
432         let name_c = CString::new(name).expect("CString::new failed");
433         let ierr = unsafe { bind_ceed::CeedOperatorSetName(self.ptr, name_c.as_ptr()) };
434         self.check_error(ierr)
435     }
436 
437     pub fn apply(&self, input: &Vector, output: &mut Vector) -> crate::Result<i32> {
438         let ierr = unsafe {
439             bind_ceed::CeedOperatorApply(
440                 self.ptr,
441                 input.ptr,
442                 output.ptr,
443                 bind_ceed::CEED_REQUEST_IMMEDIATE,
444             )
445         };
446         self.check_error(ierr)
447     }
448 
449     pub fn apply_add(&self, input: &Vector, output: &mut Vector) -> crate::Result<i32> {
450         let ierr = unsafe {
451             bind_ceed::CeedOperatorApplyAdd(
452                 self.ptr,
453                 input.ptr,
454                 output.ptr,
455                 bind_ceed::CEED_REQUEST_IMMEDIATE,
456             )
457         };
458         self.check_error(ierr)
459     }
460 
461     pub fn linear_assemble_diagonal(&self, assembled: &mut Vector) -> crate::Result<i32> {
462         let ierr = unsafe {
463             bind_ceed::CeedOperatorLinearAssembleDiagonal(
464                 self.ptr,
465                 assembled.ptr,
466                 bind_ceed::CEED_REQUEST_IMMEDIATE,
467             )
468         };
469         self.check_error(ierr)
470     }
471 
472     pub fn linear_assemble_add_diagonal(&self, assembled: &mut Vector) -> crate::Result<i32> {
473         let ierr = unsafe {
474             bind_ceed::CeedOperatorLinearAssembleAddDiagonal(
475                 self.ptr,
476                 assembled.ptr,
477                 bind_ceed::CEED_REQUEST_IMMEDIATE,
478             )
479         };
480         self.check_error(ierr)
481     }
482 
483     pub fn linear_assemble_point_block_diagonal(
484         &self,
485         assembled: &mut Vector,
486     ) -> crate::Result<i32> {
487         let ierr = unsafe {
488             bind_ceed::CeedOperatorLinearAssemblePointBlockDiagonal(
489                 self.ptr,
490                 assembled.ptr,
491                 bind_ceed::CEED_REQUEST_IMMEDIATE,
492             )
493         };
494         self.check_error(ierr)
495     }
496 
497     pub fn linear_assemble_add_point_block_diagonal(
498         &self,
499         assembled: &mut Vector,
500     ) -> crate::Result<i32> {
501         let ierr = unsafe {
502             bind_ceed::CeedOperatorLinearAssembleAddPointBlockDiagonal(
503                 self.ptr,
504                 assembled.ptr,
505                 bind_ceed::CEED_REQUEST_IMMEDIATE,
506             )
507         };
508         self.check_error(ierr)
509     }
510 }
511 
512 // -----------------------------------------------------------------------------
513 // Operator
514 // -----------------------------------------------------------------------------
515 impl<'a> Operator<'a> {
516     // Constructor
517     pub fn create<'b>(
518         ceed: &crate::Ceed,
519         qf: impl Into<QFunctionOpt<'b>>,
520         dqf: impl Into<QFunctionOpt<'b>>,
521         dqfT: impl Into<QFunctionOpt<'b>>,
522     ) -> crate::Result<Self> {
523         let mut ptr = std::ptr::null_mut();
524         let ierr = unsafe {
525             bind_ceed::CeedOperatorCreate(
526                 ceed.ptr,
527                 qf.into().to_raw(),
528                 dqf.into().to_raw(),
529                 dqfT.into().to_raw(),
530                 &mut ptr,
531             )
532         };
533         ceed.check_error(ierr)?;
534         Ok(Self {
535             op_core: OperatorCore {
536                 ptr,
537                 _lifeline: PhantomData,
538             },
539         })
540     }
541 
542     fn from_raw(ptr: bind_ceed::CeedOperator) -> crate::Result<Self> {
543         Ok(Self {
544             op_core: OperatorCore {
545                 ptr,
546                 _lifeline: PhantomData,
547             },
548         })
549     }
550 
551     /// Set name for Operator printing
552     ///
553     /// * 'name' - Name to set
554     ///
555     /// ```
556     /// # use libceed::prelude::*;
557     /// # fn main() -> libceed::Result<()> {
558     /// # let ceed = libceed::Ceed::default_init();
559     /// let qf = ceed.q_function_interior_by_name("Mass1DBuild")?;
560     ///
561     /// // Operator field arguments
562     /// let ne = 3;
563     /// let q = 4 as usize;
564     /// let mut ind: Vec<i32> = vec![0; 2 * ne];
565     /// for i in 0..ne {
566     ///     ind[2 * i + 0] = i as i32;
567     ///     ind[2 * i + 1] = (i + 1) as i32;
568     /// }
569     /// let r = ceed.elem_restriction(ne, 2, 1, 1, ne + 1, MemType::Host, &ind)?;
570     /// let strides: [i32; 3] = [1, q as i32, q as i32];
571     /// let rq = ceed.strided_elem_restriction(ne, q, 1, q * ne, strides)?;
572     ///
573     /// let b = ceed.basis_tensor_H1_Lagrange(1, 1, 2, q, QuadMode::Gauss)?;
574     ///
575     /// // Operator fields
576     /// let op = ceed
577     ///     .operator(&qf, QFunctionOpt::None, QFunctionOpt::None)?
578     ///     .name("mass")?
579     ///     .field("dx", &r, &b, VectorOpt::Active)?
580     ///     .field("weights", ElemRestrictionOpt::None, &b, VectorOpt::None)?
581     ///     .field("qdata", &rq, BasisOpt::None, VectorOpt::Active)?;
582     /// # Ok(())
583     /// # }
584     /// ```
585     #[allow(unused_mut)]
586     pub fn name(mut self, name: &str) -> crate::Result<Self> {
587         self.op_core.name(name)?;
588         Ok(self)
589     }
590 
591     /// Apply Operator to a vector
592     ///
593     /// * `input`  - Input Vector
594     /// * `output` - Output Vector
595     ///
596     /// ```
597     /// # use libceed::prelude::*;
598     /// # fn main() -> libceed::Result<()> {
599     /// # let ceed = libceed::Ceed::default_init();
600     /// let ne = 4;
601     /// let p = 3;
602     /// let q = 4;
603     /// let ndofs = p * ne - ne + 1;
604     ///
605     /// // Vectors
606     /// let x = ceed.vector_from_slice(&[-1., -0.5, 0.0, 0.5, 1.0])?;
607     /// let mut qdata = ceed.vector(ne * q)?;
608     /// qdata.set_value(0.0);
609     /// let u = ceed.vector_from_slice(&vec![1.0; ndofs])?;
610     /// let mut v = ceed.vector(ndofs)?;
611     /// v.set_value(0.0);
612     ///
613     /// // Restrictions
614     /// let mut indx: Vec<i32> = vec![0; 2 * ne];
615     /// for i in 0..ne {
616     ///     indx[2 * i + 0] = i as i32;
617     ///     indx[2 * i + 1] = (i + 1) as i32;
618     /// }
619     /// let rx = ceed.elem_restriction(ne, 2, 1, 1, ne + 1, MemType::Host, &indx)?;
620     /// let mut indu: Vec<i32> = vec![0; p * ne];
621     /// for i in 0..ne {
622     ///     indu[p * i + 0] = i as i32;
623     ///     indu[p * i + 1] = (i + 1) as i32;
624     ///     indu[p * i + 2] = (i + 2) as i32;
625     /// }
626     /// let ru = ceed.elem_restriction(ne, 3, 1, 1, ndofs, MemType::Host, &indu)?;
627     /// let strides: [i32; 3] = [1, q as i32, q as i32];
628     /// let rq = ceed.strided_elem_restriction(ne, q, 1, q * ne, strides)?;
629     ///
630     /// // Bases
631     /// let bx = ceed.basis_tensor_H1_Lagrange(1, 1, 2, q, QuadMode::Gauss)?;
632     /// let bu = ceed.basis_tensor_H1_Lagrange(1, 1, p, q, QuadMode::Gauss)?;
633     ///
634     /// // Build quadrature data
635     /// let qf_build = ceed.q_function_interior_by_name("Mass1DBuild")?;
636     /// ceed.operator(&qf_build, QFunctionOpt::None, QFunctionOpt::None)?
637     ///     .field("dx", &rx, &bx, VectorOpt::Active)?
638     ///     .field("weights", ElemRestrictionOpt::None, &bx, VectorOpt::None)?
639     ///     .field("qdata", &rq, BasisOpt::None, VectorOpt::Active)?
640     ///     .apply(&x, &mut qdata)?;
641     ///
642     /// // Mass operator
643     /// let qf_mass = ceed.q_function_interior_by_name("MassApply")?;
644     /// let op_mass = ceed
645     ///     .operator(&qf_mass, QFunctionOpt::None, QFunctionOpt::None)?
646     ///     .field("u", &ru, &bu, VectorOpt::Active)?
647     ///     .field("qdata", &rq, BasisOpt::None, &qdata)?
648     ///     .field("v", &ru, &bu, VectorOpt::Active)?;
649     ///
650     /// v.set_value(0.0);
651     /// op_mass.apply(&u, &mut v)?;
652     ///
653     /// // Check
654     /// let sum: Scalar = v.view()?.iter().sum();
655     /// assert!(
656     ///     (sum - 2.0).abs() < 10.0 * libceed::EPSILON,
657     ///     "Incorrect interval length computed"
658     /// );
659     /// # Ok(())
660     /// # }
661     /// ```
662     pub fn apply(&self, input: &Vector, output: &mut Vector) -> crate::Result<i32> {
663         self.op_core.apply(input, output)
664     }
665 
666     /// Apply Operator to a vector and add result to output vector
667     ///
668     /// * `input`  - Input Vector
669     /// * `output` - Output Vector
670     ///
671     /// ```
672     /// # use libceed::prelude::*;
673     /// # fn main() -> libceed::Result<()> {
674     /// # let ceed = libceed::Ceed::default_init();
675     /// let ne = 4;
676     /// let p = 3;
677     /// let q = 4;
678     /// let ndofs = p * ne - ne + 1;
679     ///
680     /// // Vectors
681     /// let x = ceed.vector_from_slice(&[-1., -0.5, 0.0, 0.5, 1.0])?;
682     /// let mut qdata = ceed.vector(ne * q)?;
683     /// qdata.set_value(0.0);
684     /// let u = ceed.vector_from_slice(&vec![1.0; ndofs])?;
685     /// let mut v = ceed.vector(ndofs)?;
686     ///
687     /// // Restrictions
688     /// let mut indx: Vec<i32> = vec![0; 2 * ne];
689     /// for i in 0..ne {
690     ///     indx[2 * i + 0] = i as i32;
691     ///     indx[2 * i + 1] = (i + 1) as i32;
692     /// }
693     /// let rx = ceed.elem_restriction(ne, 2, 1, 1, ne + 1, MemType::Host, &indx)?;
694     /// let mut indu: Vec<i32> = vec![0; p * ne];
695     /// for i in 0..ne {
696     ///     indu[p * i + 0] = i as i32;
697     ///     indu[p * i + 1] = (i + 1) as i32;
698     ///     indu[p * i + 2] = (i + 2) as i32;
699     /// }
700     /// let ru = ceed.elem_restriction(ne, 3, 1, 1, ndofs, MemType::Host, &indu)?;
701     /// let strides: [i32; 3] = [1, q as i32, q as i32];
702     /// let rq = ceed.strided_elem_restriction(ne, q, 1, q * ne, strides)?;
703     ///
704     /// // Bases
705     /// let bx = ceed.basis_tensor_H1_Lagrange(1, 1, 2, q, QuadMode::Gauss)?;
706     /// let bu = ceed.basis_tensor_H1_Lagrange(1, 1, p, q, QuadMode::Gauss)?;
707     ///
708     /// // Build quadrature data
709     /// let qf_build = ceed.q_function_interior_by_name("Mass1DBuild")?;
710     /// ceed.operator(&qf_build, QFunctionOpt::None, QFunctionOpt::None)?
711     ///     .field("dx", &rx, &bx, VectorOpt::Active)?
712     ///     .field("weights", ElemRestrictionOpt::None, &bx, VectorOpt::None)?
713     ///     .field("qdata", &rq, BasisOpt::None, VectorOpt::Active)?
714     ///     .apply(&x, &mut qdata)?;
715     ///
716     /// // Mass operator
717     /// let qf_mass = ceed.q_function_interior_by_name("MassApply")?;
718     /// let op_mass = ceed
719     ///     .operator(&qf_mass, QFunctionOpt::None, QFunctionOpt::None)?
720     ///     .field("u", &ru, &bu, VectorOpt::Active)?
721     ///     .field("qdata", &rq, BasisOpt::None, &qdata)?
722     ///     .field("v", &ru, &bu, VectorOpt::Active)?;
723     ///
724     /// v.set_value(1.0);
725     /// op_mass.apply_add(&u, &mut v)?;
726     ///
727     /// // Check
728     /// let sum: Scalar = v.view()?.iter().sum();
729     /// assert!(
730     ///     (sum - (2.0 + ndofs as Scalar)).abs() < 10.0 * libceed::EPSILON,
731     ///     "Incorrect interval length computed and added"
732     /// );
733     /// # Ok(())
734     /// # }
735     /// ```
736     pub fn apply_add(&self, input: &Vector, output: &mut Vector) -> crate::Result<i32> {
737         self.op_core.apply_add(input, output)
738     }
739 
740     /// Provide a field to a Operator for use by its QFunction
741     ///
742     /// * `fieldname` - Name of the field (to be matched with the name used by
743     ///                   the QFunction)
744     /// * `r`         - ElemRestriction
745     /// * `b`         - Basis in which the field resides or `BasisOpt::None`
746     /// * `v`         - Vector to be used by Operator or `VectorOpt::Active` if field
747     ///                   is active or `VectorOpt::None` if using `Weight` with the
748     ///                   QFunction
749     ///
750     ///
751     /// ```
752     /// # use libceed::prelude::*;
753     /// # fn main() -> libceed::Result<()> {
754     /// # let ceed = libceed::Ceed::default_init();
755     /// let qf = ceed.q_function_interior_by_name("Mass1DBuild")?;
756     /// let mut op = ceed.operator(&qf, QFunctionOpt::None, QFunctionOpt::None)?;
757     ///
758     /// // Operator field arguments
759     /// let ne = 3;
760     /// let q = 4;
761     /// let mut ind: Vec<i32> = vec![0; 2 * ne];
762     /// for i in 0..ne {
763     ///     ind[2 * i + 0] = i as i32;
764     ///     ind[2 * i + 1] = (i + 1) as i32;
765     /// }
766     /// let r = ceed.elem_restriction(ne, 2, 1, 1, ne + 1, MemType::Host, &ind)?;
767     ///
768     /// let b = ceed.basis_tensor_H1_Lagrange(1, 1, 2, q, QuadMode::Gauss)?;
769     ///
770     /// // Operator field
771     /// op = op.field("dx", &r, &b, VectorOpt::Active)?;
772     /// # Ok(())
773     /// # }
774     /// ```
775     #[allow(unused_mut)]
776     pub fn field<'b>(
777         mut self,
778         fieldname: &str,
779         r: impl Into<ElemRestrictionOpt<'b>>,
780         b: impl Into<BasisOpt<'b>>,
781         v: impl Into<VectorOpt<'b>>,
782     ) -> crate::Result<Self> {
783         let fieldname = CString::new(fieldname).expect("CString::new failed");
784         let fieldname = fieldname.as_ptr() as *const i8;
785         let ierr = unsafe {
786             bind_ceed::CeedOperatorSetField(
787                 self.op_core.ptr,
788                 fieldname,
789                 r.into().to_raw(),
790                 b.into().to_raw(),
791                 v.into().to_raw(),
792             )
793         };
794         self.op_core.check_error(ierr)?;
795         Ok(self)
796     }
797 
798     /// Get a slice of Operator inputs
799     ///
800     /// ```
801     /// # use libceed::prelude::*;
802     /// # fn main() -> libceed::Result<()> {
803     /// # let ceed = libceed::Ceed::default_init();
804     /// let qf = ceed.q_function_interior_by_name("Mass1DBuild")?;
805     ///
806     /// // Operator field arguments
807     /// let ne = 3;
808     /// let q = 4 as usize;
809     /// let mut ind: Vec<i32> = vec![0; 2 * ne];
810     /// for i in 0..ne {
811     ///     ind[2 * i + 0] = i as i32;
812     ///     ind[2 * i + 1] = (i + 1) as i32;
813     /// }
814     /// let r = ceed.elem_restriction(ne, 2, 1, 1, ne + 1, MemType::Host, &ind)?;
815     /// let strides: [i32; 3] = [1, q as i32, q as i32];
816     /// let rq = ceed.strided_elem_restriction(ne, q, 1, q * ne, strides)?;
817     ///
818     /// let b = ceed.basis_tensor_H1_Lagrange(1, 1, 2, q, QuadMode::Gauss)?;
819     ///
820     /// // Operator fields
821     /// let op = ceed
822     ///     .operator(&qf, QFunctionOpt::None, QFunctionOpt::None)?
823     ///     .field("dx", &r, &b, VectorOpt::Active)?
824     ///     .field("weights", ElemRestrictionOpt::None, &b, VectorOpt::None)?
825     ///     .field("qdata", &rq, BasisOpt::None, VectorOpt::Active)?;
826     ///
827     /// let inputs = op.inputs()?;
828     ///
829     /// assert_eq!(inputs.len(), 2, "Incorrect inputs array");
830     /// # Ok(())
831     /// # }
832     /// ```
833     pub fn inputs(&self) -> crate::Result<&[crate::OperatorField]> {
834         // Get array of raw C pointers for inputs
835         let mut num_inputs = 0;
836         let mut inputs_ptr = std::ptr::null_mut();
837         let ierr = unsafe {
838             bind_ceed::CeedOperatorGetFields(
839                 self.op_core.ptr,
840                 &mut num_inputs,
841                 &mut inputs_ptr,
842                 std::ptr::null_mut() as *mut bind_ceed::CeedInt,
843                 std::ptr::null_mut() as *mut *mut bind_ceed::CeedOperatorField,
844             )
845         };
846         self.op_core.check_error(ierr)?;
847         // Convert raw C pointers to fixed length slice
848         let inputs_slice = unsafe {
849             std::slice::from_raw_parts(
850                 inputs_ptr as *const crate::OperatorField,
851                 num_inputs as usize,
852             )
853         };
854         Ok(inputs_slice)
855     }
856 
857     /// Get a slice of Operator outputs
858     ///
859     /// ```
860     /// # use libceed::prelude::*;
861     /// # fn main() -> libceed::Result<()> {
862     /// # let ceed = libceed::Ceed::default_init();
863     /// let qf = ceed.q_function_interior_by_name("Mass1DBuild")?;
864     ///
865     /// // Operator field arguments
866     /// let ne = 3;
867     /// let q = 4 as usize;
868     /// let mut ind: Vec<i32> = vec![0; 2 * ne];
869     /// for i in 0..ne {
870     ///     ind[2 * i + 0] = i as i32;
871     ///     ind[2 * i + 1] = (i + 1) as i32;
872     /// }
873     /// let r = ceed.elem_restriction(ne, 2, 1, 1, ne + 1, MemType::Host, &ind)?;
874     /// let strides: [i32; 3] = [1, q as i32, q as i32];
875     /// let rq = ceed.strided_elem_restriction(ne, q, 1, q * ne, strides)?;
876     ///
877     /// let b = ceed.basis_tensor_H1_Lagrange(1, 1, 2, q, QuadMode::Gauss)?;
878     ///
879     /// // Operator fields
880     /// let op = ceed
881     ///     .operator(&qf, QFunctionOpt::None, QFunctionOpt::None)?
882     ///     .field("dx", &r, &b, VectorOpt::Active)?
883     ///     .field("weights", ElemRestrictionOpt::None, &b, VectorOpt::None)?
884     ///     .field("qdata", &rq, BasisOpt::None, VectorOpt::Active)?;
885     ///
886     /// let outputs = op.outputs()?;
887     ///
888     /// assert_eq!(outputs.len(), 1, "Incorrect outputs array");
889     /// # Ok(())
890     /// # }
891     /// ```
892     pub fn outputs(&self) -> crate::Result<&[crate::OperatorField]> {
893         // Get array of raw C pointers for outputs
894         let mut num_outputs = 0;
895         let mut outputs_ptr = std::ptr::null_mut();
896         let ierr = unsafe {
897             bind_ceed::CeedOperatorGetFields(
898                 self.op_core.ptr,
899                 std::ptr::null_mut() as *mut bind_ceed::CeedInt,
900                 std::ptr::null_mut() as *mut *mut bind_ceed::CeedOperatorField,
901                 &mut num_outputs,
902                 &mut outputs_ptr,
903             )
904         };
905         self.op_core.check_error(ierr)?;
906         // Convert raw C pointers to fixed length slice
907         let outputs_slice = unsafe {
908             std::slice::from_raw_parts(
909                 outputs_ptr as *const crate::OperatorField,
910                 num_outputs as usize,
911             )
912         };
913         Ok(outputs_slice)
914     }
915 
916     /// Check if Operator is setup correctly
917     ///
918     /// ```
919     /// # use libceed::prelude::*;
920     /// # fn main() -> libceed::Result<()> {
921     /// # let ceed = libceed::Ceed::default_init();
922     /// let ne = 4;
923     /// let p = 3;
924     /// let q = 4;
925     /// let ndofs = p * ne - ne + 1;
926     ///
927     /// // Restrictions
928     /// let mut indx: Vec<i32> = vec![0; 2 * ne];
929     /// for i in 0..ne {
930     ///     indx[2 * i + 0] = i as i32;
931     ///     indx[2 * i + 1] = (i + 1) as i32;
932     /// }
933     /// let rx = ceed.elem_restriction(ne, 2, 1, 1, ne + 1, MemType::Host, &indx)?;
934     /// let strides: [i32; 3] = [1, q as i32, q as i32];
935     /// let rq = ceed.strided_elem_restriction(ne, q, 1, q * ne, strides)?;
936     ///
937     /// // Bases
938     /// let bx = ceed.basis_tensor_H1_Lagrange(1, 1, 2, q, QuadMode::Gauss)?;
939     ///
940     /// // Build quadrature data
941     /// let qf_build = ceed.q_function_interior_by_name("Mass1DBuild")?;
942     /// let op_build = ceed
943     ///     .operator(&qf_build, QFunctionOpt::None, QFunctionOpt::None)?
944     ///     .field("dx", &rx, &bx, VectorOpt::Active)?
945     ///     .field("weights", ElemRestrictionOpt::None, &bx, VectorOpt::None)?
946     ///     .field("qdata", &rq, BasisOpt::None, VectorOpt::Active)?;
947     ///
948     /// // Check
949     /// assert!(op_build.check().is_ok(), "Operator setup incorrectly");
950     /// # Ok(())
951     /// # }
952     /// ```
953     pub fn check(self) -> crate::Result<Self> {
954         self.op_core.check()?;
955         Ok(self)
956     }
957 
958     /// Get the number of elements associated with an Operator
959     ///
960     ///
961     /// ```
962     /// # use libceed::prelude::*;
963     /// # fn main() -> libceed::Result<()> {
964     /// # let ceed = libceed::Ceed::default_init();
965     /// let qf = ceed.q_function_interior_by_name("Mass1DBuild")?;
966     /// let mut op = ceed.operator(&qf, QFunctionOpt::None, QFunctionOpt::None)?;
967     ///
968     /// // Operator field arguments
969     /// let ne = 3;
970     /// let q = 4;
971     /// let mut ind: Vec<i32> = vec![0; 2 * ne];
972     /// for i in 0..ne {
973     ///     ind[2 * i + 0] = i as i32;
974     ///     ind[2 * i + 1] = (i + 1) as i32;
975     /// }
976     /// let r = ceed.elem_restriction(ne, 2, 1, 1, ne + 1, MemType::Host, &ind)?;
977     ///
978     /// let b = ceed.basis_tensor_H1_Lagrange(1, 1, 2, q, QuadMode::Gauss)?;
979     ///
980     /// // Operator field
981     /// op = op.field("dx", &r, &b, VectorOpt::Active)?;
982     ///
983     /// // Check number of elements
984     /// assert_eq!(op.num_elements(), ne, "Incorrect number of elements");
985     /// # Ok(())
986     /// # }
987     /// ```
988     pub fn num_elements(&self) -> usize {
989         let mut nelem = 0;
990         unsafe { bind_ceed::CeedOperatorGetNumElements(self.op_core.ptr, &mut nelem) };
991         usize::try_from(nelem).unwrap()
992     }
993 
994     /// Get the number of quadrature points associated with each element of
995     ///   an Operator
996     ///
997     ///
998     /// ```
999     /// # use libceed::prelude::*;
1000     /// # fn main() -> libceed::Result<()> {
1001     /// # let ceed = libceed::Ceed::default_init();
1002     /// let qf = ceed.q_function_interior_by_name("Mass1DBuild")?;
1003     /// let mut op = ceed.operator(&qf, QFunctionOpt::None, QFunctionOpt::None)?;
1004     ///
1005     /// // Operator field arguments
1006     /// let ne = 3;
1007     /// let q = 4;
1008     /// let mut ind: Vec<i32> = vec![0; 2 * ne];
1009     /// for i in 0..ne {
1010     ///     ind[2 * i + 0] = i as i32;
1011     ///     ind[2 * i + 1] = (i + 1) as i32;
1012     /// }
1013     /// let r = ceed.elem_restriction(ne, 2, 1, 1, ne + 1, MemType::Host, &ind)?;
1014     ///
1015     /// let b = ceed.basis_tensor_H1_Lagrange(1, 1, 2, q, QuadMode::Gauss)?;
1016     ///
1017     /// // Operator field
1018     /// op = op.field("dx", &r, &b, VectorOpt::Active)?;
1019     ///
1020     /// // Check number of quadrature points
1021     /// assert_eq!(
1022     ///     op.num_quadrature_points(),
1023     ///     q,
1024     ///     "Incorrect number of quadrature points"
1025     /// );
1026     /// # Ok(())
1027     /// # }
1028     /// ```
1029     pub fn num_quadrature_points(&self) -> usize {
1030         let mut Q = 0;
1031         unsafe { bind_ceed::CeedOperatorGetNumQuadraturePoints(self.op_core.ptr, &mut Q) };
1032         usize::try_from(Q).unwrap()
1033     }
1034 
1035     /// Assemble the diagonal of a square linear Operator
1036     ///
1037     /// This overwrites a Vector with the diagonal of a linear Operator.
1038     ///
1039     /// Note: Currently only non-composite Operators with a single field and
1040     ///      composite Operators with single field sub-operators are supported.
1041     ///
1042     /// * `op`        - Operator to assemble QFunction
1043     /// * `assembled` - Vector to store assembled Operator diagonal
1044     ///
1045     /// ```
1046     /// # use libceed::prelude::*;
1047     /// # fn main() -> libceed::Result<()> {
1048     /// # let ceed = libceed::Ceed::default_init();
1049     /// let ne = 4;
1050     /// let p = 3;
1051     /// let q = 4;
1052     /// let ndofs = p * ne - ne + 1;
1053     ///
1054     /// // Vectors
1055     /// let x = ceed.vector_from_slice(&[-1., -0.5, 0.0, 0.5, 1.0])?;
1056     /// let mut qdata = ceed.vector(ne * q)?;
1057     /// qdata.set_value(0.0);
1058     /// let mut u = ceed.vector(ndofs)?;
1059     /// u.set_value(1.0);
1060     /// let mut v = ceed.vector(ndofs)?;
1061     /// v.set_value(0.0);
1062     ///
1063     /// // Restrictions
1064     /// let mut indx: Vec<i32> = vec![0; 2 * ne];
1065     /// for i in 0..ne {
1066     ///     indx[2 * i + 0] = i as i32;
1067     ///     indx[2 * i + 1] = (i + 1) as i32;
1068     /// }
1069     /// let rx = ceed.elem_restriction(ne, 2, 1, 1, ne + 1, MemType::Host, &indx)?;
1070     /// let mut indu: Vec<i32> = vec![0; p * ne];
1071     /// for i in 0..ne {
1072     ///     indu[p * i + 0] = (2 * i) as i32;
1073     ///     indu[p * i + 1] = (2 * i + 1) as i32;
1074     ///     indu[p * i + 2] = (2 * i + 2) as i32;
1075     /// }
1076     /// let ru = ceed.elem_restriction(ne, p, 1, 1, ndofs, MemType::Host, &indu)?;
1077     /// let strides: [i32; 3] = [1, q as i32, q as i32];
1078     /// let rq = ceed.strided_elem_restriction(ne, q, 1, q * ne, strides)?;
1079     ///
1080     /// // Bases
1081     /// let bx = ceed.basis_tensor_H1_Lagrange(1, 1, 2, q, QuadMode::Gauss)?;
1082     /// let bu = ceed.basis_tensor_H1_Lagrange(1, 1, p, q, QuadMode::Gauss)?;
1083     ///
1084     /// // Build quadrature data
1085     /// let qf_build = ceed.q_function_interior_by_name("Mass1DBuild")?;
1086     /// ceed.operator(&qf_build, QFunctionOpt::None, QFunctionOpt::None)?
1087     ///     .field("dx", &rx, &bx, VectorOpt::Active)?
1088     ///     .field("weights", ElemRestrictionOpt::None, &bx, VectorOpt::None)?
1089     ///     .field("qdata", &rq, BasisOpt::None, VectorOpt::Active)?
1090     ///     .apply(&x, &mut qdata)?;
1091     ///
1092     /// // Mass operator
1093     /// let qf_mass = ceed.q_function_interior_by_name("MassApply")?;
1094     /// let op_mass = ceed
1095     ///     .operator(&qf_mass, QFunctionOpt::None, QFunctionOpt::None)?
1096     ///     .field("u", &ru, &bu, VectorOpt::Active)?
1097     ///     .field("qdata", &rq, BasisOpt::None, &qdata)?
1098     ///     .field("v", &ru, &bu, VectorOpt::Active)?;
1099     ///
1100     /// // Diagonal
1101     /// let mut diag = ceed.vector(ndofs)?;
1102     /// diag.set_value(0.0);
1103     /// op_mass.linear_assemble_diagonal(&mut diag)?;
1104     ///
1105     /// // Manual diagonal computation
1106     /// let mut true_diag = ceed.vector(ndofs)?;
1107     /// true_diag.set_value(0.0)?;
1108     /// for i in 0..ndofs {
1109     ///     u.set_value(0.0);
1110     ///     {
1111     ///         let mut u_array = u.view_mut()?;
1112     ///         u_array[i] = 1.;
1113     ///     }
1114     ///
1115     ///     op_mass.apply(&u, &mut v)?;
1116     ///
1117     ///     {
1118     ///         let v_array = v.view()?;
1119     ///         let mut true_array = true_diag.view_mut()?;
1120     ///         true_array[i] = v_array[i];
1121     ///     }
1122     /// }
1123     ///
1124     /// // Check
1125     /// diag.view()?
1126     ///     .iter()
1127     ///     .zip(true_diag.view()?.iter())
1128     ///     .for_each(|(computed, actual)| {
1129     ///         assert!(
1130     ///             (*computed - *actual).abs() < 10.0 * libceed::EPSILON,
1131     ///             "Diagonal entry incorrect"
1132     ///         );
1133     ///     });
1134     /// # Ok(())
1135     /// # }
1136     /// ```
1137     pub fn linear_assemble_diagonal(&self, assembled: &mut Vector) -> crate::Result<i32> {
1138         self.op_core.linear_assemble_diagonal(assembled)
1139     }
1140 
1141     /// Assemble the diagonal of a square linear Operator
1142     ///
1143     /// This sums into a Vector with the diagonal of a linear Operator.
1144     ///
1145     /// Note: Currently only non-composite Operators with a single field and
1146     ///      composite Operators with single field sub-operators are supported.
1147     ///
1148     /// * `op`        - Operator to assemble QFunction
1149     /// * `assembled` - Vector to store assembled Operator diagonal
1150     ///
1151     ///
1152     /// ```
1153     /// # use libceed::prelude::*;
1154     /// # fn main() -> libceed::Result<()> {
1155     /// # let ceed = libceed::Ceed::default_init();
1156     /// let ne = 4;
1157     /// let p = 3;
1158     /// let q = 4;
1159     /// let ndofs = p * ne - ne + 1;
1160     ///
1161     /// // Vectors
1162     /// let x = ceed.vector_from_slice(&[-1., -0.5, 0.0, 0.5, 1.0])?;
1163     /// let mut qdata = ceed.vector(ne * q)?;
1164     /// qdata.set_value(0.0);
1165     /// let mut u = ceed.vector(ndofs)?;
1166     /// u.set_value(1.0);
1167     /// let mut v = ceed.vector(ndofs)?;
1168     /// v.set_value(0.0);
1169     ///
1170     /// // Restrictions
1171     /// let mut indx: Vec<i32> = vec![0; 2 * ne];
1172     /// for i in 0..ne {
1173     ///     indx[2 * i + 0] = i as i32;
1174     ///     indx[2 * i + 1] = (i + 1) as i32;
1175     /// }
1176     /// let rx = ceed.elem_restriction(ne, 2, 1, 1, ne + 1, MemType::Host, &indx)?;
1177     /// let mut indu: Vec<i32> = vec![0; p * ne];
1178     /// for i in 0..ne {
1179     ///     indu[p * i + 0] = (2 * i) as i32;
1180     ///     indu[p * i + 1] = (2 * i + 1) as i32;
1181     ///     indu[p * i + 2] = (2 * i + 2) as i32;
1182     /// }
1183     /// let ru = ceed.elem_restriction(ne, p, 1, 1, ndofs, MemType::Host, &indu)?;
1184     /// let strides: [i32; 3] = [1, q as i32, q as i32];
1185     /// let rq = ceed.strided_elem_restriction(ne, q, 1, q * ne, strides)?;
1186     ///
1187     /// // Bases
1188     /// let bx = ceed.basis_tensor_H1_Lagrange(1, 1, 2, q, QuadMode::Gauss)?;
1189     /// let bu = ceed.basis_tensor_H1_Lagrange(1, 1, p, q, QuadMode::Gauss)?;
1190     ///
1191     /// // Build quadrature data
1192     /// let qf_build = ceed.q_function_interior_by_name("Mass1DBuild")?;
1193     /// ceed.operator(&qf_build, QFunctionOpt::None, QFunctionOpt::None)?
1194     ///     .field("dx", &rx, &bx, VectorOpt::Active)?
1195     ///     .field("weights", ElemRestrictionOpt::None, &bx, VectorOpt::None)?
1196     ///     .field("qdata", &rq, BasisOpt::None, VectorOpt::Active)?
1197     ///     .apply(&x, &mut qdata)?;
1198     ///
1199     /// // Mass operator
1200     /// let qf_mass = ceed.q_function_interior_by_name("MassApply")?;
1201     /// let op_mass = ceed
1202     ///     .operator(&qf_mass, QFunctionOpt::None, QFunctionOpt::None)?
1203     ///     .field("u", &ru, &bu, VectorOpt::Active)?
1204     ///     .field("qdata", &rq, BasisOpt::None, &qdata)?
1205     ///     .field("v", &ru, &bu, VectorOpt::Active)?;
1206     ///
1207     /// // Diagonal
1208     /// let mut diag = ceed.vector(ndofs)?;
1209     /// diag.set_value(1.0);
1210     /// op_mass.linear_assemble_add_diagonal(&mut diag)?;
1211     ///
1212     /// // Manual diagonal computation
1213     /// let mut true_diag = ceed.vector(ndofs)?;
1214     /// true_diag.set_value(0.0)?;
1215     /// for i in 0..ndofs {
1216     ///     u.set_value(0.0);
1217     ///     {
1218     ///         let mut u_array = u.view_mut()?;
1219     ///         u_array[i] = 1.;
1220     ///     }
1221     ///
1222     ///     op_mass.apply(&u, &mut v)?;
1223     ///
1224     ///     {
1225     ///         let v_array = v.view()?;
1226     ///         let mut true_array = true_diag.view_mut()?;
1227     ///         true_array[i] = v_array[i] + 1.0;
1228     ///     }
1229     /// }
1230     ///
1231     /// // Check
1232     /// diag.view()?
1233     ///     .iter()
1234     ///     .zip(true_diag.view()?.iter())
1235     ///     .for_each(|(computed, actual)| {
1236     ///         assert!(
1237     ///             (*computed - *actual).abs() < 10.0 * libceed::EPSILON,
1238     ///             "Diagonal entry incorrect"
1239     ///         );
1240     ///     });
1241     /// # Ok(())
1242     /// # }
1243     /// ```
1244     pub fn linear_assemble_add_diagonal(&self, assembled: &mut Vector) -> crate::Result<i32> {
1245         self.op_core.linear_assemble_add_diagonal(assembled)
1246     }
1247 
1248     /// Assemble the point block diagonal of a square linear Operator
1249     ///
1250     /// This overwrites a Vector with the point block diagonal of a linear
1251     /// Operator.
1252     ///
1253     /// Note: Currently only non-composite Operators with a single field and
1254     ///      composite Operators with single field sub-operators are supported.
1255     ///
1256     /// * `op`        - Operator to assemble QFunction
1257     /// * `assembled` - Vector to store assembled CeedOperator point block
1258     ///                   diagonal, provided in row-major form with an
1259     ///                   `ncomp * ncomp` block at each node. The dimensions of
1260     ///                   this vector are derived from the active vector for
1261     ///                   the CeedOperator. The array has shape
1262     ///                   `[nodes, component out, component in]`.
1263     ///
1264     /// ```
1265     /// # use libceed::prelude::*;
1266     /// # fn main() -> libceed::Result<()> {
1267     /// # let ceed = libceed::Ceed::default_init();
1268     /// let ne = 4;
1269     /// let p = 3;
1270     /// let q = 4;
1271     /// let ncomp = 2;
1272     /// let ndofs = p * ne - ne + 1;
1273     ///
1274     /// // Vectors
1275     /// let x = ceed.vector_from_slice(&[-1., -0.5, 0.0, 0.5, 1.0])?;
1276     /// let mut qdata = ceed.vector(ne * q)?;
1277     /// qdata.set_value(0.0);
1278     /// let mut u = ceed.vector(ncomp * ndofs)?;
1279     /// u.set_value(1.0);
1280     /// let mut v = ceed.vector(ncomp * ndofs)?;
1281     /// v.set_value(0.0);
1282     ///
1283     /// // Restrictions
1284     /// let mut indx: Vec<i32> = vec![0; 2 * ne];
1285     /// for i in 0..ne {
1286     ///     indx[2 * i + 0] = i as i32;
1287     ///     indx[2 * i + 1] = (i + 1) as i32;
1288     /// }
1289     /// let rx = ceed.elem_restriction(ne, 2, 1, 1, ne + 1, MemType::Host, &indx)?;
1290     /// let mut indu: Vec<i32> = vec![0; p * ne];
1291     /// for i in 0..ne {
1292     ///     indu[p * i + 0] = (2 * i) as i32;
1293     ///     indu[p * i + 1] = (2 * i + 1) as i32;
1294     ///     indu[p * i + 2] = (2 * i + 2) as i32;
1295     /// }
1296     /// let ru = ceed.elem_restriction(ne, p, ncomp, ndofs, ncomp * ndofs, MemType::Host, &indu)?;
1297     /// let strides: [i32; 3] = [1, q as i32, q as i32];
1298     /// let rq = ceed.strided_elem_restriction(ne, q, 1, q * ne, strides)?;
1299     ///
1300     /// // Bases
1301     /// let bx = ceed.basis_tensor_H1_Lagrange(1, 1, 2, q, QuadMode::Gauss)?;
1302     /// let bu = ceed.basis_tensor_H1_Lagrange(1, ncomp, p, q, QuadMode::Gauss)?;
1303     ///
1304     /// // Build quadrature data
1305     /// let qf_build = ceed.q_function_interior_by_name("Mass1DBuild")?;
1306     /// ceed.operator(&qf_build, QFunctionOpt::None, QFunctionOpt::None)?
1307     ///     .field("dx", &rx, &bx, VectorOpt::Active)?
1308     ///     .field("weights", ElemRestrictionOpt::None, &bx, VectorOpt::None)?
1309     ///     .field("qdata", &rq, BasisOpt::None, VectorOpt::Active)?
1310     ///     .apply(&x, &mut qdata)?;
1311     ///
1312     /// // Mass operator
1313     /// let mut mass_2_comp = |[u, qdata, ..]: QFunctionInputs, [v, ..]: QFunctionOutputs| {
1314     ///     // Number of quadrature points
1315     ///     let q = qdata.len();
1316     ///
1317     ///     // Iterate over quadrature points
1318     ///     for i in 0..q {
1319     ///         v[i + 0 * q] = u[i + 1 * q] * qdata[i];
1320     ///         v[i + 1 * q] = u[i + 0 * q] * qdata[i];
1321     ///     }
1322     ///
1323     ///     // Return clean error code
1324     ///     0
1325     /// };
1326     ///
1327     /// let qf_mass = ceed
1328     ///     .q_function_interior(1, Box::new(mass_2_comp))?
1329     ///     .input("u", 2, EvalMode::Interp)?
1330     ///     .input("qdata", 1, EvalMode::None)?
1331     ///     .output("v", 2, EvalMode::Interp)?;
1332     ///
1333     /// let op_mass = ceed
1334     ///     .operator(&qf_mass, QFunctionOpt::None, QFunctionOpt::None)?
1335     ///     .field("u", &ru, &bu, VectorOpt::Active)?
1336     ///     .field("qdata", &rq, BasisOpt::None, &qdata)?
1337     ///     .field("v", &ru, &bu, VectorOpt::Active)?;
1338     ///
1339     /// // Diagonal
1340     /// let mut diag = ceed.vector(ncomp * ncomp * ndofs)?;
1341     /// diag.set_value(0.0);
1342     /// op_mass.linear_assemble_point_block_diagonal(&mut diag)?;
1343     ///
1344     /// // Manual diagonal computation
1345     /// let mut true_diag = ceed.vector(ncomp * ncomp * ndofs)?;
1346     /// true_diag.set_value(0.0)?;
1347     /// for i in 0..ndofs {
1348     ///     for j in 0..ncomp {
1349     ///         u.set_value(0.0);
1350     ///         {
1351     ///             let mut u_array = u.view_mut()?;
1352     ///             u_array[i + j * ndofs] = 1.;
1353     ///         }
1354     ///
1355     ///         op_mass.apply(&u, &mut v)?;
1356     ///
1357     ///         {
1358     ///             let v_array = v.view()?;
1359     ///             let mut true_array = true_diag.view_mut()?;
1360     ///             for k in 0..ncomp {
1361     ///                 true_array[i * ncomp * ncomp + k * ncomp + j] = v_array[i + k * ndofs];
1362     ///             }
1363     ///         }
1364     ///     }
1365     /// }
1366     ///
1367     /// // Check
1368     /// diag.view()?
1369     ///     .iter()
1370     ///     .zip(true_diag.view()?.iter())
1371     ///     .for_each(|(computed, actual)| {
1372     ///         assert!(
1373     ///             (*computed - *actual).abs() < 10.0 * libceed::EPSILON,
1374     ///             "Diagonal entry incorrect"
1375     ///         );
1376     ///     });
1377     /// # Ok(())
1378     /// # }
1379     /// ```
1380     pub fn linear_assemble_point_block_diagonal(
1381         &self,
1382         assembled: &mut Vector,
1383     ) -> crate::Result<i32> {
1384         self.op_core.linear_assemble_point_block_diagonal(assembled)
1385     }
1386 
1387     /// Assemble the point block diagonal of a square linear Operator
1388     ///
1389     /// This sums into a Vector with the point block diagonal of a linear
1390     /// Operator.
1391     ///
1392     /// Note: Currently only non-composite Operators with a single field and
1393     ///      composite Operators with single field sub-operators are supported.
1394     ///
1395     /// * `op`        -     Operator to assemble QFunction
1396     /// * `assembled` - Vector to store assembled CeedOperator point block
1397     ///                   diagonal, provided in row-major form with an
1398     ///                   `ncomp * ncomp` block at each node. The dimensions of
1399     ///                   this vector are derived from the active vector for
1400     ///                   the CeedOperator. The array has shape
1401     ///                   `[nodes, component out, component in]`.
1402     ///
1403     /// ```
1404     /// # use libceed::prelude::*;
1405     /// # fn main() -> libceed::Result<()> {
1406     /// # let ceed = libceed::Ceed::default_init();
1407     /// let ne = 4;
1408     /// let p = 3;
1409     /// let q = 4;
1410     /// let ncomp = 2;
1411     /// let ndofs = p * ne - ne + 1;
1412     ///
1413     /// // Vectors
1414     /// let x = ceed.vector_from_slice(&[-1., -0.5, 0.0, 0.5, 1.0])?;
1415     /// let mut qdata = ceed.vector(ne * q)?;
1416     /// qdata.set_value(0.0);
1417     /// let mut u = ceed.vector(ncomp * ndofs)?;
1418     /// u.set_value(1.0);
1419     /// let mut v = ceed.vector(ncomp * ndofs)?;
1420     /// v.set_value(0.0);
1421     ///
1422     /// // Restrictions
1423     /// let mut indx: Vec<i32> = vec![0; 2 * ne];
1424     /// for i in 0..ne {
1425     ///     indx[2 * i + 0] = i as i32;
1426     ///     indx[2 * i + 1] = (i + 1) as i32;
1427     /// }
1428     /// let rx = ceed.elem_restriction(ne, 2, 1, 1, ne + 1, MemType::Host, &indx)?;
1429     /// let mut indu: Vec<i32> = vec![0; p * ne];
1430     /// for i in 0..ne {
1431     ///     indu[p * i + 0] = (2 * i) as i32;
1432     ///     indu[p * i + 1] = (2 * i + 1) as i32;
1433     ///     indu[p * i + 2] = (2 * i + 2) as i32;
1434     /// }
1435     /// let ru = ceed.elem_restriction(ne, p, ncomp, ndofs, ncomp * ndofs, MemType::Host, &indu)?;
1436     /// let strides: [i32; 3] = [1, q as i32, q as i32];
1437     /// let rq = ceed.strided_elem_restriction(ne, q, 1, q * ne, strides)?;
1438     ///
1439     /// // Bases
1440     /// let bx = ceed.basis_tensor_H1_Lagrange(1, 1, 2, q, QuadMode::Gauss)?;
1441     /// let bu = ceed.basis_tensor_H1_Lagrange(1, ncomp, p, q, QuadMode::Gauss)?;
1442     ///
1443     /// // Build quadrature data
1444     /// let qf_build = ceed.q_function_interior_by_name("Mass1DBuild")?;
1445     /// ceed.operator(&qf_build, QFunctionOpt::None, QFunctionOpt::None)?
1446     ///     .field("dx", &rx, &bx, VectorOpt::Active)?
1447     ///     .field("weights", ElemRestrictionOpt::None, &bx, VectorOpt::None)?
1448     ///     .field("qdata", &rq, BasisOpt::None, VectorOpt::Active)?
1449     ///     .apply(&x, &mut qdata)?;
1450     ///
1451     /// // Mass operator
1452     /// let mut mass_2_comp = |[u, qdata, ..]: QFunctionInputs, [v, ..]: QFunctionOutputs| {
1453     ///     // Number of quadrature points
1454     ///     let q = qdata.len();
1455     ///
1456     ///     // Iterate over quadrature points
1457     ///     for i in 0..q {
1458     ///         v[i + 0 * q] = u[i + 1 * q] * qdata[i];
1459     ///         v[i + 1 * q] = u[i + 0 * q] * qdata[i];
1460     ///     }
1461     ///
1462     ///     // Return clean error code
1463     ///     0
1464     /// };
1465     ///
1466     /// let qf_mass = ceed
1467     ///     .q_function_interior(1, Box::new(mass_2_comp))?
1468     ///     .input("u", 2, EvalMode::Interp)?
1469     ///     .input("qdata", 1, EvalMode::None)?
1470     ///     .output("v", 2, EvalMode::Interp)?;
1471     ///
1472     /// let op_mass = ceed
1473     ///     .operator(&qf_mass, QFunctionOpt::None, QFunctionOpt::None)?
1474     ///     .field("u", &ru, &bu, VectorOpt::Active)?
1475     ///     .field("qdata", &rq, BasisOpt::None, &qdata)?
1476     ///     .field("v", &ru, &bu, VectorOpt::Active)?;
1477     ///
1478     /// // Diagonal
1479     /// let mut diag = ceed.vector(ncomp * ncomp * ndofs)?;
1480     /// diag.set_value(1.0);
1481     /// op_mass.linear_assemble_add_point_block_diagonal(&mut diag)?;
1482     ///
1483     /// // Manual diagonal computation
1484     /// let mut true_diag = ceed.vector(ncomp * ncomp * ndofs)?;
1485     /// true_diag.set_value(0.0)?;
1486     /// for i in 0..ndofs {
1487     ///     for j in 0..ncomp {
1488     ///         u.set_value(0.0);
1489     ///         {
1490     ///             let mut u_array = u.view_mut()?;
1491     ///             u_array[i + j * ndofs] = 1.;
1492     ///         }
1493     ///
1494     ///         op_mass.apply(&u, &mut v)?;
1495     ///
1496     ///         {
1497     ///             let v_array = v.view()?;
1498     ///             let mut true_array = true_diag.view_mut()?;
1499     ///             for k in 0..ncomp {
1500     ///                 true_array[i * ncomp * ncomp + k * ncomp + j] = v_array[i + k * ndofs];
1501     ///             }
1502     ///         }
1503     ///     }
1504     /// }
1505     ///
1506     /// // Check
1507     /// diag.view()?
1508     ///     .iter()
1509     ///     .zip(true_diag.view()?.iter())
1510     ///     .for_each(|(computed, actual)| {
1511     ///         assert!(
1512     ///             (*computed - 1.0 - *actual).abs() < 10.0 * libceed::EPSILON,
1513     ///             "Diagonal entry incorrect"
1514     ///         );
1515     ///     });
1516     /// # Ok(())
1517     /// # }
1518     /// ```
1519     pub fn linear_assemble_add_point_block_diagonal(
1520         &self,
1521         assembled: &mut Vector,
1522     ) -> crate::Result<i32> {
1523         self.op_core
1524             .linear_assemble_add_point_block_diagonal(assembled)
1525     }
1526 
1527     /// Create a multigrid coarse Operator and level transfer Operators for a
1528     ///   given Operator, creating the prolongation basis from the fine and
1529     ///   coarse grid interpolation
1530     ///
1531     /// * `p_mult_fine`  - Lvector multiplicity in parallel gather/scatter
1532     /// * `rstr_coarse`  - Coarse grid restriction
1533     /// * `basis_coarse` - Coarse grid active vector basis
1534     ///
1535     /// ```
1536     /// # use libceed::prelude::*;
1537     /// # fn main() -> libceed::Result<()> {
1538     /// # let ceed = libceed::Ceed::default_init();
1539     /// let ne = 15;
1540     /// let p_coarse = 3;
1541     /// let p_fine = 5;
1542     /// let q = 6;
1543     /// let ndofs_coarse = p_coarse * ne - ne + 1;
1544     /// let ndofs_fine = p_fine * ne - ne + 1;
1545     ///
1546     /// // Vectors
1547     /// let x_array = (0..ne + 1)
1548     ///     .map(|i| 2.0 * i as Scalar / ne as Scalar - 1.0)
1549     ///     .collect::<Vec<Scalar>>();
1550     /// let x = ceed.vector_from_slice(&x_array)?;
1551     /// let mut qdata = ceed.vector(ne * q)?;
1552     /// qdata.set_value(0.0);
1553     /// let mut u_coarse = ceed.vector(ndofs_coarse)?;
1554     /// u_coarse.set_value(1.0);
1555     /// let mut u_fine = ceed.vector(ndofs_fine)?;
1556     /// u_fine.set_value(1.0);
1557     /// let mut v_coarse = ceed.vector(ndofs_coarse)?;
1558     /// v_coarse.set_value(0.0);
1559     /// let mut v_fine = ceed.vector(ndofs_fine)?;
1560     /// v_fine.set_value(0.0);
1561     /// let mut multiplicity = ceed.vector(ndofs_fine)?;
1562     /// multiplicity.set_value(1.0);
1563     ///
1564     /// // Restrictions
1565     /// let strides: [i32; 3] = [1, q as i32, q as i32];
1566     /// let rq = ceed.strided_elem_restriction(ne, q, 1, q * ne, strides)?;
1567     ///
1568     /// let mut indx: Vec<i32> = vec![0; 2 * ne];
1569     /// for i in 0..ne {
1570     ///     indx[2 * i + 0] = i as i32;
1571     ///     indx[2 * i + 1] = (i + 1) as i32;
1572     /// }
1573     /// let rx = ceed.elem_restriction(ne, 2, 1, 1, ne + 1, MemType::Host, &indx)?;
1574     ///
1575     /// let mut indu_coarse: Vec<i32> = vec![0; p_coarse * ne];
1576     /// for i in 0..ne {
1577     ///     for j in 0..p_coarse {
1578     ///         indu_coarse[p_coarse * i + j] = (i + j) as i32;
1579     ///     }
1580     /// }
1581     /// let ru_coarse = ceed.elem_restriction(
1582     ///     ne,
1583     ///     p_coarse,
1584     ///     1,
1585     ///     1,
1586     ///     ndofs_coarse,
1587     ///     MemType::Host,
1588     ///     &indu_coarse,
1589     /// )?;
1590     ///
1591     /// let mut indu_fine: Vec<i32> = vec![0; p_fine * ne];
1592     /// for i in 0..ne {
1593     ///     for j in 0..p_fine {
1594     ///         indu_fine[p_fine * i + j] = (i + j) as i32;
1595     ///     }
1596     /// }
1597     /// let ru_fine = ceed.elem_restriction(ne, p_fine, 1, 1, ndofs_fine, MemType::Host, &indu_fine)?;
1598     ///
1599     /// // Bases
1600     /// let bx = ceed.basis_tensor_H1_Lagrange(1, 1, 2, q, QuadMode::Gauss)?;
1601     /// let bu_coarse = ceed.basis_tensor_H1_Lagrange(1, 1, p_coarse, q, QuadMode::Gauss)?;
1602     /// let bu_fine = ceed.basis_tensor_H1_Lagrange(1, 1, p_fine, q, QuadMode::Gauss)?;
1603     ///
1604     /// // Build quadrature data
1605     /// let qf_build = ceed.q_function_interior_by_name("Mass1DBuild")?;
1606     /// ceed.operator(&qf_build, QFunctionOpt::None, QFunctionOpt::None)?
1607     ///     .field("dx", &rx, &bx, VectorOpt::Active)?
1608     ///     .field("weights", ElemRestrictionOpt::None, &bx, VectorOpt::None)?
1609     ///     .field("qdata", &rq, BasisOpt::None, VectorOpt::Active)?
1610     ///     .apply(&x, &mut qdata)?;
1611     ///
1612     /// // Mass operator
1613     /// let qf_mass = ceed.q_function_interior_by_name("MassApply")?;
1614     /// let op_mass_fine = ceed
1615     ///     .operator(&qf_mass, QFunctionOpt::None, QFunctionOpt::None)?
1616     ///     .field("u", &ru_fine, &bu_fine, VectorOpt::Active)?
1617     ///     .field("qdata", &rq, BasisOpt::None, &qdata)?
1618     ///     .field("v", &ru_fine, &bu_fine, VectorOpt::Active)?;
1619     ///
1620     /// // Multigrid setup
1621     /// let (op_mass_coarse, op_prolong, op_restrict) =
1622     ///     op_mass_fine.create_multigrid_level(&multiplicity, &ru_coarse, &bu_coarse)?;
1623     ///
1624     /// // Coarse problem
1625     /// u_coarse.set_value(1.0);
1626     /// op_mass_coarse.apply(&u_coarse, &mut v_coarse)?;
1627     ///
1628     /// // Check
1629     /// let sum: Scalar = v_coarse.view()?.iter().sum();
1630     /// assert!(
1631     ///     (sum - 2.0).abs() < 50.0 * libceed::EPSILON,
1632     ///     "Incorrect interval length computed"
1633     /// );
1634     ///
1635     /// // Prolong
1636     /// op_prolong.apply(&u_coarse, &mut u_fine)?;
1637     ///
1638     /// // Fine problem
1639     /// op_mass_fine.apply(&u_fine, &mut v_fine)?;
1640     ///
1641     /// // Check
1642     /// let sum: Scalar = v_fine.view()?.iter().sum();
1643     /// assert!(
1644     ///     (sum - 2.0).abs() < 50.0 * libceed::EPSILON,
1645     ///     "Incorrect interval length computed"
1646     /// );
1647     ///
1648     /// // Restrict
1649     /// op_restrict.apply(&v_fine, &mut v_coarse)?;
1650     ///
1651     /// // Check
1652     /// let sum: Scalar = v_coarse.view()?.iter().sum();
1653     /// assert!(
1654     ///     (sum - 2.0).abs() < 50.0 * libceed::EPSILON,
1655     ///     "Incorrect interval length computed"
1656     /// );
1657     /// # Ok(())
1658     /// # }
1659     /// ```
1660     pub fn create_multigrid_level<'b>(
1661         &self,
1662         p_mult_fine: &Vector,
1663         rstr_coarse: &ElemRestriction,
1664         basis_coarse: &Basis,
1665     ) -> crate::Result<(Operator<'b>, Operator<'b>, Operator<'b>)> {
1666         let mut ptr_coarse = std::ptr::null_mut();
1667         let mut ptr_prolong = std::ptr::null_mut();
1668         let mut ptr_restrict = std::ptr::null_mut();
1669         let ierr = unsafe {
1670             bind_ceed::CeedOperatorMultigridLevelCreate(
1671                 self.op_core.ptr,
1672                 p_mult_fine.ptr,
1673                 rstr_coarse.ptr,
1674                 basis_coarse.ptr,
1675                 &mut ptr_coarse,
1676                 &mut ptr_prolong,
1677                 &mut ptr_restrict,
1678             )
1679         };
1680         self.op_core.check_error(ierr)?;
1681         let op_coarse = Operator::from_raw(ptr_coarse)?;
1682         let op_prolong = Operator::from_raw(ptr_prolong)?;
1683         let op_restrict = Operator::from_raw(ptr_restrict)?;
1684         Ok((op_coarse, op_prolong, op_restrict))
1685     }
1686 
1687     /// Create a multigrid coarse Operator and level transfer Operators for a
1688     ///   given Operator with a tensor basis for the active basis
1689     ///
1690     /// * `p_mult_fine`   - Lvector multiplicity in parallel gather/scatter
1691     /// * `rstr_coarse`   - Coarse grid restriction
1692     /// * `basis_coarse`  - Coarse grid active vector basis
1693     /// * `interp_c_to_f` - Matrix for coarse to fine
1694     ///
1695     /// ```
1696     /// # use libceed::prelude::*;
1697     /// # fn main() -> libceed::Result<()> {
1698     /// # let ceed = libceed::Ceed::default_init();
1699     /// let ne = 15;
1700     /// let p_coarse = 3;
1701     /// let p_fine = 5;
1702     /// let q = 6;
1703     /// let ndofs_coarse = p_coarse * ne - ne + 1;
1704     /// let ndofs_fine = p_fine * ne - ne + 1;
1705     ///
1706     /// // Vectors
1707     /// let x_array = (0..ne + 1)
1708     ///     .map(|i| 2.0 * i as Scalar / ne as Scalar - 1.0)
1709     ///     .collect::<Vec<Scalar>>();
1710     /// let x = ceed.vector_from_slice(&x_array)?;
1711     /// let mut qdata = ceed.vector(ne * q)?;
1712     /// qdata.set_value(0.0);
1713     /// let mut u_coarse = ceed.vector(ndofs_coarse)?;
1714     /// u_coarse.set_value(1.0);
1715     /// let mut u_fine = ceed.vector(ndofs_fine)?;
1716     /// u_fine.set_value(1.0);
1717     /// let mut v_coarse = ceed.vector(ndofs_coarse)?;
1718     /// v_coarse.set_value(0.0);
1719     /// let mut v_fine = ceed.vector(ndofs_fine)?;
1720     /// v_fine.set_value(0.0);
1721     /// let mut multiplicity = ceed.vector(ndofs_fine)?;
1722     /// multiplicity.set_value(1.0);
1723     ///
1724     /// // Restrictions
1725     /// let strides: [i32; 3] = [1, q as i32, q as i32];
1726     /// let rq = ceed.strided_elem_restriction(ne, q, 1, q * ne, strides)?;
1727     ///
1728     /// let mut indx: Vec<i32> = vec![0; 2 * ne];
1729     /// for i in 0..ne {
1730     ///     indx[2 * i + 0] = i as i32;
1731     ///     indx[2 * i + 1] = (i + 1) as i32;
1732     /// }
1733     /// let rx = ceed.elem_restriction(ne, 2, 1, 1, ne + 1, MemType::Host, &indx)?;
1734     ///
1735     /// let mut indu_coarse: Vec<i32> = vec![0; p_coarse * ne];
1736     /// for i in 0..ne {
1737     ///     for j in 0..p_coarse {
1738     ///         indu_coarse[p_coarse * i + j] = (i + j) as i32;
1739     ///     }
1740     /// }
1741     /// let ru_coarse = ceed.elem_restriction(
1742     ///     ne,
1743     ///     p_coarse,
1744     ///     1,
1745     ///     1,
1746     ///     ndofs_coarse,
1747     ///     MemType::Host,
1748     ///     &indu_coarse,
1749     /// )?;
1750     ///
1751     /// let mut indu_fine: Vec<i32> = vec![0; p_fine * ne];
1752     /// for i in 0..ne {
1753     ///     for j in 0..p_fine {
1754     ///         indu_fine[p_fine * i + j] = (i + j) as i32;
1755     ///     }
1756     /// }
1757     /// let ru_fine = ceed.elem_restriction(ne, p_fine, 1, 1, ndofs_fine, MemType::Host, &indu_fine)?;
1758     ///
1759     /// // Bases
1760     /// let bx = ceed.basis_tensor_H1_Lagrange(1, 1, 2, q, QuadMode::Gauss)?;
1761     /// let bu_coarse = ceed.basis_tensor_H1_Lagrange(1, 1, p_coarse, q, QuadMode::Gauss)?;
1762     /// let bu_fine = ceed.basis_tensor_H1_Lagrange(1, 1, p_fine, q, QuadMode::Gauss)?;
1763     ///
1764     /// // Build quadrature data
1765     /// let qf_build = ceed.q_function_interior_by_name("Mass1DBuild")?;
1766     /// ceed.operator(&qf_build, QFunctionOpt::None, QFunctionOpt::None)?
1767     ///     .field("dx", &rx, &bx, VectorOpt::Active)?
1768     ///     .field("weights", ElemRestrictionOpt::None, &bx, VectorOpt::None)?
1769     ///     .field("qdata", &rq, BasisOpt::None, VectorOpt::Active)?
1770     ///     .apply(&x, &mut qdata)?;
1771     ///
1772     /// // Mass operator
1773     /// let qf_mass = ceed.q_function_interior_by_name("MassApply")?;
1774     /// let op_mass_fine = ceed
1775     ///     .operator(&qf_mass, QFunctionOpt::None, QFunctionOpt::None)?
1776     ///     .field("u", &ru_fine, &bu_fine, VectorOpt::Active)?
1777     ///     .field("qdata", &rq, BasisOpt::None, &qdata)?
1778     ///     .field("v", &ru_fine, &bu_fine, VectorOpt::Active)?;
1779     ///
1780     /// // Multigrid setup
1781     /// let mut interp_c_to_f: Vec<Scalar> = vec![0.; p_coarse * p_fine];
1782     /// {
1783     ///     let mut coarse = ceed.vector(p_coarse)?;
1784     ///     let mut fine = ceed.vector(p_fine)?;
1785     ///     let basis_c_to_f =
1786     ///         ceed.basis_tensor_H1_Lagrange(1, 1, p_coarse, p_fine, QuadMode::GaussLobatto)?;
1787     ///     for i in 0..p_coarse {
1788     ///         coarse.set_value(0.0);
1789     ///         {
1790     ///             let mut array = coarse.view_mut()?;
1791     ///             array[i] = 1.;
1792     ///         }
1793     ///         basis_c_to_f.apply(
1794     ///             1,
1795     ///             TransposeMode::NoTranspose,
1796     ///             EvalMode::Interp,
1797     ///             &coarse,
1798     ///             &mut fine,
1799     ///         )?;
1800     ///         let array = fine.view()?;
1801     ///         for j in 0..p_fine {
1802     ///             interp_c_to_f[j * p_coarse + i] = array[j];
1803     ///         }
1804     ///     }
1805     /// }
1806     /// let (op_mass_coarse, op_prolong, op_restrict) = op_mass_fine.create_multigrid_level_tensor_H1(
1807     ///     &multiplicity,
1808     ///     &ru_coarse,
1809     ///     &bu_coarse,
1810     ///     &interp_c_to_f,
1811     /// )?;
1812     ///
1813     /// // Coarse problem
1814     /// u_coarse.set_value(1.0);
1815     /// op_mass_coarse.apply(&u_coarse, &mut v_coarse)?;
1816     ///
1817     /// // Check
1818     /// let sum: Scalar = v_coarse.view()?.iter().sum();
1819     /// assert!(
1820     ///     (sum - 2.0).abs() < 10.0 * libceed::EPSILON,
1821     ///     "Incorrect interval length computed"
1822     /// );
1823     ///
1824     /// // Prolong
1825     /// op_prolong.apply(&u_coarse, &mut u_fine)?;
1826     ///
1827     /// // Fine problem
1828     /// op_mass_fine.apply(&u_fine, &mut v_fine)?;
1829     ///
1830     /// // Check
1831     /// let sum: Scalar = v_fine.view()?.iter().sum();
1832     /// assert!(
1833     ///     (sum - 2.0).abs() < 10.0 * libceed::EPSILON,
1834     ///     "Incorrect interval length computed"
1835     /// );
1836     ///
1837     /// // Restrict
1838     /// op_restrict.apply(&v_fine, &mut v_coarse)?;
1839     ///
1840     /// // Check
1841     /// let sum: Scalar = v_coarse.view()?.iter().sum();
1842     /// assert!(
1843     ///     (sum - 2.0).abs() < 10.0 * libceed::EPSILON,
1844     ///     "Incorrect interval length computed"
1845     /// );
1846     /// # Ok(())
1847     /// # }
1848     /// ```
1849     pub fn create_multigrid_level_tensor_H1<'b>(
1850         &self,
1851         p_mult_fine: &Vector,
1852         rstr_coarse: &ElemRestriction,
1853         basis_coarse: &Basis,
1854         interpCtoF: &Vec<Scalar>,
1855     ) -> crate::Result<(Operator<'b>, Operator<'b>, Operator<'b>)> {
1856         let mut ptr_coarse = std::ptr::null_mut();
1857         let mut ptr_prolong = std::ptr::null_mut();
1858         let mut ptr_restrict = std::ptr::null_mut();
1859         let ierr = unsafe {
1860             bind_ceed::CeedOperatorMultigridLevelCreateTensorH1(
1861                 self.op_core.ptr,
1862                 p_mult_fine.ptr,
1863                 rstr_coarse.ptr,
1864                 basis_coarse.ptr,
1865                 interpCtoF.as_ptr(),
1866                 &mut ptr_coarse,
1867                 &mut ptr_prolong,
1868                 &mut ptr_restrict,
1869             )
1870         };
1871         self.op_core.check_error(ierr)?;
1872         let op_coarse = Operator::from_raw(ptr_coarse)?;
1873         let op_prolong = Operator::from_raw(ptr_prolong)?;
1874         let op_restrict = Operator::from_raw(ptr_restrict)?;
1875         Ok((op_coarse, op_prolong, op_restrict))
1876     }
1877 
1878     /// Create a multigrid coarse Operator and level transfer Operators for a
1879     ///   given Operator with a non-tensor basis for the active basis
1880     ///
1881     /// * `p_mult_fine`   - Lvector multiplicity in parallel gather/scatter
1882     /// * `rstr_coarse`   - Coarse grid restriction
1883     /// * `basis_coarse`  - Coarse grid active vector basis
1884     /// * `interp_c_to_f` - Matrix for coarse to fine
1885     ///
1886     /// ```
1887     /// # use libceed::prelude::*;
1888     /// # fn main() -> libceed::Result<()> {
1889     /// # let ceed = libceed::Ceed::default_init();
1890     /// let ne = 15;
1891     /// let p_coarse = 3;
1892     /// let p_fine = 5;
1893     /// let q = 6;
1894     /// let ndofs_coarse = p_coarse * ne - ne + 1;
1895     /// let ndofs_fine = p_fine * ne - ne + 1;
1896     ///
1897     /// // Vectors
1898     /// let x_array = (0..ne + 1)
1899     ///     .map(|i| 2.0 * i as Scalar / ne as Scalar - 1.0)
1900     ///     .collect::<Vec<Scalar>>();
1901     /// let x = ceed.vector_from_slice(&x_array)?;
1902     /// let mut qdata = ceed.vector(ne * q)?;
1903     /// qdata.set_value(0.0);
1904     /// let mut u_coarse = ceed.vector(ndofs_coarse)?;
1905     /// u_coarse.set_value(1.0);
1906     /// let mut u_fine = ceed.vector(ndofs_fine)?;
1907     /// u_fine.set_value(1.0);
1908     /// let mut v_coarse = ceed.vector(ndofs_coarse)?;
1909     /// v_coarse.set_value(0.0);
1910     /// let mut v_fine = ceed.vector(ndofs_fine)?;
1911     /// v_fine.set_value(0.0);
1912     /// let mut multiplicity = ceed.vector(ndofs_fine)?;
1913     /// multiplicity.set_value(1.0);
1914     ///
1915     /// // Restrictions
1916     /// let strides: [i32; 3] = [1, q as i32, q as i32];
1917     /// let rq = ceed.strided_elem_restriction(ne, q, 1, q * ne, strides)?;
1918     ///
1919     /// let mut indx: Vec<i32> = vec![0; 2 * ne];
1920     /// for i in 0..ne {
1921     ///     indx[2 * i + 0] = i as i32;
1922     ///     indx[2 * i + 1] = (i + 1) as i32;
1923     /// }
1924     /// let rx = ceed.elem_restriction(ne, 2, 1, 1, ne + 1, MemType::Host, &indx)?;
1925     ///
1926     /// let mut indu_coarse: Vec<i32> = vec![0; p_coarse * ne];
1927     /// for i in 0..ne {
1928     ///     for j in 0..p_coarse {
1929     ///         indu_coarse[p_coarse * i + j] = (i + j) as i32;
1930     ///     }
1931     /// }
1932     /// let ru_coarse = ceed.elem_restriction(
1933     ///     ne,
1934     ///     p_coarse,
1935     ///     1,
1936     ///     1,
1937     ///     ndofs_coarse,
1938     ///     MemType::Host,
1939     ///     &indu_coarse,
1940     /// )?;
1941     ///
1942     /// let mut indu_fine: Vec<i32> = vec![0; p_fine * ne];
1943     /// for i in 0..ne {
1944     ///     for j in 0..p_fine {
1945     ///         indu_fine[p_fine * i + j] = (i + j) as i32;
1946     ///     }
1947     /// }
1948     /// let ru_fine = ceed.elem_restriction(ne, p_fine, 1, 1, ndofs_fine, MemType::Host, &indu_fine)?;
1949     ///
1950     /// // Bases
1951     /// let bx = ceed.basis_tensor_H1_Lagrange(1, 1, 2, q, QuadMode::Gauss)?;
1952     /// let bu_coarse = ceed.basis_tensor_H1_Lagrange(1, 1, p_coarse, q, QuadMode::Gauss)?;
1953     /// let bu_fine = ceed.basis_tensor_H1_Lagrange(1, 1, p_fine, q, QuadMode::Gauss)?;
1954     ///
1955     /// // Build quadrature data
1956     /// let qf_build = ceed.q_function_interior_by_name("Mass1DBuild")?;
1957     /// ceed.operator(&qf_build, QFunctionOpt::None, QFunctionOpt::None)?
1958     ///     .field("dx", &rx, &bx, VectorOpt::Active)?
1959     ///     .field("weights", ElemRestrictionOpt::None, &bx, VectorOpt::None)?
1960     ///     .field("qdata", &rq, BasisOpt::None, VectorOpt::Active)?
1961     ///     .apply(&x, &mut qdata)?;
1962     ///
1963     /// // Mass operator
1964     /// let qf_mass = ceed.q_function_interior_by_name("MassApply")?;
1965     /// let op_mass_fine = ceed
1966     ///     .operator(&qf_mass, QFunctionOpt::None, QFunctionOpt::None)?
1967     ///     .field("u", &ru_fine, &bu_fine, VectorOpt::Active)?
1968     ///     .field("qdata", &rq, BasisOpt::None, &qdata)?
1969     ///     .field("v", &ru_fine, &bu_fine, VectorOpt::Active)?;
1970     ///
1971     /// // Multigrid setup
1972     /// let mut interp_c_to_f: Vec<Scalar> = vec![0.; p_coarse * p_fine];
1973     /// {
1974     ///     let mut coarse = ceed.vector(p_coarse)?;
1975     ///     let mut fine = ceed.vector(p_fine)?;
1976     ///     let basis_c_to_f =
1977     ///         ceed.basis_tensor_H1_Lagrange(1, 1, p_coarse, p_fine, QuadMode::GaussLobatto)?;
1978     ///     for i in 0..p_coarse {
1979     ///         coarse.set_value(0.0);
1980     ///         {
1981     ///             let mut array = coarse.view_mut()?;
1982     ///             array[i] = 1.;
1983     ///         }
1984     ///         basis_c_to_f.apply(
1985     ///             1,
1986     ///             TransposeMode::NoTranspose,
1987     ///             EvalMode::Interp,
1988     ///             &coarse,
1989     ///             &mut fine,
1990     ///         )?;
1991     ///         let array = fine.view()?;
1992     ///         for j in 0..p_fine {
1993     ///             interp_c_to_f[j * p_coarse + i] = array[j];
1994     ///         }
1995     ///     }
1996     /// }
1997     /// let (op_mass_coarse, op_prolong, op_restrict) = op_mass_fine.create_multigrid_level_H1(
1998     ///     &multiplicity,
1999     ///     &ru_coarse,
2000     ///     &bu_coarse,
2001     ///     &interp_c_to_f,
2002     /// )?;
2003     ///
2004     /// // Coarse problem
2005     /// u_coarse.set_value(1.0);
2006     /// op_mass_coarse.apply(&u_coarse, &mut v_coarse)?;
2007     ///
2008     /// // Check
2009     /// let sum: Scalar = v_coarse.view()?.iter().sum();
2010     /// assert!(
2011     ///     (sum - 2.0).abs() < 10.0 * libceed::EPSILON,
2012     ///     "Incorrect interval length computed"
2013     /// );
2014     ///
2015     /// // Prolong
2016     /// op_prolong.apply(&u_coarse, &mut u_fine)?;
2017     ///
2018     /// // Fine problem
2019     /// op_mass_fine.apply(&u_fine, &mut v_fine)?;
2020     ///
2021     /// // Check
2022     /// let sum: Scalar = v_fine.view()?.iter().sum();
2023     /// assert!(
2024     ///     (sum - 2.0).abs() < 10.0 * libceed::EPSILON,
2025     ///     "Incorrect interval length computed"
2026     /// );
2027     ///
2028     /// // Restrict
2029     /// op_restrict.apply(&v_fine, &mut v_coarse)?;
2030     ///
2031     /// // Check
2032     /// let sum: Scalar = v_coarse.view()?.iter().sum();
2033     /// assert!(
2034     ///     (sum - 2.0).abs() < 10.0 * libceed::EPSILON,
2035     ///     "Incorrect interval length computed"
2036     /// );
2037     /// # Ok(())
2038     /// # }
2039     /// ```
2040     pub fn create_multigrid_level_H1<'b>(
2041         &self,
2042         p_mult_fine: &Vector,
2043         rstr_coarse: &ElemRestriction,
2044         basis_coarse: &Basis,
2045         interpCtoF: &[Scalar],
2046     ) -> crate::Result<(Operator<'b>, Operator<'b>, Operator<'b>)> {
2047         let mut ptr_coarse = std::ptr::null_mut();
2048         let mut ptr_prolong = std::ptr::null_mut();
2049         let mut ptr_restrict = std::ptr::null_mut();
2050         let ierr = unsafe {
2051             bind_ceed::CeedOperatorMultigridLevelCreateH1(
2052                 self.op_core.ptr,
2053                 p_mult_fine.ptr,
2054                 rstr_coarse.ptr,
2055                 basis_coarse.ptr,
2056                 interpCtoF.as_ptr(),
2057                 &mut ptr_coarse,
2058                 &mut ptr_prolong,
2059                 &mut ptr_restrict,
2060             )
2061         };
2062         self.op_core.check_error(ierr)?;
2063         let op_coarse = Operator::from_raw(ptr_coarse)?;
2064         let op_prolong = Operator::from_raw(ptr_prolong)?;
2065         let op_restrict = Operator::from_raw(ptr_restrict)?;
2066         Ok((op_coarse, op_prolong, op_restrict))
2067     }
2068 }
2069 
2070 // -----------------------------------------------------------------------------
2071 // Composite Operator
2072 // -----------------------------------------------------------------------------
2073 impl<'a> CompositeOperator<'a> {
2074     // Constructor
2075     pub fn create(ceed: &crate::Ceed) -> crate::Result<Self> {
2076         let mut ptr = std::ptr::null_mut();
2077         let ierr = unsafe { bind_ceed::CeedCompositeOperatorCreate(ceed.ptr, &mut ptr) };
2078         ceed.check_error(ierr)?;
2079         Ok(Self {
2080             op_core: OperatorCore {
2081                 ptr,
2082                 _lifeline: PhantomData,
2083             },
2084         })
2085     }
2086 
2087     /// Set name for CompositeOperator printing
2088     ///
2089     /// * 'name' - Name to set
2090     ///
2091     /// ```
2092     /// # use libceed::prelude::*;
2093     /// # fn main() -> libceed::Result<()> {
2094     /// # let ceed = libceed::Ceed::default_init();
2095     ///
2096     /// // Sub operator field arguments
2097     /// let ne = 3;
2098     /// let q = 4 as usize;
2099     /// let mut ind: Vec<i32> = vec![0; 2 * ne];
2100     /// for i in 0..ne {
2101     ///     ind[2 * i + 0] = i as i32;
2102     ///     ind[2 * i + 1] = (i + 1) as i32;
2103     /// }
2104     /// let r = ceed.elem_restriction(ne, 2, 1, 1, ne + 1, MemType::Host, &ind)?;
2105     /// let strides: [i32; 3] = [1, q as i32, q as i32];
2106     /// let rq = ceed.strided_elem_restriction(ne, q, 1, q * ne, strides)?;
2107     ///
2108     /// let b = ceed.basis_tensor_H1_Lagrange(1, 1, 2, q, QuadMode::Gauss)?;
2109     ///
2110     /// let qdata_mass = ceed.vector(q * ne)?;
2111     /// let qdata_diff = ceed.vector(q * ne)?;
2112     ///
2113     /// let qf_mass = ceed.q_function_interior_by_name("MassApply")?;
2114     /// let op_mass = ceed
2115     ///     .operator(&qf_mass, QFunctionOpt::None, QFunctionOpt::None)?
2116     ///     .name("Mass term")?
2117     ///     .field("u", &r, &b, VectorOpt::Active)?
2118     ///     .field("qdata", &rq, BasisOpt::None, &qdata_mass)?
2119     ///     .field("v", &r, &b, VectorOpt::Active)?;
2120     ///
2121     /// let qf_diff = ceed.q_function_interior_by_name("Poisson1DApply")?;
2122     /// let op_diff = ceed
2123     ///     .operator(&qf_diff, QFunctionOpt::None, QFunctionOpt::None)?
2124     ///     .name("Poisson term")?
2125     ///     .field("du", &r, &b, VectorOpt::Active)?
2126     ///     .field("qdata", &rq, BasisOpt::None, &qdata_diff)?
2127     ///     .field("dv", &r, &b, VectorOpt::Active)?;
2128     ///
2129     /// let op = ceed
2130     ///     .composite_operator()?
2131     ///     .name("Screened Poisson")?
2132     ///     .sub_operator(&op_mass)?
2133     ///     .sub_operator(&op_diff)?;
2134     /// # Ok(())
2135     /// # }
2136     /// ```
2137     #[allow(unused_mut)]
2138     pub fn name(mut self, name: &str) -> crate::Result<Self> {
2139         self.op_core.name(name)?;
2140         Ok(self)
2141     }
2142 
2143     /// Apply Operator to a vector
2144     ///
2145     /// * `input`  - Inpuht Vector
2146     /// * `output` - Output Vector
2147     ///
2148     /// ```
2149     /// # use libceed::prelude::*;
2150     /// # fn main() -> libceed::Result<()> {
2151     /// # let ceed = libceed::Ceed::default_init();
2152     /// let ne = 4;
2153     /// let p = 3;
2154     /// let q = 4;
2155     /// let ndofs = p * ne - ne + 1;
2156     ///
2157     /// // Vectors
2158     /// let x = ceed.vector_from_slice(&[-1., -0.5, 0.0, 0.5, 1.0])?;
2159     /// let mut qdata_mass = ceed.vector(ne * q)?;
2160     /// qdata_mass.set_value(0.0);
2161     /// let mut qdata_diff = ceed.vector(ne * q)?;
2162     /// qdata_diff.set_value(0.0);
2163     /// let mut u = ceed.vector(ndofs)?;
2164     /// u.set_value(1.0);
2165     /// let mut v = ceed.vector(ndofs)?;
2166     /// v.set_value(0.0);
2167     ///
2168     /// // Restrictions
2169     /// let mut indx: Vec<i32> = vec![0; 2 * ne];
2170     /// for i in 0..ne {
2171     ///     indx[2 * i + 0] = i as i32;
2172     ///     indx[2 * i + 1] = (i + 1) as i32;
2173     /// }
2174     /// let rx = ceed.elem_restriction(ne, 2, 1, 1, ne + 1, MemType::Host, &indx)?;
2175     /// let mut indu: Vec<i32> = vec![0; p * ne];
2176     /// for i in 0..ne {
2177     ///     indu[p * i + 0] = i as i32;
2178     ///     indu[p * i + 1] = (i + 1) as i32;
2179     ///     indu[p * i + 2] = (i + 2) as i32;
2180     /// }
2181     /// let ru = ceed.elem_restriction(ne, 3, 1, 1, ndofs, MemType::Host, &indu)?;
2182     /// let strides: [i32; 3] = [1, q as i32, q as i32];
2183     /// let rq = ceed.strided_elem_restriction(ne, q, 1, q * ne, strides)?;
2184     ///
2185     /// // Bases
2186     /// let bx = ceed.basis_tensor_H1_Lagrange(1, 1, 2, q, QuadMode::Gauss)?;
2187     /// let bu = ceed.basis_tensor_H1_Lagrange(1, 1, p, q, QuadMode::Gauss)?;
2188     ///
2189     /// // Build quadrature data
2190     /// let qf_build_mass = ceed.q_function_interior_by_name("Mass1DBuild")?;
2191     /// ceed.operator(&qf_build_mass, QFunctionOpt::None, QFunctionOpt::None)?
2192     ///     .field("dx", &rx, &bx, VectorOpt::Active)?
2193     ///     .field("weights", ElemRestrictionOpt::None, &bx, VectorOpt::None)?
2194     ///     .field("qdata", &rq, BasisOpt::None, VectorOpt::Active)?
2195     ///     .apply(&x, &mut qdata_mass)?;
2196     ///
2197     /// let qf_build_diff = ceed.q_function_interior_by_name("Poisson1DBuild")?;
2198     /// ceed.operator(&qf_build_diff, QFunctionOpt::None, QFunctionOpt::None)?
2199     ///     .field("dx", &rx, &bx, VectorOpt::Active)?
2200     ///     .field("weights", ElemRestrictionOpt::None, &bx, VectorOpt::None)?
2201     ///     .field("qdata", &rq, BasisOpt::None, VectorOpt::Active)?
2202     ///     .apply(&x, &mut qdata_diff)?;
2203     ///
2204     /// // Application operator
2205     /// let qf_mass = ceed.q_function_interior_by_name("MassApply")?;
2206     /// let op_mass = ceed
2207     ///     .operator(&qf_mass, QFunctionOpt::None, QFunctionOpt::None)?
2208     ///     .field("u", &ru, &bu, VectorOpt::Active)?
2209     ///     .field("qdata", &rq, BasisOpt::None, &qdata_mass)?
2210     ///     .field("v", &ru, &bu, VectorOpt::Active)?;
2211     ///
2212     /// let qf_diff = ceed.q_function_interior_by_name("Poisson1DApply")?;
2213     /// let op_diff = ceed
2214     ///     .operator(&qf_diff, QFunctionOpt::None, QFunctionOpt::None)?
2215     ///     .field("du", &ru, &bu, VectorOpt::Active)?
2216     ///     .field("qdata", &rq, BasisOpt::None, &qdata_diff)?
2217     ///     .field("dv", &ru, &bu, VectorOpt::Active)?;
2218     ///
2219     /// let op_composite = ceed
2220     ///     .composite_operator()?
2221     ///     .sub_operator(&op_mass)?
2222     ///     .sub_operator(&op_diff)?;
2223     ///
2224     /// v.set_value(0.0);
2225     /// op_composite.apply(&u, &mut v)?;
2226     ///
2227     /// // Check
2228     /// let sum: Scalar = v.view()?.iter().sum();
2229     /// assert!(
2230     ///     (sum - 2.0).abs() < 10.0 * libceed::EPSILON,
2231     ///     "Incorrect interval length computed"
2232     /// );
2233     /// # Ok(())
2234     /// # }
2235     /// ```
2236     pub fn apply(&self, input: &Vector, output: &mut Vector) -> crate::Result<i32> {
2237         self.op_core.apply(input, output)
2238     }
2239 
2240     /// Apply Operator to a vector and add result to output vector
2241     ///
2242     /// * `input`  - Input Vector
2243     /// * `output` - Output Vector
2244     ///
2245     /// ```
2246     /// # use libceed::prelude::*;
2247     /// # fn main() -> libceed::Result<()> {
2248     /// # let ceed = libceed::Ceed::default_init();
2249     /// let ne = 4;
2250     /// let p = 3;
2251     /// let q = 4;
2252     /// let ndofs = p * ne - ne + 1;
2253     ///
2254     /// // Vectors
2255     /// let x = ceed.vector_from_slice(&[-1., -0.5, 0.0, 0.5, 1.0])?;
2256     /// let mut qdata_mass = ceed.vector(ne * q)?;
2257     /// qdata_mass.set_value(0.0);
2258     /// let mut qdata_diff = ceed.vector(ne * q)?;
2259     /// qdata_diff.set_value(0.0);
2260     /// let mut u = ceed.vector(ndofs)?;
2261     /// u.set_value(1.0);
2262     /// let mut v = ceed.vector(ndofs)?;
2263     /// v.set_value(0.0);
2264     ///
2265     /// // Restrictions
2266     /// let mut indx: Vec<i32> = vec![0; 2 * ne];
2267     /// for i in 0..ne {
2268     ///     indx[2 * i + 0] = i as i32;
2269     ///     indx[2 * i + 1] = (i + 1) as i32;
2270     /// }
2271     /// let rx = ceed.elem_restriction(ne, 2, 1, 1, ne + 1, MemType::Host, &indx)?;
2272     /// let mut indu: Vec<i32> = vec![0; p * ne];
2273     /// for i in 0..ne {
2274     ///     indu[p * i + 0] = i as i32;
2275     ///     indu[p * i + 1] = (i + 1) as i32;
2276     ///     indu[p * i + 2] = (i + 2) as i32;
2277     /// }
2278     /// let ru = ceed.elem_restriction(ne, 3, 1, 1, ndofs, MemType::Host, &indu)?;
2279     /// let strides: [i32; 3] = [1, q as i32, q as i32];
2280     /// let rq = ceed.strided_elem_restriction(ne, q, 1, q * ne, strides)?;
2281     ///
2282     /// // Bases
2283     /// let bx = ceed.basis_tensor_H1_Lagrange(1, 1, 2, q, QuadMode::Gauss)?;
2284     /// let bu = ceed.basis_tensor_H1_Lagrange(1, 1, p, q, QuadMode::Gauss)?;
2285     ///
2286     /// // Build quadrature data
2287     /// let qf_build_mass = ceed.q_function_interior_by_name("Mass1DBuild")?;
2288     /// ceed.operator(&qf_build_mass, QFunctionOpt::None, QFunctionOpt::None)?
2289     ///     .field("dx", &rx, &bx, VectorOpt::Active)?
2290     ///     .field("weights", ElemRestrictionOpt::None, &bx, VectorOpt::None)?
2291     ///     .field("qdata", &rq, BasisOpt::None, VectorOpt::Active)?
2292     ///     .apply(&x, &mut qdata_mass)?;
2293     ///
2294     /// let qf_build_diff = ceed.q_function_interior_by_name("Poisson1DBuild")?;
2295     /// ceed.operator(&qf_build_diff, QFunctionOpt::None, QFunctionOpt::None)?
2296     ///     .field("dx", &rx, &bx, VectorOpt::Active)?
2297     ///     .field("weights", ElemRestrictionOpt::None, &bx, VectorOpt::None)?
2298     ///     .field("qdata", &rq, BasisOpt::None, VectorOpt::Active)?
2299     ///     .apply(&x, &mut qdata_diff)?;
2300     ///
2301     /// // Application operator
2302     /// let qf_mass = ceed.q_function_interior_by_name("MassApply")?;
2303     /// let op_mass = ceed
2304     ///     .operator(&qf_mass, QFunctionOpt::None, QFunctionOpt::None)?
2305     ///     .field("u", &ru, &bu, VectorOpt::Active)?
2306     ///     .field("qdata", &rq, BasisOpt::None, &qdata_mass)?
2307     ///     .field("v", &ru, &bu, VectorOpt::Active)?;
2308     ///
2309     /// let qf_diff = ceed.q_function_interior_by_name("Poisson1DApply")?;
2310     /// let op_diff = ceed
2311     ///     .operator(&qf_diff, QFunctionOpt::None, QFunctionOpt::None)?
2312     ///     .field("du", &ru, &bu, VectorOpt::Active)?
2313     ///     .field("qdata", &rq, BasisOpt::None, &qdata_diff)?
2314     ///     .field("dv", &ru, &bu, VectorOpt::Active)?;
2315     ///
2316     /// let op_composite = ceed
2317     ///     .composite_operator()?
2318     ///     .sub_operator(&op_mass)?
2319     ///     .sub_operator(&op_diff)?;
2320     ///
2321     /// v.set_value(1.0);
2322     /// op_composite.apply_add(&u, &mut v)?;
2323     ///
2324     /// // Check
2325     /// let sum: Scalar = v.view()?.iter().sum();
2326     /// assert!(
2327     ///     (sum - (2.0 + ndofs as Scalar)).abs() < 10.0 * libceed::EPSILON,
2328     ///     "Incorrect interval length computed"
2329     /// );
2330     /// # Ok(())
2331     /// # }
2332     /// ```
2333     pub fn apply_add(&self, input: &Vector, output: &mut Vector) -> crate::Result<i32> {
2334         self.op_core.apply_add(input, output)
2335     }
2336 
2337     /// Add a sub-Operator to a Composite Operator
2338     ///
2339     /// * `subop` - Sub-Operator
2340     ///
2341     /// ```
2342     /// # use libceed::prelude::*;
2343     /// # fn main() -> libceed::Result<()> {
2344     /// # let ceed = libceed::Ceed::default_init();
2345     /// let mut op = ceed.composite_operator()?;
2346     ///
2347     /// let qf_mass = ceed.q_function_interior_by_name("MassApply")?;
2348     /// let op_mass = ceed.operator(&qf_mass, QFunctionOpt::None, QFunctionOpt::None)?;
2349     /// op = op.sub_operator(&op_mass)?;
2350     ///
2351     /// let qf_diff = ceed.q_function_interior_by_name("Poisson1DApply")?;
2352     /// let op_diff = ceed.operator(&qf_diff, QFunctionOpt::None, QFunctionOpt::None)?;
2353     /// op = op.sub_operator(&op_diff)?;
2354     /// # Ok(())
2355     /// # }
2356     /// ```
2357     #[allow(unused_mut)]
2358     pub fn sub_operator(mut self, subop: &Operator) -> crate::Result<Self> {
2359         let ierr =
2360             unsafe { bind_ceed::CeedCompositeOperatorAddSub(self.op_core.ptr, subop.op_core.ptr) };
2361         self.op_core.check_error(ierr)?;
2362         Ok(self)
2363     }
2364 
2365     /// Check if CompositeOperator is setup correctly
2366     ///
2367     /// ```
2368     /// # use libceed::prelude::*;
2369     /// # fn main() -> libceed::Result<()> {
2370     /// # let ceed = libceed::Ceed::default_init();
2371     /// let ne = 4;
2372     /// let p = 3;
2373     /// let q = 4;
2374     /// let ndofs = p * ne - ne + 1;
2375     ///
2376     /// // Restrictions
2377     /// let mut indx: Vec<i32> = vec![0; 2 * ne];
2378     /// for i in 0..ne {
2379     ///     indx[2 * i + 0] = i as i32;
2380     ///     indx[2 * i + 1] = (i + 1) as i32;
2381     /// }
2382     /// let rx = ceed.elem_restriction(ne, 2, 1, 1, ne + 1, MemType::Host, &indx)?;
2383     /// let strides: [i32; 3] = [1, q as i32, q as i32];
2384     /// let rq = ceed.strided_elem_restriction(ne, q, 1, q * ne, strides)?;
2385     ///
2386     /// // Bases
2387     /// let bx = ceed.basis_tensor_H1_Lagrange(1, 1, 2, q, QuadMode::Gauss)?;
2388     ///
2389     /// // Build quadrature data
2390     /// let qf_build_mass = ceed.q_function_interior_by_name("Mass1DBuild")?;
2391     /// let op_build_mass = ceed
2392     ///     .operator(&qf_build_mass, QFunctionOpt::None, QFunctionOpt::None)?
2393     ///     .field("dx", &rx, &bx, VectorOpt::Active)?
2394     ///     .field("weights", ElemRestrictionOpt::None, &bx, VectorOpt::None)?
2395     ///     .field("qdata", &rq, BasisOpt::None, VectorOpt::Active)?;
2396     ///
2397     /// let qf_build_diff = ceed.q_function_interior_by_name("Poisson1DBuild")?;
2398     /// let op_build_diff = ceed
2399     ///     .operator(&qf_build_diff, QFunctionOpt::None, QFunctionOpt::None)?
2400     ///     .field("dx", &rx, &bx, VectorOpt::Active)?
2401     ///     .field("weights", ElemRestrictionOpt::None, &bx, VectorOpt::None)?
2402     ///     .field("qdata", &rq, BasisOpt::None, VectorOpt::Active)?;
2403     ///
2404     /// let op_build = ceed
2405     ///     .composite_operator()?
2406     ///     .sub_operator(&op_build_mass)?
2407     ///     .sub_operator(&op_build_diff)?;
2408     ///
2409     /// // Check
2410     /// assert!(op_build.check().is_ok(), "Operator setup incorrectly");
2411     /// # Ok(())
2412     /// # }
2413     /// ```
2414     pub fn check(self) -> crate::Result<Self> {
2415         self.op_core.check()?;
2416         Ok(self)
2417     }
2418 
2419     /// Assemble the diagonal of a square linear Operator
2420     ///
2421     /// This overwrites a Vector with the diagonal of a linear Operator.
2422     ///
2423     /// Note: Currently only non-composite Operators with a single field and
2424     ///      composite Operators with single field sub-operators are supported.
2425     ///
2426     /// * `op`        - Operator to assemble QFunction
2427     /// * `assembled` - Vector to store assembled Operator diagonal
2428     pub fn linear_assemble_diagonal(&self, assembled: &mut Vector) -> crate::Result<i32> {
2429         self.op_core.linear_assemble_add_diagonal(assembled)
2430     }
2431 
2432     /// Assemble the point block diagonal of a square linear Operator
2433     ///
2434     /// This overwrites a Vector with the point block diagonal of a linear
2435     ///   Operator.
2436     ///
2437     /// Note: Currently only non-composite Operators with a single field and
2438     ///      composite Operators with single field sub-operators are supported.
2439     ///
2440     /// * `op`        - Operator to assemble QFunction
2441     /// * `assembled` - Vector to store assembled Operator diagonal
2442     pub fn linear_assemble_add_diagonal(&self, assembled: &mut Vector) -> crate::Result<i32> {
2443         self.op_core.linear_assemble_add_diagonal(assembled)
2444     }
2445 
2446     /// Assemble the diagonal of a square linear Operator
2447     ///
2448     /// This overwrites a Vector with the diagonal of a linear Operator.
2449     ///
2450     /// Note: Currently only non-composite Operators with a single field and
2451     ///      composite Operators with single field sub-operators are supported.
2452     ///
2453     /// * `op`        - Operator to assemble QFunction
2454     /// * `assembled` - Vector to store assembled CeedOperator point block
2455     ///                   diagonal, provided in row-major form with an
2456     ///                   `ncomp * ncomp` block at each node. The dimensions of
2457     ///                   this vector are derived from the active vector for
2458     ///                   the CeedOperator. The array has shape
2459     ///                   `[nodes, component out, component in]`.
2460     pub fn linear_assemble_point_block_diagonal(
2461         &self,
2462         assembled: &mut Vector,
2463     ) -> crate::Result<i32> {
2464         self.op_core.linear_assemble_add_diagonal(assembled)
2465     }
2466 
2467     /// Assemble the diagonal of a square linear Operator
2468     ///
2469     /// This sums into a Vector with the diagonal of a linear Operator.
2470     ///
2471     /// Note: Currently only non-composite Operators with a single field and
2472     ///      composite Operators with single field sub-operators are supported.
2473     ///
2474     /// * `op`        - Operator to assemble QFunction
2475     /// * `assembled` - Vector to store assembled CeedOperator point block
2476     ///                   diagonal, provided in row-major form with an
2477     ///                   `ncomp * ncomp` block at each node. The dimensions of
2478     ///                   this vector are derived from the active vector for
2479     ///                   the CeedOperator. The array has shape
2480     ///                   `[nodes, component out, component in]`.
2481     pub fn linear_assemble_add_point_block_diagonal(
2482         &self,
2483         assembled: &mut Vector,
2484     ) -> crate::Result<i32> {
2485         self.op_core.linear_assemble_add_diagonal(assembled)
2486     }
2487 }
2488 
2489 // -----------------------------------------------------------------------------
2490