|
Thyra Package Browser (Single Doxygen Collection) Version of the Day
|
00001 // @HEADER 00002 // *********************************************************************** 00003 // 00004 // Thyra: Interfaces and Support for Abstract Numerical Algorithms 00005 // Copyright (2004) Sandia Corporation 00006 // 00007 // Under terms of Contract DE-AC04-94AL85000, there is a non-exclusive 00008 // license for use of this work by or on behalf of the U.S. Government. 00009 // 00010 // This library is free software; you can redistribute it and/or modify 00011 // it under the terms of the GNU Lesser General Public License as 00012 // published by the Free Software Foundation; either version 2.1 of the 00013 // License, or (at your option) any later version. 00014 // 00015 // This library is distributed in the hope that it will be useful, but 00016 // WITHOUT ANY WARRANTY; without even the implied warranty of 00017 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU 00018 // Lesser General Public License for more details. 00019 // 00020 // You should have received a copy of the GNU Lesser General Public 00021 // License along with this library; if not, write to the Free Software 00022 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 00023 // USA 00024 // Questions? Contact Michael A. Heroux (maherou@sandia.gov) 00025 // 00026 // *********************************************************************** 00027 // @HEADER 00028 00029 #ifndef THYRA_SPMD_MULTI_VECTOR_BASE_DEF_HPP 00030 #define THYRA_SPMD_MULTI_VECTOR_BASE_DEF_HPP 00031 00032 00033 #include "Thyra_SpmdMultiVectorBase_decl.hpp" 00034 #include "Thyra_MultiVectorDefaultBase.hpp" 00035 #include "Thyra_MultiVectorAdapterBase.hpp" 00036 #include "Thyra_SpmdVectorSpaceDefaultBase.hpp" 00037 #include "Thyra_DetachedMultiVectorView.hpp" 00038 #include "Thyra_apply_op_helper.hpp" 00039 #include "RTOpPack_SPMD_apply_op.hpp" 00040 #include "RTOp_parallel_helpers.h" 00041 #include "Teuchos_Workspace.hpp" 00042 #include "Teuchos_dyn_cast.hpp" 00043 #include "Teuchos_Time.hpp" 00044 #include "Teuchos_CommHelpers.hpp" 00045 00046 00047 // Define to see some timing output! 00048 //#define THYRA_SPMD_MULTI_VECTOR_BASE_PRINT_TIMES 00049 00050 00051 namespace Thyra { 00052 00053 00054 template<class Scalar> 00055 SpmdMultiVectorBase<Scalar>::SpmdMultiVectorBase() 00056 :in_applyOp_(false), 00057 globalDim_(0), 00058 localOffset_(-1), 00059 localSubDim_(0), 00060 numCols_(0) 00061 {} 00062 00063 00064 // Overridden public functions from MultiVectorAdapterBase 00065 00066 00067 template<class Scalar> 00068 RCP< const ScalarProdVectorSpaceBase<Scalar> > 00069 SpmdMultiVectorBase<Scalar>::rangeScalarProdVecSpc() const 00070 { 00071 return Teuchos::rcp_dynamic_cast<const ScalarProdVectorSpaceBase<Scalar> >( 00072 spmdSpace(), true 00073 ); 00074 } 00075 00076 00077 // Protected functions overridden from MultiVectorBase 00078 00079 00080 template<class Scalar> 00081 void SpmdMultiVectorBase<Scalar>::mvMultiReductApplyOpImpl( 00082 const RTOpPack::RTOpT<Scalar> &pri_op, 00083 const ArrayView<const Ptr<const MultiVectorBase<Scalar> > > &multi_vecs, 00084 const ArrayView<const Ptr<MultiVectorBase<Scalar> > > &targ_multi_vecs, 00085 const ArrayView<const Ptr<RTOpPack::ReductTarget> > &reduct_objs, 00086 const Ordinal pri_global_offset_in 00087 ) const 00088 { 00089 00090 using Teuchos::dyn_cast; 00091 using Teuchos::Workspace; 00092 00093 Teuchos::WorkspaceStore* wss = Teuchos::get_default_workspace_store().get(); 00094 00095 const Ordinal numCols = this->domain()->dim(); 00096 const SpmdVectorSpaceBase<Scalar> &spmdSpc = *spmdSpace(); 00097 00098 #ifdef TEUCHOS_DEBUG 00099 TEST_FOR_EXCEPTION( 00100 in_applyOp_, std::invalid_argument, 00101 "SpmdMultiVectorBase<>::mvMultiReductApplyOpImpl(...): Error, this method is" 00102 " being entered recursively which is a clear sign that one of the methods" 00103 " acquireDetachedView(...), releaseDetachedView(...) or commitDetachedView(...)" 00104 " was not implemented properly!" 00105 ); 00106 apply_op_validate_input( 00107 "SpmdMultiVectorBase<>::mvMultiReductApplyOpImpl(...)", *this->domain(), 00108 *this->range(), pri_op, multi_vecs, targ_multi_vecs, reduct_objs, 00109 pri_global_offset_in); 00110 #endif 00111 00112 // Flag that we are in applyOp() 00113 in_applyOp_ = true; 00114 00115 // First see if this is a locally replicated vector in which case 00116 // we treat this as a local operation only. 00117 const bool locallyReplicated = (localSubDim_ == globalDim_); 00118 00119 const Range1D local_rng(localOffset_, localOffset_+localSubDim_-1); 00120 const Range1D col_rng(0, numCols-1); 00121 00122 // Create sub-vector views of all of the *participating* local data 00123 Workspace<RTOpPack::ConstSubMultiVectorView<Scalar> > 00124 sub_multi_vecs(wss,multi_vecs.size()); 00125 Workspace<RTOpPack::SubMultiVectorView<Scalar> > 00126 targ_sub_multi_vecs(wss,targ_multi_vecs.size()); 00127 for(int k = 0; k < multi_vecs.size(); ++k ) { 00128 multi_vecs[k]->acquireDetachedView(local_rng, col_rng, &sub_multi_vecs[k]); 00129 sub_multi_vecs[k].setGlobalOffset(localOffset_+pri_global_offset_in); 00130 } 00131 for(int k = 0; k < targ_multi_vecs.size(); ++k ) { 00132 targ_multi_vecs[k]->acquireDetachedView(local_rng, col_rng, &targ_sub_multi_vecs[k]); 00133 targ_sub_multi_vecs[k].setGlobalOffset(localOffset_+pri_global_offset_in); 00134 } 00135 Workspace<RTOpPack::ReductTarget*> reduct_objs_ptr(wss, reduct_objs.size()); 00136 for (int k = 0; k < reduct_objs.size(); ++k) { 00137 reduct_objs_ptr[k] = &*reduct_objs[k]; 00138 } 00139 00140 // Apply the RTOp operator object (all processors must participate) 00141 RTOpPack::SPMD_apply_op( 00142 locallyReplicated ? NULL : spmdSpc.getComm().get(), // comm 00143 pri_op, // op 00144 col_rng.size(), // num_cols 00145 multi_vecs.size(), // multi_vecs.size() 00146 multi_vecs.size() ? &sub_multi_vecs[0] : NULL, // sub_multi_vecs 00147 targ_multi_vecs.size(), // targ_multi_vecs.size() 00148 targ_multi_vecs.size() ? &targ_sub_multi_vecs[0] : NULL, // targ_sub_multi_vecs 00149 reduct_objs.size() ? &reduct_objs_ptr[0] : 0 // reduct_objs 00150 ); 00151 00152 // Free and commit the local data 00153 for(int k = 0; k < multi_vecs.size(); ++k ) { 00154 sub_multi_vecs[k].setGlobalOffset(local_rng.lbound()); 00155 multi_vecs[k]->releaseDetachedView( &sub_multi_vecs[k] ); 00156 } 00157 for(int k = 0; k < targ_multi_vecs.size(); ++k ) { 00158 targ_sub_multi_vecs[k].setGlobalOffset(local_rng.lbound()); 00159 targ_multi_vecs[k]->commitDetachedView( &targ_sub_multi_vecs[k] ); 00160 } 00161 00162 // Flag that we are leaving applyOp() 00163 in_applyOp_ = false; 00164 00165 } 00166 00167 00168 template<class Scalar> 00169 void SpmdMultiVectorBase<Scalar>::acquireDetachedMultiVectorViewImpl( 00170 const Range1D &rowRng_in, 00171 const Range1D &colRng_in, 00172 RTOpPack::ConstSubMultiVectorView<Scalar> *sub_mv 00173 ) const 00174 { 00175 using Teuchos::outArg; 00176 const Range1D rowRng = validateRowRange(rowRng_in); 00177 const Range1D colRng = validateColRange(colRng_in); 00178 if( rowRng.lbound() < localOffset_ || localOffset_+localSubDim_-1 < rowRng.ubound() ) { 00179 // rng consists of off-processor elements so use the default implementation! 00180 MultiVectorDefaultBase<Scalar>::acquireDetachedMultiVectorViewImpl( 00181 rowRng_in,colRng_in,sub_mv 00182 ); 00183 return; 00184 } 00185 ArrayRCP<const Scalar> localValues; 00186 Ordinal leadingDim = 0; 00187 this->getLocalData(outArg(localValues), outArg(leadingDim)); 00188 sub_mv->initialize( 00189 rowRng.lbound(), // globalOffset 00190 rowRng.size(), // subDim 00191 colRng.lbound(), // colOffset 00192 colRng.size(), // numSubCols 00193 localValues 00194 +(rowRng.lbound()-localOffset_) 00195 +colRng.lbound()*leadingDim, // values 00196 leadingDim // leadingDim 00197 ); 00198 } 00199 00200 00201 template<class Scalar> 00202 void SpmdMultiVectorBase<Scalar>::releaseDetachedMultiVectorViewImpl( 00203 RTOpPack::ConstSubMultiVectorView<Scalar>* sub_mv 00204 ) const 00205 { 00206 if( 00207 sub_mv->globalOffset() < localOffset_ 00208 || 00209 localOffset_+localSubDim_ < sub_mv->globalOffset()+sub_mv->subDim() 00210 ) 00211 { 00212 // Let the default implementation handle it! 00213 MultiVectorDefaultBase<Scalar>::releaseDetachedMultiVectorViewImpl(sub_mv); 00214 return; 00215 } 00216 sub_mv->uninitialize(); 00217 } 00218 00219 00220 template<class Scalar> 00221 void SpmdMultiVectorBase<Scalar>::acquireNonconstDetachedMultiVectorViewImpl( 00222 const Range1D &rowRng_in, 00223 const Range1D &colRng_in, 00224 RTOpPack::SubMultiVectorView<Scalar> *sub_mv 00225 ) 00226 { 00227 using Teuchos::outArg; 00228 const Range1D rowRng = validateRowRange(rowRng_in); 00229 const Range1D colRng = validateColRange(colRng_in); 00230 if( 00231 rowRng.lbound() < localOffset_ 00232 || 00233 localOffset_+localSubDim_-1 < rowRng.ubound() 00234 ) 00235 { 00236 // rng consists of off-processor elements so use the default implementation! 00237 MultiVectorDefaultBase<Scalar>::acquireNonconstDetachedMultiVectorViewImpl( 00238 rowRng_in, colRng_in, sub_mv 00239 ); 00240 return; 00241 } 00242 ArrayRCP<Scalar> localValues; 00243 Ordinal leadingDim = 0; 00244 this->getNonconstLocalData(outArg(localValues), outArg(leadingDim)); 00245 sub_mv->initialize( 00246 rowRng.lbound() // globalOffset 00247 ,rowRng.size() // subDim 00248 ,colRng.lbound() // colOffset 00249 ,colRng.size() // numSubCols 00250 ,localValues 00251 +(rowRng.lbound()-localOffset_) 00252 +colRng.lbound()*leadingDim // values 00253 ,leadingDim // leadingDim 00254 ); 00255 } 00256 00257 00258 template<class Scalar> 00259 void SpmdMultiVectorBase<Scalar>::commitNonconstDetachedMultiVectorViewImpl( 00260 RTOpPack::SubMultiVectorView<Scalar>* sub_mv 00261 ) 00262 { 00263 if( 00264 sub_mv->globalOffset() < localOffset_ 00265 || 00266 localOffset_+localSubDim_ < sub_mv->globalOffset()+sub_mv->subDim() 00267 ) 00268 { 00269 // Let the default implementation handle it! 00270 MultiVectorDefaultBase<Scalar>::commitNonconstDetachedMultiVectorViewImpl(sub_mv); 00271 return; 00272 } 00273 sub_mv->uninitialize(); 00274 } 00275 00276 00277 // Protected functions overridden from MultiVectorAdapterBase 00278 00279 00280 template<class Scalar> 00281 void SpmdMultiVectorBase<Scalar>::euclideanApply( 00282 const EOpTransp M_trans, 00283 const MultiVectorBase<Scalar> &X, 00284 const Ptr<MultiVectorBase<Scalar> > &Y, 00285 const Scalar alpha, 00286 const Scalar beta 00287 ) const 00288 { 00289 typedef Teuchos::ScalarTraits<Scalar> ST; 00290 using Teuchos::Workspace; 00291 Teuchos::WorkspaceStore* wss = Teuchos::get_default_workspace_store().get(); 00292 00293 #ifdef THYRA_SPMD_MULTI_VECTOR_BASE_PRINT_TIMES 00294 Teuchos::Time timerTotal("dummy",true); 00295 Teuchos::Time timer("dummy"); 00296 #endif 00297 00298 // 00299 // This function performs one of two operations. 00300 // 00301 // The first operation (M_trans == NOTRANS) is: 00302 // 00303 // Y = beta * Y + alpha * M * X 00304 // 00305 // where Y and M have compatible (distributed?) range vector 00306 // spaces and X is a locally replicated serial multi-vector. This 00307 // operation does not require any global communication. 00308 // 00309 // The second operation (M_trans == TRANS) is: 00310 // 00311 // Y = beta * Y + alpha * M' * X 00312 // 00313 // where M and X have compatible (distributed?) range vector spaces 00314 // and Y is a locally replicated serial multi-vector. This operation 00315 // requires a local reduction. 00316 // 00317 00318 // 00319 // Get spaces and validate compatibility 00320 // 00321 00322 // Get the SpmdVectorSpace 00323 const SpmdVectorSpaceBase<Scalar> &spmdSpc = *this->spmdSpace(); 00324 00325 // Get the Spmd communicator 00326 const RCP<const Teuchos::Comm<Ordinal> > 00327 comm = spmdSpc.getComm(); 00328 #ifdef TEUCHOS_DEBUG 00329 const VectorSpaceBase<Scalar> 00330 &Y_range = *Y->range(), 00331 &X_range = *X.range(); 00332 // std::cout << "SpmdMultiVectorBase<Scalar>::apply(...): comm = " << comm << std::endl; 00333 TEST_FOR_EXCEPTION( 00334 ( globalDim_ > localSubDim_ ) && comm.get()==NULL, std::logic_error 00335 ,"SpmdMultiVectorBase<Scalar>::apply(...MultiVectorBase<Scalar>...): Error!" 00336 ); 00337 // ToDo: Write a good general validation function that I can call that will replace 00338 // all of these TEST_FOR_EXCEPTION(...) uses 00339 00340 TEST_FOR_EXCEPTION( 00341 real_trans(M_trans)==NOTRANS && !spmdSpc.isCompatible(Y_range), Exceptions::IncompatibleVectorSpaces 00342 ,"SpmdMultiVectorBase<Scalar>::apply(...MultiVectorBase<Scalar>...): Error!" 00343 ); 00344 TEST_FOR_EXCEPTION( 00345 real_trans(M_trans)==TRANS && !spmdSpc.isCompatible(X_range), Exceptions::IncompatibleVectorSpaces 00346 ,"SpmdMultiVectorBase<Scalar>::apply(...MultiVectorBase<Scalar>...): Error!" 00347 ); 00348 #endif 00349 00350 // 00351 // Get explicit (local) views of Y, M and X 00352 // 00353 00354 #ifdef THYRA_SPMD_MULTI_VECTOR_BASE_PRINT_TIMES 00355 timer.start(); 00356 #endif 00357 00358 DetachedMultiVectorView<Scalar> 00359 Y_local( 00360 *Y, 00361 real_trans(M_trans)==NOTRANS ? Range1D(localOffset_,localOffset_+localSubDim_-1) : Range1D(), 00362 Range1D() 00363 ); 00364 ConstDetachedMultiVectorView<Scalar> 00365 M_local( 00366 *this, 00367 Range1D(localOffset_,localOffset_+localSubDim_-1), 00368 Range1D() 00369 ); 00370 ConstDetachedMultiVectorView<Scalar> 00371 X_local( 00372 X 00373 ,real_trans(M_trans)==NOTRANS ? Range1D() : Range1D(localOffset_,localOffset_+localSubDim_-1) 00374 ,Range1D() 00375 ); 00376 #ifdef THYRA_SPMD_MULTI_VECTOR_BASE_PRINT_TIMES 00377 timer.stop(); 00378 std::cout << "\nSpmdMultiVectorBase<Scalar>::apply(...): Time for getting view = " << timer.totalElapsedTime() << " seconds\n"; 00379 #endif 00380 #ifdef TEUCHOS_DEBUG 00381 TEST_FOR_EXCEPTION( 00382 real_trans(M_trans)==NOTRANS && ( M_local.numSubCols() != X_local.subDim() || X_local.numSubCols() != Y_local.numSubCols() ) 00383 , Exceptions::IncompatibleVectorSpaces 00384 ,"SpmdMultiVectorBase<Scalar>::apply(...MultiVectorBase<Scalar>...): Error!" 00385 ); 00386 TEST_FOR_EXCEPTION( 00387 real_trans(M_trans)==TRANS && ( M_local.subDim() != X_local.subDim() || X_local.numSubCols() != Y_local.numSubCols() ) 00388 , Exceptions::IncompatibleVectorSpaces 00389 ,"SpmdMultiVectorBase<Scalar>::apply(...MultiVectorBase<Scalar>...): Error!" 00390 ); 00391 #endif 00392 00393 // 00394 // If nonlocal (i.e. M_trans==TRANS) then create temporary storage 00395 // for: 00396 // 00397 // Y_local_tmp = alpha * M(local) * X(local) : on nonroot processes 00398 // 00399 // or 00400 // 00401 // Y_local_tmp = beta*Y_local + alpha * M(local) * X(local) : on root process (localOffset_==0) 00402 // 00403 // and set 00404 // 00405 // localBeta = ( localOffset_ == 0 ? beta : 0.0 ) 00406 // 00407 // Above, we choose localBeta such that we will only perform 00408 // Y_local = beta * Y_local + ... on one process (the root 00409 // process where localOffset_==0x). Then, when we add up Y_local 00410 // on all of the processors and we will get the correct result. 00411 // 00412 // If strictly local (i.e. M_trans == NOTRANS) then set: 00413 // 00414 // Y_local_tmp = Y_local 00415 // localBeta = beta 00416 // 00417 00418 #ifdef THYRA_SPMD_MULTI_VECTOR_BASE_PRINT_TIMES 00419 timer.start(); 00420 #endif 00421 00422 Workspace<Scalar> Y_local_tmp_store(wss, Y_local.subDim()*Y_local.numSubCols(), false); 00423 RTOpPack::SubMultiVectorView<Scalar> Y_local_tmp; 00424 Scalar localBeta; 00425 if( real_trans(M_trans) == TRANS && globalDim_ > localSubDim_ ) { 00426 // Nonlocal 00427 Y_local_tmp.initialize( 00428 0, Y_local.subDim(), 00429 0, Y_local.numSubCols(), 00430 Teuchos::arcpFromArrayView(Y_local_tmp_store()), 00431 Y_local.subDim() // leadingDim == subDim (columns are adjacent) 00432 ); 00433 if( localOffset_ == 0 ) { 00434 // Root process: Must copy Y_local into Y_local_tmp 00435 for( int j = 0; j < Y_local.numSubCols(); ++j ) { 00436 Scalar *Y_local_j = Y_local.values() + Y_local.leadingDim()*j; 00437 std::copy( Y_local_j, Y_local_j + Y_local.subDim(), Y_local_tmp.values() + Y_local_tmp.leadingDim()*j ); 00438 } 00439 localBeta = beta; 00440 } 00441 else { 00442 // Not the root process 00443 localBeta = 0.0; 00444 } 00445 } 00446 else { 00447 // Local 00448 Y_local_tmp = Y_local.smv(); // Shallow copy only! 00449 localBeta = beta; 00450 } 00451 00452 #ifdef THYRA_SPMD_MULTI_VECTOR_BASE_PRINT_TIMES 00453 timer.stop(); 00454 std::cout << "\nSpmdMultiVectorBase<Scalar>::apply(...): Time for setting up Y_local_tmp and localBeta = " << timer.totalElapsedTime() << " seconds\n"; 00455 #endif 00456 00457 // 00458 // Perform the local multiplication: 00459 // 00460 // Y(local) = localBeta * Y(local) + alpha * op(M(local)) * X(local) 00461 // 00462 // or in BLAS lingo: 00463 // 00464 // C = beta * C + alpha * op(A) * op(B) 00465 // 00466 00467 #ifdef THYRA_SPMD_MULTI_VECTOR_BASE_PRINT_TIMES 00468 timer.start(); 00469 #endif 00470 Teuchos::ETransp t_transp; 00471 if(ST::isComplex) { 00472 switch(M_trans) { 00473 case NOTRANS: t_transp = Teuchos::NO_TRANS; break; 00474 case TRANS: t_transp = Teuchos::TRANS; break; 00475 case CONJTRANS: t_transp = Teuchos::CONJ_TRANS; break; 00476 default: TEST_FOR_EXCEPT(true); 00477 } 00478 } 00479 else { 00480 switch(real_trans(M_trans)) { 00481 case NOTRANS: t_transp = Teuchos::NO_TRANS; break; 00482 case TRANS: t_transp = Teuchos::TRANS; break; 00483 default: TEST_FOR_EXCEPT(true); 00484 } 00485 } 00486 if (M_local.numSubCols() > 0) { 00487 blas_.GEMM( 00488 t_transp // TRANSA 00489 ,Teuchos::NO_TRANS // TRANSB 00490 ,Y_local.subDim() // M 00491 ,Y_local.numSubCols() // N 00492 ,real_trans(M_trans)==NOTRANS ? M_local.numSubCols() : M_local.subDim() // K 00493 ,alpha // ALPHA 00494 ,const_cast<Scalar*>(M_local.values()) // A 00495 ,M_local.leadingDim() // LDA 00496 ,const_cast<Scalar*>(X_local.values()) // B 00497 ,X_local.leadingDim() // LDB 00498 ,localBeta // BETA 00499 ,Y_local_tmp.values().get() // C 00500 ,Y_local_tmp.leadingDim() // LDC 00501 ); 00502 } 00503 else { 00504 std::fill( Y_local_tmp.values().begin(), Y_local_tmp.values().end(), 00505 ST::zero() ); 00506 } 00507 #ifdef THYRA_SPMD_MULTI_VECTOR_BASE_PRINT_TIMES 00508 timer.stop(); 00509 std::cout 00510 << "\nSpmdMultiVectorBase<Scalar>::apply(...): Time for GEMM = " 00511 << timer.totalElapsedTime() << " seconds\n"; 00512 #endif 00513 00514 if( comm.get() ) { 00515 00516 // 00517 // Perform the global reduction of Y_local_tmp back into Y_local 00518 // 00519 00520 if( real_trans(M_trans)==TRANS && globalDim_ > localSubDim_ ) { 00521 // Contiguous buffer for final reduction 00522 Workspace<Scalar> Y_local_final_buff(wss,Y_local.subDim()*Y_local.numSubCols(),false); 00523 // Perform the reduction 00524 Teuchos::reduceAll<Ordinal,Scalar>( 00525 *comm,Teuchos::REDUCE_SUM,Y_local_final_buff.size(),Y_local_tmp.values().get(), 00526 &Y_local_final_buff[0] 00527 ); 00528 // Load Y_local_final_buff back into Y_local 00529 const Scalar *Y_local_final_buff_ptr = &Y_local_final_buff[0]; 00530 for( int j = 0; j < Y_local.numSubCols(); ++j ) { 00531 Scalar *Y_local_ptr = Y_local.values() + Y_local.leadingDim()*j; 00532 for( int i = 0; i < Y_local.subDim(); ++i ) { 00533 (*Y_local_ptr++) = (*Y_local_final_buff_ptr++); 00534 } 00535 } 00536 } 00537 } 00538 else { 00539 00540 // When you get here the view Y_local will be committed back to Y 00541 // in the destructor to Y_local 00542 00543 } 00544 00545 #ifdef THYRA_SPMD_MULTI_VECTOR_BASE_PRINT_TIMES 00546 timer.stop(); 00547 std::cout 00548 << "\nSpmdMultiVectorBase<Scalar>::apply(...): Total time = " 00549 << timerTotal.totalElapsedTime() << " seconds\n"; 00550 #endif 00551 00552 } 00553 00554 00555 // Protected functions for subclasses to call 00556 00557 00558 template<class Scalar> 00559 void SpmdMultiVectorBase<Scalar>::updateSpmdSpace() 00560 { 00561 if(globalDim_ == 0) { 00562 const SpmdVectorSpaceBase<Scalar> *l_spmdSpace = this->spmdSpace().get(); 00563 if(l_spmdSpace) { 00564 globalDim_ = l_spmdSpace->dim(); 00565 localOffset_ = l_spmdSpace->localOffset(); 00566 localSubDim_ = l_spmdSpace->localSubDim(); 00567 numCols_ = this->domain()->dim(); 00568 } 00569 else { 00570 globalDim_ = 0; 00571 localOffset_ = -1; 00572 localSubDim_ = 0; 00573 numCols_ = 0; 00574 } 00575 } 00576 } 00577 00578 00579 template<class Scalar> 00580 Range1D SpmdMultiVectorBase<Scalar>::validateRowRange( const Range1D &rowRng_in ) const 00581 { 00582 const Range1D rowRng = Teuchos::full_range(rowRng_in,0,globalDim_-1); 00583 #ifdef TEUCHOS_DEBUG 00584 TEST_FOR_EXCEPTION( 00585 !( 0 <= rowRng.lbound() && rowRng.ubound() < globalDim_ ), std::invalid_argument 00586 ,"SpmdMultiVectorBase<Scalar>::validateRowRange(rowRng): Error, the range rowRng = [" 00587 <<rowRng.lbound()<<","<<rowRng.ubound()<<"] is not " 00588 "in the range [0,"<<(globalDim_-1)<<"]!" 00589 ); 00590 #endif 00591 return rowRng; 00592 } 00593 00594 00595 template<class Scalar> 00596 Range1D SpmdMultiVectorBase<Scalar>::validateColRange( const Range1D &colRng_in ) const 00597 { 00598 const Range1D colRng = Teuchos::full_range(colRng_in,0,numCols_-1); 00599 #ifdef TEUCHOS_DEBUG 00600 TEST_FOR_EXCEPTION( 00601 !(0 <= colRng.lbound() && colRng.ubound() < numCols_), std::invalid_argument 00602 ,"SpmdMultiVectorBase<Scalar>::validateColRange(colRng): Error, the range colRng = [" 00603 <<colRng.lbound()<<","<<colRng.ubound()<<"] is not " 00604 "in the range [0,"<<(numCols_-1)<<"]!" 00605 ); 00606 #endif 00607 return colRng; 00608 } 00609 00610 00611 } // end namespace Thyra 00612 00613 00614 #endif // THYRA_SPMD_MULTI_VECTOR_BASE_DEF_HPP
1.7.4