All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Xpetra_TpetraCrsMatrix_decl.hpp
Go to the documentation of this file.
1 // @HEADER
2 //
3 // ***********************************************************************
4 //
5 // Xpetra: A linear algebra interface package
6 // Copyright 2012 Sandia Corporation
7 //
8 // Under the terms of Contract DE-AC04-94AL85000 with Sandia Corporation,
9 // the U.S. Government retains certain rights in this software.
10 //
11 // Redistribution and use in source and binary forms, with or without
12 // modification, are permitted provided that the following conditions are
13 // met:
14 //
15 // 1. Redistributions of source code must retain the above copyright
16 // notice, this list of conditions and the following disclaimer.
17 //
18 // 2. Redistributions in binary form must reproduce the above copyright
19 // notice, this list of conditions and the following disclaimer in the
20 // documentation and/or other materials provided with the distribution.
21 //
22 // 3. Neither the name of the Corporation nor the names of the
23 // contributors may be used to endorse or promote products derived from
24 // this software without specific prior written permission.
25 //
26 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
27 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
28 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
29 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
30 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
31 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
32 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
33 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
34 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
35 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
36 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
37 //
38 // Questions? Contact
39 // Jonathan Hu (jhu@sandia.gov)
40 // Andrey Prokopenko (aprokop@sandia.gov)
41 // Ray Tuminaro (rstumin@sandia.gov)
42 //
43 // ***********************************************************************
44 //
45 // @HEADER
46 #ifndef XPETRA_TPETRACRSMATRIX_DECL_HPP
47 #define XPETRA_TPETRACRSMATRIX_DECL_HPP
48 
49 /* this file is automatically generated - do not edit (see scripts/tpetra.py) */
50 
51 // FIXME (mfh 03 Sep 2014) The above is probably not true anymore.
52 // Furthermore, I don't think anyone maintains the scripts.
53 // Feel free to correct this comment if I'm wrong.
54 
56 
57 #include "Tpetra_CrsMatrix.hpp"
58 #include "Tpetra_replaceDiagonalCrsMatrix.hpp"
59 
60 #include "Xpetra_CrsMatrix.hpp"
61 #include "Xpetra_TpetraMap.hpp"
62 #include "Xpetra_TpetraMultiVector.hpp"
63 #include "Xpetra_TpetraVector.hpp"
65 #include "Xpetra_Exceptions.hpp"
66 
67 namespace Xpetra {
68 
69  // TODO: move that elsewhere
70  // template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
71  // const Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> toTpetraCrsMatrix(const Xpetra::DistObject<char, LocalOrdinal, GlobalOrdinal, Node> &);
72  //
73 
74 
75  template<class Scalar = CrsMatrix<>::scalar_type,
76  class LocalOrdinal =
77  typename CrsMatrix<Scalar>::local_ordinal_type,
78  class GlobalOrdinal =
79  typename CrsMatrix<Scalar, LocalOrdinal>::global_ordinal_type,
80  class Node =
81  typename CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal>::node_type>
83  : public CrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node> //, public TpetraRowMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>
84  {
85 
86  // The following typedef are used by the XPETRA_DYNAMIC_CAST() macro.
91 
92  // The following typedefs are used by the Kokkos interface
93 #ifdef HAVE_XPETRA_KOKKOS_REFACTOR
94 #ifdef HAVE_XPETRA_TPETRA
96 #endif
97 #endif
98 
99  public:
100 
102 
103 
105  TpetraCrsMatrix(const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &rowMap, size_t maxNumEntriesPerRow, ProfileType pftype=DynamicProfile, const Teuchos::RCP< Teuchos::ParameterList > &params=Teuchos::null);
106 
108  TpetraCrsMatrix(const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &rowMap, const ArrayRCP< const size_t > &NumEntriesPerRowToAlloc, ProfileType pftype=DynamicProfile, const Teuchos::RCP< Teuchos::ParameterList > &params=Teuchos::null);
109 
111  TpetraCrsMatrix(const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &rowMap, const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &colMap, size_t maxNumEntriesPerRow, ProfileType pftype=DynamicProfile, const Teuchos::RCP< Teuchos::ParameterList > &params=Teuchos::null);
112 
115 
118 
122  const Teuchos::RCP<const Map<LocalOrdinal,GlobalOrdinal,Node> >& domainMap = Teuchos::null,
123  const Teuchos::RCP<const Map<LocalOrdinal,GlobalOrdinal,Node> >& rangeMap = Teuchos::null,
124  const Teuchos::RCP<Teuchos::ParameterList>& params = Teuchos::null);
125 
129  const Teuchos::RCP<const Map<LocalOrdinal,GlobalOrdinal,Node> >& domainMap = Teuchos::null,
130  const Teuchos::RCP<const Map<LocalOrdinal,GlobalOrdinal,Node> >& rangeMap = Teuchos::null,
131  const Teuchos::RCP<Teuchos::ParameterList>& params = Teuchos::null);
132 
135  const Import<LocalOrdinal,GlobalOrdinal,Node> & RowImporter,
136  const Teuchos::RCP<const Import<LocalOrdinal,GlobalOrdinal,Node> > DomainImporter,
137  const Teuchos::RCP<const Map<LocalOrdinal,GlobalOrdinal,Node> >& domainMap,
138  const Teuchos::RCP<const Map<LocalOrdinal,GlobalOrdinal,Node> >& rangeMap,
140 
143  const Export<LocalOrdinal,GlobalOrdinal,Node> & RowExporter,
144  const Teuchos::RCP<const Export<LocalOrdinal,GlobalOrdinal,Node> > DomainExporter,
145  const Teuchos::RCP<const Map<LocalOrdinal,GlobalOrdinal,Node> >& domainMap,
146  const Teuchos::RCP<const Map<LocalOrdinal,GlobalOrdinal,Node> >& rangeMap,
148 
149 #ifdef HAVE_XPETRA_KOKKOS_REFACTOR
150 #ifdef HAVE_XPETRA_TPETRA
173  const local_matrix_type& lclMatrix,
174  const Teuchos::RCP<Teuchos::ParameterList>& params = null);
175 
178  const local_matrix_type& lclMatrix,
181  const Teuchos::RCP<const Map<LocalOrdinal,GlobalOrdinal,Node> >& domainMap = Teuchos::null,
182  const Teuchos::RCP<const Map<LocalOrdinal,GlobalOrdinal,Node> >& rangeMap = Teuchos::null,
183  const Teuchos::RCP<Teuchos::ParameterList>& params = null);
184 #endif
185 #endif
186 
188  virtual ~TpetraCrsMatrix();
189 
191 
193 
194 
196  void insertGlobalValues(GlobalOrdinal globalRow, const ArrayView< const GlobalOrdinal > &cols, const ArrayView< const Scalar > &vals);
197 
199  void insertLocalValues(LocalOrdinal localRow, const ArrayView< const LocalOrdinal > &cols, const ArrayView< const Scalar > &vals);
200 
202  void replaceGlobalValues(GlobalOrdinal globalRow, const ArrayView< const GlobalOrdinal > &cols, const ArrayView< const Scalar > &vals);
203 
205  void
206  replaceLocalValues (LocalOrdinal localRow,
207  const ArrayView<const LocalOrdinal> &cols,
208  const ArrayView<const Scalar> &vals);
209 
211  void setAllToScalar(const Scalar &alpha);
212 
214  void scale(const Scalar &alpha);
215 
217  //** \warning This is an expert-only routine and should not be called from user code. */
218  void allocateAllValues(size_t numNonZeros,ArrayRCP<size_t> & rowptr, ArrayRCP<LocalOrdinal> & colind, ArrayRCP<Scalar> & values)
219  ;
220 
222  void setAllValues(const ArrayRCP<size_t> & rowptr, const ArrayRCP<LocalOrdinal> & colind, const ArrayRCP<Scalar> & values)
223  ;
224 
227  ;
228 
229  bool haveGlobalConstants() const
230  ;
231 
233 
235 
236 
238  void resumeFill(const RCP< ParameterList > &params=null);
239 
241  void fillComplete(const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &domainMap, const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &rangeMap, const RCP< ParameterList > &params=null);
242 
244  void fillComplete(const RCP< ParameterList > &params=null);
245 
246 
249 
252  const RCP<const Map<LocalOrdinal,GlobalOrdinal,Node> > & rangeMap,
253  const RCP<const Import<LocalOrdinal,GlobalOrdinal,Node> > &importer=Teuchos::null,
254  const RCP<const Export<LocalOrdinal,GlobalOrdinal,Node> > &exporter=Teuchos::null,
255  const RCP<ParameterList> &params=Teuchos::null);
256 
258 
260 
261 
264 
267 
270 
273 
276 
278  size_t getNodeNumRows() const;
279 
281  size_t getNodeNumCols() const;
282 
285 
287  size_t getNodeNumEntries() const;
288 
290  size_t getNumEntriesInLocalRow(LocalOrdinal localRow) const;
291 
293  size_t getGlobalMaxNumRowEntries() const;
294 
296  size_t getNodeMaxNumRowEntries() const;
297 
299  bool isLocallyIndexed() const;
300 
302  bool isGloballyIndexed() const;
303 
305  bool isFillComplete() const;
306 
308  bool isFillActive() const;
309 
312 
314  bool supportsRowViews() const;
315 
317  void getLocalRowCopy(LocalOrdinal LocalRow, const ArrayView< LocalOrdinal > &Indices, const ArrayView< Scalar > &Values, size_t &NumEntries) const;
318 
320  void getGlobalRowView(GlobalOrdinal GlobalRow, ArrayView< const GlobalOrdinal > &indices, ArrayView< const Scalar > &values) const;
321 
323  void getGlobalRowCopy(GlobalOrdinal GlobalRow, const ArrayView< GlobalOrdinal > &indices, const ArrayView< Scalar > &values, size_t &numEntries) const;
324 
326  void getLocalRowView(LocalOrdinal LocalRow, ArrayView< const LocalOrdinal > &indices, ArrayView< const Scalar > &values) const;
327 
329 
331 
332 
335 
338 
341 
343 
345 
346 
348  std::string description() const;
349 
352 
354 
356 
357  void setObjectLabel( const std::string &objectLabel );
359 
360 
361 
363  TpetraCrsMatrix(const TpetraCrsMatrix& matrix);
364 
367 
369  void getLocalDiagOffsets(Teuchos::ArrayRCP<size_t> &offsets) const;
370 
373 
376 
379 
382 
384  //{@
385 
388 
392 
396 
400 
404 
406 
407  // @}
408 #ifdef XPETRA_ENABLE_DEPRECATED_CODE
409  template<class Node2>
411  clone (const RCP<Node2> &node2) const
412  {
414  }
415 #endif
416 
418 
419 
421  bool hasMatrix() const;
422 
424  TpetraCrsMatrix(const Teuchos::RCP<Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> > &mtx);
425 
428 
431 
432 #ifdef HAVE_XPETRA_KOKKOS_REFACTOR
433 #ifdef HAVE_XPETRA_TPETRA
434  local_matrix_type getLocalMatrix () const {
436  return getTpetra_CrsMatrixNonConst()->getLocalMatrix();
437  }
438 
439  void setAllValues (const typename local_matrix_type::row_map_type& ptr,
440  const typename local_matrix_type::StaticCrsGraphType::entries_type::non_const_type& ind,
441  const typename local_matrix_type::values_type& val) {
442  getTpetra_CrsMatrixNonConst()->setAllValues(ptr,ind,val);
443  }
444 #endif
445 #endif
446 
448 
449  private:
451  }; // TpetraCrsMatrix class
452 
453 #ifdef HAVE_XPETRA_EPETRA
454 
455 #if ((defined(EPETRA_HAVE_OMP) && (!defined(HAVE_TPETRA_INST_OPENMP) || !defined(HAVE_TPETRA_INST_INT_INT))) || \
456  (!defined(EPETRA_HAVE_OMP) && (!defined(HAVE_TPETRA_INST_SERIAL) || !defined(HAVE_TPETRA_INST_INT_INT))))
457 
458  // specialization of TpetraCrsMatrix for GO=LO=int
459  template<class Scalar>
460  class TpetraCrsMatrix<Scalar,int,int,EpetraNode>
461  : public CrsMatrix<Scalar,int,int,EpetraNode> //, public TpetraRowMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>
462  {
463  typedef int LocalOrdinal;
464  typedef int GlobalOrdinal;
465  typedef EpetraNode Node;
466  // The following typedef are used by the XPETRA_DYNAMIC_CAST() macro.
471 
472  // The following typedefs are used by the Kokkos interface
473 #ifdef HAVE_XPETRA_KOKKOS_REFACTOR
474 #ifdef HAVE_XPETRA_TPETRA
476 #endif
477 #endif
478 
479  public:
480 
482 
483 
485  TpetraCrsMatrix(const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &rowMap, size_t maxNumEntriesPerRow, ProfileType pftype=DynamicProfile, const Teuchos::RCP< Teuchos::ParameterList > &params=Teuchos::null) {
487  }
488 
492  }
493 
495  TpetraCrsMatrix(const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &rowMap, const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &colMap, size_t maxNumEntriesPerRow, ProfileType pftype=DynamicProfile, const Teuchos::RCP< Teuchos::ParameterList > &params=Teuchos::null) {
497  }
498 
502  }
503 
507  }
508 
509 
510 
514  const Teuchos::RCP<const Map<LocalOrdinal,GlobalOrdinal,Node> >& domainMap = Teuchos::null,
515  const Teuchos::RCP<const Map<LocalOrdinal,GlobalOrdinal,Node> >& rangeMap = Teuchos::null,
516  const Teuchos::RCP<Teuchos::ParameterList>& params = Teuchos::null) {
518  }
519 
523  const Teuchos::RCP<const Map<LocalOrdinal,GlobalOrdinal,Node> >& domainMap = Teuchos::null,
524  const Teuchos::RCP<const Map<LocalOrdinal,GlobalOrdinal,Node> >& rangeMap = Teuchos::null,
525  const Teuchos::RCP<Teuchos::ParameterList>& params = Teuchos::null) {
527  }
528 
531  const Import<LocalOrdinal,GlobalOrdinal,Node> & RowImporter,
532  const Teuchos::RCP<const Import<LocalOrdinal,GlobalOrdinal,Node> > DomainImporter,
533  const Teuchos::RCP<const Map<LocalOrdinal,GlobalOrdinal,Node> >& domainMap,
534  const Teuchos::RCP<const Map<LocalOrdinal,GlobalOrdinal,Node> >& rangeMap,
536  {
538  }
539 
542  const Export<LocalOrdinal,GlobalOrdinal,Node> & RowExporter,
543  const Teuchos::RCP<const Export<LocalOrdinal,GlobalOrdinal,Node> > DomainExporter,
544  const Teuchos::RCP<const Map<LocalOrdinal,GlobalOrdinal,Node> >& domainMap,
545  const Teuchos::RCP<const Map<LocalOrdinal,GlobalOrdinal,Node> >& rangeMap,
547  {
549  }
550 
551 #ifdef HAVE_XPETRA_KOKKOS_REFACTOR
552 #ifdef HAVE_XPETRA_TPETRA
575  const local_matrix_type& lclMatrix,
576  const Teuchos::RCP<Teuchos::ParameterList>& params = null) {
577  XPETRA_TPETRA_ETI_EXCEPTION( typeid(TpetraCrsMatrix<Scalar,int,int,EpetraNode>).name() , "TpetraCrsMatrix<int,int>", "int", typeid(EpetraNode).name() );
578  }
579 
582  const local_matrix_type& lclMatrix,
583  const Teuchos::RCP<const Map<LocalOrdinal,GlobalOrdinal,Node> >& rowMap,
584  const Teuchos::RCP<const Map<LocalOrdinal,GlobalOrdinal,Node> >& colMap,
585  const Teuchos::RCP<const Map<LocalOrdinal,GlobalOrdinal,Node> >& domainMap = Teuchos::null,
586  const Teuchos::RCP<const Map<LocalOrdinal,GlobalOrdinal,Node> >& rangeMap = Teuchos::null,
587  const Teuchos::RCP<Teuchos::ParameterList>& params = null) {
588  XPETRA_TPETRA_ETI_EXCEPTION( typeid(TpetraCrsMatrix<Scalar,int,int,EpetraNode>).name() , "TpetraCrsMatrix<int,int>", "int", typeid(EpetraNode).name() );
589  }
590 #endif
591 #endif
592 
594  virtual ~TpetraCrsMatrix() { }
595 
597 
599 
600 
603 
606 
609 
611  void
613  const ArrayView<const LocalOrdinal> &cols,
614  const ArrayView<const Scalar> &vals)
615  { }
616 
618  void setAllToScalar(const Scalar &alpha) { }
619 
621  void scale(const Scalar &alpha) { }
622 
624  //** \warning This is an expert-only routine and should not be called from user code. */
625  void allocateAllValues(size_t numNonZeros,ArrayRCP<size_t> & rowptr, ArrayRCP<LocalOrdinal> & colind, ArrayRCP<Scalar> & values) { }
626 
628  void setAllValues(const ArrayRCP<size_t> & rowptr, const ArrayRCP<LocalOrdinal> & colind, const ArrayRCP<Scalar> & values) { }
629 
632 
633  bool haveGlobalConstants() const { return false;}
634 
636 
638 
639 
641  void resumeFill(const RCP< ParameterList > &params=null) { }
642 
644  void fillComplete(const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &domainMap, const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &rangeMap, const RCP< ParameterList > &params=null) { }
645 
647  void fillComplete(const RCP< ParameterList > &params=null) { }
648 
649 
652 
655  const RCP<const Map<LocalOrdinal,GlobalOrdinal,Node> > & rangeMap,
656  const RCP<const Import<LocalOrdinal,GlobalOrdinal,Node> > &importer=Teuchos::null,
657  const RCP<const Export<LocalOrdinal,GlobalOrdinal,Node> > &exporter=Teuchos::null,
658  const RCP<ParameterList> &params=Teuchos::null) { }
659 
661 
663 
664 
666  const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getRowMap() const { return Teuchos::null; }
667 
669  const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getColMap() const { return Teuchos::null; }
670 
673 
675  global_size_t getGlobalNumRows() const { return 0; }
676 
678  global_size_t getGlobalNumCols() const { return 0; }
679 
681  size_t getNodeNumRows() const { return 0; }
682 
684  size_t getNodeNumCols() const { return 0; }
685 
687  global_size_t getGlobalNumEntries() const { return 0; }
688 
690  size_t getNodeNumEntries() const { return 0; }
691 
693  size_t getNumEntriesInLocalRow(LocalOrdinal localRow) const { return 0; }
694 
696  size_t getGlobalMaxNumRowEntries() const { return 0; }
697 
699  size_t getNodeMaxNumRowEntries() const { return 0; }
700 
702  bool isLocallyIndexed() const { return false; }
703 
705  bool isGloballyIndexed() const { return false; }
706 
708  bool isFillComplete() const { return false; }
709 
711  bool isFillActive() const { return false; }
712 
715 
717  bool supportsRowViews() const { return false; }
718 
720  void getLocalRowCopy(LocalOrdinal LocalRow, const ArrayView< LocalOrdinal > &Indices, const ArrayView< Scalar > &Values, size_t &NumEntries) const { }
721 
724 
726  void getGlobalRowCopy(GlobalOrdinal GlobalRow, const ArrayView< GlobalOrdinal > &indices, const ArrayView< Scalar > &values, size_t &numEntries) const { }
727 
730 
732 
734 
735 
738 
740  const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getDomainMap() const { return Teuchos::null; }
741 
743  const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getRangeMap() const { return Teuchos::null; }
744 
746 
748 
749 
751  std::string description() const { return std::string(""); }
752 
755 
757 
759 
760  void setObjectLabel( const std::string &objectLabel ) { }
762 
763 
766 
769 
772 
775 
778 
781 
784 
786  //{@
787 
790 
794 
798 
802 
806 
808 
809  // @}
810 #ifdef XPETRA_ENABLE_DEPRECATED_CODE
811  template<class Node2>
813  clone (const RCP<Node2> &node2) const { return Teuchos::null; }
814 #endif
815 
817 
819  bool hasMatrix() const { return false; }
820 
822  TpetraCrsMatrix(const Teuchos::RCP<Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> > &mtx) {
824  }
825 
828 
831 
832 #ifdef HAVE_XPETRA_KOKKOS_REFACTOR
833 #ifdef HAVE_XPETRA_TPETRA
834  local_matrix_type getLocalMatrix () const {
836  TEUCHOS_UNREACHABLE_RETURN(local_matrix_type());
837  }
838 
839  void setAllValues (const typename local_matrix_type::row_map_type& ptr,
840  const typename local_matrix_type::StaticCrsGraphType::entries_type::non_const_type& ind,
841  const typename local_matrix_type::values_type& val) { }
842 #endif
843 #endif
844 
846  }; // TpetraCrsMatrix class (specialization for GO=int, NO=EpetraNode)
847 #endif
848 
849 #if ((defined(EPETRA_HAVE_OMP) && (!defined(HAVE_TPETRA_INST_OPENMP) || !defined(HAVE_TPETRA_INST_INT_LONG_LONG))) || \
850  (!defined(EPETRA_HAVE_OMP) && (!defined(HAVE_TPETRA_INST_SERIAL) || !defined(HAVE_TPETRA_INST_INT_LONG_LONG))))
851 
852  // specialization of TpetraCrsMatrix for GO=long long, NO=EpetraNode
853  template<class Scalar>
854  class TpetraCrsMatrix<Scalar,int,long long,EpetraNode>
855  : public CrsMatrix<Scalar,int,long long,EpetraNode> //, public TpetraRowMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>
856  {
857  typedef int LocalOrdinal;
858  typedef long long GlobalOrdinal;
859  typedef EpetraNode Node;
860  // The following typedef are used by the XPETRA_DYNAMIC_CAST() macro.
865 
866  // The following typedefs are used by the Kokkos interface
867 #ifdef HAVE_XPETRA_KOKKOS_REFACTOR
868 #ifdef HAVE_XPETRA_TPETRA
870 #endif
871 #endif
872 
873  public:
874 
876 
877 
879  TpetraCrsMatrix(const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &rowMap, size_t maxNumEntriesPerRow, ProfileType pftype=DynamicProfile, const Teuchos::RCP< Teuchos::ParameterList > &params=Teuchos::null) {
881  }
882 
886  }
887 
889  TpetraCrsMatrix(const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &rowMap, const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &colMap, size_t maxNumEntriesPerRow, ProfileType pftype=DynamicProfile, const Teuchos::RCP< Teuchos::ParameterList > &params=Teuchos::null) {
891  }
892 
896  }
897 
901  }
902 
903 
904 
908  const Teuchos::RCP<const Map<LocalOrdinal,GlobalOrdinal,Node> >& domainMap = Teuchos::null,
909  const Teuchos::RCP<const Map<LocalOrdinal,GlobalOrdinal,Node> >& rangeMap = Teuchos::null,
910  const Teuchos::RCP<Teuchos::ParameterList>& params = Teuchos::null) {
912  }
913 
917  const Teuchos::RCP<const Map<LocalOrdinal,GlobalOrdinal,Node> >& domainMap = Teuchos::null,
918  const Teuchos::RCP<const Map<LocalOrdinal,GlobalOrdinal,Node> >& rangeMap = Teuchos::null,
919  const Teuchos::RCP<Teuchos::ParameterList>& params = Teuchos::null) {
921  }
922 
925  const Import<LocalOrdinal,GlobalOrdinal,Node> & RowImporter,
926  const Teuchos::RCP<const Import<LocalOrdinal,GlobalOrdinal,Node> > DomainImporter,
927  const Teuchos::RCP<const Map<LocalOrdinal,GlobalOrdinal,Node> >& domainMap,
928  const Teuchos::RCP<const Map<LocalOrdinal,GlobalOrdinal,Node> >& rangeMap,
930  {
932  }
933 
936  const Export<LocalOrdinal,GlobalOrdinal,Node> & RowExporter,
937  const Teuchos::RCP<const Export<LocalOrdinal,GlobalOrdinal,Node> > DomainExporter,
938  const Teuchos::RCP<const Map<LocalOrdinal,GlobalOrdinal,Node> >& domainMap,
939  const Teuchos::RCP<const Map<LocalOrdinal,GlobalOrdinal,Node> >& rangeMap,
941  {
943  }
944 
945 #ifdef HAVE_XPETRA_KOKKOS_REFACTOR
946 #ifdef HAVE_XPETRA_TPETRA
969  const local_matrix_type& lclMatrix,
970  const Teuchos::RCP<Teuchos::ParameterList>& params = null) {
972  }
973 
976  const local_matrix_type& lclMatrix,
977  const Teuchos::RCP<const Map<LocalOrdinal,GlobalOrdinal,Node> >& rowMap,
978  const Teuchos::RCP<const Map<LocalOrdinal,GlobalOrdinal,Node> >& colMap,
979  const Teuchos::RCP<const Map<LocalOrdinal,GlobalOrdinal,Node> >& domainMap = Teuchos::null,
980  const Teuchos::RCP<const Map<LocalOrdinal,GlobalOrdinal,Node> >& rangeMap = Teuchos::null,
981  const Teuchos::RCP<Teuchos::ParameterList>& params = null) {
982  XPETRA_TPETRA_ETI_EXCEPTION( typeid(TpetraCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,EpetraNode>).name() , typeid(TpetraCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,EpetraNode>).name(), "long long", typeid(EpetraNode).name() );
983  }
984 #endif
985 #endif
986 
988  virtual ~TpetraCrsMatrix() { }
989 
991 
993 
994 
997 
1000 
1003 
1005  void
1007  const ArrayView<const LocalOrdinal> &cols,
1008  const ArrayView<const Scalar> &vals)
1009  { }
1010 
1012  void setAllToScalar(const Scalar &alpha) { }
1013 
1015  void scale(const Scalar &alpha) { }
1016 
1018  //** \warning This is an expert-only routine and should not be called from user code. */
1019  void allocateAllValues(size_t numNonZeros,ArrayRCP<size_t> & rowptr, ArrayRCP<LocalOrdinal> & colind, ArrayRCP<Scalar> & values) { }
1020 
1022  void setAllValues(const ArrayRCP<size_t> & rowptr, const ArrayRCP<LocalOrdinal> & colind, const ArrayRCP<Scalar> & values) { }
1023 
1026 
1027  bool haveGlobalConstants() const { return false;}
1028 
1030 
1032 
1033 
1035  void resumeFill(const RCP< ParameterList > &params=null) { }
1036 
1038  void fillComplete(const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &domainMap, const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &rangeMap, const RCP< ParameterList > &params=null) { }
1039 
1041  void fillComplete(const RCP< ParameterList > &params=null) { }
1042 
1043 
1046 
1049  const RCP<const Map<LocalOrdinal,GlobalOrdinal,Node> > & rangeMap,
1050  const RCP<const Import<LocalOrdinal,GlobalOrdinal,Node> > &importer=Teuchos::null,
1051  const RCP<const Export<LocalOrdinal,GlobalOrdinal,Node> > &exporter=Teuchos::null,
1052  const RCP<ParameterList> &params=Teuchos::null) { }
1053 
1055 
1057 
1058 
1060  const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getRowMap() const { return Teuchos::null; }
1061 
1063  const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getColMap() const { return Teuchos::null; }
1064 
1067 
1069  global_size_t getGlobalNumRows() const { return 0; }
1070 
1072  global_size_t getGlobalNumCols() const { return 0; }
1073 
1075  size_t getNodeNumRows() const { return 0; }
1076 
1078  size_t getNodeNumCols() const { return 0; }
1079 
1081  global_size_t getGlobalNumEntries() const { return 0; }
1082 
1084  size_t getNodeNumEntries() const { return 0; }
1085 
1087  size_t getNumEntriesInLocalRow(LocalOrdinal localRow) const { return 0; }
1088 
1090  size_t getGlobalMaxNumRowEntries() const { return 0; }
1091 
1093  size_t getNodeMaxNumRowEntries() const { return 0; }
1094 
1096  bool isLocallyIndexed() const { return false; }
1097 
1099  bool isGloballyIndexed() const { return false; }
1100 
1102  bool isFillComplete() const { return false; }
1103 
1105  bool isFillActive() const { return false; }
1106 
1109 
1111  bool supportsRowViews() const { return false; }
1112 
1114  void getLocalRowCopy(LocalOrdinal LocalRow, const ArrayView< LocalOrdinal > &Indices, const ArrayView< Scalar > &Values, size_t &NumEntries) const { }
1115 
1118 
1120  void getGlobalRowCopy(GlobalOrdinal GlobalRow, const ArrayView< GlobalOrdinal > &indices, const ArrayView< Scalar > &values, size_t &numEntries) const { }
1121 
1124 
1126 
1128 
1129 
1132 
1134  const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getDomainMap() const { return Teuchos::null; }
1135 
1137  const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getRangeMap() const { return Teuchos::null; }
1138 
1140 
1142 
1143 
1145  std::string description() const { return std::string(""); }
1146 
1149 
1151 
1153 
1154  void setObjectLabel( const std::string &objectLabel ) { }
1156 
1159 
1162 
1165 
1168 
1171 
1174 
1177 
1179  //{@
1180 
1183 
1187 
1191 
1195 
1199 
1201 
1202  // @}
1203 #ifdef XPETRA_ENABLE_DEPRECATED_CODE
1204  template<class Node2>
1206  clone (const RCP<Node2> &node2) const { return Teuchos::null; }
1207 #endif
1208 
1210 
1212  bool hasMatrix() const { return false; }
1213 
1215  TpetraCrsMatrix(const Teuchos::RCP<Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> > &mtx) {
1217  }
1218 
1221 
1224 
1225 #ifdef HAVE_XPETRA_KOKKOS_REFACTOR
1226 #ifdef HAVE_XPETRA_TPETRA
1227  local_matrix_type getLocalMatrix () const {
1229  TEUCHOS_UNREACHABLE_RETURN(local_matrix_type());
1230  }
1231 
1232  void setAllValues (const typename local_matrix_type::row_map_type& ptr,
1233  const typename local_matrix_type::StaticCrsGraphType::entries_type::non_const_type& ind,
1234  const typename local_matrix_type::values_type& val) { }
1235 #endif
1236 #endif
1237 
1239  }; // TpetraCrsMatrix class (specialization for GO=long long, NO=EpetraNode)
1240 #endif
1241 
1242 #endif // HAVE_XPETRA_EPETRA
1243 
1244 } // Xpetra namespace
1245 
1246 #define XPETRA_TPETRACRSMATRIX_SHORT
1247 #endif // XPETRA_TPETRACRSMATRIX_DECL_HPP
void doImport(const DistObject< char, LocalOrdinal, GlobalOrdinal, Node > &source, const Import< LocalOrdinal, GlobalOrdinal, Node > &importer, CombineMode CM)
Import.
void getLocalRowCopy(LocalOrdinal LocalRow, const ArrayView< LocalOrdinal > &Indices, const ArrayView< Scalar > &Values, size_t &NumEntries) const
Extract a list of entries in a specified local row of the matrix. Put into storage allocated by calli...
size_t getGlobalMaxNumRowEntries() const
Returns the maximum number of entries across all rows/columns on all nodes.
global_size_t getGlobalNumRows() const
Number of global elements in the row map of this matrix.
RCP< Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > getTpetra_CrsMatrixNonConst() const
Get the underlying Tpetra matrix.
TpetraCrsMatrix(const Teuchos::RCP< const CrsGraph< LocalOrdinal, GlobalOrdinal, Node > > &graph, const Teuchos::RCP< Teuchos::ParameterList > &params=Teuchos::null)
Constructor specifying a previously constructed graph.
std::string description() const
A simple one-line description of this object.
TpetraCrsMatrix(const Teuchos::RCP< const CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > &sourceMatrix, const Export< LocalOrdinal, GlobalOrdinal, Node > &exporter, const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &domainMap=Teuchos::null, const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &rangeMap=Teuchos::null, const Teuchos::RCP< Teuchos::ParameterList > &params=Teuchos::null)
Constructor for a fused export.
Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getMap() const
Implements DistObject interface.
RCP< Map< LocalOrdinal, GlobalOrdinal, Node2 > > clone(const Map< LocalOrdinal, GlobalOrdinal, Node1 > &map, const RCP< Node2 > &node2)
size_t getNumEntriesInLocalRow(LocalOrdinal localRow) const
Returns the current number of entries on this node in the specified local row.
void fillComplete(const RCP< ParameterList > &params=null)
Signal that data entry is complete.
TpetraCrsMatrix(const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &rowMap, const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &colMap, size_t maxNumEntriesPerRow, ProfileType pftype=DynamicProfile, const Teuchos::RCP< Teuchos::ParameterList > &params=Teuchos::null)
Constructor specifying column Map and fixed number of entries for each row.
void doImport(const DistObject< char, LocalOrdinal, GlobalOrdinal, Node > &source, const Import< LocalOrdinal, GlobalOrdinal, Node > &importer, CombineMode CM)
Import.
void setAllToScalar(const Scalar &alpha)
Set all matrix entries equal to scalarThis.
TpetraExport< LocalOrdinal, GlobalOrdinal, Node > TpetraExportClass
TpetraVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > TpetraVectorClass
const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getDomainMap() const
Returns the Map associated with the domain of this operator. This will be null until fillComplete() i...
void removeEmptyProcessesInPlace(const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &newMap)
TpetraImport< LocalOrdinal, GlobalOrdinal, Node > TpetraImportClass
void replaceDomainMapAndImporter(const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &newDomainMap, Teuchos::RCP< const Import< LocalOrdinal, GlobalOrdinal, Node > > &newImporter)
Replaces the current domainMap and importer with the user-specified objects.
void getLocalDiagCopy(Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &diag) const
Get a copy of the diagonal entries owned by this node, with local row indices.
TpetraCrsMatrix(const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &rowMap, size_t maxNumEntriesPerRow, ProfileType pftype=DynamicProfile, const Teuchos::RCP< Teuchos::ParameterList > &params=Teuchos::null)
Constructor specifying fixed number of entries for each row.
void replaceLocalValues(LocalOrdinal localRow, const ArrayView< const LocalOrdinal > &cols, const ArrayView< const Scalar > &vals)
Replace matrix entries, using local IDs.
global_size_t getGlobalNumCols() const
Number of global columns in the matrix.
void doExport(const DistObject< char, LocalOrdinal, GlobalOrdinal, Node > &dest, const Import< LocalOrdinal, GlobalOrdinal, Node > &importer, CombineMode CM)
Export.
void insertLocalValues(LocalOrdinal localRow, const ArrayView< const LocalOrdinal > &cols, const ArrayView< const Scalar > &vals)
Insert matrix entries, using local IDs.
TpetraExport< LocalOrdinal, GlobalOrdinal, Node > TpetraExportClass
const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getDomainMap() const
Returns the Map associated with the domain of this operator. This will be null until fillComplete() i...
void replaceGlobalValues(GlobalOrdinal globalRow, const ArrayView< const GlobalOrdinal > &cols, const ArrayView< const Scalar > &vals)
Replace matrix entries, using global IDs.
TpetraCrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > TpetraCrsMatrixClass
TpetraCrsMatrix(const Teuchos::RCP< Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > &mtx)
TpetraCrsMatrix constructor to wrap a Tpetra::CrsMatrix object.
bool isLocallyIndexed() const
If matrix indices are in the local range, this function returns true. Otherwise, this function return...
global_size_t getGlobalNumEntries() const
Returns the global number of entries in this matrix.
size_t getNodeNumRows() const
Returns the number of matrix rows owned on the calling node.
void getLocalDiagOffsets(Teuchos::ArrayRCP< size_t > &offsets) const
Get offsets of the diagonal entries in the matrix.
void getLocalRowView(LocalOrdinal LocalRow, ArrayView< const LocalOrdinal > &indices, ArrayView< const Scalar > &values) const
Extract a const, non-persisting view of local indices in a specified row of the matrix.
TpetraCrsMatrix(const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &rowMap, const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &colMap, const ArrayRCP< const size_t > &NumEntriesPerRowToAlloc, ProfileType pftype=DynamicProfile, const Teuchos::RCP< Teuchos::ParameterList > &params=Teuchos::null)
Constructor specifying column Map and number of entries in each row.
void getLocalDiagCopy(Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &diag) const
Get a copy of the diagonal entries owned by this node, with local row idices.
RCP< Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > mtx_
void apply(const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &X, MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &Y, Teuchos::ETransp mode=Teuchos::NO_TRANS, Scalar alpha=ScalarTraits< Scalar >::one(), Scalar beta=ScalarTraits< Scalar >::zero()) const
Computes the sparse matrix-multivector multiplication.
void setAllValues(const ArrayRCP< size_t > &rowptr, const ArrayRCP< LocalOrdinal > &colind, const ArrayRCP< Scalar > &values)
Sets the 1D pointer arrays of the graph.
void rightScale(const Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &x)
Right scale operator with given vector values.
void allocateAllValues(size_t numNonZeros, ArrayRCP< size_t > &rowptr, ArrayRCP< LocalOrdinal > &colind, ArrayRCP< Scalar > &values)
Allocates and returns ArrayRCPs of the Crs arrays — This is an Xpetra-only routine.
void replaceLocalValues(LocalOrdinal localRow, const ArrayView< const LocalOrdinal > &cols, const ArrayView< const Scalar > &vals)
Replace matrix entries, using local IDs.
void resumeFill(const RCP< ParameterList > &params=null)
void replaceLocalValues(LocalOrdinal localRow, const ArrayView< const LocalOrdinal > &cols, const ArrayView< const Scalar > &vals)
Replace matrix entries, using local IDs.
bool isLocallyIndexed() const
If matrix indices are in the local range, this function returns true. Otherwise, this function return...
size_t getNodeNumEntries() const
Returns the local number of entries in this matrix.
TpetraCrsMatrix(const Teuchos::RCP< const CrsGraph< LocalOrdinal, GlobalOrdinal, Node > > &graph, const Teuchos::RCP< Teuchos::ParameterList > &params=Teuchos::null)
Constructor specifying a previously constructed graph.
TpetraCrsMatrix(const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &rowMap, size_t maxNumEntriesPerRow, ProfileType pftype=DynamicProfile, const Teuchos::RCP< Teuchos::ParameterList > &params=Teuchos::null)
Constructor specifying fixed number of entries for each row.
void getLocalRowView(LocalOrdinal LocalRow, ArrayView< const LocalOrdinal > &indices, ArrayView< const Scalar > &values) const
Extract a const, non-persisting view of local indices in a specified row of the matrix.
size_t getNodeNumEntries() const
Returns the local number of entries in this matrix.
void setObjectLabel(const std::string &objectLabel)
bool haveGlobalConstants() const
Returns true if globalConstants have been computed; false otherwise.
void setAllValues(const ArrayRCP< size_t > &rowptr, const ArrayRCP< LocalOrdinal > &colind, const ArrayRCP< Scalar > &values)
Sets the 1D pointer arrays of the graph.
void getGlobalRowCopy(GlobalOrdinal GlobalRow, const ArrayView< GlobalOrdinal > &indices, const ArrayView< Scalar > &values, size_t &numEntries) const
Extract a list of entries in a specified global row of this matrix. Put into pre-allocated storage...
void removeEmptyProcessesInPlace(const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &newMap)
const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getRowMap() const
Returns the Map that describes the row distribution in this matrix.
global_size_t getGlobalNumRows() const
Number of global elements in the row map of this matrix.
RCP< const CrsGraph< LocalOrdinal, GlobalOrdinal, Node > > getCrsGraph() const
Returns the CrsGraph associated with this matrix.
void doExport(const DistObject< char, LocalOrdinal, GlobalOrdinal, Node > &dest, const Export< LocalOrdinal, GlobalOrdinal, Node > &exporter, CombineMode CM)
Export (using an Importer).
void scale(const Scalar &alpha)
Scale the current values of a matrix, this = alpha*this.
#define XPETRA_TPETRA_ETI_EXCEPTION(cl, obj, go, node)
void fillComplete(const RCP< ParameterList > &params=null)
Signal that data entry is complete.
void getLocalDiagCopy(Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &diag, const Teuchos::ArrayView< const size_t > &offsets) const
Get a copy of the diagonal entries owned by this node, with local row indices.
ScalarTraits< Scalar >::magnitudeType getFrobeniusNorm() const
Returns the Frobenius norm of the matrix.
size_t getNodeNumCols() const
Returns the number of columns connected to the locally owned rows of this matrix. ...
void leftScale(const Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &x)
Left scale operator with given vector values.
std::string description() const
A simple one-line description of this object.
void getAllValues(ArrayRCP< const size_t > &rowptr, ArrayRCP< const LocalOrdinal > &colind, ArrayRCP< const Scalar > &values) const
Gets the 1D pointer arrays of the graph.
void replaceDomainMapAndImporter(const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &newDomainMap, Teuchos::RCP< const Import< LocalOrdinal, GlobalOrdinal, Node > > &newImporter)
Replaces the current domainMap and importer with the user-specified objects.
void insertLocalValues(LocalOrdinal localRow, const ArrayView< const LocalOrdinal > &cols, const ArrayView< const Scalar > &vals)
Insert matrix entries, using local IDs.
size_t getNodeMaxNumRowEntries() const
Returns the maximum number of entries across all rows/columns on this node.
bool isFillComplete() const
Returns true if the matrix is in compute mode, i.e. if fillComplete() has been called.
bool isLocallyIndexed() const
If matrix indices are in the local range, this function returns true. Otherwise, this function return...
void setAllToScalar(const Scalar &alpha)
Set all matrix entries equal to scalarThis.
void insertLocalValues(LocalOrdinal localRow, const ArrayView< const LocalOrdinal > &cols, const ArrayView< const Scalar > &vals)
Insert matrix entries, using local IDs.
TpetraVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > TpetraVectorClass
const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getRowMap() const
Returns the Map that describes the row distribution in this matrix.
void getLocalDiagOffsets(Teuchos::ArrayRCP< size_t > &offsets) const
Get offsets of the diagonal entries in the matrix.
std::string description() const
A simple one-line description of this object.
bool isFillComplete() const
Returns true if the matrix is in compute mode, i.e. if fillComplete() has been called.
TpetraCrsMatrix(const Teuchos::RCP< const CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > &sourceMatrix, const Import< LocalOrdinal, GlobalOrdinal, Node > &importer, const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &domainMap=Teuchos::null, const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &rangeMap=Teuchos::null, const Teuchos::RCP< Teuchos::ParameterList > &params=Teuchos::null)
Constructor for a fused import.
TpetraImport< LocalOrdinal, GlobalOrdinal, Node > TpetraImportClass
const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getDomainMap() const
Returns the Map associated with the domain of this operator. This will be null until fillComplete() i...
TpetraCrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > TpetraCrsMatrixClass
bool haveGlobalConstants() const
Returns true if globalConstants have been computed; false otherwise.
void fillComplete(const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &domainMap, const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &rangeMap, const RCP< ParameterList > &params=null)
Signal that data entry is complete, specifying domain and range maps.
RCP< const Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > getTpetra_CrsMatrix() const
Get the underlying Tpetra matrix.
RCP< const CrsGraph< LocalOrdinal, GlobalOrdinal, Node > > getCrsGraph() const
Returns the CrsGraph associated with this matrix.
Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getMap() const
Implements DistObject interface.
void getLocalDiagOffsets(Teuchos::ArrayRCP< size_t > &offsets) const
Get offsets of the diagonal entries in the matrix.
bool isGloballyIndexed() const
If matrix indices are in the global range, this function returns true. Otherwise, this function retur...
TpetraCrsMatrix(const Teuchos::RCP< const CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > &sourceMatrix, const Export< LocalOrdinal, GlobalOrdinal, Node > &RowExporter, const Teuchos::RCP< const Export< LocalOrdinal, GlobalOrdinal, Node > > DomainExporter, const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &domainMap, const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &rangeMap, const Teuchos::RCP< Teuchos::ParameterList > &params)
Constructor for a fused export.
TpetraCrsMatrix(const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &rowMap, const ArrayRCP< const size_t > &NumEntriesPerRowToAlloc, ProfileType pftype=DynamicProfile, const Teuchos::RCP< Teuchos::ParameterList > &params=Teuchos::null)
Constructor specifying (possibly different) number of entries in each row.
const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getColMap() const
Returns the Map that describes the column distribution in this matrix.
global_size_t getGlobalNumRows() const
Number of global elements in the row map of this matrix.
void doExport(const DistObject< char, LocalOrdinal, GlobalOrdinal, Node > &dest, const Import< LocalOrdinal, GlobalOrdinal, Node > &importer, CombineMode CM)
Export.
void getLocalRowCopy(LocalOrdinal LocalRow, const ArrayView< LocalOrdinal > &Indices, const ArrayView< Scalar > &Values, size_t &NumEntries) const
Extract a list of entries in a specified local row of the matrix. Put into storage allocated by calli...
TpetraCrsMatrix(const Teuchos::RCP< Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > &mtx)
TpetraCrsMatrix constructor to wrap a Tpetra::CrsMatrix object.
void doImport(const DistObject< char, LocalOrdinal, GlobalOrdinal, Node > &source, const Export< LocalOrdinal, GlobalOrdinal, Node > &exporter, CombineMode CM)
Import (using an Exporter).
void apply(const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &X, MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &Y, Teuchos::ETransp mode=Teuchos::NO_TRANS, Scalar alpha=ScalarTraits< Scalar >::one(), Scalar beta=ScalarTraits< Scalar >::zero()) const
Computes the sparse matrix-multivector multiplication.
void removeEmptyProcessesInPlace(const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &newMap)
TpetraCrsMatrix(const Teuchos::RCP< const CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > &sourceMatrix, const Import< LocalOrdinal, GlobalOrdinal, Node > &importer, const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &domainMap=Teuchos::null, const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &rangeMap=Teuchos::null, const Teuchos::RCP< Teuchos::ParameterList > &params=Teuchos::null)
Constructor for a fused import.
void describe(Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel=Teuchos::Describable::verbLevel_default) const
Print the object with some verbosity level to an FancyOStream object.
TpetraCrsMatrix(const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &rowMap, const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &colMap, const ArrayRCP< const size_t > &NumEntriesPerRowToAlloc, ProfileType pftype=DynamicProfile, const Teuchos::RCP< Teuchos::ParameterList > &params=Teuchos::null)
Constructor specifying column Map and number of entries in each row.
void expertStaticFillComplete(const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &domainMap, const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &rangeMap, const RCP< const Import< LocalOrdinal, GlobalOrdinal, Node > > &importer=Teuchos::null, const RCP< const Export< LocalOrdinal, GlobalOrdinal, Node > > &exporter=Teuchos::null, const RCP< ParameterList > &params=Teuchos::null)
Expert static fill complete.
void allocateAllValues(size_t numNonZeros, ArrayRCP< size_t > &rowptr, ArrayRCP< LocalOrdinal > &colind, ArrayRCP< Scalar > &values)
Allocates and returns ArrayRCPs of the Crs arrays — This is an Xpetra-only routine.
RCP< const CrsGraph< LocalOrdinal, GlobalOrdinal, Node > > getCrsGraph() const
Returns the CrsGraph associated with this matrix.
size_t getNodeNumCols() const
Returns the number of columns connected to the locally owned rows of this matrix. ...
void expertStaticFillComplete(const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &domainMap, const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &rangeMap, const RCP< const Import< LocalOrdinal, GlobalOrdinal, Node > > &importer=Teuchos::null, const RCP< const Export< LocalOrdinal, GlobalOrdinal, Node > > &exporter=Teuchos::null, const RCP< ParameterList > &params=Teuchos::null)
Expert static fill complete.
void getGlobalRowView(GlobalOrdinal GlobalRow, ArrayView< const GlobalOrdinal > &indices, ArrayView< const Scalar > &values) const
Extract a const, non-persisting view of global indices in a specified row of the matrix.
TpetraVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > TpetraVectorClass
void resumeFill(const RCP< ParameterList > &params=null)
void describe(Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel=Teuchos::Describable::verbLevel_default) const
Print the object with some verbosity level to an FancyOStream object.
void rightScale(const Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &x)
Right scale operator with given vector values.
void scale(const Scalar &alpha)
Scale the current values of a matrix, this = alpha*this.
void replaceDomainMapAndImporter(const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &newDomainMap, Teuchos::RCP< const Import< LocalOrdinal, GlobalOrdinal, Node > > &newImporter)
Replaces the current domainMap and importer with the user-specified objects.
size_t getGlobalMaxNumRowEntries() const
Returns the maximum number of entries across all rows/columns on all nodes.
void getGlobalRowView(GlobalOrdinal GlobalRow, ArrayView< const GlobalOrdinal > &indices, ArrayView< const Scalar > &values) const
Extract a const, non-persisting view of global indices in a specified row of the matrix.
TpetraCrsMatrix(const Teuchos::RCP< const CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > &sourceMatrix, const Export< LocalOrdinal, GlobalOrdinal, Node > &RowExporter, const Teuchos::RCP< const Export< LocalOrdinal, GlobalOrdinal, Node > > DomainExporter, const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &domainMap, const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &rangeMap, const Teuchos::RCP< Teuchos::ParameterList > &params)
Constructor for a fused export.
size_t global_size_t
Global size_t object.
void doExport(const DistObject< char, LocalOrdinal, GlobalOrdinal, Node > &dest, const Export< LocalOrdinal, GlobalOrdinal, Node > &exporter, CombineMode CM)
Export (using an Importer).
void insertGlobalValues(GlobalOrdinal globalRow, const ArrayView< const GlobalOrdinal > &cols, const ArrayView< const Scalar > &vals)
Insert matrix entries, using global IDs.
global_size_t getGlobalNumCols() const
Number of global columns in the matrix.
static const EVerbosityLevel verbLevel_default
TpetraImport< LocalOrdinal, GlobalOrdinal, Node > TpetraImportClass
virtual ~TpetraCrsMatrix()
Destructor.
const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getRangeMap() const
Returns the Map associated with the range of this operator, which must be compatible with Y...
TpetraCrsMatrix(const Teuchos::RCP< const CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > &sourceMatrix, const Export< LocalOrdinal, GlobalOrdinal, Node > &exporter, const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &domainMap=Teuchos::null, const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &rangeMap=Teuchos::null, const Teuchos::RCP< Teuchos::ParameterList > &params=Teuchos::null)
Constructor for a fused export.
static magnitudeType magnitude(T a)
bool hasMatrix() const
Does this have an underlying matrix.
void getLocalRowCopy(LocalOrdinal LocalRow, const ArrayView< LocalOrdinal > &Indices, const ArrayView< Scalar > &Values, size_t &NumEntries) const
Extract a list of entries in a specified local row of the matrix. Put into storage allocated by calli...
bool isGloballyIndexed() const
If matrix indices are in the global range, this function returns true. Otherwise, this function retur...
TpetraCrsMatrix(const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &rowMap, const ArrayRCP< const size_t > &NumEntriesPerRowToAlloc, ProfileType pftype=DynamicProfile, const Teuchos::RCP< Teuchos::ParameterList > &params=Teuchos::null)
Constructor specifying (possibly different) number of entries in each row.
size_t getNodeMaxNumRowEntries() const
Returns the maximum number of entries across all rows/columns on this node.
void getGlobalRowView(GlobalOrdinal GlobalRow, ArrayView< const GlobalOrdinal > &indices, ArrayView< const Scalar > &values) const
Extract a const, non-persisting view of global indices in a specified row of the matrix.
void setAllValues(const ArrayRCP< size_t > &rowptr, const ArrayRCP< LocalOrdinal > &colind, const ArrayRCP< Scalar > &values)
Sets the 1D pointer arrays of the graph.
bool isFillComplete() const
Returns true if the matrix is in compute mode, i.e. if fillComplete() has been called.
void describe(Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel=Teuchos::Describable::verbLevel_default) const
Print the object with some verbosity level to an FancyOStream object.
ScalarTraits< Scalar >::magnitudeType getFrobeniusNorm() const
Returns the Frobenius norm of the matrix.
#define TEUCHOS_UNREACHABLE_RETURN(dummyReturnVal)
void replaceGlobalValues(GlobalOrdinal globalRow, const ArrayView< const GlobalOrdinal > &cols, const ArrayView< const Scalar > &vals)
Replace matrix entries, using global IDs.
TpetraCrsMatrix(const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &rowMap, const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &colMap, size_t maxNumEntriesPerRow, ProfileType pftype=DynamicProfile, const Teuchos::RCP< Teuchos::ParameterList > &params=Teuchos::null)
Constructor specifying column Map and fixed number of entries for each row.
void rightScale(const Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &x)
Right scale operator with given vector values.
RCP< const Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > getTpetra_CrsMatrix() const
Get the underlying Tpetra matrix.
void getAllValues(ArrayRCP< const size_t > &rowptr, ArrayRCP< const LocalOrdinal > &colind, ArrayRCP< const Scalar > &values) const
Gets the 1D pointer arrays of the graph.
bool hasMatrix() const
Does this have an underlying matrix.
global_size_t getGlobalNumEntries() const
Returns the global number of entries in this matrix.
bool haveGlobalConstants() const
Returns true if globalConstants have been computed; false otherwise.
global_size_t getGlobalNumCols() const
Number of global columns in the matrix.
bool supportsRowViews() const
Returns true if getLocalRowView() and getGlobalRowView() are valid for this class.
global_size_t getGlobalNumEntries() const
Returns the global number of entries in this matrix.
size_t getNodeNumEntries() const
Returns the local number of entries in this matrix.
void scale(const Scalar &alpha)
Scale the current values of a matrix, this = alpha*this.
TpetraCrsMatrix(const Teuchos::RCP< const CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > &sourceMatrix, const Import< LocalOrdinal, GlobalOrdinal, Node > &RowImporter, const Teuchos::RCP< const Import< LocalOrdinal, GlobalOrdinal, Node > > DomainImporter, const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &domainMap, const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &rangeMap, const Teuchos::RCP< Teuchos::ParameterList > &params)
Constructor for a fused import.
RCP< Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > getTpetra_CrsMatrixNonConst() const
Get the underlying Tpetra matrix.
void insertGlobalValues(GlobalOrdinal globalRow, const ArrayView< const GlobalOrdinal > &cols, const ArrayView< const Scalar > &vals)
Insert matrix entries, using global IDs.
void replaceDiag(const Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &diag)
Replace the diagonal entries of the matrix.
const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getRangeMap() const
Returns the Map associated with the range of this operator, which must be compatible with Y...
void expertStaticFillComplete(const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &domainMap, const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &rangeMap, const RCP< const Import< LocalOrdinal, GlobalOrdinal, Node > > &importer=Teuchos::null, const RCP< const Export< LocalOrdinal, GlobalOrdinal, Node > > &exporter=Teuchos::null, const RCP< ParameterList > &params=Teuchos::null)
Expert static fill complete.
bool hasMatrix() const
Does this have an underlying matrix.
void replaceDiag(const Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &diag)
Replace the diagonal entries of the matrix.
const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getRowMap() const
Returns the Map that describes the row distribution in this matrix.
bool supportsRowViews() const
Returns true if getLocalRowView() and getGlobalRowView() are valid for this class.
CombineMode
Xpetra::Combine Mode enumerable type.
TpetraCrsMatrix(const TpetraCrsMatrix &matrix)
Deep copy constructor.
void apply(const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &X, MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &Y, Teuchos::ETransp mode=Teuchos::NO_TRANS, Scalar alpha=ScalarTraits< Scalar >::one(), Scalar beta=ScalarTraits< Scalar >::zero()) const
Computes the sparse matrix-multivector multiplication.
void leftScale(const Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &x)
Left scale operator with given vector values.
RCP< const Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > getTpetra_CrsMatrix() const
Get the underlying Tpetra matrix.
void fillComplete(const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &domainMap, const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &rangeMap, const RCP< ParameterList > &params=null)
Signal that data entry is complete, specifying domain and range maps.
const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getRangeMap() const
Returns the Map associated with the range of this operator, which must be compatible with Y...
TpetraCrsMatrix(const Teuchos::RCP< const CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > &sourceMatrix, const Import< LocalOrdinal, GlobalOrdinal, Node > &RowImporter, const Teuchos::RCP< const Import< LocalOrdinal, GlobalOrdinal, Node > > DomainImporter, const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &domainMap, const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &rangeMap, const Teuchos::RCP< Teuchos::ParameterList > &params)
Constructor for a fused import.
void getLocalDiagCopy(Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &diag) const
Get a copy of the diagonal entries owned by this node, with local row idices.
void getLocalRowView(LocalOrdinal LocalRow, ArrayView< const LocalOrdinal > &indices, ArrayView< const Scalar > &values) const
Extract a const, non-persisting view of local indices in a specified row of the matrix.
Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getMap() const
Implements DistObject interface.
void doExport(const DistObject< char, LocalOrdinal, GlobalOrdinal, Node > &dest, const Import< LocalOrdinal, GlobalOrdinal, Node > &importer, CombineMode CM)
Export.
const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getColMap() const
Returns the Map that describes the column distribution in this matrix.
void getAllValues(ArrayRCP< const size_t > &rowptr, ArrayRCP< const LocalOrdinal > &colind, ArrayRCP< const Scalar > &values) const
Gets the 1D pointer arrays of the graph.
void replaceGlobalValues(GlobalOrdinal globalRow, const ArrayView< const GlobalOrdinal > &cols, const ArrayView< const Scalar > &vals)
Replace matrix entries, using global IDs.
void allocateAllValues(size_t numNonZeros, ArrayRCP< size_t > &rowptr, ArrayRCP< LocalOrdinal > &colind, ArrayRCP< Scalar > &values)
Allocates and returns ArrayRCPs of the Crs arrays — This is an Xpetra-only routine.
const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getColMap() const
Returns the Map that describes the column distribution in this matrix.
void leftScale(const Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &x)
Left scale operator with given vector values.
void doImport(const DistObject< char, LocalOrdinal, GlobalOrdinal, Node > &source, const Import< LocalOrdinal, GlobalOrdinal, Node > &importer, CombineMode CM)
Import.
size_t getNodeNumRows() const
Returns the number of matrix rows owned on the calling node.
size_t getNodeNumRows() const
Returns the number of matrix rows owned on the calling node.
bool isGloballyIndexed() const
If matrix indices are in the global range, this function returns true. Otherwise, this function retur...
void replaceDiag(const Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &diag)
Replace the diagonal entries of the matrix.
bool isFillActive() const
Returns true if the matrix is in edit mode.
TpetraCrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > TpetraCrsMatrixClass
size_t getNumEntriesInLocalRow(LocalOrdinal localRow) const
Returns the current number of entries on this node in the specified local row.
void fillComplete(const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &domainMap, const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &rangeMap, const RCP< ParameterList > &params=null)
Signal that data entry is complete, specifying domain and range maps.
size_t getGlobalMaxNumRowEntries() const
Returns the maximum number of entries across all rows/columns on all nodes.
ScalarTraits< Scalar >::magnitudeType getFrobeniusNorm() const
Returns the Frobenius norm of the matrix.
bool isFillActive() const
Returns true if the matrix is in edit mode.
TpetraCrsMatrix(const TpetraCrsMatrix &matrix)
Deep copy constructor.
RCP< Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > getTpetra_CrsMatrixNonConst() const
Get the underlying Tpetra matrix.
void insertGlobalValues(GlobalOrdinal globalRow, const ArrayView< const GlobalOrdinal > &cols, const ArrayView< const Scalar > &vals)
Insert matrix entries, using global IDs.
TpetraCrsMatrix(const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &rowMap, size_t maxNumEntriesPerRow, ProfileType pftype=DynamicProfile, const Teuchos::RCP< Teuchos::ParameterList > &params=Teuchos::null)
Constructor specifying fixed number of entries for each row.
void getGlobalRowCopy(GlobalOrdinal GlobalRow, const ArrayView< GlobalOrdinal > &indices, const ArrayView< Scalar > &values, size_t &numEntries) const
Extract a list of entries in a specified global row of this matrix. Put into pre-allocated storage...
size_t getNodeMaxNumRowEntries() const
Returns the maximum number of entries across all rows/columns on this node.
void getGlobalRowCopy(GlobalOrdinal GlobalRow, const ArrayView< GlobalOrdinal > &indices, const ArrayView< Scalar > &values, size_t &numEntries) const
Extract a list of entries in a specified global row of this matrix. Put into pre-allocated storage...
void setAllToScalar(const Scalar &alpha)
Set all matrix entries equal to scalarThis.
bool isFillActive() const
Returns true if the matrix is in edit mode.
bool supportsRowViews() const
Returns true if getLocalRowView() and getGlobalRowView() are valid for this class.
TpetraExport< LocalOrdinal, GlobalOrdinal, Node > TpetraExportClass
size_t getNumEntriesInLocalRow(LocalOrdinal localRow) const
Returns the current number of entries on this node in the specified local row.
void doImport(const DistObject< char, LocalOrdinal, GlobalOrdinal, Node > &source, const Export< LocalOrdinal, GlobalOrdinal, Node > &exporter, CombineMode CM)
Import (using an Exporter).
void getLocalDiagCopy(Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &diag, const Teuchos::ArrayView< const size_t > &offsets) const
Get a copy of the diagonal entries owned by this node, with local row indices.
size_t getNodeNumCols() const
Returns the number of columns connected to the locally owned rows of this matrix. ...