MFEM  v4.3.0
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
tesla.cpp
Go to the documentation of this file.
1 // Copyright (c) 2010-2021, 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 // Tesla Miniapp: Simple Magnetostatics Simulation Code
14 // -----------------------------------------------------
15 //
16 // This miniapp solves a simple 3D magnetostatic problem.
17 //
18 // Curl 1/mu Curl A = J + Curl mu0/mu M
19 //
20 // The permeability function is that of the vacuum with an optional diamagnetic
21 // or paramagnetic spherical shell. The optional current density takes the form
22 // of a user defined ring of current. The optional magnetization consists of a
23 // cylindrical bar of constant magnetization.
24 //
25 // The boundary conditions either apply a user selected uniform magnetic flux
26 // density or a surface current flowing between user defined surfaces.
27 //
28 // We discretize the vector potential with H(Curl) finite elements. The magnetic
29 // flux B is discretized with H(Div) finite elements.
30 //
31 // Compile with: make tesla
32 //
33 // Sample runs:
34 //
35 // A cylindrical bar magnet in a metal sphere:
36 // mpirun -np 4 tesla -bm '0 -0.5 0 0 0.5 0 0.2 1'
37 //
38 // A spherical shell of paramagnetic material in a uniform B field:
39 // mpirun -np 4 tesla -ubbc '0 0 1' -ms '0 0 0 0.2 0.4 10'
40 //
41 // A ring of current in a metal sphere:
42 // mpirun -np 4 tesla -cr '0 0 -0.2 0 0 0.2 0.2 0.4 1'
43 //
44 // A Halbach array of permanent magnets:
45 // mpirun -np 4 tesla -m ../../data/beam-hex.mesh -rs 2
46 // -ha '1 0.1 0.3 7 0.9 0.7 0 1 12'
47 //
48 // An example demonstrating the use of surface currents:
49 // mpirun -np 4 tesla -m square-angled-pipe.mesh
50 // -kbcs '3' -vbcs '1 2' -vbcv '-0.5 0.5'
51 //
52 // An example combining the paramagnetic shell, permanent magnet,
53 // and current ring:
54 // mpirun -np 4 tesla -m ../../data/inline-hex.mesh
55 // -ms '0.5 0.5 0.5 0.4 0.45 20'
56 // -bm '0.5 0.5 0.3 0.5 0.5 0.7 0.1 1'
57 // -cr '0.5 0.5 0.45 0.5 0.5 0.55 0.2 0.3 1'
58 //
59 // By default the sources and fields are all zero:
60 // mpirun -np 4 tesla
61 
62 #include "tesla_solver.hpp"
63 #include <fstream>
64 #include <iostream>
65 
66 using namespace std;
67 using namespace mfem;
68 using namespace mfem::electromagnetics;
69 
70 // Permeability Function
72 
73 static Vector pw_mu_(0); // Piecewise permeability values
74 static Vector pw_mu_inv_(0); // Piecewise inverse permeability values
75 static Vector ms_params_(0); // Center, Inner and Outer Radii, and
76 // Permeability of magnetic shell
77 double magnetic_shell(const Vector &);
78 double magnetic_shell_inv(const Vector & x) { return 1.0/magnetic_shell(x); }
79 
80 // Current Density Function
81 static Vector cr_params_(0); // Axis Start, Axis End, Inner Ring Radius,
82 // Outer Ring Radius, and Total Current
83 // of current ring (annulus)
84 void current_ring(const Vector &, Vector &);
85 
86 // Magnetization
87 static Vector bm_params_(0); // Axis Start, Axis End, Bar Radius,
88 // and Magnetic Field Magnitude
89 void bar_magnet(const Vector &, Vector &);
90 
91 static Vector ha_params_(0); // Bounding box,
92 // axis index (0->'x', 1->'y', 2->'z'),
93 // rotation axis index
94 // and number of segments
95 void halbach_array(const Vector &, Vector &);
96 
97 // A Field Boundary Condition for B = (Bx,By,Bz)
98 static Vector b_uniform_(0);
99 void a_bc_uniform(const Vector &, Vector&);
100 
101 // Phi_M Boundary Condition for H = (0,0,1)
102 double phi_m_bc_uniform(const Vector &x);
103 
104 // Prints the program's logo to the given output stream
105 void display_banner(ostream & os);
106 
107 int main(int argc, char *argv[])
108 {
109  MPI_Session mpi(argc, argv);
110 
111  if ( mpi.Root() ) { display_banner(cout); }
112 
113  // Parse command-line options.
114  const char *mesh_file = "../../data/ball-nurbs.mesh";
115  int order = 1;
116  int maxit = 100;
117  int serial_ref_levels = 0;
118  int parallel_ref_levels = 0;
119  bool visualization = true;
120  bool visit = true;
121 
122  Array<int> kbcs;
123  Array<int> vbcs;
124 
125  Vector vbcv;
126 
127  OptionsParser args(argc, argv);
128  args.AddOption(&mesh_file, "-m", "--mesh",
129  "Mesh file to use.");
130  args.AddOption(&order, "-o", "--order",
131  "Finite element order (polynomial degree).");
132  args.AddOption(&serial_ref_levels, "-rs", "--serial-ref-levels",
133  "Number of serial refinement levels.");
134  args.AddOption(&parallel_ref_levels, "-rp", "--parallel-ref-levels",
135  "Number of parallel refinement levels.");
136  args.AddOption(&b_uniform_, "-ubbc", "--uniform-b-bc",
137  "Specify if the three components of the constant magnetic flux density");
138  args.AddOption(&pw_mu_, "-pwm", "--piecewise-mu",
139  "Piecewise values of Permeability");
140  args.AddOption(&ms_params_, "-ms", "--magnetic-shell-params",
141  "Center, Inner Radius, Outer Radius, and Permeability of Magnetic Shell");
142  args.AddOption(&cr_params_, "-cr", "--current-ring-params",
143  "Axis End Points, Inner Radius, Outer Radius and Total Current of Annulus");
144  args.AddOption(&bm_params_, "-bm", "--bar-magnet-params",
145  "Axis End Points, Radius, and Magnetic Field of Cylindrical Magnet");
146  args.AddOption(&ha_params_, "-ha", "--halbach-array-params",
147  "Bounding Box Corners and Number of Segments");
148  args.AddOption(&kbcs, "-kbcs", "--surface-current-bc",
149  "Surfaces for the Surface Current (K) Boundary Condition");
150  args.AddOption(&vbcs, "-vbcs", "--voltage-bc-surf",
151  "Voltage Boundary Condition Surfaces (to drive K)");
152  args.AddOption(&vbcv, "-vbcv", "--voltage-bc-vals",
153  "Voltage Boundary Condition Values (to drive K)");
154  args.AddOption(&maxit, "-maxit", "--max-amr-iterations",
155  "Max number of iterations in the main AMR loop.");
156  args.AddOption(&visualization, "-vis", "--visualization", "-no-vis",
157  "--no-visualization",
158  "Enable or disable GLVis visualization.");
159  args.AddOption(&visit, "-visit", "--visit", "-no-visit", "--no-visit",
160  "Enable or disable VisIt visualization.");
161  args.Parse();
162  if (!args.Good())
163  {
164  if (mpi.Root())
165  {
166  args.PrintUsage(cout);
167  }
168  return 1;
169  }
170  if (mpi.Root())
171  {
172  args.PrintOptions(cout);
173  }
174 
175  // Read the (serial) mesh from the given mesh file on all processors. We
176  // can handle triangular, quadrilateral, tetrahedral, hexahedral, surface
177  // and volume meshes with the same code.
178  Mesh *mesh = new Mesh(mesh_file, 1, 1);
179 
180  if (mpi.Root())
181  {
182  cout << "Starting initialization." << endl;
183  }
184 
185  // Project a NURBS mesh to a piecewise-quadratic curved mesh
186  if (mesh->NURBSext)
187  {
188  mesh->UniformRefinement();
189  if (serial_ref_levels > 0) { serial_ref_levels--; }
190 
191  mesh->SetCurvature(2);
192  }
193 
194  // Ensure that quad and hex meshes are treated as non-conforming.
195  mesh->EnsureNCMesh();
196 
197  // Refine the serial mesh on all processors to increase the resolution. In
198  // this example we do 'ref_levels' of uniform refinement.
199  for (int l = 0; l < serial_ref_levels; l++)
200  {
201  mesh->UniformRefinement();
202  }
203 
204  // Define a parallel mesh by a partitioning of the serial mesh. Refine
205  // this mesh further in parallel to increase the resolution. Once the
206  // parallel mesh is defined, the serial mesh can be deleted.
207  ParMesh pmesh(MPI_COMM_WORLD, *mesh);
208  delete mesh;
209 
210  // Refine this mesh in parallel to increase the resolution.
211  int par_ref_levels = parallel_ref_levels;
212  for (int l = 0; l < par_ref_levels; l++)
213  {
214  pmesh.UniformRefinement();
215  }
216  // Make sure tet-only meshes are marked for local refinement.
217  pmesh.Finalize(true);
218 
219  // If values for Voltage BCs were not set issue a warning and exit
220  if ( ( vbcs.Size() > 0 && kbcs.Size() == 0 ) ||
221  ( kbcs.Size() > 0 && vbcs.Size() == 0 ) ||
222  ( vbcv.Size() < vbcs.Size() ) )
223  {
224  if ( mpi.Root() )
225  {
226  cout << "The surface current (K) boundary condition requires "
227  << "surface current boundary condition surfaces (with -kbcs), "
228  << "voltage boundary condition surface (with -vbcs), "
229  << "and voltage boundary condition values (with -vbcv)."
230  << endl;
231  }
232  return 3;
233  }
234 
235  // Create a coefficient describing the magnetic permeability
237 
238  // Create the Magnetostatic solver
239  TeslaSolver Tesla(pmesh, order, kbcs, vbcs, vbcv, *muInvCoef,
240  (b_uniform_.Size() > 0 ) ? a_bc_uniform : NULL,
241  (cr_params_.Size() > 0 ) ? current_ring : NULL,
242  (bm_params_.Size() > 0 ) ? bar_magnet :
243  (ha_params_.Size() > 0 ) ? halbach_array : NULL);
244 
245  // Initialize GLVis visualization
246  if (visualization)
247  {
248  Tesla.InitializeGLVis();
249  }
250 
251  // Initialize VisIt visualization
252  VisItDataCollection visit_dc("Tesla-AMR-Parallel", &pmesh);
253 
254  if ( visit )
255  {
256  Tesla.RegisterVisItFields(visit_dc);
257  }
258  if (mpi.Root()) { cout << "Initialization done." << endl; }
259 
260  // The main AMR loop. In each iteration we solve the problem on the current
261  // mesh, visualize the solution, estimate the error on all elements, refine
262  // the worst elements and update all objects to work with the new mesh. We
263  // refine until the maximum number of dofs in the Nedelec finite element
264  // space reaches 10 million.
265  const int max_dofs = 10000000;
266  for (int it = 1; it <= maxit; it++)
267  {
268  if (mpi.Root())
269  {
270  cout << "\nAMR Iteration " << it << endl;
271  }
272 
273  // Display the current number of DoFs in each finite element space
274  Tesla.PrintSizes();
275 
276  // Assemble all forms
277  Tesla.Assemble();
278 
279  // Solve the system and compute any auxiliary fields
280  Tesla.Solve();
281 
282  // Determine the current size of the linear system
283  int prob_size = Tesla.GetProblemSize();
284 
285  // Write fields to disk for VisIt
286  if ( visit )
287  {
288  Tesla.WriteVisItFields(it);
289  }
290 
291  // Send the solution by socket to a GLVis server.
292  if (visualization)
293  {
294  Tesla.DisplayToGLVis();
295  }
296 
297  if (mpi.Root())
298  {
299  cout << "AMR iteration " << it << " complete." << endl;
300  }
301 
302  // Check stopping criteria
303  if (prob_size > max_dofs)
304  {
305  if (mpi.Root())
306  {
307  cout << "Reached maximum number of dofs, exiting..." << endl;
308  }
309  break;
310  }
311  if ( it == maxit )
312  {
313  break;
314  }
315 
316  // Wait for user input. Ask every 10th iteration.
317  char c = 'c';
318  if (mpi.Root() && (it % 10 == 0))
319  {
320  cout << "press (q)uit or (c)ontinue --> " << flush;
321  cin >> c;
322  }
323  MPI_Bcast(&c, 1, MPI_CHAR, 0, MPI_COMM_WORLD);
324 
325  if (c != 'c')
326  {
327  break;
328  }
329 
330  // Estimate element errors using the Zienkiewicz-Zhu error estimator.
331  Vector errors(pmesh.GetNE());
332  Tesla.GetErrorEstimates(errors);
333 
334  double local_max_err = errors.Max();
335  double global_max_err;
336  MPI_Allreduce(&local_max_err, &global_max_err, 1,
337  MPI_DOUBLE, MPI_MAX, pmesh.GetComm());
338 
339  // Refine the elements whose error is larger than a fraction of the
340  // maximum element error.
341  const double frac = 0.5;
342  double threshold = frac * global_max_err;
343  if (mpi.Root()) { cout << "Refining ..." << endl; }
344  pmesh.RefineByError(errors, threshold);
345 
346  // Update the magnetostatic solver to reflect the new state of the mesh.
347  Tesla.Update();
348 
349  if (pmesh.Nonconforming() && mpi.WorldSize() > 1)
350  {
351  if (mpi.Root()) { cout << "Rebalancing ..." << endl; }
352  pmesh.Rebalance();
353 
354  // Update again after rebalancing
355  Tesla.Update();
356  }
357  }
358 
359  delete muInvCoef;
360 
361  return 0;
362 }
363 
364 // Print the Volta ascii logo to the given ostream
365 void display_banner(ostream & os)
366 {
367  os << " ___________ __ " << endl
368  << " \\__ ___/___ _____| | _____ " << endl
369  << " | |_/ __ \\ / ___/ | \\__ \\ " << endl
370  << " | |\\ ___/ \\___ \\| |__/ __ \\_ " << endl
371  << " |____| \\___ >____ >____(____ / " << endl
372  << " \\/ \\/ \\/ " << endl << flush;
373 }
374 
375 // The Permeability is a required coefficient which may be defined in
376 // various ways so we'll determine the appropriate coefficient type here.
377 Coefficient *
379 {
380  Coefficient * coef = NULL;
381 
382  if ( ms_params_.Size() > 0 )
383  {
385  }
386  else if ( pw_mu_.Size() > 0 )
387  {
388  pw_mu_inv_.SetSize(pw_mu_.Size());
389  for (int i = 0; i < pw_mu_.Size(); i++)
390  {
391  MFEM_ASSERT( pw_mu_[i] > 0.0, "permeability values must be positive" );
392  pw_mu_inv_[i] = 1.0/pw_mu_[i];
393  }
394  coef = new PWConstCoefficient(pw_mu_inv_);
395  }
396  else
397  {
398  coef = new ConstantCoefficient(1.0/mu0_);
399  }
400 
401  return coef;
402 }
403 
404 // A spherical shell with constant permeability. The sphere has inner
405 // and outer radii, center, and relative permeability specified on the
406 // command line and stored in ms_params_.
407 double magnetic_shell(const Vector &x)
408 {
409  double r2 = 0.0;
410 
411  for (int i = 0; i < x.Size(); i++)
412  {
413  r2 += (x(i) - ms_params_(i))*(x(i) - ms_params_(i));
414  }
415 
416  if ( sqrt(r2) >= ms_params_(x.Size()) &&
417  sqrt(r2) <= ms_params_(x.Size()+1) )
418  {
419  return mu0_*ms_params_(x.Size()+2);
420  }
421  return mu0_;
422 }
423 
424 // An annular ring of current density. The ring has two axis end
425 // points, inner and outer radii, and a constant current in Amperes.
426 void current_ring(const Vector &x, Vector &j)
427 {
428  MFEM_ASSERT(x.Size() == 3, "current_ring source requires 3D space.");
429 
430  j.SetSize(x.Size());
431  j = 0.0;
432 
433  Vector a(x.Size()); // Normalized Axis vector
434  Vector xu(x.Size()); // x vector relative to the axis end-point
435  Vector ju(x.Size()); // Unit vector in direction of current
436 
437  xu = x;
438 
439  for (int i=0; i<x.Size(); i++)
440  {
441  xu[i] -= cr_params_[i];
442  a[i] = cr_params_[x.Size()+i] - cr_params_[i];
443  }
444 
445  double h = a.Norml2();
446 
447  if ( h == 0.0 )
448  {
449  return;
450  }
451 
452  double ra = cr_params_[2*x.Size()+0];
453  double rb = cr_params_[2*x.Size()+1];
454  if ( ra > rb )
455  {
456  double rc = ra;
457  ra = rb;
458  rb = rc;
459  }
460  double xa = xu*a;
461 
462  if ( h > 0.0 )
463  {
464  xu.Add(-xa/(h*h),a);
465  }
466 
467  double xp = xu.Norml2();
468 
469  if ( xa >= 0.0 && xa <= h*h && xp >= ra && xp <= rb )
470  {
471  ju(0) = a(1) * xu(2) - a(2) * xu(1);
472  ju(1) = a(2) * xu(0) - a(0) * xu(2);
473  ju(2) = a(0) * xu(1) - a(1) * xu(0);
474  ju /= h;
475 
476  j.Add(cr_params_[2*x.Size()+2]/(h*(rb-ra)),ju);
477  }
478 }
479 
480 // A Cylindrical Rod of constant magnetization. The cylinder has two
481 // axis end points, a radius, and a constant magnetic field oriented
482 // along the axis.
483 void bar_magnet(const Vector &x, Vector &m)
484 {
485  m.SetSize(x.Size());
486  m = 0.0;
487 
488  Vector a(x.Size()); // Normalized Axis vector
489  Vector xu(x.Size()); // x vector relative to the axis end-point
490 
491  xu = x;
492 
493  for (int i=0; i<x.Size(); i++)
494  {
495  xu[i] -= bm_params_[i];
496  a[i] = bm_params_[x.Size()+i] - bm_params_[i];
497  }
498 
499  double h = a.Norml2();
500 
501  if ( h == 0.0 )
502  {
503  return;
504  }
505 
506  double r = bm_params_[2*x.Size()];
507  double xa = xu*a;
508 
509  if ( h > 0.0 )
510  {
511  xu.Add(-xa/(h*h),a);
512  }
513 
514  double xp = xu.Norml2();
515 
516  if ( xa >= 0.0 && xa <= h*h && xp <= r )
517  {
518  m.Add(bm_params_[2*x.Size()+1]/h,a);
519  }
520 }
521 
522 // A Square Rod of rotating magnetized segments. The rod is defined
523 // by a bounding box and a number of segments. The magnetization in
524 // each segment is constant and follows a rotating pattern.
525 void halbach_array(const Vector &x, Vector &m)
526 {
527  m.SetSize(x.Size());
528  m = 0.0;
529 
530  // Check Bounding Box
531  if ( x[0] < ha_params_[0] || x[0] > ha_params_[3] ||
532  x[1] < ha_params_[1] || x[1] > ha_params_[4] ||
533  x[2] < ha_params_[2] || x[2] > ha_params_[5] )
534  {
535  return;
536  }
537 
538  int ai = (int)ha_params_[6];
539  int ri = (int)ha_params_[7];
540  int n = (int)ha_params_[8];
541 
542  int i = (int)n * (x[ai] - ha_params_[ai]) /
543  (ha_params_[ai+3] - ha_params_[ai]);
544 
545  m[(ri + 1 + (i % 2)) % 3] = pow(-1.0,i/2);
546 }
547 
548 // To produce a uniform magnetic flux the vector potential can be set
549 // to ( By z, Bz x, Bx y).
550 void a_bc_uniform(const Vector & x, Vector & a)
551 {
552  a.SetSize(3);
553  a(0) = b_uniform_(1) * x(2);
554  a(1) = b_uniform_(2) * x(0);
555  a(2) = b_uniform_(0) * x(1);
556 }
557 
558 // To produce a uniform magnetic field the scalar potential can be set
559 // to -z (or -y in 2D).
560 double phi_m_bc_uniform(const Vector &x)
561 {
562  return -x(x.Size()-1);
563 }
int WorldSize() const
Return MPI_COMM_WORLD&#39;s size.
int Size() const
Return the logical size of the array.
Definition: array.hpp:134
A coefficient that is constant across space and time.
Definition: coefficient.hpp:78
void SetSize(int s)
Resize the vector to size s.
Definition: vector.hpp:513
int Size() const
Returns the size of the vector.
Definition: vector.hpp:190
void bar_magnet(const Vector &, Vector &)
Definition: tesla.cpp:483
int GetNE() const
Returns number of elements.
Definition: mesh.hpp:846
void GetErrorEstimates(Vector &errors)
double phi_m_bc_uniform(const Vector &x)
Definition: tesla.cpp:560
void Rebalance()
Definition: pmesh.cpp:3602
bool Nonconforming() const
Definition: mesh.hpp:1367
A simple convenience class that calls MPI_Init() at construction and MPI_Finalize() at destruction...
void Parse()
Parse the command-line options. Note that this function expects all the options provided through the ...
Definition: optparser.cpp:150
double magnetic_shell_inv(const Vector &x)
Definition: tesla.cpp:78
void EnsureNCMesh(bool simplices_nonconforming=false)
Definition: mesh.cpp:8800
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
Definition: mesh.cpp:9143
void current_ring(const Vector &, Vector &)
Definition: tesla.cpp:426
Data collection with VisIt I/O routines.
virtual void SetCurvature(int order, bool discont=false, int space_dim=-1, int ordering=1)
Definition: mesh.cpp:4882
bool Root() const
Return true if WorldRank() == 0.
void PrintUsage(std::ostream &out) const
Print the usage message.
Definition: optparser.cpp:457
void RegisterVisItFields(VisItDataCollection &visit_dc)
MPI_Comm GetComm() const
Definition: pmesh.hpp:276
Base class Coefficients that optionally depend on space and time. These are used by the BilinearFormI...
Definition: coefficient.hpp:39
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
bool RefineByError(const Array< double > &elem_error, double threshold, int nonconforming=-1, int nc_limit=0)
Definition: mesh.cpp:8869
void a_bc_uniform(const Vector &, Vector &)
Definition: tesla.cpp:550
double a
Definition: lissajous.cpp:41
NURBSExtension * NURBSext
Optional NURBS mesh extension.
Definition: mesh.hpp:206
A piecewise constant coefficient with the constants keyed off the element attribute numbers...
Definition: coefficient.hpp:94
Vector & Add(const double a, const Vector &Va)
(*this) += a * Va
Definition: vector.cpp:243
Coefficient * SetupInvPermeabilityCoefficient()
Definition: tesla.cpp:378
void PrintOptions(std::ostream &out) const
Print the options.
Definition: optparser.cpp:327
double magnetic_shell(const Vector &)
Definition: maxwell.cpp:400
virtual void Finalize(bool refine=false, bool fix_orientation=false)
Finalize the construction of a general Mesh.
Definition: pmesh.cpp:1510
A general function coefficient.
Vector data type.
Definition: vector.hpp:60
void display_banner(ostream &os)
Definition: joule.cpp:780
Class for parallel meshes.
Definition: pmesh.hpp:32
int main()
bool Good() const
Return true if the command line options were parsed successfully.
Definition: optparser.hpp:150
void halbach_array(const Vector &, Vector &)
Definition: tesla.cpp:525