/*! \file ui.cxx * \brief user interface */ #include "io.h" #include "ioutils.h" #include "logging.h" #include "stf.h" static void configure_logging(const Options &opt) { // VR's verbosity goes from 0 (less verbose) to at least 2 (very verbose) // We thus consider 2 as "trace", 1 as "debug" and so forth constexpr int NUM_LOG_LEVELS = 5; constexpr int VR_MIN_LOG_LEVEL = 2; int as_log_level = -opt.iverbose + VR_MIN_LOG_LEVEL; vr::LogLevel level = static_cast(std::min(std::max(as_log_level, 0), NUM_LOG_LEVELS - 1)); vr::init_logging(level); } static void log_options_summary(const Options &opt) { LOG(info) << "CONFIG INFO SUMMARY --------------------------"; const char *search_type; switch(opt.partsearchtype) { case PSTALL: search_type = "Searching all particles regardless of type"; break; case PSTDARK: search_type = "Searching dark matter particles"; break; case PSTGAS: search_type = "Searching gas particles"; break; case PSTSTAR: search_type = "Searching star particles"; break; case PSTBH: search_type = "Searching BH particles?! Really? Are there enough?"; break; default: search_type = "Unknown search type"; } LOG(info) << search_type; if (opt.fofbgtype == FOF6D) { LOG(info) << "Field objects found with 3d FOF"; } else if (opt.fofbgtype == FOF6D) { LOG(info) <<"Field objects found with initial 3d FOF then use single dispersion to find 6d FOFs"; } else if (opt.fofbgtype == FOF6DADAPTIVE) { LOG(info) << "Field objects found with initial 3d FOF then use adaptive dispersion to find 6d FOFs"; } else if (opt.fofbgtype == FOFSTNOSUBSET) { LOG(info) << "Field objects found with initial 3d FOF then use phase-space stream finding algorithm"; } if (opt.fofbgtype <= FOF6D && opt.iKeepFOF) { LOG(info) << "Field objects found with initial 3d FOF then use dispersion to find 6d FOFs and 3dFOF objects kept as base objects (consider them inter halo stellar mass or ICL for stellar only searches) "; } if (opt.iBaryonSearch == 1) { LOG(info) << "Baryons treated separately for substructure search"; } else if (opt.iBaryonSearch == 2) { LOG(info) << "Baryons treated separately for substructure search and also FOF field search"; } if (opt.iSingleHalo) { LOG(info) << "Field objects NOT searched for, assuming single Halo and subsearch using mean field first step"; } LOG(info) << "Allowed potential to kinetic ratio when unbinding particles: " << opt.uinfo.Eratio; if (opt.HaloMinSize != opt.MinSize) { LOG(info) << "Field objects (aka Halos) have different minimum required size than substructures: " << opt.HaloMinSize << " vs " << opt.MinSize; } LOG(info) << "Units: L=" << opt.lengthinputconversion << ", M=" << opt.massinputconversion << ", V=" << opt.velocityinputconversion << ", G=" << opt.G; if (opt.ibinaryout) { LOG(info) << "Binary output"; } if (opt.iseparatefiles) { LOG(info) << "Separate files output"; } if (opt.icomoveunit) { LOG(info) << "Converting properties into comoving, little h units"; } if (opt.iextrahalooutput) { LOG(info) << "Calculating extra halo properties"; } if (opt.iextragasoutput) { LOG(info) << "Calculating extra gas properties"; } if (opt.iextrastaroutput) { LOG(info) << "Calculating extra star properties"; } if (opt.iaperturecalc) { std::ostringstream os; os << "Calculating quantities in apertures: "; for (auto apname: opt.aperture_names_kpc) { os << apname << ","; } LOG(info) << os.str(); } if (opt.iprofilecalc) { std::ostringstream os; os << "Calculating radial profiles: "; for (auto binedge: opt.profile_bin_edges) { os << binedge << ","; } LOG(info) << os.str(); } if (opt.iextendedoutput) { LOG(info) << "Writing extended output for particle extraction from input files"; } if (opt.iSphericalOverdensityPartList) { LOG(info) << "Writing list of particles in SO regions of halos. "; } if (opt.iHaloCoreSearch) { LOG(info) << "Searching for 6dfof cores so as to disentangle mergers"; } if (opt.iHaloCoreSearch && opt.iAdaptiveCoreLinking) { LOG(info) << "With adaptive linking lengths"; } if (opt.iHaloCoreSearch && opt.iPhaseCoreGrowth) { LOG(info) << "With with phase-space tensor core assignment"; } const char *input_type; switch (opt.inputtype) { case IOGADGET: input_type = "Gadget"; break; case IOTIPSY: input_type = "Tipsy"; break; case IORAMSES: input_type = "RAMSES"; break; #ifdef USEHDF case IOHDF: input_type = "HDF"; break; #endif #ifdef USEXDR case IONCHILADA: input_type = "NCHILADA"; break; #endif default: input_type = "Unknown"; break; } LOG(info) << input_type << " file particle input"; #ifdef USEMPI LOG(info) << "MPI input communication requires allocation of temporary particle buffer of size " << vr::memory_amount(opt.mpiparticletotbufsize); #endif LOG(info) << "----------------------------------------------"; } ///routine to get arguments from command line /// \todo alter interface as now need to be able to specify only smdata file (no grid, res, normalized res) and functionality to specify eps, background fof search, etc void GetArgs(int argc, char *argv[], Options &opt) { #if defined(USEMPI) && defined(USEPARALLELHDF) opt.mpinprocswritesize=NProcs; #endif int option; int NumArgs = 0; int configflag=0; while ((option = getopt(argc, argv, ":C:I:i:s:Z:o:G:S:B:t:")) != EOF) { switch(option) { case 'C': opt.pname = optarg; configflag=1; NumArgs += 2; break; case 'I': opt.inputtype = atoi(optarg); NumArgs += 2; break; case 'i': opt.fname = optarg; NumArgs += 2; break; case 's': opt.num_files = atoi(optarg); NumArgs += 2; break; case 'Z': opt.nsnapread = atoi(optarg); NumArgs += 2; break; case 'o': opt.outname = optarg; NumArgs += 2; break; case 'G': opt.gnsphblocks = atoi(optarg); NumArgs += 2; break; case 'S': opt.gnstarblocks = atoi(optarg); NumArgs += 2; break; case 'B': opt.gnbhblocks = atoi(optarg); NumArgs += 2; break; case 't': opt.ramsessnapname = optarg; NumArgs += 2; break; case '?': usage(); } } if(configflag){ LOG_RANK0(info) << "Reading config file " << opt.pname; GetParamFile(opt); } else { LOG_RANK0(warning) << "No configuration file passed, using default values"; } #ifdef USEMPI MPI_Barrier(MPI_COMM_WORLD); #endif ConfigCheck(opt); } ///Outputs the usage to stdout void usage(void) { #ifndef USEMPI int ThisTask =0, NProcs =1; #endif Options opt; if (ThisTask==0) { cerr<<"USAGE:\n"; cerr<<"\n"; cerr<<"-C "<"< "<"<"<"<"<"<"<"<sample). within the repository. Only the keywords listed here will be used, all other words/characters are ignored. Note that if misspell a keyword it will not be used. One can check the options used by examining foo.configuration, where foo is your base output filename \ref Options.outname. This file lists all the options. Since this file is always written DO NOT name your input configuration file foo.configuration. There are numerous key words that can be passed to change the values stored in \ref Options. Here we list them \section outconfig Output related. See \ref outputs for types of output (and \ref io.cxx for implementation of the code) \arg \e Output Output base name. Overrides the name passed with the command line argument \e -o . Only implemented for completeness. \ref Options.outname \n \arg \e Write_group_array_file When producing output also produce a file which lists for every particle the group they belong to. Can be used with \b tipsy format or to tag every particle. \ref Options.iwritefof \arg \e Output_den A filename for storing the intermediate step of calculating local densities. This is particularly useful if the code is not compiled with \b STRUCDEN & \b HALOONLYDEN (see \ref STF-makeflags). \ref Options.smname \n \section searchconfig Parameters related to search type. See \ref io.cxx (and related ios like \ref gadgetio.cxx), \ref search.cxx, \ref fofalgo.h for extra details \arg \e Particle_search_type A integer describing what types of particles are searched \ref Options.partsearchtype. (See \ref SEARCHTYPES) > \n - \b 1 \e all particles are searched - \b 2 \e DarkMatter particles (which are typically defined as type 1,2,3 for gadget) are searched - \b 3 \e Star particles (which are typically defined as type 4 for gadget) are searched - \b 4 \e Gas particles (which are typically defined as type 0 for gadget) are searched \n \arg \e Baryon_searchflag 0/1/2 indicating gas/stellar search done separately from DM search \ref Options.iBaryonSearch. \n - \b 1 field search run as normal and then substructure search for baryons run using baryons identified in field search. - \b 2 field search also altered to treat baryons differently, allowing only DM particles to be used as head links (ie link dm-dm, dm-baryon, but not baryon-baryon nor baryon-dm). - \b 0 is do nothing special for baryon particles. \n \arg \e Singlehalo_search 0/1 flag indicates that no field search is going to be run and the entire volume will be treated as a background region. Useful if searching for substructures in non-cosmological simulations. But can also be co-opted for other searches using different outlier criteria and FOF algorithms. \ref Options.iSingleHalo \n \section localdensityconfig Parameters related to local density estimator. See \ref localfield.cxx, \ref bgfield.cxx & \ref localbgcomp.cxx for more details \arg \e Nsearch_velocity number of velocity neighbours used to calculate velocity density, adjust \ref Options.Nvel (suggested value is 32) \n \arg \e Nsearch_physical number of physical neighbours searched for Nv to calculate velocity density \ref Options.Nsearch (suggested value is 256) \n \arg \e Cell_fraction fraction of a halo contained in a subvolume used to characterize the background \ref Options.Ncellfac \n \arg \e Grid_type integer describing type of grid used to decompose volume for substructure search \ref Options.gridtype (see \ref GRIDTYPES) \n - \b 1 \e standard physical shannon entropy, balanced KD tree volume decomposition into cells - \b 2 \e phase phase-space shannon entropy, balanced KD tree volume decomposition into cells - \b 3 \e simple simple physical balanced KD tree decomposition of volume into cells \section fofconfig Parameters related to FOF search See \ref search.cxx \ref fofalgo.h for more details \subsection fofsubconfig Configuration for substructure search \arg \e Search_for_substructure By default field objects are searched for internal substructures but can disable this by setting this to 0 \n \arg \e Keep_FOF if field 6DFOF search is done, allows to keep structures found in 3DFOF (can be interpreted as the inter halo stellar mass when only stellar search is used).\n \arg \e FoF_search_type There are several substructure FOF criteria implemented (see \ref FOFTYPES for more types and \ref fofalgo.h for implementation) \n - \b 1 \e standard phase-space based, well tested VELOCIraptor criterion. \arg \e Outlier_threshold threshold of sigma level of outliers to be searched (default is 2.5) \ref Options.ellthreshold \n \arg \e Significance_level minimum significance level of group (default is 1) \ref Options.siglevel \n \arg \e Velocity_ratio speed ratio used in linking particles \ref Options.Vratio \n \arg \e Velocity_opening_angle angle between velocities when linking (in units of \f$ \pi \f$) \ref Options.thetaopen \n \arg \e Physical_linking_length physical linking length used in fof (if cosmological gadget file then assumed to be in units of inter particle spacing \ref gadgetio.cxx, if loading in a single halo then based on average interparticle spacing calculated in \ref haloproperties.cxx). To be deprecated and replaced by Substructure_physical_linking_length \ref Options.ellphys \n \arg \e Substructure_physical_linking_length physical linking length used in fof (if cosmological gadget file then assumed to be in units of inter particle spacing \ref gadgetio.cxx, if loading in a single halo then based on average interparticle spacing calculated in \ref haloproperties.cxx) \ref Options.ellphys \n \arg \e Minimum_size Minimum group (substructure) size \ref Options.MinSize \n \arg \e Iterative_searchflag 1/0 use interactive substructure search which is designed to first identify spatially compact candidate outlier regions and then relaxes the criteria to find the more diffuse (in phase-space) regions associate with these candidate structures \ref Options.iiterflag \n \arg \e Iterative_threshold_factor factor multiplied with \ref Options.ellthreshold when using iterative method and identifying outlier regions associated with the initial candidate list of spatially compact outlier groups. Typical values are \f$ \sim 1 \f$ \ref Options.ellfac \n \arg \e Iterative_linking_length_factor factor multiplied with \ref Options.ellphys when using iterative method and identifying outlier regions associated with the initial candidate list of spatially compact outlier groups. Typical values are \f$ \sim 2 \f$ \ref Options.ellxfac \n \arg \e Iterative_Vratio_length_factor factor multiplied with \ref Options.Vratio when using iterative method and identifying outlier regions associated with the initial candidate list of spatially compact outlier groups. Typical values are \f$ \sim 1 \f$ \ref Options.vfac \n \arg \e Iterative_ThetaOp_length_factor factor multiplied with \ref Options.thetaopen when using iterative method and identifying outlier regions associated with the initial candidate list of spatially compact outlier groups. Typical values are \f$ \sim 1 \f$ \ref Options.thetafac \n \arg \e CMrefadjustsubsearch_flag 1/0 flag indicating whether particles are moved to the rough CM velocity frame of the background before substructures are searched for. \ref Options.icmrefadjust \n \subsection foffieldconfig Configuration for field search \arg \e FoF_Field_search_type There are several FOF criteria implemented to search for so-called field objects (see \ref FOFTYPES for more types and \ref fofalgo.h for implementation) \ref Options.fofbgtype \n - \b 5 \e standard 3D FOF based algorithm - \b 4 \e standard 3D FOF based algorithm FOLLOWED by 6D FOF search using the velocity scale defined by the largest halo on ONLY particles in 3DFOF groups - \b 3 \e standard 3D FOF based algorithm FOLLOWED by 6D FOF search using the velocity scale for each 3DFOF group \arg \e Minimum_halo_size Allows field objects (or so-called halos) to require a different minimum size (typically would be <= \ref Options.MinSize. Default is -1 which sets it to \ref Options.MinSize) \ref Options.HaloMinSize \n \arg \e Halo_linking_length_factor allows one to use different physical linking lengths between field objects and substructures. (Typically for 3DFOF searches of dark matter haloes, set to value such that this times \ref Options.ellphys = 0.2 the interparticle spacing when examining cosmological simulations ) \ref Options.ellhalophysfac \n \arg \e Halo_velocity_linking_length_factor allows one to use different velocity linking lengths between field objects and substructures when using 6D FOF searches. (Since in such cases the general idea is to use the local velocity dispersion to define a scale, \f$ \geq5 \f$ times this value seems to correctly scale searches) \ref Options.ellhalovelfac \n \arg \e Halo_6D_linking_length_factor allows one to use different linking lengths between 3DFOF and 6DFOF field search. Typically values are \f$ \sim 1 \f$ \n \arg \e Halo_6D_vel_linking_length_factor scaling applied to dispersions used in 6DFOF field search. Typical values are \f$ \geq 1.25 \f$ \n \arg \e Halo_core_search 0/1/2 flag allows one to explicitly search for large 6D FOF cores that are indicative of a recent major merger. Since substructure is defined on the scale of the maximum cell size and major mergers typically result two or more phase-space dense regions that are \e larger than the cell size used in reasonable substructure searches, one can identify them using this search. The overall goal is to treat these objects differently than a substructure. However, if 2 is set, then smaller core is treated as substruture and all particles within the FOF envelop are assigned to the cores based on their phase-space distance to core particles \ref Options.iHaloCoreSearch \n \arg \e Use_adaptive_core_search 0/1 flag allows one to run complex adaptive phase-space search for large 6D FOF cores and then use these linking lengths to separate mergers. 0 is simple high density dispersively cold cores with v scale adaptive, 1 is adaptive with simple dispersions in both x and v. \arg \e Use_phase_tensor_core_growth 0/1 flag allows one to run complex phase-space growth of merger remnants (6D FOF cores found). 0 is assignment with simple x and v dispersion to nearest core particle, 1 is phase-space tensor distance assignemnt to CM of core. \arg \e Halo_core_ellx_fac scaling applied to linking length when identifying merger remnants. Typically values are \f$ \sim0.5 \f$ \ref Options.halocorexfac \arg \e Halo_core_ellv_fac scaling applied to local dispersion to define the velocity scale used to identify merger remnants. Typically values are \f$ \sim1 \f$ \ref Options.halocorevfac \arg \e Halo_core_ncellfac used to determine the minimum number of particles a merger remnants is composed of, specifically \f$ N_{\rm min}= f_{\rm ncell}* N_{\rm S} \f$. Typically values are 0.005 \arg \e Halo_core_adaptive_sigma_fac used when running fully adaptive core search, specifies the width of the physical linking length in configuration space dispersion (think of this as how many sigma to include). Typically values are \f$ \sim2 \f$. This has been tested on hydrodynamnical simulations to separate galaxy mergers. \ref Options.halocoresigmafac \arg \e Halo_core_num_loops allows the core search to iterate, shrinking the velocity linking length to used till the number of cores identified decreases or this limit is reached. Allows apative search with larger linking length to be robust. Typically values are \f$ \sim5 \f$ with loops only running often twice. \ref Options.halocorenumloops. \arg \e Halo_core_loop_ellx_fac Factor by which configuration linking length is decreased when running loops for core search. Typically values are \f$ \sim0.75 \f$. \ref Options.halocorexfaciter \arg \e Halo_core_loop_ellv_fac Factor by which velocity linking length is decreased when running loops for core search. Typically values are \f$ \sim0.75 \f$. \ref Options.halocorevfaciter \arg \e Halo_core_loop_elln_fac Factor by which min group size is changed when running loops for core search. Typically values are \f$ \sim1.25 \f$. \ref Options.halocorenumfaciter \arg \e Halo_core_phase_significance Significance a core must be in terms of phase-space distance scaled by dispersions (sigma). \section unbindconfig Unbinding Parameters See \ref unbinding and \ref unbind.cxx for more details \arg \e Unbind_flag whether or not substructures should be passed through an unbinding routine. \ref Options.uinfo & \ref UnbindInfo.unbindflag \n \arg \e Allowed_kinetic_potential_ratio ratio of kinetic to potential energy at which a particle is still considered bound, ie: particle is still bound if \f$ \alpha T+W<0 \f$, so \f$ \alpha=1 \f$ would be standard unbinding and \f$ \alpha<1 \f$ allows one to identify unbound tidal debris. Given that VELOCIraptor was designed to identify tidal streams, it makes little sense to have this set to 1 unless explicitly required. Note that the code still separates particles into bound and unbound. Typical values of \f$ \alpha\geq 0.2 \f$ seems to minimize the number of false positives in tidal debris while still identifying completely unbound tidal debris. \ref Options.uinfo & \ref UnbindInfo.Eratio \n \arg \e Min_bound_mass_frac Designed to demand a substructure still have a minimum amount of self-bound mass. \ref Options.uinfo & \ref UnbindInfo.minEfrac \n \arg \e Bound_halos 0/1/2 flag to make sure field objects such as haloes are self bound before (use 1) and also after (use 2) substructures have been identified and extracted from the halo. Demanding boundness after substructure search can have interesting consequences as it is possible that a multiple merger will appear as a single FOF halo, however all with all the cores removed, the FOF halo is actually an unbound structure. \ref Options.iBoundHalos \n \arg \e Keep_background_potential 1/0 flag When determining whether a structure is self-bound, the approach taken is to treat the candidate structure in isolation. Then determine the velocity reference frame to determine the kinetic energy of each particle and remove them. However, it is possible one wishes to keep the background particles when determining the potential, that is once one starts unbinding, don't treat the candidate structure in isolation but in a background sea. When finding tidal debris, it is useful to keep the background. \ref Options.uinfo & \ref UnbindInfo.bgpot \n \arg \e Kinetic_reference_frame_type specify kinetic frame when determining whether particle is bound. Default is to use the centre-of-mass velocity frame (0) but can also use region around minimum of the potential (1). \ref Options.uinfo & \ref UnbindInfo.cmvelreftype \n \arg \e Min_npot_ref Set the minimum number of particles used to calculate the velocity of the minimum of the potential (10). \ref Options.uinfo & \ref UnbindInfo.Npotref \n \arg \e Frac_pot_ref Set the fraction of particles used to calculate the velocity of the minimum of the potential (0.1). \ref Options.uinfo & \ref UnbindInfo.fracpotref \n \arg \e Unbinding_type Set the unbinding criteria, either just remove particles deemeed "unbound", that is those with \f$ \alpha T+W>0\f$, choosing \ref UPART. Or with \ref USYSANDPART removes "unbound" particles till system also has a true bound fraction > \ref UnbindInfo.minEfrac. \arg \e Softening_length Set the (simple plummer) gravitational softening length. \ref UnbindInfo.eps \section cosmoconfig Units & Cosmology \subsection unitconfig Units \arg \e Length_unit change the input length unit \ref Options.L \n \arg \e Velocity_unit change the input velocity unit \ref Options.V \n \arg \e Mass_unit change the input mass unit \ref Options.M \n \arg \e Gravity change the gravity unit. should be set such that \f$ v^2=GM/r \f$. \ref Options.G \n \arg \e Hubble_unit change the value of Hubble expansion (from normal 100 km/s/Mpc to units used, that is velocity unit/length unit) \ref Options.G \n \arg \e Mass_value if not mass is stored using the option \b NOMASS (see \ref STF-makeflags) then this is the mass of the particles \ref Options.MassValue \n \arg \e Length_unit_to_kpc Specify the conversion factor from the output unit to kpc \ref Options.lengthtokpc \n \arg \e Velocity_unit_to_kms Specify the conversion factor from the output unit to kms \ref Options.velocitytokms \n \arg \e Mass_unit_to_solarmass Specify the conversion factor from the output unit to solar masses \ref Options.masstosolarmass \n \subsection cosmologyconfig Cosmology \arg \e Period if not defined in the input data one can pass the period in the input units of the data. This is not always necessary, for instance the gadget snapshot format has this in the header. \ref Options.p \n \arg \e Scale_factor scale factor time if not defined in the input data. This is not always necessary, for instance the gadget snapshot format has this in the header. \ref Options.a \n \arg \e h_val the "little h" value often used in cosmological simulation when defining expansion rates, distance scales, etc. This is not always necessary, for instance the gadget snapshot format has this in the header. \ref Options.h \n \arg \e Omega_m matter density in units of the critical density used in cosmological simulations. This is not always necessary, for instance the gadget snapshot format has this in the header. \ref Options.Omega_m \n \arg \e Omega_Lambda energy density of the cosmological constant (or can technically be used for a dark energy fluid) in units of the critical density used in cosmological simulations. This is not always necessary, for instance the gadget snapshot format has this in the header. \ref Options.Omega_Lambda \n \arg \e Omega_cdm dark matter density in units of the critical density used in cosmological simulations. This is not always necessary since for simple pure DM cosmological simulations this is easily determined by the total mass in the simulation. However, for multiple resolution or non-standard DM models, this should be provided. \ref Options.Omega_cdm \n \arg \e Omega_b baryon density in units of the critical density used in cosmological simulations. This is not always necessary since for simple pure adiabatic cosmological simulations this is easily determined by the total gas mass in the simulation. \ref Options.Omega_b \n \arg \e Omega_r radiation density in units of the critical density used in cosmological simulations. Typically 0 as negligible in most simulations. \ref Options.Omega_r \n \arg \e Omega_nu neutrino density in units of the critical density used in cosmological simulations. Typically 0 as negligible in most simulations and may need alteration for consistent evolution across cosmic time. \ref Options.Omega_nu \n \arg \e Omega_DE dark energy density in units of the critical density used in cosmological simulations. Typically 0 unless what evolution based on w_of_DE. \ref Options.Omega_de \n \arg \e w_of_DE equation of state of the dark energy fluid, \f$ w=\frac{p}{\rho} \f$. This is not necessary unless one is using a cosmological simulation with a \f$ w\neq -1 \f$ \ref Options.w_de \n \arg \e Virial_density virial overdensity in units of the background matter density used in cosmological simulations. This is not absolutely necessary for gadget format as for this input, the Bryan & Norman 1998 virial density is calculated based on a LCDM cosmology as \ref Options.virlevel is set to -1 and if left, a virial density is calculated. If set and a gadget input is used, this overrides the Bryan & Norman calculation. \n \arg \e Critical_density critical density in input units used in cosmological simulations. This is not absolutely necessary for gadget format as it can be calculated from the quantities in the header for typical GR cosmologies. \ref Options.rhobg \n \section otherconfigs Other configuration options \arg \e Effective_Resolution If running a multiple resolution cosmological zoom simulation, simple method of scaling the linking length by using the period, ie: \f$ p/N_{\rm eff} \f$ \ref Options.Neff \n \arg \e Snapshot_value If halo ids need to be offset to some starting value based on the snapshot of the output, which is useful for some halo merger tree codes, one can specific a snapshot number, and all halo ids will be listed as internal haloid + \f$ sn\times10^{12}\f$. \ref Options.snapshotvalue \n \arg \e Verbose 2/1/0 flag indicating how talkative the code is (2 very verbose, 1 verbose, 0 quiet). \ref Options.iverbose \n \section propconfigs Property calculation options \arg \e Inclusive_halo_mass 1/0 flag indicating whether inclusive masses are calculated for field objects. \ref Options.iInclusiveHalo \n \arg \e Extensive_halo_properties_output 1/0 flag indicating whether to calculate/output even more halo properties. \ref Options.iextrahalooutput \n \arg \e Extended_output 1/0 flag indicating whether produce extended output for quick particle extraction from input catalog of particles in structures \ref Options.iextendedoutput \n \arg \e Iterate_cm_flag 1/0 flag indicating whether to use shrinking spheres to calculate the center of mass and velocity. \ref Options.iIterateCM \n \arg \e Sort_by_binding_energy 1/0 flag indicating whether to sort by particle binding energy or by potential (if 0). \ref Options.iSortByBindingEnergy \n \section ioconfigs I/O options \arg \e Cosmological_input 1/0 indicating that input simulation is cosmological or not. With cosmological input, a variety of length/velocity scales are set to determine such things as the virial overdensity, linking length. \ref Options.icosmologicalin \n \arg \e Input_chunk_size Amount of information to read from input file in one go (100000). \ref Options.inputbufsize \n \arg \e Write_group_array_file 0/1 flag indicating whether write a single large tipsy style group assignment file is written. \ref Options.iwritefof \n \arg \e Separate_output_files 1/0 flag indicating whether separate files are written for field and subhalo groups. \ref Options.iseparatefiles \n \arg \e Binary_output 3/2/1/0 flag indicating whether output is hdf, binary or ascii. \ref Options.ibinaryout, \ref OUTADIOS, \ref OUTHDF, \ref OUTBINARY, \ref OUTASCII \n \arg \e Comoving_units 1/0 flag indicating whether the properties output is in physical or comoving little h units. \ref Options.icomoveunit \n \section inputflags input flags related to varies input formats \arg \e NSPH_extra_blocks If gadget snapshot is loaded one can specific the number of extra SPH blocks are read/in the file. \ref Options.gnsphblocks \n \arg \e NStar_extra_blocks If gadget snapshot is loaded one can specific the number of extra Star blocks are read/in the file. \ref Options.gnstarblocks \n \arg \e NBH_extra_blocks If gadget snapshot is loaded one can specific the number of extra Black hole blocks are read/in the file. \ref Options.gnbhblocks \n \arg \e HDF_name_convention HDF dataset naming convection. See \ref hdfitems.h for what naming conventions are available and what names exist. Currently have \ref HDFNUMNAMETYPES. \ref Options.ihdfnameconvention \n \arg \e Input_includes_dm_particle If dm particle specific information is in the input file. \ref Options.iusedmparticles \n \arg \e Input_includes_gas_particle If gas particle specific information is in the input file. \ref Options.iusegasparticles \n \arg \e Input_includes_star_particle If star particle specific information is in the input file. \ref Options.iusestarparticles \n \arg \e Input_includes_bh_particle If bh/sink particle specific information is in the input file. \ref Options.iusesinkparticles \n \arg \e Input_includes_wind_particle If wind particle specific information is in the input file. \ref Options.iusewindparticles \n \arg \e Input_includes_tracer_particle If tracer particle specific information is in the input file. \ref Options.iusetracerparticles \n \arg \e Input_includes_star_particle If star particle specific information is in the input file. \ref Options.iusestarparticles \n \section mpiconfigs MPI specific options \arg \e MPI_part_allocation_fac Memory allocated in mpi mode to store particles is (1+factor)* the memory need for the initial mpi decomposition. This factor should be >0 and is mean to allow a little room for particles to be exchanged between mpi threads withouth having to require new memory allocations and copying of data. \ref Options.mpipartfac \n \arg \e MPI_particle_total_buf_size Total memory size in bytes used to store particles in temporary buffer such that particles are sent to non-reading mpi processes in one communication round in chunks of size buffer_size/NProcs/sizeof(Particle). \ref Options.mpiparticlebufsize \n */ inline void ConfigMessage(char *c) { #ifndef USEMPI int ThisTask = 0; #endif if (ThisTask==0) cerr< ExtraFieldCalculationsAllowedOptions() { vector list; list.push_back(("averagemassweighted")); list.push_back(("average")); list.push_back(("logaveragemassweighted")); list.push_back(("logaverage")); list.push_back(("totalmassweighted")); list.push_back(("total")); list.push_back(("stdmassweighted")); list.push_back(("std")); list.push_back(("logstdmassweighted")); list.push_back(("logstd")); list.push_back(("minmassweighted")); list.push_back(("min")); list.push_back(("maxmassweighted")); list.push_back(("max")); list.push_back(("aperture_total")); list.push_back(("aperture_average")); return list; } inline map ExtraFieldCalculationsStringtoIntFlag() { map list; list[string("averagemassweighted")] = CALCAVERAGEMASSWEIGHT; list[string("average")] = CALCAVERAGE; list[string("logaveragemassweighted")] = CALCLOGAVERAGEMASSWEIGHT; list[string("logaverage")] = CALCLOGAVERAGE; list[string("totalmassweighted")] = CALCTOTALMASSWEIGHT; list[string("total")] = CALCTOTAL; list[string("stdmassweighted")] = CALCSTDMASSWEIGHT; list[string("std")] = CALCSTD; list[string("logstdmassweighted")] = CALCLOGSTDMASSWEIGHT; list[string("logstd")] = CALCLOGSTD; list[string("minmassweighted")] = CALCMINMASSWEIGHT; list[string("min")] = CALCMIN; list[string("maxmassweighted")] = CALCMAXMASSWEIGHT; list[string("max")] = CALCMAX; list[string("aperture_total")] = CALCQUANTITYAPERTURETOTAL; list[string("aperture_average")] = CALCQUANTITYAPERTUREAVERAGE; list[string("unkown")] = 0; return list; } inline map ExtraFieldCalculationsIntFlagToString() { map list; list[CALCAVERAGEMASSWEIGHT] = string("averagemassweighted");; list[CALCAVERAGE] = string("average"); list[CALCLOGAVERAGEMASSWEIGHT] = string("logaveragemassweighted"); list[CALCLOGAVERAGE] = string("logaverage"); list[CALCTOTALMASSWEIGHT] = string("totalmassweighted"); list[CALCTOTAL] = string("total"); list[CALCSTDMASSWEIGHT] = string("stdmassweighted"); list[CALCSTD] = string("std"); list[CALCLOGSTDMASSWEIGHT] = string("logstdmassweighted"); list[CALCLOGSTD] = string("logstd"); list[CALCMINMASSWEIGHT] = string("minmassweighted"); list[CALCMIN] = string("min"); list[CALCMAXMASSWEIGHT] = string("maxmassweighted"); list[CALCMAX] = string("max"); list[CALCQUANTITYAPERTURETOTAL] = string("aperture_total"); list[CALCQUANTITYAPERTUREAVERAGE] = string("aperture_average"); list[0] = string("unkown"); return list; } inline void ExtraFieldCheck(string configentryname, vector &names, vector&output_names, vector &calctypes, vector &indices, vector &conversions, vector &units, vector &pairindices, vector &names_aperture, vector &output_names_aperture, vector &calctypes_aperture, vector &indices_aperture, vector &conversions_aperture, vector &units_aperture ) { if (names.size()==0) return; set unique_int; set unique_name, outputset; string outputfieldname; unsigned int entryindex, calctype; string scalctype; vector stdfuncs, avefuncs; stdfuncs.push_back(CALCSTD);avefuncs.push_back(CALCAVERAGE); stdfuncs.push_back(CALCSTDMASSWEIGHT);avefuncs.push_back(CALCAVERAGEMASSWEIGHT); stdfuncs.push_back(CALCLOGSTD);avefuncs.push_back(CALCLOGAVERAGE); stdfuncs.push_back(CALCLOGSTDMASSWEIGHT);avefuncs.push_back(CALCLOGAVERAGEMASSWEIGHT); vector newnames(names), newunits(units); vector newcalctypes(calctypes); vector newindices(indices); vector newconversions(conversions); vector allowedlist = ExtraFieldCalculationsAllowedOptions(); string list; for (auto &x:allowedlist) list+= x+string(", "); map calcstringtoint = ExtraFieldCalculationsStringtoIntFlag(); map calcinttostring = ExtraFieldCalculationsIntFlagToString(); //check for unacceptable entries for (auto i=0;i0) { names.push_back(newnames[i]); indices.push_back(newindices[i]); units.push_back(newunits[i]); calctypes.push_back(newcalctypes[i]); conversions.push_back(newconversions[i]); } else { names_aperture.push_back(newnames[i]); indices_aperture.push_back(newindices[i]); units_aperture.push_back(newunits[i]); calctypes_aperture.push_back(newcalctypes[i]); conversions_aperture.push_back(newconversions[i]); } } pairindices.resize(names.size()); for (auto i=0;i &names, vector &values) { set unique; if (num ==0) return; int iduplicates = 0; for (auto &val:names) { if(unique.count(val)) { LOG_RANK0(warning) << "Duplicate entry(ies) found in " << configentryname; LOG_RANK0(warning) << "Removing duplicate entries"; iduplicates += 1; } else unique.insert(val); } if (iduplicates) { names = vector(unique.begin(),unique.end()); values.clear(); for (auto &val:names) values.push_back(stof(val)); num = values.size(); } } inline vector StripIndexOffName( vector names1, vector names2, vector indices1, vector indices2 ) { vector result; for (auto iextra=0;iextra calcstringtoint = ExtraFieldCalculationsStringtoIntFlag(); map calcinttostring = ExtraFieldCalculationsIntFlagToString(); if (!FileExists(opt.pname)){ if (ThisTask==0) LOG(error) << "Config file: " << opt.pname << " does not exist or can't be read, terminating"; #ifdef USEMPI MPI_Abort(MPI_COMM_WORLD,9); #else exit(9); #endif } paramfile.open(opt.pname, ios::in); std::string::size_type j; //first find output name, determine the number of valid entries in if (paramfile.is_open()) { while (paramfile.good()){ getline(paramfile,line); //if line is not commented out or empty if (line[0]!='#'&&line.length()!=0) { j = line.find(sep); if (j != std::string::npos) { //clean up string tag=line.substr(0,j); strcpy(buff, tag.c_str()); pbuff=strtok(buff," "); strcpy(tbuff, pbuff); val=line.substr(j+1,line.length()-(j+1)); strcpy(buff, val.c_str()); pbuff=strtok(buff," "); if (pbuff==NULL) continue; strcpy(vbuff, pbuff); if (strcmp(tbuff, "Output")==0){ opt.outname=new char[1024]; strcpy(opt.outname,vbuff); break; } } } } #ifndef SWIFTINTERFACE if (opt.outname==NULL) { if (ThisTask==0) cerr<<"No output name given, terminating"<0); else if (strcmp(tbuff, "MPI_zcurve_mesh_decomposition_min_num_cells_per_dim")==0) opt.minnumcellperdim = atoi(vbuff); ///OpenMP related else if (strcmp(tbuff, "OMP_run_fof")==0) opt.iopenmpfof = atoi(vbuff); else if (strcmp(tbuff, "OMP_fof_region_size")==0) opt.openmpfofsize = atoi(vbuff); else if (strcmp(tbuff, "Gas_internal_property_names")==0) { pos=0; dataline=string(vbuff); while ((pos = dataline.find(delimiter)) != string::npos) { token = dataline.substr(0, pos); opt.gas_internalprop_names.push_back(token); dataline.erase(0, pos + delimiter.length()); } if (opt.gas_internalprop_index.size() == 0) { opt.gas_internalprop_index.resize(opt.gas_internalprop_names.size(),0); } if (opt.gas_internalprop_function.size() == 0) { opt.gas_internalprop_function.resize(opt.gas_internalprop_names.size(),CALCAVERAGEMASSWEIGHT); } if (opt.gas_internalprop_input_output_unit_conversion_factors.size() == 0) { opt.gas_internalprop_input_output_unit_conversion_factors.resize( opt.gas_internalprop_names.size(),1.0); } if (opt.gas_internalprop_output_units.size() == 0) { opt.gas_internalprop_output_units.resize(opt.gas_internalprop_names.size(),string("unknownunits")); } } else if (strcmp(tbuff, "Gas_internal_property_index_in_file")==0) { pos=0; dataline=string(vbuff); vector tempvec; while ((pos = dataline.find(delimiter)) != string::npos) { token = dataline.substr(0, pos); tempvec.push_back(stoi(token)); dataline.erase(0, pos + delimiter.length()); } opt.gas_internalprop_index = tempvec; } else if (strcmp(tbuff, "Gas_internal_property_calculation_type")==0) { pos=0; dataline=string(vbuff); vector tempvec; while ((pos = dataline.find(delimiter)) != string::npos) { token = dataline.substr(0, pos); tempvec.push_back(calcstringtoint[token]); dataline.erase(0, pos + delimiter.length()); } opt.gas_internalprop_function = tempvec; } else if (strcmp(tbuff, "Gas_internal_property_output_units")==0) { pos=0; dataline=string(vbuff); vector tempvec; while ((pos = dataline.find(delimiter)) != string::npos) { token = dataline.substr(0, pos); tempvec.push_back(token); dataline.erase(0, pos + delimiter.length()); } opt.gas_internalprop_output_units = tempvec; } else if (strcmp(tbuff, "Gas_internal_property_input_output_unit_conversion_factors")==0) { pos=0; dataline=string(vbuff); vector tempvec; while ((pos = dataline.find(delimiter)) != string::npos) { token = dataline.substr(0, pos); tempvec.push_back(stof(token)); dataline.erase(0, pos + delimiter.length()); } opt.gas_internalprop_input_output_unit_conversion_factors = tempvec; } else if (strcmp(tbuff, "Tcut_halogas")==0) { opt.temp_max_cut = atof(vbuff); } else if (strcmp(tbuff, "Gas_temperature_input_unit_conversion_to_output_unit")==0) { opt.temp_input_output_unit_conversion_factor = atof(vbuff); } else if (strcmp(tbuff, "Gas_hot_overdensity_normalisation")==0) { opt.hot_gas_overdensity_normalisation = atof(vbuff); // if this is define, then make sure the SO calculations are check to contain this value, otherwise modify it to contain them. auto iter = std::find(opt.SOthresholds_values_crit.begin(), opt.SOthresholds_values_crit.end(), opt.hot_gas_overdensity_normalisation); if(iter == opt.SOthresholds_values_crit.end()){ opt.SOthresholds_values_crit.push_back(opt.hot_gas_overdensity_normalisation); opt.SOnum = opt.SOnum + 1; } else if (strcmp(tbuff, "Overdensity_values_in_critical_density")==0) { pos=0; dataline=string(vbuff); while ((pos = dataline.find(delimiter)) != string::npos) { token = dataline.substr(0, pos); opt.SOthresholds_names_crit.push_back(token); opt.SOthresholds_values_crit.push_back(stof(token)); dataline.erase(0, pos + delimiter.length()); } } } else if (strcmp(tbuff, "Aperture_values_normalised_to_overdensity")==0) { pos=0; dataline=string(vbuff); vector tempvec; while ((pos = dataline.find(delimiter)) != string::npos) { token = dataline.substr(0, pos); tempvec.push_back(stof(token)); dataline.erase(0, pos + delimiter.length()); } opt.aperture_hotgas_normalised_to_overdensity = tempvec; } else if (strcmp(tbuff, "Star_internal_property_names")==0) { pos=0; dataline=string(vbuff); while ((pos = dataline.find(delimiter)) != string::npos) { token = dataline.substr(0, pos); opt.star_internalprop_names.push_back(token); dataline.erase(0, pos + delimiter.length()); } if (opt.star_internalprop_index.size() == 0) { opt.star_internalprop_index.resize(opt.star_internalprop_names.size(),0); } if (opt.gas_internalprop_function.size() == 0) { opt.star_internalprop_function.resize(opt.star_internalprop_names.size(),CALCAVERAGEMASSWEIGHT); } if (opt.star_internalprop_input_output_unit_conversion_factors.size() == 0) { opt.star_internalprop_input_output_unit_conversion_factors.resize( opt.star_internalprop_names.size(),1.0); } if (opt.star_internalprop_output_units.size() == 0) { opt.star_internalprop_output_units.resize(opt.star_internalprop_names.size(),string("unknownunits")); } } else if (strcmp(tbuff, "Star_internal_property_index_in_file")==0) { pos=0; dataline=string(vbuff); vector tempvec; while ((pos = dataline.find(delimiter)) != string::npos) { token = dataline.substr(0, pos); tempvec.push_back(stoi(token)); dataline.erase(0, pos + delimiter.length()); } opt.star_internalprop_index = tempvec; } else if (strcmp(tbuff, "Star_internal_property_calculation_type")==0) { pos=0; dataline=string(vbuff); vector tempvec; while ((pos = dataline.find(delimiter)) != string::npos) { token = dataline.substr(0, pos); tempvec.push_back(calcstringtoint[token]); dataline.erase(0, pos + delimiter.length()); } opt.star_internalprop_function = tempvec; } else if (strcmp(tbuff, "Star_internal_property_output_units")==0) { pos=0; dataline=string(vbuff); vector tempvec; while ((pos = dataline.find(delimiter)) != string::npos) { token = dataline.substr(0, pos); tempvec.push_back(token); dataline.erase(0, pos + delimiter.length()); } opt.star_internalprop_output_units = tempvec; } else if (strcmp(tbuff, "Star_internal_property_input_output_unit_conversion_factors")==0) { pos=0; dataline=string(vbuff); vector tempvec; while ((pos = dataline.find(delimiter)) != string::npos) { token = dataline.substr(0, pos); tempvec.push_back(stof(token)); dataline.erase(0, pos + delimiter.length()); } opt.star_internalprop_input_output_unit_conversion_factors = tempvec; } else if (strcmp(tbuff, "BH_internal_property_names")==0) { pos=0; dataline=string(vbuff); while ((pos = dataline.find(delimiter)) != string::npos) { token = dataline.substr(0, pos); opt.bh_internalprop_names.push_back(token); dataline.erase(0, pos + delimiter.length()); } if (opt.bh_internalprop_index.size() == 0) { opt.bh_internalprop_index.resize(opt.bh_internalprop_names.size(),0); } if (opt.bh_internalprop_function.size() == 0) { opt.bh_internalprop_function.resize(opt.bh_internalprop_names.size(),CALCAVERAGEMASSWEIGHT); } if (opt.bh_internalprop_input_output_unit_conversion_factors.size() == 0) { opt.bh_internalprop_input_output_unit_conversion_factors.resize( opt.bh_internalprop_names.size(),1.0); } if (opt.bh_internalprop_output_units.size() == 0) { opt.bh_internalprop_output_units.resize(opt.bh_internalprop_names.size(),string("unknownunits")); } } else if (strcmp(tbuff, "BH_internal_property_index_in_file")==0) { pos=0; dataline=string(vbuff); vector tempvec; while ((pos = dataline.find(delimiter)) != string::npos) { token = dataline.substr(0, pos); tempvec.push_back(stoi(token)); dataline.erase(0, pos + delimiter.length()); } opt.bh_internalprop_index = tempvec; } else if (strcmp(tbuff, "BH_internal_property_calculation_type")==0) { pos=0; dataline=string(vbuff); vector tempvec; while ((pos = dataline.find(delimiter)) != string::npos) { token = dataline.substr(0, pos); tempvec.push_back(calcstringtoint[token]); dataline.erase(0, pos + delimiter.length()); } opt.bh_internalprop_function = tempvec; } else if (strcmp(tbuff, "BH_internal_property_output_units")==0) { pos=0; dataline=string(vbuff); vector tempvec; while ((pos = dataline.find(delimiter)) != string::npos) { token = dataline.substr(0, pos); tempvec.push_back(token); dataline.erase(0, pos + delimiter.length()); } opt.bh_internalprop_output_units = tempvec; } else if (strcmp(tbuff, "BH_internal_property_input_output_unit_conversion_factors")==0) { pos=0; dataline=string(vbuff); vector tempvec; while ((pos = dataline.find(delimiter)) != string::npos) { token = dataline.substr(0, pos); tempvec.push_back(stof(token)); dataline.erase(0, pos + delimiter.length()); } opt.bh_internalprop_input_output_unit_conversion_factors = tempvec; } else if (strcmp(tbuff, "Extra_DM_internal_property_names")==0) { pos=0; dataline=string(vbuff); while ((pos = dataline.find(delimiter)) != string::npos) { token = dataline.substr(0, pos); opt.extra_dm_internalprop_names.push_back(token); dataline.erase(0, pos + delimiter.length()); } if (opt.extra_dm_internalprop_index.size() == 0) { opt.extra_dm_internalprop_index.resize(opt.extra_dm_internalprop_names.size(),0); } if (opt.extra_dm_internalprop_function.size() == 0) { opt.extra_dm_internalprop_function.resize(opt.extra_dm_internalprop_names.size(),CALCAVERAGEMASSWEIGHT); } if (opt.extra_dm_internalprop_input_output_unit_conversion_factors.size() == 0) { opt.extra_dm_internalprop_input_output_unit_conversion_factors.resize( opt.extra_dm_internalprop_names.size(),1.0); } if (opt.extra_dm_internalprop_output_units.size() == 0) { opt.extra_dm_internalprop_output_units.resize(opt.extra_dm_internalprop_names.size(),string("unknownunits")); } } else if (strcmp(tbuff, "Extra_DM_internal_property_index_in_file")==0) { pos=0; dataline=string(vbuff); vector tempvec; while ((pos = dataline.find(delimiter)) != string::npos) { token = dataline.substr(0, pos); tempvec.push_back(stoi(token)); dataline.erase(0, pos + delimiter.length()); } opt.extra_dm_internalprop_index = tempvec; } else if (strcmp(tbuff, "Extra_DM_internal_property_calculation_type")==0) { pos=0; dataline=string(vbuff); vector tempvec; while ((pos = dataline.find(delimiter)) != string::npos) { token = dataline.substr(0, pos); tempvec.push_back(calcstringtoint[token]); dataline.erase(0, pos + delimiter.length()); } opt.extra_dm_internalprop_function = tempvec; } else if (strcmp(tbuff, "Extra_DM_internal_property_output_units")==0) { pos=0; dataline=string(vbuff); vector tempvec; while ((pos = dataline.find(delimiter)) != string::npos) { token = dataline.substr(0, pos); tempvec.push_back(token); dataline.erase(0, pos + delimiter.length()); } opt.extra_dm_internalprop_output_units = tempvec; } else if (strcmp(tbuff, "Extra_DM_internal_property_input_output_unit_conversion_factors")==0) { pos=0; dataline=string(vbuff); vector tempvec; while ((pos = dataline.find(delimiter)) != string::npos) { token = dataline.substr(0, pos); tempvec.push_back(stof(token)); dataline.erase(0, pos + delimiter.length()); } opt.extra_dm_internalprop_input_output_unit_conversion_factors = tempvec; } else if (strcmp(tbuff, "Gas_chemistry_names")==0) { pos=0; dataline=string(vbuff); while ((pos = dataline.find(delimiter)) != string::npos) { token = dataline.substr(0, pos); opt.gas_chem_names.push_back(token); dataline.erase(0, pos + delimiter.length()); } if (opt.gas_chem_index.size() == 0) { opt.gas_chem_index.resize(opt.gas_chem_names.size(),0); } if (opt.gas_chem_function.size() == 0) { opt.gas_chem_function.resize(opt.gas_chem_names.size(),CALCAVERAGEMASSWEIGHT); } if (opt.gas_chem_input_output_unit_conversion_factors.size() == 0) { opt.gas_chem_input_output_unit_conversion_factors.resize( opt.gas_chem_names.size(),1.0); } if (opt.gas_chem_output_units.size() == 0) { opt.gas_chem_output_units.resize(opt.gas_chem_names.size(),string("unknownunits")); } } else if (strcmp(tbuff, "Gas_chemistry_index_in_file")==0) { pos=0; dataline=string(vbuff); vector tempvec; while ((pos = dataline.find(delimiter)) != string::npos) { token = dataline.substr(0, pos); tempvec.push_back(stoi(token)); dataline.erase(0, pos + delimiter.length()); } opt.gas_chem_index = tempvec; } else if (strcmp(tbuff, "Gas_chemistry_calculation_type")==0) { pos=0; dataline=string(vbuff); vector tempvec; while ((pos = dataline.find(delimiter)) != string::npos) { token = dataline.substr(0, pos); tempvec.push_back(calcstringtoint[token]); dataline.erase(0, pos + delimiter.length()); } opt.gas_chem_function = tempvec; } else if (strcmp(tbuff, "Gas_chemistry_output_units")==0) { pos=0; dataline=string(vbuff); vector tempvec; while ((pos = dataline.find(delimiter)) != string::npos) { token = dataline.substr(0, pos); tempvec.push_back(token); dataline.erase(0, pos + delimiter.length()); } opt.gas_chem_output_units = tempvec; } else if (strcmp(tbuff, "Gas_chemistry_input_output_unit_conversion_factors")==0) { pos=0; dataline=string(vbuff); vector tempvec; while ((pos = dataline.find(delimiter)) != string::npos) { token = dataline.substr(0, pos); tempvec.push_back(stof(token)); dataline.erase(0, pos + delimiter.length()); } opt.gas_chem_input_output_unit_conversion_factors = tempvec; } else if (strcmp(tbuff, "Star_chemistry_names")==0) { pos=0; dataline=string(vbuff); while ((pos = dataline.find(delimiter)) != string::npos) { token = dataline.substr(0, pos); opt.star_chem_names.push_back(token); dataline.erase(0, pos + delimiter.length()); } if (opt.star_chem_index.size() == 0) { opt.star_chem_index.resize(opt.star_chem_names.size(),0); } if (opt.star_chem_function.size() == 0) { opt.star_chem_function.resize(opt.star_chem_names.size(),CALCAVERAGEMASSWEIGHT); } if (opt.star_chem_input_output_unit_conversion_factors.size() == 0) { opt.star_chem_input_output_unit_conversion_factors.resize( opt.star_chem_names.size(),1.0); } if (opt.star_chem_output_units.size() == 0) { opt.star_chem_output_units.resize(opt.star_chem_names.size(),string("unknownunits")); } } else if (strcmp(tbuff, "Star_chemistry_index_in_file")==0) { pos=0; dataline=string(vbuff); vector tempvec; while ((pos = dataline.find(delimiter)) != string::npos) { token = dataline.substr(0, pos); tempvec.push_back(stoi(token)); dataline.erase(0, pos + delimiter.length()); } opt.star_chem_index = tempvec; } else if (strcmp(tbuff, "Star_chemistry_calculation_type")==0) { pos=0; dataline=string(vbuff); vector tempvec; while ((pos = dataline.find(delimiter)) != string::npos) { token = dataline.substr(0, pos); tempvec.push_back(calcstringtoint[token]); dataline.erase(0, pos + delimiter.length()); } opt.star_chem_function = tempvec; } else if (strcmp(tbuff, "Star_chemistry_output_units")==0) { pos=0; dataline=string(vbuff); vector tempvec; while ((pos = dataline.find(delimiter)) != string::npos) { token = dataline.substr(0, pos); tempvec.push_back(token); dataline.erase(0, pos + delimiter.length()); } opt.star_chem_output_units = tempvec; } else if (strcmp(tbuff, "Star_chemistry_input_output_unit_conversion_factors")==0) { pos=0; dataline=string(vbuff); vector tempvec; while ((pos = dataline.find(delimiter)) != string::npos) { token = dataline.substr(0, pos); tempvec.push_back(stof(token)); dataline.erase(0, pos + delimiter.length()); } opt.star_chem_input_output_unit_conversion_factors = tempvec; } else if (strcmp(tbuff, "BH_chemistry_names")==0) { pos=0; dataline=string(vbuff); while ((pos = dataline.find(delimiter)) != string::npos) { token = dataline.substr(0, pos); opt.bh_chem_names.push_back(token); dataline.erase(0, pos + delimiter.length()); } if (opt.bh_chem_index.size() == 0) { opt.bh_chem_index.resize(opt.bh_chem_names.size(),0); } if (opt.bh_chem_function.size() == 0) { opt.bh_chem_function.resize(opt.bh_chem_names.size(),CALCAVERAGEMASSWEIGHT); } if (opt.bh_chem_input_output_unit_conversion_factors.size() == 0) { opt.bh_chem_input_output_unit_conversion_factors.resize( opt.bh_chem_names.size(),1.0); } if (opt.bh_chem_output_units.size() == 0) { opt.bh_chem_output_units.resize(opt.bh_chem_names.size(),string("unknownunits")); } } else if (strcmp(tbuff, "BH_chemistry_index_in_file")==0) { pos=0; dataline=string(vbuff); vector tempvec; while ((pos = dataline.find(delimiter)) != string::npos) { token = dataline.substr(0, pos); tempvec.push_back(stoi(token)); dataline.erase(0, pos + delimiter.length()); } opt.bh_chem_index = tempvec; } else if (strcmp(tbuff, "BH_chemistry_calculation_type")==0) { pos=0; dataline=string(vbuff); vector tempvec; while ((pos = dataline.find(delimiter)) != string::npos) { token = dataline.substr(0, pos); tempvec.push_back(calcstringtoint[token]); dataline.erase(0, pos + delimiter.length()); } opt.bh_chem_function = tempvec; } else if (strcmp(tbuff, "BH_chemistry_output_units")==0) { pos=0; dataline=string(vbuff); vector tempvec; while ((pos = dataline.find(delimiter)) != string::npos) { token = dataline.substr(0, pos); tempvec.push_back(token); dataline.erase(0, pos + delimiter.length()); } opt.bh_chem_output_units = tempvec; } else if (strcmp(tbuff, "BH_chemistry_input_output_unit_conversion_factors")==0) { pos=0; dataline=string(vbuff); vector tempvec; while ((pos = dataline.find(delimiter)) != string::npos) { token = dataline.substr(0, pos); tempvec.push_back(stof(token)); dataline.erase(0, pos + delimiter.length()); } opt.bh_chem_input_output_unit_conversion_factors = tempvec; } else if (strcmp(tbuff, "Gas_chemistry_production_names")==0) { pos=0; dataline=string(vbuff); while ((pos = dataline.find(delimiter)) != string::npos) { token = dataline.substr(0, pos); opt.gas_chemproduction_names.push_back(token); dataline.erase(0, pos + delimiter.length()); } if (opt.gas_chemproduction_index.size() == 0) { opt.gas_chemproduction_index.resize(opt.gas_chemproduction_names.size(),0); } if (opt.gas_chemproduction_function.size() == 0) { opt.gas_chemproduction_function.resize(opt.gas_chemproduction_names.size(),CALCAVERAGEMASSWEIGHT); } if (opt.gas_chemproduction_input_output_unit_conversion_factors.size() == 0) { opt.gas_chemproduction_input_output_unit_conversion_factors.resize( opt.gas_chemproduction_names.size(),1.0); } if (opt.gas_chemproduction_output_units.size() == 0) { opt.gas_chemproduction_output_units.resize(opt.gas_chemproduction_names.size(),string("unknownunits")); } } else if (strcmp(tbuff, "Gas_chemistry_production_index_in_file")==0) { pos=0; dataline=string(vbuff); vector tempvec; while ((pos = dataline.find(delimiter)) != string::npos) { token = dataline.substr(0, pos); tempvec.push_back(stoi(token)); dataline.erase(0, pos + delimiter.length()); } opt.gas_chemproduction_index = tempvec; } else if (strcmp(tbuff, "Gas_chemistry_prodution_calculation_type")==0) { pos=0; dataline=string(vbuff); vector tempvec; while ((pos = dataline.find(delimiter)) != string::npos) { token = dataline.substr(0, pos); tempvec.push_back(calcstringtoint[token]); dataline.erase(0, pos + delimiter.length()); } opt.gas_chemproduction_function = tempvec; } else if (strcmp(tbuff, "Gas_chemistry_production_output_units")==0) { pos=0; dataline=string(vbuff); vector tempvec; while ((pos = dataline.find(delimiter)) != string::npos) { token = dataline.substr(0, pos); tempvec.push_back(token); dataline.erase(0, pos + delimiter.length()); } opt.gas_chemproduction_output_units = tempvec; } else if (strcmp(tbuff, "Gas_chemistry_production_input_output_unit_conversion_factors")==0) { pos=0; dataline=string(vbuff); vector tempvec; while ((pos = dataline.find(delimiter)) != string::npos) { token = dataline.substr(0, pos); tempvec.push_back(stof(token)); dataline.erase(0, pos + delimiter.length()); } opt.gas_chemproduction_input_output_unit_conversion_factors = tempvec; } else if (strcmp(tbuff, "Star_chemistry_production_names")==0) { pos=0; dataline=string(vbuff); while ((pos = dataline.find(delimiter)) != string::npos) { token = dataline.substr(0, pos); opt.star_chemproduction_names.push_back(token); dataline.erase(0, pos + delimiter.length()); } if (opt.star_chemproduction_index.size() == 0) { opt.star_chemproduction_index.resize(opt.star_chemproduction_names.size(),0); } if (opt.star_chemproduction_function.size() == 0) { opt.star_chemproduction_function.resize(opt.star_chemproduction_names.size(),CALCAVERAGEMASSWEIGHT); } if (opt.star_chemproduction_input_output_unit_conversion_factors.size() == 0) { opt.star_chemproduction_input_output_unit_conversion_factors.resize( opt.star_chemproduction_names.size(),1.0); } if (opt.star_chemproduction_output_units.size() == 0) { opt.star_chemproduction_output_units.resize(opt.star_chemproduction_names.size(),string("unknownunits")); } } else if (strcmp(tbuff, "Star_chemistry_production_index_in_file")==0) { pos=0; dataline=string(vbuff); vector tempvec; while ((pos = dataline.find(delimiter)) != string::npos) { token = dataline.substr(0, pos); tempvec.push_back(stoi(token)); dataline.erase(0, pos + delimiter.length()); } opt.star_chemproduction_index = tempvec; } else if (strcmp(tbuff, "Star_chemistry_prodution_calculation_type")==0) { pos=0; dataline=string(vbuff); vector tempvec; while ((pos = dataline.find(delimiter)) != string::npos) { token = dataline.substr(0, pos); tempvec.push_back(calcstringtoint[token]); dataline.erase(0, pos + delimiter.length()); } opt.star_chemproduction_function = tempvec; } else if (strcmp(tbuff, "Star_chemistry_production_output_units")==0) { pos=0; dataline=string(vbuff); vector tempvec; while ((pos = dataline.find(delimiter)) != string::npos) { token = dataline.substr(0, pos); tempvec.push_back(token); dataline.erase(0, pos + delimiter.length()); } opt.star_chemproduction_output_units = tempvec; } else if (strcmp(tbuff, "Star_chemistry_production_input_output_unit_conversion_factors")==0) { pos=0; dataline=string(vbuff); vector tempvec; while ((pos = dataline.find(delimiter)) != string::npos) { token = dataline.substr(0, pos); tempvec.push_back(stof(token)); dataline.erase(0, pos + delimiter.length()); } opt.star_chemproduction_input_output_unit_conversion_factors = tempvec; } else if (strcmp(tbuff, "BH_chemistry_production_names")==0) { pos=0; dataline=string(vbuff); while ((pos = dataline.find(delimiter)) != string::npos) { token = dataline.substr(0, pos); opt.bh_chemproduction_names.push_back(token); dataline.erase(0, pos + delimiter.length()); } if (opt.bh_chemproduction_index.size() == 0) { opt.bh_chemproduction_index.resize(opt.bh_chemproduction_names.size(),0); } if (opt.bh_chemproduction_function.size() == 0) { opt.bh_chemproduction_function.resize(opt.bh_chemproduction_names.size(),CALCAVERAGEMASSWEIGHT); } if (opt.bh_chemproduction_input_output_unit_conversion_factors.size() == 0) { opt.bh_chemproduction_input_output_unit_conversion_factors.resize( opt.bh_chemproduction_names.size(),1.0); } if (opt.bh_chemproduction_output_units.size() == 0) { opt.bh_chemproduction_output_units.resize(opt.bh_chemproduction_names.size(),string("unknownunits")); } } else if (strcmp(tbuff, "BH_chemistry_production_index_in_file")==0) { pos=0; dataline=string(vbuff); vector tempvec; while ((pos = dataline.find(delimiter)) != string::npos) { token = dataline.substr(0, pos); tempvec.push_back(stoi(token)); dataline.erase(0, pos + delimiter.length()); } opt.bh_chemproduction_index = tempvec; } else if (strcmp(tbuff, "BH_chemistry_prodution_calculation_type")==0) { pos=0; dataline=string(vbuff); vector tempvec; while ((pos = dataline.find(delimiter)) != string::npos) { token = dataline.substr(0, pos); tempvec.push_back(calcstringtoint[token]); dataline.erase(0, pos + delimiter.length()); } opt.bh_chemproduction_function = tempvec; } else if (strcmp(tbuff, "BH_chemistry_production_output_units")==0) { pos=0; dataline=string(vbuff); vector tempvec; while ((pos = dataline.find(delimiter)) != string::npos) { token = dataline.substr(0, pos); tempvec.push_back(token); dataline.erase(0, pos + delimiter.length()); } opt.bh_chemproduction_output_units = tempvec; } else if (strcmp(tbuff, "BH_chemistry_production_input_output_unit_conversion_factors")==0) { pos=0; dataline=string(vbuff); vector tempvec; while ((pos = dataline.find(delimiter)) != string::npos) { token = dataline.substr(0, pos); tempvec.push_back(stof(token)); dataline.erase(0, pos + delimiter.length()); } opt.bh_chemproduction_input_output_unit_conversion_factors = tempvec; } //output related else if (strcmp(tbuff, "Separate_output_files")==0) opt.iseparatefiles = atoi(vbuff); else if (strcmp(tbuff, "Binary_output")==0) opt.ibinaryout = atoi(vbuff); else if (strcmp(tbuff, "Comoving_units")==0) opt.icomoveunit = atoi(vbuff); else if (strcmp(tbuff, "Extended_output")==0) opt.iextendedoutput = atoi(vbuff); else if (strcmp(tbuff, "Spherical_overdensity_halo_particle_list_output")==0) opt.iSphericalOverdensityPartList = atoi(vbuff); else if (strcmp(tbuff, "Sort_by_binding_energy")==0) opt.iSortByBindingEnergy = atoi(vbuff); else if (strcmp(tbuff, "SUBFIND_like_output")==0) opt.isubfindoutput = atoi(vbuff); else if (strcmp(tbuff, "No_particle_ID_list_output")==0) opt.inoidoutput = atoi(vbuff); //gadget io related to extra info for sph, stars, bhs, else if (strcmp(tbuff, "NSPH_extra_blocks")==0) opt.gnsphblocks = atoi(vbuff); else if (strcmp(tbuff, "NStar_extra_blocks")==0) opt.gnstarblocks = atoi(vbuff); else if (strcmp(tbuff, "NBH_extra_blocks")==0) opt.gnbhblocks = atoi(vbuff); //input related to info for stars, bhs, winds/tracers, etc else if (strcmp(tbuff, "HDF_name_convention")==0) opt.ihdfnameconvention = atoi(vbuff); else if (strcmp(tbuff, "Input_includes_dm_particle")==0) opt.iusedmparticles = atoi(vbuff); else if (strcmp(tbuff, "Input_includes_gas_particle")==0) opt.iusegasparticles = atoi(vbuff); else if (strcmp(tbuff, "Input_includes_star_particle")==0) opt.iusestarparticles = atoi(vbuff); else if (strcmp(tbuff, "Input_includes_bh_particle")==0) opt.iusesinkparticles = atoi(vbuff); else if (strcmp(tbuff, "Input_includes_wind_particle")==0) opt.iusewindparticles = atoi(vbuff); else if (strcmp(tbuff, "Input_includes_tracer_particle")==0) opt.iusetracerparticles = atoi(vbuff); else if (strcmp(tbuff, "Input_includes_extradm_particle")==0) opt.iuseextradarkparticles = atoi(vbuff); } } } paramfile.close(); //cfgfile.close(); } configure_logging(opt); } void NOMASSCheck(Options &opt) { #ifdef NOMASS if (opt.MassValue<=0) { LOG_RANK0(error) << "Code compiled to not store mass per particle"; LOG_RANK0(error) << "Valid value of particle mass must be extracted from input (HDF5, Gadget) or set in the config file through Mass_value"; LOG_RANK0(error) << "Currently value <=0. Either input did not contain a valid value and update to config file needed"; ConfigExit(); } #endif } void ConfigCheck(Options &opt) { #ifndef USEMPI int ThisTask =0; #endif //if code is on the fly, no point in checking fname // input type, number of input files, etc if (opt.iontheflyfinding == false) { if (opt.fname==NULL||opt.outname==NULL){ ConfigExit("Must provide input and output file names"); } if (opt.inputtype == IOHDF && opt.ihdfnameconvention == -1) { ConfigExit("HDF input passed but the naming convention is not set in the config file. Please set HDF_name_convetion"); } if (opt.num_files<1){ ConfigExit("Invalid number of input files (<1)"); } if (opt.inputbufsize<1){ ConfigExit("Invalid read buf size (<1)"); } } if (opt.iBaryonSearch && !(opt.partsearchtype==PSTALL || opt.partsearchtype==PSTDARK)) { ConfigExit("Conflict in config file: both gas/star/etc particle type search AND the separate baryonic (gas,star,etc) search flag are on. Check config"); } if (opt.iBoundHalos && opt.iKeepFOF) { ConfigExit("Conflict in config file: Asking for Bound Field objects but also asking to keep the 3DFOF/then run 6DFOF. This is incompatible. Check config"); } if (opt.HaloMinSize==-1) opt.HaloMinSize=opt.MinSize; if (opt.lengthtokpc<=0){ ConfigExit("Invalid unit conversion, length unit to kpc is <=0 or was not set. Update config file"); } //convert reference apertures opt.lengthtokpc30pow2 /= opt.lengthtokpc*opt.lengthtokpc; opt.lengthtokpc50pow2 /= opt.lengthtokpc*opt.lengthtokpc; if (opt.velocitytokms<=0){ ConfigExit("Invalid unit conversion, velocity unit to km/s is <=0 or was not set. Update config file"); } if (opt.masstosolarmass<=0){ ConfigExit("Invalid unit conversion, mass unit to solar mass is <=0 or was not set. Update config file"); } #if defined(GASON) || defined(STARON) || defined(BHON) if (opt.SFRtosolarmassperyear<=0){ ConfigExit("Invalid unit conversion, SFR unit to solar mass per year is <=0 or was not set. Update config file"); } if (opt.metallicitytosolar<=0){ ConfigExit("Invalid unit conversion, metallicity unit to solar is <=0 or was not set. Update config file"); } if (opt.stellaragetoyrs<=0){ ConfigExit("Invalid unit conversion, stellar age unit to year is <=0 or was not set. Update config file"); } #endif #ifndef NOMASS opt.MassValue = 1.0; #endif //check gravity and hubble unit double gravity, hubunit; gravity = CalcGravitationalConstant(opt); hubunit = CalcHubbleUnit(opt); if (opt.G<=0) opt.G = gravity; else { auto diff = fabs(gravity-opt.G)/opt.G; if (diff>1e-2) { LOG_RANK0(warning) << "Configuration provides gravitational constant that differs from the default by more than 1%"; LOG_RANK0(warning) << "Expecation: " << gravity; LOG_RANK0(warning) << "Value passed: " << opt.G; } } if (opt.H<=0) opt.H = hubunit; else { auto diff = fabs(hubunit-opt.H)/opt.H; if (diff>1e-2) { LOG_RANK0(warning) << "Configuration provides hubble units that differs from the default by more than 1%"; LOG_RANK0(warning) << "Expecation: " << hubunit; LOG_RANK0(warning) << "Value passed: " << opt.H; } } #ifdef USEMPI if (opt.minnumcellperdim<8){ LOG_RANK0(warning) << "MPI mesh too coarse, minimum number of cells per dimension from which to produce z-curve decomposition is 8. Resetting to 8"; opt.minnumcellperdim = 8; } if (opt.mpiparticletotbufsize<(long int)(sizeof(Particle)*NProcs) && opt.mpiparticletotbufsize!=-1){ LOG_RANK0(error) << "Invalid input particle buffer send size, minimum input buffer size given particle byte size " << vr::memory_amount(sizeof(Particle)) << " and have " << NProcs << " MPI processes is " << vr::memory_amount(sizeof(Particle) * NProcs); ConfigExit(); } //if total buffer size is -1 then calculate individual buffer size based on default mpi size else if (opt.mpiparticletotbufsize==-1) { opt.mpiparticletotbufsize=MPIPartBufSize*NProcs*sizeof(Particle); opt.mpiparticlebufsize=MPIPartBufSize; } else { opt.mpiparticlebufsize=opt.mpiparticletotbufsize/NProcs/sizeof(Particle); } if (opt.mpipartfac<0){ ConfigExit("Invalid MPI particle allocation factor, must be >0"); } else if (opt.mpipartfac>1){ LOG_RANK0(warning) << "MPI Particle allocation factor is high (>1)"; } if (opt.mpinprocswritesize<1){ #ifdef USEPARALLELHDF LOG_RANK0(warning) << "Number of MPI task writing collectively < 1. Setting to 1"; opt.mpinprocswritesize = 1; #endif } if (opt.mpinprocswritesize>NProcs){ #ifdef USEPARALLELHDF LOG_RANK0(warning) << "Number of MPI task writing collectively > NProcs. Setting to NProcs"; opt.mpinprocswritesize = NProcs; #endif } #endif #ifdef USEOPENMP if (opt.iopenmpfof == 1 && opt.openmpfofsize < ompfofsearchnum){ LOG_RANK0(warning) << "OpenMP FOF search region is small, resetting to minimum of " << ompfofsearchnum; opt.openmpfofsize = ompfofsearchnum; } #endif #ifndef USEHDF if (opt.ibinaryout==OUTHDF){ ConfigExit("Code not compiled with HDF output enabled. Recompile with this enabled or change Binary_output"); } if (opt.isubfindoutput) { ConfigExit("Code not compiled with HDF output enabled. Recompile with this enabled to produce subfind like output or turn off subfind like output"); } #endif #ifndef USEADIOS if (opt.ibinaryout==OUTADIOS){ ConfigExit("Code not compiled with ADIOS output enabled. Recompile with this enabled or change Binary_output"); } #endif if (opt.iaperturecalc>0) { if (opt.aperturenum != opt.aperture_values_kpc.size()) { ConfigExit("Aperture calculations requested but mismatch between number stated and values provided. Check config"); } if (opt.aperturenum == 0) { LOG_RANK0(warning) << "Aperture calculations requested but number of apertures is zero"; } if (opt.apertureprojnum != opt.aperture_proj_values_kpc.size()) { ConfigExit("Projected aperture calculations requested but mismatch between number stated and values provided. Check config"); } if (opt.apertureprojnum == 0) { LOG_RANK0(warning) << "Aperture calculations requested but number of projected apertures is zero"; } if (opt.aperturenum>0) for (auto &x:opt.aperture_values_kpc) x/=opt.lengthtokpc; if (opt.apertureprojnum>0) for (auto &x:opt.aperture_proj_values_kpc) x/=opt.lengthtokpc; } if (opt.SOnum>0) { if (opt.SOnum != opt.SOthresholds_values_crit.size()) { ConfigExit("SO calculations requested but mismatch between number stated and values provided. Check config"); } } if (opt.iprofilecalc>0) { if (opt.profileminsize < 0) { LOG_RANK0(warning) << "Radial profile calculations limited to objects of < 0 size! Ignoring and setting to 0"; opt.profileminsize = 0; } if (opt.profileminFOFsize < 0) { LOG_RANK0(warning) << "Radial profile calculations limited to FOF objects of < 0 size! Ignoring and setting to 0"; opt.profileminFOFsize = 0; } if (opt.profilenbins != opt.profile_bin_edges.size()) { ConfigExit("Radial profile calculations requested but mismatch between number of edges stated and number provided. Check config"); } if (opt.profilenbins == 0) { ConfigExit("Radial profile calculations requested but number of bin edges is zero. Check config"); } if (opt.iprofilebintype == PROFILERBINTYPELOG) { for (auto i=0;i POTAPPROXMETHODRAND) { ConfigExit("In approximate potential but using invalid method for sampling particles. Use 0 for Tree and 1 for Rand. Check config."); } } set uniqueval; set outputset; string configentryname, outputfieldname, mainname; unsigned int entryindex, calctype, iduplicates; //clean up aperture list and spherical overdensity list to remove duplicates configentryname = "Overdensity_values_in_critical_density"; ListDuplicateEntryCheck(configentryname, opt.SOnum, opt.SOthresholds_names_crit, opt.SOthresholds_values_crit); configentryname = "Aperture_values_in_kpc"; ListDuplicateEntryCheck(configentryname, opt.aperturenum, opt.aperture_names_kpc, opt.aperture_values_kpc); configentryname = "Projected_aperture_values_in_kpc"; ListDuplicateEntryCheck(configentryname, opt.apertureprojnum, opt.aperture_proj_names_kpc, opt.aperture_proj_values_kpc); //set output field names for extra properties configentryname = "Gas_internal_property_names"; ExtraFieldCheck(configentryname, opt.gas_internalprop_names, opt.gas_internalprop_output_names, opt.gas_internalprop_function, opt.gas_internalprop_index, opt.gas_internalprop_input_output_unit_conversion_factors, opt.gas_internalprop_output_units, opt.gas_internalprop_index_paired_calc, opt.gas_internalprop_names_aperture, opt.gas_internalprop_output_names_aperture, opt.gas_internalprop_function_aperture, opt.gas_internalprop_index_aperture, opt.gas_internalprop_input_output_unit_conversion_factors_aperture, opt.gas_internalprop_output_units_aperture ); configentryname = "Gas_chemistry_names"; ExtraFieldCheck(configentryname, opt.gas_chem_names, opt.gas_chem_output_names, opt.gas_chem_function, opt.gas_chem_index, opt.gas_chem_input_output_unit_conversion_factors, opt.gas_chem_output_units, opt.gas_chem_index_paired_calc, opt.gas_chem_names_aperture, opt.gas_chem_output_names_aperture, opt.gas_chem_function_aperture, opt.gas_chem_index_aperture, opt.gas_chem_input_output_unit_conversion_factors_aperture, opt.gas_chem_output_units_aperture ); configentryname = "Gas_chemistry_production_names"; ExtraFieldCheck(configentryname, opt.gas_chemproduction_names, opt.gas_chemproduction_output_names, opt.gas_chemproduction_function, opt.gas_chemproduction_index, opt.gas_chemproduction_input_output_unit_conversion_factors, opt.gas_chemproduction_output_units, opt.gas_chemproduction_index_paired_calc, opt.gas_chemproduction_names_aperture, opt.gas_chemproduction_output_names_aperture, opt.gas_chemproduction_function_aperture, opt.gas_chemproduction_index_aperture, opt.gas_chemproduction_input_output_unit_conversion_factors_aperture, opt.gas_chemproduction_output_units_aperture ); opt.gas_extraprop_aperture_calc = (opt.gas_internalprop_names_aperture.size() + opt.gas_chem_names_aperture.size()+opt.gas_chemproduction_names_aperture.size()>0); configentryname = "Star_internal_property_names"; ExtraFieldCheck(configentryname, opt.star_internalprop_names, opt.star_internalprop_output_names, opt.star_internalprop_function, opt.star_internalprop_index, opt.star_internalprop_input_output_unit_conversion_factors, opt.star_internalprop_output_units, opt.star_internalprop_index_paired_calc, opt.star_internalprop_names_aperture, opt.star_internalprop_output_names_aperture, opt.star_internalprop_function_aperture, opt.star_internalprop_index_aperture, opt.star_internalprop_input_output_unit_conversion_factors_aperture, opt.star_internalprop_output_units_aperture ); configentryname = "Star_chemistry_names"; ExtraFieldCheck(configentryname, opt.star_chem_names, opt.star_chem_output_names, opt.star_chem_function, opt.star_chem_index, opt.star_chem_input_output_unit_conversion_factors, opt.star_chem_output_units, opt.star_chem_index_paired_calc, opt.star_chem_names_aperture, opt.star_chem_output_names_aperture, opt.star_chem_function_aperture, opt.star_chem_index_aperture, opt.star_chem_input_output_unit_conversion_factors_aperture, opt.star_chem_output_units_aperture ); configentryname = "Star_chemistry_production_names"; ExtraFieldCheck(configentryname, opt.star_chemproduction_names, opt.star_chemproduction_output_names, opt.star_chemproduction_function, opt.star_chemproduction_index, opt.star_chemproduction_input_output_unit_conversion_factors, opt.star_chemproduction_output_units, opt.star_chemproduction_index_paired_calc, opt.star_chemproduction_names_aperture, opt.star_chemproduction_output_names_aperture, opt.star_chemproduction_function_aperture, opt.star_chemproduction_index_aperture, opt.star_chemproduction_input_output_unit_conversion_factors_aperture, opt.star_chemproduction_output_units_aperture ); opt.star_extraprop_aperture_calc = (opt.star_internalprop_names_aperture.size() + opt.star_chem_names_aperture.size()+opt.star_chemproduction_names_aperture.size()>0); configentryname = "BH_internal_property_names"; ExtraFieldCheck(configentryname, opt.bh_internalprop_names, opt.bh_internalprop_output_names, opt.bh_internalprop_function, opt.bh_internalprop_index, opt.bh_internalprop_input_output_unit_conversion_factors, opt.bh_internalprop_output_units, opt.bh_internalprop_index_paired_calc, opt.bh_internalprop_names_aperture, opt.bh_internalprop_output_names_aperture, opt.bh_internalprop_function_aperture, opt.bh_internalprop_index_aperture, opt.bh_internalprop_input_output_unit_conversion_factors_aperture, opt.bh_internalprop_output_units_aperture ); configentryname = "BH_chemistry_names"; ExtraFieldCheck(configentryname, opt.bh_chem_names, opt.bh_chem_output_names, opt.bh_chem_function, opt.bh_chem_index, opt.bh_chem_input_output_unit_conversion_factors, opt.bh_chem_output_units, opt.bh_chem_index_paired_calc, opt.bh_chem_names_aperture, opt.bh_chem_output_names_aperture, opt.bh_chem_function_aperture, opt.bh_chem_index_aperture, opt.bh_chem_input_output_unit_conversion_factors_aperture, opt.bh_chem_output_units_aperture ); configentryname = "BH_chemistry_production_names"; ExtraFieldCheck(configentryname, opt.bh_chemproduction_names, opt.bh_chemproduction_output_names, opt.bh_chemproduction_function, opt.bh_chemproduction_index, opt.bh_chemproduction_input_output_unit_conversion_factors, opt.bh_chemproduction_output_units, opt.bh_chemproduction_index_paired_calc, opt.bh_chemproduction_names_aperture, opt.bh_chemproduction_output_names_aperture, opt.bh_chemproduction_function_aperture, opt.bh_chemproduction_index_aperture, opt.bh_chemproduction_input_output_unit_conversion_factors_aperture, opt.bh_chemproduction_output_units_aperture ); opt.bh_extraprop_aperture_calc = (opt.bh_internalprop_names_aperture.size() + opt.bh_chem_names_aperture.size()+opt.bh_chemproduction_names_aperture.size()>0); configentryname = "Extra_DM_internal_property_names"; ExtraFieldCheck(configentryname, opt.extra_dm_internalprop_names, opt.extra_dm_internalprop_output_names, opt.extra_dm_internalprop_function, opt.extra_dm_internalprop_index, opt.extra_dm_internalprop_input_output_unit_conversion_factors, opt.extra_dm_internalprop_output_units, opt.extra_dm_internalprop_index_paired_calc, opt.extra_dm_internalprop_names_aperture, opt.extra_dm_internalprop_output_names_aperture, opt.extra_dm_internalprop_function_aperture, opt.extra_dm_internalprop_index_aperture, opt.extra_dm_internalprop_input_output_unit_conversion_factors_aperture, opt.extra_dm_internalprop_output_units_aperture ); opt.extra_dm_extraprop_aperture_calc = (opt.extra_dm_internalprop_names_aperture.size() >0); if (opt.gas_extraprop_aperture_calc && opt.iaperturecalc == 0){ LOG_RANK0(error) << "Requesting extra gas properties to be calculated in apertures but apertures not set"; LOG_RANK0(error) << "Enable aperture calculations and provide apertures or change the calculations requested for extra properties. Check config"; ConfigExit(); } if (opt.star_extraprop_aperture_calc && opt.iaperturecalc == 0){ LOG_RANK0(error) << "Requesting extra star properties to be calculated in apertures but apertures not set"; LOG_RANK0(error) << "Enable aperture calculations and provide apertures or change the calculations requested for extra properties. Check config"; ConfigExit(); } if (opt.bh_extraprop_aperture_calc && opt.iaperturecalc == 0){ LOG_RANK0(error) << "Requesting extra BH properties to be calculated in apertures but apertures not set"; LOG_RANK0(error) << "Enable aperture calculations and provide apertures or change the calculations requested for extra properties. Check config"; ConfigExit(); } if (opt.extra_dm_extraprop_aperture_calc && opt.iaperturecalc == 0){ LOG_RANK0(error) << "Requesting extra dm properties to be calculated in apertures but apertures not set"; LOG_RANK0(error) << "Enable aperture calculations and provide apertures or change the calculations requested for extra properties. Check config"; ConfigExit(); } //set halo 3d fof linking length if necessary if (opt.ellhalo3dxfac == -1) { opt.ellhalo3dxfac = opt.ellhalophysfac * opt.ellphys; } else { opt.ellhalophysfac = opt.ellhalo3dxfac / opt.ellphys; } if (ThisTask==0) { log_options_summary(opt); } //store the git hash opt.git_sha1 = velociraptor::git_sha1(); #ifdef USEMPI MPI_Barrier(MPI_COMM_WORLD); #endif } ///define the configuration constructor so all input options are stored. ConfigInfo::ConfigInfo(Options &opt){ string datastring; //if compiler is super old and does not have at least std 11 implementation to_string does not exist #ifndef OLDCCOMPILER //Add version AddEntry("Git_revision", opt.git_sha1); //general search operations AddEntry("Particle_search_type", opt.partsearchtype); AddEntry("FoF_search_type", opt.foftype); AddEntry("FoF_Field_search_type", opt.fofbgtype); AddEntry("Search_for_substructure", opt.iSubSearch); AddEntry("Keep_FOF", opt.iKeepFOF); AddEntry("Iterative_searchflag", opt.iiterflag); AddEntry("Baryon_searchflag", opt.iBaryonSearch); AddEntry("CMrefadjustsubsearch_flag", opt.icmrefadjust); AddEntry("Halo_core_search", opt.iHaloCoreSearch); AddEntry("Use_adaptive_core_search", opt.iAdaptiveCoreLinking); AddEntry("Use_phase_tensor_core_growth", opt.iPhaseCoreGrowth); //local field parameters AddEntry("Local_velocity_density_approximate_calculation", opt.iLocalVelDenApproxCalcFlag); AddEntry("Cell_fraction", opt.Ncellfac); AddEntry("Grid_type", opt.gridtype); AddEntry("Nsearch_velocity", opt.Nvel); AddEntry("Nsearch_physical", opt.Nsearch); //substructure search parameters AddEntry("Outlier_threshold", opt.ellthreshold); AddEntry("Significance_level", opt.siglevel); AddEntry("Velocity_ratio", opt.Vratio); AddEntry("Velocity_opening_angle", opt.thetaopen); ///\todo this configuration option will be deprecated. Replaced by Substructure_physical_linking_length //AddEntry("Physical_linking_length", opt.ellphys); AddEntry("Substructure_physical_linking_length", opt.ellphys); AddEntry("Velocity_linking_length", opt.ellvel); AddEntry("Minimum_size", opt.MinSize); //field object specific searches AddEntry("Minimum_halo_size", opt.HaloMinSize); ///\todo this configuration option will be deprecated. Replaced by Halo_3D_physical_linking_length //AddEntry("Halo_linking_length_factor", opt.ellhalophysfac); AddEntry("Halo_3D_linking_length", opt.ellhalo3dxfac); AddEntry("Halo_velocity_linking_length_factor", opt.ellhalovelfac); //specific to 6DFOF field search AddEntry("Halo_6D_linking_length_factor", opt.ellhalo6dxfac); AddEntry("Halo_6D_vel_linking_length_factor", opt.ellhalo6dvfac); //specific search for 6d fof core searches AddEntry("Halo_core_ellx_fac", opt.halocorexfac); AddEntry("Halo_core_ellv_fac", opt.halocorevfac); AddEntry("Halo_core_ncellfac", opt.halocorenfac); AddEntry("Halo_core_adaptive_sigma_fac", opt.halocoresigmafac); AddEntry("Halo_core_num_loops", opt.halocorenumloops); AddEntry("Halo_core_loop_ellx_fac", opt.halocorexfaciter); AddEntry("Halo_core_loop_ellv_fac", opt.halocorevfaciter); AddEntry("Halo_core_loop_elln_fac", opt.halocorenumfaciter); AddEntry("Halo_core_phase_significance", opt.halocorephasedistsig); //for merging structures together AddEntry("Structure_phase_merge_dist", opt.coresubmergemindist); AddEntry("Apply_phase_merge_to_host", opt.icoresubmergewithbg); //for changing factors used in iterative search AddEntry("Iterative_threshold_factor", opt.ellfac); AddEntry("Iterative_linking_length_factor", opt.ellxfac); AddEntry("Iterative_Vratio_factor", opt.vfac); AddEntry("Iterative_ThetaOp_factor", opt.thetafac); //for changing effective resolution when rescaling linking lengh #ifdef HIGHRES AddEntry("Effective_resolution", opt.Neff); #endif //for changing effective resolution when rescaling linking lengh AddEntry("Singlehalo_search", opt.iSingleHalo); //units, cosmology AddEntry("Length_unit", opt.lengthinputconversion); AddEntry("Velocity_unit", opt.velocityinputconversion); AddEntry("Mass_unit", opt.massinputconversion); AddEntry("Length_input_unit_conversion_to_output_unit", opt.lengthinputconversion); AddEntry("Velocity_input_unit_conversion_to_output_unit", opt.velocityinputconversion); AddEntry("Mass_input_unit_conversion_to_output_unit", opt.massinputconversion); AddEntry("Star_formation_rate_input_unit_conversion_to_output_unit", opt.SFRinputconversion); AddEntry("Metallicity_input_unit_conversion_to_output_unit", opt.metallicityinputconversion); AddEntry("Stellar_age_input_is_cosmological_scalefactor", opt.istellaragescalefactor); AddEntry("Hubble_unit", opt.H); AddEntry("Gravity", opt.G); AddEntry("Mass_value", opt.MassValue); AddEntry("Length_unit_to_kpc", opt.lengthtokpc); AddEntry("Velocity_to_kms", opt.velocitytokms); AddEntry("Mass_to_solarmass", opt.masstosolarmass); AddEntry("Star_formation_rate_to_solarmassperyear", opt.SFRtosolarmassperyear); AddEntry("Metallicity_to_solarmetallicity", opt.metallicitytosolar); AddEntry("Stellar_age_to_yr", opt.stellaragetoyrs); // simulation/cosmology info AddEntry("Period", opt.p); AddEntry("Scale_factor", opt.a); AddEntry("h_val", opt.h); AddEntry("Omega_m", opt.Omega_m); AddEntry("Omega_Lambda", opt.Omega_Lambda); AddEntry("Critical_density", opt.rhobg); AddEntry("Virial_density", opt.virlevel); AddEntry("Omega_cdm", opt.Omega_cdm); AddEntry("Omega_b", opt.Omega_b); AddEntry("Omega_r", opt.Omega_r); AddEntry("Omega_nu", opt.Omega_nu); AddEntry("Omega_k", opt.Omega_k); AddEntry("Omega_DE", opt.Omega_de); AddEntry("w_of_DE", opt.w_de); //unbinding AddEntry("Unbind_flag", opt.uinfo.unbindflag); AddEntry("Unbinding_type", opt.uinfo.unbindtype); AddEntry("Bound_halos", opt.iBoundHalos); AddEntry("Allowed_kinetic_potential_ratio", opt.uinfo.Eratio); AddEntry("Min_bound_mass_frac", opt.uinfo.minEfrac); AddEntry("Keep_background_potential", opt.uinfo.bgpot); AddEntry("Kinetic_reference_frame_type", opt.uinfo.cmvelreftype); AddEntry("Min_npot_ref", opt.uinfo.Npotref); AddEntry("Frac_pot_ref", opt.uinfo.fracpotref); AddEntry("Unbinding_max_unbound_removal_fraction_per_iteration", opt.uinfo.maxunbindfrac); AddEntry("Unbinding_max_unbound_fraction", opt.uinfo.maxunboundfracforiterativeunbind); AddEntry("Unbinding_max_unbound_fraction_allowed", opt.uinfo.maxallowedunboundfrac); AddEntry("Softening_length", opt.uinfo.eps); AddEntry("Approximate_potential_calculation", opt.uinfo.iapproxpot); AddEntry("Approximate_potential_calculation_particle_number_fraction", opt.uinfo.approxpotnumfrac); AddEntry("Approximate_potential_calculation_min_particle", opt.uinfo.approxpotminnum); AddEntry("Approximate_potential_calculation_method", opt.uinfo.approxpotmethod); //property related AddEntry("Inclusive_halo_masses", opt.iInclusiveHalo); AddEntry("Extensive_halo_properties_output", opt.iextrahalooutput); AddEntry("Extensive_gas_properties_output", opt.iextragasoutput); AddEntry("Extensive_star_properties_output", opt.iextrastaroutput); AddEntry("Extensive_interloper_properties_output", opt.iextrainterloperoutput); AddEntry("Iterate_cm_flag", opt.iIterateCM); AddEntry("Sort_by_binding_energy", opt.iSortByBindingEnergy); AddEntry("Reference_frame_for_properties", opt.iPropertyReferencePosition); AddEntry("Particle_type_for_reference_frames", opt.ParticleTypeForRefenceFrame); AddEntry("Calculate_aperture_quantities", opt.iaperturecalc); AddEntry("Number_of_apertures", opt.aperturenum); AddEntry("Aperture_values_in_kpc", opt.aperture_values_kpc); AddEntry("Number_of_projected_apertures", opt.apertureprojnum); AddEntry("Projected_aperture_values_in_kpc", opt.aperture_proj_values_kpc); AddEntry("Calculate_radial_profiles", opt.iprofilecalc); if(opt.iprofilecalc) { AddEntry("Radial_profile_min_FOF_size", opt.profileminFOFsize); AddEntry("Radial_profile_min_size", opt.profileminsize); AddEntry("Number_of_radial_profile_bin_edges", opt.profilenbins); AddEntry("Radial_profile_norm", opt.iprofilenorm); AddEntry("Radial_profile_bin_edges", opt.profile_bin_edges); } AddEntry("Number_of_overdensities", opt.SOnum); AddEntry("Overdensity_values_in_critical_density", opt.SOthresholds_values_crit); AddEntry("Spherical_overdenisty_calculation_limited_to_structure_types", (opt.SphericalOverdensitySeachMaxStructLevel-HALOSTYPE)/HALOCORESTYPE); //try removing index that is now stored in vector name; name = StripIndexOffName(opt.gas_internalprop_names, opt.gas_internalprop_names_aperture, opt.gas_internalprop_index, opt.gas_internalprop_index_aperture); AddEntry("Gas_internal_property_names", name); name = StripIndexOffName(opt.gas_chem_names, opt.gas_chem_names_aperture, opt.gas_chem_index, opt.gas_chem_index_aperture); AddEntry("Gas_chemistry_names", name); name = StripIndexOffName(opt.gas_chemproduction_names, opt.gas_chemproduction_names_aperture, opt.gas_chemproduction_index, opt.gas_chemproduction_index_aperture); AddEntry("Gas_chemistry_production_names", name); name = StripIndexOffName(opt.star_internalprop_names, opt.star_internalprop_names_aperture, opt.star_internalprop_index, opt.star_internalprop_index_aperture); AddEntry("Star_internal_property_names", name); name = StripIndexOffName(opt.star_chem_names, opt.star_chem_names_aperture, opt.star_chem_index, opt.star_chem_index_aperture); AddEntry("Star_chemistry_names", name); name = StripIndexOffName(opt.star_chemproduction_names, opt.star_chemproduction_names_aperture, opt.star_chemproduction_index, opt.star_chemproduction_index_aperture); AddEntry("Star_chemistry_production_names", name); name = StripIndexOffName(opt.bh_internalprop_names, opt.bh_internalprop_names_aperture, opt.bh_internalprop_index, opt.bh_internalprop_index_aperture); AddEntry("BH_internal_property_names", name); name = StripIndexOffName(opt.bh_chem_names, opt.bh_chem_names_aperture, opt.bh_chem_index, opt.bh_chem_index_aperture); AddEntry("BH_chemistry_names", name); name = StripIndexOffName(opt.bh_chemproduction_names, opt.bh_chemproduction_names_aperture, opt.bh_chemproduction_index, opt.bh_chemproduction_index_aperture); AddEntry("BH_chemistry_production_names", name); name = StripIndexOffName(opt.extra_dm_internalprop_names, opt.extra_dm_internalprop_names_aperture, opt.extra_dm_internalprop_index, opt.extra_dm_internalprop_index_aperture); AddEntry("Extra_DM_internal_property_names", name); map calcinttostring = ExtraFieldCalculationsIntFlagToString(); vector funcname; for (auto &x:opt.gas_internalprop_function) funcname.push_back(calcinttostring[x]); for (auto &x:opt.gas_internalprop_function_aperture) funcname.push_back(calcinttostring[x]); AddEntry("Gas_internal_property_calculation_type", funcname); funcname.clear(); for (auto &x:opt.gas_chem_function) funcname.push_back(calcinttostring[x]); for (auto &x:opt.gas_chem_function_aperture) funcname.push_back(calcinttostring[x]); AddEntry("Gas_chemistry_calculation_type", funcname); funcname.clear(); for (auto &x:opt.gas_chemproduction_function) funcname.push_back(calcinttostring[x]); for (auto &x:opt.gas_chemproduction_function_aperture) funcname.push_back(calcinttostring[x]); AddEntry("Gas_chemistry_production_calculation_type", funcname); funcname.clear(); for (auto &x:opt.star_internalprop_function) funcname.push_back(calcinttostring[x]); for (auto &x:opt.star_internalprop_function_aperture) funcname.push_back(calcinttostring[x]); AddEntry("Star_internal_property_calculation_type", funcname); funcname.clear(); for (auto &x:opt.star_chem_function) funcname.push_back(calcinttostring[x]); for (auto &x:opt.star_chem_function_aperture) funcname.push_back(calcinttostring[x]); AddEntry("Star_chemistry_calculation_type", funcname); funcname.clear(); for (auto &x:opt.star_chemproduction_function) funcname.push_back(calcinttostring[x]); for (auto &x:opt.star_chemproduction_function_aperture) funcname.push_back(calcinttostring[x]); AddEntry("Star_chemistry_production_calculation_type", funcname); funcname.clear(); for (auto &x:opt.bh_internalprop_function) funcname.push_back(calcinttostring[x]); for (auto &x:opt.bh_internalprop_function_aperture) funcname.push_back(calcinttostring[x]); AddEntry("BH_internal_property_calculation_type", funcname); funcname.clear(); for (auto &x:opt.bh_chem_function) funcname.push_back(calcinttostring[x]); for (auto &x:opt.bh_chem_function_aperture) funcname.push_back(calcinttostring[x]); AddEntry("BH_chemistry_calculation_type", funcname); funcname.clear(); for (auto &x:opt.bh_chemproduction_function) funcname.push_back(calcinttostring[x]); for (auto &x:opt.bh_chemproduction_function_aperture) funcname.push_back(calcinttostring[x]); AddEntry("BH_chemistry_production_calculation_type", funcname); funcname.clear(); for (auto &x:opt.extra_dm_internalprop_function) funcname.push_back(calcinttostring[x]); for (auto &x:opt.extra_dm_internalprop_function_aperture) funcname.push_back(calcinttostring[x]); AddEntry("Extra_DM_internal_property_calculation_type", funcname); funcname.clear(); AddEntry("Gas_internal_property_index", opt.gas_internalprop_index, opt.gas_internalprop_index_aperture); AddEntry("Gas_chemistry_index", opt.gas_chem_index, opt.gas_chem_index_aperture); AddEntry("Gas_chemistry_production_index", opt.gas_chemproduction_index, opt.gas_chemproduction_index_aperture); AddEntry("Star_internal_property_index", opt.star_internalprop_index, opt.star_internalprop_index_aperture); AddEntry("Star_chemistry_index", opt.star_chem_index, opt.star_chem_index_aperture); AddEntry("Star_chemistry_production_index", opt.star_chemproduction_index, opt.star_chemproduction_index_aperture); AddEntry("BH_internal_property_index", opt.bh_internalprop_index, opt.bh_internalprop_index_aperture); AddEntry("BH_chemistry_index", opt.bh_chem_index, opt.bh_chem_index_aperture); AddEntry("BH_chemistry_production_index", opt.bh_chemproduction_index, opt.bh_chemproduction_index_aperture); AddEntry("Extra_DM_internal_property_index", opt.extra_dm_internalprop_index, opt.extra_dm_internalprop_index_aperture); AddEntry("Gas_internal_property_output_units", opt.gas_internalprop_output_units, opt.gas_internalprop_output_units_aperture); AddEntry("Gas_chemistry_output_units", opt.gas_chem_output_units, opt.gas_chem_output_units_aperture); AddEntry("Gas_chemistry_production_output_units", opt.gas_chemproduction_output_units, opt.gas_chemproduction_output_units_aperture); AddEntry("Tcut_halogas", opt.temp_max_cut); AddEntry("Gas_temperature_input_unit_conversion_to_output_unit", opt.temp_input_output_unit_conversion_factor); AddEntry("Gas_hot_overdensity_normalisation", opt.hot_gas_overdensity_normalisation); AddEntry("Aperture_values_normalised_to_overdensity", opt.aperture_hotgas_normalised_to_overdensity); AddEntry("Star_internal_property_output_units", opt.star_internalprop_output_units, opt.star_internalprop_output_units_aperture); AddEntry("Star_chemistry_output_units", opt.star_chem_output_units, opt.star_chem_output_units_aperture); AddEntry("Star_chemistry_production_output_units", opt.star_chemproduction_output_units, opt.star_chemproduction_output_units_aperture); AddEntry("BH_internal_property_output_units", opt.bh_internalprop_output_units, opt.bh_internalprop_output_units_aperture); AddEntry("BH_chemistry_output_units", opt.bh_chem_output_units, opt.bh_chem_output_units_aperture); AddEntry("BH_chemistry_production_output_units", opt.bh_chemproduction_output_units, opt.bh_chemproduction_output_units_aperture); AddEntry("Extra_DM_internal_property_output_units", opt.extra_dm_internalprop_output_units, opt.extra_dm_internalprop_output_units_aperture); AddEntry("Gas_internal_property_input_output_unit_conversion_factors", opt.gas_internalprop_input_output_unit_conversion_factors, opt.gas_internalprop_input_output_unit_conversion_factors_aperture); AddEntry("Gas_chemistry_input_output_unit_conversion_factors", opt.gas_chem_input_output_unit_conversion_factors, opt.gas_chem_input_output_unit_conversion_factors_aperture); AddEntry("Gas_chemistry_production_input_output_unit_conversion_factors", opt.gas_chemproduction_input_output_unit_conversion_factors, opt.gas_chemproduction_input_output_unit_conversion_factors_aperture); AddEntry("Star_internal_property_input_output_unit_conversion_factors", opt.star_internalprop_input_output_unit_conversion_factors, opt.star_internalprop_input_output_unit_conversion_factors_aperture); AddEntry("Star_chemistry_input_output_unit_conversion_factors", opt.star_chem_input_output_unit_conversion_factors, opt.star_chem_input_output_unit_conversion_factors_aperture); AddEntry("Star_chemistry_production_input_output_unit_conversion_factors", opt.star_chemproduction_input_output_unit_conversion_factors, opt.star_chemproduction_input_output_unit_conversion_factors_aperture); AddEntry("BH_internal_property_input_output_unit_conversion_factors", opt.bh_internalprop_input_output_unit_conversion_factors, opt.bh_internalprop_input_output_unit_conversion_factors_aperture); AddEntry("BH_chemistry_input_output_unit_conversion_factors", opt.bh_chem_input_output_unit_conversion_factors, opt.bh_chem_input_output_unit_conversion_factors_aperture); AddEntry("BH_chemistry_production_input_output_unit_conversion_factors", opt.bh_chemproduction_input_output_unit_conversion_factors, opt.bh_chemproduction_input_output_unit_conversion_factors_aperture); AddEntry("Extra_DM_internal_property_input_output_unit_conversion_factors", opt.extra_dm_internalprop_input_output_unit_conversion_factors, opt.extra_dm_internalprop_input_output_unit_conversion_factors_aperture); //other options AddEntry("Verbose", opt.iverbose); AddEntry("Write_group_array_file",opt.iwritefof); AddEntry("Snapshot_value",opt.snapshotvalue); //io related AddEntry("Cosmological_input",opt.icosmologicalin); AddEntry("Input_chunk_size",opt.inputbufsize); AddEntry("MPI_particle_total_buf_size",opt.mpiparticletotbufsize); AddEntry("Separate_output_files", opt.iseparatefiles); AddEntry("Binary_output", opt.ibinaryout); AddEntry("Comoving_units", opt.icomoveunit); AddEntry("Extended_output", opt.iextendedoutput); //HDF io related info AddEntry("HDF_name_convention", opt.ihdfnameconvention); AddEntry("Input_includes_dm_particle", opt.iusedmparticles); AddEntry("Input_includes_gas_particle",opt.iusegasparticles); AddEntry("Input_includes_star_particle", opt.iusestarparticles); AddEntry("Input_includes_bh_particle", opt.iusesinkparticles); AddEntry("Input_includes_extradm_particle", opt.iuseextradarkparticles); AddEntry("Input_includes_wind_particle", opt.iusewindparticles); AddEntry("Input_includes_tracer_particle", opt.iusetracerparticles); //gadget io related to extra info for sph, stars, bhs, AddEntry("NSPH_extra_blocks", opt.gnsphblocks); AddEntry("NStar_extra_blocks", opt.gnstarblocks); AddEntry("NBH_extra_blocks", opt.gnbhblocks); //mpi related configuration AddEntry("MPI_part_allocation_fac", opt.mpipartfac); #endif AddEntry("#Compilation Info"); #ifdef USEMPI AddEntry("#USEMPI"); #endif #ifdef USEOPENMP AddEntry("#USEOPENMP"); #endif #ifdef GASON AddEntry("#USEGAS"); #endif #ifdef STARON AddEntry("#USESTAR"); #endif #ifdef BHON AddEntry("#USEBH"); #endif #ifdef EXTRADMON AddEntry("#USEEXTRADMPROPERTIES"); #endif #ifdef HIGHRES AddEntry("#ZOOMSIM"); #endif #ifdef SWIFTINTERFACE AddEntry("#SWIFTINTERFACE"); #endif }