Nstrophy/src/navier-stokes.c

650 lines
14 KiB
C
Raw Normal View History

2018-01-11 22:48:14 +00:00
#include "navier-stokes.h"
2022-05-27 16:09:17 -04:00
#include "io.h"
2018-01-11 22:48:14 +00:00
#include <math.h>
2022-05-18 09:57:06 +02:00
#include <stdlib.h>
2023-04-12 19:05:01 -04:00
#include <string.h>
2018-01-11 22:48:14 +00:00
2022-05-18 09:57:06 +02:00
// compute solution as a function of time
int uk(
int K1,
int K2,
int N1,
int N2,
2023-04-16 00:21:42 -04:00
uint64_t nsteps,
2022-05-18 09:57:06 +02:00
double nu,
double delta,
2022-05-25 11:12:02 -04:00
double L,
2022-05-26 15:05:30 -04:00
_Complex double* u0,
_Complex double* g,
2023-04-05 20:33:38 -04:00
bool irreversible,
2023-04-16 00:21:42 -04:00
uint64_t print_freq,
uint64_t starting_time,
2022-05-27 16:09:17 -04:00
unsigned int nthreads,
FILE* savefile
2022-05-18 09:57:06 +02:00
){
_Complex double* u;
_Complex double* tmp1;
_Complex double* tmp2;
_Complex double* tmp3;
2023-04-16 00:21:42 -04:00
uint64_t t;
2022-05-18 09:57:06 +02:00
fft_vect fft1;
fft_vect fft2;
fft_vect ifft;
int kx,ky;
2022-05-18 23:52:01 +02:00
ns_init_tmps(&u, &tmp1, &tmp2, &tmp3, &fft1, &fft2, &ifft, K1, K2, N1, N2, nthreads);
2022-05-26 15:05:30 -04:00
// copy initial condition
copy_u(u, u0, K1, K2);
2022-05-18 09:57:06 +02:00
2022-05-19 17:51:45 +02:00
// print column headers
printf("# 1:i 2:t ");
t=3;
for(kx=-K1;kx<=K1;kx++){
for (ky=-K2;ky<=K2;ky++){
2023-04-16 00:21:42 -04:00
printf(" %6lu:(%4d,%4d)r ",t,kx,ky);
2022-05-19 17:51:45 +02:00
t++;
2023-04-16 00:21:42 -04:00
printf(" %6lu:(%4d,%4d)i ",t,kx,ky);
2022-05-19 17:51:45 +02:00
t++;
}
}
2022-05-18 09:57:06 +02:00
// iterate
2023-04-17 14:45:21 -04:00
for(t=starting_time;nsteps==0 || t<starting_time+nsteps;t++){
2023-04-05 20:33:38 -04:00
ns_step(u, K1, K2, N1, N2, nu, delta, L, g, fft1, fft2, ifft, tmp1, tmp2, tmp3, irreversible);
2022-05-18 09:57:06 +02:00
if(t%print_freq==0){
2023-04-16 00:21:42 -04:00
fprintf(stderr,"%lu % .8e ",t,t*delta);
printf("%8lu % .15e ",t,t*delta);
2022-05-18 09:57:06 +02:00
for(kx=-K1;kx<=K1;kx++){
for (ky=-K2;ky<=K2;ky++){
if (kx*kx+ky*ky<=1){
2023-04-11 18:45:45 -04:00
fprintf(stderr,"% .8e % .8e ",__real__ getval_sym(u,kx,ky,K2), __imag__ getval_sym(u, kx,ky,K2));
2022-05-18 09:57:06 +02:00
}
2023-04-11 18:45:45 -04:00
printf("% .15e % .15e ",__real__ getval_sym(u, kx,ky,K2), __imag__ getval_sym(u, kx,ky,K2));
2022-05-18 09:57:06 +02:00
}
}
fprintf(stderr,"\n");
printf("\n");
}
}
2022-05-27 16:09:17 -04:00
// save final entry to savefile
2023-04-14 15:01:52 -04:00
write_u_bin(u, K1, K2, savefile);
2022-05-27 16:09:17 -04:00
2022-05-18 09:57:06 +02:00
ns_free_tmps(u, tmp1, tmp2, tmp3, fft1, fft2, ifft);
return(0);
}
// compute energy, enstrophy, alpha as a function of time in the I-NS equation
int eea(
2022-05-18 09:57:06 +02:00
int K1,
int K2,
int N1,
int N2,
2023-04-16 00:21:42 -04:00
uint64_t nsteps,
2022-05-18 09:57:06 +02:00
double nu,
double delta,
2022-05-25 11:12:02 -04:00
double L,
2022-05-26 15:05:30 -04:00
_Complex double* u0,
_Complex double* g,
2023-04-05 20:33:38 -04:00
bool irreversible,
2023-04-16 00:21:42 -04:00
uint64_t print_freq,
uint64_t starting_time,
2022-05-27 16:09:17 -04:00
unsigned int nthreads,
2023-04-12 19:05:01 -04:00
FILE* savefile,
// for interrupt recovery
char* cmd_string,
char* params_string,
char* savefile_string
2022-05-18 09:57:06 +02:00
){
_Complex double* u;
_Complex double* tmp1;
_Complex double* tmp2;
_Complex double* tmp3;
double alpha, energy, enstrophy;
double avg_e,avg_a,avg_en;
2023-04-12 14:37:02 -04:00
// index
2023-04-16 00:21:42 -04:00
uint64_t t;
2022-05-18 09:57:06 +02:00
fft_vect fft1;
fft_vect fft2;
fft_vect ifft;
2022-05-18 23:52:01 +02:00
ns_init_tmps(&u, &tmp1, &tmp2, &tmp3, &fft1, &fft2, &ifft, K1, K2, N1, N2, nthreads);
2022-05-26 15:05:30 -04:00
// copy initial condition
copy_u(u, u0, K1, K2);
2022-05-18 09:57:06 +02:00
// init running average
avg_e=0;
avg_a=0;
avg_en=0;
2022-05-18 09:57:06 +02:00
// special first case when starting_time is not a multiple of print_freq
2023-04-16 00:21:42 -04:00
uint64_t first_box = print_freq - (starting_time % print_freq);
2023-04-12 14:37:02 -04:00
2022-05-18 09:57:06 +02:00
// iterate
2023-04-17 14:45:21 -04:00
for(t=starting_time;nsteps==0 || t<starting_time+nsteps;t++){
2023-04-05 20:33:38 -04:00
ns_step(u, K1, K2, N1, N2, nu, delta, L, g, fft1, fft2, ifft, tmp1, tmp2, tmp3, irreversible);
energy=compute_energy(u, K1, K2);
alpha=compute_alpha(u, K1, K2, g, L);
enstrophy=compute_enstrophy(u, K1, K2, L);
2022-05-18 09:57:06 +02:00
// running average
// reset averages
if(t % print_freq == 1){
avg_e=0;
avg_a=0;
avg_en=0;
2022-05-18 09:57:06 +02:00
}
2023-04-12 14:37:02 -04:00
// compute average
// different computationin first box if starting_time is not a multiple of print_freq
if(t < starting_time + first_box){
avg_e+=energy/first_box;
avg_a+=alpha/first_box;
avg_en+=enstrophy/first_box;
} else {
avg_e+=energy/print_freq;
avg_a+=alpha/print_freq;
avg_en+=enstrophy/print_freq;
2023-04-12 14:37:02 -04:00
}
if(t>starting_time && t%print_freq==0){
2023-04-16 00:21:42 -04:00
fprintf(stderr,"%lu % .8e % .8e % .8e % .8e % .8e % .8e % .8e\n",t,t*delta, avg_a, avg_e, avg_en, alpha, energy, enstrophy);
printf("%8lu % .15e % .15e % .15e % .15e % .15e % .15e % .15e\n",t,t*delta, avg_a, avg_e, avg_en, alpha, energy, enstrophy);
2022-05-18 09:57:06 +02:00
}
2023-04-12 15:23:35 -04:00
// catch abort signal
if (g_abort){
// print u to stderr if no savefile
if (savefile==NULL){
savefile=stderr;
}
2023-04-14 15:01:52 -04:00
2023-04-12 19:05:01 -04:00
fprintf(savefile,"# Interrupted computation. Resume with\n");
// command to resume
fprintf(savefile,"#! ");
fprintf(savefile, cmd_string);
// params
// allocate buffer for params
2023-04-14 15:11:55 -04:00
char* params=calloc(sizeof(char), strlen(params_string)+1);
2023-04-12 19:05:01 -04:00
strcpy(params, params_string);
remove_entry(params, "starting_time");
remove_entry(params, "init");
remove_entry(params, "nsteps");
2023-04-18 00:48:54 -04:00
fprintf(savefile," -p \"%s;starting_time=%lu;nsteps=%lu;init=file:%s\"", params, t+1, (nsteps == 0 ? 0 : nsteps-t-1), savefile_string);
2023-04-12 19:05:01 -04:00
free(params);
fprintf(savefile," energy\n");
2023-04-12 15:23:35 -04:00
break;
}
2022-05-18 09:57:06 +02:00
}
2023-04-12 15:23:35 -04:00
// save final entry to savefile
2023-04-14 15:01:52 -04:00
if(savefile==stderr || savefile==stdout){
write_u(u, K1, K2, savefile);
} else {
write_u_bin(u, K1, K2, savefile);
}
2023-04-12 15:23:35 -04:00
2022-05-18 09:57:06 +02:00
ns_free_tmps(u, tmp1, tmp2, tmp3, fft1, fft2, ifft);
return(0);
}
2022-05-18 22:22:42 +02:00
// 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,
2023-04-16 00:21:42 -04:00
uint64_t nsteps,
2022-05-18 22:22:42 +02:00
double nu,
double delta,
2022-05-25 11:12:02 -04:00
double L,
uint64_t starting_time,
2022-05-26 15:05:30 -04:00
_Complex double* u0,
_Complex double* g,
2023-04-05 20:33:38 -04:00
bool irreversible,
2022-05-27 16:09:17 -04:00
unsigned int nthreads,
FILE* savefile
2022-05-18 22:22:42 +02:00
){
_Complex double* u;
_Complex double* tmp1;
_Complex double* tmp2;
_Complex double* tmp3;
2023-04-16 00:21:42 -04:00
uint64_t t;
2022-05-18 22:22:42 +02:00
fft_vect fft1;
fft_vect fft2;
fft_vect ifft;
2022-05-18 23:52:01 +02:00
ns_init_tmps(&u, &tmp1, &tmp2, &tmp3, &fft1, &fft2, &ifft, K1, K2, N1, N2, nthreads);
2022-05-26 15:05:30 -04:00
// copy initial condition
copy_u(u, u0, K1, K2);
2022-05-18 22:22:42 +02:00
// iterate
for(t=starting_time;nsteps==0 || t<starting_time+nsteps;t++){
2023-04-05 20:33:38 -04:00
ns_step(u, K1, K2, N1, N2, nu, delta, L, g, fft1, fft2, ifft, tmp1, tmp2, tmp3, irreversible);
2022-05-18 22:22:42 +02:00
}
2022-05-27 16:09:17 -04:00
// save final entry to savefile
write_u(u, K1, K2, savefile);
2022-05-18 22:22:42 +02:00
ns_free_tmps(u, tmp1, tmp2, tmp3, fft1, fft2, ifft);
return(0);
}
2022-05-18 09:57:06 +02:00
// 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,
2022-05-18 23:52:01 +02:00
int N2,
unsigned int nthreads
2022-05-18 09:57:06 +02:00
){
// velocity field
2023-04-11 18:45:45 -04:00
*u=calloc(sizeof(_Complex double),(K1+1)*(2*K2+1));
2022-05-18 09:57:06 +02:00
// allocate tmp vectors for computation
2023-04-11 18:45:45 -04:00
*tmp1=calloc(sizeof(_Complex double),(K1+1)*(2*K2+1));
*tmp2=calloc(sizeof(_Complex double),(K1+1)*(2*K2+1));
*tmp3=calloc(sizeof(_Complex double),(K1+1)*(2*K2+1));
2022-05-18 09:57:06 +02:00
2022-05-18 23:52:01 +02:00
// init threads
fftw_init_threads();
fftw_plan_with_nthreads(nthreads);
2022-05-18 09:57:06 +02:00
// prepare vectors for fft
fft1->fft=fftw_malloc(sizeof(fftw_complex)*N1*N2);
fft1->fft_plan=fftw_plan_dft_2d(N1,N2, fft1->fft, fft1->fft, FFTW_FORWARD, FFTW_MEASURE);
fft2->fft=fftw_malloc(sizeof(fftw_complex)*N1*N2);
fft2->fft_plan=fftw_plan_dft_2d(N1,N2, fft2->fft, fft2->fft, FFTW_FORWARD, FFTW_MEASURE);
ifft->fft=fftw_malloc(sizeof(fftw_complex)*N1*N2);
ifft->fft_plan=fftw_plan_dft_2d(N1,N2, ifft->fft, ifft->fft, FFTW_BACKWARD, FFTW_MEASURE);
return 0;
}
// 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
){
// free memory
fftw_destroy_plan(fft1.fft_plan);
fftw_destroy_plan(fft2.fft_plan);
fftw_destroy_plan(ifft.fft_plan);
fftw_free(fft1.fft);
fftw_free(fft2.fft);
fftw_free(ifft.fft);
2022-05-18 23:52:01 +02:00
fftw_cleanup_threads();
2022-05-18 09:57:06 +02:00
free(tmp3);
free(tmp2);
free(tmp1);
free(u);
return 0;
}
2022-05-26 15:05:30 -04:00
// copy u0 to u
int copy_u(
2022-05-18 09:57:06 +02:00
_Complex double* u,
2022-05-26 15:05:30 -04:00
_Complex double* u0,
2022-05-18 09:57:06 +02:00
int K1,
int K2
){
2022-05-26 15:05:30 -04:00
int i;
2022-05-18 09:57:06 +02:00
2023-04-11 18:45:45 -04:00
for(i=0;i<(K1+1)*(2*K2+1);i++){
2022-05-26 15:05:30 -04:00
u[i]=u0[i];
2022-05-18 09:57:06 +02:00
}
return 0;
}
2023-04-05 20:33:38 -04:00
// next time step
int ns_step(
2022-05-18 09:57:06 +02:00
_Complex double* u,
int K1,
int K2,
int N1,
int N2,
double nu,
double delta,
2022-05-25 11:12:02 -04:00
double L,
_Complex double* g,
2022-05-18 09:57:06 +02:00
fft_vect fft1,
fft_vect fft2,
fft_vect ifft,
_Complex double* tmp1,
_Complex double* tmp2,
2023-04-05 20:33:38 -04:00
_Complex double* tmp3,
bool irreversible
2022-05-18 09:57:06 +02:00
){
2018-01-11 22:48:14 +00:00
int kx,ky;
// k1
2023-04-05 20:33:38 -04:00
ns_rhs(tmp1, u, K1, K2, N1, N2, nu, L, g, fft1, fft2, ifft, irreversible);
2018-01-11 22:48:14 +00:00
// add to output
2023-04-11 18:45:45 -04:00
for(kx=0;kx<=K1;kx++){
2022-05-18 09:57:06 +02:00
for(ky=-K2;ky<=K2;ky++){
2023-04-11 18:45:45 -04:00
tmp3[klookup_sym(kx,ky,K2)]=u[klookup_sym(kx,ky,K2)]+delta/6*tmp1[klookup_sym(kx,ky,K2)];
2018-01-11 22:48:14 +00:00
}
}
// u+h*k1/2
2023-04-11 18:45:45 -04:00
for(kx=0;kx<=K1;kx++){
2022-05-18 09:57:06 +02:00
for(ky=-K2;ky<=K2;ky++){
2023-04-11 18:45:45 -04:00
tmp2[klookup_sym(kx,ky,K2)]=u[klookup_sym(kx,ky,K2)]+delta/2*tmp1[klookup_sym(kx,ky,K2)];
2018-01-11 22:48:14 +00:00
}
}
// k2
2023-04-05 20:33:38 -04:00
ns_rhs(tmp1, tmp2, K1, K2, N1, N2, nu, L, g, fft1, fft2, ifft, irreversible);
2018-01-11 22:48:14 +00:00
// add to output
2023-04-11 18:45:45 -04:00
for(kx=0;kx<=K1;kx++){
2022-05-18 09:57:06 +02:00
for(ky=-K2;ky<=K2;ky++){
2023-04-11 18:45:45 -04:00
tmp3[klookup_sym(kx,ky,K2)]+=delta/3*tmp1[klookup_sym(kx,ky,K2)];
2018-01-11 22:48:14 +00:00
}
}
// u+h*k2/2
2023-04-11 18:45:45 -04:00
for(kx=0;kx<=K1;kx++){
2022-05-18 09:57:06 +02:00
for(ky=-K2;ky<=K2;ky++){
2023-04-11 18:45:45 -04:00
tmp2[klookup_sym(kx,ky,K2)]=u[klookup_sym(kx,ky,K2)]+delta/2*tmp1[klookup_sym(kx,ky,K2)];
2018-01-11 22:48:14 +00:00
}
}
// k3
2023-04-05 20:33:38 -04:00
ns_rhs(tmp1, tmp2, K1, K2, N1, N2, nu, L, g, fft1, fft2, ifft, irreversible);
2018-01-11 22:48:14 +00:00
// add to output
2023-04-11 18:45:45 -04:00
for(kx=0;kx<=K1;kx++){
2022-05-18 09:57:06 +02:00
for(ky=-K2;ky<=K2;ky++){
2023-04-11 18:45:45 -04:00
tmp3[klookup_sym(kx,ky,K2)]+=delta/3*tmp1[klookup_sym(kx,ky,K2)];
2018-01-11 22:48:14 +00:00
}
}
// u+h*k3
2023-04-11 18:45:45 -04:00
for(kx=0;kx<=K1;kx++){
2022-05-18 09:57:06 +02:00
for(ky=-K2;ky<=K2;ky++){
2023-04-11 18:45:45 -04:00
tmp2[klookup_sym(kx,ky,K2)]=u[klookup_sym(kx,ky,K2)]+delta*tmp1[klookup_sym(kx,ky,K2)];
2018-01-11 22:48:14 +00:00
}
}
// k4
2023-04-05 20:33:38 -04:00
ns_rhs(tmp1, tmp2, K1, K2, N1, N2, nu, L, g, fft1, fft2, ifft, irreversible);
2018-01-11 22:48:14 +00:00
// add to output
2023-04-11 18:45:45 -04:00
for(kx=0;kx<=K1;kx++){
2022-05-18 09:57:06 +02:00
for(ky=-K2;ky<=K2;ky++){
2023-04-11 18:45:45 -04:00
u[klookup_sym(kx,ky,K2)]=tmp3[klookup_sym(kx,ky,K2)]+delta/6*tmp1[klookup_sym(kx,ky,K2)];
2023-04-06 11:28:18 -04:00
}
}
2018-01-11 22:48:14 +00:00
return(0);
}
2023-04-05 20:33:38 -04:00
// right side of Irreversible/Reversible Navier-Stokes equation
int ns_rhs(
2022-05-18 09:57:06 +02:00
_Complex double* out,
_Complex double* u,
int K1,
int K2,
int N1,
int N2,
double nu,
2022-05-25 11:12:02 -04:00
double L,
_Complex double* g,
2022-05-18 09:57:06 +02:00
fft_vect fft1,
fft_vect fft2,
2023-04-05 20:33:38 -04:00
fft_vect ifft,
bool irreversible
2022-05-18 09:57:06 +02:00
){
2018-01-11 22:48:14 +00:00
int kx,ky;
2022-05-18 09:57:06 +02:00
int i;
2023-04-05 20:33:38 -04:00
double alpha;
2018-01-11 22:48:14 +00:00
// compute convolution term
ns_T(u,K1,K2,N1,N2,fft1,fft2,ifft);
2023-04-05 20:33:38 -04:00
if (irreversible) {
alpha=nu;
} else {
alpha=compute_alpha(u,K1,K2,g,L);
2023-04-05 20:33:38 -04:00
}
2022-05-18 22:22:42 +02:00
2023-04-11 18:45:45 -04:00
for(i=0; i<(K1+1)*(2*K2+1); i++){
out[i]=0;
}
2023-04-11 18:45:45 -04:00
for(kx=0;kx<=K1;kx++){
for(ky=-K2;ky<=K2;ky++){
if(kx!=0 || ky!=0){
2023-04-11 18:45:45 -04:00
out[klookup_sym(kx,ky,K2)]=-4*M_PI*M_PI/L/L*alpha*(kx*kx+ky*ky)*u[klookup_sym(kx,ky,K2)]+g[klookup_sym(kx,ky,K2)]+4*M_PI*M_PI/L/L/sqrt(kx*kx+ky*ky)*ifft.fft[klookup(kx,ky,N1,N2)];
}
}
}
return(0);
}
// convolution term in right side of convolution equation
int ns_T(
_Complex double* u,
int K1,
int K2,
int N1,
int N2,
fft_vect fft1,
fft_vect fft2,
fft_vect ifft
){
int kx,ky;
int i;
2018-01-12 19:20:59 +00:00
// F(px/|p|*u)*F(qy*|q|*u)
2018-01-11 22:48:14 +00:00
// init to 0
2022-05-18 09:57:06 +02:00
for(i=0; i<N1*N2; i++){
fft1.fft[i]=0;
fft2.fft[i]=0;
ifft.fft[i]=0;
2018-01-11 22:48:14 +00:00
}
// fill modes
2022-05-18 09:57:06 +02:00
for(kx=-K1;kx<=K1;kx++){
for(ky=-K2;ky<=K2;ky++){
2018-01-11 22:48:14 +00:00
if(kx!=0 || ky!=0){
2023-04-11 18:45:45 -04:00
fft1.fft[klookup(kx,ky,N1,N2)]=kx/sqrt(kx*kx+ky*ky)*getval_sym(u, kx,ky,K2)/N1;
fft2.fft[klookup(kx,ky,N1,N2)]=ky*sqrt(kx*kx+ky*ky)*getval_sym(u, kx,ky,K2)/N2;
2018-01-11 22:48:14 +00:00
}
}
}
2018-01-12 19:20:59 +00:00
2018-01-11 22:48:14 +00:00
// fft
2022-05-18 09:57:06 +02:00
fftw_execute(fft1.fft_plan);
fftw_execute(fft2.fft_plan);
// write to ifft
for(i=0;i<N1*N2;i++){
ifft.fft[i]=fft1.fft[i]*fft2.fft[i];
2018-01-11 22:48:14 +00:00
}
2018-01-12 19:20:59 +00:00
// F(py/|p|*u)*F(qx*|q|*u)
2018-01-11 22:48:14 +00:00
// init to 0
2022-05-18 09:57:06 +02:00
for(i=0; i<N1*N2; i++){
fft1.fft[i]=0;
fft2.fft[i]=0;
2018-01-11 22:48:14 +00:00
}
// fill modes
2022-05-18 09:57:06 +02:00
for(kx=-K1;kx<=K1;kx++){
for(ky=-K2;ky<=K2;ky++){
2018-01-11 22:48:14 +00:00
if(kx!=0 || ky!=0){
2023-04-11 18:45:45 -04:00
fft1.fft[klookup(kx,ky,N1,N2)]=ky/sqrt(kx*kx+ky*ky)*getval_sym(u, kx,ky,K2)/N1;
fft2.fft[klookup(kx,ky,N1,N2)]=kx*sqrt(kx*kx+ky*ky)*getval_sym(u, kx,ky,K2)/N2;
2018-01-11 22:48:14 +00:00
}
}
}
2018-01-12 19:20:59 +00:00
2018-01-11 22:48:14 +00:00
// fft
2022-05-18 09:57:06 +02:00
fftw_execute(fft1.fft_plan);
fftw_execute(fft2.fft_plan);
// write to ifft
for(i=0;i<N1*N2;i++){
ifft.fft[i]=ifft.fft[i]-fft1.fft[i]*fft2.fft[i];
2018-01-11 22:48:14 +00:00
}
2018-01-12 19:20:59 +00:00
2018-01-11 22:48:14 +00:00
// inverse fft
2022-05-18 09:57:06 +02:00
fftw_execute(ifft.fft_plan);
2022-05-17 14:31:22 +02:00
2018-01-11 22:48:14 +00:00
return(0);
}
2022-05-18 22:22:42 +02:00
// convolution term in right side of convolution equation, computed without fourier transform
int ns_T_nofft(
_Complex double* out,
_Complex double* u,
int K1,
int K2,
int N1,
int N2
){
int kx,ky;
int px,py;
int qx,qy;
// loop over K's (needs N1>=2*K1+1 and N2>=2*K2+1)
if (N1<2*K1+1 || N2<2*K2+1){
fprintf(stderr,"error: N1 and N2 need t be >= 2*K1+1 and 2*K2+1 respectively\n");
2022-05-18 22:22:42 +02:00
return(-1);
}
for(kx=-K1;kx<=K1;kx++){
for(ky=-K2;ky<=K2;ky++){
2022-05-18 22:22:42 +02:00
// init
out[klookup(kx,ky,N1,N2)]=0.;
for(px=-K1;px<=K1;px++){
for(py=-K2;py<=K2;py++){
qx=kx-px;
qy=ky-py;
// cutoff in q
if(qx>=-K1 && qx<=K1 && qy>=-K2 && qy<=K2 && qx*qx+qy*qy>0 && px*px+py*py>0){
2023-04-11 18:45:45 -04:00
out[klookup(kx,ky,N1,N2)]+=(-qx*py+qy*px)*sqrt(qx*qx+qy*qy)/sqrt(px*px+py*py)*getval_sym(u, px,py,K2)*getval_sym(u, qx,qy,K2);
2022-05-18 22:22:42 +02:00
}
}
}
}
}
return 0;
}
2018-01-11 22:48:14 +00:00
2023-04-05 20:33:38 -04:00
// compute alpha
double compute_alpha(
_Complex double* u,
int K1,
int K2,
_Complex double* g,
double L
){
_Complex double num=0;
double denom=0;
int kx,ky;
num=0.;
denom=0.;
for(kx=-K1;kx<=K1;kx++){
for(ky=-K2;ky<=K2;ky++){
num+=L*L/4/M_PI/M_PI*(kx*kx+ky*ky)*getval_sym(g, kx,ky,K2)*conj(getval_sym(u, kx,ky,K2));
denom+=__real__ (kx*kx+ky*ky)*(kx*kx+ky*ky)*getval_sym(u, kx,ky,K2)*conj(getval_sym(u, kx,ky,K2));
}
}
return __real__ num/denom;
}
// compute energy
double compute_energy(
_Complex double* u,
int K1,
int K2
){
int kx,ky;
double out=0.;
for(kx=-K1;kx<=K1;kx++){
for (ky=-K2;ky<=K2;ky++){
2023-04-11 18:45:45 -04:00
out+=__real__ (getval_sym(u, kx,ky,K2)*conj(getval_sym(u, kx,ky,K2)));
}
}
return out;
}
// compute enstrophy
double compute_enstrophy(
_Complex double* u,
int K1,
int K2,
double L
){
int kx,ky;
double out=0.;
for(kx=-K1;kx<=K1;kx++){
for (ky=-K2;ky<=K2;ky++){
2023-04-11 18:45:45 -04:00
out+=__real__ (4*M_PI*M_PI/L/L*(kx*kx+ky*ky)*getval_sym(u, kx,ky,K2)*conj(getval_sym(u, kx,ky,K2)));
}
}
return out;
2018-01-11 22:48:14 +00:00
}
2022-05-18 09:57:06 +02:00
// get index for kx,ky in array of size S
int klookup(
int kx,
int ky,
int S1,
int S2
){
return (kx>=0 ? kx : S1+kx)*S2 + (ky>=0 ? ky : S2+ky);
}
2023-04-11 18:45:45 -04:00
// 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
){
if (kx<0) {
fprintf(stderr, "bug!: attempting to access a symmetrized vector at kx<0\n Contact Ian at ian.jauslin@rutgers.edu!\n");
exit(-1);
}
return kx*(2*K2+1) + (ky>=0 ? ky : (2*K2+1)+ky);
}
// 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
){
return (kx>=0 ? u[klookup_sym(kx,ky,K2)] : conj(u[klookup_sym(-kx,-ky,K2)]));
}