File: view.rst

package info (click to toggle)
xtensor 0.25.0-2
  • links: PTS, VCS
  • area: main
  • in suites: forky, sid, trixie
  • size: 6,476 kB
  • sloc: cpp: 65,302; makefile: 202; python: 171; javascript: 8
file content (421 lines) | stat: -rw-r--r-- 14,644 bytes parent folder | download
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
.. Copyright (c) 2016, Johan Mabille, Sylvain Corlay and Wolf Vollprecht

   Distributed under the terms of the BSD 3-Clause License.

   The full license is in the file LICENSE, distributed with this software.

.. _view-description:

Views
=====

Views are used to adapt the shape of an :cpp:type:`xt::xexpression` without changing it, nor copying it. Views are
convenient tools for assigning parts of an expression: since they do not copy the underlying expression,
assigning to the view actually assigns to the underlying expression. *xtensor* provides many kinds of views.

Sliced views
------------

Sliced views consist of the combination of the :cpp:type:`xt::xexpression` to adapt, and a list of ``slice`` that specify how
the shape must be adapted. Sliced views are implemented by the :cpp:type:`xt::xview` class. Objects of this type should not be
instantiated directly, but though the :cpp:func:`xt::view` helper function.

Slices can be specified in the following ways:

- selection in a dimension by specifying an index (unsigned integer)
- :cpp:func:`xt::range(min, max) <xt::range>`, a slice representing the interval [min, max)
- :cpp:func:`xt::range(min, max, step) <xt::range>`, a slice representing the stepped interval [min, max)
- :cpp:func:`xt::all`, a slice representing all the elements of a dimension
- :cpp:func:`xt::newaxis`, a slice representing an additional dimension of length one
- :cpp:func:`xt::keep(i0, i1, i2, ...) <xt::keep>` a slice selecting non-contiguous indices to keep on the underlying expression
- :cpp:func:`xt::drop(i0, i1, i2, ...) <xt::drop>` a slice selecting non-contiguous indices to drop on the underlying expression

.. code::

    #include <vector>
    #include <xtensor/xarray.hpp>
    #include <xtensor/xview.hpp>

    std::vector<size_t> shape = {3, 2, 4};
    xt::xarray<int> a(shape);

    // View with same number of dimensions
    auto v1 = xt::view(a, xt::range(1, 3), xt::all(), xt::range(1, 3));
    // => v1.shape() = { 2, 2, 2 }
    // => v1(0, 0, 0) = a(1, 0, 1)
    // => v1(1, 1, 1) = a(2, 1, 2)

    // View reducing the number of dimensions
    auto v2 = xt::view(a, 1, xt::all(), xt::range(0, 4, 2));
    // => v2.shape() = { 2, 2 }
    // => v2(0, 0) = a(1, 0, 0)
    // => v2(1, 1) = a(1, 1, 2)

    // View increasing the number of dimensions
    auto v3 = xt::view(a, xt::all(), xt::all(), xt::newaxis(), xt::all());
    // => v3.shape() = { 3, 2, 1, 4 }
    // => v3(0, 0, 0, 0) = a(0, 0, 0)

    // View with non contiguous slices
    auto v4 = xt::view(a, xt::drop(0), xt::all(), xt::keep(0, 3));
    // => v4.shape() = { 2, 2, 2 }
    // => v4(0, 0, 0) = a(1, 0, 0)
    // => v4(1, 1, 1) = a(2, 1, 3)

    // View built with negative index
    auto v5 = xt::view(a, -2, xt::all(), xt::range(0, 4, 2));
    // => v5 == v2

The range function supports the placeholder ``_`` syntax:

.. code::

    #include <xtensor/xarray.hpp>
    #include <xtensor/xview.hpp>

    using namespace xt::placeholders;  // required for ``_`` to work

    auto a = xt::xarray<int>::from_shape({3, 2, 4});
    auto v1 = xt::view(a, xt::range(_, 2), xt::all(), xt::range(1, _));
    // The previous line is equivalent to
    auto v2 = xt::view(a, xt::range(0, 2), xt::all(), xt::range(1, 4));

:cpp:type:`xt::xview` does not perform a copy of the underlying expression.
This means if you modify an element of the :cpp:type:`xt::xview`,
you are actually also altering the underlying expression.

.. code::

    #include <vector>
    #include <xtensor/xarray.hpp>
    #include <xtensor/xview.hpp>

    std::vector<size_t> shape = {3, 2, 4};
    xt::xarray<int> a(shape, 0);

    auto v1 = xt::view(a, 1, xt::all(), xt::range(1, 3));
    v1(0, 0) = 1;
    // => a(1, 0, 1) = 1

The convenient methods :cpp:func:`xt::row` and :cpp:func:`xt::col` are available for 2-D expressions:

.. code::

    #include <vector>
    #include <xtensor/xtensor.hpp>
    #include <xtensor/xview.hpp>

    xt::xtensor<double, 2> a = {{1, 2}, {3, 4}};
    auto r = xt::row(a, 0);
    // => r = {1, 2}
    auto c = xt::col(a, -1);
    // => c = { 2, 4 }

Strided views
-------------

While the :cpp:func:`xt::view` is a compile-time static expression, xtensor also contains a dynamic
strided view in ``xstrided_view.hpp``.
The strided view and the slice vector allow to dynamically push_back slices, so when the dimension
is unknown at compile time, the slice vector can be built dynamically at runtime.
Note that the slice vector is actually a type-alias for a ``std::vector`` of a ``variant`` for
all the slice types.
The strided view does not support the slices returned by the :cpp:func:`xt::keep` and
:cpp:func:`xt::drop` functions.

.. code::

    #include <xtensor/xarray.hpp>
    #include <xtensor/xstrided_view.hpp>

    auto a = xt::xarray<int>::from_shape({3, 2, 3, 4, 5});

    xt::xstrided_slice_vector sv({xt::range(0, 1), xt::newaxis()});
    sv.push_back(1);
    sv.push_back(xt::all());

    auto v1 = xt::strided_view(a, sv);
    // v1 has the same behavior as the static view

    // Equivalent but shorter
    auto v2 = xt::strided_view(a, { xt::range(0, 1), xt::newaxis(), 1, xt::all() });
    // v2 == v1

    // ILLEGAL:
    auto v2 = xt::strided_view(a, { xt::all(), xt::all(), xt::all(), xt::keep(0, 3), xt::drop(1, 4) });
    // xt::drop and xt::keep are not supported with strided views

Since ``xtensor 0.16.3``, a new range syntax can be used with strided views:

.. code::

    #include <xtensor/xarray.hpp>
    #include <xtensor/xstrided_view.hpp>

    using namespace xt::placeholders;

    auto a = xt::xarray<int>::from_shape({3, 2, 3, 4, 5});
    auto v1 = xt::strided_view(a, {_r|0|1, 1, _r|_|2, _r|_|_|-1});
    // The previous line is equivalent to
    auto v2 = xt::strided_view(a, {xt::range(0, 1), 1, xt::range(_, 2), xt::range(_, _, -1)});

The :cpp:type:`xt::xstrided_view` is very efficient on contigous memory
(e.g. :cpp:type:`xt::xtensor` or :cpp:type:`xt::xarray`) but less efficient on\
:cpp:type:`xt::xexpression`s.

Transposed views
----------------

*xtensor* provides a lazy transposed view on any expression, whose layout is either row-major order or column major order.
Trying to build a transposed view on a expression with a dynamic layout throws an exception.

.. code::

    #include <xtensor/xarray.hpp>
    #include <xtensor/xstrided_view.hpp>

    xt::xarray<int> a = { {0, 1, 2}, {3, 4, 5} };
    auto tr = xt::transpose(a);
    // tr == { {0, 3}, {1, 4}, {2, 5} }

    xt::xarray<int, layout_type::dynamic> b = { {0, 1, 2}, {3, 4, 5} };
    auto tr2 = xt::transpose(b);
    // => throw transpose_error

Like the strided view, the transposed view is built upon the :cpp:type:`xt::xstrided_view`.

Flatten views
-------------

It is sometimes useful to have a one-dimensional view of all the elements of an expression.
*xtensor* provides two functions for that, :cpp:func:`xt::ravel` and :cpp:func:`xt::flatten`.
The former one lets you specify the order used to read the elements while the latter one
uses the layout of the expression.

.. code::

    #include <xtensor/xarray.hpp>
    #include <xtensor/xstrided_view.hpp>

    xt::xarray<int> a = { {0, 1, 2}, {3, 4, 5} };
    auto flc = xt::ravel<layout_type::column_major>(a);
    std::cout << flc << std::endl;
    // => prints { 0, 3, 1, 4, 2, 5 }

    auto fl = xt::flatten(a);
    std::cout << fl << std::endl;
    // => prints { 0, 1, 2, 3, 4, 5 }

Like the strided view and the transposed view, the flatten view is built upon the :cpp:type:`xt::xstrided_view`.

Reshape views
-------------

The reshape view allows to handle an expression as if it was given a new shape, however no additional memory allocation occurs,
the original expression keeps its shape. Like any view, the underlying expression is not copied, thus assigning a value through
the view modifies the underlying expression.

.. code::

    #include <xtensor/xarray.hpp>
    #include <xtensor/xstrided_view.hpp>

    auto a = xt::xarray<int>::from_shape({3, 2, 4});
    auto v = xt::reshape_view(a, { 4, 2, 3 });
    // a(0, 0, 3) == v(0, 1, 0)
    // a(0, 1, 0) == v(0, 1, 1)

    v(0, 2, 0) = 4;
    // a(0, 1, 2) == 4

Like the strided view and the transposed view, the reshape view is built upon the :cpp:type:`xt::xstrided_view`.

Dynamic views
-------------

The dynamic view is like the strided view, but with support of the slices returned by the
:cpp:func:`xt::keep` and :cpp:func:`xt::drop` functions.
However, this support has a cost and the dynamic view is slower than the strided view, even when no
keeping or dropping of a slice is involved.

.. code::

    #include <xtensor/xarray.hpp>
    #include <xtensor/xdynamic_view.hp>

    auto a = xt::xarray<int>::from_shape({3, 2, 3, 4, 5});
    xt::xdynamic_slice_vector sv({xt::range(0, 1), xt::newaxis()});
    sv.push_back(1);
    sv.push_back(xt::all());
    sv.push_back(xt::keep(0, 2, 3));
    sv.push_back(xt::drop(1, 2, 4));

    auto v1 = xt::dynamic_view(a, sv});

    // Equivalent but shorter
    auto v2 = xt::dynamic_view(a, { xt::range(0, 1), xt::newaxis(), 1, xt::all(), xt::keep(0, 2, 3), xt::drop(1, 2, 4) });
    // v2 == v1

Index views
-----------

Index views are one-dimensional views of an :cpp:type:`xt::xexpression`, containing the elements
whose positions are specified by a list of indices.
Like for sliced views, the elements of the underlying :cpp:type:`xt::xexpression` are not copied.
Index views should be built with the :cpp:func:`xt::index_view` helper function.

.. code::

    #include <xtensor/xarray.hpp>
    #include <xtensor/xindex_view.hpp>

    xt::xarray<double> a = {{1, 5, 3}, {4, 5, 6}};
    auto b = xt::index_view(a, {{0,0}, {1, 0}, {0, 1}});
    // => b = { 1, 4, 5 }
    b += 100;
    // => a = {{101, 5, 3}, {104, 105, 6}}

The type used for representing indices can be any 1-D container providing an ``std::vector``-like API.
The same stands for the type of the list of indices:

.. code::

    #include <xtensor/xarray.hpp>
    #include <xtensor/xindex_view.hpp>

    xt::xarray<double> a = {{1, 5, 3}, {4, 5, 6}};
    using index_type = std::array<std::size_t, 2>;
    std::vector<index_type> indices = {{0, 0}, {1, 0}, {0, 1}};
    auto b = xt::index_view(a, indices);
    // => b = { 1, 4, 5 }
    b += 100;
    // => a = {{101, 5, 3}, {104, 105, 6}}

Filter views
------------

Filters are one-dimensional views holding elements of an :cpp:type:`xt::xexpression` that verify a given condition.
Like for other views, the elements of the underlying :cpp:type:`xt::xexpression` are not copied.
Filters should be built with the :cpp:func:`xt::filter` helper function.

.. code::

    #include <xtensor/xarray.hpp>
    #include <xtensor/xindex_view.hpp>

    xt::xarray<double> a = {{1, 5, 3}, {4, 5, 6}};
    auto v = xt::filter(a, a >= 5);
    // => v = { 5, 5, 6 }
    v += 100;
    // => a = {{1, 105, 3}, {4, 105, 106}}

Filtration
----------

Sometimes, the only thing you want to do with a filter is to assign it a scalar.
Though this can be done as shown in the previous section, this is not the *optimal* way to do it.
*xtensor* provides a specially optimized mechanism for that, called filtration.
A filtration IS NOT an :cpp:type:`xt::xexpression`, the only methods it provides are scalar and
computed scalar assignments.

.. code::

    #include <xtensor/xarray.hpp>
    #include <xtensor/xindex_view.hpp>

    xt::xarray<double> a = {{1, 5, 3}, {4, 5, 6}};
    filtration(a, a >= 5) += 100;
    // => a = {{1, 105, 3}, {4, 105, 106}}

Masked view
-----------

Masked views are multidimensional views that apply a mask on an :cpp:type:`xt::xexpression`.

.. code::

    #include <xtensor/xarray.hpp>
    #include <xtensor/xmasked_view.hpp>

    xt::xarray<double> a = {{1, 5, 3}, {4, 5, 6}};
    xt::xarray<bool> mask = {{true, false, false}, {false, true, false}};

    auto m = xt::masked_view(a, mask);
    // => m = {{1, masked, masked}, {masked, 5, masked}}

    m += 100;
    // => a = {{101, 5, 3}, {4, 105, 6}}

Broadcasting views
------------------

Another type of view provided by *xtensor* is *broadcasting view*.
Such a view broadcasts an expression to the specified shape.
As long as the view is not assigned to an array, no memory allocation or copy occurs.
Broadcasting views should be built with the :cpp:func:`xt::broadcast` helper function.

.. code::

    #include <vector>
    #include <xtensor/xarray.hpp>
    #include <xtensor/xbroadcast.hpp>

    std::vector<size_t> s1 = { 2, 3 };
    std::vector<size_t> s2 = { 3, 2, 3 };

    xt::xarray<int> a1(s1);
    auto bv = xt::broadcast(a1, s2);
    // => bv(0, 0, 0) = bv(1, 0, 0) = bv(2, 0, 0) = a(0, 0)

Complex views
-------------

In the case of a tensor containing complex numbers, *xtensor* provides views returning
:cpp:type:`xt::xexpression` corresponding to the real and imaginary parts of the complex numbers.
Like for other views, the elements of the underlying :cpp:type:`xt::xexpression` are not copied.

Functions :cpp:func:`xt::real` and :cpp:func:`xt::imag` respectively return views on the real and
imaginary part of a complex expression.
The returned value is an expression holding a closure on the passed argument.

- The constness and value category (rvalue / lvalue) of :cpp:func:`xt::real(a) <xt::real>` is the same
  as that of ``a``.
  Hence, if ``a`` is a non-const lvalue, :cpp:func:`xt::real(a) <xt::real>` is an non-const lvalue
  reference, to which one can assign a real expression.
- If ``a`` has complex values, the same holds for :cpp:func:`xt::imag(a) <xt::imag>`.
  The constness and value category of :cpp:func:`xt::imag(a) <xt::imag>` is the same as that of ``a``.
- If ``a`` has real values, :cpp:func:`xt::imag(a) <xt::imag>` returns
  :cpp:func:`xt::zeros(a.shape()) <xt::zeros>`.

.. code::

    #include <complex>
    #include <xtensor/xarray.hpp>
    #include <xtensor/xcomplex.hpp>

    using namespace std::complex_literals;

    xarray<std::complex<double>> e =
        {{1.0       , 1.0 + 1.0i},
         {1.0 - 1.0i, 1.0       }};

    real(e) = zeros<double>({2, 2});
    // => e = {{0.0, 0.0 + 1.0i}, {0.0 - 1.0i, 0.0}};

Assigning to a view
-------------------

When assigning an expression ``rhs`` to a container such as :cpp:type:`xt::xarray`, the container
is resized so its shape is the same as the one of ``rhs``.
However, since views *cannot be resized*, when assigning an expression to a view, broadcasting rules are applied:

.. code::

    #include <xtensor/xarray.hpp>
    #include <xtensor/xview.hpp>

    xarray<double> a = {{0., 1., 2.}, {3., 4., 5.}};
    double b = 1.2;
    auto tr = view(a, 0, all());
    tr = b;
    // => a = {{1.2, 1.2, 1.2}, {3., 4., 5.}}