xref: /libCEED/rust/libceed/src/qfunction.rs (revision 08778c6feee59ca4815589c6a74184f38d46595f)
19df49d7eSJed Brown // Copyright (c) 2017, Lawrence Livermore National Security, LLC. Produced at
29df49d7eSJed Brown // the Lawrence Livermore National Laboratory. LLNL-CODE-734707. All Rights
39df49d7eSJed Brown // reserved. See files LICENSE and NOTICE for details.
49df49d7eSJed Brown //
59df49d7eSJed Brown // This file is part of CEED, a collection of benchmarks, miniapps, software
69df49d7eSJed Brown // libraries and APIs for efficient high-order finite element and spectral
79df49d7eSJed Brown // element discretizations for exascale applications. For more information and
89df49d7eSJed Brown // source code availability see http://github.com/ceed.
99df49d7eSJed Brown //
109df49d7eSJed Brown // The CEED research is supported by the Exascale Computing Project 17-SC-20-SC,
119df49d7eSJed Brown // a collaborative effort of two U.S. Department of Energy organizations (Office
129df49d7eSJed Brown // of Science and the National Nuclear Security Administration) responsible for
139df49d7eSJed Brown // the planning and preparation of a capable exascale ecosystem, including
149df49d7eSJed Brown // software, applications, hardware, advanced system engineering and early
159df49d7eSJed Brown // testbed platforms, in support of the nation's exascale computing imperative
169df49d7eSJed Brown 
179df49d7eSJed Brown //! A Ceed QFunction represents the spatial terms of the point-wise functions
189df49d7eSJed Brown //! describing the physics at the quadrature points.
199df49d7eSJed Brown 
209df49d7eSJed Brown use std::pin::Pin;
219df49d7eSJed Brown 
229df49d7eSJed Brown use crate::prelude::*;
239df49d7eSJed Brown 
2480a9ef05SNatalie Beams pub type QFunctionInputs<'a> = [&'a [crate::Scalar]; MAX_QFUNCTION_FIELDS];
2580a9ef05SNatalie Beams pub type QFunctionOutputs<'a> = [&'a mut [crate::Scalar]; MAX_QFUNCTION_FIELDS];
269df49d7eSJed Brown 
279df49d7eSJed Brown // -----------------------------------------------------------------------------
28*08778c6fSJeremy L Thompson // CeedQFunction Field context wrapper
29*08778c6fSJeremy L Thompson // -----------------------------------------------------------------------------
30*08778c6fSJeremy L Thompson #[derive(Debug)]
31*08778c6fSJeremy L Thompson pub struct QFunctionField<'a> {
32*08778c6fSJeremy L Thompson     ptr: bind_ceed::CeedQFunctionField,
33*08778c6fSJeremy L Thompson     _lifeline: PhantomData<&'a ()>,
34*08778c6fSJeremy L Thompson }
35*08778c6fSJeremy L Thompson 
36*08778c6fSJeremy L Thompson // -----------------------------------------------------------------------------
37*08778c6fSJeremy L Thompson // Implementations
38*08778c6fSJeremy L Thompson // -----------------------------------------------------------------------------
39*08778c6fSJeremy L Thompson impl<'a> QFunctionField<'a> {
40*08778c6fSJeremy L Thompson     /// Get the name of a QFunctionField
41*08778c6fSJeremy L Thompson     ///
42*08778c6fSJeremy L Thompson     /// ```
43*08778c6fSJeremy L Thompson     /// # use libceed::prelude::*;
44*08778c6fSJeremy L Thompson     /// # fn main() -> Result<(), libceed::CeedError> {
45*08778c6fSJeremy L Thompson     /// # let ceed = libceed::Ceed::default_init();
46*08778c6fSJeremy L Thompson     /// const Q: usize = 8;
47*08778c6fSJeremy L Thompson     /// let qf = ceed.q_function_interior_by_name("Mass2DBuild")?;
48*08778c6fSJeremy L Thompson     ///
49*08778c6fSJeremy L Thompson     /// let inputs = qf.inputs()?;
50*08778c6fSJeremy L Thompson     ///
51*08778c6fSJeremy L Thompson     /// assert_eq!(inputs[0].name(), "dx", "Incorrect input name");
52*08778c6fSJeremy L Thompson     /// assert_eq!(inputs[1].name(), "weights", "Incorrect input name");
53*08778c6fSJeremy L Thompson     /// # Ok(())
54*08778c6fSJeremy L Thompson     /// # }
55*08778c6fSJeremy L Thompson     /// ```
56*08778c6fSJeremy L Thompson     pub fn name(&self) -> &str {
57*08778c6fSJeremy L Thompson         let mut name_ptr: *mut std::os::raw::c_char = std::ptr::null_mut();
58*08778c6fSJeremy L Thompson         unsafe {
59*08778c6fSJeremy L Thompson             bind_ceed::CeedQFunctionFieldGetName(self.ptr, &mut name_ptr);
60*08778c6fSJeremy L Thompson         }
61*08778c6fSJeremy L Thompson         unsafe { CStr::from_ptr(name_ptr) }.to_str().unwrap()
62*08778c6fSJeremy L Thompson     }
63*08778c6fSJeremy L Thompson 
64*08778c6fSJeremy L Thompson     /// Get the size of a QFunctionField
65*08778c6fSJeremy L Thompson     ///
66*08778c6fSJeremy L Thompson     /// ```
67*08778c6fSJeremy L Thompson     /// # use libceed::prelude::*;
68*08778c6fSJeremy L Thompson     /// # fn main() -> Result<(), libceed::CeedError> {
69*08778c6fSJeremy L Thompson     /// # let ceed = libceed::Ceed::default_init();
70*08778c6fSJeremy L Thompson     /// const Q: usize = 8;
71*08778c6fSJeremy L Thompson     /// let qf = ceed.q_function_interior_by_name("Mass2DBuild")?;
72*08778c6fSJeremy L Thompson     ///
73*08778c6fSJeremy L Thompson     /// let inputs = qf.inputs()?;
74*08778c6fSJeremy L Thompson     ///
75*08778c6fSJeremy L Thompson     /// assert_eq!(inputs[0].size(), 4, "Incorrect input size");
76*08778c6fSJeremy L Thompson     /// assert_eq!(inputs[1].size(), 1, "Incorrect input size");
77*08778c6fSJeremy L Thompson     /// # Ok(())
78*08778c6fSJeremy L Thompson     /// # }
79*08778c6fSJeremy L Thompson     /// ```
80*08778c6fSJeremy L Thompson     pub fn size(&self) -> usize {
81*08778c6fSJeremy L Thompson         let mut size = 0;
82*08778c6fSJeremy L Thompson         unsafe {
83*08778c6fSJeremy L Thompson             bind_ceed::CeedQFunctionFieldGetSize(self.ptr, &mut size);
84*08778c6fSJeremy L Thompson         }
85*08778c6fSJeremy L Thompson         usize::try_from(size).unwrap()
86*08778c6fSJeremy L Thompson     }
87*08778c6fSJeremy L Thompson 
88*08778c6fSJeremy L Thompson     /// Get the evaluation mode of a QFunctionField
89*08778c6fSJeremy L Thompson     ///
90*08778c6fSJeremy L Thompson     /// ```
91*08778c6fSJeremy L Thompson     /// # use libceed::prelude::*;
92*08778c6fSJeremy L Thompson     /// # fn main() -> Result<(), libceed::CeedError> {
93*08778c6fSJeremy L Thompson     /// # let ceed = libceed::Ceed::default_init();
94*08778c6fSJeremy L Thompson     /// const Q: usize = 8;
95*08778c6fSJeremy L Thompson     /// let qf = ceed.q_function_interior_by_name("Mass2DBuild")?;
96*08778c6fSJeremy L Thompson     ///
97*08778c6fSJeremy L Thompson     /// let inputs = qf.inputs()?;
98*08778c6fSJeremy L Thompson     ///
99*08778c6fSJeremy L Thompson     /// assert_eq!(
100*08778c6fSJeremy L Thompson     ///     inputs[0].eval_mode(),
101*08778c6fSJeremy L Thompson     ///     EvalMode::Grad,
102*08778c6fSJeremy L Thompson     ///     "Incorrect input evaluation mode"
103*08778c6fSJeremy L Thompson     /// );
104*08778c6fSJeremy L Thompson     /// assert_eq!(
105*08778c6fSJeremy L Thompson     ///     inputs[1].eval_mode(),
106*08778c6fSJeremy L Thompson     ///     EvalMode::Weight,
107*08778c6fSJeremy L Thompson     ///     "Incorrect input evaluation mode"
108*08778c6fSJeremy L Thompson     /// );
109*08778c6fSJeremy L Thompson     /// # Ok(())
110*08778c6fSJeremy L Thompson     /// # }
111*08778c6fSJeremy L Thompson     /// ```
112*08778c6fSJeremy L Thompson     pub fn eval_mode(&self) -> crate::EvalMode {
113*08778c6fSJeremy L Thompson         let mut mode = 0;
114*08778c6fSJeremy L Thompson         unsafe {
115*08778c6fSJeremy L Thompson             bind_ceed::CeedQFunctionFieldGetEvalMode(self.ptr, &mut mode);
116*08778c6fSJeremy L Thompson         }
117*08778c6fSJeremy L Thompson         crate::EvalMode::from_u32(mode as u32)
118*08778c6fSJeremy L Thompson     }
119*08778c6fSJeremy L Thompson }
120*08778c6fSJeremy L Thompson 
121*08778c6fSJeremy L Thompson // -----------------------------------------------------------------------------
1229df49d7eSJed Brown // CeedQFunction option
1239df49d7eSJed Brown // -----------------------------------------------------------------------------
1249df49d7eSJed Brown pub enum QFunctionOpt<'a> {
1259df49d7eSJed Brown     SomeQFunction(&'a QFunction<'a>),
1269df49d7eSJed Brown     SomeQFunctionByName(&'a QFunctionByName<'a>),
1279df49d7eSJed Brown     None,
1289df49d7eSJed Brown }
1299df49d7eSJed Brown 
1309df49d7eSJed Brown /// Construct a QFunctionOpt reference from a QFunction reference
1319df49d7eSJed Brown impl<'a> From<&'a QFunction<'_>> for QFunctionOpt<'a> {
1329df49d7eSJed Brown     fn from(qfunc: &'a QFunction) -> Self {
1339df49d7eSJed Brown         debug_assert!(qfunc.qf_core.ptr != unsafe { bind_ceed::CEED_QFUNCTION_NONE });
1349df49d7eSJed Brown         Self::SomeQFunction(qfunc)
1359df49d7eSJed Brown     }
1369df49d7eSJed Brown }
1379df49d7eSJed Brown 
1389df49d7eSJed Brown /// Construct a QFunctionOpt reference from a QFunction by Name reference
1399df49d7eSJed Brown impl<'a> From<&'a QFunctionByName<'_>> for QFunctionOpt<'a> {
1409df49d7eSJed Brown     fn from(qfunc: &'a QFunctionByName) -> Self {
1419df49d7eSJed Brown         debug_assert!(qfunc.qf_core.ptr != unsafe { bind_ceed::CEED_QFUNCTION_NONE });
1429df49d7eSJed Brown         Self::SomeQFunctionByName(qfunc)
1439df49d7eSJed Brown     }
1449df49d7eSJed Brown }
1459df49d7eSJed Brown 
1469df49d7eSJed Brown impl<'a> QFunctionOpt<'a> {
1479df49d7eSJed Brown     /// Transform a Rust libCEED QFunctionOpt into C libCEED CeedQFunction
1489df49d7eSJed Brown     pub(crate) fn to_raw(self) -> bind_ceed::CeedQFunction {
1499df49d7eSJed Brown         match self {
1509df49d7eSJed Brown             Self::SomeQFunction(qfunc) => qfunc.qf_core.ptr,
1519df49d7eSJed Brown             Self::SomeQFunctionByName(qfunc) => qfunc.qf_core.ptr,
1529df49d7eSJed Brown             Self::None => unsafe { bind_ceed::CEED_QFUNCTION_NONE },
1539df49d7eSJed Brown         }
1549df49d7eSJed Brown     }
155e03682afSJeremy L Thompson 
156e03682afSJeremy L Thompson     /// Check if a QFunctionOpt is Some
157e03682afSJeremy L Thompson     ///
158e03682afSJeremy L Thompson     /// ```
159e03682afSJeremy L Thompson     /// # use libceed::prelude::*;
160e03682afSJeremy L Thompson     /// # fn main() -> Result<(), libceed::CeedError> {
161e03682afSJeremy L Thompson     /// # let ceed = libceed::Ceed::default_init();
162e03682afSJeremy L Thompson     /// let mut user_f = |[u, weights, ..]: QFunctionInputs, [v, ..]: QFunctionOutputs| {
163e03682afSJeremy L Thompson     ///     // Iterate over quadrature points
164e03682afSJeremy L Thompson     ///     v.iter_mut()
165e03682afSJeremy L Thompson     ///         .zip(u.iter().zip(weights.iter()))
166e03682afSJeremy L Thompson     ///         .for_each(|(v, (u, w))| *v = u * w);
167e03682afSJeremy L Thompson     ///
168e03682afSJeremy L Thompson     ///     // Return clean error code
169e03682afSJeremy L Thompson     ///     0
170e03682afSJeremy L Thompson     /// };
171e03682afSJeremy L Thompson     ///
172e03682afSJeremy L Thompson     /// let qf = ceed
173e03682afSJeremy L Thompson     ///     .q_function_interior(1, Box::new(user_f))?
174e03682afSJeremy L Thompson     ///     .input("u", 1, EvalMode::Interp)?
175e03682afSJeremy L Thompson     ///     .input("weights", 1, EvalMode::Weight)?
176e03682afSJeremy L Thompson     ///     .output("v", 1, EvalMode::Interp)?;
177e03682afSJeremy L Thompson     /// let qf_opt = QFunctionOpt::from(&qf);
178e03682afSJeremy L Thompson     /// assert!(qf_opt.is_some(), "Incorrect QFunctionOpt");
179e03682afSJeremy L Thompson     ///
180e03682afSJeremy L Thompson     /// let qf = ceed.q_function_interior_by_name("Mass1DBuild")?;
181e03682afSJeremy L Thompson     /// let qf_opt = QFunctionOpt::from(&qf);
182e03682afSJeremy L Thompson     /// assert!(qf_opt.is_some(), "Incorrect QFunctionOpt");
183e03682afSJeremy L Thompson     ///
184e03682afSJeremy L Thompson     /// let qf_opt = QFunctionOpt::None;
185e03682afSJeremy L Thompson     /// assert!(!qf_opt.is_some(), "Incorrect QFunctionOpt");
186e03682afSJeremy L Thompson     /// # Ok(())
187e03682afSJeremy L Thompson     /// # }
188e03682afSJeremy L Thompson     /// ```
189e03682afSJeremy L Thompson     pub fn is_some(&self) -> bool {
190e03682afSJeremy L Thompson         match self {
191e03682afSJeremy L Thompson             Self::SomeQFunction(_) => true,
192e03682afSJeremy L Thompson             Self::SomeQFunctionByName(_) => true,
193e03682afSJeremy L Thompson             Self::None => false,
194e03682afSJeremy L Thompson         }
195e03682afSJeremy L Thompson     }
196e03682afSJeremy L Thompson 
197e03682afSJeremy L Thompson     /// Check if a QFunctionOpt is SomeQFunction
198e03682afSJeremy L Thompson     ///
199e03682afSJeremy L Thompson     /// ```
200e03682afSJeremy L Thompson     /// # use libceed::prelude::*;
201e03682afSJeremy L Thompson     /// # fn main() -> Result<(), libceed::CeedError> {
202e03682afSJeremy L Thompson     /// # let ceed = libceed::Ceed::default_init();
203e03682afSJeremy L Thompson     /// let mut user_f = |[u, weights, ..]: QFunctionInputs, [v, ..]: QFunctionOutputs| {
204e03682afSJeremy L Thompson     ///     // Iterate over quadrature points
205e03682afSJeremy L Thompson     ///     v.iter_mut()
206e03682afSJeremy L Thompson     ///         .zip(u.iter().zip(weights.iter()))
207e03682afSJeremy L Thompson     ///         .for_each(|(v, (u, w))| *v = u * w);
208e03682afSJeremy L Thompson     ///
209e03682afSJeremy L Thompson     ///     // Return clean error code
210e03682afSJeremy L Thompson     ///     0
211e03682afSJeremy L Thompson     /// };
212e03682afSJeremy L Thompson     ///
213e03682afSJeremy L Thompson     /// let qf = ceed
214e03682afSJeremy L Thompson     ///     .q_function_interior(1, Box::new(user_f))?
215e03682afSJeremy L Thompson     ///     .input("u", 1, EvalMode::Interp)?
216e03682afSJeremy L Thompson     ///     .input("weights", 1, EvalMode::Weight)?
217e03682afSJeremy L Thompson     ///     .output("v", 1, EvalMode::Interp)?;
218e03682afSJeremy L Thompson     /// let qf_opt = QFunctionOpt::from(&qf);
219e03682afSJeremy L Thompson     /// assert!(qf_opt.is_some_q_function(), "Incorrect QFunctionOpt");
220e03682afSJeremy L Thompson     ///
221e03682afSJeremy L Thompson     /// let qf = ceed.q_function_interior_by_name("Mass1DBuild")?;
222e03682afSJeremy L Thompson     /// let qf_opt = QFunctionOpt::from(&qf);
223e03682afSJeremy L Thompson     /// assert!(!qf_opt.is_some_q_function(), "Incorrect QFunctionOpt");
224e03682afSJeremy L Thompson     ///
225e03682afSJeremy L Thompson     /// let qf_opt = QFunctionOpt::None;
226e03682afSJeremy L Thompson     /// assert!(!qf_opt.is_some_q_function(), "Incorrect QFunctionOpt");
227e03682afSJeremy L Thompson     /// # Ok(())
228e03682afSJeremy L Thompson     /// # }
229e03682afSJeremy L Thompson     /// ```
230e03682afSJeremy L Thompson     pub fn is_some_q_function(&self) -> bool {
231e03682afSJeremy L Thompson         match self {
232e03682afSJeremy L Thompson             Self::SomeQFunction(_) => true,
233e03682afSJeremy L Thompson             Self::SomeQFunctionByName(_) => false,
234e03682afSJeremy L Thompson             Self::None => false,
235e03682afSJeremy L Thompson         }
236e03682afSJeremy L Thompson     }
237e03682afSJeremy L Thompson 
238e03682afSJeremy L Thompson     /// Check if a QFunctionOpt is SomeQFunctionByName
239e03682afSJeremy L Thompson     ///
240e03682afSJeremy L Thompson     /// ```
241e03682afSJeremy L Thompson     /// # use libceed::prelude::*;
242e03682afSJeremy L Thompson     /// # fn main() -> Result<(), libceed::CeedError> {
243e03682afSJeremy L Thompson     /// # let ceed = libceed::Ceed::default_init();
244e03682afSJeremy L Thompson     /// let mut user_f = |[u, weights, ..]: QFunctionInputs, [v, ..]: QFunctionOutputs| {
245e03682afSJeremy L Thompson     ///     // Iterate over quadrature points
246e03682afSJeremy L Thompson     ///     v.iter_mut()
247e03682afSJeremy L Thompson     ///         .zip(u.iter().zip(weights.iter()))
248e03682afSJeremy L Thompson     ///         .for_each(|(v, (u, w))| *v = u * w);
249e03682afSJeremy L Thompson     ///
250e03682afSJeremy L Thompson     ///     // Return clean error code
251e03682afSJeremy L Thompson     ///     0
252e03682afSJeremy L Thompson     /// };
253e03682afSJeremy L Thompson     ///
254e03682afSJeremy L Thompson     /// let qf = ceed
255e03682afSJeremy L Thompson     ///     .q_function_interior(1, Box::new(user_f))?
256e03682afSJeremy L Thompson     ///     .input("u", 1, EvalMode::Interp)?
257e03682afSJeremy L Thompson     ///     .input("weights", 1, EvalMode::Weight)?
258e03682afSJeremy L Thompson     ///     .output("v", 1, EvalMode::Interp)?;
259e03682afSJeremy L Thompson     /// let qf_opt = QFunctionOpt::from(&qf);
260e03682afSJeremy L Thompson     /// assert!(
261e03682afSJeremy L Thompson     ///     !qf_opt.is_some_q_function_by_name(),
262e03682afSJeremy L Thompson     ///     "Incorrect QFunctionOpt"
263e03682afSJeremy L Thompson     /// );
264e03682afSJeremy L Thompson     ///
265e03682afSJeremy L Thompson     /// let qf = ceed.q_function_interior_by_name("Mass1DBuild")?;
266e03682afSJeremy L Thompson     /// let qf_opt = QFunctionOpt::from(&qf);
267e03682afSJeremy L Thompson     /// assert!(
268e03682afSJeremy L Thompson     ///     qf_opt.is_some_q_function_by_name(),
269e03682afSJeremy L Thompson     ///     "Incorrect QFunctionOpt"
270e03682afSJeremy L Thompson     /// );
271e03682afSJeremy L Thompson     ///
272e03682afSJeremy L Thompson     /// let qf_opt = QFunctionOpt::None;
273e03682afSJeremy L Thompson     /// assert!(
274e03682afSJeremy L Thompson     ///     !qf_opt.is_some_q_function_by_name(),
275e03682afSJeremy L Thompson     ///     "Incorrect QFunctionOpt"
276e03682afSJeremy L Thompson     /// );
277e03682afSJeremy L Thompson     /// # Ok(())
278e03682afSJeremy L Thompson     /// # }
279e03682afSJeremy L Thompson     /// ```
280e03682afSJeremy L Thompson     pub fn is_some_q_function_by_name(&self) -> bool {
281e03682afSJeremy L Thompson         match self {
282e03682afSJeremy L Thompson             Self::SomeQFunction(_) => false,
283e03682afSJeremy L Thompson             Self::SomeQFunctionByName(_) => true,
284e03682afSJeremy L Thompson             Self::None => false,
285e03682afSJeremy L Thompson         }
286e03682afSJeremy L Thompson     }
287e03682afSJeremy L Thompson 
288e03682afSJeremy L Thompson     /// Check if a QFunctionOpt is None
289e03682afSJeremy L Thompson     ///
290e03682afSJeremy L Thompson     /// ```
291e03682afSJeremy L Thompson     /// # use libceed::prelude::*;
292e03682afSJeremy L Thompson     /// # fn main() -> Result<(), libceed::CeedError> {
293e03682afSJeremy L Thompson     /// # let ceed = libceed::Ceed::default_init();
294e03682afSJeremy L Thompson     /// let mut user_f = |[u, weights, ..]: QFunctionInputs, [v, ..]: QFunctionOutputs| {
295e03682afSJeremy L Thompson     ///     // Iterate over quadrature points
296e03682afSJeremy L Thompson     ///     v.iter_mut()
297e03682afSJeremy L Thompson     ///         .zip(u.iter().zip(weights.iter()))
298e03682afSJeremy L Thompson     ///         .for_each(|(v, (u, w))| *v = u * w);
299e03682afSJeremy L Thompson     ///
300e03682afSJeremy L Thompson     ///     // Return clean error code
301e03682afSJeremy L Thompson     ///     0
302e03682afSJeremy L Thompson     /// };
303e03682afSJeremy L Thompson     ///
304e03682afSJeremy L Thompson     /// let qf = ceed
305e03682afSJeremy L Thompson     ///     .q_function_interior(1, Box::new(user_f))?
306e03682afSJeremy L Thompson     ///     .input("u", 1, EvalMode::Interp)?
307e03682afSJeremy L Thompson     ///     .input("weights", 1, EvalMode::Weight)?
308e03682afSJeremy L Thompson     ///     .output("v", 1, EvalMode::Interp)?;
309e03682afSJeremy L Thompson     /// let qf_opt = QFunctionOpt::from(&qf);
310e03682afSJeremy L Thompson     /// assert!(!qf_opt.is_none(), "Incorrect QFunctionOpt");
311e03682afSJeremy L Thompson     ///
312e03682afSJeremy L Thompson     /// let qf = ceed.q_function_interior_by_name("Mass1DBuild")?;
313e03682afSJeremy L Thompson     /// let qf_opt = QFunctionOpt::from(&qf);
314e03682afSJeremy L Thompson     /// assert!(!qf_opt.is_none(), "Incorrect QFunctionOpt");
315e03682afSJeremy L Thompson     ///
316e03682afSJeremy L Thompson     /// let qf_opt = QFunctionOpt::None;
317e03682afSJeremy L Thompson     /// assert!(qf_opt.is_none(), "Incorrect QFunctionOpt");
318e03682afSJeremy L Thompson     /// # Ok(())
319e03682afSJeremy L Thompson     /// # }
320e03682afSJeremy L Thompson     /// ```
321e03682afSJeremy L Thompson     pub fn is_none(&self) -> bool {
322e03682afSJeremy L Thompson         match self {
323e03682afSJeremy L Thompson             Self::SomeQFunction(_) => false,
324e03682afSJeremy L Thompson             Self::SomeQFunctionByName(_) => false,
325e03682afSJeremy L Thompson             Self::None => true,
326e03682afSJeremy L Thompson         }
327e03682afSJeremy L Thompson     }
3289df49d7eSJed Brown }
3299df49d7eSJed Brown 
3309df49d7eSJed Brown // -----------------------------------------------------------------------------
3319df49d7eSJed Brown // CeedQFunction context wrapper
3329df49d7eSJed Brown // -----------------------------------------------------------------------------
333c68be7a2SJeremy L Thompson #[derive(Debug)]
3349df49d7eSJed Brown pub(crate) struct QFunctionCore<'a> {
3359df49d7eSJed Brown     ptr: bind_ceed::CeedQFunction,
3361142270cSJeremy L Thompson     _lifeline: PhantomData<&'a ()>,
3379df49d7eSJed Brown }
3389df49d7eSJed Brown 
3399df49d7eSJed Brown struct QFunctionTrampolineData {
3409df49d7eSJed Brown     number_inputs: usize,
3419df49d7eSJed Brown     number_outputs: usize,
3429df49d7eSJed Brown     input_sizes: [usize; MAX_QFUNCTION_FIELDS],
3439df49d7eSJed Brown     output_sizes: [usize; MAX_QFUNCTION_FIELDS],
3449df49d7eSJed Brown     user_f: Box<QFunctionUserClosure>,
3459df49d7eSJed Brown }
3469df49d7eSJed Brown 
3479df49d7eSJed Brown pub struct QFunction<'a> {
3489df49d7eSJed Brown     qf_core: QFunctionCore<'a>,
3499df49d7eSJed Brown     qf_ctx_ptr: bind_ceed::CeedQFunctionContext,
3509df49d7eSJed Brown     trampoline_data: Pin<Box<QFunctionTrampolineData>>,
3519df49d7eSJed Brown }
3529df49d7eSJed Brown 
353c68be7a2SJeremy L Thompson #[derive(Debug)]
3549df49d7eSJed Brown pub struct QFunctionByName<'a> {
3559df49d7eSJed Brown     qf_core: QFunctionCore<'a>,
3569df49d7eSJed Brown }
3579df49d7eSJed Brown 
3589df49d7eSJed Brown // -----------------------------------------------------------------------------
3599df49d7eSJed Brown // Destructor
3609df49d7eSJed Brown // -----------------------------------------------------------------------------
3619df49d7eSJed Brown impl<'a> Drop for QFunctionCore<'a> {
3629df49d7eSJed Brown     fn drop(&mut self) {
3639df49d7eSJed Brown         unsafe {
3649df49d7eSJed Brown             if self.ptr != bind_ceed::CEED_QFUNCTION_NONE {
3659df49d7eSJed Brown                 bind_ceed::CeedQFunctionDestroy(&mut self.ptr);
3669df49d7eSJed Brown             }
3679df49d7eSJed Brown         }
3689df49d7eSJed Brown     }
3699df49d7eSJed Brown }
3709df49d7eSJed Brown 
3719df49d7eSJed Brown impl<'a> Drop for QFunction<'a> {
3729df49d7eSJed Brown     fn drop(&mut self) {
3739df49d7eSJed Brown         unsafe {
3749df49d7eSJed Brown             bind_ceed::CeedQFunctionContextDestroy(&mut self.qf_ctx_ptr);
3759df49d7eSJed Brown         }
3769df49d7eSJed Brown     }
3779df49d7eSJed Brown }
3789df49d7eSJed Brown 
3799df49d7eSJed Brown // -----------------------------------------------------------------------------
3809df49d7eSJed Brown // Display
3819df49d7eSJed Brown // -----------------------------------------------------------------------------
3829df49d7eSJed Brown impl<'a> fmt::Display for QFunctionCore<'a> {
3839df49d7eSJed Brown     fn fmt(&self, f: &mut fmt::Formatter) -> fmt::Result {
3849df49d7eSJed Brown         let mut ptr = std::ptr::null_mut();
3859df49d7eSJed Brown         let mut sizeloc = crate::MAX_BUFFER_LENGTH;
3869df49d7eSJed Brown         let cstring = unsafe {
3879df49d7eSJed Brown             let file = bind_ceed::open_memstream(&mut ptr, &mut sizeloc);
3889df49d7eSJed Brown             bind_ceed::CeedQFunctionView(self.ptr, file);
3899df49d7eSJed Brown             bind_ceed::fclose(file);
3909df49d7eSJed Brown             CString::from_raw(ptr)
3919df49d7eSJed Brown         };
3929df49d7eSJed Brown         cstring.to_string_lossy().fmt(f)
3939df49d7eSJed Brown     }
3949df49d7eSJed Brown }
3959df49d7eSJed Brown /// View a QFunction
3969df49d7eSJed Brown ///
3979df49d7eSJed Brown /// ```
3989df49d7eSJed Brown /// # use libceed::prelude::*;
399c68be7a2SJeremy L Thompson /// # fn main() -> Result<(), libceed::CeedError> {
4009df49d7eSJed Brown /// # let ceed = libceed::Ceed::default_init();
4019df49d7eSJed Brown /// let mut user_f = |[u, weights, ..]: QFunctionInputs, [v, ..]: QFunctionOutputs| {
4029df49d7eSJed Brown ///     // Iterate over quadrature points
4039df49d7eSJed Brown ///     v.iter_mut()
4049df49d7eSJed Brown ///         .zip(u.iter().zip(weights.iter()))
4059df49d7eSJed Brown ///         .for_each(|(v, (u, w))| *v = u * w);
4069df49d7eSJed Brown ///
4079df49d7eSJed Brown ///     // Return clean error code
4089df49d7eSJed Brown ///     0
4099df49d7eSJed Brown /// };
4109df49d7eSJed Brown ///
4119df49d7eSJed Brown /// let qf = ceed
412c68be7a2SJeremy L Thompson ///     .q_function_interior(1, Box::new(user_f))?
413c68be7a2SJeremy L Thompson ///     .input("u", 1, EvalMode::Interp)?
414c68be7a2SJeremy L Thompson ///     .input("weights", 1, EvalMode::Weight)?
415c68be7a2SJeremy L Thompson ///     .output("v", 1, EvalMode::Interp)?;
4169df49d7eSJed Brown ///
4179df49d7eSJed Brown /// println!("{}", qf);
418c68be7a2SJeremy L Thompson /// # Ok(())
419c68be7a2SJeremy L Thompson /// # }
4209df49d7eSJed Brown /// ```
4219df49d7eSJed Brown impl<'a> fmt::Display for QFunction<'a> {
4229df49d7eSJed Brown     fn fmt(&self, f: &mut fmt::Formatter) -> fmt::Result {
4239df49d7eSJed Brown         self.qf_core.fmt(f)
4249df49d7eSJed Brown     }
4259df49d7eSJed Brown }
4269df49d7eSJed Brown 
4279df49d7eSJed Brown /// View a QFunction by Name
4289df49d7eSJed Brown ///
4299df49d7eSJed Brown /// ```
4309df49d7eSJed Brown /// # use libceed::prelude::*;
431c68be7a2SJeremy L Thompson /// # fn main() -> Result<(), libceed::CeedError> {
4329df49d7eSJed Brown /// # let ceed = libceed::Ceed::default_init();
433c68be7a2SJeremy L Thompson /// let qf = ceed.q_function_interior_by_name("Mass1DBuild")?;
4349df49d7eSJed Brown /// println!("{}", qf);
435c68be7a2SJeremy L Thompson /// # Ok(())
436c68be7a2SJeremy L Thompson /// # }
4379df49d7eSJed Brown /// ```
4389df49d7eSJed Brown impl<'a> fmt::Display for QFunctionByName<'a> {
4399df49d7eSJed Brown     fn fmt(&self, f: &mut fmt::Formatter) -> fmt::Result {
4409df49d7eSJed Brown         self.qf_core.fmt(f)
4419df49d7eSJed Brown     }
4429df49d7eSJed Brown }
4439df49d7eSJed Brown 
4449df49d7eSJed Brown // -----------------------------------------------------------------------------
4459df49d7eSJed Brown // Core functionality
4469df49d7eSJed Brown // -----------------------------------------------------------------------------
4479df49d7eSJed Brown impl<'a> QFunctionCore<'a> {
4481142270cSJeremy L Thompson     // Error handling
4491142270cSJeremy L Thompson     #[doc(hidden)]
4501142270cSJeremy L Thompson     fn check_error(&self, ierr: i32) -> crate::Result<i32> {
4511142270cSJeremy L Thompson         let mut ptr = std::ptr::null_mut();
4521142270cSJeremy L Thompson         unsafe {
4531142270cSJeremy L Thompson             bind_ceed::CeedQFunctionGetCeed(self.ptr, &mut ptr);
4541142270cSJeremy L Thompson         }
4551142270cSJeremy L Thompson         crate::check_error(ptr, ierr)
4561142270cSJeremy L Thompson     }
4571142270cSJeremy L Thompson 
4589df49d7eSJed Brown     // Common implementation
4599df49d7eSJed Brown     pub fn apply(&self, Q: usize, u: &[Vector], v: &[Vector]) -> crate::Result<i32> {
4609df49d7eSJed Brown         let mut u_c = [std::ptr::null_mut(); MAX_QFUNCTION_FIELDS];
4619df49d7eSJed Brown         for i in 0..std::cmp::min(MAX_QFUNCTION_FIELDS, u.len()) {
4629df49d7eSJed Brown             u_c[i] = u[i].ptr;
4639df49d7eSJed Brown         }
4649df49d7eSJed Brown         let mut v_c = [std::ptr::null_mut(); MAX_QFUNCTION_FIELDS];
4659df49d7eSJed Brown         for i in 0..std::cmp::min(MAX_QFUNCTION_FIELDS, v.len()) {
4669df49d7eSJed Brown             v_c[i] = v[i].ptr;
4679df49d7eSJed Brown         }
4689df49d7eSJed Brown         let Q = i32::try_from(Q).unwrap();
4699df49d7eSJed Brown         let ierr = unsafe {
4709df49d7eSJed Brown             bind_ceed::CeedQFunctionApply(self.ptr, Q, u_c.as_mut_ptr(), v_c.as_mut_ptr())
4719df49d7eSJed Brown         };
4721142270cSJeremy L Thompson         self.check_error(ierr)
4739df49d7eSJed Brown     }
474*08778c6fSJeremy L Thompson 
475*08778c6fSJeremy L Thompson     pub fn inputs(&self) -> crate::Result<&[crate::QFunctionField]> {
476*08778c6fSJeremy L Thompson         // Get array of raw C pointers for inputs
477*08778c6fSJeremy L Thompson         let mut num_inputs = 0;
478*08778c6fSJeremy L Thompson         let mut inputs_ptr = std::ptr::null_mut();
479*08778c6fSJeremy L Thompson         let ierr = unsafe {
480*08778c6fSJeremy L Thompson             bind_ceed::CeedQFunctionGetFields(
481*08778c6fSJeremy L Thompson                 self.ptr,
482*08778c6fSJeremy L Thompson                 &mut num_inputs,
483*08778c6fSJeremy L Thompson                 &mut inputs_ptr,
484*08778c6fSJeremy L Thompson                 std::ptr::null_mut() as *mut bind_ceed::CeedInt,
485*08778c6fSJeremy L Thompson                 std::ptr::null_mut() as *mut *mut bind_ceed::CeedQFunctionField,
486*08778c6fSJeremy L Thompson             )
487*08778c6fSJeremy L Thompson         };
488*08778c6fSJeremy L Thompson         self.check_error(ierr)?;
489*08778c6fSJeremy L Thompson         // Convert raw C pointers to fixed length slice
490*08778c6fSJeremy L Thompson         let inputs_slice = unsafe {
491*08778c6fSJeremy L Thompson             std::slice::from_raw_parts(
492*08778c6fSJeremy L Thompson                 inputs_ptr as *const crate::QFunctionField,
493*08778c6fSJeremy L Thompson                 num_inputs as usize,
494*08778c6fSJeremy L Thompson             )
495*08778c6fSJeremy L Thompson         };
496*08778c6fSJeremy L Thompson         Ok(inputs_slice)
497*08778c6fSJeremy L Thompson     }
498*08778c6fSJeremy L Thompson 
499*08778c6fSJeremy L Thompson     pub fn outputs(&self) -> crate::Result<&[crate::QFunctionField]> {
500*08778c6fSJeremy L Thompson         // Get array of raw C pointers for outputs
501*08778c6fSJeremy L Thompson         let mut num_outputs = 0;
502*08778c6fSJeremy L Thompson         let mut outputs_ptr = std::ptr::null_mut();
503*08778c6fSJeremy L Thompson         let ierr = unsafe {
504*08778c6fSJeremy L Thompson             bind_ceed::CeedQFunctionGetFields(
505*08778c6fSJeremy L Thompson                 self.ptr,
506*08778c6fSJeremy L Thompson                 std::ptr::null_mut() as *mut bind_ceed::CeedInt,
507*08778c6fSJeremy L Thompson                 std::ptr::null_mut() as *mut *mut bind_ceed::CeedQFunctionField,
508*08778c6fSJeremy L Thompson                 &mut num_outputs,
509*08778c6fSJeremy L Thompson                 &mut outputs_ptr,
510*08778c6fSJeremy L Thompson             )
511*08778c6fSJeremy L Thompson         };
512*08778c6fSJeremy L Thompson         self.check_error(ierr)?;
513*08778c6fSJeremy L Thompson         // Convert raw C pointers to fixed length slice
514*08778c6fSJeremy L Thompson         let outputs_slice = unsafe {
515*08778c6fSJeremy L Thompson             std::slice::from_raw_parts(
516*08778c6fSJeremy L Thompson                 outputs_ptr as *const crate::QFunctionField,
517*08778c6fSJeremy L Thompson                 num_outputs as usize,
518*08778c6fSJeremy L Thompson             )
519*08778c6fSJeremy L Thompson         };
520*08778c6fSJeremy L Thompson         Ok(outputs_slice)
521*08778c6fSJeremy L Thompson     }
5229df49d7eSJed Brown }
5239df49d7eSJed Brown 
5249df49d7eSJed Brown // -----------------------------------------------------------------------------
5259df49d7eSJed Brown // User QFunction Closure
5269df49d7eSJed Brown // -----------------------------------------------------------------------------
52780a9ef05SNatalie Beams pub type QFunctionUserClosure = dyn FnMut(
52880a9ef05SNatalie Beams     [&[crate::Scalar]; MAX_QFUNCTION_FIELDS],
52980a9ef05SNatalie Beams     [&mut [crate::Scalar]; MAX_QFUNCTION_FIELDS],
53080a9ef05SNatalie Beams ) -> i32;
5319df49d7eSJed Brown 
5329df49d7eSJed Brown macro_rules! mut_max_fields {
5339df49d7eSJed Brown     ($e:expr) => {
5349df49d7eSJed Brown         [
5359df49d7eSJed Brown             $e, $e, $e, $e, $e, $e, $e, $e, $e, $e, $e, $e, $e, $e, $e, $e,
5369df49d7eSJed Brown         ]
5379df49d7eSJed Brown     };
5389df49d7eSJed Brown }
5399df49d7eSJed Brown unsafe extern "C" fn trampoline(
5409df49d7eSJed Brown     ctx: *mut ::std::os::raw::c_void,
5419df49d7eSJed Brown     q: bind_ceed::CeedInt,
5429df49d7eSJed Brown     inputs: *const *const bind_ceed::CeedScalar,
5439df49d7eSJed Brown     outputs: *const *mut bind_ceed::CeedScalar,
5449df49d7eSJed Brown ) -> ::std::os::raw::c_int {
5459df49d7eSJed Brown     let trampoline_data: Pin<&mut QFunctionTrampolineData> = std::mem::transmute(ctx);
5469df49d7eSJed Brown 
5479df49d7eSJed Brown     // Inputs
5489df49d7eSJed Brown     let inputs_slice: &[*const bind_ceed::CeedScalar] =
5499df49d7eSJed Brown         std::slice::from_raw_parts(inputs, MAX_QFUNCTION_FIELDS);
55080a9ef05SNatalie Beams     let mut inputs_array: [&[crate::Scalar]; MAX_QFUNCTION_FIELDS] = [&[0.0]; MAX_QFUNCTION_FIELDS];
5519df49d7eSJed Brown     inputs_slice
5529df49d7eSJed Brown         .iter()
5539df49d7eSJed Brown         .enumerate()
5549df49d7eSJed Brown         .map(|(i, &x)| {
55580a9ef05SNatalie Beams             std::slice::from_raw_parts(x, trampoline_data.input_sizes[i] * q as usize)
55680a9ef05SNatalie Beams                 as &[crate::Scalar]
5579df49d7eSJed Brown         })
5589df49d7eSJed Brown         .zip(inputs_array.iter_mut())
5599df49d7eSJed Brown         .for_each(|(x, a)| *a = x);
5609df49d7eSJed Brown 
5619df49d7eSJed Brown     // Outputs
5629df49d7eSJed Brown     let outputs_slice: &[*mut bind_ceed::CeedScalar] =
5639df49d7eSJed Brown         std::slice::from_raw_parts(outputs, MAX_QFUNCTION_FIELDS);
56480a9ef05SNatalie Beams     let mut outputs_array: [&mut [crate::Scalar]; MAX_QFUNCTION_FIELDS] =
56580a9ef05SNatalie Beams         mut_max_fields!(&mut [0.0]);
5669df49d7eSJed Brown     outputs_slice
5679df49d7eSJed Brown         .iter()
5689df49d7eSJed Brown         .enumerate()
5699df49d7eSJed Brown         .map(|(i, &x)| {
5709df49d7eSJed Brown             std::slice::from_raw_parts_mut(x, trampoline_data.output_sizes[i] * q as usize)
57180a9ef05SNatalie Beams                 as &mut [crate::Scalar]
5729df49d7eSJed Brown         })
5739df49d7eSJed Brown         .zip(outputs_array.iter_mut())
5749df49d7eSJed Brown         .for_each(|(x, a)| *a = x);
5759df49d7eSJed Brown 
5769df49d7eSJed Brown     // User closure
5779df49d7eSJed Brown     (trampoline_data.get_unchecked_mut().user_f)(inputs_array, outputs_array)
5789df49d7eSJed Brown }
5799df49d7eSJed Brown 
5809df49d7eSJed Brown // -----------------------------------------------------------------------------
5819df49d7eSJed Brown // QFunction
5829df49d7eSJed Brown // -----------------------------------------------------------------------------
5839df49d7eSJed Brown impl<'a> QFunction<'a> {
5849df49d7eSJed Brown     // Constructor
5859df49d7eSJed Brown     pub fn create(
5869df49d7eSJed Brown         ceed: &'a crate::Ceed,
5879df49d7eSJed Brown         vlength: usize,
5889df49d7eSJed Brown         user_f: Box<QFunctionUserClosure>,
5899df49d7eSJed Brown     ) -> crate::Result<Self> {
5909df49d7eSJed Brown         let source_c = CString::new("").expect("CString::new failed");
5919df49d7eSJed Brown         let mut ptr = std::ptr::null_mut();
5929df49d7eSJed Brown 
5939df49d7eSJed Brown         // Context for closure
5949df49d7eSJed Brown         let number_inputs = 0;
5959df49d7eSJed Brown         let number_outputs = 0;
5969df49d7eSJed Brown         let input_sizes = [0; MAX_QFUNCTION_FIELDS];
5979df49d7eSJed Brown         let output_sizes = [0; MAX_QFUNCTION_FIELDS];
5989df49d7eSJed Brown         let trampoline_data = unsafe {
5999df49d7eSJed Brown             Pin::new_unchecked(Box::new(QFunctionTrampolineData {
6009df49d7eSJed Brown                 number_inputs,
6019df49d7eSJed Brown                 number_outputs,
6029df49d7eSJed Brown                 input_sizes,
6039df49d7eSJed Brown                 output_sizes,
6049df49d7eSJed Brown                 user_f,
6059df49d7eSJed Brown             }))
6069df49d7eSJed Brown         };
6079df49d7eSJed Brown 
6089df49d7eSJed Brown         // Create QFunction
6099df49d7eSJed Brown         let vlength = i32::try_from(vlength).unwrap();
6109df49d7eSJed Brown         let mut ierr = unsafe {
6119df49d7eSJed Brown             bind_ceed::CeedQFunctionCreateInterior(
6129df49d7eSJed Brown                 ceed.ptr,
6139df49d7eSJed Brown                 vlength,
6149df49d7eSJed Brown                 Some(trampoline),
6159df49d7eSJed Brown                 source_c.as_ptr(),
6169df49d7eSJed Brown                 &mut ptr,
6179df49d7eSJed Brown             )
6189df49d7eSJed Brown         };
6199df49d7eSJed Brown         ceed.check_error(ierr)?;
6209df49d7eSJed Brown 
6219df49d7eSJed Brown         // Set closure
6229df49d7eSJed Brown         let mut qf_ctx_ptr = std::ptr::null_mut();
6239df49d7eSJed Brown         ierr = unsafe { bind_ceed::CeedQFunctionContextCreate(ceed.ptr, &mut qf_ctx_ptr) };
6249df49d7eSJed Brown         ceed.check_error(ierr)?;
6259df49d7eSJed Brown         ierr = unsafe {
6269df49d7eSJed Brown             bind_ceed::CeedQFunctionContextSetData(
6279df49d7eSJed Brown                 qf_ctx_ptr,
6289df49d7eSJed Brown                 crate::MemType::Host as bind_ceed::CeedMemType,
6299df49d7eSJed Brown                 crate::CopyMode::UsePointer as bind_ceed::CeedCopyMode,
6309df49d7eSJed Brown                 std::mem::size_of::<QFunctionTrampolineData>() as u64,
6319df49d7eSJed Brown                 std::mem::transmute(trampoline_data.as_ref()),
6329df49d7eSJed Brown             )
6339df49d7eSJed Brown         };
6349df49d7eSJed Brown         ceed.check_error(ierr)?;
6359df49d7eSJed Brown         ierr = unsafe { bind_ceed::CeedQFunctionSetContext(ptr, qf_ctx_ptr) };
6369df49d7eSJed Brown         ceed.check_error(ierr)?;
6379df49d7eSJed Brown         Ok(Self {
6381142270cSJeremy L Thompson             qf_core: QFunctionCore {
6391142270cSJeremy L Thompson                 ptr,
6401142270cSJeremy L Thompson                 _lifeline: PhantomData,
6411142270cSJeremy L Thompson             },
6429df49d7eSJed Brown             qf_ctx_ptr,
6439df49d7eSJed Brown             trampoline_data,
6449df49d7eSJed Brown         })
6459df49d7eSJed Brown     }
6469df49d7eSJed Brown 
6479df49d7eSJed Brown     /// Apply the action of a QFunction
6489df49d7eSJed Brown     ///
6499df49d7eSJed Brown     /// * `Q`      - The number of quadrature points
6509df49d7eSJed Brown     /// * `input`  - Array of input Vectors
6519df49d7eSJed Brown     /// * `output` - Array of output Vectors
6529df49d7eSJed Brown     ///
6539df49d7eSJed Brown     /// ```
6549df49d7eSJed Brown     /// # use libceed::prelude::*;
655c68be7a2SJeremy L Thompson     /// # fn main() -> Result<(), libceed::CeedError> {
6569df49d7eSJed Brown     /// # let ceed = libceed::Ceed::default_init();
6579df49d7eSJed Brown     /// let mut user_f = |[u, weights, ..]: QFunctionInputs, [v, ..]: QFunctionOutputs| {
6589df49d7eSJed Brown     ///     // Iterate over quadrature points
6599df49d7eSJed Brown     ///     v.iter_mut()
6609df49d7eSJed Brown     ///         .zip(u.iter().zip(weights.iter()))
6619df49d7eSJed Brown     ///         .for_each(|(v, (u, w))| *v = u * w);
6629df49d7eSJed Brown     ///
6639df49d7eSJed Brown     ///     // Return clean error code
6649df49d7eSJed Brown     ///     0
6659df49d7eSJed Brown     /// };
6669df49d7eSJed Brown     ///
6679df49d7eSJed Brown     /// let qf = ceed
668c68be7a2SJeremy L Thompson     ///     .q_function_interior(1, Box::new(user_f))?
669c68be7a2SJeremy L Thompson     ///     .input("u", 1, EvalMode::Interp)?
670c68be7a2SJeremy L Thompson     ///     .input("weights", 1, EvalMode::Weight)?
671c68be7a2SJeremy L Thompson     ///     .output("v", 1, EvalMode::Interp)?;
6729df49d7eSJed Brown     ///
6739df49d7eSJed Brown     /// const Q: usize = 8;
6749df49d7eSJed Brown     /// let mut w = [0.; Q];
6759df49d7eSJed Brown     /// let mut u = [0.; Q];
6769df49d7eSJed Brown     /// let mut v = [0.; Q];
6779df49d7eSJed Brown     ///
6789df49d7eSJed Brown     /// for i in 0..Q {
67980a9ef05SNatalie Beams     ///     let x = 2. * (i as Scalar) / ((Q as Scalar) - 1.) - 1.;
6809df49d7eSJed Brown     ///     u[i] = 2. + 3. * x + 5. * x * x;
6819df49d7eSJed Brown     ///     w[i] = 1. - x * x;
6829df49d7eSJed Brown     ///     v[i] = u[i] * w[i];
6839df49d7eSJed Brown     /// }
6849df49d7eSJed Brown     ///
685c68be7a2SJeremy L Thompson     /// let uu = ceed.vector_from_slice(&u)?;
686c68be7a2SJeremy L Thompson     /// let ww = ceed.vector_from_slice(&w)?;
687c68be7a2SJeremy L Thompson     /// let mut vv = ceed.vector(Q)?;
6889df49d7eSJed Brown     /// vv.set_value(0.0);
6899df49d7eSJed Brown     /// {
6909df49d7eSJed Brown     ///     let input = vec![uu, ww];
6919df49d7eSJed Brown     ///     let mut output = vec![vv];
692c68be7a2SJeremy L Thompson     ///     qf.apply(Q, &input, &output)?;
6939df49d7eSJed Brown     ///     vv = output.remove(0);
6949df49d7eSJed Brown     /// }
6959df49d7eSJed Brown     ///
6969df49d7eSJed Brown     /// vv.view()
6979df49d7eSJed Brown     ///     .iter()
6989df49d7eSJed Brown     ///     .zip(v.iter())
6999df49d7eSJed Brown     ///     .for_each(|(computed, actual)| {
7009df49d7eSJed Brown     ///         assert_eq!(
7019df49d7eSJed Brown     ///             *computed, *actual,
7029df49d7eSJed Brown     ///             "Incorrect value in QFunction application"
7039df49d7eSJed Brown     ///         );
7049df49d7eSJed Brown     ///     });
705c68be7a2SJeremy L Thompson     /// # Ok(())
706c68be7a2SJeremy L Thompson     /// # }
7079df49d7eSJed Brown     /// ```
7089df49d7eSJed Brown     pub fn apply(&self, Q: usize, u: &[Vector], v: &[Vector]) -> crate::Result<i32> {
7099df49d7eSJed Brown         self.qf_core.apply(Q, u, v)
7109df49d7eSJed Brown     }
7119df49d7eSJed Brown 
7129df49d7eSJed Brown     /// Add a QFunction input
7139df49d7eSJed Brown     ///
7149df49d7eSJed Brown     /// * `fieldname` - Name of QFunction field
7159df49d7eSJed Brown     /// * `size`      - Size of QFunction field, `(ncomp * dim)` for `Grad` or
7169df49d7eSJed Brown     ///                   `(ncomp * 1)` for `None`, `Interp`, and `Weight`
7179df49d7eSJed Brown     /// * `emode`     - `EvalMode::None` to use values directly, `EvalMode::Interp`
7189df49d7eSJed Brown     ///                   to use interpolated values, `EvalMode::Grad` to use
7199df49d7eSJed Brown     ///                   gradients, `EvalMode::Weight` to use quadrature weights
7209df49d7eSJed Brown     ///
7219df49d7eSJed Brown     /// ```
7229df49d7eSJed Brown     /// # use libceed::prelude::*;
723c68be7a2SJeremy L Thompson     /// # fn main() -> Result<(), libceed::CeedError> {
7249df49d7eSJed Brown     /// # let ceed = libceed::Ceed::default_init();
7259df49d7eSJed Brown     /// let mut user_f = |[u, weights, ..]: QFunctionInputs, [v, ..]: QFunctionOutputs| {
7269df49d7eSJed Brown     ///     // Iterate over quadrature points
7279df49d7eSJed Brown     ///     v.iter_mut()
7289df49d7eSJed Brown     ///         .zip(u.iter().zip(weights.iter()))
7299df49d7eSJed Brown     ///         .for_each(|(v, (u, w))| *v = u * w);
7309df49d7eSJed Brown     ///
7319df49d7eSJed Brown     ///     // Return clean error code
7329df49d7eSJed Brown     ///     0
7339df49d7eSJed Brown     /// };
7349df49d7eSJed Brown     ///
735c68be7a2SJeremy L Thompson     /// let mut qf = ceed.q_function_interior(1, Box::new(user_f))?;
7369df49d7eSJed Brown     ///
737c68be7a2SJeremy L Thompson     /// qf = qf.input("u", 1, EvalMode::Interp)?;
738c68be7a2SJeremy L Thompson     /// qf = qf.input("weights", 1, EvalMode::Weight)?;
739c68be7a2SJeremy L Thompson     /// # Ok(())
740c68be7a2SJeremy L Thompson     /// # }
7419df49d7eSJed Brown     /// ```
7429df49d7eSJed Brown     pub fn input(
7439df49d7eSJed Brown         mut self,
7449df49d7eSJed Brown         fieldname: &str,
7459df49d7eSJed Brown         size: usize,
7469df49d7eSJed Brown         emode: crate::EvalMode,
7479df49d7eSJed Brown     ) -> crate::Result<Self> {
7489df49d7eSJed Brown         let name_c = CString::new(fieldname).expect("CString::new failed");
7499df49d7eSJed Brown         let idx = self.trampoline_data.number_inputs;
7509df49d7eSJed Brown         self.trampoline_data.input_sizes[idx] = size;
7519df49d7eSJed Brown         self.trampoline_data.number_inputs += 1;
7529df49d7eSJed Brown         let (size, emode) = (
7539df49d7eSJed Brown             i32::try_from(size).unwrap(),
7549df49d7eSJed Brown             emode as bind_ceed::CeedEvalMode,
7559df49d7eSJed Brown         );
7569df49d7eSJed Brown         let ierr = unsafe {
7579df49d7eSJed Brown             bind_ceed::CeedQFunctionAddInput(self.qf_core.ptr, name_c.as_ptr(), size, emode)
7589df49d7eSJed Brown         };
7591142270cSJeremy L Thompson         self.qf_core.check_error(ierr)?;
7609df49d7eSJed Brown         Ok(self)
7619df49d7eSJed Brown     }
7629df49d7eSJed Brown 
7639df49d7eSJed Brown     /// Add a QFunction output
7649df49d7eSJed Brown     ///
7659df49d7eSJed Brown     /// * `fieldname` - Name of QFunction field
7669df49d7eSJed Brown     /// * `size`      - Size of QFunction field, `(ncomp * dim)` for `Grad` or
7679df49d7eSJed Brown     ///                   `(ncomp * 1)` for `None` and `Interp`
7689df49d7eSJed Brown     /// * `emode`     - `EvalMode::None` to use values directly, `EvalMode::Interp`
7699df49d7eSJed Brown     ///                   to use interpolated values, `EvalMode::Grad` to use
7709df49d7eSJed Brown     ///                   gradients
7719df49d7eSJed Brown     ///
7729df49d7eSJed Brown     /// ```
7739df49d7eSJed Brown     /// # use libceed::prelude::*;
774c68be7a2SJeremy L Thompson     /// # fn main() -> Result<(), libceed::CeedError> {
7759df49d7eSJed Brown     /// # let ceed = libceed::Ceed::default_init();
7769df49d7eSJed Brown     /// let mut user_f = |[u, weights, ..]: QFunctionInputs, [v, ..]: QFunctionOutputs| {
7779df49d7eSJed Brown     ///     // Iterate over quadrature points
7789df49d7eSJed Brown     ///     v.iter_mut()
7799df49d7eSJed Brown     ///         .zip(u.iter().zip(weights.iter()))
7809df49d7eSJed Brown     ///         .for_each(|(v, (u, w))| *v = u * w);
7819df49d7eSJed Brown     ///
7829df49d7eSJed Brown     ///     // Return clean error code
7839df49d7eSJed Brown     ///     0
7849df49d7eSJed Brown     /// };
7859df49d7eSJed Brown     ///
786c68be7a2SJeremy L Thompson     /// let mut qf = ceed.q_function_interior(1, Box::new(user_f))?;
7879df49d7eSJed Brown     ///
788c68be7a2SJeremy L Thompson     /// qf.output("v", 1, EvalMode::Interp)?;
789c68be7a2SJeremy L Thompson     /// # Ok(())
790c68be7a2SJeremy L Thompson     /// # }
7919df49d7eSJed Brown     /// ```
7929df49d7eSJed Brown     pub fn output(
7939df49d7eSJed Brown         mut self,
7949df49d7eSJed Brown         fieldname: &str,
7959df49d7eSJed Brown         size: usize,
7969df49d7eSJed Brown         emode: crate::EvalMode,
7979df49d7eSJed Brown     ) -> crate::Result<Self> {
7989df49d7eSJed Brown         let name_c = CString::new(fieldname).expect("CString::new failed");
7999df49d7eSJed Brown         let idx = self.trampoline_data.number_outputs;
8009df49d7eSJed Brown         self.trampoline_data.output_sizes[idx] = size;
8019df49d7eSJed Brown         self.trampoline_data.number_outputs += 1;
8029df49d7eSJed Brown         let (size, emode) = (
8039df49d7eSJed Brown             i32::try_from(size).unwrap(),
8049df49d7eSJed Brown             emode as bind_ceed::CeedEvalMode,
8059df49d7eSJed Brown         );
8069df49d7eSJed Brown         let ierr = unsafe {
8079df49d7eSJed Brown             bind_ceed::CeedQFunctionAddOutput(self.qf_core.ptr, name_c.as_ptr(), size, emode)
8089df49d7eSJed Brown         };
8091142270cSJeremy L Thompson         self.qf_core.check_error(ierr)?;
8109df49d7eSJed Brown         Ok(self)
8119df49d7eSJed Brown     }
812*08778c6fSJeremy L Thompson 
813*08778c6fSJeremy L Thompson     /// Get a slice of QFunction inputs
814*08778c6fSJeremy L Thompson     ///
815*08778c6fSJeremy L Thompson     /// ```
816*08778c6fSJeremy L Thompson     /// # use libceed::prelude::*;
817*08778c6fSJeremy L Thompson     /// # fn main() -> Result<(), libceed::CeedError> {
818*08778c6fSJeremy L Thompson     /// # let ceed = libceed::Ceed::default_init();
819*08778c6fSJeremy L Thompson     /// let mut user_f = |[u, weights, ..]: QFunctionInputs, [v, ..]: QFunctionOutputs| {
820*08778c6fSJeremy L Thompson     ///     // Iterate over quadrature points
821*08778c6fSJeremy L Thompson     ///     v.iter_mut()
822*08778c6fSJeremy L Thompson     ///         .zip(u.iter().zip(weights.iter()))
823*08778c6fSJeremy L Thompson     ///         .for_each(|(v, (u, w))| *v = u * w);
824*08778c6fSJeremy L Thompson     ///
825*08778c6fSJeremy L Thompson     ///     // Return clean error code
826*08778c6fSJeremy L Thompson     ///     0
827*08778c6fSJeremy L Thompson     /// };
828*08778c6fSJeremy L Thompson     ///
829*08778c6fSJeremy L Thompson     /// let mut qf = ceed
830*08778c6fSJeremy L Thompson     ///     .q_function_interior(1, Box::new(user_f))?
831*08778c6fSJeremy L Thompson     ///     .input("u", 1, EvalMode::Interp)?
832*08778c6fSJeremy L Thompson     ///     .input("weights", 1, EvalMode::Weight)?;
833*08778c6fSJeremy L Thompson     ///
834*08778c6fSJeremy L Thompson     /// let inputs = qf.inputs()?;
835*08778c6fSJeremy L Thompson     ///
836*08778c6fSJeremy L Thompson     /// assert_eq!(inputs.len(), 2, "Incorrect inputs array");
837*08778c6fSJeremy L Thompson     /// # Ok(())
838*08778c6fSJeremy L Thompson     /// # }
839*08778c6fSJeremy L Thompson     /// ```
840*08778c6fSJeremy L Thompson     pub fn inputs(&self) -> crate::Result<&[crate::QFunctionField]> {
841*08778c6fSJeremy L Thompson         self.qf_core.inputs()
842*08778c6fSJeremy L Thompson     }
843*08778c6fSJeremy L Thompson 
844*08778c6fSJeremy L Thompson     /// Get a slice of QFunction outputs
845*08778c6fSJeremy L Thompson     ///
846*08778c6fSJeremy L Thompson     /// ```
847*08778c6fSJeremy L Thompson     /// # use libceed::prelude::*;
848*08778c6fSJeremy L Thompson     /// # fn main() -> Result<(), libceed::CeedError> {
849*08778c6fSJeremy L Thompson     /// # let ceed = libceed::Ceed::default_init();
850*08778c6fSJeremy L Thompson     /// let mut user_f = |[u, weights, ..]: QFunctionInputs, [v, ..]: QFunctionOutputs| {
851*08778c6fSJeremy L Thompson     ///     // Iterate over quadrature points
852*08778c6fSJeremy L Thompson     ///     v.iter_mut()
853*08778c6fSJeremy L Thompson     ///         .zip(u.iter().zip(weights.iter()))
854*08778c6fSJeremy L Thompson     ///         .for_each(|(v, (u, w))| *v = u * w);
855*08778c6fSJeremy L Thompson     ///
856*08778c6fSJeremy L Thompson     ///     // Return clean error code
857*08778c6fSJeremy L Thompson     ///     0
858*08778c6fSJeremy L Thompson     /// };
859*08778c6fSJeremy L Thompson     ///
860*08778c6fSJeremy L Thompson     /// let mut qf = ceed
861*08778c6fSJeremy L Thompson     ///     .q_function_interior(1, Box::new(user_f))?
862*08778c6fSJeremy L Thompson     ///     .output("v", 1, EvalMode::Interp)?;
863*08778c6fSJeremy L Thompson     ///
864*08778c6fSJeremy L Thompson     /// let outputs = qf.outputs()?;
865*08778c6fSJeremy L Thompson     ///
866*08778c6fSJeremy L Thompson     /// assert_eq!(outputs.len(), 1, "Incorrect outputs array");
867*08778c6fSJeremy L Thompson     /// # Ok(())
868*08778c6fSJeremy L Thompson     /// # }
869*08778c6fSJeremy L Thompson     /// ```
870*08778c6fSJeremy L Thompson     pub fn outputs(&self) -> crate::Result<&[crate::QFunctionField]> {
871*08778c6fSJeremy L Thompson         self.qf_core.outputs()
872*08778c6fSJeremy L Thompson     }
8739df49d7eSJed Brown }
8749df49d7eSJed Brown 
8759df49d7eSJed Brown // -----------------------------------------------------------------------------
8769df49d7eSJed Brown // QFunction
8779df49d7eSJed Brown // -----------------------------------------------------------------------------
8789df49d7eSJed Brown impl<'a> QFunctionByName<'a> {
8799df49d7eSJed Brown     // Constructor
8809df49d7eSJed Brown     pub fn create(ceed: &'a crate::Ceed, name: &str) -> crate::Result<Self> {
8819df49d7eSJed Brown         let name_c = CString::new(name).expect("CString::new failed");
8829df49d7eSJed Brown         let mut ptr = std::ptr::null_mut();
8839df49d7eSJed Brown         let ierr = unsafe {
8849df49d7eSJed Brown             bind_ceed::CeedQFunctionCreateInteriorByName(ceed.ptr, name_c.as_ptr(), &mut ptr)
8859df49d7eSJed Brown         };
8869df49d7eSJed Brown         ceed.check_error(ierr)?;
8879df49d7eSJed Brown         Ok(Self {
8881142270cSJeremy L Thompson             qf_core: QFunctionCore {
8891142270cSJeremy L Thompson                 ptr,
8901142270cSJeremy L Thompson                 _lifeline: PhantomData,
8911142270cSJeremy L Thompson             },
8929df49d7eSJed Brown         })
8939df49d7eSJed Brown     }
8949df49d7eSJed Brown 
8959df49d7eSJed Brown     /// Apply the action of a QFunction
8969df49d7eSJed Brown     ///
8979df49d7eSJed Brown     /// * `Q`      - The number of quadrature points
8989df49d7eSJed Brown     /// * `input`  - Array of input Vectors
8999df49d7eSJed Brown     /// * `output` - Array of output Vectors
9009df49d7eSJed Brown     ///
9019df49d7eSJed Brown     /// ```
9029df49d7eSJed Brown     /// # use libceed::prelude::*;
903c68be7a2SJeremy L Thompson     /// # fn main() -> Result<(), libceed::CeedError> {
9049df49d7eSJed Brown     /// # let ceed = libceed::Ceed::default_init();
9059df49d7eSJed Brown     /// const Q: usize = 8;
906c68be7a2SJeremy L Thompson     /// let qf_build = ceed.q_function_interior_by_name("Mass1DBuild")?;
907c68be7a2SJeremy L Thompson     /// let qf_mass = ceed.q_function_interior_by_name("MassApply")?;
9089df49d7eSJed Brown     ///
9099df49d7eSJed Brown     /// let mut j = [0.; Q];
9109df49d7eSJed Brown     /// let mut w = [0.; Q];
9119df49d7eSJed Brown     /// let mut u = [0.; Q];
9129df49d7eSJed Brown     /// let mut v = [0.; Q];
9139df49d7eSJed Brown     ///
9149df49d7eSJed Brown     /// for i in 0..Q {
91580a9ef05SNatalie Beams     ///     let x = 2. * (i as Scalar) / ((Q as Scalar) - 1.) - 1.;
9169df49d7eSJed Brown     ///     j[i] = 1.;
9179df49d7eSJed Brown     ///     w[i] = 1. - x * x;
9189df49d7eSJed Brown     ///     u[i] = 2. + 3. * x + 5. * x * x;
9199df49d7eSJed Brown     ///     v[i] = w[i] * u[i];
9209df49d7eSJed Brown     /// }
9219df49d7eSJed Brown     ///
922c68be7a2SJeremy L Thompson     /// let jj = ceed.vector_from_slice(&j)?;
923c68be7a2SJeremy L Thompson     /// let ww = ceed.vector_from_slice(&w)?;
924c68be7a2SJeremy L Thompson     /// let uu = ceed.vector_from_slice(&u)?;
925c68be7a2SJeremy L Thompson     /// let mut vv = ceed.vector(Q)?;
9269df49d7eSJed Brown     /// vv.set_value(0.0);
927c68be7a2SJeremy L Thompson     /// let mut qdata = ceed.vector(Q)?;
9289df49d7eSJed Brown     /// qdata.set_value(0.0);
9299df49d7eSJed Brown     ///
9309df49d7eSJed Brown     /// {
9319df49d7eSJed Brown     ///     let mut input = vec![jj, ww];
9329df49d7eSJed Brown     ///     let mut output = vec![qdata];
933c68be7a2SJeremy L Thompson     ///     qf_build.apply(Q, &input, &output)?;
9349df49d7eSJed Brown     ///     qdata = output.remove(0);
9359df49d7eSJed Brown     /// }
9369df49d7eSJed Brown     ///
9379df49d7eSJed Brown     /// {
9389df49d7eSJed Brown     ///     let mut input = vec![qdata, uu];
9399df49d7eSJed Brown     ///     let mut output = vec![vv];
940c68be7a2SJeremy L Thompson     ///     qf_mass.apply(Q, &input, &output)?;
9419df49d7eSJed Brown     ///     vv = output.remove(0);
9429df49d7eSJed Brown     /// }
9439df49d7eSJed Brown     ///
9449df49d7eSJed Brown     /// vv.view()
9459df49d7eSJed Brown     ///     .iter()
9469df49d7eSJed Brown     ///     .zip(v.iter())
9479df49d7eSJed Brown     ///     .for_each(|(computed, actual)| {
9489df49d7eSJed Brown     ///         assert_eq!(
9499df49d7eSJed Brown     ///             *computed, *actual,
9509df49d7eSJed Brown     ///             "Incorrect value in QFunction application"
9519df49d7eSJed Brown     ///         );
9529df49d7eSJed Brown     ///     });
953c68be7a2SJeremy L Thompson     /// # Ok(())
954c68be7a2SJeremy L Thompson     /// # }
9559df49d7eSJed Brown     /// ```
9569df49d7eSJed Brown     pub fn apply(&self, Q: usize, u: &[Vector], v: &[Vector]) -> crate::Result<i32> {
9579df49d7eSJed Brown         self.qf_core.apply(Q, u, v)
9589df49d7eSJed Brown     }
959*08778c6fSJeremy L Thompson 
960*08778c6fSJeremy L Thompson     /// Get a slice of QFunction inputs
961*08778c6fSJeremy L Thompson     ///
962*08778c6fSJeremy L Thompson     /// ```
963*08778c6fSJeremy L Thompson     /// # use libceed::prelude::*;
964*08778c6fSJeremy L Thompson     /// # fn main() -> Result<(), libceed::CeedError> {
965*08778c6fSJeremy L Thompson     /// # let ceed = libceed::Ceed::default_init();
966*08778c6fSJeremy L Thompson     /// const Q: usize = 8;
967*08778c6fSJeremy L Thompson     /// let qf = ceed.q_function_interior_by_name("Mass1DBuild")?;
968*08778c6fSJeremy L Thompson     ///
969*08778c6fSJeremy L Thompson     /// let inputs = qf.inputs()?;
970*08778c6fSJeremy L Thompson     ///
971*08778c6fSJeremy L Thompson     /// assert_eq!(inputs.len(), 2, "Incorrect inputs array");
972*08778c6fSJeremy L Thompson     /// # Ok(())
973*08778c6fSJeremy L Thompson     /// # }
974*08778c6fSJeremy L Thompson     /// ```
975*08778c6fSJeremy L Thompson     pub fn inputs(&self) -> crate::Result<&[crate::QFunctionField]> {
976*08778c6fSJeremy L Thompson         self.qf_core.inputs()
977*08778c6fSJeremy L Thompson     }
978*08778c6fSJeremy L Thompson 
979*08778c6fSJeremy L Thompson     /// Get a slice of QFunction outputs
980*08778c6fSJeremy L Thompson     ///
981*08778c6fSJeremy L Thompson     /// ```
982*08778c6fSJeremy L Thompson     /// # use libceed::prelude::*;
983*08778c6fSJeremy L Thompson     /// # fn main() -> Result<(), libceed::CeedError> {
984*08778c6fSJeremy L Thompson     /// # let ceed = libceed::Ceed::default_init();
985*08778c6fSJeremy L Thompson     /// const Q: usize = 8;
986*08778c6fSJeremy L Thompson     /// let qf = ceed.q_function_interior_by_name("Mass1DBuild")?;
987*08778c6fSJeremy L Thompson     ///
988*08778c6fSJeremy L Thompson     /// let outputs = qf.outputs()?;
989*08778c6fSJeremy L Thompson     ///
990*08778c6fSJeremy L Thompson     /// assert_eq!(outputs.len(), 1, "Incorrect outputs array");
991*08778c6fSJeremy L Thompson     /// # Ok(())
992*08778c6fSJeremy L Thompson     /// # }
993*08778c6fSJeremy L Thompson     /// ```
994*08778c6fSJeremy L Thompson     pub fn outputs(&self) -> crate::Result<&[crate::QFunctionField]> {
995*08778c6fSJeremy L Thompson         self.qf_core.outputs()
996*08778c6fSJeremy L Thompson     }
9979df49d7eSJed Brown }
9989df49d7eSJed Brown 
9999df49d7eSJed Brown // -----------------------------------------------------------------------------
1000