MFEM  v4.2.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-2020, 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  unm1.SetSize(vfes_truevsize);
53  unm1 = 0.0;
54  unm2.SetSize(vfes_truevsize);
55  unm2 = 0.0;
56  fn.SetSize(vfes_truevsize);
57  Nun.SetSize(vfes_truevsize);
58  Nun = 0.0;
59  Nunm1.SetSize(vfes_truevsize);
60  Nunm1 = 0.0;
61  Nunm2.SetSize(vfes_truevsize);
62  Nunm2 = 0.0;
63  Fext.SetSize(vfes_truevsize);
64  FText.SetSize(vfes_truevsize);
65  Lext.SetSize(vfes_truevsize);
66  resu.SetSize(vfes_truevsize);
67 
68  tmp1.SetSize(vfes_truevsize);
69 
70  pn.SetSize(pfes_truevsize);
71  pn = 0.0;
72  resp.SetSize(pfes_truevsize);
73  resp = 0.0;
74  FText_bdr.SetSize(pfes_truevsize);
75  g_bdr.SetSize(pfes_truevsize);
76 
78  un_gf = 0.0;
79 
85 
86  pn_gf.SetSpace(pfes);
87  pn_gf = 0.0;
88  resp_gf.SetSpace(pfes);
89 
90  cur_step = 0;
91 
92  PrintInfo();
93 }
94 
95 void NavierSolver::Setup(double dt)
96 {
97  if (verbose && pmesh->GetMyRank() == 0)
98  {
99  mfem::out << "Setup" << std::endl;
100  if (partial_assembly)
101  {
102  mfem::out << "Using Partial Assembly" << std::endl;
103  }
104  else
105  {
106  mfem::out << "Using Full Assembly" << std::endl;
107  }
108  }
109 
110  sw_setup.Start();
111 
113  pfec_lor = new H1_FECollection(1);
115 
118 
119  Array<int> empty;
120 
121  // GLL integration rule (Numerical Integration)
123  const IntegrationRule &ir_ni = rules_ni.Get(vfes->GetFE(0)->GetGeomType(),
124  2 * order - 1);
125 
126  nlcoeff.constant = -1.0;
127  N = new ParNonlinearForm(vfes);
129  if (partial_assembly)
130  {
132  N->Setup();
133  }
134 
137  if (numerical_integ)
138  {
139  mv_blfi->SetIntRule(&ir_ni);
140  }
141  Mv_form->AddDomainIntegrator(mv_blfi);
142  if (partial_assembly)
143  {
145  }
146  Mv_form->Assemble();
147  Mv_form->FormSystemMatrix(empty, Mv);
148 
151  if (numerical_integ)
152  {
153  // blfi->SetIntRule(&ir_ni);
154  }
155  Sp_form->AddDomainIntegrator(sp_blfi);
156  if (partial_assembly)
157  {
159  }
160  Sp_form->Assemble();
162 
165  if (partial_assembly)
166  {
168  }
169  D_form->Assemble();
170  D_form->FormRectangularSystemMatrix(empty, empty, D);
171 
174  if (partial_assembly)
175  {
177  }
178  G_form->Assemble();
179  G_form->FormRectangularSystemMatrix(empty, empty, G);
180 
182  H_bdfcoeff.constant = 1.0 / dt;
183  H_form = new ParBilinearForm(vfes);
186  if (partial_assembly)
187  {
189  }
190  H_form->Assemble();
192 
196  *FText_gfcoeff),
197  vel_ess_attr);
198 
200  for (auto &vel_dbc : vel_dbcs)
201  {
203  *vel_dbc.coeff),
204  vel_dbc.attr);
205  }
206 
207  f_form = new ParLinearForm(vfes);
208  for (auto &accel_term : accel_terms)
209  {
210  auto *vdlfi = new VectorDomainLFIntegrator(*accel_term.coeff);
211  // @TODO: This order should always be the same as the nonlinear forms one!
212  // const IntegrationRule &ir = IntRules.Get(vfes->GetFE(0)->GetGeomType(),
213  // 4 * order);
214  // vdlfi->SetIntRule(&ir);
215  f_form->AddDomainIntegrator(vdlfi);
216  }
217 
218  if (partial_assembly)
219  {
220  Vector diag_pa(vfes->GetTrueVSize());
221  Mv_form->AssembleDiagonal(diag_pa);
222  MvInvPC = new OperatorJacobiSmoother(diag_pa, empty);
223  }
224  else
225  {
227  dynamic_cast<HypreSmoother *>(MvInvPC)->SetType(HypreSmoother::Jacobi, 1);
228  }
229  MvInv = new CGSolver(MPI_COMM_WORLD);
230  MvInv->iterative_mode = false;
231  MvInv->SetOperator(*Mv);
234  MvInv->SetRelTol(1e-12);
235  MvInv->SetMaxIter(200);
236 
237  if (partial_assembly)
238  {
246  SpInvPC->Mult(resp, pn);
247  SpInvOrthoPC = new OrthoSolver();
249  }
250  else
251  {
254  SpInvOrthoPC = new OrthoSolver();
256  }
257  SpInv = new CGSolver(MPI_COMM_WORLD);
258  SpInv->iterative_mode = true;
259  SpInv->SetOperator(*Sp);
260  if (pres_dbcs.empty())
261  {
263  }
264  else
265  {
267  }
270  SpInv->SetMaxIter(200);
271 
272  if (partial_assembly)
273  {
274  Vector diag_pa(vfes->GetTrueVSize());
275  H_form->AssembleDiagonal(diag_pa);
277  }
278  else
279  {
281  dynamic_cast<HypreSmoother *>(HInvPC)->SetType(HypreSmoother::Jacobi, 1);
282  }
283  HInv = new CGSolver(MPI_COMM_WORLD);
284  HInv->iterative_mode = true;
285  HInv->SetOperator(*H);
289  HInv->SetMaxIter(200);
290 
292 
293  sw_setup.Stop();
294 }
295 
296 void NavierSolver::Step(double &time, double dt, int cur_step)
297 {
298  sw_step.Start();
299 
300  time += dt;
301 
302  // Set current time for velocity dirichlet boundary conditions.
303  for (auto &vel_dbc : vel_dbcs)
304  {
305  vel_dbc.coeff->SetTime(time);
306  }
307 
308  // Set current time for pressure dirichlet boundary conditions.
309  for (auto &pres_dbc : pres_dbcs)
310  {
311  pres_dbc.coeff->SetTime(time);
312  }
313 
315 
316  if (cur_step <= 2)
317  {
318  H_bdfcoeff.constant = bd0 / dt;
319  H_form->Update();
320  H_form->Assemble();
322 
323  if (partial_assembly)
324  {
325  HInv->SetOperator(*H);
326  delete HInvPC;
327  Vector diag_pa(vfes->GetTrueVSize());
328  H_form->AssembleDiagonal(diag_pa);
331  }
332  else
333  {
334  HInv->SetOperator(*H);
335  }
336  }
337 
338  // Extrapolated f^{n+1}.
339  for (auto &accel_term : accel_terms)
340  {
341  accel_term.coeff->SetTime(time);
342  }
343 
344  f_form->Assemble();
346 
347  // Nonlinear extrapolated terms.
348  sw_extrap.Start();
349 
350  N->Mult(un, Nun);
351  Nun.Add(1.0, fn);
352 
353  {
354  const auto d_Nun = Nun.Read();
355  const auto d_Nunm1 = Nunm1.Read();
356  const auto d_Nunm2 = Nunm2.Read();
357  auto d_Fext = Fext.Write();
358  const auto ab1_ = ab1;
359  const auto ab2_ = ab2;
360  const auto ab3_ = ab3;
361  MFEM_FORALL(i, Fext.Size(),
362  d_Fext[i] = ab1_ * d_Nun[i] +
363  ab2_ * d_Nunm1[i] +
364  ab3_ * d_Nunm2[i];);
365  }
366 
367  // Rotate the solutions from previous time steps.
368  Nunm2 = Nunm1;
369  Nunm1 = Nun;
370 
371  // Fext = M^{-1} (F(u^{n}) + f^{n+1})
372  MvInv->Mult(Fext, tmp1);
375  Fext.Set(1.0, tmp1);
376 
377  // Compute BDF terms.
378  {
379  const double bd1idt = -bd1 / dt;
380  const double bd2idt = -bd2 / dt;
381  const double bd3idt = -bd3 / dt;
382  const auto d_un = un.Read();
383  const auto d_unm1 = unm1.Read();
384  const auto d_unm2 = unm2.Read();
385  auto d_Fext = Fext.ReadWrite();
386  MFEM_FORALL(i, Fext.Size(),
387  d_Fext[i] += bd1idt * d_un[i] +
388  bd2idt * d_unm1[i] +
389  bd3idt * d_unm2[i];);
390  }
391 
392  sw_extrap.Stop();
393 
394  // Pressure Poisson.
395  sw_curlcurl.Start();
396  {
397  const auto d_un = un.Read();
398  const auto d_unm1 = unm1.Read();
399  const auto d_unm2 = unm2.Read();
400  auto d_Lext = Lext.Write();
401  const auto ab1_ = ab1;
402  const auto ab2_ = ab2;
403  const auto ab3_ = ab3;
404  MFEM_FORALL(i, Lext.Size(),
405  d_Lext[i] = ab1_ * d_un[i] +
406  ab2_ * d_unm1[i] +
407  ab3_ * d_unm2[i];);
408  }
409 
411  if (pmesh->Dimension() == 2)
412  {
415  }
416  else
417  {
420  }
421 
423  Lext *= kin_vis;
424 
425  sw_curlcurl.Stop();
426 
427  // \tilde{F} = F - \nu CurlCurl(u)
428  FText.Set(-1.0, Lext);
429  FText.Add(1.0, Fext);
430 
431  // p_r = \nabla \cdot FText
432  D->Mult(FText, resp);
433  resp.Neg();
434 
435  // Add boundary terms.
439 
440  g_bdr_form->Assemble();
442  resp.Add(1.0, FText_bdr);
443  resp.Add(-bd0 / dt, g_bdr);
444 
445  if (pres_dbcs.empty())
446  {
448  }
449 
450  for (auto &pres_dbc : pres_dbcs)
451  {
452  pn_gf.ProjectBdrCoefficient(*pres_dbc.coeff, pres_dbc.attr);
453  }
454 
456 
457  Vector X1, B1;
458  if (partial_assembly)
459  {
460  auto *SpC = Sp.As<ConstrainedOperator>();
461  EliminateRHS(*Sp_form, *SpC, pres_ess_tdof, pn_gf, resp_gf, X1, B1, 1);
462  }
463  else
464  {
466  }
467  sw_spsolve.Start();
468  SpInv->Mult(B1, X1);
469  sw_spsolve.Stop();
473 
474  // If the boundary conditions on the pressure are pure Neumann remove the
475  // nullspace by removing the mean of the pressure solution. This is also
476  // ensured by the OrthoSolver wrapper for the preconditioner which removes
477  // the nullspace after every application.
478  if (pres_dbcs.empty())
479  {
480  MeanZero(pn_gf);
481  }
482 
484 
485  // Project velocity.
486  G->Mult(pn, resu);
487  resu.Neg();
488  Mv->Mult(Fext, tmp1);
489  resu.Add(1.0, tmp1);
490 
491  for (auto &vel_dbc : vel_dbcs)
492  {
493  un_gf.ProjectBdrCoefficient(*vel_dbc.coeff, vel_dbc.attr);
494  }
495 
497 
498  // Rotate solutions from previous time steps.
499  unm2 = unm1;
500  unm1 = un;
501 
502  Vector X2, B2;
503  if (partial_assembly)
504  {
505  auto *HC = H.As<ConstrainedOperator>();
506  EliminateRHS(*H_form, *HC, vel_ess_tdof, un_gf, resu_gf, X2, B2, 1);
507  }
508  else
509  {
511  }
512  sw_hsolve.Start();
513  HInv->Mult(B2, X2);
514  sw_hsolve.Stop();
518 
520 
521  sw_step.Stop();
522 
523  if (verbose && pmesh->GetMyRank() == 0)
524  {
525  mfem::out << std::setw(7) << "" << std::setw(3) << "It" << std::setw(8)
526  << "Resid" << std::setw(12) << "Reltol"
527  << "\n";
528  // If numerical integration is active, there is no solve (thus no
529  // iterations), on the inverse velocity mass application.
530  if (!numerical_integ)
531  {
532  mfem::out << std::setw(5) << "MVIN " << std::setw(5) << std::fixed
533  << iter_mvsolve << " " << std::setw(3)
534  << std::setprecision(2) << std::scientific << res_mvsolve
535  << " " << 1e-12 << "\n";
536  }
537  mfem::out << std::setw(5) << "PRES " << std::setw(5) << std::fixed
538  << iter_spsolve << " " << std::setw(3) << std::setprecision(2)
539  << std::scientific << res_spsolve << " " << rtol_spsolve
540  << "\n";
541  mfem::out << std::setw(5) << "HELM " << std::setw(5) << std::fixed
542  << iter_hsolve << " " << std::setw(3) << std::setprecision(2)
543  << std::scientific << res_hsolve << " " << rtol_hsolve
544  << "\n";
545  mfem::out << std::setprecision(8);
546  mfem::out << std::fixed;
547  }
548 }
549 
551 {
552  // Make sure not to recompute the inner product linear form every
553  // application.
554  if (mass_lf == nullptr)
555  {
556  onecoeff.constant = 1.0;
557  mass_lf = new ParLinearForm(v.ParFESpace());
559  mass_lf->Assemble();
560 
561  ParGridFunction one_gf(v.ParFESpace());
563 
564  volume = mass_lf->operator()(one_gf);
565  }
566 
567  double integ = mass_lf->operator()(v);
568 
569  v -= integ / volume;
570 }
571 
573  ConstrainedOperator &constrainedA,
574  const Array<int> &ess_tdof_list,
575  Vector &x,
576  Vector &b,
577  Vector &X,
578  Vector &B,
579  int copy_interior)
580 {
581  const Operator *Po = A.GetOutputProlongation();
582  const Operator *Pi = A.GetProlongation();
583  const Operator *Ri = A.GetRestriction();
584  A.InitTVectors(Po, Ri, Pi, x, b, X, B);
585  if (!copy_interior)
586  {
587  X.SetSubVectorComplement(ess_tdof_list, 0.0);
588  }
589  constrainedA.EliminateRHS(X, B);
590 }
591 
593 {
594  double loc_sum = v.Sum();
595  double global_sum = 0.0;
596  int loc_size = v.Size();
597  int global_size = 0;
598 
599  MPI_Allreduce(&loc_sum, &global_sum, 1, MPI_DOUBLE, MPI_SUM, MPI_COMM_WORLD);
600  MPI_Allreduce(&loc_size, &global_size, 1, MPI_INT, MPI_SUM, MPI_COMM_WORLD);
601 
602  v -= global_sum / static_cast<double>(global_size);
603 }
604 
606 {
607  FiniteElementSpace *fes = u.FESpace();
608 
609  // AccumulateAndCountZones.
610  Array<int> zones_per_vdof;
611  zones_per_vdof.SetSize(fes->GetVSize());
612  zones_per_vdof = 0;
613 
614  cu = 0.0;
615 
616  // Local interpolation.
617  int elndofs;
618  Array<int> vdofs;
619  Vector vals;
620  Vector loc_data;
621  int vdim = fes->GetVDim();
622  DenseMatrix grad_hat;
623  DenseMatrix dshape;
624  DenseMatrix grad;
625  Vector curl;
626 
627  for (int e = 0; e < fes->GetNE(); ++e)
628  {
629  fes->GetElementVDofs(e, vdofs);
630  u.GetSubVector(vdofs, loc_data);
631  vals.SetSize(vdofs.Size());
633  const FiniteElement *el = fes->GetFE(e);
634  elndofs = el->GetDof();
635  int dim = el->GetDim();
636  dshape.SetSize(elndofs, dim);
637 
638  for (int dof = 0; dof < elndofs; ++dof)
639  {
640  // Project.
641  const IntegrationPoint &ip = el->GetNodes().IntPoint(dof);
642  tr->SetIntPoint(&ip);
643 
644  // Eval and GetVectorGradientHat.
645  el->CalcDShape(tr->GetIntPoint(), dshape);
646  grad_hat.SetSize(vdim, dim);
647  DenseMatrix loc_data_mat(loc_data.GetData(), elndofs, vdim);
648  MultAtB(loc_data_mat, dshape, grad_hat);
649 
650  const DenseMatrix &Jinv = tr->InverseJacobian();
651  grad.SetSize(grad_hat.Height(), Jinv.Width());
652  Mult(grad_hat, Jinv, grad);
653 
654  curl.SetSize(3);
655  curl(0) = grad(2, 1) - grad(1, 2);
656  curl(1) = grad(0, 2) - grad(2, 0);
657  curl(2) = grad(1, 0) - grad(0, 1);
658 
659  for (int j = 0; j < curl.Size(); ++j)
660  {
661  vals(elndofs * j + dof) = curl(j);
662  }
663  }
664 
665  // Accumulate values in all dofs, count the zones.
666  for (int j = 0; j < vdofs.Size(); j++)
667  {
668  int ldof = vdofs[j];
669  cu(ldof) += vals[j];
670  zones_per_vdof[ldof]++;
671  }
672  }
673 
674  // Communication
675 
676  // Count the zones globally.
677  GroupCommunicator &gcomm = u.ParFESpace()->GroupComm();
678  gcomm.Reduce<int>(zones_per_vdof, GroupCommunicator::Sum);
679  gcomm.Bcast(zones_per_vdof);
680 
681  // Accumulate for all vdofs.
682  gcomm.Reduce<double>(cu.GetData(), GroupCommunicator::Sum);
683  gcomm.Bcast<double>(cu.GetData());
684 
685  // Compute means.
686  for (int i = 0; i < cu.Size(); i++)
687  {
688  const int nz = zones_per_vdof[i];
689  if (nz)
690  {
691  cu(i) /= nz;
692  }
693  }
694 }
695 
697  ParGridFunction &cu,
698  bool assume_scalar)
699 {
700  FiniteElementSpace *fes = u.FESpace();
701 
702  // AccumulateAndCountZones.
703  Array<int> zones_per_vdof;
704  zones_per_vdof.SetSize(fes->GetVSize());
705  zones_per_vdof = 0;
706 
707  cu = 0.0;
708 
709  // Local interpolation.
710  int elndofs;
711  Array<int> vdofs;
712  Vector vals;
713  Vector loc_data;
714  int vdim = fes->GetVDim();
715  DenseMatrix grad_hat;
716  DenseMatrix dshape;
717  DenseMatrix grad;
718  Vector curl;
719 
720  for (int e = 0; e < fes->GetNE(); ++e)
721  {
722  fes->GetElementVDofs(e, vdofs);
723  u.GetSubVector(vdofs, loc_data);
724  vals.SetSize(vdofs.Size());
726  const FiniteElement *el = fes->GetFE(e);
727  elndofs = el->GetDof();
728  int dim = el->GetDim();
729  dshape.SetSize(elndofs, dim);
730 
731  for (int dof = 0; dof < elndofs; ++dof)
732  {
733  // Project.
734  const IntegrationPoint &ip = el->GetNodes().IntPoint(dof);
735  tr->SetIntPoint(&ip);
736 
737  // Eval and GetVectorGradientHat.
738  el->CalcDShape(tr->GetIntPoint(), dshape);
739  grad_hat.SetSize(vdim, dim);
740  DenseMatrix loc_data_mat(loc_data.GetData(), elndofs, vdim);
741  MultAtB(loc_data_mat, dshape, grad_hat);
742 
743  const DenseMatrix &Jinv = tr->InverseJacobian();
744  grad.SetSize(grad_hat.Height(), Jinv.Width());
745  Mult(grad_hat, Jinv, grad);
746 
747  if (assume_scalar)
748  {
749  curl.SetSize(2);
750  curl(0) = grad(0, 1);
751  curl(1) = -grad(0, 0);
752  }
753  else
754  {
755  curl.SetSize(2);
756  curl(0) = grad(1, 0) - grad(0, 1);
757  curl(1) = 0.0;
758  }
759 
760  for (int j = 0; j < curl.Size(); ++j)
761  {
762  vals(elndofs * j + dof) = curl(j);
763  }
764  }
765 
766  // Accumulate values in all dofs, count the zones.
767  for (int j = 0; j < vdofs.Size(); j++)
768  {
769  int ldof = vdofs[j];
770  cu(ldof) += vals[j];
771  zones_per_vdof[ldof]++;
772  }
773  }
774 
775  // Communication.
776 
777  // Count the zones globally.
778  GroupCommunicator &gcomm = u.ParFESpace()->GroupComm();
779  gcomm.Reduce<int>(zones_per_vdof, GroupCommunicator::Sum);
780  gcomm.Bcast(zones_per_vdof);
781 
782  // Accumulate for all vdofs.
783  gcomm.Reduce<double>(cu.GetData(), GroupCommunicator::Sum);
784  gcomm.Bcast<double>(cu.GetData());
785 
786  // Compute means.
787  for (int i = 0; i < cu.Size(); i++)
788  {
789  const int nz = zones_per_vdof[i];
790  if (nz)
791  {
792  cu(i) /= nz;
793  }
794  }
795 }
796 
798 {
800  FiniteElementSpace *fes = u.FESpace();
801  int vdim = fes->GetVDim();
802 
803  Vector ux, uy, uz;
804  Vector ur, us, ut;
805  double cflx = 0.0;
806  double cfly = 0.0;
807  double cflz = 0.0;
808  double cflm = 0.0;
809  double cflmax = 0.0;
810 
811  for (int e = 0; e < fes->GetNE(); ++e)
812  {
813  const FiniteElement *fe = fes->GetFE(e);
814  const IntegrationRule &ir = IntRules.Get(fe->GetGeomType(),
815  fe->GetOrder());
817 
818  u.GetValues(e, ir, ux, 1);
819  ur.SetSize(ux.Size());
820  u.GetValues(e, ir, uy, 2);
821  us.SetSize(uy.Size());
822  if (vdim == 3)
823  {
824  u.GetValues(e, ir, uz, 3);
825  ut.SetSize(uz.Size());
826  }
827 
828  double hmin = pmesh->GetElementSize(e, 1) / (double) fes->GetOrder(0);
829 
830  for (int i = 0; i < ir.GetNPoints(); ++i)
831  {
832  const IntegrationPoint &ip = ir.IntPoint(i);
833  tr->SetIntPoint(&ip);
834  const DenseMatrix &invJ = tr->InverseJacobian();
835  const double detJinv = 1.0 / tr->Jacobian().Det();
836 
837  if (vdim == 2)
838  {
839  ur(i) = (ux(i) * invJ(0, 0) + uy(i) * invJ(1, 0)) * detJinv;
840  us(i) = (ux(i) * invJ(0, 1) + uy(i) * invJ(1, 1)) * detJinv;
841  }
842  else if (vdim == 3)
843  {
844  ur(i) = (ux(i) * invJ(0, 0) + uy(i) * invJ(1, 0)
845  + uz(i) * invJ(2, 0))
846  * detJinv;
847  us(i) = (ux(i) * invJ(0, 1) + uy(i) * invJ(1, 1)
848  + uz(i) * invJ(2, 1))
849  * detJinv;
850  ut(i) = (ux(i) * invJ(0, 2) + uy(i) * invJ(1, 2)
851  + uz(i) * invJ(2, 2))
852  * detJinv;
853  }
854 
855  cflx = fabs(dt * ux(i) / hmin);
856  cfly = fabs(dt * uy(i) / hmin);
857  if (vdim == 3)
858  {
859  cflz = fabs(dt * uz(i) / hmin);
860  }
861  cflm = cflx + cfly + cflz;
862  cflmax = fmax(cflmax, cflm);
863  }
864  }
865 
866  double cflmax_global = 0.0;
867  MPI_Allreduce(&cflmax,
868  &cflmax_global,
869  1,
870  MPI_DOUBLE,
871  MPI_MAX,
872  pmesh->GetComm());
873 
874  return cflmax_global;
875 }
876 
878 {
879  vel_dbcs.emplace_back(attr, coeff);
880 
881  if (verbose && pmesh->GetMyRank() == 0)
882  {
883  mfem::out << "Adding Velocity Dirichlet BC to attributes ";
884  for (int i = 0; i < attr.Size(); ++i)
885  {
886  if (attr[i] == 1)
887  {
888  mfem::out << i << " ";
889  }
890  }
891  mfem::out << std::endl;
892  }
893 
894  for (int i = 0; i < attr.Size(); ++i)
895  {
896  MFEM_ASSERT((vel_ess_attr[i] && attr[i]) == 0,
897  "Duplicate boundary definition deteceted.");
898  if (attr[i] == 1)
899  {
900  vel_ess_attr[i] = 1;
901  }
902  }
903 }
904 
906 {
908 }
909 
911 {
912  pres_dbcs.emplace_back(attr, coeff);
913 
914  if (verbose && pmesh->GetMyRank() == 0)
915  {
916  mfem::out << "Adding Pressure Dirichlet BC to attributes ";
917  for (int i = 0; i < attr.Size(); ++i)
918  {
919  if (attr[i] == 1)
920  {
921  mfem::out << i << " ";
922  }
923  }
924  mfem::out << std::endl;
925  }
926 
927  for (int i = 0; i < attr.Size(); ++i)
928  {
929  MFEM_ASSERT((pres_ess_attr[i] && attr[i]) == 0,
930  "Duplicate boundary definition deteceted.");
931  if (attr[i] == 1)
932  {
933  pres_ess_attr[i] = 1;
934  }
935  }
936 }
937 
939 {
941 }
942 
944 {
945  accel_terms.emplace_back(attr, coeff);
946 
947  if (verbose && pmesh->GetMyRank() == 0)
948  {
949  mfem::out << "Adding Acceleration term to attributes ";
950  for (int i = 0; i < attr.Size(); ++i)
951  {
952  if (attr[i] == 1)
953  {
954  mfem::out << i << " ";
955  }
956  }
957  mfem::out << std::endl;
958  }
959 }
960 
962 {
964 }
965 
967 {
968  if (step == 0)
969  {
970  bd0 = 1.0;
971  bd1 = -1.0;
972  bd2 = 0.0;
973  bd3 = 0.0;
974  ab1 = 1.0;
975  ab2 = 0.0;
976  ab3 = 0.0;
977  }
978  else if (step == 1)
979  {
980  bd0 = 3.0 / 2.0;
981  bd1 = -4.0 / 2.0;
982  bd2 = 1.0 / 2.0;
983  bd3 = 0.0;
984  ab1 = 2.0;
985  ab2 = -1.0;
986  ab3 = 0.0;
987  }
988  else if (step == 2)
989  {
990  bd0 = 11.0 / 6.0;
991  bd1 = -18.0 / 6.0;
992  bd2 = 9.0 / 6.0;
993  bd3 = -2.0 / 6.0;
994  ab1 = 3.0;
995  ab2 = -3.0;
996  ab3 = 1.0;
997  }
998 }
999 
1001 {
1002  double my_rt[6], rt_max[6];
1003 
1004  my_rt[0] = sw_setup.RealTime();
1005  my_rt[1] = sw_step.RealTime();
1006  my_rt[2] = sw_extrap.RealTime();
1007  my_rt[3] = sw_curlcurl.RealTime();
1008  my_rt[4] = sw_spsolve.RealTime();
1009  my_rt[5] = sw_hsolve.RealTime();
1010 
1011  MPI_Reduce(my_rt, rt_max, 6, MPI_DOUBLE, MPI_MAX, 0, pmesh->GetComm());
1012 
1013  if (pmesh->GetMyRank() == 0)
1014  {
1015  mfem::out << std::setw(10) << "SETUP" << std::setw(10) << "STEP"
1016  << std::setw(10) << "EXTRAP" << std::setw(10) << "CURLCURL"
1017  << std::setw(10) << "PSOLVE" << std::setw(10) << "HSOLVE"
1018  << "\n";
1019 
1020  mfem::out << std::setprecision(3) << std::setw(10) << my_rt[0]
1021  << std::setw(10) << my_rt[1] << std::setw(10) << my_rt[2]
1022  << std::setw(10) << my_rt[3] << std::setw(10) << my_rt[4]
1023  << std::setw(10) << my_rt[5] << "\n";
1024 
1025  mfem::out << std::setprecision(3) << std::setw(10) << " " << std::setw(10)
1026  << my_rt[1] / my_rt[1] << std::setw(10) << my_rt[2] / my_rt[1]
1027  << std::setw(10) << my_rt[3] / my_rt[1] << std::setw(10)
1028  << my_rt[4] / my_rt[1] << std::setw(10) << my_rt[5] / my_rt[1]
1029  << "\n";
1030 
1031  mfem::out << std::setprecision(8);
1032  }
1033 }
1034 
1036 {
1037  int fes_size0 = vfes->GlobalVSize();
1038  int fes_size1 = pfes->GlobalVSize();
1039 
1040  if (pmesh->GetMyRank() == 0)
1041  {
1042  mfem::out << "NAVIER version: " << NAVIER_VERSION << std::endl
1043  << "MFEM version: " << MFEM_VERSION << std::endl
1044  << "MFEM GIT: " << MFEM_GIT_STRING << std::endl
1045  << "Velocity #DOFs: " << fes_size0 << std::endl
1046  << "Pressure #DOFs: " << fes_size1 << std::endl;
1047  }
1048 }
1049 
1051 {
1052  delete FText_gfcoeff;
1053  delete g_bdr_form;
1054  delete FText_bdr_form;
1055  delete mass_lf;
1056  delete Mv_form;
1057  delete N;
1058  delete Sp_form;
1059  delete D_form;
1060  delete G_form;
1061  delete HInvPC;
1062  delete HInv;
1063  delete H_form;
1064  delete SpInv;
1065  delete MvInvPC;
1066  delete Sp_form_lor;
1067  delete SpInvOrthoPC;
1068  delete SpInvPC;
1069  delete f_form;
1070  delete pfes_lor;
1071  delete pfec_lor;
1072  delete pmesh_lor;
1073  delete MvInv;
1074  delete vfec;
1075  delete pfec;
1076  delete vfes;
1077  delete pfes;
1078 }
int GetNPoints() const
Returns the number of the points in the integration rule.
Definition: intrules.hpp:245
Abstract class for all finite elements.
Definition: fe.hpp:235
void PrintTimingData()
Print timing summary of the solving routine.
int Size() const
Return the logical size of the array.
Definition: array.hpp:124
void Orthogonalize(Vector &v)
Remove mean from a Vector.
Class for domain integration L(v) := (f, v)
Definition: lininteg.hpp:93
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:400
virtual void GetEssentialTrueDofs(const Array< int > &bdr_attr_is_ess, Array< int > &ess_tdof_list, int component=-1)
Definition: pfespace.cpp:775
Conjugate gradient method.
Definition: solvers.hpp:258
int GetDim() const
Returns the reference space dimension for the finite element.
Definition: fe.hpp:309
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:96
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.
void Bcast(T *ldata, int layout) const
Broadcast within each group where the master is the root.
int GetNumIterations() const
Definition: solvers.hpp:101
const IntegrationRule & Get(int GeomType, int Order)
Returns an integration rule for given GeomType and Order.
Definition: intrules.cpp:915
ParMesh * GetParMesh() const
Definition: pfespace.hpp:243
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 informations about the Navier version.
virtual void Mult(const Vector &b, Vector &x) const
Operator application: y=A(x).
Definition: solvers.cpp:535
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:459
void GetElementVDofs(int i, Array< int > &vdofs) const
Returns indexes of degrees of freedom in array dofs for i&#39;th element.
Definition: fespace.cpp:173
virtual void Setup()
Setup the NonlinearForm.
double Det() const
Definition: densemat.cpp:451
ConstantCoefficient nlcoeff
void Mult(const Table &A, const Table &B, Table &C)
C = A * B (as boolean matrices)
Definition: table.cpp:476
int Width() const
Get the width (size of input) of the Operator. Synonym with NumCols().
Definition: operator.hpp:71
void AddDomainIntegrator(NonlinearFormIntegrator *nlfi)
Adds new Domain Integrator.
void Assemble()
Assembles the linear form i.e. sums over all domain/bdr integrators.
Definition: linearform.cpp:79
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:85
void GetSubVector(const Array< int > &dofs, Vector &elemvect) const
Extract entries listed in dofs to the output Vector elemvect.
Definition: vector.cpp:495
int GetOrder() const
Returns the order of the finite element. In the case of anisotropic orders, returns the maximum order...
Definition: fe.hpp:319
Container class for integration rules.
Definition: intrules.hpp:309
Data type dense matrix using column-major storage.
Definition: densemat.hpp:23
int Size() const
Returns the size of the vector.
Definition: vector.hpp:160
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:111
double(const Vector &x, double t) ScalarFuncT
void AssembleDiagonal(Vector &diag) const
Assemble the diagonal of the bilinear form into diag.
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:132
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)
Definition: pgridfunc.cpp:474
bool iterative_mode
If true, use the second argument of Mult() as an initial guess.
Definition: operator.hpp:638
const IntegrationPoint & GetIntPoint()
Get a const reference to the currently set integration point. This will return NULL if no integration...
Definition: eltrans.hpp:90
double * Write(bool on_dev=true)
Shortcut for mfem::Write(vec.GetMemory(), vec.Size(), on_dev).
Definition: vector.hpp:380
double * GetData() const
Return a pointer to the beginning of the Vector data.
Definition: vector.hpp:169
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:103
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:255
void Stop()
Stop the stopwatch.
Definition: tic_toc.cpp:416
bool partial_assembly
Enable/disable partial assembly of forms.
void Step(double &time, double dt, int cur_step)
Compute solution at the next time step t+dt.
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.hpp:312
virtual const Operator * GetProlongation() const
Prolongation operator from linear algebra (linear system) vectors, to input vectors for the operator...
Definition: operator.hpp:105
int GetNE() const
Returns number of elements in the mesh.
Definition: fespace.hpp:427
The BoomerAMG solver in hypre.
Definition: hypre.hpp:1079
IntegrationPoint & IntPoint(int i)
Returns a reference to the i-th integration point.
Definition: intrules.hpp:248
void ComputeCurl3D(ParGridFunction &u, ParGridFunction &cu)
Compute for .
Class for parallel linear form.
Definition: plinearform.hpp:26
Communicator performing operations within groups defined by a GroupTopology with arbitrary-size data ...
void SetPrintLevel(int print_lvl)
Definition: solvers.cpp:70
FiniteElementCollection * vfec
Velocity finite element collection.
double * ReadWrite(bool on_dev=true)
Shortcut for mfem::ReadWrite(vec.GetMemory(), vec.Size(), on_dev).
Definition: vector.hpp:388
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:164
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:65
double ComputeCFL(ParGridFunction &u, double dt)
Compute CFL.
Class for boundary integration .
Definition: lininteg.hpp:172
Jacobi smoothing for a given bilinear form (no matrix necessary).
Definition: solvers.hpp:128
double b
Definition: lissajous.cpp:42
void SetPrintLevel(int print_level)
Definition: hypre.hpp:1119
const DenseMatrix & Jacobian()
Return the Jacobian matrix of the transformation at the currently set IntegrationPoint, using the method SetIntPoint().
Definition: eltrans.hpp:111
void SetMaxIter(int max_it)
Definition: solvers.hpp:98
virtual void Mult(const HypreParVector &b, HypreParVector &x) const
Solve the linear system Ax=b.
Definition: hypre.cpp:2506
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.hpp:382
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:592
void GetValues(int i, const IntegrationRule &ir, Vector &vals, int vdim=1) const
Definition: gridfunc.cpp:455
int GetVDim() const
Returns vector dimension.
Definition: fespace.hpp:389
virtual const Operator * GetRestriction() const
Restriction operator from input vectors for the operator to linear algebra (linear system) vectors...
Definition: operator.hpp:108
FiniteElementSpace * FESpace()
Definition: gridfunc.hpp:595
void AddBoundaryIntegrator(LinearFormIntegrator *lfi)
Adds new Boundary Integrator. Assumes ownership of lfi.
Definition: linearform.cpp:53
int Dimension() const
Definition: mesh.hpp:788
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::NONE.
void Start()
Clear the elapsed time and start the stopwatch.
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:656
A general vector function coefficient.
HYPRE_Int GlobalVSize() const
Definition: pfespace.hpp:249
void SetAssemblyLevel(AssemblyLevel assembly_level)
Set the desired assembly level. The default is AssemblyLevel::LEGACYFULL.
double GetElementSize(int i, int type=0)
Get the size of the i-th element relative to the perfect reference element.
Definition: mesh.cpp:74
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:460
void AddDomainIntegrator(LinearFormIntegrator *lfi)
Adds new Domain Integrator. Assumes ownership of lfi.
Definition: linearform.cpp:39
Abstract base class BilinearFormIntegrator.
Definition: bilininteg.hpp:31
Array< int > bdr_attributes
A list of all unique boundary attributes used by the Mesh.
Definition: mesh.hpp:201
int GetMyRank() const
Definition: pmesh.hpp:238
void SetRelTol(double rtol)
Definition: solvers.hpp:96
ParFiniteElementSpace * vfes
Velocity finite element space.
MPI_Comm GetComm() const
Definition: pmesh.hpp:236
void ProjectBdrCoefficient(Coefficient *coeff[], VectorCoefficient *vcoeff, Array< int > &attr)
Definition: pgridfunc.cpp:582
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Definition: fespace.hpp:87
int GetDof() const
Returns the number of degrees of freedom in the finite element.
Definition: fe.hpp:315
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 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:654
GroupCommunicator & GroupComm()
Return a reference to the internal GroupCommunicator (on VDofs)
Definition: pfespace.hpp:298
int GetOrder(int i) const
Returns the order of the i&#39;th finite element.
Definition: fespace.cpp:96
Vector & Set(const double a, const Vector &x)
(*this) = a * x
Definition: vector.cpp:228
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.
Vector & Add(const double a, const Vector &Va)
(*this) += a * Va
Definition: vector.cpp:213
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
const double * Read(bool on_dev=true) const
Shortcut for mfem::Read(vec.GetMemory(), vec.Size(), on_dev).
Definition: vector.hpp:372
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:421
virtual const SparseMatrix * GetRestrictionMatrix() const
Get the R matrix which restricts a local dof vector to true dof vector.
Definition: pfespace.hpp:339
VectorGridFunctionCoefficient * FText_gfcoeff
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:721
void AddDomainIntegrator(BilinearFormIntegrator *bfi)
Adds new Domain Integrator. Assumes ownership of bfi.
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:1798
HypreParVector * GetTrueDofs() const
Returns the true dofs in a new HypreParVector.
Definition: pgridfunc.cpp:143
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:37
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:272
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:2650
Vector data type.
Definition: vector.hpp:51
Vector coefficient defined by a vector GridFunction.
virtual void SetPreconditioner(Solver &pr)
This should be called before SetOperator.
Definition: solvers.cpp:91
Arbitrary order H1-conforming (continuous) finite elements.
Definition: fe_coll.hpp:159
void AddPresDirichletBC(Coefficient *coeff, Array< int > &attr)
Add a Dirichlet boundary condition to the pressure field.
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:90
Square Operator for imposing essential boundary conditions using only the action, Mult()...
Definition: operator.hpp:800
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:181
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:378
Solver wrapper which orthogonalizes the input and output vector.
double Sum() const
Return the sum of the vector entries.
Definition: vector.cpp:849
void ParallelAssemble(Vector &tv)
Assemble the vector on the true dofs, i.e. P^t v.
Definition: plinearform.cpp:46
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
double f(const Vector &p)
void Neg()
(*this) = -(*this)
Definition: vector.cpp:253
ParFiniteElementSpace * ParFESpace() const
Definition: pgridfunc.hpp:108