MFEM  v4.4.0
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
constraints.cpp
Go to the documentation of this file.
1 // Copyright (c) 2010-2022, 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 "constraints.hpp"
13 
14 #include "../fem/fespace.hpp"
15 #include "../fem/pfespace.hpp"
16 
17 #include <set>
18 
19 namespace mfem
20 {
21 
22 Eliminator::Eliminator(const SparseMatrix& B, const Array<int>& lagrange_tdofs_,
23  const Array<int>& primary_tdofs_,
24  const Array<int>& secondary_tdofs_)
25  :
26  lagrange_tdofs(lagrange_tdofs_),
27  primary_tdofs(primary_tdofs_),
28  secondary_tdofs(secondary_tdofs_)
29 {
30  MFEM_VERIFY(lagrange_tdofs.Size() == secondary_tdofs.Size(),
31  "Dof sizes don't match!");
32 
33  Bp.SetSize(lagrange_tdofs.Size(), primary_tdofs.Size());
34  B.GetSubMatrix(lagrange_tdofs, primary_tdofs, Bp);
35 
36  Bs.SetSize(lagrange_tdofs.Size(), secondary_tdofs.Size());
37  B.GetSubMatrix(lagrange_tdofs, secondary_tdofs, Bs);
38  BsT.Transpose(Bs);
39 
40  ipiv.SetSize(Bs.Height());
41  Bsinverse.data = Bs.HostReadWrite();
42  Bsinverse.ipiv = ipiv.HostReadWrite();
43  Bsinverse.Factor(Bs.Height());
44 
45  ipivT.SetSize(Bs.Height());
46  BsTinverse.data = BsT.HostReadWrite();
47  BsTinverse.ipiv = ipivT.HostReadWrite();
48  BsTinverse.Factor(Bs.Height());
49 }
50 
51 void Eliminator::Eliminate(const Vector& vin, Vector& vout) const
52 {
53  Bp.Mult(vin, vout);
54  Bsinverse.Solve(Bs.Height(), 1, vout);
55  vout *= -1.0;
56 }
57 
58 void Eliminator::EliminateTranspose(const Vector& vin, Vector& vout) const
59 {
60  Vector work(vin);
61  BsTinverse.Solve(Bs.Height(), 1, work);
62  Bp.MultTranspose(work, vout);
63  vout *= -1.0;
64 }
65 
66 void Eliminator::LagrangeSecondary(const Vector& vin, Vector& vout) const
67 {
68  vout = vin;
69  Bsinverse.Solve(Bs.Height(), 1, vout);
70 }
71 
73  Vector& vout) const
74 {
75  vout = vin;
76  BsTinverse.Solve(Bs.Height(), 1, vout);
77 }
78 
80 {
81  mat.SetSize(Bp.Height(), Bp.Width());
82  mat = Bp;
83  Bsinverse.Solve(Bs.Height(), Bp.Width(), mat.GetData());
84  mat *= -1.0;
85 }
86 
88  Array<Eliminator*>& eliminators_)
89  :
90  Operator(A.Height()),
91  Aop(A),
92  eliminators(eliminators_)
93 {
94 }
95 
96 void EliminationProjection::Mult(const Vector& vin, Vector& vout) const
97 {
98  MFEM_ASSERT(vin.Size() == width, "Wrong vector size!");
99  MFEM_ASSERT(vout.Size() == height, "Wrong vector size!");
100 
101  vout = vin;
102 
103  for (int k = 0; k < eliminators.Size(); ++k)
104  {
105  Eliminator* elim = eliminators[k];
106  Vector subvec_in;
107  Vector subvec_out(elim->SecondaryDofs().Size());
108  vin.GetSubVector(elim->PrimaryDofs(), subvec_in);
109  elim->Eliminate(subvec_in, subvec_out);
110  vout.SetSubVector(elim->SecondaryDofs(), subvec_out);
111  }
112 }
113 
115 {
116  MFEM_ASSERT(vin.Size() == height, "Wrong vector size!");
117  MFEM_ASSERT(vout.Size() == width, "Wrong vector size!");
118 
119  vout = vin;
120 
121  for (int k = 0; k < eliminators.Size(); ++k)
122  {
123  Eliminator* elim = eliminators[k];
124  Vector subvec_in;
125  Vector subvec_out(elim->PrimaryDofs().Size());
126  vin.GetSubVector(elim->SecondaryDofs(), subvec_in);
127  elim->EliminateTranspose(subvec_in, subvec_out);
128  vout.AddElementVector(elim->PrimaryDofs(), subvec_out);
129  vout.SetSubVector(elim->SecondaryDofs(), 0.0);
130  }
131 }
132 
134 {
135  SparseMatrix * mat = new SparseMatrix(height, width);
136 
137  for (int i = 0; i < height; ++i)
138  {
139  mat->Add(i, i, 1.0);
140  }
141 
142  for (int k = 0; k < eliminators.Size(); ++k)
143  {
144  Eliminator* elim = eliminators[k];
145  DenseMatrix mat_k;
146  elim->ExplicitAssembly(mat_k);
147  for (int iz = 0; iz < elim->SecondaryDofs().Size(); ++iz)
148  {
149  int i = elim->SecondaryDofs()[iz];
150  for (int jz = 0; jz < elim->PrimaryDofs().Size(); ++jz)
151  {
152  int j = elim->PrimaryDofs()[jz];
153  mat->Add(i, j, mat_k(iz, jz));
154  }
155  mat->Set(i, i, 0.0);
156  }
157  }
158 
159  mat->Finalize();
160  return mat;
161 }
162 
163 void EliminationProjection::BuildGTilde(const Vector& r, Vector& rtilde) const
164 {
165  MFEM_ASSERT(rtilde.Size() == Aop.Height(), "Sizes don't match!");
166 
167  rtilde = 0.0;
168  for (int k = 0; k < eliminators.Size(); ++k)
169  {
170  Eliminator* elim = eliminators[k];
171  Vector subr;
172  r.GetSubVector(elim->LagrangeDofs(), subr);
173  Vector bsinvr(subr.Size());
174  elim->LagrangeSecondary(subr, bsinvr);
175  rtilde.AddElementVector(elim->SecondaryDofs(), bsinvr);
176  }
177 }
178 
180  const Vector& disprhs, const Vector& disp, Vector& lagrangem) const
181 {
182  lagrangem = 0.0;
183  MFEM_ASSERT(disp.Size() == Aop.Height(), "Sizes don't match!");
184 
185  Vector fullrhs(Aop.Height());
186  Aop.Mult(disp, fullrhs);
187  fullrhs -= disprhs;
188  fullrhs *= -1.0;
189  for (int k = 0; k < eliminators.Size(); ++k)
190  {
191  Eliminator* elim = eliminators[k];
192  Vector localsec;
193  fullrhs.GetSubVector(elim->SecondaryDofs(), localsec);
194  Vector locallagrange(localsec.Size());
195  elim->LagrangeSecondaryTranspose(localsec, locallagrange);
196  lagrangem.AddElementVector(elim->LagrangeDofs(), locallagrange);
197  }
198 }
199 
200 #ifdef MFEM_USE_MPI
201 
203 {
204  delete h_explicit_operator;
205  for (auto elim : eliminators)
206  {
207  delete elim;
208  }
209  delete projector;
210  delete prec;
211  delete krylov;
212 }
213 
215 {
216  SparseMatrix * explicit_projector = projector->AssembleExact();
217  HypreParMatrix * h_explicit_projector =
219  hA.GetRowStarts(), explicit_projector);
220  h_explicit_projector->CopyRowStarts();
221  h_explicit_projector->CopyColStarts();
222 
223  h_explicit_operator = RAP(&hA, h_explicit_projector);
224  // next line because of square projector
228 
229  delete explicit_projector;
230  delete h_explicit_projector;
231 }
232 
234  Array<int>& primary_dofs,
235  Array<int>& secondary_dofs)
236  :
237  ConstrainedSolver(A.GetComm(), A, B),
238  hA(A),
239  krylov(nullptr),
240  prec(nullptr)
241 {
242  MFEM_VERIFY(secondary_dofs.Size() == B.Height(),
243  "Wrong number of dofs for elimination!");
244  Array<int> lagrange_dofs(secondary_dofs.Size());
245  for (int i = 0; i < lagrange_dofs.Size(); ++i)
246  {
247  lagrange_dofs[i] = i;
248  }
249  eliminators.Append(new Eliminator(B, lagrange_dofs, primary_dofs,
250  secondary_dofs));
253 }
254 
256  Array<int>& constraint_rowstarts)
257  :
258  ConstrainedSolver(A.GetComm(), A, B),
259  hA(A),
260  krylov(nullptr),
261  prec(nullptr)
262 {
263  if (!B.Empty())
264  {
265  int * I = B.GetI();
266  int * J = B.GetJ();
267  double * data = B.GetData();
268 
269  for (int k = 0; k < constraint_rowstarts.Size() - 1; ++k)
270  {
271  int constraint_size = constraint_rowstarts[k + 1] -
272  constraint_rowstarts[k];
273  Array<int> lagrange_dofs(constraint_size);
274  Array<int> primary_dofs;
275  Array<int> secondary_dofs(constraint_size);
276  secondary_dofs = -1;
277  // loop through rows, identify one secondary dof for each row
278  for (int i = constraint_rowstarts[k]; i < constraint_rowstarts[k + 1]; ++i)
279  {
280  lagrange_dofs[i - constraint_rowstarts[k]] = i;
281  for (int jptr = I[i]; jptr < I[i + 1]; ++jptr)
282  {
283  int j = J[jptr];
284  double val = data[jptr];
285  if (std::abs(val) > 1.e-12 && secondary_dofs.Find(j) == -1)
286  {
287  secondary_dofs[i - constraint_rowstarts[k]] = j;
288  break;
289  }
290  }
291  }
292  // loop through rows again, assigning non-secondary dofs as primary
293  for (int i = constraint_rowstarts[k]; i < constraint_rowstarts[k + 1]; ++i)
294  {
295  MFEM_ASSERT(secondary_dofs[i - constraint_rowstarts[k]] >= 0,
296  "Secondary dofs don't match rows!");
297  for (int jptr = I[i]; jptr < I[i + 1]; ++jptr)
298  {
299  int j = J[jptr];
300  if (secondary_dofs.Find(j) == -1)
301  {
302  primary_dofs.Append(j);
303  }
304  }
305  }
306  primary_dofs.Sort();
307  primary_dofs.Unique();
308  eliminators.Append(new Eliminator(B, lagrange_dofs, primary_dofs,
309  secondary_dofs));
310  }
311  }
314 }
315 
316 void EliminationSolver::Mult(const Vector& rhs, Vector& sol) const
317 {
318  if (!prec)
319  {
321  }
322  if (!krylov)
323  {
324  krylov = BuildKrylov();
327  }
332 
333  Vector rtilde(rhs.Size());
334  if (constraint_rhs.Size() > 0)
335  {
337  }
338  else
339  {
340  rtilde = 0.0;
341  }
342  Vector temprhs(rhs);
343  hA.Mult(-1.0, rtilde, 1.0, temprhs);
344 
345  Vector reducedrhs(rhs.Size());
346  projector->MultTranspose(temprhs, reducedrhs);
347  Vector reducedsol(rhs.Size());
348  reducedsol = 0.0;
349  krylov->Mult(reducedrhs, reducedsol);
353 
354  projector->Mult(reducedsol, sol);
356  sol += rtilde;
357 }
358 
359 void PenaltyConstrainedSolver::Initialize(HypreParMatrix& A, HypreParMatrix& B)
360 {
361  HypreParMatrix * hBT = B.Transpose();
362  HypreParMatrix * hBTB = ParMult(hBT, &B, true);
363  // this matrix doesn't get cleanly deleted?
364  // (hypre comm pkg)
365  (*hBTB) *= penalty;
366  penalized_mat = ParAdd(&A, hBTB);
367  delete hBTB;
368  delete hBT;
369 }
370 
372  HypreParMatrix& A, SparseMatrix& B, double penalty_)
373  :
374  ConstrainedSolver(A.GetComm(), A, B),
375  penalty(penalty_),
376  constraintB(B),
377  krylov(nullptr),
378  prec(nullptr)
379 {
380  int rank, size;
381  MPI_Comm_rank(A.GetComm(), &rank);
382  MPI_Comm_size(A.GetComm(), &size);
383 
384  int constraint_running_total = 0;
385  int local_constraints = B.Height();
386  MPI_Scan(&local_constraints, &constraint_running_total, 1, MPI_INT,
387  MPI_SUM, A.GetComm());
388  int global_constraints = 0;
389  if (rank == size - 1) { global_constraints = constraint_running_total; }
390  MPI_Bcast(&global_constraints, 1, MPI_INT, size - 1, A.GetComm());
391 
392  HYPRE_BigInt glob_num_rows = global_constraints;
393  HYPRE_BigInt glob_num_cols = A.N();
394  HYPRE_BigInt row_starts[2] = { constraint_running_total - local_constraints,
395  constraint_running_total
396  };
397  HYPRE_BigInt col_starts[2] = { A.ColPart()[0], A.ColPart()[1] };
398  HypreParMatrix hB(A.GetComm(), glob_num_rows, glob_num_cols,
399  row_starts, col_starts, &B);
400  hB.CopyRowStarts();
401  hB.CopyColStarts();
402  Initialize(A, hB);
403 }
404 
406  HypreParMatrix& A, HypreParMatrix& B, double penalty_)
407  :
408  ConstrainedSolver(A.GetComm(), A, B),
409  penalty(penalty_),
410  constraintB(B),
411  krylov(nullptr),
412  prec(nullptr)
413 {
414  Initialize(A, B);
415 }
416 
418 {
419  delete penalized_mat;
420  delete prec;
421  delete krylov;
422 }
423 
425 {
426  if (!prec)
427  {
429  }
430  if (!krylov)
431  {
432  krylov = BuildKrylov();
435  }
436 
437  // form penalized right-hand side
438  Vector penalized_rhs(b);
439  if (constraint_rhs.Size() > 0)
440  {
441  Vector temp(x.Size());
443  temp *= penalty;
444  penalized_rhs += temp;
445  }
446 
447  // actually solve
452  krylov->Mult(penalized_rhs, x);
456 
458  if (constraint_rhs.Size() > 0)
459  {
461  }
463 }
464 
465 #endif
466 
467 /// because IdentityOperator isn't a Solver
468 class IdentitySolver : public Solver
469 {
470 public:
471  IdentitySolver(int size) : Solver(size) { }
472  void Mult(const Vector& x, Vector& y) const { y = x; }
473  void SetOperator(const Operator& op) { }
474 };
475 
476 void SchurConstrainedSolver::Initialize()
477 {
478  offsets[0] = 0;
479  offsets[1] = A.Height();
480  offsets[2] = A.Height() + B.Height();
481 
482  block_op = new BlockOperator(offsets);
483  block_op->SetBlock(0, 0, &A);
484  block_op->SetBlock(1, 0, &B);
485  tr_B = new TransposeOperator(&B);
486  block_op->SetBlock(0, 1, tr_B);
487 
488  block_pc = new BlockDiagonalPreconditioner(block_op->RowOffsets()),
489  rel_tol = 1.e-6;
490 }
491 
492 #ifdef MFEM_USE_MPI
494  Operator& A_, Operator& B_,
495  Solver& primal_pc_)
496  :
497  ConstrainedSolver(comm, A_, B_),
498  offsets(3),
499  primal_pc(&primal_pc_),
500  dual_pc(nullptr)
501 {
502  Initialize();
504  dual_pc = new IdentitySolver(block_op->RowOffsets()[2] -
505  block_op->RowOffsets()[1]);
508 }
509 #endif
510 
512  Solver& primal_pc_)
513  :
514  ConstrainedSolver(A_, B_),
515  offsets(3),
516  primal_pc(&primal_pc_),
517  dual_pc(nullptr)
518 {
519  Initialize();
521  dual_pc = new IdentitySolver(block_op->RowOffsets()[2] -
522  block_op->RowOffsets()[1]);
525 }
526 
527 #ifdef MFEM_USE_MPI
528 // protected constructor
530  Operator& B_)
531  :
532  ConstrainedSolver(comm, A_, B_),
533  offsets(3),
534  primal_pc(nullptr),
535  dual_pc(nullptr)
536 {
537  Initialize();
538 }
539 #endif
540 
541 // protected constructor
543  :
544  ConstrainedSolver(A_, B_),
545  offsets(3),
546  primal_pc(nullptr),
547  dual_pc(nullptr)
548 {
549  Initialize();
550 }
551 
553 {
554  delete block_op;
555  delete tr_B;
556  delete block_pc;
557  delete dual_pc;
558 }
559 
561  Vector& y) const
562 {
563  GMRESSolver * gmres;
564 #ifdef MFEM_USE_MPI
565  if (GetComm() != MPI_COMM_NULL)
566  {
567  gmres = new GMRESSolver(GetComm());
568  }
569  else
570 #endif
571  {
572  gmres = new GMRESSolver;
573  }
574  gmres->SetOperator(*block_op);
575  gmres->SetRelTol(rel_tol);
576  gmres->SetAbsTol(abs_tol);
577  gmres->SetMaxIter(max_iter);
579  gmres->SetPreconditioner(
580  const_cast<BlockDiagonalPreconditioner&>(*block_pc));
581 
582  gmres->Mult(x, y);
583  final_iter = gmres->GetNumIterations();
584  delete gmres;
585 }
586 
587 #ifdef MFEM_USE_MPI
589  HypreParMatrix& hA_,
590  HypreParMatrix& hB_,
591  int dimension,
592  bool reorder)
593  :
594  SchurConstrainedSolver(comm, hA_, hB_),
595  hA(hA_),
596  hB(hB_)
597 {
598  auto h_primal_pc = new HypreBoomerAMG(hA);
599  h_primal_pc->SetPrintLevel(0);
600  if (dimension > 0)
601  {
602  h_primal_pc->SetSystemsOptions(dimension, reorder);
603  }
604  primal_pc = h_primal_pc;
605 
606  HypreParMatrix * scaledB = new HypreParMatrix(hB);
607  Vector diagA;
608  hA.GetDiag(diagA);
609  HypreParMatrix * scaledBT = scaledB->Transpose();
610  scaledBT->InvScaleRows(diagA);
611  schur_mat = ParMult(scaledB, scaledBT);
612  schur_mat->CopyRowStarts();
613  schur_mat->CopyColStarts();
614  auto h_dual_pc = new HypreBoomerAMG(*schur_mat);
615  h_dual_pc->SetPrintLevel(0);
616  dual_pc = h_dual_pc;
617  delete scaledB;
618  delete scaledBT;
619 
622 }
623 
625 {
626  delete schur_mat;
627  delete primal_pc;
628 }
629 #endif
630 
631 void ConstrainedSolver::Initialize()
632 {
633  height = A.Height() + B.Height();
634  width = A.Width() + B.Height();
635 
636  workb.SetSize(A.Height());
637  workx.SetSize(A.Height());
639  constraint_rhs = 0.0;
641 }
642 
643 #ifdef MFEM_USE_MPI
645  :
646  IterativeSolver(comm), A(A_), B(B_)
647 {
648  Initialize();
649 }
650 #endif
651 
653  :
654  A(A_), B(B_)
655 {
656  Initialize();
657 }
658 
660 {
661  MFEM_VERIFY(r.Size() == multiplier_sol.Size(), "Vector is wrong size!");
662  constraint_rhs = r;
663 }
664 
665 void ConstrainedSolver::Mult(const Vector& f, Vector &x) const
666 {
667  Vector pworkb(A.Height() + B.Height());
668  Vector pworkx(A.Height() + B.Height());
669  pworkb = 0.0;
670  pworkx = 0.0;
671  for (int i = 0; i < f.Size(); ++i)
672  {
673  pworkb(i) = f(i);
674  pworkx(i) = x(i);
675  }
676  for (int i = 0; i < B.Height(); ++i)
677  {
678  pworkb(f.Size() + i) = constraint_rhs(i);
679  }
680 
681  LagrangeSystemMult(pworkb, pworkx);
682 
683  for (int i = 0; i < f.Size(); ++i)
684  {
685  x(i) = pworkx(i);
686  }
687  for (int i = 0; i < B.Height(); ++i)
688  {
689  multiplier_sol(i) = pworkx(f.Size() + i);
690  }
691 }
692 
694  Vector& x_and_lambda) const
695 {
696  workb.MakeRef(const_cast<Vector&>(f_and_r), 0);
697  workx.MakeRef(x_and_lambda, 0);
698  Vector ref_constraint_rhs(f_and_r.GetData() + A.Height(), B.Height());
699  constraint_rhs = ref_constraint_rhs;
700  Mult(workb, workx);
701  Vector ref_constraint_sol(x_and_lambda.GetData() + A.Height(), B.Height());
702  GetMultiplierSolution(ref_constraint_sol);
703 }
704 
705 /* Helper routine to reduce code duplication - given a node (which MFEM
706  sometimes calls a "dof"), this returns what normal people call a dof but
707  which MFEM sometimes calls a "vdof" - note that MFEM's naming conventions
708  regarding this are not entirely consistent. In parallel, this always
709  returns the "truedof" in parallel numbering. */
711  int node, bool parallel, int d=0)
712 {
713 #ifdef MFEM_USE_MPI
714  if (parallel)
715  {
716  ParFiniteElementSpace* pfespace =
717  dynamic_cast<ParFiniteElementSpace*>(&fespace);
718  if (pfespace)
719  {
720  const int vdof = pfespace->DofToVDof(node, d);
721  return pfespace->GetLocalTDofNumber(vdof);
722  }
723  else
724  {
725  MFEM_ABORT("Asked for parallel form of serial object!");
726  return -1;
727  }
728  }
729  else
730 #endif
731  {
732  return fespace.DofToVDof(node, d);
733  }
734 }
735 
737  Array<int>& constrained_att,
738  Array<int>& constraint_rowstarts,
739  bool parallel)
740 {
741  int dim = fespace.GetVDim();
742 
743  // dof_constraint maps a dof (column of the constraint matrix) to
744  // a block-constraint
745  // the indexing is by tdof, but a single tdof uniquely identifies a node
746  // so we only store one tdof independent of dimension
747  std::map<int, int> dof_bconstraint;
748  // constraints[j] is a map from attribute to row number,
749  // the j itself is the index of a block-constraint
750  std::vector<std::map<int, int> > constraints;
751  int n_bconstraints = 0;
752  int n_rows = 0;
753  for (int att : constrained_att)
754  {
755  // identify tdofs on constrained boundary
756  std::set<int> constrained_tdofs;
757  for (int i = 0; i < fespace.GetNBE(); ++i)
758  {
759  if (fespace.GetBdrAttribute(i) == att)
760  {
761  Array<int> nodes;
762  // get nodes on boundary (MFEM sometimes calls these dofs, what
763  // we call dofs it calls vdofs)
764  fespace.GetBdrElementDofs(i, nodes);
765  for (auto k : nodes)
766  {
767  // get the (local) dof number corresponding to
768  // the x-coordinate dof for node k
769  int tdof = CanonicalNodeNumber(fespace, k, parallel);
770  if (tdof >= 0) { constrained_tdofs.insert(tdof); }
771  }
772  }
773  }
774  // fill in the maps identifying which constraints (rows) correspond to
775  // which tdofs
776  for (auto k : constrained_tdofs)
777  {
778  auto it = dof_bconstraint.find(k);
779  if (it == dof_bconstraint.end())
780  {
781  // build new block constraint
782  dof_bconstraint[k] = n_bconstraints++;
783  constraints.emplace_back();
784  constraints.back()[att] = n_rows++;
785  }
786  else
787  {
788  // add tdof to existing block constraint
789  constraints[it->second][att] = n_rows++;
790  }
791  }
792  }
793 
794  // reorder so block-constraints eliminated together are grouped together in
795  // adjacent rows
796  {
797  std::map<int, int> reorder_rows;
798  int new_row = 0;
799  constraint_rowstarts.DeleteAll();
800  constraint_rowstarts.Append(0);
801  for (auto& it : dof_bconstraint)
802  {
803  int bconstraint_index = it.second;
804  bool nconstraint = false;
805  for (auto& att_it : constraints[bconstraint_index])
806  {
807  auto rrit = reorder_rows.find(att_it.second);
808  if (rrit == reorder_rows.end())
809  {
810  nconstraint = true;
811  reorder_rows[att_it.second] = new_row++;
812  }
813  }
814  if (nconstraint) { constraint_rowstarts.Append(new_row); }
815  }
816  MFEM_VERIFY(new_row == n_rows, "Remapping failed!");
817  for (auto& constraint_map : constraints)
818  {
819  for (auto& it : constraint_map)
820  {
821  it.second = reorder_rows[it.second];
822  }
823  }
824  }
825 
826  SparseMatrix * mout = new SparseMatrix(n_rows, fespace.GetTrueVSize());
827 
828  // fill in constraint matrix with normal vector information
829  Vector nor(dim);
830  // how many times we have seen a node (key is truek)
831  std::map<int, int> node_visits;
832  for (int i = 0; i < fespace.GetNBE(); ++i)
833  {
834  int att = fespace.GetBdrAttribute(i);
835  if (constrained_att.FindSorted(att) != -1)
836  {
838  const FiniteElement * fe = fespace.GetBE(i);
839  const IntegrationRule& nodes = fe->GetNodes();
840 
841  Array<int> dofs;
842  fespace.GetBdrElementDofs(i, dofs);
843  MFEM_VERIFY(dofs.Size() == nodes.Size(),
844  "Something wrong in finite element space!");
845 
846  for (int j = 0; j < dofs.Size(); ++j)
847  {
848  Tr->SetIntPoint(&nodes[j]);
849  // the normal returned in the next line is scaled by h, which is
850  // probably what we want in most applications
851  CalcOrtho(Tr->Jacobian(), nor);
852 
853  int k = dofs[j];
854  int truek = CanonicalNodeNumber(fespace, k, parallel);
855  if (truek >= 0)
856  {
857  auto nv_it = node_visits.find(truek);
858  if (nv_it == node_visits.end())
859  {
860  node_visits[truek] = 1;
861  }
862  else
863  {
864  node_visits[truek]++;
865  }
866  int visits = node_visits[truek];
867  int bconstraint = dof_bconstraint[truek];
868  int row = constraints[bconstraint][att];
869  for (int d = 0; d < dim; ++d)
870  {
871  int inner_truek = CanonicalNodeNumber(fespace, k,
872  parallel, d);
873  if (visits == 1)
874  {
875  mout->Add(row, inner_truek, nor[d]);
876  }
877  else
878  {
879  mout->SetColPtr(row);
880  const double pv = mout->SearchRow(inner_truek);
881  const double scaling = ((double) (visits - 1)) /
882  ((double) visits);
883  // incremental average, based on how many times
884  // this node has been visited
885  mout->Set(row, inner_truek,
886  scaling * pv + (1.0 / visits) * nor[d]);
887  }
888 
889  }
890  }
891  }
892  }
893  }
894  mout->Finalize();
895 
896  return mout;
897 }
898 
899 #ifdef MFEM_USE_MPI
901  Array<int>& constrained_att,
902  Array<int>& constraint_rowstarts)
903 {
904  return BuildNormalConstraints(fespace, constrained_att,
905  constraint_rowstarts, true);
906 }
907 #endif
908 
909 }
Abstract class for all finite elements.
Definition: fe_base.hpp:235
HypreParMatrix & hA
T * HostReadWrite()
Shortcut for mfem::ReadWrite(a.GetMemory(), a.Size(), false).
Definition: array.hpp:324
int Size() const
Return the logical size of the array.
Definition: array.hpp:138
void Eliminate(const Vector &in, Vector &out) const
Definition: constraints.cpp:51
MPI_Comm GetComm() const
MPI communicator.
Definition: hypre.hpp:528
Class for an integration rule - an Array of IntegrationPoint.
Definition: intrules.hpp:90
void Mult(const Vector &x, Vector &y) const override
Solve for given .
int DofToVDof(int dof, int vd, int ndofs=-1) const
Definition: fespace.cpp:251
void Unique()
Removes duplicities from a sorted array. This requires operator== to be defined for T...
Definition: array.hpp:256
int GetNumIterations() const
Definition: solvers.hpp:246
MPI_Comm GetComm() const
Return the associated MPI communicator, or MPI_COMM_NULL if no communicator is set.
Definition: solvers.hpp:264
HypreParMatrix * RAP(const HypreParMatrix *A, const HypreParMatrix *P)
Returns the matrix P^t * A * P.
Definition: hypre.cpp:2673
Array< int > & RowOffsets()
Return the row offsets for block starts.
SparseMatrix * ParBuildNormalConstraints(ParFiniteElementSpace &fespace, Array< int > &constrained_att, Array< int > &constraint_rowstarts)
Parallel wrapper for BuildNormalConstraints.
Solve constrained system by solving original mixed system; see ConstrainedSolver. ...
EliminationSolver(HypreParMatrix &A, SparseMatrix &B, Array< int > &primary_dofs, Array< int > &secondary_dofs)
Constructor, with explicit splitting into primary/secondary dofs.
void Mult(const Vector &x, Vector &y) const override
Operator application: y=A(x).
Definition: constraints.cpp:96
virtual void Finalize(int skip_zeros=1)
Finalize the matrix initialization, switching the storage format from LIL to CSR. ...
Definition: sparsemat.hpp:515
SparseMatrix * BuildNormalConstraints(FiniteElementSpace &fespace, Array< int > &constrained_att, Array< int > &constraint_rowstarts, bool parallel)
Build a matrix constraining normal components to zero.
void SetSize(int s)
Resize the vector to size s.
Definition: vector.hpp:521
double & SearchRow(const int col)
Perform a fast search for an entry in the &quot;current row&quot;. See SetColPtr().
Definition: sparsemat.hpp:815
virtual DofTransformation * GetBdrElementDofs(int bel, Array< int > &dofs) const
Returns indices of degrees of freedom for boundary element &#39;bel&#39;.
Definition: fespace.cpp:2814
bool Empty() const
Check if the SparseMatrix is empty.
Definition: sparsemat.hpp:206
SchurConstrainedHypreSolver(MPI_Comm comm, HypreParMatrix &hA_, HypreParMatrix &hB_, int dimension=0, bool reorder=false)
void Mult(const Table &A, const Table &B, Table &C)
C = A * B (as boolean matrices)
Definition: table.cpp:472
int Width() const
Get the width (size of input) of the Operator. Synonym with NumCols().
Definition: operator.hpp:72
int CanonicalNodeNumber(FiniteElementSpace &fespace, int node, bool parallel, int d=0)
int * GetJ()
Return the array J.
Definition: sparsemat.hpp:214
EliminationProjection(const Operator &A, Array< Eliminator * > &eliminators)
Definition: constraints.cpp:87
An abstract class to solve the constrained system subject to the constraint .
Definition: constraints.hpp:53
void SetIntPoint(const IntegrationPoint *ip)
Set the integration point ip that weights and Jacobians will be evaluated at.
Definition: eltrans.hpp:93
void GetSubVector(const Array< int > &dofs, Vector &elemvect) const
Extract entries listed in dofs to the output Vector elemvect.
Definition: vector.cpp:534
void Initialize(HypreParMatrix &A, HypreParMatrix &B)
Data type dense matrix using column-major storage.
Definition: densemat.hpp:23
int Size() const
Returns the size of the vector.
Definition: vector.hpp:199
int GetBdrAttribute(int i) const
Definition: fespace.hpp:637
int * GetI()
Return the array I.
Definition: sparsemat.hpp:209
virtual void Mult(const Vector &x, Vector &y) const =0
Operator application: y=A(x).
Abstract parallel finite element space.
Definition: pfespace.hpp:28
double * HostReadWrite()
Shortcut for mfem::ReadWrite(GetMemory(), TotalSize(), false).
Definition: densemat.hpp:405
HypreParMatrix * h_explicit_operator
Eliminator(const SparseMatrix &B, const Array< int > &lagrange_dofs, const Array< int > &primary_tdofs, const Array< int > &secondary_tdofs)
Definition: constraints.cpp:22
bool GetConverged() const
Definition: solvers.hpp:247
int GetLocalTDofNumber(int ldof) const
Definition: pfespace.cpp:1062
void MultTranspose(const Vector &x, Vector &y) const override
Action of the transpose operator: y=A^t(x). The default behavior in class Operator is to generate an ...
double * GetData() const
Returns the matrix data array.
Definition: densemat.hpp:115
double * GetData() const
Return a pointer to the beginning of the Vector data.
Definition: vector.hpp:208
void CalcOrtho(const DenseMatrix &J, Vector &n)
Definition: densemat.cpp:2288
HypreParMatrix * ParMult(const HypreParMatrix *A, const HypreParMatrix *B, bool own_matrix)
Definition: hypre.cpp:2652
virtual Solver * BuildPreconditioner() const =0
Build preconditioner for penalized system.
void GetSubMatrix(const Array< int > &rows, const Array< int > &cols, DenseMatrix &subm) const
Definition: sparsemat.cpp:2748
void DeleteAll()
Delete the whole array.
Definition: array.hpp:846
Array< Eliminator * > eliminators
Operator & GetBlock(int i, int j)
Return a reference to block i,j.
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
ConstrainedSolver(MPI_Comm comm, Operator &A_, Operator &B_)
double GetFinalNorm() const
Definition: solvers.hpp:248
double * GetData()
Return the element data, i.e. the array A.
Definition: sparsemat.hpp:219
virtual void Mult(const Vector &f, Vector &x) const override
Solve for given .
HYPRE_BigInt N() const
Returns the global number of columns.
Definition: hypre.hpp:579
void EliminateZeroRows()
If a row contains only zeros, set its diagonal to 1.
Definition: hypre.hpp:752
virtual void SetOperator(const Operator &op) override
Also calls SetOperator for the preconditioner if there is one.
Definition: solvers.cpp:179
void LagrangeSecondary(const Vector &in, Vector &out) const
Maps Lagrange multipliers to secondary dofs, applies .
Definition: constraints.cpp:66
HYPRE_BigInt * GetRowStarts() const
Return the parallel row partitioning array.
Definition: hypre.hpp:639
void Set(const int i, const int j, const double val)
Definition: sparsemat.cpp:2630
virtual void Mult(const Vector &b, Vector &x) const
Operator application: y=A(x).
Definition: solvers.cpp:971
The BoomerAMG solver in hypre.
Definition: hypre.hpp:1443
void MultTranspose(const double *x, double *y) const
Multiply a vector with the transpose matrix.
Definition: densemat.cpp:188
void BuildGTilde(const Vector &g, Vector &gtilde) const
virtual void SetPrintLevel(int print_lvl)
Legacy method to set the level of verbosity of the solver output.
Definition: solvers.cpp:71
double f(const Vector &xvec)
Definition: lor_mms.hpp:32
int GetNBE() const
Returns number of boundary elements in the mesh.
Definition: fespace.hpp:599
Data type sparse matrix.
Definition: sparsemat.hpp:46
int Height() const
Get the height (size of output) of the Operator. Synonym with NumRows().
Definition: operator.hpp:66
int Append(const T &el)
Append element &#39;el&#39; to array, resize if necessary.
Definition: array.hpp:751
Perform elimination of a single constraint.
virtual IterativeSolver * BuildKrylov() const =0
Select krylov solver for penalized system.
double b
Definition: lissajous.cpp:42
const DenseMatrix & Jacobian()
Return the Jacobian matrix of the transformation at the currently set IntegrationPoint, using the method SetIntPoint().
Definition: eltrans.hpp:119
const Array< int > & PrimaryDofs() const
void SetMaxIter(int max_it)
Definition: solvers.hpp:200
void RecoverMultiplier(const Vector &primalrhs, const Vector &primalvars, Vector &lm) const
virtual void SetConstraintRHS(const Vector &r)
Set the right-hand side r for the constraint B x = r.
const IntegrationRule & GetNodes() const
Get a const reference to the nodes of the element.
Definition: fe_base.hpp:390
ElementTransformation * GetBdrElementTransformation(int i) const
Returns ElementTransformation for the i-th boundary element.
Definition: fespace.hpp:632
HypreParMatrix * Transpose() const
Returns the transpose of *this.
Definition: hypre.cpp:1607
virtual int GetTrueVSize() const
Return the number of vector true (conforming) dofs.
Definition: fespace.hpp:566
void ExplicitAssembly(DenseMatrix &mat) const
Return explicitly assembled in mat.
Definition: constraints.cpp:79
void GetDiag(Vector &diag) const
Get the local diagonal of the matrix.
Definition: hypre.cpp:1480
void Sort()
Sorts the array in ascending order. This requires operator&lt; to be defined for T.
Definition: array.hpp:248
void Solve(int m, int n, double *X) const
Definition: densemat.cpp:2999
int GetVDim() const
Returns vector dimension.
Definition: fespace.hpp:557
void AddElementVector(const Array< int > &dofs, const Vector &elemvect)
Add elements of the elemvect Vector to the entries listed in dofs. Negative dof values cause the -dof...
Definition: vector.cpp:626
void EliminateTranspose(const Vector &in, Vector &out) const
Transpose of Eliminate(), applies .
Definition: constraints.cpp:58
void BuildExplicitOperator()
Internal utility routine; assembles eliminated matrix explicitly.
void Add(const int i, const int j, const double val)
Definition: sparsemat.cpp:2649
void SetAbsTol(double atol)
Definition: solvers.hpp:199
PrintLevel print_options
Output behavior for the iterative solver.
Definition: solvers.hpp:137
void SetRelTol(double rtol)
Definition: solvers.hpp:198
PenaltyConstrainedSolver(HypreParMatrix &A, SparseMatrix &B, double penalty_)
Abstract base class for iterative solver.
Definition: solvers.hpp:66
void Transpose()
(*this) = (*this)^t
Definition: densemat.cpp:1359
void GetMultiplierSolution(Vector &lambda) const
Return the Lagrange multiplier solution in lambda.
Definition: constraints.hpp:76
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Definition: fespace.hpp:88
void SetColPtr(const int row) const
Initialize the SparseMatrix for fast access to the entries of the given row which becomes the &quot;curren...
Definition: sparsemat.hpp:761
int Find(const T &el) const
Return the first index where &#39;el&#39; is found; return -1 if not found.
Definition: array.hpp:810
void SetSize(int nsize)
Change the logical size of the array, keep existing entries.
Definition: array.hpp:679
int max_iter
Limit for the number of iterations the solver is allowed to do.
Definition: solvers.hpp:151
HYPRE_Int HYPRE_BigInt
virtual void LagrangeSystemMult(const Vector &f_and_r, Vector &x_and_lambda) const
Solve for (x, lambda) given (f, r)
GMRES method.
Definition: solvers.hpp:497
SchurConstrainedSolver(MPI_Comm comm, Operator &A_, Operator &B_, Solver &primal_pc_)
HypreParMatrix * ParAdd(const HypreParMatrix *A, const HypreParMatrix *B)
Returns the matrix A + B.
Definition: hypre.cpp:2613
virtual Solver * BuildPreconditioner() const =0
Build preconditioner for eliminated system.
int height
Dimension of the output / number of rows in the matrix.
Definition: operator.hpp:27
EliminationProjection * projector
double rel_tol
Relative tolerance.
Definition: solvers.hpp:154
void Mult(const Vector &x, Vector &y) const override
Solve for given .
HYPRE_BigInt GetGlobalNumRows() const
Return the global number of rows.
Definition: hypre.hpp:629
HYPRE_BigInt * ColPart()
Returns the column partitioning.
Definition: hypre.hpp:567
void InvScaleRows(const Vector &s)
Scale the local row i by 1./s(i)
Definition: hypre.cpp:2020
IterativeSolver * krylov
int dim
Definition: ex24.cpp:53
bool Factor(int m, double TOL=0.0)
Compute the LU factorization of the current matrix.
Definition: densemat.cpp:2848
virtual void LagrangeSystemMult(const Vector &x, Vector &y) const override
Solve for (x, lambda) given (f, r)
virtual void SetOperator(const Operator &op)=0
Set/update the solver for the given operator.
virtual IterativeSolver * BuildKrylov() const =0
Select krylov solver for eliminated system.
Vector data type.
Definition: vector.hpp:60
void Mult(const double *x, double *y) const
Matrix vector multiplication.
Definition: densemat.cpp:160
const Array< int > & LagrangeDofs() const
void MakeRef(Vector &base, int offset, int size)
Reset the Vector to be a reference to a sub-vector of base.
Definition: vector.hpp:585
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:1703
virtual void SetPreconditioner(Solver &pr)
This should be called before SetOperator.
Definition: solvers.cpp:173
BlockDiagonalPreconditioner * block_pc
Base class for solvers.
Definition: operator.hpp:651
const FiniteElement * GetBE(int i) const
Returns pointer to the FiniteElement in the FiniteElementCollection associated with i&#39;th boundary fac...
Definition: fespace.cpp:3102
SparseMatrix * AssembleExact() const
Assemble this projector as a (processor-local) SparseMatrix.
void SetSize(int s)
Change the size of the DenseMatrix to s x s.
Definition: densemat.hpp:105
Abstract operator.
Definition: operator.hpp:24
Wrapper for hypre&#39;s ParCSR matrix class.
Definition: hypre.hpp:337
TransposeOperator * tr_B
void LagrangeSecondaryTranspose(const Vector &in, Vector &out) const
Transpose of LagrangeSecondary()
Definition: constraints.cpp:72
int width
Dimension of the input / number of columns in the matrix.
Definition: operator.hpp:28
void SetBlock(int iRow, int iCol, Operator *op, double c=1.0)
Add a block op in the block-entry (iblock, jblock).
const Array< int > & SecondaryDofs() const
double abs_tol
Absolute tolerance.
Definition: solvers.hpp:157
void SetDiagonalBlock(int iblock, Operator *op)
Add a square block op in the block-entry (iblock, iblock).
double * data
Definition: densemat.hpp:532