GADGET-4
grid.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 CREATE_GRID
15 
16 #include <gsl/gsl_rng.h>
17 #include <math.h>
18 #include <mpi.h>
19 #include <stdlib.h>
20 
21 #include "../data/allvars.h"
22 #include "../data/dtypes.h"
23 #include "../data/intposconvert.h"
24 #include "../data/mymalloc.h"
25 #include "../logs/timer.h"
26 #include "../main/simulation.h"
27 #include "../mpi_utils/mpi_utils.h"
28 #include "../ngenic/ngenic.h"
29 #include "../system/system.h"
30 
31 void ngenic::create_grid(void)
32 {
33  long long gridSize = All.GridSize;
34  long long partTotal = gridSize * gridSize * gridSize;
35  long long partPerTask = partTotal / NTask;
36 
37  Sp->RegionLen = All.BoxSize;
38  Sp->FacCoordToInt = pow(2.0, BITS_FOR_POSITIONS) / Sp->RegionLen;
39  Sp->FacIntToCoord = Sp->RegionLen / pow(2.0, BITS_FOR_POSITIONS);
40 
42 
43  double masstot = All.Omega0 * 3 * All.Hubble * All.Hubble / (8 * M_PI * All.G) * All.BoxSize * All.BoxSize * All.BoxSize;
44 
45  double m = masstot / (partTotal);
46 
47  for(int i = 0; i < NTYPES; i++)
48  All.MassTable[i] = 0.;
49 
50  All.MassTable[1] = m;
51 
52  Sp->NumGas = 0;
53  Sp->NumPart = partPerTask;
54 
55  if(ThisTask == NTask - 1)
56  {
57  Sp->NumPart = partTotal - Sp->NumPart * (NTask - 1);
58  }
59 
60  int max_load, max_sphload;
61  MPI_Allreduce(&Sp->NumPart, &max_load, 1, MPI_INT, MPI_MAX, Communicator);
62  MPI_Allreduce(&Sp->NumGas, &max_sphload, 1, MPI_INT, MPI_MAX, Communicator);
63 
64 #ifdef GENERATE_GAS_IN_ICS
65  Sp->TotNumGas = partTotal;
66  Sp->TotNumPart = 2 * partTotal;
67  max_sphload = max_load;
68  max_load *= 2;
69 #else
70  Sp->TotNumPart = partTotal;
71  Sp->TotNumGas = 0;
72 #endif
73 
74  Sp->MaxPart = max_load / (1.0 - 2 * ALLOC_TOLERANCE);
75  Sp->MaxPartSph = max_sphload / (1.0 - 2 * ALLOC_TOLERANCE);
76 
77  Sp->allocate_memory();
78 
79  for(int i = 0; i < Sp->NumPart; i++)
80  {
81  long long ipcell = ThisTask * partPerTask + i;
82  int x = ipcell / (All.GridSize * All.GridSize);
83  int xr = ipcell % (All.GridSize * All.GridSize);
84  int y = xr / All.GridSize;
85  int z = xr % All.GridSize;
86 
87  double xyz[3];
88  xyz[0] = x * All.BoxSize / All.GridSize;
89  xyz[1] = y * All.BoxSize / All.GridSize;
90  xyz[2] = z * All.BoxSize / All.GridSize;
91 
92  Sp->pos_to_intpos(xyz, Sp->P[i].IntPos);
93 
94  Sp->P[i].Vel[0] = 0.;
95  Sp->P[i].Vel[1] = 0.;
96  Sp->P[i].Vel[2] = 0.;
97 
98  Sp->P[i].ID.set(ipcell + 1);
99 
100  Sp->P[i].setType(1);
101  }
102 
103  mpi_printf("NGENIC: generated grid of size %d\n", All.GridSize);
104 }
105 
106 #endif
global_data_all_processes All
Definition: main.cc:40
#define ALLOC_TOLERANCE
Definition: constants.h:74
#define NTYPES
Definition: constants.h:308
#define M_PI
Definition: constants.h:56
#define BITS_FOR_POSITIONS
Definition: dtypes.h:37
expr pow(half base, half exp)
Definition: half.hpp:2803
double MassTable[NTYPES]
Definition: allvars.h:269