MFEM  v4.5.1
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
navier_solver.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 "navier_solver.hpp"
13 #include "../../general/forall.hpp"
14 #include <fstream>
15 #include <iomanip>
16 
17 using namespace mfem;
18 using namespace navier;
19 
21 {
23  for (int i = 0; i < bffis->Size(); ++i)
24  {
25  dst->AddDomainIntegrator((*bffis)[i]);
26  }
27 }
28 
29 NavierSolver::NavierSolver(ParMesh *mesh, int order, double kin_vis)
30  : pmesh(mesh), order(order), kin_vis(kin_vis),
31  gll_rules(0, Quadrature1D::GaussLobatto)
32 {
33  vfec = new H1_FECollection(order, pmesh->Dimension());
34  pfec = new H1_FECollection(order);
36  pfes = new ParFiniteElementSpace(pmesh, pfec);
37 
38  // Check if fully periodic mesh
39  if (!(pmesh->bdr_attributes.Size() == 0))
40  {
42  vel_ess_attr = 0;
43 
45  pres_ess_attr = 0;
46  }
47 
48  int vfes_truevsize = vfes->GetTrueVSize();
49  int pfes_truevsize = pfes->GetTrueVSize();
50 
51  un.SetSize(vfes_truevsize);
52  un = 0.0;
53  un_next.SetSize(vfes_truevsize);
54  un_next = 0.0;
55  unm1.SetSize(vfes_truevsize);
56  unm1 = 0.0;
57  unm2.SetSize(vfes_truevsize);
58  unm2 = 0.0;
59  fn.SetSize(vfes_truevsize);
60  Nun.SetSize(vfes_truevsize);
61  Nun = 0.0;
62  Nunm1.SetSize(vfes_truevsize);
63  Nunm1 = 0.0;
64  Nunm2.SetSize(vfes_truevsize);
65  Nunm2 = 0.0;
66  Fext.SetSize(vfes_truevsize);
67  FText.SetSize(vfes_truevsize);
68  Lext.SetSize(vfes_truevsize);
69  resu.SetSize(vfes_truevsize);
70 
71  tmp1.SetSize(vfes_truevsize);
72 
73  pn.SetSize(pfes_truevsize);
74  pn = 0.0;
75  resp.SetSize(pfes_truevsize);
76  resp = 0.0;
77  FText_bdr.SetSize(pfes_truevsize);
78  g_bdr.SetSize(pfes_truevsize);
79 
81  un_gf = 0.0;
83  un_next_gf = 0.0;
84 
90 
91  pn_gf.SetSpace(pfes);
92  pn_gf = 0.0;
93  resp_gf.SetSpace(pfes);
94 
95  cur_step = 0;
96 
97  PrintInfo();
98 }
99 
100 void NavierSolver::Setup(double dt)
101 {
102  if (verbose && pmesh->GetMyRank() == 0)
103  {
104  mfem::out << "Setup" << std::endl;
105  if (partial_assembly)
106  {
107  mfem::out << "Using Partial Assembly" << std::endl;
108  }
109  else
110  {
111  mfem::out << "Using Full Assembly" << std::endl;
112  }
113  }
114 
115  sw_setup.Start();
116 
119 
120  Array<int> empty;
121 
122  // GLL integration rule (Numerical Integration)
123  const IntegrationRule &ir_ni = gll_rules.Get(vfes->GetFE(0)->GetGeomType(),
124  2 * order - 1);
125 
126  nlcoeff.constant = -1.0;
127  N = new ParNonlinearForm(vfes);
128  auto *nlc_nlfi = new VectorConvectionNLFIntegrator(nlcoeff);
129  if (numerical_integ)
130  {
131  nlc_nlfi->SetIntRule(&ir_ni);
132  }
133  N->AddDomainIntegrator(nlc_nlfi);
134  if (partial_assembly)
135  {
137  N->Setup();
138  }
139 
141  auto *mv_blfi = new VectorMassIntegrator;
142  if (numerical_integ)
143  {
144  mv_blfi->SetIntRule(&ir_ni);
145  }
146  Mv_form->AddDomainIntegrator(mv_blfi);
147  if (partial_assembly)
148  {
150  }
151  Mv_form->Assemble();
152  Mv_form->FormSystemMatrix(empty, Mv);
153 
155  auto *sp_blfi = new DiffusionIntegrator;
156  if (numerical_integ)
157  {
158  sp_blfi->SetIntRule(&ir_ni);
159  }
160  Sp_form->AddDomainIntegrator(sp_blfi);
161  if (partial_assembly)
162  {
164  }
165  Sp_form->Assemble();
167 
169  auto *vd_mblfi = new VectorDivergenceIntegrator();
170  if (numerical_integ)
171  {
172  vd_mblfi->SetIntRule(&ir_ni);
173  }
174  D_form->AddDomainIntegrator(vd_mblfi);
175  if (partial_assembly)
176  {
178  }
179  D_form->Assemble();
180  D_form->FormRectangularSystemMatrix(empty, empty, D);
181 
183  auto *g_mblfi = new GradientIntegrator();
184  if (numerical_integ)
185  {
186  g_mblfi->SetIntRule(&ir_ni);
187  }
188  G_form->AddDomainIntegrator(g_mblfi);
189  if (partial_assembly)
190  {
192  }
193  G_form->Assemble();
194  G_form->FormRectangularSystemMatrix(empty, empty, G);
195 
197  H_bdfcoeff.constant = 1.0 / dt;
198  H_form = new ParBilinearForm(vfes);
199  auto *hmv_blfi = new VectorMassIntegrator(H_bdfcoeff);
200  auto *hdv_blfi = new VectorDiffusionIntegrator(H_lincoeff);
201  if (numerical_integ)
202  {
203  hmv_blfi->SetIntRule(&ir_ni);
204  hdv_blfi->SetIntRule(&ir_ni);
205  }
206  H_form->AddDomainIntegrator(hmv_blfi);
207  H_form->AddDomainIntegrator(hdv_blfi);
208  if (partial_assembly)
209  {
211  }
212  H_form->Assemble();
214 
217  auto *ftext_bnlfi = new BoundaryNormalLFIntegrator(*FText_gfcoeff);
218  if (numerical_integ)
219  {
220  ftext_bnlfi->SetIntRule(&ir_ni);
221  }
223 
225  for (auto &vel_dbc : vel_dbcs)
226  {
227  auto *gbdr_bnlfi = new BoundaryNormalLFIntegrator(*vel_dbc.coeff);
228  if (numerical_integ)
229  {
230  gbdr_bnlfi->SetIntRule(&ir_ni);
231  }
232  g_bdr_form->AddBoundaryIntegrator(gbdr_bnlfi, vel_dbc.attr);
233  }
234 
235  f_form = new ParLinearForm(vfes);
236  for (auto &accel_term : accel_terms)
237  {
238  auto *vdlfi = new VectorDomainLFIntegrator(*accel_term.coeff);
239  // @TODO: This order should always be the same as the nonlinear forms one!
240  // const IntegrationRule &ir = IntRules.Get(vfes->GetFE(0)->GetGeomType(),
241  // 4 * order);
242  // vdlfi->SetIntRule(&ir);
243  if (numerical_integ)
244  {
245  vdlfi->SetIntRule(&ir_ni);
246  }
247  f_form->AddDomainIntegrator(vdlfi);
248  }
249 
250  if (partial_assembly)
251  {
252  Vector diag_pa(vfes->GetTrueVSize());
253  Mv_form->AssembleDiagonal(diag_pa);
254  MvInvPC = new OperatorJacobiSmoother(diag_pa, empty);
255  }
256  else
257  {
259  dynamic_cast<HypreSmoother *>(MvInvPC)->SetType(HypreSmoother::Jacobi, 1);
260  }
261  MvInv = new CGSolver(vfes->GetComm());
262  MvInv->iterative_mode = false;
263  MvInv->SetOperator(*Mv);
266  MvInv->SetRelTol(1e-12);
267  MvInv->SetMaxIter(200);
268 
269  if (partial_assembly)
270  {
274  SpInvPC->Mult(resp, pn);
277  }
278  else
279  {
284  }
285  SpInv = new CGSolver(vfes->GetComm());
286  SpInv->iterative_mode = true;
287  SpInv->SetOperator(*Sp);
288  if (pres_dbcs.empty())
289  {
291  }
292  else
293  {
295  }
298  SpInv->SetMaxIter(200);
299 
300  if (partial_assembly)
301  {
302  Vector diag_pa(vfes->GetTrueVSize());
303  H_form->AssembleDiagonal(diag_pa);
305  }
306  else
307  {
309  dynamic_cast<HypreSmoother *>(HInvPC)->SetType(HypreSmoother::Jacobi, 1);
310  }
311  HInv = new CGSolver(vfes->GetComm());
312  HInv->iterative_mode = true;
313  HInv->SetOperator(*H);
317  HInv->SetMaxIter(200);
318 
319  // If the initial condition was set, it has to be aligned with dependent
320  // Vectors and GridFunctions
322  un_next = un;
324 
325  // Set initial time step in the history array
326  dthist[0] = dt;
327 
328  if (filter_alpha != 0.0)
329  {
331  pmesh->Dimension());
333  vfec_filter,
334  pmesh->Dimension());
335 
336  un_NM1_gf.SetSpace(vfes_filter);
337  un_NM1_gf = 0.0;
338 
340  un_filtered_gf = 0.0;
341  }
342 
343  sw_setup.Stop();
344 }
345 
347 {
348  // Rotate values in time step history
349  dthist[2] = dthist[1];
350  dthist[1] = dthist[0];
351  dthist[0] = dt;
352 
353  // Rotate values in nonlinear extrapolation history
354  Nunm2 = Nunm1;
355  Nunm1 = Nun;
356 
357  // Rotate values in solution history
358  unm2 = unm1;
359  unm1 = un;
360 
361  // Update the current solution and corresponding GridFunction
363  un = un_next;
365 }
366 
367 void NavierSolver::Step(double &time, double dt, int current_step,
368  bool provisional)
369 {
370  sw_step.Start();
371 
372  SetTimeIntegrationCoefficients(current_step);
373 
374  // Set current time for velocity Dirichlet boundary conditions.
375  for (auto &vel_dbc : vel_dbcs)
376  {
377  vel_dbc.coeff->SetTime(time + dt);
378  }
379 
380  // Set current time for pressure Dirichlet boundary conditions.
381  for (auto &pres_dbc : pres_dbcs)
382  {
383  pres_dbc.coeff->SetTime(time + dt);
384  }
385 
386  H_bdfcoeff.constant = bd0 / dt;
387  H_form->Update();
388  H_form->Assemble();
390 
391  HInv->SetOperator(*H);
392  if (partial_assembly)
393  {
394  delete HInvPC;
395  Vector diag_pa(vfes->GetTrueVSize());
396  H_form->AssembleDiagonal(diag_pa);
399  }
400 
401  // Extrapolated f^{n+1}.
402  for (auto &accel_term : accel_terms)
403  {
404  accel_term.coeff->SetTime(time + dt);
405  }
406 
407  f_form->Assemble();
409 
410  // Nonlinear extrapolated terms.
411  sw_extrap.Start();
412 
413  N->Mult(un, Nun);
414  N->Mult(unm1, Nunm1);
415  N->Mult(unm2, Nunm2);
416 
417  {
418  const auto d_Nun = Nun.Read();
419  const auto d_Nunm1 = Nunm1.Read();
420  const auto d_Nunm2 = Nunm2.Read();
421  auto d_Fext = Fext.Write();
422  const auto ab1_ = ab1;
423  const auto ab2_ = ab2;
424  const auto ab3_ = ab3;
425  MFEM_FORALL(i, Fext.Size(),
426  d_Fext[i] = ab1_ * d_Nun[i] +
427  ab2_ * d_Nunm1[i] +
428  ab3_ * d_Nunm2[i];);
429  }
430 
431  Fext.Add(1.0, fn);
432 
433  // Fext = M^{-1} (F(u^{n}) + f^{n+1})
434  MvInv->Mult(Fext, tmp1);
437  Fext.Set(1.0, tmp1);
438 
439  // Compute BDF terms.
440  {
441  const double bd1idt = -bd1 / dt;
442  const double bd2idt = -bd2 / dt;
443  const double bd3idt = -bd3 / dt;
444  const auto d_un = un.Read();
445  const auto d_unm1 = unm1.Read();
446  const auto d_unm2 = unm2.Read();
447  auto d_Fext = Fext.ReadWrite();
448  MFEM_FORALL(i, Fext.Size(),
449  d_Fext[i] += bd1idt * d_un[i] +
450  bd2idt * d_unm1[i] +
451  bd3idt * d_unm2[i];);
452  }
453 
454  sw_extrap.Stop();
455 
456  // Pressure Poisson.
457  sw_curlcurl.Start();
458  {
459  const auto d_un = un.Read();
460  const auto d_unm1 = unm1.Read();
461  const auto d_unm2 = unm2.Read();
462  auto d_Lext = Lext.Write();
463  const auto ab1_ = ab1;
464  const auto ab2_ = ab2;
465  const auto ab3_ = ab3;
466  MFEM_FORALL(i, Lext.Size(),
467  d_Lext[i] = ab1_ * d_un[i] +
468  ab2_ * d_unm1[i] +
469  ab3_ * d_unm2[i];);
470  }
471 
473  if (pmesh->Dimension() == 2)
474  {
477  }
478  else
479  {
482  }
483 
485  Lext *= kin_vis;
486 
487  sw_curlcurl.Stop();
488 
489  // \tilde{F} = F - \nu CurlCurl(u)
490  FText.Set(-1.0, Lext);
491  FText.Add(1.0, Fext);
492 
493  // p_r = \nabla \cdot FText
494  D->Mult(FText, resp);
495  resp.Neg();
496 
497  // Add boundary terms.
501 
502  g_bdr_form->Assemble();
504  resp.Add(1.0, FText_bdr);
505  resp.Add(-bd0 / dt, g_bdr);
506 
507  if (pres_dbcs.empty())
508  {
510  }
511 
512  for (auto &pres_dbc : pres_dbcs)
513  {
514  pn_gf.ProjectBdrCoefficient(*pres_dbc.coeff, pres_dbc.attr);
515  }
516 
518 
519  Vector X1, B1;
520  if (partial_assembly)
521  {
522  auto *SpC = Sp.As<ConstrainedOperator>();
523  EliminateRHS(*Sp_form, *SpC, pres_ess_tdof, pn_gf, resp_gf, X1, B1, 1);
524  }
525  else
526  {
528  }
529  sw_spsolve.Start();
530  SpInv->Mult(B1, X1);
531  sw_spsolve.Stop();
535 
536  // If the boundary conditions on the pressure are pure Neumann remove the
537  // nullspace by removing the mean of the pressure solution. This is also
538  // ensured by the OrthoSolver wrapper for the preconditioner which removes
539  // the nullspace after every application.
540  if (pres_dbcs.empty())
541  {
542  MeanZero(pn_gf);
543  }
544 
546 
547  // Project velocity.
548  G->Mult(pn, resu);
549  resu.Neg();
550  Mv->Mult(Fext, tmp1);
551  resu.Add(1.0, tmp1);
552 
553  // un_next_gf = un_gf;
554 
555  for (auto &vel_dbc : vel_dbcs)
556  {
558  }
559 
561 
562  Vector X2, B2;
563  if (partial_assembly)
564  {
565  auto *HC = H.As<ConstrainedOperator>();
566  EliminateRHS(*H_form, *HC, vel_ess_tdof, un_next_gf, resu_gf, X2, B2, 1);
567  }
568  else
569  {
571  }
572  sw_hsolve.Start();
573  HInv->Mult(B2, X2);
574  sw_hsolve.Stop();
578 
580 
581  // If the current time step is not provisional, accept the computed solution
582  // and update the time step history by default.
583  if (!provisional)
584  {
586  time += dt;
587  }
588 
589  if (filter_alpha != 0.0)
590  {
593  const auto d_un_filtered_gf = un_filtered_gf.Read();
594  auto d_un_gf = un_gf.ReadWrite();
595  const auto filter_alpha_ = filter_alpha;
596  MFEM_FORALL(i,
597  un_gf.Size(),
598  d_un_gf[i] = (1.0 - filter_alpha_) * d_un_gf[i]
599  + filter_alpha_ * d_un_filtered_gf[i];);
600  }
601 
602  sw_step.Stop();
603 
604  if (verbose && pmesh->GetMyRank() == 0)
605  {
606  mfem::out << std::setw(7) << "" << std::setw(3) << "It" << std::setw(8)
607  << "Resid" << std::setw(12) << "Reltol"
608  << "\n";
609  // If numerical integration is active, there is no solve (thus no
610  // iterations), on the inverse velocity mass application.
611  if (!numerical_integ)
612  {
613  mfem::out << std::setw(5) << "MVIN " << std::setw(5) << std::fixed
614  << iter_mvsolve << " " << std::setw(3)
615  << std::setprecision(2) << std::scientific << res_mvsolve
616  << " " << 1e-12 << "\n";
617  }
618  mfem::out << std::setw(5) << "PRES " << std::setw(5) << std::fixed
619  << iter_spsolve << " " << std::setw(3) << std::setprecision(2)
620  << std::scientific << res_spsolve << " " << rtol_spsolve
621  << "\n";
622  mfem::out << std::setw(5) << "HELM " << std::setw(5) << std::fixed
623  << iter_hsolve << " " << std::setw(3) << std::setprecision(2)
624  << std::scientific << res_hsolve << " " << rtol_hsolve
625  << "\n";
626  mfem::out << std::setprecision(8);
627  mfem::out << std::fixed;
628  }
629 }
630 
632 {
633  // Make sure not to recompute the inner product linear form every
634  // application.
635  if (mass_lf == nullptr)
636  {
637  onecoeff.constant = 1.0;
638  mass_lf = new ParLinearForm(v.ParFESpace());
639  auto *dlfi = new DomainLFIntegrator(onecoeff);
640  if (numerical_integ)
641  {
642  const IntegrationRule &ir_ni = gll_rules.Get(vfes->GetFE(0)->GetGeomType(),
643  2 * order - 1);
644  dlfi->SetIntRule(&ir_ni);
645  }
647  mass_lf->Assemble();
648 
649  ParGridFunction one_gf(v.ParFESpace());
651 
652  volume = mass_lf->operator()(one_gf);
653  }
654 
655  double integ = mass_lf->operator()(v);
656 
657  v -= integ / volume;
658 }
659 
661  ConstrainedOperator &constrainedA,
662  const Array<int> &ess_tdof_list,
663  Vector &x,
664  Vector &b,
665  Vector &X,
666  Vector &B,
667  int copy_interior)
668 {
669  const Operator *Po = A.GetOutputProlongation();
670  const Operator *Pi = A.GetProlongation();
671  const Operator *Ri = A.GetRestriction();
672  A.InitTVectors(Po, Ri, Pi, x, b, X, B);
673  if (!copy_interior)
674  {
675  X.SetSubVectorComplement(ess_tdof_list, 0.0);
676  }
677  constrainedA.EliminateRHS(X, B);
678 }
679 
681 {
682  double loc_sum = v.Sum();
683  double global_sum = 0.0;
684  int loc_size = v.Size();
685  int global_size = 0;
686 
687  MPI_Allreduce(&loc_sum, &global_sum, 1, MPI_DOUBLE, MPI_SUM, pfes->GetComm());
688  MPI_Allreduce(&loc_size, &global_size, 1, MPI_INT, MPI_SUM, pfes->GetComm());
689 
690  v -= global_sum / static_cast<double>(global_size);
691 }
692 
694 {
695  FiniteElementSpace *fes = u.FESpace();
696 
697  // AccumulateAndCountZones.
698  Array<int> zones_per_vdof;
699  zones_per_vdof.SetSize(fes->GetVSize());
700  zones_per_vdof = 0;
701 
702  cu = 0.0;
703 
704  // Local interpolation.
705  int elndofs;
706  Array<int> vdofs;
707  Vector vals;
708  Vector loc_data;
709  int vdim = fes->GetVDim();
710  DenseMatrix grad_hat;
711  DenseMatrix dshape;
712  DenseMatrix grad;
713  Vector curl;
714 
715  for (int e = 0; e < fes->GetNE(); ++e)
716  {
717  fes->GetElementVDofs(e, vdofs);
718  u.GetSubVector(vdofs, loc_data);
719  vals.SetSize(vdofs.Size());
721  const FiniteElement *el = fes->GetFE(e);
722  elndofs = el->GetDof();
723  int dim = el->GetDim();
724  dshape.SetSize(elndofs, dim);
725 
726  for (int dof = 0; dof < elndofs; ++dof)
727  {
728  // Project.
729  const IntegrationPoint &ip = el->GetNodes().IntPoint(dof);
730  tr->SetIntPoint(&ip);
731 
732  // Eval and GetVectorGradientHat.
733  el->CalcDShape(tr->GetIntPoint(), dshape);
734  grad_hat.SetSize(vdim, dim);
735  DenseMatrix loc_data_mat(loc_data.GetData(), elndofs, vdim);
736  MultAtB(loc_data_mat, dshape, grad_hat);
737 
738  const DenseMatrix &Jinv = tr->InverseJacobian();
739  grad.SetSize(grad_hat.Height(), Jinv.Width());
740  Mult(grad_hat, Jinv, grad);
741 
742  curl.SetSize(3);
743  curl(0) = grad(2, 1) - grad(1, 2);
744  curl(1) = grad(0, 2) - grad(2, 0);
745  curl(2) = grad(1, 0) - grad(0, 1);
746 
747  for (int j = 0; j < curl.Size(); ++j)
748  {
749  vals(elndofs * j + dof) = curl(j);
750  }
751  }
752 
753  // Accumulate values in all dofs, count the zones.
754  for (int j = 0; j < vdofs.Size(); j++)
755  {
756  int ldof = vdofs[j];
757  cu(ldof) += vals[j];
758  zones_per_vdof[ldof]++;
759  }
760  }
761 
762  // Communication
763 
764  // Count the zones globally.
765  GroupCommunicator &gcomm = u.ParFESpace()->GroupComm();
766  gcomm.Reduce<int>(zones_per_vdof, GroupCommunicator::Sum);
767  gcomm.Bcast(zones_per_vdof);
768 
769  // Accumulate for all vdofs.
770  gcomm.Reduce<double>(cu.GetData(), GroupCommunicator::Sum);
771  gcomm.Bcast<double>(cu.GetData());
772 
773  // Compute means.
774  for (int i = 0; i < cu.Size(); i++)
775  {
776  const int nz = zones_per_vdof[i];
777  if (nz)
778  {
779  cu(i) /= nz;
780  }
781  }
782 }
783 
785  ParGridFunction &cu,
786  bool assume_scalar)
787 {
788  FiniteElementSpace *fes = u.FESpace();
789 
790  // AccumulateAndCountZones.
791  Array<int> zones_per_vdof;
792  zones_per_vdof.SetSize(fes->GetVSize());
793  zones_per_vdof = 0;
794 
795  cu = 0.0;
796 
797  // Local interpolation.
798  int elndofs;
799  Array<int> vdofs;
800  Vector vals;
801  Vector loc_data;
802  int vdim = fes->GetVDim();
803  DenseMatrix grad_hat;
804  DenseMatrix dshape;
805  DenseMatrix grad;
806  Vector curl;
807 
808  for (int e = 0; e < fes->GetNE(); ++e)
809  {
810  fes->GetElementVDofs(e, vdofs);
811  u.GetSubVector(vdofs, loc_data);
812  vals.SetSize(vdofs.Size());
814  const FiniteElement *el = fes->GetFE(e);
815  elndofs = el->GetDof();
816  int dim = el->GetDim();
817  dshape.SetSize(elndofs, dim);
818 
819  for (int dof = 0; dof < elndofs; ++dof)
820  {
821  // Project.
822  const IntegrationPoint &ip = el->GetNodes().IntPoint(dof);
823  tr->SetIntPoint(&ip);
824 
825  // Eval and GetVectorGradientHat.
826  el->CalcDShape(tr->GetIntPoint(), dshape);
827  grad_hat.SetSize(vdim, dim);
828  DenseMatrix loc_data_mat(loc_data.GetData(), elndofs, vdim);
829  MultAtB(loc_data_mat, dshape, grad_hat);
830 
831  const DenseMatrix &Jinv = tr->InverseJacobian();
832  grad.SetSize(grad_hat.Height(), Jinv.Width());
833  Mult(grad_hat, Jinv, grad);
834 
835  if (assume_scalar)
836  {
837  curl.SetSize(2);
838  curl(0) = grad(0, 1);
839  curl(1) = -grad(0, 0);
840  }
841  else
842  {
843  curl.SetSize(2);
844  curl(0) = grad(1, 0) - grad(0, 1);
845  curl(1) = 0.0;
846  }
847 
848  for (int j = 0; j < curl.Size(); ++j)
849  {
850  vals(elndofs * j + dof) = curl(j);
851  }
852  }
853 
854  // Accumulate values in all dofs, count the zones.
855  for (int j = 0; j < vdofs.Size(); j++)
856  {
857  int ldof = vdofs[j];
858  cu(ldof) += vals[j];
859  zones_per_vdof[ldof]++;
860  }
861  }
862 
863  // Communication.
864 
865  // Count the zones globally.
866  GroupCommunicator &gcomm = u.ParFESpace()->GroupComm();
867  gcomm.Reduce<int>(zones_per_vdof, GroupCommunicator::Sum);
868  gcomm.Bcast(zones_per_vdof);
869 
870  // Accumulate for all vdofs.
871  gcomm.Reduce<double>(cu.GetData(), GroupCommunicator::Sum);
872  gcomm.Bcast<double>(cu.GetData());
873 
874  // Compute means.
875  for (int i = 0; i < cu.Size(); i++)
876  {
877  const int nz = zones_per_vdof[i];
878  if (nz)
879  {
880  cu(i) /= nz;
881  }
882  }
883 }
884 
886 {
887  ParMesh *pmesh_u = u.ParFESpace()->GetParMesh();
888  FiniteElementSpace *fes = u.FESpace();
889  int vdim = fes->GetVDim();
890 
891  Vector ux, uy, uz;
892  Vector ur, us, ut;
893  double cflx = 0.0;
894  double cfly = 0.0;
895  double cflz = 0.0;
896  double cflm = 0.0;
897  double cflmax = 0.0;
898 
899  for (int e = 0; e < fes->GetNE(); ++e)
900  {
901  const FiniteElement *fe = fes->GetFE(e);
902  const IntegrationRule &ir = IntRules.Get(fe->GetGeomType(),
903  fe->GetOrder());
905 
906  u.GetValues(e, ir, ux, 1);
907  ur.SetSize(ux.Size());
908  u.GetValues(e, ir, uy, 2);
909  us.SetSize(uy.Size());
910  if (vdim == 3)
911  {
912  u.GetValues(e, ir, uz, 3);
913  ut.SetSize(uz.Size());
914  }
915 
916  double hmin = pmesh_u->GetElementSize(e, 1) /
917  (double) fes->GetElementOrder(0);
918 
919  for (int i = 0; i < ir.GetNPoints(); ++i)
920  {
921  const IntegrationPoint &ip = ir.IntPoint(i);
922  tr->SetIntPoint(&ip);
923  const DenseMatrix &invJ = tr->InverseJacobian();
924  const double detJinv = 1.0 / tr->Jacobian().Det();
925 
926  if (vdim == 2)
927  {
928  ur(i) = (ux(i) * invJ(0, 0) + uy(i) * invJ(1, 0)) * detJinv;
929  us(i) = (ux(i) * invJ(0, 1) + uy(i) * invJ(1, 1)) * detJinv;
930  }
931  else if (vdim == 3)
932  {
933  ur(i) = (ux(i) * invJ(0, 0) + uy(i) * invJ(1, 0)
934  + uz(i) * invJ(2, 0))
935  * detJinv;
936  us(i) = (ux(i) * invJ(0, 1) + uy(i) * invJ(1, 1)
937  + uz(i) * invJ(2, 1))
938  * detJinv;
939  ut(i) = (ux(i) * invJ(0, 2) + uy(i) * invJ(1, 2)
940  + uz(i) * invJ(2, 2))
941  * detJinv;
942  }
943 
944  cflx = fabs(dt * ux(i) / hmin);
945  cfly = fabs(dt * uy(i) / hmin);
946  if (vdim == 3)
947  {
948  cflz = fabs(dt * uz(i) / hmin);
949  }
950  cflm = cflx + cfly + cflz;
951  cflmax = fmax(cflmax, cflm);
952  }
953  }
954 
955  double cflmax_global = 0.0;
956  MPI_Allreduce(&cflmax,
957  &cflmax_global,
958  1,
959  MPI_DOUBLE,
960  MPI_MAX,
961  pmesh_u->GetComm());
962 
963  return cflmax_global;
964 }
965 
967 {
968  vel_dbcs.emplace_back(attr, coeff);
969 
970  if (verbose && pmesh->GetMyRank() == 0)
971  {
972  mfem::out << "Adding Velocity Dirichlet BC to attributes ";
973  for (int i = 0; i < attr.Size(); ++i)
974  {
975  if (attr[i] == 1)
976  {
977  mfem::out << i << " ";
978  }
979  }
980  mfem::out << std::endl;
981  }
982 
983  for (int i = 0; i < attr.Size(); ++i)
984  {
985  MFEM_ASSERT((vel_ess_attr[i] && attr[i]) == 0,
986  "Duplicate boundary definition deteceted.");
987  if (attr[i] == 1)
988  {
989  vel_ess_attr[i] = 1;
990  }
991  }
992 }
993 
995 {
997 }
998 
1000 {
1001  pres_dbcs.emplace_back(attr, coeff);
1002 
1003  if (verbose && pmesh->GetMyRank() == 0)
1004  {
1005  mfem::out << "Adding Pressure Dirichlet BC to attributes ";
1006  for (int i = 0; i < attr.Size(); ++i)
1007  {
1008  if (attr[i] == 1)
1009  {
1010  mfem::out << i << " ";
1011  }
1012  }
1013  mfem::out << std::endl;
1014  }
1015 
1016  for (int i = 0; i < attr.Size(); ++i)
1017  {
1018  MFEM_ASSERT((pres_ess_attr[i] && attr[i]) == 0,
1019  "Duplicate boundary definition deteceted.");
1020  if (attr[i] == 1)
1021  {
1022  pres_ess_attr[i] = 1;
1023  }
1024  }
1025 }
1026 
1028 {
1030 }
1031 
1033 {
1034  accel_terms.emplace_back(attr, coeff);
1035 
1036  if (verbose && pmesh->GetMyRank() == 0)
1037  {
1038  mfem::out << "Adding Acceleration term to attributes ";
1039  for (int i = 0; i < attr.Size(); ++i)
1040  {
1041  if (attr[i] == 1)
1042  {
1043  mfem::out << i << " ";
1044  }
1045  }
1046  mfem::out << std::endl;
1047  }
1048 }
1049 
1051 {
1053 }
1054 
1056 {
1057  // Maximum BDF order to use at current time step
1058  // step + 1 <= order <= max_bdf_order
1059  int bdf_order = std::min(step + 1, max_bdf_order);
1060 
1061  // Ratio of time step history at dt(t_{n}) - dt(t_{n-1})
1062  double rho1 = 0.0;
1063 
1064  // Ratio of time step history at dt(t_{n-1}) - dt(t_{n-2})
1065  double rho2 = 0.0;
1066 
1067  rho1 = dthist[0] / dthist[1];
1068 
1069  if (bdf_order == 3)
1070  {
1071  rho2 = dthist[1] / dthist[2];
1072  }
1073 
1074  if (step == 0 && bdf_order == 1)
1075  {
1076  bd0 = 1.0;
1077  bd1 = -1.0;
1078  bd2 = 0.0;
1079  bd3 = 0.0;
1080  ab1 = 1.0;
1081  ab2 = 0.0;
1082  ab3 = 0.0;
1083  }
1084  else if (step >= 1 && bdf_order == 2)
1085  {
1086  bd0 = (1.0 + 2.0 * rho1) / (1.0 + rho1);
1087  bd1 = -(1.0 + rho1);
1088  bd2 = pow(rho1, 2.0) / (1.0 + rho1);
1089  bd3 = 0.0;
1090  ab1 = 1.0 + rho1;
1091  ab2 = -rho1;
1092  ab3 = 0.0;
1093  }
1094  else if (step >= 2 && bdf_order == 3)
1095  {
1096  bd0 = 1.0 + rho1 / (1.0 + rho1)
1097  + (rho2 * rho1) / (1.0 + rho2 * (1 + rho1));
1098  bd1 = -1.0 - rho1 - (rho2 * rho1 * (1.0 + rho1)) / (1.0 + rho2);
1099  bd2 = pow(rho1, 2.0) * (rho2 + 1.0 / (1.0 + rho1));
1100  bd3 = -(pow(rho2, 3.0) * pow(rho1, 2.0) * (1.0 + rho1))
1101  / ((1.0 + rho2) * (1.0 + rho2 + rho2 * rho1));
1102  ab1 = ((1.0 + rho1) * (1.0 + rho2 * (1.0 + rho1))) / (1.0 + rho2);
1103  ab2 = -rho1 * (1.0 + rho2 * (1.0 + rho1));
1104  ab3 = (pow(rho2, 2.0) * rho1 * (1.0 + rho1)) / (1.0 + rho2);
1105  }
1106 }
1107 
1109 {
1110  double my_rt[6], rt_max[6];
1111 
1112  my_rt[0] = sw_setup.RealTime();
1113  my_rt[1] = sw_step.RealTime();
1114  my_rt[2] = sw_extrap.RealTime();
1115  my_rt[3] = sw_curlcurl.RealTime();
1116  my_rt[4] = sw_spsolve.RealTime();
1117  my_rt[5] = sw_hsolve.RealTime();
1118 
1119  MPI_Reduce(my_rt, rt_max, 6, MPI_DOUBLE, MPI_MAX, 0, pmesh->GetComm());
1120 
1121  if (pmesh->GetMyRank() == 0)
1122  {
1123  mfem::out << std::setw(10) << "SETUP" << std::setw(10) << "STEP"
1124  << std::setw(10) << "EXTRAP" << std::setw(10) << "CURLCURL"
1125  << std::setw(10) << "PSOLVE" << std::setw(10) << "HSOLVE"
1126  << "\n";
1127 
1128  mfem::out << std::setprecision(3) << std::setw(10) << my_rt[0]
1129  << std::setw(10) << my_rt[1] << std::setw(10) << my_rt[2]
1130  << std::setw(10) << my_rt[3] << std::setw(10) << my_rt[4]
1131  << std::setw(10) << my_rt[5] << "\n";
1132 
1133  mfem::out << std::setprecision(3) << std::setw(10) << " " << std::setw(10)
1134  << my_rt[1] / my_rt[1] << std::setw(10) << my_rt[2] / my_rt[1]
1135  << std::setw(10) << my_rt[3] / my_rt[1] << std::setw(10)
1136  << my_rt[4] / my_rt[1] << std::setw(10) << my_rt[5] / my_rt[1]
1137  << "\n";
1138 
1139  mfem::out << std::setprecision(8);
1140  }
1141 }
1142 
1144 {
1145  int fes_size0 = vfes->GlobalVSize();
1146  int fes_size1 = pfes->GlobalVSize();
1147 
1148  if (pmesh->GetMyRank() == 0)
1149  {
1150  mfem::out << "NAVIER version: " << NAVIER_VERSION << std::endl
1151  << "MFEM version: " << MFEM_VERSION << std::endl
1152  << "MFEM GIT: " << MFEM_GIT_STRING << std::endl
1153  << "Velocity #DOFs: " << fes_size0 << std::endl
1154  << "Pressure #DOFs: " << fes_size1 << std::endl;
1155  }
1156 }
1157 
1159 {
1160  delete FText_gfcoeff;
1161  delete g_bdr_form;
1162  delete FText_bdr_form;
1163  delete mass_lf;
1164  delete Mv_form;
1165  delete N;
1166  delete Sp_form;
1167  delete D_form;
1168  delete G_form;
1169  delete HInvPC;
1170  delete HInv;
1171  delete H_form;
1172  delete SpInv;
1173  delete MvInvPC;
1174  delete SpInvOrthoPC;
1175  delete SpInvPC;
1176  delete lor;
1177  delete f_form;
1178  delete MvInv;
1179  delete vfec;
1180  delete pfec;
1181  delete vfes;
1182  delete pfes;
1183  delete vfec_filter;
1184  delete vfes_filter;
1185 }
int GetNPoints() const
Returns the number of the points in the integration rule.
Definition: intrules.hpp:247
Abstract class for all finite elements.
Definition: fe_base.hpp:235
FiniteElementCollection * vfec_filter
void PrintTimingData()
Print timing summary of the solving routine.
int Size() const
Return the logical size of the array.
Definition: array.hpp:138
void Orthogonalize(Vector &v)
Remove mean from a Vector.
Class for domain integration L(v) := (f, v)
Definition: lininteg.hpp:108
void Reduce(T *ldata, void(*Op)(OpData< T >)) const
Reduce within each group where the master is the root.
int GetVSize() const
Return the number of vector dofs, i.e. GetNDofs() x GetVDim().
Definition: fespace.hpp:587
virtual void GetEssentialTrueDofs(const Array< int > &bdr_attr_is_ess, Array< int > &ess_tdof_list, int component=-1)
Definition: pfespace.cpp:1032
Conjugate gradient method.
Definition: solvers.hpp:465
int GetDim() const
Returns the reference space dimension for the finite element.
Definition: fe_base.hpp:311
FiniteElementCollection * pfec
Pressure finite element collection.
OpType * As() const
Return the Operator pointer statically cast to a specified OpType. Similar to the method Get()...
Definition: handle.hpp:104
Create and assemble a low-order refined version of a ParBilinearForm.
Definition: lor.hpp:165
Class for an integration rule - an Array of IntegrationPoint.
Definition: intrules.hpp:90
void EliminateRHS(Operator &A, ConstrainedOperator &constrainedA, const Array< int > &ess_tdof_list, Vector &x, Vector &b, Vector &X, Vector &B, int copy_interior=0)
Eliminate essential BCs in an Operator and apply to RHS.
virtual void AssembleDiagonal(Vector &diag) const
Assemble the diagonal of the bilinear form into diag. Note that diag is a true-dof Vector...
void Bcast(T *ldata, int layout) const
Broadcast within each group where the master is the root.
int GetNumIterations() const
Definition: solvers.hpp:246
const IntegrationRule & Get(int GeomType, int Order)
Returns an integration rule for given GeomType and Order.
Definition: intrules.cpp:923
ParMesh * GetParMesh() const
Definition: pfespace.hpp:277
Base class for vector Coefficients that optionally depend on time and space.
virtual void FormLinearSystem(const Array< int > &ess_tdof_list, Vector &x, Vector &b, OperatorHandle &A, Vector &X, Vector &B, int copy_interior=0)
Form the linear system A X = B, corresponding to this bilinear form and the linear form b(...
void PrintInfo()
Print information about the Navier version.
virtual void Mult(const Vector &b, Vector &x) const
Operator application: y=A(x).
Definition: solvers.cpp:711
Array< BilinearFormIntegrator * > * GetDBFI()
Access all the integrators added with AddDomainIntegrator().
void SetSize(int s)
Resize the vector to size s.
Definition: vector.hpp:513
virtual void Setup()
Setup the NonlinearForm: based on the current AssemblyLevel and the current mesh, optionally...
double Det() const
Definition: densemat.cpp:449
ConstantCoefficient nlcoeff
void Mult(const Table &A, const Table &B, Table &C)
C = A * B (as boolean matrices)
Definition: table.cpp:475
int Width() const
Get the width (size of input) of the Operator. Synonym with NumCols().
Definition: operator.hpp:73
void AddDomainIntegrator(NonlinearFormIntegrator *nlfi)
Adds new Domain Integrator.
void(const Vector &x, double t, Vector &u) VecFuncT
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
int GetOrder() const
Returns the order of the finite element. In the case of anisotropic orders, returns the maximum order...
Definition: fe_base.hpp:327
Data type dense matrix using column-major storage.
Definition: densemat.hpp:23
ParGridFunction un_filtered_gf
int Size() const
Returns the size of the vector.
Definition: vector.hpp:200
Parallel non-linear operator on the true dofs.
virtual const Operator * GetOutputProlongation() const
Prolongation operator from linear algebra (linear system) vectors, to output vectors for the operator...
Definition: operator.hpp:122
int GetElementOrder(int i) const
Returns the order of the i&#39;th finite element.
Definition: fespace.cpp:178
virtual void SetIntRule(const IntegrationRule *ir)
Prescribe a fixed IntegrationRule to use (when ir != NULL) or let the integrator choose (when ir == N...
Definition: nonlininteg.hpp:43
double(const Vector &x, double t) ScalarFuncT
ParLinearForm * FText_bdr_form
std::vector< PresDirichletBC_T > pres_dbcs
const DenseMatrix & InverseJacobian()
Return the inverse of the Jacobian matrix of the transformation at the currently set IntegrationPoint...
Definition: eltrans.hpp:145
bool numerical_integ
Enable/disable numerical integration rules of forms.
virtual void Mult(const Vector &x, Vector &y) const =0
Operator application: y=A(x).
void SetAssemblyLevel(AssemblyLevel assembly_level)
Set the desired assembly level.
double RealTime()
Definition: tic_toc.cpp:426
virtual void SetSpace(FiniteElementSpace *f)
Associate a new FiniteElementSpace with the ParGridFunction.
Definition: pgridfunc.cpp:87
Abstract parallel finite element space.
Definition: pfespace.hpp:28
ConstantCoefficient H_lincoeff
virtual void ProjectCoefficient(Coefficient &coeff)
Project coeff Coefficient to this GridFunction. The projection computation depends on the choice of t...
Definition: pgridfunc.cpp:525
bool iterative_mode
If true, use the second argument of Mult() as an initial guess.
Definition: operator.hpp:659
A class container for 1D quadrature type constants.
Definition: intrules.hpp:289
const IntegrationPoint & GetIntPoint()
Get a const reference to the currently set integration point. This will return NULL if no integration...
Definition: eltrans.hpp:98
double * GetData() const
Return a pointer to the beginning of the Vector data.
Definition: vector.hpp:209
NavierSolver(ParMesh *mesh, int order, double kin_vis)
Initialize data structures, set FE space order and kinematic viscosity.
HypreParMatrix & GetAssembledMatrix() const
Return the assembled LOR operator as a HypreParMatrix.
Definition: lor.cpp:522
ParMixedBilinearForm * G_form
ParLORDiscretization * lor
virtual void Update(FiniteElementSpace *nfes=NULL)
Update the FiniteElementSpace and delete all data associated with the old one.
double GetFinalNorm() const
Definition: solvers.hpp:248
void Setup(double dt)
Initialize forms, solvers and preconditioners.
void MeanZero(ParGridFunction &v)
Remove the mean from a ParGridFunction.
virtual int GetTrueVSize() const
Return the number of local vector true dofs.
Definition: pfespace.hpp:289
void Stop()
Stop the stopwatch.
Definition: tic_toc.cpp:416
bool partial_assembly
Enable/disable partial assembly of forms.
void AddDomainIntegrator(BilinearFormIntegrator *bfi)
Adds a domain integrator. Assumes ownership of bfi.
Geometry::Type GetGeomType() const
Returns the Geometry::Type of the reference element.
Definition: fe_base.hpp:320
virtual const Operator * GetProlongation() const
Prolongation operator from linear algebra (linear system) vectors, to input vectors for the operator...
Definition: operator.hpp:116
int GetNE() const
Returns number of elements in the mesh.
Definition: fespace.hpp:614
The BoomerAMG solver in hypre.
Definition: hypre.hpp:1579
IntegrationPoint & IntPoint(int i)
Returns a reference to the i-th integration point.
Definition: intrules.hpp:250
void ComputeCurl3D(ParGridFunction &u, ParGridFunction &cu)
Compute for .
Class for parallel linear form.
Definition: plinearform.hpp:26
HYPRE_BigInt GlobalVSize() const
Definition: pfespace.hpp:283
Communicator performing operations within groups defined by a GroupTopology with arbitrary-size data ...
DofTransformation * GetElementVDofs(int i, Array< int > &vdofs) const
Returns indexes of degrees of freedom in array dofs for i&#39;th element.
Definition: fespace.cpp:281
virtual void SetPrintLevel(int print_lvl)
Legacy method to set the level of verbosity of the solver output.
Definition: solvers.cpp:71
FiniteElementCollection * vfec
Velocity finite element collection.
MPI_Comm GetComm() const
Definition: pfespace.hpp:273
double f(const Vector &xvec)
Definition: lor_mms.hpp:32
double GetElementSize(ElementTransformation *T, int type=0)
Definition: mesh.cpp:76
virtual const FiniteElement * GetFE(int i) const
Definition: pfespace.cpp:535
static void Sum(OpData< T >)
Reduce operation Sum, instantiated for int and double.
virtual void SetFromTrueDofs(const Vector &tv)
Set the GridFunction from the given true-dof vector.
Definition: pgridfunc.hpp:169
int Height() const
Get the height (size of output) of the Operator. Synonym with NumRows().
Definition: operator.hpp:67
double ComputeCFL(ParGridFunction &u, double dt)
Compute CFL.
Class for boundary integration .
Definition: lininteg.hpp:210
virtual double * Write(bool on_dev=true)
Shortcut for mfem::Write(vec.GetMemory(), vec.Size(), on_dev).
Definition: vector.hpp:457
Jacobi smoothing for a given bilinear form (no matrix necessary).
Definition: solvers.hpp:274
void Assemble()
Assembles the ParLinearForm i.e. sums over all domain/bdr integrators.
Definition: plinearform.cpp:46
double b
Definition: lissajous.cpp:42
void SetPrintLevel(int print_level)
Definition: hypre.hpp:1659
const DenseMatrix & Jacobian()
Return the Jacobian matrix of the transformation at the currently set IntegrationPoint, using the method SetIntPoint().
Definition: eltrans.hpp:119
void SetMaxIter(int max_it)
Definition: solvers.hpp:200
virtual void Mult(const HypreParVector &b, HypreParVector &x) const
Solve the linear system Ax=b.
Definition: hypre.cpp:3906
ConstantCoefficient onecoeff
T Max() const
Find the maximal element in the array, using the comparison operator &lt; for class T.
Definition: array.cpp:68
void Assemble(int skip_zeros=1)
void Assemble(int skip_zeros=1)
Assemble the local matrix.
virtual void Mult(const Vector &x, Vector &y) const
Evaluate the action of the NonlinearForm.
const IntegrationRule & GetNodes() const
Get a const reference to the nodes of the element.
Definition: fe_base.hpp:390
std::vector< double > dthist
void SetTimeIntegrationCoefficients(int step)
Update the EXTk/BDF time integration coefficient.
virtual void FormRectangularSystemMatrix(const Array< int > &trial_tdof_list, const Array< int > &test_tdof_list, OperatorHandle &A)
Return in A a parallel (on truedofs) version of this operator.
void AddVelDirichletBC(VectorCoefficient *coeff, Array< int > &attr)
Add a Dirichlet boundary condition to the velocity field.
Parallel smoothers in hypre.
Definition: hypre.hpp:962
void GetValues(int i, const IntegrationRule &ir, Vector &vals, int vdim=1) const
Definition: gridfunc.cpp:513
int GetVDim() const
Returns vector dimension.
Definition: fespace.hpp:581
virtual const Operator * GetRestriction() const
Restriction operator from input vectors for the operator to linear algebra (linear system) vectors...
Definition: operator.hpp:119
FiniteElementSpace * FESpace()
Definition: gridfunc.hpp:668
void AddBoundaryIntegrator(LinearFormIntegrator *lfi)
Adds new Boundary Integrator. Assumes ownership of lfi.
Definition: linearform.cpp:72
int Dimension() const
Definition: mesh.hpp:1006
void ComputeCurl2D(ParGridFunction &u, ParGridFunction &cu, bool assume_scalar=false)
Compute for .
void SetAssemblyLevel(AssemblyLevel assembly_level)
Set the desired assembly level. The default is AssemblyLevel::LEGACY.
void Start()
Start the stopwatch. The elapsed time is not cleared.
Definition: tic_toc.cpp:411
void SetSubVectorComplement(const Array< int > &dofs, const double val)
Set all vector entries NOT in the dofs Array to the given val.
Definition: vector.cpp:708
A general vector function coefficient.
void SetAssemblyLevel(AssemblyLevel assembly_level)
Set the desired assembly level. The default is AssemblyLevel::LEGACY.
ElementTransformation * GetElementTransformation(int i) const
Returns ElementTransformation for the i-th element.
Definition: fespace.hpp:647
void AddDomainIntegrator(LinearFormIntegrator *lfi)
Adds new Domain Integrator. Assumes ownership of lfi.
Definition: linearform.cpp:41
Array< int > bdr_attributes
A list of all unique boundary attributes used by the Mesh.
Definition: mesh.hpp:270
int GetMyRank() const
Definition: pmesh.hpp:353
void SetRelTol(double rtol)
Definition: solvers.hpp:198
ParFiniteElementSpace * vfes
Velocity finite element space.
ComplexDenseMatrix * MultAtB(const ComplexDenseMatrix &A, const ComplexDenseMatrix &B)
Multiply the complex conjugate transpose of a matrix A with a matrix B. A^H*B.
MPI_Comm GetComm() const
Definition: pmesh.hpp:351
void ProjectBdrCoefficient(Coefficient *coeff[], VectorCoefficient *vcoeff, Array< int > &attr)
Definition: pgridfunc.cpp:633
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Definition: fespace.hpp:96
int GetDof() const
Returns the number of degrees of freedom in the finite element.
Definition: fe_base.hpp:323
ParLinearForm * mass_lf
Linear form to compute the mass matrix in various subroutines.
Base class Coefficients that optionally depend on space and time. These are used by the BilinearFormI...
Definition: coefficient.hpp:41
void ProjectGridFunction(const GridFunction &src)
Project the src GridFunction to this GridFunction, both of which must be on the same mesh...
Definition: gridfunc.cpp:1855
void InitTVectors(const Operator *Po, const Operator *Ri, const Operator *Pi, Vector &x, Vector &b, Vector &X, Vector &B) const
Initializes memory for true vectors of linear system.
Definition: operator.cpp:22
void SetSize(int nsize)
Change the logical size of the array, keep existing entries.
Definition: array.hpp:679
GroupCommunicator & GroupComm()
Return a reference to the internal GroupCommunicator (on VDofs)
Definition: pfespace.hpp:338
Vector & Set(const double a, const Vector &x)
(*this) = a * x
Definition: vector.cpp:267
std::vector< AccelTerm_T > accel_terms
virtual void RecoverFEMSolution(const Vector &X, const Vector &b, Vector &x)
bool verbose
Enable/disable verbose output.
virtual void FormSystemMatrix(const Array< int > &ess_tdof_list, OperatorHandle &A)
Form the linear system matrix A, see FormLinearSystem() for details.
Solver wrapper which orthogonalizes the input and output vector.
Definition: solvers.hpp:1135
virtual double * ReadWrite(bool on_dev=true)
Shortcut for mfem::ReadWrite(vec.GetMemory(), vec.Size(), on_dev).
Definition: vector.hpp:465
Vector & Add(const double a, const Vector &Va)
(*this) += a * Va
Definition: vector.cpp:252
double kin_vis
Kinematic viscosity (dimensionless).
Class for integration point with weight.
Definition: intrules.hpp:25
Class for parallel bilinear form using different test and trial FE spaces.
std::vector< VelDirichletBC_T > vel_dbcs
void EliminateRHS(const Vector &x, Vector &b) const
Eliminate &quot;essential boundary condition&quot; values specified in x from the given right-hand side b...
Definition: operator.cpp:457
virtual const SparseMatrix * GetRestrictionMatrix() const
Get the R matrix which restricts a local dof vector to true dof vector.
Definition: pfespace.hpp:389
VectorGridFunctionCoefficient * FText_gfcoeff
void UpdateTimestepHistory(double dt)
Rotate entries in the time step and solution history arrays.
ParFiniteElementSpace * pfes
Pressure finite element space.
int dim
Definition: ex24.cpp:53
void MultTranspose(const Vector &x, Vector &y) const
Multiply a vector with the transposed matrix. y = At * x.
Definition: sparsemat.cpp:904
void AddDomainIntegrator(BilinearFormIntegrator *bfi)
Adds new Domain Integrator. Assumes ownership of bfi.
virtual const FiniteElement * GetFE(int i) const
Returns pointer to the FiniteElement in the FiniteElementCollection associated with i&#39;th element in t...
Definition: fespace.cpp:2781
HypreParVector * GetTrueDofs() const
Returns the true dofs in a new HypreParVector.
Definition: pgridfunc.cpp:143
Class for parallel bilinear form.
ParMixedBilinearForm * D_form
virtual void SetOperator(const Operator &op)
Also calls SetOperator for the preconditioner if there is one.
Definition: solvers.hpp:479
A general function coefficient.
Vector data type.
Definition: vector.hpp:60
void SetSolver(Solver &s)
Set the solver used by the OrthoSolver.
Definition: solvers.cpp:3424
Vector coefficient defined by a vector GridFunction.
virtual void SetPreconditioner(Solver &pr)
This should be called before SetOperator.
Definition: solvers.cpp:173
Arbitrary order H1-conforming (continuous) finite elements.
Definition: fe_coll.hpp:220
void Step(double &time, double dt, int cur_step, bool provisional=false)
Compute solution at the next time step t+dt.
ParFiniteElementSpace * vfes_filter
void AddPresDirichletBC(Coefficient *coeff, Array< int > &attr)
Add a Dirichlet boundary condition to the pressure field.
double u(const Vector &xvec)
Definition: lor_mms.hpp:24
Class for parallel grid function.
Definition: pgridfunc.hpp:32
void SetSize(int s)
Change the size of the DenseMatrix to s x s.
Definition: densemat.hpp:105
Square Operator for imposing essential boundary conditions using only the action, Mult()...
Definition: operator.hpp:845
OutStream out(std::cout)
Global stream used by the library for standard output. Initially it uses the same std::streambuf as s...
Definition: globals.hpp:66
Abstract operator.
Definition: operator.hpp:24
Wrapper for hypre&#39;s ParCSR matrix class.
Definition: hypre.hpp:343
virtual void CalcDShape(const IntegrationPoint &ip, DenseMatrix &dshape) const =0
Evaluate the gradients of all shape functions of a scalar finite element in reference space at the gi...
Class for parallel meshes.
Definition: pmesh.hpp:32
IntegrationRules IntRules(0, Quadrature1D::GaussLegendre)
A global object with all integration rules (defined in intrules.cpp)
Definition: intrules.hpp:379
double Sum() const
Return the sum of the vector entries.
Definition: vector.cpp:920
virtual const double * Read(bool on_dev=true) const
Shortcut for mfem::Read(vec.GetMemory(), vec.Size(), on_dev).
Definition: vector.hpp:449
void ParallelAssemble(Vector &tv)
Assemble the vector on the true dofs, i.e. P^t v.
Definition: plinearform.cpp:95
ParMesh * pmesh
The parallel mesh.
void AddAccelTerm(VectorCoefficient *coeff, Array< int > &attr)
Add an acceleration term to the RHS of the equation.
int order
The order of the velocity and pressure space.
ConstantCoefficient H_bdfcoeff
void Neg()
(*this) = -(*this)
Definition: vector.cpp:305
ParFiniteElementSpace * ParFESpace() const
Definition: pgridfunc.hpp:113