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 2 18ded9b81dSJeremy L Thompson // 19ded9b81dSJeremy L Thompson // This example illustrates a simple usage of libCEED to compute the surface 20ded9b81dSJeremy L Thompson // area of a 3D body using matrix-free application of a diffusion operator. 21ded9b81dSJeremy L Thompson // Arbitrary mesh and solution degrees in 1D, 2D and 3D are supported from the 22ded9b81dSJeremy L Thompson // same code. 23ded9b81dSJeremy L Thompson // 24ded9b81dSJeremy L Thompson // The example has no dependencies, and is designed to be self-contained. For 25ded9b81dSJeremy L Thompson // additional examples that use external discretization libraries (MFEM, PETSc, 26ded9b81dSJeremy L Thompson // etc.) see the subdirectories in libceed/examples. 27ded9b81dSJeremy L Thompson // 28ded9b81dSJeremy L Thompson // All libCEED objects use a Ceed device object constructed based on a command 29ded9b81dSJeremy L Thompson // line argument (-ceed). 30ded9b81dSJeremy L Thompson 31ded9b81dSJeremy L Thompson use libceed::{prelude::*, Ceed}; 32ded9b81dSJeremy L Thompson use mesh; 33ded9b81dSJeremy L Thompson use structopt::StructOpt; 34ded9b81dSJeremy L Thompson 35ded9b81dSJeremy L Thompson mod opt; 36ded9b81dSJeremy L Thompson mod transform; 37ded9b81dSJeremy L Thompson 38ded9b81dSJeremy L Thompson // ---------------------------------------------------------------------------- 39954a6033SJeremy L Thompson // Example 2 40ded9b81dSJeremy L Thompson // ---------------------------------------------------------------------------- 41954a6033SJeremy L Thompson #[cfg(not(tarpaulin_include))] 4289d15d5fSJeremy L Thompson fn main() -> libceed::Result<()> { 43ded9b81dSJeremy L Thompson let options = opt::Opt::from_args(); 44ded9b81dSJeremy L Thompson example_2(options) 45ded9b81dSJeremy L Thompson } 46ded9b81dSJeremy L Thompson 4789d15d5fSJeremy L Thompson fn example_2(options: opt::Opt) -> libceed::Result<()> { 48ded9b81dSJeremy L Thompson // Process command line arguments 49ded9b81dSJeremy L Thompson let opt::Opt { 50ded9b81dSJeremy L Thompson ceed_spec, 51ded9b81dSJeremy L Thompson dim, 52ded9b81dSJeremy L Thompson mesh_degree, 53ded9b81dSJeremy L Thompson solution_degree, 54ded9b81dSJeremy L Thompson num_qpts, 55ded9b81dSJeremy L Thompson problem_size_requested, 56ded9b81dSJeremy L Thompson test, 57ded9b81dSJeremy L Thompson quiet, 58ded9b81dSJeremy L Thompson gallery, 59ded9b81dSJeremy L Thompson } = options; 60ded9b81dSJeremy L Thompson assert!(dim >= 1 && dim <= 3); 61ded9b81dSJeremy L Thompson assert!(mesh_degree >= 1); 62ded9b81dSJeremy L Thompson assert!(solution_degree >= 1); 63ded9b81dSJeremy L Thompson assert!(num_qpts >= 1); 64ded9b81dSJeremy L Thompson let ncomp_x = dim; 65ded9b81dSJeremy L Thompson let problem_size: i64; 66ded9b81dSJeremy L Thompson if problem_size_requested < 0 { 67ded9b81dSJeremy L Thompson problem_size = if test { 68ded9b81dSJeremy L Thompson 16 * 16 * (dim * dim) as i64 69ded9b81dSJeremy L Thompson } else { 70ded9b81dSJeremy L Thompson 256 * 1024 71ded9b81dSJeremy L Thompson }; 72ded9b81dSJeremy L Thompson } else { 73ded9b81dSJeremy L Thompson problem_size = problem_size_requested; 74ded9b81dSJeremy L Thompson } 75ded9b81dSJeremy L Thompson 76ded9b81dSJeremy L Thompson // Summary output 77ded9b81dSJeremy L Thompson if !quiet { 78ded9b81dSJeremy L Thompson println!("Selected options: [command line option] : <current value>"); 79ded9b81dSJeremy L Thompson println!(" Ceed specification [-c] : {}", ceed_spec); 80ded9b81dSJeremy L Thompson println!(" Mesh dimension [-d] : {}", dim); 81ded9b81dSJeremy L Thompson println!(" Mesh degree [-m] : {}", mesh_degree); 82ded9b81dSJeremy L Thompson println!(" Solution degree [-p] : {}", solution_degree); 83ded9b81dSJeremy L Thompson println!(" Num. 1D quadr. pts [-q] : {}", num_qpts); 84ded9b81dSJeremy L Thompson println!(" Approx. # unknowns [-s] : {}", problem_size); 85ded9b81dSJeremy L Thompson println!( 86ded9b81dSJeremy L Thompson " QFunction source [-g] : {}", 87ded9b81dSJeremy L Thompson if gallery { "gallery" } else { "user closure" } 88ded9b81dSJeremy L Thompson ); 89ded9b81dSJeremy L Thompson } 90ded9b81dSJeremy L Thompson 91ded9b81dSJeremy L Thompson // Initalize ceed context 92ded9b81dSJeremy L Thompson let ceed = Ceed::init(&ceed_spec); 93ded9b81dSJeremy L Thompson 94ded9b81dSJeremy L Thompson // Mesh and solution bases 95*51d03428SJeremy L Thompson let basis_mesh = 96*51d03428SJeremy L Thompson ceed.basis_tensor_H1_Lagrange(dim, ncomp_x, mesh_degree + 1, num_qpts, QuadMode::Gauss)?; 97*51d03428SJeremy L Thompson let basis_solution = 98*51d03428SJeremy L Thompson ceed.basis_tensor_H1_Lagrange(dim, 1, solution_degree + 1, num_qpts, QuadMode::Gauss)?; 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 } 110ded9b81dSJeremy L Thompson print!("\n"); 111ded9b81dSJeremy L Thompson } 112ded9b81dSJeremy L Thompson 113ded9b81dSJeremy L Thompson // Build ElemRestriction objects describing the mesh and solution discrete 114ded9b81dSJeremy L Thompson // representations 115ded9b81dSJeremy L Thompson let (restr_mesh, _) = 116e15f9bd0SJeremy L Thompson mesh::build_cartesian_restriction(&ceed, dim, num_xyz, mesh_degree, ncomp_x, num_qpts)?; 117ded9b81dSJeremy L Thompson let (_, restr_qdata) = mesh::build_cartesian_restriction( 118ded9b81dSJeremy L Thompson &ceed, 119ded9b81dSJeremy L Thompson dim, 120ded9b81dSJeremy L Thompson num_xyz, 121ded9b81dSJeremy L Thompson solution_degree, 122ded9b81dSJeremy L Thompson dim * (dim + 1) / 2, 123ded9b81dSJeremy L Thompson num_qpts, 124e15f9bd0SJeremy L Thompson )?; 125ded9b81dSJeremy L Thompson 126ded9b81dSJeremy L Thompson let (restr_solution, _) = 127e15f9bd0SJeremy L Thompson mesh::build_cartesian_restriction(&ceed, dim, num_xyz, solution_degree, 1, num_qpts)?; 128ded9b81dSJeremy L Thompson let mesh_size = restr_mesh.lvector_size(); 129ded9b81dSJeremy L Thompson let solution_size = restr_solution.lvector_size(); 130ded9b81dSJeremy L Thompson if !quiet { 131ded9b81dSJeremy L Thompson println!("Number of mesh nodes : {}", mesh_size / dim); 132ded9b81dSJeremy L Thompson println!("Number of solution nodes : {}", solution_size); 133ded9b81dSJeremy L Thompson } 134ded9b81dSJeremy L Thompson 135ded9b81dSJeremy L Thompson // Create a Vector with the mesh coordinates 136e15f9bd0SJeremy L Thompson let mut mesh_coords = mesh::cartesian_mesh_coords(&ceed, dim, num_xyz, mesh_degree, mesh_size)?; 137ded9b81dSJeremy L Thompson 138ded9b81dSJeremy L Thompson // Apply a transformation to the mesh coordinates 139e78171edSJeremy L Thompson let exact_area = transform::transform_mesh_coordinates(dim, &mut mesh_coords)?; 140ded9b81dSJeremy L Thompson 141ded9b81dSJeremy L Thompson // QFunction that builds the quadrature data for the diff operator 142ded9b81dSJeremy L Thompson // -- QFunction from user closure 143ded9b81dSJeremy L Thompson let build_diff = move |[jacobian, weights, ..]: QFunctionInputs, 144ded9b81dSJeremy L Thompson [qdata, ..]: QFunctionOutputs| { 145ded9b81dSJeremy L Thompson // Build quadrature data 146ded9b81dSJeremy L Thompson match dim { 147ded9b81dSJeremy L Thompson 1 => qdata 148ded9b81dSJeremy L Thompson .iter_mut() 149ded9b81dSJeremy L Thompson .zip(jacobian.iter().zip(weights.iter())) 150ded9b81dSJeremy L Thompson .for_each(|(qdata, (j, weight))| *qdata = weight / j), 151ded9b81dSJeremy L Thompson 2 => { 152ded9b81dSJeremy L Thompson let q = qdata.len() / 3; 153ded9b81dSJeremy L Thompson for i in 0..q { 154ded9b81dSJeremy L Thompson let j11 = jacobian[i + q * 0]; 155ded9b81dSJeremy L Thompson let j21 = jacobian[i + q * 1]; 156ded9b81dSJeremy L Thompson let j12 = jacobian[i + q * 2]; 157ded9b81dSJeremy L Thompson let j22 = jacobian[i + q * 3]; 158ded9b81dSJeremy L Thompson let qw = weights[i] / (j11 * j22 - j21 * j12); 159ded9b81dSJeremy L Thompson qdata[i + q * 0] = qw * (j12 * j12 + j22 * j22); 160ded9b81dSJeremy L Thompson qdata[i + q * 1] = qw * (j11 * j11 + j21 * j21); 161ded9b81dSJeremy L Thompson qdata[i + q * 2] = -qw * (j11 * j12 + j21 * j22); 162ded9b81dSJeremy L Thompson } 163ded9b81dSJeremy L Thompson } 164ded9b81dSJeremy L Thompson 3 => { 165ded9b81dSJeremy L Thompson let q = qdata.len() / 6; 166ded9b81dSJeremy L Thompson for i in 0..q { 167ded9b81dSJeremy L Thompson let mut a = [0.0; 9]; 168ded9b81dSJeremy L Thompson for j in 0..3 { 169ded9b81dSJeremy L Thompson for k in 0..3 { 170ded9b81dSJeremy L Thompson a[k * 3 + j] = jacobian[i + q * ((j + 1) % 3 + 3 * ((k + 1) % 3))] 171ded9b81dSJeremy L Thompson * jacobian[i + q * ((j + 2) % 3 + 3 * ((k + 2) % 3))] 172ded9b81dSJeremy L Thompson - jacobian[i + q * ((j + 1) % 3 + 3 * ((k + 2) % 3))] 173ded9b81dSJeremy L Thompson * jacobian[i + q * ((j + 2) % 3 + 3 * ((k + 1) % 3))]; 174ded9b81dSJeremy L Thompson } 175ded9b81dSJeremy L Thompson } 176ded9b81dSJeremy L Thompson let qw = weights[i] 177ded9b81dSJeremy L Thompson / (jacobian[i + q * 0] * a[0 * 3 + 0] 178ded9b81dSJeremy L Thompson + jacobian[i + q * 1] * a[1 * 3 + 1] 179ded9b81dSJeremy L Thompson + jacobian[i + q * 2] * a[2 * 3 + 2]); 180ded9b81dSJeremy L Thompson qdata[i + q * 0] = qw 181ded9b81dSJeremy L Thompson * (a[0 * 3 + 0] * a[0 * 3 + 0] 182ded9b81dSJeremy L Thompson + a[0 * 3 + 1] * a[0 * 3 + 1] 183ded9b81dSJeremy L Thompson + a[0 * 3 + 2] * a[0 * 3 + 2]); 184ded9b81dSJeremy L Thompson qdata[i + q * 1] = qw 185ded9b81dSJeremy L Thompson * (a[1 * 3 + 0] * a[1 * 3 + 0] 186ded9b81dSJeremy L Thompson + a[1 * 3 + 1] * a[1 * 3 + 1] 187ded9b81dSJeremy L Thompson + a[1 * 3 + 2] * a[1 * 3 + 2]); 188ded9b81dSJeremy L Thompson qdata[i + q * 2] = qw 189ded9b81dSJeremy L Thompson * (a[2 * 3 + 0] * a[2 * 3 + 0] 190ded9b81dSJeremy L Thompson + a[2 * 3 + 1] * a[2 * 3 + 1] 191ded9b81dSJeremy L Thompson + a[2 * 3 + 2] * a[2 * 3 + 2]); 192ded9b81dSJeremy L Thompson qdata[i + q * 3] = qw 193ded9b81dSJeremy L Thompson * (a[1 * 3 + 0] * a[2 * 3 + 0] 194ded9b81dSJeremy L Thompson + a[1 * 3 + 1] * a[2 * 3 + 1] 195ded9b81dSJeremy L Thompson + a[1 * 3 + 2] * a[2 * 3 + 2]); 196ded9b81dSJeremy L Thompson qdata[i + q * 4] = qw 197ded9b81dSJeremy L Thompson * (a[0 * 3 + 0] * a[2 * 3 + 0] 198ded9b81dSJeremy L Thompson + a[0 * 3 + 1] * a[2 * 3 + 1] 199ded9b81dSJeremy L Thompson + a[0 * 3 + 2] * a[2 * 3 + 2]); 200ded9b81dSJeremy L Thompson qdata[i + q * 5] = qw 201ded9b81dSJeremy L Thompson * (a[0 * 3 + 0] * a[1 * 3 + 0] 202ded9b81dSJeremy L Thompson + a[0 * 3 + 1] * a[1 * 3 + 1] 203ded9b81dSJeremy L Thompson + a[0 * 3 + 2] * a[1 * 3 + 2]); 204ded9b81dSJeremy L Thompson } 205ded9b81dSJeremy L Thompson } 206ded9b81dSJeremy L Thompson _ => unreachable!(), 207ded9b81dSJeremy L Thompson }; 208ded9b81dSJeremy L Thompson 209ded9b81dSJeremy L Thompson // Return clean error code 210ded9b81dSJeremy L Thompson 0 211ded9b81dSJeremy L Thompson }; 212ded9b81dSJeremy L Thompson let qf_build_closure = ceed 213e15f9bd0SJeremy L Thompson .q_function_interior(1, Box::new(build_diff))? 214e15f9bd0SJeremy L Thompson .input("dx", ncomp_x * dim, EvalMode::Grad)? 215e15f9bd0SJeremy L Thompson .input("weights", 1, EvalMode::Weight)? 216e15f9bd0SJeremy L Thompson .output("qdata", dim * (dim + 1) / 2, EvalMode::None)?; 217ded9b81dSJeremy L Thompson // -- QFunction from gallery 218ded9b81dSJeremy L Thompson let qf_build_named = { 219ded9b81dSJeremy L Thompson let name = format!("Poisson{}DBuild", dim); 220e15f9bd0SJeremy L Thompson ceed.q_function_interior_by_name(&name)? 221ded9b81dSJeremy L Thompson }; 222ded9b81dSJeremy L Thompson // -- QFunction for use with Operator 223ded9b81dSJeremy L Thompson let qf_build = if gallery { 224ded9b81dSJeremy L Thompson QFunctionOpt::SomeQFunctionByName(&qf_build_named) 225ded9b81dSJeremy L Thompson } else { 226ded9b81dSJeremy L Thompson QFunctionOpt::SomeQFunction(&qf_build_closure) 227ded9b81dSJeremy L Thompson }; 228ded9b81dSJeremy L Thompson 229ded9b81dSJeremy L Thompson // Operator that build the quadrature data for the diff operator 230ded9b81dSJeremy L Thompson let op_build = ceed 231e15f9bd0SJeremy L Thompson .operator(qf_build, QFunctionOpt::None, QFunctionOpt::None)? 232e15f9bd0SJeremy L Thompson .field("dx", &restr_mesh, &basis_mesh, VectorOpt::Active)? 233ded9b81dSJeremy L Thompson .field( 234ded9b81dSJeremy L Thompson "weights", 235ded9b81dSJeremy L Thompson ElemRestrictionOpt::None, 236ded9b81dSJeremy L Thompson &basis_mesh, 237ded9b81dSJeremy L Thompson VectorOpt::None, 238e15f9bd0SJeremy L Thompson )? 239ded9b81dSJeremy L Thompson .field( 240ded9b81dSJeremy L Thompson "qdata", 241ded9b81dSJeremy L Thompson &restr_qdata, 242ded9b81dSJeremy L Thompson BasisOpt::Collocated, 243ded9b81dSJeremy L Thompson VectorOpt::Active, 244e15f9bd0SJeremy L Thompson )?; 245ded9b81dSJeremy L Thompson 246ded9b81dSJeremy L Thompson // Compute the quadrature data for the diff operator 247ded9b81dSJeremy L Thompson let elem_qpts = num_qpts.pow(dim as u32); 248ded9b81dSJeremy L Thompson let num_elem: usize = num_xyz.iter().take(dim).product(); 249e15f9bd0SJeremy L Thompson let mut qdata = ceed.vector(num_elem * elem_qpts * dim * (dim + 1) / 2)?; 250e15f9bd0SJeremy L Thompson op_build.apply(&mesh_coords, &mut qdata)?; 251ded9b81dSJeremy L Thompson 252ded9b81dSJeremy L Thompson // QFunction that applies the diff operator 253ded9b81dSJeremy L Thompson // -- QFunction from user closure 254ded9b81dSJeremy L Thompson let apply_diff = move |[ug, qdata, ..]: QFunctionInputs, [vg, ..]: QFunctionOutputs| { 255ded9b81dSJeremy L Thompson // Apply diffusion operator 256ded9b81dSJeremy L Thompson match dim { 257ded9b81dSJeremy L Thompson 1 => vg 258ded9b81dSJeremy L Thompson .iter_mut() 259ded9b81dSJeremy L Thompson .zip(ug.iter().zip(qdata.iter())) 260ded9b81dSJeremy L Thompson .for_each(|(vg, (ug, w))| *vg = ug * w), 261ded9b81dSJeremy L Thompson 2 => { 262ded9b81dSJeremy L Thompson let q = qdata.len() / 3; 263ded9b81dSJeremy L Thompson for i in 0..q { 264ded9b81dSJeremy L Thompson let du = [ug[i + q * 0], ug[i + q * 1]]; 265ded9b81dSJeremy L Thompson let dxdxdxdx_t = [ 266ded9b81dSJeremy L Thompson [qdata[i + 0 * q], qdata[i + 2 * q]], 267ded9b81dSJeremy L Thompson [qdata[i + 2 * q], qdata[i + 1 * q]], 268ded9b81dSJeremy L Thompson ]; 269ded9b81dSJeremy L Thompson for j in 0..2 { 270ded9b81dSJeremy L Thompson vg[i + j * q] = du[0] * dxdxdxdx_t[0][j] + du[1] * dxdxdxdx_t[1][j]; 271ded9b81dSJeremy L Thompson } 272ded9b81dSJeremy L Thompson } 273ded9b81dSJeremy L Thompson } 274ded9b81dSJeremy L Thompson 3 => { 275ded9b81dSJeremy L Thompson let q = qdata.len() / 6; 276ded9b81dSJeremy L Thompson for i in 0..q { 277ded9b81dSJeremy L Thompson let du = [ug[i + q * 0], ug[i + q * 1], ug[i + q * 2]]; 278ded9b81dSJeremy L Thompson let dxdxdxdx_t = [ 279ded9b81dSJeremy L Thompson [qdata[i + 0 * q], qdata[i + 5 * q], qdata[i + 4 * q]], 280ded9b81dSJeremy L Thompson [qdata[i + 5 * q], qdata[i + 1 * q], qdata[i + 3 * q]], 281ded9b81dSJeremy L Thompson [qdata[i + 4 * q], qdata[i + 3 * q], qdata[i + 2 * q]], 282ded9b81dSJeremy L Thompson ]; 283ded9b81dSJeremy L Thompson for j in 0..3 { 284ded9b81dSJeremy L Thompson vg[i + j * q] = du[0] * dxdxdxdx_t[0][j] 285ded9b81dSJeremy L Thompson + du[1] * dxdxdxdx_t[1][j] 286ded9b81dSJeremy L Thompson + du[2] * dxdxdxdx_t[2][j]; 287ded9b81dSJeremy L Thompson } 288ded9b81dSJeremy L Thompson } 289ded9b81dSJeremy L Thompson } 290ded9b81dSJeremy L Thompson _ => unreachable!(), 291ded9b81dSJeremy L Thompson }; 292ded9b81dSJeremy L Thompson 293ded9b81dSJeremy L Thompson // Return clean error code 294ded9b81dSJeremy L Thompson 0 295ded9b81dSJeremy L Thompson }; 296ded9b81dSJeremy L Thompson let qf_diff_closure = ceed 297e15f9bd0SJeremy L Thompson .q_function_interior(1, Box::new(apply_diff))? 298e15f9bd0SJeremy L Thompson .input("du", dim, EvalMode::Grad)? 299e15f9bd0SJeremy L Thompson .input("qdata", dim * (dim + 1) / 2, EvalMode::None)? 300e15f9bd0SJeremy L Thompson .output("dv", dim, EvalMode::Grad)?; 301ded9b81dSJeremy L Thompson // -- QFunction from gallery 302ded9b81dSJeremy L Thompson let qf_diff_named = { 303ded9b81dSJeremy L Thompson let name = format!("Poisson{}DApply", dim); 304e15f9bd0SJeremy L Thompson ceed.q_function_interior_by_name(&name)? 305ded9b81dSJeremy L Thompson }; 306ded9b81dSJeremy L Thompson // -- QFunction for use with Operator 307ded9b81dSJeremy L Thompson let qf_diff = if gallery { 308ded9b81dSJeremy L Thompson QFunctionOpt::SomeQFunctionByName(&qf_diff_named) 309ded9b81dSJeremy L Thompson } else { 310ded9b81dSJeremy L Thompson QFunctionOpt::SomeQFunction(&qf_diff_closure) 311ded9b81dSJeremy L Thompson }; 312ded9b81dSJeremy L Thompson 313ded9b81dSJeremy L Thompson // Diff Operator 314ded9b81dSJeremy L Thompson let op_diff = ceed 315e15f9bd0SJeremy L Thompson .operator(qf_diff, QFunctionOpt::None, QFunctionOpt::None)? 316e15f9bd0SJeremy L Thompson .field("du", &restr_solution, &basis_solution, VectorOpt::Active)? 317e15f9bd0SJeremy L Thompson .field("qdata", &restr_qdata, BasisOpt::Collocated, &qdata)? 318e15f9bd0SJeremy L Thompson .field("dv", &restr_solution, &basis_solution, VectorOpt::Active)?; 319ded9b81dSJeremy L Thompson 320ded9b81dSJeremy L Thompson // Solution vectors 321e15f9bd0SJeremy L Thompson let mut u = ceed.vector(solution_size)?; 322e15f9bd0SJeremy L Thompson let mut v = ceed.vector(solution_size)?; 323ded9b81dSJeremy L Thompson 324ded9b81dSJeremy L Thompson // Initialize u with sum of node coordinates 325e78171edSJeremy L Thompson let coords = mesh_coords.view()?; 326e78171edSJeremy L Thompson u.view_mut()?.iter_mut().enumerate().for_each(|(i, u)| { 327ded9b81dSJeremy L Thompson *u = (0..dim).map(|d| coords[i + d * solution_size]).sum(); 328ded9b81dSJeremy L Thompson }); 329ded9b81dSJeremy L Thompson 330ded9b81dSJeremy L Thompson // Apply the diff operator 331e15f9bd0SJeremy L Thompson op_diff.apply(&u, &mut v)?; 332ded9b81dSJeremy L Thompson 333ded9b81dSJeremy L Thompson // Compute the mesh surface area 334e78171edSJeremy L Thompson let area: Scalar = v.view()?.iter().map(|v| (*v).abs()).sum(); 335ded9b81dSJeremy L Thompson 336ded9b81dSJeremy L Thompson // Output results 337ded9b81dSJeremy L Thompson if !quiet { 338ded9b81dSJeremy L Thompson println!("Exact mesh surface area : {:.12}", exact_area); 339ded9b81dSJeremy L Thompson println!("Computed mesh surface_area : {:.12}", area); 340ded9b81dSJeremy L Thompson println!("Surface area error : {:.12e}", area - exact_area); 341ded9b81dSJeremy L Thompson } 342ded9b81dSJeremy L Thompson let tolerance = match dim { 34380a9ef05SNatalie Beams 1 => 10000.0 * libceed::EPSILON, 344ded9b81dSJeremy L Thompson _ => 1E-1, 345ded9b81dSJeremy L Thompson }; 346ded9b81dSJeremy L Thompson let error = (area - exact_area).abs(); 347ded9b81dSJeremy L Thompson if error > tolerance { 348ded9b81dSJeremy L Thompson println!("Volume error too large: {:.12e}", error); 3492ba8e59cSJeremy L Thompson return Err(libceed::Error { 350e15f9bd0SJeremy L Thompson message: format!( 351ded9b81dSJeremy L Thompson "Volume error too large - expected: {:.12e}, actual: {:.12e}", 352ded9b81dSJeremy L Thompson tolerance, error 353e15f9bd0SJeremy L Thompson ), 354e15f9bd0SJeremy L Thompson }); 355ded9b81dSJeremy L Thompson } 356ded9b81dSJeremy L Thompson Ok(()) 357ded9b81dSJeremy L Thompson } 358ded9b81dSJeremy L Thompson 359ded9b81dSJeremy L Thompson // ---------------------------------------------------------------------------- 360ded9b81dSJeremy L Thompson // Tests 361ded9b81dSJeremy L Thompson // ---------------------------------------------------------------------------- 362ded9b81dSJeremy L Thompson #[cfg(test)] 363ded9b81dSJeremy L Thompson mod tests { 364ded9b81dSJeremy L Thompson use super::*; 365ded9b81dSJeremy L Thompson 366ded9b81dSJeremy L Thompson #[test] 367ded9b81dSJeremy L Thompson fn example_2_1d() { 368ded9b81dSJeremy L Thompson let options = opt::Opt { 369ded9b81dSJeremy L Thompson ceed_spec: "/cpu/self/ref/serial".to_string(), 370ded9b81dSJeremy L Thompson dim: 1, 371ded9b81dSJeremy L Thompson mesh_degree: 4, 372ded9b81dSJeremy L Thompson solution_degree: 4, 373ded9b81dSJeremy L Thompson num_qpts: 6, 374ded9b81dSJeremy L Thompson problem_size_requested: -1, 375ded9b81dSJeremy L Thompson test: true, 376ded9b81dSJeremy L Thompson quiet: true, 377ded9b81dSJeremy L Thompson gallery: false, 378ded9b81dSJeremy L Thompson }; 379ded9b81dSJeremy L Thompson assert!(example_2(options).is_ok()); 380ded9b81dSJeremy L Thompson } 381ded9b81dSJeremy L Thompson 382ded9b81dSJeremy L Thompson #[test] 383ded9b81dSJeremy L Thompson fn example_2_2d() { 384ded9b81dSJeremy L Thompson let options = opt::Opt { 385ded9b81dSJeremy L Thompson ceed_spec: "/cpu/self/ref/serial".to_string(), 386ded9b81dSJeremy L Thompson dim: 2, 387ded9b81dSJeremy L Thompson mesh_degree: 4, 388ded9b81dSJeremy L Thompson solution_degree: 4, 389ded9b81dSJeremy L Thompson num_qpts: 6, 390ded9b81dSJeremy L Thompson problem_size_requested: -1, 391ded9b81dSJeremy L Thompson test: true, 392ded9b81dSJeremy L Thompson quiet: true, 393ded9b81dSJeremy L Thompson gallery: false, 394ded9b81dSJeremy L Thompson }; 395ded9b81dSJeremy L Thompson assert!(example_2(options).is_ok()); 396ded9b81dSJeremy L Thompson } 397ded9b81dSJeremy L Thompson 398ded9b81dSJeremy L Thompson #[test] 399ded9b81dSJeremy L Thompson fn example_2_3d() { 400ded9b81dSJeremy L Thompson let options = opt::Opt { 401ded9b81dSJeremy L Thompson ceed_spec: "/cpu/self/ref/serial".to_string(), 402ded9b81dSJeremy L Thompson dim: 3, 403ded9b81dSJeremy L Thompson mesh_degree: 4, 404ded9b81dSJeremy L Thompson solution_degree: 4, 405ded9b81dSJeremy L Thompson num_qpts: 6, 406ded9b81dSJeremy L Thompson problem_size_requested: -1, 407ded9b81dSJeremy L Thompson test: true, 408954a6033SJeremy L Thompson quiet: false, 409ded9b81dSJeremy L Thompson gallery: false, 410ded9b81dSJeremy L Thompson }; 411ded9b81dSJeremy L Thompson assert!(example_2(options).is_ok()); 412ded9b81dSJeremy L Thompson } 413ded9b81dSJeremy L Thompson 414ded9b81dSJeremy L Thompson #[test] 415ded9b81dSJeremy L Thompson fn example_2_1d_gallery() { 416ded9b81dSJeremy L Thompson let options = opt::Opt { 417ded9b81dSJeremy L Thompson ceed_spec: "/cpu/self/ref/serial".to_string(), 418ded9b81dSJeremy L Thompson dim: 1, 419ded9b81dSJeremy L Thompson mesh_degree: 4, 420ded9b81dSJeremy L Thompson solution_degree: 4, 421ded9b81dSJeremy L Thompson num_qpts: 6, 422ded9b81dSJeremy L Thompson problem_size_requested: -1, 423ded9b81dSJeremy L Thompson test: true, 424ded9b81dSJeremy L Thompson quiet: true, 425ded9b81dSJeremy L Thompson gallery: true, 426ded9b81dSJeremy L Thompson }; 427ded9b81dSJeremy L Thompson assert!(example_2(options).is_ok()); 428ded9b81dSJeremy L Thompson } 429ded9b81dSJeremy L Thompson 430ded9b81dSJeremy L Thompson #[test] 431ded9b81dSJeremy L Thompson fn example_2_2d_gallery() { 432ded9b81dSJeremy L Thompson let options = opt::Opt { 433ded9b81dSJeremy L Thompson ceed_spec: "/cpu/self/ref/serial".to_string(), 434ded9b81dSJeremy L Thompson dim: 2, 435ded9b81dSJeremy L Thompson mesh_degree: 4, 436ded9b81dSJeremy L Thompson solution_degree: 4, 437ded9b81dSJeremy L Thompson num_qpts: 6, 438ded9b81dSJeremy L Thompson problem_size_requested: -1, 439ded9b81dSJeremy L Thompson test: true, 440ded9b81dSJeremy L Thompson quiet: true, 441ded9b81dSJeremy L Thompson gallery: true, 442ded9b81dSJeremy L Thompson }; 443ded9b81dSJeremy L Thompson assert!(example_2(options).is_ok()); 444ded9b81dSJeremy L Thompson } 445ded9b81dSJeremy L Thompson 446ded9b81dSJeremy L Thompson #[test] 447ded9b81dSJeremy L Thompson fn example_2_3d_gallery() { 448ded9b81dSJeremy L Thompson let options = opt::Opt { 449ded9b81dSJeremy L Thompson ceed_spec: "/cpu/self/ref/serial".to_string(), 450ded9b81dSJeremy L Thompson dim: 3, 451ded9b81dSJeremy L Thompson mesh_degree: 4, 452ded9b81dSJeremy L Thompson solution_degree: 4, 453ded9b81dSJeremy L Thompson num_qpts: 6, 454ded9b81dSJeremy L Thompson problem_size_requested: -1, 455ded9b81dSJeremy L Thompson test: true, 456ded9b81dSJeremy L Thompson quiet: true, 457ded9b81dSJeremy L Thompson gallery: true, 458ded9b81dSJeremy L Thompson }; 459ded9b81dSJeremy L Thompson assert!(example_2(options).is_ok()); 460ded9b81dSJeremy L Thompson } 461ded9b81dSJeremy L Thompson } 462ded9b81dSJeremy L Thompson 463ded9b81dSJeremy L Thompson // ---------------------------------------------------------------------------- 464