vtk-m/vtkm/worklet/zfp/ZFPEncode1.h

112 lines
3.1 KiB
C
Raw Normal View History

2018-12-17 15:22:23 +00:00
//============================================================================
// Copyright (c) Kitware, Inc.
// All rights reserved.
// See LICENSE.txt for details.
2019-04-15 23:24:21 +00:00
//
2018-12-17 15:22:23 +00:00
// This software is distributed WITHOUT ANY WARRANTY; without even
// the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR
// PURPOSE. See the above copyright notice for more information.
//============================================================================
2018-12-07 19:02:51 +00:00
#ifndef vtk_m_worklet_zfp_encode1_h
#define vtk_m_worklet_zfp_encode1_h
#include <vtkm/Types.h>
#include <vtkm/internal/ExportMacros.h>
#include <vtkm/worklet/WorkletMapField.h>
#include <vtkm/worklet/zfp/ZFPBlockWriter.h>
#include <vtkm/worklet/zfp/ZFPEncode.h>
#include <vtkm/worklet/zfp/ZFPFunctions.h>
#include <vtkm/worklet/zfp/ZFPStructs.h>
#include <vtkm/worklet/zfp/ZFPTypeInfo.h>
namespace vtkm
{
namespace worklet
{
namespace zfp
{
template <typename Scalar, typename PortalType>
VTKM_EXEC inline void GatherPartial1(Scalar* q,
const PortalType& scalars,
vtkm::Id offset,
int nx,
int sx)
{
vtkm::Id x;
2018-12-07 21:02:16 +00:00
for (x = 0; x < nx; x++, offset += sx)
2018-12-07 19:02:51 +00:00
q[x] = scalars.Get(offset);
2018-12-18 04:04:34 +00:00
PadBlock(q, vtkm::UInt32(nx), 1);
2018-12-07 19:02:51 +00:00
}
template <typename Scalar, typename PortalType>
VTKM_EXEC inline void Gather1(Scalar* fblock, const PortalType& scalars, vtkm::Id offset, int sx)
{
vtkm::Id counter = 0;
for (vtkm::Id x = 0; x < 4; x++, offset += sx)
{
fblock[counter] = scalars.Get(offset);
counter++;
}
}
struct Encode1 : public vtkm::worklet::WorkletMapField
{
protected:
vtkm::Id Dims; // field dims
vtkm::Id PaddedDims; // dims padded to a multiple of zfp block size
vtkm::Id ZFPDims; // zfp block dims
vtkm::UInt32 MaxBits; // bits per zfp block
public:
Encode1(const vtkm::Id dims, const vtkm::Id paddedDims, const vtkm::UInt32 maxbits)
: Dims(dims)
, PaddedDims(paddedDims)
, MaxBits(maxbits)
{
ZFPDims = PaddedDims / 4;
}
using ControlSignature = void(FieldIn, WholeArrayIn, AtomicArrayInOut bitstream);
2018-12-07 19:02:51 +00:00
template <class InputScalarPortal, typename BitstreamPortal>
VTKM_EXEC void operator()(const vtkm::Id blockIdx,
const InputScalarPortal& scalars,
BitstreamPortal& stream) const
{
using Scalar = typename InputScalarPortal::ValueType;
vtkm::Id zfpBlock;
zfpBlock = blockIdx % ZFPDims;
vtkm::Id logicalStart = zfpBlock * vtkm::Id(4);
constexpr vtkm::Int32 BlockSize = 4;
Scalar fblock[BlockSize];
bool partial = false;
if (logicalStart + 4 > Dims)
partial = true;
if (partial)
{
const vtkm::Int32 nx =
logicalStart + 4 > Dims ? vtkm::Int32(Dims - logicalStart) : vtkm::Int32(4);
GatherPartial1(fblock, scalars, logicalStart, nx, 1);
}
else
{
Gather1(fblock, scalars, logicalStart, 1);
}
zfp::ZFPBlockEncoder<BlockSize, Scalar, BitstreamPortal> encoder;
2018-12-18 04:04:34 +00:00
encoder.encode(fblock, static_cast<vtkm::Int32>(MaxBits), vtkm::UInt32(blockIdx), stream);
2018-12-07 19:02:51 +00:00
}
};
}
}
}
#endif