/*
densities.hh - This file is part of MUSIC -
a code to generate multi-scale initial conditions
for cosmological simulations
Copyright (C) 2010 Oliver Hahn
This program is free software: you can redistribute it and/or modify
it under the terms of the GNU General Public License as published by
the Free Software Foundation, either version 3 of the License, or
(at your option) any later version.
This program is distributed in the hope that it will be useful,
but WITHOUT ANY WARRANTY; without even the implied warranty of
MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
GNU General Public License for more details.
You should have received a copy of the GNU General Public License
along with this program. If not, see .
*/
#ifndef __DENSITIES_HH
#define __DENSITIES_HH
#ifdef SINGLE_PRECISION
#ifdef SINGLETHREAD_FFTW
#include
#else
#include
#endif
#else
#ifdef SINGLETHREAD_FFTW
#include
#else
#include
#endif
#endif
#include
#include "config_file.hh"
#include "random.hh"
#include "cosmology.hh"
#include "transfer_function.hh"
#include "general.hh"
void GenerateDensityHierarchy( config_file& cf, transfer_function *ptf, tf_type type,
refinement_hierarchy& refh, rand_gen& rand, grid_hierarchy& delta, bool bdeconvolve, bool smooth );
void GenerateDensityUnigrid( config_file& cf, transfer_function *ptf, tf_type type,
refinement_hierarchy& refh, rand_gen& rand, grid_hierarchy& delta, bool kspace, bool deconvolve, bool smooth );
void normalize_density( grid_hierarchy& delta );
/*!
* @class DensityGrid
* @brief provides infrastructure for computing the initial density field
*
* This class provides access and data management member functions that
* are used when computing the initial density field by convolution with
* transfer functions.
*/
template< typename real_t >
class DensityGrid
{
public:
int nx_; //!< number of grid cells in x-direction
int ny_; //!< number of grid cells in y-direction
int nz_; //!< number of grid cells in z-direction
int nzp_; //!< number of cells in memory (z-dir), used for Nyquist padding
//! the actual data container in the form of a 1D array
std::vector< real_t > data_;
//! constructor
/*! constructs an instance given the dimensions of the density field
* @param nx the number of cells in x
* @param ny the number of cells in y
* @param nz the number of cells in z
*/
DensityGrid( unsigned nx, unsigned ny, unsigned nz )
: nx_(nx), ny_(ny), nz_(nz)
{
data_.assign(nx_*ny_*2*(nz_/2+1),0.0);
nzp_ = 2*(nz_/2+1);
}
//! copy constructor
explicit DensityGrid( const DensityGrid & g )
: nx_(g.nx_), ny_(g.ny_), nz_(g.nz_), nzp_(g.nzp_)
{
data_ = g.data_;
}
//!destructor
~DensityGrid()
{ }
//! clears the density object
/*! sets all dimensions to zero and frees the memory
*/
void clear( void )
{
nx_ = ny_ = nz_ = nzp_ = 0;
data_.clear();
std::vector().swap(data_);
}
//! query the 3D array sizes of the density object
/*! returns the size of the 3D density object along a specified dimension
* @param i the dimension for which size is to be returned
* @returns array size along dimension i
*/
int size( int i )
{
if(i==0) return nx_;
if(i==1) return ny_;
return ny_;
}
//! zeroes the density object
/*! sets all values to 0.0
*/
void zero( void )
{
data_.assign(data_.size(),0.0);
}
//! assigns the contents of another DensityGrid to this
DensityGrid& operator=( const DensityGrid& g )
{
nx_ = g.nx_;
ny_ = g.ny_;
nz_ = g.nz_;
nzp_= g.nzp_;
data_ = g.data_;
return *this;
}
//! 3D index based data access operator
inline real_t& operator()( int i, int j, int k )
{ return data_[(i*ny_+j)*nzp_+k]; }
//! 3D index based const data access operator
inline const real_t& operator()( int i, int j, int k ) const
{ return data_[(i*ny_+j)*nzp_+k]; }
//! recover the pointer to the 1D data array
inline real_t * get_data_ptr( void )
{ return &data_[0]; }
//! fills the density field with random number values
/*! given a pointer to a random_numbers object, fills the field with random values
* @param prc pointer to a random_numbers object
* @param variance the variance of the random numbers (the values returned by prc are multiplied by this)
* @param i0 x-offset (shift) in cells of the density field with respect to the random number field
* @param j0 y-offset (shift) in cells of the density field with respect to the random number field
* @param k0 z-offset (shift) in cells of the density field with respect to the random number field
* @param setzero boolean, if true, the global mean will be subtracted
*/
void fill_rand( /*const*/ random_numbers* prc, real_t variance, int i0, int j0, int k0, bool setzero=false )
{
double sum = 0.0;
#pragma omp parallel for reduction(+:sum)
for( int i=0; i
void copy( array3& v )
{
for( int ix=0; ix<(int)nx_; ++ix )
for( int iy=0; iy<(int)ny_; ++iy )
for( int iz=0; iz<(int)nz_; ++iz )
v(ix,iy,iz) = (*this)(ix,iy,iz);
}
//! adds the data from another field with 3D index-based access operator
template< class array3 >
void copy_add( array3& v )
{
for( int ix=0; ix<(int)nx_; ++ix )
for( int iy=0; iy<(int)ny_; ++iy )
for( int iz=0; iz<(int)nz_; ++iz )
v(ix,iy,iz) += (*this)(ix,iy,iz);
}
//! subtract the mean of the field
/*! subtracting the total mean implies that a restriction does not change the mass
*/
double subtract_mean( void )
{
double sum = 0.0;
unsigned count;
for( int i=0; i=oxsub+lxsub)
|| (iy=oysub+lysub)
|| (iz=ozsub+lzsub) )
(*this)(ix,iy,iz) = 0.0;
}
}
//! sets the field to zero outside of a rectangular box
void zero_but_subgrid( int oxsub, int oysub, int ozsub, int lxsub, int lysub, int lzsub )
{
for( int ix=0; ix=oxsub+lxsub)
|| (iy=oysub+lysub)
|| (iz=ozsub+lzsub) )
(*this)(ix,iy,iz) = 0.0;
}
}
//! sets the field to zero except on the boundary of a rectangular box
void zero_but_subgrid_bnd( unsigned oxsub, unsigned oysub, unsigned ozsub, unsigned lxsub, unsigned lysub, unsigned lzsub )
{
//... correct offsets for padding (not needed for top grid)
// zero the subgrid
for( int ix=oxsub; ix<(int)(oxsub+lxsub); ++ix )
for( int iy=oysub; iy<(int)(oysub+lysub); ++iy )
for( int iz=ozsub; iz<(int)(ozsub+lzsub); ++iz )
(*this)(ix,iy,iz) = 0.0;
oxsub -= lxsub/2;
oysub -= lysub/2;
ozsub -= lzsub/2;
lxsub *= 2;
lysub *= 2;
lzsub *= 2;
// zero the outside, except the boundary
//#pragma parallel nowait
{
for( int ix=0; ix
class PaddedDensitySubGrid : public DensityGrid
{
public:
using DensityGrid::nx_;
using DensityGrid::ny_;
using DensityGrid::nz_;
using DensityGrid::data_;
using DensityGrid::fill_rand;
using DensityGrid::get_data_ptr;
int ox_, oy_, oz_;
public:
PaddedDensitySubGrid( int ox, int oy, int oz, unsigned nx, unsigned ny, unsigned nz)
: DensityGrid(2*nx,2*ny,2*nz), ox_(ox), oy_(oy), oz_(oz)
{
//if( 2*ox-(int)nx/4 < 0 || 2*oy-(int)ny/4 < 0 || 2*oz-(int)nz/4 < 0 )
// throw std::runtime_error("subgrid extends across top grid");
//.. the size in top grid cells is nx/2,ny/2,nz/2, so padding starts at ox-nx/4...
//.. loop over relevant part of top grid and copy down grid values
/*for( int ix=ox-nx/4,ixs=0; ix& o )
: DensityGrid( o ), ox_(o.ox_), oy_(o.oy_), oz_(o.oz_)
{ }
void zero_padded_subgrid( unsigned oxsub, unsigned oysub, unsigned ozsub, unsigned lxsub, unsigned lysub, unsigned lzsub )
{
//... correct offsets for padding (not needed for top grid)
oxsub += nx_/4;
oysub += ny_/4;
ozsub += nz_/4;
//... correct offsets for padding (not needed for top grid)
for( int ix=oxsub; ix<(int)(oxsub+lxsub); ++ix )
for( int iy=oysub; iy<(int)(oysub+lysub); ++iy )
for( int iz=ozsub; iz<(int)(ozsub+lzsub); ++iz )
(*this)(ix,iy,iz) = 0.0;
}
void zero_subgrid( unsigned oxsub, unsigned oysub, unsigned ozsub, unsigned lxsub, unsigned lysub, unsigned lzsub )
{
//... correct offsets for padding (not needed for top grid)
oxsub += nx_/4;
oysub += ny_/4;
ozsub += nz_/4;
for( int ix=oxsub; ix<(int)(oxsub+lxsub); ++ix )
for( int iy=oysub; iy<(int)(oysub+lysub); ++iy )
for( int iz=ozsub; iz<(int)(ozsub+lzsub); ++iz )
(*this)(ix,iy,iz) = 0.0;
}
void zero_but_subgrid_bnd( unsigned oxsub, unsigned oysub, unsigned ozsub, unsigned lxsub, unsigned lysub, unsigned lzsub )
{
//... correct offsets for padding (not needed for top grid)
oxsub += nx_/4;
oysub += ny_/4;
ozsub += nz_/4;
//lxsub += nx_/4;
//lysub += nx_/4;
//lzsub += nx_/4;
// zero the subgrid
for( int ix=oxsub; ix<(int)(oxsub+lxsub); ++ix )
for( int iy=oysub; iy<(int)(oysub+lysub); ++iy )
for( int iz=ozsub; iz<(int)(ozsub+lzsub); ++iz )
(*this)(ix,iy,iz) = 0.0;
oxsub -= lxsub/2;
oysub -= lysub/2;
ozsub -= lzsub/2;
lxsub *= 2;
lysub *= 2;
lzsub *= 2;
// zero the outside, except the boundary
//#pragma parallel block nowait
{
for( int ix=0; ix " << nx_ << std::endl;
for( int ix=oxsub+lxsub; ix
void upload_bnd_add( array3& top )
{
for( int ix=0; ix=3*nx_/4) || (iy=3*ny_/4) || (iz=3*nz_/4) )
{
int ic,jc,kc;
ic = ox_+(ix-nx_/4)/2;
jc = oy_+(iy-ny_/4)/2;
kc = oz_+(iz-nz_/4)/2;
if( ic >=0 && ic < (int)top.size(0) && jc >= 0 && jc < (int)top.size(1) && kc >= 0 && kc < (int)top.size(2) )
{
top(ic,jc,kc) += 0.125* ((*this)(ix,iy,iz)+(*this)(ix+1,iy,iz)+(*this)(ix,iy+1,iz)+(*this)(ix,iy,iz+1)
+(*this)(ix+1,iy+1,iz)+(*this)(ix+1,iy,iz+1)+(*this)(ix,iy+1,iz+1)+(*this)(ix+1,iy+1,iz+1) );
}
}
}
}
void constrain( const DensityGrid& top )
{
int ox=ox_-nx_/8, oy=oy_-ny_/8, oz=oz_-nz_/8;
if( ox < 0 || oy < 0 || oz < 0 ){
std::cerr << "offset = " << ox << ", " << oy << ", " << oz << std::endl;
std::cerr << "nx = " << nx_ << ", " << ny_ << ", " << nz_ << std::endl;
throw std::runtime_error("Subgrid extends beyond top grid. Increase levelmin or reduce size of refinement region!");
}
//...top grid is not padded
//.. fine grid will be multiplied by sqrt(8) later to account for increase in variance
//.. so, we need to divide the coarse grid value by it
real_t fac = 1.0/sqrt(8);
for( int ixs=0,ix=ox;ixs& top )
{
///... top grid is padded too...
//.. ox is shifted by nx_/4, padded overlap starts at -nx_/8
int ox=ox_+top.nx_/4-nx_/8, oy=oy_+top.ny_/4-ny_/8, oz=oz_+top.nz_/4-nz_/8;
if( ox < 0 || oy < 0 || oz < 0 )
throw std::runtime_error("subgrid extends beyond top grid");
double fac = 1.0/sqrt(8.0);
for( int ixs=0,ix=ox;ixs
void upload_bnd( unsigned nbnd, array3& v )
{
for( int ix=nx_/4-2*nbnd,ixu=-nbnd; ix<3*nx_/4+2*nbnd; ix+=2,++ixu)
for( int iy=ny_/4-2*nbnd,iyu=-nbnd; iy<3*ny_/4+2*nbnd; iy+=2,++iyu )
for( int iz=nz_/4-2*nbnd,izu=-nbnd; iz<3*nz_/4+2*nbnd; iz+=2,++izu )
{
if( (ix=3*nx_/4)
|| (iy=3*ny_/4)
|| (iz=3*nz_/4) )
{
v(ixu+ox_,iyu+oy_,izu+oz_)
= 1.0;
}
}
}
template< class array3 >
void copy_unpad( array3& v )
{
for( int ix=nx_/4,ixu=0; ix<3*nx_/4; ++ix,++ixu )
for( int iy=ny_/4,iyu=0; iy<3*ny_/4; ++iy,++iyu )
for( int iz=nz_/4,izu=0; iz<3*nz_/4; ++iz,++izu )
v(ixu,iyu,izu) = (*this)(ix,iy,iz);
}
template< class array3 >
void copy_add_unpad( array3& v )
{
for( int ix=nx_/4,ixu=0; ix<3*nx_/4; ++ix,++ixu )
for( int iy=ny_/4,iyu=0; iy<3*ny_/4; ++iy,++iyu )
for( int iz=nz_/4,izu=0; iz<3*nz_/4; ++iz,++izu )
v(ixu,iyu,izu) += (*this)(ix,iy,iz);
}
template< class array3 >
void copy_subtract_unpad( array3& v )
{
for( int ix=nx_/4,ixu=0; ix<3*nx_/4; ++ix,++ixu )
for( int iy=ny_/4,iyu=0; iy<3*ny_/4; ++iy,++iyu )
for( int iz=nz_/4,izu=0; iz<3*nz_/4; ++iz,++izu )
v(ixu,iyu,izu) -= (*this)(ix,iy,iz);
}
double oct_mean( int i, int j, int k )
{
return 0.125*((*this)(i,j,k)+(*this)(i+1,j,k)+(*this)(i,j+1,k)+(*this)(i,j,k+1)
+(*this)(i+1,j+1,k)+(*this)(i+1,j,k+1)+(*this)(i,j+1,k+1)+(*this)(i+1,j+1,k+1));
}
void set_oct_mean( int i, int j, int k, double val )
{
(*this)(i,j,k) = val;
(*this)(i+1,j,k) = val;
(*this)(i,j+1,k) = val;
(*this)(i,j,k+1) = val;
(*this)(i+1,j+1,k) = val;
(*this)(i+1,j,k+1) = val;
(*this)(i,j+1,k+1) = val;
(*this)(i+1,j+1,k+1) = val;
}
void add_oct_val( int i, int j, int k, double val )
{
(*this)(i,j,k) += val;
(*this)(i+1,j,k) += val;
(*this)(i,j+1,k) += val;
(*this)(i,j,k+1) += val;
(*this)(i+1,j+1,k) += val;
(*this)(i+1,j,k+1) += val;
(*this)(i,j+1,k+1) += val;
(*this)(i+1,j+1,k+1) += val;
}
void subtract_boundary_oct_mean( void )
{
#pragma omp parallel for
for( int ix=0; ix
inline void enforce_coarse_mean( M& v, M& V )
{
double outersum =0.0, innersum = 0.0;
unsigned
nx = v.size(0)/2,
ny = v.size(1)/2,
nz = v.size(2)/2,
ox = v.offset(0),
oy = v.offset(1),
oz = v.offset(2);
unsigned innercount = 0, outercount = 0;
for( unsigned i=0; i= ox && i < ox+nx) &&
(j >= oy && j < oy+ny) &&
(k >= oz && k < oz+nz ))
{
innersum += V(i,j,k);
++innercount;
}
else
{
outersum += V(i,j,k);
++outercount;
}
}
innersum /= innercount;
outersum /= outercount;
std::cerr << "***\n";
double dcorr = innersum * innercount / outercount;
for( unsigned i=0; i= ox && i < ox+nx) &&
(j >= oy && j < oy+ny) &&
(k >= oz && k < oz+nz )))
V(i,j,k) -= outersum + dcorr;//innersum;
}
}
template< typename M >
inline void enforce_mean( M& v, M& V )
{
int
nx = v.size(0)/2,
ny = v.size(1)/2,
nz = v.size(2)/2,
ox = v.offset(0),
oy = v.offset(1),
oz = v.offset(2);
double coarsemean = 0.0, finemean = 0.0;
unsigned count = 0;
#pragma omp parallel for reduction(+:coarsemean,finemean,count)
for( int i=0; i