localization_mod sourceΒΆ

  1module localization_mod
  2  ! MODULE localization_mod (prefix='loc' category='2. B and R matrices')
  3  !
  4  !:Purpose:  Master module for the computation of localized 3D gridpoint
  5  !           amplitude fields for each ensemble member from a given (1D)
  6  !           control vector
  7  !
  8  use midasMpi_mod
  9  use utilities_mod
 10  use localizationSpectral_mod
 11  use horizontalCoord_mod
 12  use verticalCoord_mod
 13  use ensembleStatevector_mod
 14  implicit none
 15  save
 16  private
 17
 18  ! public structure
 19  public :: struct_loc
 20  ! public procedures
 21  public :: loc_setup, loc_Lsqrt, loc_LsqrtAd, loc_finalize
 22  public :: loc_reducetompilocal, loc_reducetompilocal_r4
 23  public :: loc_expandtompiglobal, loc_expandtompiglobal_r4
 24
 25  type :: struct_loc
 26     logical             :: initialized = .false.
 27     type(struct_lsp), pointer :: lsp => null()
 28     integer             :: nEns
 29     integer             :: nEnsOverDimension
 30     integer             :: cvDim
 31     character(len=64)   :: locType
 32     type(struct_vco), pointer :: vco => null()
 33     type(struct_hco), pointer :: hco => null()
 34  end type struct_loc
 35
 36  logical, parameter :: verbose = .false.
 37
 38CONTAINS
 39
 40  !--------------------------------------------------------------------------
 41  ! loc_setup
 42  !--------------------------------------------------------------------------
 43  subroutine loc_setup(loc, cvDim_out, hco_loc, vco_loc, nEns, vertLocation, ntrunc, locType, &
 44                       locMode, horizLengthScale1, horizLengthScale2, vertLengthScale)
 45    implicit none
 46
 47    ! Arguments:
 48    type(struct_loc),          intent(inout) :: loc
 49    integer,                   intent(out)   :: cvDim_out
 50    type(struct_hco), pointer, intent(in)    :: hco_loc
 51    type(struct_vco), pointer, intent(in)    :: vco_loc
 52    integer,                   intent(in)    :: nEns
 53    integer,                   intent(in)    :: nTrunc
 54    real(8),                   intent(in)    :: vertLocation(:)
 55    real(8),                   intent(in)    :: horizLengthScale1
 56    real(8),                   intent(in)    :: horizLengthScale2
 57    real(8),                   intent(in)    :: vertLengthScale
 58    character(len=*),          intent(in)    :: locType
 59    character(len=*),          intent(in)    :: locMode
 60
 61    ! Locals:
 62    integer :: nEnsOverDimension, nLev
 63
 64    if (verbose) write(*,*) 'Entering loc_Setup'
 65
 66    !
 67    !- 2.  Setup
 68    !
 69    loc%locType = trim(locType)
 70    loc%hco => hco_loc
 71    loc%vco => vco_loc
 72
 73    if ( loc%vco%Vcode == 5002 .or. loc%vco%Vcode == 5005 ) then
 74      if (loc%vco%nLev_M > 0) then
 75        nLev = loc%vco%nLev_M
 76      else
 77        nLev = loc%vco%nLev_T
 78      end if
 79    else if ( loc%vco%nLev_depth > 0 ) then
 80      nLev = loc%vco%nLev_depth
 81    else
 82      nLev = 1
 83    end if
 84
 85    select case (trim(loc%locType))
 86    case('spectral')
 87       if (mmpi_myid == 0) write(*,*)
 88       if (mmpi_myid == 0) write(*,*) 'loc_setup: LocType = ', trim(locType)
 89       call lsp_setup(hco_loc, nEns, nLev, vertLocation, ntrunc, locType,          & ! IN
 90                      locMode, horizLengthScale1, horizLengthScale2, vertLengthScale, & ! IN
 91                      cvDim_out, loc%lsp, nEnsOverDimension)                            ! OUT
 92    case default
 93       write(*,*)
 94       write(*,*) 'locType = ', trim(locType)
 95       call utl_abort('loc_setup: unknown locType')
 96    end select
 97
 98    loc%cvDim = cvDim_out
 99    loc%nEnsOverDimension = nEnsOverDimension 
100
101    !
102    !- 3.  Ending
103    !
104    loc%initialized = .true.
105
106  end subroutine loc_setup
107
108  !--------------------------------------------------------------------------
109  ! loc_Lsqrt
110  !--------------------------------------------------------------------------
111  subroutine loc_Lsqrt(loc, controlVector, ensAmplitude, stepIndex)
112    implicit none
113
114    ! Arguments:
115    type(struct_loc), intent(in)    :: loc
116    integer,          intent(in)    :: stepIndex
117    real(8),          intent(in)    :: controlVector(:)
118    type(struct_ens), intent(inout) :: ensAmplitude
119
120    if (verbose) write(*,*) 'Entering loc_Lsqrt'
121
122    select case (trim(loc%locType))
123    case('spectral')
124      call lsp_Lsqrt(loc%lsp, controlVector, & ! IN
125                     ensAmplitude,           & ! OUT
126                     stepIndex)                ! IN
127    case default
128       call utl_abort('loc_Lsqrt: unknown locType')
129    end select
130
131  end subroutine loc_Lsqrt
132
133  !--------------------------------------------------------------------------
134  ! loc_LsqrtAd
135  !--------------------------------------------------------------------------
136  subroutine loc_LsqrtAd(loc, ensAmplitude, controlVector, stepIndex)
137    implicit none
138
139    ! Arguments:
140    type(struct_loc), intent(in)    :: loc
141    integer,          intent(in)    :: stepIndex
142    real(8),          intent(out)   :: controlVector(:)
143    type(struct_ens), intent(inout) :: ensAmplitude
144
145    if (verbose) write(*,*) 'Entering loc_LsqrtAd'
146
147    select case (trim(loc%locType))
148    case('spectral')
149       call lsp_LsqrtAd(loc%lsp,       & ! IN
150                        ensAmplitude,  & ! INOUT
151                        controlVector, & ! OUT
152                        stepIndex )      ! IN
153    case default
154       call utl_abort('loc_LsqrtAd: unknown locType')
155    end select
156
157  end subroutine loc_LsqrtAd
158
159  !--------------------------------------------------------------------------
160  ! loc_finalize
161  !--------------------------------------------------------------------------
162  subroutine loc_finalize(loc)
163    implicit none
164
165    ! Arguments:
166    type(struct_loc), intent(inout) :: loc
167
168    if (verbose) write(*,*) 'Entering loc_finalize'
169
170    select case (trim(loc%locType))
171    case('spectral')
172      call lsp_finalize(loc%lsp)
173    case default
174      call utl_abort('loc_finalize: unknown locType')
175    end select
176
177  end subroutine loc_finalize
178
179  !--------------------------------------------------------------------------
180  ! loc_reduceToMPILocal
181  !--------------------------------------------------------------------------
182  subroutine loc_reduceToMPILocal(loc,cv_mpilocal,cv_mpiglobal)
183    implicit none
184
185    ! Arguments:
186    type(struct_loc), intent(in)  :: loc
187    real(8),          intent(out) :: cv_mpilocal(:)
188    real(8),          intent(in)  :: cv_mpiglobal(:)
189
190    if (verbose) write(*,*) 'Entering loc_reduceToMPILocal'
191
192    select case (trim(loc%locType))
193    case('spectral')
194      call lsp_reduceToMPILocal(loc%lsp,      & ! IN
195                                cv_mpilocal, & ! OUT
196                                cv_mpiglobal)  ! IN
197    case default
198      call utl_abort('loc_reduceToMPILocal: unknown locType')
199    end select
200
201  end subroutine loc_reduceToMPILocal
202 
203  !--------------------------------------------------------------------------
204  ! loc_reduceToMPILocal_r4
205  !--------------------------------------------------------------------------
206  subroutine loc_reduceToMPILocal_r4(loc,cv_mpilocal,cv_mpiglobal)
207    implicit none
208
209    ! Arguments:
210    type(struct_loc), intent(in)  :: loc
211    real(4),          intent(out) :: cv_mpilocal(:)
212    real(4),          intent(in)  :: cv_mpiglobal(:)
213
214    if (verbose) write(*,*) 'Entering loc_reduceToMPILocal_r4'
215
216    select case (trim(loc%locType))
217    case('spectral')
218      call lsp_reduceToMPILocal_r4(loc%lsp,       & ! IN
219                                   cv_mpilocal,  & ! OUT
220                                   cv_mpiglobal)   ! IN
221    case default
222      call utl_abort('loc_reduceToMPILocal_r4: unknown locType')
223    end select
224
225  end subroutine loc_reduceToMPILocal_r4
226 
227  !--------------------------------------------------------------------------
228  ! loc_expandToMPIGlobal
229  !--------------------------------------------------------------------------
230  subroutine loc_expandToMPIGlobal(loc,cv_mpilocal,cv_mpiglobal)
231    implicit none
232
233    ! Arguments:
234    type(struct_loc), intent(in)  :: loc
235    real(8),          intent(in)  :: cv_mpilocal(:)
236    real(8),          intent(out) :: cv_mpiglobal(:)
237
238    if (verbose) write(*,*) 'Entering loc_expandToMPIGlobal'
239    
240    select case (trim(loc%locType))
241    case('spectral')
242      call lsp_expandToMPIGlobal(loc%lsp, cv_mpilocal,  & ! IN
243                                 cv_mpiglobal)           ! OUT
244    case default
245      call utl_abort('loc_expandToMPIGlobal: unknown locType')
246    end select
247
248  end subroutine loc_expandToMPIGlobal
249
250  !--------------------------------------------------------------------------
251  ! loc_expandToMPIGlobal_r4
252  !--------------------------------------------------------------------------
253  subroutine loc_expandToMPIGlobal_r4(loc,cv_mpilocal,cv_mpiglobal)
254    implicit none
255
256    ! Arguments:
257    type(struct_loc), intent(in)  :: loc
258    real(4),          intent(in)  :: cv_mpilocal(:)
259    real(4),          intent(out) :: cv_mpiglobal(:)
260
261    if (verbose) write(*,*) 'Entering loc_expandToMPIGlobal_r4'
262    
263    select case (trim(loc%locType))
264    case('spectral')
265      call lsp_expandToMPIGlobal_r4(loc%lsp, cv_mpilocal, & ! IN
266                                    cv_mpiglobal)          ! OUT
267    case default
268      call utl_abort('loc_expandToMPIGlobal_r4: unknown locType')
269    end select
270
271  end subroutine  loc_expandToMPIGlobal_r4
272
273end module localization_mod