16#ifndef dealii_solver_idr_h
17#define dealii_solver_idr_h
44 namespace SolverIDRImplementation
50 template <
typename VectorType>
78 operator()(
const unsigned int i,
const VectorType &temp);
89 std::vector<typename VectorMemory<VectorType>::Pointer>
data;
116template <
class VectorType = Vector<
double>>
132 const unsigned int s;
157 template <
typename MatrixType,
typename PreconditionerType>
161 const VectorType & b,
162 const PreconditionerType &preconditioner);
172 const VectorType & x,
173 const VectorType & r,
174 const VectorType & d)
const;
191 namespace SolverIDRImplementation
193 template <
class VectorType>
202 template <
class VectorType>
214 template <
class VectorType>
217 const VectorType & temp)
220 if (data[i] ==
nullptr)
223 data[i]->reinit(temp);
232template <
class VectorType>
235 const AdditionalData & data)
237 , additional_data(data)
242template <
class VectorType>
245 , additional_data(data)
250template <
class VectorType>
255 const VectorType &)
const
260template <
class VectorType>
261template <
typename MatrixType,
typename PreconditionerType>
265 const VectorType & b,
266 const PreconditionerType &preconditioner)
271 unsigned int step = 0;
273 const unsigned int s = additional_data.s;
282 VectorType &r = *r_pointer;
283 VectorType &v = *v_pointer;
284 VectorType &vhat = *vhat_pointer;
285 VectorType &uhat = *uhat_pointer;
286 VectorType &ghat = *ghat_pointer;
290 vhat.reinit(x,
true);
291 uhat.reinit(x,
true);
292 ghat.reinit(x,
true);
296 r.sadd(-1.0, 1.0, b);
299 double res = r.l2_norm();
300 iteration_state = this->iteration_status(step, res, x);
313 std::normal_distribution<> normal_distribution(0.0, 1.0);
314 for (
unsigned int i = 0; i < s; ++i)
316 VectorType &tmp_g = G(i, x);
317 VectorType &tmp_u = U(i, x);
326 VectorType &tmp_q = Q(i, x);
329 for (
auto indx : tmp_q.locally_owned_elements())
330 tmp_q(indx) = normal_distribution(rng);
336 for (
unsigned int j = 0; j < i; ++j)
339 v *= (v * tmp_q) / (tmp_q * tmp_q);
344 tmp_q *= 1.0 / tmp_q.l2_norm();
351 bool early_exit =
false;
360 for (
unsigned int i = 0; i < s; ++i)
364 for (
unsigned int k = 0; k < s; ++k)
371 std::vector<unsigned int> indices;
373 for (
unsigned int i = k; i < s; ++i, ++j)
375 indices.push_back(i);
378 Mk.extract_submatrix_from(M, indices, indices);
382 Mk_inv.vmult(gamma, phik);
389 for (
unsigned int i = k; i < s; ++i, ++j)
390 v.add(-1.0 *
gamma(j), G[i]);
391 preconditioner.vmult(vhat, v);
396 for (
unsigned int i = k; i < s; ++i, ++j)
397 uhat.add(
gamma(j), U[i]);
404 for (
unsigned int i = 0; i < k; ++i)
406 double alpha = (Q[i] * ghat) / M(i, i);
407 ghat.add(-alpha, G[i]);
408 uhat.add(-alpha, U[i]);
414 for (
unsigned int i = k; i < s; ++i)
415 M(i, k) = Q[i] * G[k];
420 double beta = phi(k) / M(k, k);
421 r.add(-1.0 * beta, G[k]);
424 print_vectors(step, x, r, U[k]);
430 iteration_state = this->iteration_status(step, res, x);
440 for (
unsigned int i = 0; i < k + 1; ++i)
442 for (
unsigned int i = k + 1; i < s; ++i)
443 phi(i) -= beta * M(i, k);
447 if (early_exit ==
true)
451 preconditioner.vmult(vhat, r);
454 omega = (v * r) / (v * v);
456 r.add(-1.0 * omega, v);
459 print_vectors(step, x, r, vhat);
463 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())
AdditionalData additional_data
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_NAMESPACE_CLOSE
#define Assert(cond, exc)
#define AssertIndexRange(index, range)
static ::ExceptionBase & ExcNotInitialized()
#define AssertThrow(cond, exc)
long double gamma(const unsigned int n)
AdditionalData(const unsigned int s=2)