GADGET-4
grav_forcetest.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#ifdef FORCETEST
15
16#include <math.h>
17#include <mpi.h>
18#include <stdio.h>
19#include <stdlib.h>
20#include <string.h>
21
22#include "../data/allvars.h"
23#include "../data/dtypes.h"
24#include "../data/intposconvert.h"
25#include "../data/mymalloc.h"
26#include "../domain/domain.h"
27#include "../gravity/ewald.h"
28#include "../gravity/grav_forcetest.h"
29#include "../gravtree/gravtree.h"
30#include "../logs/logs.h"
31#include "../logs/timer.h"
32#include "../main/main.h"
33#include "../main/simulation.h"
34#include "../mpi_utils/generic_comm.h"
35#include "../mpi_utils/mpi_utils.h"
36#include "../ngbtree/ngbtree.h"
37#include "../pm/pm.h"
38#include "../system/system.h"
39#include "../time_integration/timestep.h"
40
41/* This function computes the gravitational forces for all active particles.
42 * A new tree is constructed, if the number of force computations since
43 * it's last construction exceeds some fraction of the total
44 * particle number, otherwise tree nodes are dynamically updated if needed.
45 */
46
47/* FORCETEST_TESTFORCELAW=1 for special test to check force law for TreePM
48 * FORCETEST_TESTFORCELAW=2 for special test to check force law for TreePM+PLACEHIGHRESREGION
49 */
50
51/* local data structure for collecting particle/cell data that is sent to other processors if needed */
52struct frctest_in : data_in_generic
53{
54 MyIntPosType IntPos[3];
55 unsigned char Type;
56#if NSOFTCLASSES > 1
57 unsigned char SofteningClass;
58#endif
59};
60
61/* local data structure that holds results acquired on remote processors */
62struct frctest_out
63{
64 double Acc[3];
65 double Pot;
66 double DistToID1;
67#if defined(PMGRID) && defined(PERIODIC) && !defined(TREEPM_NOTIMESPLIT)
68 double AccShortRange[3];
69 double PotShortRange;
70#ifdef PLACEHIGHRESREGION
71 double AccVeryShortRange[3];
72 double PotVeryShortRange;
73#endif
74#endif
75};
76
78
79class frctest_comm : public my_comm
80{
81 public:
82 // need to implement a constructor that calls the constructor of the base class
83 frctest_comm(domain<simparticles> *dptr, gravtree<simparticles> *tptr, simparticles *pptr, gravtest *gtptr)
84 : my_comm(dptr, tptr, pptr)
85 {
86 }
87
88 using my_comm::D;
89 using my_comm::Thread;
90 using my_comm::Tp;
91 using my_comm::Tree;
92
93 void particle2in(frctest_in *in, int i) override
94 {
95 for(int j = 0; j < 3; j++)
96 in->IntPos[j] = Tp->P[i].IntPos[j];
97
98 in->Type = Tp->P[i].getType();
99#if NSOFTCLASSES > 1
100 in->SofteningClass = Tp->P[i].getSofteningClass();
101#endif
102 }
103
104 void out2particle(frctest_out *out, int i, int mode) override
105 {
106 if(mode == MODE_LOCAL_PARTICLES) /* initial store */
107 {
108 Tp->P[i].GravAccelDirect[0] = out->Acc[0];
109 Tp->P[i].GravAccelDirect[1] = out->Acc[1];
110 Tp->P[i].GravAccelDirect[2] = out->Acc[2];
111 Tp->P[i].PotentialDirect = out->Pot;
112
113 Tp->P[i].DistToID1 = out->DistToID1;
114#if defined(PMGRID) && defined(PERIODIC) && !defined(TREEPM_NOTIMESPLIT)
115 Tp->P[i].GravAccelShortRange[0] = out->AccShortRange[0];
116 Tp->P[i].GravAccelShortRange[1] = out->AccShortRange[1];
117 Tp->P[i].GravAccelShortRange[2] = out->AccShortRange[2];
118 Tp->P[i].PotentialShortRange = out->PotShortRange;
119#ifdef PLACEHIGHRESREGION
120 Tp->P[i].GravAccelVeryShortRange[0] = out->AccVeryShortRange[0];
121 Tp->P[i].GravAccelVeryShortRange[1] = out->AccVeryShortRange[1];
122 Tp->P[i].GravAccelVeryShortRange[2] = out->AccVeryShortRange[2];
123 Tp->P[i].PotentialVeryShortRange = out->PotVeryShortRange;
124#endif
125#endif
126 }
127 else /* combine */
128 {
129 Tp->P[i].GravAccelDirect[0] += out->Acc[0];
130 Tp->P[i].GravAccelDirect[1] += out->Acc[1];
131 Tp->P[i].GravAccelDirect[2] += out->Acc[2];
132 Tp->P[i].PotentialDirect += out->Pot;
133 if(out->DistToID1 > 0)
134 Tp->P[i].DistToID1 = out->DistToID1;
135#if defined(PMGRID) && defined(PERIODIC) && !defined(TREEPM_NOTIMESPLIT)
136 Tp->P[i].GravAccelShortRange[0] += out->AccShortRange[0];
137 Tp->P[i].GravAccelShortRange[1] += out->AccShortRange[1];
138 Tp->P[i].GravAccelShortRange[2] += out->AccShortRange[2];
139 Tp->P[i].PotentialShortRange += out->PotShortRange;
140#ifdef PLACEHIGHRESREGION
141 Tp->P[i].GravAccelVeryShortRange[0] += out->AccVeryShortRange[0];
142 Tp->P[i].GravAccelVeryShortRange[1] += out->AccVeryShortRange[1];
143 Tp->P[i].GravAccelVeryShortRange[2] += out->AccVeryShortRange[2];
144 Tp->P[i].PotentialVeryShortRange += out->PotVeryShortRange;
145#endif
146#endif
147 }
148 }
149
150 int evaluate(int target, int mode, int thread_id, int action, frctest_in *in, int numnodes, node_info *firstnode,
151 frctest_out &out) override
152 {
153 /* make sure that the particle is exported to all other tasks exactly once */
154 if(mode == MODE_LOCAL_PARTICLES)
155 {
156 for(int n = 0; n < this->D->NTopleaves; n++)
157 {
158 int task = D->TaskOfLeaf[n];
159
160 if(task == D->ThisTask)
161 continue;
162
163 if(Thread.Exportflag[task] == target)
164 continue;
165
166 int no = n + this->Tree->MaxPart + this->Tree->MaxNodes; /* a pseudo node for this task */
167
168 this->Tree->tree_export_node_threads(no, target, &Thread);
169 }
170 }
171
172 MyIntPosType *intpos = in->IntPos;
173
174 out.Pot = 0;
175
176#if defined(PMGRID) && defined(PERIODIC) && !defined(TREEPM_NOTIMESPLIT)
177 out.PotShortRange = 0;
178#ifdef PLACEHIGHRESREGION
179 out.PotVeryShortRange = 0;
180#endif
181#endif
182 for(int i = 0; i < 3; i++)
183 {
184 out.Acc[i] = 0;
185#if defined(PMGRID) && defined(PERIODIC) && !defined(TREEPM_NOTIMESPLIT)
186 out.AccShortRange[i] = 0;
187#ifdef PLACEHIGHRESREGION
188 out.AccVeryShortRange[i] = 0;
189#endif
190#endif
191 }
192
193 double disttoid1 = 0;
194
195 for(int idx = 0; idx < Tp->nsource; idx++)
196 {
197 int j = Tp->indexlist[idx];
198
199 double hmax;
200#if NSOFTCLASSES > 1
201 double h_i = All.ForceSoftening[in->SofteningClass];
202 double h_j = All.ForceSoftening[Tp->P[j].getSofteningClass()];
203
204 if(h_j > h_i)
205 hmax = h_j;
206 else
207 hmax = h_i;
208#else
209 hmax = All.ForceSoftening[0];
210#endif
211
212 double dxyz[3];
213 Tp->nearest_image_intpos_to_pos(Tp->P[j].IntPos, intpos, dxyz); /* converts the integer distance to floating point */
214
215 double r2 = dxyz[0] * dxyz[0] + dxyz[1] * dxyz[1] + dxyz[2] * dxyz[2];
216
217 double mass = Tp->P[j].getMass();
218
219 /* now evaluate the multipole moment */
220
221 double r = sqrt(r2);
222
223 if(Tp->P[j].ID.get() == 1)
224 disttoid1 = r;
225
226 /* we compute 3 different forces:
227 * (1) The correct direct summation force, if needed with Ewald correction: ftrue
228 * In the case of PM:
229 * (2) The short range direct summation force with only the erfc cut-off (this is what the tree can at best deliver): fsr
230 * (3) The expected PM force based on the long-range part of the Ewald sum. This is equal to ftrue - fsr - fsfr_periodic_images
231 * */
232
233 double wp_newton, fac_newton;
234
235 if(r > 0)
236 {
237 fac_newton = mass / (r2 * r);
238 wp_newton = -mass / r;
239 }
240 else
241 {
242 fac_newton = 0;
243 wp_newton = 0;
244 }
245
246 double fac, wp;
247
248 if(r >= hmax)
249 {
250 fac = fac_newton;
251 wp = wp_newton;
252 }
253 else
254 {
255 double h_inv = 1.0 / hmax;
256 double h3_inv = h_inv * h_inv * h_inv;
257 double u = r * h_inv;
258
259 if(u < 0.5)
260 {
261 double u2 = u * u;
262 fac = mass * h3_inv * (SOFTFAC1 + u2 * (SOFTFAC2 * u + SOFTFAC3));
263 wp = mass * h_inv * (SOFTFAC4 + u2 * (SOFTFAC5 + u2 * (SOFTFAC6 * u + SOFTFAC7)));
264 }
265 else
266 {
267 double u2 = u * u, u3 = u2 * u;
268 fac = mass * h3_inv * (SOFTFAC8 + SOFTFAC9 * u + SOFTFAC10 * u2 + SOFTFAC11 * u3 + SOFTFAC12 / u3);
269 wp = mass * h_inv * (SOFTFAC13 + SOFTFAC14 / u + u2 * (SOFTFAC1 + u * (SOFTFAC15 + u * (SOFTFAC16 + SOFTFAC17 * u))));
270 }
271 }
272
273 // The Newtonian force is: fac * dxyz
274 // The Newtonian potential is: wp
275
276#if defined(PMGRID) && defined(PERIODIC) && !defined(TREEPM_NOTIMESPLIT)
277 {
278 double asmth = Tp->Asmth[0];
279 double u = 0.5 / asmth * r;
280
281 double factor_force = (erfc(u) + 2.0 * u / sqrt(M_PI) * exp(-u * u) - 1.0);
282 double factor_pot = erfc(u);
283
284 double facs = fac + fac_newton * factor_force;
285 double wps = wp + (r > 0 ? wp_newton * (factor_pot - 1.0) : mass / (asmth * sqrt(M_PI)));
286
287 double acc_short_x = dxyz[0] * facs;
288 double acc_short_y = dxyz[1] * facs;
289 double acc_short_z = dxyz[2] * facs;
290
291#ifndef GRAVITY_TALLBOX
292 double alpha = 0.5 / asmth;
293 double pot_short =
294 wps + mass * M_PI / (alpha * alpha * All.BoxSize * All.BoxSize * All.BoxSize) * (LONG_X * LONG_Y * LONG_Z);
295#else
296 double pot_short = wps; /* the constant potential term is here computed as part of the long-range force and not the
297 short-range force, unlike in the ordinary periodic case */
298#endif
299 out.AccShortRange[0] += acc_short_x;
300 out.AccShortRange[1] += acc_short_y;
301 out.AccShortRange[2] += acc_short_z;
302 out.PotShortRange += pot_short;
303
304#ifdef PLACEHIGHRESREGION
305 if(Tp->check_high_res_point_location(Tp->P[j].IntPos) == FLAG_INSIDE &&
306 Tp->check_high_res_point_location(intpos) == FLAG_INSIDE)
307 {
308 double asmth = Tp->Asmth[1];
309 double u = 0.5 / asmth * r;
310
311 double factor_force = (erfc(u) + 2.0 * u / sqrt(M_PI) * exp(-u * u) - 1.0);
312 double factor_pot = erfc(u);
313
314 double facs = fac + fac_newton * factor_force;
315 double wps = wp + (r > 0 ? wp_newton * (factor_pot - 1.0) : mass / (asmth * sqrt(M_PI)));
316
317 double alpha = 0.5 / asmth;
318
319 acc_short_x = dxyz[0] * facs;
320 acc_short_y = dxyz[1] * facs;
321 acc_short_z = dxyz[2] * facs;
322
323 pot_short = wps + mass * M_PI / (alpha * alpha * All.BoxSize * All.BoxSize * All.BoxSize) * (LONG_X * LONG_Y * LONG_Z);
324 }
325
326 out.AccVeryShortRange[0] += acc_short_x;
327 out.AccVeryShortRange[1] += acc_short_y;
328 out.AccVeryShortRange[2] += acc_short_z;
329 out.PotVeryShortRange += pot_short;
330#endif
331 }
332#endif
333
334 // the direct force is the ordinary Newtonian force
335 out.Acc[0] += fac * dxyz[0];
336 out.Acc[1] += fac * dxyz[1];
337 out.Acc[2] += fac * dxyz[2];
338 out.Pot += wp;
339
340#ifdef PERIODIC
341 // and in the periodic case, we add in the correction potential/force
342 ewald_data ew;
343 Ewald.ewald_gridlookup(Tp->P[j].IntPos, intpos, ewald::POINTMASS, ew);
344
345 out.Pot += mass * ew.D0phi;
346 out.Acc[0] += mass * ew.D1phi[0];
347 out.Acc[1] += mass * ew.D1phi[1];
348 out.Acc[2] += mass * ew.D1phi[2];
349#endif
350 }
351
352 out.DistToID1 = disttoid1;
353
354 return 0;
355 }
356};
357
358void gravtest::gravity_forcetest(int timebin)
359{
360 TIMER_START(CPU_FORCETEST);
361
362 int *TargetList = (int *)Mem.mymalloc("TargetList", Sp->NumPart * sizeof(int));
363 int nloc = 0;
364
365 particle_data *P = Sp->P;
366
367 // create a random selection of target particles for which we compute direct summation forces
368 for(int idx = 0; idx < Sp->TimeBinsGravity.NActiveParticles; idx++)
369 {
370 int target = Sp->TimeBinsGravity.ActiveParticleList[idx];
371
372 if(target < 0)
373 continue;
374
375#ifdef FORCETEST_FIXEDPARTICLESET
376 if(All.NumCurrentTiStep == 0)
377 {
378 if(get_random_number() < FORCETEST)
379 P[target].SelectedFlag = true;
380 else
381 P[target].SelectedFlag = false;
382 }
383 if(P[target].SelectedFlag)
384 TargetList[nloc++] = target;
385#else
386 if(get_random_number() < FORCETEST)
387 TargetList[nloc++] = target;
388#endif
389 }
390
391 long long ntot;
392 sumup_large_ints(1, &nloc, &ntot, D->Communicator);
393
394 /* we pull put a separate list of the particles with non-zero masses to accelerate some of our special tests
395 * where there are potentially many particles with zero masses
396 */
397
398 Sp->nsource = 0;
399 Sp->indexlist = (int *)Mem.mymalloc("indexlist", Sp->NumPart * sizeof(int));
400
401#ifdef HIERARCHICAL_GRAVITY
402 int numidx = Sp->TimeBinsGravity.NActiveParticles;
403#else
404 int numidx = Sp->NumPart;
405#endif
406
407 for(int idx = 0; idx < numidx; idx++)
408 {
409#ifdef HIERARCHICAL_GRAVITY
410 int target = Sp->TimeBinsGravity.ActiveParticleList[idx];
411
412 if(target < 0)
413 continue;
414#else
415 int target = idx;
416#endif
417
418 if(P[target].getMass() == 0)
419 continue;
420
421 Sp->indexlist[Sp->nsource++] = target;
422 }
423
424 D->mpi_printf("FORCETEST: Testing forces for %lld particles out of %lld active ones.\n", ntot,
426
427 frctest_comm commpattern(this->D, this->GravTree, this->Sp, this);
428
429 double t0 = Logs.second();
430
431 commpattern.execute(nloc, TargetList, MODE_DEFAULT);
432
433 double t1 = Logs.second();
434 double maxt = Logs.timediff(t0, t1);
435
436 D->mpi_printf("FORCETEST: Testing forces took %g sec.\n", maxt);
437
438 Mem.myfree(Sp->indexlist);
439 Sp->indexlist = NULL;
440
441 /* now add things for comoving integration */
442
444 {
445#ifndef PERIODIC
446 double fac1 = 0.5 * All.Hubble * All.Hubble * All.Omega0 / All.G;
447
448 for(int idx = 0; idx < nloc; idx++)
449 {
450 int i = TargetList[idx];
451
452 double pos[3];
453 Sp->intpos_to_pos(Sp->P[i].IntPos, pos);
454
455 for(int j = 0; j < 3; j++)
456 Sp->P[i].GravAccelDirect[j] += fac1 * pos[j];
457 }
458#endif
459 }
460
461 /* muliply by G */
462 for(int idx = 0; idx < nloc; idx++)
463 {
464 int i = TargetList[idx];
465
466 for(int j = 0; j < 3; j++)
467 {
468 Sp->P[i].GravAccelDirect[j] *= All.G;
469#if defined(PMGRID) && defined(PERIODIC) && !defined(TREEPM_NOTIMESPLIT)
470 Sp->P[i].GravAccelShortRange[j] *= All.G;
471#ifdef PLACEHIGHRESREGION
472 Sp->P[i].GravAccelVeryShortRange[j] *= All.G;
473#endif
474#endif
475#if(FORCETEST_TESTFORCELAW == 3)
476 Sp->P[i].GravAccelDirectTest[j] *= All.G;
477#endif
478 }
479
480#if NSOFTCLASSES > 1
481 double selfpot = Sp->P[i].getMass() / (All.ForceSoftening[Sp->P[i].getSofteningClass()] / 2.8);
482#else
483 double selfpot = Sp->P[i].getMass() / (All.ForceSoftening[0] / 2.8);
484#endif
485
486 Sp->P[i].PotentialDirect += selfpot; /* remove self-potential */
487 Sp->P[i].PotentialDirect *= All.G;
488
489#if defined(PMGRID) && defined(PERIODIC) && !defined(TREEPM_NOTIMESPLIT)
490 Sp->P[i].PotentialShortRange += selfpot; /* remove self-potential */
491 Sp->P[i].PotentialShortRange *= All.G;
492
493#ifdef PLACEHIGHRESREGION
494 Sp->P[i].PotentialVeryShortRange += selfpot; /* remove self-potential */
495 Sp->P[i].PotentialVeryShortRange *= All.G;
496#endif
497#endif
498#if(FORCETEST_TESTFORCELAW == 3)
499 Sp->P[i].PotentialDirectTest += selfpot; /* remove self-potential */
500 Sp->P[i].PotentialDirectTest *= All.G;
501#endif
502 }
503
504 /* Finally, the following factor allows a computation of cosmological simulation
505 with vacuum energy in physical coordinates */
506
508 {
509 double fac1 = All.OmegaLambda * All.Hubble * All.Hubble;
510
511 for(int idx = 0; idx < nloc; idx++)
512 {
513 int i = TargetList[idx];
514
515 double pos[3];
516 Sp->intpos_to_pos(Sp->P[i].IntPos, pos);
517
518 for(int j = 0; j < 3; j++)
519 Sp->P[i].GravAccelDirect[j] += fac1 * pos[j];
520 }
521 }
522
523 /* now output the forces to a file */
524
525 int *nloc_tab = (int *)Mem.mymalloc("nloc_tab", D->NTask * sizeof(int));
526 MPI_Allgather(&nloc, 1, MPI_INT, nloc_tab, 1, MPI_INT, D->Communicator);
527
528 for(int nthis = 0; nthis < D->NTask; nthis++)
529 {
530 if(nloc_tab[nthis] > 0)
531 {
532 if(nthis == D->ThisTask)
533 {
534 char buf[MAXLEN_PATH_EXTRA];
535 sprintf(buf, "%s%s", All.OutputDir, "forcetest.txt");
536
537 if(!(Logs.FdForceTest = fopen(buf, "a")))
538 Terminate("error in opening file '%s'\n", buf);
539
540 for(int idx = 0; idx < nloc; idx++)
541 {
542 int i = TargetList[idx];
543
544 double pos[3];
545 Sp->intpos_to_pos(Sp->P[i].IntPos, pos);
546
547#if defined(PMGRID) && defined(PERIODIC) && !defined(TREEPM_NOTIMESPLIT)
548
549#ifdef PLACEHIGHRESREGION
550 int flaginside = Sp->check_high_res_point_location(P[i].IntPos);
551 fprintf(
552 Logs.FdForceTest,
553 "%d %d %lld %g %17.12g %17.12g %17.12g %17.12g %12.12g %17.12g %17.12g %17.12g %17.12g %17.12g %17.12g "
554 "%17.12g %17.12g %17.12g %17.12g %17.12g %17.12g %17.12g %17.12g %17.12g "
555 "%17.12g %17.12g %17.12g %17.12g %17.12g %17.12g %17.12g %17.12g %17.12g\n",
556 P[i].getType(), flaginside, (long long)P[i].ID.get(), All.Time, pos[0], pos[1], pos[2], P[i].DistToID1,
557 P[i].GravAccelDirect[0], P[i].GravAccelDirect[1], P[i].GravAccelDirect[2], P[i].GravAccelShortRange[0],
558 P[i].GravAccelShortRange[1], P[i].GravAccelShortRange[2], P[i].GravAccelVeryShortRange[0],
559 P[i].GravAccelVeryShortRange[1], P[i].GravAccelVeryShortRange[2], P[i].GravAccel[0], P[i].GravAccel[1],
560 P[i].GravAccel[2], P[i].GravPM[0], P[i].GravPM[1], P[i].GravPM[2], P[i].GravAccelHPM[0], P[i].GravAccelHPM[1],
561 P[i].GravAccelHPM[2], P[i].PotentialDirect, P[i].PotentialShortRange, P[i].PotentialVeryShortRange,
562 P[i].Potential, P[i].PM_Potential, P[i].PotentialHPM, Sp->Asmth[1]);
563#else
564 fprintf(
565 Logs.FdForceTest,
566 "%d %d %lld %g %17.12g %17.12g %17.12g %17.12g %17.12g %17.12g %17.12g %17.12g %17.12g %17.12g %17.12g "
567 "%17.12g "
568 "%17.12g %17.12g %17.12g %17.12g %17.12g %17.12g %17.12g %17.12g\n",
569 P[i].getType(), timebin, (long long)P[i].ID.get(), All.Time, pos[0], pos[1], pos[2], P[i].DistToID1,
570 P[i].GravAccelDirect[0], P[i].GravAccelDirect[1], P[i].GravAccelDirect[2], P[i].GravAccelShortRange[0],
571 P[i].GravAccelShortRange[1], P[i].GravAccelShortRange[2], P[i].GravAccel[0], P[i].GravAccel[1],
572 P[i].GravAccel[2], P[i].GravPM[0], P[i].GravPM[1], P[i].GravPM[2], P[i].PotentialDirect,
573 P[i].PotentialShortRange, P[i].Potential, P[i].PM_Potential);
574#endif
575#else
576 fprintf(Logs.FdForceTest,
577 "%d %d %lld %g %17.12g %17.12g %17.12g %17.12g %17.12g %17.12g %17.12g %17.12g %17.12g %17.12g %17.12g "
578 "%17.12g\n",
579 P[i].getType(), timebin, (long long)P[i].ID.get(), All.Time, pos[0], pos[1], pos[2], P[i].DistToID1,
580 P[i].GravAccelDirect[0], P[i].GravAccelDirect[1], P[i].GravAccelDirect[2], Sp->P[i].GravAccel[0],
581 Sp->P[i].GravAccel[1], Sp->P[i].GravAccel[2], P[i].Potential, P[i].PotentialDirect);
582#endif
583 }
584
585 fclose(Logs.FdForceTest);
586 }
587
588 MPI_Barrier(D->Communicator);
589 }
590 }
591 Mem.myfree(nloc_tab);
592
593 /* Now the force computation is finished */
594 if(D->ThisTask == 0)
595 {
596 double costtotal = Sp->NumPart * ntot;
597
598 fprintf(Logs.FdTimings, "DIRECT Nf= %lld step=%d timebin=%d part/sec=%g | %g ia/part=%g\n\n", ntot, All.NumCurrentTiStep,
599 timebin, ((double)ntot) / (Sp->NTask * maxt + 1.0e-20), ntot / ((maxt + 1.0e-20) * Sp->NTask),
600 ((double)(costtotal)) / (ntot + 1.0e-20));
601
603 }
604
605 Mem.myfree(TargetList);
606
607 TIMER_STOP(CPU_FORCETEST);
608}
609
610#ifdef FORCETEST_TESTFORCELAW /* in this option we assume NSOFTCLASSES >= 2 for FORCETEST_TESTFORCELAW=1, and NSOFTCLASSES >= 3 for \
611 FORCETEST_TESTFORCELAW=2*/
612void sim::gravity_forcetest_testforcelaw(void)
613{
614 int Ncycles = 40;
615 double xyz[6], eps;
616
619
620 for(int cycle = 0; cycle < Ncycles; cycle++)
621 {
622 Domain.mpi_printf("\nTEST-FORCE-LAW: cycle=%d|%d ----------------------------------\n\n", cycle, Ncycles);
623
624 double epsloc = 0, xyzloc[6] = {0, 0, 0, 0, 0, 0};
625
626 /* set particle with ID=1 to new random coordinate in box */
627 for(int n = 0; n < Sp.NumPart; n++)
628 {
629 Sp.P[n].setType(1);
630 Sp.P[n].setSofteningClass(1);
631
632 if(Sp.P[n].ID.get() == 1)
633 {
634 xyzloc[0] = All.BoxSize / LONG_X * get_random_number();
635 xyzloc[1] = All.BoxSize / LONG_Y * get_random_number();
636 xyzloc[2] = All.BoxSize / LONG_Z * get_random_number();
637
638 for(int i = 0; i < 3; i++)
639 xyzloc[3 + i] = xyzloc[i];
640
641#if defined(PLACEHIGHRESREGION) && (FORCETEST_TESTFORCELAW == 2)
642 if(get_random_number() < 0.5)
643 {
644 xyzloc[3] = All.BoxSize / LONG_X * get_random_number();
645 xyzloc[4] = All.BoxSize / LONG_Y * get_random_number();
646 xyzloc[5] = All.BoxSize / LONG_Z * get_random_number();
647 }
648#endif
649 Sp.pos_to_intpos(&xyzloc[3], Sp.P[n].IntPos);
650
651 epsloc = All.ForceSoftening[Sp.P[n].getSofteningClass()];
652 }
653 }
654
655 MPI_Allreduce(xyzloc, xyz, 6, MPI_DOUBLE, MPI_SUM, Communicator);
656 MPI_Allreduce(&epsloc, &eps, 1, MPI_DOUBLE, MPI_SUM, Communicator);
657
658 double rmin = 0.01 * eps;
659 double rmax = sqrt(pow(0.5 * All.BoxSize / LONG_X, 2) + pow(0.5 * All.BoxSize / LONG_Y, 2) + pow(0.5 * All.BoxSize / LONG_Z, 2));
660
661 for(int n = 0; n < Sp.NumPart; n++)
662 {
663 if(Sp.P[n].ID.get() != 1)
664 {
665 double r = exp(log(rmin) + (log(rmax) - log(rmin)) * get_random_number());
666 double theta = acos(2 * get_random_number() - 1);
667 double phi = 2 * M_PI * get_random_number();
668
669#if defined(PLACEHIGHRESREGION) && (FORCETEST_TESTFORCELAW == 2)
670 if(get_random_number() < 0.5)
671 {
672 r = exp(log(rmin) + (log(rmax / 100) - log(rmin)) * get_random_number());
673 Sp.P[n].setType(2);
674 }
675 else
676 {
677 Sp.P[n].setType(3);
678 }
679#endif
680
681 double dx = r * sin(theta) * cos(phi);
682 double dy = r * sin(theta) * sin(phi);
683 double dz = r * cos(theta);
684
685 double pos[3];
686 pos[0] = xyz[0] + dx;
687 pos[1] = xyz[1] + dy;
688 pos[2] = xyz[2] + dz;
689
690#if defined(PLACEHIGHRESREGION) && (FORCETEST_TESTFORCELAW == 2)
691 if(Sp.P[n].getType() == 3)
692 {
693 pos[0] = xyz[3] + dx;
694 pos[1] = xyz[4] + dy;
695 pos[2] = xyz[5] + dz;
696 }
697#endif
698 Sp.pos_to_intpos(pos, Sp.P[n].IntPos);
700 }
701 }
702
703 Domain.domain_free();
704 Domain.domain_decomposition(STANDARD); /* do domain decomposition if needed */
705
706 /* allocate space for gravity accelerations */
707
709 for(int timebin = All.HighestSynchronizedTimeBin; timebin >= 0; timebin--)
710 {
711 for(int i = Sp.TimeBinsGravity.FirstInTimeBin[timebin]; i >= 0; i = Sp.TimeBinsGravity.NextInTimeBin[i])
713 }
714
715#if defined(PMGRID) && defined(PERIODIC) && !defined(TREEPM_NOTIMESPLIT)
717#endif
718
720 }
721
722 endrun();
723}
724#endif
725
726#endif
global_data_all_processes All
Definition: main.cc:40
MyIDType get(void) const
Definition: idstorage.h:87
Definition: domain.h:31
@ POINTMASS
Definition: ewald.h:67
void ewald_gridlookup(const MyIntPosType *p_intpos, const MyIntPosType *target_intpos, enum interpolate_options flag, ewald_data &fper)
Definition: ewald.cc:355
void gravity_forcetest(int timebin)
void pos_to_intpos(T *posdiff, MyIntPosType *intpos)
void intpos_to_pos(MyIntPosType *intpos, T *posdiff)
void constrain_intpos(MyIntPosType *pos)
Definition: intposconvert.h:68
double timediff(double t0, double t1)
Definition: logs.cc:488
FILE * FdTimings
Definition: logs.h:47
double second(void)
Definition: logs.cc:471
int NTask
Definition: setcomm.h:32
MPI_Comm Communicator
Definition: setcomm.h:31
void compute_grav_accelerations(int timebin)
This routine computes the gravitational accelerations for all active particles.
Definition: gravity.cc:41
domain< simparticles > Domain
Definition: simulation.h:58
void gravity_long_range_force(void)
sph NgbTree
Definition: simulation.h:68
void endrun(void)
This function aborts the simulations.
Definition: begrun.cc:395
simparticles Sp
Definition: simulation.h:56
TimeBinData TimeBinsGravity
Definition: simparticles.h:205
particle_data * P
Definition: simparticles.h:54
void mark_active_timebins(void)
Definition: predict.cc:155
void treefree(void)
Definition: tree.cc:1232
#define SOFTFAC3
Definition: constants.h:391
#define SOFTFAC13
Definition: constants.h:401
#define FLAG_INSIDE
Definition: constants.h:30
#define SOFTFAC16
Definition: constants.h:404
#define SOFTFAC4
Definition: constants.h:392
#define SOFTFAC14
Definition: constants.h:402
#define SOFTFAC10
Definition: constants.h:398
#define SOFTFAC11
Definition: constants.h:399
#define SOFTFAC8
Definition: constants.h:396
#define MAXLEN_PATH_EXTRA
Definition: constants.h:301
#define SOFTFAC2
Definition: constants.h:390
#define SOFTFAC15
Definition: constants.h:403
#define MODE_DEFAULT
Definition: constants.h:23
#define SOFTFAC1
Definition: constants.h:389
#define MODE_LOCAL_PARTICLES
Definition: constants.h:21
#define SOFTFAC12
Definition: constants.h:400
#define M_PI
Definition: constants.h:56
#define SOFTFAC9
Definition: constants.h:397
#define SOFTFAC5
Definition: constants.h:393
#define SOFTFAC6
Definition: constants.h:394
#define SOFTFAC17
Definition: constants.h:405
#define SOFTFAC7
Definition: constants.h:395
@ STANDARD
Definition: domain.h:23
#define LONG_Y
Definition: dtypes.h:369
uint32_t MyIntPosType
Definition: dtypes.h:35
#define LONG_X
Definition: dtypes.h:361
#define LONG_Z
Definition: dtypes.h:377
ewald Ewald
Definition: main.cc:42
logs Logs
Definition: main.cc:43
#define TIMER_START(counter)
Starts the timer counter.
Definition: logs.h:197
#define TIMER_STOP(counter)
Stops the timer counter.
Definition: logs.h:220
#define Terminate(...)
Definition: macros.h:15
void sumup_large_ints(int n, int *src, long long *res, MPI_Comm comm)
memory Mem
Definition: main.cc:44
expr acos(half arg)
Definition: half.hpp:2844
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 log(half arg)
Definition: half.hpp:2745
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
int NActiveParticles
Definition: timestep.h:18
int * ActiveParticleList
Definition: timestep.h:20
int FirstInTimeBin[TIMEBINS]
Definition: timestep.h:23
long long GlobalNActiveParticles
Definition: timestep.h:19
int * NextInTimeBin
Definition: timestep.h:25
vector< MyReal > D1phi
Definition: gravtree.h:186
MyReal D0phi
Definition: gravtree.h:185
char OutputDir[MAXLEN_PATH]
Definition: allvars.h:272
double ForceSoftening[NSOFTCLASSES+NSOFTCLASSES_HYDRO+2]
Definition: allvars.h:258
Definition: tree.h:77
void setSofteningClass(unsigned char softclass)
MyDouble getMass(void)
unsigned char getSofteningClass(void)
MyIDStorage ID
Definition: particle_data.h:70
unsigned char getType(void)
vector< MyFloat > GravAccel
Definition: particle_data.h:55
void setType(unsigned char type)
MyIntPosType IntPos[3]
Definition: particle_data.h:53
void myflush(FILE *fstream)
Definition: system.cc:125
double get_random_number(void)
Definition: system.cc:49