PEXSI
 All Classes Namespaces Files Functions Variables Typedefs Pages
pselinv.hpp
Go to the documentation of this file.
1 /*
2  Copyright (c) 2012 The Regents of the University of California,
3  through Lawrence Berkeley National Laboratory.
4 
5 Authors: Lin Lin and Mathias Jacquelin
6 
7 This file is part of PEXSI. All rights reserved.
8 
9 Redistribution and use in source and binary forms, with or without
10 modification, are permitted provided that the following conditions are met:
11 
12 (1) Redistributions of source code must retain the above copyright notice, this
13 list of conditions and the following disclaimer.
14 (2) Redistributions in binary form must reproduce the above copyright notice,
15 this list of conditions and the following disclaimer in the documentation
16 and/or other materials provided with the distribution.
17 (3) Neither the name of the University of California, Lawrence Berkeley
18 National Laboratory, U.S. Dept. of Energy nor the names of its contributors may
19 be used to endorse or promote products derived from this software without
20 specific prior written permission.
21 
22 THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS "AS IS" AND
23 ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE IMPLIED
24 WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE ARE
25 DISCLAIMED. IN NO EVENT SHALL THE COPYRIGHT OWNER OR CONTRIBUTORS BE LIABLE FOR
26 ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL, EXEMPLARY, OR CONSEQUENTIAL DAMAGES
27 (INCLUDING, BUT NOT LIMITED TO, PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES;
28 LOSS OF USE, DATA, OR PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON
29 ANY THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT
30 (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
31 SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
32 
33 You are under no obligation whatsoever to provide any bug fixes, patches, or
34 upgrades to the features, functionality or performance of the source code
35 ("Enhancements") to anyone; however, if you choose to make your Enhancements
36 available either publicly, or directly to Lawrence Berkeley National
37 Laboratory, without imposing a separate written license agreement for such
38 Enhancements, then you hereby grant the following license: a non-exclusive,
39 royalty-free perpetual license to install, use, modify, prepare derivative
40 works, incorporate into other computer software, distribute, and sublicense
41 such enhancements or derivative works thereof, in binary and source code form.
42  */
46 #ifndef _PEXSI_PSELINV_HPP_
47 #define _PEXSI_PSELINV_HPP_
48 
49 // *********************************************************************
50 // Common utilities
51 // *********************************************************************
52 
53 #include "pexsi/environment.hpp"
54 #include "pexsi/NumVec.hpp"
55 #include "pexsi/NumMat.hpp"
56 #include "pexsi/sparse_matrix.hpp"
57 
59 #include "pexsi/mpi_interf.hpp"
60 #include "pexsi/utility.hpp"
61 #include "pexsi/blas.hpp"
62 #include "pexsi/lapack.hpp"
63 
64 #include "pexsi/TreeBcast.hpp"
65 
66 
67 #include <set>
68 
69 
70 //#define IDX_TO_TAG(lidx,tag) (SELINV_TAG_COUNT*(lidx)+(tag))
71 #define IDX_TO_TAG(lidx,tag,max) ((SELINV_TAG_COUNT*(lidx%(max+1))+(tag)))
72 #define IDX_TO_TAG2(sidx,lidx,tag) (SELINV_TAG_COUNT*(sidx)+(tag))
73 #define TAG_TO_IDX(tag,typetag) (((tag)-(typetag))/SELINV_TAG_COUNT)
74 
75 
76 //#define LIST_BARRIER
77 //#define ALL_BARRIER
78 
79 namespace PEXSI{
80 
81 enum MSGTYPE {LSIZE=0,LROWSIZE,USIZE,UCOLSIZE,LCONTENT,LROWCONTENT,UCONTENT,UCOLCONTENT,MSGCOUNT};
82 
83 struct ULComparator {
84  IntNumVec & lookup;
85  ULComparator(IntNumVec & v):lookup(v){}
86  bool operator() (int i,int j) { return ((lookup)(i)<(lookup)(j));}
87 };
88 
89 
90 
91 typedef std::vector<bool> bitMask;
92 typedef std::map<bitMask , std::vector<Int> > bitMaskSet;
93 
94 
95 
96 
97 
105 
106  // Member functions to setup the default value
108 };
109 
110 
116  std::string ColPerm;
117  std::string RowPerm;
118  Int Symmetric;
119 
120  // Member functions to setup the default value
121  FactorizationOptions(): Symmetric(1) {}
122 };
123 
124 
125 
126 
127 
128 /**********************************************************************
129  * Basic PSelInv data structure
130  **********************************************************************/
131 
140 struct GridType{
141  // Data
142  MPI_Comm comm;
143  MPI_Comm rowComm;
144  MPI_Comm colComm;
145  Int mpirank;
146  Int mpisize;
147  Int numProcRow;
148  Int numProcCol;
149 
150  // Member function
151  GridType( MPI_Comm Bcomm, int nprow, int npcol );
152  ~GridType();
153 };
154 
178  IntNumVec perm;
179  IntNumVec permInv;
180  IntNumVec perm_r;
181  IntNumVec permInv_r;
182  IntNumVec superIdx;
183  IntNumVec superPtr;
184  IntNumVec etree;
185 };
186 
187 
192 template<typename T>
193 struct LBlock{
194  // Variables
196  Int blockIdx;
197 
199  Int numRow;
200 
202  Int numCol;
203 
206 
207 
210 
211  // Member functions;
212  LBlock() {
213  blockIdx = -1; numRow = 0; numCol =0;
214  nzval.Resize(0,0);
215  }
216  ~LBlock() {}
217  LBlock& operator = (const LBlock& LB) {
218  blockIdx = LB.blockIdx;
219  numRow = LB.numRow;
220  numCol = LB.numCol;
221  rows = LB.rows;
222  nzval = LB.nzval;
223  return *this;
224  }
225  friend std::ostream& operator<<(std::ostream& out, const LBlock& vec) // output
226  {
227  out << "(" << vec.blockIdx << ", " << vec.numRow << ", " << vec.numCol <<std::endl<< "rows " << vec.rows <<std::endl<< "nzval " <<std::endl<< vec.nzval << ")";
228  return out;
229  }
230 
231 
232 };
233 
244 template<typename T>
245 struct UBlock{
246  // Variables
248  Int blockIdx;
249 
251  Int numRow;
252 
254  Int numCol;
255 
258 
261 
262  // Member functions;
263  UBlock() {
264  blockIdx = -1; numRow = 0; numCol =0;
265  }
266  ~UBlock() {}
267  UBlock& operator = (const UBlock& UB) {
268  blockIdx = UB.blockIdx;
269  numRow = UB.numRow;
270  numCol = UB.numCol;
271  cols = UB.cols;
272  nzval = UB.nzval;
273  return *this;
274  }
275 
276  friend std::ostream& operator<<(std::ostream& out, const UBlock& vec) // output
277  {
278  out << "(" << vec.blockIdx << ", " << vec.numRow << ", " << vec.numCol <<std::endl<< "cols " << vec.cols <<std::endl<< "nzval " <<std::endl<< vec.nzval << ")";
279  return out;
280  }
281 };
282 
283 // *********************************************************************
284 // SuperLU style utility functions
285 //
286 // The SuperLU style macros are defined here as inline functions
287 // so that the code is more portable.
288 // *********************************************************************
289 
291 inline Int MYPROC( const GridType* g )
292 { return g->mpirank; }
293 
295 inline Int MYROW( const GridType* g )
296 { return g->mpirank / g->numProcCol; }
297 
299 inline Int MYCOL( const GridType* g )
300 { return g->mpirank % g->numProcCol; }
301 
304 inline Int PROW( Int bnum, const GridType* g )
305 { return bnum % g->numProcRow; }
306 
309 inline Int PCOL( Int bnum, const GridType* g )
310 { return bnum % g->numProcCol; }
311 
314 inline Int PNUM( Int i, Int j, const GridType* g )
315 { return (i%g->numProcRow) * g->numProcCol + j%g->numProcCol; }
316 
319 inline Int LBi( Int bnum, const GridType* g )
320 { return bnum / g->numProcRow; }
321 
324 inline Int LBj( Int bnum, const GridType* g)
325 { return bnum / g->numProcCol; }
326 
329 inline Int GBi( Int iLocal, const GridType* g )
330 { return iLocal * g->numProcRow + MYROW( g ); }
331 
334 inline Int GBj( Int jLocal, const GridType* g )
335 { return jLocal * g->numProcCol + MYCOL( g ); }
336 
339 inline Int CEILING( Int a, Int b )
340 { return (a%b) ? ( a/b + 1 ) : ( a/b ); }
341 
343 inline Int BlockIdx( Int i, const SuperNodeType *s )
344 { return s->superIdx[i]; }
345 
348 inline Int FirstBlockCol( Int bnum, const SuperNodeType *s )
349 { return s->superPtr[bnum]; }
350 
351 
355 inline Int FirstBlockRow( Int bnum, const SuperNodeType *s )
356 { return s->superPtr[bnum]; }
357 
358 
360 inline Int SuperSize( Int bnum, const SuperNodeType *s )
361 { return s->superPtr[bnum+1] - s->superPtr[bnum]; }
362 
364 inline Int NumSuper( const SuperNodeType *s )
365 { return s->superPtr.m() - 1; }
366 
369 inline Int NumCol( const SuperNodeType *s )
370 { return s->superIdx.m(); }
371 
372 
373 // *********************************************************************
374 // Serialize / Deserialize
375 // *********************************************************************
376 
377 // L part
378 
394 namespace LBlockMask{
395 enum {
396  BLOCKIDX,
397  NUMROW,
398  NUMCOL,
399  ROWS,
400  NZVAL,
401  TOTAL_NUMBER
402 };
403 }
404 
405 
406 template<typename T>
407 Int inline serialize(LBlock<T>& val, std::ostream& os, const std::vector<Int>& mask){
408  if(mask[LBlockMask::BLOCKIDX]==1) serialize(val.blockIdx, os, mask);
409  if(mask[LBlockMask::NUMROW ]==1) serialize(val.numRow, os, mask);
410  if(mask[LBlockMask::NUMCOL ]==1) serialize(val.numCol, os, mask);
411  if(mask[LBlockMask::ROWS ]==1) serialize(val.rows, os, mask);
412  if(mask[LBlockMask::NZVAL ]==1) serialize(val.nzval, os, mask);
413  return 0;
414 }
415 
416 template<typename T>
417 Int inline deserialize(LBlock<T>& val, std::istream& is, const std::vector<Int>& mask){
418  if(mask[LBlockMask::BLOCKIDX]==1) deserialize(val.blockIdx, is, mask);
419  if(mask[LBlockMask::NUMROW ]==1) deserialize(val.numRow, is, mask);
420  if(mask[LBlockMask::NUMCOL ]==1) deserialize(val.numCol, is, mask);
421  if(mask[LBlockMask::ROWS ]==1) deserialize(val.rows, is, mask);
422  if(mask[LBlockMask::NZVAL ]==1) deserialize(val.nzval, is, mask);
423  return 0;
424 }
425 
426 // U part
442 namespace UBlockMask{
443 enum {
444  BLOCKIDX,
445  NUMROW,
446  NUMCOL,
447  COLS,
448  NZVAL,
449  TOTAL_NUMBER
450 };
451 }
452 
453 template<typename T>
454 Int inline serialize(UBlock<T>& val, std::ostream& os, const std::vector<Int>& mask){
455  Int i = 0;
456  if(mask[i]==1) serialize(val.blockIdx, os, mask); i++;
457  if(mask[i]==1) serialize(val.numRow, os, mask); i++;
458  if(mask[i]==1) serialize(val.numCol, os, mask); i++;
459  if(mask[i]==1) serialize(val.cols, os, mask); i++;
460  if(mask[i]==1) serialize(val.nzval, os, mask); i++;
461  return 0;
462 }
463 
464 template<typename T>
465 Int inline deserialize(UBlock<T>& val, std::istream& is, const std::vector<Int>& mask){
466  Int i = 0;
467  if(mask[i]==1) deserialize(val.blockIdx, is, mask); i++;
468  if(mask[i]==1) deserialize(val.numRow, is, mask); i++;
469  if(mask[i]==1) deserialize(val.numCol, is, mask); i++;
470  if(mask[i]==1) deserialize(val.cols, is, mask); i++;
471  if(mask[i]==1) deserialize(val.nzval, is, mask); i++;
472  return 0;
473 }
474 
475 /**********************************************************************
476  * Main data structure in PSelInv: PMatrix
477  **********************************************************************/
478 
540 
541 template<typename T>
542 class PMatrix{
543 public:
546 
547  static PMatrix<T> * Create(const GridType * pGridType, const SuperNodeType * pSuper, const PSelInvOptions * pSelInvOpt , const FactorizationOptions * pFactOpt);
548  static PMatrix<T> * Create(const FactorizationOptions * pFactOpt);
549 
550 public:
551  // This is the tag used for mpi communication for selinv
552 
553  enum{
554  SELINV_TAG_U_SIZE,
555  SELINV_TAG_U_CONTENT,
556  SELINV_TAG_L_SIZE,
557  SELINV_TAG_L_CONTENT,
558  SELINV_TAG_L_REDUCE,
559  SELINV_TAG_D_SIZE,
560  SELINV_TAG_D_CONTENT,
561  SELINV_TAG_D_REDUCE,
562  SELINV_TAG_L_SIZE_CD,
563  SELINV_TAG_L_CONTENT_CD,
564  SELINV_TAG_COUNT
565  };
566 
567 
568  std::vector<std::vector<Int> > ColBlockIdx_;
569  std::vector<std::vector<Int> > RowBlockIdx_;
570 protected:
571  // *********************************************************************
572  // Variables
573  // *********************************************************************
574  // Data variables
575 
576  const GridType* grid_;
577 
578  const SuperNodeType* super_;
579 
580  const PSelInvOptions * options_;
581  const FactorizationOptions * optionsFact_;
582 
583  Int limIndex_;
584  Int maxTag_;
585 
586  std::vector<std::vector<LBlock<T> > > L_;
587  std::vector<std::vector<UBlock<T> > > U_;
588 
589  std::vector<std::vector<Int> > workingSet_;
590 
591 
592  // Communication variables
593  BolNumMat isSendToBelow_;
594  BolNumMat isSendToRight_;
595  BolNumVec isSendToDiagonal_;
596  BolNumMat isSendToCrossDiagonal_;
597 
598  BolNumMat isRecvFromBelow_;
599  BolNumVec isRecvFromAbove_;
600  BolNumVec isRecvFromLeft_;
601  BolNumMat isRecvFromCrossDiagonal_;
602 
603 
604 #ifdef NEW_BCAST
605  std::vector<TreeBcast2<T> *> fwdToBelowTree2_;
606  std::vector<TreeBcast2<T> *> fwdToRightTree2_;
607 #endif
608  std::vector<TreeBcast *> fwdToBelowTree_;
609  std::vector<TreeBcast *> fwdToRightTree_;
610  std::vector<TreeReduce<T> *> redToLeftTree_;
611  std::vector<TreeReduce<T> *> redToAboveTree_;
612 
613 
614 
616  NumMat<T> LUpdateBuf;
617  NumMat<T> DiagBuf;
618  std::vector<Int> RowLocalPtr;
619  std::vector<Int> BlockIdxLocal;
620  std::vector<char> SstrLcolSend;
621  std::vector<char> SstrUrowSend;
622  std::vector<char> SstrLcolRecv;
623  std::vector<char> SstrUrowRecv;
624  Int SizeSstrLcolSend;
625  Int SizeSstrUrowSend;
626  Int SizeSstrLcolRecv;
627  Int SizeSstrUrowRecv;
628  Int Index;
629  Int Rank;
630  Int isReady;
631 
632 
634  SizeSstrLcolSend(0),
635  SizeSstrUrowSend(0),
636  SizeSstrLcolRecv(0),
637  SizeSstrUrowRecv(0),
638  Index(0),
639  Rank(0),
640  isReady(0){}
641 
642  SuperNodeBufferType(Int &pIndex) :
643  SizeSstrLcolSend(0),
644  SizeSstrUrowSend(0),
645  SizeSstrLcolRecv(0),
646  SizeSstrUrowRecv(0),
647  Index(pIndex),
648  Rank(0),
649  isReady(0) {}
650 
651  };
652 
653 
655  inline void SelInvIntra_P2p(Int lidx,Int & rank);
656 
658  inline void SelInv_lookup_indexes(SuperNodeBufferType & snode, std::vector<LBlock<T> > & LcolRecv, std::vector<UBlock<T> > & UrowRecv, NumMat<T> & AinvBuf,NumMat<T> & UBuf);
659 
661  inline void GetWorkSet(std::vector<Int> & snodeEtree, std::vector<std::vector<Int> > & WSet);
662 
664  inline void UnpackData(SuperNodeBufferType & snode, std::vector<LBlock<T> > & LcolRecv, std::vector<UBlock<T> > & UrowRecv);
665 
667  inline void ComputeDiagUpdate(SuperNodeBufferType & snode);
668 
670  inline void SendRecvCD_UpdateU(std::vector<SuperNodeBufferType > & arrSuperNodes, Int stepSuper);
671 
672 public:
673  // *********************************************************************
674  // Public member functions
675  // *********************************************************************
676 
677  PMatrix();
678 
679  PMatrix( const GridType* g, const SuperNodeType* s, const PEXSI::PSelInvOptions * o, const PEXSI::FactorizationOptions * oFact );
680 
681  void deallocate();
682 
683 
684  virtual ~PMatrix();
685  PMatrix( const PMatrix & C);
686  PMatrix & operator = ( const PMatrix & C);
687 
688  void Setup( const GridType* g, const SuperNodeType* s, const PEXSI::PSelInvOptions * o, const PEXSI::FactorizationOptions * oFact );
689 
690  Int NumCol() const { return super_ -> superIdx.m(); }
691 
692  Int NumSuper() const { return super_ ->superPtr.m() - 1; }
693 
695  Int NumLocalBlockCol() const { return CEILING( this->NumSuper(), grid_->numProcCol ); }
696 
698  Int NumLocalBlockRow() const { return CEILING( this->NumSuper(), grid_->numProcRow); }
699 
700 
701  std::vector< std::vector<Int> > & ColBlockIdx() { return ColBlockIdx_; }
702  std::vector< std::vector<Int> > & RowBlockIdx() { return RowBlockIdx_; }
703  std::vector<Int> & ColBlockIdx(Int jLocal) { return ColBlockIdx_[jLocal]; }
704  std::vector<Int> & RowBlockIdx(Int iLocal) { return RowBlockIdx_[iLocal]; }
705 
706 
709  Int NumBlockL( Int jLocal ) const { return L_[jLocal].size(); }
710 
713  Int NumBlockU( Int iLocal ) const { return U_[iLocal].size(); }
714 
716  const GridType* Grid() const { return grid_; }
717 
720  const SuperNodeType* SuperNode() const { return super_; }
721 
724  std::vector<LBlock<T> >& L( Int jLocal ) { return L_[jLocal]; }
725 
728  std::vector<UBlock<T> >& U( Int iLocal ) { return U_[iLocal]; }
729 
732  std::vector<std::vector<int> >& WorkingSet( ) { return workingSet_; }
733  //void WorkingSet(const std::vector<std::vector<int> > * pWSet,const std::vector<std::vector<int> > * pWRanks) { pWSet = &workingSet_; pWRanks = &workingRanks_; }
734 
737  Int CountSendToRight(Int ksup) { Int count= std::count (isSendToRight_.VecData(ksup), isSendToRight_.VecData(ksup) + grid_->numProcCol, true); return (isSendToRight_(MYCOL(grid_),ksup)?count-1:count); }
738 
741  Int CountSendToBelow(Int ksup) { Int count= std::count (isSendToBelow_.VecData(ksup), isSendToBelow_.VecData(ksup) + grid_->numProcRow, true); return (isSendToBelow_(MYROW(grid_),ksup)?count-1:count); }
742 
745  Int CountRecvFromBelow(Int ksup) { Int count= std::count (isRecvFromBelow_.VecData(ksup), isRecvFromBelow_.VecData(ksup) + grid_->numProcRow, true); return (isRecvFromBelow_(MYROW(grid_),ksup)?count-1:count); }
746 
749  Int CountSendToCrossDiagonal(Int ksup) { Int count= std::count (isSendToCrossDiagonal_.VecData(ksup), isSendToCrossDiagonal_.VecData(ksup) + grid_->numProcCol, true); return ((isSendToCrossDiagonal_(MYCOL(grid_),ksup) && MYROW(grid_)==PROW(ksup,grid_))?count-1:count); }
750 
753  Int CountRecvFromCrossDiagonal(Int ksup) { Int count= std::count (isRecvFromCrossDiagonal_.VecData(ksup), isRecvFromCrossDiagonal_.VecData(ksup) + grid_->numProcRow, true); return ((isRecvFromCrossDiagonal_(MYROW(grid_),ksup) && MYCOL(grid_)==PCOL(ksup,grid_))?count-1:count); }
754 
755 
756 
757 
760  void GetEtree(std::vector<Int> & etree_supno );
761 
762 
767  virtual void ConstructCommunicationPattern( );
768 
769 
774 
775 
776 
805  virtual void PreSelInv( );
806 
931  virtual void SelInv( );
932 
934  void SelInv_P2p( );
935 
936 
943  void GetDiagonal( NumVec<T>& diag );
944  void GetColumn ( Int colIdx, NumVec<T>& col );
945 
946 
953 
964  const DistSparseMatrix<T>& A,
965  DistSparseMatrix<T>& B );
966 
967 
978  const DistSparseMatrix<T>& A,
979  DistSparseMatrix<T>& B );
980 
981 
984  Int NnzLocal();
985 
988  LongInt Nnz();
989 
993  void GetNegativeInertia ( Real& inertia );
994 
995  // inline int IdxToTag(Int lidx, Int tag) { return SELINV_TAG_COUNT*(lidx)+(tag);}
996 
997  void DumpLU();
998 
999  void CopyLU( const PMatrix & C);
1000  inline int IdxToTag(Int lidx, Int tag) { return SELINV_TAG_COUNT*(lidx)+(tag);}
1001 };
1002 
1003 template<typename T> class PMatrixUnsym;
1004 
1005 
1006 } // namespace PEXSI
1007 
1008 
1009 #include "pexsi/pselinv_impl.hpp"
1010 
1011 
1012 #endif //_PEXSI_PSELINV_HPP_
virtual void PreSelInv()
PreSelInv prepares the structure in L_ and U_ so that SelInv only involves matrix-matrix multiplicati...
Definition: pselinv_impl.hpp:3594
Environmental variables.
A thin interface for passing parameters to set the Factorization options.
Definition: pselinv.hpp:115
Int blockIdx
Block index (supernodal index)
Definition: pselinv.hpp:248
void PMatrixToDistSparseMatrix(DistSparseMatrix< T > &A)
PMatrixToDistSparseMatrix converts the PMatrix into a distributed compressed sparse column matrix for...
Definition: pselinv_impl.hpp:4083
std::vector< LBlock< T > > & L(Int jLocal)
L returns the vector of nonzero L blocks for the local block column jLocal.
Definition: pselinv.hpp:724
SuperNodeType describes mapping between supernode and column, the permutation information, and potentially the elimination tree (not implemented here).
Definition: pselinv.hpp:177
Interface with SuperLU_Dist (version 3.0 and later)
void GetNegativeInertia(Real &inertia)
GetNegativeInertia computes the negative inertia of a PMatrix. This can be used to estimate e...
Int PROW(Int bnum, const GridType *g)
PROW returns the processor row that the bnum-th block (supernode) belongs to.
Definition: pselinv.hpp:304
Thin interface to LAPACK.
Int numRow
Number of nonzero rows.
Definition: pselinv.hpp:199
void UnpackData(SuperNodeBufferType &snode, std::vector< LBlock< T > > &LcolRecv, std::vector< UBlock< T > > &UrowRecv)
UnpackData.
Definition: pselinv_impl.hpp:1060
Int maxPipelineDepth
The maximum pipeline depth.
Definition: pselinv.hpp:104
Int MYPROC(const GridType *g)
MYPROC returns the current processor rank.
Definition: pselinv.hpp:291
Int NumLocalBlockCol() const
NumLocalBlockCol returns the total number of block columns.
Definition: pselinv.hpp:695
Int PCOL(Int bnum, const GridType *g)
PCOL returns the processor column that the bnum-th block (supernode) belongs to.
Definition: pselinv.hpp:309
Int LBi(Int bnum, const GridType *g)
LBi returns the local block number on the processor at processor row PROW( bnum, g )...
Definition: pselinv.hpp:319
Int SuperSize(Int bnum, const SuperNodeType *s)
SuperSize returns the size of the block bnum.
Definition: pselinv.hpp:360
Int LBj(Int bnum, const GridType *g)
LBj returns the local block number on the processor at processor column PCOL( bnum, g ).
Definition: pselinv.hpp:324
Int GBj(Int jLocal, const GridType *g)
GBj returns the global block number from a local block number in the column direction.
Definition: pselinv.hpp:334
Int numCol
Number of nonzero columns.
Definition: pselinv.hpp:202
static PMatrix< T > * Create(const GridType *pGridType, const SuperNodeType *pSuper, const PSelInvOptions *pSelInvOpt, const FactorizationOptions *pFactOpt)
Create is a factory method which returns a pointer either to a new PMatrix or to a new PMatrixUnsym d...
Definition: pselinv_impl.hpp:5057
Int NumCol(const SuperNodeType *s)
NumCol returns the total number of columns for a supernodal partiiton.
Definition: pselinv.hpp:369
Definition: pselinv.hpp:615
Int CountRecvFromCrossDiagonal(Int ksup)
CountRecvFromCrossDiagonal returns the number of cross diagonal processors with which current process...
Definition: pselinv.hpp:753
void GetDiagonal(NumVec< T > &diag)
GetDiagonal extracts the diagonal elements of the PMatrix.
Definition: pselinv_impl.hpp:3961
std::vector< std::vector< int > > & WorkingSet()
WorkingSet returns the ordered list of supernodes which could be done in parallel.
Definition: pselinv.hpp:732
Implementation of the parallel SelInv.
void GetWorkSet(std::vector< Int > &snodeEtree, std::vector< std::vector< Int > > &WSet)
GetWorkSet.
Definition: pselinv_impl.hpp:1227
Int FirstBlockRow(Int bnum, const SuperNodeType *s)
FirstBlockRow returns the first column of a block bnum. Note: the functionality of FirstBlockRow is e...
Definition: pselinv.hpp:355
Definition: pselinv.hpp:83
void GetEtree(std::vector< Int > &etree_supno)
GetEtree computes the supernodal elimination tree to be used later in the pipelined selected inversio...
Definition: pselinv_impl.hpp:1166
const GridType * Grid() const
Grid returns the GridType structure of the current PMatrix.
Definition: pselinv.hpp:716
LBlock stores a nonzero block in the lower triangular part or the diagonal part in PSelInv...
Definition: pselinv.hpp:193
Int PNUM(Int i, Int j, const GridType *g)
PNUM returns the processor rank that the bnum-th block (supernode) belongs to.
Definition: pselinv.hpp:314
virtual void ConstructCommunicationPattern()
ConstructCommunicationPattern constructs the communication pattern to be used later in the selected i...
Definition: pselinv_impl.hpp:2603
Sparse matrix and Distributed sparse matrix in compressed column format.
void ConstructCommunicationPattern_P2p()
ConstructCommunicationPattern_P2p constructs the communication pattern to be used later in the select...
Definition: pselinv_impl.hpp:2611
Int GBi(Int iLocal, const GridType *g)
GBi returns the global block number from a local block number in the row direction.
Definition: pselinv.hpp:329
GridType is the PSelInv way of defining the grid.
Definition: pselinv.hpp:140
void SelInvIntra_P2p(Int lidx, Int &rank)
SelInvIntra_P2p.
Definition: pselinv_impl.hpp:1465
Numerical matrix.
Interface with MPI to facilitate communication.
Int MYCOL(const GridType *g)
MYCOL returns my processor column.
Definition: pselinv.hpp:299
Int CountSendToCrossDiagonal(Int ksup)
CountSendToCrossDiagonal returns the number of cross diagonal processors with which current processor...
Definition: pselinv.hpp:749
Int NumLocalBlockRow() const
NumLocalBlockRow returns the total number of block rows.
Definition: pselinv.hpp:698
Various utility subroutines.
const SuperNodeType * SuperNode() const
SuperNode returns the supernodal partition of the current PMatrix.
Definition: pselinv.hpp:720
void PMatrixToDistSparseMatrix_OLD(const DistSparseMatrix< T > &A, DistSparseMatrix< T > &B)
PMatrixToDistSparseMatrix_OLD converts the PMatrix into a distributed compressed sparse column matrix...
Definition: pselinv_impl.hpp:4355
IntNumVec rows
Dimension numRow * 1, index (0-based) for the number of nonzero rows.
Definition: pselinv.hpp:205
Int numCol
Number of nonzero columns.
Definition: pselinv.hpp:254
UBlock stores a nonzero block in the upper triangular part in PSelInv.
Definition: pselinv.hpp:245
std::vector< UBlock< T > > & U(Int iLocal)
U returns the vector of nonzero U blocks for the local block row iLocal.
Definition: pselinv.hpp:728
Int FirstBlockCol(Int bnum, const SuperNodeType *s)
FirstBlockCol returns the first column of a block bnum.
Definition: pselinv.hpp:348
Int NnzLocal()
NnzLocal computes the number of nonzero elements (L and U) saved locally.
Definition: pselinv_impl.hpp:4965
Int blockIdx
Block index (supernodal index)
Definition: pselinv.hpp:196
Thin interface to BLAS.
PMatrix contains the main data structure and the computational routine for the parallel selected inve...
Definition: pselinv.hpp:542
void SendRecvCD_UpdateU(std::vector< SuperNodeBufferType > &arrSuperNodes, Int stepSuper)
SendRecvCD_UpdateU.
Definition: pselinv_impl.hpp:828
Int BlockIdx(Int i, const SuperNodeType *s)
BlockIdx returns the block index of a column i.
Definition: pselinv.hpp:343
Int numRow
Number of nonzero rows.
Definition: pselinv.hpp:251
Int CountSendToRight(Int ksup)
CountSendToRight returns the number of processors to the right of current processor with which it has...
Definition: pselinv.hpp:737
Int CountSendToBelow(Int ksup)
CountSendToBelow returns the number of processors below current processor with which it has to commun...
Definition: pselinv.hpp:741
Int CEILING(Int a, Int b)
CEILING is used for computing the storage space for local number of blocks.
Definition: pselinv.hpp:339
NumMat< T > nzval
Dimension numRow * numCol, nonzero elements.
Definition: pselinv.hpp:260
Int NumSuper(const SuperNodeType *s)
NumSuper returns the total number of supernodes.
Definition: pselinv.hpp:364
Numerical vector.
void SelInv_P2p()
Point-to-point version of the selected inversion.
Definition: pselinv_impl.hpp:3537
Int CountRecvFromBelow(Int ksup)
CountRecvFromBelow returns the number of processors below the current processor from which it receive...
Definition: pselinv.hpp:745
Int NumBlockL(Int jLocal) const
NumBlockL returns the number of nonzero L blocks for the local block column jLocal.
Definition: pselinv.hpp:709
NumMat< T > nzval
Dimension numRow * numCol, nonzero elements.
Definition: pselinv.hpp:209
A thin interface for passing parameters to set the PSelInv options.
Definition: pselinv.hpp:102
IntNumVec cols
Dimension numCol * 1, index (0-based) for the number of nonzero rows.
Definition: pselinv.hpp:257
Int NumBlockU(Int iLocal) const
NumBlockU returns the number of nonzero U blocks for the local block row iLocal.
Definition: pselinv.hpp:713
virtual void SelInv()
SelInv is the main function for the selected inversion.
Definition: pselinv_impl.hpp:3484
LongInt Nnz()
Nnz computes the total number of nonzero elements in the PMatrix.
Definition: pselinv_impl.hpp:4990
void SelInv_lookup_indexes(SuperNodeBufferType &snode, std::vector< LBlock< T > > &LcolRecv, std::vector< UBlock< T > > &UrowRecv, NumMat< T > &AinvBuf, NumMat< T > &UBuf)
SelInv_lookup_indexes.
Definition: pselinv_impl.hpp:471
void ComputeDiagUpdate(SuperNodeBufferType &snode)
ComputeDiagUpdate.
Definition: pselinv_impl.hpp:1126
DistSparseMatrix describes a Sparse matrix in the compressed sparse column format (CSC) and distribut...
Definition: sparse_matrix.hpp:91
PMatrixUnsym contains the main data structure and the computational routine for the parallel selected...
Definition: pselinv.hpp:1003
Int MYROW(const GridType *g)
MYROW returns my processor row.
Definition: pselinv.hpp:295