forked from ARVE-Research/makesoil
-
Notifications
You must be signed in to change notification settings - Fork 0
/
Copy pathsoilcalc.f90
386 lines (251 loc) · 9.09 KB
/
soilcalc.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
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
program soilcalc
! use Makefile
use netcdf
use parametersmod, only : i1,sp
use soilpropertiesmod, only : omcf,soildata,soilproperties
! use simplesoilmod, only : i1,i2,sp,soildata,layerinfo,simplesoil
! use pedotransfermod, only : omcf
implicit none
character(200) :: soilfile
! character(200) :: landfracfile
character(200) :: outfile
integer :: status
integer :: ncid
integer :: dimid
integer :: varid
integer :: xlen
integer :: ylen
integer :: nl
integer :: x
integer :: y
integer :: l
integer(i1), allocatable, dimension(:,:) :: usda
real(sp), allocatable, dimension(:,:,:) :: sand
real(sp), allocatable, dimension(:,:,:) :: silt
real(sp), allocatable, dimension(:,:,:) :: clay
real(sp), allocatable, dimension(:,:,:) :: cfvo
real(sp), allocatable, dimension(:,:,:) :: soc
real(sp), allocatable, dimension(:,:,:) :: bulk
real(sp), allocatable, dimension(:,:,:) :: Tsat
real(sp), allocatable, dimension(:,:,:) :: T33
real(sp), allocatable, dimension(:,:,:) :: T1500
real(sp), allocatable, dimension(:,:,:) :: whc
real(sp), allocatable, dimension(:,:,:) :: Ksat
real(sp), allocatable, dimension(:) :: zpos
real(sp), allocatable, dimension(:) :: dz
real(sp), allocatable, dimension(:,:) :: datacheck
type(soildata) :: soil
real(sp), parameter :: rmissing = -9999.
real(sp), dimension(2) :: actual_range
real(sp) :: scale
! -------------------------------------------------------
call getarg(1,soilfile)
status = nf90_open(soilfile,nf90_nowrite,ncid)
if (status /= nf90_noerr) call handle_err(status)
status = nf90_inq_dimid(ncid,'lon',dimid)
if (status /= nf90_noerr) call handle_err(status)
status = nf90_inquire_dimension(ncid,dimid,len=xlen)
if (status /= nf90_noerr) call handle_err(status)
status = nf90_inq_dimid(ncid,'lat',dimid)
if (status /= nf90_noerr) call handle_err(status)
status = nf90_inquire_dimension(ncid,dimid,len=ylen)
if (status /= nf90_noerr) call handle_err(status)
status = nf90_inq_dimid(ncid,'depth',dimid)
if (status /= nf90_noerr) call handle_err(status)
status = nf90_inquire_dimension(ncid,dimid,len=nl)
if (status /= nf90_noerr) call handle_err(status)
allocate(zpos(nl))
allocate(dz(nl))
allocate(soil%layer(nl))
status = nf90_inq_varid(ncid,'depth',varid)
if (status /= nf90_noerr) call handle_err(status)
status = nf90_get_var(ncid,varid,zpos)
if (status /= nf90_noerr) call handle_err(status)
zpos = abs(zpos)
status = nf90_inq_varid(ncid,'dz',varid)
if (status /= nf90_noerr) call handle_err(status)
status = nf90_get_var(ncid,varid,dz)
if (status /= nf90_noerr) call handle_err(status)
soil%layer%zpos = zpos
soil%layer%dz = dz
do l = 1,nl
write(0,*)l,zpos(l),dz(l)
end do
allocate(usda(xlen,ylen))
allocate(sand(xlen,ylen,nl))
allocate(silt(xlen,ylen,nl))
allocate(clay(xlen,ylen,nl))
allocate(cfvo(xlen,ylen,nl))
allocate(soc(xlen,ylen,nl))
allocate(datacheck(5,nl))
! ---------
! read input soil spatial data
status = nf90_inq_varid(ncid,'USDA',varid)
if (status /= nf90_noerr) call handle_err(status)
status = nf90_get_var(ncid,varid,usda)
if (status /= nf90_noerr) call handle_err(status)
call getvar('sand',sand)
call getvar('silt',silt)
call getvar('clay',clay)
call getvar('cfvo',cfvo)
call getvar('soc',soc)
! call getvar('bdod',bulk)
status = nf90_close(ncid)
if (status /= nf90_noerr) call handle_err(status)
! ---------
allocate(bulk(xlen,ylen,nl))
allocate(Tsat(xlen,ylen,nl))
allocate(T33(xlen,ylen,nl))
allocate(T1500(xlen,ylen,nl))
allocate(whc(xlen,ylen,nl))
allocate(Ksat(xlen,ylen,nl))
bulk = rmissing
Tsat = rmissing
T33 = rmissing
T1500 = rmissing
whc = rmissing
Ksat = rmissing
write(0,*)'calculating'
where (sand < 0.) sand = rmissing
do y = 1,ylen
do x = 1,xlen
soil%usda = usda(x,y)
if (soil%usda <= 3) cycle ! skip all cells with no soil classification
do l = 1,nl
datacheck(:,l) = [sand(x,y,l),silt(x,y,l),clay(x,y,l),cfvo(x,y,l),soc(x,y,l)]
end do
if (any(datacheck == rmissing)) cycle ! skip cells with missing data
soil%layer%sand = sand(x,y,:)
soil%layer%silt = silt(x,y,:)
soil%layer%clay = clay(x,y,:)
soil%layer%cfvo = cfvo(x,y,:)
soil%layer%orgm = soc(x,y,:) * omcf
do l = 1,nl
if (soil%layer(l)%sand + soil%layer(l)%silt + soil%layer(l)%clay > 1.) then
scale = 1. / (soil%layer(l)%sand + soil%layer(l)%silt + soil%layer(l)%clay)
soil%layer(l)%sand = soil%layer(l)%sand * scale
soil%layer(l)%silt = soil%layer(l)%silt * scale
soil%layer(l)%clay = soil%layer(l)%clay * scale
end if
end do
! write(0,*)x,y,datacheck
! call simplesoil(soil)
call soilproperties(soil)
bulk(x,y,:) = soil%layer%bulk
Tsat(x,y,:) = soil%layer%Tsat
T33(x,y,:) = soil%layer%T33
T1500(x,y,:) = soil%layer%T1500
whc(x,y,:) = soil%layer%whc
Ksat(x,y,:) = soil%layer%Ksat
if (any(soil%layer%whc <= 0)) then
write(0,*)x,y
do l = 1,nl
write(0,*)l,soil%layer(l)%sand,soil%layer(l)%clay,soil%layer(l)%whc
end do
end if
end do
end do
! ---------
write(0,*)'writing'
call getarg(2,outfile)
status = nf90_open(outfile,nf90_write,ncid)
if (status /= nf90_noerr) call handle_err(status)
! --
! bulk density
status = nf90_inq_varid(ncid,'bulk',varid)
if (status /= nf90_noerr) call handle_err(status)
status = nf90_put_var(ncid,varid,bulk)
if (status /= nf90_noerr) call handle_err(status)
actual_range = [minval(bulk,mask=bulk/=rmissing),maxval(bulk,mask=bulk/=rmissing)]
status = nf90_put_att(ncid,varid,'actual_range',actual_range)
if (status /= nf90_noerr) call handle_err(status)
write(0,*)'Bulk density range: ',actual_range
! --
! Tsat
status = nf90_inq_varid(ncid,'Tsat',varid)
if (status /= nf90_noerr) call handle_err(status)
status = nf90_put_var(ncid,varid,Tsat)
if (status /= nf90_noerr) call handle_err(status)
actual_range = [minval(Tsat,mask=Tsat/=rmissing),maxval(Tsat,mask=Tsat/=rmissing)]
status = nf90_put_att(ncid,varid,'actual_range',actual_range)
if (status /= nf90_noerr) call handle_err(status)
write(0,*)'Tsat range: ',actual_range
! --
! T33
status = nf90_inq_varid(ncid,'T33',varid)
if (status /= nf90_noerr) call handle_err(status)
status = nf90_put_var(ncid,varid,T33)
if (status /= nf90_noerr) call handle_err(status)
actual_range = [minval(T33,mask=T33/=rmissing),maxval(T33,mask=T33/=rmissing)]
status = nf90_put_att(ncid,varid,'actual_range',actual_range)
if (status /= nf90_noerr) call handle_err(status)
write(0,*)'T33 range: ',actual_range
! --
! T1500
status = nf90_inq_varid(ncid,'T1500',varid)
if (status /= nf90_noerr) call handle_err(status)
status = nf90_put_var(ncid,varid,T1500)
if (status /= nf90_noerr) call handle_err(status)
actual_range = [minval(T1500,mask=T1500/=rmissing),maxval(T1500,mask=T1500/=rmissing)]
status = nf90_put_att(ncid,varid,'actual_range',actual_range)
if (status /= nf90_noerr) call handle_err(status)
write(0,*)'T1500 range: ',actual_range
! --
! water holding capacity
status = nf90_inq_varid(ncid,'whc',varid)
if (status /= nf90_noerr) call handle_err(status)
status = nf90_put_var(ncid,varid,whc)
if (status /= nf90_noerr) call handle_err(status)
actual_range = [minval(whc,mask=whc/=rmissing),maxval(whc,mask=whc/=rmissing)]
status = nf90_put_att(ncid,varid,'actual_range',actual_range)
if (status /= nf90_noerr) call handle_err(status)
write(0,*)'whc range: ',actual_range
! --
! Ksat
status = nf90_inq_varid(ncid,'Ksat',varid)
if (status /= nf90_noerr) call handle_err(status)
actual_range = [minval(Ksat,mask=Ksat/=rmissing),maxval(Ksat,mask=Ksat/=rmissing)]
write(0,*)'Ksat range: ',actual_range
status = nf90_put_att(ncid,varid,'actual_range',actual_range)
if (status /= nf90_noerr) call handle_err(status)
status = nf90_put_var(ncid,varid,ksat)
if (status /= nf90_noerr) call handle_err(status)
! -------------------------------------------------------
contains
subroutine getvar(name,var)
use parametersmod, only : i2,sp
implicit none
character(*), intent(in) :: name
real(sp), dimension(:,:,:), intent(out) :: var
real(sp) :: scale_factor
integer(i2) :: missing_value
integer(i2), allocatable, dimension(:,:,:) :: ivar
integer :: xlen,ylen,nl
! ----
xlen = size(var,dim=1)
ylen = size(var,dim=2)
nl = size(var,dim=3)
allocate(ivar(xlen,ylen,nl))
status = nf90_inq_varid(ncid,name,varid)
if (status /= nf90_noerr) call handle_err(status)
status = nf90_get_var(ncid,varid,ivar)
if (status /= nf90_noerr) call handle_err(status)
status = nf90_get_att(ncid,varid,'missing_value',missing_value)
if (status /= nf90_noerr) call handle_err(status)
status = nf90_get_att(ncid,varid,'scale_factor',scale_factor)
if (status /= nf90_noerr) call handle_err(status)
var = rmissing
where (ivar /= missing_value) var = real(ivar) * scale_factor
write(0,*)trim(name),': ',minval(var,mask=var/=rmissing),maxval(var,mask=var/=rmissing)
end subroutine getvar
! ---------
subroutine handle_err(status)
! Internal subroutine - checks error status after each netcdf call,
! prints out text message each time an error code is returned.
integer, intent (in) :: status
if(status /= nf90_noerr) then
write(0,*)'NetCDF error: ',trim(nf90_strerror(status))
stop
end if
end subroutine handle_err
end program soilcalc