Skip to content
New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

Add an abstract base class ReadVector purely for vector access. #15197

Merged
merged 9 commits into from
Jul 2, 2023
5 changes: 5 additions & 0 deletions doc/news/changes/incompatibilities/20230701DavidWells
Original file line number Diff line number Diff line change
@@ -0,0 +1,5 @@
Changed: Most of the member functions of FEValues now have different template
parameters. As a result, some function calls which relied on implicit
conversions to ArrayView now require explicit conversions instead.
<br>
(David Wells, 2023/07/01)
4 changes: 4 additions & 0 deletions doc/news/changes/major/20230701DavidWells
Original file line number Diff line number Diff line change
@@ -0,0 +1,4 @@
New: All vector classes in deal.II now inherit from ReadVector, which provides
some common read operations.
<br>
(David Wells, 2023/07/01)
18 changes: 10 additions & 8 deletions include/deal.II/dofs/dof_accessor.h
Original file line number Diff line number Diff line change
Expand Up @@ -24,6 +24,8 @@

#include <deal.II/grid/tria_accessor.h>

#include <deal.II/lac/read_vector.h>

#include <boost/container/small_vector.hpp>

#include <set>
Expand Down Expand Up @@ -1589,11 +1591,11 @@ class DoFCellAccessor : public DoFAccessor<dimension_,
* caller to assure that the types of the numbers stored in input and output
* vectors are compatible and with similar accuracy.
*/
template <class InputVector, typename ForwardIterator>
template <typename Number, typename ForwardIterator>
void
get_dof_values(const InputVector &values,
ForwardIterator local_values_begin,
ForwardIterator local_values_end) const;
get_dof_values(const ReadVector<Number> &values,
ForwardIterator local_values_begin,
ForwardIterator local_values_end) const;

/**
* Collect the values of the given vector restricted to the dofs of this cell
Expand Down Expand Up @@ -1683,12 +1685,12 @@ class DoFCellAccessor : public DoFAccessor<dimension_,
* interpolation is presently only provided for cells by the finite element
* classes.
*/
template <class InputVector, typename number>
template <typename Number>
void
get_interpolated_dof_values(
const InputVector & values,
Vector<number> & interpolated_values,
const types::fe_index fe_index = numbers::invalid_fe_index) const;
const ReadVector<Number> &values,
Vector<Number> & interpolated_values,
const types::fe_index fe_index = numbers::invalid_fe_index) const;

/**
* This function is the counterpart to get_interpolated_dof_values(): you
Expand Down
24 changes: 15 additions & 9 deletions include/deal.II/dofs/dof_accessor.templates.h
Original file line number Diff line number Diff line change
Expand Up @@ -2524,12 +2524,12 @@ DoFCellAccessor<dimension_, space_dimension_, level_dof_access>::get_dof_values(


template <int dimension_, int space_dimension_, bool level_dof_access>
template <class InputVector, typename ForwardIterator>
template <typename Number, typename ForwardIterator>
inline void
DoFCellAccessor<dimension_, space_dimension_, level_dof_access>::get_dof_values(
const InputVector &values,
ForwardIterator local_values_begin,
ForwardIterator local_values_end) const
const ReadVector<Number> &values,
ForwardIterator local_values_begin,
ForwardIterator local_values_end) const
{
(void)local_values_end;
Assert(this->is_artificial() == false,
Expand All @@ -2547,11 +2547,17 @@ DoFCellAccessor<dimension_, space_dimension_, level_dof_access>::get_dof_values(
dof_indices(this->get_fe().n_dofs_per_cell());
internal::DoFAccessorImplementation::get_cell_dof_indices(
*this, dof_indices, this->active_fe_index());
dealii::internal::DoFAccessorImplementation::Implementation::
extract_subvector_to(values,
dof_indices.data(),
dof_indices.data() + dof_indices.size(),
local_values_begin);

boost::container::small_vector<Number, 27> values_temp(local_values_end -
local_values_begin);
auto view = make_array_view(values_temp.begin(), values_temp.end());
values.extract_subvector_to(make_array_view(dof_indices.begin(),
dof_indices.end()),
view);
using view_type = std::remove_reference_t<decltype(*local_values_begin)>;
ArrayView<view_type> values_view2(&*local_values_begin,
local_values_end - local_values_begin);
std::copy(values_temp.begin(), values_temp.end(), values_view2.begin());
}


Expand Down