Loading [MathJax]/extensions/TeX/AMSmath.js
ESPResSo
Extensible Simulation Package for Research on Soft Matter Systems
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages Concepts
boundary.hpp
Go to the documentation of this file.
1/*
2 * Copyright (C) 2020-2023 The ESPResSo project
3 *
4 * This file is part of ESPResSo.
5 *
6 * ESPResSo is free software: you can redistribute it and/or modify
7 * it under the terms of the GNU General Public License as published by
8 * the Free Software Foundation, either version 3 of the License, or
9 * (at your option) any later version.
10 *
11 * ESPResSo is distributed in the hope that it will be useful,
12 * but WITHOUT ANY WARRANTY; without even the implied warranty of
13 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
14 * GNU General Public License for more details.
15 *
16 * You should have received a copy of the GNU General Public License
17 * along with this program. If not, see <http://www.gnu.org/licenses/>.
18 */
19
20#pragma once
21
22#include "types_conversion.hpp"
23
25
26#include <utils/Vector.hpp>
27
28#include <cassert>
29#include <cstddef>
30#include <iterator>
31#include <vector>
32
33namespace walberla {
34
35inline std::vector<Utils::Vector3d>
36fill_3D_vector_array(std::vector<double> const &vec_flat,
37 Utils::Vector3i const &grid_size) {
38 auto const n_grid_points =
39 static_cast<std::size_t>(Utils::product(grid_size));
40 assert(vec_flat.size() == 3u * n_grid_points or vec_flat.size() == 3u);
41 std::vector<Utils::Vector3d> output_vector;
42 output_vector.reserve(3u * n_grid_points);
43
44 auto const vec_begin = std::begin(vec_flat);
45 auto const vec_end = std::end(vec_flat);
46 if (vec_flat.size() == 3u) {
47 auto const uniform_vector = Utils::Vector3d(vec_begin, vec_end);
48 output_vector.assign(n_grid_points, uniform_vector);
49 } else {
50 output_vector.reserve(n_grid_points);
51 for (auto it = vec_begin; it < vec_end; it += 3u) {
52 output_vector.emplace_back(it, it + 3u);
53 }
54 }
55
56 return output_vector;
57}
58
59inline std::vector<double>
60fill_3D_scalar_array(std::vector<double> const &vec_flat,
61 Utils::Vector3i const &grid_size) {
62 auto const n_grid_points =
63 static_cast<std::size_t>(Utils::product(grid_size));
64 assert(vec_flat.size() == n_grid_points or vec_flat.size() == 1u);
65 std::vector<double> output_vector;
66 output_vector.reserve(n_grid_points);
67
68 auto const vec_begin = std::begin(vec_flat);
69 auto const vec_end = std::end(vec_flat);
70 if (vec_flat.size() == 1u) {
71 auto const uniform_value = vec_flat[0];
72 output_vector.assign(n_grid_points, uniform_value);
73 } else {
74 output_vector.assign(vec_begin, vec_end);
75 }
76
77 return output_vector;
78}
79
80template <class BoundaryModel, class DataType>
81void set_boundary_from_grid(BoundaryModel &boundary,
82 LatticeWalberla const &lattice,
83 std::vector<int> const &raster_flat,
84 std::vector<DataType> const &data_flat) {
85
86 auto const &conv = es2walberla<DataType, typename BoundaryModel::value_type>;
87 auto const grid_size = lattice.get_grid_dimensions();
88 auto const gl = static_cast<int>(lattice.get_ghost_layers());
89 assert(raster_flat.size() ==
90 static_cast<std::size_t>(Utils::product(grid_size)));
91 auto const n_y = static_cast<std::size_t>(grid_size[1]);
92 auto const n_z = static_cast<std::size_t>(grid_size[2]);
93
94 for (auto &block : *lattice.get_blocks()) {
95 auto const [size_i, size_j, size_k] = boundary.block_dims(block);
96 auto const offset = lattice.get_block_corner(block, true);
97 // Get field data which knows about the indices
98 // In the loop, i,j,k are in block-local coordinates
99 for (int i = -gl; i < size_i + gl; ++i) {
100 for (int j = -gl; j < size_j + gl; ++j) {
101 for (int k = -gl; k < size_k + gl; ++k) {
102 auto const node = offset + Utils::Vector3i{{i, j, k}};
103 auto const idx = (node + grid_size) % grid_size;
104 auto const index = static_cast<std::size_t>(idx[0]) * n_y * n_z +
105 static_cast<std::size_t>(idx[1]) * n_z +
106 static_cast<std::size_t>(idx[2]);
107 if (raster_flat[index]) {
108 auto const &value = data_flat[index];
109 auto const bc = BlockAndCell{&block, Cell(i, j, k)};
110 boundary.set_node_value_at_boundary(node, conv(value), bc);
111 }
112 }
113 }
114 }
115 }
116}
117
118} // namespace walberla
Vector implementation and trait types for boost qvm interoperability.
Definition Cell.hpp:96
Class that runs and controls the BlockForest in waLBerla.
Utils::Vector3i get_block_corner(IBlock const &block, bool lower) const
static double * block(double *p, std::size_t index, std::size_t size)
Definition elc.cpp:172
T product(Vector< T, N > const &v)
Definition Vector.hpp:375
VectorXd< 3 > Vector3d
Definition Vector.hpp:165
\file PackInfoPdfDoublePrecision.cpp \author pystencils
std::vector< double > fill_3D_scalar_array(std::vector< double > const &vec_flat, Utils::Vector3i const &grid_size)
Definition boundary.hpp:60
void set_boundary_from_grid(BoundaryModel &boundary, LatticeWalberla const &lattice, std::vector< int > const &raster_flat, std::vector< DataType > const &data_flat)
Definition boundary.hpp:81
std::vector< Utils::Vector3d > fill_3D_vector_array(std::vector< double > const &vec_flat, Utils::Vector3i const &grid_size)
Definition boundary.hpp:36