25#define mycout QwMessage
27#define myendl QwLog::endl
53template <
class value_t =
float,
unsigned int value_n = 1>
67 const std::vector<coord_t>& max,
68 const std::vector<coord_t>& step) {
133 mycerr <<
"Dimensions of the grid should be strictly positive!" <<
myendl;
143 std::vector<coord_t>(
fNDim,max),
144 std::vector<coord_t>(
fNDim,step));
148 const std::vector<coord_t>& max,
149 const std::vector<coord_t>& step) {
153 if (min.size() !=
fNDim || min.size() !=
fNDim || step.size() !=
fNDim)
return;
156 for (
size_t i = 0; i <
fMin.size(); i++) {
159 fSize[i] =
static_cast<unsigned int>(int_part) + (frac_part > 0.5 ? 1 : 0);
164 for (
unsigned int i = 0; i < value_n; i++) {
167 }
catch (std::bad_alloc& e) {
168 mycerr <<
"QwInterpolator could not allocate memory for values!" <<
myendl;
185 {
return fWrap.at(dim); }
188 {
fWrap.at(dim) = wrap; }
190 if (wrap.size() !=
fNDim)
return;
196 {
return fRedux.at(dim); }
199 {
fRedux.at(dim) = redux; }
201 for (
unsigned int dim = 0; dim <
fNDim; dim++)
205 if (redux.size() !=
fNDim)
return;
220 unsigned int* cover =
new unsigned int[
fSize[dim]];
221 unsigned int* total =
new unsigned int[
fSize[dim]];
222 for (
unsigned int i = 0; i <
fSize[dim]; i++) {
227 value_t value[value_n];
228 for (
unsigned int linear_index = 0; linear_index <
fMaximumEntries; linear_index++) {
231 if (
Get(linear_index,value) && value[0] != 0)
236 for (
size_t i = 0; i <
fSize[dim]; i++) {
239 mycout <<
"bin " << i <<
", coord " << coord[dim] <<
": "
240 << double(cover[i]) / double(total[i]) * 100 <<
"%"<<
myendl;
257 if (value_n != 1)
return false;
258 if (
fNDim != 1)
return false;
259 return Set(&coord, &value);
263 if (value_n != 1)
return false;
264 return Set(coord, &value);
271 bool written =
false;
272 unsigned int linear_index;
273 for (
unsigned int dim = 0; dim <
fNDim; dim++) {
275 if (
fWrap[dim] == 0)
continue;
280 for (
size_t wrap = 0; wrap <
fWrap[dim]; wrap++) {
284 status &=
Set(linear_index, value);
288 status &=
Set(linear_index, value);
297 status &=
Set(linear_index, value);
303 bool Set(
const unsigned int linear_index,
const value_t& value) {
304 if (value_n != 1)
return false;
305 return Set(linear_index, &value);
308 bool Set(
const unsigned int linear_index,
const value_t* value) {
309 if (!
Check(linear_index))
return false;
310 for (
unsigned int i = 0; i < value_n; i++)
311 fValues[i][linear_index] = value[i];
317 bool Set(
const unsigned int* cell_index,
const value_t& value) {
318 if (value_n != 1)
return false;
319 return Set(cell_index, &value);
322 bool Set(
const unsigned int* cell_index,
const value_t* value) {
323 if (!
Check(cell_index))
return false;
324 return Set(
Index(cell_index),value);
333 if (value_n != 1)
return false;
334 if (
fNDim != 1)
return false;
336 if (
GetValue(&coord, &value))
return value;
341 if (value_n != 1)
return 0;
343 if (
GetValue(coord, &value))
return value;
348 if (value_n != 1)
return false;
353 for (
unsigned int i = 0; i < value_n; i++) value[i] = 0.0;
354 if (!
Check(coord))
return false;
355 value_t result[value_n];
359 if (!
Linear(coord, result))
return false;
364 default:
return false;
366 for (
unsigned int i = 0; i < value_n; i++) value[i] = result[i];
375 bool WriteText(std::ostream& stream)
const;
378 std::ofstream file(filename.c_str());
379 if (! file.is_open())
return false;
390 bool ReadText(std::istream& stream);
393 std::ifstream file(filename.c_str());
394 if (! file.is_open())
return false;
415 unsigned int Index(
const unsigned int* cell_index)
const;
417 unsigned int Index(
const unsigned int* cell_index,
const unsigned int offset)
const;
420 void Cell(
const coord_t coord,
unsigned int& cell_index,
double& cell_local,
const unsigned int dim)
const;
422 void Cell(
const coord_t* coord,
unsigned int* cell_index,
double* cell_local)
const;
424 void Cell(
const coord_t* coord,
unsigned int* cell_index)
const;
426 void Cell(
unsigned int linear_index,
unsigned int* cell_index)
const;
429 void Coord(
const unsigned int* cell_index,
coord_t* coord)
const;
431 void Coord(
const unsigned int linear_index,
coord_t* coord)
const;
440 for (
unsigned int dim = 0; dim <
fNDim; dim++)
452 bool Check(
const unsigned int* cell_index)
const;
454 bool Check(
const unsigned int linear_index)
const;
457 value_t
Get(
const unsigned int* cell_index)
const {
458 if (value_n != 1)
return 0;
463 value_t
Get(
const unsigned int index)
const {
467 bool Get(
const unsigned int index, value_t* value)
const {
468 for (
unsigned int i = 0; i < value_n; i++)
482template <
class value_t,
unsigned int value_n>
485 value_t* value)
const
490 for (
unsigned int i = 0; i < value_n; i++) value[i] = 0.0;
493 for (
unsigned int offset = 0; offset < (1U <<
fNDim); offset++) {
494 value_t neighbor[value_n];
498 for (
unsigned int dim = 0; dim <
fNDim; dim++)
501 for (
unsigned int i = 0; i < value_n; i++)
502 value[i] += fac * neighbor[i];
512template <
class value_t,
unsigned int value_n>
515 value_t* value)
const
527template <
class value_t,
unsigned int value_n>
530 for (
unsigned int dim = 0; dim <
fNDim; dim++)
531 if (
fWrap[dim] == 0 && (coord[dim] <
fMin[dim] || coord[dim] >
fMax[dim]))
542template <
class value_t,
unsigned int value_n>
545 for (
unsigned int dim = 0; dim <
fNDim; dim++)
546 if (cell_index[dim] >=
fSize[dim])
557template <
class value_t,
unsigned int value_n>
571template <
class value_t,
unsigned int value_n>
573 const unsigned int* cell_index)
const
575 unsigned int linear_index = 0;
577 for (
unsigned int dim = 0; dim <
fNDim; dim++)
579 linear_index +=
fExtent[dim] * cell_index[dim];
589template <
class value_t,
unsigned int value_n>
591 const unsigned int* cell_index,
592 const unsigned int pattern)
const
594 unsigned int linear_index = 0;
596 for (
unsigned int dim = 0; dim <
fNDim; dim++) {
598 unsigned int offset = (pattern >> dim) & 0x1;
599 linear_index +=
fExtent[dim] * (cell_index[dim] + offset);
611template <
class value_t,
unsigned int value_n>
614 unsigned int &cell_index,
616 const unsigned int dim)
const
619 double norm_coord = (coord -
fMin[dim]) /
fStep[dim];
621 double frac_part, int_part;
622 frac_part = modf(norm_coord, &int_part);
623 cell_local = frac_part;
624 cell_index =
static_cast<int>(int_part);
635template <
class value_t,
unsigned int value_n>
638 unsigned int* cell_index,
639 double* cell_local)
const
642 for (
unsigned int dim = 0; dim <
fNDim; dim++)
643 Cell(coord[dim], cell_index[dim], cell_local[dim], dim);
651template <
class value_t,
unsigned int value_n>
654 unsigned int* cell_index)
const
665template <
class value_t,
unsigned int value_n>
667 unsigned int linear_index,
668 unsigned int* cell_index)
const
671 for (
int dim =
fNDim-1; dim >= 0; dim--) {
672 cell_index[dim] = linear_index /
fExtent[dim];
673 linear_index -= cell_index[dim] *
fExtent[dim];
682template <
class value_t,
unsigned int value_n>
684 const unsigned int* cell_index,
687 for (
unsigned int dim = 0; dim <
fNDim; dim++)
688 coord[dim] =
fMin[dim] + cell_index[dim] *
fStep[dim];
696template <
class value_t,
unsigned int value_n>
698 const unsigned int linear_index,
709template <
class value_t,
unsigned int value_n>
713 mycout <<
"Writing text stream: ";
715 stream <<
fNDim <<
"\t" << value_n << std::endl;
717 for (
unsigned int dim = 0; dim <
fNDim; dim++)
718 stream << dim <<
"\t" <<
fMin[dim] <<
"\t" <<
fMax[dim] <<
"\t" <<
fStep[dim] << std::endl;
720 stream <<
fValues[0].size() << std::endl;
721 unsigned int entries =
fValues[0].size();
722 for (
unsigned int index = 0; index < entries; index++) {
724 for (
unsigned int i = 0; i < value_n; i++) {
725 stream <<
fValues[i][index] <<
"\t";
729 if (index % (entries / 10) == 0)
731 if (index % (entries / 10) != 0
732 && index % (entries / 40) == 0)
735 stream <<
"end" << std::endl;
745template <
class value_t,
unsigned int value_n>
749 mycout <<
"Reading text stream: ";
752 stream >>
fNDim >> n;
753 if (n != value_n)
return false;
756 for (
unsigned int dim = 0; dim <
fNDim; dim++)
760 unsigned int entries;
762 for (
unsigned int index = 0; index < entries; index++) {
764 for (
unsigned int i = 0; i < value_n; i++) {
768 if (index % (entries / 10) == 0)
770 if (index % (entries / 10) != 0
771 && index % (entries / 40) == 0)
779 if (end ==
"end")
return true;
794template <
class value_t,
unsigned int value_n>
797 std::ofstream file(filename.c_str(), std::ios::binary);
798 if (! file.is_open())
return false;
800 mycout <<
"Writing binary file: ";
802 uint32_t n = value_n;
803 uint32_t size =
sizeof(value_t);
804 file.write(
reinterpret_cast<const char*
>(&n),
sizeof(n));
805 file.write(
reinterpret_cast<const char*
>(&size),
sizeof(size));
807 uint32_t ndim =
fNDim;
808 file.write(
reinterpret_cast<const char*
>(&ndim),
sizeof(ndim));
809 for (
unsigned int dim = 0; dim <
fNDim; dim++) {
810 file.write(
reinterpret_cast<const char*
>(&
fMin[dim]),
sizeof(
fMin[dim]));
811 file.write(
reinterpret_cast<const char*
>(&
fMax[dim]),
sizeof(
fMax[dim]));
812 file.write(
reinterpret_cast<const char*
>(&
fStep[dim]),
sizeof(
fStep[dim]));
814 uint32_t entries =
fValues[0].size();
815 file.write(
reinterpret_cast<const char*
>(&entries),
sizeof(entries));
816 for (
unsigned int index = 0; index < entries; index++) {
818 for (
unsigned int i = 0; i < value_n; i++) {
819 value_t value =
fValues[i][index];
820 file.write(
reinterpret_cast<const char*
>(&value),
sizeof(value));
823 if (index % (entries / 10) == 0)
825 if (index % (entries / 10) != 0
826 && index % (entries / 40) == 0)
839template <
class value_t,
unsigned int value_n>
842 std::ifstream file(filename.c_str(), std::ios::binary);
843 if (! file.is_open())
return false;
845 mycout <<
"Reading binary file: ";
847 file.seekg(0, std::ios::end);
849 file.seekg(0, std::ios::beg);
851 file.read(
reinterpret_cast<char*
>(&n),
sizeof(n));
852 file.read(
reinterpret_cast<char*
>(&size),
sizeof(size));
853 if (n != value_n)
return false;
854 if (size !=
sizeof(value_t))
return false;
857 file.read(
reinterpret_cast<char*
>(&ndim),
sizeof(ndim));
860 for (
unsigned int dim = 0; dim <
fNDim; dim++) {
861 file.read(
reinterpret_cast<char*
>(&
fMin[dim]),
sizeof(
fMin[dim]));
862 file.read(
reinterpret_cast<char*
>(&
fMax[dim]),
sizeof(
fMax[dim]));
863 file.read(
reinterpret_cast<char*
>(&
fStep[dim]),
sizeof(
fStep[dim]));
867 uint32_t maximum_entries;
868 file.read(
reinterpret_cast<char*
>(&maximum_entries),
sizeof(maximum_entries));
870 int value_size =
sizeof(value_t);
871 unsigned int entries =
fValues[0].size();
872 for (
unsigned int index = 0; index < entries; index++) {
874 for (
unsigned int i = 0; i < value_n; i++) {
875 file.read((
char*)(&
fValues[i][index]),value_size);
878 if (index % (entries / 10) == 0)
880 if (index % (entries / 10) != 0
881 && index % (entries / 40) == 0)
A logfile class, based on an identical class in the Hermes analyzer.
EQwInterpolationMethod
Allowed interpolation methods.
@ kInterpolationMethodUnknown
void SetMinimumMaximumStep(const std::vector< coord_t > &min, const std::vector< coord_t > &max, const std::vector< coord_t > &step)
Set minimum, maximum, and step size to different values.
virtual ~QwInterpolator()
Destructor.
std::vector< size_t > fWrap
Wrap around this coordinate.
int GetDataReductionFactor(const unsigned int dim) const
Get data reduction factor.
bool Set(const unsigned int *cell_index, const value_t *value)
Set a set of values at a grid point (false if out of bounds)
unsigned int * f_cell_index
Pre-allocated cell index.
bool Set(const unsigned int linear_index, const value_t &value)
Set a single value at a linearized index (false if not possible)
void SetDimensions(const unsigned int ndim)
Set the number of coordinate dimensions and resize vectors.
coord_t GetMinimum(const unsigned int dim) const
Get minimum in dimension.
std::vector< size_t > fExtent
void SetDataReductionFactor(const unsigned int dim, const unsigned int redux)
Set data reduction factor.
void SetDataReductionFactor(const unsigned int redux)
bool WriteTextScreen() const
Write the grid to screen.
std::vector< value_t > fValues[value_n]
Table with pointers to arrays of values.
bool GetValue(const coord_t *coord, double *value) const
Get the interpolated value at coordinate (zero if out of bounds)
void PrintCoverage(const unsigned int dim)
Print coverage map for all bins in one dimension.
bool Set(const unsigned int *cell_index, const value_t &value)
Set a single value at a grid point (false if out of bounds)
bool ReadBinaryFile(const std::string &filename)
Read the grid values from binary file.
unsigned int GetWrapCoordinate(const unsigned int dim) const
Get wrapping coordinate.
QwInterpolator(const std::string &filename)
Constructor with file name.
bool Set(const coord_t *coord, const value_t &value)
Set a single value at a coordinate (false if not possible)
coord_t GetStepSize(const unsigned int dim) const
Get minimum in dimension.
void Coord(const unsigned int *cell_index, coord_t *coord) const
Return the coordinates based on the cell index (unchecked)
bool Get(const unsigned int index, value_t *value) const
Get a vector value by linearized index (unchecked)
void SetWrapCoordinate(const unsigned int dim, const size_t wrap=1)
Set wrapping coordinate.
EQwInterpolationMethod fInterpolationMethod
Interpolation method.
unsigned int GetCurrentEntries() const
Get the current number of entries.
unsigned int GetMaximumEntries() const
Get the maximum number of entries.
std::vector< coord_t > fMin
Minimum in each dimension.
unsigned int fNDim
Number of dimensions in coordinates.
QwInterpolator(const unsigned int ndim=1)
Constructor with number of dimensions.
bool Set(const coord_t &coord, const value_t &value)
Set a single value at a coordinate (false if not possible)
bool WriteBinaryFile(const std::string &filename) const
Write the grid values to binary file.
unsigned int fMaximumEntries
Maximum number of values.
unsigned int fCurrentEntries
Number of values read in.
std::vector< coord_t > fStep
Step size in each dimension.
bool Check(const coord_t *coord) const
Check for boundaries with coordinate argument.
std::vector< size_t > fSize
Number of points in each dimension.
void SetMinimumMaximumStep(const coord_t min, const coord_t max, const coord_t step)
Set minimum, maximum, and step size to single values.
value_t Get(const unsigned int *cell_index) const
Get a single value by cell index (unchecked)
void SetInterpolationMethod(const EQwInterpolationMethod method)
Set the interpolation method.
value_t Get(const unsigned int index) const
Get a single value by linearized index (unchecked)
coord_t GetMaximum(const unsigned int dim) const
Get maximum in dimension.
EQwInterpolationMethod GetInterpolationMethod() const
Get the interpolation method.
unsigned int Index(const coord_t *coord) const
Return the linearized index based on the point coordinates (unchecked)
bool InBounds(const coord_t *coord) const
Return true if the coordinate is in bounds.
bool WriteTextFile(const std::string &filename) const
Write the grid to text file.
bool ReadText(std::istream &stream)
Read the grid from text.
void Cell(const coord_t coord, unsigned int &cell_index, double &cell_local, const unsigned int dim) const
Return the cell index and local coordinates in one dimension (unchecked)
void Nearest(const coord_t *coord, unsigned int *cell_index) const
Return the cell index closest to the coordinate (could be above) (unchecked)
std::vector< size_t > fRedux
Data reduction factor.
bool Set(const coord_t *coord, const value_t *value)
Set a set of values at a coordinate (false if not possible)
bool Set(const unsigned int linear_index, const value_t *value)
Set a set of values at a linearized index (false if not possible)
bool WriteText(std::ostream &stream) const
Write the grid as text.
bool ReadTextFile(const std::string &filename)
Read the grid from text file.
void SetWrapCoordinate(const std::vector< size_t > &wrap)
void SetDataReductionFactor(const std::vector< unsigned int > &redux)
value_t GetValue(const coord_t &coord) const
Get the interpolated value at coordinate (zero if out of bounds)
bool Linear(const coord_t *coord, value_t *value) const
Linear interpolation (unchecked)
bool NearestNeighbor(const coord_t *coord, value_t *value) const
Nearest-neighbor (unchecked)
std::vector< coord_t > fMax
Maximum in each dimension.
bool GetValue(const coord_t *coord, double &value) const
Get the interpolated value at coordinate (zero if out of bounds)
value_t GetValue(const coord_t *coord) const
Get the interpolated value at coordinate (zero if out of bounds)
QwInterpolator(const std::vector< coord_t > &min, const std::vector< coord_t > &max, const std::vector< coord_t > &step)
Constructor with minimum, maximum, and step size.
double * f_cell_local
Pre-allocated local coordinates.
static std::ostream & flush(std::ostream &)
Flush the streams.