oomph::SingularNavierStokesSolutionElement< WRAPPED_NAVIER_STOKES_ELEMENT > Class Template Reference

#include <navier_stokes_elements_with_singularity.h>

+ Inheritance diagram for oomph::SingularNavierStokesSolutionElement< WRAPPED_NAVIER_STOKES_ELEMENT >:

Public Types

typedef Vector< double >(* NavierStokesVelocitySingularFctPt) (const Vector< double > &x)
 Function pointer to the velocity singular function: More...
 
typedef Vector< Vector< double > >(* NavierStokesGradVelocitySingularFctPt) (const Vector< double > &x)
 Function pointer to the gradient of the velocity singular function: More...
 
typedef double(* NavierStokesPressureSingularFctPt) (const Vector< double > &x)
 Function pointer to the pressure singular function: More...
 
typedef Vector< double >(* NavierStokesGradPressureSingularFctPt) (const Vector< double > &x)
 Function pointer to the gradient of the pressure singular function: More...
 

Public Member Functions

 SingularNavierStokesSolutionElement ()
 Constructor. More...
 
boolsingular_function_satisfies_stokes_equation ()
 
double c () const
 Find the value of the unknown amplitude C. More...
 
void set_c (const double &value)
 Find the value of the unknown amplitude C. More...
 
void pin_c ()
 Pin the value of the unknown amplitude C. More...
 
void unpin_c ()
 Unpin the value of the unknown amplitude C. More...
 
void set_wrapped_navier_stokes_element_pt (WRAPPED_NAVIER_STOKES_ELEMENT *wrapped_navier_stokes_el_pt, const Vector< double > &s, unsigned *direction_pt)
 
NavierStokesVelocitySingularFctPtvelocity_singular_fct_pt ()
 Access function to pointer to velocity singular function. More...
 
NavierStokesGradVelocitySingularFctPtgrad_velocity_singular_fct_pt ()
 Access function to pointer to gradient of velocity singular function. More...
 
NavierStokesPressureSingularFctPtpressure_singular_fct_pt ()
 Access function to pointer to pressure singular function. More...
 
NavierStokesGradPressureSingularFctPtgrad_pressure_singular_fct_pt ()
 Access function to pointer to gradient of pressure singular function. More...
 
Vector< doublevelocity_singular_function (const Vector< double > &x) const
 Evaluate velocity singular function at Eulerian position x. More...
 
Vector< Vector< double > > grad_velocity_singular_function (const Vector< double > &x) const
 
double pressure_singular_function (const Vector< double > &x) const
 Evaluate pressure singular function at Eulerian position x. More...
 
Vector< doublegrad_pressure_singular_function (const Vector< double > &x) const
 Evaluate gradient of pressure singular function at Eulerian position x. More...
 
Vector< doubleu_bar (const Vector< double > &x)
 
Vector< Vector< double > > grad_u_bar (const Vector< double > &x)
 
double p_bar (const Vector< double > &x)
 
Vector< doublegrad_p_bar (const Vector< double > &x)
 
void fill_in_contribution_to_residuals (Vector< double > &residual)
 Compute residual. More...
 
void fill_in_contribution_to_jacobian (Vector< double > &residual, DenseMatrix< double > &jacobian)
 
- Public Member Functions inherited from oomph::GeneralisedElement
 GeneralisedElement ()
 Constructor: Initialise all pointers and all values to zero. More...
 
virtual ~GeneralisedElement ()
 Virtual destructor to clean up any memory allocated by the object. More...
 
 GeneralisedElement (const GeneralisedElement &)=delete
 Broken copy constructor. More...
 
void operator= (const GeneralisedElement &)=delete
 Broken assignment operator. More...
 
Data *& internal_data_pt (const unsigned &i)
 Return a pointer to i-th internal data object. More...
 
Data *const & internal_data_pt (const unsigned &i) const
 Return a pointer to i-th internal data object (const version) More...
 
Data *& external_data_pt (const unsigned &i)
 Return a pointer to i-th external data object. More...
 
Data *const & external_data_pt (const unsigned &i) const
 Return a pointer to i-th external data object (const version) More...
 
unsigned long eqn_number (const unsigned &ieqn_local) const
 
int local_eqn_number (const unsigned long &ieqn_global) const
 
unsigned add_external_data (Data *const &data_pt, const bool &fd=true)
 
bool external_data_fd (const unsigned &i) const
 
void exclude_external_data_fd (const unsigned &i)
 
void include_external_data_fd (const unsigned &i)
 
void flush_external_data ()
 Flush all external data. More...
 
void flush_external_data (Data *const &data_pt)
 Flush the object addressed by data_pt from the external data array. More...
 
unsigned ninternal_data () const
 Return the number of internal data objects. More...
 
unsigned nexternal_data () const
 Return the number of external data objects. More...
 
unsigned ndof () const
 Return the number of equations/dofs in the element. More...
 
void dof_vector (const unsigned &t, Vector< double > &dof)
 Return the vector of dof values at time level t. More...
 
void dof_pt_vector (Vector< double * > &dof_pt)
 Return the vector of pointers to dof values. More...
 
void set_internal_data_time_stepper (const unsigned &i, TimeStepper *const &time_stepper_pt, const bool &preserve_existing_data)
 
void assign_internal_eqn_numbers (unsigned long &global_number, Vector< double * > &Dof_pt)
 
void describe_dofs (std::ostream &out, const std::string &current_string) const
 
virtual void describe_local_dofs (std::ostream &out, const std::string &current_string) const
 
void add_internal_value_pt_to_map (std::map< unsigned, double * > &map_of_value_pt)
 
virtual void assign_local_eqn_numbers (const bool &store_local_dof_pt)
 
virtual void complete_setup_of_dependencies ()
 
virtual void get_residuals (Vector< double > &residuals)
 
virtual void get_jacobian (Vector< double > &residuals, DenseMatrix< double > &jacobian)
 
virtual void get_mass_matrix (Vector< double > &residuals, DenseMatrix< double > &mass_matrix)
 
virtual void get_jacobian_and_mass_matrix (Vector< double > &residuals, DenseMatrix< double > &jacobian, DenseMatrix< double > &mass_matrix)
 
virtual void get_dresiduals_dparameter (double *const &parameter_pt, Vector< double > &dres_dparam)
 
virtual void get_djacobian_dparameter (double *const &parameter_pt, Vector< double > &dres_dparam, DenseMatrix< double > &djac_dparam)
 
virtual void get_djacobian_and_dmass_matrix_dparameter (double *const &parameter_pt, Vector< double > &dres_dparam, DenseMatrix< double > &djac_dparam, DenseMatrix< double > &dmass_matrix_dparam)
 
virtual void get_hessian_vector_products (Vector< double > const &Y, DenseMatrix< double > const &C, DenseMatrix< double > &product)
 
virtual void get_inner_products (Vector< std::pair< unsigned, unsigned >> const &history_index, Vector< double > &inner_product)
 
virtual void get_inner_product_vectors (Vector< unsigned > const &history_index, Vector< Vector< double >> &inner_product_vector)
 
virtual unsigned self_test ()
 
virtual void compute_norm (Vector< double > &norm)
 
virtual void compute_norm (double &norm)
 
virtual unsigned ndof_types () const
 
virtual void get_dof_numbers_for_unknowns (std::list< std::pair< unsigned long, unsigned >> &dof_lookup_list) const
 

Private Member Functions

void fill_in_generic_contribution_to_residuals (Vector< double > &residual, DenseMatrix< double > &jacobian, const unsigned &flag)
 Compute local residual, and, if flag=1, local jacobian matrix. More...
 

Private Attributes

WRAPPED_NAVIER_STOKES_ELEMENT * Wrapped_navier_stokes_el_pt
 Pointer to wrapped Navier-Stokes element. More...
 
NavierStokesVelocitySingularFctPt Velocity_singular_fct_pt
 Pointer to velocity singular function. More...
 
NavierStokesGradVelocitySingularFctPt Grad_velocity_singular_fct_pt
 
NavierStokesPressureSingularFctPt Pressure_singular_fct_pt
 Pointer to pressure singular function. More...
 
NavierStokesGradPressureSingularFctPt Grad_pressure_singular_fct_pt
 Pointer to gradient of pressure singular function. More...
 
Vector< doubleS_in_wrapped_navier_stokes_element
 Local coordinates of singulariity in wrapped Navier-Stokes element. More...
 
unsignedDirection_pt
 Direction of the derivative used for the residual of the element. More...
 
bool Singular_function_satisfies_stokes_equation
 

Additional Inherited Members

- Static Public Attributes inherited from oomph::GeneralisedElement
static bool Suppress_warning_about_repeated_internal_data
 
static bool Suppress_warning_about_repeated_external_data = true
 
static double Default_fd_jacobian_step = 1.0e-8
 
- Protected Member Functions inherited from oomph::GeneralisedElement
unsigned add_internal_data (Data *const &data_pt, const bool &fd=true)
 
bool internal_data_fd (const unsigned &i) const
 
void exclude_internal_data_fd (const unsigned &i)
 
void include_internal_data_fd (const unsigned &i)
 
void clear_global_eqn_numbers ()
 
void add_global_eqn_numbers (std::deque< unsigned long > const &global_eqn_numbers, std::deque< double * > const &global_dof_pt)
 
virtual void assign_internal_and_external_local_eqn_numbers (const bool &store_local_dof_pt)
 
virtual void assign_all_generic_local_eqn_numbers (const bool &store_local_dof_pt)
 
virtual void assign_additional_local_eqn_numbers ()
 
int internal_local_eqn (const unsigned &i, const unsigned &j) const
 
int external_local_eqn (const unsigned &i, const unsigned &j)
 
void fill_in_jacobian_from_internal_by_fd (Vector< double > &residuals, DenseMatrix< double > &jacobian, const bool &fd_all_data=false)
 
void fill_in_jacobian_from_internal_by_fd (DenseMatrix< double > &jacobian, const bool &fd_all_data=false)
 
void fill_in_jacobian_from_external_by_fd (Vector< double > &residuals, DenseMatrix< double > &jacobian, const bool &fd_all_data=false)
 
void fill_in_jacobian_from_external_by_fd (DenseMatrix< double > &jacobian, const bool &fd_all_data=false)
 
virtual void update_before_internal_fd ()
 
virtual void reset_after_internal_fd ()
 
virtual void update_in_internal_fd (const unsigned &i)
 
virtual void reset_in_internal_fd (const unsigned &i)
 
virtual void update_before_external_fd ()
 
virtual void reset_after_external_fd ()
 
virtual void update_in_external_fd (const unsigned &i)
 
virtual void reset_in_external_fd (const unsigned &i)
 
virtual void fill_in_contribution_to_mass_matrix (Vector< double > &residuals, DenseMatrix< double > &mass_matrix)
 
virtual void fill_in_contribution_to_jacobian_and_mass_matrix (Vector< double > &residuals, DenseMatrix< double > &jacobian, DenseMatrix< double > &mass_matrix)
 
virtual void fill_in_contribution_to_dresiduals_dparameter (double *const &parameter_pt, Vector< double > &dres_dparam)
 
virtual void fill_in_contribution_to_djacobian_dparameter (double *const &parameter_pt, Vector< double > &dres_dparam, DenseMatrix< double > &djac_dparam)
 
virtual void fill_in_contribution_to_djacobian_and_dmass_matrix_dparameter (double *const &parameter_pt, Vector< double > &dres_dparam, DenseMatrix< double > &djac_dparam, DenseMatrix< double > &dmass_matrix_dparam)
 
virtual void fill_in_contribution_to_hessian_vector_products (Vector< double > const &Y, DenseMatrix< double > const &C, DenseMatrix< double > &product)
 
virtual void fill_in_contribution_to_inner_products (Vector< std::pair< unsigned, unsigned >> const &history_index, Vector< double > &inner_product)
 
virtual void fill_in_contribution_to_inner_product_vectors (Vector< unsigned > const &history_index, Vector< Vector< double >> &inner_product_vector)
 
- Static Protected Attributes inherited from oomph::GeneralisedElement
static DenseMatrix< doubleDummy_matrix
 
static std::deque< double * > Dof_pt_deque
 

Detailed Description

template<class WRAPPED_NAVIER_STOKES_ELEMENT>
class oomph::SingularNavierStokesSolutionElement< WRAPPED_NAVIER_STOKES_ELEMENT >

We consider a singularity in the solution at the point O.

This class defines the singular functions.

velocity_singularity:

u_bar = C*velocity_singular_function

pressure_singularity:

p_bar = C*pressure_singular_function

and their gradients

The class also defines the function that computes the residual associated to C which is:

R_C = \frac{\partial p_FE}{\partial x_i} (O)

and thus regularises the FE solution by setting the pressure gradient in the coordinate direction x_i to zero. If the amplitude of the singular solution is known, pin C.

Member Typedef Documentation

◆ NavierStokesGradPressureSingularFctPt

template<class WRAPPED_NAVIER_STOKES_ELEMENT >
typedef Vector<double>(* oomph::SingularNavierStokesSolutionElement< WRAPPED_NAVIER_STOKES_ELEMENT >::NavierStokesGradPressureSingularFctPt) (const Vector< double > &x)

Function pointer to the gradient of the pressure singular function:

◆ NavierStokesGradVelocitySingularFctPt

template<class WRAPPED_NAVIER_STOKES_ELEMENT >
typedef Vector<Vector<double> >(* oomph::SingularNavierStokesSolutionElement< WRAPPED_NAVIER_STOKES_ELEMENT >::NavierStokesGradVelocitySingularFctPt) (const Vector< double > &x)

Function pointer to the gradient of the velocity singular function:

◆ NavierStokesPressureSingularFctPt

template<class WRAPPED_NAVIER_STOKES_ELEMENT >
typedef double(* oomph::SingularNavierStokesSolutionElement< WRAPPED_NAVIER_STOKES_ELEMENT >::NavierStokesPressureSingularFctPt) (const Vector< double > &x)

Function pointer to the pressure singular function:

◆ NavierStokesVelocitySingularFctPt

template<class WRAPPED_NAVIER_STOKES_ELEMENT >
typedef Vector<double>(* oomph::SingularNavierStokesSolutionElement< WRAPPED_NAVIER_STOKES_ELEMENT >::NavierStokesVelocitySingularFctPt) (const Vector< double > &x)

Function pointer to the velocity singular function:

Constructor & Destructor Documentation

◆ SingularNavierStokesSolutionElement()

template<class WRAPPED_NAVIER_STOKES_ELEMENT >
oomph::SingularNavierStokesSolutionElement< WRAPPED_NAVIER_STOKES_ELEMENT >::SingularNavierStokesSolutionElement ( )
inline

Constructor.

86  {
87  // Initialise Function pointer to velocity singular function to NULL
89 
90  // Initialise Function pointer to gradient of velocity singular
91  // function to NULL
93 
94  // Initialise Function pointer to pressure singular function to NULL
96 
97  // Initialise Function pointer to gradient of pressure singular
98  // function to NULL
100 
101  // Initalise pointer to the wrapped Navier-Stokes element which will be
102  // used to compute the residual and which includes the point O
104 
105  // Initialise the pointer to the direction of the derivative used
106  // for the residual of this element
107  Direction_pt=0;
108 
109  // Safe assumption: Singular fct does not satisfy Stokes eqn
111 
112  // Create a single item of internal Data, storing one unknown which
113  // represents the unknown C.
114  add_internal_data(new Data(1));
115  }
unsigned add_internal_data(Data *const &data_pt, const bool &fd=true)
Definition: elements.cc:62
NavierStokesPressureSingularFctPt Pressure_singular_fct_pt
Pointer to pressure singular function.
Definition: navier_stokes_elements_with_singularity.h:453
unsigned * Direction_pt
Direction of the derivative used for the residual of the element.
Definition: navier_stokes_elements_with_singularity.h:462
NavierStokesGradPressureSingularFctPt Grad_pressure_singular_fct_pt
Pointer to gradient of pressure singular function.
Definition: navier_stokes_elements_with_singularity.h:456
NavierStokesVelocitySingularFctPt Velocity_singular_fct_pt
Pointer to velocity singular function.
Definition: navier_stokes_elements_with_singularity.h:446
NavierStokesGradVelocitySingularFctPt Grad_velocity_singular_fct_pt
Definition: navier_stokes_elements_with_singularity.h:450
bool Singular_function_satisfies_stokes_equation
Definition: navier_stokes_elements_with_singularity.h:465
WRAPPED_NAVIER_STOKES_ELEMENT * Wrapped_navier_stokes_el_pt
Pointer to wrapped Navier-Stokes element.
Definition: navier_stokes_elements_with_singularity.h:443

References oomph::GeneralisedElement::add_internal_data(), oomph::SingularNavierStokesSolutionElement< WRAPPED_NAVIER_STOKES_ELEMENT >::Direction_pt, oomph::SingularNavierStokesSolutionElement< WRAPPED_NAVIER_STOKES_ELEMENT >::Grad_pressure_singular_fct_pt, oomph::SingularNavierStokesSolutionElement< WRAPPED_NAVIER_STOKES_ELEMENT >::Grad_velocity_singular_fct_pt, oomph::SingularNavierStokesSolutionElement< WRAPPED_NAVIER_STOKES_ELEMENT >::Pressure_singular_fct_pt, oomph::SingularNavierStokesSolutionElement< WRAPPED_NAVIER_STOKES_ELEMENT >::Singular_function_satisfies_stokes_equation, oomph::SingularNavierStokesSolutionElement< WRAPPED_NAVIER_STOKES_ELEMENT >::Velocity_singular_fct_pt, and oomph::SingularNavierStokesSolutionElement< WRAPPED_NAVIER_STOKES_ELEMENT >::Wrapped_navier_stokes_el_pt.

Member Function Documentation

◆ c()

template<class WRAPPED_NAVIER_STOKES_ELEMENT >
double oomph::SingularNavierStokesSolutionElement< WRAPPED_NAVIER_STOKES_ELEMENT >::c ( ) const
inline

◆ fill_in_contribution_to_jacobian()

template<class WRAPPED_NAVIER_STOKES_ELEMENT >
void oomph::SingularNavierStokesSolutionElement< WRAPPED_NAVIER_STOKES_ELEMENT >::fill_in_contribution_to_jacobian ( Vector< double > &  residuals,
DenseMatrix< double > &  jacobian 
)
inlinevirtual

Add the elemental contribution to the jacobian matrix. and the residuals vector. Note that this function will NOT initialise the residuals vector or the jacobian matrix. It must be called after the residuals vector and jacobian matrix have been initialised to zero. The default is to use finite differences to calculate the jacobian

Reimplemented from oomph::GeneralisedElement.

384  {
385  fill_in_generic_contribution_to_residuals(residual,jacobian,1);
386  }
void fill_in_generic_contribution_to_residuals(Vector< double > &residual, DenseMatrix< double > &jacobian, const unsigned &flag)
Compute local residual, and, if flag=1, local jacobian matrix.
Definition: navier_stokes_elements_with_singularity.h:391

References oomph::SingularNavierStokesSolutionElement< WRAPPED_NAVIER_STOKES_ELEMENT >::fill_in_generic_contribution_to_residuals().

◆ fill_in_contribution_to_residuals()

template<class WRAPPED_NAVIER_STOKES_ELEMENT >
void oomph::SingularNavierStokesSolutionElement< WRAPPED_NAVIER_STOKES_ELEMENT >::fill_in_contribution_to_residuals ( Vector< double > &  residual)
inlinevirtual

◆ fill_in_generic_contribution_to_residuals()

template<class WRAPPED_NAVIER_STOKES_ELEMENT >
void oomph::SingularNavierStokesSolutionElement< WRAPPED_NAVIER_STOKES_ELEMENT >::fill_in_generic_contribution_to_residuals ( Vector< double > &  residual,
DenseMatrix< double > &  jacobian,
const unsigned flag 
)
inlineprivate

Compute local residual, and, if flag=1, local jacobian matrix.

394  {
395  // Get the local eqn number of our one-and-only
396  // unknown
398  if (eqn_number>=0)
399  {
400  residual[eqn_number] = Wrapped_navier_stokes_el_pt->dpdx_fe_only
402 
403  // Do we want the Jacobian too?
404  if (flag)
405  {
406  // Find the number of pressure dofs in the wrapped Navier-Stokes
407  // element pointed by
408  // the SingularNavierStokesSolutionElement class
409  unsigned n_pres = Wrapped_navier_stokes_el_pt->npres_nst();
410 
411  // Find the dimension of the problem
412  unsigned cached_dim=Wrapped_navier_stokes_el_pt->dim();
413 
414  // Set up memory for the pressure shape functions and their derivatives
415  Shape psip(n_pres), testp(n_pres);
416  DShape dpsipdx(n_pres,cached_dim), dtestpdx(n_pres,cached_dim);
417 
418  // Compute the pressure shape functions and their derivatives
419  // at the local coordinate S_in_wrapped_navier_stokes_element
420  // (Test fcts not really needed but nobody's got around to writing
421  // a fct that only picks out the basis fcts.
422  Wrapped_navier_stokes_el_pt->dpshape_and_dptest_eulerian_nst
423  (S_in_wrapped_navier_stokes_element,psip,dpsipdx,testp,dtestpdx);
424 
425  // Derivs
426  for (unsigned j=0;j<n_pres;j++)
427  {
428  // Unknown
429  int local_unknown = Wrapped_navier_stokes_el_pt->p_local_eqn(j);
430 
431  // If not pinned
432  if (local_unknown >= 0)
433  {
434  // Add the contribution of the node to the local jacobian
435  jacobian(eqn_number,local_unknown) = dpsipdx(j,*Direction_pt);
436  }
437  }
438  }
439  }
440  }
unsigned long eqn_number(const unsigned &ieqn_local) const
Definition: elements.h:704
int internal_local_eqn(const unsigned &i, const unsigned &j) const
Definition: elements.h:267
Vector< double > S_in_wrapped_navier_stokes_element
Local coordinates of singulariity in wrapped Navier-Stokes element.
Definition: navier_stokes_elements_with_singularity.h:459
std::ptrdiff_t j
Definition: tut_arithmetic_redux_minmax.cpp:2

References oomph::SingularNavierStokesSolutionElement< WRAPPED_NAVIER_STOKES_ELEMENT >::Direction_pt, oomph::GeneralisedElement::eqn_number(), oomph::GeneralisedElement::internal_local_eqn(), j, oomph::SingularNavierStokesSolutionElement< WRAPPED_NAVIER_STOKES_ELEMENT >::S_in_wrapped_navier_stokes_element, and oomph::SingularNavierStokesSolutionElement< WRAPPED_NAVIER_STOKES_ELEMENT >::Wrapped_navier_stokes_el_pt.

Referenced by oomph::SingularNavierStokesSolutionElement< WRAPPED_NAVIER_STOKES_ELEMENT >::fill_in_contribution_to_jacobian(), and oomph::SingularNavierStokesSolutionElement< WRAPPED_NAVIER_STOKES_ELEMENT >::fill_in_contribution_to_residuals().

◆ grad_p_bar()

template<class WRAPPED_NAVIER_STOKES_ELEMENT >
Vector<double> oomph::SingularNavierStokesSolutionElement< WRAPPED_NAVIER_STOKES_ELEMENT >::grad_p_bar ( const Vector< double > &  x)
inline

Evaluate gradient of pressure singular function (including its amplitude): grad_p_bar = C * grad_pressure_singular

354  {
355  // Find the value of C
356  double c=internal_data_pt(0)->value(0);
357 
358  // Initialise the gradient of pressure
359  Vector<double> grad_p=grad_pressure_singular_function(x);
360 
361  // Find the dimension of the problem
362  unsigned cached_dim=Wrapped_navier_stokes_el_pt->dim();
363 
364  // Multiply the components of the gradient of pressure by the unknown C
365  for (unsigned d=0;d<cached_dim;d++)
366  {
367  grad_p[d] *= c;
368  }
369 
370  // Value of grad_p_bar at the position x
371  return grad_p;
372  }
Vector< double > grad_pressure_singular_function(const Vector< double > &x) const
Evaluate gradient of pressure singular function at Eulerian position x.
Definition: navier_stokes_elements_with_singularity.h:267
double c() const
Find the value of the unknown amplitude C.
Definition: navier_stokes_elements_with_singularity.h:125
list x
Definition: plotDoE.py:28

References oomph::SingularNavierStokesSolutionElement< WRAPPED_NAVIER_STOKES_ELEMENT >::c(), oomph::SingularNavierStokesSolutionElement< WRAPPED_NAVIER_STOKES_ELEMENT >::grad_pressure_singular_function(), oomph::GeneralisedElement::internal_data_pt(), oomph::Data::value(), oomph::SingularNavierStokesSolutionElement< WRAPPED_NAVIER_STOKES_ELEMENT >::Wrapped_navier_stokes_el_pt, and plotDoE::x.

◆ grad_pressure_singular_fct_pt()

template<class WRAPPED_NAVIER_STOKES_ELEMENT >
NavierStokesGradPressureSingularFctPt& oomph::SingularNavierStokesSolutionElement< WRAPPED_NAVIER_STOKES_ELEMENT >::grad_pressure_singular_fct_pt ( )
inline

Access function to pointer to gradient of pressure singular function.

References oomph::SingularNavierStokesSolutionElement< WRAPPED_NAVIER_STOKES_ELEMENT >::Grad_pressure_singular_fct_pt.

◆ grad_pressure_singular_function()

template<class WRAPPED_NAVIER_STOKES_ELEMENT >
Vector<double> oomph::SingularNavierStokesSolutionElement< WRAPPED_NAVIER_STOKES_ELEMENT >::grad_pressure_singular_function ( const Vector< double > &  x) const
inline

Evaluate gradient of pressure singular function at Eulerian position x.

268  {
269 #ifdef PARANOID
271  {
272  std::stringstream error_stream;
273  error_stream
274  << "Pointer to gradient of pressure singular function "
275  << "hasn't been defined!"
276  << std::endl;
277  throw OomphLibError(
278  error_stream.str(),
281  }
282 #endif
283 
284  // Evaluate gradient of pressure singular function
285  return (*Grad_pressure_singular_fct_pt)(x);
286  }
#define OOMPH_EXCEPTION_LOCATION
Definition: oomph_definitions.h:61
#define OOMPH_CURRENT_FUNCTION
Definition: oomph_definitions.h:86

References oomph::SingularNavierStokesSolutionElement< WRAPPED_NAVIER_STOKES_ELEMENT >::Grad_pressure_singular_fct_pt, OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION, and plotDoE::x.

Referenced by oomph::SingularNavierStokesSolutionElement< WRAPPED_NAVIER_STOKES_ELEMENT >::grad_p_bar().

◆ grad_u_bar()

template<class WRAPPED_NAVIER_STOKES_ELEMENT >
Vector<Vector<double> > oomph::SingularNavierStokesSolutionElement< WRAPPED_NAVIER_STOKES_ELEMENT >::grad_u_bar ( const Vector< double > &  x)
inline

Evaluate gradient of velocity singular function (including its amplitude): grad_u_bar = C * grad_velocity_singular; grad[i][j] = du_i/dx_j

316  {
317  // Find the value of C
318  double c=internal_data_pt(0)->value(0);
319 
320  // Initialise the gradient of velocity vector
321  Vector<Vector<double> > grad_u=grad_velocity_singular_function(x);
322 
323  // Find the dimension of the problem
324  unsigned cached_dim=Wrapped_navier_stokes_el_pt->dim();
325 
326  // Multiply the components of the gradient of velocity by the unknown C
327  for (unsigned d=0;d<cached_dim;d++)
328  {
329  for (unsigned i=0;i<cached_dim;i++)
330  {
331  grad_u[d][i] *= c;
332  }
333  }
334 
335  // Value of grad_u_bar at the position x
336  return grad_u;
337  }
int i
Definition: BiCGSTAB_step_by_step.cpp:9
Vector< Vector< double > > grad_velocity_singular_function(const Vector< double > &x) const
Definition: navier_stokes_elements_with_singularity.h:223

References oomph::SingularNavierStokesSolutionElement< WRAPPED_NAVIER_STOKES_ELEMENT >::c(), oomph::SingularNavierStokesSolutionElement< WRAPPED_NAVIER_STOKES_ELEMENT >::grad_velocity_singular_function(), i, oomph::GeneralisedElement::internal_data_pt(), oomph::Data::value(), oomph::SingularNavierStokesSolutionElement< WRAPPED_NAVIER_STOKES_ELEMENT >::Wrapped_navier_stokes_el_pt, and plotDoE::x.

◆ grad_velocity_singular_fct_pt()

template<class WRAPPED_NAVIER_STOKES_ELEMENT >
NavierStokesGradVelocitySingularFctPt& oomph::SingularNavierStokesSolutionElement< WRAPPED_NAVIER_STOKES_ELEMENT >::grad_velocity_singular_fct_pt ( )
inline

Access function to pointer to gradient of velocity singular function.

References oomph::SingularNavierStokesSolutionElement< WRAPPED_NAVIER_STOKES_ELEMENT >::Grad_velocity_singular_fct_pt.

◆ grad_velocity_singular_function()

template<class WRAPPED_NAVIER_STOKES_ELEMENT >
Vector<Vector<double> > oomph::SingularNavierStokesSolutionElement< WRAPPED_NAVIER_STOKES_ELEMENT >::grad_velocity_singular_function ( const Vector< double > &  x) const
inline

Evaluate gradient of velocity singular function at Eulerian position x. grad[i][j] = du_i/dx_j

225  {
226 #ifdef PARANOID
228  {
229  std::stringstream error_stream;
230  error_stream
231  << "Pointer to gradient of velocity singular function "
232  << "hasn't been defined!"
233  << std::endl;
234  throw OomphLibError(
235  error_stream.str(),
238  }
239 #endif
240 
241  // Evaluate gradient of velocity singular function
242  return (*Grad_velocity_singular_fct_pt)(x);
243  }

References oomph::SingularNavierStokesSolutionElement< WRAPPED_NAVIER_STOKES_ELEMENT >::Grad_velocity_singular_fct_pt, OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION, and plotDoE::x.

Referenced by oomph::SingularNavierStokesSolutionElement< WRAPPED_NAVIER_STOKES_ELEMENT >::grad_u_bar().

◆ p_bar()

template<class WRAPPED_NAVIER_STOKES_ELEMENT >
double oomph::SingularNavierStokesSolutionElement< WRAPPED_NAVIER_STOKES_ELEMENT >::p_bar ( const Vector< double > &  x)
inline

Evaluate pressure singular function (including its amplitude): p_bar = C * pressure_singular

342  {
343  // Find the value of C
344  double c=internal_data_pt(0)->value(0);
345 
346  // Value of p_bar at the position x
348  }
double pressure_singular_function(const Vector< double > &x) const
Evaluate pressure singular function at Eulerian position x.
Definition: navier_stokes_elements_with_singularity.h:246

References oomph::SingularNavierStokesSolutionElement< WRAPPED_NAVIER_STOKES_ELEMENT >::c(), oomph::GeneralisedElement::internal_data_pt(), oomph::SingularNavierStokesSolutionElement< WRAPPED_NAVIER_STOKES_ELEMENT >::pressure_singular_function(), oomph::Data::value(), and plotDoE::x.

◆ pin_c()

template<class WRAPPED_NAVIER_STOKES_ELEMENT >
void oomph::SingularNavierStokesSolutionElement< WRAPPED_NAVIER_STOKES_ELEMENT >::pin_c ( )
inline

Pin the value of the unknown amplitude C.

139  {
140  return internal_data_pt(0)->pin(0);
141  }
void pin(const unsigned &i)
Pin the i-th stored variable.
Definition: nodes.h:385

References oomph::GeneralisedElement::internal_data_pt(), and oomph::Data::pin().

◆ pressure_singular_fct_pt()

template<class WRAPPED_NAVIER_STOKES_ELEMENT >
NavierStokesPressureSingularFctPt& oomph::SingularNavierStokesSolutionElement< WRAPPED_NAVIER_STOKES_ELEMENT >::pressure_singular_fct_pt ( )
inline

Access function to pointer to pressure singular function.

194  {return Pressure_singular_fct_pt;}

References oomph::SingularNavierStokesSolutionElement< WRAPPED_NAVIER_STOKES_ELEMENT >::Pressure_singular_fct_pt.

◆ pressure_singular_function()

template<class WRAPPED_NAVIER_STOKES_ELEMENT >
double oomph::SingularNavierStokesSolutionElement< WRAPPED_NAVIER_STOKES_ELEMENT >::pressure_singular_function ( const Vector< double > &  x) const
inline

Evaluate pressure singular function at Eulerian position x.

247  {
248 #ifdef PARANOID
250  {
251  std::stringstream error_stream;
252  error_stream
253  << "Pointer to pressure singular function hasn't been defined!"
254  << std::endl;
255  throw OomphLibError(
256  error_stream.str(),
259  }
260 #endif
261 
262  // Evaluate pressure singular function
263  return (*Pressure_singular_fct_pt)(x);
264  }

References OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION, oomph::SingularNavierStokesSolutionElement< WRAPPED_NAVIER_STOKES_ELEMENT >::Pressure_singular_fct_pt, and plotDoE::x.

Referenced by oomph::SingularNavierStokesSolutionElement< WRAPPED_NAVIER_STOKES_ELEMENT >::p_bar().

◆ set_c()

template<class WRAPPED_NAVIER_STOKES_ELEMENT >
void oomph::SingularNavierStokesSolutionElement< WRAPPED_NAVIER_STOKES_ELEMENT >::set_c ( const double value)
inline

Find the value of the unknown amplitude C.

133  {
135  }
void set_value(const unsigned &i, const double &value_)
Definition: nodes.h:271
squared absolute value
Definition: GlobalFunctions.h:87

References oomph::GeneralisedElement::internal_data_pt(), oomph::Data::set_value(), and Eigen::value.

◆ set_wrapped_navier_stokes_element_pt()

template<class WRAPPED_NAVIER_STOKES_ELEMENT >
void oomph::SingularNavierStokesSolutionElement< WRAPPED_NAVIER_STOKES_ELEMENT >::set_wrapped_navier_stokes_element_pt ( WRAPPED_NAVIER_STOKES_ELEMENT *  wrapped_navier_stokes_el_pt,
const Vector< double > &  s,
unsigned direction_pt 
)
inline

Set pointer to associated wrapped Navier-Stokes element which contains the singularity (at local coordinate s). Also specify the direction in which the slope of the FE part of the pressure is set to zero. (Could also set a velocity derivative to zero but this needs to be done with a separate function. Write it if you need it...)

158  {
159  // Assign the pointer to the variable Wrapped_navier_stokes_el_pt
160  Wrapped_navier_stokes_el_pt = wrapped_navier_stokes_el_pt;
161 
162  // Find number of nodes in the element
163  unsigned nnod=wrapped_navier_stokes_el_pt->nnode();
164 
165  // Loop over the nodes of the element
166  for (unsigned j=0;j<nnod;j++)
167  {
168  // Add the node as external data in the
169  // SingularNavierStokesSolutionElement class. Note that this
170  // assumes that the pressure is stored at the nodes (Taylor Hood type
171  // NSt elements, which is assumed elsewhere too...)
173  }
174 
175  // Assign the pointer to the local coordinate at which the residual
176  // will be computed
178 
179  // Assign the pointer to the direction at which the derivative used
180  // in the residual will be computed
181  Direction_pt = direction_pt;
182  }
unsigned add_external_data(Data *const &data_pt, const bool &fd=true)
Definition: elements.cc:307
RealScalar s
Definition: level1_cplx_impl.h:130

References oomph::GeneralisedElement::add_external_data(), oomph::SingularNavierStokesSolutionElement< WRAPPED_NAVIER_STOKES_ELEMENT >::Direction_pt, j, s, oomph::SingularNavierStokesSolutionElement< WRAPPED_NAVIER_STOKES_ELEMENT >::S_in_wrapped_navier_stokes_element, and oomph::SingularNavierStokesSolutionElement< WRAPPED_NAVIER_STOKES_ELEMENT >::Wrapped_navier_stokes_el_pt.

◆ singular_function_satisfies_stokes_equation()

template<class WRAPPED_NAVIER_STOKES_ELEMENT >
bool& oomph::SingularNavierStokesSolutionElement< WRAPPED_NAVIER_STOKES_ELEMENT >::singular_function_satisfies_stokes_equation ( )
inline

Assert that singular function satisfies the Stokes equations by setting this to true or false.

120  {
122  }

References oomph::SingularNavierStokesSolutionElement< WRAPPED_NAVIER_STOKES_ELEMENT >::Singular_function_satisfies_stokes_equation.

◆ u_bar()

template<class WRAPPED_NAVIER_STOKES_ELEMENT >
Vector<double> oomph::SingularNavierStokesSolutionElement< WRAPPED_NAVIER_STOKES_ELEMENT >::u_bar ( const Vector< double > &  x)
inline

Evaluate velocity singular function (including its amplitude): u_bar = C * velocity_singular

291  {
292  // Find the value of C
293  double c=internal_data_pt(0)->value(0);
294 
295  // Initialise the velocity vector
296  Vector<double> u=velocity_singular_function(x);
297 
298  // Find the dimension of the problem
299  unsigned cached_dim=Wrapped_navier_stokes_el_pt->dim();
300 
301  // Multiply the components of the velocity vector by the unknown C
302  for (unsigned d=0;d<cached_dim;d++)
303  {
304  u[d] *= c;
305  }
306 
307  // Value of u_bar at the position x
308  return u;
309  }
Vector< double > velocity_singular_function(const Vector< double > &x) const
Evaluate velocity singular function at Eulerian position x.
Definition: navier_stokes_elements_with_singularity.h:201

References oomph::SingularNavierStokesSolutionElement< WRAPPED_NAVIER_STOKES_ELEMENT >::c(), oomph::GeneralisedElement::internal_data_pt(), oomph::Data::value(), oomph::SingularNavierStokesSolutionElement< WRAPPED_NAVIER_STOKES_ELEMENT >::velocity_singular_function(), oomph::SingularNavierStokesSolutionElement< WRAPPED_NAVIER_STOKES_ELEMENT >::Wrapped_navier_stokes_el_pt, and plotDoE::x.

◆ unpin_c()

template<class WRAPPED_NAVIER_STOKES_ELEMENT >
void oomph::SingularNavierStokesSolutionElement< WRAPPED_NAVIER_STOKES_ELEMENT >::unpin_c ( )
inline

Unpin the value of the unknown amplitude C.

145  {
146  return internal_data_pt(0)->unpin(0);
147  }
void unpin(const unsigned &i)
Unpin the i-th stored variable.
Definition: nodes.h:391

References oomph::GeneralisedElement::internal_data_pt(), and oomph::Data::unpin().

◆ velocity_singular_fct_pt()

template<class WRAPPED_NAVIER_STOKES_ELEMENT >
NavierStokesVelocitySingularFctPt& oomph::SingularNavierStokesSolutionElement< WRAPPED_NAVIER_STOKES_ELEMENT >::velocity_singular_fct_pt ( )
inline

Access function to pointer to velocity singular function.

186  {return Velocity_singular_fct_pt;}

References oomph::SingularNavierStokesSolutionElement< WRAPPED_NAVIER_STOKES_ELEMENT >::Velocity_singular_fct_pt.

◆ velocity_singular_function()

template<class WRAPPED_NAVIER_STOKES_ELEMENT >
Vector<double> oomph::SingularNavierStokesSolutionElement< WRAPPED_NAVIER_STOKES_ELEMENT >::velocity_singular_function ( const Vector< double > &  x) const
inline

Evaluate velocity singular function at Eulerian position x.

202  {
203 #ifdef PARANOID
205  {
206  std::stringstream error_stream;
207  error_stream
208  << "Pointer to velocity singular function hasn't been defined!"
209  << std::endl;
210  throw OomphLibError(
211  error_stream.str(),
214  }
215 #endif
216 
217  // Evaluate velocity singular function
218  return (*Velocity_singular_fct_pt)(x);
219  }

References OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION, oomph::SingularNavierStokesSolutionElement< WRAPPED_NAVIER_STOKES_ELEMENT >::Velocity_singular_fct_pt, and plotDoE::x.

Referenced by oomph::SingularNavierStokesSolutionElement< WRAPPED_NAVIER_STOKES_ELEMENT >::u_bar().

Member Data Documentation

◆ Direction_pt

◆ Grad_pressure_singular_fct_pt

◆ Grad_velocity_singular_fct_pt

◆ Pressure_singular_fct_pt

◆ S_in_wrapped_navier_stokes_element

template<class WRAPPED_NAVIER_STOKES_ELEMENT >
Vector<double> oomph::SingularNavierStokesSolutionElement< WRAPPED_NAVIER_STOKES_ELEMENT >::S_in_wrapped_navier_stokes_element
private

◆ Singular_function_satisfies_stokes_equation

◆ Velocity_singular_fct_pt

◆ Wrapped_navier_stokes_el_pt


The documentation for this class was generated from the following file: