MFEM  v4.6.0
Finite element discretization library
acoustics.cpp
Go to the documentation of this file.
1 // Copyright (c) 2010-2023, 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 // MFEM Ultraweak DPG example for acoustics (Helmholtz)
13 //
14 // Compile with: make acoustics
15 //
16 // Sample runs
17 //
18 // acoustics -ref 5 -o 1 -rnum 1.0
19 // acoustics -m ../../data/inline-tri.mesh -ref 4 -o 2 -sc -rnum 3.0
20 // acoustics -m ../../data/amr-quad.mesh -ref 3 -o 3 -sc -rnum 4.5 -prob 1
21 // acoustics -m ../../data/inline-quad.mesh -ref 2 -o 4 -sc -rnum 11.5 -prob 1
22 // acoustics -m ../../data/inline-hex.mesh -ref 2 -o 2 -sc -rnum 1.0
23 
24 // Description:
25 // This example code demonstrates the use of MFEM to define and solve
26 // the "ultraweak" (UW) DPG formulation for the Helmholtz problem
27 
28 // - Δ p - ω² p = f̃ , in Ω
29 // p = p₀, on ∂Ω
30 
31 // It solves two kinds of problems
32 // a) f̃ = 0 and p₀ is a plane wave
33 // b) A manufactured solution problem where p_exact is a gaussian beam
34 // This example computes and prints out convergence rates for the L² error.
35 
36 // The DPG UW deals with the First Order System
37 // ∇ p + i ω u = 0, in Ω
38 // ∇⋅u + i ω p = f, in Ω (1)
39 // p = p_0, in ∂Ω
40 // where f:=f̃/(i ω)
41 
42 // Ultraweak-DPG is obtained by integration by parts of both equations and the
43 // introduction of trace unknowns on the mesh skeleton
44 //
45 // p ∈ L²(Ω), u ∈ (L²(Ω))ᵈⁱᵐ
46 // p̂ ∈ H^1/2(Ω), û ∈ H^-1/2(Ω)
47 // -(p, ∇⋅v) + i ω (u , v) + < p̂, v⋅n> = 0, ∀ v ∈ H(div,Ω)
48 // -(u , ∇ q) + i ω (p , q) + < û, q > = (f,q) ∀ q ∈ H¹(Ω)
49 // p̂ = p₀ on ∂Ω
50 
51 // Note:
52 // p̂ := p, û := u on the mesh skeleton
53 
54 // For more information see https://doi.org/10.1016/j.camwa.2017.06.044
55 
56 // -------------------------------------------------------------
57 // | | p | u | p̂ | û | RHS |
58 // -------------------------------------------------------------
59 // | v | -(p, ∇⋅v) | i ω (u,v) | < p̂, v⋅n> | | |
60 // | | | | | | |
61 // | q | i ω (p,q) |-(u , ∇ q) | | < û,q > | (f,q) |
62 
63 // where (q,v) ∈ H¹(Ω) × H(div,Ω)
64 
65 // Here we use the "Adjoint Graph" norm on the test space i.e.,
66 // ||(q,v)||²ᵥ = ||A^*(q,v)||² + ||(q,v)||² where A is the
67 // acoustics operator defined by (1)
68 
69 #include "mfem.hpp"
70 #include "util/complexweakform.hpp"
71 #include "../common/mfem-common.hpp"
72 #include <fstream>
73 #include <iostream>
74 
75 using namespace std;
76 using namespace mfem;
77 using namespace mfem::common;
78 
79 complex<double> acoustics_solution(const Vector & X);
80 void acoustics_solution_grad(const Vector & X,vector<complex<double>> &dp);
81 complex<double> acoustics_solution_laplacian(const Vector & X);
82 
83 double p_exact_r(const Vector &x);
84 double p_exact_i(const Vector &x);
85 void u_exact_r(const Vector &x, Vector & u);
86 void u_exact_i(const Vector &x, Vector & u);
87 double rhs_func_r(const Vector &x);
88 double rhs_func_i(const Vector &x);
89 void gradp_exact_r(const Vector &x, Vector &gradu);
90 void gradp_exact_i(const Vector &x, Vector &gradu);
91 double divu_exact_r(const Vector &x);
92 double divu_exact_i(const Vector &x);
93 double d2_exact_r(const Vector &x);
94 double d2_exact_i(const Vector &x);
95 double hatp_exact_r(const Vector & X);
96 double hatp_exact_i(const Vector & X);
97 void hatu_exact_r(const Vector & X, Vector & hatu);
98 void hatu_exact_i(const Vector & X, Vector & hatu);
99 
100 int dim;
101 double omega;
102 
104 {
107 };
108 
110 
111 int main(int argc, char *argv[])
112 {
113  const char *mesh_file = "../../data/inline-quad.mesh";
114  int order = 1;
115  int delta_order = 1;
116  bool visualization = true;
117  double rnum=1.0;
118  int ref = 0;
119  bool static_cond = false;
120  int iprob = 0;
121 
122  OptionsParser args(argc, argv);
123  args.AddOption(&mesh_file, "-m", "--mesh",
124  "Mesh file to use.");
125  args.AddOption(&order, "-o", "--order",
126  "Finite element order (polynomial degree)");
127  args.AddOption(&visualization, "-vis", "--visualization", "-no-vis",
128  "--no-visualization",
129  "Enable or disable GLVis visualization.");
130  args.AddOption(&rnum, "-rnum", "--number_of_wavelengths",
131  "Number of wavelengths");
132  args.AddOption(&iprob, "-prob", "--problem", "Problem case"
133  " 0: plane wave, 1: Gaussian beam");
134  args.AddOption(&delta_order, "-do", "--delta_order",
135  "Order enrichment for DPG test space.");
136  args.AddOption(&ref, "-ref", "--refinements",
137  "Number of serial refinements.");
138  args.AddOption(&static_cond, "-sc", "--static-condensation", "-no-sc",
139  "--no-static-condensation", "Enable static condensation.");
140  args.Parse();
141  if (!args.Good())
142  {
143  args.PrintUsage(cout);
144  return 1;
145  }
146  args.PrintOptions(cout);
147 
148  if (iprob > 1) { iprob = 0; }
149  prob = (prob_type)iprob;
150 
151  omega = 2.*M_PI*rnum;
152 
153  Mesh mesh(mesh_file, 1, 1);
154  dim = mesh.Dimension();
155  MFEM_VERIFY(dim > 1, "Dimension = 1 is not supported in this example");
156 
157  // Define spaces
158  enum TrialSpace
159  {
160  p_space = 0,
161  u_space = 1,
162  hatp_space = 2,
163  hatu_space = 3
164  };
165  enum TestSpace
166  {
167  q_space = 0,
168  v_space = 1
169  };
170 
171  // L2 space for p
172  FiniteElementCollection *p_fec = new L2_FECollection(order-1,dim);
173  FiniteElementSpace *p_fes = new FiniteElementSpace(&mesh,p_fec);
174 
175  // Vector L2 space for u
176  FiniteElementCollection *u_fec = new L2_FECollection(order-1,dim);
177  FiniteElementSpace *u_fes = new FiniteElementSpace(&mesh,u_fec, dim);
178 
179  // H^1/2 space for p̂
180  FiniteElementCollection * hatp_fec = new H1_Trace_FECollection(order,dim);
181  FiniteElementSpace *hatp_fes = new FiniteElementSpace(&mesh,hatp_fec);
182 
183  // H^-1/2 space for û
184  FiniteElementCollection * hatu_fec = new RT_Trace_FECollection(order-1,dim);
185  FiniteElementSpace *hatu_fes = new FiniteElementSpace(&mesh,hatu_fec);
186 
187  // testspace fe collections
188  int test_order = order+delta_order;
189  FiniteElementCollection * q_fec = new H1_FECollection(test_order, dim);
190  FiniteElementCollection * v_fec = new RT_FECollection(test_order-1, dim);
191 
192  // Coefficients
193  ConstantCoefficient one(1.0);
194  ConstantCoefficient zero(0.0);
195  Vector vec0(dim); vec0 = 0.;
196  VectorConstantCoefficient vzero(vec0);
197  ConstantCoefficient negone(-1.0);
200  ConstantCoefficient negomeg(-omega);
201 
204 
205  trial_fes.Append(p_fes);
206  trial_fes.Append(u_fes);
207  trial_fes.Append(hatp_fes);
208  trial_fes.Append(hatu_fes);
209 
210  test_fec.Append(q_fec);
211  test_fec.Append(v_fec);
212 
213  ComplexDPGWeakForm * a = new ComplexDPGWeakForm(trial_fes,test_fec);
214 
215  // i ω (p,q)
216  a->AddTrialIntegrator(nullptr,new MixedScalarMassIntegrator(omeg),
217  TrialSpace::p_space,TestSpace::q_space);
218  // -(u , ∇ q)
219  a->AddTrialIntegrator(new TransposeIntegrator(new GradientIntegrator(negone)),
220  nullptr,TrialSpace::u_space,TestSpace::q_space);
221  // -(p, ∇⋅v)
222  a->AddTrialIntegrator(new MixedScalarWeakGradientIntegrator(one),nullptr,
223  TrialSpace::p_space,TestSpace::v_space);
224  // i ω (u,v)
225  a->AddTrialIntegrator(nullptr,
227  TrialSpace::u_space,TestSpace::v_space);
228  // < p̂, v⋅n>
229  a->AddTrialIntegrator(new NormalTraceIntegrator,nullptr,
230  TrialSpace::hatp_space,TestSpace::v_space);
231  // < û,q >
232  a->AddTrialIntegrator(new TraceIntegrator,nullptr,
233  TrialSpace::hatu_space,TestSpace::q_space);
234 
235  // test space integrators (Adjoint graph norm)
236  // (∇q,∇δq)
237  a->AddTestIntegrator(new DiffusionIntegrator(one),nullptr,
238  TestSpace::q_space, TestSpace::q_space);
239  // (q,δq)
240  a->AddTestIntegrator(new MassIntegrator(one),nullptr,
241  TestSpace::q_space, TestSpace::q_space);
242  // (∇⋅v,∇⋅δv)
243  a->AddTestIntegrator(new DivDivIntegrator(one),nullptr,
244  TestSpace::v_space, TestSpace::v_space);
245  // (v,δv)
246  a->AddTestIntegrator(new VectorFEMassIntegrator(one),nullptr,
247  TestSpace::v_space, TestSpace::v_space);
248  // -i ω (∇q,δv)
249  a->AddTestIntegrator(nullptr,new MixedVectorGradientIntegrator(negomeg),
250  TestSpace::q_space, TestSpace::v_space);
251  // i ω (v,∇ δq)
252  a->AddTestIntegrator(nullptr,new MixedVectorWeakDivergenceIntegrator(negomeg),
253  TestSpace::v_space, TestSpace::q_space);
254  // ω^2 (v,δv)
255  a->AddTestIntegrator(new VectorFEMassIntegrator(omeg2),nullptr,
256  TestSpace::v_space, TestSpace::v_space);
257  // - i ω (∇⋅v,δq)
258  a->AddTestIntegrator(nullptr,new VectorFEDivergenceIntegrator(negomeg),
259  TestSpace::v_space, TestSpace::q_space);
260  // i ω (q,∇⋅v)
261  a->AddTestIntegrator(nullptr,new MixedScalarWeakGradientIntegrator(negomeg),
262  TestSpace::q_space, TestSpace::v_space);
263  // ω^2 (q,δq)
264  a->AddTestIntegrator(new MassIntegrator(omeg2),nullptr,
265  TestSpace::q_space, TestSpace::q_space);
266 
267  // RHS
271  {
272  a->AddDomainLFIntegrator(new DomainLFIntegrator(f_rhs_r),
273  new DomainLFIntegrator(f_rhs_i), TestSpace::q_space);
274  }
275 
280 
281  socketstream p_out_r;
282  socketstream p_out_i;
283 
284  double err0 = 0.;
285  int dof0 = 0; // init to suppress gcc warning
286 
287  std::cout << "\n Ref |"
288  << " Dofs |"
289  << " ω |"
290  << " L2 Error |"
291  << " Rate |"
292  << " PCG it |" << endl;
293  std::cout << std::string(60,'-')
294  << endl;
295 
296  for (int it = 0; it<=ref; it++)
297  {
298  if (static_cond) { a->EnableStaticCondensation(); }
299  a->Assemble();
300 
301  Array<int> ess_tdof_list;
302  Array<int> ess_bdr;
303  if (mesh.bdr_attributes.Size())
304  {
305  ess_bdr.SetSize(mesh.bdr_attributes.Max());
306  ess_bdr = 1;
307  hatp_fes->GetEssentialTrueDofs(ess_bdr, ess_tdof_list);
308  }
309 
310  // shift the ess_tdofs
311  for (int j = 0; j < ess_tdof_list.Size(); j++)
312  {
313  ess_tdof_list[j] += p_fes->GetTrueVSize() + u_fes->GetTrueVSize();
314  }
315 
316  Array<int> offsets(5);
317  offsets[0] = 0;
318  offsets[1] = p_fes->GetVSize();
319  offsets[2] = u_fes->GetVSize();
320  offsets[3] = hatp_fes->GetVSize();
321  offsets[4] = hatu_fes->GetVSize();
322  offsets.PartialSum();
323 
324  Vector x(2*offsets.Last());
325  x = 0.;
326 
327  GridFunction hatp_gf_r(hatp_fes, x, offsets[2]);
328  GridFunction hatp_gf_i(hatp_fes, x, offsets.Last()+ offsets[2]);
329  hatp_gf_r.ProjectBdrCoefficient(hatpex_r, ess_bdr);
330  hatp_gf_i.ProjectBdrCoefficient(hatpex_i, ess_bdr);
331 
332  OperatorPtr Ah;
333  Vector X,B;
334  a->FormLinearSystem(ess_tdof_list,x,Ah, X,B);
335 
336  // Setup operator and preconditioner
337  ComplexOperator * Ahc = Ah.As<ComplexOperator>();
338 
339  BlockMatrix * A_r = dynamic_cast<BlockMatrix *>(&Ahc->real());
340  BlockMatrix * A_i = dynamic_cast<BlockMatrix *>(&Ahc->imag());
341 
342  int num_blocks = A_r->NumRowBlocks();
343  Array<int> tdof_offsets(2*num_blocks+1);
344  tdof_offsets[0] = 0;
345  int k = (static_cond) ? 2 : 0;
346  for (int i=0; i<num_blocks; i++)
347  {
348  tdof_offsets[i+1] = trial_fes[i+k]->GetTrueVSize();
349  tdof_offsets[num_blocks+i+1] = trial_fes[i+k]->GetTrueVSize();
350  }
351  tdof_offsets.PartialSum();
352 
353  BlockOperator A(tdof_offsets);
354  for (int i = 0; i<num_blocks; i++)
355  {
356  for (int j = 0; j<num_blocks; j++)
357  {
358  A.SetBlock(i,j,&A_r->GetBlock(i,j));
359  A.SetBlock(i,j+num_blocks,&A_i->GetBlock(i,j), -1.0);
360  A.SetBlock(i+num_blocks,j+num_blocks,&A_r->GetBlock(i,j));
361  A.SetBlock(i+num_blocks,j,&A_i->GetBlock(i,j));
362  }
363  }
364 
365  BlockDiagonalPreconditioner M(tdof_offsets);
366  M.owns_blocks = 1;
367  for (int i = 0; i<num_blocks; i++)
368  {
369  M.SetDiagonalBlock(i, new GSSmoother((SparseMatrix&)A_r->GetBlock(i,i)));
370  M.SetDiagonalBlock(num_blocks+i, new GSSmoother((SparseMatrix&)A_r->GetBlock(i,
371  i)));
372  }
373 
374  CGSolver cg;
375  cg.SetRelTol(1e-10);
376  cg.SetMaxIter(2000);
377  cg.SetPrintLevel(0);
378  cg.SetPreconditioner(M);
379  cg.SetOperator(A);
380  cg.Mult(B, X);
381 
382  a->RecoverFEMSolution(X,x);
383 
384  GridFunction p_r(p_fes, x, 0);
385  GridFunction p_i(p_fes, x, offsets.Last());
386 
387  GridFunction u_r(u_fes, x, offsets[1]);
388  GridFunction u_i(u_fes, x, offsets.Last() + offsets[1]);
389 
392 
395 
396  int dofs = 0;
397  for (int i = 0; i<trial_fes.Size(); i++)
398  {
399  dofs += trial_fes[i]->GetTrueVSize();
400  }
401 
402  double p_err_r = p_r.ComputeL2Error(p_ex_r);
403  double p_err_i = p_i.ComputeL2Error(p_ex_i);
404  double u_err_r = u_r.ComputeL2Error(u_ex_r);
405  double u_err_i = u_i.ComputeL2Error(u_ex_i);
406 
407  double L2Error = sqrt(p_err_r*p_err_r + p_err_i*p_err_i
408  +u_err_r*u_err_r + u_err_i*u_err_i);
409 
410  double rate_err = (it) ? dim*log(err0/L2Error)/log((double)dof0/dofs) : 0.0;
411 
412  err0 = L2Error;
413  dof0 = dofs;
414 
415  std::ios oldState(nullptr);
416  oldState.copyfmt(std::cout);
417  std::cout << std::right << std::setw(5) << it << " | "
418  << std::setw(10) << dof0 << " | "
419  << std::setprecision(1) << std::fixed
420  << std::setw(4) << 2*rnum << " π | "
421  << std::setprecision(3)
422  << std::setw(10) << std::scientific << err0 << " | "
423  << std::setprecision(2)
424  << std::setw(6) << std::fixed << rate_err << " | "
425  << std::setw(6) << std::fixed << cg.GetNumIterations() << " | "
426  << std::endl;
427  std::cout.copyfmt(oldState);
428 
429  if (visualization)
430  {
431  const char * keys = (it == 0 && dim == 2) ? "jRcml\n" : nullptr;
432  char vishost[] = "localhost";
433  int visport = 19916;
434  VisualizeField(p_out_r,vishost, visport, p_r,
435  "Numerical presure (real part)", 0, 0, 500, 500, keys);
436  VisualizeField(p_out_i,vishost, visport, p_i,
437  "Numerical presure (imaginary part)", 501, 0, 500, 500, keys);
438  }
439 
440  if (it == ref)
441  {
442  break;
443  }
444 
445  mesh.UniformRefinement();
446  for (int i =0; i<trial_fes.Size(); i++)
447  {
448  trial_fes[i]->Update(false);
449  }
450  a->Update();
451  }
452 
453  delete a;
454  delete q_fec;
455  delete v_fec;
456  delete hatp_fes;
457  delete hatp_fec;
458  delete hatu_fes;
459  delete hatu_fec;
460  delete u_fec;
461  delete p_fec;
462  delete u_fes;
463  delete p_fes;
464 
465  return 0;
466 }
467 
468 double p_exact_r(const Vector &x)
469 {
470  return acoustics_solution(x).real();
471 }
472 
473 double p_exact_i(const Vector &x)
474 {
475  return acoustics_solution(x).imag();
476 }
477 
478 double hatp_exact_r(const Vector & X)
479 {
480  return p_exact_r(X);
481 }
482 
483 double hatp_exact_i(const Vector & X)
484 {
485  return p_exact_i(X);
486 }
487 
488 void gradp_exact_r(const Vector &x, Vector &grad_r)
489 {
490  grad_r.SetSize(x.Size());
491  vector<complex<double>> grad;
492  acoustics_solution_grad(x,grad);
493  for (unsigned i = 0; i < grad.size(); i++)
494  {
495  grad_r[i] = grad[i].real();
496  }
497 }
498 
499 void gradp_exact_i(const Vector &x, Vector &grad_i)
500 {
501  grad_i.SetSize(x.Size());
502  vector<complex<double>> grad;
503  acoustics_solution_grad(x,grad);
504  for (unsigned i = 0; i < grad.size(); i++)
505  {
506  grad_i[i] = grad[i].imag();
507  }
508 }
509 
510 double d2_exact_r(const Vector &x)
511 {
512 
513  return acoustics_solution_laplacian(x).real();
514 }
515 
516 double d2_exact_i(const Vector &x)
517 {
518  return acoustics_solution_laplacian(x).imag();
519 }
520 
521 // u = - ∇ p / (i ω )
522 // = i (∇ p_r + i * ∇ p_i) / ω
523 // = - ∇ p_i / ω + i ∇ p_r / ω
524 void u_exact_r(const Vector &x, Vector & u)
525 {
526  gradp_exact_i(x,u);
527  u *= -1./omega;
528 }
529 
530 void u_exact_i(const Vector &x, Vector & u)
531 {
532  gradp_exact_r(x,u);
533  u *= 1./omega;
534 }
535 
536 void hatu_exact_r(const Vector & X, Vector & hatu)
537 {
538  u_exact_r(X,hatu);
539 }
540 void hatu_exact_i(const Vector & X, Vector & hatu)
541 {
542  u_exact_i(X,hatu);
543 }
544 
545 // ∇⋅u = i Δ p / ω
546 // = i (Δ p_r + i * Δ p_i) / ω
547 // = - Δ p_i / ω + i Δ p_r / ω
548 double divu_exact_r(const Vector &x)
549 {
550  return -d2_exact_i(x)/omega;
551 }
552 
553 double divu_exact_i(const Vector &x)
554 {
555  return d2_exact_r(x)/omega;
556 }
557 
558 // f = ∇⋅u + i ω p
559 // f_r = ∇⋅u_r - ω p_i
560 double rhs_func_r(const Vector &x)
561 {
562  double p = p_exact_i(x);
563  double divu = divu_exact_r(x);
564  return divu - omega * p;
565 }
566 
567 // f_i = ∇⋅u_i + ω p_r
568 double rhs_func_i(const Vector &x)
569 {
570  double p = p_exact_r(x);
571  double divu = divu_exact_i(x);
572  return divu + omega * p;
573 }
574 
575 complex<double> acoustics_solution(const Vector & X)
576 {
577  complex<double> zi = complex<double>(0., 1.);
578  switch (prob)
579  {
580  case plane_wave:
581  {
582  double beta = omega/std::sqrt((double)X.Size());
583  complex<double> alpha = beta * zi * X.Sum();
584  return exp(alpha);
585  }
586  break;
587  default:
588  {
589  double rk = omega;
590  double degrees = 45;
591  double alpha = (180+degrees) * M_PI/180.;
592  double sina = sin(alpha);
593  double cosa = cos(alpha);
594  // shift the origin
595  double xprim=X(0) + 0.1;
596  double yprim=X(1) + 0.1;
597 
598  double x = xprim*sina - yprim*cosa;
599  double y = xprim*cosa + yprim*sina;
600  //wavelength
601  double rl = 2.*M_PI/rk;
602  // beam waist radius
603  double w0 = 0.05;
604  // function w
605  double fact = rl/M_PI/(w0*w0);
606  double aux = 1. + (fact*y)*(fact*y);
607  double w = w0*sqrt(aux);
608  double phi0 = atan(fact*y);
609  double r = y + 1./y/(fact*fact);
610 
611  // pressure
612  complex<double> ze = - x*x/(w*w) - zi*rk*y - zi * M_PI * x * x/rl/r +
613  zi*phi0/2.;
614 
615  double pf = pow(2.0/M_PI/(w*w),0.25);
616  return pf*exp(ze);
617  }
618  break;
619  }
620 }
621 
622 void acoustics_solution_grad(const Vector & X, vector<complex<double>> & dp)
623 {
624  dp.resize(X.Size());
625  complex<double> zi = complex<double>(0., 1.);
626  switch (prob)
627  {
628  case plane_wave:
629  {
630  double beta = omega/std::sqrt((double)X.Size());
631  complex<double> alpha = beta * zi * X.Sum();
632  complex<double> p = exp(alpha);
633  for (int i = 0; i<X.Size(); i++)
634  {
635  dp[i] = zi * beta * p;
636  }
637  }
638  break;
639  default:
640  {
641  double rk = omega;
642  double degrees = 45;
643  double alpha = (180+degrees) * M_PI/180.;
644  double sina = sin(alpha);
645  double cosa = cos(alpha);
646  // shift the origin
647  double xprim=X(0) + 0.1;
648  double yprim=X(1) + 0.1;
649 
650  double x = xprim*sina - yprim*cosa;
651  double y = xprim*cosa + yprim*sina;
652  double dxdxprim = sina, dxdyprim = -cosa;
653  double dydxprim = cosa, dydyprim = sina;
654  //wavelength
655  double rl = 2.*M_PI/rk;
656 
657  // beam waist radius
658  double w0 = 0.05;
659 
660  // function w
661  double fact = rl/M_PI/(w0*w0);
662  double aux = 1. + (fact*y)*(fact*y);
663 
664  double w = w0*sqrt(aux);
665  double dwdy = w0*fact*fact*y/sqrt(aux);
666 
667  double phi0 = atan(fact*y);
668  double dphi0dy = cos(phi0)*cos(phi0)*fact;
669 
670  double r = y + 1./y/(fact*fact);
671  double drdy = 1. - 1./(y*y)/(fact*fact);
672 
673  // pressure
674  complex<double> ze = - x*x/(w*w) - zi*rk*y - zi * M_PI * x * x/rl/r +
675  zi*phi0/2.;
676 
677  complex<double> zdedx = -2.*x/(w*w) - 2.*zi*M_PI*x/rl/r;
678  complex<double> zdedy = 2.*x*x/(w*w*w)*dwdy - zi*rk + zi*M_PI*x*x/rl/
679  (r*r)*drdy + zi*dphi0dy/2.;
680 
681  double pf = pow(2.0/M_PI/(w*w),0.25);
682  double dpfdy = -pow(2./M_PI/(w*w),-0.75)/M_PI/(w*w*w)*dwdy;
683 
684  complex<double> zp = pf*exp(ze);
685  complex<double> zdpdx = zp*zdedx;
686  complex<double> zdpdy = dpfdy*exp(ze)+zp*zdedy;
687 
688  dp[0] = (zdpdx*dxdxprim + zdpdy*dydxprim);
689  dp[1] = (zdpdx*dxdyprim + zdpdy*dydyprim);
690  if (dim == 3) { dp[2] = 0.0; }
691  }
692  break;
693  }
694 }
695 
696 complex<double> acoustics_solution_laplacian(const Vector & X)
697 {
698  complex<double> zi = complex<double>(0., 1.);
699  switch (prob)
700  {
701  case plane_wave:
702  {
703  double beta = omega/std::sqrt((double)X.Size());
704  complex<double> alpha = beta * zi * X.Sum();
705  complex<double> p = exp(alpha);
706  return dim * beta * beta * p;
707  }
708  break;
709  default:
710  {
711  double rk = omega;
712  double degrees = 45;
713  double alpha = (180+degrees) * M_PI/180.;
714  double sina = sin(alpha);
715  double cosa = cos(alpha);
716  // shift the origin
717  double xprim=X(0) + 0.1;
718  double yprim=X(1) + 0.1;
719 
720  double x = xprim*sina - yprim*cosa;
721  double y = xprim*cosa + yprim*sina;
722  double dxdxprim = sina, dxdyprim = -cosa;
723  double dydxprim = cosa, dydyprim = sina;
724  //wavelength
725  double rl = 2.*M_PI/rk;
726 
727  // beam waist radius
728  double w0 = 0.05;
729 
730  // function w
731  double fact = rl/M_PI/(w0*w0);
732  double aux = 1. + (fact*y)*(fact*y);
733 
734  double w = w0*sqrt(aux);
735  double dwdy = w0*fact*fact*y/sqrt(aux);
736  double d2wdydy = w0*fact*fact*(1. - (fact*y)*(fact*y)/aux)/sqrt(aux);
737 
738  double phi0 = atan(fact*y);
739  double dphi0dy = cos(phi0)*cos(phi0)*fact;
740  double d2phi0dydy = -2.*cos(phi0)*sin(phi0)*fact*dphi0dy;
741 
742  double r = y + 1./y/(fact*fact);
743  double drdy = 1. - 1./(y*y)/(fact*fact);
744  double d2rdydy = 2./(y*y*y)/(fact*fact);
745 
746  // pressure
747  complex<double> ze = - x*x/(w*w) - zi*rk*y - zi * M_PI * x * x/rl/r +
748  zi*phi0/2.;
749 
750  complex<double> zdedx = -2.*x/(w*w) - 2.*zi*M_PI*x/rl/r;
751  complex<double> zdedy = 2.*x*x/(w*w*w)*dwdy - zi*rk + zi*M_PI*x*x/rl/
752  (r*r)*drdy + zi*dphi0dy/2.;
753  complex<double> zd2edxdx = -2./(w*w) - 2.*zi*M_PI/rl/r;
754  complex<double> zd2edxdy = 4.*x/(w*w*w)*dwdy + 2.*zi*M_PI*x/rl/(r*r)*drdy;
755  complex<double> zd2edydx = zd2edxdy;
756  complex<double> zd2edydy = -6.*x*x/(w*w*w*w)*dwdy*dwdy + 2.*x*x/
757  (w*w*w)*d2wdydy - 2.*zi*M_PI*x*x/rl/(r*r*r)*drdy*drdy
758  + zi*M_PI*x*x/rl/(r*r)*d2rdydy + zi/2.*d2phi0dydy;
759 
760  double pf = pow(2.0/M_PI/(w*w),0.25);
761  double dpfdy = -pow(2./M_PI/(w*w),-0.75)/M_PI/(w*w*w)*dwdy;
762  double d2pfdydy = -1./M_PI*pow(2./M_PI,-0.75)*(-1.5*pow(w,-2.5)
763  *dwdy*dwdy + pow(w,-1.5)*d2wdydy);
764 
765  complex<double> zp = pf*exp(ze);
766  complex<double> zdpdx = zp*zdedx;
767  complex<double> zdpdy = dpfdy*exp(ze)+zp*zdedy;
768  complex<double> zd2pdxdx = zdpdx*zdedx + zp*zd2edxdx;
769  complex<double> zd2pdxdy = zdpdy*zdedx + zp*zd2edxdy;
770  complex<double> zd2pdydx = dpfdy*exp(ze)*zdedx + zdpdx*zdedy + zp*zd2edydx;
771  complex<double> zd2pdydy = d2pfdydy*exp(ze) + dpfdy*exp(
772  ze)*zdedy + zdpdy*zdedy + zp*zd2edydy;
773 
774 
775  return (zd2pdxdx*dxdxprim + zd2pdydx*dydxprim)*dxdxprim
776  + (zd2pdxdy*dxdxprim + zd2pdydy*dydxprim)*dydxprim
777  + (zd2pdxdx*dxdyprim + zd2pdydx*dydyprim)*dxdyprim
778  + (zd2pdxdy*dxdyprim + zd2pdydy*dydyprim)*dydyprim;
779  }
780  break;
781  }
782 }
virtual double ComputeL2Error(Coefficient *exsol[], const IntegrationRule *irs[]=NULL, const Array< int > *elems=NULL) const
Definition: gridfunc.cpp:2786
double p_exact_i(const Vector &x)
Definition: acoustics.cpp:473
Class for domain integration L(v) := (f, v)
Definition: lininteg.hpp:108
int visport
Conjugate gradient method.
Definition: solvers.hpp:493
virtual Operator & real()
Real or imaginary part accessor methods.
Class for grid function - Vector with associated FE space.
Definition: gridfunc.hpp:30
A coefficient that is constant across space and time.
Definition: coefficient.hpp:84
void PrintOptions(std::ostream &out) const
Print the options.
Definition: optparser.cpp:331
int Dimension() const
Dimension of the reference space used within the elements.
Definition: mesh.hpp:1020
void u_exact_r(const Vector &x, Vector &u)
Definition: acoustics.cpp:524
void SetSize(int s)
Resize the vector to size s.
Definition: vector.hpp:517
Mimic the action of a complex operator using two real operators.
void hatu_exact_r(const Vector &X, Vector &hatu)
Definition: acoustics.cpp:536
void PrintUsage(std::ostream &out) const
Print the usage message.
Definition: optparser.cpp:462
Vector coefficient that is constant in space and time.
Pointer to an Operator of a specified type.
Definition: handle.hpp:33
virtual void Mult(const Vector &b, Vector &x) const
Operator application: y=A(x).
Definition: solvers.cpp:718
int Size() const
Returns the size of the vector.
Definition: vector.hpp:197
T Max() const
Find the maximal element in the array, using the comparison operator < for class T.
Definition: array.cpp:68
virtual void GetEssentialTrueDofs(const Array< int > &bdr_attr_is_ess, Array< int > &ess_tdof_list, int component=-1)
Get a list of essential true dofs, ess_tdof_list, corresponding to the boundary attributes marked in ...
Definition: fespace.cpp:587
double d2_exact_i(const Vector &x)
Definition: acoustics.cpp:516
bool Good() const
Return true if the command line options were parsed successfully.
Definition: optparser.hpp:159
Vector beta
int main(int argc, char *argv[])
Definition: acoustics.cpp:111
STL namespace.
double hatp_exact_r(const Vector &X)
Definition: acoustics.cpp:478
double hatp_exact_i(const Vector &X)
Definition: acoustics.cpp:483
(Q div u, div v) for RT elements
void gradp_exact_r(const Vector &x, Vector &gradu)
Definition: acoustics.cpp:488
prob_type prob
Definition: acoustics.cpp:109
Data type for Gauss-Seidel smoother of sparse matrix.
double divu_exact_i(const Vector &x)
Definition: acoustics.cpp:553
int NumRowBlocks() const
Return the number of row blocks.
Definition: blockmatrix.hpp:46
virtual void SetPrintLevel(int print_lvl)
Legacy method to set the level of verbosity of the solver output.
Definition: solvers.cpp:71
Class representing the DPG weak formulation for complex valued systems (see the class DPGWeakForm)...
void Parse()
Parse the command-line options. Note that this function expects all the options provided through the ...
Definition: optparser.cpp:151
A class to handle Block diagonal preconditioners in a matrix-free implementation. ...
char vishost[]
Data type sparse matrix.
Definition: sparsemat.hpp:50
int Append(const T &el)
Append element &#39;el&#39; to array, resize if necessary.
Definition: array.hpp:759
int GetNumIterations() const
Returns the number of iterations taken during the last call to Mult()
Definition: solvers.hpp:250
void UniformRefinement(int i, const DSTable &, int *, int *, int *)
Definition: mesh.cpp:10232
Arbitrary order "H^{1/2}-conforming" trace finite elements defined on the interface between mesh elem...
Definition: fe_coll.hpp:319
void SetMaxIter(int max_it)
Definition: solvers.hpp:201
double Sum() const
Return the sum of the vector entries.
Definition: vector.cpp:1286
complex< double > acoustics_solution(const Vector &X)
Definition: acoustics.cpp:575
double rhs_func_i(const Vector &x)
Definition: acoustics.cpp:568
Arbitrary order H(div)-conforming Raviart-Thomas finite elements.
Definition: fe_coll.hpp:380
A general vector function coefficient.
void acoustics_solution_grad(const Vector &X, vector< complex< double >> &dp)
Definition: acoustics.cpp:622
virtual int GetTrueVSize() const
Return the number of vector true (conforming) dofs.
Definition: fespace.hpp:712
Array< int > bdr_attributes
A list of all unique boundary attributes used by the Mesh.
Definition: mesh.hpp:275
void SetRelTol(double rtol)
Definition: solvers.hpp:199
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Definition: fespace.hpp:219
Collection of finite elements from the same family in multiple dimensions. This class is used to matc...
Definition: fe_coll.hpp:26
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
void SetSize(int nsize)
Change the logical size of the array, keep existing entries.
Definition: array.hpp:687
void PartialSum()
Fill the entries of the array with the cumulative sum of the entries.
Definition: array.cpp:103
virtual Operator & imag()
double omega
Definition: acoustics.cpp:101
double p_exact_r(const Vector &x)
Definition: acoustics.cpp:468
prob_type
Definition: acoustics.cpp:103
double a
Definition: lissajous.cpp:41
OpType * As() const
Return the Operator pointer statically cast to a specified OpType. Similar to the method Get()...
Definition: handle.hpp:104
void hatu_exact_i(const Vector &X, Vector &hatu)
Definition: acoustics.cpp:540
int GetVSize() const
Return the number of vector dofs, i.e. GetNDofs() x GetVDim().
Definition: fespace.hpp:709
double divu_exact_r(const Vector &x)
Definition: acoustics.cpp:548
SparseMatrix & GetBlock(int i, int j)
Return a reference to block (i,j). Reference may be invalid if Aij(i,j) == NULL.
Definition: blockmatrix.cpp:82
complex< double > acoustics_solution_laplacian(const Vector &X)
Definition: acoustics.cpp:696
void u_exact_i(const Vector &x, Vector &u)
Definition: acoustics.cpp:530
int Size() const
Return the logical size of the array.
Definition: array.hpp:141
T & Last()
Return the last element in the array.
Definition: array.hpp:792
double rhs_func_r(const Vector &x)
Definition: acoustics.cpp:560
virtual void SetOperator(const Operator &op)
Also calls SetOperator for the preconditioner if there is one.
Definition: solvers.hpp:507
const double alpha
Definition: ex15.cpp:369
A general function coefficient.
Vector data type.
Definition: vector.hpp:58
double d2_exact_r(const Vector &x)
Definition: acoustics.cpp:510
Arbitrary order "H^{-1/2}-conforming" face finite elements defined on the interface between mesh elem...
Definition: fe_coll.hpp:434
virtual void SetPreconditioner(Solver &pr)
This should be called before SetOperator.
Definition: solvers.cpp:173
Arbitrary order H1-conforming (continuous) finite elements.
Definition: fe_coll.hpp:259
double u(const Vector &xvec)
Definition: lor_mms.hpp:22
prob_type
Definition: ex25.cpp:146
A class to handle Block systems in a matrix-free implementation.
void gradp_exact_i(const Vector &x, Vector &gradu)
Definition: acoustics.cpp:499
int dim
Definition: acoustics.cpp:100
void SetBlock(int iRow, int iCol, Operator *op, double c=1.0)
Add a block op in the block-entry (iblock, jblock).
void ProjectBdrCoefficient(Coefficient &coeff, Array< int > &attr)
Project a Coefficient on the GridFunction, modifying only DOFs on the boundary associated with the bo...
Definition: gridfunc.hpp:468
void SetDiagonalBlock(int iblock, Operator *op)
Add a square block op in the block-entry (iblock, iblock).
Arbitrary order "L2-conforming" discontinuous finite elements.
Definition: fe_coll.hpp:327