Skip to content

Inviscid Taylor Green Vortex

Man-Long (Norman) Wong edited this page Sep 23, 2019 · 11 revisions

Tutorial 03: Inviscid Taylor-Green Vortex

Tutorial Goal

In this tutorial, the setup of an unsteady simulation in a three-dimensional domain is demonstrated. This is an inviscid, single-species problem. The following capabilities of HAMeRS will be shown:

  • Sixth-order shock capturing method WCNS6-LD
  • Fourth-order SSP-RK54 time integration method
  • Periodic boundary conditions

Problem Description

This is a non-dimensional single-species problem of a decaying vortex named after G.I. Taylor and A.E. Green. The initial conditions of the problem is:

The ratio of specific heats is 5/3. The problem has a size of [0, 2π)3. This is an incompressible problem as the mean pressure is chosen to be very large.

Initial Conditions

The file containing the initial conditions of this problem (TaylorGreenVortex3D.cpp) can be found in problems/Euler/initial_conditions and the initial conditions should be set up by ln -sf <absolute path to TaylorGreenVortex3D.cpp> src/apps/Euler/EulerInitialConditions.cpp. The code has to be re-compiled after the link to the actual initial condition file is set.

Input File Configurations

The configurations of the input file are discussed in this section. Only settings of several important blocks are discussed. The first thing to choose in the input file is whether it is an Euler (inviscid) or Navier-Stokes (diffusive and viscous) application. Since our problem is inviscid, we first set the application type to be Euler:

Application = "Euler"

The next thing to set is the block Euler:

Euler
{
    // Name of project
    project_name = "3D Taylor-Green vortex"

    // Number of species
    num_species = 1

    // Flow model to use
    flow_model = "SINGLE_SPECIES"

    Flow_model
    {
        // Equation of state to use
        equation_of_state = "IDEAL_GAS"

        Equation_of_state_mixing_rules
        {
            species_gamma = 1.666666666666667
            species_R = 1.0
        }
    }

    // Convective flux reconstructor to use
    convective_flux_reconstructor = "WCNS6_LD_HLLC_HLL"

    Convective_flux_reconstructor{}

    Boundary_data{}
}

The settings for Main is very similar to those of tutorial 1 except that dim is set to 3 in the block.

The following settings are used for CartesianGeometry:

CartesianGeometry
{
    // Lower and upper indices of computational domain
    domain_boxes = [(0, 0, 0), (31, 31, 31)]
    x_lo         = 0.0, 0.0, 0.0                                           // Lower end of computational domain
    x_up         = 6.283185307179586, 6.283185307179586, 6.283185307179586 // Upper end of computational domain

    // Periodic dimension. A non-zero value indicates that the direction is periodic
    periodic_dimension = 1, 1, 1
}

Since uniform grid is used in this tutorial, the block ExtendedTagAndInitialize is empty and max_levels is set to 1 in PatchHierarchy:

PatchHierarchy
{
    // Maximum number of levels in hierarchy
    max_levels = 1

    ratio_to_coarser
    {}

    largest_patch_size
    {
        level_0 = 1000, 1000, 1000
    }

    smallest_patch_size
    {
       level_0 = 4, 4, 4
    }
}

The following settings are used for RungeKuttaLevelIntegrator and TimeRefinementIntegrator:

RungeKuttaLevelIntegrator
{
    cfl                      = 0.6e0 // Max cfl factor used in problem
    cfl_init                 = 0.6e0 // Initial cfl factor
    lag_dt_computation       = FALSE
    use_ghosts_to_compute_dt = TRUE

   // Weights of Runge-Kutta scheme
   RungeKuttaWeights
   {
      number_steps   = 5

      alpha_0        = 1.0
      alpha_1        = 0.444370493651235, 0.555629506348765
      alpha_2        = 0.620101851488403, 0.0, 0.379898148511597
      alpha_3        = 0.178079954393132, 0.0, 0.0, 0.821920045606868
      alpha_4        = 0.0, 0.0, 0.517231671970585, 0.096059710526147, 0.386708617503269

      beta_0         = 0.391752226571890
      beta_1         = 0.0, 0.368410593050371
      beta_2         = 0.0, 0.0, 0.251891774271694
      beta_3         = 0.0, 0.0, 0.0, 0.544974750228521
      beta_4         = 0.0, 0.0, 0.0, 0.063692468666290, 0.226007483236906

      gamma_0        = 0.146811876084787
      gamma_1        = 0.0, 0.248482909444976
      gamma_2        = 0.0, 0.0, 0.104258830331981
      gamma_3        = 0.0, 0.0, 0.0, 0.210746432235061
      gamma_4        = 0.0, 0.0, 0.0, 0.063692468666290, 0.226007483236906
   }
}

TimeRefinementIntegrator
{
    start_time           = 0.0e0   // Initial simulation time
    end_time             = 5.0e0   // Final simulation time
    grow_dt              = 1.0e0   // Growth factor for timesteps
    max_integrator_steps = 100000  // Max number of simulation timesteps
}

The RungeKuttaWeights inside RungeKuttaLevelIntegrator are specified to replace the default TVD-RK3 scheme with SSP-RK54 scheme.

The settings for GriddingAlgorithm and TimerManager are very similar to those of tutorial 1 and are not discussed here.

Running HAMeRS

To run the simulation with one core, first put the input file inside a directory named tutorial_03 under HAMeRS. Then, execute the built main executable inside build/src/exec/main with the input file:

../build/exec/main <input filename>

To run the simulation with multiple cores, you can try mpirun/mpiexec/srun, depending on the MPI library used for the compilation of HAMeRS. For example, if mpirun is used:

mpirun -np <number of processors> ../build/src/exec/main <input filename>

Visualization Using VisIt

The data can be visualized by opening dumps.visit inside the viz folder with VisIt. The figures below show the magnitude of vorticity field at different times generated from VisIt: