OpenCMISS-Iron Internal API Documentation
|
This module handles all Stokes fluid routines. More...
Functions/Subroutines | |
subroutine, public | stokes_equationssetsolutionmethodset (EQUATIONS_SET, SOLUTION_METHOD, ERR, ERROR,) |
Sets/changes the solution method for a Stokes flow equation type of an fluid mechanics equations set class. More... | |
subroutine, public | stokes_equationssetspecificationset (equationsSet, specification, err, error,) |
Sets the equation specification for a Stokes flow equation of a fluid mechanics equations set. More... | |
subroutine, public | stokes_equations_set_setup (EQUATIONS_SET, EQUATIONS_SET_SETUP, ERR, ERROR,) |
Sets up the standard Stokes fluid setup. More... | |
subroutine, public | stokes_problemspecificationset (problem, problemSpecification, err, error,) |
Sets the problem specification for a Stokes fluid problem. More... | |
subroutine, public | stokes_problem_setup (PROBLEM, PROBLEM_SETUP, ERR, ERROR,) |
Sets up the Stokes problem. More... | |
subroutine, public | stokes_finite_element_calculate (EQUATIONS_SET, ELEMENT_NUMBER, ERR, ERROR,) |
Calculates the element stiffness matrices and RHS for a Stokes fluid finite element equations set. More... | |
subroutine, public | stokes_post_solve (CONTROL_LOOP, SOLVER, ERR, ERROR,) |
Sets up the Stokes problem post solve. More... | |
subroutine, public | stokes_pre_solve (CONTROL_LOOP, SOLVER, ERR, ERROR,) |
Sets up the Stokes problem pre solve. More... | |
subroutine | stokes_pre_solve_update_boundary_conditions (CONTROL_LOOP, SOLVER, ERR, ERROR,) |
Update boundary conditions for Stokes flow pre solve. More... | |
subroutine | stokes_pre_solve_ale_update_mesh (CONTROL_LOOP, SOLVER, ERR, ERROR,) |
Update mesh velocity and move mesh for ALE Stokes problem. More... | |
subroutine | stokes_pre_solve_ale_update_parameters (CONTROL_LOOP, SOLVER, ERR, ERROR,) |
Update mesh parameters for three component Laplace problem. More... | |
subroutine | stokes_post_solve_output_data (CONTROL_LOOP, SOLVER, err, error,) |
Output data post solve. More... | |
subroutine, public | stokes_boundaryconditionsanalyticcalculate (EQUATIONS_SET, BOUNDARY_CONDITIONS, ERR, ERROR,) |
Calculates the analytic solution and sets the boundary conditions for an analytic problem. More... | |
subroutine, public | stokes_equation_analytic_functions (VALUE, X, MU_PARAM, RHO_PARAM, CURRENT_TIME, VARIABLE_TYPE, GLOBAL_DERIV_INDEX, ANALYTIC_FUNCTION_TYPE, NUMBER_OF_DIMENSIONS, NUMBER_OF_COMPONENTS, COMPONENT_IDX, ERR, ERROR,) |
Calculates the various analytic solutions given X and time, can be called from within analytic calculate or elsewhere if needed. More... | |
This module handles all Stokes fluid routines.
subroutine, public stokes_equations_routines::stokes_boundaryconditionsanalyticcalculate | ( | type(equations_set_type), pointer | EQUATIONS_SET, |
type(boundary_conditions_type), pointer | BOUNDARY_CONDITIONS, | ||
integer(intg), intent(out) | ERR, | ||
type(varying_string), intent(out) | ERROR | ||
) |
Calculates the analytic solution and sets the boundary conditions for an analytic problem.
[out] | err | The error code |
[out] | error | The error string |
Definition at line 3211 of file Stokes_equations_routines.f90.
References boundary_conditions_routines::boundary_condition_fixed, base_routines::enters(), equations_set_constants::equations_set_stokes_equation_three_dim_1, equations_set_constants::equations_set_stokes_equation_three_dim_2, equations_set_constants::equations_set_stokes_equation_three_dim_3, equations_set_constants::equations_set_stokes_equation_three_dim_4, equations_set_constants::equations_set_stokes_equation_three_dim_5, equations_set_constants::equations_set_stokes_equation_two_dim_1, equations_set_constants::equations_set_stokes_equation_two_dim_2, equations_set_constants::equations_set_stokes_equation_two_dim_3, equations_set_constants::equations_set_stokes_equation_two_dim_4, equations_set_constants::equations_set_stokes_equation_two_dim_5, base_routines::exits(), constants::no_part_deriv, and stokes_equation_analytic_functions().
Referenced by fluid_mechanics_routines::fluidmechanics_boundaryconditionsanalyticcalculate().
subroutine, public stokes_equations_routines::stokes_equation_analytic_functions | ( | real(dp), intent(out) | VALUE, |
real(dp), dimension(3), intent(in) | X, | ||
real(dp) | MU_PARAM, | ||
real(dp) | RHO_PARAM, | ||
real(dp), intent(in) | CURRENT_TIME, | ||
integer(intg) | VARIABLE_TYPE, | ||
integer(intg) | GLOBAL_DERIV_INDEX, | ||
integer(intg) | ANALYTIC_FUNCTION_TYPE, | ||
integer(intg), intent(in) | NUMBER_OF_DIMENSIONS, | ||
integer(intg), intent(in) | NUMBER_OF_COMPONENTS, | ||
integer(intg), intent(in) | COMPONENT_IDX, | ||
integer(intg), intent(out) | ERR, | ||
type(varying_string), intent(out) | ERROR | ||
) |
Calculates the various analytic solutions given X and time, can be called from within analytic calculate or elsewhere if needed.
[out] | err | The error code |
[out] | error | The error string |
Definition at line 3673 of file Stokes_equations_routines.f90.
References base_routines::enters(), equations_set_constants::equations_set_stokes_equation_three_dim_1, equations_set_constants::equations_set_stokes_equation_three_dim_2, equations_set_constants::equations_set_stokes_equation_three_dim_3, equations_set_constants::equations_set_stokes_equation_three_dim_4, equations_set_constants::equations_set_stokes_equation_three_dim_5, equations_set_constants::equations_set_stokes_equation_two_dim_1, equations_set_constants::equations_set_stokes_equation_two_dim_2, equations_set_constants::equations_set_stokes_equation_two_dim_3, equations_set_constants::equations_set_stokes_equation_two_dim_4, equations_set_constants::equations_set_stokes_equation_two_dim_5, base_routines::exits(), constants::global_deriv_s1, constants::global_deriv_s1_s2, constants::global_deriv_s2, constants::no_global_deriv, and constants::pi.
Referenced by stokes_boundaryconditionsanalyticcalculate(), and stokes_pre_solve_update_boundary_conditions().
subroutine, public stokes_equations_routines::stokes_equations_set_setup | ( | type(equations_set_type), pointer | EQUATIONS_SET, |
type(equations_set_setup_type), intent(inout) | EQUATIONS_SET_SETUP, | ||
integer(intg), intent(out) | ERR, | ||
type(varying_string), intent(out) | ERROR | ||
) |
Sets up the standard Stokes fluid setup.
equations_set | A pointer to the equations set to setup | |
[in,out] | equations_set_setup | The equations set setup information |
[out] | err | The error code |
[out] | error | The error string |
Definition at line 230 of file Stokes_equations_routines.f90.
References distributed_matrix_vector::distributed_matrix_block_storage_type, distributed_matrix_vector::distributed_matrix_compressed_row_storage_type, distributed_matrix_vector::distributed_matrix_diagonal_storage_type, base_routines::enters(), equations_routines::equations_create_finish(), equations_routines::equations_create_start(), equations_set_constants::equations_first_order_dynamic, equations_set_constants::equations_linear, equations_routines::equations_linearity_type_set(), equations_routines::equations_lumped_matrices, equations_mapping_routines::equations_mapping_create_finish(), equations_mapping_routines::equations_mapping_create_start(), equations_mapping_routines::equations_mapping_dynamic_variable_type_set(), equations_mapping_routines::equations_mapping_rhs_variable_type_set(), equations_matrices_routines::equations_matrices_create_finish(), equations_matrices_routines::equations_matrices_create_start(), equations_matrices_routines::equations_matrices_dynamic_lumping_type_set(), equations_matrices_routines::equations_matrices_dynamic_storage_type_set(), equations_matrices_routines::equations_matrices_full_matrices, equations_matrices_routines::equations_matrices_linear_storage_type_set(), equations_matrices_routines::equations_matrices_sparse_matrices, equations_matrices_routines::equations_matrix_diagonal_structure, equations_matrices_routines::equations_matrix_fem_structure, equations_matrices_routines::equations_matrix_lumped, equations_matrices_routines::equations_matrix_unlumped, equations_set_constants::equations_set_ale_stokes_subtype, equations_routines::equations_set_equations_get(), equations_set_constants::equations_set_fem_solution_method, equations_set_constants::equations_set_laplace_stokes_subtype, equations_set_constants::equations_set_pgm_stokes_subtype, equations_set_constants::equations_set_setup_analytic_type, equations_set_constants::equations_set_setup_dependent_type, equations_set_constants::equations_set_setup_equations_type, equations_set_constants::equations_set_setup_finish_action, equations_set_constants::equations_set_setup_geometry_type, equations_set_constants::equations_set_setup_independent_type, equations_set_constants::equations_set_setup_initial_type, equations_set_constants::equations_set_setup_materials_type, equations_set_constants::equations_set_setup_source_type, equations_set_constants::equations_set_setup_start_action, equations_set_constants::equations_set_static_stokes_subtype, equations_set_constants::equations_set_stokes_equation_three_dim_1, equations_set_constants::equations_set_stokes_equation_three_dim_2, equations_set_constants::equations_set_stokes_equation_three_dim_3, equations_set_constants::equations_set_stokes_equation_three_dim_4, equations_set_constants::equations_set_stokes_equation_three_dim_5, equations_set_constants::equations_set_stokes_equation_two_dim_1, equations_set_constants::equations_set_stokes_equation_two_dim_2, equations_set_constants::equations_set_stokes_equation_two_dim_3, equations_set_constants::equations_set_stokes_equation_two_dim_4, equations_set_constants::equations_set_stokes_equation_two_dim_5, equations_set_constants::equations_set_transient_stokes_subtype, equations_set_constants::equations_static, equations_routines::equations_time_dependence_type_set(), equations_mapping_routines::equationsmapping_linearmatricesnumberset(), equations_mapping_routines::equationsmapping_linearmatricesvariabletypesset(), equations_matrices_routines::equationsmatrices_dynamicstructuretypeset(), equations_matrices_routines::equationsmatrices_linearstructuretypeset(), base_routines::exits(), matrix_vector::matrix_block_storage_type, matrix_vector::matrix_compressed_row_storage_type, and stokes_equationssetsolutionmethodset().
Referenced by fluid_mechanics_routines::fluid_mechanics_equations_set_setup().
subroutine, public stokes_equations_routines::stokes_equationssetsolutionmethodset | ( | type(equations_set_type), pointer | EQUATIONS_SET, |
integer(intg), intent(in) | SOLUTION_METHOD, | ||
integer(intg), intent(out) | ERR, | ||
type(varying_string), intent(out) | ERROR | ||
) |
Sets/changes the solution method for a Stokes flow equation type of an fluid mechanics equations set class.
equations_set | A pointer to the equations set to set the solution method for | |
[in] | solution_method | The solution method to set |
[out] | err | The error code |
[out] | error | The error string |
Definition at line 108 of file Stokes_equations_routines.f90.
References base_routines::enters(), equations_set_constants::equations_set_ale_stokes_subtype, equations_set_constants::equations_set_bem_solution_method, equations_set_constants::equations_set_fd_solution_method, equations_set_constants::equations_set_fem_solution_method, equations_set_constants::equations_set_fv_solution_method, equations_set_constants::equations_set_gfem_solution_method, equations_set_constants::equations_set_gfv_solution_method, equations_set_constants::equations_set_laplace_stokes_subtype, equations_set_constants::equations_set_pgm_stokes_subtype, equations_set_constants::equations_set_static_stokes_subtype, equations_set_constants::equations_set_transient_stokes_subtype, and base_routines::exits().
Referenced by fluid_mechanics_routines::fluidmechanics_equationssetsolutionmethodset(), and stokes_equations_set_setup().
subroutine, public stokes_equations_routines::stokes_equationssetspecificationset | ( | type(equations_set_type), pointer | equationsSet, |
integer(intg), dimension(:), intent(in) | specification, | ||
integer(intg), intent(out) | err, | ||
type(varying_string), intent(out) | error | ||
) |
Sets the equation specification for a Stokes flow equation of a fluid mechanics equations set.
equationsset | A pointer to the equations set to set the specification for | |
[in] | specification | The equations set specification to set |
[out] | err | The error code |
[out] | error | The error string |
Definition at line 172 of file Stokes_equations_routines.f90.
References base_routines::enters(), equations_set_constants::equations_set_ale_stokes_subtype, equations_set_constants::equations_set_fluid_mechanics_class, equations_set_constants::equations_set_laplace_stokes_subtype, equations_set_constants::equations_set_optimised_stokes_subtype, equations_set_constants::equations_set_pgm_stokes_subtype, equations_set_constants::equations_set_static_stokes_subtype, equations_set_constants::equations_set_stokes_equation_type, equations_set_constants::equations_set_transient_stokes_subtype, base_routines::errors(), and base_routines::exits().
Referenced by fluid_mechanics_routines::fluidmechanics_equationssetspecificationset().
subroutine, public stokes_equations_routines::stokes_finite_element_calculate | ( | type(equations_set_type), pointer | EQUATIONS_SET, |
integer(intg), intent(in) | ELEMENT_NUMBER, | ||
integer(intg), intent(out) | ERR, | ||
type(varying_string), intent(out) | ERROR | ||
) |
Calculates the element stiffness matrices and RHS for a Stokes fluid finite element equations set.
equations_set | A pointer to the equations set to perform the finite element calculations on | |
[in] | element_number | The element number to calculate |
[out] | err | The error code |
[out] | error | The error string |
Definition at line 1328 of file Stokes_equations_routines.f90.
References basis_routines::basis_default_quadrature_scheme, base_routines::enters(), equations_set_constants::equations_set_ale_stokes_subtype, equations_set_constants::equations_set_laplace_stokes_subtype, equations_set_constants::equations_set_pgm_stokes_subtype, equations_set_constants::equations_set_static_stokes_subtype, equations_set_constants::equations_set_stokes_equation_three_dim_1, equations_set_constants::equations_set_stokes_equation_three_dim_2, equations_set_constants::equations_set_stokes_equation_three_dim_3, equations_set_constants::equations_set_stokes_equation_three_dim_4, equations_set_constants::equations_set_stokes_equation_three_dim_5, equations_set_constants::equations_set_stokes_equation_two_dim_1, equations_set_constants::equations_set_stokes_equation_two_dim_2, equations_set_constants::equations_set_stokes_equation_two_dim_3, equations_set_constants::equations_set_stokes_equation_two_dim_4, equations_set_constants::equations_set_stokes_equation_two_dim_5, equations_set_constants::equations_set_transient_stokes_subtype, base_routines::exits(), constants::first_part_deriv, constants::no_part_deriv, constants::partial_derivative_first_derivative_map, constants::pi, problem_constants::problem_ale_stokes_subtype, problem_constants::problem_pgm_stokes_subtype, and kinds::ptr.
Referenced by fluid_mechanics_routines::fluid_mechanics_finite_element_calculate().
subroutine, public stokes_equations_routines::stokes_post_solve | ( | type(control_loop_type), pointer | CONTROL_LOOP, |
type(solver_type), pointer | SOLVER, | ||
integer(intg), intent(out) | ERR, | ||
type(varying_string), intent(out) | ERROR | ||
) |
Sets up the Stokes problem post solve.
control_loop | A pointer to the control loop to solve. | |
solver | A pointer to the solver | |
[out] | err | The error code |
[out] | error | The error string |
Definition at line 1934 of file Stokes_equations_routines.f90.
References base_routines::enters(), base_routines::exits(), base_routines::general_output_type, problem_constants::problem_ale_stokes_subtype, problem_constants::problem_laplace_stokes_subtype, problem_constants::problem_pgm_stokes_subtype, problem_constants::problem_static_stokes_subtype, problem_constants::problem_transient_stokes_subtype, solver_routines::solver_dynamic_type, solver_routines::solver_linear_type, solver_routines::solvers_solver_get(), and stokes_post_solve_output_data().
Referenced by fluid_mechanics_routines::fluid_mechanics_post_solve().
|
private |
Output data post solve.
control_loop | A pointer to the control loop to solve. | |
solver | A pointer to the solver | |
[out] | err | The error code |
[out] | error | The error string |
Definition at line 3078 of file Stokes_equations_routines.f90.
References analytic_analysis_routines::analyticanalysis_output(), control_loop_routines::control_loop_current_times_get(), control_loop_routines::control_loop_progress_output, base_routines::enters(), equations_set_constants::equations_set_navier_stokes_equation_three_dim_1, equations_set_constants::equations_set_navier_stokes_equation_three_dim_4, equations_set_constants::equations_set_navier_stokes_equation_three_dim_5, equations_set_constants::equations_set_navier_stokes_equation_two_dim_4, equations_set_constants::equations_set_navier_stokes_equation_two_dim_5, base_routines::exits(), field_io_routines::field_io_elements_export(), field_io_routines::field_io_nodes_export(), base_routines::general_output_type, problem_constants::problem_ale_stokes_subtype, problem_constants::problem_laplace_stokes_subtype, problem_constants::problem_pgm_stokes_subtype, problem_constants::problem_static_stokes_subtype, problem_constants::problem_transient_stokes_subtype, and solver_routines::solver_progress_output.
Referenced by stokes_post_solve().
subroutine, public stokes_equations_routines::stokes_pre_solve | ( | type(control_loop_type), pointer | CONTROL_LOOP, |
type(solver_type), pointer | SOLVER, | ||
integer(intg), intent(out) | ERR, | ||
type(varying_string), intent(out) | ERROR | ||
) |
Sets up the Stokes problem pre solve.
control_loop | A pointer to the control loop to solve. | |
solver | A pointer to the solver | |
[out] | err | The error code |
[out] | error | The error string |
Definition at line 2004 of file Stokes_equations_routines.f90.
References base_routines::enters(), base_routines::exits(), base_routines::general_output_type, problem_constants::problem_ale_stokes_subtype, problem_constants::problem_laplace_stokes_subtype, problem_constants::problem_pgm_stokes_subtype, problem_constants::problem_static_stokes_subtype, problem_constants::problem_transient_stokes_subtype, solver_routines::solver_dynamic_type, solver_routines::solver_linear_type, solver_routines::solvers_solver_get(), stokes_pre_solve_ale_update_mesh(), stokes_pre_solve_ale_update_parameters(), and stokes_pre_solve_update_boundary_conditions().
Referenced by fluid_mechanics_routines::fluid_mechanics_pre_solve().
|
private |
Update mesh velocity and move mesh for ALE Stokes problem.
control_loop | A pointer to the control loop to solve. | |
solver | A pointer to the solvers | |
[out] | err | The error code |
[out] | error | The error string |
Definition at line 2677 of file Stokes_equations_routines.f90.
References control_loop_routines::control_loop_current_times_get(), base_routines::enters(), base_routines::exits(), fluid_mechanics_io_routines::fluid_mechanics_io_read_data(), problem_constants::problem_ale_stokes_subtype, problem_constants::problem_laplace_stokes_subtype, problem_constants::problem_pgm_stokes_subtype, problem_constants::problem_static_stokes_subtype, problem_constants::problem_transient_stokes_subtype, solver_routines::solver_dynamic_type, solver_routines::solver_linear_type, and solver_routines::solvers_solver_get().
Referenced by stokes_pre_solve().
|
private |
Update mesh parameters for three component Laplace problem.
control_loop | A pointer to the control loop to solve. | |
solver | A pointer to the solver | |
[out] | err | The error code |
[out] | error | The error string |
Definition at line 2949 of file Stokes_equations_routines.f90.
References control_loop_routines::control_loop_current_times_get(), base_routines::enters(), base_routines::exits(), problem_constants::problem_ale_stokes_subtype, problem_constants::problem_laplace_stokes_subtype, problem_constants::problem_static_stokes_subtype, problem_constants::problem_transient_stokes_subtype, solver_routines::solver_dynamic_type, and solver_routines::solver_linear_type.
Referenced by stokes_pre_solve().
|
private |
Update boundary conditions for Stokes flow pre solve.
control_loop | A pointer to the control loop to solve. | |
solver | A pointer to the solver | |
[out] | err | The error code |
[out] | error | The error string |
Definition at line 2092 of file Stokes_equations_routines.f90.
References boundary_conditions_routines::boundary_condition_fixed, boundary_conditions_routines::boundary_condition_fixed_inlet, boundary_conditions_routines::boundary_condition_moved_wall, boundary_conditions_routines::boundary_conditions_variable_get(), control_loop_routines::control_loop_current_times_get(), base_routines::enters(), equations_set_constants::equations_set_stokes_equation_three_dim_1, equations_set_constants::equations_set_stokes_equation_three_dim_4, equations_set_constants::equations_set_stokes_equation_three_dim_5, equations_set_constants::equations_set_stokes_equation_two_dim_4, equations_set_constants::equations_set_stokes_equation_two_dim_5, base_routines::exits(), fluid_mechanics_io_routines::fluid_mechanics_io_read_boundary_conditions(), base_routines::general_output_type, constants::no_part_deriv, problem_constants::problem_ale_stokes_subtype, problem_constants::problem_laplace_stokes_subtype, problem_constants::problem_pgm_stokes_subtype, problem_constants::problem_static_stokes_subtype, problem_constants::problem_transient_stokes_subtype, solver_routines::solver_dynamic_type, solver_routines::solver_linear_type, and stokes_equation_analytic_functions().
Referenced by stokes_pre_solve().
subroutine, public stokes_equations_routines::stokes_problem_setup | ( | type(problem_type), pointer | PROBLEM, |
type(problem_setup_type), intent(inout) | PROBLEM_SETUP, | ||
integer(intg), intent(out) | ERR, | ||
type(varying_string), intent(out) | ERROR | ||
) |
Sets up the Stokes problem.
problem | A pointer to the problem set to setup a Stokes fluid on. | |
[in,out] | problem_setup | The problem setup information |
[out] | err | The error code |
[out] | error | The error string |
Definition at line 981 of file Stokes_equations_routines.f90.
References control_loop_routines::control_loop_create_finish(), control_loop_routines::control_loop_create_start(), control_loop_routines::control_loop_node, control_loop_routines::control_loop_solvers_get(), control_loop_routines::control_loop_type_set(), base_routines::enters(), base_routines::exits(), problem_constants::problem_ale_stokes_subtype, problem_constants::problem_control_time_loop_type, problem_constants::problem_laplace_stokes_subtype, problem_constants::problem_pgm_stokes_subtype, problem_constants::problem_setup_control_type, problem_constants::problem_setup_finish_action, problem_constants::problem_setup_initial_type, problem_constants::problem_setup_solver_equations_type, problem_constants::problem_setup_solvers_type, problem_constants::problem_setup_start_action, problem_constants::problem_static_stokes_subtype, problem_constants::problem_transient_stokes_subtype, solver_routines::solver_cmiss_library, solver_routines::solver_dynamic_crank_nicolson_scheme, solver_routines::solver_dynamic_degree_set(), solver_routines::solver_dynamic_first_degree, solver_routines::solver_dynamic_first_order, solver_routines::solver_dynamic_order_set(), solver_routines::solver_dynamic_scheme_set(), solver_routines::solver_dynamic_type, solver_routines::solver_equations_create_finish(), solver_routines::solver_equations_create_start(), problem_constants::solver_equations_first_order_dynamic, problem_constants::solver_equations_linear, solver_routines::solver_equations_linearity_type_set(), solver_routines::solver_equations_sparsity_type_set(), problem_constants::solver_equations_static, solver_routines::solver_equations_time_dependence_type_set(), solver_routines::solver_library_type_set(), solver_routines::solver_linear_type, solver_routines::solver_petsc_library, solver_routines::solver_solver_equations_get(), solver_routines::solver_sparse_matrices, solver_routines::solver_type_set(), solver_routines::solvers_create_finish(), solver_routines::solvers_create_start(), solver_routines::solvers_number_set(), and solver_routines::solvers_solver_get().
Referenced by fluid_mechanics_routines::fluid_mechanics_problem_setup().
subroutine, public stokes_equations_routines::stokes_problemspecificationset | ( | type(problem_type), pointer | problem, |
integer(intg), dimension(:), intent(in) | problemSpecification, | ||
integer(intg), intent(out) | err, | ||
type(varying_string), intent(out) | error | ||
) |
Sets the problem specification for a Stokes fluid problem.
problem | A pointer to the problem to set the problem specification for | |
[in] | problemspecification | The problem specification to set |
[out] | err | The error code |
[out] | error | The error string |
Definition at line 925 of file Stokes_equations_routines.f90.
References base_routines::enters(), base_routines::exits(), problem_constants::problem_ale_stokes_subtype, problem_constants::problem_fluid_mechanics_class, problem_constants::problem_laplace_stokes_subtype, problem_constants::problem_optimised_stokes_subtype, problem_constants::problem_pgm_stokes_subtype, problem_constants::problem_static_stokes_subtype, problem_constants::problem_stokes_equation_type, and problem_constants::problem_transient_stokes_subtype.
Referenced by fluid_mechanics_routines::fluidmechanics_problemspecificationset().