GADGET-4
density.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/intposconvert.h"
23 #include "../data/mymalloc.h"
24 #include "../logs/logs.h"
25 #include "../logs/timer.h"
26 #include "../main/simulation.h"
27 #include "../mpi_utils/mpi_utils.h"
28 #include "../ngbtree/ngbtree.h"
29 #include "../sort/cxxsort.h"
30 #include "../sph/kernel.h"
31 #include "../sph/sph.h"
32 #include "../system/system.h"
33 
39 /* This function checks whether there is a spatial overlap between the (rectangular) enclosing box
40  * of the particles contained in a node, and the search region.
41  */
42 inline int sph::sph_density_evaluate_particle_node_opening_criterion(pinfo &pdat, ngbnode *nop)
43 {
44  if(nop->level <= LEVEL_ALWAYS_OPEN) // always open the root node (note: full node length does not fit in the integer type)
45  return NODE_OPEN;
46 
47  if(nop->Ti_Current != All.Ti_Current)
48  nop->drift_node(All.Ti_Current, Tp);
49 
50  MyIntPosType left[3], right[3];
51 
52  left[0] = Tp->nearest_image_intpos_to_intpos_X(nop->center_offset_min[0] + nop->center[0], pdat.search_min[0]);
53  right[0] = Tp->nearest_image_intpos_to_intpos_X(nop->center_offset_max[0] + nop->center[0], pdat.search_min[0]);
54 
55  /* check whether we can stop walking along this branch */
56  if(left[0] > pdat.search_range[0] && right[0] > left[0])
57  return NODE_DISCARD;
58 
59  left[1] = Tp->nearest_image_intpos_to_intpos_Y(nop->center_offset_min[1] + nop->center[1], pdat.search_min[1]);
60  right[1] = Tp->nearest_image_intpos_to_intpos_Y(nop->center_offset_max[1] + nop->center[1], pdat.search_min[1]);
61 
62  /* check whether we can stop walking along this branch */
63  if(left[1] > pdat.search_range[1] && right[1] > left[1])
64  return NODE_DISCARD;
65 
66  left[2] = Tp->nearest_image_intpos_to_intpos_Z(nop->center_offset_min[2] + nop->center[2], pdat.search_min[2]);
67  right[2] = Tp->nearest_image_intpos_to_intpos_Z(nop->center_offset_max[2] + nop->center[2], pdat.search_min[2]);
68 
69  /* check whether we can stop walking along this branch */
70  if(left[2] > pdat.search_range[2] && right[2] > left[2])
71  return NODE_DISCARD;
72 
73  return NODE_OPEN;
74 }
75 
76 /* Check whether the potential neighbor referenced by p/p_type/shmrank is inside the smoothing length, and if yes
77  * add it to the interaction list built up for this particle.
78  */
79 inline void sph::sph_density_check_particle_particle_interaction(pinfo &pdat, int p, int p_type, unsigned char shmrank)
80 {
81 #ifdef PRESERVE_SHMEM_BINARY_INVARIANCE
82  if(skip_actual_force_computation)
83  return;
84 #endif
85 
86  if(p_type == NODE_TYPE_LOCAL_PARTICLE) /* local particle */
87  {
88  particle_data *P = get_Pp(p, shmrank);
89  sph_particle_data *SphP = get_SphPp(p, shmrank);
90 
91  if(P->getType() > 0)
92  return;
93 
94  if(P->get_Ti_Current() != All.Ti_Current)
95  Tp->drift_particle(P, SphP, All.Ti_Current); // this function avoids race conditions
96 
97  double posdiff[3];
98  Tp->nearest_image_intpos_to_pos(P->IntPos, pdat.searchcenter, posdiff); /* converts the integer distance to floating point */
99 
100  if(posdiff[0] * posdiff[0] + posdiff[1] * posdiff[1] + posdiff[2] * posdiff[2] > pdat.hsml2)
101  return;
102 
103  if(pdat.numngb >= MAX_NGBS)
104  Terminate("pdat.numngb >= MAX_NGBS");
105 
106  int n = pdat.numngb++;
107 
108  Ngbdensdat[n].IntPos = P->IntPos;
109  Ngbdensdat[n].VelPred = SphP->VelPred;
110  Ngbdensdat[n].Mass = P->getMass();
111 #ifdef PRESSURE_ENTROPY_SPH
112  Ngbdensdat[n].EntropyToInvGammaPred = SphP->EntropyToInvGammaPred;
113 #endif
114 #ifdef TIMEDEP_ART_VISC
115  Ngbdensdat[n].Csnd = SphP->Csnd;
116 #endif
117  }
118  else if(p_type == NODE_TYPE_FETCHED_PARTICLE)
119  {
120  foreign_sphpoint_data *foreignpoint = get_foreignpointsp(p - EndOfForeignNodes, shmrank);
121 
122  /* converts the integer distance to floating point */
123  double posdiff[3];
124  Tp->nearest_image_intpos_to_pos(foreignpoint->IntPos, pdat.searchcenter, posdiff);
125 
126  if(posdiff[0] * posdiff[0] + posdiff[1] * posdiff[1] + posdiff[2] * posdiff[2] > pdat.hsml2)
127  return;
128 
129  if(pdat.numngb >= MAX_NGBS)
130  Terminate("pdat.numngb >= MAX_NGBS");
131 
132  int n = pdat.numngb++;
133 
134  Ngbdensdat[n].IntPos = foreignpoint->IntPos;
135  Ngbdensdat[n].VelPred = foreignpoint->SphCore.VelPred;
136  Ngbdensdat[n].Mass = foreignpoint->Mass;
137 #ifdef PRESSURE_ENTROPY_SPH
138  Ngbdensdat[n].EntropyToInvGammaPred = foreignpoint->SphCore.EntropyToInvGammaPred;
139 #endif
140 #ifdef TIMEDEP_ART_VISC
141  Ngbdensdat[n].Csnd = foreignpoint->SphCore.Csnd;
142 #endif
143  }
144  else
145  Terminate("unexpected");
146 }
147 
148 /* Continues to walk the tree for the particle referenced in pdat by opening a node.
149  */
150 inline void sph::sph_density_open_node(pinfo &pdat, ngbnode *nop, int mintopleafnode, int committed)
151 {
152  /* open node */
153  int p = nop->nextnode;
154  unsigned char shmrank = nop->nextnode_shmrank;
155 
156  while(p != nop->sibling || (shmrank != nop->sibling_shmrank && nop->sibling >= MaxPart + D->NTopnodes))
157  {
158  if(p < 0)
159  Terminate(
160  "p=%d < 0 nop->sibling=%d nop->nextnode=%d shmrank=%d nop->sibling_shmrank=%d nop->foreigntask=%d "
161  "first_nontoplevelnode=%d",
162  p, nop->sibling, nop->nextnode, shmrank, nop->sibling_shmrank, nop->OriginTask, MaxPart + D->NTopnodes);
163 
164  int next;
165  unsigned char next_shmrank;
166  char type;
167 
168  if(p < MaxPart) /* a local particle */
169  {
170  /* note: here shmrank cannot change */
171  next = get_nextnodep(shmrank)[p];
172  next_shmrank = shmrank;
174  }
175  else if(p < MaxPart + MaxNodes) /* an internal node */
176  {
177  ngbnode *nop = get_nodep(p, shmrank);
178  next = nop->sibling;
179  next_shmrank = nop->sibling_shmrank;
180  type = NODE_TYPE_LOCAL_NODE;
181  }
182  else if(p >= ImportedNodeOffset && p < EndOfTreePoints) /* an imported Treepoint particle */
183  {
184  Terminate("not expected for SPH");
185  }
186  else if(p >= EndOfTreePoints && p < EndOfForeignNodes) /* an imported tree node */
187  {
188  ngbnode *nop = get_nodep(p, shmrank);
189  next = nop->sibling;
190  next_shmrank = nop->sibling_shmrank;
191  type = NODE_TYPE_FETCHED_NODE;
192  }
193  else if(p >= EndOfForeignNodes) /* an imported particle below an imported tree node */
194  {
195  foreign_sphpoint_data *foreignpoint = get_foreignpointsp(p - EndOfForeignNodes, shmrank);
196 
197  next = foreignpoint->Nextnode;
198  next_shmrank = foreignpoint->Nextnode_shmrank;
200  }
201  else
202  {
203  /* a pseudo point */
204  Terminate(
205  "should not happen: p=%d MaxPart=%d MaxNodes=%d ImportedNodeOffset=%d EndOfTreePoints=%d EndOfForeignNodes=%d "
206  "shmrank=%d nop->nextnode=%d nop->cannot_be_opened_locally=%d nop->not_empty=%d nop-TopNodes=%d",
208  (int)nop->cannot_be_opened_locally, (int)nop->not_empty, (int)(nop - TopNodes));
209  }
210 
211  sph_density_interact(pdat, p, type, shmrank, mintopleafnode, committed);
212 
213  p = next;
214  shmrank = next_shmrank;
215  }
216 }
217 
218 /* Take care of SPH density interaction between the particle referenced in pdat, and the node
219  * referenced through no/shmrank. The node can either be a normal node or an imported node from another
220  * shared memory machine, and the node can either be on the present MPI rank, or from another MPI rank in the
221  * local shared memory machine.
222  */
223 inline void sph::sph_density_interact(pinfo &pdat, int no, char no_type, unsigned char shmrank, int mintopleafnode, int committed)
224 {
225  if(no_type <= NODE_TYPE_FETCHED_PARTICLE) // we are interacting with a particle
226  {
227  sph_density_check_particle_particle_interaction(pdat, no, no_type, shmrank);
228  }
229  else // we are interacting with a node
230  {
231  ngbnode *nop = get_nodep(no, shmrank);
232 
233  if(nop->not_empty == 0)
234  return;
235 
236  if(no < MaxPart + MaxNodes) // we have a top-levelnode
237  if(nop->nextnode >= MaxPart + MaxNodes) // if the next node is not a top-level, we have a leaf node
238  mintopleafnode = no;
239 
240  int openflag = sph_density_evaluate_particle_node_opening_criterion(pdat, nop);
241 
242  if(openflag == NODE_OPEN) /* we need to open it */
243  {
244  if(nop->cannot_be_opened_locally.load(std::memory_order_acquire))
245  {
246  // are we in the same shared memory node?
248  {
249  Terminate("this should not happen any more");
250  }
251  else
252  {
253  tree_add_to_fetch_stack(nop, no, shmrank); // will only add unique copies
254 
255  tree_add_to_work_stack(pdat.target, no, shmrank, mintopleafnode);
256  }
257  }
258  else
259  {
260  int min_buffer_space =
262 
263  if(min_buffer_space >= committed + 8 * TREE_NUM_BEFORE_NODESPLIT)
264  sph_density_open_node(pdat, nop, mintopleafnode, committed + 8 * TREE_NUM_BEFORE_NODESPLIT);
265  else
266  tree_add_to_work_stack(pdat.target, no, shmrank, mintopleafnode);
267  }
268  }
269  }
270 }
271 
272 /* Internal driver routine to compute densities for the particles with indices listed in the targetlist
273  * array.
274  */
275 void sph::densities_determine(int ntarget, int *targetlist)
276 {
277  Ngbdensdat = (ngbdata_density *)Mem.mymalloc("Ngbdensdat", MAX_NGBS * sizeof(ngbdata_density));
278 
279  NumOnWorkStack = 0;
283 
284  WorkStack = (workstack_data *)Mem.mymalloc("WorkStack", AllocWorkStackBaseHigh * sizeof(workstack_data));
285 
286  for(int i = 0; i < ntarget; i++)
287  {
288  int target = targetlist[i];
289 
290  clear_density_result(&Tp->SphP[target]);
291 
292  WorkStack[NumOnWorkStack].Target = target;
295  WorkStack[NumOnWorkStack].MinTopLeafNode = MaxPart + D->NTopnodes;
296  NumOnWorkStack++;
297  }
298 
299 #ifdef PRESERVE_SHMEM_BINARY_INVARIANCE
300  workstack_data *WorkStackBak = (workstack_data *)Mem.mymalloc("WorkStackBak", NumOnWorkStack * sizeof(workstack_data));
301  int NumOnWorkStackBak = NumOnWorkStack;
302  memcpy(WorkStackBak, WorkStack, NumOnWorkStack * sizeof(workstack_data));
303 #endif
304 
305  // set a default size of the fetch stack equal to half the work stack (this may still be somewhat too large)
307  StackToFetch = (fetch_data *)Mem.mymalloc_movable(&StackToFetch, "StackToFetch", MaxOnFetchStack * sizeof(fetch_data));
308 
309 #ifdef PRESERVE_SHMEM_BINARY_INVARIANCE
310  for(int rep = 0; rep < 2; rep++)
311  {
312  if(rep == 0)
313  {
314  skip_actual_force_computation = true;
315  }
316  else
317  {
318  skip_actual_force_computation = false;
319  NumOnWorkStack = NumOnWorkStackBak;
320  memcpy(WorkStack, WorkStackBak, NumOnWorkStack * sizeof(workstack_data));
321  }
322 #endif
323 
324  while(NumOnWorkStack > 0) // repeat until we are out of work
325  {
326  NewOnWorkStack = 0; // gives the new entries
327  NumOnFetchStack = 0;
329 
330  TIMER_START(CPU_DENSWALK);
331 
332  int item = 0;
333 
334  while(item < NumOnWorkStack)
335  {
336  int committed = 8 * TREE_NUM_BEFORE_NODESPLIT;
337  int min_buffer_space =
339  if(min_buffer_space >= committed)
340  {
341  int target = WorkStack[item].Target;
342  int no = WorkStack[item].Node;
343  int shmrank = WorkStack[item].ShmRank;
344  int mintopleaf = WorkStack[item].MinTopLeafNode;
345  item++;
346 
347  pinfo pdat;
348  get_pinfo(target, pdat);
349 
350  if(no == MaxPart)
351  {
352  // we have a pristine particle that's processed for the first time
353  sph_density_interact(pdat, no, NODE_TYPE_LOCAL_NODE, shmrank, mintopleaf, committed);
354  }
355  else
356  {
357  // we have a node that we previously could not open
358  ngbnode *nop = get_nodep(no, shmrank);
359 
360  if(nop->cannot_be_opened_locally)
361  {
362  Terminate("item=%d: no=%d now we should be able to open it!", item, no);
363  }
364  else
365  sph_density_open_node(pdat, nop, mintopleaf, committed);
366  }
367 
368  density_evaluate_kernel(pdat);
369  }
370  else
371  break;
372  }
373 
374  if(item == 0 && NumOnWorkStack > 0)
375  Terminate("Can't even process a single particle");
376 
377  TIMER_STOP(CPU_DENSWALK);
378 
379  TIMER_START(CPU_DENSFETCH);
380 
382 
383  TIMER_STOP(CPU_DENSFETCH);
384 
385  /* now reorder the workstack such that we are first going to do residual pristine particles, and then
386  * imported nodes that hang below the first leaf nodes */
388  memmove(WorkStack, WorkStack + item, NumOnWorkStack * sizeof(workstack_data));
389 
390  /* now let's sort such that we can go deep on top-level node branches, allowing us to clear them out eventually */
392 
393  max_ncycles++;
394  }
395 #ifdef PRESERVE_SHMEM_BINARY_INVARIANCE
396  }
397 #endif
398 
399  Mem.myfree(StackToFetch);
400 #ifdef PRESERVE_SHMEM_BINARY_INVARIANCE
401  Mem.myfree(WorkStackBak);
402 #endif
403  Mem.myfree(WorkStack);
404  Mem.myfree(Ngbdensdat);
405 }
406 
407 /* This first makes sure that all active SPH particles are drifted to the current time,
408  * and then calls the SPH density computation for them.
409  */
411 {
413  {
414  /* now drift the active hydro particles if not done already */
415  for(int i = 0; i < Tp->TimeBinsHydro.NActiveParticles; i++)
416  {
419 
420  Tp->drift_particle(P, SphP, All.Ti_Current);
421 #ifdef TIMEDEP_ART_VISC
422  SphP->DivVelOld = SphP->DivVel;
423 #endif
424  }
425 
426  /* compute density (and updates pressure) */
428  }
429 }
430 
431 /* Compute the SPH densities for all particles listed in the index-array. The routine finds a suitable
432  * smoothing length by a bisection algorithm.
433  */
434 void sph::density(int *list, int ntarget)
435 {
436  TIMER_STORE;
437  TIMER_START(CPU_DENSITY);
438 
439  D->mpi_printf("SPH-DENSITY: Begin density calculation. (presently allocated=%g MB)\n", Mem.getAllocatedBytesInMB());
440  D->mpi_printf("SPH-DENSITY: Ndensities=%llu (task zero has: NumGas=%d, Ndensities=%d)\n", Tp->TimeBinsHydro.GlobalNActiveParticles,
441  Tp->NumGas, ntarget);
442 
443  double ta = Logs.second();
444 
445  /* Create list of targets. We do this here to simplify the treatment later on */
446  int *targetList = (int *)Mem.mymalloc("TargetList", Tp->NumGas * sizeof(int));
447  MyFloat *Left = (MyFloat *)Mem.mymalloc("Left", Tp->NumGas * sizeof(MyFloat));
448  MyFloat *Right = (MyFloat *)Mem.mymalloc("Right", Tp->NumGas * sizeof(MyFloat));
449 
450  int ndensities = ntarget;
451 
452  for(int i = 0; i < ntarget; i++)
453  {
454  int target = list[i];
455  targetList[i] = target;
456  Left[target] = Right[target] = 0.0;
457  }
458 
459  int iter = 0;
460 
461  // let's grab at most half the still available memory for imported points and nodes
462  int nspace = (0.5 * Mem.FreeBytes) / (sizeof(ngbnode) + 8 * sizeof(foreign_sphpoint_data));
463 
464  MaxForeignNodes = nspace;
465  MaxForeignPoints = 8 * nspace;
466  NumForeignNodes = 0;
467  NumForeignPoints = 0;
468 
471 
472  /* the following two arrays hold imported tree nodes and imported points to augment the local tree */
473  Foreign_Nodes = (ngbnode *)Mem.mymalloc_movable(&Foreign_Nodes, "Foreign_Nodes", MaxForeignNodes * sizeof(ngbnode));
474  Foreign_Points = (foreign_sphpoint_data *)Mem.mymalloc_movable(&Foreign_Points, "Foreign_Points",
476 
478 
479  max_ncycles = 0;
480 
482 
483  do
484  {
485  double t0 = Logs.second();
486 
487  /* now do the primary work with this call */
488 
489  densities_determine(ndensities, targetList);
490 
491  /* do final operations on results */
492  int npleft = 0;
493 
494  for(int i = 0; i < ndensities; i++)
495  {
496  int target = targetList[i];
497 
498  if(target >= 0)
499  {
500  if(Tp->P[target].getType() != 0)
501  Terminate("P[target].getType() != 0");
502 
503  sph_particle_data *SphP = Tp->SphP;
504  if(SphP[target].Density > 0)
505  {
506 #ifdef WENDLAND_BIAS_CORRECTION
507  SphP[target].Density -= get_density_bias(SphP[target].Hsml, Tp->P[target].getMass(), All.DesNumNgb);
508 #endif
509  SphP[target].DhsmlDensityFactor *= SphP[target].Hsml / (NUMDIMS * SphP[target].Density);
510  if(SphP[target].DhsmlDensityFactor >
511  -0.9) /* note: this would be -1 if only a single particle at zero lag is found */
512  SphP[target].DhsmlDensityFactor = 1 / (1 + SphP[target].DhsmlDensityFactor);
513  else
514  SphP[target].DhsmlDensityFactor = 1;
515 
516 #ifndef IMPROVED_VELOCITY_GRADIENTS
517  SphP[target].CurlVel = sqrt(SphP[target].Rot[0] * SphP[target].Rot[0] + SphP[target].Rot[1] * SphP[target].Rot[1] +
518  SphP[target].Rot[2] * SphP[target].Rot[2]) /
519  SphP[target].Density;
520 
521  SphP[target].DivVel /= SphP[target].Density;
522 #else
523  SphP[target].set_velocity_gradients();
524 #endif
525  SphP[target].DtHsml = (1.0 / NUMDIMS) * SphP[target].DivVel * SphP[target].Hsml;
526  SphP[target].DtDensity = -SphP[target].DivVel * SphP[target].Density;
527 
528 #ifndef PRESSURE_ENTROPY_SPH
529  SphP[target].set_thermodynamic_variables();
530 #endif
531  }
532 
533 #ifdef PRESSURE_ENTROPY_SPH
534  if(SphP[target].EntropyToInvGammaPred > 0 && SphP[target].PressureSphDensity > 0)
535  {
536  SphP[target].DhsmlDerivedDensityFactor *=
537  SphP[target].Hsml / (NUMDIMS * SphP[target].Density * SphP[target].EntropyToInvGammaPred);
538  SphP[target].DhsmlDerivedDensityFactor *= -SphP[target].DhsmlDensityFactor;
539  SphP[target].PressureSphDensity /= SphP[target].EntropyToInvGammaPred;
540 #ifdef WENDLAND_BIAS_CORRECTION /* Dehnen & Aly 2012, eq (18), (19) */
541  SphP[target].PressureSphDensity -= get_density_bias(SphP[target].Hsml, Tp->P[target].getMass(), All.DesNumNgb);
542 #endif
543  SphP[target].DtPressureSphDensity = -SphP[target].DivVel * SphP[target].PressureSphDensity;
544  SphP[target].set_thermodynamic_variables();
545  }
546  else
547  {
548  SphP[target].DhsmlDerivedDensityFactor = 0;
549  SphP[target].EntropyToInvGammaPred = 0;
550  SphP[target].PressureSphDensity = 0;
551  }
552 
553 #endif
554 #ifdef TIMEDEP_ART_VISC
555  double dt = (Tp->P[target].getTimeBinHydro() ? (((integertime)1) << Tp->P[target].getTimeBinHydro()) : 0) *
557  double dtime = All.cf_atime * dt / All.cf_atime_hubble_a;
558  SphP[target].set_viscosity_coefficient(dtime);
559 #endif
560 #ifdef ADAPTIVE_HYDRO_SOFTENING
561  Tp->P[target].setSofteningClass(Tp->get_softeningtype_for_hydro_particle(target));
562 #endif
563  /* now check whether we had enough neighbours */
564  double desnumngb = All.DesNumNgb;
565  double desnumngbdev = All.MaxNumNgbDeviation;
566 
567  double hfac = 1;
568  for(int i = 0; i < NUMDIMS; i++)
569  {
570  hfac *= SphP[target].Hsml;
571  }
572 
573  SphP[target].NumNgb = NORM_COEFF * hfac * SphP[target].Density / Tp->P[target].getMass();
574  if(SphP[target].NumNgb < (desnumngb - desnumngbdev) || (SphP[target].NumNgb > (desnumngb + desnumngbdev)))
575  {
576  if(Left[target] > 0 && Right[target] > 0)
577  if((Right[target] - Left[target]) < 1.0e-3 * Left[target])
578  {
579  /* this one should be ok */
580  continue;
581  }
582 
583  /* need to redo this particle */
584  targetList[npleft++] = target;
585 
586  if(SphP[target].NumNgb < (desnumngb - desnumngbdev))
587  Left[target] = std::max<double>(SphP[target].Hsml, Left[target]);
588  else
589  {
590  if(Right[target] != 0)
591  {
592  if(SphP[target].Hsml < Right[target])
593  Right[target] = SphP[target].Hsml;
594  }
595  else
596  Right[target] = SphP[target].Hsml;
597  }
598 
599  if(iter >= MAXITER - 10)
600  {
601  double pos[3];
602  Tp->intpos_to_pos(Tp->P[target].IntPos, pos); /* converts the integer coordinates to floating point */
603 
604  printf("target=%d Hsml=%g task=%d ID=%llu Left=%g Right=%g Ngbs=%g Right-Left=%g\n pos=(%g|%g|%g)\n", target,
605  SphP[target].Hsml, D->ThisTask, (unsigned long long)Tp->P[target].ID.get(), Left[target], Right[target],
606  SphP[target].NumNgb, Right[target] - Left[target], pos[0], pos[1], pos[2]);
607  myflush(stdout);
608  }
609 
610  if(Right[target] > 0 && Left[target] > 0)
611  SphP[target].Hsml = pow(0.5 * (pow(Left[target], 3) + pow(Right[target], 3)), 1.0 / 3);
612  else
613  {
614  if(Right[target] == 0 && Left[target] == 0)
615  Terminate("Right[i] == 0 && Left[i] == 0 SphP[i].Hsml=%g\n", SphP[target].Hsml);
616 
617  if(Right[target] == 0 && Left[target] > 0)
618  {
619  if(Tp->P[target].getType() == 0 && fabs(SphP[target].NumNgb - desnumngb) < 0.5 * desnumngb)
620  {
621  double fac = 1 - (SphP[target].NumNgb - desnumngb) / (NUMDIMS * SphP[target].NumNgb) *
622  SphP[target].DhsmlDensityFactor;
623 
624  if(fac < 1.26)
625  SphP[target].Hsml *= fac;
626  else
627  SphP[target].Hsml *= 1.26;
628  }
629  else
630  SphP[target].Hsml *= 1.26;
631  }
632 
633  if(Right[target] > 0 && Left[target] == 0)
634  {
635  if(Tp->P[target].getType() == 0 && fabs(SphP[target].NumNgb - desnumngb) < 0.5 * desnumngb && iter < 4)
636  {
637  double fac = 1 - (SphP[target].NumNgb - desnumngb) / (NUMDIMS * SphP[target].NumNgb) *
638  SphP[target].DhsmlDensityFactor;
639 
640  if(fac > 1 / 1.26)
641  SphP[target].Hsml *= fac;
642  else
643  SphP[target].Hsml /= 1.26;
644  }
645  else
646  SphP[target].Hsml /= 1.26;
647  }
648  }
649  }
650  }
651  }
652 
653  ndensities = npleft;
654 
655  double t1 = Logs.second();
656 
657  if(npleft > 0)
658  {
659  iter++;
660 
661  D->mpi_printf("SPH-DENSITY: ngb iteration %4d: took %8.3f , need to repeat for %012lld local particles.\n", iter,
662  Logs.timediff(t0, t1), npleft);
663 
664  if(iter > MAXITER)
665  Terminate("failed to converge in neighbour iteration in density()\n");
666  }
667  else
668  D->mpi_printf("SPH-DENSITY: ngb iteration %4d: took %8.3f\n", ++iter, Logs.timediff(t0, t1));
669  }
670  while(ndensities > 0);
671 
672  TIMER_START(CPU_DENSIMBALANCE);
673 
674  MPI_Allreduce(MPI_IN_PLACE, &max_ncycles, 1, MPI_INT, MPI_MAX, D->Communicator);
675 
676  TIMER_STOP(CPU_DENSIMBALANCE);
677 
679 
680  /* free temporary buffers */
681 
682  Mem.myfree(Foreign_Points);
683  Mem.myfree(Foreign_Nodes);
684 
685  Mem.myfree(Right);
686  Mem.myfree(Left);
687  Mem.myfree(targetList);
688 
689  double tb = Logs.second();
690 
691  TIMER_STOPSTART(CPU_DENSITY, CPU_LOGS);
692 
693  D->mpi_printf("SPH-DENSITY: density computation done. took %8.3f\n", Logs.timediff(ta, tb));
694 
695  struct detailed_timings
696  {
697  double tree, wait, fetch, all;
698  double numnodes;
700  double fillfacFgnNodes, fillfacFgnPoints;
701  };
702  detailed_timings timer, tisum, timax;
703 
704  timer.tree = TIMER_DIFF(CPU_DENSWALK);
705  timer.wait = TIMER_DIFF(CPU_DENSIMBALANCE);
706  timer.fetch = TIMER_DIFF(CPU_DENSFETCH);
707  timer.all = timer.tree + timer.wait + timer.fetch + TIMER_DIFF(CPU_DENSITY);
708  timer.numnodes = NumNodes;
709  timer.NumForeignNodes = NumForeignNodes;
710  timer.NumForeignPoints = NumForeignPoints;
711  timer.fillfacFgnNodes = NumForeignNodes / ((double)MaxForeignNodes);
712  timer.fillfacFgnPoints = NumForeignPoints / ((double)MaxForeignPoints);
713 
714  MPI_Reduce((double *)&timer, (double *)&tisum, (int)(sizeof(detailed_timings) / sizeof(double)), MPI_DOUBLE, MPI_SUM, 0,
715  D->Communicator);
716  MPI_Reduce((double *)&timer, (double *)&timax, (int)(sizeof(detailed_timings) / sizeof(double)), MPI_DOUBLE, MPI_MAX, 0,
717  D->Communicator);
718 
720 
721  if(D->ThisTask == 0)
722  {
723  fprintf(Logs.FdDensity, "Nf=%9lld highest active timebin=%d total-Nf=%lld\n", Tp->TimeBinsHydro.GlobalNActiveParticles,
725  fprintf(Logs.FdDensity, " work-load balance: %g part/sec: raw=%g, effective=%g\n",
726  timax.tree / ((tisum.tree + 1e-20) / D->NTask), Tp->TimeBinsGravity.GlobalNActiveParticles / (tisum.tree + 1.0e-20),
727  Tp->TimeBinsGravity.GlobalNActiveParticles / ((timax.tree + 1.0e-20) * D->NTask));
728  fprintf(Logs.FdDensity,
729  " maximum number of nodes: %g, filled: %g NumForeignNodes: max=%g avg=%g fill=%g NumForeignPoints: max=%g avg=%g "
730  "fill=%g cycles=%d\n",
731  timax.numnodes, timax.numnodes / MaxNodes, timax.NumForeignNodes, tisum.NumForeignNodes / D->NTask,
732  timax.fillfacFgnNodes, timax.NumForeignPoints, tisum.NumForeignPoints / D->NTask, timax.fillfacFgnPoints, max_ncycles);
733  fprintf(Logs.FdDensity, " avg times: <all>=%g <tree>=%g <wait>=%g <fetch>=%g sec\n", tisum.all / D->NTask,
734  tisum.tree / D->NTask, tisum.wait / D->NTask, tisum.fetch / D->NTask);
736  }
737 
738  TIMER_STOP(CPU_LOGS);
739 }
740 
741 #ifdef EXPLICIT_VECTORIZATION
742 
743 /* Main SPH compute kernel, in a version that is explicitely vectorized. Several neighbours are
744  * processed through one vector. The calculation should be semantically equivalent to the standard
745  * looped version without explicit vector instructions.
746  */
747 void sph::density_evaluate_kernel(pinfo &pdat)
748 {
749  particle_data *targetP = &Tp->P[pdat.target];
750  sph_particle_data *targetSphP = &Tp->SphP[pdat.target];
751 
752  double shinv, shinv3, shinv4;
753  kernel_hinv(targetSphP->Hsml, &shinv, &shinv3, &shinv4);
754 
755  Vec4d hinv(shinv);
756  Vec4d hinv3(shinv3);
757  Vec4d hinv4(shinv4);
758 
759  Vec4d dwnorm(NORM * shinv3);
760  Vec4d dwknorm(NORM * shinv4);
761 
762  Vec4d v_i[NUMDIMS];
763  for(int i = 0; i < NUMDIMS; i++)
764  {
765  v_i[i] = targetSphP->VelPred[i];
766  }
767 
768  Vec4d cs_i(targetSphP->Csnd);
769  const int vector_length = 4;
770  const int array_length = (pdat.numngb + vector_length - 1) & (-vector_length);
771 
772  for(int n = pdat.numngb; n < array_length; n++) /* fill up neighbour array so that sensible data is accessed */
773  Ngbdensdat[n] = Ngbdensdat[0];
774 
775  Vec4d decayVel_loc(0);
776 
777  for(int n = 0; n < array_length; n += vector_length)
778  {
779  struct ngbdata_density *ngb0 = &Ngbdensdat[n + 0];
780  struct ngbdata_density *ngb1 = &Ngbdensdat[n + 1];
781  struct ngbdata_density *ngb2 = &Ngbdensdat[n + 2];
782  struct ngbdata_density *ngb3 = &Ngbdensdat[n + 3];
783 
784  Vec4d dpos[NUMDIMS];
785 #if defined(LONG_X_BITS) || defined(LONG_Y_BITS) || defined(LONG_Z_BITS)
786  double posdiff[array_length][3];
787  for(int i = 0; i < 4; i++)
788  {
789  Tp->nearest_image_intpos_to_pos(targetP->IntPos, Ngbdensdat[n + i].IntPos, &(posdiff[i][0]));
790  }
791 
792  for(int i = 0; i < NUMDIMS; i++)
793  {
794  dpos[i] = Vec4d(posdiff[0][i], posdiff[1][i], posdiff[2][i], posdiff[3][i]);
795  }
796 #else
797  for(int i = 0; i < NUMDIMS; i++)
798  {
799  dpos[i] = Tp->nearest_image_intpos_to_doublepos_vectorial(targetP->IntPos[i], ngb0->IntPos[i], ngb1->IntPos[i],
800  ngb2->IntPos[i], ngb3->IntPos[i]);
801  }
802 #endif
803 
804  Vec4d v_j[NUMDIMS];
805  for(int i = 0; i < NUMDIMS; i++)
806  {
807  v_j[i] = Vec4d(ngb0->VelPred[i], ngb1->VelPred[i], ngb2->VelPred[i], ngb3->VelPred[i]);
808  }
809 
810  Vec4d mass_j(ngb0->Mass, ngb1->Mass, ngb2->Mass, ngb3->Mass);
811  Vec4d r2(0);
812 
813  for(int i = 0; i < NUMDIMS; i++)
814  {
815  r2 += dpos[i] * dpos[i];
816  }
817 
818  Vec4d r = sqrt(r2);
819 
820  Vec4d u = r * hinv;
821 
822  /* now calculate the kernel */
823  Vec4d wk, dwk;
824  kernel_main_vector(u, dwnorm, dwknorm, &wk, &dwk);
825 
826  if(n + vector_length > pdat.numngb) /* we have excess elements */
827  {
828  mass_j.cutoff(vector_length - (array_length - pdat.numngb));
829  wk.cutoff(vector_length - (array_length - pdat.numngb));
830  }
831 
832  Vec4d mj_wk = mass_j * wk;
833 
834  targetSphP->Density += horizontal_add(mj_wk);
835 
836 #ifdef PRESSURE_ENTROPY_SPH
837  Vec4d entr_j(ngb0->EntropyToInvGammaPred, ngb1->EntropyToInvGammaPred, ngb2->EntropyToInvGammaPred, ngb3->EntropyToInvGammaPred);
838 
839  targetSphP->PressureSphDensity += horizontal_add(mj_wk * entr_j);
840 
841  targetSphP->DhsmlDerivedDensityFactor += horizontal_add(-mass_j * entr_j * (NUMDIMS * hinv * wk + u * dwk));
842 #endif
843 
844  targetSphP->DhsmlDensityFactor += horizontal_add(-mass_j * (NUMDIMS * hinv * wk + u * dwk));
845 
846  Vec4db decision_r_gt_0 = (r > 0);
847 
848  r = select(decision_r_gt_0, r, 1.0); /* note, for r=0, we have dwk=0 */
849 
850  Vec4d mj_dwk_r = mass_j * dwk / r;
851 
852  Vec4d dv[NUMDIMS];
853  for(int i = 0; i < NUMDIMS; i++)
854  {
855  dv[i] = v_i[i] - v_j[i];
856  }
857 
858 #ifndef IMPROVED_VELOCITY_GRADIENTS
859  Vec4d dpos_times_dv(0);
860  for(int i = 0; i < NUMDIMS; i++)
861  dpos_times_dv += dpos[i] * dv[i];
862 
863  targetSphP->DivVel += horizontal_add(-mj_dwk_r * dpos_times_dv);
864 
865 #ifdef TWODIMS
866  targetSphP->Rot[2] += horizontal_add(mj_dwk_r * (dpos[1] * dv[0] - dpos[0] * dv[1]));
867 #endif
868 #ifdef THREEDIMS
869  targetSphP->Rot[0] += horizontal_add(mj_dwk_r * (dpos[2] * dv[1] - dpos[1] * dv[2]));
870  targetSphP->Rot[1] += horizontal_add(mj_dwk_r * (dpos[0] * dv[2] - dpos[2] * dv[0]));
871  targetSphP->Rot[2] += horizontal_add(mj_dwk_r * (dpos[1] * dv[0] - dpos[0] * dv[1]));
872 #endif
873 #else
874  for(int i = 0; i < NUMDIMS; i++)
875  {
876  for(int j = 0; j < NUMDIMS; j++)
877  {
878  targetSphP->dvel[i][j] -= horizontal_add(mj_dwk_r * dv[i] * dpos[j]);
879  }
880  }
881  targetSphP->dpos.dx_dx -= horizontal_add(mj_dwk_r * dpos[0] * dpos[0]);
882  targetSphP->dpos.dx_dy -= horizontal_add(mj_dwk_r * dpos[0] * dpos[1]);
883  targetSphP->dpos.dx_dz -= horizontal_add(mj_dwk_r * dpos[0] * dpos[2]);
884  targetSphP->dpos.dy_dy -= horizontal_add(mj_dwk_r * dpos[1] * dpos[1]);
885  targetSphP->dpos.dy_dz -= horizontal_add(mj_dwk_r * dpos[1] * dpos[2]);
886  targetSphP->dpos.dz_dz -= horizontal_add(mj_dwk_r * dpos[2] * dpos[2]);
887 #endif
888 #ifdef TIMEDEP_ART_VISC
889  Vec4d vdotr2 = 0;
890  for(int i = 0; i < NUMDIMS; i++)
891  {
892  vdotr2 += dpos[i] * dv[i];
893  }
894 
896  vdotr2 += All.cf_atime2_hubble_a * r2;
897 
898  Vec4d mu_ij = vdotr2 / (All.cf_afac3 * All.cf_atime * r);
899 
900  Vec4d cs_j(ngb0->Csnd, ngb1->Csnd, ngb2->Csnd, ngb3->Csnd);
901  Vec4d cs_sum = cs_i + cs_j;
902 
903  Vec4d decay_vel_2 = cs_sum - mu_ij;
904 
905  Vec4db decision = (vdotr2 > 0);
906 
907  Vec4d decay_vel = select(decision, cs_sum, decay_vel_2);
908 
909  Vec4d fac_decay_vel = select(decision_r_gt_0, 1, 0);
910 
911  decay_vel = decay_vel * fac_decay_vel;
912 
913  decayVel_loc = max(decayVel_loc, decay_vel);
914 
915 #endif
916  }
917 
918 #ifdef TIMEDEP_ART_VISC
919  for(int i = 0; i < vector_length; i++)
920  {
921  if(decayVel_loc[i] > targetSphP->decayVel)
922  targetSphP->decayVel = decayVel_loc[i];
923  }
924 #endif
925 }
926 
927 #else
928 
929 /* Main SPH compute kernel. This function carries out the SPH computations for the neighbouring particle
930  * data stored in the Ngbdensdat[] array. The results are added to the particle referenced through the pdat
931  * structure.
932  */
933 void sph::density_evaluate_kernel(pinfo &pdat)
934 {
935  particle_data *targetP = &Tp->P[pdat.target];
936  sph_particle_data *targetSphP = &Tp->SphP[pdat.target];
937 
938  kernel_density kernel;
939 
940  kernel_hinv(targetSphP->Hsml, &kernel.hinv, &kernel.hinv3, &kernel.hinv4);
941 
942  for(int n = 0; n < pdat.numngb; n++)
943  {
944  struct ngbdata_density *ngb = &Ngbdensdat[n];
945 
946  /* note: in periodic case, closest image will be found through integer wrap around */
947 
948  double posdiff[3];
949  Tp->nearest_image_intpos_to_pos(targetP->IntPos, ngb->IntPos, posdiff); /* converts the integer distance to floating point */
950 
951  for(int i = 0; i < NUMDIMS; i++)
952  {
953  kernel.dpos[i] = posdiff[i];
954  }
955 
956  double r2 = 0;
957 
958  for(int i = 0; i < NUMDIMS; i++)
959  {
960  r2 += kernel.dpos[i] * kernel.dpos[i];
961  }
962 
963  kernel.r = sqrt(r2);
964 
965  double u = kernel.r * kernel.hinv;
966 
967  kernel_main(u, kernel.hinv3, kernel.hinv4, &kernel.wk, &kernel.dwk, COMPUTE_WK_AND_DWK);
968 
969  double mass_j = ngb->Mass;
970  kernel.mj_wk = (mass_j * kernel.wk);
971 
972  targetSphP->Density += kernel.mj_wk;
973 
974 #ifdef PRESSURE_ENTROPY_SPH
975  targetSphP->PressureSphDensity += kernel.mj_wk * ngb->EntropyToInvGammaPred;
976  targetSphP->DhsmlDerivedDensityFactor +=
977  (-mass_j * ngb->EntropyToInvGammaPred * (NUMDIMS * kernel.hinv * kernel.wk + u * kernel.dwk));
978 
979 #endif
980 
981  targetSphP->DhsmlDensityFactor += (-mass_j * (NUMDIMS * kernel.hinv * kernel.wk + u * kernel.dwk));
982 
983  if(kernel.r > 0)
984  {
985  kernel.mj_dwk_r = mass_j * kernel.dwk / kernel.r;
986 
987  for(int i = 0; i < NUMDIMS; i++)
988  {
989  kernel.dv[i] = targetSphP->VelPred[i] - ngb->VelPred[i];
990  }
991 
992 #ifndef IMPROVED_VELOCITY_GRADIENTS
993  double dpos_times_dv = 0;
994  for(int i = 0; i < NUMDIMS; i++)
995  dpos_times_dv += kernel.dpos[i] * kernel.dv[i];
996 
997  targetSphP->DivVel += (-kernel.mj_dwk_r * (dpos_times_dv));
998 #ifdef TWODIMS
999  targetSphP->Rot[2] += (kernel.mj_dwk_r * (kernel.dpos[1] * kernel.dv[0] - kernel.dpos[0] * kernel.dv[1]));
1000 #endif
1001 #ifdef THREEDIMS
1002  targetSphP->Rot[0] += (kernel.mj_dwk_r * (kernel.dpos[2] * kernel.dv[1] - kernel.dpos[1] * kernel.dv[2]));
1003  targetSphP->Rot[1] += (kernel.mj_dwk_r * (kernel.dpos[0] * kernel.dv[2] - kernel.dpos[2] * kernel.dv[0]));
1004  targetSphP->Rot[2] += (kernel.mj_dwk_r * (kernel.dpos[1] * kernel.dv[0] - kernel.dpos[0] * kernel.dv[1]));
1005 #endif
1006 #else
1007  for(int i = 0; i < NUMDIMS; i++)
1008  {
1009  for(int j = 0; j < NUMDIMS; j++)
1010  {
1011  targetSphP->dvel[i][j] -= kernel.mj_dwk_r * kernel.dv[i] * kernel.dpos[j];
1012  }
1013  }
1014 
1015  targetSphP->dpos.dx_dx -= kernel.mj_dwk_r * kernel.dpos[0] * kernel.dpos[0];
1016  targetSphP->dpos.dx_dy -= kernel.mj_dwk_r * kernel.dpos[0] * kernel.dpos[1];
1017  targetSphP->dpos.dx_dz -= kernel.mj_dwk_r * kernel.dpos[0] * kernel.dpos[2];
1018  targetSphP->dpos.dy_dy -= kernel.mj_dwk_r * kernel.dpos[1] * kernel.dpos[1];
1019  targetSphP->dpos.dy_dz -= kernel.mj_dwk_r * kernel.dpos[1] * kernel.dpos[2];
1020  targetSphP->dpos.dz_dz -= kernel.mj_dwk_r * kernel.dpos[2] * kernel.dpos[2];
1021 
1022 #endif
1023 #ifdef TIMEDEP_ART_VISC
1024  double vdotr2 = 0;
1025  for(int i = 0; i < NUMDIMS; i++)
1026  {
1027  vdotr2 += kernel.dpos[i] * kernel.dv[i];
1028  }
1029 
1031  vdotr2 += All.cf_atime2_hubble_a * r2;
1032 
1033  double mu_ij = vdotr2 / (All.cf_afac3 * All.cf_atime * kernel.r);
1034  double decay_vel;
1035 
1036  if(vdotr2 < 0)
1037  decay_vel = targetSphP->Csnd + ngb->Csnd - mu_ij;
1038  else
1039  decay_vel = targetSphP->Csnd + ngb->Csnd;
1040 
1041  if(decay_vel > targetSphP->decayVel)
1042  targetSphP->decayVel = decay_vel;
1043 #endif
1044  }
1045  }
1046 }
1047 #endif
1048 
1049 /* this routine clears the fields in the SphP particle structure that are additively computed by the SPH density loop
1050  * by summing over neighbours
1051  */
1052 inline void sph::clear_density_result(sph_particle_data *SphP)
1053 {
1054  SphP->Density = 0;
1055  SphP->DhsmlDensityFactor = 0;
1056  SphP->DivVel = 0;
1057 
1058  for(int k = 0; k < 3; k++)
1059  SphP->Rot[k] = 0;
1060 
1061 #ifdef PRESSURE_ENTROPY_SPH
1062  SphP->PressureSphDensity = 0;
1063  SphP->DhsmlDerivedDensityFactor = 0;
1064 #endif
1065 #ifdef IMPROVED_VELOCITY_GRADIENTS
1066  SphP->dpos = {0};
1067  for(int i = 0; i < NUMDIMS; i++)
1068  {
1069  for(int j = 0; j < NUMDIMS; j++)
1070  {
1071  SphP->dvel[i][j] = 0;
1072  }
1073  }
1074 #endif
1075 #ifdef TIMEDEP_ART_VISC
1076  SphP->decayVel = 0;
1077 #endif
1078 }
global_data_all_processes All
Definition: main.cc:40
MyIDType get(void) const
Definition: idstorage.h:87
int NTopnodes
Definition: domain.h:52
void nearest_image_intpos_to_pos(const MyIntPosType *const a, const MyIntPosType *const b, T *posdiff)
MyIntPosType nearest_image_intpos_to_intpos_Y(const MyIntPosType a, const MyIntPosType b)
MyIntPosType nearest_image_intpos_to_intpos_Z(const MyIntPosType a, const MyIntPosType b)
void intpos_to_pos(MyIntPosType *intpos, T *posdiff)
MyIntPosType nearest_image_intpos_to_intpos_X(const MyIntPosType a, const MyIntPosType b)
double timediff(double t0, double t1)
Definition: logs.cc:488
double second(void)
Definition: logs.cc:471
FILE * FdDensity
Definition: logs.h:48
double getAllocatedBytesInMB(void)
Definition: mymalloc.h:68
size_t FreeBytes
Definition: mymalloc.h:48
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
int * GetNodeIDForSimulCommRank
int Island_ThisTask
int drift_particle(particle_data *P, sph_particle_data *SphP, integertime time1, bool ignore_light_cone=false)
This function drifts a particle i to time1.
Definition: predict.cc:313
TimeBinData TimeBinsGravity
Definition: simparticles.h:205
particle_data * P
Definition: simparticles.h:54
sph_particle_data * SphP
Definition: simparticles.h:59
TimeBinData TimeBinsHydro
Definition: simparticles.h:204
void compute_densities(void)
Definition: density.cc:410
void density(int *targetlist, int ntarget)
Definition: density.cc:434
static bool compare_workstack(const workstack_data &a, const workstack_data &b)
Definition: tree.h:197
sph_particle_data * get_SphPp(int n, unsigned char shmrank)
Definition: tree.h:292
void tree_add_to_fetch_stack(ngbnode *nop, int nodetoopen, unsigned char shmrank)
Definition: tree.h:207
void tree_add_to_work_stack(int target, int no, unsigned char shmrank, int mintopleafnode)
Definition: tree.h:232
#define NORM_COEFF
Definition: constants.h:376
#define NUMDIMS
Definition: constants.h:369
#define MAXITER
Definition: constants.h:305
int integertime
Definition: constants.h:331
double mycxxsort(T *begin, T *end, Tcomp comp)
Definition: cxxsort.h:39
float MyFloat
Definition: dtypes.h:86
uint32_t MyIntPosType
Definition: dtypes.h:35
#define LEVEL_ALWAYS_OPEN
Definition: dtypes.h:383
#define TREE_MIN_WORKSTACK_SIZE
Definition: gravtree.h:26
#define NODE_DISCARD
Definition: gravtree.h:37
#define NODE_TYPE_FETCHED_NODE
Definition: gravtree.h:33
#define TREE_EXPECTED_CYCLES
Definition: gravtree.h:27
#define NODE_TYPE_FETCHED_PARTICLE
Definition: gravtree.h:31
#define NODE_OPEN
Definition: gravtree.h:36
#define NODE_TYPE_LOCAL_NODE
Definition: gravtree.h:32
#define NODE_TYPE_LOCAL_PARTICLE
Definition: gravtree.h:29
#define COMPUTE_WK_AND_DWK
Definition: kernel.h:109
#define NORM
Definition: kernel.h:47
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 TIMER_STORE
Copies the current value of CPU times to a stored variable, such that differences with respect to thi...
Definition: logs.h:257
#define TIMER_DIFF(counter)
Returns amount elapsed for the timer since last save with TIMER_STORE.
Definition: logs.h:250
#define TIMER_STOPSTART(stop, start)
Stops the timer 'stop' and starts the timer 'start'.
Definition: logs.h:231
#define Terminate(...)
Definition: macros.h:19
shmem Shmem
Definition: main.cc:45
memory Mem
Definition: main.cc:44
half fabs(half arg)
Definition: half.hpp:2627
expr pow(half base, half exp)
Definition: half.hpp:2803
expr sqrt(half arg)
Definition: half.hpp:2777
#define MAX_NGBS
Definition: sph.h:18
int NActiveParticles
Definition: timestep.h:18
int * ActiveParticleList
Definition: timestep.h:20
long long GlobalNActiveParticles
Definition: timestep.h:19
unsigned char level
Definition: tree.h:64
unsigned char nextnode_shmrank
Definition: tree.h:66
std::atomic< unsigned char > cannot_be_opened_locally
Definition: tree.h:70
unsigned char not_empty
Definition: tree.h:73
int nextnode
Definition: tree.h:57
vector< MyIntPosType > center
Definition: tree.h:51
int OriginTask
Definition: tree.h:61
unsigned char sibling_shmrank
Definition: tree.h:65
int sibling
Definition: tree.h:53
unsigned char Nextnode_shmrank
Definition: ngbtree.h:111
sph_particle_data_hydrocore SphCore
Definition: ngbtree.h:105
MyIntPosType IntPos[3]
Definition: ngbtree.h:108
integertime Ti_Current
Definition: allvars.h:188
double mj_wk
Definition: kernel.h:22
double dpos[NUMDIMS]
Definition: kernel.h:17
double wk
Definition: kernel.h:20
double hinv
Definition: kernel.h:21
double dv[NUMDIMS]
Definition: kernel.h:19
double r
Definition: kernel.h:18
double mj_dwk_r
Definition: kernel.h:22
double hinv4
Definition: kernel.h:21
double hinv3
Definition: kernel.h:21
double dwk
Definition: kernel.h:20
std::atomic< integertime > Ti_Current
Definition: ngbtree.h:56
MySignedIntPosType center_offset_max[3]
Definition: ngbtree.h:48
MySignedIntPosType center_offset_min[3]
Definition: ngbtree.h:47
void drift_node(integertime time1, simparticles *Sp)
Definition: ngbtree.h:58
void setSofteningClass(unsigned char softclass)
integertime get_Ti_Current(void)
MyDouble getMass(void)
unsigned char getTimeBinHydro(void)
MyIDStorage ID
Definition: particle_data.h:70
unsigned char getType(void)
MyIntPosType IntPos[3]
Definition: particle_data.h:53
void set_thermodynamic_variables(void)
void myflush(FILE *fstream)
Definition: system.cc:125
#define TREE_NUM_BEFORE_NODESPLIT
Definition: tree.h:16