pml_fourier_decomposed_helmholtz_flux_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-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 // Header file for elements that are used to apply prescribed flux
27 // boundary conditions to the Fourier decomposed Helmholtz equations
28 #ifndef OOMPH_PML_FOURIER_DECOMPOSED_HELMHOLTZ_FLUX_ELEMENTS_HEADER
29 #define OOMPH_PML_FOURIER_DECOMPOSED_HELMHOLTZ_FLUX_ELEMENTS_HEADER
30 
31 
32 // Config header generated by autoconfig
33 #ifdef HAVE_CONFIG_H
34 #include <oomph-lib-config.h>
35 #endif
36 
37 #include "math.h"
38 #include <complex>
39 
40 // oomph-lib includes
41 #include "../generic/Qelements.h"
42 
43 
44 namespace oomph
45 {
46  //======================================================================
47  /// A class for elements that allow the imposition of an
48  /// applied flux on the boundaries of Fourier decomposed Helmholtz elements.
49  /// The element geometry is obtained from the FaceGeometry<ELEMENT>
50  /// policy class.
51  //======================================================================
52  template<class ELEMENT>
54  : public virtual FaceGeometry<ELEMENT>,
55  public virtual FaceElement
56  {
57  public:
58  /// Function pointer to the prescribed-flux function fct(x,f(x)) --
59  /// x is a Vector and the flux is a complex
60 
62  const Vector<double>& x, std::complex<double>& flux);
63 
64  /// Constructor, takes the pointer to the "bulk" element and the
65  /// index of the face to which the element is attached.
67  const int& face_index);
68 
69  /// Broken empty constructor
71  {
72  throw OomphLibError("Don't call empty constructor for "
73  "PMLFourierDecomposedHelmholtzFluxElement",
74  OOMPH_CURRENT_FUNCTION,
75  OOMPH_EXCEPTION_LOCATION);
76  }
77 
78  /// Broken copy constructor
80  const PMLFourierDecomposedHelmholtzFluxElement& dummy) = delete;
81 
82  /// Broken assignment operator
83  // Commented out broken assignment operator because this can lead to a
84  // conflict warning when used in the virtual inheritence hierarchy.
85  // Essentially the compiler doesn't realise that two separate
86  // implementations of the broken function are the same and so, quite
87  // rightly, it shouts.
88  /*void operator=(const PMLFourierDecomposedHelmholtzFluxElement&) =
89  * delete;*/
90 
91 
92  /// Access function for the prescribed-flux function pointer
94  {
95  return Flux_fct_pt;
96  }
97 
98 
99  /// Add the element's contribution to its residual vector
101  {
102  // Call the generic residuals function with flag set to 0
103  // using a dummy matrix argument
105  residuals, GeneralisedElement::Dummy_matrix, 0);
106  }
107 
108 
109  /// Add the element's contribution to its residual vector and its
110  /// Jacobian matrix
112  DenseMatrix<double>& jacobian)
113  {
114  // Call the generic routine with the flag set to 1
116  residuals, jacobian, 1);
117  }
118 
119  /// Output function -- forward to broken version in FiniteElement
120  /// until somebody decides what exactly they want to plot here...
121  void output(std::ostream& outfile)
122  {
123  FiniteElement::output(outfile);
124  }
125 
126  /// Output function -- forward to broken version in FiniteElement
127  /// until somebody decides what exactly they want to plot here...
128  void output(std::ostream& outfile, const unsigned& n_plot)
129  {
130  FiniteElement::output(outfile, n_plot);
131  }
132 
133 
134  /// C-style output function -- forward to broken version in FiniteElement
135  /// until somebody decides what exactly they want to plot here...
136  void output(FILE* file_pt)
137  {
138  FiniteElement::output(file_pt);
139  }
140 
141  /// C-style output function -- forward to broken version in
142  /// FiniteElement until somebody decides what exactly they want to plot
143  /// here...
144  void output(FILE* file_pt, const unsigned& n_plot)
145  {
146  FiniteElement::output(file_pt, n_plot);
147  }
148 
149 
150  /// Return the index at which the unknown value
151  /// is stored. (Real/imag part gives real index of real/imag part).
152  virtual inline std::complex<unsigned> u_index_pml_fourier_decomposed_helmholtz()
153  const
154  {
155  return std::complex<unsigned>(
158  }
159 
160 
161  protected:
162  /// Function to compute the shape and test functions and to return
163  /// the Jacobian of mapping between local and global (Eulerian)
164  /// coordinates
165  inline double shape_and_test(const Vector<double>& s,
166  Shape& psi,
167  Shape& test) const
168  {
169  // Find number of nodes
170  unsigned n_node = nnode();
171 
172  // Get the shape functions
173  shape(s, psi);
174 
175  // Set the test functions to be the same as the shape functions
176  for (unsigned i = 0; i < n_node; i++)
177  {
178  test[i] = psi[i];
179  }
180 
181  // Return the value of the jacobian
182  return J_eulerian(s);
183  }
184 
185 
186  /// Function to compute the shape and test functions and to return
187  /// the Jacobian of mapping between local and global (Eulerian)
188  /// coordinates
189  inline double shape_and_test_at_knot(const unsigned& ipt,
190  Shape& psi,
191  Shape& test) const
192  {
193  // Find number of nodes
194  unsigned n_node = nnode();
195 
196  // Get the shape functions
197  shape_at_knot(ipt, psi);
198 
199  // Set the test functions to be the same as the shape functions
200  for (unsigned i = 0; i < n_node; i++)
201  {
202  test[i] = psi[i];
203  }
204 
205  // Return the value of the jacobian
206  return J_eulerian_at_knot(ipt);
207  }
208 
209 
210  /// Function to calculate the prescribed flux at a given spatial
211  /// position
212  void get_flux(const Vector<double>& x, std::complex<double>& flux)
213  {
214  // If the function pointer is zero return zero
215  if (Flux_fct_pt == 0)
216  {
217  flux = std::complex<double>(0.0, 0.0);
218  }
219  // Otherwise call the function
220  else
221  {
222  (*Flux_fct_pt)(x, flux);
223  }
224  }
225 
226 
227  /// The index at which the real and imag part of the
228  /// unknown is stored at the nodes
230 
231 
232  /// Add the element's contribution to its residual vector.
233  /// flag=1(or 0): do (or don't) compute the contribution to the
234  /// Jacobian as well.
236  Vector<double>& residuals,
237  DenseMatrix<double>& jacobian,
238  const unsigned& flag);
239 
240 
241  /// Function pointer to the (global) prescribed-flux function
243  };
244 
245  /// ///////////////////////////////////////////////////////////////////
246  /// ///////////////////////////////////////////////////////////////////
247  /// ///////////////////////////////////////////////////////////////////
248 
249 
250  //===========================================================================
251  /// Constructor, takes the pointer to the "bulk" element, the
252  /// index of the fixed local coordinate and its value represented
253  /// by an integer (+/- 1), indicating that the face is located
254  /// at the max. or min. value of the "fixed" local coordinate
255  /// in the bulk element.
256  //===========================================================================
257  template<class ELEMENT>
260  const int& face_index)
261  : FaceGeometry<ELEMENT>(), FaceElement()
262  {
263  // Let the bulk element build the FaceElement, i.e. setup the pointers
264  // to its nodes (by referring to the appropriate nodes in the bulk
265  // element), etc.
266  bulk_el_pt->build_face_element(face_index, this);
267 
268  // Initialise the prescribed-flux function pointer to zero
269  Flux_fct_pt = 0;
270 
271  // Initialise index at which real and imag unknowns are stored
272  U_index_pml_fourier_decomposed_helmholtz = std::complex<unsigned>(0, 1);
273 
274  // Now read out indices from bulk element
276  dynamic_cast<PMLFourierDecomposedHelmholtzEquations*>(bulk_el_pt);
277  // If the cast has failed die
278  if (eqn_pt == 0)
279  {
280  std::string error_string = "Bulk element must inherit from "
281  "PMLFourierDecomposedHelmholtzEquations.";
282  throw OomphLibError(
283  error_string, OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
284  }
285  else
286  {
287  // Read the index from the (cast) bulk element
290  }
291  }
292 
293 
294  //===========================================================================
295  /// Compute the element's residual vector and the (zero) Jacobian matrix.
296  //===========================================================================
297  template<class ELEMENT>
300  Vector<double>& residuals,
301  DenseMatrix<double>& jacobian,
302  const unsigned& flag)
303  {
304  // Find out how many nodes there are3
305  const unsigned n_node = nnode();
306 
307  // Set up memory for the shape and test functions
308  Shape psif(n_node), testf(n_node);
309 
310  // Set the value of Nintpt
311  const unsigned n_intpt = integral_pt()->nweight();
312 
313  // Set the Vector to hold local coordinates
314  Vector<double> s(1);
315 
316  // Integers to hold the local equation and unknown numbers
317  int local_eqn_real = 0, local_eqn_imag = 0;
318 
319  // Loop over the integration points
320  //--------------------------------
321  for (unsigned ipt = 0; ipt < n_intpt; ipt++)
322  {
323  // Assign values of s
324  for (unsigned i = 0; i < 1; i++)
325  {
326  s[i] = integral_pt()->knot(ipt, i);
327  }
328 
329  // Get the integral weight
330  double w = integral_pt()->weight(ipt);
331 
332  // Find the shape and test functions and return the Jacobian
333  // of the mapping
334  double J = shape_and_test(s, psif, testf);
335 
336  // Premultiply the weights and the Jacobian
337  double W = w * J;
338 
339  // Need to find position to feed into flux function, initialise to zero
340  Vector<double> interpolated_x(2, 0.0);
341 
342  // Calculate coordinates
343  for (unsigned l = 0; l < n_node; l++)
344  {
345  // Loop over coordinates
346  for (unsigned i = 0; i < 2; i++)
347  {
348  interpolated_x[i] += nodal_position(l, i) * psif[l];
349  }
350  }
351 
352  // first component
353  double r = interpolated_x[0];
354 
355  // Get the imposed flux
356  std::complex<double> flux(0.0, 0.0);
357  get_flux(interpolated_x, flux);
358 
359  // Now add to the appropriate equations
360  // Loop over the test functions
361  for (unsigned l = 0; l < n_node; l++)
362  {
363  local_eqn_real =
364  nodal_local_eqn(l, U_index_pml_fourier_decomposed_helmholtz.real());
365 
366  /*IF it's not a boundary condition*/
367  if (local_eqn_real >= 0)
368  {
369  // Add the prescribed flux terms
370  residuals[local_eqn_real] -= flux.real() * testf[l] * r * W;
371 
372  // Imposed traction doesn't depend upon the solution,
373  // --> the Jacobian is always zero, so no Jacobian
374  // terms are required
375  }
376  local_eqn_imag =
377  nodal_local_eqn(l, U_index_pml_fourier_decomposed_helmholtz.imag());
378 
379  /*IF it's not a boundary condition*/
380  if (local_eqn_imag >= 0)
381  {
382  // Add the prescribed flux terms
383  residuals[local_eqn_imag] -= flux.imag() * testf[l] * r * W;
384 
385  // Imposed traction doesn't depend upon the solution,
386  // --> the Jacobian is always zero, so no Jacobian
387  // terms are required
388  }
389  }
390  }
391  }
392 
393 
394  /// //////////////////////////////////////////////////////////////////
395  /// //////////////////////////////////////////////////////////////////
396  /// //////////////////////////////////////////////////////////////////
397 
398 
399  //======================================================================
400  /// A class for elements that allow postprocessing of the
401  /// results -- currently computes radiated power over domain
402  /// boundaries.
403  /// The element geometry is obtained from the FaceGeometry<ELEMENT>
404  /// policy class.
405  //======================================================================
406  template<class ELEMENT>
408  : public virtual FaceGeometry<ELEMENT>,
409  public virtual FaceElement
410  {
411  public:
412  /// Constructor, takes the pointer to the "bulk" element and the
413  /// index of the face to which the element is attached.
415  FiniteElement* const& bulk_el_pt, const int& face_index);
416 
417  /// Broken empty constructor
419  {
420  throw OomphLibError("Don't call empty constructor for "
421  "PMLFourierDecomposedHelmholtzPowerMonitorElement",
422  OOMPH_CURRENT_FUNCTION,
423  OOMPH_EXCEPTION_LOCATION);
424  }
425 
426  /// Broken copy constructor
429 
430  /// Broken assignment operator
431  /*void operator=(const PMLFourierDecomposedHelmholtzPowerMonitorElement&) =
432  * delete;*/
433 
434 
435  /// Specify the value of nodal zeta from the face geometry
436  /// The "global" intrinsic coordinate of the element when
437  /// viewed as part of a geometric object should be given by
438  /// the FaceElement representation, by default (needed to break
439  /// indeterminacy if bulk element is SolidElement)
440  double zeta_nodal(const unsigned& n,
441  const unsigned& k,
442  const unsigned& i) const
443  {
444  return FaceElement::zeta_nodal(n, k, i);
445  }
446 
447 
448  /// Output function -- forward to broken version in FiniteElement
449  /// until somebody decides what exactly they want to plot here...
450  void output(std::ostream& outfile)
451  {
452  FiniteElement::output(outfile);
453  }
454 
455  /// Output function -- forward to broken version in FiniteElement
456  /// until somebody decides what exactly they want to plot here...
457  void output(std::ostream& outfile, const unsigned& n_plot)
458  {
459  FiniteElement::output(outfile, n_plot);
460  }
461 
462  /// C-style output function -- forward to broken version in FiniteElement
463  /// until somebody decides what exactly they want to plot here...
464  void output(FILE* file_pt)
465  {
466  FiniteElement::output(file_pt);
467  }
468 
469  /// C-style output function -- forward to broken version in
470  /// FiniteElement until somebody decides what exactly they want to plot
471  /// here...
472  void output(FILE* file_pt, const unsigned& n_plot)
473  {
474  FiniteElement::output(file_pt, n_plot);
475  }
476 
477  /// Return the index at which the real/imag unknown value
478  /// is stored.
479  virtual inline std::complex<unsigned> u_index_pml_fourier_decomposed_helmholtz()
480  const
481  {
482  return std::complex<unsigned>(
485  }
486 
487  /// Compute the element's contribution to the time-averaged
488  /// radiated power over the artificial boundary.
489  /// NOTE: This may give the wrong result
490  /// if the constitutive parameters genuinely vary!
492  {
493  // Dummy output file
494  std::ofstream outfile;
495  return global_power_contribution(outfile);
496  }
497 
498  /// Compute the element's contribution to the time-averaged
499  /// radiated power over the artificial boundary. Also output the
500  /// power density as a fct of the zenith angle in the specified
501  /// output file if it's open. NOTE: This may give the wrong result
502  /// if the constitutive parameters genuinely vary!
503  double global_power_contribution(std::ofstream& outfile)
504  {
505  // pointer to the corresponding bulk element
506  ELEMENT* bulk_elem_pt = dynamic_cast<ELEMENT*>(this->bulk_element_pt());
507 
508  // Number of nodes in bulk element
509  unsigned nnode_bulk = bulk_elem_pt->nnode();
510  const unsigned n_node_local = this->nnode();
511 
512  // get the dim of the bulk and local nodes
513  const unsigned bulk_dim = bulk_elem_pt->dim();
514  const unsigned local_dim = this->node_pt(0)->ndim();
515 
516  // Set up memory for the shape and test functions
517  Shape psi(n_node_local);
518 
519  // Set up memory for the shape functions
520  Shape psi_bulk(nnode_bulk);
521  DShape dpsi_bulk_dx(nnode_bulk, bulk_dim);
522 
523  // Set up memory for the outer unit normal
524  Vector<double> unit_normal(bulk_dim);
525 
526  // Set the value of n_intpt
527  const unsigned n_intpt = integral_pt()->nweight();
528 
529  // Set the Vector to hold local coordinates
530  Vector<double> s(local_dim - 1);
531  double power = 0.0;
532 
533  // Output?
534  if (outfile.is_open())
535  {
536  outfile << "ZONE\n";
537  }
538 
539  // Loop over the integration points
540  //--------------------------------
541  for (unsigned ipt = 0; ipt < n_intpt; ipt++)
542  {
543  // Assign values of s
544  for (unsigned i = 0; i < (local_dim - 1); i++)
545  {
546  s[i] = integral_pt()->knot(ipt, i);
547  }
548  // get the outer_unit_ext vector
549  this->outer_unit_normal(s, unit_normal);
550 
551  // Get the integral weight
552  double w = integral_pt()->weight(ipt);
553 
554  // Get jacobian of mapping
555  double J = J_eulerian(s);
556 
557  // Premultiply the weights and the Jacobian
558  double W = w * J;
559 
560  // Get local coordinates in bulk element by copy construction
562 
563  // Call the derivatives of the shape functions
564  // in the bulk -- must do this via s because this point
565  // is not an integration point the bulk element!
566  (void)bulk_elem_pt->dshape_eulerian(s_bulk, psi_bulk, dpsi_bulk_dx);
567  this->shape(s, psi);
568 
569  // Derivs of Eulerian coordinates w.r.t. local coordinates
570  std::complex<double> dphi_dn(0.0, 0.0);
571  Vector<std::complex<double>> interpolated_dphidx(bulk_dim);
572  std::complex<double> interpolated_phi(0.0, 0.0);
573  Vector<double> x(bulk_dim);
574 
575  // Calculate function value and derivatives:
576  //-----------------------------------------
577  // Loop over nodes
578  for (unsigned l = 0; l < nnode_bulk; l++)
579  {
580  // Get the nodal value of the helmholtz unknown
581  const std::complex<double> phi_value(
582  bulk_elem_pt->nodal_value(
583  l,
584  bulk_elem_pt->u_index_pml_fourier_decomposed_helmholtz().real()),
585  bulk_elem_pt->nodal_value(
586  l,
587  bulk_elem_pt->u_index_pml_fourier_decomposed_helmholtz().imag()));
588 
589  // Loop over directions
590  for (unsigned i = 0; i < bulk_dim; i++)
591  {
592  interpolated_dphidx[i] += phi_value * dpsi_bulk_dx(l, i);
593  }
594  } // End of loop over the bulk_nodes
595 
596  for (unsigned l = 0; l < n_node_local; l++)
597  {
598  // Get the nodal value of the Helmholtz unknown
599  const std::complex<double> phi_value(
602 
603  interpolated_phi += phi_value * psi(l);
604  }
605 
606  // define dphi_dn
607  for (unsigned i = 0; i < bulk_dim; i++)
608  {
609  dphi_dn += interpolated_dphidx[i] * unit_normal[i];
610  }
611 
612  // Power density
613  double integrand = (interpolated_phi.real() * dphi_dn.imag() -
614  interpolated_phi.imag() * dphi_dn.real());
615 
616  interpolated_x(s, x);
617  double theta = atan2(x[0], x[1]);
618  // Output?
619  if (outfile.is_open())
620  {
621  outfile << x[0] << " " << x[1] << " " << theta << " " << integrand
622  << "\n";
623  }
624 
625  // ...add to integral
626  power += MathematicalConstants::Pi * x[0] * integrand * W;
627  }
628 
629  return power;
630  }
631 
632  protected:
633  /// Function to compute the test functions and to return
634  /// the Jacobian of mapping between local and global (Eulerian)
635  /// coordinates
636  inline double shape_and_test(const Vector<double>& s,
637  Shape& psi,
638  Shape& test) const
639  {
640  // Get the shape functions
641  shape(s, test);
642 
643  unsigned nnod = nnode();
644  for (unsigned i = 0; i < nnod; i++)
645  {
646  psi[i] = test[i];
647  }
648 
649  // Return the value of the jacobian
650  return J_eulerian(s);
651  }
652 
653  /// Function to compute the shape, test functions and derivates
654  /// and to return
655  /// the Jacobian of mapping between local and global (Eulerian)
656  /// coordinates
658  Shape& psi,
659  Shape& test,
660  DShape& dpsi_ds,
661  DShape& dtest_ds) const
662  {
663  // Find number of nodes
664  unsigned n_node = nnode();
665 
666  // Get the shape functions
667  dshape_local(s, psi, dpsi_ds);
668 
669  // Set the test functions to be the same as the shape functions
670  for (unsigned i = 0; i < n_node; i++)
671  {
672  for (unsigned j = 0; j < 1; j++)
673  {
674  test[i] = psi[i];
675  dtest_ds(i, j) = dpsi_ds(i, j);
676  }
677  }
678  // Return the value of the jacobian
679  return J_eulerian(s);
680  }
681 
682  /// The index at which the real and imag part of the unknown
683  /// is stored at the nodes
685  };
686 
687 
688  /// ///////////////////////////////////////////////////////////////////
689  /// ///////////////////////////////////////////////////////////////////
690  /// ///////////////////////////////////////////////////////////////////
691 
692 
693  //===========================================================================
694  /// Constructor, takes the pointer to the "bulk" element and the face index.
695  //===========================================================================
696  template<class ELEMENT>
699  FiniteElement* const& bulk_el_pt, const int& face_index)
700  : FaceGeometry<ELEMENT>(), FaceElement()
701  {
702  // Let the bulk element build the FaceElement, i.e. setup the pointers
703  // to its nodes (by referring to the appropriate nodes in the bulk
704  // element), etc.
705  bulk_el_pt->build_face_element(face_index, this);
706 
707  // Set up U_index_pml_fourier_decomposedhelmholtz.
708  U_index_pml_fourier_decomposed_helmholtz = std::complex<unsigned>(0, 1);
709 
710  // Cast to the appropriate PMLFourierDecomposedHelmholtzEquation
711  // so that we can find the index at which the variable is stored
712  // We assume that the dimension of the full problem is the same
713  // as the dimension of the node, if this is not the case you will have
714  // to write custom elements, sorry
716  dynamic_cast<PMLFourierDecomposedHelmholtzEquations*>(bulk_el_pt);
717  if (eqn_pt == 0)
718  {
719  std::string error_string = "Bulk element must inherit from "
720  "PMLFourierDecomposedHelmholtzEquations.";
721  throw OomphLibError(
722  error_string, OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
723  }
724  // Otherwise read out the value
725  else
726  {
727  // Read the index from the (cast) bulk element
730  }
731  }
732 
733 } // namespace oomph
734 
735 #endif
static char t char * s
Definition: cfortran.h:568
cstr elem_len * i
Definition: cfortran.h:603
A Class for the derivatives of shape functions The class design is essentially the same as Shape,...
Definition: shape.h:278
FaceElements are elements that coincide with the faces of higher-dimensional "bulk" elements....
Definition: elements.h:4342
int & face_index()
Index of the face (a number that uniquely identifies the face in the element)
Definition: elements.h:4630
void outer_unit_normal(const Vector< double > &s, Vector< double > &unit_normal) const
Compute outer unit normal at the specified local coordinate.
Definition: elements.cc:6036
double zeta_nodal(const unsigned &n, const unsigned &k, const unsigned &i) const
In a FaceElement, the "global" intrinsic coordinate of the element along the boundary,...
Definition: elements.h:4501
Vector< double > local_coordinate_in_bulk(const Vector< double > &s) const
Return vector of local coordinates in bulk element, given the local coordinates in this FaceElement.
Definition: elements.cc:6383
FiniteElement *& bulk_element_pt()
Pointer to higher-dimensional "bulk" element.
Definition: elements.h:4739
double interpolated_x(const Vector< double > &s, const unsigned &i) const
Return FE interpolated coordinate x[i] at local coordinate s. Overloaded to get information from bulk...
Definition: elements.h:4532
double J_eulerian(const Vector< double > &s) const
Return the Jacobian of mapping from local to global coordinates at local position s....
Definition: elements.cc:5272
double J_eulerian_at_knot(const unsigned &ipt) const
Return the Jacobian of the mapping from local to global coordinates at the ipt-th integration point O...
Definition: elements.cc:5358
//////////////////////////////////////////////////////////////////// ////////////////////////////////...
Definition: elements.h:5002
A general Finite Element class.
Definition: elements.h:1317
Node *& node_pt(const unsigned &n)
Return a pointer to the local node n.
Definition: elements.h:2179
double nodal_value(const unsigned &n, const unsigned &i) const
Return the i-th value stored at local node n. Produces suitably interpolated values for hanging nodes...
Definition: elements.h:2597
virtual void output(std::ostream &outfile)
Output the element data — typically the values at the nodes in a format suitable for post-processing.
Definition: elements.h:3054
virtual void shape(const Vector< double > &s, Shape &psi) const =0
Calculate the geometric shape functions at local coordinate s. This function must be overloaded for e...
unsigned nnode() const
Return the number of nodes.
Definition: elements.h:2214
Integral *const & integral_pt() const
Return the pointer to the integration scheme (const version)
Definition: elements.h:1967
virtual void dshape_local(const Vector< double > &s, Shape &psi, DShape &dpsids) const
Function to compute the geometric shape functions and derivatives w.r.t. local coordinates at local c...
Definition: elements.h:1985
virtual void build_face_element(const int &face_index, FaceElement *face_element_pt)
Function for building a lower dimensional FaceElement on the specified face of the FiniteElement....
Definition: elements.cc:5162
virtual void shape_at_knot(const unsigned &ipt, Shape &psi) const
Return the geometric shape function at the ipt-th integration point.
Definition: elements.cc:3250
static DenseMatrix< double > Dummy_matrix
Empty dense matrix used as a dummy argument to combined residual and jacobian functions in the case w...
Definition: elements.h:227
virtual double knot(const unsigned &i, const unsigned &j) const =0
Return local coordinate s[j] of i-th integration point.
virtual unsigned nweight() const =0
Return the number of integration points of the scheme.
virtual double weight(const unsigned &i) const =0
Return weight of i-th integration point.
unsigned ndim() const
Return (Eulerian) spatial dimension of the node.
Definition: nodes.h:1054
An OomphLibError object which should be thrown when an run-time error is encountered....
///////////////////////////////////////////////////////////////////// ///////////////////////////////...
virtual std::complex< unsigned > u_index_pml_fourier_decomposed_helmholtz() const
Broken assignment operator.
A class for elements that allow the imposition of an applied flux on the boundaries of Fourier decomp...
void output(FILE *file_pt, const unsigned &n_plot)
C-style output function – forward to broken version in FiniteElement until somebody decides what exac...
virtual std::complex< unsigned > u_index_pml_fourier_decomposed_helmholtz() const
Return the index at which the unknown value is stored. (Real/imag part gives real index of real/imag ...
void(* PMLFourierDecomposedHelmholtzPrescribedFluxFctPt)(const Vector< double > &x, std::complex< double > &flux)
Function pointer to the prescribed-flux function fct(x,f(x)) – x is a Vector and the flux is a comple...
void fill_in_contribution_to_residuals(Vector< double > &residuals)
Add the element's contribution to its residual vector.
void output(FILE *file_pt)
C-style output function – forward to broken version in FiniteElement until somebody decides what exac...
std::complex< unsigned > U_index_pml_fourier_decomposed_helmholtz
The index at which the real and imag part of the unknown is stored at the nodes.
double shape_and_test_at_knot(const unsigned &ipt, Shape &psi, Shape &test) const
Function to compute the shape and test functions and to return the Jacobian of mapping between local ...
PMLFourierDecomposedHelmholtzFluxElement(const PMLFourierDecomposedHelmholtzFluxElement &dummy)=delete
Broken copy constructor.
virtual void fill_in_generic_residual_contribution_pml_fourier_decomposed_helmholtz_flux(Vector< double > &residuals, DenseMatrix< double > &jacobian, const unsigned &flag)
Add the element's contribution to its residual vector. flag=1(or 0): do (or don't) compute the contri...
void output(std::ostream &outfile)
Output function – forward to broken version in FiniteElement until somebody decides what exactly they...
double shape_and_test(const Vector< double > &s, Shape &psi, Shape &test) const
Function to compute the shape and test functions and to return the Jacobian of mapping between local ...
void fill_in_contribution_to_jacobian(Vector< double > &residuals, DenseMatrix< double > &jacobian)
Add the element's contribution to its residual vector and its Jacobian matrix.
void get_flux(const Vector< double > &x, std::complex< double > &flux)
Function to calculate the prescribed flux at a given spatial position.
PMLFourierDecomposedHelmholtzPrescribedFluxFctPt Flux_fct_pt
Function pointer to the (global) prescribed-flux function.
PMLFourierDecomposedHelmholtzPrescribedFluxFctPt & flux_fct_pt()
Broken assignment operator.
void output(std::ostream &outfile, const unsigned &n_plot)
Output function – forward to broken version in FiniteElement until somebody decides what exactly they...
////////////////////////////////////////////////////////////////// //////////////////////////////////...
std::complex< unsigned > U_index_pml_fourier_decomposed_helmholtz
The index at which the real and imag part of the unknown is stored at the nodes.
void output(FILE *file_pt)
C-style output function – forward to broken version in FiniteElement until somebody decides what exac...
double d_shape_and_test_local(const Vector< double > &s, Shape &psi, Shape &test, DShape &dpsi_ds, DShape &dtest_ds) const
Function to compute the shape, test functions and derivates and to return the Jacobian of mapping bet...
double zeta_nodal(const unsigned &n, const unsigned &k, const unsigned &i) const
Broken assignment operator.
void output(FILE *file_pt, const unsigned &n_plot)
C-style output function – forward to broken version in FiniteElement until somebody decides what exac...
void output(std::ostream &outfile, const unsigned &n_plot)
Output function – forward to broken version in FiniteElement until somebody decides what exactly they...
double global_power_contribution(std::ofstream &outfile)
Compute the element's contribution to the time-averaged radiated power over the artificial boundary....
double shape_and_test(const Vector< double > &s, Shape &psi, Shape &test) const
Function to compute the test functions and to return the Jacobian of mapping between local and global...
virtual std::complex< unsigned > u_index_pml_fourier_decomposed_helmholtz() const
Return the index at which the real/imag unknown value is stored.
void output(std::ostream &outfile)
Output function – forward to broken version in FiniteElement until somebody decides what exactly they...
double global_power_contribution()
Compute the element's contribution to the time-averaged radiated power over the artificial boundary....
PMLFourierDecomposedHelmholtzPowerMonitorElement(const PMLFourierDecomposedHelmholtzPowerMonitorElement &dummy)=delete
Broken copy constructor.
A Class for shape functions. In simple cases, the shape functions have only one index that can be tho...
Definition: shape.h:76
std::string string(const unsigned &i)
Return the i-th string or "" if the relevant string hasn't been defined.
const double Pi
50 digits from maple
//////////////////////////////////////////////////////////////////// ////////////////////////////////...