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)
102 mfem::out <<
"Using Partial Assembly" << std::endl;
106 mfem::out <<
"Using Full Assembly" << std::endl;
305 vel_dbc.coeff->SetTime(time);
311 pres_dbc.coeff->SetTime(time);
341 accel_term.coeff->SetTime(time);
358 const auto ab1_ =
ab1;
359 const auto ab2_ =
ab2;
360 const auto ab3_ =
ab3;
362 d_Fext[i] = ab1_ * d_Nun[i] +
379 const double bd1idt = -
bd1 / dt;
380 const double bd2idt = -
bd2 / dt;
381 const double bd3idt = -
bd3 / dt;
382 const auto d_un =
un.
Read();
387 d_Fext[i] += bd1idt * d_un[i] +
389 bd3idt * d_unm2[i];);
397 const auto d_un =
un.
Read();
401 const auto ab1_ =
ab1;
402 const auto ab2_ =
ab2;
403 const auto ab3_ =
ab3;
405 d_Lext[i] = ab1_ * d_un[i] +
445 if (pres_dbcs.empty())
450 for (
auto &pres_dbc : pres_dbcs)
478 if (pres_dbcs.empty())
491 for (
auto &vel_dbc : vel_dbcs)
525 mfem::out << std::setw(7) <<
"" << std::setw(3) <<
"It" << std::setw(8)
526 <<
"Resid" << std::setw(12) <<
"Reltol"
532 mfem::out << std::setw(5) <<
"MVIN " << std::setw(5) << std::fixed
534 << std::setprecision(2) << std::scientific <<
res_mvsolve
535 <<
" " << 1e-12 <<
"\n";
537 mfem::out << std::setw(5) <<
"PRES " << std::setw(5) << std::fixed
538 <<
iter_spsolve <<
" " << std::setw(3) << std::setprecision(2)
541 mfem::out << std::setw(5) <<
"HELM " << std::setw(5) << std::fixed
542 <<
iter_hsolve <<
" " << std::setw(3) << std::setprecision(2)
567 double integ =
mass_lf->operator()(v);
594 double loc_sum = v.
Sum();
595 double global_sum = 0.0;
596 int loc_size = v.
Size();
599 MPI_Allreduce(&loc_sum, &global_sum, 1, MPI_DOUBLE, MPI_SUM, MPI_COMM_WORLD);
600 MPI_Allreduce(&loc_size, &global_size, 1, MPI_INT, MPI_SUM, MPI_COMM_WORLD);
602 v -= global_sum /
static_cast<double>(global_size);
627 for (
int e = 0; e < fes->
GetNE(); ++e)
638 for (
int dof = 0; dof < elndofs; ++dof)
648 MultAtB(loc_data_mat, dshape, grad_hat);
652 Mult(grad_hat, Jinv, grad);
655 curl(0) = grad(2, 1) - grad(1, 2);
656 curl(1) = grad(0, 2) - grad(2, 0);
657 curl(2) = grad(1, 0) - grad(0, 1);
659 for (
int j = 0; j < curl.
Size(); ++j)
661 vals(elndofs * j + dof) = curl(j);
666 for (
int j = 0; j < vdofs.
Size(); j++)
670 zones_per_vdof[ldof]++;
679 gcomm.
Bcast(zones_per_vdof);
686 for (
int i = 0; i < cu.
Size(); i++)
688 const int nz = zones_per_vdof[i];
720 for (
int e = 0; e < fes->
GetNE(); ++e)
731 for (
int dof = 0; dof < elndofs; ++dof)
741 MultAtB(loc_data_mat, dshape, grad_hat);
745 Mult(grad_hat, Jinv, grad);
750 curl(0) = grad(0, 1);
751 curl(1) = -grad(0, 0);
756 curl(0) = grad(1, 0) - grad(0, 1);
760 for (
int j = 0; j < curl.
Size(); ++j)
762 vals(elndofs * j + dof) = curl(j);
767 for (
int j = 0; j < vdofs.
Size(); j++)
771 zones_per_vdof[ldof]++;
780 gcomm.
Bcast(zones_per_vdof);
787 for (
int i = 0; i < cu.
Size(); i++)
789 const int nz = zones_per_vdof[i];
811 for (
int e = 0; e < fes->
GetNE(); ++e)
839 ur(i) = (ux(i) * invJ(0, 0) + uy(i) * invJ(1, 0)) * detJinv;
840 us(i) = (ux(i) * invJ(0, 1) + uy(i) * invJ(1, 1)) * detJinv;
844 ur(i) = (ux(i) * invJ(0, 0) + uy(i) * invJ(1, 0)
845 + uz(i) * invJ(2, 0))
847 us(i) = (ux(i) * invJ(0, 1) + uy(i) * invJ(1, 1)
848 + uz(i) * invJ(2, 1))
850 ut(i) = (ux(i) * invJ(0, 2) + uy(i) * invJ(1, 2)
851 + uz(i) * invJ(2, 2))
855 cflx = fabs(dt * ux(i) / hmin);
856 cfly = fabs(dt * uy(i) / hmin);
859 cflz = fabs(dt * uz(i) / hmin);
861 cflm = cflx + cfly + cflz;
862 cflmax = fmax(cflmax, cflm);
866 double cflmax_global = 0.0;
867 MPI_Allreduce(&cflmax,
874 return cflmax_global;
883 mfem::out <<
"Adding Velocity Dirichlet BC to attributes ";
884 for (
int i = 0; i < attr.
Size(); ++i)
894 for (
int i = 0; i < attr.
Size(); ++i)
897 "Duplicate boundary definition deteceted.");
916 mfem::out <<
"Adding Pressure Dirichlet BC to attributes ";
917 for (
int i = 0; i < attr.
Size(); ++i)
927 for (
int i = 0; i < attr.
Size(); ++i)
930 "Duplicate boundary definition deteceted.");
949 mfem::out <<
"Adding Acceleration term to attributes ";
950 for (
int i = 0; i < attr.
Size(); ++i)
1002 double my_rt[6], rt_max[6];
1011 MPI_Reduce(my_rt, rt_max, 6, MPI_DOUBLE, MPI_MAX, 0,
pmesh->
GetComm());
1015 mfem::out << std::setw(10) <<
"SETUP" << std::setw(10) <<
"STEP"
1016 << std::setw(10) <<
"EXTRAP" << std::setw(10) <<
"CURLCURL"
1017 << std::setw(10) <<
"PSOLVE" << std::setw(10) <<
"HSOLVE"
1020 mfem::out << std::setprecision(3) << std::setw(10) << my_rt[0]
1021 << std::setw(10) << my_rt[1] << std::setw(10) << my_rt[2]
1022 << std::setw(10) << my_rt[3] << std::setw(10) << my_rt[4]
1023 << std::setw(10) << my_rt[5] <<
"\n";
1025 mfem::out << std::setprecision(3) << std::setw(10) <<
" " << std::setw(10)
1026 << my_rt[1] / my_rt[1] << std::setw(10) << my_rt[2] / my_rt[1]
1027 << std::setw(10) << my_rt[3] / my_rt[1] << std::setw(10)
1028 << my_rt[4] / my_rt[1] << std::setw(10) << my_rt[5] / my_rt[1]
1042 mfem::out <<
"NAVIER version: " << NAVIER_VERSION << std::endl
1043 <<
"MFEM version: " << MFEM_VERSION << std::endl
1044 <<
"MFEM GIT: " << MFEM_GIT_STRING << std::endl
1045 <<
"Velocity #DOFs: " << fes_size0 << std::endl
1046 <<
"Pressure #DOFs: " << fes_size1 << std::endl;
int GetNPoints() const
Returns the number of the points in the integration rule.
ParBilinearForm * Mv_form
Abstract class for all finite elements.
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()...
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 informations about the Navier version.
virtual void Mult(const Vector &b, Vector &x) const
Operator application: y=A(x).
FiniteElementCollection * pfec_lor
void SetSize(int s)
Resize the vector to size s.
void GetElementVDofs(int i, Array< int > &vdofs) const
Returns indexes of degrees of freedom in array dofs for i'th element.
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...
Container class for integration rules.
Data type dense matrix using column-major storage.
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...
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)
bool iterative_mode
If true, use the second argument of Mult() as an initial guess.
ParLinearForm * g_bdr_form
double * Write(bool on_dev=true)
Shortcut for mfem::Write(vec.GetMemory(), vec.Size(), on_dev).
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.
ParMixedBilinearForm * G_form
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.
void Step(double &time, double dt, int cur_step)
Compute solution at the next time step t+dt.
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 .
Communicator performing operations within groups defined by a GroupTopology with arbitrary-size data ...
void SetPrintLevel(int print_lvl)
FiniteElementCollection * vfec
Velocity finite element collection.
double * ReadWrite(bool on_dev=true)
Shortcut for mfem::ReadWrite(vec.GetMemory(), vec.Size(), on_dev).
ParFiniteElementSpace * pfes_lor
static void Sum(OpData< T >)
Reduce operation Sum, instantiated for int and double.
virtual void SetFromTrueDofs(const Vector &tv)
Set the GridFunction from the given true-dof vector.
virtual void SetOperator(const Operator &op)
Set/update the solver for the given operator.
int Height() const
Get the height (size of output) of the Operator. Synonym with NumRows().
double ComputeCFL(ParGridFunction &u, double dt)
Compute CFL.
Class for boundary integration .
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
T Max() const
Find the maximal element in the array, using the comparison operator < for class T.
const IntegrationRule & GetNodes() const
Get a const reference to the nodes of the element.
ParBilinearForm * Sp_form
void SetTimeIntegrationCoefficients(int step)
Update the EXTk/BDF time integration coefficient.
ParBilinearForm * Sp_form_lor
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()
Clear the elapsed time and start the stopwatch.
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.
HYPRE_Int GlobalVSize() const
double GetElementSize(int i, int type=0)
Get the size of the i-th element relative to the perfect reference element.
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.
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 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)
int GetOrder(int i) const
Returns the order of the i'th finite element.
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.
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
const double * Read(bool on_dev=true) const
Shortcut for mfem::Read(vec.GetMemory(), vec.Size(), on_dev).
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
ParFiniteElementSpace * pfes
Pressure finite element space.
void MultTranspose(const Vector &x, Vector &y) const
Multiply a vector with the transposed matrix. y = At * x.
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 MultAtB(const DenseMatrix &A, const DenseMatrix &B, DenseMatrix &AtB)
Multiply the transpose of a matrix A with a matrix B: At*B.
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.
ParGridFunction curlcurlu_gf
void AddPresDirichletBC(Coefficient *coeff, Array< int > &attr)
Add a Dirichlet boundary condition to the pressure field.
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)
Solver wrapper which orthogonalizes the input and output vector.
double Sum() const
Return the sum of the vector entries.
ParMesh * pmesh
The parallel mesh.
void AddAccelTerm(VectorCoefficient *coeff, Array< int > &attr)
Add an acceleration term to the RHS of the equation.
int order
The order of the velocity and pressure space.
ConstantCoefficient H_bdfcoeff
double f(const Vector &p)
void Neg()
(*this) = -(*this)
ParFiniteElementSpace * ParFESpace() const