source: trunk/LMDZ.COMMON/libf/evolution/geometry.F90 @ 4152

Last change on this file since 4152 was 4152, checked in by jbclement, 2 weeks ago

PEM:

  • Fix outputs "diagevo.nc" for 3D data.
  • Fix sign in computing exchanges due to adsorption/ice table and in balancing H2O flux from/into atmosphere.
  • Few cleanings.

JBC

File size: 12.4 KB
Line 
1MODULE geometry
2!-----------------------------------------------------------------------
3! NAME
4!     geometry
5!
6! DESCRIPTION
7!     Grid metrics/parameters and basic tools to convert between
8!     different grids.
9!
10! AUTHORS & DATE
11!     JB Clement, 12/2025
12!
13! NOTES
14!
15!-----------------------------------------------------------------------
16
17! DEPENDENCIES
18! ------------
19use numerics, only: dp, di, k4
20
21! DECLARATION
22! -----------
23implicit none
24
25! PARAMETERS
26! ----------
27integer(di),                         protected :: nlon               ! Number of longitudes
28integer(di),                         protected :: nlat               ! Number of latitudes
29integer(di),                         protected :: ngrid              ! Number of grid points
30integer(di),                         protected :: nlayer             ! Number of atmospheric layers
31integer(di),                         protected :: nsoil_PCM          ! Number of soil layers in the PCM
32integer(di),                         parameter :: nsoil = 68_di      ! Number of soil layers in the PEM
33integer(di),                         protected :: nslope             ! Number of slopes
34integer(di),                         protected :: nday               ! Number of sols in a year
35real(dp), dimension(:), allocatable, protected :: longitudes         ! Longitudes
36real(dp), dimension(:), allocatable, protected :: latitudes          ! Latitudes
37real(dp), dimension(:), allocatable, protected :: cell_area          ! Cell area
38real(dp),                            protected :: total_surface      ! Total surface of the grid/planet
39logical(k4),                         protected :: dim_init = .false. ! Flag to true if dimensions are well initialized
40
41contains
42!+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
43
44!=======================================================================
45SUBROUTINE ini_geometry()
46!-----------------------------------------------------------------------
47! NAME
48!     ini_geometry
49!
50! DESCRIPTION
51!     Initialize the parameters of module 'geometry'.
52!
53! AUTHORS & DATE
54!     JB Clement, 12/2025
55!
56! NOTES
57!
58!-----------------------------------------------------------------------
59
60! DEPENDENCIES
61! ------------
62use stoppage,   only: stop_clean
63use io_netcdf,  only: startfi_name, xios_day_name1, open_nc, close_nc, get_dim_nc, get_var_nc
64use display,    only: print_msg, LVL_NFO
65use utility,    only: int2str
66
67! DECLARATION
68! -----------
69implicit none
70
71! LOCAL VARIABLES
72! ---------------
73logical(k4) :: here, found
74
75! CODE
76! ----
77! Reading from "startfi.nc" file
78inquire(file = startfi_name,exist = here)
79if (.not. here) call stop_clean(__FILE__,__LINE__,'cannot find required file "'//startfi_name//'"!',1)
80call open_nc(startfi_name,'read')
81call get_dim_nc('physical_points',ngrid)
82call get_dim_nc('subsurface_layers',nsoil_PCM)
83call get_dim_nc('nlayer',nlayer)
84call get_dim_nc('nslope',nslope,found)
85if (.not. found) nslope = 1
86allocate(longitudes(ngrid))
87allocate(latitudes(ngrid))
88allocate(cell_area(ngrid))
89call get_var_nc('longitude',longitudes)
90call get_var_nc('latitude',latitudes)
91call get_var_nc('area',cell_area)
92call close_nc(startfi_name)
93
94! Reading from XIOS daily output file
95inquire(file = xios_day_name1,exist = here)
96if (.not. here) call stop_clean(__FILE__,__LINE__,'cannot find required file "'//xios_day_name1//'"!',1)
97call open_nc(xios_day_name1,'read')
98call get_dim_nc('lon',nlon)
99call get_dim_nc('lat',nlat)
100call get_dim_nc('time_counter',nday)
101call close_nc(xios_day_name1)
102
103! Total surface
104total_surface = sum(cell_area)
105
106! State that dimentions are well initialized
107dim_init = .true.
108call print_msg('> Reading dimensions',LVL_NFO)
109call print_msg('nlon      = '//int2str(nlon),LVL_NFO)
110call print_msg('nlat      = '//int2str(nlat),LVL_NFO)
111call print_msg('ngrid     = '//int2str(ngrid),LVL_NFO)
112call print_msg('nslope    = '//int2str(nslope),LVL_NFO)
113call print_msg('nlayer    = '//int2str(nlayer),LVL_NFO)
114call print_msg('nsoil_PCM = '//int2str(nsoil_PCM),LVL_NFO)
115call print_msg('nsoil     = '//int2str(nsoil),LVL_NFO)
116call print_msg('nday      = '//int2str(nday),LVL_NFO)
117
118! Check consistency of grids
119if (ngrid /= 2 + nlon*(nlat - 2)) call stop_clean(__FILE__,__LINE__,'geometry inconsistency between nlon/nlat and ngrid!',1)
120if (nsoil_PCM > nsoil) call stop_clean(__FILE__,__LINE__,'geometry inconsistency between nsoil and nsoil_PCM! nsoil must be greater than nsoil_PCM.',1)
121
122END SUBROUTINE ini_geometry
123!=======================================================================
124
125!=======================================================================
126SUBROUTINE end_geometry()
127!-----------------------------------------------------------------------
128! NAME
129!     end_geometry
130!
131! DESCRIPTION
132!     Deallocate geometry arrays.
133!
134! AUTHORS & DATE
135!     JB Clement, 12/2025
136!
137! NOTES
138!
139!-----------------------------------------------------------------------
140
141! DECLARATION
142! -----------
143implicit none
144
145! CODE
146! ----
147if (allocated(longitudes)) deallocate(longitudes)
148if (allocated(latitudes)) deallocate(latitudes)
149if (allocated(cell_area)) deallocate(cell_area)
150
151END SUBROUTINE end_geometry
152!=======================================================================
153
154!=======================================================================
155SUBROUTINE lonlat2vect(v_ll,v_vect)
156!-----------------------------------------------------------------------
157! NAME
158!     lonlat2vect
159!
160! DESCRIPTION
161!     Convert data from lon x lat grid (where values at the poles are
162!     duplicated) into a vector.
163!
164! AUTHORS & DATE
165!     JB Clement, 12/2025
166!
167! NOTES
168!     The longitudes -180 and +180 are not duplicated like in the PCM
169!     dynamics.
170!-----------------------------------------------------------------------
171
172! DEPENDENCIES
173! ------------
174use stoppage, only: stop_clean
175
176! DECLARATION
177! -----------
178implicit none
179
180! ARGUMENTS
181! ---------
182real(dp), dimension(:,:), intent(in)  :: v_ll
183real(dp), dimension(:),   intent(out) :: v_vect
184
185! LOCAL VARIABLES
186! ---------------
187integer(di) :: ilon, ilat, n_lon, n_lat, n_grd
188
189! CODE
190! ----
191n_grd = size(v_vect)
192n_lon = size(v_ll,1)
193n_lat = size(v_ll,2)
194if (n_grd == 1) then ! 1D case
195    if (n_lon /= 1 .or. n_lat /= 1) call stop_clean(__FILE__,__LINE__,'problem of dimensions!',1)
196    v_vect(1) = v_ll(1,1)
197else ! 3D
198    ! Check
199    if (n_grd /= n_lon*(n_lat - 2) + 2) call stop_clean(__FILE__,__LINE__,'problem of dimensions!',1)
200
201    ! Treatment of the poles
202    v_vect(1) = v_ll(1,1)
203    v_vect(n_grd) = v_ll(1,n_lat)
204
205    ! Treatment of regular points
206    do ilat = 2,n_lat - 1
207        do ilon = 1,n_lon
208            v_vect(1 + ilon + (ilat - 2)*n_lon) = v_ll(ilon,ilat)
209        end do
210    end do
211end if
212
213END SUBROUTINE lonlat2vect
214!=======================================================================
215
216!=======================================================================
217SUBROUTINE vect2lonlat(v_vect,v_ll)
218!-----------------------------------------------------------------------
219! NAME
220!     vect2lonlat
221!
222! DESCRIPTION
223!     Convert data from a vector into lon x lat grid (where values
224!     at the poles are duplicated).
225!
226! AUTHORS & DATE
227!     JB Clement, 12/2025
228!
229! NOTES
230!     The longitudes -180 and +180 are not duplicated like in the PCM
231!     dynamics.
232!-----------------------------------------------------------------------
233
234! DEPENDENCIES
235! ------------
236use stoppage, only: stop_clean
237
238! DECLARATION
239! -----------
240implicit none
241
242! ARGUMENTS
243! ---------
244real(dp), dimension(:),   intent(in)  :: v_vect
245real(dp), dimension(:,:), intent(out) :: v_ll
246
247! LOCAL VARIABLES
248! ---------------
249integer(di) :: ilon, ilat, n_lon, n_lat, n_grd
250
251! CODE
252! ----
253n_grd = size(v_vect)
254n_lon = size(v_ll,1)
255n_lat = size(v_ll,2)
256if (n_grd == 1) then ! 1D case
257    if (n_lon /= 1 .or. n_lat /= 1) call stop_clean(__FILE__,__LINE__,'problem of dimensions!',1)
258    v_ll(1,1) = v_vect(1)
259else ! 3D
260    ! Check
261    if (n_grd /= n_lon*(n_lat - 2) + 2) call stop_clean(__FILE__,__LINE__,'problem of dimensions!',1)
262
263    ! Treatment of the poles
264    v_ll(:,1) = v_vect(1)
265    v_ll(:,n_lat) = v_vect(n_grd)
266
267    ! Treatment of regular points
268    do ilat = 2,n_lat - 1
269        do ilon = 1,n_lon
270            v_ll(ilon,ilat) = v_vect(1 + ilon + (ilat - 2)*n_lon)
271        end do
272    end do
273end if
274
275END SUBROUTINE vect2lonlat
276!=======================================================================
277
278!=======================================================================
279SUBROUTINE dyngrd2vect(v_dyn,v_vect,extensive)
280!-----------------------------------------------------------------------
281! NAME
282!     dyngrd2vect
283!
284! DESCRIPTION
285!     Convert data from dynamical lon x lat grid (where poles are
286!     duplicated and longitude -180°/+180° is duplicated) into a vector.
287!
288! AUTHORS & DATE
289!     JB Clement, 12/2025
290!
291! NOTES
292!     The longitudes -180 and +180 are duplicated (PCM dynamics).
293!     For extensive variables, the values at the poles must be averaged
294!     over the number of longitudes for computation on the dynamical grid.
295!-----------------------------------------------------------------------
296
297! DEPENDENCIES
298! ------------
299use stoppage, only: stop_clean
300
301! DECLARATION
302! -----------
303implicit none
304
305! ARGUMENTS
306! ---------
307real(dp), dimension(:,:), intent(in)           :: v_dyn
308logical(k4),              intent(in), optional :: extensive
309real(dp), dimension(:),   intent(out)          :: v_vect
310
311! LOCAL VARIABLES
312! ---------------
313integer(di) :: ilon, ilat, n_lon, n_lat, n_grd
314real(dp)    :: factor
315
316! CODE
317! ----
318n_grd = size(v_vect)
319n_lon = size(v_dyn,1)
320n_lat = size(v_dyn,2)
321if (n_grd == 1) then ! 1D case
322    if (n_lon /= 1 .or. n_lat /= 1) call stop_clean(__FILE__,__LINE__,'problem of dimensions!',1)
323    v_vect(1) = v_dyn(1,1)
324else ! 3D
325    ! Check
326    if (n_grd /= (n_lon - 1)*(n_lat - 2) + 2) call stop_clean(__FILE__,__LINE__,'problem of dimensions!',1)
327
328    ! Treatment of the poles
329    factor = 1._dp
330    if (present(extensive)) then
331        ! Correction for polar values of extensive variable when mapped onto the physical grid
332        if (extensive) factor = real(n_lon - 1,dp)
333    end if
334    v_vect(1) = v_dyn(1,1)*factor
335    v_vect(n_grd) = v_dyn(1,n_lat)*factor
336
337    ! Treatment of regular points
338    do ilat = 2,n_lat - 1
339        do ilon = 1,n_lon - 1
340            v_vect(1 + ilon + (ilat - 2)*(n_lon - 1)) = v_dyn(ilon,ilat)
341        end do
342    end do
343end if
344
345END SUBROUTINE dyngrd2vect
346!=======================================================================
347
348!=======================================================================
349SUBROUTINE vect2dyngrd(v_vect,v_dyn,extensive)
350!-----------------------------------------------------------------------
351! NAME
352!     vect2dyngrd
353!
354! DESCRIPTION
355!     Convert data from a vector into lon x lat grid (where poles are
356!     duplicated and longitude -180°/+180° is duplicated).
357!
358! AUTHORS & DATE
359!     JB Clement, 12/2025
360!
361! NOTES
362!     The longitudes -180 and +180 are duplicated (PCM dynamics).
363!     For extensive variables, the values at the poles must be averaged
364!     over the number of longitudes to be computed/outputted on the
365!     dynamical grid.
366!-----------------------------------------------------------------------
367
368! DEPENDENCIES
369! ------------
370use stoppage, only: stop_clean
371
372! DECLARATION
373! -----------
374implicit none
375
376! ARGUMENTS
377! ---------
378real(dp), dimension(:),   intent(in)           :: v_vect
379logical(k4),              intent(in), optional :: extensive
380real(dp), dimension(:,:), intent(out)          :: v_dyn
381
382! LOCAL VARIABLES
383! ---------------
384integer(di) :: ilon, ilat, n_lon, n_lat, n_grd
385real(dp)    :: factor
386
387! CODE
388! ----
389n_grd = size(v_vect)
390n_lon = size(v_dyn,1)
391n_lat = size(v_dyn,2)
392if (n_grd == 1) then ! 1D case
393    if (n_lon /= 1 .or. n_lat /= 1) call stop_clean(__FILE__,__LINE__,'problem of dimensions!',1)
394    v_dyn(1,1) = v_vect(1)
395else ! 3D
396    ! Check
397    if (n_grd /= (n_lon - 1)*(n_lat - 2) + 2) call stop_clean(__FILE__,__LINE__,'problem of dimensions!',1)
398
399    ! Treatment of the poles
400    factor = 1._dp
401    if (present(extensive)) then
402        ! Correction for polar values of extensive variable when mapped onto the dynamical grid
403        if (extensive) factor = 1._dp/real(n_lon - 1,dp)
404    end if
405    v_dyn(:,1) = v_vect(1)*factor
406    v_dyn(:,n_lat) = v_vect(n_grd)*factor
407
408    ! Treatment of regular points
409    do ilat = 2,n_lat - 1
410        do ilon = 1,n_lon - 1
411            v_dyn(ilon,ilat) = v_vect(1 + ilon + (ilat - 2)*(n_lon - 1))
412        end do
413        v_dyn(n_lon,ilat) =  v_dyn(1,ilat) ! Treatment of the duplicated longitude
414    end do
415end if
416
417END SUBROUTINE vect2dyngrd
418!=======================================================================
419
420END MODULE geometry
Note: See TracBrowser for help on using the repository browser.