refineable_poisson_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 
32 #ifndef OOMPH_REFINEABLE_POISSON_ELEMENTS_HEADER
33 #define OOMPH_REFINEABLE_POISSON_ELEMENTS_HEADER
34 
35 // Config header generated by autoconfig
36 #ifdef HAVE_CONFIG_H
37  #include <oomph-lib-config.h>
38 #endif
39 
40 
41 //oomph-lib headers
42 #include "../generic/refineable_quad_element.h"
43 #include "../generic/refineable_brick_element.h"
44 #include "../generic/hp_refineable_elements.h"
45 #include "../generic/error_estimator.h"
46 #include "poisson_elements.h"
47 
48 namespace oomph
49 {
50 
51 ///////////////////////////////////////////////////////////////////////////
52 ///////////////////////////////////////////////////////////////////////////
53 
54 
55 
56 //======================================================================
57 /// Refineable version of Poisson equations
58 ///
59 ///
60 //======================================================================
61 template <unsigned DIM>
62 class RefineablePoissonEquations : public virtual PoissonEquations<DIM>,
63  public virtual RefineableElement,
64  public virtual ElementWithZ2ErrorEstimator
65 {
66  public:
67 
68  /// \short Constructor, simply call other constructors
71  { }
72 
73  /// Broken copy constructor
75  {
76  BrokenCopy::broken_copy("RefineablePoissonEquations");
77  }
78 
79  /// Broken assignment operator
81  {
82  BrokenCopy::broken_assign("RefineablePoissonEquations");
83  }
84 
85  /// Number of 'flux' terms for Z2 error estimation
86  unsigned num_Z2_flux_terms() {return DIM;}
87 
88  /// Get 'flux' for Z2 error recovery: Standard flux.from Poisson equations
90  {this->get_flux(s,flux);}
91 
92  /// Get error against and norm of exact flux
94  std::ostream &outfile,
96  double& error, double& norm);
97 
98 /// \short Get the function value u in Vector.
99 /// Note: Given the generality of the interface (this function
100 /// is usually called from black-box documentation or interpolation routines),
101 /// the values Vector sets its own size in here.
103  {
104  // Set size of Vector: u
105  values.resize(1);
106 
107  //Find number of nodes
108  unsigned n_node = nnode();
109 
110  //Local shape function
111  Shape psi(n_node);
112 
113  //Find values of shape function
114  shape(s,psi);
115 
116  //Initialise value of u
117  values[0] = 0.0;
118 
119  //Find the index at which the poisson unknown is stored
120  unsigned u_nodal_index = this->u_index_poisson();
121 
122  //Loop over the local nodes and sum up the values
123  for(unsigned l=0;l<n_node;l++)
124  {
125  values[0] += this->nodal_value(l,u_nodal_index)*psi[l];
126  }
127  }
128 
129 
130  /// \short Get the function value u in Vector.
131  /// Note: Given the generality of the interface (this function
132  /// is usually called from black-box documentation or interpolation routines),
133  /// the values Vector sets its own size in here.
134  void get_interpolated_values(const unsigned& t, const Vector<double>&s,
135  Vector<double>& values)
136  {
137  if (t!=0)
138  {
139  std::string error_message =
140  "Time-dependent version of get_interpolated_values() ";
141  error_message += "not implemented for this element \n";
142  throw
143  OomphLibError(error_message,
144  "RefineablePoissonEquations::get_interpolated_values()",
145  OOMPH_EXCEPTION_LOCATION);
146  }
147  else
148  {
149  //Make sure that we call this particular object's steady
150  //get_interpolated_values (it could get overloaded lower down)
152  }
153  }
154 
155 
156  /// Further build: Copy source function pointer from father element
158  {
159  this->Source_fct_pt=dynamic_cast<RefineablePoissonEquations<DIM>*>(
160  this->father_element_pt())->source_fct_pt();
161  }
162 
163 
164  private:
165 
166 
167 /// \short Add element's contribution to elemental residual vector and/or
168 /// Jacobian matrix
169 /// flag=1: compute both
170 /// flag=0: compute only residual vector
172  Vector<double> &residuals, DenseMatrix<double> &jacobian,
173  const unsigned& flag);
174 
175  /// \short Compute derivatives of elemental residual vector with respect
176  /// to nodal coordinates. Overwrites default implementation in
177  /// FiniteElement base class.
178  /// dresidual_dnodal_coordinates(l,i,j) = d res(l) / dX_{ij}
180  dresidual_dnodal_coordinates);
181 
182 };
183 
184 
185 //======================================================================
186 /// Refineable version of 2D QPoissonElement elements
187 ///
188 ///
189 //======================================================================
190 template <unsigned DIM, unsigned NNODE_1D>
192  public QPoissonElement<DIM,NNODE_1D>,
193  public virtual RefineablePoissonEquations<DIM>,
194  public virtual RefineableQElement<DIM>
195 {
196  public:
197 
198  /// \short Constructor, simply call the other constructors
202  RefineableQElement<DIM>(),
203  QPoissonElement<DIM,NNODE_1D>()
204  {}
205 
206 
207  /// Broken copy constructor
209  dummy)
210  {
211  BrokenCopy::broken_copy("RefineableQuadPoissonElement");
212  }
213 
214  /// Broken assignment operator
216  {
217  BrokenCopy::broken_assign("RefineableQuadPoissonElement");
218  }
219 
220  /// Number of continuously interpolated values: 1
221  unsigned ncont_interpolated_values() const {return 1;}
222 
223  /// \short Number of vertex nodes in the element
224  unsigned nvertex_node() const
226 
227  /// \short Pointer to the j-th vertex node in the element
228  Node* vertex_node_pt(const unsigned& j) const
230 
231  /// Rebuild from sons: empty
232  void rebuild_from_sons(Mesh* &mesh_pt) {}
233 
234  /// \short Order of recovery shape functions for Z2 error estimation:
235  /// Same order as shape functions.
236  unsigned nrecovery_order() {return (NNODE_1D-1);}
237 
238  /// \short Perform additional hanging node procedures for variables
239  /// that are not interpolated by all nodes. Empty.
241 
242 };
243 
244 
245 //======================================================================
246 /// p-refineable version of 2D QPoissonElement elements
247 //======================================================================
248 template<unsigned DIM>
250  public virtual RefineablePoissonEquations<DIM>,
251  public virtual PRefineableQElement<DIM>
252 {
253  public:
254 
255  /// \short Constructor, simply call the other constructors
258  PRefineableQElement<DIM>(),
259  QPoissonElement<DIM,2>()
260  {
261  // Set integration scheme
262  // (To avoid memory leaks in pre-build and p-refine where new
263  // integration schemes are created)
265  }
266 
267  /// Destructor (to avoid memory leaks)
269  {
270  delete this->integral_pt();
271  }
272 
273 
274  /// Broken copy constructor
276  {
277  BrokenCopy::broken_copy("PRefineableQPoissonElement");
278  }
279 
280  /// Broken assignment operator
282  {
283  BrokenCopy::broken_assign("PRefineableQPoissonElement");
284  }
285 
286  void further_build();
287 
288  /// Number of continuously interpolated values: 1
289  unsigned ncont_interpolated_values() const {return 1;}
290 
291  /// \short Number of vertex nodes in the element
292  unsigned nvertex_node() const
294 
295  /// \short Pointer to the j-th vertex node in the element
296  Node* vertex_node_pt(const unsigned& j) const
298 
299  /// \short Order of recovery shape functions for Z2 error estimation:
300  /// - Same order as shape functions.
301  //unsigned nrecovery_order()
302  // {
303  // if(this->nnode_1d() < 4) {return (this->nnode_1d()-1);}
304  // else {return 3;}
305  // }
306  /// - Constant recovery order, since recovery order of the first element
307  /// is used for the whole mesh.
308  unsigned nrecovery_order() {return 3;}
309 
310  void compute_energy_error(std::ostream &outfile,
312  double& error, double& norm);
313 };
314 
315 
316 
317 ////////////////////////////////////////////////////////////////////////
318 ////////////////////////////////////////////////////////////////////////
319 ////////////////////////////////////////////////////////////////////////
320 
321 
322 
323 //=======================================================================
324 /// Face geometry for the RefineableQuadPoissonElement elements: The spatial
325 /// dimension of the face elements is one lower than that of the
326 /// bulk element but they have the same number of points
327 /// along their 1D edges.
328 //=======================================================================
329 template<unsigned DIM, unsigned NNODE_1D>
331  public virtual QElement<DIM-1,NNODE_1D>
332 {
333 
334  public:
335 
336  /// \short Constructor: Call the constructor for the
337  /// appropriate lower-dimensional QElement
338  FaceGeometry() : QElement<DIM-1,NNODE_1D>() {}
339 
340 };
341 
342 }
343 
344 #endif
345 
unsigned nrecovery_order()
Order of recovery shape functions for Z2 error estimation: Same order as shape functions.
void broken_copy(const std::string &class_name)
Issue error message and terminate execution.
unsigned nrecovery_order()
Order of recovery shape functions for Z2 error estimation:
Base class for finite elements that can compute the quantities that are required for the Z2 error est...
void rebuild_from_sons(Mesh *&mesh_pt)
Rebuild from sons: empty.
virtual RefineableElement * father_element_pt() const
Return a pointer to the father element.
virtual void get_dresidual_dnodal_coordinates(RankThreeTensor< double > &dresidual_dnodal_coordinates)
Compute derivatives of elemental residual vector with respect to nodal coordinates. Overwrites default implementation in FiniteElement base class. dresidual_dnodal_coordinates(l,i,j) = d res(l) / dX_{ij}.
Node * vertex_node_pt(const unsigned &j) const
Pointer to the j-th vertex node in the element.
virtual void set_integration_scheme(Integral *const &integral_pt)
Set the spatial integration scheme.
Definition: elements.cc:3150
unsigned nvertex_node() const
Number of vertex nodes in the element.
char t
Definition: cfortran.h:572
p-refineable version of 2D QPoissonElement elements
Nodes are derived from Data, but, in addition, have a definite (Eulerian) position in a space of a gi...
Definition: nodes.h:852
RefineablePoissonEquations()
Constructor, simply call other constructors.
RefineableQPoissonElement()
Constructor, simply call the other constructors.
Node * vertex_node_pt(const unsigned &j) const
Pointer to the j-th vertex node in the element.
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...
PRefineableQPoissonElement(const PRefineableQPoissonElement< DIM > &dummy)
Broken copy constructor.
void operator=(const RefineableQPoissonElement< DIM, NNODE_1D > &)
Broken assignment operator.
void(* SteadyExactSolutionFctPt)(const Vector< double > &, Vector< double > &)
Function pointer for function that computes vector-valued steady "exact solution" as ...
Definition: elements.h:1723
virtual unsigned u_index_poisson() const
Return the index at which the unknown value is stored. The default value, 0, is appropriate for singl...
A Rank 3 Tensor class.
Definition: matrices.h:1337
Integral *const & integral_pt() const
Return the pointer to the integration scheme (const version)
Definition: elements.h:1908
unsigned nvertex_node() const
Number of vertex nodes in the element.
virtual unsigned nvertex_node() const
Definition: elements.h:2374
static char t char * s
Definition: cfortran.h:572
RefineableQPoissonElement(const RefineableQPoissonElement< DIM, NNODE_1D > &dummy)
Broken copy constructor.
unsigned ncont_interpolated_values() const
Number of continuously interpolated values: 1.
~PRefineableQPoissonElement()
Destructor (to avoid memory leaks)
RefineablePoissonEquations(const RefineablePoissonEquations< DIM > &dummy)
Broken copy constructor.
unsigned num_Z2_flux_terms()
Number of &#39;flux&#39; terms for Z2 error estimation.
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 get_flux(const Vector< double > &s, Vector< double > &flux) const
Get flux: flux[i] = du/dx_i.
void operator=(const PRefineableQPoissonElement< DIM > &)
Broken assignment operator.
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 broken_assign(const std::string &class_name)
Issue error message and terminate execution.
PoissonSourceFctPt & source_fct_pt()
Access function: Pointer to source function.
void operator=(const RefineablePoissonEquations< DIM > &)
Broken assignment operator.
void further_setup_hanging_nodes()
Perform additional hanging node procedures for variables that are not interpolated by all nodes...
PRefineableQPoissonElement()
Constructor, simply call the other constructors.
std::string string(const unsigned &i)
Return the i-th string or "" if the relevant string hasn&#39;t been defined.
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
FaceGeometry()
Constructor: Call the constructor for the appropriate lower-dimensional QElement. ...
void fill_in_generic_residual_contribution_poisson(Vector< double > &residuals, DenseMatrix< double > &jacobian, const unsigned &flag)
Add element&#39;s contribution to elemental residual vector and/or Jacobian matrix flag=1: compute both f...
unsigned ncont_interpolated_values() const
Number of continuously interpolated values: 1.
void get_Z2_flux(const Vector< double > &s, Vector< double > &flux)
Get &#39;flux&#39; for Z2 error recovery: Standard flux.from Poisson equations.
unsigned nnode() const
Return the number of nodes.
Definition: elements.h:2146
void compute_exact_Z2_error(std::ostream &outfile, FiniteElement::SteadyExactSolutionFctPt exact_flux_pt, double &error, double &norm)
Get error against and norm of exact flux.
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...
void further_build()
Further build: Copy source function pointer from father element.
A general mesh class.
Definition: mesh.h:74
PoissonSourceFctPt Source_fct_pt
Pointer to source function: