1*9ba83ac0SJeremy L Thompson // Copyright (c) 2017-2026, Lawrence Livermore National Security, LLC and other CEED contributors.
23d8e8822SJeremy L Thompson // All Rights Reserved. See the top-level LICENSE and NOTICE files for details.
3ded9b81dSJeremy L Thompson //
43d8e8822SJeremy L Thompson // SPDX-License-Identifier: BSD-2-Clause
5ded9b81dSJeremy L Thompson //
63d8e8822SJeremy L Thompson // This file is part of CEED: http://github.com/ceed
73d8e8822SJeremy L Thompson //
8ded9b81dSJeremy L Thompson // libCEED Example 1
9ded9b81dSJeremy L Thompson //
10ded9b81dSJeremy L Thompson // This example illustrates a simple usage of libCEED to compute the volume of a
11ded9b81dSJeremy L Thompson // 3D body using matrix-free application of a mass operator. Arbitrary mesh and
12ded9b81dSJeremy L Thompson // solution orders in 1D, 2D and 3D are supported from the same code.
13ded9b81dSJeremy L Thompson //
14ded9b81dSJeremy L Thompson // The example has no dependencies, and is designed to be self-contained. For
15ded9b81dSJeremy L Thompson // additional examples that use external discretization libraries (MFEM, PETSc,
16ded9b81dSJeremy L Thompson // etc.) see the subdirectories in libceed/examples.
17ded9b81dSJeremy L Thompson //
18ded9b81dSJeremy L Thompson // All libCEED objects use a Ceed device object constructed based on a command
19ded9b81dSJeremy L Thompson // line argument (-ceed).
20ded9b81dSJeremy L Thompson
21947f93aaSJed Brown use clap::Parser;
22eb07d68fSJeremy L Thompson use libceed::{
23eb07d68fSJeremy L Thompson BasisOpt, Ceed, ElemRestrictionOpt, QFunctionInputs, QFunctionOpt, QFunctionOutputs, VectorOpt,
24eb07d68fSJeremy L Thompson };
25ded9b81dSJeremy L Thompson mod opt;
26ded9b81dSJeremy L Thompson mod transform;
27ded9b81dSJeremy L Thompson
28ded9b81dSJeremy L Thompson // ----------------------------------------------------------------------------
29ded9b81dSJeremy L Thompson // Example 1
30ded9b81dSJeremy L Thompson // ----------------------------------------------------------------------------
main() -> libceed::Result<()>3189d15d5fSJeremy L Thompson fn main() -> libceed::Result<()> {
32947f93aaSJed Brown let options = opt::Opt::parse();
33ded9b81dSJeremy L Thompson example_1(options)
34ded9b81dSJeremy L Thompson }
35ded9b81dSJeremy L Thompson
3678c2cefaSJeremy L Thompson #[allow(clippy::erasing_op)]
3778c2cefaSJeremy L Thompson #[allow(clippy::identity_op)]
example_1(options: opt::Opt) -> libceed::Result<()>3889d15d5fSJeremy L Thompson fn example_1(options: opt::Opt) -> libceed::Result<()> {
39ded9b81dSJeremy L Thompson // Process command line arguments
40ded9b81dSJeremy L Thompson let opt::Opt {
41ded9b81dSJeremy L Thompson ceed_spec,
42ded9b81dSJeremy L Thompson dim,
43ded9b81dSJeremy L Thompson mesh_degree,
44ded9b81dSJeremy L Thompson solution_degree,
45ded9b81dSJeremy L Thompson num_qpts,
46ded9b81dSJeremy L Thompson problem_size_requested,
47ded9b81dSJeremy L Thompson test,
48ded9b81dSJeremy L Thompson quiet,
49ded9b81dSJeremy L Thompson gallery,
50ded9b81dSJeremy L Thompson } = options;
51bf55b007SJeremy L Thompson assert!((1..=3).contains(&dim));
52ded9b81dSJeremy L Thompson assert!(mesh_degree >= 1);
53ded9b81dSJeremy L Thompson assert!(solution_degree >= 1);
54ded9b81dSJeremy L Thompson assert!(num_qpts >= 1);
55ded9b81dSJeremy L Thompson let ncomp_x = dim;
56bf55b007SJeremy L Thompson let problem_size: i64 = if problem_size_requested < 0 {
57bf55b007SJeremy L Thompson if test {
58bf55b007SJeremy L Thompson 8 * 16
59ded9b81dSJeremy L Thompson } else {
60bf55b007SJeremy L Thompson 256 * 1024
61ded9b81dSJeremy L Thompson }
62bf55b007SJeremy L Thompson } else {
63bf55b007SJeremy L Thompson problem_size_requested
64bf55b007SJeremy L Thompson };
65ded9b81dSJeremy L Thompson
66ded9b81dSJeremy L Thompson // Summary output
67ded9b81dSJeremy L Thompson if !quiet {
68ded9b81dSJeremy L Thompson println!("Selected options: [command line option] : <current value>");
69ded9b81dSJeremy L Thompson println!(" Ceed specification [-c] : {}", ceed_spec);
70ded9b81dSJeremy L Thompson println!(" Mesh dimension [-d] : {}", dim);
71ded9b81dSJeremy L Thompson println!(" Mesh degree [-m] : {}", mesh_degree);
72ded9b81dSJeremy L Thompson println!(" Solution degree [-p] : {}", solution_degree);
73ded9b81dSJeremy L Thompson println!(" Num. 1D quadr. pts [-q] : {}", num_qpts);
74ded9b81dSJeremy L Thompson println!(" Approx. # unknowns [-s] : {}", problem_size);
75ded9b81dSJeremy L Thompson println!(
76ded9b81dSJeremy L Thompson " QFunction source [-g] : {}",
77ded9b81dSJeremy L Thompson if gallery { "gallery" } else { "user closure" }
78ded9b81dSJeremy L Thompson );
79ded9b81dSJeremy L Thompson }
80ded9b81dSJeremy L Thompson
81ded9b81dSJeremy L Thompson // Initalize ceed context
82ded9b81dSJeremy L Thompson let ceed = Ceed::init(&ceed_spec);
83ded9b81dSJeremy L Thompson
84ded9b81dSJeremy L Thompson // Mesh and solution bases
85eb07d68fSJeremy L Thompson let basis_mesh = ceed.basis_tensor_H1_Lagrange(
86eb07d68fSJeremy L Thompson dim,
87eb07d68fSJeremy L Thompson ncomp_x,
88eb07d68fSJeremy L Thompson mesh_degree + 1,
89eb07d68fSJeremy L Thompson num_qpts,
90eb07d68fSJeremy L Thompson libceed::QuadMode::Gauss,
91eb07d68fSJeremy L Thompson )?;
92eb07d68fSJeremy L Thompson let basis_solution = ceed.basis_tensor_H1_Lagrange(
93eb07d68fSJeremy L Thompson dim,
94eb07d68fSJeremy L Thompson 1,
95eb07d68fSJeremy L Thompson solution_degree + 1,
96eb07d68fSJeremy L Thompson num_qpts,
97eb07d68fSJeremy L Thompson libceed::QuadMode::Gauss,
98eb07d68fSJeremy L Thompson )?;
99ded9b81dSJeremy L Thompson
100ded9b81dSJeremy L Thompson // Determine mesh size from approximate problem size
101ded9b81dSJeremy L Thompson let num_xyz = mesh::cartesian_mesh_size(dim, solution_degree, problem_size);
102ded9b81dSJeremy L Thompson if !quiet {
103ded9b81dSJeremy L Thompson print!("\nMesh size : nx = {}", num_xyz[0]);
104ded9b81dSJeremy L Thompson if dim > 1 {
105ded9b81dSJeremy L Thompson print!(", ny = {}", num_xyz[1]);
106ded9b81dSJeremy L Thompson }
107ded9b81dSJeremy L Thompson if dim > 2 {
108ded9b81dSJeremy L Thompson print!(", nz = {}", num_xyz[2]);
109ded9b81dSJeremy L Thompson }
110bf55b007SJeremy L Thompson println!();
111ded9b81dSJeremy L Thompson }
112ded9b81dSJeremy L Thompson
113ded9b81dSJeremy L Thompson // Build ElemRestriction objects describing the mesh and solution discrete
114ded9b81dSJeremy L Thompson // representations
115edb2538eSJeremy L Thompson let (rstr_mesh, _) =
116e15f9bd0SJeremy L Thompson mesh::build_cartesian_restriction(&ceed, dim, num_xyz, mesh_degree, ncomp_x, num_qpts)?;
117edb2538eSJeremy L Thompson let (rstr_solution, rstr_qdata) =
118e15f9bd0SJeremy L Thompson mesh::build_cartesian_restriction(&ceed, dim, num_xyz, solution_degree, 1, num_qpts)?;
119edb2538eSJeremy L Thompson let mesh_size = rstr_mesh.lvector_size();
120edb2538eSJeremy L Thompson let solution_size = rstr_solution.lvector_size();
121ded9b81dSJeremy L Thompson if !quiet {
122ded9b81dSJeremy L Thompson println!("Number of mesh nodes : {}", mesh_size / dim);
123ded9b81dSJeremy L Thompson println!("Number of solution nodes : {}", solution_size);
124ded9b81dSJeremy L Thompson }
125ded9b81dSJeremy L Thompson
126ded9b81dSJeremy L Thompson // Create a Vector with the mesh coordinates
127e15f9bd0SJeremy L Thompson let mut mesh_coords = mesh::cartesian_mesh_coords(&ceed, dim, num_xyz, mesh_degree, mesh_size)?;
128ded9b81dSJeremy L Thompson
129ded9b81dSJeremy L Thompson // Apply a transformation to the mesh coordinates
130e78171edSJeremy L Thompson let exact_volume = transform::transform_mesh_coordinates(dim, mesh_size, &mut mesh_coords)?;
131ded9b81dSJeremy L Thompson
132ded9b81dSJeremy L Thompson // QFunction that builds the quadrature data for the mass operator
133ded9b81dSJeremy L Thompson // -- QFunction from user closure
134ded9b81dSJeremy L Thompson let build_mass = move |[jacobian, weights, ..]: QFunctionInputs,
135ded9b81dSJeremy L Thompson [qdata, ..]: QFunctionOutputs| {
136ded9b81dSJeremy L Thompson // Build quadrature data
137ded9b81dSJeremy L Thompson match dim {
138ded9b81dSJeremy L Thompson 1 => qdata
139ded9b81dSJeremy L Thompson .iter_mut()
140ded9b81dSJeremy L Thompson .zip(jacobian.iter().zip(weights.iter()))
141ded9b81dSJeremy L Thompson .for_each(|(qdata, (j, weight))| *qdata = j * weight),
142ded9b81dSJeremy L Thompson 2 => {
143ded9b81dSJeremy L Thompson let q = qdata.len();
144ded9b81dSJeremy L Thompson qdata.iter_mut().zip(weights.iter()).enumerate().for_each(
145ded9b81dSJeremy L Thompson |(i, (qdata, weight))| {
146ded9b81dSJeremy L Thompson *qdata = (jacobian[i + q * 0] * jacobian[i + q * 3]
147ded9b81dSJeremy L Thompson - jacobian[i + q * 1] * jacobian[i + q * 2])
148ded9b81dSJeremy L Thompson * weight
149ded9b81dSJeremy L Thompson },
150ded9b81dSJeremy L Thompson );
151ded9b81dSJeremy L Thompson }
152ded9b81dSJeremy L Thompson 3 => {
153ded9b81dSJeremy L Thompson let q = qdata.len();
154ded9b81dSJeremy L Thompson qdata.iter_mut().zip(weights.iter()).enumerate().for_each(
155ded9b81dSJeremy L Thompson |(i, (qdata, weight))| {
156ded9b81dSJeremy L Thompson *qdata = (jacobian[i + q * 0]
157ded9b81dSJeremy L Thompson * (jacobian[i + q * 4] * jacobian[i + q * 8]
158ded9b81dSJeremy L Thompson - jacobian[i + q * 5] * jacobian[i + q * 7])
159ded9b81dSJeremy L Thompson - jacobian[i + q * 1]
160ded9b81dSJeremy L Thompson * (jacobian[i + q * 3] * jacobian[i + q * 8]
161ded9b81dSJeremy L Thompson - jacobian[i + q * 5] * jacobian[i + q * 6])
162ded9b81dSJeremy L Thompson + jacobian[i + q * 2]
163ded9b81dSJeremy L Thompson * (jacobian[i + q * 3] * jacobian[i + q * 7]
164ded9b81dSJeremy L Thompson - jacobian[i + q * 4] * jacobian[i + q * 6]))
165ded9b81dSJeremy L Thompson * weight
166ded9b81dSJeremy L Thompson },
167ded9b81dSJeremy L Thompson );
168ded9b81dSJeremy L Thompson }
169ded9b81dSJeremy L Thompson _ => unreachable!(),
170ded9b81dSJeremy L Thompson };
171ded9b81dSJeremy L Thompson
172ded9b81dSJeremy L Thompson // Return clean error code
173ded9b81dSJeremy L Thompson 0
174ded9b81dSJeremy L Thompson };
175ded9b81dSJeremy L Thompson let qf_build_closure = ceed
176e15f9bd0SJeremy L Thompson .q_function_interior(1, Box::new(build_mass))?
177eb07d68fSJeremy L Thompson .input("dx", ncomp_x * dim, libceed::EvalMode::Grad)?
178eb07d68fSJeremy L Thompson .input("weights", 1, libceed::EvalMode::Weight)?
179eb07d68fSJeremy L Thompson .output("qdata", 1, libceed::EvalMode::None)?;
180ded9b81dSJeremy L Thompson // -- QFunction from gallery
181ded9b81dSJeremy L Thompson let qf_build_named = {
182ded9b81dSJeremy L Thompson let name = format!("Mass{}DBuild", dim);
183e15f9bd0SJeremy L Thompson ceed.q_function_interior_by_name(&name)?
184ded9b81dSJeremy L Thompson };
185ded9b81dSJeremy L Thompson // -- QFunction for use with Operator
186ded9b81dSJeremy L Thompson let qf_build = if gallery {
187ded9b81dSJeremy L Thompson QFunctionOpt::SomeQFunctionByName(&qf_build_named)
188ded9b81dSJeremy L Thompson } else {
189ded9b81dSJeremy L Thompson QFunctionOpt::SomeQFunction(&qf_build_closure)
190ded9b81dSJeremy L Thompson };
191ded9b81dSJeremy L Thompson
192ded9b81dSJeremy L Thompson // Operator that build the quadrature data for the mass operator
193ded9b81dSJeremy L Thompson let op_build = ceed
194e15f9bd0SJeremy L Thompson .operator(qf_build, QFunctionOpt::None, QFunctionOpt::None)?
195ea6b5821SJeremy L Thompson .name("build qdata")?
196edb2538eSJeremy L Thompson .field("dx", &rstr_mesh, &basis_mesh, VectorOpt::Active)?
197ded9b81dSJeremy L Thompson .field(
198ded9b81dSJeremy L Thompson "weights",
199ded9b81dSJeremy L Thompson ElemRestrictionOpt::None,
200ded9b81dSJeremy L Thompson &basis_mesh,
201ded9b81dSJeremy L Thompson VectorOpt::None,
202e15f9bd0SJeremy L Thompson )?
203edb2538eSJeremy L Thompson .field("qdata", &rstr_qdata, BasisOpt::None, VectorOpt::Active)?
2046f97ff0aSJeremy L Thompson .check()?;
205ded9b81dSJeremy L Thompson
206ded9b81dSJeremy L Thompson // Compute the quadrature data for the mass operator
207ded9b81dSJeremy L Thompson let elem_qpts = num_qpts.pow(dim as u32);
208ded9b81dSJeremy L Thompson let num_elem: usize = num_xyz.iter().take(dim).product();
209e15f9bd0SJeremy L Thompson let mut qdata = ceed.vector(num_elem * elem_qpts)?;
210e15f9bd0SJeremy L Thompson op_build.apply(&mesh_coords, &mut qdata)?;
211ded9b81dSJeremy L Thompson
212ded9b81dSJeremy L Thompson // QFunction that applies the mass operator
213ded9b81dSJeremy L Thompson // -- QFunction from user closure
214ded9b81dSJeremy L Thompson let apply_mass = |[u, qdata, ..]: QFunctionInputs, [v, ..]: QFunctionOutputs| {
215ded9b81dSJeremy L Thompson // Apply mass operator
216ded9b81dSJeremy L Thompson v.iter_mut()
217ded9b81dSJeremy L Thompson .zip(u.iter().zip(qdata.iter()))
218ded9b81dSJeremy L Thompson .for_each(|(v, (u, w))| *v = u * w);
219ded9b81dSJeremy L Thompson // Return clean error code
220ded9b81dSJeremy L Thompson 0
221ded9b81dSJeremy L Thompson };
222ded9b81dSJeremy L Thompson let qf_mass_closure = ceed
223e15f9bd0SJeremy L Thompson .q_function_interior(1, Box::new(apply_mass))?
224eb07d68fSJeremy L Thompson .input("u", 1, libceed::EvalMode::Interp)?
225eb07d68fSJeremy L Thompson .input("qdata", 1, libceed::EvalMode::None)?
226eb07d68fSJeremy L Thompson .output("v", 1, libceed::EvalMode::Interp)?;
227ded9b81dSJeremy L Thompson // -- QFunction from gallery
228e15f9bd0SJeremy L Thompson let qf_mass_named = ceed.q_function_interior_by_name("MassApply")?;
229ded9b81dSJeremy L Thompson // -- QFunction for use with Operator
230ded9b81dSJeremy L Thompson let qf_mass = if gallery {
231ded9b81dSJeremy L Thompson QFunctionOpt::SomeQFunctionByName(&qf_mass_named)
232ded9b81dSJeremy L Thompson } else {
233ded9b81dSJeremy L Thompson QFunctionOpt::SomeQFunction(&qf_mass_closure)
234ded9b81dSJeremy L Thompson };
235ded9b81dSJeremy L Thompson
236ded9b81dSJeremy L Thompson // Mass Operator
237ded9b81dSJeremy L Thompson let op_mass = ceed
238e15f9bd0SJeremy L Thompson .operator(qf_mass, QFunctionOpt::None, QFunctionOpt::None)?
239ea6b5821SJeremy L Thompson .name("mass")?
240edb2538eSJeremy L Thompson .field("u", &rstr_solution, &basis_solution, VectorOpt::Active)?
241edb2538eSJeremy L Thompson .field("qdata", &rstr_qdata, BasisOpt::None, &qdata)?
242edb2538eSJeremy L Thompson .field("v", &rstr_solution, &basis_solution, VectorOpt::Active)?
2436f97ff0aSJeremy L Thompson .check()?;
244ded9b81dSJeremy L Thompson
245ded9b81dSJeremy L Thompson // Solution vectors
246e15f9bd0SJeremy L Thompson let u = ceed.vector_from_slice(&vec![1.0; solution_size])?;
247e15f9bd0SJeremy L Thompson let mut v = ceed.vector(solution_size)?;
248ded9b81dSJeremy L Thompson
249ded9b81dSJeremy L Thompson // Apply the mass operator
250e15f9bd0SJeremy L Thompson op_mass.apply(&u, &mut v)?;
251ded9b81dSJeremy L Thompson
252ded9b81dSJeremy L Thompson // Compute the mesh volume
253eb07d68fSJeremy L Thompson let volume: libceed::Scalar = v.view()?.iter().sum();
254ded9b81dSJeremy L Thompson
255ded9b81dSJeremy L Thompson // Output results
256ded9b81dSJeremy L Thompson if !quiet {
257ded9b81dSJeremy L Thompson println!("Exact mesh volume : {:.12}", exact_volume);
258ded9b81dSJeremy L Thompson println!("Computed mesh volume : {:.12}", volume);
259ded9b81dSJeremy L Thompson println!(
260ded9b81dSJeremy L Thompson "Volume error : {:.12e}",
261ded9b81dSJeremy L Thompson volume - exact_volume
262ded9b81dSJeremy L Thompson );
263ded9b81dSJeremy L Thompson }
264ded9b81dSJeremy L Thompson let tolerance = match dim {
2652247a93fSRezgar Shakeri 1 => 500.0 * libceed::EPSILON,
266ded9b81dSJeremy L Thompson _ => 1E-5,
267ded9b81dSJeremy L Thompson };
268ded9b81dSJeremy L Thompson let error = (volume - exact_volume).abs();
269ded9b81dSJeremy L Thompson if error > tolerance {
270ded9b81dSJeremy L Thompson println!("Volume error too large: {:.12e}", error);
2712ba8e59cSJeremy L Thompson return Err(libceed::Error {
272e15f9bd0SJeremy L Thompson message: format!(
273ded9b81dSJeremy L Thompson "Volume error too large - expected: {:.12e}, actual: {:.12e}",
274ded9b81dSJeremy L Thompson tolerance, error
275e15f9bd0SJeremy L Thompson ),
276e15f9bd0SJeremy L Thompson });
277ded9b81dSJeremy L Thompson }
278ded9b81dSJeremy L Thompson Ok(())
279ded9b81dSJeremy L Thompson }
280ded9b81dSJeremy L Thompson
281ded9b81dSJeremy L Thompson // ----------------------------------------------------------------------------
282ded9b81dSJeremy L Thompson // Tests
283ded9b81dSJeremy L Thompson // ----------------------------------------------------------------------------
284ded9b81dSJeremy L Thompson #[cfg(test)]
285ded9b81dSJeremy L Thompson mod tests {
286ded9b81dSJeremy L Thompson use super::*;
287ded9b81dSJeremy L Thompson
288ded9b81dSJeremy L Thompson #[test]
example_1_1d()289ded9b81dSJeremy L Thompson fn example_1_1d() {
290ded9b81dSJeremy L Thompson let options = opt::Opt {
291ded9b81dSJeremy L Thompson ceed_spec: "/cpu/self/ref/serial".to_string(),
292ded9b81dSJeremy L Thompson dim: 1,
293ded9b81dSJeremy L Thompson mesh_degree: 4,
294ded9b81dSJeremy L Thompson solution_degree: 4,
295ded9b81dSJeremy L Thompson num_qpts: 6,
296ded9b81dSJeremy L Thompson problem_size_requested: -1,
297ded9b81dSJeremy L Thompson test: true,
298ded9b81dSJeremy L Thompson quiet: true,
299ded9b81dSJeremy L Thompson gallery: false,
300ded9b81dSJeremy L Thompson };
301ded9b81dSJeremy L Thompson assert!(example_1(options).is_ok());
302ded9b81dSJeremy L Thompson }
303ded9b81dSJeremy L Thompson
304ded9b81dSJeremy L Thompson #[test]
example_1_2d()305ded9b81dSJeremy L Thompson fn example_1_2d() {
306ded9b81dSJeremy L Thompson let options = opt::Opt {
307ded9b81dSJeremy L Thompson ceed_spec: "/cpu/self/ref/serial".to_string(),
308ded9b81dSJeremy L Thompson dim: 2,
309ded9b81dSJeremy L Thompson mesh_degree: 4,
310ded9b81dSJeremy L Thompson solution_degree: 4,
311ded9b81dSJeremy L Thompson num_qpts: 6,
312ded9b81dSJeremy L Thompson problem_size_requested: -1,
313ded9b81dSJeremy L Thompson test: true,
314ded9b81dSJeremy L Thompson quiet: true,
315ded9b81dSJeremy L Thompson gallery: false,
316ded9b81dSJeremy L Thompson };
317ded9b81dSJeremy L Thompson assert!(example_1(options).is_ok());
318ded9b81dSJeremy L Thompson }
319ded9b81dSJeremy L Thompson
320ded9b81dSJeremy L Thompson #[test]
example_1_3d()321ded9b81dSJeremy L Thompson fn example_1_3d() {
322ded9b81dSJeremy L Thompson let options = opt::Opt {
323ded9b81dSJeremy L Thompson ceed_spec: "/cpu/self/ref/serial".to_string(),
324ded9b81dSJeremy L Thompson dim: 3,
325ded9b81dSJeremy L Thompson mesh_degree: 4,
326ded9b81dSJeremy L Thompson solution_degree: 4,
327ded9b81dSJeremy L Thompson num_qpts: 6,
328ded9b81dSJeremy L Thompson problem_size_requested: -1,
329ded9b81dSJeremy L Thompson test: true,
330954a6033SJeremy L Thompson quiet: false,
331ded9b81dSJeremy L Thompson gallery: false,
332ded9b81dSJeremy L Thompson };
333ded9b81dSJeremy L Thompson assert!(example_1(options).is_ok());
334ded9b81dSJeremy L Thompson }
335ded9b81dSJeremy L Thompson
336ded9b81dSJeremy L Thompson #[test]
example_1_1d_gallery()337ded9b81dSJeremy L Thompson fn example_1_1d_gallery() {
338ded9b81dSJeremy L Thompson let options = opt::Opt {
339ded9b81dSJeremy L Thompson ceed_spec: "/cpu/self/ref/serial".to_string(),
340ded9b81dSJeremy L Thompson dim: 1,
341ded9b81dSJeremy L Thompson mesh_degree: 4,
342ded9b81dSJeremy L Thompson solution_degree: 4,
343ded9b81dSJeremy L Thompson num_qpts: 6,
344ded9b81dSJeremy L Thompson problem_size_requested: -1,
345ded9b81dSJeremy L Thompson test: true,
346ded9b81dSJeremy L Thompson quiet: true,
347ded9b81dSJeremy L Thompson gallery: true,
348ded9b81dSJeremy L Thompson };
349ded9b81dSJeremy L Thompson assert!(example_1(options).is_ok());
350ded9b81dSJeremy L Thompson }
351ded9b81dSJeremy L Thompson
352ded9b81dSJeremy L Thompson #[test]
example_1_2d_gallery()353ded9b81dSJeremy L Thompson fn example_1_2d_gallery() {
354ded9b81dSJeremy L Thompson let options = opt::Opt {
355ded9b81dSJeremy L Thompson ceed_spec: "/cpu/self/ref/serial".to_string(),
356ded9b81dSJeremy L Thompson dim: 2,
357ded9b81dSJeremy L Thompson mesh_degree: 4,
358ded9b81dSJeremy L Thompson solution_degree: 4,
359ded9b81dSJeremy L Thompson num_qpts: 6,
360ded9b81dSJeremy L Thompson problem_size_requested: -1,
361ded9b81dSJeremy L Thompson test: true,
362ded9b81dSJeremy L Thompson quiet: true,
363ded9b81dSJeremy L Thompson gallery: true,
364ded9b81dSJeremy L Thompson };
365ded9b81dSJeremy L Thompson assert!(example_1(options).is_ok());
366ded9b81dSJeremy L Thompson }
367ded9b81dSJeremy L Thompson
368ded9b81dSJeremy L Thompson #[test]
example_1_3d_gallery()369ded9b81dSJeremy L Thompson fn example_1_3d_gallery() {
370ded9b81dSJeremy L Thompson let options = opt::Opt {
371ded9b81dSJeremy L Thompson ceed_spec: "/cpu/self/ref/serial".to_string(),
372ded9b81dSJeremy L Thompson dim: 3,
373ded9b81dSJeremy L Thompson mesh_degree: 4,
374ded9b81dSJeremy L Thompson solution_degree: 4,
375ded9b81dSJeremy L Thompson num_qpts: 6,
376ded9b81dSJeremy L Thompson problem_size_requested: -1,
377ded9b81dSJeremy L Thompson test: true,
378ded9b81dSJeremy L Thompson quiet: true,
379ded9b81dSJeremy L Thompson gallery: true,
380ded9b81dSJeremy L Thompson };
381ded9b81dSJeremy L Thompson assert!(example_1(options).is_ok());
382ded9b81dSJeremy L Thompson }
383ded9b81dSJeremy L Thompson }
384ded9b81dSJeremy L Thompson
385ded9b81dSJeremy L Thompson // ----------------------------------------------------------------------------
386