ANS-1-1
|
In the next class, we define the main problem at hand. Here, we implement the top-level logic of solving a time dependent FSI problem in a monolithic ALE framework.
The initial framework of our program is based on the step-22 tutorial program, which explains best how to deal with vector-valued problems in deal.II. However, we extend that program by several additional elements: i) additional non-linearity in the fluid (convection term) -> requires non-linear solution algorithm ii) non-linear structure problem that is fully coupled to the fluid -> second source of non-linearities due to the transformation iii) implementation of a Newton-like method to solve the non-linear problem iv) Solution of a fully time-dependent problem -> implementation of various time-stepping schemes
To construct the ALE mapping for the fluid mesh motion, we solve an additional partial differential equation that is given by the biharmonic equation. This kind of equation will be split into two equations (Lit. Ciarlet), to avoid H^2 conforming finite elements.
All equations are written in a common global system that is referred to as a monolithic solution algorithm.
The discretization of the continuous problem is organized as follows:
The program is organized as follows. First, we set up runtime parameters and the system as done in other deal.II tutorial steps. Then, we assemble the system matrix (Jacobian of Newton's method) and system right hand side (residual of Newton's method) for the non-linear system. Two functions for the boundary values are provided because we are only supposed to apply boundary values in the first Newton step. In the subsequent Newton steps all Dirichlet values have to be equal zero. Afterwards, the routines for solving the linear system and the Newton iteration are self-explaining. The following function is standard in deal.II tutorial steps: writing the solutions to graphical output. The last three functions provide the framework to compute functional values of interest. For the given fluid-structure interaction problem, we compute the displacement in the x- and y-directions of the structure at a certain point. We are also interested in the observation of the drag- and lift evaluations, which are achieved by line-integration over faces.
Definition at line 914 of file step-fsi.cc.
FSIALEProblem< dim >::FSIALEProblem | ( | const unsigned int | degree | ) |
The constructor of this class is comparable to other tutorials steps, e.g., step-22, and step-31. We are going to use the following finite element discretization: Q_2^c for the fluid, Q_2^c for the structure, P_1^dc for the pressure, and Q_2^c for the additional displacement.
Definition at line 990 of file step-fsi.cc.
FSIALEProblem< dim >::~FSIALEProblem | ( | ) |
Definition at line 1005 of file step-fsi.cc.
void FSIALEProblem< dim >::run | ( | ) |
As usual in deal.II programs, we have to call the run method. It handles the output stream to the terminal. Second, we define some output skip that is necessary (and really useful) to avoid to much printing of solutions. For large time dependent problems it is sufficient to print only each tenth solution. Third, we perform the time stepping scheme of the solution process.
Definition at line 2969 of file step-fsi.cc.
Referenced by main().
void FSIALEProblem< dim >::set_runtime_parameters | ( | ) | [private] |
In this method, we set up runtime parameters that could also come from a paramter file. We propose three different configurations FSI 1, FSI 2, and FSI 3. The reader is invited to change these values to obtain other results.
Definition at line 1016 of file step-fsi.cc.
References Assert, GridIn< dim, spacedim >::attach_triangulation(), ExcInternalError(), and GridIn< dim, spacedim >::read_ucd().
void FSIALEProblem< dim >::setup_system | ( | ) | [private] |
Definition at line 1103 of file step-fsi.cc.
References BlockCompressedSimpleSparsityPattern::block(), BlockCompressedSimpleSparsityPattern::collect_sizes(), DoFRenumbering::component_wise(), DoFTools::count_dofs_per_block(), DoFRenumbering::Cuthill_McKee(), DoFTools::make_hanging_node_constraints(), DoFTools::make_sparsity_pattern(), and CompressedSimpleSparsityPattern::reinit().
void FSIALEProblem< dim >::assemble_system_matrix | ( | ) | [private] |
In this function, we assemble the Jacobian matrix for the Newton iteration. The fluid and the structure equations are computed on different sub-domains in the mesh and ask for the corresponding material ids. The fluid equations are defined on mesh cells with the material id == 0 and the structure equations on cells with the material id == 1.
To compensate the well-known problem in fluid dynamics on the outflow boundary, we also add some correction term on the outflow boundary. This relation is known as `do-nothing' condition. In the inner loops of the local_cell_matrix, the time dependent equations are discretized with a finite difference scheme. Quasi-stationary processes (FSI 1) can be computed by the BE scheme. The other two schemes are useful for non-stationary computations (FSI 2 and FSI 3).
Assembling of the inner most loop is treated with help of the fe.system_to_component_index(j).first function from the library. Using this function makes the assembling process much faster than running over all local degrees of freedom.
Definition at line 1259 of file step-fsi.cc.
References DoFHandler< dim, spacedim >::begin_active(), NSEALE::get_accelaration_term_LinAll(), ALETransformations::get_F_Inverse_LinU(), FEFaceValues< dim, spacedim >::get_function_grads(), FEValues< dim, spacedim >::get_function_grads(), FEFaceValues< dim, spacedim >::get_function_values(), FEValues< dim, spacedim >::get_function_values(), NSEALE::get_stress_fluid_ALE_2nd_term_LinAll_short(), FEFaceValues< dim, spacedim >::JxW(), FEValues< dim, spacedim >::JxW(), FEFaceValues< dim, spacedim >::normal_vector(), FEValues< dim, spacedim >::reinit(), FEFaceValues< dim, spacedim >::reinit(), QGauss< dim >::size(), update_gradients, update_JxW_values, update_normal_vectors, update_quadrature_points, and update_values.
void FSIALEProblem< dim >::assemble_system_rhs | ( | ) | [private] |
In this function we assemble the semi-linear form of the right hand side of Newton's method (its residual). The framework is in principal the same as for the system matrix.
Definition at line 1838 of file step-fsi.cc.
References DoFHandler< dim, spacedim >::begin_active(), Tensor< rank_, dim, Number >::clear(), FEValues< dim, spacedim >::get_function_grads(), FEFaceValues< dim, spacedim >::get_function_grads(), FEValues< dim, spacedim >::get_function_values(), FEFaceValues< dim, spacedim >::get_function_values(), FEFaceValues< dim, spacedim >::JxW(), FEValues< dim, spacedim >::JxW(), FEFaceValues< dim, spacedim >::normal_vector(), FEValues< dim, spacedim >::reinit(), FEFaceValues< dim, spacedim >::reinit(), QGauss< dim >::size(), update_gradients, update_JxW_values, update_normal_vectors, update_quadrature_points, update_values, and value.
void FSIALEProblem< dim >::set_initial_bc | ( | const double | time | ) | [private] |
In this function, we impose boundary conditions for the whole system. The fluid inflow is prescribed by a parabolic profile. The usual structure displacement shall be fixed at all outer boundaries. Consequently our formulation of the mixed biharmonic equation requires no Dirichlet zero values for the second displacement variable $w$ (see the standard literature to elasticity and Ciarlet). The pressure variable is not subjected to any Dirichlet boundary conditions and is left free in this method. Please note, that the interface between fluid and structure has no physical boundary due to our formulation. Interface conditions are automatically fulfilled: that is one major advantage of the `monolithic' formulation.
Definition at line 2457 of file step-fsi.cc.
References VectorTools::interpolate_boundary_values().
void FSIALEProblem< dim >::set_newton_bc | ( | ) | [private] |
This function applies boundary conditions to the Newton iteration steps. For all variables that have Dirichlet conditions on some (or all) parts of the outer boundary, we now apply zero-Dirichlet conditions.
Definition at line 2517 of file step-fsi.cc.
References VectorTools::interpolate_boundary_values().
void FSIALEProblem< dim >::solve | ( | ) | [private] |
In this function, we solve the linear systems inside the nonlinear Newton iteration. We only use a direct solver from UMFPACK. The reason is twofold: First, the focus of this implementation is more on time-dependent problems. Hence, a huge amount of spatial degrees of freedom is not our primal goal. For this, a direct solver is an adequate tool. Second, the devolpement of an iterative solver based on the GMRES scheme for instance, is difficult to derive. Only a few results are known in the literature so far. However, Baerbel Janssen and the author have already had a first try for our implementation at hand (that is working, of coarse) but suffers from a good preconditioner for the GMRES scheme (Lit. B. Janssen, T. Wick, ECCOMAS 2010). This lack will be resolved in upcoming work but we invite everybody to collaborate with us if he/she has a resonable idea.
Definition at line 2678 of file step-fsi.cc.
References SparseDirectUMFPACK::factorize(), and SparseDirectUMFPACK::vmult().
void FSIALEProblem< dim >::newton_iteration | ( | const double | time | ) | [private] |
This function performs the Newton iteration to solve the non-linear system of equations. First, we declare some standard parameters of the solution method. Addionally, we also implement an easy line search algorithm.
Definition at line 2563 of file step-fsi.cc.
References Timer::reset(), Timer::start(), and Timer::stop().
void FSIALEProblem< dim >::output_results | ( | const unsigned int | refinement_cycle, |
const BlockVector< double > | solution | ||
) | const [private] |
Definition at line 2714 of file step-fsi.cc.
References DataOut< dim, DH >::add_data_vector(), DataOut< dim, DH >::attach_dof_handler(), DataOut< dim, DH >::build_patches(), DataComponentInterpretation::component_is_part_of_vector, DataComponentInterpretation::component_is_scalar, Utilities::int_to_string(), and DataOut< dim, DH >::write_vtk().
double FSIALEProblem< dim >::compute_point_value | ( | Point< dim > | p, |
const unsigned int | component | ||
) | const [private] |
With help of this function, we extract point values for a certain component from our discrete solution. We use it to gain the displacements of the structure in the x- and y-directions.
Definition at line 2770 of file step-fsi.cc.
References VectorTools::point_value().
void FSIALEProblem< dim >::compute_drag_lift_fsi_fluid_tensor | ( | ) | [private] |
Now, we arrive at the function that is responsible to compute the line integrals for the drag and the lift. Note, that by a proper transformation via the Gauss theorem, the both quantities could also be achieved by domain integral computation. Nevertheless, we choose the line integration because deal.II provides all routines for face value evaluation.
Definition at line 2791 of file step-fsi.cc.
References DoFHandler< dim, spacedim >::begin_active(), Tensor< rank_, dim, Number >::clear(), FEFaceValues< dim, spacedim >::get_function_grads(), FEFaceValues< dim, spacedim >::get_function_values(), FEFaceValues< dim, spacedim >::JxW(), FEFaceValues< dim, spacedim >::normal_vector(), FEFaceValues< dim, spacedim >::reinit(), update_gradients, update_JxW_values, update_normal_vectors, and update_values.
void FSIALEProblem< dim >::compute_functional_values | ( | ) | [private] |
Here, we compute the four quantities of interest: the x and y-displacements of the structure, the drag, and the lift.
Definition at line 2943 of file step-fsi.cc.
const unsigned int FSIALEProblem< dim >::degree [private] |
Definition at line 943 of file step-fsi.cc.
Triangulation<dim> FSIALEProblem< dim >::triangulation [private] |
Definition at line 945 of file step-fsi.cc.
FESystem<dim> FSIALEProblem< dim >::fe [private] |
Definition at line 946 of file step-fsi.cc.
DoFHandler<dim> FSIALEProblem< dim >::dof_handler [private] |
Definition at line 947 of file step-fsi.cc.
ConstraintMatrix FSIALEProblem< dim >::constraints [private] |
Definition at line 949 of file step-fsi.cc.
BlockSparsityPattern FSIALEProblem< dim >::sparsity_pattern [private] |
Definition at line 951 of file step-fsi.cc.
BlockSparseMatrix<double> FSIALEProblem< dim >::system_matrix [private] |
Definition at line 952 of file step-fsi.cc.
BlockVector<double> FSIALEProblem< dim >::solution [private] |
Definition at line 954 of file step-fsi.cc.
BlockVector<double> FSIALEProblem< dim >::newton_update [private] |
Definition at line 954 of file step-fsi.cc.
BlockVector<double> FSIALEProblem< dim >::old_timestep_solution [private] |
Definition at line 954 of file step-fsi.cc.
BlockVector<double> FSIALEProblem< dim >::system_rhs [private] |
Definition at line 955 of file step-fsi.cc.
TimerOutput FSIALEProblem< dim >::timer [private] |
Definition at line 957 of file step-fsi.cc.
unsigned int FSIALEProblem< dim >::timestep_number [private] |
Definition at line 960 of file step-fsi.cc.
unsigned int FSIALEProblem< dim >::max_no_timesteps [private] |
Definition at line 961 of file step-fsi.cc.
double FSIALEProblem< dim >::timestep [private] |
Definition at line 962 of file step-fsi.cc.
double FSIALEProblem< dim >::theta [private] |
Definition at line 962 of file step-fsi.cc.
double FSIALEProblem< dim >::time [private] |
Definition at line 962 of file step-fsi.cc.
std::string FSIALEProblem< dim >::time_stepping_scheme [private] |
Definition at line 963 of file step-fsi.cc.
double FSIALEProblem< dim >::density_fluid [private] |
Definition at line 966 of file step-fsi.cc.
double FSIALEProblem< dim >::viscosity [private] |
Definition at line 966 of file step-fsi.cc.
double FSIALEProblem< dim >::density_structure [private] |
Definition at line 969 of file step-fsi.cc.
double FSIALEProblem< dim >::lame_coefficient_mu [private] |
Definition at line 970 of file step-fsi.cc.
double FSIALEProblem< dim >::lame_coefficient_lambda [private] |
Definition at line 970 of file step-fsi.cc.
double FSIALEProblem< dim >::poisson_ratio_nu [private] |
Definition at line 970 of file step-fsi.cc.
double FSIALEProblem< dim >::cell_diameter [private] |
Definition at line 973 of file step-fsi.cc.
double FSIALEProblem< dim >::alpha_u [private] |
Definition at line 974 of file step-fsi.cc.
double FSIALEProblem< dim >::alpha_w [private] |
Definition at line 974 of file step-fsi.cc.