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
7ded9b81dSJeremy L Thompson
8eb07d68fSJeremy L Thompson use libceed::{Ceed, ElemRestriction, Vector};
978c2cefaSJeremy L Thompson use std::convert::TryInto;
10ded9b81dSJeremy L Thompson
11ded9b81dSJeremy L Thompson // ----------------------------------------------------------------------------
12ded9b81dSJeremy L Thompson // Determine problem size in each dimension from size and dimenison
13ded9b81dSJeremy L Thompson // ----------------------------------------------------------------------------
cartesian_mesh_size(dim: usize, solution_degree: usize, problem_size: i64) -> [usize; 3]14ded9b81dSJeremy L Thompson pub fn cartesian_mesh_size(dim: usize, solution_degree: usize, problem_size: i64) -> [usize; 3] {
15ded9b81dSJeremy L Thompson // Use the approximate formula:
16ded9b81dSJeremy L Thompson // prob_size ~ num_elem * degree^dim
17ded9b81dSJeremy L Thompson let mut num_elem = problem_size / solution_degree.pow(dim as u32) as i64;
18ded9b81dSJeremy L Thompson let mut s = 0; // find s: num_elem / 2 < 2^s <= num_elem
19ded9b81dSJeremy L Thompson while num_elem > 1 {
20ded9b81dSJeremy L Thompson num_elem /= 2;
21ded9b81dSJeremy L Thompson s += 1;
22ded9b81dSJeremy L Thompson }
23ded9b81dSJeremy L Thompson
24ded9b81dSJeremy L Thompson // Size per dimension
25ded9b81dSJeremy L Thompson let mut r = s % dim;
2678c2cefaSJeremy L Thompson let xyz: [usize; 3] = (0..3)
2778c2cefaSJeremy L Thompson .map(|_| -> usize {
28ded9b81dSJeremy L Thompson let mut sd = s / dim;
29ded9b81dSJeremy L Thompson if r > 0 {
30ded9b81dSJeremy L Thompson sd += 1;
31ded9b81dSJeremy L Thompson r -= 1;
32ded9b81dSJeremy L Thompson }
3378c2cefaSJeremy L Thompson 1 << sd
3478c2cefaSJeremy L Thompson })
3578c2cefaSJeremy L Thompson .collect::<Vec<usize>>()
3678c2cefaSJeremy L Thompson .try_into()
3778c2cefaSJeremy L Thompson .unwrap();
3878c2cefaSJeremy L Thompson xyz
39ded9b81dSJeremy L Thompson }
40ded9b81dSJeremy L Thompson
41ded9b81dSJeremy L Thompson // ----------------------------------------------------------------------------
42ded9b81dSJeremy L Thompson // Build element restriction objects for the mesh
43ded9b81dSJeremy L Thompson // ----------------------------------------------------------------------------
build_cartesian_restriction( ceed: &Ceed, dim: usize, num_xyz: [usize; 3], degree: usize, num_comp: usize, num_qpts: usize, ) -> libceed::Result<(ElemRestriction, ElemRestriction)>44ded9b81dSJeremy L Thompson pub fn build_cartesian_restriction(
45ded9b81dSJeremy L Thompson ceed: &Ceed,
46ded9b81dSJeremy L Thompson dim: usize,
47ded9b81dSJeremy L Thompson num_xyz: [usize; 3],
48ded9b81dSJeremy L Thompson degree: usize,
49ded9b81dSJeremy L Thompson num_comp: usize,
50ded9b81dSJeremy L Thompson num_qpts: usize,
5189d15d5fSJeremy L Thompson ) -> libceed::Result<(ElemRestriction, ElemRestriction)> {
52ded9b81dSJeremy L Thompson let p = degree + 1;
53ded9b81dSJeremy L Thompson let num_nodes = p.pow(dim as u32); // number of nodes per element
54ded9b81dSJeremy L Thompson let elem_qpts = num_qpts.pow(dim as u32); // number of quadrature pts per element
55ded9b81dSJeremy L Thompson
56ded9b81dSJeremy L Thompson // Problem dimensions
57ded9b81dSJeremy L Thompson let mut num_d = [0; 3];
58ded9b81dSJeremy L Thompson let mut num_elem = 1;
59ded9b81dSJeremy L Thompson let mut scalar_size = 1;
60ded9b81dSJeremy L Thompson for d in 0..dim {
61ded9b81dSJeremy L Thompson num_elem *= num_xyz[d];
62ded9b81dSJeremy L Thompson num_d[d] = num_xyz[d] * (p - 1) + 1;
63ded9b81dSJeremy L Thompson scalar_size *= num_d[d];
64ded9b81dSJeremy L Thompson }
65ded9b81dSJeremy L Thompson
66ded9b81dSJeremy L Thompson // elem: 0 1 n-1
67ded9b81dSJeremy L Thompson // |---*-...-*---|---*-...-*---|- ... -|--...--|
68ded9b81dSJeremy L Thompson // nodes: 0 1 p-1 p p+1 2*p n*p
69ded9b81dSJeremy L Thompson let mut elem_nodes = vec![0; num_elem * num_nodes];
70ded9b81dSJeremy L Thompson for e in 0..num_elem {
71ded9b81dSJeremy L Thompson let mut e_xyz = [1; 3];
72ded9b81dSJeremy L Thompson let mut re = e;
73ded9b81dSJeremy L Thompson for d in 0..dim {
74ded9b81dSJeremy L Thompson e_xyz[d] = re % num_xyz[d];
75ded9b81dSJeremy L Thompson re /= num_xyz[d];
76ded9b81dSJeremy L Thompson }
77ded9b81dSJeremy L Thompson let loc_offset = e * num_nodes;
78ded9b81dSJeremy L Thompson for loc_nodes in 0..num_nodes {
79ded9b81dSJeremy L Thompson let mut global_nodes = 0;
80ded9b81dSJeremy L Thompson let mut global_nodes_stride = 1;
81ded9b81dSJeremy L Thompson let mut r_nodes = loc_nodes;
82ded9b81dSJeremy L Thompson for d in 0..dim {
83ded9b81dSJeremy L Thompson global_nodes += (e_xyz[d] * (p - 1) + r_nodes % p) * global_nodes_stride;
84ded9b81dSJeremy L Thompson global_nodes_stride *= num_d[d];
85ded9b81dSJeremy L Thompson r_nodes /= p;
86ded9b81dSJeremy L Thompson }
87ded9b81dSJeremy L Thompson elem_nodes[loc_offset + loc_nodes] = global_nodes as i32;
88ded9b81dSJeremy L Thompson }
89ded9b81dSJeremy L Thompson }
90ded9b81dSJeremy L Thompson
91ded9b81dSJeremy L Thompson // Mesh/solution data restriction
92edb2538eSJeremy L Thompson let rstr = ceed.elem_restriction(
93ded9b81dSJeremy L Thompson num_elem,
94ded9b81dSJeremy L Thompson num_nodes,
95ded9b81dSJeremy L Thompson num_comp,
96ded9b81dSJeremy L Thompson scalar_size,
97ded9b81dSJeremy L Thompson num_comp * scalar_size,
98eb07d68fSJeremy L Thompson libceed::MemType::Host,
99ded9b81dSJeremy L Thompson &elem_nodes,
100e15f9bd0SJeremy L Thompson )?;
101ded9b81dSJeremy L Thompson
10278c2cefaSJeremy L Thompson // Quadrature data restriction
103edb2538eSJeremy L Thompson let rstr_qdata = ceed.strided_elem_restriction(
104ded9b81dSJeremy L Thompson num_elem,
105ded9b81dSJeremy L Thompson elem_qpts,
106ded9b81dSJeremy L Thompson num_comp,
107ded9b81dSJeremy L Thompson num_comp * elem_qpts * num_elem,
108eb07d68fSJeremy L Thompson libceed::CEED_STRIDES_BACKEND,
109e15f9bd0SJeremy L Thompson )?;
110edb2538eSJeremy L Thompson Ok((rstr, rstr_qdata))
111ded9b81dSJeremy L Thompson }
112ded9b81dSJeremy L Thompson
113ded9b81dSJeremy L Thompson // ----------------------------------------------------------------------------
114ded9b81dSJeremy L Thompson // Set mesh coordinates
115ded9b81dSJeremy L Thompson // ----------------------------------------------------------------------------
cartesian_mesh_coords( ceed: &Ceed, dim: usize, num_xyz: [usize; 3], mesh_degree: usize, mesh_size: usize, ) -> libceed::Result<Vector>116ded9b81dSJeremy L Thompson pub fn cartesian_mesh_coords(
117ded9b81dSJeremy L Thompson ceed: &Ceed,
118ded9b81dSJeremy L Thompson dim: usize,
119ded9b81dSJeremy L Thompson num_xyz: [usize; 3],
120ded9b81dSJeremy L Thompson mesh_degree: usize,
121ded9b81dSJeremy L Thompson mesh_size: usize,
12289d15d5fSJeremy L Thompson ) -> libceed::Result<Vector> {
123ded9b81dSJeremy L Thompson let p = mesh_degree + 1;
124ded9b81dSJeremy L Thompson let mut num_d = [0; 3];
125ded9b81dSJeremy L Thompson let mut scalar_size = 1;
126ded9b81dSJeremy L Thompson for d in 0..dim {
127ded9b81dSJeremy L Thompson num_d[d] = num_xyz[d] * (p - 1) + 1;
128ded9b81dSJeremy L Thompson scalar_size *= num_d[d];
129ded9b81dSJeremy L Thompson }
130ded9b81dSJeremy L Thompson
131ded9b81dSJeremy L Thompson // Lobatto points
132eb07d68fSJeremy L Thompson let lobatto_basis =
133eb07d68fSJeremy L Thompson ceed.basis_tensor_H1_Lagrange(1, 1, 2, p, libceed::QuadMode::GaussLobatto)?;
134e15f9bd0SJeremy L Thompson let nodes_corners = ceed.vector_from_slice(&[0.0, 1.0])?;
135e15f9bd0SJeremy L Thompson let mut nodes_full = ceed.vector(p)?;
136ded9b81dSJeremy L Thompson lobatto_basis.apply(
137ded9b81dSJeremy L Thompson 1,
138eb07d68fSJeremy L Thompson libceed::TransposeMode::NoTranspose,
139eb07d68fSJeremy L Thompson libceed::EvalMode::Interp,
140ded9b81dSJeremy L Thompson &nodes_corners,
141ded9b81dSJeremy L Thompson &mut nodes_full,
142e15f9bd0SJeremy L Thompson )?;
143ded9b81dSJeremy L Thompson
144ded9b81dSJeremy L Thompson // Coordinates for mesh
145e15f9bd0SJeremy L Thompson let mut mesh_coords = ceed.vector(mesh_size)?;
1469c774eddSJeremy L Thompson mesh_coords.set_value(0.0)?;
147ded9b81dSJeremy L Thompson {
148e78171edSJeremy L Thompson let mut coords = mesh_coords.view_mut()?;
149e78171edSJeremy L Thompson let nodes = nodes_full.view()?;
150ded9b81dSJeremy L Thompson for gs_nodes in 0..scalar_size {
151ded9b81dSJeremy L Thompson let mut r_nodes = gs_nodes;
152ded9b81dSJeremy L Thompson for d in 0..dim {
153ded9b81dSJeremy L Thompson let d_1d = r_nodes % num_d[d];
154eb07d68fSJeremy L Thompson coords[gs_nodes + scalar_size * d] = ((d_1d / (p - 1)) as libceed::Scalar
155eb07d68fSJeremy L Thompson + nodes[d_1d % (p - 1)])
156eb07d68fSJeremy L Thompson / num_xyz[d] as libceed::Scalar;
157ded9b81dSJeremy L Thompson r_nodes /= num_d[d];
158ded9b81dSJeremy L Thompson }
159ded9b81dSJeremy L Thompson }
160ded9b81dSJeremy L Thompson }
161e15f9bd0SJeremy L Thompson Ok(mesh_coords)
162ded9b81dSJeremy L Thompson }
163ded9b81dSJeremy L Thompson
164ded9b81dSJeremy L Thompson // ----------------------------------------------------------------------------
165