43 #include "Thyra_DefaultSpmdVectorSpace.hpp"
44 #include "Thyra_DefaultSpmdMultiVector.hpp"
45 #include "Thyra_DefaultSpmdVector.hpp"
46 #include "Thyra_DetachedVectorView.hpp"
47 #include "Thyra_DetachedMultiVectorView.hpp"
48 #include "Thyra_VectorSpaceFactoryBase.hpp"
49 #include "Thyra_ProductVectorSpaceBase.hpp"
50 #include "Teuchos_Assert.hpp"
51 #include "Teuchos_dyn_cast.hpp"
53 #include "Teuchos_DefaultSerialComm.hpp"
55 #include "Teuchos_DefaultMpiComm.hpp"
58 #include "Epetra_Map.h"
59 #include "Epetra_Comm.h"
60 #include "Epetra_SerialComm.h"
62 # include "Epetra_MpiComm.h"
64 #include "Epetra_MultiVector.h"
65 #include "Epetra_Vector.h"
76 unwrapSingleProductVectorSpace(
81 using Teuchos::rcp_dynamic_cast;
83 const RCP<const ProductVectorSpaceBase<double> > pvs =
84 rcp_dynamic_cast<
const ProductVectorSpaceBase<double> >(vs_in);
87 return pvs->getBlock(0);
105 using Teuchos::rcp_dynamic_cast;
106 using Teuchos::set_extra_data;
110 if( serialEpetraComm.
get() ) {
113 set_extra_data( serialEpetraComm,
"serialEpetraComm", Teuchos::inOutArg(serialComm) );
120 mpiEpetraComm = rcp_dynamic_cast<
const Epetra_MpiComm>(epetraComm);
121 if( mpiEpetraComm.
get() ) {
123 rawMpiComm = Teuchos::opaqueWrapper(mpiEpetraComm->
Comm());
124 set_extra_data( mpiEpetraComm,
"mpiEpetraComm", Teuchos::inOutArg(rawMpiComm) );
126 mpiComm =
rcp(
new Teuchos::MpiComm<Ordinal>(rawMpiComm));
133 return Teuchos::null;
147 !epetra_map.
get(), std::invalid_argument,
148 "create_VectorSpace::initialize(...): Error!" );
149 #endif // TEUCHOS_DEBUG
154 Teuchos::set_extra_data(epetra_map,
"epetra_map", inoutArg(comm));
158 defaultSpmdVectorSpace<double>(
159 comm, localSubDim, epetra_map->NumGlobalElements64(),
170 Teuchos::set_extra_data( epetra_map,
"epetra_map", inoutArg(vs) );
182 using Teuchos::rcp_dynamic_cast;
188 return parentSpace->smallVecSpcFcty()->createVecSpc(dim);
209 return Teuchos::null;
218 Teuchos::arcp(localValues,0,epetra_v->Map().NumMyElements(),
false),
222 Teuchos::set_extra_data<RCP<Epetra_Vector> >( epetra_v,
"Epetra_Vector",
223 Teuchos::inOutArg(v) );
245 return Teuchos::null;
254 Teuchos::arcp(localValues,0,epetra_v->Map().NumMyElements(),
false),
258 Teuchos::set_extra_data<RCP<const Epetra_Vector> >( epetra_v,
"Epetra_Vector",
259 Teuchos::inOutArg(v) );
271 using Teuchos::rcp_dynamic_cast;
277 unwrapSingleProductVectorSpace(range_in),
282 unwrapSingleProductVectorSpace(domain_in),
289 if (!epetra_mv.
get() )
290 return Teuchos::null;
297 double *localValues;
int leadingDim;
298 if( epetra_mv->ConstantStride() ) {
299 epetra_mv->ExtractView( &localValues, &leadingDim );
310 Teuchos::arcp(localValues,0,leadingDim*epetra_mv->NumVectors(),
false),
314 Teuchos::set_extra_data<RCP<Epetra_MultiVector> >(
315 epetra_mv,
"Epetra_MultiVector", Teuchos::inOutArg(mv) );
327 using Teuchos::rcp_dynamic_cast;
333 unwrapSingleProductVectorSpace(range_in),
338 unwrapSingleProductVectorSpace(domain_in),
345 if (!epetra_mv.
get())
346 return Teuchos::null;
353 double *localValues;
int leadingDim;
354 if( epetra_mv->ConstantStride() ) {
355 epetra_mv->ExtractView( &localValues, &leadingDim );
366 Teuchos::arcp(localValues,0,leadingDim*epetra_mv->NumVectors(),
false),
370 Teuchos::set_extra_data<RCP<const Epetra_MultiVector> >(
371 epetra_mv,
"Epetra_MultiVector", Teuchos::inOutArg(mv) );
381 using Teuchos::ptrFromRef;
382 using Teuchos::ptr_dynamic_cast;
385 using Teuchos::MpiComm;
391 ptr_dynamic_cast<
const SerialComm<Ordinal> >(comm);
398 ptr_dynamic_cast<
const MpiComm<Ordinal> >(comm);
402 "SPMD std::vector space has a communicator that is "
403 "neither a serial comm nor an MPI comm");
414 TEUCHOS_TEST_FOR_EXCEPTION(
is_null(serialComm), std::runtime_error,
415 "SPMD std::vector space has a communicator that is "
416 "neither a serial comm nor an MPI comm");
422 TEUCHOS_TEST_FOR_EXCEPTION(
is_null(epetraComm), std::runtime_error,
423 "null communicator created");
436 using Teuchos::rcpFromRef;
437 using Teuchos::rcpFromPtr;
438 using Teuchos::rcp_dynamic_cast;
439 using Teuchos::ptrFromRef;
440 using Teuchos::ptr_dynamic_cast;
451 "Error, the concrete VectorSpaceBase object of type "
453 " SpmdVectorSpaceBase or the ProductVectorSpaceBase interfaces!" );
455 const int numBlocks = (
nonnull(prod_vs) ? prod_vs->numBlocks() : 1);
461 for (
int block_i = 0; block_i < numBlocks; ++block_i) {
464 prod_vs->getBlock(block_i),
true);
469 spmd_vs_blocks.
push_back(rcpFromPtr(spmd_vs));
474 int myLocalElements = 0;
475 for (
int block_i = 0; block_i < numBlocks; ++block_i) {
476 myLocalElements += spmd_vs_blocks[block_i]->localSubDim();
483 #ifndef EPETRA_NO_32BIT_GLOBAL_INDICES
488 for (
int block_i = 0; block_i < numBlocks; ++block_i) {
490 const int lowGIDInBlock = spmd_vs_i->localOffset();
491 const int numLocalElementsInBlock = spmd_vs_i->localSubDim();
492 for (
int i=0; i < numLocalElementsInBlock; ++i, ++count) {
493 myGIDs[count] = blockOffset + lowGIDInBlock + i;
495 blockOffset += spmd_vs_i->dim();
498 #ifndef EPETRA_NO_32BIT_GLOBAL_INDICES
499 const int globalDim = vs_in.
dim();
501 const long long globalDim = vs_in.
dim();
513 const RCP<
const VectorSpaceBase<double>>& vs,
514 const RCP<const Epetra_Comm>& comm)
520 const Ptr<const RCP<const Epetra_Map> >
521 epetra_map_ptr = Teuchos::get_optional_extra_data<RCP<const Epetra_Map> >(
526 return *epetra_map_ptr;
531 "Error! No RCP Epetra_Map attached to the input vector space RCP, "
532 "and the input comm RCP is null.\n");
543 using Teuchos::get_optional_extra_data;
548 "Thyra::get_Epetra_Vector(map,v)", *epetra_vs, *v->space() );
555 epetra_v_ptr = get_optional_extra_data<RCP<Epetra_Vector> >(
560 return *epetra_v_ptr;
590 ,
Range1D(localOffset,localOffset+localSubDim-1)
601 ,const_cast<double*>(emvv->values())
609 Teuchos::set_extra_data( emvv,
"emvv", Teuchos::inOutArg(epetra_v),
610 Teuchos::PRE_DESTROY );
619 const RCP<VectorBase<double> > &v,
620 const RCP<const Epetra_Map>& map
627 const Ptr<const RCP<Epetra_Vector> >
628 epetra_v_ptr = Teuchos::get_optional_extra_data<RCP<Epetra_Vector> >(
633 return *epetra_v_ptr;
638 "Error! No RCP Epetra_Vector attached to the input vector RCP, "
639 "and the input map RCP is null.\n");
650 using Teuchos::get_optional_extra_data;
655 "Thyra::get_Epetra_Vector(map,v)", *epetra_vs, *v->space() );
662 epetra_v_ptr = get_optional_extra_data<RCP<const Epetra_Vector> >(
667 return *epetra_v_ptr;
669 epetra_nonconst_v_ptr = get_optional_extra_data<RCP<Epetra_Vector> >(
673 if(!
is_null(epetra_nonconst_v_ptr))
674 return *epetra_nonconst_v_ptr;
689 ,
Range1D(localOffset,localOffset+localSubDim-1)
700 ,const_cast<double*>(evv->values())
708 Teuchos::set_extra_data( evv,
"evv", Teuchos::inOutArg(epetra_v),
709 Teuchos::PRE_DESTROY );
718 const RCP<
const VectorBase<double> > &v,
719 const RCP<const Epetra_Map>& map
726 const Ptr<const RCP<const Epetra_Vector> >
727 epetra_v_ptr = Teuchos::get_optional_extra_data<RCP<const Epetra_Vector> >(
732 return *epetra_v_ptr;
737 "Error! No RCP to Epetra_Vector attached to the input vector RCP, "
738 "and the input map RCP is null.\n");
749 using Teuchos::get_optional_extra_data;
754 "Thyra::get_Epetra_MultiVector(map,mv)", *epetra_vs, *mv->range() );
761 epetra_mv_ptr = get_optional_extra_data<RCP<Epetra_MultiVector> >(
762 mv,
"Epetra_MultiVector");
766 return *epetra_mv_ptr;
793 ,
Range1D(localOffset,localOffset+localSubDim-1)
803 ,const_cast<double*>(emmvv->values())
814 Teuchos::set_extra_data( emmvv,
"emmvv", Teuchos::inOutArg(epetra_mv),
815 Teuchos::PRE_DESTROY );
817 Teuchos::set_extra_data( mv,
"mv", Teuchos::inOutArg(epetra_mv) );
826 const RCP<MultiVectorBase<double> > &mv,
827 const RCP<const Epetra_Map>& map
834 const Ptr<const RCP<Epetra_MultiVector> >
835 epetra_mv_ptr = Teuchos::get_optional_extra_data<RCP<Epetra_MultiVector> >(
836 mv,
"Epetra_MultiVector");
840 return *epetra_mv_ptr;
845 "Error! No RCP to Epetra_MultiVector attached to the input vector RCP, "
846 "and the input map RCP is null.\n");
857 using Teuchos::get_optional_extra_data;
864 "Thyra::get_Epetra_MultiVector(map,mv)",
865 *epetra_vs, *mv->range() );
874 = get_optional_extra_data<RCP<const Epetra_MultiVector> >(
875 mv,
"Epetra_MultiVector" );
879 return *epetra_mv_ptr;
896 ,
Range1D(localOffset,localOffset+localSubDim-1)
907 ,const_cast<double*>(emvv->values())
918 Teuchos::set_extra_data( emvv,
"emvv", Teuchos::inOutArg(epetra_mv),
919 Teuchos::PRE_DESTROY );
921 Teuchos::set_extra_data( mv,
"mv", Teuchos::inOutArg(epetra_mv) );
931 const RCP<
const MultiVectorBase<double> > &mv,
932 const RCP<const Epetra_Map>& map
939 const Ptr<const RCP<const Epetra_MultiVector> >
940 epetra_mv_ptr = Teuchos::get_optional_extra_data<RCP<const Epetra_MultiVector> >(
941 mv,
"Epetra_MultiVector");
945 return *epetra_mv_ptr;
950 "Error! No RCP to Epetra_MultiVector attached to the input vector RCP, "
951 "and the input map RCP is null.\n");
962 using Teuchos::rcpWithEmbeddedObj;
963 using Teuchos::ptrFromRef;
964 using Teuchos::ptr_dynamic_cast;
965 using Teuchos::outArg;
973 mvSpmdMv->getNonconstLocalData(outArg(mvData), outArg(mvLeadingDim));
976 return rcpWithEmbeddedObj(
978 ::
View, map, mvData.getRawPtr(), mvLeadingDim, mv.
domain()->dim()
991 using Teuchos::rcpWithEmbeddedObj;
992 using Teuchos::ptrFromRef;
993 using Teuchos::ptr_dynamic_cast;
994 using Teuchos::outArg;
1002 mvSpmdMv->getLocalData(outArg(mvData), outArg(mvLeadingDim));
1005 return rcpWithEmbeddedObj(
1007 ::
View, map, const_cast<double*>(mvData.getRawPtr()), mvLeadingDim, mv.
domain()->dim()
RCP< Epetra_MultiVector > get_Epetra_MultiVector(const Epetra_Map &map, const RCP< MultiVectorBase< double > > &mv)
Get a non-const Epetra_MultiVector view from a non-const MultiVectorBase object if possible...
bool DistributedGlobal() const
#define THYRA_ASSERT_VEC_SPACES(FUNC_NAME, VS1, VS2)
This is a very useful macro that should be used to validate that two vector spaces are compatible...
Create an explicit non-mutable (const) view of a MultiVectorBase object.
bool is_null(const boost::shared_ptr< T > &p)
Base interface class for SPMD multi-vectors.
RCP< const Epetra_Map > get_Epetra_Map(const VectorSpaceBase< double > &vs, const RCP< const Epetra_Comm > &comm)
Get (or create) an Epetra_Map object given an VectorSpaceBase object an optionally an extra Epetra_Co...
RCP< const VectorSpaceBase< double > > create_VectorSpace(const RCP< const Epetra_Map > &epetra_map)
Create an VectorSpaceBase object given an Epetra_Map object.
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
T_To & dyn_cast(T_From &from)
Create an explicit non-mutable (const) view of a VectorBase object.
TEUCHOSCORE_LIB_DLL_EXPORT std::string demangleName(const std::string &mangledName)
Abstract interface for objects that represent a space for vectors.
RCP< MultiVectorBase< double > > create_MultiVector(const RCP< Epetra_MultiVector > &epetra_mv, const RCP< const VectorSpaceBase< double > > &range=Teuchos::null, const RCP< const VectorSpaceBase< double > > &domain=Teuchos::null)
Create a non-const MultiVectorBase object from a non-const Epetra_MultiVector object.
Create an explicit mutable (non-const) view of a MultiVectorBase object.
int NumMyElements() const
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
Teuchos::Ordinal Ordinal
Type for the dimension of a vector space. `*.
Interface for a collection of column vectors called a multi-vector.
RCP< VectorBase< double > > create_Vector(const RCP< Epetra_Vector > &epetra_v, const RCP< const VectorSpaceBase< double > > &space=Teuchos::null)
Create a non-const VectorBase object from a non-const Epetra_Vector object.
Create an explicit mutable (non-const) view of a VectorBase object.
RCP< Epetra_Vector > get_Epetra_Vector(const Epetra_Map &map, const RCP< VectorBase< double > > &v)
Get a non-const Epetra_Vector view from a non-const VectorBase object if possible.
Abstract interface for finite-dimensional dense vectors.
const Epetra_Comm & Comm() const
const RCP< T > & assert_not_null() const
Efficient concrete implementation subclass for SPMD vectors.
RCP< const Epetra_Comm > get_Epetra_Comm(const Teuchos::Comm< Ordinal > &comm)
Get (or create) and Epetra_Comm given a Teuchos::Comm object.
virtual RCP< const VectorSpaceBase< Scalar > > domain() const =0
Return a smart pointer for the domain space for this operator.
bool nonnull(const boost::shared_ptr< T > &p)
void push_back(const value_type &x)
virtual Ordinal localOffset() const =0
Returns the offset for the local sub-vector stored on this process.
TypeTo as(const TypeFrom &t)
Efficient concrete implementation subclass for SPMD multi-vectors.
#define TEUCHOS_ASSERT_EQUALITY(val1, val2)
int ExtractView(double **V) const
RCP< const VectorSpaceBase< double > > create_LocallyReplicatedVectorSpace(const RCP< const VectorSpaceBase< double > > &parentSpace, const int dim)
Create a VectorSpaceBase object that creates locally replicated vector objects.
RCP< const Teuchos::Comm< Ordinal > > create_Comm(const RCP< const Epetra_Comm > &epetraComm)
Given an Epetra_Comm object, return an equivalent Teuchos::Comm object.
Base abstract VectorSpaceBase class for all SPMD-based vector spaces.
virtual Ordinal dim() const =0
Return the dimension of the vector space.
#define TEUCHOS_TEST_FOR_EXCEPT(throw_exception_test)
virtual Ordinal localSubDim() const =0
Returns the number of local elements stored on this process.