unstructured_two_d_helmholtz.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
27 // perfectly matched layer treatment for the exterior boundaries
28 
29 #include<fenv.h>
30 
31 #include "math.h"
32 #include <complex>
33 
34 // Generic routines
35 #include "generic.h"
36 
37 // The pml Helmholtz equations
38 #include "pml_helmholtz.h"
39 
40 // The meshes needed in the PML constructions
41 #include "meshes/triangle_mesh.h"
42 #include "meshes/rectangular_quadmesh.h"
43 
44 using namespace oomph;
45 using namespace std;
46 
47 /// //////////////////////////////////////////////////////////////////
48 /// //////////////////////////////////////////////////////////////////
49 /// //////////////////////////////////////////////////////////////////
50 
51 //===== start_of_namespace=============================================
52 /// Namespace for the Helmholtz problem parameters
53 //=====================================================================
55 {
56 
57  /// Wavenumber (also known as k), k=omega/c
58  double Wavenumber = sqrt(50.0);
59 
60  /// Square of the wavenumber (also known as k^2)
62 
63 } // end of namespace
64 
65 
66 /// //////////////////////////////////////////////////////////////////
67 /// //////////////////////////////////////////////////////////////////
68 /// //////////////////////////////////////////////////////////////////
69 
70 //========= start_of_problem_class=====================================
71 /// Problem class to demonstrate use of perfectly matched layers
72 /// for Helmholtz problems.
73 //=====================================================================
74 template<class ELEMENT>
75 class PMLProblem : public Problem
76 {
77 
78 public:
79 
80  /// Constructor
81  PMLProblem();
82 
83  /// Destructor (empty)
85 
86  /// Update the problem specs before solve (empty)
88 
89  /// Update the problem specs after solve (empty)
91 
92  /// Doc the solution. DocInfo object stores flags/labels for where the
93  /// output gets written to
94  void doc_solution(DocInfo& doc_info);
95 
96  /// Create PML meshes
97  void create_pml_meshes();
98 
99  // Apply boundary conditions
100  void apply_boundary_conditions();
101 
102 #ifdef ADAPTIVE
103 
104  /// Actions before adapt: Wipe the PML meshes
105  void actions_before_adapt();
106 
107  /// Actions after adapt: Rebuild the PML meshes
108  void actions_after_adapt();
109 
110 #endif
111 
112 
113 #ifdef ADAPTIVE
114 
115 private:
116 
117  /// Pointer to the refineable "bulk" mesh
118  RefineableTriangleMesh<ELEMENT>* Bulk_mesh_pt;
119 
120 #else
121 
122 private:
123 
124  /// Pointer to the "bulk" mesh
125  TriangleMesh<ELEMENT>* Bulk_mesh_pt;
126 
127 #endif
128 
129 
130  /// Pointer to the right PML mesh
132 
133  /// Pointer to the top PML mesh
135 
136  /// Pointer to the left PML mesh
138 
139  /// Pointer to the bottom PML mesh
141 
142  /// Pointer to the top right corner PML mesh
144 
145  /// Pointer to the top left corner PML mesh
147 
148  /// Pointer to the bottom right corner PML mesh
150 
151  /// Pointer to the bottom left corner PML mesh
153 
154  /// Trace file
155  ofstream Trace_file;
156 
157 }; // end of problem class
158 
159 
160 
161 //=======start_of_constructor=============================================
162 /// Constructor for Helmholtz problem
163 //========================================================================
164 template<class ELEMENT>
166 {
167 
168  // Open trace file
169  Trace_file.open("RESLT/trace.dat");
170 
171  // Create circle representing inner boundary
172  double a=0.2;
173  double x_c=0.0;
174  double y_c=0.0;
175  Circle* inner_circle_pt=new Circle(x_c,y_c,a);
176 
177  // Outer boundary
178  //---------------
179  TriangleMeshClosedCurve* outer_boundary_pt=0;
180 
181  unsigned n_segments = 20;
182  Vector<TriangleMeshCurveSection*> outer_boundary_line_pt(4);
183 
184  // Each polyline only has three vertices, provide storage for their
185  // coordinates
186  Vector<Vector<double> > vertex_coord(2);
187  for(unsigned i=0;i<2;i++)
188  {
189  vertex_coord[i].resize(2);
190  }
191 
192  // First polyline
193  vertex_coord[0][0]=-2.0;
194  vertex_coord[0][1]=-2.0;
195  vertex_coord[1][0]=-2.0;
196  vertex_coord[1][1]=2.0;
197 
198  // Build the 1st boundary polyline
199  unsigned boundary_id=2;
200  outer_boundary_line_pt[0] = new TriangleMeshPolyLine(vertex_coord,
201  boundary_id);
202 
203  // Second boundary polyline
204  vertex_coord[0][0]=-2.0;
205  vertex_coord[0][1]=2.0;
206  vertex_coord[1][0]=2.0;
207  vertex_coord[1][1]=2.0;
208 
209  // Build the 2nd boundary polyline
210  boundary_id=3;
211  outer_boundary_line_pt[1] = new TriangleMeshPolyLine(vertex_coord,
212  boundary_id);
213 
214  // Third boundary polyline
215  vertex_coord[0][0]=2.0;
216  vertex_coord[0][1]=2.0;
217  vertex_coord[1][0]=2.0;
218  vertex_coord[1][1]=-2.0;
219 
220  // Build the 3rd boundary polyline
221  boundary_id=4;
222  outer_boundary_line_pt[2] = new TriangleMeshPolyLine(vertex_coord,
223  boundary_id);
224 
225  // Fourth boundary polyline
226  vertex_coord[0][0]=2.0;
227  vertex_coord[0][1]=-2.0;
228  vertex_coord[1][0]=-2.0;
229  vertex_coord[1][1]=-2.0;
230 
231  // Build the 4th boundary polyline
232  boundary_id=5;
233  outer_boundary_line_pt[3] = new TriangleMeshPolyLine(vertex_coord,
234  boundary_id);
235 
236  // Create the triangle mesh polygon for outer boundary
237  outer_boundary_pt = new TriangleMeshPolygon(outer_boundary_line_pt);
238 
239  // Inner circular boundary
240  //------------------------
241  Vector<TriangleMeshCurveSection*> inner_boundary_line_pt(2);
242 
243  // The intrinsic coordinates for the beginning and end of the curve
244  double s_start = 0.0;
245  double s_end = MathematicalConstants::Pi;
246  boundary_id = 0;
247  inner_boundary_line_pt[0]=
248  new TriangleMeshCurviLine(inner_circle_pt,
249  s_start,
250  s_end,
251  n_segments,
252  boundary_id);
253 
254  // The intrinsic coordinates for the beginning and end of the curve
255  s_start = MathematicalConstants::Pi;
256  s_end = 2.0*MathematicalConstants::Pi;
257  boundary_id = 1;
258  inner_boundary_line_pt[1]=
259  new TriangleMeshCurviLine(inner_circle_pt,
260  s_start,
261  s_end,
262  n_segments,
263  boundary_id);
264 
265 
266  // Combine to hole
267  //----------------
268  Vector<TriangleMeshClosedCurve*> hole_pt(1);
269  Vector<double> hole_coords(2);
270  hole_coords[0]=0.0;
271  hole_coords[1]=0.0;
272  hole_pt[0]=new TriangleMeshClosedCurve(inner_boundary_line_pt,
273  hole_coords);
274 
275 
276  // Use the TriangleMeshParameters object for helping on the manage
277  // of the TriangleMesh parameters. The only parameter that needs to take
278  // is the outer boundary.
279  TriangleMeshParameters triangle_mesh_parameters(outer_boundary_pt);
280 
281  // Specify the closed curve using the TriangleMeshParameters object
282  triangle_mesh_parameters.internal_closed_curve_pt() = hole_pt;
283 
284  // Target element size in bulk mesh
285  triangle_mesh_parameters.element_area() = 0.1;
286 
287 #ifdef ADAPTIVE
288 
289  // Build adaptive "bulk" mesh
290  Bulk_mesh_pt=new RefineableTriangleMesh<ELEMENT>(triangle_mesh_parameters);
291 
292  // Create/set error estimator
293  Bulk_mesh_pt->spatial_error_estimator_pt()=new Z2ErrorEstimator;
294 
295  // Choose error tolerances to force some uniform refinement
296  Bulk_mesh_pt->min_permitted_error()=0.00004;
297  Bulk_mesh_pt->max_permitted_error()=0.0001;
298 
299 #else
300 
301  // Build "bulk" mesh
302  Bulk_mesh_pt=new TriangleMesh<ELEMENT>(triangle_mesh_parameters);
303 
304 #endif
305 
306  // Create the main triangular mesh
307  add_sub_mesh(Bulk_mesh_pt);
308 
309  // Create PML meshes and add them to the global mesh
310  create_pml_meshes();
311 
312  // Build the entire mesh from its submeshes
313  build_global_mesh();
314 
315  // Let's have a look where the boundaries are
316  this->mesh_pt()->output("global_mesh.dat");
317  this->mesh_pt()->output_boundaries("global_mesh_boundary.dat");
318 
319  // Complete the build of all elements so they are fully functional
320  unsigned n_element = this->mesh_pt()->nelement();
321  for(unsigned e=0;e<n_element;e++)
322  {
323  // Upcast from GeneralisedElement to Helmholtz bulk element
324  PMLHelmholtzEquations<2> *el_pt =
325  dynamic_cast<PMLHelmholtzEquations<2>*>(mesh_pt()->element_pt(e));
326 
327  //Set the k_squared double pointer
328  el_pt->k_squared_pt() = &GlobalParameters::K_squared;
329  }
330 
331  // Apply boundary conditions
332  apply_boundary_conditions();
333 
334  // Setup equation numbering scheme
335  cout <<"Number of equations: " << assign_eqn_numbers() << std::endl;
336 
337 } // end of constructor
338 
339 
340 #ifdef ADAPTIVE
341 
342 //=====================start_of_actions_before_adapt======================
343 /// Actions before adapt: Wipe the mesh of face elements
344 //========================================================================
345 template<class ELEMENT>
347 {
348  // Before adapting the added PML meshes must be removed
349  // as they are not refineable and are to be rebuilt from the
350  // newly refined triangular mesh
351  delete PML_right_mesh_pt;
352  PML_right_mesh_pt=0;
353  delete PML_top_mesh_pt;
354  PML_top_mesh_pt=0;
355  delete PML_left_mesh_pt;
356  PML_left_mesh_pt=0;
357  delete PML_bottom_mesh_pt;
358  PML_bottom_mesh_pt=0;
359  delete PML_top_right_mesh_pt;
360  PML_top_right_mesh_pt=0;
361  delete PML_top_left_mesh_pt;
362  PML_top_left_mesh_pt=0;
363  delete PML_bottom_right_mesh_pt;
364  PML_bottom_right_mesh_pt=0;
365  delete PML_bottom_left_mesh_pt;
366  PML_bottom_left_mesh_pt=0;
367 
368  // Rebuild the Problem's global mesh from its various sub-meshes
369  // but first flush all its submeshes
370  flush_sub_meshes();
371 
372  // Then add the triangular mesh back
373  add_sub_mesh(Bulk_mesh_pt);
374 
375  // Rebuild the global mesh such that it now stores
376  // the triangular mesh only
377  rebuild_global_mesh();
378 
379 }// end of actions_before_adapt
380 
381 
382 //=====================start_of_actions_after_adapt=======================
383 /// Actions after adapt: Rebuild the face element meshes
384 //========================================================================
385 template<class ELEMENT>
387 {
388 
389  // Build PML meshes and add them to the global mesh
390  create_pml_meshes();
391 
392  // Build the entire mesh from its submeshes
393  rebuild_global_mesh();
394 
395  // Complete the build of all elements so they are fully functional
396 
397  // Loop over the entire mesh elements to set up element-specific
398  // things that cannot be handled by constructor
399  unsigned n_element = this->mesh_pt()->nelement();
400 
401  for(unsigned e=0;e<n_element;e++)
402  {
403  // Upcast from GeneralisedElement to PMLHelmholtz bulk element
404  PMLHelmholtzEquations<2> *el_pt =
405  dynamic_cast<PMLHelmholtzEquations<2>*>(mesh_pt()->element_pt(e));
406 
407  //Set the frequency function pointer
408  el_pt->k_squared_pt() = &GlobalParameters::K_squared;
409  }
410 
411  // Re-apply boundary conditions
412  apply_boundary_conditions();
413 
414 }// end of actions_after_adapt
415 
416 #endif
417 
418 //==================start_of_apply_boundary_conditions====================
419 /// Apply boundary conditions
420 //========================================================================
421 template<class ELEMENT>
423 {
424 
425  // Boundary conditions are set on the surface of the circle
426  // as a constant nonzero Dirichlet boundary condition
427  unsigned n_bound = Bulk_mesh_pt->nboundary();
428 
429  for(unsigned b=0;b<n_bound;b++)
430  {
431  unsigned n_node = Bulk_mesh_pt->nboundary_node(b);
432  for (unsigned n=0;n<n_node;n++)
433  {
434  if ((b==0) || (b==1))
435  {
436  Node* nod_pt=Bulk_mesh_pt->boundary_node_pt(b,n);
437  nod_pt->pin(0);
438  nod_pt->pin(1);
439 
440  nod_pt->set_value(0,0.1);
441  nod_pt->set_value(1,0.0);
442  }
443  }
444  }
445 
446 }// end of apply_boundary_conditions
447 
448 
449 //=====================start_of_doc=======================================
450 /// Doc the solution: doc_info contains labels/output directory etc.
451 //========================================================================
452 template<class ELEMENT>
453 void PMLProblem<ELEMENT>::doc_solution(DocInfo& doc_info)
454 {
455 
456  ofstream some_file,some_file2;
457  char filename[100];
458 
459  // Number of plot points
460  unsigned npts;
461  npts=5;
462 
463  // Output solution
464  //-----------------
465  sprintf(filename,"%s/soln%i.dat",doc_info.directory().c_str(),
466  doc_info.number());
467  some_file.open(filename);
468  Bulk_mesh_pt->output(some_file,npts);
469  some_file.close();
470 
471  // Output coarse solution
472  //-----------------------
473  sprintf(filename,"%s/coarse_soln%i.dat",doc_info.directory().c_str(),
474  doc_info.number());
475  some_file.open(filename);
476  unsigned npts_coarse=2;
477  Bulk_mesh_pt->output(some_file,npts_coarse);
478  some_file.close();
479 
480 
481  // Output solution within pml domains
482  //-----------------------------------
483  sprintf(filename,"%s/pml_soln%i.dat",doc_info.directory().c_str(),
484  doc_info.number());
485  some_file.open(filename);
486  PML_top_mesh_pt->output(some_file,npts);
487  PML_right_mesh_pt->output(some_file,npts);
488  PML_bottom_mesh_pt->output(some_file,npts);
489  PML_left_mesh_pt->output(some_file,npts);
490  PML_top_right_mesh_pt->output(some_file,npts);
491  PML_bottom_right_mesh_pt->output(some_file,npts);
492  PML_top_left_mesh_pt->output(some_file,npts);
493  PML_bottom_left_mesh_pt->output(some_file,npts);
494  some_file.close();
495 
496 
497 
498 
499  // Write norm of solution to trace file
500  double norm=0.0;
501  Bulk_mesh_pt->compute_norm(norm);
502  Trace_file << norm << std::endl;
503 
504  // //Do animation of Helmholtz solution
505  // //-----------------------------------
506  // unsigned nstep=40;
507  // for (unsigned i=0;i<nstep;i++)
508  // {
509  // sprintf(filename,"%s/helmholtz_animation%i_frame%i.dat",
510  // doc_info.directory().c_str(),
511  // doc_info.number(),i);
512  // some_file.open(filename);
513  // double phi=2.0*MathematicalConstants::Pi*double(i)/double(nstep-1);
514  // unsigned nelem=Bulk_mesh_pt->nelement();
515  // for (unsigned e=0;e<nelem;e++)
516  // {
517  // ELEMENT* el_pt=dynamic_cast<ELEMENT*>(
518  // Bulk_mesh_pt->element_pt(e));
519  // el_pt->output_real(some_file,phi,npts);
520  // }
521  // some_file.close();
522  // }
523 
524 } // end of doc
525 
526 //============start_of_create_pml_meshes======================================
527 /// Create PML meshes and add them to the problem's sub-meshes
528 //============================================================================
529 template<class ELEMENT>
531 {
532 
533  // Bulk mesh left boundary id
534  unsigned int left_boundary_id = 2;
535 
536  // Bulk mesh top boundary id
537  unsigned int top_boundary_id = 3;
538 
539  // Bulk mesh right boundary id
540  unsigned int right_boundary_id = 4;
541 
542  // Bulk mesh bottom boundary id
543  unsigned int bottom_boundary_id = 5;
544 
545  // PML width in elements for the right layer
546  unsigned n_x_right_pml = 3;
547 
548  // PML width in elements for the top layer
549  unsigned n_y_top_pml = 3;
550 
551  // PML width in elements for the left layer
552  unsigned n_x_left_pml = 3;
553 
554  // PML width in elements for the left layer
555  unsigned n_y_bottom_pml = 3;
556 
557  // Outer physical length of the PML layers
558  // defaults to 0.2, so 10% of the size of the
559  // physical domain
560  double width_x_right_pml = 0.2;
561  double width_y_top_pml = 0.2;
562  double width_x_left_pml = 0.2;
563  double width_y_bottom_pml = 0.2;
564 
565  // Build the PML meshes based on the new adapted triangular mesh
566  PML_right_mesh_pt =
567  TwoDimensionalPMLHelper::create_right_pml_mesh
568  <PMLLayerElement<ELEMENT> >
569  (Bulk_mesh_pt,right_boundary_id,
570  n_x_right_pml, width_x_right_pml);
571  PML_top_mesh_pt =
572  TwoDimensionalPMLHelper::create_top_pml_mesh
573  <PMLLayerElement<ELEMENT> >
574  (Bulk_mesh_pt, top_boundary_id,
575  n_y_top_pml, width_y_top_pml);
576  PML_left_mesh_pt =
577  TwoDimensionalPMLHelper::create_left_pml_mesh
578  <PMLLayerElement<ELEMENT> >
579  (Bulk_mesh_pt, left_boundary_id,
580  n_x_left_pml, width_x_left_pml);
581  PML_bottom_mesh_pt=
582  TwoDimensionalPMLHelper::create_bottom_pml_mesh
583  <PMLLayerElement<ELEMENT> >
584  (Bulk_mesh_pt, bottom_boundary_id,
585  n_y_bottom_pml, width_y_bottom_pml);
586 
587  // Add submeshes to the global mesh
588  add_sub_mesh(PML_right_mesh_pt);
589  add_sub_mesh(PML_top_mesh_pt);
590  add_sub_mesh(PML_left_mesh_pt);
591  add_sub_mesh(PML_bottom_mesh_pt);
592 
593  // Rebuild corner PML meshes
594  PML_top_right_mesh_pt =
595  TwoDimensionalPMLHelper::create_top_right_pml_mesh
596  <PMLLayerElement<ELEMENT> >
597  (PML_right_mesh_pt, PML_top_mesh_pt,
598  Bulk_mesh_pt, right_boundary_id);
599 
600  PML_bottom_right_mesh_pt =
601  TwoDimensionalPMLHelper::create_bottom_right_pml_mesh
602  <PMLLayerElement<ELEMENT> >
603  (PML_right_mesh_pt, PML_bottom_mesh_pt,
604  Bulk_mesh_pt, right_boundary_id);
605 
606  PML_top_left_mesh_pt =
607  TwoDimensionalPMLHelper::create_top_left_pml_mesh
608  <PMLLayerElement<ELEMENT> >
609  (PML_left_mesh_pt, PML_top_mesh_pt,
610  Bulk_mesh_pt, left_boundary_id);
611 
612  PML_bottom_left_mesh_pt =
613  TwoDimensionalPMLHelper::create_bottom_left_pml_mesh
614  <PMLLayerElement<ELEMENT> >
615  (PML_left_mesh_pt, PML_bottom_mesh_pt,
616  Bulk_mesh_pt, left_boundary_id);
617 
618  // Add submeshes to the global mesh
619  add_sub_mesh(PML_top_right_mesh_pt);
620  add_sub_mesh(PML_bottom_right_mesh_pt);
621  add_sub_mesh(PML_top_left_mesh_pt);
622  add_sub_mesh(PML_bottom_left_mesh_pt);
623 
624 } // end of create_pml_meshes
625 
626 
627 
628 //==========start_of_main=================================================
629 /// Solve 2D Helmholtz problem
630 //========================================================================
631 int main(int argc, char **argv)
632 {
633  //Set up the problem
634  //------------------
635 
636 #ifdef ADAPTIVE
637 
638  // Set up the problem with projectable 2D six-node elements from the
639  // TPMLHelmholtzElement family.
640  PMLProblem<ProjectablePMLHelmholtzElement
641  <TPMLHelmholtzElement<2,3> > > problem;
642 
643  // Set up the problem with 2D ten-node elements from the
644  // TPMLHelmholtzElement family.
645  // PMLProblem<ProjectablePMLHelmholtzElement
646  // <TPMLHelmholtzElement<2,4> > > problem;
647 
648 #else
649 
650  // Set up the problem with 2D six-node elements from the
651  // TPMLHelmholtzElement family.
653 
654  // Set up the problem with 2D ten-node elements from the
655  // TPMLHelmholtzElement family.
656  // PMLProblem<TPMLHelmholtzElement<2,4> > problem;
657 
658 #endif
659 
660  // Create label for output
661  //------------------------
662  DocInfo doc_info;
663 
664  // Set output directory
665  doc_info.set_directory("RESLT");
666 
667 
668 #ifdef ADAPTIVE
669 
670  // Max. number of adaptations
671  unsigned max_adapt=1;
672 
673  // Solve the problem with the adaptive Newton's method, allowing
674  // up to max_adapt mesh adaptations after every solve.
675  problem.newton_solve(max_adapt);
676 
677 #else
678 
679  // Solve the problem with Newton's method
680  problem.newton_solve();
681 
682 #endif
683 
684  //Output solution
685  problem.doc_solution(doc_info);
686 
687 } //end of main
////////////////////////////////////////////////////////////////// //////////////////////////////////...
Mesh * PML_bottom_left_mesh_pt
Pointer to the bottom left corner PML mesh.
void actions_before_newton_solve()
Update the problem specs before solve (empty)
void apply_boundary_conditions()
Apply boundary conditions.
void actions_before_adapt()
Actions before adapt: Wipe the PML meshes.
~PMLProblem()
Destructor (empty)
Mesh * PML_bottom_right_mesh_pt
Pointer to the bottom right corner PML mesh.
Mesh * PML_bottom_mesh_pt
Pointer to the bottom PML mesh.
ofstream Trace_file
Trace file.
TriangleMesh< ELEMENT > * Bulk_mesh_pt
Pointer to the "bulk" mesh.
void actions_after_adapt()
Actions after adapt: Rebuild the PML meshes.
Mesh * PML_top_mesh_pt
Pointer to the top PML mesh.
Mesh * PML_top_right_mesh_pt
Pointer to the top right corner PML mesh.
void create_pml_meshes()
Create PML meshes.
Mesh * PML_right_mesh_pt
Pointer to the right PML mesh.
Mesh * PML_top_left_mesh_pt
Pointer to the top left corner PML mesh.
Mesh * PML_left_mesh_pt
Pointer to the left PML mesh.
void actions_after_newton_solve()
Update the problem specs after solve (empty)
RefineableTriangleMesh< ELEMENT > * Bulk_mesh_pt
Pointer to the refineable "bulk" mesh.
void doc_solution(DocInfo &doc_info)
Doc the solution. DocInfo object stores flags/labels for where the output gets written to.
////////////////////////////////////////////////////////////////// //////////////////////////////////...
double Wavenumber
Wavenumber (also known as k), k=omega/c.
double K_squared
Square of the wavenumber (also known as k^2)
int main(int argc, char **argv)
Solve 2D Helmholtz problem.