Advanced Multi-Physics (AMP)
On-Line Documentation
NativePetscVectorData.h
Go to the documentation of this file.
1#ifndef included_AMP_NativePetscVectorData
2#define included_AMP_NativePetscVectorData
3
4#include "AMP/utils/AMP_MPI.h"
5#include "AMP/vectors/data/GhostDataHelper.hpp"
6#include "AMP/vectors/data/VectorData.h"
7
8#include "petscvec.h"
9
10
11namespace AMP::LinearAlgebra {
12
13
25class NativePetscVectorData : public GhostDataHelper<PetscScalar>
26{
27public:
33 explicit NativePetscVectorData( Vec v, bool deleteable, AMP_MPI comm = AMP_MPI() );
34
37
38 // Overrides to makeConsistent ensure that the internal Petsc Vector
39 // is in a valid state by calling assemble/resetArray as needed
40 void makeConsistent( ScatterType t ) override;
41 void makeConsistent() override;
42
43 std::string VectorDataName() const override { return "NativePetscVector"; }
44 size_t numberOfDataBlocks() const override;
45 size_t sizeOfDataBlock( size_t i ) const override;
46 void putRawData( const void *, const typeID & ) override;
47 void getRawData( void *, const typeID & ) const override;
48
49 void setValuesByLocalID( size_t, const size_t *, const void *, const typeID & ) override;
50 void addValuesByLocalID( size_t, const size_t *, const void *, const typeID & ) override;
51 void getValuesByLocalID( size_t, const size_t *, void *, const typeID & ) const override;
52
53 // Return the id of the data
54 uint64_t getDataID() const override
55 {
56 return reinterpret_cast<uint64_t>( getRawDataBlockAsVoid( 0 ) );
57 }
58 void *getRawDataBlockAsVoid( size_t i ) override;
59 const void *getRawDataBlockAsVoid( size_t i ) const override;
60 size_t sizeofDataBlockType( size_t ) const override { return sizeof( PetscScalar ); }
61 typeID getType( size_t ) const override { return getTypeID<PetscScalar>(); }
62 void swapData( VectorData &rhs ) override;
63
66 std::shared_ptr<VectorData> cloneData( const std::string &name = "" ) const override;
67
68 void assemble() override;
69
71 Vec &getVec() { return d_petscVec; }
72
74 const Vec &getVec() const { return d_petscVec; }
75
76protected:
77 void resetArray();
78 void resetArray() const;
79
80private:
82 bool d_bDeleteMe = false;
83 Vec d_petscVec = nullptr;
84 mutable PetscScalar *d_pArray = nullptr; // mutable so that we can cache the value
85 mutable const PetscScalar *d_pArrayRead = nullptr;
86};
87
88
89} // namespace AMP::LinearAlgebra
90
91#endif
Provides C++ wrapper around MPI routines.
Definition AMP_MPI.h:63
std::shared_ptr< VectorData > cloneData(const std::string &name="") const override
Clone the data.
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.
std::string VectorDataName() const override
Get the type name.
void makeConsistent() override
Update shared values on entire communicator.
void assemble() override
This method is used to implement the assemble interface of PETSc.
void makeConsistent(ScatterType t) override
Update shared values on entire communicator.
void getValuesByLocalID(size_t, const size_t *, void *, const typeID &) const override
Get local values in the vector by their global offset.
const void * getRawDataBlockAsVoid(size_t i) const override
Return a pointer to a particular block of memory in the vector.
const Vec & getVec() const
Get the PETSc vector.
NativePetscVectorData(Vec v, bool deleteable, AMP_MPI comm=AMP_MPI())
Construct a wrapper for a PETSc Vec from a set of parameters.
typeID getType(size_t) const override
Return the typeid of the given block.
size_t sizeOfDataBlock(size_t i) const override
Number of elements in a data block.
uint64_t getDataID() const override
A unique id for the underlying data allocation.
void addValuesByLocalID(size_t, const size_t *, const void *, const typeID &) override
Add values to vector entities by their local offset.
void getRawData(void *, const typeID &) const override
Copy data out of this vector.
void setValuesByLocalID(size_t, const size_t *, const void *, const typeID &) override
Set values in the vector by their local offset.
void putRawData(const void *, const typeID &) override
Copy data into this vector.
size_t numberOfDataBlocks() const override
Number of blocks of contiguous data in the Vector.
void swapData(VectorData &rhs) override
Swap the data with another VectorData object.
A class used to hold vector data.
Definition VectorData.h:38
ScatterType
Flag to choose algorithm for makeConsistent.
Definition Variable.h:21
Class to store type info.
Definition typeid.h:210
PETSc vector.



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