source: palm/trunk/SOURCE/spectrum.f90 @ 1801

Last change on this file since 1801 was 1787, checked in by raasch, 9 years ago

last commit documented

  • Property svn:keywords set to Id
File size: 20.6 KB
Line 
1!> @file spectrum.f90
2!--------------------------------------------------------------------------------!
3! This file is part of PALM.
4!
5! PALM is free software: you can redistribute it and/or modify it under the terms
6! of the GNU General Public License as published by the Free Software Foundation,
7! either version 3 of the License, or (at your option) any later version.
8!
9! PALM is distributed in the hope that it will be useful, but WITHOUT ANY
10! WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR
11! A PARTICULAR PURPOSE.  See the GNU General Public License for more details.
12!
13! You should have received a copy of the GNU General Public License along with
14! PALM. If not, see <http://www.gnu.org/licenses/>.
15!
16! Copyright 1997-2014 Leibniz Universitaet Hannover
17!--------------------------------------------------------------------------------!
18!
19! Current revisions:
20! -----------------
21!
22!
23! Former revisions:
24! -----------------
25! $Id: spectrum.f90 1787 2016-03-08 06:57:00Z raasch $
26!
27! 1786 2016-03-08 05:49:27Z raasch
28! routine is modularized, filename renamed from calc_spectra to spectrum,
29! privious data module spectrum moved from modules.f90 to here,
30! cpp-direktives for spectra removed, immediate return if no spectra levels are
31! given
32!
33! 1682 2015-10-07 23:56:08Z knoop
34! Code annotations made doxygen readable
35!
36! 1575 2015-03-27 09:56:27Z raasch
37! adjustments for psolver-queries
38!
39! 1511 2014-12-16 15:54:16Z suehring
40! Bugfix concerning spectra normalization
41!
42! 1431 2014-07-15 14:47:17Z suehring
43! Wavenumber-integrated spectra coincide with respective variance.
44!
45! 1342 2014-03-26 17:04:47Z kanani
46! REAL constants defined as wp-kinds
47!
48! 1324 2014-03-21 09:13:16Z suehring
49! Bugfix: nzb_x, nzb_yd, nyn_x, nyn_x, nzt_x, nzt_yd belong to transpose_indices
50!
51! 1320 2014-03-20 08:40:49Z raasch
52! ONLY-attribute added to USE-statements,
53! kind-parameters added to all INTEGER and REAL declaration statements,
54! kinds are defined in new module kinds,
55! revision history before 2012 removed,
56! comment fields (!:) to be used for variable explanations added to
57! all variable declaration statements
58!
59! 1318 2014-03-17 13:35:16Z raasch
60! module interfaces removed
61!
62! 1216 2013-08-26 09:31:42Z raasch
63! resorting of array moved to separate routine resort_for_zx,
64! one argument removed from the transpose_..d routines
65!
66! 1120 2013-04-05 15:11:35Z raasch
67! bugfix: calls of fft_x|y replaced by fft_x|y_1d
68!
69! 1036 2012-10-22 13:43:42Z raasch
70! code put under GPL (PALM 3.9)
71!
72! 1003 2012-09-14 14:35:53Z raasch
73! adjustment of array tend for cases with unequal subdomain sizes removed
74!
75! Revision 1.1  2001/01/05 15:08:07  raasch
76! Initial revision
77!
78!
79! Description:
80! ------------
81!> Calculate horizontal spectra along x and y.
82!> ATTENTION: 1d-decomposition along y still needs improvement, because in that
83!>            case the gridpoint number along z still depends on the PE number
84!>            because transpose_xz has to be used (and possibly also
85!>            transpose_zyd needs modification).
86!------------------------------------------------------------------------------!
87 MODULE spectrum
88
89    USE kinds
90
91    PRIVATE
92
93    CHARACTER (LEN=6),  DIMENSION(1:5) ::  header_char = (/ 'PS(u) ', 'PS(v) ',&
94                                           'PS(w) ', 'PS(pt)', 'PS(q) ' /)
95    CHARACTER (LEN=2),  DIMENSION(10)  ::  spectra_direction = 'x'
96    CHARACTER (LEN=10), DIMENSION(10)  ::  data_output_sp  = ' '
97    CHARACTER (LEN=25), DIMENSION(1:5) ::  utext_char =                    &
98                                           (/ '-power spectrum of u     ', &
99                                              '-power spectrum of v     ', &
100                                              '-power spectrum of w     ', &
101                                              '-power spectrum of ^1185 ', &
102                                              '-power spectrum of q     ' /)
103    CHARACTER (LEN=39), DIMENSION(1:5) ::  ytext_char =                        &
104                                 (/ 'k ^2236 ^2566^2569<u(k) in m>2s>->2    ', &
105                                    'k ^2236 ^2566^2569<v(k) in m>2s>->2    ', &
106                                    'k ^2236 ^2566^2569<w(k) in m>2s>->2    ', &
107                                    'k ^2236 ^2566^2569<^1185(k) in m>2s>->2', &
108                                    'k ^2236 ^2566^2569<q(k) in m>2s>->2    ' /)
109
110    INTEGER(iwp) ::  klist_x = 0, klist_y = 0, n_sp_x = 0, n_sp_y = 0
111
112    INTEGER(iwp) ::  comp_spectra_level(100) = 999999,                   &
113                     lstyles(100) = (/ 0, 7, 3, 10, 1, 4, 9, 2, 6, 8,    &
114                                       0, 7, 3, 10, 1, 4, 9, 2, 6, 8,    &
115                                       0, 7, 3, 10, 1, 4, 9, 2, 6, 8,    &
116                                       0, 7, 3, 10, 1, 4, 9, 2, 6, 8,    &
117                                       0, 7, 3, 10, 1, 4, 9, 2, 6, 8,    &
118                                       0, 7, 3, 10, 1, 4, 9, 2, 6, 8,    &
119                                       0, 7, 3, 10, 1, 4, 9, 2, 6, 8,    &
120                                       0, 7, 3, 10, 1, 4, 9, 2, 6, 8,    &
121                                       0, 7, 3, 10, 1, 4, 9, 2, 6, 8,    &
122                                       0, 7, 3, 10, 1, 4, 9, 2, 6, 8 /), &
123                     plot_spectra_level(100) = 999999
124
125    REAL(wp)    ::  time_to_start_sp = 0.0_wp
126
127    PUBLIC  comp_spectra_level, data_output_sp, header_char, klist_x, klist_y, &
128            lstyles, n_sp_x, n_sp_y, plot_spectra_level, spectra_direction,    &
129            utext_char, ytext_char
130
131    SAVE
132
133    INTERFACE calc_spectra
134       MODULE PROCEDURE calc_spectra
135    END INTERFACE calc_spectra
136
137    INTERFACE preprocess_spectra
138       MODULE PROCEDURE preprocess_spectra
139    END INTERFACE preprocess_spectra
140
141    INTERFACE calc_spectra_x
142       MODULE PROCEDURE calc_spectra_x
143    END INTERFACE calc_spectra_x
144
145    INTERFACE calc_spectra_y
146       MODULE PROCEDURE calc_spectra_y
147    END INTERFACE calc_spectra_y
148
149    PUBLIC calc_spectra
150
151
152 CONTAINS
153
154    SUBROUTINE calc_spectra
155
156       USE arrays_3d,                                                          &
157           ONLY:  d, tend
158
159       USE control_parameters,                                                 &
160           ONLY:  average_count_sp, bc_lr_cyc, bc_ns_cyc, message_string, psolver
161
162       USE cpulog,                                                             &
163           ONLY:  cpu_log, log_point
164
165       USE fft_xy,                                                             &
166           ONLY:  fft_init
167
168       USE indices,                                                            &
169           ONLY:  nxl, nxr, nyn, nys, nzb, nzt
170
171       USE kinds
172
173       USE pegrid,                                                             &
174           ONLY:  myid, pdims
175
176       IMPLICIT NONE
177
178       INTEGER(iwp) ::  m  !<
179       INTEGER(iwp) ::  pr !<
180
181
182!
183!--    Check if user gave any levels for spectra to be calculated
184       IF ( comp_spectra_level(1) == 999999 )  RETURN
185
186       CALL cpu_log( log_point(30), 'calc_spectra', 'start' )
187
188!
189!--    Initialize ffts
190       CALL fft_init
191
192!
193!--    Reallocate array d in required size
194       IF ( psolver(1:9) == 'multigrid' )  THEN
195          DEALLOCATE( d )
196          ALLOCATE( d(nzb+1:nzt,nys:nyn,nxl:nxr) )
197       ENDIF
198
199       m = 1
200       DO WHILE ( data_output_sp(m) /= ' '  .AND.  m <= 10 )
201!
202!--       Transposition from z --> x  ( y --> x in case of a 1d-decomposition
203!--       along x)
204          IF ( INDEX( spectra_direction(m), 'x' ) /= 0 )  THEN
205
206!
207!--          Calculation of spectra works for cyclic boundary conditions only
208             IF ( .NOT. bc_lr_cyc )  THEN
209
210                message_string = 'non-cyclic lateral boundaries along x do'//  &
211                                 ' not & allow calculation of spectra along x'
212                CALL message( 'calc_spectra', 'PA0160', 1, 2, 0, 6, 0 )
213             ENDIF
214
215             CALL preprocess_spectra( m, pr )
216
217#if defined( __parallel )
218             IF ( pdims(2) /= 1 )  THEN
219                CALL resort_for_zx( d, tend )
220                CALL transpose_zx( tend, d )
221             ELSE
222                CALL transpose_yxd( d, d )
223             ENDIF
224             CALL calc_spectra_x( d, pr, m )
225#else
226             message_string = 'sorry, calculation of spectra in non paral' //  &
227                              'lel mode& is still not realized'
228             CALL message( 'calc_spectra', 'PA0161', 1, 2, 0, 6, 0 )
229#endif
230
231          ENDIF
232
233!
234!--       Transposition from z --> y (d is rearranged only in case of a
235!--       1d-decomposition along x)
236          IF ( INDEX( spectra_direction(m), 'y' ) /= 0 )  THEN
237
238!
239!--          Calculation of spectra works for cyclic boundary conditions only
240             IF ( .NOT. bc_ns_cyc )  THEN
241                IF ( myid == 0 )  THEN
242                   message_string = 'non-cyclic lateral boundaries along y' // &
243                                    ' do not & allow calculation of spectr' // &
244                                    'a along y'
245                   CALL message( 'calc_spectra', 'PA0162', 1, 2, 0, 6, 0 )
246                ENDIF
247                CALL local_stop
248             ENDIF
249
250             CALL preprocess_spectra( m, pr )
251
252#if defined( __parallel )
253             CALL transpose_zyd( d, d )
254             CALL calc_spectra_y( d, pr, m )
255#else
256             message_string = 'sorry, calculation of spectra in non paral' //  &
257                              'lel mode& is still not realized'
258             CALL message( 'calc_spectra', 'PA0161', 1, 2, 0, 6, 0 )
259#endif
260
261          ENDIF
262
263!
264!--       Increase counter for next spectrum
265          m = m + 1
266         
267       ENDDO
268
269!
270!--    Increase counter for averaging process in routine plot_spectra
271       average_count_sp = average_count_sp + 1
272
273       CALL cpu_log( log_point(30), 'calc_spectra', 'stop' )
274
275    END SUBROUTINE calc_spectra
276
277
278!------------------------------------------------------------------------------!
279! Description:
280! ------------
281!> @todo Missing subroutine description.
282!------------------------------------------------------------------------------!
283    SUBROUTINE preprocess_spectra( m, pr )
284
285       USE arrays_3d,                                                          &
286           ONLY:  d, pt, q, u, v, w
287
288       USE indices,                                                            &
289           ONLY:  ngp_2dh, nxl, nxr, nyn, nys, nzb, nzt
290
291       USE kinds
292
293#if defined( __lc )
294       USE MPI
295#else
296       INCLUDE "mpif.h"
297#endif
298       USE pegrid,                                                             &
299           ONLY:  collective_wait, comm2d, ierr
300
301       USE statistics,                                                         &
302           ONLY:  hom, var_d
303
304
305       IMPLICIT NONE
306
307       INTEGER(iwp) :: i  !<
308       INTEGER(iwp) :: j  !<
309       INTEGER(iwp) :: k  !<
310       INTEGER(iwp) :: m  !<
311       INTEGER(iwp) :: pr !<
312
313       REAL(wp), DIMENSION(nzb:nzt+1) :: var_d_l
314
315       SELECT CASE ( TRIM( data_output_sp(m) ) )
316         
317       CASE ( 'u' )
318          pr = 1
319          d(nzb+1:nzt,nys:nyn,nxl:nxr) = u(nzb+1:nzt,nys:nyn,nxl:nxr)
320       
321       CASE ( 'v' )
322          pr = 2
323          d(nzb+1:nzt,nys:nyn,nxl:nxr) = v(nzb+1:nzt,nys:nyn,nxl:nxr)
324       
325       CASE ( 'w' )
326          pr = 3
327          d(nzb+1:nzt,nys:nyn,nxl:nxr) = w(nzb+1:nzt,nys:nyn,nxl:nxr)
328       
329       CASE ( 'pt' )
330          pr = 4
331          d(nzb+1:nzt,nys:nyn,nxl:nxr) = pt(nzb+1:nzt,nys:nyn,nxl:nxr)
332       
333       CASE ( 'q' )
334          pr = 41
335          d(nzb+1:nzt,nys:nyn,nxl:nxr) = q(nzb+1:nzt,nys:nyn,nxl:nxr)
336       
337       CASE DEFAULT
338!
339!--       The DEFAULT case is reached either if the parameter data_output_sp(m)
340!--       contains a wrong character string or if the user has coded a special
341!--       case in the user interface. There, the subroutine user_spectra
342!--       checks which of these two conditions applies.
343          CALL user_spectra( 'preprocess', m, pr )
344         
345       END SELECT
346
347!
348!--    Subtract horizontal mean from the array, for which spectra have to be
349!--    calculated
350       var_d_l(:) = 0.0_wp
351       DO  i = nxl, nxr
352          DO  j = nys, nyn
353             DO  k = nzb+1, nzt
354                d(k,j,i)   = d(k,j,i) - hom(k,1,pr,0)
355                var_d_l(k) = var_d_l(k) + d(k,j,i) * d(k,j,i)
356             ENDDO
357          ENDDO
358       ENDDO
359!
360!--    Compute total variance from local variances
361       var_d(:) = 0.0_wp
362#if defined( __parallel ) 
363       IF ( collective_wait )  CALL MPI_BARRIER( comm2d, ierr )
364       CALL MPI_ALLREDUCE( var_d_l(0), var_d(0), nzt+1-nzb, MPI_REAL, MPI_SUM, &
365                           comm2d, ierr )
366#else
367       var_d(:) = var_d_l(:)
368#endif
369       var_d(:) = var_d(:) / ngp_2dh(0)
370
371    END SUBROUTINE preprocess_spectra
372
373
374!------------------------------------------------------------------------------!
375! Description:
376! ------------
377!> @todo Missing subroutine description.
378!------------------------------------------------------------------------------!
379    SUBROUTINE calc_spectra_x( ddd, pr, m )
380
381       USE control_parameters,                                                 &
382           ONLY:  fft_method
383
384       USE fft_xy,                                                             &
385           ONLY:  fft_x_1d
386
387       USE grid_variables,                                                     &
388           ONLY:  dx
389
390       USE indices,                                                            &
391           ONLY:  nx, ny
392
393       USE kinds
394
395#if defined( __lc )
396       USE MPI
397#else
398       INCLUDE "mpif.h"
399#endif
400       USE pegrid,                                                             &
401           ONLY:  comm2d, ierr, myid
402
403       USE statistics,                                                         &
404           ONLY:  spectrum_x, var_d
405
406       USE transpose_indices,                                                  &
407           ONLY:  nyn_x, nys_x, nzb_x, nzt_x
408
409
410       IMPLICIT NONE
411
412       INTEGER(iwp) ::  i         !<
413       INTEGER(iwp) ::  ishape(1) !<
414       INTEGER(iwp) ::  j         !<
415       INTEGER(iwp) ::  k         !<
416       INTEGER(iwp) ::  m         !<
417       INTEGER(iwp) ::  n         !<
418       INTEGER(iwp) ::  pr        !<
419
420       REAL(wp) ::  exponent     !<
421       REAL(wp) ::  sum_spec_dum !< wavenumber-integrated spectrum
422   
423       REAL(wp), DIMENSION(0:nx) ::  work !<
424   
425       REAL(wp), DIMENSION(0:nx/2) ::  sums_spectra_l !<
426   
427       REAL(wp), DIMENSION(0:nx/2,100) ::  sums_spectra !<
428   
429       REAL(wp), DIMENSION(0:nx,nys_x:nyn_x,nzb_x:nzt_x) ::  ddd !<
430
431!
432!--    Exponent for geometric average
433       exponent = 1.0_wp / ( ny + 1.0_wp )
434
435!
436!--    Loop over all levels defined by the user
437       n = 1
438       DO WHILE ( comp_spectra_level(n) /= 999999  .AND.  n <= 100 )
439
440          k = comp_spectra_level(n)
441
442!
443!--       Calculate FFT only if the corresponding level is situated on this PE
444          IF ( k >= nzb_x  .AND.  k <= nzt_x )  THEN
445         
446             DO  j = nys_x, nyn_x
447
448                work = ddd(0:nx,j,k)
449                CALL fft_x_1d( work, 'forward' )
450
451                ddd(0,j,k) = dx * work(0)**2
452                DO  i = 1, nx/2
453                   ddd(i,j,k) = dx * ( work(i)**2 + work(nx+1-i)**2 )
454                ENDDO
455
456             ENDDO
457
458!
459!--          Local sum and geometric average of these spectra
460!--          (WARNING: no global sum should be performed, because floating
461!--          point overflow may occur)
462             DO  i = 0, nx/2
463
464                sums_spectra_l(i) = 1.0_wp
465                DO  j = nys_x, nyn_x
466                   sums_spectra_l(i) = sums_spectra_l(i) * ddd(i,j,k)**exponent
467                ENDDO
468
469             ENDDO
470         
471          ELSE
472
473             sums_spectra_l = 1.0_wp
474
475          ENDIF
476
477!
478!--       Global sum of spectra on PE0 (from where they are written on file)
479          sums_spectra(:,n) = 0.0_wp
480#if defined( __parallel )   
481          CALL MPI_BARRIER( comm2d, ierr )  ! Necessary?
482          CALL MPI_REDUCE( sums_spectra_l(0), sums_spectra(0,n), nx/2+1,       &
483                           MPI_REAL, MPI_PROD, 0, comm2d, ierr )
484#else
485          sums_spectra(:,n) = sums_spectra_l
486#endif
487!
488!--       Normalize spectra by variance
489          sum_spec_dum = SUM( sums_spectra(:,n) )
490          IF ( sum_spec_dum /= 0.0_wp )  THEN
491             sums_spectra(:,n) = sums_spectra(:,n) * var_d(k) / sum_spec_dum
492          ENDIF
493          n = n + 1
494
495       ENDDO
496       n = n - 1
497
498       IF ( myid == 0 )  THEN
499!
500!--       Sum of spectra for later averaging (see routine data_output_spectra)
501          DO  i = 1, nx/2
502             DO k = 1, n
503                spectrum_x(i,k,m) = spectrum_x(i,k,m) + sums_spectra(i,k)
504             ENDDO
505          ENDDO
506
507       ENDIF
508!
509!--    n_sp_x is needed by data_output_spectra_x
510       n_sp_x = n
511
512    END SUBROUTINE calc_spectra_x
513
514
515!------------------------------------------------------------------------------!
516! Description:
517! ------------
518!> @todo Missing subroutine description.
519!------------------------------------------------------------------------------!
520    SUBROUTINE calc_spectra_y( ddd, pr, m )
521
522       USE control_parameters,                                                 &
523           ONLY:  fft_method
524
525       USE fft_xy,                                                             &
526           ONLY:  fft_y_1d
527
528       USE grid_variables,                                                     &
529           ONLY:  dy
530
531       USE indices,                                                            &
532           ONLY:  nx, ny
533
534       USE kinds
535
536#if defined( __lc )
537       USE MPI
538#else
539       INCLUDE "mpif.h"
540#endif
541       USE pegrid,                                                             &
542           ONLY:  comm2d, ierr, myid
543
544       USE statistics,                                                         &
545           ONLY:  spectrum_y, var_d
546
547       USE transpose_indices,                                                  &
548           ONLY:  nxl_yd, nxr_yd, nzb_yd, nzt_yd
549
550
551       IMPLICIT NONE
552
553       INTEGER(iwp) ::  i         !<
554       INTEGER(iwp) ::  j         !<
555       INTEGER(iwp) ::  jshape(1) !<
556       INTEGER(iwp) ::  k         !<
557       INTEGER(iwp) ::  m         !<
558       INTEGER(iwp) ::  n         !<
559       INTEGER(iwp) ::  pr        !<
560
561       REAL(wp) ::  exponent !<
562       REAL(wp) ::  sum_spec_dum !< wavenumber-integrated spectrum
563   
564       REAL(wp), DIMENSION(0:ny) ::  work !<
565   
566       REAL(wp), DIMENSION(0:ny/2) ::  sums_spectra_l !<
567   
568       REAL(wp), DIMENSION(0:ny/2,100) ::  sums_spectra !<
569   
570       REAL(wp), DIMENSION(0:ny,nxl_yd:nxr_yd,nzb_yd:nzt_yd) :: ddd !<
571
572
573!
574!--    Exponent for geometric average
575       exponent = 1.0_wp / ( nx + 1.0_wp )
576
577!
578!--    Loop over all levels defined by the user
579       n = 1
580       DO WHILE ( comp_spectra_level(n) /= 999999  .AND.  n <= 100 )
581
582          k = comp_spectra_level(n)
583
584!
585!--       Calculate FFT only if the corresponding level is situated on this PE
586          IF ( k >= nzb_yd  .AND.  k <= nzt_yd )  THEN
587         
588             DO  i = nxl_yd, nxr_yd
589
590                work = ddd(0:ny,i,k)
591                CALL fft_y_1d( work, 'forward' )
592
593                ddd(0,i,k) = dy * work(0)**2
594                DO  j = 1, ny/2
595                   ddd(j,i,k) = dy * ( work(j)**2 + work(ny+1-j)**2 )
596                ENDDO
597
598             ENDDO
599
600!
601!--          Local sum and geometric average of these spectra
602!--          (WARNING: no global sum should be performed, because floating
603!--          point overflow may occur)
604             DO  j = 0, ny/2
605
606                sums_spectra_l(j) = 1.0_wp
607                DO  i = nxl_yd, nxr_yd
608                   sums_spectra_l(j) = sums_spectra_l(j) * ddd(j,i,k)**exponent
609                ENDDO
610
611             ENDDO
612         
613          ELSE
614
615             sums_spectra_l = 1.0_wp
616
617          ENDIF
618
619!
620!--       Global sum of spectra on PE0 (from where they are written on file)
621          sums_spectra(:,n) = 0.0_wp
622#if defined( __parallel )   
623          CALL MPI_BARRIER( comm2d, ierr )  ! Necessary?
624          CALL MPI_REDUCE( sums_spectra_l(0), sums_spectra(0,n), ny/2+1,       &
625                           MPI_REAL, MPI_PROD, 0, comm2d, ierr )
626#else
627          sums_spectra(:,n) = sums_spectra_l
628#endif
629!
630!--       Normalize spectra by variance
631          sum_spec_dum = SUM( sums_spectra(:,n) )
632          IF ( SUM(sums_spectra(:,n)) /= 0.0_wp )  THEN
633             sums_spectra(:,n) = sums_spectra(:,n) *                           &
634                                 var_d(k) / SUM(sums_spectra(:,n))
635          ENDIF
636          n = n + 1
637
638       ENDDO
639       n = n - 1
640
641
642       IF ( myid == 0 )  THEN
643!
644!--       Sum of spectra for later averaging (see routine data_output_spectra)
645          DO  j = 1, ny/2
646             DO k = 1, n
647                spectrum_y(j,k,m) = spectrum_y(j,k,m) + sums_spectra(j,k)
648             ENDDO
649          ENDDO
650
651       ENDIF
652
653!
654!--    n_sp_y is needed by data_output_spectra_y
655       n_sp_y = n
656
657    END SUBROUTINE calc_spectra_y
658
659 END MODULE spectrum
Note: See TracBrowser for help on using the repository browser.