ExaDG
Loading...
Searching...
No Matches
interpolation.h
1/* ______________________________________________________________________
2 *
3 * ExaDG - High-Order Discontinuous Galerkin for the Exa-Scale
4 *
5 * Copyright (C) 2021 by the ExaDG authors
6 *
7 * This program is free software: you can redistribute it and/or modify
8 * it under the terms of the GNU General Public License as published by
9 * the Free Software Foundation, either version 3 of the License, or
10 * (at your option) any later version.
11 *
12 * This program is distributed in the hope that it will be useful,
13 * but WITHOUT ANY WARRANTY; without even the implied warranty of
14 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
15 * GNU General Public License for more details.
16 *
17 * You should have received a copy of the GNU General Public License
18 * along with this program. If not, see <https://www.gnu.org/licenses/>.
19 * ______________________________________________________________________
20 */
21
22#ifndef DEAL_SPECTRUM_INTERPOLATION
23#define DEAL_SPECTRUM_INTERPOLATION
24
25// std
26#include <mpi.h>
27#include <stdlib.h>
28
29// deal.II
30#include <deal.II/base/aligned_vector.h>
31#include <deal.II/base/point.h>
32#include <deal.II/base/quadrature_lib.h>
33#include <deal.II/base/utilities.h>
34#include <deal.II/fe/fe_dgq.h>
35#include <deal.II/fe/fe_values.h>
36#include <deal.II/matrix_free/evaluation_kernels.h>
37
38// ExaDG
39#include <exadg/postprocessor/spectral_analysis/setup.h>
40
41
42
43namespace dealspectrum
44{
50{
51public:
52 MPI_Comm const comm;
53 // reference to DEAL.SPECTRUM setup
54 Setup & s;
55 // is initialized?
56 bool initialized;
57 // source vector (values to be interpolated)
58 double * src = NULL;
59 // destination vector (interpolated values)
60 double * dst = NULL;
61 // shape function (gauss lobatto to equidistant)
62 dealii::AlignedVector<double> shape_values;
63 // number of cells in each direction
64 int cells;
65 // dimensions
66 int DIM;
67 // process local cell range on sfc
68 int start;
69 int end;
70 // points in each direction in cell & points per cell for source...
71 unsigned int points_source;
72 int dofs_source;
73 // ... and for target
74 unsigned int points_target;
75 int dofs_target;
76
81 Interpolator(MPI_Comm const & comm, Setup & s) : comm(comm), s(s), initialized(false)
82 {
83 }
84
85 virtual ~Interpolator()
86 {
87 // not initialized -> nothing to clean up
88 if(not initialized)
89 return;
90
91 delete[] src;
92 delete[] dst;
93 }
94
95 template<typename MAPPING>
96 void
97 init(MAPPING & MAP)
98 {
99 int start, end;
100 MAP.getLocalRange(start, end);
101
102 init(start, end);
103 }
104
105 void
106 init(int start, int end)
107 {
108 // check if already initialized
109 if(this->initialized)
110 return;
111 this->initialized = true;
112
113 // get cell range...
114 this->start = start;
115 this->end = end;
116 // ...number of local cells
117 cells = end - start;
118
119 // get dimensions
120 DIM = s.dim;
121
122 // number of Gauss-Lobatto points:
123 points_source = s.points_src;
124 // ... of equidistant points:
125 points_target = s.points_dst;
126
127 // number of Gauss-Lobatto points per cell and ...
128 dofs_source = dealii::Utilities::pow(points_source, DIM);
129 // ...number of equidistant points per cell
130 dofs_target = dealii::Utilities::pow(points_target, DIM);
131
132 // allocate memory for source (only needed if values are read by IO)...
133 src = new double[cells * dofs_source * DIM];
134 // ... and target
135 dst = new double[cells * dofs_target * DIM];
136
137 // fill shape values
138 fill_shape_values(shape_values);
139 }
140
147 std::vector<dealii::Point<1>>
148 get_equidistant_inner(unsigned int const n_points)
149 {
150 std::vector<dealii::Point<1>> points(n_points);
151
152 for(unsigned int i = 0; i < n_points; i++)
153 points[i][0] = (i + 0.5) / n_points;
154
155 return points;
156 }
157
163 void
164 fill_shape_values(dealii::AlignedVector<double> & shape_values)
165 {
166 // determine coefficients with deal.II functions
167 dealii::FE_DGQ<1> fe(points_source - 1);
168 dealii::FullMatrix<double> matrix(points_target, points_source);
169 dealii::FE_DGQArbitraryNodes<1>(dealii::Quadrature<1>(get_equidistant_inner(points_target)))
170 .get_interpolation_matrix(fe, matrix);
171
172 // ... and convert to linearized format
173 shape_values.resize(points_source * points_target);
174 for(unsigned int i = 0; i < points_source; ++i)
175 for(unsigned int q = 0; q < points_target; ++q)
176 shape_values[i * points_target + q] = matrix(q, i);
177 }
178
183 void
185 {
186 double const * src = this->src;
187 interpolate(src);
188 }
189
197 void
198 interpolate(double const *& src)
199 {
200 // allocate dst- and src-vector
201 dealii::AlignedVector<double> temp1(MAX(dofs_source, dofs_target));
202 dealii::AlignedVector<double> temp2(MAX(dofs_source, dofs_target));
203
204 // get start point of arrays
205 double const * src_ = src;
206 double * dst_ = dst;
207 // loop over all cells
208 for(int c = 0; c < cells; c++)
209 {
210 // loop over all velocity directions
211 for(int d = 0; d < DIM; d++)
212 {
213 // perform interpolation
214 if(DIM == 2)
215 {
216 // ... 2D
217 dealii::internal::
218 EvaluatorTensorProduct<dealii::internal::evaluate_general, 2, 0, 0, double>
219 eval_val(shape_values, shape_values, shape_values, points_source, points_target);
220 eval_val.template values<0, true, false>(src_, temp2.begin());
221 eval_val.template values<1, true, false>(temp2.begin(), temp1.begin());
222 }
223 else
224 {
225 // ... 3D
226 dealii::internal::
227 EvaluatorTensorProduct<dealii::internal::evaluate_general, 3, 0, 0, double>
228 eval_val(shape_values, shape_values, shape_values, points_source, points_target);
229 eval_val.template values<0, true, false>(src_, temp1.begin());
230 eval_val.template values<1, true, false>(temp1.begin(), temp2.begin());
231 eval_val.template values<2, true, false>(temp2.begin(), temp1.begin());
232 }
233 // write interpolated data permuted into destination vector
234 for(int i = 0; i < dofs_target; i++)
235 dst_[i * DIM + d] = temp1[i];
236 // go to next dimension
237 src_ += dofs_source;
238 }
239 // go to next cell
240 dst_ += dofs_target * DIM;
241 }
242 }
243
249 void
251 {
253 }
254
263 void
264 interpolate_simple(double *& src)
265 {
266 double * src_ = src;
267 double * dst_ = dst;
268
269 // loop over all cells
270 for(int c = 0; c < cells; c++)
271 {
272 // loop over all velocity directions
273 for(int d = 0; d < DIM; d++)
274 {
275 // write data permuted into destination vector
276 for(int i = 0; i < dofs_target; i++)
277 dst_[i * DIM + d] = src_[i];
278 // go to next dimension
279 src_ += dofs_source;
280 }
281 // go to next cell
282 dst_ += dofs_target * DIM;
283 }
284 }
285
291 void
292 deserialize(char const * filename)
293 {
294 deserialize(filename, src);
295 }
296
303 void
304 deserialize(char const * filename, double *& src)
305 {
306 io(0, filename, src);
307 }
308
314 void
315 serialize(char const * filename)
316 {
317 double const * src = this->src;
318 serialize(filename, src);
319 }
320
327 void
328 serialize(char const * filename, double const *& src)
329 {
330 io(1, filename, src);
331 }
332
340 template<typename VEC>
341 void
342 io(int type, char const * filename, VEC & src)
343 {
344 // dofs to read/write per field
345 unsigned long int dofs = cells * dofs_source;
346 // local displacement in file (in bytes)
347 MPI_Offset disp =
348 8 * sizeof(int) + static_cast<unsigned long int>(start) * dofs_source * sizeof(double) * DIM;
349
350 // create view
351 MPI_Datatype stype;
352 MPI_Type_contiguous(dofs, MPI_DOUBLE, &stype);
353 MPI_Type_commit(&stype);
354
355 // ooen file ...
356 MPI_File fh;
357 if(type == 0)
358 MPI_File_open(comm, filename, MPI_MODE_RDONLY, MPI_INFO_NULL, &fh);
359 else
360 MPI_File_open(comm, filename, MPI_MODE_WRONLY, MPI_INFO_NULL, &fh);
361
362 // ... set view
363 MPI_File_set_view(fh, disp, MPI_DOUBLE, MPI_DOUBLE, "native", MPI_INFO_NULL);
364
365 if(type == 0)
366 // ... read file
367 MPI_File_read_all(fh, (void *)src, dofs * DIM, MPI_DOUBLE, MPI_STATUSES_IGNORE);
368 else
369 // ... write file
370 MPI_File_write_all(fh, src, dofs * DIM, MPI_DOUBLE, MPI_STATUSES_IGNORE);
371
372 // ... close file
373 MPI_File_close(&fh);
374
375 // clean up
376 MPI_Type_free(&stype);
377 }
378};
379
380} // namespace dealspectrum
381
382#endif
Definition interpolation.h:50
std::vector< dealii::Point< 1 > > get_equidistant_inner(unsigned int const n_points)
Definition interpolation.h:148
void interpolate(double const *&src)
Definition interpolation.h:198
void deserialize(char const *filename, double *&src)
Definition interpolation.h:304
Interpolator(MPI_Comm const &comm, Setup &s)
Definition interpolation.h:81
void serialize(char const *filename, double const *&src)
Definition interpolation.h:328
void interpolate()
Definition interpolation.h:184
void serialize(char const *filename)
Definition interpolation.h:315
void interpolate_simple()
Definition interpolation.h:250
void deserialize(char const *filename)
Definition interpolation.h:292
void interpolate_simple(double *&src)
Definition interpolation.h:264
void io(int type, char const *filename, VEC &src)
Definition interpolation.h:342
void fill_shape_values(dealii::AlignedVector< double > &shape_values)
Definition interpolation.h:164
Definition setup.h:56