MFEM  v4.3.0
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
div_free_solver.cpp
Go to the documentation of this file.
1 // Copyright (c) 2010-2021, Lawrence Livermore National Security, LLC. Produced
2 // at the Lawrence Livermore National Laboratory. All Rights reserved. See files
3 // LICENSE and NOTICE for details. LLNL-CODE-806117.
4 //
5 // This file is part of the MFEM library. For more information and source code
6 // availability visit https://mfem.org.
7 //
8 // MFEM is free software; you can redistribute it and/or modify it under the
9 // terms of the BSD-3 license. We welcome feedback and contributions, see file
10 // CONTRIBUTING.md for details.
11 
12 #include "div_free_solver.hpp"
13 
14 using namespace std;
15 using namespace mfem;
16 using namespace blocksolvers;
17 
18 void SetOptions(IterativeSolver& solver, const IterSolveParameters& param)
19 {
20  solver.SetPrintLevel(param.print_level);
21  solver.SetMaxIter(param.max_iter);
22  solver.SetAbsTol(param.abs_tol);
23  solver.SetRelTol(param.rel_tol);
24 }
25 
27  const HypreParMatrix& P)
28 {
29  OperatorPtr R(Rt.Transpose());
30  OperatorPtr RA(ParMult(R.As<HypreParMatrix>(), &A));
31  return ParMult(RA.As<HypreParMatrix>(), &P, true);
32 }
33 
34 void GetRowColumnsRef(const SparseMatrix& A, int row, Array<int>& cols)
35 {
36  cols.MakeRef(const_cast<int*>(A.GetRowColumns(row)), A.RowSize(row));
37 }
38 
40 {
41  int* I = new int[fes.GetNE()+1];
42  copy_n(fes.GetElementToDofTable().GetI(), fes.GetNE()+1, I);
43  Array<int> J(new int[I[fes.GetNE()]], I[fes.GetNE()]);
44  copy_n(fes.GetElementToDofTable().GetJ(), J.Size(), J.begin());
45  fes.AdjustVDofs(J);
46  double* D = new double[J.Size()];
47  fill_n(D, J.Size(), 1.0);
48  return SparseMatrix(I, J, D, fes.GetNE(), fes.GetVSize());
49 }
50 
51 DFSSpaces::DFSSpaces(int order, int num_refine, ParMesh *mesh,
52  const Array<int>& ess_attr, const DFSParameters& param)
53  : hdiv_fec_(order, mesh->Dimension()), l2_fec_(order, mesh->Dimension()),
54  l2_0_fec_(0, mesh->Dimension()), ess_bdr_attr_(ess_attr), level_(0)
55 {
56  if (mesh->GetElement(0)->GetType() == Element::TETRAHEDRON && order)
57  {
58  mfem_error("DFSDataCollector: High order spaces on tetrahedra are not supported");
59  }
60 
61  data_.param = param;
62 
63  if (mesh->Dimension() == 3)
64  {
65  hcurl_fec_.reset(new ND_FECollection(order+1, mesh->Dimension()));
66  }
67  else
68  {
69  hcurl_fec_.reset(new H1_FECollection(order+1, mesh->Dimension()));
70  }
71 
72  all_bdr_attr_.SetSize(ess_attr.Size(), 1);
73  hdiv_fes_.reset(new ParFiniteElementSpace(mesh, &hdiv_fec_));
74  l2_fes_.reset(new ParFiniteElementSpace(mesh, &l2_fec_));
75  coarse_hdiv_fes_.reset(new ParFiniteElementSpace(*hdiv_fes_));
76  coarse_l2_fes_.reset(new ParFiniteElementSpace(*l2_fes_));
77  l2_0_fes_.reset(new ParFiniteElementSpace(mesh, &l2_0_fec_));
78  l2_0_fes_->SetUpdateOperatorType(Operator::MFEM_SPARSEMAT);
79  el_l2dof_.reserve(num_refine+1);
80  el_l2dof_.push_back(ElemToDof(*coarse_l2_fes_));
81 
82  data_.agg_hdivdof.resize(num_refine);
83  data_.agg_l2dof.resize(num_refine);
84  data_.P_hdiv.resize(num_refine, OperatorPtr(Operator::Hypre_ParCSR));
85  data_.P_l2.resize(num_refine, OperatorPtr(Operator::Hypre_ParCSR));
86  data_.Q_l2.resize(num_refine);
87  hdiv_fes_->GetEssentialTrueDofs(ess_attr, data_.coarsest_ess_hdivdofs);
88  data_.C.resize(num_refine+1);
89 
90  hcurl_fes_.reset(new ParFiniteElementSpace(mesh, hcurl_fec_.get()));
91  coarse_hcurl_fes_.reset(new ParFiniteElementSpace(*hcurl_fes_));
92  data_.P_hcurl.resize(num_refine, OperatorPtr(Operator::Hypre_ParCSR));
93 }
94 
96  const SparseMatrix& agg_elem,
97  const SparseMatrix& elem_dof,
98  const HypreParMatrix& dof_truedof,
99  Array<HYPRE_Int>& agg_starts)
100 {
101  OperatorPtr agg_dof(Mult(agg_elem, elem_dof));
102  SparseMatrix& agg_dof_ref = *agg_dof.As<SparseMatrix>();
103  OperatorPtr agg_tdof(dof_truedof.LeftDiagMult(agg_dof_ref, agg_starts));
104  OperatorPtr agg_tdof_T(agg_tdof.As<HypreParMatrix>()->Transpose());
105  SparseMatrix tdof_agg, is_shared;
106  HYPRE_Int* trash;
107  agg_tdof_T.As<HypreParMatrix>()->GetDiag(tdof_agg);
108  agg_tdof_T.As<HypreParMatrix>()->GetOffd(is_shared, trash);
109 
110  int * I = new int [tdof_agg.NumRows()+1]();
111  int * J = new int[tdof_agg.NumNonZeroElems()];
112 
113  Array<int> is_bdr;
114  FiniteElementSpace::ListToMarker(bdr_truedofs, tdof_agg.NumRows(), is_bdr);
115 
116  int counter = 0;
117  for (int i = 0; i < tdof_agg.NumRows(); ++i)
118  {
119  bool agg_bdr = is_bdr[i] || is_shared.RowSize(i) || tdof_agg.RowSize(i)>1;
120  if (agg_bdr) { I[i+1] = I[i]; continue; }
121  I[i+1] = I[i] + 1;
122  J[counter++] = tdof_agg.GetRowColumns(i)[0];
123  }
124 
125  double * D = new double[I[tdof_agg.NumRows()]];
126  std::fill_n(D, I[tdof_agg.NumRows()], 1.0);
127 
128  SparseMatrix intdof_agg(I, J, D, tdof_agg.NumRows(), tdof_agg.NumCols());
129  return Transpose(intdof_agg);
130 }
131 
132 void DFSSpaces::MakeDofRelationTables(int level)
133 {
134  Array<HYPRE_Int> agg_starts(Array<HYPRE_Int>(l2_0_fes_->GetDofOffsets(), 2));
135  auto& elem_agg = (const SparseMatrix&)*l2_0_fes_->GetUpdateOperator();
136  OperatorPtr agg_elem(Transpose(elem_agg));
137  SparseMatrix& agg_el = *agg_elem.As<SparseMatrix>();
138 
139  el_l2dof_.push_back(ElemToDof(*l2_fes_));
140  data_.agg_l2dof[level].Reset(Mult(agg_el, el_l2dof_[level+1]));
141 
142  Array<int> bdr_tdofs;
143  hdiv_fes_->GetEssentialTrueDofs(all_bdr_attr_, bdr_tdofs);
144  auto tmp = AggToInteriorDof(bdr_tdofs, agg_el, ElemToDof(*hdiv_fes_),
145  *hdiv_fes_->Dof_TrueDof_Matrix(), agg_starts);
146  data_.agg_hdivdof[level].Reset(tmp);
147 }
148 
150 {
151  auto GetP = [this](OperatorPtr& P, unique_ptr<ParFiniteElementSpace>& cfes,
152  ParFiniteElementSpace& fes, bool remove_zero)
153  {
154  fes.Update();
155  fes.GetTrueTransferOperator(*cfes, P);
156  if (remove_zero)
157  {
158  P.As<HypreParMatrix>()->DropSmallEntries(1e-16);
159  }
160  (level_ < (int)data_.P_l2.size()-1) ? cfes->Update() : cfes.reset();
161  };
162 
163  GetP(data_.P_hdiv[level_], coarse_hdiv_fes_, *hdiv_fes_, true);
164  GetP(data_.P_l2[level_], coarse_l2_fes_, *l2_fes_, false);
165  MakeDofRelationTables(level_);
166 
167  GetP(data_.P_hcurl[level_], coarse_hcurl_fes_, *hcurl_fes_, true);
168 
169  Vector trash1(hcurl_fes_->GetVSize()), trash2(hdiv_fes_->GetVSize());
170  ParDiscreteLinearOperator curl(hcurl_fes_.get(), hdiv_fes_.get());
172  curl.Assemble();
173  curl.EliminateTrialDofs(ess_bdr_attr_, trash1, trash2);
174  curl.Finalize();
175  data_.C[level_+1].Reset(curl.ParallelAssemble());
176 
177  ++level_;
178 
179  if (level_ == (int)data_.P_l2.size()) { DataFinalize(); }
180 }
181 
182 void DFSSpaces::DataFinalize()
183 {
184  ParBilinearForm mass(l2_fes_.get());
186  mass.Assemble();
187  mass.Finalize();
188  OperatorPtr W(mass.LoseMat());
189 
190  SparseMatrix P_l2;
191  for (int l = (int)data_.P_l2.size()-1; l >= 0; --l)
192  {
193  data_.P_l2[l].As<HypreParMatrix>()->GetDiag(P_l2);
194  OperatorPtr PT_l2(Transpose(P_l2));
195  auto PTW = Mult(*PT_l2.As<SparseMatrix>(), *W.As<SparseMatrix>());
196  auto cW = Mult(*PTW, P_l2);
197  auto cW_inv = new SymDirectSubBlockSolver(*cW, el_l2dof_[l]);
198  data_.Q_l2[l].Reset(new ProductOperator(cW_inv, PTW, true, true));
199  W.Reset(cW);
200  }
201 
202  l2_0_fes_.reset();
203 }
204 
206  : Solver(B.NumRows()), BBT_solver_(B.GetComm())
207 {
208  OperatorPtr BT(B.Transpose());
209  BBT_.Reset(ParMult(&B, BT.As<HypreParMatrix>()));
210  BBT_.As<HypreParMatrix>()->CopyColStarts();
211 
212  BBT_prec_.Reset(new HypreBoomerAMG(*BBT_.As<HypreParMatrix>()));
213  BBT_prec_.As<HypreBoomerAMG>()->SetPrintLevel(0);
214 
215  SetOptions(BBT_solver_, param);
216  BBT_solver_.SetOperator(*BBT_);
217  BBT_solver_.SetPreconditioner(*BBT_prec_.As<HypreBoomerAMG>());
218 }
219 
221  : Solver(M.NumRows()+B.NumRows()), local_system_(height), offset_(M.NumRows())
222 {
223  local_system_.CopyMN(M, 0, 0);
224  local_system_.CopyMN(B, offset_, 0);
225  local_system_.CopyMNt(B, 0, offset_);
226 
227  local_system_.SetRow(offset_, 0.0);
228  local_system_.SetCol(offset_, 0.0);
229  local_system_(offset_, offset_) = -1.0;
230  local_solver_.SetOperator(local_system_);
231 }
232 
233 void LocalSolver::Mult(const Vector &x, Vector &y) const
234 {
235  const double x0 = x[offset_];
236  const_cast<Vector&>(x)[offset_] = 0.0;
237 
238  y.SetSize(local_system_.NumRows());
239  local_solver_.Mult(x, y);
240 
241  const_cast<Vector&>(x)[offset_] = x0;
242 }
243 
245  const HypreParMatrix& B,
246  const SparseMatrix& agg_hdivdof,
247  const SparseMatrix& agg_l2dof,
248  const HypreParMatrix& P_l2,
249  const HypreParMatrix& Q_l2)
250  : Solver(M.NumRows() + B.NumRows()), agg_hdivdof_(agg_hdivdof),
251  agg_l2dof_(agg_l2dof), solvers_loc_(agg_l2dof.NumRows())
252 {
253  coarse_l2_projector_.Reset(new ProductOperator(&P_l2, &Q_l2, false, false));
254 
255  offsets_loc_.SetSize(3, 0);
256  offsets_.SetSize(3, 0);
257  offsets_[1] = M.NumRows();
258  offsets_[2] = M.NumRows() + B.NumRows();
259 
260  SparseMatrix M_diag, B_diag;
261  M.GetDiag(M_diag);
262  B.GetDiag(B_diag);
263 
264  DenseMatrix B_loc, M_loc;
265 
266  for (int agg = 0; agg < (int)solvers_loc_.size(); agg++)
267  {
268  GetRowColumnsRef(agg_hdivdof_, agg, hdivdofs_loc_);
269  GetRowColumnsRef(agg_l2dof_, agg, l2dofs_loc_);
270  M_loc.SetSize(hdivdofs_loc_.Size(), hdivdofs_loc_.Size());
271  B_loc.SetSize(l2dofs_loc_.Size(), hdivdofs_loc_.Size());
272  M_diag.GetSubMatrix(hdivdofs_loc_, hdivdofs_loc_, M_loc);
273  B_diag.GetSubMatrix(l2dofs_loc_, hdivdofs_loc_, B_loc);
274  solvers_loc_[agg].Reset(new LocalSolver(M_loc, B_loc));
275  }
276 }
277 
278 void SaddleSchwarzSmoother::Mult(const Vector & x, Vector & y) const
279 {
280  y.SetSize(offsets_[2]);
281  y = 0.0;
282 
283  BlockVector blk_y(y.GetData(), offsets_);
284  BlockVector Pi_x(offsets_); // aggregate-wise average free projection of x
285  static_cast<Vector&>(Pi_x) = x;
286 
287  // Right hand side: F_l = F - W_l P_l2[l] (W_{l+1})^{-1} P_l2[l]^T F
288  // This ensures the existence of solutions to the local problems
289  Vector coarse_l2_projection(Pi_x.BlockSize(1));
290  coarse_l2_projector_->MultTranspose(Pi_x.GetBlock(1), coarse_l2_projection);
291 
292  Pi_x.GetBlock(1) -= coarse_l2_projection;
293 
294  for (int agg = 0; agg < (int)solvers_loc_.size(); agg++)
295  {
296  GetRowColumnsRef(agg_hdivdof_, agg, hdivdofs_loc_);
297  GetRowColumnsRef(agg_l2dof_, agg, l2dofs_loc_);
298 
299  offsets_loc_[1] = hdivdofs_loc_.Size();
300  offsets_loc_[2] = offsets_loc_[1]+l2dofs_loc_.Size();
301 
302  BlockVector rhs_loc(offsets_loc_), sol_loc(offsets_loc_);
303  Pi_x.GetBlock(0).GetSubVector(hdivdofs_loc_, rhs_loc.GetBlock(0));
304  Pi_x.GetBlock(1).GetSubVector(l2dofs_loc_, rhs_loc.GetBlock(1));
305 
306  solvers_loc_[agg]->Mult(rhs_loc, sol_loc);
307 
308  blk_y.GetBlock(0).AddElementVector(hdivdofs_loc_, sol_loc.GetBlock(0));
309  blk_y.GetBlock(1).AddElementVector(l2dofs_loc_, sol_loc.GetBlock(1));
310  }
311 
312  coarse_l2_projector_->Mult(blk_y.GetBlock(1), coarse_l2_projection);
313  blk_y.GetBlock(1) -= coarse_l2_projection;
314 }
315 
317  IterSolveParameters param)
318  : DarcySolver(M.NumRows(), B.NumRows()), op_(offsets_), prec_(offsets_),
319  BT_(B.Transpose()), solver_(M.GetComm())
320 {
321  op_.SetBlock(0,0, &M);
322  op_.SetBlock(0,1, BT_.As<HypreParMatrix>());
323  op_.SetBlock(1,0, &B);
324 
325  Vector Md;
326  M.GetDiag(Md);
327  BT_.As<HypreParMatrix>()->InvScaleRows(Md);
328  S_.Reset(ParMult(&B, BT_.As<HypreParMatrix>()));
329  BT_.As<HypreParMatrix>()->ScaleRows(Md);
330 
331  prec_.SetDiagonalBlock(0, new HypreDiagScale(M));
332  prec_.SetDiagonalBlock(1, new HypreBoomerAMG(*S_.As<HypreParMatrix>()));
333  static_cast<HypreBoomerAMG&>(prec_.GetDiagonalBlock(1)).SetPrintLevel(0);
334  prec_.owns_blocks = true;
335 
336  SetOptions(solver_, param);
337  solver_.SetOperator(op_);
338  solver_.SetPreconditioner(prec_);
339 }
340 
341 void BDPMinresSolver::Mult(const Vector & x, Vector & y) const
342 {
343  solver_.Mult(x, y);
344  for (int dof : ess_zero_dofs_) { y[dof] = 0.0; }
345 }
346 
348  const DFSData& data)
349  : DarcySolver(M.NumRows(), B.NumRows()), data_(data), param_(data.param),
350  BT_(B.Transpose()), BBT_solver_(B, param_.BBT_solve_param),
351  ops_offsets_(data.P_l2.size()+1), ops_(ops_offsets_.size()),
352  blk_Ps_(ops_.Size()-1), smoothers_(ops_.Size())
353 {
354  ops_offsets_.back().MakeRef(DarcySolver::offsets_);
355  ops_.Last() = new BlockOperator(ops_offsets_.back());
356  ops_.Last()->SetBlock(0, 0, const_cast<HypreParMatrix*>(&M));
357  ops_.Last()->SetBlock(1, 0, const_cast<HypreParMatrix*>(&B));
358  ops_.Last()->SetBlock(0, 1, BT_.Ptr());
359 
360  for (int l = data.P_l2.size(); l >= 0; --l)
361  {
362  auto& M_f = static_cast<HypreParMatrix&>(ops_[l]->GetBlock(0, 0));
363  auto& B_f = static_cast<HypreParMatrix&>(ops_[l]->GetBlock(1, 0));
364 
365  if (l == 0)
366  {
367  SparseMatrix M_f_diag, B_f_diag;
368  M_f.GetDiag(M_f_diag);
369  B_f.GetDiag(B_f_diag);
370  for (int dof : data.coarsest_ess_hdivdofs)
371  {
372  M_f_diag.EliminateRowCol(dof);
373  B_f_diag.EliminateCol(dof);
374  }
375 
376  const IterSolveParameters& param = param_.coarse_solve_param;
377  auto coarse_solver = new BDPMinresSolver(M_f, B_f, param);
378  if (ops_.Size() > 1)
379  {
380  coarse_solver->SetEssZeroDofs(data.coarsest_ess_hdivdofs);
381  }
382  smoothers_[l] = coarse_solver;
383  continue;
384  }
385 
386  HypreParMatrix& P_hdiv_l = *data.P_hdiv[l-1].As<HypreParMatrix>();
387  HypreParMatrix& P_l2_l = *data.P_l2[l-1].As<HypreParMatrix>();
388  SparseMatrix& agg_hdivdof_l = *data.agg_hdivdof[l-1].As<SparseMatrix>();
389  SparseMatrix& agg_l2dof_l = *data.agg_l2dof[l-1].As<SparseMatrix>();
390  HypreParMatrix& Q_l2_l = *data.Q_l2[l-1].As<HypreParMatrix>();
391  HypreParMatrix* C_l = data.C[l].As<HypreParMatrix>();
392 
393  auto S0 = new SaddleSchwarzSmoother(M_f, B_f, agg_hdivdof_l,
394  agg_l2dof_l, P_l2_l, Q_l2_l);
395  if (param_.coupled_solve)
396  {
397  auto S1 = new BlockDiagonalPreconditioner(ops_offsets_[l]);
398  S1->SetDiagonalBlock(0, new AuxSpaceSmoother(M_f, C_l));
399  S1->owns_blocks = true;
400  smoothers_[l] = new ProductSolver(ops_[l], S0, S1, false, true, true);
401  }
402  else
403  {
404  smoothers_[l] = S0;
405  }
406 
407  HypreParMatrix* M_c = TwoStepsRAP(P_hdiv_l, M_f, P_hdiv_l);
408  HypreParMatrix* B_c = TwoStepsRAP(P_l2_l, B_f, P_hdiv_l);
409 
410  ops_offsets_[l-1].SetSize(3, 0);
411  ops_offsets_[l-1][1] = M_c->NumRows();
412  ops_offsets_[l-1][2] = M_c->NumRows() + B_c->NumRows();
413 
414  blk_Ps_[l-1] = new BlockOperator(ops_offsets_[l], ops_offsets_[l-1]);
415  blk_Ps_[l-1]->SetBlock(0, 0, &P_hdiv_l);
416  blk_Ps_[l-1]->SetBlock(1, 1, &P_l2_l);
417 
418  ops_[l-1] = new BlockOperator(ops_offsets_[l-1]);
419  ops_[l-1]->SetBlock(0, 0, M_c);
420  ops_[l-1]->SetBlock(1, 0, B_c);
421  ops_[l-1]->SetBlock(0, 1, B_c->Transpose());
422  ops_[l-1]->owns_blocks = true;
423  }
424 
425  Array<bool> own_ops(ops_.Size());
426  Array<bool> own_smoothers(smoothers_.Size());
427  Array<bool> own_Ps(blk_Ps_.Size());
428  own_ops = true;
429  own_smoothers = true;
430  own_Ps = true;
431 
432  if (data_.P_l2.size() == 0) { return; }
433 
434  if (param_.coupled_solve)
435  {
436  solver_.Reset(new GMRESSolver(B.GetComm()));
437  solver_.As<GMRESSolver>()->SetOperator(*(ops_.Last()));
438  prec_.Reset(new Multigrid(ops_, smoothers_, blk_Ps_,
439  own_ops, own_smoothers, own_Ps));
440  }
441  else
442  {
443  Array<HypreParMatrix*> ops(data_.P_hcurl.size()+1);
444  Array<Solver*> smoothers(ops.Size());
445  Array<HypreParMatrix*> Ps(data_.P_hcurl.size());
446  own_Ps = false;
447 
448  HypreParMatrix& C_finest = *data.C.back().As<HypreParMatrix>();
449  ops.Last() = TwoStepsRAP(C_finest, M, C_finest);
450  ops.Last()->EliminateZeroRows();
451  ops.Last()->DropSmallEntries(1e-14);
452 
453  solver_.Reset(new CGSolver(B.GetComm()));
454  solver_.As<CGSolver>()->SetOperator(*ops.Last());
455  smoothers.Last() = new HypreSmoother(*ops.Last());
456  static_cast<HypreSmoother*>(smoothers.Last())->SetOperatorSymmetry(true);
457 
458  for (int l = Ps.Size()-1; l >= 0; --l)
459  {
460  Ps[l] = data_.P_hcurl[l].As<HypreParMatrix>();
461  ops[l] = TwoStepsRAP(*Ps[l], *ops[l+1], *Ps[l]);
462  ops[l]->DropSmallEntries(1e-14);
463  smoothers[l] = new HypreSmoother(*ops[l]);
464  static_cast<HypreSmoother*>(smoothers[l])->SetOperatorSymmetry(true);
465  }
466 
467  prec_.Reset(new Multigrid(ops, smoothers, Ps, own_ops, own_smoothers, own_Ps));
468  }
469 
470  solver_.As<IterativeSolver>()->SetPreconditioner(*prec_.As<Solver>());
471  SetOptions(*solver_.As<IterativeSolver>(), param_);
472 }
473 
475 {
476  if (param_.coupled_solve) { return; }
477  for (int i = 0; i < ops_.Size(); ++i)
478  {
479  delete ops_[i];
480  delete smoothers_[i];
481  if (i == ops_.Size() - 1) { break; }
482  delete blk_Ps_[i];
483  }
484 }
485 
486 void DivFreeSolver::SolveParticular(const Vector& rhs, Vector& sol) const
487 {
488  std::vector<Vector> rhss(smoothers_.Size());
489  std::vector<Vector> sols(smoothers_.Size());
490 
491  rhss.back().SetDataAndSize(const_cast<Vector&>(rhs), rhs.Size());
492  sols.back().SetDataAndSize(sol, sol.Size());
493 
494  for (int l = blk_Ps_.Size()-1; l >= 0; --l)
495  {
496  rhss[l].SetSize(blk_Ps_[l]->NumCols());
497  sols[l].SetSize(blk_Ps_[l]->NumCols());
498 
499  sols[l] = 0.0;
500  rhss[l] = 0.0;
501 
502  blk_Ps_[l]->MultTranspose(rhss[l+1], rhss[l]);
503  }
504 
505  for (int l = 0; l < smoothers_.Size(); ++l)
506  {
507  smoothers_[l]->Mult(rhss[l], sols[l]);
508  }
509 
510  for (int l = 0; l < blk_Ps_.Size(); ++l)
511  {
512  Vector P_sol(blk_Ps_[l]->NumRows());
513  blk_Ps_[l]->Mult(sols[l], P_sol);
514  sols[l+1] += P_sol;
515  }
516 }
517 
518 void DivFreeSolver::SolveDivFree(const Vector &rhs, Vector& sol) const
519 {
520  Vector rhs_divfree(data_.C.back()->NumCols());
521  data_.C.back()->MultTranspose(rhs, rhs_divfree);
522 
523  Vector potential_divfree(rhs_divfree.Size());
524  potential_divfree = 0.0;
525  solver_->Mult(rhs_divfree, potential_divfree);
526 
527  data_.C.back()->Mult(potential_divfree, sol);
528 }
529 
530 void DivFreeSolver::SolvePotential(const Vector& rhs, Vector& sol) const
531 {
532  Vector rhs_p(BT_->NumCols());
533  BT_->MultTranspose(rhs, rhs_p);
534  BBT_solver_.Mult(rhs_p, sol);
535 }
536 
537 void DivFreeSolver::Mult(const Vector & x, Vector & y) const
538 {
539  MFEM_VERIFY(x.Size() == offsets_[2], "MLDivFreeSolver: x size is invalid");
540  MFEM_VERIFY(y.Size() == offsets_[2], "MLDivFreeSolver: y size is invalid");
541 
542  if (ops_.Size() == 1) { smoothers_[0]->Mult(x, y); return; }
543 
544  BlockVector blk_y(y, offsets_);
545 
546  BlockVector resid(offsets_);
547  ops_.Last()->Mult(y, resid);
548  add(1.0, x, -1.0, resid, resid);
549 
550  BlockVector correction(offsets_);
551  correction = 0.0;
552 
553  if (param_.coupled_solve)
554  {
555  solver_->Mult(resid, correction);
556  y += correction;
557  }
558  else
559  {
560  StopWatch ch;
561  ch.Start();
562 
563  SolveParticular(resid, correction);
564  blk_y += correction;
565 
566  if (param_.verbose)
567  {
568  cout << "Particular solution found in " << ch.RealTime() << "s.\n";
569  }
570 
571  ch.Clear();
572  ch.Start();
573 
574  ops_.Last()->Mult(y, resid);
575  add(1.0, x, -1.0, resid, resid);
576 
577  SolveDivFree(resid.GetBlock(0), correction.GetBlock(0));
578  blk_y.GetBlock(0) += correction.GetBlock(0);
579 
580  if (param_.verbose)
581  {
582  cout << "Divergence free solution found in " << ch.RealTime() << "s.\n";
583  }
584 
585  ch.Clear();
586  ch.Start();
587 
588  auto M = dynamic_cast<HypreParMatrix&>(ops_.Last()->GetBlock(0, 0));
589  M.Mult(-1.0, correction.GetBlock(0), 1.0, resid.GetBlock(0));
590  SolvePotential(resid.GetBlock(0), correction.GetBlock(1));
591  blk_y.GetBlock(1) += correction.GetBlock(1);
592 
593  if (param_.verbose)
594  {
595  cout << "Scalar potential found in " << ch.RealTime() << "s.\n";
596  }
597  }
598 }
599 
601 {
602  if (ops_.Size() == 1)
603  {
604  return static_cast<BDPMinresSolver*>(smoothers_[0])->GetNumIterations();
605  }
606  return solver_.As<IterativeSolver>()->GetNumIterations();
607 }
int Size() const
Return the logical size of the array.
Definition: array.hpp:134
Multigrid solver class.
Definition: multigrid.hpp:25
int RowSize(const int i) const
Returns the number of elements in row i.
Definition: sparsemat.cpp:311
int GetVSize() const
Return the number of vector dofs, i.e. GetNDofs() x GetVDim().
Definition: fespace.hpp:540
Conjugate gradient method.
Definition: solvers.hpp:316
MPI_Comm GetComm() const
MPI communicator.
Definition: hypre.hpp:472
OpType * As() const
Return the Operator pointer statically cast to a specified OpType. Similar to the method Get()...
Definition: handle.hpp:99
int * GetJ()
Definition: table.hpp:114
void EliminateCol(int col, DiagonalPolicy dpolicy=DIAG_ZERO)
Eliminates the column col from the matrix.
Definition: sparsemat.cpp:1595
void SetCol(int c, const double *col)
Definition: densemat.cpp:1803
A class to handle Vectors in a block fashion.
Definition: blockvector.hpp:30
void SetRow(int r, const double *row)
Definition: densemat.cpp:1788
DivFreeSolver(const HypreParMatrix &M, const HypreParMatrix &B, const DFSData &data)
void SetSize(int s)
Resize the vector to size s.
Definition: vector.hpp:513
int * GetRowColumns(const int row)
Return a pointer to the column indices in a row.
Definition: sparsemat.cpp:358
void Mult(const Table &A, const Table &B, Table &C)
C = A * B (as boolean matrices)
Definition: table.cpp:476
std::vector< OperatorPtr > Q_l2
Pointer to an Operator of a specified type.
Definition: handle.hpp:33
Abstract solver class for Darcy&#39;s flow.
Data type dense matrix using column-major storage.
Definition: densemat.hpp:23
int Size() const
Returns the size of the vector.
Definition: vector.hpp:190
virtual void Mult(const Vector &x, Vector &y) const
Operator application: y=A(x).
virtual void Mult(const Vector &x, Vector &y) const =0
Operator application: y=A(x).
double RealTime()
Definition: tic_toc.cpp:426
Abstract parallel finite element space.
Definition: pfespace.hpp:28
SaddleSchwarzSmoother(const HypreParMatrix &M, const HypreParMatrix &B, const SparseMatrix &agg_hdivdof, const SparseMatrix &agg_l2dof, const HypreParMatrix &P_l2, const HypreParMatrix &Q_l2)
std::vector< OperatorPtr > P_hcurl
double * GetData() const
Return a pointer to the beginning of the Vector data.
Definition: vector.hpp:199
HypreParMatrix * ParMult(const HypreParMatrix *A, const HypreParMatrix *B, bool own_matrix)
Definition: hypre.cpp:2464
void GetSubMatrix(const Array< int > &rows, const Array< int > &cols, DenseMatrix &subm) const
Definition: sparsemat.cpp:2676
void add(const Vector &v1, const Vector &v2, Vector &v)
Definition: vector.cpp:291
void AddDomainInterpolator(DiscreteInterpolator *di)
Adds a domain interpolator. Assumes ownership of di.
virtual void MultTranspose(const Vector &x, Vector &y) const
Action of the transpose operator: y=A^t(x). The default behavior in class Operator is to generate an ...
Definition: operator.hpp:93
IterSolveParameters coarse_solve_param
virtual void Mult(const Vector &x, Vector &y) const
Operator application: y=A(x).
LocalSolver(const DenseMatrix &M, const DenseMatrix &B)
Operator & GetDiagonalBlock(int iblock)
Return a reference to block i,i.
virtual void Mult(const Vector &x, Vector &y) const
Operator application: y=A(x).
ID for class SparseMatrix.
Definition: operator.hpp:255
int GetNE() const
Returns number of elements in the mesh.
Definition: fespace.hpp:567
The BoomerAMG solver in hypre.
Definition: hypre.hpp:1387
Data for the divergence free solver.
virtual void Mult(const Vector &b, Vector &x) const
Operator application: y=A(x).
Definition: solvers.cpp:1453
General product operator: x -&gt; (A*B)(x) = A(B(x)).
Definition: operator.hpp:736
virtual void SetPreconditioner(Solver &pr)
This should be called before SetOperator.
Definition: solvers.hpp:439
void SetPrintLevel(int print_lvl)
Definition: solvers.cpp:71
A class to handle Block diagonal preconditioners in a matrix-free implementation. ...
Data type sparse matrix.
Definition: sparsemat.hpp:41
Wrapper for the block-diagonal-preconditioned MINRES defined in ex5p.cpp.
Jacobi preconditioner in hypre.
Definition: hypre.hpp:1234
SparseMatrix ElemToDof(const ParFiniteElementSpace &fes)
void mfem_error(const char *msg)
Function called when an error is encountered. Used by the macros MFEM_ABORT, MFEM_ASSERT, MFEM_VERIFY.
Definition: error.cpp:153
Timing object.
Definition: tic_toc.hpp:34
void SetMaxIter(int max_it)
Definition: solvers.hpp:100
std::vector< OperatorPtr > P_hdiv
virtual void SetOperator(const Operator &op)
Set/update the solver for the given operator.
HypreParMatrix * Transpose() const
Returns the transpose of *this.
Definition: hypre.cpp:1471
const Element * GetElement(int i) const
Definition: mesh.hpp:942
Parallel smoothers in hypre.
Definition: hypre.hpp:840
virtual void SetOperator(const Operator &op)
Set/update the solver for the given operator.
Definition: densemat.cpp:3268
void GetDiag(Vector &diag) const
Get the local diagonal of the matrix.
Definition: hypre.cpp:1344
virtual void SetOperator(const Operator &op)
Also calls SetOperator for the preconditioner if there is one.
Definition: solvers.cpp:1439
int Dimension() const
Definition: mesh.hpp:911
Operator * Ptr() const
Access the underlying Operator pointer.
Definition: handle.hpp:82
void Start()
Clear the elapsed time and start the stopwatch.
Definition: tic_toc.cpp:411
BDPMinresSolver(HypreParMatrix &M, HypreParMatrix &B, IterSolveParameters param)
void SetAbsTol(double atol)
Definition: solvers.hpp:99
void SetRelTol(double rtol)
Definition: solvers.hpp:98
Abstract base class for iterative solver.
Definition: solvers.hpp:66
int NumRows() const
Get the number of rows (size of output) of the Operator. Synonym with Height().
Definition: operator.hpp:69
OperatorHandle OperatorPtr
Add an alternative name for OperatorHandle – OperatorPtr.
Definition: handle.hpp:207
void SetSize(int nsize)
Change the logical size of the array, keep existing entries.
Definition: array.hpp:674
void EliminateRowCol(int rc, const double sol, Vector &rhs, DiagonalPolicy dpolicy=DIAG_ONE)
Eliminate row rc and column rc and modify the rhs using sol.
Definition: sparsemat.cpp:1708
void Transpose(const Table &A, Table &At, int ncols_A_)
Transpose a Table.
Definition: table.cpp:414
int NumCols() const
Get the number of columns (size of input) of the Operator. Synonym with Width().
Definition: operator.hpp:75
std::vector< OperatorPtr > agg_hdivdof
std::vector< OperatorPtr > agg_l2dof
GMRES method.
Definition: solvers.hpp:348
Parameters for the divergence free solver.
void Mult(const double *x, double *y) const
Matrix vector multiplication with the inverse of dense matrix.
Definition: densemat.cpp:3275
std::vector< OperatorPtr > C
void GetDiag(Vector &d) const
Returns the Diagonal of A.
Definition: sparsemat.cpp:549
void CopyMNt(const DenseMatrix &A, int row_offset, int col_offset)
Copy matrix A^t to the location in *this at row_offset, col_offset.
Definition: densemat.cpp:1558
void SetOptions(IterativeSolver &solver, const IterSolveParameters &param)
HypreParMatrix * TwoStepsRAP(const HypreParMatrix &Rt, const HypreParMatrix &A, const HypreParMatrix &P)
Block diagonal solver for symmetric A, each block is inverted by direct solver.
void AddDomainIntegrator(BilinearFormIntegrator *bfi)
Adds new Domain Integrator. Assumes ownership of bfi.
virtual void Mult(const Vector &x, Vector &y) const
Operator application: y=A(x).
void GetRowColumnsRef(const SparseMatrix &A, int row, Array< int > &cols)
Class for parallel bilinear form.
virtual void Mult(const Vector &x, Vector &y) const
Operator application: y=A(x).
static void ListToMarker(const Array< int > &list, int marker_size, Array< int > &marker, int mark_val=-1)
Convert an array of indices (list) to a Boolean marker array where all indices in the list are marked...
Definition: fespace.cpp:559
void MakeRef(T *, int)
Make this Array a reference to a pointer.
Definition: array.hpp:859
virtual void SetOperator(const Operator &op)
Also calls SetOperator for the preconditioner if there is one.
Definition: solvers.hpp:330
Arbitrary order H(curl)-conforming Nedelec finite elements.
Definition: fe_coll.hpp:416
Vector data type.
Definition: vector.hpp:60
ID for class HypreParMatrix.
Definition: operator.hpp:256
HYPRE_Int Mult(HypreParVector &x, HypreParVector &y, double alpha=1.0, double beta=0.0) const
Computes y = alpha * A * x + beta * y.
Definition: hypre.cpp:1529
virtual void SetPreconditioner(Solver &pr)
This should be called before SetOperator.
Definition: solvers.cpp:92
Parameters for iterative solver.
void DropSmallEntries(double tol)
Wrapper for hypre_ParCSRMatrixDropSmallEntries in different versions of hypre. Drop off-diagonal entr...
Definition: hypre.cpp:2006
int * GetI()
Definition: table.hpp:113
Arbitrary order H1-conforming (continuous) finite elements.
Definition: fe_coll.hpp:216
void CopyMN(const DenseMatrix &A, int m, int n, int Aro, int Aco)
Copy the m x n submatrix of A at row/col offsets Aro/Aco to *this.
Definition: densemat.cpp:1532
const Table & GetElementToDofTable() const
Return a reference to the internal Table that stores the lists of scalar dofs, for each mesh element...
Definition: fespace.hpp:700
HypreParMatrix * LeftDiagMult(const SparseMatrix &D, HYPRE_BigInt *row_starts=NULL) const
Multiply the HypreParMatrix on the left by a block-diagonal parallel matrix D and return the result a...
Definition: hypre.cpp:1705
Base class for solvers.
Definition: operator.hpp:648
void SetSize(int s)
Change the size of the DenseMatrix to s x s.
Definition: densemat.hpp:105
SparseMatrix * AggToInteriorDof(const Array< int > &bdr_truedofs, const SparseMatrix &agg_elem, const SparseMatrix &elem_dof, const HypreParMatrix &dof_truedof, Array< HYPRE_Int > &agg_starts)
Solver for local problems in SaddleSchwarzSmoother.
Wrapper for hypre&#39;s ParCSR matrix class.
Definition: hypre.hpp:277
A class to handle Block systems in a matrix-free implementation.
BBTSolver(const HypreParMatrix &B, IterSolveParameters param)
Class for parallel meshes.
Definition: pmesh.hpp:32
virtual Type GetType() const =0
Returns element&#39;s type.
std::vector< OperatorPtr > P_l2
static void AdjustVDofs(Array< int > &vdofs)
Definition: fespace.cpp:249
void Clear()
Clear the elapsed time on the stopwatch and restart it if it&#39;s running.
Definition: tic_toc.cpp:406
void SetBlock(int iRow, int iCol, Operator *op, double c=1.0)
Add a block op in the block-entry (iblock, jblock).
void Reset(OpType *A, bool own_A=true)
Reset the OperatorHandle to the given OpType pointer, A.
Definition: handle.hpp:140
void SetDiagonalBlock(int iblock, Operator *op)
Add a square block op in the block-entry (iblock, iblock).
Vector & GetBlock(int i)
Get the i-th vector in the block.
Definition: blockvector.hpp:90