FMS  2025.04
Flexible Modeling System
get_grid_version.F90
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 !> @defgroup get_grid_version_mod get_grid_version_mod
19 !> @ingroup data_override
20 !> @brief get_grid implementations and helper routines for @ref data_override_mod
21 
22 !> @addtogroup get_grid_version_mod
23 !> @{
24 module get_grid_version_mod
25 use constants_mod, only: deg_to_rad
26 use platform_mod, only: r4_kind, r8_kind, fms_path_len
27 use mpp_mod, only : mpp_error,fatal,note, mpp_min, mpp_max
28 use mpp_domains_mod, only : domain2d, operator(.NE.),operator(.EQ.)
29 use mpp_domains_mod, only : mpp_get_compute_domain, mpp_get_global_domain
30 use fms2_io_mod, only : fmsnetcdfdomainfile_t, fmsnetcdffile_t, open_file, close_file, &
31  variable_exists, read_data, get_variable_size, get_variable_num_dimensions
32 use mosaic2_mod, only : get_mosaic_tile_grid
33 
34 implicit none
35 
37  module procedure get_grid_version_1_r4
38  module procedure get_grid_version_1_r8
39 end interface get_grid_version_1
40 
42  module procedure get_grid_version_2_r4
43  module procedure get_grid_version_2_r8
44 end interface get_grid_version_2
45 
46 contains
47 
48 !> Get lon and lat of three model (target) grids from grid_spec.nc
49 subroutine check_grid_sizes(domain_name, Domain, nlon, nlat)
50 character(len=12), intent(in) :: domain_name
51 type (domain2d), intent(in) :: Domain
52 integer, intent(in) :: nlon, nlat
53 
54 character(len=184) :: error_message
55 integer :: xsize, ysize
56 
57 call mpp_get_global_domain(domain, xsize=xsize, ysize=ysize)
58 if(nlon .NE. xsize .OR. nlat .NE. ysize) then
59  error_message = 'Error in data_override_init. Size of grid as specified by '// &
60  ' does not conform to that specified by grid_spec.nc.'// &
61  ' From : by From grid_spec.nc: by '
62  error_message( 59: 70) = domain_name
63  error_message(130:141) = domain_name
64  write(error_message(143:146),'(i4)') xsize
65  write(error_message(150:153),'(i4)') ysize
66  write(error_message(174:177),'(i4)') nlon
67  write(error_message(181:184),'(i4)') nlat
68  call mpp_error(fatal,error_message)
69 endif
70 end subroutine check_grid_sizes
71 
72 #include "get_grid_version_r4.fh"
73 #include "get_grid_version_r8.fh"
74 
75 end module get_grid_version_mod
76 !> @}
77 ! close documentation grouping
Close a netcdf or domain file opened with open_file or open_virtual_file.
Definition: fms2_io.F90:165
Opens a given netcdf or domain file.
Definition: fms2_io.F90:121
Read data from a defined field in a file.
Definition: fms2_io.F90:291
subroutine check_grid_sizes(domain_name, Domain, nlon, nlat)
Get lon and lat of three model (target) grids from grid_spec.nc.
subroutine, public get_mosaic_tile_grid(grid_file, fileobj, domain, tile_count)
Gets the name of a mosaic tile grid file.
Definition: mosaic2.F90:428
These routines retrieve the axis specifications associated with the compute domains....
These routines retrieve the axis specifications associated with the global domains....
The domain2D type contains all the necessary information to define the global, compute and data domai...
Error handler.
Definition: mpp.F90:381
Reduction operations. Find the max of scalar a from the PEs in pelist result is also automatically br...
Definition: mpp.F90:537
Reduction operations. Find the min of scalar a from the PEs in pelist result is also automatically br...
Definition: mpp.F90:559