New URL for NEMO forge!   http://forge.nemo-ocean.eu

Since March 2022 along with NEMO 4.2 release, the code development moved to a self-hosted GitLab.
This present forge is now archived and remained online for history.
interp_nearest.f90 in utils/tools/SIREN/src – NEMO

source: utils/tools/SIREN/src/interp_nearest.f90 @ 13369

Last change on this file since 13369 was 13369, checked in by jpaul, 4 years ago

update: cf changelog inside documentation

File size: 11.7 KB
Line 
1!----------------------------------------------------------------------
2! NEMO system team, System and Interface for oceanic RElocable Nesting
3!----------------------------------------------------------------------
4!
5! DESCRIPTION:
6!> @brief
7!> This module manage nearest interpolation on regular grid.
8!>
9!> @details
10!> to compute nearest interpolation:<br/>
11!> @code
12!> CALL interp_nearest_fill(dd_value, dd_fill, id_detect, id_rho, ld_even [,ld_discont] )
13!> @endcode
14!>    - dd_value is 2D array of variable value
15!>    - dd_fill is the FillValue of variable
16!>    - id_detect is 2D array of point to be interpolated (see interp module)
17!>    - id_rho  is array of refinment factor
18!>    - ld_even indicates even refinment or not
19!>    - ld_discont indicates longitudinal discontinuity (-180°/180°, 0°/360°) or not
20!>
21!> @author
22!> J.Paul
23!>
24!> @date September, 2014 - Initial version
25!>
26!> @note Software governed by the CeCILL licence     (NEMOGCM/NEMO_CeCILL.txt)
27!----------------------------------------------------------------------
28MODULE interp_nearest
29
30   USE netcdf                          ! nf90 library
31   USE global                          ! global variable
32   USE kind                            ! F90 kind parameter
33   USE logger                          ! log file manager
34   USE fct                             ! basic useful function
35
36   IMPLICIT NONE
37   ! NOTE_avoid_public_variables_if_possible
38
39   ! type and variable
40
41   ! function and subroutine
42   PUBLIC :: interp_nearest_fill  !< compute interpolation using nearest method
43
44   PRIVATE :: interp_nearest__2D       !< compute binearest interpolation on 2D gid
45   PRIVATE :: interp_nearest__1D       !< compute   nearest interpolation on 1D gid
46   PRIVATE :: interp_nearest__2D_fill  !< fill value using binearest interpolation
47   PRIVATE :: interp_nearest__1D_fill  !< fill value using   nearest interpolation
48
49CONTAINS
50   !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
51   SUBROUTINE interp_nearest_fill(dd_value, id_detect, id_rho)
52   !-------------------------------------------------------------------
53   !> @brief
54   !> This subroutine compute horizontal nearest interpolation on 4D array of value.
55   !>
56   !> @author J.Paul
57   !> @date September, 2014 - Initial Version
58   !>
59   !> @param[inout] dd_value  2D array of variable value
60   !> @param[inout] id_detect 2D array of point to be interpolated
61   !> @param[in]    id_rho    array of refinment factor
62   !-------------------------------------------------------------------
63
64      IMPLICIT NONE
65
66      ! Argument
67      REAL(dp)        , DIMENSION(:,:,:,:), INTENT(INOUT) :: dd_value
68      INTEGER(I4)     , DIMENSION(:,:,:)  , INTENT(INOUT) :: id_detect
69      INTEGER(I4)     , DIMENSION(:)      , INTENT(IN   ) :: id_rho
70
71      ! local variable
72      INTEGER(i4), DIMENSION(4)                  :: il_shape
73
74      INTEGER(I4), DIMENSION(:,:,:), ALLOCATABLE :: il_detect
75
76      ! loop indices
77      INTEGER(i4) :: ji
78      INTEGER(i4) :: jj
79      INTEGER(i4) :: jk
80      INTEGER(i4) :: jl
81      !----------------------------------------------------------------
82
83      il_shape(:)=SHAPE(dd_value)
84
85      ALLOCATE(il_detect(il_shape(1),il_shape(2),il_shape(3)))
86      DO jl=1,il_shape(4)
87         il_detect(:,:,:)=id_detect(:,:,:)
88         ! loop on vertical level
89         DO jk=1,il_shape(3)
90
91            ! I-J plan
92            CALL interp_nearest__2D(dd_value(:,:,jk,jl),&
93            &                       il_detect(:,:,jk),  &
94            &                       id_rho(jp_I), id_rho(jp_J) )
95            IF( ANY(il_detect(:,:,jk)==1) )THEN
96               ! I direction
97               DO jj=1,il_shape(2)
98                  CALL interp_nearest__1D( dd_value(:,jj,jk,jl),&
99                  &                        il_detect(:,jj,jk),  &
100                  &                        id_rho(jp_I) )
101               ENDDO
102               IF( ALL(il_detect(:,:,jk)==0) )THEN
103                  CYCLE
104               ELSE
105                  ! J direction
106                  DO ji=1,il_shape(1)
107                     CALL interp_nearest__1D( dd_value(ji,:,jk,jl),&
108                     &                        il_detect(ji,:,jk),  &
109                     &                        id_rho(jp_J) )
110                  ENDDO
111               ENDIF
112            ENDIF
113
114         ENDDO
115      ENDDO
116
117      id_detect(:,:,:)=il_detect(:,:,:)
118      DEALLOCATE(il_detect)
119
120   END SUBROUTINE interp_nearest_fill
121   !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
122   SUBROUTINE interp_nearest__2D(dd_value, id_detect, id_rhoi, id_rhoj)
123   !-------------------------------------------------------------------
124   !> @brief
125   !> This subroutine compute nearest interpolation on 2D array of value.
126   !>
127   !> @author J.Paul
128   !> @date September, 2014 - Initial Version
129   !>
130   !> @param[inout] dd_value  2D array of variable value
131   !> @param[inout] id_detect 2D array of point to be interpolated
132   !> @param[in] id_rhoi      refinment factor in i-direction
133   !> @param[in] id_rhoj      refinment factor in j-direction
134   !> @param[in] id_rhok      refinment factor in k-direction
135   !-------------------------------------------------------------------
136
137      IMPLICIT NONE
138
139      ! Argument
140      REAL(dp)        , DIMENSION(:,:), INTENT(INOUT) :: dd_value
141      INTEGER(I4)     , DIMENSION(:,:), INTENT(INOUT) :: id_detect
142      INTEGER(I4)                     , INTENT(IN   ) :: id_rhoi
143      INTEGER(I4)                     , INTENT(IN   ) :: id_rhoj
144
145      ! local variable
146      INTEGER(i4), DIMENSION(2)                :: il_shape
147
148      ! loop indices
149      INTEGER(i4) :: ji
150      INTEGER(i4) :: jj
151
152      !----------------------------------------------------------------
153
154      IF( ANY(id_detect(:,:)==1) )THEN
155
156         il_shape(:)=SHAPE(dd_value)
157
158         DO jj=1,il_shape(2)-1,id_rhoj
159            DO ji=1,il_shape(1)-1,id_rhoi
160
161               ! check if point to be interpolated
162               IF( ALL(id_detect(ji:ji+id_rhoi,   &
163               &                 jj:jj+id_rhoj)==0) ) CYCLE
164
165               ! compute value on detetected point
166               CALL interp_nearest__2D_fill(dd_value( ji:ji+id_rhoi,   &
167               &                                      jj:jj+id_rhoj ), &
168               &                            id_detect(ji:ji+id_rhoi,   &
169               &                                      jj:jj+id_rhoj ) )
170
171            ENDDO
172         ENDDO
173
174      ENDIF
175
176   END SUBROUTINE interp_nearest__2D
177   !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
178   SUBROUTINE interp_nearest__1D(dd_value, id_detect, id_rhoi)
179   !-------------------------------------------------------------------
180   !> @brief
181   !> This subroutine compute nearest interpolation on 1D array of value.
182   !>
183   !> @author J.Paul
184   !> @date September, 2014 - Initial Version
185   !>
186   !> @param[inout] dd_value  1D array of variable value
187   !> @param[inout] id_detect 1D array of point to be interpolated
188   !> @param[in]    id_rhoi   refinment factor
189   !-------------------------------------------------------------------
190
191      IMPLICIT NONE
192
193      ! Argument
194      REAL(dp)        , DIMENSION(:), INTENT(INOUT) :: dd_value
195      INTEGER(I4)     , DIMENSION(:), INTENT(INOUT) :: id_detect
196      INTEGER(I4)                   , INTENT(IN   ) :: id_rhoi
197
198      ! local variable
199      INTEGER(i4), DIMENSION(1)              :: il_shape
200
201      ! loop indices
202      INTEGER(i4) :: ji
203      !----------------------------------------------------------------
204
205      IF( ANY(id_detect(:)==1) )THEN
206         il_shape(:)=SHAPE(dd_value)
207
208         DO ji=1,il_shape(1)-1,id_rhoi
209
210            ! check if point to be interpolated
211            IF( ALL(id_detect(ji:ji+id_rhoi)==0) ) CYCLE
212
213            ! compute value on detetected point
214            CALL interp_nearest__1D_fill( dd_value( ji:ji+id_rhoi ), &
215            &                             id_detect(ji:ji+id_rhoi ) )
216
217         ENDDO
218
219      ENDIF
220
221   END SUBROUTINE interp_nearest__1D
222   !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
223   SUBROUTINE interp_nearest__2D_fill(dd_value, id_detect)
224   !-------------------------------------------------------------------
225   !> @brief
226   !> This subroutine compute nearest interpolation of a 2D array of value.
227   !>
228   !> @author J.Paul
229   !> @date September, 2014 - Initial Version
230   !>
231   !> @param[inout] dd_value  2D array of mixed grid value
232   !> @param[inout] id_detect 2D array of point to be interpolated
233   !-------------------------------------------------------------------
234
235      IMPLICIT NONE
236
237      ! Argument
238      REAL(dp)   , DIMENSION(:,:)  , INTENT(INOUT) :: dd_value
239      INTEGER(i4), DIMENSION(:,:)  , INTENT(INOUT) :: id_detect
240
241      ! local variable
242      INTEGER(i4), DIMENSION(2) :: il_shape
243
244      INTEGER(i4) :: il_i1
245      INTEGER(i4) :: il_i2
246      INTEGER(i4) :: il_j1
247      INTEGER(i4) :: il_j2
248
249      INTEGER(i4) :: il_half1
250      INTEGER(i4) :: il_half2
251
252      ! loop indices
253      INTEGER(i4) :: ji
254      INTEGER(i4) :: jj
255      !----------------------------------------------------------------
256
257      il_shape(:)=SHAPE(dd_value(:,:))
258
259      il_i1=1
260      il_i2=il_shape(1)
261
262      il_j1=1
263      il_j2=il_shape(2)
264
265      il_half1=CEILING(il_shape(1)*0.5)
266      il_half2=CEILING(il_shape(2)*0.5)
267
268      DO jj=1,il_half2
269
270         DO ji=1,il_half1
271
272            ! lower left point
273            IF(id_detect(ji,jj)==1)THEN
274
275               dd_value( ji,jj)=dd_value(il_i1,il_j1)
276               id_detect(ji,jj)=0
277
278            ENDIF
279
280            ! lower right point
281            IF(id_detect(il_shape(1)-ji+1,jj)==1)THEN
282
283               dd_value( il_shape(1)-ji+1,jj)=dd_value(il_i2,il_j1)
284               id_detect(il_shape(1)-ji+1,jj)=0
285
286            ENDIF
287
288            ! upper left point
289            IF(id_detect(ji,il_shape(2)-jj+1)==1)THEN
290
291               dd_value( ji,il_shape(2)-jj+1)=dd_value(il_i1,il_j2)
292               id_detect(ji,il_shape(2)-jj+1)=0
293
294            ENDIF
295
296            ! upper right point
297            IF(id_detect(il_shape(1)-ji+1,il_shape(2)-jj+1)==1)THEN
298
299               dd_value( il_shape(1)-ji+1,il_shape(2)-jj+1)=dd_value(il_i2,il_j2)
300               id_detect(il_shape(1)-ji+1,il_shape(2)-jj+1)=0
301
302            ENDIF
303
304         ENDDO
305
306      ENDDO
307
308   END SUBROUTINE interp_nearest__2D_fill
309   !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
310   SUBROUTINE interp_nearest__1D_fill(dd_value, id_detect)
311   !-------------------------------------------------------------------
312   !> @brief
313   !> This subroutine compute nearest interpolation of a 1D array of value.
314   !>
315   !> @author J.Paul
316   !> @date September, 2014 - Initial Version
317   !>
318   !> @param[inout] dd_value  1D array of mixed grid value
319   !> @param[inout] id_detect 1D array of point to be interpolated
320   !-------------------------------------------------------------------
321
322      IMPLICIT NONE
323
324      ! Argument
325      REAL(dp)   , DIMENSION(:), INTENT(INOUT) :: dd_value
326      INTEGER(i4), DIMENSION(:), INTENT(INOUT) :: id_detect
327
328      ! local variable
329      INTEGER(i4), DIMENSION(1) :: il_shape
330
331      INTEGER(i4) :: il_i1
332      INTEGER(i4) :: il_i2
333
334      INTEGER(i4) :: il_half1
335
336      ! loop indices
337      INTEGER(i4) :: ji
338      !----------------------------------------------------------------
339
340      il_shape(:)=SHAPE(dd_value)
341
342      il_i1=1
343      il_i2=il_shape(1)
344
345      il_half1=CEILING(il_shape(1)*0.5)
346
347      DO ji=1,il_half1
348
349         ! lower left point
350         IF(id_detect(ji)==1)THEN
351
352            dd_value( ji)=dd_value(il_i1)
353            id_detect(ji)=0
354
355         ENDIF
356
357         ! lower right point
358         IF(id_detect(il_shape(1)-ji+1)==1)THEN
359
360            dd_value( il_shape(1)-ji+1)=dd_value(il_i2)
361            id_detect(il_shape(1)-ji+1)=0
362
363         ENDIF
364
365      ENDDO
366
367   END SUBROUTINE interp_nearest__1D_fill
368   !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
369END MODULE interp_nearest
Note: See TracBrowser for help on using the repository browser.