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