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