Coverage Report

Created: 2026-03-14 06:47

next uncovered line (L), next uncovered region (R), next uncovered branch (B)
/rust/registry/src/index.crates.io-1949cf8c6b5b557f/nalgebra-0.34.1/src/base/edition.rs
Line
Count
Source
1
use num::{One, Zero};
2
use std::cmp;
3
#[cfg(any(feature = "std", feature = "alloc"))]
4
use std::iter::ExactSizeIterator;
5
use std::ptr;
6
7
use crate::base::allocator::{Allocator, Reallocator};
8
use crate::base::constraint::{DimEq, SameNumberOfColumns, SameNumberOfRows, ShapeConstraint};
9
#[cfg(any(feature = "std", feature = "alloc"))]
10
use crate::base::dimension::Dyn;
11
use crate::base::dimension::{Const, Dim, DimAdd, DimDiff, DimMin, DimMinimum, DimSub, DimSum, U1};
12
use crate::base::storage::{RawStorage, RawStorageMut, ReshapableStorage};
13
use crate::base::{DefaultAllocator, Matrix, OMatrix, RowVector, Scalar, Vector};
14
use crate::{Storage, UninitMatrix};
15
use std::mem::MaybeUninit;
16
17
/// # Triangular matrix extraction
18
impl<T: Scalar + Zero, R: Dim, C: Dim, S: Storage<T, R, C>> Matrix<T, R, C, S> {
19
    /// Extracts the upper triangular part of this matrix (including the diagonal).
20
    #[inline]
21
    #[must_use]
22
    pub fn upper_triangle(&self) -> OMatrix<T, R, C>
23
    where
24
        DefaultAllocator: Allocator<R, C>,
25
    {
26
        let mut res = self.clone_owned();
27
        res.fill_lower_triangle(T::zero(), 1);
28
29
        res
30
    }
31
32
    /// Extracts the lower triangular part of this matrix (including the diagonal).
33
    #[inline]
34
    #[must_use]
35
    pub fn lower_triangle(&self) -> OMatrix<T, R, C>
36
    where
37
        DefaultAllocator: Allocator<R, C>,
38
    {
39
        let mut res = self.clone_owned();
40
        res.fill_upper_triangle(T::zero(), 1);
41
42
        res
43
    }
44
}
45
46
/// # Rows and columns extraction
47
impl<T: Scalar, R: Dim, C: Dim, S: Storage<T, R, C>> Matrix<T, R, C, S> {
48
    /// Creates a new matrix by extracting the given set of rows from `self`.
49
    #[cfg(any(feature = "std", feature = "alloc"))]
50
    #[must_use]
51
    pub fn select_rows<'a, I>(&self, irows: I) -> OMatrix<T, Dyn, C>
52
    where
53
        I: IntoIterator<Item = &'a usize>,
54
        I::IntoIter: ExactSizeIterator + Clone,
55
        DefaultAllocator: Allocator<Dyn, C>,
56
    {
57
        let irows = irows.into_iter();
58
        let ncols = self.shape_generic().1;
59
        let mut res = Matrix::uninit(Dyn(irows.len()), ncols);
60
61
        // First, check that all the indices from irows are valid.
62
        // This will allow us to use unchecked access in the inner loop.
63
        for i in irows.clone() {
64
            assert!(*i < self.nrows(), "Row index out of bounds.")
65
        }
66
67
        for j in 0..ncols.value() {
68
            // TODO: use unchecked column indexing
69
            let mut res = res.column_mut(j);
70
            let src = self.column(j);
71
72
            for (destination, source) in irows.clone().enumerate() {
73
                // Safety: all indices are in range.
74
                unsafe {
75
                    *res.vget_unchecked_mut(destination) =
76
                        MaybeUninit::new(src.vget_unchecked(*source).clone());
77
                }
78
            }
79
        }
80
81
        // Safety: res is now fully initialized.
82
        unsafe { res.assume_init() }
83
    }
84
85
    /// Creates a new matrix by extracting the given set of columns from `self`.
86
    #[cfg(any(feature = "std", feature = "alloc"))]
87
    #[must_use]
88
    pub fn select_columns<'a, I>(&self, icols: I) -> OMatrix<T, R, Dyn>
89
    where
90
        I: IntoIterator<Item = &'a usize>,
91
        I::IntoIter: ExactSizeIterator,
92
        DefaultAllocator: Allocator<R, Dyn>,
93
    {
94
        let icols = icols.into_iter();
95
        let nrows = self.shape_generic().0;
96
        let mut res = Matrix::uninit(nrows, Dyn(icols.len()));
97
98
        for (destination, source) in icols.enumerate() {
99
            // NOTE: this is basically a copy_frow but wrapping the values insnide of MaybeUninit.
100
            res.column_mut(destination)
101
                .zip_apply(&self.column(*source), |out, e| *out = MaybeUninit::new(e));
102
        }
103
104
        // Safety: res is now fully initialized.
105
        unsafe { res.assume_init() }
106
    }
107
}
108
109
/// # Set rows, columns, and diagonal
110
impl<T: Scalar, R: Dim, C: Dim, S: RawStorageMut<T, R, C>> Matrix<T, R, C, S> {
111
    /// Fills the diagonal of this matrix with the content of the given vector.
112
    #[inline]
113
    pub fn set_diagonal<R2: Dim, S2>(&mut self, diag: &Vector<T, R2, S2>)
114
    where
115
        R: DimMin<C>,
116
        S2: RawStorage<T, R2>,
117
        ShapeConstraint: DimEq<DimMinimum<R, C>, R2>,
118
    {
119
        let (nrows, ncols) = self.shape();
120
        let min_nrows_ncols = cmp::min(nrows, ncols);
121
        assert_eq!(diag.len(), min_nrows_ncols, "Mismatched dimensions.");
122
123
        for i in 0..min_nrows_ncols {
124
            unsafe { *self.get_unchecked_mut((i, i)) = diag.vget_unchecked(i).clone() }
125
        }
126
    }
127
128
    /// Fills the diagonal of this matrix with the content of the given iterator.
129
    ///
130
    /// This will fill as many diagonal elements as the iterator yields, up to the
131
    /// minimum of the number of rows and columns of `self`, and starting with the
132
    /// diagonal element at index (0, 0).
133
    #[inline]
134
0
    pub fn set_partial_diagonal(&mut self, diag: impl Iterator<Item = T>) {
135
0
        let (nrows, ncols) = self.shape();
136
0
        let min_nrows_ncols = cmp::min(nrows, ncols);
137
138
0
        for (i, val) in diag.enumerate().take(min_nrows_ncols) {
139
0
            unsafe { *self.get_unchecked_mut((i, i)) = val }
140
        }
141
0
    }
142
143
    /// Fills the selected row of this matrix with the content of the given vector.
144
    #[inline]
145
    pub fn set_row<C2: Dim, S2>(&mut self, i: usize, row: &RowVector<T, C2, S2>)
146
    where
147
        S2: RawStorage<T, U1, C2>,
148
        ShapeConstraint: SameNumberOfColumns<C, C2>,
149
    {
150
        self.row_mut(i).copy_from(row);
151
    }
152
153
    /// Fills the selected column of this matrix with the content of the given vector.
154
    #[inline]
155
    pub fn set_column<R2: Dim, S2>(&mut self, i: usize, column: &Vector<T, R2, S2>)
156
    where
157
        S2: RawStorage<T, R2, U1>,
158
        ShapeConstraint: SameNumberOfRows<R, R2>,
159
    {
160
        self.column_mut(i).copy_from(column);
161
    }
162
}
163
164
/// # In-place filling
165
impl<T, R: Dim, C: Dim, S: RawStorageMut<T, R, C>> Matrix<T, R, C, S> {
166
    /// Sets all the elements of this matrix to the value returned by the closure.
167
    #[inline]
168
    pub fn fill_with(&mut self, val: impl Fn() -> T) {
169
        for e in self.iter_mut() {
170
            *e = val()
171
        }
172
    }
173
174
    /// Sets all the elements of this matrix to `val`.
175
    #[inline]
176
    pub fn fill(&mut self, val: T)
177
    where
178
        T: Scalar,
179
    {
180
        for e in self.iter_mut() {
181
            *e = val.clone()
182
        }
183
    }
184
185
    /// Fills `self` with the identity matrix.
186
    #[inline]
187
    pub fn fill_with_identity(&mut self)
188
    where
189
        T: Scalar + Zero + One,
190
    {
191
        self.fill(T::zero());
192
        self.fill_diagonal(T::one());
193
    }
194
195
    /// Sets all the diagonal elements of this matrix to `val`.
196
    #[inline]
197
    pub fn fill_diagonal(&mut self, val: T)
198
    where
199
        T: Scalar,
200
    {
201
        let (nrows, ncols) = self.shape();
202
        let n = cmp::min(nrows, ncols);
203
204
        for i in 0..n {
205
            unsafe { *self.get_unchecked_mut((i, i)) = val.clone() }
206
        }
207
    }
208
209
    /// Sets all the elements of the selected row to `val`.
210
    #[inline]
211
    pub fn fill_row(&mut self, i: usize, val: T)
212
    where
213
        T: Scalar,
214
    {
215
        assert!(i < self.nrows(), "Row index out of bounds.");
216
        for j in 0..self.ncols() {
217
            unsafe { *self.get_unchecked_mut((i, j)) = val.clone() }
218
        }
219
    }
220
221
    /// Sets all the elements of the selected column to `val`.
222
    #[inline]
223
    pub fn fill_column(&mut self, j: usize, val: T)
224
    where
225
        T: Scalar,
226
    {
227
        assert!(j < self.ncols(), "Row index out of bounds.");
228
        for i in 0..self.nrows() {
229
            unsafe { *self.get_unchecked_mut((i, j)) = val.clone() }
230
        }
231
    }
232
233
    /// Sets all the elements of the lower-triangular part of this matrix to `val`.
234
    ///
235
    /// The parameter `shift` allows some subdiagonals to be left untouched:
236
    /// * If `shift = 0` then the diagonal is overwritten as well.
237
    /// * If `shift = 1` then the diagonal is left untouched.
238
    /// * If `shift > 1`, then the diagonal and the first `shift - 1` subdiagonals are left
239
    ///   untouched.
240
    #[inline]
241
0
    pub fn fill_lower_triangle(&mut self, val: T, shift: usize)
242
0
    where
243
0
        T: Scalar,
244
    {
245
0
        for j in 0..self.ncols() {
246
0
            for i in (j + shift)..self.nrows() {
247
0
                unsafe { *self.get_unchecked_mut((i, j)) = val.clone() }
248
            }
249
        }
250
0
    }
251
252
    /// Sets all the elements of the upper-triangular part of this matrix to `val`.
253
    ///
254
    /// The parameter `shift` allows some superdiagonals to be left untouched:
255
    /// * If `shift = 0` then the diagonal is overwritten as well.
256
    /// * If `shift = 1` then the diagonal is left untouched.
257
    /// * If `shift > 1`, then the diagonal and the first `shift - 1` superdiagonals are left
258
    ///   untouched.
259
    #[inline]
260
    pub fn fill_upper_triangle(&mut self, val: T, shift: usize)
261
    where
262
        T: Scalar,
263
    {
264
        for j in shift..self.ncols() {
265
            // TODO: is there a more efficient way to avoid the min ?
266
            // (necessary for rectangular matrices)
267
            for i in 0..cmp::min(j + 1 - shift, self.nrows()) {
268
                unsafe { *self.get_unchecked_mut((i, j)) = val.clone() }
269
            }
270
        }
271
    }
272
}
273
274
impl<T: Scalar, D: Dim, S: RawStorageMut<T, D, D>> Matrix<T, D, D, S> {
275
    /// Copies the upper-triangle of this matrix to its lower-triangular part.
276
    ///
277
    /// This makes the matrix symmetric. Panics if the matrix is not square.
278
    pub fn fill_lower_triangle_with_upper_triangle(&mut self) {
279
        assert!(self.is_square(), "The input matrix should be square.");
280
281
        let dim = self.nrows();
282
        for j in 0..dim {
283
            for i in j + 1..dim {
284
                unsafe {
285
                    *self.get_unchecked_mut((i, j)) = self.get_unchecked((j, i)).clone();
286
                }
287
            }
288
        }
289
    }
290
291
    /// Copies the upper-triangle of this matrix to its upper-triangular part.
292
    ///
293
    /// This makes the matrix symmetric. Panics if the matrix is not square.
294
    pub fn fill_upper_triangle_with_lower_triangle(&mut self) {
295
        assert!(self.is_square(), "The input matrix should be square.");
296
297
        for j in 1..self.ncols() {
298
            for i in 0..j {
299
                unsafe {
300
                    *self.get_unchecked_mut((i, j)) = self.get_unchecked((j, i)).clone();
301
                }
302
            }
303
        }
304
    }
305
}
306
307
/// # In-place swapping
308
impl<T: Scalar, R: Dim, C: Dim, S: RawStorageMut<T, R, C>> Matrix<T, R, C, S> {
309
    /// Swaps two rows in-place.
310
    #[inline]
311
0
    pub fn swap_rows(&mut self, irow1: usize, irow2: usize) {
312
0
        assert!(irow1 < self.nrows() && irow2 < self.nrows());
313
314
0
        if irow1 != irow2 {
315
            // TODO: optimize that.
316
0
            for i in 0..self.ncols() {
317
0
                unsafe { self.swap_unchecked((irow1, i), (irow2, i)) }
318
            }
319
0
        }
320
        // Otherwise do nothing.
321
0
    }
322
323
    /// Swaps two columns in-place.
324
    #[inline]
325
    pub fn swap_columns(&mut self, icol1: usize, icol2: usize) {
326
        assert!(icol1 < self.ncols() && icol2 < self.ncols());
327
328
        if icol1 != icol2 {
329
            // TODO: optimize that.
330
            for i in 0..self.nrows() {
331
                unsafe { self.swap_unchecked((i, icol1), (i, icol2)) }
332
            }
333
        }
334
        // Otherwise do nothing.
335
    }
336
}
337
338
/*
339
 *
340
 * TODO: specialize all the following for slices.
341
 *
342
 */
343
/// # Rows and columns removal
344
impl<T: Scalar, R: Dim, C: Dim, S: Storage<T, R, C>> Matrix<T, R, C, S> {
345
    /*
346
     *
347
     * Column removal.
348
     *
349
     */
350
    /// Removes the `i`-th column from this matrix.
351
    #[inline]
352
    pub fn remove_column(self, i: usize) -> OMatrix<T, R, DimDiff<C, U1>>
353
    where
354
        C: DimSub<U1>,
355
        DefaultAllocator: Reallocator<T, R, C, R, DimDiff<C, U1>>,
356
    {
357
        self.remove_fixed_columns::<1>(i)
358
    }
359
360
    /// Removes all columns in `indices`   
361
    #[cfg(any(feature = "std", feature = "alloc"))]
362
    pub fn remove_columns_at(self, indices: &[usize]) -> OMatrix<T, R, Dyn>
363
    where
364
        C: DimSub<Dyn, Output = Dyn>,
365
        DefaultAllocator: Reallocator<T, R, C, R, Dyn>,
366
    {
367
        let mut m = self.into_owned();
368
        let (nrows, ncols) = m.shape_generic();
369
        let mut offset: usize = 0;
370
        let mut target: usize = 0;
371
        while offset + target < ncols.value() {
372
            if indices.contains(&(target + offset)) {
373
                // Safety: the resulting pointer is within range.
374
                let col_ptr = unsafe { m.data.ptr_mut().add((target + offset) * nrows.value()) };
375
                // Drop every element in the column we are about to overwrite.
376
                // We use the a similar technique as in `Vec::truncate`.
377
                let s = ptr::slice_from_raw_parts_mut(col_ptr, nrows.value());
378
                // Safety: we drop the column in-place, which is OK because we will overwrite these
379
                //         entries later in the loop, or discard them with the `reallocate_copy`
380
                //         afterwards.
381
                unsafe { ptr::drop_in_place(s) };
382
383
                offset += 1;
384
            } else {
385
                unsafe {
386
                    let ptr_source = m.data.ptr().add((target + offset) * nrows.value());
387
                    let ptr_target = m.data.ptr_mut().add(target * nrows.value());
388
389
                    // Copy the data, overwriting what we dropped.
390
                    ptr::copy(ptr_source, ptr_target, nrows.value());
391
                    target += 1;
392
                }
393
            }
394
        }
395
396
        // Safety: The new size is smaller than the old size, so
397
        //         DefaultAllocator::reallocate_copy will initialize
398
        //         every element of the new matrix which can then
399
        //         be assumed to be initialized.
400
        unsafe {
401
            let new_data = DefaultAllocator::reallocate_copy(
402
                nrows,
403
                ncols.sub(Dyn::from_usize(offset)),
404
                m.data,
405
            );
406
407
            Matrix::from_data(new_data).assume_init()
408
        }
409
    }
410
411
    /// Removes all rows in `indices`   
412
    #[cfg(any(feature = "std", feature = "alloc"))]
413
    pub fn remove_rows_at(self, indices: &[usize]) -> OMatrix<T, Dyn, C>
414
    where
415
        R: DimSub<Dyn, Output = Dyn>,
416
        DefaultAllocator: Reallocator<T, R, C, Dyn, C>,
417
    {
418
        let mut m = self.into_owned();
419
        let (nrows, ncols) = m.shape_generic();
420
        let mut offset: usize = 0;
421
        let mut target: usize = 0;
422
        while offset + target < nrows.value() * ncols.value() {
423
            if indices.contains(&((target + offset) % nrows.value())) {
424
                // Safety: the resulting pointer is within range.
425
                unsafe {
426
                    let elt_ptr = m.data.ptr_mut().add(target + offset);
427
                    // Safety: we drop the component in-place, which is OK because we will overwrite these
428
                    //         entries later in the loop, or discard them with the `reallocate_copy`
429
                    //         afterwards.
430
                    ptr::drop_in_place(elt_ptr)
431
                };
432
                offset += 1;
433
            } else {
434
                unsafe {
435
                    let ptr_source = m.data.ptr().add(target + offset);
436
                    let ptr_target = m.data.ptr_mut().add(target);
437
438
                    // Copy the data, overwriting what we dropped in the previous iterations.
439
                    ptr::copy(ptr_source, ptr_target, 1);
440
                    target += 1;
441
                }
442
            }
443
        }
444
445
        // Safety: The new size is smaller than the old size, so
446
        //         DefaultAllocator::reallocate_copy will initialize
447
        //         every element of the new matrix which can then
448
        //         be assumed to be initialized.
449
        unsafe {
450
            let new_data = DefaultAllocator::reallocate_copy(
451
                nrows.sub(Dyn::from_usize(offset / ncols.value())),
452
                ncols,
453
                m.data,
454
            );
455
456
            Matrix::from_data(new_data).assume_init()
457
        }
458
    }
459
460
    /// Removes `D::dim()` consecutive columns from this matrix, starting with the `i`-th
461
    /// (included).
462
    #[inline]
463
    pub fn remove_fixed_columns<const D: usize>(
464
        self,
465
        i: usize,
466
    ) -> OMatrix<T, R, DimDiff<C, Const<D>>>
467
    where
468
        C: DimSub<Const<D>>,
469
        DefaultAllocator: Reallocator<T, R, C, R, DimDiff<C, Const<D>>>,
470
    {
471
        self.remove_columns_generic(i, Const::<D>)
472
    }
473
474
    /// Removes `n` consecutive columns from this matrix, starting with the `i`-th (included).
475
    #[inline]
476
    #[cfg(any(feature = "std", feature = "alloc"))]
477
    pub fn remove_columns(self, i: usize, n: usize) -> OMatrix<T, R, Dyn>
478
    where
479
        C: DimSub<Dyn, Output = Dyn>,
480
        DefaultAllocator: Reallocator<T, R, C, R, Dyn>,
481
    {
482
        self.remove_columns_generic(i, Dyn(n))
483
    }
484
485
    /// Removes `nremove.value()` columns from this matrix, starting with the `i`-th (included).
486
    ///
487
    /// This is the generic implementation of `.remove_columns(...)` and
488
    /// `.remove_fixed_columns(...)` which have nicer API interfaces.
489
    #[inline]
490
    pub fn remove_columns_generic<D>(self, i: usize, nremove: D) -> OMatrix<T, R, DimDiff<C, D>>
491
    where
492
        D: Dim,
493
        C: DimSub<D>,
494
        DefaultAllocator: Reallocator<T, R, C, R, DimDiff<C, D>>,
495
    {
496
        let mut m = self.into_owned();
497
        let (nrows, ncols) = m.shape_generic();
498
        assert!(
499
            i + nremove.value() <= ncols.value(),
500
            "Column index out of range."
501
        );
502
503
        let need_column_shifts = nremove.value() != 0 && i + nremove.value() < ncols.value();
504
        if need_column_shifts {
505
            // The first `deleted_i * nrows` are left untouched.
506
            let copied_value_start = i + nremove.value();
507
508
            unsafe {
509
                let ptr_in = m.data.ptr().add(copied_value_start * nrows.value());
510
                let ptr_out = m.data.ptr_mut().add(i * nrows.value());
511
512
                // Drop all the elements of the columns we are about to overwrite.
513
                // We use the a similar technique as in `Vec::truncate`.
514
                let s = ptr::slice_from_raw_parts_mut(ptr_out, nremove.value() * nrows.value());
515
                // Safety: we drop the column in-place, which is OK because we will overwrite these
516
                //         entries with `ptr::copy` afterward.
517
                ptr::drop_in_place(s);
518
519
                ptr::copy(
520
                    ptr_in,
521
                    ptr_out,
522
                    (ncols.value() - copied_value_start) * nrows.value(),
523
                );
524
            }
525
        } else {
526
            // All the columns to remove are at the end of the buffer. Drop them.
527
            unsafe {
528
                let ptr = m.data.ptr_mut().add(i * nrows.value());
529
                let s = ptr::slice_from_raw_parts_mut(ptr, nremove.value() * nrows.value());
530
                ptr::drop_in_place(s)
531
            };
532
        }
533
534
        // Safety: The new size is smaller than the old size, so
535
        //         DefaultAllocator::reallocate_copy will initialize
536
        //         every element of the new matrix which can then
537
        //         be assumed to be initialized.
538
        unsafe {
539
            let new_data = DefaultAllocator::reallocate_copy(nrows, ncols.sub(nremove), m.data);
540
            Matrix::from_data(new_data).assume_init()
541
        }
542
    }
543
544
    /*
545
     *
546
     * Row removal.
547
     *
548
     */
549
    /// Removes the `i`-th row from this matrix.
550
    #[inline]
551
    pub fn remove_row(self, i: usize) -> OMatrix<T, DimDiff<R, U1>, C>
552
    where
553
        R: DimSub<U1>,
554
        DefaultAllocator: Reallocator<T, R, C, DimDiff<R, U1>, C>,
555
    {
556
        self.remove_fixed_rows::<1>(i)
557
    }
558
559
    /// Removes `D::dim()` consecutive rows from this matrix, starting with the `i`-th (included).
560
    #[inline]
561
    pub fn remove_fixed_rows<const D: usize>(self, i: usize) -> OMatrix<T, DimDiff<R, Const<D>>, C>
562
    where
563
        R: DimSub<Const<D>>,
564
        DefaultAllocator: Reallocator<T, R, C, DimDiff<R, Const<D>>, C>,
565
    {
566
        self.remove_rows_generic(i, Const::<D>)
567
    }
568
569
    /// Removes `n` consecutive rows from this matrix, starting with the `i`-th (included).
570
    #[inline]
571
    #[cfg(any(feature = "std", feature = "alloc"))]
572
    pub fn remove_rows(self, i: usize, n: usize) -> OMatrix<T, Dyn, C>
573
    where
574
        R: DimSub<Dyn, Output = Dyn>,
575
        DefaultAllocator: Reallocator<T, R, C, Dyn, C>,
576
    {
577
        self.remove_rows_generic(i, Dyn(n))
578
    }
579
580
    /// Removes `nremove.value()` rows from this matrix, starting with the `i`-th (included).
581
    ///
582
    /// This is the generic implementation of `.remove_rows(...)` and `.remove_fixed_rows(...)`
583
    /// which have nicer API interfaces.
584
    #[inline]
585
    pub fn remove_rows_generic<D>(self, i: usize, nremove: D) -> OMatrix<T, DimDiff<R, D>, C>
586
    where
587
        D: Dim,
588
        R: DimSub<D>,
589
        DefaultAllocator: Reallocator<T, R, C, DimDiff<R, D>, C>,
590
    {
591
        let mut m = self.into_owned();
592
        let (nrows, ncols) = m.shape_generic();
593
        assert!(
594
            i + nremove.value() <= nrows.value(),
595
            "Row index out of range."
596
        );
597
598
        if nremove.value() != 0 {
599
            unsafe {
600
                compress_rows(
601
                    m.as_mut_slice(),
602
                    nrows.value(),
603
                    ncols.value(),
604
                    i,
605
                    nremove.value(),
606
                );
607
            }
608
        }
609
610
        // Safety: The new size is smaller than the old size, so
611
        //         DefaultAllocator::reallocate_copy will initialize
612
        //         every element of the new matrix which can then
613
        //         be assumed to be initialized.
614
        unsafe {
615
            let new_data = DefaultAllocator::reallocate_copy(nrows.sub(nremove), ncols, m.data);
616
            Matrix::from_data(new_data).assume_init()
617
        }
618
    }
619
}
620
621
/// # Rows and columns insertion
622
impl<T: Scalar, R: Dim, C: Dim, S: Storage<T, R, C>> Matrix<T, R, C, S> {
623
    /*
624
     *
625
     * Columns insertion.
626
     *
627
     */
628
    /// Inserts a column filled with `val` at the `i-th` position.
629
    #[inline]
630
    pub fn insert_column(self, i: usize, val: T) -> OMatrix<T, R, DimSum<C, U1>>
631
    where
632
        C: DimAdd<U1>,
633
        DefaultAllocator: Reallocator<T, R, C, R, DimSum<C, U1>>,
634
    {
635
        self.insert_fixed_columns::<1>(i, val)
636
    }
637
638
    /// Inserts `D` columns filled with `val` starting at the `i-th` position.
639
    #[inline]
640
    pub fn insert_fixed_columns<const D: usize>(
641
        self,
642
        i: usize,
643
        val: T,
644
    ) -> OMatrix<T, R, DimSum<C, Const<D>>>
645
    where
646
        C: DimAdd<Const<D>>,
647
        DefaultAllocator: Reallocator<T, R, C, R, DimSum<C, Const<D>>>,
648
    {
649
        let mut res = unsafe { self.insert_columns_generic_uninitialized(i, Const::<D>) };
650
        res.fixed_columns_mut::<D>(i)
651
            .fill_with(|| MaybeUninit::new(val.clone()));
652
653
        // Safety: the result is now fully initialized. The added columns have
654
        //         been initialized by the `fill_with` above, and the rest have
655
        //         been initialized by `insert_columns_generic_uninitialized`.
656
        unsafe { res.assume_init() }
657
    }
658
659
    /// Inserts `n` columns filled with `val` starting at the `i-th` position.
660
    #[inline]
661
    #[cfg(any(feature = "std", feature = "alloc"))]
662
    pub fn insert_columns(self, i: usize, n: usize, val: T) -> OMatrix<T, R, Dyn>
663
    where
664
        C: DimAdd<Dyn, Output = Dyn>,
665
        DefaultAllocator: Reallocator<T, R, C, R, Dyn>,
666
    {
667
        let mut res = unsafe { self.insert_columns_generic_uninitialized(i, Dyn(n)) };
668
        res.columns_mut(i, n)
669
            .fill_with(|| MaybeUninit::new(val.clone()));
670
671
        // Safety: the result is now fully initialized. The added columns have
672
        //         been initialized by the `fill_with` above, and the rest have
673
        //         been initialized by `insert_columns_generic_uninitialized`.
674
        unsafe { res.assume_init() }
675
    }
676
677
    /// Inserts `ninsert.value()` columns starting at the `i-th` place of this matrix.
678
    ///
679
    /// # Safety
680
    /// The output matrix has all its elements initialized except for the the components of the
681
    /// added columns.
682
    #[inline]
683
    pub unsafe fn insert_columns_generic_uninitialized<D>(
684
        self,
685
        i: usize,
686
        ninsert: D,
687
    ) -> UninitMatrix<T, R, DimSum<C, D>>
688
    where
689
        D: Dim,
690
        C: DimAdd<D>,
691
        DefaultAllocator: Reallocator<T, R, C, R, DimSum<C, D>>,
692
    {
693
        unsafe {
694
            let m = self.into_owned();
695
            let (nrows, ncols) = m.shape_generic();
696
            let mut res = Matrix::from_data(DefaultAllocator::reallocate_copy(
697
                nrows,
698
                ncols.add(ninsert),
699
                m.data,
700
            ));
701
702
            assert!(i <= ncols.value(), "Column insertion index out of range.");
703
704
            if ninsert.value() != 0 && i != ncols.value() {
705
                let ptr_in = res.data.ptr().add(i * nrows.value());
706
                let ptr_out = res
707
                    .data
708
                    .ptr_mut()
709
                    .add((i + ninsert.value()) * nrows.value());
710
711
                ptr::copy(ptr_in, ptr_out, (ncols.value() - i) * nrows.value())
712
            }
713
714
            res
715
        }
716
    }
717
718
    /*
719
     *
720
     * Rows insertion.
721
     *
722
     */
723
    /// Inserts a row filled with `val` at the `i-th` position.
724
    #[inline]
725
    pub fn insert_row(self, i: usize, val: T) -> OMatrix<T, DimSum<R, U1>, C>
726
    where
727
        R: DimAdd<U1>,
728
        DefaultAllocator: Reallocator<T, R, C, DimSum<R, U1>, C>,
729
    {
730
        self.insert_fixed_rows::<1>(i, val)
731
    }
732
733
    /// Inserts `D::dim()` rows filled with `val` starting at the `i-th` position.
734
    #[inline]
735
    pub fn insert_fixed_rows<const D: usize>(
736
        self,
737
        i: usize,
738
        val: T,
739
    ) -> OMatrix<T, DimSum<R, Const<D>>, C>
740
    where
741
        R: DimAdd<Const<D>>,
742
        DefaultAllocator: Reallocator<T, R, C, DimSum<R, Const<D>>, C>,
743
    {
744
        let mut res = unsafe { self.insert_rows_generic_uninitialized(i, Const::<D>) };
745
        res.fixed_rows_mut::<D>(i)
746
            .fill_with(|| MaybeUninit::new(val.clone()));
747
748
        // Safety: the result is now fully initialized. The added rows have
749
        //         been initialized by the `fill_with` above, and the rest have
750
        //         been initialized by `insert_rows_generic_uninitialized`.
751
        unsafe { res.assume_init() }
752
    }
753
754
    /// Inserts `n` rows filled with `val` starting at the `i-th` position.
755
    #[inline]
756
    #[cfg(any(feature = "std", feature = "alloc"))]
757
    pub fn insert_rows(self, i: usize, n: usize, val: T) -> OMatrix<T, Dyn, C>
758
    where
759
        R: DimAdd<Dyn, Output = Dyn>,
760
        DefaultAllocator: Reallocator<T, R, C, Dyn, C>,
761
    {
762
        let mut res = unsafe { self.insert_rows_generic_uninitialized(i, Dyn(n)) };
763
        res.rows_mut(i, n)
764
            .fill_with(|| MaybeUninit::new(val.clone()));
765
766
        // Safety: the result is now fully initialized. The added rows have
767
        //         been initialized by the `fill_with` above, and the rest have
768
        //         been initialized by `insert_rows_generic_uninitialized`.
769
        unsafe { res.assume_init() }
770
    }
771
772
    /// Inserts `ninsert.value()` rows at the `i-th` place of this matrix.
773
    ///
774
    /// # Safety
775
    /// The added rows values are not initialized.
776
    /// This is the generic implementation of `.insert_rows(...)` and
777
    /// `.insert_fixed_rows(...)` which have nicer API interfaces.
778
    #[inline]
779
    pub unsafe fn insert_rows_generic_uninitialized<D>(
780
        self,
781
        i: usize,
782
        ninsert: D,
783
    ) -> UninitMatrix<T, DimSum<R, D>, C>
784
    where
785
        D: Dim,
786
        R: DimAdd<D>,
787
        DefaultAllocator: Reallocator<T, R, C, DimSum<R, D>, C>,
788
    {
789
        unsafe {
790
            let m = self.into_owned();
791
            let (nrows, ncols) = m.shape_generic();
792
            let mut res = Matrix::from_data(DefaultAllocator::reallocate_copy(
793
                nrows.add(ninsert),
794
                ncols,
795
                m.data,
796
            ));
797
798
            assert!(i <= nrows.value(), "Row insertion index out of range.");
799
800
            if ninsert.value() != 0 {
801
                extend_rows(
802
                    res.as_mut_slice(),
803
                    nrows.value(),
804
                    ncols.value(),
805
                    i,
806
                    ninsert.value(),
807
                );
808
            }
809
810
            res
811
        }
812
    }
813
}
814
815
/// # Resizing and reshaping
816
impl<T: Scalar, R: Dim, C: Dim, S: Storage<T, R, C>> Matrix<T, R, C, S> {
817
    /// Resizes this matrix so that it contains `new_nrows` rows and `new_ncols` columns.
818
    ///
819
    /// The values are copied such that `self[(i, j)] == result[(i, j)]`. If the result has more
820
    /// rows and/or columns than `self`, then the extra rows or columns are filled with `val`.
821
    #[cfg(any(feature = "std", feature = "alloc"))]
822
    pub fn resize(self, new_nrows: usize, new_ncols: usize, val: T) -> OMatrix<T, Dyn, Dyn>
823
    where
824
        DefaultAllocator: Reallocator<T, R, C, Dyn, Dyn>,
825
    {
826
        self.resize_generic(Dyn(new_nrows), Dyn(new_ncols), val)
827
    }
828
829
    /// Resizes this matrix vertically, i.e., so that it contains `new_nrows` rows while keeping the same number of columns.
830
    ///
831
    /// The values are copied such that `self[(i, j)] == result[(i, j)]`. If the result has more
832
    /// rows than `self`, then the extra rows are filled with `val`.
833
    #[cfg(any(feature = "std", feature = "alloc"))]
834
    pub fn resize_vertically(self, new_nrows: usize, val: T) -> OMatrix<T, Dyn, C>
835
    where
836
        DefaultAllocator: Reallocator<T, R, C, Dyn, C>,
837
    {
838
        let ncols = self.shape_generic().1;
839
        self.resize_generic(Dyn(new_nrows), ncols, val)
840
    }
841
842
    /// Resizes this matrix horizontally, i.e., so that it contains `new_ncolumns` columns while keeping the same number of columns.
843
    ///
844
    /// The values are copied such that `self[(i, j)] == result[(i, j)]`. If the result has more
845
    /// columns than `self`, then the extra columns are filled with `val`.
846
    #[cfg(any(feature = "std", feature = "alloc"))]
847
    pub fn resize_horizontally(self, new_ncols: usize, val: T) -> OMatrix<T, R, Dyn>
848
    where
849
        DefaultAllocator: Reallocator<T, R, C, R, Dyn>,
850
    {
851
        let nrows = self.shape_generic().0;
852
        self.resize_generic(nrows, Dyn(new_ncols), val)
853
    }
854
855
    /// Resizes this matrix so that it contains `R2::value()` rows and `C2::value()` columns.
856
    ///
857
    /// The values are copied such that `self[(i, j)] == result[(i, j)]`. If the result has more
858
    /// rows and/or columns than `self`, then the extra rows or columns are filled with `val`.
859
    pub fn fixed_resize<const R2: usize, const C2: usize>(
860
        self,
861
        val: T,
862
    ) -> OMatrix<T, Const<R2>, Const<C2>>
863
    where
864
        DefaultAllocator: Reallocator<T, R, C, Const<R2>, Const<C2>>,
865
    {
866
        self.resize_generic(Const::<R2>, Const::<C2>, val)
867
    }
868
869
    /// Resizes `self` such that it has dimensions `new_nrows × new_ncols`.
870
    ///
871
    /// The values are copied such that `self[(i, j)] == result[(i, j)]`. If the result has more
872
    /// rows and/or columns than `self`, then the extra rows or columns are filled with `val`.
873
    #[inline]
874
    pub fn resize_generic<R2: Dim, C2: Dim>(
875
        self,
876
        new_nrows: R2,
877
        new_ncols: C2,
878
        val: T,
879
    ) -> OMatrix<T, R2, C2>
880
    where
881
        DefaultAllocator: Reallocator<T, R, C, R2, C2>,
882
    {
883
        let (nrows, ncols) = self.shape();
884
        let mut data = self.into_owned();
885
886
        if new_nrows.value() == nrows {
887
            if new_ncols.value() < ncols {
888
                unsafe {
889
                    let num_cols_to_delete = ncols - new_ncols.value();
890
                    let col_ptr = data.data.ptr_mut().add(new_ncols.value() * nrows);
891
                    let s = ptr::slice_from_raw_parts_mut(col_ptr, num_cols_to_delete * nrows);
892
                    // Safety: drop the elements of the deleted columns.
893
                    //         these are the elements that will be truncated
894
                    //         by the `reallocate_copy` afterward.
895
                    ptr::drop_in_place(s)
896
                };
897
            }
898
899
            let res = unsafe { DefaultAllocator::reallocate_copy(new_nrows, new_ncols, data.data) };
900
            let mut res = Matrix::from_data(res);
901
902
            if new_ncols.value() > ncols {
903
                res.columns_range_mut(ncols..)
904
                    .fill_with(|| MaybeUninit::new(val.clone()));
905
            }
906
907
            // Safety: the result is now fully initialized by `reallocate_copy` and
908
            //         `fill_with` (if the output has more columns than the input).
909
            unsafe { res.assume_init() }
910
        } else {
911
            let mut res;
912
913
            unsafe {
914
                if new_nrows.value() < nrows {
915
                    compress_rows(
916
                        data.as_mut_slice(),
917
                        nrows,
918
                        ncols,
919
                        new_nrows.value(),
920
                        nrows - new_nrows.value(),
921
                    );
922
                    res = Matrix::from_data(DefaultAllocator::reallocate_copy(
923
                        new_nrows, new_ncols, data.data,
924
                    ));
925
                } else {
926
                    res = Matrix::from_data(DefaultAllocator::reallocate_copy(
927
                        new_nrows, new_ncols, data.data,
928
                    ));
929
                    extend_rows(
930
                        res.as_mut_slice(),
931
                        nrows,
932
                        new_ncols.value(),
933
                        nrows,
934
                        new_nrows.value() - nrows,
935
                    );
936
                }
937
            }
938
939
            if new_ncols.value() > ncols {
940
                res.columns_range_mut(ncols..)
941
                    .fill_with(|| MaybeUninit::new(val.clone()));
942
            }
943
944
            if new_nrows.value() > nrows {
945
                res.view_range_mut(nrows.., ..cmp::min(ncols, new_ncols.value()))
946
                    .fill_with(|| MaybeUninit::new(val.clone()));
947
            }
948
949
            // Safety: the result is now fully initialized by `reallocate_copy` and
950
            //         `fill_with` (whenever applicable).
951
            unsafe { res.assume_init() }
952
        }
953
    }
954
955
    /// Reshapes `self` such that it has dimensions `new_nrows × new_ncols`.
956
    ///
957
    /// This will reinterpret `self` as if it is a matrix with `new_nrows` rows and `new_ncols`
958
    /// columns. The arrangements of the component in the output matrix are the same as what
959
    /// would be obtained by `Matrix::from_slice_generic(self.as_slice(), new_nrows, new_ncols)`.
960
    ///
961
    /// If `self` is a dynamically-sized matrix, then its components are neither copied nor moved.
962
    /// If `self` is staticyll-sized, then a copy may happen in some situations.
963
    /// This function will panic if the given dimensions are such that the number of elements of
964
    /// the input matrix are not equal to the number of elements of the output matrix.
965
    ///
966
    /// # Examples
967
    ///
968
    /// ```
969
    /// # use nalgebra::{Matrix3x2, Matrix2x3, DMatrix, Const, Dyn};
970
    ///
971
    /// let m1 = Matrix2x3::new(
972
    ///     1.1, 1.2, 1.3,
973
    ///     2.1, 2.2, 2.3
974
    /// );
975
    /// let m2 = Matrix3x2::new(
976
    ///     1.1, 2.2,
977
    ///     2.1, 1.3,
978
    ///     1.2, 2.3
979
    /// );
980
    /// let reshaped = m1.reshape_generic(Const::<3>, Const::<2>);
981
    /// assert_eq!(reshaped, m2);
982
    ///
983
    /// let dm1 = DMatrix::from_row_slice(
984
    ///     4,
985
    ///     3,
986
    ///     &[
987
    ///         1.0, 0.0, 0.0,
988
    ///         0.0, 0.0, 1.0,
989
    ///         0.0, 0.0, 0.0,
990
    ///         0.0, 1.0, 0.0
991
    ///     ],
992
    /// );
993
    /// let dm2 = DMatrix::from_row_slice(
994
    ///     6,
995
    ///     2,
996
    ///     &[
997
    ///         1.0, 0.0,
998
    ///         0.0, 1.0,
999
    ///         0.0, 0.0,
1000
    ///         0.0, 1.0,
1001
    ///         0.0, 0.0,
1002
    ///         0.0, 0.0,
1003
    ///     ],
1004
    /// );
1005
    /// let reshaped = dm1.reshape_generic(Dyn(6), Dyn(2));
1006
    /// assert_eq!(reshaped, dm2);
1007
    /// ```
1008
    pub fn reshape_generic<R2, C2>(
1009
        self,
1010
        new_nrows: R2,
1011
        new_ncols: C2,
1012
    ) -> Matrix<T, R2, C2, S::Output>
1013
    where
1014
        R2: Dim,
1015
        C2: Dim,
1016
        S: ReshapableStorage<T, R, C, R2, C2>,
1017
    {
1018
        let data = self.data.reshape_generic(new_nrows, new_ncols);
1019
        Matrix::from_data(data)
1020
    }
1021
}
1022
1023
/// # In-place resizing
1024
#[cfg(any(feature = "std", feature = "alloc"))]
1025
impl<T: Scalar> OMatrix<T, Dyn, Dyn> {
1026
    /// Resizes this matrix in-place.
1027
    ///
1028
    /// The values are copied such that `self[(i, j)] == result[(i, j)]`. If the result has more
1029
    /// rows and/or columns than `self`, then the extra rows or columns are filled with `val`.
1030
    ///
1031
    /// Defined only for owned fully-dynamic matrices, i.e., `DMatrix`.
1032
    pub fn resize_mut(&mut self, new_nrows: usize, new_ncols: usize, val: T)
1033
    where
1034
        DefaultAllocator: Reallocator<T, Dyn, Dyn, Dyn, Dyn>,
1035
    {
1036
        // TODO: avoid the clone.
1037
        *self = self.clone().resize(new_nrows, new_ncols, val);
1038
    }
1039
}
1040
1041
#[cfg(any(feature = "std", feature = "alloc"))]
1042
impl<T: Scalar, C: Dim> OMatrix<T, Dyn, C>
1043
where
1044
    DefaultAllocator: Allocator<Dyn, C>,
1045
{
1046
    /// Changes the number of rows of this matrix in-place.
1047
    ///
1048
    /// The values are copied such that `self[(i, j)] == result[(i, j)]`. If the result has more
1049
    /// rows than `self`, then the extra rows are filled with `val`.
1050
    ///
1051
    /// Defined only for owned matrices with a dynamic number of rows (for example, `DVector`).
1052
    #[cfg(any(feature = "std", feature = "alloc"))]
1053
    pub fn resize_vertically_mut(&mut self, new_nrows: usize, val: T)
1054
    where
1055
        DefaultAllocator: Reallocator<T, Dyn, C, Dyn, C>,
1056
    {
1057
        // TODO: avoid the clone.
1058
        *self = self.clone().resize_vertically(new_nrows, val);
1059
    }
1060
}
1061
1062
#[cfg(any(feature = "std", feature = "alloc"))]
1063
impl<T: Scalar, R: Dim> OMatrix<T, R, Dyn>
1064
where
1065
    DefaultAllocator: Allocator<R, Dyn>,
1066
{
1067
    /// Changes the number of column of this matrix in-place.
1068
    ///
1069
    /// The values are copied such that `self[(i, j)] == result[(i, j)]`. If the result has more
1070
    /// columns than `self`, then the extra columns are filled with `val`.
1071
    ///
1072
    /// Defined only for owned matrices with a dynamic number of columns (for example, `DVector`).
1073
    #[cfg(any(feature = "std", feature = "alloc"))]
1074
    pub fn resize_horizontally_mut(&mut self, new_ncols: usize, val: T)
1075
    where
1076
        DefaultAllocator: Reallocator<T, R, Dyn, R, Dyn>,
1077
    {
1078
        // TODO: avoid the clone.
1079
        *self = self.clone().resize_horizontally(new_ncols, val);
1080
    }
1081
}
1082
1083
// Move the elements of `data` in such a way that the matrix with
1084
// the rows `[i, i + nremove[` deleted is represented in a contiguous
1085
// way in `data` after this method completes.
1086
// Every deleted element are manually dropped by this method.
1087
unsafe fn compress_rows<T: Scalar>(
1088
    data: &mut [T],
1089
    nrows: usize,
1090
    ncols: usize,
1091
    i: usize,
1092
    nremove: usize,
1093
) {
1094
    unsafe {
1095
        let new_nrows = nrows - nremove;
1096
1097
        if nremove == 0 {
1098
            return; // Nothing to remove or drop.
1099
        }
1100
1101
        if new_nrows == 0 || ncols == 0 {
1102
            // The output matrix is empty, drop everything.
1103
            ptr::drop_in_place(data);
1104
            return;
1105
        }
1106
1107
        // Safety: because `nremove != 0`, the pointers given to `ptr::copy`
1108
        //         won’t alias.
1109
        let ptr_in = data.as_ptr();
1110
        let ptr_out = data.as_mut_ptr();
1111
1112
        let mut curr_i = i;
1113
1114
        for k in 0..ncols - 1 {
1115
            // Safety: we drop the row elements in-place because we will overwrite these
1116
            //         entries later with the `ptr::copy`.
1117
            let s = ptr::slice_from_raw_parts_mut(ptr_out.add(curr_i), nremove);
1118
            ptr::drop_in_place(s);
1119
            ptr::copy(
1120
                ptr_in.add(curr_i + (k + 1) * nremove),
1121
                ptr_out.add(curr_i),
1122
                new_nrows,
1123
            );
1124
1125
            curr_i += new_nrows;
1126
        }
1127
1128
        /*
1129
         * Deal with the last column from which less values have to be copied.
1130
         */
1131
        // Safety: we drop the row elements in-place because we will overwrite these
1132
        //         entries later with the `ptr::copy`.
1133
        let s = ptr::slice_from_raw_parts_mut(ptr_out.add(curr_i), nremove);
1134
        ptr::drop_in_place(s);
1135
        let remaining_len = nrows - i - nremove;
1136
        ptr::copy(
1137
            ptr_in.add(nrows * ncols - remaining_len),
1138
            ptr_out.add(curr_i),
1139
            remaining_len,
1140
        );
1141
    }
1142
}
1143
1144
// Moves entries of a matrix buffer to make place for `ninsert` empty rows starting at the `i-th` row index.
1145
// The `data` buffer is assumed to contained at least `(nrows + ninsert) * ncols` elements.
1146
unsafe fn extend_rows<T>(data: &mut [T], nrows: usize, ncols: usize, i: usize, ninsert: usize) {
1147
    unsafe {
1148
        let new_nrows = nrows + ninsert;
1149
1150
        if new_nrows == 0 || ncols == 0 {
1151
            return; // Nothing to do as the output matrix is empty.
1152
        }
1153
1154
        let ptr_in = data.as_ptr();
1155
        let ptr_out = data.as_mut_ptr();
1156
1157
        let remaining_len = nrows - i;
1158
        let mut curr_i = new_nrows * ncols - remaining_len;
1159
1160
        // Deal with the last column from which less values have to be copied.
1161
        ptr::copy(
1162
            ptr_in.add(nrows * ncols - remaining_len),
1163
            ptr_out.add(curr_i),
1164
            remaining_len,
1165
        );
1166
1167
        for k in (0..ncols - 1).rev() {
1168
            curr_i -= new_nrows;
1169
1170
            ptr::copy(ptr_in.add(k * nrows + i), ptr_out.add(curr_i), nrows);
1171
        }
1172
    }
1173
}
1174
1175
/// Extend the number of columns of the `Matrix` with elements from
1176
/// a given iterator.
1177
#[cfg(any(feature = "std", feature = "alloc"))]
1178
impl<T, R, S> Extend<T> for Matrix<T, R, Dyn, S>
1179
where
1180
    T: Scalar,
1181
    R: Dim,
1182
    S: Extend<T>,
1183
{
1184
    /// Extend the number of columns of the `Matrix` with elements
1185
    /// from the given iterator.
1186
    ///
1187
    /// # Example
1188
    /// ```
1189
    /// # use nalgebra::{DMatrix, Dyn, Matrix, OMatrix, Matrix3};
1190
    ///
1191
    /// let data = vec![0, 1, 2,      // column 1
1192
    ///                 3, 4, 5];     // column 2
1193
    ///
1194
    /// let mut matrix = DMatrix::from_vec(3, 2, data);
1195
    ///
1196
    /// matrix.extend(vec![6, 7, 8]); // column 3
1197
    ///
1198
    /// assert!(matrix.eq(&Matrix3::new(0, 3, 6,
1199
    ///                                 1, 4, 7,
1200
    ///                                 2, 5, 8)));
1201
    /// ```
1202
    ///
1203
    /// # Panics
1204
    /// This function panics if the number of elements yielded by the
1205
    /// given iterator is not a multiple of the number of rows of the
1206
    /// `Matrix`.
1207
    ///
1208
    /// ```should_panic
1209
    /// # use nalgebra::{DMatrix, Dyn, OMatrix};
1210
    /// let data = vec![0, 1, 2,  // column 1
1211
    ///                 3, 4, 5]; // column 2
1212
    ///
1213
    /// let mut matrix = DMatrix::from_vec(3, 2, data);
1214
    ///
1215
    /// // The following panics because the vec length is not a multiple of 3.
1216
    /// matrix.extend(vec![6, 7, 8, 9]);
1217
    /// ```
1218
    fn extend<I: IntoIterator<Item = T>>(&mut self, iter: I) {
1219
        self.data.extend(iter);
1220
    }
1221
}
1222
1223
/// Extend the number of rows of the `Vector` with elements from
1224
/// a given iterator.
1225
#[cfg(any(feature = "std", feature = "alloc"))]
1226
impl<T, S> Extend<T> for Matrix<T, Dyn, U1, S>
1227
where
1228
    T: Scalar,
1229
    S: Extend<T>,
1230
{
1231
    /// Extend the number of rows of a `Vector` with elements
1232
    /// from the given iterator.
1233
    ///
1234
    /// # Example
1235
    /// ```
1236
    /// # use nalgebra::DVector;
1237
    /// let mut vector = DVector::from_vec(vec![0, 1, 2]);
1238
    /// vector.extend(vec![3, 4, 5]);
1239
    /// assert!(vector.eq(&DVector::from_vec(vec![0, 1, 2, 3, 4, 5])));
1240
    /// ```
1241
    fn extend<I: IntoIterator<Item = T>>(&mut self, iter: I) {
1242
        self.data.extend(iter);
1243
    }
1244
}
1245
1246
#[cfg(any(feature = "std", feature = "alloc"))]
1247
impl<T, R, S, RV, SV> Extend<Vector<T, RV, SV>> for Matrix<T, R, Dyn, S>
1248
where
1249
    T: Scalar,
1250
    R: Dim,
1251
    S: Extend<Vector<T, RV, SV>>,
1252
    RV: Dim,
1253
    SV: RawStorage<T, RV>,
1254
    ShapeConstraint: SameNumberOfRows<R, RV>,
1255
{
1256
    /// Extends the number of columns of a `Matrix` with `Vector`s
1257
    /// from a given iterator.
1258
    ///
1259
    /// # Example
1260
    /// ```
1261
    /// # use nalgebra::{DMatrix, Vector3, Matrix3x4};
1262
    ///
1263
    /// let data = vec![0, 1, 2,          // column 1
1264
    ///                 3, 4, 5];         // column 2
1265
    ///
1266
    /// let mut matrix = DMatrix::from_vec(3, 2, data);
1267
    ///
1268
    /// matrix.extend(
1269
    ///   vec![Vector3::new(6,  7,  8),   // column 3
1270
    ///        Vector3::new(9, 10, 11)]); // column 4
1271
    ///
1272
    /// assert!(matrix.eq(&Matrix3x4::new(0, 3, 6,  9,
1273
    ///                                   1, 4, 7, 10,
1274
    ///                                   2, 5, 8, 11)));
1275
    /// ```
1276
    ///
1277
    /// # Panics
1278
    /// This function panics if the dimension of each `Vector` yielded
1279
    /// by the given iterator is not equal to the number of rows of
1280
    /// this `Matrix`.
1281
    ///
1282
    /// ```should_panic
1283
    /// # use nalgebra::{DMatrix, Vector2, Matrix3x4};
1284
    /// let mut matrix =
1285
    ///   DMatrix::from_vec(3, 2,
1286
    ///                     vec![0, 1, 2,   // column 1
1287
    ///                          3, 4, 5]); // column 2
1288
    ///
1289
    /// // The following panics because this matrix can only be extended with 3-dimensional vectors.
1290
    /// matrix.extend(
1291
    ///   vec![Vector2::new(6,  7)]); // too few dimensions!
1292
    /// ```
1293
    ///
1294
    /// ```should_panic
1295
    /// # use nalgebra::{DMatrix, Vector4, Matrix3x4};
1296
    /// let mut matrix =
1297
    ///   DMatrix::from_vec(3, 2,
1298
    ///                     vec![0, 1, 2,   // column 1
1299
    ///                          3, 4, 5]); // column 2
1300
    ///
1301
    /// // The following panics because this matrix can only be extended with 3-dimensional vectors.
1302
    /// matrix.extend(
1303
    ///   vec![Vector4::new(6, 7, 8, 9)]); // too few dimensions!
1304
    /// ```
1305
    fn extend<I: IntoIterator<Item = Vector<T, RV, SV>>>(&mut self, iter: I) {
1306
        self.data.extend(iter);
1307
    }
1308
}