MFEM  v4.0
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
joule.cpp
Go to the documentation of this file.
1 // Copyright (c) 2010, Lawrence Livermore National Security, LLC. Produced at
2 // the Lawrence Livermore National Laboratory. LLNL-CODE-443211. All Rights
3 // reserved. See file COPYRIGHT for details.
4 //
5 // This file is part of the MFEM library. For more information and source code
6 // availability see http://mfem.org.
7 //
8 // MFEM is free software; you can redistribute it and/or modify it under the
9 // terms of the GNU Lesser General Public License (as published by the Free
10 // Software Foundation) version 2.1 dated February 1999.
11 //
12 // -----------------------------------------------------
13 // Joule Miniapp: Transient Magnetics and Joule Heating
14 // -----------------------------------------------------
15 //
16 // This miniapp solves a time dependent eddy current problem, resulting in Joule
17 // heating.
18 //
19 // This version has electrostatic potential, Phi, which is a source term in the
20 // EM diffusion equation. The potential itself is driven by essential BC's
21 //
22 // Div sigma Grad Phi = 0
23 // sigma E = Curl B/mu - sigma grad Phi
24 // dB/dt = - Curl E
25 // F = -k Grad T
26 // c dT/dt = -Div(F) + sigma E.E,
27 //
28 // where B is the magnetic flux, E is the electric field, T is the temperature,
29 // F is the thermal flux, sigma is electrical conductivity, mu is the magnetic
30 // permeability, and alpha is the thermal diffusivity. The geometry of the
31 // domain is assumed to be as follows:
32 //
33 // boundary attribute 3
34 // +---------------------+
35 // boundary --->| | boundary
36 // attribute 1 | | attribute 2
37 // (driven) +---------------------+
38 //
39 // The voltage BC condition is essential BC on attribute 1 (front) and 2 (rear)
40 // given by function p_bc() at bottom of this file.
41 //
42 // The E-field boundary condition specifies the essential BC (n cross E) on
43 // attribute 1 (front) and 2 (rear) given by function edot_bc at bottom of this
44 // file. The E-field can be set on attribute 3 also.
45 //
46 // The thermal boundary condition for the flux F is the natural BC on attribute 1
47 // (front) and 2 (rear). This means that dT/dt = 0 on the boundaries, and the
48 // initial T = 0.
49 //
50 // See Section 3 for how the material propertied are assigned to mesh
51 // attributes, this needs to be changed for different applications.
52 //
53 // See Section 5 for how the boundary conditions are assigned to mesh
54 // attributes, this needs to be changed for different applications.
55 //
56 // This code supports a simple version of AMR, all elements containing material
57 // attribute 1 are (optionally) refined.
58 //
59 // Compile with: make joule
60 //
61 // Sample runs:
62 //
63 // mpirun -np 8 joule -m cylinder-hex.mesh -p rod
64 // mpirun -np 8 joule -m cylinder-tet.mesh -sc 1 -amr 1 -p rod
65 // mpirun -np 8 joule -m cylinder-hex-q2.gen -s 22 -dt 0.1 -tf 240.0 -p rod
66 //
67 // Options:
68 //
69 // -m [string] the mesh file name
70 // -o [int] the order of the basis
71 // -rs [int] number of times to serially refine the mesh
72 // -rp [int] number of times to refine the mesh in parallel
73 // -s [int] time integrator 1=Backward Euler, 2=SDIRK2, 3=SDIRK3,
74 // 22=Midpoint, 23=SDIRK23, 34=SDIRK34
75 // -tf [double] the final time
76 // -dt [double] time step
77 // -mu [double] the magnetic permeability
78 // -cnd [double] the electrical conductivity
79 // -f [double] the frequency of the applied EM BC
80 // -vis [int] GLVis -vis = true -no-vis = false
81 // -vs [int] visualization step
82 // -k [string] base file name for output file
83 // -print [int] print solution (gridfunctions) to disk 0 = no, 1 = yes
84 // -amr [int] 0 = no amr, 1 = amr
85 // -sc [int] 0 = no static condensation, 1 = use static condensation
86 // -p [string] specify the problem to run, "rod", "coil", or "test"
87 //
88 //
89 // NOTE: Example meshes for this miniapp are the included cylinder/rod meshes:
90 // cylinder-hex.mesh, cylinder-tet.mesh, cylinder-hex-q2.gen,
91 // cylinder-tet-q2.gen, as well as the coil.gen mesh which can be
92 // downloaded from github.com/mfem/data (its size is 21MB). Note that the
93 // meshes with the "gen" extension require MFEM to be built with NetCDF.
94 //
95 // NOTE: This code is set up to solve two example problems, 1) a straight metal
96 // rod surrounded by air, 2) a metal rod surrounded by a metal coil all
97 // surrounded by air. To specify problem (1) use the command line options
98 // "-p rod -m cylinder-hex-q2.gen", to specify problem (2) use the
99 // command line options "-p coil -m coil.gen". Problem (1) has two
100 // materials and problem (2) has three materials, and the BC's are
101 // different.
102 //
103 // NOTE: We write out, optionally, grid functions for P, E, B, W, F, and
104 // T. These can be visualized using "glvis -np 4 -m mesh.mesh -g E",
105 // assuming we used 4 processors.
106 
107 #include "joule_solver.hpp"
108 #include <memory>
109 #include <iostream>
110 #include <fstream>
111 
112 using namespace std;
113 using namespace mfem;
114 using namespace mfem::electromagnetics;
115 
116 void display_banner(ostream & os);
117 
118 static double mj_ = 0.0;
119 static double sj_ = 0.0;
120 static double wj_ = 0.0;
121 
122 int main(int argc, char *argv[])
123 {
124  // 1. Initialize MPI.
125  MPI_Session mpi(argc, argv);
126  int myid = mpi.WorldRank();
127 
128  // print the cool banner
129  if (mpi.Root()) { display_banner(cout); }
130 
131  // 2. Parse command-line options.
132  const char *mesh_file = "cylinder-hex.mesh";
133  int ser_ref_levels = 0;
134  int par_ref_levels = 0;
135  int order = 2;
136  int ode_solver_type = 1;
137  double t_final = 100.0;
138  double dt = 0.5;
139  double mu = 1.0;
140  double sigma = 2.0*M_PI*10;
141  double Tcapacity = 1.0;
142  double Tconductivity = 0.01;
143  double freq = 1.0/60.0;
144  bool visualization = true;
145  bool visit = true;
146  int vis_steps = 1;
147  int gfprint = 0;
148  const char *basename = "Joule";
149  int amr = 0;
150  int debug = 0;
151  const char *problem = "rod";
152 
153  OptionsParser args(argc, argv);
154  args.AddOption(&mesh_file, "-m", "--mesh",
155  "Mesh file to use.");
156  args.AddOption(&ser_ref_levels, "-rs", "--refine-serial",
157  "Number of times to refine the mesh uniformly in serial.");
158  args.AddOption(&par_ref_levels, "-rp", "--refine-parallel",
159  "Number of times to refine the mesh uniformly in parallel.");
160  args.AddOption(&order, "-o", "--order",
161  "Order (degree) of the finite elements.");
162  args.AddOption(&ode_solver_type, "-s", "--ode-solver",
163  "ODE solver: 1 - Backward Euler, 2 - SDIRK2, 3 - SDIRK3\n\t."
164  "\t 22 - Mid-Point, 23 - SDIRK23, 34 - SDIRK34.");
165  args.AddOption(&t_final, "-tf", "--t-final",
166  "Final time; start time is 0.");
167  args.AddOption(&dt, "-dt", "--time-step",
168  "Time step.");
169  args.AddOption(&mu, "-mu", "--permeability",
170  "Magnetic permeability coefficient.");
171  args.AddOption(&sigma, "-cnd", "--sigma",
172  "Conductivity coefficient.");
173  args.AddOption(&freq, "-f", "--frequency",
174  "Frequency of oscillation.");
175  args.AddOption(&visualization, "-vis", "--visualization", "-no-vis",
176  "--no-visualization",
177  "Enable or disable GLVis visualization.");
178  args.AddOption(&visit, "-visit", "--visit", "-no-visit", "--no-visit",
179  "Enable or disable VisIt visualization.");
180  args.AddOption(&vis_steps, "-vs", "--visualization-steps",
181  "Visualize every n-th timestep.");
182  args.AddOption(&basename, "-k", "--outputfilename",
183  "Name of the visit dump files");
184  args.AddOption(&gfprint, "-print", "--print",
185  "Print results (grid functions) to disk.");
186  args.AddOption(&amr, "-amr", "--amr",
187  "Enable AMR");
188  args.AddOption(&STATIC_COND, "-sc", "--static-condensation",
189  "Enable static condensation");
190  args.AddOption(&debug, "-debug", "--debug",
191  "Print matrices and vectors to disk");
192  args.AddOption(&SOLVER_PRINT_LEVEL, "-hl", "--hypre-print-level",
193  "Hypre print level");
194  args.AddOption(&problem, "-p", "--problem",
195  "Name of problem to run");
196  args.Parse();
197  if (!args.Good())
198  {
199  if (mpi.Root())
200  {
201  args.PrintUsage(cout);
202  }
203  return 1;
204  }
205  if (mpi.Root())
206  {
207  args.PrintOptions(cout);
208  }
209 
210  mj_ = mu;
211  sj_ = sigma;
212  wj_ = 2.0*M_PI*freq;
213 
214  if (mpi.Root())
215  {
216  cout << "\nSkin depth sqrt(2.0/(wj*mj*sj)) = " << sqrt(2.0/(wj_*mj_*sj_))
217  << "\nSkin depth sqrt(2.0*dt/(mj*sj)) = " << sqrt(2.0*dt/(mj_*sj_))
218  << endl;
219  }
220 
221  // 3. Here material properties are assigned to mesh attributes. This code is
222  // not general, it is assumed the mesh has 3 regions each with a different
223  // integer attribute: 1, 2 or 3.
224  //
225  // The coil problem has three regions: 1) coil, 2) air, 3) the rod.
226  // The rod problem has two regions: 1) rod, 2) air.
227  //
228  // We can use the same material maps for both problems.
229 
230  std::map<int, double> sigmaMap, InvTcondMap, TcapMap, InvTcapMap;
231  double sigmaAir;
232  double TcondAir;
233  double TcapAir;
234  if (strcmp(problem,"rod")==0 || strcmp(problem,"coil")==0)
235  {
236  sigmaAir = 1.0e-6 * sigma;
237  TcondAir = 1.0e6 * Tconductivity;
238  TcapAir = 1.0 * Tcapacity;
239  }
240  else
241  {
242  cerr << "Problem " << problem << " not recognized\n";
243  mfem_error();
244  }
245 
246  if (strcmp(problem,"rod")==0 || strcmp(problem,"coil")==0)
247  {
248 
249  sigmaMap.insert(pair<int, double>(1, sigma));
250  sigmaMap.insert(pair<int, double>(2, sigmaAir));
251  sigmaMap.insert(pair<int, double>(3, sigmaAir));
252 
253  InvTcondMap.insert(pair<int, double>(1, 1.0/Tconductivity));
254  InvTcondMap.insert(pair<int, double>(2, 1.0/TcondAir));
255  InvTcondMap.insert(pair<int, double>(3, 1.0/TcondAir));
256 
257  TcapMap.insert(pair<int, double>(1, Tcapacity));
258  TcapMap.insert(pair<int, double>(2, TcapAir));
259  TcapMap.insert(pair<int, double>(3, TcapAir));
260 
261  InvTcapMap.insert(pair<int, double>(1, 1.0/Tcapacity));
262  InvTcapMap.insert(pair<int, double>(2, 1.0/TcapAir));
263  InvTcapMap.insert(pair<int, double>(3, 1.0/TcapAir));
264  }
265  else
266  {
267  cerr << "Problem " << problem << " not recognized\n";
268  mfem_error();
269  }
270 
271  // 4. Read the serial mesh from the given mesh file on all processors. We can
272  // handle triangular, quadrilateral, tetrahedral and hexahedral meshes
273  // with the same code.
274  Mesh *mesh;
275  mesh = new Mesh(mesh_file, 1, 1);
276  int dim = mesh->Dimension();
277 
278  // 5. Assign the boundary conditions
279  Array<int> ess_bdr(mesh->bdr_attributes.Max());
280  Array<int> thermal_ess_bdr(mesh->bdr_attributes.Max());
281  Array<int> poisson_ess_bdr(mesh->bdr_attributes.Max());
282  if (strcmp(problem,"coil")==0)
283  {
284  // BEGIN CODE FOR THE COIL PROBLEM
285  // For the coil in a box problem we have surfaces 1) coil end (+),
286  // 2) coil end (-), 3) five sides of box, 4) side of box with coil BC
287 
288  ess_bdr = 0;
289  ess_bdr[0] = 1; // boundary attribute 4 (index 3) is fixed
290  ess_bdr[1] = 1; // boundary attribute 4 (index 3) is fixed
291  ess_bdr[2] = 1; // boundary attribute 4 (index 3) is fixed
292  ess_bdr[3] = 1; // boundary attribute 4 (index 3) is fixed
293 
294  // Same as above, but this is for the thermal operator for HDiv
295  // formulation the essential BC is the flux
296 
297  thermal_ess_bdr = 0;
298  thermal_ess_bdr[2] = 1; // boundary attribute 4 (index 3) is fixed
299 
300  // Same as above, but this is for the poisson eq for H1 formulation the
301  // essential BC is the value of Phi
302 
303  poisson_ess_bdr = 0;
304  poisson_ess_bdr[0] = 1; // boundary attribute 1 (index 0) is fixed
305  poisson_ess_bdr[1] = 1; // boundary attribute 2 (index 1) is fixed
306  // END CODE FOR THE COIL PROBLEM
307  }
308  else if (strcmp(problem,"rod")==0)
309  {
310  // BEGIN CODE FOR THE STRAIGHT ROD PROBLEM
311  // the boundary conditions below are for the straight rod problem
312 
313  ess_bdr = 0;
314  ess_bdr[0] = 1; // boundary attribute 1 (index 0) is fixed (front)
315  ess_bdr[1] = 1; // boundary attribute 2 (index 1) is fixed (rear)
316  ess_bdr[2] = 1; // boundary attribute 3 (index 2) is fixed (outer)
317 
318  // Same as above, but this is for the thermal operator. For HDiv
319  // formulation the essential BC is the flux, which is zero on the front
320  // and sides. Note the Natural BC is T = 0 on the outer surface.
321 
322  thermal_ess_bdr = 0;
323  thermal_ess_bdr[0] = 1; // boundary attribute 1 (index 0) is fixed (front)
324  thermal_ess_bdr[1] = 1; // boundary attribute 2 (index 1) is fixed (rear)
325 
326  // Same as above, but this is for the poisson eq for H1 formulation the
327  // essential BC is the value of Phi
328 
329  poisson_ess_bdr = 0;
330  poisson_ess_bdr[0] = 1; // boundary attribute 1 (index 0) is fixed (front)
331  poisson_ess_bdr[1] = 1; // boundary attribute 2 (index 1) is fixed (back)
332  // END CODE FOR THE STRAIGHT ROD PROBLEM
333  }
334  else
335  {
336  cerr << "Problem " << problem << " not recognized\n";
337  mfem_error();
338  }
339 
340  // The following is required for mesh refinement
341  mesh->EnsureNCMesh();
342 
343  // 6. Define the ODE solver used for time integration. Several implicit
344  // methods are available, including singly diagonal implicit Runge-Kutta
345  // (SDIRK).
346  ODESolver *ode_solver;
347  switch (ode_solver_type)
348  {
349  // Implicit L-stable methods
350  case 1: ode_solver = new BackwardEulerSolver; break;
351  case 2: ode_solver = new SDIRK23Solver(2); break;
352  case 3: ode_solver = new SDIRK33Solver; break;
353  // Implicit A-stable methods (not L-stable)
354  case 22: ode_solver = new ImplicitMidpointSolver; break;
355  case 23: ode_solver = new SDIRK23Solver; break;
356  case 34: ode_solver = new SDIRK34Solver; break;
357  default:
358  if (mpi.Root())
359  {
360  cout << "Unknown ODE solver type: " << ode_solver_type << '\n';
361  }
362  delete mesh;
363  return 3;
364  }
365 
366  // 7. Refine the mesh in serial to increase the resolution. In this example
367  // we do 'ser_ref_levels' of uniform refinement, where 'ser_ref_levels' is
368  // a command-line parameter.
369  for (int lev = 0; lev < ser_ref_levels; lev++)
370  {
371  mesh->UniformRefinement();
372  }
373 
374  // 8. Define a parallel mesh by a partitioning of the serial mesh. Refine
375  // this mesh further in parallel to increase the resolution. Once the
376  // parallel mesh is defined, the serial mesh can be deleted.
377  ParMesh *pmesh = new ParMesh(MPI_COMM_WORLD, *mesh);
378  delete mesh;
379  for (int lev = 0; lev < par_ref_levels; lev++)
380  {
381  pmesh->UniformRefinement();
382  }
383  // Make sure tet-only meshes are marked for local refinement.
384  pmesh->Finalize(true);
385 
386  // 9. Apply non-uniform non-conforming mesh refinement to the mesh. The
387  // whole metal region is refined once, before the start of the time loop,
388  // i.e. this is not based on any error estimator.
389  if (amr == 1)
390  {
391  Array<int> ref_list;
392  int numElems = pmesh->GetNE();
393  for (int ielem = 0; ielem < numElems; ielem++)
394  {
395  int thisAtt = pmesh->GetAttribute(ielem);
396  if (thisAtt == 1)
397  {
398  ref_list.Append(ielem);
399  }
400  }
401 
402  pmesh->GeneralRefinement(ref_list);
403 
404  ref_list.DeleteAll();
405  }
406 
407  // 10. Reorient the mesh. Must be done after refinement but before definition
408  // of higher order Nedelec spaces
409  pmesh->ReorientTetMesh();
410 
411  // 11. Rebalance the mesh. Since the mesh was adaptively refined in a
412  // non-uniform way it will be computationally unbalanced.
413  if (pmesh->Nonconforming())
414  {
415  pmesh->Rebalance();
416  }
417 
418  // 12. Define the parallel finite element spaces. We use:
419  //
420  // H(curl) for electric field,
421  // H(div) for magnetic flux,
422  // H(div) for thermal flux,
423  // H(grad)/H1 for electrostatic potential,
424  // L2 for temperature
425 
426  // L2 contains discontinuous "cell-center" finite elements, type 2 is
427  // "positive"
428  L2_FECollection L2FEC(order-1, dim);
429 
430  // ND contains Nedelec "edge-centered" vector finite elements with continuous
431  // tangential component.
432  ND_FECollection HCurlFEC(order, dim);
433 
434  // RT contains Raviart-Thomas "face-centered" vector finite elements with
435  // continuous normal component.
436  RT_FECollection HDivFEC(order-1, dim);
437 
438  // H1 contains continuous "node-centered" Lagrange finite elements.
439  H1_FECollection HGradFEC(order, dim);
440 
441  ParFiniteElementSpace L2FESpace(pmesh, &L2FEC);
442  ParFiniteElementSpace HCurlFESpace(pmesh, &HCurlFEC);
443  ParFiniteElementSpace HDivFESpace(pmesh, &HDivFEC);
444  ParFiniteElementSpace HGradFESpace(pmesh, &HGradFEC);
445 
446  // The terminology is TrueVSize is the unique (non-redundant) number of dofs
447  HYPRE_Int glob_size_l2 = L2FESpace.GlobalTrueVSize();
448  HYPRE_Int glob_size_nd = HCurlFESpace.GlobalTrueVSize();
449  HYPRE_Int glob_size_rt = HDivFESpace.GlobalTrueVSize();
450  HYPRE_Int glob_size_h1 = HGradFESpace.GlobalTrueVSize();
451 
452  if (mpi.Root())
453  {
454  cout << "Number of Temperature Flux unknowns: " << glob_size_rt << endl;
455  cout << "Number of Temperature unknowns: " << glob_size_l2 << endl;
456  cout << "Number of Electric Field unknowns: " << glob_size_nd << endl;
457  cout << "Number of Magnetic Field unknowns: " << glob_size_rt << endl;
458  cout << "Number of Electrostatic unknowns: " << glob_size_h1 << endl;
459  }
460 
461  int Vsize_l2 = L2FESpace.GetVSize();
462  int Vsize_nd = HCurlFESpace.GetVSize();
463  int Vsize_rt = HDivFESpace.GetVSize();
464  int Vsize_h1 = HGradFESpace.GetVSize();
465 
466  // the big BlockVector stores the fields as
467  // 0 Temperature
468  // 1 Temperature Flux
469  // 2 P field
470  // 3 E field
471  // 4 B field
472  // 5 Joule Heating
473 
474  Array<int> true_offset(7);
475  true_offset[0] = 0;
476  true_offset[1] = true_offset[0] + Vsize_l2;
477  true_offset[2] = true_offset[1] + Vsize_rt;
478  true_offset[3] = true_offset[2] + Vsize_h1;
479  true_offset[4] = true_offset[3] + Vsize_nd;
480  true_offset[5] = true_offset[4] + Vsize_rt;
481  true_offset[6] = true_offset[5] + Vsize_l2;
482 
483  // The BlockVector is a large contiguous chunk of memory for storing required
484  // data for the hypre vectors, in this case: the temperature L2, the T-flux
485  // HDiv, the E-field HCurl, and the B-field HDiv, and scalar potential P.
486  BlockVector F(true_offset);
487 
488  // grid functions E, B, T, F, P, and w which is the Joule heating
489  ParGridFunction E_gf, B_gf, T_gf, F_gf, w_gf, P_gf;
490  T_gf.MakeRef(&L2FESpace,F, true_offset[0]);
491  F_gf.MakeRef(&HDivFESpace,F, true_offset[1]);
492  P_gf.MakeRef(&HGradFESpace,F,true_offset[2]);
493  E_gf.MakeRef(&HCurlFESpace,F,true_offset[3]);
494  B_gf.MakeRef(&HDivFESpace,F, true_offset[4]);
495  w_gf.MakeRef(&L2FESpace,F, true_offset[5]);
496 
497  // 13. Get the boundary conditions, set up the exact solution grid functions
498  // These VectorCoefficients have an Eval function. Note that e_exact and
499  // b_exact in this case are exact analytical solutions, taking a 3-vector
500  // point as input and returning a 3-vector field
503  FunctionCoefficient T_exact(t_exact);
504  E_exact.SetTime(0.0);
505  B_exact.SetTime(0.0);
506 
507  // 14. Initialize the Diffusion operator, the GLVis visualization and print
508  // the initial energies.
509  MagneticDiffusionEOperator oper(true_offset[6], L2FESpace, HCurlFESpace,
510  HDivFESpace, HGradFESpace,
511  ess_bdr, thermal_ess_bdr, poisson_ess_bdr,
512  mu, sigmaMap, TcapMap, InvTcapMap,
513  InvTcondMap);
514 
515  // This function initializes all the fields to zero or some provided IC
516  oper.Init(F);
517 
518  socketstream vis_T, vis_E, vis_B, vis_w, vis_P;
519  char vishost[] = "localhost";
520  int visport = 19916;
521  if (visualization)
522  {
523  // Make sure all ranks have sent their 'v' solution before initiating
524  // another set of GLVis connections (one from each rank):
525  MPI_Barrier(pmesh->GetComm());
526 
527  vis_T.precision(8);
528  vis_E.precision(8);
529  vis_B.precision(8);
530  vis_P.precision(8);
531  vis_w.precision(8);
532 
533  int Wx = 0, Wy = 0; // window position
534  int Ww = 350, Wh = 350; // window size
535  int offx = Ww+10, offy = Wh+45; // window offsets
536 
537  miniapps::VisualizeField(vis_P, vishost, visport,
538  P_gf, "Electric Potential (Phi)", Wx, Wy, Ww, Wh);
539  Wx += offx;
540 
541  miniapps::VisualizeField(vis_E, vishost, visport,
542  E_gf, "Electric Field (E)", Wx, Wy, Ww, Wh);
543  Wx += offx;
544 
545  miniapps::VisualizeField(vis_B, vishost, visport,
546  B_gf, "Magnetic Field (B)", Wx, Wy, Ww, Wh);
547  Wx = 0;
548  Wy += offy;
549 
550  miniapps::VisualizeField(vis_w, vishost, visport,
551  w_gf, "Joule Heating", Wx, Wy, Ww, Wh);
552 
553  Wx += offx;
554 
555  miniapps::VisualizeField(vis_T, vishost, visport,
556  T_gf, "Temperature", Wx, Wy, Ww, Wh);
557  }
558  // VisIt visualization
559  VisItDataCollection visit_dc(basename, pmesh);
560  if ( visit )
561  {
562  visit_dc.RegisterField("E", &E_gf);
563  visit_dc.RegisterField("B", &B_gf);
564  visit_dc.RegisterField("T", &T_gf);
565  visit_dc.RegisterField("w", &w_gf);
566  visit_dc.RegisterField("Phi", &P_gf);
567  visit_dc.RegisterField("F", &F_gf);
568 
569  visit_dc.SetCycle(0);
570  visit_dc.SetTime(0.0);
571  visit_dc.Save();
572  }
573 
574  E_exact.SetTime(0.0);
575  B_exact.SetTime(0.0);
576 
577  // 15. Perform time-integration (looping over the time iterations, ti, with a
578  // time-step dt). The object oper is the MagneticDiffusionOperator which
579  // has a Mult() method and an ImplicitSolve() method which are used by
580  // the time integrators.
581  ode_solver->Init(oper);
582  double t = 0.0;
583 
584  bool last_step = false;
585  for (int ti = 1; !last_step; ti++)
586  {
587  if (t + dt >= t_final - dt/2)
588  {
589  last_step = true;
590  }
591 
592  // F is the vector of dofs, t is the current time, and dt is the time step
593  // to advance.
594  ode_solver->Step(F, t, dt);
595 
596  if (debug == 1)
597  {
598  oper.Debug(basename,t);
599  }
600 
601  if (gfprint == 1)
602  {
603  ostringstream T_name, E_name, B_name, F_name, w_name, P_name, mesh_name;
604  T_name << basename << "_" << setfill('0') << setw(6) << t << "_"
605  << "T." << setfill('0') << setw(6) << myid;
606  E_name << basename << "_" << setfill('0') << setw(6) << t << "_"
607  << "E." << setfill('0') << setw(6) << myid;
608  B_name << basename << "_" << setfill('0') << setw(6) << t << "_"
609  << "B." << setfill('0') << setw(6) << myid;
610  F_name << basename << "_" << setfill('0') << setw(6) << t << "_"
611  << "F." << setfill('0') << setw(6) << myid;
612  w_name << basename << "_" << setfill('0') << setw(6) << t << "_"
613  << "w." << setfill('0') << setw(6) << myid;
614  P_name << basename << "_" << setfill('0') << setw(6) << t << "_"
615  << "P." << setfill('0') << setw(6) << myid;
616  mesh_name << basename << "_" << setfill('0') << setw(6) << t << "_"
617  << "mesh." << setfill('0') << setw(6) << myid;
618 
619  ofstream mesh_ofs(mesh_name.str().c_str());
620  mesh_ofs.precision(8);
621  pmesh->Print(mesh_ofs);
622  mesh_ofs.close();
623 
624  ofstream T_ofs(T_name.str().c_str());
625  T_ofs.precision(8);
626  T_gf.Save(T_ofs);
627  T_ofs.close();
628 
629  ofstream E_ofs(E_name.str().c_str());
630  E_ofs.precision(8);
631  E_gf.Save(E_ofs);
632  E_ofs.close();
633 
634  ofstream B_ofs(B_name.str().c_str());
635  B_ofs.precision(8);
636  B_gf.Save(B_ofs);
637  B_ofs.close();
638 
639  ofstream F_ofs(F_name.str().c_str());
640  F_ofs.precision(8);
641  F_gf.Save(B_ofs);
642  F_ofs.close();
643 
644  ofstream P_ofs(P_name.str().c_str());
645  P_ofs.precision(8);
646  P_gf.Save(P_ofs);
647  P_ofs.close();
648 
649  ofstream w_ofs(w_name.str().c_str());
650  w_ofs.precision(8);
651  w_gf.Save(w_ofs);
652  w_ofs.close();
653  }
654 
655  if (last_step || (ti % vis_steps) == 0)
656  {
657  double el = oper.ElectricLosses(E_gf);
658 
659  if (mpi.Root())
660  {
661  cout << fixed;
662  cout << "step " << setw(6) << ti << ",\tt = " << setw(6)
663  << setprecision(3) << t
664  << ",\tdot(E, J) = " << setprecision(8) << el << endl;
665  }
666 
667  // Make sure all ranks have sent their 'v' solution before initiating
668  // another set of GLVis connections (one from each rank):
669  MPI_Barrier(pmesh->GetComm());
670 
671  if (visualization)
672  {
673  int Wx = 0, Wy = 0; // window position
674  int Ww = 350, Wh = 350; // window size
675  int offx = Ww+10, offy = Wh+45; // window offsets
676 
677  miniapps::VisualizeField(vis_P, vishost, visport,
678  P_gf, "Electric Potential (Phi)", Wx, Wy, Ww, Wh);
679  Wx += offx;
680 
681  miniapps::VisualizeField(vis_E, vishost, visport,
682  E_gf, "Electric Field (E)", Wx, Wy, Ww, Wh);
683  Wx += offx;
684 
685  miniapps::VisualizeField(vis_B, vishost, visport,
686  B_gf, "Magnetic Field (B)", Wx, Wy, Ww, Wh);
687 
688  Wx = 0;
689  Wy += offy;
690 
691  miniapps::VisualizeField(vis_w, vishost, visport,
692  w_gf, "Joule Heating", Wx, Wy, Ww, Wh);
693 
694  Wx += offx;
695 
696  miniapps::VisualizeField(vis_T, vishost, visport,
697  T_gf, "Temperature", Wx, Wy, Ww, Wh);
698  }
699 
700  if (visit)
701  {
702  visit_dc.SetCycle(ti);
703  visit_dc.SetTime(t);
704  visit_dc.Save();
705  }
706  }
707  }
708  if (visualization)
709  {
710  vis_T.close();
711  vis_E.close();
712  vis_B.close();
713  vis_w.close();
714  vis_P.close();
715  }
716 
717  // 16. Free the used memory.
718  delete ode_solver;
719  delete pmesh;
720 
721  return 0;
722 }
723 
724 namespace mfem
725 {
726 
727 namespace electromagnetics
728 {
729 
730 void edot_bc(const Vector &x, Vector &E)
731 {
732  E = 0.0;
733 }
734 
735 void e_exact(const Vector &x, double t, Vector &E)
736 {
737  E[0] = 0.0;
738  E[1] = 0.0;
739  E[2] = 0.0;
740 }
741 
742 void b_exact(const Vector &x, double t, Vector &B)
743 {
744  B[0] = 0.0;
745  B[1] = 0.0;
746  B[2] = 0.0;
747 }
748 
749 double t_exact(Vector &x)
750 {
751  double T = 0.0;
752  return T;
753 }
754 
755 double p_bc(const Vector &x, double t)
756 {
757  // the value
758  double T;
759  if (x[2] < 0.0)
760  {
761  T = 1.0;
762  }
763  else
764  {
765  T = -1.0;
766  }
767 
768  return T*cos(wj_ * t);
769 }
770 
771 } // namespace electromagnetics
772 
773 } // namespace mfem
774 
775 void display_banner(ostream & os)
776 {
777  os << " ____. .__ " << endl
778  << " | | ____ __ __| | ____ " << endl
779  << " | |/ _ \\| | \\ | _/ __ \\ " << endl
780  << "/\\__| ( <_> ) | / |_\\ ___/ " << endl
781  << "\\________|\\____/|____/|____/\\___ >" << endl
782  << " \\/ " << endl
783  << flush;
784 }
int GetVSize() const
Return the number of vector dofs, i.e. GetNDofs() x GetVDim().
Definition: fespace.hpp:347
void SetCycle(int c)
Set time cycle (for time-dependent simulations)
A class to handle Vectors in a block fashion.
Definition: blockvector.hpp:30
virtual void ReorientTetMesh()
See the remarks for the serial version in mesh.hpp.
Definition: pmesh.cpp:2402
int offx
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]...
void e_exact(const Vector &x, double t, Vector &E)
Definition: joule.cpp:735
int Wx
double t_exact(Vector &x)
Definition: joule.cpp:749
virtual void Init(TimeDependentOperator &f)
Associate a TimeDependentOperator with the ODE solver.
Definition: ode.cpp:18
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:676
virtual void Save()
Save the collection and a VisIt root file.
virtual void Save(std::ostream &out) const
Definition: pgridfunc.cpp:488
Abstract parallel finite element space.
Definition: pfespace.hpp:28
int main(int argc, char *argv[])
Definition: ex1.cpp:58
Backward Euler ODE solver. L-stable.
Definition: ode.hpp:210
void DeleteAll()
Delete whole array.
Definition: array.hpp:785
void Debug(const char *basefilename, double time)
virtual void MakeRef(FiniteElementSpace *f, double *v)
Make the ParGridFunction reference external data on a new FiniteElementSpace.
Definition: pgridfunc.cpp:101
double p_bc(const Vector &x, double t)
Definition: joule.cpp:755
void Rebalance()
Load balance the mesh. NC meshes only.
Definition: pmesh.cpp:3085
bool Nonconforming() const
Definition: mesh.hpp:1136
A simple convenience class that calls MPI_Init() at construction and MPI_Finalize() at destruction...
int dim
Definition: ex3.cpp:48
int Append(const T &el)
Append element to array, resize if necessary.
Definition: array.hpp:690
void mfem_error(const char *msg)
Function called when an error is encountered. Used by the macros MFEM_ABORT, MFEM_ASSERT, MFEM_VERIFY.
Definition: error.cpp:146
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
Definition: mesh.cpp:7591
Data collection with VisIt I/O routines.
T Max() const
Find the maximal element in the array, using the comparison operator &lt; for class T.
Definition: array.cpp:68
HYPRE_Int GlobalTrueVSize() const
Definition: pfespace.hpp:251
bool Root() const
Return true if WorldRank() == 0.
virtual void Print(std::ostream &out=mfem::out) const
Definition: pmesh.cpp:3863
int Dimension() const
Definition: mesh.hpp:713
void PrintUsage(std::ostream &out) const
Definition: optparser.cpp:434
void SetTime(double t)
Set physical time (for time-dependent simulations)
Arbitrary order H(div)-conforming Raviart-Thomas finite elements.
Definition: fe_coll.hpp:180
double ElectricLosses(ParGridFunction &E_gf) const
Array< int > bdr_attributes
A list of all unique boundary attributes used by the Mesh.
Definition: mesh.hpp:179
int WorldRank() const
Return MPI_COMM_WORLD&#39;s rank.
MPI_Comm GetComm() const
Definition: pmesh.hpp:223
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)
Definition: optparser.hpp:76
int problem
Definition: ex15.cpp:57
void SetTime(double t)
virtual void RegisterField(const std::string &field_name, GridFunction *gf)
Add a grid function to the collection and update the root file.
Implicit midpoint method. A-stable, not L-stable.
Definition: ode.hpp:223
int Wh
void PrintOptions(std::ostream &out) const
Definition: optparser.cpp:304
void E_exact(const Vector &, Vector &)
Definition: ex3.cpp:212
virtual void Finalize(bool refine=false, bool fix_orientation=false)
Finalize the construction of a general Mesh.
Definition: pmesh.cpp:1281
void EnsureNCMesh(bool triangles_nonconforming=false)
Definition: mesh.cpp:7253
void edot_bc(const Vector &x, Vector &E)
Definition: joule.cpp:730
class for C-function coefficient
Arbitrary order H(curl)-conforming Nedelec finite elements.
Definition: fe_coll.hpp:240
Vector data type.
Definition: vector.hpp:48
int offy
int Wy
Arbitrary order H1-conforming (continuous) finite elements.
Definition: fe_coll.hpp:83
int Ww
void b_exact(const Vector &x, double t, Vector &B)
Definition: joule.cpp:742
Class for parallel grid function.
Definition: pgridfunc.hpp:32
void display_banner(ostream &os)
Definition: joule.cpp:775
void GeneralRefinement(const Array< Refinement > &refinements, int nonconforming=-1, int nc_limit=0)
Definition: mesh.cpp:7185
double freq
Definition: ex3.cpp:47
Class for parallel meshes.
Definition: pmesh.hpp:32
int GetAttribute(int i) const
Return the attribute of element i.
Definition: mesh.hpp:970
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)
Definition: fem_extras.cpp:94
Arbitrary order &quot;L2-conforming&quot; discontinuous finite elements.
Definition: fe_coll.hpp:134
double sigma(const Vector &x)
Definition: maxwell.cpp:102
bool Good() const
Definition: optparser.hpp:122