16 #ifndef dealii_constrained_linear_operator_h 17 #define dealii_constrained_linear_operator_h 19 #include <deal.II/lac/affine_constraints.h> 20 #include <deal.II/lac/linear_operator.h> 21 #include <deal.II/lac/packaged_operation.h> 24 DEAL_II_NAMESPACE_OPEN
62 template <
typename Range,
typename Domain,
typename Payload>
70 return_op.
vmult_add = [&constraints](Range &v,
const Domain &u) {
72 ::
ExcMessage(
"The domain and range vectors must be different " 73 "storage locations"));
79 const auto &locally_owned_elements = v.locally_owned_elements();
80 for (
const auto &line : constraints.
get_lines())
82 const auto i = line.index;
83 if (locally_owned_elements.is_element(i))
86 const auto &entries = line.entries;
89 const auto pos = entries[j].first;
90 v(i) += u(pos) * entries[j].second;
98 return_op.
Tvmult_add = [&constraints](Domain &v,
const Range &u) {
100 ::
ExcMessage(
"The domain and range vectors must be different " 101 "storage locations"));
107 const auto &locally_owned_elements = v.locally_owned_elements();
108 for (
const auto &line : constraints.
get_lines())
110 const auto i = line.index;
112 if (locally_owned_elements.is_element(i))
117 const auto &entries = line.entries;
120 const auto pos = entries[j].first;
121 if (locally_owned_elements.is_element(pos))
122 v(pos) += u(i) * entries[j].second;
130 const auto vmult_add = return_op.
vmult_add;
131 return_op.
vmult = [vmult_add](Range &v,
const Domain &u) {
138 return_op.
Tvmult = [Tvmult_add](Domain &v,
const Range &u) {
158 template <
typename Range,
typename Domain,
typename Payload>
166 return_op.
vmult_add = [&constraints](Range &v,
const Domain &u) {
167 const auto &locally_owned_elements = v.locally_owned_elements();
168 for (
const auto &line : constraints.
get_lines())
170 const auto i = line.index;
171 if (locally_owned_elements.is_element(i))
180 return_op.
Tvmult_add = [&constraints](Domain &v,
const Range &u) {
181 const auto &locally_owned_elements = v.locally_owned_elements();
182 for (
const auto &line : constraints.
get_lines())
184 const auto i = line.index;
185 if (locally_owned_elements.is_element(i))
195 const auto vmult_add = return_op.
vmult_add;
196 return_op.
vmult = [vmult_add](Range &v,
const Domain &u) {
203 return_op.
Tvmult = [Tvmult_add](Domain &v,
const Range &u) {
249 template <
typename Range,
typename Domain,
typename Payload>
258 return Ct * linop * C + Id_c;
296 template <
typename Range,
typename Domain,
typename Payload>
301 const Range & right_hand_side)
307 return_comp.
apply_add = [&constraints, &linop, &right_hand_side](Range &v) {
313 linop.reinit_domain_vector(*k,
false);
316 v += Ct * (right_hand_side - linop * *k);
320 const auto apply_add = return_comp.
apply_add;
321 return_comp.
apply = [apply_add](Range &v) {
331 DEAL_II_NAMESPACE_CLOSE
std::function< void(Range &v)> apply
LinearOperator< Range, Domain, Payload > distribute_constraints_linear_operator(const AffineConstraints< typename Range::value_type > &constraints, const LinearOperator< Range, Domain, Payload > &exemplar)
std::function< void(Range &v, bool omit_zeroing_entries)> reinit_vector
std::function< void(Range &v, bool omit_zeroing_entries)> reinit_range_vector
std::function< void(Domain &v, const Range &u)> Tvmult_add
std::function< void(Domain &v, const Range &u)> Tvmult
const LineRange get_lines() const
LinearOperator< Range, Domain, Payload > constrained_linear_operator(const AffineConstraints< typename Range::value_type > &constraints, const LinearOperator< Range, Domain, Payload > &linop)
LinearOperator< Domain, Range, Payload > transpose_operator(const LinearOperator< Range, Domain, Payload > &op)
static ::ExceptionBase & ExcMessage(std::string arg1)
#define Assert(cond, exc)
std::function< void(Range &v, const Domain &u)> vmult
PackagedOperation< Range > constrained_right_hand_side(const AffineConstraints< typename Range::value_type > &constraints, const LinearOperator< Range, Domain, Payload > &linop, const Range &right_hand_side)
std::function< void(Range &v)> apply_add
LinearOperator< Range, Domain, Payload > project_to_constrained_linear_operator(const AffineConstraints< typename Range::value_type > &constraints, const LinearOperator< Range, Domain, Payload > &exemplar)
unsigned int global_dof_index
std::function< void(Range &v, const Domain &u)> vmult_add
void distribute(VectorType &vec) const
static bool equal(const T *p1, const T *p2)