SLIP  1.4
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
HyperVolume.hpp
Go to the documentation of this file.
1 /*
2  * Copyright(c):
3  * Signal Image and Communications (SIC) Department
4  * http://www.sic.sp2mi.univ-poitiers.fr/
5  * - University of Poitiers, France http://www.univ-poitiers.fr
6  * - XLIM Institute UMR CNRS 7252 http://www.xlim.fr/
7  *
8  * and
9  *
10  * D2 Fluid, Thermic and Combustion
11  * - University of Poitiers, France http://www.univ-poitiers.fr
12  * - PPRIME Institute - UPR CNRS 3346 http://www.pprime.fr
13  * - ISAE-ENSMA http://www.ensma.fr
14  *
15  * Contributor(s):
16  * The SLIP team,
17  * Benoit Tremblais <tremblais_AT_sic.univ-poitiers.fr>,
18  * Laurent David <laurent.david_AT_lea.univ-poitiers.fr>,
19  * Ludovic Chatellier <ludovic.chatellier_AT_univ-poitiers.fr>,
20  * Lionel Thomas <lionel.thomas_AT_univ-poitiers.fr>,
21  * Denis Arrivault <arrivault_AT_sic.univ-poitiers.fr>,
22  * Julien Dombre <julien.dombre_AT_univ-poitiers.fr>.
23  *
24  * Description:
25  * The Simple Library of Image Processing (SLIP) is a new image processing
26  * library. It is written in the C++ language following as much as possible
27  * the ISO/ANSI C++ standard. It is consequently compatible with any system
28  * satisfying the ANSI C++ complience. It works on different Unix , Linux ,
29  * Mircrosoft Windows and Mac OS X plateforms. SLIP is a research library that
30  * was created by the Signal, Image and Communications (SIC) departement of
31  * the XLIM, UMR 7252 CNRS Institute in collaboration with the Fluids, Thermic
32  * and Combustion departement of the P', UPR 3346 CNRS Institute of the
33  * University of Poitiers.
34  *
35  * The SLIP Library source code has been registered to the APP (French Agency
36  * for the Protection of Programs) by the University of Poitiers and CNRS,
37  * under registration number IDDN.FR.001.300034.000.S.P.2010.000.21000.
38 
39  * http://www.sic.sp2mi.univ-poitiers.fr/slip/
40  *
41  * This software is governed by the CeCILL-C license under French law and
42  * abiding by the rules of distribution of free software. You can use,
43  * modify and/ or redistribute the software under the terms of the CeCILL-C
44  * license as circulated by CEA, CNRS and INRIA at the following URL
45  * http://www.cecill.info.
46  * As a counterpart to the access to the source code and rights to copy,
47  * modify and redistribute granted by the license, users are provided only
48  * with a limited warranty and the software's author, the holder of the
49  * economic rights, and the successive licensors have only limited
50  * liability.
51  *
52  * In this respect, the user's attention is drawn to the risks associated
53  * with loading, using, modifying and/or developing or reproducing the
54  * software by the user in light of its specific status of free software,
55  * that may mean that it is complicated to manipulate, and that also
56  * therefore means that it is reserved for developers and experienced
57  * professionals having in-depth computer knowledge. Users are therefore
58  * encouraged to load and test the software's suitability as regards their
59  * requirements in conditions enabling the security of their systems and/or
60  * data to be ensured and, more generally, to use and operate it in the
61  * same conditions as regards security.
62  *
63  * The fact that you are presently reading this means that you have had
64  * knowledge of the CeCILL-C license and that you accept its terms.
65  */
66 
74 #ifndef SLIP_HYPER_VOLUME_HPP
75 #define SLIP_HYPER_VOLUME_HPP
76 
77 #include <iostream>
78 #include <iomanip>
79 #include <fstream>
80 #include <iterator>
81 #include <cassert>
82 #include <numeric>
83 #include <algorithm>
84 #include <cmath>
85 #include <vector>
86 #include <string>
87 #include <cstddef>
88 
89 #include "Matrix4d.hpp"
90 #include "stride_iterator.hpp"
91 #include "apply.hpp"
92 #include "iterator4d_box.hpp"
93 #include "iterator4d_range.hpp"
94 
95 #include <boost/serialization/access.hpp>
96 #include <boost/serialization/split_member.hpp>
97 #include <boost/serialization/version.hpp>
98 
99 
100 namespace slip
101 {
102 
103 template <typename T>
105 
106 template <typename T>
107 class Matrix4d;
108 
109 template <typename T>
110 std::ostream& operator<<(std::ostream & out, const slip::HyperVolume<T>& a);
111 
112 template<typename T>
113 bool operator==(const slip::HyperVolume<T>& x,
114  const slip::HyperVolume<T>& y);
115 
116 template<typename T>
117 bool operator!=(const slip::HyperVolume<T>& x,
118  const slip::HyperVolume<T>& y);
119 
120 template<typename T>
121 bool operator<(const slip::HyperVolume<T>& x,
122  const slip::HyperVolume<T>& y);
123 
124 template<typename T>
125 bool operator>(const slip::HyperVolume<T>& x,
126  const slip::HyperVolume<T>& y);
127 
128 template<typename T>
129 bool operator<=(const slip::HyperVolume<T>& x,
130  const slip::HyperVolume<T>& y);
131 
132 template<typename T>
133 bool operator>=(const slip::HyperVolume<T>& x,
134  const slip::HyperVolume<T>& y);
135 
161 template <typename T>
162 class HyperVolume
163 {
164 public :
165 
166  typedef T value_type;
167  typedef HyperVolume<T> self;
168 
169  typedef value_type* pointer;
170  typedef const value_type* const_pointer;
172  typedef const value_type& const_reference;
173 
174  typedef ptrdiff_t difference_type;
175  typedef std::size_t size_type;
176 
177  typedef pointer iterator;
179 
180  typedef std::reverse_iterator<iterator> reverse_iterator;
181  typedef std::reverse_iterator<const_iterator> const_reverse_iterator;
182 
183 
184  //slab, slice, row and col iterator
193 
202 
203  typedef std::reverse_iterator<slab_iterator> reverse_slab_iterator;
204  typedef std::reverse_iterator<const_slab_iterator> const_reverse_slab_iterator;
205  typedef std::reverse_iterator<slice_iterator> reverse_slice_iterator;
206  typedef std::reverse_iterator<const_slice_iterator> const_reverse_slice_iterator;
207  typedef std::reverse_iterator<iterator> reverse_row_iterator;
208  typedef std::reverse_iterator<const_iterator> const_reverse_row_iterator;
209  typedef std::reverse_iterator<col_iterator> reverse_col_iterator;
210  typedef std::reverse_iterator<const_col_iterator> const_reverse_col_iterator;
211 
212  typedef std::reverse_iterator<slab_range_iterator> reverse_slab_range_iterator;
213  typedef std::reverse_iterator<const_slab_range_iterator> const_reverse_slab_range_iterator;
214  typedef std::reverse_iterator<slice_range_iterator> reverse_slice_range_iterator;
215  typedef std::reverse_iterator<const_slice_range_iterator> const_reverse_slice_range_iterator;
216  typedef std::reverse_iterator<row_range_iterator> reverse_row_range_iterator;
217  typedef std::reverse_iterator<const_row_range_iterator> const_reverse_row_range_iterator;
218  typedef std::reverse_iterator<col_range_iterator> reverse_col_range_iterator;
219  typedef std::reverse_iterator<const_col_range_iterator> const_reverse_col_range_iterator;
220 
221  //iterator 4d
226 
227  typedef std::reverse_iterator<iterator4d> reverse_iterator4d;
228  typedef std::reverse_iterator<const_iterator4d> const_reverse_iterator4d;
229  typedef std::reverse_iterator<iterator4d_range> reverse_iterator4d_range;
230  typedef std::reverse_iterator<const_iterator4d_range> const_reverse_iterator4d_range;
231 
232  //default iterator of the container
235 
236  //Range
238 
239  static const std::size_t DIM = 4;
240 public:
245 
249  HyperVolume();
250 
260  HyperVolume(const std::size_t d1,
261  const std::size_t d2,
262  const std::size_t d3,
263  const std::size_t d4);
264 
273  HyperVolume(const std::size_t d1,
274  const std::size_t d2,
275  const std::size_t d3,
276  const std::size_t d4,
277  const T& val);
286  HyperVolume(const std::size_t d1,
287  const std::size_t d2,
288  const std::size_t d3,
289  const std::size_t d4,
290  const T* val);
291 
304  template<typename InputIterator>
306  const size_type d2,
307  const size_type d3,
308  const size_type d4,
309  InputIterator first,
310  InputIterator last):
311  matrix_(new slip::Matrix4d<T>(d1,d2,d3,d4,first,last))
312  {}
313 
317  HyperVolume(const HyperVolume<T>& rhs);
318 
319 
323  ~HyperVolume();
324 
325 
337  void resize(std::size_t d1,
338  std::size_t d2,
339  std::size_t d3,
340  std::size_t d4,
341  const T& val = T());
342 
343  //********************************************************************
348 
349  //****************************************************************************
350  // One dimensional iterators
351  //****************************************************************************
352 
353 
354 
355  //----------------------Global iterators------------------------------
356 
376  const_iterator begin() const;
377 
397  iterator begin();
398 
418  iterator end();
419 
439  const_iterator end() const;
440 
461 
482 
503 
524 
527  //****************************************************************************
528  // One dimensional stride iterators
529  //****************************************************************************
530 
535  //--------------------One dimensional slab iterators----------------------
536 
560  slab_iterator slab_begin(const size_type slice, const size_type row, const size_type col);
561 
585  const_slab_iterator slab_begin(const size_type slice, const size_type row, const size_type col) const;
586 
610  slab_iterator slab_end(const size_type slice, const size_type row, const size_type col);
611 
635  const_slab_iterator slab_end(const size_type slice, const size_type row, const size_type col) const;
636 
637 
662  reverse_slab_iterator slab_rbegin(const size_type slice, const size_type row, const size_type col);
663 
688  const_reverse_slab_iterator slab_rbegin(const size_type slice, const size_type row, const size_type col) const;
689 
714  reverse_slab_iterator slab_rend(const size_type slice, const size_type row, const size_type col);
715 
716 
741  const_reverse_slab_iterator slab_rend(const size_type slice, const size_type row, const size_type col) const;
742 
749  //--------------------One dimensional slice iterators----------------------
750 
751 
775  slice_iterator slice_begin(const size_type slab, const size_type row, const size_type col);
776 
800  const_slice_iterator slice_begin(const size_type slab, const size_type row, const size_type col) const;
801 
825  slice_iterator slice_end(const size_type slab, const size_type row, const size_type col);
826 
850  const_slice_iterator slice_end(const size_type slab, const size_type row, const size_type col) const;
851 
852 
877  reverse_slice_iterator slice_rbegin(const size_type slab, const size_type row, const size_type col);
878 
903  const_reverse_slice_iterator slice_rbegin(const size_type slab, const size_type row, const size_type col) const;
904 
929  reverse_slice_iterator slice_rend(const size_type slab, const size_type row, const size_type col);
930 
931 
956  const_reverse_slice_iterator slice_rend(const size_type slab, const size_type row, const size_type col) const;
957 
964  //-------------------row iterators----------
965 
979  row_iterator row_begin(const size_type slab, const size_type slice,
980  const size_type row);
981 
995  const_row_iterator row_begin(const size_type slab, const size_type slice,
996  const size_type row) const;
997 
998 
1012  row_iterator row_end(const size_type slab, const size_type slice,
1013  const size_type row);
1014 
1015 
1029  const_row_iterator row_end(const size_type slab, const size_type slice,
1030  const size_type row) const;
1031 
1032 
1046  reverse_row_iterator row_rbegin(const size_type slab, const size_type slice,
1047  const size_type row);
1048 
1062  const_reverse_row_iterator row_rbegin(const size_type slab, const size_type slice,
1063  const size_type row) const;
1064 
1065 
1079  reverse_row_iterator row_rend(const size_type slab, const size_type slice,
1080  const size_type row);
1081 
1082 
1096  const_reverse_row_iterator row_rend(const size_type slab, const size_type slice,
1097  const size_type row) const;
1098 
1104  //-------------------col iterators----------
1105 
1119  col_iterator col_begin(const size_type slab, const size_type slice,
1120  const size_type col);
1121 
1122 
1136  const_col_iterator col_begin(const size_type slab, const size_type slice,
1137  const size_type col) const;
1138 
1152  col_iterator col_end(const size_type slab, const size_type slice,
1153  const size_type col);
1154 
1155 
1169  const_col_iterator col_end(const size_type slab, const size_type slice,
1170  const size_type col) const;
1171 
1172 
1186  reverse_col_iterator col_rbegin(const size_type slab, const size_type slice,
1187  const size_type col);
1188 
1189 
1203  const_reverse_col_iterator col_rbegin(const size_type slab, const size_type slice,
1204  const size_type col) const;
1205 
1219  reverse_col_iterator col_rend(const size_type slab, const size_type slice,
1220  const size_type col);
1221 
1222 
1236  const_reverse_col_iterator col_rend(const size_type slab, const size_type slice,
1237  const size_type col) const;
1238 
1241  //****************************************************************************
1242  // One dimensional range iterators
1243  //****************************************************************************
1244 
1249  //------------------------slab range iterators -----------------------
1250 
1278  slab_range_iterator slab_begin(const size_type slice,const size_type row,const size_type col,
1279  const slip::Range<int>& range);
1280 
1281 
1309  slab_range_iterator slab_end(const size_type slice,const size_type row,const size_type col,
1310  const slip::Range<int>& range);
1311 
1312 
1339  const_slab_range_iterator slab_begin(const size_type slice,const size_type row,const size_type col,
1340  const slip::Range<int>& range) const;
1341 
1368  const_slab_range_iterator slab_end(const size_type slice,const size_type row,const size_type col,
1369  const slip::Range<int>& range) const;
1370 
1371 
1400  reverse_slab_range_iterator slab_rbegin(const size_type slice,const size_type row,const size_type col,
1401  const slip::Range<int>& range);
1402 
1403 
1432  reverse_slab_range_iterator slab_rend(const size_type slice,const size_type row,const size_type col,
1433  const slip::Range<int>& range);
1434 
1435 
1465  const slip::Range<int>& range) const;
1466 
1496  const slip::Range<int>& range) const;
1497 
1498 
1506  //------------------------slice range iterators -----------------------
1507 
1535  slice_range_iterator slice_begin(const size_type slab,const size_type row,const size_type col,
1536  const slip::Range<int>& range);
1537 
1538 
1566  slice_range_iterator slice_end(const size_type slab,const size_type row,const size_type col,
1567  const slip::Range<int>& range);
1568 
1569 
1596  const_slice_range_iterator slice_begin(const size_type slab,const size_type row,const size_type col,
1597  const slip::Range<int>& range) const;
1598 
1625  const_slice_range_iterator slice_end(const size_type slab,const size_type row,const size_type col,
1626  const slip::Range<int>& range) const;
1627 
1628 
1658  const slip::Range<int>& range);
1659 
1660 
1689  reverse_slice_range_iterator slice_rend(const size_type slab,const size_type row,const size_type col,
1690  const slip::Range<int>& range);
1691 
1692 
1722  const slip::Range<int>& range) const;
1723 
1753  const slip::Range<int>& range) const;
1754 
1755 
1762  //------------------------row range iterators -----------------------
1763 
1791  row_range_iterator row_begin(const size_type slab,const size_type slice,const size_type row,
1792  const slip::Range<int>& range);
1793 
1821  row_range_iterator row_end(const size_type slab,const size_type slice,const size_type row,
1822  const slip::Range<int>& range);
1823 
1824 
1852  const_row_range_iterator row_begin(const size_type slab,const size_type slice,const size_type row,
1853  const slip::Range<int>& range) const;
1854 
1855 
1882  const_row_range_iterator row_end(const size_type slab,const size_type slice,const size_type row,
1883  const slip::Range<int>& range) const;
1884 
1911  reverse_row_range_iterator row_rbegin(const size_type slab,const size_type slice,const size_type row,
1912  const slip::Range<int>& range);
1913 
1914 
1942  reverse_row_range_iterator row_rend(const size_type slab,const size_type slice,const size_type row,
1943  const slip::Range<int>& range);
1944 
1945 
1946 
1972  const_reverse_row_range_iterator row_rbegin(const size_type slab,const size_type slice,const size_type row,
1973  const slip::Range<int>& range) const;
1974 
1975 
2003  const_reverse_row_range_iterator row_rend(const size_type slab,const size_type slice,const size_type row,
2004  const slip::Range<int>& range) const;
2005 
2011  //------------------------col range iterators -----------------------
2012 
2040  col_range_iterator col_begin(const size_type slab,const size_type slice,const size_type col,
2041  const slip::Range<int>& range);
2042 
2071  col_range_iterator col_end(const size_type slab,const size_type slice,const size_type col,
2072  const slip::Range<int>& range);
2073 
2074 
2102  const_col_range_iterator col_begin(const size_type slab,const size_type slice,const size_type col,
2103  const slip::Range<int>& range) const;
2104 
2133  const_col_range_iterator col_end(const size_type slab,const size_type slice,const size_type col,
2134  const slip::Range<int>& range) const;
2135 
2162  reverse_col_range_iterator col_rbegin(const size_type slab,const size_type slice,const size_type col,
2163  const slip::Range<int>& range);
2164 
2192  reverse_col_range_iterator col_rend(const size_type slab,const size_type slice,const size_type col,
2193  const slip::Range<int>& range);
2194 
2221  const_reverse_col_range_iterator col_rbegin(const size_type slab,const size_type slice,const size_type col,
2222  const slip::Range<int>& range) const;
2223 
2250  const_reverse_col_range_iterator col_rend(const size_type slab,const size_type slice,const size_type col,
2251  const slip::Range<int>& range) const;
2252 
2255  //****************************************************************************
2256  // Four dimensional iterators
2257  //****************************************************************************
2258 
2263 
2264  //------------------------ Global iterators------------------------------------
2265 
2283 
2284 
2302 
2303 
2321 
2322 
2340 
2358 
2376 
2394 
2395 
2413 
2420 
2421  //------------------------ Box iterators------------------------------------
2422 
2444 
2445 
2468 
2469 
2492 
2493 
2516 
2517 
2518 
2542 
2566 
2590 
2591 
2615 
2616 
2624 
2625  //------------------------ Range iterators------------------------------------
2626 
2654  iterator4d_range first_front_upper_left(const range & slab_range, const range & slice_range,
2655  const range & row_range, const range & col_range);
2656 
2684  iterator4d_range last_back_bottom_right(const range & slab_range, const range & slice_range,
2685  const range & row_range, const range & col_range);
2686 
2687 
2715  const_iterator4d_range first_front_upper_left(const range & slab_range, const range & slice_range,
2716  const range & row_range, const range & col_range) const;
2717 
2718 
2746  const_iterator4d_range last_back_bottom_right(const range & slab_range, const range & slice_range,
2747  const range & row_range, const range & col_range) const;
2748 
2777  reverse_iterator4d_range rfirst_front_upper_left(const range & slab_range, const range & slice_range,
2778  const range & row_range, const range & col_range);
2779 
2808  reverse_iterator4d_range rlast_back_bottom_right(const range & slab_range, const range & slice_range,
2809  const range & row_range, const range & col_range);
2810 
2839  const_reverse_iterator4d_range rfirst_front_upper_left(const range & slab_range, const range & slice_range,
2840  const range & row_range, const range & col_range) const;
2841 
2870  const_reverse_iterator4d_range rlast_back_bottom_right(const range & slab_range, const range & slice_range,
2871  const range & row_range, const range & col_range) const;
2872 
2873 
2876  //********************************************************************
2877 
2887  friend std::ostream& operator<< <>(std::ostream & out, const HyperVolume<T>& a);
2888 
2889 
2904  self& operator=(const HyperVolume<T> & rhs);
2905 
2906 
2907 
2915  self& operator=(const T& value);
2916 
2917 
2923  void fill(const T& value)
2924  {
2925  std::fill_n(this->begin(),this->size(),value);
2926  }
2927 
2934  void fill(const T* value)
2935  {
2936  std::copy(value,value + this->size(), this->begin());
2937  }
2938 
2947  template<typename InputIterator>
2948  void fill(InputIterator first,
2949  InputIterator last)
2950  {
2951  std::copy(first,last, this->begin());
2952  }
2967  friend bool operator== <>(const HyperVolume<T>& x,
2968  const HyperVolume<T>& y);
2969 
2976  friend bool operator!= <>(const HyperVolume<T>& x,
2977  const HyperVolume<T>& y);
2978 
2985  friend bool operator< <>(const HyperVolume<T>& x,
2986  const HyperVolume<T>& y);
2987 
2994  friend bool operator> <>(const HyperVolume<T>& x,
2995  const HyperVolume<T>& y);
2996 
3003  friend bool operator<= <>(const HyperVolume<T>& x,
3004  const HyperVolume<T>& y);
3005 
3012  friend bool operator>= <>(const HyperVolume<T>& x,
3013  const HyperVolume<T>& y);
3014 
3015 
3025 
3026  T*** operator[](const size_type l);
3027 
3028  const T** const* operator[](const size_type l) const;
3029 
3030 
3046  reference operator()(const size_type l,
3047  const size_type k,
3048  const size_type i,
3049  const size_type j);
3050 
3067  const size_type k,
3068  const size_type i,
3069  const size_type j) const;
3070 
3077  std::string name() const;
3078 
3079 
3084  size_type dim1() const;
3085 
3090  size_type slabs() const;
3091 
3096  size_type dim2() const;
3097 
3102  size_type slices() const;
3103 
3108  size_type dim3() const;
3109 
3114  size_type rows() const;
3115 
3120  size_type dim4() const;
3121 
3126  size_type cols() const;
3127 
3132  size_type columns() const;
3133 
3137  size_type size() const;
3138 
3139 
3143  size_type max_size() const;
3144 
3145 
3149  bool empty()const;
3150 
3160  void swap(HyperVolume& M);
3161 
3171  self& operator+=(const T& val);
3172  self& operator-=(const T& val);
3173  self& operator*=(const T& val);
3174  self& operator/=(const T& val);
3175 
3176 
3177  self operator-() const;
3178 
3179 
3180 
3181  self& operator+=(const self& rhs);
3182  self& operator-=(const self& rhs);
3183  self& operator*=(const self& rhs);
3184  self& operator/=(const self& rhs);
3185 
3186 
3198  T& min() const;
3199 
3200 
3206  T& max() const;
3207 
3212  T sum() const;
3213 
3214 
3221  HyperVolume<T>& apply(T (*fun)(T));
3222 
3229  HyperVolume<T>& apply(T (*fun)(const T&));
3230 
3233 private:
3234  Matrix4d<T>* matrix_;
3235 
3236 private:
3238  template<class Archive>
3239  void save(Archive & ar, const unsigned int version) const
3240  {
3241  ar & this->matrix_;
3242  }
3243  template<class Archive>
3244  void load(Archive & ar, const unsigned int version)
3245  {
3246  ar & this->matrix_;
3247  }
3248  BOOST_SERIALIZATION_SPLIT_MEMBER();
3249 };
3250 
3271 
3272 }//slip::
3273 
3274 
3275 namespace slip{
3286 template<typename T>
3287 HyperVolume<T> operator+(const HyperVolume<T>& M1,
3288  const HyperVolume<T>& M2);
3289 
3296 template<typename T>
3297 HyperVolume<T> operator+(const HyperVolume<T>& M1,
3298  const T& val);
3299 
3306 template<typename T>
3307 HyperVolume<T> operator+(const T& val,
3308  const HyperVolume<T>& M1);
3309 
3320 template<typename T>
3321 HyperVolume<T> operator-(const HyperVolume<T>& M1,
3322  const HyperVolume<T>& M2);
3323 
3330 template<typename T>
3331 HyperVolume<T> operator-(const HyperVolume<T>& M1,
3332  const T& val);
3333 
3340 template<typename T>
3341 HyperVolume<T> operator-(const T& val,
3342  const HyperVolume<T>& M1);
3343 
3354 template<typename T>
3355 HyperVolume<T> operator*(const HyperVolume<T>& M1,
3356  const HyperVolume<T>& M2);
3357 
3364 template<typename T>
3365 HyperVolume<T> operator*(const HyperVolume<T>& M1,
3366  const T& val);
3367 
3374 template<typename T>
3375 HyperVolume<T> operator*(const T& val,
3376  const HyperVolume<T>& M1);
3377 
3388 template<typename T>
3389 HyperVolume<T> operator/(const HyperVolume<T>& M1,
3390  const HyperVolume<T>& M2);
3391 
3398 template<typename T>
3399 HyperVolume<T> operator/(const HyperVolume<T>& M1,
3400  const T& val);
3401 
3402 
3409 template<typename T>
3410 T& min(const HyperVolume<T>& M1);
3411 
3418 template<typename T>
3419 T& max(const HyperVolume<T>& M1);
3420 
3427 template<typename T>
3428 HyperVolume<T> abs(const HyperVolume<T>& V);
3429 
3436 template<typename T>
3437 HyperVolume<T> sqrt(const HyperVolume<T>& V);
3438 
3445 template<typename T>
3446 HyperVolume<T> cos(const HyperVolume<T>& V);
3447 
3454 template<typename T>
3455 HyperVolume<T> acos(const HyperVolume<T>& V);
3456 
3463 template<typename T>
3464 HyperVolume<T> sin(const HyperVolume<T>& V);
3465 
3472 template<typename T>
3473 HyperVolume<T> asin(const HyperVolume<T>& V);
3474 
3475 
3482 template<typename T>
3483 HyperVolume<T> tan(const HyperVolume<T>& V);
3484 
3491 template<typename T>
3492 HyperVolume<T> atan(const HyperVolume<T>& V);
3493 
3500 template<typename T>
3501 HyperVolume<T> exp(const HyperVolume<T>& V);
3502 
3509 template<typename T>
3510 HyperVolume<T> log(const HyperVolume<T>& V);
3511 
3518 template<typename T>
3519 HyperVolume<T> cosh(const HyperVolume<T>& V);
3520 
3527 template<typename T>
3528 HyperVolume<T> sinh(const HyperVolume<T>& V);
3529 
3536 template<typename T>
3537 HyperVolume<T> tanh(const HyperVolume<T>& V);
3538 
3545 template<typename T>
3546 HyperVolume<T> log10(const HyperVolume<T>& V);
3547 
3548 
3549 }//slip::
3550 
3551 
3552 
3553 namespace slip
3554 {
3556 // Constructors
3558 template<typename T>
3559 inline
3561 matrix_(new slip::Matrix4d<T>())
3562 {}
3563 
3564 template<typename T>
3565 inline
3567  const typename HyperVolume<T>::size_type d2,
3568  const typename HyperVolume<T>::size_type d3,
3569  const typename HyperVolume<T>::size_type d4):
3570  matrix_(new slip::Matrix4d<T>(d1,d2,d3,d4))
3571  {}
3572 
3573 template<typename T>
3574 inline
3576  const typename HyperVolume<T>::size_type d2,
3577  const typename HyperVolume<T>::size_type d3,
3578  const typename HyperVolume<T>::size_type d4,
3579  const T& val):
3580  matrix_(new slip::Matrix4d<T>(d1,d2,d3,d4,val))
3581  {}
3582 
3583 template<typename T>
3584 inline
3586  const typename HyperVolume<T>::size_type d2,
3587  const typename HyperVolume<T>::size_type d3,
3588  const typename HyperVolume<T>::size_type d4,
3589  const T* val):
3590  matrix_(new slip::Matrix4d<T>(d1,d2,d3,d4,val))
3591  {}
3592 
3593 
3594 template<typename T>
3595 inline
3597 matrix_(new slip::Matrix4d<T>((*rhs.matrix_)))
3598 {}
3599 
3600 template<typename T>
3601 inline
3603 {
3604  delete matrix_;
3605 }
3606 
3608 
3609 
3611 // Assignment operators
3613 template<typename T>
3614 inline
3616 {
3617  if(this != &rhs)
3618  {
3619  *matrix_ = *(rhs.matrix_);
3620  }
3621  return *this;
3622 }
3623 
3624 template<typename T>
3625 inline
3627 {
3628  std::fill_n((*matrix_)[0][0][0],matrix_->size(),value);
3629  return *this;
3630 }
3631 
3632 template<typename T>
3633 inline
3635  const typename HyperVolume<T>::size_type d2,
3636  const typename HyperVolume<T>::size_type d3,
3637  const typename HyperVolume<T>::size_type d4,
3638  const T& val)
3639  {
3640  matrix_->resize(d1,d2,d3,d4,val);
3641  }
3643 
3644 
3646 // Iterators
3648 
3649 //****************************************************************************
3650 // One dimensional iterators
3651 //****************************************************************************
3652 
3653 
3654 
3655 //----------------------Global iterators------------------------------
3656 
3657 template<typename T>
3658 inline
3660 {
3661  return matrix_->begin();
3662 }
3663 
3664 template<typename T>
3665 inline
3667 {
3668  return matrix_->end();
3669 }
3670 
3671 template<typename T>
3672 inline
3674 {
3675  Matrix4d<T> const * tp(matrix_);
3676  return tp->begin();
3677 }
3678 
3679 template<typename T>
3680 inline
3682 {
3683  Matrix4d<T> const * tp(matrix_);
3684  return tp->end();
3685 }
3686 
3687 
3688 template<typename T>
3689 inline
3691 {
3692  return typename HyperVolume<T>::reverse_iterator(this->end());
3693 }
3694 
3695 template<typename T>
3696 inline
3698 {
3699  return typename HyperVolume<T>::reverse_iterator(this->begin());
3700 }
3701 
3702 template<typename T>
3703 inline
3705 {
3706  return typename HyperVolume<T>::const_reverse_iterator(this->end());
3707 }
3708 
3709 template<typename T>
3710 inline
3712 {
3713  return typename HyperVolume<T>::const_reverse_iterator(this->begin());
3714 }
3715 
3716 //--------------------Constant slice global iterators----------------------
3717 
3718 //--------------------One dimensional slab iterators----------------------
3719 
3720 
3721 template<typename T>
3722 inline
3724  const typename HyperVolume<T>::size_type row, const typename HyperVolume<T>::size_type col){
3725  return matrix_->slab_begin(slice,row,col);
3726 }
3727 
3728 
3729 template<typename T>
3730 inline
3731 typename HyperVolume<T>::const_slab_iterator HyperVolume<T>::slab_begin(const typename HyperVolume<T>::size_type slice,
3732  const typename HyperVolume<T>::size_type row, const typename HyperVolume<T>::size_type col) const{
3733  Matrix4d<T> const * tp(matrix_);
3734  return tp->slab_begin(slice,row,col);
3735 }
3736 
3737 template<typename T>
3738 inline
3739 typename HyperVolume<T>::slab_iterator HyperVolume<T>::slab_end(const typename HyperVolume<T>::size_type slice,
3740  const typename HyperVolume<T>::size_type row, const typename HyperVolume<T>::size_type col){
3741  return matrix_->slab_end(slice,row,col);
3742 }
3743 
3744 
3745 template<typename T>
3746 inline
3747 typename HyperVolume<T>::const_slab_iterator HyperVolume<T>::slab_end(const typename HyperVolume<T>::size_type slice,
3748  const typename HyperVolume<T>::size_type row, const typename HyperVolume<T>::size_type col) const{
3749  Matrix4d<T> const * tp(matrix_);
3750  return tp->slab_end(slice,row,col);
3751 }
3752 
3753 
3754 template<typename T>
3755 inline
3756 typename HyperVolume<T>::reverse_slab_iterator HyperVolume<T>::slab_rbegin(const typename HyperVolume<T>::size_type slice,
3757  const typename HyperVolume<T>::size_type row, const typename HyperVolume<T>::size_type col){
3758  return matrix_->slab_rbegin(slice,row,col);
3759 }
3760 
3761 template<typename T>
3762 inline
3763 typename HyperVolume<T>::const_reverse_slab_iterator HyperVolume<T>::slab_rbegin(const typename HyperVolume<T>::size_type slice,
3764  const typename HyperVolume<T>::size_type row, const typename HyperVolume<T>::size_type col) const{
3765  Matrix4d<T> const * tp(matrix_);
3766  return tp->slab_rbegin(slice,row,col);
3767 }
3768 
3769 
3770 template<typename T>
3771 inline
3772 typename HyperVolume<T>::reverse_slab_iterator HyperVolume<T>::slab_rend(const typename HyperVolume<T>::size_type slice,
3773  const typename HyperVolume<T>::size_type row, const typename HyperVolume<T>::size_type col){
3774  return matrix_->slab_rend(slice,row,col);
3775 }
3776 
3777 template<typename T>
3778 inline
3779 typename HyperVolume<T>::const_reverse_slab_iterator HyperVolume<T>::slab_rend(const typename HyperVolume<T>::size_type slice,
3780  const typename HyperVolume<T>::size_type row, const typename HyperVolume<T>::size_type col) const{
3781  Matrix4d<T> const * tp(matrix_);
3782  return tp->slab_rend(slice,row,col);
3783 }
3784 
3785 
3786 //-------------------- One dimensional slice iterators----------------------
3787 
3788 template<typename T>
3789 inline
3790 typename HyperVolume<T>::slice_iterator
3791 HyperVolume<T>::slice_begin(const typename HyperVolume<T>::size_type slab, const typename HyperVolume<T>::size_type row,
3792  const typename HyperVolume<T>::size_type col)
3793  {
3794  return matrix_->slice_begin(slab,row,col);
3795  }
3796 
3797 template<typename T>
3798 inline
3799 typename HyperVolume<T>::const_slice_iterator
3800 HyperVolume<T>::slice_begin(const typename HyperVolume<T>::size_type slab, const typename HyperVolume<T>::size_type row,
3801  const typename HyperVolume<T>::size_type col) const
3802  {
3803  Matrix4d<T> const * tp(matrix_);
3804  return tp->slice_begin(slab,row,col);
3805  }
3806 
3807 template<typename T>
3808 inline
3809 typename HyperVolume<T>::slice_iterator
3810 HyperVolume<T>::slice_end(const typename HyperVolume<T>::size_type slab, const typename HyperVolume<T>::size_type row,
3811  const typename HyperVolume<T>::size_type col)
3812  {
3813  return matrix_->slice_end(slab,row,col);
3814  }
3815 
3816 template<typename T>
3817 inline
3818 typename HyperVolume<T>::const_slice_iterator
3819 HyperVolume<T>::slice_end(const typename HyperVolume<T>::size_type slab, const typename HyperVolume<T>::size_type row,
3820  const typename HyperVolume<T>::size_type col) const
3821  {
3822  Matrix4d<T> const * tp(matrix_);
3823  return tp->slice_end(slab,row,col);
3824  }
3825 
3826 template<typename T>
3827 inline
3828 typename HyperVolume<T>::reverse_slice_iterator
3829 HyperVolume<T>::slice_rbegin(const typename HyperVolume<T>::size_type slab, const typename HyperVolume<T>::size_type row,
3830  const typename HyperVolume<T>::size_type col)
3831  {
3832  return matrix_->slice_rbegin(slab,row,col);
3833  }
3834 
3835 template<typename T>
3836 inline
3837 typename HyperVolume<T>::const_reverse_slice_iterator
3838 HyperVolume<T>::slice_rbegin(const typename HyperVolume<T>::size_type slab, const typename HyperVolume<T>::size_type row,
3839  const typename HyperVolume<T>::size_type col) const
3840  {
3841  Matrix4d<T> const * tp(matrix_);
3842  return tp->slice_rbegin(slab,row,col);
3843  }
3844 
3845 template<typename T>
3846 inline
3847 typename HyperVolume<T>::reverse_slice_iterator
3848 HyperVolume<T>::slice_rend(const typename HyperVolume<T>::size_type slab, const typename HyperVolume<T>::size_type row,
3849  const typename HyperVolume<T>::size_type col)
3850  {
3851  return matrix_->slice_rend(slab,row,col);
3852  }
3853 
3854 template<typename T>
3855 inline
3856 typename HyperVolume<T>::const_reverse_slice_iterator
3857 HyperVolume<T>::slice_rend(const typename HyperVolume<T>::size_type slab, const typename HyperVolume<T>::size_type row,
3858  const typename HyperVolume<T>::size_type col) const
3859  {
3860  Matrix4d<T> const * tp(matrix_);
3861  return tp->slice_rend(slab,row,col);
3862  }
3863 
3864 //--------------------One dimensional row iterators----------------------
3865 
3866 template<typename T>
3867 inline
3868 typename HyperVolume<T>::row_iterator
3869 HyperVolume<T>::row_begin(const typename HyperVolume<T>::size_type slab, const typename HyperVolume<T>::size_type slice,
3870  const typename HyperVolume<T>::size_type row)
3871  {
3872  return matrix_->row_begin(slab,slice,row);
3873  }
3874 
3875 template<typename T>
3876 inline
3877 typename HyperVolume<T>::const_row_iterator
3878 HyperVolume<T>::row_begin(const typename HyperVolume<T>::size_type slab, const typename HyperVolume<T>::size_type slice,
3879  const typename HyperVolume<T>::size_type row) const
3880  {
3881  Matrix4d<T> const * tp(matrix_);
3882  return tp->row_begin(slab,slice,row);
3883  }
3884 
3885 template<typename T>
3886 inline
3887 typename HyperVolume<T>::row_iterator
3888 HyperVolume<T>::row_end(const typename HyperVolume<T>::size_type slab, const typename HyperVolume<T>::size_type slice,
3889  const typename HyperVolume<T>::size_type row)
3890  {
3891  return matrix_->row_end(slab,slice,row);
3892  }
3893 
3894 template<typename T>
3895 inline
3896 typename HyperVolume<T>::const_row_iterator
3897 HyperVolume<T>::row_end(const typename HyperVolume<T>::size_type slab, const typename HyperVolume<T>::size_type slice,
3898  const typename HyperVolume<T>::size_type row) const
3899  {
3900  Matrix4d<T> const * tp(matrix_);
3901  return tp->row_end(slab,slice,row);
3902  }
3903 
3904 
3905 template<typename T>
3906 inline
3907 typename HyperVolume<T>::reverse_row_iterator
3908 HyperVolume<T>::row_rbegin(const typename HyperVolume<T>::size_type slab, const typename HyperVolume<T>::size_type slice,
3909  const typename HyperVolume<T>::size_type row)
3910  {
3911  return matrix_->row_rbegin(slab,slice,row);
3912  }
3913 
3914 template<typename T>
3915 inline
3916 typename HyperVolume<T>::const_reverse_row_iterator
3917 HyperVolume<T>::row_rbegin(const typename HyperVolume<T>::size_type slab, const typename HyperVolume<T>::size_type slice,
3918  const typename HyperVolume<T>::size_type row) const
3919  {
3920  Matrix4d<T> const * tp(matrix_);
3921  return tp->row_rbegin(slab,slice,row);
3922  }
3923 
3924 
3925 template<typename T>
3926 inline
3927 typename HyperVolume<T>::reverse_row_iterator
3928 HyperVolume<T>::row_rend(const typename HyperVolume<T>::size_type slab, const typename HyperVolume<T>::size_type slice,
3929  const typename HyperVolume<T>::size_type row)
3930  {
3931  return matrix_->row_rend(slab,slice,row);
3932  }
3933 
3934 template<typename T>
3935 inline
3936 typename HyperVolume<T>::const_reverse_row_iterator
3937 HyperVolume<T>::row_rend(const typename HyperVolume<T>::size_type slab, const typename HyperVolume<T>::size_type slice,
3938  const typename HyperVolume<T>::size_type row) const
3939  {
3940  Matrix4d<T> const * tp(matrix_);
3941  return tp->row_rend(slab,slice,row);
3942  }
3943 
3944 //--------------------One dimensional col iterators----------------------
3945 
3946 template<typename T>
3947 inline
3948 typename HyperVolume<T>::col_iterator
3949 HyperVolume<T>::col_begin(const typename HyperVolume<T>::size_type slab, const typename HyperVolume<T>::size_type slice,
3950  const typename HyperVolume<T>::size_type col)
3951  {
3952  return matrix_->col_begin(slab,slice,col);
3953  }
3954 
3955 template<typename T>
3956 inline
3957 typename HyperVolume<T>::const_col_iterator
3958 HyperVolume<T>::col_begin(const typename HyperVolume<T>::size_type slab, const typename HyperVolume<T>::size_type slice,
3959  const typename HyperVolume<T>::size_type col) const
3960  {
3961  Matrix4d<T> const * tp(matrix_);
3962  return tp->col_begin(slab,slice,col);
3963  }
3964 
3965 
3966 template<typename T>
3967 inline
3968 typename HyperVolume<T>::col_iterator
3969 HyperVolume<T>::col_end(const typename HyperVolume<T>::size_type slab, const typename HyperVolume<T>::size_type slice,
3970  const typename HyperVolume<T>::size_type col)
3971  {
3972  return matrix_->col_end(slab,slice,col);
3973  }
3974 
3975 template<typename T>
3976 inline
3977 typename HyperVolume<T>::const_col_iterator
3978 HyperVolume<T>::col_end(const typename HyperVolume<T>::size_type slab, const typename HyperVolume<T>::size_type slice,
3979  const typename HyperVolume<T>::size_type col) const
3980  {
3981  Matrix4d<T> const * tp(matrix_);
3982  return tp->col_end(slab,slice,col);
3983  }
3984 
3985 template<typename T>
3986 inline
3987 typename HyperVolume<T>::reverse_col_iterator
3988 HyperVolume<T>::col_rbegin(const typename HyperVolume<T>::size_type slab, const typename HyperVolume<T>::size_type slice,
3989  const typename HyperVolume<T>::size_type col)
3990  {
3991  return matrix_->col_rbegin(slab,slice,col);
3992  }
3993 
3994 template<typename T>
3995 inline
3996 typename HyperVolume<T>::const_reverse_col_iterator
3997 HyperVolume<T>::col_rbegin(const typename HyperVolume<T>::size_type slab, const typename HyperVolume<T>::size_type slice,
3998  const typename HyperVolume<T>::size_type col) const
3999  {
4000  Matrix4d<T> const * tp(matrix_);
4001  return tp->col_rbegin(slab,slice,col);
4002  }
4003 
4004 
4005 template<typename T>
4006 inline
4007 typename HyperVolume<T>::reverse_col_iterator
4008 HyperVolume<T>::col_rend(const typename HyperVolume<T>::size_type slab, const typename HyperVolume<T>::size_type slice,
4009  const typename HyperVolume<T>::size_type col)
4010  {
4011  return matrix_->col_rend(slab,slice,col);
4012  }
4013 
4014 template<typename T>
4015 inline
4016 typename HyperVolume<T>::const_reverse_col_iterator
4017 HyperVolume<T>::col_rend(const typename HyperVolume<T>::size_type slab, const typename HyperVolume<T>::size_type slice,
4018  const typename HyperVolume<T>::size_type col) const
4019  {
4020  Matrix4d<T> const * tp(matrix_);
4021  return tp->col_rend(slab,slice,col);
4022  }
4023 
4024 //--------------------slab range iterators----------------------
4025 
4026 template<typename T>
4027 inline
4028 typename HyperVolume<T>::slab_range_iterator HyperVolume<T>::slab_begin(const typename HyperVolume<T>::size_type slice,
4029  const typename HyperVolume<T>::size_type row,const typename HyperVolume<T>::size_type col,
4030  const slip::Range<int>& range){
4031  return matrix_->slab_begin(slice,row,col,range);
4032 }
4033 
4034 
4035 template<typename T>
4036 inline
4037 typename HyperVolume<T>::slab_range_iterator HyperVolume<T>::slab_end(const typename HyperVolume<T>::size_type slice,
4038  const typename HyperVolume<T>::size_type row,const typename HyperVolume<T>::size_type col,
4039  const slip::Range<int>& range){
4040  return matrix_->slab_end(slice,row,col,range);
4041 }
4042 
4043 
4044 template<typename T>
4045 inline
4046 typename HyperVolume<T>::const_slab_range_iterator HyperVolume<T>::slab_begin(const typename HyperVolume<T>::size_type slice,
4047  const typename HyperVolume<T>::size_type row,const typename HyperVolume<T>::size_type col,
4048  const slip::Range<int>& range) const{
4049  Matrix4d<T> const * tp(matrix_);
4050  return tp->slab_begin(slice,row,col,range);
4051 }
4052 
4053 template<typename T>
4054 inline
4055 typename HyperVolume<T>::const_slab_range_iterator HyperVolume<T>::slab_end(const typename HyperVolume<T>::size_type slice,
4056  const typename HyperVolume<T>::size_type row,const typename HyperVolume<T>::size_type col,
4057  const slip::Range<int>& range) const{
4058  Matrix4d<T> const * tp(matrix_);
4059  return tp->slab_end(slice,row,col,range);
4060 }
4061 
4062 template<typename T>
4063 inline
4064 typename HyperVolume<T>::reverse_slab_range_iterator HyperVolume<T>::slab_rbegin(const typename HyperVolume<T>::size_type slice,
4065  const typename HyperVolume<T>::size_type row,const typename HyperVolume<T>::size_type col,
4066  const slip::Range<int>& range){
4067  return matrix_->slab_rbegin(slice,row,col,range);
4068 }
4069 
4070 
4071 template<typename T>
4072 inline
4073 typename HyperVolume<T>::reverse_slab_range_iterator HyperVolume<T>::slab_rend(const HyperVolume<T>::size_type slice,
4074  const HyperVolume<T>::size_type row, const HyperVolume<T>::size_type col,
4075  const slip::Range<int>& range){
4076  return matrix_->slab_rend(slice,row,col,range);
4077 }
4078 
4079 
4080 template<typename T>
4081 inline
4082 typename HyperVolume<T>::const_reverse_slab_range_iterator HyperVolume<T>::slab_rbegin(const typename HyperVolume<T>::size_type slice,
4083  const typename HyperVolume<T>::size_type row,const typename HyperVolume<T>::size_type col,
4084  const slip::Range<int>& range) const{
4085  Matrix4d<T> const * tp(matrix_);
4086  return tp->slab_rbegin(slice,row,col,range);
4087 }
4088 
4089 template<typename T>
4090 inline
4091 typename HyperVolume<T>::const_reverse_slab_range_iterator HyperVolume<T>::slab_rend(const typename HyperVolume<T>::size_type slice,
4092  const typename HyperVolume<T>::size_type row,const typename HyperVolume<T>::size_type col,
4093  const slip::Range<int>& range) const{
4094  Matrix4d<T> const * tp(matrix_);
4095  return tp->slab_rend(slice,row,col,range);
4096 }
4097 
4098 
4099 //--------------------Constant slice range iterators----------------------
4100 
4101 template<typename T>
4102 inline
4103 typename HyperVolume<T>::slice_range_iterator
4104 HyperVolume<T>::slice_begin(const typename HyperVolume<T>::size_type slab,
4105  const typename HyperVolume<T>::size_type row, const typename HyperVolume<T>::size_type col,
4106  const slip::Range<int>& range)
4107  {
4108  return matrix_->slice_begin(slab,row,col,range);
4109  }
4110 
4111 template<typename T>
4112 inline
4113 typename HyperVolume<T>::const_slice_range_iterator
4114 HyperVolume<T>::slice_begin(const typename HyperVolume<T>::size_type slab,
4115  const typename HyperVolume<T>::size_type row, const typename HyperVolume<T>::size_type col,
4116  const slip::Range<int>& range) const
4117  {
4118  Matrix4d<T> const * tp(matrix_);
4119  return tp->slice_begin(slab,row,col,range);
4120  }
4121 
4122 template<typename T>
4123 inline
4124 typename HyperVolume<T>::slice_range_iterator
4125 HyperVolume<T>::slice_end(const typename HyperVolume<T>::size_type slab,
4126  const typename HyperVolume<T>::size_type row, const typename HyperVolume<T>::size_type col,
4127  const slip::Range<int>& range)
4128  {
4129  return matrix_->slice_end(slab,row,col,range);
4130  }
4131 
4132 template<typename T>
4133 inline
4134 typename HyperVolume<T>::const_slice_range_iterator
4135 HyperVolume<T>::slice_end(const typename HyperVolume<T>::size_type slab,
4136  const typename HyperVolume<T>::size_type row, const typename HyperVolume<T>::size_type col,
4137  const slip::Range<int>& range) const
4138  {
4139  Matrix4d<T> const * tp(matrix_);
4140  return tp->slice_end(slab,row,col,range);
4141  }
4142 
4143 template<typename T>
4144 inline
4145 typename HyperVolume<T>::reverse_slice_range_iterator
4146 HyperVolume<T>::slice_rbegin(const typename HyperVolume<T>::size_type slab,
4147  const typename HyperVolume<T>::size_type row, const typename HyperVolume<T>::size_type col,
4148  const slip::Range<int>& range)
4149  {
4150  return matrix_->slice_rbegin(slab,row,col,range);
4151  }
4152 
4153 template<typename T>
4154 inline
4155 typename HyperVolume<T>::const_reverse_slice_range_iterator
4156 HyperVolume<T>::slice_rbegin(const typename HyperVolume<T>::size_type slab,
4157  const typename HyperVolume<T>::size_type row, const typename HyperVolume<T>::size_type col,
4158  const slip::Range<int>& range) const
4159  {
4160  Matrix4d<T> const * tp(matrix_);
4161  return tp->slice_rbegin(slab,row,col,range);
4162  }
4163 
4164 template<typename T>
4165 inline
4166 typename HyperVolume<T>::reverse_slice_range_iterator
4167 HyperVolume<T>::slice_rend(const typename HyperVolume<T>::size_type slab,
4168  const typename HyperVolume<T>::size_type row, const typename HyperVolume<T>::size_type col,
4169  const slip::Range<int>& range)
4170  {
4171  return matrix_->slice_rend(slab,row,col,range);
4172  }
4173 
4174 template<typename T>
4175 inline
4176 typename HyperVolume<T>::const_reverse_slice_range_iterator
4177 HyperVolume<T>::slice_rend(const typename HyperVolume<T>::size_type slab,
4178  const typename HyperVolume<T>::size_type row, const typename HyperVolume<T>::size_type col,
4179  const slip::Range<int>& range) const
4180  {
4181  Matrix4d<T> const * tp(matrix_);
4182  return tp->slice_rend(slab,row,col,range);
4183  }
4184 
4185 //--------------------Constant row range iterators----------------------
4186 
4187 template<typename T>
4188 inline
4189 typename HyperVolume<T>::row_range_iterator
4190 HyperVolume<T>::row_begin(const typename HyperVolume<T>::size_type slab,
4191  const typename HyperVolume<T>::size_type slice, const typename HyperVolume<T>::size_type row,const slip::Range<int>& range)
4192  {
4193  return matrix_->row_begin(slab,slice,row,range);
4194  }
4195 
4196 template<typename T>
4197 inline
4198 typename HyperVolume<T>::const_row_range_iterator
4199 HyperVolume<T>::row_begin(const typename HyperVolume<T>::size_type slab,
4200  const typename HyperVolume<T>::size_type slice,const typename HyperVolume<T>::size_type row,
4201  const slip::Range<int>& range) const
4202  {
4203  Matrix4d<T> const * tp(matrix_);
4204  return tp->row_begin(slab,slice,row,range);
4205  }
4206 
4207 template<typename T>
4208 inline
4209 typename HyperVolume<T>::row_range_iterator
4210 HyperVolume<T>::row_end(const typename HyperVolume<T>::size_type slab,
4211  const typename HyperVolume<T>::size_type slice,const typename HyperVolume<T>::size_type row,
4212  const slip::Range<int>& range)
4213  {
4214  return matrix_->row_end(slab,slice,row,range);
4215  }
4216 
4217 template<typename T>
4218 inline
4219 typename HyperVolume<T>::const_row_range_iterator
4220 HyperVolume<T>::row_end(const typename HyperVolume<T>::size_type slab,
4221  const typename HyperVolume<T>::size_type slice,const typename HyperVolume<T>::size_type row,
4222  const slip::Range<int>& range) const
4223  {
4224  Matrix4d<T> const * tp(matrix_);
4225  return tp->row_end(slab,slice,row,range);
4226  }
4227 
4228 template<typename T>
4229 inline
4230 typename HyperVolume<T>::reverse_row_range_iterator
4231 HyperVolume<T>::row_rbegin(const typename HyperVolume<T>::size_type slab,
4232  const typename HyperVolume<T>::size_type slice,const typename HyperVolume<T>::size_type row,
4233  const slip::Range<int>& range)
4234  {
4235  return matrix_->row_rbegin(slab,slice,row,range);
4236  }
4237 
4238 template<typename T>
4239 inline
4240 typename HyperVolume<T>::const_reverse_row_range_iterator
4241 HyperVolume<T>::row_rbegin(const typename HyperVolume<T>::size_type slab,
4242  const typename HyperVolume<T>::size_type slice,const typename HyperVolume<T>::size_type row,
4243  const slip::Range<int>& range) const
4244  {
4245  Matrix4d<T> const * tp(matrix_);
4246  return tp->row_rbegin(slab,slice,row,range);
4247  }
4248 
4249 template<typename T>
4250 inline
4251 typename HyperVolume<T>::reverse_row_range_iterator
4252 HyperVolume<T>::row_rend(const typename HyperVolume<T>::size_type slab,
4253  const typename HyperVolume<T>::size_type slice,const typename HyperVolume<T>::size_type row,
4254  const slip::Range<int>& range)
4255  {
4256  return matrix_->row_rend(slab,slice,row,range);
4257  }
4258 
4259 template<typename T>
4260 inline
4261 typename HyperVolume<T>::const_reverse_row_range_iterator
4262 HyperVolume<T>::row_rend(const typename HyperVolume<T>::size_type slab,
4263  const typename HyperVolume<T>::size_type slice,const typename HyperVolume<T>::size_type row,
4264  const slip::Range<int>& range) const
4265  {
4266  Matrix4d<T> const * tp(matrix_);
4267  return tp->row_rend(slab,slice,row,range);
4268  }
4269 
4270 //--------------------Constant col range iterators----------------------
4271 
4272 template<typename T>
4273 inline
4274 typename HyperVolume<T>::col_range_iterator
4275 HyperVolume<T>::col_begin(const typename HyperVolume<T>::size_type slab,
4276  const typename HyperVolume<T>::size_type slice,const typename HyperVolume<T>::size_type col,
4277  const slip::Range<int>& range)
4278  {
4279  return matrix_->col_begin(slab,slice,col,range);
4280  }
4281 
4282 template<typename T>
4283 inline
4284 typename HyperVolume<T>::const_col_range_iterator
4285 HyperVolume<T>::col_begin(const typename HyperVolume<T>::size_type slab,
4286  const typename HyperVolume<T>::size_type slice,const typename HyperVolume<T>::size_type col,
4287  const slip::Range<int>& range) const
4288  {
4289  Matrix4d<T> const * tp(matrix_);
4290  return tp->col_begin(slab,slice,col,range);
4291  }
4292 
4293 template<typename T>
4294 inline
4295 typename HyperVolume<T>::col_range_iterator
4296 HyperVolume<T>::col_end(const typename HyperVolume<T>::size_type slab,
4297  const typename HyperVolume<T>::size_type slice,const typename HyperVolume<T>::size_type col,
4298  const slip::Range<int>& range)
4299  {
4300  return matrix_->col_end(slab,slice,col,range);
4301  }
4302 
4303 template<typename T>
4304 inline
4305 typename HyperVolume<T>::const_col_range_iterator
4306 HyperVolume<T>::col_end(const typename HyperVolume<T>::size_type slab,
4307  const typename HyperVolume<T>::size_type slice,const typename HyperVolume<T>::size_type col,
4308  const slip::Range<int>& range) const
4309  {
4310  Matrix4d<T> const * tp(matrix_);
4311  return tp->col_end(slab,slice,col,range);
4312  }
4313 
4314 template<typename T>
4315 inline
4316 typename HyperVolume<T>::reverse_col_range_iterator
4317 HyperVolume<T>::col_rbegin(const typename HyperVolume<T>::size_type slab,
4318  const typename HyperVolume<T>::size_type slice,const typename HyperVolume<T>::size_type col,
4319  const slip::Range<int>& range)
4320  {
4321  return matrix_->col_rbegin(slab,slice,col,range);
4322  }
4323 
4324 template<typename T>
4325 inline
4326 typename HyperVolume<T>::const_reverse_col_range_iterator
4327 HyperVolume<T>::col_rbegin(const typename HyperVolume<T>::size_type slab,
4328  const typename HyperVolume<T>::size_type slice,const typename HyperVolume<T>::size_type col,
4329  const slip::Range<int>& range) const
4330  {
4331  Matrix4d<T> const * tp(matrix_);
4332  return tp->col_rbegin(slab,slice,col,range);
4333  }
4334 
4335 template<typename T>
4336 inline
4337 typename HyperVolume<T>::reverse_col_range_iterator
4338 HyperVolume<T>::col_rend(const typename HyperVolume<T>::size_type slab,
4339  const typename HyperVolume<T>::size_type slice,const typename HyperVolume<T>::size_type col,
4340  const slip::Range<int>& range)
4341  {
4342  return matrix_->col_rend(slab,slice,col,range);
4343  }
4344 
4345 template<typename T>
4346 inline
4347 typename HyperVolume<T>::const_reverse_col_range_iterator
4348 HyperVolume<T>::col_rend(const typename HyperVolume<T>::size_type slab,
4349  const typename HyperVolume<T>::size_type slice,const typename HyperVolume<T>::size_type col,
4350  const slip::Range<int>& range) const
4351  {
4352  Matrix4d<T> const * tp(matrix_);
4353  return tp->col_rend(slab,slice,col,range);
4354  }
4355 
4356 
4357 //****************************************************************************
4358 // Four dimensional iterators
4359 //****************************************************************************
4360 
4361 //------------------------ Global iterators------------------------------------
4362 
4363 template<typename T>
4364 inline
4366 {
4367  return matrix_->first_front_upper_left();
4368 }
4369 
4370 template<typename T>
4371 inline
4373 {
4374  Matrix4d<T> const * tp(matrix_);
4375  return tp->first_front_upper_left();
4376 }
4377 
4378 
4379 template<typename T>
4380 inline
4382 {
4383  return matrix_->last_back_bottom_right();
4384 }
4385 
4386 template<typename T>
4387 inline
4389 {
4390  Matrix4d<T> const * tp(matrix_);
4391  return tp->last_back_bottom_right();
4392 }
4393 
4394 template<typename T>
4395 inline
4398 {
4399  return matrix_->rlast_back_bottom_right();
4400 }
4401 
4402 template<typename T>
4403 inline
4406 {
4407  Matrix4d<T> const * tp(matrix_);
4408  return tp->rlast_back_bottom_right();
4409 }
4410 
4411 template<typename T>
4412 inline
4415 {
4416  return matrix_->rfirst_front_upper_left();
4417 }
4418 
4419 template<typename T>
4420 inline
4423 {
4424  Matrix4d<T> const * tp(matrix_);
4425  return tp->rfirst_front_upper_left();
4426 }
4427 
4428 //------------------------ Box iterators------------------------------------
4429 
4430 template<typename T>
4431 inline
4433 {
4434  return matrix_->first_front_upper_left(box);
4435 }
4436 
4437 template<typename T>
4438 inline
4440 {
4441  Matrix4d<T> const * tp(matrix_);
4442  return tp->first_front_upper_left(box);
4443 }
4444 
4445 
4446 template<typename T>
4447 inline
4450 {
4451  return matrix_->last_back_bottom_right(box);
4452 }
4453 
4454 template<typename T>
4455 inline
4458 {
4459  Matrix4d<T> const * tp(matrix_);
4460  return tp->last_back_bottom_right(box);
4461 }
4462 
4463 template<typename T>
4464 inline
4467 {
4468  return matrix_->rlast_back_bottom_right(box);
4469 }
4470 
4471 template<typename T>
4472 inline
4475 {
4476  Matrix4d<T> const * tp(matrix_);
4477  return tp->rlast_back_bottom_right(box);
4478 }
4479 
4480 template<typename T>
4481 inline
4484 {
4485  return matrix_->rfirst_front_upper_left(box);
4486 }
4487 
4488 template<typename T>
4489 inline
4492 {
4493  Matrix4d<T> const * tp(matrix_);
4494  return tp->rfirst_front_upper_left(box);
4495 }
4496 
4497 //------------------------ Range iterators------------------------------------
4498 
4499 template<typename T>
4500 inline
4503  const typename HyperVolume<T>::range& slice_range,
4504  const typename HyperVolume<T>::range& row_range,
4505  const typename HyperVolume<T>::range& col_range)
4506  {
4507  return matrix_->first_front_upper_left(slab_range,slice_range,row_range,col_range);
4508  }
4509 
4510 template<typename T>
4511 inline
4512 typename HyperVolume<T>::iterator4d_range
4513 HyperVolume<T>::last_back_bottom_right(const typename HyperVolume<T>::range& slab_range,
4514  const typename HyperVolume<T>::range& slice_range,
4515  const typename HyperVolume<T>::range& row_range,
4516  const typename HyperVolume<T>::range& col_range)
4517  {
4518  return matrix_->last_back_bottom_right(slab_range,slice_range,row_range,col_range);
4519  }
4520 
4521 
4522 template<typename T>
4523 inline
4524 typename HyperVolume<T>::const_iterator4d_range
4525 HyperVolume<T>::first_front_upper_left(const typename HyperVolume<T>::range& slab_range,
4526  const typename HyperVolume<T>::range& slice_range,
4527  const typename HyperVolume<T>::range& row_range,
4528  const typename HyperVolume<T>::range& col_range) const
4529  {
4530  Matrix4d<T> const * tp(matrix_);
4531  return tp->first_front_upper_left(slab_range,slice_range,row_range,col_range);
4532  }
4533 
4534 
4535 template<typename T>
4536 inline
4537 typename HyperVolume<T>::const_iterator4d_range
4538 HyperVolume<T>::last_back_bottom_right(const typename HyperVolume<T>::range& slab_range,
4539  const typename HyperVolume<T>::range& slice_range,
4540  const typename HyperVolume<T>::range& row_range,
4541  const typename HyperVolume<T>::range& col_range) const
4542  {
4543  Matrix4d<T> const * tp(matrix_);
4544  return tp->last_back_bottom_right(slab_range,slice_range,row_range,col_range);
4545  }
4546 
4547 template<typename T>
4548 inline
4549 typename HyperVolume<T>::reverse_iterator4d_range
4550 HyperVolume<T>::rfirst_front_upper_left(const typename HyperVolume<T>::range& slab_range,
4551  const typename HyperVolume<T>::range& slice_range,
4552  const typename HyperVolume<T>::range& row_range,
4553  const typename HyperVolume<T>::range& col_range)
4554  {
4555  return matrix_->rfirst_front_upper_left(slab_range,slice_range,row_range,col_range);
4556  }
4557 
4558 template<typename T>
4559 inline
4560 typename HyperVolume<T>::const_reverse_iterator4d_range
4561 HyperVolume<T>::rfirst_front_upper_left(const typename HyperVolume<T>::range& slab_range,
4562  const typename HyperVolume<T>::range& slice_range,
4563  const typename HyperVolume<T>::range& row_range,
4564  const typename HyperVolume<T>::range& col_range) const
4565  {
4566  Matrix4d<T> const * tp(matrix_);
4567  return tp->rfirst_front_upper_left(slab_range,slice_range,row_range,col_range);
4568  }
4569 
4570 template<typename T>
4571 inline
4572 typename HyperVolume<T>::reverse_iterator4d_range
4573 HyperVolume<T>::rlast_back_bottom_right(const typename HyperVolume<T>::range& slab_range,
4574  const typename HyperVolume<T>::range& slice_range,
4575  const typename HyperVolume<T>::range& row_range,
4576  const typename HyperVolume<T>::range& col_range)
4577  {
4578  return matrix_->rlast_back_bottom_right(slab_range,slice_range,row_range,col_range);
4579  }
4580 
4581 template<typename T>
4582 inline
4583 typename HyperVolume<T>::const_reverse_iterator4d_range
4584 HyperVolume<T>::rlast_back_bottom_right(const typename HyperVolume<T>::range& slab_range,
4585  const typename HyperVolume<T>::range& slice_range,
4586  const typename HyperVolume<T>::range& row_range,
4587  const typename HyperVolume<T>::range& col_range) const
4588  {
4589  Matrix4d<T> const * tp(matrix_);
4590  return tp->rlast_back_bottom_right(slab_range,slice_range,row_range,col_range);
4591  }
4592 
4593 
4595 
4596 
4597 template <typename T>
4598 inline
4599 std::ostream& operator<<(std::ostream & out, const HyperVolume<T>& a)
4600 {
4601  out<<*(a.matrix_);
4602  return out;
4603 }
4604 
4605 
4607 
4608 
4610 // Elements access operators
4612 template<typename T>
4613 inline
4614 T***
4615 HyperVolume<T>::operator[](const typename HyperVolume<T>::size_type l)
4616 {
4617  return (*matrix_)[l];
4618 }
4619 
4620 template<typename T>
4621 inline
4622 const T** const*
4623 HyperVolume<T>::operator[](const typename HyperVolume<T>::size_type l) const
4624 {
4625  Matrix4d<T> const * tp(matrix_);
4626  return tp->operator[](l);
4627 }
4628 
4629 template<typename T>
4630 inline
4631 typename HyperVolume<T>::reference
4632 HyperVolume<T>::operator()(const typename HyperVolume<T>::size_type l,
4633  const typename HyperVolume<T>::size_type k,
4634  const typename HyperVolume<T>::size_type i,
4635  const typename HyperVolume<T>::size_type j)
4636 {
4637  return (*matrix_)[l][k][i][j];
4638 }
4639 
4640 template<typename T>
4641 inline
4642 typename HyperVolume<T>::const_reference
4643 HyperVolume<T>::operator()(const typename HyperVolume<T>::size_type l,
4644  const typename HyperVolume<T>::size_type k,
4645  const typename HyperVolume<T>::size_type i,
4646  const typename HyperVolume<T>::size_type j) const
4647 {
4648  return (*matrix_)[l][k][i][j];
4649 }
4650 
4652 
4653 template<typename T>
4654 inline
4655 std::string
4656 HyperVolume<T>::name() const {return "HyperVolume";}
4657 
4658 template<typename T>
4659 inline
4661 HyperVolume<T>::dim1() const {return matrix_->dim1();}
4662 
4663 template<typename T>
4664 inline
4666 HyperVolume<T>::slabs() const {return matrix_->dim1();}
4667 
4668 template<typename T>
4669 inline
4671 HyperVolume<T>::dim2() const {return matrix_->dim2();}
4672 
4673 template<typename T>
4674 inline
4676 HyperVolume<T>::slices() const {return matrix_->dim2();}
4677 
4678 template<typename T>
4679 inline
4681 HyperVolume<T>::dim3() const {return matrix_->dim3();}
4682 
4683 template<typename T>
4684 inline
4686 HyperVolume<T>::rows() const {return matrix_->dim3();}
4687 
4688 template<typename T>
4689 inline
4691 HyperVolume<T>::dim4() const {return matrix_->dim4();}
4692 
4693 template<typename T>
4694 inline
4696 HyperVolume<T>::cols() const {return matrix_->dim4();}
4697 
4698 template<typename T>
4699 inline
4701 HyperVolume<T>::columns() const {return matrix_->dim4();}
4702 
4703 
4704 template<typename T>
4705 inline
4707 HyperVolume<T>::size() const {return matrix_->size();}
4708 
4709 template<typename T>
4710 inline
4712 HyperVolume<T>::max_size() const {return matrix_->max_size();}
4713 
4714 template<typename T>
4715 inline
4716 bool HyperVolume<T>::empty()const {return matrix_->empty();}
4717 
4718 template<typename T>
4719 inline
4721 {
4722  matrix_->swap(*(M.matrix_));
4723 }
4724 
4725 
4727 /* @{ */
4728  template<typename T>
4729 inline
4731  const HyperVolume<T>& y)
4732  {
4733  return ( x.size() == y.size()
4734  && std::equal(x.begin(),x.end(),y.begin()));
4735  }
4736 
4737  template<typename T>
4738  inline
4740  const HyperVolume<T>& y)
4741  {
4742  return !(x == y);
4743  }
4744  /* @} */
4745 
4746 
4748  /* @{ */
4749  template<typename T>
4750  inline
4751  bool operator<(const HyperVolume<T>& x,
4752  const HyperVolume<T>& y)
4753  {
4754  return std::lexicographical_compare(x.begin(), x.end(),
4755  y.begin(), y.end());
4756  }
4757 
4758 
4759  template<typename T>
4760  inline
4761  bool operator>(const HyperVolume<T>& x,
4762  const HyperVolume<T>& y)
4763  {
4764  return (y < x);
4765  }
4766 
4767  template<typename T>
4768  inline
4769  bool operator<=(const HyperVolume<T>& x,
4770  const HyperVolume<T>& y)
4771  {
4772  return !(y < x);
4773  }
4774 
4775  template<typename T>
4776  inline
4778  const HyperVolume<T>& y)
4779  {
4780  return !(x < y);
4781  }
4782  /* @} */
4784 
4785 
4787  // arithmetic and mathematic operators
4789  template<typename T>
4790  inline
4792  {
4793  std::transform(matrix_->begin(),matrix_->end(),matrix_->begin(),std::bind2nd(std::plus<T>(),val));
4794  return *this;
4795  }
4796 
4797  template<typename T>
4798  inline
4800  {
4801  std::transform(matrix_->begin(),matrix_->end(),matrix_->begin(),std::bind2nd(std::minus<T>(),val));
4802  return *this;
4803  }
4804 
4805  template<typename T>
4806  inline
4808  {
4809  std::transform(matrix_->begin(),matrix_->end(),matrix_->begin(),std::bind2nd(std::multiplies<T>(),val));
4810  return *this;
4811  }
4812 
4813  template<typename T>
4814  inline
4816  {
4817  std::transform(matrix_->begin(),matrix_->end(),matrix_->begin(),std::bind2nd(std::divides<T>(),val));
4818  return *this;
4819  }
4820 
4821  template<typename T>
4822  inline
4824  {
4825  HyperVolume<T> tmp(*this);
4826  std::transform(matrix_->begin(),matrix_->end(),tmp.begin(),std::negate<T>());
4827  return tmp;
4828  }
4829 
4830  template<typename T>
4831  inline
4833  {
4834  assert(this->dim1() == rhs.dim1());
4835  assert(this->dim2() == rhs.dim2());
4836  assert(this->dim3() == rhs.dim3());
4837  assert(this->dim4() == rhs.dim4());
4838 
4839  std::transform(matrix_->begin(),matrix_->end(),(*(rhs.matrix_)).begin(),matrix_->begin(),std::plus<T>());
4840  return *this;
4841  }
4842 
4843  template<typename T>
4844  inline
4846  {
4847  assert(this->dim1() == rhs.dim1());
4848  assert(this->dim2() == rhs.dim2());
4849  assert(this->dim3() == rhs.dim3());
4850  assert(this->dim4() == rhs.dim4());
4851  std::transform(matrix_->begin(),matrix_->end(),(*(rhs.matrix_)).begin(),matrix_->begin(),std::minus<T>());
4852  return *this;
4853  }
4854 
4855  template<typename T>
4856  inline
4858  {
4859  assert(this->dim1() == rhs.dim1());
4860  assert(this->dim2() == rhs.dim2());
4861  assert(this->dim3() == rhs.dim3());
4862  assert(this->dim4() == rhs.dim4());
4863  std::transform(matrix_->begin(),matrix_->end(),(*(rhs.matrix_)).begin(),matrix_->begin(),std::multiplies<T>());
4864  return *this;
4865  }
4866 
4867  template<typename T>
4868  inline
4870  {
4871  assert(this->size() == rhs.size());
4872  std::transform(matrix_->begin(),matrix_->end(),(*(rhs.matrix_)).begin(),matrix_->begin(),std::divides<T>());
4873  return *this;
4874  }
4875 
4876 
4877 
4878  template<typename T>
4879  inline
4881  {
4882  assert(matrix_->size() != 0);
4883  return *std::min_element(matrix_->begin(),matrix_->end());
4884  }
4885 
4886  template<typename T>
4887  inline
4889  {
4890  assert(matrix_->size() != 0);
4891  return *std::max_element(matrix_->begin(),matrix_->end());
4892  }
4893 
4894  template<typename T>
4895  inline
4897  {
4898  assert(matrix_->size() != 0);
4899  return std::accumulate(matrix_->begin(),matrix_->end(),T());
4900  }
4901 
4902 
4903  template<typename T>
4904  inline
4906  {
4907  slip::apply(this->begin(),this->end(),this->begin(),fun);
4908  return *this;
4909  }
4910 
4911  template<typename T>
4912  inline
4914  {
4915  slip::apply(this->begin(),this->end(),this->begin(),fun);
4916  return *this;
4917  }
4918 
4919 
4920 
4921 
4922  template<typename T>
4923  inline
4925  const HyperVolume<T>& M2)
4926  {
4927  assert(M1.dim1() == M2.dim1());
4928  assert(M1.dim2() == M2.dim2());
4929  assert(M1.dim3() == M2.dim3());
4930  assert(M1.dim4() == M2.dim4());
4931 
4932  HyperVolume<T> tmp(M1.dim1(),M1.dim2(),M1.dim3(),M1.dim4());
4933  std::transform(M1.begin(),M1.end(),M2.begin(),tmp.begin(),std::plus<T>());
4934  return tmp;
4935  }
4936 
4937  template<typename T>
4938  inline
4940  const T& val)
4941  {
4942  HyperVolume<T> tmp(M1);
4943  tmp+=val;
4944  return tmp;
4945  }
4946 
4947  template<typename T>
4948  inline
4950  const HyperVolume<T>& M1)
4951  {
4952  return M1 + val;
4953  }
4954 
4955 
4956  template<typename T>
4957  inline
4959  const HyperVolume<T>& M2)
4960  {
4961  assert(M1.dim1() == M2.dim1());
4962  assert(M1.dim2() == M2.dim2());
4963  assert(M1.dim3() == M2.dim3());
4964  assert(M1.dim4() == M2.dim4());
4965 
4966  HyperVolume<T> tmp(M1.dim1(),M1.dim2(),M1.dim3(),M1.dim4());
4967  std::transform(M1.begin(),M1.end(),M2.begin(),tmp.begin(),std::minus<T>());
4968  return tmp;
4969  }
4970 
4971  template<typename T>
4972  inline
4974  const T& val)
4975  {
4976  HyperVolume<T> tmp(M1);
4977  tmp-=val;
4978  return tmp;
4979  }
4980 
4981  template<typename T>
4982  inline
4984  const HyperVolume<T>& M1)
4985  {
4986  return -(M1 - val);
4987  }
4988 
4989  template<typename T>
4990  inline
4992  const HyperVolume<T>& M2)
4993  {
4994  assert(M1.dim1() == M2.dim1());
4995  assert(M1.dim2() == M2.dim2());
4996  assert(M1.dim3() == M2.dim3());
4997  assert(M1.dim4() == M2.dim4());
4998 
4999  HyperVolume<T> tmp(M1.dim1(),M1.dim2(),M1.dim3(),M1.dim4());
5000  std::transform(M1.begin(),M1.end(),M2.begin(),tmp.begin(),std::multiplies<T>());
5001  return tmp;
5002  }
5003 
5004  template<typename T>
5005  inline
5007  const T& val)
5008  {
5009  HyperVolume<T> tmp(M1);
5010  tmp*=val;
5011  return tmp;
5012  }
5013 
5014  template<typename T>
5015  inline
5017  const HyperVolume<T>& M1)
5018  {
5019  return M1 * val;
5020  }
5021 
5022  template<typename T>
5023  inline
5025  const HyperVolume<T>& M2)
5026  {
5027  assert(M1.dim1() == M2.dim1());
5028  assert(M1.dim2() == M2.dim2());
5029  assert(M1.dim3() == M2.dim3());
5030  assert(M1.dim4() == M2.dim4());
5031 
5032  HyperVolume<T> tmp(M1.dim1(),M1.dim2(),M1.dim3(),M1.dim4());
5033  std::transform(M1.begin(),M1.end(),M2.begin(),tmp.begin(),std::divides<T>());
5034  return tmp;
5035  }
5036 
5037  template<typename T>
5038  inline
5040  const T& val)
5041  {
5042  HyperVolume<T> tmp(M1);
5043  tmp/=val;
5044  return tmp;
5045  }
5046 
5047 
5048  template<typename T>
5049  inline
5050  T& min(const HyperVolume<T>& M1)
5051  {
5052  return M1.min();
5053  }
5054 
5055  template<typename T>
5056  inline
5057  T& max(const HyperVolume<T>& M1)
5058  {
5059  return M1.max();
5060  }
5061 
5062  template<typename T>
5063  inline
5065  {
5066 
5067  HyperVolume<T> tmp(M.dim1(),M.dim2(),M.dim3(),M.dim4());
5068  slip::apply(M.begin(),M.end(),tmp.begin(),std::abs);
5069  return tmp;
5070  }
5071 
5072  template<typename T>
5073  inline
5075  {
5076  HyperVolume<T> tmp(M.dim1(),M.dim2(),M.dim3(),M.dim4());
5077  slip::apply(M.begin(),M.end(),tmp.begin(),std::sqrt);
5078  return tmp;
5079  }
5080 
5081  template<typename T>
5082  inline
5084  {
5085  HyperVolume<T> tmp(M.dim1(),M.dim2(),M.dim3(),M.dim4());
5086  slip::apply(M.begin(),M.end(),tmp.begin(),std::cos);
5087  return tmp;
5088  }
5089 
5090  template<typename T>
5091  inline
5093  {
5094  HyperVolume<T> tmp(M.dim1(),M.dim2(),M.dim3(),M.dim4());
5095  slip::apply(M.begin(),M.end(),tmp.begin(),std::acos);
5096  return tmp;
5097  }
5098 
5099  template<typename T>
5100  inline
5102  {
5103  HyperVolume<T> tmp(M.dim1(),M.dim2(),M.dim3(),M.dim4());
5104  slip::apply(M.begin(),M.end(),tmp.begin(),std::sin);
5105  return tmp;
5106  }
5107 
5108  template<typename T>
5109  inline
5111  {
5112  HyperVolume<T> tmp(M.dim1(),M.dim2(),M.dim3(),M.dim4());
5113  slip::apply(M.begin(),M.end(),tmp.begin(),std::asin);
5114  return tmp;
5115  }
5116 
5117  template<typename T>
5118  inline
5120  {
5121  HyperVolume<T> tmp(M.dim1(),M.dim2(),M.dim3(),M.dim4());
5122  slip::apply(M.begin(),M.end(),tmp.begin(),std::tan);
5123  return tmp;
5124  }
5125 
5126  template<typename T>
5127  inline
5129  {
5130  HyperVolume<T> tmp(M.dim1(),M.dim2(),M.dim3(),M.dim4());
5131  slip::apply(M.begin(),M.end(),tmp.begin(),std::atan);
5132  return tmp;
5133  }
5134 
5135  template<typename T>
5136  inline
5138  {
5139  HyperVolume<T> tmp(M.dim1(),M.dim2(),M.dim3(),M.dim4());
5140  slip::apply(M.begin(),M.end(),tmp.begin(),std::exp);
5141  return tmp;
5142  }
5143 
5144  template<typename T>
5145  inline
5147  {
5148  HyperVolume<T> tmp(M.dim1(),M.dim2(),M.dim3(),M.dim4());
5149  slip::apply(M.begin(),M.end(),tmp.begin(),std::log);
5150  return tmp;
5151  }
5152 
5153  template<typename T>
5154  inline
5156  {
5157  HyperVolume<T> tmp(M.dim1(),M.dim2(),M.dim3(),M.dim4());
5158  slip::apply(M.begin(),M.end(),tmp.begin(),std::cosh);
5159  return tmp;
5160  }
5161 
5162  template<typename T>
5163  inline
5165  {
5166  HyperVolume<T> tmp(M.dim1(),M.dim2(),M.dim3(),M.dim4());
5167  slip::apply(M.begin(),M.end(),tmp.begin(),std::sinh);
5168  return tmp;
5169  }
5170 
5171  template<typename T>
5172  inline
5174  {
5175  HyperVolume<T> tmp(M.dim1(),M.dim2(),M.dim3(),M.dim4());
5176  slip::apply(M.begin(),M.end(),tmp.begin(),std::tanh);
5177  return tmp;
5178  }
5179 
5180  template<typename T>
5181  inline
5183  {
5184  HyperVolume<T> tmp(M.dim1(),M.dim2(),M.dim3(),M.dim4());
5185  slip::apply(M.begin(),M.end(),tmp.begin(),std::log10);
5186  return tmp;
5187  }
5188 
5189 
5190 
5191 }//slip::
5192 
5193 #endif //SLIP_HYPER_VOLUME_HPP
slip::stride_iterator< pointer > slice_iterator
std::reverse_iterator< col_range_iterator > reverse_col_range_iterator
T & max() const
Returns the max element of the HyperVolume according to the operator <.
HyperVolume< T > tanh(const HyperVolume< T > &M)
bool operator!=(const Array< T > &x, const Array< T > &y)
Definition: Array.hpp:1448
This is some iterator to iterate a 4d container into a Box area defined by the subscripts of the 4d c...
iterator4d first_front_upper_left()
Returns a read/write iterator4d that points to the first element of the Matrix4d. It points to the fi...
Definition: Matrix4d.hpp:4372
slip::stride_iterator< const_pointer > const_row_range_iterator
iterator4d last_back_bottom_right()
Returns a read/write iterator4d that points to the past the end element of the HyperVolume. It points to past the end element of the last back bottom right element of the HyperVolume.
reverse_iterator4d rfirst_front_upper_left()
Returns a read/write reverse iterator4d. It points to the last back bottom right element of the Matri...
Definition: Matrix4d.hpp:4421
T *** operator[](const size_type l)
reverse_iterator rend()
Returns a read/write reverse iterator that points to one before the first element in the HyperVolume...
slice_iterator slice_end(const size_type slab, const size_type row, const size_type col)
Returns a read/write iterator that points to the one past the end element of the line (slab...
iterator4d last_back_bottom_right()
Returns a read/write iterator4d that points to the past the end element of the Matrix4d. It points to past the end element of the last back bottom right element of the Matrix4d.
Definition: Matrix4d.hpp:4388
reverse_row_iterator row_rbegin(const size_type slab, const size_type slice, const size_type row)
Returns a read/write reverse iterator that points to the last element of the row row of a given slab ...
slip::Range< int > range
HyperVolume(const size_type d1, const size_type d2, const size_type d3, const size_type d4, InputIterator first, InputIterator last)
Contructs a HyperVolume from a range.
col_iterator col_end(const size_type slab, const size_type slice, const size_type col)
Returns a read/write iterator that points to the past-the-end element of the column column of a given...
std::reverse_iterator< const_row_range_iterator > const_reverse_row_range_iterator
const value_type * const_pointer
std::reverse_iterator< slab_range_iterator > reverse_slab_range_iterator
size_type slices() const
Returns the number of slices (second dimension size) in the HyperVolume.
size_type dim2() const
Returns the number of slices (second dimension size) in the HyperVolume.
size_type max_size() const
Returns the maximal size (number of elements) in the HyperVolume.
col_iterator col_begin(const size_type slab, const size_type slice, const size_type col)
Returns a read/write iterator that points to the first element of the column column of a given slab s...
slip::HyperVolume< unsigned long > HyperVolume_ul
unsigned long alias
iterator4d default_iterator
HyperVolume< T > & apply(T(*fun)(T))
Applys the one-parameter C-function fun to each element of the HyperVolume.
This is some iterator to iterate a 4d container into two Range defined by the indices and strides of ...
reverse_iterator4d rfirst_front_upper_left()
Returns a read/write reverse iterator4d. It points to the last back bottom right element of the Hyper...
slip::stride_iterator< const_slab_iterator > const_slab_range_iterator
std::reverse_iterator< slice_iterator > reverse_slice_iterator
HyperVolume< T > cosh(const HyperVolume< T > &M)
Numerical matrix4d class. This container statisfies the RandomAccessContainer concepts of the STL exc...
slip::HyperVolume< float > HyperVolume_f
float alias
std::reverse_iterator< iterator4d_range > reverse_iterator4d_range
T & min(const GrayscaleImage< T > &M1)
Returns the min element of a GrayscaleImage.
row_iterator row_end(const size_type slab, const size_type slice, const size_type row)
Returns a read/write iterator that points to the past-the-end element of the row row of a given slab ...
slip::Array4d< T >::const_iterator4d const_iterator4d
This is some iterator to iterate a 4d container into a Box area defined by the subscripts of the 4d c...
Color< T > operator+(const Color< T > &V1, const Color< T > &V2)
Definition: Color.hpp:602
const_iterator4d const_default_iterator
std::reverse_iterator< row_range_iterator > reverse_row_range_iterator
const value_type & const_reference
T sum() const
Returns the sum of the elements of the HyperVolume.
size_type cols() const
Returns the number of columns (fourth dimension size) in the HyperVolume.
slice_iterator slice_begin(const size_type slab, const size_type row, const size_type col)
Returns a read/write iterator that points to the first element of the line (slab,row,col) through the slices in the HyperVolume. Iteration is done in ordinary element order (increasing slice number).
HyperVolume< T > atan(const HyperVolume< T > &M)
std::reverse_iterator< const_slab_range_iterator > const_reverse_slab_range_iterator
slip::stride_iterator< const_pointer > const_slice_iterator
slip::stride_iterator< const_pointer > const_slab_iterator
size_type rows() const
Returns the number of rows (third dimension size) in the HyperVolume.
friend class boost::serialization::access
value_type & reference
self & operator*=(const T &val)
std::reverse_iterator< const_iterator4d_range > const_reverse_iterator4d_range
HyperVolume< T > exp(const HyperVolume< T > &M)
std::reverse_iterator< col_iterator > reverse_col_iterator
size_type dim4() const
Returns the number of columns (fourth dimension size) in the HyperVolume.
Provides a class to manipulate Matrix4d.
slip::stride_iterator< pointer > col_iterator
std::reverse_iterator< const_iterator4d > const_reverse_iterator4d
HyperVolume< T > abs(const HyperVolume< T > &M)
bool operator>(const Array< T > &x, const Array< T > &y)
Definition: Array.hpp:1469
slip::stride_iterator< slice_iterator > slice_range_iterator
reverse_slice_iterator slice_rbegin(const size_type slab, const size_type row, const size_type col)
Returns a read/write iterator that points to the last element of the line (slab,row,col) through the slices in the HyperVolume. Iteration is done in reverse element order (decreasing slice number).
size_type slabs() const
Returns the number of slabs (first dimension size) in the HyperVolume.
HyperVolume< T > sin(const HyperVolume< T > &M)
std::reverse_iterator< iterator > reverse_row_iterator
Numerical container class This is a four-dimensional dynamic and generic container. This container satisfies the BidirectionnalContainer concepts of the STL. It is also an 4d extension of the RandomAccessContainer concept. That is to say the bracket element access is replaced by the triple bracket element access. Data are stored using a Matrix4d class. It extends the interface of Matrix4d adding image read/write operations. These operations are done using the ImageMagick library.
self & operator/=(const T &val)
HyperVolume< T > cos(const HyperVolume< T > &M)
void fill(const T *value)
Fills the container range [begin(),begin()+size()) with a copy of the value array.
reverse_slab_iterator slab_rbegin(const size_type slice, const size_type row, const size_type col)
Returns a read/write iterator that points to the last element of the line (slice,row,col) through the slabs in the HyperVolume. Iteration is done in reverse element order (decreasing slab number).
const_pointer const_iterator
void swap(HyperVolume &M)
Swaps data with another HyperVolume.
slip::stride_iterator< const_slice_iterator > const_slice_range_iterator
HyperVolume()
Constructs a HyperVolume.
size_type columns() const
Returns the number of columns (fourth dimension size) in the HyperVolume.
std::reverse_iterator< iterator > reverse_iterator
slip::HyperVolume< short > HyperVolume_s
short alias
bool empty() const
Returns true if the HyperVolume is empty. (Thus size() == 0)
slip::stride_iterator< slab_iterator > slab_range_iterator
row_iterator row_begin(const size_type slab, const size_type slice, const size_type row)
Returns a read/write iterator that points to the first element of the row row of a given slab slab an...
slip::Array4d< T >::iterator4d iterator4d
const_pointer const_row_iterator
std::reverse_iterator< const_col_range_iterator > const_reverse_col_range_iterator
self & operator-=(const T &val)
iterator4d first_front_upper_left()
Returns a read/write iterator4d that points to the first element of the HyperVolume. It points to the first front upper left element of the HyperVolume.
const_iterator begin() const
Returns a read-only (constant) iterator that points to the first element in the Matrix4d. Iteration is done in ordinary element order.
Definition: Matrix4d.hpp:3668
~HyperVolume()
Destructor of the HyperVolume.
void copy(_II first, _II last, _OI output_first)
Copy algorithm optimized for slip iterators.
Definition: copy_ext.hpp:177
reverse_slice_iterator slice_rend(const size_type slab, const size_type row, const size_type col)
Returns a read/write iterator that points to the one before the first element of the line (slab...
slip::Array4d< T >::iterator4d_range iterator4d_range
size_type dim1() const
Returns the number of slabs (first dimension size) in the HyperVolume.
void fill(const T &value)
Fills the container range [begin(),begin()+size()) with copies of value.
slab_iterator slab_end(const size_type slice, const size_type row, const size_type col)
Returns a read/write iterator that points to the one past the end element of the line (slice...
reverse_row_iterator row_rend(const size_type slab, const size_type slice, const size_type row)
Returns a read/write reverse iterator that points to the first element of the row row of a given slab...
iterator end()
Returns a read/write iterator that points one past the last element in the Matrix4d. Iteration is done in ordinary element order.
Definition: Matrix4d.hpp:3691
reverse_col_iterator col_rend(const size_type slab, const size_type slice, const size_type col)
Returns a read/write reverse iterator that points to the first element of the column column of a give...
reference operator()(const size_type l, const size_type k, const size_type i, const size_type j)
Subscript access to the data contained in the HyperVolume.
HyperVolume< T > sqrt(const HyperVolume< T > &M)
slip::stride_iterator< const_pointer > const_col_iterator
HyperVolume< T > sinh(const HyperVolume< T > &M)
slip::HyperVolume< long > HyperVolume_l
long alias
slip::HyperVolume< unsigned char > HyperVolume_uc
unsigned char alias
std::size_t size_type
ptrdiff_t difference_type
void resize(std::size_t d1, std::size_t d2, std::size_t d3, std::size_t d4, const T &val=T())
Resizes a HyperVolume.
reverse_slab_iterator slab_rend(const size_type slice, const size_type row, const size_type col)
Returns a read/write iterator that points to the one before the first element of the line (slice...
HyperVolume< T > log(const HyperVolume< T > &M)
std::reverse_iterator< const_iterator > const_reverse_row_iterator
std::string name() const
Returns the name of the class.
Some const iterator to iterate a 4d container into two Range defined by the indices and strides of th...
slip::stride_iterator< col_iterator > col_range_iterator
HyperVolume< T > acos(const HyperVolume< T > &M)
std::reverse_iterator< slice_range_iterator > reverse_slice_range_iterator
std::reverse_iterator< slab_iterator > reverse_slab_iterator
slip::stride_iterator< pointer > slab_iterator
std::reverse_iterator< const_slab_iterator > const_reverse_slab_iterator
slab_iterator slab_begin(const size_type slice, const size_type row, const size_type col)
Returns a read/write iterator that points to the first element of the line (slice,row,col) through the slabs in the HyperVolume. Iteration is done in ordinary element order (increasing slab number).
HyperVolume< T > log10(const HyperVolume< T > &M)
slip::HyperVolume< unsigned short > HyperVolume_us
unsigned long alias
reverse_iterator4d rlast_back_bottom_right()
Returns a read/write reverse iterator4d. It points to past the first front upper left element of the ...
value_type * pointer
const_iterator begin() const
Returns a read-only (constant) iterator that points to the first element in the HyperVolume. Iteration is done in ordinary element order.
slip::HyperVolume< char > HyperVolume_c
char alias
T & min() const
Returns the min element of the HyperVolume according to the operator <.
std::reverse_iterator< const_slice_iterator > const_reverse_slice_iterator
std::reverse_iterator< const_col_iterator > const_reverse_col_iterator
slip::Array4d< T >::const_iterator4d_range const_iterator4d_range
slip::stride_iterator< pointer > row_range_iterator
Provides some algorithms to apply C like functions to ranges.
Provides a class to manipulate iterator4d within a slip::Range. It is used to iterate throw 4d contai...
Provides a class to iterate a 1d range according to a step.
self & operator=(const HyperVolume< T > &rhs)
Assign a HyperVolume.
HyperVolume< T > tan(const HyperVolume< T > &M)
size_type dim3() const
Returns the number of rows (third dimension size) in the HyperVolume.
self operator-() const
reverse_col_iterator col_rbegin(const size_type slab, const size_type slice, const size_type col)
Returns a read/write reverse iterator that points to the last element of the column column of a given...
bool operator==(const Array< T > &x, const Array< T > &y)
Definition: Array.hpp:1439
HyperVolume< T > asin(const HyperVolume< T > &M)
slip::HyperVolume< unsigned int > HyperVolume_ui
unsigned int alias
Provides a class to manipulate iterator4d within a slip::Box4d. It is used to iterate throw 4d contai...
Color< T > operator*(const Color< T > &V1, const Color< T > &V2)
Definition: Color.hpp:658
self & operator+=(const T &val)
Add val to each element of the HyperVolume.
slip::HyperVolume< double > HyperVolume_d
double alias
reverse_iterator4d rlast_back_bottom_right()
Returns a read/write reverse iterator4d. It points to past the first front upper left element of the ...
Definition: Matrix4d.hpp:4404
std::reverse_iterator< const_iterator > const_reverse_iterator
void apply(InputIterator first, InputIterator last, OutputIterator result, typename std::iterator_traits< OutputIterator >::value_type(*function)(typename std::iterator_traits< InputIterator >::value_type))
Applies a C-function to each element of a range.
Definition: apply.hpp:128
T & max(const GrayscaleImage< T > &M1)
Returns the max element of a GrayscaleImage.
Color< T > operator/(const Color< T > &V1, const Color< T > &V2)
Definition: Color.hpp:686
iterator end()
Returns a read/write iterator that points one past the last element in the HyperVolume. Iteration is done in ordinary element order.
std::reverse_iterator< const_slice_range_iterator > const_reverse_slice_range_iterator
std::reverse_iterator< iterator4d > reverse_iterator4d
bool operator>=(const Array< T > &x, const Array< T > &y)
Definition: Array.hpp:1485
slip::stride_iterator< const_col_iterator > const_col_range_iterator
Color< T > operator-(const Color< T > &V1, const Color< T > &V2)
Definition: Color.hpp:630
static const std::size_t DIM
void fill(InputIterator first, InputIterator last)
Fills the container range [begin(),begin()+size()) with a copy of the range [first,last)
reverse_iterator rbegin()
Returns a read/write reverse iterator that points to the last element in the HyperVolume. Iteration is done in reverse element order.
size_type size() const
Returns the number of elements in the HyperVolume.
slip::HyperVolume< int > HyperVolume_i
int alias