MFEM  v4.4.0
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 {
32  vfec = new H1_FECollection(order, pmesh->Dimension());
33  pfec = new H1_FECollection(order);
35  pfes = new ParFiniteElementSpace(pmesh, pfec);
36 
37  // Check if fully periodic mesh
38  if (!(pmesh->bdr_attributes.Size() == 0))
39  {
41  vel_ess_attr = 0;
42 
44  pres_ess_attr = 0;
45  }
46 
47  int vfes_truevsize = vfes->GetTrueVSize();
48  int pfes_truevsize = pfes->GetTrueVSize();
49 
50  un.SetSize(vfes_truevsize);
51  un = 0.0;
52  un_next.SetSize(vfes_truevsize);
53  un_next = 0.0;
54  unm1.SetSize(vfes_truevsize);
55  unm1 = 0.0;
56  unm2.SetSize(vfes_truevsize);
57  unm2 = 0.0;
58  fn.SetSize(vfes_truevsize);
59  Nun.SetSize(vfes_truevsize);
60  Nun = 0.0;
61  Nunm1.SetSize(vfes_truevsize);
62  Nunm1 = 0.0;
63  Nunm2.SetSize(vfes_truevsize);
64  Nunm2 = 0.0;
65  Fext.SetSize(vfes_truevsize);
66  FText.SetSize(vfes_truevsize);
67  Lext.SetSize(vfes_truevsize);
68  resu.SetSize(vfes_truevsize);
69 
70  tmp1.SetSize(vfes_truevsize);
71 
72  pn.SetSize(pfes_truevsize);
73  pn = 0.0;
74  resp.SetSize(pfes_truevsize);
75  resp = 0.0;
76  FText_bdr.SetSize(pfes_truevsize);
77  g_bdr.SetSize(pfes_truevsize);
78 
80  un_gf = 0.0;
82  un_next_gf = 0.0;
83 
89 
90  pn_gf.SetSpace(pfes);
91  pn_gf = 0.0;
92  resp_gf.SetSpace(pfes);
93 
94  cur_step = 0;
95 
96  PrintInfo();
97 }
98 
99 void NavierSolver::Setup(double dt)
100 {
101  if (verbose && pmesh->GetMyRank() == 0)
102  {
103  mfem::out << "Setup" << std::endl;
104  if (partial_assembly)
105  {
106  mfem::out << "Using Partial Assembly" << std::endl;
107  }
108  else
109  {
110  mfem::out << "Using Full Assembly" << std::endl;
111  }
112  }
113 
114  sw_setup.Start();
115 
116  pmesh_lor = new ParMesh(
118  pfec_lor = new H1_FECollection(1);
120 
123 
124  Array<int> empty;
125 
126  // GLL integration rule (Numerical Integration)
128  const IntegrationRule &ir_ni = rules_ni.Get(vfes->GetFE(0)->GetGeomType(),
129  2 * order - 1);
130 
131  nlcoeff.constant = -1.0;
132  N = new ParNonlinearForm(vfes);
134  if (partial_assembly)
135  {
137  N->Setup();
138  }
139 
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 
156  if (numerical_integ)
157  {
158  // blfi->SetIntRule(&ir_ni);
159  }
160  Sp_form->AddDomainIntegrator(sp_blfi);
161  if (partial_assembly)
162  {
164  }
165  Sp_form->Assemble();
167 
170  if (partial_assembly)
171  {
173  }
174  D_form->Assemble();
175  D_form->FormRectangularSystemMatrix(empty, empty, D);
176 
179  if (partial_assembly)
180  {
182  }
183  G_form->Assemble();
184  G_form->FormRectangularSystemMatrix(empty, empty, G);
185 
187  H_bdfcoeff.constant = 1.0 / dt;
188  H_form = new ParBilinearForm(vfes);
191  if (partial_assembly)
192  {
194  }
195  H_form->Assemble();
197 
201  *FText_gfcoeff),
202  vel_ess_attr);
203 
205  for (auto &vel_dbc : vel_dbcs)
206  {
208  *vel_dbc.coeff),
209  vel_dbc.attr);
210  }
211 
212  f_form = new ParLinearForm(vfes);
213  for (auto &accel_term : accel_terms)
214  {
215  auto *vdlfi = new VectorDomainLFIntegrator(*accel_term.coeff);
216  // @TODO: This order should always be the same as the nonlinear forms one!
217  // const IntegrationRule &ir = IntRules.Get(vfes->GetFE(0)->GetGeomType(),
218  // 4 * order);
219  // vdlfi->SetIntRule(&ir);
220  f_form->AddDomainIntegrator(vdlfi);
221  }
222 
223  if (partial_assembly)
224  {
225  Vector diag_pa(vfes->GetTrueVSize());
226  Mv_form->AssembleDiagonal(diag_pa);
227  MvInvPC = new OperatorJacobiSmoother(diag_pa, empty);
228  }
229  else
230  {
232  dynamic_cast<HypreSmoother *>(MvInvPC)->SetType(HypreSmoother::Jacobi, 1);
233  }
234  MvInv = new CGSolver(vfes->GetComm());
235  MvInv->iterative_mode = false;
236  MvInv->SetOperator(*Mv);
239  MvInv->SetRelTol(1e-12);
240  MvInv->SetMaxIter(200);
241 
242  if (partial_assembly)
243  {
251  SpInvPC->Mult(resp, pn);
254  }
255  else
256  {
261  }
262  SpInv = new CGSolver(vfes->GetComm());
263  SpInv->iterative_mode = true;
264  SpInv->SetOperator(*Sp);
265  if (pres_dbcs.empty())
266  {
268  }
269  else
270  {
272  }
275  SpInv->SetMaxIter(200);
276 
277  if (partial_assembly)
278  {
279  Vector diag_pa(vfes->GetTrueVSize());
280  H_form->AssembleDiagonal(diag_pa);
282  }
283  else
284  {
286  dynamic_cast<HypreSmoother *>(HInvPC)->SetType(HypreSmoother::Jacobi, 1);
287  }
288  HInv = new CGSolver(vfes->GetComm());
289  HInv->iterative_mode = true;
290  HInv->SetOperator(*H);
294  HInv->SetMaxIter(200);
295 
296  // If the initial condition was set, it has to be aligned with dependent
297  // Vectors and GridFunctions
299  un_next = un;
301 
302  // Set initial time step in the history array
303  dthist[0] = dt;
304 
305  if (filter_alpha != 0.0)
306  {
308  pmesh->Dimension());
310  vfec_filter,
311  pmesh->Dimension());
312 
313  un_NM1_gf.SetSpace(vfes_filter);
314  un_NM1_gf = 0.0;
315 
317  un_filtered_gf = 0.0;
318  }
319 
320  sw_setup.Stop();
321 }
322 
324 {
325  // Rotate values in time step history
326  dthist[2] = dthist[1];
327  dthist[1] = dthist[0];
328  dthist[0] = dt;
329 
330  // Rotate values in nonlinear extrapolation history
331  Nunm2 = Nunm1;
332  Nunm1 = Nun;
333 
334  // Rotate values in solution history
335  unm2 = unm1;
336  unm1 = un;
337 
338  // Update the current solution and corresponding GridFunction
340  un = un_next;
342 }
343 
344 void NavierSolver::Step(double &time, double dt, int current_step,
345  bool provisional)
346 {
347  sw_step.Start();
348 
349  SetTimeIntegrationCoefficients(current_step);
350 
351  // Set current time for velocity Dirichlet boundary conditions.
352  for (auto &vel_dbc : vel_dbcs)
353  {
354  vel_dbc.coeff->SetTime(time + dt);
355  }
356 
357  // Set current time for pressure Dirichlet boundary conditions.
358  for (auto &pres_dbc : pres_dbcs)
359  {
360  pres_dbc.coeff->SetTime(time + dt);
361  }
362 
363  H_bdfcoeff.constant = bd0 / dt;
364  H_form->Update();
365  H_form->Assemble();
367 
368  HInv->SetOperator(*H);
369  if (partial_assembly)
370  {
371  delete HInvPC;
372  Vector diag_pa(vfes->GetTrueVSize());
373  H_form->AssembleDiagonal(diag_pa);
376  }
377 
378  // Extrapolated f^{n+1}.
379  for (auto &accel_term : accel_terms)
380  {
381  accel_term.coeff->SetTime(time + dt);
382  }
383 
384  f_form->Assemble();
386 
387  // Nonlinear extrapolated terms.
388  sw_extrap.Start();
389 
390  N->Mult(un, Nun);
391  Nun.Add(1.0, fn);
392 
393  {
394  const auto d_Nun = Nun.Read();
395  const auto d_Nunm1 = Nunm1.Read();
396  const auto d_Nunm2 = Nunm2.Read();
397  auto d_Fext = Fext.Write();
398  const auto ab1_ = ab1;
399  const auto ab2_ = ab2;
400  const auto ab3_ = ab3;
401  MFEM_FORALL(i, Fext.Size(),
402  d_Fext[i] = ab1_ * d_Nun[i] +
403  ab2_ * d_Nunm1[i] +
404  ab3_ * d_Nunm2[i];);
405  }
406 
407  // Fext = M^{-1} (F(u^{n}) + f^{n+1})
408  MvInv->Mult(Fext, tmp1);
411  Fext.Set(1.0, tmp1);
412 
413  // Compute BDF terms.
414  {
415  const double bd1idt = -bd1 / dt;
416  const double bd2idt = -bd2 / dt;
417  const double bd3idt = -bd3 / dt;
418  const auto d_un = un.Read();
419  const auto d_unm1 = unm1.Read();
420  const auto d_unm2 = unm2.Read();
421  auto d_Fext = Fext.ReadWrite();
422  MFEM_FORALL(i, Fext.Size(),
423  d_Fext[i] += bd1idt * d_un[i] +
424  bd2idt * d_unm1[i] +
425  bd3idt * d_unm2[i];);
426  }
427 
428  sw_extrap.Stop();
429 
430  // Pressure Poisson.
431  sw_curlcurl.Start();
432  {
433  const auto d_un = un.Read();
434  const auto d_unm1 = unm1.Read();
435  const auto d_unm2 = unm2.Read();
436  auto d_Lext = Lext.Write();
437  const auto ab1_ = ab1;
438  const auto ab2_ = ab2;
439  const auto ab3_ = ab3;
440  MFEM_FORALL(i, Lext.Size(),
441  d_Lext[i] = ab1_ * d_un[i] +
442  ab2_ * d_unm1[i] +
443  ab3_ * d_unm2[i];);
444  }
445 
447  if (pmesh->Dimension() == 2)
448  {
451  }
452  else
453  {
456  }
457 
459  Lext *= kin_vis;
460 
461  sw_curlcurl.Stop();
462 
463  // \tilde{F} = F - \nu CurlCurl(u)
464  FText.Set(-1.0, Lext);
465  FText.Add(1.0, Fext);
466 
467  // p_r = \nabla \cdot FText
468  D->Mult(FText, resp);
469  resp.Neg();
470 
471  // Add boundary terms.
475 
476  g_bdr_form->Assemble();
478  resp.Add(1.0, FText_bdr);
479  resp.Add(-bd0 / dt, g_bdr);
480 
481  if (pres_dbcs.empty())
482  {
484  }
485 
486  for (auto &pres_dbc : pres_dbcs)
487  {
488  pn_gf.ProjectBdrCoefficient(*pres_dbc.coeff, pres_dbc.attr);
489  }
490 
492 
493  Vector X1, B1;
494  if (partial_assembly)
495  {
496  auto *SpC = Sp.As<ConstrainedOperator>();
497  EliminateRHS(*Sp_form, *SpC, pres_ess_tdof, pn_gf, resp_gf, X1, B1, 1);
498  }
499  else
500  {
502  }
503  sw_spsolve.Start();
504  SpInv->Mult(B1, X1);
505  sw_spsolve.Stop();
509 
510  // If the boundary conditions on the pressure are pure Neumann remove the
511  // nullspace by removing the mean of the pressure solution. This is also
512  // ensured by the OrthoSolver wrapper for the preconditioner which removes
513  // the nullspace after every application.
514  if (pres_dbcs.empty())
515  {
516  MeanZero(pn_gf);
517  }
518 
520 
521  // Project velocity.
522  G->Mult(pn, resu);
523  resu.Neg();
524  Mv->Mult(Fext, tmp1);
525  resu.Add(1.0, tmp1);
526 
527  // un_next_gf = un_gf;
528 
529  for (auto &vel_dbc : vel_dbcs)
530  {
532  }
533 
535 
536  Vector X2, B2;
537  if (partial_assembly)
538  {
539  auto *HC = H.As<ConstrainedOperator>();
540  EliminateRHS(*H_form, *HC, vel_ess_tdof, un_next_gf, resu_gf, X2, B2, 1);
541  }
542  else
543  {
545  }
546  sw_hsolve.Start();
547  HInv->Mult(B2, X2);
548  sw_hsolve.Stop();
552 
554 
555  // If the current time step is not provisional, accept the computed solution
556  // and update the time step history by default.
557  if (!provisional)
558  {
560  time += dt;
561  }
562 
563  if (filter_alpha != 0.0)
564  {
567  const auto d_un_filtered_gf = un_filtered_gf.Read();
568  auto d_un_gf = un_gf.ReadWrite();
569  const auto filter_alpha_ = filter_alpha;
570  MFEM_FORALL(i,
571  un_gf.Size(),
572  d_un_gf[i] = (1.0 - filter_alpha_) * d_un_gf[i]
573  + filter_alpha_ * d_un_filtered_gf[i];);
574  }
575 
576  sw_step.Stop();
577 
578  if (verbose && pmesh->GetMyRank() == 0)
579  {
580  mfem::out << std::setw(7) << "" << std::setw(3) << "It" << std::setw(8)
581  << "Resid" << std::setw(12) << "Reltol"
582  << "\n";
583  // If numerical integration is active, there is no solve (thus no
584  // iterations), on the inverse velocity mass application.
585  if (!numerical_integ)
586  {
587  mfem::out << std::setw(5) << "MVIN " << std::setw(5) << std::fixed
588  << iter_mvsolve << " " << std::setw(3)
589  << std::setprecision(2) << std::scientific << res_mvsolve
590  << " " << 1e-12 << "\n";
591  }
592  mfem::out << std::setw(5) << "PRES " << std::setw(5) << std::fixed
593  << iter_spsolve << " " << std::setw(3) << std::setprecision(2)
594  << std::scientific << res_spsolve << " " << rtol_spsolve
595  << "\n";
596  mfem::out << std::setw(5) << "HELM " << std::setw(5) << std::fixed
597  << iter_hsolve << " " << std::setw(3) << std::setprecision(2)
598  << std::scientific << res_hsolve << " " << rtol_hsolve
599  << "\n";
600  mfem::out << std::setprecision(8);
601  mfem::out << std::fixed;
602  }
603 }
604 
606 {
607  // Make sure not to recompute the inner product linear form every
608  // application.
609  if (mass_lf == nullptr)
610  {
611  onecoeff.constant = 1.0;
612  mass_lf = new ParLinearForm(v.ParFESpace());
614  mass_lf->Assemble();
615 
616  ParGridFunction one_gf(v.ParFESpace());
618 
619  volume = mass_lf->operator()(one_gf);
620  }
621 
622  double integ = mass_lf->operator()(v);
623 
624  v -= integ / volume;
625 }
626 
628  ConstrainedOperator &constrainedA,
629  const Array<int> &ess_tdof_list,
630  Vector &x,
631  Vector &b,
632  Vector &X,
633  Vector &B,
634  int copy_interior)
635 {
636  const Operator *Po = A.GetOutputProlongation();
637  const Operator *Pi = A.GetProlongation();
638  const Operator *Ri = A.GetRestriction();
639  A.InitTVectors(Po, Ri, Pi, x, b, X, B);
640  if (!copy_interior)
641  {
642  X.SetSubVectorComplement(ess_tdof_list, 0.0);
643  }
644  constrainedA.EliminateRHS(X, B);
645 }
646 
648 {
649  double loc_sum = v.Sum();
650  double global_sum = 0.0;
651  int loc_size = v.Size();
652  int global_size = 0;
653 
654  MPI_Allreduce(&loc_sum, &global_sum, 1, MPI_DOUBLE, MPI_SUM, pfes->GetComm());
655  MPI_Allreduce(&loc_size, &global_size, 1, MPI_INT, MPI_SUM, pfes->GetComm());
656 
657  v -= global_sum / static_cast<double>(global_size);
658 }
659 
661 {
662  FiniteElementSpace *fes = u.FESpace();
663 
664  // AccumulateAndCountZones.
665  Array<int> zones_per_vdof;
666  zones_per_vdof.SetSize(fes->GetVSize());
667  zones_per_vdof = 0;
668 
669  cu = 0.0;
670 
671  // Local interpolation.
672  int elndofs;
673  Array<int> vdofs;
674  Vector vals;
675  Vector loc_data;
676  int vdim = fes->GetVDim();
677  DenseMatrix grad_hat;
678  DenseMatrix dshape;
679  DenseMatrix grad;
680  Vector curl;
681 
682  for (int e = 0; e < fes->GetNE(); ++e)
683  {
684  fes->GetElementVDofs(e, vdofs);
685  u.GetSubVector(vdofs, loc_data);
686  vals.SetSize(vdofs.Size());
688  const FiniteElement *el = fes->GetFE(e);
689  elndofs = el->GetDof();
690  int dim = el->GetDim();
691  dshape.SetSize(elndofs, dim);
692 
693  for (int dof = 0; dof < elndofs; ++dof)
694  {
695  // Project.
696  const IntegrationPoint &ip = el->GetNodes().IntPoint(dof);
697  tr->SetIntPoint(&ip);
698 
699  // Eval and GetVectorGradientHat.
700  el->CalcDShape(tr->GetIntPoint(), dshape);
701  grad_hat.SetSize(vdim, dim);
702  DenseMatrix loc_data_mat(loc_data.GetData(), elndofs, vdim);
703  MultAtB(loc_data_mat, dshape, grad_hat);
704 
705  const DenseMatrix &Jinv = tr->InverseJacobian();
706  grad.SetSize(grad_hat.Height(), Jinv.Width());
707  Mult(grad_hat, Jinv, grad);
708 
709  curl.SetSize(3);
710  curl(0) = grad(2, 1) - grad(1, 2);
711  curl(1) = grad(0, 2) - grad(2, 0);
712  curl(2) = grad(1, 0) - grad(0, 1);
713 
714  for (int j = 0; j < curl.Size(); ++j)
715  {
716  vals(elndofs * j + dof) = curl(j);
717  }
718  }
719 
720  // Accumulate values in all dofs, count the zones.
721  for (int j = 0; j < vdofs.Size(); j++)
722  {
723  int ldof = vdofs[j];
724  cu(ldof) += vals[j];
725  zones_per_vdof[ldof]++;
726  }
727  }
728 
729  // Communication
730 
731  // Count the zones globally.
732  GroupCommunicator &gcomm = u.ParFESpace()->GroupComm();
733  gcomm.Reduce<int>(zones_per_vdof, GroupCommunicator::Sum);
734  gcomm.Bcast(zones_per_vdof);
735 
736  // Accumulate for all vdofs.
737  gcomm.Reduce<double>(cu.GetData(), GroupCommunicator::Sum);
738  gcomm.Bcast<double>(cu.GetData());
739 
740  // Compute means.
741  for (int i = 0; i < cu.Size(); i++)
742  {
743  const int nz = zones_per_vdof[i];
744  if (nz)
745  {
746  cu(i) /= nz;
747  }
748  }
749 }
750 
752  ParGridFunction &cu,
753  bool assume_scalar)
754 {
755  FiniteElementSpace *fes = u.FESpace();
756 
757  // AccumulateAndCountZones.
758  Array<int> zones_per_vdof;
759  zones_per_vdof.SetSize(fes->GetVSize());
760  zones_per_vdof = 0;
761 
762  cu = 0.0;
763 
764  // Local interpolation.
765  int elndofs;
766  Array<int> vdofs;
767  Vector vals;
768  Vector loc_data;
769  int vdim = fes->GetVDim();
770  DenseMatrix grad_hat;
771  DenseMatrix dshape;
772  DenseMatrix grad;
773  Vector curl;
774 
775  for (int e = 0; e < fes->GetNE(); ++e)
776  {
777  fes->GetElementVDofs(e, vdofs);
778  u.GetSubVector(vdofs, loc_data);
779  vals.SetSize(vdofs.Size());
781  const FiniteElement *el = fes->GetFE(e);
782  elndofs = el->GetDof();
783  int dim = el->GetDim();
784  dshape.SetSize(elndofs, dim);
785 
786  for (int dof = 0; dof < elndofs; ++dof)
787  {
788  // Project.
789  const IntegrationPoint &ip = el->GetNodes().IntPoint(dof);
790  tr->SetIntPoint(&ip);
791 
792  // Eval and GetVectorGradientHat.
793  el->CalcDShape(tr->GetIntPoint(), dshape);
794  grad_hat.SetSize(vdim, dim);
795  DenseMatrix loc_data_mat(loc_data.GetData(), elndofs, vdim);
796  MultAtB(loc_data_mat, dshape, grad_hat);
797 
798  const DenseMatrix &Jinv = tr->InverseJacobian();
799  grad.SetSize(grad_hat.Height(), Jinv.Width());
800  Mult(grad_hat, Jinv, grad);
801 
802  if (assume_scalar)
803  {
804  curl.SetSize(2);
805  curl(0) = grad(0, 1);
806  curl(1) = -grad(0, 0);
807  }
808  else
809  {
810  curl.SetSize(2);
811  curl(0) = grad(1, 0) - grad(0, 1);
812  curl(1) = 0.0;
813  }
814 
815  for (int j = 0; j < curl.Size(); ++j)
816  {
817  vals(elndofs * j + dof) = curl(j);
818  }
819  }
820 
821  // Accumulate values in all dofs, count the zones.
822  for (int j = 0; j < vdofs.Size(); j++)
823  {
824  int ldof = vdofs[j];
825  cu(ldof) += vals[j];
826  zones_per_vdof[ldof]++;
827  }
828  }
829 
830  // Communication.
831 
832  // Count the zones globally.
833  GroupCommunicator &gcomm = u.ParFESpace()->GroupComm();
834  gcomm.Reduce<int>(zones_per_vdof, GroupCommunicator::Sum);
835  gcomm.Bcast(zones_per_vdof);
836 
837  // Accumulate for all vdofs.
838  gcomm.Reduce<double>(cu.GetData(), GroupCommunicator::Sum);
839  gcomm.Bcast<double>(cu.GetData());
840 
841  // Compute means.
842  for (int i = 0; i < cu.Size(); i++)
843  {
844  const int nz = zones_per_vdof[i];
845  if (nz)
846  {
847  cu(i) /= nz;
848  }
849  }
850 }
851 
853 {
854  ParMesh *pmesh_u = u.ParFESpace()->GetParMesh();
855  FiniteElementSpace *fes = u.FESpace();
856  int vdim = fes->GetVDim();
857 
858  Vector ux, uy, uz;
859  Vector ur, us, ut;
860  double cflx = 0.0;
861  double cfly = 0.0;
862  double cflz = 0.0;
863  double cflm = 0.0;
864  double cflmax = 0.0;
865 
866  for (int e = 0; e < fes->GetNE(); ++e)
867  {
868  const FiniteElement *fe = fes->GetFE(e);
869  const IntegrationRule &ir = IntRules.Get(fe->GetGeomType(),
870  fe->GetOrder());
872 
873  u.GetValues(e, ir, ux, 1);
874  ur.SetSize(ux.Size());
875  u.GetValues(e, ir, uy, 2);
876  us.SetSize(uy.Size());
877  if (vdim == 3)
878  {
879  u.GetValues(e, ir, uz, 3);
880  ut.SetSize(uz.Size());
881  }
882 
883  double hmin = pmesh_u->GetElementSize(e, 1) /
884  (double) fes->GetElementOrder(0);
885 
886  for (int i = 0; i < ir.GetNPoints(); ++i)
887  {
888  const IntegrationPoint &ip = ir.IntPoint(i);
889  tr->SetIntPoint(&ip);
890  const DenseMatrix &invJ = tr->InverseJacobian();
891  const double detJinv = 1.0 / tr->Jacobian().Det();
892 
893  if (vdim == 2)
894  {
895  ur(i) = (ux(i) * invJ(0, 0) + uy(i) * invJ(1, 0)) * detJinv;
896  us(i) = (ux(i) * invJ(0, 1) + uy(i) * invJ(1, 1)) * detJinv;
897  }
898  else if (vdim == 3)
899  {
900  ur(i) = (ux(i) * invJ(0, 0) + uy(i) * invJ(1, 0)
901  + uz(i) * invJ(2, 0))
902  * detJinv;
903  us(i) = (ux(i) * invJ(0, 1) + uy(i) * invJ(1, 1)
904  + uz(i) * invJ(2, 1))
905  * detJinv;
906  ut(i) = (ux(i) * invJ(0, 2) + uy(i) * invJ(1, 2)
907  + uz(i) * invJ(2, 2))
908  * detJinv;
909  }
910 
911  cflx = fabs(dt * ux(i) / hmin);
912  cfly = fabs(dt * uy(i) / hmin);
913  if (vdim == 3)
914  {
915  cflz = fabs(dt * uz(i) / hmin);
916  }
917  cflm = cflx + cfly + cflz;
918  cflmax = fmax(cflmax, cflm);
919  }
920  }
921 
922  double cflmax_global = 0.0;
923  MPI_Allreduce(&cflmax,
924  &cflmax_global,
925  1,
926  MPI_DOUBLE,
927  MPI_MAX,
928  pmesh_u->GetComm());
929 
930  return cflmax_global;
931 }
932 
934 {
935  vel_dbcs.emplace_back(attr, coeff);
936 
937  if (verbose && pmesh->GetMyRank() == 0)
938  {
939  mfem::out << "Adding Velocity Dirichlet BC to attributes ";
940  for (int i = 0; i < attr.Size(); ++i)
941  {
942  if (attr[i] == 1)
943  {
944  mfem::out << i << " ";
945  }
946  }
947  mfem::out << std::endl;
948  }
949 
950  for (int i = 0; i < attr.Size(); ++i)
951  {
952  MFEM_ASSERT((vel_ess_attr[i] && attr[i]) == 0,
953  "Duplicate boundary definition deteceted.");
954  if (attr[i] == 1)
955  {
956  vel_ess_attr[i] = 1;
957  }
958  }
959 }
960 
962 {
964 }
965 
967 {
968  pres_dbcs.emplace_back(attr, coeff);
969 
970  if (verbose && pmesh->GetMyRank() == 0)
971  {
972  mfem::out << "Adding Pressure 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((pres_ess_attr[i] && attr[i]) == 0,
986  "Duplicate boundary definition deteceted.");
987  if (attr[i] == 1)
988  {
989  pres_ess_attr[i] = 1;
990  }
991  }
992 }
993 
995 {
997 }
998 
1000 {
1001  accel_terms.emplace_back(attr, coeff);
1002 
1003  if (verbose && pmesh->GetMyRank() == 0)
1004  {
1005  mfem::out << "Adding Acceleration term 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 
1018 {
1020 }
1021 
1023 {
1024  // Maximum BDF order to use at current time step
1025  // step + 1 <= order <= max_bdf_order
1026  int bdf_order = std::min(step + 1, max_bdf_order);
1027 
1028  // Ratio of time step history at dt(t_{n}) - dt(t_{n-1})
1029  double rho1 = 0.0;
1030 
1031  // Ratio of time step history at dt(t_{n-1}) - dt(t_{n-2})
1032  double rho2 = 0.0;
1033 
1034  rho1 = dthist[0] / dthist[1];
1035 
1036  if (bdf_order == 3)
1037  {
1038  rho2 = dthist[1] / dthist[2];
1039  }
1040 
1041  if (step == 0 && bdf_order == 1)
1042  {
1043  bd0 = 1.0;
1044  bd1 = -1.0;
1045  bd2 = 0.0;
1046  bd3 = 0.0;
1047  ab1 = 1.0;
1048  ab2 = 0.0;
1049  ab3 = 0.0;
1050  }
1051  else if (step >= 1 && bdf_order == 2)
1052  {
1053  bd0 = (1.0 + 2.0 * rho1) / (1.0 + rho1);
1054  bd1 = -(1.0 + rho1);
1055  bd2 = pow(rho1, 2.0) / (1.0 + rho1);
1056  bd3 = 0.0;
1057  ab1 = 1.0 + rho1;
1058  ab2 = -rho1;
1059  ab3 = 0.0;
1060  }
1061  else if (step >= 2 && bdf_order == 3)
1062  {
1063  bd0 = 1.0 + rho1 / (1.0 + rho1)
1064  + (rho2 * rho1) / (1.0 + rho2 * (1 + rho1));
1065  bd1 = -1.0 - rho1 - (rho2 * rho1 * (1.0 + rho1)) / (1.0 + rho2);
1066  bd2 = pow(rho1, 2.0) * (rho2 + 1.0 / (1.0 + rho1));
1067  bd3 = -(pow(rho2, 3.0) * pow(rho1, 2.0) * (1.0 + rho1))
1068  / ((1.0 + rho2) * (1.0 + rho2 + rho2 * rho1));
1069  ab1 = ((1.0 + rho1) * (1.0 + rho2 * (1.0 + rho1))) / (1.0 + rho2);
1070  ab2 = -rho1 * (1.0 + rho2 * (1.0 + rho1));
1071  ab3 = (pow(rho2, 2.0) * rho1 * (1.0 + rho1)) / (1.0 + rho2);
1072  }
1073 }
1074 
1076 {
1077  double my_rt[6], rt_max[6];
1078 
1079  my_rt[0] = sw_setup.RealTime();
1080  my_rt[1] = sw_step.RealTime();
1081  my_rt[2] = sw_extrap.RealTime();
1082  my_rt[3] = sw_curlcurl.RealTime();
1083  my_rt[4] = sw_spsolve.RealTime();
1084  my_rt[5] = sw_hsolve.RealTime();
1085 
1086  MPI_Reduce(my_rt, rt_max, 6, MPI_DOUBLE, MPI_MAX, 0, pmesh->GetComm());
1087 
1088  if (pmesh->GetMyRank() == 0)
1089  {
1090  mfem::out << std::setw(10) << "SETUP" << std::setw(10) << "STEP"
1091  << std::setw(10) << "EXTRAP" << std::setw(10) << "CURLCURL"
1092  << std::setw(10) << "PSOLVE" << std::setw(10) << "HSOLVE"
1093  << "\n";
1094 
1095  mfem::out << std::setprecision(3) << std::setw(10) << my_rt[0]
1096  << std::setw(10) << my_rt[1] << std::setw(10) << my_rt[2]
1097  << std::setw(10) << my_rt[3] << std::setw(10) << my_rt[4]
1098  << std::setw(10) << my_rt[5] << "\n";
1099 
1100  mfem::out << std::setprecision(3) << std::setw(10) << " " << std::setw(10)
1101  << my_rt[1] / my_rt[1] << std::setw(10) << my_rt[2] / my_rt[1]
1102  << std::setw(10) << my_rt[3] / my_rt[1] << std::setw(10)
1103  << my_rt[4] / my_rt[1] << std::setw(10) << my_rt[5] / my_rt[1]
1104  << "\n";
1105 
1106  mfem::out << std::setprecision(8);
1107  }
1108 }
1109 
1111 {
1112  int fes_size0 = vfes->GlobalVSize();
1113  int fes_size1 = pfes->GlobalVSize();
1114 
1115  if (pmesh->GetMyRank() == 0)
1116  {
1117  mfem::out << "NAVIER version: " << NAVIER_VERSION << std::endl
1118  << "MFEM version: " << MFEM_VERSION << std::endl
1119  << "MFEM GIT: " << MFEM_GIT_STRING << std::endl
1120  << "Velocity #DOFs: " << fes_size0 << std::endl
1121  << "Pressure #DOFs: " << fes_size1 << std::endl;
1122  }
1123 }
1124 
1126 {
1127  delete FText_gfcoeff;
1128  delete g_bdr_form;
1129  delete FText_bdr_form;
1130  delete mass_lf;
1131  delete Mv_form;
1132  delete N;
1133  delete Sp_form;
1134  delete D_form;
1135  delete G_form;
1136  delete HInvPC;
1137  delete HInv;
1138  delete H_form;
1139  delete SpInv;
1140  delete MvInvPC;
1141  delete Sp_form_lor;
1142  delete SpInvOrthoPC;
1143  delete SpInvPC;
1144  delete f_form;
1145  delete pfes_lor;
1146  delete pfec_lor;
1147  delete pmesh_lor;
1148  delete MvInv;
1149  delete vfec;
1150  delete pfec;
1151  delete vfes;
1152  delete pfes;
1153  delete vfec_filter;
1154  delete vfes_filter;
1155 }
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:97
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:563
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
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
static ParMesh MakeRefined(ParMesh &orig_mesh, int ref_factor, int ref_type)
Create a uniformly refined (by any factor) version of orig_mesh.
Definition: pmesh.cpp:1366
Array< BilinearFormIntegrator * > * GetDBFI()
Access all the integrators added with AddDomainIntegrator().
FiniteElementCollection * pfec_lor
void SetSize(int s)
Resize the vector to size s.
Definition: vector.hpp:521
virtual void Setup()
Setup the NonlinearForm: based on the current AssemblyLevel and the current mesh, optionally...
double Det() const
Definition: densemat.cpp:436
ConstantCoefficient nlcoeff
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
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:534
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
Container class for integration rules.
Definition: intrules.hpp:311
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:199
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:121
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:655
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:208
NavierSolver(ParMesh *mesh, int order, double kin_vis)
Initialize data structures, set FE space order and kinematic viscosity.
ParMixedBilinearForm * G_form
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:115
int GetNE() const
Returns number of elements in the mesh.
Definition: fespace.hpp:590
The BoomerAMG solver in hypre.
Definition: hypre.hpp:1443
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
ParFiniteElementSpace * pfes_lor
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
virtual void SetOperator(const Operator &op)
Set/update the solver for the given operator.
int Height() const
Get the height (size of output) of the Operator. Synonym with NumRows().
Definition: operator.hpp:66
double ComputeCFL(ParGridFunction &u, double dt)
Compute CFL.
Class for boundary integration .
Definition: lininteg.hpp:178
virtual double * Write(bool on_dev=true)
Shortcut for mfem::Write(vec.GetMemory(), vec.Size(), on_dev).
Definition: vector.hpp:446
Jacobi smoothing for a given bilinear form (no matrix necessary).
Definition: solvers.hpp:274
double b
Definition: lissajous.cpp:42
void SetPrintLevel(int print_level)
Definition: hypre.hpp:1523
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:3618
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.
ParBilinearForm * Sp_form_lor
void AddVelDirichletBC(VectorCoefficient *coeff, Array< int > &attr)
Add a Dirichlet boundary condition to the velocity field.
Parallel smoothers in hypre.
Definition: hypre.hpp:896
void GetValues(int i, const IntegrationRule &ir, Vector &vals, int vdim=1) const
Definition: gridfunc.cpp:510
int GetVDim() const
Returns vector dimension.
Definition: fespace.hpp:557
virtual const Operator * GetRestriction() const
Restriction operator from input vectors for the operator to linear algebra (linear system) vectors...
Definition: operator.hpp:118
FiniteElementSpace * FESpace()
Definition: gridfunc.hpp:652
void AddBoundaryIntegrator(LinearFormIntegrator *lfi)
Adds new Boundary Integrator. Assumes ownership of lfi.
Definition: linearform.cpp:70
int Dimension() const
Definition: mesh.hpp:999
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:695
FDualNumber< tbase > pow(const FDualNumber< tbase > &a, const FDualNumber< tbase > &b)
pow([dual number],[dual number])
Definition: fdual.hpp:543
A general vector function coefficient.
void SetAssemblyLevel(AssemblyLevel assembly_level)
Set the desired assembly level. The default is AssemblyLevel::LEGACY.
void UseExternalIntegrators()
Indicate that integrators are not owned by the BilinearForm.
ElementTransformation * GetElementTransformation(int i) const
Returns ElementTransformation for the i-th element.
Definition: fespace.hpp:623
void AddDomainIntegrator(LinearFormIntegrator *lfi)
Adds new Domain Integrator. Assumes ownership of lfi.
Definition: linearform.cpp:39
Abstract base class BilinearFormIntegrator.
Definition: bilininteg.hpp:34
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:290
void SetRelTol(double rtol)
Definition: solvers.hpp:198
ParFiniteElementSpace * vfes
Velocity finite element space.
MPI_Comm GetComm() const
Definition: pmesh.hpp:288
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:88
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:39
void ProjectGridFunction(const GridFunction &src)
Project the src GridFunction to this GridFunction, both of which must be on the same mesh...
Definition: gridfunc.cpp:1851
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.
virtual double * ReadWrite(bool on_dev=true)
Shortcut for mfem::ReadWrite(vec.GetMemory(), vec.Size(), on_dev).
Definition: vector.hpp:454
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:813
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:2783
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.
void MultAtB(const DenseMatrix &A, const DenseMatrix &B, DenseMatrix &AtB)
Multiply the transpose of a matrix A with a matrix B: At*B.
Definition: densemat.cpp:2633
Vector data type.
Definition: vector.hpp:60
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:216
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:841
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:337
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
Solver wrapper which orthogonalizes the input and output vector.
double Sum() const
Return the sum of the vector entries.
Definition: vector.cpp:900
virtual const double * Read(bool on_dev=true) const
Shortcut for mfem::Read(vec.GetMemory(), vec.Size(), on_dev).
Definition: vector.hpp:438
void ParallelAssemble(Vector &tv)
Assemble the vector on the true dofs, i.e. P^t v.
Definition: plinearform.cpp:87
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:292
ParFiniteElementSpace * ParFESpace() const
Definition: pgridfunc.hpp:113