MFEM v4.7.0
Finite element discretization library
Loading...
Searching...
No Matches
multidomain_rt.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// This miniapp is a variant of the multidomain miniapp which aims to extend
13// the demonstration given therein to PDEs involving H(div) finite elements.
14//
15// A 3D domain comprised of an outer box with a cylinder shaped inside is used.
16//
17// A pressure wave diffusion equation is described on the outer box domain
18//
19// dp/dt = ∇(κ∇•p) in outer box
20// n•p = n•p_wall on outside wall
21// ∇•p = 0 on inside (cylinder) wall
22//
23// with pressure p and coefficient κ (non-physical in this example). In this
24// context the pressure is a vector quantity equal to the force per unit area
25// exerted on an elastic material with negligible shear strength.
26//
27// A convection-diffusion equation is described inside the cylinder domain
28//
29// dp/dt = ∇(κ∇•p) - α∇(v•p) in inner cylinder
30// n•p = n•p_wall on cylinder wall (obtained from
31// pressure equation)
32// ∇•p = 0 else
33//
34// with pressure p, coefficients κ, α and prescribed velocity profile v.
35//
36// To couple the solutions of both equations, a segregated solve with one way
37// coupling approach is used. The pressure equation of the outer box is solved
38// from the timestep p_box(t) to p_box(t+dt). Then for the convection-diffusion
39// equation p_wall is set to p_box(t+dt) and the equation is solved for p(t+dt)
40// which results in a first-order one way coupling. It is important to note
41// that when using Raviart-Thomas basis functions, as in this example, only the
42// normal component of p is communicated between the two regions.
43
44#include "mfem.hpp"
45#include <fstream>
46#include <iostream>
47
48using namespace mfem;
49
50// Prescribed velocity profile for the convection-diffusion equation inside the
51// cylinder. The profile is constructed s.t. it approximates a no-slip (v=0)
52// directly at the cylinder wall boundary.
53void velocity_profile(const Vector &c, Vector &q)
54{
55 real_t A = 1.0;
56 real_t x = c(0);
57 real_t y = c(1);
58 real_t r = sqrt(pow(x, 2.0) + pow(y, 2.0));
59
60 q(0) = 0.0;
61 q(1) = 0.0;
62
63 if (std::abs(r) >= 0.25 - 1e-8)
64 {
65 q(2) = 0.0;
66 }
67 else
68 {
69 q(2) = A * exp(-(pow(x, 2.0) / 2.0 + pow(y, 2.0) / 2.0));
70 }
71}
72
73void square_xy(const Vector &p, Vector &v)
74{
75 v.SetSize(3);
76
77 v[0] = 2.0 * p[0];
78 v[1] = 2.0 * p[1];
79 v[2] = 0.0;
80}
81
82/**
83 * @brief Convection-diffusion time dependent operator
84 *
85 * dp/dt = ∇(κ∇•p) - α∇(v•p)
86 *
87 * Can also be used to create a diffusion or convection only operator by setting
88 * α or κ to zero.
89 */
90class ConvectionDiffusionTDO : public TimeDependentOperator
91{
92public:
93 /**
94 * @brief Construct a new convection-diffusion time dependent operator.
95 *
96 * @param fes The ParFiniteElementSpace the solution is defined on
97 * @param ess_tdofs All essential true dofs in the Raviart-Thomas space
98 * @param alpha The convection coefficient
99 * @param kappa The diffusion coefficient
100 */
101 ConvectionDiffusionTDO(ParFiniteElementSpace &fes,
102 Array<int> ess_tdofs,
103 real_t alpha = 1.0,
104 real_t kappa = 1.0e-1)
105 : TimeDependentOperator(fes.GetTrueVSize()),
106 Mform(&fes),
107 Kform(&fes),
108 bform(&fes),
109 ess_tdofs_(ess_tdofs),
110 M_solver(fes.GetComm())
111 {
112 d = new ConstantCoefficient(-kappa);
115
117
119 Mform.Assemble(0);
120 Mform.Finalize();
121
124 Kform.Assemble(0);
125
126 Array<int> empty;
127 Kform.FormSystemMatrix(empty, K);
128 Mform.FormSystemMatrix(ess_tdofs_, M);
129
130 bform.Assemble();
131 b = bform.ParallelAssemble();
132
133 M_solver.iterative_mode = false;
134 M_solver.SetRelTol(1e-8);
135 M_solver.SetAbsTol(0.0);
136 M_solver.SetMaxIter(100);
137 M_solver.SetPrintLevel(0);
139 M_solver.SetPreconditioner(M_prec);
140 M_solver.SetOperator(*M);
141
142 t1.SetSize(height);
143 t2.SetSize(height);
144 }
145
146 void Mult(const Vector &u, Vector &du_dt) const override
147 {
148 K->Mult(u, t1);
149 t1.Add(1.0, *b);
150 M_solver.Mult(t1, du_dt);
151 du_dt.SetSubVector(ess_tdofs_, 0.0);
152 }
153
154 ~ConvectionDiffusionTDO()
155 {
156 delete aq;
157 delete q;
158 delete d;
159 delete b;
160 }
161
162 /// Mass form
163 ParBilinearForm Mform;
164
165 /// Stiffness form. Might include diffusion, convection or both.
166 ParBilinearForm Kform;
167
168 /// Mass opeperator
170
171 /// Stiffness opeperator. Might include diffusion, convection or both.
173
174 /// RHS form
175 ParLinearForm bform;
176
177 /// RHS vector
178 Vector *b = nullptr;
179
180 /// Velocity coefficient
181 VectorCoefficient *q = nullptr;
182
183 /// alpha * Velocity coefficient
184 VectorCoefficient *aq = nullptr;
185
186 /// Diffusion coefficient
187 Coefficient *d = nullptr;
188
189 /// Essential true dof array. Relevant for eliminating boundary conditions
190 /// when using a Raviart-Thomas space.
191 Array<int> ess_tdofs_;
192
193 real_t current_dt = -1.0;
194
195 /// Mass matrix solver
196 CGSolver M_solver;
197
198 /// Mass matrix preconditioner
199 HypreSmoother M_prec;
200
201 /// Auxiliary vectors
202 mutable Vector t1, t2;
203};
204
205int main(int argc, char *argv[])
206{
207 Mpi::Init();
208 Hypre::Init();
209 int num_procs = Mpi::WorldSize();
210 int myid = Mpi::WorldRank();
211
212 int order = 1;
213 real_t t_final = 5.0;
214 real_t dt = 1.0e-5;
215 bool visualization = true;
216 int vis_steps = 10;
217
218 OptionsParser args(argc, argv);
219 args.AddOption(&order, "-o", "--order",
220 "Finite element order (polynomial degree).");
221 args.AddOption(&t_final, "-tf", "--t-final",
222 "Final time; start time is 0.");
223 args.AddOption(&dt, "-dt", "--time-step",
224 "Time step.");
225 args.AddOption(&visualization, "-vis", "--visualization", "-no-vis",
226 "--no-visualization",
227 "Enable or disable GLVis visualization.");
228 args.AddOption(&vis_steps, "-vs", "--visualization-steps",
229 "Visualize every n-th timestep.");
230 args.ParseCheck();
231
232 Mesh *serial_mesh = new Mesh("multidomain-hex.mesh");
233 ParMesh parent_mesh = ParMesh(MPI_COMM_WORLD, *serial_mesh);
234 delete serial_mesh;
235
236 parent_mesh.UniformRefinement();
237
238 RT_FECollection fec(order, parent_mesh.Dimension());
239
240 // Create the sub-domains and accompanying Finite Element spaces from
241 // corresponding attributes. This specific mesh has two domain attributes and
242 // 9 boundary attributes.
243 Array<int> cylinder_domain_attributes(1);
244 cylinder_domain_attributes[0] = 1;
245
246 auto cylinder_submesh =
247 ParSubMesh::CreateFromDomain(parent_mesh, cylinder_domain_attributes);
248
249 ParFiniteElementSpace fes_cylinder(&cylinder_submesh, &fec);
250
251 Array<int> inflow_attributes(cylinder_submesh.bdr_attributes.Max());
252 inflow_attributes = 0;
253 inflow_attributes[7] = 1;
254
255 Array<int> inner_cylinder_wall_attributes(
256 cylinder_submesh.bdr_attributes.Max());
257 inner_cylinder_wall_attributes = 0;
258 inner_cylinder_wall_attributes[8] = 1;
259
260 // For the convection-diffusion equation inside the cylinder domain, the
261 // inflow surface and outer wall are treated as Dirichlet boundary
262 // conditions.
263 Array<int> inflow_tdofs, interface_tdofs, ess_tdofs;
264 fes_cylinder.GetEssentialTrueDofs(inflow_attributes, inflow_tdofs);
265 fes_cylinder.GetEssentialTrueDofs(inner_cylinder_wall_attributes,
266 interface_tdofs);
267 ess_tdofs.Append(inflow_tdofs);
268 ess_tdofs.Append(interface_tdofs);
269 ess_tdofs.Sort();
270 ess_tdofs.Unique();
271 ConvectionDiffusionTDO cd_tdo(fes_cylinder, ess_tdofs);
272
273 ParGridFunction pressure_cylinder_gf(&fes_cylinder);
274 pressure_cylinder_gf = 0.0;
275
276 Vector pressure_cylinder;
277 pressure_cylinder_gf.GetTrueDofs(pressure_cylinder);
278
279 RK3SSPSolver cd_ode_solver;
280 cd_ode_solver.Init(cd_tdo);
281
282 Array<int> outer_domain_attributes(1);
283 outer_domain_attributes[0] = 2;
284
285 auto block_submesh = ParSubMesh::CreateFromDomain(parent_mesh,
286 outer_domain_attributes);
287
288 ParFiniteElementSpace fes_block(&block_submesh, &fec);
289
290 Array<int> block_wall_attributes(block_submesh.bdr_attributes.Max());
291 block_wall_attributes = 1;
292 block_wall_attributes[8] = 0;
293
294 Array<int> outer_cylinder_wall_attributes(
295 block_submesh.bdr_attributes.Max());
296 outer_cylinder_wall_attributes = 0;
297 outer_cylinder_wall_attributes[8] = 1;
298
299 fes_block.GetEssentialTrueDofs(block_wall_attributes, ess_tdofs);
300
301 ConvectionDiffusionTDO d_tdo(fes_block, ess_tdofs, 0.0, 1.0);
302
303 ParGridFunction pressure_block_gf(&fes_block);
304 pressure_block_gf = 0.0;
305
307 pressure_block_gf.ProjectBdrCoefficientNormal(one,
308 block_wall_attributes);
309
310 Vector pressure_block;
311 pressure_block_gf.GetTrueDofs(pressure_block);
312
313 RK3SSPSolver d_ode_solver;
314 d_ode_solver.Init(d_tdo);
315
316 Array<int> cylinder_surface_attributes(1);
317 cylinder_surface_attributes[0] = 9;
318
319 auto cylinder_surface_submesh =
320 ParSubMesh::CreateFromBoundary(parent_mesh, cylinder_surface_attributes);
321
322 char vishost[] = "localhost";
323 int visport = 19916;
324 socketstream cyl_sol_sock;
325 if (visualization)
326 {
327 cyl_sol_sock.open(vishost, visport);
328 cyl_sol_sock << "parallel " << num_procs << " " << myid << "\n";
329 cyl_sol_sock.precision(8);
330 cyl_sol_sock << "solution\n" << cylinder_submesh
331 << pressure_cylinder_gf
332 << "window_title \"Time step: " << 0 << "\""
333 << "keys cvv\n autoscale off\n valuerange 0 1.414\n"
334 << "pause\n" << std::flush;
335 }
336 socketstream block_sol_sock;
337 if (visualization)
338 {
339 block_sol_sock.open(vishost, visport);
340 block_sol_sock << "parallel " << num_procs << " " << myid << "\n";
341 block_sol_sock.precision(8);
342 block_sol_sock << "solution\n" << block_submesh << pressure_block_gf
343 << "window_title \"Time step: " << 0 << "\""
344 << "window_geometry 400 0 400 350\n"
345 << "keys cvv\n autoscale off\n valuerange 0 1.414\n"
346 << "pause\n" << std::flush;
347 }
348
349 // Create the transfer map needed in the time integration loop
350 auto pressure_block_to_cylinder_map = ParSubMesh::CreateTransferMap(
351 pressure_block_gf,
352 pressure_cylinder_gf);
353
354 real_t t = 0.0;
355 bool last_step = false;
356 for (int ti = 1; !last_step; ti++)
357 {
358 if (t + dt >= t_final - dt/2)
359 {
360 last_step = true;
361 }
362
363 // Advance the diffusion equation on the outer block to the next time step
364 d_ode_solver.Step(pressure_block, t, dt);
365 {
366 // Transfer the solution from the inner surface of the outer block to
367 // the cylinder outer surface to act as a boundary condition.
368 pressure_block_gf.SetFromTrueDofs(pressure_block);
369
370 pressure_block_to_cylinder_map.Transfer(pressure_block_gf,
371 pressure_cylinder_gf);
372
373 pressure_cylinder_gf.GetTrueDofs(pressure_cylinder);
374 }
375 // Advance the convection-diffusion equation on the outer block to the
376 // next time step
377 cd_ode_solver.Step(pressure_cylinder, t, dt);
378
379 if (last_step || (ti % vis_steps) == 0)
380 {
381 if (myid == 0)
382 {
383 out << "step " << ti << ", t = " << t << std::endl;
384 }
385
386 pressure_cylinder_gf.SetFromTrueDofs(pressure_cylinder);
387 pressure_block_gf.SetFromTrueDofs(pressure_block);
388
389 if (visualization)
390 {
391 cyl_sol_sock << "parallel " << num_procs << " " << myid << "\n";
392 cyl_sol_sock << "solution\n" << cylinder_submesh
393 << pressure_cylinder_gf
394 << "window_title \"Time step: " << ti << "\""
395 << std::flush;
396 block_sol_sock << "parallel " << num_procs << " " << myid << "\n";
397 block_sol_sock << "solution\n" << block_submesh
398 << pressure_block_gf
399 << "window_title \"Time step: " << ti << "\""
400 << std::flush;
401 }
402 }
403 }
404
405 return 0;
406}
void Sort()
Sorts the array in ascending order. This requires operator< to be defined for T.
Definition array.hpp:261
int Append(const T &el)
Append element 'el' to array, resize if necessary.
Definition array.hpp:769
void Unique()
Removes duplicities from a sorted array. This requires operator== to be defined for T.
Definition array.hpp:269
virtual void Finalize(int skip_zeros=1)
Finalizes the matrix initialization if the AssemblyLevel is AssemblyLevel::LEGACY....
void AddDomainIntegrator(BilinearFormIntegrator *bfi)
Adds new Domain Integrator. Assumes ownership of bfi.
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
Base class Coefficients that optionally depend on space and time. These are used by the BilinearFormI...
A coefficient that is constant across space and time.
for Raviart-Thomas elements
void ProjectBdrCoefficientNormal(VectorCoefficient &vcoeff, const Array< int > &bdr_attr)
Parallel smoothers in hypre.
Definition hypre.hpp:1020
void SetType(HypreSmoother::Type type, int relax_times=1)
Set the relaxation type and number of sweeps.
Definition hypre.cpp:3550
static void Init()
Initialize hypre by calling HYPRE_Init() and set default options. After calling Hypre::Init(),...
Definition hypre.hpp:74
void SetRelTol(real_t rtol)
Definition solvers.hpp:209
virtual void SetPreconditioner(Solver &pr)
This should be called before SetOperator.
Definition solvers.cpp:173
virtual void SetPrintLevel(int print_lvl)
Legacy method to set the level of verbosity of the solver output.
Definition solvers.cpp:71
void SetMaxIter(int max_it)
Definition solvers.hpp:211
void SetAbsTol(real_t atol)
Definition solvers.hpp:210
Mesh data type.
Definition mesh.hpp:56
int Dimension() const
Dimension of the reference space used within the elements.
Definition mesh.hpp:1160
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 int WorldSize()
Return the size of 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).
Pointer to an Operator of a specified type.
Definition handle.hpp:34
int height
Dimension of the output / number of rows in the matrix.
Definition operator.hpp:27
virtual void Mult(const Vector &x, Vector &y) const =0
Operator application: y=A(x).
void ParseCheck(std::ostream &out=mfem::out)
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
Class for parallel bilinear form.
void Assemble(int skip_zeros=1)
Assemble the local matrix.
virtual void FormSystemMatrix(const Array< int > &ess_tdof_list, OperatorHandle &A)
Form the linear system matrix A, see FormLinearSystem() for details.
Abstract parallel finite element space.
Definition pfespace.hpp:29
void GetEssentialTrueDofs(const Array< int > &bdr_attr_is_ess, Array< int > &ess_tdof_list, int component=-1) const override
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.
Class for parallel linear form.
void ParallelAssemble(Vector &tv)
Assemble the vector on the true dofs, i.e. P^t v.
void Assemble()
Assembles the ParLinearForm i.e. sums over all domain/bdr integrators.
Class for parallel meshes.
Definition pmesh.hpp:34
static ParTransferMap CreateTransferMap(const ParGridFunction &src, const ParGridFunction &dst)
Create a Transfer Map object.
static ParSubMesh CreateFromBoundary(const ParMesh &parent, Array< int > &boundary_attributes)
Create a surface ParSubMesh from it's parent.
Definition psubmesh.cpp:32
static ParSubMesh CreateFromDomain(const ParMesh &parent, Array< int > &domain_attributes)
Create a domain ParSubMesh from it's parent.
Definition psubmesh.cpp:26
Third-order, strong stability preserving (SSP) Runge-Kutta method.
Definition ode.hpp:151
void Init(TimeDependentOperator &f_) override
Associate a TimeDependentOperator with the ODE solver.
Definition ode.cpp:68
void Step(Vector &x, real_t &t, real_t &dt) override
Perform a time step from time t [in] to time t [out] based on the requested step size dt [in].
Definition ode.cpp:76
Arbitrary order H(div)-conforming Raviart-Thomas finite elements.
Definition fe_coll.hpp:386
Vector coefficient defined as a product of scalar and vector coefficients.
bool iterative_mode
If true, use the second argument of Mult() as an initial guess.
Definition operator.hpp:686
Base abstract class for first order time dependent operators.
Definition operator.hpp:317
Base class for vector Coefficients that optionally depend on time and space.
A general vector function coefficient.
Vector data type.
Definition vector.hpp:80
void SetSubVector(const Array< int > &dofs, const real_t value)
Set the entries listed in dofs to the given value.
Definition vector.cpp:604
void SetSize(int s)
Resize the vector to size s.
Definition vector.hpp:538
Vector & Add(const real_t a, const Vector &Va)
(*this) += a * Va
Definition vector.cpp:247
int open(const char hostname[], int port)
Open the socket stream on 'port' at 'hostname'.
const real_t alpha
Definition ex15.cpp:369
real_t kappa
Definition ex24.cpp:54
int main()
void square_xy(const Vector &p, Vector &v)
void velocity_profile(const Vector &c, Vector &q)
const int visport
real_t u(const Vector &xvec)
Definition lor_mms.hpp:22
OutStream out(std::cout)
Global stream used by the library for standard output. Initially it uses the same std::streambuf as s...
Definition globals.hpp:66
float real_t
Definition config.hpp:43
const char vishost[]
real_t p(const Vector &x, real_t t)
RefCoord t[3]