linearised_navier_stokes_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-2023 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 linearised axisymmetric Navier-Stokes elements
27 
28 #ifndef OOMPH_LINEARISED_NAVIER_STOKES_ELEMENTS_HEADER
29 #define OOMPH_LINEARISED_NAVIER_STOKES_ELEMENTS_HEADER
30 
31 // Config header generated by autoconfig
32 #ifdef HAVE_CONFIG_H
33 #include <oomph-lib-config.h>
34 #endif
35 
36 // oomph-lib includes
37 #include "../generic/Qelements.h"
38 #include "../generic/fsi.h"
39 
41 
42 namespace oomph
43 {
44 #define DIM 2
45 
46 
47  //=======================================================================
48  /// A class for elements that solve the linearised version of the
49  /// unsteady Navier--Stokes equations in cylindrical polar coordinates,
50  /// where we have Fourier-decomposed in the azimuthal direction so that
51  /// the theta-dependance is replaced by an azimuthal mode number.
52  //=======================================================================
54  {
55  private:
56  /// Static "magic" number that indicates that the pressure is not
57  /// stored at a node
59 
60  /// Static default value for the physical constants
61  /// (all initialised to zero)
63 
64  /// Static default value for the physical ratios (all initialised to one)
66 
67  protected:
68  // Physical constants
69  // ------------------
70 
71  /// Pointer to the viscosity ratio (relative to the
72  /// viscosity used in the definition of the Reynolds number)
74 
75  /// Pointer to the density ratio (relative to the
76  /// density used in the definition of the Reynolds number)
78 
79  /// Pointer to global Reynolds number
80  double* Re_pt;
81 
82  /// Pointer to global Reynolds number x Strouhal number (=Womersley)
83  double* ReSt_pt;
84 
85  /// Pointer to eigenvalue
86  double* Lambda_pt;
87 
88  /// Pointer to frequency
89  double* Omega_pt;
90 
91  /// Pointer to the normalisation element
94 
95 
96  /// Index of datum where eigenvalue is stored
98 
100 
101  /// Pointer to base flow solution (velocity components) function
102  void (*Base_flow_u_fct_pt)(const double& time,
103  const Vector<double>& x,
104  Vector<double>& result);
105 
106  /// Pointer to derivatives of base flow solution velocity
107  /// components w.r.t. global coordinates (r and z) function
108  void (*Base_flow_dudx_fct_pt)(const double& time,
109  const Vector<double>& x,
110  DenseMatrix<double>& result);
111 
112  /// Boolean flag to indicate if ALE formulation is disabled when
113  /// the time-derivatives are computed. Only set to true if you're sure
114  /// that the mesh is stationary.
116 
117  /// Access function for the local equation number
118  /// information for the i-th component of the pressure.
119  /// p_local_eqn[n,i] = local equation number or < 0 if pinned.
120  virtual int p_local_eqn(const unsigned& n, const unsigned& i) = 0;
121 
122  /// Compute the shape functions and their derivatives
123  /// w.r.t. global coordinates at local coordinate s.
124  /// Return Jacobian of mapping between local and global coordinates.
126  const Vector<double>& s,
127  Shape& psi,
128  DShape& dpsidx,
129  Shape& test,
130  DShape& dtestdx) const = 0;
131 
132  /// Compute the shape functions and their derivatives
133  /// w.r.t. global coordinates at the ipt-th integration point.
134  /// Return Jacobian of mapping between local and global coordinates.
136  const unsigned& ipt,
137  Shape& psi,
138  DShape& dpsidx,
139  Shape& test,
140  DShape& dtestdx) const = 0;
141 
142  /// Compute the pressure shape functions at local coordinate s
144  Shape& psi) const = 0;
145 
146  /// Compute the pressure shape and test functions at local coordinate s
148  Shape& psi,
149  Shape& test) const = 0;
150 
151  /// Calculate the velocity components of the base flow solution
152  /// at a given time and Eulerian position
153  virtual void get_base_flow_u(const double& time,
154  const unsigned& ipt,
155  const Vector<double>& x,
156  Vector<double>& result) const
157  {
158  // If the function pointer is zero return zero
159  if (Base_flow_u_fct_pt == 0)
160  {
161  // Loop over velocity components and set base flow solution to zero
162  for (unsigned i = 0; i < DIM; i++)
163  {
164  result[i] = 0.0;
165  }
166  }
167  // Otherwise call the function
168  else
169  {
170  (*Base_flow_u_fct_pt)(time, x, result);
171  }
172  }
173 
174  /// Calculate the derivatives of the velocity components of the
175  /// base flow solution w.r.t. global coordinates (r and z) at a given
176  /// time and Eulerian position
177  virtual void get_base_flow_dudx(const double& time,
178  const unsigned& ipt,
179  const Vector<double>& x,
180  DenseMatrix<double>& result) const
181  {
182  // If the function pointer is zero return zero
183  if (Base_flow_dudx_fct_pt == 0)
184  {
185  // Loop over velocity components
186  for (unsigned i = 0; i < DIM; i++)
187  {
188  // Loop over coordinate directions and set to zero
189  for (unsigned j = 0; j < DIM; j++)
190  {
191  result(i, j) = 0.0;
192  }
193  }
194  }
195  // Otherwise call the function
196  else
197  {
198  (*Base_flow_dudx_fct_pt)(time, x, result);
199  }
200  }
201 
202 
203  inline int eigenvalue_local_eqn(const unsigned& i)
204  {
205  return this->external_local_eqn(this->Data_number_of_eigenvalue,
206  this->Index_of_eigenvalue + i);
207  }
208 
209 
210  /// Compute the residuals for the Navier-Stokes equations;
211  /// flag=1(or 0): do (or don't) compute the Jacobian as well.
213  Vector<double>& residuals,
214  DenseMatrix<double>& jacobian,
215  DenseMatrix<double>& mass_matrix,
216  unsigned flag);
217 
218  public:
219  /// Constructor: NULL the base flow solution and the
220  /// derivatives of the base flow function
223  {
224  // Set all the physical parameter pointers to the default value of zero
227 
230 
231  // Set to sensible defaults
234 
235  // Set the physical ratios to the default value of one
238 
239  // Null out normalisation
241  }
242 
243  /// Vector to decide whether the stress-divergence form is used or not.
244  // N.B. This needs to be public so that the intel compiler gets things
245  // correct. Somehow the access function messes things up when going to
246  // refineable navier--stokes
248 
249  // Access functions for the physical constants
250  // -------------------------------------------
251 
252  /// Reynolds number
253  const double& re() const
254  {
255  return *Re_pt;
256  }
257 
258  /// Product of Reynolds and Strouhal number (=Womersley number)
259  const double& re_st() const
260  {
261  return *ReSt_pt;
262  }
263 
264  const double& lambda() const
265  {
266  return *Lambda_pt;
267  }
268 
269  const double& omega() const
270  {
271  return *Omega_pt;
272  }
273 
274  /// Pointer to Reynolds number
275  double*& re_pt()
276  {
277  return Re_pt;
278  }
279 
280  /// Pointer to product of Reynolds and Strouhal number (=Womersley number)
281  double*& re_st_pt()
282  {
283  return ReSt_pt;
284  }
285 
286  /// Pointer to lambda
287  double*& lambda_pt()
288  {
289  return Lambda_pt;
290  }
291 
292  /// Pointer to frequency
293  double*& omega_pt()
294  {
295  return Omega_pt;
296  }
297 
298  /// Pointer to normalisation element
300  {
302  }
303 
304  /// the boolean flag check_nodal_data is set to false.
307  normalisation_el_pt)
308  {
309  // Set the normalisation element
310  Normalisation_element_pt = normalisation_el_pt;
311 
312  // Add eigenvalue unknown as external data to this element
314  this->add_external_data(normalisation_el_pt->eigenvalue_data_pt());
315 
316  // Which value corresponds to the eigenvalue
317  Index_of_eigenvalue = normalisation_el_pt->index_of_eigenvalue();
318 
319  // Now set the pointers to the eigenvalues
320  Lambda_pt = normalisation_el_pt->eigenvalue_data_pt()->value_pt(
322  Omega_pt = normalisation_el_pt->eigenvalue_data_pt()->value_pt(
323  Index_of_eigenvalue + 1);
324  }
325 
326 
327  /// Viscosity ratio for element: element's viscosity relative
328  /// to the viscosity used in the definition of the Reynolds number
329  const double& viscosity_ratio() const
330  {
331  return *Viscosity_Ratio_pt;
332  }
333 
334  /// Pointer to the viscosity ratio
336  {
337  return Viscosity_Ratio_pt;
338  }
339 
340  /// Density ratio for element: element's density relative
341  /// to the viscosity used in the definition of the Reynolds number
342  const double& density_ratio() const
343  {
344  return *Density_Ratio_pt;
345  }
346 
347  /// Pointer to the density ratio
348  double*& density_ratio_pt()
349  {
350  return Density_Ratio_pt;
351  }
352 
353  /// Access function for the base flow solution pointer
354  void (*&base_flow_u_fct_pt())(const double& time,
355  const Vector<double>& x,
356  Vector<double>& f)
357  {
358  return Base_flow_u_fct_pt;
359  }
360 
361  /// Access function for the derivatives of the base flow
362  /// w.r.t. global coordinates solution pointer
363  void (*&base_flow_dudx_fct_pt())(const double& time,
364  const Vector<double>& x,
366  {
367  return Base_flow_dudx_fct_pt;
368  }
369 
370  /// Return the number of pressure degrees of freedom
371  /// associated with a single pressure component in the element
372  virtual unsigned npres_linearised_nst() const = 0;
373 
374  /// Return the index at which the i-th unknown velocity
375  /// component is stored. The default value, i, is appropriate for
376  /// single-physics problems. In derived multi-physics elements, this
377  /// function should be overloaded to reflect the chosen storage scheme.
378  /// Note that these equations require that the unknowns are always
379  /// stored at the same indices at each node.
380  virtual inline unsigned u_index_linearised_nst(const unsigned& i) const
381  {
382  return i;
383  }
384 
385  /// Return the i-th component of du/dt at local node n.
386  /// Uses suitably interpolated value for hanging nodes.
387  double du_dt_linearised_nst(const unsigned& n, const unsigned& i) const
388  {
389  // Get the data's timestepper
391 
392  // Initialise dudt
393  double dudt = 0.0;
394 
395  // Loop over the timesteps, if there is a non-steady timestepper
396  if (!time_stepper_pt->is_steady())
397  {
398  // Get the index at which the velocity is stored
399  const unsigned u_nodal_index = u_index_linearised_nst(i);
400 
401  // Determine number of timsteps (past & present)
402  const unsigned n_time = time_stepper_pt->ntstorage();
403 
404  // Add the contributions to the time derivative
405  for (unsigned t = 0; t < n_time; t++)
406  {
407  dudt +=
408  time_stepper_pt->weight(1, t) * nodal_value(t, n, u_nodal_index);
409  }
410  }
411 
412  return dudt;
413  }
414 
415  /// Disable ALE, i.e. assert the mesh is not moving -- you do this
416  /// at your own risk!
417  void disable_ALE()
418  {
419  ALE_is_disabled = true;
420  }
421 
422  /// (Re-)enable ALE, i.e. take possible mesh motion into account
423  /// when evaluating the time-derivative. Note: By default, ALE is
424  /// enabled, at the expense of possibly creating unnecessary work
425  /// in problems where the mesh is, in fact, stationary.
426  void enable_ALE()
427  {
428  ALE_is_disabled = false;
429  }
430 
431  /// Return the i-th pressure value at local pressure "node" n_p.
432  /// Uses suitably interpolated value for hanging nodes.
433  virtual double p_linearised_nst(const unsigned& n_p,
434  const unsigned& i) const = 0;
435 
436  /// Pin the real or imaginary part of the problem
437  /// Input integer 0 for real 1 for imaginary
438  virtual void pin_real_or_imag(const unsigned& real) = 0;
439  virtual void unpin_real_or_imag(const unsigned& real) = 0;
440 
441  /// Pin the normalisation dofs
443 
444  /// Which nodal value represents the pressure?
445  // N.B. This function has return type "int" (rather than "unsigned"
446  // as in the u_index case) so that we can return the "magic" number
447  // "Pressure_not_stored_at_node" ( = -100 )
448  virtual inline int p_index_linearised_nst(const unsigned& i) const
449  {
451  }
452 
453  /// Strain-rate tensor: \f$ e_{ij} \f$
454  /// where \f$ i,j = r,z,\theta \f$ (in that order)
455  void strain_rate(const Vector<double>& s,
457  const unsigned& real) const;
458 
459  /// Output function: r, z, U^C, U^S, V^C, V^S, W^C, W^S, P^C, P^S
460  /// in tecplot format. Default number of plot points
461  void output(std::ostream& outfile)
462  {
463  const unsigned nplot = 5;
464  output(outfile, nplot);
465  }
466 
467  /// Output function: r, z, U^C, U^S, V^C, V^S, W^C, W^S, P^C, P^S
468  /// in tecplot format. nplot points in each coordinate direction
469  void output(std::ostream& outfile, const unsigned& nplot);
470 
471  /// Output function: r, z, U^C, U^S, V^C, V^S, W^C, W^S, P^C, P^S
472  /// in tecplot format. Default number of plot points
473  void output(FILE* file_pt)
474  {
475  const unsigned nplot = 5;
476  output(file_pt, nplot);
477  }
478 
479  /// Output function: r, z, U^C, U^S, V^C, V^S, W^C, W^S, P^C, P^S
480  /// in tecplot format. nplot points in each coordinate direction
481  void output(FILE* file_pt, const unsigned& nplot);
482 
483  /// Output function: r, z, U^C, U^S, V^C, V^S, W^C, W^S,
484  /// in tecplot format. nplot points in each coordinate direction
485  /// at timestep t (t=0: present; t>0: previous timestep)
486  void output_veloc(std::ostream& outfile,
487  const unsigned& nplot,
488  const unsigned& t);
489 
490  /// Compute the element's residual Vector
492  {
493  // Call the generic residuals function with flag set to 0
494  // and using a dummy matrix argument
496  residuals,
499  0);
500  }
501 
502  /// Compute the element's residual Vector and the jacobian matrix.
503  /// Virtual function can be overloaded by hanging-node version.
504  /*void fill_in_contribution_to_jacobian(Vector<double> &residuals,
505  DenseMatrix<double> &jacobian)
506  {
507  // Call the generic routine with the flag set to 1
508  fill_in_generic_residual_contribution_linearised_nst(
509  residuals,jacobian,GeneralisedElement::Dummy_matrix,1);
510  }*/
511 
512  /// Add the element's contribution to its residuals vector,
513  /// jacobian matrix and mass matrix
514  /*void fill_in_contribution_to_jacobian_and_mass_matrix(
515  Vector<double> &residuals, DenseMatrix<double> &jacobian,
516  DenseMatrix<double> &mass_matrix)
517  {
518  // Call the generic routine with the flag set to 2
519  fill_in_generic_residual_contribution_linearised_nst(
520  residuals,jacobian,mass_matrix,2);
521  }*/
522 
523  /// Return the i-th component of the FE interpolated velocity
524  /// u[i] at local coordinate s
526  const unsigned& i) const
527  {
528  // Determine number of nodes in the element
529  const unsigned n_node = nnode();
530 
531  // Provide storage for local shape functions
532  Shape psi(n_node);
533 
534  // Find values of shape functions
535  shape(s, psi);
536 
537  // Get the index at which the velocity is stored
538  const unsigned u_nodal_index = u_index_linearised_nst(i);
539 
540  // Initialise value of u
541  double interpolated_u = 0.0;
542 
543  // Loop over the local nodes and sum
544  for (unsigned l = 0; l < n_node; l++)
545  {
546  interpolated_u += nodal_value(l, u_nodal_index) * psi[l];
547  }
548 
549  return (interpolated_u);
550  }
551 
552  /// Return the i-th component of the FE interpolated pressure
553  /// p[i] at local coordinate s
555  const unsigned& i) const
556  {
557  // Determine number of pressure nodes in the element
558  const unsigned n_pressure_nodes = npres_linearised_nst();
559 
560  // Provide storage for local shape functions
561  Shape psi(n_pressure_nodes);
562 
563  // Find values of shape functions
564  pshape_linearised_nst(s, psi);
565 
566  // Initialise value of p
567  double interpolated_p = 0.0;
568 
569  // Loop over the local nodes and sum
570  for (unsigned l = 0; l < n_pressure_nodes; l++)
571  {
572  // N.B. The pure virtual function p_linearised_nst(...)
573  // automatically calculates the index at which the pressure value
574  // is stored, so we don't need to worry about this here
575  interpolated_p += p_linearised_nst(l, i) * psi[l];
576  }
577 
578  return (interpolated_p);
579  }
580 
581  }; // End of LinearisedNavierStokesEquations class definition
582 
583 
584  /// ///////////////////////////////////////////////////////////////////////////
585  /// ///////////////////////////////////////////////////////////////////////////
586  /// ///////////////////////////////////////////////////////////////////////////
587 
588 
589  //=======================================================================
590  /// Crouzeix-Raviart elements are Navier-Stokes elements with quadratic
591  /// interpolation for velocities and positions, but a discontinuous
592  /// linear pressure interpolation
593  //=======================================================================
595  : public virtual QElement<2, 3>,
596  public virtual LinearisedNavierStokesEquations
597  {
598  private:
599  /// Static array of ints to hold required number of variables at nodes
600  static const unsigned Initial_Nvalue[];
601 
602  protected:
603  /// Internal indices that indicate at which internal data the
604  /// pressure values are stored. We note that there are two pressure
605  /// values, corresponding to the functions P^C(r,z,t) and P^S(r,z,t)
606  /// which multiply the cosine and sine terms respectively.
608 
609  /// Velocity shape and test functions and their derivatives
610  /// w.r.t. global coordinates at local coordinate s (taken from geometry).
611  /// Return Jacobian of mapping between local and global coordinates.
613  const Vector<double>& s,
614  Shape& psi,
615  DShape& dpsidx,
616  Shape& test,
617  DShape& dtestdx) const;
618 
619  /// Velocity shape and test functions and their derivatives
620  /// w.r.t. global coordinates at the ipt-th integation point
621  /// (taken from geometry).
622  /// Return Jacobian of mapping between local and global coordinates.
624  const unsigned& ipt,
625  Shape& psi,
626  DShape& dpsidx,
627  Shape& test,
628  DShape& dtestdx) const;
629 
630  /// Compute the pressure shape functions at local coordinate s
631  inline void pshape_linearised_nst(const Vector<double>& s,
632  Shape& psi) const;
633 
634  /// Compute the pressure shape and test functions at local coordinate s
635  inline void pshape_linearised_nst(const Vector<double>& s,
636  Shape& psi,
637  Shape& test) const;
638 
639  public:
640  /// Constructor: there are three internal values for each
641  /// of the two pressure components
643  : QElement<2, 3>(),
646  {
647  // Loop over the two pressure components
648  // and two normalisation constraints
649  for (unsigned i = 0; i < 4; i++)
650  {
651  // Allocate and add one internal data object for each of the two
652  // pressure components that store the three pressure values
654  this->add_internal_data(new Data(3));
655  }
656  }
657 
658  /// Return number of values (pinned or dofs) required at local node n
659  virtual unsigned required_nvalue(const unsigned& n) const;
660 
661  /// Return the pressure value i at internal dof i_internal
662  /// (Discontinous pressure interpolation -- no need to cater for
663  /// hanging nodes)
664  double p_linearised_nst(const unsigned& i_internal, const unsigned& i) const
665  {
667  ->value(i_internal);
668  }
669 
670  // Pin the normalisation dofs
672  {
673  for (unsigned i = 2; i < 4; i++)
674  {
675  this->internal_data_pt(P_linearised_nst_internal_index[i])->pin_all();
676  }
677  }
678 
679  void pin_real_or_imag(const unsigned& real_index)
680  {
681  unsigned n_node = this->nnode();
682  for (unsigned n = 0; n < n_node; n++)
683  {
684  Node* nod_pt = this->node_pt(n);
685 
686  for (unsigned i = 0; i < DIM; ++i)
687  {
688  // Provided it's not constrained then pin it
689  if (!nod_pt->is_constrained(i))
690  {
691  this->node_pt(n)->pin(2 * i + real_index);
692  }
693  }
694  }
695 
696  // Similarly for the pressure
697  this->internal_data_pt(P_linearised_nst_internal_index[real_index])
698  ->pin_all();
699  }
700 
701  void unpin_real_or_imag(const unsigned& real_index)
702  {
703  unsigned n_node = this->nnode();
704  for (unsigned n = 0; n < n_node; n++)
705  {
706  Node* nod_pt = this->node_pt(n);
707 
708  for (unsigned i = 0; i < DIM; ++i)
709  {
710  // Provided it's not constrained then unpin it
711  if (!nod_pt->is_constrained(i))
712  {
713  nod_pt->unpin(2 * i + real_index);
714  }
715  }
716  }
717 
718  // Similarly for the pressure
719  this->internal_data_pt(P_linearised_nst_internal_index[real_index])
720  ->unpin_all();
721  }
722 
724  {
725  unsigned n_node = this->nnode();
726  for (unsigned n = 0; n < n_node; n++)
727  {
728  Node* nod_pt = this->node_pt(n);
729 
730  // Transfer the eigenfunctions to the normalisation constraints
731  for (unsigned i = 0; i < DIM; ++i)
732  {
733  for (unsigned j = 0; j < 2; ++j)
734  {
735  nod_pt->set_value(2 * (DIM + i) + j, nod_pt->value(2 * i + j));
736  }
737  }
738  }
739 
740  // Similarly for the pressure
741  for (unsigned i = 0; i < 2; ++i)
742  {
743  Data* local_data_pt =
744  this->internal_data_pt(P_linearised_nst_internal_index[i]);
745  Data* norm_local_data_pt =
746  this->internal_data_pt(P_linearised_nst_internal_index[2 + i]);
747  for (unsigned j = 0; j < 3; j++)
748  {
749  norm_local_data_pt->set_value(j, local_data_pt->value(j));
750  }
751  }
752  }
753 
754 
755  /// Return number of pressure values corresponding to a
756  /// single pressure component
757  unsigned npres_linearised_nst() const
758  {
759  return 3;
760  }
761 
762  /// Fix both components of the internal pressure degrees
763  /// of freedom p_dof to pvalue
764  void fix_pressure(const unsigned& p_dof, const double& pvalue)
765  {
766  // Loop over the two pressure components
767  for (unsigned i = 0; i < 2; i++)
768  {
769  this->internal_data_pt(P_linearised_nst_internal_index[i])->pin(p_dof);
771  ->set_value(p_dof, pvalue);
772  }
773  }
774 
775  /// Overload the access function for the i-th component of the
776  /// pressure's local equation numbers
777  inline int p_local_eqn(const unsigned& n, const unsigned& i)
778  {
780  }
781 
782  /// Redirect output to NavierStokesEquations output
783  void output(std::ostream& outfile)
784  {
786  }
787 
788  /// Redirect output to NavierStokesEquations output
789  void output(std::ostream& outfile, const unsigned& n_plot)
790  {
792  }
793 
794  /// Redirect output to NavierStokesEquations output
795  void output(FILE* file_pt)
796  {
798  }
799 
800  /// Redirect output to NavierStokesEquations output
801  void output(FILE* file_pt, const unsigned& n_plot)
802  {
804  }
805 
806  /// The number of "dof-blocks" that degrees of freedom in this
807  /// element are sub-divided into: Velocity and pressure.
808  unsigned ndof_types() const
809  {
810  return 2 * (DIM + 1);
811  }
812 
813  }; // End of LinearisedQCrouzeixRaviartElement class definition
814 
815 
816  // Inline functions
817  // ----------------
818 
819  //=======================================================================
820  /// Derivatives of the shape functions and test functions w.r.t.
821  /// global (Eulerian) coordinates at local coordinate s.
822  /// Return Jacobian of mapping between local and global coordinates.
823  //=======================================================================
826  Shape& psi,
827  DShape& dpsidx,
828  Shape& test,
829  DShape& dtestdx) const
830  {
831  // Call the geometrical shape functions and derivatives
832  const double J = this->dshape_eulerian(s, psi, dpsidx);
833  // The test functions are equal to the shape functions
834  test = psi;
835  dtestdx = dpsidx;
836  // Return the Jacobian
837  return J;
838  }
839 
840  //=======================================================================
841  /// Derivatives of the shape functions and test functions w.r.t.
842  /// global (Eulerian) coordinates at the ipt-th integration point.
843  /// Return Jacobian of mapping between local and global coordinates.
844  //=======================================================================
847  Shape& psi,
848  DShape& dpsidx,
849  Shape& test,
850  DShape& dtestdx) const
851  {
852  // Call the geometrical shape functions and derivatives
853  const double J = this->dshape_eulerian_at_knot(ipt, psi, dpsidx);
854 
855  // Loop over the test functions and derivatives and set them
856  // equal to the shape functions
857  test = psi;
858  dtestdx = dpsidx;
859  // Return the Jacobian
860  return J;
861  }
862 
863  //=======================================================================
864  /// Pressure shape functions
865  //=======================================================================
867  const Vector<double>& s, Shape& psi) const
868  {
869  psi[0] = 1.0;
870  psi[1] = s[0];
871  psi[2] = s[1];
872  }
873 
874  //=======================================================================
875  /// Define the pressure shape and test functions
876  //=======================================================================
878  const Vector<double>& s, Shape& psi, Shape& test) const
879  {
880  // Call the pressure shape functions
881  pshape_linearised_nst(s, psi);
882 
883  // Loop over the test functions and set them equal to the shape functions
884  for (unsigned i = 0; i < 3; i++)
885  {
886  test[i] = psi[i];
887  }
888  }
889 
890  //=======================================================================
891  /// Face geometry of the linearised axisym Crouzeix-Raviart elements
892  //=======================================================================
893  template<>
895  : public virtual QElement<1, 3>
896  {
897  public:
898  FaceGeometry() : QElement<1, 3>() {}
899  };
900 
901  //=======================================================================
902  /// Face geometry of face geometry of the linearised axisymmetric
903  /// Crouzeix Raviart elements
904  //=======================================================================
905  template<>
907  : public virtual PointElement
908  {
909  public:
911  };
912 
913 
914  /// ///////////////////////////////////////////////////////////////////////////
915  /// ///////////////////////////////////////////////////////////////////////////
916  /// ///////////////////////////////////////////////////////////////////////////
917 
918 
919  //=======================================================================
920  /// Taylor--Hood elements are Navier--Stokes elements with quadratic
921  /// interpolation for velocities and positions and continuous linear
922  /// pressure interpolation
923  //=======================================================================
925  : public virtual QElement<2, 3>,
926  public virtual LinearisedNavierStokesEquations
927  {
928  private:
929  /// Static array of ints to hold number of variables at node
930  static const unsigned Initial_Nvalue[];
931 
932  protected:
933  /// Static array of ints to hold conversion from pressure
934  /// node numbers to actual node numbers
935  static const unsigned Pconv[];
936 
937  /// Velocity shape and test functions and their derivatives
938  /// w.r.t. global coordinates at local coordinate s (taken from geometry).
939  /// Return Jacobian of mapping between local and global coordinates.
941  const Vector<double>& s,
942  Shape& psi,
943  DShape& dpsidx,
944  Shape& test,
945  DShape& dtestdx) const;
946 
947  /// Velocity shape and test functions and their derivatives
948  /// w.r.t. global coordinates the ipt-th integation point
949  /// (taken from geometry).
950  /// Return Jacobian of mapping between local and global coordinates.
952  const unsigned& ipt,
953  Shape& psi,
954  DShape& dpsidx,
955  Shape& test,
956  DShape& dtestdx) const;
957 
958  /// Compute the pressure shape functions at local coordinate s
959  inline void pshape_linearised_nst(const Vector<double>& s,
960  Shape& psi) const;
961 
962  /// Compute the pressure shape and test functions at local coordinte s
963  inline void pshape_linearised_nst(const Vector<double>& s,
964  Shape& psi,
965  Shape& test) const;
966 
967  public:
968  /// Constructor, no internal data points
971  {
972  }
973 
974  /// Number of values (pinned or dofs) required at node n. Can
975  /// be overwritten for hanging node version
976  inline virtual unsigned required_nvalue(const unsigned& n) const
977  {
978  return Initial_Nvalue[n];
979  }
980 
981  /// Which nodal value represents the pressure? Overload version
982  /// in base class which returns static int "Pressure_not_stored_at_node"
983  virtual int p_index_linearised_nst(const unsigned& i) const
984  {
985  return (2 * DIM + i);
986  }
987 
988  /// Access function for the i-th component of pressure
989  /// at local pressure node n_p (const version).
990  double p_linearised_nst(const unsigned& n_p, const unsigned& i) const
991  {
992  return nodal_value(Pconv[n_p], p_index_linearised_nst(i));
993  }
994 
995  // Pin the normalisation dofs
997  {
998  throw OomphLibError("This is not implemented yet\n",
999  OOMPH_CURRENT_FUNCTION,
1000  OOMPH_EXCEPTION_LOCATION);
1001  }
1002 
1003 
1004  virtual void pin_real_or_imag(const unsigned& real)
1005  {
1006  throw OomphLibError("This is not implemented yet\n",
1007  OOMPH_CURRENT_FUNCTION,
1008  OOMPH_EXCEPTION_LOCATION);
1009  }
1010 
1011  virtual void unpin_real_or_imag(const unsigned& real)
1012  {
1013  throw OomphLibError("This is not implemented yet\n",
1014  OOMPH_CURRENT_FUNCTION,
1015  OOMPH_EXCEPTION_LOCATION);
1016  }
1017 
1018 
1019  /// Return number of pressure values corresponding to a
1020  /// single pressure component
1021  unsigned npres_linearised_nst() const
1022  {
1023  return 4;
1024  }
1025 
1026  /// Fix both components of the pressure at local pressure
1027  /// node n_p to pvalue
1028  void fix_pressure(const unsigned& n_p, const double& pvalue)
1029  {
1030  // Loop over the two pressure components
1031  for (unsigned i = 0; i < 2; i++)
1032  {
1033  this->node_pt(Pconv[n_p])->pin(p_index_linearised_nst(i));
1034  this->node_pt(Pconv[n_p])->set_value(p_index_linearised_nst(i), pvalue);
1035  }
1036  }
1037 
1038  /// Overload the access function for the i-th component of the
1039  /// pressure's local equation numbers
1040  inline int p_local_eqn(const unsigned& n, const unsigned& i)
1041  {
1043  }
1044 
1045  /// Redirect output to NavierStokesEquations output
1046  void output(std::ostream& outfile)
1047  {
1049  }
1050 
1051  /// Redirect output to NavierStokesEquations output
1052  void output(std::ostream& outfile, const unsigned& n_plot)
1053  {
1054  LinearisedNavierStokesEquations::output(outfile, n_plot);
1055  }
1056 
1057  /// Redirect output to NavierStokesEquations output
1058  void output(FILE* file_pt)
1059  {
1061  }
1062 
1063  /// Redirect output to NavierStokesEquations output
1064  void output(FILE* file_pt, const unsigned& n_plot)
1065  {
1066  LinearisedNavierStokesEquations::output(file_pt, n_plot);
1067  }
1068 
1069  /// Returns the number of "dof-blocks" that degrees of freedom
1070  /// in this element are sub-divided into: Velocity and pressure.
1071  unsigned ndof_types() const
1072  {
1073  return 8;
1074  }
1075 
1076  }; // End of LinearisedQTaylorHoodElement class definition
1077 
1078 
1079  // Inline functions
1080  // ----------------
1081 
1082  //=======================================================================
1083  /// Derivatives of the shape functions and test functions w.r.t
1084  /// global (Eulerian) coordinates at local coordinate s.
1085  /// Return Jacobian of mapping between local and global coordinates.
1086  //=======================================================================
1089  Shape& psi,
1090  DShape& dpsidx,
1091  Shape& test,
1092  DShape& dtestdx) const
1093  {
1094  // Call the geometrical shape functions and derivatives
1095  const double J = this->dshape_eulerian(s, psi, dpsidx);
1096 
1097  test = psi;
1098  dtestdx = dpsidx;
1099 
1100  // Return the Jacobian
1101  return J;
1102  }
1103 
1104  //=======================================================================
1105  /// Derivatives of the shape functions and test functions w.r.t
1106  /// global (Eulerian) coordinates at the ipt-th integration point.
1107  /// Return Jacobian of mapping between local and global coordinates.
1108  //=======================================================================
1111  Shape& psi,
1112  DShape& dpsidx,
1113  Shape& test,
1114  DShape& dtestdx) const
1115  {
1116  // Call the geometrical shape functions and derivatives
1117  const double J = this->dshape_eulerian_at_knot(ipt, psi, dpsidx);
1118 
1119  test = psi;
1120  dtestdx = dpsidx;
1121 
1122  // Return the Jacobian
1123  return J;
1124  }
1125 
1126  //=======================================================================
1127  /// Pressure shape functions
1128  //=======================================================================
1130  const Vector<double>& s, Shape& psi) const
1131  {
1132  // Allocate local storage for the pressure shape functions
1133  double psi1[2], psi2[2];
1134 
1135  // Call the one-dimensional shape functions
1136  OneDimLagrange::shape<2>(s[0], psi1);
1137  OneDimLagrange::shape<2>(s[1], psi2);
1138 
1139  // Now let's loop over the nodal points in the element
1140  // s1 is the "r" coordinate, s2 the "z"
1141  for (unsigned i = 0; i < 2; i++)
1142  {
1143  for (unsigned j = 0; j < 2; j++)
1144  {
1145  // Multiply the two 1D functions together to get the 2D function
1146  psi[2 * i + j] = psi2[i] * psi1[j];
1147  }
1148  }
1149  }
1150 
1151  //=======================================================================
1152  /// Pressure shape and test functions
1153  //=======================================================================
1155  const Vector<double>& s, Shape& psi, Shape& test) const
1156  {
1157  // Call the pressure shape functions
1158  pshape_linearised_nst(s, psi);
1159 
1160  // Loop over the test functions and set them equal to the shape functions
1161  for (unsigned i = 0; i < 4; i++)
1162  {
1163  test[i] = psi[i];
1164  }
1165  }
1166 
1167  //=======================================================================
1168  /// Face geometry of the linearised axisymmetric Taylor Hood elements
1169  //=======================================================================
1170  template<>
1172  : public virtual QElement<1, 3>
1173  {
1174  public:
1175  FaceGeometry() : QElement<1, 3>() {}
1176  };
1177 
1178  //=======================================================================
1179  /// Face geometry of the face geometry of the linearised
1180  /// axisymmetric Taylor Hood elements
1181  //=======================================================================
1182  template<>
1184  : public virtual PointElement
1185  {
1186  public:
1188  };
1189 
1190 
1191 } // namespace oomph
1192 
1193 #endif
static char t char * s
Definition: cfortran.h:568
cstr elem_len * i
Definition: cfortran.h:603
char t
Definition: cfortran.h:568
A Class for the derivatives of shape functions The class design is essentially the same as Shape,...
Definition: shape.h:278
A class that represents a collection of data; each Data object may contain many different individual ...
Definition: nodes.h:86
void pin_all()
Pin all the stored variables.
Definition: nodes.h:397
void pin(const unsigned &i)
Pin the i-th stored variable.
Definition: nodes.h:385
void unpin(const unsigned &i)
Unpin the i-th stored variable.
Definition: nodes.h:391
double * value_pt(const unsigned &i) const
Return the pointer to the i-the stored value. Typically this is required when direct access to the st...
Definition: nodes.h:324
void unpin_all()
Unpin all the stored variables.
Definition: nodes.h:407
TimeStepper *& time_stepper_pt()
Return the pointer to the timestepper.
Definition: nodes.h:238
void set_value(const unsigned &i, const double &value_)
Set the i-th stored data value to specified value. The only reason that we require an explicit set fu...
Definition: nodes.h:271
double value(const unsigned &i) const
Return i-th stored value. This function is not virtual so that it can be inlined. This means that if ...
Definition: nodes.h:293
bool is_constrained(const unsigned &i)
Test whether the i-th variable is constrained (1: true; 0: false).
Definition: nodes.h:472
//////////////////////////////////////////////////////////////////// ////////////////////////////////...
Definition: elements.h:4998
A general Finite Element class.
Definition: elements.h:1313
Node *& node_pt(const unsigned &n)
Return a pointer to the local node n.
Definition: elements.h:2175
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:2593
virtual double dshape_eulerian_at_knot(const unsigned &ipt, Shape &psi, DShape &dpsidx) const
Return the geometric shape functions and also first derivatives w.r.t. global coordinates at the ipt-...
Definition: elements.cc:3325
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...
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:1432
unsigned nnode() const
Return the number of nodes.
Definition: elements.h:2210
double dshape_eulerian(const Vector< double > &s, Shape &psi, DShape &dpsidx) const
Compute the geometric shape functions and also first derivatives w.r.t. global coordinates at local c...
Definition: elements.cc:3298
unsigned add_internal_data(Data *const &data_pt, const bool &fd=true)
Add a (pointer to an) internal data object to the element and return the index required to obtain it ...
Definition: elements.cc:62
Data *& internal_data_pt(const unsigned &i)
Return a pointer to i-th internal data object.
Definition: elements.h:622
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
int internal_local_eqn(const unsigned &i, const unsigned &j) const
Return the local equation number corresponding to the j-th value stored at the i-th internal data.
Definition: elements.h:267
int external_local_eqn(const unsigned &i, const unsigned &j)
Return the local equation number corresponding to the j-th value stored at the i-th external data.
Definition: elements.h:311
unsigned add_external_data(Data *const &data_pt, const bool &fd=true)
Add a (pointer to an) external data object to the element and return its index (i....
Definition: elements.cc:307
TimeStepper *& time_stepper_pt()
Access function for pointer to time stepper: Null if object is not time-dependent.
Definition: geom_objects.h:192
A class that is used to implement the constraint that the eigenfunction has a particular normalisatio...
Data * eigenvalue_data_pt()
Access to Data that contains the traded pressure.
unsigned index_of_eigenvalue()
Return the index of Data object at which the traded pressure is stored.
A class for elements that solve the linearised version of the unsteady Navier–Stokes equations in cyl...
LinearisedNavierStokesEquations()
Constructor: NULL the base flow solution and the derivatives of the base flow function.
virtual unsigned u_index_linearised_nst(const unsigned &i) const
Return the index at which the i-th unknown velocity component is stored. The default value,...
double *& re_pt()
Pointer to Reynolds number.
virtual void pshape_linearised_nst(const Vector< double > &s, Shape &psi, Shape &test) const =0
Compute the pressure shape and test functions at local coordinate s.
LinearisedNavierStokesEigenfunctionNormalisationElement * normalisation_element_pt()
Pointer to normalisation element.
virtual void pin_real_or_imag(const unsigned &real)=0
Pin the real or imaginary part of the problem Input integer 0 for real 1 for imaginary.
void fill_in_contribution_to_residuals(Vector< double > &residuals)
Compute the element's residual Vector.
static Vector< double > Gamma
Vector to decide whether the stress-divergence form is used or not.
virtual unsigned npres_linearised_nst() const =0
Return the number of pressure degrees of freedom associated with a single pressure component in the e...
void output(FILE *file_pt)
Output function: r, z, U^C, U^S, V^C, V^S, W^C, W^S, P^C, P^S in tecplot format. Default number of pl...
void(* Base_flow_dudx_fct_pt)(const double &time, const Vector< double > &x, DenseMatrix< double > &result)
Pointer to derivatives of base flow solution velocity components w.r.t. global coordinates (r and z) ...
unsigned Data_number_of_eigenvalue
Index of datum where eigenvalue is stored.
virtual int p_local_eqn(const unsigned &n, const unsigned &i)=0
Access function for the local equation number information for the i-th component of the pressure....
LinearisedNavierStokesEigenfunctionNormalisationElement * Normalisation_element_pt
Pointer to the normalisation element.
const double & re_st() const
Product of Reynolds and Strouhal number (=Womersley number)
void strain_rate(const Vector< double > &s, DenseMatrix< double > &strain_rate, const unsigned &real) const
Strain-rate tensor: where (in that order)
void output_veloc(std::ostream &outfile, const unsigned &nplot, const unsigned &t)
Output function: r, z, U^C, U^S, V^C, V^S, W^C, W^S, in tecplot format. nplot points in each coordina...
void(*&)(const double &time, const Vector< double > &x, DenseMatrix< double > &f) base_flow_dudx_fct_pt()
Access function for the derivatives of the base flow w.r.t. global coordinates solution pointer.
virtual void pshape_linearised_nst(const Vector< double > &s, Shape &psi) const =0
Compute the pressure shape functions at local coordinate s.
bool ALE_is_disabled
Boolean flag to indicate if ALE formulation is disabled when the time-derivatives are computed....
double du_dt_linearised_nst(const unsigned &n, const unsigned &i) const
Return the i-th component of du/dt at local node n. Uses suitably interpolated value for hanging node...
void set_eigenfunction_normalisation_element(LinearisedNavierStokesEigenfunctionNormalisationElement *const &normalisation_el_pt)
the boolean flag check_nodal_data is set to false.
virtual void get_base_flow_u(const double &time, const unsigned &ipt, const Vector< double > &x, Vector< double > &result) const
Calculate the velocity components of the base flow solution at a given time and Eulerian position.
virtual double p_linearised_nst(const unsigned &n_p, const unsigned &i) const =0
Return the i-th pressure value at local pressure "node" n_p. Uses suitably interpolated value for han...
static int Pressure_not_stored_at_node
Static "magic" number that indicates that the pressure is not stored at a node.
virtual void pin_pressure_normalisation_dofs()=0
Pin the normalisation dofs.
virtual double dshape_and_dtest_eulerian_at_knot_linearised_nst(const unsigned &ipt, Shape &psi, DShape &dpsidx, Shape &test, DShape &dtestdx) const =0
Compute the shape functions and their derivatives w.r.t. global coordinates at the ipt-th integration...
void output(std::ostream &outfile)
Output function: r, z, U^C, U^S, V^C, V^S, W^C, W^S, P^C, P^S in tecplot format. Default number of pl...
double * Viscosity_Ratio_pt
Pointer to the viscosity ratio (relative to the viscosity used in the definition of the Reynolds numb...
void disable_ALE()
Disable ALE, i.e. assert the mesh is not moving – you do this at your own risk!
void(*&)(const double &time, const Vector< double > &x, Vector< double > &f) base_flow_u_fct_pt()
Access function for the base flow solution pointer.
void(* Base_flow_u_fct_pt)(const double &time, const Vector< double > &x, Vector< double > &result)
Pointer to base flow solution (velocity components) function.
virtual void fill_in_generic_residual_contribution_linearised_nst(Vector< double > &residuals, DenseMatrix< double > &jacobian, DenseMatrix< double > &mass_matrix, unsigned flag)
Compute the residuals for the Navier-Stokes equations; flag=1(or 0): do (or don't) compute the Jacobi...
void enable_ALE()
(Re-)enable ALE, i.e. take possible mesh motion into account when evaluating the time-derivative....
double interpolated_p_linearised_nst(const Vector< double > &s, const unsigned &i) const
Return the i-th component of the FE interpolated pressure p[i] at local coordinate s.
double * ReSt_pt
Pointer to global Reynolds number x Strouhal number (=Womersley)
virtual void get_base_flow_dudx(const double &time, const unsigned &ipt, const Vector< double > &x, DenseMatrix< double > &result) const
Calculate the derivatives of the velocity components of the base flow solution w.r....
double * Re_pt
Pointer to global Reynolds number.
double *& viscosity_ratio_pt()
Pointer to the viscosity ratio.
double *& re_st_pt()
Pointer to product of Reynolds and Strouhal number (=Womersley number)
double *& density_ratio_pt()
Pointer to the density ratio.
double * Density_Ratio_pt
Pointer to the density ratio (relative to the density used in the definition of the Reynolds number)
static double Default_Physical_Constant_Value
Static default value for the physical constants (all initialised to zero)
const double & viscosity_ratio() const
Viscosity ratio for element: element's viscosity relative to the viscosity used in the definition of ...
double interpolated_u_linearised_nst(const Vector< double > &s, const unsigned &i) const
Compute the element's residual Vector and the jacobian matrix. Virtual function can be overloaded by ...
virtual void unpin_real_or_imag(const unsigned &real)=0
virtual int p_index_linearised_nst(const unsigned &i) const
Which nodal value represents the pressure?
static double Default_Physical_Ratio_Value
Static default value for the physical ratios (all initialised to one)
const double & density_ratio() const
Density ratio for element: element's density relative to the viscosity used in the definition of the ...
virtual double dshape_and_dtest_eulerian_linearised_nst(const Vector< double > &s, Shape &psi, DShape &dpsidx, Shape &test, DShape &dtestdx) const =0
Compute the shape functions and their derivatives w.r.t. global coordinates at local coordinate s....
/////////////////////////////////////////////////////////////////////////// /////////////////////////...
void fix_pressure(const unsigned &p_dof, const double &pvalue)
Fix both components of the internal pressure degrees of freedom p_dof to pvalue.
void output(std::ostream &outfile)
Redirect output to NavierStokesEquations output.
double p_linearised_nst(const unsigned &i_internal, const unsigned &i) const
Return the pressure value i at internal dof i_internal (Discontinous pressure interpolation – no need...
void output(FILE *file_pt)
Redirect output to NavierStokesEquations output.
LinearisedQCrouzeixRaviartElement()
Constructor: there are three internal values for each of the two pressure components.
void output(FILE *file_pt, const unsigned &n_plot)
Redirect output to NavierStokesEquations output.
double dshape_and_dtest_eulerian_linearised_nst(const Vector< double > &s, Shape &psi, DShape &dpsidx, Shape &test, DShape &dtestdx) const
Velocity shape and test functions and their derivatives w.r.t. global coordinates at local coordinate...
static const unsigned Initial_Nvalue[]
Static array of ints to hold required number of variables at nodes.
void pin_real_or_imag(const unsigned &real_index)
Pin the real or imaginary part of the problem Input integer 0 for real 1 for imaginary.
unsigned npres_linearised_nst() const
Return number of pressure values corresponding to a single pressure component.
Vector< unsigned > P_linearised_nst_internal_index
Internal indices that indicate at which internal data the pressure values are stored....
void pshape_linearised_nst(const Vector< double > &s, Shape &psi) const
Compute the pressure shape functions at local coordinate s.
int p_local_eqn(const unsigned &n, const unsigned &i)
Overload the access function for the i-th component of the pressure's local equation numbers.
void output(std::ostream &outfile, const unsigned &n_plot)
Redirect output to NavierStokesEquations output.
void pin_pressure_normalisation_dofs()
Pin the normalisation dofs.
double dshape_and_dtest_eulerian_at_knot_linearised_nst(const unsigned &ipt, Shape &psi, DShape &dpsidx, Shape &test, DShape &dtestdx) const
Velocity shape and test functions and their derivatives w.r.t. global coordinates at the ipt-th integ...
unsigned ndof_types() const
The number of "dof-blocks" that degrees of freedom in this element are sub-divided into: Velocity and...
virtual unsigned required_nvalue(const unsigned &n) const
Return number of values (pinned or dofs) required at local node n.
/////////////////////////////////////////////////////////////////////////// /////////////////////////...
double dshape_and_dtest_eulerian_at_knot_linearised_nst(const unsigned &ipt, Shape &psi, DShape &dpsidx, Shape &test, DShape &dtestdx) const
Velocity shape and test functions and their derivatives w.r.t. global coordinates the ipt-th integati...
void fix_pressure(const unsigned &n_p, const double &pvalue)
Fix both components of the pressure at local pressure node n_p to pvalue.
void output(std::ostream &outfile)
Redirect output to NavierStokesEquations output.
void pin_pressure_normalisation_dofs()
Pin the normalisation dofs.
void output(FILE *file_pt, const unsigned &n_plot)
Redirect output to NavierStokesEquations output.
unsigned ndof_types() const
Returns the number of "dof-blocks" that degrees of freedom in this element are sub-divided into: Velo...
void output(std::ostream &outfile, const unsigned &n_plot)
Redirect output to NavierStokesEquations output.
unsigned npres_linearised_nst() const
Return number of pressure values corresponding to a single pressure component.
static const unsigned Pconv[]
Static array of ints to hold conversion from pressure node numbers to actual node numbers.
virtual void pin_real_or_imag(const unsigned &real)
Pin the real or imaginary part of the problem Input integer 0 for real 1 for imaginary.
static const unsigned Initial_Nvalue[]
Static array of ints to hold number of variables at node.
void output(FILE *file_pt)
Redirect output to NavierStokesEquations output.
void pshape_linearised_nst(const Vector< double > &s, Shape &psi) const
Compute the pressure shape functions at local coordinate s.
virtual int p_index_linearised_nst(const unsigned &i) const
Which nodal value represents the pressure? Overload version in base class which returns static int "P...
double p_linearised_nst(const unsigned &n_p, const unsigned &i) const
Access function for the i-th component of pressure at local pressure node n_p (const version).
LinearisedQTaylorHoodElement()
Constructor, no internal data points.
virtual unsigned required_nvalue(const unsigned &n) const
Number of values (pinned or dofs) required at node n. Can be overwritten for hanging node version.
int p_local_eqn(const unsigned &n, const unsigned &i)
Overload the access function for the i-th component of the pressure's local equation numbers.
double dshape_and_dtest_eulerian_linearised_nst(const Vector< double > &s, Shape &psi, DShape &dpsidx, Shape &test, DShape &dtestdx) const
Velocity shape and test functions and their derivatives w.r.t. global coordinates at local coordinate...
virtual void unpin_real_or_imag(const unsigned &real)
Nodes are derived from Data, but, in addition, have a definite (Eulerian) position in a space of a gi...
Definition: nodes.h:906
double value(const unsigned &i) const
Return i-th value (dofs or pinned) at this node either directly or via hanging node representation....
Definition: nodes.cc:2408
An OomphLibError object which should be thrown when an run-time error is encountered....
/////////////////////////////////////////////////////////////////////// /////////////////////////////...
Definition: elements.h:3439
/////////////////////////////////////////////////////////////////////// /////////////////////////////...
Definition: Qelements.h:459
A Class for shape functions. In simple cases, the shape functions have only one index that can be tho...
Definition: shape.h:76
////////////////////////////////////////////////////////////////////// //////////////////////////////...
Definition: timesteppers.h:231
unsigned ntstorage() const
Return the number of doubles required to represent history (one for steady)
Definition: timesteppers.h:601
virtual double weight(const unsigned &i, const unsigned &j) const
Access function for j-th weight for the i-th derivative.
Definition: timesteppers.h:594
bool is_steady() const
Flag to indicate if a timestepper has been made steady (possibly temporarily to switch off time-depen...
Definition: timesteppers.h:389
void shape< 2 >(const double &s, double *Psi)
1D shape functions specialised to linear order (2 Nodes)
Definition: shape.h:608
//////////////////////////////////////////////////////////////////// ////////////////////////////////...