Introduction
Introduction Statistics Contact Development Disclaimer Help
tadded output file time parameter - ns2dfd - 2D finite difference Navier Stokes…
git clone git://src.adamsgaard.dk/ns2dfd
Log
Files
Refs
LICENSE
---
commit 3eba95c436a2e58872c5a56f133fa08df65a214a
parent 08a48407f92cd5c0c603c8e5207ab92b8a5df5a5
Author: Anders Damsgaard <[email protected]>
Date: Sun, 2 Mar 2014 21:12:11 +0100
added output file time parameter
Diffstat:
M ns2dfd.py | 17 ++++++++++++++---
M src/file_io.c | 6 ++++--
M src/file_io.h | 4 ++--
M src/main.c | 34 +++++++++++++++++++++++------…
M src/utility.c | 13 ++++++++++++-
5 files changed, 57 insertions(+), 17 deletions(-)
---
diff --git a/ns2dfd.py b/ns2dfd.py
t@@ -16,6 +16,7 @@ class fluid:
self.init_grid()
self.current_time()
self.end_time()
+ self.file_time()
self.safety_factor()
self.max_iterations()
self.tolerance_criteria()
t@@ -61,6 +62,14 @@ class fluid:
'''
self.t_end = numpy.asarray(t_end)
+ def file_time(self, t_file = 0.1):
+ '''
+ Set the simulation output file interval.
+ :param t_file: The time when to stop the simulation.
+ :type t_file: float
+ '''
+ self.t_file = numpy.asarray(t_file)
+
def safety_factor(self, tau = 0.5):
'''
Define the safety factor for the time step size control. Default value…
t@@ -150,9 +159,10 @@ class fluid:
print('Input file: ' + targetfile)
fh = open(targetfile, 'rb')
- self.t = numpy.fromfile(fh, dtype=numpy.float64, count=1)
- self.t_end = numpy.fromfile(fh, dtype=numpy.float64, count=1)
- self.tau = numpy.fromfile(fh, dtype=numpy.float64, count=1)
+ self.t = numpy.fromfile(fh, dtype=numpy.float64, count=1)
+ self.t_end = numpy.fromfile(fh, dtype=numpy.float64, count=1)
+ self.t_file = numpy.fromfile(fh, dtype=numpy.float64, count=1)
+ self.tau = numpy.fromfile(fh, dtype=numpy.float64, count=1)
self.itermax = numpy.fromfile(fh, dtype=numpy.int32, count=1)
self.epsilon = numpy.fromfile(fh, dtype=numpy.float64, count=1)
t@@ -209,6 +219,7 @@ class fluid:
fh.write(self.t.astype(numpy.float64))
fh.write(self.t_end.astype(numpy.float64))
+ fh.write(self.t_file.astype(numpy.float64))
fh.write(self.tau.astype(numpy.float64))
fh.write(self.itermax.astype(numpy.int32))
diff --git a/src/file_io.c b/src/file_io.c
t@@ -3,7 +3,7 @@
/* Read the variable values from a file on disk */
int read_file(char *path,
- double *t, double *t_end, double *tau, int *itermax,
+ double *t, double *t_end, double *t_file, double *tau, int *itermax,
double *epsilon, double *omega, double *gamma,
double *gx, double *gy, double *re,
int *w_left, int *w_right, int *w_top, int *w_bottom,
t@@ -21,6 +21,7 @@ int read_file(char *path,
fread(t, sizeof(double), 1, fp);
fread(t_end, sizeof(double), 1, fp);
+ fread(t_file, sizeof(double), 1, fp);
fread(tau, sizeof(double), 1, fp);
fread(itermax, sizeof(int), 1, fp);
t@@ -72,7 +73,7 @@ int read_file(char *path,
/* Write the variable values to a file on disk */
int write_file(char *path,
- double *t, double *t_end, double *tau, int *itermax,
+ double *t, double *t_end, double *t_file, double *tau, int *itermax,
double *epsilon, double *omega, double *gamma,
double *gx, double *gy, double *re,
int *w_left, int *w_right, int *w_top, int *w_bottom,
t@@ -90,6 +91,7 @@ int write_file(char *path,
fwrite(t, sizeof(double), 1, fp);
fwrite(t_end, sizeof(double), 1, fp);
+ fwrite(t_file, sizeof(double), 1, fp);
fwrite(tau, sizeof(double), 1, fp);
fwrite(itermax, sizeof(int), 1, fp);
diff --git a/src/file_io.h b/src/file_io.h
t@@ -2,7 +2,7 @@
#define FILE_IO_H_
int read_file(char *path,
- double *t, double *t_end, double *tau, int *itermax,
+ double *t, double *t_end, double *t_file, double *tau, int *itermax,
double *epsilon, double *omega, double *gamma,
double *gx, double *gy, double *re,
int *w_left, int *w_right, int *w_top, int *w_bottom,
t@@ -11,7 +11,7 @@ int read_file(char *path,
double ***P, double ***U, double ***V);
int write_file(char *path,
- double *t, double *t_end, double *tau, int *itermax,
+ double *t, double *t_end, double *t_file, double *tau, int *itermax,
double *epsilon, double *omega, double *gamma,
double *gx, double *gy, double *re,
int *w_left, int *w_right, int *w_top, int *w_bottom,
diff --git a/src/main.c b/src/main.c
t@@ -2,6 +2,7 @@
#include <stdlib.h>
#include <ctype.h>
#include <unistd.h>
+#include <string.h>
#include "file_io.h"
#include "utility.h"
t@@ -9,7 +10,7 @@
int main(int argc, char** argv)
{
- double t, t_end, tau;
+ double t, t_file, t_end, tau;
int itermax;
double epsilon, omega, gamma;
double gx, gy, re;
t@@ -20,6 +21,10 @@ int main(int argc, char** argv)
double dt;
long n = 0;
+ int nfile = 0;
+ double t_file_elapsed = 0.0;
+ char filename[50];
+ char *simulation_id;
int c;
while ((c = getopt(argc, argv, "hv")) != -1)
t@@ -51,7 +56,7 @@ int main(int argc, char** argv)
}
if (optind == argc - 1) {
- read_file(argv[optind], &t, &t_end, &tau, &itermax,
+ read_file(argv[optind], &t, &t_end, &t_file, &tau, &itermax,
&epsilon, &omega, &gamma,
&gx, &gy, &re, &w_left, &w_right, &w_top, &w_bottom,
&dx, &dy, &nx, &ny, &P, &U, &V);
t@@ -60,20 +65,31 @@ int main(int argc, char** argv)
return 1;
}
- /*write_file("unnamed2.dat", &t, &t_end, &tau, &itermax,
- &epsilon, &omega, &gamma,
- &gx, &gy, &re, &w_left, &w_right, &w_top, &w_bottom,
- &dx, &dy, &nx, &ny, &P, &U, &V);*/
+ simulation_id = argv[optind];
+ char* dot = strchr(simulation_id, '.');
+ dot = '\0';
- /*while (t < t_end) {*/
+ printf("%s\n", simulation_id);
+
+ while (t < t_end) {
dt = select_time_step(tau, re, dx, dy, nx, ny, U, V);
printf("dt = %f\n", dt);
- /*t += dt;
+ if (t_file_elapsed >= t_file) {
+ write_file("unnamed2.dat", &t, &t_end, &t_file, &tau, &itermax,
+ &epsilon, &omega, &gamma,
+ &gx, &gy, &re, &w_left, &w_right, &w_top, &w_bottom,
+ &dx, &dy, &nx, &ny, &P, &U, &V);
+ t_file_elapsed = 0.0;
+ }
+
+ t += dt;
n++;
- }*/
+ t_file_elapsed += dt;
+ break;
+ }
free_memory(P, U, V, nx);
return 0;
diff --git a/src/utility.c b/src/utility.c
t@@ -68,6 +68,7 @@ double select_time_step(double tau, double re, double dx, do…
{
double t_diff, t_adv_u, t_adv_v;
double u_max, v_max;
+ double dt;
u_max = max_abs_value(U, nx, ny);
v_max = max_abs_value(V, nx, ny);
t@@ -75,5 +76,15 @@ double select_time_step(double tau, double re, double dx, d…
t_adv_u = dx/(u_max + 1.0e-12);
t_adv_v = dx/(v_max + 1.0e-12);
- return fmin(t_diff, fmin(t_adv_u, t_adv_v));
+ dt = fmin(t_diff, fmin(t_adv_u, t_adv_v));
+
+ if (dt < 1.0e-12) {
+ fprintf(stderr, "select_time_step: Error, the time step is too small: "
+ "%f s.\n"
+ "tau = %f, Re = %f, u_max = %f, v_max = %f\n"
+ "t_diff = %f s, t_adv_u = %f s, t_adv_v = %f s\n",
+ dt, tau, re, u_max, v_max, t_diff, t_adv_u, t_adv_v);
+ exit(1);
+ }
+ return dt;
}
You are viewing proxied material from mx1.adamsgaard.dk. The copyright of proxied material belongs to its original authors. Any comments or complaints in relation to proxied material should be directed to the original authors of the content concerned. Please see the disclaimer for more details.