16using namespace navier;
21 for (
int i = 0; i < bffis->
Size(); ++i)
28 : pmesh(mesh), order(order), kin_vis(kin_vis),
105 mfem::out <<
"Using Partial Assembly" << std::endl;
109 mfem::out <<
"Using Full Assembly" << std::endl;
129 nlc_nlfi->SetIntRule(&ir_ni);
170 vd_mblfi->SetIntRule(&ir_ni);
184 g_mblfi->SetIntRule(&ir_ni);
201 hmv_blfi->SetIntRule(&ir_ni);
202 hdv_blfi->SetIntRule(&ir_ni);
218 ftext_bnlfi->SetIntRule(&ir_ni);
228 gbdr_bnlfi->SetIntRule(&ir_ni);
243 vdlfi->SetIntRule(&ir_ni);
375 vel_dbc.coeff->SetTime(time + dt);
381 pres_dbc.coeff->SetTime(time + dt);
402 accel_term.coeff->SetTime(time + dt);
420 const auto ab1_ =
ab1;
421 const auto ab2_ =
ab2;
422 const auto ab3_ =
ab3;
425 d_Fext[i] = ab1_ * d_Nun[i] +
444 const auto d_un =
un.
Read();
450 d_Fext[i] += bd1idt * d_un[i] +
461 const auto d_un =
un.
Read();
465 const auto ab1_ =
ab1;
466 const auto ab2_ =
ab2;
467 const auto ab3_ =
ab3;
470 d_Lext[i] = ab1_ * d_un[i] +
602 d_un_gf[i] = (1.0 - filter_alpha_) * d_un_gf[i]
603 + filter_alpha_ * d_un_filtered_gf[i];
611 mfem::out << std::setw(7) <<
"" << std::setw(3) <<
"It" << std::setw(8)
612 <<
"Resid" << std::setw(12) <<
"Reltol"
618 mfem::out << std::setw(5) <<
"MVIN " << std::setw(5) << std::fixed
620 << std::setprecision(2) << std::scientific <<
res_mvsolve
623 mfem::out << std::setw(5) <<
"PRES " << std::setw(5) << std::fixed
624 <<
iter_spsolve <<
" " << std::setw(3) << std::setprecision(2)
627 mfem::out << std::setw(5) <<
"HELM " << std::setw(5) << std::fixed
628 <<
iter_hsolve <<
" " << std::setw(3) << std::setprecision(2)
649 dlfi->SetIntRule(&ir_ni);
689 int loc_size = v.
Size();
694 MPI_Allreduce(&loc_size, &global_size, 1, MPI_INT, MPI_SUM,
pfes->
GetComm());
696 v -= global_sum /
static_cast<real_t>(global_size);
721 for (
int e = 0; e < fes->
GetNE(); ++e)
724 u.GetSubVector(vdofs, loc_data);
732 for (
int dof = 0; dof < elndofs; ++dof)
736 tr->SetIntPoint(&ip);
742 MultAtB(loc_data_mat, dshape, grad_hat);
746 Mult(grad_hat, Jinv, grad);
749 curl(0) = grad(2, 1) - grad(1, 2);
750 curl(1) = grad(0, 2) - grad(2, 0);
751 curl(2) = grad(1, 0) - grad(0, 1);
753 for (
int j = 0; j < curl.
Size(); ++j)
755 vals(elndofs * j + dof) = curl(j);
760 for (
int j = 0; j < vdofs.
Size(); j++)
764 zones_per_vdof[ldof]++;
773 gcomm.
Bcast(zones_per_vdof);
780 for (
int i = 0; i < cu.
Size(); i++)
782 const int nz = zones_per_vdof[i];
814 for (
int e = 0; e < fes->
GetNE(); ++e)
817 u.GetSubVector(vdofs, loc_data);
825 for (
int dof = 0; dof < elndofs; ++dof)
829 tr->SetIntPoint(&ip);
835 MultAtB(loc_data_mat, dshape, grad_hat);
839 Mult(grad_hat, Jinv, grad);
844 curl(0) = grad(0, 1);
845 curl(1) = -grad(0, 0);
850 curl(0) = grad(1, 0) - grad(0, 1);
854 for (
int j = 0; j < curl.
Size(); ++j)
856 vals(elndofs * j + dof) = curl(j);
861 for (
int j = 0; j < vdofs.
Size(); j++)
865 zones_per_vdof[ldof]++;
874 gcomm.
Bcast(zones_per_vdof);
881 for (
int i = 0; i < cu.
Size(); i++)
883 const int nz = zones_per_vdof[i];
893 ParMesh *pmesh_u =
u.ParFESpace()->GetParMesh();
905 for (
int e = 0; e < fes->
GetNE(); ++e)
912 u.GetValues(e, ir, ux, 1);
914 u.GetValues(e, ir, uy, 2);
918 u.GetValues(e, ir, uz, 3);
928 tr->SetIntPoint(&ip);
930 const real_t detJinv = 1.0 / tr->Jacobian().Det();
934 ur(i) = (ux(i) * invJ(0, 0) + uy(i) * invJ(1, 0)) * detJinv;
935 us(i) = (ux(i) * invJ(0, 1) + uy(i) * invJ(1, 1)) * detJinv;
939 ur(i) = (ux(i) * invJ(0, 0) + uy(i) * invJ(1, 0)
940 + uz(i) * invJ(2, 0))
942 us(i) = (ux(i) * invJ(0, 1) + uy(i) * invJ(1, 1)
943 + uz(i) * invJ(2, 1))
945 ut(i) = (ux(i) * invJ(0, 2) + uy(i) * invJ(1, 2)
946 + uz(i) * invJ(2, 2))
950 cflx = fabs(dt * ux(i) / hmin);
951 cfly = fabs(dt * uy(i) / hmin);
954 cflz = fabs(dt * uz(i) / hmin);
956 cflm = cflx + cfly + cflz;
957 cflmax = fmax(cflmax, cflm);
961 real_t cflmax_global = 0.0;
962 MPI_Allreduce(&cflmax,
969 return cflmax_global;
978 mfem::out <<
"Adding Velocity Dirichlet BC to attributes ";
979 for (
int i = 0; i < attr.
Size(); ++i)
989 for (
int i = 0; i < attr.
Size(); ++i)
992 "Duplicate boundary definition deteceted.");
1011 mfem::out <<
"Adding Pressure Dirichlet BC to attributes ";
1012 for (
int i = 0; i < attr.
Size(); ++i)
1022 for (
int i = 0; i < attr.
Size(); ++i)
1025 "Duplicate boundary definition deteceted.");
1044 mfem::out <<
"Adding Acceleration term to attributes ";
1045 for (
int i = 0; i < attr.
Size(); ++i)
1080 if (step == 0 && bdf_order == 1)
1090 else if (step >= 1 && bdf_order == 2)
1092 bd0 = (1.0 + 2.0 * rho1) / (1.0 + rho1);
1093 bd1 = -(1.0 + rho1);
1094 bd2 = pow(rho1, 2.0) / (1.0 + rho1);
1100 else if (step >= 2 && bdf_order == 3)
1102 bd0 = 1.0 + rho1 / (1.0 + rho1)
1103 + (rho2 * rho1) / (1.0 + rho2 * (1 + rho1));
1104 bd1 = -1.0 - rho1 - (rho2 * rho1 * (1.0 + rho1)) / (1.0 + rho2);
1105 bd2 = pow(rho1, 2.0) * (rho2 + 1.0 / (1.0 + rho1));
1106 bd3 = -(pow(rho2, 3.0) * pow(rho1, 2.0) * (1.0 + rho1))
1107 / ((1.0 + rho2) * (1.0 + rho2 + rho2 * rho1));
1108 ab1 = ((1.0 + rho1) * (1.0 + rho2 * (1.0 + rho1))) / (1.0 + rho2);
1109 ab2 = -rho1 * (1.0 + rho2 * (1.0 + rho1));
1110 ab3 = (pow(rho2, 2.0) * rho1 * (1.0 + rho1)) / (1.0 + rho2);
1116 real_t my_rt[6], rt_max[6];
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;
T Max() const
Find the maximal element in the array, using the comparison operator < for class T.
void SetSize(int nsize)
Change the logical size of the array, keep existing entries.
int Size() const
Return the logical size of the array.
Class for boundary integration .
Conjugate gradient method.
void Mult(const Vector &b, Vector &x) const override
Iterative solution of the linear system using the Conjugate Gradient method.
void SetOperator(const Operator &op) override
Set/update the solver for the given operator.
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(),...
void EliminateRHS(const Vector &x, Vector &b) const
Eliminate "essential boundary condition" values specified in x from the given right-hand side b.
Data type dense matrix using column-major storage.
void SetSize(int s)
Change the size of the DenseMatrix to s x s.
Class for domain integration .
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
ElementTransformation * GetElementTransformation(int i) const
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 ...
virtual const FiniteElement * GetFE(int i) const
Returns pointer to the FiniteElement in the FiniteElementCollection associated with i'th element in t...
int GetNE() const
Returns number of elements in the mesh.
int GetElementOrder(int i) const
Returns the order of the i'th finite element.
int GetVSize() const
Return the number of vector dofs, i.e. GetNDofs() x GetVDim().
int GetVDim() const
Returns the vector dimension of the finite element space.
Abstract class for all finite elements.
int GetOrder() const
Returns the order of the finite element. In the case of anisotropic orders, returns the maximum order...
int GetDim() const
Returns the reference space dimension for the finite element.
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...
Geometry::Type GetGeomType() const
Returns the Geometry::Type of the reference element.
int GetDof() const
Returns the number of degrees of freedom in the finite element.
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.
The BoomerAMG solver in hypre.
void SetPrintLevel(int print_level)
virtual void Mult(const HypreParVector &b, HypreParVector &x) const
Solve the linear system Ax=b.
Wrapper for hypre's ParCSR matrix class.
Parallel smoothers in hypre.
Class for integration point with weight.
Class for an integration rule - an Array of IntegrationPoint.
int GetNPoints() const
Returns the number of the points in the integration rule.
IntegrationPoint & IntPoint(int i)
Returns a reference to the i-th integration point.
const IntegrationRule & Get(int GeomType, int Order)
Returns an integration rule for given GeomType and Order.
virtual void SetIntRule(const IntegrationRule *ir)
Prescribe a fixed IntegrationRule to use, or set to null to let the integrator choose an appropriate ...
real_t GetFinalNorm() const
Returns the final residual norm after termination of the solver during the last call to Mult().
void SetRelTol(real_t rtol)
virtual void SetPreconditioner(Solver &pr)
This should be called before SetOperator.
int GetNumIterations() const
Returns the number of iterations taken during the last call to Mult()
virtual void SetPrintLevel(int print_lvl)
Legacy method to set the level of verbosity of the solver output.
void SetMaxIter(int max_it)
Array< int > bdr_attributes
A list of all unique boundary attributes used by the Mesh.
Geometry::Type GetTypicalElementGeometry() const
If the local mesh is not empty, return GetElementGeometry(0); otherwise, return a typical Geometry pr...
int Dimension() const
Dimension of the reference space used within the elements.
real_t GetElementSize(int i, int type=0)
Get the size of the i-th element relative to the perfect reference element.
OpType * As() const
Return the Operator pointer statically cast to a specified OpType. Similar to the method Get().
Jacobi smoothing for a given bilinear form (no matrix necessary).
int Height() const
Get the height (size of output) of the Operator. Synonym with NumRows().
virtual const Operator * GetRestriction() const
Restriction operator from input vectors for the operator to linear algebra (linear system) vectors....
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().
virtual const Operator * GetOutputProlongation() const
Prolongation operator from linear algebra (linear system) vectors, to output vectors for the operator...
virtual const Operator * GetProlongation() const
Prolongation operator from linear algebra (linear system) vectors, to input vectors for the operator....
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.
Solver wrapper which orthogonalizes the input and output vector.
void SetSolver(Solver &s)
Set the solver used by the OrthoSolver.
Abstract parallel finite element space.
void GetEssentialTrueDofs(const Array< int > &bdr_attr_is_ess, Array< int > &ess_tdof_list, int component=-1) const override
HYPRE_BigInt GlobalVSize() const
int GetTrueVSize() const override
Return the number of local vector true dofs.
const SparseMatrix * GetRestrictionMatrix() const override
Get the R matrix which restricts a local dof vector to true dof vector.
Class for parallel grid function.
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.
Create and assemble a low-order refined version of a ParBilinearForm.
HypreParMatrix & GetAssembledMatrix() const
Return the assembled LOR operator as a HypreParMatrix.
Class for parallel meshes.
A class container for 1D quadrature type constants.
bool iterative_mode
If true, use the second argument of Mult() as an initial guess.
void MultTranspose(const Vector &x, Vector &y) const override
Multiply a vector with the transposed matrix. y = At * x.
double RealTime()
Return the number of real seconds elapsed since the stopwatch was started.
void Start()
Start the stopwatch. The elapsed time is not cleared.
void Stop()
Stop the stopwatch.
Base class for vector Coefficients that optionally depend on time and space.
A general vector function coefficient.
Vector coefficient defined by a vector GridFunction.
virtual const real_t * Read(bool on_dev=true) const
Shortcut for mfem::Read(vec.GetMemory(), vec.Size(), on_dev).
void Neg()
(*this) = -(*this)
virtual real_t * ReadWrite(bool on_dev=true)
Shortcut for mfem::ReadWrite(vec.GetMemory(), vec.Size(), on_dev).
Vector & Set(const real_t a, const Vector &x)
(*this) = a * x
int Size() const
Returns the size of the vector.
real_t Sum() const
Return the sum of the vector entries.
void SetSize(int s)
Resize the vector to size s.
real_t * GetData() const
Return a pointer to the beginning of the Vector data.
void SetSubVectorComplement(const Array< int > &dofs, const real_t val)
Set all vector entries NOT in the dofs Array to the given val.
virtual real_t * Write(bool on_dev=true)
Shortcut for mfem::Write(vec.GetMemory(), vec.Size(), on_dev).
Vector & Add(const real_t a, const Vector &Va)
(*this) += a * Va
void Setup(real_t dt)
Initialize forms, solvers and preconditioners.
ParBilinearForm * Sp_form
ParFiniteElementSpace * vfes_filter
Array< int > pres_ess_attr
real_t kin_vis
Kinematic viscosity (dimensionless).
real_t ComputeCFL(ParGridFunction &u, real_t dt)
Compute CFL.
ParGridFunction un_filtered_gf
ParLORDiscretization * lor
bool verbose
Enable/disable verbose output.
void AddVelDirichletBC(VectorCoefficient *coeff, Array< int > &attr)
Add a Dirichlet boundary condition to the velocity field.
OrthoSolver * SpInvOrthoPC
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
ParGridFunction curlcurlu_gf
int order
The order of the velocity and pressure space.
bool partial_assembly
Enable/disable partial assembly of forms.
IntegrationRules gll_rules
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
Array< int > vel_ess_attr
ConstantCoefficient onecoeff
ParGridFunction un_next_gf
ParLinearForm * FText_bdr_form
ParFiniteElementSpace * pfes
Pressure finite element space.
ConstantCoefficient nlcoeff
ParBilinearForm * Mv_form
void PrintInfo()
Print information about the Navier version.
std::vector< AccelTerm_T > accel_terms
void MeanZero(ParGridFunction &v)
Remove the mean from a ParGridFunction.
Array< int > pres_ess_tdof
void ComputeCurl2D(ParGridFunction &u, ParGridFunction &cu, bool assume_scalar=false)
Compute for .
ParLinearForm * g_bdr_form
Array< int > vel_ess_tdof
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.
ParGridFunction un_NM1_gf
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
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)
void Mult(const Table &A, const Table &B, Table &C)
C = A * B (as boolean matrices)
OutStream out(std::cout)
Global stream used by the library for standard output. Initially it uses the same std::streambuf as s...
ComplexDenseMatrix * MultAtB(const ComplexDenseMatrix &A, const ComplexDenseMatrix &B)
Multiply the complex conjugate transpose of a matrix A with a matrix B. A^H*B.
std::function< real_t(const Vector &)> f(real_t mass_coeff)
void forall(int N, lambda &&body)
IntegrationRules IntRules(0, Quadrature1D::GaussLegendre)
A global object with all integration rules (defined in intrules.cpp)
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.