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_FOURIER_DECOMPOSED_HELMHOLTZ_ELEMENTS_HEADER
28 #define OOMPH_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 
46 namespace oomph
47 {
48  //========================================================================
50  //========================================================================
51  namespace Legendre_functions_helper
52  {
54  extern double factorial(const unsigned& l);
55 
57  extern double plgndr1(const unsigned& n, const double& x);
58 
60  extern double plgndr2(const unsigned& l,
61  const unsigned& m,
62  const double& x);
63 
64  } // namespace Legendre_functions_helper
65 
66 
70 
71 
72  //=============================================================
82  //=============================================================
84  {
85  public:
89  const Vector<double>& x, std::complex<double>& f);
90 
91 
95  {
96  }
97 
100  const FourierDecomposedHelmholtzEquations& dummy) = delete;
101 
103  // Commented out broken assignment operator because this can lead to a
104  // conflict warning when used in the virtual inheritence hierarchy.
105  // Essentially the compiler doesn't realise that two separate
106  // implementations of the broken function are the same and so, quite
107  // rightly, it shouts.
108  /*void operator=(const FourierDecomposedHelmholtzEquations&) = delete;*/
109 
110 
114  virtual inline std::complex<unsigned> u_index_fourier_decomposed_helmholtz()
115  const
116  {
117  return std::complex<unsigned>(0, 1);
118  }
119 
120 
122  double*& k_squared_pt()
123  {
124  return K_squared_pt;
125  }
126 
128  double k_squared()
129  {
130  if (K_squared_pt == 0)
131  {
132  return 0.0;
133  }
134  else
135  {
136  return *K_squared_pt;
137  }
138  }
139 
142  {
143  return N_fourier_pt;
144  }
145 
148  {
149  if (N_fourier_pt == 0)
150  {
151  return 0;
152  }
153  else
154  {
155  return *N_fourier_pt;
156  }
157  }
158 
160  void output(std::ostream& outfile)
161  {
162  const unsigned n_plot = 5;
163  output(outfile, n_plot);
164  }
165 
168  void output(std::ostream& outfile, const unsigned& n_plot);
169 
175  void output_real(std::ostream& outfile,
176  const double& phi,
177  const unsigned& n_plot);
178 
180  void output(FILE* file_pt)
181  {
182  const unsigned n_plot = 5;
183  output(file_pt, n_plot);
184  }
185 
188  void output(FILE* file_pt, const unsigned& n_plot);
189 
192  void output_fct(std::ostream& outfile,
193  const unsigned& n_plot,
195 
198  virtual void output_fct(
199  std::ostream& outfile,
200  const unsigned& n_plot,
201  const double& time,
203  {
204  throw OomphLibError("There is no time-dependent output_fct() for "
205  "FourierDecomposedHelmholtz elements ",
208  }
209 
210 
216  void output_real_fct(std::ostream& outfile,
217  const double& phi,
218  const unsigned& n_plot,
220 
221 
223  void compute_error(std::ostream& outfile,
225  double& error,
226  double& norm);
227 
228 
230  void compute_error(std::ostream& outfile,
232  const double& time,
233  double& error,
234  double& norm)
235  {
236  throw OomphLibError("There is no time-dependent compute_error() for "
237  "FourierDecomposedHelmholtz elements",
240  }
241 
243  void compute_norm(double& norm);
244 
247  {
248  return Source_fct_pt;
249  }
250 
253  {
254  return Source_fct_pt;
255  }
256 
262  const unsigned& ipt,
263  const Vector<double>& x,
264  std::complex<double>& source) const
265  {
266  // If no source function has been set, return zero
267  if (Source_fct_pt == 0)
268  {
269  source = std::complex<double>(0.0, 0.0);
270  }
271  else
272  {
273  // Get source strength
274  (*Source_fct_pt)(x, source);
275  }
276  }
277 
278 
280  void get_flux(const Vector<double>& s,
281  Vector<std::complex<double>>& flux) const
282  {
283  // Find out how many nodes there are in the element
284  const unsigned n_node = nnode();
285 
286  // Set up memory for the shape and test functions
287  Shape psi(n_node);
288  DShape dpsidx(n_node, 2);
289 
290  // Call the derivatives of the shape and test functions
291  dshape_eulerian(s, psi, dpsidx);
292 
293  // Initialise to zero
294  const std::complex<double> zero(0.0, 0.0);
295  for (unsigned j = 0; j < 2; j++)
296  {
297  flux[j] = zero;
298  }
299 
300  // Loop over nodes
301  for (unsigned l = 0; l < n_node; l++)
302  {
303  // Cache the complex value of the unknown
304  const std::complex<double> u_value(
307 
308  // Loop over derivative directions
309  for (unsigned j = 0; j < 2; j++)
310  {
311  flux[j] += u_value * dpsidx(l, j);
312  }
313  }
314  }
315 
316 
319  {
320  // Call the generic residuals function with flag set to 0
321  // using a dummy matrix argument
323  residuals, GeneralisedElement::Dummy_matrix, 0);
324  }
325 
326 
330  DenseMatrix<double>& jacobian)
331  {
332  // Call the generic routine with the flag set to 1
334  residuals, jacobian, 1);
335  }
336 
337 
340  inline std::complex<double> interpolated_u_fourier_decomposed_helmholtz(
341  const Vector<double>& s) const
342  {
343  // Find number of nodes
344  const unsigned n_node = nnode();
345 
346  // Local shape function
347  Shape psi(n_node);
348 
349  // Find values of shape function
350  shape(s, psi);
351 
352  // Initialise value of u
353  std::complex<double> interpolated_u(0.0, 0.0);
354 
355  // Get the index at which the helmholtz unknown is stored
356  const unsigned u_nodal_index_real =
358  const unsigned u_nodal_index_imag =
360 
361  // Loop over the local nodes and sum
362  for (unsigned l = 0; l < n_node; l++)
363  {
364  // Make a temporary complex number from the stored data
365  const std::complex<double> u_value(
366  this->nodal_value(l, u_nodal_index_real),
367  this->nodal_value(l, u_nodal_index_imag));
368  // Add to the interpolated value
369  interpolated_u += u_value * psi[l];
370  }
371  return interpolated_u;
372  }
373 
374 
376  unsigned self_test();
377 
378 
379  protected:
383  const Vector<double>& s,
384  Shape& psi,
385  DShape& dpsidx,
386  Shape& test,
387  DShape& dtestdx) const = 0;
388 
389 
393  const unsigned& ipt,
394  Shape& psi,
395  DShape& dpsidx,
396  Shape& test,
397  DShape& dtestdx) const = 0;
398 
402  Vector<double>& residuals,
403  DenseMatrix<double>& jacobian,
404  const unsigned& flag);
405 
408 
410  double* K_squared_pt;
411 
414  };
415 
416 
420 
421 
422  //======================================================================
426  //======================================================================
427  template<unsigned NNODE_1D>
429  : public virtual QElement<2, NNODE_1D>,
431  {
432  private:
435  static const unsigned Initial_Nvalue;
436 
437  public:
442  {
443  }
444 
447  const QFourierDecomposedHelmholtzElement<NNODE_1D>& dummy) = delete;
448 
450  /*void operator=(const QFourierDecomposedHelmholtzElement<NNODE_1D>&) =
451  * delete;*/
452 
453 
456  inline unsigned required_nvalue(const unsigned& n) const
457  {
458  return Initial_Nvalue;
459  }
460 
462  void output(std::ostream& outfile)
463  {
465  }
466 
469  void output(std::ostream& outfile, const unsigned& n_plot)
470  {
472  }
473 
479  void output_real(std::ostream& outfile,
480  const double& phi,
481  const unsigned& n_plot)
482  {
484  }
485 
487  void output(FILE* file_pt)
488  {
490  }
491 
494  void output(FILE* file_pt, const unsigned& n_plot)
495  {
497  }
498 
501  void output_fct(std::ostream& outfile,
502  const unsigned& n_plot,
504  {
506  outfile, n_plot, exact_soln_pt);
507  }
508 
514  void output_real_fct(std::ostream& outfile,
515  const double& phi,
516  const unsigned& n_plot,
518  {
520  outfile, phi, n_plot, exact_soln_pt);
521  }
522 
523 
527  void output_fct(std::ostream& outfile,
528  const unsigned& n_plot,
529  const double& time,
531  {
533  outfile, n_plot, time, exact_soln_pt);
534  }
535 
536  protected:
540  const Vector<double>& s,
541  Shape& psi,
542  DShape& dpsidx,
543  Shape& test,
544  DShape& dtestdx) const;
545 
546 
550  const unsigned& ipt,
551  Shape& psi,
552  DShape& dpsidx,
553  Shape& test,
554  DShape& dtestdx) const;
555  };
556 
557 
558  // Inline functions:
559 
560 
561  //======================================================================
566  //======================================================================
567  template<unsigned NNODE_1D>
570  const Vector<double>& s,
571  Shape& psi,
572  DShape& dpsidx,
573  Shape& test,
574  DShape& dtestdx) const
575  {
576  // Call the geometrical shape functions and derivatives
577  const double J = this->dshape_eulerian(s, psi, dpsidx);
578 
579  // Set the test functions equal to the shape functions
580  test = psi;
581  dtestdx = dpsidx;
582 
583  // Return the jacobian
584  return J;
585  }
586 
587 
588  //======================================================================
593  //======================================================================
594  template<unsigned NNODE_1D>
597  const unsigned& ipt,
598  Shape& psi,
599  DShape& dpsidx,
600  Shape& test,
601  DShape& dtestdx) const
602  {
603  // Call the geometrical shape functions and derivatives
604  const double J = this->dshape_eulerian_at_knot(ipt, psi, dpsidx);
605 
606  // Set the pointers of the test functions
607  test = psi;
608  dtestdx = dpsidx;
609 
610  // Return the jacobian
611  return J;
612  }
613 
617 
618 
619  //=======================================================================
624  //=======================================================================
625  template<unsigned NNODE_1D>
627  : public virtual QElement<1, NNODE_1D>
628  {
629  public:
632  FaceGeometry() : QElement<1, NNODE_1D>() {}
633  };
634 
635 
639 
640 
641  //==========================================================
643  //==========================================================
644  template<class FOURIER_DECOMPOSED_HELMHOLTZ_ELEMENT>
646  : public virtual ProjectableElement<FOURIER_DECOMPOSED_HELMHOLTZ_ELEMENT>
647  {
648  public:
652 
657  {
658 #ifdef PARANOID
659  if (fld > 1)
660  {
661  std::stringstream error_stream;
662  error_stream << "Fourier decomposed Helmholtz elements only store 2 "
663  "fields so fld = "
664  << fld << " is illegal \n";
665  throw OomphLibError(
666  error_stream.str(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
667  }
668 #endif
669 
670  // Create the vector
671  unsigned nnod = this->nnode();
672  Vector<std::pair<Data*, unsigned>> data_values(nnod);
673 
674  // Loop over all nodes
675  for (unsigned j = 0; j < nnod; j++)
676  {
677  // Add the data value associated field: The node itself
678  data_values[j] = std::make_pair(this->node_pt(j), fld);
679  }
680 
681  // Return the vector
682  return data_values;
683  }
684 
687  {
688  return 2;
689  }
690 
693  unsigned nhistory_values_for_projection(const unsigned& fld)
694  {
695 #ifdef PARANOID
696  if (fld > 1)
697  {
698  std::stringstream error_stream;
699  error_stream << "Helmholtz elements only store two fields so fld = "
700  << fld << " is illegal\n";
701  throw OomphLibError(
702  error_stream.str(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
703  }
704 #endif
705  return this->node_pt(0)->ntstorage();
706  }
707 
711  {
712  return this->node_pt(0)->position_time_stepper_pt()->ntstorage();
713  }
714 
717  double jacobian_and_shape_of_field(const unsigned& fld,
718  const Vector<double>& s,
719  Shape& psi)
720  {
721 #ifdef PARANOID
722  if (fld > 1)
723  {
724  std::stringstream error_stream;
725  error_stream << "Helmholtz elements only store two fields so fld = "
726  << fld << " is illegal.\n";
727  throw OomphLibError(
728  error_stream.str(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
729  }
730 #endif
731  unsigned n_dim = this->dim();
732  unsigned n_node = this->nnode();
733  Shape test(n_node);
734  DShape dpsidx(n_node, n_dim), dtestdx(n_node, n_dim);
735  double J = this->dshape_and_dtest_eulerian_fourier_decomposed_helmholtz(
736  s, psi, dpsidx, test, dtestdx);
737  return J;
738  }
739 
740 
743  double get_field(const unsigned& t,
744  const unsigned& fld,
745  const Vector<double>& s)
746  {
747 #ifdef PARANOID
748  if (fld > 1)
749  {
750  std::stringstream error_stream;
751  error_stream << "Helmholtz elements only store two fields so fld = "
752  << fld << " is illegal\n";
753  throw OomphLibError(
754  error_stream.str(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
755  }
756 #endif
757  // Find the index at which the variable is stored
758  std::complex<unsigned> complex_u_nodal_index =
759  this->u_index_fourier_decomposed_helmholtz();
760  unsigned u_nodal_index = 0;
761  if (fld == 0)
762  {
763  u_nodal_index = complex_u_nodal_index.real();
764  }
765  else
766  {
767  u_nodal_index = complex_u_nodal_index.imag();
768  }
769 
770 
771  // Local shape function
772  unsigned n_node = this->nnode();
773  Shape psi(n_node);
774 
775  // Find values of shape function
776  this->shape(s, psi);
777 
778  // Initialise value of u
779  double interpolated_u = 0.0;
780 
781  // Sum over the local nodes
782  for (unsigned l = 0; l < n_node; l++)
783  {
784  interpolated_u += this->nodal_value(t, l, u_nodal_index) * psi[l];
785  }
786  return interpolated_u;
787  }
788 
789 
791  unsigned nvalue_of_field(const unsigned& fld)
792  {
793 #ifdef PARANOID
794  if (fld > 1)
795  {
796  std::stringstream error_stream;
797  error_stream << "Helmholtz elements only store two fields so fld = "
798  << fld << " is illegal\n";
799  throw OomphLibError(
800  error_stream.str(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
801  }
802 #endif
803  return this->nnode();
804  }
805 
806 
808  int local_equation(const unsigned& fld, const unsigned& j)
809  {
810 #ifdef PARANOID
811  if (fld > 1)
812  {
813  std::stringstream error_stream;
814  error_stream << "Helmholtz elements only store two fields so fld = "
815  << fld << " is illegal\n";
816  throw OomphLibError(
817  error_stream.str(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
818  }
819 #endif
820  std::complex<unsigned> complex_u_nodal_index =
821  this->u_index_fourier_decomposed_helmholtz();
822  unsigned u_nodal_index = 0;
823  if (fld == 0)
824  {
825  u_nodal_index = complex_u_nodal_index.real();
826  }
827  else
828  {
829  u_nodal_index = complex_u_nodal_index.imag();
830  }
831  return this->nodal_local_eqn(j, u_nodal_index);
832  }
833 
834 
837  void output(std::ostream& outfile, const unsigned& nplot)
838  {
840  }
841  };
842 
843 
844  //=======================================================================
847  //=======================================================================
848  template<class ELEMENT>
850  : public virtual FaceGeometry<ELEMENT>
851  {
852  public:
853  FaceGeometry() : FaceGeometry<ELEMENT>() {}
854  };
855 
856 
857  //=======================================================================
860  //=======================================================================
861  template<class ELEMENT>
864  : public virtual FaceGeometry<FaceGeometry<ELEMENT>>
865  {
866  public:
868  };
869 
870 
871 } // namespace oomph
872 
873 #endif
AnnoyingScalar imag(const AnnoyingScalar &)
Definition: AnnoyingScalar.h:132
const unsigned n
Definition: CG3DPackingUnitTest.cpp:11
JacobiRotation< float > J
Definition: Jacobi_makeJacobi.cpp:3
Definition: shape.h:278
unsigned ntstorage() const
Definition: nodes.cc:879
FaceGeometry()
Definition: fourier_decomposed_helmholtz_elements.h:853
FaceGeometry()
Definition: fourier_decomposed_helmholtz_elements.h:632
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
Definition: fourier_decomposed_helmholtz_elements.h:84
FourierDecomposedHelmholtzSourceFctPt Source_fct_pt
Pointer to source function:
Definition: fourier_decomposed_helmholtz_elements.h:407
void output_real_fct(std::ostream &outfile, const double &phi, const unsigned &n_plot, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt)
Definition: fourier_decomposed_helmholtz_elements.cc:574
virtual void fill_in_generic_residual_contribution_fourier_decomposed_helmholtz(Vector< double > &residuals, DenseMatrix< double > &jacobian, const unsigned &flag)
Definition: fourier_decomposed_helmholtz_elements.cc:183
int * N_fourier_pt
Pointer to Fourier wave number.
Definition: fourier_decomposed_helmholtz_elements.h:413
virtual void output_fct(std::ostream &outfile, const unsigned &n_plot, const double &time, FiniteElement::UnsteadyExactSolutionFctPt exact_soln_pt)
Definition: fourier_decomposed_helmholtz_elements.h:198
double * K_squared_pt
Pointer to square of wavenumber.
Definition: fourier_decomposed_helmholtz_elements.h:410
unsigned self_test()
Self-test: Return 0 for OK.
Definition: fourier_decomposed_helmholtz_elements.cc:373
int *& fourier_wavenumber_pt()
Get pointer to Fourier wavenumber.
Definition: fourier_decomposed_helmholtz_elements.h:141
void(* FourierDecomposedHelmholtzSourceFctPt)(const Vector< double > &x, std::complex< double > &f)
Definition: fourier_decomposed_helmholtz_elements.h:88
void output_real(std::ostream &outfile, const double &phi, const unsigned &n_plot)
Definition: fourier_decomposed_helmholtz_elements.cc:441
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: fourier_decomposed_helmholtz_elements.h:280
void fill_in_contribution_to_jacobian(Vector< double > &residuals, DenseMatrix< double > &jacobian)
Definition: fourier_decomposed_helmholtz_elements.h:329
void compute_norm(double &norm)
Compute norm of fe solution.
Definition: fourier_decomposed_helmholtz_elements.cc:706
double k_squared()
Get the square of wavenumber.
Definition: fourier_decomposed_helmholtz_elements.h:128
double *& k_squared_pt()
Get pointer to square of wavenumber.
Definition: fourier_decomposed_helmholtz_elements.h:122
void output_fct(std::ostream &outfile, const unsigned &n_plot, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt)
Definition: fourier_decomposed_helmholtz_elements.cc:520
virtual double dshape_and_dtest_eulerian_fourier_decomposed_helmholtz(const Vector< double > &s, Shape &psi, DShape &dpsidx, Shape &test, DShape &dtestdx) const =0
virtual double dshape_and_dtest_eulerian_at_knot_fourier_decomposed_helmholtz(const unsigned &ipt, Shape &psi, DShape &dpsidx, Shape &test, DShape &dtestdx) const =0
FourierDecomposedHelmholtzSourceFctPt source_fct_pt() const
Access function: Pointer to source function. Const version.
Definition: fourier_decomposed_helmholtz_elements.h:252
FourierDecomposedHelmholtzSourceFctPt & source_fct_pt()
Access function: Pointer to source function.
Definition: fourier_decomposed_helmholtz_elements.h:246
FourierDecomposedHelmholtzEquations(const FourierDecomposedHelmholtzEquations &dummy)=delete
Broken copy constructor.
void output(std::ostream &outfile)
Output with default number of plot points.
Definition: fourier_decomposed_helmholtz_elements.h:160
FourierDecomposedHelmholtzEquations()
Constructor.
Definition: fourier_decomposed_helmholtz_elements.h:93
std::complex< double > interpolated_u_fourier_decomposed_helmholtz(const Vector< double > &s) const
Definition: fourier_decomposed_helmholtz_elements.h:340
void compute_error(std::ostream &outfile, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt, double &error, double &norm)
Get error against and norm of exact solution.
Definition: fourier_decomposed_helmholtz_elements.cc:626
virtual std::complex< unsigned > u_index_fourier_decomposed_helmholtz() const
Broken assignment operator.
Definition: fourier_decomposed_helmholtz_elements.h:114
virtual void get_source_fourier_decomposed_helmholtz(const unsigned &ipt, const Vector< double > &x, std::complex< double > &source) const
Definition: fourier_decomposed_helmholtz_elements.h:261
void output(FILE *file_pt)
C_style output with default number of plot points.
Definition: fourier_decomposed_helmholtz_elements.h:180
void fill_in_contribution_to_residuals(Vector< double > &residuals)
Add the element's contribution to its residual vector (wrapper)
Definition: fourier_decomposed_helmholtz_elements.h:318
void compute_error(std::ostream &outfile, FiniteElement::UnsteadyExactSolutionFctPt exact_soln_pt, const double &time, double &error, double &norm)
Dummy, time dependent error checker.
Definition: fourier_decomposed_helmholtz_elements.h:230
int fourier_wavenumber()
Get the Fourier wavenumber.
Definition: fourier_decomposed_helmholtz_elements.h:147
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
Definition: projection.h:183
Fourier decomposed Helmholtz upgraded to become projectable.
Definition: fourier_decomposed_helmholtz_elements.h:647
double get_field(const unsigned &t, const unsigned &fld, const Vector< double > &s)
Definition: fourier_decomposed_helmholtz_elements.h:743
unsigned nfields_for_projection()
Number of fields to be projected: 2 (real and imag part)
Definition: fourier_decomposed_helmholtz_elements.h:686
ProjectableFourierDecomposedHelmholtzElement()
Definition: fourier_decomposed_helmholtz_elements.h:651
double jacobian_and_shape_of_field(const unsigned &fld, const Vector< double > &s, Shape &psi)
Definition: fourier_decomposed_helmholtz_elements.h:717
Vector< std::pair< Data *, unsigned > > data_values_of_field(const unsigned &fld)
Definition: fourier_decomposed_helmholtz_elements.h:656
int local_equation(const unsigned &fld, const unsigned &j)
Return local equation number of value j in field fld.
Definition: fourier_decomposed_helmholtz_elements.h:808
unsigned nhistory_values_for_coordinate_projection()
Definition: fourier_decomposed_helmholtz_elements.h:710
unsigned nvalue_of_field(const unsigned &fld)
Return number of values in field fld: One per node.
Definition: fourier_decomposed_helmholtz_elements.h:791
void output(std::ostream &outfile, const unsigned &nplot)
Definition: fourier_decomposed_helmholtz_elements.h:837
unsigned nhistory_values_for_projection(const unsigned &fld)
Definition: fourier_decomposed_helmholtz_elements.h:693
Definition: Qelements.h:459
Definition: fourier_decomposed_helmholtz_elements.h:431
void output_real(std::ostream &outfile, const double &phi, const unsigned &n_plot)
Definition: fourier_decomposed_helmholtz_elements.h:479
void output_fct(std::ostream &outfile, const unsigned &n_plot, const double &time, FiniteElement::UnsteadyExactSolutionFctPt exact_soln_pt)
Definition: fourier_decomposed_helmholtz_elements.h:527
void output(FILE *file_pt)
C-style output function: r,z,u.
Definition: fourier_decomposed_helmholtz_elements.h:487
void output_fct(std::ostream &outfile, const unsigned &n_plot, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt)
Definition: fourier_decomposed_helmholtz_elements.h:501
void output(std::ostream &outfile)
Output function: r,z,u.
Definition: fourier_decomposed_helmholtz_elements.h:462
void output_real_fct(std::ostream &outfile, const double &phi, const unsigned &n_plot, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt)
Definition: fourier_decomposed_helmholtz_elements.h:514
static const unsigned Initial_Nvalue
Definition: fourier_decomposed_helmholtz_elements.h:435
QFourierDecomposedHelmholtzElement()
Definition: fourier_decomposed_helmholtz_elements.h:440
void output(FILE *file_pt, const unsigned &n_plot)
Definition: fourier_decomposed_helmholtz_elements.h:494
double dshape_and_dtest_eulerian_at_knot_fourier_decomposed_helmholtz(const unsigned &ipt, Shape &psi, DShape &dpsidx, Shape &test, DShape &dtestdx) const
Definition: fourier_decomposed_helmholtz_elements.h:596
unsigned required_nvalue(const unsigned &n) const
Broken assignment operator.
Definition: fourier_decomposed_helmholtz_elements.h:456
void output(std::ostream &outfile, const unsigned &n_plot)
Definition: fourier_decomposed_helmholtz_elements.h:469
double dshape_and_dtest_eulerian_fourier_decomposed_helmholtz(const Vector< double > &s, Shape &psi, DShape &dpsidx, Shape &test, DShape &dtestdx) const
Definition: fourier_decomposed_helmholtz_elements.h:569
QFourierDecomposedHelmholtzElement(const QFourierDecomposedHelmholtzElement< NNODE_1D > &dummy)=delete
Broken copy constructor.
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
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
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