HyPar  1.0
Finite-Difference Hyperbolic-Parabolic PDE Solver on Cartesian Grids
2D Navier-Stokes Equations - Laminar Flow over Flat Plate

Location: hypar/Examples/2D/NavierStokes2D/FlatPlateLaminar_PETSc_Implicit (This directory contains all the input files needed to run this case. If there is a Run.m, run it in MATLAB to quickly set up, run, and visualize the example).

Governing equations: 2D Euler Equations (navierstokes2d.h)

Reference:

  • Hirsch, "Numerical Computation of Internal & External Flows", Volume 1 (Fundamentals of Computational Fluid Dynamics), 2nd Edition, Elsevier, Section 12.3.2 (page 618-625).

Domain: \(-0.25 \le x \le 1\), \(0 \le y \le 0.25\)

Boundary conditions:

  • Symmetry BC on \(y=0, -0.25 \le x < 0\) (imposed through "slip-wall" _SLIP_WALL_ with 0 wall velocity).
  • No-slip wall BC on \(y=0, 0 \le x \le 1\) (_NOSLIP_WALL_ with 0 wall velocity).
  • Subsonic outflow on \(y=0.25\) (_SUBSONIC_OUTFLOW_) with pressure \(p=1/\gamma\).
  • Subsonic inflow on \(x=0\) (_SUBSONIC_INFLOW_) with density \(\rho=1\), and velocity \((u,v) = (0.3,0)\).
  • Subsonic outflow on \(x=1\) (_SUBSONIC_OUTFLOW_) with pressure \(p=1/\gamma\).

Initial solution: Uniform flow with \(\rho=1, u=0.3, v=0, p=1/\gamma\).

Other parameters:

Note: Pressure is taken as \(1/\gamma\) in the above so that the freestream speed of sound is 1.

Numerical method:

This is a steady state problem - the solution converges to a steady laminar flow over a flat plate, and the residuals decrease with time. The skin friction coefficient is given by

\begin{equation} c_f = \frac {0.664} {\sqrt{Re_x}}, Re_x = \frac {\rho u x} {\mu} \end{equation}

(Blasius solution).

Input files required:

.petscrc

# See PETSc documentation for more details (https://petsc.org/release/overview/).
# Note that if the following are specified in this file, the corresponding inputs in solver.inp are *ignored*.
# + "-ts_dt" (time step size): ignores "dt" in solver.inp
# + "-ts_max_steps" (maximum number of time iterations): ignores "n_iter" in solver.inp
# + "-ts_max_time" (final simulation time): ignores "n_iter" X "dt" in solver.inp
# If these are not specified here, then the values in solver.inp are used.
# Use PETSc time-integration
-use-petscts
# Time integration scheme type - ARK
-ts_type arkimex
-ts_arkimex_type 2e
-ts_arkimex_fully_implicit
# no time-step adaptivity
-ts_adapt_type none
# Nonlinear solver (SNES) type
-snes_type newtonls
# Set relative tolerance
-snes_rtol 1e-4
# Set absolute tolerance
-snes_atol 1e-4
# Set step size tolerance
-snes_stol 1e-16
# Linear solver (KSP) type
-ksp_type gmres
# Set relative tolerance
-ksp_rtol 1e-4
# Set absolute tolerance
-ksp_atol 1e-4
# use a preconditioner for solving the system
-with_pc
# preconditioner type - SOR
-pc_type sor
-pc_sor_omega 1.0
-pc_sor_its 1
# apply right preconditioner
-ksp_pc_side RIGHT

solver.inp

begin
ndims 2
nvars 4
size 226 502
iproc 2 3
ghost 3
n_iter 400
restart_iter 0
time_scheme rk
time_scheme_type 44
hyp_space_scheme weno5
hyp_interp_type components
par_space_type nonconservative-2stage
par_space_scheme 4
dt 0.05
screen_op_iter 1
file_op_iter 10
ip_file_type binary
op_file_format binary
op_overwrite yes
model navierstokes2d
end

boundary.inp

5
slip-wall 1 1 -0.25 0.00 0 0
0.0 0.0
noslip-wall 1 1 0.00 2.00 0 0
0.0 0.0
subsonic-outflow 1 -1 -0.25 2.00 0 0
0.7142857143
subsonic-inflow 0 1 0 0 0.00 0.25
1.0 0.3 0.0
subsonic-outflow 0 -1 0 0 0.00 0.25
0.7142857143

physics.inp

begin
gamma 1.4
upwinding roe
Pr 0.72
Minf 0.3
Re 100000
end

weno.inp (optional)

begin
mapped 0
borges 0
yc 0
no_limiting 1
epsilon 0.000001
p 2.0
rc 0.3
xi 0.001
end

To generate initial.inp (initial solution), compile and run the following code in the run directory.

#include <stdio.h>
#include <stdlib.h>
#include <math.h>
#include <string.h>
double GAMMA = 1.4;
int main()
{
int NI,NJ,ndims;
char ip_file_type[50];
strcpy(ip_file_type,"ascii");
FILE *in;
printf("Reading file \"solver.inp\"...\n");
in= fopen("solver.inp","r");
if (!in) printf("Error: Input file \"solver.inp\" not found. Default values will be used.\n");
else {
char word[500];
fscanf(in,"%s",word);
if (!strcmp(word, "begin")){
while (strcmp(word, "end")){
fscanf(in,"%s",word);
if (!strcmp(word, "ndims")) fscanf(in,"%d",&ndims);
else if (!strcmp(word, "size")) {
fscanf(in,"%d",&NI);
fscanf(in,"%d",&NJ);
} else if (!strcmp(word, "ip_file_type")) fscanf(in,"%s",ip_file_type);
}
} else printf("Error: Illegal format in solver.inp. Crash and burn!\n");
}
fclose(in);
if (ndims != 2) {
printf("ndims is not 2 in solver.inp. this code is to generate 2D solution\n");
return(0);
}
printf("Grid:\t\t\t%d x %d\n",NI,NJ);
double xmin = -0.25;
double xmax = 1.0;
double ymin = 0.0;
double ymax = 0.25;
double Lx = xmax - xmin;
double Ly = ymax - ymin;
int i,j;
double dx = Lx / ((double)(NI-1));
double dy = Ly / ((double)(NJ-1));
double *x, *y, *U;
FILE *out;
x = (double*) calloc (NI , sizeof(double));
y = (double*) calloc (NJ , sizeof(double));
U = (double*) calloc (4*NI*NJ, sizeof(double));
for (i = 0; i < NI; i++){
for (j = 0; j < NJ; j++){
x[i] = xmin + i*dx;
y[j] = ymin + j*dy;
int p = i + NI*j;
double rho, u, v, P;
rho = 1.0;
P = 1.0/1.4;
u = 0.3;
v = 0.0;
U[4*p+0] = rho;
U[4*p+1] = rho*u;
U[4*p+2] = rho*v;
U[4*p+3] = P/(GAMMA-1.0) + 0.5*rho*(u*u+v*v);
}
}
if (!strcmp(ip_file_type,"ascii")) {
printf("Writing ASCII initial solution file initial.inp\n");
out = fopen("initial.inp","w");
for (i = 0; i < NI; i++) fprintf(out,"%1.16E ",x[i]);
fprintf(out,"\n");
for (j = 0; j < NJ; j++) fprintf(out,"%1.16E ",y[j]);
fprintf(out,"\n");
for (j = 0; j < NJ; j++) {
for (i = 0; i < NI; i++) {
int p = i + NI*j;
fprintf(out,"%1.16E ",U[4*p+0]);
}
}
fprintf(out,"\n");
for (j = 0; j < NJ; j++) {
for (i = 0; i < NI; i++) {
int p = i + NI*j;
fprintf(out,"%1.16E ",U[4*p+1]);
}
}
fprintf(out,"\n");
for (j = 0; j < NJ; j++) {
for (i = 0; i < NI; i++) {
int p = i + NI*j;
fprintf(out,"%1.16E ",U[4*p+2]);
}
}
fprintf(out,"\n");
for (j = 0; j < NJ; j++) {
for (i = 0; i < NI; i++) {
int p = i + NI*j;
fprintf(out,"%1.16E ",U[4*p+3]);
}
}
fprintf(out,"\n");
fclose(out);
} else if ((!strcmp(ip_file_type,"binary")) || (!strcmp(ip_file_type,"bin"))) {
printf("Writing binary initial solution file initial.inp\n");
out = fopen("initial.inp","wb");
fwrite(x,sizeof(double),NI,out);
fwrite(y,sizeof(double),NJ,out);
fwrite(U,sizeof(double),4*NI*NJ,out);
fclose(out);
}
free(x);
free(y);
free(U);
return(0);
}

Output:

Note that iproc is set to

  2 3

in solver.inp (i.e., 2 processor along x, and 3 processor along y). Thus, this example should be run with 6 MPI ranks (or change iproc).

After running the code, there should be one output file op.bin, since HyPar::op_overwrite is set to yes in solver.inp.

HyPar::op_file_format is set to binary in solver.inp, and thus, all the files are written out in the binary format, see WriteBinary(). The binary file contains the conserved variables \(\left(\rho, \rho u, \rho v, e\right)\). The following two codes are available to convert the binary output file:

  • hypar/Extras/BinaryToTecplot.c - convert binary output file to Tecplot file (works only for 2D and 3D).
  • hypar/Extras/BinaryToText.c - convert binary output file to an ASCII text file (to visualize in, for example, MATLAB).

The following plot was obtained by converting the binary file to the Tecplot format, and using VisIt to plot it (it shows the density and velocity):

Solution_2DNavStokFlatPlate.png

The following plot shows a magnified view of the boundary layer:

Solution_2DNavStokFlatPlateMagnified.png

The following file computes the skin friction as a function of the Reynolds number and writes to to a text file SkinFriction.dat with 4 columns: Reynolds number, computed skin friction coefficient, exact skin friction coefficient ( \(0.664/\sqrt{Re_x}\)), and normal velocity gradient on the plate surface ( \(\left.\partial u/\partial y\right|_{y=0}\)). Compile and run it in the run directory.

/*
This code calculates the skin friction
as a function of Reynolds number for
the flow over a flat plate.
*/
#include <stdio.h>
#include <stdlib.h>
#include <string.h>
#include <math.h>
int main()
{
FILE *out, *in, *inputs;
char filename[50], op_file_format[50];
double Re_inf, M_inf;
inputs = fopen("solver.inp","r");
if (!inputs) {
fprintf(stderr,"Error: File \"solver.inp\" not found.\n");
return(1);
} else {
char word[100];
fscanf(inputs,"%s",word);
if (!strcmp(word, "begin")){
while (strcmp(word, "end")){
fscanf(inputs,"%s",word);
if (!strcmp(word, "op_file_format")) fscanf(inputs,"%s" ,op_file_format);
}
}
fclose(inputs);
}
if (strcmp(op_file_format,"binary") && strcmp(op_file_format,"bin")) {
printf("Error: solution output needs to be in binary files.\n");
return(0);
}
inputs = fopen("physics.inp","r");
if (!inputs) {
fprintf(stderr,"Error: File \"physics.inp\" not found.\n");
return(1);
} else {
char word[100];
fscanf(inputs,"%s",word);
if (!strcmp(word, "begin")){
while (strcmp(word, "end")){
fscanf(inputs,"%s",word);
if (!strcmp(word, "Re")) fscanf(inputs,"%lf",&Re_inf);
else if (!strcmp(word, "Minf")) fscanf(inputs,"%lf",&M_inf);
}
}
fclose(inputs);
}
strcpy(filename,"op.bin");
in = fopen(filename,"rb");
out = fopen("SkinFriction.dat","w");
if (!in) {
printf("File op.bin found. Exiting.\n");
return(0);
} else {
printf("Reading file %s.\n",filename);
int ndims, nvars, dims[3];
double *U,*x,*y,*z;
fread(&ndims,sizeof(int),1,in);
fread(&nvars,sizeof(int),1,in);
if (ndims != 2) {
printf("Error: ndims in %s not equal to 2!\n",filename);
return(0);
}
if (nvars != 4) {
printf("Error: nvars in %s not equal to 4!\n",filename);
return(0);
}
fread(dims,sizeof(int),ndims,in);
printf("Dimensions: %d x %d\n",dims[0],dims[1]);
int NI = dims[0];
int NJ = dims[1];
int sizeu = NI*NJ*nvars;
x = (double*) calloc (NI ,sizeof(double));
y = (double*) calloc (NJ ,sizeof(double));
U = (double*) calloc (sizeu,sizeof(double));
fread(x,sizeof(double),NI ,in);
fread(y,sizeof(double),NJ ,in);
fread(U,sizeof(double),sizeu,in);
int i, j, q;
j = 0;
for (i = 1; i < NI-1; i++) {
int q1,q2,q3;
double rho1, u1, v1, e1, p1;
double rho2, u2, v2, e2, p2;
double rho3, u3, v3, e3, p3;
/* calculate du_dy at the wall */
q1 = i + NI*j;
rho1 = U[nvars*q1+0];
u1 = U[nvars*q1+1] / rho1;
v1 = U[nvars*q1+2] / rho1;
e1 = U[nvars*q1+3];
p1 = 0.4 * (e1 - 0.5*rho1*(u1*u1+v1*v1));
q2 = i + NI*(j+1);
rho2 = U[nvars*q2+0];
u2 = U[nvars*q2+1] / rho2;
v2 = U[nvars*q2+2] / rho2;
e2 = U[nvars*q2+3];
p2 = 0.4 * (e2 - 0.5*rho2*(u2*u2+v2*v2));
q3 = i + NI*(j+2);
rho3 = U[nvars*q3+0];
u3 = U[nvars*q3+1] / rho3;
v3 = U[nvars*q3+2] / rho3;
e3 = U[nvars*q3+3];
p3 = 0.4 * (e3 - 0.5*rho3*(u3*u3+v3*v3));
double u0 = -u1;
double u_wall = 0.5*u0 + 0.5*u1;
double u_flow1 = 0.5*u1 + 0.5*u2;
double u_flow2 = 0.5*u2 + 0.5*u3;
double dy = y[j+1] - y[j];
double du_dy = (-3.0*u_wall + 4.0*u_flow1 - u_flow2) / (2.0*dy);
/* Calculate Re_x */
double Re_x = Re_inf * x[i];
/* Calculate exact solution */
double exact_cf = 0.664 / sqrt(Re_x);
/* Calculate the computed cf */
double cf = (M_inf/Re_inf) * (2.0/(M_inf*M_inf)) * du_dy;
if (x[i] > 0) fprintf(out,"%1.16E %1.16E %1.16E %1.16E\n",Re_x,cf,exact_cf,du_dy);
}
free(U);
free(x);
free(y);
fclose(in);
fclose(out);
}
return(0);
}

The following figure showing the exact and computed skin friction coefficients was obtained by plotting SkinFriction.dat:

Solution_2DNavStokFlatPlateSkinFrictionPETSc.png

The file function_counts.dat reports the computational expense (in terms of the number of function counts):

400
22324
22324
22324
2331
2331
1131
19993

The numbers are, respectively,

Expected screen output:

HyPar - Parallel (MPI) version with 6 processes
Compiled with PETSc time integration.
Reading solver inputs from file "solver.inp".
No. of dimensions : 2
No. of variables : 4
Domain size : 226 502
Processes along each dimension : 2 3
No. of ghosts pts : 3
No. of iter. : 400
Restart iteration : 0
Time integration scheme : PETSc
Spatial discretization scheme (hyperbolic) : weno5
Split hyperbolic flux term? : no
Interpolation type for hyperbolic term : components
Spatial discretization type (parabolic ) : nonconservative-2stage
Spatial discretization scheme (parabolic ) : 4
Time Step : 5.000000E-02
Check for conservation : no
Screen output iterations : 1
File output iterations : 10
Initial solution file type : binary
Initial solution read mode : serial
Solution file write mode : serial
Solution file format : binary
Overwrite solution file : yes
Physical model : navierstokes2d
Partitioning domain.
Allocating data arrays.
Reading array from binary file initial.inp (Serial mode).
Volume integral of the initial solution:
0: 3.1451541361714314E-01
1: 9.4354624085191963E-02
2: 0.0000000000000000E+00
3: 5.7578786078649591E-01
Reading boundary conditions from "boundary.inp".
Boundary slip-wall: Along dimension 1 and face +1
Boundary noslip-wall: Along dimension 1 and face +1
Boundary subsonic-outflow: Along dimension 1 and face -1
Boundary subsonic-inflow: Along dimension 0 and face +1
Boundary subsonic-outflow: Along dimension 0 and face -1
5 boundary condition(s) read.
Initializing solvers.
Reading WENO parameters from weno.inp.
Initializing physics. Model = "navierstokes2d"
Reading physical model inputs from file "physics.inp".
Setting up PETSc time integration...
Implicit-Explicit time-integration:-
Hyperbolic term: Implicit
Parabolic term: Implicit
Source term: Implicit
SolvePETSc(): Problem type is nonlinear.
** Starting PETSc time integration **
Writing solution file op.bin.
Iteration: 1 Time: 5.000E-02 Max CFL: 1.011E+02 Max Diff. No.: -1.000E+00
Iteration: 2 Time: 1.000E-01 Max CFL: 1.013E+02 Max Diff. No.: -1.000E+00
Iteration: 3 Time: 1.500E-01 Max CFL: 1.014E+02 Max Diff. No.: -1.000E+00
Iteration: 4 Time: 2.000E-01 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 5 Time: 2.500E-01 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 6 Time: 3.000E-01 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 7 Time: 3.500E-01 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 8 Time: 4.000E-01 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 9 Time: 4.500E-01 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 10 Time: 5.000E-01 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Writing solution file op.bin.
Iteration: 11 Time: 5.500E-01 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 12 Time: 6.000E-01 Max CFL: 1.014E+02 Max Diff. No.: -1.000E+00
Iteration: 13 Time: 6.500E-01 Max CFL: 1.014E+02 Max Diff. No.: -1.000E+00
Iteration: 14 Time: 7.000E-01 Max CFL: 1.014E+02 Max Diff. No.: -1.000E+00
Iteration: 15 Time: 7.500E-01 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 16 Time: 8.000E-01 Max CFL: 1.014E+02 Max Diff. No.: -1.000E+00
Iteration: 17 Time: 8.500E-01 Max CFL: 1.014E+02 Max Diff. No.: -1.000E+00
Iteration: 18 Time: 9.000E-01 Max CFL: 1.014E+02 Max Diff. No.: -1.000E+00
Iteration: 19 Time: 9.500E-01 Max CFL: 1.014E+02 Max Diff. No.: -1.000E+00
Iteration: 20 Time: 1.000E+00 Max CFL: 1.014E+02 Max Diff. No.: -1.000E+00
Writing solution file op.bin.
Iteration: 21 Time: 1.050E+00 Max CFL: 1.014E+02 Max Diff. No.: -1.000E+00
Iteration: 22 Time: 1.100E+00 Max CFL: 1.014E+02 Max Diff. No.: -1.000E+00
Iteration: 23 Time: 1.150E+00 Max CFL: 1.014E+02 Max Diff. No.: -1.000E+00
Iteration: 24 Time: 1.200E+00 Max CFL: 1.014E+02 Max Diff. No.: -1.000E+00
Iteration: 25 Time: 1.250E+00 Max CFL: 1.014E+02 Max Diff. No.: -1.000E+00
Iteration: 26 Time: 1.300E+00 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 27 Time: 1.350E+00 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 28 Time: 1.400E+00 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 29 Time: 1.450E+00 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 30 Time: 1.500E+00 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Writing solution file op.bin.
Iteration: 31 Time: 1.550E+00 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 32 Time: 1.600E+00 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 33 Time: 1.650E+00 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 34 Time: 1.700E+00 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 35 Time: 1.750E+00 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 36 Time: 1.800E+00 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 37 Time: 1.850E+00 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 38 Time: 1.900E+00 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 39 Time: 1.950E+00 Max CFL: 1.014E+02 Max Diff. No.: -1.000E+00
Iteration: 40 Time: 2.000E+00 Max CFL: 1.014E+02 Max Diff. No.: -1.000E+00
Writing solution file op.bin.
Iteration: 41 Time: 2.050E+00 Max CFL: 1.014E+02 Max Diff. No.: -1.000E+00
Iteration: 42 Time: 2.100E+00 Max CFL: 1.014E+02 Max Diff. No.: -1.000E+00
Iteration: 43 Time: 2.150E+00 Max CFL: 1.014E+02 Max Diff. No.: -1.000E+00
Iteration: 44 Time: 2.200E+00 Max CFL: 1.014E+02 Max Diff. No.: -1.000E+00
Iteration: 45 Time: 2.250E+00 Max CFL: 1.014E+02 Max Diff. No.: -1.000E+00
Iteration: 46 Time: 2.300E+00 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 47 Time: 2.350E+00 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 48 Time: 2.400E+00 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 49 Time: 2.450E+00 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 50 Time: 2.500E+00 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Writing solution file op.bin.
Iteration: 51 Time: 2.550E+00 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 52 Time: 2.600E+00 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 53 Time: 2.650E+00 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 54 Time: 2.700E+00 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 55 Time: 2.750E+00 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 56 Time: 2.800E+00 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 57 Time: 2.850E+00 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 58 Time: 2.900E+00 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 59 Time: 2.950E+00 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 60 Time: 3.000E+00 Max CFL: 1.014E+02 Max Diff. No.: -1.000E+00
Writing solution file op.bin.
Iteration: 61 Time: 3.050E+00 Max CFL: 1.014E+02 Max Diff. No.: -1.000E+00
Iteration: 62 Time: 3.100E+00 Max CFL: 1.014E+02 Max Diff. No.: -1.000E+00
Iteration: 63 Time: 3.150E+00 Max CFL: 1.014E+02 Max Diff. No.: -1.000E+00
Iteration: 64 Time: 3.200E+00 Max CFL: 1.014E+02 Max Diff. No.: -1.000E+00
Iteration: 65 Time: 3.250E+00 Max CFL: 1.014E+02 Max Diff. No.: -1.000E+00
Iteration: 66 Time: 3.300E+00 Max CFL: 1.014E+02 Max Diff. No.: -1.000E+00
Iteration: 67 Time: 3.350E+00 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 68 Time: 3.400E+00 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 69 Time: 3.450E+00 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 70 Time: 3.500E+00 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Writing solution file op.bin.
Iteration: 71 Time: 3.550E+00 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 72 Time: 3.600E+00 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 73 Time: 3.650E+00 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 74 Time: 3.700E+00 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 75 Time: 3.750E+00 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 76 Time: 3.800E+00 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 77 Time: 3.850E+00 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 78 Time: 3.900E+00 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 79 Time: 3.950E+00 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 80 Time: 4.000E+00 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Writing solution file op.bin.
Iteration: 81 Time: 4.050E+00 Max CFL: 1.014E+02 Max Diff. No.: -1.000E+00
Iteration: 82 Time: 4.100E+00 Max CFL: 1.014E+02 Max Diff. No.: -1.000E+00
Iteration: 83 Time: 4.150E+00 Max CFL: 1.014E+02 Max Diff. No.: -1.000E+00
Iteration: 84 Time: 4.200E+00 Max CFL: 1.014E+02 Max Diff. No.: -1.000E+00
Iteration: 85 Time: 4.250E+00 Max CFL: 1.014E+02 Max Diff. No.: -1.000E+00
Iteration: 86 Time: 4.300E+00 Max CFL: 1.014E+02 Max Diff. No.: -1.000E+00
Iteration: 87 Time: 4.350E+00 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 88 Time: 4.400E+00 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 89 Time: 4.450E+00 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 90 Time: 4.500E+00 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Writing solution file op.bin.
Iteration: 91 Time: 4.550E+00 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 92 Time: 4.600E+00 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 93 Time: 4.650E+00 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 94 Time: 4.700E+00 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 95 Time: 4.750E+00 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 96 Time: 4.800E+00 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 97 Time: 4.850E+00 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 98 Time: 4.900E+00 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 99 Time: 4.950E+00 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 100 Time: 5.000E+00 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Writing solution file op.bin.
Iteration: 101 Time: 5.050E+00 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 102 Time: 5.100E+00 Max CFL: 1.014E+02 Max Diff. No.: -1.000E+00
Iteration: 103 Time: 5.150E+00 Max CFL: 1.014E+02 Max Diff. No.: -1.000E+00
Iteration: 104 Time: 5.200E+00 Max CFL: 1.014E+02 Max Diff. No.: -1.000E+00
Iteration: 105 Time: 5.250E+00 Max CFL: 1.014E+02 Max Diff. No.: -1.000E+00
Iteration: 106 Time: 5.300E+00 Max CFL: 1.014E+02 Max Diff. No.: -1.000E+00
Iteration: 107 Time: 5.350E+00 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 108 Time: 5.400E+00 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 109 Time: 5.450E+00 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 110 Time: 5.500E+00 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Writing solution file op.bin.
Iteration: 111 Time: 5.550E+00 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 112 Time: 5.600E+00 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 113 Time: 5.650E+00 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 114 Time: 5.700E+00 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 115 Time: 5.750E+00 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 116 Time: 5.800E+00 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 117 Time: 5.850E+00 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 118 Time: 5.900E+00 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 119 Time: 5.950E+00 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 120 Time: 6.000E+00 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Writing solution file op.bin.
Iteration: 121 Time: 6.050E+00 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 122 Time: 6.100E+00 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 123 Time: 6.150E+00 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 124 Time: 6.200E+00 Max CFL: 1.014E+02 Max Diff. No.: -1.000E+00
Iteration: 125 Time: 6.250E+00 Max CFL: 1.014E+02 Max Diff. No.: -1.000E+00
Iteration: 126 Time: 6.300E+00 Max CFL: 1.014E+02 Max Diff. No.: -1.000E+00
Iteration: 127 Time: 6.350E+00 Max CFL: 1.014E+02 Max Diff. No.: -1.000E+00
Iteration: 128 Time: 6.400E+00 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 129 Time: 6.450E+00 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 130 Time: 6.500E+00 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Writing solution file op.bin.
Iteration: 131 Time: 6.550E+00 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 132 Time: 6.600E+00 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 133 Time: 6.650E+00 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 134 Time: 6.700E+00 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 135 Time: 6.750E+00 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 136 Time: 6.800E+00 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 137 Time: 6.850E+00 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 138 Time: 6.900E+00 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 139 Time: 6.950E+00 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 140 Time: 7.000E+00 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Writing solution file op.bin.
Iteration: 141 Time: 7.050E+00 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 142 Time: 7.100E+00 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 143 Time: 7.150E+00 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 144 Time: 7.200E+00 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 145 Time: 7.250E+00 Max CFL: 1.014E+02 Max Diff. No.: -1.000E+00
Iteration: 146 Time: 7.300E+00 Max CFL: 1.014E+02 Max Diff. No.: -1.000E+00
Iteration: 147 Time: 7.350E+00 Max CFL: 1.014E+02 Max Diff. No.: -1.000E+00
Iteration: 148 Time: 7.400E+00 Max CFL: 1.014E+02 Max Diff. No.: -1.000E+00
Iteration: 149 Time: 7.450E+00 Max CFL: 1.014E+02 Max Diff. No.: -1.000E+00
Iteration: 150 Time: 7.500E+00 Max CFL: 1.014E+02 Max Diff. No.: -1.000E+00
Writing solution file op.bin.
Iteration: 151 Time: 7.550E+00 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 152 Time: 7.600E+00 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 153 Time: 7.650E+00 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 154 Time: 7.700E+00 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 155 Time: 7.750E+00 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 156 Time: 7.800E+00 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 157 Time: 7.850E+00 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 158 Time: 7.900E+00 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 159 Time: 7.950E+00 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 160 Time: 8.000E+00 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Writing solution file op.bin.
Iteration: 161 Time: 8.050E+00 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 162 Time: 8.100E+00 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 163 Time: 8.150E+00 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 164 Time: 8.200E+00 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 165 Time: 8.250E+00 Max CFL: 1.014E+02 Max Diff. No.: -1.000E+00
Iteration: 166 Time: 8.300E+00 Max CFL: 1.014E+02 Max Diff. No.: -1.000E+00
Iteration: 167 Time: 8.350E+00 Max CFL: 1.014E+02 Max Diff. No.: -1.000E+00
Iteration: 168 Time: 8.400E+00 Max CFL: 1.014E+02 Max Diff. No.: -1.000E+00
Iteration: 169 Time: 8.450E+00 Max CFL: 1.014E+02 Max Diff. No.: -1.000E+00
Iteration: 170 Time: 8.500E+00 Max CFL: 1.014E+02 Max Diff. No.: -1.000E+00
Writing solution file op.bin.
Iteration: 171 Time: 8.550E+00 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 172 Time: 8.600E+00 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 173 Time: 8.650E+00 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 174 Time: 8.700E+00 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 175 Time: 8.750E+00 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 176 Time: 8.800E+00 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 177 Time: 8.850E+00 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 178 Time: 8.900E+00 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 179 Time: 8.950E+00 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 180 Time: 9.000E+00 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Writing solution file op.bin.
Iteration: 181 Time: 9.050E+00 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 182 Time: 9.100E+00 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 183 Time: 9.150E+00 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 184 Time: 9.200E+00 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 185 Time: 9.250E+00 Max CFL: 1.014E+02 Max Diff. No.: -1.000E+00
Iteration: 186 Time: 9.300E+00 Max CFL: 1.014E+02 Max Diff. No.: -1.000E+00
Iteration: 187 Time: 9.350E+00 Max CFL: 1.014E+02 Max Diff. No.: -1.000E+00
Iteration: 188 Time: 9.400E+00 Max CFL: 1.014E+02 Max Diff. No.: -1.000E+00
Iteration: 189 Time: 9.450E+00 Max CFL: 1.014E+02 Max Diff. No.: -1.000E+00
Iteration: 190 Time: 9.500E+00 Max CFL: 1.014E+02 Max Diff. No.: -1.000E+00
Writing solution file op.bin.
Iteration: 191 Time: 9.550E+00 Max CFL: 1.014E+02 Max Diff. No.: -1.000E+00
Iteration: 192 Time: 9.600E+00 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 193 Time: 9.650E+00 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 194 Time: 9.700E+00 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 195 Time: 9.750E+00 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 196 Time: 9.800E+00 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 197 Time: 9.850E+00 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 198 Time: 9.900E+00 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 199 Time: 9.950E+00 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 200 Time: 1.000E+01 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Writing solution file op.bin.
Iteration: 201 Time: 1.005E+01 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 202 Time: 1.010E+01 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 203 Time: 1.015E+01 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 204 Time: 1.020E+01 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 205 Time: 1.025E+01 Max CFL: 1.014E+02 Max Diff. No.: -1.000E+00
Iteration: 206 Time: 1.030E+01 Max CFL: 1.014E+02 Max Diff. No.: -1.000E+00
Iteration: 207 Time: 1.035E+01 Max CFL: 1.014E+02 Max Diff. No.: -1.000E+00
Iteration: 208 Time: 1.040E+01 Max CFL: 1.014E+02 Max Diff. No.: -1.000E+00
Iteration: 209 Time: 1.045E+01 Max CFL: 1.014E+02 Max Diff. No.: -1.000E+00
Iteration: 210 Time: 1.050E+01 Max CFL: 1.014E+02 Max Diff. No.: -1.000E+00
Writing solution file op.bin.
Iteration: 211 Time: 1.055E+01 Max CFL: 1.014E+02 Max Diff. No.: -1.000E+00
Iteration: 212 Time: 1.060E+01 Max CFL: 1.014E+02 Max Diff. No.: -1.000E+00
Iteration: 213 Time: 1.065E+01 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 214 Time: 1.070E+01 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 215 Time: 1.075E+01 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 216 Time: 1.080E+01 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 217 Time: 1.085E+01 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 218 Time: 1.090E+01 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 219 Time: 1.095E+01 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 220 Time: 1.100E+01 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Writing solution file op.bin.
Iteration: 221 Time: 1.105E+01 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 222 Time: 1.110E+01 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 223 Time: 1.115E+01 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 224 Time: 1.120E+01 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 225 Time: 1.125E+01 Max CFL: 1.014E+02 Max Diff. No.: -1.000E+00
Iteration: 226 Time: 1.130E+01 Max CFL: 1.014E+02 Max Diff. No.: -1.000E+00
Iteration: 227 Time: 1.135E+01 Max CFL: 1.014E+02 Max Diff. No.: -1.000E+00
Iteration: 228 Time: 1.140E+01 Max CFL: 1.014E+02 Max Diff. No.: -1.000E+00
Iteration: 229 Time: 1.145E+01 Max CFL: 1.014E+02 Max Diff. No.: -1.000E+00
Iteration: 230 Time: 1.150E+01 Max CFL: 1.014E+02 Max Diff. No.: -1.000E+00
Writing solution file op.bin.
Iteration: 231 Time: 1.155E+01 Max CFL: 1.014E+02 Max Diff. No.: -1.000E+00
Iteration: 232 Time: 1.160E+01 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 233 Time: 1.165E+01 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 234 Time: 1.170E+01 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 235 Time: 1.175E+01 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 236 Time: 1.180E+01 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 237 Time: 1.185E+01 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 238 Time: 1.190E+01 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 239 Time: 1.195E+01 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 240 Time: 1.200E+01 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Writing solution file op.bin.
Iteration: 241 Time: 1.205E+01 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 242 Time: 1.210E+01 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 243 Time: 1.215E+01 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 244 Time: 1.220E+01 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 245 Time: 1.225E+01 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 246 Time: 1.230E+01 Max CFL: 1.014E+02 Max Diff. No.: -1.000E+00
Iteration: 247 Time: 1.235E+01 Max CFL: 1.014E+02 Max Diff. No.: -1.000E+00
Iteration: 248 Time: 1.240E+01 Max CFL: 1.014E+02 Max Diff. No.: -1.000E+00
Iteration: 249 Time: 1.245E+01 Max CFL: 1.014E+02 Max Diff. No.: -1.000E+00
Iteration: 250 Time: 1.250E+01 Max CFL: 1.014E+02 Max Diff. No.: -1.000E+00
Writing solution file op.bin.
Iteration: 251 Time: 1.255E+01 Max CFL: 1.014E+02 Max Diff. No.: -1.000E+00
Iteration: 252 Time: 1.260E+01 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 253 Time: 1.265E+01 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 254 Time: 1.270E+01 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 255 Time: 1.275E+01 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 256 Time: 1.280E+01 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 257 Time: 1.285E+01 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 258 Time: 1.290E+01 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 259 Time: 1.295E+01 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 260 Time: 1.300E+01 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Writing solution file op.bin.
Iteration: 261 Time: 1.305E+01 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 262 Time: 1.310E+01 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 263 Time: 1.315E+01 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 264 Time: 1.320E+01 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 265 Time: 1.325E+01 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 266 Time: 1.330E+01 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 267 Time: 1.335E+01 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 268 Time: 1.340E+01 Max CFL: 1.014E+02 Max Diff. No.: -1.000E+00
Iteration: 269 Time: 1.345E+01 Max CFL: 1.014E+02 Max Diff. No.: -1.000E+00
Iteration: 270 Time: 1.350E+01 Max CFL: 1.014E+02 Max Diff. No.: -1.000E+00
Writing solution file op.bin.
Iteration: 271 Time: 1.355E+01 Max CFL: 1.014E+02 Max Diff. No.: -1.000E+00
Iteration: 272 Time: 1.360E+01 Max CFL: 1.014E+02 Max Diff. No.: -1.000E+00
Iteration: 273 Time: 1.365E+01 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 274 Time: 1.370E+01 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 275 Time: 1.375E+01 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 276 Time: 1.380E+01 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 277 Time: 1.385E+01 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 278 Time: 1.390E+01 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 279 Time: 1.395E+01 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 280 Time: 1.400E+01 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Writing solution file op.bin.
Iteration: 281 Time: 1.405E+01 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 282 Time: 1.410E+01 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 283 Time: 1.415E+01 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 284 Time: 1.420E+01 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 285 Time: 1.425E+01 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 286 Time: 1.430E+01 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 287 Time: 1.435E+01 Max CFL: 1.014E+02 Max Diff. No.: -1.000E+00
Iteration: 288 Time: 1.440E+01 Max CFL: 1.014E+02 Max Diff. No.: -1.000E+00
Iteration: 289 Time: 1.445E+01 Max CFL: 1.014E+02 Max Diff. No.: -1.000E+00
Iteration: 290 Time: 1.450E+01 Max CFL: 1.014E+02 Max Diff. No.: -1.000E+00
Writing solution file op.bin.
Iteration: 291 Time: 1.455E+01 Max CFL: 1.014E+02 Max Diff. No.: -1.000E+00
Iteration: 292 Time: 1.460E+01 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 293 Time: 1.465E+01 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 294 Time: 1.470E+01 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 295 Time: 1.475E+01 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 296 Time: 1.480E+01 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 297 Time: 1.485E+01 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 298 Time: 1.490E+01 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 299 Time: 1.495E+01 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 300 Time: 1.500E+01 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Writing solution file op.bin.
Iteration: 301 Time: 1.505E+01 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 302 Time: 1.510E+01 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 303 Time: 1.515E+01 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 304 Time: 1.520E+01 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 305 Time: 1.525E+01 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 306 Time: 1.530E+01 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 307 Time: 1.535E+01 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 308 Time: 1.540E+01 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 309 Time: 1.545E+01 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 310 Time: 1.550E+01 Max CFL: 1.014E+02 Max Diff. No.: -1.000E+00
Writing solution file op.bin.
Iteration: 311 Time: 1.555E+01 Max CFL: 1.014E+02 Max Diff. No.: -1.000E+00
Iteration: 312 Time: 1.560E+01 Max CFL: 1.014E+02 Max Diff. No.: -1.000E+00
Iteration: 313 Time: 1.565E+01 Max CFL: 1.014E+02 Max Diff. No.: -1.000E+00
Iteration: 314 Time: 1.570E+01 Max CFL: 1.014E+02 Max Diff. No.: -1.000E+00
Iteration: 315 Time: 1.575E+01 Max CFL: 1.014E+02 Max Diff. No.: -1.000E+00
Iteration: 316 Time: 1.580E+01 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 317 Time: 1.585E+01 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 318 Time: 1.590E+01 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 319 Time: 1.595E+01 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 320 Time: 1.600E+01 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Writing solution file op.bin.
Iteration: 321 Time: 1.605E+01 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 322 Time: 1.610E+01 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 323 Time: 1.615E+01 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 324 Time: 1.620E+01 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 325 Time: 1.625E+01 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 326 Time: 1.630E+01 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 327 Time: 1.635E+01 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 328 Time: 1.640E+01 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 329 Time: 1.645E+01 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 330 Time: 1.650E+01 Max CFL: 1.014E+02 Max Diff. No.: -1.000E+00
Writing solution file op.bin.
Iteration: 331 Time: 1.655E+01 Max CFL: 1.014E+02 Max Diff. No.: -1.000E+00
Iteration: 332 Time: 1.660E+01 Max CFL: 1.014E+02 Max Diff. No.: -1.000E+00
Iteration: 333 Time: 1.665E+01 Max CFL: 1.014E+02 Max Diff. No.: -1.000E+00
Iteration: 334 Time: 1.670E+01 Max CFL: 1.014E+02 Max Diff. No.: -1.000E+00
Iteration: 335 Time: 1.675E+01 Max CFL: 1.014E+02 Max Diff. No.: -1.000E+00
Iteration: 336 Time: 1.680E+01 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 337 Time: 1.685E+01 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 338 Time: 1.690E+01 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 339 Time: 1.695E+01 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 340 Time: 1.700E+01 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Writing solution file op.bin.
Iteration: 341 Time: 1.705E+01 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 342 Time: 1.710E+01 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 343 Time: 1.715E+01 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 344 Time: 1.720E+01 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 345 Time: 1.725E+01 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 346 Time: 1.730E+01 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 347 Time: 1.735E+01 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 348 Time: 1.740E+01 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 349 Time: 1.745E+01 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 350 Time: 1.750E+01 Max CFL: 1.014E+02 Max Diff. No.: -1.000E+00
Writing solution file op.bin.
Iteration: 351 Time: 1.755E+01 Max CFL: 1.014E+02 Max Diff. No.: -1.000E+00
Iteration: 352 Time: 1.760E+01 Max CFL: 1.014E+02 Max Diff. No.: -1.000E+00
Iteration: 353 Time: 1.765E+01 Max CFL: 1.014E+02 Max Diff. No.: -1.000E+00
Iteration: 354 Time: 1.770E+01 Max CFL: 1.014E+02 Max Diff. No.: -1.000E+00
Iteration: 355 Time: 1.775E+01 Max CFL: 1.014E+02 Max Diff. No.: -1.000E+00
Iteration: 356 Time: 1.780E+01 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 357 Time: 1.785E+01 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 358 Time: 1.790E+01 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 359 Time: 1.795E+01 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 360 Time: 1.800E+01 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Writing solution file op.bin.
Iteration: 361 Time: 1.805E+01 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 362 Time: 1.810E+01 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 363 Time: 1.815E+01 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 364 Time: 1.820E+01 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 365 Time: 1.825E+01 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 366 Time: 1.830E+01 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 367 Time: 1.835E+01 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 368 Time: 1.840E+01 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 369 Time: 1.845E+01 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 370 Time: 1.850E+01 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Writing solution file op.bin.
Iteration: 371 Time: 1.855E+01 Max CFL: 1.014E+02 Max Diff. No.: -1.000E+00
Iteration: 372 Time: 1.860E+01 Max CFL: 1.014E+02 Max Diff. No.: -1.000E+00
Iteration: 373 Time: 1.865E+01 Max CFL: 1.014E+02 Max Diff. No.: -1.000E+00
Iteration: 374 Time: 1.870E+01 Max CFL: 1.014E+02 Max Diff. No.: -1.000E+00
Iteration: 375 Time: 1.875E+01 Max CFL: 1.014E+02 Max Diff. No.: -1.000E+00
Iteration: 376 Time: 1.880E+01 Max CFL: 1.014E+02 Max Diff. No.: -1.000E+00
Iteration: 377 Time: 1.885E+01 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 378 Time: 1.890E+01 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 379 Time: 1.895E+01 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 380 Time: 1.900E+01 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Writing solution file op.bin.
Iteration: 381 Time: 1.905E+01 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 382 Time: 1.910E+01 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 383 Time: 1.915E+01 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 384 Time: 1.920E+01 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 385 Time: 1.925E+01 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 386 Time: 1.930E+01 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 387 Time: 1.935E+01 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 388 Time: 1.940E+01 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 389 Time: 1.945E+01 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 390 Time: 1.950E+01 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Writing solution file op.bin.
Iteration: 391 Time: 1.955E+01 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 392 Time: 1.960E+01 Max CFL: 1.014E+02 Max Diff. No.: -1.000E+00
Iteration: 393 Time: 1.965E+01 Max CFL: 1.014E+02 Max Diff. No.: -1.000E+00
Iteration: 394 Time: 1.970E+01 Max CFL: 1.014E+02 Max Diff. No.: -1.000E+00
Iteration: 395 Time: 1.975E+01 Max CFL: 1.014E+02 Max Diff. No.: -1.000E+00
Iteration: 396 Time: 1.980E+01 Max CFL: 1.014E+02 Max Diff. No.: -1.000E+00
Iteration: 397 Time: 1.985E+01 Max CFL: 1.014E+02 Max Diff. No.: -1.000E+00
Iteration: 398 Time: 1.990E+01 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 399 Time: 1.995E+01 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Iteration: 400 Time: 2.000E+01 Max CFL: 1.015E+02 Max Diff. No.: -1.000E+00
Writing solution file op.bin.
** Completed PETSc time integration **
Computed errors:
L1 Error : 0.0000000000000000E+00
L2 Error : 0.0000000000000000E+00
Linfinity Error : 0.0000000000000000E+00
Conservation Errors:
0.0000000000000000E+00
0.0000000000000000E+00
0.0000000000000000E+00
0.0000000000000000E+00
Solver runtime (in seconds): 1.9543866630000000E+03
Total runtime (in seconds): 1.9544529339999999E+03
Deallocating arrays.
Finished.