MP-Gadget  5.0.1.dev1-76bc7d4726-dirty
Functions
save.c File Reference
#include <math.h>
#include <stdlib.h>
#include <stdio.h>
#include <stddef.h>
#include <string.h>
#include <mpi.h>
#include "allvars.h"
#include "proto.h"
#include <bigfile-mpi.h>
#include <libgadget/types.h>
#include <libgadget/physconst.h>
#include <libgadget/utils.h>
#include <libgadget/cosmology.h>
#include <libgadget/walltime.h>
#include <libgadget/utils/mymalloc.h>
Include dependency graph for save.c:

Go to the source code of this file.

Functions

void _bigfile_utils_create_block_from_c_array (BigFile *bf, void *baseptr, const char *name, const char *dtype, size_t dims[], ptrdiff_t elsize, int NumFiles, int NumWriters, MPI_Comm comm)
 
static void saveblock (BigFile *bf, void *baseptr, int ptype, const char *bname, const char *dtype, int items_per_particle, const int NumPart, ptrdiff_t elsize, int NumFiles, int NumWriters)
 
void write_particle_data (IDGenerator *idgen, const int Type, BigFile *bf, const uint64_t FirstID, const int SavePrePos, int NumFiles, int NumWriters, struct ic_part_data *curICP)
 
void compute_mass (double *mass, int64_t TotNumPartCDM, int64_t TotNumPartGas, int64_t TotNuPart, double nufrac, const double BoxSize, Cosmology *CP, const struct genic_config GenicConfig)
 
void saveheader (BigFile *bf, int64_t TotNumPartCDM, int64_t TotNumPartGas, int64_t TotNuPart, double nufrac, const double BoxSize, Cosmology *CP, const struct genic_config GenicConfig)
 

Function Documentation

◆ _bigfile_utils_create_block_from_c_array()

void _bigfile_utils_create_block_from_c_array ( BigFile *  bf,
void *  baseptr,
const char *  name,
const char *  dtype,
size_t  dims[],
ptrdiff_t  elsize,
int  NumFiles,
int  NumWriters,
MPI_Comm  comm 
)

Definition at line 18 of file save.c.

19 {
20  BigBlock block;
21  BigArray array;
22  BigBlockPtr ptr;
23  ptrdiff_t strides[2];
24  int64_t TotNumPart;
25  MPI_Allreduce(&dims[0], &TotNumPart, 1, MPI_INT64, MPI_SUM, comm);
26 
27  strides[1] = dtype_itemsize(dtype);
28  strides[0] = elsize;
29 
30  big_array_init(&array, baseptr, dtype, 2, dims, strides);
31 
32  if(0 != big_file_mpi_create_block(bf, &block, name, dtype, dims[1], NumFiles, TotNumPart, comm)) {
33  endrun(0, "%s:%s\n", big_file_get_error_message(), name);
34  }
35 
36  if(0 != big_block_seek(&block, &ptr, 0)) {
37  endrun(0, "Failed to seek:%s\n", big_file_get_error_message());
38  }
39 
40  if(0 != big_block_mpi_write(&block, &ptr, &array, NumWriters, comm)) {
41  endrun(0, "Failed to write :%s\n", big_file_get_error_message());
42  }
43 
44  if(0 != big_block_mpi_close(&block, comm)) {
45  endrun(0, "%s:%s\n", big_file_get_error_message(), name);
46  }
47 }
const char * name
Definition: densitykernel.c:93
void endrun(int where, const char *fmt,...)
Definition: endrun.c:147
#define MPI_INT64
Definition: system.h:12
int TotNumPart
Definition: test_exchange.c:24

References endrun(), MPI_INT64, name, and TotNumPart.

Referenced by save_transfer(), and saveblock().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ compute_mass()

void compute_mass ( double *  mass,
int64_t  TotNumPartCDM,
int64_t  TotNumPartGas,
int64_t  TotNuPart,
double  nufrac,
const double  BoxSize,
Cosmology CP,
const struct genic_config  GenicConfig 
)

Definition at line 90 of file save.c.

91 {
92  double UnitTime_in_s = GenicConfig.units.UnitLength_in_cm / GenicConfig.units.UnitVelocity_in_cm_per_s;
93  double Grav = GRAVITY / pow(GenicConfig.units.UnitLength_in_cm, 3) * GenicConfig.units.UnitMass_in_g * pow(UnitTime_in_s, 2);
94  const double OmegatoMass = 3 * CP->Hubble * CP->Hubble / (8 * M_PI * Grav) * pow(BoxSize, 3);
95  double OmegaCDM = CP->Omega0;
96  mass[0] = mass[2] = mass[3] = mass[4] = mass[5] = 0;
97  if (TotNumPartGas > 0) {
98  mass[0] = (CP->OmegaBaryon) * 3 * CP->Hubble * CP->Hubble / (8 * M_PI * Grav) * pow(BoxSize, 3) / TotNumPartGas;
99  OmegaCDM -= CP->OmegaBaryon;
100  }
101  if(CP->MNu[0] + CP->MNu[1] + CP->MNu[2] > 0) {
102  double OmegaNu = get_omega_nu(&CP->ONu, 1);
103  OmegaCDM -= OmegaNu;
104  if(TotNuPart > 0) {
105  mass[2] = nufrac * OmegaNu * OmegatoMass / TotNuPart;
106  }
107  }
108  mass[1] = OmegaCDM * OmegatoMass / TotNumPartCDM;
109 }
double get_omega_nu(const _omega_nu *const omnu, const double a)
#define GRAVITY
Definition: physconst.h:5
static Cosmology * CP
Definition: power.c:27
double Omega0
Definition: cosmology.h:10
double OmegaBaryon
Definition: cosmology.h:19
double MNu[3]
Definition: cosmology.h:25
_omega_nu ONu
Definition: cosmology.h:24
double Hubble
Definition: cosmology.h:21
double UnitMass_in_g
Definition: unitsystem.h:8
double UnitLength_in_cm
Definition: unitsystem.h:10
double UnitVelocity_in_cm_per_s
Definition: unitsystem.h:9
struct UnitSystem units
Definition: allvars.h:35

References CP, get_omega_nu(), GRAVITY, Cosmology::Hubble, Cosmology::MNu, Cosmology::Omega0, Cosmology::OmegaBaryon, Cosmology::ONu, UnitSystem::UnitLength_in_cm, UnitSystem::UnitMass_in_g, genic_config::units, and UnitSystem::UnitVelocity_in_cm_per_s.

Referenced by main(), and saveheader().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ saveblock()

static void saveblock ( BigFile *  bf,
void *  baseptr,
int  ptype,
const char *  bname,
const char *  dtype,
int  items_per_particle,
const int  NumPart,
ptrdiff_t  elsize,
int  NumFiles,
int  NumWriters 
)
static

Definition at line 49 of file save.c.

49  {
50  size_t dims[2];
51  char name[128];
52  snprintf(name, 128, "%d/%s", ptype, bname);
53 
54  dims[0] = NumPart;
55  dims[1] = items_per_particle;
56  _bigfile_utils_create_block_from_c_array(bf, baseptr, name, dtype, dims, elsize, NumFiles, NumWriters, MPI_COMM_WORLD);
57 }
void _bigfile_utils_create_block_from_c_array(BigFile *bf, void *baseptr, const char *name, const char *dtype, size_t dims[], ptrdiff_t elsize, int NumFiles, int NumWriters, MPI_Comm comm)
Definition: save.c:18
static enum TransferType ptype
Definition: zeldovich.c:146

References _bigfile_utils_create_block_from_c_array(), name, and ptype.

Referenced by write_particle_data().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ saveheader()

void saveheader ( BigFile *  bf,
int64_t  TotNumPartCDM,
int64_t  TotNumPartGas,
int64_t  TotNuPart,
double  nufrac,
const double  BoxSize,
Cosmology CP,
const struct genic_config  GenicConfig 
)

Definition at line 112 of file save.c.

112  {
113  BigBlock bheader;
114  if(0 != big_file_mpi_create_block(bf, &bheader, "Header", NULL, 0, 0, 0, MPI_COMM_WORLD)) {
115  endrun(0, "failed to create block %s:%s\n", "Header",
116  big_file_get_error_message());
117  }
118 
119  int64_t totnumpart[6] = {TotNumPartGas, TotNumPartCDM, TotNuPart, 0, 0, 0};
120  double mass[6] = {0};
121  compute_mass(mass, TotNumPartCDM, TotNumPartGas, TotNuPart, nufrac, BoxSize, CP, GenicConfig);
122 
123  double redshift = 1.0 / GenicConfig.TimeIC - 1.;
124 
125  int rt =(0 != big_block_set_attr(&bheader, "TotNumPart", totnumpart, "i8", 6)) ||
126  (0 != big_block_set_attr(&bheader, "MassTable", mass, "f8", 6)) ||
127  (big_block_set_attr(&bheader, "Time", &GenicConfig.TimeIC, "f8", 1)) ||
128  (big_block_set_attr(&bheader, "Redshift", &redshift, "f8", 1)) ||
129  (big_block_set_attr(&bheader, "BoxSize", &BoxSize, "f8", 1)) ||
130  (big_block_set_attr(&bheader, "UsePeculiarVelocity", &GenicConfig.UsePeculiarVelocity, "i4", 1)) ||
131  (big_block_set_attr(&bheader, "Omega0", &CP->Omega0, "f8", 1)) ||
132  (big_block_set_attr(&bheader, "FractionNuInParticles", &nufrac, "f8", 1)) ||
133  (big_block_set_attr(&bheader, "OmegaBaryon", &CP->OmegaBaryon, "f8", 1)) ||
134  (big_block_set_attr(&bheader, "OmegaLambda", &CP->OmegaLambda, "f8", 1)) ||
135  (big_block_set_attr(&bheader, "UnitLength_in_cm", &GenicConfig.units.UnitLength_in_cm, "f8", 1)) ||
136  (big_block_set_attr(&bheader, "UnitMass_in_g", &GenicConfig.units.UnitMass_in_g, "f8", 1)) ||
137  (big_block_set_attr(&bheader, "UnitVelocity_in_cm_per_s", &GenicConfig.units.UnitVelocity_in_cm_per_s, "f8", 1)) ||
138  (big_block_set_attr(&bheader, "HubbleParam", &CP->HubbleParam, "f8", 1)) ||
139  (big_block_set_attr(&bheader, "InvertPhase", &GenicConfig.InvertPhase, "i4", 1)) ||
140  (big_block_set_attr(&bheader, "Seed", &GenicConfig.Seed, "i8", 1)) ||
141  (big_block_set_attr(&bheader, "UnitaryAmplitude", &GenicConfig.UnitaryAmplitude, "i4", 1));
142  if(rt) {
143  endrun(0, "failed to create attr %s",
144  big_file_get_error_message());
145  }
146 
147  if(0 != big_block_mpi_close(&bheader, MPI_COMM_WORLD)) {
148  endrun(0, "failed to close block %s:%s", "Header",
149  big_file_get_error_message());
150  }
151 }
void compute_mass(double *mass, int64_t TotNumPartCDM, int64_t TotNumPartGas, int64_t TotNuPart, double nufrac, const double BoxSize, Cosmology *CP, const struct genic_config GenicConfig)
Definition: save.c:90
double HubbleParam
Definition: cosmology.h:20
double OmegaLambda
Definition: cosmology.h:14
int UnitaryAmplitude
Definition: allvars.h:23
int Seed
Definition: allvars.h:22
int InvertPhase
Definition: allvars.h:24
int UsePeculiarVelocity
Definition: allvars.h:39
double TimeIC
Definition: allvars.h:38

References compute_mass(), CP, endrun(), Cosmology::HubbleParam, genic_config::InvertPhase, Cosmology::Omega0, Cosmology::OmegaBaryon, Cosmology::OmegaLambda, genic_config::Seed, genic_config::TimeIC, genic_config::UnitaryAmplitude, UnitSystem::UnitLength_in_cm, UnitSystem::UnitMass_in_g, genic_config::units, UnitSystem::UnitVelocity_in_cm_per_s, and genic_config::UsePeculiarVelocity.

Referenced by main().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ write_particle_data()

void write_particle_data ( IDGenerator idgen,
const int  Type,
BigFile *  bf,
const uint64_t  FirstID,
const int  SavePrePos,
int  NumFiles,
int  NumWriters,
struct ic_part_data curICP 
)

Definition at line 61 of file save.c.

68 {
69  /* Write particles */
70  if(SavePrePos)
71  saveblock(bf, &curICP[0].PrePos, Type, "PrePosition", "f8", 3, idgen->NumPart, sizeof(curICP[0]), NumFiles, NumWriters);
72  saveblock(bf, &curICP[0].Density, Type, "ICDensity", "f4", 1, idgen->NumPart, sizeof(curICP[0]), NumFiles, NumWriters);
73  saveblock(bf, &curICP[0].Pos, Type, "Position", "f8", 3, idgen->NumPart, sizeof(curICP[0]), NumFiles, NumWriters);
74  saveblock(bf, &curICP[0].Vel, Type, "Velocity", "f4", 3, idgen->NumPart, sizeof(curICP[0]), NumFiles, NumWriters);
75  /*Generate and write IDs*/
76  uint64_t * ids = (uint64_t *) mymalloc("IDs", idgen->NumPart * sizeof(uint64_t));
77  memset(ids, 0, idgen->NumPart * sizeof(uint64_t));
78  int i;
79  #pragma omp parallel for
80  for(i = 0; i < idgen->NumPart; i++)
81  {
82  ids[i] = idgen_create_id_from_index(idgen, i) + FirstID;
83  }
84  saveblock(bf, ids, Type, "ID", "u8", 1, idgen->NumPart, sizeof(uint64_t), NumFiles, NumWriters);
85  myfree(ids);
86  walltime_measure("/Write");
87 }
static struct ic_part_data * curICP
Definition: glass.c:177
#define mymalloc(name, size)
Definition: mymalloc.h:15
#define myfree(x)
Definition: mymalloc.h:19
uint64_t idgen_create_id_from_index(IDGenerator *idgen, int index)
Definition: zeldovich.c:69
static void saveblock(BigFile *bf, void *baseptr, int ptype, const char *bname, const char *dtype, int items_per_particle, const int NumPart, ptrdiff_t elsize, int NumFiles, int NumWriters)
Definition: save.c:49
int NumPart
Definition: proto.h:14
#define walltime_measure(name)
Definition: walltime.h:8

References curICP, idgen_create_id_from_index(), myfree, mymalloc, IDGenerator::NumPart, saveblock(), and walltime_measure.

Referenced by main().

Here is the call graph for this function:
Here is the caller graph for this function: