16#ifndef dealii_solver_idr_h
17#define dealii_solver_idr_h
44 namespace SolverIDRImplementation
50 template <
typename VectorType>
77 operator()(
const unsigned int i,
const VectorType &temp);
88 std::vector<typename VectorMemory<VectorType>::Pointer>
data;
115template <
class VectorType = Vector<
double>>
131 const unsigned int s;
156 template <
typename MatrixType,
typename PreconditionerType>
160 const VectorType &
b,
161 const PreconditionerType &preconditioner);
171 const VectorType & x,
172 const VectorType & r,
173 const VectorType &
d)
const;
190 namespace SolverIDRImplementation
192 template <
class VectorType>
201 template <
class VectorType>
213 template <
class VectorType>
216 const VectorType & temp)
219 if (data[i] ==
nullptr)
222 data[i]->reinit(temp);
231template <
class VectorType>
234 const AdditionalData & data)
236 , additional_data(data)
241template <
class VectorType>
244 , additional_data(data)
249template <
class VectorType>
254 const VectorType &)
const
259template <
class VectorType>
260template <
typename MatrixType,
typename PreconditionerType>
264 const VectorType &
b,
265 const PreconditionerType &preconditioner)
270 unsigned int step = 0;
272 const unsigned int s = additional_data.s;
281 VectorType &r = *r_pointer;
282 VectorType &v = *v_pointer;
283 VectorType &vhat = *vhat_pointer;
284 VectorType &uhat = *uhat_pointer;
285 VectorType &ghat = *ghat_pointer;
289 vhat.reinit(x,
true);
290 uhat.reinit(x,
true);
291 ghat.reinit(x,
true);
295 r.sadd(-1.0, 1.0,
b);
298 double res = r.l2_norm();
299 iteration_state = this->iteration_status(step, res, x);
312 std::normal_distribution<> normal_distribution(0.0, 1.0);
313 for (
unsigned int i = 0; i < s; ++i)
315 VectorType &tmp_g = G(i, x);
316 VectorType &tmp_u =
U(i, x);
325 VectorType &tmp_q = Q(i, x);
328 for (
auto indx : tmp_q.locally_owned_elements())
329 tmp_q(indx) = normal_distribution(rng);
335 for (
unsigned int j = 0; j < i; ++j)
338 v *= (v * tmp_q) / (tmp_q * tmp_q);
343 tmp_q *= 1.0 / tmp_q.l2_norm();
350 bool early_exit =
false;
359 for (
unsigned int i = 0; i < s; ++i)
363 for (
unsigned int k = 0; k < s; ++k)
370 std::vector<unsigned int> indices;
372 for (
unsigned int i = k; i < s; ++i, ++j)
374 indices.push_back(i);
377 Mk.extract_submatrix_from(M, indices, indices);
381 Mk_inv.vmult(
gamma, phik);
388 for (
unsigned int i = k; i < s; ++i, ++j)
389 v.add(-1.0 *
gamma(j), G[i]);
390 preconditioner.vmult(vhat, v);
395 for (
unsigned int i = k; i < s; ++i, ++j)
403 for (
unsigned int i = 0; i < k; ++i)
405 double alpha = (Q[i] * ghat) / M(i, i);
406 ghat.add(-alpha, G[i]);
407 uhat.add(-alpha,
U[i]);
413 for (
unsigned int i = k; i < s; ++i)
414 M(i, k) = Q[i] * G[k];
419 double beta = phi(k) / M(k, k);
420 r.add(-1.0 * beta, G[k]);
423 print_vectors(step, x, r,
U[k]);
429 iteration_state = this->iteration_status(step, res, x);
439 for (
unsigned int i = 0; i < k + 1; ++i)
441 for (
unsigned int i = k + 1; i < s; ++i)
442 phi(i) -= beta * M(i, k);
446 if (early_exit ==
true)
450 preconditioner.vmult(vhat, r);
453 omega = (v * r) / (v * v);
455 r.add(-1.0 * omega, v);
458 print_vectors(step, x, r, vhat);
462 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)
SymmetricTensor< 2, dim, Number > b(const Tensor< 2, dim, Number > &F)
SymmetricTensor< 2, dim, Number > d(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
long double gamma(const unsigned int n)
AdditionalData(const unsigned int s=2)