///////////////////////////////////////////////////////////////////// ///////////////////////////////////////////////////////////////////// More...
#include <womersley_elements.h>
Public Types | |
typedef double(* | PrescribedVolumeFluxFctPt) (const double &time) |
Function pointer to fct that prescribes volume flux q=fct(t) – mainly used for validation purposes. More... | |
Public Member Functions | |
WomersleyImpedanceTubeBase (const double &length, PrescribedVolumeFluxFctPt prescribed_volume_flux_fct_pt) | |
Constructor: Specify length of tube and pointer to function that specifies the prescribed volume flux. Outlet pressure is set to zero. More... | |
WomersleyImpedanceTubeBase (const double &length, Mesh *navier_stokes_outflow_mesh_pt) | |
Constructor: Specify length of tube and the pointer to the mesh of either NavierStokesImpedanceTractionElements or NavierStokesFluxControlElements that are attached to the outflow cross-section of a (higher-dimensional) Navier Stokes mesh and provide the inflow into the ImpedanceTube. Outlet pressure is set to zero. More... | |
double & | p_out () |
Access fct to outlet pressure. More... | |
virtual Mesh * | build_mesh_and_apply_boundary_conditions (TimeStepper *time_stepper_pt)=0 |
Pure virtual function in which the user of a derived class must create the mesh of WomersleyElements (of the type specified by the class's template argument) and apply the boundary conditions. The Womersley elements use the timestepper specified as the input argument. More... | |
void | setup () |
Set up the Womersley tubes so that a subsequent call to get_response(...) computes the inlet pressure for the currently prescribed instantaneous flow rate. Steady version! More... | |
void | setup (double *re_st_pt, const double &dt, const double &q_initial, TimeStepper *time_stepper_pt=0) |
Set up the Womersley tubes so that a subsequent call to get_response(...) computes the inlet pressure for the currently prescribed instantaneous flow rate, assuming that at all previous times the tube conveyed steady, fully-developed flow with flowrate q_initial. dt specifies the timestep for the subsequent time integration. Specify: Womersley number, (constant) timestep, the initial volume flux (from which the subsequent impulsive start is performed) and, optionally the pointer to the timestepper to be used in the Womersley elements (defaults to BDF<2>). More... | |
WomersleyProblem< ELEMENT, DIM > * | womersley_problem_pt () |
Access to underlying Womersley problem. More... | |
void | shift_time_values (const double &dt) |
Shift history values to allow coputation of next timestep. Note: When used with a full Navier-Stokes problem this function must be called in actions_before_implicit_timestep() More... | |
double | total_volume_flux_into_impedance_tube () |
Compute total current volume flux into the "impedance tube" that provides the flow resistance (flux is either obtained from the function that specifies it externally or by by adding up the flux through all NavierStokesImpedanceTractionElements in the mesh pointed to by the Navier_stokes_outflow_mesh_pt. More... | |
void | get_response (double &p_in, double &dp_in_dq) |
Compute inlet pressure, p_in, required to achieve the currently imposed, instantaneous volume flux q prescribed by total_volume_flux_into_impedance_tube(), and its derivative, dp_in/dq. More... | |
Public Member Functions inherited from oomph::TemplateFreeWomersleyImpedanceTubeBase | |
TemplateFreeWomersleyImpedanceTubeBase () | |
Empty constructor. More... | |
virtual | ~TemplateFreeWomersleyImpedanceTubeBase () |
Empty virtual destructor. More... | |
Protected Member Functions | |
void | precompute_aux_integrals () |
Precompute auxiliary integrals required for the computation of the Jacobian in the NavierStokesImpedanceTractionElement. Also pass the pointer to the pre-computed integrals to the elements in the Navier_stokes_outflow_mesh_pt so they can refer to it. More... | |
Protected Attributes | |
double | Length |
Length of the tube. More... | |
double | Dp_in_dq |
Derivative of inflow pressure w.r.t. instantaenous volume flux (Note: Can be pre-computed) More... | |
double * | Current_volume_flux_pt |
Pointer to double that specifies the currently imposed instantaneous volume flux into the impedance tube. This is used to communicate with the Womersley elements which require access to the flux via a pointer to a double. More... | |
WomersleyProblem< ELEMENT, DIM > * | Womersley_problem_pt |
Pointer to Womersley problem that determines the pressure gradient along the tube. More... | |
double | P_out |
Outlet pressure. More... | |
PrescribedVolumeFluxFctPt | Prescribed_volume_flux_fct_pt |
Pointer to function that specifies the prescribed volume flux. More... | |
Mesh * | Navier_stokes_outflow_mesh_pt |
Pointer to the mesh of NavierStokesImpedanceTractionElements that are attached to the outflow cross-section of the higher-dimensional Navier Stokes mesh and provide the inflow into the Impedance tube. More... | |
std::map< unsigned, double > * | Aux_integral_pt |
Pointer to auxiliary integral, containing the derivative of the total volume flux through the outflow boundary of the (higher-dimensional) Navier-Stokes mesh w.r.t. to the discrete (global) (velocity) degrees of freedom. More... | |
Private Attributes | |
bool | Using_flux_control_elements |
Additional Inherited Members | |
Static Public Attributes inherited from oomph::TemplateFreeWomersleyImpedanceTubeBase | |
static double | Zero = 0.0 |
Zero! More... | |
///////////////////////////////////////////////////////////////////// /////////////////////////////////////////////////////////////////////
Base class for Womersley impedance tube. Allows the computation of the inlet pressure p_in into a uniform tube of specified length that is assumed to convey fully-developed, but time-dependent flow with a presribed instantaneous flow rate, q. Also computes the derivative dp_in/dq required when this is used to determine impedance-type outlet boundary conditions in a Navier-Stokes computation.
Definition at line 1069 of file womersley_elements.h.
typedef double(* oomph::WomersleyImpedanceTubeBase< ELEMENT, DIM >::PrescribedVolumeFluxFctPt) (const double &time) |
Function pointer to fct that prescribes volume flux q=fct(t) – mainly used for validation purposes.
Definition at line 1075 of file womersley_elements.h.
|
inline |
Constructor: Specify length of tube and pointer to function that specifies the prescribed volume flux. Outlet pressure is set to zero.
Definition at line 1080 of file womersley_elements.h.
References oomph::WomersleyImpedanceTubeBase< ELEMENT, DIM >::Aux_integral_pt, and oomph::WomersleyImpedanceTubeBase< ELEMENT, DIM >::Current_volume_flux_pt.
|
inline |
Constructor: Specify length of tube and the pointer to the mesh of either NavierStokesImpedanceTractionElements or NavierStokesFluxControlElements that are attached to the outflow cross-section of a (higher-dimensional) Navier Stokes mesh and provide the inflow into the ImpedanceTube. Outlet pressure is set to zero.
Definition at line 1103 of file womersley_elements.h.
References oomph::WomersleyImpedanceTubeBase< ELEMENT, DIM >::Aux_integral_pt, oomph::WomersleyImpedanceTubeBase< ELEMENT, DIM >::Current_volume_flux_pt, e, oomph::Mesh::element_pt(), oomph::Mesh::nelement(), oomph::NavierStokesImpedanceTractionElementBase::set_external_data_from_navier_stokes_outflow_mesh(), and oomph::WomersleyImpedanceTubeBase< ELEMENT, DIM >::Using_flux_control_elements.
|
pure virtual |
Pure virtual function in which the user of a derived class must create the mesh of WomersleyElements (of the type specified by the class's template argument) and apply the boundary conditions. The Womersley elements use the timestepper specified as the input argument.
Implemented in oomph::WomersleyOutflowImpedanceTube< ELEMENT, DIM >.
Referenced by oomph::WomersleyImpedanceTubeBase< ELEMENT, DIM >::setup().
|
inlinevirtual |
Compute inlet pressure, p_in, required to achieve the currently imposed, instantaneous volume flux q prescribed by total_volume_flux_into_impedance_tube(), and its derivative, dp_in/dq.
Implements oomph::TemplateFreeWomersleyImpedanceTubeBase.
Definition at line 1375 of file womersley_elements.h.
References oomph::WomersleyImpedanceTubeBase< ELEMENT, DIM >::Current_volume_flux_pt, oomph::WomersleyImpedanceTubeBase< ELEMENT, DIM >::Dp_in_dq, oomph::WomersleyImpedanceTubeBase< ELEMENT, DIM >::Length, oomph::WomersleyImpedanceTubeBase< ELEMENT, DIM >::P_out, oomph::WomersleyImpedanceTubeBase< ELEMENT, DIM >::total_volume_flux_into_impedance_tube(), and oomph::WomersleyImpedanceTubeBase< ELEMENT, DIM >::Womersley_problem_pt.
|
inline |
Access fct to outlet pressure.
Definition at line 1170 of file womersley_elements.h.
References oomph::WomersleyImpedanceTubeBase< ELEMENT, DIM >::P_out.
|
inlineprotected |
Precompute auxiliary integrals required for the computation of the Jacobian in the NavierStokesImpedanceTractionElement. Also pass the pointer to the pre-computed integrals to the elements in the Navier_stokes_outflow_mesh_pt so they can refer to it.
Definition at line 1400 of file womersley_elements.h.
References oomph::NavierStokesImpedanceTractionElementBase::add_element_contribution_to_aux_integral(), oomph::WomersleyImpedanceTubeBase< ELEMENT, DIM >::Aux_integral_pt, e, oomph::Mesh::element_pt(), oomph::WomersleyImpedanceTubeBase< ELEMENT, DIM >::Navier_stokes_outflow_mesh_pt, oomph::Mesh::nelement(), oomph::NavierStokesImpedanceTractionElementBase::set_aux_integral_pt(), and oomph::NavierStokesImpedanceTractionElementBase::set_impedance_tube_pt().
Referenced by oomph::WomersleyImpedanceTubeBase< ELEMENT, DIM >::setup().
|
inline |
Set up the Womersley tubes so that a subsequent call to get_response(...) computes the inlet pressure for the currently prescribed instantaneous flow rate. Steady version!
Definition at line 1187 of file womersley_elements.h.
References oomph::Mesh::Default_TimeStepper, and oomph::TemplateFreeWomersleyImpedanceTubeBase::Zero.
|
inline |
Set up the Womersley tubes so that a subsequent call to get_response(...) computes the inlet pressure for the currently prescribed instantaneous flow rate, assuming that at all previous times the tube conveyed steady, fully-developed flow with flowrate q_initial. dt specifies the timestep for the subsequent time integration. Specify: Womersley number, (constant) timestep, the initial volume flux (from which the subsequent impulsive start is performed) and, optionally the pointer to the timestepper to be used in the Womersley elements (defaults to BDF<2>).
By default, we do want to suppress the output from the Newton solver
Definition at line 1207 of file womersley_elements.h.
References oomph::WomersleyImpedanceTubeBase< ELEMENT, DIM >::build_mesh_and_apply_boundary_conditions(), oomph::WomersleyImpedanceTubeBase< ELEMENT, DIM >::Current_volume_flux_pt, oomph::WomersleyImpedanceTubeBase< ELEMENT, DIM >::Dp_in_dq, oomph::WomersleyImpedanceTubeBase< ELEMENT, DIM >::Length, oomph::WomersleyImpedanceTubeBase< ELEMENT, DIM >::Navier_stokes_outflow_mesh_pt, oomph::oomph_info, oomph::WomersleyImpedanceTubeBase< ELEMENT, DIM >::precompute_aux_integrals(), oomph::WomersleyImpedanceTubeBase< ELEMENT, DIM >::Using_flux_control_elements, and oomph::WomersleyImpedanceTubeBase< ELEMENT, DIM >::Womersley_problem_pt.
|
inline |
Shift history values to allow coputation of next timestep. Note: When used with a full Navier-Stokes problem this function must be called in actions_before_implicit_timestep()
Definition at line 1311 of file womersley_elements.h.
References i, and oomph::WomersleyImpedanceTubeBase< ELEMENT, DIM >::Womersley_problem_pt.
|
inline |
Compute total current volume flux into the "impedance tube" that provides the flow resistance (flux is either obtained from the function that specifies it externally or by by adding up the flux through all NavierStokesImpedanceTractionElements in the mesh pointed to by the Navier_stokes_outflow_mesh_pt.
Definition at line 1336 of file womersley_elements.h.
References e, oomph::Mesh::element_pt(), oomph::WomersleyImpedanceTubeBase< ELEMENT, DIM >::Navier_stokes_outflow_mesh_pt, oomph::Mesh::nelement(), oomph::WomersleyImpedanceTubeBase< ELEMENT, DIM >::Prescribed_volume_flux_fct_pt, oomph::WomersleyImpedanceTubeBase< ELEMENT, DIM >::Using_flux_control_elements, and oomph::WomersleyImpedanceTubeBase< ELEMENT, DIM >::Womersley_problem_pt.
Referenced by oomph::WomersleyImpedanceTubeBase< ELEMENT, DIM >::get_response().
|
inline |
Access to underlying Womersley problem.
Definition at line 1302 of file womersley_elements.h.
References oomph::WomersleyImpedanceTubeBase< ELEMENT, DIM >::Womersley_problem_pt.
|
protected |
Pointer to auxiliary integral, containing the derivative of the total volume flux through the outflow boundary of the (higher-dimensional) Navier-Stokes mesh w.r.t. to the discrete (global) (velocity) degrees of freedom.
Definition at line 1462 of file womersley_elements.h.
Referenced by oomph::WomersleyImpedanceTubeBase< ELEMENT, DIM >::precompute_aux_integrals(), and oomph::WomersleyImpedanceTubeBase< ELEMENT, DIM >::WomersleyImpedanceTubeBase().
|
protected |
Pointer to double that specifies the currently imposed instantaneous volume flux into the impedance tube. This is used to communicate with the Womersley elements which require access to the flux via a pointer to a double.
Definition at line 1441 of file womersley_elements.h.
Referenced by oomph::WomersleyImpedanceTubeBase< ELEMENT, DIM >::get_response(), oomph::WomersleyImpedanceTubeBase< ELEMENT, DIM >::setup(), and oomph::WomersleyImpedanceTubeBase< ELEMENT, DIM >::WomersleyImpedanceTubeBase().
|
protected |
Derivative of inflow pressure w.r.t. instantaenous volume flux (Note: Can be pre-computed)
Definition at line 1435 of file womersley_elements.h.
Referenced by oomph::WomersleyImpedanceTubeBase< ELEMENT, DIM >::get_response(), and oomph::WomersleyImpedanceTubeBase< ELEMENT, DIM >::setup().
|
protected |
Length of the tube.
Definition at line 1431 of file womersley_elements.h.
Referenced by oomph::WomersleyImpedanceTubeBase< ELEMENT, DIM >::get_response(), and oomph::WomersleyImpedanceTubeBase< ELEMENT, DIM >::setup().
|
protected |
Pointer to the mesh of NavierStokesImpedanceTractionElements that are attached to the outflow cross-section of the higher-dimensional Navier Stokes mesh and provide the inflow into the Impedance tube.
Definition at line 1456 of file womersley_elements.h.
Referenced by oomph::WomersleyOutflowImpedanceTube< ELEMENT, DIM >::build_mesh_and_apply_boundary_conditions(), oomph::WomersleyImpedanceTubeBase< ELEMENT, DIM >::precompute_aux_integrals(), oomph::WomersleyImpedanceTubeBase< ELEMENT, DIM >::setup(), and oomph::WomersleyImpedanceTubeBase< ELEMENT, DIM >::total_volume_flux_into_impedance_tube().
|
protected |
Outlet pressure.
Definition at line 1448 of file womersley_elements.h.
Referenced by oomph::WomersleyImpedanceTubeBase< ELEMENT, DIM >::get_response(), and oomph::WomersleyImpedanceTubeBase< ELEMENT, DIM >::p_out().
|
protected |
Pointer to function that specifies the prescribed volume flux.
Definition at line 1451 of file womersley_elements.h.
Referenced by oomph::WomersleyImpedanceTubeBase< ELEMENT, DIM >::total_volume_flux_into_impedance_tube().
|
private |
Definition at line 1467 of file womersley_elements.h.
Referenced by oomph::WomersleyImpedanceTubeBase< ELEMENT, DIM >::setup(), oomph::WomersleyImpedanceTubeBase< ELEMENT, DIM >::total_volume_flux_into_impedance_tube(), and oomph::WomersleyImpedanceTubeBase< ELEMENT, DIM >::WomersleyImpedanceTubeBase().
|
protected |
Pointer to Womersley problem that determines the pressure gradient along the tube.
Definition at line 1445 of file womersley_elements.h.
Referenced by oomph::WomersleyImpedanceTubeBase< ELEMENT, DIM >::get_response(), oomph::WomersleyImpedanceTubeBase< ELEMENT, DIM >::setup(), oomph::WomersleyImpedanceTubeBase< ELEMENT, DIM >::shift_time_values(), oomph::WomersleyImpedanceTubeBase< ELEMENT, DIM >::total_volume_flux_into_impedance_tube(), and oomph::WomersleyImpedanceTubeBase< ELEMENT, DIM >::womersley_problem_pt().