20 using namespace navier;
22 struct s_NavierContext
24 int element_subdivisions = 1;
26 double kinvis = 1.0 / 1600.0;
27 double t_final = 10 * 1e-3;
31 bool visualization =
false;
32 bool checkres =
false;
41 u(0) = sin(xi) * cos(yi) * cos(zi);
42 u(1) = -cos(xi) * sin(yi) * cos(zi);
46 class QuantitiesOfInterest
49 QuantitiesOfInterest(
ParMesh *pmesh)
54 onecoeff.constant = 1.0;
60 one_gf.ProjectCoefficient(onecoeff);
62 volume = mass_lf->operator()(one_gf);
73 for (
int i = 0; i < fes->
GetNE(); i++)
76 double intorder = 2 * fe->
GetOrder();
90 double vel2 = velx(j) * velx(j) + vely(j) * vely(j)
97 double global_integral = 0.0;
105 return 0.5 * global_integral / volume;
108 ~QuantitiesOfInterest() {
delete mass_lf; };
145 for (
int e = 0; e < fes->
GetNE(); ++e)
157 for (
int dof = 0; dof < elndofs; ++dof)
168 MultAtB(loc_data_mat, dshape, grad_hat);
172 Mult(grad_hat, Jinv, grad);
174 double q_val = 0.5 * (
sq(grad(0, 0)) +
sq(grad(1, 1)) +
sq(grad(2, 2)))
175 + grad(0, 1) * grad(1, 0) + grad(0, 2) * grad(2, 0)
176 + grad(1, 2) * grad(2, 1);
182 for (
int j = 0; j < dofs.
Size(); j++)
186 zones_per_vdof[ldof]++;
195 gcomm.
Bcast(zones_per_vdof);
202 for (
int i = 0; i < q.
Size(); i++)
204 const int nz = zones_per_vdof[i];
212 int main(
int argc,
char *argv[])
219 "--element-subdivisions",
220 "Number of 1d uniform subdivisions for each element.");
224 "Order (degree) of the finite elements.");
225 args.
AddOption(&
ctx.dt,
"-dt",
"--time-step",
"Time step.");
226 args.
AddOption(&
ctx.t_final,
"-tf",
"--final-time",
"Final time.");
232 "Enable partial assembly.");
238 "Enable numerical integration rules.");
243 "--no-visualization",
244 "Enable or disable GLVis visualization.");
251 "Enable or disable checking of the result. Returns -1 on failure.");
266 Mesh *orig_mesh =
new Mesh(
"../../data/periodic-cube.mesh");
268 ctx.element_subdivisions,
276 int nel = mesh->
GetNE();
279 mfem::out <<
"Number of elements: " << nel << std::endl;
282 auto *pmesh =
new ParMesh(MPI_COMM_WORLD, *mesh);
297 double t_final =
ctx.t_final;
298 bool last_step =
false;
300 flowsolver.
Setup(dt);
310 QuantitiesOfInterest kin_energy(pmesh);
324 double u_inf_loc = u_gf->
Normlinf();
325 double p_inf_loc = p_gf->
Normlinf();
328 double ke = kin_energy.ComputeKineticEnergy(*u_gf);
335 int nel1d = std::round(pow(nel, 1.0 / 3.0));
337 printf(
"%11s %11s %11s %11s %11s\n",
"Time",
"dt",
"u_inf",
"p_inf",
"ke");
338 printf(
"%.5E %.5E %.5E %.5E %.5E\n", t, dt, u_inf, p_inf, ke);
340 f = fopen(fname.c_str(),
"w");
341 fprintf(f,
"3D Taylor Green Vortex\n");
342 fprintf(f,
"order = %d\n",
ctx.order);
343 fprintf(f,
"grid = %d x %d x %d\n", nel1d, nel1d, nel1d);
344 fprintf(f,
"dofs per component = %d\n", ngridpts);
345 fprintf(f,
"=================================================\n");
346 fprintf(f,
" time kinetic energy\n");
347 fprintf(f,
"%20.16e %20.16e\n", t, ke);
352 for (
int step = 0; !last_step; ++step)
354 if (t + dt >= t_final - dt / 2)
359 flowsolver.
Step(t, dt, step);
361 if ((step + 1) % 100 == 0 || last_step)
370 double u_inf_loc = u_gf->
Normlinf();
371 double p_inf_loc = p_gf->
Normlinf();
374 double ke = kin_energy.ComputeKineticEnergy(*u_gf);
377 printf(
"%.5E %.5E %.5E %.5E %.5E\n", t, dt, u_inf, p_inf, ke);
378 fprintf(f,
"%20.16e %20.16e\n", t, ke);
390 double ke_expected = 1.25e-1;
391 if (fabs(ke - ke_expected) > tol)
395 mfem::out <<
"Result has a larger error than expected."
int GetNPoints() const
Returns the number of the points in the integration rule.
Abstract class for all finite elements.
void PrintTimingData()
Print timing summary of the solving routine.
int Size() const
Return the logical size of the array.
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().
int GetDim() const
Returns the reference space dimension for the finite element.
Class for an integration rule - an Array of IntegrationPoint.
Class for grid function - Vector with associated FE space.
void Bcast(T *ldata, int layout) const
Broadcast within each group where the master is the root.
void vel_tgv(const Vector &x, double t, Vector &u)
void SetCycle(int c)
Set time cycle (for time-dependent simulations)
const IntegrationRule & Get(int GeomType, int Order)
Returns an integration rule for given GeomType and Order.
void SetDataFormat(VTKFormat fmt)
A coefficient that is constant across space and time.
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.
std::string to_string(int i)
Convert an integer to an std::string.
Helper class for ParaView visualization data.
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 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...
Data type dense matrix using column-major storage.
int Size() const
Returns the size of the vector.
int GetNE() const
Returns number of elements.
Abstract parallel finite element space.
double Normlinf() const
Returns the l_infinity norm of the vector.
virtual void ProjectCoefficient(Coefficient &coeff)
int main(int argc, char *argv[])
double * GetData() const
Return a pointer to the beginning of the Vector data.
Nodes: x_i = i/(n-1), i=0,...,n-1.
void Setup(double dt)
Initialize forms, solvers and preconditioners.
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.
int GetNE() const
Returns number of elements in the mesh.
IntegrationPoint & IntPoint(int i)
Returns a reference to the i-th integration point.
void ComputeCurl3D(ParGridFunction &u, ParGridFunction &cu)
Compute for .
A simple convenience class that calls MPI_Init() at construction and MPI_Finalize() at destruction...
Communicator performing operations within groups defined by a GroupTopology with arbitrary-size data ...
virtual void RegisterField(const std::string &field_name, GridFunction *gf)
Add a grid function to the collection.
static void Sum(OpData< T >)
Reduce operation Sum, instantiated for int and double.
void Parse()
Parse the command-line options. Note that this function expects all the options provided through the ...
int Height() const
Get the height (size of output) of the Operator. Synonym with NumRows().
struct s_NavierContext ctx
ParGridFunction * GetCurrentPressure()
Return a pointer to the current pressure ParGridFunction.
const IntegrationRule & GetNodes() const
Get a const reference to the nodes of the element.
bool Root() const
Return true if WorldRank() == 0.
void EnablePA(bool pa)
Enable partial assembly for every operator.
void SetHighOrderOutput(bool high_order_output_)
void GetValues(int i, const IntegrationRule &ir, Vector &vals, int vdim=1) const
int GetVDim() const
Returns vector dimension.
FiniteElementSpace * FESpace()
void PrintUsage(std::ostream &out) const
Print the usage message.
void SetTime(double t)
Set physical time (for time-dependent simulations)
A general vector function coefficient.
HYPRE_Int GlobalVSize() const
ElementTransformation * GetElementTransformation(int i) const
Returns ElementTransformation for the i-th element.
double GlobalLpNorm(const double p, double loc_norm, MPI_Comm comm)
Compute a global Lp norm from the local Lp norms computed by each processor.
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.
void ComputeQCriterion(ParGridFunction &u, ParGridFunction &q)
void AddOption(bool *var, const char *enable_short_name, const char *enable_long_name, const char *disable_short_name, const char *disable_long_name, const char *description, bool required=false)
Add a boolean option and set 'var' to receive the value. Enable/disable tags are used to set the bool...
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)
ParGridFunction * GetCurrentVelocity()
Return a pointer to the current velocity ParGridFunction.
Class for integration point with weight.
const FiniteElement * GetFE(int i) const
Returns pointer to the FiniteElement in the FiniteElementCollection associated with i'th element in t...
void PrintOptions(std::ostream &out) const
Print the options.
double infinity()
Define a shortcut for std::numeric_limits<double>::infinity()
void SetLevelsOfDetail(int levels_of_detail_)
void MultAtB(const DenseMatrix &A, const DenseMatrix &B, DenseMatrix &AtB)
Multiply the transpose of a matrix A with a matrix B: At*B.
void GetNodes(Vector &node_coord) const
Arbitrary order H1-conforming (continuous) finite elements.
virtual void Save() override
Class for parallel grid function.
void SetSize(int s)
Change the size of the DenseMatrix to s x s.
OutStream out(std::cout)
Global stream used by the library for standard output. Initially it uses the same std::streambuf as s...
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.
Transient incompressible Navier Stokes solver in a split scheme formulation.
IntegrationRules IntRules(0, Quadrature1D::GaussLegendre)
A global object with all integration rules (defined in intrules.cpp)
double f(const Vector &p)
ParFiniteElementSpace * ParFESpace() const
bool Good() const
Return true if the command line options were parsed successfully.