12 #include "gadgetconfig.h"
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"
27 #include "../main/simulation.h"
28 #include "../mpi_utils/shared_mem_handler.h"
29 #include "../sort/cxxsort.h"
30 #include "../system/system.h"
71 mpi_printf(
"EWALD: initialize Ewald correction...\n");
83 int recomputeflag = 0;
88 if((fd = fopen(buf,
"r")))
90 mpi_printf(
"\nEWALD: reading Ewald tables from file `%s'\n", buf);
95 #ifndef GRAVITY_TALLBOX
98 int ewaldtype = GRAVITY_TALLBOX + 1;
103 mpi_printf(
"\nEWALD: something's wrong with this table file. Discarding it.\n");
122 mpi_printf(
"\nEWALD: No usable Ewald tables in file `%s' found. Recomputing them...\n", buf);
126 int size = (
ENX + 1) * (
ENY + 1) * (
ENZ + 1);
131 for(
int n = first; n < first + count; n++)
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));
139 if(((n - first) % (count / 20)) == 0)
141 printf(
"%4.1f percent done\n", (n - first) / (count / 100.0));
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;
152 #ifndef GRAVITY_TALLBOX
157 #if(HIGHEST_NEEDEDORDER_EWALD_DPHI + EWALD_TAYLOR_ORDER) >= 4
160 #if(HIGHEST_NEEDEDORDER_EWALD_DPHI + EWALD_TAYLOR_ORDER) >= 5
163 #if(HIGHEST_NEEDEDORDER_EWALD_DPHI + EWALD_TAYLOR_ORDER) >= 6
166 #if(HIGHEST_NEEDEDORDER_EWALD_DPHI + EWALD_TAYLOR_ORDER) >= 7
174 switch(GRAVITY_TALLBOX)
177 ewdp->D1phi[0] = force[2];
178 ewdp->D1phi[1] = force[0];
179 ewdp->D1phi[2] = force[1];
183 ewdp->D1phi[0] = force[0];
184 ewdp->D1phi[1] = force[2];
185 ewdp->D1phi[2] = force[1];
189 ewdp->D1phi[0] = force[0];
190 ewdp->D1phi[1] = force[1];
191 ewdp->D1phi[2] = force[2];
197 int *recvcnts = (
int *)
Mem.mymalloc(
"recvcnts",
NTask *
sizeof(
int));
198 int *recvoffs = (
int *)
Mem.mymalloc(
"recvoffs",
NTask *
sizeof(
int));
200 for(
int i = 0; i <
NTask; i++)
208 MPI_Allgatherv(MPI_IN_PLACE, size *
sizeof(
ewald_data), MPI_BYTE, Ewd, recvcnts, recvoffs, MPI_BYTE,
Communicator);
210 Mem.myfree(recvoffs);
211 Mem.myfree(recvcnts);
213 mpi_printf(
"\nEWALD: writing Ewald tables to file `%s'\n", buf);
217 if((fd = fopen(buf,
"w")))
224 #ifndef GRAVITY_TALLBOX
235 Terminate(
"can't write to file '%s'\n", buf);
250 for(
int i = 0; i <=
ENX; i++)
251 for(
int j = 0; j <=
ENY; j++)
252 for(
int k = 0; k <=
ENZ; k++)
260 #if(HIGHEST_NEEDEDORDER_EWALD_DPHI + EWALD_TAYLOR_ORDER) >= 4
263 #if(HIGHEST_NEEDEDORDER_EWALD_DPHI + EWALD_TAYLOR_ORDER) >= 5
266 #if(HIGHEST_NEEDEDORDER_EWALD_DPHI + EWALD_TAYLOR_ORDER) >= 6
269 #if(HIGHEST_NEEDEDORDER_EWALD_DPHI + EWALD_TAYLOR_ORDER) >= 7
274 mpi_printf(
"EWALD: Initialization of periodic boundaries finished.\n");
276 ewald_is_initialized = 1;
300 test_interpolation_accuracy();
313 MyIntPosType temppos[3] = {p_intpos[0] - target_intpos[0], p_intpos[1] - target_intpos[1], p_intpos[2] - target_intpos[2]};
316 gridpos[0] = (temppos[0] + halflen) & ewaldmask;
317 gridpos[1] = (temppos[1] + halflen) & ewaldmask;
318 gridpos[2] = (temppos[2] + halflen) & ewaldmask;
327 int signx = 1, signy = 1, signz = 1;
334 else if(i ==
EN && gridpos[0] < temppos[0])
342 else if(j ==
EN && gridpos[1] < temppos[1])
350 else if(k ==
EN && gridpos[2] < temppos[2])
353 fper = Ewd[ewd_offset(i, j, k)];
357 fper.
D1phi[0] *= signx;
358 fper.
D1phi[1] *= signy;
359 fper.
D1phi[2] *= signz;
369 fper.
D3phi[
dXYZ] *= signx * signy * signz;
376 #if EWALD_TAYLOR_ORDER == 3
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;
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);
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;
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;
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;
415 fper.D5phi[
rXXYZZ] *= signy;
416 fper.D5phi[
rXXYYZ] *= signz;
417 fper.D5phi[
rXYYZZ] *= signx;
419 fper.D5phi[
rXXXYZ] *= signx * signy * signz;
420 fper.D5phi[
rXYYYZ] *= signx * signy * signz;
421 fper.D5phi[
rXYZZZ] *= signx * signy * signz;
423 fper.
D2phi += fper.
D3phi * off + 0.5 * ((fper.D4phi * off) * off) + (1.0 / 6) * (((fper.D5phi * off) * off) * off);
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;
446 fper.
D3phi += fper.D4phi * off + 0.5 * ((fper.D5phi * off) * off) + (1.0 / 6) * (((fper.D6phi * off) * off) * off);
449 #if(HIGHEST_NEEDEDORDER_EWALD_DPHI + EWALD_TAYLOR_ORDER) >= 7
454 fper.D7phi[
tXXXXXYZ] *= signx * signy * signz;
461 fper.D7phi[
tXXXYYYZ] *= signx * signy * signz;
463 fper.D7phi[
tXXXYZZZ] *= signx * signy * signz;
472 fper.D7phi[
tXYYYYYZ] *= signx * signy * signz;
474 fper.D7phi[
tXYYYZZZ] *= signx * signy * signz;
476 fper.D7phi[
tXYZZZZZ] *= signx * signy * signz;
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);
497 #ifndef GRAVITY_TALLBOX
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;
516 fper.
D2phi += fper.
D3phi * off + 0.5 * ((fper.D4phi * off) * off);
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;
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;
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;
538 fper.D5phi[
rXXYZZ] *= signy;
539 fper.D5phi[
rXXYYZ] *= signz;
540 fper.D5phi[
rXYYZZ] *= signx;
542 fper.D5phi[
rXXXYZ] *= signx * signy * signz;
543 fper.D5phi[
rXYYYZ] *= signx * signy * signz;
544 fper.D5phi[
rXYZZZ] *= signx * signy * signz;
546 fper.
D3phi += fper.D4phi * off + 0.5 * ((fper.D5phi * off) * off);
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;
569 fper.D4phi += fper.D5phi * off + 0.5 * ((fper.D6phi * off) * off);
572 #if(HIGHEST_NEEDEDORDER_EWALD_DPHI + EWALD_TAYLOR_ORDER) >= 7
577 fper.D7phi[
tXXXXXYZ] *= signx * signy * signz;
584 fper.D7phi[
tXXXYYYZ] *= signx * signy * signz;
586 fper.D7phi[
tXXXYZZZ] *= signx * signy * signz;
595 fper.D7phi[
tXYYYYYZ] *= signx * signy * signz;
597 fper.D7phi[
tXYYYZZZ] *= signx * signy * signz;
599 fper.D7phi[
tXYZZZZZ] *= signx * signy * signz;
610 fper.D5phi += fper.D6phi * off + 0.5 * ((fper.D7phi * off) * off);
625 static int printed = 0;
629 #ifndef GRAVITY_TALLBOX
632 double alpha = 2.0 / leff;
633 double alpha2 = alpha * alpha;
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);
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);
645 mpi_printf(
"EWALD: D0 table: qxmax=%d qymax=%d qzmax=%d nxmax=%d nymax=%d nzmax=%d\n", qxmax, qymax, qzmax, nxmax, nymax,
650 for(
int nx = -qxmax; nx <= qxmax; nx++)
651 for(
int ny = -qymax; ny <= qymax; ny++)
652 for(
int nz = -qzmax; nz <= qzmax; nz++)
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);
658 double r2 = dx * dx + dy * dy + dz * dz;
661 double rinv = (r > 0) ? 1.0 / r : 0.0;
665 if(nx != 0 || ny != 0 || nz != 0)
667 g0 = -
erfc(alpha * r) * rinv;
680 if((alpha * r) < 0.5)
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);
688 g0 =
erf(alpha * r) * rinv;
695 for(
int nx = -nxmax; nx <= nxmax; nx++)
696 for(
int ny = -nymax; ny <= nymax; ny++)
697 for(
int nz = -nzmax; nz <= nzmax; nz++)
699 if(nx != 0 || ny != 0 || nz != 0)
704 double k2 = kx * kx + ky * ky + kz * kz;
705 double kdotx = (x * kx + y * ky + z * kz);
716 double leff =
sqrt(BOXX * BOXY);
717 double alpha = 2.0 / leff;
719 int qxmax = (int)(8.0 / (BOXX * alpha) + 0.5);
720 int qymax = (int)(8.0 / (BOXY * alpha) + 0.5);
722 int nxmax = (int)(2.0 * alpha * BOXX + 0.5);
723 int nymax = (int)(2.0 * alpha * BOXY + 0.5);
727 mpi_printf(
"EWALD: D0 table: qxmax=%d qymax=%d nxmax=%d nymax=%d\n", qxmax, qymax, nxmax, nymax);
731 for(
int nx = -qxmax; nx <= qxmax; nx++)
732 for(
int ny = -qymax; ny <= qymax; ny++)
734 double dx = x - nx * BOXX;
735 double dy = y - ny * BOXY;
736 double r =
sqrt(dx * dx + dy * dy + z * z);
738 double rinv = (r > 0) ? 1.0 / r : 0.0;
742 if(nx != 0 || ny != 0)
744 g0 = -
erfc(alpha * r) * rinv;
750 if((alpha * r) < 0.5)
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);
758 g0 =
erf(alpha * r) * rinv;
765 double alpha2 = alpha * alpha;
767 for(
int nx = -nxmax; nx <= nxmax; nx++)
768 for(
int ny = -nymax; ny <= nymax; ny++)
770 if(nx != 0 || ny != 0)
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;
779 double ex =
exp(-k * z);
781 D0 += -
M_PI / (BOXX * BOXY) * (
erfc(k / (2 * alpha) + alpha * z) / ex + ex *
erfc(k / (2 * alpha) - alpha * z)) / k;
785 double ex =
exp(k * z);
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;
809 static int printed = 0;
813 #ifndef GRAVITY_TALLBOX
816 double alpha = 2.0 / leff;
817 double alpha2 = alpha * alpha;
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);
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);
829 mpi_printf(
"EWALD: D1 table: qxmax=%d qymax=%d qzmax=%d nxmax=%d nymax=%d nzmax=%d\n", qxmax, qymax, qzmax, nxmax, nymax,
834 for(
int nx = -qxmax; nx <= qxmax; nx++)
835 for(
int ny = -qymax; ny <= qymax; ny++)
836 for(
int nz = -qzmax; nz <= qzmax; nz++)
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);
844 double r2 = dx * dx + dy * dy + dz * dz;
847 double rinv = (r > 0) ? 1.0 / r : 0.0;
848 double r2inv = rinv * rinv;
849 double r3inv = r2inv * rinv;
853 if(nx != 0 || ny != 0 || nz != 0)
855 g1 = (
erfc(alpha * r) + 2.0 * alpha * r /
sqrt(
M_PI) *
exp(-alpha2 * r2)) * r3inv;
877 if((alpha * r) < 0.5)
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);
885 g1 = (-
erf(alpha * r) + 2.0 * alpha * r /
sqrt(
M_PI) *
exp(-alpha2 * r2)) * r3inv;
892 for(
int nx = -nxmax; nx <= nxmax; nx++)
893 for(
int ny = -nymax; ny <= nymax; ny++)
894 for(
int nz = -nzmax; nz <= nzmax; nz++)
899 double k2 = kx * kx + ky * ky + kz * kz;
903 double kdotx = (x * kx + y * ky + z * kz);
913 double leff =
sqrt(BOXX * BOXY);
914 double alpha = 2.0 / leff;
915 double alpha2 = alpha * alpha;
917 int qxmax = (int)(8.0 / (BOXX * alpha) + 0.5);
918 int qymax = (int)(8.0 / (BOXY * alpha) + 0.5);
920 int nxmax = (int)(2.0 * alpha * BOXX + 0.5);
921 int nymax = (int)(2.0 * alpha * BOXY + 0.5);
925 mpi_printf(
"EWALD: D1 table: qxmax=%d qymax=%d nxmax=%d nymax=%d\n", qxmax, qymax, nxmax, nymax);
929 for(
int nx = -qxmax; nx <= qxmax; nx++)
930 for(
int ny = -qymax; ny <= qymax; ny++)
932 double dx = x - nx * BOXX;
933 double dy = y - ny * BOXY;
938 double r2 = dx * dx + dy * dy + dz * dz;
941 double rinv = (r > 0) ? 1.0 / r : 0.0;
942 double r2inv = rinv * rinv;
943 double r3inv = r2inv * rinv;
947 if(nx != 0 || ny != 0)
949 g1 = (
erfc(alpha * r) + 2.0 * alpha * r /
sqrt(
M_PI) *
exp(-alpha2 * r2)) * r3inv;
955 if((alpha * r) < 0.5)
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);
963 g1 = (-
erf(alpha * r) + 2.0 * alpha * r /
sqrt(
M_PI) *
exp(-alpha2 * r2)) * r3inv;
970 for(
int nx = -nxmax; nx <= nxmax; nx++)
971 for(
int ny = -nymax; ny <= nymax; ny++)
973 if(nx != 0 || ny != 0)
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;
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;
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));
990 D1[2] += 2.0 *
M_PI / (BOXX * BOXY) *
erf(alpha * z);
998 static int printed = 0;
1003 double alpha = 2.0 / leff;
1004 double alpha2 = alpha * alpha;
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);
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);
1016 mpi_printf(
"EWALD: D2 table: qxmax=%d qymax=%d qzmax=%d nxmax=%d nymax=%d nzmax=%d\n", qxmax, qymax, qzmax, nxmax, nymax,
1021 for(
int nx = -qxmax; nx <= qxmax; nx++)
1022 for(
int ny = -qymax; ny <= qymax; ny++)
1023 for(
int nz = -qzmax; nz <= qzmax; nz++)
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);
1031 double r2 = dx * dx + dy * dy + dz * dz;
1032 double r =
sqrt(r2);
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;
1041 if(nx != 0 || ny != 0 || nz != 0)
1043 g1 = (
erfc(alpha * r) + 2.0 * alpha * r /
sqrt(
M_PI) *
exp(-alpha2 * r2)) * r3inv;
1045 g2 = -(3.0 *
erfc(alpha * r) + (6.0 * alpha * r + 4.0 *
pow(alpha * r, 3)) /
sqrt(
M_PI) *
exp(-alpha2 * r2)) * r5inv;
1067 if((alpha * r) < 0.5)
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);
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);
1079 g1 = (-
erf(alpha * r) + 2.0 * alpha * r /
sqrt(
M_PI) *
exp(-alpha2 * r2)) * r3inv;
1081 g2 = (3.0 *
erf(alpha * r) - (6.0 * alpha * r + 4.0 *
pow(alpha * r, 3)) /
sqrt(
M_PI) *
exp(-alpha2 * r2)) * r5inv;
1085 D2 += g2 * (dxyz % dxyz);
1091 for(
int nx = -nxmax; nx <= nxmax; nx++)
1092 for(
int ny = -nymax; ny <= nymax; ny++)
1093 for(
int nz = -nzmax; nz <= nzmax; nz++)
1095 if(nx != 0 || ny != 0 || nz != 0)
1100 double k2 = kx * kx + ky * ky + kz * kz;
1102 double kdotx = (x * kx + y * ky + z * kz);
1107 D2 += (val * kxyz) % kxyz;
1116 static int printed = 0;
1118 #ifdef GRAVITY_TALLBOX
1119 Terminate(
"GRAVITY_TALLBOX is not implemented for MULTIPOLE_ORDER >= 3");
1125 double alpha = 2.0 / leff;
1126 double alpha2 = alpha * alpha;
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);
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);
1138 mpi_printf(
"EWALD: D3 table: qxmax=%d qymax=%d qzmax=%d nxmax=%d nymax=%d nzmax=%d\n", qxmax, qymax, qzmax, nxmax, nymax,
1143 for(
int nx = -qxmax; nx <= qxmax; nx++)
1144 for(
int ny = -qymax; ny <= qymax; ny++)
1145 for(
int nz = -qzmax; nz <= qzmax; nz++)
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);
1153 double r2 = dx * dx + dy * dy + dz * dz;
1154 double r =
sqrt(r2);
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;
1165 if(nx != 0 || ny != 0 || nz != 0)
1167 g2 = -(3.0 *
erfc(alpha * r) + (6.0 * alpha * r + 4.0 *
pow(alpha * r, 3)) /
sqrt(
M_PI) *
exp(-alpha2 * r2)) * r5inv;
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)) *
1192 if((alpha * r) < 0.5)
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);
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);
1204 g2 = (3.0 *
erf(alpha * r) - (6.0 * alpha * r + 4.0 *
pow(alpha * r, 3)) /
sqrt(
M_PI) *
exp(-alpha2 * r2)) * r5inv;
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)) *
1218 for(
int nx = -nxmax; nx <= nxmax; nx++)
1219 for(
int ny = -nymax; ny <= nymax; ny++)
1220 for(
int nz = -nzmax; nz <= nzmax; nz++)
1222 if(nx != 0 || ny != 0 || nz != 0)
1227 double k2 = kx * kx + ky * ky + kz * kz;
1229 double kdotx = (x * kx + y * ky + z * kz);
1234 D3 += (val * kxyz) % (kxyz % kxyz);
1243 static int printed = 0;
1245 #ifdef GRAVITY_TALLBOX
1246 Terminate(
"GRAVITY_TALLBOX is not implemented for MULTIPOLE_ORDER >= 4");
1252 double alpha = 2.0 / leff;
1253 double alpha2 = alpha * alpha;
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);
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);
1265 mpi_printf(
"EWALD: D4 table: qxmax=%d qymax=%d qzmax=%d nxmax=%d nymax=%d nzmax=%d\n", qxmax, qymax, qzmax, nxmax, nymax,
1270 for(
int nx = -qxmax; nx <= qxmax; nx++)
1271 for(
int ny = -qymax; ny <= qymax; ny++)
1272 for(
int nz = -qzmax; nz <= qzmax; nz++)
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);
1280 double r2 = dx * dx + dy * dy + dz * dz;
1281 double r =
sqrt(r2);
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;
1293 if(nx != 0 || ny != 0 || nz != 0)
1295 g2 = -(3.0 *
erfc(alpha * r) + (6.0 * alpha * r + 4.0 *
pow(alpha * r, 3)) /
sqrt(
M_PI) *
exp(-alpha2 * r2)) * r5inv;
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)) *
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)) /
1326 if((alpha * r) < 0.5)
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);
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);
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);
1342 g2 = (3.0 *
erf(alpha * r) - (6.0 * alpha * r + 4.0 *
pow(alpha * r, 3)) /
sqrt(
M_PI) *
exp(-alpha2 * r2)) * r5inv;
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)) *
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)) /
1359 setup_D4(
ADD, D4, dxyz, aux2, aux3, aux4, g2, g3, g4);
1362 for(
int nx = -nxmax; nx <= nxmax; nx++)
1363 for(
int ny = -nymax; ny <= nymax; ny++)
1364 for(
int nz = -nzmax; nz <= nzmax; nz++)
1366 if(nx != 0 || ny != 0 || nz != 0)
1371 double k2 = kx * kx + ky * ky + kz * kz;
1373 double kdotx = (x * kx + y * ky + z * kz);
1378 D4 += (val * kxyz) % ((kxyz % (kxyz % kxyz)));
1387 static int printed = 0;
1389 #ifdef GRAVITY_TALLBOX
1390 Terminate(
"GRAVITY_TALLBOX is not implemented for MULTIPOLE_ORDER >= 4");
1396 double alpha = 2.0 / leff;
1397 double alpha2 = alpha * alpha;
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);
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);
1409 mpi_printf(
"EWALD: D5 table: qxmax=%d qymax=%d qzmax=%d nxmax=%d nymax=%d nzmax=%d\n", qxmax, qymax, qzmax, nxmax, nymax,
1414 for(
int nx = -qxmax; nx <= qxmax; nx++)
1415 for(
int ny = -qymax; ny <= qymax; ny++)
1416 for(
int nz = -qzmax; nz <= qzmax; nz++)
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);
1424 double r2 = dx * dx + dy * dy + dz * dz;
1425 double r =
sqrt(r2);
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;
1438 if(nx != 0 || ny != 0 || nz != 0)
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)) *
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)) /
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)) /
1474 if((alpha * r) < 0.5)
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);
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);
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);
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)) *
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)) /
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)) /
1510 setup_D5(
ADD, D5, dxyz, aux3, aux4, aux5, g3, g4, g5);
1513 for(
int nx = -nxmax; nx <= nxmax; nx++)
1514 for(
int ny = -nymax; ny <= nymax; ny++)
1515 for(
int nz = -nzmax; nz <= nzmax; nz++)
1520 double k2 = kx * kx + ky * ky + kz * kz;
1524 double kdotx = (x * kx + y * ky + z * kz);
1529 D5 += (val * kxyz) % (kxyz % ((kxyz % (kxyz % kxyz))));
1538 static int printed = 0;
1540 #ifdef GRAVITY_TALLBOX
1541 Terminate(
"GRAVITY_TALLBOX is not implemented for MULTIPOLE_ORDER >= 4");
1547 double alpha = 2.0 / leff;
1548 double alpha2 = alpha * alpha;
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);
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);
1560 mpi_printf(
"EWALD: D6 table: qxmax=%d qymax=%d qzmax=%d nxmax=%d nymax=%d nzmax=%d\n", qxmax, qymax, qzmax, nxmax, nymax,
1565 for(
int nx = -qxmax; nx <= qxmax; nx++)
1566 for(
int ny = -qymax; ny <= qymax; ny++)
1567 for(
int nz = -qzmax; nz <= qzmax; nz++)
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);
1575 double r2 = dx * dx + dy * dy + dz * dz;
1576 double r =
sqrt(r2);
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;
1588 double g3, g4, g5, g6;
1590 if(nx != 0 || ny != 0 || nz != 0)
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)) *
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)) /
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)) /
1606 g6 = (-10395.0 *
erfc(alpha * r) -
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)) /
1633 if((alpha * r) < 0.5)
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);
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);
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);
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);
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)) *
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)) /
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)) /
1667 g6 = (10395.0 *
erf(alpha * r) -
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)) /
1679 for(
int nx = -nxmax; nx <= nxmax; nx++)
1680 for(
int ny = -nymax; ny <= nymax; ny++)
1681 for(
int nz = -nzmax; nz <= nzmax; nz++)
1686 double k2 = kx * kx + ky * ky + kz * kz;
1690 double kdotx = (x * kx + y * ky + z * kz);
1695 D6 += (val * kxyz) % (kxyz % (kxyz % ((kxyz % (kxyz % kxyz)))));
1704 static int printed = 0;
1706 #ifdef GRAVITY_TALLBOX
1707 Terminate(
"GRAVITY_TALLBOX is not implemented for MULTIPOLE_ORDER >= 4");
1713 double alpha = 2.0 / leff;
1714 double alpha2 = alpha * alpha;
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);
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);
1726 mpi_printf(
"EWALD: D7 table: qxmax=%d qymax=%d qzmax=%d nxmax=%d nymax=%d nzmax=%d\n", qxmax, qymax, qzmax, nxmax, nymax,
1731 for(
int nx = -qxmax; nx <= qxmax; nx++)
1732 for(
int ny = -qymax; ny <= qymax; ny++)
1733 for(
int nz = -qzmax; nz <= qzmax; nz++)
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);
1741 double r2 = dx * dx + dy * dy + dz * dz;
1742 double r =
sqrt(r2);
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;
1755 double g4, g5, g6, g7;
1757 if(nx != 0 || ny != 0 || nz != 0)
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)) /
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)) /
1769 g6 = (-10395.0 *
erfc(alpha * r) -
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)) /
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)) /
1804 if((alpha * r) < 0.5)
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);
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);
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);
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);
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)) /
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)) /
1834 g6 = (10395.0 *
erf(alpha * r) -
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)) /
1841 g7 = (-135135.0 *
erf(alpha * r) +
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)) /
1854 for(
int nx = -nxmax; nx <= nxmax; nx++)
1855 for(
int ny = -nymax; ny <= nymax; ny++)
1856 for(
int nz = -nzmax; nz <= nzmax; nz++)
1861 double k2 = kx * kx + ky * ky + kz * kz;
1865 double kdotx = (x * kx + y * ky + z * kz);
1870 D7 += (val * kxyz) % (kxyz % (kxyz % (kxyz % (kxyz % (kxyz % kxyz)))));
global_data_all_processes All
symtensor3< double > ewald_D3(double x, double y, double z)
double ewald_D0(double x, double y, double z)
This function computes the potential correction term by means of Ewald summation.
symtensor4< double > ewald_D4(double x, double y, double z)
void ewald_init(void)
This function initializes tables with the correction force and the correction potential due to the pe...
symtensor6< double > ewald_D6(double x, double y, double z)
symtensor7< double > ewald_D7(double x, double y, double z)
symtensor5< double > ewald_D5(double x, double y, double z)
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...
symtensor2< double > ewald_D2(double x, double y, double z)
void ewald_gridlookup(const MyIntPosType *p_intpos, const MyIntPosType *target_intpos, enum interpolate_options flag, ewald_data &fper)
void nearest_image_intpos_to_pos(const MyIntPosType *const a, const MyIntPosType *const b, T *posdiff)
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,...)
void ** SharedMemBaseAddr
#define BITS_FOR_POSITIONS
#define EN
Base dimension of cubical Ewald lookup table, for one octant.
#define HIGHEST_NEEDEDORDER_EWALD_DPHI
#define EWALD_TAYLOR_ORDER
expr pow(half base, half exp)
symtensor3< MyReal > D3phi
symtensor2< MyReal > D2phi
void setup_D6(enum setup_options opt, symtensor6< T > &D6, vector< T > &dxyz, TypeGfac g3, TypeGfac g4, TypeGfac g5, TypeGfac g6)
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)
void setup_D3(enum setup_options opt, symtensor3< T > &D3, vector< T > &dxyz, symtensor2< T > &aux2, symtensor3< T > &aux3, TypeGfac g2, TypeGfac g3)
void setup_D7(enum setup_options opt, symtensor7< T > &D7, vector< T > &dxyz, TypeGfac g4, TypeGfac g5, TypeGfac g6, TypeGfac g7)
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)
void subdivide_evenly(long long N, int pieces, int index_bin, long long *first, int *count)
void myflush(FILE *fstream)