refineable_pml_helmholtz_elements.h
Go to the documentation of this file.
1 // LIC// ====================================================================
2 // LIC// This file forms part of oomph-lib, the object-oriented,
3 // LIC// multi-physics finite-element library, available
4 // LIC// at http://www.oomph-lib.org.
5 // LIC//
6 // LIC// Copyright (C) 2006-2022 Matthias Heil and Andrew Hazel
7 // LIC//
8 // LIC// This library is free software; you can redistribute it and/or
9 // LIC// modify it under the terms of the GNU Lesser General Public
10 // LIC// License as published by the Free Software Foundation; either
11 // LIC// version 2.1 of the License, or (at your option) any later version.
12 // LIC//
13 // LIC// This library is distributed in the hope that it will be useful,
14 // LIC// but WITHOUT ANY WARRANTY; without even the implied warranty of
15 // LIC// MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
16 // LIC// Lesser General Public License for more details.
17 // LIC//
18 // LIC// You should have received a copy of the GNU Lesser General Public
19 // LIC// License along with this library; if not, write to the Free Software
20 // LIC// Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA
21 // LIC// 02110-1301 USA.
22 // LIC//
23 // LIC// The authors may be contacted at oomph-lib@maths.man.ac.uk.
24 // LIC//
25 // LIC//====================================================================
26 // Header file for refineable QPMLHelmholtzElement elements
27 
28 #ifndef OOMPH_REFINEABLE_PML_HELMHOLTZ_ELEMENTS_HEADER
29 #define OOMPH_REFINEABLE_PML_HELMHOLTZ_ELEMENTS_HEADER
30 
31 // Config header generated by autoconfig
32 #ifdef HAVE_CONFIG_H
33 #include <oomph-lib-config.h>
34 #endif
35 
36 
37 // oomph-lib headers
38 #include "../generic/refineable_quad_element.h"
39 #include "../generic/refineable_brick_element.h"
40 #include "../generic/hp_refineable_elements.h"
41 #include "../generic/error_estimator.h"
42 #include "pml_helmholtz_elements.h"
43 
44 namespace oomph
45 {
48 
49 
50  //======================================================================
54  //======================================================================
55  template<unsigned DIM>
57  : public virtual PMLHelmholtzEquations<DIM>,
58  public virtual RefineableElement,
59  public virtual ElementWithZ2ErrorEstimator
60  {
61  public:
67  {
68  }
69 
72  const RefineablePMLHelmholtzEquations<DIM>& dummy) = delete;
73 
75  // Commented out broken assignment operator because this can lead to a
76  // conflict warning when used in the virtual inheritence hierarchy.
77  // Essentially the compiler doesn't realise that two separate
78  // implementations of the broken function are the same and so, quite
79  // rightly, it shouts.
80  /*void operator=(const RefineablePMLHelmholtzEquations<DIM>&) = delete;*/
81 
83  unsigned num_Z2_flux_terms()
84  {
85  return 2 * DIM;
86  }
87 
91  {
92  Vector<std::complex<double>> complex_flux(DIM);
93  this->get_flux(s, complex_flux);
94  unsigned count = 0;
95  for (unsigned i = 0; i < DIM; i++)
96  {
97  flux[count] = complex_flux[i].real();
98  count++;
99  flux[count] = complex_flux[i].imag();
100  count++;
101  }
102  }
103 
104 
110  Vector<double>& values)
111  {
112  // Set size of Vector: u
113  values.resize(2);
114 
115  // Find number of nodes
116  unsigned n_node = nnode();
117 
118  // Local shape function
119  Shape psi(n_node);
120 
121  // Find values of shape function
122  shape(s, psi);
123 
124  // Initialise value of u
125  values[0] = 0.0;
126  values[1] = 0.0;
127 
128  // Find the index at which the pml_helmholtz unknown is stored
129  std::complex<unsigned> u_nodal_index = this->u_index_helmholtz();
130 
131  // Loop over the local nodes and sum up the values
132  for (unsigned l = 0; l < n_node; l++)
133  {
134  values[0] += this->nodal_value(l, u_nodal_index.real()) * psi[l];
135  values[1] += this->nodal_value(l, u_nodal_index.imag()) * psi[l];
136  }
137  }
138 
139 
144  void get_interpolated_values(const unsigned& t,
145  const Vector<double>& s,
146  Vector<double>& values)
147  {
148  if (t != 0)
149  {
150  std::string error_message =
151  "Time-dependent version of get_interpolated_values() ";
152  error_message += "not implemented for this element \n";
153  throw OomphLibError(
154  error_message,
155  "RefineablePMLHelmholtzEquations::get_interpolated_values()",
157  }
158  else
159  {
160  // Make sure that we call this particular object's steady
161  // get_interpolated_values (it could get overloaded lower down)
163  values);
164  }
165  }
166 
167 
170  {
172  this->father_element_pt())
173  ->source_fct_pt();
174  }
175 
176 
177  private:
183  Vector<double>& residuals,
184  DenseMatrix<double>& jacobian,
185  const unsigned& flag);
186  };
187 
188 
189  //======================================================================
193  //======================================================================
194  template<unsigned DIM, unsigned NNODE_1D>
196  : public QPMLHelmholtzElement<DIM, NNODE_1D>,
197  public virtual RefineablePMLHelmholtzEquations<DIM>,
198  public virtual RefineableQElement<DIM>
199  {
200  public:
203  : RefineableElement(),
206  QPMLHelmholtzElement<DIM, NNODE_1D>()
207  {
208  }
209 
210 
213  const RefineableQPMLHelmholtzElement<DIM, NNODE_1D>& dummy) = delete;
214 
216  /*void operator=(const RefineableQPMLHelmholtzElement<DIM,NNODE_1D>&) =
217  * delete;*/
218 
220  unsigned ncont_interpolated_values() const
221  {
222  return 2;
223  }
224 
226  unsigned nvertex_node() const
227  {
229  }
230 
232  Node* vertex_node_pt(const unsigned& j) const
233  {
235  }
236 
238  void rebuild_from_sons(Mesh*& mesh_pt) {}
239 
242  unsigned nrecovery_order()
243  {
244  return (NNODE_1D - 1);
245  }
246 
250  };
251 
252 
256 
257 
258  //=======================================================================
264  //=======================================================================
265  template<unsigned DIM, unsigned NNODE_1D>
267  : public virtual QElement<DIM - 1, NNODE_1D>
268  {
269  public:
272  FaceGeometry() : QElement<DIM - 1, NNODE_1D>() {}
273  };
274 
275  //=======================================================================
278  //=======================================================================
279  template<unsigned NNODE_1D>
281  : public virtual RefineableQPMLHelmholtzElement<2, NNODE_1D>
282  {
283  public:
287  };
288 
289 } // namespace oomph
290 
291 #endif
int i
Definition: BiCGSTAB_step_by_step.cpp:9
Definition: error_estimator.h:79
FaceGeometry()
Definition: refineable_pml_helmholtz_elements.h:272
Definition: elements.h:4998
double nodal_value(const unsigned &n, const unsigned &i) const
Definition: elements.h:2593
virtual unsigned nvertex_node() const
Definition: elements.h:2491
virtual void shape(const Vector< double > &s, Shape &psi) const =0
unsigned nnode() const
Return the number of nodes.
Definition: elements.h:2210
virtual Node * vertex_node_pt(const unsigned &j) const
Definition: elements.h:2500
Definition: mesh.h:67
Definition: nodes.h:906
Definition: oomph_definitions.h:222
Definition: pml_helmholtz_elements.h:62
PMLHelmholtzSourceFctPt Source_fct_pt
Pointer to source function:
Definition: pml_helmholtz_elements.h:727
void get_flux(const Vector< double > &s, Vector< std::complex< double >> &flux) const
Get flux: flux[i] = du/dx_i for real and imag part.
Definition: pml_helmholtz_elements.h:446
PMLHelmholtzSourceFctPt & source_fct_pt()
Access function: Pointer to source function.
Definition: pml_helmholtz_elements.h:400
virtual std::complex< unsigned > u_index_helmholtz() const
Broken assignment operator.
Definition: pml_helmholtz_elements.h:91
PMLLayerElement()
Definition: refineable_pml_helmholtz_elements.h:286
Definition: pml_meshes.h:48
Definition: Qelements.h:459
Definition: pml_helmholtz_elements.h:754
Definition: refineable_elements.h:97
virtual RefineableElement * father_element_pt() const
Return a pointer to the father element.
Definition: refineable_elements.h:539
Definition: refineable_pml_helmholtz_elements.h:60
void get_interpolated_values(const unsigned &t, const Vector< double > &s, Vector< double > &values)
Definition: refineable_pml_helmholtz_elements.h:144
void further_build()
Further build: Copy source function pointer from father element.
Definition: refineable_pml_helmholtz_elements.h:169
void get_interpolated_values(const Vector< double > &s, Vector< double > &values)
Definition: refineable_pml_helmholtz_elements.h:109
unsigned num_Z2_flux_terms()
Broken assignment operator.
Definition: refineable_pml_helmholtz_elements.h:83
void get_Z2_flux(const Vector< double > &s, Vector< double > &flux)
Definition: refineable_pml_helmholtz_elements.h:90
void fill_in_generic_residual_contribution_helmholtz(Vector< double > &residuals, DenseMatrix< double > &jacobian, const unsigned &flag)
Definition: refineable_pml_helmholtz_elements.cc:41
RefineablePMLHelmholtzEquations()
Constructor, simply call other constructors.
Definition: refineable_pml_helmholtz_elements.h:63
RefineablePMLHelmholtzEquations(const RefineablePMLHelmholtzEquations< DIM > &dummy)=delete
Broken copy constructor.
Definition: Qelements.h:2259
Definition: refineable_pml_helmholtz_elements.h:199
unsigned nrecovery_order()
Definition: refineable_pml_helmholtz_elements.h:242
unsigned ncont_interpolated_values() const
Broken assignment operator.
Definition: refineable_pml_helmholtz_elements.h:220
Node * vertex_node_pt(const unsigned &j) const
Pointer to the j-th vertex node in the element.
Definition: refineable_pml_helmholtz_elements.h:232
unsigned nvertex_node() const
Number of vertex nodes in the element.
Definition: refineable_pml_helmholtz_elements.h:226
RefineableQPMLHelmholtzElement()
Constructor, simply call the other constructors.
Definition: refineable_pml_helmholtz_elements.h:202
void rebuild_from_sons(Mesh *&mesh_pt)
Rebuild from sons: empty.
Definition: refineable_pml_helmholtz_elements.h:238
RefineableQPMLHelmholtzElement(const RefineableQPMLHelmholtzElement< DIM, NNODE_1D > &dummy)=delete
Broken copy constructor.
void further_setup_hanging_nodes()
Definition: refineable_pml_helmholtz_elements.h:249
Definition: shape.h:76
RealScalar s
Definition: level1_cplx_impl.h:130
#define DIM
Definition: linearised_navier_stokes_elements.h:44
void flux(const double &time, const Vector< double > &x, double &flux)
Get flux applied along boundary x=0.
Definition: pretend_melt.cc:59
std::string string(const unsigned &i)
Definition: oomph_definitions.cc:286
DRAIG: Change all instances of (SPATIAL_DIM) to (DIM-1).
Definition: AnisotropicHookean.h:10
t
Definition: plotPSD.py:36
#define OOMPH_EXCEPTION_LOCATION
Definition: oomph_definitions.h:61
std::ptrdiff_t j
Definition: tut_arithmetic_redux_minmax.cpp:2