oomph::Hankel_functions_for_helmholtz_problem Namespace Reference

Functions

void Hankel_first (const unsigned &n, const double &x, Vector< std::complex< double >> &h, Vector< std::complex< double >> &hp)
 
void CHankel_first (const unsigned &n, const std::complex< double > &x, Vector< std::complex< double >> &h, Vector< std::complex< double >> &hp)
 

Detailed Description

//////////////////////////////////////////////////////////////////// //////////////////////////////////////////////////////////////////// //////////////////////////////////////////////////////////////////// //////////////////////////////////////////////////////////////////// Namespace to provide Hankel function of the first kind and various orders – needed for Helmholtz computations.

Function Documentation

◆ CHankel_first()

void oomph::Hankel_functions_for_helmholtz_problem::CHankel_first ( const unsigned n,
const std::complex< double > &  x,
Vector< std::complex< double >> &  h,
Vector< std::complex< double >> &  hp 
)

Compute Hankel function of the first kind of orders 0...n and its derivates at coordinate x. The function returns the vector then its derivative (complex version). This functionality is only required in the computation of the solution for the complex- shifted Laplacian preconditioner.

99  {
100  // Set the highest order actually calculated
101  int n_actual = 0;
102 
103  // Create a vector for the Bessel function of the 1st kind
104  Vector<std::complex<double>> jn(n + 1);
105 
106  // Create a vector for the Bessel function of the 2nd kind
107  Vector<std::complex<double>> yn(n + 1);
108 
109  // Create a vector for the Bessel function (1st kind) derivative
110  Vector<std::complex<double>> jnp(n + 1);
111 
112  // Create a vector for the Bessel function (2nd kind) derivative
113  Vector<std::complex<double>> ynp(n + 1);
114 
115  // Call the (complex) Bessel function to calculate the solution
117  int(n), x, n_actual, &jn[0], &yn[0], &jnp[0], &ynp[0]);
118 
119 #ifdef PARANOID
120  // Tell the user if they tried to calculate higher order terms
121  if (n_actual != int(n))
122  {
123  // Create an output stream
124  std::ostringstream error_message_stream;
125 
126  // Create the error message
127  error_message_stream << "CRBond_Bessel::cbessjyna() only computed "
128  << n_actual << " rather than " << n
129  << " Bessel functions.\n";
130 
131  // Throw the error message
132  throw OomphLibError(error_message_stream.str(),
135  }
136 #endif
137 
138  // Loop over the number of terms requested (only the first entry
139  // has *actually* been calculated)
140  for (unsigned i = 0; i < n; i++)
141  {
142  // Set the entries in the Hankel function vector
143  h[i] = std::complex<double>(jn[i].real() - yn[i].imag(),
144  jn[i].imag() + yn[i].real());
145 
146  // Set the entries in the Hankel function derivative vector
147  hp[i] = std::complex<double>(jnp[i].real() - ynp[i].imag(),
148  jnp[i].imag() + ynp[i].real());
149  }
150  } // End of Hankel_first
AnnoyingScalar imag(const AnnoyingScalar &)
Definition: AnnoyingScalar.h:132
int i
Definition: BiCGSTAB_step_by_step.cpp:9
const unsigned n
Definition: CG3DPackingUnitTest.cpp:11
float real
Definition: datatypes.h:10
int cbessjyna(int n, complex< double > z, int &nm, complex< double > *cj, complex< double > *cy, complex< double > *cjp, complex< double > *cyp)
Definition: crbond_bessel.cc:1858
list x
Definition: plotDoE.py:28
#define OOMPH_EXCEPTION_LOCATION
Definition: oomph_definitions.h:61
#define OOMPH_CURRENT_FUNCTION
Definition: oomph_definitions.h:86

References CRBond_Bessel::cbessjyna(), i, imag(), n, OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION, and plotDoE::x.

◆ Hankel_first()

void oomph::Hankel_functions_for_helmholtz_problem::Hankel_first ( const unsigned n,
const double x,
Vector< std::complex< double >> &  h,
Vector< std::complex< double >> &  hp 
)

Compute Hankel function of the first kind of orders 0...n and its derivates at coordinate x. The function returns the vector then its derivative.

66  {
67  int n_actual = 0;
68  Vector<double> jn(n + 1), yn(n + 1), jnp(n + 1), ynp(n + 1);
70  int(n), x, n_actual, &jn[0], &yn[0], &jnp[0], &ynp[0]);
71 #ifdef PARANOID
72  if (n_actual != int(n))
73  {
74  std::ostringstream error_stream;
75  error_stream << "CRBond_Bessel::bessjyna() only computed " << n_actual
76  << " rather than " << n << " Bessel functions.\n";
77  throw OomphLibError(
79  }
80 #endif
81  for (unsigned i = 0; i < n; i++)
82  {
83  h[i] = std::complex<double>(jn[i], yn[i]);
84  hp[i] = std::complex<double>(jnp[i], ynp[i]);
85  }
86  } // End of Hankel_first
int bessjyna(int n, double x, int &nm, double *jn, double *yn, double *jnp, double *ynp)
Definition: crbond_bessel.cc:852

References CRBond_Bessel::bessjyna(), i, n, OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION, and plotDoE::x.

Referenced by GlobalParameters::get_exact_u(), GlobalParameters::get_exact_u_bessel(), and oomph::HelmholtzDtNMesh< ELEMENT >::setup_gamma().