Cabana 0.8.0-dev
 
Loading...
Searching...
No Matches
Cabana_Grid_GlobalMesh.hpp
Go to the documentation of this file.
1/****************************************************************************
2 * Copyright (c) 2018-2023 by the Cabana authors *
3 * All rights reserved. *
4 * *
5 * This file is part of the Cabana library. Cabana is distributed under a *
6 * BSD 3-clause license. For the licensing terms see the LICENSE file in *
7 * the top-level directory. *
8 * *
9 * SPDX-License-Identifier: BSD-3-Clause *
10 ****************************************************************************/
11
16#ifndef CABANA_GRID_GLOBALMESH_HPP
17#define CABANA_GRID_GLOBALMESH_HPP
18
19#include <Cabana_Grid_Types.hpp>
20
21#include <array>
22#include <cmath>
23#include <limits>
24#include <memory>
25#include <stdexcept>
26#include <type_traits>
27#include <vector>
28
29namespace Cabana
30{
31namespace Grid
32{
33//---------------------------------------------------------------------------//
34// Forward declaration of global mesh.
35template <class MeshType>
36class GlobalMesh;
37
38//---------------------------------------------------------------------------//
47template <class MeshType>
49{
50 public:
52 using mesh_type = MeshType;
53
55 using scalar_type = typename mesh_type::scalar_type;
56
58 static constexpr std::size_t num_space_dim = mesh_type::num_space_dim;
59
62 const std::array<scalar_type, num_space_dim>& global_low_corner,
63 const std::array<scalar_type, num_space_dim>& global_high_corner,
64 const scalar_type cell_size )
65 : _global_low_corner( global_low_corner )
66 , _global_high_corner( global_high_corner )
67 {
68 // Check that the domain is evenly divisible by the cell size in each
69 // dimension within round-off error.
70 for ( std::size_t d = 0; d < num_space_dim; ++d )
71 {
72 _cell_size[d] = cell_size;
73 scalar_type ext = globalNumCell( d ) * _cell_size[d];
74 if ( std::abs( ext - extent( d ) ) >
75 scalar_type( 100.0 ) *
76 std::numeric_limits<scalar_type>::epsilon() )
77 throw std::logic_error(
78 "Cabana::Grid::GlobalMesh: Extent not evenly divisible by "
79 "uniform cell size" );
80 }
81 }
82
85 const std::array<scalar_type, num_space_dim>& global_low_corner,
86 const std::array<scalar_type, num_space_dim>& global_high_corner,
87 const std::array<scalar_type, num_space_dim>& cell_size )
88 : _global_low_corner( global_low_corner )
89 , _global_high_corner( global_high_corner )
90 , _cell_size( cell_size )
91 {
92 // Check that the domain is evenly divisible by the cell size in each
93 // dimension within round-off error.
94 for ( std::size_t d = 0; d < num_space_dim; ++d )
95 {
96 scalar_type ext = globalNumCell( d ) * _cell_size[d];
97 if ( std::abs( ext - extent( d ) ) >
98 scalar_type( 100.0 ) *
99 std::numeric_limits<scalar_type>::epsilon() )
100 throw std::logic_error(
101 "Cabana::Grid::GlocalMesh: Extent not evenly divisible by "
102 "uniform cell size" );
103 }
104 }
105
108 const std::array<scalar_type, num_space_dim>& global_low_corner,
109 const std::array<scalar_type, num_space_dim>& global_high_corner,
110 const std::array<int, num_space_dim>& global_num_cell )
111 : _global_low_corner( global_low_corner )
112 , _global_high_corner( global_high_corner )
113 {
114 // Compute the cell size in each dimension.
115 for ( std::size_t d = 0; d < num_space_dim; ++d )
116 _cell_size[d] = ( _global_high_corner[d] - _global_low_corner[d] ) /
117 global_num_cell[d];
118
119 // Check that the domain is evenly divisible by the cell size in each
120 // dimension within round-off error and that we got the expected
121 // number of cells.
122 for ( std::size_t d = 0; d < num_space_dim; ++d )
123 {
124 scalar_type ext = globalNumCell( d ) * _cell_size[d];
125 if ( std::abs( ext - extent( d ) ) >
126 scalar_type( 100.0 ) *
127 std::numeric_limits<scalar_type>::epsilon() )
128 throw std::logic_error(
129 "Cabana::Grid::GlocalMesh: Extent not evenly divisible by "
130 "uniform cell size" );
131 if ( globalNumCell( d ) != global_num_cell[d] )
132 throw std::logic_error( "Cabana::Grid::GlocalMesh: Global "
133 "number of cells mismatch" );
134 }
135 }
136
137 // GLOBAL MESH INTERFACE
138
141 scalar_type lowCorner( const std::size_t dim ) const
142 {
143 return _global_low_corner[dim];
144 }
145
148 scalar_type highCorner( const std::size_t dim ) const
149 {
150 return _global_high_corner[dim];
151 }
152
155 scalar_type extent( const std::size_t dim ) const
156 {
157 return highCorner( dim ) - lowCorner( dim );
158 }
159
162 int globalNumCell( const std::size_t dim ) const
163 {
164 return std::rint( extent( dim ) / _cell_size[dim] );
165 }
166
167 // UNIFORM MESH SPECIFIC
168
171 scalar_type cellSize( const std::size_t dim ) const
172 {
173 return _cell_size[dim];
174 }
175
176 private:
177 std::array<scalar_type, num_space_dim> _global_low_corner;
178 std::array<scalar_type, num_space_dim> _global_high_corner;
179 std::array<scalar_type, num_space_dim> _cell_size;
180};
181
192template <class Scalar, std::size_t NumSpaceDim>
193std::shared_ptr<GlobalMesh<UniformMesh<Scalar, NumSpaceDim>>>
195 const std::array<Scalar, NumSpaceDim>& global_low_corner,
196 const std::array<Scalar, NumSpaceDim>& global_high_corner,
197 const Scalar cell_size )
198{
199 return std::make_shared<GlobalMesh<UniformMesh<Scalar, NumSpaceDim>>>(
200 global_low_corner, global_high_corner, cell_size );
201}
202
213template <class Scalar, std::size_t NumSpaceDim>
214std::shared_ptr<GlobalMesh<UniformMesh<Scalar, NumSpaceDim>>>
216 const std::array<Scalar, NumSpaceDim>& global_low_corner,
217 const std::array<Scalar, NumSpaceDim>& global_high_corner,
218 const std::array<Scalar, NumSpaceDim>& cell_size )
219{
220 return std::make_shared<GlobalMesh<UniformMesh<Scalar, NumSpaceDim>>>(
221 global_low_corner, global_high_corner, cell_size );
222}
223
234template <class Scalar, std::size_t NumSpaceDim>
235std::shared_ptr<GlobalMesh<UniformMesh<Scalar, NumSpaceDim>>>
237 const std::array<Scalar, NumSpaceDim>& global_low_corner,
238 const std::array<Scalar, NumSpaceDim>& global_high_corner,
239 const std::array<int, NumSpaceDim>& global_num_cell )
240{
241 return std::make_shared<GlobalMesh<UniformMesh<Scalar, NumSpaceDim>>>(
242 global_low_corner, global_high_corner, global_num_cell );
243}
244
245//---------------------------------------------------------------------------//
256template <class Scalar, std::size_t NumSpaceDim>
257std::shared_ptr<GlobalMesh<SparseMesh<Scalar, NumSpaceDim>>>
259 const std::array<Scalar, NumSpaceDim>& global_low_corner,
260 const std::array<Scalar, NumSpaceDim>& global_high_corner,
261 const Scalar cell_size )
262{
263 return std::make_shared<GlobalMesh<SparseMesh<Scalar, NumSpaceDim>>>(
264 global_low_corner, global_high_corner, cell_size );
265}
266
277template <class Scalar, std::size_t NumSpaceDim>
278std::shared_ptr<GlobalMesh<SparseMesh<Scalar, NumSpaceDim>>>
280 const std::array<Scalar, NumSpaceDim>& global_low_corner,
281 const std::array<Scalar, NumSpaceDim>& global_high_corner,
282 const std::array<Scalar, NumSpaceDim>& cell_size )
283{
284 return std::make_shared<GlobalMesh<SparseMesh<Scalar, NumSpaceDim>>>(
285 global_low_corner, global_high_corner, cell_size );
286}
287
298template <class Scalar, std::size_t NumSpaceDim>
299std::shared_ptr<GlobalMesh<SparseMesh<Scalar, NumSpaceDim>>>
301 const std::array<Scalar, NumSpaceDim>& global_low_corner,
302 const std::array<Scalar, NumSpaceDim>& global_high_corner,
303 const std::array<int, NumSpaceDim>& global_num_cell )
304{
305 return std::make_shared<GlobalMesh<SparseMesh<Scalar, NumSpaceDim>>>(
306 global_low_corner, global_high_corner, global_num_cell );
307}
308
309//---------------------------------------------------------------------------//
319template <class Scalar>
320class GlobalMesh<NonUniformMesh<Scalar, 3>>
321{
322 public:
325
327 using scalar_type = Scalar;
328
330 static constexpr std::size_t num_space_dim = 3;
331
334 GlobalMesh( const std::vector<Scalar>& i_edges,
335 const std::vector<Scalar>& j_edges,
336 const std::vector<Scalar>& k_edges )
337 : _edges( { i_edges, j_edges, k_edges } )
338 {
339 }
340
341 // GLOBAL MESH INTERFACE
342
345 Scalar lowCorner( const std::size_t dim ) const
346 {
347 return _edges[dim].front();
348 }
349
352 Scalar highCorner( const std::size_t dim ) const
353 {
354 return _edges[dim].back();
355 }
356
359 Scalar extent( const std::size_t dim ) const
360 {
361 return highCorner( dim ) - lowCorner( dim );
362 }
363
366 int globalNumCell( const std::size_t dim ) const
367 {
368 return _edges[dim].size() - 1;
369 }
370
371 // NON-UNIFORM MESH SPECIFIC
372
375 const std::vector<Scalar>& nonUniformEdge( const std::size_t dim ) const
376 {
377 return _edges[dim];
378 }
379
380 private:
381 std::array<std::vector<Scalar>, 3> _edges;
382};
383
390template <class Scalar>
391std::shared_ptr<GlobalMesh<NonUniformMesh<Scalar, 3>>>
392createNonUniformGlobalMesh( const std::vector<Scalar>& i_edges,
393 const std::vector<Scalar>& j_edges,
394 const std::vector<Scalar>& k_edges )
395{
396 return std::make_shared<GlobalMesh<NonUniformMesh<Scalar, 3>>>(
397 i_edges, j_edges, k_edges );
398}
399
400//---------------------------------------------------------------------------//
410template <class Scalar>
411class GlobalMesh<NonUniformMesh<Scalar, 2>>
412{
413 public:
416
418 using scalar_type = Scalar;
419
421 static constexpr std::size_t num_space_dim = 2;
422
425 GlobalMesh( const std::vector<Scalar>& i_edges,
426 const std::vector<Scalar>& j_edges )
427 : _edges( { i_edges, j_edges } )
428 {
429 }
430
431 // GLOBAL MESH INTERFACE
432
435 Scalar lowCorner( const std::size_t dim ) const
436 {
437 return _edges[dim].front();
438 }
439
442 Scalar highCorner( const std::size_t dim ) const
443 {
444 return _edges[dim].back();
445 }
446
449 Scalar extent( const std::size_t dim ) const
450 {
451 return highCorner( dim ) - lowCorner( dim );
452 }
453
456 int globalNumCell( const std::size_t dim ) const
457 {
458 return _edges[dim].size() - 1;
459 }
460
461 // NON-UNIFORM MESH SPECIFIC
462
465 const std::vector<Scalar>& nonUniformEdge( const std::size_t dim ) const
466 {
467 return _edges[dim];
468 }
469
470 private:
471 std::array<std::vector<Scalar>, 2> _edges;
472};
473
480template <class Scalar>
481std::shared_ptr<GlobalMesh<NonUniformMesh<Scalar, 2>>>
482createNonUniformGlobalMesh( const std::vector<Scalar>& i_edges,
483 const std::vector<Scalar>& j_edges )
484{
485 return std::make_shared<GlobalMesh<NonUniformMesh<Scalar, 2>>>( i_edges,
486 j_edges );
487}
488
489//---------------------------------------------------------------------------//
490
491} // namespace Grid
492} // namespace Cabana
493
494#endif // end CABANA_GRID_GLOBALMESH_HPP
std::shared_ptr< GlobalMesh< SparseMesh< Scalar, NumSpaceDim > > > createSparseGlobalMesh(const std::array< Scalar, NumSpaceDim > &global_low_corner, const std::array< Scalar, NumSpaceDim > &global_high_corner, const Scalar cell_size)
Create sparse mesh with uniform cell size.
Definition Cabana_Grid_GlobalMesh.hpp:258
std::shared_ptr< GlobalMesh< NonUniformMesh< Scalar, 3 > > > createNonUniformGlobalMesh(const std::vector< Scalar > &i_edges, const std::vector< Scalar > &j_edges, const std::vector< Scalar > &k_edges)
Create a non-uniform 3D mesh.
Definition Cabana_Grid_GlobalMesh.hpp:392
std::shared_ptr< GlobalMesh< UniformMesh< Scalar, NumSpaceDim > > > createUniformGlobalMesh(const std::array< Scalar, NumSpaceDim > &global_low_corner, const std::array< Scalar, NumSpaceDim > &global_high_corner, const Scalar cell_size)
Create uniform mesh with uniform cell size.
Definition Cabana_Grid_GlobalMesh.hpp:194
Grid type tags.
NonUniformMesh< Scalar, 2 > mesh_type
Mesh type.
Definition Cabana_Grid_GlobalMesh.hpp:415
Scalar extent(const std::size_t dim) const
Get the extent of a given dimension.
Definition Cabana_Grid_GlobalMesh.hpp:449
Scalar lowCorner(const std::size_t dim) const
Get the global low corner of the mesh.
Definition Cabana_Grid_GlobalMesh.hpp:435
static constexpr std::size_t num_space_dim
Spatial dimension.
Definition Cabana_Grid_GlobalMesh.hpp:421
int globalNumCell(const std::size_t dim) const
Get the global number of cells in a given dimension.
Definition Cabana_Grid_GlobalMesh.hpp:456
Scalar highCorner(const std::size_t dim) const
Get the global high corner of the mesh.
Definition Cabana_Grid_GlobalMesh.hpp:442
Scalar scalar_type
Scalar type.
Definition Cabana_Grid_GlobalMesh.hpp:418
const std::vector< Scalar > & nonUniformEdge(const std::size_t dim) const
Get the edge array in a given dimension.
Definition Cabana_Grid_GlobalMesh.hpp:465
GlobalMesh(const std::vector< Scalar > &i_edges, const std::vector< Scalar > &j_edges)
2D constructor.
Definition Cabana_Grid_GlobalMesh.hpp:425
static constexpr std::size_t num_space_dim
Spatial dimension.
Definition Cabana_Grid_GlobalMesh.hpp:330
Scalar lowCorner(const std::size_t dim) const
Get the global low corner of the mesh.
Definition Cabana_Grid_GlobalMesh.hpp:345
NonUniformMesh< Scalar, 3 > mesh_type
Mesh type.
Definition Cabana_Grid_GlobalMesh.hpp:324
Scalar scalar_type
Scalar type.
Definition Cabana_Grid_GlobalMesh.hpp:327
int globalNumCell(const std::size_t dim) const
Get the global number of cells in a given dimension.
Definition Cabana_Grid_GlobalMesh.hpp:366
Scalar highCorner(const std::size_t dim) const
Get the global high corner of the mesh.
Definition Cabana_Grid_GlobalMesh.hpp:352
Scalar extent(const std::size_t dim) const
Get the extent of a given dimension.
Definition Cabana_Grid_GlobalMesh.hpp:359
GlobalMesh(const std::vector< Scalar > &i_edges, const std::vector< Scalar > &j_edges, const std::vector< Scalar > &k_edges)
3D constructor.
Definition Cabana_Grid_GlobalMesh.hpp:334
const std::vector< Scalar > & nonUniformEdge(const std::size_t dim) const
Get the edge array in a given dimension.
Definition Cabana_Grid_GlobalMesh.hpp:375
Global mesh partial specialization for uniform mesh.
Definition Cabana_Grid_GlobalMesh.hpp:49
scalar_type lowCorner(const std::size_t dim) const
Get the global low corner of the mesh.
Definition Cabana_Grid_GlobalMesh.hpp:141
static constexpr std::size_t num_space_dim
Spatial dimension.
Definition Cabana_Grid_GlobalMesh.hpp:58
GlobalMesh(const std::array< scalar_type, num_space_dim > &global_low_corner, const std::array< scalar_type, num_space_dim > &global_high_corner, const std::array< int, num_space_dim > &global_num_cell)
Number of global cells constructor.
Definition Cabana_Grid_GlobalMesh.hpp:107
scalar_type extent(const std::size_t dim) const
Get the extent of a given dimension.
Definition Cabana_Grid_GlobalMesh.hpp:155
scalar_type highCorner(const std::size_t dim) const
Get the global high corner of the mesh.
Definition Cabana_Grid_GlobalMesh.hpp:148
MeshType mesh_type
Mesh type.
Definition Cabana_Grid_GlobalMesh.hpp:52
GlobalMesh(const std::array< scalar_type, num_space_dim > &global_low_corner, const std::array< scalar_type, num_space_dim > &global_high_corner, const scalar_type cell_size)
Cell size constructor where all cell dimensions are the same.
Definition Cabana_Grid_GlobalMesh.hpp:61
scalar_type cellSize(const std::size_t dim) const
Get the uniform cell size in a given dimension.
Definition Cabana_Grid_GlobalMesh.hpp:171
GlobalMesh(const std::array< scalar_type, num_space_dim > &global_low_corner, const std::array< scalar_type, num_space_dim > &global_high_corner, const std::array< scalar_type, num_space_dim > &cell_size)
Cell size constructor - each cell dimension can be different.
Definition Cabana_Grid_GlobalMesh.hpp:84
typename mesh_type::scalar_type scalar_type
Scalar type.
Definition Cabana_Grid_GlobalMesh.hpp:55
int globalNumCell(const std::size_t dim) const
Get the global number of cells in a given dimension.
Definition Cabana_Grid_GlobalMesh.hpp:162
Core: particle data structures and algorithms.
Definition Cabana_AoSoA.hpp:36
Non-uniform mesh tag.
Definition Cabana_Grid_Types.hpp:240