MFEM  v4.5.1
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  else
323  {
325  }
326  if (!krylov)
327  {
328  krylov = BuildKrylov();
331  }
336 
337  Vector rtilde(rhs.Size());
338  if (constraint_rhs.Size() > 0)
339  {
341  }
342  else
343  {
344  rtilde = 0.0;
345  }
346  Vector temprhs(rhs);
347  hA.Mult(-1.0, rtilde, 1.0, temprhs);
348 
349  Vector reducedrhs(rhs.Size());
350  projector->MultTranspose(temprhs, reducedrhs);
351  Vector reducedsol(rhs.Size());
352  reducedsol = 0.0;
353  krylov->Mult(reducedrhs, reducedsol);
357 
358  projector->Mult(reducedsol, sol);
360  sol += rtilde;
361 }
362 
363 void PenaltyConstrainedSolver::Initialize(HypreParMatrix& A, HypreParMatrix& B,
364  HypreParMatrix& D)
365 {
366  HypreParMatrix * hBTB = RAP(&D, &B);
367  // this matrix doesn't get cleanly deleted?
368  // (hypre comm pkg)
369  penalized_mat = ParAdd(&A, hBTB);
370  delete hBTB;
371 }
372 
374  HypreParMatrix& A, SparseMatrix& B, double penalty_)
375  :
376  ConstrainedSolver(A.GetComm(), A, B),
377  penalty(B.Height()),
378  constraintB(B),
379  krylov(nullptr),
380  prec(nullptr)
381 {
382  int rank, size;
383  MPI_Comm_rank(A.GetComm(), &rank);
384  MPI_Comm_size(A.GetComm(), &size);
385 
386  int constraint_running_total = 0;
387  int local_constraints = B.Height();
388  MPI_Scan(&local_constraints, &constraint_running_total, 1, MPI_INT,
389  MPI_SUM, A.GetComm());
390  int global_constraints = 0;
391  if (rank == size - 1) { global_constraints = constraint_running_total; }
392  MPI_Bcast(&global_constraints, 1, MPI_INT, size - 1, A.GetComm());
393 
394  HYPRE_BigInt glob_num_rows = global_constraints;
395  HYPRE_BigInt glob_num_cols = A.N();
396  HYPRE_BigInt row_starts[2] = { constraint_running_total - local_constraints,
397  constraint_running_total
398  };
399  HYPRE_BigInt col_starts[2] = { A.ColPart()[0], A.ColPart()[1] };
400  HypreParMatrix hB(A.GetComm(), glob_num_rows, glob_num_cols,
401  row_starts, col_starts, &B);
402  hB.CopyRowStarts();
403  hB.CopyColStarts();
404  penalty=penalty_;
406  HypreParMatrix hD(hB.GetComm(), hB.M(), hB.RowPart(), &D);
407  hD.CopyRowStarts();
408  hD.CopyColStarts();
409  Initialize(A, hB, hD);
410 }
411 
413  HypreParMatrix& A, HypreParMatrix& B, double penalty_)
414  :
415  ConstrainedSolver(A.GetComm(), A, B),
416  penalty(B.Height()),
417  constraintB(B),
418  krylov(nullptr),
419  prec(nullptr)
420 {
421  penalty=penalty_;
423  HypreParMatrix hD(B.GetComm(), B.M(), B.RowPart(), &D);
424  hD.CopyRowStarts();
425  hD.CopyColStarts();
426  Initialize(A, B, hD);
427 }
428 
430  HypreParMatrix& A, HypreParMatrix& B, Vector& penalty_)
431  :
432  ConstrainedSolver(A.GetComm(), A, B),
433  penalty(penalty_),
434  constraintB(B),
435  krylov(nullptr),
436  prec(nullptr)
437 {
438  SparseMatrix D(penalty_);
439  HypreParMatrix hD(B.GetComm(), B.M(), B.RowPart(), &D);
440  hD.CopyRowStarts();
441  hD.CopyColStarts();
442  Initialize(A, B, hD);
443 }
444 
446 {
447  delete penalized_mat;
448  delete prec;
449  delete krylov;
450 }
451 
453 {
454  if (!prec)
455  {
457  }
458  else
459  {
461  }
462  if (!krylov)
463  {
464  krylov = BuildKrylov();
467  }
468 
469  // form penalized right-hand side
470  Vector penalized_rhs(b);
471  if (constraint_rhs.Size() > 0)
472  {
473  Vector temp_rhs(constraint_rhs.Size());
475  D.Mult(constraint_rhs, temp_rhs);
476  Vector temp(x.Size());
477  constraintB.MultTranspose(temp_rhs, temp);
478  penalized_rhs += temp;
479  }
480 
481  // actually solve
486  krylov->Mult(penalized_rhs, x);
490 
492  if (constraint_rhs.Size() > 0)
493  {
495  }
497 }
498 
499 #endif
500 
501 /// because IdentityOperator isn't a Solver
502 class IdentitySolver : public Solver
503 {
504 public:
505  IdentitySolver(int size) : Solver(size) { }
506  void Mult(const Vector& x, Vector& y) const { y = x; }
507  void SetOperator(const Operator& op) { }
508 };
509 
510 void SchurConstrainedSolver::Initialize()
511 {
512  offsets[0] = 0;
513  offsets[1] = A.Height();
514  offsets[2] = A.Height() + B.Height();
515 
516  block_op = new BlockOperator(offsets);
517  block_op->SetBlock(0, 0, &A);
518  block_op->SetBlock(1, 0, &B);
519  tr_B = new TransposeOperator(&B);
520  block_op->SetBlock(0, 1, tr_B);
521 
522  block_pc = new BlockDiagonalPreconditioner(block_op->RowOffsets()),
523  rel_tol = 1.e-6;
524 }
525 
526 #ifdef MFEM_USE_MPI
528  Operator& A_, Operator& B_,
529  Solver& primal_pc_)
530  :
531  ConstrainedSolver(comm, A_, B_),
532  offsets(3),
533  primal_pc(&primal_pc_),
534  dual_pc(nullptr)
535 {
536  Initialize();
538  dual_pc = new IdentitySolver(block_op->RowOffsets()[2] -
539  block_op->RowOffsets()[1]);
542 }
543 #endif
544 
546  Solver& primal_pc_)
547  :
548  ConstrainedSolver(A_, B_),
549  offsets(3),
550  primal_pc(&primal_pc_),
551  dual_pc(nullptr)
552 {
553  Initialize();
555  dual_pc = new IdentitySolver(block_op->RowOffsets()[2] -
556  block_op->RowOffsets()[1]);
559 }
560 
561 #ifdef MFEM_USE_MPI
562 // protected constructor
564  Operator& B_)
565  :
566  ConstrainedSolver(comm, A_, B_),
567  offsets(3),
568  primal_pc(nullptr),
569  dual_pc(nullptr)
570 {
571  Initialize();
572 }
573 #endif
574 
575 // protected constructor
577  :
578  ConstrainedSolver(A_, B_),
579  offsets(3),
580  primal_pc(nullptr),
581  dual_pc(nullptr)
582 {
583  Initialize();
584 }
585 
587 {
588  delete block_op;
589  delete tr_B;
590  delete block_pc;
591  delete dual_pc;
592 }
593 
595  Vector& y) const
596 {
597  GMRESSolver * gmres;
598 #ifdef MFEM_USE_MPI
599  if (GetComm() != MPI_COMM_NULL)
600  {
601  gmres = new GMRESSolver(GetComm());
602  }
603  else
604 #endif
605  {
606  gmres = new GMRESSolver;
607  }
608  gmres->SetOperator(*block_op);
609  gmres->SetRelTol(rel_tol);
610  gmres->SetAbsTol(abs_tol);
611  gmres->SetMaxIter(max_iter);
613  gmres->SetPreconditioner(
614  const_cast<BlockDiagonalPreconditioner&>(*block_pc));
615 
616  gmres->Mult(x, y);
617  final_iter = gmres->GetNumIterations();
618  converged = gmres->GetConverged();
619  final_norm = gmres->GetFinalNorm();
620  delete gmres;
621 }
622 
623 #ifdef MFEM_USE_MPI
625  HypreParMatrix& hA_,
626  HypreParMatrix& hB_,
627  Solver * prec,
628  int dimension,
629  bool reorder)
630  :
631  SchurConstrainedSolver(comm, hA_, hB_),
632  hA(hA_),
633  hB(hB_)
634 {
635  if (prec == nullptr)
636  {
637  auto h_primal_pc = new HypreBoomerAMG(hA);
638  h_primal_pc->SetPrintLevel(0);
639  if (dimension > 0)
640  {
641  h_primal_pc->SetSystemsOptions(dimension, reorder);
642  }
643  primal_pc = h_primal_pc;
644  }
645  else
646  {
647  primal_pc = prec;
648  }
649 
650  HypreParMatrix * scaledB = new HypreParMatrix(hB);
651  Vector diagA;
652  hA.GetDiag(diagA);
653  HypreParMatrix * scaledBT = scaledB->Transpose();
654  scaledBT->InvScaleRows(diagA);
655  schur_mat = ParMult(scaledB, scaledBT);
656  schur_mat->CopyRowStarts();
657  schur_mat->CopyColStarts();
658  auto h_dual_pc = new HypreBoomerAMG(*schur_mat);
659  h_dual_pc->SetPrintLevel(0);
660  dual_pc = h_dual_pc;
661  delete scaledB;
662  delete scaledBT;
663 
666 }
667 
669 {
670  delete schur_mat;
671  delete primal_pc;
672 }
673 #endif
674 
675 void ConstrainedSolver::Initialize()
676 {
677  height = A.Height() + B.Height();
678  width = A.Width() + B.Height();
679 
680  workb.SetSize(A.Height());
681  workx.SetSize(A.Height());
683  constraint_rhs = 0.0;
685 }
686 
687 #ifdef MFEM_USE_MPI
689  :
690  IterativeSolver(comm), A(A_), B(B_)
691 {
692  Initialize();
693 }
694 #endif
695 
697  :
698  A(A_), B(B_)
699 {
700  Initialize();
701 }
702 
704 {
705  MFEM_VERIFY(r.Size() == multiplier_sol.Size(), "Vector is wrong size!");
706  constraint_rhs = r;
707 }
708 
709 void ConstrainedSolver::Mult(const Vector& f, Vector &x) const
710 {
711  Vector pworkb(A.Height() + B.Height());
712  Vector pworkx(A.Height() + B.Height());
713  pworkb = 0.0;
714  pworkx = 0.0;
715  for (int i = 0; i < f.Size(); ++i)
716  {
717  pworkb(i) = f(i);
718  pworkx(i) = x(i);
719  }
720  for (int i = 0; i < B.Height(); ++i)
721  {
722  pworkb(f.Size() + i) = constraint_rhs(i);
723  }
724 
725  LagrangeSystemMult(pworkb, pworkx);
726 
727  for (int i = 0; i < f.Size(); ++i)
728  {
729  x(i) = pworkx(i);
730  }
731  for (int i = 0; i < B.Height(); ++i)
732  {
733  multiplier_sol(i) = pworkx(f.Size() + i);
734  }
735 }
736 
738  Vector& x_and_lambda) const
739 {
740  workb.MakeRef(const_cast<Vector&>(f_and_r), 0);
741  workx.MakeRef(x_and_lambda, 0);
742  Vector ref_constraint_rhs(f_and_r.GetData() + A.Height(), B.Height());
743  constraint_rhs = ref_constraint_rhs;
744  Mult(workb, workx);
745  Vector ref_constraint_sol(x_and_lambda.GetData() + A.Height(), B.Height());
746  GetMultiplierSolution(ref_constraint_sol);
747 }
748 
749 /* Helper routine to reduce code duplication - given a node (which MFEM
750  sometimes calls a "dof"), this returns what normal people call a dof but
751  which MFEM sometimes calls a "vdof" - note that MFEM's naming conventions
752  regarding this are not entirely consistent. In parallel, this always
753  returns the "truedof" in parallel numbering. */
755  int node, bool parallel, int d=0)
756 {
757 #ifdef MFEM_USE_MPI
758  if (parallel)
759  {
760  ParFiniteElementSpace* pfespace =
761  dynamic_cast<ParFiniteElementSpace*>(&fespace);
762  if (pfespace)
763  {
764  const int vdof = pfespace->DofToVDof(node, d);
765  return pfespace->GetLocalTDofNumber(vdof);
766  }
767  else
768  {
769  MFEM_ABORT("Asked for parallel form of serial object!");
770  return -1;
771  }
772  }
773  else
774 #endif
775  {
776  return fespace.DofToVDof(node, d);
777  }
778 }
779 
781  Array<int>& constrained_att,
782  Array<int>& constraint_rowstarts,
783  bool parallel)
784 {
785  int dim = fespace.GetVDim();
786 
787  // dof_constraint maps a dof (column of the constraint matrix) to
788  // a block-constraint
789  // the indexing is by tdof, but a single tdof uniquely identifies a node
790  // so we only store one tdof independent of dimension
791  std::map<int, int> dof_bconstraint;
792  // constraints[j] is a map from attribute to row number,
793  // the j itself is the index of a block-constraint
794  std::vector<std::map<int, int> > constraints;
795  int n_bconstraints = 0;
796  int n_rows = 0;
797  for (int att : constrained_att)
798  {
799  // identify tdofs on constrained boundary
800  std::set<int> constrained_tdofs;
801  for (int i = 0; i < fespace.GetNBE(); ++i)
802  {
803  if (fespace.GetBdrAttribute(i) == att)
804  {
805  Array<int> nodes;
806  // get nodes on boundary (MFEM sometimes calls these dofs, what
807  // we call dofs it calls vdofs)
808  fespace.GetBdrElementDofs(i, nodes);
809  for (auto k : nodes)
810  {
811  // get the (local) dof number corresponding to
812  // the x-coordinate dof for node k
813  int tdof = CanonicalNodeNumber(fespace, k, parallel);
814  if (tdof >= 0) { constrained_tdofs.insert(tdof); }
815  }
816  }
817  }
818  // fill in the maps identifying which constraints (rows) correspond to
819  // which tdofs
820  for (auto k : constrained_tdofs)
821  {
822  auto it = dof_bconstraint.find(k);
823  if (it == dof_bconstraint.end())
824  {
825  // build new block constraint
826  dof_bconstraint[k] = n_bconstraints++;
827  constraints.emplace_back();
828  constraints.back()[att] = n_rows++;
829  }
830  else
831  {
832  // add tdof to existing block constraint
833  constraints[it->second][att] = n_rows++;
834  }
835  }
836  }
837 
838  // reorder so block-constraints eliminated together are grouped together in
839  // adjacent rows
840  {
841  std::map<int, int> reorder_rows;
842  int new_row = 0;
843  constraint_rowstarts.DeleteAll();
844  constraint_rowstarts.Append(0);
845  for (auto& it : dof_bconstraint)
846  {
847  int bconstraint_index = it.second;
848  bool nconstraint = false;
849  for (auto& att_it : constraints[bconstraint_index])
850  {
851  auto rrit = reorder_rows.find(att_it.second);
852  if (rrit == reorder_rows.end())
853  {
854  nconstraint = true;
855  reorder_rows[att_it.second] = new_row++;
856  }
857  }
858  if (nconstraint) { constraint_rowstarts.Append(new_row); }
859  }
860  MFEM_VERIFY(new_row == n_rows, "Remapping failed!");
861  for (auto& constraint_map : constraints)
862  {
863  for (auto& it : constraint_map)
864  {
865  it.second = reorder_rows[it.second];
866  }
867  }
868  }
869 
870  SparseMatrix * mout = new SparseMatrix(n_rows, fespace.GetTrueVSize());
871 
872  // fill in constraint matrix with normal vector information
873  Vector nor(dim);
874  // how many times we have seen a node (key is truek)
875  std::map<int, int> node_visits;
876  for (int i = 0; i < fespace.GetNBE(); ++i)
877  {
878  int att = fespace.GetBdrAttribute(i);
879  if (constrained_att.FindSorted(att) != -1)
880  {
882  const FiniteElement * fe = fespace.GetBE(i);
883  const IntegrationRule& nodes = fe->GetNodes();
884 
885  Array<int> dofs;
886  fespace.GetBdrElementDofs(i, dofs);
887  MFEM_VERIFY(dofs.Size() == nodes.Size(),
888  "Something wrong in finite element space!");
889 
890  for (int j = 0; j < dofs.Size(); ++j)
891  {
892  Tr->SetIntPoint(&nodes[j]);
893  // the normal returned in the next line is scaled by h, which is
894  // probably what we want in most applications
895  CalcOrtho(Tr->Jacobian(), nor);
896 
897  int k = dofs[j];
898  int truek = CanonicalNodeNumber(fespace, k, parallel);
899  if (truek >= 0)
900  {
901  auto nv_it = node_visits.find(truek);
902  if (nv_it == node_visits.end())
903  {
904  node_visits[truek] = 1;
905  }
906  else
907  {
908  node_visits[truek]++;
909  }
910  int visits = node_visits[truek];
911  int bconstraint = dof_bconstraint[truek];
912  int row = constraints[bconstraint][att];
913  for (int d = 0; d < dim; ++d)
914  {
915  int inner_truek = CanonicalNodeNumber(fespace, k,
916  parallel, d);
917  if (visits == 1)
918  {
919  mout->Add(row, inner_truek, nor[d]);
920  }
921  else
922  {
923  mout->SetColPtr(row);
924  const double pv = mout->SearchRow(inner_truek);
925  const double scaling = ((double) (visits - 1)) /
926  ((double) visits);
927  // incremental average, based on how many times
928  // this node has been visited
929  mout->Set(row, inner_truek,
930  scaling * pv + (1.0 / visits) * nor[d]);
931  }
932 
933  }
934  }
935  }
936  }
937  }
938  mout->Finalize();
939 
940  return mout;
941 }
942 
943 #ifdef MFEM_USE_MPI
945  Array<int>& constrained_att,
946  Array<int>& constraint_rowstarts)
947 {
948  return BuildNormalConstraints(fespace, constrained_att,
949  constraint_rowstarts, true);
950 }
951 #endif
952 
953 }
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:534
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
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
void Initialize(HypreParMatrix &A, HypreParMatrix &B, HypreParMatrix &D)
virtual void Finalize(int skip_zeros=1)
Finalize the matrix initialization, switching the storage format from LIL to CSR. ...
Definition: sparsemat.hpp:545
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:513
double & SearchRow(const int col)
Perform a fast search for an entry in the &quot;current row&quot;. See SetColPtr().
Definition: sparsemat.hpp:845
virtual DofTransformation * GetBdrElementDofs(int bel, Array< int > &dofs) const
Returns indices of degrees of freedom for boundary element &#39;bel&#39;.
Definition: fespace.cpp:2812
bool Empty() const
Check if the SparseMatrix is empty.
Definition: sparsemat.hpp:219
void Mult(const Table &A, const Table &B, Table &C)
C = A * B (as boolean matrices)
Definition: table.cpp:475
SchurConstrainedHypreSolver(MPI_Comm comm, HypreParMatrix &hA_, HypreParMatrix &hB_, Solver *prec=nullptr, int dimension=0, bool reorder=false)
int Width() const
Get the width (size of input) of the Operator. Synonym with NumCols().
Definition: operator.hpp:73
int CanonicalNodeNumber(FiniteElementSpace &fespace, int node, bool parallel, int d=0)
int * GetJ()
Return the array J.
Definition: sparsemat.hpp:227
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:547
Data type dense matrix using column-major storage.
Definition: densemat.hpp:23
int Size() const
Returns the size of the vector.
Definition: vector.hpp:200
int GetBdrAttribute(int i) const
Definition: fespace.hpp:661
int * GetI()
Return the array I.
Definition: sparsemat.hpp:222
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:463
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:209
void CalcOrtho(const DenseMatrix &J, Vector &n)
Definition: densemat.cpp:2654
HypreParMatrix * ParMult(const HypreParMatrix *A, const HypreParMatrix *B, bool own_matrix)
Definition: hypre.cpp:2849
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:2882
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:94
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:232
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:585
HYPRE_BigInt M() const
Returns the global number of rows.
Definition: hypre.hpp:583
void EliminateZeroRows()
If a row contains only zeros, set its diagonal to 1.
Definition: hypre.hpp:790
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:645
void Set(const int i, const int j, const double val)
Definition: sparsemat.cpp:2764
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:1579
void MultTranspose(const double *x, double *y) const
Multiply a vector with the transpose matrix.
Definition: densemat.cpp:201
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:623
Data type sparse matrix.
Definition: sparsemat.hpp:50
int Height() const
Get the height (size of output) of the Operator. Synonym with NumRows().
Definition: operator.hpp:67
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:656
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:590
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
virtual void Solve(int m, int n, double *X) const
Definition: densemat.cpp:3381
double * data
Definition: densemat.hpp:599
int GetVDim() const
Returns vector dimension.
Definition: fespace.hpp:581
void RAP(const DenseMatrix &A, const DenseMatrix &P, DenseMatrix &RAP)
Definition: densemat.cpp:3213
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:639
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:2783
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:1381
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:96
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:791
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:2809
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:635
HYPRE_BigInt * ColPart()
Returns the column partitioning.
Definition: hypre.hpp:573
void InvScaleRows(const Vector &s)
Scale the local row i by 1./s(i)
Definition: hypre.cpp:2052
HYPRE_BigInt * RowPart()
Returns the row partitioning.
Definition: hypre.hpp:569
IterativeSolver * krylov
int dim
Definition: ex24.cpp:53
virtual bool Factor(int m, double TOL=0.0)
Compute the LU factorization of the current matrix.
Definition: densemat.cpp:3230
constexpr int dimension
This example only works in 3D. Kernels for 2D are not implemented.
Definition: hooke.cpp:45
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:173
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:577
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:1732
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:655
const FiniteElement * GetBE(int i) const
Returns pointer to the FiniteElement in the FiniteElementCollection associated with i&#39;th boundary fac...
Definition: fespace.cpp:3100
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:343
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).