#endif
+namespace mg
+{
+ /**
+ * A structure with boost::signal objects for optional processing in a
+ * non-nested multigrid solver.
+ *
+ * Similarly to mg::Signals, each signal is called twice: once before and once
+ * after the action is performed. The two calls only differ in the boolean argument @p before, which is true the first time and false the second.
+ *
+ */
+ struct SignalsNonNested
+ {
+ /**
+ * This signal is triggered before and after the call to actual evaluation
+ * function inside RemotePointEvaluation::evaluate_and_process() during
+ * prolongation.
+ */
+ boost::signals2::signal<void(const bool before)> evaluation_prolongation;
+
+ /**
+ * This signal is triggered before and after the call to actual evaluation
+ * function inside RemotePointEvaluation::process_and_evaluate() during
+ * restriction.
+ */
+ boost::signals2::signal<void(const bool before)> evaluation_restriction;
+
+ /**
+ * This signal is triggered before and after the call to
+ * RemotePointEvaluation::evaluate_and_process() used in
+ * MGTwoLevelTransferNonNested::prolongate_and_add(). The difference with the @p evaluation
+ * signal is that also the communication phase is included.
+ */
+ boost::signals2::signal<void(const bool before)> evaluate_and_process;
+
+ /**
+ * This signal is triggered before and after the call to
+ * RemotePointEvaluation::process_and_evaluate() used in
+ * MGTwoLevelTransferNonNested::restrict_and_add(). Similarly to the @p evaluate_and_process signal,
+ * also the communication phase is included.
+ */
+ boost::signals2::signal<void(const bool before)> process_and_evaluate;
+ };
+} // namespace mg
/**
* Global coarsening utility functions.
private:
using VectorizedArrayType = VectorizedArray<Number, 1>;
+ mg::SignalsNonNested signals_non_nested;
+
public:
/**
* AdditionalData structure that can be used to tweak parameters
std::size_t
memory_consumption() const override;
+ /**
+ * Connect a function to mg::SignalsNonNested::evaluation_restriction.
+ *
+ */
+ boost::signals2::connection
+ connect_evaluation_prolongation(const std::function<void(const bool)> &slot);
+
+ /**
+ * Connect a function to mg::SignalsNonNested::evaluation_restriction.
+ *
+ */
+ boost::signals2::connection
+ connect_evaluation_restriction(const std::function<void(const bool)> &slot);
+
+ /**
+ * Connect a function to mg::SignalsNonNested::evaluate_and_process.
+ *
+ */
+ boost::signals2::connection
+ connect_evaluate_and_process(const std::function<void(const bool)> &slot);
+
+ /**
+ * Connect a function to mg::SignalsNonNested::process_and_evaluate.
+ *
+ */
+ boost::signals2::connection
+ connect_process_and_evaluate(const std::function<void(const bool)> &slot);
+
protected:
AdditionalData additional_data;
/**
std::vector<value_type> buffer;
const auto evaluation_function = [&](auto &values, const auto &cell_data) {
+ this->signals_non_nested.evaluation_prolongation(true);
std::vector<Number> solution_values;
FEPointEvaluation<n_components, dim, dim, Number> evaluator(*mapping_info,
values[q + cell_data.reference_point_ptrs[cell]] =
evaluator.get_value(q);
}
+ this->signals_non_nested.evaluation_prolongation(false);
};
+ this->signals_non_nested.evaluate_and_process(true);
rpe.template evaluate_and_process<value_type>(evaluation_point_results,
buffer,
evaluation_function);
+ this->signals_non_nested.evaluate_and_process(false);
// Weight operator in case some points are owned by multiple cells.
if (rpe.is_map_unique() == false)
const auto evaluation_function = [&](const auto &values,
const auto &cell_data) {
+ this->signals_non_nested.evaluation_restriction(true);
std::vector<Number> solution_values;
FEPointEvaluation<n_components, dim, dim, Number> evaluator(*mapping_info,
*fe_coarse);
solution_values.size(),
true);
}
+ this->signals_non_nested.evaluation_restriction(false);
};
+ this->signals_non_nested.process_and_evaluate(true);
rpe.template process_and_evaluate<value_type>(evaluation_point_results,
buffer,
evaluation_function);
+ this->signals_non_nested.process_and_evaluate(false);
}
}
+
+template <int dim, typename Number>
+boost::signals2::connection
+MGTwoLevelTransferNonNested<dim, LinearAlgebra::distributed::Vector<Number>>::
+ connect_evaluation_prolongation(const std::function<void(const bool)> &slot)
+{
+ return signals_non_nested.evaluation_prolongation.connect(slot);
+}
+
+
+
+template <int dim, typename Number>
+boost::signals2::connection
+MGTwoLevelTransferNonNested<dim, LinearAlgebra::distributed::Vector<Number>>::
+ connect_evaluation_restriction(const std::function<void(const bool)> &slot)
+{
+ return signals_non_nested.evaluation_restriction.connect(slot);
+}
+
+
+
+template <int dim, typename Number>
+boost::signals2::connection
+MGTwoLevelTransferNonNested<dim, LinearAlgebra::distributed::Vector<Number>>::
+ connect_evaluate_and_process(const std::function<void(const bool)> &slot)
+{
+ return signals_non_nested.evaluate_and_process.connect(slot);
+}
+
+
+
+template <int dim, typename Number>
+boost::signals2::connection
+MGTwoLevelTransferNonNested<dim, LinearAlgebra::distributed::Vector<Number>>::
+ connect_process_and_evaluate(const std::function<void(const bool)> &slot)
+{
+ return signals_non_nested.process_and_evaluate.connect(slot);
+}
+
+
DEAL_II_NAMESPACE_CLOSE
#endif