pml_mapping_functions.h
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// Version 1.0; svn revision $LastChangedRevision: 1097 $
7 //LIC//
8 //LIC// $LastChangedDate: 2015-12-17 11:53:17 +0000 (Thu, 17 Dec 2015) $
9 //LIC//
10 //LIC// Copyright (C) 2006-2016 Matthias Heil and Andrew Hazel
11 //LIC//
12 //LIC// This library is free software; you can redistribute it and/or
13 //LIC// modify it under the terms of the GNU Lesser General Public
14 //LIC// License as published by the Free Software Foundation; either
15 //LIC// version 2.1 of the License, or (at your option) any later version.
16 //LIC//
17 //LIC// This library is distributed in the hope that it will be useful,
18 //LIC// but WITHOUT ANY WARRANTY; without even the implied warranty of
19 //LIC// MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
20 //LIC// Lesser General Public License for more details.
21 //LIC//
22 //LIC// You should have received a copy of the GNU Lesser General Public
23 //LIC// License along with this library; if not, write to the Free Software
24 //LIC// Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA
25 //LIC// 02110-1301 USA.
26 //LIC//
27 //LIC// The authors may be contacted at oomph-lib@maths.man.ac.uk.
28 //LIC//
29 //LIC//====================================================================
30 #ifndef OOMPH_PML_MAPPING_FUNCTIONS_HEADER
31 #define OOMPH_PML_MAPPING_FUNCTIONS_HEADER
32 
33 // Config header generated by autoconfig
34 #ifdef HAVE_CONFIG_H
35  #include <oomph-lib-config.h>
36 #endif
37 
38 #include "complex_matrices.h"
39 #include "oomph_utilities.h"
40 
41 namespace oomph
42 {
43 
44 //=======================================================================
45 /// Class to hold the mapping function (gamma) for the Pml which defines
46 /// how the coordinates are transformed in the Pml. This class holds
47 /// the one dimensional or uniaxial case which is the most common
48 //=======================================================================
50 {
51 
52 public:
53 
54  /// Default constructor (empty)
56 
57  /// \short Pure virtual to return Pml mapping gamma, where gamma is the
58  /// \f$d\tilde x / d x\f$ as function of \f$\nu\f$ where \f$\nu = x - h\f$ where h is
59  /// the vector from the origin to the start of the Pml
60  virtual std::complex<double> dtransformed_nu_dnu(const double& nu,
61  const double& delta,
62  const double& wavenumber_squared,
63  const double& alpha_shift=0.0) = 0;
64 
65  //
66  virtual std::complex<double> transformed_nu(const double& nu,
67  const double& delta,
68  const double& wavenumber_squared,
69  const double& alpha_shift=0.0) = 0;
70 
71 };
72 
73 //=======================================================================
74 /// A mapping function propsed by Bermudez et al, appears to be the best
75 /// for the Helmholtz equations and so this will be the default mapping
76 /// (see definition of PmlHelmholtzEquations)
77 //=======================================================================
79 {
80 
81 public:
82 
83  /// Default constructor (empty)
85 
86  /// \short Overwrite the pure Pml mapping coefficient function to return the
87  /// coeffcients proposed by Bermudez et al
88  std::complex<double> dtransformed_nu_dnu(const double& nu,
89  const double& delta,
90  const double& k,
91  const double& alpha_shift=0.0)
92  {
93  return 1.0 + MathematicalConstants::I / k * (1.0/std::fabs(delta - nu));
94  }
95 
96  //
97  std::complex<double> transformed_nu(const double& nu,
98  const double& delta,
99  const double& k,
100  const double& alpha_shift=0.0)
101  {
102  return nu - MathematicalConstants::I/k * log(1.0 - std::fabs(nu/delta));
103  }
104 
105 };
106 
107 //=======================================================================
108 /// A mapping function proposed by Bermudez et al, similar to the one above
109 /// but is continuous across the inner Pml boundary
110 /// appears to be the best for TimeHarmonicLinearElasticity
111 /// and so this will be the default mapping
112 //=======================================================================
114 {
115 
116 public:
117 
118  /// Default constructor (empty)
120 
121  /// \short Overwrite the pure Pml mapping coefficient function to return the
122  /// coeffcients proposed by Bermudez et al
123  std::complex<double> dtransformed_nu_dnu(const double& nu,
124  const double& delta,
125  const double& k,
126  const double& alpha_shift=0.0)
127  {
128  /// return \f$\gamma=1 + (i/k)(1/|outer_boundary - x|-1/|pml width|)\f$
129  return 1.0 + MathematicalConstants::I / k
130  *( 1.0/std::fabs(delta-nu) - 1.0/std::fabs(delta));
131  }
132 
133  //
134  std::complex<double> transformed_nu(const double& nu,
135  const double& delta,
136  const double& k,
137  const double& alpha_shift=0.0)
138  {
139  return nu - MathematicalConstants::I/k
140  *( log(1.0-std::fabs(nu/delta)) - nu/std::fabs(delta) );
141  }
142 
143 };
144 
145 //=======================================================================
146 /// A slight modification to the mapping function propsed by Bermudez et al.
147 /// The transformation is independent of the scale (thickness) of the PML.
148 /// See Jonathan Deakin's PhD thesis (University of Manchester) for more
149 /// information.
150 //=======================================================================
152 {
153 
154 public:
155 
156  /// Default constructor (empty)
158 
159  /// \short Overwrite the pure Pml mapping coefficient function to return the
160  /// coeffcients proposed by Bermudez et al
161  std::complex<double> dtransformed_nu_dnu(const double& nu,
162  const double& delta,
163  const double& k,
164  const double& alpha_shift=0.0)
165  {
166  return MathematicalConstants::I / k * (1.0/std::fabs(delta - nu));
167  }
168 
169  //
170  std::complex<double> transformed_nu(const double& nu,
171  const double& delta,
172  const double& k,
173  const double& alpha_shift=0.0)
174  {
175  return MathematicalConstants::I/k * log(1.0 - std::fabs(nu/delta));
176  }
177 
178 };
179 
180 //=======================================================================
181 /// Class to hold the mapping function (gamma) for the PML which defines
182 /// how the coordinates are transformed in the PML. This PML mapping
183 /// aims to transform the solution into a straight line, and requires a
184 /// lot of information from the element. Returns a mapping Jacobian which
185 /// I haven't quite worked out yet
186 //=======================================================================
188 {
189 
190 public:
191 
192  /// Default constructor (empty)
194 
195  /// \short Pure virtual to return PML mapping Jacobian
196  virtual void get_mapping_jacobian(
197  const Vector<double>& x,
198  const Vector<double>& x_inner,
199  const Vector<double>& x_outer,
200  const Vector<double>& dx_inner_dacross,
201  const Vector<double>& dp_dacross,
202  const double& k,
203  std::complex<double>& tnu,
204  std::complex<double>& dtnu_dnu,
205  std::complex<double>& dtnu_dacross,
206  const double& alpha=0.0
207  ) = 0;
208 };
209 
210 //=======================================================================
211 /// Class to hold the mapping function (gamma) for the PML which defines
212 /// how the coordinates are transformed in the PML. This PML mapping
213 /// aims to transform the solution into a straight line, and requires a
214 /// lot of information from the element. Returns a mapping Jacobian which
215 /// I haven't quite worked out yet
216 //=======================================================================
219 {
220 
221 public:
222 
223  /// Default constructor (empty)
225 
226  /// \short Search along the line 0 to nun for pole (or closest point to)
227  virtual void pole_line_search(
228  const Vector<double>& x_inner,
229  const Vector<double>& p,
230  const double& k,
231  double& nun,
232  std::complex<double>& tnu,
233  const double& alpha=0.0) = 0;
234 
235  /// Make a Newton step towards a pole (pole is when du/dtnu=0s)
236  virtual bool newton_step_to_pole(
237  const Vector<double>& x_inner,
238  const Vector<double>& p,
239  const Vector<double>& dx_inner_dacross,
240  const Vector<double>& dp_dacross,
241  double& nu,
242  double& zeta,
243  std::complex<double>& tnu,
244  std::complex<double>& alpha,
245  std::complex<double>& beta) = 0;
246 
247  /// Can this be taken out?
248  /// Make a Newton step towards a pole (pole is when du/dtnu=0s)
249  virtual void set_initial_guess(
250  const Vector<double>& x_inner,
251  const Vector<double>& p,
252  const double& nun,
253  const std::complex<double>& tnu) = 0;
254 };
255 
256 }
257 
258 #endif
std::complex< double > transformed_nu(const double &nu, const double &delta, const double &k, const double &alpha_shift=0.0)
std::complex< double > dtransformed_nu_dnu(const double &nu, const double &delta, const double &k, const double &alpha_shift=0.0)
Overwrite the pure Pml mapping coefficient function to return the coeffcients proposed by Bermudez et...
virtual std::complex< double > dtransformed_nu_dnu(const double &nu, const double &delta, const double &wavenumber_squared, const double &alpha_shift=0.0)=0
Pure virtual to return Pml mapping gamma, where gamma is the as function of where where h is the v...
std::complex< double > transformed_nu(const double &nu, const double &delta, const double &k, const double &alpha_shift=0.0)
std::complex< double > dtransformed_nu_dnu(const double &nu, const double &delta, const double &k, const double &alpha_shift=0.0)
Overwrite the pure Pml mapping coefficient function to return the coeffcients proposed by Bermudez et...
virtual std::complex< double > transformed_nu(const double &nu, const double &delta, const double &wavenumber_squared, const double &alpha_shift=0.0)=0
UniaxialPMLMapping()
Default constructor (empty)
std::complex< double > dtransformed_nu_dnu(const double &nu, const double &delta, const double &k, const double &alpha_shift=0.0)
Overwrite the pure Pml mapping coefficient function to return the coeffcients proposed by Bermudez et...
BermudezPMLMapping()
Default constructor (empty)
TangentiallyDiscontinuousConformalPMLMapping()
Default constructor (empty)
ScaleFreeBermudezPMLMapping()
Default constructor (empty)
const std::complex< double > I(0.0, 1.0)
The imaginary unit.
std::complex< double > transformed_nu(const double &nu, const double &delta, const double &k, const double &alpha_shift=0.0)
TangentiallyVaryingConformalPMLMapping()
Default constructor (empty)
C1BermudezPMLMapping()
Default constructor (empty)