COMBINATORIAL_BLAS 1.6
 
Loading...
Searching...
No Matches
SpDCCols.cpp
Go to the documentation of this file.
1/****************************************************************/
2/* Parallel Combinatorial BLAS Library (for Graph Computations) */
3/* version 1.6 -------------------------------------------------*/
4/* date: 6/15/2017 ---------------------------------------------*/
5/* authors: Ariful Azad, Aydin Buluc --------------------------*/
6/****************************************************************/
7/*
8 Copyright (c) 2010-2017, The Regents of the University of California
9
10 Permission is hereby granted, free of charge, to any person obtaining a copy
11 of this software and associated documentation files (the "Software"), to deal
12 in the Software without restriction, including without limitation the rights
13 to use, copy, modify, merge, publish, distribute, sublicense, and/or sell
14 copies of the Software, and to permit persons to whom the Software is
15 furnished to do so, subject to the following conditions:
16
17 The above copyright notice and this permission notice shall be included in
18 all copies or substantial portions of the Software.
19
20 THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS OR
21 IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
22 FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL THE
23 AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER
24 LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING FROM,
25 OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER DEALINGS IN
26 THE SOFTWARE.
27 */
28
29
30#include "SpDCCols.h"
31#include "Deleter.h"
32#include <algorithm>
33#include <functional>
34#include <vector>
35#include <climits>
36#include <iomanip>
37#include <cassert>
38
39namespace combblas {
40
41/****************************************************************************/
42/********************* PUBLIC CONSTRUCTORS/DESTRUCTORS **********************/
43/****************************************************************************/
44
45template <class IT, class NT>
46const IT SpDCCols<IT,NT>::esscount = static_cast<IT>(4);
47
48
49template <class IT, class NT>
50SpDCCols<IT,NT>::SpDCCols():dcsc(NULL), m(0), n(0), nnz(0), splits(0){
51}
52
53// Allocate all the space necessary
54template <class IT, class NT>
56:m(nRow), n(nCol), nnz(size), splits(0)
57{
58 if(nnz > 0)
59 dcsc = new Dcsc<IT,NT>(nnz, nzc);
60 else
61 dcsc = NULL;
62}
63
64template <class IT, class NT>
66{
67 if(nnz > 0)
68 {
69 if(dcsc != NULL)
70 {
71 if(splits > 0)
72 {
73 for(int i=0; i<splits; ++i)
74 delete dcscarr[i];
75 delete [] dcscarr;
76 }
77 else
78 {
79 delete dcsc;
80 }
81 }
82 }
83}
84
85// Copy constructor (constructs a new object. i.e. this is NEVER called on an existing object)
86// Derived's copy constructor can safely call Base's default constructor as base has no data members
87template <class IT, class NT>
89: m(rhs.m), n(rhs.n), nnz(rhs.nnz), splits(rhs.splits)
90{
91 if(splits > 0)
92 {
93 for(int i=0; i<splits; ++i)
94 CopyDcsc(rhs.dcscarr[i]);
95 }
96 else
97 {
98 CopyDcsc(rhs.dcsc);
99 }
100}
101
108template <class IT, class NT>
110: m(rhs.m), n(rhs.n), nnz(rhs.nnz), splits(0)
111{
112
113 if(nnz == 0) // m by n matrix of complete zeros
114 {
115 if(transpose) std::swap(m,n);
116 dcsc = NULL;
117 }
118 else
119 {
120 if(transpose)
121 {
122 std::swap(m,n);
123 IT localnzc = 1;
124 for(IT i=1; i< rhs.nnz; ++i)
125 {
126 if(rhs.rowindex(i) != rhs.rowindex(i-1))
127 {
128 ++localnzc;
129 }
130 }
131 dcsc = new Dcsc<IT,NT>(rhs.nnz,localnzc);
132 dcsc->jc[0] = rhs.rowindex(0);
133 dcsc->cp[0] = 0;
134
135 for(IT i=0; i<rhs.nnz; ++i)
136 {
137 dcsc->ir[i] = rhs.colindex(i); // copy rhs.jc to ir since this transpose=true
138 dcsc->numx[i] = rhs.numvalue(i);
139 }
140
141 IT jspos = 1;
142 for(IT i=1; i<rhs.nnz; ++i)
143 {
144 if(rhs.rowindex(i) != dcsc->jc[jspos-1])
145 {
146 dcsc->jc[jspos] = rhs.rowindex(i); // copy rhs.ir to jc since this transpose=true
147 dcsc->cp[jspos++] = i;
148 }
149 }
150 dcsc->cp[jspos] = rhs.nnz;
151 }
152 else
153 {
154 IT localnzc = 1;
155 for(IT i=1; i<rhs.nnz; ++i)
156 {
157 if(rhs.colindex(i) != rhs.colindex(i-1))
158 {
159 ++localnzc;
160 }
161 }
162 dcsc = new Dcsc<IT,NT>(rhs.nnz,localnzc);
163 dcsc->jc[0] = rhs.colindex(0);
164 dcsc->cp[0] = 0;
165
166 for(IT i=0; i<rhs.nnz; ++i)
167 {
168 dcsc->ir[i] = rhs.rowindex(i); // copy rhs.ir to ir since this transpose=false
169 dcsc->numx[i] = rhs.numvalue(i);
170 }
171
172 IT jspos = 1;
173 for(IT i=1; i<rhs.nnz; ++i)
174 {
175 if(rhs.colindex(i) != dcsc->jc[jspos-1])
176 {
177 dcsc->jc[jspos] = rhs.colindex(i); // copy rhs.jc to jc since this transpose=true
178 dcsc->cp[jspos++] = i;
179 }
180 }
181 dcsc->cp[jspos] = rhs.nnz;
182 }
183 }
184}
185
186
187
188
197template <class IT, class NT>
198SpDCCols<IT,NT>::SpDCCols(IT nRow, IT nCol, IT nTuples, const std::tuple<IT, IT, NT>* tuples, bool transpose)
199: m(nRow), n(nCol), nnz(nTuples), splits(0)
200{
201
202 if(nnz == 0) // m by n matrix of complete zeros
203 {
204 dcsc = NULL;
205 }
206 else
207 {
208 int totThreads=1;
209#ifdef _OPENMP
210#pragma omp parallel
211 {
212 totThreads = omp_get_num_threads();
213 }
214#endif
215
216 std::vector <IT> tstart(totThreads);
217 std::vector <IT> tend(totThreads);
218 std::vector <IT> tdisp(totThreads+1);
219
220 // extra memory, but replaces an O(nnz) loop by an O(nzc) loop
221 IT* temp_jc = new IT[nTuples];
222 IT* temp_cp = new IT[nTuples];
223
224#ifdef _OPENMP
225#pragma omp parallel
226#endif
227 {
228 int threadID = 0;
229#ifdef _OPENMP
230 threadID = omp_get_thread_num();
231#endif
232 IT start = threadID * (nTuples / totThreads);
233 IT end = (threadID + 1) * (nTuples / totThreads);
234 if(threadID == (totThreads-1)) end = nTuples;
235 IT curpos=start;
236 if(end>start) // no work for the current thread
237 {
238 temp_jc[start] = std::get<1>(tuples[start]);
239 temp_cp[start] = start;
240 for (IT i = start+1; i < end; ++i)
241 {
242 if(std::get<1>(tuples[i]) != temp_jc[curpos] )
243 {
244 temp_jc[++curpos] = std::get<1>(tuples[i]);
245 temp_cp[curpos] = i;
246 }
247 }
248 }
249
250 tstart[threadID] = start;
251 if(end>start) tend[threadID] = curpos+1;
252 else tend[threadID] = end; // start=end
253 }
254
255
256 // serial part
257 for(int t=totThreads-1; t>0; --t)
258 {
259 if(tend[t] > tstart[t] && tend[t-1] > tstart[t-1])
260 {
261 if(temp_jc[tstart[t]] == temp_jc[tend[t-1]-1])
262 {
263 tstart[t] ++;
264 }
265 }
266 }
267
268 tdisp[0] = 0;
269 for(int t=0; t<totThreads; ++t)
270 {
271 tdisp[t+1] = tdisp[t] + tend[t] - tstart[t];
272 }
273
276
277#ifdef _OPENMP
278#pragma omp parallel
279#endif
280 {
281 int threadID = 0;
282#ifdef _OPENMP
283 threadID = omp_get_thread_num();
284#endif
285 std::copy(temp_jc + tstart[threadID], temp_jc + tend[threadID], dcsc->jc + tdisp[threadID]);
286 std::copy(temp_cp + tstart[threadID], temp_cp + tend[threadID], dcsc->cp + tdisp[threadID]);
287 }
288 dcsc->cp[localnzc] = nTuples;
289
290 delete [] temp_jc;
291 delete [] temp_cp;
292
293#ifdef _OPENMP
294#pragma omp parallel for schedule (static)
295#endif
296 for(IT i=0; i<nTuples; ++i)
297 {
298 dcsc->ir[i] = std::get<0>(tuples[i]);
299 dcsc->numx[i] = std::get<2>(tuples[i]);
300 }
301 }
302
303 if(transpose) Transpose(); // this is not efficient, think to improve later. We included this parameter anyway to make this constructor different from another constracttor when the fourth argument is passed as 0.
304}
305
306
307
308/*
309template <class IT, class NT>
310SpDCCols<IT,NT>::SpDCCols(IT nRow, IT nCol, IT nTuples, const tuple<IT, IT, NT>* tuples)
311: m(nRow), n(nCol), nnz(nTuples), splits(0)
312{
313
314 if(nnz == 0) // m by n matrix of complete zeros
315 {
316 dcsc = NULL;
317 }
318 else
319 {
320 IT localnzc = 1;
321#pragma omp parallel for schedule (static) default(shared) reduction(+:localnzc)
322 for(IT i=1; i<nTuples; ++i) // not scaling well, try my own version
323 {
324 if(std::get<1>(tuples[i]) != std::get<1>(tuples[i-1]))
325 {
326 ++localnzc;
327 }
328 }
329
330 dcsc = new Dcsc<IT,NT>(nTuples,localnzc);
331 dcsc->jc[0] = std::get<1>(tuples[0]);
332 dcsc->cp[0] = 0;
333
334#pragma omp parallel for schedule (static)
335 for(IT i=0; i<nTuples; ++i)
336 {
337 dcsc->ir[i] = std::get<0>(tuples[i]);
338 dcsc->numx[i] = std::get<2>(tuples[i]);
339 }
340
341 IT jspos = 1;
342 for(IT i=1; i<nTuples; ++i) // now this loop
343 {
344 if(std::get<1>(tuples[i]) != dcsc->jc[jspos-1])
345 {
346 dcsc->jc[jspos] = std::get<1>(tuples[i]);
347 dcsc->cp[jspos++] = i;
348 }
349 }
350 dcsc->cp[jspos] = nTuples;
351 }
352}
353*/
354
355/****************************************************************************/
356/************************** PUBLIC OPERATORS ********************************/
357/****************************************************************************/
358
364template <class IT, class NT>
366{
367 // this pointer stores the address of the class instance
368 // check for self assignment using address comparison
369 if(this != &rhs)
370 {
371 if(dcsc != NULL && nnz > 0)
372 {
373 delete dcsc;
374 }
375 if(rhs.dcsc != NULL)
376 {
377 dcsc = new Dcsc<IT,NT>(*(rhs.dcsc));
378 nnz = rhs.nnz;
379 }
380 else
381 {
382 dcsc = NULL;
383 nnz = 0;
384 }
385
386 m = rhs.m;
387 n = rhs.n;
388 splits = rhs.splits;
389 }
390 return *this;
391}
392
393template <class IT, class NT>
395{
396 // this pointer stores the address of the class instance
397 // check for self assignment using address comparison
398 if(this != &rhs)
399 {
400 if(m == rhs.m && n == rhs.n)
401 {
402 if(rhs.nnz == 0)
403 {
404 return *this;
405 }
406 else if(nnz == 0)
407 {
408 dcsc = new Dcsc<IT,NT>(*(rhs.dcsc));
409 nnz = dcsc->nz;
410 }
411 else
412 {
413 (*dcsc) += (*(rhs.dcsc));
414 nnz = dcsc->nz;
415 }
416 }
417 else
418 {
419 std::cout<< "Not addable: " << m << "!=" << rhs.m << " or " << n << "!=" << rhs.n <<std::endl;
420 }
421 }
422 else
423 {
424 std::cout<< "Missing feature (A+A): Use multiply with 2 instead !"<<std::endl;
425 }
426 return *this;
427}
428
429template <class IT, class NT>
430template <typename _UnaryOperation, typename GlobalIT>
432{
433 if(nnz > 0)
434 {
436 if (inPlace)
437 {
438 nnz = dcsc->nz;
439
440 if(nnz == 0)
441 {
442 delete dcsc;
443 dcsc = NULL;
444 }
445 return NULL;
446 }
447 else
448 {
449 // wrap the new pruned Dcsc into a new SpDCCols
451 retcols->dcsc = ret;
452 retcols->nnz = retcols->dcsc->nz;
453 retcols->n = n;
454 retcols->m = m;
455 return retcols;
456 }
457 }
458 else
459 {
460 if (inPlace)
461 {
462 return NULL;
463 }
464 else
465 {
467 retcols->dcsc = NULL;
468 retcols->nnz = 0;
469 retcols->n = n;
470 retcols->m = m;
471 return retcols;
472 }
473 }
474}
475
476template <class IT, class NT>
477template <typename _UnaryOperation>
479{
480 if(nnz > 0)
481 {
482 Dcsc<IT,NT>* ret = dcsc->Prune (__unary_op, inPlace);
483 if (inPlace)
484 {
485 nnz = dcsc->nz;
486
487 if(nnz == 0)
488 {
489 delete dcsc;
490 dcsc = NULL;
491 }
492 return NULL;
493 }
494 else
495 {
496 // wrap the new pruned Dcsc into a new SpDCCols
498 retcols->dcsc = ret;
499 retcols->nnz = retcols->dcsc->nz;
500 retcols->n = n;
501 retcols->m = m;
502 return retcols;
503 }
504 }
505 else
506 {
507 if (inPlace)
508 {
509 return NULL;
510 }
511 else
512 {
514 retcols->dcsc = NULL;
515 retcols->nnz = 0;
516 retcols->n = n;
517 retcols->m = m;
518 return retcols;
519 }
520 }
521}
522
523
524
525template <class IT, class NT>
526template <typename _BinaryOperation>
528{
529 if(nnz > 0)
530 {
531 Dcsc<IT,NT>* ret = dcsc->PruneColumn (pvals, __binary_op, inPlace);
532 if (inPlace)
533 {
534 nnz = dcsc->nz;
535
536 if(nnz == 0)
537 {
538 delete dcsc;
539 dcsc = NULL;
540 }
541 return NULL;
542 }
543 else
544 {
545 // wrap the new pruned Dcsc into a new SpDCCols
547 retcols->dcsc = ret;
548 retcols->nnz = retcols->dcsc->nz;
549 retcols->n = n;
550 retcols->m = m;
551 return retcols;
552 }
553 }
554 else
555 {
556 if (inPlace)
557 {
558 return NULL;
559 }
560 else
561 {
563 retcols->dcsc = NULL;
564 retcols->nnz = 0;
565 retcols->n = n;
566 retcols->m = m;
567 return retcols;
568 }
569 }
570}
571
572template <class IT, class NT>
573void SpDCCols<IT,NT>::PruneColumnByIndex(const std::vector<IT>& ci)
574{
575 if (nnz > 0)
576 {
577 dcsc->PruneColumnByIndex(ci);
578 nnz = dcsc->nz;
579 }
580}
581
582template <class IT, class NT>
583template <typename _BinaryOperation>
585{
586 if(nnz > 0)
587 {
588 Dcsc<IT,NT>* ret = dcsc->PruneColumn (pinds, pvals, __binary_op, inPlace);
589 if (inPlace)
590 {
591 nnz = dcsc->nz;
592
593 if(nnz == 0)
594 {
595 delete dcsc;
596 dcsc = NULL;
597 }
598 return NULL;
599 }
600 else
601 {
602 // wrap the new pruned Dcsc into a new SpDCCols
604 retcols->dcsc = ret;
605 retcols->nnz = retcols->dcsc->nz;
606 retcols->n = n;
607 retcols->m = m;
608 return retcols;
609 }
610 }
611 else
612 {
613 if (inPlace)
614 {
615 return NULL;
616 }
617 else
618 {
620 retcols->dcsc = NULL;
621 retcols->nnz = 0;
622 retcols->n = n;
623 retcols->m = m;
624 return retcols;
625 }
626 }
627}
628
629
630template <class IT, class NT>
632{
633 if(this != &rhs)
634 {
635 if(m == rhs.m && n == rhs.n)
636 {
637 if(rhs.nnz == 0)
638 {
639 return;
640 }
641 else if (rhs.nnz != 0 && nnz != 0)
642 {
643 dcsc->SetDifference (*(rhs.dcsc));
644 nnz = dcsc->nz;
645 if(nnz == 0 )
646 dcsc = NULL;
647 }
648 }
649 else
650 {
651 std::cout<< "Matrices do not conform for A - B !"<<std::endl;
652 }
653 }
654 else
655 {
656 std::cout<< "Missing feature (A .* A): Use Square_EWise() instead !"<<std::endl;
657 }
658}
659
660// Aydin (June 2021): Make the exclude case of this call SetDifference above instead
661template <class IT, class NT>
663{
664 if(this != &rhs)
665 {
666 if(m == rhs.m && n == rhs.n)
667 {
668 if(rhs.nnz == 0)
669 {
670 if(exclude) // then we don't exclude anything
671 {
672 return;
673 }
674 else // A .* zeros() is zeros()
675 {
676 *this = SpDCCols<IT,NT>(0,m,n,0); // completely reset the matrix
677 }
678 }
679 else if (rhs.nnz != 0 && nnz != 0)
680 {
681 dcsc->EWiseMult (*(rhs.dcsc), exclude);
682 nnz = dcsc->nz;
683 if(nnz == 0 )
684 dcsc = NULL;
685 }
686 }
687 else
688 {
689 std::cout<< "Matrices do not conform for A .* op(B) !"<<std::endl;
690 }
691 }
692 else
693 {
694 std::cout<< "Missing feature (A .* A): Use Square_EWise() instead !"<<std::endl;
695 }
696}
697
701template <class IT, class NT>
703{
704 if(m == m_scaler && n == n_scaler)
705 {
706 if(nnz > 0)
707 dcsc->EWiseScale (scaler);
708 }
709 else
710 {
711 std::cout<< "Matrices do not conform for EWiseScale !"<<std::endl;
712 }
713}
714
715
716/****************************************************************************/
717/********************* PUBLIC MEMBER FUNCTIONS ******************************/
718/****************************************************************************/
719
720template <class IT, class NT>
722{
723 m = _m;
724 n = _n;
725 nnz = _nz;
726
727 if(nnz > 0)
728 dcsc = new Dcsc<IT,NT>(_cp, _jc, _ir, _numx, _nz, _nzc, false); // memory not owned by DCSC
729 else
730 dcsc = NULL;
731}
732
733template <class IT, class NT>
734void SpDCCols<IT,NT>::CreateImpl(const std::vector<IT> & essentials)
735{
736 assert(essentials.size() == esscount);
737 nnz = essentials[0];
738 m = essentials[1];
739 n = essentials[2];
740
741 if(nnz > 0)
742 dcsc = new Dcsc<IT,NT>(nnz,essentials[3]);
743 else
744 dcsc = NULL;
745}
746
747template <class IT, class NT>
748void SpDCCols<IT,NT>::CreateImpl(IT size, IT nRow, IT nCol, std::tuple<IT, IT, NT> * mytuples)
749{
751 tuples.SortColBased();
752
753#ifdef DEBUG
754 std::pair<IT,IT> rlim = tuples.RowLimits();
755 std::pair<IT,IT> clim = tuples.ColLimits();
756
757 std::ofstream oput;
758 std::stringstream ss;
759 std::string rank;
760 int myrank;
762 ss << myrank;
763 ss >> rank;
764 std::string ofilename = "Read";
765 ofilename += rank;
766 oput.open(ofilename.c_str(), std::ios_base::app );
767 oput << "Creating of dimensions " << nRow << "-by-" << nCol << " of size: " << size <<
768 " with row range (" << rlim.first << "," << rlim.second << ") and column range (" << clim.first << "," << clim.second << ")" << std::endl;
769 if(tuples.getnnz() > 0)
770 {
771 IT minfr = joker::get<0>(tuples.front());
772 IT minto = joker::get<1>(tuples.front());
773 IT maxfr = joker::get<0>(tuples.back());
774 IT maxto = joker::get<1>(tuples.back());
775
776 oput << "Min: " << minfr << ", " << minto << "; Max: " << maxfr << ", " << maxto << std::endl;
777 }
778 oput.close();
779#endif
780
781 SpDCCols<IT,NT> object(tuples, false);
782 *this = object;
783}
784
785
786template <class IT, class NT>
787std::vector<IT> SpDCCols<IT,NT>::GetEssentials() const
788{
789 std::vector<IT> essentials(esscount);
790 essentials[0] = nnz;
791 essentials[1] = m;
792 essentials[2] = n;
793 essentials[3] = (nnz > 0) ? dcsc->nzc : 0;
794 return essentials;
795}
796
797template <class IT, class NT>
798template <typename NNT>
800{
802 if(nnz > 0)
803 convert = new Dcsc<IT,NNT>(*dcsc);
804 else
805 convert = NULL;
806
807 return SpDCCols<IT,NNT>(m, n, convert);
808}
809
810
811template <class IT, class NT>
812template <typename NIT, typename NNT>
814{
816 if(nnz > 0)
817 convert = new Dcsc<NIT,NNT>(*dcsc);
818 else
819 convert = NULL;
820
821 return SpDCCols<NIT,NNT>(m, n, convert);
822}
823
824
825template <class IT, class NT>
827{
828 Arr<IT,NT> arr(3,1);
829
830 if(nnz > 0)
831 {
832 arr.indarrs[0] = LocArr<IT,IT>(dcsc->cp, dcsc->nzc+1);
833 arr.indarrs[1] = LocArr<IT,IT>(dcsc->jc, dcsc->nzc);
834 arr.indarrs[2] = LocArr<IT,IT>(dcsc->ir, dcsc->nz);
835 arr.numarrs[0] = LocArr<NT,IT>(dcsc->numx, dcsc->nz);
836 }
837 else
838 {
839 arr.indarrs[0] = LocArr<IT,IT>(NULL, 0);
840 arr.indarrs[1] = LocArr<IT,IT>(NULL, 0);
841 arr.indarrs[2] = LocArr<IT,IT>(NULL, 0);
842 arr.numarrs[0] = LocArr<NT,IT>(NULL, 0);
843
844 }
845 return arr;
846}
847
853template <class IT, class NT>
855{
856 if(nnz > 0)
857 {
859 Atuples.SortRowBased();
860
861 // destruction of (*this) is handled by the assignment operator
862 *this = SpDCCols<IT,NT>(Atuples,true);
863 }
864 else
865 {
866 *this = SpDCCols<IT,NT>(0, n, m, 0);
867 }
868}
869
870
876template <class IT, class NT>
878{
880 Atuples.SortRowBased();
881
882 return SpDCCols<IT,NT>(Atuples,true);
883}
884
890template <class IT, class NT>
892{
894 Atuples.SortRowBased();
895
896 return new SpDCCols<IT,NT>(Atuples,true);
897}
898
905template <class IT, class NT>
907{
908 IT cut = n/2;
909 if(cut == 0)
910 {
911 std::cout<< "Matrix is too small to be splitted" << std::endl;
912 return;
913 }
914
917
918 if(nnz != 0)
919 {
920 dcsc->Split(Adcsc, Bdcsc, cut);
921 }
922
924 partB = SpDCCols<IT,NT> (m, n-cut, Bdcsc);
925
926 // handle destruction through assignment operator
927 *this = SpDCCols<IT, NT>();
928}
929
935template <class IT, class NT>
937{
938 if(parts < 2)
939 {
940 matrices.emplace_back(*this);
941 }
942 else
943 {
944 std::vector<IT> cuts(parts-1);
945 for(int i=0; i< (parts-1); ++i)
946 {
947 cuts[i] = (i+1) * (n/parts);
948 }
949 if(n < parts)
950 {
951 std::cout<< "Matrix is too small to be splitted" << std::endl;
952 return;
953 }
954 std::vector< Dcsc<IT,NT> * > dcscs(parts, NULL);
955
956 if(nnz != 0)
957 {
958 dcsc->ColSplit(dcscs, cuts);
959 }
960
961 for(int i=0; i< (parts-1); ++i)
962 {
964 matrices.emplace_back(matrix);
965 }
967 matrices.emplace_back(matrix);
968 }
969 *this = SpDCCols<IT, NT>(); // handle destruction through assignment operator
970}
971
977template <class IT, class NT>
979{
980 if(parts < 2)
981 {
982 matrices.emplace_back(new SpDCCols<IT,NT>(*this));
983 }
984 else
985 {
986 std::vector<IT> cuts(parts-1);
987 for(int i=0; i< (parts-1); ++i)
988 {
989 cuts[i] = (i+1) * (n/parts);
990 }
991 if(n < parts)
992 {
993 std::cout<< "Matrix is too small to be splitted" << std::endl;
994 return;
995 }
996 std::vector< Dcsc<IT,NT> * > dcscs(parts, NULL);
997
998 if(nnz != 0)
999 {
1000 dcsc->ColSplit(dcscs, cuts);
1001 }
1002
1003 for(int i=0; i< (parts-1); ++i)
1004 {
1006 matrices.emplace_back(matrix);
1007 }
1009 matrices.emplace_back(matrix);
1010 }
1011 *this = SpDCCols<IT, NT>(); // handle destruction through assignment operator
1012}
1013
1019template <class IT, class NT>
1020void SpDCCols<IT,NT>::ColSplit(std::vector<IT> & cutSizes, std::vector< SpDCCols<IT,NT> > & matrices)
1021{
1022 IT totn = 0;
1023 int parts = cutSizes.size();
1024 for(int i = 0; i < parts; i++) totn += cutSizes[i];
1025 if(parts < 2){
1026 matrices.emplace_back(*this);
1027 }
1028 else if(totn != n){
1029 std::cout << "Cut sizes are not appropriate" << std::endl;
1030 return;
1031 }
1032 else{
1033 std::vector<IT> cuts(parts-1);
1034 cuts[0] = cutSizes[0];
1035 for(int i = 1; i < parts-1; i++){
1036 cuts[i] = cuts[i-1] + cutSizes[i];
1037 }
1038 std::vector< Dcsc<IT,NT> * > dcscs(parts, NULL);
1039
1040 if(nnz != 0){
1041 dcsc->ColSplit(dcscs, cuts);
1042 }
1043
1044 for(int i=0; i< parts; ++i){
1046 matrices.emplace_back(matrix);
1047 }
1048 }
1049 *this = SpDCCols<IT, NT>();
1050}
1051
1058template <class IT, class NT>
1059void SpDCCols<IT,NT>::ColSplit(std::vector<IT> & cutSizes, std::vector< SpDCCols<IT,NT>* > & matrices)
1060{
1061 IT totn = 0;
1062 int parts = cutSizes.size();
1063 for(int i = 0; i < parts; i++) totn += cutSizes[i];
1064 if(parts < 2){
1065 matrices.emplace_back(new SpDCCols<IT,NT>(*this));
1066 }
1067 else if(totn != n){
1068 std::cout << "Cut sizes are not appropriate" << std::endl;
1069 return;
1070 }
1071 else{
1072 std::vector<IT> cuts(parts-1);
1073 cuts[0] = cutSizes[0];
1074 for(int i = 1; i < parts-1; i++){
1075 cuts[i] = cuts[i-1] + cutSizes[i];
1076 }
1077 std::vector< Dcsc<IT,NT> * > dcscs(parts, NULL);
1078
1079 if(nnz != 0){
1080 dcsc->ColSplit(dcscs, cuts);
1081 }
1082
1083 for(int i=0; i< parts; ++i){
1085 matrices.emplace_back(matrix);
1086 }
1087 }
1088 *this = SpDCCols<IT, NT>();
1089}
1090
1095template <class IT, class NT>
1097{
1098 std::vector< SpDCCols<IT,NT> * > nonempties;
1099 std::vector< Dcsc<IT,NT> * > dcscs;
1100 std::vector< IT > offsets;
1101 IT runningoffset = 0;
1102
1103 for(size_t i=0; i< matrices.size(); ++i)
1104 {
1105 if(matrices[i].nnz != 0)
1106 {
1107 nonempties.push_back(&(matrices[i]));
1108 dcscs.push_back(matrices[i].dcsc);
1109 offsets.push_back(runningoffset);
1110 }
1111 runningoffset += matrices[i].n;
1112 }
1113
1114 if(nonempties.size() < 1)
1115 {
1116#ifdef DEBUG
1117 std::cout << "Nothing to ColConcatenate" << std::endl;
1118#endif
1119 n = runningoffset;
1120 }/*
1121 else if(nonempties.size() < 2)
1122 {
1123 *this = *(nonempties[0]);
1124 n = runningoffset;
1125 }*/
1126 else // nonempties.size() > 1
1127 {
1128 Dcsc<IT,NT> * Cdcsc = new Dcsc<IT,NT>();
1129 Cdcsc->ColConcatenate(dcscs, offsets);
1131 }
1132
1133 // destruct parameters
1134 for(size_t i=0; i< matrices.size(); ++i)
1135 {
1137 }
1138}
1139
1144template <class IT, class NT>
1146{
1147 std::vector< SpDCCols<IT,NT> * > nonempties;
1148 std::vector< Dcsc<IT,NT> * > dcscs;
1149 std::vector< IT > offsets;
1150 IT runningoffset = 0;
1151
1152 for(size_t i=0; i< matrices.size(); ++i)
1153 {
1154 if(matrices[i]->nnz != 0)
1155 {
1156 nonempties.push_back(matrices[i]);
1157 dcscs.push_back(matrices[i]->dcsc);
1158 offsets.push_back(runningoffset);
1159 }
1160 runningoffset += matrices[i]->n;
1161 }
1162
1163 if(nonempties.size() < 1)
1164 {
1165#ifdef DEBUG
1166 std::cout << "Nothing to ColConcatenate" << std::endl;
1167#endif
1168 n = runningoffset;
1169 }/*
1170 else if(nonempties.size() < 2)
1171 {
1172 *this = *(nonempties[0]);
1173 n = runningoffset;
1174 }*/
1175 else // nonempties.size() > 1
1176 {
1177 Dcsc<IT,NT> * Cdcsc = new Dcsc<IT,NT>();
1178 Cdcsc->ColConcatenate(dcscs, offsets);
1180 }
1181
1182 // destruct parameters
1183 for(size_t i=0; i< matrices.size(); ++i)
1184 {
1185 delete matrices[i];
1186 }
1187}
1188
1189
1194template <class IT, class NT>
1196{
1197 assert( partA.m == partB.m );
1198
1199 Dcsc<IT,NT> * Cdcsc = new Dcsc<IT,NT>();
1200
1201 if(partA.nnz == 0 && partB.nnz == 0)
1202 {
1203 Cdcsc = NULL;
1204 }
1205 else if(partA.nnz == 0)
1206 {
1207 Cdcsc = new Dcsc<IT,NT>(*(partB.dcsc));
1208 std::transform(Cdcsc->jc, Cdcsc->jc + Cdcsc->nzc, Cdcsc->jc, std::bind2nd(std::plus<IT>(), partA.n));
1209 }
1210 else if(partB.nnz == 0)
1211 {
1212 Cdcsc = new Dcsc<IT,NT>(*(partA.dcsc));
1213 }
1214 else
1215 {
1216 Cdcsc->Merge(partA.dcsc, partB.dcsc, partA.n);
1217 }
1218 *this = SpDCCols<IT,NT> (partA.m, partA.n + partB.n, Cdcsc);
1219
1222}
1223
1230template <class IT, class NT>
1231template <class SR>
1233{
1234 if(A.isZero() || B.isZero())
1235 {
1236 return -1; // no need to do anything
1237 }
1239 SpHelper::SpIntersect(*(A.dcsc), *(B.dcsc), cols, rows, isect1, isect2, itr1, itr2);
1240
1241 IT kisect = static_cast<IT>(itr1-isect1); // size of the intersection ((itr1-isect1) == (itr2-isect2))
1242 if(kisect == 0)
1243 {
1245 return -1;
1246 }
1247
1249 IT cnz = SpHelper::SpCartesian< SR > (*(A.dcsc), *(B.dcsc), kisect, isect1, isect2, multstack);
1251
1252 IT mdim = A.m;
1253 IT ndim = B.m; // since B has already been transposed
1254 if(isZero())
1255 {
1256 dcsc = new Dcsc<IT,NT>(multstack, mdim, ndim, cnz);
1257 }
1258 else
1259 {
1260 dcsc->AddAndAssign(multstack, mdim, ndim, cnz);
1261 }
1262 nnz = dcsc->nz;
1263
1264 delete [] multstack;
1265 return 1;
1266}
1267
1274template <class IT, class NT>
1275template <typename SR>
1277{
1278 if(A.isZero() || B.isZero())
1279 {
1280 return -1; // no need to do anything
1281 }
1283 int cnz = SpHelper::SpColByCol< SR > (*(A.dcsc), *(B.dcsc), A.n, multstack);
1284
1285 IT mdim = A.m;
1286 IT ndim = B.n;
1287 if(isZero())
1288 {
1289 dcsc = new Dcsc<IT,NT>(multstack, mdim, ndim, cnz);
1290 }
1291 else
1292 {
1293 dcsc->AddAndAssign(multstack, mdim, ndim, cnz);
1294 }
1295 nnz = dcsc->nz;
1296
1297 delete [] multstack;
1298 return 1;
1299}
1300
1301
1302template <class IT, class NT>
1303template <typename SR>
1305{
1306 std::cout << "PlusEq_AtXBn function has not been implemented yet !" << std::endl;
1307 return 0;
1308}
1309
1310template <class IT, class NT>
1311template <typename SR>
1313{
1314 std::cout << "PlusEq_AtXBt function has not been implemented yet !" << std::endl;
1315 return 0;
1316}
1317
1318
1319template <class IT, class NT>
1321{
1322 IT * itr = std::find(dcsc->jc, dcsc->jc + dcsc->nzc, ci);
1323 if(itr != dcsc->jc + dcsc->nzc)
1324 {
1325 IT irbeg = dcsc->cp[itr - dcsc->jc];
1326 IT irend = dcsc->cp[itr - dcsc->jc + 1];
1327
1328 IT * ele = std::find(dcsc->ir + irbeg, dcsc->ir + irend, ri);
1329 if(ele != dcsc->ir + irend)
1330 {
1331 SpDCCols<IT,NT> SingEleMat(1, 1, 1, 1); // 1-by-1 matrix with 1 nonzero
1332 *(SingEleMat.dcsc->numx) = dcsc->numx[ele - dcsc->ir];
1333 *(SingEleMat.dcsc->ir) = *ele;
1334 *(SingEleMat.dcsc->jc) = *itr;
1335 (SingEleMat.dcsc->cp)[0] = 0;
1336 (SingEleMat.dcsc->cp)[1] = 1;
1337 return SingEleMat;
1338 }
1339 else
1340 {
1341 return SpDCCols<IT,NT>(); // 0-by-0 empty matrix
1342 }
1343 }
1344 else
1345 {
1346 return SpDCCols<IT,NT>(); // 0-by-0 empty matrix
1347 }
1348}
1349
1354template <class IT, class NT>
1355SpDCCols<IT,NT> SpDCCols<IT,NT>::operator() (const std::vector<IT> & ri, const std::vector<IT> & ci) const
1356{
1357 typedef PlusTimesSRing<NT,NT> PT;
1358
1359 IT rsize = ri.size();
1360 IT csize = ci.size();
1361
1362 if(rsize == 0 && csize == 0)
1363 {
1364 // return an m x n matrix of complete zeros
1365 // since we don't know whether columns or rows are indexed
1366 return SpDCCols<IT,NT> (0, m, n, 0);
1367 }
1368 else if(rsize == 0)
1369 {
1370 return ColIndex(ci);
1371 }
1372 else if(csize == 0)
1373 {
1374 SpDCCols<IT,NT> LeftMatrix(rsize, rsize, this->m, ri, true);
1375 return LeftMatrix.OrdColByCol< PT >(*this);
1376 }
1377 else // this handles the (rsize=1 && csize=1) case as well
1378 {
1379 SpDCCols<IT,NT> LeftMatrix(rsize, rsize, this->m, ri, true);
1380 SpDCCols<IT,NT> RightMatrix(csize, this->n, csize, ci, false);
1381 return LeftMatrix.OrdColByCol< PT >( OrdColByCol< PT >(RightMatrix) );
1382 }
1383}
1384
1385template <class IT, class NT>
1386std::ofstream & SpDCCols<IT,NT>::put(std::ofstream & outfile) const
1387{
1388 if(nnz == 0)
1389 {
1390 outfile << "Matrix doesn't have any nonzeros" <<std::endl;
1391 return outfile;
1392 }
1393 SpTuples<IT,NT> tuples(*this);
1394 outfile << tuples << std::endl;
1395 return outfile;
1396}
1397
1398
1399template <class IT, class NT>
1400std::ifstream & SpDCCols<IT,NT>::get(std::ifstream & infile)
1401{
1402 std::cout << "Getting... SpDCCols" << std::endl;
1403 IT m, n, nnz;
1404 infile >> m >> n >> nnz;
1405 SpTuples<IT,NT> tuples(nnz, m, n);
1406 infile >> tuples;
1407 tuples.SortColBased();
1408
1409 SpDCCols<IT,NT> object(tuples, false);
1410 *this = object;
1411 return infile;
1412}
1413
1414
1415template<class IT, class NT>
1416void SpDCCols<IT,NT>::PrintInfo(std::ofstream & out) const
1417{
1418 out << "m: " << m ;
1419 out << ", n: " << n ;
1420 out << ", nnz: "<< nnz ;
1421
1422 if(splits > 0)
1423 {
1424 out << ", local splits: " << splits << std::endl;
1425 }
1426 else
1427 {
1428 if(dcsc != NULL)
1429 {
1430 out << ", nzc: "<< dcsc->nzc << std::endl;
1431 }
1432 else
1433 {
1434 out <<", nzc: "<< 0 << std::endl;
1435 }
1436 }
1437}
1438
1439template<class IT, class NT>
1441{
1442 std::cout << "m: " << m ;
1443 std::cout << ", n: " << n ;
1444 std::cout << ", nnz: "<< nnz ;
1445
1446 if(splits > 0)
1447 {
1448 std::cout << ", local splits: " << splits << std::endl;
1449 }
1450 else
1451 {
1452 if(dcsc != NULL)
1453 {
1454 std::cout << ", nzc: "<< dcsc->nzc << std::endl;
1455 }
1456 else
1457 {
1458 std::cout <<", nzc: "<< 0 << std::endl;
1459 }
1460
1461 if(m < PRINT_LIMIT && n < PRINT_LIMIT) // small enough to print
1462 {
1463 std::string ** A = SpHelper::allocate2D<std::string>(m,n);
1464 for(IT i=0; i< m; ++i)
1465 for(IT j=0; j<n; ++j)
1466 A[i][j] = "-";
1467 if(dcsc != NULL)
1468 {
1469 for(IT i=0; i< dcsc->nzc; ++i)
1470 {
1471 for(IT j = dcsc->cp[i]; j<dcsc->cp[i+1]; ++j)
1472 {
1473 IT colid = dcsc->jc[i];
1474 IT rowid = dcsc->ir[j];
1475 A[rowid][colid] = std::to_string(dcsc->numx[j]);
1476 }
1477 }
1478 }
1479 for(IT i=0; i< m; ++i)
1480 {
1481 for(IT j=0; j<n; ++j)
1482 {
1483 std::cout << A[i][j];
1484 std::cout << "\t";
1485 }
1486 std::cout << std::endl;
1487 }
1489 }
1490 }
1491}
1492
1493
1494/****************************************************************************/
1495/********************* PRIVATE CONSTRUCTORS/DESTRUCTORS *********************/
1496/****************************************************************************/
1497
1499template <class IT, class NT>
1501:dcsc(mydcsc), m(nRow), n(nCol), splits(0)
1502{
1503 if (mydcsc == NULL)
1504 nnz = 0;
1505 else
1506 nnz = mydcsc->nz;
1507}
1508
1510template <class IT, class NT>
1511SpDCCols<IT,NT>::SpDCCols (IT size, IT nRow, IT nCol, const std::vector<IT> & indices, bool isRow)
1512:m(nRow), n(nCol), nnz(size), splits(0)
1513{
1514 if(size > 0)
1515 dcsc = new Dcsc<IT,NT>(size,indices,isRow);
1516 else
1517 dcsc = NULL;
1518}
1519
1520
1521/****************************************************************************/
1522/************************* PRIVATE MEMBER FUNCTIONS *************************/
1523/****************************************************************************/
1524
1525template <class IT, class NT>
1526inline void SpDCCols<IT,NT>::CopyDcsc(Dcsc<IT,NT> * source)
1527{
1528 // source dcsc will be NULL if number of nonzeros is zero
1529 if(source != NULL)
1530 dcsc = new Dcsc<IT,NT>(*source);
1531 else
1532 dcsc = NULL;
1533}
1534
1541template <class IT, class NT>
1542SpDCCols<IT,NT> SpDCCols<IT,NT>::ColIndex(const std::vector<IT> & ci) const
1543{
1544 IT csize = ci.size();
1545 if(nnz == 0) // nothing to index
1546 {
1547 return SpDCCols<IT,NT>(0, m, csize, 0);
1548 }
1549 else if(ci.empty())
1550 {
1551 return SpDCCols<IT,NT>(0, m,0, 0);
1552 }
1553
1554 // First pass for estimation
1555 IT estsize = 0;
1556 IT estnzc = 0;
1557 for(IT i=0, j=0; i< dcsc->nzc && j < csize;)
1558 {
1559 if((dcsc->jc)[i] < ci[j])
1560 {
1561 ++i;
1562 }
1563 else if ((dcsc->jc)[i] > ci[j])
1564 {
1565 ++j;
1566 }
1567 else
1568 {
1569 estsize += (dcsc->cp)[i+1] - (dcsc->cp)[i];
1570 ++estnzc;
1571 ++i;
1572 ++j;
1573 }
1574 }
1575
1577 if(estnzc == 0)
1578 {
1579 return SubA; // no need to run the second pass
1580 }
1581 SubA.dcsc->cp[0] = 0;
1582 IT cnzc = 0;
1583 IT cnz = 0;
1584 for(IT i=0, j=0; i < dcsc->nzc && j < csize;)
1585 {
1586 if((dcsc->jc)[i] < ci[j])
1587 {
1588 ++i;
1589 }
1590 else if ((dcsc->jc)[i] > ci[j]) // an empty column for the output
1591 {
1592 ++j;
1593 }
1594 else
1595 {
1596 IT columncount = (dcsc->cp)[i+1] - (dcsc->cp)[i];
1597 SubA.dcsc->jc[cnzc++] = j;
1598 SubA.dcsc->cp[cnzc] = SubA.dcsc->cp[cnzc-1] + columncount;
1599 std::copy(dcsc->ir + dcsc->cp[i], dcsc->ir + dcsc->cp[i+1], SubA.dcsc->ir + cnz);
1600 std::copy(dcsc->numx + dcsc->cp[i], dcsc->numx + dcsc->cp[i+1], SubA.dcsc->numx + cnz);
1601 cnz += columncount;
1602 ++i;
1603 ++j;
1604 }
1605 }
1606 return SubA;
1607}
1608
1609template <class IT, class NT>
1610template <typename SR, typename NTR>
1611SpDCCols< IT, typename promote_trait<NT,NTR>::T_promote > SpDCCols<IT,NT>::OrdOutProdMult(const SpDCCols<IT,NTR> & rhs) const
1612{
1613 typedef typename promote_trait<NT,NTR>::T_promote T_promote;
1614
1615 if(isZero() || rhs.isZero())
1616 {
1617 return SpDCCols< IT, T_promote > (0, m, rhs.n, 0); // return an empty matrix
1618 }
1619 SpDCCols<IT,NTR> Btrans = rhs.TransposeConst();
1620
1623
1624 IT kisect = static_cast<IT>(itr1-isect1); // size of the intersection ((itr1-isect1) == (itr2-isect2))
1625 if(kisect == 0)
1626 {
1628 return SpDCCols< IT, T_promote > (0, m, rhs.n, 0);
1629 }
1631 IT cnz = SpHelper::SpCartesian< SR > (*dcsc, *(Btrans.dcsc), kisect, isect1, isect2, multstack);
1633
1635 if(cnz > 0)
1636 {
1638 delete [] multstack;
1639 }
1640 return SpDCCols< IT,T_promote > (m, rhs.n, mydcsc);
1641}
1642
1643
1644template <class IT, class NT>
1645template <typename SR, typename NTR>
1646SpDCCols< IT, typename promote_trait<NT,NTR>::T_promote > SpDCCols<IT,NT>::OrdColByCol(const SpDCCols<IT,NTR> & rhs) const
1647{
1648 typedef typename promote_trait<NT,NTR>::T_promote T_promote;
1649
1650 if(isZero() || rhs.isZero())
1651 {
1652 return SpDCCols<IT, T_promote> (0, m, rhs.n, 0); // return an empty matrix
1653 }
1655 IT cnz = SpHelper::SpColByCol< SR > (*dcsc, *(rhs.dcsc), n, multstack);
1656
1658 if(cnz > 0)
1659 {
1661 delete [] multstack;
1662 }
1663 return SpDCCols< IT,T_promote > (m, rhs.n, mydcsc);
1664}
1665
1666}
int64_t IT
double NT
void convert(string ifname, string ofname, string sort="revsize")
Definition test.cpp:53
SpDCCols< IT, NT > & operator=(const SpDCCols< IT, NT > &rhs)
Definition SpDCCols.cpp:365
friend class SpDCCols
Definition SpDCCols.h:369
void PrintInfo() const
SpDCCols< IT, NT > * PruneI(_UnaryOperation __unary_op, bool inPlace, GlobalIT rowOffset, GlobalIT colOffset)
Definition SpDCCols.cpp:431
std::ofstream & put(std::ofstream &outfile) const
std::ifstream & get(std::ifstream &infile)
void ColSplit(int parts, std::vector< SpDCCols< IT, NT > > &matrices)
Definition SpDCCols.cpp:936
int PlusEq_AnXBn(const SpDCCols< IT, NT > &A, const SpDCCols< IT, NT > &B)
void PruneColumnByIndex(const std::vector< IT > &ci)
Definition SpDCCols.cpp:573
void ColConcatenate(std::vector< SpDCCols< IT, NT > > &matrices)
SpDCCols< IT, NT > * Prune(_UnaryOperation __unary_op, bool inPlace)
Definition SpDCCols.cpp:478
Dcsc< IT, NT > * dcsc
Definition SpDCCols.h:358
static const IT esscount
Definition SpDCCols.h:303
std::vector< IT > GetEssentials() const
Definition SpDCCols.cpp:787
void Merge(SpDCCols< IT, NT > &partA, SpDCCols< IT, NT > &partB)
int PlusEq_AnXBt(const SpDCCols< IT, NT > &A, const SpDCCols< IT, NT > &B)
SpDCCols< IT, NT > operator()(IT ri, IT ci) const
SpDCCols< IT, NT > & operator+=(const SpDCCols< IT, NT > &rhs)
Definition SpDCCols.cpp:394
SpDCCols< IT, NT > TransposeConst() const
Const version, doesn't touch the existing object.
Definition SpDCCols.cpp:877
void SetDifference(const SpDCCols< IT, NT > &rhs)
Definition SpDCCols.cpp:631
void Transpose()
Mutator version, replaces the calling object.
Definition SpDCCols.cpp:854
int PlusEq_AtXBt(const SpDCCols< IT, NT > &A, const SpDCCols< IT, NT > &B)
void EWiseScale(NT **scaler, IT m_scaler, IT n_scaler)
Definition SpDCCols.cpp:702
friend SpDCCols< IU, typename promote_trait< NU1, NU2 >::T_promote > EWiseMult(const SpDCCols< IU, NU1 > &A, const SpDCCols< IU, NU2 > &B, bool exclude)
Definition Friends.h:1011
void CreateImpl(const std::vector< IT > &essentials)
Definition SpDCCols.cpp:734
SpDCCols< IT, NT > * TransposeConstPtr() const
Definition SpDCCols.cpp:891
Arr< IT, NT > GetArrays() const
Definition SpDCCols.cpp:826
SpDCCols< IT, NT > * PruneColumn(NT *pvals, _BinaryOperation __binary_op, bool inPlace)
Definition SpDCCols.cpp:527
int PlusEq_AtXBn(const SpDCCols< IT, NT > &A, const SpDCCols< IT, NT > &B)
void Split(SpDCCols< IT, NT > &partA, SpDCCols< IT, NT > &partB)
Definition SpDCCols.cpp:906
static void SpIntersect(const Dcsc< IT, NT1 > &Adcsc, const Dcsc< IT, NT2 > &Bdcsc, Isect< IT > *&cols, Isect< IT > *&rows, Isect< IT > *&isect1, Isect< IT > *&isect2, Isect< IT > *&itr1, Isect< IT > *&itr2)
Definition SpHelper.h:346
static void deallocate2D(T **array, I m)
Definition SpHelper.h:249
int size
Definition common.h:20
int rank
#define PRINT_LIMIT
Definition SpDefs.h:63
void DeleteAll(A arr1)
Definition Deleter.h:48
double A