CosmoBolognaLib
Free Software C++/Python libraries for cosmological calculations
Field3D.cpp
Go to the documentation of this file.
1 /*******************************************************************
2  * Copyright (C) 2010 by Federico Marulli and Alfonso Veropalumbo *
3  * federico.marulli3@unibo.it *
4  * *
5  * This program is free software; you can redistribute it and/or *
6  * modify it under the terms of the GNU General Public License as *
7  * published by the Free Software Foundation; either version 2 of *
8  * the License, or (at your option) any later version. *
9  * *
10  * This program is distributed in the hope that it will be useful, *
11  * but WITHOUT ANY WARRANTY; without even the implied warranty of *
12  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
13  * GNU General Public License for more details. *
14  * *
15  * You should have received a copy of the GNU General Public *
16  * License along with this program; if not, write to the Free *
17  * Software Foundation, Inc., *
18  * 59 Temple Place - Suite 330, Boston, MA 02111-1307, USA. *
19  *******************************************************************/
20 
34 #include "Field3D.h"
35 
36 using namespace std;
37 
38 using namespace cbl;
39 using namespace data;
40 
41 
42 // ============================================================================
43 
44 
45 cbl::data::Field3D::Field3D (const double deltaR, const double minX, const double maxX, const double minY, const double maxY, const double minZ, const double maxZ)
46 {
47  set_parameters(deltaR, minX, maxX, minY, maxY, minZ, maxZ);
48 }
49 
50 
51 // ============================================================================
52 
53 
54 cbl::data::Field3D::Field3D (const int nx, const int ny, const int nz, const double minX, const double maxX, const double minY, const double maxY, const double minZ, const double maxZ)
55 {
56  set_parameters(nx, ny, nz, minX, maxX, minY, maxY, minZ, maxZ);
57 }
58 
59 
60 // ============================================================================
61 
62 
63 void cbl::data::Field3D::set_parameters (const double deltaR, const double minX, const double maxX, const double minY, const double maxY, const double minZ, const double maxZ)
64 {
65  m_MinX = minX;
66  m_MaxX = maxX;
67  m_MinY = minY;
68  m_MaxY = maxY;
69  m_MinZ = minZ;
70  m_MaxZ = maxZ;
71 
72  double DeltaX = m_MaxX-m_MinX;
73  double DeltaY = m_MaxY-m_MinY;
74  double DeltaZ = m_MaxZ-m_MinZ;
75 
76  m_nX = DeltaX/deltaR;
77  m_nX = (m_nX%2==0) ? m_nX : m_nX+1;
78  m_deltaX = DeltaX/m_nX;
79 
80  m_nY = DeltaY/deltaR;
81  m_nY = (m_nY%2==0) ? m_nY : m_nY+1;
82  m_deltaY = DeltaY/m_nY;
83 
84  m_nZ = DeltaZ/deltaR;
85  m_nZ = (m_nZ%2==0) ? m_nZ : m_nZ+1;
86  m_deltaZ = DeltaZ/m_nZ;
87 
88  m_Volume = DeltaX*DeltaY*DeltaZ;
89 
90  m_nZF = m_nZ*0.5+1;
91 
92  m_nCells = m_nX*m_nY*m_nZ;
93  m_nCells_Fourier = m_nX*m_nY*m_nZF;
94 
95 // Set coordinates of the grid
96 
97  double xfact = 2*par::pi/DeltaX;
98  double yfact = 2*par::pi/DeltaY;
99  double zfact = 2*par::pi/DeltaZ;
100 
101  for(int i=0; i<m_nX; i++){
102  m_X.push_back(m_MinX+(i+0.5)*m_deltaX);
103  m_kX.push_back( (i<=m_nX/2) ? xfact*i : -(m_nX-i)*xfact);
104  }
105 
106  for(int i=0; i<m_nY; i++){
107  m_Y.push_back(m_MinY+(i+0.5)*m_deltaY);
108  m_kY.push_back( (i<=m_nY/2) ? yfact*i : -(m_nY-i)*yfact);
109  }
110 
111  for(int i=0; i<m_nZ; i++)
112  m_Z.push_back(m_MinZ+(i+0.5)*m_deltaZ);
113 
114  for(int i=0; i<m_nZF; i++)
115  m_kZ.push_back(zfact*i);
116 
117 }
118 
119 
120 // ============================================================================
121 
122 
123 void cbl::data::Field3D::set_parameters (const int nx, const int ny, const int nz, const double minX, const double maxX, const double minY, const double maxY, const double minZ, const double maxZ)
124 {
125  m_MinX = minX;
126  m_MaxX = maxX;
127  m_MinY = minY;
128  m_MaxY = maxY;
129  m_MinZ = minZ;
130  m_MaxZ = maxZ;
131 
132  m_nX = (nx%2==0) ? nx : nx+1;
133  m_nY = (ny%2==0) ? ny : ny+1;
134  m_nZ = (nz%2==0) ? nz : nz+1;
135 
136  m_nZF = m_nZ*0.5+1;
137 
138  double DeltaX = m_MaxX-m_MinX;
139  double DeltaY = m_MaxY-m_MinY;
140  double DeltaZ = m_MaxZ-m_MinZ;
141 
142  m_deltaX = DeltaX/m_nX;
143  m_deltaY = DeltaY/m_nY;
144  m_deltaZ = DeltaZ/m_nZ;
145 
146  m_Volume = DeltaX*DeltaY*DeltaZ;
147 
148  m_nZF = m_nZ*0.5+1;
149 
150  m_nCells = m_nX*m_nY*m_nZ;
151  m_nCells_Fourier = m_nX*m_nY*m_nZF;
152 
153 // Set coordinates of the grid
154 
155  double xfact = 2*par::pi/DeltaX;
156  double yfact = 2*par::pi/DeltaY;
157  double zfact = 2*par::pi/DeltaZ;
158 
159  for(int i=0; i<m_nX; i++){
160  m_X.push_back(m_MinX+(i+0.5)*m_deltaX);
161  m_kX.push_back( (i<=m_nX/2) ? xfact*i : -(m_nX-i)*xfact);
162  }
163 
164  for(int i=0; i<m_nY; i++){
165  m_Y.push_back(m_MinY+(i+0.5)*m_deltaY);
166  m_kY.push_back( (i<=m_nY/2) ? yfact*i : -(m_nY-i)*yfact);
167  }
168 
169  for(int i=0; i<m_nZ; i++)
170  m_Z.push_back(m_MinZ+(i+0.5)*m_deltaZ);
171 
172  for(int i=0; i<m_nZF; i++)
173  m_kZ.push_back(zfact*i);
174 }
175 
176 
177 // ============================================================================
178 
179 
180 cbl::data::ScalarField3D::ScalarField3D (const double deltaR, const double minX, const double maxX, const double minY, const double maxY, const double minZ, const double maxZ) : Field3D(deltaR, minX, maxX, minY, maxY, minZ, maxZ)
181 {
182  m_field = fftw_alloc_real(m_nCells);
183  m_field_FourierSpace = fftw_alloc_complex(m_nCells_Fourier);
184 
185  for(int i=0;i<m_nCells;i++)
186  m_field[i]=0;
187 
188  for(int i=0;i<m_nCells_Fourier;i++){
189  m_field_FourierSpace[i][0] = 0;
190  m_field_FourierSpace[i][1] = 0;
191  }
192 }
193 
194 
195 // ============================================================================
196 
197 
198 cbl::data::ScalarField3D::ScalarField3D (const int nx, const int ny, const int nz, const double minX, const double maxX, const double minY, const double maxY, const double minZ, const double maxZ) : Field3D(nx, ny, nz, minX, maxX, minY, maxY, minZ, maxZ)
199 {
200  m_field = fftw_alloc_real(m_nCells);
201  m_field_FourierSpace = fftw_alloc_complex(m_nCells_Fourier);
202 
203  for (int i=0; i<m_nCells; i++)
204  m_field[i] = 0;
205 
206  for (int i=0; i<m_nCells_Fourier; i++) {
207  m_field_FourierSpace[i][0] = 0;
208  m_field_FourierSpace[i][1] = 0;
209  }
210 
211 }
212 
213 
214 // ============================================================================
215 
216 
218 {
219 
220  for (int i=0; i<m_nCells_Fourier; i++) {
221  m_field_FourierSpace[i][0] = 0;
222  m_field_FourierSpace[i][1] = 0;
223  }
224 
225  fftw_plan real2complex;
226  real2complex = fftw_plan_dft_r2c_3d(m_nX, m_nY, m_nZ, m_field, m_field_FourierSpace, FFTW_ESTIMATE);
227  fftw_execute(real2complex);
228  fftw_destroy_plan(real2complex);
229 
230  for(int i=0;i<m_nCells_Fourier;i++){
231  m_field_FourierSpace[i][0] = m_field_FourierSpace[i][0]/m_nCells;
232  m_field_FourierSpace[i][1] = m_field_FourierSpace[i][1]/m_nCells;
233  }
234 
235 }
236 
237 
238 // ============================================================================
239 
240 
242 {
243  for (int i=0; i<m_nCells; i++)
244  m_field[i] = 0;
245 
246  fftw_plan complex2real;
247  complex2real = fftw_plan_dft_c2r_3d(m_nX, m_nY, m_nZ, m_field_FourierSpace, m_field, FFTW_ESTIMATE);
248  fftw_execute(complex2real);
249  fftw_destroy_plan(complex2real);
250 
251 }
252 
253 
254 // ============================================================================
255 
256 
258 {
259  FourierTransformField();
260 
261  double xfact = 2*par::pi/(m_nX*m_deltaX);
262  double yfact = 2*par::pi/(m_nY*m_deltaY);
263  double zfact = 2*par::pi/(m_nZ*m_deltaZ);
264 
265  double ks2 = kernel_size*kernel_size;
266 
267  for (int i=0; i<m_nX; i++) {
268  double kx = (i<=m_nX/2) ? xfact*i : -(m_nX-i)*xfact;
269  for (int j=0; j<m_nY; j++) {
270  double ky = (j<=m_nY/2) ? yfact*j : -(m_nY-j)*yfact;
271  for (int k=0; k<m_nZF; k++) {
272  double kz = zfact*k;
273  double kk2 = kx*kx+ky*ky+kz*kz;
274  long int index = inds_to_index_Fourier(i,j,k);
275  m_field_FourierSpace[index][0] = m_field_FourierSpace[index][0]*exp(-0.5*kk2*ks2);
276  m_field_FourierSpace[index][1] = m_field_FourierSpace[index][1]*exp(-0.5*kk2*ks2);
277  }
278  }
279  }
280 
281  FourierAntiTransformField();
282 
283 }
284 
285 
286 // ============================================================================
287 
288 
289 void cbl::data::ScalarField3D::set_ScalarField (const double value, const int i, const int j, const int k, const bool add)
290 {
291  m_field[inds_to_index(i,j,k)] = (add) ? m_field[inds_to_index(i,j,k)] + value : value;
292 }
293 
294 
295 // ============================================================================
296 
297 
298 void cbl::data::ScalarField3D::set_ScalarField_FourierSpace_real (const double value, const int i, const int j, const int k, const bool add)
299 {
300  m_field_FourierSpace[inds_to_index_Fourier(i,j,k)][0] = (add) ? m_field_FourierSpace[inds_to_index_Fourier(i,j,k)][0] + value: value;
301 }
302 
303 
304 // ============================================================================
305 
306 
307 void cbl::data::ScalarField3D::set_ScalarField_FourierSpace_complex (const double value, const int i, const int j, const int k, const bool add)
308 {
309  m_field_FourierSpace[inds_to_index_Fourier(i,j,k)][1] = (add) ? m_field_FourierSpace[inds_to_index_Fourier(i,j,k)][1] + value: value;
310 }
311 
312 
313 // ============================================================================
314 
315 
316 double cbl::data::ScalarField3D::ScalarField (const int i, const int j, const int k) const
317 {
318  return m_field[inds_to_index(i,j,k)];
319 }
320 
321 
322 // ============================================================================
323 
324 
325 double cbl::data::ScalarField3D::ScalarField_FourierSpace_real (const int i, const int j, const int k) const
326 {
327  return m_field_FourierSpace[inds_to_index_Fourier(i,j,k)][0];
328 }
329 
330 
331 // ============================================================================
332 
333 
334 double cbl::data::ScalarField3D::ScalarField_FourierSpace_complex (const int i, const int j, const int k) const
335 {
336  return m_field_FourierSpace[inds_to_index_Fourier(i,j,k)][1];
337 }
338 
339 
340 // ============================================================================
341 
342 
343 double cbl::data::ScalarField3D::ScalarField (const std::vector<double> pos) const
344 {
345  int i = min(int((pos[0]-m_MinX)/m_deltaX), m_nX-1);
346  int j = min(int((pos[1]-m_MinY)/m_deltaY), m_nY-1);
347  int k = min(int((pos[2]-m_MinZ)/m_deltaZ), m_nZ-1);
348 
349  return m_field[inds_to_index(i,j,k)];
350 }
351 
352 
353 // ============================================================================
354 
355 
357 {
358  vector<double> vv(m_field, m_field+m_nCells);
359  return vv;
360 }
361 
362 
363 // ============================================================================
364 
365 
367 {
368  for(int i=0;i<m_nCells;i++)
369  m_field[i]=0;
370 
371  for(int i=0;i<m_nCells_Fourier;i++){
372  m_field_FourierSpace[i][0] = 0;
373  m_field_FourierSpace[i][1] = 0;
374  }
375 }
376 
377 // ============================================================================
378 
379 
380 cbl::data::VectorField3D::VectorField3D (const double deltaR, const double minX, const double maxX, const double minY, const double maxY, const double minZ, const double maxZ) : Field3D(deltaR, minX, maxX, minY, maxY, minZ, maxZ)
381 {
382  m_field.resize(3);
383 
384  m_field_FourierSpace.resize(3);
385 
386 
387  m_field[0] = fftw_alloc_real(m_nCells);
388  m_field[1] = fftw_alloc_real(m_nCells);
389  m_field[2] = fftw_alloc_real(m_nCells);
390 
391  m_field_FourierSpace[0] = fftw_alloc_complex(m_nCells_Fourier);
392  m_field_FourierSpace[1] = fftw_alloc_complex(m_nCells_Fourier);
393  m_field_FourierSpace[2] = fftw_alloc_complex(m_nCells_Fourier);
394 
395  for(int i=0;i<m_nCells;i++){
396  m_field[0][i]=0;
397  m_field[1][i]=0;
398  m_field[2][i]=0;
399  }
400  for(int i=0;i<m_nCells_Fourier;i++){
401  m_field_FourierSpace[0][i][0] = 0;
402  m_field_FourierSpace[0][i][1] = 0;
403  m_field_FourierSpace[1][i][0] = 0;
404  m_field_FourierSpace[1][i][1] = 0;
405  m_field_FourierSpace[2][i][0] = 0;
406  m_field_FourierSpace[2][i][1] = 0;
407  }
408 }
409 
410 
411 // ============================================================================
412 
413 
414 cbl::data::VectorField3D::VectorField3D (const int nx, const int ny, const int nz, const double minX, const double maxX, const double minY, const double maxY, const double minZ, const double maxZ) : Field3D(nx, ny, nz, minX, maxX, minY, maxY, minZ, maxZ)
415 {
416  m_field.resize(3);
417 
418  m_field_FourierSpace.resize(3);
419 
420 
421  m_field[0] = fftw_alloc_real(m_nCells);
422  m_field[1] = fftw_alloc_real(m_nCells);
423  m_field[2] = fftw_alloc_real(m_nCells);
424 
425  m_field_FourierSpace[0] = fftw_alloc_complex(m_nCells_Fourier);
426  m_field_FourierSpace[1] = fftw_alloc_complex(m_nCells_Fourier);
427  m_field_FourierSpace[2] = fftw_alloc_complex(m_nCells_Fourier);
428 
429  for(int i=0;i<m_nCells;i++){
430  m_field[0][i]=0;
431  m_field[1][i]=0;
432  m_field[2][i]=0;
433  }
434  for(int i=0;i<m_nCells_Fourier;i++){
435  m_field_FourierSpace[0][i][0] = 0;
436  m_field_FourierSpace[0][i][1] = 0;
437  m_field_FourierSpace[1][i][0] = 0;
438  m_field_FourierSpace[1][i][1] = 0;
439  m_field_FourierSpace[2][i][0] = 0;
440  m_field_FourierSpace[2][i][1] = 0;
441  }
442 
443 }
444 
445 
446 // ============================================================================
447 
448 
450 {
451 
452  for(int i=0;i<m_nCells_Fourier;i++){
453  m_field_FourierSpace[0][i][0] = 0;
454  m_field_FourierSpace[0][i][1] = 0;
455  m_field_FourierSpace[1][i][0] = 0;
456  m_field_FourierSpace[1][i][1] = 0;
457  m_field_FourierSpace[2][i][0] = 0;
458  m_field_FourierSpace[2][i][1] = 0;
459  }
460 
461  fftw_plan real2complex;
462 
463  real2complex = fftw_plan_dft_r2c_3d(m_nX, m_nY, m_nZ, m_field[0], m_field_FourierSpace[0], FFTW_ESTIMATE);
464  fftw_execute(real2complex);
465  fftw_destroy_plan(real2complex);
466 
467  real2complex = fftw_plan_dft_r2c_3d(m_nX, m_nY, m_nZ, m_field[1], m_field_FourierSpace[1], FFTW_ESTIMATE);
468  fftw_execute(real2complex);
469  fftw_destroy_plan(real2complex);
470 
471  real2complex = fftw_plan_dft_r2c_3d(m_nX, m_nY, m_nZ, m_field[2], m_field_FourierSpace[2], FFTW_ESTIMATE);
472  fftw_execute(real2complex);
473  fftw_destroy_plan(real2complex);
474 
475  for(int i=0;i<m_nCells_Fourier;i++){
476  m_field_FourierSpace[0][i][0] = m_field_FourierSpace[0][i][0]/m_nCells;
477  m_field_FourierSpace[0][i][1] = m_field_FourierSpace[0][i][1]/m_nCells;
478  m_field_FourierSpace[1][i][0] = m_field_FourierSpace[1][i][0]/m_nCells;
479  m_field_FourierSpace[1][i][1] = m_field_FourierSpace[1][i][1]/m_nCells;
480  m_field_FourierSpace[2][i][0] = m_field_FourierSpace[2][i][0]/m_nCells;
481  m_field_FourierSpace[2][i][1] = m_field_FourierSpace[2][i][1]/m_nCells;
482  }
483 
484 }
485 
486 
487 // ============================================================================
488 
489 
491 {
492  for (int i=0; i<m_nCells; i++) {
493  m_field[0][i] = 0;
494  m_field[1][i] = 0;
495  m_field[2][i] = 0;
496  }
497 
498  fftw_plan complex2real;
499 
500  complex2real = fftw_plan_dft_c2r_3d(m_nX, m_nY, m_nZ, m_field_FourierSpace[0], m_field[0], FFTW_ESTIMATE);
501  fftw_execute(complex2real);
502  fftw_destroy_plan(complex2real);
503 
504  complex2real = fftw_plan_dft_c2r_3d(m_nX, m_nY, m_nZ, m_field_FourierSpace[1], m_field[1], FFTW_ESTIMATE);
505  fftw_execute(complex2real);
506  fftw_destroy_plan(complex2real);
507 
508  complex2real = fftw_plan_dft_c2r_3d(m_nX, m_nY, m_nZ, m_field_FourierSpace[2], m_field[2], FFTW_ESTIMATE);
509  fftw_execute(complex2real);
510  fftw_destroy_plan(complex2real);
511 
512 }
513 
514 
515 // ============================================================================
516 
517 
518 void cbl::data::VectorField3D::set_VectorField (const vector<double> value, const int i, const int j, const int k, const bool add)
519 {
520  m_field[0][inds_to_index(i,j,k)] = (add) ? m_field[0][inds_to_index(i,j,k)]+value[0]: value[0];
521  m_field[1][inds_to_index(i,j,k)] = (add) ? m_field[1][inds_to_index(i,j,k)]+value[1]: value[1];
522  m_field[2][inds_to_index(i,j,k)] = (add) ? m_field[2][inds_to_index(i,j,k)]+value[2]: value[2];
523 }
524 
525 
526 // ============================================================================
527 
528 
529 void cbl::data::VectorField3D::set_VectorField_FourierSpace_real (const std::vector<double> value, const int i, const int j, const int k, const bool add)
530 {
531  m_field_FourierSpace[0][inds_to_index_Fourier(i,j,k)][0] = (add) ? m_field_FourierSpace[0][inds_to_index(i,j,k)][0]+value[0] : value[0];
532  m_field_FourierSpace[1][inds_to_index_Fourier(i,j,k)][0] = (add) ? m_field_FourierSpace[1][inds_to_index(i,j,k)][0]+value[1] : value[1];
533  m_field_FourierSpace[2][inds_to_index_Fourier(i,j,k)][0] = (add) ? m_field_FourierSpace[2][inds_to_index(i,j,k)][0]+value[2] : value[2];
534 }
535 
536 
537 // ============================================================================
538 
539 
540 void cbl::data::VectorField3D::set_VectorField_FourierSpace_complex(const vector<double> value, const int i, const int j, const int k, const bool add)
541 {
542  m_field_FourierSpace[0][inds_to_index_Fourier(i,j,k)][1] = (add) ? m_field_FourierSpace[0][inds_to_index(i,j,k)][1]+value[0] : value[0];
543  m_field_FourierSpace[1][inds_to_index_Fourier(i,j,k)][1] = (add) ? m_field_FourierSpace[1][inds_to_index(i,j,k)][1]+value[1] : value[1];
544  m_field_FourierSpace[2][inds_to_index_Fourier(i,j,k)][1] = (add) ? m_field_FourierSpace[2][inds_to_index(i,j,k)][1]+value[2] : value[2];
545 }
546 
547 
548 // ============================================================================
549 
550 
551 vector<double> cbl::data::VectorField3D::VectorField(const int i, const int j, const int k) const
552 {
553  return {m_field[0][inds_to_index(i,j,k)], m_field[1][inds_to_index(i,j,k)], m_field[2][inds_to_index(i,j,k)]};
554 }
555 
556 
557 // ============================================================================
558 
559 
560 vector<double> cbl::data::VectorField3D::VectorField_FourierSpace_real(const int i, const int j, const int k) const
561 {
562  return {m_field_FourierSpace[0][inds_to_index_Fourier(i,j,k)][0], m_field_FourierSpace[1][inds_to_index_Fourier(i,j,k)][0], m_field_FourierSpace[2][inds_to_index_Fourier(i,j,k)][0]};
563 }
564 
565 
566 // ============================================================================
567 
568 
569 vector<double> cbl::data::VectorField3D::VectorField_FourierSpace_complex(const int i, const int j, const int k) const
570 {
571  return {m_field_FourierSpace[0][inds_to_index_Fourier(i,j,k)][1], m_field_FourierSpace[1][inds_to_index_Fourier(i,j,k)][1], m_field_FourierSpace[2][inds_to_index_Fourier(i,j,k)][1]};
572 }
573 
574 
575 // ============================================================================
576 
577 
578 std::vector<double> cbl::data::VectorField3D::VectorField (const std::vector<double> pos) const
579 {
580  int i = min(int((pos[0]-m_MinX)/m_deltaX), m_nX-1);
581  int j = min(int((pos[1]-m_MinY)/m_deltaY), m_nY-1);
582  int k = min(int((pos[2]-m_MinZ)/m_deltaZ), m_nZ-1);
583 
584  return {m_field[0][inds_to_index(i,j,k)], m_field[1][inds_to_index(i,j,k)], m_field[2][inds_to_index(i,j,k)]};
585 }
586 
587 
588 // ============================================================================
589 
590 
592 {
593  for(int i=0;i<m_nCells;i++){
594  m_field[0][i]=0;
595  m_field[1][i]=0;
596  m_field[2][i]=0;
597  }
598  for(int i=0;i<m_nCells_Fourier;i++){
599  m_field_FourierSpace[0][i][0] = 0;
600  m_field_FourierSpace[0][i][1] = 0;
601  m_field_FourierSpace[1][i][0] = 0;
602  m_field_FourierSpace[1][i][1] = 0;
603  m_field_FourierSpace[2][i][0] = 0;
604  m_field_FourierSpace[2][i][1] = 0;
605  }
606 }
The class field3D.
The class Field3D.
Definition: Field3D.h:51
int m_nCells
number of cells
Definition: Field3D.h:68
Field3D()=default
default constructor
void set_parameters(const double deltaR, const double minX, const double maxX, const double minY, const double maxY, const double minZ, const double maxZ)
set the parameters
Definition: Field3D.cpp:63
virtual std::vector< double > ScalarField_FourierSpace_real() const
get the value of the scalar field, Fourier space, real part
Definition: Field3D.h:478
virtual std::vector< double > ScalarField_FourierSpace_complex() const
get the value of the scalar field, Fourier space, complex part
Definition: Field3D.h:486
int m_nCells_Fourier
number of cells, Fourier space
Definition: Field3D.h:71
ScalarField3D()=default
default constructor
std::vector< double > ScalarField() const
get the values of the scalar field
Definition: Field3D.cpp:356
void set_ScalarField_FourierSpace_complex(const double value, const int i, const int j, const int k, const bool add=0)
set the value of the scalar field in Fourier space, real part
Definition: Field3D.cpp:307
void FourierTransformField()
perform the Fourier transform on the field
Definition: Field3D.cpp:217
void FourierAntiTransformField()
perform the anti-Fourier transform on the field
Definition: Field3D.cpp:241
double * m_field
scalar field
Definition: Field3D.h:653
void GaussianConvolutionField(const double kernel_size)
perform a smoothing of the field with a gaussian kernel
Definition: Field3D.cpp:257
void reset()
set to 0 the fields
Definition: Field3D.cpp:366
void set_ScalarField_FourierSpace_real(const double value, const int i, const int j, const int k, const bool add=0)
set the value of the scalar field in Fourier space, real part
Definition: Field3D.cpp:298
void set_ScalarField(const double value, const int i, const int j, const int k, const bool add=0)
set the value of the scalar field
Definition: Field3D.cpp:289
fftw_complex * m_field_FourierSpace
fourier transform of the scalar field
Definition: Field3D.h:656
void reset()
set to 0 the fields
Definition: Field3D.cpp:591
void FourierTransformField()
perform the anti-Fourier transform on the field
Definition: Field3D.cpp:449
void set_VectorField_FourierSpace_real(const std::vector< double > value, const int i, const int j, const int k, const bool add=0)
set the value of the vector field, Fourier space, real part
Definition: Field3D.cpp:529
void set_VectorField_FourierSpace_complex(const std::vector< double > value, const int i, const int j, const int k, const bool add=0)
set the value of the vector field, Fourier space, complex part
Definition: Field3D.cpp:540
std::vector< double * > m_field
vector field
Definition: Field3D.h:875
std::vector< fftw_complex * > m_field_FourierSpace
vector field in fourier space
Definition: Field3D.h:878
std::vector< double > VectorField_FourierSpace_complex(const int i, const int j, const int k) const
get the value of the vector field, Fourier space, complex part
Definition: Field3D.cpp:569
void set_VectorField(const std::vector< double > value, const int i, const int j, const int k, const bool add=0)
set the value of the vectorr field
Definition: Field3D.cpp:518
std::vector< double > VectorField_FourierSpace_real(const int i, const int j, const int k) const
get the value of the vector field, Fourier space, real part
Definition: Field3D.cpp:560
VectorField3D()=default
default constructor
std::vector< double > VectorField(const int i, const int j, const int k) const
get the value of the vector field
Definition: Field3D.cpp:551
void FourierAntiTransformField()
perform the anti-Fourier transform on the field
Definition: Field3D.cpp:490
static const double pi
: the ratio of a circle's circumference to its diameter
Definition: Constants.h:199
The global namespace of the CosmoBolognaLib
Definition: CAMB.h:38