Fixed CIC errors. Added support for CosmoPower in _cosmotool.pyx

This commit is contained in:
Guilhem Lavaux 2014-06-03 09:51:26 +02:00
parent 0a28eb3e04
commit 302ed9a912
8 changed files with 230 additions and 110 deletions

View file

@ -11,9 +11,12 @@ SET(CosmoTool_SRCS
miniargs.cpp
growthFactor.cpp
cosmopower.cpp
cic.cpp
)
IF (ENABLE_OPENMP)
set(CosmoTool_SRCS ${CosmoTool_SRCS} cic.cpp)
ENDIF (ENABLE_OPENMP)
IF(FOUND_NETCDF3)
SET(CosmoTool_SRCS ${CosmoTool_SRCS} yorick_nc3.cpp)
ELSE(FOUND_NETCDF3)

View file

@ -32,7 +32,7 @@ same conditions as regards security.
The fact that you are presently reading this means that you have had
knowledge of the CeCILL license and that you accept its terms.
+*/
#include <omp.h>
#include <assert.h>
#include <math.h>
#include <inttypes.h>
@ -60,89 +60,29 @@ void CICFilter::resetMesh()
void CICFilter::putParticles(CICParticles *particles, uint32_t N)
{
#if 0
uint32_t numCorners = 1 << NUMDIMS;
int threadUsed = omp_get_max_threads();
double *threadedDensity[threadUsed];
bool threadActivated[threadUsed];
uint32_t tUsedMin[threadUsed], tUsedMax[threadUsed];
for (uint32_t i = 0; i < N; i++)
{
Coordinates xyz;
int32_t ixyz[NUMDIMS];
int32_t rxyz[NUMDIMS];
CICType alpha[NUMDIMS];
CICType beta[NUMDIMS];
for (int j = 0; j < NUMDIMS; j++)
{
xyz[j] = (particles[i].coords[j] / spatialLen * szGrid);
ixyz[j] = (int32_t)floor(xyz[j] - 0.5);
beta[j] = xyz[j] - ixyz[j] - 0.5;
alpha[j] = 1 - beta[j];
if (ixyz[j] < 0)
ixyz[j] = szGrid-1;
}
for (int t = 0; t < threadUsed; t++)
{
threadedDensity[t] = new double[totalSize];
std::fill(threadedDensity[t], threadedDensity[t]+totalSize, 0);
}
CICType tot_mass = 0;
for (int j = 0; j < numCorners; j++)
{
CICType rel_mass = 1;
uint32_t idx = 0;
uint32_t mul = 1;
uint32_t mul2 = 1;
std::fill(threadActivated, threadActivated+threadUsed, false);
std::fill(tUsedMin, tUsedMin+threadUsed, totalSize);
std::fill(tUsedMax, tUsedMax+threadUsed, 0);
for (int k = 0; k < NUMDIMS; k++)
{
uint32_t ipos = ((j & mul2) != 0);
if (ipos == 1)
{
rel_mass *= beta[k];
}
else
{
rel_mass *= alpha[k];
}
rxyz[k] = ixyz[k] + ipos;
if (rxyz[k] >= szGrid)
idx += (rxyz[k] - szGrid) * mul;
else
idx += rxyz[k] * mul;
mul2 *= 2;
mul *= szGrid;
}
assert(rel_mass > 0);
assert(rel_mass < 1);
assert(idx < totalSize);
densityGrid[idx] += rel_mass * particles[i].mass;
tot_mass += rel_mass;
}
assert(tot_mass < 1.1);
assert(tot_mass > 0.9);
}
#endif
#if 0
for (uint32_t i = 0; i < N; i++)
{
Coordinates xyz;
int32_t ixyz[NUMDIMS];
for (int j = 0; j < NUMDIMS; j++)
{
xyz[j] = (particles[i].coords[j] / spatialLen * szGrid);
ixyz[j] = (int32_t)round(xyz[j] - 0.5);
if (ixyz[j] < 0)
ixyz[j] = szGrid-1;
else if (ixyz[j] >= szGrid)
ixyz[j] = 0;
}
uint32_t idx = ixyz[0] + ixyz[1] * szGrid + ixyz[2] * szGrid * szGrid;
densityGrid[idx] += particles[i].mass;
}
#endif
#pragma omp parallel
{
int thisThread = omp_get_thread_num();
double *dg = threadedDensity[thisThread];
threadActivated[thisThread] = true;
#pragma omp for schedule(static)
for (uint32_t i = 0; i < N; i++)
{
CICType x, y, z;
@ -193,44 +133,57 @@ void CICFilter::putParticles(CICParticles *particles, uint32_t N)
// 000
idx = ix + (iy + iz * szGrid) * szGrid;
densityGrid[idx] +=
dg[idx] +=
mass * beta_x * beta_y * beta_z;
// 100
idx = ix2 + (iy + iz * szGrid) * szGrid;
densityGrid[idx] +=
dg[idx] +=
mass * alpha_x * beta_y * beta_z;
// 010
idx = ix + (iy2 + iz * szGrid) * szGrid;
densityGrid[idx] +=
dg[idx] +=
mass * beta_x * alpha_y * beta_z;
// 110
idx = ix2 + (iy2 + iz * szGrid) * szGrid;
densityGrid[idx] +=
dg[idx] +=
mass * alpha_x * alpha_y * beta_z;
// 001
idx = ix + (iy + iz2 * szGrid) * szGrid;
densityGrid[idx] +=
dg[idx] +=
mass * beta_x * beta_y * alpha_z;
// 101
idx = ix2 + (iy + iz2 * szGrid) * szGrid;
densityGrid[idx] +=
dg[idx] +=
mass * alpha_x * beta_y * alpha_z;
// 011
idx = ix + (iy2 + iz2 * szGrid) * szGrid;
densityGrid[idx] +=
dg[idx] +=
mass * beta_x * alpha_y * alpha_z;
// 111
idx = ix2 + (iy2 + iz2 * szGrid) * szGrid;
densityGrid[idx] +=
dg[idx] +=
mass * alpha_x * alpha_y * alpha_z;
tUsedMin[thisThread] = std::min(tUsedMin[thisThread], idx);
tUsedMax[thisThread] = std::max(tUsedMax[thisThread], idx);
}
}
for (int t = 0; t < threadUsed; t++)
{
if (!threadActivated[t])
continue;
#pragma omp parallel for schedule(static)
for (long p = tUsedMin[t]; p < tUsedMax[t]; p++)
densityGrid[p] += threadedDensity[t][p];
}
}
void CICFilter::getDensityField(CICType*& field, uint32_t& res)
@ -238,3 +191,4 @@ void CICFilter::getDensityField(CICType*& field, uint32_t& res)
field = densityGrid;
res = totalSize;
}