GADGET-4
init_entropy.cc
Go to the documentation of this file.
1/*******************************************************************************
2 * \copyright This file is part of the GADGET4 N-body/SPH code developed
3 * \copyright by Volker Springel. Copyright (C) 2014-2020 by Volker Springel
4 * \copyright (vspringel@mpa-garching.mpg.de) and all contributing authors.
5 *******************************************************************************/
6
12#include "gadgetconfig.h"
13
14#include <math.h>
15#include <mpi.h>
16#include <stdio.h>
17#include <stdlib.h>
18#include <string.h>
19
20#include "../data/allvars.h"
21#include "../data/dtypes.h"
22#include "../data/mymalloc.h"
23#include "../domain/domain.h"
24#include "../io/io.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/peano.h"
30#include "../sph/kernel.h"
31#include "../system/system.h"
32#include "../time_integration/timestep.h"
33
34#ifdef PRESSURE_ENTROPY_SPH
35
36#define MAX_ITER_ENTROPY 100
37#define ENTROPY_TOLERANCE 1.0e-5
38
43void sph::init_entropy(void)
44{
46 TIMER_START(CPU_DENSITY);
47
48 D->mpi_printf("SPH-INIT-ENTROPY: Begin entropy calculation. (presently allocated=%g MB)\n", Mem.getAllocatedBytesInMB());
49
50 /* Create list of targets. We do this here to simplify the treatment later on */
51 int *targetList = (int *)Mem.mymalloc("targetlist", Tp->NumGas * sizeof(int));
52
53 int ndensities = 0;
54
55 for(int i = 0; i < Tp->TimeBinsHydro.NActiveParticles; i++)
56 {
57 int target = Tp->TimeBinsHydro.ActiveParticleList[i];
58 if(target < 0 || target >= Tp->NumGas)
59 Terminate("target=%d i=%d\n", target, i);
60 targetList[ndensities++] = target;
61 }
62
63 int iter = 0;
64
65 // let's grab at most half the still available memory for imported points and nodes
66 int nspace = (0.33 * Mem.FreeBytes) / (sizeof(ngbnode) + 8 * sizeof(foreign_sphpoint_data));
67
68 MaxForeignNodes = nspace;
69 MaxForeignPoints = 8 * nspace;
72
75
76 /* the following two arrays hold imported tree nodes and imported points to augment the local tree */
77 Foreign_Nodes = (ngbnode *)Mem.mymalloc_movable(&Foreign_Nodes, "Foreign_Nodes", MaxForeignNodes * sizeof(ngbnode));
78 Foreign_Points = (foreign_sphpoint_data *)Mem.mymalloc_movable(&Foreign_Points, "Foreign_Points",
80
81
82
83 double tstart = Logs.second();
84
85 int global_left_particles = 0;
86
87 MPI_Allreduce(&ndensities, &global_left_particles, 1, MPI_INT, MPI_SUM, D->Communicator);
88
89 do
90 {
91 double t0 = Logs.second();
92
93 /* Since EntropyToInvGammaPred of remote particles can change, we have to import the particles in every iteration */
94
95 MaxForeignNodes = nspace;
96 MaxForeignPoints = 8 * nspace;
99
102
104
105 max_ncycles = 0;
106
108
109 /* now do the primary work with this call */
110 densities_determine(ndensities, targetList);
111
112 MPI_Allreduce(MPI_IN_PLACE, &max_ncycles, 1, MPI_INT, MPI_MAX, D->Communicator);
113
115
116 /* do final operations on results */
117
118 double entropy_old;
119
120 int npleft = 0;
121
122 for(int i = 0; i < ndensities; i++)
123 {
124 int target = targetList[i];
125 if(target >= 0)
126 {
127 if(Tp->P[target].getType() != 0)
128 Terminate("P[target].getType() != 0");
129
130 sph_particle_data *SphP = Tp->SphP;
131
132 if(SphP[target].EntropyToInvGammaPred > 0 && SphP[target].Density > 0)
133 {
134 entropy_old = SphP[target].Entropy;
135 SphP[target].PressureSphDensity /= SphP[target].EntropyToInvGammaPred;
136 SphP[target].Entropy =
137 GAMMA_MINUS1 * SphP[target].EntropyPred / pow(SphP[target].PressureSphDensity * All.cf_a3inv, GAMMA_MINUS1);
138 SphP[target].EntropyToInvGammaPred = pow(SphP[target].Entropy, 1.0 / GAMMA);
139 }
140 else
141 {
142 entropy_old = SphP[target].Entropy;
143 SphP[target].PressureSphDensity = 0;
144 SphP[target].Entropy = 0;
145 SphP[target].EntropyToInvGammaPred = 0;
146 }
147 /* entropy has not converged yet */
148 if(fabs(entropy_old - SphP[target].Entropy) > ENTROPY_TOLERANCE * entropy_old)
149 targetList[npleft++] = target;
150 }
151 }
152
153 ndensities = npleft;
154
155 MPI_Allreduce(&ndensities, &global_left_particles, 1, MPI_INT, MPI_SUM, D->Communicator);
156
157 double t1 = Logs.second();
158
159 if(npleft > 0)
160 {
161 iter++;
162
163 D->mpi_printf("SPH-INIT-ENTROPY: ngb iteration %4d: took %8.3f , need to repeat for %012lld local particles.\n", iter,
164 Logs.timediff(t0, t1), npleft);
165
166 if(iter > MAXITER)
167 Terminate("failed to converge in neighbour iteration in density()\n");
168 }
169 else
170 D->mpi_printf("SPH-INIT-ENTROPY: ngb iteration %4d: took %8.3f\n", ++iter, Logs.timediff(t0, t1));
171 }
172 while(global_left_particles > 0);
173
174
175 TIMER_STOP(CPU_DENSITY);
176
177
178 /* free temporary buffers */
179
180 Mem.myfree(Foreign_Points);
181 Mem.myfree(Foreign_Nodes);
182
183 Mem.myfree(targetList);
184
185 double tb = Logs.second();
186
187 D->mpi_printf("SPH-INIT-ENTROPY: entropy calculation is done. took: %8.3f\n", Logs.timediff(tstart, tb));
188}
189
194void sph::setup_entropy_to_invgamma(void)
195{
197
198 /* Initialize entropy and entropy^invgamma with a fist guess coming from standard SPH density estimate. */
199 /* EntropyPred is untouched since it contains the internal energies needed for the iterative process; it */
200 /* will be set in init.c to the correct value. */
201 for(int i = 0; i < Tp->NumGas; i++)
202 {
204 Tp->SphP[i].EntropyToInvGammaPred = pow(Tp->SphP[i].Entropy, 1.0 / GAMMA);
205 }
206
207 init_entropy();
208}
209#endif
global_data_all_processes All
Definition: main.cc:40
double timediff(double t0, double t1)
Definition: logs.cc:488
double second(void)
Definition: logs.cc:471
double getAllocatedBytesInMB(void)
Definition: mymalloc.h:71
size_t FreeBytes
Definition: mymalloc.h:48
particle_data * P
Definition: simparticles.h:54
sph_particle_data * SphP
Definition: simparticles.h:59
TimeBinData TimeBinsHydro
Definition: simparticles.h:204
#define GAMMA_MINUS1
Definition: constants.h:110
#define GAMMA
Definition: constants.h:99
#define MAXITER
Definition: constants.h:305
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 Terminate(...)
Definition: macros.h:15
memory Mem
Definition: main.cc:44
half fabs(half arg)
Definition: half.hpp:2627
expr pow(half base, half exp)
Definition: half.hpp:2803
int NActiveParticles
Definition: timestep.h:18
int * ActiveParticleList
Definition: timestep.h:20
void set_cosmo_factors_for_current_time(void)
Definition: allvars.cc:20
unsigned char getType(void)