Stokhos Package Browser (Single Doxygen Collection)  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Tpetra_TsqrAdaptor_MP_Vector.hpp
Go to the documentation of this file.
1 // @HEADER
2 // *****************************************************************************
3 // Tpetra: Templated Linear Algebra Services Package
4 //
5 // Copyright 2008 NTESS and the Tpetra contributors.
6 // SPDX-License-Identifier: BSD-3-Clause
7 // *****************************************************************************
8 // @HEADER
9 
10 #ifndef TPETRA_TSQR_ADAPTOR_MP_VECTOR_HPP
11 #define TPETRA_TSQR_ADAPTOR_MP_VECTOR_HPP
12 
13 #include <Tpetra_ConfigDefs.hpp> // HAVE_TPETRA_TSQR, etc.
14 
15 #ifdef HAVE_TPETRA_TSQR
16 
18 
19 # include "Tsqr_NodeTsqrFactory.hpp" // create intranode TSQR object
20 # include "Tsqr.hpp" // full (internode + intranode) TSQR
21 # include "Tsqr_DistTsqr.hpp" // internode TSQR
22 // Subclass of TSQR::MessengerBase, implemented using Teuchos
23 // communicator template helper functions
24 # include "Tsqr_TeuchosMessenger.hpp"
25 # include "Tpetra_MultiVector.hpp"
27 # include <stdexcept>
28 
29 // Base TsqrAdator template we will specialize
30 # include "Tpetra_TsqrAdaptor.hpp"
31 
32 namespace Tpetra {
33 
40  template <class Storage, class LO, class GO, class Node>
41  class TsqrAdaptor< Tpetra::MultiVector< Sacado::MP::Vector<Storage>,
42  LO, GO, Node > > :
44  public:
45  typedef Tpetra::MultiVector< Sacado::MP::Vector<Storage>, LO, GO, Node > MV;
46  typedef typename MV::scalar_type mp_scalar_type;
47 
48  // For Sacado::MP::Vector< Storage<Ordinal,Scalar,Device> > this is Scalar
50  typedef typename mp_scalar_type::ordinal_type mp_ordinal_type;
51  typedef typename MV::local_ordinal_type ordinal_type;
53  typedef typename Teuchos::ScalarTraits<scalar_type>::magnitudeType magnitude_type;
54 
55  private:
56  using node_tsqr_factory_type =
57  TSQR::NodeTsqrFactory<scalar_type, ordinal_type,
58  typename MV::device_type>;
59  using node_tsqr_type = TSQR::NodeTsqr<ordinal_type, scalar_type>;
60  using dist_tsqr_type = TSQR::DistTsqr<ordinal_type, scalar_type>;
61  using tsqr_type = TSQR::Tsqr<ordinal_type, scalar_type>;
62 
63  public:
70  TsqrAdaptor (const Teuchos::RCP<Teuchos::ParameterList>& plist) :
71  nodeTsqr_ (node_tsqr_factory_type::getNodeTsqr ()),
72  distTsqr_ (new dist_tsqr_type),
73  tsqr_ (new tsqr_type (nodeTsqr_, distTsqr_)),
74  ready_ (false)
75  {
76  setParameterList (plist);
77  }
78 
80  TsqrAdaptor () :
81  nodeTsqr_ (node_tsqr_factory_type::getNodeTsqr ()),
82  distTsqr_ (new dist_tsqr_type),
83  tsqr_ (new tsqr_type (nodeTsqr_, distTsqr_)),
84  ready_ (false)
85  {
86  setParameterList (Teuchos::null);
87  }
88 
90  getValidParameters () const
91  {
92  using Teuchos::RCP;
93  using Teuchos::rcp;
95  using Teuchos::parameterList;
96 
97  if (defaultParams_.is_null()) {
98  RCP<ParameterList> params = parameterList ("TSQR implementation");
99  params->set ("NodeTsqr", *(nodeTsqr_->getValidParameters ()));
100  params->set ("DistTsqr", *(distTsqr_->getValidParameters ()));
101  defaultParams_ = params;
102  }
103  return defaultParams_;
104  }
105 
106  void
107  setParameterList (const Teuchos::RCP<Teuchos::ParameterList>& plist)
108  {
110  using Teuchos::parameterList;
111  using Teuchos::RCP;
112  using Teuchos::sublist;
113 
114  RCP<ParameterList> params = plist.is_null() ?
115  parameterList (*getValidParameters ()) : plist;
116  nodeTsqr_->setParameterList (sublist (params, "NodeTsqr"));
117  distTsqr_->setParameterList (sublist (params, "DistTsqr"));
118 
119  this->setMyParamList (params);
120  }
121 
143  void
144  factorExplicit (MV& A,
145  MV& Q,
146  dense_matrix_type& R,
147  const bool forceNonnegativeDiagonal=false)
148  {
149  prepareTsqr (Q); // Finish initializing TSQR.
150 
151  ordinal_type numRows;
152  ordinal_type numCols;
153  ordinal_type LDA;
154  ordinal_type LDQ;
155  scalar_type* A_ptr;
156  scalar_type* Q_ptr;
157 
158  getNonConstView (numRows, numCols, A_ptr, LDA, A);
159  getNonConstView (numRows, numCols, Q_ptr, LDQ, Q);
160  const bool contiguousCacheBlocks = false;
161  tsqr_->factorExplicitRaw (numRows, numCols, A_ptr, LDA,
162  Q_ptr, LDQ, R.values (), R.stride (),
163  contiguousCacheBlocks,
164  forceNonnegativeDiagonal);
165  }
166 
197  int
198  revealRank (MV& Q,
199  dense_matrix_type& R,
200  const magnitude_type& tol)
201  {
202  prepareTsqr (Q); // Finish initializing TSQR.
203 
204  // FIXME (mfh 18 Oct 2010) Check Teuchos::Comm<int> object in Q
205  // to make sure it is the same communicator as the one we are
206  // using in our dist_tsqr_type implementation.
207 
208  ordinal_type numRows;
209  ordinal_type numCols;
210  scalar_type* Q_ptr;
211  ordinal_type LDQ;
212  getNonConstView (numRows, numCols, Q_ptr, LDQ, Q);
213  const bool contiguousCacheBlocks = false;
214  return tsqr_->revealRankRaw (numRows, numCols, Q_ptr, LDQ,
215  R.values (), R.stride (), tol,
216  contiguousCacheBlocks);
217  }
218 
219  private:
222 
225 
228 
230  mutable Teuchos::RCP<const Teuchos::ParameterList> defaultParams_;
231 
233  bool ready_;
234 
255  void
256  prepareTsqr (const MV& mv)
257  {
258  if (! ready_) {
259  prepareDistTsqr (mv);
260  ready_ = true;
261  }
262  }
263 
270  void
271  prepareDistTsqr (const MV& mv)
272  {
273  using Teuchos::RCP;
274  using Teuchos::rcp_implicit_cast;
275  typedef TSQR::TeuchosMessenger<scalar_type> mess_type;
276  typedef TSQR::MessengerBase<scalar_type> base_mess_type;
277 
278  RCP<const Teuchos::Comm<int> > comm = mv.getMap()->getComm();
279  RCP<mess_type> mess (new mess_type (comm));
280  RCP<base_mess_type> messBase = rcp_implicit_cast<base_mess_type> (mess);
281  distTsqr_->init (messBase);
282  }
283 
290  static void
291  getNonConstView (ordinal_type& numRows,
292  ordinal_type& numCols,
293  scalar_type*& A_ptr,
294  ordinal_type& LDA,
295  const MV& A)
296  {
297  // FIXME (mfh 25 Oct 2010) We should be able to run TSQR even if
298  // storage of A uses nonconstant stride internally. We would
299  // have to copy and pack into a matrix with constant stride, and
300  // then unpack on exit. For now we choose just to raise an
301  // exception.
303  (! A.isConstantStride (), std::invalid_argument,
304  "TSQR does not currently support Tpetra::MultiVector "
305  "inputs that do not have constant stride.");
306 
307  // FIXME (mfh 16 Jan 2016) When I got here, I found strides[0]
308  // instead of strides[1] for the stride. I don't think this is
309  // right. However, I don't know about these Stokhos scalar
310  // types so I'll just do what was here.
311  //
312  // STOKHOS' TYPES ARE NOT TESTED WITH TSQR REGULARLY SO IT IS
313  // POSSIBLE THAT THE ORIGINAL CODE WAS WRONG.
314 
315  typedef typename MV::dual_view_type view_type;
316  typedef typename view_type::t_dev::array_type flat_array_type;
317 
318  // Reinterpret the data as a longer array of the base scalar
319  // type. TSQR currently forbids MultiVector input with
320  // nonconstant stride, so we need not worry about that here.
321 
322  flat_array_type flat_mv = A.getLocalViewDevice();
323 
324  numRows = static_cast<ordinal_type> (flat_mv.extent(0));
325  numCols = static_cast<ordinal_type> (flat_mv.extent(1));
326  A_ptr = flat_mv.data ();
327 
328  ordinal_type strides[2];
329  flat_mv.stride (strides);
330  LDA = strides[0];
331  }
332  };
333 
334 } // namespace Tpetra
335 
336 #endif // HAVE_TPETRA_TSQR
337 
338 #endif // TPETRA_TSQR_ADAPTOR_MP_VECTOR_HPP
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
bool is_null(const RCP< T > &p)
Tpetra::KokkosClassic::DefaultNode::DefaultNodeType Node