tDebugging memory - ns2dfd - 2D finite difference Navier Stokes solver for flui… | |
git clone git://src.adamsgaard.dk/ns2dfd | |
Log | |
Files | |
Refs | |
LICENSE | |
--- | |
commit 428ee066696d71e182801d30915c3ee16ea06e9b | |
parent 21beace0920e40917cfb89632dc3bb0387b2370a | |
Author: Anders Damsgaard <[email protected]> | |
Date: Sun, 2 Mar 2014 11:33:22 +0100 | |
Debugging memory | |
Diffstat: | |
M src/Makefile | 13 +++++++++++++ | |
M src/file_io.c | 64 +++++++++++++++++++++++++++++… | |
M src/main.c | 27 +++++++++++++++++++++++++++ | |
M src/utility.c | 49 +++++++++++++++++++----------… | |
M src/utility.h | 8 +++++--- | |
5 files changed, 140 insertions(+), 21 deletions(-) | |
--- | |
diff --git a/src/Makefile b/src/Makefile | |
t@@ -0,0 +1,13 @@ | |
+CFLAGS=-Wall -pedantic -g | |
+ | |
+BIN=../ns2dfd | |
+ | |
+$(BIN): main.o file_io.o utility.o | |
+ $(CC) $(CFLAGS) $^ -o $@ | |
+ | |
+main.o: file_io.h utility.h | |
+file_io.o: utility.h | |
+ | |
+clean: | |
+ $(RM) *.o | |
+ $(RM) $(BIN) | |
diff --git a/src/file_io.c b/src/file_io.c | |
t@@ -0,0 +1,64 @@ | |
+#include <stdio.h> | |
+#include "utility.h" | |
+ | |
+/* Read the variable values from a file on disk */ | |
+int read_file(char *path, | |
+ double *t, double *t_end, 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, | |
+ double *dx, double *dy, | |
+ int *nx, int *ny, | |
+ double **P, double **U, double **V) | |
+{ | |
+ int i, j; | |
+ double tmp; | |
+ FILE *fp; | |
+ if ((fp = fopen(path, "rb")) == NULL) { | |
+ fprintf(stderr, "read_file: Could not open file %s\n", path); | |
+ return 1; | |
+ } else { | |
+ fread(t, sizeof(double), 1, fp); | |
+ fread(t_end, sizeof(double), 1, fp); | |
+ fread(tau, sizeof(double), 1, fp); | |
+ | |
+ fread(itermax, sizeof(int), 1, fp); | |
+ fread(epsilon, sizeof(double), 1, fp); | |
+ fread(omega, sizeof(double), 1, fp); | |
+ fread(gamma, sizeof(double), 1, fp); | |
+ | |
+ fread(gx, sizeof(double), 1, fp); | |
+ fread(gy, sizeof(double), 1, fp); | |
+ fread(re, sizeof(double), 1, fp); | |
+ | |
+ fread(w_left, sizeof(int), 1, fp); | |
+ fread(w_right, sizeof(int), 1, fp); | |
+ fread(w_top, sizeof(int), 1, fp); | |
+ fread(w_bottom, sizeof(int), 1, fp); | |
+ | |
+ fread(dx, sizeof(double), 1, fp); | |
+ fread(dy, sizeof(double), 1, fp); | |
+ fread(nx, sizeof(int), 1, fp); | |
+ fread(ny, sizeof(int), 1, fp); | |
+ | |
+ allocate_memory(P, U, V, *nx, *ny); | |
+ | |
+ for (j=0; j<*ny+1; j++) | |
+ for (i=0; i<*nx+1; i++) { | |
+ fread(&tmp, sizeof(double), 1, fp); | |
+ P[i][j] = tmp; | |
+ } | |
+ /*fread(P[i][j], sizeof(double), 1, fp);*/ | |
+ | |
+ for (j=0; j<*ny+2; j++) | |
+ for (i=0; i<*nx+2; i++) | |
+ fread(&U[i][j], sizeof(double), 1, fp); | |
+ | |
+ for (j=0; j<*ny+2; j++) | |
+ for (i=0; i<*nx+2; i++) | |
+ fread(&V[i][j], sizeof(double), 1, fp); | |
+ | |
+ fclose(fp); | |
+ } | |
+ return 0; | |
+} | |
diff --git a/src/main.c b/src/main.c | |
t@@ -0,0 +1,27 @@ | |
+#include <stdio.h> | |
+#include "file_io.h" | |
+#include "utility.h" | |
+ | |
+int main(int argc, char** argv) | |
+{ | |
+ | |
+ double t, t_end, tau; | |
+ int itermax; | |
+ double epsilon, omega, gamma; | |
+ double gx, gy, re; | |
+ int w_left, w_right, w_top, w_bottom; | |
+ double dx, dy; | |
+ int nx, ny; | |
+ double **P, **U, **V; | |
+ | |
+ read_file("unnamed.dat", &t, &t_end, &tau, &itermax, &epsilon, &omega, &ga… | |
+ &gx, &gy, &re, &w_left, &w_right, &w_top, &w_bottom, | |
+ &dx, &dy, &nx, &ny, P, U, V); | |
+ | |
+ printf("omega = %f\n", omega); | |
+ printf("P[0][0] = %f\n", P[0][0]); | |
+ printf("V[%d][%d] = %f\n", nx, ny, V[nx][ny]); | |
+ | |
+ free_memory(P, U, V, nx); | |
+ return 0; | |
+} | |
diff --git a/src/utility.c b/src/utility.c | |
t@@ -1,36 +1,49 @@ | |
#include <stdio.h> | |
#include <stdlib.h> | |
-/* Linear index from two-dimensional index */ | |
-unsigned int idx(int x, int y, int nx) | |
+int allocate_matrix(double ***M, int nx, int ny) | |
{ | |
- return x + y*nx; | |
+ int i; | |
+ *M = (double**)malloc(sizeof(double*)*nx); | |
+ if (*M == NULL) | |
+ return 1; | |
+ | |
+ for (i=0; i<nx; i++) { | |
+ (*M)[i] = (double*)malloc(sizeof(double)*ny); | |
+ if ((*M)[i] == NULL) | |
+ return 1; | |
+ } | |
+ return 0; | |
+} | |
+ | |
+void free_matrix(double ***M, int nx) | |
+{ | |
+ int i; | |
+ for (i=0; i<nx; i++) | |
+ free((*M)[i]); | |
+ free(*M); | |
} | |
-/* Allocate memory for 2d arrays */ | |
-int allocate_memory(double *p, double *u, double *v, int nx, int ny) | |
+int allocate_memory(double **P, double **U, double **V, int nx, int ny) | |
{ | |
- p = malloc(sizeof(double)*(nx+1)*(ny+1)); | |
- if (p == NULL) { | |
- fprintf(stderr, "allocate_memory: Could not allocate memory for p\n"); | |
+ if (allocate_matrix(&P, nx+1, ny+1) != 0) { | |
+ fprintf(stderr, "allocate_memory: Could not allocate memory for P\n"); | |
return 1; | |
} | |
- u = malloc(sizeof(double)*(nx+2)*(ny+2)); | |
- if (u == NULL) { | |
- fprintf(stderr, "allocate_memory: Could not allocate memory for u\n"); | |
+ if (allocate_matrix(&U, nx+2, ny+2) != 0) { | |
+ fprintf(stderr, "allocate_memory: Could not allocate memory for U\n"); | |
return 1; | |
} | |
- v = malloc(sizeof(double)*(nx+2)*(ny+2)); | |
- if (v == NULL) { | |
- fprintf(stderr, "allocate_memory: Could not allocate memory for v\n"); | |
+ if (allocate_matrix(&V, nx+2, ny+2) != 0) { | |
+ fprintf(stderr, "allocate_memory: Could not allocate memory for V\n"); | |
return 1; | |
} | |
return 0; | |
} | |
-void free_memory(double *p, double *u, double *v) | |
+void free_memory(double **P, double **U, double **V, int nx) | |
{ | |
- free(p); | |
- free(u); | |
- free(v); | |
+ free_matrix(&P, nx); | |
+ free_matrix(&U, nx); | |
+ free_matrix(&V, nx); | |
} | |
diff --git a/src/utility.h b/src/utility.h | |
t@@ -1,8 +1,10 @@ | |
#ifndef UTILITY_H_ | |
#define UTILITY_H_ | |
-unsigned int idx(int x, int y, int nx); | |
-int allocate_memory(double *p, double *u, double *v, int nx, int ny); | |
-void free_memory(double *p, double *u, double *v); | |
+int allocate_matrix(double ***M, int nx, int ny); | |
+void free_matrix(double ***M, int nx); | |
+ | |
+int allocate_memory(double **P, double **U, double **V, int nx, int ny); | |
+void free_memory(double **P, double **U, double **V, int nx); | |
#endif |