MONC
Functions/Subroutines | Variables
mcclatchey_profiles Module Reference

Functions/Subroutines

subroutine read_mcclatchey_profiles (current_state, mcc)
 
subroutine check_mcc_status (status)
 Will check for a NetCDF status of the McClatchey profile file and write to log_log error any decoded statuses. More...
 
subroutine read_mcc_dimensions (ncid, plev_dim)
 Reads the dimensions from the NetCDF file. More...
 
subroutine read_mcc_variables (filename, ncid, plev_dim, MCC_PROFILE_KEY, temp_mcc_3d, plev)
 
subroutine read_single_mcc_variable (ncid, key, data1d, data3d)
 Reads a single variable out of a NetCDF file. More...
 
subroutine calculate_radiation_levels (current_state, mcc)
 

Variables

character(len= *), parameter pressure_key = "plev"
 NetCDF data plev key in McClatchey profile. More...
 
character(len= *), parameter temperature_profile_key = "t"
 NetCDF data for temperature. More...
 
character(len= *), parameter vapour_profile_key = "q"
 NetCDF data for vapour. More...
 
character(len= *), parameter ozone_profile_key = "o3"
 NetCDF data for ozone. More...
 

Function/Subroutine Documentation

◆ calculate_radiation_levels()

subroutine mcclatchey_profiles::calculate_radiation_levels ( type(model_state_type), intent(inout), target  current_state,
type(str_mcc_profiles), intent(inout)  mcc 
)

Definition at line 235 of file mcclatchey_profiles.F90.

236  ! This routine calculates the total number of levels to be
237  ! used in the radiation calcs (socrates). This includes the
238  ! number of levels in the MONC testcase plus all the
239  ! McClattchey levels which fill in the top of atmosphere.
240  ! The number of levels is dependent on the McClattchey profile
241  ! selected in the MONC configuration
242  type(model_state_type), target, intent(inout) :: current_state
243 
244  type(str_mcc_profiles), intent(inout) :: mcc
245 
246  real(kind=default_precision), parameter :: &
247  fracdp = 0.75, &
248  ! fraction of to dp which the 1st rad leevl can be from top
249  mindp = 2000., &
250  ! Minimum values dp in Pa
251  maxdp = 10000.
252  ! Maximum dp we will allow from top
253  real(kind=default_precision) :: pmindiff, pcutoff
254  real(kind=default_precision), allocatable :: prefn_loc(:) !(current_state%local_grid%size(Z_INDEX))
255 
256  integer :: k
257 
258  allocate(prefn_loc(current_state%local_grid%size(z_index)))
259 
260  prefn_loc(:) = current_state%global_grid%configuration%vertical%prefn(:)
261 
262  pmindiff = min(max(fracdp* &
263  (prefn_loc(current_state%local_grid%size(z_index)-1) - &
264  prefn_loc(current_state%local_grid%size(z_index))), &
265  mindp), maxdp)
266  pcutoff = prefn_loc(current_state%local_grid%size(z_index)) - pmindiff
267 
268  mcc%cut = 0
269  do k=1, mcc%levs
270  if (mcc%p_n(k) .lt. pcutoff) mcc%cut = k
271  enddo
272 
273  mcc%irad_levs = current_state%local_grid%size(z_index) + mcc%cut
274  mcc%n_levels = mcc%irad_levs - 1
275 
276  deallocate(prefn_loc)
277 

◆ check_mcc_status()

subroutine mcclatchey_profiles::check_mcc_status ( integer, intent(in)  status)

Will check for a NetCDF status of the McClatchey profile file and write to log_log error any decoded statuses.

Parameters
statusThe NetCDF status flag

Definition at line 144 of file mcclatchey_profiles.F90.

145  integer, intent(in) :: status
146 
147  if (status /= nf90_noerr) then
148  call log_log(log_error, "NetCDF returned error code of "//trim(nf90_strerror(status)))
149  end if
Here is the call graph for this function:
Here is the caller graph for this function:

◆ read_mcc_dimensions()

subroutine mcclatchey_profiles::read_mcc_dimensions ( integer, intent(in)  ncid,
integer, intent(out)  plev_dim 
)

Reads the dimensions from the NetCDF file.

Parameters
ncidThe NetCDF file id
time_dimNumber of elements in the time dimension

Definition at line 158 of file mcclatchey_profiles.F90.

159  integer, intent(in) :: ncid
160  integer, intent(out) :: plev_dim
161  integer :: plev_dimid
162 
163  call check_mcc_status(nf90_inq_dimid(ncid, pressure_key, plev_dimid))
164 
165  call check_mcc_status(nf90_inquire_dimension(ncid, plev_dimid, len=plev_dim))
166 
Here is the call graph for this function:
Here is the caller graph for this function:

◆ read_mcc_variables()

subroutine mcclatchey_profiles::read_mcc_variables ( character(*), intent(in)  filename,
integer, intent(in)  ncid,
integer, intent(in)  plev_dim,
character(*), intent(in)  MCC_PROFILE_KEY,
real(kind=default_precision), dimension(:,:,:), intent(inout), allocatable  temp_mcc_3d,
real(kind=default_precision), dimension(:), intent(inout), optional, allocatable  plev 
)

Definition at line 169 of file mcclatchey_profiles.F90.

171 
172  character(*), intent(in) :: filename
173  character(*), intent(in) :: MCC_PROFILE_KEY
174 
175  real(kind=default_precision), dimension(:), allocatable, intent(inout), optional :: plev
176  real(kind=default_precision), dimension(:,:,:), allocatable, intent(inout) :: temp_mcc_3d
177 
178  integer, intent(in) :: ncid, plev_dim
179  integer :: status, variable_id
180 
181  ! Do some checking on the variable contents so that we can deal with different
182  ! variable names or missing variables
183 
184  ! pressure levels, only done once with the first call...
185  if (present(plev)) then
186  status=nf90_inq_varid(ncid, pressure_key, variable_id)
187  if (status==nf90_noerr)then
188  allocate(plev(plev_dim))
189  call read_single_mcc_variable(ncid, pressure_key, data1d=plev)
190  else
191  call log_log(log_error, "No recognized pressure variable found in"//trim(filename))
192  end if
193  endif
194 
195  status=nf90_inq_varid(ncid, mcc_profile_key, variable_id)
196  if (status==nf90_noerr)then
197  call read_single_mcc_variable(ncid, mcc_profile_key, data3d=temp_mcc_3d)
198  else
199  call log_log(log_error, "No recognized pressure variable found in"//trim(filename))
200  end if
201 
Here is the call graph for this function:
Here is the caller graph for this function:

◆ read_mcclatchey_profiles()

subroutine mcclatchey_profiles::read_mcclatchey_profiles ( type(model_state_type), intent(inout), target  current_state,
type(str_mcc_profiles), intent(inout)  mcc 
)

Definition at line 24 of file mcclatchey_profiles.F90.

25  ! This routine is the top-level routine to read the
26  ! McClatchey (mcc) data from the file. It is assumed the data is
27  ! in netcdf format. In Socrates, the mcc data is in CDL format
28  ! (netcdf ASCII), so it needs to be converted using ncgen
29  !
30  ! This routine also derives centre grid values for the McClatchey profiles
31  ! It is assumed that the staggering of the grid is the same as MONC, but flipped
32  !
33  type(model_state_type), target, intent(inout) :: current_state
34 
35  type(str_mcc_profiles), intent(inout) :: mcc
36 
37  real(kind=default_precision), allocatable :: weight_upper, weight_lower
38 
39  real(kind=default_precision), allocatable :: temp_mcc_3d(:,:,:)
40 
41  integer, parameter :: MAX_FILE_LEN=200
42  character(MAX_FILE_LEN) :: mcc_input_file
43 
44  integer :: ncid, plev_dim, k
45  ! Read McClatchey profile from the global or user config. Seperate
46  ! profile files required from temperature, vapour and ozone. Pressure
47  ! is in all files. After each get_string, there is check status for each
48  ! file.
49 
50  ! read pressure and temperature
51  mcc_input_file = &
52  options_get_string(current_state%options_database, "mcc_temperature_profile")
53  if (log_get_logging_level() .ge. log_debug) then
55  "Reading in McClatchey pressure and temperature from:"//trim(mcc_input_file) )
56  endif
57  call check_mcc_status(nf90_open(path = trim(mcc_input_file), mode = nf90_nowrite, ncid = ncid))
58  call read_mcc_dimensions(ncid, plev_dim)
59  ! set number of McClatchey levels
60  mcc%levs = plev_dim
61  ! check if temp_mcc_3d allocated and allocate if needed
62  if(.not. allocated(temp_mcc_3d)) then
63  allocate(temp_mcc_3d(1,1,plev_dim))
64  endif
65  call read_mcc_variables(trim(mcc_input_file), ncid, plev_dim, &
66  temperature_profile_key, temp_mcc_3d, mcc%p_level)
67  allocate(mcc%t_level(plev_dim))
68  mcc%t_level(:)=temp_mcc_3d(1,1,:)
69  call check_mcc_status(nf90_close(ncid))
70 
71  ! work out the centre grid values of pressure (linear interpolation)
72  allocate(mcc%p_n(plev_dim))
73  do k = 2, mcc%levs
74  mcc%p_n(k)=(0.5*(mcc%p_level(k-1)+(mcc%p_level(k))))
75  enddo
76  mcc%p_n(1)=mcc%p_level(1)+(0.5*(mcc%p_n(2)-mcc%p_level(1)))
77  ! work out the centre grid values of using (log(p) weighting but linear interpolation t)
78  allocate(mcc%t_n(plev_dim))
79  do k = 2, mcc%levs
80  weight_upper = log(mcc%p_level(k))-log(mcc%p_n(k-1))
81  weight_lower = log(mcc%p_level(k))-log(mcc%p_level(k-1))
82  mcc%t_n(k) = (((weight_upper*mcc%t_level(k))) + &
83  (weight_lower*mcc%t_level(k-1)))/(weight_upper+weight_lower)
84  enddo
85  mcc%t_n(1)=mcc%t_level(1)+(0.5*(mcc%t_n(2)-mcc%t_level(1)))
86 
87  ! read mcc vapour field
88  mcc_input_file = &
89  options_get_string(current_state%options_database, "mcc_vapour_profile")
90  if (log_get_logging_level() .ge. log_debug) then
92  "Reading in McClatchey vapour profile from:"//trim(mcc_input_file) )
93  endif
94  call check_mcc_status(nf90_open(path = trim(mcc_input_file), mode = nf90_nowrite, ncid = ncid))
95  call read_mcc_variables(trim(mcc_input_file), ncid, plev_dim, vapour_profile_key, temp_mcc_3d)
96  allocate(mcc%q_level(plev_dim))
97  mcc%q_level(:)=temp_mcc_3d(1,1,:)
98  call check_mcc_status(nf90_close(ncid))
99 
100  ! work out the centre grid values of using (log(p) weighting but log interpolation q and o3)
101  ! change to log of mixing ratio
102  allocate(mcc%q_n(plev_dim))
103  do k = 2, mcc%levs
104  weight_upper = log(mcc%p_level(k))-log(mcc%p_n(k-1))
105  weight_lower = log(mcc%p_level(k))-log(mcc%p_level(k-1))
106  mcc%q_n(k) = (((weight_upper*log(mcc%q_level(k)))) + &
107  (weight_lower*log(mcc%q_level(k-1))))/(weight_upper+weight_lower)
108  mcc%q_n(k) = exp(mcc%q_n(k))
109  enddo
110  mcc%q_n(1)=mcc%q_level(1)+(0.5*((mcc%q_n(2)-mcc%q_level(1))))
111 
112  ! read mcc ozone field
113  mcc_input_file = &
114  options_get_string(current_state%options_database, "mcc_ozone_profile")
115  if (log_get_logging_level() .ge. log_debug) then
116  call log_master_log(log_debug, &
117  "Reading in McClatchey ozone profile from:"//trim(mcc_input_file) )
118  endif
119  call check_mcc_status(nf90_open(path = trim(mcc_input_file), mode = nf90_nowrite, ncid = ncid))
120  call read_mcc_variables(trim(mcc_input_file), ncid, plev_dim, ozone_profile_key, temp_mcc_3d)
121  allocate(mcc%o3_level(plev_dim))
122  mcc%o3_level(:)=temp_mcc_3d(1,1,:)
123  call check_mcc_status(nf90_close(ncid))
124 
125  ! work out the centre grid values of using (log(p) weighting but log interpolation q and o3)
126  ! change to log of mixing ratio
127  allocate(mcc%o3_n(plev_dim))
128  do k = 2, mcc%levs
129  weight_upper = log(mcc%p_level(k))-log(mcc%p_n(k-1))
130  weight_lower = log(mcc%p_level(k))-log(mcc%p_level(k-1))
131  mcc%o3_n(k) = (((weight_upper*log(mcc%o3_level(k)))) + &
132  (weight_lower*log(mcc%o3_level(k-1))))/(weight_upper+weight_lower)
133  mcc%o3_n(k) = exp(mcc%o3_n(k))
134  enddo
135  mcc%o3_n(1)=mcc%o3_level(1)+(0.5*(mcc%o3_n(2)-mcc%o3_level(1)))
136 
137  deallocate(temp_mcc_3d)
138 
Here is the call graph for this function:

◆ read_single_mcc_variable()

subroutine mcclatchey_profiles::read_single_mcc_variable ( integer, intent(in)  ncid,
character(len=*), intent(in)  key,
real(kind=default_precision), dimension(:), intent(inout), optional  data1d,
real(kind=default_precision), dimension(:,:,:), intent(inout), optional  data3d 
)

Reads a single variable out of a NetCDF file.

Parameters
ncidThe NetCDF file id
keyThe variable key (name) to access
data1dOptional one dimensional data to read into
data3dOptional three dimensional data to read into

Definition at line 211 of file mcclatchey_profiles.F90.

212  integer, intent(in) :: ncid
213  character(len=*), intent(in) :: key
214  real(kind=default_precision), dimension(:), intent(inout), optional :: data1d
215  real(kind=default_precision), dimension(:,:,:), intent(inout), optional :: data3d
216 
217  integer :: variable_id
218  real(kind=default_precision), dimension(:,:,:), allocatable :: sdata
219 
220  call check_mcc_status(nf90_inq_varid(ncid, key, variable_id))
221 
222  if (.not. present(data1d) .and. .not. present(data3d)) return
223 
224  if (present(data1d)) then
225  call check_mcc_status(nf90_get_var(ncid, variable_id, data1d))
226  else
227  ! 3D will reshape the data to take account of the column-row major C-F transposition
228  allocate(sdata(size(data3d,1),size(data3d,2), size(data3d,3)))
229  call check_mcc_status(nf90_get_var(ncid, variable_id, sdata))
230  data3d(:,:,:)=reshape(sdata(:,:,:),(/size(data3d,1),size(data3d,2),size(data3d,3)/))
231  deallocate(sdata)
232  end if
Here is the call graph for this function:
Here is the caller graph for this function:

Variable Documentation

◆ ozone_profile_key

character(len=*), parameter mcclatchey_profiles::ozone_profile_key = "o3"

NetCDF data for ozone.

Definition at line 16 of file mcclatchey_profiles.F90.

◆ pressure_key

character(len=*), parameter mcclatchey_profiles::pressure_key = "plev"

NetCDF data plev key in McClatchey profile.

Definition at line 16 of file mcclatchey_profiles.F90.

16  character(len=*), parameter :: &
17  PRESSURE_KEY = "plev" , & !< NetCDF data plev key in McClatchey profile
18  temperature_profile_key = "t" , &
19  vapour_profile_key = "q" , &
20  ozone_profile_key = "o3"

◆ temperature_profile_key

character(len=*), parameter mcclatchey_profiles::temperature_profile_key = "t"

NetCDF data for temperature.

Definition at line 16 of file mcclatchey_profiles.F90.

◆ vapour_profile_key

character(len=*), parameter mcclatchey_profiles::vapour_profile_key = "q"

NetCDF data for vapour.

Definition at line 16 of file mcclatchey_profiles.F90.

logging_mod::log_error
integer, parameter, public log_error
Only log ERROR messages.
Definition: logging.F90:11
logging_mod::log_log
subroutine, public log_log(level, message, str)
Logs a message at the specified level. If the level is above the current level then the message is ig...
Definition: logging.F90:75
logging_mod::log_get_logging_level
integer function, public log_get_logging_level()
Retrieves the current logging level.
Definition: logging.F90:122
logging_mod::log_debug
integer, parameter, public log_debug
Log DEBUG, INFO, WARNING and ERROR messages.
Definition: logging.F90:14
optionsdatabase_mod::options_get_string
character(len=string_length) function, public options_get_string(options_database, key, index)
Retrieves a string value from the database that matches the provided key.
Definition: optionsdatabase.F90:280
logging_mod::log_master_log
subroutine, public log_master_log(level, message)
Will log just from the master process.
Definition: logging.F90:47
datadefn_mod::default_precision
integer, parameter, public default_precision
MPI communication type which we use for the prognostic and calculation data.
Definition: datadefn.F90:17