15#ifndef dealii_solver_idr_h
16#define dealii_solver_idr_h
46 namespace SolverIDRImplementation
52 template <
typename VectorType>
80 operator()(
const unsigned int i,
const VectorType &temp);
91 std::vector<typename VectorMemory<VectorType>::Pointer>
data;
118template <
typename VectorType = Vector<
double>>
135 const unsigned int s;
160 template <
typename MatrixType,
typename PreconditionerType>
164 void solve(const MatrixType &A,
167 const PreconditionerType &preconditioner);
176 print_vectors(const
unsigned int step,
179 const VectorType &d) const;
196 namespace SolverIDRImplementation
198 template <
typename VectorType>
207 template <
typename VectorType>
219 template <
typename VectorType>
222 const VectorType &temp)
225 if (
data[i] ==
nullptr)
228 data[i]->reinit(temp,
true);
235 template <
typename VectorType,
236 std::enable_if_t<!IsBlockVector<VectorType>::value, VectorType>
239 n_blocks(
const VectorType &)
246 template <
typename VectorType,
247 std::enable_if_t<IsBlockVector<VectorType>::value, VectorType> * =
250 n_blocks(
const VectorType &vector)
252 return vector.n_blocks();
257 template <
typename VectorType,
258 std::enable_if_t<!IsBlockVector<VectorType>::value, VectorType>
261 block(VectorType &vector,
const unsigned int b)
270 template <
typename VectorType,
271 std::enable_if_t<IsBlockVector<VectorType>::value, VectorType> * =
273 typename VectorType::BlockType &
274 block(VectorType &vector,
const unsigned int b)
276 return vector.block(b);
284template <
typename VectorType>
288 const AdditionalData &
data)
290 , additional_data(
data)
295template <
typename VectorType>
299 , additional_data(
data)
304template <
typename VectorType>
309 const VectorType &)
const
314template <
typename VectorType>
316template <
typename MatrixType,
typename PreconditionerType>
323 const PreconditionerType &preconditioner)
328 unsigned int step = 0;
330 const unsigned int s = additional_data.
s;
343 uhat.reinit(x,
true);
347 r.sadd(-1.0, 1.0, b);
349 using value_type =
typename VectorType::value_type;
353 real_type res = r.l2_norm();
354 iteration_state = this->iteration_status(step, res, x);
367 std::normal_distribution<> normal_distribution(0.0, 1.0);
368 for (
unsigned int i = 0; i < s; ++i)
382 for (
unsigned int b = 0;
383 b < internal::SolverIDRImplementation::n_blocks(tmp_q);
385 for (
auto indx :
internal::SolverIDRImplementation::block(tmp_q,
b)
386 .locally_owned_elements())
387 internal::SolverIDRImplementation::block(tmp_q,
b)(indx) =
388 normal_distribution(rng);
394 for (
unsigned int j = 0; j < i; ++j)
397 v *= (v * tmp_q) / (tmp_q * tmp_q);
402 tmp_q *= 1.0 / tmp_q.l2_norm();
407 value_type omega = 1.;
409 bool early_exit =
false;
418 for (
unsigned int i = 0; i < s; ++i)
422 for (
unsigned int k = 0; k < s; ++k)
429 std::vector<unsigned int> indices;
431 for (
unsigned int i = k; i < s; ++i, ++j)
433 indices.push_back(i);
436 Mk.extract_submatrix_from(M, indices, indices);
440 Mk_inv.vmult(gamma, phik);
447 for (
unsigned int i = k, j = 0; i < s; ++i, ++j)
448 v.add(-
gamma(j), G[i]);
451 preconditioner.vmult(uhat, v);
455 uhat.sadd(omega,
gamma(0), U[k]);
456 for (
unsigned int i = k + 1, j = 1; i < s; ++i, ++j)
457 uhat.add(
gamma(j), U[i]);
468 value_type alpha = Q[0] * G[k] / M(0, 0);
469 for (
unsigned int i = 1; i < k; ++i)
471 const value_type alpha_old = alpha;
472 alpha = G[k].add_and_dot(-alpha, G[i - 1], Q[i]) / M(i, i);
476 uhat.add(-alpha_old, U[i - 1], -alpha, U[i]);
478 M(k, k) = G[k].add_and_dot(-alpha, G[k - 1], Q[k]);
480 uhat.add(-alpha, U[k - 1]);
483 M(k, k) = G[k] * Q[k];
488 for (
unsigned int i = k + 1; i < s; ++i)
489 M(i, k) = Q[i] * G[k];
493 const value_type beta = phi(k) / M(k, k);
497 print_vectors(step, x, r, U[k]);
503 iteration_state = this->iteration_status(step, res, x);
513 for (
unsigned int i = 0; i < k + 1; ++i)
515 for (
unsigned int i = k + 1; i < s; ++i)
516 phi(i) -= beta * M(i, k);
520 if (early_exit ==
true)
524 preconditioner.vmult(uhat, r);
527 omega = (v * r) / (v * v);
529 res =
std::sqrt(r.add_and_dot(-1.0 * omega, v, r));
532 print_vectors(step, x, r, uhat);
535 iteration_state = this->iteration_status(step, res, x);
@ iterate
Continue iteration.
@ success
Stop iteration, goal reached.
virtual void print_vectors(const unsigned int step, const VectorType &x, const VectorType &r, const VectorType &d) const
SolverIDR(SolverControl &cn, const AdditionalData &data=AdditionalData())
void solve(const MatrixType &A, VectorType &x, const VectorType &b, const PreconditionerType &preconditioner)
SolverIDR(SolverControl &cn, VectorMemory< VectorType > &mem, const AdditionalData &data=AdditionalData())
virtual ~SolverIDR() override=default
TmpVectors(const unsigned int s_param, VectorMemory< VectorType > &vmem)
VectorType & operator()(const unsigned int i, const VectorType &temp)
VectorType & operator[](const unsigned int i) const
VectorMemory< VectorType > & mem
std::vector< typename VectorMemory< VectorType >::Pointer > data
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_CXX20_REQUIRES(condition)
#define DEAL_II_NAMESPACE_CLOSE
#define Assert(cond, exc)
#define AssertDimension(dim1, dim2)
#define AssertIndexRange(index, range)
static ::ExceptionBase & ExcNotInitialized()
#define AssertThrow(cond, exc)
std::vector< index_type > data
Tpetra::Vector< Number, LO, GO, NodeType< MemorySpace > > VectorType
SymmetricTensor< 2, dim, Number > b(const Tensor< 2, dim, Number > &F)
long double gamma(const unsigned int n)
::VectorizedArray< Number, width > sqrt(const ::VectorizedArray< Number, width > &)
AdditionalData(const unsigned int s=2)