-
Notifications
You must be signed in to change notification settings - Fork 0
/
PMP2mod_power.f90
417 lines (396 loc) · 15 KB
/
PMP2mod_power.f90
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
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
!-------------------------------------------------
!
! Routines to estimate Power Spectrum
!
!-------------------------------------------------
Module Power
use Tools
use fft5
Real*4, Allocatable, Dimension(:,:) :: Pk,dNharm,kharm
Integer*4 :: iCIC=1
Contains
!----------------------------------------------------
! power spectrum P(k)
!
! iFlag =0 - only real space
! =1 - real + 3 projections P0 + 3 projections P2
!
SUBROUTINE GetPower(iFlag)
!
!----------------------------------------------------
use Density
Integer*4 :: moment
Real*4 :: kNyq
Character*120 :: Name
Mem_current = Memory(INT(21*NGRID,8))
Allocate(Pk(NGRID,0:6),dNharm(NGRID,0:6),kharm(NGRID,0:6))
moment = INT(100.*(1./AEXPN-1.)+0.5) ! make integer out of Z_moment
Pk(:,:) = 0.
dNharm(:,:) = 0.
kharm(:,:) = 0.
boxsize = Box
kNyq = 2.*3.1415926/Box*NGRID/2.
If(iFlag==0)Then ! --- only real space Pk
write(Name,'(2(a,i4.4),3(a,i3.3))')'PowerDM.R.',moment, &
'.',Nrealization,'.dat'
OPEN(18,FILE=TRIM(Name),STATUS='UNKNOWN')
write(18,'(a)')HEADER
write(18,'(a,f7.4,a,i4,a,f8.3,a,i4,a,f8.2)')' Aexpn =',AEXPN,' Step=',ISTEP, &
' Redshift= ',1./AEXPN-1.,' Ngrid= ',Ngrid,' kNyq= ',kNyq
CALL POWERfft5(0,0)
WRITE (18,'(3x,a,T9,a,T19,a,T30,a,T43,a)') 'bin','k/Mpch', &
'N','Preal'
DO I=1,NGRID/2
WRITE (18,'(i6,f9.5,es11.4,7es13.5)')I,kharm(I,0),dNharm(I,0),Pk(I,0)
EndDo
Else ! --- Real + redshift space
!-------------- change to smaller grid
Mem_current = Memory(-1_8*NGRID*NGRID*NGRID)
Deallocate (FI)
NGRID_old = NGRID ! store old value of NGRID
NGRID = NGRID/2
Mem_current = Memory(1_8*NGRID*NGRID*NGRID)
ALLOCATE(FI(NGRID,NGRID,NGRID))
write(Name,'(2(a,i4.4),3(a,i3.3))')'PowerDM.Z.',moment,'.',Nrealization, &
'.V',INT(sigv),'.d',INT(DensThr),'.dat'
OPEN(18,FILE=TRIM(Name),STATUS='UNKNOWN')
write(18,'(a)')HEADER
write(18,'(a,f7.4,a,i4,a,f8.3,a,i4,a,f8.2)')' Aexpn =',AEXPN,' Step=',ISTEP, &
' Redshift= ',1./AEXPN-1.,' Ngrid= ',Ngrid,' kNyq= ',kNyq
Do iSwitch = 0, 3
Call DENSITrsd(iSwitch,NGRID_old) ! find density for all particles in PMfiles
CALL POWERfft5(0,iSwitch) ! Power spectrum monopole
If(iSwitch/=0)Then
Call DENSITrsd(iSwitch,NGRID_old) ! find density for all particles in PMfiles
CALL POWERfft5(1,iSwitch) ! Power spectrum quadrupole
End If
EndDo
WRITE (18,'(3x,a,T9,a,T19,a,T30,a,T43,a,T56,a)') 'bin','k/Mpch', &
'N','Preal','P0_redshift','P2_redshift'
DO I=1,NGRID/2
P2 = (Pk(I,1)+ Pk(I,3)+Pk(I,5))/3.
P0 = (Pk(I,2)+ Pk(I,4)+Pk(I,6))/3.
!P2 = (Pk(I,1)+ Pk(I,3))/2.
!P0 = (Pk(I,2)+ Pk(I,4))/2.
WRITE (18,110) I,kharm(I,0),dNharm(I,0),Pk(I,0),P0,P2 !(Pk(I,j),j=0,6)
EndDo
110 format(i6,f9.5,es11.4,7es13.5)
Mem_current = Memory(-1_8*NGRID*NGRID*NGRID)
Deallocate (FI)
NGRID = NGRID_old ! restore the old value of NGRID
Mem_current = Memory(1_8*NGRID*NGRID*NGRID)
ALLOCATE(FI(NGRID,NGRID,NGRID))
End If
Mem_current = Memory(-INT(21*NGRID,8))
DeAllocate(Pk,dNharm,kharm)
close(18)
end SUBROUTINE GetPower
!----------------------------------------------------
! power spectrum P(k) for tracers
!
! iFlag =0 - only real space
! =1 - real + 3 projections P0 + 3 projections P2
!
SUBROUTINE GetPowerGal(iFlag)
!
!----------------------------------------------------
use Density
Integer*4 :: moment
Real*4 :: kNyq
Character*120 :: Name
Mem_current = Memory(INT(21*NGRID,8))
Allocate(Pk(NGRID,0:6),dNharm(NGRID,0:6),kharm(NGRID,0:6))
moment = INT(100.*(1./AEXPN-1.)+0.5) ! make integer out of Z_moment
Pk(:,:) = 0.
dNharm(:,:) = 0.
kharm(:,:) = 0.
boxsize = Box
kNyq = 2.*3.1415926/Box*NGRID/2.
If(iFlag==0)Then ! --- only real space Pk
write(Name,'(2(a,i4.4),3(a,i3.3))')'PowerGal.R.',moment, &
'.',Nrealization,'.dat'
OPEN(18,FILE=TRIM(Name),STATUS='UNKNOWN')
write(18,'(a)')HEADER
write(18,'(a,f7.4,a,i4,a,f8.3,a,i4,a,f8.2,a,i10)')' Aexpn =',AEXPN,' Step=',ISTEP, &
' Redshift= ',1./AEXPN-1.,' Ngrid= ',Ngrid,' kNyq= ',kNyq
write(18,'(a,i10,a,f8.3,a,es13.4,a,f8.3)'),' Ngalaxies=',Ngalaxies, &
' DensThresh= ',BiasPars(1), &
' Normalize = ',BiasPars(2), &
' Slope = ',BiasPars(3)
CALL POWERfft5(0,0)
WRITE (18,'(3x,a,T9,a,T19,a,T30,a,T43,a)') 'bin','k/Mpch', &
'N','Preal'
DO I=1,NGRID/2
WRITE (18,'(i6,f9.5,es11.4,7es13.5)')I,kharm(I,0),dNharm(I,0),Pk(I,0)
EndDo
Else ! --- Real + redshift space
!-------------- change to smaller grid
Mem_current = Memory(-1_8*NGRID*NGRID*NGRID)
Deallocate (FI)
NGRID_old = NGRID ! store old value of NGRID
NGRID = NGRID/2
Mem_current = Memory(1_8*NGRID*NGRID*NGRID)
ALLOCATE(FI(NGRID,NGRID,NGRID))
write(Name,'(2(a,i4.4),3(a,i3.3))')'PowerGal.Z.',moment,'.',Nrealization, &
'.dat'
OPEN(18,FILE=TRIM(Name),STATUS='UNKNOWN')
write(18,'(a)')HEADER
write(18,'(a,f7.4,a,i4,a,f8.3,a,i4,a,f8.2,a,i10)')' Aexpn =',AEXPN,' Step=',ISTEP, &
' Redshift= ',1./AEXPN-1.,' Ngrid= ',Ngrid,' kNyq= ',kNyq
write(18,'(a,i10,a,f8.3,a,es13.4,a,f8.3)'),' Ngalaxies=',Ngalaxies, &
' DensThresh= ',BiasPars(1), &
' Normalize = ',BiasPars(2), &
' Slope = ',BiasPars(3)
Do iSwitch = 0, 3
Call DensGal(iSwitch) ! find density for all particles in PMfiles
CALL POWERfft5(0,iSwitch) ! Power spectrum monopole
If(iSwitch/=0)Then
Call DensGal(iSwitch) ! find density for all particles in PMfiles
CALL POWERfft5(1,iSwitch) ! Power spectrum quadrupole
End If
EndDo
WRITE (18,'(3x,a,T9,a,T19,a,T30,a,T43,a,T56,a)') 'bin','k/Mpch', &
'N','Preal','P0_redshift','P2_redshift'
DO I=1,NGRID/2
P2 = (Pk(I,1)+ Pk(I,3)+Pk(I,5))/3.
P0 = (Pk(I,2)+ Pk(I,4)+Pk(I,6))/3.
!P2 = (Pk(I,1)+ Pk(I,3))/2.
!P0 = (Pk(I,2)+ Pk(I,4))/2.
WRITE (18,110) I,kharm(I,0),dNharm(I,0),Pk(I,0),P0,P2 !(Pk(I,j),j=0,6)
EndDo
110 format(i6,f9.5,es11.4,7es13.5)
Mem_current = Memory(-1_8*NGRID*NGRID*NGRID)
Deallocate (FI)
NGRID = NGRID_old ! restore the old value of NGRID
Mem_current = Memory(1_8*NGRID*NGRID*NGRID)
ALLOCATE(FI(NGRID,NGRID,NGRID))
End If
Mem_current = Memory(-INT(21*NGRID,8))
DeAllocate(Pk,dNharm,kharm)
close(18)
end SUBROUTINE GetPowerGal
!----------------------------------------------------
! power spectrum P(k)
! for given density field FI
SUBROUTINE POWERfft5(iPol,iSwitch)
!----------------------------------------------------
PARAMETER ( Pi=3.141592653)
Real*8, Allocatable,DIMENSION(:) :: DPOWER,dk
Real*8, Allocatable,DIMENSION(:) :: iVolume
Real*8, Allocatable,DIMENSION(:,:) :: DPW,dkk
Real*8, Allocatable,DIMENSION(:,:) :: iVol
integer*4, parameter :: Nlensav = 16384
integer*4, parameter :: Nlenwrk = 16384
real*8,SAVE :: work(1:Nlenwrk)
real*8,SAVE :: wsave(1:Nlensav)
real*8 :: r(Nlenwrk)
REAL*8 :: SIGMA,k_Ny,PkCorrec,Pkk,dens
Integer*4 :: OMP_GET_MAX_THREADS,OMP_GET_THREAD_NUM
!$OMP THREADPRIVATE(work,wsave)
Call Timing(4,-1)
lensav = Ngrid+int(log(real(Ngrid,kind = 4))/log(2.0E+00))+4
lenwrk = Ngrid
NPOWER =NGRID
boxsize = Box
iThreads = OMP_GET_MAX_THREADS()
ALLOCATE( DPOWER(NPOWER),iVolume(NPOWER),dk(NPOWER))
ALLOCATE( DPW(NPOWER,iThreads),iVol(NPOWER,iThreads), &
dkk(NPOWER,iThreads))
k_Ny = 2.*Pi/boxsize*NGrid/2.
dens = Nparticles/boxsize**3
NSPEC = NGRID/2
If(Boxsize.le.0..or.Boxsize.gt.1.e+4)Then
write (*,*) ' Wrong Box size=',Boxsize,' STOP'
Return
EndIf
SIGMA =0.
Scale = 2.*Pi/Boxsize
DO M =1,NGRID
iVolume(M) =0
DPOWER(M) =0.
dk(M) =0.
Do i=1,iThreads
iVol(M,i) =0
DPW(M,i) =0.
dkk(M,i) =0.
EndDo
ENDDO
call rfft1i ( Ngrid, wsave, lensav, ier ) ! Initialize FFT
inc = 1
lenr = Ngrid
write(*,*) ' time =',seconds()/60., ' go to first fft'
!$OMP PARALLEL DO DEFAULT(SHARED) copyin(wsave,work) &
!$OMP PRIVATE ( k,j,i ,r, ier)
Do k=1,NGRID ! ------ x-direction
Do j=1,NGRID
Do i=1,NGRID
r(i) = FI(i,j,k)
EndDo
call rfft1f ( Ngrid, inc, r, lenr, wsave, lensav, work, lenwrk, ier )
Do i=1,NGRID
FI(i,j,k) = r(i)
EndDo
EndDo
EndDo
Do k=1,NGRID
Do j=1,NGRID
FI(1,j,k) = FI(1,j,k)*sqrt(2.)
FI(NGRID,j,k) = FI(NGRID,j,k)*sqrt(2.)
end Do
end Do
write(*,*) ' time =',seconds()/60., ' go to second fft'
!$OMP PARALLEL DO DEFAULT(SHARED) copyin(wsave,work) &
!$OMP PRIVATE (k, j,i ,r, ier )
Do k=1,NGRID ! ------ y-direction
Do i=1,NGRID
Do j=1,NGRID
r(j) = FI(i,j,k)
EndDo
call rfft1f ( Ngrid, inc, r, lenr, wsave, lensav, work, lenwrk, ier )
Do j=1,NGRID
FI(i,j,k) = r(j)
EndDo
EndDo
EndDo
Do k=1,NGRID
Do i=1,NGRID
FI(i,1,k) = FI(i,1,k)*sqrt(2.)
FI(i,NGRID,k) = FI(i,NGRID,k)*sqrt(2.)
end Do
end Do
write(*,*) ' time =',seconds()/60., ' go to transposition'
!$OMP PARALLEL DO DEFAULT(SHARED) &
!$OMP PRIVATE ( k,j,i,aa)
DO J=1,Ngrid
DO K=1,Ngrid-1
DO I=K+1,Ngrid
aa = FI(I,J,K)
FI(I,J,K) =FI(K,J,I)
FI(K,J,I) =aa
ENDDO
ENDDO
ENDDO
write(*,*) ' time =',seconds()/60., ' go to third fft'
!$OMP PARALLEL DO DEFAULT(SHARED) copyin(wsave,work) &
!$OMP PRIVATE ( k,j,i ,r, ier)
Do j=1,NGRID ! ------ z-direction
Do i=1,NGRID
Do k=1,NGRID
r(k) = FI(k,i,j)
EndDo
call rfft1f ( Ngrid, inc, r, lenr, wsave, lensav, work, lenwrk, ier )
Do k=1,NGRID
FI(k,i,j) = r(k)
EndDo
EndDo
EndDo
Do j=1,NGRID
Do i=1,NGRID
FI(1,i,j) = FI(1,i,j)*sqrt(2.)
FI(NGRID,i,j) = FI(NGRID,i,j)*sqrt(2.)
end Do
end Do
write(*,*) ' time =',seconds()/60.,' Normalize '
!$OMP PARALLEL DO DEFAULT(SHARED) &
!$OMP PRIVATE ( k,j,i)
Do k=1,NGRID
Do j=1,NGRID
Do i=1,NGRID
FI(i,j,k) = FI(i,j,k)**2/8.
end Do
end Do
end Do
FI(1,1,1) =0.
! kStep =NGRID/iThreads
! If(kstep*iThreads.ne.NGRID)STOP '-- Wrong number of iThreads'
write (*,*) ' Power .... threads =',iThreads
!$OMP PARALLEL DO DEFAULT(SHARED) &
!$OMP PRIVATE ( k,j,i,MI,MJ,MK,AMP,RK,MP,iOMP,MP2,d1,d2,w1,w2,cosmu2,P2) &
!$OMP REDUCTION(+:SIGMA)
Do K=1,NGRID
iOMP = OMP_GET_THREAD_NUM()+1
!iOMP = 1
MK = (K/2)**2
!write(*,*) K,K/2,MK
DO J=1,NGRID
MJ = (J/2)**2
DO I=1,NGRID
MI = (I/2)**2
AMP = FI(I,J,K)
SIGMA = SIGMA +AMP
RK = MI +MJ+MK ! k**2
! Power Spectrum
If(RK.lt.0.5)cycle ! do not take k =0 component
RK =Sqrt(RK)
If(iCIC ==1)Then
MP =max(INT(RK +0.0001),1) ! left boundary
MP2 = min(MP+1,NGRID) ! right boundary
d1 = MP
d2 = d1 + 1.
w1 = MP - RK +1.
w2 = 1. - w1
If(iPol==0)Then
DPW(MP,iOMP) =DPW(MP,iOMP) +AMP*w1
DPW(MP2,iOMP) =DPW(MP2,iOMP) +AMP*w2
Else
cosmu2 = MI/RK**2 ! square of cos()
P2 = 1.5*cosmu2 -0.5 ! Lagrange Pol_2
DPW(MP,iOMP) =DPW(MP,iOMP) +AMP*w1*P2
DPW(MP2,iOMP) =DPW(MP2,iOMP) +AMP*w2*P2
EndIf ! iPol
iVol(MP,iOMP) =iVol(MP,iOMP)+w1
dkk(MP,iOMP) =dkk(MP,iOMP) +d1*w1
iVol(MP2,iOMP) =iVol(MP2,iOMP)+w2
dkk(MP2,iOMP) =dkk(MP2,iOMP) +d2*w2
Else ! NGP in k-space
MP = max(INT(RK +0.5),1) ! middle bin
If(iPol==0)Then
DPW(MP,iOMP) =DPW(MP,iOMP) +AMP
Else
cosmu2 = MI/RK**2 ! square of cos()
P2 = 1.5*cosmu2 -0.5 ! Lagrange Pol_2
DPW(MP,iOMP) =DPW(MP,iOMP) +AMP*P2
EndIf
iVol(MP,iOMP) =iVol(MP,iOMP)+1.
dkk(MP,iOMP) =dkk(MP,iOMP) +RK
End If
ENDDO
ENDDO
ENDDO
Do MP=1,NGRID
Do iP=1, iThreads
DPOWER(MP) =DPOWER(MP) +DPW(MP,iP)
iVolume(MP) =iVolume(MP)+iVol(MP,iP)
dk(MP) =dk(MP) +dkk(MP,iP)
EndDo
EndDo
If(iPol == 1)DPOWER(:) = 2.5*DPOWER(:) !! (2L+1)/2 for quadrupole
WRITE (*,100) SQRT(SIGMA),2.*Pi*(float(Nparticles))**0.3333/boxsize/2., &
2.*Pi/boxsize*NGRID/2.
!WRITE (17,100) SQRT(SIGMA),2.*Pi*(float(Nparticles))**0.3333/boxsize/2., &
! 2.*Pi/boxsize*NGRID/2.
100 FORMAT(8x,'RMS DRho/Rho =',F7.3, &
/8x,'k_Nyquist for Particles =',f6.2, &
/8x,'k_Nyquist for Grid =',f6.2, &
/4X,'n',2X,'k(h/Mpc)',1X, &
'Harmonics ','Power Spectrum(h^3) Pk_corrected')
110 FORMAT(I6,F11.6,1p,g12.4,3X,5G14.5)
ii = 2*iSwitch-iPol
!write(*,*) ' iPol/sw =',iPol,iSwitch,ii
DO I=1,NGRID/2
IF(iVolume(I).GT.0)THEN
wk = dk(I)/iVolume(I)*Scale
Pkk = Boxsize**3*DPOWER(I)/iVolume(I)
PkCorrec = Pkk/(1.-0.6667*(sin(Pi*wk/k_Ny/2.))**2)-1./dens
!if(I< 20)WRITE (*,110) I,wk,iVolume(I),Pkk,PkCorrec
!WRITE (17,110) I,wk,iVolume(I),Pkk,PkCorrec
kharm(I,ii) = wk
dNharm(I,ii) = iVolume(I)
Pk(I,ii) = Pkk !!!PkCorrec
ENDIF
ENDDO
DEALLOCATE( DPOWER,iVolume,dk)
DEALLOCATE( DPW,iVol,dkk)
Call Timing(4,1)
END SUBROUTINE POWERfft5
end Module Power