13 #include "../../general/forall.hpp" 18 using namespace navier;
23 for (
int i = 0; i < bffis->
Size(); ++i)
30 : pmesh(mesh), order(order), kin_vis(kin_vis),
107 mfem::out <<
"Using Partial Assembly" << std::endl;
111 mfem::out <<
"Using Full Assembly" << std::endl;
131 nlc_nlfi->SetIntRule(&ir_ni);
172 vd_mblfi->SetIntRule(&ir_ni);
186 g_mblfi->SetIntRule(&ir_ni);
203 hmv_blfi->SetIntRule(&ir_ni);
204 hdv_blfi->SetIntRule(&ir_ni);
220 ftext_bnlfi->SetIntRule(&ir_ni);
230 gbdr_bnlfi->SetIntRule(&ir_ni);
245 vdlfi->SetIntRule(&ir_ni);
377 vel_dbc.coeff->SetTime(time + dt);
383 pres_dbc.coeff->SetTime(time + dt);
404 accel_term.coeff->SetTime(time + dt);
422 const auto ab1_ =
ab1;
423 const auto ab2_ =
ab2;
424 const auto ab3_ =
ab3;
427 d_Fext[i] = ab1_ * d_Nun[i] +
443 const double bd1idt = -
bd1 / dt;
444 const double bd2idt = -
bd2 / dt;
445 const double bd3idt = -
bd3 / dt;
446 const auto d_un =
un.
Read();
452 d_Fext[i] += bd1idt * d_un[i] +
463 const auto d_un =
un.
Read();
467 const auto ab1_ =
ab1;
468 const auto ab2_ =
ab2;
469 const auto ab3_ =
ab3;
472 d_Lext[i] = ab1_ * d_un[i] +
604 d_un_gf[i] = (1.0 - filter_alpha_) * d_un_gf[i]
605 + filter_alpha_ * d_un_filtered_gf[i];
613 mfem::out << std::setw(7) <<
"" << std::setw(3) <<
"It" << std::setw(8)
614 <<
"Resid" << std::setw(12) <<
"Reltol" 620 mfem::out << std::setw(5) <<
"MVIN " << std::setw(5) << std::fixed
622 << std::setprecision(2) << std::scientific <<
res_mvsolve 623 <<
" " << 1e-12 <<
"\n";
625 mfem::out << std::setw(5) <<
"PRES " << std::setw(5) << std::fixed
626 <<
iter_spsolve <<
" " << std::setw(3) << std::setprecision(2)
629 mfem::out << std::setw(5) <<
"HELM " << std::setw(5) << std::fixed
630 <<
iter_hsolve <<
" " << std::setw(3) << std::setprecision(2)
651 dlfi->SetIntRule(&ir_ni);
662 double integ =
mass_lf->operator()(v);
689 double loc_sum = v.
Sum();
690 double global_sum = 0.0;
691 int loc_size = v.
Size();
694 MPI_Allreduce(&loc_sum, &global_sum, 1, MPI_DOUBLE, MPI_SUM,
pfes->
GetComm());
695 MPI_Allreduce(&loc_size, &global_size, 1, MPI_INT, MPI_SUM,
pfes->
GetComm());
697 v -= global_sum /
static_cast<double>(global_size);
722 for (
int e = 0; e < fes->
GetNE(); ++e)
725 u.GetSubVector(vdofs, loc_data);
733 for (
int dof = 0; dof < elndofs; ++dof)
743 MultAtB(loc_data_mat, dshape, grad_hat);
747 Mult(grad_hat, Jinv, grad);
750 curl(0) = grad(2, 1) - grad(1, 2);
751 curl(1) = grad(0, 2) - grad(2, 0);
752 curl(2) = grad(1, 0) - grad(0, 1);
754 for (
int j = 0; j < curl.
Size(); ++j)
756 vals(elndofs * j + dof) = curl(j);
761 for (
int j = 0; j < vdofs.
Size(); j++)
765 zones_per_vdof[ldof]++;
774 gcomm.
Bcast(zones_per_vdof);
781 for (
int i = 0; i < cu.
Size(); i++)
783 const int nz = zones_per_vdof[i];
815 for (
int e = 0; e < fes->
GetNE(); ++e)
818 u.GetSubVector(vdofs, loc_data);
826 for (
int dof = 0; dof < elndofs; ++dof)
836 MultAtB(loc_data_mat, dshape, grad_hat);
840 Mult(grad_hat, Jinv, grad);
845 curl(0) = grad(0, 1);
846 curl(1) = -grad(0, 0);
851 curl(0) = grad(1, 0) - grad(0, 1);
855 for (
int j = 0; j < curl.
Size(); ++j)
857 vals(elndofs * j + dof) = curl(j);
862 for (
int j = 0; j < vdofs.
Size(); j++)
866 zones_per_vdof[ldof]++;
875 gcomm.
Bcast(zones_per_vdof);
882 for (
int i = 0; i < cu.
Size(); i++)
884 const int nz = zones_per_vdof[i];
894 ParMesh *pmesh_u =
u.ParFESpace()->GetParMesh();
906 for (
int e = 0; e < fes->
GetNE(); ++e)
913 u.GetValues(e, ir, ux, 1);
915 u.GetValues(e, ir, uy, 2);
919 u.GetValues(e, ir, uz, 3);
935 ur(i) = (ux(i) * invJ(0, 0) + uy(i) * invJ(1, 0)) * detJinv;
936 us(i) = (ux(i) * invJ(0, 1) + uy(i) * invJ(1, 1)) * detJinv;
940 ur(i) = (ux(i) * invJ(0, 0) + uy(i) * invJ(1, 0)
941 + uz(i) * invJ(2, 0))
943 us(i) = (ux(i) * invJ(0, 1) + uy(i) * invJ(1, 1)
944 + uz(i) * invJ(2, 1))
946 ut(i) = (ux(i) * invJ(0, 2) + uy(i) * invJ(1, 2)
947 + uz(i) * invJ(2, 2))
951 cflx = fabs(dt * ux(i) / hmin);
952 cfly = fabs(dt * uy(i) / hmin);
955 cflz = fabs(dt * uz(i) / hmin);
957 cflm = cflx + cfly + cflz;
958 cflmax = fmax(cflmax, cflm);
962 double cflmax_global = 0.0;
963 MPI_Allreduce(&cflmax,
970 return cflmax_global;
979 mfem::out <<
"Adding Velocity Dirichlet BC to attributes ";
980 for (
int i = 0; i < attr.
Size(); ++i)
990 for (
int i = 0; i < attr.
Size(); ++i)
993 "Duplicate boundary definition deteceted.");
1012 mfem::out <<
"Adding Pressure Dirichlet BC to attributes ";
1013 for (
int i = 0; i < attr.
Size(); ++i)
1023 for (
int i = 0; i < attr.
Size(); ++i)
1026 "Duplicate boundary definition deteceted.");
1045 mfem::out <<
"Adding Acceleration term to attributes ";
1046 for (
int i = 0; i < attr.
Size(); ++i)
1081 if (step == 0 && bdf_order == 1)
1091 else if (step >= 1 && bdf_order == 2)
1093 bd0 = (1.0 + 2.0 * rho1) / (1.0 + rho1);
1094 bd1 = -(1.0 + rho1);
1095 bd2 = pow(rho1, 2.0) / (1.0 + rho1);
1101 else if (step >= 2 && bdf_order == 3)
1103 bd0 = 1.0 + rho1 / (1.0 + rho1)
1104 + (rho2 * rho1) / (1.0 + rho2 * (1 + rho1));
1105 bd1 = -1.0 - rho1 - (rho2 * rho1 * (1.0 + rho1)) / (1.0 + rho2);
1106 bd2 = pow(rho1, 2.0) * (rho2 + 1.0 / (1.0 + rho1));
1107 bd3 = -(pow(rho2, 3.0) * pow(rho1, 2.0) * (1.0 + rho1))
1108 / ((1.0 + rho2) * (1.0 + rho2 + rho2 * rho1));
1109 ab1 = ((1.0 + rho1) * (1.0 + rho2 * (1.0 + rho1))) / (1.0 + rho2);
1110 ab2 = -rho1 * (1.0 + rho2 * (1.0 + rho1));
1111 ab3 = (pow(rho2, 2.0) * rho1 * (1.0 + rho1)) / (1.0 + rho2);
1117 double my_rt[6], rt_max[6];
1126 MPI_Reduce(my_rt, rt_max, 6, MPI_DOUBLE, MPI_MAX, 0,
pmesh->
GetComm());
1130 mfem::out << std::setw(10) <<
"SETUP" << std::setw(10) <<
"STEP" 1131 << std::setw(10) <<
"EXTRAP" << std::setw(10) <<
"CURLCURL" 1132 << std::setw(10) <<
"PSOLVE" << std::setw(10) <<
"HSOLVE" 1135 mfem::out << std::setprecision(3) << std::setw(10) << my_rt[0]
1136 << std::setw(10) << my_rt[1] << std::setw(10) << my_rt[2]
1137 << std::setw(10) << my_rt[3] << std::setw(10) << my_rt[4]
1138 << std::setw(10) << my_rt[5] <<
"\n";
1140 mfem::out << std::setprecision(3) << std::setw(10) <<
" " << std::setw(10)
1141 << my_rt[1] / my_rt[1] << std::setw(10) << my_rt[2] / my_rt[1]
1142 << std::setw(10) << my_rt[3] / my_rt[1] << std::setw(10)
1143 << my_rt[4] / my_rt[1] << std::setw(10) << my_rt[5] / my_rt[1]
1157 mfem::out <<
"NAVIER version: " << NAVIER_VERSION << std::endl
1158 <<
"MFEM version: " << MFEM_VERSION << std::endl
1159 <<
"MFEM GIT: " << MFEM_GIT_STRING << std::endl
1160 <<
"Velocity #DOFs: " << fes_size0 << std::endl
1161 <<
"Pressure #DOFs: " << fes_size1 << std::endl;
ParBilinearForm * Mv_form
Abstract class for all finite elements.
FiniteElementCollection * vfec_filter
Array< int > pres_ess_attr
void PrintTimingData()
Print timing summary of the solving routine.
void Orthogonalize(Vector &v)
Remove mean from a Vector.
Class for domain integration L(v) := (f, v)
virtual void GetEssentialTrueDofs(const Array< int > &bdr_attr_is_ess, Array< int > &ess_tdof_list, int component=-1)
Conjugate gradient method.
FiniteElementCollection * pfec
Pressure finite element collection.
int GetNPoints() const
Returns the number of the points in the integration rule.
Create and assemble a low-order refined version of a ParBilinearForm.
Class for an integration rule - an Array of IntegrationPoint.
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.
const IntegrationRule & Get(int GeomType, int Order)
Returns an integration rule for given GeomType and Order.
Base class for vector Coefficients that optionally depend on time and space.
void PrintInfo()
Print information about the Navier version.
int Dimension() const
Dimension of the reference space used within the elements.
void Bcast(T *ldata, int layout) const
Broadcast within each group where the master is the root.
void SetSize(int s)
Resize the vector to size s.
ConstantCoefficient nlcoeff
void Mult(const Table &A, const Table &B, Table &C)
C = A * B (as boolean matrices)
virtual const Operator * GetRestriction() const
Restriction operator from input vectors for the operator to linear algebra (linear system) vectors...
void(const Vector &x, double t, Vector &u) VecFuncT
int Width() const
Get the width (size of input) of the Operator. Synonym with NumCols().
virtual void Mult(const Vector &b, Vector &x) const
Operator application: y=A(x).
int Size() const
Returns the size of the vector.
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.
void vel_dbc(const Vector &x, double t, Vector &u)
Fluid data.
T Max() const
Find the maximal element in the array, using the comparison operator < for class T.
Data type dense matrix using column-major storage.
ParGridFunction un_filtered_gf
virtual const double * Read(bool on_dev=true) const
Shortcut for mfem::Read(vec.GetMemory(), vec.Size(), on_dev).
double(const Vector &x, double t) ScalarFuncT
ParLinearForm * FText_bdr_form
std::vector< PresDirichletBC_T > pres_dbcs
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).
double RealTime()
Return the number of real seconds elapsed since the stopwatch was started.
virtual void SetSpace(FiniteElementSpace *f)
Associate a new FiniteElementSpace with the ParGridFunction.
Abstract parallel finite element space.
ConstantCoefficient H_lincoeff
virtual void ProjectCoefficient(Coefficient &coeff)
Project coeff Coefficient to this GridFunction. The projection computation depends on the choice of t...
bool iterative_mode
If true, use the second argument of Mult() as an initial guess.
A class container for 1D quadrature type constants.
ParLinearForm * g_bdr_form
NavierSolver(ParMesh *mesh, int order, double kin_vis)
Initialize data structures, set FE space order and kinematic viscosity.
virtual const SparseMatrix * GetRestrictionMatrix() const
Get the R matrix which restricts a local dof vector to true dof vector.
ParMixedBilinearForm * G_form
ParLORDiscretization * lor
virtual const FiniteElement * GetFE(int i) const
Returns pointer to the FiniteElement in the FiniteElementCollection associated with i'th element in t...
IntegrationRules IntRules(0, Quadrature1D::GaussLegendre)
A global object with all integration rules (defined in intrules.cpp)
double GetFinalNorm() const
Returns the final residual norm after termination of the solver during the last call to Mult()...
void Setup(double dt)
Initialize forms, solvers and preconditioners.
void MeanZero(ParGridFunction &v)
Remove the mean from a ParGridFunction.
Array< int > vel_ess_attr
std::function< double(const Vector &)> f(double mass_coeff)
void Stop()
Stop the stopwatch.
Array< int > vel_ess_tdof
bool partial_assembly
Enable/disable partial assembly of forms.
Geometry::Type GetGeomType() const
Returns the Geometry::Type of the reference element.
The BoomerAMG solver in hypre.
IntegrationPoint & IntPoint(int i)
Returns a reference to the i-th integration point.
void ComputeCurl3D(ParGridFunction &u, ParGridFunction &cu)
Compute for .
Communicator performing operations within groups defined by a GroupTopology with arbitrary-size data ...
virtual void SetPrintLevel(int print_lvl)
Legacy method to set the level of verbosity of the solver output.
FiniteElementCollection * vfec
Velocity finite element collection.
double GetElementSize(ElementTransformation *T, int type=0)
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.
double ComputeCFL(ParGridFunction &u, double dt)
Compute CFL.
Class for boundary integration .
virtual double * Write(bool on_dev=true)
Shortcut for mfem::Write(vec.GetMemory(), vec.Size(), on_dev).
Jacobi smoothing for a given bilinear form (no matrix necessary).
int GetNumIterations() const
Returns the number of iterations taken during the last call to Mult()
void SetPrintLevel(int print_level)
void SetMaxIter(int max_it)
ConstantCoefficient onecoeff
ParGridFunction un_next_gf
double Sum() const
Return the sum of the vector entries.
IntegrationRules gll_rules
std::vector< double > dthist
ParBilinearForm * Sp_form
void SetTimeIntegrationCoefficients(int step)
Update the EXTk/BDF time integration coefficient.
void AddVelDirichletBC(VectorCoefficient *coeff, Array< int > &attr)
Add a Dirichlet boundary condition to the velocity field.
void CopyDBFIntegrators(ParBilinearForm *src, ParBilinearForm *dst)
Parallel smoothers in hypre.
void ComputeCurl2D(ParGridFunction &u, ParGridFunction &cu, bool assume_scalar=false)
Compute for .
int GetNE() const
Returns number of elements in the mesh.
virtual const Operator * GetProlongation() const
Prolongation operator from linear algebra (linear system) vectors, to input vectors for the operator...
void Start()
Start the stopwatch. The elapsed time is not cleared.
void SetSubVectorComplement(const Array< int > &dofs, const double val)
Set all vector entries NOT in the dofs Array to the given val.
double * GetData() const
Return a pointer to the beginning of the Vector data.
int GetVDim() const
Returns vector dimension.
A general vector function coefficient.
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 ...
HypreParMatrix & GetAssembledMatrix() const
Return the assembled LOR operator as a HypreParMatrix.
virtual void MultTranspose(const Vector &x, Vector &y) const
Multiply a vector with the transposed matrix. y = At * x.
Array< int > bdr_attributes
A list of all unique boundary attributes used by the Mesh.
void SetRelTol(double rtol)
void EliminateRHS(const Vector &x, Vector &b) const
Eliminate "essential boundary condition" values specified in x from the given right-hand side b...
ParFiniteElementSpace * vfes
Velocity finite element space.
int GetDim() const
Returns the reference space dimension for the finite element.
ComplexDenseMatrix * MultAtB(const ComplexDenseMatrix &A, const ComplexDenseMatrix &B)
Multiply the complex conjugate transpose of a matrix A with a matrix B. A^H*B.
virtual const FiniteElement * GetFE(int i) const
void forall(int N, lambda &&body)
void ProjectBdrCoefficient(Coefficient *coeff[], VectorCoefficient *vcoeff, Array< int > &attr)
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
HYPRE_BigInt GlobalVSize() const
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...
OutStream out(std::cout)
Global stream used by the library for standard output. Initially it uses the same std::streambuf as s...
void ProjectGridFunction(const GridFunction &src)
Project the src GridFunction to this GridFunction, both of which must be on the same mesh...
ParFiniteElementSpace * ParFESpace() const
void SetSize(int nsize)
Change the logical size of the array, keep existing entries.
int GetDof() const
Returns the number of degrees of freedom in the finite element.
int Height() const
Get the height (size of output) of the Operator. Synonym with NumRows().
Vector & Set(const double a, const Vector &x)
(*this) = a * x
std::vector< AccelTerm_T > accel_terms
OrthoSolver * SpInvOrthoPC
virtual int GetTrueVSize() const
Return the number of local vector true dofs.
bool verbose
Enable/disable verbose output.
Solver wrapper which orthogonalizes the input and output vector.
virtual double * ReadWrite(bool on_dev=true)
Shortcut for mfem::ReadWrite(vec.GetMemory(), vec.Size(), on_dev).
Vector & Add(const double a, const Vector &Va)
(*this) += a * Va
double kin_vis
Kinematic viscosity (dimensionless).
Class for integration point with weight.
Array< int > pres_ess_tdof
int GetElementOrder(int i) const
Returns the order of the i'th finite element.
std::vector< VelDirichletBC_T > vel_dbcs
OpType * As() const
Return the Operator pointer statically cast to a specified OpType. Similar to the method Get()...
VectorGridFunctionCoefficient * FText_gfcoeff
void UpdateTimestepHistory(double dt)
Rotate entries in the time step and solution history arrays.
ParFiniteElementSpace * pfes
Pressure finite element space.
ParGridFunction un_NM1_gf
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().
int Size() const
Return the logical size of the array.
ParMixedBilinearForm * D_form
virtual void SetOperator(const Operator &op)
Also calls SetOperator for the preconditioner if there is one.
A general function coefficient.
HypreParVector * GetTrueDofs() const
Returns the true dofs in a new HypreParVector.
void SetSolver(Solver &s)
Set the solver used by the OrthoSolver.
Vector coefficient defined by a vector GridFunction.
virtual void SetPreconditioner(Solver &pr)
This should be called before SetOperator.
virtual const Operator * GetOutputProlongation() const
Prolongation operator from linear algebra (linear system) vectors, to output vectors for the operator...
Arbitrary order H1-conforming (continuous) finite elements.
virtual void Mult(const HypreParVector &b, HypreParVector &x) const
Solve the linear system Ax=b.
void Step(double &time, double dt, int cur_step, bool provisional=false)
Compute solution at the next time step t+dt.
ParFiniteElementSpace * vfes_filter
ParGridFunction curlcurlu_gf
void AddPresDirichletBC(Coefficient *coeff, Array< int > &attr)
Add a Dirichlet boundary condition to the pressure field.
double u(const Vector &xvec)
Class for parallel grid function.
void SetSize(int s)
Change the size of the DenseMatrix to s x s.
Square Operator for imposing essential boundary conditions using only the action, Mult()...
Wrapper for hypre's ParCSR matrix class.
const IntegrationRule & GetNodes() const
Get a const reference to the nodes of the element.
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.
int GetOrder() const
Returns the order of the finite element. In the case of anisotropic orders, returns the maximum order...
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)
ElementTransformation * GetElementTransformation(int i) const
Returns ElementTransformation for the i-th element.