Skip to content

Commit 3ca87a2

Browse files
authored
Merge pull request #4430 from roystgnr/petsc_serial_ghosted
Use is_effectively_ghosted() in PetscVector
2 parents 62e5b34 + 16c3fb3 commit 3ca87a2

3 files changed

Lines changed: 31 additions & 22 deletions

File tree

include/numerics/numeric_vector.h

Lines changed: 13 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -168,6 +168,19 @@ class NumericVector : public ReferenceCountedObject<NumericVector<T>>,
168168
bool is_effectively_serial() const
169169
{ return (this->n_processors()==1) || (_type == SERIAL); }
170170

171+
/**
172+
* \returns Whether the vector is effectively ghosted, i.e. whether
173+
* vector operations require handling of ghosted coefficients.
174+
*
175+
* This is true for GHOSTED vectors in parallel, but in serial
176+
* subclasses may return false here for GHOSTED vectors since every
177+
* ghosted vector is effectively serial and thus may be able to take
178+
* advantage of serial-specific code paths.
179+
*/
180+
181+
bool is_effectively_ghosted() const
182+
{ return (this->n_processors()!=1) && (_type == GHOSTED); }
183+
171184
/**
172185
* \returns The type (SERIAL, PARALLEL, GHOSTED) of the vector.
173186
*

include/numerics/petsc_vector.h

Lines changed: 6 additions & 10 deletions
Original file line numberDiff line numberDiff line change
@@ -916,7 +916,7 @@ void PetscVector<T>::close ()
916916

917917
VecAssemblyBeginEnd(this->comm(), _vec);
918918

919-
if (this->type() == GHOSTED)
919+
if (this->is_effectively_ghosted())
920920
VecGhostUpdateBeginEnd(this->comm(), _vec, INSERT_VALUES, SCATTER_FORWARD);
921921

922922
this->_is_closed = true;
@@ -961,7 +961,7 @@ void PetscVector<T>::zero ()
961961

962962
PetscScalar z=0.;
963963

964-
if (this->type() != GHOSTED)
964+
if (!this->is_effectively_ghosted())
965965
LibmeshPetscCall(VecSet (_vec, z));
966966
else
967967
{
@@ -1145,10 +1145,8 @@ T PetscVector<T>::operator() (const numeric_index_type i) const
11451145
const numeric_index_type local_index = this->map_global_to_local_index(i);
11461146

11471147
#ifndef NDEBUG
1148-
if (this->type() == GHOSTED)
1149-
{
1150-
libmesh_assert_less (local_index, _local_size);
1151-
}
1148+
if (this->is_effectively_ghosted())
1149+
libmesh_assert_less (local_index, _local_size);
11521150
#endif
11531151

11541152
return static_cast<T>(_read_only_values[local_index]);
@@ -1169,10 +1167,8 @@ void PetscVector<T>::get(const std::vector<numeric_index_type> & index,
11691167
{
11701168
const numeric_index_type local_index = this->map_global_to_local_index(index[i]);
11711169
#ifndef NDEBUG
1172-
if (this->type() == GHOSTED)
1173-
{
1174-
libmesh_assert_less (local_index, _local_size);
1175-
}
1170+
if (this->is_effectively_ghosted())
1171+
libmesh_assert_less (local_index, _local_size);
11761172
#endif
11771173
values[i] = static_cast<T>(_read_only_values[local_index]);
11781174
}

src/numerics/petsc_vector.C

Lines changed: 12 additions & 12 deletions
Original file line numberDiff line numberDiff line change
@@ -338,7 +338,7 @@ void PetscVector<T>::add (const T a_in, const NumericVector<T> & v_in)
338338
libmesh_assert(this->comm().verify(
339339
static_cast<typename std::underlying_type<ParallelType>::type>(this->type())));
340340

341-
if (this->type() == GHOSTED)
341+
if (this->is_effectively_ghosted())
342342
VecGhostUpdateBeginEnd(this->comm(), _vec, INSERT_VALUES, SCATTER_FORWARD);
343343

344344
this->_is_closed = true;
@@ -379,7 +379,7 @@ void PetscVector<T>::scale (const T factor_in)
379379
libmesh_assert(this->comm().verify(
380380
static_cast<typename std::underlying_type<ParallelType>::type>(this->type())));
381381

382-
if (this->type() == GHOSTED)
382+
if (this->is_effectively_ghosted())
383383
VecGhostUpdateBeginEnd(this->comm(), _vec, INSERT_VALUES, SCATTER_FORWARD);
384384
}
385385

@@ -421,7 +421,7 @@ void PetscVector<T>::abs()
421421
libmesh_assert(this->comm().verify(
422422
static_cast<typename std::underlying_type<ParallelType>::type>(this->type())));
423423

424-
if (this->type() == GHOSTED)
424+
if (this->is_effectively_ghosted())
425425
VecGhostUpdateBeginEnd(this->comm(), _vec, INSERT_VALUES, SCATTER_FORWARD);
426426
}
427427

@@ -486,7 +486,7 @@ PetscVector<T>::operator = (const T s_in)
486486
libmesh_assert(this->comm().verify(
487487
static_cast<typename std::underlying_type<ParallelType>::type>(this->type())));
488488

489-
if (this->type() == GHOSTED)
489+
if (this->is_effectively_ghosted())
490490
VecGhostUpdateBeginEnd(this->comm(), _vec, INSERT_VALUES, SCATTER_FORWARD);
491491
}
492492

@@ -569,7 +569,7 @@ PetscVector<T>::operator = (const PetscVector<T> & v)
569569
libmesh_assert(this->comm().verify(
570570
static_cast<typename std::underlying_type<ParallelType>::type>(this->type())));
571571

572-
if (this->type() == GHOSTED)
572+
if (this->is_effectively_ghosted())
573573
VecGhostUpdateBeginEnd(this->comm(), _vec, INSERT_VALUES, SCATTER_FORWARD);
574574

575575
this->_is_closed = true;
@@ -610,7 +610,7 @@ PetscVector<T>::operator = (const std::vector<T> & v)
610610
}
611611

612612
// Make sure ghost dofs are up to date
613-
if (this->type() == GHOSTED)
613+
if (this->is_effectively_ghosted())
614614
this->close();
615615

616616
return *this;
@@ -969,7 +969,7 @@ void PetscVector<T>::pointwise_mult (const NumericVector<T> & vec1,
969969
libmesh_assert(this->comm().verify(
970970
static_cast<typename std::underlying_type<ParallelType>::type>(this->type())));
971971

972-
if (this->type() == GHOSTED)
972+
if (this->is_effectively_ghosted())
973973
VecGhostUpdateBeginEnd(this->comm(), _vec, INSERT_VALUES, SCATTER_FORWARD);
974974

975975
this->_is_closed = true;
@@ -993,7 +993,7 @@ void PetscVector<T>::pointwise_divide (const NumericVector<T> & vec1,
993993
libmesh_assert(this->comm().verify(
994994
static_cast<typename std::underlying_type<ParallelType>::type>(this->type())));
995995

996-
if (this->type() == GHOSTED)
996+
if (this->is_effectively_ghosted())
997997
VecGhostUpdateBeginEnd(this->comm(), _vec, INSERT_VALUES, SCATTER_FORWARD);
998998

999999
this->_is_closed = true;
@@ -1056,7 +1056,7 @@ void PetscVector<T>::create_subvector(NumericVector<T> & subvector,
10561056
parallel_object_only();
10571057

10581058
libmesh_error_msg_if(
1059-
subvector.type() == GHOSTED && !subvector.is_effectively_serial(),
1059+
subvector.is_effectively_ghosted(),
10601060
"We do not support scattering parallel information to ghosts for subvectors");
10611061

10621062
this->_restore_array();
@@ -1184,7 +1184,7 @@ void PetscVector<T>::_get_array(bool read_only) const
11841184
std::scoped_lock lock(_petsc_get_restore_array_mutex);
11851185
if (!_array_is_present.load(std::memory_order_relaxed))
11861186
{
1187-
if (this->type() != GHOSTED)
1187+
if (!this->is_effectively_ghosted())
11881188
{
11891189
if (read_only)
11901190
{
@@ -1243,7 +1243,7 @@ void PetscVector<T>::_restore_array() const
12431243
std::scoped_lock lock(_petsc_get_restore_array_mutex);
12441244
if (_array_is_present.load(std::memory_order_relaxed))
12451245
{
1246-
if (this->type() != GHOSTED)
1246+
if (!this->is_effectively_ghosted())
12471247
{
12481248
if (_values_read_only)
12491249
LibmeshPetscCall(VecRestoreArrayRead (_vec, &_read_only_values));
@@ -1307,7 +1307,7 @@ PetscVector<T>::restore_subvector(std::unique_ptr<NumericVector<T>> subvector,
13071307
Vec subvec = petsc_subvector->vec();
13081308
LibmeshPetscCall(VecRestoreSubVector(_vec, parent_is, &subvec));
13091309

1310-
if (this->type() == GHOSTED)
1310+
if (this->is_effectively_ghosted())
13111311
VecGhostUpdateBeginEnd(this->comm(), _vec, INSERT_VALUES, SCATTER_FORWARD);
13121312

13131313
this->_is_closed = true;

0 commit comments

Comments
 (0)