GEOSX
Public Types | Public Member Functions | Protected Attributes | List of all members
geosx::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 geosx::PetscVector:
Inheritance graph
[legend]

Public Types

using Vec = struct _p_Vec *
 Alias for PETSc vector struct pointer.
 

Public Member Functions

const Vecunwrapped () const
 Returns a const pointer to the underlying Vec. More...
 
Vecunwrapped ()
 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 createWithLocalSize (localIndex const localSize, MPI_Comm const &comm) override
 Create a vector based on local number of elements. More...
 
virtual void createWithGlobalSize (globalIndex const globalSize, MPI_Comm const &comm) override
 Create a vector based on global number of elements. More...
 
virtual void create (arrayView1d< real64 const > const &localValues, MPI_Comm const &comm) override
 Construct parallel vector from a local array. More...
 
virtual void open () override
 Open the vector for modifying entries.
 
virtual void close () override
 Assemble vector. More...
 
virtual void reset () override
 Reset the matrix to default state.
 
virtual void set (globalIndex const globalRow, real64 const value) override
 Set vector value. More...
 
virtual void add (globalIndex const globalRow, real64 const value) override
 Add into vector value. More...
 
virtual void set (globalIndex const *globalIndices, real64 const *values, localIndex size) override
 Set vector values. More...
 
virtual void add (globalIndex const *globalIndices, real64 const *values, localIndex size) override
 Add vector values. More...
 
virtual void set (arraySlice1d< globalIndex const > const &globalIndices, arraySlice1d< real64 const > const &values) override
 Set vector values using array1d. More...
 
virtual void add (arraySlice1d< globalIndex const > const &globalIndices, arraySlice1d< real64 const > const &values) override
 Add into vector values using array1d. More...
 
virtual void set (real64 const value) override
 Set all elements to a constant value. More...
 
virtual void zero () override
 Set vector elements to zero.
 
virtual void rand (unsigned const seed=1984) override
 Set vector elements to random entries. More...
 
virtual void scale (real64 const scalingFactor) override
 Multiply all elements by scalingFactor. 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 componentwise 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 real64 get (globalIndex const globalRow) const override
 Get a value by index. More...
 
void get (arraySlice1d< globalIndex const > const &globalIndices, arraySlice1d< real64 > const &values) const override
 Get a sequence of values by index. 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 localIndex getLocalRowID (globalIndex const globalRow) const override
 Map a global row index to local row index. More...
 
virtual globalIndex getGlobalRowID (localIndex const localRow) const override
 Map a local row index to global row index. More...
 
virtual real64 const * extractLocalVector () const override
 Extract a view of the local portion of the array. More...
 
virtual real64extractLocalVector () override
 Extract a view of the local portion of the array. More...
 
virtual MPI_Comm getComm () const override
 Get the communicator used by this vector. More...
 

Protected Attributes

Vec m_vec
 

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]

geosx::PetscVector::PetscVector ( PetscVector const &  src)

Copy constructor.

Parameters
srcPetscVector to be copied.

◆ PetscVector() [2/2]

geosx::PetscVector::PetscVector ( PetscVector &&  src)
noexcept

Move constructor.

Parameters
srcPetscVector to move from

Member Function Documentation

◆ add() [1/3]

virtual void geosx::PetscVector::add ( globalIndex const  globalRow,
real64 const  value 
)
overridevirtual

Add into vector value.

Parameters
globalRowglobal row
valueValues to add in given row

Add into vector value at given row.

Implements geosx::VectorBase< PetscVector >.

◆ add() [2/3]

virtual void geosx::PetscVector::add ( globalIndex const *  globalIndices,
real64 const *  values,
localIndex  size 
)
overridevirtual

Add vector values.

Parameters
globalIndicesglobal row indices
valuesvalues to add in given rows
sizenumber of elements

Add vector values at given elements.

Implements geosx::VectorBase< PetscVector >.

◆ add() [3/3]

virtual void geosx::PetscVector::add ( arraySlice1d< globalIndex const > const &  globalIndices,
arraySlice1d< real64 const > const &  values 
)
overridevirtual

Add into vector values using array1d.

Parameters
globalIndicesglobal rows indices
valuesvalues to add in given rows

Add into vector values at given rows.

Implements geosx::VectorBase< PetscVector >.

◆ axpby()

virtual void geosx::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 geosx::VectorBase< PetscVector >.

◆ axpy()

virtual void geosx::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 geosx::VectorBase< PetscVector >.

◆ close()

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

Assemble vector.

Performs parallel communication to scatter assembled entries to appropriate locations

Implements geosx::VectorBase< PetscVector >.

◆ copy()

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

Update vector y as y = x.

Parameters
xvector to copy

Implements geosx::VectorBase< PetscVector >.

◆ create()

virtual void geosx::PetscVector::create ( arrayView1d< real64 const > const &  localValues,
MPI_Comm const &  comm 
)
overridevirtual

Construct parallel vector from a local array.

Parameters
localValueslocal data to put into vector
commMPI communicator to use

Implements geosx::VectorBase< PetscVector >.

◆ created()

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

Query vector creation status.

Returns
true if vector has been created

Implements geosx::VectorBase< PetscVector >.

◆ createWithGlobalSize()

virtual void geosx::PetscVector::createWithGlobalSize ( globalIndex const  globalSize,
MPI_Comm const &  comm 
)
overridevirtual

Create a vector based on global number of elements.

Parameters
globalSizeGlobal number of elements
commMPI communicator to use

Create a vector based on global number of elements. Every processors gets the same number of local elements except proc 0, which gets any remainder elements as well if the split can't be done evenly.

Implements geosx::VectorBase< PetscVector >.

◆ createWithLocalSize()

virtual void geosx::PetscVector::createWithLocalSize ( 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. For specifying a global size and having automatic partitioning, see createWithGlobalSize().

Implements geosx::VectorBase< PetscVector >.

◆ dot()

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

Dot product with the vector vec.

Parameters
vecvector to dot-product with
Returns
dot product

Implements geosx::VectorBase< PetscVector >.

◆ extractLocalVector() [1/2]

virtual real64 const* geosx::PetscVector::extractLocalVector ( ) const
overridevirtual

Extract a view of the local portion of the array.

Returns
pointer to local vector data

Implements geosx::VectorBase< PetscVector >.

◆ extractLocalVector() [2/2]

virtual real64* geosx::PetscVector::extractLocalVector ( )
overridevirtual

Extract a view of the local portion of the array.

Returns
pointer to local vector data

Implements geosx::VectorBase< PetscVector >.

◆ get() [1/2]

virtual real64 geosx::PetscVector::get ( globalIndex const  globalRow) const
overridevirtual

Get a value by index.

Parameters
globalRowglobal row index
Returns
value at global index globalRow

Implements geosx::VectorBase< PetscVector >.

◆ get() [2/2]

void geosx::PetscVector::get ( arraySlice1d< globalIndex const > const &  globalIndices,
arraySlice1d< real64 > const &  values 
) const
overridevirtual

Get a sequence of values by index.

Parameters
[in]globalIndicesarray of global row indices
[out]valuesarray of vector values

Implements geosx::VectorBase< PetscVector >.

◆ getComm()

virtual MPI_Comm geosx::PetscVector::getComm ( ) const
overridevirtual

Get the communicator used by this vector.

Returns
the MPI communicator

Implements geosx::VectorBase< PetscVector >.

◆ getGlobalRowID()

virtual globalIndex geosx::PetscVector::getGlobalRowID ( localIndex const  localRow) const
overridevirtual

Map a local row index to global row index.

Parameters
[in]localRowthe local row index
Returns
global row index corresponding to localRow

Implements geosx::VectorBase< PetscVector >.

◆ getLocalRowID()

virtual localIndex geosx::PetscVector::getLocalRowID ( globalIndex const  globalRow) const
overridevirtual

Map a global row index to local row index.

Parameters
[in]globalRowthe global row index
Returns
global row index corresponding to globalRow

Implements geosx::VectorBase< PetscVector >.

◆ globalSize()

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

Returns the global of the vector.

Returns
the global size

Implements geosx::VectorBase< PetscVector >.

◆ ilower()

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

Get lower bound of local partition.

Returns
index of the first global row owned by this processor

Implements geosx::VectorBase< PetscVector >.

◆ iupper()

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

Get upper bound of local partition.

Returns
next index after last global row owned by that processor
Note
The intention is for [ilower; iupper) to be used as a half-open index range

Implements geosx::VectorBase< PetscVector >.

◆ localSize()

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

Returns the local size of the vector.

Returns
the local size (on this processor)

Implements geosx::VectorBase< PetscVector >.

◆ norm1()

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

1-norm of the vector.

Returns
the 1-norm value

Implements geosx::VectorBase< PetscVector >.

◆ norm2()

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

2-norm of the vector.

Returns
the 2-norm value

Implements geosx::VectorBase< PetscVector >.

◆ normInf()

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

Infinity-norm of the vector.

Returns
the inf-norm value

Implements geosx::VectorBase< PetscVector >.

◆ operator=() [1/2]

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

Copy assignment.

Parameters
srcPetscVector to be copied.
Returns
the new vector.

◆ operator=() [2/2]

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

Move assignment.

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

◆ pointwiseProduct()

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

Compute the componentwise multiplication y = v * x.

Parameters
xfirst vector (input)
ysecond vector (output)

Implements geosx::VectorBase< PetscVector >.

◆ print()

virtual void geosx::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 geosx::VectorBase< PetscVector >.

◆ rand()

virtual void geosx::PetscVector::rand ( unsigned const  seed = 1984)
overridevirtual

Set vector elements to random entries.

Parameters
seedthe random number seed to use

Implements geosx::VectorBase< PetscVector >.

◆ reciprocal()

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

Replace vector elements by their reciprocals.

Note
No guarding is done against division by zero.

Implements geosx::VectorBase< PetscVector >.

◆ scale()

virtual void geosx::PetscVector::scale ( real64 const  scalingFactor)
overridevirtual

Multiply all elements by scalingFactor.

Parameters
scalingFactorscaling Factor

Implements geosx::VectorBase< PetscVector >.

◆ set() [1/4]

virtual void geosx::PetscVector::set ( globalIndex const  globalRow,
real64 const  value 
)
overridevirtual

Set vector value.

Parameters
globalRowglobal row index
valueValue to add at given row

Set vector value at given element.

Implements geosx::VectorBase< PetscVector >.

◆ set() [2/4]

virtual void geosx::PetscVector::set ( globalIndex const *  globalIndices,
real64 const *  values,
localIndex  size 
)
overridevirtual

Set vector values.

Parameters
globalIndicesglobal row indices
valuesValues to add in given rows
sizeNumber of elements

Set vector values at given elements.

Implements geosx::VectorBase< PetscVector >.

◆ set() [3/4]

virtual void geosx::PetscVector::set ( arraySlice1d< globalIndex const > const &  globalIndices,
arraySlice1d< real64 const > const &  values 
)
overridevirtual

Set vector values using array1d.

Parameters
globalIndicesglobal row indices
valuesvalues to add in given rows

Set vector values at given elements.

Implements geosx::VectorBase< PetscVector >.

◆ set() [4/4]

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

Set all elements to a constant value.

Parameters
valuevalue to set vector elements to

Implements geosx::VectorBase< PetscVector >.

◆ unwrapped() [1/2]

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

Returns a const pointer to the underlying Vec.

Returns
the const pointer to the underlying Vec.

◆ unwrapped() [2/2]

Vec& geosx::PetscVector::unwrapped ( )

Returns a non-const pointer to the underlying Vec.

Returns
the non-const pointer to the underlying Vec.

◆ write()

virtual void geosx::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 geosx::VectorBase< PetscVector >.

Member Data Documentation

◆ m_vec

Vec geosx::PetscVector::m_vec
protected

Pointer to underlying PETSc Vec

Definition at line 250 of file PetscVector.hpp.


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