FMS  2025.04
Flexible Modeling System
column_diagnostics.inc
1 !***********************************************************************
2 !* Apache License 2.0
3 !*
4 !* This file is part of the GFDL Flexible Modeling System (FMS).
5 !*
6 !* Licensed under the Apache License, Version 2.0 (the "License");
7 !* you may not use this file except in compliance with the License.
8 !* You may obtain a copy of the License at
9 !*
10 !* http://www.apache.org/licenses/LICENSE-2.0
11 !*
12 !* FMS is distributed in the hope that it will be useful, but WITHOUT
13 !* WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied;
14 !* without even the implied warranty of MERCHANTABILITY or FITNESS FOR A
15 !* PARTICULAR PURPOSE. See the License for the specific language
16 !* governing permissions and limitations under the License.
17 !***********************************************************************
18 !> @file
19 
20 !> @brief initialize_diagnostic_columns returns the (i, j, lat, lon) coord-
21 !! inates of any diagnostic columns that are located on the current
22 !! processor.
23 subroutine initialize_diagnostic_columns_ &
24  (module, num_diag_pts_latlon, num_diag_pts_ij, &
25  global_i , global_j , global_lat_latlon, &
26  global_lon_latlon, lonb_in, latb_in, &
27  do_column_diagnostics, &
28  diag_lon, diag_lat, diag_i, diag_j, diag_units)
29 
30 !---------------------------------------------------------------------
31 ! initialize_diagnostic_columns returns the (i, j, lat, lon) coord-
32 ! inates of any diagnostic columns that are located on the current
33 ! processor.
34 !----------------------------------------------------------------------
35 
36 !---------------------------------------------------------------------
37 character(len=*), intent(in) :: module !< module calling this subroutine
38 integer, intent(in) :: num_diag_pts_latlon !< number of diagnostic columns specified
39  !! by lat-lon coordinates
40 integer, intent(in) :: num_diag_pts_ij !< number of diagnostic columns specified
41  !! by global (i,j) coordinates
42 integer, dimension(:), intent(in) :: global_i !< specified global i coordinates
43 integer, dimension(:), intent(in) :: global_j !< specified global j coordinates
44 real(FMS_CD_KIND_), dimension(:), intent(in) :: global_lat_latlon !< specified global lat coordinates
45 real(FMS_CD_KIND_), dimension(:), intent(in) :: global_lon_latlon !< specified global lon coordinates
46 real(FMS_CD_KIND_), dimension(:,:), intent(in) :: lonb_in, latb_in
47 logical, dimension(:,:), intent(out) :: do_column_diagnostics !< is a diagnostic column in this jrow ?
48 integer, dimension(:), intent(inout) :: diag_i !< processor i indices of diagnstic columns
49 integer, dimension(:), intent(inout) :: diag_j !< processor j indices of diagnstic columns
50 real(FMS_CD_KIND_), dimension(:), intent(out) :: diag_lat !< latitudes of diagnostic columns [ degrees ]
51 real(FMS_CD_KIND_), dimension(:), intent(out) :: diag_lon !< longitudes of diagnostic columns [ degrees ]
52 integer, dimension(:), intent(out) :: diag_units !< unit number for each diagnostic column
53 !---------------------------------------------------------------------
54 
55 !---------------------------------------------------------------------
56 ! intent(in) variables:
57 !
58 ! module module calling this subroutine
59 ! num_diag_pts_latlon number of diagnostic columns specified
60 ! by lat-lon coordinates
61 ! num_diag_pts_ij number of diagnostic columns specified
62 ! by global (i,j) coordinates
63 ! global_i specified global i coordinates
64 ! global_j specified global j coordinates
65 ! global_lat_latlon specified global lat coordinates
66 ! global_lon_latlon specified global lon coordinates
67 !
68 ! intent(out) variables:
69 !
70 ! do_column_diagnostics is a diagnostic column in this jrow ?
71 ! diag_i processor i indices of diagnstic columns
72 ! diag_j processor j indices of diagnstic columns
73 ! diag_lat latitudes of diagnostic columns
74 ! [ degrees ]
75 ! diag_lon longitudes of diagnostic columns
76 ! [ degrees ]
77 ! diag_units unit number for each diagnostic column
78 !
79 !---------------------------------------------------------------------
80 
81 !--------------------------------------------------------------------
82 ! local variables:
83 
84  real(FMS_CD_KIND_), dimension(size(diag_i,1)) :: global_lat !< latitudes for all diagnostic columns [ degrees ]
85  real(FMS_CD_KIND_), dimension(size(diag_i,1)) :: global_lon !< longitudes for all diagnostic columns [ degrees ]
86  real(FMS_CD_KIND_), dimension(size(latb_in,1)-1, size(latb_in,2)-1) :: &
87  distance, distance_x, distance_y, &
88  distance_x2, distance2
89  real(FMS_CD_KIND_), dimension(size(latb_in,1), size(latb_in,2)) :: latb_deg
90  real(FMS_CD_KIND_), dimension(size(lonb_in,1), size(lonb_in,2)) :: lonb_deg
91  real(FMS_CD_KIND_) :: dellat, dellon
92  real(FMS_CD_KIND_) :: latb_max, latb_min, lonb_max, lonb_min
93 
94  integer :: num_diag_pts !< total number of diagnostic columns
95  integer :: i !< do loop indices
96  integer :: j !< do loop indices
97  integer :: nn !< do loop indices
98  real(FMS_CD_KIND_) :: ref_lat
99  real(FMS_CD_KIND_) :: current_distance
100  character(len=8) :: char !< character string for diaganostic column index
101  character(len=FMS_FILE_LEN) :: filename !< filename for output file for diagnostic column
102  logical :: allow_ij_input
103  logical :: open_file
104  integer :: io
105 
106  integer, parameter :: lkind=fms_cd_kind_
107  real(FMS_CD_KIND_) :: tmp
108 !--------------------------------------------------------------------
109 ! local variables:
110 !
111 ! global_lat latitudes for all diagnostic columns [ degrees ]
112 ! global_lon longitudes for all diagnostic columns
113 ! [ degrees ]
114 ! num_diag_pts total number of diagnostic columns
115 ! i, j, nn do loop indices
116 ! char character string for diaganostic column index
117 ! filename filename for output file for diagnostic column
118 !
119 !---------------------------------------------------------------------
120 
121  if (.not. module_is_initialized) call column_diagnostics_init
122 
123 !--------------------------------------------------------------------
124 ! save the input lat and lon fields. define the delta of latitude
125 ! and longitude.
126 !--------------------------------------------------------------------
127  latb_deg = real( real(latb_in,r8_kind)*radian, fms_cd_kind_) !< unit conversion in r8_kind
128  lonb_deg = real( real(lonb_in,r8_kind)*radian, fms_cd_kind_ ) !< unit conversion in r8_kind
129  dellat = latb_in(1,2) - latb_in(1,1)
130  dellon = lonb_in(2,1) - lonb_in(1,1)
131  latb_max = maxval(latb_deg(:,:))
132  latb_min = minval(latb_deg(:,:))
133  lonb_max = maxval(lonb_deg(:,:))
134  lonb_min = minval(lonb_deg(:,:))
135  if (lonb_min < 10.0_lkind .or. lonb_max > 350.0_lkind) then
136  lonb_min = 0.0_lkind
137  lonb_max = 360.0_lkind
138  endif
139 
140  allow_ij_input = .true.
141  ref_lat = latb_in(1,1)
142  do i =2,size(latb_in,1)
143  if (latb_in(i,1) /= ref_lat) then
144  allow_ij_input = .false.
145  exit
146  endif
147  end do
148 
149  if ( .not. allow_ij_input .and. num_diag_pts_ij /= 0) then
150  call error_mesg ('column_diagnostics_mod', &
151  'cannot specify column diagnostics column with (i,j) &
152  &coordinates when using cubed sphere -- must specify &
153  & lat/lon coordinates', fatal)
154  endif
155 
156 !----------------------------------------------------------------------
157 
158 !----------------------------------------------------------------------
159 ! initialize column_diagnostics flag and diag unit numbers. define
160 ! total number of diagnostic columns.
161 !----------------------------------------------------------------------
162  do_column_diagnostics = .false.
163  diag_units(:) = -1
164  diag_i(:) = -99
165  diag_j(:) = -99
166  diag_lat(:) = -999.0_lkind
167  diag_lon(:) = -999.0_lkind
168  num_diag_pts = size(diag_i(:))
169 
170 !--------------------------------------------------------------------
171 ! define an array of lat-lon values for all diagnostic columns.
172 !--------------------------------------------------------------------
173  do nn = 1, num_diag_pts_latlon
174  global_lat(nn) = global_lat_latlon(nn)
175  global_lon(nn) = global_lon_latlon(nn)
176  end do
177 
178  do nn = 1, num_diag_pts_ij
179  tmp = (-0.5_lkind*acos(-1.0_lkind) + 0.5_lkind*dellat) + real(global_j(nn)-1,fms_cd_kind_)*dellat
180  global_lat(nn+num_diag_pts_latlon) = real( real(tmp,r8_kind)*radian, fms_cd_kind_ )
181  tmp = 0.5_lkind*dellon + real(global_i(nn)-1,fms_cd_kind_)*dellon
182  global_lon(nn+num_diag_pts_latlon) = real( real(tmp,r8_kind)*radian, fms_cd_kind_ )
183  end do
184 
185 !----------------------------------------------------------------------
186 ! loop over all diagnostic points to check for their presence on
187 ! this processor.
188 !----------------------------------------------------------------------
189  do nn=1,num_diag_pts
190  open_file = .false.
191 
192 !----------------------------------------------------------------------
193 ! verify that the values of lat and lon are valid.
194 !----------------------------------------------------------------------
195  if (global_lon(nn) >= 0.0_lkind .and. &
196  global_lon(nn) <= 360.0_lkind) then
197  else
198  call error_mesg ('column_diagnostics_mod', &
199  ' invalid longitude', fatal)
200  endif
201  if (global_lat(nn) >= -90.0_lkind .and. &
202  global_lat(nn) <= 90.0_lkind) then
203  else
204  call error_mesg ('column_diagnostics_mod', &
205  ' invalid latitude', fatal)
206  endif
207 
208 !--------------------------------------------------------------------
209 ! if the desired diagnostics column is within the current
210 ! processor's domain, define the total and coordinate distances from
211 ! each of the processor's grid points to the diagnostics point.
212 !--------------------------------------------------------------------
213 
214  if (global_lat(nn) .ge. latb_min .and. &
215  global_lat(nn) .le. latb_max) then
216  if (global_lon(nn) .ge. lonb_min .and.&
217  global_lon(nn) .le. lonb_max) then
218  do j=1,size(latb_deg,2) - 1
219  do i=1,size(lonb_deg,1) - 1
220  distance_y(i,j) = abs(global_lat(nn) - latb_deg(i,j))
221  distance_x(i,j) = abs(global_lon(nn) - lonb_deg(i,j))
222  distance_x2(i,j) = abs((global_lon(nn)-360.0_lkind) - lonb_deg(i,j))
223  distance(i,j) = (global_lat(nn)-latb_deg(i,j))**2 + (global_lon(nn)-lonb_deg(i,j))**2
224  distance2(i,j) = (global_lat(nn)-latb_deg(i,j))**2 + ((global_lon(nn)-360.0_lkind) - lonb_deg(i,j))**2
225  end do
226  end do
227 
228 !--------------------------------------------------------------------
229 ! find the grid point on the processor that is within the specified
230 ! critical distance and also closest to the requested diagnostics
231 ! column. save the (i,j) coordinates and (lon,lat) of this model
232 ! grid point. set a flag indicating that a disgnostics file should
233 ! be opened on this processor for this diagnostic point.
234 !--------------------------------------------------------------------
235  current_distance = distance(1,1)
236  do j=1,size(latb_deg,2) - 1
237  do i=1,size(lonb_deg,1) - 1
238  if (distance_x(i,j) <= real(crit_xdistance,fms_cd_kind_) .and. &
239  distance_y(i,j) <= real(crit_ydistance,fms_cd_kind_)) then
240  if (distance(i,j) < current_distance) then
241  current_distance = distance(i,j)
242  do_column_diagnostics(i,j) = .true.
243  diag_j(nn) = j
244  diag_i(nn) = i
245  diag_lon(nn) = lonb_deg(i,j)
246  diag_lat(nn) = latb_deg(i,j)
247  open_file = .true.
248  endif
249  endif
250 
251 !---------------------------------------------------------------------
252 ! check needed because of the 0.0 / 360.0 longitude periodicity.
253 !---------------------------------------------------------------------
254  if (distance_x2(i,j)<= real(crit_xdistance,fms_cd_kind_) .and. &
255  distance_y(i,j) <= real(crit_ydistance,fms_cd_kind_)) then
256  if (distance2(i,j) < current_distance) then
257  current_distance = distance2(i,j)
258  do_column_diagnostics(i,j) = .true.
259  diag_j(nn) = j
260  diag_i(nn) = i
261  diag_lon(nn) = lonb_deg(i,j)
262  diag_lat(nn) = latb_deg(i,j)
263  open_file = .true.
264  endif
265  endif
266  end do
267  end do
268 
269 !--------------------------------------------------------------------
270 ! if the point has been found on this processor, open a diagnostics
271 ! file.
272 !--------------------------------------------------------------------
273  if (open_file) then
274  write (char, '(i2)') nn
275  filename = trim(module) // '_point' // &
276  trim(adjustl(char)) // '.out'
277  if(mpp_npes() > 10000) then
278  write( filename,'(a,i6.6)' )trim(filename)//'.', mpp_pe()-mpp_root_pe()
279  else
280  write( filename,'(a,i4.4)' )trim(filename)//'.', mpp_pe()-mpp_root_pe()
281  endif
282  open(newunit=diag_units(nn), file=trim(filename), action='WRITE', position='rewind', iostat=io)
283  if(io/=0) call error_mesg ('column_diagnostics_mod', 'Error in opening file '//trim(filename), fatal)
284  endif ! (open_file)
285  endif
286  endif
287  end do
288 
289 !---------------------------------------------------------------------
290 
291 
292 end subroutine initialize_diagnostic_columns_
293 
294 
295 
296 
297 !####################################################################
298 !> @brief column_diagnostics_header writes out information concerning
299 !! time and location of following data into the column_diagnostics
300 !! output file.
301 subroutine column_diagnostics_header_ &
302  (module, diag_unit, time, nn, diag_lon, &
303  diag_lat, diag_i, diag_j)
304 
305 !--------------------------------------------------------------------
306 ! column_diagnostics_header writes out information concerning
307 ! time and location of following data into the column_diagnostics
308 ! output file.
309 !--------------------------------------------------------------------
310 
311 !--------------------------------------------------------------------
312 character(len=*), intent(in) :: module !< module name calling this subroutine
313 type(time_type), intent(in) :: Time !< current model time [ time_type ]
314 integer, intent(in) :: diag_unit !< unit number for column_diagnostics output
315 integer, intent(in) :: nn !< index of diagnostic column currently active
316 real(FMS_CD_KIND_), dimension(:), intent(in) :: diag_lon !< longitude of current diagnostic column [ degrees ]
317 real(FMS_CD_KIND_), dimension(:), intent(in) :: diag_lat !< latitude of current diagnostic column [ degrees ]
318 integer, dimension(:), intent(in) :: diag_i !< i coordinate of current diagnostic column
319 integer, dimension(:), intent(in) :: diag_j !< j coordinate of current diagnostic column
320 
321 !--------------------------------------------------------------------
322 ! intent(in) variables
323 !
324 ! module module name calling this subroutine
325 ! Time current model time [ time_type ]
326 ! diag_unit unit number for column_diagnostics output
327 ! nn index of diagnostic column currently active
328 ! diag_lon longitude of current diagnostic column [ degrees ]
329 ! diag_lat latitude of current diagnostic column [ degrees ]
330 ! diag_i i coordinate of current diagnostic column
331 ! diag_j j coordinate of current diagnostic column
332 !
333 !---------------------------------------------------------------------
334 
335 
336 !--------------------------------------------------------------------
337 ! local variables:
338 
339  integer :: year !< integers defining the current time
340  integer :: month !< integers defining the current time
341  integer :: day !< integers defining the current time
342  integer :: hour !< integers defining the current time
343  integer :: minute !< integers defining the current time
344  integer :: second !< integers defining the current time
345  character(len=9) :: mon !< character string for the current month
346  character(len=64) :: header !< title for the output
347 
348 !--------------------------------------------------------------------
349 ! local variables:
350 !
351 ! year, month, day, hour, minute, seconds
352 ! integers defining the current time
353 ! mon character string for the current month
354 ! header title for the output
355 !
356 !--------------------------------------------------------------------
357 
358  if (.not. module_is_initialized) call column_diagnostics_init
359 
360 !--------------------------------------------------------------------
361 ! convert the time type to a date and time for printing. convert
362 ! month to a character string.
363 !--------------------------------------------------------------------
364  call get_date (time, year, month, day, hour, minute, second)
365  mon = month_name(month)
366 
367 !---------------------------------------------------------------------
368 ! write timestamp and column location information to the diagnostic
369 ! columns output unit.
370 !---------------------------------------------------------------------
371  write (diag_unit,'(a)') ' '
372  write (diag_unit,'(a)') ' '
373  write (diag_unit,'(a)') &
374  '======================================================'
375  write (diag_unit,'(a)') ' '
376  header = ' PRINTING ' // module // ' DIAGNOSTICS'
377  write (diag_unit,'(a)') header
378  write (diag_unit,'(a)') ' '
379  write (diag_unit,'(a, i6,2x, a,i4,i4,i4,i4)') ' time stamp:', &
380  year, trim(mon), day, &
381  hour, minute, second
382  write (diag_unit,'(a, i4)') &
383  ' DIAGNOSTIC POINT COORDINATES, point #', nn
384  write (diag_unit,'(a)') ' '
385  write (diag_unit,'(a,f8.3,a,f8.3)') ' longitude = ', &
386  diag_lon(nn), ' latitude = ', diag_lat(nn)
387  write (diag_unit,'(a, i6, a,i6,a,i6)') &
388  ' on processor # ', mpp_pe(), &
389  ' : processor i =', diag_i(nn), &
390  ' , processor j =', diag_j(nn)
391  write (diag_unit,'(a)') ' '
392 
393 !---------------------------------------------------------------------
394 
395 end subroutine column_diagnostics_header_
396 !@}
integer function mpp_npes()
Returns processor count for current pelist.
Definition: mpp_util.inc:420
integer function mpp_pe()
Returns processor ID.
Definition: mpp_util.inc:406