File: grib_f90_head.f90

package info (click to toggle)
eccodes 2.43.0-1
  • links: PTS, VCS
  • area: main
  • in suites: forky, sid
  • size: 420,752 kB
  • sloc: cpp: 153,709; sh: 20,505; ansic: 12,673; f90: 6,854; python: 3,011; perl: 2,744; javascript: 1,427; yacc: 854; fortran: 468; lex: 359; makefile: 141
file content (157 lines) | stat: -rw-r--r-- 7,264 bytes parent folder | download
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
! (C) Copyright 2005- ECMWF.
!
! This software is licensed under the terms of the Apache Licence Version 2.0
! which can be obtained at http://www.apache.org/licenses/LICENSE-2.0.
!
! In applying this licence, ECMWF does not waive the privileges and immunities granted to it by
! virtue of its status as an intergovernmental organisation nor does it submit to any jurisdiction.
!
!
!> Module grib_api
!>
!> The grib_api module provides the Fortran 90 interface of the GRIB API.
module grib_api

  implicit none
  include "grib_kinds.h"
  include "grib_api_constants.h"
  include "grib_api_externals.h"
  include "grib_api_visibility.h"

  real(8), parameter, public :: GRIB_MISSING_DOUBLE = -1.D+100
  integer(4), parameter, public :: GRIB_MISSING_LONG = 2147483647

  !> Create a new message in memory from an integer or character array containing the coded message.
  !>
  !> The message can be accessed through its gribid and it will be available\n
  !> until @ref grib_release is called. A reference to the original coded\n
  !> message is kept in the new message structure.
  !>
  !> In case of error, if the status parameter (optional) is not given, the program will
  !> exit with an error message.\n Otherwise the error message can be
  !> gathered with @ref grib_get_error_string.
  !>
  !>
  !> \b Examples: \ref grib_copy_message.f90 "grib_copy_message.f90"
  !>
  !> @param gribid      id of the grib loaded in memory
  !> @param message     array containing the coded message
  !> @param status      GRIB_SUCCESS if OK, integer value on error
  interface grib_new_from_message
    module procedure grib_new_from_message_int4
    module procedure grib_new_from_message_char
  end interface grib_new_from_message
  

  !> Create a message pointing to an character array containing the coded message.
  !>
  !> The message can be accessed through its gribid and it will be available\n
  !> until @ref grib_release is called or (attention) the character array is deallocated!
  !>
  !> In case of error, if the status parameter (optional) is not given, the program will
  !> exit with an error message.\n Otherwise the error message can be
  !> gathered with @ref grib_get_error_string.
  !>
  !> @param gribid      id of the grib loaded in memory
  !> @param message     array containing the coded message
  !> @param status      GRIB_SUCCESS if OK, integer value on error  
  interface grib_new_from_message_no_copy
    module procedure grib_new_from_message_no_copy_int4
    module procedure grib_new_from_message_no_copy_char
  end interface grib_new_from_message_no_copy


  !> Get a value of specified index from an array key.
  !>
  !> Given a gribid and key name as input a value corresponding to the given index
  !> is returned. The index is zero based i.e. the first element has
  !> zero index, the second element index one and so on.
  !> If the parameter index is an array all the values corresponding to the indexes
  !> list is returned.
  !> The gribid references to a grib message loaded in memory.
  !>
  !> In case of error, if the status parameter (optional) is not given, the program will
  !> exit with an error message.\n Otherwise the error message can be
  !> gathered with @ref grib_get_error_string.
  !>
  !> \b Examples: \ref grib_nearest.f90 "grib_nearest.f90"
  !>
  !> @see grib_new_from_file, grib_release, grib_get
  !>
  !> @param[in] gribid    id of the GRIB loaded in memory
  !> @param[in] key       key name
  !> @param[in] index     index can be a scalar or array of integer(4)
  !> @param[out] value    value can be a scalar or array of integer(4),real(4),real(8)
  !> @param[out] status   GRIB_SUCCESS if OK, integer value on error
  interface grib_get_element
    module procedure grib_get_real4_element, &
      grib_get_real8_element, &
      grib_get_real4_elements, &
      grib_get_real8_elements
  end interface grib_get_element

  !> Find the nearest point/points of a given latitude/longitude point.
  !>
  !> The value in the nearest point (or the four nearest points) is returned as well as the
  !> zero based index (which can be used in @ref grib_get_element)
  !> and its distance from the given point using the following
  !> formula radius * acos( sin(lat1)*sin(lat2)+cos(lat1)*cos(lat2)*cos(lon1-lon2) ).
  !>
  !> If the is_lsm flag is .true. the input field gribid is considered as
  !> a land sea mask and the nearest land point is returned.\n
  !> The nearest land point among the four neighbours is:
  !> - the nearest point with land sea mask value >= 0.5
  !> - the nearest without any other condition if all the four have land sea mask value <0.5.
  !>
  !> Arrays (real(8)) of latitude/longitude can be provided to find with one call
  !> the values,indexes and distances for all the lat/lon points listed in the arrays.
  !>
  !> If a single latitude/longitude point is provided and outlat,outlon,value,distance,index
  !> are defined as arrays with four elements the lat/lon coordinates and values, distances
  !> and indexes of the four nearest points are returned.
  !>
  !> In case of error, if the status parameter (optional) is not given, the program will
  !> exit with an error message.\n Otherwise the error message can be
  !> gathered with @ref grib_get_error_string.
  !>
  !> \b Examples: \ref grib_nearest.f90 "grib_nearest.f90"
  !>
  !> @param[in] gribid     id of the grib loaded in memory
  !> @param[in] is_lsm      .true. if the nearest land point is required otherwise .false.
  !> @param[in] inlat       latitude of the point in degrees
  !> @param[in] inlon       longitudes of the point in degrees
  !> @param[out] outlat     latitude of the nearest point in degrees
  !> @param[out] outlon     longitude of the nearest point in degrees
  !> @param[out] distance   distance between the given point and its nearest (in km)
  !> @param[out] index      zero based index
  !> @param[out] value      value of the field in the nearest point
  !> @param[out] status     GRIB_SUCCESS if OK, integer value on error
  interface grib_find_nearest
    module procedure grib_find_nearest_single, &
      grib_find_nearest_four_single, &
      grib_find_nearest_multiple
  end interface grib_find_nearest

  !> Get latitude/longitude and data values.
  !>
  !> Latitudes, longitudes, data values arrays are returned.
  !> They must be properly allocated by the caller and their required
  !> dimension can be obtained with \ref grib_get_size or by getting (with \ref grib_get)
  !> the value of the integer key "numberOfPoints".
  !>
  !> In case of error, if the status parameter (optional) is not given, the program will
  !> exit with an error message.\n Otherwise the error message can be
  !> gathered with @ref grib_get_error_string.
  !>
  !> \b Examples: \ref grib_get_data.f90 "grib_get_data.f90"
  !>
  !> @param[in] gribid       id of the grib loaded in memory
  !> @param[out] lats        latitudes array with dimension "size"
  !> @param[out] lons        longitudes array with dimension "size"
  !> @param[out] values      data values array with dimension "size"
  !> @param[out] status      GRIB_SUCCESS if OK, integer value on error
  interface grib_get_data
    module procedure grib_get_data_real4, &
      grib_get_data_real8
  end interface grib_get_data