GADGET-4
begrun.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 "compiler-command-line-args.h"
13#include "gadgetconfig.h"
14
15#include <hdf5.h>
16#include <math.h>
17#include <mpi.h>
18#include <stdio.h>
19#include <stdlib.h>
20#include <string.h>
21#include <sys/stat.h>
22#include <sys/types.h>
23#include <unistd.h>
24
25#include "../cooling_sfr/cooling.h"
26#include "../data/allvars.h"
27#include "../data/dtypes.h"
28#include "../data/mymalloc.h"
29#include "../data/simparticles.h"
30#include "../fmm/fmm.h"
31#include "../fof/fof.h"
32#include "../gitversion/version.h"
33#include "../gravity/ewald.h"
34#include "../gravtree/gravtree.h"
35#include "../io/hdf5_util.h"
36#include "../io/io.h"
37#include "../io/parameters.h"
38#include "../lightcone/lightcone.h"
39#include "../logs/logs.h"
40#include "../logs/timer.h"
41#include "../main/main.h"
42#include "../main/simulation.h"
43#include "../mergertree/mergertree.h"
44#include "../mpi_utils/mpi_utils.h"
45#include "../mpi_utils/shared_mem_handler.h"
46#include "../ngbtree/ngbtree.h"
47#include "../pm/pm.h"
48#include "../system/system.h"
49#include "../time_integration/driftfac.h"
50#include "../time_integration/timestep.h"
51
58void sim::hello(void)
59{
61 "\n ___ __ ____ ___ ____ ____ __\n / __) /__\\ ( _ \\ / __)( ___)(_ _)___ /. |\n"
62 "( (_-. /(__)\\ )(_) )( (_-. )__) )( (___)(_ _)\n \\___/(__)(__)(____/ \\___/(____) (__) (_)\n\n");
63
64 mpi_printf("This is Gadget, version %s.\nGit commit %s, %s\n\n", GADGET_VERSION, GIT_COMMIT, GIT_DATE);
65
66 mpi_printf("Code was compiled with the following compiler and flags:\n%s\n\n\n", compiler_flags);
67
68 mpi_printf("Code was compiled with the following settings:\n");
69
70 if(ThisTask == 0)
71 {
73 }
74
75 mpi_printf("\n\nRunning on %d MPI tasks.\n\n", NTask);
76
77#ifdef EXPLICIT_VECTORIZATION
78
79 int instrset = instrset_detect();
80 mpi_printf("CPU supports instruction sets: ");
81 if(instrset >= 1)
82 mpi_printf("SSE ");
83 if(instrset >= 2)
84 mpi_printf("SSE2 ");
85 if(instrset >= 3)
86 mpi_printf("SSE3 ");
87 if(instrset >= 4)
88 mpi_printf("SSSE3 ");
89 if(instrset >= 5)
90 mpi_printf("SSE4.1 ");
91 if(instrset >= 6)
92 mpi_printf("SSE4.2 ");
93 if(instrset >= 7)
94 mpi_printf("AVX ");
95 if(instrset >= 8)
96 mpi_printf("AVX2 ");
97 if(instrset >= 9)
98 mpi_printf("AVX512F ");
99 if(instrset >= 10)
100 mpi_printf("AVX512VL AVX512BW AVX512DQ ");
101 mpi_printf("\n\n");
102
103 if(instrset < 7)
104 {
106 "You compiled with explicit vectorization in terms of AVX instructions, but this CPU does not support AVX or higher.\n\n");
107 endrun();
108 }
109#endif
110
111 mpi_printf("\n");
112 mpi_printf("BEGRUN: Size of particle structure %4d [bytes]\n", (int)sizeof(particle_data));
113 mpi_printf("BEGRUN: Size of sph particle structure %4d [bytes]\n", (int)sizeof(sph_particle_data));
114 mpi_printf("BEGRUN: Size of gravity tree node %4d [bytes]\n", (int)sizeof(gravnode));
115 mpi_printf("BEGRUN: Size of neighbour tree node %4d [bytes]\n", (int)sizeof(ngbnode));
116 mpi_printf("BEGRUN: Size of subfind auxiliary data %4d [bytes]\n", (int)sizeof(subfind_data));
117
118 mpi_printf("\n");
119}
120
130void sim::begrun1(const char *parameterFile)
131{
133
134 int errorFlag = All.read_parameter_file(parameterFile); /* ... read in parameters for this run on task 0*/
135
136 if(ThisTask == 0)
137 {
138 int n = strlen(All.OutputDir);
139 if(n > 0)
140 if(All.OutputDir[n - 1] != '/')
141 strcat(All.OutputDir, "/");
142
143 mkdir(All.OutputDir, 02755);
144 }
145
146 /* now communicate the relevant parameters to the other processes, *including* the shared memory handler */
147 /* this also tells the shared memory handler how much memory it may allocate */
148 MPI_Bcast(All.get_data_ptr(), All.get_data_size(), MPI_BYTE, 0, MPI_COMM_WORLD);
149
150#ifdef HOST_MEMORY_REPORTING
152#endif
153
155
156 MPI_Bcast(&errorFlag, 1, MPI_INT, 0, Communicator);
157
158 if(errorFlag)
159 {
161 {
162 char c = 0;
163 // need to send this flag to our shared memory rank so that it also ends itself
164 MPI_Send(&c, 1, MPI_BYTE, Shmem.MyShmRankInGlobal, TAG_KEY, MPI_COMM_WORLD);
165 }
166
167 mpi_printf("\nWe stop because of an error in the parameterfile.\n\n");
168 MPI_Finalize();
169 exit(0);
170 }
171
172 if(All.OutputListOn)
174 else
176
177 All.write_used_parameters(All.OutputDir, "parameters-usedvalues");
178
180
181#ifdef ENABLE_HEALTHTEST
182 healthtest();
183#endif
184
185 set_units();
186
189
191
192#ifdef LIGHTCONE_PARTICLES
193 LightCone.lightcone_init_geometry(All.LightConeDefinitionFile);
194#endif
195
196 /* disable automatic error printing of HDF5 - we check for errors ourselves */
197 H5Eset_auto(H5E_DEFAULT, NULL, NULL);
198
199#ifdef DEBUG
201#endif
202
203#ifdef PMGRID
204 GravTree.short_range_init();
205#endif
206
209
210#ifdef COOLING
213 CoolSfr.InitCool();
214#endif
215
216#ifdef STARFORMATION
217 CoolSfr.init_clouds();
218#endif
219
220#if((!defined(PMGRID) || (defined(PMGRID) && defined(TREEPM_NOTIMESPLIT))) && defined(SELFGRAVITY) && defined(PERIODIC)) || \
221 defined(FORCETEST)
223#endif
224
226
227#ifdef DEBUG_SYMTENSORS
229#endif
230
231#ifdef PMGRID
234 {
235#ifdef PERIODIC
236 PM.pm_init_periodic(&Sp);
237#ifdef PLACEHIGHRESREGION
238 PM.pm_init_nonperiodic(&Sp);
239#endif
240#else
241 PM.pm_init_nonperiodic(&Sp);
242#endif
243 }
244#endif
245
247
249
250#ifdef DEBUG_MD5
251 Logs.log_debug_md5("START");
252#endif
253}
254
262void sim::begrun2(void)
263{
264 char contfname[MAXLEN_PATH_EXTRA];
265 sprintf(contfname, "%scont", All.OutputDir);
266 unlink(contfname);
267
270
272 {
274
276 {
277 // We actually have multiple shared memory nodes in which we set aside one MPI rank for shared memory communication.
278 // Tell the shared memory node to initialize its drift table as well, so that it can be used for particle preditions on the
279 // fly
280 if(Shmem.Island_ThisTask == 0)
281 {
282 // update All on shared memory handler, to be sure that be get the correct times
283 MPI_Send(All.get_data_ptr(), All.get_data_size(), MPI_BYTE, Shmem.MyShmRankInGlobal, TAG_DRIFT_INIT, MPI_COMM_WORLD);
284 }
285 }
286 }
287
288#ifdef LIGHTCONE
289#ifdef LIGHTCONE_PARTICLES
290 if(LightCone.lightcone_init_times())
291 endrun();
292#endif
293#ifdef LIGHTCONE_MASSMAPS
294 LightCone.lightcone_init_massmaps();
295 if(LightCone.lightcone_massmap_report_boundaries())
296 endrun();
297#endif
298 double linklength = 0;
299
300#ifdef FOF
301 fof<simparticles> FoF{Communicator, &Sp, &Domain};
302 linklength = FoF.fof_get_comoving_linking_length();
303#endif
304
305 LightCone.lightcone_init_intposconverter(linklength);
306#endif
307
309 All.Ti_nextoutput = find_next_outputtime(All.Ti_Current + 100);
310 else if(All.RestartFlag == RST_RESUME)
311 All.Ti_nextoutput = find_next_outputtime(All.Ti_Current + 1);
312 else
313 All.Ti_nextoutput = find_next_outputtime(All.Ti_Current);
314
316
317#ifdef REDUCE_FLUSH
318 All.FlushLast = Logs.CPUThisRun;
319#endif
320
321 // update All on shared memory handler, just to allow it to access its elements if needed
323 MPI_Send(All.get_data_ptr(), All.get_data_size(), MPI_BYTE, Shmem.MyShmRankInGlobal, TAG_ALL_UPDATE, MPI_COMM_WORLD);
324
325#if defined(FORCETEST) && defined(FORCETEST_TESTFORCELAW)
326 gravity_forcetest_testforcelaw();
327#endif
328}
329
336{
340
343 else
345
350
352 {
353 /* check whether the supplied value of All.Hubble makes sense */
354 if(All.HubbleParam != 1.0)
355 {
356 double Hubble_expected = HUBBLE * All.UnitTime_in_s;
357 if(fabs(Hubble_expected - All.Hubble) > 1.0e-3 * All.Hubble)
358 Terminate(
359 "You have supplied All.Hubble=%g/All.HubbleParam=%g. For this choice, we would expect All.Hubble=%g. We better stop\n",
360 All.Hubble, All.HubbleParam, Hubble_expected);
361 }
362 else
363 {
364 double Hubble_expected = 0.7 * HUBBLE * All.UnitTime_in_s;
365 if(All.Hubble < 0.5 * Hubble_expected || All.Hubble > 2.0 * Hubble_expected)
366 Terminate(
367 "You have supplied All.Hubble=%g/All.HubbleParam=%g. For this choice, we would expect All.Hubble somewhere in the "
368 "ballpark of %g. We better stop\n",
369 All.Hubble, All.HubbleParam, Hubble_expected);
370 }
371 }
372
373 mpi_printf("BEGRUN: Hubble (internal units) = %g\n", All.Hubble);
374 mpi_printf("BEGRUN: h = %g\n", All.HubbleParam);
375 mpi_printf("BEGRUN: G (internal units) = %g\n", All.G);
376 mpi_printf("BEGRUN: UnitMass_in_g = %g\n", All.UnitMass_in_g);
377 mpi_printf("BEGRUN: UnitLenth_in_cm = %g\n", All.UnitLength_in_cm);
378 mpi_printf("BEGRUN: UnitTime_in_s = %g\n", All.UnitTime_in_s);
379 mpi_printf("BEGRUN: UnitVelocity_in_cm_per_s = %g\n", All.UnitVelocity_in_cm_per_s);
380 mpi_printf("BEGRUN: UnitDensity_in_cgs = %g\n", All.UnitDensity_in_cgs);
381 mpi_printf("BEGRUN: UnitEnergy_in_cgs = %g\n", All.UnitEnergy_in_cgs);
382 mpi_printf("\n");
383
384#ifdef STARFORMATION
385 CoolSfr.set_units_sfr();
386#endif
387}
388
395void sim::endrun(void)
396{
397 mpi_printf("endrun called, calling MPI_Finalize()\nbye!\n\n");
398 fflush(stdout);
399
401 {
402 char c = 0;
403 // need to send this flag to our shared memory rank so that it also ends itself
404 MPI_Send(&c, 1, MPI_BYTE, Shmem.MyShmRankInGlobal, TAG_KEY, MPI_COMM_WORLD);
405 }
406
407 /* The hdf5 library will sometimes register an atexit() handler that calls its error handler.
408 * This is set to my_hdf_error_handler, which calls MPI_Abort.
409 * Calling MPI_Abort after MPI_Finalize is not allowed.
410 * Hence unset the HDF error handler here*/
411 H5Eset_auto(H5E_DEFAULT, NULL, NULL);
412
413 MPI_Finalize();
414 exit(0);
415}
global_data_all_processes All
Definition: main.cc:40
void init_drift_table(void)
Definition: driftfac.cc:26
void ewald_init(void)
This function initializes tables with the correction force and the correction potential due to the pe...
Definition: ewald.cc:69
double CPUThisRun
Definition: logs.h:43
void open_logfiles(void)
Open files for logging.
Definition: logs.cc:35
void log_debug_md5(const char *msg)
void check_maxmemsize_setting(int maxmemsize)
Definition: mymalloc.cc:684
void mymalloc_init(int maxmemsize, enum restart_options restartflag)
Initialize memory manager.
Definition: mymalloc.cc:58
void write_used_parameters(const char *dirname, const char *fname)
Definition: parameters.cc:195
int read_parameter_file(const char *fname)
This function parses the parameter file.
Definition: parameters.cc:57
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 Island_ThisTask
int Island_NTask
int MyShmRankInGlobal
domain< simparticles > Domain
Definition: simulation.h:58
void begrun1(const char *parameterFile)
This function performs the initial set-up of the simulation.
Definition: begrun.cc:130
gwalk GravTree
Definition: simulation.h:65
void hello(void)
Definition: begrun.cc:58
void set_units(void)
Computes conversion factors between internal code units and the cgs-system.
Definition: begrun.cc:335
void begrun2(void)
This function does late setup, after the IC file has been loaded but before run() is called.
Definition: begrun.cc:262
void healthtest(void)
Definition: healthtest.cc:36
void endrun(void)
This function aborts the simulations.
Definition: begrun.cc:395
simparticles Sp
Definition: simulation.h:56
TimeBinData TimeBinsGravity
Definition: simparticles.h:205
TimeBinData TimeBinsHydro
Definition: simparticles.h:204
#define GADGET_VERSION
Definition: constants.h:15
#define HUBBLE
Definition: constants.h:125
#define GRAVITY
Definition: constants.h:118
#define MAXLEN_PATH_EXTRA
Definition: constants.h:301
#define SEC_PER_YEAR
Definition: constants.h:128
#define SEC_PER_MEGAYEAR
Definition: constants.h:127
@ RST_STARTFROMSNAP
Definition: dtypes.h:315
@ RST_POWERSPEC
Definition: dtypes.h:317
@ RST_BEGIN
Definition: dtypes.h:313
@ RST_RESUME
Definition: dtypes.h:314
ewald Ewald
Definition: main.cc:42
void my_create_HDF5_special_integer_types(void)
Definition: hdf5_util.cc:135
void my_create_HDF5_halfprec_handler(void)
Definition: hdf5_util.cc:121
logs Logs
Definition: main.cc:43
#define Terminate(...)
Definition: macros.h:15
shmem Shmem
Definition: main.cc:45
driftfac Driftfac
Definition: main.cc:41
void output_compile_time_options(void)
int instrset_detect(void)
void my_mpi_types_init(void)
Definition: mpi_types.cc:72
#define TAG_ALL_UPDATE
Definition: mpi_utils.h:100
#define TAG_KEY
Definition: mpi_utils.h:29
#define TAG_DRIFT_INIT
Definition: mpi_utils.h:99
memory Mem
Definition: main.cc:44
half fabs(half arg)
Definition: half.hpp:2627
expr pow(half base, half exp)
Definition: half.hpp:2803
void timebins_init(const char *name, int *MaxPart)
Definition: timestep.cc:439
double UnitVelocity_in_cm_per_s
Definition: allvars.h:119
enum restart_options RestartFlag
Definition: allvars.h:68
void read_outputlist(char *fname)
This function reads a table with a list of desired output times.
Definition: allvars.cc:211
void register_parameters(void)
Definition: allvars.cc:55
integertime Ti_Current
Definition: allvars.h:188
char OutputListFilename[MAXLEN_PATH]
Definition: allvars.h:272
char OutputDir[MAXLEN_PATH]
Definition: allvars.h:272
integertime Ti_nextoutput
Definition: allvars.h:189
void some_parameter_checks(void)
Definition: allvars.cc:255
double GravityConstantInternal
Definition: allvars.h:128
void set_cosmo_factors_for_current_time(void)
Definition: allvars.cc:20
size_t get_data_size(void)
Definition: allvars.h:362
char * get_data_ptr(void)
Definition: allvars.h:360
void symtensor_test(void)
void init_rng(int thistask)
Definition: system.cc:42
void enable_core_dumps_and_fpu_exceptions(void)
const char * GIT_COMMIT
const char * GIT_DATE