HyPar  1.0
Finite-Difference Hyperbolic-Parabolic PDE Solver on Cartesian Grids
Steady, incompressible, viscous flow around a cylinder

Location: hypar/Examples/3D/NavierStokes3D/2D_Cylinder/Steady_Viscous_Incompressible

Governing equations: 3D Navier-Stokes Equations (navierstokes3d.h)

Domain: The domain consists of a fine uniform grid around the cylinder defined by [-2,6] X [-2,2], and a stretched grid beyond this zone. Note: This is a 2D flow simulated using a 3-dimensional setup by taking the length of the domain along z to be very small and with only 3 grid points (the domain size along z must be smaller than the cylinder length).

Geometry: A cylinder of radius 1.0 centered at (0,0) (hypar/Examples/STLGeometries/cylinder.stl)

The following images shows the grid and the cylinder:

Domain3D_Cylinder.png
Domain2D_Cylinder.png

Boundary conditions:

Reference:

  • Taneda, S., "Experimental Investigation of the Wakes behind Cylinders and Plates at Low Reynolds Numbers," Journal of the Physical Society of Japan, Vol. 11, 302–307, 1956.
  • Dennis, S. C. R., Chang, G.-Z., "Numerical solutions for steady flow past a circular cylinder at Reynolds numbers up to 100", Journal of Fluid Mechanics, 42 (3), 1970, pp. 471-489.

Initial solution: \(\rho=1, u=0.1, v=w=0, p=1/\gamma\) everywhere in the domain.

Other parameters (all dimensional quantities are in SI units):

Numerical Method:

Input files required:

These files are all located in: hypar/Examples/3D/NavierStokes3D/2D_Cylinder/Steady_Viscous_Incompressible/

solver.inp

begin
ndims 3
nvars 5
size 200 96 3
iproc 2 2 1
ghost 3
n_iter 50000
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.04
screen_op_iter 200
file_op_iter 1000
ip_file_type binary
op_file_format binary
op_overwrite yes
model navierstokes3d
immersed_body cylinder.stl
end

boundary.inp

6
subsonic-inflow 0 1 0 0 -9999.0 9999.0 -9999.0 9999.0
1.0 0.1 0.0 0.0
subsonic-outflow 0 -1 0 0 -9999.0 9999.0 -9999.0 9999.0
0.7142857142857143
subsonic-ambivalent 1 1 -9999.0 9999.0 0 0 -9999.0 9999.0
1.0 0.1 0.0 0.0 0.7142857142857143
subsonic-ambivalent 1 -1 -9999.0 9999.0 0 0 -9999.0 9999.0
1.0 0.1 0.0 0.0 0.7142857142857143
periodic 2 1 -9999.0 9999.0 -9999.0 9999.0 0 0
periodic 2 -1 -9999.0 9999.0 -9999.0 9999.0 0 0

physics.inp : The following file specifies a Reynolds number of 10 (corresponding to \(Re_D\) of 20). To try other Reynolds numbers, change it here.

begin
gamma 1.4
upwinding roe
Pr 0.72
Minf 0.1
Re 10
end

cylinder.stl : the filename "cylinder.stl" must match the input for immersed_body in solver.inp.
Located at hypar/Examples/STLGeometries/cylinder.stl

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>
int main()
{
const double GAMMA = 1.4;
int NI,NJ,NK,ndims;
char ip_file_type[50];
FILE *in;
strcpy(ip_file_type,"ascii");
printf("Reading file \"solver.inp\"...\n");
in = fopen("solver.inp","r");
if (!in) {
printf("Error: Input file \"solver.inp\" not found.\n");
return(0);
} 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);
fscanf(in,"%d",&NK);
} 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");
return(0);
}
}
fclose(in);
if (ndims != 3) {
printf("Error: ndims is not 3 in solver.inp.\n");
return(0);
}
printf("Grid: %3d x %3d x %3d.\n",NI,NJ,NK);
double Lx = 8.0;
double Ly = 4.0;
double sf_x1 = 1.15;
double sf_x2 = 1.20;
double sf_y = 1.4;
int i,j,k;
double u = 0.1,
v = 0.0,
w = 0.0,
rho = 1.0,
P = 1.0/GAMMA;
double *x, *y, *z, *U;
double dx, dy, dz;
FILE *out;
x = (double*) calloc (NI , sizeof(double));
y = (double*) calloc (NJ , sizeof(double));
z = (double*) calloc (NK , sizeof(double));
U = (double*) calloc (5*NI*NJ*NK, sizeof(double));
x[NI/8] = -Lx/4;
x[7*NI/8] = 3*Lx/4;
dx = Lx / ((double)(6*NI/8));
for (i = NI/8 ; i < 7*NI/8; i++) x[i] = x[NI/8] + dx * (i-NI/8);
for (i = 7*NI/8; i < NI ; i++) x[i] = x[i-1] + sf_x2 * (x[i-1] - x[i-2]);
for (i = NI/8-1; i >= 0 ; i--) x[i] = x[i+1] - sf_x1 * (x[i+2] - x[i+1]);
y[NJ/8] = -Ly/2;
y[7*NJ/8] = Ly/2;
dy = Ly / ((double)(6*NJ/8));
for (j = NJ/8 ; j < 7*NJ/8; j++) y[j] = y[NJ/8] + dy * (j-NJ/8);
for (j = 7*NJ/8; j < NJ ; j++) y[j] = y[j-1] + sf_y * (y[j-1] - y[j-2]);
for (j = NJ/8-1; j >= 0 ; j--) y[j] = y[j+1] - sf_y * (y[j+2] - y[j+1]);
dz = 0.5 * (dx +dy);
for (k = 0; k < NK; k++) z[k] = (k-NK/2) * dz;
for (i = 0; i < NI; i++){
for (j = 0; j < NJ; j++){
for (k = 0; k < NK; k++){
int p = i + NI*j + NI*NJ*k;
U[5*p+0] = rho;
U[5*p+1] = rho*u;
U[5*p+2] = rho*v;
U[5*p+3] = rho*w;
U[5*p+4] = P/(GAMMA-1.0) + 0.5*rho*(u*u+v*v+w*w);
}
}
}
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 (k = 0; k < NK; k++) fprintf(out,"%1.16E ",z[k]);
fprintf(out,"\n");
for (k = 0; k < NK; k++) {
for (j = 0; j < NJ; j++) {
for (i = 0; i < NI; i++) {
int p = i + NK*j + NI*NJ*k;
fprintf(out,"%1.16E ",U[5*p+0]);
}
}
}
fprintf(out,"\n");
for (k = 0; k < NK; k++) {
for (j = 0; j < NJ; j++) {
for (i = 0; i < NI; i++) {
int p = i + NK*j + NI*NJ*k;
fprintf(out,"%1.16E ",U[5*p+1]);
}
}
}
fprintf(out,"\n");
for (k = 0; k < NK; k++) {
for (j = 0; j < NJ; j++) {
for (i = 0; i < NI; i++) {
int p = i + NK*j + NI*NJ*k;
fprintf(out,"%1.16E ",U[5*p+2]);
}
}
}
fprintf(out,"\n");
for (k = 0; k < NK; k++) {
for (j = 0; j < NJ; j++) {
for (i = 0; i < NI; i++) {
int p = i + NK*j + NI*NJ*k;
fprintf(out,"%1.16E ",U[5*p+3]);
}
}
}
fprintf(out,"\n");
for (k = 0; k < NK; k++) {
for (j = 0; j < NJ; j++) {
for (i = 0; i < NI; i++) {
int p = i + NK*j + NI*NJ*k;
fprintf(out,"%1.16E ",U[5*p+4]);
}
}
}
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(z,sizeof(double),NK,out);
fwrite(U,sizeof(double),5*NI*NJ*NK,out);
fclose(out);
}
free(x);
free(y);
free(z);
free(U);
return(0);
}

Output:

Note that iproc is set to

  2 2 1

in solver.inp (i.e., 2 processors along x, 2 processors along y, and 1 processor along z). Thus, this example should be run with 4 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.
  • hypar/Extras/BinaryToText.c - convert binary output file to an ASCII text file (to visualize in, for example, MATLAB).

The file Extras/ExtractSlice.c can be used to extract a slice perpendicular to any dimension at a specified location along that dimension. The extracted slice is written out in the same binary format as the original solutions files (with the same names op_xxxxx.bin) in a subdirectory called slices (Note: make the subdirectory called slices before running this code). The codes Extras/BinaryToTecplot.c and Extras/BinaryToText.c can then be used (in the slices subdirectory) to convert the binary slice solution file to Tecplot or plain text files. Note that it needs the relevant solver.inp that can be created as follows:

  • Copy the original solver.inp to the slices subdirectory.
  • In solver.inp, set ndims as 2, and remove the component of size and iproc corresponding to the dimension being eliminated while extracting the slices (in this case, it is z or the 3rd component).

The following figure shows the flow (density and streamlines) at \(Re_D=20\):

Solution_3DNavStokCylinder_ReD020.png

The parameter Re can be changed to 20 in physics.inp to run this simulation for \(Re_D=40\), and following figure shows the solution:

Solution_3DNavStokCylinder_ReD040.png

The following plot shows the wake length ( \(L/D\)) as a function of the Reynolds number ( \(Re_D\)) for the computed solutions and experimental results reported in the reference above:

Solution_3DNavStokCylinder_WL.png

In addition to the main solution, the code also writes out a file with the aerodynamic forces on the immersed body. This file is called surface.dat (if HyPar::op_overwrite is "yes") or surface_nnnnn.dat (if HyPar::op_overwrite is "no", "nnnnn" is a numerical index) (in this example, the file surface.dat is written out). This is an ASCII file in the Tecplot format, where the immersed body and the forces on it are represented using the "FETRIANGLE" type. The following image shows the surface pressure on the cylinder (front-view):

IBSurface_3DNavStokCylinder.png

Since this is a 2D simulation, the value of the surface pressure on the end-surfaces of the cylinder (zmin and zmax) are not physically relevant.

Expected screen output (for \(Re_D = 20\)):

HyPar - Parallel (MPI) version with 4 processes
Compiled with PETSc time integration.
Reading solver inputs from file "solver.inp".
No. of dimensions : 3
No. of variables : 5
Domain size : 200 96 3
Processes along each dimension : 2 2 1
No. of ghosts pts : 3
No. of iter. : 50000
Restart iteration : 0
Time integration scheme : rk (44)
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 : 4.000000E-02
Check for conservation : no
Screen output iterations : 200
File output iterations : 1000
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 : navierstokes3d
Immersed Body : cylinder.stl
Partitioning domain.
Allocating data arrays.
Reading array from binary file initial.inp (Serial mode).
Volume integral of the initial solution:
0: 2.5654236261034362E+02
1: 2.5654236261037145E+01
2: 0.0000000000000000E+00
3: 0.0000000000000000E+00
4: 4.5939407361726899E+02
Reading boundary conditions from "boundary.inp".
Boundary subsonic-inflow: Along dimension 0 and face +1
Boundary subsonic-outflow: Along dimension 0 and face -1
Boundary subsonic-ambivalent: Along dimension 1 and face +1
Boundary subsonic-ambivalent: Along dimension 1 and face -1
Boundary periodic: Along dimension 2 and face +1
Boundary periodic: Along dimension 2 and face -1
6 boundary condition(s) read.
Immersed body read from cylinder.stl:
Number of facets: 628
Bounding box: [-1.0,1.0] X [-1.0,1.0] X [-5.0,5.0]
Number of grid points inside immersed body: 3288 ( 5.7%).
Number of immersed boundary points : 882 ( 1.5%).
Immersed body simulation mode : 2d (xy).
Initializing solvers.
Warning: File weno.inp not found. Using default parameters for WENO5/CRWENO5/HCWENO5 scheme.
Initializing physics. Model = "navierstokes3d"
Reading physical model inputs from file "physics.inp".
Writing solution file iblank.bin.
Setting up time integration.
Solving in time (from 0 to 50000 iterations)
Writing immersed body surface data file surface.dat.
Writing solution file op.bin.
Iteration: 200 Time: 8.000E+00 Max CFL: 8.612E-01 Max Diff. No.: -1.000E+00 Norm: 2.8716E-04
Iteration: 400 Time: 1.600E+01 Max CFL: 8.540E-01 Max Diff. No.: -1.000E+00 Norm: 6.3423E-05
Iteration: 600 Time: 2.400E+01 Max CFL: 8.509E-01 Max Diff. No.: -1.000E+00 Norm: 5.6642E-05
Iteration: 800 Time: 3.200E+01 Max CFL: 8.511E-01 Max Diff. No.: -1.000E+00 Norm: 5.8059E-05
Iteration: 1000 Time: 4.000E+01 Max CFL: 8.455E-01 Max Diff. No.: -1.000E+00 Norm: 7.1796E-05
Writing immersed body surface data file surface.dat.
Writing solution file op.bin.
Iteration: 1200 Time: 4.800E+01 Max CFL: 8.452E-01 Max Diff. No.: -1.000E+00 Norm: 4.7222E-05
Iteration: 1400 Time: 5.600E+01 Max CFL: 8.450E-01 Max Diff. No.: -1.000E+00 Norm: 2.3224E-05
Iteration: 1600 Time: 6.400E+01 Max CFL: 8.443E-01 Max Diff. No.: -1.000E+00 Norm: 1.5946E-05
Iteration: 1800 Time: 7.200E+01 Max CFL: 8.437E-01 Max Diff. No.: -1.000E+00 Norm: 1.6644E-05
Iteration: 2000 Time: 8.000E+01 Max CFL: 8.433E-01 Max Diff. No.: -1.000E+00 Norm: 1.2999E-05
Writing immersed body surface data file surface.dat.
Writing solution file op.bin.
Iteration: 2200 Time: 8.800E+01 Max CFL: 8.421E-01 Max Diff. No.: -1.000E+00 Norm: 1.4633E-05
Iteration: 2400 Time: 9.600E+01 Max CFL: 8.417E-01 Max Diff. No.: -1.000E+00 Norm: 8.0353E-06
Iteration: 2600 Time: 1.040E+02 Max CFL: 8.416E-01 Max Diff. No.: -1.000E+00 Norm: 9.5690E-06
Iteration: 2800 Time: 1.120E+02 Max CFL: 8.413E-01 Max Diff. No.: -1.000E+00 Norm: 7.0334E-06
Iteration: 3000 Time: 1.200E+02 Max CFL: 8.407E-01 Max Diff. No.: -1.000E+00 Norm: 8.0227E-06
Writing immersed body surface data file surface.dat.
Writing solution file op.bin.
Iteration: 3200 Time: 1.280E+02 Max CFL: 8.406E-01 Max Diff. No.: -1.000E+00 Norm: 6.5859E-06
Iteration: 3400 Time: 1.360E+02 Max CFL: 8.406E-01 Max Diff. No.: -1.000E+00 Norm: 4.6185E-06
Iteration: 3600 Time: 1.440E+02 Max CFL: 8.406E-01 Max Diff. No.: -1.000E+00 Norm: 4.5334E-06
Iteration: 3800 Time: 1.520E+02 Max CFL: 8.406E-01 Max Diff. No.: -1.000E+00 Norm: 6.1795E-06
Iteration: 4000 Time: 1.600E+02 Max CFL: 8.401E-01 Max Diff. No.: -1.000E+00 Norm: 4.9332E-06
Writing immersed body surface data file surface.dat.
Writing solution file op.bin.
Iteration: 4200 Time: 1.680E+02 Max CFL: 8.399E-01 Max Diff. No.: -1.000E+00 Norm: 4.1792E-06
Iteration: 4400 Time: 1.760E+02 Max CFL: 8.401E-01 Max Diff. No.: -1.000E+00 Norm: 2.9561E-06
Iteration: 4600 Time: 1.840E+02 Max CFL: 8.401E-01 Max Diff. No.: -1.000E+00 Norm: 3.4669E-06
Iteration: 4800 Time: 1.920E+02 Max CFL: 8.399E-01 Max Diff. No.: -1.000E+00 Norm: 2.9673E-06
Iteration: 5000 Time: 2.000E+02 Max CFL: 8.396E-01 Max Diff. No.: -1.000E+00 Norm: 3.0482E-06
Writing immersed body surface data file surface.dat.
Writing solution file op.bin.
Iteration: 5200 Time: 2.080E+02 Max CFL: 8.393E-01 Max Diff. No.: -1.000E+00 Norm: 3.8329E-06
Iteration: 5400 Time: 2.160E+02 Max CFL: 8.392E-01 Max Diff. No.: -1.000E+00 Norm: 3.3113E-06
Iteration: 5600 Time: 2.240E+02 Max CFL: 8.394E-01 Max Diff. No.: -1.000E+00 Norm: 3.0367E-06
Iteration: 5800 Time: 2.320E+02 Max CFL: 8.395E-01 Max Diff. No.: -1.000E+00 Norm: 2.0487E-06
Iteration: 6000 Time: 2.400E+02 Max CFL: 8.393E-01 Max Diff. No.: -1.000E+00 Norm: 2.1751E-06
Writing immersed body surface data file surface.dat.
Writing solution file op.bin.
Iteration: 6200 Time: 2.480E+02 Max CFL: 8.392E-01 Max Diff. No.: -1.000E+00 Norm: 2.8692E-06
Iteration: 6400 Time: 2.560E+02 Max CFL: 8.391E-01 Max Diff. No.: -1.000E+00 Norm: 1.3672E-06
Iteration: 6600 Time: 2.640E+02 Max CFL: 8.392E-01 Max Diff. No.: -1.000E+00 Norm: 2.3427E-06
Iteration: 6800 Time: 2.720E+02 Max CFL: 8.392E-01 Max Diff. No.: -1.000E+00 Norm: 3.3279E-06
Iteration: 7000 Time: 2.800E+02 Max CFL: 8.391E-01 Max Diff. No.: -1.000E+00 Norm: 1.4366E-06
Writing immersed body surface data file surface.dat.
Writing solution file op.bin.
Iteration: 7200 Time: 2.880E+02 Max CFL: 8.389E-01 Max Diff. No.: -1.000E+00 Norm: 2.6672E-06
Iteration: 7400 Time: 2.960E+02 Max CFL: 8.389E-01 Max Diff. No.: -1.000E+00 Norm: 2.5985E-06
Iteration: 7600 Time: 3.040E+02 Max CFL: 8.389E-01 Max Diff. No.: -1.000E+00 Norm: 1.2827E-06
Iteration: 7800 Time: 3.120E+02 Max CFL: 8.389E-01 Max Diff. No.: -1.000E+00 Norm: 1.8008E-06
Iteration: 8000 Time: 3.200E+02 Max CFL: 8.388E-01 Max Diff. No.: -1.000E+00 Norm: 1.4720E-06
Writing immersed body surface data file surface.dat.
Writing solution file op.bin.
Iteration: 8200 Time: 3.280E+02 Max CFL: 8.387E-01 Max Diff. No.: -1.000E+00 Norm: 1.4423E-06
Iteration: 8400 Time: 3.360E+02 Max CFL: 8.386E-01 Max Diff. No.: -1.000E+00 Norm: 1.4104E-06
Iteration: 8600 Time: 3.440E+02 Max CFL: 8.388E-01 Max Diff. No.: -1.000E+00 Norm: 1.0364E-06
Iteration: 8800 Time: 3.520E+02 Max CFL: 8.388E-01 Max Diff. No.: -1.000E+00 Norm: 1.3879E-06
Iteration: 9000 Time: 3.600E+02 Max CFL: 8.388E-01 Max Diff. No.: -1.000E+00 Norm: 1.2611E-06
Writing immersed body surface data file surface.dat.
Writing solution file op.bin.
Iteration: 9200 Time: 3.680E+02 Max CFL: 8.386E-01 Max Diff. No.: -1.000E+00 Norm: 1.0633E-06
Iteration: 9400 Time: 3.760E+02 Max CFL: 8.385E-01 Max Diff. No.: -1.000E+00 Norm: 1.1508E-06
Iteration: 9600 Time: 3.840E+02 Max CFL: 8.385E-01 Max Diff. No.: -1.000E+00 Norm: 1.3081E-06
Iteration: 9800 Time: 3.920E+02 Max CFL: 8.386E-01 Max Diff. No.: -1.000E+00 Norm: 1.0999E-06
Iteration: 10000 Time: 4.000E+02 Max CFL: 8.386E-01 Max Diff. No.: -1.000E+00 Norm: 9.2836E-07
Writing immersed body surface data file surface.dat.
Writing solution file op.bin.
Iteration: 10200 Time: 4.080E+02 Max CFL: 8.386E-01 Max Diff. No.: -1.000E+00 Norm: 1.5148E-06
Iteration: 10400 Time: 4.160E+02 Max CFL: 8.385E-01 Max Diff. No.: -1.000E+00 Norm: 1.6311E-06
Iteration: 10600 Time: 4.240E+02 Max CFL: 8.384E-01 Max Diff. No.: -1.000E+00 Norm: 5.6694E-07
Iteration: 10800 Time: 4.320E+02 Max CFL: 8.385E-01 Max Diff. No.: -1.000E+00 Norm: 1.9044E-06
Iteration: 11000 Time: 4.400E+02 Max CFL: 8.385E-01 Max Diff. No.: -1.000E+00 Norm: 1.7657E-06
Writing immersed body surface data file surface.dat.
Writing solution file op.bin.
Iteration: 11200 Time: 4.480E+02 Max CFL: 8.384E-01 Max Diff. No.: -1.000E+00 Norm: 5.4494E-07
Iteration: 11400 Time: 4.560E+02 Max CFL: 8.384E-01 Max Diff. No.: -1.000E+00 Norm: 1.4810E-06
Iteration: 11600 Time: 4.640E+02 Max CFL: 8.384E-01 Max Diff. No.: -1.000E+00 Norm: 1.0132E-06
Iteration: 11800 Time: 4.720E+02 Max CFL: 8.384E-01 Max Diff. No.: -1.000E+00 Norm: 6.2476E-07
Iteration: 12000 Time: 4.800E+02 Max CFL: 8.384E-01 Max Diff. No.: -1.000E+00 Norm: 7.9453E-07
Writing immersed body surface data file surface.dat.
Writing solution file op.bin.
Iteration: 12200 Time: 4.880E+02 Max CFL: 8.383E-01 Max Diff. No.: -1.000E+00 Norm: 8.3038E-07
Iteration: 12400 Time: 4.960E+02 Max CFL: 8.383E-01 Max Diff. No.: -1.000E+00 Norm: 8.4823E-07
Iteration: 12600 Time: 5.040E+02 Max CFL: 8.383E-01 Max Diff. No.: -1.000E+00 Norm: 5.4890E-07
Iteration: 12800 Time: 5.120E+02 Max CFL: 8.383E-01 Max Diff. No.: -1.000E+00 Norm: 7.0400E-07
Iteration: 13000 Time: 5.200E+02 Max CFL: 8.384E-01 Max Diff. No.: -1.000E+00 Norm: 7.9581E-07
Writing immersed body surface data file surface.dat.
Writing solution file op.bin.
Iteration: 13200 Time: 5.280E+02 Max CFL: 8.383E-01 Max Diff. No.: -1.000E+00 Norm: 6.0061E-07
Iteration: 13400 Time: 5.360E+02 Max CFL: 8.382E-01 Max Diff. No.: -1.000E+00 Norm: 5.4059E-07
Iteration: 13600 Time: 5.440E+02 Max CFL: 8.382E-01 Max Diff. No.: -1.000E+00 Norm: 7.1623E-07
Iteration: 13800 Time: 5.520E+02 Max CFL: 8.382E-01 Max Diff. No.: -1.000E+00 Norm: 9.9992E-07
Iteration: 14000 Time: 5.600E+02 Max CFL: 8.383E-01 Max Diff. No.: -1.000E+00 Norm: 5.7588E-07
Writing immersed body surface data file surface.dat.
Writing solution file op.bin.
Iteration: 14200 Time: 5.680E+02 Max CFL: 8.382E-01 Max Diff. No.: -1.000E+00 Norm: 8.3110E-07
Iteration: 14400 Time: 5.760E+02 Max CFL: 8.382E-01 Max Diff. No.: -1.000E+00 Norm: 1.2286E-06
Iteration: 14600 Time: 5.840E+02 Max CFL: 8.382E-01 Max Diff. No.: -1.000E+00 Norm: 5.1600E-07
Iteration: 14800 Time: 5.920E+02 Max CFL: 8.382E-01 Max Diff. No.: -1.000E+00 Norm: 6.8904E-07
Iteration: 15000 Time: 6.000E+02 Max CFL: 8.381E-01 Max Diff. No.: -1.000E+00 Norm: 1.0995E-06
Writing immersed body surface data file surface.dat.
Writing solution file op.bin.
Iteration: 15200 Time: 6.080E+02 Max CFL: 8.381E-01 Max Diff. No.: -1.000E+00 Norm: 6.5784E-07
Iteration: 15400 Time: 6.160E+02 Max CFL: 8.381E-01 Max Diff. No.: -1.000E+00 Norm: 3.5429E-07
Iteration: 15600 Time: 6.240E+02 Max CFL: 8.381E-01 Max Diff. No.: -1.000E+00 Norm: 8.2808E-07
Iteration: 15800 Time: 6.320E+02 Max CFL: 8.381E-01 Max Diff. No.: -1.000E+00 Norm: 6.9343E-07
Iteration: 16000 Time: 6.400E+02 Max CFL: 8.381E-01 Max Diff. No.: -1.000E+00 Norm: 3.9166E-07
Writing immersed body surface data file surface.dat.
Writing solution file op.bin.
Iteration: 16200 Time: 6.480E+02 Max CFL: 8.381E-01 Max Diff. No.: -1.000E+00 Norm: 5.8707E-07
Iteration: 16400 Time: 6.560E+02 Max CFL: 8.380E-01 Max Diff. No.: -1.000E+00 Norm: 6.2642E-07
Iteration: 16600 Time: 6.640E+02 Max CFL: 8.380E-01 Max Diff. No.: -1.000E+00 Norm: 4.5526E-07
Iteration: 16800 Time: 6.720E+02 Max CFL: 8.381E-01 Max Diff. No.: -1.000E+00 Norm: 2.5147E-07
Iteration: 17000 Time: 6.800E+02 Max CFL: 8.381E-01 Max Diff. No.: -1.000E+00 Norm: 5.0853E-07
Writing immersed body surface data file surface.dat.
Writing solution file op.bin.
Iteration: 17200 Time: 6.880E+02 Max CFL: 8.381E-01 Max Diff. No.: -1.000E+00 Norm: 6.0967E-07
Iteration: 17400 Time: 6.960E+02 Max CFL: 8.380E-01 Max Diff. No.: -1.000E+00 Norm: 4.0915E-07
Iteration: 17600 Time: 7.040E+02 Max CFL: 8.380E-01 Max Diff. No.: -1.000E+00 Norm: 4.0966E-07
Iteration: 17800 Time: 7.120E+02 Max CFL: 8.380E-01 Max Diff. No.: -1.000E+00 Norm: 7.5952E-07
Iteration: 18000 Time: 7.200E+02 Max CFL: 8.380E-01 Max Diff. No.: -1.000E+00 Norm: 5.8736E-07
Writing immersed body surface data file surface.dat.
Writing solution file op.bin.
Iteration: 18200 Time: 7.280E+02 Max CFL: 8.380E-01 Max Diff. No.: -1.000E+00 Norm: 2.1899E-07
Iteration: 18400 Time: 7.360E+02 Max CFL: 8.380E-01 Max Diff. No.: -1.000E+00 Norm: 7.5336E-07
Iteration: 18600 Time: 7.440E+02 Max CFL: 8.380E-01 Max Diff. No.: -1.000E+00 Norm: 6.5657E-07
Iteration: 18800 Time: 7.520E+02 Max CFL: 8.380E-01 Max Diff. No.: -1.000E+00 Norm: 1.7811E-07
Iteration: 19000 Time: 7.600E+02 Max CFL: 8.380E-01 Max Diff. No.: -1.000E+00 Norm: 6.2469E-07
Writing immersed body surface data file surface.dat.
Writing solution file op.bin.
Iteration: 19200 Time: 7.680E+02 Max CFL: 8.379E-01 Max Diff. No.: -1.000E+00 Norm: 7.0291E-07
Iteration: 19400 Time: 7.760E+02 Max CFL: 8.379E-01 Max Diff. No.: -1.000E+00 Norm: 3.1273E-07
Iteration: 19600 Time: 7.840E+02 Max CFL: 8.380E-01 Max Diff. No.: -1.000E+00 Norm: 3.6727E-07
Iteration: 19800 Time: 7.920E+02 Max CFL: 8.380E-01 Max Diff. No.: -1.000E+00 Norm: 6.0523E-07
Iteration: 20000 Time: 8.000E+02 Max CFL: 8.380E-01 Max Diff. No.: -1.000E+00 Norm: 4.3255E-07
Writing immersed body surface data file surface.dat.
Writing solution file op.bin.
Iteration: 20200 Time: 8.080E+02 Max CFL: 8.379E-01 Max Diff. No.: -1.000E+00 Norm: 2.5776E-07
Iteration: 20400 Time: 8.160E+02 Max CFL: 8.379E-01 Max Diff. No.: -1.000E+00 Norm: 4.3176E-07
Iteration: 20600 Time: 8.240E+02 Max CFL: 8.379E-01 Max Diff. No.: -1.000E+00 Norm: 4.9758E-07
Iteration: 20800 Time: 8.320E+02 Max CFL: 8.379E-01 Max Diff. No.: -1.000E+00 Norm: 3.1633E-07
Iteration: 21000 Time: 8.400E+02 Max CFL: 8.379E-01 Max Diff. No.: -1.000E+00 Norm: 2.1786E-07
Writing immersed body surface data file surface.dat.
Writing solution file op.bin.
Iteration: 21200 Time: 8.480E+02 Max CFL: 8.379E-01 Max Diff. No.: -1.000E+00 Norm: 5.2598E-07
Iteration: 21400 Time: 8.560E+02 Max CFL: 8.379E-01 Max Diff. No.: -1.000E+00 Norm: 4.8853E-07
Iteration: 21600 Time: 8.640E+02 Max CFL: 8.379E-01 Max Diff. No.: -1.000E+00 Norm: 1.8436E-07
Iteration: 21800 Time: 8.720E+02 Max CFL: 8.379E-01 Max Diff. No.: -1.000E+00 Norm: 4.9345E-07
Iteration: 22000 Time: 8.800E+02 Max CFL: 8.379E-01 Max Diff. No.: -1.000E+00 Norm: 5.8239E-07
Writing immersed body surface data file surface.dat.
Writing solution file op.bin.
Iteration: 22200 Time: 8.880E+02 Max CFL: 8.379E-01 Max Diff. No.: -1.000E+00 Norm: 2.2815E-07
Iteration: 22400 Time: 8.960E+02 Max CFL: 8.379E-01 Max Diff. No.: -1.000E+00 Norm: 3.4285E-07
Iteration: 22600 Time: 9.040E+02 Max CFL: 8.379E-01 Max Diff. No.: -1.000E+00 Norm: 5.8436E-07
Iteration: 22800 Time: 9.120E+02 Max CFL: 8.379E-01 Max Diff. No.: -1.000E+00 Norm: 3.5835E-07
Iteration: 23000 Time: 9.200E+02 Max CFL: 8.378E-01 Max Diff. No.: -1.000E+00 Norm: 2.4104E-07
Writing immersed body surface data file surface.dat.
Writing solution file op.bin.
Iteration: 23200 Time: 9.280E+02 Max CFL: 8.378E-01 Max Diff. No.: -1.000E+00 Norm: 5.1939E-07
Iteration: 23400 Time: 9.360E+02 Max CFL: 8.378E-01 Max Diff. No.: -1.000E+00 Norm: 4.4233E-07
Iteration: 23600 Time: 9.440E+02 Max CFL: 8.378E-01 Max Diff. No.: -1.000E+00 Norm: 1.5591E-07
Iteration: 23800 Time: 9.520E+02 Max CFL: 8.378E-01 Max Diff. No.: -1.000E+00 Norm: 3.2582E-07
Iteration: 24000 Time: 9.600E+02 Max CFL: 8.378E-01 Max Diff. No.: -1.000E+00 Norm: 4.5036E-07
Writing immersed body surface data file surface.dat.
Writing solution file op.bin.
Iteration: 24200 Time: 9.680E+02 Max CFL: 8.378E-01 Max Diff. No.: -1.000E+00 Norm: 2.9978E-07
Iteration: 24400 Time: 9.760E+02 Max CFL: 8.378E-01 Max Diff. No.: -1.000E+00 Norm: 1.9134E-07
Iteration: 24600 Time: 9.840E+02 Max CFL: 8.378E-01 Max Diff. No.: -1.000E+00 Norm: 4.2376E-07
Iteration: 24800 Time: 9.920E+02 Max CFL: 8.378E-01 Max Diff. No.: -1.000E+00 Norm: 4.3198E-07
Iteration: 25000 Time: 1.000E+03 Max CFL: 8.378E-01 Max Diff. No.: -1.000E+00 Norm: 1.5994E-07
Writing immersed body surface data file surface.dat.
Writing solution file op.bin.
Iteration: 25200 Time: 1.008E+03 Max CFL: 8.378E-01 Max Diff. No.: -1.000E+00 Norm: 3.1032E-07
Iteration: 25400 Time: 1.016E+03 Max CFL: 8.378E-01 Max Diff. No.: -1.000E+00 Norm: 4.8410E-07
Iteration: 25600 Time: 1.024E+03 Max CFL: 8.378E-01 Max Diff. No.: -1.000E+00 Norm: 2.8539E-07
Iteration: 25800 Time: 1.032E+03 Max CFL: 8.378E-01 Max Diff. No.: -1.000E+00 Norm: 2.1036E-07
Iteration: 26000 Time: 1.040E+03 Max CFL: 8.378E-01 Max Diff. No.: -1.000E+00 Norm: 4.7019E-07
Writing immersed body surface data file surface.dat.
Writing solution file op.bin.
Iteration: 26200 Time: 1.048E+03 Max CFL: 8.378E-01 Max Diff. No.: -1.000E+00 Norm: 3.8263E-07
Iteration: 26400 Time: 1.056E+03 Max CFL: 8.378E-01 Max Diff. No.: -1.000E+00 Norm: 9.6516E-08
Iteration: 26600 Time: 1.064E+03 Max CFL: 8.378E-01 Max Diff. No.: -1.000E+00 Norm: 3.6475E-07
Iteration: 26800 Time: 1.072E+03 Max CFL: 8.378E-01 Max Diff. No.: -1.000E+00 Norm: 4.1589E-07
Iteration: 27000 Time: 1.080E+03 Max CFL: 8.378E-01 Max Diff. No.: -1.000E+00 Norm: 1.9194E-07
Writing immersed body surface data file surface.dat.
Writing solution file op.bin.
Iteration: 27200 Time: 1.088E+03 Max CFL: 8.378E-01 Max Diff. No.: -1.000E+00 Norm: 2.5096E-07
Iteration: 27400 Time: 1.096E+03 Max CFL: 8.378E-01 Max Diff. No.: -1.000E+00 Norm: 4.0283E-07
Iteration: 27600 Time: 1.104E+03 Max CFL: 8.378E-01 Max Diff. No.: -1.000E+00 Norm: 2.9946E-07
Iteration: 27800 Time: 1.112E+03 Max CFL: 8.378E-01 Max Diff. No.: -1.000E+00 Norm: 1.0476E-07
Iteration: 28000 Time: 1.120E+03 Max CFL: 8.378E-01 Max Diff. No.: -1.000E+00 Norm: 3.1693E-07
Writing immersed body surface data file surface.dat.
Writing solution file op.bin.
Iteration: 28200 Time: 1.128E+03 Max CFL: 8.378E-01 Max Diff. No.: -1.000E+00 Norm: 3.7616E-07
Iteration: 28400 Time: 1.136E+03 Max CFL: 8.378E-01 Max Diff. No.: -1.000E+00 Norm: 1.8538E-07
Iteration: 28600 Time: 1.144E+03 Max CFL: 8.377E-01 Max Diff. No.: -1.000E+00 Norm: 2.2525E-07
Iteration: 28800 Time: 1.152E+03 Max CFL: 8.377E-01 Max Diff. No.: -1.000E+00 Norm: 4.0184E-07
Iteration: 29000 Time: 1.160E+03 Max CFL: 8.377E-01 Max Diff. No.: -1.000E+00 Norm: 2.9509E-07
Writing immersed body surface data file surface.dat.
Writing solution file op.bin.
Iteration: 29200 Time: 1.168E+03 Max CFL: 8.378E-01 Max Diff. No.: -1.000E+00 Norm: 9.5927E-08
Iteration: 29400 Time: 1.176E+03 Max CFL: 8.378E-01 Max Diff. No.: -1.000E+00 Norm: 3.4870E-07
Iteration: 29600 Time: 1.184E+03 Max CFL: 8.378E-01 Max Diff. No.: -1.000E+00 Norm: 3.6288E-07
Iteration: 29800 Time: 1.192E+03 Max CFL: 8.377E-01 Max Diff. No.: -1.000E+00 Norm: 1.3083E-07
Iteration: 30000 Time: 1.200E+03 Max CFL: 8.377E-01 Max Diff. No.: -1.000E+00 Norm: 2.6300E-07
Writing immersed body surface data file surface.dat.
Writing solution file op.bin.
Iteration: 30200 Time: 1.208E+03 Max CFL: 8.377E-01 Max Diff. No.: -1.000E+00 Norm: 3.8056E-07
Iteration: 30400 Time: 1.216E+03 Max CFL: 8.377E-01 Max Diff. No.: -1.000E+00 Norm: 2.2423E-07
Iteration: 30600 Time: 1.224E+03 Max CFL: 8.377E-01 Max Diff. No.: -1.000E+00 Norm: 1.3543E-07
Iteration: 30800 Time: 1.232E+03 Max CFL: 8.378E-01 Max Diff. No.: -1.000E+00 Norm: 3.2728E-07
Iteration: 31000 Time: 1.240E+03 Max CFL: 8.377E-01 Max Diff. No.: -1.000E+00 Norm: 2.9447E-07
Writing immersed body surface data file surface.dat.
Writing solution file op.bin.
Iteration: 31200 Time: 1.248E+03 Max CFL: 8.377E-01 Max Diff. No.: -1.000E+00 Norm: 1.1169E-07
Iteration: 31400 Time: 1.256E+03 Max CFL: 8.377E-01 Max Diff. No.: -1.000E+00 Norm: 2.4691E-07
Iteration: 31600 Time: 1.264E+03 Max CFL: 8.377E-01 Max Diff. No.: -1.000E+00 Norm: 3.3174E-07
Iteration: 31800 Time: 1.272E+03 Max CFL: 8.377E-01 Max Diff. No.: -1.000E+00 Norm: 1.9605E-07
Iteration: 32000 Time: 1.280E+03 Max CFL: 8.377E-01 Max Diff. No.: -1.000E+00 Norm: 1.2799E-07
Writing immersed body surface data file surface.dat.
Writing solution file op.bin.
Iteration: 32200 Time: 1.288E+03 Max CFL: 8.377E-01 Max Diff. No.: -1.000E+00 Norm: 3.1298E-07
Iteration: 32400 Time: 1.296E+03 Max CFL: 8.377E-01 Max Diff. No.: -1.000E+00 Norm: 2.8125E-07
Iteration: 32600 Time: 1.304E+03 Max CFL: 8.377E-01 Max Diff. No.: -1.000E+00 Norm: 9.1322E-08
Iteration: 32800 Time: 1.312E+03 Max CFL: 8.377E-01 Max Diff. No.: -1.000E+00 Norm: 2.6033E-07
Iteration: 33000 Time: 1.320E+03 Max CFL: 8.377E-01 Max Diff. No.: -1.000E+00 Norm: 3.2770E-07
Writing immersed body surface data file surface.dat.
Writing solution file op.bin.
Iteration: 33200 Time: 1.328E+03 Max CFL: 8.377E-01 Max Diff. No.: -1.000E+00 Norm: 1.5822E-07
Iteration: 33400 Time: 1.336E+03 Max CFL: 8.377E-01 Max Diff. No.: -1.000E+00 Norm: 1.6008E-07
Iteration: 33600 Time: 1.344E+03 Max CFL: 8.377E-01 Max Diff. No.: -1.000E+00 Norm: 3.1519E-07
Iteration: 33800 Time: 1.352E+03 Max CFL: 8.377E-01 Max Diff. No.: -1.000E+00 Norm: 2.3733E-07
Iteration: 34000 Time: 1.360E+03 Max CFL: 8.377E-01 Max Diff. No.: -1.000E+00 Norm: 9.1190E-08
Writing immersed body surface data file surface.dat.
Writing solution file op.bin.
Iteration: 34200 Time: 1.368E+03 Max CFL: 8.377E-01 Max Diff. No.: -1.000E+00 Norm: 2.6325E-07
Iteration: 34400 Time: 1.376E+03 Max CFL: 8.377E-01 Max Diff. No.: -1.000E+00 Norm: 2.8337E-07
Iteration: 34600 Time: 1.384E+03 Max CFL: 8.377E-01 Max Diff. No.: -1.000E+00 Norm: 1.2001E-07
Iteration: 34800 Time: 1.392E+03 Max CFL: 8.377E-01 Max Diff. No.: -1.000E+00 Norm: 1.6551E-07
Iteration: 35000 Time: 1.400E+03 Max CFL: 8.377E-01 Max Diff. No.: -1.000E+00 Norm: 2.8223E-07
Writing immersed body surface data file surface.dat.
Writing solution file op.bin.
Iteration: 35200 Time: 1.408E+03 Max CFL: 8.377E-01 Max Diff. No.: -1.000E+00 Norm: 2.0503E-07
Iteration: 35400 Time: 1.416E+03 Max CFL: 8.377E-01 Max Diff. No.: -1.000E+00 Norm: 8.8727E-08
Iteration: 35600 Time: 1.424E+03 Max CFL: 8.377E-01 Max Diff. No.: -1.000E+00 Norm: 2.4851E-07
Iteration: 35800 Time: 1.432E+03 Max CFL: 8.377E-01 Max Diff. No.: -1.000E+00 Norm: 2.6411E-07
Iteration: 36000 Time: 1.440E+03 Max CFL: 8.377E-01 Max Diff. No.: -1.000E+00 Norm: 1.0088E-07
Writing immersed body surface data file surface.dat.
Writing solution file op.bin.
Iteration: 36200 Time: 1.448E+03 Max CFL: 8.377E-01 Max Diff. No.: -1.000E+00 Norm: 1.7408E-07
Iteration: 36400 Time: 1.456E+03 Max CFL: 8.377E-01 Max Diff. No.: -1.000E+00 Norm: 2.7814E-07
Iteration: 36600 Time: 1.464E+03 Max CFL: 8.377E-01 Max Diff. No.: -1.000E+00 Norm: 1.7851E-07
Iteration: 36800 Time: 1.472E+03 Max CFL: 8.377E-01 Max Diff. No.: -1.000E+00 Norm: 9.9721E-08
Iteration: 37000 Time: 1.480E+03 Max CFL: 8.377E-01 Max Diff. No.: -1.000E+00 Norm: 2.5491E-07
Writing immersed body surface data file surface.dat.
Writing solution file op.bin.
Iteration: 37200 Time: 1.488E+03 Max CFL: 8.377E-01 Max Diff. No.: -1.000E+00 Norm: 2.3551E-07
Iteration: 37400 Time: 1.496E+03 Max CFL: 8.377E-01 Max Diff. No.: -1.000E+00 Norm: 7.2067E-08
Iteration: 37600 Time: 1.504E+03 Max CFL: 8.377E-01 Max Diff. No.: -1.000E+00 Norm: 1.8881E-07
Iteration: 37800 Time: 1.512E+03 Max CFL: 8.377E-01 Max Diff. No.: -1.000E+00 Norm: 2.5326E-07
Iteration: 38000 Time: 1.520E+03 Max CFL: 8.377E-01 Max Diff. No.: -1.000E+00 Norm: 1.4126E-07
Writing immersed body surface data file surface.dat.
Writing solution file op.bin.
Iteration: 38200 Time: 1.528E+03 Max CFL: 8.377E-01 Max Diff. No.: -1.000E+00 Norm: 1.1245E-07
Iteration: 38400 Time: 1.536E+03 Max CFL: 8.377E-01 Max Diff. No.: -1.000E+00 Norm: 2.3730E-07
Iteration: 38600 Time: 1.544E+03 Max CFL: 8.377E-01 Max Diff. No.: -1.000E+00 Norm: 2.0423E-07
Iteration: 38800 Time: 1.552E+03 Max CFL: 8.377E-01 Max Diff. No.: -1.000E+00 Norm: 6.2945E-08
Iteration: 39000 Time: 1.560E+03 Max CFL: 8.377E-01 Max Diff. No.: -1.000E+00 Norm: 1.8383E-07
Writing immersed body surface data file surface.dat.
Writing solution file op.bin.
Iteration: 39200 Time: 1.568E+03 Max CFL: 8.377E-01 Max Diff. No.: -1.000E+00 Norm: 2.3538E-07
Iteration: 39400 Time: 1.576E+03 Max CFL: 8.377E-01 Max Diff. No.: -1.000E+00 Norm: 1.2224E-07
Iteration: 39600 Time: 1.584E+03 Max CFL: 8.377E-01 Max Diff. No.: -1.000E+00 Norm: 1.1692E-07
Iteration: 39800 Time: 1.592E+03 Max CFL: 8.377E-01 Max Diff. No.: -1.000E+00 Norm: 2.3277E-07
Iteration: 40000 Time: 1.600E+03 Max CFL: 8.377E-01 Max Diff. No.: -1.000E+00 Norm: 1.8441E-07
Writing immersed body surface data file surface.dat.
Writing solution file op.bin.
Iteration: 40200 Time: 1.608E+03 Max CFL: 8.377E-01 Max Diff. No.: -1.000E+00 Norm: 5.5882E-08
Iteration: 40400 Time: 1.616E+03 Max CFL: 8.377E-01 Max Diff. No.: -1.000E+00 Norm: 1.9215E-07
Iteration: 40600 Time: 1.624E+03 Max CFL: 8.377E-01 Max Diff. No.: -1.000E+00 Norm: 2.1808E-07
Iteration: 40800 Time: 1.632E+03 Max CFL: 8.377E-01 Max Diff. No.: -1.000E+00 Norm: 9.3465E-08
Iteration: 41000 Time: 1.640E+03 Max CFL: 8.377E-01 Max Diff. No.: -1.000E+00 Norm: 1.3166E-07
Writing immersed body surface data file surface.dat.
Writing solution file op.bin.
Iteration: 41200 Time: 1.648E+03 Max CFL: 8.377E-01 Max Diff. No.: -1.000E+00 Norm: 2.1996E-07
Iteration: 41400 Time: 1.656E+03 Max CFL: 8.377E-01 Max Diff. No.: -1.000E+00 Norm: 1.5302E-07
Iteration: 41600 Time: 1.664E+03 Max CFL: 8.377E-01 Max Diff. No.: -1.000E+00 Norm: 6.3539E-08
Iteration: 41800 Time: 1.672E+03 Max CFL: 8.377E-01 Max Diff. No.: -1.000E+00 Norm: 1.8667E-07
Iteration: 42000 Time: 1.680E+03 Max CFL: 8.377E-01 Max Diff. No.: -1.000E+00 Norm: 1.9292E-07
Writing immersed body surface data file surface.dat.
Writing solution file op.bin.
Iteration: 42200 Time: 1.688E+03 Max CFL: 8.377E-01 Max Diff. No.: -1.000E+00 Norm: 7.5612E-08
Iteration: 42400 Time: 1.696E+03 Max CFL: 8.377E-01 Max Diff. No.: -1.000E+00 Norm: 1.3301E-07
Iteration: 42600 Time: 1.704E+03 Max CFL: 8.377E-01 Max Diff. No.: -1.000E+00 Norm: 2.0518E-07
Iteration: 42800 Time: 1.712E+03 Max CFL: 8.377E-01 Max Diff. No.: -1.000E+00 Norm: 1.3333E-07
Iteration: 43000 Time: 1.720E+03 Max CFL: 8.377E-01 Max Diff. No.: -1.000E+00 Norm: 6.7674E-08
Writing immersed body surface data file surface.dat.
Writing solution file op.bin.
Iteration: 43200 Time: 1.728E+03 Max CFL: 8.377E-01 Max Diff. No.: -1.000E+00 Norm: 1.8463E-07
Iteration: 43400 Time: 1.736E+03 Max CFL: 8.377E-01 Max Diff. No.: -1.000E+00 Norm: 1.7728E-07
Iteration: 43600 Time: 1.744E+03 Max CFL: 8.377E-01 Max Diff. No.: -1.000E+00 Norm: 5.9122E-08
Iteration: 43800 Time: 1.752E+03 Max CFL: 8.377E-01 Max Diff. No.: -1.000E+00 Norm: 1.4089E-07
Iteration: 44000 Time: 1.760E+03 Max CFL: 8.377E-01 Max Diff. No.: -1.000E+00 Norm: 1.9474E-07
Writing immersed body surface data file surface.dat.
Writing solution file op.bin.
Iteration: 44200 Time: 1.768E+03 Max CFL: 8.377E-01 Max Diff. No.: -1.000E+00 Norm: 1.0951E-07
Iteration: 44400 Time: 1.776E+03 Max CFL: 8.377E-01 Max Diff. No.: -1.000E+00 Norm: 8.0062E-08
Iteration: 44600 Time: 1.784E+03 Max CFL: 8.377E-01 Max Diff. No.: -1.000E+00 Norm: 1.8043E-07
Iteration: 44800 Time: 1.792E+03 Max CFL: 8.377E-01 Max Diff. No.: -1.000E+00 Norm: 1.5413E-07
Iteration: 45000 Time: 1.800E+03 Max CFL: 8.377E-01 Max Diff. No.: -1.000E+00 Norm: 4.9563E-08
Writing immersed body surface data file surface.dat.
Writing solution file op.bin.
Iteration: 45200 Time: 1.808E+03 Max CFL: 8.377E-01 Max Diff. No.: -1.000E+00 Norm: 1.4241E-07
Iteration: 45400 Time: 1.816E+03 Max CFL: 8.377E-01 Max Diff. No.: -1.000E+00 Norm: 1.7652E-07
Iteration: 45600 Time: 1.824E+03 Max CFL: 8.377E-01 Max Diff. No.: -1.000E+00 Norm: 8.9007E-08
Iteration: 45800 Time: 1.832E+03 Max CFL: 8.377E-01 Max Diff. No.: -1.000E+00 Norm: 8.6250E-08
Iteration: 46000 Time: 1.840E+03 Max CFL: 8.377E-01 Max Diff. No.: -1.000E+00 Norm: 1.7093E-07
Writing immersed body surface data file surface.dat.
Writing solution file op.bin.
Iteration: 46200 Time: 1.848E+03 Max CFL: 8.377E-01 Max Diff. No.: -1.000E+00 Norm: 1.3614E-07
Iteration: 46400 Time: 1.856E+03 Max CFL: 8.377E-01 Max Diff. No.: -1.000E+00 Norm: 4.5794E-08
Iteration: 46600 Time: 1.864E+03 Max CFL: 8.377E-01 Max Diff. No.: -1.000E+00 Norm: 1.4233E-07
Iteration: 46800 Time: 1.872E+03 Max CFL: 8.377E-01 Max Diff. No.: -1.000E+00 Norm: 1.6384E-07
Iteration: 47000 Time: 1.880E+03 Max CFL: 8.377E-01 Max Diff. No.: -1.000E+00 Norm: 7.1933E-08
Writing immersed body surface data file surface.dat.
Writing solution file op.bin.
Iteration: 47200 Time: 1.888E+03 Max CFL: 8.377E-01 Max Diff. No.: -1.000E+00 Norm: 9.4149E-08
Iteration: 47400 Time: 1.896E+03 Max CFL: 8.377E-01 Max Diff. No.: -1.000E+00 Norm: 1.6509E-07
Iteration: 47600 Time: 1.904E+03 Max CFL: 8.377E-01 Max Diff. No.: -1.000E+00 Norm: 1.1752E-07
Iteration: 47800 Time: 1.912E+03 Max CFL: 8.377E-01 Max Diff. No.: -1.000E+00 Norm: 4.8625E-08
Iteration: 48000 Time: 1.920E+03 Max CFL: 8.377E-01 Max Diff. No.: -1.000E+00 Norm: 1.4292E-07
Writing immersed body surface data file surface.dat.
Writing solution file op.bin.
Iteration: 48200 Time: 1.928E+03 Max CFL: 8.377E-01 Max Diff. No.: -1.000E+00 Norm: 1.4748E-07
Iteration: 48400 Time: 1.936E+03 Max CFL: 8.377E-01 Max Diff. No.: -1.000E+00 Norm: 5.5165E-08
Iteration: 48600 Time: 1.944E+03 Max CFL: 8.377E-01 Max Diff. No.: -1.000E+00 Norm: 1.0027E-07
Iteration: 48800 Time: 1.952E+03 Max CFL: 8.377E-01 Max Diff. No.: -1.000E+00 Norm: 1.5391E-07
Iteration: 49000 Time: 1.960E+03 Max CFL: 8.377E-01 Max Diff. No.: -1.000E+00 Norm: 9.8802E-08
Writing immersed body surface data file surface.dat.
Writing solution file op.bin.
Iteration: 49200 Time: 1.968E+03 Max CFL: 8.377E-01 Max Diff. No.: -1.000E+00 Norm: 5.3834E-08
Iteration: 49400 Time: 1.976E+03 Max CFL: 8.377E-01 Max Diff. No.: -1.000E+00 Norm: 1.3837E-07
Iteration: 49600 Time: 1.984E+03 Max CFL: 8.377E-01 Max Diff. No.: -1.000E+00 Norm: 1.3210E-07
Iteration: 49800 Time: 1.992E+03 Max CFL: 8.377E-01 Max Diff. No.: -1.000E+00 Norm: 4.4076E-08
Iteration: 50000 Time: 2.000E+03 Max CFL: 8.377E-01 Max Diff. No.: -1.000E+00 Norm: 1.0302E-07
Writing immersed body surface data file surface.dat.
Writing solution file op.bin.
Completed time integration (Final time: 2000.000000).
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
0.0000000000000000E+00
Solver runtime (in seconds): 4.3481479369000001E+04
Total runtime (in seconds): 4.3481578939999999E+04
Deallocating arrays.
Finished.