diff --git a/.arclint b/.arclint
index 73c2265..cae0567 100644
--- a/.arclint
+++ b/.arclint
@@ -1,18 +1,24 @@
{
"linters": {
"python": {
"type": "flake8",
"include": [ "(\\.py$)", "(^hbm$)" ],
"exclude": ["(^third-party/)", "(legacy/)"],
"severity.rules": {
"(^E)": "error",
"(^W)": "warning",
"(^F)": "advice"
}
},
+ "format": {
+ "type": "clang-format",
+ "include": [ "(\\.cc$)", "(\\.hh$)" ],
+ "exclude": "(^third-party/)"
+ },
+
"merges": {
"type": "merge-conflict"
}
}
}
diff --git a/src/core/array.hh b/src/core/array.hh
index 7d48317..727e24a 100644
--- a/src/core/array.hh
+++ b/src/core/array.hh
@@ -1,160 +1,160 @@
/**
* @file
* @section LICENSE
*
* Copyright (©) 2016-19 EPFL (École Polytechnique Fédérale de Lausanne),
* Laboratory (LSMS - Laboratoire de Simulation en Mécanique des Solides)
*
* This program is free software: you can redistribute it and/or modify
* it under the terms of the GNU Affero General Public License as published
* by the Free Software Foundation, either version 3 of the License, or
* (at your option) any later version.
*
* This program 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 Affero General Public License for more details.
*
* You should have received a copy of the GNU Affero General Public License
* along with this program. If not, see .
*
*/
/* -------------------------------------------------------------------------- */
#ifndef __ARRAY__HH__
#define __ARRAY__HH__
/* -------------------------------------------------------------------------- */
#include "tamaas.hh"
#include
#include
#include
#include
/* -------------------------------------------------------------------------- */
namespace tamaas {
/**
* @brief Generic storage class with wrapping capacities
*/
template
class Array final {
public:
/// Default
Array() = default;
/// Empty array of given size
Array(UInt size) : Array() { this->resize(size); }
/// Copy constructor (deep)
Array(const Array& v) : Array() {
this->resize(v._size);
thrust::copy(v._data, v._data + v._size, this->_data);
}
/// Move constructor (transfers data ownership)
- Array(Array&& v) : Array() {
+ Array(Array&& v) noexcept : Array() {
_data = std::exchange(v._data, nullptr);
_size = std::exchange(v._size, 0);
wrapped = std::exchange(v.wrapped, false);
}
/// Wrap array on data
Array(T* data, UInt size) : Array() { this->wrapMemory(data, size); }
/// Destructor
~Array() {
- if (wrapped == false)
+ if (not wrapped)
_alloc.deallocate(_data, _size);
}
public:
/// Copy operator
Array& operator=(const Array& v) {
this->wrapped = false;
if (v.size() != this->size())
this->resize(v.size());
thrust::copy(v._data, v._data + _size, _data);
return *this;
}
/// Move operator
- Array& operator=(Array&& v) {
+ Array& operator=(Array&& v) noexcept {
if (this != &v) {
- if (!wrapped)
+ if (not wrapped)
_alloc.deallocate(_data, _size);
_data = std::exchange(v._data, nullptr);
_size = std::exchange(v._size, 0);
wrapped = std::exchange(v.wrapped, false);
}
return *this;
}
void wrap(Array & other) {
this->wrapMemory(other._data, other._size);
}
void wrap(const Array & other) {
this->wrapMemory(other._data, other._size);
}
/// Wrap a memory pointer
void wrapMemory(T* data, UInt size) {
this->_data = data;
this->_size = size;
this->wrapped = true;
}
/// Assign a sigle value to the array
void assign(UInt size, const T& value) {
this->resize(size);
thrust::fill(_data, _data + size, value);
}
void setWrapped(bool w) { wrapped = w; }
/// Data pointer access (const)
const T* data() const { return _data; }
/// Data pointer access (non-const)
T* data() { return _data; }
/// Data pointer setter
void setData(T* new_ptr) { _data = new_ptr; }
/// Resize array
void resize(UInt size) {
- if (wrapped == true)
+ if (wrapped)
TAMAAS_EXCEPTION("cannot resize wrapped array");
// Erase array
if (size == 0) {
_alloc.deallocate(_data, _size);
this->_data = nullptr;
this->_size = 0;
return;
}
// Do nothing
if (size == this->_size)
return;
// Allocate new data
_alloc.deallocate(_data, _size);
this->_data = _alloc.allocate(size);
this->_size = size;
}
/// Access operator
inline T& operator[](UInt i) {
TAMAAS_ASSERT(i < _size, "memory overflow");
return _data[i];
}
/// Access operator (const)
inline const T& operator[](UInt i) const {
TAMAAS_ASSERT(i < _size, "memory overflow");
return _data[i];
}
/// Get size of array
inline UInt size() const { return _size; }
private:
T* _data = nullptr;
UInt _size = 0;
bool wrapped = false;
DEFAULT_ALLOCATOR(T) _alloc = DEFAULT_ALLOCATOR(T)();
};
} // namespace tamaas
/* -------------------------------------------------------------------------- */
#endif /* __ARRAY__HH__ */
diff --git a/src/core/fftw_allocator.hh b/src/core/fftw_allocator.hh
index 92e6554..135fa90 100644
--- a/src/core/fftw_allocator.hh
+++ b/src/core/fftw_allocator.hh
@@ -1,50 +1,52 @@
/**
* @file
* @section LICENSE
*
* Copyright (©) 2016-19 EPFL (École Polytechnique Fédérale de Lausanne),
* Laboratory (LSMS - Laboratoire de Simulation en Mécanique des Solides)
*
* This program is free software: you can redistribute it and/or modify
* it under the terms of the GNU Affero General Public License as published
* by the Free Software Foundation, either version 3 of the License, or
* (at your option) any later version.
*
* This program 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 Affero General Public License for more details.
*
* You should have received a copy of the GNU Affero General Public License
* along with this program. If not, see .
*
*/
/* -------------------------------------------------------------------------- */
#ifndef __FFTW_ALLOCATOR_HH__
#define __FFTW_ALLOCATOR_HH__
/* -------------------------------------------------------------------------- */
#include
#include
/* -------------------------------------------------------------------------- */
namespace tamaas {
/// Class allocating
/// [SIMD](http://www.fftw.org/fftw3_doc/SIMD-alignment-and-fftw_005fmalloc.html#SIMD-alignment-and-fftw_005fmalloc)
/// aligned memory
template
struct FFTWAllocator {
/// Allocate memory
- T* allocate(std::size_t n) {
+ T* allocate(std::size_t n) noexcept {
T* p = nullptr;
p = (T*)fftw_malloc(sizeof(T) * n);
return p;
}
/// Free memory
- void deallocate(T* p, __attribute__((unused)) std::size_t n) { fftw_free(p); }
+ void deallocate(T* p, __attribute__((unused)) std::size_t n) noexcept {
+ fftw_free(p);
+ }
};
} // namespace tamaas
#endif // __FFTW_ALLOCATOR_HH__
diff --git a/src/core/grid.cpp b/src/core/grid.cpp
index 74940ae..9224743 100644
--- a/src/core/grid.cpp
+++ b/src/core/grid.cpp
@@ -1,121 +1,122 @@
/**
* @file
* @section LICENSE
*
* Copyright (©) 2016-19 EPFL (École Polytechnique Fédérale de Lausanne),
* Laboratory (LSMS - Laboratoire de Simulation en Mécanique des Solides)
*
* This program is free software: you can redistribute it and/or modify
* it under the terms of the GNU Affero General Public License as published
* by the Free Software Foundation, either version 3 of the License, or
* (at your option) any later version.
*
* This program 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 Affero General Public License for more details.
*
* You should have received a copy of the GNU Affero General Public License
* along with this program. If not, see .
*
*/
/* -------------------------------------------------------------------------- */
#include "grid.hh"
#include "tamaas.hh"
#include
#include
#include
/* -------------------------------------------------------------------------- */
namespace tamaas {
/* -------------------------------------------------------------------------- */
template
Grid::Grid() : GridBase() {
this->n.fill(0);
this->strides.fill(1);
this->nb_components = 1;
}
/* -------------------------------------------------------------------------- */
template
Grid::Grid(const std::array& n, UInt nb_components,
value_type* data)
: GridBase() {
this->n = n;
this->nb_components = nb_components;
UInt size = this->computeSize();
this->data.wrapMemory(data, size);
this->computeStrides();
}
template
Grid::Grid(const Grid& o)
: GridBase(o), n(o.n), strides(o.strides) {}
template
-Grid::Grid(Grid&& o)
- : GridBase(o), n(std::move(o.n)), strides(std::move(o.strides)) {}
+Grid::Grid(Grid&& o) noexcept
+ : GridBase(std::move(o)), n(std::move(o.n)),
+ strides(std::move(o.strides)) {}
/* -------------------------------------------------------------------------- */
template
void Grid::resize(const std::array& n) {
this->resize(n.begin(), n.end());
}
/* -------------------------------------------------------------------------- */
template
void Grid::resize(const std::vector& n) {
TAMAAS_ASSERT(n.size() == dim, "Shape vector not matching grid dimensions");
this->resize(n.begin(), n.end());
}
/* -------------------------------------------------------------------------- */
template
void Grid::resize(std::initializer_list n) {
TAMAAS_ASSERT(n.size() == dim,
"Shape initializer list not matching grid dimensions");
this->resize(std::begin(n), std::end(n));
}
/* -------------------------------------------------------------------------- */
template
void Grid::computeStrides() {
std::copy(n.begin() + 1, n.end(), strides.rbegin() + 2);
strides[dim] = 1;
strides[dim - 1] = this->nb_components;
std::partial_sum(strides.rbegin(), strides.rend(), strides.rbegin(),
std::multiplies());
}
/* -------------------------------------------------------------------------- */
template
void Grid::printself(std::ostream& str) const {
str << "Grid(" << dim << ", " << this->nb_components << ") {";
for (auto& val : *this) {
str << val << ", ";
}
str << "\b\b}";
}
/* -------------------------------------------------------------------------- */
/// Class instanciation
#define GRID_INSTANCIATE_TYPE(type) \
template class Grid; \
template class Grid; \
template class Grid
GRID_INSTANCIATE_TYPE(Real);
GRID_INSTANCIATE_TYPE(UInt);
GRID_INSTANCIATE_TYPE(Complex);
GRID_INSTANCIATE_TYPE(Int);
GRID_INSTANCIATE_TYPE(bool);
#undef GRID_INSTANCIATE_TYPE
} // namespace tamaas
diff --git a/src/core/grid.hh b/src/core/grid.hh
index 8eb56bf..bfa30c3 100644
--- a/src/core/grid.hh
+++ b/src/core/grid.hh
@@ -1,222 +1,222 @@
/**
* @file
* @section LICENSE
*
* Copyright (©) 2016-19 EPFL (École Polytechnique Fédérale de Lausanne),
* Laboratory (LSMS - Laboratoire de Simulation en Mécanique des Solides)
*
* This program is free software: you can redistribute it and/or modify
* it under the terms of the GNU Affero General Public License as published
* by the Free Software Foundation, either version 3 of the License, or
* (at your option) any later version.
*
* This program 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 Affero General Public License for more details.
*
* You should have received a copy of the GNU Affero General Public License
* along with this program. If not, see .
*
*/
/* -------------------------------------------------------------------------- */
#ifndef __GRID_HH__
#define __GRID_HH__
/* -------------------------------------------------------------------------- */
#include "grid_base.hh"
#include "tamaas.hh"
#include
#include
#include
#include
/* -------------------------------------------------------------------------- */
namespace tamaas {
/**
* @brief Multi-dimensional & multi-component array class
*
* This class is a container for multi-component data stored on a multi-
* dimensional grid.
*
* The access function is the parenthesis operator. For a grid of dimension d,
* the operator takes d+1 arguments: the first d arguments are the position on
* the grid and the last one is the component of the stored data.
*
* It is also possible to use the access operator with only one argument, it is
* then considering the grid as a flat array, accessing the given cell of the
* array.
*/
template
class Grid : public GridBase {
template
using is_valid_index = fold_trait;
public:
/* ------------------------------------------------------------------------ */
/* Types */
/* ------------------------------------------------------------------------ */
using value_type = T;
using reference = value_type&;
static constexpr UInt dimension = dim;
/* ------------------------------------------------------------------------ */
/* Constructors */
/* ------------------------------------------------------------------------ */
public:
/// Constructor by default (empty array)
Grid();
/// Constructor from iterators
template
Grid(RandomAccessIterator begin, RandomAccessIterator end,
UInt nb_components);
/// Constructor
Grid(const std::array& n, UInt nb_components)
: Grid(std::begin(n), std::end(n), nb_components) {}
/// Wrap on given data
Grid(const std::array& n, UInt nb_components, value_type* data);
/// Constructor with vectors
Grid(const std::vector& n, UInt nb_components)
: Grid(std::begin(n), std::end(n), nb_components) {}
/// Constructor with initializer list
Grid(const std::initializer_list& n, UInt nb_components)
: Grid(std::begin(n), std::end(n), nb_components) {}
/// Copy constructor
Grid(const Grid& o);
/// Move constructor (transfers data ownership)
- Grid(Grid&& o);
+ Grid(Grid&& o) noexcept;
/// Destructor
~Grid() override = default;
private:
/// Init from standard container
template
void init(const Container& n, UInt nb_components);
public:
/* ------------------------------------------------------------------------ */
/* Common operations */
/* ------------------------------------------------------------------------ */
/// Resize array
void resize(const std::array& n);
/// Resize array (from std::vector)
void resize(const std::vector& n);
/// Resize array (from initializer list)
void resize(std::initializer_list n);
/// Resize array (from iterators)
template
void resize(ForwardIt begin, ForwardIt end);
/// Compute size
inline UInt computeSize() const;
/// Get grid dimension
inline UInt getDimension() const override { return dim; }
/// Compute strides
virtual void computeStrides();
/// Print
virtual void printself(std::ostream& str) const;
/// Get sizes
const std::array& sizes() const { return n; }
/// Get strides
const std::array& getStrides() const { return this->strides; }
/* ------------------------------------------------------------------------ */
/* Access operators (these cannot be moved outside the class) */
/* ------------------------------------------------------------------------ */
/// Variadic access operator (non-const)
template
inline std::enable_if_t::value, T&>
operator()(T1... args) {
/// Checking access integrity
constexpr UInt nargs = sizeof...(T1);
static_assert(nargs == dim + 1 || nargs == 1 || nargs == dim,
"number of arguments in operator() does not match dimension");
constexpr UInt start = (nargs == 1) ? dim : 0;
UInt offset = unpackOffset(0, start, args...);
return this->data[offset];
}
/// Variadic access operator
template
inline std::enable_if_t::value, const T&>
operator()(T1... args) const {
/// Checking access integrity
constexpr UInt nargs = sizeof...(T1);
static_assert(nargs == dim + 1 || nargs == 1 || nargs == dim,
"number of arguments in operator() does not match dimension");
constexpr UInt start = (nargs == 1) ? dim : 0;
UInt offset = unpackOffset(0, start, args...);
return this->data[offset];
}
/// Tuple index access operator
template
inline T& operator()(std::array tuple);
template
inline const T& operator()(std::array tuple) const;
private:
/// Unpacking the arguments of variadic ()
template
inline UInt unpackOffset(UInt offset, UInt index_pos, UInt index,
T1... rest) const;
/// End case for recursion
template
inline UInt unpackOffset(UInt offset, UInt index_pos, UInt index) const;
/// Computing offset for a tuple index
template
inline UInt computeOffset(std::array tuple) const;
/* ------------------------------------------------------------------------ */
/* Move/Copy operators */
/* ------------------------------------------------------------------------ */
public:
using GridBase::operator=;
// = operator
Grid& operator=(const Grid& other);
// = operator (not const input, otherwise gcc is confused)
Grid& operator=(Grid& other);
// = operator (move)
- Grid& operator=(Grid&& other);
+ Grid& operator=(Grid&& other) noexcept;
// Copy data from another grid
template
void copy(const Grid& other);
// Move data from another grid
template
- void move(Grid&& other);
+ void move(Grid&& other) noexcept;
// Wrap memory (non-const)
void wrap(Grid& other) {
GridBase::wrap(other);
std::copy(other.n.begin(), other.n.end(), n.begin());
this->computeStrides();
}
// Wrap memory
void wrap(const Grid& other) {
GridBase::wrap(other);
std::copy(other.n.begin(), other.n.end(), n.begin());
this->computeStrides();
}
/* ------------------------------------------------------------------------ */
/* Member variables */
/* ------------------------------------------------------------------------ */
protected:
std::array n; ///< shape of grid: size per dimension
std::array strides; ///< strides for access
};
} // namespace tamaas
/* -------------------------------------------------------------------------- */
/* Inline/template function definitions */
/* -------------------------------------------------------------------------- */
#include "grid_tmpl.hh"
/* -------------------------------------------------------------------------- */
#endif // __GRID_HH__
diff --git a/src/core/grid_base.hh b/src/core/grid_base.hh
index 9cf1d59..b0398d6 100644
--- a/src/core/grid_base.hh
+++ b/src/core/grid_base.hh
@@ -1,331 +1,331 @@
/**
* @file
* @section LICENSE
*
* Copyright (©) 2016-19 EPFL (École Polytechnique Fédérale de Lausanne),
* Laboratory (LSMS - Laboratoire de Simulation en Mécanique des Solides)
*
* This program is free software: you can redistribute it and/or modify
* it under the terms of the GNU Affero General Public License as published
* by the Free Software Foundation, either version 3 of the License, or
* (at your option) any later version.
*
* This program 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 Affero General Public License for more details.
*
* You should have received a copy of the GNU Affero General Public License
* along with this program. If not, see .
*
*/
/* -------------------------------------------------------------------------- */
#ifndef __GRID_BASE_HH__
#define __GRID_BASE_HH__
/* -------------------------------------------------------------------------- */
#include "array.hh"
#include "iterator.hh"
#include "loop.hh"
#include "tamaas.hh"
#include "static_types.hh"
#include
#include
#include
#include
/* -------------------------------------------------------------------------- */
namespace tamaas {
/**
* @brief Dimension agnostic grid with components stored per points
*/
template
class GridBase {
public:
/// Constructor by default
GridBase() = default;
/// Copy constructor
GridBase(const GridBase& o) { this->copy(o); }
/// Move constructor (transfers data ownership)
- GridBase(GridBase&& o)
+ GridBase(GridBase&& o) noexcept
: data(std::move(o.data)),
nb_components(std::exchange(o.nb_components, 1)) {}
/// Destructor
virtual ~GridBase() = default;
/* ------------------------------------------------------------------------ */
/* Iterator class */
/* ------------------------------------------------------------------------ */
public:
using iterator = iterator_::iterator;
using const_iterator = iterator_::iterator;
using value_type = T;
using reference = value_type&;
public:
/// Get grid dimension
virtual UInt getDimension() const { return 1; }
/// Get internal data pointer (const)
const T* getInternalData() const { return this->data.data(); }
/// Get internal data pointer (non-const)
T* getInternalData() { return this->data.data(); }
/// Get number of components
UInt getNbComponents() const { return nb_components; }
/// Set number of components
void setNbComponents(UInt n) { nb_components = n; }
/// Get total size
virtual UInt dataSize() const { return this->data.size(); }
/// Get number of points
UInt getNbPoints() const {
return this->dataSize() / this->getNbComponents();
}
/// Set components
void uniformSetComponents(const GridBase& vec);
/// Resize
void resize(UInt size) {
this->data.resize(size);
this->data.assign(size, T(0.));
}
/* ------------------------------------------------------------------------ */
/* Iterator methods */
/* ------------------------------------------------------------------------ */
virtual iterator begin(UInt n = 1) {
return iterator(this->getInternalData(), n);
}
virtual iterator end(UInt n = 1) {
return iterator(this->getInternalData() + this->dataSize(), n);
}
virtual const_iterator begin(UInt n = 1) const {
return const_iterator(this->getInternalData(), n);
}
virtual const_iterator end(UInt n = 1) const {
return const_iterator(this->getInternalData() + this->dataSize(), n);
}
/* ------------------------------------------------------------------------ */
/* Operators */
/* ------------------------------------------------------------------------ */
inline T& operator()(UInt i) { return this->data[i]; }
inline const T& operator()(UInt i) const { return this->data[i]; }
#define VEC_OPERATOR(op) \
template \
GridBase& operator op(const GridBase& other)
VEC_OPERATOR(+=);
VEC_OPERATOR(*=);
VEC_OPERATOR(-=);
VEC_OPERATOR(/=);
#undef VEC_OPERATOR
template
T dot(const GridBase& other) const;
#define SCALAR_OPERATOR(op) GridBase& operator op(const T& e)
SCALAR_OPERATOR(+=);
SCALAR_OPERATOR(*=);
SCALAR_OPERATOR(-=);
SCALAR_OPERATOR(/=);
SCALAR_OPERATOR(=);
#undef SCALAR_OPERATOR
/// Broadcast operators
#define BROADCAST_OPERATOR(op) \
template \
GridBase& operator op(const StaticArray& b)
BROADCAST_OPERATOR(+=);
BROADCAST_OPERATOR(-=);
BROADCAST_OPERATOR(*=);
BROADCAST_OPERATOR(/=);
BROADCAST_OPERATOR(=);
#undef BROADCAST_OPERATOR
/* ------------------------------------------------------------------------ */
/* Min/max */
/* ------------------------------------------------------------------------ */
T min() const;
T max() const;
T sum() const;
T mean() const { return sum() / static_cast(dataSize()); }
T var() const;
/* ------------------------------------------------------------------------ */
/* Move/Copy */
/* ------------------------------------------------------------------------ */
public:
GridBase& operator=(const GridBase& o) {
this->copy(o);
return *this;
}
GridBase& operator=(GridBase& o) {
this->copy(o);
return *this;
}
- GridBase& operator=(GridBase&& o) {
+ GridBase& operator=(GridBase&& o) noexcept {
this->move(std::move(o));
return *this;
}
template
void copy(const GridBase& other) {
if (other.dataSize() != this->dataSize())
this->resize(other.dataSize());
thrust::copy(other.begin(), other.end(), this->begin());
nb_components = other.nb_components;
}
template
- void move(GridBase&& other) {
+ void move(GridBase&& other) noexcept {
data = std::move(other.data);
nb_components = std::exchange(other.nb_components, 1);
}
GridBase& wrap(GridBase& o) {
data.wrap(o.data);
nb_components = o.nb_components;
return *this;
}
GridBase& wrap(const GridBase& o) {
data.wrap(o.data);
nb_components = o.nb_components;
return *this;
}
protected:
Array data;
UInt nb_components = 1;
};
/* -------------------------------------------------------------------------- */
template
void GridBase::uniformSetComponents(const GridBase& vec) {
if (!(vec.dataSize() == this->nb_components))
TAMAAS_EXCEPTION("Cannot set grid field with values of vector");
auto begin_it(begin());
auto end_it(end()); // silencing nvcc warnings
for (auto it = begin(); it < end_it; ++it) {
UInt i = it - begin_it;
*it = vec(i % this->nb_components);
}
}
/* --------------------------------------------------------------------------
*/
#define SCALAR_OPERATOR_IMPL(op) \
template \
inline GridBase& GridBase::operator op(const T& e) { \
Loop::loop([e] CUDA_LAMBDA(T& val) { val op e; }, *this); \
return *this; \
}
SCALAR_OPERATOR_IMPL(+=);
SCALAR_OPERATOR_IMPL(*=);
SCALAR_OPERATOR_IMPL(-=);
SCALAR_OPERATOR_IMPL(/=);
SCALAR_OPERATOR_IMPL(=);
#undef SCALAR_OPERATOR_IMPL
/* -------------------------------------------------------------------------- */
#define VEC_OPERATOR_IMPL(op) \
template \
template \
inline GridBase& GridBase::operator op(const GridBase& other) { \
TAMAAS_ASSERT(other.dataSize() == this->dataSize(), \
"surface size does not match"); \
Loop::loop([] CUDA_LAMBDA(T& x, const T1& y) { x op y; }, *this, other); \
return *this; \
}
VEC_OPERATOR_IMPL(+=);
VEC_OPERATOR_IMPL(-=);
VEC_OPERATOR_IMPL(*=);
VEC_OPERATOR_IMPL(/=);
#undef VEC_OPERATOR
/* -------------------------------------------------------------------------- */
#define BROADCAST_OPERATOR_IMPL(op) \
template \
template \
GridBase& GridBase::operator op(const StaticArray& b) { \
TAMAAS_ASSERT(this->dataSize() % size == 0, \
"Broadcast operator cannot broadcast" \
<< size << " on array of size " << this->dataSize()); \
Tensor::type, size> a(b); \
Loop::loop([&a] CUDA_LAMBDA(VectorProxy g) { g op a; }, \
range>(*this)); \
return *this; \
}
BROADCAST_OPERATOR_IMPL(+=);
BROADCAST_OPERATOR_IMPL(-=);
BROADCAST_OPERATOR_IMPL(*=);
BROADCAST_OPERATOR_IMPL(/=);
BROADCAST_OPERATOR_IMPL(=);
#undef BROADCAST_OPERATOR_IMPL
/* -------------------------------------------------------------------------- */
template
inline T GridBase::min() const {
// const auto id = [] CUDA_LAMBDA (const T&x){return x;};
return Loop::reduce([] CUDA_LAMBDA(const T& x) { return x; },
*this);
}
/* -------------------------------------------------------------------------- */
template
inline T GridBase::max() const {
return Loop::reduce([] CUDA_LAMBDA(const T& x) { return x; },
*this);
}
/* -------------------------------------------------------------------------- */
template
inline T GridBase::sum() const {
return Loop::reduce([] CUDA_LAMBDA(const T& x) { return x; },
*this);
}
/* -------------------------------------------------------------------------- */
template
inline T GridBase::var() const {
const T mu = mean();
const T var = Loop::reduce(
[mu] CUDA_LAMBDA(const T& x) { return (x - mu) * (x - mu); }, *this);
return var / (dataSize() - 1);
}
/* -------------------------------------------------------------------------- */
template
template
T GridBase::dot(const GridBase& other) const {
return Loop::reduce(
[] CUDA_LAMBDA(const T& x, const T1& y) { return x * y; }, *this, other);
}
} // namespace tamaas
#endif // __GRID_BASE_HH__
diff --git a/src/core/grid_hermitian.hh b/src/core/grid_hermitian.hh
index 41e947b..b914a05 100644
--- a/src/core/grid_hermitian.hh
+++ b/src/core/grid_hermitian.hh
@@ -1,163 +1,163 @@
/**
* @file
* @section LICENSE
*
* Copyright (©) 2016-19 EPFL (École Polytechnique Fédérale de Lausanne),
* Laboratory (LSMS - Laboratoire de Simulation en Mécanique des Solides)
*
* This program is free software: you can redistribute it and/or modify
* it under the terms of the GNU Affero General Public License as published
* by the Free Software Foundation, either version 3 of the License, or
* (at your option) any later version.
*
* This program 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 Affero General Public License for more details.
*
* You should have received a copy of the GNU Affero General Public License
* along with this program. If not, see .
*
*/
/* -------------------------------------------------------------------------- */
#ifndef __GRID_HERMITIAN_HH__
#define __GRID_HERMITIAN_HH__
/* -------------------------------------------------------------------------- */
#include "grid.hh"
#include "tamaas.hh"
#include
#include
/* -------------------------------------------------------------------------- */
namespace tamaas {
static complex dummy;
/**
* @brief Multi-dimensional, multi-component herimitian array
*
* This class represents an array of hermitian data, meaning it has dimensions
* of: n1 * n2 * n3 * ... * (nx / 2 + 1)
*
* However, it acts as a fully dimensioned array, returning a dummy reference
* for data outside its real dimension, allowing one to write full (and
* inefficient) loops without really worrying about the reduced dimension.
*
* It would however be better to just use the true dimensions of the surface
* for all intents and purposes, as it saves computation time.
*/
template
class GridHermitian : public Grid, dim> {
public:
using value_type = complex;
/* ------------------------------------------------------------------------ */
/* Constructors */
/* ------------------------------------------------------------------------ */
public:
GridHermitian() = default;
GridHermitian(const GridHermitian& o) = default;
- GridHermitian(GridHermitian&& o) = default;
+ GridHermitian(GridHermitian&& o) noexcept = default;
using Grid::Grid;
using Grid, dim>::operator=;
/* ------------------------------------------------------------------------ */
/* Access operators */
/* ------------------------------------------------------------------------ */
template
inline complex& operator()(T1... args);
template
inline const complex& operator()(T1... args) const;
static std::array
hermitianDimensions(std::array n) {
n[dim - 1] /= 2;
n[dim - 1] += 1;
return n;
}
static std::vector
hermitianDimensions(std::vector n) {
TAMAAS_ASSERT(n.size() == dim,
"Requested size dimension does not match class dimension");
n[dim - 1] /= 2;
n[dim - 1] += 1;
return n;
}
private:
template
inline void packTuple(UInt* t, UInt i, T1... args) const;
template
inline void packTuple(UInt* t, UInt i) const;
};
/* -------------------------------------------------------------------------- */
/* Inline function definitions */
/* -------------------------------------------------------------------------- */
template
template
inline void GridHermitian::packTuple(UInt* t, UInt i,
T1... args) const {
*t = i;
packTuple(t + 1, args...);
}
template
template
inline void GridHermitian::packTuple(UInt* t, UInt i) const {
*t = i;
}
template
template
inline complex& GridHermitian::operator()(T1... args) {
constexpr UInt nargs = sizeof...(T1);
static_assert(nargs <= dim + 1, "Too many arguments");
UInt last_index = std::get(std::forward_as_tuple(args...));
if (nargs != 1 && last_index >= this->n[dim - 1]) {
std::array tuple = {{0}};
packTuple(tuple.data(), args...);
for (UInt i = 0; i < dim; i++) {
if (tuple[i] && i != dim - 1)
tuple[i] = this->n[i] - tuple[i];
else if (tuple[i] && i == dim - 1)
tuple[i] = 2 * (this->n[i] - 1) - tuple[i];
}
dummy = conj(this->Grid, dim>::operator()(tuple));
return dummy;
}
else
return Grid, dim>::operator()(args...);
}
template
template
inline const complex& GridHermitian::operator()(T1... args) const {
constexpr UInt nargs = sizeof...(T1);
static_assert(nargs <= dim + 1, "Too many arguments");
UInt last_index = std::get(std::forward_as_tuple(args...));
if (nargs != 1 && last_index >= this->n[dim - 1]) {
std::array tuple = {{0}};
packTuple(tuple.data(), args...);
for (UInt i = 0; i < dim; i++) {
if (tuple[i] && i != dim - 1)
tuple[i] = this->n[i] - tuple[i];
else if (tuple[i] && i == dim - 1)
tuple[i] = 2 * (this->n[i] - 1) - tuple[i];
}
dummy = conj(this->Grid, dim>::operator()(tuple));
return dummy;
}
else
return Grid, dim>::operator()(args...);
}
} // namespace tamaas
#endif // __GRID_HERMITIAN_HH__
diff --git a/src/core/grid_tmpl.hh b/src/core/grid_tmpl.hh
index c928c76..b18c0cf 100644
--- a/src/core/grid_tmpl.hh
+++ b/src/core/grid_tmpl.hh
@@ -1,159 +1,159 @@
/**
* @file
* @section LICENSE
*
* Copyright (©) 2016-19 EPFL (École Polytechnique Fédérale de Lausanne),
* Laboratory (LSMS - Laboratoire de Simulation en Mécanique des Solides)
*
* This program is free software: you can redistribute it and/or modify
* it under the terms of the GNU Affero General Public License as published
* by the Free Software Foundation, either version 3 of the License, or
* (at your option) any later version.
*
* This program 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 Affero General Public License for more details.
*
* You should have received a copy of the GNU Affero General Public License
* along with this program. If not, see .
*
*/
/* -------------------------------------------------------------------------- */
#ifndef __GRID_TMPL_HH__
#define __GRID_TMPL_HH__
/* -------------------------------------------------------------------------- */
#include "grid.hh"
#include "tamaas.hh"
namespace tamaas {
template
template
Grid::Grid(RandomAccessIterator begin, RandomAccessIterator end,
UInt nb_components)
: GridBase() {
if (end - begin != dim) {
TAMAAS_EXCEPTION("Provided sizes (" << end - begin
<< ") for grid do not match dimension ("
<< dim << ")");
}
this->nb_components = nb_components;
this->resize(begin, end);
}
/* -------------------------------------------------------------------------- */
template
inline UInt Grid::computeSize() const {
UInt size = 1;
for (UInt i = 0; i < dim; i++)
size *= n[i];
size *= this->nb_components;
return size;
}
/* -------------------------------------------------------------------------- */
template
template
void Grid::resize(ForwardIt begin, ForwardIt end) {
std::copy(begin, end, this->n.begin());
UInt size = this->computeSize();
GridBase::resize(size);
this->computeStrides();
}
/* -------------------------------------------------------------------------- */
template
template
inline UInt Grid::unpackOffset(UInt offset, UInt index_pos, UInt index,
T1... rest) const {
offset += index * strides[index_pos];
return unpackOffset(offset, index_pos + 1, rest...); // tail-rec bb
}
template
template
inline UInt Grid::unpackOffset(UInt offset, UInt index_pos,
UInt index) const {
return offset + index * strides[index_pos];
}
/* -------------------------------------------------------------------------- */
template
template
inline UInt Grid::computeOffset(std::array tuple) const {
static_assert(tdim == dim or tdim == dim + 1, "Tuple dimension is invalid");
return std::inner_product(tuple.begin(), tuple.end(), strides.begin(), 0);
}
template
template
inline T& Grid::operator()(std::array tuple) {
UInt offset = computeOffset(tuple);
return this->data[offset];
}
template
template
inline const T& Grid::operator()(std::array tuple) const {
UInt offset = computeOffset(tuple);
return this->data[offset];
}
/* -------------------------------------------------------------------------- */
template
Grid& Grid::operator=(const Grid& other) {
this->copy(other);
return *this;
}
template
Grid& Grid::operator=(Grid& other) {
this->copy(other);
return *this;
}
template
-Grid& Grid::operator=(Grid&& other) {
+Grid& Grid::operator=(Grid&& other) noexcept {
this->move(std::move(other));
return *this;
}
/* -------------------------------------------------------------------------- */
template
template
void Grid::copy(const Grid& other) {
GridBase::copy(other);
this->n = other.n;
this->strides = other.strides;
}
template
template
-void Grid::move(Grid&& other) {
+void Grid::move(Grid&& other) noexcept {
GridBase::move(std::move(other));
this->n = std::move(other.n);
this->strides = std::move(other.strides);
}
/* -------------------------------------------------------------------------- */
/* Stream operator */
/* -------------------------------------------------------------------------- */
template
inline std::ostream& operator<<(std::ostream& stream,
const Grid& _this) {
_this.printself(stream);
return stream;
}
} // namespace tamaas
/* -------------------------------------------------------------------------- */
#endif // __GRID_TMPL_HH__
diff --git a/src/core/grid_view.hh b/src/core/grid_view.hh
index 5589ecd..cd34ac4 100644
--- a/src/core/grid_view.hh
+++ b/src/core/grid_view.hh
@@ -1,159 +1,159 @@
/**
* @file
* @section LICENSE
*
* Copyright (©) 2016-19 EPFL (École Polytechnique Fédérale de Lausanne),
* Laboratory (LSMS - Laboratoire de Simulation en Mécanique des Solides)
*
* This program is free software: you can redistribute it and/or modify
* it under the terms of the GNU Affero General Public License as published
* by the Free Software Foundation, either version 3 of the License, or
* (at your option) any later version.
*
* This program 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 Affero General Public License for more details.
*
* You should have received a copy of the GNU Affero General Public License
* along with this program. If not, see .
*
*/
/* -------------------------------------------------------------------------- */
#ifndef __GRID_VIEW_HH__
#define __GRID_VIEW_HH__
/* -------------------------------------------------------------------------- */
#include "grid.hh"
#include "tamaas.hh"
#include
/* -------------------------------------------------------------------------- */
namespace tamaas {
/* -------------------------------------------------------------------------- */
/**
* @brief View type on grid
* This is a view on a *contiguous* chunk of data defined by a grid
*/
template class Base, typename T, UInt base_dim,
UInt dim>
class GridView : public Base {
public:
using iterator = typename Base::iterator;
using const_iterator = typename Base::const_iterator;
using value_type = typename Base::value_type;
using referecence = typename Base::reference;
public:
/// Constructor
GridView(GridBase::value_type>& grid_base,
const std::vector& multi_index, Int component = -1);
/// Move constructor
- GridView(GridView&& o)
+ GridView(GridView&& o) noexcept
: Base(std::forward>(o)),
grid(std::exchange(o.grid, nullptr)) {}
/// Destructor
~GridView() override = default;
UInt dataSize() const override { return this->computeSize(); }
/// Iterators
iterator begin(UInt /*n*/ = 1) override {
return iterator(this->getInternalData(), this->strides.back());
}
iterator end(UInt /*n*/ = 1) override {
return iterator(this->getInternalData() +
this->dataSize() * this->strides.back(),
this->strides.back());
}
const_iterator begin(UInt /*n*/ = 1) const override {
return const_iterator(this->getInternalData(), this->strides.back());
}
const_iterator end(UInt /*n*/ = 1) const override {
return const_iterator(this->getInternalData() +
this->dataSize() * this->strides.back(),
this->strides.back());
}
protected:
Base* grid;
};
/* -------------------------------------------------------------------------- */
template class Base, typename T, UInt base_dim,
typename... Args>
GridView