/*
Copyright 2013--2018 James E. McClure, Virginia Polytechnic & State University
Copyright Equnior ASA
This file is part of the Open Porous Media project (OPM).
OPM 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.
OPM 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 OPM. If not, see .
*/
/*
Copyright 2013--2018 James E. McClure, Virginia Polytechnic & State University
Copyright Equnior ASA
This file is part of the Open Porous Media project (OPM).
OPM 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.
OPM 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 OPM. If not, see .
*/
#ifndef included_ArrayClass_hpp
#define included_ArrayClass_hpp
#include "common/Array.h"
#include "common/FunctionTable.h"
#include "common/FunctionTable.hpp"
#include "common/Utilities.h"
#include
#include
#include
#include
#include
#include
/********************************************************
* External instantiations *
********************************************************/
extern template class Array;
extern template class Array;
extern template class Array;
extern template class Array;
extern template class Array;
extern template class Array;
extern template class Array;
extern template class Array;
extern template class Array;
extern template class Array;
extern template class Array;
/********************************************************
* Macros to help instantiate functions *
********************************************************/
// clang-format off
#define instantiateArrayConstructors( TYPE ) \
template Array::Array(); \
template Array::~Array(); \
template Array::Array( const ArraySize & ); \
template Array::Array( size_t ); \
template Array::Array( size_t, size_t ); \
template Array::Array( size_t, size_t, size_t ); \
template Array::Array( size_t, size_t, size_t, size_t ); \
template Array::Array( size_t, size_t, size_t, size_t, size_t ); \
template Array::Array( const std::vector &, const TYPE * ); \
template Array::Array( std::initializer_list ); \
template Array::Array( std::initializer_list> ); \
template Array::Array( const Array & ); \
template Array::Array( Array && ); \
template void Array::reshape( ArraySize const& ); \
template void Array::squeeze(); \
template std::unique_ptr> \
Array::constView(ArraySize const&, std::shared_ptr const&); \
template void Array::viewRaw( ArraySize const&, TYPE*, bool, bool ); \
template void Array::view2(ArraySize const&, std::shared_ptr ); \
template Array &Array::operator=( const Array & ); \
template Array &Array::operator=( Array && );
// clang-format on
/********************************************************
* Constructors *
********************************************************/
template
Array::Array()
: d_isCopyable(true), d_isFixedSize(false), d_data(nullptr) {}
template
Array::Array(const ArraySize &N)
: d_isCopyable(true), d_isFixedSize(false) {
allocate(N);
}
template
Array::Array(size_t N)
: d_isCopyable(true), d_isFixedSize(false) {
allocate(ArraySize(N));
}
template
Array::Array(size_t N_rows, size_t N_cols)
: d_isCopyable(true), d_isFixedSize(false) {
allocate(ArraySize(N_rows, N_cols));
}
template
Array::Array(size_t N1, size_t N2, size_t N3)
: d_isCopyable(true), d_isFixedSize(false) {
allocate(ArraySize(N1, N2, N3));
}
template
Array::Array(size_t N1, size_t N2, size_t N3, size_t N4)
: d_isCopyable(true), d_isFixedSize(false) {
allocate(ArraySize(N1, N2, N3, N4));
}
template
Array::Array(size_t N1, size_t N2, size_t N3, size_t N4,
size_t N5)
: d_isCopyable(true), d_isFixedSize(false) {
allocate(ArraySize(N1, N2, N3, N4, N5));
}
template
Array::Array(const std::vector &N,
const TYPE *data)
: d_isCopyable(true), d_isFixedSize(false) {
allocate(N);
if (data)
copy(data);
}
template
Array::Array(std::string str)
: d_isCopyable(true), d_isFixedSize(false) {
allocate(0);
if ((int)std::count(str.begin(), str.end(), ' ') == (int)str.length()) {
// Empty string
return;
}
// Remove unnecessary whitespace
while (str.front() == ' ')
str.erase(0, 1);
while (str.back() == ' ')
str.resize(str.length() - 1);
while (str.find(',') != std::string::npos)
str[str.find(',')] = ' ';
while (str.find(" ") != std::string::npos)
str.replace(str.find(" "), 2, " ");
// Check if the string is of the format [...]
if (str.front() == '[' && str.back() == ']') {
str.erase(0, 1);
str.resize(str.length() - 1);
*this = Array(str);
return;
}
// Check if we are dealing with a 2D array
if (str.find(';') != std::string::npos) {
size_t i1 = 0;
size_t i2 = str.find(';');
std::vector x;
while (i2 > i1) {
auto tmp = str.substr(i1, i2 - i1);
x.emplace_back(Array(tmp));
i1 = i2 + 1;
i2 = str.find(';', i1 + 1);
if (i2 == std::string::npos)
i2 = str.length();
}
for (auto &y : x)
y.reshape({1, y.length()});
*this = cat(x, 0);
return;
}
// Begin parsing the array constructor
size_t i1 = 0;
size_t i2 = str.find(' ');
std::vector data;
while (i2 > i1) {
auto tmp = str.substr(i1, i2 - i1);
int type = std::count(tmp.begin(), tmp.end(), ':');
if (type == 0) {
data.push_back(std::stod(tmp));
} else if (type == 1) {
size_t k = tmp.find(':');
TYPE x1 = std::stod(tmp.substr(0, k));
TYPE x2 = std::stod(tmp.substr(k + 1));
double tol = 1e-8 * (x2 - x1);
for (TYPE x = x1; x <= x2 + tol; x += 1)
data.push_back(x);
} else if (type == 2) {
size_t k1 = tmp.find(':');
size_t k2 = tmp.find(':', k1 + 1);
TYPE x1 = std::stod(tmp.substr(0, k1));
TYPE dx = std::stod(tmp.substr(k1 + 1, k2 - k1 - 1));
TYPE x2 = std::stod(tmp.substr(k2 + 1));
double tol = 1e-8 * (x2 - x1);
for (TYPE x = x1; x <= x2 + tol; x += dx)
data.push_back(x);
} else {
throw std::logic_error("Failed to parse string constructor: " +
str);
}
i1 = i2;
i2 = str.find(' ', i1 + 1);
if (i2 == std::string::npos)
i2 = str.length();
}
allocate(data.size());
if (std::is_same::value) {
for (size_t i = 0; i < data.size(); i++)
d_data[i] = data[i];
} else {
copy(data.data());
}
}
template
Array::Array(std::initializer_list x)
: d_isCopyable(true), d_isFixedSize(false) {
allocate({x.size()});
auto it = x.begin();
for (size_t i = 0; i < x.size(); ++i, ++it)
d_data[i] = *it;
}
template
Array::Array(
std::initializer_list> x)
: d_isCopyable(true), d_isFixedSize(false) {
size_t Nx = x.size();
size_t Ny = 0;
for (const auto y : x)
Ny = std::max(Ny, y.size());
allocate({Nx, Ny});
auto itx = x.begin();
for (size_t i = 0; i < x.size(); ++i, ++itx) {
auto ity = itx->begin();
for (size_t j = 0; j < itx->size(); ++j, ++ity) {
d_data[i + j * Nx] = *ity;
}
}
}
template
void Array::allocate(const ArraySize &N) {
if (d_isFixedSize)
throw std::logic_error("Array cannot be resized");
d_size = N;
auto length = d_size.length();
d_data = nullptr;
if (length > 0) {
try {
d_data = new TYPE[length];
} catch (...) {
throw std::logic_error("Failed to allocate array");
}
}
d_ptr.reset(d_data, [](TYPE *p) { delete[] p; });
}
template
Array::Array(const Array &rhs)
: d_isCopyable(true), d_isFixedSize(false) {
if (!rhs.d_isCopyable)
throw std::logic_error("Array cannot be copied");
allocate(rhs.size());
copy(rhs.d_data);
}
template
Array::Array(Array &&rhs)
: d_isCopyable(rhs.d_isCopyable), d_isFixedSize(rhs.d_isFixedSize),
d_size(rhs.d_size), d_data(rhs.d_data), d_ptr(std::move(rhs.d_ptr)) {
rhs.d_size = ArraySize();
rhs.d_data = nullptr;
rhs.d_ptr = nullptr;
}
template
Array &Array::
operator=(const Array &rhs) {
if (this == &rhs)
return *this;
if (!rhs.d_isCopyable)
throw std::logic_error("Array cannot be copied");
allocate(rhs.size());
copy(rhs.d_data);
return *this;
}
template
Array &Array::
operator=(Array &&rhs) {
if (this == &rhs)
return *this;
d_isCopyable = rhs.d_isCopyable;
d_isFixedSize = rhs.d_isFixedSize;
d_size = rhs.d_size;
d_data = rhs.d_data;
d_ptr = std::move(rhs.d_ptr);
rhs.d_size = ArraySize();
rhs.d_data = nullptr;
rhs.d_ptr = nullptr;
return *this;
}
template
Array &Array::
operator=(const std::vector &rhs) {
allocate(ArraySize(rhs.size()));
for (size_t i = 0; i < rhs.size(); i++)
d_data[i] = rhs[i];
return *this;
}
template
Array::~Array() {}
template
void Array::clear() {
d_isCopyable = true;
d_isFixedSize = false;
d_size = ArraySize();
d_ptr.reset();
d_data = nullptr;
}
/********************************************************
* Copy/move values from one array to another (resize) *
********************************************************/
template
static inline void moveValues(const ArraySize &N1, const ArraySize &N2,
TYPE *data1, TYPE *data2) {
for (size_t i5 = 0; i5 < std::min(N1[4], N2[4]); i5++) {
for (size_t i4 = 0; i4 < std::min(N1[3], N2[3]); i4++) {
for (size_t i3 = 0; i3 < std::min(N1[2], N2[2]); i3++) {
for (size_t i2 = 0; i2 < std::min(N1[1], N2[1]); i2++) {
for (size_t i1 = 0; i1 < std::min(N1[0], N2[0]); i1++) {
size_t index1 = N1.index(i1, i2, i3, i4, i5);
size_t index2 = N2.index(i1, i2, i3, i4, i5);
data2[index2] = std::move(data1[index1]);
}
}
}
}
}
}
template
static inline void copyValues(const ArraySize &N1, const ArraySize &N2,
const TYPE *data1, TYPE *data2) {
for (size_t i5 = 0; i5 < std::min(N1[4], N2[4]); i5++) {
for (size_t i4 = 0; i4 < std::min(N1[3], N2[3]); i4++) {
for (size_t i3 = 0; i3 < std::min(N1[2], N2[2]); i3++) {
for (size_t i2 = 0; i2 < std::min(N1[1], N2[1]); i2++) {
for (size_t i1 = 0; i1 < std::min(N1[0], N2[0]); i1++) {
size_t index1 = N1.index(i1, i2, i3, i4, i5);
size_t index2 = N2.index(i1, i2, i3, i4, i5);
data2[index2] = data1[index1];
}
}
}
}
}
}
/********************************************************
* Resize the array *
********************************************************/
template
void Array::resize(const ArraySize &N) {
// Check if the array actually changed size
bool equal = true;
for (size_t i = 0; i < ArraySize::maxDim(); i++)
equal = equal && N[i] == d_size[i];
if (equal) {
d_size = N;
return;
}
// Store the old data
auto N0 = d_size;
auto data0 = d_ptr;
// Allocate new data
allocate(N);
// Copy the old values
if (N.length() > 0 && d_size.length() > 0) {
if (data0.use_count() <= 1) {
// We own the data, use std:move
moveValues(N0, N, data0.get(), d_data);
} else if (std::is_copy_constructible::value) {
// We do not own the data, copy
copyValues(N0, N, data0.get(), d_data);
} else {
throw std::logic_error("No copy constructor");
}
}
}
template
void Array::resizeDim(int dim, size_t N,
const TYPE &value) {
if (dim < 0 || dim > d_size.ndim())
throw std::out_of_range("Invalid dimension");
size_t N0 = d_size[dim];
auto size = d_size;
size.resize(dim, N);
resize(size);
size_t n1 = 1, n2 = 1;
for (int d = 0; d < dim; d++)
n1 *= size[d];
for (size_t d = dim + 1; d < size.ndim(); d++)
n2 *= size[d];
for (size_t k = 0; k < n2; k++) {
for (size_t j = N0; j < N; j++) {
for (size_t i = 0; i < n1; i++) {
d_data[i + j * n1 + k * n1 * N] = value;
}
}
}
}
/********************************************************
* Reshape/squeeze the array *
********************************************************/
template
void Array::reshape(const ArraySize &N) {
if (N.length() != d_size.length())
throw std::logic_error(
"reshape is not allowed to change the array size");
d_size = N;
}
template
void Array::squeeze() {
d_size.squeeze();
}
/********************************************************
* Shift/permute the array *
********************************************************/
template
Array Array::shiftDim(int N) const {
if (N > 0)
N = N % d_size.ndim();
if (N == 0) {
// No shift required
return *this;
} else if (N > 0) {
// Shift to the left and wrap
std::vector index(d_size.ndim());
size_t i = 0;
for (size_t j = N; j < index.size(); j++, i++)
index[i] = j;
for (size_t j = 0; i < index.size(); j++, i++)
index[i] = j;
return permute(index);
} else {
// Shift to the right (padding with singletons)
N = -N;
ASSERT(d_size.ndim() + N < (int)ArraySize::maxDim());
size_t dims[10] = {1};
for (int i = 0; i < ndim(); i++)
dims[N + i] = d_size[i];
auto y = *this;
y.reshape(ArraySize(ndim() + N, dims));
return y;
}
}
template
Array
Array::permute(const std::vector &index) const {
// Check the permutation
ASSERT((int)index.size() == ndim());
for (int i = 0; i < ndim(); i++) {
ASSERT(index[i] < ndim());
for (int j = 0; j < i; j++)
ASSERT(index[i] != index[j]);
}
// Create the new Array
size_t dims[5] = {1u, 1u, 1u, 1u, 1u};
for (size_t i = 0; i < index.size(); i++)
dims[i] = d_size[index[i]];
Array y(ArraySize(ndim(), dims));
y.fill(-1);
ASSERT(y.length() == this->length());
// Fill the data
size_t N[5] = {1u, 1u, 1u, 1u, 1u};
for (int i = 0; i < ndim(); i++) {
std::array ijk = {0, 0, 0, 0, 0};
ijk[index[i]] = 1;
N[i] = d_size.index(ijk);
}
size_t tmp = (dims[0] - 1) * N[0] + (dims[1] - 1) * N[1] +
(dims[2] - 1) * N[2] + (dims[3] - 1) * N[3] +
(dims[4] - 1) * N[4] + 1;
ASSERT(tmp == length());
for (size_t i4 = 0; i4 < dims[4]; i4++) {
for (size_t i3 = 0; i3 < dims[3]; i3++) {
for (size_t i2 = 0; i2 < dims[2]; i2++) {
for (size_t i1 = 0; i1 < dims[1]; i1++) {
for (size_t i0 = 0; i0 < dims[0]; i0++) {
size_t index2 = i0 * N[0] + i1 * N[1] + i2 * N[2] +
i3 * N[3] + i4 * N[4];
y(i0, i1, i2, i3, i4) = d_data[index2];
}
}
}
}
}
return y;
}
/********************************************************
* Subset the array *
********************************************************/
// Helper function to check subset indices
template
inline void Array::checkSubsetIndex(
const std::vector> &range) const {
bool test = (int)range.size() == d_size.ndim();
for (size_t d = 0; d < range.size(); d++)
test = test && range[d].j <= d_size[d];
if (!test)
throw std::logic_error("indices for subset are invalid");
}
template
std::vector>
Array::convert(const std::vector &index) const {
std::vector> range(d_size.ndim());
if (index.size() % 2 != 0 ||
static_cast(index.size() / 2) < d_size.ndim())
throw std::logic_error("indices for subset are invalid");
for (int d = 0; d < d_size.ndim(); d++)
range[d] = Range(index[2 * d + 0], index[2 * d + 1]);
return range;
}
// Helper function to return dimensions for the subset array
template
void Array::getSubsetArrays(
const std::vector> &index, std::array &first,
std::array &last, std::array &inc,
std::array &N) {
first.fill(0);
last.fill(0);
inc.fill(1);
N.fill(1);
size_t ndim = index.size();
for (size_t d = 0; d < ndim; d++) {
first[d] = index[d].i;
last[d] = index[d].j;
inc[d] = index[d].k;
N[d] = (last[d] - first[d] + inc[d]) / inc[d];
}
}
template
Array Array::subset(
const std::vector> &index) const {
// Get the subset indicies
checkSubsetIndex(index);
std::array first, last, inc, N1;
getSubsetArrays(index, first, last, inc, N1);
ArraySize S1(d_size.ndim(), N1.data());
// Create the new array
Array subset_array(S1);
// Fill the new array
TYPE *subset_data = subset_array.data();
for (size_t i4 = first[4], k1 = 0; i4 <= last[4]; i4 += inc[4]) {
for (size_t i3 = first[3]; i3 <= last[3]; i3 += inc[3]) {
for (size_t i2 = first[2]; i2 <= last[2]; i2 += inc[2]) {
for (size_t i1 = first[1]; i1 <= last[1]; i1 += inc[1]) {
for (size_t i0 = first[0]; i0 <= last[0];
i0 += inc[0], k1++) {
size_t k2 = d_size.index(i0, i1, i2, i3, i4);
subset_data[k1] = d_data[k2];
}
}
}
}
}
return subset_array;
}
template
Array
Array::subset(const std::vector &index) const {
auto range = convert(index);
return subset(range);
}
template
void Array::copySubset(
const std::vector> &index,
const Array &subset) {
// Get the subset indices
checkSubsetIndex(index);
std::array first, last, inc, N1;
getSubsetArrays(index, first, last, inc, N1);
// Copy the sub-array
const TYPE *src_data = subset.data();
for (size_t i4 = first[4], k1 = 0; i4 <= last[4]; i4 += inc[4]) {
for (size_t i3 = first[3]; i3 <= last[3]; i3 += inc[3]) {
for (size_t i2 = first[2]; i2 <= last[2]; i2 += inc[2]) {
for (size_t i1 = first[1]; i1 <= last[1]; i1 += inc[1]) {
for (size_t i0 = first[0]; i0 <= last[0];
i0 += inc[0], k1++) {
size_t k2 = d_size.index(i0, i1, i2, i3, i4);
d_data[k2] = src_data[k1];
}
}
}
}
}
}
template
void Array::addSubset(
const std::vector> &index,
const Array &subset) {
// Get the subset indices
checkSubsetIndex(index);
std::array first, last, inc, N1;
getSubsetArrays(index, first, last, inc, N1);
// add the sub-array
for (size_t i4 = first[4], k1 = 0; i4 <= last[4]; i4 += inc[4]) {
for (size_t i3 = first[3]; i3 <= last[3]; i3 += inc[3]) {
for (size_t i2 = first[2]; i2 <= last[2]; i2 += inc[2]) {
for (size_t i1 = first[1]; i1 <= last[1]; i1 += inc[1]) {
for (size_t i0 = first[0]; i0 <= last[0];
i0 += inc[0], k1++) {
size_t k2 = d_size.index(i0, i1, i2, i3, i4);
d_data[k2] += subset.d_data[k1];
}
}
}
}
}
}
template
void Array::copySubset(
const std::vector &index,
const Array &subset) {
auto range = convert(index);
copySubset(range, subset);
}
template
void Array::addSubset(
const std::vector &index,
const Array &subset) {
auto range = convert(index);
addSubset(range, subset);
}
/********************************************************
* Operator overloading *
********************************************************/
template
bool Array::operator==(const Array &rhs) const {
if (this == &rhs)
return true;
if (d_size != rhs.d_size)
return false;
bool match = true;
for (size_t i = 0; i < d_size.length(); i++)
match = match && d_data[i] == rhs.d_data[i];
return match;
}
/********************************************************
* Get a view of an C array *
********************************************************/
template
std::unique_ptr>
Array::view(const ArraySize &N,
std::shared_ptr data) {
auto array = std::make_unique>();
array->d_size = N;
array->d_ptr = data;
array->d_data = array->d_ptr.get();
return array;
}
template
std::unique_ptr>
Array::constView(
const ArraySize &N, std::shared_ptr const &data) {
auto array = std::make_unique>();
array->d_size = N;
array->d_ptr = std::const_pointer_cast(data);
array->d_data = array->d_ptr.get();
return array;
}
template
void Array::view2(Array &src) {
view2(src.size(), src.getPtr());
d_data = src.d_data;
}
template
void Array::view2(const ArraySize &N,
std::shared_ptr data) {
d_size = N;
d_ptr = data;
d_data = d_ptr.get();
}
template
void Array::viewRaw(const ArraySize &N, TYPE *data,
bool isCopyable, bool isFixedSize) {
d_isCopyable = isCopyable;
d_isFixedSize = isFixedSize;
d_ptr.reset();
d_size = N;
d_data = data;
}
/********************************************************
* Basic functions *
********************************************************/
template
void Array::swap(Array &other) {
// check that dimensions match
if (d_size != other.d_size)
throw std::logic_error("length of arrays do not match");
// swap the data
std::swap(d_data, other.d_data);
std::swap(d_ptr, other.d_ptr);
}
template
void Array::pow(
const Array &baseArray, const TYPE &exp) {
// not insisting on the shapes being the same
// but insisting on the total size being the same
if (d_size.length() != baseArray.length())
throw std::logic_error("length of arrays do not match");
const auto base_data = baseArray.data();
for (size_t i = 0; i < d_size.length(); i++)
d_data[i] = std::pow(base_data[i], exp);
}
/********************************************************
* Replicate the array *
********************************************************/
template
Array
Array::repmat(const std::vector &N_rep) const {
std::vector N2(d_size.begin(), d_size.end());
if (N2.size() < N_rep.size())
N2.resize(N_rep.size(), 1);
std::array N1, Nr;
N1.fill(1);
Nr.fill(1);
for (size_t d = 0; d < N_rep.size(); d++) {
N1[d] = d_size[d];
Nr[d] = N_rep[d];
N2[d] *= N_rep[d];
}
Array y(N2);
TYPE *y2 = y.data();
for (size_t i4 = 0, index = 0; i4 < N1[4]; i4++) {
for (size_t j4 = 0; j4 < Nr[4]; j4++) {
for (size_t i3 = 0; i3 < N1[3]; i3++) {
for (size_t j3 = 0; j3 < Nr[3]; j3++) {
for (size_t i2 = 0; i2 < N1[2]; i2++) {
for (size_t j2 = 0; j2 < Nr[2]; j2++) {
for (size_t i1 = 0; i1 < N1[1]; i1++) {
for (size_t j1 = 0; j1 < Nr[1]; j1++) {
for (size_t i0 = 0; i0 < N1[0]; i0++) {
size_t k =
d_size.index(i0, i1, i2, i3, i4);
TYPE x = d_data[k];
for (size_t j0 = 0; j0 < Nr[0];
j0++, index++)
y2[index] = x;
}
}
}
}
}
}
}
}
}
return y;
}
/********************************************************
* Simple math operations *
********************************************************/
template
bool Array::NaNs() const {
bool test = false;
for (size_t i = 0; i < d_size.length(); i++)
test = test || d_data[i] != d_data[i];
return test;
}
template
TYPE Array::mean(void) const {
TYPE x = sum() / d_size.length();
return x;
}
template
Array Array::min(int dir) const {
auto size_ans = d_size;
size_ans.resize(dir, 1);
Array ans(size_ans);
size_t N1 = 1, N2 = 1, N3 = 1;
for (int d = 0; d < std::min(dir, d_size.ndim()); d++)
N1 *= d_size[d];
N2 = d_size[dir];
for (size_t d = dir + 1; d < d_size.ndim(); d++)
N3 *= d_size[d];
TYPE *data2 = ans.d_data;
for (size_t i3 = 0; i3 < N3; i3++) {
for (size_t i1 = 0; i1 < N1; i1++) {
TYPE x = d_data[i1 + i3 * N1 * N2];
for (size_t i2 = 0; i2 < N2; i2++)
x = std::min(x, d_data[i1 + i2 * N1 + i3 * N1 * N2]);
data2[i1 + i3 * N1] = x;
}
}
return ans;
}
template
Array Array::max(int dir) const {
auto size_ans = d_size;
size_ans.resize(dir, 1);
Array ans(size_ans);
size_t N1 = 1, N2 = 1, N3 = 1;
for (int d = 0; d < std::min(dir, d_size.ndim()); d++)
N1 *= d_size[d];
N2 = d_size[dir];
DISABLE_WARNINGS // Suppress false array subscript is above array bounds
for (size_t d = dir + 1; d < d_size.ndim(); d++) N3 *= d_size[d];
ENABLE_WARNINGS // Enable warnings
TYPE *data2 = ans.d_data;
for (size_t i3 = 0; i3 < N3; i3++) {
for (size_t i1 = 0; i1 < N1; i1++) {
TYPE x = d_data[i1 + i3 * N1 * N2];
for (size_t i2 = 0; i2 < N2; i2++)
x = std::max(x, d_data[i1 + i2 * N1 + i3 * N1 * N2]);
data2[i1 + i3 * N1] = x;
}
}
return ans;
}
template
Array Array::sum(int dir) const {
auto size_ans = d_size;
size_ans.resize(dir, 1);
Array ans(size_ans);
size_t N1 = 1, N2 = 1, N3 = 1;
for (int d = 0; d < std::min(dir, d_size.ndim()); d++)
N1 *= d_size[d];
N2 = d_size[dir];
DISABLE_WARNINGS
for (size_t d = dir + 1; d < d_size.ndim(); d++)
N3 *= d_size[d];
ENABLE_WARNINGS
TYPE *data2 = ans.d_data;
for (size_t i3 = 0; i3 < N3; i3++) {
for (size_t i1 = 0; i1 < N1; i1++) {
TYPE x = 0;
for (size_t i2 = 0; i2 < N2; i2++)
x += d_data[i1 + i2 * N1 + i3 * N1 * N2];
data2[i1 + i3 * N1] = x;
}
}
return ans;
}
template
TYPE Array::min(
const std::vector> &range) const {
// Get the subset indicies
checkSubsetIndex(range);
std::array first, last, inc, N1;
getSubsetArrays(range, first, last, inc, N1);
TYPE x = std::numeric_limits::max();
for (size_t i4 = first[4]; i4 <= last[4]; i4 += inc[4]) {
for (size_t i3 = first[3]; i3 <= last[3]; i3 += inc[3]) {
for (size_t i2 = first[2]; i2 <= last[2]; i2 += inc[2]) {
for (size_t i1 = first[1]; i1 <= last[1]; i1 += inc[1]) {
for (size_t i0 = first[0]; i0 <= last[0]; i0 += inc[0]) {
size_t k1 = d_size.index(i0, i1, i2, i3, i4);
x = std::min(x, d_data[k1]);
}
}
}
}
}
return x;
}
template
TYPE Array::max(
const std::vector> &range) const {
// Get the subset indicies
checkSubsetIndex(range);
std::array first, last, inc, N1;
getSubsetArrays(range, first, last, inc, N1);
TYPE x = std::numeric_limits::min();
for (size_t i4 = first[4]; i4 <= last[4]; i4 += inc[4]) {
for (size_t i3 = first[3]; i3 <= last[3]; i3 += inc[3]) {
for (size_t i2 = first[2]; i2 <= last[2]; i2 += inc[2]) {
for (size_t i1 = first[1]; i1 <= last[1]; i1 += inc[1]) {
for (size_t i0 = first[0]; i0 <= last[0]; i0 += inc[0]) {
size_t k1 = d_size.index(i0, i1, i2, i3, i4);
x = std::max(x, d_data[k1]);
}
}
}
}
}
return x;
}
template
TYPE Array::sum(
const std::vector> &range) const {
// Get the subset indicies
checkSubsetIndex(range);
std::array first, last, inc, N1;
getSubsetArrays(range, first, last, inc, N1);
TYPE x = 0;
for (size_t i4 = first[4]; i4 <= last[4]; i4 += inc[4]) {
for (size_t i3 = first[3]; i3 <= last[3]; i3 += inc[3]) {
for (size_t i2 = first[2]; i2 <= last[2]; i2 += inc[2]) {
for (size_t i1 = first[1]; i1 <= last[1]; i1 += inc[1]) {
for (size_t i0 = first[0]; i0 <= last[0]; i0 += inc[0]) {
size_t k1 = d_size.index(i0, i1, i2, i3, i4);
x += d_data[k1];
}
}
}
}
}
return x;
}
template
TYPE Array::mean(
const std::vector> &range) const {
// Get the subset indicies
checkSubsetIndex(range);
std::array first, last, inc, N1;
getSubsetArrays(range, first, last, inc, N1);
size_t n = 1;
for (auto &d : N1)
n *= d;
TYPE x = sum(range) / n;
return x;
}
template
TYPE Array::min(const std::vector &index) const {
auto range = convert(index);
return min(range);
}
template
TYPE Array::max(const std::vector &index) const {
auto range = convert(index);
return max(range);
}
template
TYPE Array::sum(const std::vector &index) const {
auto range = convert(index);
return sum(range);
}
template
TYPE Array::mean(const std::vector &index) const {
auto range = convert(index);
return mean(range);
}
/********************************************************
* Find all elements that match the given operation *
********************************************************/
template
std::vector Array::find(
const TYPE &value,
std::function compare) const {
std::vector result;
result.reserve(d_size.length());
for (size_t i = 0; i < d_size.length(); i++) {
if (compare(d_data[i], value))
result.push_back(i);
}
return result;
}
/********************************************************
* Print an array to an output stream *
********************************************************/
template
void Array::print(std::ostream &os,
const std::string &name,
const std::string &prefix) const {
if (d_size.ndim() == 1) {
for (size_t i = 0; i < d_size[0]; i++)
os << prefix << name << "[" << i << "] = " << d_data[i]
<< std::endl;
} else if (d_size.ndim() == 2) {
os << prefix << name << ":" << std::endl;
for (size_t i = 0; i < d_size[0]; i++) {
for (size_t j = 0; j < d_size[1]; j++)
os << prefix << " " << operator()(i, j);
os << std::endl;
}
} else {
throw std::logic_error("Not programmed for this dimension");
}
}
/********************************************************
* Reverse dimensions (transpose) *
********************************************************/
template
Array Array