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;
426 d_Fext[i] = ab1_ * d_Nun[i] +
441 const double bd1idt = -
bd1 / dt;
442 const double bd2idt = -
bd2 / dt;
443 const double bd3idt = -
bd3 / dt;
444 const auto d_un =
un.
Read();
449 d_Fext[i] += bd1idt * d_un[i] +
451 bd3idt * d_unm2[i];);
459 const auto d_un =
un.
Read();
463 const auto ab1_ =
ab1;
464 const auto ab2_ =
ab2;
465 const auto ab3_ =
ab3;
467 d_Lext[i] = ab1_ * d_un[i] +
507 if (pres_dbcs.empty())
512 for (
auto &pres_dbc : pres_dbcs)
540 if (pres_dbcs.empty())
598 d_un_gf[i] = (1.0 - filter_alpha_) * d_un_gf[i]
599 + filter_alpha_ * d_un_filtered_gf[i];);
606 mfem::out << std::setw(7) <<
"" << std::setw(3) <<
"It" << std::setw(8)
607 <<
"Resid" << std::setw(12) <<
"Reltol"
613 mfem::out << std::setw(5) <<
"MVIN " << std::setw(5) << std::fixed
615 << std::setprecision(2) << std::scientific <<
res_mvsolve
616 <<
" " << 1e-12 <<
"\n";
618 mfem::out << std::setw(5) <<
"PRES " << std::setw(5) << std::fixed
619 <<
iter_spsolve <<
" " << std::setw(3) << std::setprecision(2)
622 mfem::out << std::setw(5) <<
"HELM " << std::setw(5) << std::fixed
623 <<
iter_hsolve <<
" " << std::setw(3) << std::setprecision(2)
644 dlfi->SetIntRule(&ir_ni);
655 double integ =
mass_lf->operator()(v);
682 double loc_sum = v.
Sum();
683 double global_sum = 0.0;
684 int loc_size = v.
Size();
687 MPI_Allreduce(&loc_sum, &global_sum, 1, MPI_DOUBLE, MPI_SUM,
pfes->
GetComm());
688 MPI_Allreduce(&loc_size, &global_size, 1, MPI_INT, MPI_SUM,
pfes->
GetComm());
690 v -= global_sum /
static_cast<double>(global_size);
715 for (
int e = 0; e < fes->
GetNE(); ++e)
726 for (
int dof = 0; dof < elndofs; ++dof)
736 MultAtB(loc_data_mat, dshape, grad_hat);
740 Mult(grad_hat, Jinv, grad);
743 curl(0) = grad(2, 1) - grad(1, 2);
744 curl(1) = grad(0, 2) - grad(2, 0);
745 curl(2) = grad(1, 0) - grad(0, 1);
747 for (
int j = 0; j < curl.
Size(); ++j)
749 vals(elndofs * j + dof) = curl(j);
754 for (
int j = 0; j < vdofs.
Size(); j++)
758 zones_per_vdof[ldof]++;
767 gcomm.
Bcast(zones_per_vdof);
774 for (
int i = 0; i < cu.
Size(); i++)
776 const int nz = zones_per_vdof[i];
808 for (
int e = 0; e < fes->
GetNE(); ++e)
819 for (
int dof = 0; dof < elndofs; ++dof)
829 MultAtB(loc_data_mat, dshape, grad_hat);
833 Mult(grad_hat, Jinv, grad);
838 curl(0) = grad(0, 1);
839 curl(1) = -grad(0, 0);
844 curl(0) = grad(1, 0) - grad(0, 1);
848 for (
int j = 0; j < curl.
Size(); ++j)
850 vals(elndofs * j + dof) = curl(j);
855 for (
int j = 0; j < vdofs.
Size(); j++)
859 zones_per_vdof[ldof]++;
868 gcomm.
Bcast(zones_per_vdof);
875 for (
int i = 0; i < cu.
Size(); i++)
877 const int nz = zones_per_vdof[i];
899 for (
int e = 0; e < fes->
GetNE(); ++e)
928 ur(i) = (ux(i) * invJ(0, 0) + uy(i) * invJ(1, 0)) * detJinv;
929 us(i) = (ux(i) * invJ(0, 1) + uy(i) * invJ(1, 1)) * detJinv;
933 ur(i) = (ux(i) * invJ(0, 0) + uy(i) * invJ(1, 0)
934 + uz(i) * invJ(2, 0))
936 us(i) = (ux(i) * invJ(0, 1) + uy(i) * invJ(1, 1)
937 + uz(i) * invJ(2, 1))
939 ut(i) = (ux(i) * invJ(0, 2) + uy(i) * invJ(1, 2)
940 + uz(i) * invJ(2, 2))
944 cflx = fabs(dt * ux(i) / hmin);
945 cfly = fabs(dt * uy(i) / hmin);
948 cflz = fabs(dt * uz(i) / hmin);
950 cflm = cflx + cfly + cflz;
951 cflmax = fmax(cflmax, cflm);
955 double cflmax_global = 0.0;
956 MPI_Allreduce(&cflmax,
963 return cflmax_global;
972 mfem::out <<
"Adding Velocity Dirichlet BC to attributes ";
973 for (
int i = 0; i < attr.
Size(); ++i)
983 for (
int i = 0; i < attr.
Size(); ++i)
986 "Duplicate boundary definition deteceted.");
1005 mfem::out <<
"Adding Pressure Dirichlet BC to attributes ";
1006 for (
int i = 0; i < attr.
Size(); ++i)
1016 for (
int i = 0; i < attr.
Size(); ++i)
1019 "Duplicate boundary definition deteceted.");
1038 mfem::out <<
"Adding Acceleration term to attributes ";
1039 for (
int i = 0; i < attr.
Size(); ++i)
1074 if (step == 0 && bdf_order == 1)
1084 else if (step >= 1 && bdf_order == 2)
1086 bd0 = (1.0 + 2.0 * rho1) / (1.0 + rho1);
1087 bd1 = -(1.0 + rho1);
1088 bd2 = pow(rho1, 2.0) / (1.0 + rho1);
1094 else if (step >= 2 && bdf_order == 3)
1096 bd0 = 1.0 + rho1 / (1.0 + rho1)
1097 + (rho2 * rho1) / (1.0 + rho2 * (1 + rho1));
1098 bd1 = -1.0 - rho1 - (rho2 * rho1 * (1.0 + rho1)) / (1.0 + rho2);
1099 bd2 = pow(rho1, 2.0) * (rho2 + 1.0 / (1.0 + rho1));
1100 bd3 = -(pow(rho2, 3.0) * pow(rho1, 2.0) * (1.0 + rho1))
1101 / ((1.0 + rho2) * (1.0 + rho2 + rho2 * rho1));
1102 ab1 = ((1.0 + rho1) * (1.0 + rho2 * (1.0 + rho1))) / (1.0 + rho2);
1103 ab2 = -rho1 * (1.0 + rho2 * (1.0 + rho1));
1104 ab3 = (pow(rho2, 2.0) * rho1 * (1.0 + rho1)) / (1.0 + rho2);
1110 double my_rt[6], rt_max[6];
1119 MPI_Reduce(my_rt, rt_max, 6, MPI_DOUBLE, MPI_MAX, 0,
pmesh->
GetComm());
1123 mfem::out << std::setw(10) <<
"SETUP" << std::setw(10) <<
"STEP"
1124 << std::setw(10) <<
"EXTRAP" << std::setw(10) <<
"CURLCURL"
1125 << std::setw(10) <<
"PSOLVE" << std::setw(10) <<
"HSOLVE"
1128 mfem::out << std::setprecision(3) << std::setw(10) << my_rt[0]
1129 << std::setw(10) << my_rt[1] << std::setw(10) << my_rt[2]
1130 << std::setw(10) << my_rt[3] << std::setw(10) << my_rt[4]
1131 << std::setw(10) << my_rt[5] <<
"\n";
1133 mfem::out << std::setprecision(3) << std::setw(10) <<
" " << std::setw(10)
1134 << my_rt[1] / my_rt[1] << std::setw(10) << my_rt[2] / my_rt[1]
1135 << std::setw(10) << my_rt[3] / my_rt[1] << std::setw(10)
1136 << my_rt[4] / my_rt[1] << std::setw(10) << my_rt[5] / my_rt[1]
1150 mfem::out <<
"NAVIER version: " << NAVIER_VERSION << std::endl
1151 <<
"MFEM version: " << MFEM_VERSION << std::endl
1152 <<
"MFEM GIT: " << MFEM_GIT_STRING << std::endl
1153 <<
"Velocity #DOFs: " << fes_size0 << std::endl
1154 <<
"Pressure #DOFs: " << fes_size1 << std::endl;
int GetNPoints() const
Returns the number of the points in the integration rule.
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.
int Size() const
Return the logical size of the array.
void Orthogonalize(Vector &v)
Remove mean from a Vector.
Class for domain integration L(v) := (f, v)
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().
virtual void GetEssentialTrueDofs(const Array< int > &bdr_attr_is_ess, Array< int > &ess_tdof_list, int component=-1)
Conjugate gradient method.
int GetDim() const
Returns the reference space dimension for the finite element.
FiniteElementCollection * pfec
Pressure finite element collection.
OpType * As() const
Return the Operator pointer statically cast to a specified OpType. Similar to the method Get()...
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.
void Bcast(T *ldata, int layout) const
Broadcast within each group where the master is the root.
int GetNumIterations() const
const IntegrationRule & Get(int GeomType, int Order)
Returns an integration rule for given GeomType and Order.
ParMesh * GetParMesh() const
Base class for vector Coefficients that optionally depend on time and space.
void PrintInfo()
Print information about the Navier version.
virtual void Mult(const Vector &b, Vector &x) const
Operator application: y=A(x).
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)
int Width() const
Get the width (size of input) of the Operator. Synonym with NumCols().
void(const Vector &x, double t, Vector &u) VecFuncT
void GetSubVector(const Array< int > &dofs, Vector &elemvect) const
Extract entries listed in dofs to the output Vector elemvect.
int GetOrder() const
Returns the order of the finite element. In the case of anisotropic orders, returns the maximum order...
void vel_dbc(const Vector &x, double t, Vector &u)
Fluid data.
Data type dense matrix using column-major storage.
ParGridFunction un_filtered_gf
int Size() const
Returns the size of the vector.
virtual const Operator * GetOutputProlongation() const
Prolongation operator from linear algebra (linear system) vectors, to output vectors for the operator...
int GetElementOrder(int i) const
Returns the order of the i'th finite element.
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).
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
double * GetData() const
Return a pointer to the beginning of the Vector data.
NavierSolver(ParMesh *mesh, int order, double kin_vis)
Initialize data structures, set FE space order and kinematic viscosity.
HypreParMatrix & GetAssembledMatrix() const
Return the assembled LOR operator as a HypreParMatrix.
ParMixedBilinearForm * G_form
ParLORDiscretization * lor
double GetFinalNorm() const
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.
Array< int > vel_ess_attr
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.
virtual const Operator * GetProlongation() const
Prolongation operator from linear algebra (linear system) vectors, to input vectors for the operator...
int GetNE() const
Returns number of elements in the mesh.
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 .
HYPRE_BigInt GlobalVSize() const
Communicator performing operations within groups defined by a GroupTopology with arbitrary-size data ...
DofTransformation * GetElementVDofs(int i, Array< int > &vdofs) const
Returns indexes of degrees of freedom in array dofs for i'th element.
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 f(const Vector &xvec)
double GetElementSize(ElementTransformation *T, int type=0)
virtual const FiniteElement * GetFE(int i) const
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.
int Height() const
Get the height (size of output) of the Operator. Synonym with NumRows().
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).
void SetPrintLevel(int print_level)
void SetMaxIter(int max_it)
virtual void Mult(const HypreParVector &b, HypreParVector &x) const
Solve the linear system Ax=b.
ConstantCoefficient onecoeff
ParGridFunction un_next_gf
T Max() const
Find the maximal element in the array, using the comparison operator < for class T.
IntegrationRules gll_rules
const IntegrationRule & GetNodes() const
Get a const reference to the nodes of the element.
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 GetValues(int i, const IntegrationRule &ir, Vector &vals, int vdim=1) const
int GetVDim() const
Returns vector dimension.
virtual const Operator * GetRestriction() const
Restriction operator from input vectors for the operator to linear algebra (linear system) vectors...
FiniteElementSpace * FESpace()
void ComputeCurl2D(ParGridFunction &u, ParGridFunction &cu, bool assume_scalar=false)
Compute for .
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.
A general vector function coefficient.
ElementTransformation * GetElementTransformation(int i) const
Returns ElementTransformation for the i-th element.
Array< int > bdr_attributes
A list of all unique boundary attributes used by the Mesh.
void SetRelTol(double rtol)
ParFiniteElementSpace * vfes
Velocity finite element space.
ComplexDenseMatrix * MultAtB(const ComplexDenseMatrix &A, const ComplexDenseMatrix &B)
Multiply the complex conjugate transpose of a matrix A with a matrix B. A^H*B.
void ProjectBdrCoefficient(Coefficient *coeff[], VectorCoefficient *vcoeff, Array< int > &attr)
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
int GetDof() const
Returns the number of degrees of freedom in the finite element.
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...
void ProjectGridFunction(const GridFunction &src)
Project the src GridFunction to this GridFunction, both of which must be on the same mesh...
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 SetSize(int nsize)
Change the logical size of the array, keep existing entries.
GroupCommunicator & GroupComm()
Return a reference to the internal GroupCommunicator (on VDofs)
Vector & Set(const double a, const Vector &x)
(*this) = a * x
std::vector< AccelTerm_T > accel_terms
OrthoSolver * SpInvOrthoPC
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
std::vector< VelDirichletBC_T > vel_dbcs
void EliminateRHS(const Vector &x, Vector &b) const
Eliminate "essential boundary condition" values specified in x from the given right-hand side b...
virtual const SparseMatrix * GetRestrictionMatrix() const
Get the R matrix which restricts a local dof vector to true dof vector.
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 MultTranspose(const Vector &x, Vector &y) const
Multiply a vector with the transposed matrix. y = At * x.
virtual const FiniteElement * GetFE(int i) const
Returns pointer to the FiniteElement in the FiniteElementCollection associated with i'th element in t...
HypreParVector * GetTrueDofs() const
Returns the true dofs in a new HypreParVector.
ParMixedBilinearForm * D_form
virtual void SetOperator(const Operator &op)
Also calls SetOperator for the preconditioner if there is one.
A general function coefficient.
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.
Arbitrary order H1-conforming (continuous) finite elements.
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()...
OutStream out(std::cout)
Global stream used by the library for standard output. Initially it uses the same std::streambuf as s...
Wrapper for hypre's ParCSR matrix class.
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.
IntegrationRules IntRules(0, Quadrature1D::GaussLegendre)
A global object with all integration rules (defined in intrules.cpp)
double Sum() const
Return the sum of the vector entries.
virtual const double * Read(bool on_dev=true) const
Shortcut for mfem::Read(vec.GetMemory(), vec.Size(), on_dev).
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)
ParFiniteElementSpace * ParFESpace() const