libMesh
|
This class provides a nice interface to PETSc's Vec object. More...
#include <petsc_vector.h>
Public Member Functions | |
PetscVector (const Parallel::Communicator &comm_in, const ParallelType type=AUTOMATIC) | |
Dummy-Constructor. More... | |
PetscVector (const Parallel::Communicator &comm_in, const numeric_index_type n, const ParallelType type=AUTOMATIC) | |
Constructor. More... | |
PetscVector (const Parallel::Communicator &comm_in, const numeric_index_type n, const numeric_index_type n_local, const ParallelType type=AUTOMATIC) | |
Constructor. More... | |
PetscVector (const Parallel::Communicator &comm_in, const numeric_index_type N, const numeric_index_type n_local, const std::vector< numeric_index_type > &ghost, const ParallelType type=AUTOMATIC) | |
Constructor. More... | |
PetscVector (Vec v, const Parallel::Communicator &comm_in) | |
Constructor. More... | |
PetscVector< T > & | operator= (const PetscVector< T > &v) |
Copy assignment operator. More... | |
PetscVector (PetscVector &&)=delete | |
This class manages a C-style struct (Vec) manually, so we don't want to allow any automatic copy/move functions to be generated, and we can't default the destructor. More... | |
PetscVector (const PetscVector &)=delete | |
PetscVector & | operator= (PetscVector &&)=delete |
virtual | ~PetscVector () |
virtual void | close () override |
Calls the NumericVector's internal assembly routines, ensuring that the values are consistent across processors. More... | |
virtual void | clear () noexcept override |
clear() is called from the destructor, so it should not throw. More... | |
virtual void | zero () override |
Set all entries to zero. More... | |
virtual std::unique_ptr< NumericVector< T > > | zero_clone () const override |
virtual std::unique_ptr< NumericVector< T > > | clone () const override |
virtual void | init (const numeric_index_type N, const numeric_index_type n_local, const bool fast=false, const ParallelType type=AUTOMATIC) override |
Change the dimension of the vector to n . More... | |
virtual void | init (const numeric_index_type N, const bool fast=false, const ParallelType type=AUTOMATIC) override |
Call init() with n_local = N. More... | |
virtual void | init (const numeric_index_type N, const numeric_index_type n_local, const std::vector< numeric_index_type > &ghost, const bool fast=false, const ParallelType=AUTOMATIC) override |
Create a vector that holds the local indices plus those specified in the ghost argument. More... | |
virtual void | init (const NumericVector< T > &other, const bool fast=false) override |
Creates a vector that has the same dimension and storage type as other , including ghost dofs. More... | |
virtual NumericVector< T > & | operator= (const T s) override |
Sets all entries of the vector to the value s . More... | |
virtual NumericVector< T > & | operator= (const NumericVector< T > &v) override |
This looks like a copy assignment operator, but note that, unlike normal copy assignment operators, it is pure virtual. More... | |
virtual NumericVector< T > & | operator= (const std::vector< T > &v) override |
Sets (*this)(i) = v(i) for each entry of the vector. More... | |
virtual Real | min () const override |
virtual Real | max () const override |
virtual T | sum () const override |
virtual Real | l1_norm () const override |
virtual Real | l2_norm () const override |
virtual Real | linfty_norm () const override |
virtual numeric_index_type | size () const override |
virtual numeric_index_type | local_size () const override |
virtual numeric_index_type | first_local_index () const override |
virtual numeric_index_type | last_local_index () const override |
numeric_index_type | map_global_to_local_index (const numeric_index_type i) const |
virtual T | operator() (const numeric_index_type i) const override |
virtual void | get (const std::vector< numeric_index_type > &index, T *values) const override |
Access multiple components at once. More... | |
PetscScalar * | get_array () |
Get read/write access to the raw PETSc Vector data array. More... | |
const PetscScalar * | get_array_read () const |
Get read only access to the raw PETSc Vector data array. More... | |
void | restore_array () |
Restore the data array. More... | |
virtual NumericVector< T > & | operator+= (const NumericVector< T > &v) override |
Adds v to *this, \( \vec{u} \leftarrow \vec{u} + \vec{v} \). More... | |
virtual NumericVector< T > & | operator-= (const NumericVector< T > &v) override |
Subtracts v from *this, \( \vec{u} \leftarrow \vec{u} - \vec{v} \). More... | |
virtual void | reciprocal () override |
Computes the component-wise reciprocal, \( u_i \leftarrow \frac{1}{u_i} \, \forall i\). More... | |
virtual void | conjugate () override |
Negates the imaginary component of each entry in the vector. More... | |
virtual void | set (const numeric_index_type i, const T value) override |
Sets v(i) = value . More... | |
virtual void | add (const numeric_index_type i, const T value) override |
Adds value to the vector entry specified by i . More... | |
virtual void | add (const T s) override |
Adds s to each entry of the vector, \( u_i \leftarrow u_i + s \). More... | |
virtual void | add (const NumericVector< T > &v) override |
Adds v to this , \( \vec{u} \leftarrow \vec{u} + \vec{v} \). More... | |
virtual void | add (const T a, const NumericVector< T > &v) override |
Vector addition with a scalar multiple, \( \vec{u} \leftarrow \vec{u} + a\vec{v} \). More... | |
virtual void | add_vector (const T *v, const std::vector< numeric_index_type > &dof_indices) override |
Computes \( \vec{u} \leftarrow \vec{u} + \vec{v} \), where v is a pointer and each dof_indices [i] specifies where to add value v [i]. More... | |
virtual void | add_vector (const NumericVector< T > &v, const SparseMatrix< T > &A) override |
Computes \( \vec{u} \leftarrow \vec{u} + A \vec{v} \), i.e. More... | |
virtual void | add_vector_transpose (const NumericVector< T > &v, const SparseMatrix< T > &A) override |
Computes \( \vec{u} \leftarrow \vec{u} + A^T \vec{v} \), i.e. More... | |
void | add_vector_conjugate_transpose (const NumericVector< T > &v, const SparseMatrix< T > &A) |
\( U \leftarrow U + A^H v \). More... | |
virtual void | insert (const T *v, const std::vector< numeric_index_type > &dof_indices) override |
Inserts the entries of v in *this at the locations specified by v . More... | |
virtual void | scale (const T factor) override |
Scale each element of the vector by the given factor . More... | |
virtual NumericVector< T > & | operator*= (const NumericVector< T > &v) override |
Computes the component-wise multiplication of this vector's entries by another's, \( u_i \leftarrow u_i v_i \, \forall i\). More... | |
virtual NumericVector< T > & | operator/= (const NumericVector< T > &v) override |
Computes the component-wise division of this vector's entries by another's, \( u_i \leftarrow \frac{u_i}{v_i} \, \forall i\). More... | |
virtual void | abs () override |
Sets \( u_i \leftarrow |u_i| \) for each entry in the vector. More... | |
virtual T | dot (const NumericVector< T > &v) const override |
T | indefinite_dot (const NumericVector< T > &v) const |
virtual void | localize (std::vector< T > &v_local) const override |
Creates a copy of the global vector in the local vector v_local . More... | |
virtual void | localize (NumericVector< T > &v_local) const override |
Same, but fills a NumericVector<T> instead of a std::vector . More... | |
virtual void | localize (NumericVector< T > &v_local, const std::vector< numeric_index_type > &send_list) const override |
Creates a local vector v_local containing only information relevant to this processor, as defined by the send_list . More... | |
virtual void | localize (std::vector< T > &v_local, const std::vector< numeric_index_type > &indices) const override |
Fill in the local std::vector "v_local" with the global indices given in "indices". More... | |
virtual void | localize (const numeric_index_type first_local_idx, const numeric_index_type last_local_idx, const std::vector< numeric_index_type > &send_list) override |
Updates a local vector with selected values from neighboring processors, as defined by send_list . More... | |
virtual void | localize_to_one (std::vector< T > &v_local, const processor_id_type proc_id=0) const override |
Creates a local copy of the global vector in v_local only on processor proc_id . More... | |
virtual void | pointwise_mult (const NumericVector< T > &vec1, const NumericVector< T > &vec2) override |
Computes \( u_i \leftarrow u_i v_i \) (summation not implied) i.e. More... | |
virtual void | pointwise_divide (const NumericVector< T > &vec1, const NumericVector< T > &vec2) override |
Computes \( u_i \leftarrow \frac{v_{1,i}}{v_{2,i}} \) (summation not implied) i.e. More... | |
virtual void | print_matlab (const std::string &name="") const override |
Print the contents of the vector in Matlab's sparse matrix format. More... | |
virtual void | create_subvector (NumericVector< T > &subvector, const std::vector< numeric_index_type > &rows, bool supplying_global_rows=true) const override |
Fills in subvector from this vector using the indices in rows . More... | |
virtual void | swap (NumericVector< T > &v) override |
Swaps the contents of this with v . More... | |
virtual std::size_t | max_allowed_id () const override |
Vec | vec () |
Vec | vec () const |
virtual std::unique_ptr< NumericVector< T > > | get_subvector (const std::vector< numeric_index_type > &rows) override |
Creates a view into this vector using the indices in rows . More... | |
virtual void | restore_subvector (std::unique_ptr< NumericVector< T >> subvector, const std::vector< numeric_index_type > &rows) override |
Restores a view into this vector using the indices in rows . More... | |
template<> | |
void | localize_to_one (std::vector< Real > &v_local, const processor_id_type timpi_mpi_var(pid)) const |
template<> | |
void | localize_to_one (std::vector< Complex > &v_local, const processor_id_type pid) const |
virtual bool | initialized () const |
ParallelType | type () const |
ParallelType & | type () |
void | set_type (ParallelType t) |
Allow the user to change the ParallelType of the NumericVector under some circumstances. More... | |
virtual bool | closed () const |
virtual Real | subset_l1_norm (const std::set< numeric_index_type > &indices) const |
virtual Real | subset_l2_norm (const std::set< numeric_index_type > &indices) const |
virtual Real | subset_linfty_norm (const std::set< numeric_index_type > &indices) const |
Real | l2_norm_diff (const NumericVector< T > &other_vec) const |
Real | l1_norm_diff (const NumericVector< T > &other_vec) const |
virtual T | el (const numeric_index_type i) const |
void | get (const std::vector< numeric_index_type > &index, std::vector< T > &values) const |
Access multiple components at once. More... | |
NumericVector< T > & | operator*= (const T a) |
Scales the vector by a , \( \vec{u} \leftarrow a\vec{u} \). More... | |
NumericVector< T > & | operator/= (const T a) |
Scales the vector by 1/a , \( \vec{u} \leftarrow \frac{1}{a}\vec{u} \). More... | |
void | add_vector (const std::vector< T > &v, const std::vector< numeric_index_type > &dof_indices) |
Computes \( \vec{u} \leftarrow \vec{u} + \vec{v} \), where v is a std::vector and each dof_indices [i] specifies where to add value v [i]. More... | |
void | add_vector (const NumericVector< T > &v, const std::vector< numeric_index_type > &dof_indices) |
Computes \( \vec{u} \leftarrow \vec{u} + \vec{v} \), where v is a NumericVector and each dof_indices [i] specifies where to add value v(i) . More... | |
void | add_vector (const DenseVector< T > &v, const std::vector< numeric_index_type > &dof_indices) |
Computes \( \vec{u} \leftarrow \vec{u} + \vec{v} \), where v is a DenseVector and each dof_indices [i] specifies where to add value v(i) . More... | |
void | add_vector (const NumericVector< T > &v, const ShellMatrix< T > &A) |
Computes \( \vec{u} \leftarrow \vec{u} + A \vec{v} \), i.e. More... | |
void | insert (const std::vector< T > &v, const std::vector< numeric_index_type > &dof_indices) |
Inserts the entries of v in *this at the locations specified by v . More... | |
void | insert (const NumericVector< T > &v, const std::vector< numeric_index_type > &dof_indices) |
Inserts the entries of v in *this at the locations specified by v . More... | |
void | insert (const DenseVector< T > &v, const std::vector< numeric_index_type > &dof_indices) |
Inserts the entries of v in *this at the locations specified by v . More... | |
void | insert (const DenseSubVector< T > &v, const std::vector< numeric_index_type > &dof_indices) |
Inserts the entries of v in *this at the locations specified by v . More... | |
virtual int | compare (const NumericVector< T > &other_vector, const Real threshold=TOLERANCE) const |
virtual int | local_relative_compare (const NumericVector< T > &other_vector, const Real threshold=TOLERANCE) const |
virtual int | global_relative_compare (const NumericVector< T > &other_vector, const Real threshold=TOLERANCE) const |
virtual void | print (std::ostream &os=libMesh::out) const |
Prints the local contents of the vector, by default to libMesh::out. More... | |
template<> | |
void | print (std::ostream &os) const |
virtual void | print_global (std::ostream &os=libMesh::out) const |
Prints the global contents of the vector, by default to libMesh::out. More... | |
template<> | |
void | print_global (std::ostream &os) const |
virtual void | read_matlab (const std::string &filename) |
Read the contents of the vector from the Matlab-script format used by PETSc. More... | |
bool | readable () const |
bool | compatible (const NumericVector< T > &v) const |
const Parallel::Communicator & | comm () const |
processor_id_type | n_processors () const |
processor_id_type | processor_id () const |
Static Public Member Functions | |
static std::unique_ptr< NumericVector< T > > | build (const Parallel::Communicator &comm, SolverPackage solver_package=libMesh::default_solver_package(), ParallelType parallel_type=AUTOMATIC) |
Builds a NumericVector on the processors in communicator comm using the linear solver package specified by solver_package . More... | |
static std::string | get_info () |
Gets a string containing the reference information. More... | |
static void | print_info (std::ostream &out_stream=libMesh::out) |
Prints the reference information, by default to libMesh::out . More... | |
static unsigned int | n_objects () |
Prints the number of outstanding (created, but not yet destroyed) objects. More... | |
static void | enable_print_counter_info () |
Methods to enable/disable the reference counter output from print_info() More... | |
static void | disable_print_counter_info () |
Protected Types | |
typedef std::map< std::string, std::pair< unsigned int, unsigned int > > | Counts |
Data structure to log the information. More... | |
Protected Member Functions | |
void | increment_constructor_count (const std::string &name) noexcept |
Increments the construction counter. More... | |
void | increment_destructor_count (const std::string &name) noexcept |
Increments the destruction counter. More... | |
Protected Attributes | |
bool | _is_closed |
Flag which tracks whether the vector's values are consistent on all processors after insertion or addition of values has occurred on some or all processors. More... | |
bool | _is_initialized |
true once init() has been called. More... | |
ParallelType | _type |
Type of vector. More... | |
std::mutex | _numeric_vector_mutex |
Mutex for performing thread-safe operations. More... | |
const Parallel::Communicator & | _communicator |
Static Protected Attributes | |
static Counts | _counts |
Actually holds the data. More... | |
static Threads::atomic< unsigned int > | _n_objects |
The number of objects. More... | |
static Threads::spin_mutex | _mutex |
Mutual exclusion object to enable thread-safe reference counting. More... | |
static bool | _enable_print_counter = true |
Flag to control whether reference count information is printed when print_info is called. More... | |
Private Types | |
typedef std::unordered_map< numeric_index_type, numeric_index_type > | GlobalToLocalMap |
Type for map that maps global to local ghost cells. More... | |
Private Member Functions | |
void | _get_array (bool read_only) const |
Queries the array (and the local form if the vector is ghosted) from PETSc. More... | |
void | _restore_array () const |
Restores the array (and the local form if the vector is ghosted) to PETSc. More... | |
template<NormType N> | |
Real | norm () const |
Private Attributes | |
Vec | _vec |
Actual PETSc vector datatype to hold vector entries. More... | |
std::atomic< bool > | _array_is_present |
If true , the actual PETSc array of the values of the vector is currently accessible. More... | |
bool | _array_is_present |
numeric_index_type | _first |
First local index. More... | |
numeric_index_type | _last |
Last local index. More... | |
numeric_index_type | _local_size |
Size of the local values from _get_array() More... | |
Vec | _local_form |
PETSc vector datatype to hold the local form of a ghosted vector. More... | |
const PetscScalar * | _read_only_values |
Pointer to the actual PETSc array of the values of the vector. More... | |
PetscScalar * | _values |
Pointer to the actual PETSc array of the values of the vector. More... | |
std::mutex | _petsc_get_restore_array_mutex |
Mutex for _get_array and _restore_array. More... | |
GlobalToLocalMap | _global_to_local_map |
Map that maps global to local ghost cells (will be empty if not in ghost cell mode). More... | |
bool | _destroy_vec_on_exit |
This boolean value should only be set to false for the constructor which takes a PETSc Vec object. More... | |
bool | _values_manually_retrieved |
Whether or not the data array has been manually retrieved using get_array() or get_array_read() More... | |
bool | _values_read_only |
Whether or not the data array is for read only access. More... | |
This class provides a nice interface to PETSc's Vec object.
All overridden virtual functions are documented in numeric_vector.h.
Definition at line 73 of file petsc_vector.h.
|
protectedinherited |
Data structure to log the information.
The log is identified by the class name.
Definition at line 119 of file reference_counter.h.
|
private |
Type for map that maps global to local ghost cells.
Definition at line 444 of file petsc_vector.h.
|
inlineexplicit |
Dummy-Constructor.
Dimension=0
Definition at line 485 of file petsc_vector.h.
|
inlineexplicit |
Constructor.
Set dimension to n
and initialize all elements with zero.
Definition at line 507 of file petsc_vector.h.
|
inline |
Constructor.
Set local dimension to n_local
, the global dimension to n
, and initialize all elements with zero.
Definition at line 531 of file petsc_vector.h.
|
inline |
Constructor.
Set local dimension to n_local
, the global dimension to n
, but additionally reserve memory for the indices specified by the ghost
argument.
Definition at line 556 of file petsc_vector.h.
|
inline |
Constructor.
Creates a PetscVector assuming you already have a valid PETSc Vec object. In this case, v
is NOT destroyed by the PetscVector constructor when this object goes out of scope. This allows ownership of v
to remain with the original creator, and to simply provide additional functionality with the PetscVector.
Definition at line 584 of file petsc_vector.h.
|
delete |
This class manages a C-style struct (Vec) manually, so we don't want to allow any automatic copy/move functions to be generated, and we can't default the destructor.
|
delete |
|
inlinevirtual |
Definition at line 682 of file petsc_vector.h.
|
private |
Queries the array (and the local form if the vector is ghosted) from PETSc.
read_only | Whether or not a read only copy of the raw data |
Definition at line 1182 of file petsc_vector.C.
|
private |
Restores the array (and the local form if the vector is ghosted) to PETSc.
Definition at line 1256 of file petsc_vector.C.
Referenced by libMesh::PetscVector< libMesh::Number >::add(), libMesh::PetscVector< libMesh::Number >::create_subvector(), libMesh::PetscVector< libMesh::Number >::init(), and libMesh::PetscVector< libMesh::Number >::operator=().
|
overridevirtual |
Sets \( u_i \leftarrow |u_i| \) for each entry in the vector.
Implements libMesh::NumericVector< T >.
Definition at line 411 of file petsc_vector.C.
|
overridevirtual |
Adds value
to the vector entry specified by i
.
Note that library implementations of this method are thread safe, e.g. we will lock _numeric_vector_mutex
before writing to the vector
Implements libMesh::NumericVector< T >.
Definition at line 165 of file petsc_vector.C.
Referenced by main().
|
overridevirtual |
Adds s
to each entry of the vector, \( u_i \leftarrow u_i + s \).
Implements libMesh::NumericVector< T >.
Definition at line 294 of file petsc_vector.C.
|
overridevirtual |
Adds v
to this
, \( \vec{u} \leftarrow \vec{u} + \vec{v} \).
Equivalent to calling operator+=()
.
Implements libMesh::NumericVector< T >.
Definition at line 305 of file petsc_vector.C.
|
overridevirtual |
Vector addition with a scalar multiple, \( \vec{u} \leftarrow \vec{u} + a\vec{v} \).
Equivalent to calling operator+=()
.
Implements libMesh::NumericVector< T >.
Definition at line 315 of file petsc_vector.C.
|
overridevirtual |
Computes \( \vec{u} \leftarrow \vec{u} + \vec{v} \), where v
is a pointer and each dof_indices
[i] specifies where to add value v
[i].
This should be overridden in subclasses for efficiency. Note that library implementations of this method are thread safe
Reimplemented from libMesh::NumericVector< T >.
Definition at line 182 of file petsc_vector.C.
|
overridevirtual |
Computes \( \vec{u} \leftarrow \vec{u} + A \vec{v} \), i.e.
adds the product of a SparseMatrix
A
and a NumericVector
v
to this
.
Implements libMesh::NumericVector< T >.
Definition at line 204 of file petsc_vector.C.
|
inlineinherited |
Computes \( \vec{u} \leftarrow \vec{u} + \vec{v} \), where v
is a std::vector and each dof_indices
[i] specifies where to add value v
[i].
This method is guaranteed to be thread safe as long as library implementations of NumericVector
are used
Definition at line 959 of file numeric_vector.h.
|
inherited |
Computes \( \vec{u} \leftarrow \vec{u} + \vec{v} \), where v
is a NumericVector and each dof_indices
[i] specifies where to add value v(i)
.
This method is guaranteed to be thread-safe as long as library implementations of NumericVector
are used
Definition at line 660 of file numeric_vector.C.
|
inlineinherited |
Computes \( \vec{u} \leftarrow \vec{u} + \vec{v} \), where v
is a DenseVector and each dof_indices
[i] specifies where to add value v(i)
.
This method is guaranteed to be thread safe as long as library implementations of NumericVector
are used
Definition at line 971 of file numeric_vector.h.
|
inherited |
Computes \( \vec{u} \leftarrow \vec{u} + A \vec{v} \), i.e.
adds the product of a ShellMatrix
A
and a NumericVector
v
to this
.
Definition at line 674 of file numeric_vector.C.
void libMesh::PetscVector< T >::add_vector_conjugate_transpose | ( | const NumericVector< T > & | v, |
const SparseMatrix< T > & | A | ||
) |
\( U \leftarrow U + A^H v \).
Adds the product of the conjugate-transpose of SparseMatrix
A
and a NumericVector
v
to this
.
Definition at line 260 of file petsc_vector.C.
|
overridevirtual |
Computes \( \vec{u} \leftarrow \vec{u} + A^T \vec{v} \), i.e.
adds the product of the transpose of a SparseMatrix
A
and a NumericVector
v
to this
.
Implements libMesh::NumericVector< T >.
Definition at line 232 of file petsc_vector.C.
|
staticinherited |
Builds a NumericVector
on the processors in communicator comm
using the linear solver package specified by solver_package
.
Definition at line 61 of file numeric_vector.C.
Referenced by libMesh::ExactSolution::_compute_error(), libMesh::UniformRefinementEstimator::_estimate_error(), libMesh::System::add_vector(), libMesh::TransientRBConstruction::allocate_data_structures(), libMesh::RBConstruction::allocate_data_structures(), libMesh::TransientRBConstruction::assemble_affine_expansion(), libMesh::AdvectionSystem::assemble_claw_rhs(), libMesh::EquationSystems::build_parallel_elemental_solution_vector(), libMesh::EquationSystems::build_parallel_solution_vector(), libMesh::System::calculate_norm(), libMesh::RBConstruction::compute_Fq_representor_innerprods(), libMesh::RBConstruction::compute_output_dual_innerprods(), libMesh::RBConstruction::compute_residual_dual_norm_slow(), libMesh::DiagonalMatrix< T >::DiagonalMatrix(), libMesh::AdjointRefinementEstimator::estimate_error(), libMesh::ExactErrorEstimator::estimate_error(), libMesh::CondensedEigenSystem::get_eigenpair(), libMesh::StaticCondensation::init(), main(), libMesh::TransientRBConstruction::mass_matrix_scaled_matvec(), libMesh::DofMap::process_mesh_constraint_rows(), libMesh::InterMeshProjection::project_system_vectors(), libMesh::RBEvaluation::read_in_vectors_from_multiple_files(), libMesh::TransientRBConstruction::read_riesz_representors_from_files(), libMesh::RBConstruction::read_riesz_representors_from_files(), OverlappingAlgebraicGhostingTest::run_ghosting_test(), OverlappingCouplingGhostingTest::run_sparsity_pattern_test(), libMesh::TransientRBConstruction::set_error_temporal_data(), libMesh::ClawSystem::solve_conservation_law(), ConstraintOperatorTest::test1DCoarseningOperator(), MeshfunctionDFEM::test_mesh_function_dfem(), MeshfunctionDFEM::test_mesh_function_dfem_grad(), MeshFunctionTest::test_p_level(), DiagonalMatrixTest::testNumerics(), SystemsTest::testProjectCubeWithMeshFunction(), libMesh::RBConstruction::train_reduced_basis_with_POD(), libMesh::MeshFunctionSolutionTransfer::transfer(), libMesh::TransientRBConstruction::truth_assembly(), libMesh::RBConstruction::truth_assembly(), libMesh::TransientRBConstruction::update_RB_initial_condition_all_N(), libMesh::TransientRBConstruction::update_RB_system_matrices(), libMesh::RBConstruction::update_RB_system_matrices(), libMesh::TransientRBConstruction::update_residual_terms(), and libMesh::RBConstruction::update_residual_terms().
|
inlineoverridevirtualnoexcept |
clear() is called from the destructor, so it should not throw.
Reimplemented from libMesh::NumericVector< T >.
Definition at line 887 of file petsc_vector.h.
|
inlineoverridevirtual |
Implements libMesh::NumericVector< T >.
Definition at line 951 of file petsc_vector.h.
Referenced by PetscVectorTest::testPetscOperations().
|
inlineoverridevirtual |
Calls the NumericVector's internal assembly routines, ensuring that the values are consistent across processors.
Implements libMesh::NumericVector< T >.
Definition at line 869 of file petsc_vector.h.
Referenced by libMesh::__libmesh_petsc_diff_solver_monitor(), libMesh::SlepcEigenSolver< libMesh::Number >::get_eigenpair(), libMesh::PetscVector< libMesh::Number >::localize(), libMesh::PetscLinearSolver< Number >::solve_base(), PetscVectorTest::testGetArray(), and PetscVectorTest::testPetscOperations().
|
inlinevirtualinherited |
true
if the vector is closed and ready for computation, false otherwise. Definition at line 182 of file numeric_vector.h.
Referenced by libMesh::System::add_vector(), libMesh::PetscVector< libMesh::Number >::localize(), libMesh::DofMap::max_constraint_error(), libMesh::EigenSparseVector< T >::operator=(), libMesh::LaspackVector< T >::operator=(), and libMesh::PetscVector< libMesh::Number >::operator=().
|
inlineinherited |
Parallel::Communicator
object used by this mesh. Definition at line 97 of file parallel_object.h.
References libMesh::ParallelObject::_communicator.
Referenced by libMesh::__libmesh_petsc_diff_solver_jacobian(), libMesh::__libmesh_petsc_diff_solver_monitor(), libMesh::__libmesh_petsc_diff_solver_residual(), libMesh::__libmesh_tao_equality_constraints(), libMesh::__libmesh_tao_equality_constraints_jacobian(), libMesh::__libmesh_tao_gradient(), libMesh::__libmesh_tao_hessian(), libMesh::__libmesh_tao_inequality_constraints(), libMesh::__libmesh_tao_inequality_constraints_jacobian(), libMesh::__libmesh_tao_objective(), libMesh::MeshRefinement::_coarsen_elements(), libMesh::ExactSolution::_compute_error(), libMesh::UniformRefinementEstimator::_estimate_error(), libMesh::Partitioner::_find_global_index_by_pid_map(), libMesh::BoundaryInfo::_find_id_maps(), libMesh::PetscLinearSolver< Number >::_petsc_shell_matrix_get_diagonal(), libMesh::SlepcEigenSolver< libMesh::Number >::_petsc_shell_matrix_get_diagonal(), libMesh::PetscLinearSolver< Number >::_petsc_shell_matrix_mult(), libMesh::SlepcEigenSolver< libMesh::Number >::_petsc_shell_matrix_mult(), libMesh::PetscLinearSolver< Number >::_petsc_shell_matrix_mult_add(), libMesh::MeshRefinement::_refine_elements(), libMesh::MeshRefinement::_smooth_flags(), libMesh::DofMap::add_constraints_to_send_list(), add_cube_convex_hull_to_mesh(), libMesh::PetscDMWrapper::add_dofs_helper(), libMesh::PetscDMWrapper::add_dofs_to_section(), libMesh::TransientRBConstruction::add_IC_to_RB_space(), libMesh::RBEIMEvaluation::add_interpolation_data(), libMesh::CondensedEigenSystem::add_matrices(), libMesh::EigenSystem::add_matrices(), libMesh::System::add_matrix(), libMesh::RBConstruction::add_scaled_matrix_and_vector(), libMesh::System::add_variable(), libMesh::System::add_variables(), libMesh::System::add_vector(), libMesh::MeshTools::Modification::all_tri(), libMesh::LaplaceMeshSmoother::allgather_graph(), libMesh::DofMap::allgather_recursive_constraints(), libMesh::TransientRBConstruction::allocate_data_structures(), libMesh::RBConstruction::allocate_data_structures(), libMesh::TransientRBConstruction::assemble_affine_expansion(), libMesh::AdvectionSystem::assemble_claw_rhs(), libMesh::FEMSystem::assemble_qoi(), libMesh::Nemesis_IO::assert_symmetric_cmaps(), libMesh::MeshCommunication::assign_global_indices(), libMesh::Partitioner::assign_partitioning(), libMesh::MeshTools::Generation::build_extrusion(), libMesh::Partitioner::build_graph(), libMesh::BoundaryInfo::build_node_list_from_side_list(), libMesh::EquationSystems::build_parallel_elemental_solution_vector(), libMesh::EquationSystems::build_parallel_solution_vector(), libMesh::PetscDMWrapper::build_section(), libMesh::PetscDMWrapper::build_sf(), libMesh::MeshBase::cache_elem_data(), libMesh::System::calculate_norm(), libMesh::DofMap::check_dirichlet_bcid_consistency(), libMesh::RBConstruction::compute_Fq_representor_innerprods(), libMesh::RBConstruction::compute_max_error_bound(), libMesh::Nemesis_IO_Helper::compute_num_global_elem_blocks(), libMesh::Nemesis_IO_Helper::compute_num_global_nodesets(), libMesh::Nemesis_IO_Helper::compute_num_global_sidesets(), libMesh::RBConstruction::compute_output_dual_innerprods(), libMesh::RBConstruction::compute_residual_dual_norm_slow(), libMesh::RBSCMConstruction::compute_SCM_bounds_on_training_set(), libMesh::DofMap::computed_sparsity_already(), libMesh::Problem_Interface::computeF(), libMesh::Problem_Interface::computeJacobian(), libMesh::Problem_Interface::computePreconditioner(), libMesh::ContinuationSystem::ContinuationSystem(), libMesh::MeshBase::copy_constraint_rows(), libMesh::ExodusII_IO::copy_elemental_solution(), libMesh::ExodusII_IO::copy_nodal_solution(), libMesh::ExodusII_IO::copy_scalar_solution(), libMesh::CondensedEigenSystem::copy_super_to_sub(), libMesh::MeshTools::correct_node_proc_ids(), libMesh::MeshTools::create_bounding_box(), libMesh::DofMap::create_dof_constraints(), libMesh::MeshTools::create_nodal_bounding_box(), libMesh::MeshRefinement::create_parent_error_vector(), libMesh::MeshTools::create_processor_bounding_box(), libMesh::MeshTools::create_subdomain_bounding_box(), libMesh::PetscMatrix< T >::create_submatrix_nosort(), create_wrapped_function(), libMesh::MeshCommunication::delete_remote_elements(), libMesh::RBEIMEvaluation::distribute_bfs(), DMlibMeshFunction(), DMlibMeshJacobian(), DMlibMeshSetSystem_libMesh(), DMVariableBounds_libMesh(), libMesh::DTKSolutionTransfer::DTKSolutionTransfer(), libMesh::MeshRefinement::eliminate_unrefined_patches(), libMesh::RBEIMConstruction::enrich_eim_approximation_on_interiors(), libMesh::RBEIMConstruction::enrich_eim_approximation_on_nodes(), libMesh::RBEIMConstruction::enrich_eim_approximation_on_sides(), libMesh::TransientRBConstruction::enrich_RB_space(), libMesh::EpetraVector< T >::EpetraVector(), AssembleOptimization::equality_constraints(), libMesh::PatchRecoveryErrorEstimator::estimate_error(), libMesh::WeightedPatchRecoveryErrorEstimator::estimate_error(), libMesh::AdjointRefinementEstimator::estimate_error(), libMesh::ExactErrorEstimator::estimate_error(), libMesh::MeshRefinement::flag_elements_by_elem_fraction(), libMesh::MeshRefinement::flag_elements_by_error_fraction(), libMesh::MeshRefinement::flag_elements_by_error_tolerance(), libMesh::MeshRefinement::flag_elements_by_mean_stddev(), libMesh::MeshRefinement::flag_elements_by_nelem_target(), libMesh::RBEIMEvaluation::gather_bfs(), libMesh::DofMap::gather_constraints(), libMesh::MeshfreeInterpolation::gather_remote_data(), libMesh::CondensedEigenSystem::get_eigenpair(), libMesh::RBEIMEvaluation::get_eim_basis_function_node_value(), libMesh::RBEIMEvaluation::get_eim_basis_function_side_value(), libMesh::RBEIMEvaluation::get_eim_basis_function_value(), libMesh::MeshBase::get_info(), libMesh::System::get_info(), libMesh::DofMap::get_info(), libMesh::RBEIMEvaluation::get_interior_basis_functions_as_vecs(), libMesh::ImplicitSystem::get_linear_solver(), libMesh::RBEIMConstruction::get_max_abs_value(), libMesh::RBEIMConstruction::get_node_max_abs_value(), libMesh::RBEIMEvaluation::get_parametrized_function_node_value(), libMesh::RBEIMEvaluation::get_parametrized_function_side_value(), libMesh::RBEIMEvaluation::get_parametrized_function_value(), libMesh::RBEIMConstruction::get_random_point(), AssembleOptimization::inequality_constraints(), AssembleOptimization::inequality_constraints_jacobian(), libMesh::LocationMap< T >::init(), libMesh::TimeSolver::init(), libMesh::StaticCondensation::init(), libMesh::SystemSubsetBySubdomain::init(), libMesh::PetscDMWrapper::init_and_attach_petscdm(), libMesh::AdvectionSystem::init_data(), libMesh::ClawSystem::init_data(), libMesh::PetscDMWrapper::init_petscdm(), libMesh::ExodusII_IO_Helper::initialize(), libMesh::OptimizationSystem::initialize_equality_constraints_storage(), libMesh::OptimizationSystem::initialize_inequality_constraints_storage(), libMesh::RBEIMConstruction::initialize_parametrized_functions_in_training_set(), libMesh::RBEIMConstruction::inner_product(), integrate_function(), libMesh::MeshTools::libmesh_assert_consistent_distributed(), libMesh::MeshTools::libmesh_assert_consistent_distributed_nodes(), libMesh::MeshTools::libmesh_assert_contiguous_dof_ids(), libMesh::MeshTools::libmesh_assert_equal_connectivity(), libMesh::MeshTools::libmesh_assert_equal_points(), libMesh::MeshTools::libmesh_assert_parallel_consistent_new_node_procids(), libMesh::MeshTools::libmesh_assert_parallel_consistent_procids< Elem >(), libMesh::MeshTools::libmesh_assert_parallel_consistent_procids< Node >(), libMesh::MeshTools::libmesh_assert_topology_consistent_procids< Node >(), libMesh::MeshTools::libmesh_assert_valid_boundary_ids(), libMesh::MeshTools::libmesh_assert_valid_constraint_rows(), libMesh::MeshTools::libmesh_assert_valid_dof_ids(), libMesh::MeshTools::libmesh_assert_valid_neighbors(), libMesh::DistributedMesh::libmesh_assert_valid_parallel_flags(), libMesh::DistributedMesh::libmesh_assert_valid_parallel_object_ids(), libMesh::DistributedMesh::libmesh_assert_valid_parallel_p_levels(), libMesh::MeshTools::libmesh_assert_valid_refinement_flags(), libMesh::MeshTools::libmesh_assert_valid_unique_ids(), libMesh::libmesh_petsc_linesearch_shellfunc(), libMesh::libmesh_petsc_preconditioner_apply(), libMesh::libmesh_petsc_recalculate_monitor(), libMesh::libmesh_petsc_snes_fd_residual(), libMesh::libmesh_petsc_snes_jacobian(), libMesh::libmesh_petsc_snes_mffd_interface(), libMesh::libmesh_petsc_snes_mffd_residual(), libMesh::libmesh_petsc_snes_postcheck(), libMesh::libmesh_petsc_snes_precheck(), libMesh::libmesh_petsc_snes_residual(), libMesh::libmesh_petsc_snes_residual_helper(), libMesh::MeshRefinement::limit_level_mismatch_at_edge(), libMesh::MeshRefinement::limit_level_mismatch_at_node(), libMesh::MeshRefinement::limit_overrefined_boundary(), libMesh::MeshRefinement::limit_underrefined_boundary(), libMesh::LinearImplicitSystem::LinearImplicitSystem(), main(), libMesh::MeshRefinement::make_coarsening_compatible(), libMesh::MeshCommunication::make_elems_parallel_consistent(), libMesh::MeshRefinement::make_flags_parallel_consistent(), libMesh::MeshCommunication::make_new_node_proc_ids_parallel_consistent(), libMesh::MeshCommunication::make_new_nodes_parallel_consistent(), libMesh::MeshCommunication::make_node_bcids_parallel_consistent(), libMesh::MeshCommunication::make_node_ids_parallel_consistent(), libMesh::MeshCommunication::make_node_proc_ids_parallel_consistent(), libMesh::MeshCommunication::make_node_unique_ids_parallel_consistent(), libMesh::MeshCommunication::make_nodes_parallel_consistent(), libMesh::MeshCommunication::make_p_levels_parallel_consistent(), libMesh::MeshRefinement::make_refinement_compatible(), libMesh::TransientRBConstruction::mass_matrix_scaled_matvec(), libMesh::FEMSystem::mesh_position_set(), libMesh::TriangulatorInterface::MeshedHole::MeshedHole(), LinearElasticityWithContact::move_mesh(), libMesh::DistributedMesh::n_active_elem(), libMesh::MeshTools::n_active_levels(), libMesh::BoundaryInfo::n_boundary_conds(), libMesh::MeshTools::n_connected_components(), libMesh::DofMap::n_constrained_dofs(), libMesh::MeshBase::n_constraint_rows(), libMesh::DofMap::n_dofs(), libMesh::DofMap::n_dofs_per_processor(), libMesh::BoundaryInfo::n_edge_conds(), libMesh::CondensedEigenSystem::n_global_non_condensed_dofs(), libMesh::MeshTools::n_levels(), MixedOrderTest::n_neighbor_links(), libMesh::BoundaryInfo::n_nodeset_conds(), libMesh::SparsityPattern::Build::n_nonzeros(), libMesh::MeshTools::n_p_levels(), libMesh::BoundaryInfo::n_shellface_conds(), libMesh::RBEIMEvaluation::node_distribute_bfs(), libMesh::RBEIMEvaluation::node_gather_bfs(), libMesh::RBEIMConstruction::node_inner_product(), libMesh::MeshBase::operator==(), libMesh::DistributedMesh::parallel_max_elem_id(), libMesh::DistributedMesh::parallel_max_node_id(), libMesh::ReplicatedMesh::parallel_max_unique_id(), libMesh::DistributedMesh::parallel_max_unique_id(), libMesh::DistributedMesh::parallel_n_elem(), libMesh::DistributedMesh::parallel_n_nodes(), libMesh::SparsityPattern::Build::parallel_sync(), libMesh::BoundaryInfo::parallel_sync_node_ids(), libMesh::BoundaryInfo::parallel_sync_side_ids(), libMesh::MeshTools::paranoid_n_levels(), libMesh::Partitioner::partition(), libMesh::Partitioner::partition_unpartitioned_elements(), libMesh::petsc_auto_fieldsplit(), libMesh::System::point_gradient(), libMesh::System::point_hessian(), libMesh::System::point_value(), libMesh::MeshBase::prepare_for_use(), libMesh::MeshBase::print_constraint_rows(), libMesh::DofMap::print_dof_constraints(), libMesh::DofMap::process_mesh_constraint_rows(), libMesh::Partitioner::processor_pairs_to_interface_nodes(), libMesh::InterMeshProjection::project_system_vectors(), FEMParameters::read(), libMesh::Nemesis_IO::read(), libMesh::XdrIO::read(), libMesh::EquationSystems::read(), libMesh::ExodusII_IO::read_header(), libMesh::CheckpointIO::read_header(), libMesh::XdrIO::read_header(), libMesh::System::read_header(), libMesh::RBEIMEvaluation::read_in_interior_basis_functions(), libMesh::RBEIMEvaluation::read_in_node_basis_functions(), libMesh::RBEIMEvaluation::read_in_side_basis_functions(), libMesh::RBEvaluation::read_in_vectors_from_multiple_files(), libMesh::System::read_legacy_data(), libMesh::TransientRBConstruction::read_riesz_representors_from_files(), libMesh::RBConstruction::read_riesz_representors_from_files(), libMesh::System::read_SCALAR_dofs(), libMesh::XdrIO::read_serialized_bc_names(), libMesh::XdrIO::read_serialized_bcs_helper(), libMesh::System::read_serialized_blocked_dof_objects(), libMesh::XdrIO::read_serialized_connectivity(), libMesh::XdrIO::read_serialized_nodes(), libMesh::XdrIO::read_serialized_nodesets(), libMesh::XdrIO::read_serialized_subdomain_names(), libMesh::System::read_serialized_vector(), libMesh::Nemesis_IO_Helper::read_var_names_impl(), libMesh::MeshBase::recalculate_n_partitions(), libMesh::MeshRefinement::refine_and_coarsen_elements(), libMesh::SimplexRefiner::refine_via_edges(), libMesh::StaticCondensationDofMap::reinit(), libMesh::DistributedMesh::renumber_dof_objects(), libMesh::DistributedMesh::renumber_nodes_and_elements(), LinearElasticityWithContact::residual_and_jacobian(), OverlappingAlgebraicGhostingTest::run_ghosting_test(), OverlappingCouplingGhostingTest::run_sparsity_pattern_test(), scale_mesh_and_plot(), libMesh::DofMap::scatter_constraints(), libMesh::CheckpointIO::select_split_config(), libMesh::GenericProjector< FFunctor, GFunctor, FValue, ProjectionAction >::send_and_insert_dof_values(), libMesh::TransientRBConstruction::set_error_temporal_data(), libMesh::Partitioner::set_interface_node_processor_ids_BFS(), libMesh::Partitioner::set_interface_node_processor_ids_linear(), libMesh::Partitioner::set_interface_node_processor_ids_petscpartitioner(), libMesh::Partitioner::set_node_processor_ids(), libMesh::DofMap::set_nonlocal_dof_objects(), libMesh::Partitioner::set_parent_processor_ids(), libMesh::PetscDMWrapper::set_point_range_in_section(), libMesh::PetscDiffSolver::setup_petsc_data(), libMesh::RBEIMEvaluation::side_distribute_bfs(), libMesh::RBEIMEvaluation::side_gather_bfs(), libMesh::RBEIMConstruction::side_inner_product(), libMesh::Partitioner::single_partition(), libMesh::LaplaceMeshSmoother::smooth(), libMesh::VariationalMeshSmoother::smooth(), libMesh::ClawSystem::solve_conservation_law(), libMesh::split_mesh(), libMesh::RBEIMConstruction::store_eim_solutions_for_training_set(), libMesh::MeshBase::subdomain_ids(), libMesh::BoundaryInfo::sync(), ConstraintOperatorTest::test1DCoarseningNewNodes(), ConstraintOperatorTest::test1DCoarseningOperator(), libMesh::MeshRefinement::test_level_one(), MeshfunctionDFEM::test_mesh_function_dfem(), MeshfunctionDFEM::test_mesh_function_dfem_grad(), MeshFunctionTest::test_p_level(), libMesh::MeshRefinement::test_unflagged(), DofMapTest::testBadElemFECombo(), SystemsTest::testBlockRestrictedVarNDofs(), BoundaryInfoTest::testBoundaryOnChildrenErrors(), VolumeTest::testC0PolygonMethods(), VolumeTest::testC0PolyhedronMethods(), ConstraintOperatorTest::testCoreform(), ConnectedComponentsTest::testEdge(), MeshInputTest::testExodusIGASidesets(), MeshTriangulationTest::testFoundCenters(), PointLocatorTest::testLocator(), BoundaryInfoTest::testMesh(), PointLocatorTest::testPlanar(), MeshTriangulationTest::testPoly2TriRefinementBase(), SystemsTest::testProjectCubeWithMeshFunction(), BoundaryInfoTest::testRenumber(), CheckpointIOTest::testSplitter(), MeshInputTest::testTetgenIO(), MeshTriangulationTest::testTriangulatorInterp(), MeshTriangulationTest::testTriangulatorMeshedHoles(), MeshTriangulationTest::testTriangulatorRoundHole(), libMesh::MeshTools::total_weight(), libMesh::RBConstruction::train_reduced_basis_with_POD(), libMesh::MeshFunctionSolutionTransfer::transfer(), libMesh::MeshfreeSolutionTransfer::transfer(), libMesh::Poly2TriTriangulator::triangulate(), libMesh::TransientRBConstruction::truth_assembly(), libMesh::RBConstruction::truth_assembly(), libMesh::MeshRefinement::uniformly_coarsen(), update_current_local_solution(), libMesh::TransientRBConstruction::update_RB_initial_condition_all_N(), libMesh::TransientRBConstruction::update_RB_system_matrices(), libMesh::RBConstruction::update_RB_system_matrices(), libMesh::TransientRBConstruction::update_residual_terms(), libMesh::RBConstruction::update_residual_terms(), libMesh::MeshTools::volume(), libMesh::STLIO::write(), libMesh::NameBasedIO::write(), libMesh::XdrIO::write(), libMesh::VTKIO::write_nodal_data(), libMesh::RBEIMEvaluation::write_out_interior_basis_functions(), libMesh::RBEIMEvaluation::write_out_node_basis_functions(), libMesh::RBEIMEvaluation::write_out_side_basis_functions(), libMesh::RBEvaluation::write_out_vectors(), libMesh::TransientRBConstruction::write_riesz_representors_to_files(), libMesh::RBConstruction::write_riesz_representors_to_files(), libMesh::System::write_SCALAR_dofs(), libMesh::XdrIO::write_serialized_bcs_helper(), libMesh::System::write_serialized_blocked_dof_objects(), libMesh::XdrIO::write_serialized_connectivity(), libMesh::XdrIO::write_serialized_nodes(), libMesh::XdrIO::write_serialized_nodesets(), libMesh::RBDataSerialization::RBEvaluationSerialization::write_to_file(), libMesh::RBDataSerialization::TransientRBEvaluationSerialization::write_to_file(), libMesh::RBDataSerialization::RBEIMEvaluationSerialization::write_to_file(), and libMesh::RBDataSerialization::RBSCMEvaluationSerialization::write_to_file().
|
virtualinherited |
-1
when this
is equivalent to other_vector
(up to the given threshold
), or the first index where abs
(a[i]-b[i]) exceeds the threshold. Definition at line 153 of file numeric_vector.C.
|
inherited |
v
are able to be read for global operations with one-another Definition at line 692 of file numeric_vector.C.
|
overridevirtual |
Negates the imaginary component of each entry in the vector.
Implements libMesh::NumericVector< T >.
Definition at line 153 of file petsc_vector.C.
|
overridevirtual |
Fills in subvector
from this vector using the indices in rows
.
Similar to the create_submatrix()
routine for the SparseMatrix class, it is currently only implemented for PetscVectors. The boolean parameter communicates whether the supplied vector of rows corresponds to all the rows that should be used in the subvector's index set, e.g. whether the rows correspond to the global collective. If the rows supplied are only the local indices, then the boolean parameter should be set to false
Reimplemented from libMesh::NumericVector< T >.
Definition at line 1086 of file petsc_vector.C.
|
staticinherited |
Definition at line 100 of file reference_counter.C.
References libMesh::ReferenceCounter::_enable_print_counter.
|
overridevirtual |
v
.Uses the complex-conjugate of v
in the complex-valued case.
Implements libMesh::NumericVector< T >.
Definition at line 426 of file petsc_vector.C.
|
inlinevirtualinherited |
|
staticinherited |
Methods to enable/disable the reference counter output from print_info()
Definition at line 94 of file reference_counter.C.
References libMesh::ReferenceCounter::_enable_print_counter.
|
inlineoverridevirtual |
0
. Implements libMesh::NumericVector< T >.
Definition at line 997 of file petsc_vector.h.
Referenced by PetscVectorTest::testPetscOperations().
|
inlineoverridevirtual |
Access multiple components at once.
values
will not be reallocated; it should already have enough space. The default implementation calls operator()
for each index, but some implementations may supply faster methods here.
Reimplemented from libMesh::NumericVector< T >.
Definition at line 1119 of file petsc_vector.h.
|
inlineinherited |
Access multiple components at once.
values
will be resized, if necessary, and filled. The default implementation calls operator()
for each index, but some implementations may supply faster methods here.
Definition at line 944 of file numeric_vector.h.
|
inline |
Get read/write access to the raw PETSc Vector data array.
Definition at line 1142 of file petsc_vector.h.
Referenced by PetscVectorTest::testGetArray().
|
inline |
Get read only access to the raw PETSc Vector data array.
Definition at line 1153 of file petsc_vector.h.
Referenced by PetscVectorTest::testGetArray().
|
staticinherited |
Gets a string containing the reference information.
Definition at line 47 of file reference_counter.C.
References libMesh::ReferenceCounter::_counts, and libMesh::Quality::name().
Referenced by libMesh::ReferenceCounter::print_info().
|
overridevirtual |
Creates a view into this vector using the indices in rows
.
Similar to the create_submatrix()
routine for the SparseMatrix class, it is currently only implemented for PetscVectors.
Reimplemented from libMesh::NumericVector< T >.
Definition at line 1295 of file petsc_vector.C.
|
virtualinherited |
-1
when this
is equivalent to other_vector
(up to the given global relative threshold
), or the first index where abs
(a[i]-b[i])/max_j(a[j],b[j]) exceeds the threshold. Definition at line 211 of file numeric_vector.C.
|
inlineprotectednoexceptinherited |
Increments the construction counter.
Should be called in the constructor of any derived class that will be reference counted.
Definition at line 183 of file reference_counter.h.
References libMesh::err, libMesh::BasicOStreamProxy< charT, traits >::get(), libMesh::Quality::name(), and libMesh::Threads::spin_mtx.
Referenced by libMesh::ReferenceCountedObject< RBParametrized >::ReferenceCountedObject().
|
inlineprotectednoexceptinherited |
Increments the destruction counter.
Should be called in the destructor of any derived class that will be reference counted.
Definition at line 207 of file reference_counter.h.
References libMesh::err, libMesh::BasicOStreamProxy< charT, traits >::get(), libMesh::Quality::name(), and libMesh::Threads::spin_mtx.
Referenced by libMesh::ReferenceCountedObject< RBParametrized >::~ReferenceCountedObject().
T libMesh::PetscVector< T >::indefinite_dot | ( | const NumericVector< T > & | v | ) | const |
v
.Definition at line 447 of file petsc_vector.C.
|
inlineoverridevirtual |
Change the dimension of the vector to n
.
The reserved memory for this vector remains unchanged if possible. If n==0
, all memory is freed. Therefore, if you want to resize the vector and release the memory not needed, you have to first call init(0)
and then init(n)
. This behaviour is analogous to that of the STL containers.
On fast==false
, the vector is filled by zeros.
Implements libMesh::NumericVector< T >.
Definition at line 691 of file petsc_vector.h.
Referenced by libMesh::PetscVector< libMesh::Number >::clone(), libMesh::PetscVector< libMesh::Number >::localize(), and libMesh::PetscVector< libMesh::Number >::zero_clone().
|
inlineoverridevirtual |
Call init()
with n_local = N.
Implements libMesh::NumericVector< T >.
Definition at line 755 of file petsc_vector.h.
|
overridevirtual |
Create a vector that holds the local indices plus those specified in the ghost
argument.
Implements libMesh::NumericVector< T >.
|
inlineoverridevirtual |
Creates a vector that has the same dimension and storage type as other
, including ghost dofs.
Implements libMesh::NumericVector< T >.
Definition at line 831 of file petsc_vector.h.
|
inlinevirtualinherited |
true
if the vector has been initialized, false otherwise. Definition at line 149 of file numeric_vector.h.
Referenced by libMesh::System::add_vector(), libMesh::ImplicitSystem::assemble(), libMesh::PetscVector< libMesh::Number >::create_subvector(), and libMesh::PetscVector< libMesh::Number >::init().
|
overridevirtual |
Inserts the entries of v
in *this at the locations specified by v
.
Note that library implementations of this method are thread safe
Reimplemented from libMesh::NumericVector< T >.
Definition at line 349 of file petsc_vector.C.
|
inlineinherited |
Inserts the entries of v
in *this at the locations specified by v
.
This method is guaranteed to be thread-safe as long as library implementations of NumericVector
are used
Definition at line 983 of file numeric_vector.h.
|
inherited |
Inserts the entries of v
in *this at the locations specified by v
.
This method is guaranteed to be thread-safe as long as library implementations of NumericVector
are used
Definition at line 140 of file numeric_vector.C.
|
inlineinherited |
Inserts the entries of v
in *this at the locations specified by v
.
This method is guaranteed to be thread-safe as long as library implementations of NumericVector
are used
Definition at line 995 of file numeric_vector.h.
|
inlineinherited |
Inserts the entries of v
in *this at the locations specified by v
.
This method is guaranteed to be thread-safe as long as library implementations of NumericVector
are used
Definition at line 1007 of file numeric_vector.h.
|
overridevirtual |
Implements libMesh::NumericVector< T >.
Definition at line 74 of file petsc_vector.C.
|
inherited |
this
. Definition at line 632 of file numeric_vector.C.
|
overridevirtual |
Implements libMesh::NumericVector< T >.
Definition at line 79 of file petsc_vector.C.
|
inherited |
this
. Definition at line 616 of file numeric_vector.C.
|
inlineoverridevirtual |
size()
. Implements libMesh::NumericVector< T >.
Definition at line 1019 of file petsc_vector.h.
Referenced by PetscVectorTest::testPetscOperations().
|
overridevirtual |
Implements libMesh::NumericVector< T >.
Definition at line 84 of file petsc_vector.C.
|
virtualinherited |
-1
when this
is equivalent to other_vector
, (up to the given local relative threshold
), or the first index where abs
(a[i]-b[i])/max(a[i],b[i]) exceeds the threshold. Definition at line 181 of file numeric_vector.C.
|
inlineoverridevirtual |
index_stop
- index_start
. Implements libMesh::NumericVector< T >.
Definition at line 982 of file petsc_vector.h.
Referenced by libMesh::PetscVector< libMesh::Number >::localize(), and libMesh::PetscVector< libMesh::Number >::operator=().
|
overridevirtual |
Creates a copy of the global vector in the local vector v_local
.
Implements libMesh::NumericVector< T >.
Definition at line 787 of file petsc_vector.C.
Referenced by libMesh::PetscVector< libMesh::Number >::localize().
|
overridevirtual |
Same, but fills a NumericVector<T>
instead of a std::vector
.
Implements libMesh::NumericVector< T >.
Definition at line 578 of file petsc_vector.C.
|
overridevirtual |
Creates a local vector v_local
containing only information relevant to this processor, as defined by the send_list
.
Implements libMesh::NumericVector< T >.
Definition at line 620 of file petsc_vector.C.
|
overridevirtual |
Fill in the local std::vector "v_local" with the global indices given in "indices".
Example:
* On 4 procs *this = {a, b, c, d, e, f, g, h, i} is a parallel vector. * On each proc, the indices arrays are set up as: * proc0, indices = {1,2,4,5} * proc1, indices = {2,5,6,8} * proc2, indices = {2,3,6,7} * proc3, indices = {0,1,2,3} * * After calling this version of localize, the v_local vectors are: * proc0, v_local = {b,c,e,f} * proc1, v_local = {c,f,g,i} * proc2, v_local = {c,d,g,h} * proc3, v_local = {a,b,c,d} *
This function is useful in parallel I/O routines, when you have a parallel vector of solution values which you want to write a subset of.
Implements libMesh::NumericVector< T >.
Definition at line 680 of file petsc_vector.C.
|
overridevirtual |
Updates a local vector with selected values from neighboring processors, as defined by send_list
.
Implements libMesh::NumericVector< T >.
Definition at line 729 of file petsc_vector.C.
|
overridevirtual |
Creates a local copy of the global vector in v_local
only on processor proc_id
.
By default the data is sent to processor 0. This method is useful for outputting data from one processor.
Implements libMesh::NumericVector< T >.
void libMesh::PetscVector< Real >::localize_to_one | ( | std::vector< Real > & | v_local, |
const processor_id_type | timpi_mpi_varpid | ||
) | const |
Definition at line 821 of file petsc_vector.C.
void libMesh::PetscVector< Complex >::localize_to_one | ( | std::vector< Complex > & | v_local, |
const processor_id_type | pid | ||
) | const |
Definition at line 904 of file petsc_vector.C.
|
inline |
i
.If the index is not a ghost cell, this is done by subtraction the number of the first local index. If it is a ghost cell, it has to be looked up in the map.
Definition at line 1041 of file petsc_vector.h.
|
inlineoverridevirtual |
Implements libMesh::NumericVector< T >.
Definition at line 1192 of file petsc_vector.h.
|
inlineoverridevirtual |
Although we define numeric_index_type, and it is typically required that the NumericVector's underlying implementation's indices have size equal to sizeof(numeric_index_type), these two types can still have different signedness, which affects the maximum number of values which can be stored in the NumericVector.
Implements libMesh::NumericVector< T >.
Definition at line 1238 of file petsc_vector.h.
|
inlineoverridevirtual |
Implements libMesh::NumericVector< T >.
Definition at line 1173 of file petsc_vector.h.
|
inlinestaticinherited |
Prints the number of outstanding (created, but not yet destroyed) objects.
Definition at line 85 of file reference_counter.h.
References libMesh::ReferenceCounter::_n_objects.
Referenced by libMesh::LibMeshInit::~LibMeshInit().
|
inlineinherited |
Definition at line 103 of file parallel_object.h.
References libMesh::ParallelObject::_communicator, libMesh::libmesh_assert(), and TIMPI::Communicator::size().
Referenced by libMesh::Partitioner::_find_global_index_by_pid_map(), libMesh::BoundaryInfo::_find_id_maps(), libMesh::DofMap::add_constraints_to_send_list(), libMesh::PetscDMWrapper::add_dofs_to_section(), libMesh::DistributedMesh::add_elem(), libMesh::DofMap::add_neighbors_to_send_list(), libMesh::DistributedMesh::add_node(), libMesh::System::add_vector(), libMesh::LaplaceMeshSmoother::allgather_graph(), libMesh::DofMap::allgather_recursive_constraints(), libMesh::FEMSystem::assembly(), libMesh::Nemesis_IO::assert_symmetric_cmaps(), libMesh::Partitioner::assign_partitioning(), libMesh::AztecLinearSolver< T >::AztecLinearSolver(), libMesh::Partitioner::build_graph(), libMesh::EquationSystems::build_parallel_elemental_solution_vector(), libMesh::DistributedMesh::clear(), libMesh::DistributedMesh::clear_elems(), libMesh::Nemesis_IO_Helper::compute_border_node_ids(), libMesh::Nemesis_IO_Helper::construct_nemesis_filename(), libMesh::UnstructuredMesh::copy_nodes_and_elements(), libMesh::ExodusII_IO::copy_scalar_solution(), libMesh::Nemesis_IO::copy_scalar_solution(), libMesh::UnstructuredMesh::create_pid_mesh(), libMesh::MeshTools::create_processor_bounding_box(), libMesh::DofMap::distribute_dofs(), libMesh::DofMap::distribute_scalar_dofs(), libMesh::DistributedMesh::DistributedMesh(), libMesh::EnsightIO::EnsightIO(), libMesh::RBEIMEvaluation::gather_bfs(), libMesh::MeshBase::get_info(), libMesh::StaticCondensation::init(), libMesh::SystemSubsetBySubdomain::init(), libMesh::PetscDMWrapper::init_petscdm(), libMesh::Nemesis_IO_Helper::initialize(), libMesh::ExodusII_IO_Helper::initialize(), libMesh::DistributedMesh::insert_elem(), libMesh::MeshTools::libmesh_assert_contiguous_dof_ids(), libMesh::MeshTools::libmesh_assert_parallel_consistent_new_node_procids(), libMesh::MeshTools::libmesh_assert_parallel_consistent_procids< Elem >(), libMesh::MeshTools::libmesh_assert_parallel_consistent_procids< Node >(), libMesh::MeshTools::libmesh_assert_topology_consistent_procids< Node >(), libMesh::MeshTools::libmesh_assert_valid_boundary_ids(), libMesh::MeshTools::libmesh_assert_valid_dof_ids(), libMesh::MeshTools::libmesh_assert_valid_neighbors(), libMesh::MeshTools::libmesh_assert_valid_refinement_flags(), libMesh::DofMap::local_variable_indices(), libMesh::MeshRefinement::make_coarsening_compatible(), libMesh::MeshBase::n_active_elem_on_proc(), libMesh::DofMap::n_dofs_per_processor(), libMesh::MeshBase::n_elem_on_proc(), libMesh::MeshBase::n_nodes_on_proc(), libMesh::RBEIMEvaluation::node_gather_bfs(), libMesh::Partitioner::partition(), libMesh::MeshBase::partition(), libMesh::Partitioner::partition_unpartitioned_elements(), libMesh::System::point_gradient(), libMesh::System::point_hessian(), libMesh::System::point_value(), libMesh::DofMap::prepare_send_list(), libMesh::MeshBase::print_constraint_rows(), libMesh::DofMap::print_dof_constraints(), libMesh::NameBasedIO::read(), libMesh::Nemesis_IO::read(), libMesh::CheckpointIO::read(), libMesh::CheckpointIO::read_connectivity(), libMesh::XdrIO::read_header(), libMesh::CheckpointIO::read_nodes(), libMesh::System::read_parallel_data(), libMesh::System::read_SCALAR_dofs(), libMesh::System::read_serialized_blocked_dof_objects(), libMesh::System::read_serialized_vector(), libMesh::DistributedMesh::renumber_dof_objects(), libMesh::Partitioner::repartition(), OverlappingFunctorTest::run_partitioner_test(), libMesh::DofMap::scatter_constraints(), libMesh::DistributedMesh::set_next_unique_id(), libMesh::DofMap::set_nonlocal_dof_objects(), libMesh::PetscDMWrapper::set_point_range_in_section(), WriteVecAndScalar::setupTests(), libMesh::RBEIMEvaluation::side_gather_bfs(), DistributedMeshTest::testRemoteElemError(), CheckpointIOTest::testSplitter(), libMesh::MeshRefinement::uniformly_coarsen(), libMesh::DistributedMesh::update_parallel_id_counts(), libMesh::GMVIO::write_binary(), libMesh::GMVIO::write_discontinuous_gmv(), libMesh::ExodusII_IO_Helper::write_nodal_coordinates(), libMesh::VTKIO::write_nodal_data(), libMesh::ExodusII_IO::write_nodal_data(), libMesh::System::write_parallel_data(), libMesh::System::write_SCALAR_dofs(), libMesh::XdrIO::write_serialized_bcs_helper(), libMesh::System::write_serialized_blocked_dof_objects(), libMesh::XdrIO::write_serialized_connectivity(), libMesh::XdrIO::write_serialized_nodes(), and libMesh::XdrIO::write_serialized_nodesets().
Definition at line 60 of file petsc_vector.C.
|
inlineoverridevirtual |
Implements libMesh::NumericVector< T >.
Definition at line 1099 of file petsc_vector.h.
|
overridevirtual |
Computes the component-wise multiplication of this vector's entries by another's, \( u_i \leftarrow u_i v_i \, \forall i\).
Implements libMesh::NumericVector< T >.
Definition at line 385 of file petsc_vector.C.
|
inlineinherited |
Scales the vector by a
, \( \vec{u} \leftarrow a\vec{u} \).
Equivalent to u.scale(a)
Definition at line 420 of file numeric_vector.h.
|
overridevirtual |
Adds v
to *this, \( \vec{u} \leftarrow \vec{u} + \vec{v} \).
Equivalent to u.add(1, v)
.
Implements libMesh::NumericVector< T >.
Definition at line 93 of file petsc_vector.C.
|
overridevirtual |
Subtracts v
from *this, \( \vec{u} \leftarrow \vec{u} - \vec{v} \).
Equivalent to u.add
(-1, v).
Implements libMesh::NumericVector< T >.
Definition at line 109 of file petsc_vector.C.
|
overridevirtual |
Computes the component-wise division of this vector's entries by another's, \( u_i \leftarrow \frac{u_i}{v_i} \, \forall i\).
Implements libMesh::NumericVector< T >.
Definition at line 398 of file petsc_vector.C.
|
inlineinherited |
Scales the vector by 1/a
, \( \vec{u} \leftarrow \frac{1}{a}\vec{u} \).
Equivalent to u.scale
(1./a)
Definition at line 438 of file numeric_vector.h.
PetscVector< T > & libMesh::PetscVector< T >::operator= | ( | const PetscVector< T > & | v | ) |
Copy assignment operator.
Calls VecCopy after performing various checks.
Definition at line 512 of file petsc_vector.C.
|
delete |
|
overridevirtual |
Sets all entries of the vector to the value s
.
Implements libMesh::NumericVector< T >.
Definition at line 470 of file petsc_vector.C.
|
overridevirtual |
This looks like a copy assignment operator, but note that, unlike normal copy assignment operators, it is pure virtual.
This function should be overridden in derived classes so that they can be copied correctly via references to the base class. This design usually isn't a good idea in general, but in this context it works because we usually don't have a mix of different kinds of NumericVectors active in the library at a single time.
Implements libMesh::NumericVector< T >.
Definition at line 496 of file petsc_vector.C.
|
overridevirtual |
Sets (*this)(i) = v(i) for each entry of the vector.
Case 1: The vector is the same size of The global vector. Only add the local components.
Case 2: The vector is the same size as our local piece. Insert directly to the local piece.
Implements libMesh::NumericVector< T >.
Definition at line 540 of file petsc_vector.C.
|
overridevirtual |
Computes \( u_i \leftarrow \frac{v_{1,i}}{v_{2,i}} \) (summation not implied) i.e.
the pointwise (component-wise) division of vec1
and vec2
, and stores the result in *this
.
Implements libMesh::NumericVector< T >.
Definition at line 1012 of file petsc_vector.C.
|
overridevirtual |
Computes \( u_i \leftarrow u_i v_i \) (summation not implied) i.e.
the pointwise (component-wise) product of vec1
and vec2
, and stores the result in *this
.
Implements libMesh::NumericVector< T >.
Definition at line 989 of file petsc_vector.C.
|
inlinevirtualinherited |
Prints the local contents of the vector, by default to libMesh::out.
Definition at line 1040 of file numeric_vector.h.
|
inlineinherited |
Definition at line 1022 of file numeric_vector.h.
|
inlinevirtualinherited |
Prints the global contents of the vector, by default to libMesh::out.
Definition at line 1077 of file numeric_vector.h.
|
inlineinherited |
Definition at line 1055 of file numeric_vector.h.
|
staticinherited |
Prints the reference information, by default to libMesh::out
.
Definition at line 81 of file reference_counter.C.
References libMesh::ReferenceCounter::_enable_print_counter, and libMesh::ReferenceCounter::get_info().
Referenced by libMesh::LibMeshInit::~LibMeshInit().
|
overridevirtual |
Print the contents of the vector in Matlab's sparse matrix format.
Optionally prints the vector to the file named name
. If name
is not specified it is dumped to the screen.
Reimplemented from libMesh::NumericVector< T >.
Definition at line 1035 of file petsc_vector.C.
|
inlineinherited |
Definition at line 114 of file parallel_object.h.
References libMesh::ParallelObject::_communicator, and TIMPI::Communicator::rank().
Referenced by libMesh::BoundaryInfo::_find_id_maps(), libMesh::PetscDMWrapper::add_dofs_to_section(), libMesh::DistributedMesh::add_elem(), libMesh::BoundaryInfo::add_elements(), libMesh::DofMap::add_neighbors_to_send_list(), libMesh::DistributedMesh::add_node(), libMesh::MeshTools::Modification::all_tri(), libMesh::DofMap::allgather_recursive_constraints(), libMesh::FEMSystem::assembly(), libMesh::Nemesis_IO::assert_symmetric_cmaps(), libMesh::Partitioner::assign_partitioning(), libMesh::Nemesis_IO_Helper::build_element_and_node_maps(), libMesh::Partitioner::build_graph(), libMesh::InfElemBuilder::build_inf_elem(), libMesh::BoundaryInfo::build_node_list_from_side_list(), libMesh::EquationSystems::build_parallel_elemental_solution_vector(), libMesh::EquationSystems::build_parallel_solution_vector(), libMesh::MeshFunction::check_found_elem(), libMesh::DistributedMesh::clear(), libMesh::DistributedMesh::clear_elems(), libMesh::ExodusII_IO_Helper::close(), libMesh::Nemesis_IO_Helper::compute_border_node_ids(), libMesh::Nemesis_IO_Helper::compute_communication_map_parameters(), libMesh::Nemesis_IO_Helper::compute_internal_and_border_elems_and_internal_nodes(), libMesh::RBConstruction::compute_max_error_bound(), libMesh::Nemesis_IO_Helper::compute_node_communication_maps(), libMesh::Nemesis_IO_Helper::compute_num_global_elem_blocks(), libMesh::Nemesis_IO_Helper::compute_num_global_nodesets(), libMesh::Nemesis_IO_Helper::compute_num_global_sidesets(), libMesh::Nemesis_IO_Helper::construct_nemesis_filename(), libMesh::ExodusII_IO::copy_elemental_solution(), libMesh::ExodusII_IO::copy_nodal_solution(), libMesh::ExodusII_IO::copy_scalar_solution(), libMesh::Nemesis_IO::copy_scalar_solution(), libMesh::MeshTools::correct_node_proc_ids(), libMesh::ExodusII_IO_Helper::create(), libMesh::DistributedMesh::delete_elem(), libMesh::MeshCommunication::delete_remote_elements(), libMesh::DofMap::distribute_dofs(), libMesh::DofMap::distribute_scalar_dofs(), libMesh::DistributedMesh::DistributedMesh(), libMesh::DofMapBase::end_dof(), libMesh::DofMapBase::end_old_dof(), libMesh::EnsightIO::EnsightIO(), libMesh::GenericProjector< FFunctor, GFunctor, FValue, ProjectionAction >::SubFunctor::find_dofs_to_send(), libMesh::UnstructuredMesh::find_neighbors(), libMesh::DofMapBase::first_dof(), libMesh::DofMapBase::first_old_dof(), libMesh::RBEIMEvaluation::gather_bfs(), libMesh::Nemesis_IO_Helper::get_cmap_params(), libMesh::Nemesis_IO_Helper::get_eb_info_global(), libMesh::Nemesis_IO_Helper::get_elem_cmap(), libMesh::Nemesis_IO_Helper::get_elem_map(), libMesh::MeshBase::get_info(), libMesh::DofMap::get_info(), libMesh::Nemesis_IO_Helper::get_init_global(), libMesh::Nemesis_IO_Helper::get_init_info(), libMesh::RBEIMEvaluation::get_interior_basis_functions_as_vecs(), libMesh::Nemesis_IO_Helper::get_loadbal_param(), libMesh::DofMap::get_local_constraints(), libMesh::MeshBase::get_local_constraints(), libMesh::Nemesis_IO_Helper::get_node_cmap(), libMesh::Nemesis_IO_Helper::get_node_map(), libMesh::Nemesis_IO_Helper::get_ns_param_global(), libMesh::Nemesis_IO_Helper::get_ss_param_global(), libMesh::SparsityPattern::Build::handle_vi_vj(), libMesh::LaplaceMeshSmoother::init(), libMesh::SystemSubsetBySubdomain::init(), HeatSystem::init_data(), libMesh::ExodusII_IO_Helper::initialize(), libMesh::ExodusII_IO_Helper::initialize_element_variables(), libMesh::ExodusII_IO_Helper::initialize_global_variables(), libMesh::ExodusII_IO_Helper::initialize_nodal_variables(), libMesh::DistributedMesh::insert_elem(), libMesh::DofMap::is_evaluable(), libMesh::SparsityPattern::Build::join(), libMesh::TransientRBEvaluation::legacy_write_offline_data_to_files(), libMesh::RBSCMEvaluation::legacy_write_offline_data_to_files(), libMesh::RBEvaluation::legacy_write_offline_data_to_files(), libMesh::MeshTools::libmesh_assert_consistent_distributed(), libMesh::MeshTools::libmesh_assert_consistent_distributed_nodes(), libMesh::MeshTools::libmesh_assert_contiguous_dof_ids(), libMesh::MeshTools::libmesh_assert_parallel_consistent_procids< Elem >(), libMesh::MeshTools::libmesh_assert_valid_neighbors(), libMesh::DistributedMesh::libmesh_assert_valid_parallel_object_ids(), libMesh::DofMap::local_variable_indices(), main(), libMesh::MeshRefinement::make_coarsening_compatible(), AugmentSparsityOnInterface::mesh_reinit(), libMesh::TriangulatorInterface::MeshedHole::MeshedHole(), libMesh::MeshBase::n_active_local_elem(), libMesh::BoundaryInfo::n_boundary_conds(), libMesh::MeshTools::n_connected_components(), libMesh::MeshBase::n_constraint_rows(), libMesh::BoundaryInfo::n_edge_conds(), libMesh::DofMapBase::n_local_dofs(), libMesh::MeshBase::n_local_elem(), libMesh::MeshBase::n_local_nodes(), libMesh::BoundaryInfo::n_nodeset_conds(), libMesh::BoundaryInfo::n_shellface_conds(), libMesh::RBEIMEvaluation::node_gather_bfs(), libMesh::DistributedMesh::own_node(), libMesh::BoundaryInfo::parallel_sync_node_ids(), libMesh::BoundaryInfo::parallel_sync_side_ids(), libMesh::System::point_gradient(), libMesh::System::point_hessian(), libMesh::System::point_value(), libMesh::MeshBase::print_constraint_rows(), libMesh::DofMap::print_dof_constraints(), libMesh::DofMap::process_mesh_constraint_rows(), libMesh::Nemesis_IO_Helper::put_cmap_params(), libMesh::Nemesis_IO_Helper::put_elem_cmap(), libMesh::Nemesis_IO_Helper::put_elem_map(), libMesh::Nemesis_IO_Helper::put_loadbal_param(), libMesh::Nemesis_IO_Helper::put_node_cmap(), libMesh::Nemesis_IO_Helper::put_node_map(), libMesh::NameBasedIO::read(), libMesh::Nemesis_IO::read(), libMesh::XdrIO::read(), libMesh::CheckpointIO::read(), libMesh::EquationSystems::read(), libMesh::ExodusII_IO_Helper::read_elem_num_map(), libMesh::ExodusII_IO_Helper::read_global_values(), libMesh::ExodusII_IO::read_header(), libMesh::CheckpointIO::read_header(), libMesh::XdrIO::read_header(), libMesh::System::read_header(), libMesh::System::read_legacy_data(), libMesh::DynaIO::read_mesh(), libMesh::ExodusII_IO_Helper::read_node_num_map(), libMesh::System::read_parallel_data(), libMesh::TransientRBConstruction::read_riesz_representors_from_files(), libMesh::RBConstruction::read_riesz_representors_from_files(), libMesh::System::read_SCALAR_dofs(), libMesh::XdrIO::read_serialized_bc_names(), libMesh::XdrIO::read_serialized_bcs_helper(), libMesh::System::read_serialized_blocked_dof_objects(), libMesh::XdrIO::read_serialized_connectivity(), libMesh::System::read_serialized_data(), libMesh::XdrIO::read_serialized_nodes(), libMesh::XdrIO::read_serialized_nodesets(), libMesh::XdrIO::read_serialized_subdomain_names(), libMesh::System::read_serialized_vector(), libMesh::System::read_serialized_vectors(), libMesh::Nemesis_IO_Helper::read_var_names_impl(), libMesh::SimplexRefiner::refine_via_edges(), libMesh::StaticCondensationDofMap::reinit(), libMesh::DistributedMesh::renumber_dof_objects(), libMesh::DistributedMesh::renumber_nodes_and_elements(), libMesh::DofMap::scatter_constraints(), libMesh::CheckpointIO::select_split_config(), libMesh::DistributedMesh::set_next_unique_id(), libMesh::DofMap::set_nonlocal_dof_objects(), libMesh::PetscDMWrapper::set_point_range_in_section(), libMesh::RBEIMEvaluation::side_gather_bfs(), ExodusTest< elem_type >::test_read_gold(), ExodusTest< elem_type >::test_write(), MeshInputTest::testAbaqusRead(), MeshInputTest::testBadGmsh(), MeshInputTest::testCopyElementSolutionImpl(), MeshInputTest::testCopyElementVectorImpl(), MeshInputTest::testCopyNodalSolutionImpl(), DefaultCouplingTest::testCoupling(), PointNeighborCouplingTest::testCoupling(), MeshInputTest::testDynaFileMappings(), MeshInputTest::testDynaNoSplines(), MeshInputTest::testDynaReadElem(), MeshInputTest::testDynaReadPatch(), MeshInputTest::testExodusFileMappings(), MeshInputTest::testExodusIGASidesets(), MeshInputTest::testExodusWriteElementDataFromDiscontinuousNodalData(), MeshInputTest::testGmshBCIDOverlap(), MeshInputTest::testGoodGmsh(), MeshInputTest::testGoodSTL(), MeshInputTest::testGoodSTLBinary(), MeshInputTest::testLowOrderEdgeBlocks(), SystemsTest::testProjectMatrix3D(), BoundaryInfoTest::testShellFaceConstraints(), MeshInputTest::testSingleElementImpl(), WriteVecAndScalar::testSolution(), CheckpointIOTest::testSplitter(), MeshInputTest::testTetgenIO(), libMesh::MeshTools::total_weight(), libMesh::NetGenMeshInterface::triangulate(), libMesh::MeshRefinement::uniformly_coarsen(), libMesh::DistributedMesh::update_parallel_id_counts(), libMesh::DTKAdapter::update_variable_values(), libMesh::MeshTools::volume(), libMesh::STLIO::write(), libMesh::NameBasedIO::write(), libMesh::XdrIO::write(), libMesh::CheckpointIO::write(), libMesh::EquationSystems::write(), libMesh::GMVIO::write_discontinuous_gmv(), libMesh::ExodusII_IO::write_element_data(), libMesh::ExodusII_IO_Helper::write_element_values(), libMesh::ExodusII_IO_Helper::write_element_values_element_major(), libMesh::ExodusII_IO_Helper::write_elements(), libMesh::ExodusII_IO_Helper::write_elemset_data(), libMesh::ExodusII_IO_Helper::write_elemsets(), libMesh::ExodusII_IO::write_global_data(), libMesh::ExodusII_IO_Helper::write_global_values(), libMesh::System::write_header(), libMesh::ExodusII_IO::write_information_records(), libMesh::ExodusII_IO_Helper::write_information_records(), libMesh::ExodusII_IO_Helper::write_nodal_coordinates(), libMesh::UCDIO::write_nodal_data(), libMesh::VTKIO::write_nodal_data(), libMesh::ExodusII_IO::write_nodal_data(), libMesh::ExodusII_IO::write_nodal_data_common(), libMesh::ExodusII_IO::write_nodal_data_discontinuous(), libMesh::ExodusII_IO_Helper::write_nodal_values(), libMesh::ExodusII_IO_Helper::write_nodeset_data(), libMesh::Nemesis_IO_Helper::write_nodesets(), libMesh::ExodusII_IO_Helper::write_nodesets(), libMesh::RBEIMEvaluation::write_out_interior_basis_functions(), libMesh::RBEIMEvaluation::write_out_node_basis_functions(), libMesh::RBEIMEvaluation::write_out_side_basis_functions(), write_output_solvedata(), libMesh::System::write_parallel_data(), libMesh::RBConstruction::write_riesz_representors_to_files(), libMesh::System::write_SCALAR_dofs(), libMesh::XdrIO::write_serialized_bc_names(), libMesh::XdrIO::write_serialized_bcs_helper(), libMesh::System::write_serialized_blocked_dof_objects(), libMesh::XdrIO::write_serialized_connectivity(), libMesh::System::write_serialized_data(), libMesh::XdrIO::write_serialized_nodes(), libMesh::XdrIO::write_serialized_nodesets(), libMesh::XdrIO::write_serialized_subdomain_names(), libMesh::System::write_serialized_vector(), libMesh::System::write_serialized_vectors(), libMesh::ExodusII_IO_Helper::write_sideset_data(), libMesh::Nemesis_IO_Helper::write_sidesets(), libMesh::ExodusII_IO_Helper::write_sidesets(), libMesh::ExodusII_IO::write_timestep(), libMesh::ExodusII_IO_Helper::write_timestep(), and libMesh::ExodusII_IO::write_timestep_discontinuous().
|
virtualinherited |
Read the contents of the vector from the Matlab-script format used by PETSc.
If the size of the matrix in filename
appears consistent with the existing partitioning of this
then the existing parallel decomposition will be retained. If not, then this
will be initialized with the size from the file, linearly partitioned onto the number of processors available.
Definition at line 297 of file numeric_vector.C.
|
inherited |
Definition at line 685 of file numeric_vector.C.
Referenced by libMesh::NumericVector< Number >::add_vector(), libMesh::NumericVector< Number >::compatible(), and libMesh::NumericVector< Number >::insert().
|
overridevirtual |
Computes the component-wise reciprocal, \( u_i \leftarrow \frac{1}{u_i} \, \forall i\).
Implements libMesh::NumericVector< T >.
Definition at line 141 of file petsc_vector.C.
|
inline |
Restore the data array.
Definition at line 1163 of file petsc_vector.h.
Referenced by PetscVectorTest::testGetArray().
|
overridevirtual |
Restores a view into this vector using the indices in rows
.
This is currently only implemented for PetscVectors.
Reimplemented from libMesh::NumericVector< T >.
Definition at line 1315 of file petsc_vector.C.
|
overridevirtual |
Scale each element of the vector by the given factor
.
Implements libMesh::NumericVector< T >.
Definition at line 368 of file petsc_vector.C.
|
overridevirtual |
Sets v(i) = value
.
Note that library implementations of this method are thread safe, e.g. we will lock _numeric_vector_mutex
before writing to the vector
Implements libMesh::NumericVector< T >.
Definition at line 124 of file petsc_vector.C.
Referenced by PetscVectorTest::testGetArray(), and PetscVectorTest::testPetscOperations().
|
inherited |
Allow the user to change the ParallelType of the NumericVector under some circumstances.
If the NumericVector has not been initialized yet, then it is generally safe to change the ParallelType. otherwise, if the NumericVector has already been initialized with a specific type, we cannot change it without doing some extra copying/reinitialization work, and we currently throw an error if this is requested.
Definition at line 97 of file numeric_vector.C.
Referenced by libMesh::System::add_vector().
|
inlineoverridevirtual |
Implements libMesh::NumericVector< T >.
Definition at line 964 of file petsc_vector.h.
Referenced by libMesh::PetscVector< libMesh::Number >::add(), libMesh::PetscVector< libMesh::Number >::localize(), and libMesh::PetscVector< libMesh::Number >::operator=().
|
virtualinherited |
Definition at line 559 of file numeric_vector.C.
Referenced by libMesh::System::discrete_var_norm().
|
virtualinherited |
Definition at line 576 of file numeric_vector.C.
Referenced by libMesh::System::discrete_var_norm().
|
virtualinherited |
Definition at line 593 of file numeric_vector.C.
Referenced by libMesh::System::discrete_var_norm().
|
overridevirtual |
Implements libMesh::NumericVector< T >.
Definition at line 44 of file petsc_vector.C.
|
inlineoverridevirtual |
Swaps the contents of this with v
.
There should be enough indirection in subclasses to make this an O(1) header-swap operation.
Reimplemented from libMesh::NumericVector< T >.
Definition at line 1211 of file petsc_vector.h.
Referenced by libMesh::__libmesh_petsc_diff_solver_residual(), libMesh::__libmesh_tao_equality_constraints(), libMesh::__libmesh_tao_equality_constraints_jacobian(), libMesh::__libmesh_tao_gradient(), libMesh::__libmesh_tao_hessian(), libMesh::__libmesh_tao_inequality_constraints(), libMesh::__libmesh_tao_inequality_constraints_jacobian(), libMesh::__libmesh_tao_objective(), DMlibMeshFunction(), libMesh::libmesh_petsc_snes_fd_residual(), libMesh::libmesh_petsc_snes_mffd_residual(), libMesh::libmesh_petsc_snes_precheck(), libMesh::libmesh_petsc_snes_residual(), libMesh::libmesh_petsc_snes_residual_helper(), and update_current_local_solution().
|
inlineinherited |
Definition at line 154 of file numeric_vector.h.
Referenced by libMesh::System::add_vector(), libMesh::MeshFunction::check_found_elem(), libMesh::DofMap::constrain_element_residual(), libMesh::PetscVector< libMesh::Number >::create_subvector(), libMesh::EpetraVector< T >::EpetraVector(), libMesh::DofMap::heterogeneously_constrain_element_jacobian_and_residual(), libMesh::DofMap::heterogeneously_constrain_element_residual(), libMesh::DistributedVector< T >::init(), libMesh::LaspackVector< T >::init(), libMesh::EigenSparseVector< T >::init(), libMesh::EpetraVector< T >::init(), libMesh::PetscVector< libMesh::Number >::localize(), libMesh::System::project_vector(), and libMesh::System::read_serialized_vector().
|
inlineinherited |
Definition at line 164 of file numeric_vector.h.
|
inline |
Definition at line 345 of file petsc_vector.h.
Referenced by libMesh::PetscVector< libMesh::Number >::add(), libMesh::PetscPreconditioner< T >::apply(), libMesh::SlepcEigenSolver< libMesh::Number >::attach_deflation_space(), DMCreateGlobalVector_libMesh(), libMesh::PetscShellMatrix< T >::get_diagonal(), libMesh::PetscMatrix< T >::get_diagonal(), libMesh::TaoOptimizationSolver< T >::get_dual_variables(), libMesh::SlepcEigenSolver< libMesh::Number >::get_eigenpair(), main(), libMesh::PetscVector< libMesh::Number >::pointwise_divide(), libMesh::PetscVector< libMesh::Number >::pointwise_mult(), libMesh::PetscDiffSolver::solve(), libMesh::PetscLinearSolver< Number >::solve_base(), libMesh::PetscShellMatrix< T >::vector_mult(), and libMesh::PetscShellMatrix< T >::vector_mult_add().
|
inline |
Definition at line 347 of file petsc_vector.h.
|
inlineoverridevirtual |
Set all entries to zero.
Equivalent to v
= 0, but more obvious and faster.
Implements libMesh::NumericVector< T >.
Definition at line 912 of file petsc_vector.h.
Referenced by libMesh::__libmesh_tao_equality_constraints(), libMesh::__libmesh_tao_gradient(), libMesh::__libmesh_tao_inequality_constraints(), libMesh::libmesh_petsc_snes_fd_residual(), libMesh::libmesh_petsc_snes_mffd_residual(), libMesh::libmesh_petsc_snes_residual(), and libMesh::TaoOptimizationSolver< T >::solve().
|
inlineoverridevirtual |
Implements libMesh::NumericVector< T >.
Definition at line 939 of file petsc_vector.h.
|
mutableprivate |
If true
, the actual PETSc array of the values of the vector is currently accessible.
That means that the members _local_form
and _values
are valid.
Definition at line 369 of file petsc_vector.h.
Referenced by libMesh::PetscVector< libMesh::Number >::swap().
|
mutableprivate |
Definition at line 371 of file petsc_vector.h.
|
protectedinherited |
Definition at line 120 of file parallel_object.h.
Referenced by libMesh::EquationSystems::build_parallel_elemental_solution_vector(), libMesh::EquationSystems::build_parallel_solution_vector(), libMesh::StaticCondensation::close(), libMesh::ParallelObject::comm(), libMesh::CondensedEigenSystem::initialize_condensed_matrices(), libMesh::ParallelObject::n_processors(), libMesh::ParallelObject::operator=(), libMesh::ParallelObject::processor_id(), libMesh::BoundaryInfo::regenerate_id_sets(), and libMesh::StaticCondensationDofMap::reinit().
|
staticprotectedinherited |
Actually holds the data.
Definition at line 124 of file reference_counter.h.
Referenced by libMesh::ReferenceCounter::get_info().
|
private |
This boolean value should only be set to false for the constructor which takes a PETSc Vec object.
Definition at line 456 of file petsc_vector.h.
Referenced by libMesh::PetscVector< libMesh::Number >::swap().
|
staticprotectedinherited |
Flag to control whether reference count information is printed when print_info is called.
Definition at line 143 of file reference_counter.h.
Referenced by libMesh::ReferenceCounter::disable_print_counter_info(), libMesh::ReferenceCounter::enable_print_counter_info(), and libMesh::ReferenceCounter::print_info().
|
mutableprivate |
|
private |
Map that maps global to local ghost cells (will be empty if not in ghost cell mode).
Definition at line 450 of file petsc_vector.h.
Referenced by libMesh::PetscVector< libMesh::Number >::init(), and libMesh::PetscVector< libMesh::Number >::swap().
|
protectedinherited |
Flag which tracks whether the vector's values are consistent on all processors after insertion or addition of values has occurred on some or all processors.
Definition at line 831 of file numeric_vector.h.
Referenced by libMesh::NumericVector< Number >::closed(), libMesh::PetscVector< libMesh::Number >::create_subvector(), libMesh::EpetraVector< T >::EpetraVector(), libMesh::DistributedVector< T >::localize(), libMesh::DistributedVector< T >::operator=(), and libMesh::DistributedVector< T >::swap().
|
protectedinherited |
true
once init() has been called.
Definition at line 836 of file numeric_vector.h.
Referenced by libMesh::PetscVector< libMesh::Number >::create_subvector(), libMesh::EpetraVector< T >::EpetraVector(), libMesh::NumericVector< Number >::initialized(), libMesh::DistributedVector< T >::localize(), libMesh::DistributedVector< T >::operator=(), and libMesh::DistributedVector< T >::swap().
|
mutableprivate |
|
mutableprivate |
PETSc vector datatype to hold the local form of a ghosted vector.
The contents of this field are only valid if the vector is ghosted and _array_is_present
is true
.
Definition at line 398 of file petsc_vector.h.
Referenced by libMesh::PetscVector< libMesh::Number >::swap().
|
mutableprivate |
Size of the local values from _get_array()
Definition at line 391 of file petsc_vector.h.
|
staticprotectedinherited |
Mutual exclusion object to enable thread-safe reference counting.
Definition at line 137 of file reference_counter.h.
|
staticprotectedinherited |
The number of objects.
Print the reference count information when the number returns to 0.
Definition at line 132 of file reference_counter.h.
Referenced by libMesh::ReferenceCounter::n_objects(), libMesh::ReferenceCounter::ReferenceCounter(), and libMesh::ReferenceCounter::~ReferenceCounter().
|
protectedinherited |
Mutex for performing thread-safe operations.
Definition at line 846 of file numeric_vector.h.
|
mutableprivate |
Mutex for _get_array and _restore_array.
This is part of the object to keep down thread contention when reading from multiple PetscVectors simultaneously
Definition at line 425 of file petsc_vector.h.
|
mutableprivate |
Pointer to the actual PETSc array of the values of the vector.
This pointer is only valid if _array_is_present
is true
. We're using PETSc's VecGetArrayRead() function, which requires a constant PetscScalar *, but _get_array and _restore_array are const member functions, so _values also needs to be mutable (otherwise it is a "const PetscScalar * const" in that context).
Definition at line 408 of file petsc_vector.h.
|
protectedinherited |
Type of vector.
Definition at line 841 of file numeric_vector.h.
Referenced by libMesh::PetscVector< libMesh::Number >::create_subvector(), libMesh::DistributedVector< T >::DistributedVector(), libMesh::EigenSparseVector< T >::EigenSparseVector(), libMesh::EpetraVector< T >::EpetraVector(), libMesh::PetscVector< libMesh::Number >::init(), libMesh::LaspackVector< T >::LaspackVector(), and libMesh::NumericVector< Number >::type().
|
mutableprivate |
Pointer to the actual PETSc array of the values of the vector.
This pointer is only valid if _array_is_present
is true
. We're using PETSc's VecGetArrayRead() function, which requires a constant PetscScalar *, but _get_array and _restore_array are const member functions, so _values also needs to be mutable (otherwise it is a "const PetscScalar * const" in that context).
Definition at line 418 of file petsc_vector.h.
Referenced by libMesh::PetscVector< libMesh::Number >::swap().
|
mutableprivate |
Whether or not the data array has been manually retrieved using get_array() or get_array_read()
Definition at line 462 of file petsc_vector.h.
|
mutableprivate |
Whether or not the data array is for read only access.
Definition at line 467 of file petsc_vector.h.
|
private |
Actual PETSc vector datatype to hold vector entries.
Definition at line 360 of file petsc_vector.h.
Referenced by libMesh::PetscVector< libMesh::Number >::add_vector(), libMesh::PetscVector< libMesh::Number >::add_vector_conjugate_transpose(), libMesh::PetscVector< libMesh::Number >::add_vector_transpose(), libMesh::PetscVector< libMesh::Number >::create_subvector(), libMesh::PetscVector< libMesh::Number >::dot(), libMesh::PetscVector< libMesh::Number >::indefinite_dot(), libMesh::PetscVector< libMesh::Number >::init(), libMesh::PetscVector< libMesh::Number >::localize(), libMesh::PetscVector< libMesh::Number >::operator*=(), libMesh::PetscVector< libMesh::Number >::operator/=(), libMesh::PetscVector< libMesh::Number >::operator=(), libMesh::PetscVector< libMesh::Number >::swap(), and libMesh::PetscVector< libMesh::Number >::vec().