Skip to content

Commit

Permalink
Merge pull request #13076 from peterrum/gather_evaluate_mixed
Browse files Browse the repository at this point in the history
MatrixFree::gather_evalute(): allow mixed numbers
  • Loading branch information
kronbichler committed Dec 15, 2021
2 parents 9a90f85 + 0178acc commit 0acb31b
Show file tree
Hide file tree
Showing 5 changed files with 318 additions and 39 deletions.
30 changes: 14 additions & 16 deletions include/deal.II/matrix_free/evaluation_kernels.h
Original file line number Diff line number Diff line change
Expand Up @@ -3163,9 +3163,9 @@ namespace internal
const unsigned int n_components,
const EvaluationFlags::EvaluationFlags evaluation_flag,
typename Processor::Number2_ * global_vector_ptr,
const std::vector<ArrayView<const typename Processor::Number_>> *sm_ptr,
const EvaluationData & fe_eval,
typename Processor::VectorizedArrayType_ * temp1)
const std::vector<ArrayView<const typename Processor::Number2_>> *sm_ptr,
const EvaluationData & fe_eval,
typename Processor::VectorizedArrayType_ * temp1)
{
constexpr int dim = Processor::dim_;
constexpr int fe_degree = Processor::fe_degree_;
Expand Down Expand Up @@ -3530,7 +3530,7 @@ namespace internal
dof_info
.dof_indices_contiguous_sm[dof_access_index]
[cell * n_lanes + v];
vector_ptrs[v] = const_cast<Number *>(
vector_ptrs[v] = const_cast<Number2_ *>(
sm_ptr->operator[](temp.first).data() +
temp.second + index_offset);
}
Expand All @@ -3555,7 +3555,7 @@ namespace internal
dof_info
.dof_indices_contiguous_sm[dof_access_index]
[cells[v]];
vector_ptrs[v] = const_cast<Number *>(
vector_ptrs[v] = const_cast<Number2_ *>(
sm_ptr->operator[](temp.first).data() +
temp.second + index_offset);
}
Expand Down Expand Up @@ -3682,18 +3682,17 @@ namespace internal



template <int dim,
typename Number,
typename VectorizedArrayType,
typename Number2 = Number>
template <int dim, typename Number2, typename VectorizedArrayType>
struct FEFaceEvaluationImplGatherEvaluateSelector
{
using Number = typename VectorizedArrayType::value_type;

template <int fe_degree, int n_q_points_1d>
static bool
run(const unsigned int n_components,
const EvaluationFlags::EvaluationFlags evaluation_flag,
const Number2 * src_ptr,
const std::vector<ArrayView<const Number>> * sm_ptr,
const std::vector<ArrayView<const Number2>> * sm_ptr,
FEEvaluationData<dim, VectorizedArrayType, true> &fe_eval)
{
Assert(fe_degree > -1, ExcInternalError());
Expand Down Expand Up @@ -3806,7 +3805,7 @@ namespace internal
supports(
const EvaluationFlags::EvaluationFlags evaluation_flag,
const MatrixFreeFunctions::ShapeInfo<VectorizedArrayType> &shape_info,
const Number * vector_ptr,
const Number2 * vector_ptr,
MatrixFreeFunctions::DoFInfo::IndexStorageVariants storage)
{
const unsigned int fe_degree = shape_info.data.front().fe_degree;
Expand Down Expand Up @@ -3884,7 +3883,7 @@ namespace internal
hermite_grad(T0 & temp_1,
T0 & temp_2,
const T1 &src_ptr_1,
const T2 &src_ptr_2,
const T1 &src_ptr_2,
const T2 &grad_weight)
{
// case 3a)
Expand All @@ -3904,12 +3903,11 @@ namespace internal



template <int dim,
typename Number,
typename VectorizedArrayType,
typename Number2 = Number>
template <int dim, typename Number2, typename VectorizedArrayType>
struct FEFaceEvaluationImplIntegrateScatterSelector
{
using Number = typename VectorizedArrayType::value_type;

template <int fe_degree, int n_q_points_1d>
static bool
run(const unsigned int n_components,
Expand Down
50 changes: 28 additions & 22 deletions include/deal.II/matrix_free/fe_evaluation.h
Original file line number Diff line number Diff line change
Expand Up @@ -8056,35 +8056,38 @@ FEFaceEvaluation<dim,
this->data->data.front().n_q_points_1d) &&
internal::FEFaceEvaluationImplGatherEvaluateSelector<
dim,
Number,
VectorizedArrayType>::supports(evaluation_flag,
*this->data,
internal::get_beginning<Number>(
input_vector),
this->dof_info->index_storage_variants
[this->dof_access_index][this->cell]))
typename VectorType::value_type,
VectorizedArrayType>::
supports(evaluation_flag,
*this->data,
internal::get_beginning<typename VectorType::value_type>(
input_vector),
this->dof_info->index_storage_variants[this->dof_access_index]
[this->cell]))
{
if (fe_degree > -1)
{
internal::FEFaceEvaluationImplGatherEvaluateSelector<
dim,
Number,
typename VectorType::value_type,
VectorizedArrayType>::template run<fe_degree,
n_q_points_1d>(
n_components,
evaluation_flag,
internal::get_beginning<Number>(input_vector),
internal::get_beginning<typename VectorType::value_type>(
input_vector),
shared_vector_data,
*this);
}
else
{
internal::FEFaceEvaluationGatherFactory<
dim,
Number,
typename VectorType::value_type,
VectorizedArrayType>::evaluate(n_components,
evaluation_flag,
internal::get_beginning<Number>(
internal::get_beginning<
typename VectorType::value_type>(
input_vector),
shared_vector_data,
*this);
Expand Down Expand Up @@ -8170,35 +8173,38 @@ FEFaceEvaluation<dim,
this->data->data.front().n_q_points_1d) &&
internal::FEFaceEvaluationImplGatherEvaluateSelector<
dim,
Number,
VectorizedArrayType>::supports(integration_flag,
*this->data,
internal::get_beginning<Number>(
destination),
this->dof_info->index_storage_variants
[this->dof_access_index][this->cell]))
typename VectorType::value_type,
VectorizedArrayType>::
supports(integration_flag,
*this->data,
internal::get_beginning<typename VectorType::value_type>(
destination),
this->dof_info->index_storage_variants[this->dof_access_index]
[this->cell]))
{
if (fe_degree > -1)
{
internal::FEFaceEvaluationImplIntegrateScatterSelector<
dim,
Number,
typename VectorType::value_type,
VectorizedArrayType>::template run<fe_degree,
n_q_points_1d>(
n_components,
integration_flag,
internal::get_beginning<Number>(destination),
internal::get_beginning<typename VectorType::value_type>(
destination),
shared_vector_data,
*this);
}
else
{
internal::FEFaceEvaluationGatherFactory<
dim,
Number,
typename VectorType::value_type,
VectorizedArrayType>::integrate(n_components,
integration_flag,
internal::get_beginning<Number>(
internal::get_beginning<
typename VectorType::value_type>(
destination),
shared_vector_data,
*this);
Expand Down
7 changes: 6 additions & 1 deletion source/matrix_free/evaluation_template_factory.inst.in
Original file line number Diff line number Diff line change
Expand Up @@ -25,7 +25,12 @@ for (deal_II_dimension : DIMENSIONS;

template struct dealii::internal::FEFaceEvaluationGatherFactory<
deal_II_dimension,
deal_II_scalar_vectorized::value_type,
double,
deal_II_scalar_vectorized>;

template struct dealii::internal::FEFaceEvaluationGatherFactory<
deal_II_dimension,
float,
deal_II_scalar_vectorized>;

// inverse mass
Expand Down

0 comments on commit 0acb31b

Please sign in to comment.