42 #ifndef __Epetra_TsqrAdaptor_hpp
43 #define __Epetra_TsqrAdaptor_hpp
59 #include <Tpetra_ConfigDefs.hpp>
61 #if defined(HAVE_TPETRA_EPETRA) && defined(HAVE_TPETRA_TSQR)
63 #include <Kokkos_DefaultNode.hpp>
64 #include <Tsqr_NodeTsqrFactory.hpp>
66 #include <Tsqr_DistTsqr.hpp>
67 #include <Epetra_Comm.h>
71 #include <Epetra_MultiVector.h>
72 #include <Teuchos_ParameterListAcceptorDefaultBase.hpp>
102 class TsqrAdaptor :
public Teuchos::ParameterListAcceptorDefaultBase {
104 typedef Epetra_MultiVector MV;
112 typedef double scalar_type;
120 typedef int ordinal_type;
136 typedef Teuchos::SerialDenseMatrix<ordinal_type, scalar_type> dense_matrix_type;
143 typedef double magnitude_type;
146 typedef TSQR::MatView<ordinal_type, scalar_type> matview_type;
147 typedef TSQR::NodeTsqrFactory<node_type, scalar_type, ordinal_type> node_tsqr_factory_type;
150 typedef node_tsqr_factory_type::node_tsqr_type node_tsqr_type;
151 typedef TSQR::DistTsqr<ordinal_type, scalar_type> dist_tsqr_type;
152 typedef TSQR::Tsqr<ordinal_type, scalar_type, node_tsqr_type> tsqr_type;
161 TsqrAdaptor (
const Teuchos::RCP<Teuchos::ParameterList>& plist) :
162 nodeTsqr_ (new node_tsqr_type),
163 distTsqr_ (new dist_tsqr_type),
164 tsqr_ (new tsqr_type (nodeTsqr_, distTsqr_)),
167 setParameterList (plist);
172 nodeTsqr_ (new node_tsqr_type),
173 distTsqr_ (new dist_tsqr_type),
174 tsqr_ (new tsqr_type (nodeTsqr_, distTsqr_)),
177 setParameterList (Teuchos::null);
180 Teuchos::RCP<const Teuchos::ParameterList>
181 getValidParameters ()
const
185 using Teuchos::ParameterList;
186 using Teuchos::parameterList;
188 if (defaultParams_.is_null()) {
189 RCP<ParameterList> params = parameterList (
"TSQR implementation");
190 params->set (
"NodeTsqr", *(nodeTsqr_->getValidParameters ()));
191 params->set (
"DistTsqr", *(distTsqr_->getValidParameters ()));
192 defaultParams_ = params;
194 return defaultParams_;
198 setParameterList (
const Teuchos::RCP<Teuchos::ParameterList>& plist)
200 using Teuchos::ParameterList;
201 using Teuchos::parameterList;
203 using Teuchos::sublist;
205 RCP<ParameterList> params = plist.is_null() ?
206 parameterList (*getValidParameters ()) : plist;
207 nodeTsqr_->setParameterList (sublist (params,
"NodeTsqr"));
208 distTsqr_->setParameterList (sublist (params,
"DistTsqr"));
210 this->setMyParamList (params);
235 factorExplicit (MV& A,
237 dense_matrix_type& R,
238 const bool forceNonnegativeDiagonal=
false)
242 scalar_type*
const A_ptr = A.Values ();
243 scalar_type*
const Q_ptr = Q.Values ();
244 scalar_type*
const R_ptr = R.values ();
245 const ordinal_type numRows = A.MyLength ();
246 const ordinal_type numCols = A.NumVectors ();
247 const ordinal_type lda = A.Stride ();
248 const ordinal_type ldq = Q.Stride ();
249 const ordinal_type ldr = R.stride ();
251 const bool contiguousCacheBlocks =
false;
252 tsqr_->factorExplicitRaw (numRows, numCols, A_ptr, lda,
253 Q_ptr, ldq, R_ptr, ldr,
254 contiguousCacheBlocks,
255 forceNonnegativeDiagonal);
290 dense_matrix_type& R,
291 const magnitude_type& tol)
293 TEUCHOS_TEST_FOR_EXCEPTION
294 (! Q.ConstantStride (), std::invalid_argument,
"TsqrAdaptor::"
295 "revealRank: Input MultiVector Q must have constant stride.");
300 return tsqr_->revealRankRaw (Q.MyLength (), Q.NumVectors (),
301 Q.Values (), Q.Stride (),
302 R.values (), R.stride (), tol,
false);
307 Teuchos::RCP<node_tsqr_type> nodeTsqr_;
310 Teuchos::RCP<dist_tsqr_type> distTsqr_;
313 Teuchos::RCP<tsqr_type> tsqr_;
316 mutable Teuchos::RCP<const Teuchos::ParameterList> defaultParams_;
340 prepareTsqr (
const MV& mv)
343 prepareDistTsqr (mv);
344 prepareNodeTsqr (mv);
353 prepareNodeTsqr (
const MV& mv)
358 Teuchos::ParameterList plist;
359 Teuchos::RCP<node_type> node (
new node_type (plist));
360 node_tsqr_factory_type::prepareNodeTsqr (nodeTsqr_, node);
370 prepareDistTsqr (
const MV& mv)
374 using TSQR::Epetra::makeTsqrMessenger;
375 typedef TSQR::MessengerBase<scalar_type> base_mess_type;
380 RCP<const Epetra_Comm> comm = rcp (mv.Comm().Clone());
381 RCP<base_mess_type> messBase = makeTsqrMessenger<scalar_type> (comm);
382 distTsqr_->init (messBase);
388 #endif // defined(HAVE_TPETRA_EPETRA) && defined(HAVE_TPETRA_TSQR)
390 #endif // __Epetra_TsqrAdaptor_hpp
Function for wrapping Epetra_Comm in a communicator wrapper that TSQR can use.
::Kokkos::Compat::KokkosDeviceWrapperNode< execution_space > node_type
Default value of Node template parameter.