unstructured_scattering.cc
Go to the documentation of this file.
1 //LIC// ====================================================================
2 //LIC// This file forms part of oomph-lib, the object-oriented,
3 //LIC// multi-physics finite-element library, available
4 //LIC// at http://www.oomph-lib.org.
5 //LIC//
6 //LIC// Copyright (C) 2006-2024 Matthias Heil and Andrew Hazel
7 //LIC//
8 //LIC// This library is free software; you can redistribute it and/or
9 //LIC// modify it under the terms of the GNU Lesser General Public
10 //LIC// License as published by the Free Software Foundation; either
11 //LIC// version 2.1 of the License, or (at your option) any later version.
12 //LIC//
13 //LIC// This library is distributed in the hope that it will be useful,
14 //LIC// but WITHOUT ANY WARRANTY; without even the implied warranty of
15 //LIC// MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
16 //LIC// Lesser General Public License for more details.
17 //LIC//
18 //LIC// You should have received a copy of the GNU Lesser General Public
19 //LIC// License along with this library; if not, write to the Free Software
20 //LIC// Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA
21 //LIC// 02110-1301 USA.
22 //LIC//
23 //LIC// The authors may be contacted at oomph-lib@maths.man.ac.uk.
24 //LIC//
25 //LIC//====================================================================
26 //Driver for a specific 2D Helmholtz problem with flux boundary conditions
27 //uses two separate meshes for bulk and surface mesh
28 
29 #include "math.h"
30 #include <complex>
31 
32 
33 //Generic routines
34 #include "generic.h"
35 
36 // The Helmholtz equations
37 #include "helmholtz.h"
38 
39 // The mesh
40 #include "meshes/triangle_mesh.h"
41 
42 // Get the Bessel functions
43 #include "oomph_crbond_bessel.h"
44 
45 using namespace oomph;
46 using namespace std;
47 
48 
49 
50 /// //////////////////////////////////////////////////////////////////
51 /// //////////////////////////////////////////////////////////////////
52 /// //////////////////////////////////////////////////////////////////
53 
54 //===== start_of_namespace=============================================
55 /// Namespace for the Helmholtz problem parameters
56 //=====================================================================
57 namespace GlobalParameters
58 {
59 
60  /// Square of the wavenumber
61  double K_squared=10.0;
62 
63  /// Number of terms used in the computation
64  /// of the exact solution
65  unsigned N_fourier=10;
66 
67  /// Flag to choose the Dirichlet to Neumann BC
68  /// or ABC BC
69  bool DtN_BC=false;
70 
71  /// Flag to choose wich order to use
72  // in the ABCs BC: 1 for ABC 1st order...
73  unsigned ABC_order=3;
74 
75  /// Radius of outer boundary (must be a circle!)
76  double Outer_radius=1.5;
77 
78  /// Imaginary unit
79  std::complex<double> I(0.0,1.0);
80 
81  /// Exact solution for scattered field
82  /// (vector returns real and impaginary parts).
83  void get_exact_u(const Vector<double>& x, Vector<double>& u)
84  {
85  // Switch to polar coordinates
86  double r;
87  r=sqrt(x[0]*x[0]+x[1]*x[1]);
88  double theta;
89  theta=atan2(x[1],x[0]);
90 
91  // Argument for Bessel/Hankel functions
92  double rr=sqrt(K_squared)*r;
93 
94  // Evaluate Bessel/Hankel functions
95  complex <double > u_ex(0.0,0.0);
96  Vector<double> jn(N_fourier+1), yn(N_fourier+1),
97  jnp(N_fourier+1), ynp(N_fourier+1);
98  Vector<double> jn_a(N_fourier+1),yn_a(N_fourier+1),
99  jnp_a(N_fourier+1), ynp_a(N_fourier+1);
100  Vector<complex<double> > h(N_fourier+1),h_a(N_fourier+1),
101  hp(N_fourier+1), hp_a(N_fourier+1);
102 
103  // We want to compute N_fourier terms but the function
104  // may return fewer than that.
105  int n_actual=0;
106  CRBond_Bessel::bessjyna(N_fourier,sqrt(K_squared),n_actual,
107  &jn_a[0],&yn_a[0],
108  &jnp_a[0],&ynp_a[0]);
109 
110  // Shout if things went wrong
111 #ifdef PARANOID
112  if (n_actual!=int(N_fourier))
113  {
114  std::ostringstream error_stream;
115  error_stream << "CRBond_Bessel::bessjyna() only computed "
116  << n_actual << " rather than " << N_fourier
117  << " Bessel functions.\n";
118  throw OomphLibError(error_stream.str(),
119  OOMPH_CURRENT_FUNCTION,
120  OOMPH_EXCEPTION_LOCATION);
121  }
122 #endif
123 
124  // Evaluate Hankel at actual radius
125  Hankel_functions_for_helmholtz_problem::Hankel_first(N_fourier,rr,h,hp);
126 
127  // Evaluate Hankel at inner (unit) radius
128  Hankel_functions_for_helmholtz_problem::Hankel_first(N_fourier
129  ,sqrt(K_squared),
130  h_a,hp_a);
131 
132  // Compute the sum: Separate the computation of the negative
133  // and positive terms
134  // ALH: The construction with the static cast to a double is to avoid
135  // a floating point exception when running unoptimised on my machine
136  for (unsigned i=0;i<N_fourier;i++)
137  {
138  u_ex-=pow(I,i)*h[i]*jnp_a[i]*pow(exp(I*theta),static_cast<double>(i))/hp_a[i];
139  }
140  for (unsigned i=1;i<N_fourier;i++)
141  {
142  u_ex-=pow(I,i)*h[i]*jnp_a[i]*pow(exp(-I*theta),static_cast<double>(i))/hp_a[i];
143  }
144 
145  // Get the real & imaginary part of the result
146  u[0]=real(u_ex);
147  u[1]=imag(u_ex);
148 
149  }// end of get_exact_u
150 
151 
152 
153  /// Flux (normal derivative) on the unit disk
154  /// for a planar incoming wave
155  void prescribed_incoming_flux(const Vector<double>& x,
156  complex<double>& flux)
157  {
158  // Switch to polar coordinates
159  double r;
160  r=sqrt(x[0]*x[0]+x[1]*x[1]);
161  double theta;
162  theta=atan2(x[1],x[0]);
163 
164  // Argument of the Bessel/Hankel fcts
165  double rr=sqrt(K_squared)*r;
166 
167  // Compute Bessel/Hankel functions
168  Vector<double> jn(N_fourier+1), yn(N_fourier+1),
169  jnp(N_fourier+1), ynp(N_fourier+1);
170 
171  // We want to compute N_fourier terms but the function
172  // may return fewer than that.
173  int n_actual=0;
174  CRBond_Bessel::bessjyna(N_fourier,rr,n_actual,&jn[0],&yn[0],
175  &jnp[0],&ynp[0]);
176 
177  // Shout if things went wrong...
178 #ifdef PARANOID
179  if (n_actual!=int(N_fourier))
180  {
181  std::ostringstream error_stream;
182  error_stream << "CRBond_Bessel::bessjyna() only computed "
183  << n_actual << " rather than " << N_fourier
184  << " Bessel functions.\n";
185  throw OomphLibError(error_stream.str(),
186  OOMPH_CURRENT_FUNCTION,
187  OOMPH_EXCEPTION_LOCATION);
188  }
189 #endif
190 
191  // Compute the sum: Separate the computation of the negative and
192  // positive terms
193  flux=std::complex<double>(0.0,0.0);
194  for (unsigned i=0;i<N_fourier;i++)
195  {
196  flux+=pow(I,i)*(sqrt(K_squared))*pow(exp(I*theta),i)*jnp[i];
197  }
198  for (unsigned i=1;i<N_fourier;i++)
199  {
200  flux+=pow(I,i)*(sqrt(K_squared))*pow(exp(-I*theta),i)*jnp[i];
201  }
202 
203 
204  }// end of prescribed_incoming_flux
205 
206 } // end of namespace
207 
208 
209 
210 
211 /// //////////////////////////////////////////////////////////////////
212 /// //////////////////////////////////////////////////////////////////
213 /// //////////////////////////////////////////////////////////////////
214 
215 
216 
217 //========= start_of_problem_class=====================================
218 /// Problem class to compute scattering of planar wave from unit disk
219 //=====================================================================
220 template<class ELEMENT>
221 class ScatteringProblem : public Problem
222 {
223 
224 public:
225 
226  /// Constructor
228 
229  /// Destructor (empty)
231 
232  /// Doc the solution. DocInfo object stores flags/labels for where the
233  /// output gets written to
234  void doc_solution(DocInfo& doc_info);
235 
236  /// Update the problem specs before solve (empty)
238 
239  /// Update the problem specs after solve (empty)
241 
242  /// Recompute gamma integral before checking Newton residuals
244  {
246  {
247  Helmholtz_outer_boundary_mesh_pt->setup_gamma();
248  }
249  }
250 
251  /// Actions before adapt: Wipe the mesh of prescribed flux elements
253 
254  /// Actions after adapt: Rebuild the mesh of prescribed flux elements
256 
257  /// Create BC elements on boundary b of the Mesh pointed
258  /// to by bulk_mesh_pt and add them to the specified survace Mesh
260  const unsigned &b, Mesh* const &bulk_mesh_pt,
261  Mesh* const & helmholtz_outer_boundary_mesh_pt);
262 
263  /// Create Helmholtz flux elements on boundary b of the Mesh pointed
264  /// to by bulk_mesh_pt and add them to the specified surface Mesh
265  void create_flux_elements(const unsigned &b, Mesh* const &bulk_mesh_pt,
266  Mesh* const & helmholtz_inner_boundary_mesh_pt);
267 
268  /// Delete boundary face elements and wipe the surface mesh
269  void delete_face_elements( Mesh* const & boundary_mesh_pt);
270 
271  /// Set pointer to prescribed-flux function for all
272  /// elements in the surface mesh on the surface of the unit disk
274 
275  /// Set up boundary condition elements on outer boundary
277 
278 #ifdef ADAPTIVE
279 
280  /// Pointer to the "bulk" mesh
281  RefineableTriangleMesh<ELEMENT>* Bulk_mesh_pt;
282 
283 #else
284 
285  /// Pointer to the "bulk" mesh
286  TriangleMesh<ELEMENT>* Bulk_mesh_pt;
287 
288 #endif
289 
290  /// Pointer to mesh containing the DtN (or ABC) boundary
291  /// condition elements
292  HelmholtzDtNMesh<ELEMENT>* Helmholtz_outer_boundary_mesh_pt;
293 
294  /// Pointer to the mesh containing
295  /// the Helmholtz inner boundary condition elements
296  Mesh* Helmholtz_inner_boundary_mesh_pt;
297 
298  /// Trace file
299  ofstream Trace_file;
300 
301 }; // end of problem class
302 
303 
304 
305 //=======start_of_constructor=============================================
306 /// Constructor for Helmholtz problem
307 //========================================================================
308 template<class ELEMENT>
311 {
312 
313  // Open trace file
314  Trace_file.open("RESLT/trace.dat");
315 
316  // Setup "bulk" mesh
317 
318  // Inner radius
319  double a=1.0;
320 
321  // Thickness of annular computational domain
322  double h=0.5;
323 
324  // Set outer radius
326 
327  // Create circles representing inner and outer boundary
328  double x_c=0.0;
329  double y_c=0.0;
330  Circle* inner_circle_pt=new Circle(x_c,y_c,a);
331  Circle* outer_circle_pt=new Circle(x_c,y_c,
333 
334  // Outer boundary
335  //---------------
336  TriangleMeshClosedCurve* outer_boundary_pt=0;
337 
338  unsigned n_segments=40;
339  Vector<TriangleMeshCurveSection*> outer_boundary_line_pt(2);
340 
341  // The intrinsic coordinates for the beginning and end of the curve
342  double s_start = 0.0;
343  double s_end = MathematicalConstants::Pi;
344  unsigned boundary_id = 2;
345  outer_boundary_line_pt[0]=
346  new TriangleMeshCurviLine(outer_circle_pt,
347  s_start,
348  s_end,
349  n_segments,
350  boundary_id);
351 
352  // The intrinsic coordinates for the beginning and end of the curve
353  s_start = MathematicalConstants::Pi;
354  s_end = 2.0*MathematicalConstants::Pi;
355  boundary_id = 3;
356  outer_boundary_line_pt[1]=
357  new TriangleMeshCurviLine(outer_circle_pt,
358  s_start,
359  s_end,
360  n_segments,
361  boundary_id);
362 
363  // Create closed curve for outer boundary
364  outer_boundary_pt=new TriangleMeshClosedCurve(outer_boundary_line_pt);
365 
366  // Inner circular boundary
367  //------------------------
368  Vector<TriangleMeshCurveSection*> inner_boundary_line_pt(2);
369 
370  // The intrinsic coordinates for the beginning and end of the curve
371  s_start = 0.0;
372  s_end = MathematicalConstants::Pi;
373  boundary_id = 0;
374  inner_boundary_line_pt[0]=
375  new TriangleMeshCurviLine(inner_circle_pt,
376  s_start,
377  s_end,
378  n_segments,
379  boundary_id);
380 
381  // The intrinsic coordinates for the beginning and end of the curve
382  s_start = MathematicalConstants::Pi;
383  s_end = 2.0*MathematicalConstants::Pi;
384  boundary_id = 1;
385  inner_boundary_line_pt[1]=
386  new TriangleMeshCurviLine(inner_circle_pt,
387  s_start,
388  s_end,
389  n_segments,
390  boundary_id);
391 
392 
393  // Combine to hole
394  //----------------
395  Vector<TriangleMeshClosedCurve*> hole_pt(1);
396  Vector<double> hole_coords(2);
397  hole_coords[0]=0.0;
398  hole_coords[1]=0.0;
399  hole_pt[0]=new TriangleMeshClosedCurve(inner_boundary_line_pt,
400  hole_coords);
401 
402 
403 
404  // Use the TriangleMeshParameters object for helping on the manage of the
405  // TriangleMesh parameters. The only parameter that needs to take is the
406  // outer boundary.
407  TriangleMeshParameters triangle_mesh_parameters(outer_boundary_pt);
408 
409  // Specify the closed curve using the TriangleMeshParameters object
410  triangle_mesh_parameters.internal_closed_curve_pt() = hole_pt;
411 
412 #ifdef ADAPTIVE
413 
414  // Build "bulk" mesh
415  Bulk_mesh_pt=new RefineableTriangleMesh<ELEMENT>(triangle_mesh_parameters);
416 
417  // Create/set error estimator
418  Bulk_mesh_pt->spatial_error_estimator_pt()=new Z2ErrorEstimator;
419 
420  // Choose error tolerances to force some uniform refinement
421  Bulk_mesh_pt->min_permitted_error()=0.00004;
422  Bulk_mesh_pt->max_permitted_error()=0.0001;
423 
424 #else
425 
426  // Build "bulk" mesh
427  Bulk_mesh_pt=new TriangleMesh<ELEMENT>(triangle_mesh_parameters);
428 
429 #endif
430 
431  // Pointer to mesh containing the Helmholtz outer boundary condition
432  // elements. Specify outer radius and number of Fourier terms to be
433  // used in gamma integral
434  Helmholtz_outer_boundary_mesh_pt =
435  new HelmholtzDtNMesh<ELEMENT>(a+h,GlobalParameters::N_fourier);
436 
437  // Create outer boundary elements from all elements that are
438  // adjacent to the outer boundary , but add them to a separate mesh.
439  create_outer_bc_elements(2,Bulk_mesh_pt,Helmholtz_outer_boundary_mesh_pt);
440  create_outer_bc_elements(3,Bulk_mesh_pt,Helmholtz_outer_boundary_mesh_pt);
441 
442  // Pointer to mesh containing the Helmholtz inner boundary condition
443  // elements. Specify outer radius
444  Helmholtz_inner_boundary_mesh_pt = new Mesh;
445 
446  // Create prescribed-flux elements from all elements that are
447  // adjacent to the inner boundary , but add them to a separate mesh.
448  create_flux_elements(0,Bulk_mesh_pt,Helmholtz_inner_boundary_mesh_pt);
449  create_flux_elements(1,Bulk_mesh_pt,Helmholtz_inner_boundary_mesh_pt);
450 
451  // Add the several sub meshes to the problem
452  add_sub_mesh(Bulk_mesh_pt);
453  add_sub_mesh(Helmholtz_outer_boundary_mesh_pt);
454  add_sub_mesh(Helmholtz_inner_boundary_mesh_pt);
455 
456  // Build the Problem's global mesh from its various sub-meshes
457  build_global_mesh();
458 
459  // Complete the build of all elements so they are fully functional
460 
461  // Loop over the Helmholtz bulk elements to set up element-specific
462  // things that cannot be handled by constructor: Pass pointer to
463  // wave number squared
464  unsigned n_element = Bulk_mesh_pt->nelement();
465  for(unsigned e=0;e<n_element;e++)
466  {
467  // Upcast from GeneralisedElement to Helmholtz bulk element
468  ELEMENT *el_pt = dynamic_cast<ELEMENT*>(Bulk_mesh_pt->element_pt(e));
469 
470  //Set the k_squared pointer
471  el_pt->k_squared_pt() = &GlobalParameters::K_squared;
472  }
473 
474  // Set up elements on outer boundary
475  setup_outer_boundary();
476 
477  // Set pointer to prescribed flux function for flux elements
478  set_prescribed_incoming_flux_pt();
479 
480  // Setup equation numbering scheme
481  cout <<"Number of equations: " << assign_eqn_numbers() << std::endl;
482 
483 } // end of constructor
484 
485 //=====================start_of_actions_before_adapt======================
486 /// Actions before adapt: Wipe the mesh of face elements
487 //========================================================================
488 template<class ELEMENT>
490 {
491  // Kill the flux elements and wipe the boundary meshs
492  delete_face_elements(Helmholtz_outer_boundary_mesh_pt);
493  delete_face_elements(Helmholtz_inner_boundary_mesh_pt);
494 
495  // Rebuild the Problem's global mesh from its various sub-meshes
496  rebuild_global_mesh();
497 
498 }// end of actions_before_adapt
499 
500 
501 //=====================start_of_actions_after_adapt=======================
502 /// Actions after adapt: Rebuild the face element meshes
503 //========================================================================
504 template<class ELEMENT>
506 {
507 
508 
509  // Complete the build of all elements so they are fully functional
510 
511  // Loop over the Helmholtz bulk elements to set up element-specific
512  // things that cannot be handled by constructor: Pass pointer to
513  // wave number squared
514  unsigned n_element = Bulk_mesh_pt->nelement();
515  for(unsigned e=0;e<n_element;e++)
516  {
517  // Upcast from GeneralisedElement to Helmholtz bulk element
518  ELEMENT *el_pt = dynamic_cast<ELEMENT*>(Bulk_mesh_pt->element_pt(e));
519 
520  //Set the k_squared pointer
521  el_pt->k_squared_pt() = &GlobalParameters::K_squared;
522  }
523 
524  // Create prescribed-flux elements and BC elements
525  // from all elements that are adjacent to the boundaries and add them to
526  // Helmholtz_boundary_meshes
527  create_flux_elements(0,Bulk_mesh_pt,Helmholtz_inner_boundary_mesh_pt);
528  create_flux_elements(1,Bulk_mesh_pt,Helmholtz_inner_boundary_mesh_pt);
529  create_outer_bc_elements(2,Bulk_mesh_pt,Helmholtz_outer_boundary_mesh_pt);
530  create_outer_bc_elements(3,Bulk_mesh_pt,Helmholtz_outer_boundary_mesh_pt);
531 
532  // Rebuild the Problem's global mesh from its various sub-meshes
533  rebuild_global_mesh();
534 
535  // Set pointer to prescribed flux function and DtN mesh
536  setup_outer_boundary();
537  set_prescribed_incoming_flux_pt();
538 
539 
540 }// end of actions_after_adapt
541 
542 
543 //==================start_of_setup_outer_boundary=========================
544 /// Set pointers for elements on outer boundary
545 //========================================================================
546 template<class ELEMENT>
548 {
549  // Loop over the flux elements to pass pointer to DtN
550  // BC for the outer boundary
551  unsigned n_element=Helmholtz_outer_boundary_mesh_pt->nelement();
552  for(unsigned e=0;e<n_element;e++)
553  {
554  // Dirichlet to Neumann BC
556  {
557  // Upcast from GeneralisedElement to Helmholtz flux element
558  HelmholtzDtNBoundaryElement<ELEMENT> *el_pt =
559  dynamic_cast< HelmholtzDtNBoundaryElement<ELEMENT>*>(
560  Helmholtz_outer_boundary_mesh_pt->element_pt(e));
561 
562  // Set pointer to the mesh that contains all the boundary condition
563  // elements on this boundary
564  el_pt->set_outer_boundary_mesh_pt(Helmholtz_outer_boundary_mesh_pt);
565  }
566  // ABCs BC
567  else
568  {
569  // Upcast from GeneralisedElement to appropriate type
570  HelmholtzAbsorbingBCElement<ELEMENT> *el_pt =
571  dynamic_cast< HelmholtzAbsorbingBCElement<ELEMENT>*>(
572  Helmholtz_outer_boundary_mesh_pt->element_pt(e));
573 
574  // Set pointer to outer radius of artificial boundary
575  el_pt->outer_radius_pt()=&GlobalParameters::Outer_radius;
576 
577  // Set order of absorbing boundary condition
578  el_pt->abc_order_pt()=&GlobalParameters::ABC_order;
579  }
580  }
581 }
582 
583 
584 //==================start_of_set_prescribed_incoming_flux_pt==============
585 /// Set pointer to prescribed incoming-flux function for all
586 /// elements in the inner boundary
587 //========================================================================
588 template<class ELEMENT>
590 {
591  // Loop over the flux elements to pass pointer to prescribed
592  // flux function for the inner boundary
593  unsigned n_element=Helmholtz_inner_boundary_mesh_pt->nelement();
594  for(unsigned e=0;e<n_element;e++)
595  {
596  // Upcast from GeneralisedElement to Helmholtz flux element
597  HelmholtzFluxElement<ELEMENT> *el_pt =
598  dynamic_cast< HelmholtzFluxElement<ELEMENT>*>(
599  Helmholtz_inner_boundary_mesh_pt->element_pt(e));
600 
601  // Set the pointer to the prescribed flux function
602  el_pt->flux_fct_pt() =
604  }
605 
606 }// end of set prescribed_incoming_flux pt
607 
608 //=====================start_of_doc=======================================
609 /// Doc the solution: doc_info contains labels/output directory etc.
610 //========================================================================
611 template<class ELEMENT>
613  doc_info)
614 {
615 
616  ofstream some_file,some_file2;
617  char filename[100];
618 
619  // Number of plot points
620  unsigned npts;
621  npts=5;
622 
623  // Total radiated power
624  double power=0.0;
625 
626  // Compute/output the radiated power
627  //----------------------------------
628  sprintf(filename,"%s/power%i.dat",doc_info.directory().c_str(),
629  doc_info.number());
630  some_file.open(filename);
631 
632  // Accumulate contribution from elements
633  unsigned nn_element=Helmholtz_outer_boundary_mesh_pt->nelement();
634  for(unsigned e=0;e<nn_element;e++)
635  {
636  HelmholtzBCElementBase<ELEMENT> *el_pt =
637  dynamic_cast< HelmholtzBCElementBase<ELEMENT>*>(
638  Helmholtz_outer_boundary_mesh_pt->element_pt(e));
639  power += el_pt->global_power_contribution(some_file);
640  }
641  some_file.close();
642  oomph_info << "Total radiated power: " << power << std::endl;
643 
644  // Output solution
645  //-----------------
646  sprintf(filename,"%s/soln%i.dat",doc_info.directory().c_str(),
647  doc_info.number());
648  some_file.open(filename);
649  Bulk_mesh_pt->output(some_file,npts);
650  some_file.close();
651 
652  // Output exact solution
653  //----------------------
654  sprintf(filename,"%s/exact_soln%i.dat",doc_info.directory().c_str(),
655  doc_info.number());
656  some_file.open(filename);
657  Bulk_mesh_pt->output_fct(some_file,npts,GlobalParameters::get_exact_u);
658  some_file.close();
659 
660  double error,norm;
661  sprintf(filename,"%s/error%i.dat",doc_info.directory().c_str(),
662  doc_info.number());
663  some_file.open(filename);
664  Bulk_mesh_pt->compute_error(some_file,GlobalParameters::get_exact_u,
665  error,norm);
666  some_file.close();
667 
668  // Doc L2 error and norm of solution
669  oomph_info << "\nNorm of error : " << sqrt(error) << std::endl;
670  oomph_info << "Norm of solution: " << sqrt(norm) << std::endl << std::endl;
671 
672 
673  // Write power to trace file
674  Trace_file << power << std::endl;
675 
676  // Do animation of Helmholtz solution
677  //-----------------------------------
678  unsigned nstep=40;
679  for (unsigned i=0;i<nstep;i++)
680  {
681  sprintf(filename,"%s/helmholtz_animation%i_frame%i.dat",
682  doc_info.directory().c_str(),
683  doc_info.number(),i);
684  some_file.open(filename);
685  sprintf(filename,"%s/exact_helmholtz_animation%i_frame%i.dat",
686  doc_info.directory().c_str(),
687  doc_info.number(),i);
688  some_file2.open(filename);
689  double phi=2.0*MathematicalConstants::Pi*double(i)/double(nstep-1);
690  unsigned nelem=Bulk_mesh_pt->nelement();
691  for (unsigned e=0;e<nelem;e++)
692  {
693  ELEMENT* el_pt=dynamic_cast<ELEMENT*>(
694  Bulk_mesh_pt->element_pt(e));
695  el_pt->output_real(some_file,phi,npts);
696  el_pt->output_real_fct(some_file2,phi,npts,
698  }
699  some_file.close();
700  some_file2.close();
701  }
702 
703 } // end of doc
704 
705 //============start_of_create_flux_elements==================================
706 /// Create Helmholtz inner Flux Elements on the b-th boundary of
707 /// the Mesh object pointed to by bulk_mesh_pt and add the elements
708 /// to the Mesh object pointed to by helmholtz_inner_boundary_mesh_pt
709 //============================================================================
710 template<class ELEMENT>
712 create_flux_elements(const unsigned &b, Mesh* const &bulk_mesh_pt,
713  Mesh* const & helmholtz_inner_boundary_mesh_pt)
714 {
715  // Loop over the bulk elements adjacent to boundary b
716  unsigned n_element = bulk_mesh_pt->nboundary_element(b);
717  for(unsigned e=0;e<n_element;e++)
718  {
719  // Get pointer to the bulk element that is adjacent to boundary b
720  ELEMENT* bulk_elem_pt = dynamic_cast<ELEMENT*>(
721  bulk_mesh_pt->boundary_element_pt(b,e));
722 
723  //Find the index of the face of element e along boundary b
724  int face_index = bulk_mesh_pt->face_index_at_boundary(b,e);
725 
726  // Build the corresponding prescribed incoming-flux element
727  HelmholtzFluxElement<ELEMENT>* flux_element_pt = new
728  HelmholtzFluxElement<ELEMENT>(bulk_elem_pt,face_index);
729 
730  //Add the prescribed incoming-flux element to the surface mesh
731  helmholtz_inner_boundary_mesh_pt->add_element_pt(flux_element_pt);
732 
733  } //end of loop over bulk elements adjacent to boundary b
734 
735 } // end of create_flux_elements
736 
737 
738 
739 //============start_of_create_outer_bc_elements==============================
740 /// Create outer BC elements on the b-th boundary of
741 /// the Mesh object pointed to by bulk_mesh_pt and add the elements
742 /// to the Mesh object pointed to by helmholtz_outer_boundary_mesh_pt.
743 //===========================================================================
744 template<class ELEMENT>
746 create_outer_bc_elements(const unsigned &b, Mesh* const &bulk_mesh_pt,
747  Mesh* const & helmholtz_outer_boundary_mesh_pt)
748 {
749  // Loop over the bulk elements adjacent to boundary b?
750  unsigned n_element = bulk_mesh_pt->nboundary_element(b);
751  for(unsigned e=0;e<n_element;e++)
752  {
753  // Get pointer to the bulk element that is adjacent to boundary b
754  ELEMENT* bulk_elem_pt = dynamic_cast<ELEMENT*>(
755  bulk_mesh_pt->boundary_element_pt(b,e));
756 
757  //Find the index of the face of element e along boundary b
758  int face_index = bulk_mesh_pt->face_index_at_boundary(b,e);
759 
760  // Build the corresponding outer flux element
761 
762  // Dirichlet to Neumann boundary conditon
764  {
765  HelmholtzDtNBoundaryElement<ELEMENT>* flux_element_pt = new
766  HelmholtzDtNBoundaryElement<ELEMENT>(bulk_elem_pt,face_index);
767 
768  //Add the flux boundary element to the helmholtz_outer_boundary_mesh
769  helmholtz_outer_boundary_mesh_pt->add_element_pt(flux_element_pt);
770  }
771  // ABCs BC
772  else
773  {
774  HelmholtzAbsorbingBCElement<ELEMENT>* flux_element_pt = new
775  HelmholtzAbsorbingBCElement<ELEMENT>(bulk_elem_pt,face_index);
776 
777  //Add the flux boundary element to the helmholtz_outer_boundary_mesh
778  helmholtz_outer_boundary_mesh_pt->add_element_pt(flux_element_pt);
779  }
780  } //end of loop over bulk elements adjacent to boundary b
781 } // end of create_outer_bc_elements
782 
783 
784 //============start_of_delete_face_elements================
785 /// Delete face elements and wipe the boundary mesh
786 //==========================================================
787 template<class ELEMENT>
789 delete_face_elements(Mesh* const & boundary_mesh_pt)
790 {
791  // Loop over the surface elements
792  unsigned n_element = boundary_mesh_pt->nelement();
793  for(unsigned e=0;e<n_element;e++)
794  {
795  // Kill surface element
796  delete boundary_mesh_pt->element_pt(e);
797  }
798 
799  // Wipe the mesh
800  boundary_mesh_pt->flush_element_and_node_storage();
801 
802 } // end of delete_outer_face_elements
803 
804 
805 
806 //==========start_of_main=================================================
807 /// Solve 2D Helmholtz problem for scattering of a planar wave from a
808 /// unit disk
809 //========================================================================
810 int main(int argc, char **argv)
811 {
812  // Store command line arguments
813  CommandLineArgs::setup(argc,argv);
814 
815  // Define case to be run
816  unsigned i_case=0;
817  CommandLineArgs::specify_command_line_flag("--case",&i_case);
818 
819  // Parse command line
820  CommandLineArgs::parse_and_assign();
821 
822  // Doc what has actually been specified on the command line
823  CommandLineArgs::doc_specified_flags();
824 
825  // Now set flags accordingly
826  switch(i_case)
827  {
828  case 0:
830  break;
831 
832  case 1:
835  break;
836 
837  case 2:
840  break;
841 
842  case 3:
845  break;
846  }
847 
848 
849  //Set up the problem
850  //------------------
851 
852 #ifdef ADAPTIVE
853 
854  //Set up the problem with 2D nine-node elements from the
855  //QHelmholtzElement family.
857  problem;
858 
859 #else
860 
861  // Set up the problem with 2D six-node elements from the
862  // THelmholtzElement family.
864 
865 
866 #endif
867 
868  // Create label for output
869  //------------------------
870  DocInfo doc_info;
871 
872  // Set output directory
873  doc_info.set_directory("RESLT");
874 
875 
876 #ifdef ADAPTIVE
877 
878  // Max. number of adaptations
879  unsigned max_adapt=1;
880 
881  // Solve the problem with Newton's method, allowing
882  // up to max_adapt mesh adaptations after every solve.
883  problem.newton_solve(max_adapt);
884 
885 #else
886 
887  // Solve the problem with Newton's method
888  problem.newton_solve();
889 
890 #endif
891 
892  //Output solution
893  problem.doc_solution(doc_info);
894 
895 } //end of main
896 
897 
////////////////////////////////////////////////////////////////// //////////////////////////////////...
Definition: scattering.cc:223
~ScatteringProblem()
Destructor (empty)
void actions_after_adapt()
Actions after adapt: Rebuild the mesh of prescribed flux elements.
void create_outer_bc_elements(const unsigned &b, Mesh *const &bulk_mesh_pt, Mesh *const &helmholtz_outer_boundary_mesh_pt)
Create BC elements on boundary b of the Mesh pointed to by bulk_mesh_pt and add them to the specified...
ofstream Trace_file
Trace file.
TriangleMesh< ELEMENT > * Bulk_mesh_pt
Pointer to the "bulk" mesh.
void setup_outer_boundary()
Set up boundary condition elements on outer boundary.
void actions_before_newton_solve()
Update the problem specs before solve (empty)
ScatteringProblem()
Constructor.
void delete_face_elements(Mesh *const &boundary_mesh_pt)
Delete boundary face elements and wipe the surface mesh.
void actions_after_newton_solve()
Update the problem specs after solve (empty)
void actions_before_adapt()
Actions before adapt: Wipe the mesh of prescribed flux elements.
void actions_before_newton_convergence_check()
Recompute gamma integral before checking Newton residuals.
void set_prescribed_incoming_flux_pt()
Set pointer to prescribed-flux function for all elements in the surface mesh on the surface of the un...
RefineableTriangleMesh< ELEMENT > * Bulk_mesh_pt
Pointer to the "bulk" mesh.
void create_flux_elements(const unsigned &b, Mesh *const &bulk_mesh_pt, Mesh *const &helmholtz_inner_boundary_mesh_pt)
Create Helmholtz flux elements on boundary b of the Mesh pointed to by bulk_mesh_pt and add them to t...
void doc_solution(DocInfo &doc_info)
Doc the solution. DocInfo object stores flags/labels for where the output gets written to.
Namespace for "global" problem parameters.
Definition: barrel.cc:45
void prescribed_incoming_flux(const Vector< double > &x, complex< double > &flux)
Flux (normal derivative) on the unit disk for a planar incoming wave.
Definition: scattering.cc:156
unsigned ABC_order
Flag to choose wich order to use.
Definition: scattering.cc:76
double Outer_radius
Radius of outer boundary (must be a circle!)
Definition: scattering.cc:79
bool DtN_BC
Flag to choose the Dirichlet to Neumann BC or ABC BC.
Definition: scattering.cc:72
std::complex< double > I(0.0, 1.0)
Imaginary unit.
double K_squared
Square of the wavenumber.
Definition: scattering.cc:64
void get_exact_u(const Vector< double > &x, Vector< double > &u)
Exact solution for scattered field (vector returns real and impaginary parts).
Definition: scattering.cc:86
unsigned N_fourier
Number of terms used in the computation of the exact solution.
Definition: scattering.cc:68
int main(int argc, char **argv)
Solve 2D Helmholtz problem for scattering of a planar wave from a unit disk.