xref: /libCEED/examples/rust/ex1-volume/src/main.rs (revision e15f9bd09af0280c89b79924fa9af7dd2e3e30be)
1ded9b81dSJeremy L Thompson // Copyright (c) 2017-2021, Lawrence Livermore National Security, LLC.
2ded9b81dSJeremy L Thompson // Produced at the Lawrence Livermore National Laboratory. LLNL-CODE-734707.
3ded9b81dSJeremy L Thompson // All Rights reserved. See files LICENSE and NOTICE for details.
4ded9b81dSJeremy L Thompson //
5ded9b81dSJeremy L Thompson // This file is part of CEED, a collection of benchmarks, miniapps, software
6ded9b81dSJeremy L Thompson // libraries and APIs for efficient high-order finite element and spectral
7ded9b81dSJeremy L Thompson // element discretizations for exascale applications. For more information and
8ded9b81dSJeremy L Thompson // source code availability see http://github.com/ceed.
9ded9b81dSJeremy L Thompson //
10ded9b81dSJeremy L Thompson // The CEED research is supported by the Exascale Computing Project 17-SC-20-SC,
11ded9b81dSJeremy L Thompson // a collaborative effort of two U.S. Department of Energy organizations (Office
12ded9b81dSJeremy L Thompson // of Science and the National Nuclear Security Administration) responsible for
13ded9b81dSJeremy L Thompson // the planning and preparation of a capable exascale ecosystem, including
14ded9b81dSJeremy L Thompson // software, applications, hardware, advanced system engineering and early
15ded9b81dSJeremy L Thompson // testbed platforms, in support of the nation's exascale computing imperative.
16ded9b81dSJeremy L Thompson 
17ded9b81dSJeremy L Thompson //                             libCEED Example 1
18ded9b81dSJeremy L Thompson //
19ded9b81dSJeremy L Thompson // This example illustrates a simple usage of libCEED to compute the volume of a
20ded9b81dSJeremy L Thompson // 3D body using matrix-free application of a mass operator.  Arbitrary mesh and
21ded9b81dSJeremy L Thompson // solution orders in 1D, 2D and 3D are supported from the same code.
22ded9b81dSJeremy L Thompson //
23ded9b81dSJeremy L Thompson // The example has no dependencies, and is designed to be self-contained. For
24ded9b81dSJeremy L Thompson // additional examples that use external discretization libraries (MFEM, PETSc,
25ded9b81dSJeremy L Thompson // etc.) see the subdirectories in libceed/examples.
26ded9b81dSJeremy L Thompson //
27ded9b81dSJeremy L Thompson // All libCEED objects use a Ceed device object constructed based on a command
28ded9b81dSJeremy L Thompson // line argument (-ceed).
29ded9b81dSJeremy L Thompson 
30ded9b81dSJeremy L Thompson use libceed::{prelude::*, Ceed};
31ded9b81dSJeremy L Thompson use mesh;
32ded9b81dSJeremy L Thompson use structopt::StructOpt;
33ded9b81dSJeremy L Thompson 
34ded9b81dSJeremy L Thompson mod opt;
35ded9b81dSJeremy L Thompson mod transform;
36ded9b81dSJeremy L Thompson 
37ded9b81dSJeremy L Thompson // ----------------------------------------------------------------------------
38ded9b81dSJeremy L Thompson // Example 1
39ded9b81dSJeremy L Thompson // ----------------------------------------------------------------------------
40954a6033SJeremy L Thompson #[cfg(not(tarpaulin_include))]
41*e15f9bd0SJeremy L Thompson fn main() -> Result<(), libceed::CeedError> {
42ded9b81dSJeremy L Thompson     let options = opt::Opt::from_args();
43ded9b81dSJeremy L Thompson     example_1(options)
44ded9b81dSJeremy L Thompson }
45ded9b81dSJeremy L Thompson 
46*e15f9bd0SJeremy L Thompson fn example_1(options: opt::Opt) -> Result<(), libceed::CeedError> {
47ded9b81dSJeremy L Thompson     // Process command line arguments
48ded9b81dSJeremy L Thompson     let opt::Opt {
49ded9b81dSJeremy L Thompson         ceed_spec,
50ded9b81dSJeremy L Thompson         dim,
51ded9b81dSJeremy L Thompson         mesh_degree,
52ded9b81dSJeremy L Thompson         solution_degree,
53ded9b81dSJeremy L Thompson         num_qpts,
54ded9b81dSJeremy L Thompson         problem_size_requested,
55ded9b81dSJeremy L Thompson         test,
56ded9b81dSJeremy L Thompson         quiet,
57ded9b81dSJeremy L Thompson         gallery,
58ded9b81dSJeremy L Thompson     } = options;
59ded9b81dSJeremy L Thompson     assert!(dim >= 1 && dim <= 3);
60ded9b81dSJeremy L Thompson     assert!(mesh_degree >= 1);
61ded9b81dSJeremy L Thompson     assert!(solution_degree >= 1);
62ded9b81dSJeremy L Thompson     assert!(num_qpts >= 1);
63ded9b81dSJeremy L Thompson     let ncomp_x = dim;
64ded9b81dSJeremy L Thompson     let problem_size: i64;
65ded9b81dSJeremy L Thompson     if problem_size_requested < 0 {
66ded9b81dSJeremy L Thompson         problem_size = if test { 8 * 16 } else { 256 * 1024 };
67ded9b81dSJeremy L Thompson     } else {
68ded9b81dSJeremy L Thompson         problem_size = problem_size_requested;
69ded9b81dSJeremy L Thompson     }
70ded9b81dSJeremy L Thompson 
71ded9b81dSJeremy L Thompson     // Summary output
72ded9b81dSJeremy L Thompson     if !quiet {
73ded9b81dSJeremy L Thompson         println!("Selected options: [command line option] : <current value>");
74ded9b81dSJeremy L Thompson         println!("    Ceed specification [-c] : {}", ceed_spec);
75ded9b81dSJeremy L Thompson         println!("    Mesh dimension     [-d] : {}", dim);
76ded9b81dSJeremy L Thompson         println!("    Mesh degree        [-m] : {}", mesh_degree);
77ded9b81dSJeremy L Thompson         println!("    Solution degree    [-p] : {}", solution_degree);
78ded9b81dSJeremy L Thompson         println!("    Num. 1D quadr. pts [-q] : {}", num_qpts);
79ded9b81dSJeremy L Thompson         println!("    Approx. # unknowns [-s] : {}", problem_size);
80ded9b81dSJeremy L Thompson         println!(
81ded9b81dSJeremy L Thompson             "    QFunction source   [-g] : {}",
82ded9b81dSJeremy L Thompson             if gallery { "gallery" } else { "user closure" }
83ded9b81dSJeremy L Thompson         );
84ded9b81dSJeremy L Thompson     }
85ded9b81dSJeremy L Thompson 
86ded9b81dSJeremy L Thompson     // Initalize ceed context
87ded9b81dSJeremy L Thompson     let ceed = Ceed::init(&ceed_spec);
88ded9b81dSJeremy L Thompson 
89ded9b81dSJeremy L Thompson     // Mesh and solution bases
90*e15f9bd0SJeremy L Thompson     let basis_mesh = ceed
91*e15f9bd0SJeremy L Thompson         .basis_tensor_H1_Lagrange(dim, ncomp_x, mesh_degree + 1, num_qpts, QuadMode::Gauss)
92*e15f9bd0SJeremy L Thompson         .unwrap();
93*e15f9bd0SJeremy L Thompson     let basis_solution = ceed
94*e15f9bd0SJeremy L Thompson         .basis_tensor_H1_Lagrange(dim, 1, solution_degree + 1, num_qpts, QuadMode::Gauss)
95*e15f9bd0SJeremy L Thompson         .unwrap();
96ded9b81dSJeremy L Thompson 
97ded9b81dSJeremy L Thompson     // Determine mesh size from approximate problem size
98ded9b81dSJeremy L Thompson     let num_xyz = mesh::cartesian_mesh_size(dim, solution_degree, problem_size);
99ded9b81dSJeremy L Thompson     if !quiet {
100ded9b81dSJeremy L Thompson         print!("\nMesh size                   : nx = {}", num_xyz[0]);
101ded9b81dSJeremy L Thompson         if dim > 1 {
102ded9b81dSJeremy L Thompson             print!(", ny = {}", num_xyz[1]);
103ded9b81dSJeremy L Thompson         }
104ded9b81dSJeremy L Thompson         if dim > 2 {
105ded9b81dSJeremy L Thompson             print!(", nz = {}", num_xyz[2]);
106ded9b81dSJeremy L Thompson         }
107ded9b81dSJeremy L Thompson         print!("\n");
108ded9b81dSJeremy L Thompson     }
109ded9b81dSJeremy L Thompson 
110ded9b81dSJeremy L Thompson     // Build ElemRestriction objects describing the mesh and solution discrete
111ded9b81dSJeremy L Thompson     // representations
112ded9b81dSJeremy L Thompson     let (restr_mesh, _) =
113*e15f9bd0SJeremy L Thompson         mesh::build_cartesian_restriction(&ceed, dim, num_xyz, mesh_degree, ncomp_x, num_qpts)?;
114ded9b81dSJeremy L Thompson     let (restr_solution, restr_qdata) =
115*e15f9bd0SJeremy L Thompson         mesh::build_cartesian_restriction(&ceed, dim, num_xyz, solution_degree, 1, num_qpts)?;
116ded9b81dSJeremy L Thompson     let mesh_size = restr_mesh.lvector_size();
117ded9b81dSJeremy L Thompson     let solution_size = restr_solution.lvector_size();
118ded9b81dSJeremy L Thompson     if !quiet {
119ded9b81dSJeremy L Thompson         println!("Number of mesh nodes        : {}", mesh_size / dim);
120ded9b81dSJeremy L Thompson         println!("Number of solution nodes    : {}", solution_size);
121ded9b81dSJeremy L Thompson     }
122ded9b81dSJeremy L Thompson 
123ded9b81dSJeremy L Thompson     // Create a Vector with the mesh coordinates
124*e15f9bd0SJeremy L Thompson     let mut mesh_coords = mesh::cartesian_mesh_coords(&ceed, dim, num_xyz, mesh_degree, mesh_size)?;
125ded9b81dSJeremy L Thompson 
126ded9b81dSJeremy L Thompson     // Apply a transformation to the mesh coordinates
127ded9b81dSJeremy L Thompson     let exact_volume = transform::transform_mesh_coordinates(dim, mesh_size, &mut mesh_coords);
128ded9b81dSJeremy L Thompson 
129ded9b81dSJeremy L Thompson     // QFunction that builds the quadrature data for the mass operator
130ded9b81dSJeremy L Thompson     // -- QFunction from user closure
131ded9b81dSJeremy L Thompson     let build_mass = move |[jacobian, weights, ..]: QFunctionInputs,
132ded9b81dSJeremy L Thompson                            [qdata, ..]: QFunctionOutputs| {
133ded9b81dSJeremy L Thompson         // Build quadrature data
134ded9b81dSJeremy L Thompson         match dim {
135ded9b81dSJeremy L Thompson             1 => qdata
136ded9b81dSJeremy L Thompson                 .iter_mut()
137ded9b81dSJeremy L Thompson                 .zip(jacobian.iter().zip(weights.iter()))
138ded9b81dSJeremy L Thompson                 .for_each(|(qdata, (j, weight))| *qdata = j * weight),
139ded9b81dSJeremy L Thompson             2 => {
140ded9b81dSJeremy L Thompson                 let q = qdata.len();
141ded9b81dSJeremy L Thompson                 qdata.iter_mut().zip(weights.iter()).enumerate().for_each(
142ded9b81dSJeremy L Thompson                     |(i, (qdata, weight))| {
143ded9b81dSJeremy L Thompson                         *qdata = (jacobian[i + q * 0] * jacobian[i + q * 3]
144ded9b81dSJeremy L Thompson                             - jacobian[i + q * 1] * jacobian[i + q * 2])
145ded9b81dSJeremy L Thompson                             * weight
146ded9b81dSJeremy L Thompson                     },
147ded9b81dSJeremy L Thompson                 );
148ded9b81dSJeremy L Thompson             }
149ded9b81dSJeremy L Thompson             3 => {
150ded9b81dSJeremy L Thompson                 let q = qdata.len();
151ded9b81dSJeremy L Thompson                 qdata.iter_mut().zip(weights.iter()).enumerate().for_each(
152ded9b81dSJeremy L Thompson                     |(i, (qdata, weight))| {
153ded9b81dSJeremy L Thompson                         *qdata = (jacobian[i + q * 0]
154ded9b81dSJeremy L Thompson                             * (jacobian[i + q * 4] * jacobian[i + q * 8]
155ded9b81dSJeremy L Thompson                                 - jacobian[i + q * 5] * jacobian[i + q * 7])
156ded9b81dSJeremy L Thompson                             - jacobian[i + q * 1]
157ded9b81dSJeremy L Thompson                                 * (jacobian[i + q * 3] * jacobian[i + q * 8]
158ded9b81dSJeremy L Thompson                                     - jacobian[i + q * 5] * jacobian[i + q * 6])
159ded9b81dSJeremy L Thompson                             + jacobian[i + q * 2]
160ded9b81dSJeremy L Thompson                                 * (jacobian[i + q * 3] * jacobian[i + q * 7]
161ded9b81dSJeremy L Thompson                                     - jacobian[i + q * 4] * jacobian[i + q * 6]))
162ded9b81dSJeremy L Thompson                             * weight
163ded9b81dSJeremy L Thompson                     },
164ded9b81dSJeremy L Thompson                 );
165ded9b81dSJeremy L Thompson             }
166ded9b81dSJeremy L Thompson             _ => unreachable!(),
167ded9b81dSJeremy L Thompson         };
168ded9b81dSJeremy L Thompson 
169ded9b81dSJeremy L Thompson         // Return clean error code
170ded9b81dSJeremy L Thompson         0
171ded9b81dSJeremy L Thompson     };
172ded9b81dSJeremy L Thompson     let qf_build_closure = ceed
173*e15f9bd0SJeremy L Thompson         .q_function_interior(1, Box::new(build_mass))?
174*e15f9bd0SJeremy L Thompson         .input("dx", ncomp_x * dim, EvalMode::Grad)?
175*e15f9bd0SJeremy L Thompson         .input("weights", 1, EvalMode::Weight)?
176*e15f9bd0SJeremy L Thompson         .output("qdata", 1, EvalMode::None)?;
177ded9b81dSJeremy L Thompson     // -- QFunction from gallery
178ded9b81dSJeremy L Thompson     let qf_build_named = {
179ded9b81dSJeremy L Thompson         let name = format!("Mass{}DBuild", dim);
180*e15f9bd0SJeremy L Thompson         ceed.q_function_interior_by_name(&name)?
181ded9b81dSJeremy L Thompson     };
182ded9b81dSJeremy L Thompson     // -- QFunction for use with Operator
183ded9b81dSJeremy L Thompson     let qf_build = if gallery {
184ded9b81dSJeremy L Thompson         QFunctionOpt::SomeQFunctionByName(&qf_build_named)
185ded9b81dSJeremy L Thompson     } else {
186ded9b81dSJeremy L Thompson         QFunctionOpt::SomeQFunction(&qf_build_closure)
187ded9b81dSJeremy L Thompson     };
188ded9b81dSJeremy L Thompson 
189ded9b81dSJeremy L Thompson     // Operator that build the quadrature data for the mass operator
190ded9b81dSJeremy L Thompson     let op_build = ceed
191*e15f9bd0SJeremy L Thompson         .operator(qf_build, QFunctionOpt::None, QFunctionOpt::None)?
192*e15f9bd0SJeremy L Thompson         .field("dx", &restr_mesh, &basis_mesh, VectorOpt::Active)?
193ded9b81dSJeremy L Thompson         .field(
194ded9b81dSJeremy L Thompson             "weights",
195ded9b81dSJeremy L Thompson             ElemRestrictionOpt::None,
196ded9b81dSJeremy L Thompson             &basis_mesh,
197ded9b81dSJeremy L Thompson             VectorOpt::None,
198*e15f9bd0SJeremy L Thompson         )?
199ded9b81dSJeremy L Thompson         .field(
200ded9b81dSJeremy L Thompson             "qdata",
201ded9b81dSJeremy L Thompson             &restr_qdata,
202ded9b81dSJeremy L Thompson             BasisOpt::Collocated,
203ded9b81dSJeremy L Thompson             VectorOpt::Active,
204*e15f9bd0SJeremy L Thompson         )?;
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();
209*e15f9bd0SJeremy L Thompson     let mut qdata = ceed.vector(num_elem * elem_qpts)?;
210*e15f9bd0SJeremy 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
223*e15f9bd0SJeremy L Thompson         .q_function_interior(1, Box::new(apply_mass))?
224*e15f9bd0SJeremy L Thompson         .input("u", 1, EvalMode::Interp)?
225*e15f9bd0SJeremy L Thompson         .input("qdata", 1, EvalMode::None)?
226*e15f9bd0SJeremy L Thompson         .output("v", 1, EvalMode::Interp)?;
227ded9b81dSJeremy L Thompson     // -- QFunction from gallery
228*e15f9bd0SJeremy 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
238*e15f9bd0SJeremy L Thompson         .operator(qf_mass, QFunctionOpt::None, QFunctionOpt::None)?
239*e15f9bd0SJeremy L Thompson         .field("u", &restr_solution, &basis_solution, VectorOpt::Active)?
240*e15f9bd0SJeremy L Thompson         .field("qdata", &restr_qdata, BasisOpt::Collocated, &qdata)?
241*e15f9bd0SJeremy L Thompson         .field("v", &restr_solution, &basis_solution, VectorOpt::Active)?;
242ded9b81dSJeremy L Thompson 
243ded9b81dSJeremy L Thompson     // Solution vectors
244*e15f9bd0SJeremy L Thompson     let u = ceed.vector_from_slice(&vec![1.0; solution_size])?;
245*e15f9bd0SJeremy L Thompson     let mut v = ceed.vector(solution_size)?;
246ded9b81dSJeremy L Thompson 
247ded9b81dSJeremy L Thompson     // Apply the mass operator
248*e15f9bd0SJeremy L Thompson     op_mass.apply(&u, &mut v)?;
249ded9b81dSJeremy L Thompson 
250ded9b81dSJeremy L Thompson     // Compute the mesh volume
251ded9b81dSJeremy L Thompson     let volume: f64 = v.view().iter().sum();
252ded9b81dSJeremy L Thompson 
253ded9b81dSJeremy L Thompson     // Output results
254ded9b81dSJeremy L Thompson     if !quiet {
255ded9b81dSJeremy L Thompson         println!("Exact mesh volume           : {:.12}", exact_volume);
256ded9b81dSJeremy L Thompson         println!("Computed mesh volume        : {:.12}", volume);
257ded9b81dSJeremy L Thompson         println!(
258ded9b81dSJeremy L Thompson             "Volume error                : {:.12e}",
259ded9b81dSJeremy L Thompson             volume - exact_volume
260ded9b81dSJeremy L Thompson         );
261ded9b81dSJeremy L Thompson     }
262ded9b81dSJeremy L Thompson     let tolerance = match dim {
263ded9b81dSJeremy L Thompson         1 => 1E-14,
264ded9b81dSJeremy L Thompson         2 => 1E-7,
265ded9b81dSJeremy L Thompson         _ => 1E-5,
266ded9b81dSJeremy L Thompson     };
267ded9b81dSJeremy L Thompson     let error = (volume - exact_volume).abs();
268ded9b81dSJeremy L Thompson     if error > tolerance {
269ded9b81dSJeremy L Thompson         println!("Volume error too large: {:.12e}", error);
270*e15f9bd0SJeremy L Thompson         return Err(libceed::CeedError {
271*e15f9bd0SJeremy L Thompson             message: format!(
272ded9b81dSJeremy L Thompson                 "Volume error too large - expected: {:.12e}, actual: {:.12e}",
273ded9b81dSJeremy L Thompson                 tolerance, error
274*e15f9bd0SJeremy L Thompson             ),
275*e15f9bd0SJeremy L Thompson         });
276ded9b81dSJeremy L Thompson     }
277ded9b81dSJeremy L Thompson     Ok(())
278ded9b81dSJeremy L Thompson }
279ded9b81dSJeremy L Thompson 
280ded9b81dSJeremy L Thompson // ----------------------------------------------------------------------------
281ded9b81dSJeremy L Thompson // Tests
282ded9b81dSJeremy L Thompson // ----------------------------------------------------------------------------
283ded9b81dSJeremy L Thompson #[cfg(test)]
284ded9b81dSJeremy L Thompson mod tests {
285ded9b81dSJeremy L Thompson     use super::*;
286ded9b81dSJeremy L Thompson 
287ded9b81dSJeremy L Thompson     #[test]
288ded9b81dSJeremy L Thompson     fn example_1_1d() {
289ded9b81dSJeremy L Thompson         let options = opt::Opt {
290ded9b81dSJeremy L Thompson             ceed_spec: "/cpu/self/ref/serial".to_string(),
291ded9b81dSJeremy L Thompson             dim: 1,
292ded9b81dSJeremy L Thompson             mesh_degree: 4,
293ded9b81dSJeremy L Thompson             solution_degree: 4,
294ded9b81dSJeremy L Thompson             num_qpts: 6,
295ded9b81dSJeremy L Thompson             problem_size_requested: -1,
296ded9b81dSJeremy L Thompson             test: true,
297ded9b81dSJeremy L Thompson             quiet: true,
298ded9b81dSJeremy L Thompson             gallery: false,
299ded9b81dSJeremy L Thompson         };
300ded9b81dSJeremy L Thompson         assert!(example_1(options).is_ok());
301ded9b81dSJeremy L Thompson     }
302ded9b81dSJeremy L Thompson 
303ded9b81dSJeremy L Thompson     #[test]
304ded9b81dSJeremy L Thompson     fn example_1_2d() {
305ded9b81dSJeremy L Thompson         let options = opt::Opt {
306ded9b81dSJeremy L Thompson             ceed_spec: "/cpu/self/ref/serial".to_string(),
307ded9b81dSJeremy L Thompson             dim: 2,
308ded9b81dSJeremy L Thompson             mesh_degree: 4,
309ded9b81dSJeremy L Thompson             solution_degree: 4,
310ded9b81dSJeremy L Thompson             num_qpts: 6,
311ded9b81dSJeremy L Thompson             problem_size_requested: -1,
312ded9b81dSJeremy L Thompson             test: true,
313ded9b81dSJeremy L Thompson             quiet: true,
314ded9b81dSJeremy L Thompson             gallery: false,
315ded9b81dSJeremy L Thompson         };
316ded9b81dSJeremy L Thompson         assert!(example_1(options).is_ok());
317ded9b81dSJeremy L Thompson     }
318ded9b81dSJeremy L Thompson 
319ded9b81dSJeremy L Thompson     #[test]
320ded9b81dSJeremy L Thompson     fn example_1_3d() {
321ded9b81dSJeremy L Thompson         let options = opt::Opt {
322ded9b81dSJeremy L Thompson             ceed_spec: "/cpu/self/ref/serial".to_string(),
323ded9b81dSJeremy L Thompson             dim: 3,
324ded9b81dSJeremy L Thompson             mesh_degree: 4,
325ded9b81dSJeremy L Thompson             solution_degree: 4,
326ded9b81dSJeremy L Thompson             num_qpts: 6,
327ded9b81dSJeremy L Thompson             problem_size_requested: -1,
328ded9b81dSJeremy L Thompson             test: true,
329954a6033SJeremy L Thompson             quiet: false,
330ded9b81dSJeremy L Thompson             gallery: false,
331ded9b81dSJeremy L Thompson         };
332ded9b81dSJeremy L Thompson         assert!(example_1(options).is_ok());
333ded9b81dSJeremy L Thompson     }
334ded9b81dSJeremy L Thompson 
335ded9b81dSJeremy L Thompson     #[test]
336ded9b81dSJeremy L Thompson     fn example_1_1d_gallery() {
337ded9b81dSJeremy L Thompson         let options = opt::Opt {
338ded9b81dSJeremy L Thompson             ceed_spec: "/cpu/self/ref/serial".to_string(),
339ded9b81dSJeremy L Thompson             dim: 1,
340ded9b81dSJeremy L Thompson             mesh_degree: 4,
341ded9b81dSJeremy L Thompson             solution_degree: 4,
342ded9b81dSJeremy L Thompson             num_qpts: 6,
343ded9b81dSJeremy L Thompson             problem_size_requested: -1,
344ded9b81dSJeremy L Thompson             test: true,
345ded9b81dSJeremy L Thompson             quiet: true,
346ded9b81dSJeremy L Thompson             gallery: true,
347ded9b81dSJeremy L Thompson         };
348ded9b81dSJeremy L Thompson         assert!(example_1(options).is_ok());
349ded9b81dSJeremy L Thompson     }
350ded9b81dSJeremy L Thompson 
351ded9b81dSJeremy L Thompson     #[test]
352ded9b81dSJeremy L Thompson     fn example_1_2d_gallery() {
353ded9b81dSJeremy L Thompson         let options = opt::Opt {
354ded9b81dSJeremy L Thompson             ceed_spec: "/cpu/self/ref/serial".to_string(),
355ded9b81dSJeremy L Thompson             dim: 2,
356ded9b81dSJeremy L Thompson             mesh_degree: 4,
357ded9b81dSJeremy L Thompson             solution_degree: 4,
358ded9b81dSJeremy L Thompson             num_qpts: 6,
359ded9b81dSJeremy L Thompson             problem_size_requested: -1,
360ded9b81dSJeremy L Thompson             test: true,
361ded9b81dSJeremy L Thompson             quiet: true,
362ded9b81dSJeremy L Thompson             gallery: true,
363ded9b81dSJeremy L Thompson         };
364ded9b81dSJeremy L Thompson         assert!(example_1(options).is_ok());
365ded9b81dSJeremy L Thompson     }
366ded9b81dSJeremy L Thompson 
367ded9b81dSJeremy L Thompson     #[test]
368ded9b81dSJeremy L Thompson     fn example_1_3d_gallery() {
369ded9b81dSJeremy L Thompson         let options = opt::Opt {
370ded9b81dSJeremy L Thompson             ceed_spec: "/cpu/self/ref/serial".to_string(),
371ded9b81dSJeremy L Thompson             dim: 3,
372ded9b81dSJeremy L Thompson             mesh_degree: 4,
373ded9b81dSJeremy L Thompson             solution_degree: 4,
374ded9b81dSJeremy L Thompson             num_qpts: 6,
375ded9b81dSJeremy L Thompson             problem_size_requested: -1,
376ded9b81dSJeremy L Thompson             test: true,
377ded9b81dSJeremy L Thompson             quiet: true,
378ded9b81dSJeremy L Thompson             gallery: true,
379ded9b81dSJeremy L Thompson         };
380ded9b81dSJeremy L Thompson         assert!(example_1(options).is_ok());
381ded9b81dSJeremy L Thompson     }
382ded9b81dSJeremy L Thompson }
383ded9b81dSJeremy L Thompson 
384ded9b81dSJeremy L Thompson // ----------------------------------------------------------------------------
385