45 #include "libmesh/libmesh.h"
46 #include "libmesh/mesh.h"
47 #include "libmesh/mesh_generation.h"
48 #include "libmesh/exodusII_io.h"
49 #include "libmesh/equation_systems.h"
50 #include "libmesh/fe.h"
51 #include "libmesh/quadrature_gauss.h"
52 #include "libmesh/dof_map.h"
53 #include "libmesh/sparse_matrix.h"
54 #include "libmesh/numeric_vector.h"
55 #include "libmesh/dense_matrix.h"
56 #include "libmesh/dense_vector.h"
57 #include "libmesh/linear_implicit_system.h"
58 #include "libmesh/transient_system.h"
59 #include "libmesh/perf_log.h"
60 #include "libmesh/boundary_info.h"
61 #include "libmesh/utility.h"
62 #include "libmesh/dirichlet_boundaries.h"
63 #include "libmesh/zero_function.h"
64 #include "libmesh/const_function.h"
65 #include "libmesh/parsed_function.h"
66 #include "libmesh/enum_solver_package.h"
72 #include "libmesh/dense_submatrix.h"
73 #include "libmesh/dense_subvector.h"
76 #include "libmesh/elem.h"
84 const std::string & system_name);
92 int main (
int argc,
char** argv)
99 "--enable-petsc, --enable-trilinos, or --enable-eigen");
102 libmesh_example_requires(2 <= LIBMESH_DIM,
"2D support");
105 #ifndef LIBMESH_ENABLE_DIRICHLET
106 libmesh_example_requires(
false,
"--enable-dirichlet");
143 system.add_variable (
"vel_x",
SECOND);
144 system.add_variable (
"vel_y",
SECOND);
149 system.add_variable (
"p",
FIRST);
161 equation_systems.
init ();
167 PerfLog perf_log(
"Systems Example 2");
176 navier_stokes_system.time = 0.0;
177 const unsigned int n_timesteps = 15;
181 const unsigned int n_nonlinear_steps = 15;
182 const Real nonlinear_tolerance = 1.e-5;
186 equation_systems.
parameters.
set<
unsigned int>(
"linear solver maximum iterations") = 250;
199 std::unique_ptr<NumericVector<Number>>
200 last_nonlinear_soln (navier_stokes_system.solution->clone());
205 for (
unsigned int t_step=1; t_step<=n_timesteps; ++t_step)
209 navier_stokes_system.time += dt;
215 << navier_stokes_system.time
222 *navier_stokes_system.
old_local_solution = *navier_stokes_system.current_local_solution;
226 const Real initial_linear_solver_tol = 1.e-6;
227 equation_systems.
parameters.
set<
Real> (
"linear solver tolerance") = initial_linear_solver_tol;
230 bool converged =
false;
233 for (
unsigned int l=0; l<n_nonlinear_steps; ++l)
236 last_nonlinear_soln->zero();
237 last_nonlinear_soln->add(*navier_stokes_system.solution);
240 perf_log.
push(
"linear solve");
241 equation_systems.
get_system(
"Navier-Stokes").solve();
242 perf_log.
pop(
"linear solve");
246 last_nonlinear_soln->add (-1., *navier_stokes_system.solution);
249 last_nonlinear_soln->close();
252 const Real norm_delta = last_nonlinear_soln->l2_norm();
255 const unsigned int n_linear_iterations = navier_stokes_system.n_linear_iterations();
258 const Real final_linear_residual = navier_stokes_system.final_linear_residual();
273 if (n_linear_iterations == 0 &&
274 (navier_stokes_system.final_linear_residual() >= nonlinear_tolerance || l==0))
276 Real old_linear_solver_tolerance = equation_systems.
parameters.
get<
Real> (
"linear solver tolerance");
277 equation_systems.
parameters.
set<
Real> (
"linear solver tolerance") = 1.e-3 * old_linear_solver_tolerance;
284 << n_linear_iterations
285 <<
", final residual: "
286 << final_linear_residual
287 <<
" Nonlinear convergence: ||u - u_old|| = "
294 if ((norm_delta < nonlinear_tolerance) &&
295 (navier_stokes_system.final_linear_residual() < nonlinear_tolerance))
309 Real new_linear_solver_tolerance = std::min(Utility::pow<2>(final_linear_residual), initial_linear_solver_tol);
310 equation_systems.
parameters.
set<
Real> (
"linear solver tolerance") = new_linear_solver_tolerance;
315 libmesh_error_msg(
"Error: Newton iterations failed to converge!");
317 #ifdef LIBMESH_HAVE_EXODUS_API
319 const unsigned int write_interval = 1;
321 if ((t_step+1)%write_interval == 0)
323 exo_io.write_timestep(
"out.e",
326 navier_stokes_system.time);
328 #endif // #ifdef LIBMESH_HAVE_EXODUS_API
343 const std::string & libmesh_dbg_var(system_name))
347 libmesh_assert_equal_to (system_name,
"Navier-Stokes");
361 const unsigned int u_var = navier_stokes_system.variable_number (
"vel_x");
362 const unsigned int v_var = navier_stokes_system.variable_number (
"vel_y");
363 const unsigned int p_var = navier_stokes_system.variable_number (
"p");
367 FEType fe_vel_type = navier_stokes_system.variable_type(u_var);
370 FEType fe_pres_type = navier_stokes_system.variable_type(p_var);
385 fe_vel->attach_quadrature_rule (&qrule);
386 fe_pres->attach_quadrature_rule (&qrule);
392 const std::vector<Real> & JxW = fe_vel->get_JxW();
395 const std::vector<std::vector<Real>> & phi = fe_vel->get_phi();
399 const std::vector<std::vector<RealGradient>> & dphi = fe_vel->get_dphi();
403 const std::vector<std::vector<Real>> & psi = fe_pres->get_phi();
412 const DofMap & dof_map = navier_stokes_system.get_dof_map();
422 Kuu(Ke), Kuv(Ke), Kup(Ke),
423 Kvu(Ke), Kvv(Ke), Kvp(Ke),
424 Kpu(Ke), Kpv(Ke), Kpp(Ke);
434 std::vector<dof_id_type> dof_indices;
435 std::vector<dof_id_type> dof_indices_u;
436 std::vector<dof_id_type> dof_indices_v;
437 std::vector<dof_id_type> dof_indices_p;
450 const Real theta = 1.;
458 const bool pin_pressure = es.
parameters.
get<
bool>(
"pin_pressure");
476 const unsigned int n_dofs = dof_indices.size();
477 const unsigned int n_u_dofs = dof_indices_u.size();
478 const unsigned int n_v_dofs = dof_indices_v.size();
479 const unsigned int n_p_dofs = dof_indices_p.size();
485 fe_vel->reinit (elem);
486 fe_pres->reinit (elem);
494 Ke.
resize (n_dofs, n_dofs);
510 Kuu.reposition (u_var*n_u_dofs, u_var*n_u_dofs, n_u_dofs, n_u_dofs);
511 Kuv.reposition (u_var*n_u_dofs, v_var*n_u_dofs, n_u_dofs, n_v_dofs);
512 Kup.reposition (u_var*n_u_dofs, p_var*n_u_dofs, n_u_dofs, n_p_dofs);
514 Kvu.reposition (v_var*n_v_dofs, u_var*n_v_dofs, n_v_dofs, n_u_dofs);
515 Kvv.reposition (v_var*n_v_dofs, v_var*n_v_dofs, n_v_dofs, n_v_dofs);
516 Kvp.reposition (v_var*n_v_dofs, p_var*n_v_dofs, n_v_dofs, n_p_dofs);
518 Kpu.reposition (p_var*n_u_dofs, u_var*n_u_dofs, n_p_dofs, n_u_dofs);
519 Kpv.reposition (p_var*n_u_dofs, v_var*n_u_dofs, n_p_dofs, n_v_dofs);
520 Kpp.
reposition (p_var*n_u_dofs, p_var*n_u_dofs, n_p_dofs, n_p_dofs);
522 Fu.reposition (u_var*n_u_dofs, n_u_dofs);
523 Fv.reposition (v_var*n_u_dofs, n_v_dofs);
532 for (
unsigned int qp=0; qp<qrule.n_points(); qp++)
535 Number u = 0., u_old = 0.;
536 Number v = 0., v_old = 0.;
543 for (
unsigned int l=0; l<n_u_dofs; l++)
546 u_old += phi[l][qp]*navier_stokes_system.
old_solution (dof_indices_u[l]);
547 v_old += phi[l][qp]*navier_stokes_system.
old_solution (dof_indices_v[l]);
552 u += phi[l][qp]*navier_stokes_system.current_solution (dof_indices_u[l]);
553 v += phi[l][qp]*navier_stokes_system.current_solution (dof_indices_v[l]);
554 grad_u.
add_scaled (dphi[l][qp],navier_stokes_system.current_solution (dof_indices_u[l]));
555 grad_v.
add_scaled (dphi[l][qp],navier_stokes_system.current_solution (dof_indices_v[l]));
559 for (
unsigned int l=0; l<n_p_dofs; l++)
560 p_old += psi[l][qp]*navier_stokes_system.
old_solution (dof_indices_p[l]);
566 const Number u_x = grad_u(0);
567 const Number u_y = grad_u(1);
568 const Number v_x = grad_v(0);
569 const Number v_y = grad_v(1);
574 for (
unsigned int i=0; i<n_u_dofs; i++)
576 Fu(i) += JxW[qp]*(u_old*phi[i][qp] -
577 (1.-theta)*dt*(U_old*grad_u_old)*phi[i][qp] +
578 (1.-theta)*dt*p_old*dphi[i][qp](0) -
579 (1.-theta)*dt*nu*(grad_u_old*dphi[i][qp]) +
580 theta*dt*(U*grad_u)*phi[i][qp]);
583 Fv(i) += JxW[qp]*(v_old*phi[i][qp] -
584 (1.-theta)*dt*(U_old*grad_v_old)*phi[i][qp] +
585 (1.-theta)*dt*p_old*dphi[i][qp](1) -
586 (1.-theta)*dt*nu*(grad_v_old*dphi[i][qp]) +
587 theta*dt*(U*grad_v)*phi[i][qp]);
594 for (
unsigned int j=0; j<n_u_dofs; j++)
596 Kuu(i,j) += JxW[qp]*(phi[i][qp]*phi[j][qp] +
597 theta*dt*nu*(dphi[i][qp]*dphi[j][qp]) +
598 theta*dt*(U*dphi[j][qp])*phi[i][qp] +
599 theta*dt*u_x*phi[i][qp]*phi[j][qp]);
601 Kuv(i,j) += JxW[qp]*theta*dt*u_y*phi[i][qp]*phi[j][qp];
603 Kvv(i,j) += JxW[qp]*(phi[i][qp]*phi[j][qp] +
604 theta*dt*nu*(dphi[i][qp]*dphi[j][qp]) +
605 theta*dt*(U*dphi[j][qp])*phi[i][qp] +
606 theta*dt*v_y*phi[i][qp]*phi[j][qp]);
608 Kvu(i,j) += JxW[qp]*theta*dt*v_x*phi[i][qp]*phi[j][qp];
612 for (
unsigned int j=0; j<n_p_dofs; j++)
614 Kup(i,j) += JxW[qp]*(-theta*dt*psi[j][qp]*dphi[i][qp](0));
615 Kvp(i,j) += JxW[qp]*(-theta*dt*psi[j][qp]*dphi[i][qp](1));
623 for (
unsigned int i=0; i<n_p_dofs; i++)
624 for (
unsigned int j=0; j<n_u_dofs; j++)
626 Kpu(i,j) += JxW[qp]*psi[i][qp]*dphi[j][qp](0);
627 Kpv(i,j) += JxW[qp]*psi[i][qp]*dphi[j][qp](1);
642 const Real penalty = 1.e10;
643 const unsigned int pressure_node = 0;
644 const Real p_value = 0.0;
645 for (
auto c : elem->node_index_range())
646 if (elem->node_id(c) == pressure_node)
649 Fp(c) += penalty*p_value;
662 navier_stokes_system.matrix->add_matrix (Ke, dof_indices);
663 navier_stokes_system.rhs->add_vector (Fe, dof_indices);
676 system.get_equation_systems().parameters.set<
bool>(
"pin_pressure") =
true;
678 #ifdef LIBMESH_ENABLE_DIRICHLET
680 u_var = system.variable_number(
"vel_x"),
681 v_var = system.variable_number(
"vel_y");
684 DofMap & dof_map = system.get_dof_map();
688 std::set<boundary_id_type> boundary_ids;
689 boundary_ids.insert(0);
690 boundary_ids.insert(1);
691 boundary_ids.insert(3);
693 std::vector<unsigned int> variables;
694 variables.push_back(u_var);
695 variables.push_back(v_var);
703 std::set<boundary_id_type> boundary_ids;
704 boundary_ids.insert(2);
706 std::vector<unsigned int> variables;
707 variables.push_back(u_var);
715 std::set<boundary_id_type> boundary_ids;
716 boundary_ids.insert(2);
718 std::vector<unsigned int> variables;
719 variables.push_back(v_var);
725 #endif // LIBMESH_ENABLE_DIRICHLET
734 system.get_equation_systems().parameters.set<
bool>(
"pin_pressure") =
false;
736 #ifdef LIBMESH_ENABLE_DIRICHLET
738 u_var = system.variable_number(
"vel_x"),
739 v_var = system.variable_number(
"vel_y");
742 DofMap & dof_map = system.get_dof_map();
746 std::set<boundary_id_type> boundary_ids;
747 boundary_ids.insert(0);
749 std::vector<unsigned int> variables;
750 variables.push_back(u_var);
751 variables.push_back(v_var);
759 std::set<boundary_id_type> boundary_ids;
760 boundary_ids.insert(3);
762 std::vector<unsigned int> variables;
763 variables.push_back(u_var);
771 std::set<boundary_id_type> boundary_ids;
772 boundary_ids.insert(2);
774 std::vector<unsigned int> variables;
775 variables.push_back(u_var);
778 std::vector<std::string> additional_vars;
779 additional_vars.push_back(
"k");
780 std::vector<Number> initial_vals;
781 initial_vals.push_back(1.);
791 std::set<boundary_id_type> boundary_ids;
792 boundary_ids.insert(2);
794 std::vector<unsigned int> variables;
795 variables.push_back(v_var);
798 std::vector<std::string> additional_vars;
799 additional_vars.push_back(
"k");
800 std::vector<Number> initial_vals;
801 initial_vals.push_back(1.);
813 #endif // LIBMESH_ENABLE_DIRICHLET
822 system.get_equation_systems().parameters.set<
bool>(
"pin_pressure") =
false;
824 #ifdef LIBMESH_ENABLE_DIRICHLET
826 u_var = system.variable_number(
"vel_x"),
827 v_var = system.variable_number(
"vel_y");
830 DofMap & dof_map = system.get_dof_map();
834 std::set<boundary_id_type> boundary_ids;
835 boundary_ids.insert(0);
836 boundary_ids.insert(2);
838 std::vector<unsigned int> variables;
839 variables.push_back(u_var);
840 variables.push_back(v_var);
848 std::set<boundary_id_type> boundary_ids;
849 boundary_ids.insert(3);
851 std::vector<unsigned int> variables;
852 variables.push_back(u_var);
860 std::set<boundary_id_type> boundary_ids;
861 boundary_ids.insert(3);
863 std::vector<unsigned int> variables;
864 variables.push_back(v_var);
870 #endif // LIBMESH_ENABLE_DIRICHLET