IFPACK  Development
Ifpack_SupportGraph.h
1 /*@HEADER
2 // ***********************************************************************
3 //
4 // Ifpack: Object-Oriented Algebraic Preconditioner Package
5 // Copyright (2002) Sandia Corporation
6 //
7 // Under terms of Contract DE-AC04-94AL85000, there is a non-exclusive
8 // license for use of this work by or on behalf of the U.S. Government.
9 //
10 // Redistribution and use in source and binary forms, with or without
11 // modification, are permitted provided that the following conditions are
12 // met:
13 //
14 // 1. Redistributions of source code must retain the above copyright
15 // notice, this list of conditions and the following disclaimer.
16 //
17 // 2. Redistributions in binary form must reproduce the above copyright
18 // notice, this list of conditions and the following disclaimer in the
19 // documentation and/or other materials provided with the distribution.
20 //
21 // 3. Neither the name of the Corporation nor the names of the
22 // contributors may be used to endorse or promote products derived from
23 // this software without specific prior written permission.
24 //
25 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
26 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
27 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
28 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
29 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
30 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
31 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
32 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
33 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
34 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
35 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
36 //
37 // Questions? Contact Michael A. Heroux (maherou@sandia.gov)
38 //
39 // ***********************************************************************
40 //@HEADER
41 */
42 
43 #ifndef IFPACK_SUPPORTGRAPH_H
44 #define IFPACK_SUPPORTGRAPH_H
45 
46 #include "Ifpack_ConfigDefs.h"
47 #include "Ifpack_Condest.h"
48 #include "Ifpack_Preconditioner.h"
49 #include "Ifpack_Amesos.h"
50 #include "Ifpack_Condest.h"
51 #include "Epetra_Map.h"
52 #include "Epetra_Comm.h"
53 #include "Epetra_Time.h"
54 #include "Epetra_Vector.h"
55 #include "Epetra_MultiVector.h"
56 #include "Epetra_LinearProblem.h"
57 #include "Epetra_RowMatrix.h"
58 #include "Epetra_CrsMatrix.h"
59 
60 #include "Teuchos_ParameterList.hpp"
61 #include "Teuchos_RefCountPtr.hpp"
62 
63 #include <boost/graph/adjacency_list.hpp>
64 #include <boost/graph/kruskal_min_spanning_tree.hpp>
65 #include <boost/graph/prim_minimum_spanning_tree.hpp>
66 #include <boost/config.hpp>
67 
68 using Teuchos::RefCountPtr;
69 using Teuchos::rcp;
70 typedef std::pair<int, int> E;
71 using namespace boost;
72 
73 typedef adjacency_list < vecS, vecS, undirectedS,
74  no_property, property < edge_weight_t, double > > Graph;
75 typedef graph_traits < Graph >::edge_descriptor Edge;
76 typedef graph_traits < Graph >::vertex_descriptor Vertex;
77 
78 
79 
80 template<typename T=Ifpack_Amesos> class Ifpack_SupportGraph :
81 public virtual Ifpack_Preconditioner
82 {
83 
84  public:
85 
87 
90 
92 
93 
95 
105  virtual int SetUseTranspose(bool UseTranspose_in);
106 
108 
109 
111 
113 
121  virtual int Apply(const Epetra_MultiVector& X, Epetra_MultiVector& Y) const;
122 
124 
135  virtual int ApplyInverse(const Epetra_MultiVector& X, Epetra_MultiVector& Y) const;
136 
138  virtual double NormInf() const {return(0.0);}
139 
141 
142 
144 
146  virtual const char * Label() const;
147 
149  virtual bool UseTranspose() const {return(UseTranspose_);};
150 
152  virtual bool HasNormInf() const {return(false);};
153 
155  virtual const Epetra_Comm & Comm() const {return(Matrix_->Comm());};
156 
158  virtual const Epetra_Map & OperatorDomainMap() const {return(Matrix_->OperatorDomainMap());};
159 
161  virtual const Epetra_Map & OperatorRangeMap() const {return(Matrix_->OperatorRangeMap());};
162 
164 
165 
167 
169  virtual bool IsInitialized() const
170  {
171  return(IsInitialized_);
172  }
173 
175  virtual bool IsComputed() const
176  {
177  return(IsComputed_);
178  }
179 
181 
193  virtual int SetParameters(Teuchos::ParameterList& List);
194 
196 
199  virtual int Initialize();
201 
204  virtual int Compute();
205 
207 
208 
210 
211 
213 
216  virtual double Condest(const Ifpack_CondestType CT = Ifpack_Cheap,
217  const int MaxIters = 1550,
218  const double Tol = 1e-9,
219  Epetra_RowMatrix* Matrix_in = 0);
220 
222  virtual double Condest() const
223  {
224  return(Condest_);
225  }
226 
228  virtual const Epetra_RowMatrix& Matrix() const
229  {
230  return(*Matrix_);
231  }
232 
234  virtual std::ostream& Print(std::ostream&) const;
235 
237  virtual int NumInitialize() const
238  {
239  return(NumInitialize_);
240  }
241 
243  virtual int NumCompute() const
244  {
245  return(NumCompute_);
246  }
247 
249  virtual int NumApplyInverse() const
250  {
251  return(NumApplyInverse_);
252  }
253 
255  virtual double InitializeTime() const
256  {
257  return(InitializeTime_);
258  }
259 
261  virtual double ComputeTime() const
262  {
263  return(ComputeTime_);
264  }
265 
267  virtual double ApplyInverseTime() const
268  {
269  return(ApplyInverseTime_);
270  }
271 
273  virtual double InitializeFlops() const
274  {
275  return(InitializeFlops_);
276  }
277 
279  virtual double ComputeFlops() const
280  {
281  return(ComputeFlops_);
282  }
283 
285  virtual double ApplyInverseFlops() const
286  {
287  return(ApplyInverseFlops_);
288  }
289 
290 
292 
293  protected:
294 
296  int FindSupport();
297 
299  Teuchos::RefCountPtr<const Epetra_RowMatrix> Matrix_;
300 
302  Teuchos::RefCountPtr<Epetra_CrsMatrix> Support_;
303 
305  std::string Label_;
306 
309 
312 
315 
317  Teuchos::ParameterList List_;
318 
320  double Condest_;
321 
324 
327 
329  mutable int NumApplyInverse_;
330 
333 
335  double ComputeTime_;
336 
338  mutable double ApplyInverseTime_;
339 
342 
345 
347  mutable double ApplyInverseFlops_;
348 
350  Teuchos::RefCountPtr<Epetra_Time> Time_;
351 
353  Teuchos::RefCountPtr<T> Inverse_;
354 
357 
359  double DiagPertRel_;
360 
362  double DiagPertAbs_;
363 
365  double KeepDiag_;
366 
369 
370 }; // class Ifpack_SupportGraph<T>
371 
372 
373 
374 //==============================================================================
375 template<typename T>
377 Matrix_(rcp(Matrix_in,false)),
378  IsInitialized_(false),
379  IsComputed_(false),
380  UseTranspose_(false),
381  Condest_(-1.0),
382  NumInitialize_(0),
383  NumCompute_(0),
384  NumApplyInverse_(0),
385  InitializeTime_(0.0),
386  ComputeTime_(0.0),
387  ApplyInverseTime_(0.0),
388  InitializeFlops_(0.0),
389  ComputeFlops_(0.0),
390  ApplyInverseFlops_(0.0),
391  NumForests_(1),
392  DiagPertRel_(1.0),
393  DiagPertAbs_(0.0),
394  KeepDiag_(1.0),
395  Randomize_(0)
396 {
397 
398  Teuchos::ParameterList List_in;
399  SetParameters(List_in);
400 }
401 //==============================================================================
402 template<typename T>
404 {
405 
406  // Extract matrix dimensions
407  long long rows = (*Matrix_).NumGlobalRows64();
408  long long cols = (*Matrix_).NumGlobalCols64();
409  int num_edges = ((*Matrix_).NumMyNonzeros() - (*Matrix_).NumMyDiagonals())/2;
410  std::cout << "global num rows " << rows << std::endl;
411 
412  // Assert square matrix
413  IFPACK_CHK_ERR((rows == cols));
414 
415  if(rows > std::numeric_limits<int>::max())
416  {
417  std::cerr << "Ifpack_SupportGraph<T>::FindSupport: global num rows won't fit an int. " << rows << std::endl;
418  IFPACK_CHK_ERR(1);
419  }
420 
421  // Rename for clarity
422  int num_verts = (int) rows;
423 
424  // Create data structures for the BGL code and temp data structures for extraction
425  E *edge_array = new E[num_edges];
426  double *weights = new double[num_edges];
427 
428  int num_entries;
429  int max_num_entries = (*Matrix_).MaxNumEntries();
430  double *values = new double[max_num_entries];
431  int *indices = new int[max_num_entries];
432 
433  double * diagonal = new double[num_verts];
434 
435 
436  for(int i = 0; i < max_num_entries; i++)
437  {
438  values[i]=0;
439  indices[i]=0;
440  }
441 
442  // Extract from the epetra matrix keeping only one edge per pair (assume symmetric)
443  int k = 0;
444  for(int i = 0; i < num_verts; i++)
445  {
446  (*Matrix_).ExtractMyRowCopy(i,max_num_entries,num_entries,values,indices);
447 
448  for(int j = 0; j < num_entries; j++)
449  {
450 
451  if(i == indices[j])
452  {
453  diagonal[i] = values[j];
454  // Diagonal pertubation, only if requested
455  if (DiagPertRel_)
456  diagonal[i] *= DiagPertRel_;
457  if (DiagPertAbs_)
458  diagonal[i] += DiagPertAbs_;
459  }
460 
461  if(i < indices[j])
462  {
463  edge_array[k] = E(i,indices[j]);
464  weights[k] = values[j];
465  if (Randomize_)
466  {
467  // Add small random pertubation.
468  weights[k] *= (1.0 + 1e-8 * drand48());
469  }
470 
471  k++;
472  }
473  }
474  }
475 
476  // Create BGL graph
477  Graph g(edge_array, edge_array + num_edges, weights, num_verts);
478 
479 
480  property_map < Graph, edge_weight_t >::type weight = get(edge_weight, g);
481 
482  std::vector < Edge > spanning_tree;
483 
484  // Run Kruskal, actually maximal weight ST since edges are negative
485  kruskal_minimum_spanning_tree(g, std::back_inserter(spanning_tree));
486 
487 
488  std::vector<int> NumNz(num_verts,1);
489 
490  //Find the degree of all the vertices
491  for (std::vector < Edge >::iterator ei = spanning_tree.begin();
492  ei != spanning_tree.end(); ++ei)
493  {
494  NumNz[source(*ei,g)] = NumNz[source(*ei,g)] + 1;
495  NumNz[target(*ei,g)] = NumNz[target(*ei,g)] + 1;
496  }
497 
498 
499  // Create an stl vector of stl vectors to hold indices and values (neighbour edges)
500  std::vector< std::vector< int > > Indices(num_verts);
501  // TODO: Optimize for performance, may use arrays instead of vectors
502  //std::vector<int> Indices[num_verts];
503  //std::vector<double> Values[num_verts];
504 
505  std::vector< std::vector< double > > Values(num_verts);
506 
507  for(int i = 0; i < num_verts; i++)
508  {
509  std::vector<int> temp(NumNz[i],0);
510  std::vector<double> temp2(NumNz[i],0);
511  Indices[i] = temp;
512  Values[i] = temp2;
513  }
514 
515  int *l = new int[num_verts];
516  for(int i = 0; i < num_verts; i++)
517  {
518  Indices[i][0] = i;
519  l[i] = 1;
520  }
521 
522  // Add each spanning forest (tree) to the support graph and
523  // remove it from original graph
524  for(int i = 0; i < NumForests_; i++)
525  {
526  if(i > 0)
527  {
528  spanning_tree.clear();
529  kruskal_minimum_spanning_tree(g,std::back_inserter(spanning_tree));
530  for(std::vector < Edge >::iterator ei = spanning_tree.begin();
531  ei != spanning_tree.end(); ++ei)
532  {
533  NumNz[source(*ei,g)] = NumNz[source(*ei,g)] + 1;
534  NumNz[target(*ei,g)] = NumNz[target(*ei,g)] + 1;
535  }
536  for(int i = 0; i < num_verts; i++)
537  {
538  Indices[i].resize(NumNz[i]);
539  Values[i].resize(NumNz[i]);
540  }
541  }
542 
543  for (std::vector < Edge >::iterator ei = spanning_tree.begin();
544  ei != spanning_tree.end(); ++ei)
545  {
546  // Assume standard Laplacian with constant row-sum.
547  // Edge weights are negative, so subtract to make diagonal positive
548  Indices[source(*ei,g)][0] = source(*ei,g);
549  Values[source(*ei,g)][0] = Values[source(*ei,g)][0] - weight[*ei];
550  Indices[target(*ei,g)][0] = target(*ei,g);
551  Values[target(*ei,g)][0] = Values[target(*ei,g)][0] - weight[*ei];
552 
553  Indices[source(*ei,g)][l[source(*ei,g)]] = target(*ei,g);
554  Values[source(*ei,g)][l[source(*ei,g)]] = weight[*ei];
555  l[source(*ei,g)] = l[source(*ei,g)] + 1;
556 
557  Indices[target(*ei,g)][l[target(*ei,g)]] = source(*ei,g);
558  Values[target(*ei,g)][l[target(*ei,g)]] = weight[*ei];
559  l[target(*ei,g)] = l[target(*ei,g)] + 1;
560 
561  remove_edge(*ei,g);
562  }
563 
564  }
565 
566 
567  // Set diagonal to weighted average of Laplacian preconditioner
568  // and the original matrix
569 
570  // First compute the "diagonal surplus" (in the original input matrix)
571  // If input is a (pure, Dirichlet) graph Laplacian , this will be 0
572  Epetra_Vector ones(Matrix_->OperatorDomainMap());
573  Epetra_Vector surplus(Matrix_->OperatorRangeMap());
574 
575  ones.PutScalar(1.0);
576  Matrix_->Multiply(false, ones, surplus);
577 
578  for(int i = 0; i < num_verts; i++)
579  {
580  Values[i][0] += surplus[i];
581  Values[i][0] = KeepDiag_*diagonal[i] +
582  (1.-KeepDiag_) * Values[i][0];
583  }
584 
585  // Create the CrsMatrix for the support graph
586  Support_ = rcp(new Epetra_CrsMatrix(Copy, Matrix().RowMatrixRowMap(),l, false));
587 
588 #ifndef EPETRA_NO_64BIT_GLOBAL_INDICES
589  if((*Matrix_).RowMatrixRowMap().GlobalIndicesLongLong())
590  {
591  // Fill in the matrix with the stl vectors for each row
592  for(int i = 0; i < num_verts; i++)
593  {
594  std::vector<long long> IndicesLL(l[i]);
595  for(int k = 0; k < l[i]; ++k)
596  IndicesLL[k] = Indices[i][k];
597 
598  (*Support_).InsertGlobalValues(i,l[i],&Values[i][0],&IndicesLL[0]);
599  }
600  }
601  else
602 #endif
603 #ifndef EPETRA_NO_32BIT_GLOBAL_INDICES
604  if((*Matrix_).RowMatrixRowMap().GlobalIndicesInt())
605  {
606  // Fill in the matrix with the stl vectors for each row
607  for(int i = 0; i < num_verts; i++)
608  {
609  (*Support_).InsertGlobalValues(i,l[i],&Values[i][0],&Indices[i][0]);
610  }
611  }
612  else
613 #endif
614  throw "Ifpack_SupportGraph::FindSupport: GlobalIndices unknown.";;
615 
616  (*Support_).FillComplete();
617 
618  delete edge_array;
619  delete weights;
620  delete values;
621  delete indices;
622  delete l;
623  delete diagonal;
624 
625  return(0);
626 }
627 //==============================================================================
628 template<typename T>
629 int Ifpack_SupportGraph<T>::SetParameters(Teuchos::ParameterList& List_in)
630 {
631  List_ = List_in;
632  NumForests_ = List_in.get("MST: forest number", NumForests_);
633  KeepDiag_ = List_in.get("MST: keep diagonal", KeepDiag_);
634  Randomize_ = List_in.get("MST: randomize", Randomize_);
635  // Diagonal pertubation parameters have weird names to be compatible with rest of Ifpack!
636  DiagPertRel_ = List_in.get("fact: relative threshold", DiagPertRel_);
637  DiagPertAbs_ = List_in.get("fact: absolute threshold", DiagPertAbs_);
638 
639  return(0);
640 }
641 //==============================================================================
642 template<typename T>
644 {
645  IsInitialized_ = false;
646  IsComputed_ = false;
647 
648 
649  if (Time_ == Teuchos::null)
650  {
651  Time_ = Teuchos::rcp( new Epetra_Time(Comm()) );
652  }
653 
654 
655  Time_->ResetStartTime();
656 
657  FindSupport();
658 
659  Inverse_ = Teuchos::rcp(new T(Support_.get()));
660 
661  IFPACK_CHK_ERR(Inverse_->Initialize());
662 
663  IsInitialized_ = true;
664  ++NumInitialize_;
665  InitializeTime_ += Time_->ElapsedTime();
666 
667  return(0);
668 
669 }
670 //==============================================================================
671 template<typename T>
673 {
674  if (IsInitialized() == false)
675  IFPACK_CHK_ERR(Initialize());
676 
677  Time_->ResetStartTime();
678  IsComputed_ = false;
679  Condest_ = -1.0;
680 
681  IFPACK_CHK_ERR(Inverse_->Compute());
682 
683  IsComputed_ = true;
684  ++NumCompute_;
685  ComputeTime_ += Time_->ElapsedTime();
686 
687 
688  return(0);
689 }
690 //==============================================================================
691 template<typename T>
693 {
694  // store the flag -- it will be set in Initialize() if Inverse_ does not
695  // exist.
696  UseTranspose_ = UseTranspose_in;
697 
698  // If Inverse_ exists, pass it right now.
699  if (Inverse_!=Teuchos::null)
700  IFPACK_CHK_ERR(Inverse_->SetUseTranspose(UseTranspose_in));
701 
702  return(0);
703 }
704 //==============================================================================
705 template<typename T>
708 {
709  IFPACK_CHK_ERR(Matrix_->Apply(X,Y));
710  return(0);
711 }
712 //==============================================================================
713 template<typename T>
714 const char * Ifpack_SupportGraph<T>::Label() const
715 {
716  return(Label_.c_str());
717 }
718 //==============================================================================
719 template<typename T>
722 {
723  if (!IsComputed())
724  IFPACK_CHK_ERR(-3);
725 
726  Time_->ResetStartTime();
727 
728 
729  Inverse_->ApplyInverse(X,Y);
730 
731  ++NumApplyInverse_;
732  ApplyInverseTime_ += Time_->ElapsedTime();
733 
734  return(0);
735 }
736 //==============================================================================
737 template<typename T>
738 std::ostream& Ifpack_SupportGraph<T>::
739 Print(std::ostream& os) const
740 {
741  os << "================================================================================" << std::endl;
742  os << "Ifpack_SupportGraph: " << Label () << std::endl << std::endl;
743  os << "Condition number estimate = " << Condest() << std::endl;
744  os << "Global number of rows = " << Matrix_->NumGlobalRows64() << std::endl;
745  os << "Number of edges in support graph = " << (Support_->NumGlobalNonzeros64()-Support_->NumGlobalDiagonals64())/2 << std::endl;
746  os << "Fraction of off diagonals of support graph/off diagonals of original = "
747  << ((double)Support_->NumGlobalNonzeros64()-Support_->NumGlobalDiagonals64())/(Matrix_->NumGlobalNonzeros64()-Matrix_->NumGlobalDiagonals64());
748  os << std::endl;
749  os << "Phase # calls Total Time (s) Total MFlops MFlops/s" << std::endl;
750  os << "----- ------- -------------- ------------ --------" << std::endl;
751  os << "Initialize() " << std::setw(10) << NumInitialize_
752  << " " << std::setw(15) << InitializeTime_
753  << " 0.0 0.0" << std::endl;
754  os << "Compute() " << std::setw(10) << NumCompute_
755  << " " << std::setw(22) << ComputeTime_
756  << " " << std::setw(15) << 1.0e-6 * ComputeFlops_;
757  if (ComputeTime_ != 0.0)
758  os << " " << std::setw(15) << 1.0e-6 * ComputeFlops_ / ComputeTime_ << std::endl;
759  else
760  os << " " << std::setw(15) << 0.0 << std::endl;
761  os << "ApplyInverse() " << std::setw(10) << NumApplyInverse_
762  << " " << std::setw(22) << ApplyInverseTime_
763  << " " << std::setw(15) << 1.0e-6 * ApplyInverseFlops_;
764  if (ApplyInverseTime_ != 0.0)
765  os << " " << std::setw(15) << 1.0e-6 * ApplyInverseFlops_ / ApplyInverseTime_ << std::endl;
766  else
767  os << " " << std::setw(15) << 0.0 << std::endl;
768 
769  os << std::endl << std::endl;
770  os << "Now calling the underlying preconditioner's print()" << std::endl;
771 
772  Inverse_->Print(os);
773 
774  return os;
775 }
776 //==============================================================================
777 template<typename T>
779 Condest(const Ifpack_CondestType CT, const int MaxIters,
780  const double Tol, Epetra_RowMatrix* Matrix_in)
781 {
782  if (!IsComputed()) // cannot compute right now
783  {
784  return(-1.0);
785  }
786 
787  Condest_ = Ifpack_Condest(*this, CT, MaxIters, Tol, Matrix_in);
788 
789  return(Condest_);
790 }
791 
792 #endif // IFPACK_SUPPORTGRAPH_H
Ifpack_SupportGraph::NumApplyInverse_
int NumApplyInverse_
Contains the number of successful call to ApplyInverse().
Definition: Ifpack_SupportGraph.h:329
Ifpack_SupportGraph::IsComputed_
bool IsComputed_
If true, the preconditioner has been successfully computed.
Definition: Ifpack_SupportGraph.h:311
Ifpack_SupportGraph::KeepDiag_
double KeepDiag_
Contains the option to keep the diagonal of original matrix, or weighted average.
Definition: Ifpack_SupportGraph.h:365
Ifpack_SupportGraph::ApplyInverseFlops
virtual double ApplyInverseFlops() const
Returns the total number of flops to apply the preconditioner.
Definition: Ifpack_SupportGraph.h:285
Epetra_MultiVector::PutScalar
int PutScalar(double ScalarConstant)
Ifpack_SupportGraph::Comm
virtual const Epetra_Comm & Comm() const
Returns a pointer to the Epetra_Comm communicator associated with this operator.
Definition: Ifpack_SupportGraph.h:155
Ifpack_SupportGraph::Time_
Teuchos::RefCountPtr< Epetra_Time > Time_
Object used for timing purposes.
Definition: Ifpack_SupportGraph.h:350
Ifpack_SupportGraph::Print
virtual std::ostream & Print(std::ostream &) const
Prints on ostream basic information about this object.
Definition: Ifpack_SupportGraph.h:739
Ifpack_SupportGraph::SetUseTranspose
virtual int SetUseTranspose(bool UseTranspose_in)
If set true, transpose of this operator will be applied (not implemented).
Definition: Ifpack_SupportGraph.h:692
Ifpack_SupportGraph::Condest
virtual double Condest() const
Returns the computed condition number.
Definition: Ifpack_SupportGraph.h:222
Ifpack_SupportGraph::ComputeFlops
virtual double ComputeFlops() const
Returns the total number of flops to compute the preconditioner.
Definition: Ifpack_SupportGraph.h:279
Ifpack_SupportGraph::InitializeTime_
double InitializeTime_
Contains the time for all successful calls to Initialize().
Definition: Ifpack_SupportGraph.h:332
Ifpack_SupportGraph::NumInitialize_
int NumInitialize_
Contains the number of successful calls to Initialize().
Definition: Ifpack_SupportGraph.h:323
Copy
Copy
Ifpack_SupportGraph::Matrix_
Teuchos::RefCountPtr< const Epetra_RowMatrix > Matrix_
Pointers to the matrix to be preconditioned.
Definition: Ifpack_SupportGraph.h:299
Ifpack_SupportGraph::UseTranspose_
bool UseTranspose_
If true, solve with the transpose (not supported by all solvers).
Definition: Ifpack_SupportGraph.h:314
Ifpack_SupportGraph::OperatorDomainMap
virtual const Epetra_Map & OperatorDomainMap() const
Returns the Epetra_Map object associated with the domain of this operator.
Definition: Ifpack_SupportGraph.h:158
Ifpack_SupportGraph::NumCompute_
int NumCompute_
Contains the number of successful call to Compute().
Definition: Ifpack_SupportGraph.h:326
Epetra_Comm
Ifpack_SupportGraph::DiagPertAbs_
double DiagPertAbs_
Absolute diagonal pertubation.
Definition: Ifpack_SupportGraph.h:362
Ifpack_SupportGraph::Support_
Teuchos::RefCountPtr< Epetra_CrsMatrix > Support_
Pointers to the matrix of the support graph.
Definition: Ifpack_SupportGraph.h:302
Ifpack_SupportGraph::Matrix
virtual const Epetra_RowMatrix & Matrix() const
Returns a const reference to the internally stored matrix.
Definition: Ifpack_SupportGraph.h:228
Ifpack_SupportGraph::OperatorRangeMap
virtual const Epetra_Map & OperatorRangeMap() const
Returns the Epetra_Map object associated with the range of this operator.
Definition: Ifpack_SupportGraph.h:161
Ifpack_SupportGraph::ApplyInverse
virtual int ApplyInverse(const Epetra_MultiVector &X, Epetra_MultiVector &Y) const
Applies the preconditioner to X, returns the result in Y.
Definition: Ifpack_SupportGraph.h:721
Epetra_CrsMatrix
Ifpack_SupportGraph::Compute
virtual int Compute()
Computes the preconditioners.
Definition: Ifpack_SupportGraph.h:672
Ifpack_SupportGraph::ApplyInverseTime_
double ApplyInverseTime_
Contains the time for all successful calls to ApplyInverse().
Definition: Ifpack_SupportGraph.h:338
Ifpack_SupportGraph::UseTranspose
virtual bool UseTranspose() const
Returns the current UseTranspose setting.
Definition: Ifpack_SupportGraph.h:149
Ifpack_SupportGraph::Initialize
virtual int Initialize()
Initialize the preconditioner.
Definition: Ifpack_SupportGraph.h:643
Ifpack_SupportGraph::IsInitialized
virtual bool IsInitialized() const
Returns true if the preconditioner has been successfully initialized.
Definition: Ifpack_SupportGraph.h:169
Ifpack_SupportGraph::SetParameters
virtual int SetParameters(Teuchos::ParameterList &List)
Sets all the parameters for the preconditioner.
Definition: Ifpack_SupportGraph.h:629
Epetra_RowMatrix
Ifpack_SupportGraph::ApplyInverseTime
virtual double ApplyInverseTime() const
Returns the total time spent in ApplyInverse().
Definition: Ifpack_SupportGraph.h:267
Ifpack_SupportGraph::ApplyInverseFlops_
double ApplyInverseFlops_
Contain sthe number of flops for ApplyInverse().
Definition: Ifpack_SupportGraph.h:347
Ifpack_SupportGraph::NumApplyInverse
virtual int NumApplyInverse() const
Returns the number of calls to ApplyInverse().
Definition: Ifpack_SupportGraph.h:249
Ifpack_SupportGraph::ComputeTime_
double ComputeTime_
Contains the time for all successful calls to Compute().
Definition: Ifpack_SupportGraph.h:335
Ifpack_SupportGraph::Inverse_
Teuchos::RefCountPtr< T > Inverse_
Pointer to the local solver.
Definition: Ifpack_SupportGraph.h:353
Ifpack_SupportGraph::NumCompute
virtual int NumCompute() const
Returns the number of calls to Compute().
Definition: Ifpack_SupportGraph.h:243
Ifpack_SupportGraph::FindSupport
int FindSupport()
Compute the support graph.
Definition: Ifpack_SupportGraph.h:403
Epetra_Vector
Ifpack_SupportGraph::Randomize_
int Randomize_
Option to add random pertubation to edge weights, to get random spanning trees.
Definition: Ifpack_SupportGraph.h:368
Ifpack_SupportGraph::Label
virtual const char * Label() const
Returns a character string describing the operator.
Definition: Ifpack_SupportGraph.h:714
Ifpack_SupportGraph::InitializeFlops
virtual double InitializeFlops() const
Returns the number of flops in the initialization phase.
Definition: Ifpack_SupportGraph.h:273
Ifpack_SupportGraph::NumInitialize
virtual int NumInitialize() const
Returns the number of calls to Initialize().
Definition: Ifpack_SupportGraph.h:237
Ifpack_SupportGraph::Label_
std::string Label_
Contains the label of this object.
Definition: Ifpack_SupportGraph.h:305
Epetra_Time
Ifpack_SupportGraph::Ifpack_SupportGraph
Ifpack_SupportGraph(Epetra_RowMatrix *Matrix_in)
Constructor.
Definition: Ifpack_SupportGraph.h:376
Ifpack_SupportGraph::IsComputed
virtual bool IsComputed() const
Returns true if the preconditioner has been successfully computed.
Definition: Ifpack_SupportGraph.h:175
Epetra_MultiVector
Ifpack_SupportGraph::Condest_
double Condest_
Contains the estimated condition number.
Definition: Ifpack_SupportGraph.h:320
Ifpack_SupportGraph::DiagPertRel_
double DiagPertRel_
Relative diagonal pertubation.
Definition: Ifpack_SupportGraph.h:359
Ifpack_Preconditioner
Ifpack_Preconditioner: basic class for preconditioning in Ifpack.
Definition: Ifpack_Preconditioner.h:136
Ifpack_SupportGraph::List_
Teuchos::ParameterList List_
Stores a copy of the list given in SetParameters()
Definition: Ifpack_SupportGraph.h:317
Ifpack_SupportGraph::Apply
virtual int Apply(const Epetra_MultiVector &X, Epetra_MultiVector &Y) const
Applies the matrix to an Epetra_MultiVector.
Definition: Ifpack_SupportGraph.h:707
Ifpack_SupportGraph::HasNormInf
virtual bool HasNormInf() const
Returns true if this object can provide an approximate Inf-norm, false otherwise.
Definition: Ifpack_SupportGraph.h:152
Ifpack_SupportGraph::IsInitialized_
bool IsInitialized_
If true, the preconditioner has been successfully initialized.
Definition: Ifpack_SupportGraph.h:308
Ifpack_SupportGraph::InitializeFlops_
double InitializeFlops_
Contains the number of flops for Initialize().
Definition: Ifpack_SupportGraph.h:341
Ifpack_SupportGraph::NumForests_
int NumForests_
Contains the number of forests in the support graph.
Definition: Ifpack_SupportGraph.h:356
Epetra_Map
Ifpack_SupportGraph::ComputeFlops_
double ComputeFlops_
Contains the number of flops for Compute().
Definition: Ifpack_SupportGraph.h:344
Ifpack_SupportGraph
Definition: Ifpack_SupportGraph.h:80
Ifpack_SupportGraph::ComputeTime
virtual double ComputeTime() const
Returns the total time spent in Compute().
Definition: Ifpack_SupportGraph.h:261
Ifpack_SupportGraph::InitializeTime
virtual double InitializeTime() const
Returns the total time spent in Initialize().
Definition: Ifpack_SupportGraph.h:255
Ifpack_SupportGraph::NormInf
virtual double NormInf() const
Returns the infinity norm of the global matrix (not implemented)
Definition: Ifpack_SupportGraph.h:138