Thyra
Version of the Day
|
Base VectorSpaceBase
class for all SPMD vector spaces with contiguous local-to-global indexing.
More...
#include <Thyra_SpmdVectorSpaceDefaultBase_decl.hpp>
Public Member Functions | |
SpmdVectorSpaceDefaultBase () | |
Public Member Functions inherited from Thyra::SpmdVectorSpaceBase< Scalar > | |
virtual Teuchos::RCP< const Teuchos::Comm< Ordinal > > | getComm () const =0 |
Returns the SPMD communicator. More... | |
virtual Ordinal | localSubDim () const =0 |
Returns the number of local elements stored on this process. More... | |
Public Member Functions inherited from Thyra::VectorSpaceBase< Scalar > | |
void | scalarProds (const MultiVectorBase< Scalar > &X, const MultiVectorBase< Scalar > &Y, const ArrayView< Scalar > &scalarProds_out) const |
Return the scalar product of each column in two multi-vectors in the vector space. More... | |
virtual bool | hasInCoreView (const Range1D &rng=Range1D(), const EViewType viewType=VIEW_TYPE_DETACHED, const EStrideType strideType=STRIDE_TYPE_NONUNIT) const |
Returns true if this->acquireDetachedView(rng,...) returns a direct view of the range of data requested. More... | |
virtual RCP< const VectorSpaceBase< Scalar > > | clone () const |
Clone this object (if supported). More... | |
virtual RCP< MultiVectorBase < Scalar > > | createCachedMembersView (const RTOpPack::SubMultiVectorView< Scalar > &raw_mv, bool initialize=true) const |
Create a (possibly) cached multi-vector member that is a non-const view of raw multi-vector data. The caching mechanism must be implemented by child classes, by default this just calls the regular createMembersView . More... | |
virtual RCP< const MultiVectorBase< Scalar > > | createCachedMembersView (const RTOpPack::ConstSubMultiVectorView< Scalar > &raw_mv) const |
Create a (possibly) cached multi-vector member that is a const view of raw multi-vector data. The caching mechanism must be implemented by child classes, by default this just calls the regular createMembersView . More... | |
Public Member Functions inherited from Thyra::ScalarProdVectorSpaceBase< Scalar > | |
ScalarProdVectorSpaceBase () | |
Construct to use dot product as the default. More... | |
ScalarProdVectorSpaceBase (const RCP< const ScalarProdBase< Scalar > > &scalarProd) | |
Construct with a different scalar product. More... | |
virtual void | setScalarProd (const RCP< const ScalarProdBase< Scalar > > &scalarProd) |
Set a different scalar product. More... | |
RCP< const ScalarProdBase < Scalar > > | getScalarProd () const |
Return the current scalar product. More... | |
bool | isEuclidean () const |
Returns getScalarProd()->isEuclidean() More... | |
Scalar | scalarProd (const VectorBase< Scalar > &x, const VectorBase< Scalar > &y) const |
Returns getScalarProd()->scalarProd(x,y) More... | |
void | scalarProdsImpl (const MultiVectorBase< Scalar > &X, const MultiVectorBase< Scalar > &Y, const ArrayView< Scalar > &scalarProds_out) const |
Calls getScalarProd()->scalarProds(X,Y,scalar_prods) More... | |
Protected Member Functions | |
virtual void | updateState (const Ordinal globalDim, const bool isLocallyReplicated=false) |
This function must be called whenever the state of this changes and some internal state must be updated. More... | |
Protected Member Functions inherited from Thyra::VectorSpaceBase< Scalar > | |
virtual RCP< VectorBase< Scalar > > | createMember () const =0 |
Create a vector member from the vector space. More... | |
Protected Member Functions inherited from Thyra::VectorSpaceDefaultBase< Scalar > | |
RCP< MultiVectorBase< Scalar > > | createMembers (int numMembers) const |
RCP< VectorBase< Scalar > > | createMemberView (const RTOpPack::SubVectorView< Scalar > &raw_v) const |
RCP< const VectorBase< Scalar > > | createMemberView (const RTOpPack::ConstSubVectorView< Scalar > &raw_v) const |
RCP< MultiVectorBase< Scalar > > | createMembersView (const RTOpPack::SubMultiVectorView< Scalar > &raw_mv) const |
RCP< const MultiVectorBase < Scalar > > | createMembersView (const RTOpPack::ConstSubMultiVectorView< Scalar > &raw_mv) const |
Overridden from SpmdVectorSpaceBase | |
Ordinal | localOffset () const |
Ordinal | mapCode () const |
bool | isLocallyReplicated () const |
Returns true if vector space is locally replicated space. More... | |
Overridden form Teuchos::Describable | |
std::string | description () const |
Overridden from VectorSpaceBase | |
Ordinal | dim () const |
Returns the sum of the local number of elements on every process. More... | |
Teuchos::RCP< const VectorSpaceFactoryBase< Scalar > > | smallVecSpcFcty () const |
Returns a DefaultSpmdVectorSpaceFactory object that has been given getComm() . More... | |
bool | isCompatible (const VectorSpaceBase< Scalar > &vecSpc) const |
Checks the general compatibility of parallel (or serial on one process) Spmd-based vector spaces. More... | |
Additional Inherited Members | |
Related Functions inherited from Thyra::VectorSpaceBase< Scalar > | |
template<class Scalar > | |
RCP< const VectorSpaceBase < Scalar > > | makeHaveOwnership (const RCP< const VectorSpaceBase< Scalar > > &vs) |
Helper function that clones a VectorSpaceBase object if the RCP does not have ownership. More... | |
template<class Scalar > | |
RCP< VectorBase< Scalar > > | createMember (const RCP< const VectorSpaceBase< Scalar > > &vs, const std::string &label="") |
Create a vector member from the vector space. More... | |
template<class Scalar > | |
RCP< VectorBase< Scalar > > | createMember (const VectorSpaceBase< Scalar > &vs, const std::string &label="") |
Calls createMember(Teuchos::rcp(&vs,false)) . More... | |
template<class Scalar > | |
RCP< MultiVectorBase< Scalar > > | createMembers (const RCP< const VectorSpaceBase< Scalar > > &vs, int numMembers, const std::string &label="") |
Create a set of vector members (a MultiVectorBase ) from the vector space. More... | |
template<class Scalar > | |
RCP< MultiVectorBase< Scalar > > | createMembers (const RCP< const VectorSpaceBase< Scalar > > &vs, const RCP< const VectorSpaceBase< Scalar > > &domain, const std::string &label="") |
Create a set of vector members (a MultiVectorBase ) from the vector space. More... | |
template<class Scalar > | |
RCP< MultiVectorBase< Scalar > > | createMembers (const VectorSpaceBase< Scalar > &vs, int numMembers, const std::string &label="") |
Calls createMembers(Teuchos::rcp(&vs,false),numMembers) . More... | |
template<class Scalar > | |
RCP< VectorBase< Scalar > > | createMemberView (const RCP< const VectorSpaceBase< Scalar > > &vs, const RTOpPack::SubVectorView< Scalar > &raw_v, const std::string &label="") |
Create a vector member that is a non-const view of raw data. More... | |
template<class Scalar > | |
RCP< VectorBase< Scalar > > | createMemberView (const VectorSpaceBase< Scalar > &vs, const RTOpPack::SubVectorView< Scalar > &raw_v, const std::string &label="") |
Calls createMemberView(Teuchos::rcp(&vs,false),raw_v) . More... | |
template<class Scalar > | |
RCP< const VectorBase< Scalar > > | createMemberView (const RCP< const VectorSpaceBase< Scalar > > &vs, const RTOpPack::ConstSubVectorView< Scalar > &raw_v, const std::string &label="") |
Create a vector member that is a const view of raw data. More... | |
template<class Scalar > | |
RCP< const VectorBase< Scalar > > | createMemberView (const VectorSpaceBase< Scalar > &vs, const RTOpPack::ConstSubVectorView< Scalar > &raw_v, const std::string &label="") |
Calls createMemberView(Teuchos::rcp(&vs,false),raw_v) . More... | |
template<class Scalar > | |
RCP< MultiVectorBase< Scalar > > | createMembersView (const RCP< const VectorSpaceBase< Scalar > > &vs, const RTOpPack::SubMultiVectorView< Scalar > &raw_mv, const std::string &label="") |
Create a multi-vector member that is a non-const view of raw data. More... | |
template<class Scalar > | |
RCP< MultiVectorBase< Scalar > > | createMembersView (const VectorSpaceBase< Scalar > &vs, const RTOpPack::SubMultiVectorView< Scalar > &raw_mv, const std::string &label="") |
Calls createMembersView(Teuchos::rcp(&vs,false),raw_mv) . More... | |
template<class Scalar > | |
RCP< const MultiVectorBase < Scalar > > | createMembersView (const RCP< const VectorSpaceBase< Scalar > > &vs, const RTOpPack::ConstSubMultiVectorView< Scalar > &raw_mv, const std::string &label="") |
Create a multi-vector member that is a const view of raw data. More... | |
template<class Scalar > | |
RCP< const MultiVectorBase < Scalar > > | createMembersView (const VectorSpaceBase< Scalar > &vs, const RTOpPack::ConstSubMultiVectorView< Scalar > &raw_mv, const std::string &label="") |
Calls createMembersView(Teuchos::rcp(&vs,false),raw_mv) . More... | |
Related Functions inherited from Thyra::ScalarProdVectorSpaceBase< Scalar > | |
template<class Scalar > | |
RCP< const ScalarProdVectorSpaceBase < Scalar > > | createSmallScalarProdVectorSpaceBase (const RCP< const VectorSpaceBase< Scalar > > &vs, const Ordinal dim) |
Create a small vector space casted to ScalarProdVectorSpaceBase. More... | |
Base VectorSpaceBase
class for all SPMD vector spaces with contiguous local-to-global indexing.
See SpmdVectorSpaceBase
for details on what this class represents in an abstract way.
Notes to subclass developers:
The pure virtual methods mpiComm()
, localSubDim()
and createMember()
are the only functions that must be overridden.
If this
this is in an uninitialized state then localSubDim()
should return 0
.
It should never be necessary to override the virtual functions mapCode()
and isCompatible()
as these functions have very good and very general implementations.
If optimized implementations of multi-vectors can be supported, then the createMembers()
method should also be overridden.
This class defines a very general default implementation for smallVecSpcFcty()
that returns a DefaultSpmdVectorSpaceFactory
object which in turn can be used to create DefaultSpmdVectorSpace
objects. These are the vector space objects that are used by the domain space of multi-vectors created by createMembers()
. The class DefaultSpmdVectorSpace
creates DefaultSpmdVector
and DefaultSpmdMultiVector
objects. This implementation of smallVecSpcFcty()
is very general should be very appropriate for many different concrete implementations.
Note: It is very important that subclasses call the updateState()
function whenever the state of *this
changes in a way that might affect the behavior of any of the public member functions. For example, if a different value of localSubDim()
will be returned the next time it is called by a client, then updateState()
needs to be called by the subclass. External clients should never need to worry about this function and that is why updateState()
is declared protected.
Definition at line 61 of file Thyra_SpmdVectorSpaceDefaultBase_decl.hpp.
Thyra::SpmdVectorSpaceDefaultBase< Scalar >::SpmdVectorSpaceDefaultBase | ( | ) |
Definition at line 24 of file Thyra_SpmdVectorSpaceDefaultBase_def.hpp.
|
virtual |
This method has a default implementation which just assigns this offset based on counting up localSubDim()
on each process and then setting localOffset()
by the rank of the process. For example, if there are 5 elements in process 0 and 4 elements in process rank, then localOffset
on each of these processes will be set as: localOffset=0
on process 0, localOffset=5
on process 1, localOffset=9
on process 2 and so on.
Implements Thyra::SpmdVectorSpaceBase< Scalar >.
Definition at line 34 of file Thyra_SpmdVectorSpaceDefaultBase_def.hpp.
|
virtual |
This method takes the data getComm()
, numProc
(where numProc=this->getComm()->size()
), localOffset()
or localSubDim()
on each process and then uses it to compute a value for mapCode
(using a single global reduction if numProc > 1
) which is returned from this function.
The value returned from this default implementation of this method must not be changed or this approach breaks down. The only reason for overriding this method is for the subclass to be alerted of when this method is called but not what is returned from this method. If a subclass developer does not understand what this means then don't override this method!
The default implementation will always return return > 0
(unless this
is uninitialized) so that if this method is overridden to return return <=
then this is a flag that the underlying vector map does not satisfy the assumptions of this vector space interface and vectors that are in *this
vector space can not collaborate with other Spmd-based vector implementations.
Implements Thyra::SpmdVectorSpaceBase< Scalar >.
Definition at line 41 of file Thyra_SpmdVectorSpaceDefaultBase_def.hpp.
|
virtual |
Returns true if vector space is locally replicated space.
Implements Thyra::SpmdVectorSpaceBase< Scalar >.
Definition at line 48 of file Thyra_SpmdVectorSpaceDefaultBase_def.hpp.
|
virtual |
Reimplemented from Teuchos::Describable.
Definition at line 55 of file Thyra_SpmdVectorSpaceDefaultBase_def.hpp.
|
virtual |
Returns the sum of the local number of elements on every process.
Implements Thyra::VectorSpaceBase< Scalar >.
Definition at line 80 of file Thyra_SpmdVectorSpaceDefaultBase_def.hpp.
|
virtual |
Returns a DefaultSpmdVectorSpaceFactory
object that has been given getComm()
.
Implements Thyra::VectorSpaceBase< Scalar >.
Definition at line 88 of file Thyra_SpmdVectorSpaceDefaultBase_def.hpp.
|
virtual |
Checks the general compatibility of parallel (or serial on one process) Spmd-based vector spaces.
*this
and vecSpace
are both serial in-core vectors or if vecSpc
is of type SpmdVectorSpaceDefaultBase<Scalar>
and both *this
and vecSpc
have the same Spmd communicators and the same mapping of elements to processes.Postconditions:
return = ( this->hasInCoreView() && vecSpc->hasInCoreView() ) || ( (mpiVecSpc = dynamic_cast<const SpmdVectorSpaceDefaultBase<Scalar>*>(&vecSpc)) && this->getComm() == mpiVecSpc->getComm() && this->mapCode() == mpiVecSpc->getComm())
.
If the mapping of vector elements to processes is not as described above then this method should be overridden in a way that is specific to the vector implementation.
Implements Thyra::VectorSpaceBase< Scalar >.
Definition at line 95 of file Thyra_SpmdVectorSpaceDefaultBase_def.hpp.
|
protectedvirtual |
This function must be called whenever the state of this
changes and some internal state must be updated.
globalDim | [in] If globalDim > 0 then this determines the global dimension of the vector space. If globalDim==this->localSubDim() and this->localSubDim() is the same on every process and is > 0 then this is a locally replicated vector space. If globalDim < 0 then the global dimension is computed using a global reduction. If this->getComm()->size()==1 then this argument is ignored. |
isLocallyReplicated | [in] If isLocallyReplicated==true , then tt>globalDim must be equal to this->localSubDim() on every process but no global reduction will be done to verify this. If false , an global reduction is used to automatically determine if globalDim == this->localSubDim() on every process and therefore construct a locally replicated vector space. |
Note that calling this function may involve one or more global reductions being called if this is parallel vector space so it should only be called when needed by subclasses and RCP should be use to share a constructed VS with lots of clients.
Usually, this operation only needs to be called once for every new parallel vector space constructed and very few parallel vector spaces will be created per application usually.
Definition at line 133 of file Thyra_SpmdVectorSpaceDefaultBase_def.hpp.