phaseflow
FEM solver for the Navier-Stokes-Boussinesq equations coupled with enthalpy-based phase change
|
#include <phaseflow.h>
Public Member Functions | |
Phaseflow () | |
void | init (std::string file_path) |
void | run (const std::string parameter_file="") |
Public Attributes | |
Parameters::StructuredParameters | params |
Private Member Functions | |
void | create_coarse_grid () |
void | setup_system () |
void | assemble_system () |
Assemble the system for the Newton linearized Navier-Stokes-Boussinesq equtaions. More... | |
void | interpolate_boundary_values (Function< dim > *function, std::map< types::global_dof_index, double > &boundary_values) const |
void | apply_boundary_values_and_constraints () |
void | solve_linear_system () |
Solve the linear system. More... | |
void | step_newton () |
bool | solve_nonlinear_problem () |
void | set_time_step_size (double new_size) |
void | step_time () |
void | write_solution () |
void | append_verification_table () |
void | write_verification_table () |
Private Attributes | |
Triangulation< dim > | triangulation |
FESystem< dim, dim > | fe |
const FEValuesExtractors::Vector | velocity_extractor |
const FEValuesExtractors::Scalar | pressure_extractor |
const FEValuesExtractors::Scalar | temperature_extractor |
DoFHandler< dim > | dof_handler |
ConstraintMatrix | constraints |
SparsityPattern | sparsity_pattern |
SparseMatrix< double > | system_matrix |
Vector< double > | solution |
Vector< double > | newton_residual |
Vector< double > | newton_solution |
Vector< double > | old_solution |
Vector< double > | old_newton_solution |
Vector< double > | system_rhs |
double | time |
double | new_time |
double | time_step_size |
unsigned int | time_step_counter |
unsigned int | boundary_count |
std::vector< unsigned int > | manifold_ids |
std::vector< std::string > | manifold_descriptors |
Functions::ParsedFunction< dim > | source_function |
Functions::ParsedFunction< dim > | initial_values_function |
Functions::ParsedFunction< dim > | boundary_function |
Functions::ParsedFunction< dim > | exact_solution_function |
TableHandler | verification_table |
std::string | verification_table_file_name = "verification_table.txt" |
Phaseflow::Phaseflow< dim >::Phaseflow | ( | ) |
|
private |
|
private |
Apply the boundary conditions (strong and natural) and apply constraints (including those for hanging nodes
|
private |
Assemble the system for the Newton linearized Navier-Stokes-Boussinesq equtaions.
This implements equation (17) from Danaila et al. 2014.
This is the bouyancy force function from Danaila 2014,
$f_B() = theta Ra / (Pr Re^2)$
in terms of the Reynolds, Prandtl, and Rayleigh numbers, defined as:
$Re = \rho_{ref} V_{ref} L_ref/\mu_l, Pr = \nu_l / \alpha_l, Ra = g \beta L_{ref}^3 (T_h - T_c) / (\nu_l \alpha_l)$
Danaila's paper assumes gravity is pointing uniformly downward and treats gravity as a scalar. This implementation is generalized to a dim-dimensional gravity vector. This requires that we interpret the Rayleigh number also as vector-valued.
This models the bouyancy force as only a function of $$.
Following the implementation approach in
http://dealii.org/8.4.1/doxygen/deal.II/group__vector__valued.html
and
http://dealii.org/8.4.1/doxygen/deal.II/step_20.html#Assemblingthelinearsystem
Local parameters
lambda function for classical (linear) Boussinesq bouyancy
Analytical derivative of classical (linear) Boussinesq bouyancy
lambda functions for linear, bilinear, and trilinear operator
Organize data
Set local variables to match notation in Danaila 2014
Assemble element-wise
|
private |
void Phaseflow::Phaseflow< dim >::init | ( | std::string | file_path | ) |
|
private |
void Phaseflow::Phaseflow< dim >::run | ( | const std::string | parameter_file = "" | ) |
|
private |
|
private |
Setup the linear system objects.
|
private |
Solve the linear system.
|
private |
Iterate the Newton method to solve the nonlinear problem
|
private |
Setup and solve a Newton iteration
|
private |
Step the simulation from the current time step to the next time step.
This requires iterating through each Newton substep of the timestep, assembling and solving a linear system for each substep, until convergence.
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
These strings label types of manifolds used for exact geometry
|
private |
These ID's label manifolds used for exact geometry
|
private |
|
private |
|
private |
|
private |
|
private |
Parameters::StructuredParameters Phaseflow::Phaseflow< dim >::params |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |