1
0
Fork 0
mirror of https://github.com/cosmo-sims/monofonIC.git synced 2024-09-19 17:03:45 +02:00

first commit of the masked particle load for the upcoming Virgo consortium simulations

This commit is contained in:
Oliver Hahn 2022-02-04 00:33:54 +01:00
parent 762bef3aa2
commit 9ea23ea581
3 changed files with 263 additions and 20 deletions

View file

@ -20,7 +20,8 @@ DoFixing = no # do mode fixing à la Angulo&Pontzen (https://arxiv.o
DoInversion = no # invert phases (for paired simulations)
ParticleLoad = sc # particle load, can be 'sc' (1x), 'bcc' (2x) or 'fcc' (4x)
# (increases number of particles by given factor!), or 'glass'
# (increases number of particles by given factor!),
# or 'glass' or 'masked'
## if `ParticleLoad = glass' then specify here where to load the glass distribution from
# GlassFileName = glass128.hdf5

View file

@ -29,13 +29,20 @@ namespace particle
enum lattice
{
lattice_glass = -1,
lattice_masked = -2, // masked lattices are SC lattices with a specifiable mask leaving out particles
lattice_glass = -1, // glass: needs a glass file
lattice_sc = 0, // SC : simple cubic
lattice_bcc = 1, // BCC: body-centered cubic
lattice_fcc = 2, // FCC: face-centered cubic
lattice_rsc = 3, // RSC: refined simple cubic
};
const std::vector<std::vector<bool>> lattice_masks =
{
// mask from Richings et al. https://arxiv.org/pdf/2005.14495.pdf
{true,true,true,true,true,true,true,false},
};
const std::vector<std::vector<vec3_t<real_t>>> lattice_shifts =
{
// first shift must always be zero! (otherwise set_positions and set_velocities break)
@ -149,11 +156,40 @@ namespace particle
private:
particle::container particles_;
size_t global_num_particles_;
static constexpr int masksize_ = 2;
std::array<int, masksize_*masksize_*masksize_> particle_type_mask_;
inline int get_mask_value( const vec3_t<size_t>& global_idx_3d ) const
{
int sig = ((global_idx_3d[0]%masksize_)*masksize_
+global_idx_3d[1]%masksize_)*masksize_
+global_idx_3d[2]%masksize_;
return particle_type_mask_[sig];
}
public:
lattice_generator(lattice lattice_type, const bool b64reals, const bool b64ids, const bool bwithmasses, size_t IDoffset, const field_t &field, config_file &cf)
/**
* @brief Construct a new lattice generator object
*
* @param lattice_type Type of the lattice (currently: 0=SC, 1=BCC, 2=FCC, 3=double SC, -1=glass, -2=masked SC)
* @param lattice_index Index of the 'atom type', i.e. whether CDM, baryons, ... (currently only supports 0 or 1)
* @param b64reals Boolean whether 64bit floating point shall be used to store positions/velocities/masses
* @param b64ids Boolean whether 64bit integers should be used for IDS
* @param bwithmasses Boolean whether distinct masses for each particle shall be stored
* @param IDoffset The global ID offset to be applied to this generation of particles
* @param field Reference to the field from which particles shall be generated (used only to set dimensions)
* @param cf Reference to the config_file object
*/
lattice_generator(lattice lattice_type, int lattice_index, const bool b64reals, const bool b64ids, const bool bwithmasses, size_t IDoffset, const field_t &field, config_file &cf)
: global_num_particles_(0)
{
if (lattice_type != lattice_glass)
// initialise the particle mask with zeros (only used if lattice_type==lattice_masked)
for( auto& m : particle_type_mask_) m = 0;
if (lattice_type >= 0) // These are the Bravais lattices
{
music::wlog << "Glass ICs will currently be incorrect due to disabled ghost zone updates! ";
@ -163,6 +199,9 @@ namespace particle
const size_t overload = 1ull << std::max<int>(0, lattice_type); // 1 for sc, 2 for bcc, 4 for fcc, 8 for rsc
// allocate memory for all local particles
particles_.allocate(overload * num_p_in_load, b64reals, b64ids, bwithmasses);
// set the global number of particles for this lattice_type and lattice_index
global_num_particles_ = field.global_size() * overload;
// set particle IDs to the Lagrangian coordinate (1D encoded) with additionally the field shift encoded as well
IDoffset = IDoffset * overload * field.global_size();
@ -188,7 +227,82 @@ namespace particle
}
}
}
else if( lattice_type == lattice_masked )
{
if( field.global_size()%8 != 0 ){
music::elog << "For masked lattice type, linear field resolution must be a multiple of two." << std::endl;
abort();
}
if( lattice_index > 1 || lattice_index < 0 ){
music::elog << "For masked lattice type, lattice index must be 0 or 1." << std::endl;
abort();
}
// set the particle mask
if( cf.get_value_safe("setup","DoBaryons",false) )
{
// mask everywehere 0, except the last element
for( auto& m : particle_type_mask_) m = 0;
particle_type_mask_[7] = 1;
}else{
// mask everywehere 0
for( auto& m : particle_type_mask_) m = 0;
}
// count number of particles taking into account masking
size_t ipcount = 0;
for (size_t i = 0; i < field.rsize(0); ++i)
{
for (size_t j = 0; j < field.rsize(1); ++j)
{
for (size_t k = 0; k < field.rsize(2); ++k)
{
if( this->get_mask_value(field.get_cell_idx_3d(i,j,k)) != lattice_index ) continue;
++ipcount;
}
}
}
// set global number of particles
#if defined(USE_MPI)
MPI_Allreduce( &ipcount, &global_num_particles_, 1, MPI_UNSIGNED_LONG_LONG, MPI_SUM, MPI_COMM_WORLD );
#else
global_num_particles_ = ipcount;
#endif
// number of modes present in the field
const size_t num_p_in_load = ipcount;
// allocate memory for all local particles
particles_.allocate(num_p_in_load, b64reals, b64ids, bwithmasses);
// set particle IDs to the Lagrangian coordinate (1D encoded) with additionally the field shift encoded as well
IDoffset = IDoffset * field.global_size();
for (size_t i = 0, ipcount = 0; i < field.rsize(0); ++i)
{
for (size_t j = 0; j < field.rsize(1); ++j)
{
for (size_t k = 0; k < field.rsize(2); ++k)
{
if( this->get_mask_value(field.get_cell_idx_3d(i,j,k)) != lattice_index ) continue;
if (b64ids)
{
particles_.set_id64(ipcount, IDoffset + field.get_cell_idx_1d(i, j, k));
}
else
{
particles_.set_id32(ipcount, IDoffset + field.get_cell_idx_1d(i, j, k));
}
++ipcount;
}
}
}
}
else if( lattice_type == lattice_glass )
{
glass_ptr_ = std::make_unique<glass>( cf, field );
particles_.allocate(glass_ptr_->size(), b64reals, b64ids, false);
@ -206,14 +320,24 @@ namespace particle
}
}
}
music::ilog << "Created Particles [" << lattice_index << "] : " << global_num_particles_ << std::endl;
}
// invalidates field, phase shifted to unspecified position after return
void set_masses(const lattice lattice_type, bool is_second_lattice, const real_t munit, const bool b64reals, field_t &field, config_file &cf)
void set_masses(const lattice lattice_type, int lattice_index, const real_t munit, const bool b64reals, field_t &field, config_file &cf)
{
// works only for Bravais types
if (lattice_type >= 0)
// works only for Bravais types and masked type
assert( lattice_type>=0 || lattice_type==lattice_masked );
if (lattice_type >= 0) // Bravais lattices
{
if( lattice_index > 1 || lattice_index < 0 ){
music::elog << "For Bravais lattice type, lattice index must be 0 or 1." << std::endl;
abort();
}
const size_t overload = 1ull << std::max<int>(0, lattice_type); // 1 for sc, 2 for bcc, 4 for fcc, 8 for rsc
const size_t num_p_in_load = field.local_size();
const real_t pmeanmass = munit / real_t(field.global_size()* overload);
@ -225,7 +349,7 @@ namespace particle
for (int ishift = 0; ishift < (1 << lattice_type); ++ishift)
{
// if we are dealing with the secondary lattice, apply a global shift
if (ishift == 0 && is_second_lattice)
if (ishift == 0 && lattice_index > 0)
{
field.shift_field(second_lattice_shift[lattice_type]);
}
@ -261,7 +385,49 @@ namespace particle
music::ilog << "Particle Mass : mean/munit = " << mean_pm/munit << " ; fractional RMS = " << std_pm / mean_pm * 100.0 << "%" << std::endl;
if(std_pm / mean_pm > 0.1 ) music::wlog << "Particle mass perturbation larger than 10%, consider decreasing \n\t the starting redshift or disabling baryon decaying modes." << std::endl;
if(bmass_negative) music::elog << "Negative particle mass produced! Decrease the starting \n\t redshift or disable baryon decaying modes!" << std::endl;
} else if( lattice_type == lattice_masked ) {
if( field.global_size()%8 != 0 ){
music::elog << "For masked lattice type, linear field resolution must be a multiple of two." << std::endl;
abort();
}
if( lattice_index > 1 || lattice_index < 0 ){
music::elog << "For masked lattice type, lattice index must be 0 or 1." << std::endl;
abort();
}
const real_t pmeanmass = munit / global_num_particles_;
bool bmass_negative = false;
auto mean_pm = field.mean() * pmeanmass;
auto std_pm = field.std() * pmeanmass;
// read out values from phase shifted field and set assoc. particle's value
for (size_t i = 0, ipcount = 0; i < field.size(0); ++i)
{
for (size_t j = 0; j < field.size(1); ++j)
{
for (size_t k = 0; k < field.size(2); ++k)
{
if( this->get_mask_value(field.get_cell_idx_3d(i,j,k)) != lattice_index ) continue;
// get
const auto pmass = pmeanmass * field.relem(i, j, k);
// check for negative mass
bmass_negative |= pmass<0.0;
// set
if (b64reals) particles_.set_mass64(ipcount++, pmass);
else particles_.set_mass32(ipcount++, pmass);
}
}
}
// diagnostics
music::ilog << "Particle Mass : mean/munit = " << mean_pm/munit << " ; fractional RMS = " << std_pm / mean_pm * 100.0 << "%" << std::endl;
if(std_pm / mean_pm > 0.1 ) music::wlog << "Particle mass perturbation larger than 10%, consider decreasing \n\t the starting redshift or disabling baryon decaying modes." << std::endl;
if(bmass_negative) music::elog << "Negative particle mass produced! Decrease the starting \n\t redshift or disable baryon decaying modes!" << std::endl;
}else{
// should not happen
music::elog << "Cannot have individual particle masses for glasses!" << std::endl;
@ -270,15 +436,21 @@ namespace particle
}
// invalidates field, phase shifted to unspecified position after return
void set_positions(const lattice lattice_type, bool is_second_lattice, int idim, real_t lunit, const bool b64reals, field_t &field, config_file &cf)
void set_positions(const lattice lattice_type, int lattice_index, int idim, real_t lunit, const bool b64reals, field_t &field, config_file &cf)
{
if (lattice_type >= 0)
if (lattice_type >= 0) // Bravais lattice
{
if( lattice_index > 1 || lattice_index < 0 ){
music::elog << "For Bravais lattice type, lattice index must be 0 or 1." << std::endl;
abort();
}
const size_t num_p_in_load = field.local_size();
for (int ishift = 0; ishift < (1 << lattice_type); ++ishift)
{
// if we are dealing with the secondary lattice, apply a global shift
if (ishift == 0 && is_second_lattice)
if (ishift == 0 && lattice_index==1)
{
field.shift_field(second_lattice_shift[lattice_type]);
}
@ -296,7 +468,7 @@ namespace particle
{
for (size_t k = 0; k < field.size(2); ++k)
{
auto pos = field.template get_unit_r_shifted<real_t>(i, j, k, lattice_shifts[lattice_type][ishift] + (is_second_lattice ? second_lattice_shift[lattice_type] : vec3_t<real_t>{real_t(0.), real_t(0.), real_t(0.)}));
auto pos = field.template get_unit_r_shifted<real_t>(i, j, k, lattice_shifts[lattice_type][ishift] + (lattice_index==1 ? second_lattice_shift[lattice_type] : vec3_t<real_t>{real_t(0.), real_t(0.), real_t(0.)}));
if (b64reals)
{
particles_.set_pos64(ipcount++, idim, pos[idim] * lunit + field.relem(i, j, k));
@ -310,6 +482,38 @@ namespace particle
}
}
}
else if( lattice_type == lattice_masked )
{
if( field.global_size()%8 != 0 ){
music::elog << "For masked lattice type, linear field resolution must be a multiple of two." << std::endl;
abort();
}
if( lattice_index > 1 || lattice_index < 0 ){
music::elog << "For masked lattice type, lattice index must be 0 or 1." << std::endl;
abort();
}
for (size_t i = 0, ipcount = 0; i < field.size(0); ++i)
{
for (size_t j = 0; j < field.size(1); ++j)
{
for (size_t k = 0; k < field.size(2); ++k)
{
if( this->get_mask_value(field.get_cell_idx_3d(i,j,k)) != lattice_index ) continue;
// get position (in box units) of the current cell of 3d array 'field'
auto pos = field.template get_unit_r<real_t>(i, j, k);
// add the displacement to get the particle position
if (b64reals){
particles_.set_pos64(ipcount++, idim, pos[idim] * lunit + field.relem(i, j, k));
}else{
particles_.set_pos32(ipcount++, idim, pos[idim] * lunit + field.relem(i, j, k));
}
}
}
}
}
else
{
glass_ptr_->update_ghosts( field );
@ -330,15 +534,20 @@ namespace particle
}
}
void set_velocities(lattice lattice_type, bool is_second_lattice, int idim, const bool b64reals, field_t &field, config_file &cf)
void set_velocities(lattice lattice_type, int lattice_index, int idim, const bool b64reals, field_t &field, config_file &cf)
{
if (lattice_type >= 0)
if (lattice_type >= 0) // Bravais lattice
{
if( lattice_index > 1 || lattice_index < 0 ){
music::elog << "For Bravais lattice type, lattice index must be 0 or 1." << std::endl;
abort();
}
const size_t num_p_in_load = field.local_size();
for (int ishift = 0; ishift < (1 << lattice_type); ++ishift)
{
// if we are dealing with the secondary lattice, apply a global shift
if (ishift == 0 && is_second_lattice)
if (ishift == 0 && lattice_index==1)
{
field.shift_field(second_lattice_shift[lattice_type]);
}
@ -368,6 +577,35 @@ namespace particle
}
}
}
else if( lattice_type == lattice_masked )
{
if( field.global_size()%8 != 0 ){
music::elog << "For masked lattice type, linear field resolution must be a multiple of two." << std::endl;
abort();
}
if( lattice_index > 1 || lattice_index < 0 ){
music::elog << "For masked lattice type, lattice index must be 0 or 1." << std::endl;
abort();
}
for (size_t i = 0, ipcount = 0; i < field.size(0); ++i)
{
for (size_t j = 0; j < field.size(1); ++j)
{
for (size_t k = 0; k < field.size(2); ++k)
{
if( this->get_mask_value(field.get_cell_idx_3d(i,j,k)) != lattice_index ) continue;
if (b64reals){
particles_.set_vel64(ipcount++, idim, field.relem(i, j, k));
}else{
particles_.set_vel32(ipcount++, idim, field.relem(i, j, k));
}
}
}
}
}
else
{
glass_ptr_->update_ghosts( field );

View file

@ -116,7 +116,8 @@ int run( config_file& the_config )
: ((lattice_str=="fcc")? particle::lattice_fcc
: ((lattice_str=="rsc")? particle::lattice_rsc
: ((lattice_str=="glass")? particle::lattice_glass
: particle::lattice_sc))));
: ((lattice_str=="masked")? particle::lattice_masked
: particle::lattice_sc)))));
//--------------------------------------------------------------------------------------------------------
//! apply fixing of the complex mode amplitude following Angulo & Pontzen (2016) [https://arxiv.org/abs/1603.05253]
@ -536,8 +537,11 @@ int run( config_file& the_config )
size_t IDoffset = (this_species == cosmo_species::baryon)? ((the_output_plugin->has_64bit_ids())? 1 : 1): 0 ;
// allocate particle structure and generate particle IDs
bool secondary_lattice = (this_species == cosmo_species::baryon &&
the_output_plugin->write_species_as(this_species) == output_type::particles) ? true : false;
particle_lattice_generator_ptr =
std::make_unique<particle::lattice_generator<Grid_FFT<real_t>>>( lattice_type, the_output_plugin->has_64bit_reals(), the_output_plugin->has_64bit_ids(),
std::make_unique<particle::lattice_generator<Grid_FFT<real_t>>>( lattice_type, secondary_lattice, the_output_plugin->has_64bit_reals(), the_output_plugin->has_64bit_ids(),
bDoBaryons, IDoffset, tmp, the_config );
}
@ -545,7 +549,7 @@ int run( config_file& the_config )
if( bDoBaryons && (the_output_plugin->write_species_as( this_species ) == output_type::particles
|| the_output_plugin->write_species_as( this_species ) == output_type::field_lagrangian) )
{
bool shifted_lattice = (this_species == cosmo_species::baryon &&
bool secondary_lattice = (this_species == cosmo_species::baryon &&
the_output_plugin->write_species_as(this_species) == output_type::particles) ? true : false;
const real_t munit = the_output_plugin->mass_unit();
@ -568,7 +572,7 @@ int run( config_file& the_config )
});
if( the_output_plugin->write_species_as( this_species ) == output_type::particles ){
particle_lattice_generator_ptr->set_masses( lattice_type, shifted_lattice, 1.0, the_output_plugin->has_64bit_reals(), rho, the_config );
particle_lattice_generator_ptr->set_masses( lattice_type, secondary_lattice, 1.0, the_output_plugin->has_64bit_reals(), rho, the_config );
}else if( the_output_plugin->write_species_as( this_species ) == output_type::field_lagrangian ){
the_output_plugin->write_grid_data( rho, this_species, fluid_component::mass );
}