GADGET-4
ewald.h
Go to the documentation of this file.
1 /*******************************************************************************
2  * \copyright This file is part of the GADGET4 N-body/SPH code developed
3  * \copyright by Volker Springel. Copyright (C) 2014-2020 by Volker Springel
4  * \copyright (vspringel@mpa-garching.mpg.de) and all contributing authors.
5  *******************************************************************************/
6 
12 #ifndef EWALD_H
13 #define EWALD_H
14 
15 #include "gadgetconfig.h"
16 
17 #include <math.h>
18 #include <mpi.h>
19 #include <stdio.h>
20 #include <stdlib.h>
21 #include <string.h>
22 
23 #include "../data/allvars.h"
24 #include "../data/dtypes.h"
25 #include "../data/mymalloc.h"
26 #include "../data/symtensors.h"
27 #include "../gravity/ewaldtensors.h"
28 #include "../gravtree/gravtree.h"
29 #include "../io/io.h"
30 #include "../system/system.h"
31 
37 #if defined(GRAVITY_TALLBOX) && defined(FMM)
38 #error "FFM and GRAVITY_TALLBOX cannot be used together"
39 #endif
40 
41 #define EWLEVEL 6
42 
43 #define EN (1 << EWLEVEL)
44 
45 #define ENX (DBX * EN)
46 #define ENY (DBY * EN)
47 #define ENZ (DBZ * EN)
48 
53 class ewald : public intposconvert, public io_streamcount, public setcomm
54 {
55  public:
56  ewald(void) : setcomm("delayed init") {}
57 
58  void ewald_init(void);
59 
60  struct ewald_header
61  {
63  };
64 
66  {
70  };
71 
72  void ewald_corr(double dx, double dy, double dz, enum interpolate_options, ewald_data &fper);
73  void ewald_corr_exact(double dx, double dy, double dz, enum interpolate_options flag, ewald_data &fper);
74  void ewald_gridlookup(const MyIntPosType *p_intpos, const MyIntPosType *target_intpos, enum interpolate_options flag,
75  ewald_data &fper);
76 
77  double ewald_D0(double x, double y, double z);
78  vector<double> ewald_D1(double x, double y, double z);
79  symtensor2<double> ewald_D2(double x, double y, double z);
80  symtensor3<double> ewald_D3(double x, double y, double z);
81  symtensor4<double> ewald_D4(double x, double y, double z);
82  symtensor5<double> ewald_D5(double x, double y, double z);
83  symtensor6<double> ewald_D6(double x, double y, double z);
84  symtensor7<double> ewald_D7(double x, double y, double z);
85 
89 
90  private:
91  ewald_data *Ewd; // points to an [ENX + 1][ENY + 1][ENZ + 1] array
92 
93  inline int ewd_offset(int i, int j, int k) { return (i * (ENY + 1) + j) * (ENZ + 1) + k; }
94 
95  double Ewd_fac_intp[3];
96 
97  /*
98  * in D0phi we store the correction potential:
99  *
100  * phi = 1/x + pi/alpha^2 - sum_q (erfc(alpha |x-q|)/|x-q|) - 4pi/V sum_k exp(-k^2/(4alpha^2))/k^2 cos(k*x)
101  *
102  * in D1phi we store the correction force (first derivative of correction potential)
103  *
104  * dphi/dx_i
105  *
106  * in D2Phi we store the correction tensor (second derivatives of correction potential)
107  *
108  * d2phi/(dx_i dx_j)
109  *
110  * in D3phi we store the third order correction tensor (third derivatives of correction potential)
111  *
112  * d3phi/(dx_i dx_j dx_k)
113  */
114 
115  int ewald_is_initialized = 0;
116 
117  void test_interpolation_accuracy(void);
118 
119  void ewald_interpolate(double dx, double dy, double dz, enum interpolate_options, ewald_data &fper);
120 
121  public:
122 #if defined(EVALPOTENTIAL) && defined(FMM) && defined(PERIODIC) && !defined(PMGRID)
123  double ewald_gridlookup_origin_D0(void) { return Ewd->D0phi; }
124 #endif
125 };
126 
127 extern ewald Ewald;
128 
129 #endif
Ewald correction functionality.
Definition: ewald.h:54
symtensor3< double > ewald_D3(double x, double y, double z)
Definition: ewald.cc:1114
double ewald_D0(double x, double y, double z)
This function computes the potential correction term by means of Ewald summation.
Definition: ewald.cc:623
symtensor4< double > ewald_D4(double x, double y, double z)
Definition: ewald.cc:1241
void ewald_corr(double dx, double dy, double dz, enum interpolate_options, ewald_data &fper)
void ewald_init(void)
This function initializes tables with the correction force and the correction potential due to the pe...
Definition: ewald.cc:69
ewaldtensor6< double > ewald_P6(void)
symtensor6< double > ewald_D6(double x, double y, double z)
Definition: ewald.cc:1536
interpolate_options
Definition: ewald.h:66
@ POINTMASS
Definition: ewald.h:67
@ LONGRANGE_FORCETEST
Definition: ewald.h:69
@ MULTIPOLES
Definition: ewald.h:68
ewald(void)
Definition: ewald.h:56
symtensor7< double > ewald_D7(double x, double y, double z)
Definition: ewald.cc:1702
symtensor5< double > ewald_D5(double x, double y, double z)
Definition: ewald.cc:1385
vector< double > ewald_D1(double x, double y, double z)
This function computes the force correction term (difference between full force of infinite lattice a...
Definition: ewald.cc:807
void ewald_corr_exact(double dx, double dy, double dz, enum interpolate_options flag, ewald_data &fper)
ewaldtensor8< double > ewald_P8(void)
ewaldtensor10< double > ewald_P10(void)
symtensor2< double > ewald_D2(double x, double y, double z)
Definition: ewald.cc:996
void ewald_gridlookup(const MyIntPosType *p_intpos, const MyIntPosType *target_intpos, enum interpolate_options flag, ewald_data &fper)
Definition: ewald.cc:304
uint32_t MyIntPosType
Definition: dtypes.h:35
ewald Ewald
Definition: main.cc:42
#define ENZ
Definition: ewald.h:47
#define ENY
Definition: ewald.h:46
MyReal D0phi
Definition: gravtree.h:185