1use super::*;
31use crate::{assert, group_helpers::VecGroup};
32use core::{cell::Cell, iter::zip, ops::Range, slice::SliceIndex};
33use dyn_stack::GlobalPodBuffer;
34use group_helpers::SliceGroup;
35
36pub use permutation::{Index, SignedIndex};
37
38mod ghost {
39 pub use crate::constrained::{group_helpers::*, permutation::*, sparse::*, *};
40}
41
42mod mem {
43 #[inline]
44 pub fn fill_zero<I: bytemuck::Zeroable>(slice: &mut [I]) {
45 let len = slice.len();
46 unsafe { core::ptr::write_bytes(slice.as_mut_ptr(), 0u8, len) }
47 }
48}
49
50#[inline(always)]
51#[track_caller]
52#[doc(hidden)]
53pub unsafe fn __get_unchecked<I, R: Clone + SliceIndex<[I]>>(slice: &[I], i: R) -> &R::Output {
54 #[cfg(debug_assertions)]
55 {
56 let _ = &slice[i.clone()];
57 }
58 unsafe { slice.get_unchecked(i) }
59}
60#[inline(always)]
61#[track_caller]
62#[doc(hidden)]
63pub unsafe fn __get_unchecked_mut<I, R: Clone + SliceIndex<[I]>>(
64 slice: &mut [I],
65 i: R,
66) -> &mut R::Output {
67 #[cfg(debug_assertions)]
68 {
69 let _ = &slice[i.clone()];
70 }
71 unsafe { slice.get_unchecked_mut(i) }
72}
73
74#[inline(always)]
75#[doc(hidden)]
76pub fn windows2<I>(slice: &[I]) -> impl DoubleEndedIterator<Item = &[I; 2]> {
77 slice
78 .windows(2)
79 .map(|window| unsafe { &*(window.as_ptr() as *const [I; 2]) })
80}
81
82#[inline]
83#[doc(hidden)]
84pub const fn repeat_byte(byte: u8) -> usize {
85 union Union {
86 bytes: [u8; 32],
87 value: usize,
88 }
89
90 let data = Union { bytes: [byte; 32] };
91 unsafe { data.value }
92}
93
94pub struct SymbolicSparseColMatRef<'a, I> {
119 nrows: usize,
120 ncols: usize,
121 col_ptr: &'a [I],
122 col_nnz: Option<&'a [I]>,
123 row_ind: &'a [I],
124}
125
126pub struct SymbolicSparseRowMatRef<'a, I> {
151 nrows: usize,
152 ncols: usize,
153 row_ptr: &'a [I],
154 row_nnz: Option<&'a [I]>,
155 col_ind: &'a [I],
156}
157
158#[derive(Clone)]
174pub struct SymbolicSparseColMat<I> {
175 nrows: usize,
176 ncols: usize,
177 col_ptr: Vec<I>,
178 col_nnz: Option<Vec<I>>,
179 row_ind: Vec<I>,
180}
181
182#[derive(Clone)]
198pub struct SymbolicSparseRowMat<I> {
199 nrows: usize,
200 ncols: usize,
201 row_ptr: Vec<I>,
202 row_nnz: Option<Vec<I>>,
203 col_ind: Vec<I>,
204}
205
206impl<I> Copy for SymbolicSparseColMatRef<'_, I> {}
207impl<I> Clone for SymbolicSparseColMatRef<'_, I> {
208 #[inline]
209 fn clone(&self) -> Self {
210 *self
211 }
212}
213impl<I> Copy for SymbolicSparseRowMatRef<'_, I> {}
214impl<I> Clone for SymbolicSparseRowMatRef<'_, I> {
215 #[inline]
216 fn clone(&self) -> Self {
217 *self
218 }
219}
220
221impl<I: Index> SymbolicSparseRowMat<I> {
222 #[inline]
228 #[track_caller]
229 pub fn new_checked(
230 nrows: usize,
231 ncols: usize,
232 row_ptrs: Vec<I>,
233 nnz_per_row: Option<Vec<I>>,
234 col_indices: Vec<I>,
235 ) -> Self {
236 SymbolicSparseRowMatRef::new_checked(
237 nrows,
238 ncols,
239 &row_ptrs,
240 nnz_per_row.as_deref(),
241 &col_indices,
242 );
243
244 Self {
245 nrows,
246 ncols,
247 row_ptr: row_ptrs,
248 row_nnz: nnz_per_row,
249 col_ind: col_indices,
250 }
251 }
252
253 #[inline]
260 #[track_caller]
261 pub fn new_unsorted_checked(
262 nrows: usize,
263 ncols: usize,
264 row_ptrs: Vec<I>,
265 nnz_per_row: Option<Vec<I>>,
266 col_indices: Vec<I>,
267 ) -> Self {
268 SymbolicSparseRowMatRef::new_unsorted_checked(
269 nrows,
270 ncols,
271 &row_ptrs,
272 nnz_per_row.as_deref(),
273 &col_indices,
274 );
275
276 Self {
277 nrows,
278 ncols,
279 row_ptr: row_ptrs,
280 row_nnz: nnz_per_row,
281 col_ind: col_indices,
282 }
283 }
284
285 #[inline(always)]
291 #[track_caller]
292 pub unsafe fn new_unchecked(
293 nrows: usize,
294 ncols: usize,
295 row_ptrs: Vec<I>,
296 nnz_per_row: Option<Vec<I>>,
297 col_indices: Vec<I>,
298 ) -> Self {
299 SymbolicSparseRowMatRef::new_unchecked(
300 nrows,
301 ncols,
302 &row_ptrs,
303 nnz_per_row.as_deref(),
304 &col_indices,
305 );
306
307 Self {
308 nrows,
309 ncols,
310 row_ptr: row_ptrs,
311 row_nnz: nnz_per_row,
312 col_ind: col_indices,
313 }
314 }
315
316 #[inline]
323 pub fn into_parts(self) -> (usize, usize, Vec<I>, Option<Vec<I>>, Vec<I>) {
324 (
325 self.nrows,
326 self.ncols,
327 self.row_ptr,
328 self.row_nnz,
329 self.col_ind,
330 )
331 }
332
333 #[inline]
335 pub fn as_ref(&self) -> SymbolicSparseRowMatRef<'_, I> {
336 SymbolicSparseRowMatRef {
337 nrows: self.nrows,
338 ncols: self.ncols,
339 row_ptr: &self.row_ptr,
340 row_nnz: self.row_nnz.as_deref(),
341 col_ind: &self.col_ind,
342 }
343 }
344
345 #[inline]
347 pub fn nrows(&self) -> usize {
348 self.nrows
349 }
350 #[inline]
352 pub fn ncols(&self) -> usize {
353 self.ncols
354 }
355
356 #[inline]
361 pub fn into_transpose(self) -> SymbolicSparseColMat<I> {
362 SymbolicSparseColMat {
363 nrows: self.ncols,
364 ncols: self.nrows,
365 col_ptr: self.row_ptr,
366 col_nnz: self.row_nnz,
367 row_ind: self.col_ind,
368 }
369 }
370
371 #[inline]
376 pub fn to_owned(&self) -> Result<SymbolicSparseRowMat<I>, FaerError> {
377 self.as_ref().to_owned()
378 }
379
380 #[inline]
385 pub fn to_col_major(&self) -> Result<SymbolicSparseColMat<I>, FaerError> {
386 self.as_ref().to_col_major()
387 }
388
389 #[inline]
397 pub fn compute_nnz(&self) -> usize {
398 self.as_ref().compute_nnz()
399 }
400
401 #[inline]
403 pub fn row_ptrs(&self) -> &[I] {
404 &self.row_ptr
405 }
406
407 #[inline]
409 pub fn nnz_per_row(&self) -> Option<&[I]> {
410 self.row_nnz.as_deref()
411 }
412
413 #[inline]
415 pub fn col_indices(&self) -> &[I] {
416 &self.col_ind
417 }
418
419 #[inline]
425 #[track_caller]
426 pub fn col_indices_of_row_raw(&self, i: usize) -> &[I] {
427 &self.col_ind[self.row_range(i)]
428 }
429
430 #[inline]
436 #[track_caller]
437 pub fn col_indices_of_row(
438 &self,
439 i: usize,
440 ) -> impl '_ + ExactSizeIterator + DoubleEndedIterator<Item = usize> {
441 self.as_ref().col_indices_of_row(i)
442 }
443
444 #[inline]
450 #[track_caller]
451 pub fn row_range(&self, i: usize) -> Range<usize> {
452 self.as_ref().row_range(i)
453 }
454
455 #[inline]
461 #[track_caller]
462 pub unsafe fn row_range_unchecked(&self, i: usize) -> Range<usize> {
463 self.as_ref().row_range_unchecked(i)
464 }
465}
466
467impl<I: Index> SymbolicSparseColMat<I> {
468 #[inline]
474 #[track_caller]
475 pub fn new_checked(
476 nrows: usize,
477 ncols: usize,
478 col_ptrs: Vec<I>,
479 nnz_per_col: Option<Vec<I>>,
480 row_indices: Vec<I>,
481 ) -> Self {
482 SymbolicSparseColMatRef::new_checked(
483 nrows,
484 ncols,
485 &col_ptrs,
486 nnz_per_col.as_deref(),
487 &row_indices,
488 );
489
490 Self {
491 nrows,
492 ncols,
493 col_ptr: col_ptrs,
494 col_nnz: nnz_per_col,
495 row_ind: row_indices,
496 }
497 }
498
499 #[inline]
506 #[track_caller]
507 pub fn new_unsorted_checked(
508 nrows: usize,
509 ncols: usize,
510 col_ptrs: Vec<I>,
511 nnz_per_col: Option<Vec<I>>,
512 row_indices: Vec<I>,
513 ) -> Self {
514 SymbolicSparseColMatRef::new_unsorted_checked(
515 nrows,
516 ncols,
517 &col_ptrs,
518 nnz_per_col.as_deref(),
519 &row_indices,
520 );
521
522 Self {
523 nrows,
524 ncols,
525 col_ptr: col_ptrs,
526 col_nnz: nnz_per_col,
527 row_ind: row_indices,
528 }
529 }
530
531 #[inline(always)]
537 #[track_caller]
538 pub unsafe fn new_unchecked(
539 nrows: usize,
540 ncols: usize,
541 col_ptrs: Vec<I>,
542 nnz_per_col: Option<Vec<I>>,
543 row_indices: Vec<I>,
544 ) -> Self {
545 SymbolicSparseRowMatRef::new_unchecked(
546 nrows,
547 ncols,
548 &col_ptrs,
549 nnz_per_col.as_deref(),
550 &row_indices,
551 );
552
553 Self {
554 nrows,
555 ncols,
556 col_ptr: col_ptrs,
557 col_nnz: nnz_per_col,
558 row_ind: row_indices,
559 }
560 }
561
562 #[inline]
569 pub fn into_parts(self) -> (usize, usize, Vec<I>, Option<Vec<I>>, Vec<I>) {
570 (
571 self.nrows,
572 self.ncols,
573 self.col_ptr,
574 self.col_nnz,
575 self.row_ind,
576 )
577 }
578
579 #[inline]
581 pub fn as_ref(&self) -> SymbolicSparseColMatRef<'_, I> {
582 SymbolicSparseColMatRef {
583 nrows: self.nrows,
584 ncols: self.ncols,
585 col_ptr: &self.col_ptr,
586 col_nnz: self.col_nnz.as_deref(),
587 row_ind: &self.row_ind,
588 }
589 }
590
591 #[inline]
593 pub fn nrows(&self) -> usize {
594 self.nrows
595 }
596 #[inline]
598 pub fn ncols(&self) -> usize {
599 self.ncols
600 }
601
602 #[inline]
607 pub fn into_transpose(self) -> SymbolicSparseRowMat<I> {
608 SymbolicSparseRowMat {
609 nrows: self.ncols,
610 ncols: self.nrows,
611 row_ptr: self.col_ptr,
612 row_nnz: self.col_nnz,
613 col_ind: self.row_ind,
614 }
615 }
616
617 #[inline]
622 pub fn to_owned(&self) -> Result<SymbolicSparseColMat<I>, FaerError> {
623 self.as_ref().to_owned()
624 }
625
626 #[inline]
631 pub fn to_row_major(&self) -> Result<SymbolicSparseRowMat<I>, FaerError> {
632 self.as_ref().to_row_major()
633 }
634
635 #[inline]
643 pub fn compute_nnz(&self) -> usize {
644 self.as_ref().compute_nnz()
645 }
646
647 #[inline]
649 pub fn col_ptrs(&self) -> &[I] {
650 &self.col_ptr
651 }
652
653 #[inline]
655 pub fn nnz_per_col(&self) -> Option<&[I]> {
656 self.col_nnz.as_deref()
657 }
658
659 #[inline]
661 pub fn row_indices(&self) -> &[I] {
662 &self.row_ind
663 }
664
665 #[inline]
671 #[track_caller]
672 pub fn row_indices_of_col_raw(&self, j: usize) -> &[I] {
673 &self.row_ind[self.col_range(j)]
674 }
675
676 #[inline]
682 #[track_caller]
683 pub fn row_indices_of_col(
684 &self,
685 j: usize,
686 ) -> impl '_ + ExactSizeIterator + DoubleEndedIterator<Item = usize> {
687 self.as_ref().row_indices_of_col(j)
688 }
689
690 #[inline]
696 #[track_caller]
697 pub fn col_range(&self, j: usize) -> Range<usize> {
698 self.as_ref().col_range(j)
699 }
700
701 #[inline]
707 #[track_caller]
708 pub unsafe fn col_range_unchecked(&self, j: usize) -> Range<usize> {
709 self.as_ref().col_range_unchecked(j)
710 }
711}
712
713impl<'a, I: Index> SymbolicSparseRowMatRef<'a, I> {
714 #[inline]
720 #[track_caller]
721 pub fn new_checked(
722 nrows: usize,
723 ncols: usize,
724 row_ptrs: &'a [I],
725 nnz_per_row: Option<&'a [I]>,
726 col_indices: &'a [I],
727 ) -> Self {
728 assert!(all(
729 ncols <= I::Signed::MAX.zx(),
730 nrows <= I::Signed::MAX.zx(),
731 ));
732 assert!(row_ptrs.len() == nrows + 1);
733 for &[c, c_next] in windows2(row_ptrs) {
734 assert!(c <= c_next);
735 }
736 assert!(row_ptrs[ncols].zx() <= col_indices.len());
737
738 if let Some(nnz_per_row) = nnz_per_row {
739 for (&nnz_i, &[c, c_next]) in zip(nnz_per_row, windows2(row_ptrs)) {
740 assert!(nnz_i <= c_next - c);
741 let col_indices = &col_indices[c.zx()..c.zx() + nnz_i.zx()];
742 if !col_indices.is_empty() {
743 let mut j_prev = col_indices[0];
744 for &j in &col_indices[1..] {
745 assert!(j_prev < j);
746 j_prev = j;
747 }
748 let ncols = I::truncate(ncols);
749 assert!(j_prev < ncols);
750 }
751 }
752 } else {
753 for &[c, c_next] in windows2(row_ptrs) {
754 let col_indices = &col_indices[c.zx()..c_next.zx()];
755 if !col_indices.is_empty() {
756 let mut j_prev = col_indices[0];
757 for &j in &col_indices[1..] {
758 assert!(j_prev < j);
759 j_prev = j;
760 }
761 let ncols = I::truncate(ncols);
762 assert!(j_prev < ncols);
763 }
764 }
765 }
766
767 Self {
768 nrows,
769 ncols,
770 row_ptr: row_ptrs,
771 row_nnz: nnz_per_row,
772 col_ind: col_indices,
773 }
774 }
775
776 #[inline]
783 #[track_caller]
784 pub fn new_unsorted_checked(
785 nrows: usize,
786 ncols: usize,
787 row_ptrs: &'a [I],
788 nnz_per_row: Option<&'a [I]>,
789 col_indices: &'a [I],
790 ) -> Self {
791 assert!(all(
792 ncols <= I::Signed::MAX.zx(),
793 nrows <= I::Signed::MAX.zx(),
794 ));
795 assert!(row_ptrs.len() == nrows + 1);
796 for &[c, c_next] in windows2(row_ptrs) {
797 assert!(c <= c_next);
798 }
799 assert!(row_ptrs[ncols].zx() <= col_indices.len());
800
801 if let Some(nnz_per_row) = nnz_per_row {
802 for (&nnz_i, &[c, c_next]) in zip(nnz_per_row, windows2(row_ptrs)) {
803 assert!(nnz_i <= c_next - c);
804 for &j in &col_indices[c.zx()..c.zx() + nnz_i.zx()] {
805 assert!(j < I::truncate(ncols));
806 }
807 }
808 } else {
809 let c0 = row_ptrs[0].zx();
810 let cn = row_ptrs[ncols].zx();
811 for &j in &col_indices[c0..cn] {
812 assert!(j < I::truncate(ncols));
813 }
814 }
815
816 Self {
817 nrows,
818 ncols,
819 row_ptr: row_ptrs,
820 row_nnz: nnz_per_row,
821 col_ind: col_indices,
822 }
823 }
824
825 #[inline(always)]
831 #[track_caller]
832 pub unsafe fn new_unchecked(
833 nrows: usize,
834 ncols: usize,
835 row_ptrs: &'a [I],
836 nnz_per_row: Option<&'a [I]>,
837 col_indices: &'a [I],
838 ) -> Self {
839 assert!(all(
840 ncols <= <I::Signed as SignedIndex>::MAX.zx(),
841 nrows <= <I::Signed as SignedIndex>::MAX.zx(),
842 ));
843 assert!(row_ptrs.len() == nrows + 1);
844 assert!(row_ptrs[nrows].zx() <= col_indices.len());
845
846 Self {
847 nrows,
848 ncols,
849 row_ptr: row_ptrs,
850 row_nnz: nnz_per_row,
851 col_ind: col_indices,
852 }
853 }
854
855 #[inline]
857 pub fn nrows(&self) -> usize {
858 self.nrows
859 }
860 #[inline]
862 pub fn ncols(&self) -> usize {
863 self.ncols
864 }
865
866 #[inline]
868 pub fn transpose(self) -> SymbolicSparseColMatRef<'a, I> {
869 SymbolicSparseColMatRef {
870 nrows: self.ncols,
871 ncols: self.nrows,
872 col_ptr: self.row_ptr,
873 col_nnz: self.row_nnz,
874 row_ind: self.col_ind,
875 }
876 }
877
878 #[inline]
883 pub fn to_owned(&self) -> Result<SymbolicSparseRowMat<I>, FaerError> {
884 self.transpose()
885 .to_owned()
886 .map(SymbolicSparseColMat::into_transpose)
887 }
888
889 #[inline]
894 pub fn to_col_major(&self) -> Result<SymbolicSparseColMat<I>, FaerError> {
895 self.transpose().to_row_major().map(|m| m.into_transpose())
896 }
897
898 #[inline]
906 pub fn compute_nnz(&self) -> usize {
907 self.transpose().compute_nnz()
908 }
909
910 #[inline]
912 pub fn row_ptrs(&self) -> &'a [I] {
913 self.row_ptr
914 }
915
916 #[inline]
918 pub fn nnz_per_row(&self) -> Option<&'a [I]> {
919 self.row_nnz
920 }
921
922 #[inline]
924 pub fn col_indices(&self) -> &'a [I] {
925 self.col_ind
926 }
927
928 #[inline]
934 #[track_caller]
935 pub fn col_indices_of_row_raw(&self, i: usize) -> &'a [I] {
936 &self.col_ind[self.row_range(i)]
937 }
938
939 #[inline]
945 #[track_caller]
946 pub fn col_indices_of_row(
947 &self,
948 i: usize,
949 ) -> impl 'a + ExactSizeIterator + DoubleEndedIterator<Item = usize> {
950 self.col_indices_of_row_raw(i).iter().map(
951 #[inline(always)]
952 |&i| i.zx(),
953 )
954 }
955
956 #[inline]
962 #[track_caller]
963 pub fn row_range(&self, i: usize) -> Range<usize> {
964 let start = self.row_ptr[i].zx();
965 let end = self
966 .row_nnz
967 .map(|row_nnz| row_nnz[i].zx() + start)
968 .unwrap_or(self.row_ptr[i + 1].zx());
969
970 start..end
971 }
972
973 #[inline]
979 #[track_caller]
980 pub unsafe fn row_range_unchecked(&self, i: usize) -> Range<usize> {
981 let start = __get_unchecked(self.row_ptr, i).zx();
982 let end = self
983 .row_nnz
984 .map(|row_nnz| (__get_unchecked(row_nnz, i).zx() + start))
985 .unwrap_or(__get_unchecked(self.row_ptr, i + 1).zx());
986
987 start..end
988 }
989}
990
991impl<'a, I: Index> SymbolicSparseColMatRef<'a, I> {
992 #[inline]
998 #[track_caller]
999 pub fn new_checked(
1000 nrows: usize,
1001 ncols: usize,
1002 col_ptrs: &'a [I],
1003 nnz_per_col: Option<&'a [I]>,
1004 row_indices: &'a [I],
1005 ) -> Self {
1006 assert!(all(
1007 ncols <= I::Signed::MAX.zx(),
1008 nrows <= I::Signed::MAX.zx(),
1009 ));
1010 assert!(col_ptrs.len() == ncols + 1);
1011 for &[c, c_next] in windows2(col_ptrs) {
1012 assert!(c <= c_next);
1013 }
1014 assert!(col_ptrs[ncols].zx() <= row_indices.len());
1015
1016 if let Some(nnz_per_col) = nnz_per_col {
1017 for (&nnz_j, &[c, c_next]) in zip(nnz_per_col, windows2(col_ptrs)) {
1018 assert!(nnz_j <= c_next - c);
1019 let row_indices = &row_indices[c.zx()..c.zx() + nnz_j.zx()];
1020 if !row_indices.is_empty() {
1021 let mut i_prev = row_indices[0];
1022 for &i in &row_indices[1..] {
1023 assert!(i_prev < i);
1024 i_prev = i;
1025 }
1026 let nrows = I::truncate(nrows);
1027 assert!(i_prev < nrows);
1028 }
1029 }
1030 } else {
1031 for &[c, c_next] in windows2(col_ptrs) {
1032 let row_indices = &row_indices[c.zx()..c_next.zx()];
1033 if !row_indices.is_empty() {
1034 let mut i_prev = row_indices[0];
1035 for &i in &row_indices[1..] {
1036 assert!(i_prev < i);
1037 i_prev = i;
1038 }
1039 let nrows = I::truncate(nrows);
1040 assert!(i_prev < nrows);
1041 }
1042 }
1043 }
1044
1045 Self {
1046 nrows,
1047 ncols,
1048 col_ptr: col_ptrs,
1049 col_nnz: nnz_per_col,
1050 row_ind: row_indices,
1051 }
1052 }
1053
1054 #[inline]
1061 #[track_caller]
1062 pub fn new_unsorted_checked(
1063 nrows: usize,
1064 ncols: usize,
1065 col_ptrs: &'a [I],
1066 nnz_per_col: Option<&'a [I]>,
1067 row_indices: &'a [I],
1068 ) -> Self {
1069 assert!(all(
1070 ncols <= I::Signed::MAX.zx(),
1071 nrows <= I::Signed::MAX.zx(),
1072 ));
1073 assert!(col_ptrs.len() == ncols + 1);
1074 for &[c, c_next] in windows2(col_ptrs) {
1075 assert!(c <= c_next);
1076 }
1077 assert!(col_ptrs[ncols].zx() <= row_indices.len());
1078
1079 if let Some(nnz_per_col) = nnz_per_col {
1080 for (&nnz_j, &[c, c_next]) in zip(nnz_per_col, windows2(col_ptrs)) {
1081 assert!(nnz_j <= c_next - c);
1082 for &i in &row_indices[c.zx()..c.zx() + nnz_j.zx()] {
1083 assert!(i < I::truncate(nrows));
1084 }
1085 }
1086 } else {
1087 let c0 = col_ptrs[0].zx();
1088 let cn = col_ptrs[ncols].zx();
1089 for &i in &row_indices[c0..cn] {
1090 assert!(i < I::truncate(nrows));
1091 }
1092 }
1093
1094 Self {
1095 nrows,
1096 ncols,
1097 col_ptr: col_ptrs,
1098 col_nnz: nnz_per_col,
1099 row_ind: row_indices,
1100 }
1101 }
1102
1103 #[inline(always)]
1109 #[track_caller]
1110 pub unsafe fn new_unchecked(
1111 nrows: usize,
1112 ncols: usize,
1113 col_ptrs: &'a [I],
1114 nnz_per_col: Option<&'a [I]>,
1115 row_indices: &'a [I],
1116 ) -> Self {
1117 assert!(all(
1118 ncols <= <I::Signed as SignedIndex>::MAX.zx(),
1119 nrows <= <I::Signed as SignedIndex>::MAX.zx(),
1120 ));
1121 assert!(col_ptrs.len() == ncols + 1);
1122 assert!(col_ptrs[ncols].zx() <= row_indices.len());
1123
1124 Self {
1125 nrows,
1126 ncols,
1127 col_ptr: col_ptrs,
1128 col_nnz: nnz_per_col,
1129 row_ind: row_indices,
1130 }
1131 }
1132
1133 #[inline]
1135 pub fn nrows(&self) -> usize {
1136 self.nrows
1137 }
1138 #[inline]
1140 pub fn ncols(&self) -> usize {
1141 self.ncols
1142 }
1143
1144 #[inline]
1146 pub fn transpose(self) -> SymbolicSparseRowMatRef<'a, I> {
1147 SymbolicSparseRowMatRef {
1148 nrows: self.ncols,
1149 ncols: self.nrows,
1150 row_ptr: self.col_ptr,
1151 row_nnz: self.col_nnz,
1152 col_ind: self.row_ind,
1153 }
1154 }
1155
1156 #[inline]
1161 pub fn to_owned(&self) -> Result<SymbolicSparseColMat<I>, FaerError> {
1162 Ok(SymbolicSparseColMat {
1163 nrows: self.nrows,
1164 ncols: self.ncols,
1165 col_ptr: try_collect(self.col_ptr.iter().copied())?,
1166 col_nnz: self
1167 .col_nnz
1168 .map(|nnz| try_collect(nnz.iter().copied()))
1169 .transpose()?,
1170 row_ind: try_collect(self.row_ind.iter().copied())?,
1171 })
1172 }
1173
1174 #[inline]
1179 pub fn to_row_major(&self) -> Result<SymbolicSparseRowMat<I>, FaerError> {
1180 let mut col_ptr = try_zeroed::<I>(self.nrows + 1)?;
1181 let mut row_ind = try_zeroed::<I>(self.compute_nnz())?;
1182
1183 let mut mem = GlobalPodBuffer::try_new(StackReq::new::<I>(self.nrows))
1184 .map_err(|_| FaerError::OutOfMemory)?;
1185
1186 util::adjoint_symbolic(&mut col_ptr, &mut row_ind, *self, PodStack::new(&mut mem));
1187
1188 let transpose = unsafe {
1189 SymbolicSparseColMat::new_unchecked(self.ncols, self.nrows, col_ptr, None, row_ind)
1190 };
1191
1192 Ok(transpose.into_transpose())
1193 }
1194
1195 #[inline]
1203 pub fn compute_nnz(&self) -> usize {
1204 match self.col_nnz {
1205 Some(col_nnz) => {
1206 let mut nnz = 0usize;
1207 for &nnz_j in col_nnz {
1208 nnz += nnz_j.zx();
1210 }
1211 nnz
1212 }
1213 None => self.col_ptr[self.ncols].zx() - self.col_ptr[0].zx(),
1214 }
1215 }
1216
1217 #[inline]
1219 pub fn col_ptrs(&self) -> &'a [I] {
1220 self.col_ptr
1221 }
1222
1223 #[inline]
1225 pub fn nnz_per_col(&self) -> Option<&'a [I]> {
1226 self.col_nnz
1227 }
1228
1229 #[inline]
1231 pub fn row_indices(&self) -> &'a [I] {
1232 self.row_ind
1233 }
1234
1235 #[inline]
1241 #[track_caller]
1242 pub fn row_indices_of_col_raw(&self, j: usize) -> &'a [I] {
1243 &self.row_ind[self.col_range(j)]
1244 }
1245
1246 #[inline]
1252 #[track_caller]
1253 pub fn row_indices_of_col(
1254 &self,
1255 j: usize,
1256 ) -> impl 'a + ExactSizeIterator + DoubleEndedIterator<Item = usize> {
1257 self.row_indices_of_col_raw(j).iter().map(
1258 #[inline(always)]
1259 |&i| i.zx(),
1260 )
1261 }
1262
1263 #[inline]
1269 #[track_caller]
1270 pub fn col_range(&self, j: usize) -> Range<usize> {
1271 let start = self.col_ptr[j].zx();
1272 let end = self
1273 .col_nnz
1274 .map(|col_nnz| col_nnz[j].zx() + start)
1275 .unwrap_or(self.col_ptr[j + 1].zx());
1276
1277 start..end
1278 }
1279
1280 #[inline]
1286 #[track_caller]
1287 pub unsafe fn col_range_unchecked(&self, j: usize) -> Range<usize> {
1288 let start = __get_unchecked(self.col_ptr, j).zx();
1289 let end = self
1290 .col_nnz
1291 .map(|col_nnz| (__get_unchecked(col_nnz, j).zx() + start))
1292 .unwrap_or(__get_unchecked(self.col_ptr, j + 1).zx());
1293
1294 start..end
1295 }
1296}
1297
1298impl<I: Index> core::fmt::Debug for SymbolicSparseColMatRef<'_, I> {
1299 fn fmt(&self, f: &mut core::fmt::Formatter<'_>) -> core::fmt::Result {
1300 let mat = *self;
1301 let mut iter = (0..mat.ncols()).into_iter().flat_map(move |j| {
1302 struct Wrapper(usize, usize);
1303 impl core::fmt::Debug for Wrapper {
1304 fn fmt(&self, f: &mut std::fmt::Formatter<'_>) -> std::fmt::Result {
1305 let row = self.0;
1306 let col = self.1;
1307 write!(f, "({row}, {col}")
1308 }
1309 }
1310
1311 mat.row_indices_of_col(j).map(move |i| Wrapper(i, j))
1312 });
1313
1314 f.debug_list().entries(&mut iter).finish()
1315 }
1316}
1317
1318impl<I: Index> core::fmt::Debug for SymbolicSparseRowMatRef<'_, I> {
1319 fn fmt(&self, f: &mut core::fmt::Formatter<'_>) -> core::fmt::Result {
1320 let mat = *self;
1321 let mut iter = (0..mat.nrows()).into_iter().flat_map(move |i| {
1322 struct Wrapper(usize, usize);
1323 impl core::fmt::Debug for Wrapper {
1324 fn fmt(&self, f: &mut std::fmt::Formatter<'_>) -> std::fmt::Result {
1325 let row = self.0;
1326 let col = self.1;
1327 write!(f, "({row}, {col}")
1328 }
1329 }
1330
1331 mat.col_indices_of_row(i).map(move |j| Wrapper(i, j))
1332 });
1333
1334 f.debug_list().entries(&mut iter).finish()
1335 }
1336}
1337impl<I: Index, E: Entity> core::fmt::Debug for SparseColMatRef<'_, I, E> {
1338 fn fmt(&self, f: &mut core::fmt::Formatter<'_>) -> core::fmt::Result {
1339 let mat = *self;
1340 let mut iter = (0..mat.ncols()).into_iter().flat_map(move |j| {
1341 struct Wrapper<E>(usize, usize, E);
1342 impl<E: core::fmt::Debug> core::fmt::Debug for Wrapper<E> {
1343 fn fmt(&self, f: &mut std::fmt::Formatter<'_>) -> std::fmt::Result {
1344 let row = self.0;
1345 let col = self.1;
1346 let val = &self.2;
1347 write!(f, "({row}, {col}, {val:?})")
1348 }
1349 }
1350
1351 mat.row_indices_of_col(j)
1352 .zip(SliceGroup::<E>::new(mat.values_of_col(j)).into_ref_iter())
1353 .map(move |(i, val)| Wrapper(i, j, val.read()))
1354 });
1355
1356 f.debug_list().entries(&mut iter).finish()
1357 }
1358}
1359
1360impl<I: Index, E: Entity> core::fmt::Debug for SparseRowMatRef<'_, I, E> {
1361 fn fmt(&self, f: &mut core::fmt::Formatter<'_>) -> core::fmt::Result {
1362 let mat = *self;
1363 let mut iter = (0..mat.nrows()).into_iter().flat_map(move |i| {
1364 struct Wrapper<E>(usize, usize, E);
1365 impl<E: core::fmt::Debug> core::fmt::Debug for Wrapper<E> {
1366 fn fmt(&self, f: &mut std::fmt::Formatter<'_>) -> std::fmt::Result {
1367 let row = self.0;
1368 let col = self.1;
1369 let val = &self.2;
1370 write!(f, "({row}, {col}, {val:?})")
1371 }
1372 }
1373
1374 mat.col_indices_of_row(i)
1375 .zip(SliceGroup::<E>::new(mat.values_of_row(i)).into_ref_iter())
1376 .map(move |(j, val)| Wrapper(i, j, val.read()))
1377 });
1378
1379 f.debug_list().entries(&mut iter).finish()
1380 }
1381}
1382
1383impl<I: Index, E: Entity> core::fmt::Debug for SparseColMatMut<'_, I, E> {
1384 fn fmt(&self, f: &mut core::fmt::Formatter<'_>) -> core::fmt::Result {
1385 self.rb().fmt(f)
1386 }
1387}
1388
1389impl<I: Index, E: Entity> core::fmt::Debug for SparseColMat<I, E> {
1390 fn fmt(&self, f: &mut core::fmt::Formatter<'_>) -> core::fmt::Result {
1391 self.as_ref().fmt(f)
1392 }
1393}
1394
1395impl<I: Index> core::fmt::Debug for SymbolicSparseColMat<I> {
1396 fn fmt(&self, f: &mut core::fmt::Formatter<'_>) -> core::fmt::Result {
1397 self.as_ref().fmt(f)
1398 }
1399}
1400
1401impl<I: Index, E: Entity> core::fmt::Debug for SparseRowMatMut<'_, I, E> {
1402 fn fmt(&self, f: &mut core::fmt::Formatter<'_>) -> core::fmt::Result {
1403 self.rb().fmt(f)
1404 }
1405}
1406
1407impl<I: Index, E: Entity> core::fmt::Debug for SparseRowMat<I, E> {
1408 fn fmt(&self, f: &mut core::fmt::Formatter<'_>) -> core::fmt::Result {
1409 self.as_ref().fmt(f)
1410 }
1411}
1412
1413impl<I: Index> core::fmt::Debug for SymbolicSparseRowMat<I> {
1414 fn fmt(&self, f: &mut core::fmt::Formatter<'_>) -> core::fmt::Result {
1415 self.as_ref().fmt(f)
1416 }
1417}
1418
1419pub type SparseRowMatMut<'a, I, E> = Matrix<inner::SparseRowMatMut<'a, I, E>>;
1424
1425pub type SparseColMatMut<'a, I, E> = Matrix<inner::SparseColMatMut<'a, I, E>>;
1430
1431pub type SparseRowMatRef<'a, I, E> = Matrix<inner::SparseRowMatRef<'a, I, E>>;
1433
1434pub type SparseColMatRef<'a, I, E> = Matrix<inner::SparseColMatRef<'a, I, E>>;
1436
1437pub type SparseRowMat<I, E> = Matrix<inner::SparseRowMat<I, E>>;
1439
1440pub type SparseColMat<I, E> = Matrix<inner::SparseColMat<I, E>>;
1442
1443impl<'a, I: Index, E: Entity> SparseRowMatMut<'a, I, E> {
1444 #[inline]
1451 #[track_caller]
1452 pub fn new(
1453 symbolic: SymbolicSparseRowMatRef<'a, I>,
1454 values: GroupFor<E, &'a mut [E::Unit]>,
1455 ) -> Self {
1456 let values = SliceGroupMut::new(values);
1457 assert!(symbolic.col_indices().len() == values.len());
1458 Self {
1459 inner: inner::SparseRowMatMut { symbolic, values },
1460 }
1461 }
1462
1463 #[inline]
1468 pub fn to_owned(&self) -> Result<SparseRowMat<I, E::Canonical>, FaerError>
1469 where
1470 E: Conjugate,
1471 E::Canonical: ComplexField,
1472 {
1473 self.rb().to_owned()
1474 }
1475
1476 #[inline]
1481 pub fn to_col_major(&self) -> Result<SparseColMat<I, E::Canonical>, FaerError>
1482 where
1483 E: Conjugate,
1484 E::Canonical: ComplexField,
1485 {
1486 self.rb().to_col_major()
1487 }
1488
1489 #[inline]
1491 pub fn transpose_mut(self) -> SparseColMatMut<'a, I, E> {
1492 SparseColMatMut {
1493 inner: inner::SparseColMatMut {
1494 symbolic: SymbolicSparseColMatRef {
1495 nrows: self.inner.symbolic.ncols,
1496 ncols: self.inner.symbolic.nrows,
1497 col_ptr: self.inner.symbolic.row_ptr,
1498 col_nnz: self.inner.symbolic.row_nnz,
1499 row_ind: self.inner.symbolic.col_ind,
1500 },
1501 values: self.inner.values,
1502 },
1503 }
1504 }
1505
1506 #[inline]
1508 pub fn canonicalize_mut(self) -> (SparseRowMatMut<'a, I, E::Canonical>, Conj)
1509 where
1510 E: Conjugate,
1511 {
1512 (
1513 SparseRowMatMut {
1514 inner: inner::SparseRowMatMut {
1515 symbolic: self.inner.symbolic,
1516 values: unsafe {
1517 SliceGroupMut::<'a, E::Canonical>::new(transmute_unchecked::<
1518 GroupFor<E, &mut [UnitFor<E::Canonical>]>,
1519 GroupFor<E::Canonical, &mut [UnitFor<E::Canonical>]>,
1520 >(
1521 E::faer_map(self.inner.values.into_inner(), |slice| {
1522 let len = slice.len();
1523 core::slice::from_raw_parts_mut(
1524 slice.as_mut_ptr() as *mut UnitFor<E>
1525 as *mut UnitFor<E::Canonical>,
1526 len,
1527 )
1528 }),
1529 ))
1530 },
1531 },
1532 },
1533 if coe::is_same::<E, E::Canonical>() {
1534 Conj::No
1535 } else {
1536 Conj::Yes
1537 },
1538 )
1539 }
1540
1541 #[inline]
1543 pub fn conjugate_mut(self) -> SparseRowMatMut<'a, I, E::Conj>
1544 where
1545 E: Conjugate,
1546 {
1547 SparseRowMatMut {
1548 inner: inner::SparseRowMatMut {
1549 symbolic: self.inner.symbolic,
1550 values: unsafe {
1551 SliceGroupMut::<'a, E::Conj>::new(transmute_unchecked::<
1552 GroupFor<E, &mut [UnitFor<E::Conj>]>,
1553 GroupFor<E::Conj, &mut [UnitFor<E::Conj>]>,
1554 >(E::faer_map(
1555 self.inner.values.into_inner(),
1556 |slice| {
1557 let len = slice.len();
1558 core::slice::from_raw_parts_mut(
1559 slice.as_mut_ptr() as *mut UnitFor<E> as *mut UnitFor<E::Conj>,
1560 len,
1561 )
1562 },
1563 )))
1564 },
1565 },
1566 }
1567 }
1568
1569 #[inline]
1571 pub fn adjoint_mut(self) -> SparseColMatMut<'a, I, E::Conj>
1572 where
1573 E: Conjugate,
1574 {
1575 self.transpose_mut().conjugate_mut()
1576 }
1577
1578 #[inline]
1580 pub fn values_mut(self) -> GroupFor<E, &'a mut [E::Unit]> {
1581 self.inner.values.into_inner()
1582 }
1583
1584 #[inline]
1590 #[track_caller]
1591 pub fn values_of_row_mut(self, i: usize) -> GroupFor<E, &'a mut [E::Unit]> {
1592 let range = self.symbolic().row_range(i);
1593 self.inner.values.subslice(range).into_inner()
1594 }
1595
1596 #[inline]
1598 pub fn symbolic(&self) -> SymbolicSparseRowMatRef<'a, I> {
1599 self.inner.symbolic
1600 }
1601
1602 #[inline]
1604 pub fn into_parts(
1605 self,
1606 ) -> (
1607 SymbolicSparseRowMatRef<'a, I>,
1608 GroupFor<E, &'a mut [E::Unit]>,
1609 ) {
1610 (self.inner.symbolic, self.inner.values.into_inner())
1611 }
1612}
1613
1614impl<'a, I: Index, E: Entity> SparseColMatMut<'a, I, E> {
1615 #[inline]
1622 #[track_caller]
1623 pub fn new(
1624 symbolic: SymbolicSparseColMatRef<'a, I>,
1625 values: GroupFor<E, &'a mut [E::Unit]>,
1626 ) -> Self {
1627 let values = SliceGroupMut::new(values);
1628 assert!(symbolic.row_indices().len() == values.len());
1629 Self {
1630 inner: inner::SparseColMatMut { symbolic, values },
1631 }
1632 }
1633
1634 #[inline]
1639 pub fn to_owned(&self) -> Result<SparseColMat<I, E::Canonical>, FaerError>
1640 where
1641 E: Conjugate,
1642 E::Canonical: ComplexField,
1643 {
1644 self.rb().to_owned()
1645 }
1646
1647 #[inline]
1652 pub fn to_row_major(&self) -> Result<SparseRowMat<I, E::Canonical>, FaerError>
1653 where
1654 E: Conjugate,
1655 E::Canonical: ComplexField,
1656 {
1657 self.rb().to_row_major()
1658 }
1659
1660 #[inline]
1662 pub fn transpose_mut(self) -> SparseRowMatMut<'a, I, E> {
1663 SparseRowMatMut {
1664 inner: inner::SparseRowMatMut {
1665 symbolic: SymbolicSparseRowMatRef {
1666 nrows: self.inner.symbolic.ncols,
1667 ncols: self.inner.symbolic.nrows,
1668 row_ptr: self.inner.symbolic.col_ptr,
1669 row_nnz: self.inner.symbolic.col_nnz,
1670 col_ind: self.inner.symbolic.row_ind,
1671 },
1672 values: self.inner.values,
1673 },
1674 }
1675 }
1676
1677 #[inline]
1679 pub fn conjugate_mut(self) -> SparseColMatMut<'a, I, E::Conj>
1680 where
1681 E: Conjugate,
1682 {
1683 SparseColMatMut {
1684 inner: inner::SparseColMatMut {
1685 symbolic: self.inner.symbolic,
1686 values: unsafe {
1687 SliceGroupMut::<'a, E::Conj>::new(transmute_unchecked::<
1688 GroupFor<E, &mut [UnitFor<E::Conj>]>,
1689 GroupFor<E::Conj, &mut [UnitFor<E::Conj>]>,
1690 >(E::faer_map(
1691 self.inner.values.into_inner(),
1692 |slice| {
1693 let len = slice.len();
1694 core::slice::from_raw_parts_mut(
1695 slice.as_ptr() as *mut UnitFor<E> as *mut UnitFor<E::Conj>,
1696 len,
1697 )
1698 },
1699 )))
1700 },
1701 },
1702 }
1703 }
1704
1705 #[inline]
1707 pub fn canonicalize_mut(self) -> (SparseColMatMut<'a, I, E::Canonical>, Conj)
1708 where
1709 E: Conjugate,
1710 {
1711 (
1712 SparseColMatMut {
1713 inner: inner::SparseColMatMut {
1714 symbolic: self.inner.symbolic,
1715 values: unsafe {
1716 SliceGroupMut::<'a, E::Canonical>::new(transmute_unchecked::<
1717 GroupFor<E, &mut [UnitFor<E::Canonical>]>,
1718 GroupFor<E::Canonical, &mut [UnitFor<E::Canonical>]>,
1719 >(
1720 E::faer_map(self.inner.values.into_inner(), |slice| {
1721 let len = slice.len();
1722 core::slice::from_raw_parts_mut(
1723 slice.as_mut_ptr() as *mut UnitFor<E>
1724 as *mut UnitFor<E::Canonical>,
1725 len,
1726 )
1727 }),
1728 ))
1729 },
1730 },
1731 },
1732 if coe::is_same::<E, E::Canonical>() {
1733 Conj::No
1734 } else {
1735 Conj::Yes
1736 },
1737 )
1738 }
1739
1740 #[inline]
1742 pub fn adjoint_mut(self) -> SparseRowMatMut<'a, I, E::Conj>
1743 where
1744 E: Conjugate,
1745 {
1746 self.transpose_mut().conjugate_mut()
1747 }
1748
1749 #[inline]
1751 pub fn values_mut(self) -> GroupFor<E, &'a mut [E::Unit]> {
1752 self.inner.values.into_inner()
1753 }
1754
1755 #[inline]
1761 #[track_caller]
1762 pub fn values_of_col_mut(self, j: usize) -> GroupFor<E, &'a mut [E::Unit]> {
1763 let range = self.symbolic().col_range(j);
1764 self.inner.values.subslice(range).into_inner()
1765 }
1766
1767 #[inline]
1769 pub fn symbolic(&self) -> SymbolicSparseColMatRef<'a, I> {
1770 self.inner.symbolic
1771 }
1772
1773 #[inline]
1775 pub fn into_parts_mut(
1776 self,
1777 ) -> (
1778 SymbolicSparseColMatRef<'a, I>,
1779 GroupFor<E, &'a mut [E::Unit]>,
1780 ) {
1781 (self.inner.symbolic, self.inner.values.into_inner())
1782 }
1783}
1784
1785impl<'a, I: Index, E: Entity> SparseRowMatRef<'a, I, E> {
1786 #[inline]
1793 #[track_caller]
1794 pub fn new(
1795 symbolic: SymbolicSparseRowMatRef<'a, I>,
1796 values: GroupFor<E, &'a [E::Unit]>,
1797 ) -> Self {
1798 let values = SliceGroup::new(values);
1799 assert!(symbolic.col_indices().len() == values.len());
1800 Self {
1801 inner: inner::SparseRowMatRef { symbolic, values },
1802 }
1803 }
1804
1805 #[inline]
1807 pub fn values(self) -> GroupFor<E, &'a [E::Unit]> {
1808 self.inner.values.into_inner()
1809 }
1810
1811 #[inline]
1816 pub fn to_owned(&self) -> Result<SparseRowMat<I, E::Canonical>, FaerError>
1817 where
1818 E: Conjugate,
1819 E::Canonical: ComplexField,
1820 {
1821 self.transpose()
1822 .to_owned()
1823 .map(SparseColMat::into_transpose)
1824 }
1825
1826 #[inline]
1831 pub fn to_col_major(&self) -> Result<SparseColMat<I, E::Canonical>, FaerError>
1832 where
1833 E: Conjugate,
1834 E::Canonical: ComplexField,
1835 {
1836 self.transpose()
1837 .to_row_major()
1838 .map(SparseRowMat::into_transpose)
1839 }
1840
1841 #[inline]
1843 pub fn transpose(self) -> SparseColMatRef<'a, I, E> {
1844 SparseColMatRef {
1845 inner: inner::SparseColMatRef {
1846 symbolic: SymbolicSparseColMatRef {
1847 nrows: self.inner.symbolic.ncols,
1848 ncols: self.inner.symbolic.nrows,
1849 col_ptr: self.inner.symbolic.row_ptr,
1850 col_nnz: self.inner.symbolic.row_nnz,
1851 row_ind: self.inner.symbolic.col_ind,
1852 },
1853 values: self.inner.values,
1854 },
1855 }
1856 }
1857
1858 #[inline]
1860 pub fn conjugate(self) -> SparseRowMatRef<'a, I, E::Conj>
1861 where
1862 E: Conjugate,
1863 {
1864 SparseRowMatRef {
1865 inner: inner::SparseRowMatRef {
1866 symbolic: self.inner.symbolic,
1867 values: unsafe {
1868 SliceGroup::<'a, E::Conj>::new(transmute_unchecked::<
1869 GroupFor<E, &[UnitFor<E::Conj>]>,
1870 GroupFor<E::Conj, &[UnitFor<E::Conj>]>,
1871 >(E::faer_map(
1872 self.inner.values.into_inner(),
1873 |slice| {
1874 let len = slice.len();
1875 core::slice::from_raw_parts(
1876 slice.as_ptr() as *const UnitFor<E> as *const UnitFor<E::Conj>,
1877 len,
1878 )
1879 },
1880 )))
1881 },
1882 },
1883 }
1884 }
1885
1886 #[inline]
1888 pub fn canonicalize(self) -> (SparseRowMatRef<'a, I, E::Canonical>, Conj)
1889 where
1890 E: Conjugate,
1891 {
1892 (
1893 SparseRowMatRef {
1894 inner: inner::SparseRowMatRef {
1895 symbolic: self.inner.symbolic,
1896 values: unsafe {
1897 SliceGroup::<'a, E::Canonical>::new(transmute_unchecked::<
1898 GroupFor<E, &[UnitFor<E::Canonical>]>,
1899 GroupFor<E::Canonical, &[UnitFor<E::Canonical>]>,
1900 >(E::faer_map(
1901 self.inner.values.into_inner(),
1902 |slice| {
1903 let len = slice.len();
1904 core::slice::from_raw_parts(
1905 slice.as_ptr() as *const UnitFor<E>
1906 as *const UnitFor<E::Canonical>,
1907 len,
1908 )
1909 },
1910 )))
1911 },
1912 },
1913 },
1914 if coe::is_same::<E, E::Canonical>() {
1915 Conj::No
1916 } else {
1917 Conj::Yes
1918 },
1919 )
1920 }
1921
1922 #[inline]
1924 pub fn adjoint(self) -> SparseColMatRef<'a, I, E::Conj>
1925 where
1926 E: Conjugate,
1927 {
1928 self.transpose().conjugate()
1929 }
1930
1931 #[inline]
1937 #[track_caller]
1938 pub fn values_of_row(self, i: usize) -> GroupFor<E, &'a [E::Unit]> {
1939 self.inner.values.subslice(self.row_range(i)).into_inner()
1940 }
1941
1942 #[inline]
1944 pub fn symbolic(&self) -> SymbolicSparseRowMatRef<'a, I> {
1945 self.inner.symbolic
1946 }
1947
1948 #[inline]
1950 pub fn into_parts(self) -> (SymbolicSparseRowMatRef<'a, I>, GroupFor<E, &'a [E::Unit]>) {
1951 (self.inner.symbolic, self.inner.values.into_inner())
1952 }
1953}
1954
1955impl<'a, I: Index, E: Entity> SparseColMatRef<'a, I, E> {
1956 #[inline]
1963 #[track_caller]
1964 pub fn new(
1965 symbolic: SymbolicSparseColMatRef<'a, I>,
1966 values: GroupFor<E, &'a [E::Unit]>,
1967 ) -> Self {
1968 let values = SliceGroup::new(values);
1969 assert!(symbolic.row_indices().len() == values.len());
1970 Self {
1971 inner: inner::SparseColMatRef { symbolic, values },
1972 }
1973 }
1974
1975 #[inline]
1980 pub fn to_owned(&self) -> Result<SparseColMat<I, E::Canonical>, FaerError>
1981 where
1982 E: Conjugate,
1983 E::Canonical: ComplexField,
1984 {
1985 let symbolic = self.symbolic().to_owned()?;
1986 let mut values = VecGroup::<E::Canonical>::new();
1987
1988 values
1989 .try_reserve_exact(self.inner.values.len())
1990 .map_err(|_| FaerError::OutOfMemory)?;
1991
1992 values.resize(
1993 self.inner.values.len(),
1994 E::Canonical::faer_zero().faer_into_units(),
1995 );
1996
1997 let src = self.inner.values;
1998 let dst = values.as_slice_mut();
1999
2000 for (mut dst, src) in core::iter::zip(dst.into_mut_iter(), src.into_ref_iter()) {
2001 dst.write(src.read().canonicalize());
2002 }
2003
2004 Ok(SparseColMat {
2005 inner: inner::SparseColMat { symbolic, values },
2006 })
2007 }
2008
2009 #[inline]
2014 pub fn to_row_major(&self) -> Result<SparseRowMat<I, E::Canonical>, FaerError>
2015 where
2016 E: Conjugate,
2017 E::Canonical: ComplexField,
2018 {
2019 let mut col_ptr = try_zeroed::<I>(self.nrows + 1)?;
2020 let nnz = self.compute_nnz();
2021 let mut row_ind = try_zeroed::<I>(nnz)?;
2022 let mut values = VecGroup::<E::Canonical>::new();
2023 values
2024 .try_reserve_exact(nnz)
2025 .map_err(|_| FaerError::OutOfMemory)?;
2026 values.resize(nnz, E::Canonical::faer_zero().faer_into_units());
2027
2028 let mut mem = GlobalPodBuffer::try_new(StackReq::new::<I>(self.nrows))
2029 .map_err(|_| FaerError::OutOfMemory)?;
2030
2031 let (this, conj) = self.canonicalize();
2032
2033 if conj == Conj::No {
2034 util::transpose::<I, E::Canonical>(
2035 &mut col_ptr,
2036 &mut row_ind,
2037 values.as_slice_mut().into_inner(),
2038 this,
2039 PodStack::new(&mut mem),
2040 );
2041 } else {
2042 util::adjoint::<I, E::Canonical>(
2043 &mut col_ptr,
2044 &mut row_ind,
2045 values.as_slice_mut().into_inner(),
2046 this,
2047 PodStack::new(&mut mem),
2048 );
2049 }
2050
2051 let transpose = unsafe {
2052 SparseColMat::new(
2053 SymbolicSparseColMat::new_unchecked(self.ncols, self.nrows, col_ptr, None, row_ind),
2054 values.into_inner(),
2055 )
2056 };
2057
2058 Ok(transpose.into_transpose())
2059 }
2060
2061 #[inline]
2063 pub fn transpose(self) -> SparseRowMatRef<'a, I, E> {
2064 SparseRowMatRef {
2065 inner: inner::SparseRowMatRef {
2066 symbolic: SymbolicSparseRowMatRef {
2067 nrows: self.inner.symbolic.ncols,
2068 ncols: self.inner.symbolic.nrows,
2069 row_ptr: self.inner.symbolic.col_ptr,
2070 row_nnz: self.inner.symbolic.col_nnz,
2071 col_ind: self.inner.symbolic.row_ind,
2072 },
2073 values: self.inner.values,
2074 },
2075 }
2076 }
2077
2078 #[inline]
2080 pub fn conjugate(self) -> SparseColMatRef<'a, I, E::Conj>
2081 where
2082 E: Conjugate,
2083 {
2084 SparseColMatRef {
2085 inner: inner::SparseColMatRef {
2086 symbolic: self.inner.symbolic,
2087 values: unsafe {
2088 SliceGroup::<'a, E::Conj>::new(transmute_unchecked::<
2089 GroupFor<E, &[UnitFor<E::Conj>]>,
2090 GroupFor<E::Conj, &[UnitFor<E::Conj>]>,
2091 >(E::faer_map(
2092 self.inner.values.into_inner(),
2093 |slice| {
2094 let len = slice.len();
2095 core::slice::from_raw_parts(
2096 slice.as_ptr() as *const UnitFor<E> as *const UnitFor<E::Conj>,
2097 len,
2098 )
2099 },
2100 )))
2101 },
2102 },
2103 }
2104 }
2105
2106 #[inline]
2108 pub fn canonicalize(self) -> (SparseColMatRef<'a, I, E::Canonical>, Conj)
2109 where
2110 E: Conjugate,
2111 {
2112 (
2113 SparseColMatRef {
2114 inner: inner::SparseColMatRef {
2115 symbolic: self.inner.symbolic,
2116 values: unsafe {
2117 SliceGroup::<'a, E::Canonical>::new(transmute_unchecked::<
2118 GroupFor<E, &[UnitFor<E::Canonical>]>,
2119 GroupFor<E::Canonical, &[UnitFor<E::Canonical>]>,
2120 >(E::faer_map(
2121 self.inner.values.into_inner(),
2122 |slice| {
2123 let len = slice.len();
2124 core::slice::from_raw_parts(
2125 slice.as_ptr() as *const UnitFor<E>
2126 as *const UnitFor<E::Canonical>,
2127 len,
2128 )
2129 },
2130 )))
2131 },
2132 },
2133 },
2134 if coe::is_same::<E, E::Canonical>() {
2135 Conj::No
2136 } else {
2137 Conj::Yes
2138 },
2139 )
2140 }
2141
2142 #[inline]
2144 pub fn adjoint(self) -> SparseRowMatRef<'a, I, E::Conj>
2145 where
2146 E: Conjugate,
2147 {
2148 self.transpose().conjugate()
2149 }
2150
2151 #[inline]
2153 pub fn values(self) -> GroupFor<E, &'a [E::Unit]> {
2154 self.inner.values.into_inner()
2155 }
2156
2157 #[inline]
2163 #[track_caller]
2164 pub fn values_of_col(self, j: usize) -> GroupFor<E, &'a [E::Unit]> {
2165 self.inner.values.subslice(self.col_range(j)).into_inner()
2166 }
2167
2168 #[inline]
2170 pub fn symbolic(&self) -> SymbolicSparseColMatRef<'a, I> {
2171 self.inner.symbolic
2172 }
2173
2174 #[inline]
2176 pub fn into_parts(self) -> (SymbolicSparseColMatRef<'a, I>, GroupFor<E, &'a [E::Unit]>) {
2177 (self.inner.symbolic, self.inner.values.into_inner())
2178 }
2179}
2180
2181impl<I: Index, E: Entity> SparseColMat<I, E> {
2182 #[inline]
2189 #[track_caller]
2190 pub fn new(symbolic: SymbolicSparseColMat<I>, values: GroupFor<E, Vec<E::Unit>>) -> Self {
2191 let values = VecGroup::from_inner(values);
2192 assert!(symbolic.row_indices().len() == values.len());
2193 Self {
2194 inner: inner::SparseColMat { symbolic, values },
2195 }
2196 }
2197
2198 #[inline]
2203 pub fn to_owned(&self) -> Result<SparseColMat<I, E::Canonical>, FaerError>
2204 where
2205 E: Conjugate,
2206 E::Canonical: ComplexField,
2207 {
2208 self.as_ref().to_owned()
2209 }
2210
2211 #[inline]
2216 pub fn to_row_major(&self) -> Result<SparseRowMat<I, E::Canonical>, FaerError>
2217 where
2218 E: Conjugate,
2219 E::Canonical: ComplexField,
2220 {
2221 self.as_ref().to_row_major()
2222 }
2223
2224 #[inline]
2226 pub fn into_parts(self) -> (SymbolicSparseColMat<I>, GroupFor<E, Vec<E::Unit>>) {
2227 (self.inner.symbolic, self.inner.values.into_inner())
2228 }
2229
2230 #[inline]
2232 pub fn as_ref(&self) -> SparseColMatRef<'_, I, E> {
2233 SparseColMatRef {
2234 inner: inner::SparseColMatRef {
2235 symbolic: self.inner.symbolic.as_ref(),
2236 values: self.inner.values.as_slice(),
2237 },
2238 }
2239 }
2240
2241 #[inline]
2245 pub fn as_mut(&mut self) -> SparseColMatMut<'_, I, E> {
2246 SparseColMatMut {
2247 inner: inner::SparseColMatMut {
2248 symbolic: self.inner.symbolic.as_ref(),
2249 values: self.inner.values.as_slice_mut(),
2250 },
2251 }
2252 }
2253
2254 #[inline]
2256 pub fn values(&self) -> GroupFor<E, &'_ [E::Unit]> {
2257 self.inner.values.as_slice().into_inner()
2258 }
2259
2260 #[inline]
2262 pub fn values_mut(&mut self) -> GroupFor<E, &'_ mut [E::Unit]> {
2263 self.inner.values.as_slice_mut().into_inner()
2264 }
2265
2266 #[inline]
2271 pub fn into_transpose(self) -> SparseRowMat<I, E> {
2272 SparseRowMat {
2273 inner: inner::SparseRowMat {
2274 symbolic: SymbolicSparseRowMat {
2275 nrows: self.inner.symbolic.ncols,
2276 ncols: self.inner.symbolic.nrows,
2277 row_ptr: self.inner.symbolic.col_ptr,
2278 row_nnz: self.inner.symbolic.col_nnz,
2279 col_ind: self.inner.symbolic.row_ind,
2280 },
2281 values: self.inner.values,
2282 },
2283 }
2284 }
2285
2286 #[inline]
2288 pub fn into_conjugate(self) -> SparseColMat<I, E::Conj>
2289 where
2290 E: Conjugate,
2291 {
2292 SparseColMat {
2293 inner: inner::SparseColMat {
2294 symbolic: self.inner.symbolic,
2295 values: unsafe {
2296 VecGroup::<E::Conj>::from_inner(transmute_unchecked::<
2297 GroupFor<E, Vec<UnitFor<E::Conj>>>,
2298 GroupFor<E::Conj, Vec<UnitFor<E::Conj>>>,
2299 >(E::faer_map(
2300 self.inner.values.into_inner(),
2301 |mut slice| {
2302 let len = slice.len();
2303 let cap = slice.capacity();
2304 let ptr =
2305 slice.as_mut_ptr() as *mut UnitFor<E> as *mut UnitFor<E::Conj>;
2306
2307 Vec::from_raw_parts(ptr, len, cap)
2308 },
2309 )))
2310 },
2311 },
2312 }
2313 }
2314
2315 #[inline]
2317 pub fn into_adjoint(self) -> SparseRowMat<I, E::Conj>
2318 where
2319 E: Conjugate,
2320 {
2321 self.into_transpose().into_conjugate()
2322 }
2323}
2324
2325impl<I: Index, E: Entity> SparseRowMat<I, E> {
2326 #[inline]
2333 #[track_caller]
2334 pub fn new(symbolic: SymbolicSparseRowMat<I>, values: GroupFor<E, Vec<E::Unit>>) -> Self {
2335 let values = VecGroup::from_inner(values);
2336 assert!(symbolic.col_indices().len() == values.len());
2337 Self {
2338 inner: inner::SparseRowMat { symbolic, values },
2339 }
2340 }
2341
2342 #[inline]
2347 pub fn to_owned(&self) -> Result<SparseRowMat<I, E::Canonical>, FaerError>
2348 where
2349 E: Conjugate,
2350 E::Canonical: ComplexField,
2351 {
2352 self.as_ref().to_owned()
2353 }
2354
2355 #[inline]
2360 pub fn to_col_major(&self) -> Result<SparseColMat<I, E::Canonical>, FaerError>
2361 where
2362 E: Conjugate,
2363 E::Canonical: ComplexField,
2364 {
2365 self.as_ref().to_col_major()
2366 }
2367
2368 #[inline]
2370 pub fn into_parts(self) -> (SymbolicSparseRowMat<I>, GroupFor<E, Vec<E::Unit>>) {
2371 (self.inner.symbolic, self.inner.values.into_inner())
2372 }
2373
2374 #[inline]
2376 pub fn as_ref(&self) -> SparseRowMatRef<'_, I, E> {
2377 SparseRowMatRef {
2378 inner: inner::SparseRowMatRef {
2379 symbolic: self.inner.symbolic.as_ref(),
2380 values: self.inner.values.as_slice(),
2381 },
2382 }
2383 }
2384
2385 #[inline]
2389 pub fn as_mut(&mut self) -> SparseRowMatMut<'_, I, E> {
2390 SparseRowMatMut {
2391 inner: inner::SparseRowMatMut {
2392 symbolic: self.inner.symbolic.as_ref(),
2393 values: self.inner.values.as_slice_mut(),
2394 },
2395 }
2396 }
2397
2398 #[inline]
2400 pub fn values(&self) -> GroupFor<E, &'_ [E::Unit]> {
2401 self.inner.values.as_slice().into_inner()
2402 }
2403
2404 #[inline]
2406 pub fn values_mut(&mut self) -> GroupFor<E, &'_ mut [E::Unit]> {
2407 self.inner.values.as_slice_mut().into_inner()
2408 }
2409
2410 #[inline]
2415 pub fn into_transpose(self) -> SparseColMat<I, E> {
2416 SparseColMat {
2417 inner: inner::SparseColMat {
2418 symbolic: SymbolicSparseColMat {
2419 nrows: self.inner.symbolic.ncols,
2420 ncols: self.inner.symbolic.nrows,
2421 col_ptr: self.inner.symbolic.row_ptr,
2422 col_nnz: self.inner.symbolic.row_nnz,
2423 row_ind: self.inner.symbolic.col_ind,
2424 },
2425 values: self.inner.values,
2426 },
2427 }
2428 }
2429
2430 #[inline]
2432 pub fn into_conjugate(self) -> SparseRowMat<I, E::Conj>
2433 where
2434 E: Conjugate,
2435 {
2436 SparseRowMat {
2437 inner: inner::SparseRowMat {
2438 symbolic: self.inner.symbolic,
2439 values: unsafe {
2440 VecGroup::<E::Conj>::from_inner(transmute_unchecked::<
2441 GroupFor<E, Vec<UnitFor<E::Conj>>>,
2442 GroupFor<E::Conj, Vec<UnitFor<E::Conj>>>,
2443 >(E::faer_map(
2444 self.inner.values.into_inner(),
2445 |mut slice| {
2446 let len = slice.len();
2447 let cap = slice.capacity();
2448 let ptr =
2449 slice.as_mut_ptr() as *mut UnitFor<E> as *mut UnitFor<E::Conj>;
2450
2451 Vec::from_raw_parts(ptr, len, cap)
2452 },
2453 )))
2454 },
2455 },
2456 }
2457 }
2458
2459 #[inline]
2461 pub fn into_adjoint(self) -> SparseColMat<I, E::Conj>
2462 where
2463 E: Conjugate,
2464 {
2465 self.into_transpose().into_conjugate()
2466 }
2467}
2468
2469const _: () = {
2471 impl<'a, I: Index, E: Entity> core::ops::Deref for SparseRowMatMut<'a, I, E> {
2472 type Target = SymbolicSparseRowMatRef<'a, I>;
2473 #[inline]
2474 fn deref(&self) -> &Self::Target {
2475 &self.inner.symbolic
2476 }
2477 }
2478
2479 impl<'a, I: Index, E: Entity> core::ops::Deref for SparseColMatMut<'a, I, E> {
2480 type Target = SymbolicSparseColMatRef<'a, I>;
2481 #[inline]
2482 fn deref(&self) -> &Self::Target {
2483 &self.inner.symbolic
2484 }
2485 }
2486
2487 impl<'a, I: Index, E: Entity> core::ops::Deref for SparseRowMatRef<'a, I, E> {
2488 type Target = SymbolicSparseRowMatRef<'a, I>;
2489 #[inline]
2490 fn deref(&self) -> &Self::Target {
2491 &self.inner.symbolic
2492 }
2493 }
2494
2495 impl<'a, I: Index, E: Entity> core::ops::Deref for SparseColMatRef<'a, I, E> {
2496 type Target = SymbolicSparseColMatRef<'a, I>;
2497 #[inline]
2498 fn deref(&self) -> &Self::Target {
2499 &self.inner.symbolic
2500 }
2501 }
2502
2503 impl<I: Index, E: Entity> core::ops::Deref for SparseRowMat<I, E> {
2504 type Target = SymbolicSparseRowMat<I>;
2505 #[inline]
2506 fn deref(&self) -> &Self::Target {
2507 &self.inner.symbolic
2508 }
2509 }
2510
2511 impl<I: Index, E: Entity> core::ops::Deref for SparseColMat<I, E> {
2512 type Target = SymbolicSparseColMat<I>;
2513 #[inline]
2514 fn deref(&self) -> &Self::Target {
2515 &self.inner.symbolic
2516 }
2517 }
2518
2519 impl<'short, I: Index, E: Entity> ReborrowMut<'short> for SparseRowMatRef<'_, I, E> {
2520 type Target = SparseRowMatRef<'short, I, E>;
2521
2522 #[inline]
2523 fn rb_mut(&'short mut self) -> Self::Target {
2524 *self
2525 }
2526 }
2527
2528 impl<'short, I: Index, E: Entity> Reborrow<'short> for SparseRowMatRef<'_, I, E> {
2529 type Target = SparseRowMatRef<'short, I, E>;
2530
2531 #[inline]
2532 fn rb(&'short self) -> Self::Target {
2533 *self
2534 }
2535 }
2536
2537 impl<'a, I: Index, E: Entity> IntoConst for SparseRowMatRef<'a, I, E> {
2538 type Target = SparseRowMatRef<'a, I, E>;
2539
2540 #[inline]
2541 fn into_const(self) -> Self::Target {
2542 self
2543 }
2544 }
2545
2546 impl<'short, I: Index, E: Entity> ReborrowMut<'short> for SparseColMatRef<'_, I, E> {
2547 type Target = SparseColMatRef<'short, I, E>;
2548
2549 #[inline]
2550 fn rb_mut(&'short mut self) -> Self::Target {
2551 *self
2552 }
2553 }
2554
2555 impl<'short, I: Index, E: Entity> Reborrow<'short> for SparseColMatRef<'_, I, E> {
2556 type Target = SparseColMatRef<'short, I, E>;
2557
2558 #[inline]
2559 fn rb(&'short self) -> Self::Target {
2560 *self
2561 }
2562 }
2563
2564 impl<'a, I: Index, E: Entity> IntoConst for SparseColMatRef<'a, I, E> {
2565 type Target = SparseColMatRef<'a, I, E>;
2566
2567 #[inline]
2568 fn into_const(self) -> Self::Target {
2569 self
2570 }
2571 }
2572
2573 impl<'short, I: Index, E: Entity> ReborrowMut<'short> for SparseRowMatMut<'_, I, E> {
2574 type Target = SparseRowMatMut<'short, I, E>;
2575
2576 #[inline]
2577 fn rb_mut(&'short mut self) -> Self::Target {
2578 SparseRowMatMut {
2579 inner: inner::SparseRowMatMut {
2580 symbolic: self.inner.symbolic,
2581 values: self.inner.values.rb_mut(),
2582 },
2583 }
2584 }
2585 }
2586
2587 impl<'short, I: Index, E: Entity> Reborrow<'short> for SparseRowMatMut<'_, I, E> {
2588 type Target = SparseRowMatRef<'short, I, E>;
2589
2590 #[inline]
2591 fn rb(&'short self) -> Self::Target {
2592 SparseRowMatRef {
2593 inner: inner::SparseRowMatRef {
2594 symbolic: self.inner.symbolic,
2595 values: self.inner.values.rb(),
2596 },
2597 }
2598 }
2599 }
2600
2601 impl<'a, I: Index, E: Entity> IntoConst for SparseRowMatMut<'a, I, E> {
2602 type Target = SparseRowMatRef<'a, I, E>;
2603
2604 #[inline]
2605 fn into_const(self) -> Self::Target {
2606 SparseRowMatRef {
2607 inner: inner::SparseRowMatRef {
2608 symbolic: self.inner.symbolic,
2609 values: self.inner.values.into_const(),
2610 },
2611 }
2612 }
2613 }
2614
2615 impl<'short, I: Index, E: Entity> ReborrowMut<'short> for SparseColMatMut<'_, I, E> {
2616 type Target = SparseColMatMut<'short, I, E>;
2617
2618 #[inline]
2619 fn rb_mut(&'short mut self) -> Self::Target {
2620 SparseColMatMut {
2621 inner: inner::SparseColMatMut {
2622 symbolic: self.inner.symbolic,
2623 values: self.inner.values.rb_mut(),
2624 },
2625 }
2626 }
2627 }
2628
2629 impl<'short, I: Index, E: Entity> Reborrow<'short> for SparseColMatMut<'_, I, E> {
2630 type Target = SparseColMatRef<'short, I, E>;
2631
2632 #[inline]
2633 fn rb(&'short self) -> Self::Target {
2634 SparseColMatRef {
2635 inner: inner::SparseColMatRef {
2636 symbolic: self.inner.symbolic,
2637 values: self.inner.values.rb(),
2638 },
2639 }
2640 }
2641 }
2642
2643 impl<'a, I: Index, E: Entity> IntoConst for SparseColMatMut<'a, I, E> {
2644 type Target = SparseColMatRef<'a, I, E>;
2645
2646 #[inline]
2647 fn into_const(self) -> Self::Target {
2648 SparseColMatRef {
2649 inner: inner::SparseColMatRef {
2650 symbolic: self.inner.symbolic,
2651 values: self.inner.values.into_const(),
2652 },
2653 }
2654 }
2655 }
2656};
2657
2658#[derive(Copy, Clone, Debug, Eq, PartialEq, Ord, PartialOrd)]
2660#[non_exhaustive]
2661pub enum CreationError {
2662 Generic(FaerError),
2664 OutOfBounds {
2666 row: usize,
2668 col: usize,
2670 },
2671}
2672
2673impl From<FaerError> for CreationError {
2674 #[inline]
2675 fn from(value: FaerError) -> Self {
2676 Self::Generic(value)
2677 }
2678}
2679impl core::fmt::Display for CreationError {
2680 #[inline]
2681 fn fmt(&self, f: &mut core::fmt::Formatter<'_>) -> core::fmt::Result {
2682 core::fmt::Debug::fmt(self, f)
2683 }
2684}
2685
2686#[cfg(feature = "std")]
2687impl std::error::Error for CreationError {}
2688
2689#[inline]
2690#[track_caller]
2691fn try_zeroed<I: bytemuck::Pod>(n: usize) -> Result<alloc::vec::Vec<I>, FaerError> {
2692 let mut v = alloc::vec::Vec::new();
2693 v.try_reserve_exact(n).map_err(|_| FaerError::OutOfMemory)?;
2694 unsafe {
2695 core::ptr::write_bytes::<I>(v.as_mut_ptr(), 0u8, n);
2696 v.set_len(n);
2697 }
2698 Ok(v)
2699}
2700
2701#[inline]
2702#[track_caller]
2703fn try_collect<I: IntoIterator>(iter: I) -> Result<alloc::vec::Vec<I::Item>, FaerError> {
2704 let iter = iter.into_iter();
2705 let mut v = alloc::vec::Vec::new();
2706 v.try_reserve_exact(iter.size_hint().0)
2707 .map_err(|_| FaerError::OutOfMemory)?;
2708 v.extend(iter);
2709 Ok(v)
2710}
2711
2712#[derive(Debug, Clone)]
2716pub struct ValuesOrder<I> {
2717 argsort: Vec<usize>,
2718 all_nnz: usize,
2719 nnz: usize,
2720 __marker: PhantomData<I>,
2721}
2722
2723#[derive(Debug, Copy, Clone)]
2725pub enum FillMode {
2726 Replace,
2728 Add,
2730}
2731
2732const _: () = {
2734 const TOP_BIT: usize = 1usize << (usize::BITS - 1);
2735 const TOP_BIT_MASK: usize = TOP_BIT - 1;
2736
2737 impl<I: Index> SymbolicSparseColMat<I> {
2738 fn try_new_from_indices_impl(
2739 nrows: usize,
2740 ncols: usize,
2741 indices: impl Fn(usize) -> (I, I),
2742 all_nnz: usize,
2743 ) -> Result<(Self, ValuesOrder<I>), CreationError> {
2744 if nrows > I::Signed::MAX.zx() || ncols > I::Signed::MAX.zx() {
2745 return Err(CreationError::Generic(FaerError::IndexOverflow));
2746 }
2747
2748 if all_nnz == 0 {
2749 return Ok((
2750 Self {
2751 nrows,
2752 ncols,
2753 col_ptr: try_zeroed(ncols + 1)?,
2754 col_nnz: None,
2755 row_ind: Vec::new(),
2756 },
2757 ValuesOrder {
2758 argsort: Vec::new(),
2759 all_nnz,
2760 nnz: 0,
2761 __marker: PhantomData,
2762 },
2763 ));
2764 }
2765
2766 let mut argsort = try_collect(0..all_nnz)?;
2767 argsort.sort_unstable_by_key(|&i| {
2768 let (row, col) = indices(i);
2769 (col, row)
2770 });
2771
2772 let mut n_duplicates = 0usize;
2773 let mut current_bit = 0usize;
2774
2775 let mut prev = indices(argsort[0]);
2776 for i in 1..all_nnz {
2777 let idx = indices(argsort[i]);
2778 let same_as_prev = idx == prev;
2779 prev = idx;
2780 current_bit = ((current_bit == ((same_as_prev as usize) << (usize::BITS - 1)))
2781 as usize)
2782 << (usize::BITS - 1);
2783 argsort[i] |= current_bit;
2784
2785 n_duplicates += same_as_prev as usize;
2786 }
2787
2788 let nnz = all_nnz - n_duplicates;
2789 if nnz > I::Signed::MAX.zx() {
2790 return Err(CreationError::Generic(FaerError::IndexOverflow));
2791 }
2792
2793 let mut col_ptr = try_zeroed::<I>(ncols + 1)?;
2794 let mut row_ind = try_zeroed::<I>(nnz)?;
2795
2796 let mut original_pos = 0usize;
2797 let mut new_pos = 0usize;
2798
2799 for j in 0..ncols {
2800 let mut n_unique = 0usize;
2801
2802 while original_pos < all_nnz {
2803 let (row, col) = indices(argsort[original_pos] & TOP_BIT_MASK);
2804 if row.zx() >= nrows || col.zx() >= ncols {
2805 return Err(CreationError::OutOfBounds {
2806 row: row.zx(),
2807 col: col.zx(),
2808 });
2809 }
2810
2811 if col.zx() != j {
2812 break;
2813 }
2814
2815 row_ind[new_pos] = row;
2816
2817 n_unique += 1;
2818
2819 new_pos += 1;
2820 original_pos += 1;
2821
2822 while original_pos < all_nnz
2823 && indices(argsort[original_pos] & TOP_BIT_MASK) == (row, col)
2824 {
2825 original_pos += 1;
2826 }
2827 }
2828
2829 col_ptr[j + 1] = col_ptr[j] + I::truncate(n_unique);
2830 }
2831
2832 Ok((
2833 Self {
2834 nrows,
2835 ncols,
2836 col_ptr,
2837 col_nnz: None,
2838 row_ind,
2839 },
2840 ValuesOrder {
2841 argsort,
2842 all_nnz,
2843 nnz,
2844 __marker: PhantomData,
2845 },
2846 ))
2847 }
2848
2849 fn try_new_from_nonnegative_indices_impl(
2850 nrows: usize,
2851 ncols: usize,
2852 indices: impl Fn(usize) -> (I::Signed, I::Signed),
2853 all_nnz: usize,
2854 ) -> Result<(Self, ValuesOrder<I>), CreationError> {
2855 if nrows > I::Signed::MAX.zx() || ncols > I::Signed::MAX.zx() {
2856 return Err(CreationError::Generic(FaerError::IndexOverflow));
2857 }
2858
2859 let mut argsort = try_collect(0..all_nnz)?;
2860 argsort.sort_unstable_by_key(|&i| {
2861 let (row, col) = indices(i);
2862 let ignore = (row < I::Signed::truncate(0)) | (col < I::Signed::truncate(0));
2863 (ignore, col, row)
2864 });
2865
2866 let all_nnz = argsort.partition_point(|&i| {
2867 let (row, col) = indices(i);
2868 let ignore = (row < I::Signed::truncate(0)) | (col < I::Signed::truncate(0));
2869 !ignore
2870 });
2871
2872 if all_nnz == 0 {
2873 return Ok((
2874 Self {
2875 nrows,
2876 ncols,
2877 col_ptr: try_zeroed(ncols + 1)?,
2878 col_nnz: None,
2879 row_ind: Vec::new(),
2880 },
2881 ValuesOrder {
2882 argsort: Vec::new(),
2883 all_nnz,
2884 nnz: 0,
2885 __marker: PhantomData,
2886 },
2887 ));
2888 }
2889
2890 let mut n_duplicates = 0usize;
2891 let mut current_bit = 0usize;
2892
2893 let mut prev = indices(argsort[0]);
2894
2895 for i in 1..all_nnz {
2896 let idx = indices(argsort[i]);
2897 let same_as_prev = idx == prev;
2898 prev = idx;
2899 current_bit = ((current_bit == ((same_as_prev as usize) << (usize::BITS - 1)))
2900 as usize)
2901 << (usize::BITS - 1);
2902 argsort[i] |= current_bit;
2903
2904 n_duplicates += same_as_prev as usize;
2905 }
2906
2907 let nnz = all_nnz - n_duplicates;
2908 if nnz > I::Signed::MAX.zx() {
2909 return Err(CreationError::Generic(FaerError::IndexOverflow));
2910 }
2911
2912 let mut col_ptr = try_zeroed::<I>(ncols + 1)?;
2913 let mut row_ind = try_zeroed::<I>(nnz)?;
2914
2915 let mut original_pos = 0usize;
2916 let mut new_pos = 0usize;
2917
2918 for j in 0..ncols {
2919 let mut n_unique = 0usize;
2920
2921 while original_pos < all_nnz {
2922 let (row, col) = indices(argsort[original_pos] & TOP_BIT_MASK);
2923 if row.zx() >= nrows || col.zx() >= ncols {
2924 return Err(CreationError::OutOfBounds {
2925 row: row.zx(),
2926 col: col.zx(),
2927 });
2928 }
2929
2930 if col.zx() != j {
2931 break;
2932 }
2933
2934 row_ind[new_pos] = I::from_signed(row);
2935
2936 n_unique += 1;
2937
2938 new_pos += 1;
2939 original_pos += 1;
2940
2941 while original_pos < all_nnz
2942 && indices(argsort[original_pos] & TOP_BIT_MASK) == (row, col)
2943 {
2944 original_pos += 1;
2945 }
2946 }
2947
2948 col_ptr[j + 1] = col_ptr[j] + I::truncate(n_unique);
2949 }
2950
2951 Ok((
2952 Self {
2953 nrows,
2954 ncols,
2955 col_ptr,
2956 col_nnz: None,
2957 row_ind,
2958 },
2959 ValuesOrder {
2960 argsort,
2961 all_nnz,
2962 nnz,
2963 __marker: PhantomData,
2964 },
2965 ))
2966 }
2967
2968 #[inline]
2971 pub fn try_new_from_indices(
2972 nrows: usize,
2973 ncols: usize,
2974 indices: &[(I, I)],
2975 ) -> Result<(Self, ValuesOrder<I>), CreationError> {
2976 Self::try_new_from_indices_impl(nrows, ncols, |i| indices[i], indices.len())
2977 }
2978
2979 #[inline]
2984 pub fn try_new_from_nonnegative_indices(
2985 nrows: usize,
2986 ncols: usize,
2987 indices: &[(I::Signed, I::Signed)],
2988 ) -> Result<(Self, ValuesOrder<I>), CreationError> {
2989 Self::try_new_from_nonnegative_indices_impl(nrows, ncols, |i| indices[i], indices.len())
2990 }
2991 }
2992
2993 impl<I: Index> SymbolicSparseRowMat<I> {
2994 #[inline]
2997 pub fn try_new_from_indices(
2998 nrows: usize,
2999 ncols: usize,
3000 indices: &[(I, I)],
3001 ) -> Result<(Self, ValuesOrder<I>), CreationError> {
3002 SymbolicSparseColMat::try_new_from_indices_impl(
3003 ncols,
3004 nrows,
3005 |i| {
3006 let (row, col) = indices[i];
3007 (col, row)
3008 },
3009 indices.len(),
3010 )
3011 .map(|(m, o)| (m.into_transpose(), o))
3012 }
3013
3014 #[inline]
3019 pub fn try_new_from_nonnegative_indices(
3020 nrows: usize,
3021 ncols: usize,
3022 indices: &[(I::Signed, I::Signed)],
3023 ) -> Result<(Self, ValuesOrder<I>), CreationError> {
3024 SymbolicSparseColMat::try_new_from_nonnegative_indices_impl(
3025 ncols,
3026 nrows,
3027 |i| {
3028 let (row, col) = indices[i];
3029 (col, row)
3030 },
3031 indices.len(),
3032 )
3033 .map(|(m, o)| (m.into_transpose(), o))
3034 }
3035 }
3036
3037 impl<I: Index, E: ComplexField> SparseColMat<I, E> {
3038 #[track_caller]
3039 fn new_from_order_and_values_impl(
3040 symbolic: SymbolicSparseColMat<I>,
3041 order: &ValuesOrder<I>,
3042 all_values: impl Fn(usize) -> E,
3043 values_len: usize,
3044 ) -> Result<Self, FaerError> {
3045 {
3046 let nnz = order.argsort.len();
3047 assert!(values_len == nnz);
3048 }
3049
3050 let all_nnz = order.all_nnz;
3051
3052 let mut values = VecGroup::<E>::new();
3053 match values.try_reserve_exact(order.nnz) {
3054 Ok(()) => {}
3055 Err(_) => return Err(FaerError::OutOfMemory),
3056 };
3057
3058 let mut pos = 0usize;
3059 let mut pos_unique = usize::MAX;
3060 let mut current_bit = TOP_BIT;
3061
3062 while pos < all_nnz {
3063 let argsort_pos = order.argsort[pos];
3064 let extracted_bit = argsort_pos & TOP_BIT;
3065 let argsort_pos = argsort_pos & TOP_BIT_MASK;
3066
3067 let val = all_values(argsort_pos);
3068 if extracted_bit != current_bit {
3069 values.push(val.faer_into_units());
3070 pos_unique = pos_unique.wrapping_add(1);
3071 } else {
3072 let old_val = values.as_slice().read(pos_unique);
3073 values
3074 .as_slice_mut()
3075 .write(pos_unique, old_val.faer_add(val));
3076 }
3077
3078 current_bit = extracted_bit;
3079
3080 pos += 1;
3081 }
3082
3083 Ok(Self {
3084 inner: inner::SparseColMat { symbolic, values },
3085 })
3086 }
3087
3088 #[track_caller]
3092 pub fn new_from_order_and_values(
3093 symbolic: SymbolicSparseColMat<I>,
3094 order: &ValuesOrder<I>,
3095 values: GroupFor<E, &[E::Unit]>,
3096 ) -> Result<Self, FaerError> {
3097 let values = SliceGroup::<'_, E>::new(values);
3098 Self::new_from_order_and_values_impl(symbolic, order, |i| values.read(i), values.len())
3099 }
3100
3101 #[track_caller]
3103 pub fn try_new_from_triplets(
3104 nrows: usize,
3105 ncols: usize,
3106 triplets: &[(I, I, E)],
3107 ) -> Result<Self, CreationError> {
3108 let (symbolic, order) = SymbolicSparseColMat::try_new_from_indices_impl(
3109 nrows,
3110 ncols,
3111 |i| {
3112 let (row, col, _) = triplets[i];
3113 (row, col)
3114 },
3115 triplets.len(),
3116 )?;
3117 Ok(Self::new_from_order_and_values_impl(
3118 symbolic,
3119 &order,
3120 |i| triplets[i].2,
3121 triplets.len(),
3122 )?)
3123 }
3124
3125 #[track_caller]
3127 pub fn try_new_from_nonnegative_triplets(
3128 nrows: usize,
3129 ncols: usize,
3130 triplets: &[(I::Signed, I::Signed, E)],
3131 ) -> Result<Self, CreationError> {
3132 let (symbolic, order) =
3133 SymbolicSparseColMat::<I>::try_new_from_nonnegative_indices_impl(
3134 nrows,
3135 ncols,
3136 |i| {
3137 let (row, col, _) = triplets[i];
3138 (row, col)
3139 },
3140 triplets.len(),
3141 )?;
3142 Ok(Self::new_from_order_and_values_impl(
3143 symbolic,
3144 &order,
3145 |i| triplets[i].2,
3146 triplets.len(),
3147 )?)
3148 }
3149 }
3150
3151 impl<I: Index, E: ComplexField> SparseRowMat<I, E> {
3152 #[track_caller]
3156 pub fn new_from_order_and_values(
3157 symbolic: SymbolicSparseRowMat<I>,
3158 order: &ValuesOrder<I>,
3159 values: GroupFor<E, &[E::Unit]>,
3160 ) -> Result<Self, FaerError> {
3161 SparseColMat::new_from_order_and_values(symbolic.into_transpose(), order, values)
3162 .map(SparseColMat::into_transpose)
3163 }
3164
3165 #[track_caller]
3167 pub fn try_new_from_triplets(
3168 nrows: usize,
3169 ncols: usize,
3170 triplets: &[(I, I, E)],
3171 ) -> Result<Self, CreationError> {
3172 let (symbolic, order) = SymbolicSparseColMat::try_new_from_indices_impl(
3173 ncols,
3174 nrows,
3175 |i| {
3176 let (row, col, _) = triplets[i];
3177 (col, row)
3178 },
3179 triplets.len(),
3180 )?;
3181 Ok(SparseColMat::new_from_order_and_values_impl(
3182 symbolic,
3183 &order,
3184 |i| triplets[i].2,
3185 triplets.len(),
3186 )?
3187 .into_transpose())
3188 }
3189
3190 #[track_caller]
3192 pub fn try_new_from_nonnegative_triplets(
3193 nrows: usize,
3194 ncols: usize,
3195 triplets: &[(I::Signed, I::Signed, E)],
3196 ) -> Result<Self, CreationError> {
3197 let (symbolic, order) =
3198 SymbolicSparseColMat::<I>::try_new_from_nonnegative_indices_impl(
3199 ncols,
3200 nrows,
3201 |i| {
3202 let (row, col, _) = triplets[i];
3203 (col, row)
3204 },
3205 triplets.len(),
3206 )?;
3207 Ok(SparseColMat::new_from_order_and_values_impl(
3208 symbolic,
3209 &order,
3210 |i| triplets[i].2,
3211 triplets.len(),
3212 )?
3213 .into_transpose())
3214 }
3215 }
3216
3217 impl<I: Index, E: ComplexField> SparseColMatMut<'_, I, E> {
3218 pub fn fill_from_order_and_values(
3225 &mut self,
3226 order: &ValuesOrder<I>,
3227 values: GroupFor<E, &[E::Unit]>,
3228 mode: FillMode,
3229 ) {
3230 let values = SliceGroup::<'_, E>::new(values);
3231
3232 {
3233 let nnz = order.argsort.len();
3234 assert!(values.len() == nnz);
3235 assert!(order.nnz == self.inner.values.len());
3236 }
3237 let all_nnz = order.all_nnz;
3238 let mut dst = self.inner.values.rb_mut();
3239
3240 let mut pos = 0usize;
3241 let mut pos_unique = usize::MAX;
3242 let mut current_bit = TOP_BIT;
3243
3244 match mode {
3245 FillMode::Replace => {
3246 while pos < all_nnz {
3247 let argsort_pos = order.argsort[pos];
3248 let extracted_bit = argsort_pos & TOP_BIT;
3249 let argsort_pos = argsort_pos & TOP_BIT_MASK;
3250
3251 let val = values.read(argsort_pos);
3252 if extracted_bit != current_bit {
3253 pos_unique = pos_unique.wrapping_add(1);
3254 dst.write(pos_unique, val);
3255 } else {
3256 let old_val = dst.read(pos_unique);
3257 dst.write(pos_unique, old_val.faer_add(val));
3258 }
3259
3260 current_bit = extracted_bit;
3261
3262 pos += 1;
3263 }
3264 }
3265 FillMode::Add => {
3266 while pos < all_nnz {
3267 let argsort_pos = order.argsort[pos];
3268 let extracted_bit = argsort_pos & TOP_BIT;
3269 let argsort_pos = argsort_pos & TOP_BIT_MASK;
3270
3271 let val = values.read(argsort_pos);
3272 if extracted_bit != current_bit {
3273 pos_unique = pos_unique.wrapping_add(1);
3274 }
3275
3276 let old_val = dst.read(pos_unique);
3277 dst.write(pos_unique, old_val.faer_add(val));
3278
3279 current_bit = extracted_bit;
3280
3281 pos += 1;
3282 }
3283 }
3284 }
3285 }
3286 }
3287
3288 impl<I: Index, E: ComplexField> SparseRowMatMut<'_, I, E> {
3289 pub fn fill_from_order_and_values(
3296 &mut self,
3297 order: &ValuesOrder<I>,
3298 values: GroupFor<E, &[E::Unit]>,
3299 mode: FillMode,
3300 ) {
3301 self.rb_mut()
3302 .transpose_mut()
3303 .fill_from_order_and_values(order, values, mode);
3304 }
3305 }
3306};
3307
3308pub mod util {
3310 use super::*;
3311 use crate::{assert, debug_assert};
3312
3313 pub fn sort_indices<I: Index, E: Entity>(
3315 col_ptrs: &[I],
3316 row_indices: &mut [I],
3317 values: GroupFor<E, &mut [E::Unit]>,
3318 ) {
3319 assert!(col_ptrs.len() >= 1);
3320 let mut values = SliceGroupMut::<'_, E>::new(values);
3321
3322 let n = col_ptrs.len() - 1;
3323 for j in 0..n {
3324 let start = col_ptrs[j].zx();
3325 let end = col_ptrs[j + 1].zx();
3326
3327 unsafe {
3328 crate::sort::sort_indices(
3329 &mut row_indices[start..end],
3330 values.rb_mut().subslice(start..end),
3331 );
3332 }
3333 }
3334 }
3335
3336 #[doc(hidden)]
3337 pub unsafe fn ghost_permute_hermitian_unsorted<'n, 'out, I: Index, E: ComplexField>(
3338 new_values: SliceGroupMut<'out, E>,
3339 new_col_ptrs: &'out mut [I],
3340 new_row_indices: &'out mut [I],
3341 A: ghost::SparseColMatRef<'n, 'n, '_, I, E>,
3342 perm: ghost::PermutationRef<'n, '_, I, E>,
3343 in_side: Side,
3344 out_side: Side,
3345 sort: bool,
3346 stack: PodStack<'_>,
3347 ) -> ghost::SparseColMatMut<'n, 'n, 'out, I, E> {
3348 let N = A.ncols();
3349 let n = *A.ncols();
3350
3351 assert!(new_col_ptrs.len() == n + 1);
3353 let (_, perm_inv) = perm.into_arrays();
3354
3355 let (current_row_position, _) = stack.make_raw::<I>(n);
3356 let current_row_position = ghost::Array::from_mut(current_row_position, N);
3357
3358 mem::fill_zero(current_row_position.as_mut());
3359 let col_counts = &mut *current_row_position;
3360 match (in_side, out_side) {
3361 (Side::Lower, Side::Lower) => {
3362 for old_j in N.indices() {
3363 let new_j = perm_inv[old_j].zx();
3364 for old_i in A.row_indices_of_col(old_j) {
3365 if old_i >= old_j {
3366 let new_i = perm_inv[old_i].zx();
3367 let new_min = Ord::min(new_i, new_j);
3368 col_counts[new_min] += I::truncate(1);
3371 }
3372 }
3373 }
3374 }
3375 (Side::Lower, Side::Upper) => {
3376 for old_j in N.indices() {
3377 let new_j = perm_inv[old_j].zx();
3378 for old_i in A.row_indices_of_col(old_j) {
3379 if old_i >= old_j {
3380 let new_i = perm_inv[old_i].zx();
3381 let new_max = Ord::max(new_i, new_j);
3382 col_counts[new_max] += I::truncate(1);
3385 }
3386 }
3387 }
3388 }
3389 (Side::Upper, Side::Lower) => {
3390 for old_j in N.indices() {
3391 let new_j = perm_inv[old_j].zx();
3392 for old_i in A.row_indices_of_col(old_j) {
3393 if old_i <= old_j {
3394 let new_i = perm_inv[old_i].zx();
3395 let new_min = Ord::min(new_i, new_j);
3396 col_counts[new_min] += I::truncate(1);
3399 }
3400 }
3401 }
3402 }
3403 (Side::Upper, Side::Upper) => {
3404 for old_j in N.indices() {
3405 let new_j = perm_inv[old_j].zx();
3406 for old_i in A.row_indices_of_col(old_j) {
3407 if old_i <= old_j {
3408 let new_i = perm_inv[old_i].zx();
3409 let new_max = Ord::max(new_i, new_j);
3410 col_counts[new_max] += I::truncate(1);
3413 }
3414 }
3415 }
3416 }
3417 }
3418
3419 new_col_ptrs[0] = I::truncate(0);
3424 for (count, [ci0, ci1]) in zip(
3425 col_counts.as_mut(),
3426 windows2(Cell::as_slice_of_cells(Cell::from_mut(&mut *new_col_ptrs))),
3427 ) {
3428 let ci0 = ci0.get();
3429 ci1.set(ci0 + *count);
3430 *count = ci0;
3431 }
3432 let nnz = new_col_ptrs[n].zx();
3435 let new_row_indices = &mut new_row_indices[..nnz];
3436 let mut new_values = new_values.subslice(0..nnz);
3437
3438 ghost::Size::with(
3439 nnz,
3440 #[inline(always)]
3441 |NNZ| {
3442 let mut new_values =
3443 ghost::ArrayGroupMut::new(new_values.rb_mut().into_inner(), NNZ);
3444 let new_row_indices = ghost::Array::from_mut(new_row_indices, NNZ);
3445
3446 let conj_if = |cond: bool, x: E| {
3447 if !coe::is_same::<E, E::Real>() && cond {
3448 x.faer_conj()
3449 } else {
3450 x
3451 }
3452 };
3453
3454 match (in_side, out_side) {
3455 (Side::Lower, Side::Lower) => {
3456 for old_j in N.indices() {
3457 let new_j_ = perm_inv[old_j];
3458 let new_j = new_j_.zx();
3459
3460 for (old_i, val) in zip(
3461 A.row_indices_of_col(old_j),
3462 SliceGroup::<'_, E>::new(A.values_of_col(old_j)).into_ref_iter(),
3463 ) {
3464 if old_i >= old_j {
3465 let new_i_ = perm_inv[old_i];
3466 let new_i = new_i_.zx();
3467
3468 let new_max = Ord::max(new_i_, new_j_);
3469 let new_min = Ord::min(new_i, new_j);
3470 let current_row_pos: &mut I =
3471 &mut current_row_position[new_min];
3472 let row_pos = unsafe {
3474 ghost::Idx::new_unchecked(current_row_pos.zx(), NNZ)
3475 };
3476 *current_row_pos += I::truncate(1);
3477 new_values
3478 .write(row_pos, conj_if(new_min == new_i, val.read()));
3479 new_row_indices[row_pos] = *new_max;
3481 }
3482 }
3483 }
3484 }
3485 (Side::Lower, Side::Upper) => {
3486 for old_j in N.indices() {
3487 let new_j_ = perm_inv[old_j];
3488 let new_j = new_j_.zx();
3489
3490 for (old_i, val) in zip(
3491 A.row_indices_of_col(old_j),
3492 SliceGroup::<'_, E>::new(A.values_of_col(old_j)).into_ref_iter(),
3493 ) {
3494 if old_i >= old_j {
3495 let new_i_ = perm_inv[old_i];
3496 let new_i = new_i_.zx();
3497
3498 let new_max = Ord::max(new_i, new_j);
3499 let new_min = Ord::min(new_i_, new_j_);
3500 let current_row_pos = &mut current_row_position[new_max];
3501 let row_pos = unsafe {
3503 ghost::Idx::new_unchecked(current_row_pos.zx(), NNZ)
3504 };
3505 *current_row_pos += I::truncate(1);
3506 new_values
3507 .write(row_pos, conj_if(new_max == new_i, val.read()));
3508 new_row_indices[row_pos] = *new_min;
3510 }
3511 }
3512 }
3513 }
3514 (Side::Upper, Side::Lower) => {
3515 for old_j in N.indices() {
3516 let new_j_ = perm_inv[old_j];
3517 let new_j = new_j_.zx();
3518
3519 for (old_i, val) in zip(
3520 A.row_indices_of_col(old_j),
3521 SliceGroup::<'_, E>::new(A.values_of_col(old_j)).into_ref_iter(),
3522 ) {
3523 if old_i <= old_j {
3524 let new_i_ = perm_inv[old_i];
3525 let new_i = new_i_.zx();
3526
3527 let new_max = Ord::max(new_i_, new_j_);
3528 let new_min = Ord::min(new_i, new_j);
3529 let current_row_pos = &mut current_row_position[new_min];
3530 let row_pos = unsafe {
3532 ghost::Idx::new_unchecked(current_row_pos.zx(), NNZ)
3533 };
3534 *current_row_pos += I::truncate(1);
3535 new_values
3536 .write(row_pos, conj_if(new_min == new_i, val.read()));
3537 new_row_indices[row_pos] = *new_max;
3539 }
3540 }
3541 }
3542 }
3543 (Side::Upper, Side::Upper) => {
3544 for old_j in N.indices() {
3545 let new_j_ = perm_inv[old_j];
3546 let new_j = new_j_.zx();
3547
3548 for (old_i, val) in zip(
3549 A.row_indices_of_col(old_j),
3550 SliceGroup::<'_, E>::new(A.values_of_col(old_j)).into_ref_iter(),
3551 ) {
3552 if old_i <= old_j {
3553 let new_i_ = perm_inv[old_i];
3554 let new_i = new_i_.zx();
3555
3556 let new_max = Ord::max(new_i, new_j);
3557 let new_min = Ord::min(new_i_, new_j_);
3558 let current_row_pos = &mut current_row_position[new_max];
3559 let row_pos = unsafe {
3561 ghost::Idx::new_unchecked(current_row_pos.zx(), NNZ)
3562 };
3563 *current_row_pos += I::truncate(1);
3564 new_values
3565 .write(row_pos, conj_if(new_max == new_i, val.read()));
3566 new_row_indices[row_pos] = *new_min;
3568 }
3569 }
3570 }
3571 }
3572 }
3573 debug_assert!(current_row_position.as_ref() == &new_col_ptrs[1..]);
3574 },
3575 );
3576
3577 if sort {
3578 sort_indices::<I, E>(
3579 new_col_ptrs,
3580 new_row_indices,
3581 new_values.rb_mut().into_inner(),
3582 );
3583 }
3584
3585 unsafe {
3590 ghost::SparseColMatMut::new(
3591 SparseColMatMut::new(
3592 SymbolicSparseColMatRef::new_unchecked(
3593 n,
3594 n,
3595 new_col_ptrs,
3596 None,
3597 new_row_indices,
3598 ),
3599 new_values.into_inner(),
3600 ),
3601 N,
3602 N,
3603 )
3604 }
3605 }
3606
3607 #[doc(hidden)]
3608 pub unsafe fn ghost_permute_hermitian_unsorted_symbolic<'n, 'out, I: Index>(
3609 new_col_ptrs: &'out mut [I],
3610 new_row_indices: &'out mut [I],
3611 A: ghost::SymbolicSparseColMatRef<'n, 'n, '_, I>,
3612 perm: ghost::PermutationRef<'n, '_, I, Symbolic>,
3613 in_side: Side,
3614 out_side: Side,
3615 stack: PodStack<'_>,
3616 ) -> ghost::SymbolicSparseColMatRef<'n, 'n, 'out, I> {
3617 let old_values = &*Symbolic::materialize(A.into_inner().row_indices().len());
3618 let new_values = Symbolic::materialize(new_row_indices.len());
3619 *ghost_permute_hermitian_unsorted(
3620 SliceGroupMut::<'_, Symbolic>::new(new_values),
3621 new_col_ptrs,
3622 new_row_indices,
3623 ghost::SparseColMatRef::new(
3624 SparseColMatRef::new(A.into_inner(), old_values),
3625 A.nrows(),
3626 A.ncols(),
3627 ),
3628 perm,
3629 in_side,
3630 out_side,
3631 false,
3632 stack,
3633 )
3634 }
3635
3636 #[doc(hidden)]
3641 pub unsafe fn permute_hermitian_unsorted<'out, I: Index, E: ComplexField>(
3642 new_values: GroupFor<E, &'out mut [E::Unit]>,
3643 new_col_ptrs: &'out mut [I],
3644 new_row_indices: &'out mut [I],
3645 A: SparseColMatRef<'_, I, E>,
3646 perm: crate::permutation::PermutationRef<'_, I, E>,
3647 in_side: Side,
3648 out_side: Side,
3649 stack: PodStack<'_>,
3650 ) -> SparseColMatMut<'out, I, E> {
3651 ghost::Size::with(A.nrows(), |N| {
3652 assert!(A.nrows() == A.ncols());
3653 ghost_permute_hermitian_unsorted(
3654 SliceGroupMut::new(new_values),
3655 new_col_ptrs,
3656 new_row_indices,
3657 ghost::SparseColMatRef::new(A, N, N),
3658 ghost::PermutationRef::new(perm, N),
3659 in_side,
3660 out_side,
3661 false,
3662 stack,
3663 )
3664 .into_inner()
3665 })
3666 }
3667
3668 pub fn permute_hermitian<'out, I: Index, E: ComplexField>(
3676 new_values: GroupFor<E, &'out mut [E::Unit]>,
3677 new_col_ptrs: &'out mut [I],
3678 new_row_indices: &'out mut [I],
3679 A: SparseColMatRef<'_, I, E>,
3680 perm: crate::permutation::PermutationRef<'_, I, E>,
3681 in_side: Side,
3682 out_side: Side,
3683 stack: PodStack<'_>,
3684 ) -> SparseColMatMut<'out, I, E> {
3685 ghost::Size::with(A.nrows(), |N| {
3686 assert!(A.nrows() == A.ncols());
3687 unsafe {
3688 ghost_permute_hermitian_unsorted(
3689 SliceGroupMut::new(new_values),
3690 new_col_ptrs,
3691 new_row_indices,
3692 ghost::SparseColMatRef::new(A, N, N),
3693 ghost::PermutationRef::new(perm, N),
3694 in_side,
3695 out_side,
3696 true,
3697 stack,
3698 )
3699 }
3700 .into_inner()
3701 })
3702 }
3703
3704 #[doc(hidden)]
3705 pub fn ghost_adjoint_symbolic<'m, 'n, 'a, I: Index>(
3706 new_col_ptrs: &'a mut [I],
3707 new_row_indices: &'a mut [I],
3708 A: ghost::SymbolicSparseColMatRef<'m, 'n, '_, I>,
3709 stack: PodStack<'_>,
3710 ) -> ghost::SymbolicSparseColMatRef<'n, 'm, 'a, I> {
3711 let old_values = &*Symbolic::materialize(A.into_inner().row_indices().len());
3712 let new_values = Symbolic::materialize(new_row_indices.len());
3713 *ghost_adjoint(
3714 new_col_ptrs,
3715 new_row_indices,
3716 SliceGroupMut::<'_, Symbolic>::new(new_values),
3717 ghost::SparseColMatRef::new(
3718 SparseColMatRef::new(A.into_inner(), old_values),
3719 A.nrows(),
3720 A.ncols(),
3721 ),
3722 stack,
3723 )
3724 }
3725
3726 #[doc(hidden)]
3727 pub fn ghost_adjoint<'m, 'n, 'a, I: Index, E: ComplexField>(
3728 new_col_ptrs: &'a mut [I],
3729 new_row_indices: &'a mut [I],
3730 new_values: SliceGroupMut<'a, E>,
3731 A: ghost::SparseColMatRef<'m, 'n, '_, I, E>,
3732 stack: PodStack<'_>,
3733 ) -> ghost::SparseColMatMut<'n, 'm, 'a, I, E> {
3734 let M = A.nrows();
3735 let N = A.ncols();
3736 assert!(new_col_ptrs.len() == *M + 1);
3737
3738 let (col_count, _) = stack.make_raw::<I>(*M);
3739 let col_count = ghost::Array::from_mut(col_count, M);
3740 mem::fill_zero(col_count.as_mut());
3741
3742 for j in N.indices() {
3744 for i in A.row_indices_of_col(j) {
3745 col_count[i] += I::truncate(1);
3746 }
3747 }
3748
3749 new_col_ptrs[0] = I::truncate(0);
3750 for (j, [pj0, pj1]) in zip(
3752 M.indices(),
3753 windows2(Cell::as_slice_of_cells(Cell::from_mut(new_col_ptrs))),
3754 ) {
3755 let cj = &mut col_count[j];
3756 let pj = pj0.get();
3757 pj1.set(pj + *cj);
3759 *cj = pj;
3760 }
3761
3762 let new_row_indices = &mut new_row_indices[..new_col_ptrs[*M].zx()];
3763 let mut new_values = new_values.subslice(0..new_col_ptrs[*M].zx());
3764 let current_row_position = &mut *col_count;
3765 for j in N.indices() {
3767 let j_: ghost::Idx<'n, I> = j.truncate::<I>();
3768 for (i, val) in zip(
3769 A.row_indices_of_col(j),
3770 SliceGroup::<'_, E>::new(A.values_of_col(j)).into_ref_iter(),
3771 ) {
3772 let ci = &mut current_row_position[i];
3773
3774 unsafe {
3776 *new_row_indices.get_unchecked_mut(ci.zx()) = *j_;
3777 new_values.write_unchecked(ci.zx(), val.read().faer_conj())
3778 };
3779 *ci += I::truncate(1);
3780 }
3781 }
3782 debug_assert!(current_row_position.as_ref() == &new_col_ptrs[1..]);
3786
3787 ghost::SparseColMatMut::new(
3791 unsafe {
3792 SparseColMatMut::new(
3793 SymbolicSparseColMatRef::new_unchecked(
3794 *N,
3795 *M,
3796 new_col_ptrs,
3797 None,
3798 new_row_indices,
3799 ),
3800 new_values.into_inner(),
3801 )
3802 },
3803 N,
3804 M,
3805 )
3806 }
3807
3808 #[doc(hidden)]
3809 pub fn ghost_transpose<'m, 'n, 'a, I: Index, E: Entity>(
3810 new_col_ptrs: &'a mut [I],
3811 new_row_indices: &'a mut [I],
3812 new_values: SliceGroupMut<'a, E>,
3813 A: ghost::SparseColMatRef<'m, 'n, '_, I, E>,
3814 stack: PodStack<'_>,
3815 ) -> ghost::SparseColMatMut<'n, 'm, 'a, I, E> {
3816 let M = A.nrows();
3817 let N = A.ncols();
3818 assert!(new_col_ptrs.len() == *M + 1);
3819
3820 let (col_count, _) = stack.make_raw::<I>(*M);
3821 let col_count = ghost::Array::from_mut(col_count, M);
3822 mem::fill_zero(col_count.as_mut());
3823
3824 for j in N.indices() {
3826 for i in A.row_indices_of_col(j) {
3827 col_count[i] += I::truncate(1);
3828 }
3829 }
3830
3831 new_col_ptrs[0] = I::truncate(0);
3832 for (j, [pj0, pj1]) in zip(
3834 M.indices(),
3835 windows2(Cell::as_slice_of_cells(Cell::from_mut(new_col_ptrs))),
3836 ) {
3837 let cj = &mut col_count[j];
3838 let pj = pj0.get();
3839 pj1.set(pj + *cj);
3841 *cj = pj;
3842 }
3843
3844 let new_row_indices = &mut new_row_indices[..new_col_ptrs[*M].zx()];
3845 let mut new_values = new_values.subslice(0..new_col_ptrs[*M].zx());
3846 let current_row_position = &mut *col_count;
3847 for j in N.indices() {
3849 let j_: ghost::Idx<'n, I> = j.truncate::<I>();
3850 for (i, val) in zip(
3851 A.row_indices_of_col(j),
3852 SliceGroup::<'_, E>::new(A.values_of_col(j)).into_ref_iter(),
3853 ) {
3854 let ci = &mut current_row_position[i];
3855
3856 unsafe {
3858 *new_row_indices.get_unchecked_mut(ci.zx()) = *j_;
3859 new_values.write_unchecked(ci.zx(), val.read())
3860 };
3861 *ci += I::truncate(1);
3862 }
3863 }
3864 debug_assert!(current_row_position.as_ref() == &new_col_ptrs[1..]);
3868
3869 ghost::SparseColMatMut::new(
3873 unsafe {
3874 SparseColMatMut::new(
3875 SymbolicSparseColMatRef::new_unchecked(
3876 *N,
3877 *M,
3878 new_col_ptrs,
3879 None,
3880 new_row_indices,
3881 ),
3882 new_values.into_inner(),
3883 )
3884 },
3885 N,
3886 M,
3887 )
3888 }
3889
3890 pub fn transpose<'a, I: Index, E: Entity>(
3897 new_col_ptrs: &'a mut [I],
3898 new_row_indices: &'a mut [I],
3899 new_values: GroupFor<E, &'a mut [E::Unit]>,
3900 A: SparseColMatRef<'_, I, E>,
3901 stack: PodStack<'_>,
3902 ) -> SparseColMatMut<'a, I, E> {
3903 ghost::Size::with2(A.nrows(), A.ncols(), |M, N| {
3904 ghost_transpose(
3905 new_col_ptrs,
3906 new_row_indices,
3907 SliceGroupMut::new(new_values),
3908 ghost::SparseColMatRef::new(A, M, N),
3909 stack,
3910 )
3911 .into_inner()
3912 })
3913 }
3914
3915 pub fn adjoint<'a, I: Index, E: ComplexField>(
3922 new_col_ptrs: &'a mut [I],
3923 new_row_indices: &'a mut [I],
3924 new_values: GroupFor<E, &'a mut [E::Unit]>,
3925 A: SparseColMatRef<'_, I, E>,
3926 stack: PodStack<'_>,
3927 ) -> SparseColMatMut<'a, I, E> {
3928 ghost::Size::with2(A.nrows(), A.ncols(), |M, N| {
3929 ghost_adjoint(
3930 new_col_ptrs,
3931 new_row_indices,
3932 SliceGroupMut::new(new_values),
3933 ghost::SparseColMatRef::new(A, M, N),
3934 stack,
3935 )
3936 .into_inner()
3937 })
3938 }
3939
3940 pub fn adjoint_symbolic<'a, I: Index>(
3947 new_col_ptrs: &'a mut [I],
3948 new_row_indices: &'a mut [I],
3949 A: SymbolicSparseColMatRef<'_, I>,
3950 stack: PodStack<'_>,
3951 ) -> SymbolicSparseColMatRef<'a, I> {
3952 ghost::Size::with2(A.nrows(), A.ncols(), |M, N| {
3953 ghost_adjoint_symbolic(
3954 new_col_ptrs,
3955 new_row_indices,
3956 ghost::SymbolicSparseColMatRef::new(A, M, N),
3957 stack,
3958 )
3959 .into_inner()
3960 })
3961 }
3962}
3963
3964pub mod ops {
3966 use super::*;
3967 use crate::assert;
3968
3969 #[track_caller]
3975 pub fn binary_op<I: Index, E: Entity, LhsE: Entity, RhsE: Entity>(
3976 lhs: SparseColMatRef<'_, I, LhsE>,
3977 rhs: SparseColMatRef<'_, I, RhsE>,
3978 f: impl FnMut(LhsE, RhsE) -> E,
3979 ) -> Result<SparseColMat<I, E>, FaerError> {
3980 assert!(lhs.nrows() == rhs.nrows());
3981 assert!(lhs.ncols() == rhs.ncols());
3982 let mut f = f;
3983 let m = lhs.nrows();
3984 let n = lhs.ncols();
3985
3986 let mut col_ptrs = try_zeroed::<I>(n + 1)?;
3987
3988 let mut nnz = 0usize;
3989 for j in 0..n {
3990 let lhs = lhs.row_indices_of_col_raw(j);
3991 let rhs = rhs.row_indices_of_col_raw(j);
3992
3993 let mut lhs_pos = 0usize;
3994 let mut rhs_pos = 0usize;
3995 while lhs_pos < lhs.len() && rhs_pos < rhs.len() {
3996 let lhs = lhs[lhs_pos];
3997 let rhs = rhs[rhs_pos];
3998
3999 lhs_pos += (lhs <= rhs) as usize;
4000 rhs_pos += (rhs <= lhs) as usize;
4001 nnz += 1;
4002 }
4003 nnz += lhs.len() - lhs_pos;
4004 nnz += rhs.len() - rhs_pos;
4005 col_ptrs[j + 1] = I::truncate(nnz);
4006 }
4007
4008 if nnz > I::Signed::MAX.zx() {
4009 return Err(FaerError::IndexOverflow);
4010 }
4011
4012 let mut row_indices = try_zeroed(nnz)?;
4013 let mut values = VecGroup::<E>::new();
4014 values
4015 .try_reserve_exact(nnz)
4016 .map_err(|_| FaerError::OutOfMemory)?;
4017 values.resize(nnz, unsafe { core::mem::zeroed() });
4018
4019 let mut nnz = 0usize;
4020 for j in 0..n {
4021 let mut values = values.as_slice_mut();
4022 let lhs_values = SliceGroup::<LhsE>::new(lhs.values_of_col(j));
4023 let rhs_values = SliceGroup::<RhsE>::new(rhs.values_of_col(j));
4024 let lhs = lhs.row_indices_of_col_raw(j);
4025 let rhs = rhs.row_indices_of_col_raw(j);
4026
4027 let mut lhs_pos = 0usize;
4028 let mut rhs_pos = 0usize;
4029 while lhs_pos < lhs.len() && rhs_pos < rhs.len() {
4030 let lhs = lhs[lhs_pos];
4031 let rhs = rhs[rhs_pos];
4032
4033 match lhs.cmp(&rhs) {
4034 core::cmp::Ordering::Less => {
4035 row_indices[nnz] = lhs;
4036 values.write(
4037 nnz,
4038 f(lhs_values.read(lhs_pos), unsafe { core::mem::zeroed() }),
4039 );
4040 }
4041 core::cmp::Ordering::Equal => {
4042 row_indices[nnz] = lhs;
4043 values.write(nnz, f(lhs_values.read(lhs_pos), rhs_values.read(rhs_pos)));
4044 }
4045 core::cmp::Ordering::Greater => {
4046 row_indices[nnz] = rhs;
4047 values.write(
4048 nnz,
4049 f(unsafe { core::mem::zeroed() }, rhs_values.read(rhs_pos)),
4050 );
4051 }
4052 }
4053
4054 lhs_pos += (lhs <= rhs) as usize;
4055 rhs_pos += (rhs <= lhs) as usize;
4056 nnz += 1;
4057 }
4058 row_indices[nnz..nnz + lhs.len() - lhs_pos].copy_from_slice(&lhs[lhs_pos..]);
4059 for (mut dst, src) in values
4060 .rb_mut()
4061 .subslice(nnz..nnz + lhs.len() - lhs_pos)
4062 .into_mut_iter()
4063 .zip(lhs_values.subslice(lhs_pos..lhs.len()).into_ref_iter())
4064 {
4065 dst.write(f(src.read(), unsafe { core::mem::zeroed() }));
4066 }
4067 nnz += lhs.len() - lhs_pos;
4068
4069 row_indices[nnz..nnz + rhs.len() - rhs_pos].copy_from_slice(&rhs[rhs_pos..]);
4070 for (mut dst, src) in values
4071 .rb_mut()
4072 .subslice(nnz..nnz + rhs.len() - rhs_pos)
4073 .into_mut_iter()
4074 .zip(rhs_values.subslice(rhs_pos..rhs.len()).into_ref_iter())
4075 {
4076 dst.write(f(unsafe { core::mem::zeroed() }, src.read()));
4077 }
4078 nnz += rhs.len() - rhs_pos;
4079 }
4080
4081 Ok(SparseColMat::<I, E>::new(
4082 SymbolicSparseColMat::<I>::new_checked(m, n, col_ptrs, None, row_indices),
4083 values.into_inner(),
4084 ))
4085 }
4086
4087 #[track_caller]
4095 pub fn binary_op_assign_into<I: Index, E: Entity, SrcE: Entity>(
4096 dst: SparseColMatMut<'_, I, E>,
4097 src: SparseColMatRef<'_, I, SrcE>,
4098 f: impl FnMut(E, SrcE) -> E,
4099 ) {
4100 {
4101 assert!(dst.nrows() == src.nrows());
4102 assert!(dst.ncols() == src.ncols());
4103
4104 let n = dst.ncols();
4105 let mut dst = dst;
4106 let mut f = f;
4107 unsafe {
4108 assert!(f(core::mem::zeroed(), core::mem::zeroed()) == core::mem::zeroed());
4109 }
4110
4111 for j in 0..n {
4112 let (dst, dst_val) = dst.rb_mut().into_parts_mut();
4113
4114 let mut dst_val = SliceGroupMut::<E>::new(dst_val).subslice(dst.col_range(j));
4115 let src_val = SliceGroup::<SrcE>::new(src.values_of_col(j));
4116
4117 let dst = dst.row_indices_of_col_raw(j);
4118 let src = src.row_indices_of_col_raw(j);
4119
4120 let mut dst_pos = 0usize;
4121 let mut src_pos = 0usize;
4122
4123 while src_pos < src.len() {
4124 let src = src[src_pos];
4125
4126 if dst[dst_pos] < src {
4127 dst_val.write(
4128 dst_pos,
4129 f(dst_val.read(dst_pos), unsafe { core::mem::zeroed() }),
4130 );
4131 dst_pos += 1;
4132 continue;
4133 }
4134
4135 assert!(dst[dst_pos] == src);
4136
4137 dst_val.write(dst_pos, f(dst_val.read(dst_pos), src_val.read(src_pos)));
4138
4139 src_pos += 1;
4140 dst_pos += 1;
4141 }
4142 while dst_pos < dst.len() {
4143 dst_val.write(
4144 dst_pos,
4145 f(dst_val.read(dst_pos), unsafe { core::mem::zeroed() }),
4146 );
4147 dst_pos += 1;
4148 }
4149 }
4150 }
4151 }
4152
4153 #[track_caller]
4161 pub fn ternary_op_assign_into<I: Index, E: Entity, LhsE: Entity, RhsE: Entity>(
4162 dst: SparseColMatMut<'_, I, E>,
4163 lhs: SparseColMatRef<'_, I, LhsE>,
4164 rhs: SparseColMatRef<'_, I, RhsE>,
4165 f: impl FnMut(E, LhsE, RhsE) -> E,
4166 ) {
4167 {
4168 assert!(dst.nrows() == lhs.nrows());
4169 assert!(dst.ncols() == lhs.ncols());
4170 assert!(dst.nrows() == rhs.nrows());
4171 assert!(dst.ncols() == rhs.ncols());
4172
4173 let n = dst.ncols();
4174 let mut dst = dst;
4175 let mut f = f;
4176 unsafe {
4177 assert!(
4178 f(
4179 core::mem::zeroed(),
4180 core::mem::zeroed(),
4181 core::mem::zeroed()
4182 ) == core::mem::zeroed()
4183 );
4184 }
4185
4186 for j in 0..n {
4187 let (dst, dst_val) = dst.rb_mut().into_parts_mut();
4188
4189 let mut dst_val = SliceGroupMut::<E>::new(dst_val);
4190 let lhs_val = SliceGroup::<LhsE>::new(lhs.values_of_col(j));
4191 let rhs_val = SliceGroup::<RhsE>::new(rhs.values_of_col(j));
4192
4193 let dst = dst.row_indices_of_col_raw(j);
4194 let rhs = rhs.row_indices_of_col_raw(j);
4195 let lhs = lhs.row_indices_of_col_raw(j);
4196
4197 let mut dst_pos = 0usize;
4198 let mut lhs_pos = 0usize;
4199 let mut rhs_pos = 0usize;
4200
4201 while lhs_pos < lhs.len() && rhs_pos < rhs.len() {
4202 let lhs = lhs[lhs_pos];
4203 let rhs = rhs[rhs_pos];
4204
4205 if dst[dst_pos] < Ord::min(lhs, rhs) {
4206 dst_val.write(
4207 dst_pos,
4208 f(
4209 dst_val.read(dst_pos),
4210 unsafe { core::mem::zeroed() },
4211 unsafe { core::mem::zeroed() },
4212 ),
4213 );
4214 dst_pos += 1;
4215 continue;
4216 }
4217
4218 assert!(dst[dst_pos] == Ord::min(lhs, rhs));
4219
4220 match lhs.cmp(&rhs) {
4221 core::cmp::Ordering::Less => {
4222 dst_val.write(
4223 dst_pos,
4224 f(dst_val.read(dst_pos), lhs_val.read(lhs_pos), unsafe {
4225 core::mem::zeroed()
4226 }),
4227 );
4228 }
4229 core::cmp::Ordering::Equal => {
4230 dst_val.write(
4231 dst_pos,
4232 f(
4233 dst_val.read(dst_pos),
4234 lhs_val.read(lhs_pos),
4235 rhs_val.read(rhs_pos),
4236 ),
4237 );
4238 }
4239 core::cmp::Ordering::Greater => {
4240 dst_val.write(
4241 dst_pos,
4242 f(
4243 dst_val.read(dst_pos),
4244 unsafe { core::mem::zeroed() },
4245 rhs_val.read(rhs_pos),
4246 ),
4247 );
4248 }
4249 }
4250
4251 lhs_pos += (lhs <= rhs) as usize;
4252 rhs_pos += (rhs <= lhs) as usize;
4253 dst_pos += 1;
4254 }
4255 while lhs_pos < lhs.len() {
4256 let lhs = lhs[lhs_pos];
4257 if dst[dst_pos] < lhs {
4258 dst_val.write(
4259 dst_pos,
4260 f(
4261 dst_val.read(dst_pos),
4262 unsafe { core::mem::zeroed() },
4263 unsafe { core::mem::zeroed() },
4264 ),
4265 );
4266 dst_pos += 1;
4267 continue;
4268 }
4269 dst_val.write(
4270 dst_pos,
4271 f(dst_val.read(dst_pos), lhs_val.read(lhs_pos), unsafe {
4272 core::mem::zeroed()
4273 }),
4274 );
4275 lhs_pos += 1;
4276 dst_pos += 1;
4277 }
4278 while rhs_pos < rhs.len() {
4279 let rhs = rhs[rhs_pos];
4280 if dst[dst_pos] < rhs {
4281 dst_val.write(
4282 dst_pos,
4283 f(
4284 dst_val.read(dst_pos),
4285 unsafe { core::mem::zeroed() },
4286 unsafe { core::mem::zeroed() },
4287 ),
4288 );
4289 dst_pos += 1;
4290 continue;
4291 }
4292 dst_val.write(
4293 dst_pos,
4294 f(
4295 dst_val.read(dst_pos),
4296 unsafe { core::mem::zeroed() },
4297 rhs_val.read(rhs_pos),
4298 ),
4299 );
4300 rhs_pos += 1;
4301 dst_pos += 1;
4302 }
4303 while rhs_pos < rhs.len() {
4304 let rhs = rhs[rhs_pos];
4305 dst_pos += dst[dst_pos..].binary_search(&rhs).unwrap();
4306 dst_val.write(
4307 dst_pos,
4308 f(
4309 dst_val.read(dst_pos),
4310 unsafe { core::mem::zeroed() },
4311 rhs_val.read(rhs_pos),
4312 ),
4313 );
4314 rhs_pos += 1;
4315 }
4316 }
4317 }
4318 }
4319
4320 #[track_caller]
4325 #[inline]
4326 pub fn union_symbolic<I: Index>(
4327 lhs: SymbolicSparseColMatRef<'_, I>,
4328 rhs: SymbolicSparseColMatRef<'_, I>,
4329 ) -> Result<SymbolicSparseColMat<I>, FaerError> {
4330 Ok(binary_op(
4331 SparseColMatRef::<I, Symbolic>::new(lhs, Symbolic::materialize(lhs.compute_nnz())),
4332 SparseColMatRef::<I, Symbolic>::new(rhs, Symbolic::materialize(rhs.compute_nnz())),
4333 #[inline(always)]
4334 |_, _| Symbolic,
4335 )?
4336 .into_parts()
4337 .0)
4338 }
4339
4340 #[track_caller]
4345 #[inline]
4346 pub fn add<
4347 I: Index,
4348 E: ComplexField,
4349 LhsE: Conjugate<Canonical = E>,
4350 RhsE: Conjugate<Canonical = E>,
4351 >(
4352 lhs: SparseColMatRef<'_, I, LhsE>,
4353 rhs: SparseColMatRef<'_, I, RhsE>,
4354 ) -> Result<SparseColMat<I, E>, FaerError> {
4355 binary_op(lhs, rhs, |lhs, rhs| {
4356 lhs.canonicalize().faer_add(rhs.canonicalize())
4357 })
4358 }
4359
4360 #[track_caller]
4365 #[inline]
4366 pub fn sub<
4367 I: Index,
4368 LhsE: Conjugate<Canonical = E>,
4369 RhsE: Conjugate<Canonical = E>,
4370 E: ComplexField,
4371 >(
4372 lhs: SparseColMatRef<'_, I, LhsE>,
4373 rhs: SparseColMatRef<'_, I, RhsE>,
4374 ) -> Result<SparseColMat<I, E>, FaerError> {
4375 binary_op(lhs, rhs, |lhs, rhs| {
4376 lhs.canonicalize().faer_sub(rhs.canonicalize())
4377 })
4378 }
4379
4380 pub fn add_assign<I: Index, E: ComplexField, RhsE: Conjugate<Canonical = E>>(
4387 dst: SparseColMatMut<'_, I, E>,
4388 rhs: SparseColMatRef<'_, I, RhsE>,
4389 ) {
4390 binary_op_assign_into(dst, rhs, |dst, rhs| dst.faer_add(rhs.canonicalize()))
4391 }
4392
4393 pub fn sub_assign<I: Index, E: ComplexField, RhsE: Conjugate<Canonical = E>>(
4400 dst: SparseColMatMut<'_, I, E>,
4401 rhs: SparseColMatRef<'_, I, RhsE>,
4402 ) {
4403 binary_op_assign_into(dst, rhs, |dst, rhs| dst.faer_sub(rhs.canonicalize()))
4404 }
4405
4406 #[track_caller]
4413 #[inline]
4414 pub fn add_into<
4415 I: Index,
4416 E: ComplexField,
4417 LhsE: Conjugate<Canonical = E>,
4418 RhsE: Conjugate<Canonical = E>,
4419 >(
4420 dst: SparseColMatMut<'_, I, E>,
4421 lhs: SparseColMatRef<'_, I, LhsE>,
4422 rhs: SparseColMatRef<'_, I, RhsE>,
4423 ) {
4424 ternary_op_assign_into(dst, lhs, rhs, |_, lhs, rhs| {
4425 lhs.canonicalize().faer_add(rhs.canonicalize())
4426 })
4427 }
4428
4429 #[track_caller]
4436 #[inline]
4437 pub fn sub_into<
4438 I: Index,
4439 E: ComplexField,
4440 LhsE: Conjugate<Canonical = E>,
4441 RhsE: Conjugate<Canonical = E>,
4442 >(
4443 dst: SparseColMatMut<'_, I, E>,
4444 lhs: SparseColMatRef<'_, I, LhsE>,
4445 rhs: SparseColMatRef<'_, I, RhsE>,
4446 ) {
4447 ternary_op_assign_into(dst, lhs, rhs, |_, lhs, rhs| {
4448 lhs.canonicalize().faer_sub(rhs.canonicalize())
4449 })
4450 }
4451}
4452
4453impl<'a, I: Index, E: Entity> SparseColMatRef<'a, I, E> {
4454 #[track_caller]
4461 pub fn get(self, row: usize, col: usize) -> Option<GroupFor<E, &'a E::Unit>> {
4462 assert!(row < self.nrows());
4463 assert!(col < self.ncols());
4464
4465 let Ok(pos) = self
4466 .row_indices_of_col_raw(col)
4467 .binary_search(&I::truncate(row))
4468 else {
4469 return None;
4470 };
4471
4472 Some(E::faer_map(self.values_of_col(col), |slice| &slice[pos]))
4473 }
4474}
4475
4476impl<'a, I: Index, E: Entity> SparseColMatMut<'a, I, E> {
4477 #[track_caller]
4484 pub fn get(self, row: usize, col: usize) -> Option<GroupFor<E, &'a E::Unit>> {
4485 self.into_const().get(row, col)
4486 }
4487
4488 #[track_caller]
4495 pub fn get_mut(self, row: usize, col: usize) -> Option<GroupFor<E, &'a mut E::Unit>> {
4496 assert!(row < self.nrows());
4497 assert!(col < self.ncols());
4498
4499 let Ok(pos) = self
4500 .row_indices_of_col_raw(col)
4501 .binary_search(&I::truncate(row))
4502 else {
4503 return None;
4504 };
4505
4506 Some(E::faer_map(self.values_of_col_mut(col), |slice| {
4507 &mut slice[pos]
4508 }))
4509 }
4510}
4511
4512impl<I: Index, E: Entity> SparseColMat<I, E> {
4513 #[track_caller]
4520 pub fn get(&self, row: usize, col: usize) -> Option<GroupFor<E, &'_ E::Unit>> {
4521 self.as_ref().get(row, col)
4522 }
4523
4524 #[track_caller]
4531 pub fn get_mut(&mut self, row: usize, col: usize) -> Option<GroupFor<E, &'_ mut E::Unit>> {
4532 self.as_mut().get_mut(row, col)
4533 }
4534}
4535
4536impl<'a, I: Index, E: Entity> SparseRowMatRef<'a, I, E> {
4537 #[track_caller]
4544 pub fn get(self, row: usize, col: usize) -> Option<GroupFor<E, &'a E::Unit>> {
4545 assert!(row < self.nrows());
4546 assert!(col < self.ncols());
4547
4548 let Ok(pos) = self
4549 .col_indices_of_row_raw(row)
4550 .binary_search(&I::truncate(col))
4551 else {
4552 return None;
4553 };
4554
4555 Some(E::faer_map(self.values_of_row(row), |slice| &slice[pos]))
4556 }
4557}
4558
4559impl<'a, I: Index, E: Entity> SparseRowMatMut<'a, I, E> {
4560 #[track_caller]
4567 pub fn get_mut(self, row: usize, col: usize) -> Option<GroupFor<E, &'a mut E::Unit>> {
4568 assert!(row < self.nrows());
4569 assert!(col < self.ncols());
4570
4571 let Ok(pos) = self
4572 .col_indices_of_row_raw(row)
4573 .binary_search(&I::truncate(col))
4574 else {
4575 return None;
4576 };
4577
4578 Some(E::faer_map(self.values_of_row_mut(row), |slice| {
4579 &mut slice[pos]
4580 }))
4581 }
4582}
4583
4584impl<I: Index, E: Entity> SparseRowMat<I, E> {
4585 #[track_caller]
4592 pub fn get(&self, row: usize, col: usize) -> Option<GroupFor<E, &'_ E::Unit>> {
4593 self.as_ref().get(row, col)
4594 }
4595
4596 #[track_caller]
4603 pub fn get_mut(&mut self, row: usize, col: usize) -> Option<GroupFor<E, &'_ mut E::Unit>> {
4604 self.as_mut().get_mut(row, col)
4605 }
4606}
4607
4608impl<I: Index, E: SimpleEntity> core::ops::Index<(usize, usize)> for SparseColMatRef<'_, I, E> {
4609 type Output = E;
4610
4611 #[track_caller]
4612 fn index(&self, (row, col): (usize, usize)) -> &Self::Output {
4613 self.get(row, col).unwrap()
4614 }
4615}
4616
4617impl<I: Index, E: SimpleEntity> core::ops::Index<(usize, usize)> for SparseRowMatRef<'_, I, E> {
4618 type Output = E;
4619
4620 #[track_caller]
4621 fn index(&self, (row, col): (usize, usize)) -> &Self::Output {
4622 self.get(row, col).unwrap()
4623 }
4624}
4625
4626impl<I: Index, E: SimpleEntity> core::ops::Index<(usize, usize)> for SparseColMatMut<'_, I, E> {
4627 type Output = E;
4628
4629 #[track_caller]
4630 fn index(&self, (row, col): (usize, usize)) -> &Self::Output {
4631 self.rb().get(row, col).unwrap()
4632 }
4633}
4634
4635impl<I: Index, E: SimpleEntity> core::ops::Index<(usize, usize)> for SparseRowMatMut<'_, I, E> {
4636 type Output = E;
4637
4638 #[track_caller]
4639 fn index(&self, (row, col): (usize, usize)) -> &Self::Output {
4640 self.rb().get(row, col).unwrap()
4641 }
4642}
4643
4644impl<I: Index, E: SimpleEntity> core::ops::IndexMut<(usize, usize)> for SparseColMatMut<'_, I, E> {
4645 #[track_caller]
4646 fn index_mut(&mut self, (row, col): (usize, usize)) -> &mut Self::Output {
4647 self.rb_mut().get_mut(row, col).unwrap()
4648 }
4649}
4650
4651impl<I: Index, E: SimpleEntity> core::ops::IndexMut<(usize, usize)> for SparseRowMatMut<'_, I, E> {
4652 #[track_caller]
4653 fn index_mut(&mut self, (row, col): (usize, usize)) -> &mut Self::Output {
4654 self.rb_mut().get_mut(row, col).unwrap()
4655 }
4656}
4657
4658impl<I: Index, E: SimpleEntity> core::ops::Index<(usize, usize)> for SparseColMat<I, E> {
4659 type Output = E;
4660
4661 #[track_caller]
4662 fn index(&self, (row, col): (usize, usize)) -> &Self::Output {
4663 self.as_ref().get(row, col).unwrap()
4664 }
4665}
4666
4667impl<I: Index, E: SimpleEntity> core::ops::Index<(usize, usize)> for SparseRowMat<I, E> {
4668 type Output = E;
4669
4670 #[track_caller]
4671 fn index(&self, (row, col): (usize, usize)) -> &Self::Output {
4672 self.as_ref().get(row, col).unwrap()
4673 }
4674}
4675
4676impl<I: Index, E: SimpleEntity> core::ops::IndexMut<(usize, usize)> for SparseColMat<I, E> {
4677 #[track_caller]
4678 fn index_mut(&mut self, (row, col): (usize, usize)) -> &mut Self::Output {
4679 self.as_mut().get_mut(row, col).unwrap()
4680 }
4681}
4682
4683impl<I: Index, E: SimpleEntity> core::ops::IndexMut<(usize, usize)> for SparseRowMat<I, E> {
4684 #[track_caller]
4685 fn index_mut(&mut self, (row, col): (usize, usize)) -> &mut Self::Output {
4686 self.as_mut().get_mut(row, col).unwrap()
4687 }
4688}
4689
4690pub mod mul {
4692 use super::*;
4699 use crate::{
4700 assert,
4701 constrained::{self, Size},
4702 };
4703
4704 #[track_caller]
4710 pub fn sparse_dense_matmul<
4711 I: Index,
4712 E: ComplexField,
4713 LhsE: Conjugate<Canonical = E>,
4714 RhsE: Conjugate<Canonical = E>,
4715 >(
4716 acc: MatMut<'_, E>,
4717 lhs: SparseColMatRef<'_, I, LhsE>,
4718 rhs: MatRef<'_, RhsE>,
4719 alpha: Option<E>,
4720 beta: E,
4721 parallelism: Parallelism,
4722 ) {
4723 assert!(all(
4724 acc.nrows() == lhs.nrows(),
4725 acc.ncols() == rhs.ncols(),
4726 lhs.ncols() == rhs.nrows(),
4727 ));
4728
4729 let _ = parallelism;
4730 let m = acc.nrows();
4731 let n = acc.ncols();
4732 let k = lhs.ncols();
4733
4734 let mut acc = acc;
4735
4736 match alpha {
4737 Some(alpha) => {
4738 if alpha != E::faer_one() {
4739 zipped!(acc.rb_mut())
4740 .for_each(|unzipped!(mut dst)| dst.write(dst.read().faer_mul(alpha)))
4741 }
4742 }
4743 None => acc.fill_zero(),
4744 }
4745
4746 Size::with2(m, n, |m, n| {
4747 Size::with(k, |k| {
4748 let mut acc = constrained::MatMut::new(acc, m, n);
4749 let lhs = constrained::sparse::SparseColMatRef::new(lhs, m, k);
4750 let rhs = constrained::MatRef::new(rhs, k, n);
4751
4752 for j in n.indices() {
4753 for depth in k.indices() {
4754 let rhs_kj = rhs.read(depth, j).canonicalize().faer_mul(beta);
4755 for (i, lhs_ik) in zip(
4756 lhs.row_indices_of_col(depth),
4757 SliceGroup::<'_, LhsE>::new(lhs.values_of_col(depth)).into_ref_iter(),
4758 ) {
4759 acc.write(
4760 i,
4761 j,
4762 acc.read(i, j)
4763 .faer_add(lhs_ik.read().canonicalize().faer_mul(rhs_kj)),
4764 );
4765 }
4766 }
4767 }
4768 });
4769 });
4770 }
4771
4772 #[track_caller]
4778 pub fn dense_sparse_matmul<
4779 I: Index,
4780 E: ComplexField,
4781 LhsE: Conjugate<Canonical = E>,
4782 RhsE: Conjugate<Canonical = E>,
4783 >(
4784 acc: MatMut<'_, E>,
4785 lhs: MatRef<'_, LhsE>,
4786 rhs: SparseColMatRef<'_, I, RhsE>,
4787 alpha: Option<E>,
4788 beta: E,
4789 parallelism: Parallelism,
4790 ) {
4791 assert!(all(
4792 acc.nrows() == lhs.nrows(),
4793 acc.ncols() == rhs.ncols(),
4794 lhs.ncols() == rhs.nrows(),
4795 ));
4796
4797 let _ = parallelism;
4798 let m = acc.nrows();
4799 let n = acc.ncols();
4800 let k = lhs.ncols();
4801
4802 let mut acc = acc;
4803
4804 match alpha {
4805 Some(alpha) => {
4806 if alpha != E::faer_one() {
4807 zipped!(acc.rb_mut())
4808 .for_each(|unzipped!(mut dst)| dst.write(dst.read().faer_mul(alpha)))
4809 }
4810 }
4811 None => acc.fill_zero(),
4812 }
4813
4814 Size::with2(m, n, |m, n| {
4815 Size::with(k, |k| {
4816 let mut acc = constrained::MatMut::new(acc, m, n);
4817 let lhs = constrained::MatRef::new(lhs, m, k);
4818 let rhs = constrained::sparse::SparseColMatRef::new(rhs, k, n);
4819
4820 for i in m.indices() {
4821 for j in n.indices() {
4822 let mut acc_ij = E::faer_zero();
4823 for (depth, rhs_kj) in zip(
4824 rhs.row_indices_of_col(j),
4825 SliceGroup::<'_, RhsE>::new(rhs.values_of_col(j)).into_ref_iter(),
4826 ) {
4827 let lhs_ik = lhs.read(i, depth);
4828 acc_ij = acc_ij.faer_add(
4829 lhs_ik.canonicalize().faer_mul(rhs_kj.read().canonicalize()),
4830 );
4831 }
4832
4833 acc.write(i, j, acc.read(i, j).faer_add(beta.faer_mul(acc_ij)));
4834 }
4835 }
4836 });
4837 });
4838 }
4839}
4840
4841#[cfg(feature = "std")]
4842#[cfg_attr(docsrs, doc(cfg(feature = "std")))]
4843impl<I: Index, E: Entity> matrixcompare_core::Matrix<E> for SparseColMatRef<'_, I, E> {
4844 #[inline]
4845 fn rows(&self) -> usize {
4846 self.nrows()
4847 }
4848 #[inline]
4849 fn cols(&self) -> usize {
4850 self.ncols()
4851 }
4852 #[inline]
4853 fn access(&self) -> matrixcompare_core::Access<'_, E> {
4854 matrixcompare_core::Access::Sparse(self)
4855 }
4856}
4857
4858#[cfg(feature = "std")]
4859#[cfg_attr(docsrs, doc(cfg(feature = "std")))]
4860impl<I: Index, E: Entity> matrixcompare_core::SparseAccess<E> for SparseColMatRef<'_, I, E> {
4861 #[inline]
4862 fn nnz(&self) -> usize {
4863 self.compute_nnz()
4864 }
4865
4866 #[inline]
4867 fn fetch_triplets(&self) -> Vec<(usize, usize, E)> {
4868 let mut triplets = Vec::new();
4869 for j in 0..self.ncols() {
4870 for (i, val) in self
4871 .row_indices_of_col(j)
4872 .zip(SliceGroup::<'_, E>::new(self.values_of_col(j)).into_ref_iter())
4873 {
4874 triplets.push((i, j, val.read()))
4875 }
4876 }
4877 triplets
4878 }
4879}
4880
4881#[cfg(feature = "std")]
4882#[cfg_attr(docsrs, doc(cfg(feature = "std")))]
4883impl<I: Index, E: Entity> matrixcompare_core::Matrix<E> for SparseRowMatRef<'_, I, E> {
4884 #[inline]
4885 fn rows(&self) -> usize {
4886 self.nrows()
4887 }
4888 #[inline]
4889 fn cols(&self) -> usize {
4890 self.ncols()
4891 }
4892 #[inline]
4893 fn access(&self) -> matrixcompare_core::Access<'_, E> {
4894 matrixcompare_core::Access::Sparse(self)
4895 }
4896}
4897
4898#[cfg(feature = "std")]
4899#[cfg_attr(docsrs, doc(cfg(feature = "std")))]
4900impl<I: Index, E: Entity> matrixcompare_core::SparseAccess<E> for SparseRowMatRef<'_, I, E> {
4901 #[inline]
4902 fn nnz(&self) -> usize {
4903 self.compute_nnz()
4904 }
4905
4906 #[inline]
4907 fn fetch_triplets(&self) -> Vec<(usize, usize, E)> {
4908 let mut triplets = Vec::new();
4909 for i in 0..self.nrows() {
4910 for (j, val) in self
4911 .col_indices_of_row(i)
4912 .zip(SliceGroup::<'_, E>::new(self.values_of_row(i)).into_ref_iter())
4913 {
4914 triplets.push((i, j, val.read()))
4915 }
4916 }
4917 triplets
4918 }
4919}
4920
4921#[cfg(feature = "std")]
4922#[cfg_attr(docsrs, doc(cfg(feature = "std")))]
4923impl<I: Index, E: Entity> matrixcompare_core::Matrix<E> for SparseColMatMut<'_, I, E> {
4924 #[inline]
4925 fn rows(&self) -> usize {
4926 self.nrows()
4927 }
4928 #[inline]
4929 fn cols(&self) -> usize {
4930 self.ncols()
4931 }
4932 #[inline]
4933 fn access(&self) -> matrixcompare_core::Access<'_, E> {
4934 matrixcompare_core::Access::Sparse(self)
4935 }
4936}
4937
4938#[cfg(feature = "std")]
4939#[cfg_attr(docsrs, doc(cfg(feature = "std")))]
4940impl<I: Index, E: Entity> matrixcompare_core::SparseAccess<E> for SparseColMatMut<'_, I, E> {
4941 #[inline]
4942 fn nnz(&self) -> usize {
4943 self.compute_nnz()
4944 }
4945
4946 #[inline]
4947 fn fetch_triplets(&self) -> Vec<(usize, usize, E)> {
4948 self.rb().fetch_triplets()
4949 }
4950}
4951
4952#[cfg(feature = "std")]
4953#[cfg_attr(docsrs, doc(cfg(feature = "std")))]
4954impl<I: Index, E: Entity> matrixcompare_core::Matrix<E> for SparseColMat<I, E> {
4955 #[inline]
4956 fn rows(&self) -> usize {
4957 self.nrows()
4958 }
4959 #[inline]
4960 fn cols(&self) -> usize {
4961 self.ncols()
4962 }
4963 #[inline]
4964 fn access(&self) -> matrixcompare_core::Access<'_, E> {
4965 matrixcompare_core::Access::Sparse(self)
4966 }
4967}
4968
4969#[cfg(feature = "std")]
4970#[cfg_attr(docsrs, doc(cfg(feature = "std")))]
4971impl<I: Index, E: Entity> matrixcompare_core::SparseAccess<E> for SparseColMat<I, E> {
4972 #[inline]
4973 fn nnz(&self) -> usize {
4974 self.compute_nnz()
4975 }
4976
4977 #[inline]
4978 fn fetch_triplets(&self) -> Vec<(usize, usize, E)> {
4979 self.as_ref().fetch_triplets()
4980 }
4981}
4982
4983#[cfg(feature = "std")]
4984#[cfg_attr(docsrs, doc(cfg(feature = "std")))]
4985impl<I: Index, E: Entity> matrixcompare_core::Matrix<E> for SparseRowMatMut<'_, I, E> {
4986 #[inline]
4987 fn rows(&self) -> usize {
4988 self.nrows()
4989 }
4990 #[inline]
4991 fn cols(&self) -> usize {
4992 self.ncols()
4993 }
4994 #[inline]
4995 fn access(&self) -> matrixcompare_core::Access<'_, E> {
4996 matrixcompare_core::Access::Sparse(self)
4997 }
4998}
4999
5000#[cfg(feature = "std")]
5001#[cfg_attr(docsrs, doc(cfg(feature = "std")))]
5002impl<I: Index, E: Entity> matrixcompare_core::SparseAccess<E> for SparseRowMatMut<'_, I, E> {
5003 #[inline]
5004 fn nnz(&self) -> usize {
5005 self.compute_nnz()
5006 }
5007
5008 #[inline]
5009 fn fetch_triplets(&self) -> Vec<(usize, usize, E)> {
5010 self.rb().fetch_triplets()
5011 }
5012}
5013
5014#[cfg(feature = "std")]
5015#[cfg_attr(docsrs, doc(cfg(feature = "std")))]
5016impl<I: Index, E: Entity> matrixcompare_core::Matrix<E> for SparseRowMat<I, E> {
5017 #[inline]
5018 fn rows(&self) -> usize {
5019 self.nrows()
5020 }
5021 #[inline]
5022 fn cols(&self) -> usize {
5023 self.ncols()
5024 }
5025 #[inline]
5026 fn access(&self) -> matrixcompare_core::Access<'_, E> {
5027 matrixcompare_core::Access::Sparse(self)
5028 }
5029}
5030
5031#[cfg(feature = "std")]
5032#[cfg_attr(docsrs, doc(cfg(feature = "std")))]
5033impl<I: Index, E: Entity> matrixcompare_core::SparseAccess<E> for SparseRowMat<I, E> {
5034 #[inline]
5035 fn nnz(&self) -> usize {
5036 self.compute_nnz()
5037 }
5038
5039 #[inline]
5040 fn fetch_triplets(&self) -> Vec<(usize, usize, E)> {
5041 self.as_ref().fetch_triplets()
5042 }
5043}
5044
5045#[cfg(test)]
5046mod tests {
5047 use super::*;
5048 use crate::assert;
5049
5050 #[test]
5051 fn test_from_indices() {
5052 let nrows = 5;
5053 let ncols = 4;
5054
5055 let indices = &[(0, 0), (1, 2), (0, 0), (1, 1), (0, 1), (3, 3), (3, 3usize)];
5056 let values = &[1.0, 2.0, 3.0, 4.0, 5.0, 6.0, 7.0f64];
5057
5058 let triplets = &[
5059 (0, 0, 1.0),
5060 (1, 2, 2.0),
5061 (0, 0, 3.0),
5062 (1, 1, 4.0),
5063 (0, 1, 5.0),
5064 (3, 3, 6.0),
5065 (3, 3usize, 7.0),
5066 ];
5067
5068 {
5069 let mat = SymbolicSparseColMat::try_new_from_indices(nrows, ncols, indices);
5070 assert!(mat.is_ok());
5071
5072 let (mat, order) = mat.unwrap();
5073 assert!(mat.nrows() == nrows);
5074 assert!(mat.ncols() == ncols);
5075 assert!(mat.col_ptrs() == &[0, 1, 3, 4, 5]);
5076 assert!(mat.nnz_per_col() == None);
5077 assert!(mat.row_indices() == &[0, 0, 1, 1, 3]);
5078
5079 let mat =
5080 SparseColMat::<_, f64>::new_from_order_and_values(mat, &order, values).unwrap();
5081 assert!(mat.as_ref().values() == &[1.0 + 3.0, 5.0, 4.0, 2.0, 6.0 + 7.0]);
5082 }
5083
5084 {
5085 let mat = SparseColMat::try_new_from_triplets(nrows, ncols, triplets);
5086 assert!(mat.is_ok());
5087 let mat = mat.unwrap();
5088
5089 assert!(mat.nrows() == nrows);
5090 assert!(mat.ncols() == ncols);
5091 assert!(mat.col_ptrs() == &[0, 1, 3, 4, 5]);
5092 assert!(mat.nnz_per_col() == None);
5093 assert!(mat.row_indices() == &[0, 0, 1, 1, 3]);
5094 assert!(mat.values() == &[1.0 + 3.0, 5.0, 4.0, 2.0, 6.0 + 7.0]);
5095 }
5096
5097 {
5098 let mat = SymbolicSparseRowMat::try_new_from_indices(nrows, ncols, indices);
5099 assert!(mat.is_ok());
5100
5101 let (mat, order) = mat.unwrap();
5102 assert!(mat.nrows() == nrows);
5103 assert!(mat.ncols() == ncols);
5104 assert!(mat.row_ptrs() == &[0, 2, 4, 4, 5, 5]);
5105 assert!(mat.nnz_per_row() == None);
5106 assert!(mat.col_indices() == &[0, 1, 1, 2, 3]);
5107
5108 let mat =
5109 SparseRowMat::<_, f64>::new_from_order_and_values(mat, &order, values).unwrap();
5110 assert!(mat.values() == &[1.0 + 3.0, 5.0, 4.0, 2.0, 6.0 + 7.0]);
5111 }
5112 {
5113 let mat = SparseRowMat::try_new_from_triplets(nrows, ncols, triplets);
5114 assert!(mat.is_ok());
5115
5116 let mat = mat.unwrap();
5117 assert!(mat.nrows() == nrows);
5118 assert!(mat.ncols() == ncols);
5119 assert!(mat.row_ptrs() == &[0, 2, 4, 4, 5, 5]);
5120 assert!(mat.nnz_per_row() == None);
5121 assert!(mat.col_indices() == &[0, 1, 1, 2, 3]);
5122 assert!(mat.as_ref().values() == &[1.0 + 3.0, 5.0, 4.0, 2.0, 6.0 + 7.0]);
5123 }
5124 }
5125
5126 #[test]
5127 fn test_from_nonnegative_indices() {
5128 let nrows = 5;
5129 let ncols = 4;
5130
5131 let indices = &[
5132 (0, 0),
5133 (1, 2),
5134 (0, 0),
5135 (1, 1),
5136 (0, 1),
5137 (-1, 2),
5138 (-2, 1),
5139 (-3, -4),
5140 (3, 3),
5141 (3, 3isize),
5142 ];
5143 let values = &[
5144 1.0,
5145 2.0,
5146 3.0,
5147 4.0,
5148 5.0,
5149 f64::NAN,
5150 f64::NAN,
5151 f64::NAN,
5152 6.0,
5153 7.0f64,
5154 ];
5155
5156 let triplets = &[
5157 (0, 0, 1.0),
5158 (1, 2, 2.0),
5159 (0, 0, 3.0),
5160 (1, 1, 4.0),
5161 (0, 1, 5.0),
5162 (-1, 2, f64::NAN),
5163 (-2, 1, f64::NAN),
5164 (-3, -4, f64::NAN),
5165 (3, 3, 6.0),
5166 (3, 3isize, 7.0),
5167 ];
5168
5169 {
5170 let mat = SymbolicSparseColMat::<usize>::try_new_from_nonnegative_indices(
5171 nrows, ncols, indices,
5172 );
5173 assert!(mat.is_ok());
5174
5175 let (mat, order) = mat.unwrap();
5176 assert!(mat.nrows() == nrows);
5177 assert!(mat.ncols() == ncols);
5178 assert!(mat.col_ptrs() == &[0, 1, 3, 4, 5]);
5179 assert!(mat.nnz_per_col() == None);
5180 assert!(mat.row_indices() == &[0, 0, 1, 1, 3]);
5181
5182 let mat =
5183 SparseColMat::<_, f64>::new_from_order_and_values(mat, &order, values).unwrap();
5184 assert!(mat.as_ref().values() == &[1.0 + 3.0, 5.0, 4.0, 2.0, 6.0 + 7.0]);
5185 }
5186
5187 {
5188 let mat =
5189 SparseColMat::<usize, _>::try_new_from_nonnegative_triplets(nrows, ncols, triplets);
5190 assert!(mat.is_ok());
5191 let mat = mat.unwrap();
5192
5193 assert!(mat.nrows() == nrows);
5194 assert!(mat.ncols() == ncols);
5195 assert!(mat.col_ptrs() == &[0, 1, 3, 4, 5]);
5196 assert!(mat.nnz_per_col() == None);
5197 assert!(mat.row_indices() == &[0, 0, 1, 1, 3]);
5198 assert!(mat.values() == &[1.0 + 3.0, 5.0, 4.0, 2.0, 6.0 + 7.0]);
5199 }
5200
5201 {
5202 let mat = SymbolicSparseRowMat::<usize>::try_new_from_nonnegative_indices(
5203 nrows, ncols, indices,
5204 );
5205 assert!(mat.is_ok());
5206
5207 let (mat, order) = mat.unwrap();
5208 assert!(mat.nrows() == nrows);
5209 assert!(mat.ncols() == ncols);
5210 assert!(mat.row_ptrs() == &[0, 2, 4, 4, 5, 5]);
5211 assert!(mat.nnz_per_row() == None);
5212 assert!(mat.col_indices() == &[0, 1, 1, 2, 3]);
5213
5214 let mat =
5215 SparseRowMat::<_, f64>::new_from_order_and_values(mat, &order, values).unwrap();
5216 assert!(mat.values() == &[1.0 + 3.0, 5.0, 4.0, 2.0, 6.0 + 7.0]);
5217 }
5218 {
5219 let mat =
5220 SparseRowMat::<usize, _>::try_new_from_nonnegative_triplets(nrows, ncols, triplets);
5221 assert!(mat.is_ok());
5222
5223 let mat = mat.unwrap();
5224 assert!(mat.nrows() == nrows);
5225 assert!(mat.ncols() == ncols);
5226 assert!(mat.row_ptrs() == &[0, 2, 4, 4, 5, 5]);
5227 assert!(mat.nnz_per_row() == None);
5228 assert!(mat.col_indices() == &[0, 1, 1, 2, 3]);
5229 assert!(mat.as_ref().values() == &[1.0 + 3.0, 5.0, 4.0, 2.0, 6.0 + 7.0]);
5230 }
5231 {
5232 let order = SymbolicSparseRowMat::<usize>::try_new_from_nonnegative_indices(
5233 nrows, ncols, indices,
5234 )
5235 .unwrap()
5236 .1;
5237
5238 let new_values = &mut [f64::NAN; 5];
5239 let mut mat = SparseRowMatMut::<'_, usize, f64>::new(
5240 SymbolicSparseRowMatRef::new_checked(
5241 nrows,
5242 ncols,
5243 &[0, 2, 4, 4, 5, 5],
5244 None,
5245 &[0, 1, 1, 2, 3],
5246 ),
5247 new_values,
5248 );
5249 mat.fill_from_order_and_values(&order, values, FillMode::Replace);
5250
5251 assert!(&*new_values == &[1.0 + 3.0, 5.0, 4.0, 2.0, 6.0 + 7.0]);
5252 }
5253 }
5254
5255 #[test]
5256 fn test_from_indices_oob_row() {
5257 let nrows = 5;
5258 let ncols = 4;
5259
5260 let indices = &[
5261 (0, 0),
5262 (1, 2),
5263 (0, 0),
5264 (1, 1),
5265 (0, 1),
5266 (3, 3),
5267 (3, 3),
5268 (5, 3usize),
5269 ];
5270 let err = SymbolicSparseColMat::try_new_from_indices(nrows, ncols, indices);
5271 assert!(err.is_err());
5272 let err = err.unwrap_err();
5273 assert!(err == CreationError::OutOfBounds { row: 5, col: 3 });
5274 }
5275
5276 #[test]
5277 fn test_from_indices_oob_col() {
5278 let nrows = 5;
5279 let ncols = 4;
5280
5281 let indices = &[
5282 (0, 0),
5283 (1, 2),
5284 (0, 0),
5285 (1, 1),
5286 (0, 1),
5287 (3, 3),
5288 (3, 3),
5289 (2, 4usize),
5290 ];
5291 let err = SymbolicSparseColMat::try_new_from_indices(nrows, ncols, indices);
5292 assert!(err.is_err());
5293 let err = err.unwrap_err();
5294 assert!(err == CreationError::OutOfBounds { row: 2, col: 4 });
5295 }
5296
5297 #[test]
5298 fn test_add_intersecting() {
5299 let lhs = SparseColMat::<usize, f64>::try_new_from_triplets(
5300 5,
5301 4,
5302 &[
5303 (1, 0, 1.0),
5304 (2, 1, 2.0),
5305 (3, 2, 3.0),
5306 (0, 0, 4.0),
5307 (1, 1, 5.0),
5308 (2, 2, 6.0),
5309 (3, 3, 7.0),
5310 (2, 0, 8.0),
5311 (3, 1, 9.0),
5312 (4, 2, 10.0),
5313 (0, 2, 11.0),
5314 (1, 3, 12.0),
5315 (4, 0, 13.0),
5316 ],
5317 )
5318 .unwrap();
5319
5320 let rhs = SparseColMat::<usize, f64>::try_new_from_triplets(
5321 5,
5322 4,
5323 &[
5324 (1, 0, 10.0),
5325 (2, 1, 14.0),
5326 (3, 2, 15.0),
5327 (4, 3, 16.0),
5328 (0, 1, 17.0),
5329 (1, 2, 18.0),
5330 (2, 3, 19.0),
5331 (3, 0, 20.0),
5332 (4, 1, 21.0),
5333 (0, 3, 22.0),
5334 ],
5335 )
5336 .unwrap();
5337
5338 let sum = ops::add(lhs.as_ref(), rhs.as_ref()).unwrap();
5339 assert!(sum.compute_nnz() == lhs.compute_nnz() + rhs.compute_nnz() - 3);
5340
5341 for j in 0..4 {
5342 for i in 0..5 {
5343 assert!(sum.row_indices_of_col_raw(j)[i] == i);
5344 }
5345 }
5346
5347 for j in 0..4 {
5348 for i in 0..5 {
5349 assert!(
5350 sum[(i, j)] == lhs.get(i, j).unwrap_or(&0.0) + rhs.get(i, j).unwrap_or(&0.0)
5351 );
5352 }
5353 }
5354 }
5355
5356 #[test]
5357 fn test_add_disjoint() {
5358 let lhs = SparseColMat::<usize, f64>::try_new_from_triplets(
5359 5,
5360 4,
5361 &[
5362 (0, 0, 1.0),
5363 (1, 1, 2.0),
5364 (2, 2, 3.0),
5365 (3, 3, 4.0),
5366 (2, 0, 5.0),
5367 (3, 1, 6.0),
5368 (4, 2, 7.0),
5369 (0, 2, 8.0),
5370 (1, 3, 9.0),
5371 (4, 0, 10.0),
5372 ],
5373 )
5374 .unwrap();
5375
5376 let rhs = SparseColMat::<usize, f64>::try_new_from_triplets(
5377 5,
5378 4,
5379 &[
5380 (1, 0, 11.0),
5381 (2, 1, 12.0),
5382 (3, 2, 13.0),
5383 (4, 3, 14.0),
5384 (0, 1, 15.0),
5385 (1, 2, 16.0),
5386 (2, 3, 17.0),
5387 (3, 0, 18.0),
5388 (4, 1, 19.0),
5389 (0, 3, 20.0),
5390 ],
5391 )
5392 .unwrap();
5393
5394 let sum = ops::add(lhs.as_ref(), rhs.as_ref()).unwrap();
5395 assert!(sum.compute_nnz() == lhs.compute_nnz() + rhs.compute_nnz());
5396
5397 for j in 0..4 {
5398 for i in 0..5 {
5399 assert!(sum.row_indices_of_col_raw(j)[i] == i);
5400 }
5401 }
5402
5403 for j in 0..4 {
5404 for i in 0..5 {
5405 assert!(
5406 sum[(i, j)] == lhs.get(i, j).unwrap_or(&0.0) + rhs.get(i, j).unwrap_or(&0.0)
5407 );
5408 }
5409 }
5410 }
5411}