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

Last change on this file since 2836 was 2718, checked in by maronga, 7 years ago

deleting of deprecated files; headers updated where needed

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