-
Notifications
You must be signed in to change notification settings - Fork 0
/
BEC_Harmonic_GS_imag_time.cpp
executable file
·142 lines (117 loc) · 5.22 KB
/
BEC_Harmonic_GS_imag_time.cpp
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
#include"fftw3.h"
#include<vector>
#include<iostream>
#include"headers/bec.h"
#include<algorithm>
#include<functional>
#include<iterator>
#include<cmath>
#include<complex>
using namespace std;
inline vector<complex<double> > S3 (const double deltat, const double omega, const double g, vector<complex<double>
> data, const vector<complex<double> >& potential, fftw_complex* in, fftw_complex* out)
{
complex<double>I(0,1);
unsigned int N = data.size();
vector<complex<double> > kvec(N);
vector<complex<double> > kfac(N);
for (unsigned int i = 0; i<N; ++i) {
kvec[i]= (double(i)<double(N)/2.) ? double(i): -double(N-i);
kfac[i]= exp(-1.*deltat*omega/2.*kvec[i]*kvec[i]/2.);
}
in=reinterpret_cast<fftw_complex*>(&data.front());
fftw_plan fourier1 = fftw_plan_dft_1d(N, in, out, FFTW_FORWARD, FFTW_ESTIMATE);
fftw_execute(fourier1);
fftw_destroy_plan(fourier1);
for(unsigned int i = 0; i<N; ++i){
double real_part = out[i][0]; double imag_part = out[i][1];
complex<double> fac = kfac[i];
out[i][0] = (real(fac)*real_part-imag(fac)*imag_part)/double(N);
out[i][1] = (real(fac)*imag_part+imag(fac)*real_part)/double(N);
}
fftw_plan backfourier1 = fftw_plan_dft_1d(N, out, in, FFTW_BACKWARD, FFTW_ESTIMATE);
fftw_execute(backfourier1);
fftw_destroy_plan(backfourier1);
for(unsigned int i = 0; i<N; ++i)
data[i]*=exp(-1*g*deltat*norm(data[i]));
transform(data.begin(), data.end(), potential.begin(), data.begin(), multiplies<complex<double> >() );
in=reinterpret_cast<fftw_complex*>(&data.front());
fftw_plan fourier3 = fftw_plan_dft_1d(N, in, out, FFTW_FORWARD, FFTW_ESTIMATE);
fftw_execute(fourier3);
fftw_destroy_plan(fourier3);
for(unsigned int i = 0; i<N; ++i){
double real_part = out[i][0]; double imag_part = out[i][1];
complex<double> fac = kfac[i];
out[i][0] = (real(fac)*real_part-imag(fac)*imag_part)/double(N);
out[i][1] = (real(fac)*imag_part+imag(fac)*real_part)/double(N);
}
fftw_plan backfourier3 = fftw_plan_dft_1d(N, out, in, FFTW_BACKWARD, FFTW_ESTIMATE);
fftw_execute(backfourier3);
fftw_destroy_plan(backfourier3);
return data;
}
int main(){
const double pi = 4.*atan(1.);
const double g = 0.5;
const double L = 10.;
const double q = 40.*pi/L;
double newnorm = 0.;
double d = 0.1;
int N=1024;
const int period = 10000;
const double deltat = 2.*pi/double(period);
const double omega[4] = {1.3151861047504085, -1.1776799841887, 0.2355733213359357, 0.784513610477560};
typedef complex<double> complex_type;
typedef fftw_complex complex_c_type;
typedef vector<complex<double> >::size_type size_t;
const complex_type I(0,1);
vector<complex_type> data(N, 100.);
/*
for(size_t i = 0; i<N; ++i){
double xpos =20.*0.5-double(i)*20./double(N);
data[i]=1./sqrt(2*pi)*exp(-xpos*xpos/2.);
}
*/
vector<complex<double> > potential(N);
for(unsigned int i = 0; i<N; ++i){
double xpos =L*0.5-double(i)*L/double(N);
double sinqx = sin(q*xpos);
potential[i]=exp(-1.*deltat*((xpos*xpos/2.)+q*q*0.5*sinqx*sinqx));
}
fftw_complex* in = new fftw_complex[sizeof(fftw_complex)*N];
fftw_complex* out = new fftw_complex[sizeof(fftw_complex)*N];
double t=0;
int stepcounter=0;
for(;t<50; t+=deltat, ++stepcounter){
data=S3(deltat, 1., g, data, potential, in, out);
newnorm = L_2_norm(data.begin(), data.end(), L);
transform(data.begin(), data.end(), data.begin(), bind2nd(divides<complex<double> >(), newnorm) );
/*
data=S3(deltat, omega[3], g, data, potential, in, out);
newnorm = L_2_norm(data.begin(), data.end(), L);
transform(data.begin(), data.end(), data.begin(), bind2nd(divides<complex<double> >(), newnorm) );
data=S3(deltat, omega[2], g, data, potential, in, out);
newnorm = L_2_norm(data.begin(), data.end(), L);
transform(data.begin(), data.end(), data.begin(), bind2nd(divides<complex<double> >(), newnorm) );
data=S3(deltat, omega[1], g, data, potential, in, out);
newnorm = L_2_norm(data.begin(), data.end(), L);
transform(data.begin(), data.end(), data.begin(), bind2nd(divides<complex<double> >(), newnorm) );
data=S3(deltat, omega[0], g, data, potential, in, out);
newnorm = L_2_norm(data.begin(), data.end(), L);
transform(data.begin(), data.end(), data.begin(), bind2nd(divides<complex<double> >(), newnorm) );
data=S3(deltat, omega[1], g, data, potential, in, out);
newnorm = L_2_norm(data.begin(), data.end(), L);
transform(data.begin(), data.end(), data.begin(), bind2nd(divides<complex<double> >(), newnorm) );
data=S3(deltat, omega[2], g, data, potential, in, out);
newnorm = L_2_norm(data.begin(), data.end(), L);
transform(data.begin(), data.end(), data.begin(), bind2nd(divides<complex<double> >(), newnorm) );
data=S3(deltat, omega[3], g, data, potential, in, out);
newnorm = L_2_norm(data.begin(), data.end(), L);
transform(data.begin(), data.end(), data.begin(), bind2nd(divides<complex<double> >(), newnorm) );
*/
}
for(int i = 0; i<N;++i)
cout<<norm(data[i])<<endl;
cout<<"&"<<endl;
return 0;
}