MFEM  v4.1.0
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
ex22.cpp
Go to the documentation of this file.
1 // MFEM Example 22
2 //
3 // Compile with: make ex22
4 //
5 // Sample runs: ex22 -m ../data/inline-segment.mesh -o 3
6 // ex22 -m ../data/inline-tri.mesh -o 3
7 // ex22 -m ../data/inline-quad.mesh -o 3
8 // ex22 -m ../data/inline-quad.mesh -o 3 -p 1
9 // ex22 -m ../data/inline-quad.mesh -o 3 -p 2
10 // ex22 -m ../data/inline-tet.mesh -o 2
11 // ex22 -m ../data/inline-hex.mesh -o 2
12 // ex22 -m ../data/inline-hex.mesh -o 2 -p 1
13 // ex22 -m ../data/inline-hex.mesh -o 2 -p 2
14 // ex22 -m ../data/star.mesh -r 1 -o 2 -sigma 10.0
15 //
16 // Description: This example code demonstrates the use of MFEM to define and
17 // solve simple complex-valued linear systems. It implements three
18 // variants of a damped harmonic oscillator:
19 //
20 // 1) A scalar H1 field
21 // -Div(a Grad u) - omega^2 b u + i omega c u = 0
22 //
23 // 2) A vector H(Curl) field
24 // Curl(a Curl u) - omega^2 b u + i omega c u = 0
25 //
26 // 3) A vector H(Div) field
27 // -Grad(a Div u) - omega^2 b u + i omega c u = 0
28 //
29 // In each case the field is driven by a forced oscillation, with
30 // angular frequency omega, imposed at the boundary or a portion
31 // of the boundary.
32 //
33 // In electromagnetics, the coefficients are typically named the
34 // permeability, mu = 1/a, permittivity, epsilon = b, and
35 // conductivity, sigma = c. The user can specify these constants
36 // using either set of names.
37 //
38 // The example also demonstrates how to display a time-varying
39 // solution as a sequence of fields sent to a single GLVis socket.
40 //
41 // We recommend viewing examples 1, 3 and 4 before viewing this
42 // example.
43 
44 #include "mfem.hpp"
45 #include <fstream>
46 #include <iostream>
47 
48 using namespace std;
49 using namespace mfem;
50 
51 static double mu_ = 1.0;
52 static double epsilon_ = 1.0;
53 static double sigma_ = 20.0;
54 static double omega_ = 10.0;
55 
56 double u0_real_exact(const Vector &);
57 double u0_imag_exact(const Vector &);
58 
59 void u1_real_exact(const Vector &, Vector &);
60 void u1_imag_exact(const Vector &, Vector &);
61 
62 void u2_real_exact(const Vector &, Vector &);
63 void u2_imag_exact(const Vector &, Vector &);
64 
65 bool check_for_inline_mesh(const char * mesh_file);
66 
67 int main(int argc, char *argv[])
68 {
69  // 1. Parse command-line options.
70  const char *mesh_file = "../data/inline-quad.mesh";
71  int ref_levels = 0;
72  int order = 1;
73  int prob = 0;
74  double freq = -1.0;
75  double a_coef = 0.0;
76  bool visualization = 1;
77  bool herm_conv = true;
78  bool exact_sol = true;
79 
80  OptionsParser args(argc, argv);
81  args.AddOption(&mesh_file, "-m", "--mesh",
82  "Mesh file to use.");
83  args.AddOption(&ref_levels, "-r", "--refine",
84  "Number of times to refine the mesh uniformly.");
85  args.AddOption(&order, "-o", "--order",
86  "Finite element order (polynomial degree).");
87  args.AddOption(&prob, "-p", "--problem-type",
88  "Choose between 0: H_1, 1: H(Curl), or 2: H(Div) "
89  "damped harmonic oscillator.");
90  args.AddOption(&a_coef, "-a", "--stiffness-coef",
91  "Stiffness coefficient (spring constant or 1/mu).");
92  args.AddOption(&epsilon_, "-b", "--mass-coef",
93  "Mass coefficient (or epsilon).");
94  args.AddOption(&sigma_, "-c", "--damping-coef",
95  "Damping coefficient (or sigma).");
96  args.AddOption(&mu_, "-mu", "--permeability",
97  "Permeability of free space (or 1/(spring constant)).");
98  args.AddOption(&epsilon_, "-eps", "--permittivity",
99  "Permittivity of free space (or mass constant).");
100  args.AddOption(&sigma_, "-sigma", "--conductivity",
101  "Conductivity (or damping constant).");
102  args.AddOption(&freq, "-f", "--frequency",
103  "Frequency (in Hz).");
104  args.AddOption(&herm_conv, "-herm", "--hermitian", "-no-herm",
105  "--no-hermitian", "Use convention for Hermitian operators.");
106  args.AddOption(&visualization, "-vis", "--visualization", "-no-vis",
107  "--no-visualization",
108  "Enable or disable GLVis visualization.");
109  args.Parse();
110  if (!args.Good())
111  {
112  args.PrintUsage(cout);
113  return 1;
114  }
115  args.PrintOptions(cout);
116 
117  MFEM_VERIFY(prob >= 0 && prob <=2,
118  "Unrecognized problem type: " << prob);
119 
120  if ( a_coef != 0.0 )
121  {
122  mu_ = 1.0 / a_coef;
123  }
124  if ( freq > 0.0 )
125  {
126  omega_ = 2.0 * M_PI * freq;
127  }
128 
129  exact_sol = check_for_inline_mesh(mesh_file);
130  if (exact_sol)
131  {
132  cout << "Identified a mesh with known exact solution" << endl;
133  }
134 
136  herm_conv ? ComplexOperator::HERMITIAN : ComplexOperator::BLOCK_SYMMETRIC;
137 
138  // 2. Read the mesh from the given mesh file. We can handle triangular,
139  // quadrilateral, tetrahedral, hexahedral, surface and volume meshes
140  // with the same code.
141  Mesh *mesh = new Mesh(mesh_file, 1, 1);
142  int dim = mesh->Dimension();
143 
144  // 3. Refine the mesh to increase resolution. In this example we do
145  // 'ref_levels' of uniform refinement where the user specifies
146  // the number of levels with the '-r' option.
147  for (int l = 0; l < ref_levels; l++)
148  {
149  mesh->UniformRefinement();
150  }
151 
152  // 4. Define a finite element space on the mesh. Here we use continuous
153  // Lagrange, Nedelec, or Raviart-Thomas finite elements of the specified
154  // order.
155  if (dim == 1 && prob != 0 )
156  {
157  cout << "Switching to problem type 0, H1 basis functions, "
158  << "for 1 dimensional mesh." << endl;
159  prob = 0;
160  }
161 
162  FiniteElementCollection *fec = NULL;
163  switch (prob)
164  {
165  case 0: fec = new H1_FECollection(order, dim); break;
166  case 1: fec = new ND_FECollection(order, dim); break;
167  case 2: fec = new RT_FECollection(order - 1, dim); break;
168  default: break; // This should be unreachable
169  }
170  FiniteElementSpace *fespace = new FiniteElementSpace(mesh, fec);
171  cout << "Number of finite element unknowns: " << fespace->GetTrueVSize()
172  << endl;
173 
174  // 5. Determine the list of true (i.e. conforming) essential boundary dofs.
175  // In this example, the boundary conditions are defined based on the type
176  // of mesh and the problem type.
177  Array<int> ess_tdof_list;
178  Array<int> ess_bdr;
179  if (mesh->bdr_attributes.Size())
180  {
181  ess_bdr.SetSize(mesh->bdr_attributes.Max());
182  ess_bdr = 1;
183  fespace->GetEssentialTrueDofs(ess_bdr, ess_tdof_list);
184  }
185 
186  // 6. Set up the linear form b(.) which corresponds to the right-hand side of
187  // the FEM linear system.
188  ComplexLinearForm b(fespace, conv);
189  b.Vector::operator=(0.0);
190 
191  // 7. Define the solution vector u as a complex finite element grid function
192  // corresponding to fespace. Initialize u with initial guess of 1+0i or
193  // the exact solution if it is known.
194  ComplexGridFunction u(fespace);
195  ComplexGridFunction * u_exact = NULL;
196  if (exact_sol) { u_exact = new ComplexGridFunction(fespace); }
197 
204 
205  ConstantCoefficient zeroCoef(0.0);
206  ConstantCoefficient oneCoef(1.0);
207 
208  Vector zeroVec(dim); zeroVec = 0.0;
209  Vector oneVec(dim); oneVec = 0.0; oneVec[(prob==2)?(dim-1):0] = 1.0;
210  VectorConstantCoefficient zeroVecCoef(zeroVec);
211  VectorConstantCoefficient oneVecCoef(oneVec);
212 
213  u = 0.0;
214  switch (prob)
215  {
216  case 0:
217  if (exact_sol)
218  {
219  u.ProjectBdrCoefficient(u0_r, u0_i, ess_bdr);
220  u_exact->ProjectCoefficient(u0_r, u0_i);
221  }
222  else
223  {
224  u.ProjectBdrCoefficient(oneCoef, zeroCoef, ess_bdr);
225  }
226  break;
227  case 1:
228  if (exact_sol)
229  {
230  u.ProjectBdrCoefficientTangent(u1_r, u1_i, ess_bdr);
231  u_exact->ProjectCoefficient(u1_r, u1_i);
232  }
233  else
234  {
235  u.ProjectBdrCoefficientTangent(oneVecCoef, zeroVecCoef, ess_bdr);
236  }
237  break;
238  case 2:
239  if (exact_sol)
240  {
241  u.ProjectBdrCoefficientNormal(u2_r, u2_i, ess_bdr);
242  u_exact->ProjectCoefficient(u2_r, u2_i);
243  }
244  else
245  {
246  u.ProjectBdrCoefficientNormal(oneVecCoef, zeroVecCoef, ess_bdr);
247  }
248  break;
249  default: break; // This should be unreachable
250  }
251 
252  if (visualization && exact_sol)
253  {
254  char vishost[] = "localhost";
255  int visport = 19916;
256  socketstream sol_sock_r(vishost, visport);
257  socketstream sol_sock_i(vishost, visport);
258  sol_sock_r.precision(8);
259  sol_sock_i.precision(8);
260  sol_sock_r << "solution\n" << *mesh << u_exact->real()
261  << "window_title 'Exact: Real Part'" << flush;
262  sol_sock_i << "solution\n" << *mesh << u_exact->imag()
263  << "window_title 'Exact: Imaginary Part'" << flush;
264  }
265 
266  // 8. Set up the sesquilinear form a(.,.) on the finite element space
267  // corresponding to the damped harmonic oscillator operator of the
268  // appropriate type:
269  //
270  // 0) A scalar H1 field
271  // -Div(a Grad) - omega^2 b + i omega c
272  //
273  // 1) A vector H(Curl) field
274  // Curl(a Curl) - omega^2 b + i omega c
275  //
276  // 2) A vector H(Div) field
277  // -Grad(a Div) - omega^2 b + i omega c
278  //
279  ConstantCoefficient stiffnessCoef(1.0/mu_);
280  ConstantCoefficient massCoef(-omega_ * omega_ * epsilon_);
281  ConstantCoefficient lossCoef(omega_ * sigma_);
282  ConstantCoefficient negMassCoef(omega_ * omega_ * epsilon_);
283 
284  SesquilinearForm *a = new SesquilinearForm(fespace, conv);
285  switch (prob)
286  {
287  case 0:
288  a->AddDomainIntegrator(new DiffusionIntegrator(stiffnessCoef),
289  NULL);
290  a->AddDomainIntegrator(new MassIntegrator(massCoef),
291  new MassIntegrator(lossCoef));
292  break;
293  case 1:
294  a->AddDomainIntegrator(new CurlCurlIntegrator(stiffnessCoef),
295  NULL);
297  new VectorFEMassIntegrator(lossCoef));
298  break;
299  case 2:
300  a->AddDomainIntegrator(new DivDivIntegrator(stiffnessCoef),
301  NULL);
303  new VectorFEMassIntegrator(lossCoef));
304  break;
305  default: break; // This should be unreachable
306  }
307 
308  // 8a. Set up the bilinear form for the preconditioner corresponding to the
309  // appropriate operator
310  //
311  // 0) A scalar H1 field
312  // -Div(a Grad) - omega^2 b + omega c
313  //
314  // 1) A vector H(Curl) field
315  // Curl(a Curl) + omega^2 b + omega c
316  //
317  // 2) A vector H(Div) field
318  // -Grad(a Div) - omega^2 b + omega c
319  //
320  BilinearForm *pcOp = new BilinearForm(fespace);
321  switch (prob)
322  {
323  case 0:
324  pcOp->AddDomainIntegrator(new DiffusionIntegrator(stiffnessCoef));
325  pcOp->AddDomainIntegrator(new MassIntegrator(massCoef));
326  pcOp->AddDomainIntegrator(new MassIntegrator(lossCoef));
327  break;
328  case 1:
329  pcOp->AddDomainIntegrator(new CurlCurlIntegrator(stiffnessCoef));
330  pcOp->AddDomainIntegrator(new VectorFEMassIntegrator(negMassCoef));
331  pcOp->AddDomainIntegrator(new VectorFEMassIntegrator(lossCoef));
332  break;
333  case 2:
334  pcOp->AddDomainIntegrator(new DivDivIntegrator(stiffnessCoef));
335  pcOp->AddDomainIntegrator(new VectorFEMassIntegrator(massCoef));
336  pcOp->AddDomainIntegrator(new VectorFEMassIntegrator(lossCoef));
337  break;
338  default: break; // This should be unreachable
339  }
340 
341  // 9. Assemble the form and the corresponding linear system, applying any
342  // necessary transformations such as: assembly, eliminating boundary
343  // conditions, conforming constraints for non-conforming AMR, etc.
344  a->Assemble();
345  pcOp->Assemble();
346 
347  OperatorHandle A;
348  Vector B, U;
349 
350  a->FormLinearSystem(ess_tdof_list, u, b, A, U, B);
351  u = 0.0;
352  U = 0.0;
353 
354  OperatorHandle PCOp;
355  pcOp->FormSystemMatrix(ess_tdof_list, PCOp);
356 
357  {
358  ComplexSparseMatrix * Asp =
359  dynamic_cast<ComplexSparseMatrix*>(A.Ptr());
360 
361  cout << "Size of linear system: "
362  << 2 * Asp->real().Width() << endl << endl;
363  }
364 
365  // 10. Define and apply a GMRES solver for AU=B with a block diagonal
366  // preconditioner based on the appropriate sparse smoother.
367  {
368  Array<int> blockOffsets;
369  blockOffsets.SetSize(3);
370  blockOffsets[0] = 0;
371  blockOffsets[1] = PCOp.Ptr()->Height();
372  blockOffsets[2] = PCOp.Ptr()->Height();
373  blockOffsets.PartialSum();
374 
375  BlockDiagonalPreconditioner BDP(blockOffsets);
376 
377  Operator * pc_r = NULL;
378  Operator * pc_i = NULL;
379 
380  double s = 1.0;
381  switch (prob)
382  {
383  case 0:
384  pc_r = new DSmoother(*PCOp.As<SparseMatrix>());
385  break;
386  case 1:
387  pc_r = new GSSmoother(*PCOp.As<SparseMatrix>());
388  s = -1.0;
389  break;
390  case 2:
391  pc_r = new DSmoother(*PCOp.As<SparseMatrix>());
392  break;
393 
394  default: break; // This should be unreachable
395  }
396  pc_i = new ScaledOperator(pc_r,
397  (conv == ComplexOperator::HERMITIAN) ?
398  s:-s);
399 
400  BDP.SetDiagonalBlock(0, pc_r);
401  BDP.SetDiagonalBlock(1, pc_i);
402  BDP.owns_blocks = 1;
403 
404  GMRESSolver gmres;
405  gmres.SetPreconditioner(BDP);
406  gmres.SetOperator(*A.Ptr());
407  gmres.SetRelTol(1e-12);
408  gmres.SetMaxIter(1000);
409  gmres.SetPrintLevel(1);
410  gmres.Mult(B, U);
411  }
412 
413  // 11. Recover the solution as a finite element grid function and compute the
414  // errors if the exact solution is known.
415  a->RecoverFEMSolution(U, b, u);
416 
417  if (exact_sol)
418  {
419  double err_r = -1.0;
420  double err_i = -1.0;
421 
422  switch (prob)
423  {
424  case 0:
425  err_r = u.real().ComputeL2Error(u0_r);
426  err_i = u.imag().ComputeL2Error(u0_i);
427  break;
428  case 1:
429  err_r = u.real().ComputeL2Error(u1_r);
430  err_i = u.imag().ComputeL2Error(u1_i);
431  break;
432  case 2:
433  err_r = u.real().ComputeL2Error(u2_r);
434  err_i = u.imag().ComputeL2Error(u2_i);
435  break;
436  default: break; // This should be unreachable
437  }
438 
439  cout << endl;
440  cout << "|| Re (u_h - u) ||_{L^2} = " << err_r << endl;
441  cout << "|| Im (u_h - u) ||_{L^2} = " << err_i << endl;
442  cout << endl;
443  }
444 
445  // 12. Save the refined mesh and the solution. This output can be viewed
446  // later using GLVis: "glvis -m mesh -g sol".
447  {
448  ofstream mesh_ofs("refined.mesh");
449  mesh_ofs.precision(8);
450  mesh->Print(mesh_ofs);
451 
452  ofstream sol_r_ofs("sol_r.gf");
453  ofstream sol_i_ofs("sol_i.gf");
454  sol_r_ofs.precision(8);
455  sol_i_ofs.precision(8);
456  u.real().Save(sol_r_ofs);
457  u.imag().Save(sol_i_ofs);
458  }
459 
460  // 13. Send the solution by socket to a GLVis server.
461  if (visualization)
462  {
463  char vishost[] = "localhost";
464  int visport = 19916;
465  socketstream sol_sock_r(vishost, visport);
466  socketstream sol_sock_i(vishost, visport);
467  sol_sock_r.precision(8);
468  sol_sock_i.precision(8);
469  sol_sock_r << "solution\n" << *mesh << u.real()
470  << "window_title 'Solution: Real Part'" << flush;
471  sol_sock_i << "solution\n" << *mesh << u.imag()
472  << "window_title 'Solution: Imaginary Part'" << flush;
473  }
474  if (visualization && exact_sol)
475  {
476  *u_exact -= u;
477 
478  char vishost[] = "localhost";
479  int visport = 19916;
480  socketstream sol_sock_r(vishost, visport);
481  socketstream sol_sock_i(vishost, visport);
482  sol_sock_r.precision(8);
483  sol_sock_i.precision(8);
484  sol_sock_r << "solution\n" << *mesh << u_exact->real()
485  << "window_title 'Error: Real Part'" << flush;
486  sol_sock_i << "solution\n" << *mesh << u_exact->imag()
487  << "window_title 'Error: Imaginary Part'" << flush;
488  }
489  if (visualization)
490  {
491  GridFunction u_t(fespace);
492  u_t = u.real();
493  char vishost[] = "localhost";
494  int visport = 19916;
495  socketstream sol_sock(vishost, visport);
496  sol_sock.precision(8);
497  sol_sock << "solution\n" << *mesh << u_t
498  << "window_title 'Harmonic Solution (t = 0.0 T)'"
499  << "pause\n" << flush;
500 
501  cout << "GLVis visualization paused."
502  << " Press space (in the GLVis window) to resume it.\n";
503  int num_frames = 32;
504  int i = 0;
505  while (sol_sock)
506  {
507  double t = (double)(i % num_frames) / num_frames;
508  ostringstream oss;
509  oss << "Harmonic Solution (t = " << t << " T)";
510 
511  add(cos( 2.0 * M_PI * t), u.real(),
512  sin(-2.0 * M_PI * t), u.imag(), u_t);
513  sol_sock << "solution\n" << *mesh << u_t
514  << "window_title '" << oss.str() << "'" << flush;
515  i++;
516  }
517  }
518 
519  // 14. Free the used memory.
520  delete a;
521  delete u_exact;
522  delete pcOp;
523  delete fespace;
524  delete fec;
525  delete mesh;
526 
527  return 0;
528 }
529 
530 bool check_for_inline_mesh(const char * mesh_file)
531 {
532  string file(mesh_file);
533  size_t p0 = file.find_last_of("/");
534  string s0 = file.substr((p0==string::npos)?0:(p0+1),7);
535  return s0 == "inline-";
536 }
537 
538 complex<double> u0_exact(const Vector &x)
539 {
540  int dim = x.Size();
541  complex<double> i(0.0, 1.0);
542  complex<double> alpha = (epsilon_ * omega_ - i * sigma_);
543  complex<double> kappa = std::sqrt(mu_ * omega_* alpha);
544  return std::exp(-i * kappa * x[dim - 1]);
545 }
546 
547 double u0_real_exact(const Vector &x)
548 {
549  return u0_exact(x).real();
550 }
551 
552 double u0_imag_exact(const Vector &x)
553 {
554  return u0_exact(x).imag();
555 }
556 
557 void u1_real_exact(const Vector &x, Vector &v)
558 {
559  int dim = x.Size();
560  v.SetSize(dim); v = 0.0; v[0] = u0_real_exact(x);
561 }
562 
563 void u1_imag_exact(const Vector &x, Vector &v)
564 {
565  int dim = x.Size();
566  v.SetSize(dim); v = 0.0; v[0] = u0_imag_exact(x);
567 }
568 
569 void u2_real_exact(const Vector &x, Vector &v)
570 {
571  int dim = x.Size();
572  v.SetSize(dim); v = 0.0; v[dim-1] = u0_real_exact(x);
573 }
574 
575 void u2_imag_exact(const Vector &x, Vector &v)
576 {
577  int dim = x.Size();
578  v.SetSize(dim); v = 0.0; v[dim-1] = u0_imag_exact(x);
579 }
int Size() const
Logical size of the array.
Definition: array.hpp:124
virtual void Print(std::ostream &out=mfem::out) const
Definition: mesh.hpp:1188
OpType * As() const
Return the Operator pointer statically cast to a specified OpType. Similar to the method Get()...
Definition: handle.hpp:96
virtual void ProjectBdrCoefficientTangent(VectorCoefficient &real_coeff, VectorCoefficient &imag_coeff, Array< int > &attr)
Class for grid function - Vector with associated FE space.
Definition: gridfunc.hpp:27
Data type for scaled Jacobi-type smoother of sparse matrix.
void Assemble(int skip_zeros=1)
Assembles the form i.e. sums over all domain/bdr integrators.
Subclass constant coefficient.
Definition: coefficient.hpp:67
virtual double ComputeL2Error(Coefficient &exsol, const IntegrationRule *irs[]=NULL) const
Definition: gridfunc.hpp:318
void u2_imag_exact(const Vector &, Vector &)
Definition: ex22.cpp:575
void SetSize(int s)
Resize the vector to size s.
Definition: vector.hpp:407
virtual void ProjectCoefficient(Coefficient &real_coeff, Coefficient &imag_coeff)
Definition: complex_fem.cpp:76
Integrator for (curl u, curl v) for Nedelec elements.
int Width() const
Get the width (size of input) of the Operator. Synonym with NumCols().
Definition: operator.hpp:63
Pointer to an Operator of a specified type.
Definition: handle.hpp:33
int Size() const
Returns the size of the vector.
Definition: vector.hpp:157
virtual void GetEssentialTrueDofs(const Array< int > &bdr_attr_is_ess, Array< int > &ess_tdof_list, int component=-1)
Definition: fespace.cpp:366
complex< double > u0_exact(const Vector &x)
Definition: ex22.cpp:538
int main(int argc, char *argv[])
Definition: ex1.cpp:62
(Q div u, div v) for RT elements
void u1_imag_exact(const Vector &, Vector &)
Definition: ex22.cpp:563
void add(const Vector &v1, const Vector &v2, Vector &v)
Definition: vector.cpp:238
Data type for Gauss-Seidel smoother of sparse matrix.
virtual void ProjectBdrCoefficientNormal(VectorCoefficient &real_coeff, VectorCoefficient &imag_coeff, Array< int > &attr)
virtual void Save(std::ostream &out) const
Save the GridFunction to an output stream.
Definition: gridfunc.cpp:2731
virtual void Mult(const Vector &b, Vector &x) const
Operator application: y=A(x).
Definition: solvers.cpp:594
void SetPrintLevel(int print_lvl)
Definition: solvers.cpp:70
virtual void SetOperator(const Operator &op)
Also calls SetOperator for the preconditioner if there is one.
Definition: solvers.cpp:97
A class to handle Block diagonal preconditioners in a matrix-free implementation. ...
Data type sparse matrix.
Definition: sparsemat.hpp:40
int Height() const
Get the height (size of output) of the Operator. Synonym with NumRows().
Definition: operator.hpp:57
double b
Definition: lissajous.cpp:42
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
Definition: mesh.cpp:7982
double u0_imag_exact(const Vector &)
Definition: ex22.cpp:552
virtual void ProjectBdrCoefficient(Coefficient &real_coeff, Coefficient &imag_coeff, Array< int > &attr)
Definition: complex_fem.cpp:92
void SetMaxIter(int max_it)
Definition: solvers.hpp:65
T Max() const
Find the maximal element in the array, using the comparison operator &lt; for class T.
Definition: array.cpp:68
virtual void FormSystemMatrix(const Array< int > &ess_tdof_list, OperatorHandle &A)
Form the linear system matrix A, see FormLinearSystem() for details.
virtual int GetTrueVSize() const
Return the number of vector true (conforming) dofs.
Definition: fespace.hpp:384
GridFunction & real()
Definition: complex_fem.hpp:69
int Dimension() const
Definition: mesh.hpp:744
void PrintUsage(std::ostream &out) const
Definition: optparser.cpp:434
Operator * Ptr() const
Access the underlying Operator pointer.
Definition: handle.hpp:82
Arbitrary order H(div)-conforming Raviart-Thomas finite elements.
Definition: fe_coll.hpp:191
virtual SparseMatrix & real()
Real or imaginary part accessor methods.
Array< int > bdr_attributes
A list of all unique boundary attributes used by the Mesh.
Definition: mesh.hpp:189
void SetRelTol(double rtol)
Definition: solvers.hpp:63
Scaled Operator B: x -&gt; a A(x).
Definition: operator.hpp:544
Specialization of the ComplexOperator built from a pair of Sparse Matrices.
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Definition: fespace.hpp:87
GridFunction & imag()
Definition: complex_fem.hpp:70
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
void SetSize(int nsize)
Change logical size of the array, keep existing entries.
Definition: array.hpp:635
void PartialSum()
Partial Sum.
Definition: array.cpp:103
GMRES method.
Definition: solvers.hpp:184
virtual void RecoverFEMSolution(const Vector &X, const Vector &b, Vector &x)
double a
Definition: lissajous.cpp:41
void FormLinearSystem(const Array< int > &ess_tdof_list, Vector &x, Vector &b, OperatorHandle &A, Vector &X, Vector &B, int copy_interior=0)
void Assemble(int skip_zeros=1)
Assemble the local matrix.
int dim
Definition: ex24.cpp:43
void AddDomainIntegrator(BilinearFormIntegrator *bfi)
Adds new Domain Integrator. Assumes ownership of bfi.
void PrintOptions(std::ostream &out) const
Definition: optparser.cpp:304
void AddDomainIntegrator(BilinearFormIntegrator *bfi_real, BilinearFormIntegrator *bfi_imag)
Adds new Domain Integrator.
const double alpha
Definition: ex15.cpp:336
class for C-function coefficient
Arbitrary order H(curl)-conforming Nedelec finite elements.
Definition: fe_coll.hpp:251
double kappa
Definition: ex3.cpp:54
Vector data type.
Definition: vector.hpp:48
void u1_real_exact(const Vector &, Vector &)
Definition: ex22.cpp:557
virtual void SetPreconditioner(Solver &pr)
This should be called before SetOperator.
Definition: solvers.cpp:91
Arbitrary order H1-conforming (continuous) finite elements.
Definition: fe_coll.hpp:83
Abstract operator.
Definition: operator.hpp:24
bool check_for_inline_mesh(const char *mesh_file)
Definition: ex22.cpp:530
double freq
Definition: ex3.cpp:54
void SetDiagonalBlock(int iblock, Operator *op)
Add a square block op in the block-entry (iblock, iblock).
void u2_real_exact(const Vector &, Vector &)
Definition: ex22.cpp:569
bool Good() const
Definition: optparser.hpp:125
double u0_real_exact(const Vector &)
Definition: ex22.cpp:547