This tutorial can be run as a Jupyter notebook. Try it on our cloud resources!

Compressible Flow Simulation

Learning objectives

Welcome to the tutorial on the Compressible Flow simulation in Nektar++ framework. This tutorial aims to show how to setup and run a simulation for inviscid compressible flows (governed by Euler equations) using a well known geometry of NACA0012 airfoil.

After completetion of this tutorial you will be familiar with:

  • Configuring a session file suitable for compressible Euler equations in Nektar++
  • Solving the compressible flows explicitly and implicitly
**Note:** It is suggested to download the cases to your local machine or cluster to complete the tutorial.

Introduction

In this tutorial we walk through the user how to configure a session file for solving the compressible Euler equations using a simple flow example. We consider the 2D flow over the NACA0012 airfoil.

Governing equation

The compressible solver implemented in Nektar++ uses the conservative form of the governing equations.

For Navier-Stokes equations: \begin{equation} \label{eq::NS} \frac{\partial \mathbf{Q}}{\partial t} + \frac{\partial \mathbf{F}^c_1}{\partial x_1} + \frac{\partial \mathbf{F}^c_2}{\partial x_2} + \frac{\partial \mathbf{F}^c_3}{\partial x_3} = \frac{\partial \mathbf{F}^v_1}{\partial x_1} + \frac{\partial \mathbf{F}^v_2}{\partial x_2} + \frac{\partial \mathbf{F}^v_3}{\partial x_3} \end{equation}

For Euler equations: \begin{equation} \label{eq::Euler} \frac{\partial \mathbf{Q}}{\partial t} + \frac{\partial \mathbf{F}^c_1}{\partial x_1} + \frac{\partial \mathbf{F}^c_2}{\partial x_2} + \frac{\partial \mathbf{F}^c_3}{\partial x_3} = 0 \end{equation}

In the above $x_i$ denotes the $x$, $y$, and $z$-directions for $i=1,2,3$, respectively, $\mathbf{Q}$ is the vector of conservative variables, $\mathbf{F}^c_i = \mathbf{F}^c_i(\mathbf{Q})$ and $\mathbf{F}^v_i = \mathbf{F}^v_i(\mathbf{Q},\nabla \mathbf{Q})$ are vectors for convective flux and viscous flux, respectively \begin{equation} \label{eq::NS_vectors} \mathbf{Q} =\left( \begin{array}{c} \rho \ \rho u \ \rho v \ \rho w \ E \end{array} \right), \hspace{0.5 cm} \mathbf{F}^c_i =\left( \begin{array}{c} \rho u_i \ \rho u_1 u_i + p \delta_{1i} \ \rho u_2 u_i + p \delta_{2i} \ \rho u_3 u_i + p \delta_{3i} \ u_i (E + p) \end{array} \right), \hspace{0.5 cm} \mathbf{F}^v_i =\left( \begin{array}{c} 0 \ \tau_{1i} \ \tau_{2i} \ \tau_{3i} \ u_j\tau_{ji} + \kappa \partial T / \partial x_i \end{array} \right) \end{equation} where $\rho$ is the density, $u_i$ is the velocity component in the $x_i$-direction, $p$ is the pressure, $\delta_{ij}$ is the Kronecker delta function, and $E$ is the total energy per unit volume. As for the viscous flux vectors, $T$ is the temperature, $\kappa$ is the thermal conductivity, $\tau_{ij}$ is stress tensor component.

Problem description

In this tutorial we simulate a 2D compressible flow over NACA0012 airfoil. The simulation is non-dimesionalized so that the chord length is scaled to a unit, $L=1.0$, and the following freestream parameters are used:

\(\gamma\) \(Ma\) \(R\) \(T_{\infty}\) \(\rho_{\infty}\) \(U_{\infty}\) \(\alpha\)
\(1.4\) \(0.4\) \(\frac{1}{\gamma Ma^2}\) \(1.0\) \(1.0\) \(1.0\) \(0.0\)

The flow domain is a rectangle of sizes $[-5,10] \times [-5,5]$ as shown in the figure below.

Full domain and mesh Near wall mesh

We apply the following boundary conditions (BCs):

  • No-penetration condition on the airfoil surface;
  • Riemann boundary condition at the boundary of the other four boundaries (inflow, outflow, top, and bottom)
**Note:** For the top and bottom boundaries, a periodic boundary condition in the $y$-direction can be used instead.

Session file setup: explicit Euler

In this section we go through the settings in the session file (.xml) for explicit simulations. Please open the session file session_naca_ex.xml. First, let us complete the <EXPANSIONS> tags.

Expansion bases

Task 1: Please write the missing parts in <EXPANSIONS> tag. A completed <EXPANSIONS> tag has the following standard form:

  <EXPANSIONS>
    <E COMPOSITE="C[1]" NUMMODES="2" FIELDS="rho,rhou,rhov,E" TYPE="MODIFIED" /> <!--quad elements-->
    <E COMPOSITE="C[2]" NUMMODES="2" FIELDS="rho,rhou,rhov,E" TYPE="MODIFIED" /> <!--tri elements-->
  </EXPANSIONS>`

Here,

  • COMPOSITE refers to the ID assigned to the geometric domain (or sub-domain) defined in the geometry xml file.
  • NUMMODES is used to indentify the polynomial order of the expansion basis functions as P = NUMMODES-1.
  • TYPE denotes the basis functions. Possible choices are MODIFIED, GLL_LAGRANGE,GLL_LAGRANGE_SEM. One may note that the first basis type uses modal expansion whereas the remaining two consider nodal expansion.

Conditions setup

Now, we will complete the <CONDITIONS> tag. There are five sections in this tag:

  • SolverInfo
  • Parameters
  • Variables
  • Boundary Conditions
  • Initial Conditions

Solver info

Task 2: Please complete the properties EQTYPE, TimeIntegrationMethod, and UpwindType as shown below,

      <I PROPERTY="EQType"                VALUE="EulerCFE"            />
      <I PROPERTY="TimeIntegrationMethod" VALUE="RungeKutta2_SSP"     />
      <I PROPERTY="UpwindType"            VALUE="Roe"                 />

where

  • EQTYPE sets the kind of equations we want to solve. Here we use the Euler equations;
  • TimeIntegrationMethod sets the integration method for time advancing. Here we use the 2nd order Runge-Kutta;
  • UpwindType sets the Riemann solver numerical inviscid flux computation. Here we use the approximate Riemann solver of Roe.

Parameters

In <PARAMETERS> tag, we can define the parameters present in the flow problem as well as the parameters required by various solvers of Nektar++. Hence a parameter can either be defined in the context of the problem we are solving (dummay variable), or be predefined within Nektar++ (solver parameter). One should note that the names of the solver parameters can’t be altered.

Task 3: Assign values to the following parameters as

      <P> Ma                    = 0.4                     </P>
      <P> GasConstant           = 1.0/(Gamma*Ma^2)        </P>
      <P> AoA                   = 0.0/180.0*PI            </P>
      <P> Timestep              = 5.0e-5                  </P>
      <P> NumSteps              = 40000                   </P>    

Here,

  • Ma is the Mach number. It is not necessary but set to compute non-dimensionalized GasConstant conveniently;
  • GasConstant is the gas constant predefined in Nektar++.
  • TimeStep sets the time-step size for the time integration. One should select magnitude of TimeStep strategically to avoid the situation of CFL (Courant-Friedrichs-Lewy) number to be too large. However, in the following implicit session file, this value can be much larger!
  • NumSteps sets the number of time-steps required for the simulation.

Variables

In <VARIABLES> tag, we define the number (and name) of solution variables. We prescribed a unique ID to each variable.

Task 4: Provide a ID for the conservative variable $\rho v$ as

      <V ID="2"> rhov </V>

Boundary conditions

Boundary conditions are defined by two tags. First we define the boundary regions in the <BOUNDARYREGIONS> tag in terms of composite C[ ] entities from the <GEOMETRY> tag section of the mesh xml file. Each boundary region has a unique ID. For the current computational domain there are seven composites (see towards the bottom section of naca.xml file) marked as C[1] to C[7].

Task 5: Define the inflow region using C[4] as shown below

      <B ID="3"> C[4] </B>   <!-- inflow -->

The next tag we define is <BOUNDARYCONDITIONS> tag by which we actually specify the boundary conditions for each boundary ID specified in the above<BOUNDARYREGIONS> tag.

Task 6: Define the the conservative variables at inflow as the farfield quantities as

      <REGION REF="3">
        <D VAR="rho"    VALUE="rhoInf"     />
        <D VAR="rhou"   VALUE="rhoInf * uInf"   />
        <D VAR="rhov"   VALUE="rhoInf * vInf"   />
        <D VAR="E"      VALUE="pInf / (Gamma - 1) + 0.5 * rhoInf * (uInf * uInf + vInf * vInf)"/>
      </REGION>
           

Initial conditions

Task 7: Provide the initial conditions as the freestream conditions in the “InitialConditions” FUNCTION as

        <E VAR="rho"    VALUE="rhoInf"     />
        <E VAR="rhou"   VALUE="rhoInf * uInf"   />
        <E VAR="rhov"   VALUE="rhoInf * vInf"   />
        <E VAR="E"      VALUE="pInf / (Gamma - 1) + 0.5 * rhoInf * (uInf * uInf + vInf * vInf)"/>
                 

Running the solver

The two xml files provided in this tutorial are:

The flow solvers in Nektar++ can be run either in serial or parallel. However, for this case it is suggessted to run in parallel on your local machine. Here we use ! (the jupyter notebook syntax for bash command) to run the solver in this online interactive tutorial which used Docker imange and Binder to run Nektar++ in the background.

Task 8: Go to the next cell and click the Run icon to run the solver. In this example we use 4 physical cores to solve the problem.

**Note:** Use of high number of processors might not work in the current online tutorial. Also, parallel run might be slower than the serial run due to resource constrain.
**Note:** To run the solver in serial, remove the `mpirun -np N` part.
!mpirun -np 4 CompressibleFlowSolver naca.xml session_naca_ex-Complete.xml
Command output
gs_setup: 290 unique labels shared
   pairwise times (avg, min, max): 1.24138e-06 1.2205e-06 1.2558e-06
   crystal router                : 2.245e-06 2.2254e-06 2.27103e-06
   all reduce                    : 2.77909e-06 2.72729e-06 2.8288e-06
   used all_to_all method: pairwise
   handle bytes (avg, min, max): 3756 3092 4516
   buffer bytes (avg, min, max): 2320 1888 2816
gs_setup: 290 unique labels shared
   pairwise times (avg, min, max): 1.15705e-06 1.14245e-06 1.19517e-06
   crystal router                : 2.25196e-06 2.23555e-06 2.25799e-06
   all reduce                    : 2.64696e-06 2.60947e-06 2.66479e-06
   used all_to_all method: pairwise
   handle bytes (avg, min, max): 3756 3092 4516
   buffer bytes (avg, min, max): 2320 1888 2816
gs_setup: 290 unique labels shared
   pairwise times (avg, min, max): 1.26855e-06 1.23549e-06 1.29472e-06
   crystal router                : 2.22756e-06 2.19028e-06 2.25045e-06
   all reduce                    : 2.62342e-06 2.61217e-06 2.6295e-06
   used all_to_all method: pairwise
   handle bytes (avg, min, max): 3756 3092 4516
   buffer bytes (avg, min, max): 2320 1888 2816
gs_setup: 290 unique labels shared
   pairwise times (avg, min, max): 1.189e-06 1.1649e-06 1.21007e-06
   crystal router                : 2.22535e-06 2.20472e-06 2.24514e-06
   all reduce                    : 2.6423e-06 2.62968e-06 2.66181e-06
   used all_to_all method: pairwise
   handle bytes (avg, min, max): 3756 3092 4516
   buffer bytes (avg, min, max): 2320 1888 2816
=======================================================================
	        EquationType: EulerCFE
	        Session Name: naca
	        Spatial Dim.: 2
	  Max SEM Exp. Order: 2
	      Num. Processes: 4
	      Expansion Dim.: 2
	      Riemann Solver: Roe
	      Advection Type: 
	     Projection Type: Discontinuous Galerkin
	           Advection: explicit
	       AdvectionType: WeakDG
	           Diffusion: explicit
	           Time Step: 5e-05
	        No. of Steps: 40000
	 Checkpoints (steps): 5000
	    Integration Type: RungeKutta
=======================================================================
Initial Conditions:
  - Field rho: rhoInf
  - Field rhou: rhoInf * uInf
  - Field rhov: rhoInf * vInf
  - Field E: pInf / (Gamma - 1) + 0.5 * rhoInf * (uInf * uInf + vInf * vInf)
Writing: "naca_0.chk" (0.00510452s, XML)
Steps: 100      Time: 0.005        CPU Time: 0.349357s
Steps: 200      Time: 0.01         CPU Time: 0.349136s
Steps: 300      Time: 0.015        CPU Time: 0.349077s
Steps: 400      Time: 0.02         CPU Time: 0.349361s
Steps: 500      Time: 0.025        CPU Time: 0.349556s
Steps: 600      Time: 0.03         CPU Time: 0.349698s
Steps: 700      Time: 0.035        CPU Time: 0.349086s
Steps: 800      Time: 0.04         CPU Time: 0.349081s
Steps: 900      Time: 0.045        CPU Time: 0.350365s
Steps: 1000     Time: 0.05         CPU Time: 0.349186s
Steps: 1100     Time: 0.055        CPU Time: 0.349222s
Steps: 1200     Time: 0.06         CPU Time: 0.349043s
Steps: 1300     Time: 0.065        CPU Time: 0.349554s
Steps: 1400     Time: 0.07         CPU Time: 0.349886s
Steps: 1500     Time: 0.075        CPU Time: 0.349423s
Steps: 1600     Time: 0.08         CPU Time: 0.350352s
Steps: 1700     Time: 0.085        CPU Time: 0.349456s
Steps: 1800     Time: 0.09         CPU Time: 0.349567s
Steps: 1900     Time: 0.095        CPU Time: 0.349763s
Steps: 2000     Time: 0.1          CPU Time: 0.34946s
Steps: 2100     Time: 0.105        CPU Time: 0.35014s
Steps: 2200     Time: 0.11         CPU Time: 0.349586s
Steps: 2300     Time: 0.115        CPU Time: 0.350087s
Steps: 2400     Time: 0.12         CPU Time: 0.349962s
Steps: 2500     Time: 0.125        CPU Time: 0.349807s
Steps: 2600     Time: 0.13         CPU Time: 0.350099s
Steps: 2700     Time: 0.135        CPU Time: 0.349445s
Steps: 2800     Time: 0.14         CPU Time: 0.350617s
Steps: 2900     Time: 0.145        CPU Time: 0.350219s
Steps: 3000     Time: 0.15         CPU Time: 0.349878s
Steps: 3100     Time: 0.155        CPU Time: 0.350212s
Steps: 3200     Time: 0.16         CPU Time: 0.34964s
Steps: 3300     Time: 0.165        CPU Time: 0.349778s
Steps: 3400     Time: 0.17         CPU Time: 0.349801s
Steps: 3500     Time: 0.175        CPU Time: 0.349752s
Steps: 3600     Time: 0.18         CPU Time: 0.349635s
Steps: 3700     Time: 0.185        CPU Time: 0.349684s
Steps: 3800     Time: 0.19         CPU Time: 0.349837s
Steps: 3900     Time: 0.195        CPU Time: 0.349567s
Steps: 4000     Time: 0.2          CPU Time: 0.349773s
Steps: 4100     Time: 0.205        CPU Time: 0.34964s
Steps: 4200     Time: 0.21         CPU Time: 0.349168s
Steps: 4300     Time: 0.215        CPU Time: 0.348828s
Steps: 4400     Time: 0.22         CPU Time: 0.348461s
Steps: 4500     Time: 0.225        CPU Time: 0.348045s
Steps: 4600     Time: 0.23         CPU Time: 0.34944s
Steps: 4700     Time: 0.235        CPU Time: 0.349537s
Steps: 4800     Time: 0.24         CPU Time: 0.349781s
Steps: 4900     Time: 0.245        CPU Time: 0.34941s
Steps: 5000     Time: 0.25         CPU Time: 0.348858s
renaming "naca_1.chk" -> "naca_1.bak2.chk"
Writing: "naca_1.chk" (0.0152179s, XML)
Steps: 5100     Time: 0.255        CPU Time: 0.350648s
Steps: 5200     Time: 0.26         CPU Time: 0.349934s
Steps: 5300     Time: 0.265        CPU Time: 0.349477s
Steps: 5400     Time: 0.27         CPU Time: 0.348177s
Steps: 5500     Time: 0.275        CPU Time: 0.350254s
Steps: 5600     Time: 0.28         CPU Time: 0.348672s
Steps: 5700     Time: 0.285        CPU Time: 0.349876s
Steps: 5800     Time: 0.29         CPU Time: 0.348738s
Steps: 5900     Time: 0.295        CPU Time: 0.348365s
Steps: 6000     Time: 0.3          CPU Time: 0.34873s
Steps: 6100     Time: 0.305        CPU Time: 0.348408s
Steps: 6200     Time: 0.31         CPU Time: 0.348819s
Steps: 6300     Time: 0.315        CPU Time: 0.348528s
Steps: 6400     Time: 0.32         CPU Time: 0.348705s
Steps: 6500     Time: 0.325        CPU Time: 0.34825s
Steps: 6600     Time: 0.33         CPU Time: 0.348484s
Steps: 6700     Time: 0.335        CPU Time: 0.34927s
Steps: 6800     Time: 0.34         CPU Time: 0.349485s
Steps: 6900     Time: 0.345        CPU Time: 0.35028s
Steps: 7000     Time: 0.35         CPU Time: 0.350081s
Steps: 7100     Time: 0.355        CPU Time: 0.349542s
Steps: 7200     Time: 0.36         CPU Time: 0.350168s
Steps: 7300     Time: 0.365        CPU Time: 0.34987s
Steps: 7400     Time: 0.37         CPU Time: 0.349866s
Steps: 7500     Time: 0.375        CPU Time: 0.34902s
Steps: 7600     Time: 0.38         CPU Time: 0.349551s
Steps: 7700     Time: 0.385        CPU Time: 0.349362s
Steps: 7800     Time: 0.39         CPU Time: 0.349267s
Steps: 7900     Time: 0.395        CPU Time: 0.34897s
Steps: 8000     Time: 0.4          CPU Time: 0.348667s
Steps: 8100     Time: 0.405        CPU Time: 0.348994s
Steps: 8200     Time: 0.41         CPU Time: 0.349999s
Steps: 8300     Time: 0.415        CPU Time: 0.349803s
Steps: 8400     Time: 0.42         CPU Time: 0.349437s
Steps: 8500     Time: 0.425        CPU Time: 0.349582s
Steps: 8600     Time: 0.43         CPU Time: 0.350066s
Steps: 8700     Time: 0.435        CPU Time: 0.349356s
Steps: 8800     Time: 0.44         CPU Time: 0.349774s
Steps: 8900     Time: 0.445        CPU Time: 0.349627s
Steps: 9000     Time: 0.45         CPU Time: 0.349002s
Steps: 9100     Time: 0.455        CPU Time: 0.34986s
Steps: 9200     Time: 0.46         CPU Time: 0.350875s
Steps: 9300     Time: 0.465        CPU Time: 0.350024s
Steps: 9400     Time: 0.47         CPU Time: 0.350568s
Steps: 9500     Time: 0.475        CPU Time: 0.350803s
Steps: 9600     Time: 0.48         CPU Time: 0.350345s
Steps: 9700     Time: 0.485        CPU Time: 0.350916s
Steps: 9800     Time: 0.49         CPU Time: 0.35036s
Steps: 9900     Time: 0.495        CPU Time: 0.349677s
Steps: 10000    Time: 0.5          CPU Time: 0.350555s
renaming "naca_2.chk" -> "naca_2.bak2.chk"
Writing: "naca_2.chk" (0.0151992s, XML)
Steps: 10100    Time: 0.505        CPU Time: 0.349925s
Steps: 10200    Time: 0.51         CPU Time: 0.350372s
Steps: 10300    Time: 0.515        CPU Time: 0.351154s
Steps: 10400    Time: 0.52         CPU Time: 0.349833s
Steps: 10500    Time: 0.525        CPU Time: 0.350052s
Steps: 10600    Time: 0.53         CPU Time: 0.349149s
Steps: 10700    Time: 0.535        CPU Time: 0.350286s
Steps: 10800    Time: 0.54         CPU Time: 0.350235s
Steps: 10900    Time: 0.545        CPU Time: 0.349498s
Steps: 11000    Time: 0.55         CPU Time: 0.352918s
Steps: 11100    Time: 0.555        CPU Time: 0.350615s
Steps: 11200    Time: 0.56         CPU Time: 0.35082s
Steps: 11300    Time: 0.565        CPU Time: 0.352131s
Steps: 11400    Time: 0.57         CPU Time: 0.349375s
Steps: 11500    Time: 0.575        CPU Time: 0.352983s
Steps: 11600    Time: 0.58         CPU Time: 0.351309s
Steps: 11700    Time: 0.585        CPU Time: 0.350926s
Steps: 11800    Time: 0.59         CPU Time: 0.351226s
Steps: 11900    Time: 0.595        CPU Time: 0.350476s
Steps: 12000    Time: 0.6          CPU Time: 0.35233s
Steps: 12100    Time: 0.605        CPU Time: 0.350974s
Steps: 12200    Time: 0.61         CPU Time: 0.350453s
Steps: 12300    Time: 0.615        CPU Time: 0.35218s
Steps: 12400    Time: 0.62         CPU Time: 0.350863s
Steps: 12500    Time: 0.625        CPU Time: 0.352102s
Steps: 12600    Time: 0.63         CPU Time: 0.351867s
Steps: 12700    Time: 0.635        CPU Time: 0.350737s
Steps: 12800    Time: 0.64         CPU Time: 0.351152s
Steps: 12900    Time: 0.645        CPU Time: 0.351376s
Steps: 13000    Time: 0.65         CPU Time: 0.351333s
Steps: 13100    Time: 0.655        CPU Time: 0.351345s
Steps: 13200    Time: 0.66         CPU Time: 0.350793s
Steps: 13300    Time: 0.665        CPU Time: 0.350794s
Steps: 13400    Time: 0.67         CPU Time: 0.35157s
Steps: 13500    Time: 0.675        CPU Time: 0.349503s
Steps: 13600    Time: 0.68         CPU Time: 0.350014s
Steps: 13700    Time: 0.685        CPU Time: 0.351534s
Steps: 13800    Time: 0.69         CPU Time: 0.350154s
Steps: 13900    Time: 0.695        CPU Time: 0.355124s
Steps: 14000    Time: 0.7          CPU Time: 0.358614s
Steps: 14100    Time: 0.705        CPU Time: 0.358251s
Steps: 14200    Time: 0.71         CPU Time: 0.358253s
Steps: 14300    Time: 0.715        CPU Time: 0.357743s
Steps: 14400    Time: 0.72         CPU Time: 0.358442s
Steps: 14500    Time: 0.725        CPU Time: 0.358269s
Steps: 14600    Time: 0.73         CPU Time: 0.357616s
Steps: 14700    Time: 0.735        CPU Time: 0.357628s
Steps: 14800    Time: 0.74         CPU Time: 0.357659s
Steps: 14900    Time: 0.745        CPU Time: 0.360401s
Steps: 15000    Time: 0.75         CPU Time: 0.360811s
renaming "naca_3.chk" -> "naca_3.bak2.chk"
Writing: "naca_3.chk" (0.0161381s, XML)
Steps: 15100    Time: 0.755        CPU Time: 0.358451s
Steps: 15200    Time: 0.76         CPU Time: 0.359347s
Steps: 15300    Time: 0.765        CPU Time: 0.358909s
Steps: 15400    Time: 0.77         CPU Time: 0.358759s
Steps: 15500    Time: 0.775        CPU Time: 0.359998s
Steps: 15600    Time: 0.78         CPU Time: 0.359621s
Steps: 15700    Time: 0.785        CPU Time: 0.359316s
Steps: 15800    Time: 0.79         CPU Time: 0.358893s
Steps: 15900    Time: 0.795        CPU Time: 0.359847s
Steps: 16000    Time: 0.8          CPU Time: 0.36042s
Steps: 16100    Time: 0.805        CPU Time: 0.360746s
Steps: 16200    Time: 0.81         CPU Time: 0.36124s
Steps: 16300    Time: 0.815        CPU Time: 0.35764s
Steps: 16400    Time: 0.82         CPU Time: 0.359763s
Steps: 16500    Time: 0.825        CPU Time: 0.35989s
Steps: 16600    Time: 0.83         CPU Time: 0.359476s
Steps: 16700    Time: 0.835        CPU Time: 0.359171s
Steps: 16800    Time: 0.84         CPU Time: 0.357132s
Steps: 16900    Time: 0.845        CPU Time: 0.356754s
Steps: 17000    Time: 0.85         CPU Time: 0.356526s
Steps: 17100    Time: 0.855        CPU Time: 0.356483s
Steps: 17200    Time: 0.86         CPU Time: 0.356598s
Steps: 17300    Time: 0.865        CPU Time: 0.356972s
Steps: 17400    Time: 0.87         CPU Time: 0.35819s
Steps: 17500    Time: 0.875        CPU Time: 0.359195s
Steps: 17600    Time: 0.88         CPU Time: 0.358137s
Steps: 17700    Time: 0.885        CPU Time: 0.359437s
Steps: 17800    Time: 0.89         CPU Time: 0.357383s
Steps: 17900    Time: 0.895        CPU Time: 0.356882s
Steps: 18000    Time: 0.9          CPU Time: 0.35685s
Steps: 18100    Time: 0.905        CPU Time: 0.356162s
Steps: 18200    Time: 0.91         CPU Time: 0.35718s
Steps: 18300    Time: 0.915        CPU Time: 0.357885s
Steps: 18400    Time: 0.92         CPU Time: 0.357392s
Steps: 18500    Time: 0.925        CPU Time: 0.357194s
Steps: 18600    Time: 0.93         CPU Time: 0.35693s
Steps: 18700    Time: 0.935        CPU Time: 0.357397s
Steps: 18800    Time: 0.94         CPU Time: 0.357053s
Steps: 18900    Time: 0.945        CPU Time: 0.356886s
Steps: 19000    Time: 0.95         CPU Time: 0.356778s
Steps: 19100    Time: 0.955        CPU Time: 0.356925s
Steps: 19200    Time: 0.96         CPU Time: 0.356846s
Steps: 19300    Time: 0.965        CPU Time: 0.357189s
Steps: 19400    Time: 0.97         CPU Time: 0.356142s
Steps: 19500    Time: 0.975        CPU Time: 0.356854s
Steps: 19600    Time: 0.98         CPU Time: 0.356605s
Steps: 19700    Time: 0.985        CPU Time: 0.357102s
Steps: 19800    Time: 0.99         CPU Time: 0.356594s
Steps: 19900    Time: 0.995        CPU Time: 0.357292s
Steps: 20000    Time: 1            CPU Time: 0.356778s
renaming "naca_4.chk" -> "naca_4.bak2.chk"
Writing: "naca_4.chk" (0.016031s, XML)
Steps: 20100    Time: 1.005        CPU Time: 0.357127s
Steps: 20200    Time: 1.01         CPU Time: 0.357011s
Steps: 20300    Time: 1.015        CPU Time: 0.357135s
Steps: 20400    Time: 1.02         CPU Time: 0.356911s
Steps: 20500    Time: 1.025        CPU Time: 0.356964s
Steps: 20600    Time: 1.03         CPU Time: 0.357102s
Steps: 20700    Time: 1.035        CPU Time: 0.356695s
Steps: 20800    Time: 1.04         CPU Time: 0.35745s
Steps: 20900    Time: 1.045        CPU Time: 0.356872s
Steps: 21000    Time: 1.05         CPU Time: 0.357054s
Steps: 21100    Time: 1.055        CPU Time: 0.356836s
Steps: 21200    Time: 1.06         CPU Time: 0.356851s
Steps: 21300    Time: 1.065        CPU Time: 0.357187s
Steps: 21400    Time: 1.07         CPU Time: 0.356886s
Steps: 21500    Time: 1.075        CPU Time: 0.357206s
Steps: 21600    Time: 1.08         CPU Time: 0.35685s
Steps: 21700    Time: 1.085        CPU Time: 0.35705s
Steps: 21800    Time: 1.09         CPU Time: 0.357189s
Steps: 21900    Time: 1.095        CPU Time: 0.357059s
Steps: 22000    Time: 1.1          CPU Time: 0.35748s
Steps: 22100    Time: 1.105        CPU Time: 0.357155s
Steps: 22200    Time: 1.11         CPU Time: 0.356848s
Steps: 22300    Time: 1.115        CPU Time: 0.357042s
Steps: 22400    Time: 1.12         CPU Time: 0.35718s
Steps: 22500    Time: 1.125        CPU Time: 0.356956s
Steps: 22600    Time: 1.13         CPU Time: 0.356967s
Steps: 22700    Time: 1.135        CPU Time: 0.356425s
Steps: 22800    Time: 1.14         CPU Time: 0.356749s
Steps: 22900    Time: 1.145        CPU Time: 0.35693s
Steps: 23000    Time: 1.15         CPU Time: 0.356737s
Steps: 23100    Time: 1.155        CPU Time: 0.35758s
Steps: 23200    Time: 1.16         CPU Time: 0.356061s
Steps: 23300    Time: 1.165        CPU Time: 0.357441s
Steps: 23400    Time: 1.17         CPU Time: 0.357468s
Steps: 23500    Time: 1.175        CPU Time: 0.357s  
Steps: 23600    Time: 1.18         CPU Time: 0.357155s
Steps: 23700    Time: 1.185        CPU Time: 0.356976s
Steps: 23800    Time: 1.19         CPU Time: 0.358022s
Steps: 23900    Time: 1.195        CPU Time: 0.357293s
Steps: 24000    Time: 1.2          CPU Time: 0.358887s
Steps: 24100    Time: 1.205        CPU Time: 0.357865s
Steps: 24200    Time: 1.21         CPU Time: 0.359367s
Steps: 24300    Time: 1.215        CPU Time: 0.360107s
Steps: 24400    Time: 1.22         CPU Time: 0.359588s
Steps: 24500    Time: 1.225        CPU Time: 0.359722s
Steps: 24600    Time: 1.23         CPU Time: 0.358225s
Steps: 24700    Time: 1.235        CPU Time: 0.360664s
Steps: 24800    Time: 1.24         CPU Time: 0.359166s
Steps: 24900    Time: 1.245        CPU Time: 0.359518s
Steps: 25000    Time: 1.25         CPU Time: 0.358162s
renaming "naca_5.chk" -> "naca_5.bak2.chk"
Writing: "naca_5.chk" (0.016035s, XML)
Steps: 25100    Time: 1.255        CPU Time: 0.358484s
Steps: 25200    Time: 1.26         CPU Time: 0.359304s
Steps: 25300    Time: 1.265        CPU Time: 0.359234s
Steps: 25400    Time: 1.27         CPU Time: 0.360685s
Steps: 25500    Time: 1.275        CPU Time: 0.359273s
Steps: 25600    Time: 1.28         CPU Time: 0.360199s
Steps: 25700    Time: 1.285        CPU Time: 0.359381s
Steps: 25800    Time: 1.29         CPU Time: 0.358243s
Steps: 25900    Time: 1.295        CPU Time: 0.361418s
Steps: 26000    Time: 1.3          CPU Time: 0.359699s
Steps: 26100    Time: 1.305        CPU Time: 0.359781s
Steps: 26200    Time: 1.31         CPU Time: 0.358146s
Steps: 26300    Time: 1.315        CPU Time: 0.360716s
Steps: 26400    Time: 1.32         CPU Time: 0.358729s
Steps: 26500    Time: 1.325        CPU Time: 0.358628s
Steps: 26600    Time: 1.33         CPU Time: 0.360351s
Steps: 26700    Time: 1.335        CPU Time: 0.358981s
Steps: 26800    Time: 1.34         CPU Time: 0.358519s
Steps: 26900    Time: 1.345        CPU Time: 0.360991s
Steps: 27000    Time: 1.35         CPU Time: 0.357483s
Steps: 27100    Time: 1.355        CPU Time: 0.360148s
Steps: 27200    Time: 1.36         CPU Time: 0.357616s
Steps: 27300    Time: 1.365        CPU Time: 0.35822s
Steps: 27400    Time: 1.37         CPU Time: 0.358493s
Steps: 27500    Time: 1.375        CPU Time: 0.359454s
Steps: 27600    Time: 1.38         CPU Time: 0.358943s
Steps: 27700    Time: 1.385        CPU Time: 0.357636s
Steps: 27800    Time: 1.39         CPU Time: 0.357034s
Steps: 27900    Time: 1.395        CPU Time: 0.35806s
Steps: 28000    Time: 1.4          CPU Time: 0.357591s
Steps: 28100    Time: 1.405        CPU Time: 0.357353s
Steps: 28200    Time: 1.41         CPU Time: 0.357388s
Steps: 28300    Time: 1.415        CPU Time: 0.357321s
Steps: 28400    Time: 1.42         CPU Time: 0.36019s
Steps: 28500    Time: 1.425        CPU Time: 0.358763s
Steps: 28600    Time: 1.43         CPU Time: 0.357987s
Steps: 28700    Time: 1.435        CPU Time: 0.360043s
Steps: 28800    Time: 1.44         CPU Time: 0.359679s
Steps: 28900    Time: 1.445        CPU Time: 0.358079s
Steps: 29000    Time: 1.45         CPU Time: 0.359054s
Steps: 29100    Time: 1.455        CPU Time: 0.360043s
Steps: 29200    Time: 1.46         CPU Time: 0.359445s
Steps: 29300    Time: 1.465        CPU Time: 0.358676s
Steps: 29400    Time: 1.47         CPU Time: 0.358939s
Steps: 29500    Time: 1.475        CPU Time: 0.359122s
Steps: 29600    Time: 1.48         CPU Time: 0.359515s
Steps: 29700    Time: 1.485        CPU Time: 0.35973s
Steps: 29800    Time: 1.49         CPU Time: 0.35896s
Steps: 29900    Time: 1.495        CPU Time: 0.357453s
Steps: 30000    Time: 1.5          CPU Time: 0.359362s
renaming "naca_6.chk" -> "naca_6.bak2.chk"
Writing: "naca_6.chk" (0.0157765s, XML)
Steps: 30100    Time: 1.505        CPU Time: 0.357887s
Steps: 30200    Time: 1.51         CPU Time: 0.358418s
Steps: 30300    Time: 1.515        CPU Time: 0.358042s
Steps: 30400    Time: 1.52         CPU Time: 0.35905s
Steps: 30500    Time: 1.525        CPU Time: 0.356793s
Steps: 30600    Time: 1.53         CPU Time: 0.356277s
Steps: 30700    Time: 1.535        CPU Time: 0.356594s
Steps: 30800    Time: 1.54         CPU Time: 0.356311s
Steps: 30900    Time: 1.545        CPU Time: 0.356221s
Steps: 31000    Time: 1.55         CPU Time: 0.356733s
Steps: 31100    Time: 1.555        CPU Time: 0.356165s
Steps: 31200    Time: 1.56         CPU Time: 0.357478s
Steps: 31300    Time: 1.565        CPU Time: 0.357009s
Steps: 31400    Time: 1.57         CPU Time: 0.35622s
Steps: 31500    Time: 1.575        CPU Time: 0.356484s
Steps: 31600    Time: 1.58         CPU Time: 0.356478s
Steps: 31700    Time: 1.585        CPU Time: 0.355942s
Steps: 31800    Time: 1.59         CPU Time: 0.356791s
Steps: 31900    Time: 1.595        CPU Time: 0.356379s
Steps: 32000    Time: 1.6          CPU Time: 0.357554s
Steps: 32100    Time: 1.605        CPU Time: 0.358097s
Steps: 32200    Time: 1.61         CPU Time: 0.358666s
Steps: 32300    Time: 1.615        CPU Time: 0.357418s
Steps: 32400    Time: 1.62         CPU Time: 0.358085s
Steps: 32500    Time: 1.625        CPU Time: 0.358368s
Steps: 32600    Time: 1.63         CPU Time: 0.359139s
Steps: 32700    Time: 1.635        CPU Time: 0.358828s
Steps: 32800    Time: 1.64         CPU Time: 0.357191s
Steps: 32900    Time: 1.645        CPU Time: 0.359192s
Steps: 33000    Time: 1.65         CPU Time: 0.359737s
Steps: 33100    Time: 1.655        CPU Time: 0.358535s
Steps: 33200    Time: 1.66         CPU Time: 0.357955s
Steps: 33300    Time: 1.665        CPU Time: 0.357641s
Steps: 33400    Time: 1.67         CPU Time: 0.359588s
Steps: 33500    Time: 1.675        CPU Time: 0.358668s
Steps: 33600    Time: 1.68         CPU Time: 0.357682s
Steps: 33700    Time: 1.685        CPU Time: 0.359248s
Steps: 33800    Time: 1.69         CPU Time: 0.357286s
Steps: 33900    Time: 1.695        CPU Time: 0.357686s
Steps: 34000    Time: 1.7          CPU Time: 0.358391s
Steps: 34100    Time: 1.705        CPU Time: 0.359018s
Steps: 34200    Time: 1.71         CPU Time: 0.359797s
Steps: 34300    Time: 1.715        CPU Time: 0.357732s
Steps: 34400    Time: 1.72         CPU Time: 0.358386s
Steps: 34500    Time: 1.725        CPU Time: 0.357589s
Steps: 34600    Time: 1.73         CPU Time: 0.357943s
Steps: 34700    Time: 1.735        CPU Time: 0.358335s
Steps: 34800    Time: 1.74         CPU Time: 0.358036s
Steps: 34900    Time: 1.745        CPU Time: 0.357502s
Steps: 35000    Time: 1.75         CPU Time: 0.358375s
renaming "naca_7.chk" -> "naca_7.bak2.chk"
Writing: "naca_7.chk" (0.0162489s, XML)
Steps: 35100    Time: 1.755        CPU Time: 0.359081s
Steps: 35200    Time: 1.76         CPU Time: 0.358648s
Steps: 35300    Time: 1.765        CPU Time: 0.35904s
Steps: 35400    Time: 1.77         CPU Time: 0.358688s
Steps: 35500    Time: 1.775        CPU Time: 0.35936s
Steps: 35600    Time: 1.78         CPU Time: 0.358047s
Steps: 35700    Time: 1.785        CPU Time: 0.358759s
Steps: 35800    Time: 1.79         CPU Time: 0.358351s
Steps: 35900    Time: 1.795        CPU Time: 0.358483s
Steps: 36000    Time: 1.8          CPU Time: 0.359007s
Steps: 36100    Time: 1.805        CPU Time: 0.35696s
Steps: 36200    Time: 1.81         CPU Time: 0.358293s
Steps: 36300    Time: 1.815        CPU Time: 0.357488s
Steps: 36400    Time: 1.82         CPU Time: 0.358266s
Steps: 36500    Time: 1.825        CPU Time: 0.357542s
Steps: 36600    Time: 1.83         CPU Time: 0.358586s
Steps: 36700    Time: 1.835        CPU Time: 0.357797s
Steps: 36800    Time: 1.84         CPU Time: 0.357786s
Steps: 36900    Time: 1.845        CPU Time: 0.357564s
Steps: 37000    Time: 1.85         CPU Time: 0.357686s
Steps: 37100    Time: 1.855        CPU Time: 0.357925s
Steps: 37200    Time: 1.86         CPU Time: 0.358379s
Steps: 37300    Time: 1.865        CPU Time: 0.35782s
Steps: 37400    Time: 1.87         CPU Time: 0.357977s
Steps: 37500    Time: 1.875        CPU Time: 0.358187s
Steps: 37600    Time: 1.88         CPU Time: 0.357385s
Steps: 37700    Time: 1.885        CPU Time: 0.357549s
Steps: 37800    Time: 1.89         CPU Time: 0.357892s
Steps: 37900    Time: 1.895        CPU Time: 0.357774s
Steps: 38000    Time: 1.9          CPU Time: 0.358141s
Steps: 38100    Time: 1.905        CPU Time: 0.357687s
Steps: 38200    Time: 1.91         CPU Time: 0.357519s
Steps: 38300    Time: 1.915        CPU Time: 0.357101s
Steps: 38400    Time: 1.92         CPU Time: 0.35869s
Steps: 38500    Time: 1.925        CPU Time: 0.357695s
Steps: 38600    Time: 1.93         CPU Time: 0.357515s
Steps: 38700    Time: 1.935        CPU Time: 0.357765s
Steps: 38800    Time: 1.94         CPU Time: 0.357611s
Steps: 38900    Time: 1.945        CPU Time: 0.358236s
Steps: 39000    Time: 1.95         CPU Time: 0.357267s
Steps: 39100    Time: 1.955        CPU Time: 0.357465s
Steps: 39200    Time: 1.96         CPU Time: 0.35783s
Steps: 39300    Time: 1.965        CPU Time: 0.357324s
Steps: 39400    Time: 1.97         CPU Time: 0.357603s
Steps: 39500    Time: 1.975        CPU Time: 0.35779s
Steps: 39600    Time: 1.98         CPU Time: 0.357175s
Steps: 39700    Time: 1.985        CPU Time: 0.35737s
Steps: 39800    Time: 1.99         CPU Time: 0.357758s
Steps: 39900    Time: 1.995        CPU Time: 0.360275s
Steps: 40000    Time: 2            CPU Time: 0.35823s
renaming "naca_8.chk" -> "naca_8.bak1.chk"
Writing: "naca_8.chk" (0.0158526s, XML)
Time-integration  : 142.1s
renaming "naca.fld" -> "naca.bak4.fld"
Writing: "naca.fld" (0.0158572s, XML)
-------------------------------------------
Total Computation Time = 158s
-------------------------------------------
L 2 error (variable rho) : 12.2441
L inf error (variable rho) : 1.08473
L 2 error (variable rhou) : 12.2431
L inf error (variable rhou) : 1.40827
L 2 error (variable rhov) : 0.0787936
L inf error (variable rhov) : 0.617216
L 2 error (variable E) : 142.775
L inf error (variable E) : 12.6828

Post-processing

To visualise the flowfield first we convert the .fld file (generated by the solver) to a .vtu format file by using the FieldConvert ulitity available in the Nektar++.

Task 9: Run the following code cell to convet the .fld into a .vtu file.

!FieldConvert -f naca.xml session_naca_ex-Complete.xml naca.fld flowfield_ex.vtu
Command output
Writing: "flowfield.vtu"
Written file: flowfield.vtu
**Note:** `-f` is supplied to force output, otherwise the user would be prompted whether to overwrite an existing file.
**Note:** It is strongly recommended to use `-f` option when one is running the commands in Jupyter notebook cells.

Task 10: Run the following python script which uses pyvista combined with itkwidgets to visualise the flowfield in-browser.

import pyvista as pv

# First read the VTK file
mesh = pv.read('flowfield_ex.vtu')

# Now create an itkwidget to visualise this in-browser.
pl = pv.PlotterITK()
pl.add_mesh(mesh, scalars=mesh['u'], smooth_shading=True)
pl.show(True)
Viewer(geometries=[{'vtkClass': 'vtkPolyData', 'points': {'vtkClass': 'vtkPoints', 'name': '_points', 'numberO…
Viewer(geometries=[{'vtkClass': 'vtkPolyData', 'points': {'vtkClass': 'vtkPoints', 'name': '_points', 'numberO…
The velocity field solving by set `NUMMODES=4`

Session file setup: implicit Euler

In Nektar++ the compressible flows can also be solved implicitly. The majority of the session file is the same as the one for explicit solving while the differences only sit in the SOLVERINFO and PARAMETERS tag. Please open the session file session_naca_im.xml and complete the following steps.

Solver info

Task 11: Please complete the properties EQTYPE, AdvectionAdvancement, and TimeIntegrationMethod as shown below,

      <I PROPERTY="EQType"                VALUE="EulerImplicitCFE"/>
      <I PROPERTY="AdvectionAdvancement"  VALUE="Implicit"        />
      <I PROPERTY="TimeIntegrationMethod" VALUE="DIRKOrder2"      />

where EQTYPE is set to be the implicit type of Euler equations, AdvectionAdvancement and TimeIntegrationMethod of the implicit versions are used accordingly.

Parameters

Since we are now solving the Euler equations implicitly, the restrictions on time step is very much relexed. However, to start from the uniform freestream condition, we still use a smaller time step to keep convergence.

Task 12: Assign values to the following parameters as

      <P> Timestep              = 1.0e-4                  </P>
      <P> NumSteps              = 100                     </P>    

Running the solver: starting

Task 13: Go to the next cell and click the Run icon to run the impicit solver.

!mpirun -np 4 CompressibleFlowSolver naca.xml session_naca_im.xml
Command output
gs_setup: 290 unique labels shared
   pairwise times (avg, min, max): 1.24138e-06 1.2205e-06 1.2558e-06
   crystal router                : 2.245e-06 2.2254e-06 2.27103e-06
   all reduce                    : 2.77909e-06 2.72729e-06 2.8288e-06
   used all_to_all method: pairwise
   handle bytes (avg, min, max): 3756 3092 4516
   buffer bytes (avg, min, max): 2320 1888 2816
gs_setup: 290 unique labels shared
   pairwise times (avg, min, max): 1.15705e-06 1.14245e-06 1.19517e-06
   crystal router                : 2.25196e-06 2.23555e-06 2.25799e-06
   all reduce                    : 2.64696e-06 2.60947e-06 2.66479e-06
   used all_to_all method: pairwise
   handle bytes (avg, min, max): 3756 3092 4516
   buffer bytes (avg, min, max): 2320 1888 2816
gs_setup: 290 unique labels shared
   pairwise times (avg, min, max): 1.26855e-06 1.23549e-06 1.29472e-06
   crystal router                : 2.22756e-06 2.19028e-06 2.25045e-06
   all reduce                    : 2.62342e-06 2.61217e-06 2.6295e-06
   used all_to_all method: pairwise
   handle bytes (avg, min, max): 3756 3092 4516
   buffer bytes (avg, min, max): 2320 1888 2816
gs_setup: 290 unique labels shared
   pairwise times (avg, min, max): 1.189e-06 1.1649e-06 1.21007e-06
   crystal router                : 2.22535e-06 2.20472e-06 2.24514e-06
   all reduce                    : 2.6423e-06 2.62968e-06 2.66181e-06
   used all_to_all method: pairwise
   handle bytes (avg, min, max): 3756 3092 4516
   buffer bytes (avg, min, max): 2320 1888 2816
=======================================================================
	        EquationType: EulerCFE
	        Session Name: naca
	        Spatial Dim.: 2
	  Max SEM Exp. Order: 2
	      Num. Processes: 4
	      Expansion Dim.: 2
	      Riemann Solver: Roe
	      Advection Type: 
	     Projection Type: Discontinuous Galerkin
	           Advection: explicit
	       AdvectionType: WeakDG
	           Diffusion: explicit
	           Time Step: 5e-05
	        No. of Steps: 40000
	 Checkpoints (steps): 5000
	    Integration Type: RungeKutta
=======================================================================
Initial Conditions:
  - Field rho: rhoInf
  - Field rhou: rhoInf * uInf
  - Field rhov: rhoInf * vInf
  - Field E: pInf / (Gamma - 1) + 0.5 * rhoInf * (uInf * uInf + vInf * vInf)
Writing: "naca_0.chk" (0.00510452s, XML)
Steps: 100      Time: 0.005        CPU Time: 0.349357s
Steps: 200      Time: 0.01         CPU Time: 0.349136s
Steps: 300      Time: 0.015        CPU Time: 0.349077s
Steps: 400      Time: 0.02         CPU Time: 0.349361s
Steps: 500      Time: 0.025        CPU Time: 0.349556s
Steps: 600      Time: 0.03         CPU Time: 0.349698s
Steps: 700      Time: 0.035        CPU Time: 0.349086s
Steps: 800      Time: 0.04         CPU Time: 0.349081s
Steps: 900      Time: 0.045        CPU Time: 0.350365s
Steps: 1000     Time: 0.05         CPU Time: 0.349186s
Steps: 1100     Time: 0.055        CPU Time: 0.349222s
Steps: 1200     Time: 0.06         CPU Time: 0.349043s
Steps: 1300     Time: 0.065        CPU Time: 0.349554s
Steps: 1400     Time: 0.07         CPU Time: 0.349886s
Steps: 1500     Time: 0.075        CPU Time: 0.349423s
Steps: 1600     Time: 0.08         CPU Time: 0.350352s
Steps: 1700     Time: 0.085        CPU Time: 0.349456s
Steps: 1800     Time: 0.09         CPU Time: 0.349567s
Steps: 1900     Time: 0.095        CPU Time: 0.349763s
Steps: 2000     Time: 0.1          CPU Time: 0.34946s
Steps: 2100     Time: 0.105        CPU Time: 0.35014s
Steps: 2200     Time: 0.11         CPU Time: 0.349586s
Steps: 2300     Time: 0.115        CPU Time: 0.350087s
Steps: 2400     Time: 0.12         CPU Time: 0.349962s
Steps: 2500     Time: 0.125        CPU Time: 0.349807s
Steps: 2600     Time: 0.13         CPU Time: 0.350099s
Steps: 2700     Time: 0.135        CPU Time: 0.349445s
Steps: 2800     Time: 0.14         CPU Time: 0.350617s
Steps: 2900     Time: 0.145        CPU Time: 0.350219s
Steps: 3000     Time: 0.15         CPU Time: 0.349878s
Steps: 3100     Time: 0.155        CPU Time: 0.350212s
Steps: 3200     Time: 0.16         CPU Time: 0.34964s
Steps: 3300     Time: 0.165        CPU Time: 0.349778s
Steps: 3400     Time: 0.17         CPU Time: 0.349801s
Steps: 3500     Time: 0.175        CPU Time: 0.349752s
Steps: 3600     Time: 0.18         CPU Time: 0.349635s
Steps: 3700     Time: 0.185        CPU Time: 0.349684s
Steps: 3800     Time: 0.19         CPU Time: 0.349837s
Steps: 3900     Time: 0.195        CPU Time: 0.349567s
Steps: 4000     Time: 0.2          CPU Time: 0.349773s
Steps: 4100     Time: 0.205        CPU Time: 0.34964s
Steps: 4200     Time: 0.21         CPU Time: 0.349168s
Steps: 4300     Time: 0.215        CPU Time: 0.348828s
Steps: 4400     Time: 0.22         CPU Time: 0.348461s
Steps: 4500     Time: 0.225        CPU Time: 0.348045s
Steps: 4600     Time: 0.23         CPU Time: 0.34944s
Steps: 4700     Time: 0.235        CPU Time: 0.349537s
Steps: 4800     Time: 0.24         CPU Time: 0.349781s
Steps: 4900     Time: 0.245        CPU Time: 0.34941s
Steps: 5000     Time: 0.25         CPU Time: 0.348858s
renaming "naca_1.chk" -> "naca_1.bak2.chk"
Writing: "naca_1.chk" (0.0152179s, XML)
Steps: 5100     Time: 0.255        CPU Time: 0.350648s
Steps: 5200     Time: 0.26         CPU Time: 0.349934s
Steps: 5300     Time: 0.265        CPU Time: 0.349477s
Steps: 5400     Time: 0.27         CPU Time: 0.348177s
Steps: 5500     Time: 0.275        CPU Time: 0.350254s
Steps: 5600     Time: 0.28         CPU Time: 0.348672s
Steps: 5700     Time: 0.285        CPU Time: 0.349876s
Steps: 5800     Time: 0.29         CPU Time: 0.348738s
Steps: 5900     Time: 0.295        CPU Time: 0.348365s
Steps: 6000     Time: 0.3          CPU Time: 0.34873s
Steps: 6100     Time: 0.305        CPU Time: 0.348408s
Steps: 6200     Time: 0.31         CPU Time: 0.348819s
Steps: 6300     Time: 0.315        CPU Time: 0.348528s
Steps: 6400     Time: 0.32         CPU Time: 0.348705s
Steps: 6500     Time: 0.325        CPU Time: 0.34825s
Steps: 6600     Time: 0.33         CPU Time: 0.348484s
Steps: 6700     Time: 0.335        CPU Time: 0.34927s
Steps: 6800     Time: 0.34         CPU Time: 0.349485s
Steps: 6900     Time: 0.345        CPU Time: 0.35028s
Steps: 7000     Time: 0.35         CPU Time: 0.350081s
Steps: 7100     Time: 0.355        CPU Time: 0.349542s
Steps: 7200     Time: 0.36         CPU Time: 0.350168s
Steps: 7300     Time: 0.365        CPU Time: 0.34987s
Steps: 7400     Time: 0.37         CPU Time: 0.349866s
Steps: 7500     Time: 0.375        CPU Time: 0.34902s
Steps: 7600     Time: 0.38         CPU Time: 0.349551s
Steps: 7700     Time: 0.385        CPU Time: 0.349362s
Steps: 7800     Time: 0.39         CPU Time: 0.349267s
Steps: 7900     Time: 0.395        CPU Time: 0.34897s
Steps: 8000     Time: 0.4          CPU Time: 0.348667s
Steps: 8100     Time: 0.405        CPU Time: 0.348994s
Steps: 8200     Time: 0.41         CPU Time: 0.349999s
Steps: 8300     Time: 0.415        CPU Time: 0.349803s
Steps: 8400     Time: 0.42         CPU Time: 0.349437s
Steps: 8500     Time: 0.425        CPU Time: 0.349582s
Steps: 8600     Time: 0.43         CPU Time: 0.350066s
Steps: 8700     Time: 0.435        CPU Time: 0.349356s
Steps: 8800     Time: 0.44         CPU Time: 0.349774s
Steps: 8900     Time: 0.445        CPU Time: 0.349627s
Steps: 9000     Time: 0.45         CPU Time: 0.349002s
Steps: 9100     Time: 0.455        CPU Time: 0.34986s
Steps: 9200     Time: 0.46         CPU Time: 0.350875s
Steps: 9300     Time: 0.465        CPU Time: 0.350024s
Steps: 9400     Time: 0.47         CPU Time: 0.350568s
Steps: 9500     Time: 0.475        CPU Time: 0.350803s
Steps: 9600     Time: 0.48         CPU Time: 0.350345s
Steps: 9700     Time: 0.485        CPU Time: 0.350916s
Steps: 9800     Time: 0.49         CPU Time: 0.35036s
Steps: 9900     Time: 0.495        CPU Time: 0.349677s
Steps: 10000    Time: 0.5          CPU Time: 0.350555s
renaming "naca_2.chk" -> "naca_2.bak2.chk"
Writing: "naca_2.chk" (0.0151992s, XML)
Steps: 10100    Time: 0.505        CPU Time: 0.349925s
Steps: 10200    Time: 0.51         CPU Time: 0.350372s
Steps: 10300    Time: 0.515        CPU Time: 0.351154s
Steps: 10400    Time: 0.52         CPU Time: 0.349833s
Steps: 10500    Time: 0.525        CPU Time: 0.350052s
Steps: 10600    Time: 0.53         CPU Time: 0.349149s
Steps: 10700    Time: 0.535        CPU Time: 0.350286s
Steps: 10800    Time: 0.54         CPU Time: 0.350235s
Steps: 10900    Time: 0.545        CPU Time: 0.349498s
Steps: 11000    Time: 0.55         CPU Time: 0.352918s
Steps: 11100    Time: 0.555        CPU Time: 0.350615s
Steps: 11200    Time: 0.56         CPU Time: 0.35082s
Steps: 11300    Time: 0.565        CPU Time: 0.352131s
Steps: 11400    Time: 0.57         CPU Time: 0.349375s
Steps: 11500    Time: 0.575        CPU Time: 0.352983s
Steps: 11600    Time: 0.58         CPU Time: 0.351309s
Steps: 11700    Time: 0.585        CPU Time: 0.350926s
Steps: 11800    Time: 0.59         CPU Time: 0.351226s
Steps: 11900    Time: 0.595        CPU Time: 0.350476s
Steps: 12000    Time: 0.6          CPU Time: 0.35233s
Steps: 12100    Time: 0.605        CPU Time: 0.350974s
Steps: 12200    Time: 0.61         CPU Time: 0.350453s
Steps: 12300    Time: 0.615        CPU Time: 0.35218s
Steps: 12400    Time: 0.62         CPU Time: 0.350863s
Steps: 12500    Time: 0.625        CPU Time: 0.352102s
Steps: 12600    Time: 0.63         CPU Time: 0.351867s
Steps: 12700    Time: 0.635        CPU Time: 0.350737s
Steps: 12800    Time: 0.64         CPU Time: 0.351152s
Steps: 12900    Time: 0.645        CPU Time: 0.351376s
Steps: 13000    Time: 0.65         CPU Time: 0.351333s
Steps: 13100    Time: 0.655        CPU Time: 0.351345s
Steps: 13200    Time: 0.66         CPU Time: 0.350793s
Steps: 13300    Time: 0.665        CPU Time: 0.350794s
Steps: 13400    Time: 0.67         CPU Time: 0.35157s
Steps: 13500    Time: 0.675        CPU Time: 0.349503s
Steps: 13600    Time: 0.68         CPU Time: 0.350014s
Steps: 13700    Time: 0.685        CPU Time: 0.351534s
Steps: 13800    Time: 0.69         CPU Time: 0.350154s
Steps: 13900    Time: 0.695        CPU Time: 0.355124s
Steps: 14000    Time: 0.7          CPU Time: 0.358614s
Steps: 14100    Time: 0.705        CPU Time: 0.358251s
Steps: 14200    Time: 0.71         CPU Time: 0.358253s
Steps: 14300    Time: 0.715        CPU Time: 0.357743s
Steps: 14400    Time: 0.72         CPU Time: 0.358442s
Steps: 14500    Time: 0.725        CPU Time: 0.358269s
Steps: 14600    Time: 0.73         CPU Time: 0.357616s
Steps: 14700    Time: 0.735        CPU Time: 0.357628s
Steps: 14800    Time: 0.74         CPU Time: 0.357659s
Steps: 14900    Time: 0.745        CPU Time: 0.360401s
Steps: 15000    Time: 0.75         CPU Time: 0.360811s
renaming "naca_3.chk" -> "naca_3.bak2.chk"
Writing: "naca_3.chk" (0.0161381s, XML)
Steps: 15100    Time: 0.755        CPU Time: 0.358451s
Steps: 15200    Time: 0.76         CPU Time: 0.359347s
Steps: 15300    Time: 0.765        CPU Time: 0.358909s
Steps: 15400    Time: 0.77         CPU Time: 0.358759s
Steps: 15500    Time: 0.775        CPU Time: 0.359998s
Steps: 15600    Time: 0.78         CPU Time: 0.359621s
Steps: 15700    Time: 0.785        CPU Time: 0.359316s
Steps: 15800    Time: 0.79         CPU Time: 0.358893s
Steps: 15900    Time: 0.795        CPU Time: 0.359847s
Steps: 16000    Time: 0.8          CPU Time: 0.36042s
Steps: 16100    Time: 0.805        CPU Time: 0.360746s
Steps: 16200    Time: 0.81         CPU Time: 0.36124s
Steps: 16300    Time: 0.815        CPU Time: 0.35764s
Steps: 16400    Time: 0.82         CPU Time: 0.359763s
Steps: 16500    Time: 0.825        CPU Time: 0.35989s
Steps: 16600    Time: 0.83         CPU Time: 0.359476s
Steps: 16700    Time: 0.835        CPU Time: 0.359171s
Steps: 16800    Time: 0.84         CPU Time: 0.357132s
Steps: 16900    Time: 0.845        CPU Time: 0.356754s
Steps: 17000    Time: 0.85         CPU Time: 0.356526s
Steps: 17100    Time: 0.855        CPU Time: 0.356483s
Steps: 17200    Time: 0.86         CPU Time: 0.356598s
Steps: 17300    Time: 0.865        CPU Time: 0.356972s
Steps: 17400    Time: 0.87         CPU Time: 0.35819s
Steps: 17500    Time: 0.875        CPU Time: 0.359195s
Steps: 17600    Time: 0.88         CPU Time: 0.358137s
Steps: 17700    Time: 0.885        CPU Time: 0.359437s
Steps: 17800    Time: 0.89         CPU Time: 0.357383s
Steps: 17900    Time: 0.895        CPU Time: 0.356882s
Steps: 18000    Time: 0.9          CPU Time: 0.35685s
Steps: 18100    Time: 0.905        CPU Time: 0.356162s
Steps: 18200    Time: 0.91         CPU Time: 0.35718s
Steps: 18300    Time: 0.915        CPU Time: 0.357885s
Steps: 18400    Time: 0.92         CPU Time: 0.357392s
Steps: 18500    Time: 0.925        CPU Time: 0.357194s
Steps: 18600    Time: 0.93         CPU Time: 0.35693s
Steps: 18700    Time: 0.935        CPU Time: 0.357397s
Steps: 18800    Time: 0.94         CPU Time: 0.357053s
Steps: 18900    Time: 0.945        CPU Time: 0.356886s
Steps: 19000    Time: 0.95         CPU Time: 0.356778s
Steps: 19100    Time: 0.955        CPU Time: 0.356925s
Steps: 19200    Time: 0.96         CPU Time: 0.356846s
Steps: 19300    Time: 0.965        CPU Time: 0.357189s
Steps: 19400    Time: 0.97         CPU Time: 0.356142s
Steps: 19500    Time: 0.975        CPU Time: 0.356854s
Steps: 19600    Time: 0.98         CPU Time: 0.356605s
Steps: 19700    Time: 0.985        CPU Time: 0.357102s
Steps: 19800    Time: 0.99         CPU Time: 0.356594s
Steps: 19900    Time: 0.995        CPU Time: 0.357292s
Steps: 20000    Time: 1            CPU Time: 0.356778s
renaming "naca_4.chk" -> "naca_4.bak2.chk"
Writing: "naca_4.chk" (0.016031s, XML)
Steps: 20100    Time: 1.005        CPU Time: 0.357127s
Steps: 20200    Time: 1.01         CPU Time: 0.357011s
Steps: 20300    Time: 1.015        CPU Time: 0.357135s
Steps: 20400    Time: 1.02         CPU Time: 0.356911s
Steps: 20500    Time: 1.025        CPU Time: 0.356964s
Steps: 20600    Time: 1.03         CPU Time: 0.357102s
Steps: 20700    Time: 1.035        CPU Time: 0.356695s
Steps: 20800    Time: 1.04         CPU Time: 0.35745s
Steps: 20900    Time: 1.045        CPU Time: 0.356872s
Steps: 21000    Time: 1.05         CPU Time: 0.357054s
Steps: 21100    Time: 1.055        CPU Time: 0.356836s
Steps: 21200    Time: 1.06         CPU Time: 0.356851s
Steps: 21300    Time: 1.065        CPU Time: 0.357187s
Steps: 21400    Time: 1.07         CPU Time: 0.356886s
Steps: 21500    Time: 1.075        CPU Time: 0.357206s
Steps: 21600    Time: 1.08         CPU Time: 0.35685s
Steps: 21700    Time: 1.085        CPU Time: 0.35705s
Steps: 21800    Time: 1.09         CPU Time: 0.357189s
Steps: 21900    Time: 1.095        CPU Time: 0.357059s
Steps: 22000    Time: 1.1          CPU Time: 0.35748s
Steps: 22100    Time: 1.105        CPU Time: 0.357155s
Steps: 22200    Time: 1.11         CPU Time: 0.356848s
Steps: 22300    Time: 1.115        CPU Time: 0.357042s
Steps: 22400    Time: 1.12         CPU Time: 0.35718s
Steps: 22500    Time: 1.125        CPU Time: 0.356956s
Steps: 22600    Time: 1.13         CPU Time: 0.356967s
Steps: 22700    Time: 1.135        CPU Time: 0.356425s
Steps: 22800    Time: 1.14         CPU Time: 0.356749s
Steps: 22900    Time: 1.145        CPU Time: 0.35693s
Steps: 23000    Time: 1.15         CPU Time: 0.356737s
Steps: 23100    Time: 1.155        CPU Time: 0.35758s
Steps: 23200    Time: 1.16         CPU Time: 0.356061s
Steps: 23300    Time: 1.165        CPU Time: 0.357441s
Steps: 23400    Time: 1.17         CPU Time: 0.357468s
Steps: 23500    Time: 1.175        CPU Time: 0.357s  
Steps: 23600    Time: 1.18         CPU Time: 0.357155s
Steps: 23700    Time: 1.185        CPU Time: 0.356976s
Steps: 23800    Time: 1.19         CPU Time: 0.358022s
Steps: 23900    Time: 1.195        CPU Time: 0.357293s
Steps: 24000    Time: 1.2          CPU Time: 0.358887s
Steps: 24100    Time: 1.205        CPU Time: 0.357865s
Steps: 24200    Time: 1.21         CPU Time: 0.359367s
Steps: 24300    Time: 1.215        CPU Time: 0.360107s
Steps: 24400    Time: 1.22         CPU Time: 0.359588s
Steps: 24500    Time: 1.225        CPU Time: 0.359722s
Steps: 24600    Time: 1.23         CPU Time: 0.358225s
Steps: 24700    Time: 1.235        CPU Time: 0.360664s
Steps: 24800    Time: 1.24         CPU Time: 0.359166s
Steps: 24900    Time: 1.245        CPU Time: 0.359518s
Steps: 25000    Time: 1.25         CPU Time: 0.358162s
renaming "naca_5.chk" -> "naca_5.bak2.chk"
Writing: "naca_5.chk" (0.016035s, XML)
Steps: 25100    Time: 1.255        CPU Time: 0.358484s
Steps: 25200    Time: 1.26         CPU Time: 0.359304s
Steps: 25300    Time: 1.265        CPU Time: 0.359234s
Steps: 25400    Time: 1.27         CPU Time: 0.360685s
Steps: 25500    Time: 1.275        CPU Time: 0.359273s
Steps: 25600    Time: 1.28         CPU Time: 0.360199s
Steps: 25700    Time: 1.285        CPU Time: 0.359381s
Steps: 25800    Time: 1.29         CPU Time: 0.358243s
Steps: 25900    Time: 1.295        CPU Time: 0.361418s
Steps: 26000    Time: 1.3          CPU Time: 0.359699s
Steps: 26100    Time: 1.305        CPU Time: 0.359781s
Steps: 26200    Time: 1.31         CPU Time: 0.358146s
Steps: 26300    Time: 1.315        CPU Time: 0.360716s
Steps: 26400    Time: 1.32         CPU Time: 0.358729s
Steps: 26500    Time: 1.325        CPU Time: 0.358628s
Steps: 26600    Time: 1.33         CPU Time: 0.360351s
Steps: 26700    Time: 1.335        CPU Time: 0.358981s
Steps: 26800    Time: 1.34         CPU Time: 0.358519s
Steps: 26900    Time: 1.345        CPU Time: 0.360991s
Steps: 27000    Time: 1.35         CPU Time: 0.357483s
Steps: 27100    Time: 1.355        CPU Time: 0.360148s
Steps: 27200    Time: 1.36         CPU Time: 0.357616s
Steps: 27300    Time: 1.365        CPU Time: 0.35822s
Steps: 27400    Time: 1.37         CPU Time: 0.358493s
Steps: 27500    Time: 1.375        CPU Time: 0.359454s
Steps: 27600    Time: 1.38         CPU Time: 0.358943s
Steps: 27700    Time: 1.385        CPU Time: 0.357636s
Steps: 27800    Time: 1.39         CPU Time: 0.357034s
Steps: 27900    Time: 1.395        CPU Time: 0.35806s
Steps: 28000    Time: 1.4          CPU Time: 0.357591s
Steps: 28100    Time: 1.405        CPU Time: 0.357353s
Steps: 28200    Time: 1.41         CPU Time: 0.357388s
Steps: 28300    Time: 1.415        CPU Time: 0.357321s
Steps: 28400    Time: 1.42         CPU Time: 0.36019s
Steps: 28500    Time: 1.425        CPU Time: 0.358763s
Steps: 28600    Time: 1.43         CPU Time: 0.357987s
Steps: 28700    Time: 1.435        CPU Time: 0.360043s
Steps: 28800    Time: 1.44         CPU Time: 0.359679s
Steps: 28900    Time: 1.445        CPU Time: 0.358079s
Steps: 29000    Time: 1.45         CPU Time: 0.359054s
Steps: 29100    Time: 1.455        CPU Time: 0.360043s
Steps: 29200    Time: 1.46         CPU Time: 0.359445s
Steps: 29300    Time: 1.465        CPU Time: 0.358676s
Steps: 29400    Time: 1.47         CPU Time: 0.358939s
Steps: 29500    Time: 1.475        CPU Time: 0.359122s
Steps: 29600    Time: 1.48         CPU Time: 0.359515s
Steps: 29700    Time: 1.485        CPU Time: 0.35973s
Steps: 29800    Time: 1.49         CPU Time: 0.35896s
Steps: 29900    Time: 1.495        CPU Time: 0.357453s
Steps: 30000    Time: 1.5          CPU Time: 0.359362s
renaming "naca_6.chk" -> "naca_6.bak2.chk"
Writing: "naca_6.chk" (0.0157765s, XML)
Steps: 30100    Time: 1.505        CPU Time: 0.357887s
Steps: 30200    Time: 1.51         CPU Time: 0.358418s
Steps: 30300    Time: 1.515        CPU Time: 0.358042s
Steps: 30400    Time: 1.52         CPU Time: 0.35905s
Steps: 30500    Time: 1.525        CPU Time: 0.356793s
Steps: 30600    Time: 1.53         CPU Time: 0.356277s
Steps: 30700    Time: 1.535        CPU Time: 0.356594s
Steps: 30800    Time: 1.54         CPU Time: 0.356311s
Steps: 30900    Time: 1.545        CPU Time: 0.356221s
Steps: 31000    Time: 1.55         CPU Time: 0.356733s
Steps: 31100    Time: 1.555        CPU Time: 0.356165s
Steps: 31200    Time: 1.56         CPU Time: 0.357478s
Steps: 31300    Time: 1.565        CPU Time: 0.357009s
Steps: 31400    Time: 1.57         CPU Time: 0.35622s
Steps: 31500    Time: 1.575        CPU Time: 0.356484s
Steps: 31600    Time: 1.58         CPU Time: 0.356478s
Steps: 31700    Time: 1.585        CPU Time: 0.355942s
Steps: 31800    Time: 1.59         CPU Time: 0.356791s
Steps: 31900    Time: 1.595        CPU Time: 0.356379s
Steps: 32000    Time: 1.6          CPU Time: 0.357554s
Steps: 32100    Time: 1.605        CPU Time: 0.358097s
Steps: 32200    Time: 1.61         CPU Time: 0.358666s
Steps: 32300    Time: 1.615        CPU Time: 0.357418s
Steps: 32400    Time: 1.62         CPU Time: 0.358085s
Steps: 32500    Time: 1.625        CPU Time: 0.358368s
Steps: 32600    Time: 1.63         CPU Time: 0.359139s
Steps: 32700    Time: 1.635        CPU Time: 0.358828s
Steps: 32800    Time: 1.64         CPU Time: 0.357191s
Steps: 32900    Time: 1.645        CPU Time: 0.359192s
Steps: 33000    Time: 1.65         CPU Time: 0.359737s
Steps: 33100    Time: 1.655        CPU Time: 0.358535s
Steps: 33200    Time: 1.66         CPU Time: 0.357955s
Steps: 33300    Time: 1.665        CPU Time: 0.357641s
Steps: 33400    Time: 1.67         CPU Time: 0.359588s
Steps: 33500    Time: 1.675        CPU Time: 0.358668s
Steps: 33600    Time: 1.68         CPU Time: 0.357682s
Steps: 33700    Time: 1.685        CPU Time: 0.359248s
Steps: 33800    Time: 1.69         CPU Time: 0.357286s
Steps: 33900    Time: 1.695        CPU Time: 0.357686s
Steps: 34000    Time: 1.7          CPU Time: 0.358391s
Steps: 34100    Time: 1.705        CPU Time: 0.359018s
Steps: 34200    Time: 1.71         CPU Time: 0.359797s
Steps: 34300    Time: 1.715        CPU Time: 0.357732s
Steps: 34400    Time: 1.72         CPU Time: 0.358386s
Steps: 34500    Time: 1.725        CPU Time: 0.357589s
Steps: 34600    Time: 1.73         CPU Time: 0.357943s
Steps: 34700    Time: 1.735        CPU Time: 0.358335s
Steps: 34800    Time: 1.74         CPU Time: 0.358036s
Steps: 34900    Time: 1.745        CPU Time: 0.357502s
Steps: 35000    Time: 1.75         CPU Time: 0.358375s
renaming "naca_7.chk" -> "naca_7.bak2.chk"
Writing: "naca_7.chk" (0.0162489s, XML)
Steps: 35100    Time: 1.755        CPU Time: 0.359081s
Steps: 35200    Time: 1.76         CPU Time: 0.358648s
Steps: 35300    Time: 1.765        CPU Time: 0.35904s
Steps: 35400    Time: 1.77         CPU Time: 0.358688s
Steps: 35500    Time: 1.775        CPU Time: 0.35936s
Steps: 35600    Time: 1.78         CPU Time: 0.358047s
Steps: 35700    Time: 1.785        CPU Time: 0.358759s
Steps: 35800    Time: 1.79         CPU Time: 0.358351s
Steps: 35900    Time: 1.795        CPU Time: 0.358483s
Steps: 36000    Time: 1.8          CPU Time: 0.359007s
Steps: 36100    Time: 1.805        CPU Time: 0.35696s
Steps: 36200    Time: 1.81         CPU Time: 0.358293s
Steps: 36300    Time: 1.815        CPU Time: 0.357488s
Steps: 36400    Time: 1.82         CPU Time: 0.358266s
Steps: 36500    Time: 1.825        CPU Time: 0.357542s
Steps: 36600    Time: 1.83         CPU Time: 0.358586s
Steps: 36700    Time: 1.835        CPU Time: 0.357797s
Steps: 36800    Time: 1.84         CPU Time: 0.357786s
Steps: 36900    Time: 1.845        CPU Time: 0.357564s
Steps: 37000    Time: 1.85         CPU Time: 0.357686s
Steps: 37100    Time: 1.855        CPU Time: 0.357925s
Steps: 37200    Time: 1.86         CPU Time: 0.358379s
Steps: 37300    Time: 1.865        CPU Time: 0.35782s
Steps: 37400    Time: 1.87         CPU Time: 0.357977s
Steps: 37500    Time: 1.875        CPU Time: 0.358187s
Steps: 37600    Time: 1.88         CPU Time: 0.357385s
Steps: 37700    Time: 1.885        CPU Time: 0.357549s
Steps: 37800    Time: 1.89         CPU Time: 0.357892s
Steps: 37900    Time: 1.895        CPU Time: 0.357774s
Steps: 38000    Time: 1.9          CPU Time: 0.358141s
Steps: 38100    Time: 1.905        CPU Time: 0.357687s
Steps: 38200    Time: 1.91         CPU Time: 0.357519s
Steps: 38300    Time: 1.915        CPU Time: 0.357101s
Steps: 38400    Time: 1.92         CPU Time: 0.35869s
Steps: 38500    Time: 1.925        CPU Time: 0.357695s
Steps: 38600    Time: 1.93         CPU Time: 0.357515s
Steps: 38700    Time: 1.935        CPU Time: 0.357765s
Steps: 38800    Time: 1.94         CPU Time: 0.357611s
Steps: 38900    Time: 1.945        CPU Time: 0.358236s
Steps: 39000    Time: 1.95         CPU Time: 0.357267s
Steps: 39100    Time: 1.955        CPU Time: 0.357465s
Steps: 39200    Time: 1.96         CPU Time: 0.35783s
Steps: 39300    Time: 1.965        CPU Time: 0.357324s
Steps: 39400    Time: 1.97         CPU Time: 0.357603s
Steps: 39500    Time: 1.975        CPU Time: 0.35779s
Steps: 39600    Time: 1.98         CPU Time: 0.357175s
Steps: 39700    Time: 1.985        CPU Time: 0.35737s
Steps: 39800    Time: 1.99         CPU Time: 0.357758s
Steps: 39900    Time: 1.995        CPU Time: 0.360275s
Steps: 40000    Time: 2            CPU Time: 0.35823s
renaming "naca_8.chk" -> "naca_8.bak1.chk"
Writing: "naca_8.chk" (0.0158526s, XML)
Time-integration  : 142.1s
renaming "naca.fld" -> "naca.bak4.fld"
Writing: "naca.fld" (0.0158572s, XML)
-------------------------------------------
Total Computation Time = 158s
-------------------------------------------
L 2 error (variable rho) : 12.2441
L inf error (variable rho) : 1.08473
L 2 error (variable rhou) : 12.2431
L inf error (variable rhou) : 1.40827
L 2 error (variable rhov) : 0.0787936
L inf error (variable rhov) : 0.617216
L 2 error (variable E) : 142.775
L inf error (variable E) : 12.6828

Task 14: Use the command in the next block to rename the output naca.fld file as naca_im_start.fld

!mv naca.fld/ naca_im_start.fld/

Now we can use a much larger time step restarting from the previous result. Please follow these step to re-run the solver.

Task 15: Assign values to the following parameters as

      <P> Timestep              = 2.0e-3                  </P>
      <P> NumSteps              = 1000                    </P>    

Task 16: Set the “initialConditions” FUNCTION as the output of the starting

      <FUNCTION NAME="InitialConditions">
        <F VAR="rho,rhou,rhov,E"    FILE="naca_im_start.fld"/>
      </FUNCTION>           

Task 17: Go to the next cell and click the Run icon to restart the impicit solver.

!mpirun -np 4 CompressibleFlowSolver naca.xml session_naca_im.xml
Command output
gs_setup: 290 unique labels shared
   pairwise times (avg, min, max): 1.24732e-06 1.215e-06 1.27451e-06
   crystal router                : 2.30079e-06 2.26684e-06 2.32877e-06
   all reduce                    : 2.83667e-06 2.79583e-06 2.86503e-06
   used all_to_all method: pairwise
   handle bytes (avg, min, max): 3756 3092 4516
   buffer bytes (avg, min, max): 2320 1888 2816
gs_setup: 290 unique labels shared
   pairwise times (avg, min, max): 1.30166e-06 1.27563e-06 1.32388e-06
   crystal router                : 2.29748e-06 2.26544e-06 2.33203e-06
   all reduce                    : 2.64156e-06 2.61944e-06 2.65921e-06
   used all_to_all method: pairwise
   handle bytes (avg, min, max): 3756 3092 4516
   buffer bytes (avg, min, max): 2320 1888 2816
gs_setup: 290 unique labels shared
   pairwise times (avg, min, max): 1.20776e-06 1.17887e-06 1.22907e-06
   crystal router                : 2.32121e-06 2.30866e-06 2.33734e-06
   all reduce                    : 2.6498e-06 2.63294e-06 2.6634e-06
   used all_to_all method: pairwise
   handle bytes (avg, min, max): 3756 3092 4516
   buffer bytes (avg, min, max): 2320 1888 2816
gs_setup: 290 unique labels shared
   pairwise times (avg, min, max): 1.20704e-06 1.18352e-06 1.23838e-06
   crystal router                : 2.2938e-06 2.26917e-06 2.31992e-06
   all reduce                    : 2.65795e-06 2.64859e-06 2.66964e-06
   used all_to_all method: pairwise
   handle bytes (avg, min, max): 3756 3092 4516
   buffer bytes (avg, min, max): 2320 1888 2816
=======================================================================
	        EquationType: EulerImplicitCFE
	        Session Name: naca
	        Spatial Dim.: 2
	  Max SEM Exp. Order: 2
	      Num. Processes: 4
	      Expansion Dim.: 2
	      Riemann Solver: Roe
	      Advection Type: 
	     Projection Type: Discontinuous Galerkin
	           Advection: implicit
	       AdvectionType: WeakDG
	           Diffusion: explicit
	           Time Step: 0.002
	        No. of Steps: 1000
	 Checkpoints (steps): 500
	    Integration Type: DIRK
=======================================================================
Initial Conditions:
  - Field rho: from file naca_im_start.fld
  - Field rhou: from file naca_im_start.fld
  - Field rhov: from file naca_im_start.fld
  - Field E: from file naca_im_start.fld
renaming "naca_3.chk" -> "naca_3.bak0.chk"
Writing: "naca_3.chk" (0.0277014s, XML)
Steps: 1        Time: 0.012        CPU Time: 0.235347s
Steps: 2        Time: 0.014        CPU Time: 0.12753s
Steps: 3        Time: 0.016        CPU Time: 0.127583s
Steps: 4        Time: 0.018        CPU Time: 0.127929s
Steps: 5        Time: 0.02         CPU Time: 0.127647s
Steps: 6        Time: 0.022        CPU Time: 0.12773s
Steps: 7        Time: 0.024        CPU Time: 0.127867s
Steps: 8        Time: 0.026        CPU Time: 0.128035s
Steps: 9        Time: 0.028        CPU Time: 0.127135s
Steps: 10       Time: 0.03         CPU Time: 0.127493s
Steps: 11       Time: 0.032        CPU Time: 0.127846s
Steps: 12       Time: 0.034        CPU Time: 0.127323s
Steps: 13       Time: 0.036        CPU Time: 0.127665s
Steps: 14       Time: 0.038        CPU Time: 0.127755s
Steps: 15       Time: 0.04         CPU Time: 0.12783s
Steps: 16       Time: 0.042        CPU Time: 0.127658s
Steps: 17       Time: 0.044        CPU Time: 0.127537s
Steps: 18       Time: 0.046        CPU Time: 0.159404s
Steps: 19       Time: 0.048        CPU Time: 0.0658412s
Steps: 20       Time: 0.05         CPU Time: 0.0659865s
Steps: 21       Time: 0.052        CPU Time: 0.065949s
Steps: 22       Time: 0.054        CPU Time: 0.0659689s
Steps: 23       Time: 0.056        CPU Time: 0.0659177s
Steps: 24       Time: 0.058        CPU Time: 0.0660655s
Steps: 25       Time: 0.06         CPU Time: 0.0661115s
Steps: 26       Time: 0.062        CPU Time: 0.066029s
Steps: 27       Time: 0.064        CPU Time: 0.0660165s
Steps: 28       Time: 0.066        CPU Time: 0.0658099s
Steps: 29       Time: 0.068        CPU Time: 0.0660275s
Steps: 30       Time: 0.07         CPU Time: 0.065839s
Steps: 31       Time: 0.072        CPU Time: 0.0658516s
Steps: 32       Time: 0.074        CPU Time: 0.0661354s
Steps: 33       Time: 0.076        CPU Time: 0.0661474s
Steps: 34       Time: 0.078        CPU Time: 0.0658339s
Steps: 35       Time: 0.08         CPU Time: 0.0659601s
Steps: 36       Time: 0.082        CPU Time: 0.0658791s
Steps: 37       Time: 0.084        CPU Time: 0.0658887s
Steps: 38       Time: 0.086        CPU Time: 0.0659342s
Steps: 39       Time: 0.088        CPU Time: 0.0658331s
Steps: 40       Time: 0.09         CPU Time: 0.0662052s
Steps: 41       Time: 0.092        CPU Time: 0.0660277s
Steps: 42       Time: 0.094        CPU Time: 0.065958s
Steps: 43       Time: 0.096        CPU Time: 0.0660508s
Steps: 44       Time: 0.098        CPU Time: 0.0657424s
Steps: 45       Time: 0.1          CPU Time: 0.0660583s
Steps: 46       Time: 0.102        CPU Time: 0.0659107s
Steps: 47       Time: 0.104        CPU Time: 0.0659403s
Steps: 48       Time: 0.106        CPU Time: 0.0657553s
Steps: 49       Time: 0.108        CPU Time: 0.0661316s
Steps: 50       Time: 0.11         CPU Time: 0.0659703s
Steps: 51       Time: 0.112        CPU Time: 0.0657308s
Steps: 52       Time: 0.114        CPU Time: 0.161281s
Steps: 53       Time: 0.116        CPU Time: 0.0658117s
Steps: 54       Time: 0.118        CPU Time: 0.0656717s
Steps: 55       Time: 0.12         CPU Time: 0.0671829s
Steps: 56       Time: 0.122        CPU Time: 0.0666851s
Steps: 57       Time: 0.124        CPU Time: 0.0665443s
Steps: 58       Time: 0.126        CPU Time: 0.0666969s
Steps: 59       Time: 0.128        CPU Time: 0.0670582s
Steps: 60       Time: 0.13         CPU Time: 0.0664842s
Steps: 61       Time: 0.132        CPU Time: 0.0669602s
Steps: 62       Time: 0.134        CPU Time: 0.0669429s
Steps: 63       Time: 0.136        CPU Time: 0.0667499s
Steps: 64       Time: 0.138        CPU Time: 0.0666158s
Steps: 65       Time: 0.14         CPU Time: 0.066709s
Steps: 66       Time: 0.142        CPU Time: 0.0666988s
Steps: 67       Time: 0.144        CPU Time: 0.0669002s
Steps: 68       Time: 0.146        CPU Time: 0.0667011s
Steps: 69       Time: 0.148        CPU Time: 0.0667225s
Steps: 70       Time: 0.15         CPU Time: 0.0665036s
Steps: 71       Time: 0.152        CPU Time: 0.0670465s
Steps: 72       Time: 0.154        CPU Time: 0.0666204s
Steps: 73       Time: 0.156        CPU Time: 0.0664507s
Steps: 74       Time: 0.158        CPU Time: 0.0669407s
Steps: 75       Time: 0.16         CPU Time: 0.0668734s
Steps: 76       Time: 0.162        CPU Time: 0.0667501s
Steps: 77       Time: 0.164        CPU Time: 0.0668897s
Steps: 78       Time: 0.166        CPU Time: 0.0669701s
Steps: 79       Time: 0.168        CPU Time: 0.0665681s
Steps: 80       Time: 0.17         CPU Time: 0.0667746s
Steps: 81       Time: 0.172        CPU Time: 0.0668894s
Steps: 82       Time: 0.174        CPU Time: 0.0667435s
Steps: 83       Time: 0.176        CPU Time: 0.0666187s
Steps: 84       Time: 0.178        CPU Time: 0.0665287s
Steps: 85       Time: 0.18         CPU Time: 0.0666284s
Steps: 86       Time: 0.182        CPU Time: 0.161857s
Steps: 87       Time: 0.184        CPU Time: 0.065984s
Steps: 88       Time: 0.186        CPU Time: 0.0662417s
Steps: 89       Time: 0.188        CPU Time: 0.0663232s
Steps: 90       Time: 0.19         CPU Time: 0.0663985s
Steps: 91       Time: 0.192        CPU Time: 0.0660593s
Steps: 92       Time: 0.194        CPU Time: 0.0663288s
Steps: 93       Time: 0.196        CPU Time: 0.0659721s
Steps: 94       Time: 0.198        CPU Time: 0.0662001s
Steps: 95       Time: 0.2          CPU Time: 0.0659644s
Steps: 96       Time: 0.202        CPU Time: 0.0662228s
Steps: 97       Time: 0.204        CPU Time: 0.0660284s
Steps: 98       Time: 0.206        CPU Time: 0.0657793s
Steps: 99       Time: 0.208        CPU Time: 0.066517s
Steps: 100      Time: 0.21         CPU Time: 0.0662547s
Steps: 101      Time: 0.212        CPU Time: 0.0664966s
Steps: 102      Time: 0.214        CPU Time: 0.0662491s
Steps: 103      Time: 0.216        CPU Time: 0.0662487s
Steps: 104      Time: 0.218        CPU Time: 0.0664906s
Steps: 105      Time: 0.22         CPU Time: 0.0662842s
Steps: 106      Time: 0.222        CPU Time: 0.0664857s
Steps: 107      Time: 0.224        CPU Time: 0.066252s
Steps: 108      Time: 0.226        CPU Time: 0.0661335s
Steps: 109      Time: 0.228        CPU Time: 0.0664098s
Steps: 110      Time: 0.23         CPU Time: 0.0661458s
Steps: 111      Time: 0.232        CPU Time: 0.0662417s
Steps: 112      Time: 0.234        CPU Time: 0.0665046s
Steps: 113      Time: 0.236        CPU Time: 0.0660551s
Steps: 114      Time: 0.238        CPU Time: 0.0664189s
Steps: 115      Time: 0.24         CPU Time: 0.0964345s
Steps: 116      Time: 0.242        CPU Time: 0.12675s
Steps: 117      Time: 0.244        CPU Time: 0.126418s
Steps: 118      Time: 0.246        CPU Time: 0.223103s
Steps: 119      Time: 0.248        CPU Time: 0.126063s
Steps: 120      Time: 0.25         CPU Time: 0.126818s
Steps: 121      Time: 0.252        CPU Time: 0.127748s
Steps: 122      Time: 0.254        CPU Time: 0.127612s
Steps: 123      Time: 0.256        CPU Time: 0.12744s
Steps: 124      Time: 0.258        CPU Time: 0.126892s
Steps: 125      Time: 0.26         CPU Time: 0.12774s
Steps: 126      Time: 0.262        CPU Time: 0.127614s
Steps: 127      Time: 0.264        CPU Time: 0.127329s
Steps: 128      Time: 0.266        CPU Time: 0.127616s
Steps: 129      Time: 0.268        CPU Time: 0.127907s
Steps: 130      Time: 0.27         CPU Time: 0.127683s
Steps: 131      Time: 0.272        CPU Time: 0.129263s
Steps: 132      Time: 0.274        CPU Time: 0.127747s
Steps: 133      Time: 0.276        CPU Time: 0.128006s
Steps: 134      Time: 0.278        CPU Time: 0.127665s
Steps: 135      Time: 0.28         CPU Time: 0.225378s
Steps: 136      Time: 0.282        CPU Time: 0.12606s
Steps: 137      Time: 0.284        CPU Time: 0.0963404s
Steps: 138      Time: 0.286        CPU Time: 0.0657737s
Steps: 139      Time: 0.288        CPU Time: 0.0662037s
Steps: 140      Time: 0.29         CPU Time: 0.0659191s
Steps: 141      Time: 0.292        CPU Time: 0.0661723s
Steps: 142      Time: 0.294        CPU Time: 0.0661655s
Steps: 143      Time: 0.296        CPU Time: 0.0661125s
Steps: 144      Time: 0.298        CPU Time: 0.0660954s
Steps: 145      Time: 0.3          CPU Time: 0.12608s
Steps: 146      Time: 0.302        CPU Time: 0.126199s
Steps: 147      Time: 0.304        CPU Time: 0.126377s
Steps: 148      Time: 0.306        CPU Time: 0.12611s
Steps: 149      Time: 0.308        CPU Time: 0.126243s
Steps: 150      Time: 0.31         CPU Time: 0.126193s
Steps: 151      Time: 0.312        CPU Time: 0.0958777s
Steps: 152      Time: 0.314        CPU Time: 0.066275s
Steps: 153      Time: 0.316        CPU Time: 0.0664031s
Steps: 154      Time: 0.318        CPU Time: 0.0659992s
Steps: 155      Time: 0.32         CPU Time: 0.0660761s
Steps: 156      Time: 0.322        CPU Time: 0.0663321s
Steps: 157      Time: 0.324        CPU Time: 0.0657227s
Steps: 158      Time: 0.326        CPU Time: 0.0660302s
Steps: 159      Time: 0.328        CPU Time: 0.0660509s
Steps: 160      Time: 0.33         CPU Time: 0.164966s
Steps: 161      Time: 0.332        CPU Time: 0.0663287s
Steps: 162      Time: 0.334        CPU Time: 0.0662595s
Steps: 163      Time: 0.336        CPU Time: 0.0673066s
Steps: 164      Time: 0.338        CPU Time: 0.0669844s
Steps: 165      Time: 0.34         CPU Time: 0.0673488s
Steps: 166      Time: 0.342        CPU Time: 0.0670383s
Steps: 167      Time: 0.344        CPU Time: 0.0669789s
Steps: 168      Time: 0.346        CPU Time: 0.0669417s
Steps: 169      Time: 0.348        CPU Time: 0.0669602s
Steps: 170      Time: 0.35         CPU Time: 0.0670245s
Steps: 171      Time: 0.352        CPU Time: 0.0671581s
Steps: 172      Time: 0.354        CPU Time: 0.0671645s
Steps: 173      Time: 0.356        CPU Time: 0.0668446s
Steps: 174      Time: 0.358        CPU Time: 0.0670899s
Steps: 175      Time: 0.36         CPU Time: 0.0670895s
Steps: 176      Time: 0.362        CPU Time: 0.0671282s
Steps: 177      Time: 0.364        CPU Time: 0.0673371s
Steps: 178      Time: 0.366        CPU Time: 0.0670136s
Steps: 179      Time: 0.368        CPU Time: 0.0670415s
Steps: 180      Time: 0.37         CPU Time: 0.06712s
Steps: 181      Time: 0.372        CPU Time: 0.067017s
Steps: 182      Time: 0.374        CPU Time: 0.0671039s
Steps: 183      Time: 0.376        CPU Time: 0.0670293s
Steps: 184      Time: 0.378        CPU Time: 0.0670792s
Steps: 185      Time: 0.38         CPU Time: 0.0670638s
Steps: 186      Time: 0.382        CPU Time: 0.0671263s
Steps: 187      Time: 0.384        CPU Time: 0.0667261s
Steps: 188      Time: 0.386        CPU Time: 0.0667074s
Steps: 189      Time: 0.388        CPU Time: 0.0669614s
Steps: 190      Time: 0.39         CPU Time: 0.067143s
Steps: 191      Time: 0.392        CPU Time: 0.0672564s
Steps: 192      Time: 0.394        CPU Time: 0.0672207s
Steps: 193      Time: 0.396        CPU Time: 0.067894s
Steps: 194      Time: 0.398        CPU Time: 0.165653s
Steps: 195      Time: 0.4          CPU Time: 0.0660099s
Steps: 196      Time: 0.402        CPU Time: 0.0656378s
Steps: 197      Time: 0.404        CPU Time: 0.0655944s
Steps: 198      Time: 0.406        CPU Time: 0.066049s
Steps: 199      Time: 0.408        CPU Time: 0.0659789s
Steps: 200      Time: 0.41         CPU Time: 0.0657404s
Steps: 201      Time: 0.412        CPU Time: 0.0661461s
Steps: 202      Time: 0.414        CPU Time: 0.0659077s
Steps: 203      Time: 0.416        CPU Time: 0.0662013s
Steps: 204      Time: 0.418        CPU Time: 0.066305s
Steps: 205      Time: 0.42         CPU Time: 0.0661088s
Steps: 206      Time: 0.422        CPU Time: 0.0660281s
Steps: 207      Time: 0.424        CPU Time: 0.0657782s
Steps: 208      Time: 0.426        CPU Time: 0.0661167s
Steps: 209      Time: 0.428        CPU Time: 0.0658909s
Steps: 210      Time: 0.43         CPU Time: 0.0665848s
Steps: 211      Time: 0.432        CPU Time: 0.0663525s
Steps: 212      Time: 0.434        CPU Time: 0.0661063s
Steps: 213      Time: 0.436        CPU Time: 0.0659833s
Steps: 214      Time: 0.438        CPU Time: 0.0658159s
Steps: 215      Time: 0.44         CPU Time: 0.0659037s
Steps: 216      Time: 0.442        CPU Time: 0.0658364s
Steps: 217      Time: 0.444        CPU Time: 0.0661783s
Steps: 218      Time: 0.446        CPU Time: 0.0658743s
Steps: 219      Time: 0.448        CPU Time: 0.0672369s
Steps: 220      Time: 0.45         CPU Time: 0.0663279s
Steps: 221      Time: 0.452        CPU Time: 0.0660595s
Steps: 222      Time: 0.454        CPU Time: 0.0659003s
Steps: 223      Time: 0.456        CPU Time: 0.0659885s
Steps: 224      Time: 0.458        CPU Time: 0.0660229s
Steps: 225      Time: 0.46         CPU Time: 0.0666917s
Steps: 226      Time: 0.462        CPU Time: 0.0661015s
Steps: 227      Time: 0.464        CPU Time: 0.0661756s
Steps: 228      Time: 0.466        CPU Time: 0.164416s
Steps: 229      Time: 0.468        CPU Time: 0.0662641s
Steps: 230      Time: 0.47         CPU Time: 0.0671939s
Steps: 231      Time: 0.472        CPU Time: 0.0668863s
Steps: 232      Time: 0.474        CPU Time: 0.0671626s
Steps: 233      Time: 0.476        CPU Time: 0.0672682s
Steps: 234      Time: 0.478        CPU Time: 0.0673009s
Steps: 235      Time: 0.48         CPU Time: 0.067178s
Steps: 236      Time: 0.482        CPU Time: 0.0672038s
Steps: 237      Time: 0.484        CPU Time: 0.0670345s
Steps: 238      Time: 0.486        CPU Time: 0.0672563s
Steps: 239      Time: 0.488        CPU Time: 0.0671691s
Steps: 240      Time: 0.49         CPU Time: 0.0671196s
Steps: 241      Time: 0.492        CPU Time: 0.0669293s
Steps: 242      Time: 0.494        CPU Time: 0.0670788s
Steps: 243      Time: 0.496        CPU Time: 0.0668462s
Steps: 244      Time: 0.498        CPU Time: 0.0670426s
Steps: 245      Time: 0.5          CPU Time: 0.0672086s
Steps: 246      Time: 0.502        CPU Time: 0.0670067s
Steps: 247      Time: 0.504        CPU Time: 0.0670062s
Steps: 248      Time: 0.506        CPU Time: 0.0673925s
Steps: 249      Time: 0.508        CPU Time: 0.0672432s
Steps: 250      Time: 0.51         CPU Time: 0.0668663s
Steps: 251      Time: 0.512        CPU Time: 0.0671902s
Steps: 252      Time: 0.514        CPU Time: 0.0670785s
Steps: 253      Time: 0.516        CPU Time: 0.067241s
Steps: 254      Time: 0.518        CPU Time: 0.06717s
Steps: 255      Time: 0.52         CPU Time: 0.0669337s
Steps: 256      Time: 0.522        CPU Time: 0.0670124s
Steps: 257      Time: 0.524        CPU Time: 0.0672434s
Steps: 258      Time: 0.526        CPU Time: 0.0673057s
Steps: 259      Time: 0.528        CPU Time: 0.067068s
Steps: 260      Time: 0.53         CPU Time: 0.0671737s
Steps: 261      Time: 0.532        CPU Time: 0.0669115s
Steps: 262      Time: 0.534        CPU Time: 0.165073s
Steps: 263      Time: 0.536        CPU Time: 0.0659161s
Steps: 264      Time: 0.538        CPU Time: 0.0660812s
Steps: 265      Time: 0.54         CPU Time: 0.0659327s
Steps: 266      Time: 0.542        CPU Time: 0.0658781s
Steps: 267      Time: 0.544        CPU Time: 0.0660551s
Steps: 268      Time: 0.546        CPU Time: 0.0663633s
Steps: 269      Time: 0.548        CPU Time: 0.0660642s
Steps: 270      Time: 0.55         CPU Time: 0.0660168s
Steps: 271      Time: 0.552        CPU Time: 0.0659331s
Steps: 272      Time: 0.554        CPU Time: 0.0659796s
Steps: 273      Time: 0.556        CPU Time: 0.0659656s
Steps: 274      Time: 0.558        CPU Time: 0.0659303s
Steps: 275      Time: 0.56         CPU Time: 0.0657227s
Steps: 276      Time: 0.562        CPU Time: 0.0661793s
Steps: 277      Time: 0.564        CPU Time: 0.0660566s
Steps: 278      Time: 0.566        CPU Time: 0.0659258s
Steps: 279      Time: 0.568        CPU Time: 0.0660852s
Steps: 280      Time: 0.57         CPU Time: 0.0657204s
Steps: 281      Time: 0.572        CPU Time: 0.0658634s
Steps: 282      Time: 0.574        CPU Time: 0.0662805s
Steps: 283      Time: 0.576        CPU Time: 0.0659634s
Steps: 284      Time: 0.578        CPU Time: 0.0657077s
Steps: 285      Time: 0.58         CPU Time: 0.0661529s
Steps: 286      Time: 0.582        CPU Time: 0.0660262s
Steps: 287      Time: 0.584        CPU Time: 0.0659954s
Steps: 288      Time: 0.586        CPU Time: 0.065787s
Steps: 289      Time: 0.588        CPU Time: 0.0658868s
Steps: 290      Time: 0.59         CPU Time: 0.0659948s
Steps: 291      Time: 0.592        CPU Time: 0.0656549s
Steps: 292      Time: 0.594        CPU Time: 0.0659506s
Steps: 293      Time: 0.596        CPU Time: 0.0662117s
Steps: 294      Time: 0.598        CPU Time: 0.0657969s
Steps: 295      Time: 0.6          CPU Time: 0.0659473s
Steps: 296      Time: 0.602        CPU Time: 0.164246s
Steps: 297      Time: 0.604        CPU Time: 0.065761s
Steps: 298      Time: 0.606        CPU Time: 0.0670737s
Steps: 299      Time: 0.608        CPU Time: 0.0669452s
Steps: 300      Time: 0.61         CPU Time: 0.0669245s
Steps: 301      Time: 0.612        CPU Time: 0.0671158s
Steps: 302      Time: 0.614        CPU Time: 0.0670081s
Steps: 303      Time: 0.616        CPU Time: 0.0667708s
Steps: 304      Time: 0.618        CPU Time: 0.0669924s
Steps: 305      Time: 0.62         CPU Time: 0.0671367s
Steps: 306      Time: 0.622        CPU Time: 0.0668997s
Steps: 307      Time: 0.624        CPU Time: 0.0671113s
Steps: 308      Time: 0.626        CPU Time: 0.0670429s
Steps: 309      Time: 0.628        CPU Time: 0.0668895s
Steps: 310      Time: 0.63         CPU Time: 0.0669439s
Steps: 311      Time: 0.632        CPU Time: 0.0669533s
Steps: 312      Time: 0.634        CPU Time: 0.0667639s
Steps: 313      Time: 0.636        CPU Time: 0.0670071s
Steps: 314      Time: 0.638        CPU Time: 0.0668667s
Steps: 315      Time: 0.64         CPU Time: 0.0666594s
Steps: 316      Time: 0.642        CPU Time: 0.0666099s
Steps: 317      Time: 0.644        CPU Time: 0.0670197s
Steps: 318      Time: 0.646        CPU Time: 0.0668421s
Steps: 319      Time: 0.648        CPU Time: 0.066964s
Steps: 320      Time: 0.65         CPU Time: 0.0667984s
Steps: 321      Time: 0.652        CPU Time: 0.0665498s
Steps: 322      Time: 0.654        CPU Time: 0.0674205s
Steps: 323      Time: 0.656        CPU Time: 0.0668645s
Steps: 324      Time: 0.658        CPU Time: 0.0665917s
Steps: 325      Time: 0.66         CPU Time: 0.066948s
Steps: 326      Time: 0.662        CPU Time: 0.0667789s
Steps: 327      Time: 0.664        CPU Time: 0.0667243s
Steps: 328      Time: 0.666        CPU Time: 0.0673049s
Steps: 329      Time: 0.668        CPU Time: 0.0666142s
Steps: 330      Time: 0.67         CPU Time: 0.164701s
Steps: 331      Time: 0.672        CPU Time: 0.0654653s
Steps: 332      Time: 0.674        CPU Time: 0.0672326s
Steps: 333      Time: 0.676        CPU Time: 0.0661208s
Steps: 334      Time: 0.678        CPU Time: 0.0659609s
Steps: 335      Time: 0.68         CPU Time: 0.0663785s
Steps: 336      Time: 0.682        CPU Time: 0.0660999s
Steps: 337      Time: 0.684        CPU Time: 0.0664329s
Steps: 338      Time: 0.686        CPU Time: 0.0660343s
Steps: 339      Time: 0.688        CPU Time: 0.0665622s
Steps: 340      Time: 0.69         CPU Time: 0.066812s
Steps: 341      Time: 0.692        CPU Time: 0.0666436s
Steps: 342      Time: 0.694        CPU Time: 0.0662466s
Steps: 343      Time: 0.696        CPU Time: 0.0661274s
Steps: 344      Time: 0.698        CPU Time: 0.0664171s
Steps: 345      Time: 0.7          CPU Time: 0.0661411s
Steps: 346      Time: 0.702        CPU Time: 0.0661484s
Steps: 347      Time: 0.704        CPU Time: 0.066484s
Steps: 348      Time: 0.706        CPU Time: 0.0665721s
Steps: 349      Time: 0.708        CPU Time: 0.0660215s
Steps: 350      Time: 0.71         CPU Time: 0.0662464s
Steps: 351      Time: 0.712        CPU Time: 0.0664143s
Steps: 352      Time: 0.714        CPU Time: 0.066317s
Steps: 353      Time: 0.716        CPU Time: 0.0666861s
Steps: 354      Time: 0.718        CPU Time: 0.0663884s
Steps: 355      Time: 0.72         CPU Time: 0.0662332s
Steps: 356      Time: 0.722        CPU Time: 0.0663939s
Steps: 357      Time: 0.724        CPU Time: 0.0662873s
Steps: 358      Time: 0.726        CPU Time: 0.0659351s
Steps: 359      Time: 0.728        CPU Time: 0.0658187s
Steps: 360      Time: 0.73         CPU Time: 0.0661984s
Steps: 361      Time: 0.732        CPU Time: 0.0660216s
Steps: 362      Time: 0.734        CPU Time: 0.0661623s
Steps: 363      Time: 0.736        CPU Time: 0.0661155s
Steps: 364      Time: 0.738        CPU Time: 0.163249s
Steps: 365      Time: 0.74         CPU Time: 0.0662954s
Steps: 366      Time: 0.742        CPU Time: 0.0666751s
Steps: 367      Time: 0.744        CPU Time: 0.0665671s
Steps: 368      Time: 0.746        CPU Time: 0.0666497s
Steps: 369      Time: 0.748        CPU Time: 0.0670669s
Steps: 370      Time: 0.75         CPU Time: 0.0672601s
Steps: 371      Time: 0.752        CPU Time: 0.066606s
Steps: 372      Time: 0.754        CPU Time: 0.0668398s
Steps: 373      Time: 0.756        CPU Time: 0.0662812s
Steps: 374      Time: 0.758        CPU Time: 0.0667574s
Steps: 375      Time: 0.76         CPU Time: 0.0668228s
Steps: 376      Time: 0.762        CPU Time: 0.0665015s
Steps: 377      Time: 0.764        CPU Time: 0.0668004s
Steps: 378      Time: 0.766        CPU Time: 0.067014s
Steps: 379      Time: 0.768        CPU Time: 0.0665918s
Steps: 380      Time: 0.77         CPU Time: 0.0667517s
Steps: 381      Time: 0.772        CPU Time: 0.0666713s
Steps: 382      Time: 0.774        CPU Time: 0.0664019s
Steps: 383      Time: 0.776        CPU Time: 0.0669231s
Steps: 384      Time: 0.778        CPU Time: 0.0668759s
Steps: 385      Time: 0.78         CPU Time: 0.0667995s
Steps: 386      Time: 0.782        CPU Time: 0.0669999s
Steps: 387      Time: 0.784        CPU Time: 0.0668822s
Steps: 388      Time: 0.786        CPU Time: 0.0664235s
Steps: 389      Time: 0.788        CPU Time: 0.0668345s
Steps: 390      Time: 0.79         CPU Time: 0.0667333s
Steps: 391      Time: 0.792        CPU Time: 0.0665934s
Steps: 392      Time: 0.794        CPU Time: 0.0669674s
Steps: 393      Time: 0.796        CPU Time: 0.0672054s
Steps: 394      Time: 0.798        CPU Time: 0.0668297s
Steps: 395      Time: 0.8          CPU Time: 0.0668547s
Steps: 396      Time: 0.802        CPU Time: 0.0668874s
Steps: 397      Time: 0.804        CPU Time: 0.0662866s
Steps: 398      Time: 0.806        CPU Time: 0.16391s
Steps: 399      Time: 0.808        CPU Time: 0.0660066s
Steps: 400      Time: 0.81         CPU Time: 0.0659615s
Steps: 401      Time: 0.812        CPU Time: 0.0656353s
Steps: 402      Time: 0.814        CPU Time: 0.0658349s
Steps: 403      Time: 0.816        CPU Time: 0.0658084s
Steps: 404      Time: 0.818        CPU Time: 0.0656171s
Steps: 405      Time: 0.82         CPU Time: 0.0656387s
Steps: 406      Time: 0.822        CPU Time: 0.0661054s
Steps: 407      Time: 0.824        CPU Time: 0.0658918s
Steps: 408      Time: 0.826        CPU Time: 0.0655149s
Steps: 409      Time: 0.828        CPU Time: 0.0659655s
Steps: 410      Time: 0.83         CPU Time: 0.0657614s
Steps: 411      Time: 0.832        CPU Time: 0.0659851s
Steps: 412      Time: 0.834        CPU Time: 0.0654777s
Steps: 413      Time: 0.836        CPU Time: 0.0658393s
Steps: 414      Time: 0.838        CPU Time: 0.0657104s
Steps: 415      Time: 0.84         CPU Time: 0.0661591s
Steps: 416      Time: 0.842        CPU Time: 0.0660561s
Steps: 417      Time: 0.844        CPU Time: 0.0657913s
Steps: 418      Time: 0.846        CPU Time: 0.0655357s
Steps: 419      Time: 0.848        CPU Time: 0.0657989s
Steps: 420      Time: 0.85         CPU Time: 0.0655586s
Steps: 421      Time: 0.852        CPU Time: 0.0657657s
Steps: 422      Time: 0.854        CPU Time: 0.0662776s
Steps: 423      Time: 0.856        CPU Time: 0.0660959s
Steps: 424      Time: 0.858        CPU Time: 0.066052s
Steps: 425      Time: 0.86         CPU Time: 0.0660341s
Steps: 426      Time: 0.862        CPU Time: 0.0656546s
Steps: 427      Time: 0.864        CPU Time: 0.0657386s
Steps: 428      Time: 0.866        CPU Time: 0.0658108s
Steps: 429      Time: 0.868        CPU Time: 0.0656987s
Steps: 430      Time: 0.87         CPU Time: 0.0656422s
Steps: 431      Time: 0.872        CPU Time: 0.0657963s
Steps: 432      Time: 0.874        CPU Time: 0.163855s
Steps: 433      Time: 0.876        CPU Time: 0.0955602s
Steps: 434      Time: 0.878        CPU Time: 0.127238s
Steps: 435      Time: 0.88         CPU Time: 0.127473s
Steps: 436      Time: 0.882        CPU Time: 0.127711s
Steps: 437      Time: 0.884        CPU Time: 0.127691s
Steps: 438      Time: 0.886        CPU Time: 0.127521s
Steps: 439      Time: 0.888        CPU Time: 0.127853s
Steps: 440      Time: 0.89         CPU Time: 0.128164s
Steps: 441      Time: 0.892        CPU Time: 0.128086s
Steps: 442      Time: 0.894        CPU Time: 0.127874s
Steps: 443      Time: 0.896        CPU Time: 0.12732s
Steps: 444      Time: 0.898        CPU Time: 0.127739s
Steps: 445      Time: 0.9          CPU Time: 0.127144s
Steps: 446      Time: 0.902        CPU Time: 0.127903s
Steps: 447      Time: 0.904        CPU Time: 0.127374s
Steps: 448      Time: 0.906        CPU Time: 0.128052s
Steps: 449      Time: 0.908        CPU Time: 0.128117s
Steps: 450      Time: 0.91         CPU Time: 0.224356s
Steps: 451      Time: 0.912        CPU Time: 0.125659s
Steps: 452      Time: 0.914        CPU Time: 0.12584s
Steps: 453      Time: 0.916        CPU Time: 0.126731s
Steps: 454      Time: 0.918        CPU Time: 0.126983s
Steps: 455      Time: 0.92         CPU Time: 0.126801s
Steps: 456      Time: 0.922        CPU Time: 0.126223s
Steps: 457      Time: 0.924        CPU Time: 0.126261s
Steps: 458      Time: 0.926        CPU Time: 0.125792s
Steps: 459      Time: 0.928        CPU Time: 0.126295s
Steps: 460      Time: 0.93         CPU Time: 0.125754s
Steps: 461      Time: 0.932        CPU Time: 0.125709s
Steps: 462      Time: 0.934        CPU Time: 0.126585s
Steps: 463      Time: 0.936        CPU Time: 0.126299s
Steps: 464      Time: 0.938        CPU Time: 0.126387s
Steps: 465      Time: 0.94         CPU Time: 0.125506s
Steps: 466      Time: 0.942        CPU Time: 0.125923s
Steps: 467      Time: 0.944        CPU Time: 0.225s  
Steps: 468      Time: 0.946        CPU Time: 0.125712s
Steps: 469      Time: 0.948        CPU Time: 0.125956s
Steps: 470      Time: 0.95         CPU Time: 0.127204s
Steps: 471      Time: 0.952        CPU Time: 0.126998s
Steps: 472      Time: 0.954        CPU Time: 0.127571s
Steps: 473      Time: 0.956        CPU Time: 0.127811s
Steps: 474      Time: 0.958        CPU Time: 0.127427s
Steps: 475      Time: 0.96         CPU Time: 0.127196s
Steps: 476      Time: 0.962        CPU Time: 0.127027s
Steps: 477      Time: 0.964        CPU Time: 0.12724s
Steps: 478      Time: 0.966        CPU Time: 0.127375s
Steps: 479      Time: 0.968        CPU Time: 0.12727s
Steps: 480      Time: 0.97         CPU Time: 0.127467s
Steps: 481      Time: 0.972        CPU Time: 0.127939s
Steps: 482      Time: 0.974        CPU Time: 0.127067s
Steps: 483      Time: 0.976        CPU Time: 0.127906s
Steps: 484      Time: 0.978        CPU Time: 0.223278s
Steps: 485      Time: 0.98         CPU Time: 0.125531s
Steps: 486      Time: 0.982        CPU Time: 0.125711s
Steps: 487      Time: 0.984        CPU Time: 0.127168s
Steps: 488      Time: 0.986        CPU Time: 0.126327s
Steps: 489      Time: 0.988        CPU Time: 0.126166s
Steps: 490      Time: 0.99         CPU Time: 0.126272s
Steps: 491      Time: 0.992        CPU Time: 0.126762s
Steps: 492      Time: 0.994        CPU Time: 0.12628s
Steps: 493      Time: 0.996        CPU Time: 0.127287s
Steps: 494      Time: 0.998        CPU Time: 0.127218s
Steps: 495      Time: 1            CPU Time: 0.126625s
Steps: 496      Time: 1.002        CPU Time: 0.12612s
Steps: 497      Time: 1.004        CPU Time: 0.0961312s
Steps: 498      Time: 1.006        CPU Time: 0.0664836s
Steps: 499      Time: 1.008        CPU Time: 0.0662628s
Steps: 500      Time: 1.01         CPU Time: 0.0660716s
renaming "naca_4.chk" -> "naca_4.bak0.chk"
Writing: "naca_4.chk" (0.0154979s, XML)
Steps: 501      Time: 1.012        CPU Time: 0.0663716s
Steps: 502      Time: 1.014        CPU Time: 0.0665879s
Steps: 503      Time: 1.016        CPU Time: 0.0671289s
Steps: 504      Time: 1.018        CPU Time: 0.164694s
Steps: 505      Time: 1.02         CPU Time: 0.0655611s
Steps: 506      Time: 1.022        CPU Time: 0.065999s
Steps: 507      Time: 1.024        CPU Time: 0.0656405s
Steps: 508      Time: 1.026        CPU Time: 0.0660293s
Steps: 509      Time: 1.028        CPU Time: 0.0658532s
Steps: 510      Time: 1.03         CPU Time: 0.0658245s
Steps: 511      Time: 1.032        CPU Time: 0.0660153s
Steps: 512      Time: 1.034        CPU Time: 0.0660282s
Steps: 513      Time: 1.036        CPU Time: 0.065862s
Steps: 514      Time: 1.038        CPU Time: 0.0658197s
Steps: 515      Time: 1.04         CPU Time: 0.0660893s
Steps: 516      Time: 1.042        CPU Time: 0.0658325s
Steps: 517      Time: 1.044        CPU Time: 0.0657959s
Steps: 518      Time: 1.046        CPU Time: 0.0663729s
Steps: 519      Time: 1.048        CPU Time: 0.0664445s
Steps: 520      Time: 1.05         CPU Time: 0.0660299s
Steps: 521      Time: 1.052        CPU Time: 0.0655556s
Steps: 522      Time: 1.054        CPU Time: 0.0663083s
Steps: 523      Time: 1.056        CPU Time: 0.0661182s
Steps: 524      Time: 1.058        CPU Time: 0.0656946s
Steps: 525      Time: 1.06         CPU Time: 0.065922s
Steps: 526      Time: 1.062        CPU Time: 0.0660314s
Steps: 527      Time: 1.064        CPU Time: 0.066131s
Steps: 528      Time: 1.066        CPU Time: 0.0658351s
Steps: 529      Time: 1.068        CPU Time: 0.0660589s
Steps: 530      Time: 1.07         CPU Time: 0.065694s
Steps: 531      Time: 1.072        CPU Time: 0.0662217s
Steps: 532      Time: 1.074        CPU Time: 0.0660392s
Steps: 533      Time: 1.076        CPU Time: 0.0660537s
Steps: 534      Time: 1.078        CPU Time: 0.0660652s
Steps: 535      Time: 1.08         CPU Time: 0.0661091s
Steps: 536      Time: 1.082        CPU Time: 0.0656077s
Steps: 537      Time: 1.084        CPU Time: 0.0656937s
Steps: 538      Time: 1.086        CPU Time: 0.16402s
Steps: 539      Time: 1.088        CPU Time: 0.0668061s
Steps: 540      Time: 1.09         CPU Time: 0.067345s
Steps: 541      Time: 1.092        CPU Time: 0.0672834s
Steps: 542      Time: 1.094        CPU Time: 0.0671592s
Steps: 543      Time: 1.096        CPU Time: 0.0670947s
Steps: 544      Time: 1.098        CPU Time: 0.0672588s
Steps: 545      Time: 1.1          CPU Time: 0.0670682s
Steps: 546      Time: 1.102        CPU Time: 0.0671423s
Steps: 547      Time: 1.104        CPU Time: 0.0671936s
Steps: 548      Time: 1.106        CPU Time: 0.0668957s
Steps: 549      Time: 1.108        CPU Time: 0.0672187s
Steps: 550      Time: 1.11         CPU Time: 0.0670799s
Steps: 551      Time: 1.112        CPU Time: 0.0670895s
Steps: 552      Time: 1.114        CPU Time: 0.0669561s
Steps: 553      Time: 1.116        CPU Time: 0.0671826s
Steps: 554      Time: 1.118        CPU Time: 0.0670828s
Steps: 555      Time: 1.12         CPU Time: 0.0669988s
Steps: 556      Time: 1.122        CPU Time: 0.0671238s
Steps: 557      Time: 1.124        CPU Time: 0.0671433s
Steps: 558      Time: 1.126        CPU Time: 0.0670889s
Steps: 559      Time: 1.128        CPU Time: 0.0671682s
Steps: 560      Time: 1.13         CPU Time: 0.0669582s
Steps: 561      Time: 1.132        CPU Time: 0.0669144s
Steps: 562      Time: 1.134        CPU Time: 0.0671622s
Steps: 563      Time: 1.136        CPU Time: 0.0666737s
Steps: 564      Time: 1.138        CPU Time: 0.0667329s
Steps: 565      Time: 1.14         CPU Time: 0.0668242s
Steps: 566      Time: 1.142        CPU Time: 0.0667271s
Steps: 567      Time: 1.144        CPU Time: 0.0668105s
Steps: 568      Time: 1.146        CPU Time: 0.0668932s
Steps: 569      Time: 1.148        CPU Time: 0.0665937s
Steps: 570      Time: 1.15         CPU Time: 0.0666797s
Steps: 571      Time: 1.152        CPU Time: 0.0669004s
Steps: 572      Time: 1.154        CPU Time: 0.163916s
Steps: 573      Time: 1.156        CPU Time: 0.0664443s
Steps: 574      Time: 1.158        CPU Time: 0.0663908s
Steps: 575      Time: 1.16         CPU Time: 0.0663449s
Steps: 576      Time: 1.162        CPU Time: 0.0660608s
Steps: 577      Time: 1.164        CPU Time: 0.066124s
Steps: 578      Time: 1.166        CPU Time: 0.0656831s
Steps: 579      Time: 1.168        CPU Time: 0.0658536s
Steps: 580      Time: 1.17         CPU Time: 0.0661835s
Steps: 581      Time: 1.172        CPU Time: 0.0660261s
Steps: 582      Time: 1.174        CPU Time: 0.0659678s
Steps: 583      Time: 1.176        CPU Time: 0.0662136s
Steps: 584      Time: 1.178        CPU Time: 0.0662305s
Steps: 585      Time: 1.18         CPU Time: 0.0662689s
Steps: 586      Time: 1.182        CPU Time: 0.0659769s
Steps: 587      Time: 1.184        CPU Time: 0.0665023s
Steps: 588      Time: 1.186        CPU Time: 0.0662661s
Steps: 589      Time: 1.188        CPU Time: 0.0660097s
Steps: 590      Time: 1.19         CPU Time: 0.0662339s
Steps: 591      Time: 1.192        CPU Time: 0.0660102s
Steps: 592      Time: 1.194        CPU Time: 0.0660935s
Steps: 593      Time: 1.196        CPU Time: 0.0661942s
Steps: 594      Time: 1.198        CPU Time: 0.0659537s
Steps: 595      Time: 1.2          CPU Time: 0.0658511s
Steps: 596      Time: 1.202        CPU Time: 0.0660146s
Steps: 597      Time: 1.204        CPU Time: 0.0661159s
Steps: 598      Time: 1.206        CPU Time: 0.06605s
Steps: 599      Time: 1.208        CPU Time: 0.0662013s
Steps: 600      Time: 1.21         CPU Time: 0.0662296s
Steps: 601      Time: 1.212        CPU Time: 0.0661539s
Steps: 602      Time: 1.214        CPU Time: 0.0661509s
Steps: 603      Time: 1.216        CPU Time: 0.0663135s
Steps: 604      Time: 1.218        CPU Time: 0.0662992s
Steps: 605      Time: 1.22         CPU Time: 0.066249s
Steps: 606      Time: 1.222        CPU Time: 0.164886s
Steps: 607      Time: 1.224        CPU Time: 0.0669721s
Steps: 608      Time: 1.226        CPU Time: 0.0671921s
Steps: 609      Time: 1.228        CPU Time: 0.0673167s
Steps: 610      Time: 1.23         CPU Time: 0.0671249s
Steps: 611      Time: 1.232        CPU Time: 0.0671339s
Steps: 612      Time: 1.234        CPU Time: 0.0673609s
Steps: 613      Time: 1.236        CPU Time: 0.0669526s
Steps: 614      Time: 1.238        CPU Time: 0.0668679s
Steps: 615      Time: 1.24         CPU Time: 0.067167s
Steps: 616      Time: 1.242        CPU Time: 0.0668972s
Steps: 617      Time: 1.244        CPU Time: 0.0669675s
Steps: 618      Time: 1.246        CPU Time: 0.0673075s
Steps: 619      Time: 1.248        CPU Time: 0.067491s
Steps: 620      Time: 1.25         CPU Time: 0.0673544s
Steps: 621      Time: 1.252        CPU Time: 0.0668831s
Steps: 622      Time: 1.254        CPU Time: 0.0667971s
Steps: 623      Time: 1.256        CPU Time: 0.0673603s
Steps: 624      Time: 1.258        CPU Time: 0.0672078s
Steps: 625      Time: 1.26         CPU Time: 0.0672354s
Steps: 626      Time: 1.262        CPU Time: 0.0671908s
Steps: 627      Time: 1.264        CPU Time: 0.0674049s
Steps: 628      Time: 1.266        CPU Time: 0.0670915s
Steps: 629      Time: 1.268        CPU Time: 0.0670478s
Steps: 630      Time: 1.27         CPU Time: 0.0671269s
Steps: 631      Time: 1.272        CPU Time: 0.0670895s
Steps: 632      Time: 1.274        CPU Time: 0.06723s
Steps: 633      Time: 1.276        CPU Time: 0.0673908s
Steps: 634      Time: 1.278        CPU Time: 0.0673855s
Steps: 635      Time: 1.28         CPU Time: 0.0669654s
Steps: 636      Time: 1.282        CPU Time: 0.0670571s
Steps: 637      Time: 1.284        CPU Time: 0.067197s
Steps: 638      Time: 1.286        CPU Time: 0.0671402s
Steps: 639      Time: 1.288        CPU Time: 0.0672785s
Steps: 640      Time: 1.29         CPU Time: 0.16421s
Steps: 641      Time: 1.292        CPU Time: 0.0660765s
Steps: 642      Time: 1.294        CPU Time: 0.0661887s
Steps: 643      Time: 1.296        CPU Time: 0.0659785s
Steps: 644      Time: 1.298        CPU Time: 0.0660946s
Steps: 645      Time: 1.3          CPU Time: 0.0660811s
Steps: 646      Time: 1.302        CPU Time: 0.066113s
Steps: 647      Time: 1.304        CPU Time: 0.0660056s
Steps: 648      Time: 1.306        CPU Time: 0.0660612s
Steps: 649      Time: 1.308        CPU Time: 0.0661887s
Steps: 650      Time: 1.31         CPU Time: 0.0658185s
Steps: 651      Time: 1.312        CPU Time: 0.0662921s
Steps: 652      Time: 1.314        CPU Time: 0.0661955s
Steps: 653      Time: 1.316        CPU Time: 0.0658635s
Steps: 654      Time: 1.318        CPU Time: 0.0658637s
Steps: 655      Time: 1.32         CPU Time: 0.0660102s
Steps: 656      Time: 1.322        CPU Time: 0.0657596s
Steps: 657      Time: 1.324        CPU Time: 0.0655935s
Steps: 658      Time: 1.326        CPU Time: 0.0659933s
Steps: 659      Time: 1.328        CPU Time: 0.0658084s
Steps: 660      Time: 1.33         CPU Time: 0.0658654s
Steps: 661      Time: 1.332        CPU Time: 0.0659728s
Steps: 662      Time: 1.334        CPU Time: 0.0659779s
Steps: 663      Time: 1.336        CPU Time: 0.0657101s
Steps: 664      Time: 1.338        CPU Time: 0.0660382s
Steps: 665      Time: 1.34         CPU Time: 0.0660306s
Steps: 666      Time: 1.342        CPU Time: 0.0658044s
Steps: 667      Time: 1.344        CPU Time: 0.0661791s
Steps: 668      Time: 1.346        CPU Time: 0.065987s
Steps: 669      Time: 1.348        CPU Time: 0.0656214s
Steps: 670      Time: 1.35         CPU Time: 0.0656294s
Steps: 671      Time: 1.352        CPU Time: 0.0666082s
Steps: 672      Time: 1.354        CPU Time: 0.0656694s
Steps: 673      Time: 1.356        CPU Time: 0.0660476s
Steps: 674      Time: 1.358        CPU Time: 0.164407s
Steps: 675      Time: 1.36         CPU Time: 0.0678624s
Steps: 676      Time: 1.362        CPU Time: 0.0688172s
Steps: 677      Time: 1.364        CPU Time: 0.0685224s
Steps: 678      Time: 1.366        CPU Time: 0.0688477s
Steps: 679      Time: 1.368        CPU Time: 0.0686272s
Steps: 680      Time: 1.37         CPU Time: 0.0687987s
Steps: 681      Time: 1.372        CPU Time: 0.0685475s
Steps: 682      Time: 1.374        CPU Time: 0.06859s
Steps: 683      Time: 1.376        CPU Time: 0.0685736s
Steps: 684      Time: 1.378        CPU Time: 0.0686692s
Steps: 685      Time: 1.38         CPU Time: 0.0683725s
Steps: 686      Time: 1.382        CPU Time: 0.0687105s
Steps: 687      Time: 1.384        CPU Time: 0.0688626s
Steps: 688      Time: 1.386        CPU Time: 0.0688812s
Steps: 689      Time: 1.388        CPU Time: 0.068497s
Steps: 690      Time: 1.39         CPU Time: 0.0684464s
Steps: 691      Time: 1.392        CPU Time: 0.0686307s
Steps: 692      Time: 1.394        CPU Time: 0.0685508s
Steps: 693      Time: 1.396        CPU Time: 0.0681682s
Steps: 694      Time: 1.398        CPU Time: 0.0686284s
Steps: 695      Time: 1.4          CPU Time: 0.0687474s
Steps: 696      Time: 1.402        CPU Time: 0.068658s
Steps: 697      Time: 1.404        CPU Time: 0.0683176s
Steps: 698      Time: 1.406        CPU Time: 0.0686589s
Steps: 699      Time: 1.408        CPU Time: 0.0689339s
Steps: 700      Time: 1.41         CPU Time: 0.068528s
Steps: 701      Time: 1.412        CPU Time: 0.0683074s
Steps: 702      Time: 1.414        CPU Time: 0.06885s
Steps: 703      Time: 1.416        CPU Time: 0.0687522s
Steps: 704      Time: 1.418        CPU Time: 0.0695425s
Steps: 705      Time: 1.42         CPU Time: 0.0685287s
Steps: 706      Time: 1.422        CPU Time: 0.0686594s
Steps: 707      Time: 1.424        CPU Time: 0.0680259s
Steps: 708      Time: 1.426        CPU Time: 0.167506s
Steps: 709      Time: 1.428        CPU Time: 0.067164s
Steps: 710      Time: 1.43         CPU Time: 0.0678041s
Steps: 711      Time: 1.432        CPU Time: 0.0677023s
Steps: 712      Time: 1.434        CPU Time: 0.0672114s
Steps: 713      Time: 1.436        CPU Time: 0.067753s
Steps: 714      Time: 1.438        CPU Time: 0.0676932s
Steps: 715      Time: 1.44         CPU Time: 0.0673137s
Steps: 716      Time: 1.442        CPU Time: 0.0676422s
Steps: 717      Time: 1.444        CPU Time: 0.0680092s
Steps: 718      Time: 1.446        CPU Time: 0.0673631s
Steps: 719      Time: 1.448        CPU Time: 0.0674576s
Steps: 720      Time: 1.45         CPU Time: 0.0676187s
Steps: 721      Time: 1.452        CPU Time: 0.0671454s
Steps: 722      Time: 1.454        CPU Time: 0.0675453s
Steps: 723      Time: 1.456        CPU Time: 0.0679321s
Steps: 724      Time: 1.458        CPU Time: 0.0676568s
Steps: 725      Time: 1.46         CPU Time: 0.0671638s
Steps: 726      Time: 1.462        CPU Time: 0.0673744s
Steps: 727      Time: 1.464        CPU Time: 0.0678024s
Steps: 728      Time: 1.466        CPU Time: 0.0676144s
Steps: 729      Time: 1.468        CPU Time: 0.0671391s
Steps: 730      Time: 1.47         CPU Time: 0.0674575s
Steps: 731      Time: 1.472        CPU Time: 0.0680985s
Steps: 732      Time: 1.474        CPU Time: 0.0673991s
Steps: 733      Time: 1.476        CPU Time: 0.0677714s
Steps: 734      Time: 1.478        CPU Time: 0.0677669s
Steps: 735      Time: 1.48         CPU Time: 0.0673055s
Steps: 736      Time: 1.482        CPU Time: 0.0674726s
Steps: 737      Time: 1.484        CPU Time: 0.0674697s
Steps: 738      Time: 1.486        CPU Time: 0.0674096s
Steps: 739      Time: 1.488        CPU Time: 0.0677035s
Steps: 740      Time: 1.49         CPU Time: 0.067646s
Steps: 741      Time: 1.492        CPU Time: 0.0673846s
Steps: 742      Time: 1.494        CPU Time: 0.16731s
Steps: 743      Time: 1.496        CPU Time: 0.0680298s
Steps: 744      Time: 1.498        CPU Time: 0.0683395s
Steps: 745      Time: 1.5          CPU Time: 0.0684183s
Steps: 746      Time: 1.502        CPU Time: 0.0682138s
Steps: 747      Time: 1.504        CPU Time: 0.0684528s
Steps: 748      Time: 1.506        CPU Time: 0.0684352s
Steps: 749      Time: 1.508        CPU Time: 0.0678106s
Steps: 750      Time: 1.51         CPU Time: 0.0685358s
Steps: 751      Time: 1.512        CPU Time: 0.0682356s
Steps: 752      Time: 1.514        CPU Time: 0.0682038s
Steps: 753      Time: 1.516        CPU Time: 0.0684271s
Steps: 754      Time: 1.518        CPU Time: 0.0681748s
Steps: 755      Time: 1.52         CPU Time: 0.0681929s
Steps: 756      Time: 1.522        CPU Time: 0.0686146s
Steps: 757      Time: 1.524        CPU Time: 0.0684096s
Steps: 758      Time: 1.526        CPU Time: 0.0681887s
Steps: 759      Time: 1.528        CPU Time: 0.0686182s
Steps: 760      Time: 1.53         CPU Time: 0.0685426s
Steps: 761      Time: 1.532        CPU Time: 0.06812s
Steps: 762      Time: 1.534        CPU Time: 0.0683676s
Steps: 763      Time: 1.536        CPU Time: 0.0684913s
Steps: 764      Time: 1.538        CPU Time: 0.0681752s
Steps: 765      Time: 1.54         CPU Time: 0.0682261s
Steps: 766      Time: 1.542        CPU Time: 0.0686096s
Steps: 767      Time: 1.544        CPU Time: 0.0681351s
Steps: 768      Time: 1.546        CPU Time: 0.0680296s
Steps: 769      Time: 1.548        CPU Time: 0.0686129s
Steps: 770      Time: 1.55         CPU Time: 0.068212s
Steps: 771      Time: 1.552        CPU Time: 0.0684347s
Steps: 772      Time: 1.554        CPU Time: 0.0685436s
Steps: 773      Time: 1.556        CPU Time: 0.0684096s
Steps: 774      Time: 1.558        CPU Time: 0.0680562s
Steps: 775      Time: 1.56         CPU Time: 0.0685916s
Steps: 776      Time: 1.562        CPU Time: 0.166654s
Steps: 777      Time: 1.564        CPU Time: 0.0672426s
Steps: 778      Time: 1.566        CPU Time: 0.0678238s
Steps: 779      Time: 1.568        CPU Time: 0.0679944s
Steps: 780      Time: 1.57         CPU Time: 0.0674798s
Steps: 781      Time: 1.572        CPU Time: 0.0675413s
Steps: 782      Time: 1.574        CPU Time: 0.067597s
Steps: 783      Time: 1.576        CPU Time: 0.0676954s
Steps: 784      Time: 1.578        CPU Time: 0.0679914s
Steps: 785      Time: 1.58         CPU Time: 0.0678937s
Steps: 786      Time: 1.582        CPU Time: 0.0674724s
Steps: 787      Time: 1.584        CPU Time: 0.0676255s
Steps: 788      Time: 1.586        CPU Time: 0.0678725s
Steps: 789      Time: 1.588        CPU Time: 0.0676785s
Steps: 790      Time: 1.59         CPU Time: 0.067682s
Steps: 791      Time: 1.592        CPU Time: 0.0678758s
Steps: 792      Time: 1.594        CPU Time: 0.068015s
Steps: 793      Time: 1.596        CPU Time: 0.0677898s
Steps: 794      Time: 1.598        CPU Time: 0.0673934s
Steps: 795      Time: 1.6          CPU Time: 0.0679742s
Steps: 796      Time: 1.602        CPU Time: 0.0680615s
Steps: 797      Time: 1.604        CPU Time: 0.067451s
Steps: 798      Time: 1.606        CPU Time: 0.0677489s
Steps: 799      Time: 1.608        CPU Time: 0.0682407s
Steps: 800      Time: 1.61         CPU Time: 0.0670312s
Steps: 801      Time: 1.612        CPU Time: 0.0677704s
Steps: 802      Time: 1.614        CPU Time: 0.0678729s
Steps: 803      Time: 1.616        CPU Time: 0.067503s
Steps: 804      Time: 1.618        CPU Time: 0.0678097s
Steps: 805      Time: 1.62         CPU Time: 0.068055s
Steps: 806      Time: 1.622        CPU Time: 0.0675712s
Steps: 807      Time: 1.624        CPU Time: 0.0676915s
Steps: 808      Time: 1.626        CPU Time: 0.0682019s
Steps: 809      Time: 1.628        CPU Time: 0.0674776s
Steps: 810      Time: 1.63         CPU Time: 0.167831s
Steps: 811      Time: 1.632        CPU Time: 0.0693578s
Steps: 812      Time: 1.634        CPU Time: 0.0684064s
Steps: 813      Time: 1.636        CPU Time: 0.068945s
Steps: 814      Time: 1.638        CPU Time: 0.0687344s
Steps: 815      Time: 1.64         CPU Time: 0.0689657s
Steps: 816      Time: 1.642        CPU Time: 0.0689048s
Steps: 817      Time: 1.644        CPU Time: 0.0691113s
Steps: 818      Time: 1.646        CPU Time: 0.0680139s
Steps: 819      Time: 1.648        CPU Time: 0.0690878s
Steps: 820      Time: 1.65         CPU Time: 0.0681125s
Steps: 821      Time: 1.652        CPU Time: 0.0691723s
Steps: 822      Time: 1.654        CPU Time: 0.0691004s
Steps: 823      Time: 1.656        CPU Time: 0.068786s
Steps: 824      Time: 1.658        CPU Time: 0.0690516s
Steps: 825      Time: 1.66         CPU Time: 0.0996486s
Steps: 826      Time: 1.662        CPU Time: 0.131742s
Steps: 827      Time: 1.664        CPU Time: 0.131897s
Steps: 828      Time: 1.666        CPU Time: 0.131266s
Steps: 829      Time: 1.668        CPU Time: 0.13129s
Steps: 830      Time: 1.67         CPU Time: 0.131582s
Steps: 831      Time: 1.672        CPU Time: 0.131815s
Steps: 832      Time: 1.674        CPU Time: 0.131152s
Steps: 833      Time: 1.676        CPU Time: 0.131523s
Steps: 834      Time: 1.678        CPU Time: 0.131399s
Steps: 835      Time: 1.68         CPU Time: 0.227672s
Steps: 836      Time: 1.682        CPU Time: 0.128868s
Steps: 837      Time: 1.684        CPU Time: 0.129107s
Steps: 838      Time: 1.686        CPU Time: 0.0677898s
Steps: 839      Time: 1.688        CPU Time: 0.0677319s
Steps: 840      Time: 1.69         CPU Time: 0.0673839s
Steps: 841      Time: 1.692        CPU Time: 0.0673508s
Steps: 842      Time: 1.694        CPU Time: 0.0678902s
Steps: 843      Time: 1.696        CPU Time: 0.0675893s
Steps: 844      Time: 1.698        CPU Time: 0.0674098s
Steps: 845      Time: 1.7          CPU Time: 0.0678582s
Steps: 846      Time: 1.702        CPU Time: 0.0677416s
Steps: 847      Time: 1.704        CPU Time: 0.0671056s
Steps: 848      Time: 1.706        CPU Time: 0.0680113s
Steps: 849      Time: 1.708        CPU Time: 0.0676492s
Steps: 850      Time: 1.71         CPU Time: 0.0672749s
Steps: 851      Time: 1.712        CPU Time: 0.0676177s
Steps: 852      Time: 1.714        CPU Time: 0.0677063s
Steps: 853      Time: 1.716        CPU Time: 0.0673854s
Steps: 854      Time: 1.718        CPU Time: 0.0676281s
Steps: 855      Time: 1.72         CPU Time: 0.0676262s
Steps: 856      Time: 1.722        CPU Time: 0.0673155s
Steps: 857      Time: 1.724        CPU Time: 0.0677862s
Steps: 858      Time: 1.726        CPU Time: 0.0676s 
Steps: 859      Time: 1.728        CPU Time: 0.0672531s
Steps: 860      Time: 1.73         CPU Time: 0.0674707s
Steps: 861      Time: 1.732        CPU Time: 0.067538s
Steps: 862      Time: 1.734        CPU Time: 0.0675525s
Steps: 863      Time: 1.736        CPU Time: 0.0671715s
Steps: 864      Time: 1.738        CPU Time: 0.0676323s
Steps: 865      Time: 1.74         CPU Time: 0.0674928s
Steps: 866      Time: 1.742        CPU Time: 0.166616s
Steps: 867      Time: 1.744        CPU Time: 0.0682471s
Steps: 868      Time: 1.746        CPU Time: 0.0685549s
Steps: 869      Time: 1.748        CPU Time: 0.0686765s
Steps: 870      Time: 1.75         CPU Time: 0.0682995s
Steps: 871      Time: 1.752        CPU Time: 0.0685933s
Steps: 872      Time: 1.754        CPU Time: 0.0685311s
Steps: 873      Time: 1.756        CPU Time: 0.0680402s
Steps: 874      Time: 1.758        CPU Time: 0.0685865s
Steps: 875      Time: 1.76         CPU Time: 0.0684824s
Steps: 876      Time: 1.762        CPU Time: 0.0996198s
Steps: 877      Time: 1.764        CPU Time: 0.130556s
Steps: 878      Time: 1.766        CPU Time: 0.130935s
Steps: 879      Time: 1.768        CPU Time: 0.131179s
Steps: 880      Time: 1.77         CPU Time: 0.130465s
Steps: 881      Time: 1.772        CPU Time: 0.131097s
Steps: 882      Time: 1.774        CPU Time: 0.131005s
Steps: 883      Time: 1.776        CPU Time: 0.130917s
Steps: 884      Time: 1.778        CPU Time: 0.130545s
Steps: 885      Time: 1.78         CPU Time: 0.130389s
Steps: 886      Time: 1.782        CPU Time: 0.131147s
Steps: 887      Time: 1.784        CPU Time: 0.130698s
Steps: 888      Time: 1.786        CPU Time: 0.227754s
Steps: 889      Time: 1.788        CPU Time: 0.128832s
Steps: 890      Time: 1.79         CPU Time: 0.128988s
Steps: 891      Time: 1.792        CPU Time: 0.128795s
Steps: 892      Time: 1.794        CPU Time: 0.128638s
Steps: 893      Time: 1.796        CPU Time: 0.128447s
Steps: 894      Time: 1.798        CPU Time: 0.129042s
Steps: 895      Time: 1.8          CPU Time: 0.129106s
Steps: 896      Time: 1.802        CPU Time: 0.128738s
Steps: 897      Time: 1.804        CPU Time: 0.128866s
Steps: 898      Time: 1.806        CPU Time: 0.128616s
Steps: 899      Time: 1.808        CPU Time: 0.128902s
Steps: 900      Time: 1.81         CPU Time: 0.128857s
Steps: 901      Time: 1.812        CPU Time: 0.128644s
Steps: 902      Time: 1.814        CPU Time: 0.129003s
Steps: 903      Time: 1.816        CPU Time: 0.12871s
Steps: 904      Time: 1.818        CPU Time: 0.128568s
Steps: 905      Time: 1.82         CPU Time: 0.230286s
Steps: 906      Time: 1.822        CPU Time: 0.129017s
Steps: 907      Time: 1.824        CPU Time: 0.130353s
Steps: 908      Time: 1.826        CPU Time: 0.130529s
Steps: 909      Time: 1.828        CPU Time: 0.130696s
Steps: 910      Time: 1.83         CPU Time: 0.131122s
Steps: 911      Time: 1.832        CPU Time: 0.130437s
Steps: 912      Time: 1.834        CPU Time: 0.130845s
Steps: 913      Time: 1.836        CPU Time: 0.130337s
Steps: 914      Time: 1.838        CPU Time: 0.130608s
Steps: 915      Time: 1.84         CPU Time: 0.130831s
Steps: 916      Time: 1.842        CPU Time: 0.130421s
Steps: 917      Time: 1.844        CPU Time: 0.131105s
Steps: 918      Time: 1.846        CPU Time: 0.130709s
Steps: 919      Time: 1.848        CPU Time: 0.130742s
Steps: 920      Time: 1.85         CPU Time: 0.0997901s
Steps: 921      Time: 1.852        CPU Time: 0.0687435s
Steps: 922      Time: 1.854        CPU Time: 0.0684174s
Steps: 923      Time: 1.856        CPU Time: 0.166384s
Steps: 924      Time: 1.858        CPU Time: 0.0672753s
Steps: 925      Time: 1.86         CPU Time: 0.0678125s
Steps: 926      Time: 1.862        CPU Time: 0.0680093s
Steps: 927      Time: 1.864        CPU Time: 0.0671753s
Steps: 928      Time: 1.866        CPU Time: 0.0676646s
Steps: 929      Time: 1.868        CPU Time: 0.0678861s
Steps: 930      Time: 1.87         CPU Time: 0.0673003s
Steps: 931      Time: 1.872        CPU Time: 0.0677724s
Steps: 932      Time: 1.874        CPU Time: 0.0676601s
Steps: 933      Time: 1.876        CPU Time: 0.0673572s
Steps: 934      Time: 1.878        CPU Time: 0.0677659s
Steps: 935      Time: 1.88         CPU Time: 0.067506s
Steps: 936      Time: 1.882        CPU Time: 0.0676228s
Steps: 937      Time: 1.884        CPU Time: 0.0677223s
Steps: 938      Time: 1.886        CPU Time: 0.0676277s
Steps: 939      Time: 1.888        CPU Time: 0.067542s
Steps: 940      Time: 1.89         CPU Time: 0.0675123s
Steps: 941      Time: 1.892        CPU Time: 0.0678136s
Steps: 942      Time: 1.894        CPU Time: 0.0678701s
Steps: 943      Time: 1.896        CPU Time: 0.0673936s
Steps: 944      Time: 1.898        CPU Time: 0.0678668s
Steps: 945      Time: 1.9          CPU Time: 0.0676903s
Steps: 946      Time: 1.902        CPU Time: 0.0675852s
Steps: 947      Time: 1.904        CPU Time: 0.0674695s
Steps: 948      Time: 1.906        CPU Time: 0.0675942s
Steps: 949      Time: 1.908        CPU Time: 0.0675231s
Steps: 950      Time: 1.91         CPU Time: 0.0677544s
Steps: 951      Time: 1.912        CPU Time: 0.0676849s
Steps: 952      Time: 1.914        CPU Time: 0.0677765s
Steps: 953      Time: 1.916        CPU Time: 0.0676378s
Steps: 954      Time: 1.918        CPU Time: 0.0677551s
Steps: 955      Time: 1.92         CPU Time: 0.0677076s
Steps: 956      Time: 1.922        CPU Time: 0.0672524s
Steps: 957      Time: 1.924        CPU Time: 0.166283s
Steps: 958      Time: 1.926        CPU Time: 0.0673562s
Steps: 959      Time: 1.928        CPU Time: 0.0685373s
Steps: 960      Time: 1.93         CPU Time: 0.0684724s
Steps: 961      Time: 1.932        CPU Time: 0.0682402s
Steps: 962      Time: 1.934        CPU Time: 0.0688092s
Steps: 963      Time: 1.936        CPU Time: 0.0684279s
Steps: 964      Time: 1.938        CPU Time: 0.0681907s
Steps: 965      Time: 1.94         CPU Time: 0.0687095s
Steps: 966      Time: 1.942        CPU Time: 0.0683071s
Steps: 967      Time: 1.944        CPU Time: 0.0681998s
Steps: 968      Time: 1.946        CPU Time: 0.0688318s
Steps: 969      Time: 1.948        CPU Time: 0.0997656s
Steps: 970      Time: 1.95         CPU Time: 0.130773s
Steps: 971      Time: 1.952        CPU Time: 0.130745s
Steps: 972      Time: 1.954        CPU Time: 0.130407s
Steps: 973      Time: 1.956        CPU Time: 0.130931s
Steps: 974      Time: 1.958        CPU Time: 0.13027s
Steps: 975      Time: 1.96         CPU Time: 0.130992s
Steps: 976      Time: 1.962        CPU Time: 0.130658s
Steps: 977      Time: 1.964        CPU Time: 0.130618s
Steps: 978      Time: 1.966        CPU Time: 0.130926s
Steps: 979      Time: 1.968        CPU Time: 0.130674s
Steps: 980      Time: 1.97         CPU Time: 0.227641s
Steps: 981      Time: 1.972        CPU Time: 0.128976s
Steps: 982      Time: 1.974        CPU Time: 0.129092s
Steps: 983      Time: 1.976        CPU Time: 0.129091s
Steps: 984      Time: 1.978        CPU Time: 0.128763s
Steps: 985      Time: 1.98         CPU Time: 0.129394s
Steps: 986      Time: 1.982        CPU Time: 0.128867s
Steps: 987      Time: 1.984        CPU Time: 0.129365s
Steps: 988      Time: 1.986        CPU Time: 0.128651s
Steps: 989      Time: 1.988        CPU Time: 0.128626s
Steps: 990      Time: 1.99         CPU Time: 0.129501s
Steps: 991      Time: 1.992        CPU Time: 0.128947s
Steps: 992      Time: 1.994        CPU Time: 0.129196s
Steps: 993      Time: 1.996        CPU Time: 0.129379s
Steps: 994      Time: 1.998        CPU Time: 0.129037s
Steps: 995      Time: 2            CPU Time: 0.129178s
Steps: 996      Time: 2.002        CPU Time: 0.128706s
Steps: 997      Time: 2.004        CPU Time: 0.227751s
Steps: 998      Time: 2.006        CPU Time: 0.129278s
Steps: 999      Time: 2.008        CPU Time: 0.130401s
Steps: 1000     Time: 2.01         CPU Time: 0.130978s
renaming "naca_5.chk" -> "naca_5.bak0.chk"
Writing: "naca_5.chk" (0.0159675s, XML)
Time-integration  : 82.5745s
renaming "naca.fld" -> "naca.bak6.fld"
Writing: "naca.fld" (0.0158489s, XML)
-------------------------------------------
Total Computation Time = 83s
-------------------------------------------
L 2 error (variable rho) : 12.2441
L inf error (variable rho) : 1.08474
L 2 error (variable rhou) : 12.2431
L inf error (variable rhou) : 1.40795
L 2 error (variable rhov) : 0.0787586
L inf error (variable rhov) : 0.617043
L 2 error (variable E) : 142.775
L inf error (variable E) : 12.6824

Post-processing

To visualise the flowfield first we convert the .fld file (generated by the solver) to a .vtu format file by using the FieldConvert ulitity available in the Nektar++.

Task 18: Run the following code cell to convet the .fld into a .vtu file.

!FieldConvert -f naca.xml session_naca_im.xml naca.fld flowfield_im.vtu
Command output
Writing: "flowfield_im.vtu"
Written file: flowfield_im.vtu
import pyvista as pv

# First read the VTK file
mesh = pv.read('flowfield_im.vtu')

# Now create an itkwidget to visualise this in-browser.
pl = pv.PlotterITK()
pl.add_mesh(mesh, scalars=mesh['u'], smooth_shading=True)
pl.show(True)
Viewer(geometries=[{'vtkClass': 'vtkPolyData', 'points': {'vtkClass': 'vtkPoints', 'name': '_points', 'numberO…
Viewer(geometries=[{'vtkClass': 'vtkPolyData', 'points': {'vtkClass': 'vtkPoints', 'name': '_points', 'numberO…

Extra tasks

Task 19: A better resolution can be obtained by increasing the expansion order under the <EXPANSIONS> tag. We can set the NUMMODES=5 and check the result. Please noted that a higher order expansion may require the time step to be smaller to avoid divergence.


**Note:** To run the Navier-Stokes simulation, under the `` tag change `EulerCFE` or `EulerImplicitCFE` to be `NavierStokesCFE` or `NavierStokesImplicitCFE` for an explicit or implicit simulation, respectively . `DiffusionType`, `DiffusionAdvancement`, and `ViscosityType` need to be specified as well. Other differences include the parameters related to viscosity and heat transfer under `` tag and the type of wall boundary condition under `` tag. An example for implicitly solving Navier-Stokes equations is provided in [**session_naca_NS-Complete.xml**](session_naca_NS-Complete.xml). More information can also be found in the [tutorial](https://doc.nektar.info/tutorials/latest/cfs/CylinderSubsonic_NS/cfs-CylinderSubsonic_NS.pdf) </div> # Take-a-ways It is worthy to keep in mind the following points: - The comprssible flows can be solved either explicitly or implicitly; - The explicit solving requires a smaller time step to satisfy the CFL condition while the restriction is very much relexed for the implicit solving; # References 1. Nektar++: Spectral/hp Element Framework, Version 5.0.1, [User guide](http://www.nektar.info/downloads), January, 2021. 2. Kirby R, Sherwin S. Stabilisation of spectral/hp methods through spectral vanishing viscosity: Application to fluid mechanics modelling. Comput. Methods Appl. Mech. Engrg. 2006; 195: 3128–3144 ## End of the tutorial. </div> </div> </div>