-
Notifications
You must be signed in to change notification settings - Fork 2
/
nn_cf_net.F90
288 lines (218 loc) · 9.73 KB
/
nn_cf_net.F90
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
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
module nn_cf_net_mod
!! neural_net convection emulator
!! Module containing code pertaining only to the Neural Net functionalities of
!! the M2LiNES Convection Flux parameterisation.
!---------------------------------------------------------------------
! Libraries to use
use netcdf
use precision, only: sp
implicit none
private
!---------------------------------------------------------------------
! public interfaces
public relu, net_forward, nn_cf_net_init, nn_cf_net_finalize
!---------------------------------------------------------------------
! local/private data
! Neural Net Parameters
integer :: n_in
!! Combined length of input features
integer :: n_features_out
!! Number of output features (variables)
integer, allocatable, dimension(:) :: feature_out_sizes
!! Vector storing the length of each of the input features
integer :: n_lev
!! number of atmospheric layers (model levels) output is supplied for
! Dimension of each hidden layer
integer :: n_h1
integer :: n_h2
integer :: n_h3
integer :: n_h4
integer :: n_out
! Weights at each hidden layer of the NN
real(sp), allocatable, dimension(:,:) :: r_w1
real(sp), allocatable, dimension(:,:) :: r_w2
real(sp), allocatable, dimension(:,:) :: r_w3
real(sp), allocatable, dimension(:,:) :: r_w4
real(sp), allocatable, dimension(:,:) :: r_w5
! Biases at each hidden layer of the NN
real(sp), allocatable, dimension(:) :: r_b1
real(sp), allocatable, dimension(:) :: r_b2
real(sp), allocatable, dimension(:) :: r_b3
real(sp), allocatable, dimension(:) :: r_b4
real(sp), allocatable, dimension(:) :: r_b5
! Scale factors for inputs
real(sp), allocatable, dimension(:) :: xscale_mean
real(sp), allocatable, dimension(:) :: xscale_stnd
! vectors at each hidden layer of the NN
real(sp), allocatable, dimension(:) :: z1
real(sp), allocatable, dimension(:) :: z2
real(sp), allocatable, dimension(:) :: z3
real(sp), allocatable, dimension(:) :: z4
! Scale factors for outputs
real(sp), allocatable, dimension(:) :: yscale_mean
real(sp), allocatable, dimension(:) :: yscale_stnd
!---------------------------------------------------------------------
! Functions and Subroutines
contains
!-----------------------------------------------------------------
! Public Subroutines
subroutine relu(logits)
!! Applies ReLU to a vector.
real(sp), dimension(:), intent(inout) :: logits
!! vector to which ReLU will be applied
where (logits .lt. 0.0) logits = 0.0
end subroutine relu
subroutine net_forward(features, logits)
!! Run forward method of the Neural Net.
real(sp), dimension(:) :: features
!! Vector of input features
real(sp), dimension(:), intent(out) :: logits
!! Output vector
integer :: i
!! Loop counter
integer :: out_dim_counter
!!
integer :: out_pos, feature_size, f
!!
! Initialise hidden layer values
! TODO Is this really neccessary...???
z1 = 0.
z2 = 0.
z3 = 0.
z4 = 0.
!Normalize features
features = (features - xscale_mean) / xscale_stnd
z1 = matmul(features, r_w1) + r_b1
call relu(z1)
z2 = matmul(z1, r_w2) + r_b2
call relu(z2)
z3 = matmul(z2, r_w3) + r_b3
call relu(z3)
z4 = matmul(z3, r_w4) + r_b4
call relu(z4)
logits = matmul(z4, r_w5) + r_b5
! Apply scaling and normalisation of each output feature in logits
out_pos = 0
do f = 1, n_features_out
feature_size = feature_out_sizes(f)
logits(out_pos+1:out_pos+feature_size) = (logits(out_pos+1:out_pos+feature_size)*yscale_stnd(f)) + yscale_mean(f)
out_pos = out_pos + feature_size
end do
end subroutine
subroutine nn_cf_net_init(nn_filename, n_inputs, n_outputs, nrf)
!! Initialise the neural net
integer, intent(out) :: n_inputs, n_outputs
integer, intent(in) :: nrf
!! number of atmospheric layers in each input
! This will be the netCDF ID for the file and data variable.
integer :: ncid
integer :: in_dimid, out_dimid, single_dimid
integer :: h1_dimid, h2_dimid, h3_dimid, h4_dimid
integer :: r_w1_varid, r_w2_varid, r_b1_varid, r_b2_varid
integer :: r_w3_varid, r_w4_varid, r_b3_varid, r_b4_varid
integer :: r_w5_varid, r_b5_varid
integer :: xscale_mean_varid, xscale_stnd_varid
integer :: yscale_mean_varid, yscale_stnd_varid
character(len=136), intent(in) :: nn_filename
!! NetCDF filename from which to read model weights
!-------------allocate arrays and read data-------------------
! Open the file. NF90_NOWRITE tells netCDF we want read-only
! access
! Get the varid or dimid for each variable or dimension based
! on its name.
call check( nf90_open(trim(nn_filename),NF90_NOWRITE,ncid ))
call check( nf90_inq_dimid(ncid, 'N_in', in_dimid))
call check( nf90_inquire_dimension(ncid, in_dimid, len=n_in))
call check( nf90_inq_dimid(ncid, 'N_h1', h1_dimid))
call check( nf90_inquire_dimension(ncid, h1_dimid, len=n_h1))
call check( nf90_inq_dimid(ncid, 'N_h2', h2_dimid))
call check( nf90_inquire_dimension(ncid, h2_dimid, len=n_h2))
call check( nf90_inq_dimid(ncid, 'N_h3', h3_dimid))
call check( nf90_inquire_dimension(ncid, h3_dimid, len=n_h3))
call check( nf90_inq_dimid(ncid, 'N_h4', h4_dimid))
call check( nf90_inquire_dimension(ncid, h4_dimid, len=n_h4))
call check( nf90_inq_dimid(ncid, 'N_out', out_dimid))
call check( nf90_inquire_dimension(ncid, out_dimid, len=n_out))
call check( nf90_inq_dimid(ncid, 'N_out_dim', out_dimid))
call check( nf90_inquire_dimension(ncid, out_dimid, len=n_features_out))
print *, 'size of features to NN', n_in
print *, 'size of outputs from NN', n_out
! Allocate arrays for NN based on sizes read from nc file
allocate(r_w1(n_in,n_h1))
allocate(r_w2(n_h1,n_h2))
allocate(r_w3(n_h2,n_h3))
allocate(r_w4(n_h3,n_h4))
allocate(r_w5(n_h4,n_out))
allocate(r_b1(n_h1))
allocate(r_b2(n_h2))
allocate(r_b3(n_h3))
allocate(r_b4(n_h4))
allocate(r_b5(n_out))
allocate(z1(n_h1))
allocate(z2(n_h2))
allocate(z3(n_h3))
allocate(z4(n_h4))
allocate(xscale_mean(n_in))
allocate(xscale_stnd(n_in))
allocate(yscale_mean(n_features_out))
allocate(yscale_stnd(n_features_out))
allocate(feature_out_sizes(n_features_out))
! Read NN data in from nc file
call check( nf90_inq_varid(ncid, "w1", r_w1_varid))
call check( nf90_get_var(ncid, r_w1_varid, r_w1))
call check( nf90_inq_varid(ncid, "w2", r_w2_varid))
call check( nf90_get_var(ncid, r_w2_varid, r_w2))
call check( nf90_inq_varid(ncid, "w3", r_w3_varid))
call check( nf90_get_var(ncid, r_w3_varid, r_w3))
call check( nf90_inq_varid(ncid, "w4", r_w4_varid))
call check( nf90_get_var(ncid, r_w4_varid, r_w4))
call check( nf90_inq_varid(ncid, "w5", r_w5_varid))
call check( nf90_get_var(ncid, r_w5_varid, r_w5))
call check( nf90_inq_varid(ncid, "b1", r_b1_varid))
call check( nf90_get_var(ncid, r_b1_varid, r_b1))
call check( nf90_inq_varid(ncid, "b2", r_b2_varid))
call check( nf90_get_var(ncid, r_b2_varid, r_b2))
call check( nf90_inq_varid(ncid, "b3", r_b3_varid))
call check( nf90_get_var(ncid, r_b3_varid, r_b3))
call check( nf90_inq_varid(ncid, "b4", r_b4_varid))
call check( nf90_get_var(ncid, r_b4_varid, r_b4))
call check( nf90_inq_varid(ncid, "b5", r_b5_varid))
call check( nf90_get_var(ncid, r_b5_varid, r_b5))
call check( nf90_inq_varid(ncid,"fscale_mean", xscale_mean_varid))
call check( nf90_get_var( ncid, xscale_mean_varid,xscale_mean ))
call check( nf90_inq_varid(ncid,"fscale_stnd", xscale_stnd_varid))
call check( nf90_get_var( ncid, xscale_stnd_varid,xscale_stnd ))
call check( nf90_inq_varid(ncid,"oscale_mean", yscale_mean_varid))
call check( nf90_get_var( ncid, yscale_mean_varid,yscale_mean ))
call check( nf90_inq_varid(ncid,"oscale_stnd", yscale_stnd_varid))
call check( nf90_get_var( ncid, yscale_stnd_varid,yscale_stnd ))
! Close the nc file
call check( nf90_close(ncid))
write(*, *) 'Finished reading NN regression file.'
n_inputs = n_in
n_outputs = n_out
! Set sizes of output features based on nrf
n_lev = nrf
! Set sizes of the feature groups
feature_out_sizes(:) = n_lev
feature_out_sizes(2:3) = n_lev-1
end subroutine nn_cf_net_init
subroutine nn_cf_net_finalize()
!! Clean up NN space by deallocating any arrays.
! Deallocate arrays
deallocate(r_w1, r_w2, r_w3, r_w4, r_w5)
deallocate(r_b1, r_b2, r_b3, r_b4, r_b5)
deallocate(z1, z2, z3, z4)
deallocate(xscale_mean, xscale_stnd, yscale_mean, yscale_stnd)
deallocate(feature_out_sizes)
end subroutine nn_cf_net_finalize
subroutine check(err_status)
!! Check error status after netcdf call and print message for
!! error codes.
integer, intent(in) :: err_status
!! error status from nf90 function
if(err_status /= nf90_noerr) then
write(*, *) trim(nf90_strerror(err_status))
end if
end subroutine check
end module nn_cf_net_mod