oscillating_sphere.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-2023 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 pml Fourier-decomposed Helmholtz problem
27 
28 #include <complex>
29 #include <cmath>
30 
31 //Generic routines
32 #include "generic.h"
33 
34 // The Helmholtz equations
35 #include "pml_fourier_decomposed_helmholtz.h"
36 
37 // The mesh
38 #include "meshes/triangle_mesh.h"
39 
40 // Get the Bessel functions
41 #include "oomph_crbond_bessel.h"
42 
43 using namespace oomph;
44 using namespace std;
45 
46 
47 /// //////////////////////////////////////////////////////////////////
48 /// //////////////////////////////////////////////////////////////////
49 /// //////////////////////////////////////////////////////////////////
50 
51 
52 
53 //===== start_of_namespace=============================================
54 /// Namespace for the Fourier decomposed Helmholtz problem parameters
55 //=====================================================================
57 {
58  /// Output directory
59  string Directory="RESLT";
60 
61  /// Frequency
62  double K_squared = 10.0;
63 
64  /// Default physical PML thickness
65  double PML_thickness=4.0;
66 
67  /// Default number of elements within PMLs
68  unsigned Nel_pml=15;
69 
70  /// Target area for initial mesh
71  double Element_area = 0.1;
72 
73  /// The default Fourier wave number
74  int N_fourier=0;
75 
76  /// Number of terms in the exact solution
77  unsigned N_terms=6;
78 
79  /// Coefficients in the exact solution
80  Vector<double> Coeff(N_terms,1.0);
81 
82  /// Imaginary unit
83  std::complex<double> I(0.0,1.0);
84 
85  /// Exact solution as a Vector of size 2, containing real and imag parts
86  void get_exact_u(const Vector<double>& x, Vector<double>& u)
87  {
88  double K = sqrt(K_squared);
89 
90  // Switch to spherical coordinates
91  double R=sqrt(x[0]*x[0]+x[1]*x[1]);
92 
93  double theta;
94  theta=atan2(x[0],x[1]);
95 
96  // Argument for Bessel/Hankel functions
97  double kr= K*R;
98 
99  // Need half-order Bessel functions
100  double bessel_offset=0.5;
101 
102  // Evaluate Bessel/Hankel functions
103  Vector<double> jv(N_terms);
104  Vector<double> yv(N_terms);
105  Vector<double> djv(N_terms);
106  Vector<double> dyv(N_terms);
107  double order_max_in=double(N_terms-1)+bessel_offset;
108  double order_max_out=0;
109 
110  // This function returns vectors containing
111  // J_k(x), Y_k(x) and their derivatives
112  // up to k=order_max, with k increasing in
113  // integer increments starting with smallest
114  // positive value. So, e.g. for order_max=3.5
115  // jv[0] contains J_{1/2}(x),
116  // jv[1] contains J_{3/2}(x),
117  // jv[2] contains J_{5/2}(x),
118  // jv[3] contains J_{7/2}(x).
119  CRBond_Bessel::bessjyv(order_max_in,
120  kr,
121  order_max_out,
122  &jv[0],&yv[0],
123  &djv[0],&dyv[0]);
124 
125  // Assemble exact solution (actually no need to add terms
126  // below i=N_fourier as Legendre polynomial would be zero anyway)
127  complex<double> u_ex(0.0,0.0);
128  for(unsigned i=N_fourier;i<N_terms;i++)
129  {
130  //Associated_legendre_functions
131  double p=Legendre_functions_helper::plgndr2(i,N_fourier,
132  cos(theta));
133  // Set exact solution
134  u_ex+=Coeff[i]*sqrt(MathematicalConstants::Pi/(2.0*kr))*(jv[i]+I*yv[i])*p;
135  }
136 
137  // Get the real & imaginary part of the result
138  u[0]=u_ex.real();
139  u[1]=u_ex.imag();
140 
141  }//end of get_exact_u
142 
143 
144 
145  /// Get -du/dr (spherical r) for exact solution. Equal to prescribed
146  /// flux on inner boundary.
147  void exact_minus_dudr(const Vector<double>& x, std::complex<double>& flux)
148  {
149  double K = sqrt(K_squared);
150 
151  // Initialise flux
152  flux=std::complex<double>(0.0,0.0);
153 
154  // Switch to spherical coordinates
155  double R=sqrt(x[0]*x[0]+x[1]*x[1]);
156 
157  double theta;
158  theta=atan2(x[0],x[1]);
159 
160  // Argument for Bessel/Hankel functions
161  double kr=K*R;
162 
163  // Need half-order Bessel functions
164  double bessel_offset=0.5;
165 
166  // Evaluate Bessel/Hankel functions
167  Vector<double> jv(N_terms);
168  Vector<double> yv(N_terms);
169  Vector<double> djv(N_terms);
170  Vector<double> dyv(N_terms);
171  double order_max_in=double(N_terms-1)+bessel_offset;
172  double order_max_out=0;
173 
174  // This function returns vectors containing
175  // J_k(x), Y_k(x) and their derivatives
176  // up to k=order_max, with k increasing in
177  // integer increments starting with smallest
178  // positive value. So, e.g. for order_max=3.5
179  // jv[0] contains J_{1/2}(x),
180  // jv[1] contains J_{3/2}(x),
181  // jv[2] contains J_{5/2}(x),
182  // jv[3] contains J_{7/2}(x).
183  CRBond_Bessel::bessjyv(order_max_in,
184  kr,
185  order_max_out,
186  &jv[0],&yv[0],
187  &djv[0],&dyv[0]);
188 
189 
190  // Assemble exact solution (actually no need to add terms
191  // below i=N_fourier as Legendre polynomial would be zero anyway)
192  complex<double> u_ex(0.0,0.0);
193  for(unsigned i=N_fourier;i<N_terms;i++)
194  {
195  //Associated_legendre_functions
196  double p=Legendre_functions_helper::plgndr2(i,N_fourier,
197  cos(theta));
198  // Set flux of exact solution
199  flux-=Coeff[i]*sqrt(MathematicalConstants::Pi/(2.0*kr))*p*
200  ( K*(djv[i]+I*dyv[i]) - (0.5*(jv[i]+I*yv[i])/R) );
201  }
202  } // end of exact_normal_derivative
203 
204 
205  /// Radial position of point source
206  double R_source = 2.0;
207 
208  /// Axial position of point source
209  double Z_source = 2.0;
210 
211  /// Point source magnitude (Complex)
212  std::complex<double> Magnitude(100.0,100.0);
213 
214 } // end of namespace
215 
216 
217 /// //////////////////////////////////////////////////////////////////
218 /// //////////////////////////////////////////////////////////////////
219 /// //////////////////////////////////////////////////////////////////
220 
221 
222 namespace oomph
223 {
224 
225 //========= start_of_point_source_wrapper==============================
226 /// Class to impose point source to (wrapped) Helmholtz element
227 //=====================================================================
228 template<class ELEMENT>
229 class PMLHelmholtzPointSourceElement : public virtual ELEMENT
230 {
231 
232 public:
233 
234  /// Constructor
236  {
237  // Initialise
238  Point_source_magnitude=std::complex<double>(0.0,0.0);
239  }
240 
241  /// Destructor (empty)
243 
244  /// Set local coordinate and magnitude of point source
245  void setup(const Vector<double>& s_point_source,
246  const std::complex<double>& magnitude)
247  {
248  S_point_source=s_point_source;
249  Point_source_magnitude=magnitude;
250  }
251 
252 
253  /// Add the element's contribution to its residual vector (wrapper)
254  void fill_in_contribution_to_residuals(Vector<double> &residuals)
255  {
256  //Call the generic residuals function with flag set to 0
257  //using a dummy matrix argument
258  ELEMENT::fill_in_generic_residual_contribution_pml_fourier_decomposed_helmholtz(
259  residuals,GeneralisedElement::Dummy_matrix,0);
260 
261  // Add point source contribution
262  fill_in_point_source_contribution_to_residuals(residuals);
263  }
264 
265 
266  /// Add the element's contribution to its residual vector and
267  /// element Jacobian matrix (wrapper)
268  void fill_in_contribution_to_jacobian(Vector<double> &residuals,
269  DenseMatrix<double> &jacobian)
270  {
271  //Call the generic routine with the flag set to 1
272  ELEMENT::fill_in_generic_residual_contribution_pml_fourier_decomposed_helmholtz(residuals,
273  jacobian,1);
274 
275  // Add point source contribution
276  fill_in_point_source_contribution_to_residuals(residuals);
277  }
278 
279 
280 private:
281 
282 
283  /// Add the point source contribution to the residual vector
284  void fill_in_point_source_contribution_to_residuals(Vector<double> &residuals)
285  {
286  // No further action
287  if (S_point_source.size()==0) return;
288 
289  //Find out how many nodes there are
290  const unsigned n_node = this->nnode();
291 
292  //Set up memory for the shape/test functions
293  Shape psi(n_node);
294 
295  //Integers to store the local equation and unknown numbers
296  int local_eqn_real=0;
297  int local_eqn_imag=0;
298 
299  // Get shape/test fcts
300  this->shape(S_point_source,psi);
301 
302  // Assemble residuals
303  //--------------------------------
304 
305  // Loop over the test functions
306  for(unsigned l=0;l<n_node;l++)
307  {
308  // first, compute the real part contribution
309  //-------------------------------------------
310 
311  //Get the local equation
312  local_eqn_real =
313  this->nodal_local_eqn
314  (l,this->u_index_pml_fourier_decomposed_helmholtz().real());
315 
316  /*IF it's not a boundary condition*/
317  if(local_eqn_real >= 0)
318  {
319  residuals[local_eqn_real] -= 2.0*Point_source_magnitude.real()*psi(l);
320  }
321 
322  // Second, compute the imaginary part contribution
323  //------------------------------------------------
324 
325  //Get the local equation
326  local_eqn_imag =
327  this->nodal_local_eqn
328  (l,this->u_index_pml_fourier_decomposed_helmholtz().imag());
329 
330  /*IF it's not a boundary condition*/
331  if(local_eqn_imag >= 0)
332  {
333  // Add body force/source term and Helmholtz bit
334  residuals[local_eqn_imag] -= 2.0*Point_source_magnitude.imag()*psi(l);
335  }
336  }
337  }
338 
339  /// Local coordinates of point at which point source is applied
340  Vector<double> S_point_source;
341 
342  /// Magnitude of point source (complex!)
343  std::complex<double> Point_source_magnitude;
344 
345  };
346 
347 
348 
349 
350 //=======================================================================
351 /// Face geometry for element is the same as that for the underlying
352 /// wrapped element
353 //=======================================================================
354  template<class ELEMENT>
355  class FaceGeometry<PMLHelmholtzPointSourceElement<ELEMENT> >
356  : public virtual FaceGeometry<ELEMENT>
357  {
358  public:
359  FaceGeometry() : FaceGeometry<ELEMENT>() {}
360  };
361 
362 
363 //=======================================================================
364 /// Face geometry of the Face Geometry for element is the same as
365 /// that for the underlying wrapped element
366 //=======================================================================
367  template<class ELEMENT>
368  class FaceGeometry<FaceGeometry<PMLHelmholtzPointSourceElement<ELEMENT> > >
369  : public virtual FaceGeometry<FaceGeometry<ELEMENT> >
370  {
371  public:
372  FaceGeometry() : FaceGeometry<FaceGeometry<ELEMENT> >() {}
373  };
374 
375 
376 //=======================================================================
377 /// Policy class defining the elements to be used in the actual
378 /// PML layers.
379 //=======================================================================
380  template<class ELEMENT>
381 class PMLLayerElement<PMLHelmholtzPointSourceElement<ELEMENT> > :
382  public virtual PMLLayerElement<ELEMENT>
383 {
384 
385  public:
386 
387  /// Constructor: Call the constructor for the
388  /// appropriate Element
389  PMLLayerElement() : PMLLayerElement<ELEMENT>()
390  {}
391 
392 };
393 
394 
395 
396 //=======================================================================
397 /// Policy class defining the elements to be used in the actual
398 /// PML layers.
399 //=======================================================================
400  template<class ELEMENT>
401 class PMLLayerElement<
402  ProjectablePMLFourierDecomposedHelmholtzElement<
403  PMLHelmholtzPointSourceElement<ELEMENT> > >:
404  public virtual PMLLayerElement<ELEMENT>
405 {
406 
407  public:
408 
409  /// Constructor: Call the constructor for the
410  /// appropriate Element
411  PMLLayerElement() : PMLLayerElement<ELEMENT>()
412  {}
413 
414 };
415 
416 }
417 
418 
419 
420 
421 /// //////////////////////////////////////////////////////////////////
422 /// //////////////////////////////////////////////////////////////////
423 /// //////////////////////////////////////////////////////////////////
424 
425 //========= start_of_problem_class=====================================
426 /// Problem class
427 //=====================================================================
428 template<class ELEMENT>
430 {
431 
432 public:
433 
434  /// Constructor
436 
437  /// Destructor (empty)
439 
440  /// Update the problem specs before solve (empty)
442 
443  /// Update the problem after solve (empty)
445 
446  /// Doc the solution. DocInfo object stores flags/labels for where the
447  /// output gets written to
448  void doc_solution(DocInfo& doc_info);
449 
450  /// Create PML meshes
451  void create_pml_meshes();
452 
453  /// Create mesh of face elements that monitor the radiated power
454  void create_power_monitor_mesh();
455 
456  #ifdef ADAPTIVE
457 
458  /// Actions before adapt: Wipe the mesh of prescribed flux elements
459  void actions_before_adapt();
460 
461  /// Actions after adapt: Rebuild the mesh of prescribed flux elements
462  void actions_after_adapt();
463 
464  #endif
465 
466 
467  // Apply boundary condtions for odd Fourier wavenumber
468  void complete_problem_setup();
469 
470 private:
471 
472  /// Create flux elements on inner boundary
473  void create_flux_elements_on_inner_boundary();
474 
475  /// Delete boundary face elements and wipe the surface mesh
476  void delete_face_elements(Mesh* const & boundary_mesh_pt)
477  {
478  // Loop over the surface elements
479  unsigned n_element = boundary_mesh_pt->nelement();
480  for(unsigned e=0;e<n_element;e++)
481  {
482  // Kill surface element
483  delete boundary_mesh_pt->element_pt(e);
484  }
485 
486  // Wipe the mesh
487  boundary_mesh_pt->flush_element_and_node_storage();
488  }
489 
490  // Apply boundary condtions for odd Fourier wavenumber
491  void apply_zero_dirichlet_boundary_conditions();
492 
493  /// Pointer to mesh that stores the power monitor elements
495 
496 
497 #ifdef ADAPTIVE
498 
499  // Create point source element (only used in adaptive run)
500  void setup_point_source();
501 
502  /// Pointer to the refineable "bulk" mesh
503  RefineableTriangleMesh<ELEMENT>* Bulk_mesh_pt;
504 
505  /// Mesh as geometric object representation of bulk mesh
506  MeshAsGeomObject* Mesh_as_geom_obj_pt;
507 
508 #else
509 
510  /// Pointer to the "bulk" mesh
511  TriangleMesh<ELEMENT>* Bulk_mesh_pt;
512 
513 #endif
514 
515  /// Mesh of FaceElements that apply the flux bc on the inner boundary
517 
518  /// Pointer to the right PML mesh
520 
521  /// Pointer to the top PML mesh
523 
524  /// Pointer to the bottom PML mesh
526 
527  /// Pointer to the top right corner PML mesh
529 
530  /// Pointer to the bottom right corner PML mesh
532 
533  /// Trace file
534  ofstream Trace_file;
535 
536 }; // end of problem class
537 
538 
539 //===================start_of_create_power_monitor_mesh===================
540 /// Create BC elements on outer boundary
541 //========================================================================
542 template<class ELEMENT>
545 {
546  // Loop over outer boundaries
547  for (unsigned b=1;b<4;b++)
548  {
549  // Loop over the bulk elements adjacent to boundary b?
550  unsigned n_element = Bulk_mesh_pt->nboundary_element(b);
551  for(unsigned e=0;e<n_element;e++)
552  {
553  // Get pointer to the bulk element that is adjacent to boundary b
554  ELEMENT* bulk_elem_pt = dynamic_cast<ELEMENT*>(
555  Bulk_mesh_pt->boundary_element_pt(b,e));
556 
557  //Find the index of the face of element e along boundary b
558  int face_index = Bulk_mesh_pt->face_index_at_boundary(b,e);
559 
560  // Build the corresponding element
561  PMLFourierDecomposedHelmholtzPowerMonitorElement<ELEMENT>*
562  flux_element_pt = new
563  PMLFourierDecomposedHelmholtzPowerMonitorElement<ELEMENT>
564  (bulk_elem_pt,face_index);
565 
566  //Add the flux boundary element
567  Power_monitor_mesh_pt->add_element_pt(flux_element_pt);
568  }
569  }
570 } // end of create_power_monitor_mesh
571 
572 
573 
574 #ifdef ADAPTIVE
575 
576 //=====================start_of_actions_before_adapt======================
577 /// Actions before adapt: Wipe the mesh of face elements
578 //========================================================================
579 template<class ELEMENT>
582 {
583  // Before adapting the added PML meshes must be removed
584  // as they are not refineable and are to be rebuilt from the
585  // newly refined triangular mesh
586  delete PML_right_mesh_pt;
587  PML_right_mesh_pt=0;
588  delete PML_top_mesh_pt;
589  PML_top_mesh_pt=0;
590  delete PML_bottom_mesh_pt;
591  PML_bottom_mesh_pt=0;
592  delete PML_top_right_mesh_pt;
593  PML_top_right_mesh_pt=0;
594  delete PML_bottom_right_mesh_pt;
595  PML_bottom_right_mesh_pt=0;
596 
597  // Wipe the power monitor elements
598  delete_face_elements(Power_monitor_mesh_pt);
599 
600  // Rebuild the Problem's global mesh from its various sub-meshes
601  // but first flush all its submeshes
602  flush_sub_meshes();
603 
604  // Then add the triangular mesh back
605  add_sub_mesh(Bulk_mesh_pt);
606 
607  // Rebuild the Problem's global mesh from its various sub-meshes
608  rebuild_global_mesh();
609 
610 }// end of actions_before_adapt
611 
612 
613 
614 //=====================start_of_actions_after_adapt=======================
615 /// Actions after adapt: Rebuild the face element meshes
616 //========================================================================
617 template<class ELEMENT>
620 {
621 
622  // Build PML meshes and add them to the global mesh
623  create_pml_meshes();
624 
625  // Re-attach the power monitor elements
626  create_power_monitor_mesh();
627 
628  // Build the entire mesh from its submeshes
629  rebuild_global_mesh();
630 
631  // Complete the build of all elements
632  complete_problem_setup();
633 
634  // Setup point source
635  setup_point_source();
636 
637 }// end of actions_after_adapt
638 
639 #endif
640 
641 
642 
643 //=================start_of_complete_problem_setup==================
644 // Complete the build of all elements so that they are fully
645 // functional
646 //==================================================================
647 template<class ELEMENT>
650 {
651  // Complete the build of all elements so they are fully functional
652  unsigned n_element = this->mesh_pt()->nelement();
653  for(unsigned i=0;i<n_element;i++)
654  {
655  // Upcast from GeneralsedElement to the present element
656  PMLFourierDecomposedHelmholtzEquations *el_pt = dynamic_cast<
657  PMLFourierDecomposedHelmholtzEquations*>(
658  mesh_pt()->element_pt(i));
659 
660  if (!(el_pt==0))
661  {
662  //Set the frequency pointer
663  el_pt->k_squared_pt()=&ProblemParameters::K_squared;
664 
665  // Set pointer to Fourier wave number
666  el_pt->pml_fourier_wavenumber_pt()=&ProblemParameters::N_fourier;
667  }
668  }
669 
670  // If the Fourier wavenumber is odd, then apply zero dirichlet boundary
671  // conditions on the two straight boundaries on the symmetry line.
672  if (ProblemParameters::N_fourier % 2 == 1)
673  {
674  cout
675  << "Zero Dirichlet boundary condition has been applied on symmetry line\n";
676  cout << "due to an odd Fourier wavenumber\n" << std::endl;
677  apply_zero_dirichlet_boundary_conditions();
678  }
679 
680 } // end of complete_problem_setup
681 
682 
683 #ifdef ADAPTIVE
684 
685 //==================start_of_setup_point_source===========================
686 /// Set point source
687 //========================================================================
688 template<class ELEMENT>
691 {
692  // Create mesh as geometric object
693  delete Mesh_as_geom_obj_pt;
694  Mesh_as_geom_obj_pt= new MeshAsGeomObject(Bulk_mesh_pt);
695 
696  // Position of point source
697  Vector<double> x_point_source;
698  x_point_source.resize(2);
699  x_point_source[0]=ProblemParameters::R_source;
700  x_point_source[1]=ProblemParameters::Z_source;
701 
702  GeomObject* sub_geom_object_pt=0;
703  Vector<double> s_point_source(2);
704  Mesh_as_geom_obj_pt->locate_zeta(x_point_source,sub_geom_object_pt,
705  s_point_source);
706 
707  // Set point force
708  if(x_point_source[0]==0.0)
709  {
711  {
712  // if source on z axis, only contribution to residual comes
713  // from Fourier wavenumber zero
714  ProblemParameters::Magnitude=complex<double>(0.0,0.0);
715  }
716  }
717 
718  // Set point force
719  dynamic_cast<ELEMENT*>(sub_geom_object_pt)->
720  setup(s_point_source,ProblemParameters::Magnitude);
721 
722 }
723 
724 
725 #endif
726 
727 //=========start_of_apply_zero_dirichlet_boundary_conditions========
728 // Apply extra bounday conditions if given an odd Fourier wavenumber
729 //==================================================================
730 template<class ELEMENT>
733 {
734  // Apply zero dirichlet conditions on the bottom straight boundary
735  // and the top straight boundary located on the symmetry line.
736 
737  // Bottom straight boundary on symmetry line:
738  {
739  //Boundary id
740  unsigned b=0;
741 
742  // How many nodes are there?
743  unsigned n_node=Bulk_mesh_pt->nboundary_node(b);
744  for (unsigned n=0;n<n_node;n++)
745  {
746  // Get the node
747  Node* nod_pt=Bulk_mesh_pt->boundary_node_pt(b,n);
748 
749  // Pin the node
750  nod_pt->pin(0);
751  nod_pt->pin(1);
752 
753  // Set the node's value
754  nod_pt->set_value(0, 0.0);
755  nod_pt->set_value(1, 0.0);
756  }
757  }
758 
759 // Top straight boundary on symmetry line:
760  {
761  //Boundary id
762  unsigned b=4;
763 
764  // How many nodes are there?
765  unsigned n_node=Bulk_mesh_pt->nboundary_node(b);
766  for (unsigned n=0;n<n_node;n++)
767  {
768  // Get the node
769  Node* nod_pt=Bulk_mesh_pt->boundary_node_pt(b,n);
770 
771  // Pin the node
772  nod_pt->pin(0);
773  nod_pt->pin(1);
774 
775  // Set the node's value
776  nod_pt->set_value(0, 0.0);
777  nod_pt->set_value(1, 0.0);
778  }
779  }
780 
781 
782 } // end of apply_zero_dirichlet_boundary_conditions
783 
784 
785 //=======start_of_constructor=============================================
786 /// Constructor for Pml Fourier-decomposed Helmholtz problem
787 //========================================================================
788 template<class ELEMENT>
791 {
792  string trace_file_location = ProblemParameters::Directory + "/trace.dat";
793 
794  // Open trace file
795  Trace_file.open(trace_file_location.c_str());
796 
797  /// Setup "bulk" mesh
798 
799  // Create the circle that represents the inner boundary
800  double x_c=0.0;
801  double y_c=0.0;
802  double r_min=1.0;
803  double r_max=3.0;
804  Circle* inner_circle_pt=new Circle(x_c,y_c,r_min);
805 
806 
807  // Edges/boundary segments making up outer boundary
808  //-------------------------------------------------
809  Vector<TriangleMeshCurveSection*> outer_boundary_line_pt(6);
810 
811  // All poly boundaries are defined by two vertices
812  Vector<Vector<double> > boundary_vertices(2);
813 
814 
815  // Bottom straight boundary on symmetry line
816  //------------------------------------------
817  boundary_vertices[0].resize(2);
818  boundary_vertices[0][0]=0.0;
819  boundary_vertices[0][1]=-r_min;
820  boundary_vertices[1].resize(2);
821  boundary_vertices[1][0]=0.0;
822  boundary_vertices[1][1]=-r_max;
823 
824  unsigned boundary_id=0;
825  outer_boundary_line_pt[0]=
826  new TriangleMeshPolyLine(boundary_vertices,boundary_id);
827 
828 
829  // Bottom boundary of bulk mesh
830  //-----------------------------
831  boundary_vertices[0][0]=0.0;
832  boundary_vertices[0][1]=-r_max;
833  boundary_vertices[1][0]=r_max;;
834  boundary_vertices[1][1]=-r_max;
835 
836  boundary_id=1;
837  outer_boundary_line_pt[1]=
838  new TriangleMeshPolyLine(boundary_vertices,boundary_id);
839 
840 
841  // Right boundary of bulk mesh
842  //----------------------------
843  boundary_vertices[0][0]=r_max;
844  boundary_vertices[0][1]=-r_max;
845  boundary_vertices[1][0]=r_max;;
846  boundary_vertices[1][1]=r_max;
847 
848  boundary_id=2;
849  outer_boundary_line_pt[2]=
850  new TriangleMeshPolyLine(boundary_vertices,boundary_id);
851 
852 
853  // Top boundary of bulk mesh
854  //--------------------------
855  boundary_vertices[0][0]=r_max;
856  boundary_vertices[0][1]=r_max;
857  boundary_vertices[1][0]=0.0;
858  boundary_vertices[1][1]=r_max;
859 
860  boundary_id=3;
861  outer_boundary_line_pt[3]=
862  new TriangleMeshPolyLine(boundary_vertices,boundary_id);
863 
864 // Top straight boundary on symmetry line
865  //---------------------------------------
866  boundary_vertices[0][0]=0.0;
867  boundary_vertices[0][1]=r_max;
868  boundary_vertices[1][0]=0.0;
869  boundary_vertices[1][1]=r_min;
870 
871  boundary_id=4;
872  outer_boundary_line_pt[4]=
873  new TriangleMeshPolyLine(boundary_vertices,boundary_id);
874 
875 
876  // Inner circular boundary:
877  //-------------------------
878 
879  // Number of segments used for representing the curvilinear boundary
880  unsigned n_segments = 20;
881 
882  // The intrinsic coordinates for the beginning and end of the curve
883  double s_start = 0.5*MathematicalConstants::Pi;
884  double s_end = -0.5*MathematicalConstants::Pi;
885 
886  boundary_id = 5;
887  outer_boundary_line_pt[5]=
888  new TriangleMeshCurviLine(inner_circle_pt,
889  s_start,
890  s_end,
891  n_segments,
892  boundary_id);
893 
894 
895  // Create closed curve that defines outer boundary
896  //------------------------------------------------
897  TriangleMeshClosedCurve *outer_boundary_pt =
898  new TriangleMeshClosedCurve(outer_boundary_line_pt);
899 
900 
901  // Use the TriangleMeshParameters object for helping on the manage of the
902  // TriangleMesh parameters. The only parameter that needs to take is the
903  // outer boundary.
904  TriangleMeshParameters triangle_mesh_parameters(outer_boundary_pt);
905 
906 
907  // Specify maximum element area
908  double element_area = ProblemParameters::Element_area;
909  triangle_mesh_parameters.element_area() = element_area;
910 
911 #ifdef ADAPTIVE
912 
913  // Build "bulk" mesh
914  Bulk_mesh_pt=new RefineableTriangleMesh<ELEMENT>(triangle_mesh_parameters);
915 
916  // Add the bulk mesh to the problem
917  add_sub_mesh(Bulk_mesh_pt);
918 
919  // Initialise mesh as geom object
920  Mesh_as_geom_obj_pt=0;
921 
922  // Create/set error estimator
923  Bulk_mesh_pt->spatial_error_estimator_pt()=new Z2ErrorEstimator;
924 
925  // Choose error tolerances to force some uniform refinement
926  Bulk_mesh_pt->min_permitted_error()=0.00004;
927  Bulk_mesh_pt->max_permitted_error()=0.0001;
928 
929  // Reduce wavenumber to make effect of singularity more prominent
931 
932 #else
933 
934  // Create the bulk mesh
935  Bulk_mesh_pt= new TriangleMesh<ELEMENT>(triangle_mesh_parameters);
936 
937  // Add the bulk mesh to the problem
938  add_sub_mesh(Bulk_mesh_pt);
939 
940  // Create flux elements on inner boundary
941  Helmholtz_inner_boundary_mesh_pt=new Mesh;
942  create_flux_elements_on_inner_boundary();
943 
944  // ...and add the mesh to the problem
945  add_sub_mesh(Helmholtz_inner_boundary_mesh_pt);
946 
947 #endif
948 
949 
950  // Attach the power monitor elements
951  Power_monitor_mesh_pt=new Mesh;
952  create_power_monitor_mesh();
953 
954  // Create the pml meshes
955  create_pml_meshes();
956 
957  // Build the Problem's global mesh from its various sub-meshes
958  build_global_mesh();
959 
960  // Complete the build of all elements
961  complete_problem_setup();
962 
963  // Setup equation numbering scheme
964  cout <<"Number of equations: " << assign_eqn_numbers() << std::endl;
965 
966 } // end of constructor
967 
968 
969 
970 //===============start_of_doc=============================================
971 /// Doc the solution: doc_info contains labels/output directory etc.
972 //========================================================================
973 template<class ELEMENT>
975 doc_solution(DocInfo& doc_info)
976 {
977 
978  ofstream some_file;
979  char filename[100];
980 
981  // Number of plot points: npts x npts
982  unsigned npts=5;
983 
984 
985  // Total radiated power
986  double power=0.0;
987 
988  // Compute/output the radiated power
989  //----------------------------------
990  sprintf(filename,"%s/power%i.dat",doc_info.directory().c_str(),
991  doc_info.number());
992  some_file.open(filename);
993 
994  // Accumulate contribution from elements
995  unsigned nn_element=Power_monitor_mesh_pt->nelement();
996  for(unsigned e=0;e<nn_element;e++)
997  {
998  PMLFourierDecomposedHelmholtzPowerMonitorElement<ELEMENT> *el_pt =
999  dynamic_cast<PMLFourierDecomposedHelmholtzPowerMonitorElement
1000  <ELEMENT>*>(Power_monitor_mesh_pt->element_pt(e));
1001  power += el_pt->global_power_contribution(some_file);
1002  }
1003  some_file.close();
1004  oomph_info << "Total radiated power: " << power << std::endl;
1005 
1006 
1007  // Output solution within the bulk mesh
1008  //-------------------------------------
1009  sprintf(filename,"%s/soln%i.dat",doc_info.directory().c_str(),
1010  doc_info.number());
1011  some_file.open(filename);
1012  Bulk_mesh_pt->output(some_file,npts);
1013  some_file.close();
1014 
1015  // Output solution within pml domains
1016  //-----------------------------------
1017  sprintf(filename,"%s/pml_soln%i.dat",doc_info.directory().c_str(),
1018  doc_info.number());
1019  some_file.open(filename);
1020  PML_top_mesh_pt->output(some_file,npts);
1021  PML_right_mesh_pt->output(some_file,npts);
1022  PML_bottom_mesh_pt->output(some_file,npts);
1023  PML_top_right_mesh_pt->output(some_file,npts);
1024  PML_bottom_right_mesh_pt->output(some_file,npts);
1025  some_file.close();
1026 
1027 
1028  // Output exact solution
1029  //----------------------
1030  sprintf(filename,"%s/exact_soln%i.dat",doc_info.directory().c_str(),
1031  doc_info.number());
1032  some_file.open(filename);
1033  Bulk_mesh_pt->output_fct(some_file,npts,ProblemParameters::get_exact_u);
1034  some_file.close();
1035 
1036 
1037  // Doc error and return of the square of the L2 error
1038  //---------------------------------------------------
1039  double error,norm;
1040  sprintf(filename,"%s/error%i.dat",doc_info.directory().c_str(),
1041  doc_info.number());
1042  some_file.open(filename);
1043  Bulk_mesh_pt->compute_error(some_file,ProblemParameters::get_exact_u,
1044  error,norm);
1045  some_file.close();
1046 
1047  // Doc L2 error and norm of solution
1048  cout << "\nNorm of error : " << sqrt(error) << std::endl;
1049  cout << "Norm of solution: " << sqrt(norm) << std::endl << std::endl;
1050 
1051  sprintf(filename,"%s/int_error%i.dat",doc_info.directory().c_str(),
1052  doc_info.number());
1053  some_file.open(filename);
1054  // Doc L2 error and norm of solution
1055  some_file << "\nNorm of error : " << sqrt(error) << std::endl;
1056  some_file << "Norm of solution: " << sqrt(norm) << std::endl <<std::endl;
1057  some_file << "Relative error: " << sqrt(error)/sqrt(norm) << std::endl;
1058  some_file.close();
1059 
1060  // Write norm and radiated power of solution to trace file
1061  Bulk_mesh_pt->compute_norm(norm);
1062  Trace_file << norm << " " << power << std::endl;
1063 
1064 } // end of doc
1065 
1066 
1067 //============start_of_create_flux_elements=================
1068 /// Create flux elements on inner boundary
1069 //==========================================================
1070 template<class ELEMENT>
1073 {
1074  // Apply flux bc on inner boundary (boundary 5)
1075  unsigned b=5;
1076 
1077 // Loop over the bulk elements adjacent to boundary b
1078  unsigned n_element = Bulk_mesh_pt->nboundary_element(b);
1079  for(unsigned e=0;e<n_element;e++)
1080  {
1081  // Get pointer to the bulk element that is adjacent to boundary b
1082  ELEMENT* bulk_elem_pt = dynamic_cast<ELEMENT*>(
1083  Bulk_mesh_pt->boundary_element_pt(b,e));
1084 
1085  //Find the index of the face of element e along boundary b
1086  int face_index = Bulk_mesh_pt->face_index_at_boundary(b,e);
1087 
1088  // Build the corresponding prescribed incoming-flux element
1089  PMLFourierDecomposedHelmholtzFluxElement<ELEMENT>*
1090  flux_element_pt = new
1091  PMLFourierDecomposedHelmholtzFluxElement<ELEMENT>
1092  (bulk_elem_pt,face_index);
1093 
1094  //Add the prescribed incoming-flux element to the surface mesh
1095  Helmholtz_inner_boundary_mesh_pt->add_element_pt(flux_element_pt);
1096 
1097  // Set the pointer to the prescribed flux function
1098  flux_element_pt->flux_fct_pt() = &ProblemParameters::exact_minus_dudr;
1099 
1100  } //end of loop over bulk elements adjacent to boundary b
1101 
1102 } // end of create flux elements on inner boundary
1103 
1104 
1105 
1106 //============start_of_create_pml_meshes======================================
1107 /// Create PML meshes and add them to the problem's sub-meshes
1108 //============================================================================
1109 template<class ELEMENT>
1111 {
1112  // Bulk mesh bottom boundary id
1113  unsigned int bottom_boundary_id = 1;
1114 
1115  // Bulk mesh right boundary id
1116  unsigned int right_boundary_id = 2;
1117 
1118  // Bulk mesh top boundary id
1119  unsigned int top_boundary_id = 3;
1120 
1121  // PML width in elements for the right layer
1122  unsigned n_x_right_pml = ProblemParameters::Nel_pml;
1123 
1124  // PML width in elements for the top layer
1125  unsigned n_y_top_pml = ProblemParameters::Nel_pml;
1126 
1127  // PML width in elements for the bottom layer
1128  unsigned n_y_bottom_pml = ProblemParameters::Nel_pml;
1129 
1130 
1131  // Outer physical length of the PML layers
1132  // defaults to 4.0
1133  double width_x_right_pml = ProblemParameters::PML_thickness;
1134  double width_y_top_pml = ProblemParameters::PML_thickness;
1135  double width_y_bottom_pml = ProblemParameters::PML_thickness;
1136 
1137  // Build the PML meshes based on the new adapted triangular mesh
1138  PML_right_mesh_pt = TwoDimensionalPMLHelper::create_right_pml_mesh
1139  <PMLLayerElement<ELEMENT> >(Bulk_mesh_pt,
1140  right_boundary_id,
1141  n_x_right_pml,
1142  width_x_right_pml);
1143  PML_top_mesh_pt = TwoDimensionalPMLHelper::create_top_pml_mesh
1144  <PMLLayerElement<ELEMENT> >(Bulk_mesh_pt,
1145  top_boundary_id,
1146  n_y_top_pml,
1147  width_y_top_pml);
1148  PML_bottom_mesh_pt= TwoDimensionalPMLHelper::create_bottom_pml_mesh
1149  <PMLLayerElement<ELEMENT> >(Bulk_mesh_pt,
1150  bottom_boundary_id,
1151  n_y_bottom_pml,
1152  width_y_bottom_pml);
1153 
1154  // Add submeshes to the global mesh
1155  add_sub_mesh(PML_right_mesh_pt);
1156  add_sub_mesh(PML_top_mesh_pt);
1157  add_sub_mesh(PML_bottom_mesh_pt);
1158 
1159  // Rebuild corner PML meshes
1160  PML_top_right_mesh_pt =
1161  TwoDimensionalPMLHelper::create_top_right_pml_mesh
1162  <PMLLayerElement<ELEMENT> >(PML_right_mesh_pt,
1163  PML_top_mesh_pt,
1164  Bulk_mesh_pt,
1165  right_boundary_id);
1166 
1167  PML_bottom_right_mesh_pt =
1168  TwoDimensionalPMLHelper::create_bottom_right_pml_mesh
1169  <PMLLayerElement<ELEMENT> >(PML_right_mesh_pt,
1170  PML_bottom_mesh_pt,
1171  Bulk_mesh_pt,
1172  right_boundary_id);
1173 
1174  // Add submeshes to the global mesh
1175  add_sub_mesh(PML_top_right_mesh_pt);
1176  add_sub_mesh(PML_bottom_right_mesh_pt);
1177 
1178 } // end of create_pml_meshes
1179 
1180 
1181 //===== start_of_main=====================================================
1182 /// Driver code for Pml Fourier decomposed Helmholtz problem
1183 //========================================================================
1184 int main(int argc, char **argv)
1185 {
1186  // Store command line arguments
1187  CommandLineArgs::setup(argc,argv);
1188 
1189  // Define possible command line arguments and parse the ones that
1190  // were actually specified
1191 
1192  // Fourier wavenumber
1193  CommandLineArgs::specify_command_line_flag("--Fourier_wavenumber",
1195 
1196  // Output directory
1197  CommandLineArgs::specify_command_line_flag("--dir",
1199 
1200  // PML thickness
1201  CommandLineArgs::specify_command_line_flag("--pml_thick",
1203 
1204  // Number of elements within PMLs
1205  CommandLineArgs::specify_command_line_flag("--npml_element",
1207 
1208  // Target Element size on first mesh generation
1209  CommandLineArgs::specify_command_line_flag("--element_area",
1211  // k squared (wavenumber squared)
1212  CommandLineArgs::specify_command_line_flag("--k_squared",
1214 
1215  // Validation run?
1216  CommandLineArgs::specify_command_line_flag("--validate");
1217 
1218 
1219  // Demonstrate across a range of omega (or k) values with good mesh for a
1220  // visual test of accuracy (put in by Jonathan Deakin)
1221  CommandLineArgs::specify_command_line_flag("--demonstrate");
1222 
1223  // Parse command line
1224  CommandLineArgs::parse_and_assign();
1225 
1226  // Doc what has actually been specified on the command line
1227  CommandLineArgs::doc_specified_flags();
1228 
1229 
1230  // Create label for output
1231  DocInfo doc_info;
1232 
1233  // Set output directory
1234  doc_info.set_directory(ProblemParameters::Directory);
1235 
1236 #ifdef ADAPTIVE
1237 
1238  // Create the problem with 2D projectable six-node elements from the
1239  // TPMLFourierDecomposedHelmholtzElement family,
1240  // allowing for the imposition of a point source (via another
1241  // templated wrapper)
1243  ProjectablePMLFourierDecomposedHelmholtzElement<
1245  TPMLFourierDecomposedHelmholtzElement<3> > > > problem;
1246 
1247 #else
1248 
1249  // Create the problem with 2D six-node elements from the
1250  // TPMLFourierDecomposedHelmholtzElement family.
1252  <TPMLFourierDecomposedHelmholtzElement<3> >
1253  problem;
1254 
1255 #endif
1256 
1257  // Step number
1258  doc_info.number()=ProblemParameters::N_fourier;
1259 
1260 #ifdef ADAPTIVE
1261 
1262  // Max. number of adaptations
1263  unsigned max_adapt=4;
1264 
1265  // Validation run?
1266  if (CommandLineArgs::command_line_flag_has_been_set("--validate"))
1267  {
1268  max_adapt=1;
1269  }
1270 
1271  // Solve the problem, allowing
1272  // up to max_adapt mesh adaptations after every solve.
1273  problem.newton_solve(max_adapt);
1274 
1275 #else
1276 
1277 
1278 // Solve the problem with Newton's method
1279 problem.newton_solve();
1280 
1281 
1282 
1283 #endif
1284 
1285  //Output the solution
1286  problem.doc_solution(doc_info);
1287 
1288 
1289 } //end of main
////////////////////////////////////////////////////////////////// //////////////////////////////////...
void actions_after_newton_solve()
Update the problem after solve (empty)
Mesh * PML_right_mesh_pt
Pointer to the right PML mesh.
MeshAsGeomObject * Mesh_as_geom_obj_pt
Mesh as geometric object representation of bulk mesh.
void delete_face_elements(Mesh *const &boundary_mesh_pt)
Delete boundary face elements and wipe the surface mesh.
Mesh * PML_top_mesh_pt
Pointer to the top PML mesh.
TriangleMesh< ELEMENT > * Bulk_mesh_pt
Pointer to the "bulk" mesh.
Mesh * Helmholtz_inner_boundary_mesh_pt
Mesh of FaceElements that apply the flux bc on the inner boundary.
void setup_point_source()
Set point source.
void create_power_monitor_mesh()
Create mesh of face elements that monitor the radiated power.
Mesh * PML_bottom_right_mesh_pt
Pointer to the bottom right corner PML mesh.
Mesh * Power_monitor_mesh_pt
Pointer to mesh that stores the power monitor elements.
~PMLFourierDecomposedHelmholtzProblem()
Destructor (empty)
void actions_after_adapt()
Actions after adapt: Rebuild the mesh of prescribed flux elements.
RefineableTriangleMesh< ELEMENT > * Bulk_mesh_pt
Pointer to the refineable "bulk" mesh.
Mesh * PML_bottom_mesh_pt
Pointer to the bottom PML mesh.
void actions_before_adapt()
Actions before adapt: Wipe the mesh of prescribed flux elements.
void create_pml_meshes()
Create PML meshes.
Mesh * PML_top_right_mesh_pt
Pointer to the top right corner PML mesh.
void actions_before_newton_solve()
Update the problem specs before solve (empty)
void doc_solution(DocInfo &doc_info)
Doc the solution. DocInfo object stores flags/labels for where the output gets written to.
void create_flux_elements_on_inner_boundary()
Create flux elements on inner boundary.
Class to impose point source to (wrapped) Helmholtz element.
void fill_in_contribution_to_jacobian(Vector< double > &residuals, DenseMatrix< double > &jacobian)
Add the element's contribution to its residual vector and element Jacobian matrix (wrapper)
void fill_in_point_source_contribution_to_residuals(Vector< double > &residuals)
Add the point source contribution to the residual vector.
std::complex< double > Point_source_magnitude
Magnitude of point source (complex!)
void fill_in_contribution_to_residuals(Vector< double > &residuals)
Add the element's contribution to its residual vector (wrapper)
void setup(const Vector< double > &s_point_source, const std::complex< double > &magnitude)
Set local coordinate and magnitude of point source.
Vector< double > S_point_source
Local coordinates of point at which point source is applied.
PMLLayerElement()
Constructor: Call the constructor for the appropriate Element.
////////////////////////////////////////////////////////////////// //////////////////////////////////...
Vector< double > Coeff(N_terms, 1.0)
Coefficients in the exact solution.
double Z_source
Axial position of point source.
double R_source
Radial position of point source.
unsigned N_terms
Number of terms in the exact solution.
std::complex< double > I(0.0, 1.0)
Imaginary unit.
string Directory
Output directory.
double K_squared
Frequency.
void exact_minus_dudr(const Vector< double > &x, std::complex< double > &flux)
Get -du/dr (spherical r) for exact solution. Equal to prescribed flux on inner boundary.
int N_fourier
The default Fourier wave number.
double Element_area
Target area for initial mesh.
std::complex< double > Magnitude(100.0, 100.0)
Point source magnitude (Complex)
double PML_thickness
Default physical PML thickness.
unsigned Nel_pml
Default number of elements within PMLs.
void get_exact_u(const Vector< double > &x, Vector< double > &u)
Exact solution as a Vector of size 2, containing real and imag parts.
////////////////////////////////////////////////////////////////// //////////////////////////////////...
int main(int argc, char **argv)
Driver code for Pml Fourier decomposed Helmholtz problem.