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-2023 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_FOURIER_DECOMPOSED_HELMHOLTZ_FLUX_ELEMENTS_HEADER
29 #define OOMPH_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  "FourierDecomposedHelmholtzFluxElement",
74  OOMPH_CURRENT_FUNCTION,
75  OOMPH_EXCEPTION_LOCATION);
76  }
77 
78  /// Broken copy constructor
80  const FourierDecomposedHelmholtzFluxElement& 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 FourierDecomposedHelmholtzFluxElement&) = delete;*/
89 
90 
91  /// Access function for the prescribed-flux function pointer
93  {
94  return Flux_fct_pt;
95  }
96 
97 
98  /// Add the element's contribution to its residual vector
100  {
101  // Call the generic residuals function with flag set to 0
102  // using a dummy matrix argument
104  residuals, GeneralisedElement::Dummy_matrix, 0);
105  }
106 
107 
108  /// Add the element's contribution to its residual vector and its
109  /// Jacobian matrix
111  DenseMatrix<double>& jacobian)
112  {
113  // Call the generic routine with the flag set to 1
115  residuals, jacobian, 1);
116  }
117 
118  /// Output function -- forward to broken version in FiniteElement
119  /// until somebody decides what exactly they want to plot here...
120  void output(std::ostream& outfile)
121  {
122  FiniteElement::output(outfile);
123  }
124 
125  /// Output function -- forward to broken version in FiniteElement
126  /// until somebody decides what exactly they want to plot here...
127  void output(std::ostream& outfile, const unsigned& n_plot)
128  {
129  FiniteElement::output(outfile, n_plot);
130  }
131 
132 
133  /// C-style output function -- forward to broken version in FiniteElement
134  /// until somebody decides what exactly they want to plot here...
135  void output(FILE* file_pt)
136  {
137  FiniteElement::output(file_pt);
138  }
139 
140  /// C-style output function -- forward to broken version in
141  /// FiniteElement until somebody decides what exactly they want to plot
142  /// here...
143  void output(FILE* file_pt, const unsigned& n_plot)
144  {
145  FiniteElement::output(file_pt, n_plot);
146  }
147 
148 
149  /// Return the index at which the unknown value
150  /// is stored. (Real/imag part gives real index of real/imag part).
151  virtual inline std::complex<unsigned> u_index_fourier_decomposed_helmholtz()
152  const
153  {
154  return std::complex<unsigned>(
157  }
158 
159 
160  protected:
161  /// Function to compute the shape and test functions and to return
162  /// the Jacobian of mapping between local and global (Eulerian)
163  /// coordinates
164  inline double shape_and_test(const Vector<double>& s,
165  Shape& psi,
166  Shape& test) const
167  {
168  // Find number of nodes
169  unsigned n_node = nnode();
170 
171  // Get the shape functions
172  shape(s, psi);
173 
174  // Set the test functions to be the same as the shape functions
175  for (unsigned i = 0; i < n_node; i++)
176  {
177  test[i] = psi[i];
178  }
179 
180  // Return the value of the jacobian
181  return J_eulerian(s);
182  }
183 
184 
185  /// Function to compute the shape and test functions and to return
186  /// the Jacobian of mapping between local and global (Eulerian)
187  /// coordinates
188  inline double shape_and_test_at_knot(const unsigned& ipt,
189  Shape& psi,
190  Shape& test) const
191  {
192  // Find number of nodes
193  unsigned n_node = nnode();
194 
195  // Get the shape functions
196  shape_at_knot(ipt, psi);
197 
198  // Set the test functions to be the same as the shape functions
199  for (unsigned i = 0; i < n_node; i++)
200  {
201  test[i] = psi[i];
202  }
203 
204  // Return the value of the jacobian
205  return J_eulerian_at_knot(ipt);
206  }
207 
208 
209  /// Function to calculate the prescribed flux at a given spatial
210  /// position
211  void get_flux(const Vector<double>& x, std::complex<double>& flux)
212  {
213  // If the function pointer is zero return zero
214  if (Flux_fct_pt == 0)
215  {
216  flux = std::complex<double>(0.0, 0.0);
217  }
218  // Otherwise call the function
219  else
220  {
221  (*Flux_fct_pt)(x, flux);
222  }
223  }
224 
225 
226  /// The index at which the real and imag part of the
227  /// unknown is stored at the nodes
228  std::complex<unsigned> U_index_fourier_decomposed_helmholtz;
229 
230 
231  /// Add the element's contribution to its residual vector.
232  /// flag=1(or 0): do (or don't) compute the contribution to the
233  /// Jacobian as well.
235  Vector<double>& residuals,
236  DenseMatrix<double>& jacobian,
237  const unsigned& flag);
238 
239 
240  /// Function pointer to the (global) prescribed-flux function
242  };
243 
244  /// ///////////////////////////////////////////////////////////////////
245  /// ///////////////////////////////////////////////////////////////////
246  /// ///////////////////////////////////////////////////////////////////
247 
248 
249  //===========================================================================
250  /// Constructor, takes the pointer to the "bulk" element, the
251  /// index of the fixed local coordinate and its value represented
252  /// by an integer (+/- 1), indicating that the face is located
253  /// at the max. or min. value of the "fixed" local coordinate
254  /// in the bulk element.
255  //===========================================================================
256  template<class ELEMENT>
259  const int& face_index)
260  : FaceGeometry<ELEMENT>(), FaceElement()
261  {
262  // Let the bulk element build the FaceElement, i.e. setup the pointers
263  // to its nodes (by referring to the appropriate nodes in the bulk
264  // element), etc.
265  bulk_el_pt->build_face_element(face_index, this);
266 
267  // Initialise the prescribed-flux function pointer to zero
268  Flux_fct_pt = 0;
269 
270 
271  // Initialise index at which real and imag unknowns are stored
272  U_index_fourier_decomposed_helmholtz = std::complex<unsigned>(0, 1);
273 
274  // Now read out indices from bulk element
276  dynamic_cast<FourierDecomposedHelmholtzEquations*>(bulk_el_pt);
277  // If the cast has failed die
278  if (eqn_pt == 0)
279  {
280  std::string error_string =
281  "Bulk element must inherit from FourierDecomposedHelmholtzEquations.";
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 are
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_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_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 } // namespace oomph
394 
395 #endif
static char t char * s
Definition: cfortran.h:568
cstr elem_len * i
Definition: cfortran.h:603
FaceElements are elements that coincide with the faces of higher-dimensional "bulk" elements....
Definition: elements.h:4338
int & face_index()
Index of the face (a number that uniquely identifies the face in the element)
Definition: elements.h:4626
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:5242
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:5328
//////////////////////////////////////////////////////////////////// ////////////////////////////////...
Definition: elements.h:4998
A general Finite Element class.
Definition: elements.h:1313
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:3050
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:2210
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:5132
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:3220
//////////////////////////////////////////////////////////////////// ////////////////////////////////...
virtual std::complex< unsigned > u_index_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)
C-style output function – forward to broken version in FiniteElement until somebody decides what exac...
virtual std::complex< unsigned > u_index_fourier_decomposed_helmholtz() const
Return the index at which the unknown value is stored. (Real/imag part gives real index of real/imag ...
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 get_flux(const Vector< double > &x, std::complex< double > &flux)
Function to calculate the prescribed flux at a given spatial position.
FourierDecomposedHelmholtzFluxElement(const FourierDecomposedHelmholtzFluxElement &dummy)=delete
Broken copy constructor.
FourierDecomposedHelmholtzPrescribedFluxFctPt & 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...
void fill_in_contribution_to_residuals(Vector< double > &residuals)
Add the element's contribution to its residual vector.
virtual void fill_in_generic_residual_contribution_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...
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 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)
Output function – forward to broken version in FiniteElement until somebody decides what exactly they...
std::complex< unsigned > U_index_fourier_decomposed_helmholtz
The index at which the real and imag part of the unknown is stored at the nodes.
void(* FourierDecomposedHelmholtzPrescribedFluxFctPt)(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_jacobian(Vector< double > &residuals, DenseMatrix< double > &jacobian)
Add the element's contribution to its residual vector and its Jacobian matrix.
FourierDecomposedHelmholtzPrescribedFluxFctPt Flux_fct_pt
Function pointer to the (global) prescribed-flux function.
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
An OomphLibError object which should be thrown when an run-time error is encountered....
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.
//////////////////////////////////////////////////////////////////// ////////////////////////////////...