Advanced Multi-Physics (AMP)
On-Line Documentation
ArrayVectorData.h
Go to the documentation of this file.
1#ifndef included_AMP_ArrayVectorData
2#define included_AMP_ArrayVectorData
3
4#include <string>
5
6#include "AMP/utils/Array.h"
7#include "AMP/utils/FunctionTable.h"
8#include "AMP/utils/Memory.h"
9#include "AMP/utils/UtilityMacros.h"
10#include "AMP/vectors/data/GhostDataHelper.hpp"
11#include "AMP/vectors/data/VectorData.h"
12
13
14namespace AMP::LinearAlgebra {
15
19template<typename T, typename FUN = FunctionTable<T>, typename Allocator = AMP::HostAllocator<void>>
20class ArrayVectorData : public GhostDataHelper<T, Allocator>
21{
22private:
25 size_t d_offset = 0;
28
29public:
35 static std::shared_ptr<ArrayVectorData> create( const ArraySize &localSize );
36
44 static std::shared_ptr<ArrayVectorData>
45 create( const ArraySize &localSize, const ArraySize &blockIndex, AMP_MPI comm );
46
47 static std::shared_ptr<ArrayVectorData>
48 create( const size_t localSize, std::shared_ptr<CommunicationList> commList, T *data );
49
51 ArrayVectorData() = default;
52
54 virtual ~ArrayVectorData() {}
55
57 ArrayVectorData( const ArrayVectorData & ) = delete;
58
61 const AMP_MPI &getComm() const override { return d_comm; }
62
64 void resize( const ArraySize &localDims );
65
68
70 const Array<T, FUN, Allocator> &getArray( void ) const { return d_array; }
71
77 size_t numberOfDataBlocks() const override { return 1; }
78
83 size_t sizeOfDataBlock( size_t i = 0 ) const override { return i == 0 ? d_array.length() : 0; }
84
89 void putRawData( const void *buf, const typeID &id ) override;
90
96 void getRawData( void *buf, const typeID &id ) const override;
97
108 void setValuesByLocalID( size_t num,
109 const size_t *indices,
110 const void *vals,
111 const typeID &id ) override;
112
124 void addValuesByLocalID( size_t num,
125 const size_t *indices,
126 const void *vals,
127 const typeID &id ) override;
128
137 void getValuesByLocalID( size_t num,
138 const size_t *indices,
139 void *vals,
140 const typeID &id ) const override;
141
149 uint64_t getDataID() const override { return reinterpret_cast<uint64_t>( d_array.data() ); }
150
152 std::string VectorDataName() const override { return "ArrayVector"; }
153
154protected:
159 void *getRawDataBlockAsVoid( size_t i ) override
160 {
161 AMP_ASSERT( i == 0 );
162 return d_array.data();
163 }
164
169 const void *getRawDataBlockAsVoid( size_t i ) const override
170 {
171 AMP_ASSERT( i == 0 );
172 return d_array.data();
173 }
174
175 typeID getType( size_t ) const override { return getTypeID<T>(); }
176 size_t sizeofDataBlockType( size_t ) const override { return sizeof( T ); }
177 void swapData( VectorData & ) override;
178 std::shared_ptr<VectorData> cloneData( const std::string &name = "" ) const override;
179};
180
181} // namespace AMP::LinearAlgebra
182
183#include "AMP/vectors/data/ArrayVectorData.hpp"
184
185#endif
Provides C++ wrapper around MPI routines.
Definition AMP_MPI.h:63
Simple class to store the array dimensions.
Definition ArraySize.h:138
size_t sizeOfDataBlock(size_t i=0) const override
Number of elements in a data block.
ArrayVectorData(const ArrayVectorData &)=delete
Copy constructor.
ArrayVectorData()=default
Empty constructor.
static std::shared_ptr< ArrayVectorData > create(const size_t localSize, std::shared_ptr< CommunicationList > commList, T *data)
void setValuesByLocalID(size_t num, const size_t *indices, const void *vals, const typeID &id) override
Set values in the vector by their local offset.
void swapData(VectorData &) override
Swap the data with another VectorData object.
const Array< T, FUN, Allocator > & getArray(void) const
return a const reference to the internal data container
const AMP_MPI & getComm() const override
Return the communicator this Vector spans.
void resize(const ArraySize &localDims)
resize the ArrayVector and reset the internal data structures
void putRawData(const void *buf, const typeID &id) override
Copy data into this vector.
void getValuesByLocalID(size_t num, const size_t *indices, void *vals, const typeID &id) const override
Get local values in the vector by their lcoal offset.
static std::shared_ptr< ArrayVectorData > create(const ArraySize &localSize, const ArraySize &blockIndex, AMP_MPI comm)
Create a ArrayVector.
size_t numberOfDataBlocks() const override
Number of blocks of contiguous data in the Vector.
void * getRawDataBlockAsVoid(size_t i) override
Return a pointer to a particular block of memory in the vector.
size_t sizeofDataBlockType(size_t) const override
Return the result of sizeof(TYPE) for the given data block.
static std::shared_ptr< ArrayVectorData > create(const ArraySize &localSize)
Create a ArrayVector.
std::shared_ptr< VectorData > cloneData(const std::string &name="") const override
Clone the data.
std::string VectorDataName() const override
Get the type name.
void getRawData(void *buf, const typeID &id) const override
Copy data out of this vector.
uint64_t getDataID() const override
A unique id for the underlying data allocation.
Array< T, FUN, Allocator > & getArray(void)
return a non-const reference to the internal data container
AMP::Array< T, FUN, Allocator > d_array
void addValuesByLocalID(size_t num, const size_t *indices, const void *vals, const typeID &id) override
Add values to vector entities by their local offset.
typeID getType(size_t) const override
Return the typeid of the given block.
const void * getRawDataBlockAsVoid(size_t i) const override
Return a pointer to a particular block of memory in the vector.
A class used to hold vector data.
Definition VectorData.h:38
#define AMP_ASSERT(EXP)
Assert error.
Class to store type info.
Definition typeid.h:210



Advanced Multi-Physics (AMP)
Oak Ridge National Laboratory
Idaho National Laboratory
Los Alamos National Laboratory
This page automatically produced from the
source code by doxygen
Last updated: Tue Mar 10 2026 13:06:41.
Comments on this page