linear_wave_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-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 elements that are used to apply prescribed flux
27// boundary conditions to the LinearWave equations
28#ifndef OOMPH_WAVE_FLUX_ELEMENTS_HEADER
29#define OOMPH_WAVE_FLUX_ELEMENTS_HEADER
30
31#ifdef HAVE_CONFIG_H
32#include <oomph-lib-config.h>
33#endif
34
35// oomph-lib ncludes
36#include "generic/Qelements.h"
37
38namespace oomph
39{
40 ////////////////////////////////////////////////////////////////////////
41 ////////////////////////////////////////////////////////////////////////
42 ////////////////////////////////////////////////////////////////////////
43
44 //======================================================================
45 /// A class for elements that allow the imposition of an
46 /// applied flux on the boundaries of LinearWave elements.
47 /// The element geometry is obtained from the FaceGeometry<ELEMENT>
48 /// policy class.
49 //======================================================================
50 template<class ELEMENT>
51 class LinearWaveFluxElement : public virtual FaceGeometry<ELEMENT>,
52 public virtual FaceElement
53 {
54 public:
55 /// Function pointer to the prescribed-flux function fct(x,f(x)) --
56 /// x is a Vector!
57 typedef void (*LinearWavePrescribedFluxFctPt)(const double& time,
58 const Vector<double>& x,
59 double& flux);
60
61
62 /// Constructor, takes the pointer to the "bulk" element and the
63 /// index of the face to be created
65
66
67 /// Broken empty constructor
69 {
70 throw OomphLibError(
71 "Don't call empty constructor for LinearWaveFluxElement",
74 }
75
76 /// Broken copy constructor
78
79 /// Broken assignment operator
80 void operator=(const LinearWaveFluxElement&) = delete;
81
82
83 /// Access function for the prescribed-flux function pointer
88
89
90 /// Compute the element residual vector
92 {
93 // Call the generic residuals function with flag set to 0
94 // using the dummy matrix argument
97 }
98
99
100 /// Compute the element's residual vector and its Jacobian matrix
102 DenseMatrix<double>& jacobian)
103 {
104 // Call the generic routine with the flag set to 1
106 residuals, jacobian, 1);
107 }
108
109 /// Specify the value of nodal zeta from the face geometry
110 /// The "global" intrinsic coordinate of the element when
111 /// viewed as part of a geometric object should be given by
112 /// the FaceElement representation, by default (needed to break
113 /// indeterminacy if bulk element is SolidElement)
114 double zeta_nodal(const unsigned& n,
115 const unsigned& k,
116 const unsigned& i) const
117 {
118 return FaceElement::zeta_nodal(n, k, i);
119 }
120
121 /// Output function -- forward to broken version in FiniteElement
122 /// until somebody decides what exactly they want to plot here...
123 void output(std::ostream& outfile)
124 {
126 }
127
128 /// Output function -- forward to broken version in FiniteElement
129 /// until somebody decides what exactly they want to plot here...
130 void output(std::ostream& outfile, const unsigned& n_plot)
131 {
133 }
134
135
136 /// Output function -- forward to broken version in FiniteElement
137 /// until somebody decides what exactly they want to plot here...
142
143 /// Output function -- forward to broken version in FiniteElement
144 /// until somebody decides what exactly they want to plot here...
145 void output(FILE* file_pt, const unsigned& n_plot)
146 {
148 }
149
150
151 protected:
152 /// Function to compute the shape and test functions and to return
153 /// the Jacobian of mapping between local and global (Eulerian)
154 /// coordinates
155 inline double shape_and_test(const Vector<double>& s,
156 Shape& psi,
157 Shape& test) const
158 {
159 // Find number of nodes
160 unsigned n_node = nnode();
161
162 // Get the shape functions
163 shape(s, psi);
164
165 // Set the test functions to be the same as the shape functions
166 for (unsigned i = 0; i < n_node; i++)
167 {
168 test[i] = psi[i];
169 }
170
171 // Return the value of the jacobian
172 return J_eulerian(s);
173 }
174
175
176 /// Function to calculate the prescribed flux at a given spatial
177 /// position and at a gien time
178 void get_flux(const double& time, const Vector<double>& x, double& flux)
179 {
180 // If the function pointer is zero return zero
181 if (Flux_fct_pt == 0)
182 {
183 flux = 0.0;
184 }
185 // Otherwise call the function
186 else
187 {
188 (*Flux_fct_pt)(time, x, flux);
189 }
190 }
191
192
193 private:
194 /// Compute the element residual vector.
195 /// flag=1(or 0): do (or don't) compute the Jacobian as well.
197 Vector<double>& residuals, DenseMatrix<double>& jacobian, unsigned flag);
198
199 /// Function pointer to the (global) prescribed-flux function
201
202 /// The spatial dimension of the problem
203 unsigned Dim;
204
205 /// The index at which the unknown is stored at the nodes
207 };
208
209
210 ////////////////////////////////////////////////////////////////////////
211 ////////////////////////////////////////////////////////////////////////
212 ////////////////////////////////////////////////////////////////////////
213
214
215 //===========================================================================
216 /// Constructor, takes the pointer to the "bulk" element and the
217 /// index of the face to be created.
218 //===========================================================================
219 template<class ELEMENT>
221 FiniteElement* bulk_el_pt, const int& face_index)
222 : FaceGeometry<ELEMENT>(), FaceElement()
223 {
224 // Let the bulk element build the FaceElement, i.e. setup the pointers
225 // to its nodes (by referring to the appropriate nodes in the bulk
226 // element), etc.
228
229#ifdef PARANOID
230 {
231 // Check that the element is not a refineable 3d element
232 ELEMENT* elem_pt = dynamic_cast<ELEMENT*>(bulk_el_pt);
233
234 // If it's three-d
235 if (elem_pt->dim() == 3)
236 {
237 // Is it refineable
239 dynamic_cast<RefineableElement*>(elem_pt);
240 if (ref_el_pt != 0)
241 {
242 if (this->has_hanging_nodes())
243 {
244 throw OomphLibError("This flux element will not work correctly if "
245 "nodes are hanging\n",
248 }
249 }
250 }
251 }
252#endif
253
254 // Initialise the prescribed-flux function pointer to zero
255 Flux_fct_pt = 0;
256
257 // Extract the dimension of the problem from the dimension of
258 // the first node
259 Dim = node_pt(0)->ndim();
260
261 // Set up U_index_lin_wave. Initialise to zero, which probably won't change
262 // in most cases, oh well, the price we pay for generality
264
265 // Cast to the appropriate LinearWaveEquation so that we can
266 // find the index at which the variable is stored
267 // We assume that the dimension of the full problem is the same
268 // as the dimension of the node, if this is not the case you will have
269 // to write custom elements, sorry
270 switch (Dim)
271 {
272 // One dimensional problem
273 case 1:
274 {
276 dynamic_cast<LinearWaveEquations<1>*>(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 LinearWaveEquations.";
282 error_string +=
283 "Nodes are one dimensional, but cannot cast the bulk element to\n";
284 error_string += "LinearWaveEquations<1>\n.";
285 error_string += "If you desire this functionality, you must "
286 "implement it yourself\n";
287
288 throw OomphLibError(
290 }
291 // Otherwise read out the value
292 else
293 {
294 // Read the index from the (cast) bulk element
295 U_index_lin_wave = eqn_pt->u_index_lin_wave();
296 }
297 }
298 break;
299
300 // Two dimensional problem
301 case 2:
302 {
304 dynamic_cast<LinearWaveEquations<2>*>(bulk_el_pt);
305 // If the cast has failed die
306 if (eqn_pt == 0)
307 {
308 std::string error_string =
309 "Bulk element must inherit from LinearWaveEquations.";
310 error_string +=
311 "Nodes are two dimensional, but cannot cast the bulk element to\n";
312 error_string += "LinearWaveEquations<2>\n.";
313 error_string += "If you desire this functionality, you must "
314 "implement it yourself\n";
315
316 throw OomphLibError(
318 }
319 else
320 {
321 // Read the index from the (cast) bulk element.
322 U_index_lin_wave = eqn_pt->u_index_lin_wave();
323 }
324 }
325 break;
326
327 // Three dimensional problem
328 case 3:
329 {
331 dynamic_cast<LinearWaveEquations<3>*>(bulk_el_pt);
332 // If the cast has failed die
333 if (eqn_pt == 0)
334 {
335 std::string error_string =
336 "Bulk element must inherit from LinearWaveEquations.";
337 error_string += "Nodes are three dimensional, but cannot cast the "
338 "bulk element to\n";
339 error_string += "LinearWaveEquations<3>\n.";
340 error_string += "If you desire this functionality, you must "
341 "implement it yourself\n";
342
343 throw OomphLibError(
345 }
346 else
347 {
348 // Read the index from the (cast) bulk element.
349 U_index_lin_wave = eqn_pt->u_index_lin_wave();
350 }
351 }
352 break;
353
354 // Any other case is an error
355 default:
356 std::ostringstream error_stream;
357 error_stream << "Dimension of node is " << Dim
358 << ". It should be 1,2, or 3!" << std::endl;
359
360 throw OomphLibError(
362 break;
363 }
364 }
365
366
367 //===========================================================================
368 /// Compute the element's residual vector and the (zero) Jacobian matrix.
369 //===========================================================================
370 template<class ELEMENT>
373 Vector<double>& residuals, DenseMatrix<double>& jacobian, unsigned flag)
374 {
375 // Find out how many nodes there are
376 const unsigned n_node = nnode();
377
378 // Get continuous time from timestepper of first node
379 double time = node_pt(0)->time_stepper_pt()->time_pt()->time();
380
381 // Set up memory for the shape and test functions
383
384 // Set the value of n_intpt
385 const unsigned n_intpt = integral_pt()->nweight();
386
387 // Set the Vector to hold local coordinates
388 Vector<double> s(Dim - 1);
389
390 // Integer to store the local equation and unknowns
391 int local_eqn = 0;
392
393 // Locally cache the index at which the variable is stored
394 const unsigned u_index_lin_wave = U_index_lin_wave;
395
396 // Loop over the integration points
397 //--------------------------------
398 for (unsigned ipt = 0; ipt < n_intpt; ipt++)
399 {
400 // Assign values of s
401 for (unsigned i = 0; i < (Dim - 1); i++)
402 {
403 s[i] = integral_pt()->knot(ipt, i);
404 }
405
406 // Get the integral weight
407 double w = integral_pt()->weight(ipt);
408
409 // Find the shape and test functions and return the Jacobian
410 // of the mapping
411 double J = shape_and_test(s, psif, testf);
412
413 // Premultiply the weights and the Jacobian
414 double W = w * J;
415
416 // Need to find position to feed into flux function
418
419 // Initialise to zero
420 for (unsigned i = 0; i < Dim; i++)
421 {
422 interpolated_x[i] = 0.0;
423 }
424
425 // Calculate velocities and derivatives
426 for (unsigned l = 0; l < n_node; l++)
427 {
428 // Loop over velocity components
429 for (unsigned i = 0; i < Dim; i++)
430 {
432 }
433 }
434
435 // Get the imposed flux
436 double flux;
437 get_flux(time, interpolated_x, flux);
438
439 // Now add to the appropriate equations
440
441 // Loop over the test functions
442 for (unsigned l = 0; l < n_node; l++)
443 {
444 local_eqn = nodal_local_eqn(l, u_index_lin_wave);
445 /*IF it's not a boundary condition*/
446 if (local_eqn >= 0)
447 {
448 // Add the prescribed flux terms
449 residuals[local_eqn] -= flux * testf[l] * W;
450
451 // Imposed traction doesn't depend upon the solution,
452 // --> the Jacobian is always zero, so no Jacobian
453 // terms are required
454 }
455 }
456 }
457 }
458
459
460} // namespace oomph
461
462#endif
static char t char * s
Definition cfortran.h:568
cstr elem_len * i
Definition cfortran.h:603
void get_flux(const Vector< double > &s, Vector< double > &flux) const
Get flux: .
TimeStepper *& time_stepper_pt()
Return the pointer to the timestepper.
Definition nodes.h:238
FaceElements are elements that coincide with the faces of higher-dimensional "bulk" elements....
Definition elements.h:4342
int & face_index()
Index of the face (a number that uniquely identifies the face in the element)
Definition elements.h:4630
double zeta_nodal(const unsigned &n, const unsigned &k, const unsigned &i) const
In a FaceElement, the "global" intrinsic coordinate of the element along the boundary,...
Definition elements.h:4501
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:5273
FaceGeometry class definition: This policy class is used to allow construction of face elements that ...
Definition elements.h:5002
A general Finite Element class.
Definition elements.h:1317
Integral *const & integral_pt() const
Return the pointer to the integration scheme (const version)
Definition elements.h:1967
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:3054
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...
virtual double interpolated_x(const Vector< double > &s, const unsigned &i) const
Return FE interpolated coordinate x[i] at local coordinate s.
Definition elements.cc:3992
int nodal_local_eqn(const unsigned &n, const unsigned &i) const
Return the local equation number corresponding to the i-th value at the n-th local node.
Definition elements.h:1436
unsigned dim() const
Return the spatial dimension of the element, i.e. the number of local coordinates required to paramet...
Definition elements.h:2615
unsigned nnode() const
Return the number of nodes.
Definition elements.h:2214
double nodal_position(const unsigned &n, const unsigned &i) const
Return the i-th coordinate at local node n. If the node is hanging, the appropriate interpolation is ...
Definition elements.h:2321
Node *& node_pt(const unsigned &n)
Return a pointer to the local node n.
Definition elements.h:2179
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:5163
bool has_hanging_nodes() const
Return boolean to indicate if any of the element's nodes are geometrically hanging.
Definition elements.h:2474
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
virtual double knot(const unsigned &i, const unsigned &j) const =0
Return local coordinate s[j] of i-th integration point.
virtual unsigned nweight() const =0
Return the number of integration points of the scheme.
virtual double weight(const unsigned &i) const =0
Return weight of i-th integration point.
A class for elements that allow the imposition of an applied flux on the boundaries of LinearWave ele...
void output(FILE *file_pt, 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)
Compute the element residual vector.
void(* LinearWavePrescribedFluxFctPt)(const double &time, const Vector< double > &x, double &flux)
Function pointer to the prescribed-flux function fct(x,f(x)) – x is a Vector!
LinearWavePrescribedFluxFctPt & flux_fct_pt()
Access function for the prescribed-flux function pointer.
unsigned U_index_lin_wave
The index at which the unknown is stored at the nodes.
LinearWaveFluxElement(const LinearWaveFluxElement &dummy)=delete
Broken copy constructor.
double zeta_nodal(const unsigned &n, const unsigned &k, const unsigned &i) const
Specify the value of nodal zeta from the face geometry The "global" intrinsic coordinate of the eleme...
void fill_in_contribution_to_jacobian(Vector< double > &residuals, DenseMatrix< double > &jacobian)
Compute the element's residual vector and its Jacobian matrix.
unsigned Dim
The spatial dimension of the problem.
void output(FILE *file_pt)
Output function – forward to broken version in FiniteElement until somebody decides what exactly they...
void operator=(const LinearWaveFluxElement &)=delete
Broken assignment operator.
LinearWavePrescribedFluxFctPt Flux_fct_pt
Function pointer to the (global) prescribed-flux function.
void get_flux(const double &time, const Vector< double > &x, double &flux)
Function to calculate the prescribed flux at a given spatial position and at a gien time.
LinearWaveFluxElement()
Broken empty constructor.
void fill_in_generic_residual_contribution_lin_wave_flux(Vector< double > &residuals, DenseMatrix< double > &jacobian, unsigned flag)
Compute the element residual vector. flag=1(or 0): do (or don't) compute the Jacobian as well.
void output(std::ostream &outfile)
Output function – forward to broken version in FiniteElement until somebody decides what exactly they...
void output(std::ostream &outfile, const unsigned &n_plot)
Output function – forward to broken version in FiniteElement until somebody decides what exactly they...
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 ...
unsigned ndim() const
Return (Eulerian) spatial dimension of the node.
Definition nodes.h:1054
An OomphLibError object which should be thrown when an run-time error is encountered....
RefineableElements are FiniteElements that may be subdivided into children to provide a better local ...
A Class for shape functions. In simple cases, the shape functions have only one index that can be tho...
Definition shape.h:76
TAdvectionDiffusionReactionElement<NREAGENT,DIM,NNODE_1D> elements are isoparametric triangular DIM-d...
TAdvectionDiffusionReactionElement()
Constructor: Call constructors for TElement and AdvectionDiffusionReaction equations.
Time *const & time_pt() const
Access function for the pointer to time (const version)
double & time()
Return the current value of the continuous time.
DRAIG: Change all instances of (SPATIAL_DIM) to (DIM-1).