GEOSX
Public Member Functions | List of all members
geos::PetscVector Class Referencefinal

This class creates and provides basic support for Vec vector object type used in PETSc. More...

#include <PetscVector.hpp>

Inheritance diagram for geos::PetscVector:
Inheritance graph
[legend]

Public Member Functions

const Vec & unwrapped () const
 Returns a const pointer to the underlying Vec. More...
 
Vec & unwrapped ()
 Returns a non-const pointer to the underlying Vec. More...
 
Constructor/Destructor Methods
 PetscVector ()
 Empty vector constructor.
 
 PetscVector (PetscVector const &src)
 Copy constructor. More...
 
 PetscVector (PetscVector &&src) noexcept
 Move constructor. More...
 
PetscVectoroperator= (PetscVector const &src)
 Copy assignment. More...
 
PetscVectoroperator= (PetscVector &&src) noexcept
 Move assignment. More...
 
 ~PetscVector ()
 Destructor.
 
VectorBase interface
virtual bool created () const override
 Query vector creation status. More...
 
virtual void create (localIndex const localSize, MPI_Comm const &comm) override
 Create a vector based on local number of elements. More...
 
virtual void close () override
 Close vector for modification. More...
 
virtual void touch () override
 Notify the vector about external modification through direct data pointer. More...
 
virtual void reset () override
 Reset the vector to default state.
 
virtual void set (real64 const value) override
 Set all elements to a constant value. More...
 
virtual void rand (unsigned const seed) override
 Set vector elements to random entries. More...
 
virtual void scale (real64 const scalingFactor) override
 Multiply all elements by factor. More...
 
virtual void reciprocal () override
 Replace vector elements by their reciprocals. More...
 
virtual real64 dot (PetscVector const &vec) const override
 Dot product with the vector vec. More...
 
virtual void copy (PetscVector const &x) override
 Update vector y as y = x. More...
 
virtual void axpy (real64 const alpha, PetscVector const &x) override
 Update vector y as y = alpha*x + y. More...
 
virtual void axpby (real64 const alpha, PetscVector const &x, real64 const beta) override
 Update vector y as y = alpha*x + beta*y. More...
 
virtual void pointwiseProduct (PetscVector const &x, PetscVector &y) const override
 Compute the component-wise multiplication y = v * x. More...
 
virtual real64 norm1 () const override
 1-norm of the vector. More...
 
virtual real64 norm2 () const override
 2-norm of the vector. More...
 
virtual real64 normInf () const override
 Infinity-norm of the vector. More...
 
virtual globalIndex globalSize () const override
 Returns the global of the vector. More...
 
virtual localIndex localSize () const override
 Returns the local size of the vector. More...
 
virtual globalIndex ilower () const override
 Get lower bound of local partition. More...
 
virtual globalIndex iupper () const override
 Get upper bound of local partition. More...
 
virtual void print (std::ostream &os=std::cout) const override
 Print the vector in Trilinos format to the terminal. More...
 
virtual void write (string const &filename, LAIOutputFormat const format=LAIOutputFormat::MATRIX_MARKET) const override
 Write the vector to a file. More...
 
virtual MPI_Comm comm () const override
 Get the communicator used by this vector. More...
 
void setName (string const &name)
 Set a name for the vector (mainly used during various logging). More...
 
bool closed () const
 Query vector closed status. More...
 
bool ready () const
 Query vector ready status. More...
 
virtual arrayView1d< real64open ()
 Open the vector for modifying entries. More...
 
virtual void zero ()
 Set vector elements to zero.
 
arrayView1d< real64 const > values () const
 

Detailed Description

This class creates and provides basic support for Vec vector object type used in PETSc.

Definition at line 44 of file PetscVector.hpp.

Constructor & Destructor Documentation

◆ PetscVector() [1/2]

geos::PetscVector::PetscVector ( PetscVector const &  src)

Copy constructor.

Parameters
srcPetscVector to be copied.

◆ PetscVector() [2/2]

geos::PetscVector::PetscVector ( PetscVector &&  src)
noexcept

Move constructor.

Parameters
srcPetscVector to move from

Member Function Documentation

◆ axpby()

virtual void geos::PetscVector::axpby ( real64 const  alpha,
PetscVector const &  x,
real64 const  beta 
)
overridevirtual

Update vector y as y = alpha*x + beta*y.

Parameters
alphascaling factor for added vector
xvector to add
betascaling factor for self vector

Implements geos::VectorBase< PetscVector >.

◆ axpy()

virtual void geos::PetscVector::axpy ( real64 const  alpha,
PetscVector const &  x 
)
overridevirtual

Update vector y as y = alpha*x + y.

Parameters
alphascaling factor for added vector
xvector to add

Implements geos::VectorBase< PetscVector >.

◆ close()

virtual void geos::PetscVector::close ( )
overridevirtual

Close vector for modification.

After calling this method, the view obtained via open() should not be used.

Implements geos::VectorBase< PetscVector >.

◆ closed()

bool geos::VectorBase< VECTOR >::closed
inline

Query vector closed status.

Returns
true if vector has been opened and has not been closed since; false otherwise

Definition at line 63 of file VectorBase.hpp.

◆ comm()

virtual MPI_Comm geos::PetscVector::comm ( ) const
overridevirtual

Get the communicator used by this vector.

Returns
the MPI communicator

◆ copy()

virtual void geos::PetscVector::copy ( PetscVector const &  x)
overridevirtual

Update vector y as y = x.

Parameters
xvector to copy
Note
Unlike copy assignment operator, this method expects both vectors to be created and have identical parallel distributions, and never reallocates memory.

Implements geos::VectorBase< PetscVector >.

◆ create()

virtual void geos::PetscVector::create ( localIndex const  localSize,
MPI_Comm const &  comm 
)
overridevirtual

Create a vector based on local number of elements.

Parameters
localSizelocal number of elements
commMPI communicator to use

Create a vector based on local number of elements. Global size is the sum across processors.

Reimplemented from geos::VectorBase< PetscVector >.

◆ created()

virtual bool geos::PetscVector::created ( ) const
overridevirtual

Query vector creation status.

Returns
true if vector has been created

◆ dot()

virtual real64 geos::PetscVector::dot ( PetscVector const &  vec) const
overridevirtual

Dot product with the vector vec.

Parameters
vecvector to dot-product with
Returns
dot product

Implements geos::VectorBase< PetscVector >.

◆ globalSize()

virtual globalIndex geos::PetscVector::globalSize ( ) const
overridevirtual

Returns the global of the vector.

Returns
the global size

◆ ilower()

virtual globalIndex geos::PetscVector::ilower ( ) const
overridevirtual

Get lower bound of local partition.

Returns
index of the first global row owned by this processor

◆ iupper()

virtual globalIndex geos::PetscVector::iupper ( ) const
overridevirtual

Get upper bound of local partition.

Returns
next index after last global row owned by that processor
Note
[ v.ilower(); v.iupper() ) is a half-open index range

◆ localSize()

virtual localIndex geos::PetscVector::localSize ( ) const
overridevirtual

Returns the local size of the vector.

Returns
the local size (on this processor)

◆ norm1()

virtual real64 geos::PetscVector::norm1 ( ) const
overridevirtual

1-norm of the vector.

Returns
the 1-norm value

◆ norm2()

virtual real64 geos::PetscVector::norm2 ( ) const
overridevirtual

2-norm of the vector.

Returns
the 2-norm value

◆ normInf()

virtual real64 geos::PetscVector::normInf ( ) const
overridevirtual

Infinity-norm of the vector.

Returns
the inf-norm value

◆ open()

virtual arrayView1d< real64 > geos::VectorBase< VECTOR >::open
inline

Open the vector for modifying entries.

Returns
an array view to assemble local values into

Definition at line 128 of file VectorBase.hpp.

◆ operator=() [1/2]

PetscVector& geos::PetscVector::operator= ( PetscVector &&  src)
noexcept

Move assignment.

Parameters
srcPetscVector to be moved from.
Returns
the new vector.

◆ operator=() [2/2]

PetscVector& geos::PetscVector::operator= ( PetscVector const &  src)

Copy assignment.

Parameters
srcPetscVector to be copied.
Returns
the new vector.

◆ pointwiseProduct()

virtual void geos::PetscVector::pointwiseProduct ( PetscVector const &  x,
PetscVector y 
) const
overridevirtual

Compute the component-wise multiplication y = v * x.

Parameters
xfirst vector (input)
ysecond vector (output)

Implements geos::VectorBase< PetscVector >.

◆ print()

virtual void geos::PetscVector::print ( std::ostream &  os = std::cout) const
overridevirtual

Print the vector in Trilinos format to the terminal.

Parameters
osthe output stream to print to

Implements geos::VectorBase< PetscVector >.

◆ rand()

virtual void geos::PetscVector::rand ( unsigned const  seed)
overridevirtual

Set vector elements to random entries.

Parameters
seedthe random number seed to use

Implements geos::VectorBase< PetscVector >.

◆ ready()

bool geos::VectorBase< VECTOR >::ready
inline

Query vector ready status.

Returns
true if vector has been created and is currently closed

Definition at line 75 of file VectorBase.hpp.

◆ reciprocal()

virtual void geos::PetscVector::reciprocal ( )
overridevirtual

Replace vector elements by their reciprocals.

Note
No guarding is done against division by zero.

Implements geos::VectorBase< PetscVector >.

◆ scale()

virtual void geos::PetscVector::scale ( real64 const  factor)
overridevirtual

Multiply all elements by factor.

Parameters
factorscaling factor

Implements geos::VectorBase< PetscVector >.

◆ set()

virtual void geos::PetscVector::set ( real64 const  value)
overridevirtual

Set all elements to a constant value.

Parameters
valuevalue to set vector elements to

Implements geos::VectorBase< PetscVector >.

◆ setName()

void geos::VectorBase< VECTOR >::setName
inline

Set a name for the vector (mainly used during various logging).

Parameters
namethe name

Definition at line 112 of file VectorBase.hpp.

◆ touch()

virtual void geos::PetscVector::touch ( )
overridevirtual

Notify the vector about external modification through direct data pointer.

This method MUST be called after any changes to vector values performed through direct pointer access, so that the underlying Array object can be made aware of external changes in a specific memory space.

Implements geos::VectorBase< PetscVector >.

◆ unwrapped() [1/2]

Vec& geos::PetscVector::unwrapped ( )

Returns a non-const pointer to the underlying Vec.

Returns
the non-const pointer to the underlying Vec.

◆ unwrapped() [2/2]

const Vec& geos::PetscVector::unwrapped ( ) const

Returns a const pointer to the underlying Vec.

Returns
the const pointer to the underlying Vec.

◆ values()

arrayView1d< real64 const > geos::VectorBase< VECTOR >::values
inline
Returns
a const access view to local vector values

Definition at line 304 of file VectorBase.hpp.

◆ write()

virtual void geos::PetscVector::write ( string const &  filename,
LAIOutputFormat const  format = LAIOutputFormat::MATRIX_MARKET 
) const
overridevirtual

Write the vector to a file.

Parameters
filenamename of the output file
[in]formatoutput format

Implements geos::VectorBase< PetscVector >.


The documentation for this class was generated from the following file: