Ian Jauslin
summaryrefslogtreecommitdiff
blob: 293be0da6574280eff1b69f525ce07dd03e42433 (plain)
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
#ifndef NAVIERSTOKES_H
#define NAVIERSTOKES_H

#include <complex.h>
#include <stdbool.h>
#include <fftw3.h>

#define M_PI 3.14159265358979323846

// abort signal
extern volatile bool g_abort;

// arrays on which the ffts are performed
typedef struct fft_vects {
  fftw_complex* fft;
  fftw_plan fft_plan;
} fft_vect;

// compute u_k
int uk( int K1, int K2, int N1, int N2, unsigned int nsteps, double nu, double delta, double L, _Complex double* u0, _Complex double* g, bool irreversible, unsigned int print_freq, unsigned int starting_time, unsigned int nthreadsl, FILE* savefile);

// compute energy, enstrophy and alpha
int eea( int K1, int K2, int N1, int N2, unsigned int nsteps, double nu, double delta, double L, _Complex double* u0, _Complex double* g, bool irreversible, unsigned int print_freq, unsigned int running_avg_window, unsigned int starting_time, unsigned int nthreads, FILE* savefile);

// compute solution as a function of time, but do not print anything (useful for debugging)
int quiet( int K1, int K2, int N1, int N2, unsigned int nsteps, double nu, double delta, double L, _Complex double* u0, _Complex double* g, bool irreversible, unsigned int nthreads, FILE* savefile);


// initialize vectors for computation
int  ns_init_tmps( _Complex double **u, _Complex double ** tmp1, _Complex double **tmp2, _Complex double **tmp3, fft_vect* fft1, fft_vect *fft2, fft_vect *ifft, int K1, int K2, int N1, int N2, unsigned int nthreads);
// release vectors
int ns_free_tmps( _Complex double* u, _Complex double* tmp1, _Complex double *tmp2,_Complex double *tmp3, fft_vect fft1, fft_vect fft2, fft_vect ifft);

// copy u0 to u
int copy_u( _Complex double* u, _Complex double* u0, int K1, int K2);

// next time step for Irreversible/reversible Navier-Stokes equation
int ns_step( _Complex double* u, int K1, int K2, int N1, int N2, double nu, double delta, double L, _Complex double* g, fft_vect fft1, fft_vect fft2,fft_vect ifft, _Complex double* tmp1, _Complex double *tmp2, _Complex double *tmp3, bool irreversible);

// right side of Irreversible/reversible Navier-Stokes equation
int ns_rhs( _Complex double* out, _Complex double* u, int K1, int K2, int N1, int N2, double nu, double L, _Complex double* g, fft_vect fft1, fft_vect fft2, fft_vect ifft, bool irreversible);

// convolution term in right side of equation
int ns_T( _Complex double* u, int K1, int K2, int N1, int N2, fft_vect fft1, fft_vect fft2, fft_vect ifft);

// convolution term in right side of equation (computed without fft)
int ns_T_nofft( _Complex double* out, _Complex double* u, int K1, int K2, int N1, int N2);

// compute alpha
double compute_alpha( _Complex double* u, int K1, int K2, int N1, int N2, _Complex double* g, double L, _Complex double* T);
// compute energy
double compute_energy( _Complex double* u, int K1, int K2);
// compute enstrophy
double compute_enstrophy( _Complex double* u, int K1, int K2, double L);

// get index for kx,ky in array of size S
int klookup( int kx, int ky, int S1, int S2);
// get index for kx,ky in array of size K1,K2 in which only the terms with kx>=0 are stored
int klookup_sym( int kx, int ky, int K2);
// get u_{kx,ky} from a vector u in which only the values for kx>=0 are stored, assuming u_{-k}=u_k^*
_Complex double getval_sym( _Complex double* u, int kx, int ky, int K2);

#endif