refineable_linear_wave_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// Version 1.0; svn revision $LastChangedRevision$
7 //LIC//
8 //LIC// $LastChangedDate$
9 //LIC//
10 //LIC// Copyright (C) 2006-2016 Matthias Heil and Andrew Hazel
11 //LIC//
12 //LIC// This library is free software; you can redistribute it and/or
13 //LIC// modify it under the terms of the GNU Lesser General Public
14 //LIC// License as published by the Free Software Foundation; either
15 //LIC// version 2.1 of the License, or (at your option) any later version.
16 //LIC//
17 //LIC// This library is distributed in the hope that it will be useful,
18 //LIC// but WITHOUT ANY WARRANTY; without even the implied warranty of
19 //LIC// MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
20 //LIC// Lesser General Public License for more details.
21 //LIC//
22 //LIC// You should have received a copy of the GNU Lesser General Public
23 //LIC// License along with this library; if not, write to the Free Software
24 //LIC// Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA
25 //LIC// 02110-1301 USA.
26 //LIC//
27 //LIC// The authors may be contacted at oomph-lib@maths.man.ac.uk.
28 //LIC//
29 //LIC//====================================================================
30 //Header file for refineable QPoissonElement elements
31 #ifndef OOMPH_REFINEABLE_LINEAR_WAVE_ELEMENTS_HEADER
32 #define OOMPH_REFINEABLE_LINEAR_WAVE_ELEMENTS_HEADER
33 
34 // Config header generated by autoconfig
35 #ifdef HAVE_CONFIG_H
36  #include <oomph-lib-config.h>
37 #endif
38 
39 //oomph-lib headers
40 #include "../generic/refineable_quad_element.h"
41 #include "../generic/refineable_brick_element.h"
42 #include "../generic/error_estimator.h"
43 #include "linear_wave_elements.h"
44 
45 namespace oomph
46 {
47 
48 //======================================================================
49 /// \short Refineable version of LinearWave equations.
50 //======================================================================
51 template <unsigned DIM>
53  public virtual RefineableElement,
54  public virtual ElementWithZ2ErrorEstimator
55 {
56  public:
57 
58  /// \short Constructor
60  LinearWaveEquations<DIM>(),
63  {}
64 
65 
66  /// Broken copy constructor
68  {
69  BrokenCopy::broken_copy("RefineableLinearWaveEquations");
70  }
71 
72  /// Broken assignment operator
74  {
75  BrokenCopy::broken_assign("RefineableLinearWaveEquations");
76  }
77 
78  /// Number of 'flux' terms for Z2 error estimation
79  unsigned num_Z2_flux_terms() {return DIM;}
80 
81  /// Get 'flux' for Z2 error recovery: Standard flux.from LinearWave equations
83  {this->get_flux(s,flux);}
84 
85 
86 /// \short Get the function value u in Vector.
87 /// Note: Given the generality of the interface (this function
88 /// is usually called from black-box documentation or interpolation routines),
89 /// the values Vector sets its own size in here.
91  {
92  // Set size of Vector: u
93  values.resize(1);
94 
95  //Find number of nodes
96  unsigned n_node = nnode();
97 
98  //Find the index at which the unknown is stored
99  unsigned u_nodal_index = this->u_index_lin_wave();
100 
101  //Local shape function
102  Shape psi(n_node);
103 
104  //Find values of shape function
105  shape(s,psi);
106 
107  //Initialise value of u
108  values[0] = 0.0;
109 
110  //Loop over the local nodes and sum
111  for(unsigned l=0;l<n_node;l++)
112  {
113  values[0] += this->nodal_value(l,u_nodal_index)*psi[l];
114  }
115  }
116 
117 
118  /// \short Get the function value u in Vector.
119  /// Note: Given the generality of the interface (this function
120  /// is usually called from black-box documentation or interpolation routines),
121  /// the values Vector sets its own size in here.
122  void get_interpolated_values(const unsigned& t, const Vector<double>&s,
123  Vector<double>& values)
124  {
125  // Set size of Vector:
126  values.resize(1);
127 
128  // Initialise
129  values[0]=0.0;
130 
131  //Find out how many nodes there are
132  unsigned n_node = nnode();
133 
134  //Find the index at which the unknown is stored
135  unsigned u_nodal_index = this->u_index_lin_wave();
136 
137  // Shape functions
138  Shape psi(n_node);
139  shape(s,psi);
140 
141  //Calculate value
142  for(unsigned l=0;l<n_node;l++)
143  {
144  values[0] += this->nodal_value(t,l,u_nodal_index)*psi[l];
145  }
146  }
147 
148 
149  /// Further build: Copy source function pointer from father element
151  {
153  this->father_element_pt())->source_fct_pt();
154  }
155 
156 
157  private:
158 
159 /// \short Add element's contribution to elemental residual vector and/or
160 /// Jacobian matrix
161 /// flag=1: compute both
162 /// flag=0: compute only residual vector
164  Vector<double> &residuals, DenseMatrix<double> &jacobian,
165  unsigned flag);
166 
167 };
168 
169 
170 //======================================================================
171 /// Refineable version of 2D QLinearWaveElement elements
172 ///
173 ///
174 //======================================================================
175 template <unsigned DIM, unsigned NNODE_1D>
177 public QLinearWaveElement<DIM,NNODE_1D>,
178  public virtual RefineableLinearWaveEquations<DIM>,
179  public virtual RefineableQElement<DIM>
180 {
181  public:
182 
183  /// \short Constructor
187  RefineableQElement<DIM>(),
188  QLinearWaveElement<DIM,NNODE_1D>()
189  {}
190 
191 
192  /// Broken copy constructor
194  dummy)
195  {
196  BrokenCopy::broken_copy("RefineableQuadLinearWaveElement");
197  }
198 
199  /// Broken assignment operator
201  {
202  BrokenCopy::broken_assign("RefineableQuadLinearWaveElement");
203  }
204 
205  /// Number of continuously interpolated values: 1
206  unsigned ncont_interpolated_values() const {return 1;}
207 
208  /// \short Number of vertex nodes in the element
209  unsigned nvertex_node() const
211 
212  /// \short Pointer to the j-th vertex node in the element
213  Node* vertex_node_pt(const unsigned& j) const
215 
216  /// Rebuild from sons: empty
217  void rebuild_from_sons(Mesh* &mesh_pt) {}
218 
219  /// \short Order of recovery shape functions for Z2 error estimation:
220  /// Same order as shape functions.
221  unsigned nrecovery_order() {return (NNODE_1D-1);}
222 
223  /// \short Perform additional hanging node procedures for variables
224  /// that are not interpolated by all nodes. Empty.
226 
227 };
228 ////////////////////////////////////////////////////////////////////////
229 ////////////////////////////////////////////////////////////////////////
230 ////////////////////////////////////////////////////////////////////////
231 
232 
233 
234 //=======================================================================
235 /// Face geometry for the RefineableQuadLinearWaveElement elements: The spatial
236 /// dimension of the face elements is one lower than that of the
237 /// bulk element but they have the same number of points
238 /// along their 1D edges.
239 //=======================================================================
240 template<unsigned DIM, unsigned NNODE_1D>
242  public virtual QElement<DIM-1,NNODE_1D>
243 {
244 
245  public:
246 
247  /// \short Constructor: Call the constructor for the
248  /// appropriate lower-dimensional QElement
249  FaceGeometry() : QElement<DIM-1,NNODE_1D>() {}
250 
251 };
252 
253 }
254 
255 #endif
256 
LinearWaveSourceFctPt & source_fct_pt()
Access function: Pointer to source function.
void broken_copy(const std::string &class_name)
Issue error message and terminate execution.
Base class for finite elements that can compute the quantities that are required for the Z2 error est...
void fill_in_generic_residual_contribution_lin_wave(Vector< double > &residuals, DenseMatrix< double > &jacobian, unsigned flag)
Add element&#39;s contribution to elemental residual vector and/or Jacobian matrix flag=1: compute both f...
virtual RefineableElement * father_element_pt() const
Return a pointer to the father element.
void rebuild_from_sons(Mesh *&mesh_pt)
Rebuild from sons: empty.
void further_setup_hanging_nodes()
Perform additional hanging node procedures for variables that are not interpolated by all nodes...
void get_interpolated_values(const unsigned &t, const Vector< double > &s, Vector< double > &values)
Get the function value u in Vector. Note: Given the generality of the interface (this function is usu...
void operator=(const RefineableQLinearWaveElement< DIM, NNODE_1D > &)
Broken assignment operator.
char t
Definition: cfortran.h:572
Nodes are derived from Data, but, in addition, have a definite (Eulerian) position in a space of a gi...
Definition: nodes.h:852
unsigned nrecovery_order()
Order of recovery shape functions for Z2 error estimation: Same order as shape functions.
void operator=(const RefineableLinearWaveEquations< DIM > &)
Broken assignment operator.
RefineableQLinearWaveElement(const RefineableQLinearWaveElement< DIM, NNODE_1D > &dummy)
Broken copy constructor.
virtual unsigned u_index_lin_wave() const
Return the index at which the unknown value is stored. The default value, 0, is appropriate for singl...
LinearWaveSourceFctPt Source_fct_pt
Pointer to source function:
void further_build()
Further build: Copy source function pointer from father element.
Refineable version of LinearWave equations.
void get_flux(const Vector< double > &s, Vector< double > &flux) const
Get flux: flux[i] = du/dx_i.
RefineableLinearWaveEquations(const RefineableLinearWaveEquations< DIM > &dummy)
Broken copy constructor.
void get_interpolated_values(const Vector< double > &s, Vector< double > &values)
Get the function value u in Vector. Note: Given the generality of the interface (this function is usu...
virtual unsigned nvertex_node() const
Definition: elements.h:2374
static char t char * s
Definition: cfortran.h:572
FaceGeometry()
Constructor: Call the constructor for the appropriate lower-dimensional QElement. ...
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:2470
void broken_assign(const std::string &class_name)
Issue error message and terminate execution.
Node * vertex_node_pt(const unsigned &j) const
Pointer to the j-th vertex node in the element.
unsigned ncont_interpolated_values() const
Number of continuously interpolated values: 1.
virtual Node * vertex_node_pt(const unsigned &j) const
Pointer to the j-th vertex node in the element. Broken virtual function in "pure" finite elements...
Definition: elements.h:2383
unsigned nvertex_node() const
Number of vertex nodes in the element.
unsigned nnode() const
Return the number of nodes.
Definition: elements.h:2146
void get_Z2_flux(const Vector< double > &s, Vector< double > &flux)
Get &#39;flux&#39; for Z2 error recovery: Standard flux.from LinearWave equations.
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...
A general mesh class.
Definition: mesh.h:74
unsigned num_Z2_flux_terms()
Number of &#39;flux&#39; terms for Z2 error estimation.