Thyra  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Thyra_SpmdMultiVectorDefaultBase_def.hpp
1 // @HEADER
2 // *****************************************************************************
3 // Thyra: Interfaces and Support for Abstract Numerical Algorithms
4 //
5 // Copyright 2004 NTESS and the Thyra contributors.
6 // SPDX-License-Identifier: BSD-3-Clause
7 // *****************************************************************************
8 // @HEADER
9 
10 #ifndef THYRA_SPMD_MULTI_VECTOR_DEFAULT_BASE_DEF_HPP
11 #define THYRA_SPMD_MULTI_VECTOR_DEFAULT_BASE_DEF_HPP
12 
13 // disable clang warnings
14 #if defined (__clang__) && !defined (__INTEL_COMPILER)
15 #pragma clang system_header
16 #endif
17 
18 
19 #include "Thyra_SpmdMultiVectorDefaultBase_decl.hpp"
20 #include "Thyra_MultiVectorDefaultBase.hpp"
21 #include "Thyra_MultiVectorAdapterBase.hpp"
22 #include "Thyra_SpmdVectorSpaceDefaultBase.hpp"
23 #include "Thyra_DetachedMultiVectorView.hpp"
24 #include "Thyra_apply_op_helper.hpp"
25 #include "Thyra_SpmdLocalDataAccess.hpp"
26 #include "RTOpPack_SPMD_apply_op.hpp"
27 #include "RTOp_parallel_helpers.h"
28 #include "Teuchos_Workspace.hpp"
29 #include "Teuchos_dyn_cast.hpp"
30 #include "Teuchos_Time.hpp"
31 #include "Teuchos_CommHelpers.hpp"
32 
33 
34 // Define to see some timing output!
35 //#define THYRA_SPMD_MULTI_VECTOR_BASE_PRINT_TIMES
36 
37 
38 namespace Thyra {
39 
40 
41 // Constructors / initializers / accessors
42 
43 
44 template<class Scalar>
46  :in_applyOp_(false),
47  globalDim_(0),
48  localOffset_(-1),
49  localSubDim_(0),
50  numCols_(0)
51 {}
52 
53 
54 // Overridden public functions from MultiVectorAdapterBase
55 
56 
57 template<class Scalar>
60 {
61  return Teuchos::rcp_dynamic_cast<const ScalarProdVectorSpaceBase<Scalar> >(
62  this->spmdSpace(), true
63  );
64 }
65 
66 
67 // Overridden public functions from MultiVectorAdapterBase
68 
69 
70 template<class Scalar>
73 {
74  using Teuchos::outArg;
75  ArrayRCP<Scalar> localValues;
76  Ordinal leadingDim = -1;
77  this->getNonconstLocalData(outArg(localValues), outArg(leadingDim));
79  localOffset_, // globalOffset
80  localSubDim_,
81  0, // colOffset
82  numCols_,
83  localValues,
84  leadingDim
85  );
86 }
87 
88 
89 template<class Scalar>
92 {
93  using Teuchos::outArg;
94  ArrayRCP<const Scalar> localValues;
95  Ordinal leadingDim = -1;
96  this->getLocalData(outArg(localValues), outArg(leadingDim));
98  localOffset_, // globalOffset
99  localSubDim_,
100  0, // colOffset
101  numCols_,
102  localValues,
103  leadingDim
104  );
105 }
106 
107 
108 // Protected functions overridden from MultiVectorBase
109 
110 
111 template<class Scalar>
113  const RTOpPack::RTOpT<Scalar> &pri_op,
114  const ArrayView<const Ptr<const MultiVectorBase<Scalar> > > &multi_vecs,
115  const ArrayView<const Ptr<MultiVectorBase<Scalar> > > &targ_multi_vecs,
116  const ArrayView<const Ptr<RTOpPack::ReductTarget> > &reduct_objs,
117  const Ordinal pri_global_offset_in
118  ) const
119 {
120 
121  using Teuchos::dyn_cast;
122  using Teuchos::Workspace;
123  using Teuchos::rcpFromPtr;
124 
126 
127  const Ordinal numCols = this->domain()->dim();
128  const SpmdVectorSpaceBase<Scalar> &spmdSpc = *this->spmdSpace();
129 
130 #ifdef TEUCHOS_DEBUG
132  in_applyOp_, std::invalid_argument,
133  "SpmdMultiVectorDefaultBase<>::mvMultiReductApplyOpImpl(...): Error, this method is"
134  " being entered recursively which is a clear sign that one of the methods"
135  " acquireDetachedView(...), releaseDetachedView(...) or commitDetachedView(...)"
136  " was not implemented properly!"
137  );
139  "SpmdMultiVectorDefaultBase<>::mvMultiReductApplyOpImpl(...)", *this->domain(),
140  *this->range(), pri_op, multi_vecs, targ_multi_vecs, reduct_objs,
141  pri_global_offset_in);
142 #endif
143 
144  // Flag that we are in applyOp()
145  in_applyOp_ = true;
146 
147  // First see if this is a locally replicated vector in which case
148  // we treat this as a local operation only.
149  const bool locallyReplicated = spmdSpc.isLocallyReplicated();
150 
151  const Range1D local_rng(localOffset_, localOffset_+localSubDim_-1);
152  const Range1D col_rng(0, numCols-1);
153 
154  // Create sub-vector views of all of the *participating* local data
155  Workspace<RTOpPack::ConstSubMultiVectorView<Scalar> >
156  sub_multi_vecs(wss,multi_vecs.size());
157  Workspace<RTOpPack::SubMultiVectorView<Scalar> >
158  targ_sub_multi_vecs(wss,targ_multi_vecs.size());
159  for(int k = 0; k < multi_vecs.size(); ++k ) {
160  sub_multi_vecs[k] = getLocalSubMultiVectorView<Scalar>(rcpFromPtr(multi_vecs[k]));
161  sub_multi_vecs[k].setGlobalOffset(localOffset_+pri_global_offset_in);
162  }
163  for(int k = 0; k < targ_multi_vecs.size(); ++k ) {
164  targ_sub_multi_vecs[k] =
165  getNonconstLocalSubMultiVectorView<Scalar>(rcpFromPtr(targ_multi_vecs[k]));
166  targ_sub_multi_vecs[k].setGlobalOffset(localOffset_+pri_global_offset_in);
167  }
168  Workspace<RTOpPack::ReductTarget*> reduct_objs_ptr(wss, reduct_objs.size());
169  for (int k = 0; k < reduct_objs.size(); ++k) {
170  reduct_objs_ptr[k] = &*reduct_objs[k];
171  }
172 
173  // Apply the RTOp operator object (all processors must participate)
174  RTOpPack::SPMD_apply_op(
175  locallyReplicated ? NULL : spmdSpc.getComm().get(), // comm
176  pri_op, // op
177  col_rng.size(), // num_cols
178  sub_multi_vecs.size(), // multi_vecs.size()
179  sub_multi_vecs.getRawPtr(), // sub_multi_vecs
180  targ_sub_multi_vecs.size(), // targ_multi_vecs.size()
181  targ_sub_multi_vecs.getRawPtr(), // targ_sub_multi_vecs
182  reduct_objs_ptr.getRawPtr() // reduct_objs
183  );
184 
185  // Free and commit the local data
186  for(int k = 0; k < multi_vecs.size(); ++k ) {
187  sub_multi_vecs[k] = RTOpPack::ConstSubMultiVectorView<Scalar>();
188  }
189  for(int k = 0; k < targ_multi_vecs.size(); ++k ) {
190  targ_sub_multi_vecs[k] = RTOpPack::SubMultiVectorView<Scalar>();
191  }
192 
193  // Flag that we are leaving applyOp()
194  in_applyOp_ = false;
195 
196 }
197 
198 
199 template<class Scalar>
201  const Range1D &rowRng_in,
202  const Range1D &colRng_in,
204  ) const
205 {
206  using Teuchos::outArg;
207  const Range1D rowRng = validateRowRange(rowRng_in);
208  const Range1D colRng = validateColRange(colRng_in);
209  if( rowRng.lbound() < localOffset_ || localOffset_+localSubDim_-1 < rowRng.ubound() ) {
210  // rng consists of off-processor elements so use the default implementation!
212  rowRng_in,colRng_in,sub_mv
213  );
214  return;
215  }
216  ArrayRCP<const Scalar> localValues;
217  Ordinal leadingDim = 0;
218  this->getLocalData(outArg(localValues), outArg(leadingDim));
219  sub_mv->initialize(
220  rowRng.lbound(), // globalOffset
221  rowRng.size(), // subDim
222  colRng.lbound(), // colOffset
223  colRng.size(), // numSubCols
224  localValues
225  +(rowRng.lbound()-localOffset_)
226  +colRng.lbound()*leadingDim, // values
227  leadingDim // leadingDim
228  );
229 }
230 
231 
232 template<class Scalar>
235  ) const
236 {
237  if(
238  sub_mv->globalOffset() < localOffset_
239  ||
240  localOffset_+localSubDim_ < sub_mv->globalOffset()+sub_mv->subDim()
241  )
242  {
243  // Let the default implementation handle it!
245  return;
246  }
247  sub_mv->uninitialize();
248 }
249 
250 
251 template<class Scalar>
253  const Range1D &rowRng_in,
254  const Range1D &colRng_in,
256  )
257 {
258  using Teuchos::outArg;
259  const Range1D rowRng = validateRowRange(rowRng_in);
260  const Range1D colRng = validateColRange(colRng_in);
261  if(
262  rowRng.lbound() < localOffset_
263  ||
264  localOffset_+localSubDim_-1 < rowRng.ubound()
265  )
266  {
267  // rng consists of off-processor elements so use the default implementation!
269  rowRng_in, colRng_in, sub_mv
270  );
271  return;
272  }
273  ArrayRCP<Scalar> localValues;
274  Ordinal leadingDim = 0;
275  this->getNonconstLocalData(outArg(localValues), outArg(leadingDim));
276  sub_mv->initialize(
277  rowRng.lbound() // globalOffset
278  ,rowRng.size() // subDim
279  ,colRng.lbound() // colOffset
280  ,colRng.size() // numSubCols
281  ,localValues
282  +(rowRng.lbound()-localOffset_)
283  +colRng.lbound()*leadingDim // values
284  ,leadingDim // leadingDim
285  );
286 }
287 
288 
289 template<class Scalar>
292  )
293 {
294  if(
295  sub_mv->globalOffset() < localOffset_
296  ||
297  localOffset_+localSubDim_ < sub_mv->globalOffset()+sub_mv->subDim()
298  )
299  {
300  // Let the default implementation handle it!
302  return;
303  }
304  sub_mv->uninitialize();
305 }
306 
307 
308 // Protected functions overridden from MultiVectorAdapterBase
309 
310 
311 template<class Scalar>
313  const EOpTransp M_trans,
314  const MultiVectorBase<Scalar> &X,
315  const Ptr<MultiVectorBase<Scalar> > &Y,
316  const Scalar alpha,
317  const Scalar beta
318  ) const
319 {
321  using Teuchos::Workspace;
322  using Teuchos::rcpFromPtr;
324 
325 #ifdef THYRA_SPMD_MULTI_VECTOR_BASE_PRINT_TIMES
326  Teuchos::Time timerTotal("dummy",true);
327  Teuchos::Time timer("dummy");
328 #endif
329 
330  //
331  // This function performs one of two operations.
332  //
333  // The first operation (M_trans == NOTRANS) is:
334  //
335  // Y = beta * Y + alpha * M * X
336  //
337  // where Y and M have compatible (distributed?) range vector
338  // spaces and X is a locally replicated serial multi-vector. This
339  // operation does not require any global communication.
340  //
341  // The second operation (M_trans == TRANS) is:
342  //
343  // Y = beta * Y + alpha * M' * X
344  //
345  // where M and X have compatible (distributed?) range vector spaces
346  // and Y is a locally replicated serial multi-vector. This operation
347  // requires a local reduction.
348  //
349 
350  //
351  // Get spaces and validate compatibility
352  //
353 
354  // Get the SpmdVectorSpace
355  const SpmdVectorSpaceBase<Scalar> &spmdSpc = *this->spmdSpace();
356 
357  // Get the Spmd communicator
358  const RCP<const Teuchos::Comm<Ordinal> > comm = spmdSpc.getComm();
359  const int procRank = (nonnull(comm) ? comm->getRank() : 0 );
360 
361 #ifdef TEUCHOS_DEBUG
363  &Y_range = *Y->range(),
364  &X_range = *X.range();
365 // std::cout << "SpmdMultiVectorDefaultBase<Scalar>::apply(...): comm = " << comm << std::endl;
367  ( globalDim_ > localSubDim_ ) && is_null(comm), std::logic_error
368  ,"SpmdMultiVectorDefaultBase<Scalar>::apply(...MultiVectorBase<Scalar>...): Error!"
369  );
370  // ToDo: Write a good general validation function that I can call that will replace
371  // all of these TEUCHOS_TEST_FOR_EXCEPTION(...) uses
372 
375  ,"SpmdMultiVectorDefaultBase<Scalar>::apply(...MultiVectorBase<Scalar>...): Error!"
376  );
378  real_trans(M_trans)==TRANS && !spmdSpc.isCompatible(X_range), Exceptions::IncompatibleVectorSpaces
379  ,"SpmdMultiVectorDefaultBase<Scalar>::apply(...MultiVectorBase<Scalar>...): Error!"
380  );
381 #endif
382 
383  //
384  // Get explicit (local) views of Y, M and X
385  //
386 
387 #ifdef THYRA_SPMD_MULTI_VECTOR_BASE_PRINT_TIMES
388  timer.start();
389 #endif
390 
392  Y_local = getNonconstLocalSubMultiVectorView<Scalar>(rcpFromPtr(Y));
394  M_local = getLocalSubMultiVectorView<Scalar>(rcpFromRef(*this)),
395  X_local = getLocalSubMultiVectorView<Scalar>(rcpFromRef(X));
396 /*
397  DetachedMultiVectorView<Scalar>
398  Y_local(
399  *Y,
400  real_trans(M_trans)==NOTRANS ? Range1D(localOffset_,localOffset_+localSubDim_-1) : Range1D(),
401  Range1D()
402  );
403  ConstDetachedMultiVectorView<Scalar>
404  M_local(
405  *this,
406  Range1D(localOffset_,localOffset_+localSubDim_-1),
407  Range1D()
408  );
409  ConstDetachedMultiVectorView<Scalar>
410  X_local(
411  X
412  ,real_trans(M_trans)==NOTRANS ? Range1D() : Range1D(localOffset_,localOffset_+localSubDim_-1)
413  ,Range1D()
414  );
415 */
416 #ifdef THYRA_SPMD_MULTI_VECTOR_BASE_PRINT_TIMES
417  timer.stop();
418  std::cout << "\nSpmdMultiVectorDefaultBase<Scalar>::apply(...): Time for getting view = " << timer.totalElapsedTime() << " seconds\n";
419 #endif
420 #ifdef TEUCHOS_DEBUG
422  real_trans(M_trans)==NOTRANS && ( M_local.numSubCols() != X_local.subDim() || X_local.numSubCols() != Y_local.numSubCols() )
424  ,"SpmdMultiVectorDefaultBase<Scalar>::apply(...MultiVectorBase<Scalar>...): Error!"
425  );
427  real_trans(M_trans)==TRANS && ( M_local.subDim() != X_local.subDim() || X_local.numSubCols() != Y_local.numSubCols() )
429  ,"SpmdMultiVectorDefaultBase<Scalar>::apply(...MultiVectorBase<Scalar>...): Error!"
430  );
431 #endif
432 
433  //
434  // If nonlocal (i.e. M_trans==TRANS) then create temporary storage
435  // for:
436  //
437  // Y_local_tmp = alpha * M(local) * X(local) : on nonroot processes
438  //
439  // or
440  //
441  // Y_local_tmp = beta*Y_local + alpha * M(local) * X(local) : on root process (localOffset_==0)
442  //
443  // and set
444  //
445  // localBeta = ( localOffset_ == 0 ? beta : 0.0 )
446  //
447  // Above, we choose localBeta such that we will only perform
448  // Y_local = beta * Y_local + ... on one process (the root
449  // process where localOffset_==0x). Then, when we add up Y_local
450  // on all of the processors and we will get the correct result.
451  //
452  // If strictly local (i.e. M_trans == NOTRANS) then set:
453  //
454  // Y_local_tmp = Y_local
455  // localBeta = beta
456  //
457 
458 #ifdef THYRA_SPMD_MULTI_VECTOR_BASE_PRINT_TIMES
459  timer.start();
460 #endif
461 
462  // determine if both fields are locally replicated. If they are
463  // there is no need to do any reduction operations.
464  bool locallyReplicated = false;
465  {
466  bool locallyReplicated_this = spmdSpc.isLocallyReplicated();
467  bool locallyReplicated_x = locallyReplicated_this; // x must be compatible with "this"
468  bool locallyReplicated_y = false;
469 
470  Ptr<SpmdMultiVectorBase<Scalar> > spmd_Y = Teuchos::ptr_dynamic_cast<SpmdMultiVectorBase<Scalar> >(Y);
471  if (nonnull(spmd_Y))
472  locallyReplicated_y = spmd_Y->spmdSpace()->isLocallyReplicated();
473 
474  locallyReplicated = locallyReplicated_this && locallyReplicated_x && locallyReplicated_y;
475  }
476 
477  bool isNonLocalAdjoint =
478  (
479  real_trans(M_trans) == TRANS
480  &&
481  (globalDim_ > localSubDim_ || (nonnull(comm) && comm->getSize() > 1))
482  );
483 
484  if (locallyReplicated)
485  isNonLocalAdjoint = false;
486 
487  Workspace<Scalar> Y_local_tmp_store(wss, Y_local.subDim()*Y_local.numSubCols(), false);
489  Scalar localBeta;
490  if (isNonLocalAdjoint) {
491  // Nonlocal
492  Y_local_tmp.initialize(
493  0, Y_local.subDim(),
494  0, Y_local.numSubCols(),
495  Teuchos::arcpFromArrayView(Y_local_tmp_store()),
496  Y_local.subDim() // leadingDim == subDim (columns are adjacent)
497  );
498  if (procRank == 0) {
499  // Root process: Must copy Y_local into Y_local_tmp
500  for( int j = 0; j < Y_local.numSubCols(); ++j ) {
501  typedef typename ArrayRCP<const Scalar>::const_iterator Y_local_values_iter_t;
502  const Y_local_values_iter_t Y_local_j =
503  Y_local.values().begin() + Y_local.leadingDim()*j;
504  std::copy( Y_local_j, Y_local_j + Y_local.subDim(),
505  Y_local_tmp.values().begin() + Y_local_tmp.leadingDim()*j );
506  }
507  localBeta = beta;
508  }
509  else {
510  // Not the root process
511  localBeta = 0.0;
512  }
513  }
514  else {
515  // Local
516  Y_local_tmp = Y_local; // Shallow copy only!
517  localBeta = beta;
518  }
519 
520 #ifdef THYRA_SPMD_MULTI_VECTOR_BASE_PRINT_TIMES
521  timer.stop();
522  std::cout << "\nSpmdMultiVectorDefaultBase<Scalar>::apply(...): Time for setting up Y_local_tmp and localBeta = " << timer.totalElapsedTime() << " seconds\n";
523 #endif
524 
525  //
526  // Perform the local multiplication:
527  //
528  // Y(local) = localBeta * Y(local) + alpha * op(M(local)) * X(local)
529  //
530  // or in BLAS lingo:
531  //
532  // C = beta * C + alpha * op(A) * op(B)
533  //
534 
535 #ifdef THYRA_SPMD_MULTI_VECTOR_BASE_PRINT_TIMES
536  timer.start();
537 #endif
538  Teuchos::ETransp t_transp;
539  if(ST::isComplex) {
540  switch(M_trans) {
541  case NOTRANS: t_transp = Teuchos::NO_TRANS; break;
542  case TRANS: t_transp = Teuchos::TRANS; break;
543  case CONJTRANS: t_transp = Teuchos::CONJ_TRANS; break;
544  default: TEUCHOS_TEST_FOR_EXCEPT(true);
545  }
546  }
547  else {
548  switch(real_trans(M_trans)) {
549  case NOTRANS: t_transp = Teuchos::NO_TRANS; break;
550  case TRANS: t_transp = Teuchos::TRANS; break;
551  default: TEUCHOS_TEST_FOR_EXCEPT(true);
552  }
553  }
554  if (M_local.numSubCols() > 0) {
555  // AGS: Added std::max on ld? below, following what is done in
556  // Epetra_MultiVector Multiply use of GEMM. Allows for 0 length.
557  blas_.GEMM(
558  t_transp // TRANSA
559  ,Teuchos::NO_TRANS // TRANSB
560  ,Y_local.subDim() // M
561  ,Y_local.numSubCols() // N
562  ,real_trans(M_trans)==NOTRANS ? M_local.numSubCols() : M_local.subDim() // K
563  ,alpha // ALPHA
564  ,const_cast<Scalar*>(M_local.values().getRawPtr()) // A
565  ,std::max((int) M_local.leadingDim(),1) // LDA
566  ,const_cast<Scalar*>(X_local.values().getRawPtr()) // B
567  ,std::max((int) X_local.leadingDim(),1) // LDB
568  ,localBeta // BETA
569  ,Y_local_tmp.values().getRawPtr() // C
570  ,std::max((int) Y_local_tmp.leadingDim(),1) // LDC
571  );
572  }
573  else {
574  std::fill( Y_local_tmp.values().begin(), Y_local_tmp.values().end(),
575  ST::zero() );
576  }
577 #ifdef THYRA_SPMD_MULTI_VECTOR_BASE_PRINT_TIMES
578  timer.stop();
579  std::cout
580  << "\nSpmdMultiVectorDefaultBase<Scalar>::apply(...): Time for GEMM = "
581  << timer.totalElapsedTime() << " seconds\n";
582 #endif
583 
584  if (nonnull(comm)) {
585 
586  //
587  // Perform the global reduction of Y_local_tmp back into Y_local
588  //
589 
590  if (isNonLocalAdjoint) {
591  // Contiguous buffer for final reduction
592  Workspace<Scalar> Y_local_final_buff(wss,Y_local.subDim()*Y_local.numSubCols(),false);
593  // Perform the reduction
594  Teuchos::reduceAll<Ordinal,Scalar>(
595  *comm, Teuchos::REDUCE_SUM, Y_local_final_buff.size(), Y_local_tmp.values().getRawPtr(),
596  Y_local_final_buff.getRawPtr()
597  );
598  // Load Y_local_final_buff back into Y_local
599  typedef typename ArrayView<const Scalar>::const_iterator Y_local_final_buff_iter_t;
600  const ArrayView<const Scalar> Y_local_final_buff_av = Y_local_final_buff;
601  Y_local_final_buff_iter_t Y_local_final_buff_ptr = Y_local_final_buff_av.begin();
602  for( int j = 0; j < Y_local.numSubCols(); ++j ) {
603  typedef typename ArrayRCP<Scalar>::iterator Y_local_values_iter_t;
604  Y_local_values_iter_t Y_local_ptr =
605  Y_local.values().begin() + Y_local.leadingDim()*j;
606  for( int i = 0; i < Y_local.subDim(); ++i ) {
607  (*Y_local_ptr++) = (*Y_local_final_buff_ptr++);
608  }
609  }
610  }
611  }
612  else {
613 
614  // When you get here the view Y_local will be committed back to Y
615  // in the destructor to Y_local
616 
617  }
618 
619 #ifdef THYRA_SPMD_MULTI_VECTOR_BASE_PRINT_TIMES
620  timer.stop();
621  std::cout
622  << "\nSpmdMultiVectorDefaultBase<Scalar>::apply(...): Total time = "
623  << timerTotal.totalElapsedTime() << " seconds\n";
624 #endif
625 
626 }
627 
628 
629 // Protected functions for subclasses to call
630 
631 
632 template<class Scalar>
634 {
635  if(globalDim_ == 0) {
636  const SpmdVectorSpaceBase<Scalar> *l_spmdSpace = this->spmdSpace().get();
637  if(l_spmdSpace) {
638  globalDim_ = l_spmdSpace->dim();
639  localOffset_ = l_spmdSpace->localOffset();
640  localSubDim_ = l_spmdSpace->localSubDim();
641  numCols_ = this->domain()->dim();
642  }
643  else {
644  globalDim_ = 0;
645  localOffset_ = -1;
646  localSubDim_ = 0;
647  numCols_ = 0;
648  }
649  }
650 }
651 
652 
653 template<class Scalar>
655 {
656  const Range1D rowRng = Teuchos::full_range(rowRng_in,0,globalDim_-1);
657 #ifdef TEUCHOS_DEBUG
659  !( 0 <= rowRng.lbound() && rowRng.ubound() < globalDim_ ), std::invalid_argument
660  ,"SpmdMultiVectorDefaultBase<Scalar>::validateRowRange(rowRng): Error, the range rowRng = ["
661  <<rowRng.lbound()<<","<<rowRng.ubound()<<"] is not "
662  "in the range [0,"<<(globalDim_-1)<<"]!"
663  );
664 #endif
665  return rowRng;
666 }
667 
668 
669 template<class Scalar>
671 {
672  const Range1D colRng = Teuchos::full_range(colRng_in,0,numCols_-1);
673 #ifdef TEUCHOS_DEBUG
675  !(0 <= colRng.lbound() && colRng.ubound() < numCols_), std::invalid_argument
676  ,"SpmdMultiVectorDefaultBase<Scalar>::validateColRange(colRng): Error, the range colRng = ["
677  <<colRng.lbound()<<","<<colRng.ubound()<<"] is not "
678  "in the range [0,"<<(numCols_-1)<<"]!"
679  );
680 #endif
681  return colRng;
682 }
683 
684 
685 } // end namespace Thyra
686 
687 
688 #endif // THYRA_SPMD_MULTI_VECTOR_DEFAULT_BASE_DEF_HPP
bool is_null(const boost::shared_ptr< T > &p)
virtual RCP< const VectorSpaceBase< Scalar > > range() const =0
Return a smart pointer for the range space for this operator.
EOpTransp
Enumeration for determining how a linear operator is applied. `*.
virtual void releaseDetachedMultiVectorViewImpl(RTOpPack::ConstSubMultiVectorView< Scalar > *sub_mv) const
Thrown if vector spaces are incompatible.
void apply_op_validate_input(const std::string &func_name, const VectorSpaceBase< Scalar > &space, const RTOpPack::RTOpT< Scalar > &op, const ArrayView< const Ptr< const VectorBase< Scalar > > > &vecs, const ArrayView< const Ptr< VectorBase< Scalar > > > &targ_vecs, const Ptr< RTOpPack::ReductTarget > &reduct_obj, const Ordinal global_offset)
Validate the inputs to VectorBase::applyOp().
RTOpPack::SubMultiVectorView< Scalar > getNonconstLocalSubMultiVectorImpl()
void euclideanApply(const EOpTransp M_trans, const MultiVectorBase< Scalar > &X, const Ptr< MultiVectorBase< Scalar > > &Y, const Scalar alpha, const Scalar beta) const
Uses GEMM() and Teuchos::reduceAll() to implement.
iterator begin() const
Base interface class for SPMD multi-vectors.
void mvMultiReductApplyOpImpl(const RTOpPack::RTOpT< Scalar > &primary_op, const ArrayView< const Ptr< const MultiVectorBase< Scalar > > > &multi_vecs, const ArrayView< const Ptr< MultiVectorBase< Scalar > > > &targ_multi_vecs, const ArrayView< const Ptr< RTOpPack::ReductTarget > > &reduct_objs, const Ordinal primary_global_offset) const
virtual void acquireNonconstDetachedMultiVectorViewImpl(const Range1D &rowRng, const Range1D &colRng, RTOpPack::SubMultiVectorView< Scalar > *sub_mv)
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
Use the non-transposed operator.
const_pointer const_iterator
T_To & dyn_cast(T_From &from)
EOpTransp real_trans(EOpTransp transp)
Return NOTRANS or TRANS for real scalar valued operators and this also is used for determining struct...
Use the transposed operator with complex-conjugate clements (same as TRANS for real scalar types)...
const ArrayRCP< Scalar > values() const
T * get() const
Abstract interface for objects that represent a space for vectors.
Ordinal ubound() const
RTOpPack::ConstSubMultiVectorView< Scalar > getLocalSubMultiVectorImpl() const
Use the transposed operator.
RCP< const ScalarProdVectorSpaceBase< Scalar > > rangeScalarProdVecSpc() const
Returns spmdSpace.
void start(bool reset=false)
Teuchos::Ordinal Ordinal
Type for the dimension of a vector space. `*.
Interface for a collection of column vectors called a multi-vector.
double stop()
void initialize(Ordinal globalOffset_in, Ordinal subDim_in, Ordinal colOffset_in, Ordinal numSubCols_in, const ArrayRCP< Scalar > &values_in, Ordinal leadingDim_in)
virtual void commitNonconstDetachedMultiVectorViewImpl(RTOpPack::SubMultiVectorView< Scalar > *sub_mv)
void acquireDetachedMultiVectorViewImpl(const Range1D &rowRng, const Range1D &colRng, RTOpPack::ConstSubMultiVectorView< Scalar > *sub_mv) const
void commitNonconstDetachedMultiVectorViewImpl(RTOpPack::SubMultiVectorView< Scalar > *sub_mv)
Ordinal lbound() const
bool nonnull(const boost::shared_ptr< T > &p)
virtual Teuchos::RCP< const Teuchos::Comm< Ordinal > > getComm() const =0
Returns the SPMD communicator.
virtual Ordinal localOffset() const =0
Returns the offset for the local sub-vector stored on this process.
void releaseDetachedMultiVectorViewImpl(RTOpPack::ConstSubMultiVectorView< Scalar > *sub_mv) const
Ordinal size() const
double totalElapsedTime(bool readCurrentTime=false) const
void initialize(Ordinal globalOffset_in, Ordinal subDim_in, Ordinal colOffset_in, Ordinal numSubCols_in, const ArrayRCP< const Scalar > &values_in, Ordinal leadingDim_in)
Range1D validateColRange(const Range1D &rowCol) const
Validate and resize the column range.
virtual void acquireDetachedMultiVectorViewImpl(const Range1D &rowRng, const Range1D &colRng, RTOpPack::ConstSubMultiVectorView< Scalar > *sub_mv) const
void acquireNonconstDetachedMultiVectorViewImpl(const Range1D &rowRng, const Range1D &colRng, RTOpPack::SubMultiVectorView< Scalar > *sub_mv)
virtual bool isCompatible(const VectorSpaceBase< Scalar > &vecSpc) const =0
Compare the compatibility of two vector spaces.
Base abstract VectorSpaceBase class for all SPMD-based vector spaces.
virtual Ordinal dim() const =0
Return the dimension of the vector space.
virtual void updateSpmdSpace()
Subclasses should call whenever the structure of the VectorSpaceBase changes.
#define TEUCHOS_TEST_FOR_EXCEPT(throw_exception_test)
iterator begin() const
TEUCHOSCORE_LIB_DLL_EXPORT Teuchos::RCP< WorkspaceStore > get_default_workspace_store()
virtual bool isLocallyReplicated() const =0
Returns true if vector space is locally replicated space.
Range1D validateRowRange(const Range1D &rowRng) const
Validate and resize the row range.
virtual Ordinal localSubDim() const =0
Returns the number of local elements stored on this process.