pml_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 Pml Helmholtz equations
28 #ifndef OOMPH_PML_HELMHOLTZ_FLUX_ELEMENTS_HEADER
29 #define OOMPH_PML_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 // oomph-lib ncludes
38 #include "../generic/Qelements.h"
39 #include "math.h"
40 #include <complex>
41 
42 namespace oomph
43 {
44  //======================================================================
45  /// A class for elements that allow the post-processing of
46  /// radiated power and flux on the boundaries of PMLHelmholtz elements.
47  /// The element geometry is obtained from the FaceGeometry<ELEMENT>
48  /// policy class.
49  //======================================================================
50  template<class ELEMENT>
51  class PMLHelmholtzPowerElement : public virtual FaceGeometry<ELEMENT>,
52  public virtual FaceElement
53  {
54  public:
55  /// Constructor, takes the pointer to the "bulk" element and the
56  /// index of the face to which the element is attached.
57  PMLHelmholtzPowerElement(FiniteElement* const& bulk_el_pt,
58  const int& face_index);
59 
60  /// Broken empty constructor
62  {
63  throw OomphLibError(
64  "Don't call empty constructor for PMLHelmholtzPowerElement",
65  OOMPH_CURRENT_FUNCTION,
66  OOMPH_EXCEPTION_LOCATION);
67  }
68 
69  /// Broken copy constructor
71 
72  /// Broken assignment operator
73  // Commented out broken assignment operator because this can lead to a
74  // conflict warning when used in the virtual inheritence hierarchy.
75  // Essentially the compiler doesn't realise that two separate
76  // implementations of the broken function are the same and so, quite
77  // rightly, it shouts.
78  /*void operator=(const PMLHelmholtzPowerElement&) = delete;*/
79 
80 
81  /// Specify the value of nodal zeta from the face geometry
82  /// The "global" intrinsic coordinate of the element when
83  /// viewed as part of a geometric object should be given by
84  /// the FaceElement representation, by default (needed to break
85  /// indeterminacy if bulk element is SolidElement)
86  double zeta_nodal(const unsigned& n,
87  const unsigned& k,
88  const unsigned& i) const
89  {
90  return FaceElement::zeta_nodal(n, k, i);
91  }
92 
93 
94  /// Return the index at which the unknown value
95  /// is stored.
96  virtual inline std::complex<unsigned> u_index_helmholtz() const
97  {
98  return std::complex<unsigned>(U_index_helmholtz.real(),
99  U_index_helmholtz.imag());
100  }
101 
102 
103  /// Compute the element's contribution to the time-averaged
104  /// radiated power over the artificial boundary
106  {
107  // Dummy output file
108  std::ofstream outfile;
109  return global_power_contribution(outfile);
110  }
111 
112  /// Compute the element's contribution to the time-averaged
113  /// radiated power over the artificial boundary. Also output the
114  /// power density in the specified
115  /// output file if it's open.
116  double global_power_contribution(std::ofstream& outfile)
117  {
118  // pointer to the corresponding bulk element
119  ELEMENT* bulk_elem_pt = dynamic_cast<ELEMENT*>(this->bulk_element_pt());
120 
121  // Number of nodes in bulk element
122  unsigned nnode_bulk = bulk_elem_pt->nnode();
123  const unsigned n_node_local = nnode();
124 
125  // get the dim of the bulk and local nodes
126  const unsigned bulk_dim = bulk_elem_pt->dim();
127  const unsigned local_dim = this->dim();
128 
129  // Set up memory for the shape and test functions
130  Shape psi(n_node_local);
131 
132  // Set up memory for the shape functions
133  Shape psi_bulk(nnode_bulk);
134  DShape dpsi_bulk_dx(nnode_bulk, bulk_dim);
135 
136  // Set up memory for the outer unit normal
137  Vector<double> unit_normal(bulk_dim);
138 
139  // Set the value of n_intpt
140  const unsigned n_intpt = integral_pt()->nweight();
141 
142  // Set the Vector to hold local coordinates
143  Vector<double> s(local_dim);
144  double power = 0.0;
145 
146  // Output?
147  if (outfile.is_open())
148  {
149  outfile << "ZONE\n";
150  }
151 
152  // Loop over the integration points
153  //--------------------------------
154  for (unsigned ipt = 0; ipt < n_intpt; ipt++)
155  {
156  // Assign values of s
157  for (unsigned i = 0; i < local_dim; i++)
158  {
159  s[i] = integral_pt()->knot(ipt, i);
160  }
161  // get the outer_unit_ext vector
162  this->outer_unit_normal(s, unit_normal);
163 
164  // Get the integral weight
165  double w = integral_pt()->weight(ipt);
166 
167  // Get jacobian of mapping
168  double J = J_eulerian(s);
169 
170  // Premultiply the weights and the Jacobian
171  double W = w * J;
172 
173  // Get local coordinates in bulk element by copy construction
175 
176  // Call the derivatives of the shape functions
177  // in the bulk -- must do this via s because this point
178  // is not an integration point the bulk element!
179  (void)bulk_elem_pt->dshape_eulerian(s_bulk, psi_bulk, dpsi_bulk_dx);
180  this->shape(s, psi);
181 
182  // Derivs of Eulerian coordinates w.r.t. local coordinates
183  std::complex<double> dphi_dn(0.0, 0.0);
184  Vector<std::complex<double>> interpolated_dphidx(
185  bulk_dim, std::complex<double>(0.0, 0.0));
186  std::complex<double> interpolated_phi(0.0, 0.0);
187  Vector<double> x(bulk_dim);
188 
189  // Calculate function value and derivatives:
190  //-----------------------------------------
191  // Loop over nodes
192  for (unsigned l = 0; l < nnode_bulk; l++)
193  {
194  // Get the nodal value of the helmholtz unknown
195  const std::complex<double> phi_value(
196  bulk_elem_pt->nodal_value(l,
197  bulk_elem_pt->u_index_helmholtz().real()),
198  bulk_elem_pt->nodal_value(
199  l, bulk_elem_pt->u_index_helmholtz().imag()));
200 
201  // Loop over directions
202  for (unsigned i = 0; i < bulk_dim; i++)
203  {
204  interpolated_dphidx[i] += phi_value * dpsi_bulk_dx(l, i);
205  }
206  } // End of loop over the bulk_nodes
207 
208  for (unsigned l = 0; l < n_node_local; l++)
209  {
210  // Get the nodal value of the helmholtz unknown
211  const std::complex<double> phi_value(
212  nodal_value(l, u_index_helmholtz().real()),
213  nodal_value(l, u_index_helmholtz().imag()));
214 
215  interpolated_phi += phi_value * psi(l);
216  }
217 
218  // define dphi_dn
219  for (unsigned i = 0; i < bulk_dim; i++)
220  {
221  dphi_dn += interpolated_dphidx[i] * unit_normal[i];
222  }
223 
224  // Power density
225  double integrand = 0.5 * (interpolated_phi.real() * dphi_dn.imag() -
226  interpolated_phi.imag() * dphi_dn.real());
227 
228  // Output?
229  if (outfile.is_open())
230  {
231  interpolated_x(s, x);
232  double phi = atan2(x[1], x[0]);
233  outfile << x[0] << " " << x[1] << " " << phi << " " << integrand
234  << "\n";
235  }
236 
237  // ...add to integral
238  power += integrand * W;
239  }
240 
241  return power;
242  }
243 
244 
245  /// Compute the element's contribution to the time-averaged
246  /// flux
247  std::complex<double> global_flux_contribution()
248  {
249  // Dummy output file
250  std::ofstream outfile;
251  return global_flux_contribution(outfile);
252  }
253 
254  /// Compute the element's contribution to the integral of
255  /// the flux over the artificial boundary. Also output the
256  /// flux in the specified
257  /// output file if it's open.
258  std::complex<double> global_flux_contribution(std::ofstream& outfile)
259  {
260  // pointer to the corresponding bulk element
261  ELEMENT* bulk_elem_pt = dynamic_cast<ELEMENT*>(this->bulk_element_pt());
262 
263  // Number of nodes in bulk element
264  unsigned nnode_bulk = bulk_elem_pt->nnode();
265  const unsigned n_node_local = nnode();
266 
267  // get the dim of the bulk and local nodes
268  const unsigned bulk_dim = bulk_elem_pt->dim();
269  const unsigned local_dim = this->dim();
270 
271  // Set up memory for the shape and test functions
272  Shape psi(n_node_local);
273 
274  // Set up memory for the shape functions
275  Shape psi_bulk(nnode_bulk);
276  DShape dpsi_bulk_dx(nnode_bulk, bulk_dim);
277 
278  // Set up memory for the outer unit normal
279  Vector<double> unit_normal(bulk_dim);
280 
281  // Set the value of n_intpt
282  const unsigned n_intpt = integral_pt()->nweight();
283 
284  // Set the Vector to hold local coordinates
285  Vector<double> s(local_dim);
286  std::complex<double> flux(0.0, 0.0);
287 
288  // Output?
289  if (outfile.is_open())
290  {
291  outfile << "ZONE\n";
292  }
293 
294  // Loop over the integration points
295  //--------------------------------
296  for (unsigned ipt = 0; ipt < n_intpt; ipt++)
297  {
298  // Assign values of s
299  for (unsigned i = 0; i < local_dim; i++)
300  {
301  s[i] = integral_pt()->knot(ipt, i);
302  }
303  // get the outer_unit_ext vector
304  this->outer_unit_normal(s, unit_normal);
305 
306  // Get the integral weight
307  double w = integral_pt()->weight(ipt);
308 
309  // Get jacobian of mapping
310  double J = J_eulerian(s);
311 
312  // Premultiply the weights and the Jacobian
313  double W = w * J;
314 
315  // Get local coordinates in bulk element by copy construction
317 
318  // Call the derivatives of the shape functions
319  // in the bulk -- must do this via s because this point
320  // is not an integration point the bulk element!
321  (void)bulk_elem_pt->dshape_eulerian(s_bulk, psi_bulk, dpsi_bulk_dx);
322  this->shape(s, psi);
323 
324  // Derivs of Eulerian coordinates w.r.t. local coordinates
325  std::complex<double> dphi_dn(0.0, 0.0);
326  Vector<std::complex<double>> interpolated_dphidx(
327  bulk_dim, std::complex<double>(0.0, 0.0));
328  Vector<double> x(bulk_dim);
329 
330  // Calculate function value and derivatives:
331  //-----------------------------------------
332  // Loop over nodes
333  for (unsigned l = 0; l < nnode_bulk; l++)
334  {
335  // Get the nodal value of the helmholtz unknown
336  const std::complex<double> phi_value(
337  bulk_elem_pt->nodal_value(l,
338  bulk_elem_pt->u_index_helmholtz().real()),
339  bulk_elem_pt->nodal_value(
340  l, bulk_elem_pt->u_index_helmholtz().imag()));
341 
342  // Loop over directions
343  for (unsigned i = 0; i < bulk_dim; i++)
344  {
345  interpolated_dphidx[i] += phi_value * dpsi_bulk_dx(l, i);
346  }
347  } // End of loop over the bulk_nodes
348 
349 
350  // define dphi_dn
351  for (unsigned i = 0; i < bulk_dim; i++)
352  {
353  dphi_dn += interpolated_dphidx[i] * unit_normal[i];
354  }
355 
356  // Output?
357  if (outfile.is_open())
358  {
359  interpolated_x(s, x);
360  outfile << x[0] << " " << x[1] << " " << dphi_dn.real() << " "
361  << dphi_dn.imag() << "\n";
362  }
363 
364  // ...add to integral
365  flux += dphi_dn * W;
366  }
367 
368  return flux;
369  }
370 
371 
372  protected:
373  /// The index at which the real and imag part of the unknown is
374  /// stored at the nodes
375  std::complex<unsigned> U_index_helmholtz;
376 
377  /// The spatial dimension of the problem
378  unsigned Dim;
379  };
380 
381  /// ///////////////////////////////////////////////////////////////////
382  /// ///////////////////////////////////////////////////////////////////
383  /// ///////////////////////////////////////////////////////////////////
384 
385 
386  //===========================================================================
387  /// Constructor, takes the pointer to the "bulk" element, the
388  /// index of the fixed local coordinate and its value represented
389  /// by an integer (+/- 1), indicating that the face is located
390  /// at the max. or min. value of the "fixed" local coordinate
391  /// in the bulk element.
392  //===========================================================================
393  template<class ELEMENT>
395  FiniteElement* const& bulk_el_pt, const int& face_index)
396  : FaceGeometry<ELEMENT>(), FaceElement()
397  {
398 #ifdef PARANOID
399  {
400  // Check that the element is not a refineable 3d element
401  ELEMENT* elem_pt = dynamic_cast<ELEMENT*>(bulk_el_pt);
402  // If it's three-d
403  if (elem_pt->dim() == 3)
404  {
405  // Is it refineable
406  RefineableElement* ref_el_pt =
407  dynamic_cast<RefineableElement*>(elem_pt);
408  if (ref_el_pt != 0)
409  {
410  if (this->has_hanging_nodes())
411  {
412  throw OomphLibError("This flux element will not work correctly if "
413  "nodes are hanging\n",
414  OOMPH_CURRENT_FUNCTION,
415  OOMPH_EXCEPTION_LOCATION);
416  }
417  }
418  }
419  }
420 #endif
421 
422  // Let the bulk element build the FaceElement, i.e. setup the pointers
423  // to its nodes (by referring to the appropriate nodes in the bulk
424  // element), etc.
425  bulk_el_pt->build_face_element(face_index, this);
426 
427  // Extract the dimension of the problem from the dimension of
428  // the first node
429  Dim = this->node_pt(0)->ndim();
430 
431  // Set up U_index_helmholtz. Initialise to zero, which probably won't change
432  // in most cases, oh well, the price we pay for generality
433  U_index_helmholtz = std::complex<unsigned>(0, 1);
434 
435  // Cast to the appropriate PMLHelmholtzEquation so that we can
436  // find the index at which the variable is stored
437  // We assume that the dimension of the full problem is the same
438  // as the dimension of the node, if this is not the case you will have
439  // to write custom elements, sorry
440  switch (Dim)
441  {
442  // One dimensional problem
443  case 1:
444  {
445  PMLHelmholtzEquations<1>* eqn_pt =
446  dynamic_cast<PMLHelmholtzEquations<1>*>(bulk_el_pt);
447  // If the cast has failed die
448  if (eqn_pt == 0)
449  {
450  std::string error_string =
451  "Bulk element must inherit from PMLHelmholtzEquations.";
452  error_string +=
453  "Nodes are one dimensional, but cannot cast the bulk element to\n";
454  error_string += "PMLHelmholtzEquations<1>\n.";
455  error_string += "If you desire this functionality, you must "
456  "implement it yourself\n";
457 
458  throw OomphLibError(
459  error_string, OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
460  }
461  // Otherwise read out the value
462  else
463  {
464  // Read the index from the (cast) bulk element
466  }
467  }
468  break;
469 
470  // Two dimensional problem
471  case 2:
472  {
473  PMLHelmholtzEquations<2>* eqn_pt =
474  dynamic_cast<PMLHelmholtzEquations<2>*>(bulk_el_pt);
475  // If the cast has failed die
476  if (eqn_pt == 0)
477  {
478  std::string error_string =
479  "Bulk element must inherit from PMLHelmholtzEquations.";
480  error_string +=
481  "Nodes are two dimensional, but cannot cast the bulk element to\n";
482  error_string += "PMLHelmholtzEquations<2>\n.";
483  error_string += "If you desire this functionality, you must "
484  "implement it yourself\n";
485 
486  throw OomphLibError(
487  error_string, OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
488  }
489  else
490  {
491  // Read the index from the (cast) bulk element
493  }
494  }
495 
496  break;
497 
498  // Three dimensional problem
499  case 3:
500  {
501  PMLHelmholtzEquations<3>* eqn_pt =
502  dynamic_cast<PMLHelmholtzEquations<3>*>(bulk_el_pt);
503  // If the cast has failed die
504  if (eqn_pt == 0)
505  {
506  std::string error_string =
507  "Bulk element must inherit from PMLHelmholtzEquations.";
508  error_string += "Nodes are three dimensional, but cannot cast the "
509  "bulk element to\n";
510  error_string += "PMLHelmholtzEquations<3>\n.";
511  error_string += "If you desire this functionality, you must "
512  "implement it yourself\n";
513 
514  throw OomphLibError(
515  error_string, OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
516  }
517  else
518  {
519  // Read the index from the (cast) bulk element
521  }
522  }
523  break;
524 
525  // Any other case is an error
526  default:
527  std::ostringstream error_stream;
528  error_stream << "Dimension of node is " << Dim
529  << ". It should be 1,2, or 3!" << std::endl;
530 
531  throw OomphLibError(
532  error_stream.str(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
533  break;
534  }
535  }
536 
537 
538  /// ///////////////////////////////////////////////////////////////////////
539  /// ///////////////////////////////////////////////////////////////////////
540  /// ///////////////////////////////////////////////////////////////////////
541 
542 
543  //======================================================================
544  /// A class for elements that allow the imposition of an
545  /// applied flux on the boundaries of PMLHelmholtz elements.
546  /// The element geometry is obtained from the FaceGeometry<ELEMENT>
547  /// policy class.
548  //======================================================================
549  template<class ELEMENT>
550  class PMLHelmholtzFluxElement : public virtual FaceGeometry<ELEMENT>,
551  public virtual FaceElement
552  {
553  public:
554  /// Function pointer to the prescribed-flux function fct(x,f(x)) --
555  /// x is a Vector and the flux is a complex. NOTE THAT f(x) represents
556  /// c^2 du/dn!
558  std::complex<double>& flux);
559 
560  /// Constructor, takes the pointer to the "bulk" element and the
561  /// index of the face to which the element is attached.
562  PMLHelmholtzFluxElement(FiniteElement* const& bulk_el_pt,
563  const int& face_index);
564 
565  /// Broken empty constructor
567  {
568  throw OomphLibError(
569  "Don't call empty constructor for PMLHelmholtzFluxElement",
570  OOMPH_CURRENT_FUNCTION,
571  OOMPH_EXCEPTION_LOCATION);
572  }
573 
574  /// Broken copy constructor
576 
577  /// Broken assignment operator
578  /*void operator=(const PMLHelmholtzFluxElement&) = delete;*/
579 
580 
581  /// Access function for the prescribed-flux function pointer
583  {
584  return Flux_fct_pt;
585  }
586 
587 
588  /// Add the element's contribution to its residual vector
590  {
591  // Call the generic residuals function with flag set to 0
592  // using a dummy matrix argument
594  residuals, GeneralisedElement::Dummy_matrix, 0);
595  }
596 
597 
598  /// Add the element's contribution to its residual vector and its
599  /// Jacobian matrix
601  DenseMatrix<double>& jacobian)
602  {
603  // Call the generic routine with the flag set to 1
605  residuals, jacobian, 1);
606  }
607 
608 
609  /// Specify the value of nodal zeta from the face geometry
610  /// The "global" intrinsic coordinate of the element when
611  /// viewed as part of a geometric object should be given by
612  /// the FaceElement representation, by default (needed to break
613  /// indeterminacy if bulk element is SolidElement)
614  double zeta_nodal(const unsigned& n,
615  const unsigned& k,
616  const unsigned& i) const
617  {
618  return FaceElement::zeta_nodal(n, k, i);
619  }
620 
621 
622  /// Output function -- forward to broken version in FiniteElement
623  /// until somebody decides what exactly they want to plot here...
624  void output(std::ostream& outfile)
625  {
626  FiniteElement::output(outfile);
627  }
628 
629  /// Output function -- forward to broken version in FiniteElement
630  /// until somebody decides what exactly they want to plot here...
631  void output(std::ostream& outfile, const unsigned& n_plot)
632  {
633  FiniteElement::output(outfile, n_plot);
634  }
635 
636 
637  /// C-style output function -- forward to broken version in FiniteElement
638  /// until somebody decides what exactly they want to plot here...
639  void output(FILE* file_pt)
640  {
641  FiniteElement::output(file_pt);
642  }
643 
644  /// C-style output function -- forward to broken version in
645  /// FiniteElement until somebody decides what exactly they want to plot
646  /// here...
647  void output(FILE* file_pt, const unsigned& n_plot)
648  {
649  FiniteElement::output(file_pt, n_plot);
650  }
651 
652 
653  /// Return the index at which the unknown value
654  /// is stored.
655  virtual inline std::complex<unsigned> u_index_helmholtz() const
656  {
657  return std::complex<unsigned>(U_index_helmholtz.real(),
658  U_index_helmholtz.imag());
659  }
660 
661 
662  /// The number of "DOF types" that degrees of freedom in this element
663  /// are sub-divided into: real and imaginary part
664  unsigned ndof_types() const
665  {
666  return 2;
667  }
668 
669  /// Create a list of pairs for all unknowns in this element,
670  /// so that the first entry in each pair contains the global equation
671  /// number of the unknown, while the second one contains the number
672  /// of the "DOF type" that this unknown is associated with.
673  /// (Function can obviously only be called if the equation numbering
674  /// scheme has been set up.) Real=0; Imag=1
676  std::list<std::pair<unsigned long, unsigned>>& dof_lookup_list) const
677  {
678  // temporary pair (used to store dof lookup prior to being added to list)
679  std::pair<unsigned, unsigned> dof_lookup;
680 
681  // number of nodes
682  unsigned n_node = this->nnode();
683 
684  // loop over the nodes
685  for (unsigned n = 0; n < n_node; n++)
686  {
687  // determine local eqn number for real part
688  int local_eqn_number =
689  this->nodal_local_eqn(n, this->U_index_helmholtz.real());
690 
691  // ignore pinned values
692  if (local_eqn_number >= 0)
693  {
694  // store dof lookup in temporary pair: First entry in pair
695  // is global equation number; second entry is dof type
696  dof_lookup.first = this->eqn_number(local_eqn_number);
697  dof_lookup.second = 0;
698 
699  // add to list
700  dof_lookup_list.push_front(dof_lookup);
701  }
702 
703  // determine local eqn number for imag part
705  this->nodal_local_eqn(n, this->U_index_helmholtz.imag());
706 
707  // ignore pinned values
708  if (local_eqn_number >= 0)
709  {
710  // store dof lookup in temporary pair: First entry in pair
711  // is global equation number; second entry is dof type
712  dof_lookup.first = this->eqn_number(local_eqn_number);
713  dof_lookup.second = 1;
714 
715  // add to list
716  dof_lookup_list.push_front(dof_lookup);
717  }
718  }
719  }
720 
721  protected:
722  /// Function to compute the shape and test functions and to return
723  /// the Jacobian of mapping between local and global (Eulerian)
724  /// coordinates
725  inline double shape_and_test(const Vector<double>& s,
726  Shape& psi,
727  Shape& test) const
728  {
729  // Find number of nodes
730  unsigned n_node = nnode();
731 
732  // Get the shape functions
733  shape(s, psi);
734 
735  // Set the test functions to be the same as the shape functions
736  for (unsigned i = 0; i < n_node; i++)
737  {
738  test[i] = psi[i];
739  }
740 
741  // Return the value of the jacobian
742  return J_eulerian(s);
743  }
744 
745 
746  /// Function to compute the shape and test functions and to return
747  /// the Jacobian of mapping between local and global (Eulerian)
748  /// coordinates
749  inline double shape_and_test_at_knot(const unsigned& ipt,
750  Shape& psi,
751  Shape& test) const
752  {
753  // Find number of nodes
754  unsigned n_node = nnode();
755 
756  // Get the shape functions
757  shape_at_knot(ipt, psi);
758 
759  // Set the test functions to be the same as the shape functions
760  for (unsigned i = 0; i < n_node; i++)
761  {
762  test[i] = psi[i];
763  }
764 
765  // Return the value of the jacobian
766  return J_eulerian_at_knot(ipt);
767  }
768 
769 
770  /// Function to calculate the prescribed flux at a given spatial
771  /// position
772  void get_flux(const Vector<double>& x, std::complex<double>& flux)
773  {
774  // If the function pointer is zero return zero
775  if (Flux_fct_pt == 0)
776  {
777  flux = std::complex<double>(0.0, 0.0);
778  }
779  // Otherwise call the function
780  else
781  {
782  (*Flux_fct_pt)(x, flux);
783  }
784  }
785 
786 
787  /// The index at which the real and imag part of the unknown is
788  /// stored at the nodes
789  std::complex<unsigned> U_index_helmholtz;
790 
791 
792  /// Add the element's contribution to its residual vector.
793  /// flag=1(or 0): do (or don't) compute the contribution to the
794  /// Jacobian as well.
796  Vector<double>& residuals,
797  DenseMatrix<double>& jacobian,
798  const unsigned& flag);
799 
800 
801  /// Function pointer to the (global) prescribed-flux function
803 
804  /// The spatial dimension of the problem
805  unsigned Dim;
806  };
807 
808  /// ///////////////////////////////////////////////////////////////////
809  /// ///////////////////////////////////////////////////////////////////
810  /// ///////////////////////////////////////////////////////////////////
811 
812 
813  //===========================================================================
814  /// Constructor, takes the pointer to the "bulk" element, the
815  /// index of the fixed local coordinate and its value represented
816  /// by an integer (+/- 1), indicating that the face is located
817  /// at the max. or min. value of the "fixed" local coordinate
818  /// in the bulk element.
819  //===========================================================================
820  template<class ELEMENT>
822  FiniteElement* const& bulk_el_pt, const int& face_index)
823  : FaceGeometry<ELEMENT>(), FaceElement()
824  {
825 #ifdef PARANOID
826  {
827  // Check that the element is not a refineable 3d element
828  ELEMENT* elem_pt = dynamic_cast<ELEMENT*>(bulk_el_pt);
829  // If it's three-d
830  if (elem_pt->dim() == 3)
831  {
832  // Is it refineable
833  RefineableElement* ref_el_pt =
834  dynamic_cast<RefineableElement*>(elem_pt);
835  if (ref_el_pt != 0)
836  {
837  if (this->has_hanging_nodes())
838  {
839  throw OomphLibError("This flux element will not work correctly if "
840  "nodes are hanging\n",
841  OOMPH_CURRENT_FUNCTION,
842  OOMPH_EXCEPTION_LOCATION);
843  }
844  }
845  }
846  }
847 #endif
848 
849  // Let the bulk element build the FaceElement, i.e. setup the pointers
850  // to its nodes (by referring to the appropriate nodes in the bulk
851  // element), etc.
852  bulk_el_pt->build_face_element(face_index, this);
853 
854  // Initialise the prescribed-flux function pointer to zero
855  Flux_fct_pt = 0;
856 
857  // Extract the dimension of the problem from the dimension of
858  // the first node
859  Dim = this->node_pt(0)->ndim();
860 
861  // Set up U_index_helmholtz. Initialise to zero, which probably won't change
862  // in most cases, oh well, the price we pay for generality
863  U_index_helmholtz = std::complex<unsigned>(0, 1);
864 
865  // Cast to the appropriate PMLHelmholtzEquation so that we can
866  // find the index at which the variable is stored
867  // We assume that the dimension of the full problem is the same
868  // as the dimension of the node, if this is not the case you will have
869  // to write custom elements, sorry
870  switch (Dim)
871  {
872  // One dimensional problem
873  case 1:
874  {
875  PMLHelmholtzEquations<1>* eqn_pt =
876  dynamic_cast<PMLHelmholtzEquations<1>*>(bulk_el_pt);
877  // If the cast has failed die
878  if (eqn_pt == 0)
879  {
880  std::string error_string =
881  "Bulk element must inherit from PMLHelmholtzEquations.";
882  error_string +=
883  "Nodes are one dimensional, but cannot cast the bulk element to\n";
884  error_string += "PMLHelmholtzEquations<1>\n.";
885  error_string += "If you desire this functionality, you must "
886  "implement it yourself\n";
887 
888  throw OomphLibError(
889  error_string, OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
890  }
891  // Otherwise read out the value
892  else
893  {
894  // Read the index from the (cast) bulk element
896  }
897  }
898  break;
899 
900  // Two dimensional problem
901  case 2:
902  {
903  PMLHelmholtzEquations<2>* eqn_pt =
904  dynamic_cast<PMLHelmholtzEquations<2>*>(bulk_el_pt);
905  // If the cast has failed die
906  if (eqn_pt == 0)
907  {
908  std::string error_string =
909  "Bulk element must inherit from PMLHelmholtzEquations.";
910  error_string +=
911  "Nodes are two dimensional, but cannot cast the bulk element to\n";
912  error_string += "PMLHelmholtzEquations<2>\n.";
913  error_string += "If you desire this functionality, you must "
914  "implement it yourself\n";
915 
916  throw OomphLibError(
917  error_string, OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
918  }
919  else
920  {
921  // Read the index from the (cast) bulk element
923  }
924  }
925 
926  break;
927 
928  // Three dimensional problem
929  case 3:
930  {
931  PMLHelmholtzEquations<3>* eqn_pt =
932  dynamic_cast<PMLHelmholtzEquations<3>*>(bulk_el_pt);
933  // If the cast has failed die
934  if (eqn_pt == 0)
935  {
936  std::string error_string =
937  "Bulk element must inherit from PMLHelmholtzEquations.";
938  error_string += "Nodes are three dimensional, but cannot cast the "
939  "bulk element to\n";
940  error_string += "PMLHelmholtzEquations<3>\n.";
941  error_string += "If you desire this functionality, you must "
942  "implement it yourself\n";
943 
944  throw OomphLibError(
945  error_string, OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
946  }
947  else
948  {
949  // Read the index from the (cast) bulk element
951  }
952  }
953  break;
954 
955  // Any other case is an error
956  default:
957  std::ostringstream error_stream;
958  error_stream << "Dimension of node is " << Dim
959  << ". It should be 1,2, or 3!" << std::endl;
960 
961  throw OomphLibError(
962  error_stream.str(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
963  break;
964  }
965  }
966 
967 
968  //===========================================================================
969  /// Compute the element's residual vector and the (zero) Jacobian matrix.
970  //===========================================================================
971  template<class ELEMENT>
974  Vector<double>& residuals,
975  DenseMatrix<double>& jacobian,
976  const unsigned& flag)
977  {
978  // Find out how many nodes there are
979  const unsigned n_node = nnode();
980 
981  // Set up memory for the shape and test functions
982  Shape psif(n_node), testf(n_node);
983 
984  // Set the value of Nintpt
985  const unsigned n_intpt = integral_pt()->nweight();
986 
987  // Set the Vector to hold local coordinates
988  Vector<double> s(Dim - 1);
989 
990  // Integers to hold the local equation and unknown numbers
991  int local_eqn_real = 0, local_eqn_imag = 0;
992 
993  // Loop over the integration points
994  //--------------------------------
995  for (unsigned ipt = 0; ipt < n_intpt; ipt++)
996  {
997  // Assign values of s
998  for (unsigned i = 0; i < (Dim - 1); i++)
999  {
1000  s[i] = integral_pt()->knot(ipt, i);
1001  }
1002 
1003  // Get the integral weight
1004  double w = integral_pt()->weight(ipt);
1005 
1006  // Find the shape and test functions and return the Jacobian
1007  // of the mapping
1008  double J = shape_and_test(s, psif, testf);
1009 
1010  // Premultiply the weights and the Jacobian
1011  double W = w * J;
1012 
1013  // Need to find position to feed into flux function, initialise to zero
1014  Vector<double> interpolated_x(Dim, 0.0);
1015 
1016  // Calculate Eulerian position of integration point
1017  for (unsigned l = 0; l < n_node; l++)
1018  {
1019  for (unsigned i = 0; i < Dim; i++)
1020  {
1021  interpolated_x[i] += nodal_position(l, i) * psif[l];
1022  }
1023  }
1024 
1025  // Get the imposed flux
1026  std::complex<double> flux(0.0, 0.0);
1027  get_flux(interpolated_x, flux);
1028 
1029  // Now add to the appropriate equations
1030  // Loop over the test functions
1031  for (unsigned l = 0; l < n_node; l++)
1032  {
1033  local_eqn_real = nodal_local_eqn(l, U_index_helmholtz.real());
1034  /*IF it's not a boundary condition*/
1035  if (local_eqn_real >= 0)
1036  {
1037  // Add the prescribed flux terms
1038  residuals[local_eqn_real] -= flux.real() * testf[l] * W;
1039 
1040  // Imposed traction doesn't depend upon the solution,
1041  // --> the Jacobian is always zero, so no Jacobian
1042  // terms are required
1043  }
1044  local_eqn_imag = nodal_local_eqn(l, U_index_helmholtz.imag());
1045  /*IF it's not a boundary condition*/
1046  if (local_eqn_imag >= 0)
1047  {
1048  // Add the prescribed flux terms
1049  residuals[local_eqn_imag] -= flux.imag() * testf[l] * W;
1050 
1051  // Imposed traction doesn't depend upon the solution,
1052  // --> the Jacobian is always zero, so no Jacobian
1053  // terms are required
1054  }
1055  }
1056  }
1057  }
1058 
1059 
1060 } // namespace oomph
1061 
1062 #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...
int nodal_local_eqn(const unsigned &n, const unsigned &i) const
Return the local equation number corresponding to the i-th value at the n-th local node.
Definition: elements.h:1436
unsigned dim() const
Return the spatial dimension of the element, i.e. the number of local coordinates required to paramet...
Definition: elements.h:2615
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 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
bool has_hanging_nodes() const
Return boolean to indicate if any of the element's nodes are geometrically hanging.
Definition: elements.h:2474
unsigned long eqn_number(const unsigned &ieqn_local) const
Return the global equation number corresponding to the ieqn_local-th local equation number.
Definition: elements.h:708
int local_eqn_number(const unsigned long &ieqn_global) const
Return the local equation number corresponding to the ieqn_global-th global equation number....
Definition: elements.h:730
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....
A class for all isoparametric elements that solve the Helmholtz equations with pml capabilities....
virtual std::complex< unsigned > u_index_helmholtz() const
Broken assignment operator.
/////////////////////////////////////////////////////////////////////// /////////////////////////////...
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.
virtual std::complex< unsigned > u_index_helmholtz() const
Return the index at which the unknown value is stored.
void output(FILE *file_pt, const unsigned &n_plot)
C-style output function – forward to broken version in FiniteElement until somebody decides what exac...
std::complex< unsigned > U_index_helmholtz
The index at which the real and imag part of the unknown is stored at the nodes.
void output(std::ostream &outfile, const unsigned &n_plot)
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 get_dof_numbers_for_unknowns(std::list< std::pair< unsigned long, unsigned >> &dof_lookup_list) const
Create a list of pairs for all unknowns in this element, so that the first entry in each pair contain...
PMLHelmholtzFluxElement(const PMLHelmholtzFluxElement &dummy)=delete
Broken copy constructor.
double zeta_nodal(const unsigned &n, const unsigned &k, const unsigned &i) const
Specify the value of nodal zeta from the face geometry The "global" intrinsic coordinate of the eleme...
unsigned Dim
The spatial dimension of the problem.
unsigned ndof_types() const
The number of "DOF types" that degrees of freedom in this element are sub-divided into: real and imag...
PMLHelmholtzPrescribedFluxFctPt Flux_fct_pt
Function pointer to the (global) prescribed-flux function.
PMLHelmholtzPrescribedFluxFctPt & flux_fct_pt()
Broken assignment operator.
void get_flux(const Vector< double > &x, std::complex< double > &flux)
Function to calculate the prescribed flux at a given spatial position.
virtual void fill_in_generic_residual_contribution_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(FILE *file_pt)
C-style output function – forward to broken version in FiniteElement until somebody decides what exac...
PMLHelmholtzFluxElement()
Broken empty constructor.
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 ...
void fill_in_contribution_to_residuals(Vector< double > &residuals)
Add the element's contribution to its residual vector.
void(* PMLHelmholtzPrescribedFluxFctPt)(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 output(std::ostream &outfile)
Output function – forward to broken version in FiniteElement until somebody decides what exactly they...
A class for elements that allow the post-processing of radiated power and flux on the boundaries of P...
double global_power_contribution(std::ofstream &outfile)
Compute the element's contribution to the time-averaged radiated power over the artificial boundary....
std::complex< unsigned > U_index_helmholtz
The index at which the real and imag part of the unknown is stored at the nodes.
std::complex< double > global_flux_contribution()
Compute the element's contribution to the time-averaged flux.
double zeta_nodal(const unsigned &n, const unsigned &k, const unsigned &i) const
Broken assignment operator.
double global_power_contribution()
Compute the element's contribution to the time-averaged radiated power over the artificial boundary.
unsigned Dim
The spatial dimension of the problem.
PMLHelmholtzPowerElement()
Broken empty constructor.
PMLHelmholtzPowerElement(const PMLHelmholtzPowerElement &dummy)=delete
Broken copy constructor.
std::complex< double > global_flux_contribution(std::ofstream &outfile)
Compute the element's contribution to the integral of the flux over the artificial boundary....
virtual std::complex< unsigned > u_index_helmholtz() const
Return the index at which the unknown value is stored.
RefineableElements are FiniteElements that may be subdivided into children to provide a better local ...
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.
unsigned Dim
Dimension of zeta tuples (set by get_dim_helper) – needed because we store the scalar coordinates in ...
Definition: multi_domain.cc:60
//////////////////////////////////////////////////////////////////// ////////////////////////////////...