Ian Jauslin
summaryrefslogtreecommitdiff
blob: 3def963dca0ef447e15228174d1de62e4cc998db (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
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
#include "navier-stokes.h"
#include <math.h>

#define M_PI 3.14159265358979323846

#define CHK 1

// next time step for Irreversible Navier-Stokes equation
int ins_step(_Complex double* u, ns_params params, fft_vects vects, _Complex double* tmp1, _Complex double* tmp2, _Complex double* tmp3){
  int kx,ky;

  // k1
  ins_rhs(tmp1, u, params, vects);
  // add to output
  for(kx=-params.K;kx<=params.K;kx++){
    for(ky=-params.K;ky<=params.K;ky++){
      tmp3[KLOOKUP(kx,ky,params.S)]=u[KLOOKUP(kx,ky,params.S)]+params.h/6*tmp1[KLOOKUP(kx,ky,params.S)];
    }
  }

  // u+h*k1/2
  for(kx=-params.K;kx<=params.K;kx++){
    for(ky=-params.K;ky<=params.K;ky++){
      tmp2[KLOOKUP(kx,ky,params.S)]=u[KLOOKUP(kx,ky,params.S)]+params.h/2*tmp1[KLOOKUP(kx,ky,params.S)];
    }
  }
  // k2
  ins_rhs(tmp1, tmp2, params, vects);
  // add to output
  for(kx=-params.K;kx<=params.K;kx++){
    for(ky=-params.K;ky<=params.K;ky++){
      tmp3[KLOOKUP(kx,ky,params.S)]+=params.h/3*tmp1[KLOOKUP(kx,ky,params.S)];
    }
  }

  // u+h*k2/2
  for(kx=-params.K;kx<=params.K;kx++){
    for(ky=-params.K;ky<=params.K;ky++){
      tmp2[KLOOKUP(kx,ky,params.S)]=u[KLOOKUP(kx,ky,params.S)]+params.h/2*tmp1[KLOOKUP(kx,ky,params.S)];
    }
  }
  // k3
  ins_rhs(tmp1, tmp2, params, vects);
  // add to output
  for(kx=-params.K;kx<=params.K;kx++){
    for(ky=-params.K;ky<=params.K;ky++){
      tmp3[KLOOKUP(kx,ky,params.S)]+=params.h/3*tmp1[KLOOKUP(kx,ky,params.S)];
    }
  }

  // u+h*k3
  for(kx=-params.K;kx<=params.K;kx++){
    for(ky=-params.K;ky<=params.K;ky++){
      tmp2[KLOOKUP(kx,ky,params.S)]=u[KLOOKUP(kx,ky,params.S)]+params.h*tmp1[KLOOKUP(kx,ky,params.S)];
    }
  }
  // k4
  ins_rhs(tmp1, tmp2, params, vects);
  // add to output
  for(kx=-params.K;kx<=params.K;kx++){
    for(ky=-params.K;ky<=params.K;ky++){
      u[KLOOKUP(kx,ky,params.S)]=tmp3[KLOOKUP(kx,ky,params.S)]+params.h/6*tmp1[KLOOKUP(kx,ky,params.S)];
    }
  }

  return(0);
}

// right side of Irreversible Navier-Stokes equation
int ins_rhs(_Complex double* out, _Complex double* u, ns_params params, fft_vects vects){
  int kx,ky;

#if CHK==0
  // F(u/|p|)*F(q1*q2*u/|q|)
  // init to 0
  for(kx=0; kx<params.N*params.N; kx++){
    vects.fft1[kx]=0;
    vects.fft2[kx]=0;
    vects.invfft[kx]=0;
  }
  // fill modes
  for(kx=-params.K;kx<=params.K;kx++){
    for(ky=-params.K;ky<=params.K;ky++){
      if(kx!=0 || ky!=0){
        vects.fft1[KLOOKUP(kx,ky,params.N)]=u[KLOOKUP(kx,ky,params.S)]/sqrt(kx*kx+ky*ky);
        vects.fft2[KLOOKUP(kx,ky,params.N)]=kx*ky*u[KLOOKUP(kx,ky,params.S)]/sqrt(kx*kx+ky*ky);
      }
    }
  }
  // fft
  fftw_execute(vects.fft1_plan);
  fftw_execute(vects.fft2_plan);
  // write to invfft
  for(kx=-2*params.K;kx<=2*params.K;kx++){
    for(ky=-2*params.K;ky<=2*params.K;ky++){
      vects.invfft[KLOOKUP(kx,ky,params.N)]=vects.fft1[KLOOKUP(kx,ky,params.N)]*vects.fft2[KLOOKUP(kx,ky,params.N)];
    }
  }
  // inverse fft
  fftw_execute(vects.invfft_plan);
  
  // write out
  for(kx=0; kx<params.S*params.S; kx++){
    out[kx]=0;
  }
  for(kx=-params.K;kx<=params.K;kx++){
    for(ky=-params.K;ky<=params.K;ky++){
      if(kx!=0 || ky!=0){
        out[KLOOKUP(kx,ky,params.S)]=-4*M_PI*M_PI*params.nu*(kx*kx+ky*ky)*u[KLOOKUP(kx,ky,params.S)]+params.g[KLOOKUP(kx,ky,params.S)]+4*M_PI*M_PI/sqrt(kx*kx+ky*ky)*vects.invfft[KLOOKUP(kx,ky,params.N)]*(kx*kx-ky*ky)/params.N/params.N;
      }
    }
  }

  // F(u/|p|)*F((q1*q1-q2*q2)*u/|q|)
  // init to 0
  for(kx=0; kx<params.N*params.N; kx++){
    vects.fft2[kx]=0;
    vects.invfft[kx]=0;
  }
  // fill modes
  for(kx=-params.K;kx<=params.K;kx++){
    for(ky=-params.K;ky<=params.K;ky++){
      if(kx!=0 || ky!=0){
        vects.fft2[KLOOKUP(kx,ky,params.N)]=(kx*kx-ky*ky)*u[KLOOKUP(kx,ky,params.S)]/sqrt(kx*kx+ky*ky);
      }
    }
  }
  // fft
  fftw_execute(vects.fft2_plan);
  // write to invfft
  for(kx=-2*params.K;kx<=2*params.K;kx++){
    for(ky=-2*params.K;ky<=2*params.K;ky++){
      vects.invfft[KLOOKUP(kx,ky,params.N)]=vects.fft1[KLOOKUP(kx,ky,params.N)]*vects.fft2[KLOOKUP(kx,ky,params.N)];
    }
  }
  // inverse fft
  fftw_execute(vects.invfft_plan);
  
  // write out
  for(kx=-params.K;kx<=params.K;kx++){
    for(ky=-params.K;ky<=params.K;ky++){
      if(kx!=0 || ky!=0){
        out[KLOOKUP(kx,ky,params.S)]=out[KLOOKUP(kx,ky,params.S)]-4*M_PI*M_PI/sqrt(kx*kx+ky*ky)*vects.invfft[KLOOKUP(kx,ky,params.N)]*(kx*ky)/params.N/params.N;
      }
    }
  }
#elif CHK == 1
  // F(-p2/|p|*u)*F(q1*|q|*u)
  // init to 0
  for(kx=0; kx<params.N*params.N; kx++){
    vects.fft1[kx]=0;
    vects.fft2[kx]=0;
    vects.invfft[kx]=0;
  }
  // fill modes
  for(kx=-params.K;kx<=params.K;kx++){
    for(ky=-params.K;ky<=params.K;ky++){
      if(kx!=0 || ky!=0){
        vects.fft1[KLOOKUP(kx,ky,params.N)]=-kx/sqrt(kx*kx+ky*ky)*u[KLOOKUP(kx,ky,params.S)];
        vects.fft2[KLOOKUP(kx,ky,params.N)]=kx*sqrt(kx*kx+ky*ky)*u[KLOOKUP(kx,ky,params.S)];
      }
    }
  }
  // fft
  fftw_execute(vects.fft1_plan);
  fftw_execute(vects.fft2_plan);
  // write to invfft
  for(kx=-2*params.K;kx<=2*params.K;kx++){
    for(ky=-2*params.K;ky<=2*params.K;ky++){
      vects.invfft[KLOOKUP(kx,ky,params.N)]=vects.fft1[KLOOKUP(kx,ky,params.N)]*vects.fft2[KLOOKUP(kx,ky,params.N)] - vects.fft1[KLOOKUP(ky,kx,params.N)]*vects.fft2[KLOOKUP(ky,kx,params.N)];
    }
  }
  
  // inverse fft
  fftw_execute(vects.invfft_plan);
  
  
  // write out
  for(kx=0; kx<params.S*params.S; kx++){
    out[kx]=0;
  }
  for(kx=-params.K;kx<=params.K;kx++){
    for(ky=-params.K;ky<=params.K;ky++){
      if(kx!=0 || ky!=0){
        out[KLOOKUP(kx,ky,params.S)]=-4*M_PI*M_PI*params.nu/params.S*(kx*kx+ky*ky)*u[KLOOKUP(kx,ky,params.S)]+params.g[KLOOKUP(kx,ky,params.S)]+2*M_PI/sqrt(kx*kx+ky*ky)/params.S*vects.invfft[KLOOKUP(kx,ky,params.N)]/params.N/params.N;
      }
    }
  }

#elif CHK==2
  // F(u)*F(q1*u)
  // init to 0
  for(kx=0; kx<params.N*params.N; kx++){
    vects.fft1[kx]=0;
    vects.fft2[kx]=0;
    vects.invfft[kx]=0;
  }
  // fill modes
  for(kx=-params.K;kx<=params.K;kx++){
    for(ky=-params.K;ky<=params.K;ky++){
      // u_k
      vects.fft1[KLOOKUP(kx,ky,params.N)]=u[KLOOKUP(kx,ky,params.S)];
      // 2i\pi k_x u_k
      vects.fft2[KLOOKUP(kx,ky,params.N)]=2*M_PI*I*kx*u[KLOOKUP(kx,ky,params.S)];
    }
  }
  // fft
  fftw_execute(vects.fft1_plan);
  fftw_execute(vects.fft2_plan);
  // write to invfft
  for(kx=-2*params.K;kx<=2*params.K;kx++){
    for(ky=-2*params.K;ky<=2*params.K;ky++){
      vects.invfft[KLOOKUP(kx,ky,params.N)]=vects.fft1[KLOOKUP(kx,ky,params.N)]*vects.fft2[KLOOKUP(kx,ky,params.N)];
    }
  }

  // F(p1/p2*u)*F(q2*u)
  // init to 0
  for(kx=0; kx<params.N*params.N; kx++){
    vects.fft1[kx]=0;
    vects.fft2[kx]=0;
  }
  // fill modes
  for(kx=-params.K;kx<=params.K;kx++){
    for(ky=-params.K;ky<=params.K;ky++){
      // k_x/k_y u_k
      if(ky!=0){
	vects.fft1[KLOOKUP(kx,ky,params.N)]=kx/ky*u[KLOOKUP(kx,ky,params.S)];
      }
      // 2i\pi k_y u_k
      vects.fft2[KLOOKUP(kx,ky,params.N)]=2*M_PI*I*ky*u[KLOOKUP(kx,ky,params.S)];
    }
  }
  // fft
  fftw_execute(vects.fft1_plan);
  fftw_execute(vects.fft2_plan);
  // write to invfft
  for(kx=-2*params.K;kx<=2*params.K;kx++){
    for(ky=-2*params.K;ky<=2*params.K;ky++){
      vects.invfft[KLOOKUP(kx,ky,params.N)]+=-vects.fft1[KLOOKUP(kx,ky,params.N)]*vects.fft2[KLOOKUP(kx,ky,params.N)];
    }
  }

  // inverse fft
  fftw_execute(vects.invfft_plan);

  /*
  // check: convolution instead of fft
  for(kx=0; kx<params.S*params.S; kx++){
    out[kx]=0;
  }
  for(kx=-params.K;kx<=params.K;kx++){
    for(ky=-params.K;ky<=params.K;ky++){
      for(px=-params.K;px<=params.K;px++){
	for(py=-params.K;py<=params.K;py++){
	  if(kx-px<=params.K && kx-px>=-params.K && ky-py<=params.K && ky-py>=-params.K){
	    out[KLOOKUP(kx,ky,params.S)]+=2*M_PI*I*(u[KLOOKUP(px,py,params.S)]*(kx-px)*u[KLOOKUP(kx-px,ky-py,params.S)]-(py==0?0:px/py*u[KLOOKUP(px,py,params.S)]*(ky-py)*u[KLOOKUP(kx-px,ky-py,params.S)]));
	  }
	}
      }
      dd=(__real__ vects.invfft[KLOOKUP(kx,ky,params.N)]/params.N/params.N-__real__ out[KLOOKUP(kx,ky,params.S)])*(__real__ vects.invfft[KLOOKUP(kx,ky,params.N)]/params.N/params.N-__real__ out[KLOOKUP(kx,ky,params.S)])+(__imag__ vects.invfft[KLOOKUP(kx,ky,params.N)]/params.N/params.N-__imag__ out[KLOOKUP(kx,ky,params.S)])*(__imag__ vects.invfft[KLOOKUP(kx,ky,params.N)]/params.N/params.N-__imag__ out[KLOOKUP(kx,ky,params.S)]);
      if(dd>1e-25){
	printf("%d %d % .8e\n",kx,ky, dd);
      }
    }
  }
  */


  // write out
  for(kx=0; kx<params.S*params.S; kx++){
    out[kx]=0;
  }
  for(kx=-params.K;kx<=params.K;kx++){
    for(ky=-params.K;ky<=params.K;ky++){
      out[KLOOKUP(kx,ky,params.S)]=-4*M_PI*M_PI*params.nu*(kx*kx+ky*ky)*u[KLOOKUP(kx,ky,params.S)]+params.g[KLOOKUP(kx,ky,params.S)]+vects.invfft[KLOOKUP(kx,ky,params.N)]/params.N/params.N;
    }
  }

#endif
  return(0);
}


// compute alpha
_Complex double compute_alpha(_Complex double* u, ns_params params){
  _Complex double num=0;
  _Complex double denom=0;
  int kx,ky;

  for(kx=-params.K;kx<=params.K;kx++){
    for(ky=-params.K;ky<=params.K;ky++){
      denom+=(kx*kx+ky*ky)*(kx*kx+ky*ky)*u[KLOOKUP(kx,ky,params.S)]*conj(u[KLOOKUP(kx,ky,params.S)])*(1+(ky!=0?kx*kx/ky/ky:0));
      num+=(kx*kx+ky*ky)*u[KLOOKUP(kx,ky,params.S)]*conj(params.g[KLOOKUP(kx,ky,params.S)])*(1+(ky!=0?kx*kx/ky/ky:0));
    }
  }

  return(num/denom);
}