-
Notifications
You must be signed in to change notification settings - Fork 1
/
Copy pathmakethickness.f90
177 lines (112 loc) · 4.29 KB
/
makethickness.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
program soildepth
! gfortran -o makethickness makethickness.f90 -I/usr/local/include -L/usr/local/lib -lnetcdff
! calculate soil depth using files from the Pelletier et al. (2016) files
use iso_fortran_env
use netcdf
implicit none
integer, parameter :: i1 = int8
integer, parameter :: i2 = int16
integer, parameter :: sp = real32
character(100) :: lowlandfile ! = 'upland_valley-bottom_and_lowland_sedimentary_deposit_thickness.nc' (byte nd: -1)
character(100) :: uplandfile ! = 'upland_hill-slope_soil_thickness.nc' (float nd: -1)
character(100) :: slpfracfile ! = 'hill-slope_valley-bottom.nc' (float nd: not defined, but is -1)
character(100) :: outfile
integer(i1), allocatable, dimension(:,:) :: lowland
real(sp), allocatable, dimension(:,:) :: upland
real(sp), allocatable, dimension(:,:) :: landformf
real(sp), allocatable, dimension(:,:) :: thickness
integer :: status
integer :: ncid
integer :: dimid
integer :: varid
integer :: xlen
integer :: ylen
integer :: x
integer :: y
real(sp), dimension(2) :: actual_range
! ---------
call getarg(1,lowlandfile)
call getarg(2,uplandfile)
call getarg(3,slpfracfile)
status = nf90_open(lowlandfile,nf90_nowrite,ncid)
if (status /= nf90_noerr) call handle_err(status)
status = nf90_inq_dimid(ncid,'lon',dimid)
if (status == nf90_ebaddim) status = nf90_inq_dimid(ncid,'x',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_ebaddim) status = nf90_inq_dimid(ncid,'y',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)
! ---
allocate(lowland(xlen,ylen))
status = nf90_inq_varid(ncid,'Band1',varid)
if (status /= nf90_noerr) call handle_err(status)
status = nf90_get_var(ncid,varid,lowland)
if (status /= nf90_noerr) call handle_err(status)
status = nf90_close(ncid)
if (status /= nf90_noerr) call handle_err(status)
! ---
allocate(upland(xlen,ylen))
status = nf90_open(uplandfile,nf90_nowrite,ncid)
if (status /= nf90_noerr) call handle_err(status)
status = nf90_inq_varid(ncid,'Band1',varid)
if (status /= nf90_noerr) call handle_err(status)
status = nf90_get_var(ncid,varid,upland)
if (status /= nf90_noerr) call handle_err(status)
status = nf90_close(ncid)
if (status /= nf90_noerr) call handle_err(status)
! ---
allocate(landformf(xlen,ylen))
status = nf90_open(slpfracfile,nf90_nowrite,ncid)
if (status /= nf90_noerr) call handle_err(status)
status = nf90_inq_varid(ncid,'Band1',varid)
if (status /= nf90_noerr) call handle_err(status)
status = nf90_get_var(ncid,varid,landformf)
if (status /= nf90_noerr) call handle_err(status)
status = nf90_close(ncid)
if (status /= nf90_noerr) call handle_err(status)
! ---
allocate(thickness(xlen,ylen))
thickness = -9999.
do y = 1,ylen
do x = 1,xlen
if (lowland(x,y) < 0 .and. upland(x,y) < 0.) cycle
if (upland(x,y) < 0.) then
thickness(x,y) = real(lowland(x,y))
else if (lowland(x,y) < 0) then
thickness(x,y) = upland(x,y)
else
thickness(x,y) = upland(x,y) * landformf(x,y) + real(lowland(x,y)) * (1. - landformf(x,y))
end if
end do
end do
actual_range(1) = minval(thickness,mask=thickness /= -9999.)
actual_range(2) = maxval(thickness,mask=thickness /= -9999.)
write(0,*)'soil depth actual range',actual_range
! ---
call getarg(4,outfile)
status = nf90_open(outfile,nf90_write,ncid)
if (status /= nf90_noerr) call handle_err(status)
status = nf90_inq_varid(ncid,'thickness',varid)
if (status /= nf90_noerr) call handle_err(status)
status = nf90_put_var(ncid,varid,thickness)
if (status /= nf90_noerr) call handle_err(status)
status = nf90_put_att(ncid,varid,'actual_range',actual_range)
if (status /= nf90_noerr) call handle_err(status)
status = nf90_close(ncid)
if (status /= nf90_noerr) call handle_err(status)
! ---------
contains
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 soildepth