IBAMR
An adaptive and distributed-memory parallel implementation of the immersed boundary (IB) method
Static Public Member Functions | Protected Member Functions | List of all members
IBTK::PETScSAMRAIVectorReal Class Reference

Class PETScSAMRAIVectorReal is a class for wrapping SAMRAI::solv::SAMRAIVectorReal objects so that they may be used with the PETSc solver package. More...

#include </home/runner/work/IBAMR/IBAMR/ibtk/include/ibtk/PETScSAMRAIVectorReal.h>

Static Public Member Functions

static Vec createPETScVector (SAMRAI::tbox::Pointer< SAMRAI::solv::SAMRAIVectorReal< NDIM, PetscScalar > > samrai_vec, MPI_Comm comm=PETSC_COMM_WORLD)
 
static void destroyPETScVector (Vec petsc_vec)
 
static void getSAMRAIVector (Vec petsc_vec, SAMRAI::tbox::Pointer< SAMRAI::solv::SAMRAIVectorReal< NDIM, PetscScalar > > *samrai_vec)
 
static void restoreSAMRAIVector (Vec petsc_vec, SAMRAI::tbox::Pointer< SAMRAI::solv::SAMRAIVectorReal< NDIM, PetscScalar > > *samrai_vec)
 
static void getSAMRAIVectorRead (Vec petsc_vec, SAMRAI::tbox::Pointer< SAMRAI::solv::SAMRAIVectorReal< NDIM, PetscScalar > > *samrai_vec)
 
static void restoreSAMRAIVectorRead (Vec petsc_vec, SAMRAI::tbox::Pointer< SAMRAI::solv::SAMRAIVectorReal< NDIM, PetscScalar > > *samrai_vec)
 
static void replaceSAMRAIVector (Vec petsc_vec, SAMRAI::tbox::Pointer< SAMRAI::solv::SAMRAIVectorReal< NDIM, PetscScalar > > samrai_vec)
 

Protected Member Functions

 PETScSAMRAIVectorReal (SAMRAI::tbox::Pointer< SAMRAI::solv::SAMRAIVectorReal< NDIM, PetscScalar > > samrai_vector, bool vector_created_via_duplicate, MPI_Comm comm)
 

Detailed Description

Class PETScSAMRAIVectorReal is a class for wrapping SAMRAI::solv::SAMRAIVectorReal objects so that they may be used with the PETSc solver package.

Class PETScSAMRAIVectorReal wraps a real-valued SAMRAI vector (see SAMRAI::solv::SAMRAIVectorReal class) object so that it may be used with the PETSc solver package. A SAMRAI vector is defined as a collection of patch data components and associated operations living on some subset of levels in a structured AMR mesh hierarchy.

Observe that there are only three public member functions in this class. They are used to create and destroy PETSc vectors (i.e., Vec objects) and to obtain the SAMRAI vector associated with the PETSc vector. In particular, note that the constructor and destructor of this class are protected members. The construction and destruction of instances of this class may occur only through the static member functions that create and destroy PETSc vector objects.

Finally, we remark that PETSc allows vectors with complex-valued entries. This class and the class SAMRAI::solv::SAMRAIVectorReal assume real-values vectors, i.e., data of type double or float. The (currently unimplemented) class PETScSAMRAIVectorComplex must be used for complex data.

See also
SAMRAI::solv::SAMRAIVectorReal

Member Function Documentation

◆ createPETScVector()

Vec IBTK::PETScSAMRAIVectorReal::createPETScVector ( SAMRAI::tbox::Pointer< SAMRAI::solv::SAMRAIVectorReal< NDIM, PetscScalar > >  samrai_vec,
MPI_Comm  comm = PETSC_COMM_WORLD 
)
inlinestatic

Create and return a new PETSc vector object that wraps the SAMRAI vector object, so that the SAMRAI vector may be manipulated by PETSc routines. It is important to note that this function does not allocate storage for the vector data. Data must be allocated through the SAMRAI vector object directly.

Note
Each call to createPETScVector() should be matched with a corresponding call to destroyPETScVector().

◆ destroyPETScVector()

void IBTK::PETScSAMRAIVectorReal::destroyPETScVector ( Vec  petsc_vec)
inlinestatic

Destroy a given PETSc vector object. It is important to note that this function does not deallocate storage for the vector data. Vector data must be deallocated through the SAMRAI vector object.

Note
Each call to createPETScVector() should be matched with a corresponding call to destroyPETScVector().

◆ getSAMRAIVector()

void IBTK::PETScSAMRAIVectorReal::getSAMRAIVector ( Vec  petsc_vec,
SAMRAI::tbox::Pointer< SAMRAI::solv::SAMRAIVectorReal< NDIM, PetscScalar > > *  samrai_vec 
)
inlinestatic

Get a pointer to the SAMRAI vector object associated with the given PETSc vector object.

Note
The SAMRAI vector must be restored by calling restoreSAMRAIVector().

◆ getSAMRAIVectorRead()

void IBTK::PETScSAMRAIVectorReal::getSAMRAIVectorRead ( Vec  petsc_vec,
SAMRAI::tbox::Pointer< SAMRAI::solv::SAMRAIVectorReal< NDIM, PetscScalar > > *  samrai_vec 
)
inlinestatic

Get a pointer to the SAMRAI vector object associated with the given PETSc vector object. This vector must be treated as read only.

Note
The SAMRAI vector must be restored by calling restoreSAMRAIVectorRead().

◆ replaceSAMRAIVector()

void IBTK::PETScSAMRAIVectorReal::replaceSAMRAIVector ( Vec  petsc_vec,
SAMRAI::tbox::Pointer< SAMRAI::solv::SAMRAIVectorReal< NDIM, PetscScalar > >  samrai_vec 
)
inlinestatic

Replace the SAMRAI vector object associated with the given PETSc vector object.

◆ restoreSAMRAIVector()

void IBTK::PETScSAMRAIVectorReal::restoreSAMRAIVector ( Vec  petsc_vec,
SAMRAI::tbox::Pointer< SAMRAI::solv::SAMRAIVectorReal< NDIM, PetscScalar > > *  samrai_vec 
)
inlinestatic

Restore the SAMRAI vector object associated with the given PETSc vector object.

◆ restoreSAMRAIVectorRead()

void IBTK::PETScSAMRAIVectorReal::restoreSAMRAIVectorRead ( Vec  petsc_vec,
SAMRAI::tbox::Pointer< SAMRAI::solv::SAMRAIVectorReal< NDIM, PetscScalar > > *  samrai_vec 
)
inlinestatic

Restore the SAMRAI vector object associated with the given PETSc vector object.


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