MFEM  v4.5.1
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
navier_cht.cpp
Go to the documentation of this file.
1 // Copyright (c) 2010-2022, Lawrence Livermore National Security, LLC. Produced
2 // at the Lawrence Livermore National Laboratory. All Rights reserved. See files
3 // LICENSE and NOTICE for details. LLNL-CODE-806117.
4 //
5 // This file is part of the MFEM library. For more information and source code
6 // availability visit https://mfem.org.
7 //
8 // MFEM is free software; you can redistribute it and/or modify it under the
9 // terms of the BSD-3 license. We welcome feedback and contributions, see file
10 // CONTRIBUTING.md for details.
11 //
12 // --------------------------------------------------
13 // Overlapping Grids Miniapp: Conjugate heat transfer
14 // --------------------------------------------------
15 //
16 // This example code demonstrates use of MFEM to solve different physics in
17 // different domains using overlapping grids: A solid block with its base at a
18 // fixed temperature is cooled by incoming flow. The Fluid domain models the
19 // entire domain, minus the solid block, and the incompressible Navier-Stokes
20 // equations are solved on it:
21 //
22 // ________________________________________
23 // | |
24 // | FLUID DOMAIN |
25 // | |
26 // -->inflow | ______ | --> outflow
27 // (attr=1) | | | | (attr=2)
28 // |_______________| |_________________|
29 //
30 // Inhomogeneous Dirichlet conditions are imposed at inflow (attr=1) and
31 // homogeneous Dirichlet conditions are imposed on all surface (attr=3) except
32 // the outflow (attr=2) which has Neumann boundary conditions for velocity.
33 //
34 // In contrast to the Fluid domain, the Thermal domain includes the solid block,
35 // and the advection-diffusion equation is solved on it:
36 //
37 // dT/dt + u.grad T = kappa \nabla^2 T
38 //
39 // (attr=3)
40 // ________________________________________
41 // | |
42 // | THERMAL DOMAIN |
43 // (attr=1) | kappa1 |
44 // T=0 | ______ |
45 // | |kappa2| |
46 // |_______________|______|_________________|
47 // (attr=4) (attr=2) (attr=4)
48 // T=10
49 //
50 // Inhomogeneous boundary conditions (T=10) are imposed on the base of the solid
51 // block (attr=2) and homogeneous boundary conditions are imposed at the inflow
52 // region (attr=1). All other surfaces have Neumann condition.
53 //
54 // The one-sided coupling between the two domains is via transfer of the
55 // advection velocity (u) from fluid domain to thermal domain at each time step.
56 // mpirun -np 4 navier_cht -r1 3 -r2 2 -np1 2 -np2 2
57 
58 #include "mfem.hpp"
59 #include "navier_solver.hpp"
60 #include <fstream>
61 #include <iostream>
62 
63 using namespace std;
64 using namespace mfem;
65 using namespace navier;
66 
67 struct schwarz_common
68 {
69  // common
70  double dt = 2e-2;
71  double t_final = 250*dt;
72  // fluid
73  int fluid_order = 4;
74  double fluid_kin_vis = 0.001;
75  // solid
76  int solid_order = 4;
77  int ode_solver_type = 3;
78  double alpha = 1.0e-2;
79  double kappa = 0.5;
80 } schwarz;
81 
82 // Dirichlet conditions for velocity
83 void vel_dbc(const Vector &x, double t, Vector &u);
84 // solid conductivity
85 double kappa_fun(const Vector &x);
86 // initial condition for temperature
87 double temp_init(const Vector &x);
88 
89 class ConductionOperator : public TimeDependentOperator
90 {
91 protected:
92  ParFiniteElementSpace &fespace;
93  Array<int> ess_tdof_list; // this list remains empty for pure Neumann b.c.
94 
95  mutable ParBilinearForm *M;
96  ParBilinearForm *K;
97 
98  HypreParMatrix Mmat;
99  HypreParMatrix Kmat;
100  HypreParMatrix *T; // T = M + dt K
101  double current_dt;
102 
103  mutable CGSolver M_solver; // Krylov solver for inverting the mass matrix M
104  HypreSmoother M_prec; // Preconditioner for the mass matrix M
105 
106  CGSolver T_solver; // Implicit solver for T = M + dt K
107  HypreSmoother T_prec; // Preconditioner for the implicit solver
108 
109  double alpha, kappa, udir;
110 
111  mutable Vector z; // auxiliary vector
112 
113 public:
114  ConductionOperator(ParFiniteElementSpace &f, double alpha, double kappa,
116 
117  virtual void Mult(const Vector &u, Vector &du_dt) const;
118  /** Solve the Backward-Euler equation: k = f(u + dt*k, t), for the unknown k.
119  This is the only requirement for high-order SDIRK implicit integration.*/
120  virtual void ImplicitSolve(const double dt, const Vector &u, Vector &k);
121 
122  /// Update the diffusion BilinearForm K using the given true-dof vector `u`.
123  void SetParameters(VectorGridFunctionCoefficient adv_gf_c);
124 
125  virtual ~ConductionOperator();
126 };
127 
128 void VisualizeField(socketstream &sock, const char *vishost, int visport,
129  ParGridFunction &gf, const char *title,
130  int x = 0, int y = 0, int w = 400, int h = 400,
131  bool vec = false);
132 
133 int main(int argc, char *argv[])
134 {
135  // Initialize MPI and HYPRE.
136  Mpi::Init(argc, argv);
137  int num_procs = Mpi::WorldSize();
138  int myid = Mpi::WorldRank();
139  Hypre::Init();
140 
141  // Parse command-line options.
142  int lim_meshes = 2; // should be greater than nmeshes
143  Array <const char *> mesh_file_list(lim_meshes);
144  Array <int> np_list(lim_meshes),
145  rs_levels(lim_meshes);
146  rs_levels = 0;
147  np_list = 1;
148  bool visualization = true;
149 
150  OptionsParser args(argc, argv);
151  args.AddOption(&np_list[0], "-np1", "--np1",
152  "number of MPI ranks for mesh 1");
153  args.AddOption(&np_list[1], "-np2", "--np2",
154  "number of MPI ranks for mesh 1");
155  args.AddOption(&rs_levels[0], "-r1", "--refine-serial 1",
156  "Number of times to refine the mesh 1 uniformly in serial.");
157  args.AddOption(&rs_levels[1], "-r2", "--refine-serial 2",
158  "Number of times to refine the mesh 2 uniformly in serial.");
159  args.AddOption(&visualization, "-vis", "--visualization", "-no-vis",
160  "--no-visualization",
161  "Enable or disable GLVis visualization.");
162 
163  args.Parse();
164  if (!args.Good())
165  {
166  args.PrintUsage(cout);
167  return 1;
168  }
169  if (myid == 0)
170  {
171  args.PrintOptions(cout);
172  }
173 
174  const int nmeshes = 2;
175  mesh_file_list[0] = "fluid-cht.mesh";
176  mesh_file_list[1] = "solid-cht.mesh";
177 
178  // Setup MPI communicator for each mesh
179  MPI_Comm *comml = new MPI_Comm;
180  int color = 0;
181  int npsum = 0;
182  for (int i = 0; i < nmeshes; i++)
183  {
184  npsum += np_list[i];
185  if (myid < npsum) { color = i; break; }
186  }
187  MPI_Comm_split(MPI_COMM_WORLD, color, myid, comml);
188  int myidlocal, numproclocal;
189  MPI_Comm_rank(*comml, &myidlocal);
190  MPI_Comm_size(*comml, &numproclocal);
191 
192  Mesh *mesh = new Mesh(mesh_file_list[color], 1, 1);
193  int dim = mesh->Dimension();
194  mesh->SetCurvature(color == 0 ? schwarz.fluid_order : schwarz.solid_order);
195 
196  for (int lev = 0; lev < rs_levels[color]; lev++)
197  {
198  mesh->UniformRefinement();
199  }
200 
201 
202  if (color == 0 && myidlocal == 0)
203  {
204  std::cout << "Number of elements: " << mesh->GetNE() << std::endl;
205  }
206 
207  // Setup ParMesh based on the communicator for each mesh
208  ParMesh *pmesh;
209  pmesh = new ParMesh(*comml, *mesh);
210  delete mesh;
211 
212  // Setup pointer for FESpaces, GridFunctions, and Solvers
213  H1_FECollection *fec_s = NULL; //FECollection for solid
214  ParFiniteElementSpace *fes_s = NULL; //FESpace for solid
215  ParFiniteElementSpace *adv_fes_s = NULL; //FESpace for advection in solid
216  ParGridFunction *u_gf = NULL; //Velocity solution on both meshes
217  ParGridFunction *t_gf = NULL; //Temperature solution
218  NavierSolver *flowsolver = NULL; //Fluid solver
219  ConductionOperator *coper = NULL; //Temperature solver
220  Vector t_tdof; //Temperature true-dof vector
221 
222  double t = 0,
223  dt = schwarz.dt,
224  t_final = schwarz.t_final;
225  bool last_step = false;
226 
227  // Setup flow solver on mesh for fluid
228  if (color == 0)
229  {
230  flowsolver = new NavierSolver(pmesh, schwarz.fluid_order,
231  schwarz.fluid_kin_vis);
232  flowsolver->EnablePA(true);
233  u_gf = flowsolver->GetCurrentVelocity();
234  Vector init_vel(dim);
235  init_vel = 0.;
236  VectorConstantCoefficient u_excoeff(init_vel);
237  u_gf->ProjectCoefficient(u_excoeff);
238 
239  // Dirichlet boundary conditions for fluid
240  Array<int> attr(pmesh->bdr_attributes.Max());
241  // Inlet is attribute 1.
242  attr[0] = 1;
243  // Walls is attribute 3.
244  attr[2] = 1;
245  flowsolver->AddVelDirichletBC(vel_dbc, attr);
246 
247  flowsolver->Setup(dt);
248  u_gf = flowsolver->GetCurrentVelocity();
249  }
250 
251  // Setup temperature solver for mesh on solid
252  ODESolver *ode_solver = NULL;
253  Vector vxyz;
254  if (color == 1)
255  {
256  switch (schwarz.ode_solver_type)
257  {
258  // Implicit L-stable methods
259  case 1: ode_solver = new BackwardEulerSolver; break;
260  case 2: ode_solver = new SDIRK23Solver(2); break;
261  case 3: ode_solver = new SDIRK33Solver; break;
262  // Explicit methods
263  case 11: ode_solver = new ForwardEulerSolver; break;
264  case 12: ode_solver = new RK2Solver(0.5); break; // midpoint method
265  case 13: ode_solver = new RK3SSPSolver; break;
266  case 14: ode_solver = new RK4Solver; break;
267  case 15: ode_solver = new GeneralizedAlphaSolver(0.5); break;
268  // Implicit A-stable methods (not L-stable)
269  case 22: ode_solver = new ImplicitMidpointSolver; break;
270  case 23: ode_solver = new SDIRK23Solver; break;
271  case 24: ode_solver = new SDIRK34Solver; break;
272  default:
273  std::cout << "Unknown ODE solver type: " << schwarz.ode_solver_type << '\n';
274  delete mesh;
275  return 3;
276  }
277  fec_s = new H1_FECollection(schwarz.solid_order, dim);
278  fes_s = new ParFiniteElementSpace(pmesh, fec_s);
279  adv_fes_s = new ParFiniteElementSpace(pmesh, fec_s, 2);
280  t_gf = new ParGridFunction(fes_s);
281  u_gf = new ParGridFunction(adv_fes_s);
282 
284  t_gf->ProjectCoefficient(t_0);
285  t_gf->SetTrueVector();
286  t_gf->GetTrueDofs(t_tdof);
287 
288  // Create a list of points for the interior where the gridfunction will
289  // be interpolate from the fluid mesh
290  vxyz = *pmesh->GetNodes();
291  }
292 
293  // Setup FindPointsGSLIB. Note: we set it up with MPI_COMM_WORLD to enable
294  // communication between ParMesh for solid and fluid zones.
295  OversetFindPointsGSLIB finder(MPI_COMM_WORLD);
296  finder.Setup(*pmesh, color);
297 
298  // Tag each point to be found with the same id as the mesh
299  Array<unsigned int> color_array;
300  color_array.SetSize(vxyz.Size());
301  for (int i = 0; i < color_array.Size(); i++)
302  {
303  color_array[i] = (unsigned int)color;
304  }
305  Vector interp_vals(vxyz.Size());
306 
307  // Interpolate velocity solution on both meshes. Since the velocity solution
308  // does not exist on the temperature mesh, it just passes in a dummy
309  // gridfunction that is not used in any way on the fluid mesh.
310  finder.Interpolate(vxyz, color_array, *u_gf, interp_vals);
311 
312  // Transfer the interpolated solution to solid mesh and setup a coefficient.
314  if (color == 1)
315  {
316  *u_gf = interp_vals;
317  adv_gf_c.SetGridFunction(u_gf);
318  coper = new ConductionOperator(*fes_s, schwarz.alpha, schwarz.kappa,
319  adv_gf_c);
320  coper->SetParameters(adv_gf_c);
321  }
322 
323  // Visualize the solution.
324  char vishost[] = "localhost";
325  int visport = 19916;
326  socketstream vis_sol;
327  int Ww = 350, Wh = 350; // window size
328  int Wx = color*Ww+10, Wy = 0; // window position
329  if (visualization)
330  {
331  if (color == 0)
332  {
333  VisualizeField(vis_sol, vishost, visport, *u_gf,
334  "Velocity", Wx, Wy, Ww, Wh);
335  }
336  else
337  {
338  VisualizeField(vis_sol, vishost, visport, *t_gf,
339  "Temperature", Wx, Wy, Ww, Wh);
340  }
341  }
342 
343  if (ode_solver) { ode_solver->Init(*coper); }
344 
345  for (int step = 0; !last_step; ++step)
346  {
347  if (t + dt >= t_final - dt / 2)
348  {
349  last_step = true;
350  }
351 
352  double cfl;
353  if (flowsolver)
354  {
355  flowsolver->Step(t, dt, step);
356  cfl = flowsolver->ComputeCFL(*u_gf, dt);
357  }
358  if (ode_solver)
359  {
360  ode_solver->Step(t_tdof, t, dt);
361  t_gf->SetFromTrueDofs(t_tdof);
362  }
363  finder.Interpolate(vxyz, color_array, *u_gf, interp_vals);
364  if (color == 1)
365  {
366  *u_gf = interp_vals;
367  adv_gf_c.SetGridFunction(u_gf);
368  coper->SetParameters(adv_gf_c);
369  }
370 
371  if (visualization)
372  {
373  if (color == 0)
374  {
375  VisualizeField(vis_sol, vishost, visport, *u_gf,
376  "Velocity", Wx, Wy, Ww, Wh);
377  }
378  else
379  {
380  VisualizeField(vis_sol, vishost, visport, *t_gf,
381  "Temperature", Wx, Wy, Ww, Wh);
382  }
383  }
384 
385  if (color == 0 && myidlocal == 0)
386  {
387  printf("%11s %11s %11s\n", "Time", "dt", "CFL");
388  printf("%.5E %.5E %.5E\n", t, dt,cfl);
389  fflush(stdout);
390  }
391  }
392 
393  if (flowsolver) { flowsolver->PrintTimingData(); }
394 
395  finder.FreeData();
396  delete coper;
397  delete t_gf;
398  if (color == 1) { delete u_gf; }
399  delete adv_fes_s;
400  delete fes_s;
401  delete fec_s;
402  delete ode_solver;
403  delete flowsolver;
404  delete pmesh;
405  delete comml;
406 
407  return 0;
408 }
409 
410 ConductionOperator::ConductionOperator(ParFiniteElementSpace &f, double al,
411  double kap,
413  : TimeDependentOperator(f.GetTrueVSize(), 0.0), fespace(f), M(NULL), K(NULL),
414  T(NULL), current_dt(0.0),
415  M_solver(f.GetComm()), T_solver(f.GetComm()), udir(10), z(height)
416 {
417  const double rel_tol = 1e-8;
418 
419  Array<int> ess_bdr(f.GetParMesh()->bdr_attributes.Max());
420  // Dirichlet boundary condition on inlet and isothermal section of wall.
421  ess_bdr = 0;
422  ess_bdr[0] = 1; // inlet
423  ess_bdr[1] = 1; // homogeneous isothermal section of bottom wall
424  ess_bdr[2] = 0; // top wall
425  ess_bdr[3] = 0; // inhomogeneous isothermal section of bottom wall
426  f.GetEssentialTrueDofs(ess_bdr, ess_tdof_list);
427 
428  M = new ParBilinearForm(&fespace);
429  M->AddDomainIntegrator(new MassIntegrator());
430  M->Assemble(0); // keep sparsity pattern of M and K the same
431  M->FormSystemMatrix(ess_tdof_list, Mmat);
432 
433  M_solver.iterative_mode = false;
434  M_solver.SetRelTol(rel_tol);
435  M_solver.SetAbsTol(0.0);
436  M_solver.SetMaxIter(100);
437  M_solver.SetPrintLevel(0);
438  M_prec.SetType(HypreSmoother::Jacobi);
439  M_solver.SetPreconditioner(M_prec);
440  M_solver.SetOperator(Mmat);
441 
442  alpha = al;
443  kappa = kap;
444 
445  T_solver.iterative_mode = false;
446  T_solver.SetRelTol(rel_tol);
447  T_solver.SetAbsTol(0.0);
448  T_solver.SetMaxIter(100);
449  T_solver.SetPrintLevel(0);
450  T_solver.SetPreconditioner(T_prec);
451 
452  SetParameters(adv_gf_c);
453 }
454 
455 void ConductionOperator::Mult(const Vector &u, Vector &du_dt) const
456 {
457  // Compute:
458  // du_dt = M^{-1}*-K(u)
459  // for du_dt
460 
461  Kmat.Mult(u, z);
462  z.Neg(); // z = -z
463  K->EliminateVDofsInRHS(ess_tdof_list, u, z);
464 
465  M_solver.Mult(z, du_dt);
466  du_dt.Print();
467  du_dt.SetSubVector(ess_tdof_list, 0.0);
468 }
469 
470 void ConductionOperator::ImplicitSolve(const double dt,
471  const Vector &u, Vector &du_dt)
472 {
473  // Solve the equation:
474  // du_dt = M^{-1}*[-K(u + dt*du_dt)]
475  // for du_dt
476  if (!T)
477  {
478  T = Add(1.0, Mmat, dt, Kmat);
479  current_dt = dt;
480  T_solver.SetOperator(*T);
481  }
482  MFEM_VERIFY(dt == current_dt, ""); // SDIRK methods use the same dt
483  Kmat.Mult(u, z);
484  z.Neg();
485  K->EliminateVDofsInRHS(ess_tdof_list, u, z);
486 
487  T_solver.Mult(z, du_dt);
488  du_dt.SetSubVector(ess_tdof_list, 0.0);
489 }
490 
491 void ConductionOperator::SetParameters(VectorGridFunctionCoefficient adv_gf_c)
492 {
493  ParGridFunction u_alpha_gf(&fespace);
494  FunctionCoefficient kapfuncoef(kappa_fun);
495  u_alpha_gf.ProjectCoefficient(kapfuncoef);
496 
497  delete K;
498  K = new ParBilinearForm(&fespace);
499 
500  GridFunctionCoefficient u_coeff(&u_alpha_gf);
501 
502  K->AddDomainIntegrator(new DiffusionIntegrator(u_coeff));
503  K->AddDomainIntegrator(new MixedDirectionalDerivativeIntegrator(adv_gf_c));
504  K->Assemble(0); // keep sparsity pattern of M and K the same
505  K->FormSystemMatrix(ess_tdof_list, Kmat);
506  delete T;
507  T = NULL; // re-compute T on the next ImplicitSolve
508 }
509 
510 ConductionOperator::~ConductionOperator()
511 {
512  delete T;
513  delete M;
514  delete K;
515 }
516 
517 void VisualizeField(socketstream &sock, const char *vishost, int visport,
518  ParGridFunction &gf, const char *title,
519  int x, int y, int w, int h, bool vec)
520 {
521  gf.HostRead();
522  ParMesh &pmesh = *gf.ParFESpace()->GetParMesh();
523  MPI_Comm comm = pmesh.GetComm();
524 
525  int num_procs, myid;
526  MPI_Comm_size(comm, &num_procs);
527  MPI_Comm_rank(comm, &myid);
528 
529  bool newly_opened = false;
530  int connection_failed;
531 
532  do
533  {
534  if (myid == 0)
535  {
536  if (!sock.is_open() || !sock)
537  {
538  sock.open(vishost, visport);
539  sock.precision(8);
540  newly_opened = true;
541  }
542  sock << "solution\n";
543  }
544 
545  pmesh.PrintAsOne(sock);
546  gf.SaveAsOne(sock);
547 
548  if (myid == 0 && newly_opened)
549  {
550  const char* keys = (gf.FESpace()->GetMesh()->Dimension() == 2)
551  ? "mAcRjlmm" : "mmaaAcl";
552 
553  sock << "window_title '" << title << "'\n"
554  << "window_geometry "
555  << x << " " << y << " " << w << " " << h << "\n"
556  << "keys " << keys;
557  if ( vec ) { sock << "vvv"; }
558  sock << std::endl;
559  }
560 
561  if (myid == 0)
562  {
563  connection_failed = !sock && !newly_opened;
564  }
565  MPI_Bcast(&connection_failed, 1, MPI_INT, 0, comm);
566  }
567  while (connection_failed);
568 }
569 
570 /// Fluid data
571 // Dirichlet conditions for velocity
572 void vel_dbc(const Vector &x, double t, Vector &u)
573 {
574  double xi = x(0);
575  double yi = x(1);
576 
577  u(0) = 0.;
578  u(1) = 0.;
579  if (std::fabs(xi+2.5)<1.e-5) { u(0) = 0.25*yi*(3-yi)/(1.5*1.5); }
580 }
581 
582 /// Solid data
583 // solid conductivity
584 double kappa_fun(const Vector &x)
585 {
586  return x(1) <= 1.0 && std::fabs(x(0)) < 0.5 ? 5.: 1.0;
587 }
588 
589 // initial temperature
590 double temp_init(const Vector &x)
591 {
592  double t_init = 1.0;
593  if (x(1) < 0.5)
594  {
595  t_init = 10*(std::exp(-x(1)*x(1)));
596  }
597  if (std::fabs(x(0)) >= 0.5)
598  {
599  double dx = std::fabs(x(0))-0.5;
600  t_init *= std::exp(-10*dx*dx);
601  }
602  return t_init;
603 }
void PrintTimingData()
Print timing summary of the solving routine.
int Size() const
Return the logical size of the array.
Definition: array.hpp:138
void SetSubVector(const Array< int > &dofs, const double value)
Set the entries listed in dofs to the given value.
Definition: vector.cpp:573
virtual void GetEssentialTrueDofs(const Array< int > &bdr_attr_is_ess, Array< int > &ess_tdof_list, int component=-1)
Definition: pfespace.cpp:1032
Conjugate gradient method.
Definition: solvers.hpp:465
ParMesh * GetParMesh() const
Definition: pfespace.hpp:277
Base abstract class for first order time dependent operators.
Definition: operator.hpp:289
Vector coefficient that is constant in space and time.
void Mult(const Table &A, const Table &B, Table &C)
C = A * B (as boolean matrices)
Definition: table.cpp:475
virtual void Step(Vector &x, double &t, double &dt)=0
Perform a time step from time t [in] to time t [out] based on the requested step size dt [in]...
int Wx
Coefficient defined by a GridFunction. This coefficient is mesh dependent.
int Size() const
Returns the size of the vector.
Definition: vector.hpp:200
void Interpolate(const Vector &point_pos, Array< unsigned int > &point_id, const GridFunction &field_in, Vector &field_out, int point_pos_ordering=Ordering::byNODES)
Definition: gslib.cpp:1263
Abstract class for solving systems of ODEs: dx/dt = f(x,t)
Definition: ode.hpp:22
int GetNE() const
Returns number of elements.
Definition: mesh.hpp:923
double kappa
Definition: ex24.cpp:54
Abstract parallel finite element space.
Definition: pfespace.hpp:28
virtual void ProjectCoefficient(Coefficient &coeff)
Project coeff Coefficient to this GridFunction. The projection computation depends on the choice of t...
Definition: pgridfunc.cpp:525
Backward Euler ODE solver. L-stable.
Definition: ode.hpp:391
void Setup(double dt)
Initialize forms, solvers and preconditioners.
void Add(const DenseMatrix &A, const DenseMatrix &B, double alpha, DenseMatrix &C)
C = A + alpha*B.
Definition: densemat.cpp:2282
void SetTrueVector()
Shortcut for calling GetTrueDofs() with GetTrueVector() as argument.
Definition: gridfunc.hpp:146
virtual void SetFromTrueDofs(const Vector &tv)
Set the GridFunction from the given true-dof vector.
Definition: pgridfunc.hpp:169
void Parse()
Parse the command-line options. Note that this function expects all the options provided through the ...
Definition: optparser.cpp:151
double ComputeCFL(ParGridFunction &u, double dt)
Compute CFL.
OversetFindPointsGSLIB enables use of findpts for arbitrary number of overlapping grids...
Definition: gslib.hpp:214
constexpr char vishost[]
virtual void FreeData()
Definition: gslib.cpp:265
Mesh * GetMesh() const
Returns the mesh.
Definition: fespace.hpp:441
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
Definition: mesh.cpp:9816
constexpr int visport
virtual void SetCurvature(int order, bool discont=false, int space_dim=-1, int ordering=1)
Definition: mesh.cpp:5343
T Max() const
Find the maximal element in the array, using the comparison operator &lt; for class T.
Definition: array.cpp:68
void EnablePA(bool pa)
Enable partial assembly for every operator.
void AddVelDirichletBC(VectorCoefficient *coeff, Array< int > &attr)
Add a Dirichlet boundary condition to the velocity field.
Parallel smoothers in hypre.
Definition: hypre.hpp:962
FiniteElementSpace * FESpace()
Definition: gridfunc.hpp:668
int Dimension() const
Definition: mesh.hpp:1006
void PrintUsage(std::ostream &out) const
Print the usage message.
Definition: optparser.cpp:454
virtual const double * HostRead() const
Shortcut for mfem::Read(vec.GetMemory(), vec.Size(), false).
Definition: vector.hpp:453
The classical explicit forth-order Runge-Kutta method, RK4.
Definition: ode.hpp:162
Array< int > bdr_attributes
A list of all unique boundary attributes used by the Mesh.
Definition: mesh.hpp:270
bool is_open()
True if the socketstream is open, false otherwise.
MPI_Comm GetComm() const
Definition: pmesh.hpp:351
void Print(std::ostream &out=mfem::out, int width=8) const
Prints vector to stream out.
Definition: vector.cpp:724
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 &#39;var&#39; to receive the value. Enable/disable tags are used to set the bool...
Definition: optparser.hpp:82
void SetSize(int nsize)
Change the logical size of the array, keep existing entries.
Definition: array.hpp:679
Third-order, strong stability preserving (SSP) Runge-Kutta method.
Definition: ode.hpp:149
ParGridFunction * GetCurrentVelocity()
Return a pointer to the current velocity ParGridFunction.
void SaveAsOne(const char *fname, int precision=16) const
Definition: pgridfunc.cpp:899
Implicit midpoint method. A-stable, not L-stable.
Definition: ode.hpp:404
int Wh
int dim
Definition: ex24.cpp:53
void SetGridFunction(const GridFunction *gf)
Set the grid function for this coefficient. Also sets the Vector dimension to match that of the gf...
void PrintOptions(std::ostream &out) const
Print the options.
Definition: optparser.cpp:324
HypreParVector * GetTrueDofs() const
Returns the true dofs in a new HypreParVector.
Definition: pgridfunc.cpp:143
Class for parallel bilinear form.
int open(const char hostname[], int port)
Open the socket stream on &#39;port&#39; at &#39;hostname&#39;.
RefCoord t[3]
const double alpha
Definition: ex15.cpp:369
A general function coefficient.
Vector data type.
Definition: vector.hpp:60
void GetNodes(Vector &node_coord) const
Definition: mesh.cpp:7908
Vector coefficient defined by a vector GridFunction.
void PrintAsOne(std::ostream &out=mfem::out) const
Definition: pmesh.cpp:4914
int Wy
Arbitrary order H1-conforming (continuous) finite elements.
Definition: fe_coll.hpp:220
int Ww
void Step(double &time, double dt, int cur_step, bool provisional=false)
Compute solution at the next time step t+dt.
double u(const Vector &xvec)
Definition: lor_mms.hpp:24
Class for parallel grid function.
Definition: pgridfunc.hpp:32
The classical forward Euler method.
Definition: ode.hpp:116
Wrapper for hypre&#39;s ParCSR matrix class.
Definition: hypre.hpp:343
Class for parallel meshes.
Definition: pmesh.hpp:32
Transient incompressible Navier Stokes solver in a split scheme formulation.
void Setup(Mesh &m, const int meshid, GridFunction *gfmax=NULL, const double bb_t=0.1, const double newt_tol=1.0e-12, const int npt_max=256)
Definition: gslib.cpp:1106
virtual void Init(TimeDependentOperator &f_)
Associate a TimeDependentOperator with the ODE solver.
Definition: ode.cpp:18
int main()
double f(const Vector &p)
ParFiniteElementSpace * ParFESpace() const
Definition: pgridfunc.hpp:113
bool Good() const
Return true if the command line options were parsed successfully.
Definition: optparser.hpp:150