/*******************************************************************\ particle_generator.hh - This file is part of MUSIC2 - a code to generate initial conditions for cosmological simulations CHANGELOG (only majors, for details see repo): 10/2019 - Oliver Hahn - first implementation \*******************************************************************/ #pragma once #include #include #if defined(USE_HDF5) #include "HDF_IO.hh" #endif namespace particle { enum lattice{ lattice_glass = -1, 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> > lattice_shifts = { // first shift must always be zero! (otherwise set_positions and set_velocities break) /* SC : */ {{0.0,0.0,0.0}}, /* BCC: */ {{0.0,0.0,0.0},{0.5,0.5,0.5}}, /* FCC: */ {{0.0,0.0,0.0},{0.0,0.5,0.5},{0.5,0.0,0.5},{0.5,0.5,0.0}}, /* RSC: */ {{0.0,0.0,0.0},{0.0,0.0,0.5},{0.0,0.5,0.0},{0.0,0.5,0.5},{0.5,0.0,0.0},{0.5,0.0,0.5},{0.5,0.5,0.0},{0.5,0.5,0.5}}, }; const std::vector> second_lattice_shift = { /* SC : */ {0.5, 0.5, 0.5}, // this corresponds to CsCl lattice /* BCC: */ {0.5, 0.5, 0.0}, // is there a diatomic lattice with BCC base?!? /* FCC: */ {0.5, 0.5, 0.5}, // this corresponds to NaCl lattice // /* FCC: */ {0.25, 0.25, 0.25}, // this corresponds to Zincblende/GaAs lattice /* RSC: */ {0.25, 0.25, 0.25}, }; template void initialize_lattice( container& particles, lattice lattice_type, const bool b64reals, const bool b64ids, const size_t IDoffset, const field_t& field, config_file& cf ){ if( lattice_type != lattice_glass ) { // number of modes present in the field const size_t num_p_in_load = field.local_size(); // unless SC lattice is used, particle number is a multiple of the number of modes (=num_p_in_load): const size_t overload = 1ull<(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 ); // set particle IDs to the Lagrangian coordinate (1D encoded) with additionally the field shift encoded as well for( size_t i=0,ipcount=0; i("setup","GlassFileName"); std::vector glass_dims; HDFGetDatasetExtent( glass_fname, "/PartType1/Coordinates", glass_dims ); music::ilog << "Glass file contains " << glass_dims[0] << " particles." << std::endl; size_t ntiles = cf.get_value("setup","GlassTiles"); size_t num_p = glass_dims[0] * ntiles*ntiles*ntiles / MPI::get_size(); size_t off_p = MPI::get_rank() * num_p; particles.allocate( num_p, b64reals, b64ids ); for( size_t i=0; i void set_positions( container& particles, const lattice lattice_type, bool is_second_lattice, int idim, real_t lunit, const bool b64reals, field_t& field, config_file& cf ) { // works only for Bravais types if( lattice_type >= 0 ){ const size_t num_p_in_load = field.local_size(); for( int ishift=0; ishift<(1<0 ){ field.shift_field( lattice_shifts[lattice_type][ishift] - lattice_shifts[lattice_type][ishift-1] ); } // read out values from phase shifted field and set assoc. particle's value const auto ipcount0 = ishift * num_p_in_load; for( size_t i=0,ipcount=ipcount0; i(i,j,k,lattice_shifts[lattice_type][ishift] + (is_second_lattice? second_lattice_shift[lattice_type] : vec3_t{0.,0.,0.}) ); 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{ #if defined(USE_HDF5) std::string glass_fname = cf.get_value("setup","GlassFileName"); size_t ntiles = cf.get_value("setup","GlassTiles"); real_t lglassbox = 1.0; HDFReadGroupAttribute( glass_fname, "Header", "BoxSize", lglassbox ); std::vector glass_pos; HDFReadDataset( glass_fname, "/PartType1/Coordinates", glass_pos ); size_t np_in_file = glass_pos.size()/3; size_t num_p = np_in_file * ntiles*ntiles*ntiles / MPI::get_size(); size_t off_p = num_p * MPI::get_rank(); std::vector< std::array > glass_posr(num_p,{0.0,0.0,0.0}); std::array ng({field.n_[0],field.n_[1],field.n_[2]}); for( size_t i=0; i interp( field ); interp.domain_decompose_pos( glass_posr ); for( size_t i=0; i void set_velocities(container &particles, lattice lattice_type, bool is_second_lattice, int idim, const bool b64reals, field_t &field, config_file& cf) { // works only for Bravais types if( lattice_type >= 0 ){ const size_t num_p_in_load = field.local_size(); for( int ishift=0; ishift<(1< 0){ field.shift_field( lattice_shifts[lattice_type][ishift]-lattice_shifts[lattice_type][ishift-1] ); } // read out values from phase shifted field and set assoc. particle's value const auto ipcount0 = ishift * num_p_in_load; for( size_t i=0,ipcount=ipcount0; i("setup","GlassFileName"); size_t ntiles = cf.get_value("setup","GlassTiles"); real_t lglassbox = 1.0; HDFReadGroupAttribute( glass_fname, "Header", "BoxSize", lglassbox ); std::vector glass_pos; HDFReadDataset( glass_fname, "/PartType1/Coordinates", glass_pos ); size_t np_in_file = glass_pos.size()/3; size_t num_p = np_in_file * ntiles*ntiles*ntiles / MPI::get_size(); size_t off_p = num_p * MPI::get_rank(); std::vector< std::array > glass_posr(num_p,{0.0,0.0,0.0}); std::array ng({field.n_[0],field.n_[1],field.n_[2]}); for( size_t i=0; i interp( field ); interp.domain_decompose_pos( glass_posr ); for( size_t i=0; i