void Delta::BuildDensityGrid(Particle& pp, Vec& delta) // RS 2010/07/06: Commented! { _mydestroy(delta); // Clean up delta DensityGrid dg2(N, L); DensityGrid dg3(N, L); delta = dg1.Allocate(); bool rmask = (rmasksmooth > 1e-5); // RS: boolean shorthand // First the particles // RS 2010/07/06: CIC assign the list of particles to delta. // If we're using a mask made from smoothed randoms, don't bother trimming, // because the mask will just be a "dummy mask". if (!rmask) pp.TrimMask(*maskptr); pp.SlabDecompose(dg1); dg1.CIC(delta,pp, false); PetscPrintf(PETSC_COMM_WORLD, "Particles assigned to the grid.....\n"); // RS 2010/07/06: If rmasksmooth > 0, smooth CIC-assigned particles before // computing overdensity. Otherwise we'll get unpleasant edge effects. if (rmask) dg1.GaussSmooth(delta, rmasksmooth); // Normalization factor // RS: This equals the *average* number (per allowed cell) of "randoms" // compared to the number of particles from which we're building delta. VecPointwiseMult(delta, W); double fac; VecSum(delta, &fac); fac = nrand/fac; // Now construct delta int lo, hi; dg1.slab(lo, hi); // RS: Set up the wrappers to allow us to index delta by (x,y,z). dg1.Pull(delta); dg2.Pull(W); dg3.Pull(_W); for (int ix=lo; ix < hi; ++ix) for (int iy=0; iy < N; ++iy) for (int iz=0; iz < N; ++iz) if (dg2(ix, iy,iz) > 1.e-5) { // RS: In each cell location, if W > 0, divide by _W/fac. // Since _W is the average density of regular-grid "randoms" in // each cell (accounting for edge effects), _W/fac corresponds to // the average density of pp particles. So W/(_W/fac) = W/_W*fac, // i.e., the density contrast, which is what we want. dg1(ix, iy, iz) = (dg1(ix, iy, iz)/dg3(ix, iy, iz))*fac - 1.0; } else { // RS: If W = 0 (to within epsilon), set delta to zero. dg1(ix, iy, iz) = 0.0; } // RS: Deallocate the memory we just set up... and we'll be done. dg1.Push(delta, INSERT_VALUES, false); dg2.Push(W, INSERT_VALUES, false); dg3.Push(_W, INSERT_VALUES, false); VecSum(delta, &fac); PetscPrintf(PETSC_COMM_WORLD, "Sum (delta).....%e\n", fac); }
int main(int argc, char *args[]) { PetscErrorCode ierr; ierr=PetscInitialize(&argc,&args,(char *) 0, help); CHKERRQ(ierr); { // Get MPI rank int rank; MPI_Comm_rank(PETSC_COMM_WORLD, &rank); // Read in configuration file char configfn[200]; PetscTruth flg; // See if we need help PetscOptionsHasName("-help", &flg); if (flg) exit(0); PetscOptionsGetString("-configfn", configfn, 200, &flg); if (!flg) RAISE_ERR(99,"Specify configuration file"); // Read in parameters ShellParams pars(string(configfn), "shellradial"); /**************************** * Define seeds ****************************/ int const_seed=2147; int randompart_seed=8271; /****************************************** * Define the mask and delta ******************************************/ Shell maskss(pars.shell.xcen, pars.shell.ycen, pars.shell.zcen, pars.shell.rmin, pars.shell.rmax); Delta del1(pars.Ngrid, pars.Lbox, pars.shell.nover, pars.shell.thresh, maskss); /******************************************** * Read in pk prior *******************************************/ vector<double> kvec, pkvec; int nk; // Rank 0 job reads in pk function and broadcasts to everyone if (rank == 0) { double k1, pk1; ifstream fp(pars.pkprior.fn.c_str()); do { fp >> k1 >> pk1; if (fp) { kvec.push_back(k1); pkvec.push_back(pars.pkprior.bias*pk1 + pars.pkprior.noise); } } while (fp); fp.close(); nk = kvec.size(); } MPI_Bcast(&nk, 1, MPI_INT, 0, PETSC_COMM_WORLD); if (rank !=0) {kvec.resize(nk, 0.0); pkvec.resize(nk, 0.0);} MPI_Bcast(&kvec[0], nk, MPI_DOUBLE, 0, PETSC_COMM_WORLD); MPI_Bcast(&pkvec[0], nk, MPI_DOUBLE, 0, PETSC_COMM_WORLD); PetscPrintf(PETSC_COMM_WORLD, "%d lines read in from %s...\n",nk, pars.pkprior.fn.c_str()); // Define the potential solver PotentialSolve psolve(pars.Ngrid, pars.Lbox, pars.recon.maxit); psolve.SetupOperator(RADIAL, pars.recon.fval, pars.recon.origin); // Loop over files list<ShellParams::fn>::iterator files; for (files = pars.fnlist.begin(); files !=pars.fnlist.end(); ++files) { /* ****************************** * First we get the various options and print out useful information * ********************************/ ostringstream hdr; hdr << "# Input file is " << files->in << endl; hdr << "# Output file is " << files->out << endl; hdr << "# Ngrid=" << setw(5) << pars.Ngrid << endl; hdr << "# boxsize=" << setw(8) << fixed << setprecision(2) << pars.Lbox << endl; hdr << "# bias=" << setw(8) << setprecision(2) << pars.recon.bias << endl; hdr << "# fval=" << setw(8) << setprecision(2) << pars.recon.fval << endl; hdr << "# smooth=" << setw(8) << setprecision(2) << pars.recon.smooth << endl; hdr << "# " << setw(4) << pars.xi.Nbins << " Xi bins from " << setw(8) << setprecision(2) << pars.xi.rmin << " with spacing of " << pars.xi.dr << endl; hdr << "# " << "Correlation function smoothed with a smoothing scale of" << setw(8) << setprecision(2) << pars.xi.smooth << endl; hdr << "# " << "Survey shell centered at the origin, rmin=" << setprecision(2) << pars.shell.rmin << ", rmax=" << pars.shell.rmax << endl; hdr << "# " << "Shell center at " << setprecision(3) << pars.shell.xcen << ", " << pars.shell.ycen << ", " << pars.shell.zcen << endl; hdr << "# Mask threshold =" << pars.shell.thresh << endl; PetscPrintf(PETSC_COMM_WORLD, (hdr.str()).c_str()); /**************************************** * Read in the particle data here and slab decompose ****************************************/ Particle pp; DensityGrid dg(pars.Ngrid, pars.Lbox); pp.TPMReadSerial(files->in.c_str(), pars.Lbox, true, true); PetscPrintf(PETSC_COMM_WORLD,"Read in %i particles.....\n",pp.npart); pp.SlabDecompose(dg); // Test slab decomp bool good = dg.TestSlabDecompose(pp); if (good) {PetscSynchronizedPrintf(PETSC_COMM_WORLD,"Slab decomposition succeeded on process %i\n",rank);} else {PetscSynchronizedPrintf(PETSC_COMM_WORLD,"Slab decomposition FAILED on process %i\n",rank);} PetscSynchronizedFlush(PETSC_COMM_WORLD); if (!good) RAISE_ERR(99, "Slab decomposition failed"); /************************************************* * Generate the density field and constrained realization *************************************************/ Vec grid=PETSC_NULL; del1.BuildDensityGridSmooth(pars.recon.smooth, pp, grid); PetscPrintf(PETSC_COMM_WORLD, "Density grid computed.....\n"); del1.HoffmanRibak(grid, kvec, pkvec, const_seed); const_seed+=1; PetscPrintf(PETSC_COMM_WORLD, "Constrained realization computed.....\n"); /************************************************ * Now we solve for the potential ************************************************/ // Allocate potential solver Vec pot; pot = dg.Allocate(); if (psolve.Solve(grid, pot, pars.recon.bias)) { // If the potential calculation converged PetscPrintf(PETSC_COMM_WORLD,"Potential calculated....\n"); /************************************************ * Now we shift data and randoms ************************************************/ // Generate random particles Vec dp, qx, qy, qz; Particle pr; pr.RandomInit(pp.npart*pars.recon.nrandomfac, pars.Lbox, randompart_seed); randompart_seed +=1; pr.SlabDecompose(dg); pr.TrimMask(maskss); // Compute derivatives at data positions and shift dp = dg.Deriv(pot, 0); qx = dg.Interp3d(dp, pp); _mydestroy(dp); dp = dg.Deriv(pot, 1); qy = dg.Interp3d(dp, pp); _mydestroy(dp); dp = dg.Deriv(pot, 2); qz = dg.Interp3d(dp, pp); _mydestroy(dp); // Print some statistics double sum[3]; VecSum(qx,&sum[0]); VecSum(qy, &sum[1]); VecSum(qz, &sum[2]); for (int ii=0; ii < 3; ++ii) sum[ii] /= pp.npart; PetscPrintf(PETSC_COMM_WORLD, "Mean x,y,z displacements on particles is : %10.4f,%10.4f,%10.4f\n",sum[0],sum[1],sum[2]); VecNorm(qx,NORM_2,&sum[0]); VecNorm(qy, NORM_2,&sum[1]); VecNorm(qz, NORM_2,&sum[2]); for (int ii=0; ii < 3; ++ii) sum[ii] /= sqrt(pp.npart); PetscPrintf(PETSC_COMM_WORLD, "RMS x,y,z displacements on particles is : %10.4f,%10.4f,%10.4f\n",sum[0],sum[1],sum[2]); // Do the shift in two pieces -- first the real space shift VecAXPY(pp.px, -1.0, qx); VecAXPY(pp.py, -1.0, qy); VecAXPY(pp.pz, -1.0, qz); PetscPrintf(PETSC_COMM_WORLD, "Standard displacements complete\n"); // Now shift to remove redshift distortions - this requires a little coordinate geometry // We need q.p/|r| pp.RadialDisplace(qx, qy, qz, pars.recon.origin, -pars.recon.fval); PetscPrintf(PETSC_COMM_WORLD, "zspace displacements complete\n"); // Cleanup _mydestroy(qx); _mydestroy(qy); _mydestroy(qz); // Do the same for the randoms dp = dg.Deriv(pot, 0); qx = dg.Interp3d(dp, pr); _mydestroy(dp); dp = dg.Deriv(pot, 1); qy = dg.Interp3d(dp, pr); _mydestroy(dp); dp = dg.Deriv(pot, 2); qz = dg.Interp3d(dp, pr); _mydestroy(dp); VecSum(qx,&sum[0]); VecSum(qy, &sum[1]); VecSum(qz, &sum[2]); for (int ii=0; ii < 3; ++ii) sum[ii] /= pr.npart; PetscPrintf(PETSC_COMM_WORLD, "Mean x,y,z displacements on randoms is : %10.4f,%10.4f,%10.4f\n",sum[0],sum[1],sum[2]); VecNorm(qx,NORM_2,&sum[0]); VecNorm(qy, NORM_2,&sum[1]); VecNorm(qz, NORM_2,&sum[2]); for (int ii=0; ii < 3; ++ii) sum[ii] /= sqrt(pr.npart); PetscPrintf(PETSC_COMM_WORLD, "RMS x,y,z displacements on randoms is : %10.4f,%10.4f,%10.4f\n",sum[0],sum[1],sum[2]); VecAXPY(pr.px, -1.0, qx); VecAXPY(pr.py, -1.0, qy); VecAXPY(pr.pz, -1.0, qz); PetscPrintf(PETSC_COMM_WORLD,"Displacements calculated....\n"); // Clean up _mydestroy(qx); _mydestroy(qy); _mydestroy(qz); // Shifted data and random grid Vec gridr=PETSC_NULL; del1.BuildDensityGrid(pp, grid); del1.BuildDensityGrid(pr, gridr); VecAXPY(grid, -1.0, gridr); // Correlation fn PkStruct xi2(pars.xi.rmin, pars.xi.dr, pars.xi.Nbins); dg.XiFFT_W(grid, del1.W, pars.xi.smooth, xi2); _mydestroy(gridr); // Outputs FILE *fp; double _rvec, _xi2, _n1; PetscFOpen(PETSC_COMM_WORLD,files->out.c_str(),"w", &fp); PetscFPrintf(PETSC_COMM_WORLD, fp, (hdr.str()).c_str()); for (int ii = xi2.lo; ii < xi2.hi; ++ii) { _xi2 = xi2(ii, _rvec, _n1); if (_n1>0) PetscSynchronizedFPrintf(PETSC_COMM_WORLD, fp, "%6i %9.3f %15.8e\n",ii,_rvec,_xi2); } PetscSynchronizedFlush(PETSC_COMM_WORLD); PetscFClose(PETSC_COMM_WORLD,fp); } // Cleanup _mydestroy(grid); _mydestroy(pot); } }
Delta::Delta(int _N, double _L, int _nover, double thresh, Mask3D& _mask) // RS 2010/07/06: Comments! { // Cache values N=_N; L=_L; nover=_nover; maskptr=&_mask; rmasksmooth = 0.0; // Set up the density // RS: dg2 is some kind of local variable; dg1 is a class member. // W and _W are contents of a grid of the same size and dimensions // as dg1 and dg2. DensityGrid dg2(N, L); dg1.Init(N, L); _W = dg1.Allocate(); W = dg1.Allocate(); VecSet(W, 0.0); // Now do the randoms // RS: "Doing" the randoms in this context means "oversampling" each grid // cell into subcells, where the main cell is nover times as large as each // subcell and hence there are nover^3 subcells to a cell. nrand=0.0; double nexp, nthresh; nexp = pow((double) nover, 3); // RS: thresh is the *fraction* of the main cell's volume to be included // inside the mask for it to be counted as not masked, so nthresh is the // corresponding number of subcells. nthresh = nexp*thresh; { double tmp = 1./nover, dx, dy, dz; VecSet(_W, 0.0); // Manually set to zero Vec W1; VecDuplicate(_W, &W1); // RS: Loop over all subcell *offsets*. The loop over actual grid cells // is implicit in the calls to methods like GridInitialize below. for (int ix = 0; ix < nover; ++ix) { PetscPrintf(PETSC_COMM_WORLD,"Starting random loop %i of %i....\n",ix,nover-1); dx = ix*tmp + tmp/2.0; for (int iy = 0; iy < nover; ++iy) { dy = iy*tmp + tmp/2.0; for (int iz = 0; iz < nover; ++iz) { dz = iz*tmp + tmp/2.0; Particle rr; // Definition here is important, since it ensures destruction, before recreating // RS: Generate a bunch of particles uniformly spaced, and offset // from the low-(x,y,z) corner of each grid cell by (dx,dy,dz). rr.GridInitialize(N, dx, dy, dz, L); // RS: Remove from the list any particles which were not allowed // by the mask. Add to the total. rr.TrimMask(*maskptr); nrand += rr.npart; // rr.SlabDecompose(dg2); We don't need to slab decompose because nproc divides N by construction // RS: CIC assign remaining particles to temporary grid vector W1. dg1.CIC(W1, rr, false); // No overdensity // RS: Now add to _W the CIC contents of W1. VecAXPY(_W, 1.0, W1); } } } } // RS: At this point, nrand equals the total number of oversampled subcells // found to be allowed by the masked region. _W contains the number of // particles allowed by an oversampled mask at each grid location. PetscPrintf(PETSC_COMM_WORLD, "%f grid particles accumulated....\n", nrand); // Now set the mask // RS: For each cell, if the contents of _W exceed nthresh, set W = 1. // Use dg1 and dg2 as wrappers to index the copies of _W and W. int lo, hi; dg1.slab(lo, hi); dg1.Pull(_W); dg2.Pull(W); for (int ix=lo; ix < hi; ++ix) for (int iy=0; iy < N; ++iy) for (int iz=0; iz < N; ++iz) if (dg1(ix, iy,iz) > nthresh) { dg2(ix, iy,iz) = 1.0; } else { dg1(ix, iy, iz) = 0.0; // Good to explicitly set this to zero as well } dg1.Push(_W, INSERT_VALUES, false); dg2.Push(W, INSERT_VALUES, false); double tmp; VecSum(W, &tmp); PetscPrintf(PETSC_COMM_WORLD,"Nonzero grid points in W : %f\n",tmp); VecSum(_W, &nrand); PetscPrintf(PETSC_COMM_WORLD,"Sum of background density vector : %f\n",nrand); }
void BuildDensityGrid(int N, double L, Mask3D& mask, Particle& pp, int nover, double thresh, Vec& delta, Vec& W, double smooth) { /* Parameters : * INPUT --- * N : grid size * L : box size * mask : A 3D positive mask (assumed to be 1 or 0) * nover : Oversampling rate for defining the mean density (we use the GridInitialize routine) * This determines the displacement. * For a periodic box, nover means there will be nover**3 particles per grid. * thresh: Fractional threshold below which to remove grid point from survey * smooth: This is an optional parameter -- if set to > 0, it will smooth the density and mask * fields before computing the overdensity. * OUTPUT --- * delta : What do think this is??? * W : Window function, 1 in survey, 0 elsewhere. * Note that the output vectors are not cleared, so you had better make sure you don't leak memory. */ // Allocate two density grids DensityGrid dg1(N,L), dg2(N,L); delta = dg1.Allocate(); W = dg2.Allocate(); // First the particles pp.TrimMask(mask); pp.SlabDecompose(dg1); dg2.CIC(delta,pp, false); PetscPrintf(PETSC_COMM_WORLD, "Particles assigned to the grid.....\n"); // Now do the randoms double nrand=0.0, nexp, nthresh; nexp = pow((double) nover, 3); nthresh = nexp*thresh; { double tmp = 1./nover, dx, dy, dz; VecSet(W, 0.0); // Manually set to zero Vec W1; VecDuplicate(W, &W1); for (int ix = 0; ix < nover; ++ix) { PetscPrintf(PETSC_COMM_WORLD,"Starting random loop %i of %i....\n",ix,nover-1); dx = ix*tmp + tmp/2.0; for (int iy = 0; iy < nover; ++iy) { dy = iy*tmp + tmp/2.0; for (int iz = 0; iz < nover; ++iz) { dz = iz*tmp + tmp/2.0; Particle rr; // Definition here is important, since it ensures destruction, before recreating rr.GridInitialize(N, dx, dy, dz, L); rr.TrimMask(mask); nrand += rr.npart; // rr.SlabDecompose(dg2); We don't need to slab decompose because nproc divides N by construction dg2.CIC(W1, rr, false); // No overdensity VecAXPY(W, 1.0, W1); } } } } PetscPrintf(PETSC_COMM_WORLD, "%f grid particles accumulated....\n", nrand); PetscPrintf(PETSC_COMM_WORLD, "Expected number of particles=%f, threshold=%f\n", nexp, nthresh); // Smooth if desired if (smooth > 1.e-5) { dg1.GaussSmooth(delta, smooth); dg2.GaussSmooth(W, smooth); } // Now build up delta int lo, hi; nrand /= pp.npart; dg1.slab(lo, hi); dg1.Pull(delta); dg2.Pull(W); for (int ix=lo; ix < hi; ++ix) for (int iy=0; iy < N; ++iy) for (int iz=0; iz < N; ++iz) { if (dg2(ix, iy,iz) > nthresh) { dg1(ix, iy,iz) = (dg1(ix,iy,iz)/dg2(ix, iy, iz)) * nrand - 1.0; dg2(ix, iy,iz) = 1.0; } else { dg1(ix, iy, iz) = 0.0; dg2(ix, iy, iz) = 0.0; } } dg1.Push(delta, INSERT_VALUES, false); dg2.Push(W, INSERT_VALUES, false); VecSum(W, &nrand); PetscPrintf(PETSC_COMM_WORLD,"Nonzero grid points in W : %f\n",nrand); }