updated top-level README and version_decl for V4.5 (#1847)
[WRF.git] / var / da / da_setup_structures / da_truncate_spectra.inc
blobf6e9a399ff41ec8a8fae5e0b4db7703b1cba25bd
1 subroutine da_truncate_spectra(max_wave, nw, power_trunc, power, max_wave_trunc)
3    !-----------------------------------------------------------------------
4    ! Purpose: TBD
5    !-----------------------------------------------------------------------
7    implicit none
9    integer, intent(in)  :: max_wave         ! Smallest wave for domain.
10    integer, intent(in)  :: nw               ! Dimension of input power spectrum.
11    real,    intent(in)  :: power_trunc      ! Power truncation (fraction).
12    real*8,  intent(in)  :: power(0:nw)      ! Power spectrum.
13    integer, intent(out) :: max_wave_trunc   ! Smallest wave after truncation.
15    integer :: l                ! Loop counter.
16    real    :: truncated_power  ! Truncated power.
17    real    :: cumul_power      ! Cumulative power.
19    if (trace_use) call da_trace_entry("da_truncate_spectra")
21    truncated_power = power_trunc * sum(power(0:nw))
23    cumul_power = 0.0
24    max_wave_trunc = max_wave
25    do l = 0, nw - 1 
26       cumul_power = cumul_power + power(l)
27       if (cumul_power > truncated_power) then
28          max_wave_trunc = l - 1
29          exit
30       end if
31    end do
33    if (max_wave_trunc > max_wave) then
34       write(unit=message(1),fmt='(a)') &
35          'da_truncate_spectra: Power requested needs higher resolution.'     
36       write(unit=message(2),fmt='(a,i8)') &
37          'Maximum grid wavenumber =  ', max_wave
38       write(unit=message(3),fmt='(a,i8)') &
39          'Truncating to wavenumber = ', max_wave_trunc
40       call da_warning(__FILE__,__LINE__,message(1:3))
41       max_wave_trunc = max_wave
42    end if
44    if (trace_use) call da_trace_exit("da_truncate_spectra")
46 end subroutine da_truncate_spectra