pml_fourier_decomposed_helmholtz_elements.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// Copyright (C) 2006-2022 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 // Header file for Fourier-decomposed Helmholtz elements
27 #ifndef OOMPH_PML_FOURIER_DECOMPOSED_HELMHOLTZ_ELEMENTS_HEADER
28 #define OOMPH_PML_FOURIER_DECOMPOSED_HELMHOLTZ_ELEMENTS_HEADER
29 
30 
31 // Config header generated by autoconfig
32 #ifdef HAVE_CONFIG_H
33 #include <oomph-lib-config.h>
34 #endif
35 
36 #include "math.h"
37 #include <complex>
38 
39 
40 // OOMPH-LIB headers
41 #include "../generic/projection.h"
42 #include "../generic/nodes.h"
43 #include "../generic/Qelements.h"
44 #include "../generic/oomph_utilities.h"
45 #include "../generic/pml_meshes.h"
46 #include "../generic/projection.h"
47 #include "../generic/oomph_definitions.h"
48 
49 
50 namespace oomph
51 {
52  //========================================================================
54  //========================================================================
55  namespace Legendre_functions_helper
56  {
58  extern double factorial(const unsigned& l);
59 
61  extern double plgndr1(const unsigned& n, const double& x);
62 
64  extern double plgndr2(const unsigned& l,
65  const unsigned& m,
66  const double& x);
67 
68  } // namespace Legendre_functions_helper
69 
70 
74 
75  //=======================================================================
78  //=======================================================================
80  {
81  public:
84 
88  virtual std::complex<double> gamma(const double& nu_i,
89  const double& pml_width_i,
90  const double& k_squared) = 0;
91 
95  virtual std::complex<double> transformed_coordinate(
96  const double& nu_i,
97  const double& pml_width_i,
98  const double& pml_inner_boundary,
99  const double& k_squared) = 0;
100  };
101 
102 
103  //=======================================================================
107  //=======================================================================
110  {
111  public:
114 
117  std::complex<double> gamma(const double& nu_i,
118  const double& pml_width_i,
119  const double& k_squared)
120  {
123  return 1.0 +
124  std::complex<double>(1.0 / sqrt(k_squared), 0) *
125  std::complex<double>(0.0, 1.0 / (std::fabs(pml_width_i - nu_i)));
126  }
127 
130  std::complex<double> transformed_coordinate(
131  const double& nu_i,
132  const double& pml_width_i,
133  const double& pml_inner_boundary,
134  const double& k_squared)
135  {
137  double log_arg = 1.0 - std::fabs(nu_i / pml_width_i);
138  return std::complex<double>(pml_inner_boundary + nu_i,
139  -log(log_arg) / sqrt(k_squared));
140  }
141  };
142 
143 
147 
148 
149  //=============================================================
158  //=============================================================
160  : public virtual PMLElementBase<2>,
161  public virtual FiniteElement
162  {
163  public:
167  const Vector<double>& x, std::complex<double>& f);
168 
172  {
173  // Provide pointer to static method (Save memory)
177 
179  }
180 
181 
184  const PMLFourierDecomposedHelmholtzEquations& dummy) = delete;
185 
187  // Commented out broken assignment operator because this can lead to a
188  // conflict warning when used in the virtual inheritence hierarchy.
189  // Essentially the compiler doesn't realise that two separate
190  // implementations of the broken function are the same and so, quite
191  // rightly, it shouts.
192  /*void operator=(const PMLFourierDecomposedHelmholtzEquations&) = delete;*/
193 
197  virtual inline std::complex<unsigned> u_index_pml_fourier_decomposed_helmholtz()
198  const
199  {
200  return std::complex<unsigned>(0, 1);
201  }
202 
203 
205  double*& k_squared_pt()
206  {
207  return K_squared_pt;
208  }
209 
210 
212  double k_squared()
213  {
214 #ifdef PARANOID
215  if (K_squared_pt == 0)
216  {
217  throw OomphLibError(
218  "Please set pointer to k_squared using access fct to pointer!",
221  }
222 #endif
223  return *K_squared_pt;
224  }
225 
227  double*& alpha_pt()
228  {
229  return Alpha_pt;
230  }
231 
232 
234  double alpha()
235  {
236  return *Alpha_pt;
237  }
238 
241  {
242  return N_pml_fourier_pt;
243  }
244 
247  {
248  if (N_pml_fourier_pt == 0)
249  {
250  return 0;
251  }
252  else
253  {
254  return *N_pml_fourier_pt;
255  }
256  }
257 
258 
260  void output(std::ostream& outfile)
261  {
262  const unsigned n_plot = 5;
263  output(outfile, n_plot);
264  }
265 
268  void output(std::ostream& outfile, const unsigned& n_plot);
269 
275  void output_real(std::ostream& outfile,
276  const double& phi,
277  const unsigned& n_plot);
278 
280  void output(FILE* file_pt)
281  {
282  const unsigned n_plot = 5;
283  output(file_pt, n_plot);
284  }
285 
288  void output(FILE* file_pt, const unsigned& n_plot);
289 
292  void output_fct(std::ostream& outfile,
293  const unsigned& n_plot,
295 
298  virtual void output_fct(
299  std::ostream& outfile,
300  const unsigned& n_plot,
301  const double& time,
303  {
304  throw OomphLibError("There is no time-dependent output_fct() for "
305  "PMLFourierDecomposedHelmholtz elements ",
308  }
309 
310 
316  void output_real_fct(std::ostream& outfile,
317  const double& phi,
318  const unsigned& n_plot,
320 
321 
323  void compute_error(std::ostream& outfile,
325  double& error,
326  double& norm);
327 
328 
330  void compute_error(std::ostream& outfile,
332  const double& time,
333  double& error,
334  double& norm)
335  {
336  throw OomphLibError("There is no time-dependent compute_error() for "
337  "PMLFourierDecomposedHelmholtz elements",
340  }
341 
343  void compute_norm(double& norm);
344 
345 
348  {
349  return Source_fct_pt;
350  }
351 
352 
355  {
356  return Source_fct_pt;
357  }
358 
359 
365  const unsigned& ipt,
366  const Vector<double>& x,
367  std::complex<double>& source) const
368  {
369  // If no source function has been set, return zero
370  if (Source_fct_pt == 0)
371  {
372  source = std::complex<double>(0.0, 0.0);
373  }
374  else
375  {
376  // Get source strength
377  (*Source_fct_pt)(x, source);
378  }
379  }
380 
385  Vector<unsigned>& values_to_pin)
386  {
387  values_to_pin.resize(2);
388  for (unsigned j = 0; j < 2; j++)
389  {
390  values_to_pin[j] = j;
391  }
392  }
393 
394 
396  void get_flux(const Vector<double>& s,
397  Vector<std::complex<double>>& flux) const
398  {
399  // Find out how many nodes there are in the element
400  const unsigned n_node = nnode();
401 
402  // Set up memory for the shape and test functions
403  Shape psi(n_node);
404  DShape dpsidx(n_node, 2);
405 
406  // Call the derivatives of the shape and test functions
407  dshape_eulerian(s, psi, dpsidx);
408 
409  // Initialise to zero
410  const std::complex<double> zero(0.0, 0.0);
411  for (unsigned j = 0; j < 2; j++)
412  {
413  flux[j] = zero;
414  }
415 
416  // Loop over nodes
417  for (unsigned l = 0; l < n_node; l++)
418  {
419  // Cache the complex value of the unknown
420  const std::complex<double> u_value(
421  this->nodal_value(l,
423  this->nodal_value(l,
425 
426  // Loop over derivative directions
427  for (unsigned j = 0; j < 2; j++)
428  {
429  flux[j] += u_value * dpsidx(l, j);
430  }
431  }
432  }
433 
434 
437  {
438  // Call the generic residuals function with flag set to 0
439  // using a dummy matrix argument
441  residuals, GeneralisedElement::Dummy_matrix, 0);
442  }
443 
444 
448  DenseMatrix<double>& jacobian)
449  {
450  // Call the generic routine with the flag set to 1
452  residuals, jacobian, 1);
453  }
454 
455 
459  const Vector<double>& s) const
460  {
461  // Find number of nodes
462  const unsigned n_node = nnode();
463 
464  // Local shape function
465  Shape psi(n_node);
466 
467  // Find values of shape function
468  shape(s, psi);
469 
470  // Initialise value of u
471  std::complex<double> interpolated_u(0.0, 0.0);
472 
473  // Get the index at which the helmholtz unknown is stored
474  const unsigned u_nodal_index_real =
476  const unsigned u_nodal_index_imag =
478 
479  // Loop over the local nodes and sum
480  for (unsigned l = 0; l < n_node; l++)
481  {
482  // Make a temporary complex number from the stored data
483  const std::complex<double> u_value(
484  this->nodal_value(l, u_nodal_index_real),
485  this->nodal_value(l, u_nodal_index_imag));
486  // Add to the interpolated value
487  interpolated_u += u_value * psi[l];
488  }
489  return interpolated_u;
490  }
491 
492 
494  unsigned self_test();
495 
496 
497  protected:
503  const unsigned& ipt,
504  const Vector<double>& x,
505  Vector<std::complex<double>>& pml_laplace_factor,
506  std::complex<double>& pml_k_squared_factor)
507  {
509  Vector<double> nu(2);
510  for (unsigned k = 0; k < 2; k++)
511  {
512  nu[k] = x[k] - this->Pml_inner_boundary[k];
513  }
514 
517  Vector<double> pml_width(2);
518  for (unsigned k = 0; k < 2; k++)
519  {
520  pml_width[k] =
521  this->Pml_outer_boundary[k] - this->Pml_inner_boundary[k];
522  }
523 
524  // Declare gamma_i vectors of complex numbers for PML weights
525  Vector<std::complex<double>> pml_gamma(2);
526 
527  if (this->Pml_is_enabled)
528  {
529  // Cache k_squared to pass into mapping function
530  double k_squared_local = k_squared();
531 
532  for (unsigned k = 0; k < 2; k++)
533  {
534  // If PML is enabled in the respective direction
535  if (this->Pml_direction_active[k])
536  {
538  nu[k], pml_width[k], k_squared_local);
539  }
540  else
541  {
542  pml_gamma[k] = 1.0;
543  }
544  }
545 
550  pml_laplace_factor[0] = pml_gamma[1] / pml_gamma[0];
551  pml_laplace_factor[1] = pml_gamma[0] / pml_gamma[1];
552  pml_k_squared_factor = pml_gamma[0] * pml_gamma[1];
553  }
554  else
555  {
558  for (unsigned k = 0; k < 2; k++)
559  {
560  pml_laplace_factor[k] = std::complex<double>(1.0, 0.0);
561  }
562 
563  pml_k_squared_factor = std::complex<double>(1.0, 0.0);
564  }
565  }
566 
567 
570  void compute_complex_r(const unsigned& ipt,
571  const Vector<double>& x,
572  std::complex<double>& complex_r)
573  {
574  // Cache current position r
575  double r = x[0];
576 
582 
583  // If the complex r variable is imaginary
584  if (this->Pml_is_enabled && (this->Pml_direction_active[0]))
585  {
586  double nu = x[0] - Pml_inner_boundary[0];
587  double pml_width = Pml_outer_boundary[0] - Pml_inner_boundary[0];
588  double k_squared_local = k_squared();
589 
590  // Determine the complex r variable
591  complex_r =
593  nu, pml_width, Pml_inner_boundary[0], k_squared_local);
594  }
595  else
596  {
597  // The complex r variable is infact purely real, and
598  // is equal to x[0]
599  complex_r = std::complex<double>(r, 0.0);
600  }
601 
602  } // end of compute_complex_r
603 
606  {
608  }
609 
612  const
613  {
615  }
616 
621 
622 
626  const Vector<double>& s,
627  Shape& psi,
628  DShape& dpsidx,
629  Shape& test,
630  DShape& dtestdx) const = 0;
631 
632 
636  const unsigned& ipt,
637  Shape& psi,
638  DShape& dpsidx,
639  Shape& test,
640  DShape& dtestdx) const = 0;
641 
645  Vector<double>& residuals,
646  DenseMatrix<double>& jacobian,
647  const unsigned& flag);
648 
651 
653  double* K_squared_pt;
654 
659 
661  double* Alpha_pt;
662 
665 
668  };
669 
670 
674 
675 
676  //======================================================================
680  //======================================================================
681  template<unsigned NNODE_1D>
683  : public virtual QElement<2, NNODE_1D>,
685  {
686  private:
689  static const unsigned Initial_Nvalue;
690 
691  public:
696  {
697  }
698 
702 
704  /*void operator=(const
705  QPMLFourierDecomposedHelmholtzElement<NNODE_1D>&) = delete;*/
706 
707 
710  inline unsigned required_nvalue(const unsigned& n) const
711  {
712  return Initial_Nvalue;
713  }
714 
716  void output(std::ostream& outfile)
717  {
719  }
720 
723  void output(std::ostream& outfile, const unsigned& n_plot)
724  {
726  }
727 
733  void output_real(std::ostream& outfile,
734  const double& phi,
735  const unsigned& n_plot)
736  {
738  }
739 
741  void output(FILE* file_pt)
742  {
744  }
745 
748  void output(FILE* file_pt, const unsigned& n_plot)
749  {
751  }
752 
755  void output_fct(std::ostream& outfile,
756  const unsigned& n_plot,
758  {
760  outfile, n_plot, exact_soln_pt);
761  }
762 
768  void output_real_fct(std::ostream& outfile,
769  const double& phi,
770  const unsigned& n_plot,
772  {
774  outfile, phi, n_plot, exact_soln_pt);
775  }
776 
777 
781  void output_fct(std::ostream& outfile,
782  const unsigned& n_plot,
783  const double& time,
785  {
787  outfile, n_plot, time, exact_soln_pt);
788  }
789 
790  protected:
794  const Vector<double>& s,
795  Shape& psi,
796  DShape& dpsidx,
797  Shape& test,
798  DShape& dtestdx) const;
799 
800 
804  const unsigned& ipt,
805  Shape& psi,
806  DShape& dpsidx,
807  Shape& test,
808  DShape& dtestdx) const;
809  };
810 
811 
812  // Inline functions:
813 
814 
815  //======================================================================
820  //======================================================================
821  template<unsigned NNODE_1D>
824  const Vector<double>& s,
825  Shape& psi,
826  DShape& dpsidx,
827  Shape& test,
828  DShape& dtestdx) const
829  {
830  // Call the geometrical shape functions and derivatives
831  const double J = this->dshape_eulerian(s, psi, dpsidx);
832 
833  // Set the test functions equal to the shape functions
834  test = psi;
835  dtestdx = dpsidx;
836 
837  // Return the jacobian
838  return J;
839  }
840 
841 
842  //======================================================================
847  //======================================================================
848  template<unsigned NNODE_1D>
851  const unsigned& ipt,
852  Shape& psi,
853  DShape& dpsidx,
854  Shape& test,
855  DShape& dtestdx) const
856  {
857  // Call the geometrical shape functions and derivatives
858  const double J = this->dshape_eulerian_at_knot(ipt, psi, dpsidx);
859 
860  // Set the pointers of the test functions
861  test = psi;
862  dtestdx = dpsidx;
863 
864  // Return the jacobian
865  return J;
866  }
867 
871 
872 
873  //=======================================================================
879  //=======================================================================
880  template<unsigned NNODE_1D>
882  : public virtual QElement<1, NNODE_1D>
883  {
884  public:
887  FaceGeometry() : QElement<1, NNODE_1D>() {}
888  };
889 
890 
894 
895 
896  //==========================================================
898  //==========================================================
899  template<class FOURIER_DECOMPOSED_HELMHOLTZ_ELEMENT>
901  : public virtual ProjectableElement<FOURIER_DECOMPOSED_HELMHOLTZ_ELEMENT>
902  {
903  public:
907 
912  {
913 #ifdef PARANOID
914  if (fld > 1)
915  {
916  std::stringstream error_stream;
917  error_stream << "Fourier decomposed Helmholtz elements only store 2 "
918  "fields so fld = "
919  << fld << " is illegal \n";
920  throw OomphLibError(
921  error_stream.str(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
922  }
923 #endif
924 
925  // Create the vector
926  unsigned nnod = this->nnode();
927  Vector<std::pair<Data*, unsigned>> data_values(nnod);
928 
929  // Loop over all nodes
930  for (unsigned j = 0; j < nnod; j++)
931  {
932  // Add the data value associated field: The node itself
933  data_values[j] = std::make_pair(this->node_pt(j), fld);
934  }
935 
936  // Return the vector
937  return data_values;
938  }
939 
942  {
943  return 2;
944  }
945 
948  unsigned nhistory_values_for_projection(const unsigned& fld)
949  {
950 #ifdef PARANOID
951  if (fld > 1)
952  {
953  std::stringstream error_stream;
954  error_stream << "Helmholtz elements only store two fields so fld = "
955  << fld << " is illegal\n";
956  throw OomphLibError(
957  error_stream.str(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
958  }
959 #endif
960  return this->node_pt(0)->ntstorage();
961  }
962 
966  {
967  return this->node_pt(0)->position_time_stepper_pt()->ntstorage();
968  }
969 
972  double jacobian_and_shape_of_field(const unsigned& fld,
973  const Vector<double>& s,
974  Shape& psi)
975  {
976 #ifdef PARANOID
977  if (fld > 1)
978  {
979  std::stringstream error_stream;
980  error_stream << "Helmholtz elements only store two fields so fld = "
981  << fld << " is illegal.\n";
982  throw OomphLibError(
983  error_stream.str(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
984  }
985 #endif
986  unsigned n_dim = this->dim();
987  unsigned n_node = this->nnode();
988  Shape test(n_node);
989  DShape dpsidx(n_node, n_dim), dtestdx(n_node, n_dim);
990  double J =
991  this->dshape_and_dtest_eulerian_pml_fourier_decomposed_helmholtz(
992  s, psi, dpsidx, test, dtestdx);
993  return J;
994  }
995 
996 
999  double get_field(const unsigned& t,
1000  const unsigned& fld,
1001  const Vector<double>& s)
1002  {
1003 #ifdef PARANOID
1004  if (fld > 1)
1005  {
1006  std::stringstream error_stream;
1007  error_stream << "Helmholtz elements only store two fields so fld = "
1008  << fld << " is illegal\n";
1009  throw OomphLibError(
1010  error_stream.str(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
1011  }
1012 #endif
1013  // Find the index at which the variable is stored
1014  std::complex<unsigned> complex_u_nodal_index =
1015  this->u_index_pml_fourier_decomposed_helmholtz();
1016  unsigned u_nodal_index = 0;
1017  if (fld == 0)
1018  {
1019  u_nodal_index = complex_u_nodal_index.real();
1020  }
1021  else
1022  {
1023  u_nodal_index = complex_u_nodal_index.imag();
1024  }
1025 
1026 
1027  // Local shape function
1028  unsigned n_node = this->nnode();
1029  Shape psi(n_node);
1030 
1031  // Find values of shape function
1032  this->shape(s, psi);
1033 
1034  // Initialise value of u
1035  double interpolated_u = 0.0;
1036 
1037  // Sum over the local nodes
1038  for (unsigned l = 0; l < n_node; l++)
1039  {
1040  interpolated_u += this->nodal_value(t, l, u_nodal_index) * psi[l];
1041  }
1042  return interpolated_u;
1043  }
1044 
1045 
1047  unsigned nvalue_of_field(const unsigned& fld)
1048  {
1049 #ifdef PARANOID
1050  if (fld > 1)
1051  {
1052  std::stringstream error_stream;
1053  error_stream << "Helmholtz elements only store two fields so fld = "
1054  << fld << " is illegal\n";
1055  throw OomphLibError(
1056  error_stream.str(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
1057  }
1058 #endif
1059  return this->nnode();
1060  }
1061 
1062 
1064  int local_equation(const unsigned& fld, const unsigned& j)
1065  {
1066 #ifdef PARANOID
1067  if (fld > 1)
1068  {
1069  std::stringstream error_stream;
1070  error_stream << "Helmholtz elements only store two fields so fld = "
1071  << fld << " is illegal\n";
1072  throw OomphLibError(
1073  error_stream.str(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
1074  }
1075 #endif
1076  std::complex<unsigned> complex_u_nodal_index =
1077  this->u_index_pml_fourier_decomposed_helmholtz();
1078  unsigned u_nodal_index = 0;
1079  if (fld == 0)
1080  {
1081  u_nodal_index = complex_u_nodal_index.real();
1082  }
1083  else
1084  {
1085  u_nodal_index = complex_u_nodal_index.imag();
1086  }
1087  return this->nodal_local_eqn(j, u_nodal_index);
1088  }
1089 
1090 
1093  void output(std::ostream& outfile, const unsigned& nplot)
1094  {
1096  }
1097  };
1098 
1099 
1100  //=======================================================================
1103  //=======================================================================
1104  template<class ELEMENT>
1106  : public virtual FaceGeometry<ELEMENT>
1107  {
1108  public:
1109  FaceGeometry() : FaceGeometry<ELEMENT>() {}
1110  };
1111 
1112 
1113  //=======================================================================
1116  //=======================================================================
1117  template<class ELEMENT>
1120  : public virtual FaceGeometry<FaceGeometry<ELEMENT>>
1121  {
1122  public:
1124  };
1125 
1126 
1130 
1131  //=======================================================================
1134  //=======================================================================
1135  template<unsigned NNODE_1D>
1137  : public virtual QPMLFourierDecomposedHelmholtzElement<NNODE_1D>
1138  {
1139  public:
1143  };
1144 
1145 } // namespace oomph
1146 
1147 #endif
AnnoyingScalar imag(const AnnoyingScalar &)
Definition: AnnoyingScalar.h:132
AnnoyingScalar sqrt(const AnnoyingScalar &x)
Definition: AnnoyingScalar.h:134
const unsigned n
Definition: CG3DPackingUnitTest.cpp:11
JacobiRotation< float > J
Definition: Jacobi_makeJacobi.cpp:3
Definition: pml_fourier_decomposed_helmholtz_elements.h:110
std::complex< double > gamma(const double &nu_i, const double &pml_width_i, const double &k_squared)
Definition: pml_fourier_decomposed_helmholtz_elements.h:117
std::complex< double > transformed_coordinate(const double &nu_i, const double &pml_width_i, const double &pml_inner_boundary, const double &k_squared)
Definition: pml_fourier_decomposed_helmholtz_elements.h:130
BermudezPMLMappingAndTransformedCoordinate()
Default constructor (empty)
Definition: pml_fourier_decomposed_helmholtz_elements.h:113
Definition: shape.h:278
unsigned ntstorage() const
Definition: nodes.cc:879
FaceGeometry()
Definition: pml_fourier_decomposed_helmholtz_elements.h:1123
FaceGeometry()
Definition: pml_fourier_decomposed_helmholtz_elements.h:1109
FaceGeometry()
Definition: pml_fourier_decomposed_helmholtz_elements.h:887
Definition: elements.h:4998
Definition: elements.h:1313
Node *& node_pt(const unsigned &n)
Return a pointer to the local node n.
Definition: elements.h:2175
double nodal_value(const unsigned &n, const unsigned &i) const
Definition: elements.h:2593
virtual void shape(const Vector< double > &s, Shape &psi) const =0
int nodal_local_eqn(const unsigned &n, const unsigned &i) const
Definition: elements.h:1432
unsigned dim() const
Definition: elements.h:2611
unsigned nnode() const
Return the number of nodes.
Definition: elements.h:2210
void(* SteadyExactSolutionFctPt)(const Vector< double > &, Vector< double > &)
Definition: elements.h:1759
void(* UnsteadyExactSolutionFctPt)(const double &, const Vector< double > &, Vector< double > &)
Definition: elements.h:1765
double dshape_eulerian(const Vector< double > &s, Shape &psi, DShape &dpsidx) const
Definition: elements.cc:3298
static DenseMatrix< double > Dummy_matrix
Definition: elements.h:227
TimeStepper *& position_time_stepper_pt()
Return a pointer to the position timestepper.
Definition: nodes.h:1022
Definition: oomph_definitions.h:222
Base class for elements with pml capabilities.
Definition: pml_meshes.h:60
bool Pml_is_enabled
Boolean indicating if element is used in pml mode.
Definition: pml_meshes.h:119
std::vector< bool > Pml_direction_active
Definition: pml_meshes.h:124
Vector< double > Pml_outer_boundary
Definition: pml_meshes.h:134
Vector< double > Pml_inner_boundary
Definition: pml_meshes.h:129
Definition: pml_fourier_decomposed_helmholtz_elements.h:162
std::complex< double > interpolated_u_pml_fourier_decomposed_helmholtz(const Vector< double > &s) const
Definition: pml_fourier_decomposed_helmholtz_elements.h:458
int pml_fourier_wavenumber()
Get the Fourier wavenumber.
Definition: pml_fourier_decomposed_helmholtz_elements.h:246
virtual void fill_in_generic_residual_contribution_pml_fourier_decomposed_helmholtz(Vector< double > &residuals, DenseMatrix< double > &jacobian, const unsigned &flag)
Definition: pml_fourier_decomposed_helmholtz_elements.cc:192
PMLMappingAndTransformedCoordinate *const & pml_mapping_and_transformed_coordinate_pt() const
Return a pointer to the PML Mapping object (const version)
Definition: pml_fourier_decomposed_helmholtz_elements.h:611
void compute_error(std::ostream &outfile, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt, double &error, double &norm)
Get error against and norm of exact solution.
Definition: pml_fourier_decomposed_helmholtz_elements.cc:715
double * K_squared_pt
Pointer to k^2 (wavenumber squared)
Definition: pml_fourier_decomposed_helmholtz_elements.h:653
void output_real(std::ostream &outfile, const double &phi, const unsigned &n_plot)
Definition: pml_fourier_decomposed_helmholtz_elements.cc:529
double *& alpha_pt()
Get pointer to complex shift.
Definition: pml_fourier_decomposed_helmholtz_elements.h:227
void(* PMLFourierDecomposedHelmholtzSourceFctPt)(const Vector< double > &x, std::complex< double > &f)
Definition: pml_fourier_decomposed_helmholtz_elements.h:166
void output(std::ostream &outfile)
Output with default number of plot points.
Definition: pml_fourier_decomposed_helmholtz_elements.h:260
virtual double dshape_and_dtest_eulerian_pml_fourier_decomposed_helmholtz(const Vector< double > &s, Shape &psi, DShape &dpsidx, Shape &test, DShape &dtestdx) const =0
void compute_complex_r(const unsigned &ipt, const Vector< double > &x, std::complex< double > &complex_r)
Definition: pml_fourier_decomposed_helmholtz_elements.h:570
static double Default_Physical_Constant_Value
Static default value for the physical constants (initialised to zero)
Definition: pml_fourier_decomposed_helmholtz_elements.h:664
void compute_error(std::ostream &outfile, FiniteElement::UnsteadyExactSolutionFctPt exact_soln_pt, const double &time, double &error, double &norm)
Dummy, time dependent error checker.
Definition: pml_fourier_decomposed_helmholtz_elements.h:330
double * Alpha_pt
Pointer to wavenumber complex shift.
Definition: pml_fourier_decomposed_helmholtz_elements.h:661
unsigned self_test()
Self-test: Return 0 for OK.
Definition: pml_fourier_decomposed_helmholtz_elements.cc:460
void fill_in_contribution_to_jacobian(Vector< double > &residuals, DenseMatrix< double > &jacobian)
Definition: pml_fourier_decomposed_helmholtz_elements.h:447
virtual void output_fct(std::ostream &outfile, const unsigned &n_plot, const double &time, FiniteElement::UnsteadyExactSolutionFctPt exact_soln_pt)
Definition: pml_fourier_decomposed_helmholtz_elements.h:298
PMLFourierDecomposedHelmholtzEquations(const PMLFourierDecomposedHelmholtzEquations &dummy)=delete
Broken copy constructor.
double k_squared()
Get k squared.
Definition: pml_fourier_decomposed_helmholtz_elements.h:212
virtual void get_source_pml_fourier_decomposed_helmholtz(const unsigned &ipt, const Vector< double > &x, std::complex< double > &source) const
Definition: pml_fourier_decomposed_helmholtz_elements.h:364
int * N_pml_fourier_pt
Pointer to Fourier wave number.
Definition: pml_fourier_decomposed_helmholtz_elements.h:667
PMLFourierDecomposedHelmholtzSourceFctPt Source_fct_pt
Pointer to source function:
Definition: pml_fourier_decomposed_helmholtz_elements.h:650
PMLFourierDecomposedHelmholtzSourceFctPt & source_fct_pt()
Access function: Pointer to source function.
Definition: pml_fourier_decomposed_helmholtz_elements.h:347
PMLFourierDecomposedHelmholtzSourceFctPt source_fct_pt() const
Access function: Pointer to source function. Const version.
Definition: pml_fourier_decomposed_helmholtz_elements.h:354
int *& pml_fourier_wavenumber_pt()
Get pointer to Fourier wavenumber.
Definition: pml_fourier_decomposed_helmholtz_elements.h:240
void compute_norm(double &norm)
Compute norm of fe solution.
Definition: pml_fourier_decomposed_helmholtz_elements.cc:795
void output_fct(std::ostream &outfile, const unsigned &n_plot, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt)
Definition: pml_fourier_decomposed_helmholtz_elements.cc:609
void fill_in_contribution_to_residuals(Vector< double > &residuals)
Add the element's contribution to its residual vector (wrapper)
Definition: pml_fourier_decomposed_helmholtz_elements.h:436
void values_to_be_pinned_on_outer_pml_boundary(Vector< unsigned > &values_to_pin)
Definition: pml_fourier_decomposed_helmholtz_elements.h:384
void compute_pml_coefficients(const unsigned &ipt, const Vector< double > &x, Vector< std::complex< double >> &pml_laplace_factor, std::complex< double > &pml_k_squared_factor)
Definition: pml_fourier_decomposed_helmholtz_elements.h:502
void output(FILE *file_pt)
C_style output with default number of plot points.
Definition: pml_fourier_decomposed_helmholtz_elements.h:280
double *& k_squared_pt()
Get pointer to frequency.
Definition: pml_fourier_decomposed_helmholtz_elements.h:205
void output_real_fct(std::ostream &outfile, const double &phi, const unsigned &n_plot, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt)
Definition: pml_fourier_decomposed_helmholtz_elements.cc:663
double alpha()
Get complex shift.
Definition: pml_fourier_decomposed_helmholtz_elements.h:234
virtual double dshape_and_dtest_eulerian_at_knot_pml_fourier_decomposed_helmholtz(const unsigned &ipt, Shape &psi, DShape &dpsidx, Shape &test, DShape &dtestdx) const =0
PMLMappingAndTransformedCoordinate *& pml_mapping_and_transformed_coordinate_pt()
Return a pointer to the PML Mapping object.
Definition: pml_fourier_decomposed_helmholtz_elements.h:605
static BermudezPMLMappingAndTransformedCoordinate Default_pml_mapping_and_transformed_coordinate
Definition: pml_fourier_decomposed_helmholtz_elements.h:620
PMLFourierDecomposedHelmholtzEquations()
Constructor.
Definition: pml_fourier_decomposed_helmholtz_elements.h:170
PMLMappingAndTransformedCoordinate * Pml_mapping_and_transformed_coordinate_pt
Definition: pml_fourier_decomposed_helmholtz_elements.h:658
void get_flux(const Vector< double > &s, Vector< std::complex< double >> &flux) const
Get flux: flux[i] = du/dx_i for real and imag part.
Definition: pml_fourier_decomposed_helmholtz_elements.h:396
virtual std::complex< unsigned > u_index_pml_fourier_decomposed_helmholtz() const
Broken assignment operator.
Definition: pml_fourier_decomposed_helmholtz_elements.h:197
PMLLayerElement()
Definition: pml_fourier_decomposed_helmholtz_elements.h:1142
Definition: pml_meshes.h:48
Definition: pml_fourier_decomposed_helmholtz_elements.h:80
virtual std::complex< double > gamma(const double &nu_i, const double &pml_width_i, const double &k_squared)=0
PMLMappingAndTransformedCoordinate()
Default constructor (empty)
Definition: pml_fourier_decomposed_helmholtz_elements.h:83
virtual std::complex< double > transformed_coordinate(const double &nu_i, const double &pml_width_i, const double &pml_inner_boundary, const double &k_squared)=0
Definition: projection.h:183
Fourier decomposed Helmholtz upgraded to become projectable.
Definition: pml_fourier_decomposed_helmholtz_elements.h:902
unsigned nhistory_values_for_projection(const unsigned &fld)
Definition: pml_fourier_decomposed_helmholtz_elements.h:948
int local_equation(const unsigned &fld, const unsigned &j)
Return local equation number of value j in field fld.
Definition: pml_fourier_decomposed_helmholtz_elements.h:1064
double get_field(const unsigned &t, const unsigned &fld, const Vector< double > &s)
Definition: pml_fourier_decomposed_helmholtz_elements.h:999
ProjectablePMLFourierDecomposedHelmholtzElement()
Definition: pml_fourier_decomposed_helmholtz_elements.h:906
unsigned nfields_for_projection()
Number of fields to be projected: 2 (real and imag part)
Definition: pml_fourier_decomposed_helmholtz_elements.h:941
unsigned nhistory_values_for_coordinate_projection()
Definition: pml_fourier_decomposed_helmholtz_elements.h:965
double jacobian_and_shape_of_field(const unsigned &fld, const Vector< double > &s, Shape &psi)
Definition: pml_fourier_decomposed_helmholtz_elements.h:972
Vector< std::pair< Data *, unsigned > > data_values_of_field(const unsigned &fld)
Definition: pml_fourier_decomposed_helmholtz_elements.h:911
void output(std::ostream &outfile, const unsigned &nplot)
Definition: pml_fourier_decomposed_helmholtz_elements.h:1093
unsigned nvalue_of_field(const unsigned &fld)
Return number of values in field fld: One per node.
Definition: pml_fourier_decomposed_helmholtz_elements.h:1047
Definition: Qelements.h:459
Definition: pml_fourier_decomposed_helmholtz_elements.h:685
double dshape_and_dtest_eulerian_pml_fourier_decomposed_helmholtz(const Vector< double > &s, Shape &psi, DShape &dpsidx, Shape &test, DShape &dtestdx) const
Definition: pml_fourier_decomposed_helmholtz_elements.h:823
void output(std::ostream &outfile, const unsigned &n_plot)
Definition: pml_fourier_decomposed_helmholtz_elements.h:723
void output_fct(std::ostream &outfile, const unsigned &n_plot, const double &time, FiniteElement::UnsteadyExactSolutionFctPt exact_soln_pt)
Definition: pml_fourier_decomposed_helmholtz_elements.h:781
static const unsigned Initial_Nvalue
Definition: pml_fourier_decomposed_helmholtz_elements.h:689
void output_fct(std::ostream &outfile, const unsigned &n_plot, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt)
Definition: pml_fourier_decomposed_helmholtz_elements.h:755
void output_real_fct(std::ostream &outfile, const double &phi, const unsigned &n_plot, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt)
Definition: pml_fourier_decomposed_helmholtz_elements.h:768
QPMLFourierDecomposedHelmholtzElement(const QPMLFourierDecomposedHelmholtzElement< NNODE_1D > &dummy)=delete
Broken copy constructor.
void output(FILE *file_pt)
C-style output function: r,z,u.
Definition: pml_fourier_decomposed_helmholtz_elements.h:741
unsigned required_nvalue(const unsigned &n) const
Broken assignment operator.
Definition: pml_fourier_decomposed_helmholtz_elements.h:710
double dshape_and_dtest_eulerian_at_knot_pml_fourier_decomposed_helmholtz(const unsigned &ipt, Shape &psi, DShape &dpsidx, Shape &test, DShape &dtestdx) const
Definition: pml_fourier_decomposed_helmholtz_elements.h:850
void output(FILE *file_pt, const unsigned &n_plot)
Definition: pml_fourier_decomposed_helmholtz_elements.h:748
QPMLFourierDecomposedHelmholtzElement()
Definition: pml_fourier_decomposed_helmholtz_elements.h:694
void output_real(std::ostream &outfile, const double &phi, const unsigned &n_plot)
Definition: pml_fourier_decomposed_helmholtz_elements.h:733
void output(std::ostream &outfile)
Output function: r,z,u.
Definition: pml_fourier_decomposed_helmholtz_elements.h:716
Definition: shape.h:76
unsigned ntstorage() const
Definition: timesteppers.h:601
static int f(const TensorMap< Tensor< int, 3 > > &tensor)
Definition: cxx11_tensor_map.cpp:237
float real
Definition: datatypes.h:10
RealScalar s
Definition: level1_cplx_impl.h:130
int * m
Definition: level2_cplx_impl.h:294
char char char int int * k
Definition: level2_impl.h:374
EIGEN_STRONG_INLINE EIGEN_DEVICE_FUNC bfloat16 log(const bfloat16 &a)
Definition: BFloat16.h:618
squared absolute sa ArrayBase::abs2 DOXCOMMA MatrixBase::cwiseAbs2 sa Eigen::abs2 DOXCOMMA Eigen::pow DOXCOMMA ArrayBase::square nearest sa Eigen::floor DOXCOMMA Eigen::ceil DOXCOMMA ArrayBase::round nearest integer not less than the given sa Eigen::floor DOXCOMMA ArrayBase::ceil not a number test
Definition: GlobalFunctions.h:109
void flux(const double &time, const Vector< double > &x, double &flux)
Get flux applied along boundary x=0.
Definition: pretend_melt.cc:59
void source(const Vector< double > &x, Vector< double > &f)
Source function.
Definition: unstructured_two_d_circle.cc:46
r
Definition: UniformPSDSelfTest.py:20
Real fabs(const Real &a)
Definition: boostmultiprec.cpp:117
int error
Definition: calibrate.py:297
double factorial(const unsigned &l)
Factorial.
Definition: fourier_decomposed_helmholtz_elements.cc:40
double plgndr1(const unsigned &n, const double &x)
Legendre polynomials depending on one parameter.
Definition: fourier_decomposed_helmholtz_elements.cc:50
double plgndr2(const unsigned &l, const unsigned &m, const double &x)
Legendre polynomials depending on two parameters.
Definition: fourier_decomposed_helmholtz_elements.cc:97
DRAIG: Change all instances of (SPATIAL_DIM) to (DIM-1).
Definition: AnisotropicHookean.h:10
list x
Definition: plotDoE.py:28
t
Definition: plotPSD.py:36
Definition: indexed_view.cpp:20
#define OOMPH_EXCEPTION_LOCATION
Definition: oomph_definitions.h:61
#define OOMPH_CURRENT_FUNCTION
Definition: oomph_definitions.h:86
void output(std::ostream &outfile, const unsigned &nplot)
Overload output function.
Definition: overloaded_element_body.h:490
EIGEN_DONT_INLINE Scalar zero()
Definition: svd_common.h:232
std::ptrdiff_t j
Definition: tut_arithmetic_redux_minmax.cpp:2