GADGET-4
ewald.cc
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 #include "gadgetconfig.h"
13 
14 #include <math.h>
15 #include <mpi.h>
16 #include <stdio.h>
17 #include <stdlib.h>
18 #include <string.h>
19 
20 #include "../data/allvars.h"
21 #include "../data/dtypes.h"
22 #include "../data/mymalloc.h"
23 #include "../gravity/ewald.h"
24 #include "../gravity/ewaldtensors.h"
25 #include "../gravtree/gravtree.h"
26 #include "../io/io.h"
27 #include "../main/simulation.h"
28 #include "../mpi_utils/shared_mem_handler.h"
29 #include "../sort/cxxsort.h"
30 #include "../system/system.h"
31 
70 {
71  mpi_printf("EWALD: initialize Ewald correction...\n");
72 
76 
77  Ewd = (ewald_data *)Mem.mymalloc("Ewd", sizeof(ewald_data) * (ENX + 1) * (ENY + 1) * (ENZ + 1));
78 
79  char buf[200];
80  sprintf(buf, "ewald_table_%d-%d-%d_%d-%d-%d_precision%d-order%d.dat", LONG_X, LONG_Y, LONG_Z, ENX, ENY, ENZ, (int)sizeof(MyReal),
82 
83  int recomputeflag = 0;
84 
85  if(ThisTask == 0)
86  {
87  FILE *fd;
88  if((fd = fopen(buf, "r")))
89  {
90  mpi_printf("\nEWALD: reading Ewald tables from file `%s'\n", buf);
91 
92  ewald_header tabh;
93  my_fread(&tabh, sizeof(ewald_header), 1, fd);
94 
95 #ifndef GRAVITY_TALLBOX
96  int ewaldtype = -1;
97 #else
98  int ewaldtype = GRAVITY_TALLBOX + 1;
99 #endif
100  if(tabh.resx != ENX || tabh.resy != ENY || tabh.resz != ENZ || tabh.varsize != sizeof(MyFloat) ||
101  tabh.ewaldtype != ewaldtype)
102  {
103  mpi_printf("\nEWALD: something's wrong with this table file. Discarding it.\n");
104  recomputeflag = 1;
105  }
106  else
107  {
108  my_fread(Ewd, sizeof(ewald_data), (ENX + 1) * (ENY + 1) * (ENZ + 1), fd);
109 
110  recomputeflag = 0;
111  }
112  fclose(fd);
113  }
114  else
115  recomputeflag = 1;
116  }
117 
118  MPI_Bcast(&recomputeflag, 1, MPI_INT, 0, Communicator);
119 
120  if(recomputeflag)
121  {
122  mpi_printf("\nEWALD: No usable Ewald tables in file `%s' found. Recomputing them...\n", buf);
123 
124  /* ok, let's recompute things. Actually, we do that in parallel. */
125 
126  int size = (ENX + 1) * (ENY + 1) * (ENZ + 1);
127  int first, count;
128 
129  subdivide_evenly(size, NTask, ThisTask, &first, &count);
130 
131  for(int n = first; n < first + count; n++)
132  {
133  int i = n / ((ENY + 1) * (ENZ + 1));
134  int j = (n - i * (ENY + 1) * (ENZ + 1)) / (ENZ + 1);
135  int k = (n - i * (ENY + 1) * (ENZ + 1) - j * (ENZ + 1));
136 
137  if(ThisTask == 0)
138  {
139  if(((n - first) % (count / 20)) == 0)
140  {
141  printf("%4.1f percent done\n", (n - first) / (count / 100.0));
142  myflush(stdout);
143  }
144  }
145 
146  double xx = 0.5 * DBX * (1.0 / LONG_X) * ((double)i) / ENX;
147  double yy = 0.5 * DBY * (1.0 / LONG_Y) * ((double)j) / ENY;
148  double zz = 0.5 * DBZ * (1.0 / LONG_Z) * ((double)k) / ENZ;
149 
150  ewald_data *ewdp = Ewd + ewd_offset(i, j, k);
151 
152 #ifndef GRAVITY_TALLBOX
153  ewdp->D0phi = ewald_D0(xx, yy, zz);
154  ewdp->D1phi = ewald_D1(xx, yy, zz);
155  ewdp->D2phi = ewald_D2(xx, yy, zz);
156  ewdp->D3phi = ewald_D3(xx, yy, zz);
157 #if(HIGHEST_NEEDEDORDER_EWALD_DPHI + EWALD_TAYLOR_ORDER) >= 4
158  ewdp->D4phi = ewald_D4(xx, yy, zz);
159 #endif
160 #if(HIGHEST_NEEDEDORDER_EWALD_DPHI + EWALD_TAYLOR_ORDER) >= 5
161  ewdp->D5phi = ewald_D5(xx, yy, zz);
162 #endif
163 #if(HIGHEST_NEEDEDORDER_EWALD_DPHI + EWALD_TAYLOR_ORDER) >= 6
164  ewdp->D6phi = ewald_D6(xx, yy, zz);
165 #endif
166 #if(HIGHEST_NEEDEDORDER_EWALD_DPHI + EWALD_TAYLOR_ORDER) >= 7
167  ewdp->D7phi = ewald_D7(xx, yy, zz);
168 #endif
169 #else
170  ewdp->D0phi = ewald_D0(xx, zz, yy);
171 
172  vector<double> force = ewald_D1(yy, zz, xx);
173 
174  switch(GRAVITY_TALLBOX)
175  {
176  case 0:
177  ewdp->D1phi[0] = force[2];
178  ewdp->D1phi[1] = force[0];
179  ewdp->D1phi[2] = force[1];
180  break;
181 
182  case 1:
183  ewdp->D1phi[0] = force[0];
184  ewdp->D1phi[1] = force[2];
185  ewdp->D1phi[2] = force[1];
186  break;
187 
188  case 2:
189  ewdp->D1phi[0] = force[0];
190  ewdp->D1phi[1] = force[1];
191  ewdp->D1phi[2] = force[2];
192  break;
193  }
194 #endif
195  }
196 
197  int *recvcnts = (int *)Mem.mymalloc("recvcnts", NTask * sizeof(int));
198  int *recvoffs = (int *)Mem.mymalloc("recvoffs", NTask * sizeof(int));
199 
200  for(int i = 0; i < NTask; i++)
201  {
202  int off, cnt;
203  subdivide_evenly(size, NTask, i, &off, &cnt);
204  recvcnts[i] = cnt * sizeof(ewald_data);
205  recvoffs[i] = off * sizeof(ewald_data);
206  }
207 
208  MPI_Allgatherv(MPI_IN_PLACE, size * sizeof(ewald_data), MPI_BYTE, Ewd, recvcnts, recvoffs, MPI_BYTE, Communicator);
209 
210  Mem.myfree(recvoffs);
211  Mem.myfree(recvcnts);
212 
213  mpi_printf("\nEWALD: writing Ewald tables to file `%s'\n", buf);
214  if(ThisTask == 0)
215  {
216  FILE *fd;
217  if((fd = fopen(buf, "w")))
218  {
219  ewald_header tabh;
220  tabh.resx = ENX;
221  tabh.resy = ENY;
222  tabh.resz = ENZ;
223  tabh.varsize = sizeof(MyFloat);
224 #ifndef GRAVITY_TALLBOX
225  tabh.ewaldtype = -1;
226 #else
227  tabh.ewaldtype = GRAVITY_TALLBOX + 1;
228 #endif
229  my_fwrite(&tabh, sizeof(ewald_header), 1, fd);
230 
231  my_fwrite(Ewd, sizeof(ewald_data), (ENX + 1) * (ENY + 1) * (ENZ + 1), fd);
232  fclose(fd);
233  }
234  else
235  Terminate("can't write to file '%s'\n", buf);
236  }
237  }
238  else
239  {
240  /* here we got them from disk */
241  int len = (ENX + 1) * (ENY + 1) * (ENZ + 1) * sizeof(ewald_data);
242  MPI_Bcast(Ewd, len, MPI_BYTE, 0, Communicator);
243  }
244 
245  Ewd_fac_intp[0] = 2.0 * EN * LONG_X / All.BoxSize;
246  Ewd_fac_intp[1] = 2.0 * EN * LONG_Y / All.BoxSize;
247  Ewd_fac_intp[2] = 2.0 * EN * LONG_Z / All.BoxSize;
248 
249  /* now scale things to the boxsize that is actually used */
250  for(int i = 0; i <= ENX; i++)
251  for(int j = 0; j <= ENY; j++)
252  for(int k = 0; k <= ENZ; k++)
253  {
254  ewald_data *ewdp = Ewd + ewd_offset(i, j, k);
255 
256  ewdp->D0phi *= 1 / All.BoxSize; /* potential */
257  ewdp->D1phi *= 1 / pow(All.BoxSize, 2);
258  ewdp->D2phi *= 1 / pow(All.BoxSize, 3);
259  ewdp->D3phi *= 1 / pow(All.BoxSize, 4);
260 #if(HIGHEST_NEEDEDORDER_EWALD_DPHI + EWALD_TAYLOR_ORDER) >= 4
261  ewdp->D4phi *= 1 / pow(All.BoxSize, 5);
262 #endif
263 #if(HIGHEST_NEEDEDORDER_EWALD_DPHI + EWALD_TAYLOR_ORDER) >= 5
264  ewdp->D5phi *= 1 / pow(All.BoxSize, 6);
265 #endif
266 #if(HIGHEST_NEEDEDORDER_EWALD_DPHI + EWALD_TAYLOR_ORDER) >= 6
267  ewdp->D6phi *= 1 / pow(All.BoxSize, 7);
268 #endif
269 #if(HIGHEST_NEEDEDORDER_EWALD_DPHI + EWALD_TAYLOR_ORDER) >= 7
270  ewdp->D7phi *= 1 / pow(All.BoxSize, 8);
271 #endif
272  }
273 
274  mpi_printf("EWALD: Initialization of periodic boundaries finished.\n");
275 
276  ewald_is_initialized = 1;
277 
279  {
280  // We actually have multiple shared memory nodes in which we set aside one MPI rank for shared memory communication.
281  // In this case, move the ewaldtable to the communication rank in order to consume this memory only once on the node
282 
283  if(Shmem.Island_ThisTask == 0)
284  {
285  size_t tab_len = sizeof(ewald_data) * (ENX + 1) * (ENY + 1) * (ENZ + 1);
286 
287  MPI_Send(&tab_len, sizeof(tab_len), MPI_BYTE, Shmem.MyShmRankInGlobal, TAG_EWALD_ALLOC, MPI_COMM_WORLD);
288  MPI_Send(Ewd, tab_len, MPI_BYTE, Shmem.MyShmRankInGlobal, TAG_DMOM, MPI_COMM_WORLD);
289  }
290 
291  Mem.myfree(Ewd);
292 
293  ptrdiff_t off;
294  MPI_Bcast(&off, sizeof(ptrdiff_t), MPI_BYTE, Shmem.Island_NTask - 1, Shmem.SharedMemComm);
295 
296  Ewd = (ewald_data *)((char *)Shmem.SharedMemBaseAddr[Shmem.Island_NTask - 1] + off);
297  }
298 
299 #ifdef EWALD_TEST
300  test_interpolation_accuracy();
301 #endif
302 }
303 
304 void ewald::ewald_gridlookup(const MyIntPosType *p_intpos, const MyIntPosType *target_intpos, enum interpolate_options flag,
305  ewald_data &fper)
306 {
307  // we determine the closest available point in our Ewald look-up table
308 
309  static MyIntPosType const halflen = ((MyIntPosType)1) << ((BITS_FOR_POSITIONS - 1) - (EWLEVEL + 1));
310  static MyIntPosType const intlen = halflen << 1;
311  static MyIntPosType const ewaldmask = ~(intlen - 1);
312 
313  MyIntPosType temppos[3] = {p_intpos[0] - target_intpos[0], p_intpos[1] - target_intpos[1], p_intpos[2] - target_intpos[2]};
314 
315  MyIntPosType gridpos[3];
316  gridpos[0] = (temppos[0] + halflen) & ewaldmask;
317  gridpos[1] = (temppos[1] + halflen) & ewaldmask;
318  gridpos[2] = (temppos[2] + halflen) & ewaldmask;
319 
320  vector<double> off;
321  nearest_image_intpos_to_pos(temppos, gridpos, off.da);
322 
323  int i = (gridpos[0] >> (BITS_FOR_POSITIONS - (EWLEVEL + 1)));
324  int j = (gridpos[1] >> (BITS_FOR_POSITIONS - (EWLEVEL + 1)));
325  int k = (gridpos[2] >> (BITS_FOR_POSITIONS - (EWLEVEL + 1)));
326 
327  int signx = 1, signy = 1, signz = 1;
328 
329  if(i > EN)
330  {
331  i = 2 * EN - i;
332  signx = -1;
333  }
334  else if(i == EN && gridpos[0] < temppos[0])
335  signx = -1;
336 
337  if(j > EN)
338  {
339  j = 2 * EN - j;
340  signy = -1;
341  }
342  else if(j == EN && gridpos[1] < temppos[1])
343  signy = -1;
344 
345  if(k > EN)
346  {
347  k = 2 * EN - k;
348  signz = -1;
349  }
350  else if(k == EN && gridpos[2] < temppos[2])
351  signz = -1;
352 
353  fper = Ewd[ewd_offset(i, j, k)];
354 
355  /* change signs as needed */
356 
357  fper.D1phi[0] *= signx;
358  fper.D1phi[1] *= signy;
359  fper.D1phi[2] *= signz;
360 
361  fper.D2phi[qXY] *= signx * signy;
362  fper.D2phi[qXZ] *= signx * signz;
363  fper.D2phi[qYZ] *= signy * signz;
364 
365  fper.D3phi[dXXX] *= signx;
366  fper.D3phi[dXXY] *= signy;
367  fper.D3phi[dXXZ] *= signz;
368  fper.D3phi[dXYY] *= signx;
369  fper.D3phi[dXYZ] *= signx * signy * signz;
370  fper.D3phi[dXZZ] *= signx;
371  fper.D3phi[dYYY] *= signy;
372  fper.D3phi[dYYZ] *= signz;
373  fper.D3phi[dYZZ] *= signy;
374  fper.D3phi[dZZZ] *= signz;
375 
376 #if EWALD_TAYLOR_ORDER == 3
377 
378  fper.D4phi[sXXXY] *= signx * signy;
379  fper.D4phi[sXYYY] *= signx * signy;
380  fper.D4phi[sXXXZ] *= signx * signz;
381  fper.D4phi[sXZZZ] *= signx * signz;
382  fper.D4phi[sYYYZ] *= signy * signz;
383  fper.D4phi[sYZZZ] *= signy * signz;
384  fper.D4phi[sXXYZ] *= signy * signz;
385  fper.D4phi[sXYYZ] *= signx * signz;
386  fper.D4phi[sXYZZ] *= signx * signy;
387 
388  // now Taylor corrections
389 
390  fper.D0phi += fper.D1phi * off + 0.5 * ((fper.D2phi * off) * off) + (1.0 / 6) * (((fper.D3phi * off) * off) * off);
391  fper.D1phi += fper.D2phi * off + 0.5 * ((fper.D3phi * off) * off) + (1.0 / 6) * (((fper.D4phi * off) * off) * off);
392 
393  if(flag == POINTMASS)
394  return;
395 
396 #if(HIGHEST_NEEDEDORDER_EWALD_DPHI + EWALD_TAYLOR_ORDER) >= 5
397  fper.D5phi[rXXXXX] *= signx;
398  fper.D5phi[rYYYYY] *= signy;
399  fper.D5phi[rZZZZZ] *= signz;
400 
401  fper.D5phi[rXXXXY] *= signy;
402  fper.D5phi[rXXXXZ] *= signz;
403  fper.D5phi[rXYYYY] *= signx;
404  fper.D5phi[rXZZZZ] *= signx;
405  fper.D5phi[rYYYYZ] *= signz;
406  fper.D5phi[rYZZZZ] *= signy;
407 
408  fper.D5phi[rXXXYY] *= signx;
409  fper.D5phi[rXXXZZ] *= signx;
410  fper.D5phi[rXXYYY] *= signy;
411  fper.D5phi[rXXZZZ] *= signz;
412  fper.D5phi[rYYYZZ] *= signy;
413  fper.D5phi[rYYZZZ] *= signz;
414 
415  fper.D5phi[rXXYZZ] *= signy;
416  fper.D5phi[rXXYYZ] *= signz;
417  fper.D5phi[rXYYZZ] *= signx;
418 
419  fper.D5phi[rXXXYZ] *= signx * signy * signz;
420  fper.D5phi[rXYYYZ] *= signx * signy * signz;
421  fper.D5phi[rXYZZZ] *= signx * signy * signz;
422 
423  fper.D2phi += fper.D3phi * off + 0.5 * ((fper.D4phi * off) * off) + (1.0 / 6) * (((fper.D5phi * off) * off) * off);
424 #endif
425 
426 #if(HIGHEST_NEEDEDORDER_EWALD_DPHI + EWALD_TAYLOR_ORDER) >= 6
427  fper.D6phi[pXXXXXY] *= signx * signy;
428  fper.D6phi[pXXXXXZ] *= signx * signz;
429  fper.D6phi[pXXXXYZ] *= signy * signz;
430  fper.D6phi[pXXXYYY] *= signx * signy;
431  fper.D6phi[pXXXYYZ] *= signx * signz;
432  fper.D6phi[pXXXYZZ] *= signx * signy;
433  fper.D6phi[pXXXZZZ] *= signx * signz;
434  fper.D6phi[pXXYYYZ] *= signy * signz;
435  fper.D6phi[pXXYZZZ] *= signy * signz;
436  fper.D6phi[pXYYYYY] *= signx * signy;
437  fper.D6phi[pXYYYYZ] *= signx * signz;
438  fper.D6phi[pXYYYZZ] *= signx * signy;
439  fper.D6phi[pXYYZZZ] *= signx * signz;
440  fper.D6phi[pXYZZZZ] *= signx * signy;
441  fper.D6phi[pXZZZZZ] *= signx * signz;
442  fper.D6phi[pYYYYYZ] *= signy * signz;
443  fper.D6phi[pYYYZZZ] *= signy * signz;
444  fper.D6phi[pYZZZZZ] *= signy * signz;
445 
446  fper.D3phi += fper.D4phi * off + 0.5 * ((fper.D5phi * off) * off) + (1.0 / 6) * (((fper.D6phi * off) * off) * off);
447 #endif
448 
449 #if(HIGHEST_NEEDEDORDER_EWALD_DPHI + EWALD_TAYLOR_ORDER) >= 7
450  fper.D7phi[tXXXXXXX] *= signx;
451  fper.D7phi[tXXXXXXY] *= signy;
452  fper.D7phi[tXXXXXXZ] *= signz;
453  fper.D7phi[tXXXXXYY] *= signx;
454  fper.D7phi[tXXXXXYZ] *= signx * signy * signz;
455  fper.D7phi[tXXXXXZZ] *= signx;
456  fper.D7phi[tXXXXYYY] *= signy;
457  fper.D7phi[tXXXXYYZ] *= signz;
458  fper.D7phi[tXXXXYZZ] *= signy;
459  fper.D7phi[tXXXXZZZ] *= signz;
460  fper.D7phi[tXXXYYYY] *= signx;
461  fper.D7phi[tXXXYYYZ] *= signx * signy * signz;
462  fper.D7phi[tXXXYYZZ] *= signx;
463  fper.D7phi[tXXXYZZZ] *= signx * signy * signz;
464  fper.D7phi[tXXXZZZZ] *= signx;
465  fper.D7phi[tXXYYYYY] *= signy;
466  fper.D7phi[tXXYYYYZ] *= signz;
467  fper.D7phi[tXXYYYZZ] *= signy;
468  fper.D7phi[tXXYYZZZ] *= signz;
469  fper.D7phi[tXXYZZZZ] *= signy;
470  fper.D7phi[tXXZZZZZ] *= signz;
471  fper.D7phi[tXYYYYYY] *= signx;
472  fper.D7phi[tXYYYYYZ] *= signx * signy * signz;
473  fper.D7phi[tXYYYYZZ] *= signx;
474  fper.D7phi[tXYYYZZZ] *= signx * signy * signz;
475  fper.D7phi[tXYYZZZZ] *= signx;
476  fper.D7phi[tXYZZZZZ] *= signx * signy * signz;
477  fper.D7phi[tXZZZZZZ] *= signx;
478  fper.D7phi[tYYYYYYY] *= signy;
479  fper.D7phi[tYYYYYYZ] *= signz;
480  fper.D7phi[tYYYYYZZ] *= signy;
481  fper.D7phi[tYYYYZZZ] *= signz;
482  fper.D7phi[tYYYZZZZ] *= signy;
483  fper.D7phi[tYYZZZZZ] *= signz;
484  fper.D7phi[tYZZZZZZ] *= signy;
485  fper.D7phi[tZZZZZZZ] *= signz;
486 
487  fper.D4phi += fper.D5phi * off + 0.5 * ((fper.D6phi * off) * off) + (1.0 / 6) * (((fper.D7phi * off) * off) * off);
488  fper.D5phi += fper.D6phi * off + 0.5 * ((fper.D7phi * off) * off);
489 #endif
490 
491 #else
492 
493  // only second order Taylor expansion, i.e. EWALD_TAYLOR_ORDER==2
494 
495  // now Taylor corrections
496 
497 #ifndef GRAVITY_TALLBOX
498  fper.D0phi += fper.D1phi * off + 0.5 * ((fper.D2phi * off) * off);
499  fper.D1phi += fper.D2phi * off + 0.5 * ((fper.D3phi * off) * off);
500 #endif
501 
502  if(flag == POINTMASS)
503  return;
504 
505 #if(HIGHEST_NEEDEDORDER_EWALD_DPHI + EWALD_TAYLOR_ORDER) >= 4
506  fper.D4phi[sXXXY] *= signx * signy;
507  fper.D4phi[sXYYY] *= signx * signy;
508  fper.D4phi[sXXXZ] *= signx * signz;
509  fper.D4phi[sXZZZ] *= signx * signz;
510  fper.D4phi[sYYYZ] *= signy * signz;
511  fper.D4phi[sYZZZ] *= signy * signz;
512  fper.D4phi[sXXYZ] *= signy * signz;
513  fper.D4phi[sXYYZ] *= signx * signz;
514  fper.D4phi[sXYZZ] *= signx * signy;
515 
516  fper.D2phi += fper.D3phi * off + 0.5 * ((fper.D4phi * off) * off);
517 #endif
518 
519 #if(HIGHEST_NEEDEDORDER_EWALD_DPHI + EWALD_TAYLOR_ORDER) >= 5
520  fper.D5phi[rXXXXX] *= signx;
521  fper.D5phi[rYYYYY] *= signy;
522  fper.D5phi[rZZZZZ] *= signz;
523 
524  fper.D5phi[rXXXXY] *= signy;
525  fper.D5phi[rXXXXZ] *= signz;
526  fper.D5phi[rXYYYY] *= signx;
527  fper.D5phi[rXZZZZ] *= signx;
528  fper.D5phi[rYYYYZ] *= signz;
529  fper.D5phi[rYZZZZ] *= signy;
530 
531  fper.D5phi[rXXXYY] *= signx;
532  fper.D5phi[rXXXZZ] *= signx;
533  fper.D5phi[rXXYYY] *= signy;
534  fper.D5phi[rXXZZZ] *= signz;
535  fper.D5phi[rYYYZZ] *= signy;
536  fper.D5phi[rYYZZZ] *= signz;
537 
538  fper.D5phi[rXXYZZ] *= signy;
539  fper.D5phi[rXXYYZ] *= signz;
540  fper.D5phi[rXYYZZ] *= signx;
541 
542  fper.D5phi[rXXXYZ] *= signx * signy * signz;
543  fper.D5phi[rXYYYZ] *= signx * signy * signz;
544  fper.D5phi[rXYZZZ] *= signx * signy * signz;
545 
546  fper.D3phi += fper.D4phi * off + 0.5 * ((fper.D5phi * off) * off);
547 #endif
548 
549 #if(HIGHEST_NEEDEDORDER_EWALD_DPHI + EWALD_TAYLOR_ORDER) >= 6
550  fper.D6phi[pXXXXXY] *= signx * signy;
551  fper.D6phi[pXXXXXZ] *= signx * signz;
552  fper.D6phi[pXXXXYZ] *= signy * signz;
553  fper.D6phi[pXXXYYY] *= signx * signy;
554  fper.D6phi[pXXXYYZ] *= signx * signz;
555  fper.D6phi[pXXXYZZ] *= signx * signy;
556  fper.D6phi[pXXXZZZ] *= signx * signz;
557  fper.D6phi[pXXYYYZ] *= signy * signz;
558  fper.D6phi[pXXYZZZ] *= signy * signz;
559  fper.D6phi[pXYYYYY] *= signx * signy;
560  fper.D6phi[pXYYYYZ] *= signx * signz;
561  fper.D6phi[pXYYYZZ] *= signx * signy;
562  fper.D6phi[pXYYZZZ] *= signx * signz;
563  fper.D6phi[pXYZZZZ] *= signx * signy;
564  fper.D6phi[pXZZZZZ] *= signx * signz;
565  fper.D6phi[pYYYYYZ] *= signy * signz;
566  fper.D6phi[pYYYZZZ] *= signy * signz;
567  fper.D6phi[pYZZZZZ] *= signy * signz;
568 
569  fper.D4phi += fper.D5phi * off + 0.5 * ((fper.D6phi * off) * off);
570 #endif
571 
572 #if(HIGHEST_NEEDEDORDER_EWALD_DPHI + EWALD_TAYLOR_ORDER) >= 7
573  fper.D7phi[tXXXXXXX] *= signx;
574  fper.D7phi[tXXXXXXY] *= signy;
575  fper.D7phi[tXXXXXXZ] *= signz;
576  fper.D7phi[tXXXXXYY] *= signx;
577  fper.D7phi[tXXXXXYZ] *= signx * signy * signz;
578  fper.D7phi[tXXXXXZZ] *= signx;
579  fper.D7phi[tXXXXYYY] *= signy;
580  fper.D7phi[tXXXXYYZ] *= signz;
581  fper.D7phi[tXXXXYZZ] *= signy;
582  fper.D7phi[tXXXXZZZ] *= signz;
583  fper.D7phi[tXXXYYYY] *= signx;
584  fper.D7phi[tXXXYYYZ] *= signx * signy * signz;
585  fper.D7phi[tXXXYYZZ] *= signx;
586  fper.D7phi[tXXXYZZZ] *= signx * signy * signz;
587  fper.D7phi[tXXXZZZZ] *= signx;
588  fper.D7phi[tXXYYYYY] *= signy;
589  fper.D7phi[tXXYYYYZ] *= signz;
590  fper.D7phi[tXXYYYZZ] *= signy;
591  fper.D7phi[tXXYYZZZ] *= signz;
592  fper.D7phi[tXXYZZZZ] *= signy;
593  fper.D7phi[tXXZZZZZ] *= signz;
594  fper.D7phi[tXYYYYYY] *= signx;
595  fper.D7phi[tXYYYYYZ] *= signx * signy * signz;
596  fper.D7phi[tXYYYYZZ] *= signx;
597  fper.D7phi[tXYYYZZZ] *= signx * signy * signz;
598  fper.D7phi[tXYYZZZZ] *= signx;
599  fper.D7phi[tXYZZZZZ] *= signx * signy * signz;
600  fper.D7phi[tXZZZZZZ] *= signx;
601  fper.D7phi[tYYYYYYY] *= signy;
602  fper.D7phi[tYYYYYYZ] *= signz;
603  fper.D7phi[tYYYYYZZ] *= signy;
604  fper.D7phi[tYYYYZZZ] *= signz;
605  fper.D7phi[tYYYZZZZ] *= signy;
606  fper.D7phi[tYYZZZZZ] *= signz;
607  fper.D7phi[tYZZZZZZ] *= signy;
608  fper.D7phi[tZZZZZZZ] *= signz;
609 
610  fper.D5phi += fper.D6phi * off + 0.5 * ((fper.D7phi * off) * off);
611 #endif
612 
613 #endif
614 }
615 
623 double ewald::ewald_D0(double x, double y, double z)
624 {
625  static int printed = 0;
626 
627  double D0 = 0.0;
628 
629 #ifndef GRAVITY_TALLBOX
630 
631  double leff = pow((1.0 / LONG_X) * (1.0 / LONG_Y) * (1.0 / LONG_Z), 1.0 / 3);
632  double alpha = 2.0 / leff;
633  double alpha2 = alpha * alpha;
634 
635  int qxmax = (int)(8.0 * LONG_X / alpha + 0.5);
636  int qymax = (int)(8.0 * LONG_Y / alpha + 0.5);
637  int qzmax = (int)(8.0 * LONG_Z / alpha + 0.5);
638 
639  int nxmax = (int)(2.0 * alpha / LONG_X + 0.5);
640  int nymax = (int)(2.0 * alpha / LONG_Y + 0.5);
641  int nzmax = (int)(2.0 * alpha / LONG_Z + 0.5);
642 
643  if(printed == 0)
644  {
645  mpi_printf("EWALD: D0 table: qxmax=%d qymax=%d qzmax=%d nxmax=%d nymax=%d nzmax=%d\n", qxmax, qymax, qzmax, nxmax, nymax,
646  nzmax);
647  printed = 1;
648  }
649 
650  for(int nx = -qxmax; nx <= qxmax; nx++)
651  for(int ny = -qymax; ny <= qymax; ny++)
652  for(int nz = -qzmax; nz <= qzmax; nz++)
653  {
654  double dx = x - nx * (1.0 / LONG_X);
655  double dy = y - ny * (1.0 / LONG_Y);
656  double dz = z - nz * (1.0 / LONG_Z);
657 
658  double r2 = dx * dx + dy * dy + dz * dz;
659  double r = sqrt(r2);
660 
661  double rinv = (r > 0) ? 1.0 / r : 0.0;
662 
663  double g0;
664 
665  if(nx != 0 || ny != 0 || nz != 0)
666  {
667  g0 = -erfc(alpha * r) * rinv;
668  }
669  else
670  {
671  /* we add the 1/r term here to the (0|0|0) entry, followed by differentiation, and the limit r->0 to obtain accurate
672  * results at the origin
673  */
674 
675  /* for small r:
676  *
677  * [1- erfc(a r)]/r = 2 a/sqrt(pi) * [ 1 - (a r)^2/3 + (a r)^4 / 10 - (a r)^6 / 42 + (a r)^8 / 216 - ...]
678  */
679 
680  if((alpha * r) < 0.5)
681  {
682  g0 = 2.0 * pow(alpha, 1) / sqrt(M_PI) *
683  (1.0 - pow(alpha * r, 2) / 3.0 + pow(alpha * r, 4) / 10.0 - pow(alpha * r, 6) / 42.0 +
684  pow(alpha * r, 8) / 216.0 - pow(alpha * r, 10) / 1320.0);
685  }
686  else
687  {
688  g0 = erf(alpha * r) * rinv;
689  }
690  }
691 
692  D0 += g0;
693  }
694 
695  for(int nx = -nxmax; nx <= nxmax; nx++)
696  for(int ny = -nymax; ny <= nymax; ny++)
697  for(int nz = -nzmax; nz <= nzmax; nz++)
698  {
699  if(nx != 0 || ny != 0 || nz != 0)
700  {
701  double kx = (2.0 * M_PI * LONG_X) * nx;
702  double ky = (2.0 * M_PI * LONG_Y) * ny;
703  double kz = (2.0 * M_PI * LONG_Z) * nz;
704  double k2 = kx * kx + ky * ky + kz * kz;
705  double kdotx = (x * kx + y * ky + z * kz);
706 
707  D0 += -4.0 * M_PI * (LONG_X * LONG_Y * LONG_Z) / k2 * exp(-k2 / (4.0 * alpha2)) * cos(kdotx);
708  }
709  }
710 
711  D0 += M_PI * (LONG_X * LONG_Y * LONG_Z) / (alpha * alpha);
712 
713 #else
714  /* in the tallbox case, the third dimension, z, is assumed to be the non-periodic one */
715 
716  double leff = sqrt(BOXX * BOXY);
717  double alpha = 2.0 / leff;
718 
719  int qxmax = (int)(8.0 / (BOXX * alpha) + 0.5);
720  int qymax = (int)(8.0 / (BOXY * alpha) + 0.5);
721 
722  int nxmax = (int)(2.0 * alpha * BOXX + 0.5);
723  int nymax = (int)(2.0 * alpha * BOXY + 0.5);
724 
725  if(printed == 0)
726  {
727  mpi_printf("EWALD: D0 table: qxmax=%d qymax=%d nxmax=%d nymax=%d\n", qxmax, qymax, nxmax, nymax);
728  printed = 1;
729  }
730 
731  for(int nx = -qxmax; nx <= qxmax; nx++)
732  for(int ny = -qymax; ny <= qymax; ny++)
733  {
734  double dx = x - nx * BOXX;
735  double dy = y - ny * BOXY;
736  double r = sqrt(dx * dx + dy * dy + z * z);
737 
738  double rinv = (r > 0) ? 1.0 / r : 0.0;
739 
740  double g0;
741 
742  if(nx != 0 || ny != 0)
743  {
744  g0 = -erfc(alpha * r) * rinv;
745  }
746  else
747  {
748  /* we add the 1/r term here to the (0|0) entry */
749 
750  if((alpha * r) < 0.5)
751  {
752  g0 = 2.0 * pow(alpha, 1) / sqrt(M_PI) *
753  (1.0 - pow(alpha * r, 2) / 3.0 + pow(alpha * r, 4) / 10.0 - pow(alpha * r, 6) / 42.0 + pow(alpha * r, 8) / 216.0 -
754  pow(alpha * r, 10) / 1320.0);
755  }
756  else
757  {
758  g0 = erf(alpha * r) * rinv;
759  }
760  }
761 
762  D0 += g0;
763  }
764 
765  double alpha2 = alpha * alpha;
766 
767  for(int nx = -nxmax; nx <= nxmax; nx++)
768  for(int ny = -nymax; ny <= nymax; ny++)
769  {
770  if(nx != 0 || ny != 0)
771  {
772  double kx = (2.0 * M_PI / BOXX) * nx;
773  double ky = (2.0 * M_PI / BOXY) * ny;
774  double k2 = kx * kx + ky * ky;
775  double k = sqrt(k2);
776 
777  if(k * z > 0)
778  {
779  double ex = exp(-k * z);
780  if(ex > 0)
781  D0 += -M_PI / (BOXX * BOXY) * (erfc(k / (2 * alpha) + alpha * z) / ex + ex * erfc(k / (2 * alpha) - alpha * z)) / k;
782  }
783  else
784  {
785  double ex = exp(k * z);
786  if(ex > 0)
787  D0 += -M_PI / (BOXX * BOXY) * cos(kx * x + ky * y) *
788  (ex * erfc(k / (2 * alpha) + alpha * z) + erfc(k / (2 * alpha) - alpha * z) / ex) / k;
789  }
790  }
791  }
792 
793  D0 += 2.0 * alpha / sqrt(M_PI) + 2 * sqrt(M_PI) / (BOXX * BOXY) * (exp(-alpha2 * z * z) / alpha + sqrt(M_PI) * z * erf(alpha * z));
794 
795 #endif
796 
797  return D0;
798 }
799 
807 vector<double> ewald::ewald_D1(double x, double y, double z)
808 {
809  static int printed = 0;
810 
811  vector<double> D1 = 0.0;
812 
813 #ifndef GRAVITY_TALLBOX
814 
815  double leff = pow((1.0 / LONG_X) * (1.0 / LONG_Y) * (1.0 / LONG_Z), 1.0 / 3);
816  double alpha = 2.0 / leff;
817  double alpha2 = alpha * alpha;
818 
819  int qxmax = (int)(8.0 * LONG_X / alpha + 0.5);
820  int qymax = (int)(8.0 * LONG_Y / alpha + 0.5);
821  int qzmax = (int)(8.0 * LONG_Z / alpha + 0.5);
822 
823  int nxmax = (int)(2.0 * alpha / LONG_X + 0.5);
824  int nymax = (int)(2.0 * alpha / LONG_Y + 0.5);
825  int nzmax = (int)(2.0 * alpha / LONG_Z + 0.5);
826 
827  if(printed == 0)
828  {
829  mpi_printf("EWALD: D1 table: qxmax=%d qymax=%d qzmax=%d nxmax=%d nymax=%d nzmax=%d\n", qxmax, qymax, qzmax, nxmax, nymax,
830  nzmax);
831  printed = 1;
832  }
833 
834  for(int nx = -qxmax; nx <= qxmax; nx++)
835  for(int ny = -qymax; ny <= qymax; ny++)
836  for(int nz = -qzmax; nz <= qzmax; nz++)
837  {
838  double dx = x - nx * (1.0 / LONG_X);
839  double dy = y - ny * (1.0 / LONG_Y);
840  double dz = z - nz * (1.0 / LONG_Z);
841 
842  vector<double> dxyz(dx, dy, dz);
843 
844  double r2 = dx * dx + dy * dy + dz * dz;
845  double r = sqrt(r2);
846 
847  double rinv = (r > 0) ? 1.0 / r : 0.0;
848  double r2inv = rinv * rinv;
849  double r3inv = r2inv * rinv;
850 
851  double g1;
852 
853  if(nx != 0 || ny != 0 || nz != 0)
854  {
855  g1 = (erfc(alpha * r) + 2.0 * alpha * r / sqrt(M_PI) * exp(-alpha2 * r2)) * r3inv;
856  }
857  else
858  {
859  /* we add the 1/r term here to the (0|0|0) entry, followed by differentiation, and the limit r->0 to obtain accurate
860  * results at the origin
861  */
862 
863  /* Note, for small r:
864  *
865  * [1/- erfc(a r)]/r = 2 a/sqrt(pi) * [ 1 - (a r)^2/3 + (a r)^4 / 10 - (a r)^6 / 42 + (a r)^8 / 216 - ...]
866  *
867  * Hence for r = 0:
868  *
869  * g0 = 2 * alpha / sqrt(pi)
870  * g1 = -4/3 * alpha^3 / sqrt(pi)
871  * g2 = 8/5 * alpha^5 / sqrt(pi)
872  * g3 = -16/7 * alpha^7 / sqrt(pi)
873  * g4 = 32/9 * alpha^9 / sqrt(pi)
874  * g5 = -64/11 * alpha^11/ sqrt(pi)
875  */
876 
877  if((alpha * r) < 0.5)
878  {
879  g1 = 4.0 * pow(alpha, 3) / sqrt(M_PI) *
880  (-1.0 / 3.0 + pow(alpha * r, 2) / 5.0 - pow(alpha * r, 4) / 14.0 + pow(alpha * r, 6) / 54.0 -
881  pow(alpha * r, 8) / 264.0 + pow(alpha * r, 10) / 1560.0);
882  }
883  else
884  {
885  g1 = (-erf(alpha * r) + 2.0 * alpha * r / sqrt(M_PI) * exp(-alpha2 * r2)) * r3inv;
886  }
887  }
888 
889  D1 += g1 * dxyz;
890  }
891 
892  for(int nx = -nxmax; nx <= nxmax; nx++)
893  for(int ny = -nymax; ny <= nymax; ny++)
894  for(int nz = -nzmax; nz <= nzmax; nz++)
895  {
896  double kx = (2.0 * M_PI * LONG_X) * nx;
897  double ky = (2.0 * M_PI * LONG_Y) * ny;
898  double kz = (2.0 * M_PI * LONG_Z) * nz;
899  double k2 = kx * kx + ky * ky + kz * kz;
900 
901  if(k2 > 0)
902  {
903  double kdotx = (x * kx + y * ky + z * kz);
904  double val = 4.0 * M_PI * (LONG_X * LONG_Y * LONG_Z) / k2 * exp(-k2 / (4.0 * alpha2)) * sin(kdotx);
905  D1[0] += kx * val;
906  D1[1] += ky * val;
907  D1[2] += kz * val;
908  }
909  }
910 #else
911  /* this is the case with periodicity only in two dimensions */
912 
913  double leff = sqrt(BOXX * BOXY);
914  double alpha = 2.0 / leff;
915  double alpha2 = alpha * alpha;
916 
917  int qxmax = (int)(8.0 / (BOXX * alpha) + 0.5);
918  int qymax = (int)(8.0 / (BOXY * alpha) + 0.5);
919 
920  int nxmax = (int)(2.0 * alpha * BOXX + 0.5);
921  int nymax = (int)(2.0 * alpha * BOXY + 0.5);
922 
923  if(printed == 0)
924  {
925  mpi_printf("EWALD: D1 table: qxmax=%d qymax=%d nxmax=%d nymax=%d\n", qxmax, qymax, nxmax, nymax);
926  printed = 1;
927  }
928 
929  for(int nx = -qxmax; nx <= qxmax; nx++)
930  for(int ny = -qymax; ny <= qymax; ny++)
931  {
932  double dx = x - nx * BOXX;
933  double dy = y - ny * BOXY;
934  double dz = z;
935 
936  vector<double> dxyz(dx, dy, dz);
937 
938  double r2 = dx * dx + dy * dy + dz * dz;
939  double r = sqrt(r2);
940 
941  double rinv = (r > 0) ? 1.0 / r : 0.0;
942  double r2inv = rinv * rinv;
943  double r3inv = r2inv * rinv;
944 
945  double g1;
946 
947  if(nx != 0 || ny != 0)
948  {
949  g1 = (erfc(alpha * r) + 2.0 * alpha * r / sqrt(M_PI) * exp(-alpha2 * r2)) * r3inv;
950  }
951  else
952  {
953  /* we add the 1/r term here to the (0|0) entry */
954 
955  if((alpha * r) < 0.5)
956  {
957  g1 = 4.0 * pow(alpha, 3) / sqrt(M_PI) *
958  (-1.0 / 3.0 + pow(alpha * r, 2) / 5.0 - pow(alpha * r, 4) / 14.0 + pow(alpha * r, 6) / 54.0 -
959  pow(alpha * r, 8) / 264.0 + pow(alpha * r, 10) / 1560.0);
960  }
961  else
962  {
963  g1 = (-erf(alpha * r) + 2.0 * alpha * r / sqrt(M_PI) * exp(-alpha2 * r2)) * r3inv;
964  }
965  }
966 
967  D1 += g1 * dxyz;
968  }
969 
970  for(int nx = -nxmax; nx <= nxmax; nx++)
971  for(int ny = -nymax; ny <= nymax; ny++)
972  {
973  if(nx != 0 || ny != 0)
974  {
975  double kx = (2.0 * M_PI / BOXX) * nx;
976  double ky = (2.0 * M_PI / BOXY) * ny;
977  double k2 = kx * kx + ky * ky;
978  double k = sqrt(k2);
979 
980  double val = M_PI / (BOXX * BOXY) * sin(kx * x + ky * y) *
981  (exp(k * z) * erfc(k / (2 * alpha) + alpha * z) + exp(-k * z) * erfc(k / (2 * alpha) - alpha * z)) / k;
982 
983  D1[0] -= -kx * val;
984  D1[1] -= -ky * val;
985  D1[2] -= M_PI / (BOXX * BOXY) * cos(kx * x + ky * y) *
986  (exp(k * z) * erfc(k / (2 * alpha) + alpha * z) - exp(-k * z) * erfc(k / (2 * alpha) - alpha * z));
987  }
988  }
989 
990  D1[2] += 2.0 * M_PI / (BOXX * BOXY) * erf(alpha * z);
991 #endif
992 
993  return D1; // now in dimensionless form;
994 }
995 
996 symtensor2<double> ewald::ewald_D2(double x, double y, double z)
997 {
998  static int printed = 0;
999 
1000  symtensor2<double> D2 = 0.0;
1001 
1002  double leff = pow((1.0 / LONG_X) * (1.0 / LONG_Y) * (1.0 / LONG_Z), 1.0 / 3);
1003  double alpha = 2.0 / leff;
1004  double alpha2 = alpha * alpha;
1005 
1006  int qxmax = (int)(8.0 * LONG_X / alpha + 0.5);
1007  int qymax = (int)(8.0 * LONG_Y / alpha + 0.5);
1008  int qzmax = (int)(8.0 * LONG_Z / alpha + 0.5);
1009 
1010  int nxmax = (int)(2.0 * alpha / LONG_X + 0.5);
1011  int nymax = (int)(2.0 * alpha / LONG_Y + 0.5);
1012  int nzmax = (int)(2.0 * alpha / LONG_Z + 0.5);
1013 
1014  if(printed == 0)
1015  {
1016  mpi_printf("EWALD: D2 table: qxmax=%d qymax=%d qzmax=%d nxmax=%d nymax=%d nzmax=%d\n", qxmax, qymax, qzmax, nxmax, nymax,
1017  nzmax);
1018  printed = 1;
1019  }
1020 
1021  for(int nx = -qxmax; nx <= qxmax; nx++)
1022  for(int ny = -qymax; ny <= qymax; ny++)
1023  for(int nz = -qzmax; nz <= qzmax; nz++)
1024  {
1025  double dx = x - nx * (1.0 / LONG_X);
1026  double dy = y - ny * (1.0 / LONG_Y);
1027  double dz = z - nz * (1.0 / LONG_Z);
1028 
1029  vector<double> dxyz(dx, dy, dz);
1030 
1031  double r2 = dx * dx + dy * dy + dz * dz;
1032  double r = sqrt(r2);
1033 
1034  double rinv = (r > 0) ? 1.0 / r : 0.0;
1035  double r2inv = rinv * rinv;
1036  double r3inv = r2inv * rinv;
1037  double r5inv = r3inv * r2inv;
1038 
1039  double g1, g2;
1040 
1041  if(nx != 0 || ny != 0 || nz != 0)
1042  {
1043  g1 = (erfc(alpha * r) + 2.0 * alpha * r / sqrt(M_PI) * exp(-alpha2 * r2)) * r3inv;
1044 
1045  g2 = -(3.0 * erfc(alpha * r) + (6.0 * alpha * r + 4.0 * pow(alpha * r, 3)) / sqrt(M_PI) * exp(-alpha2 * r2)) * r5inv;
1046  }
1047  else
1048  {
1049  /* we add the 1/r term here to the (0|0|0) entry, followed by differentiation, and the limit r->0 to obtain accurate
1050  * results at the origin
1051  */
1052 
1053  /* Note, for small r:
1054  *
1055  * [1/- erfc(a r)]/r = 2 a/sqrt(pi) * [ 1 - (a r)^2/3 + (a r)^4 / 10 - (a r)^6 / 42 + (a r)^8 / 216 - ...]
1056  *
1057  * Hence for r = 0:
1058  *
1059  * g0 = 2 * alpha / sqrt(pi)
1060  * g1 = -4/3 * alpha^3 / sqrt(pi)
1061  * g2 = 8/5 * alpha^5 / sqrt(pi)
1062  * g3 = -16/7 * alpha^7 / sqrt(pi)
1063  * g4 = 32/9 * alpha^9 / sqrt(pi)
1064  * g5 = -64/11 * alpha^11/ sqrt(pi)
1065  */
1066 
1067  if((alpha * r) < 0.5)
1068  {
1069  g1 = 4.0 * pow(alpha, 3) / sqrt(M_PI) *
1070  (-1.0 / 3.0 + pow(alpha * r, 2) / 5.0 - pow(alpha * r, 4) / 14.0 + pow(alpha * r, 6) / 54.0 -
1071  pow(alpha * r, 8) / 264.0 + pow(alpha * r, 10) / 1560.0);
1072 
1073  g2 = 8.0 * pow(alpha, 5) / sqrt(M_PI) *
1074  (1.0 / 5.0 - pow(alpha * r, 2) / 7.0 + pow(alpha * r, 4) / 18.0 - pow(alpha * r, 6) / 66.0 +
1075  pow(alpha * r, 8) / 312.0 - pow(alpha * r, 10) / 1800.0);
1076  }
1077  else
1078  {
1079  g1 = (-erf(alpha * r) + 2.0 * alpha * r / sqrt(M_PI) * exp(-alpha2 * r2)) * r3inv;
1080 
1081  g2 = (3.0 * erf(alpha * r) - (6.0 * alpha * r + 4.0 * pow(alpha * r, 3)) / sqrt(M_PI) * exp(-alpha2 * r2)) * r5inv;
1082  }
1083  }
1084 
1085  D2 += g2 * (dxyz % dxyz);
1086  D2[qXX] += g1;
1087  D2[qYY] += g1;
1088  D2[qZZ] += g1;
1089  }
1090 
1091  for(int nx = -nxmax; nx <= nxmax; nx++)
1092  for(int ny = -nymax; ny <= nymax; ny++)
1093  for(int nz = -nzmax; nz <= nzmax; nz++)
1094  {
1095  if(nx != 0 || ny != 0 || nz != 0)
1096  {
1097  double kx = (2.0 * M_PI * LONG_X) * nx;
1098  double ky = (2.0 * M_PI * LONG_Y) * ny;
1099  double kz = (2.0 * M_PI * LONG_Z) * nz;
1100  double k2 = kx * kx + ky * ky + kz * kz;
1101 
1102  double kdotx = (x * kx + y * ky + z * kz);
1103  double val = 4.0 * M_PI * (LONG_X * LONG_Y * LONG_Z) / k2 * exp(-k2 / (4.0 * alpha2)) * cos(kdotx);
1104 
1105  vector<double> kxyz(kx, ky, kz);
1106 
1107  D2 += (val * kxyz) % kxyz;
1108  }
1109  }
1110 
1111  return D2;
1112 }
1113 
1114 symtensor3<double> ewald::ewald_D3(double x, double y, double z)
1115 {
1116  static int printed = 0;
1117 
1118 #ifdef GRAVITY_TALLBOX
1119  Terminate("GRAVITY_TALLBOX is not implemented for MULTIPOLE_ORDER >= 3");
1120 #endif
1121 
1122  symtensor3<double> D3 = 0.0;
1123 
1124  double leff = pow((1.0 / LONG_X) * (1.0 / LONG_Y) * (1.0 / LONG_Z), 1.0 / 3);
1125  double alpha = 2.0 / leff;
1126  double alpha2 = alpha * alpha;
1127 
1128  int qxmax = (int)(8.0 * LONG_X / alpha + 0.5);
1129  int qymax = (int)(8.0 * LONG_Y / alpha + 0.5);
1130  int qzmax = (int)(8.0 * LONG_Z / alpha + 0.5);
1131 
1132  int nxmax = (int)(2.0 * alpha / LONG_X + 0.5);
1133  int nymax = (int)(2.0 * alpha / LONG_Y + 0.5);
1134  int nzmax = (int)(2.0 * alpha / LONG_Z + 0.5);
1135 
1136  if(printed == 0)
1137  {
1138  mpi_printf("EWALD: D3 table: qxmax=%d qymax=%d qzmax=%d nxmax=%d nymax=%d nzmax=%d\n", qxmax, qymax, qzmax, nxmax, nymax,
1139  nzmax);
1140  printed = 1;
1141  }
1142 
1143  for(int nx = -qxmax; nx <= qxmax; nx++)
1144  for(int ny = -qymax; ny <= qymax; ny++)
1145  for(int nz = -qzmax; nz <= qzmax; nz++)
1146  {
1147  double dx = x - nx * (1.0 / LONG_X);
1148  double dy = y - ny * (1.0 / LONG_Y);
1149  double dz = z - nz * (1.0 / LONG_Z);
1150 
1151  vector<double> dxyz(dx, dy, dz);
1152 
1153  double r2 = dx * dx + dy * dy + dz * dz;
1154  double r = sqrt(r2);
1155 
1156  double rinv = (r > 0) ? 1.0 / r : 0.0;
1157  double r2inv = rinv * rinv;
1158  double r3inv = r2inv * rinv;
1159  double r4inv = r2inv * r2inv;
1160  double r5inv = r2inv * r3inv;
1161  double r7inv = r3inv * r4inv;
1162 
1163  double g2, g3;
1164 
1165  if(nx != 0 || ny != 0 || nz != 0)
1166  {
1167  g2 = -(3.0 * erfc(alpha * r) + (6.0 * alpha * r + 4.0 * pow(alpha * r, 3)) / sqrt(M_PI) * exp(-alpha2 * r2)) * r5inv;
1168 
1169  g3 = (15.0 * erfc(alpha * r) +
1170  (30.0 * alpha * r + 20.0 * pow(alpha * r, 3) + 8.0 * pow(alpha * r, 5)) / sqrt(M_PI) * exp(-alpha2 * r2)) *
1171  r7inv;
1172  }
1173  else
1174  {
1175  /* we add the 1/r term here to the (0|0|0) entry, followed by differentiation, and the limit r->0 to obtain accurate
1176  * results at the origin
1177  */
1178 
1179  /* Note, for small r:
1180  *
1181  * [1/- erfc(a r)]/r = 2 a/sqrt(pi) * [ 1 - (a r)^2/3 + (a r)^4 / 10 - (a r)^6 / 42 + (a r)^8 / 216 - ...]
1182  *
1183  * Hence for r = 0:
1184  *
1185  * g0 = 2 * alpha / sqrt(pi)
1186  * g1 = -4/3 * alpha^3 / sqrt(pi)
1187  * g2 = 8/5 * alpha^5 / sqrt(pi)
1188  * g3 = -16/7 * alpha^7 / sqrt(pi)
1189  * g4 = 32/9 * alpha^9 / sqrt(pi)
1190  * g5 = -64/11 * alpha^11/ sqrt(pi)
1191  */
1192  if((alpha * r) < 0.5)
1193  {
1194  g2 = 8.0 * pow(alpha, 5) / sqrt(M_PI) *
1195  (1.0 / 5.0 - pow(alpha * r, 2) / 7.0 + pow(alpha * r, 4) / 18.0 - pow(alpha * r, 6) / 66.0 +
1196  pow(alpha * r, 8) / 312.0 - pow(alpha * r, 10) / 1800.0);
1197 
1198  g3 = 16.0 * pow(alpha, 7) / sqrt(M_PI) *
1199  (-1.0 / 7.0 + pow(alpha * r, 2) / 9.0 - pow(alpha * r, 4) / 22.0 + pow(alpha * r, 6) / 78.0 -
1200  pow(alpha * r, 8) / 360.0 + pow(alpha * r, 10) / 2040.0);
1201  }
1202  else
1203  {
1204  g2 = (3.0 * erf(alpha * r) - (6.0 * alpha * r + 4.0 * pow(alpha * r, 3)) / sqrt(M_PI) * exp(-alpha2 * r2)) * r5inv;
1205 
1206  g3 = (-15.0 * erf(alpha * r) +
1207  (30.0 * alpha * r + 20.0 * pow(alpha * r, 3) + 8.0 * pow(alpha * r, 5)) / sqrt(M_PI) * exp(-alpha2 * r2)) *
1208  r7inv;
1209  }
1210  }
1211 
1212  symtensor2<double> aux2 = dxyz % dxyz;
1213  symtensor3<double> aux3;
1214 
1215  setup_D3(ADD, D3, dxyz, aux2, aux3, g2, g3);
1216  }
1217 
1218  for(int nx = -nxmax; nx <= nxmax; nx++)
1219  for(int ny = -nymax; ny <= nymax; ny++)
1220  for(int nz = -nzmax; nz <= nzmax; nz++)
1221  {
1222  if(nx != 0 || ny != 0 || nz != 0)
1223  {
1224  double kx = (2.0 * M_PI * LONG_X) * nx;
1225  double ky = (2.0 * M_PI * LONG_Y) * ny;
1226  double kz = (2.0 * M_PI * LONG_Z) * nz;
1227  double k2 = kx * kx + ky * ky + kz * kz;
1228 
1229  double kdotx = (x * kx + y * ky + z * kz);
1230  double val = -4.0 * M_PI * (LONG_X * LONG_Y * LONG_Z) / k2 * exp(-k2 / (4.0 * alpha2)) * sin(kdotx);
1231 
1232  vector<double> kxyz(kx, ky, kz);
1233 
1234  D3 += (val * kxyz) % (kxyz % kxyz);
1235  }
1236  }
1237 
1238  return D3;
1239 }
1240 
1241 symtensor4<double> ewald::ewald_D4(double x, double y, double z)
1242 {
1243  static int printed = 0;
1244 
1245 #ifdef GRAVITY_TALLBOX
1246  Terminate("GRAVITY_TALLBOX is not implemented for MULTIPOLE_ORDER >= 4");
1247 #endif
1248 
1249  symtensor4<double> D4 = 0.0;
1250 
1251  double leff = pow((1.0 / LONG_X) * (1.0 / LONG_Y) * (1.0 / LONG_Z), 1.0 / 3);
1252  double alpha = 2.0 / leff;
1253  double alpha2 = alpha * alpha;
1254 
1255  int qxmax = (int)(8.0 * LONG_X / alpha + 0.5);
1256  int qymax = (int)(8.0 * LONG_Y / alpha + 0.5);
1257  int qzmax = (int)(8.0 * LONG_Z / alpha + 0.5);
1258 
1259  int nxmax = (int)(2.0 * alpha / LONG_X + 0.5);
1260  int nymax = (int)(2.0 * alpha / LONG_Y + 0.5);
1261  int nzmax = (int)(2.0 * alpha / LONG_Z + 0.5);
1262 
1263  if(printed == 0)
1264  {
1265  mpi_printf("EWALD: D4 table: qxmax=%d qymax=%d qzmax=%d nxmax=%d nymax=%d nzmax=%d\n", qxmax, qymax, qzmax, nxmax, nymax,
1266  nzmax);
1267  printed = 1;
1268  }
1269 
1270  for(int nx = -qxmax; nx <= qxmax; nx++)
1271  for(int ny = -qymax; ny <= qymax; ny++)
1272  for(int nz = -qzmax; nz <= qzmax; nz++)
1273  {
1274  double dx = x - nx * (1.0 / LONG_X);
1275  double dy = y - ny * (1.0 / LONG_Y);
1276  double dz = z - nz * (1.0 / LONG_Z);
1277 
1278  vector<double> dxyz(dx, dy, dz);
1279 
1280  double r2 = dx * dx + dy * dy + dz * dz;
1281  double r = sqrt(r2);
1282 
1283  double rinv = (r > 0) ? 1.0 / r : 0.0;
1284  double r2inv = rinv * rinv;
1285  double r3inv = r2inv * rinv;
1286  double r4inv = r2inv * r2inv;
1287  double r5inv = r2inv * r3inv;
1288  double r7inv = r3inv * r4inv;
1289  double r9inv = r4inv * r5inv;
1290 
1291  double g2, g3, g4;
1292 
1293  if(nx != 0 || ny != 0 || nz != 0)
1294  {
1295  g2 = -(3.0 * erfc(alpha * r) + (6.0 * alpha * r + 4.0 * pow(alpha * r, 3)) / sqrt(M_PI) * exp(-alpha2 * r2)) * r5inv;
1296 
1297  g3 = (15.0 * erfc(alpha * r) +
1298  (30.0 * alpha * r + 20.0 * pow(alpha * r, 3) + 8.0 * pow(alpha * r, 5)) / sqrt(M_PI) * exp(-alpha2 * r2)) *
1299  r7inv;
1300 
1301  g4 = -(105.0 * erfc(alpha * r) +
1302  (210.0 * alpha * r + 140.0 * pow(alpha * r, 3) + 56.0 * pow(alpha * r, 5) + 16.0 * pow(alpha * r, 7)) /
1303  sqrt(M_PI) * exp(-alpha2 * r2)) *
1304  r9inv;
1305  }
1306  else
1307  {
1308  /* we add the 1/r term here to the (0|0|0) entry, followed by differentiation, and the limit r->0 to obtain accurate
1309  * results at the origin
1310  */
1311 
1312  /* Note, for small r:
1313  *
1314  * [1/- erfc(a r)]/r = 2 a/sqrt(pi) * [ 1 - (a r)^2/3 + (a r)^4 / 10 - (a r)^6 / 42 + (a r)^8 / 216 - ...]
1315  *
1316  * Hence for r = 0:
1317  *
1318  * g0 = 2 * alpha / sqrt(pi)
1319  * g1 = -4/3 * alpha^3 / sqrt(pi)
1320  * g2 = 8/5 * alpha^5 / sqrt(pi)
1321  * g3 = -16/7 * alpha^7 / sqrt(pi)
1322  * g4 = 32/9 * alpha^9 / sqrt(pi)
1323  * g5 = -64/11 * alpha^11/ sqrt(pi)
1324  */
1325 
1326  if((alpha * r) < 0.5)
1327  {
1328  g2 = 8.0 * pow(alpha, 5) / sqrt(M_PI) *
1329  (1.0 / 5.0 - pow(alpha * r, 2) / 7.0 + pow(alpha * r, 4) / 18.0 - pow(alpha * r, 6) / 66.0 +
1330  pow(alpha * r, 8) / 312.0 - pow(alpha * r, 10) / 1800.0);
1331 
1332  g3 = 16.0 * pow(alpha, 7) / sqrt(M_PI) *
1333  (-1.0 / 7.0 + pow(alpha * r, 2) / 9.0 - pow(alpha * r, 4) / 22.0 + pow(alpha * r, 6) / 78.0 -
1334  pow(alpha * r, 8) / 360.0 + pow(alpha * r, 10) / 2040.0);
1335 
1336  g4 = 32.0 * pow(alpha, 9) / sqrt(M_PI) *
1337  (1.0 / 9.0 - pow(alpha * r, 2) / 11.0 + pow(alpha * r, 4) / 26.0 - pow(alpha * r, 6) / 90.0 +
1338  pow(alpha * r, 8) / 408.0 - pow(alpha * r, 10) / 2280.0);
1339  }
1340  else
1341  {
1342  g2 = (3.0 * erf(alpha * r) - (6.0 * alpha * r + 4.0 * pow(alpha * r, 3)) / sqrt(M_PI) * exp(-alpha2 * r2)) * r5inv;
1343 
1344  g3 = (-15.0 * erf(alpha * r) +
1345  (30.0 * alpha * r + 20.0 * pow(alpha * r, 3) + 8.0 * pow(alpha * r, 5)) / sqrt(M_PI) * exp(-alpha2 * r2)) *
1346  r7inv;
1347 
1348  g4 = (105.0 * erf(alpha * r) -
1349  (210.0 * alpha * r + 140.0 * pow(alpha * r, 3) + 56.0 * pow(alpha * r, 5) + 16.0 * pow(alpha * r, 7)) /
1350  sqrt(M_PI) * exp(-alpha2 * r2)) *
1351  r9inv;
1352  }
1353  }
1354 
1355  symtensor2<double> aux2 = dxyz % dxyz;
1356  symtensor3<double> aux3 = dxyz % aux2;
1357  symtensor4<double> aux4;
1358 
1359  setup_D4(ADD, D4, dxyz, aux2, aux3, aux4, g2, g3, g4);
1360  }
1361 
1362  for(int nx = -nxmax; nx <= nxmax; nx++)
1363  for(int ny = -nymax; ny <= nymax; ny++)
1364  for(int nz = -nzmax; nz <= nzmax; nz++)
1365  {
1366  if(nx != 0 || ny != 0 || nz != 0)
1367  {
1368  double kx = (2.0 * M_PI * LONG_X) * nx;
1369  double ky = (2.0 * M_PI * LONG_Y) * ny;
1370  double kz = (2.0 * M_PI * LONG_Z) * nz;
1371  double k2 = kx * kx + ky * ky + kz * kz;
1372 
1373  double kdotx = (x * kx + y * ky + z * kz);
1374  double val = -4.0 * M_PI * (LONG_X * LONG_Y * LONG_Z) / k2 * exp(-k2 / (4.0 * alpha2)) * cos(kdotx);
1375 
1376  vector<double> kxyz(kx, ky, kz);
1377 
1378  D4 += (val * kxyz) % ((kxyz % (kxyz % kxyz)));
1379  }
1380  }
1381 
1382  return D4;
1383 }
1384 
1385 symtensor5<double> ewald::ewald_D5(double x, double y, double z)
1386 {
1387  static int printed = 0;
1388 
1389 #ifdef GRAVITY_TALLBOX
1390  Terminate("GRAVITY_TALLBOX is not implemented for MULTIPOLE_ORDER >= 4");
1391 #endif
1392 
1393  symtensor5<double> D5 = 0.0;
1394 
1395  double leff = pow((1.0 / LONG_X) * (1.0 / LONG_Y) * (1.0 / LONG_Z), 1.0 / 3);
1396  double alpha = 2.0 / leff;
1397  double alpha2 = alpha * alpha;
1398 
1399  int qxmax = (int)(8.0 * LONG_X / alpha + 0.5);
1400  int qymax = (int)(8.0 * LONG_Y / alpha + 0.5);
1401  int qzmax = (int)(8.0 * LONG_Z / alpha + 0.5);
1402 
1403  int nxmax = (int)(2.0 * alpha / LONG_X + 0.5);
1404  int nymax = (int)(2.0 * alpha / LONG_Y + 0.5);
1405  int nzmax = (int)(2.0 * alpha / LONG_Z + 0.5);
1406 
1407  if(printed == 0)
1408  {
1409  mpi_printf("EWALD: D5 table: qxmax=%d qymax=%d qzmax=%d nxmax=%d nymax=%d nzmax=%d\n", qxmax, qymax, qzmax, nxmax, nymax,
1410  nzmax);
1411  printed = 1;
1412  }
1413 
1414  for(int nx = -qxmax; nx <= qxmax; nx++)
1415  for(int ny = -qymax; ny <= qymax; ny++)
1416  for(int nz = -qzmax; nz <= qzmax; nz++)
1417  {
1418  double dx = x - nx * (1.0 / LONG_X);
1419  double dy = y - ny * (1.0 / LONG_Y);
1420  double dz = z - nz * (1.0 / LONG_Z);
1421 
1422  vector<double> dxyz(dx, dy, dz);
1423 
1424  double r2 = dx * dx + dy * dy + dz * dz;
1425  double r = sqrt(r2);
1426 
1427  double rinv = (r > 0) ? 1.0 / r : 0.0;
1428  double r2inv = rinv * rinv;
1429  double r3inv = r2inv * rinv;
1430  double r4inv = r2inv * r2inv;
1431  double r5inv = r2inv * r3inv;
1432  double r7inv = r3inv * r4inv;
1433  double r9inv = r4inv * r5inv;
1434  double r11inv = r4inv * r7inv;
1435 
1436  double g3, g4, g5;
1437 
1438  if(nx != 0 || ny != 0 || nz != 0)
1439  {
1440  g3 = (15.0 * erfc(alpha * r) +
1441  (30.0 * alpha * r + 20.0 * pow(alpha * r, 3) + 8.0 * pow(alpha * r, 5)) / sqrt(M_PI) * exp(-alpha2 * r2)) *
1442  r7inv;
1443 
1444  g4 = -(105.0 * erfc(alpha * r) +
1445  (210.0 * alpha * r + 140.0 * pow(alpha * r, 3) + 56.0 * pow(alpha * r, 5) + 16.0 * pow(alpha * r, 7)) /
1446  sqrt(M_PI) * exp(-alpha2 * r2)) *
1447  r9inv;
1448 
1449  g5 = (945.0 * erfc(alpha * r) + (1890.0 * alpha * r + 1260.0 * pow(alpha * r, 3) + 504.0 * pow(alpha * r, 5) +
1450  144.0 * pow(alpha * r, 7) + 32.0 * pow(alpha * r, 9)) /
1451  sqrt(M_PI) * exp(-alpha2 * r2)) *
1452  r11inv;
1453  }
1454  else
1455  {
1456  /* we add the 1/r term here to the (0|0|0) entry, followed by differentiation, and the limit r->0 to obtain accurate
1457  * results at the origin
1458  */
1459 
1460  /* Note, for small r:
1461  *
1462  * [1/- erfc(a r)]/r = 2 a/sqrt(pi) * [ 1 - (a r)^2/3 + (a r)^4 / 10 - (a r)^6 / 42 + (a r)^8 / 216 - ...]
1463  *
1464  * Hence for r = 0:
1465  *
1466  * g0 = 2 * alpha / sqrt(pi)
1467  * g1 = -4/3 * alpha^3 / sqrt(pi)
1468  * g2 = 8/5 * alpha^5 / sqrt(pi)
1469  * g3 = -16/7 * alpha^7 / sqrt(pi)
1470  * g4 = 32/9 * alpha^9 / sqrt(pi)
1471  * g5 = -64/11 * alpha^11/ sqrt(pi)
1472  */
1473 
1474  if((alpha * r) < 0.5)
1475  {
1476  g3 = 16.0 * pow(alpha, 7) / sqrt(M_PI) *
1477  (-1.0 / 7.0 + pow(alpha * r, 2) / 9.0 - pow(alpha * r, 4) / 22.0 + pow(alpha * r, 6) / 78.0 -
1478  pow(alpha * r, 8) / 360.0 + pow(alpha * r, 10) / 2040.0);
1479 
1480  g4 = 32.0 * pow(alpha, 9) / sqrt(M_PI) *
1481  (1.0 / 9.0 - pow(alpha * r, 2) / 11.0 + pow(alpha * r, 4) / 26.0 - pow(alpha * r, 6) / 90.0 +
1482  pow(alpha * r, 8) / 408.0 - pow(alpha * r, 10) / 2280.0);
1483 
1484  g5 = 64.0 * pow(alpha, 11) / sqrt(M_PI) *
1485  (-1.0 / 11.0 + pow(alpha * r, 2) / 13.0 - pow(alpha * r, 4) / 30.0 + pow(alpha * r, 6) / 102.0 -
1486  pow(alpha * r, 8) / 456.0 + pow(alpha * r, 10) / 2520.0);
1487  }
1488  else
1489  {
1490  g3 = (-15.0 * erf(alpha * r) +
1491  (30.0 * alpha * r + 20.0 * pow(alpha * r, 3) + 8.0 * pow(alpha * r, 5)) / sqrt(M_PI) * exp(-alpha2 * r2)) *
1492  r7inv;
1493 
1494  g4 = (105.0 * erf(alpha * r) -
1495  (210.0 * alpha * r + 140.0 * pow(alpha * r, 3) + 56.0 * pow(alpha * r, 5) + 16.0 * pow(alpha * r, 7)) /
1496  sqrt(M_PI) * exp(-alpha2 * r2)) *
1497  r9inv;
1498 
1499  g5 = (-945.0 * erf(alpha * r) + (1890.0 * alpha * r + 1260.0 * pow(alpha * r, 3) + 504.0 * pow(alpha * r, 5) +
1500  144.0 * pow(alpha * r, 7) + 32.0 * pow(alpha * r, 9)) /
1501  sqrt(M_PI) * exp(-alpha2 * r2)) *
1502  r11inv;
1503  }
1504  }
1505 
1506  symtensor3<double> aux3 = dxyz % (dxyz % dxyz);
1507  symtensor4<double> aux4 = dxyz % aux3;
1508  symtensor5<double> aux5;
1509 
1510  setup_D5(ADD, D5, dxyz, aux3, aux4, aux5, g3, g4, g5);
1511  }
1512 
1513  for(int nx = -nxmax; nx <= nxmax; nx++)
1514  for(int ny = -nymax; ny <= nymax; ny++)
1515  for(int nz = -nzmax; nz <= nzmax; nz++)
1516  {
1517  double kx = (2.0 * M_PI * LONG_X) * nx;
1518  double ky = (2.0 * M_PI * LONG_Y) * ny;
1519  double kz = (2.0 * M_PI * LONG_Z) * nz;
1520  double k2 = kx * kx + ky * ky + kz * kz;
1521 
1522  if(k2 > 0)
1523  {
1524  double kdotx = (x * kx + y * ky + z * kz);
1525  double val = 4.0 * M_PI * (LONG_X * LONG_Y * LONG_Z) / k2 * exp(-k2 / (4.0 * alpha2)) * sin(kdotx);
1526 
1527  vector<double> kxyz(kx, ky, kz);
1528 
1529  D5 += (val * kxyz) % (kxyz % ((kxyz % (kxyz % kxyz))));
1530  }
1531  }
1532 
1533  return D5;
1534 }
1535 
1536 symtensor6<double> ewald::ewald_D6(double x, double y, double z)
1537 {
1538  static int printed = 0;
1539 
1540 #ifdef GRAVITY_TALLBOX
1541  Terminate("GRAVITY_TALLBOX is not implemented for MULTIPOLE_ORDER >= 4");
1542 #endif
1543 
1544  symtensor6<double> D6 = 0.0;
1545 
1546  double leff = pow((1.0 / LONG_X) * (1.0 / LONG_Y) * (1.0 / LONG_Z), 1.0 / 3);
1547  double alpha = 2.0 / leff;
1548  double alpha2 = alpha * alpha;
1549 
1550  int qxmax = (int)(8.0 * LONG_X / alpha + 0.5);
1551  int qymax = (int)(8.0 * LONG_Y / alpha + 0.5);
1552  int qzmax = (int)(8.0 * LONG_Z / alpha + 0.5);
1553 
1554  int nxmax = (int)(2.0 * alpha / LONG_X + 0.5);
1555  int nymax = (int)(2.0 * alpha / LONG_Y + 0.5);
1556  int nzmax = (int)(2.0 * alpha / LONG_Z + 0.5);
1557 
1558  if(printed == 0)
1559  {
1560  mpi_printf("EWALD: D6 table: qxmax=%d qymax=%d qzmax=%d nxmax=%d nymax=%d nzmax=%d\n", qxmax, qymax, qzmax, nxmax, nymax,
1561  nzmax);
1562  printed = 1;
1563  }
1564 
1565  for(int nx = -qxmax; nx <= qxmax; nx++)
1566  for(int ny = -qymax; ny <= qymax; ny++)
1567  for(int nz = -qzmax; nz <= qzmax; nz++)
1568  {
1569  double dx = x - nx * (1.0 / LONG_X);
1570  double dy = y - ny * (1.0 / LONG_Y);
1571  double dz = z - nz * (1.0 / LONG_Z);
1572 
1573  vector<double> dxyz(dx, dy, dz);
1574 
1575  double r2 = dx * dx + dy * dy + dz * dz;
1576  double r = sqrt(r2);
1577 
1578  double rinv = (r > 0) ? 1.0 / r : 0.0;
1579  double r2inv = rinv * rinv;
1580  double r3inv = r2inv * rinv;
1581  double r4inv = r2inv * r2inv;
1582  double r5inv = r2inv * r3inv;
1583  double r7inv = r3inv * r4inv;
1584  double r9inv = r4inv * r5inv;
1585  double r11inv = r4inv * r7inv;
1586  double r13inv = r4inv * r9inv;
1587 
1588  double g3, g4, g5, g6;
1589 
1590  if(nx != 0 || ny != 0 || nz != 0)
1591  {
1592  g3 = (15.0 * erfc(alpha * r) +
1593  (30.0 * alpha * r + 20.0 * pow(alpha * r, 3) + 8.0 * pow(alpha * r, 5)) / sqrt(M_PI) * exp(-alpha2 * r2)) *
1594  r7inv;
1595 
1596  g4 = -(105.0 * erfc(alpha * r) +
1597  (210.0 * alpha * r + 140.0 * pow(alpha * r, 3) + 56.0 * pow(alpha * r, 5) + 16.0 * pow(alpha * r, 7)) /
1598  sqrt(M_PI) * exp(-alpha2 * r2)) *
1599  r9inv;
1600 
1601  g5 = (945.0 * erfc(alpha * r) + (1890.0 * alpha * r + 1260.0 * pow(alpha * r, 3) + 504.0 * pow(alpha * r, 5) +
1602  144.0 * pow(alpha * r, 7) + 32.0 * pow(alpha * r, 9)) /
1603  sqrt(M_PI) * exp(-alpha2 * r2)) *
1604  r11inv;
1605 
1606  g6 = (-10395.0 * erfc(alpha * r) -
1607  2.0 *
1608  (10395.0 * alpha * r + 6930.0 * pow(alpha * r, 3) + 2772.0 * pow(alpha * r, 5) + 792.0 * pow(alpha * r, 7) +
1609  176.0 * pow(alpha * r, 9) + 32.0 * pow(alpha * r, 11)) /
1610  sqrt(M_PI) * exp(-alpha2 * r2)) *
1611  r13inv;
1612  }
1613  else
1614  {
1615  /* we add the 1/r term here to the (0|0|0) entry, followed by differentiation, and the limit r->0 to obtain accurate
1616  * results at the origin
1617  */
1618 
1619  /* Note, for small r:
1620  *
1621  * [1/- erfc(a r)]/r = 2 a/sqrt(pi) * [ 1 - (a r)^2/3 + (a r)^4 / 10 - (a r)^6 / 42 + (a r)^8 / 216 - ...]
1622  *
1623  * Hence for r = 0:
1624  *
1625  * g0 = 2 * alpha / sqrt(pi)
1626  * g1 = -4/3 * alpha^3 / sqrt(pi)
1627  * g2 = 8/5 * alpha^5 / sqrt(pi)
1628  * g3 = -16/7 * alpha^7 / sqrt(pi)
1629  * g4 = 32/9 * alpha^9 / sqrt(pi)
1630  * g5 = -64/11 * alpha^11/ sqrt(pi)
1631  */
1632 
1633  if((alpha * r) < 0.5)
1634  {
1635  g3 = 16.0 * pow(alpha, 7) / sqrt(M_PI) *
1636  (-1.0 / 7.0 + pow(alpha * r, 2) / 9.0 - pow(alpha * r, 4) / 22.0 + pow(alpha * r, 6) / 78.0 -
1637  pow(alpha * r, 8) / 360.0 + pow(alpha * r, 10) / 2040.0);
1638 
1639  g4 = 32.0 * pow(alpha, 9) / sqrt(M_PI) *
1640  (1.0 / 9.0 - pow(alpha * r, 2) / 11.0 + pow(alpha * r, 4) / 26.0 - pow(alpha * r, 6) / 90.0 +
1641  pow(alpha * r, 8) / 408.0 - pow(alpha * r, 10) / 2280.0);
1642 
1643  g5 = 64.0 * pow(alpha, 11) / sqrt(M_PI) *
1644  (-1.0 / 11.0 + pow(alpha * r, 2) / 13.0 - pow(alpha * r, 4) / 30.0 + pow(alpha * r, 6) / 102.0 -
1645  pow(alpha * r, 8) / 456.0 + pow(alpha * r, 10) / 2520.0);
1646 
1647  g6 = 128.0 * pow(alpha, 13) / sqrt(M_PI) *
1648  (1.0 / 13.0 - pow(alpha * r, 2) / 15.0 + pow(alpha * r, 4) / 34.0 - pow(alpha * r, 6) / 114.0 +
1649  pow(alpha * r, 8) / 504.0 - pow(alpha * r, 10) / 2760.0);
1650  }
1651  else
1652  {
1653  g3 = (-15.0 * erf(alpha * r) +
1654  (30.0 * alpha * r + 20.0 * pow(alpha * r, 3) + 8.0 * pow(alpha * r, 5)) / sqrt(M_PI) * exp(-alpha2 * r2)) *
1655  r7inv;
1656 
1657  g4 = (105.0 * erf(alpha * r) -
1658  (210.0 * alpha * r + 140.0 * pow(alpha * r, 3) + 56.0 * pow(alpha * r, 5) + 16.0 * pow(alpha * r, 7)) /
1659  sqrt(M_PI) * exp(-alpha2 * r2)) *
1660  r9inv;
1661 
1662  g5 = (-945.0 * erf(alpha * r) + (1890.0 * alpha * r + 1260.0 * pow(alpha * r, 3) + 504.0 * pow(alpha * r, 5) +
1663  144.0 * pow(alpha * r, 7) + 32.0 * pow(alpha * r, 9)) /
1664  sqrt(M_PI) * exp(-alpha2 * r2)) *
1665  r11inv;
1666 
1667  g6 = (10395.0 * erf(alpha * r) -
1668  2.0 *
1669  (10395.0 * alpha * r + 6930.0 * pow(alpha * r, 3) + 2772.0 * pow(alpha * r, 5) +
1670  792.0 * pow(alpha * r, 7) + 176.0 * pow(alpha * r, 9) + 32.0 * pow(alpha * r, 11)) /
1671  sqrt(M_PI) * exp(-alpha2 * r2)) *
1672  r13inv;
1673  }
1674  }
1675 
1676  setup_D6(ADD, D6, dxyz, g3, g4, g5, g6);
1677  }
1678 
1679  for(int nx = -nxmax; nx <= nxmax; nx++)
1680  for(int ny = -nymax; ny <= nymax; ny++)
1681  for(int nz = -nzmax; nz <= nzmax; nz++)
1682  {
1683  double kx = (2.0 * M_PI * LONG_X) * nx;
1684  double ky = (2.0 * M_PI * LONG_Y) * ny;
1685  double kz = (2.0 * M_PI * LONG_Z) * nz;
1686  double k2 = kx * kx + ky * ky + kz * kz;
1687 
1688  if(k2 > 0)
1689  {
1690  double kdotx = (x * kx + y * ky + z * kz);
1691  double val = 4.0 * M_PI * (LONG_X * LONG_Y * LONG_Z) / k2 * exp(-k2 / (4.0 * alpha2)) * cos(kdotx);
1692 
1693  vector<double> kxyz(kx, ky, kz);
1694 
1695  D6 += (val * kxyz) % (kxyz % (kxyz % ((kxyz % (kxyz % kxyz)))));
1696  }
1697  }
1698 
1699  return D6;
1700 }
1701 
1702 symtensor7<double> ewald::ewald_D7(double x, double y, double z)
1703 {
1704  static int printed = 0;
1705 
1706 #ifdef GRAVITY_TALLBOX
1707  Terminate("GRAVITY_TALLBOX is not implemented for MULTIPOLE_ORDER >= 4");
1708 #endif
1709 
1710  symtensor7<double> D7 = 0.0;
1711 
1712  double leff = pow((1.0 / LONG_X) * (1.0 / LONG_Y) * (1.0 / LONG_Z), 1.0 / 3);
1713  double alpha = 2.0 / leff;
1714  double alpha2 = alpha * alpha;
1715 
1716  int qxmax = (int)(8.0 * LONG_X / alpha + 0.5);
1717  int qymax = (int)(8.0 * LONG_Y / alpha + 0.5);
1718  int qzmax = (int)(8.0 * LONG_Z / alpha + 0.5);
1719 
1720  int nxmax = (int)(2.0 * alpha / LONG_X + 0.5);
1721  int nymax = (int)(2.0 * alpha / LONG_Y + 0.5);
1722  int nzmax = (int)(2.0 * alpha / LONG_Z + 0.5);
1723 
1724  if(printed == 0)
1725  {
1726  mpi_printf("EWALD: D7 table: qxmax=%d qymax=%d qzmax=%d nxmax=%d nymax=%d nzmax=%d\n", qxmax, qymax, qzmax, nxmax, nymax,
1727  nzmax);
1728  printed = 1;
1729  }
1730 
1731  for(int nx = -qxmax; nx <= qxmax; nx++)
1732  for(int ny = -qymax; ny <= qymax; ny++)
1733  for(int nz = -qzmax; nz <= qzmax; nz++)
1734  {
1735  double dx = x - nx * (1.0 / LONG_X);
1736  double dy = y - ny * (1.0 / LONG_Y);
1737  double dz = z - nz * (1.0 / LONG_Z);
1738 
1739  vector<double> dxyz(dx, dy, dz);
1740 
1741  double r2 = dx * dx + dy * dy + dz * dz;
1742  double r = sqrt(r2);
1743 
1744  double rinv = (r > 0) ? 1.0 / r : 0.0;
1745  double r2inv = rinv * rinv;
1746  double r3inv = r2inv * rinv;
1747  double r4inv = r2inv * r2inv;
1748  double r5inv = r2inv * r3inv;
1749  double r7inv = r3inv * r4inv;
1750  double r9inv = r4inv * r5inv;
1751  double r11inv = r4inv * r7inv;
1752  double r13inv = r4inv * r9inv;
1753  double r15inv = r4inv * r11inv;
1754 
1755  double g4, g5, g6, g7;
1756 
1757  if(nx != 0 || ny != 0 || nz != 0)
1758  {
1759  g4 = -(105.0 * erfc(alpha * r) +
1760  (210.0 * alpha * r + 140.0 * pow(alpha * r, 3) + 56.0 * pow(alpha * r, 5) + 16.0 * pow(alpha * r, 7)) /
1761  sqrt(M_PI) * exp(-alpha2 * r2)) *
1762  r9inv;
1763 
1764  g5 = (945.0 * erfc(alpha * r) + (1890.0 * alpha * r + 1260.0 * pow(alpha * r, 3) + 504.0 * pow(alpha * r, 5) +
1765  144.0 * pow(alpha * r, 7) + 32.0 * pow(alpha * r, 9)) /
1766  sqrt(M_PI) * exp(-alpha2 * r2)) *
1767  r11inv;
1768 
1769  g6 = (-10395.0 * erfc(alpha * r) -
1770  2.0 *
1771  (10395.0 * alpha * r + 6930.0 * pow(alpha * r, 3) + 2772.0 * pow(alpha * r, 5) + 792.0 * pow(alpha * r, 7) +
1772  176.0 * pow(alpha * r, 9) + 32.0 * pow(alpha * r, 11)) /
1773  sqrt(M_PI) * exp(-alpha2 * r2)) *
1774  r13inv;
1775 
1776  g7 =
1777  (135135.0 * erfc(alpha * r) + 2.0 *
1778  (135135.0 * alpha * r + 90090.0 * pow(alpha * r, 3) + 36036.0 * pow(alpha * r, 5) +
1779  10296.0 * pow(alpha * r, 7) + 2288.0 * pow(alpha * r, 9) +
1780  416.0 * pow(alpha * r, 11) + 64.0 * pow(alpha * r, 13)) /
1781  sqrt(M_PI) * exp(-alpha2 * r2)) *
1782  r15inv;
1783  }
1784  else
1785  {
1786  /* we add the 1/r term here to the (0|0|0) entry, followed by differentiation, and the limit r->0 to obtain accurate
1787  * results at the origin
1788  */
1789 
1790  /* Note, for small r:
1791  *
1792  * [1/- erfc(a r)]/r = 2 a/sqrt(pi) * [ 1 - (a r)^2/3 + (a r)^4 / 10 - (a r)^6 / 42 + (a r)^8 / 216 - ...]
1793  *
1794  * Hence for r = 0:
1795  *
1796  * g0 = 2 * alpha / sqrt(pi)
1797  * g1 = -4/3 * alpha^3 / sqrt(pi)
1798  * g2 = 8/5 * alpha^5 / sqrt(pi)
1799  * g3 = -16/7 * alpha^7 / sqrt(pi)
1800  * g4 = 32/9 * alpha^9 / sqrt(pi)
1801  * g5 = -64/11 * alpha^11/ sqrt(pi)
1802  */
1803 
1804  if((alpha * r) < 0.5)
1805  {
1806  g4 = 32.0 * pow(alpha, 9) / sqrt(M_PI) *
1807  (1.0 / 9.0 - pow(alpha * r, 2) / 11.0 + pow(alpha * r, 4) / 26.0 - pow(alpha * r, 6) / 90.0 +
1808  pow(alpha * r, 8) / 408.0 - pow(alpha * r, 10) / 2280.0);
1809 
1810  g5 = 64.0 * pow(alpha, 11) / sqrt(M_PI) *
1811  (-1.0 / 11.0 + pow(alpha * r, 2) / 13.0 - pow(alpha * r, 4) / 30.0 + pow(alpha * r, 6) / 102.0 -
1812  pow(alpha * r, 8) / 456.0 + pow(alpha * r, 10) / 2520.0);
1813 
1814  g6 = 128.0 * pow(alpha, 13) / sqrt(M_PI) *
1815  (1.0 / 13.0 - pow(alpha * r, 2) / 15.0 + pow(alpha * r, 4) / 34.0 - pow(alpha * r, 6) / 114.0 +
1816  pow(alpha * r, 8) / 504.0 - pow(alpha * r, 10) / 2760.0);
1817 
1818  g7 = 256.0 * pow(alpha, 15) / sqrt(M_PI) *
1819  (-1.0 / 15.0 + pow(alpha * r, 2) / 17.0 - pow(alpha * r, 4) / 38.0 + pow(alpha * r, 6) / 126.0 -
1820  pow(alpha * r, 8) / 552.0 + pow(alpha * r, 10) / 3000.0);
1821  }
1822  else
1823  {
1824  g4 = (105.0 * erf(alpha * r) -
1825  (210.0 * alpha * r + 140.0 * pow(alpha * r, 3) + 56.0 * pow(alpha * r, 5) + 16.0 * pow(alpha * r, 7)) /
1826  sqrt(M_PI) * exp(-alpha2 * r2)) *
1827  r9inv;
1828 
1829  g5 = (-945.0 * erf(alpha * r) + (1890.0 * alpha * r + 1260.0 * pow(alpha * r, 3) + 504.0 * pow(alpha * r, 5) +
1830  144.0 * pow(alpha * r, 7) + 32.0 * pow(alpha * r, 9)) /
1831  sqrt(M_PI) * exp(-alpha2 * r2)) *
1832  r11inv;
1833 
1834  g6 = (10395.0 * erf(alpha * r) -
1835  2.0 *
1836  (10395.0 * alpha * r + 6930.0 * pow(alpha * r, 3) + 2772.0 * pow(alpha * r, 5) +
1837  792.0 * pow(alpha * r, 7) + 176.0 * pow(alpha * r, 9) + 32.0 * pow(alpha * r, 11)) /
1838  sqrt(M_PI) * exp(-alpha2 * r2)) *
1839  r13inv;
1840 
1841  g7 = (-135135.0 * erf(alpha * r) +
1842  2.0 *
1843  (135135.0 * alpha * r + 90090.0 * pow(alpha * r, 3) + 36036.0 * pow(alpha * r, 5) +
1844  10296.0 * pow(alpha * r, 7) + 2288.0 * pow(alpha * r, 9) + 416.0 * pow(alpha * r, 11) +
1845  64.0 * pow(alpha * r, 13)) /
1846  sqrt(M_PI) * exp(-alpha2 * r2)) *
1847  r15inv;
1848  }
1849  }
1850 
1851  setup_D7(ADD, D7, dxyz, g4, g5, g6, g7);
1852  }
1853 
1854  for(int nx = -nxmax; nx <= nxmax; nx++)
1855  for(int ny = -nymax; ny <= nymax; ny++)
1856  for(int nz = -nzmax; nz <= nzmax; nz++)
1857  {
1858  double kx = (2.0 * M_PI * LONG_X) * nx;
1859  double ky = (2.0 * M_PI * LONG_Y) * ny;
1860  double kz = (2.0 * M_PI * LONG_Z) * nz;
1861  double k2 = kx * kx + ky * ky + kz * kz;
1862 
1863  if(k2 > 0)
1864  {
1865  double kdotx = (x * kx + y * ky + z * kz);
1866  double val = -4.0 * M_PI * (LONG_X * LONG_Y * LONG_Z) / k2 * exp(-k2 / (4.0 * alpha2)) * sin(kdotx);
1867 
1868  vector<double> kxyz(kx, ky, kz);
1869 
1870  D7 += (val * kxyz) % (kxyz % (kxyz % (kxyz % (kxyz % (kxyz % kxyz)))));
1871  }
1872  }
1873 
1874  return D7;
1875 }
global_data_all_processes All
Definition: main.cc:40
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_init(void)
This function initializes tables with the correction force and the correction potential due to the pe...
Definition: ewald.cc:69
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
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
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
MyReal FacCoordToInt
Definition: intposconvert.h:38
void nearest_image_intpos_to_pos(const MyIntPosType *const a, const MyIntPosType *const b, T *posdiff)
MyReal RegionLen
Definition: intposconvert.h:39
MyReal FacIntToCoord
Definition: intposconvert.h:37
size_t my_fread(void *ptr, size_t size, size_t nmemb, FILE *stream)
A wrapper for the fread() function.
size_t my_fwrite(const void *ptr, size_t size, size_t nmemb, FILE *stream)
A wrapper for the fwrite() function.
void mpi_printf(const char *fmt,...)
Definition: setcomm.h:55
int ThisTask
Definition: setcomm.h:33
int NTask
Definition: setcomm.h:32
MPI_Comm Communicator
Definition: setcomm.h:31
MPI_Comm SharedMemComm
int Island_ThisTask
int Island_NTask
void ** SharedMemBaseAddr
int MyShmRankInGlobal
T da[3]
Definition: symtensors.h:106
#define M_PI
Definition: constants.h:56
float MyFloat
Definition: dtypes.h:86
#define LONG_Y
Definition: dtypes.h:369
#define DBY
Definition: dtypes.h:420
uint32_t MyIntPosType
Definition: dtypes.h:35
double MyReal
Definition: dtypes.h:82
#define LONG_X
Definition: dtypes.h:361
#define BITS_FOR_POSITIONS
Definition: dtypes.h:37
#define DBX
Definition: dtypes.h:419
#define LONG_Z
Definition: dtypes.h:377
#define DBZ
Definition: dtypes.h:421
#define EN
Base dimension of cubical Ewald lookup table, for one octant.
Definition: ewald.h:43
#define ENX
Definition: ewald.h:45
#define ENZ
Definition: ewald.h:47
#define ENY
Definition: ewald.h:46
#define EWLEVEL
Definition: ewald.h:41
#define HIGHEST_NEEDEDORDER_EWALD_DPHI
Definition: gravtree.h:153
#define EWALD_TAYLOR_ORDER
Definition: gravtree.h:179
#define Terminate(...)
Definition: macros.h:19
shmem Shmem
Definition: main.cc:45
#define TAG_DMOM
Definition: mpi_utils.h:30
#define TAG_EWALD_ALLOC
Definition: mpi_utils.h:22
memory Mem
Definition: main.cc:44
expr exp(half arg)
Definition: half.hpp:2724
expr cos(half arg)
Definition: half.hpp:2823
expr sin(half arg)
Definition: half.hpp:2816
expr erf(half arg)
Definition: half.hpp:2918
expr pow(half base, half exp)
Definition: half.hpp:2803
expr erfc(half arg)
Definition: half.hpp:2925
expr sqrt(half arg)
Definition: half.hpp:2777
vector< MyReal > D1phi
Definition: gravtree.h:186
symtensor3< MyReal > D3phi
Definition: gravtree.h:188
MyReal D0phi
Definition: gravtree.h:185
symtensor2< MyReal > D2phi
Definition: gravtree.h:187
#define tXXXZZZZ
#define tZZZZZZZ
#define tXYYZZZZ
#define pXYYYZZ
#define rYZZZZ
#define tYYZZZZZ
#define rXXYZZ
#define tXYYYYYZ
#define dYZZ
#define tXXYYYYZ
#define rYYZZZ
#define sXXXY
#define pXXXXYZ
#define rXXZZZ
#define dXYY
#define rXXXYY
#define tXXXXXXX
#define sXXYZ
#define pXYYYYY
#define tXXZZZZZ
#define tXYYYZZZ
#define pXYYZZZ
#define tYYYZZZZ
#define qXX
#define sXYYZ
#define pXXXYYY
#define dZZZ
#define rXXXYZ
#define tYYYYYZZ
#define dYYY
#define dXYZ
#define pXXYYYZ
#define sYYYZ
#define rXXYYY
#define tXXXYYYZ
#define tXYZZZZZ
#define rXXXXZ
#define sXZZZ
#define qYY
#define rXXXXY
#define tXXYYYYY
#define rYYYYZ
#define pXXXXXZ
#define rXYYYZ
#define tXXXXYZZ
#define tXXXYZZZ
#define pXYZZZZ
#define qYZ
#define tXXXYYZZ
#define tYYYYYYZ
#define rXYYYY
#define sXYYY
#define dXXY
#define pXXXYYZ
#define tXXXXXYZ
#define tXXXXZZZ
#define pYZZZZZ
#define dXXX
#define dXXZ
#define rYYYZZ
#define pXYYYYZ
#define dYYZ
#define pYYYZZZ
#define tXXYZZZZ
#define tXXXXXXY
#define rXXXXX
#define tXXXXXZZ
#define pXXXYZZ
#define rXXYYZ
#define tYZZZZZZ
#define tXXXXYYY
#define sXYZZ
#define qZZ
#define pXXXXXY
#define tXXXXXXZ
#define rZZZZZ
#define rXXXZZ
#define rXYYZZ
#define sXXXZ
#define tXYYYYZZ
#define sYZZZ
#define pXXYZZZ
#define tXYYYYYY
#define tXXXYYYY
#define tYYYYYYY
#define pYYYYYZ
#define qXY
#define rXYZZZ
#define tYYYYZZZ
#define tXXXXXYY
#define rYYYYY
#define tXZZZZZZ
#define pXXXZZZ
#define tXXYYZZZ
#define dXZZ
#define tXXXXYYZ
#define rXZZZZ
#define qXZ
#define tXXYYYZZ
#define pXZZZZZ
void setup_D6(enum setup_options opt, symtensor6< T > &D6, vector< T > &dxyz, TypeGfac g3, TypeGfac g4, TypeGfac g5, TypeGfac g6)
Definition: symtensors.h:1912
void setup_D5(enum setup_options opt, symtensor5< T > &D5, vector< T > &dxyz, symtensor3< T > &aux3, symtensor4< T > &aux4, symtensor5< T > &aux5, TypeGfac g3, TypeGfac g4, TypeGfac g5)
Definition: symtensors.h:1842
void setup_D3(enum setup_options opt, symtensor3< T > &D3, vector< T > &dxyz, symtensor2< T > &aux2, symtensor3< T > &aux3, TypeGfac g2, TypeGfac g3)
Definition: symtensors.h:1775
void setup_D7(enum setup_options opt, symtensor7< T > &D7, vector< T > &dxyz, TypeGfac g4, TypeGfac g5, TypeGfac g6, TypeGfac g7)
Definition: symtensors.h:2028
@ ADD
Definition: symtensors.h:1771
void setup_D4(enum setup_options opt, symtensor4< T > &D4, vector< T > &dxyz, symtensor2< T > &aux2, symtensor3< T > &aux3, symtensor4< T > &aux4, TypeGfac g2, TypeGfac g3, TypeGfac g4)
Definition: symtensors.h:1801
void subdivide_evenly(long long N, int pieces, int index_bin, long long *first, int *count)
Definition: system.cc:51
void myflush(FILE *fstream)
Definition: system.cc:125