MFEM  v4.3.0
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
ex20p.cpp
Go to the documentation of this file.
1 // MFEM Example 20 - Parallel Version
2 //
3 // Compile with: make ex20p
4 //
5 // Sample runs: mpirun -np 4 ex20p
6 // mpirun -np 4 ex20p -p 1 -o 1 -n 120 -dt 0.1
7 // mpirun -np 4 ex20p -p 1 -o 2 -n 60 -dt 0.2
8 // mpirun -np 4 ex20p -p 1 -o 3 -n 40 -dt 0.3
9 // mpirun -np 4 ex20p -p 1 -o 4 -n 30 -dt 0.4
10 //
11 // Description: This example demonstrates the use of the variable order,
12 // symplectic ODE integration algorithm. Symplectic integration
13 // algorithms are designed to conserve energy when integrating, in
14 // time, systems of ODEs which are derived from Hamiltonian
15 // systems.
16 //
17 // Hamiltonian systems define the energy of a system as a function
18 // of time (t), a set of generalized coordinates (q), and their
19 // corresponding generalized momenta (p).
20 //
21 // H(q,p,t) = T(p) + V(q,t)
22 //
23 // Hamilton's equations then specify how q and p evolve in time:
24 //
25 // dq/dt = dH/dp
26 // dp/dt = -dH/dq
27 //
28 // To use the symplectic integration classes we need to define an
29 // mfem::Operator P which evaluates the action of dH/dp, and an
30 // mfem::TimeDependentOperator F which computes -dH/dq.
31 //
32 // This example offers five simple 1D Hamiltonians:
33 // 0) Simple Harmonic Oscillator (mass on a spring)
34 // H = ( p^2 / m + q^2 / k ) / 2
35 // 1) Pendulum
36 // H = ( p^2 / m - k ( 1 - cos(q) ) ) / 2
37 // 2) Gaussian Potential Well
38 // H = ( p^2 / m ) / 2 - k exp(-q^2 / 2)
39 // 3) Quartic Potential
40 // H = ( p^2 / m + k ( 1 + q^2 ) q^2 ) / 2
41 // 4) Negative Quartic Potential
42 // H = ( p^2 / m + k ( 1 - q^2 /8 ) q^2 ) / 2
43 //
44 // In all cases these Hamiltonians are shifted by constant values
45 // so that the energy will remain positive. The mean and standard
46 // deviation of the computed energies at each time step are
47 // displayed upon completion. When run in parallel the same
48 // Hamiltonian system is evolved on each processor but starting
49 // from different initial conditions.
50 //
51 // We then use GLVis to visualize the results in a non-standard way
52 // by defining the axes to be q, p, and t rather than x, y, and z.
53 // In this space we build a ribbon-like mesh on each processor with
54 // nodes at (0,0,t) and (q,p,t). When these ribbons are bonded
55 // together on the t-axis they resemble a Rotini pasta. Finally we
56 // plot the energy as a function of time as a scalar field on this
57 // Rotini-like mesh.
58 //
59 // For a more traditional plot of the results, including q, p, and
60 // H from each processor, can be obtained by selecting the "-gp"
61 // option. This creates a collection of data files and an input
62 // deck for the GnuPlot application (not included with MFEM). To
63 // visualize these results on most linux systems type the command
64 // "gnuplot gnuplot_ex20p.inp". The data files, named
65 // "ex20p_?????.dat", should be simple enough to display with other
66 // plotting programs as well.
67 
68 #include "mfem.hpp"
69 #include <fstream>
70 #include <iostream>
71 
72 using namespace std;
73 using namespace mfem;
74 
75 // Constants used in the Hamiltonian
76 static int prob_ = 0;
77 static double m_ = 1.0;
78 static double k_ = 1.0;
79 
80 // Hamiltonian functional, see below for implementation
81 double hamiltonian(double q, double p, double t);
82 
83 class GradT : public Operator
84 {
85 public:
86  GradT() : Operator(1) {}
87  void Mult(const Vector &x, Vector &y) const { y.Set(1.0/m_, x); }
88 };
89 
90 class NegGradV : public TimeDependentOperator
91 {
92 public:
93  NegGradV() : TimeDependentOperator(1) {}
94  void Mult(const Vector &x, Vector &y) const;
95 };
96 
97 int main(int argc, char *argv[])
98 {
99  // 1. Initialize MPI.
100  int num_procs, myid;
101  MPI_Comm comm = MPI_COMM_WORLD;
102  MPI_Init(&argc, &argv);
103  MPI_Comm_size(comm, &num_procs);
104  MPI_Comm_rank(comm, &myid);
105 
106  // 2. Parse command-line options.
107  int order = 1;
108  int nsteps = 100;
109  double dt = 0.1;
110  bool visualization = true;
111  bool gnuplot = false;
112 
113  OptionsParser args(argc, argv);
114  args.AddOption(&order, "-o", "--order",
115  "Time integration order.");
116  args.AddOption(&prob_, "-p", "--problem-type",
117  "Problem Type:\n"
118  "\t 0 - Simple Harmonic Oscillator\n"
119  "\t 1 - Pendulum\n"
120  "\t 2 - Gaussian Potential Well\n"
121  "\t 3 - Quartic Potential\n"
122  "\t 4 - Negative Quartic Potential");
123  args.AddOption(&nsteps, "-n", "--number-of-steps",
124  "Number of time steps.");
125  args.AddOption(&dt, "-dt", "--time-step",
126  "Time step size.");
127  args.AddOption(&m_, "-m", "--mass",
128  "Mass.");
129  args.AddOption(&k_, "-k", "--spring-const",
130  "Spring constant.");
131  args.AddOption(&visualization, "-vis", "--visualization", "-no-vis",
132  "--no-visualization",
133  "Enable or disable GLVis visualization.");
134  args.AddOption(&gnuplot, "-gp", "--gnuplot", "-no-gp", "--no-gnuplot",
135  "Enable or disable GnuPlot visualization.");
136  args.Parse();
137  if (!args.Good())
138  {
139  if (myid == 0)
140  {
141  args.PrintUsage(cout);
142  }
143  MPI_Finalize();
144  return 1;
145  }
146  if (myid == 0)
147  {
148  args.PrintOptions(cout);
149  }
150 
151  // 3. Create and Initialize the Symplectic Integration Solver
152  SIAVSolver siaSolver(order);
153  GradT P;
154  NegGradV F;
155  siaSolver.Init(P,F);
156 
157  // 4. Set the initial conditions
158  double t = 0.0;
159  Vector q(1), p(1);
160  Vector e(nsteps+1);
161  q(0) = sin(2.0*M_PI*(double)myid/num_procs);
162  p(0) = cos(2.0*M_PI*(double)myid/num_procs);
163 
164  // 5. Prepare GnuPlot output file if needed
165  ostringstream oss;
166  ofstream ofs;
167  if (gnuplot)
168  {
169  oss << "ex20p_" << setfill('0') << setw(5) << myid << ".dat";
170  ofs.open(oss.str().c_str());
171  ofs << t << "\t" << q(0) << "\t" << p(0) << endl;
172  }
173 
174  // 6. Create a Mesh for visualization in phase space
175  int nverts = (visualization) ? 2*num_procs*(nsteps+1) : 0;
176  int nelems = (visualization) ? (nsteps * num_procs) : 0;
177  Mesh mesh(2, nverts, nelems, 0, 3);
178 
179  int *part = (visualization) ? (new int[nelems]) : NULL;
180  int v[4];
181  Vector x0(3); x0 = 0.0;
182  Vector x1(3); x1 = 0.0;
183 
184  // 7. Perform time-stepping
185  double e_mean = 0.0;
186 
187  for (int i = 0; i < nsteps; i++)
188  {
189  // 7a. Record initial state
190  if (i == 0)
191  {
192  e[0] = hamiltonian(q(0),p(0),t);
193  e_mean += e[0];
194 
195  if (visualization)
196  {
197  for (int j = 0; j < num_procs; j++)
198  {
199  mesh.AddVertex(x0);
200  x1[0] = q(0);
201  x1[1] = p(0);
202  x1[2] = 0.0;
203  mesh.AddVertex(x1);
204  }
205  }
206  }
207 
208  // 7b. Advance the state of the system
209  siaSolver.Step(q,p,t,dt);
210  e[i+1] = hamiltonian(q(0),p(0),t);
211  e_mean += e[i+1];
212 
213  // 7c. Record the state of the system
214  if (gnuplot)
215  {
216  ofs << t << "\t" << q(0) << "\t" << p(0) << "\t" << e[i+1] << endl;
217  }
218 
219  // 7d. Add results to GLVis visualization
220  if (visualization)
221  {
222  x0[2] = t;
223  for (int j = 0; j < num_procs; j++)
224  {
225  mesh.AddVertex(x0);
226  x1[0] = q(0);
227  x1[1] = p(0);
228  x1[2] = t;
229  mesh.AddVertex(x1);
230  v[0] = 2 * num_procs * i + 2 * j;
231  v[1] = 2 * num_procs * (i + 1) + 2 * j;
232  v[2] = 2 * num_procs * (i + 1) + 2 * j + 1;
233  v[3] = 2 * num_procs * i + 2 * j + 1;
234  mesh.AddQuad(v);
235  part[num_procs * i + j] = j;
236  }
237  }
238  }
239 
240  // 8. Compute and display mean and standard deviation of the energy
241  e_mean /= (nsteps + 1);
242  double e_var = 0.0;
243  for (int i = 0; i <= nsteps; i++)
244  {
245  e_var += pow(e[i] - e_mean, 2);
246  }
247  e_var /= (nsteps + 1);
248  double e_sd = sqrt(e_var);
249 
250  if (myid == 0)
251  {
252  cout << endl << "Mean and standard deviation of the energy" << endl;
253  }
254  for (int i = 0; i < num_procs; i++)
255  {
256  if (myid == i)
257  {
258  cout << myid << ": " << e_mean << "\t" << e_sd << endl;
259  }
260  MPI_Barrier(comm);
261  }
262 
263  // 9. Finalize the GnuPlot output
264  if (gnuplot)
265  {
266  ofs.close();
267  if (myid == 0)
268  {
269  ofs.open("gnuplot_ex20p.inp");
270  for (int i = 0; i < num_procs; i++)
271  {
272  ostringstream ossi;
273  ossi << "ex20p_" << setfill('0') << setw(5) << i << ".dat";
274  if (i == 0)
275  {
276  ofs << "plot";
277  }
278  ofs << " '" << ossi.str() << "' using 1:2 w l t 'q" << i << "',"
279  << " '" << ossi.str() << "' using 1:3 w l t 'p" << i << "',"
280  << " '" << ossi.str() << "' using 1:4 w l t 'H" << i << "'";
281  if (i < num_procs-1)
282  {
283  ofs << ",";
284  }
285  else
286  {
287  ofs << ";" << endl;
288  }
289  }
290  ofs.close();
291  }
292  }
293 
294  // 10. Finalize the GLVis output
295  if (visualization)
296  {
297  mesh.FinalizeQuadMesh(1);
298  ParMesh pmesh(comm, mesh, part);
299  delete [] part;
300 
301  H1_FECollection fec(order = 1, 2);
302  ParFiniteElementSpace fespace(&pmesh, &fec);
303  ParGridFunction energy(&fespace);
304  energy = 0.0;
305  for (int i = 0; i <= nsteps; i++)
306  {
307  energy[2*i+0] = e[i];
308  energy[2*i+1] = e[i];
309  }
310 
311  char vishost[] = "localhost";
312  int visport = 19916;
313  socketstream sock(vishost, visport);
314  sock.precision(8);
315  sock << "parallel " << num_procs << " " << myid << "\n"
316  << "solution\n" << pmesh << energy
317  << "window_title 'Energy in Phase Space'\n"
318  << "keys\n maac\n" << "axis_labels 'q' 'p' 't'\n"<< flush;
319  }
320 
321  MPI_Finalize();
322 }
323 
324 double hamiltonian(double q, double p, double t)
325 {
326  double h = 1.0 - 0.5 / m_ + 0.5 * p * p / m_;
327  switch (prob_)
328  {
329  case 1:
330  h += k_ * (1.0 - cos(q));
331  break;
332  case 2:
333  h += k_ * (1.0 - exp(-0.5 * q * q));
334  break;
335  case 3:
336  h += 0.5 * k_ * (1.0 + q * q) * q * q;
337  break;
338  case 4:
339  h += 0.5 * k_ * (1.0 - 0.125 * q * q) * q * q;
340  break;
341  default:
342  h += 0.5 * k_ * q * q;
343  break;
344  }
345  return h;
346 }
347 
348 void NegGradV::Mult(const Vector &x, Vector &y) const
349 {
350  switch (prob_)
351  {
352  case 1:
353  y(0) = - k_* sin(x(0));
354  break;
355  case 2:
356  y(0) = - k_ * x(0) * exp(-0.5 * x(0) * x(0));
357  break;
358  case 3:
359  y(0) = - k_ * (1.0 + 2.0 * x(0) * x(0)) * x(0);
360  break;
361  case 4:
362  y(0) = - k_ * (1.0 - 0.25 * x(0) * x(0)) * x(0);
363  break;
364  default:
365  y(0) = - k_ * x(0);
366  break;
367  };
368 }
double hamiltonian(double q, double p, double t)
Definition: ex20.cpp:259
virtual void Init(Operator &P, TimeDependentOperator &F)
Definition: ode.cpp:885
int AddQuad(int v1, int v2, int v3, int v4, int attr=1)
Definition: mesh.cpp:1324
Base abstract class for first order time dependent operators.
Definition: operator.hpp:282
void Mult(const Table &A, const Table &B, Table &C)
C = A * B (as boolean matrices)
Definition: table.cpp:476
Abstract parallel finite element space.
Definition: pfespace.hpp:28
void Step(Vector &q, Vector &p, double &t, double &dt) override
Definition: ode.cpp:964
int AddVertex(double x, double y=0.0, double z=0.0)
Definition: mesh.cpp:1263
void Parse()
Parse the command-line options. Note that this function expects all the options provided through the ...
Definition: optparser.cpp:150
constexpr char vishost[]
constexpr int visport
void PrintUsage(std::ostream &out) const
Print the usage message.
Definition: optparser.cpp:457
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
Vector & Set(const double a, const Vector &x)
(*this) = a * x
Definition: vector.cpp:258
void FinalizeQuadMesh(int generate_edges=0, int refine=0, bool fix_orientation=true)
Finalize the construction of a quadrilateral Mesh.
Definition: mesh.cpp:1585
void PrintOptions(std::ostream &out) const
Print the options.
Definition: optparser.cpp:327
RefCoord t[3]
Vector data type.
Definition: vector.hpp:60
Arbitrary order H1-conforming (continuous) finite elements.
Definition: fe_coll.hpp:216
Class for parallel grid function.
Definition: pgridfunc.hpp:32
Abstract operator.
Definition: operator.hpp:24
Class for parallel meshes.
Definition: pmesh.hpp:32
int main()
Variable order Symplectic Integration Algorithm (orders 1-4)
Definition: ode.hpp:589
bool Good() const
Return true if the command line options were parsed successfully.
Definition: optparser.hpp:150