MFEM  v3.4
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
volta.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 // Volta Miniapp: Simple Electrostatics Simulation Code
14 // -----------------------------------------------------
15 //
16 // This miniapp solves a simple 2D or 3D electrostatic problem.
17 //
18 // Div eps Grad Phi = rho
19 //
20 // The permittivity function is that of the vacuum with an optional dielectric
21 // sphere. The charge density is either zero or a user defined sphere of charge.
22 //
23 // Boundary conditions for the electric potential consist of a user defined
24 // piecewise constant potential or a potential leading to a user selected
25 // uniform electric field.
26 //
27 // We discretize the electric potential with H1 finite elements. The electric
28 // field E is discretized with Nedelec finite elements.
29 //
30 // Compile with: make volta
31 //
32 // Sample runs:
33 //
34 // A cylinder at constant voltage in a square, grounded metal pipe:
35 // mpirun -np 4 volta -m ../../data/square-disc.mesh
36 // -dbcs '1 2 3 4 5 6 7 8' -dbcv '0 0 0 0 1 1 1 1'
37 //
38 // A cylinder with a constant surface charge density in a square,
39 // grounded metal pipe:
40 // mpirun -np 4 volta -m ../../data/square-disc.mesh
41 // -nbcs '5 6 7 8' -nbcv '5e-11 5e-11 5e-11 5e-11'
42 // -dbcs '1 2 3 4'
43 //
44 // A cylindrical voltaic pile within a grounded metal sphere:
45 // mpirun -np 4 volta -dbcs 1 -vp '0 -0.5 0 0 0.5 0 0.2 1'
46 //
47 // A charged sphere, off-center, within a grounded metal sphere:
48 // mpirun -np 4 volta -dbcs 1 -cs '0.0 0.5 0.0 0.2 2.0e-11'
49 //
50 // A dielectric sphere suspended in a uniform electric field:
51 // mpirun -np 4 volta -dbcs 1 -dbcg -ds '0.0 0.0 0.0 0.2 8.0'
52 //
53 // An example using piecewise constant permittivity values
54 // mpirun -np 4 volta -m llnl.mesh -dbcs '4' -dbcv '0'
55 // -cs '8.5 8.5 17 1.57' -pwe '1 1 1 0.001'
56 //
57 // By default the sources and fields are all zero:
58 // mpirun -np 4 volta
59 
60 #include "volta_solver.hpp"
61 #include <fstream>
62 #include <iostream>
63 
64 using namespace std;
65 using namespace mfem;
66 using namespace mfem::electromagnetics;
67 
68 // Permittivity Functions
70 
71 static Vector pw_eps_(0); // Piecewise permittivity values
72 static Vector ds_params_(0); // Center, Radius, and Permittivity
73 // of dielectric sphere
74 double dielectric_sphere(const Vector &);
75 
76 // Charge Density Function
77 static Vector cs_params_(0); // Center, Radius, and Total Charge
78 // of charged sphere
79 double charged_sphere(const Vector &);
80 
81 // Polarization
82 static Vector vp_params_(0); // Axis Start, Axis End, Cylinder Radius,
83 // and Polarization Magnitude
84 void voltaic_pile(const Vector &, Vector &);
85 
86 // Phi Boundary Condition
87 static Vector e_uniform_(0);
88 double phi_bc_uniform(const Vector &);
89 
90 // Prints the program's logo to the given output stream
91 void display_banner(ostream & os);
92 
93 int main(int argc, char *argv[])
94 {
95  MPI_Session mpi(argc, argv);
96 
97  if ( mpi.Root() ) { display_banner(cout); }
98 
99  // Parse command-line options.
100  const char *mesh_file = "../../data/ball-nurbs.mesh";
101  int order = 1;
102  int maxit = 100;
103  int serial_ref_levels = 0;
104  int parallel_ref_levels = 0;
105  bool visualization = true;
106  bool visit = true;
107 
108  Array<int> dbcs;
109  Array<int> nbcs;
110 
111  Vector dbcv;
112  Vector nbcv;
113 
114  bool dbcg = false;
115 
116  OptionsParser args(argc, argv);
117  args.AddOption(&mesh_file, "-m", "--mesh",
118  "Mesh file to use.");
119  args.AddOption(&order, "-o", "--order",
120  "Finite element order (polynomial degree).");
121  args.AddOption(&serial_ref_levels, "-rs", "--serial-ref-levels",
122  "Number of serial refinement levels.");
123  args.AddOption(&parallel_ref_levels, "-rp", "--parallel-ref-levels",
124  "Number of parallel refinement levels.");
125  args.AddOption(&e_uniform_, "-uebc", "--uniform-e-bc",
126  "Specify if the three components of the constant "
127  "electric field");
128  args.AddOption(&pw_eps_, "-pwe", "--piecewise-eps",
129  "Piecewise values of Permittivity");
130  args.AddOption(&ds_params_, "-ds", "--dielectric-sphere-params",
131  "Center, Radius, and Permittivity of Dielectric Sphere");
132  args.AddOption(&cs_params_, "-cs", "--charged-sphere-params",
133  "Center, Radius, and Total Charge of Charged Sphere");
134  args.AddOption(&vp_params_, "-vp", "--voltaic-pile-params",
135  "Axis End Points, Radius, and "
136  "Polarization of Cylindrical Voltaic Pile");
137  args.AddOption(&dbcs, "-dbcs", "--dirichlet-bc-surf",
138  "Dirichlet Boundary Condition Surfaces");
139  args.AddOption(&dbcv, "-dbcv", "--dirichlet-bc-vals",
140  "Dirichlet Boundary Condition Values");
141  args.AddOption(&dbcg, "-dbcg", "--dirichlet-bc-gradient",
142  "-no-dbcg", "--no-dirichlet-bc-gradient",
143  "Dirichlet Boundary Condition Gradient (phi = -z)");
144  args.AddOption(&nbcs, "-nbcs", "--neumann-bc-surf",
145  "Neumann Boundary Condition Surfaces");
146  args.AddOption(&nbcv, "-nbcv", "--neumann-bc-vals",
147  "Neumann Boundary Condition Values");
148  args.AddOption(&maxit, "-maxit", "--max-amr-iterations",
149  "Max number of iterations in the main AMR loop.");
150  args.AddOption(&visualization, "-vis", "--visualization", "-no-vis",
151  "--no-visualization",
152  "Enable or disable GLVis visualization.");
153  args.AddOption(&visit, "-visit", "--visit", "-no-visit", "--no-visit",
154  "Enable or disable VisIt visualization.");
155  args.Parse();
156  if (!args.Good())
157  {
158  if (mpi.Root())
159  {
160  args.PrintUsage(cout);
161  }
162  return 1;
163  }
164  if (mpi.Root())
165  {
166  args.PrintOptions(cout);
167  }
168 
169  // Read the (serial) mesh from the given mesh file on all processors. We
170  // can handle triangular, quadrilateral, tetrahedral, hexahedral, surface
171  // and volume meshes with the same code.
172  Mesh *mesh = new Mesh(mesh_file, 1, 1);
173  int sdim = mesh->SpaceDimension();
174 
175  if (mpi.Root())
176  {
177  cout << "Starting initialization." << endl;
178  }
179 
180  // Project a NURBS mesh to a piecewise-quadratic curved mesh
181  if (mesh->NURBSext)
182  {
183  mesh->UniformRefinement();
184  if (serial_ref_levels > 0) { serial_ref_levels--; }
185 
186  mesh->SetCurvature(2);
187  }
188 
189  // Ensure that quad and hex meshes are treated as non-conforming.
190  mesh->EnsureNCMesh();
191 
192  // Refine the serial mesh on all processors to increase the resolution. In
193  // this example we do 'ref_levels' of uniform refinement. NURBS meshes are
194  // refined at least twice, as they are typically coarse.
195  for (int l = 0; l < serial_ref_levels; l++)
196  {
197  mesh->UniformRefinement();
198  }
199 
200  // Define a parallel mesh by a partitioning of the serial mesh. Refine
201  // this mesh further in parallel to increase the resolution. Once the
202  // parallel mesh is defined, the serial mesh can be deleted.
203  ParMesh pmesh(MPI_COMM_WORLD, *mesh);
204  delete mesh;
205 
206  // Refine this mesh in parallel to increase the resolution.
207  int par_ref_levels = parallel_ref_levels;
208  for (int l = 0; l < par_ref_levels; l++)
209  {
210  pmesh.UniformRefinement();
211  }
212 
213  // If the gradient bc was selected but the E field was not specified
214  // set a default vector value.
215  if ( dbcg && e_uniform_.Size() != sdim )
216  {
217  e_uniform_.SetSize(sdim);
218  e_uniform_ = 0.0;
219  e_uniform_(sdim-1) = 1.0;
220  }
221 
222  // If values for Dirichlet BCs were not set assume they are zero
223  if ( dbcv.Size() < dbcs.Size() && !dbcg )
224  {
225  dbcv.SetSize(dbcs.Size());
226  dbcv = 0.0;
227  }
228 
229  // If values for Neumann BCs were not set assume they are zero
230  if ( nbcv.Size() < nbcs.Size() )
231  {
232  nbcv.SetSize(nbcs.Size());
233  nbcv = 0.0;
234  }
235 
236  // Create a coefficient describing the dielectric permittivity
238 
239  // Create the Electrostatic solver
240  VoltaSolver Volta(pmesh, order, dbcs, dbcv, nbcs, nbcv, *epsCoef,
241  ( e_uniform_.Size() > 0 ) ? phi_bc_uniform : NULL,
242  ( cs_params_.Size() > 0 ) ? charged_sphere : NULL,
243  ( vp_params_.Size() > 0 ) ? voltaic_pile : NULL);
244 
245  // Initialize GLVis visualization
246  if (visualization)
247  {
248  Volta.InitializeGLVis();
249  }
250 
251  // Initialize VisIt visualization
252  VisItDataCollection visit_dc("Volta-AMR-Parallel", &pmesh);
253 
254  if ( visit )
255  {
256  Volta.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 nodal finite element space
264  // 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  Volta.PrintSizes();
275 
276  // Assemble all forms
277  Volta.Assemble();
278 
279  // Solve the system and compute any auxiliary fields
280  Volta.Solve();
281 
282  // Determine the current size of the linear system
283  int prob_size = Volta.GetProblemSize();
284 
285  // Write fields to disk for VisIt
286  if ( visit )
287  {
288  Volta.WriteVisItFields(it);
289  }
290 
291  // Send the solution by socket to a GLVis server.
292  if (visualization)
293  {
294  Volta.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  Volta.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.7;
342  double threshold = frac * global_max_err;
343  if (mpi.Root()) { cout << "Refining ..." << endl; }
344  pmesh.RefineByError(errors, threshold);
345 
346  // Update the electrostatic solver to reflect the new state of the mesh.
347  Volta.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  Volta.Update();
356  }
357  }
358 
359  delete epsCoef;
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  << " \\ Y / _ \\| |\\ __\\__ \\ " << endl
370  << " \\ ( <_> ) |_| | / __ \\_ " << endl
371  << " \\___/ \\____/|____/__| (____ / " << endl
372  << " \\/ " << endl << flush;
373 }
374 
375 // The Permittivity 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 ( ds_params_.Size() > 0 )
383  {
385  }
386  else if ( pw_eps_.Size() > 0 )
387  {
388  coef = new PWConstCoefficient(pw_eps_);
389  }
390  else
391  {
392  coef = new ConstantCoefficient(epsilon0_);
393  }
394 
395  return coef;
396 }
397 
398 // A sphere with constant permittivity. The sphere has a radius,
399 // center, and permittivity specified on the command line and stored
400 // in ds_params_.
401 double dielectric_sphere(const Vector &x)
402 {
403  double r2 = 0.0;
404 
405  for (int i=0; i<x.Size(); i++)
406  {
407  r2 += (x(i)-ds_params_(i))*(x(i)-ds_params_(i));
408  }
409 
410  if ( sqrt(r2) <= ds_params_(x.Size()) )
411  {
412  return ds_params_(x.Size()+1) * epsilon0_;
413  }
414  return epsilon0_;
415 }
416 
417 // A sphere with constant charge density. The sphere has a radius,
418 // center, and total charge specified on the command line and stored
419 // in cs_params_.
420 double charged_sphere(const Vector &x)
421 {
422  double r2 = 0.0;
423  double rho = 0.0;
424 
425  if ( cs_params_(x.Size()) > 0.0 )
426  {
427  switch ( x.Size() )
428  {
429  case 2:
430  rho = cs_params_(x.Size()+1) /
431  (M_PI * pow(cs_params_(x.Size()), 2));
432  break;
433  case 3:
434  rho = 0.75 * cs_params_(x.Size()+1) /
435  (M_PI * pow(cs_params_(x.Size()), 3));
436  break;
437  default:
438  rho = 0.0;
439  }
440  }
441 
442  for (int i=0; i<x.Size(); i++)
443  {
444  r2 += (x(i) - cs_params_(i)) * (x(i) - cs_params_(i));
445  }
446 
447  if ( sqrt(r2) <= cs_params_(x.Size()) )
448  {
449  return rho;
450  }
451  return 0.0;
452 }
453 
454 // A Cylindrical Rod of constant polarization. The cylinder has two
455 // axis end points, a radius, and a constant electric polarization oriented
456 // along the axis.
457 void voltaic_pile(const Vector &x, Vector &p)
458 {
459  p.SetSize(x.Size());
460  p = 0.0;
461 
462  Vector a(x.Size()); // Normalized Axis vector
463  Vector xu(x.Size()); // x vector relative to the axis end-point
464 
465  xu = x;
466 
467  for (int i=0; i<x.Size(); i++)
468  {
469  xu[i] -= vp_params_[i];
470  a[i] = vp_params_[x.Size()+i] - vp_params_[i];
471  }
472 
473  double h = a.Norml2();
474 
475  if ( h == 0.0 )
476  {
477  return;
478  }
479 
480  double r = vp_params_[2 * x.Size()];
481  double xa = xu * a;
482 
483  if ( h > 0.0 )
484  {
485  xu.Add(-xa / (h * h), a);
486  }
487 
488  double xp = xu.Norml2();
489 
490  if ( xa >= 0.0 && xa <= h*h && xp <= r )
491  {
492  p.Add(vp_params_[2 * x.Size() + 1] / h, a);
493  }
494 }
495 
496 // To produce a uniform electric field the potential can be set
497 // to (- Ex x - Ey y - Ez z).
498 double phi_bc_uniform(const Vector &x)
499 {
500  double phi = 0.0;
501 
502  for (int i=0; i<x.Size(); i++)
503  {
504  phi -= x(i) * e_uniform_(i);
505  }
506 
507  return phi;
508 }
double phi_bc_uniform(const Vector &)
Definition: volta.cpp:498
int WorldSize() const
Return MPI_COMM_WORLD&#39;s size.
int Size() const
Logical size of the array.
Definition: array.hpp:133
Subclass constant coefficient.
Definition: coefficient.hpp:57
void SetSize(int s)
Resize the vector to size s.
Definition: vector.hpp:328
double Norml2() const
Returns the l2 norm of the vector.
Definition: vector.cpp:666
double dielectric_sphere(const Vector &)
Definition: maxwell.cpp:381
int Size() const
Returns the size of the vector.
Definition: vector.hpp:120
int GetNE() const
Returns number of elements.
Definition: mesh.hpp:618
void RegisterVisItFields(VisItDataCollection &visit_dc)
int main(int argc, char *argv[])
Definition: ex1.cpp:45
void voltaic_pile(const Vector &, Vector &)
Definition: volta.cpp:457
void GetErrorEstimates(Vector &errors)
void Rebalance()
Load balance the mesh. NC meshes only.
Definition: pmesh.cpp:2746
bool Nonconforming() const
Definition: mesh.hpp:1014
A simple convenience class that calls MPI_Init() at construction and MPI_Finalize() at destruction...
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
Definition: mesh.cpp:6741
Data collection with VisIt I/O routines.
void SetCurvature(int order, bool discont=false, int space_dim=-1, int ordering=1)
Definition: mesh.cpp:3355
bool Root() const
Return true if WorldRank() == 0.
void PrintUsage(std::ostream &out) const
Definition: optparser.cpp:434
double charged_sphere(const Vector &)
Definition: volta.cpp:420
int SpaceDimension() const
Definition: mesh.hpp:646
MPI_Comm GetComm() const
Definition: pmesh.hpp:123
Base class Coefficient that may optionally depend on time.
Definition: coefficient.hpp:31
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:74
bool RefineByError(const Array< double > &elem_error, double threshold, int nonconforming=-1, int nc_limit=0)
Definition: mesh.cpp:6469
NURBSExtension * NURBSext
Optional NURBS mesh extension.
Definition: mesh.hpp:176
class for piecewise constant coefficient
Definition: coefficient.hpp:72
Vector & Add(const double a, const Vector &Va)
(*this) += a * Va
Definition: vector.cpp:201
void PrintOptions(std::ostream &out) const
Definition: optparser.cpp:304
void EnsureNCMesh(bool triangles_nonconforming=false)
Definition: mesh.cpp:6407
class for C-function coefficient
Vector data type.
Definition: vector.hpp:48
Coefficient * SetupPermittivityCoefficient()
Definition: volta.cpp:378
void display_banner(ostream &os)
Definition: joule.cpp:785
Class for parallel meshes.
Definition: pmesh.hpp:32
bool Good() const
Definition: optparser.hpp:120