HyPar  1.0
Finite-Difference Hyperbolic-Parabolic PDE Solver on Cartesian Grids
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros Pages
Unsteady, compressible, viscous flow around an adiabatic sphere

Location: hypar/Examples/3D/NavierStokes3D/Sphere/Unsteady_Viscous_Compressible_Adiabatic

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

Domain: The domain consists of a fine uniform grid around the sphere defined by [-2,6] X [-2,2] X [-2,2], and a stretched grid beyond this zone.

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

The following image shows the sphere:

Surface3D_Sphere.png

The following images shows the grid and the sphere:

Domain3D_Sphere1.png
Domain3D_Sphere2.png

Boundary conditions:

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

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

  • Specific heat ratio \(\gamma = 1.4\) (NavierStokes3D::gamma)
  • Freestream Mach number \(M_{\infty} = 0.8\) (NavierStokes3D::Minf)
  • Prandlt number \(Pr = 0.72\) (NavierStokes3D::Pr)
  • Reynolds number \(Re = \frac {\rho u L } {\mu} = 1000000\) (NavierStokes3D::Re) (Note: since the diameter of the sphere is 1.0, the diameter-based Reynolds number is the same as the specified Reynolds number \(Re_D = Re = 1000000\)).

Numerical Method:

Input files required:

These files are all located in: hypar/Examples/3D/NavierStokes3D/Sphere/Unsteady_Viscous_Compressible_Adiabatic/

solver.inp

begin
ndims 3
nvars 5
size 200 96 96
iproc 8 4 4
ghost 3
n_iter 25000
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.01
screen_op_iter 100
file_op_iter 250
ip_file_type binary
op_file_format binary
op_overwrite no
model navierstokes3d
immersed_body sphere.stl
end

boundary.inp

6
subsonic-inflow 0 1 0 0 -9999.0 9999.0 -9999.0 9999.0
1.0 0.8 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.8 0.0 0.0 0.7142857142857143
subsonic-ambivalent 1 -1 -9999.0 9999.0 0 0 -9999.0 9999.0
1.0 0.8 0.0 0.0 0.7142857142857143
subsonic-ambivalent 2 1 -9999.0 9999.0 -9999.0 9999.0 0 0
1.0 0.8 0.0 0.0 0.7142857142857143
subsonic-ambivalent 2 -1 -9999.0 9999.0 -9999.0 9999.0 0 0
1.0 0.8 0.0 0.0 0.7142857142857143

physics.inp : The following file specifies a Reynolds number of 1,000,000. To try other Reynolds numbers, change it here.

begin
gamma 1.4
upwinding roe
Pr 0.72
Minf 0.8
Re 1000000
end

sphere.stl : the filename "sphere.stl" must match the input for immersed_body in solver.inp.
Located at hypar/Examples/STLGeometries/sphere.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 Lz = 4.0;
double sf_x1 = 1.15;
double sf_x2 = 1.20;
double sf_y = 1.4;
double sf_z = 1.4;
int i,j,k;
double u = 0.8,
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-1));
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-1));
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]);
z[NK/8] = -Lz/2;
z[7*NK/8] = Lz/2;
dz = Lz / ((double)(6*NK/8-1));
for (k = NK/8 ; k < 7*NK/8; k++) z[k] = z[NK/8] + dz * (k-NK/8);
for (k = 7*NK/8; k < NK ; k++) z[k] = z[k-1] + sf_z * (z[k-1] - z[k-2]);
for (k = NK/8-1; k >= 0 ; k--) z[k] = z[k+1] - sf_z * (z[k+2] - z[k+1]);
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

  8 4 4

in solver.inp (i.e., 8 processors along x, 4 processors along y, and 4 processor along z). Thus, this example should be run with 128 MPI ranks (or change iproc).

After running the code, there should be 101 output files op_nnnnn.bin, since HyPar::op_overwrite is set to no 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).

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_nnnnn.dat (if HyPar::op_overwrite is "no") or surface.dat (if HyPar::op_overwrite is "yes") (in this example, the files surface_nnnnn.dat are 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 visualizations were generated using these solution files.

  • Pressure and velocity vector on a 2D x-y plane:
    Solution_3DNavStokSphere_ReD1mil_2Dflow.gif
  • Streamlines colored by velocity magnitude and sphere surface colored by temperature:
    Solution_3DNavStokSphere_ReD1mil_3Dflow.gif

Expected screen output:

HyPar - Parallel (MPI) version with 128 processes
Compiled with PETSc time integration.
Allocated simulation object(s).
Reading solver inputs from file "solver.inp".
No. of dimensions : 3
No. of variables : 5
Domain size : 200 96 96
Processes along each dimension : 8 4 4
No. of ghosts pts : 3
No. of iter. : 25000
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 : 1.000000E-02
Check for conservation : no
Screen output iterations : 100
File output iterations : 250
Initial solution file type : binary
Initial solution read mode : serial
Solution file write mode : serial
Solution file format : binary
Overwrite solution file : no
Physical model : navierstokes3d
Immersed Body : sphere.stl
Partitioning domain and allocating data arrays.
Reading array from binary file initial.inp (Serial mode).
Volume integral of the initial solution:
0: 4.6753965134569866E+04
1: 3.7403172107655490E+04
2: 0.0000000000000000E+00
3: 0.0000000000000000E+00
4: 9.8450492297650984E+04
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 subsonic-ambivalent: Along dimension 2 and face +1
Boundary subsonic-ambivalent: Along dimension 2 and face -1
6 boundary condition(s) read.
Immersed body read from sphere.stl:
Number of facets: 12320
Bounding box: [-0.5,0.5] X [-0.5,0.5] X [-0.5,0.5]
Number of grid points inside immersed body: 3619 ( 0.2%).
Number of immersed boundary points : 1988 ( 0.1%).
Immersed body simulation mode : 3d.
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_00000.bin.
Setting up time integration.
Solving in time (from 0 to 25000 iterations)
Writing immersed body surface data file surface00000.dat.
Writing solution file op_00000.bin.
Iteration: 100 Time: 1.000E+00 Max CFL: 6.607E-01 Max Diff. No.: -1.000E+00 Norm: 2.1587E-03
Iteration: 200 Time: 2.000E+00 Max CFL: 6.534E-01 Max Diff. No.: -1.000E+00 Norm: 1.8496E-03
Writing immersed body surface data file surface00001.dat.
Writing solution file op_00001.bin.
Iteration: 300 Time: 3.000E+00 Max CFL: 6.406E-01 Max Diff. No.: -1.000E+00 Norm: 1.4917E-03
Iteration: 400 Time: 4.000E+00 Max CFL: 6.355E-01 Max Diff. No.: -1.000E+00 Norm: 8.0364E-04
Iteration: 500 Time: 5.000E+00 Max CFL: 6.341E-01 Max Diff. No.: -1.000E+00 Norm: 4.1132E-04
Writing immersed body surface data file surface00002.dat.
Writing solution file op_00002.bin.
Iteration: 600 Time: 6.000E+00 Max CFL: 6.286E-01 Max Diff. No.: -1.000E+00 Norm: 3.2654E-04
Iteration: 700 Time: 7.000E+00 Max CFL: 6.233E-01 Max Diff. No.: -1.000E+00 Norm: 2.8960E-04
Writing immersed body surface data file surface00003.dat.
Writing solution file op_00003.bin.
Iteration: 800 Time: 8.000E+00 Max CFL: 6.230E-01 Max Diff. No.: -1.000E+00 Norm: 3.1168E-04
Iteration: 900 Time: 9.000E+00 Max CFL: 6.240E-01 Max Diff. No.: -1.000E+00 Norm: 3.9269E-04
Iteration: 1000 Time: 1.000E+01 Max CFL: 6.156E-01 Max Diff. No.: -1.000E+00 Norm: 4.5250E-04
Writing immersed body surface data file surface00004.dat.
Writing solution file op_00004.bin.
Iteration: 1100 Time: 1.100E+01 Max CFL: 6.161E-01 Max Diff. No.: -1.000E+00 Norm: 5.0082E-04
Iteration: 1200 Time: 1.200E+01 Max CFL: 6.183E-01 Max Diff. No.: -1.000E+00 Norm: 5.5941E-04
Writing immersed body surface data file surface00005.dat.
Writing solution file op_00005.bin.
Iteration: 1300 Time: 1.300E+01 Max CFL: 6.172E-01 Max Diff. No.: -1.000E+00 Norm: 6.4344E-04
Iteration: 1400 Time: 1.400E+01 Max CFL: 6.168E-01 Max Diff. No.: -1.000E+00 Norm: 7.4400E-04
Iteration: 1500 Time: 1.500E+01 Max CFL: 6.151E-01 Max Diff. No.: -1.000E+00 Norm: 8.2837E-04
Writing immersed body surface data file surface00006.dat.
Writing solution file op_00006.bin.
Iteration: 1600 Time: 1.600E+01 Max CFL: 6.130E-01 Max Diff. No.: -1.000E+00 Norm: 9.4428E-04
Iteration: 1700 Time: 1.700E+01 Max CFL: 6.146E-01 Max Diff. No.: -1.000E+00 Norm: 1.0079E-03
Writing immersed body surface data file surface00007.dat.
Writing solution file op_00007.bin.
Iteration: 1800 Time: 1.800E+01 Max CFL: 6.126E-01 Max Diff. No.: -1.000E+00 Norm: 1.0326E-03
Iteration: 1900 Time: 1.900E+01 Max CFL: 6.094E-01 Max Diff. No.: -1.000E+00 Norm: 9.7359E-04
Iteration: 2000 Time: 2.000E+01 Max CFL: 6.103E-01 Max Diff. No.: -1.000E+00 Norm: 8.8888E-04
Writing immersed body surface data file surface00008.dat.
Writing solution file op_00008.bin.
Iteration: 2100 Time: 2.100E+01 Max CFL: 6.111E-01 Max Diff. No.: -1.000E+00 Norm: 8.6589E-04
Iteration: 2200 Time: 2.200E+01 Max CFL: 6.105E-01 Max Diff. No.: -1.000E+00 Norm: 9.1291E-04
Writing immersed body surface data file surface00009.dat.
Writing solution file op_00009.bin.
Iteration: 2300 Time: 2.300E+01 Max CFL: 6.095E-01 Max Diff. No.: -1.000E+00 Norm: 9.9810E-04
Iteration: 2400 Time: 2.400E+01 Max CFL: 6.094E-01 Max Diff. No.: -1.000E+00 Norm: 1.0475E-03
Iteration: 2500 Time: 2.500E+01 Max CFL: 6.099E-01 Max Diff. No.: -1.000E+00 Norm: 1.0806E-03
Writing immersed body surface data file surface00010.dat.
Writing solution file op_00010.bin.
Iteration: 2600 Time: 2.600E+01 Max CFL: 6.110E-01 Max Diff. No.: -1.000E+00 Norm: 1.0470E-03
Iteration: 2700 Time: 2.700E+01 Max CFL: 6.090E-01 Max Diff. No.: -1.000E+00 Norm: 9.8278E-04
Writing immersed body surface data file surface00011.dat.
Writing solution file op_00011.bin.
Iteration: 2800 Time: 2.800E+01 Max CFL: 6.054E-01 Max Diff. No.: -1.000E+00 Norm: 9.2103E-04
Iteration: 2900 Time: 2.900E+01 Max CFL: 6.073E-01 Max Diff. No.: -1.000E+00 Norm: 9.3249E-04
Iteration: 3000 Time: 3.000E+01 Max CFL: 6.083E-01 Max Diff. No.: -1.000E+00 Norm: 9.9446E-04
Writing immersed body surface data file surface00012.dat.
Writing solution file op_00012.bin.
Iteration: 3100 Time: 3.100E+01 Max CFL: 6.068E-01 Max Diff. No.: -1.000E+00 Norm: 1.0440E-03
Iteration: 3200 Time: 3.200E+01 Max CFL: 6.067E-01 Max Diff. No.: -1.000E+00 Norm: 1.0875E-03
Writing immersed body surface data file surface00013.dat.
Writing solution file op_00013.bin.
Iteration: 3300 Time: 3.300E+01 Max CFL: 6.064E-01 Max Diff. No.: -1.000E+00 Norm: 1.1055E-03
Iteration: 3400 Time: 3.400E+01 Max CFL: 6.056E-01 Max Diff. No.: -1.000E+00 Norm: 1.0696E-03
Iteration: 3500 Time: 3.500E+01 Max CFL: 6.040E-01 Max Diff. No.: -1.000E+00 Norm: 1.0011E-03
Writing immersed body surface data file surface00014.dat.
Writing solution file op_00014.bin.
Iteration: 3600 Time: 3.600E+01 Max CFL: 6.035E-01 Max Diff. No.: -1.000E+00 Norm: 9.7561E-04
Iteration: 3700 Time: 3.700E+01 Max CFL: 6.034E-01 Max Diff. No.: -1.000E+00 Norm: 9.2655E-04
Writing immersed body surface data file surface00015.dat.
Writing solution file op_00015.bin.
Iteration: 3800 Time: 3.800E+01 Max CFL: 6.018E-01 Max Diff. No.: -1.000E+00 Norm: 8.3056E-04
Iteration: 3900 Time: 3.900E+01 Max CFL: 6.020E-01 Max Diff. No.: -1.000E+00 Norm: 8.1601E-04
Iteration: 4000 Time: 4.000E+01 Max CFL: 6.009E-01 Max Diff. No.: -1.000E+00 Norm: 8.5789E-04
Writing immersed body surface data file surface00016.dat.
Writing solution file op_00016.bin.
Iteration: 4100 Time: 4.100E+01 Max CFL: 6.017E-01 Max Diff. No.: -1.000E+00 Norm: 9.0531E-04
Iteration: 4200 Time: 4.200E+01 Max CFL: 6.018E-01 Max Diff. No.: -1.000E+00 Norm: 8.8005E-04
Writing immersed body surface data file surface00017.dat.
Writing solution file op_00017.bin.
Iteration: 4300 Time: 4.300E+01 Max CFL: 6.035E-01 Max Diff. No.: -1.000E+00 Norm: 8.8410E-04
Iteration: 4400 Time: 4.400E+01 Max CFL: 6.017E-01 Max Diff. No.: -1.000E+00 Norm: 9.6493E-04
Iteration: 4500 Time: 4.500E+01 Max CFL: 6.016E-01 Max Diff. No.: -1.000E+00 Norm: 1.0225E-03
Writing immersed body surface data file surface00018.dat.
Writing solution file op_00018.bin.
Iteration: 4600 Time: 4.600E+01 Max CFL: 6.018E-01 Max Diff. No.: -1.000E+00 Norm: 1.0776E-03
Iteration: 4700 Time: 4.700E+01 Max CFL: 6.034E-01 Max Diff. No.: -1.000E+00 Norm: 1.2139E-03
Writing immersed body surface data file surface00019.dat.
Writing solution file op_00019.bin.
Iteration: 4800 Time: 4.800E+01 Max CFL: 6.042E-01 Max Diff. No.: -1.000E+00 Norm: 1.3031E-03
Iteration: 4900 Time: 4.900E+01 Max CFL: 6.052E-01 Max Diff. No.: -1.000E+00 Norm: 1.3112E-03
Iteration: 5000 Time: 5.000E+01 Max CFL: 6.062E-01 Max Diff. No.: -1.000E+00 Norm: 1.2187E-03
Writing immersed body surface data file surface00020.dat.
Writing solution file op_00020.bin.
Iteration: 5100 Time: 5.100E+01 Max CFL: 6.034E-01 Max Diff. No.: -1.000E+00 Norm: 1.1466E-03
Iteration: 5200 Time: 5.200E+01 Max CFL: 6.030E-01 Max Diff. No.: -1.000E+00 Norm: 1.0300E-03
Writing immersed body surface data file surface00021.dat.
Writing solution file op_00021.bin.
Iteration: 5300 Time: 5.300E+01 Max CFL: 6.031E-01 Max Diff. No.: -1.000E+00 Norm: 9.4448E-04
Iteration: 5400 Time: 5.400E+01 Max CFL: 6.017E-01 Max Diff. No.: -1.000E+00 Norm: 9.2512E-04
Iteration: 5500 Time: 5.500E+01 Max CFL: 6.010E-01 Max Diff. No.: -1.000E+00 Norm: 8.7789E-04
Writing immersed body surface data file surface00022.dat.
Writing solution file op_00022.bin.
Iteration: 5600 Time: 5.600E+01 Max CFL: 6.016E-01 Max Diff. No.: -1.000E+00 Norm: 8.2529E-04
Iteration: 5700 Time: 5.700E+01 Max CFL: 6.020E-01 Max Diff. No.: -1.000E+00 Norm: 8.5311E-04
Writing immersed body surface data file surface00023.dat.
Writing solution file op_00023.bin.
Iteration: 5800 Time: 5.800E+01 Max CFL: 6.013E-01 Max Diff. No.: -1.000E+00 Norm: 9.0786E-04
Iteration: 5900 Time: 5.900E+01 Max CFL: 6.024E-01 Max Diff. No.: -1.000E+00 Norm: 9.2271E-04
Iteration: 6000 Time: 6.000E+01 Max CFL: 6.018E-01 Max Diff. No.: -1.000E+00 Norm: 8.7618E-04
Writing immersed body surface data file surface00024.dat.
Writing solution file op_00024.bin.
Iteration: 6100 Time: 6.100E+01 Max CFL: 6.001E-01 Max Diff. No.: -1.000E+00 Norm: 9.3817E-04
Iteration: 6200 Time: 6.200E+01 Max CFL: 5.982E-01 Max Diff. No.: -1.000E+00 Norm: 9.9966E-04
Writing immersed body surface data file surface00025.dat.
Writing solution file op_00025.bin.
Iteration: 6300 Time: 6.300E+01 Max CFL: 6.000E-01 Max Diff. No.: -1.000E+00 Norm: 1.0380E-03
Iteration: 6400 Time: 6.400E+01 Max CFL: 6.011E-01 Max Diff. No.: -1.000E+00 Norm: 1.0823E-03
Iteration: 6500 Time: 6.500E+01 Max CFL: 6.017E-01 Max Diff. No.: -1.000E+00 Norm: 1.1532E-03
Writing immersed body surface data file surface00026.dat.
Writing solution file op_00026.bin.
Iteration: 6600 Time: 6.600E+01 Max CFL: 6.019E-01 Max Diff. No.: -1.000E+00 Norm: 1.2094E-03
Iteration: 6700 Time: 6.700E+01 Max CFL: 6.025E-01 Max Diff. No.: -1.000E+00 Norm: 1.2183E-03
Writing immersed body surface data file surface00027.dat.
Writing solution file op_00027.bin.
Iteration: 6800 Time: 6.800E+01 Max CFL: 6.029E-01 Max Diff. No.: -1.000E+00 Norm: 1.1623E-03
Iteration: 6900 Time: 6.900E+01 Max CFL: 6.023E-01 Max Diff. No.: -1.000E+00 Norm: 1.0849E-03
Iteration: 7000 Time: 7.000E+01 Max CFL: 6.023E-01 Max Diff. No.: -1.000E+00 Norm: 1.0404E-03
Writing immersed body surface data file surface00028.dat.
Writing solution file op_00028.bin.
Iteration: 7100 Time: 7.100E+01 Max CFL: 6.019E-01 Max Diff. No.: -1.000E+00 Norm: 9.8946E-04
Iteration: 7200 Time: 7.200E+01 Max CFL: 6.020E-01 Max Diff. No.: -1.000E+00 Norm: 9.6658E-04
Writing immersed body surface data file surface00029.dat.
Writing solution file op_00029.bin.
Iteration: 7300 Time: 7.300E+01 Max CFL: 6.034E-01 Max Diff. No.: -1.000E+00 Norm: 1.0008E-03
Iteration: 7400 Time: 7.400E+01 Max CFL: 6.027E-01 Max Diff. No.: -1.000E+00 Norm: 1.0635E-03
Iteration: 7500 Time: 7.500E+01 Max CFL: 6.040E-01 Max Diff. No.: -1.000E+00 Norm: 1.1097E-03
Writing immersed body surface data file surface00030.dat.
Writing solution file op_00030.bin.
Iteration: 7600 Time: 7.600E+01 Max CFL: 6.037E-01 Max Diff. No.: -1.000E+00 Norm: 1.0747E-03
Iteration: 7700 Time: 7.700E+01 Max CFL: 6.036E-01 Max Diff. No.: -1.000E+00 Norm: 1.0728E-03
Writing immersed body surface data file surface00031.dat.
Writing solution file op_00031.bin.
Iteration: 7800 Time: 7.800E+01 Max CFL: 6.026E-01 Max Diff. No.: -1.000E+00 Norm: 1.0994E-03
Iteration: 7900 Time: 7.900E+01 Max CFL: 6.051E-01 Max Diff. No.: -1.000E+00 Norm: 1.1264E-03
Iteration: 8000 Time: 8.000E+01 Max CFL: 6.044E-01 Max Diff. No.: -1.000E+00 Norm: 1.1489E-03
Writing immersed body surface data file surface00032.dat.
Writing solution file op_00032.bin.
Iteration: 8100 Time: 8.100E+01 Max CFL: 6.034E-01 Max Diff. No.: -1.000E+00 Norm: 1.1636E-03
Iteration: 8200 Time: 8.200E+01 Max CFL: 6.052E-01 Max Diff. No.: -1.000E+00 Norm: 1.1493E-03
Writing immersed body surface data file surface00033.dat.
Writing solution file op_00033.bin.
Iteration: 8300 Time: 8.300E+01 Max CFL: 6.045E-01 Max Diff. No.: -1.000E+00 Norm: 1.0936E-03
Iteration: 8400 Time: 8.400E+01 Max CFL: 6.046E-01 Max Diff. No.: -1.000E+00 Norm: 1.0605E-03
Iteration: 8500 Time: 8.500E+01 Max CFL: 6.061E-01 Max Diff. No.: -1.000E+00 Norm: 1.0044E-03
Writing immersed body surface data file surface00034.dat.
Writing solution file op_00034.bin.
Iteration: 8600 Time: 8.600E+01 Max CFL: 6.053E-01 Max Diff. No.: -1.000E+00 Norm: 9.4573E-04
Iteration: 8700 Time: 8.700E+01 Max CFL: 6.043E-01 Max Diff. No.: -1.000E+00 Norm: 9.6967E-04
Writing immersed body surface data file surface00035.dat.
Writing solution file op_00035.bin.
Iteration: 8800 Time: 8.800E+01 Max CFL: 6.052E-01 Max Diff. No.: -1.000E+00 Norm: 1.0531E-03
Iteration: 8900 Time: 8.900E+01 Max CFL: 6.060E-01 Max Diff. No.: -1.000E+00 Norm: 1.2167E-03
Iteration: 9000 Time: 9.000E+01 Max CFL: 6.053E-01 Max Diff. No.: -1.000E+00 Norm: 1.3216E-03
Writing immersed body surface data file surface00036.dat.
Writing solution file op_00036.bin.
Iteration: 9100 Time: 9.100E+01 Max CFL: 6.062E-01 Max Diff. No.: -1.000E+00 Norm: 1.3024E-03
Iteration: 9200 Time: 9.200E+01 Max CFL: 6.074E-01 Max Diff. No.: -1.000E+00 Norm: 1.2463E-03
Writing immersed body surface data file surface00037.dat.
Writing solution file op_00037.bin.
Iteration: 9300 Time: 9.300E+01 Max CFL: 6.086E-01 Max Diff. No.: -1.000E+00 Norm: 1.1997E-03
Iteration: 9400 Time: 9.400E+01 Max CFL: 6.087E-01 Max Diff. No.: -1.000E+00 Norm: 1.2057E-03
Iteration: 9500 Time: 9.500E+01 Max CFL: 6.076E-01 Max Diff. No.: -1.000E+00 Norm: 1.2309E-03
Writing immersed body surface data file surface00038.dat.
Writing solution file op_00038.bin.
Iteration: 9600 Time: 9.600E+01 Max CFL: 6.075E-01 Max Diff. No.: -1.000E+00 Norm: 1.1424E-03
Iteration: 9700 Time: 9.700E+01 Max CFL: 6.077E-01 Max Diff. No.: -1.000E+00 Norm: 1.0772E-03
Writing immersed body surface data file surface00039.dat.
Writing solution file op_00039.bin.
Iteration: 9800 Time: 9.800E+01 Max CFL: 6.077E-01 Max Diff. No.: -1.000E+00 Norm: 1.1352E-03
Iteration: 9900 Time: 9.900E+01 Max CFL: 6.063E-01 Max Diff. No.: -1.000E+00 Norm: 1.2613E-03
Iteration: 10000 Time: 1.000E+02 Max CFL: 6.071E-01 Max Diff. No.: -1.000E+00 Norm: 1.3418E-03
Writing immersed body surface data file surface00040.dat.
Writing solution file op_00040.bin.
Iteration: 10100 Time: 1.010E+02 Max CFL: 6.071E-01 Max Diff. No.: -1.000E+00 Norm: 1.3374E-03
Iteration: 10200 Time: 1.020E+02 Max CFL: 6.067E-01 Max Diff. No.: -1.000E+00 Norm: 1.2798E-03
Writing immersed body surface data file surface00041.dat.
Writing solution file op_00041.bin.
Iteration: 10300 Time: 1.030E+02 Max CFL: 6.063E-01 Max Diff. No.: -1.000E+00 Norm: 1.1626E-03
Iteration: 10400 Time: 1.040E+02 Max CFL: 6.051E-01 Max Diff. No.: -1.000E+00 Norm: 1.0093E-03
Iteration: 10500 Time: 1.050E+02 Max CFL: 6.055E-01 Max Diff. No.: -1.000E+00 Norm: 8.2690E-04
Writing immersed body surface data file surface00042.dat.
Writing solution file op_00042.bin.
Iteration: 10600 Time: 1.060E+02 Max CFL: 6.048E-01 Max Diff. No.: -1.000E+00 Norm: 8.1764E-04
Iteration: 10700 Time: 1.070E+02 Max CFL: 6.042E-01 Max Diff. No.: -1.000E+00 Norm: 8.6238E-04
Writing immersed body surface data file surface00043.dat.
Writing solution file op_00043.bin.
Iteration: 10800 Time: 1.080E+02 Max CFL: 6.030E-01 Max Diff. No.: -1.000E+00 Norm: 8.2463E-04
Iteration: 10900 Time: 1.090E+02 Max CFL: 6.029E-01 Max Diff. No.: -1.000E+00 Norm: 8.1196E-04
Iteration: 11000 Time: 1.100E+02 Max CFL: 6.030E-01 Max Diff. No.: -1.000E+00 Norm: 8.1546E-04
Writing immersed body surface data file surface00044.dat.
Writing solution file op_00044.bin.
Iteration: 11100 Time: 1.110E+02 Max CFL: 6.027E-01 Max Diff. No.: -1.000E+00 Norm: 8.3051E-04
Iteration: 11200 Time: 1.120E+02 Max CFL: 6.022E-01 Max Diff. No.: -1.000E+00 Norm: 9.0698E-04
Writing immersed body surface data file surface00045.dat.
Writing solution file op_00045.bin.
Iteration: 11300 Time: 1.130E+02 Max CFL: 6.018E-01 Max Diff. No.: -1.000E+00 Norm: 9.8449E-04
Iteration: 11400 Time: 1.140E+02 Max CFL: 6.013E-01 Max Diff. No.: -1.000E+00 Norm: 9.8283E-04
Iteration: 11500 Time: 1.150E+02 Max CFL: 6.018E-01 Max Diff. No.: -1.000E+00 Norm: 1.0001E-03
Writing immersed body surface data file surface00046.dat.
Writing solution file op_00046.bin.
Iteration: 11600 Time: 1.160E+02 Max CFL: 6.011E-01 Max Diff. No.: -1.000E+00 Norm: 9.8092E-04
Iteration: 11700 Time: 1.170E+02 Max CFL: 6.016E-01 Max Diff. No.: -1.000E+00 Norm: 9.3357E-04
Writing immersed body surface data file surface00047.dat.
Writing solution file op_00047.bin.
Iteration: 11800 Time: 1.180E+02 Max CFL: 6.030E-01 Max Diff. No.: -1.000E+00 Norm: 9.0505E-04
Iteration: 11900 Time: 1.190E+02 Max CFL: 6.024E-01 Max Diff. No.: -1.000E+00 Norm: 9.4405E-04
Iteration: 12000 Time: 1.200E+02 Max CFL: 6.012E-01 Max Diff. No.: -1.000E+00 Norm: 9.4815E-04
Writing immersed body surface data file surface00048.dat.
Writing solution file op_00048.bin.
Iteration: 12100 Time: 1.210E+02 Max CFL: 6.001E-01 Max Diff. No.: -1.000E+00 Norm: 1.0411E-03
Iteration: 12200 Time: 1.220E+02 Max CFL: 6.006E-01 Max Diff. No.: -1.000E+00 Norm: 1.1021E-03
Writing immersed body surface data file surface00049.dat.
Writing solution file op_00049.bin.
Iteration: 12300 Time: 1.230E+02 Max CFL: 6.015E-01 Max Diff. No.: -1.000E+00 Norm: 1.1723E-03
Iteration: 12400 Time: 1.240E+02 Max CFL: 6.011E-01 Max Diff. No.: -1.000E+00 Norm: 1.2203E-03
Iteration: 12500 Time: 1.250E+02 Max CFL: 6.016E-01 Max Diff. No.: -1.000E+00 Norm: 1.1761E-03
Writing immersed body surface data file surface00050.dat.
Writing solution file op_00050.bin.
Iteration: 12600 Time: 1.260E+02 Max CFL: 6.024E-01 Max Diff. No.: -1.000E+00 Norm: 1.1207E-03
Iteration: 12700 Time: 1.270E+02 Max CFL: 6.032E-01 Max Diff. No.: -1.000E+00 Norm: 1.0327E-03
Writing immersed body surface data file surface00051.dat.
Writing solution file op_00051.bin.
Iteration: 12800 Time: 1.280E+02 Max CFL: 6.047E-01 Max Diff. No.: -1.000E+00 Norm: 1.0277E-03
Iteration: 12900 Time: 1.290E+02 Max CFL: 6.041E-01 Max Diff. No.: -1.000E+00 Norm: 1.0808E-03
Iteration: 13000 Time: 1.300E+02 Max CFL: 6.038E-01 Max Diff. No.: -1.000E+00 Norm: 1.1248E-03
Writing immersed body surface data file surface00052.dat.
Writing solution file op_00052.bin.
Iteration: 13100 Time: 1.310E+02 Max CFL: 1.546E+00 Max Diff. No.: -1.000E+00 Norm: 1.0776E-03
Iteration: 13200 Time: 1.320E+02 Max CFL: 1.477E+00 Max Diff. No.: -1.000E+00 Norm: 9.8623E-04
Writing immersed body surface data file surface00053.dat.
Writing solution file op_00053.bin.
Iteration: 13300 Time: 1.330E+02 Max CFL: 1.484E+00 Max Diff. No.: -1.000E+00 Norm: 9.5300E-04
Iteration: 13400 Time: 1.340E+02 Max CFL: 8.324E-01 Max Diff. No.: -1.000E+00 Norm: 9.8610E-04
Iteration: 13500 Time: 1.350E+02 Max CFL: 6.018E-01 Max Diff. No.: -1.000E+00 Norm: 1.0066E-03
Writing immersed body surface data file surface00054.dat.
Writing solution file op_00054.bin.
Iteration: 13600 Time: 1.360E+02 Max CFL: 6.009E-01 Max Diff. No.: -1.000E+00 Norm: 9.9106E-04
Iteration: 13700 Time: 1.370E+02 Max CFL: 6.008E-01 Max Diff. No.: -1.000E+00 Norm: 1.0551E-03
Writing immersed body surface data file surface00055.dat.
Writing solution file op_00055.bin.
Iteration: 13800 Time: 1.380E+02 Max CFL: 6.017E-01 Max Diff. No.: -1.000E+00 Norm: 1.0601E-03
Iteration: 13900 Time: 1.390E+02 Max CFL: 6.057E-01 Max Diff. No.: -1.000E+00 Norm: 1.0130E-03
Iteration: 14000 Time: 1.400E+02 Max CFL: 6.067E-01 Max Diff. No.: -1.000E+00 Norm: 1.0080E-03
Writing immersed body surface data file surface00056.dat.
Writing solution file op_00056.bin.
Iteration: 14100 Time: 1.410E+02 Max CFL: 6.062E-01 Max Diff. No.: -1.000E+00 Norm: 1.0214E-03
Iteration: 14200 Time: 1.420E+02 Max CFL: 6.064E-01 Max Diff. No.: -1.000E+00 Norm: 1.0908E-03
Writing immersed body surface data file surface00057.dat.
Writing solution file op_00057.bin.
Iteration: 14300 Time: 1.430E+02 Max CFL: 6.064E-01 Max Diff. No.: -1.000E+00 Norm: 1.1500E-03
Iteration: 14400 Time: 1.440E+02 Max CFL: 6.068E-01 Max Diff. No.: -1.000E+00 Norm: 1.2072E-03
Iteration: 14500 Time: 1.450E+02 Max CFL: 6.061E-01 Max Diff. No.: -1.000E+00 Norm: 1.2617E-03
Writing immersed body surface data file surface00058.dat.
Writing solution file op_00058.bin.
Iteration: 14600 Time: 1.460E+02 Max CFL: 6.065E-01 Max Diff. No.: -1.000E+00 Norm: 1.2841E-03
Iteration: 14700 Time: 1.470E+02 Max CFL: 6.045E-01 Max Diff. No.: -1.000E+00 Norm: 1.2285E-03
Writing immersed body surface data file surface00059.dat.
Writing solution file op_00059.bin.
Iteration: 14800 Time: 1.480E+02 Max CFL: 6.053E-01 Max Diff. No.: -1.000E+00 Norm: 1.1534E-03
Iteration: 14900 Time: 1.490E+02 Max CFL: 6.055E-01 Max Diff. No.: -1.000E+00 Norm: 1.1809E-03
Iteration: 15000 Time: 1.500E+02 Max CFL: 6.049E-01 Max Diff. No.: -1.000E+00 Norm: 1.1963E-03
Writing immersed body surface data file surface00060.dat.
Writing solution file op_00060.bin.
Iteration: 15100 Time: 1.510E+02 Max CFL: 6.064E-01 Max Diff. No.: -1.000E+00 Norm: 1.2028E-03
Iteration: 15200 Time: 1.520E+02 Max CFL: 6.058E-01 Max Diff. No.: -1.000E+00 Norm: 1.2102E-03
Writing immersed body surface data file surface00061.dat.
Writing solution file op_00061.bin.
Iteration: 15300 Time: 1.530E+02 Max CFL: 6.056E-01 Max Diff. No.: -1.000E+00 Norm: 1.1441E-03
Iteration: 15400 Time: 1.540E+02 Max CFL: 6.068E-01 Max Diff. No.: -1.000E+00 Norm: 1.1099E-03
Iteration: 15500 Time: 1.550E+02 Max CFL: 6.075E-01 Max Diff. No.: -1.000E+00 Norm: 1.0958E-03
Writing immersed body surface data file surface00062.dat.
Writing solution file op_00062.bin.
Iteration: 15600 Time: 1.560E+02 Max CFL: 6.076E-01 Max Diff. No.: -1.000E+00 Norm: 1.0242E-03
Iteration: 15700 Time: 1.570E+02 Max CFL: 6.068E-01 Max Diff. No.: -1.000E+00 Norm: 1.0342E-03
Writing immersed body surface data file surface00063.dat.
Writing solution file op_00063.bin.
Iteration: 15800 Time: 1.580E+02 Max CFL: 6.080E-01 Max Diff. No.: -1.000E+00 Norm: 1.0191E-03
Iteration: 15900 Time: 1.590E+02 Max CFL: 1.284E+00 Max Diff. No.: -1.000E+00 Norm: 1.0093E-03
Iteration: 16000 Time: 1.600E+02 Max CFL: 6.068E-01 Max Diff. No.: -1.000E+00 Norm: 1.0188E-03
Writing immersed body surface data file surface00064.dat.
Writing solution file op_00064.bin.
Iteration: 16100 Time: 1.610E+02 Max CFL: 6.051E-01 Max Diff. No.: -1.000E+00 Norm: 9.8276E-04
Iteration: 16200 Time: 1.620E+02 Max CFL: 6.059E-01 Max Diff. No.: -1.000E+00 Norm: 9.2220E-04
Writing immersed body surface data file surface00065.dat.
Writing solution file op_00065.bin.
Iteration: 16300 Time: 1.630E+02 Max CFL: 6.054E-01 Max Diff. No.: -1.000E+00 Norm: 8.9649E-04
Iteration: 16400 Time: 1.640E+02 Max CFL: 6.048E-01 Max Diff. No.: -1.000E+00 Norm: 8.9448E-04
Iteration: 16500 Time: 1.650E+02 Max CFL: 6.044E-01 Max Diff. No.: -1.000E+00 Norm: 9.1285E-04
Writing immersed body surface data file surface00066.dat.
Writing solution file op_00066.bin.
Iteration: 16600 Time: 1.660E+02 Max CFL: 6.043E-01 Max Diff. No.: -1.000E+00 Norm: 9.6216E-04
Iteration: 16700 Time: 1.670E+02 Max CFL: 6.038E-01 Max Diff. No.: -1.000E+00 Norm: 9.2653E-04
Writing immersed body surface data file surface00067.dat.
Writing solution file op_00067.bin.
Iteration: 16800 Time: 1.680E+02 Max CFL: 6.037E-01 Max Diff. No.: -1.000E+00 Norm: 8.9057E-04
Iteration: 16900 Time: 1.690E+02 Max CFL: 6.027E-01 Max Diff. No.: -1.000E+00 Norm: 8.9011E-04
Iteration: 17000 Time: 1.700E+02 Max CFL: 6.023E-01 Max Diff. No.: -1.000E+00 Norm: 9.0964E-04
Writing immersed body surface data file surface00068.dat.
Writing solution file op_00068.bin.
Iteration: 17100 Time: 1.710E+02 Max CFL: 6.022E-01 Max Diff. No.: -1.000E+00 Norm: 9.6842E-04
Iteration: 17200 Time: 1.720E+02 Max CFL: 6.029E-01 Max Diff. No.: -1.000E+00 Norm: 1.0174E-03
Writing immersed body surface data file surface00069.dat.
Writing solution file op_00069.bin.
Iteration: 17300 Time: 1.730E+02 Max CFL: 6.032E-01 Max Diff. No.: -1.000E+00 Norm: 9.6810E-04
Iteration: 17400 Time: 1.740E+02 Max CFL: 6.033E-01 Max Diff. No.: -1.000E+00 Norm: 8.8410E-04
Iteration: 17500 Time: 1.750E+02 Max CFL: 6.034E-01 Max Diff. No.: -1.000E+00 Norm: 9.2965E-04
Writing immersed body surface data file surface00070.dat.
Writing solution file op_00070.bin.
Iteration: 17600 Time: 1.760E+02 Max CFL: 6.043E-01 Max Diff. No.: -1.000E+00 Norm: 9.9215E-04
Iteration: 17700 Time: 1.770E+02 Max CFL: 6.036E-01 Max Diff. No.: -1.000E+00 Norm: 9.9665E-04
Writing immersed body surface data file surface00071.dat.
Writing solution file op_00071.bin.
Iteration: 17800 Time: 1.780E+02 Max CFL: 6.043E-01 Max Diff. No.: -1.000E+00 Norm: 9.4065E-04
Iteration: 17900 Time: 1.790E+02 Max CFL: 6.037E-01 Max Diff. No.: -1.000E+00 Norm: 1.0291E-03
Iteration: 18000 Time: 1.800E+02 Max CFL: 6.038E-01 Max Diff. No.: -1.000E+00 Norm: 1.1830E-03
Writing immersed body surface data file surface00072.dat.
Writing solution file op_00072.bin.
Iteration: 18100 Time: 1.810E+02 Max CFL: 6.049E-01 Max Diff. No.: -1.000E+00 Norm: 1.3203E-03
Iteration: 18200 Time: 1.820E+02 Max CFL: 6.046E-01 Max Diff. No.: -1.000E+00 Norm: 1.4205E-03
Writing immersed body surface data file surface00073.dat.
Writing solution file op_00073.bin.
Iteration: 18300 Time: 1.830E+02 Max CFL: 6.040E-01 Max Diff. No.: -1.000E+00 Norm: 1.4025E-03
Iteration: 18400 Time: 1.840E+02 Max CFL: 6.052E-01 Max Diff. No.: -1.000E+00 Norm: 1.2797E-03
Iteration: 18500 Time: 1.850E+02 Max CFL: 6.056E-01 Max Diff. No.: -1.000E+00 Norm: 1.1985E-03
Writing immersed body surface data file surface00074.dat.
Writing solution file op_00074.bin.
Iteration: 18600 Time: 1.860E+02 Max CFL: 6.065E-01 Max Diff. No.: -1.000E+00 Norm: 1.2232E-03
Iteration: 18700 Time: 1.870E+02 Max CFL: 6.070E-01 Max Diff. No.: -1.000E+00 Norm: 1.1732E-03
Writing immersed body surface data file surface00075.dat.
Writing solution file op_00075.bin.
Iteration: 18800 Time: 1.880E+02 Max CFL: 6.066E-01 Max Diff. No.: -1.000E+00 Norm: 1.1418E-03
Iteration: 18900 Time: 1.890E+02 Max CFL: 6.074E-01 Max Diff. No.: -1.000E+00 Norm: 1.1873E-03
Iteration: 19000 Time: 1.900E+02 Max CFL: 6.079E-01 Max Diff. No.: -1.000E+00 Norm: 1.1947E-03
Writing immersed body surface data file surface00076.dat.
Writing solution file op_00076.bin.
Iteration: 19100 Time: 1.910E+02 Max CFL: 6.084E-01 Max Diff. No.: -1.000E+00 Norm: 1.2128E-03
Iteration: 19200 Time: 1.920E+02 Max CFL: 6.090E-01 Max Diff. No.: -1.000E+00 Norm: 1.1832E-03
Writing immersed body surface data file surface00077.dat.
Writing solution file op_00077.bin.
Iteration: 19300 Time: 1.930E+02 Max CFL: 6.102E-01 Max Diff. No.: -1.000E+00 Norm: 1.1854E-03
Iteration: 19400 Time: 1.940E+02 Max CFL: 6.106E-01 Max Diff. No.: -1.000E+00 Norm: 1.2050E-03
Iteration: 19500 Time: 1.950E+02 Max CFL: 6.107E-01 Max Diff. No.: -1.000E+00 Norm: 1.2275E-03
Writing immersed body surface data file surface00078.dat.
Writing solution file op_00078.bin.
Iteration: 19600 Time: 1.960E+02 Max CFL: 6.110E-01 Max Diff. No.: -1.000E+00 Norm: 1.3046E-03
Iteration: 19700 Time: 1.970E+02 Max CFL: 6.103E-01 Max Diff. No.: -1.000E+00 Norm: 1.3311E-03
Writing immersed body surface data file surface00079.dat.
Writing solution file op_00079.bin.
Iteration: 19800 Time: 1.980E+02 Max CFL: 6.108E-01 Max Diff. No.: -1.000E+00 Norm: 1.3359E-03
Iteration: 19900 Time: 1.990E+02 Max CFL: 6.107E-01 Max Diff. No.: -1.000E+00 Norm: 1.2573E-03
Iteration: 20000 Time: 2.000E+02 Max CFL: 6.096E-01 Max Diff. No.: -1.000E+00 Norm: 1.2022E-03
Writing immersed body surface data file surface00080.dat.
Writing solution file op_00080.bin.
Iteration: 20100 Time: 2.010E+02 Max CFL: 6.101E-01 Max Diff. No.: -1.000E+00 Norm: 1.1989E-03
Iteration: 20200 Time: 2.020E+02 Max CFL: 6.091E-01 Max Diff. No.: -1.000E+00 Norm: 1.1973E-03
Writing immersed body surface data file surface00081.dat.
Writing solution file op_00081.bin.
Iteration: 20300 Time: 2.030E+02 Max CFL: 6.087E-01 Max Diff. No.: -1.000E+00 Norm: 1.2371E-03
Iteration: 20400 Time: 2.040E+02 Max CFL: 6.077E-01 Max Diff. No.: -1.000E+00 Norm: 1.1029E-03
Iteration: 20500 Time: 2.050E+02 Max CFL: 6.055E-01 Max Diff. No.: -1.000E+00 Norm: 9.9278E-04
Writing immersed body surface data file surface00082.dat.
Writing solution file op_00082.bin.
Iteration: 20600 Time: 2.060E+02 Max CFL: 6.050E-01 Max Diff. No.: -1.000E+00 Norm: 9.0469E-04
Iteration: 20700 Time: 2.070E+02 Max CFL: 6.030E-01 Max Diff. No.: -1.000E+00 Norm: 8.1020E-04
Writing immersed body surface data file surface00083.dat.
Writing solution file op_00083.bin.
Iteration: 20800 Time: 2.080E+02 Max CFL: 6.022E-01 Max Diff. No.: -1.000E+00 Norm: 8.2543E-04
Iteration: 20900 Time: 2.090E+02 Max CFL: 6.016E-01 Max Diff. No.: -1.000E+00 Norm: 8.3492E-04
Iteration: 21000 Time: 2.100E+02 Max CFL: 6.005E-01 Max Diff. No.: -1.000E+00 Norm: 8.1449E-04
Writing immersed body surface data file surface00084.dat.
Writing solution file op_00084.bin.
Iteration: 21100 Time: 2.110E+02 Max CFL: 6.001E-01 Max Diff. No.: -1.000E+00 Norm: 7.9149E-04
Iteration: 21200 Time: 2.120E+02 Max CFL: 6.000E-01 Max Diff. No.: -1.000E+00 Norm: 7.6865E-04
Writing immersed body surface data file surface00085.dat.
Writing solution file op_00085.bin.
Iteration: 21300 Time: 2.130E+02 Max CFL: 5.999E-01 Max Diff. No.: -1.000E+00 Norm: 8.8928E-04
Iteration: 21400 Time: 2.140E+02 Max CFL: 6.002E-01 Max Diff. No.: -1.000E+00 Norm: 9.9502E-04
Iteration: 21500 Time: 2.150E+02 Max CFL: 6.007E-01 Max Diff. No.: -1.000E+00 Norm: 1.0482E-03
Writing immersed body surface data file surface00086.dat.
Writing solution file op_00086.bin.
Iteration: 21600 Time: 2.160E+02 Max CFL: 6.010E-01 Max Diff. No.: -1.000E+00 Norm: 1.0408E-03
Iteration: 21700 Time: 2.170E+02 Max CFL: 6.012E-01 Max Diff. No.: -1.000E+00 Norm: 1.0477E-03
Writing immersed body surface data file surface00087.dat.
Writing solution file op_00087.bin.
Iteration: 21800 Time: 2.180E+02 Max CFL: 6.000E-01 Max Diff. No.: -1.000E+00 Norm: 1.0202E-03
Iteration: 21900 Time: 2.190E+02 Max CFL: 6.007E-01 Max Diff. No.: -1.000E+00 Norm: 9.8021E-04
Iteration: 22000 Time: 2.200E+02 Max CFL: 6.010E-01 Max Diff. No.: -1.000E+00 Norm: 9.8226E-04
Writing immersed body surface data file surface00088.dat.
Writing solution file op_00088.bin.
Iteration: 22100 Time: 2.210E+02 Max CFL: 6.020E-01 Max Diff. No.: -1.000E+00 Norm: 1.0081E-03
Iteration: 22200 Time: 2.220E+02 Max CFL: 6.026E-01 Max Diff. No.: -1.000E+00 Norm: 1.0253E-03
Writing immersed body surface data file surface00089.dat.
Writing solution file op_00089.bin.
Iteration: 22300 Time: 2.230E+02 Max CFL: 6.011E-01 Max Diff. No.: -1.000E+00 Norm: 1.1114E-03
Iteration: 22400 Time: 2.240E+02 Max CFL: 6.029E-01 Max Diff. No.: -1.000E+00 Norm: 1.1458E-03
Iteration: 22500 Time: 2.250E+02 Max CFL: 6.021E-01 Max Diff. No.: -1.000E+00 Norm: 1.1563E-03
Writing immersed body surface data file surface00090.dat.
Writing solution file op_00090.bin.
Iteration: 22600 Time: 2.260E+02 Max CFL: 6.021E-01 Max Diff. No.: -1.000E+00 Norm: 1.1697E-03
Iteration: 22700 Time: 2.270E+02 Max CFL: 6.031E-01 Max Diff. No.: -1.000E+00 Norm: 1.1329E-03
Writing immersed body surface data file surface00091.dat.
Writing solution file op_00091.bin.
Iteration: 22800 Time: 2.280E+02 Max CFL: 6.037E-01 Max Diff. No.: -1.000E+00 Norm: 1.1228E-03
Iteration: 22900 Time: 2.290E+02 Max CFL: 6.044E-01 Max Diff. No.: -1.000E+00 Norm: 1.0832E-03
Iteration: 23000 Time: 2.300E+02 Max CFL: 6.047E-01 Max Diff. No.: -1.000E+00 Norm: 1.0486E-03
Writing immersed body surface data file surface00092.dat.
Writing solution file op_00092.bin.
Iteration: 23100 Time: 2.310E+02 Max CFL: 6.045E-01 Max Diff. No.: -1.000E+00 Norm: 1.0785E-03
Iteration: 23200 Time: 2.320E+02 Max CFL: 6.061E-01 Max Diff. No.: -1.000E+00 Norm: 1.1937E-03
Writing immersed body surface data file surface00093.dat.
Writing solution file op_00093.bin.
Iteration: 23300 Time: 2.330E+02 Max CFL: 6.059E-01 Max Diff. No.: -1.000E+00 Norm: 1.2537E-03
Iteration: 23400 Time: 2.340E+02 Max CFL: 6.068E-01 Max Diff. No.: -1.000E+00 Norm: 1.2414E-03
Iteration: 23500 Time: 2.350E+02 Max CFL: 6.066E-01 Max Diff. No.: -1.000E+00 Norm: 1.2200E-03
Writing immersed body surface data file surface00094.dat.
Writing solution file op_00094.bin.
Iteration: 23600 Time: 2.360E+02 Max CFL: 6.056E-01 Max Diff. No.: -1.000E+00 Norm: 1.1886E-03
Iteration: 23700 Time: 2.370E+02 Max CFL: 6.061E-01 Max Diff. No.: -1.000E+00 Norm: 1.1737E-03
Writing immersed body surface data file surface00095.dat.
Writing solution file op_00095.bin.
Iteration: 23800 Time: 2.380E+02 Max CFL: 6.076E-01 Max Diff. No.: -1.000E+00 Norm: 1.2204E-03
Iteration: 23900 Time: 2.390E+02 Max CFL: 6.088E-01 Max Diff. No.: -1.000E+00 Norm: 1.1511E-03
Iteration: 24000 Time: 2.400E+02 Max CFL: 6.087E-01 Max Diff. No.: -1.000E+00 Norm: 1.0819E-03
Writing immersed body surface data file surface00096.dat.
Writing solution file op_00096.bin.
Iteration: 24100 Time: 2.410E+02 Max CFL: 6.074E-01 Max Diff. No.: -1.000E+00 Norm: 1.0918E-03
Iteration: 24200 Time: 2.420E+02 Max CFL: 6.064E-01 Max Diff. No.: -1.000E+00 Norm: 1.2059E-03
Writing immersed body surface data file surface00097.dat.
Writing solution file op_00097.bin.
Iteration: 24300 Time: 2.430E+02 Max CFL: 6.069E-01 Max Diff. No.: -1.000E+00 Norm: 1.2782E-03
Iteration: 24400 Time: 2.440E+02 Max CFL: 6.078E-01 Max Diff. No.: -1.000E+00 Norm: 1.2969E-03
Iteration: 24500 Time: 2.450E+02 Max CFL: 6.060E-01 Max Diff. No.: -1.000E+00 Norm: 1.1951E-03
Writing immersed body surface data file surface00098.dat.
Writing solution file op_00098.bin.
Iteration: 24600 Time: 2.460E+02 Max CFL: 6.053E-01 Max Diff. No.: -1.000E+00 Norm: 1.1339E-03
Iteration: 24700 Time: 2.470E+02 Max CFL: 6.037E-01 Max Diff. No.: -1.000E+00 Norm: 1.0123E-03
Writing immersed body surface data file surface00099.dat.
Writing solution file op_00099.bin.
Iteration: 24800 Time: 2.480E+02 Max CFL: 6.037E-01 Max Diff. No.: -1.000E+00 Norm: 9.0485E-04
Iteration: 24900 Time: 2.490E+02 Max CFL: 6.033E-01 Max Diff. No.: -1.000E+00 Norm: 9.0238E-04
Iteration: 25000 Time: 2.500E+02 Max CFL: 6.030E-01 Max Diff. No.: -1.000E+00 Norm: 9.7567E-04
Writing immersed body surface data file surface00100.dat.
Writing solution file op_00100.bin.
Completed time integration (Final time: 250.000000).
Solver runtime (in seconds): 9.8006815889999998E+03
Total runtime (in seconds): 9.8018351559999992E+03
Deallocating arrays.
Finished.