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