//////////////////////////////////////////////////////////////// //////////////////////////////////////////////////////////////// More...
#include <multi_domain_boussinesq_elements.h>
Public Member Functions | |
RefineableAdvectionDiffusionBoussinesqElement () | |
Constructor: call the underlying constructors. More... | |
void | output (std::ostream &outfile, const unsigned &nplot) |
Output function: Output x, y, theta at Nplot^DIM plot points. More... | |
void | output (std::ostream &outfile) |
Overload the standard output function with the broken default. More... | |
void | output (FILE *file_pt) |
C-style output function: Broken default. More... | |
void | output (FILE *file_pt, const unsigned &n_plot) |
C-style output function: Broken default. More... | |
void | get_wind_adv_diff (const unsigned &ipt, const Vector< double > &s, const Vector< double > &x, Vector< double > &wind) const |
Overload the wind function in the advection-diffusion equations. This provides the coupling from the Navier–Stokes equations to the advection-diffusion equations because the wind is the fluid velocity, obtained from the "external" element. More... | |
void | fill_in_contribution_to_jacobian (Vector< double > &residuals, DenseMatrix< double > &jacobian) |
Compute the element's residual vector and the Jacobian matrix. More... | |
void | identify_all_field_data_for_external_interaction (Vector< std::set< FiniteElement * >> const &external_elements_pt, std::set< std::pair< Data *, unsigned >> &paired_interaction_data) |
Overload the function that must return all field data involved in the interaction with the external (Navier Stokes) element. Only the velocity dofs in the Navier Stokes element affect the interaction with the current element. More... | |
void | fill_in_contribution_to_jacobian_and_mass_matrix (Vector< double > &residuals, DenseMatrix< double > &jacobian, DenseMatrix< double > &mass_matrix) |
Add the element's contribution to its residuals vector, jacobian matrix and mass matrix. More... | |
void | get_dwind_adv_diff_dexternal_element_data (const unsigned &ipt, const unsigned &i, Vector< double > &result, Vector< unsigned > &global_eqn_number) |
Fill in the derivatives of the wind with respect to the external unknowns. More... | |
void | fill_in_off_diagonal_block_analytic (Vector< double > &residuals, DenseMatrix< double > &jacobian) |
Compute the contribution of the external degrees of freedom (velocities) on the advection-diffusion equations. More... | |
void | get_dof_numbers_for_unknowns (std::list< std::pair< unsigned long, unsigned >> &dof_lookup_list) const |
Classify dofs for use in block preconditioner. More... | |
unsigned | ndof_types () const |
Specify number of dof types for use in block preconditioner. More... | |
//////////////////////////////////////////////////////////////// ////////////////////////////////////////////////////////////////
Build an AdvectionDiffusionElement that inherits from ElementWithExternalElement so that it can "communicate" with the a NavierStokesElement that provides its wind.
Definition at line 375 of file multi_domain_boussinesq_elements.h.
|
inline |
Constructor: call the underlying constructors.
Definition at line 381 of file multi_domain_boussinesq_elements.h.
|
inline |
Compute the element's residual vector and the Jacobian matrix.
Definition at line 476 of file multi_domain_boussinesq_elements.h.
|
inline |
Add the element's contribution to its residuals vector, jacobian matrix and mass matrix.
Definition at line 507 of file multi_domain_boussinesq_elements.h.
|
inline |
Compute the contribution of the external degrees of freedom (velocities) on the advection-diffusion equations.
Definition at line 548 of file multi_domain_boussinesq_elements.h.
Referenced by oomph::RefineableAdvectionDiffusionBoussinesqElement< AD_ELEMENT, NST_ELEMENT >::fill_in_contribution_to_jacobian().
|
inline |
Classify dofs for use in block preconditioner.
Definition at line 702 of file multi_domain_boussinesq_elements.h.
|
inline |
Fill in the derivatives of the wind with respect to the external unknowns.
Definition at line 522 of file multi_domain_boussinesq_elements.h.
Referenced by oomph::RefineableAdvectionDiffusionBoussinesqElement< AD_ELEMENT, NST_ELEMENT >::fill_in_off_diagonal_block_analytic().
|
inline |
Overload the wind function in the advection-diffusion equations. This provides the coupling from the Navier–Stokes equations to the advection-diffusion equations because the wind is the fluid velocity, obtained from the "external" element.
Definition at line 456 of file multi_domain_boussinesq_elements.h.
void oomph::RefineableAdvectionDiffusionBoussinesqElement< AD_ELEMENT, NST_ELEMENT >::identify_all_field_data_for_external_interaction | ( | Vector< std::set< FiniteElement * >> const & | external_elements_pt, |
std::set< std::pair< Data *, unsigned >> & | paired_interaction_data | ||
) |
Overload the function that must return all field data involved in the interaction with the external (Navier Stokes) element. Only the velocity dofs in the Navier Stokes element affect the interaction with the current element.
//////////////////////////////////////////////////////////////////////// ////////////////////////////////////////////////////////////////////////
Overload the function that must return all field data involved in the interaction with the external (Navier Stokes) element. Only the velocity dofs in the Navier Stokes element affect the interaction with the current element.
Definition at line 763 of file multi_domain_boussinesq_elements.h.
|
inline |
Specify number of dof types for use in block preconditioner.
Definition at line 744 of file multi_domain_boussinesq_elements.h.
|
inline |
C-style output function: Broken default.
Definition at line 440 of file multi_domain_boussinesq_elements.h.
|
inline |
C-style output function: Broken default.
Definition at line 446 of file multi_domain_boussinesq_elements.h.
|
inline |
Overload the standard output function with the broken default.
Definition at line 434 of file multi_domain_boussinesq_elements.h.
|
inline |
Output function: Output x, y, theta at Nplot^DIM plot points.
Definition at line 400 of file multi_domain_boussinesq_elements.h.