refineable_spectral_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// Copyright (C) 2006-2025 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 QSpectralPoissonElement elements
27
28#ifndef OOMPH_REFINEABLE_SPECTRAL_POISSON_ELEMENTS_HEADER
29#define OOMPH_REFINEABLE_SPECTRAL_POISSON_ELEMENTS_HEADER
30
31// Config header
32#ifdef HAVE_CONFIG_H
33#include <oomph-lib-config.h>
34#endif
35
36
37// oomph-lib headers
42
43namespace oomph
44{
45 //======================================================================
46 /// Refineable version of 2D QSpectralPoissonElement elements
47 ///
48 ///
49 //======================================================================
50 template<unsigned DIM, unsigned NNODE_1D>
52 : public QSpectralPoissonElement<DIM, NNODE_1D>,
53 public virtual RefineablePoissonEquations<DIM>,
54 public virtual RefineableQSpectralElement<DIM>
55 {
56 public:
57 /// Constructor: Pass refinement level to refineable quad element
58 /// (default 0 = root)
66
67
68 /// 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 RefineableQSpectralPoissonElement<DIM,NNODE_1D>&) =
79 * delete;*/
80
81 /// Number of continuously interpolated values: 1
83 {
84 return 1;
85 }
86
87 /// Number of vertex nodes in the element
92
93 /// Pointer to the j-th vertex node in the element
94 Node* vertex_node_pt(const unsigned& j) const
95 {
97 }
98
104
105 /// Function to describe the local dofs of the element. The ostream
106 /// specifies the output stream to which the description
107 /// is written; the string stores the currently
108 /// assembled output that is ultimately written to the
109 /// output stream by Data::describe_dofs(...); it is typically
110 /// built up incrementally as we descend through the
111 /// call hierarchy of this function when called from
112 /// Problem::describe_dofs(...)
113 void describe_local_dofs(std::ostream& out,
114 const std::string& current_string) const
115 {
117 }
118
119 /// Order of recovery shape functions for Z2 error estimation:
120 /// Same order as shape functions.
122 {
123 if (NNODE_1D < 4)
124 {
125 return (NNODE_1D - 1);
126 }
127 else
128 {
129 return 3;
130 }
131 }
132
133 /// Perform additional hanging node procedures for variables
134 /// that are not interpolated by all nodes. Empty.
136 };
137
138
139 //=======================================================================
140 /// Face geometry for the RefineableQuadPoissonElement elements: The spatial
141 /// dimension of the face elements is one lower than that of the
142 /// bulk element but they have the same number of points
143 /// along their 1D edges.
144 //=======================================================================
145 template<unsigned DIM, unsigned NNODE_1D>
147 : public virtual QSpectralElement<DIM - 1, NNODE_1D>
148 {
149 public:
150 /// Constructor: Call the constructor for the
151 /// appropriate lower-dimensional QElement
153 };
154
155} // namespace oomph
156
157
158#endif
FaceGeometry()
Constructor: Call the constructor for the appropriate lower-dimensional QElement.
FaceGeometry class definition: This policy class is used to allow construction of face elements that ...
Definition elements.h:5002
virtual unsigned nvertex_node() const
Return the number of vertex nodes in this element. Broken virtual function in "pure" finite elements.
Definition elements.h:2495
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:2504
virtual void assign_all_generic_local_eqn_numbers(const bool &store_local_dof_pt)
Assign all the local equation numbering schemes that can be applied generically for the element....
Definition elements.h:253
virtual void describe_local_dofs(std::ostream &out, const std::string &current_string) const
Function to describe the local dofs of the element. The ostream specifies the output stream to which ...
Definition elements.cc:578
Nodes are derived from Data, but, in addition, have a definite (Eulerian) position in a space of a gi...
Definition nodes.h:906
General QLegendreElement class.
QSpectralPoissonElement elements are linear/quadrilateral/brick-shaped Poisson elements with isoparam...
RefineableElements are FiniteElements that may be subdivided into children to provide a better local ...
Refineable version of Poisson equations.
A class that is used to template the refineable Q elements by dimension. It's really nothing more tha...
Definition Qelements.h:2259
A class that is used to template the refineable Q spectral elements by dimension. It's really nothing...
Refineable version of 2D QSpectralPoissonElement elements.
unsigned nvertex_node() const
Number of vertex nodes in the element.
unsigned ncont_interpolated_values() const
Broken assignment operator.
unsigned nrecovery_order()
Order of recovery shape functions for Z2 error estimation: Same order as shape functions.
void further_setup_hanging_nodes()
Perform additional hanging node procedures for variables that are not interpolated by all nodes....
RefineableQSpectralPoissonElement()
Constructor: Pass refinement level to refineable quad element (default 0 = root)
void describe_local_dofs(std::ostream &out, const std::string &current_string) const
Function to describe the local dofs of the element. The ostream specifies the output stream to which ...
void assign_all_generic_local_eqn_numbers(const bool &store_local_dof_pt)
Overloaded version of the calculation of the local equation numbers. If the boolean argument is true ...
Node * vertex_node_pt(const unsigned &j) const
Pointer to the j-th vertex node in the element.
RefineableQSpectralPoissonElement(const RefineableQSpectralPoissonElement< DIM, NNODE_1D > &dummy)=delete
Broken copy constructor.
TAdvectionDiffusionReactionElement<NREAGENT,DIM,NNODE_1D> elements are isoparametric triangular DIM-d...
DRAIG: Change all instances of (SPATIAL_DIM) to (DIM-1).