MFEM v4.7.0
Finite element discretization library
Loading...
Searching...
No Matches
navier_solver.cpp
Go to the documentation of this file.
1// Copyright (c) 2010-2024, 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"
14#include <fstream>
15#include <iomanip>
16
17using namespace mfem;
18using namespace navier;
19
21{
23 for (int i = 0; i < bffis->Size(); ++i)
24 {
25 dst->AddDomainIntegrator((*bffis)[i]);
26 }
27}
28
29NavierSolver::NavierSolver(ParMesh *mesh, int order, real_t kin_vis)
30 : pmesh(mesh), order(order), kin_vis(kin_vis),
31 gll_rules(0, Quadrature1D::GaussLobatto)
32{
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
92 pn_gf = 0.0;
94
95 cur_step = 0;
96
97 PrintInfo();
98}
99
101{
102 if (verbose && pmesh->GetMyRank() == 0)
103 {
104 mfem::out << "Setup" << std::endl;
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);
135 {
137 N->Setup();
138 }
139
141 auto *mv_blfi = new VectorMassIntegrator;
142 if (numerical_integ)
143 {
144 mv_blfi->SetIntRule(&ir_ni);
145 }
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 }
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);
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);
190 {
192 }
193 G_form->Assemble();
194 G_form->FormRectangularSystemMatrix(empty, empty, G);
195
197 H_bdfcoeff.constant = 1.0 / dt;
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);
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
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 }
248 }
249
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;
267 MvInv->SetMaxIter(200);
268
270 {
274 SpInvPC->Mult(resp, pn);
277 }
278 else
279 {
284 }
285 SpInv = new CGSolver(vfes->GetComm());
286 SpInv->iterative_mode = true;
288 if (pres_dbcs.empty())
289 {
291 }
292 else
293 {
295 }
298 SpInv->SetMaxIter(200);
299
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());
334 pmesh->Dimension());
335
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
367void NavierSolver::Step(real_t &time, real_t 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);
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.
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(Fext.Size(), [=] MFEM_HOST_DEVICE (int i)
426 {
427 d_Fext[i] = ab1_ * d_Nun[i] +
428 ab2_ * d_Nunm1[i] +
429 ab3_ * d_Nunm2[i];
430 });
431 }
432
433 Fext.Add(1.0, fn);
434
435 // Fext = M^{-1} (F(u^{n}) + f^{n+1})
436 MvInv->Mult(Fext, tmp1);
439 Fext.Set(1.0, tmp1);
440
441 // Compute BDF terms.
442 {
443 const real_t bd1idt = -bd1 / dt;
444 const real_t bd2idt = -bd2 / dt;
445 const real_t bd3idt = -bd3 / dt;
446 const auto d_un = un.Read();
447 const auto d_unm1 = unm1.Read();
448 const auto d_unm2 = unm2.Read();
449 auto d_Fext = Fext.ReadWrite();
450 mfem::forall(Fext.Size(), [=] MFEM_HOST_DEVICE (int i)
451 {
452 d_Fext[i] += bd1idt * d_un[i] +
453 bd2idt * d_unm1[i] +
454 bd3idt * d_unm2[i];
455 });
456 }
457
458 sw_extrap.Stop();
459
460 // Pressure Poisson.
462 {
463 const auto d_un = un.Read();
464 const auto d_unm1 = unm1.Read();
465 const auto d_unm2 = unm2.Read();
466 auto d_Lext = Lext.Write();
467 const auto ab1_ = ab1;
468 const auto ab2_ = ab2;
469 const auto ab3_ = ab3;
470 mfem::forall(Lext.Size(), [=] MFEM_HOST_DEVICE (int i)
471 {
472 d_Lext[i] = ab1_ * d_un[i] +
473 ab2_ * d_unm1[i] +
474 ab3_ * d_unm2[i];
475 });
476 }
477
479 if (pmesh->Dimension() == 2)
480 {
483 }
484 else
485 {
488 }
489
491 Lext *= kin_vis;
492
494
495 // \tilde{F} = F - \nu CurlCurl(u)
496 FText.Set(-1.0, Lext);
497 FText.Add(1.0, Fext);
498
499 // p_r = \nabla \cdot FText
500 D->Mult(FText, resp);
501 resp.Neg();
502
503 // Add boundary terms.
507
510 resp.Add(1.0, FText_bdr);
511 resp.Add(-bd0 / dt, g_bdr);
512
513 if (pres_dbcs.empty())
514 {
516 }
517
518 for (auto &pres_dbc : pres_dbcs)
519 {
520 pn_gf.ProjectBdrCoefficient(*pres_dbc.coeff, pres_dbc.attr);
521 }
522
524
525 Vector X1, B1;
527 {
528 auto *SpC = Sp.As<ConstrainedOperator>();
529 EliminateRHS(*Sp_form, *SpC, pres_ess_tdof, pn_gf, resp_gf, X1, B1, 1);
530 }
531 else
532 {
534 }
536 SpInv->Mult(B1, X1);
541
542 // If the boundary conditions on the pressure are pure Neumann remove the
543 // nullspace by removing the mean of the pressure solution. This is also
544 // ensured by the OrthoSolver wrapper for the preconditioner which removes
545 // the nullspace after every application.
546 if (pres_dbcs.empty())
547 {
549 }
550
552
553 // Project velocity.
554 G->Mult(pn, resu);
555 resu.Neg();
556 Mv->Mult(Fext, tmp1);
557 resu.Add(1.0, tmp1);
558
559 // un_next_gf = un_gf;
560
561 for (auto &vel_dbc : vel_dbcs)
562 {
564 }
565
567
568 Vector X2, B2;
570 {
571 auto *HC = H.As<ConstrainedOperator>();
573 }
574 else
575 {
577 }
579 HInv->Mult(B2, X2);
580 sw_hsolve.Stop();
584
586
587 // If the current time step is not provisional, accept the computed solution
588 // and update the time step history by default.
589 if (!provisional)
590 {
592 time += dt;
593 }
594
595 if (filter_alpha != 0.0)
596 {
599 const auto d_un_filtered_gf = un_filtered_gf.Read();
600 auto d_un_gf = un_gf.ReadWrite();
601 const auto filter_alpha_ = filter_alpha;
602 mfem::forall(un_gf.Size(), [=] MFEM_HOST_DEVICE (int i)
603 {
604 d_un_gf[i] = (1.0 - filter_alpha_) * d_un_gf[i]
605 + filter_alpha_ * d_un_filtered_gf[i];
606 });
607 }
608
609 sw_step.Stop();
610
611 if (verbose && pmesh->GetMyRank() == 0)
612 {
613 mfem::out << std::setw(7) << "" << std::setw(3) << "It" << std::setw(8)
614 << "Resid" << std::setw(12) << "Reltol"
615 << "\n";
616 // If numerical integration is active, there is no solve (thus no
617 // iterations), on the inverse velocity mass application.
618 if (!numerical_integ)
619 {
620 mfem::out << std::setw(5) << "MVIN " << std::setw(5) << std::fixed
621 << iter_mvsolve << " " << std::setw(3)
622 << std::setprecision(2) << std::scientific << res_mvsolve
623 << " " << rtol_mvsolve << "\n";
624 }
625 mfem::out << std::setw(5) << "PRES " << std::setw(5) << std::fixed
626 << iter_spsolve << " " << std::setw(3) << std::setprecision(2)
627 << std::scientific << res_spsolve << " " << rtol_spsolve
628 << "\n";
629 mfem::out << std::setw(5) << "HELM " << std::setw(5) << std::fixed
630 << iter_hsolve << " " << std::setw(3) << std::setprecision(2)
631 << std::scientific << res_hsolve << " " << rtol_hsolve
632 << "\n";
633 mfem::out << std::setprecision(8);
634 mfem::out << std::scientific;
635 }
636}
637
639{
640 // Make sure not to recompute the inner product linear form every
641 // application.
642 if (mass_lf == nullptr)
643 {
644 onecoeff.constant = 1.0;
646 auto *dlfi = new DomainLFIntegrator(onecoeff);
647 if (numerical_integ)
648 {
649 const IntegrationRule &ir_ni = gll_rules.Get(vfes->GetFE(0)->GetGeomType(),
650 2 * order - 1);
651 dlfi->SetIntRule(&ir_ni);
652 }
654 mass_lf->Assemble();
655
656 ParGridFunction one_gf(v.ParFESpace());
658
659 volume = mass_lf->operator()(one_gf);
660 }
661
662 real_t integ = mass_lf->operator()(v);
663
664 v -= integ / volume;
665}
666
668 ConstrainedOperator &constrainedA,
669 const Array<int> &ess_tdof_list,
670 Vector &x,
671 Vector &b,
672 Vector &X,
673 Vector &B,
674 int copy_interior)
675{
676 const Operator *Po = A.GetOutputProlongation();
677 const Operator *Pi = A.GetProlongation();
678 const Operator *Ri = A.GetRestriction();
679 A.InitTVectors(Po, Ri, Pi, x, b, X, B);
680 if (!copy_interior)
681 {
682 X.SetSubVectorComplement(ess_tdof_list, 0.0);
683 }
684 constrainedA.EliminateRHS(X, B);
685}
686
688{
689 real_t loc_sum = v.Sum();
690 real_t global_sum = 0.0;
691 int loc_size = v.Size();
692 int global_size = 0;
693
694 MPI_Allreduce(&loc_sum, &global_sum, 1, MPITypeMap<real_t>::mpi_type, MPI_SUM,
695 pfes->GetComm());
696 MPI_Allreduce(&loc_size, &global_size, 1, MPI_INT, MPI_SUM, pfes->GetComm());
697
698 v -= global_sum / static_cast<real_t>(global_size);
699}
700
702{
703 FiniteElementSpace *fes = u.FESpace();
704
705 // AccumulateAndCountZones.
706 Array<int> zones_per_vdof;
707 zones_per_vdof.SetSize(fes->GetVSize());
708 zones_per_vdof = 0;
709
710 cu = 0.0;
711
712 // Local interpolation.
713 int elndofs;
714 Array<int> vdofs;
715 Vector vals;
716 Vector loc_data;
717 int vdim = fes->GetVDim();
718 DenseMatrix grad_hat;
719 DenseMatrix dshape;
720 DenseMatrix grad;
721 Vector curl;
722
723 for (int e = 0; e < fes->GetNE(); ++e)
724 {
725 fes->GetElementVDofs(e, vdofs);
726 u.GetSubVector(vdofs, loc_data);
727 vals.SetSize(vdofs.Size());
729 const FiniteElement *el = fes->GetFE(e);
730 elndofs = el->GetDof();
731 int dim = el->GetDim();
732 dshape.SetSize(elndofs, dim);
733
734 for (int dof = 0; dof < elndofs; ++dof)
735 {
736 // Project.
737 const IntegrationPoint &ip = el->GetNodes().IntPoint(dof);
738 tr->SetIntPoint(&ip);
739
740 // Eval and GetVectorGradientHat.
741 el->CalcDShape(tr->GetIntPoint(), dshape);
742 grad_hat.SetSize(vdim, dim);
743 DenseMatrix loc_data_mat(loc_data.GetData(), elndofs, vdim);
744 MultAtB(loc_data_mat, dshape, grad_hat);
745
746 const DenseMatrix &Jinv = tr->InverseJacobian();
747 grad.SetSize(grad_hat.Height(), Jinv.Width());
748 Mult(grad_hat, Jinv, grad);
749
750 curl.SetSize(3);
751 curl(0) = grad(2, 1) - grad(1, 2);
752 curl(1) = grad(0, 2) - grad(2, 0);
753 curl(2) = grad(1, 0) - grad(0, 1);
754
755 for (int j = 0; j < curl.Size(); ++j)
756 {
757 vals(elndofs * j + dof) = curl(j);
758 }
759 }
760
761 // Accumulate values in all dofs, count the zones.
762 for (int j = 0; j < vdofs.Size(); j++)
763 {
764 int ldof = vdofs[j];
765 cu(ldof) += vals[j];
766 zones_per_vdof[ldof]++;
767 }
768 }
769
770 // Communication
771
772 // Count the zones globally.
773 GroupCommunicator &gcomm = u.ParFESpace()->GroupComm();
774 gcomm.Reduce<int>(zones_per_vdof, GroupCommunicator::Sum);
775 gcomm.Bcast(zones_per_vdof);
776
777 // Accumulate for all vdofs.
779 gcomm.Bcast<real_t>(cu.GetData());
780
781 // Compute means.
782 for (int i = 0; i < cu.Size(); i++)
783 {
784 const int nz = zones_per_vdof[i];
785 if (nz)
786 {
787 cu(i) /= nz;
788 }
789 }
790}
791
793 ParGridFunction &cu,
794 bool assume_scalar)
795{
796 FiniteElementSpace *fes = u.FESpace();
797
798 // AccumulateAndCountZones.
799 Array<int> zones_per_vdof;
800 zones_per_vdof.SetSize(fes->GetVSize());
801 zones_per_vdof = 0;
802
803 cu = 0.0;
804
805 // Local interpolation.
806 int elndofs;
807 Array<int> vdofs;
808 Vector vals;
809 Vector loc_data;
810 int vdim = fes->GetVDim();
811 DenseMatrix grad_hat;
812 DenseMatrix dshape;
813 DenseMatrix grad;
814 Vector curl;
815
816 for (int e = 0; e < fes->GetNE(); ++e)
817 {
818 fes->GetElementVDofs(e, vdofs);
819 u.GetSubVector(vdofs, loc_data);
820 vals.SetSize(vdofs.Size());
822 const FiniteElement *el = fes->GetFE(e);
823 elndofs = el->GetDof();
824 int dim = el->GetDim();
825 dshape.SetSize(elndofs, dim);
826
827 for (int dof = 0; dof < elndofs; ++dof)
828 {
829 // Project.
830 const IntegrationPoint &ip = el->GetNodes().IntPoint(dof);
831 tr->SetIntPoint(&ip);
832
833 // Eval and GetVectorGradientHat.
834 el->CalcDShape(tr->GetIntPoint(), dshape);
835 grad_hat.SetSize(vdim, dim);
836 DenseMatrix loc_data_mat(loc_data.GetData(), elndofs, vdim);
837 MultAtB(loc_data_mat, dshape, grad_hat);
838
839 const DenseMatrix &Jinv = tr->InverseJacobian();
840 grad.SetSize(grad_hat.Height(), Jinv.Width());
841 Mult(grad_hat, Jinv, grad);
842
843 if (assume_scalar)
844 {
845 curl.SetSize(2);
846 curl(0) = grad(0, 1);
847 curl(1) = -grad(0, 0);
848 }
849 else
850 {
851 curl.SetSize(2);
852 curl(0) = grad(1, 0) - grad(0, 1);
853 curl(1) = 0.0;
854 }
855
856 for (int j = 0; j < curl.Size(); ++j)
857 {
858 vals(elndofs * j + dof) = curl(j);
859 }
860 }
861
862 // Accumulate values in all dofs, count the zones.
863 for (int j = 0; j < vdofs.Size(); j++)
864 {
865 int ldof = vdofs[j];
866 cu(ldof) += vals[j];
867 zones_per_vdof[ldof]++;
868 }
869 }
870
871 // Communication.
872
873 // Count the zones globally.
874 GroupCommunicator &gcomm = u.ParFESpace()->GroupComm();
875 gcomm.Reduce<int>(zones_per_vdof, GroupCommunicator::Sum);
876 gcomm.Bcast(zones_per_vdof);
877
878 // Accumulate for all vdofs.
880 gcomm.Bcast<real_t>(cu.GetData());
881
882 // Compute means.
883 for (int i = 0; i < cu.Size(); i++)
884 {
885 const int nz = zones_per_vdof[i];
886 if (nz)
887 {
888 cu(i) /= nz;
889 }
890 }
891}
892
894{
895 ParMesh *pmesh_u = u.ParFESpace()->GetParMesh();
896 FiniteElementSpace *fes = u.FESpace();
897 int vdim = fes->GetVDim();
898
899 Vector ux, uy, uz;
900 Vector ur, us, ut;
901 real_t cflx = 0.0;
902 real_t cfly = 0.0;
903 real_t cflz = 0.0;
904 real_t cflm = 0.0;
905 real_t cflmax = 0.0;
906
907 for (int e = 0; e < fes->GetNE(); ++e)
908 {
909 const FiniteElement *fe = fes->GetFE(e);
910 const IntegrationRule &ir = IntRules.Get(fe->GetGeomType(),
911 fe->GetOrder());
913
914 u.GetValues(e, ir, ux, 1);
915 ur.SetSize(ux.Size());
916 u.GetValues(e, ir, uy, 2);
917 us.SetSize(uy.Size());
918 if (vdim == 3)
919 {
920 u.GetValues(e, ir, uz, 3);
921 ut.SetSize(uz.Size());
922 }
923
924 real_t hmin = pmesh_u->GetElementSize(e, 1) /
925 (real_t) fes->GetElementOrder(0);
926
927 for (int i = 0; i < ir.GetNPoints(); ++i)
928 {
929 const IntegrationPoint &ip = ir.IntPoint(i);
930 tr->SetIntPoint(&ip);
931 const DenseMatrix &invJ = tr->InverseJacobian();
932 const real_t detJinv = 1.0 / tr->Jacobian().Det();
933
934 if (vdim == 2)
935 {
936 ur(i) = (ux(i) * invJ(0, 0) + uy(i) * invJ(1, 0)) * detJinv;
937 us(i) = (ux(i) * invJ(0, 1) + uy(i) * invJ(1, 1)) * detJinv;
938 }
939 else if (vdim == 3)
940 {
941 ur(i) = (ux(i) * invJ(0, 0) + uy(i) * invJ(1, 0)
942 + uz(i) * invJ(2, 0))
943 * detJinv;
944 us(i) = (ux(i) * invJ(0, 1) + uy(i) * invJ(1, 1)
945 + uz(i) * invJ(2, 1))
946 * detJinv;
947 ut(i) = (ux(i) * invJ(0, 2) + uy(i) * invJ(1, 2)
948 + uz(i) * invJ(2, 2))
949 * detJinv;
950 }
951
952 cflx = fabs(dt * ux(i) / hmin);
953 cfly = fabs(dt * uy(i) / hmin);
954 if (vdim == 3)
955 {
956 cflz = fabs(dt * uz(i) / hmin);
957 }
958 cflm = cflx + cfly + cflz;
959 cflmax = fmax(cflmax, cflm);
960 }
961 }
962
963 real_t cflmax_global = 0.0;
964 MPI_Allreduce(&cflmax,
965 &cflmax_global,
966 1,
968 MPI_MAX,
969 pmesh_u->GetComm());
970
971 return cflmax_global;
972}
973
975{
976 vel_dbcs.emplace_back(attr, coeff);
977
978 if (verbose && pmesh->GetMyRank() == 0)
979 {
980 mfem::out << "Adding Velocity Dirichlet BC to attributes ";
981 for (int i = 0; i < attr.Size(); ++i)
982 {
983 if (attr[i] == 1)
984 {
985 mfem::out << i << " ";
986 }
987 }
988 mfem::out << std::endl;
989 }
990
991 for (int i = 0; i < attr.Size(); ++i)
992 {
993 MFEM_ASSERT((vel_ess_attr[i] && attr[i]) == 0,
994 "Duplicate boundary definition deteceted.");
995 if (attr[i] == 1)
996 {
997 vel_ess_attr[i] = 1;
998 }
999 }
1000}
1001
1006
1008{
1009 pres_dbcs.emplace_back(attr, coeff);
1010
1011 if (verbose && pmesh->GetMyRank() == 0)
1012 {
1013 mfem::out << "Adding Pressure Dirichlet BC to attributes ";
1014 for (int i = 0; i < attr.Size(); ++i)
1015 {
1016 if (attr[i] == 1)
1017 {
1018 mfem::out << i << " ";
1019 }
1020 }
1021 mfem::out << std::endl;
1022 }
1023
1024 for (int i = 0; i < attr.Size(); ++i)
1025 {
1026 MFEM_ASSERT((pres_ess_attr[i] && attr[i]) == 0,
1027 "Duplicate boundary definition deteceted.");
1028 if (attr[i] == 1)
1029 {
1030 pres_ess_attr[i] = 1;
1031 }
1032 }
1033}
1034
1039
1041{
1042 accel_terms.emplace_back(attr, coeff);
1043
1044 if (verbose && pmesh->GetMyRank() == 0)
1045 {
1046 mfem::out << "Adding Acceleration term to attributes ";
1047 for (int i = 0; i < attr.Size(); ++i)
1048 {
1049 if (attr[i] == 1)
1050 {
1051 mfem::out << i << " ";
1052 }
1053 }
1054 mfem::out << std::endl;
1055 }
1056}
1057
1062
1064{
1065 // Maximum BDF order to use at current time step
1066 // step + 1 <= order <= max_bdf_order
1067 int bdf_order = std::min(step + 1, max_bdf_order);
1068
1069 // Ratio of time step history at dt(t_{n}) - dt(t_{n-1})
1070 real_t rho1 = 0.0;
1071
1072 // Ratio of time step history at dt(t_{n-1}) - dt(t_{n-2})
1073 real_t rho2 = 0.0;
1074
1075 rho1 = dthist[0] / dthist[1];
1076
1077 if (bdf_order == 3)
1078 {
1079 rho2 = dthist[1] / dthist[2];
1080 }
1081
1082 if (step == 0 && bdf_order == 1)
1083 {
1084 bd0 = 1.0;
1085 bd1 = -1.0;
1086 bd2 = 0.0;
1087 bd3 = 0.0;
1088 ab1 = 1.0;
1089 ab2 = 0.0;
1090 ab3 = 0.0;
1091 }
1092 else if (step >= 1 && bdf_order == 2)
1093 {
1094 bd0 = (1.0 + 2.0 * rho1) / (1.0 + rho1);
1095 bd1 = -(1.0 + rho1);
1096 bd2 = pow(rho1, 2.0) / (1.0 + rho1);
1097 bd3 = 0.0;
1098 ab1 = 1.0 + rho1;
1099 ab2 = -rho1;
1100 ab3 = 0.0;
1101 }
1102 else if (step >= 2 && bdf_order == 3)
1103 {
1104 bd0 = 1.0 + rho1 / (1.0 + rho1)
1105 + (rho2 * rho1) / (1.0 + rho2 * (1 + rho1));
1106 bd1 = -1.0 - rho1 - (rho2 * rho1 * (1.0 + rho1)) / (1.0 + rho2);
1107 bd2 = pow(rho1, 2.0) * (rho2 + 1.0 / (1.0 + rho1));
1108 bd3 = -(pow(rho2, 3.0) * pow(rho1, 2.0) * (1.0 + rho1))
1109 / ((1.0 + rho2) * (1.0 + rho2 + rho2 * rho1));
1110 ab1 = ((1.0 + rho1) * (1.0 + rho2 * (1.0 + rho1))) / (1.0 + rho2);
1111 ab2 = -rho1 * (1.0 + rho2 * (1.0 + rho1));
1112 ab3 = (pow(rho2, 2.0) * rho1 * (1.0 + rho1)) / (1.0 + rho2);
1113 }
1114}
1115
1117{
1118 real_t my_rt[6], rt_max[6];
1119
1120 my_rt[0] = sw_setup.RealTime();
1121 my_rt[1] = sw_step.RealTime();
1122 my_rt[2] = sw_extrap.RealTime();
1123 my_rt[3] = sw_curlcurl.RealTime();
1124 my_rt[4] = sw_spsolve.RealTime();
1125 my_rt[5] = sw_hsolve.RealTime();
1126
1127 MPI_Reduce(my_rt, rt_max, 6, MPITypeMap<real_t>::mpi_type, MPI_MAX, 0,
1128 pmesh->GetComm());
1129
1130 if (pmesh->GetMyRank() == 0)
1131 {
1132 mfem::out << std::setw(10) << "SETUP" << std::setw(10) << "STEP"
1133 << std::setw(10) << "EXTRAP" << std::setw(10) << "CURLCURL"
1134 << std::setw(10) << "PSOLVE" << std::setw(10) << "HSOLVE"
1135 << "\n";
1136
1137 mfem::out << std::setprecision(3) << std::setw(10) << my_rt[0]
1138 << std::setw(10) << my_rt[1] << std::setw(10) << my_rt[2]
1139 << std::setw(10) << my_rt[3] << std::setw(10) << my_rt[4]
1140 << std::setw(10) << my_rt[5] << "\n";
1141
1142 mfem::out << std::setprecision(3) << std::setw(10) << " " << std::setw(10)
1143 << my_rt[1] / my_rt[1] << std::setw(10) << my_rt[2] / my_rt[1]
1144 << std::setw(10) << my_rt[3] / my_rt[1] << std::setw(10)
1145 << my_rt[4] / my_rt[1] << std::setw(10) << my_rt[5] / my_rt[1]
1146 << "\n";
1147
1148 mfem::out << std::setprecision(8);
1149 }
1150}
1151
1153{
1154 int fes_size0 = vfes->GlobalVSize();
1155 int fes_size1 = pfes->GlobalVSize();
1156
1157 if (pmesh->GetMyRank() == 0)
1158 {
1159 mfem::out << "NAVIER version: " << NAVIER_VERSION << std::endl
1160 << "MFEM version: " << MFEM_VERSION << std::endl
1161 << "MFEM GIT: " << MFEM_GIT_STRING << std::endl
1162 << "Velocity #DOFs: " << fes_size0 << std::endl
1163 << "Pressure #DOFs: " << fes_size1 << std::endl;
1164 }
1165}
1166
1168{
1169 delete FText_gfcoeff;
1170 delete g_bdr_form;
1171 delete FText_bdr_form;
1172 delete mass_lf;
1173 delete Mv_form;
1174 delete N;
1175 delete Sp_form;
1176 delete D_form;
1177 delete G_form;
1178 delete HInvPC;
1179 delete HInv;
1180 delete H_form;
1181 delete SpInv;
1182 delete MvInvPC;
1183 delete SpInvOrthoPC;
1184 delete SpInvPC;
1185 delete lor;
1186 delete f_form;
1187 delete MvInv;
1188 delete vfec;
1189 delete pfec;
1190 delete vfes;
1191 delete pfes;
1192 delete vfec_filter;
1193 delete vfes_filter;
1194}
T Max() const
Find the maximal element in the array, using the comparison operator < for class T.
Definition array.cpp:68
void SetSize(int nsize)
Change the logical size of the array, keep existing entries.
Definition array.hpp:697
int Size() const
Return the logical size of the array.
Definition array.hpp:144
void SetAssemblyLevel(AssemblyLevel assembly_level)
Set the desired assembly level.
void AddDomainIntegrator(BilinearFormIntegrator *bfi)
Adds new Domain Integrator. Assumes ownership of bfi.
Array< BilinearFormIntegrator * > * GetDBFI()
Access all the integrators added with AddDomainIntegrator().
Class for boundary integration .
Definition lininteg.hpp:211
Conjugate gradient method.
Definition solvers.hpp:513
virtual void SetOperator(const Operator &op)
Also calls SetOperator for the preconditioner if there is one.
Definition solvers.hpp:526
virtual void Mult(const Vector &b, Vector &x) const
Iterative solution of the linear system using the Conjugate Gradient method.
Definition solvers.cpp:718
Base class Coefficients that optionally depend on space and time. These are used by the BilinearFormI...
Square Operator for imposing essential boundary conditions using only the action, Mult(),...
Definition operator.hpp:895
void EliminateRHS(const Vector &x, Vector &b) const
Eliminate "essential boundary condition" values specified in x from the given right-hand side b.
Definition operator.cpp:559
Data type dense matrix using column-major storage.
Definition densemat.hpp:24
void SetSize(int s)
Change the size of the DenseMatrix to s x s.
Definition densemat.hpp:105
Class for domain integration .
Definition lininteg.hpp:109
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Definition fespace.hpp:220
ElementTransformation * GetElementTransformation(int i) const
Returns ElementTransformation for the i-th element.
Definition fespace.hpp:773
DofTransformation * GetElementVDofs(int i, Array< int > &vdofs) const
Returns indices of degrees of freedom for the i'th element. The returned indices are offsets into an ...
Definition fespace.cpp:287
virtual const FiniteElement * GetFE(int i) const
Returns pointer to the FiniteElement in the FiniteElementCollection associated with i'th element in t...
Definition fespace.cpp:3168
int GetNE() const
Returns number of elements in the mesh.
Definition fespace.hpp:740
int GetElementOrder(int i) const
Returns the order of the i'th finite element.
Definition fespace.cpp:176
int GetVSize() const
Return the number of vector dofs, i.e. GetNDofs() x GetVDim().
Definition fespace.hpp:713
int GetVDim() const
Returns vector dimension.
Definition fespace.hpp:706
Abstract class for all finite elements.
Definition fe_base.hpp:239
int GetOrder() const
Returns the order of the finite element. In the case of anisotropic orders, returns the maximum order...
Definition fe_base.hpp:333
int GetDim() const
Returns the reference space dimension for the finite element.
Definition fe_base.hpp:316
const IntegrationRule & GetNodes() const
Get a const reference to the nodes of the element.
Definition fe_base.hpp:395
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...
Geometry::Type GetGeomType() const
Returns the Geometry::Type of the reference element.
Definition fe_base.hpp:326
int GetDof() const
Returns the number of degrees of freedom in the finite element.
Definition fe_base.hpp:329
A general function coefficient.
void ProjectGridFunction(const GridFunction &src)
Project the src GridFunction to this GridFunction, both of which must be on the same mesh.
Communicator performing operations within groups defined by a GroupTopology with arbitrary-size data ...
void Reduce(T *ldata, void(*Op)(OpData< T >)) const
Reduce within each group where the master is the root.
static void Sum(OpData< T >)
Reduce operation Sum, instantiated for int and double.
void Bcast(T *ldata, int layout) const
Broadcast within each group where the master is the root.
Arbitrary order H1-conforming (continuous) finite elements.
Definition fe_coll.hpp:260
The BoomerAMG solver in hypre.
Definition hypre.hpp:1691
void SetPrintLevel(int print_level)
Definition hypre.hpp:1771
virtual void Mult(const HypreParVector &b, HypreParVector &x) const
Solve the linear system Ax=b.
Definition hypre.cpp:4047
Wrapper for hypre's ParCSR matrix class.
Definition hypre.hpp:388
Parallel smoothers in hypre.
Definition hypre.hpp:1020
Class for integration point with weight.
Definition intrules.hpp:35
Class for an integration rule - an Array of IntegrationPoint.
Definition intrules.hpp:100
int GetNPoints() const
Returns the number of the points in the integration rule.
Definition intrules.hpp:256
IntegrationPoint & IntPoint(int i)
Returns a reference to the i-th integration point.
Definition intrules.hpp:259
const IntegrationRule & Get(int GeomType, int Order)
Returns an integration rule for given GeomType and Order.
real_t GetFinalNorm() const
Returns the final residual norm after termination of the solver during the last call to Mult().
Definition solvers.hpp:275
void SetRelTol(real_t rtol)
Definition solvers.hpp:209
virtual void SetPreconditioner(Solver &pr)
This should be called before SetOperator.
Definition solvers.cpp:173
int GetNumIterations() const
Returns the number of iterations taken during the last call to Mult()
Definition solvers.hpp:260
virtual void SetPrintLevel(int print_lvl)
Legacy method to set the level of verbosity of the solver output.
Definition solvers.cpp:71
void SetMaxIter(int max_it)
Definition solvers.hpp:211
void AddDomainIntegrator(LinearFormIntegrator *lfi)
Adds new Domain Integrator. Assumes ownership of lfi.
void AddBoundaryIntegrator(LinearFormIntegrator *lfi)
Adds new Boundary Integrator. Assumes ownership of lfi.
Array< int > bdr_attributes
A list of all unique boundary attributes used by the Mesh.
Definition mesh.hpp:282
int Dimension() const
Dimension of the reference space used within the elements.
Definition mesh.hpp:1160
real_t GetElementSize(int i, int type=0)
Get the size of the i-th element relative to the perfect reference element.
Definition mesh.cpp:107
void Assemble(int skip_zeros=1)
void SetAssemblyLevel(AssemblyLevel assembly_level)
Set the desired assembly level. The default is AssemblyLevel::LEGACY.
void AddDomainIntegrator(BilinearFormIntegrator *bfi)
Adds a domain integrator. Assumes ownership of bfi.
virtual void SetIntRule(const IntegrationRule *ir)
Prescribe a fixed IntegrationRule to use (when ir != NULL) or let the integrator choose (when ir == N...
void AddDomainIntegrator(NonlinearFormIntegrator *nlfi)
Adds new Domain Integrator.
void SetAssemblyLevel(AssemblyLevel assembly_level)
Set the desired assembly level. The default is AssemblyLevel::LEGACY.
virtual void Setup()
Setup the NonlinearForm: based on the current AssemblyLevel and the current mesh, optionally,...
OpType * As() const
Return the Operator pointer statically cast to a specified OpType. Similar to the method Get().
Definition handle.hpp:104
Jacobi smoothing for a given bilinear form (no matrix necessary).
Definition solvers.hpp:313
Abstract operator.
Definition operator.hpp:25
int Height() const
Get the height (size of output) of the Operator. Synonym with NumRows().
Definition operator.hpp:66
virtual const Operator * GetRestriction() const
Restriction operator from input vectors for the operator to linear algebra (linear system) vectors....
Definition operator.hpp:143
virtual void Mult(const Vector &x, Vector &y) const =0
Operator application: y=A(x).
int Width() const
Get the width (size of input) of the Operator. Synonym with NumCols().
Definition operator.hpp:72
virtual const Operator * GetOutputProlongation() const
Prolongation operator from linear algebra (linear system) vectors, to output vectors for the operator...
Definition operator.hpp:147
virtual const Operator * GetProlongation() const
Prolongation operator from linear algebra (linear system) vectors, to input vectors for the operator....
Definition operator.hpp:139
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
Solver wrapper which orthogonalizes the input and output vector.
Definition solvers.hpp:1218
void SetSolver(Solver &s)
Set the solver used by the OrthoSolver.
Definition solvers.cpp:3436
Class for parallel bilinear form.
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 Assemble(int skip_zeros=1)
Assemble the local matrix.
virtual void RecoverFEMSolution(const Vector &X, const Vector &b, Vector &x)
virtual void Update(FiniteElementSpace *nfes=NULL)
Update the FiniteElementSpace and delete all data associated with the old one.
virtual void AssembleDiagonal(Vector &diag) const
Assemble the diagonal of the bilinear form into diag. Note that diag is a true-dof Vector.
virtual void FormSystemMatrix(const Array< int > &ess_tdof_list, OperatorHandle &A)
Form the linear system matrix A, see FormLinearSystem() for details.
Abstract parallel finite element space.
Definition pfespace.hpp:29
MPI_Comm GetComm() const
Definition pfespace.hpp:273
void GetEssentialTrueDofs(const Array< int > &bdr_attr_is_ess, Array< int > &ess_tdof_list, int component=-1) const override
HYPRE_BigInt GlobalVSize() const
Definition pfespace.hpp:283
int GetTrueVSize() const override
Return the number of local vector true dofs.
Definition pfespace.hpp:289
const SparseMatrix * GetRestrictionMatrix() const override
Get the R matrix which restricts a local dof vector to true dof vector.
Definition pfespace.hpp:389
const FiniteElement * GetFE(int i) const override
Definition pfespace.cpp:534
Class for parallel grid function.
Definition pgridfunc.hpp:33
HypreParVector * GetTrueDofs() const
Returns the true dofs in a new HypreParVector.
void SetFromTrueDofs(const Vector &tv) override
Set the GridFunction from the given true-dof vector.
void ProjectBdrCoefficient(Coefficient *coeff[], VectorCoefficient *vcoeff, const Array< int > &attr)
ParFiniteElementSpace * ParFESpace() const
void ProjectCoefficient(Coefficient &coeff) override
Project coeff Coefficient to this GridFunction. The projection computation depends on the choice of t...
void SetSpace(FiniteElementSpace *f) override
Associate a new FiniteElementSpace with the ParGridFunction.
Definition pgridfunc.cpp:96
Create and assemble a low-order refined version of a ParBilinearForm.
Definition lor.hpp:166
HypreParMatrix & GetAssembledMatrix() const
Return the assembled LOR operator as a HypreParMatrix.
Definition lor.cpp:522
Class for parallel linear form.
void ParallelAssemble(Vector &tv)
Assemble the vector on the true dofs, i.e. P^t v.
void Assemble()
Assembles the ParLinearForm i.e. sums over all domain/bdr integrators.
Class for parallel meshes.
Definition pmesh.hpp:34
MPI_Comm GetComm() const
Definition pmesh.hpp:402
int GetMyRank() const
Definition pmesh.hpp:404
Class for parallel bilinear form using different test and trial FE spaces.
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.
Parallel non-linear operator on the true dofs.
virtual void Mult(const Vector &x, Vector &y) const
Evaluate the action of the NonlinearForm.
A class container for 1D quadrature type constants.
Definition intrules.hpp:394
bool iterative_mode
If true, use the second argument of Mult() as an initial guess.
Definition operator.hpp:686
virtual void MultTranspose(const Vector &x, Vector &y) const
Multiply a vector with the transposed matrix. y = At * x.
double RealTime()
Return the number of real seconds elapsed since the stopwatch was started.
Definition tic_toc.cpp:432
void Start()
Start the stopwatch. The elapsed time is not cleared.
Definition tic_toc.cpp:411
void Stop()
Stop the stopwatch.
Definition tic_toc.cpp:422
Base class for vector Coefficients that optionally depend on time and space.
A general vector function coefficient.
Vector coefficient defined by a vector GridFunction.
Vector data type.
Definition vector.hpp:80
virtual const real_t * Read(bool on_dev=true) const
Shortcut for mfem::Read(vec.GetMemory(), vec.Size(), on_dev).
Definition vector.hpp:474
void Neg()
(*this) = -(*this)
Definition vector.cpp:300
virtual real_t * ReadWrite(bool on_dev=true)
Shortcut for mfem::ReadWrite(vec.GetMemory(), vec.Size(), on_dev).
Definition vector.hpp:490
Vector & Set(const real_t a, const Vector &x)
(*this) = a * x
Definition vector.cpp:262
int Size() const
Returns the size of the vector.
Definition vector.hpp:218
real_t Sum() const
Return the sum of the vector entries.
Definition vector.cpp:1283
void SetSize(int s)
Resize the vector to size s.
Definition vector.hpp:538
real_t * GetData() const
Return a pointer to the beginning of the Vector data.
Definition vector.hpp:227
void SetSubVectorComplement(const Array< int > &dofs, const real_t val)
Set all vector entries NOT in the dofs Array to the given val.
Definition vector.cpp:739
virtual real_t * Write(bool on_dev=true)
Shortcut for mfem::Write(vec.GetMemory(), vec.Size(), on_dev).
Definition vector.hpp:482
Vector & Add(const real_t a, const Vector &Va)
(*this) += a * Va
Definition vector.cpp:247
void Setup(real_t dt)
Initialize forms, solvers and preconditioners.
ParFiniteElementSpace * vfes_filter
real_t kin_vis
Kinematic viscosity (dimensionless).
real_t ComputeCFL(ParGridFunction &u, real_t dt)
Compute CFL.
ParLORDiscretization * lor
bool verbose
Enable/disable verbose output.
void AddVelDirichletBC(VectorCoefficient *coeff, Array< int > &attr)
Add a Dirichlet boundary condition to the velocity field.
FiniteElementCollection * vfec
Velocity finite element collection.
std::vector< PresDirichletBC_T > pres_dbcs
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.
VectorGridFunctionCoefficient * FText_gfcoeff
int order
The order of the velocity and pressure space.
bool partial_assembly
Enable/disable partial assembly of forms.
void PrintTimingData()
Print timing summary of the solving routine.
void AddPresDirichletBC(Coefficient *coeff, Array< int > &attr)
Add a Dirichlet boundary condition to the pressure field.
FiniteElementCollection * vfec_filter
ConstantCoefficient onecoeff
ParFiniteElementSpace * pfes
Pressure finite element space.
ConstantCoefficient nlcoeff
void PrintInfo()
Print information about the Navier version.
std::vector< AccelTerm_T > accel_terms
void MeanZero(ParGridFunction &v)
Remove the mean from a ParGridFunction.
void ComputeCurl2D(ParGridFunction &u, ParGridFunction &cu, bool assume_scalar=false)
Compute for .
bool numerical_integ
Enable/disable numerical integration rules of forms.
ParMixedBilinearForm * D_form
void UpdateTimestepHistory(real_t dt)
Rotate entries in the time step and solution history arrays.
ParLinearForm * mass_lf
Linear form to compute the mass matrix in various subroutines.
ConstantCoefficient H_lincoeff
ParMesh * pmesh
The parallel mesh.
ConstantCoefficient H_bdfcoeff
void AddAccelTerm(VectorCoefficient *coeff, Array< int > &attr)
Add an acceleration term to the RHS of the equation.
std::vector< VelDirichletBC_T > vel_dbcs
void Step(real_t &time, real_t dt, int cur_step, bool provisional=false)
Compute solution at the next time step t+dt.
ParMixedBilinearForm * G_form
FiniteElementCollection * pfec
Pressure finite element collection.
void SetTimeIntegrationCoefficients(int step)
Update the EXTk/BDF time integration coefficient.
void Orthogonalize(Vector &v)
Remove mean from a Vector.
ParFiniteElementSpace * vfes
Velocity finite element space.
void ComputeCurl3D(ParGridFunction &u, ParGridFunction &cu)
Compute for .
NavierSolver(ParMesh *mesh, int order, real_t kin_vis)
Initialize data structures, set FE space order and kinematic viscosity.
std::vector< real_t > dthist
int dim
Definition ex24.cpp:53
real_t b
Definition lissajous.cpp:42
real_t(const Vector &x, real_t t) ScalarFuncT
void(const Vector &x, real_t t, Vector &u) VecFuncT
real_t u(const Vector &xvec)
Definition lor_mms.hpp:22
void Mult(const Table &A, const Table &B, Table &C)
C = A * B (as boolean matrices)
Definition table.cpp:476
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
ComplexDenseMatrix * MultAtB(const ComplexDenseMatrix &A, const ComplexDenseMatrix &B)
Multiply the complex conjugate transpose of a matrix A with a matrix B. A^H*B.
float real_t
Definition config.hpp:43
std::function< real_t(const Vector &)> f(real_t mass_coeff)
Definition lor_mms.hpp:30
void forall(int N, lambda &&body)
Definition forall.hpp:754
IntegrationRules IntRules(0, Quadrature1D::GaussLegendre)
A global object with all integration rules (defined in intrules.cpp)
Definition intrules.hpp:486
void vel_dbc(const Vector &x, real_t t, Vector &u)
Fluid data.
void CopyDBFIntegrators(ParBilinearForm *src, ParBilinearForm *dst)
Helper struct to convert a C++ type to an MPI type.