source: palm/trunk/SOURCE/spectra_mod.f90 @ 4359

Last change on this file since 4359 was 4182, checked in by scharf, 5 years ago
  • corrected "Former revisions" section
  • minor formatting in "Former revisions" section
  • added "Author" section
  • Property svn:keywords set to Id
File size: 26.3 KB
Line 
1!> @file spectra_mod.f90
2!------------------------------------------------------------------------------!
3! This file is part of the PALM model system.
4!
5! PALM is free software: you can redistribute it and/or modify it under the
6! terms of the GNU General Public License as published by the Free Software
7! Foundation, either version 3 of the License, or (at your option) any later
8! version.
9!
10! PALM is distributed in the hope that it will be useful, but WITHOUT ANY
11! WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR
12! A PARTICULAR PURPOSE.  See the GNU General Public License for more details.
13!
14! You should have received a copy of the GNU General Public License along with
15! PALM. If not, see <http://www.gnu.org/licenses/>.
16!
17! Copyright 1997-2019 Leibniz Universitaet Hannover
18!------------------------------------------------------------------------------!
19!
20! Current revisions:
21! -----------------
22!
23!
24! Former revisions:
25! -----------------
26! $Id: spectra_mod.f90 4182 2019-08-22 15:20:23Z suehring $
27! Corrected "Former revisions" section
28!
29! 3805 2019-03-20 15:26:35Z raasch
30! unused variable removed
31!
32! 3655 2019-01-07 16:51:22Z knoop
33! Renamed output variables
34!
35! Revision 1.1  2001/01/05 15:08:07  raasch
36! Initial revision
37!
38!
39! Description:
40! ------------
41!> Calculate horizontal spectra along x and y.
42!> ATTENTION: 1d-decomposition along y still needs improvement, because in that
43!>            case the gridpoint number along z still depends on the PE number
44!>            because transpose_xz has to be used (and possibly also
45!>            transpose_zyd needs modification).
46!------------------------------------------------------------------------------!
47 MODULE spectra_mod
48
49    USE kinds
50
51    PRIVATE
52
53    CHARACTER (LEN=2),  DIMENSION(10) ::  spectra_direction = 'x'
54    CHARACTER (LEN=10), DIMENSION(10) ::  data_output_sp  = ' '
55
56    INTEGER(iwp) ::  average_count_sp = 0
57    INTEGER(iwp) ::  dosp_time_count = 0
58    INTEGER(iwp) ::  n_sp_x = 0, n_sp_y = 0
59
60    INTEGER(iwp) ::  comp_spectra_level(100) = 999999
61
62    LOGICAL ::  calculate_spectra   = .FALSE.  !< internal switch that spectra are calculated
63    LOGICAL ::  spectra_initialized = .FALSE.  !< internal switch that spectra related quantities are initialized
64
65    REAL(wp) ::  averaging_interval_sp = 9999999.9_wp  !< averaging interval for spectra output
66    REAL(wp) ::  dt_dosp = 9999999.9_wp                !< time interval for spectra output
67    REAL(wp) ::  skip_time_dosp = 9999999.9_wp         !< no output of spectra data before this interval has passed
68
69    REAL(wp), DIMENSION(:), ALLOCATABLE ::  var_d
70
71    REAL(wp), DIMENSION(:,:,:), ALLOCATABLE ::  spectrum_x, spectrum_y
72
73    SAVE
74
75    INTERFACE calc_spectra
76       MODULE PROCEDURE calc_spectra
77    END INTERFACE calc_spectra
78
79    INTERFACE preprocess_spectra
80       MODULE PROCEDURE preprocess_spectra
81    END INTERFACE preprocess_spectra
82
83    INTERFACE calc_spectra_x
84       MODULE PROCEDURE calc_spectra_x
85    END INTERFACE calc_spectra_x
86
87    INTERFACE calc_spectra_y
88       MODULE PROCEDURE calc_spectra_y
89    END INTERFACE calc_spectra_y
90
91    INTERFACE spectra_check_parameters
92       MODULE PROCEDURE spectra_check_parameters
93    END INTERFACE spectra_check_parameters
94
95    INTERFACE spectra_header
96       MODULE PROCEDURE spectra_header
97    END INTERFACE spectra_header
98
99    INTERFACE spectra_init
100       MODULE PROCEDURE spectra_init
101    END INTERFACE spectra_init
102
103    INTERFACE spectra_parin
104       MODULE PROCEDURE spectra_parin
105    END INTERFACE spectra_parin
106
107    PUBLIC average_count_sp, averaging_interval_sp, calc_spectra,              &
108           calculate_spectra, comp_spectra_level, data_output_sp,              &
109           dosp_time_count, dt_dosp, n_sp_x, n_sp_y,                           &
110           skip_time_dosp, spectra_check_parameters, spectra_direction,        &
111           spectra_header, spectra_init, spectra_parin, spectrum_x,            &
112           spectrum_y, var_d
113
114
115 CONTAINS
116
117!------------------------------------------------------------------------------!
118! Description:
119! ------------
120!> Parin for &spectra_par for calculating spectra
121!------------------------------------------------------------------------------!
122    SUBROUTINE spectra_parin
123
124       USE control_parameters,                                                 &
125           ONLY:  dt_data_output, message_string
126
127       IMPLICIT NONE
128
129       CHARACTER (LEN=80) ::  line  !< dummy string that contains the current  &
130                                    !< line of the parameter file
131
132       NAMELIST /spectra_par/  averaging_interval_sp, comp_spectra_level,      &
133                               data_output_sp, dt_dosp, skip_time_dosp,        &
134                               spectra_direction
135
136       NAMELIST /spectra_parameters/                                           &
137                               averaging_interval_sp, comp_spectra_level,      &
138                               data_output_sp, dt_dosp, skip_time_dosp,        &
139                               spectra_direction
140!
141!--    Position the namelist-file at the beginning (it was already opened in
142!--    parin), search for the namelist-group of the package and position the
143!--    file at this line.
144       line = ' '
145
146!
147!--    Try to find the spectra package
148       REWIND ( 11 )
149       line = ' '
150       DO WHILE ( INDEX( line, '&spectra_parameters' ) == 0 )
151          READ ( 11, '(A)', END=12 )  line
152       ENDDO
153       BACKSPACE ( 11 )
154
155!
156!--    Read namelist
157       READ ( 11, spectra_parameters, ERR = 10 )
158
159!
160!--    Default setting of dt_dosp here (instead of check_parameters), because
161!--    its current value is needed in init_pegrid
162       IF ( dt_dosp == 9999999.9_wp )  dt_dosp = dt_data_output
163
164!
165!--    Set general switch that spectra shall be calculated
166       calculate_spectra = .TRUE.
167
168       GOTO 14
169
170 10    BACKSPACE( 11 )
171       READ( 11 , '(A)') line
172       CALL parin_fail_message( 'spectra_parameters', line )
173!
174!--    Try to find the old namelist
175 12    REWIND ( 11 )
176       line = ' '
177       DO WHILE ( INDEX( line, '&spectra_par' ) == 0 )
178          READ ( 11, '(A)', END=14 )  line
179       ENDDO
180       BACKSPACE ( 11 )
181
182!
183!--    Read namelist
184       READ ( 11, spectra_par, ERR = 13, END = 14 )
185
186       
187       message_string = 'namelist spectra_par is deprecated and will be ' // &
188                     'removed in near future. Please use namelist ' //       &
189                     'spectra_parameters instead'
190       CALL message( 'spectra_parin', 'PA0487', 0, 1, 0, 6, 0 )
191!
192!--    Default setting of dt_dosp here (instead of check_parameters), because
193!--    its current value is needed in init_pegrid
194       IF ( dt_dosp == 9999999.9_wp )  dt_dosp = dt_data_output
195
196!
197!--    Set general switch that spectra shall be calculated
198       calculate_spectra = .TRUE.
199       
200       GOTO 14
201
202 13    BACKSPACE( 11 )
203       READ( 11 , '(A)') line
204       CALL parin_fail_message( 'spectra_par', line )
205       
206       
207 14    CONTINUE
208
209    END SUBROUTINE spectra_parin
210
211
212
213!------------------------------------------------------------------------------!
214! Description:
215! ------------
216!> Initialization of spectra related variables
217!------------------------------------------------------------------------------!
218    SUBROUTINE spectra_init
219
220       USE indices,                                                            &
221           ONLY:  nx, ny, nzb, nzt
222
223       IMPLICIT NONE
224
225       IF ( spectra_initialized )  RETURN
226
227       IF ( dt_dosp /= 9999999.9_wp )  THEN
228
229          IF ( .NOT. ALLOCATED( spectrum_x ) )  THEN
230             ALLOCATE( spectrum_x( 1:nx/2, 1:100, 1:10 ) )
231             spectrum_x = 0.0_wp
232          ENDIF
233
234          IF ( .NOT. ALLOCATED( spectrum_y ) )  THEN
235             ALLOCATE( spectrum_y( 1:ny/2, 1:100, 1:10 ) )
236             spectrum_y = 0.0_wp
237          ENDIF
238
239          ALLOCATE( var_d(nzb:nzt+1) )
240          var_d = 0.0_wp
241       ENDIF
242
243       spectra_initialized = .TRUE.
244
245    END SUBROUTINE spectra_init
246
247
248
249!------------------------------------------------------------------------------!
250! Description:
251! ------------
252!> Check spectra related quantities
253!------------------------------------------------------------------------------!
254    SUBROUTINE spectra_check_parameters
255
256       USE control_parameters,                                                 &
257           ONLY:  averaging_interval, message_string, skip_time_data_output
258
259       IMPLICIT NONE
260
261!
262!--    Check the average interval
263       IF ( averaging_interval_sp == 9999999.9_wp )  THEN
264          averaging_interval_sp = averaging_interval
265       ENDIF
266
267       IF ( averaging_interval_sp > dt_dosp )  THEN
268          WRITE( message_string, * )  'averaging_interval_sp = ',              &
269                averaging_interval_sp, ' must be <= dt_dosp = ', dt_dosp
270          CALL message( 'spectra_check_parameters', 'PA0087', 1, 2, 0, 6, 0 )
271       ENDIF
272
273!
274!--    Set the default skip time interval for data output, if necessary
275       IF ( skip_time_dosp == 9999999.9_wp )                                   &
276                                          skip_time_dosp = skip_time_data_output
277
278    END SUBROUTINE spectra_check_parameters
279
280
281
282!------------------------------------------------------------------------------!
283! Description:
284! ------------
285!> Header output for spectra
286!>
287!> @todo Output of netcdf data format and compression level
288!------------------------------------------------------------------------------!
289    SUBROUTINE spectra_header ( io )
290
291       USE control_parameters,                                                 &
292           ONLY:  dt_averaging_input_pr
293
294!       USE netcdf_interface,                                                  &
295!           ONLY:  netcdf_data_format_string, netcdf_deflate
296
297       IMPLICIT NONE
298
299!       CHARACTER (LEN=40) ::  output_format       !< internal string
300
301       INTEGER(iwp) ::  i                         !< internal counter
302       INTEGER(iwp), INTENT(IN) ::  io            !< Unit of the output file
303
304!
305!--    Spectra output
306       IF ( dt_dosp /= 9999999.9_wp )  THEN
307          WRITE ( io, 1 )
308
309!          output_format = netcdf_data_format_string
310!          IF ( netcdf_deflate == 0 )  THEN
311!             WRITE ( io, 2 )  output_format
312!          ELSE
313!             WRITE ( io, 3 )  TRIM( output_format ), netcdf_deflate
314!          ENDIF
315          WRITE ( io, 2 )  'see profiles or other quantities'
316          WRITE ( io, 4 )  dt_dosp
317          IF ( skip_time_dosp /= 0.0_wp )  WRITE ( io, 5 )  skip_time_dosp
318          WRITE ( io, 6 )  ( data_output_sp(i), i = 1,10 ),     &
319                           ( spectra_direction(i), i = 1,10 ),  &
320                           ( comp_spectra_level(i), i = 1,100 ), &
321                           averaging_interval_sp, dt_averaging_input_pr
322       ENDIF
323
324     1 FORMAT ('    Spectra:')
325     2 FORMAT ('       Output format: ',A/)
326!     3 FORMAT ('       Output format: ',A, '   compressed with level: ',I1/)
327     4 FORMAT ('       Output every ',F7.1,' s'/)
328     5 FORMAT ('       No output during initial ',F8.2,' s')
329     6 FORMAT ('       Arrays:     ', 10(A5,',')/                         &
330               '       Directions: ', 10(A5,',')/                         &
331               '       height levels  k = ', 20(I3,',')/                  &
332               '                          ', 20(I3,',')/                  &
333               '                          ', 20(I3,',')/                  &
334               '                          ', 20(I3,',')/                  &
335               '                          ', 19(I3,','),I3,'.'/           &
336               '       Time averaged over ', F7.1, ' s,' /                &
337               '       Profiles for the time averaging are taken every ', &
338                    F6.1,' s')
339
340    END SUBROUTINE spectra_header
341
342
343
344    SUBROUTINE calc_spectra
345
346       USE arrays_3d,                                                          &
347           ONLY:  d, tend
348
349       USE control_parameters,                                                 &
350           ONLY:  bc_lr_cyc, bc_ns_cyc, message_string, psolver
351
352       USE cpulog,                                                             &
353           ONLY:  cpu_log, log_point
354
355       USE fft_xy,                                                             &
356           ONLY:  fft_init
357
358       USE indices,                                                            &
359           ONLY:  nxl, nxr, nyn, nys, nzb, nzt
360
361       USE kinds
362
363       USE pegrid,                                                             &
364           ONLY:  myid, pdims
365
366       IMPLICIT NONE
367
368       INTEGER(iwp) ::  m  !<
369       INTEGER(iwp) ::  pr !<
370
371
372!
373!--    Check if user gave any levels for spectra to be calculated
374       IF ( comp_spectra_level(1) == 999999 )  RETURN
375
376       CALL cpu_log( log_point(30), 'calc_spectra', 'start' )
377
378!
379!--    Initialize spectra related quantities
380       CALL spectra_init
381
382!
383!--    Initialize ffts
384       CALL fft_init
385
386!
387!--    Reallocate array d in required size
388       IF ( psolver(1:9) == 'multigrid' )  THEN
389          DEALLOCATE( d )
390          ALLOCATE( d(nzb+1:nzt,nys:nyn,nxl:nxr) )
391       ENDIF
392
393       m = 1
394       DO WHILE ( data_output_sp(m) /= ' '  .AND.  m <= 10 )
395!
396!--       Transposition from z --> x  ( y --> x in case of a 1d-decomposition
397!--       along x)
398          IF ( INDEX( spectra_direction(m), 'x' ) /= 0 )  THEN
399
400!
401!--          Calculation of spectra works for cyclic boundary conditions only
402             IF ( .NOT. bc_lr_cyc )  THEN
403
404                message_string = 'non-cyclic lateral boundaries along x do'//  &
405                                 ' not &  allow calculation of spectra along x'
406                CALL message( 'calc_spectra', 'PA0160', 1, 2, 0, 6, 0 )
407             ENDIF
408
409             CALL preprocess_spectra( m, pr )
410
411#if defined( __parallel )
412             IF ( pdims(2) /= 1 )  THEN
413                CALL resort_for_zx( d, tend )
414                CALL transpose_zx( tend, d )
415             ELSE
416                CALL transpose_yxd( d, d )
417             ENDIF
418             CALL calc_spectra_x( d, m )
419#else
420             message_string = 'sorry, calculation of spectra in non paral' //  &
421                              'lel mode& is still not realized'
422             CALL message( 'calc_spectra', 'PA0161', 1, 2, 0, 6, 0 )
423#endif
424
425          ENDIF
426
427!
428!--       Transposition from z --> y (d is rearranged only in case of a
429!--       1d-decomposition along x)
430          IF ( INDEX( spectra_direction(m), 'y' ) /= 0 )  THEN
431
432!
433!--          Calculation of spectra works for cyclic boundary conditions only
434             IF ( .NOT. bc_ns_cyc )  THEN
435                IF ( myid == 0 )  THEN
436                   message_string = 'non-cyclic lateral boundaries along y' // &
437                                    ' do not & allow calculation of spectra' //&
438                                    ' along y'
439                   CALL message( 'calc_spectra', 'PA0162', 1, 2, 0, 6, 0 )
440                ENDIF
441                CALL local_stop
442             ENDIF
443
444             CALL preprocess_spectra( m, pr )
445
446#if defined( __parallel )
447             CALL transpose_zyd( d, d )
448             CALL calc_spectra_y( d, m )
449#else
450             message_string = 'sorry, calculation of spectra in non paral' //  &
451                              'lel mode& is still not realized'
452             CALL message( 'calc_spectra', 'PA0161', 1, 2, 0, 6, 0 )
453#endif
454
455          ENDIF
456
457!
458!--       Increase counter for next spectrum
459          m = m + 1
460         
461       ENDDO
462
463!
464!--    Increase counter for averaging process in routine plot_spectra
465       average_count_sp = average_count_sp + 1
466
467       CALL cpu_log( log_point(30), 'calc_spectra', 'stop' )
468
469    END SUBROUTINE calc_spectra
470
471
472!------------------------------------------------------------------------------!
473! Description:
474! ------------
475!> @todo Missing subroutine description.
476!------------------------------------------------------------------------------!
477    SUBROUTINE preprocess_spectra( m, pr )
478
479       USE arrays_3d,                                                          &
480           ONLY:  d, pt, q, s, u, v, w
481
482       USE indices,                                                            &
483           ONLY:  ngp_2dh, nxl, nxr, nyn, nys, nzb, nzt
484
485       USE kinds
486
487#if defined( __parallel )
488#if !defined( __mpifh )
489       USE MPI
490#endif
491#endif
492
493       USE pegrid,                                                             &
494           ONLY:  collective_wait, comm2d, ierr
495
496       USE statistics,                                                         &
497           ONLY:  hom
498
499
500       IMPLICIT NONE
501
502#if defined( __parallel )
503#if defined( __mpifh )
504       INCLUDE "mpif.h"
505#endif
506#endif
507
508       INTEGER(iwp) :: i  !<
509       INTEGER(iwp) :: j  !<
510       INTEGER(iwp) :: k  !<
511       INTEGER(iwp) :: m  !<
512       INTEGER(iwp) :: pr !<
513
514       REAL(wp), DIMENSION(nzb:nzt+1) :: var_d_l
515
516       SELECT CASE ( TRIM( data_output_sp(m) ) )
517         
518       CASE ( 'u' )
519          pr = 1
520          d(nzb+1:nzt,nys:nyn,nxl:nxr) = u(nzb+1:nzt,nys:nyn,nxl:nxr)
521       
522       CASE ( 'v' )
523          pr = 2
524          d(nzb+1:nzt,nys:nyn,nxl:nxr) = v(nzb+1:nzt,nys:nyn,nxl:nxr)
525       
526       CASE ( 'w' )
527          pr = 3
528          d(nzb+1:nzt,nys:nyn,nxl:nxr) = w(nzb+1:nzt,nys:nyn,nxl:nxr)
529       
530       CASE ( 'theta' )
531          pr = 4
532          d(nzb+1:nzt,nys:nyn,nxl:nxr) = pt(nzb+1:nzt,nys:nyn,nxl:nxr)
533       
534       CASE ( 'q' )
535          pr = 41
536          d(nzb+1:nzt,nys:nyn,nxl:nxr) = q(nzb+1:nzt,nys:nyn,nxl:nxr)
537         
538       CASE ( 's' )
539          pr = 117
540          d(nzb+1:nzt,nys:nyn,nxl:nxr) = s(nzb+1:nzt,nys:nyn,nxl:nxr)
541       
542       CASE DEFAULT
543!
544!--       The DEFAULT case is reached either if the parameter data_output_sp(m)
545!--       contains a wrong character string or if the user has coded a special
546!--       case in the user interface. There, the subroutine user_spectra
547!--       checks which of these two conditions applies.
548          CALL user_spectra( 'preprocess', m, pr )
549         
550       END SELECT
551
552!
553!--    Subtract horizontal mean from the array, for which spectra have to be
554!--    calculated. Moreover, calculate variance of the respective quantitiy,
555!--    later used for normalizing spectra output.
556       var_d_l(:) = 0.0_wp
557       DO  i = nxl, nxr
558          DO  j = nys, nyn
559             DO  k = nzb+1, nzt
560                d(k,j,i)   = d(k,j,i) - hom(k,1,pr,0)
561                var_d_l(k) = var_d_l(k) + d(k,j,i) * d(k,j,i)
562             ENDDO
563          ENDDO
564       ENDDO
565!
566!--    Compute total variance from local variances
567       var_d(:) = 0.0_wp
568#if defined( __parallel ) 
569       IF ( collective_wait )  CALL MPI_BARRIER( comm2d, ierr )
570       CALL MPI_ALLREDUCE( var_d_l(0), var_d(0), nzt+1-nzb, MPI_REAL, MPI_SUM, &
571                           comm2d, ierr )
572#else
573       var_d(:) = var_d_l(:)
574#endif
575       var_d(:) = var_d(:) / ngp_2dh(0)
576
577    END SUBROUTINE preprocess_spectra
578
579
580!------------------------------------------------------------------------------!
581! Description:
582! ------------
583!> @todo Missing subroutine description.
584!------------------------------------------------------------------------------!
585    SUBROUTINE calc_spectra_x( ddd, m )
586
587       USE fft_xy,                                                             &
588           ONLY:  fft_x_1d
589
590       USE grid_variables,                                                     &
591           ONLY:  dx
592
593       USE indices,                                                            &
594           ONLY:  nx, ny
595
596       USE kinds
597
598#if defined( __parallel )
599#if !defined( __mpifh )
600       USE MPI
601#endif
602#endif
603
604       USE pegrid,                                                             &
605           ONLY:  comm2d, ierr, myid
606
607       USE transpose_indices,                                                  &
608           ONLY:  nyn_x, nys_x, nzb_x, nzt_x
609
610
611       IMPLICIT NONE
612
613#if defined( __parallel )
614#if defined( __mpifh )
615       INCLUDE "mpif.h"
616#endif
617#endif
618
619       INTEGER(iwp) ::  i         !<
620       INTEGER(iwp) ::  j         !<
621       INTEGER(iwp) ::  k         !<
622       INTEGER(iwp) ::  m         !<
623       INTEGER(iwp) ::  n         !<
624
625       REAL(wp) ::  exponent     !<
626       REAL(wp) ::  sum_spec_dum !< wavenumber-integrated spectrum
627   
628       REAL(wp), DIMENSION(0:nx) ::  work !<
629   
630       REAL(wp), DIMENSION(0:nx/2) ::  sums_spectra_l !<
631   
632       REAL(wp), DIMENSION(0:nx/2,100) ::  sums_spectra !<
633   
634       REAL(wp), DIMENSION(0:nx,nys_x:nyn_x,nzb_x:nzt_x) ::  ddd !<
635
636!
637!--    Exponent for geometric average
638       exponent = 1.0_wp / ( ny + 1.0_wp )
639
640!
641!--    Loop over all levels defined by the user
642       n = 1
643       DO WHILE ( comp_spectra_level(n) /= 999999  .AND.  n <= 100 )
644
645          k = comp_spectra_level(n)
646
647!
648!--       Calculate FFT only if the corresponding level is situated on this PE
649          IF ( k >= nzb_x  .AND.  k <= nzt_x )  THEN
650         
651             DO  j = nys_x, nyn_x
652
653                work = ddd(0:nx,j,k)
654                CALL fft_x_1d( work, 'forward' )
655
656                ddd(0,j,k) = dx * work(0)**2
657                DO  i = 1, nx/2
658                   ddd(i,j,k) = dx * ( work(i)**2 + work(nx+1-i)**2 )
659                ENDDO
660
661             ENDDO
662
663!
664!--          Local sum and geometric average of these spectra
665!--          (WARNING: no global sum should be performed, because floating
666!--          point overflow may occur)
667             DO  i = 0, nx/2
668
669                sums_spectra_l(i) = 1.0_wp
670                DO  j = nys_x, nyn_x
671                   sums_spectra_l(i) = sums_spectra_l(i) * ddd(i,j,k)**exponent
672                ENDDO
673
674             ENDDO
675         
676          ELSE
677
678             sums_spectra_l = 1.0_wp
679
680          ENDIF
681
682!
683!--       Global sum of spectra on PE0 (from where they are written on file)
684          sums_spectra(:,n) = 0.0_wp
685#if defined( __parallel )   
686          CALL MPI_BARRIER( comm2d, ierr )  ! Necessary?
687          CALL MPI_REDUCE( sums_spectra_l(0), sums_spectra(0,n), nx/2+1,       &
688                           MPI_REAL, MPI_PROD, 0, comm2d, ierr )
689#else
690          sums_spectra(:,n) = sums_spectra_l
691#endif
692!
693!--       Normalize spectra by variance
694          sum_spec_dum = SUM( sums_spectra(1:nx/2,n) )
695
696          IF ( sum_spec_dum /= 0.0_wp )  THEN
697             sums_spectra(1:nx/2,n) = sums_spectra(1:nx/2,n) *                 &
698                                      var_d(k) / sum_spec_dum
699          ENDIF
700          n = n + 1
701
702       ENDDO
703       n = n - 1
704
705       IF ( myid == 0 )  THEN
706!
707!--       Sum of spectra for later averaging (see routine data_output_spectra)
708          DO  i = 1, nx/2
709             DO k = 1, n
710                spectrum_x(i,k,m) = spectrum_x(i,k,m) + sums_spectra(i,k)
711             ENDDO
712          ENDDO
713
714       ENDIF
715!
716!--    n_sp_x is needed by data_output_spectra_x
717       n_sp_x = n
718
719    END SUBROUTINE calc_spectra_x
720
721
722!------------------------------------------------------------------------------!
723! Description:
724! ------------
725!> @todo Missing subroutine description.
726!------------------------------------------------------------------------------!
727    SUBROUTINE calc_spectra_y( ddd, m )
728
729       USE fft_xy,                                                             &
730           ONLY:  fft_y_1d
731
732       USE grid_variables,                                                     &
733           ONLY:  dy
734
735       USE indices,                                                            &
736           ONLY:  nx, ny
737
738       USE kinds
739
740#if defined( __parallel )
741#if !defined( __mpifh )
742       USE MPI
743#endif
744#endif
745
746       USE pegrid,                                                             &
747           ONLY:  comm2d, ierr, myid
748
749       USE transpose_indices,                                                  &
750           ONLY:  nxl_yd, nxr_yd, nzb_yd, nzt_yd
751
752
753       IMPLICIT NONE
754
755#if defined( __parallel )
756#if defined( __mpifh )
757       INCLUDE "mpif.h"
758#endif
759#endif
760
761       INTEGER(iwp) ::  i         !<
762       INTEGER(iwp) ::  j         !<
763       INTEGER(iwp) ::  k         !<
764       INTEGER(iwp) ::  m         !<
765       INTEGER(iwp) ::  n         !<
766
767       REAL(wp) ::  exponent !<
768       REAL(wp) ::  sum_spec_dum !< wavenumber-integrated spectrum
769   
770       REAL(wp), DIMENSION(0:ny) ::  work !<
771   
772       REAL(wp), DIMENSION(0:ny/2) ::  sums_spectra_l !<
773   
774       REAL(wp), DIMENSION(0:ny/2,100) ::  sums_spectra !<
775   
776       REAL(wp), DIMENSION(0:ny,nxl_yd:nxr_yd,nzb_yd:nzt_yd) :: ddd !<
777
778
779!
780!--    Exponent for geometric average
781       exponent = 1.0_wp / ( nx + 1.0_wp )
782
783!
784!--    Loop over all levels defined by the user
785       n = 1
786       DO WHILE ( comp_spectra_level(n) /= 999999  .AND.  n <= 100 )
787
788          k = comp_spectra_level(n)
789
790!
791!--       Calculate FFT only if the corresponding level is situated on this PE
792          IF ( k >= nzb_yd  .AND.  k <= nzt_yd )  THEN
793         
794             DO  i = nxl_yd, nxr_yd
795
796                work = ddd(0:ny,i,k)
797                CALL fft_y_1d( work, 'forward' )
798
799                ddd(0,i,k) = dy * work(0)**2
800                DO  j = 1, ny/2
801                   ddd(j,i,k) = dy * ( work(j)**2 + work(ny+1-j)**2 )
802                ENDDO
803
804             ENDDO
805
806!
807!--          Local sum and geometric average of these spectra
808!--          (WARNING: no global sum should be performed, because floating
809!--          point overflow may occur)
810             DO  j = 0, ny/2
811
812                sums_spectra_l(j) = 1.0_wp
813                DO  i = nxl_yd, nxr_yd
814                   sums_spectra_l(j) = sums_spectra_l(j) * ddd(j,i,k)**exponent
815                ENDDO
816
817             ENDDO
818         
819          ELSE
820
821             sums_spectra_l = 1.0_wp
822
823          ENDIF
824
825!
826!--       Global sum of spectra on PE0 (from where they are written on file)
827          sums_spectra(:,n) = 0.0_wp
828#if defined( __parallel )   
829          CALL MPI_BARRIER( comm2d, ierr )  ! Necessary?
830          CALL MPI_REDUCE( sums_spectra_l(0), sums_spectra(0,n), ny/2+1,       &
831                           MPI_REAL, MPI_PROD, 0, comm2d, ierr )
832#else
833          sums_spectra(:,n) = sums_spectra_l
834#endif
835!
836!--       Normalize spectra by variance
837          sum_spec_dum = SUM( sums_spectra(1:ny/2,n) )
838          IF ( sum_spec_dum /= 0.0_wp )  THEN
839             sums_spectra(1:ny/2,n) = sums_spectra(1:ny/2,n) *                 &
840                                      var_d(k) / sum_spec_dum
841          ENDIF
842          n = n + 1
843
844       ENDDO
845       n = n - 1
846
847
848       IF ( myid == 0 )  THEN
849!
850!--       Sum of spectra for later averaging (see routine data_output_spectra)
851          DO  j = 1, ny/2
852             DO k = 1, n
853                spectrum_y(j,k,m) = spectrum_y(j,k,m) + sums_spectra(j,k)
854             ENDDO
855          ENDDO
856
857       ENDIF
858
859!
860!--    n_sp_y is needed by data_output_spectra_y
861       n_sp_y = n
862
863    END SUBROUTINE calc_spectra_y
864
865 END MODULE spectra_mod
Note: See TracBrowser for help on using the repository browser.