updated top-level README and version_decl for V4.5 (#1847)
[WRF.git] / var / da / da_spectral / da_initialize_h.inc
blob899632cef3e2e8ed05d1c419e0264d1c5fd48803
1 subroutine da_initialize_h(ni, nj, max_wavenumber, lensav, alp_size, wsave, lon, sinlon, coslon, &
2    lat, sinlat, coslat, int_wgts, alp)
4    !-----------------------------------------------------------------------
5    ! Purpose: TBD
6    !-----------------------------------------------------------------------
8    implicit none
10    integer, intent(in)  :: ni                         ! Number of longitudes.
11    integer, intent(in)  :: nj                         ! Number of latitudes.
12    integer, intent(in)  :: max_wavenumber             ! Smallest scale required (ni/2 - 1).
13    integer, intent(in)  :: lensav                     ! Size of FFTs wsave array.
14    integer, intent(in)  :: alp_size                   ! Size of ALP array.
15    real,    intent(out) :: wsave(1:lensav)            ! Primes for FFT.
16    real,    intent(out) :: lon(1:ni)                  ! Longitude (radians).
17    real,    intent(out) :: sinlon(1:ni)               ! sine(longitude).
18    real,    intent(out) :: coslon(1:ni)               ! cosine(longitude).
19    real,    intent(out) :: lat(1:nj)                  ! Latitude (radians, from south).
20    real,    intent(out) :: sinlat(1:nj)               ! sine(latitude).
21    real,    intent(out) :: coslat(1:nj)               ! cosine(latitude).
22    real,    intent(out) :: int_wgts(1:nj)             ! Legendre integration weights.
23    real,    intent(out) :: alp(1:alp_size)            ! Associated Legendre Polynomial.
25    integer :: i                          ! Loop counters.
27    if (trace_use) call da_trace_entry("da_initialize_h")
29    !----------------------------------------------------------------------------
30    ! [1] Initialize FFT coefficients.'
31    !----------------------------------------------------------------------------
33    wsave(:) = 0.0
34 #ifdef FFTPACK
35    call rfft1i(ni, wsave, lensav, ierr)
36 #else
37    call da_error(__FILE__,__LINE__,(/"Needs to be compiled with FFTPACK"/))
38 #endif
40    if (ierr /= 0) then
41      write(unit=message(1),fmt='(A,I4)') &
42         "Fourier initialization failed. ierr = ", ierr
43      call da_error(__FILE__,__LINE__,message(1:1))
44    end if
46    !----------------------------------------------------------------------------
47    ! [2] Calculate latitudes, and their sines/cosines.'
48    !---------------------------------------------------------------------------
50    if (gaussian_lats) then
51       call da_get_gausslats(nj, lat, int_wgts, sinlat, coslat)
52    else
53       call da_get_reglats(nj, lat, sinlat, coslat, int_wgts)
54    end if
56    do i = 1, ni
57       lon(i) = 2.0 * pi / real(ni) * real(i - 1)
58       sinlon(i) = sin(lon(i))
59       coslon(i) = cos(lon(i))
60    end do
62    !----------------------------------------------------------------------------
63    ! [3] Initialize Legendre coefficients.'
64    !----------------------------------------------------------------------------
66    call da_setlegpol(nj, max_wavenumber, alp_size, sinlat, coslat, alp)
68    if (trace_use) call da_trace_exit("da_initialize_h")
70 end subroutine da_initialize_h