Skip to content
New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

Netcdf ufs #48

Merged
merged 11 commits into from
Feb 9, 2023
50 changes: 29 additions & 21 deletions README.md
Original file line number Diff line number Diff line change
Expand Up @@ -7,7 +7,14 @@ Patrick Campbell, Zachary Moon, and Wei-Ting Hung

Build canopy model:

Canopy-App requires NetCDF-Fortran Libraries (i.e., -lnetcdf -lnetcdff) for NetCDF I/O Option. See included Makefile for example.
Canopy-App requires NetCDF-Fortran Libraries (i.e., -lnetcdf -lnetcdff) when using the 2D NetCDF I/O Option (i.e., infmt_opt=0).
See included Makefile for example (currently developed using netcdf-c/4.7.4-vh and netcdf-fortran/4.5.3-ff modules).

Compilation Options with or without Debug or NetCDF options:
- `DEBUG = 0(off; default) or DEBUG =1(on)`
- `NETCDF = 0(off) or NETCDF =1(on; default)`

Example: `DEBUG=1 NETCDF=1 make`

Compile, edit namelist, and run canopy model:
- `make`
Expand Down Expand Up @@ -44,8 +51,8 @@ Current Canopy-App components:
Namelist Option : `file_vars` Full name of input file (Supports either text or NetCDF format with following formats:
`.txt`, `.nc`, `.ncf`, or `.nc4`)

- See example file inputs for variables and format (`input_variables.txt`, `input_variables_1D.nc`, or `input_variables_2D.nc`).
- Canopy-App assumes the NetCDF input files are in CF-Convention; recommend using double or float for real variables.
- See example file inputs for variables and format (`input_variables.txt` or `testf000.nc`).
- Canopy-App assumes the NetCDF input files are in CF-Convention and test file is based on UFS-GFSv16; recommend using double or float for real variables.
- Canopy-App can also be run with a single point of 1D input data in a text file (e.g. `input_variables_point.txt`).

**Current Canopy-App Output:** Outputs 3D canopy winds, canopy vertical/eddy diffusivity values, and
Expand All @@ -58,21 +65,22 @@ Current Canopy-App components:

| Variable Name | Variable Description and Units |
| --------------- | ------------------------------------------------- |
| LAT | Latitude (degrees) |
| LON | Longitude (degrees; from 0-360) |
| TIME | Timestamp (days since YYYY-N-D 0:0:0) (NetCDF Only) |
| FH | Forest canopy height (m) |
| HREF | Reference height above canopy (m) |
| WS | Wind speed at HREF (m/s), e.g., 10 m |
| CLU | Forest clumping index (dimensionless) |
| LAI | Leaf area index (m2/m2) |
| VTYPE | Vegetation type (dimensionless), VIIRS Only |
| FFRAC | Forest fraction (dimensionless) |
| UST | Friction velocity (m/s) |
| CSZ | Cosine of the solar zenith angle (dimensionless) |
| Z0 | Total surface roughness length (m) |
| MOL | Monin-Obukhov Length (m) |
| FRP | Total Fire Radiative Power (MW/grid cell area) |
| lat | Latitude (degrees) |
| lon | Longitude (degrees; from 0-360) |
| time | Timestamp (days since YYYY-N-D 0:0:0) (NetCDF Only) |
| fh | Forest canopy height (m) |
| href | Reference height above canopy (m) - 10 m |
| ugrd10m | U wind at HREF (m/s), e.g., 10 m |
| vgrd10m | V wind at HREF (m/s), e.g., 10 m |
| clu | Forest clumping index (dimensionless) |
| lai | Leaf area index (m2/m2) |
| vtype | Vegetation type (dimensionless), VIIRS Only |
| ffrac | Forest fraction (dimensionless) |
| fricv | Friction velocity (m/s) |
| csz | Cosine of the solar zenith angle (dimensionless) |
| sfcr | Total surface roughness length (m) |
| mol | Monin-Obukhov Length (m) |
| frp | Total Fire Radiative Power (MW/grid cell area) |

**Table 2. Current User Namelist Options**

Expand All @@ -93,18 +101,18 @@ Current Canopy-App components:
| lambdars | Value representing influence of roughness sublayer (Massman et al., 2017) |
| dx_opt | 0=Calculation of dx resolution/distance from lon; 1=user set dx grid resolution |
| dx_set | user set real value of grid resolution (m) only if dx_opt=1 |
| flameh_opt | 0=Calculation of flame height from FRP (Byram, 1959); 1=user set flameh; 2=FRP calculation where available (active fires), otherwise user set flameh |
| flameh_opt | 0=Calculation of flame height from FRP (Byram, 1959); 1=user set flameh; 2=FRP calculation where available (active fires), elsewhere user set flameh |
| flameh_set | user set real value of flame height (m) only if flame_opt=1 |
| pai_opt | integer (0=PAI fixed from Katul et al. 2004 veg types-->default; 1=PAI Massman et al. 2017 Eq. 19 calc; 2=PAI from model LAI+WAI; 3=user set PAI value) |
| pai_set | user set real value of PAI (default=4.0; only used if pai_opt=3) |
| lu_opt | integer (0=VIIRS LU type input mapped to Massman et al., currently only option ) |
| lu_opt | integer for input model land use type (0=VIIRS 17 Cat (default) or 1=MODIS-IGBP 20 Cat (valid LU types 1-10 and 12); input mapped to Massman et al.) |
| z0_opt | integer (0=use model input or 1 = vegtype dependent z0 for first estimate) |
| lai_thresh | user set real value of LAI threshold for contiguous canopy (m2/m2) |
| frt_thresh | user set real value of forest fraction threshold for contiguous canopy |
| fch_thresh | user set real value of canopy height threshold for contiguous canopy (m) |
| rsl_opt | user set option to include more explicit stability and Roughness SubLayer (RSL) effects in calculation of wind speed at canopy top (Uc) from reference height 0 = off; 1 = on |

**Note: If modres >> flameh then some error in WAF calculation will be incurred. Suggestion is to use relative fine modres compared to average flame heights (at least <= 0.5 m) if WAF is required.
**Note: If modres >> flameh then some error in WAF calculation will be incurred. Suggestion is to use relative fine modres (at least <= 0.5 m) compared to average flame heights (e.g., ~ 1.0 m) if WAF is required.


Main Citations (further references contained within):
Expand Down
54 changes: 34 additions & 20 deletions canopy_calcs.F90
Original file line number Diff line number Diff line change
Expand Up @@ -51,25 +51,29 @@ SUBROUTINE canopy_calcs
do j=1, nlat

hcmref = variables_2d(i,j)%fh
ubzref = variables_2d(i,j)%ws
uref = variables_2d(i,j)%ugrd10m
vref = variables_2d(i,j)%vgrd10m
cluref = variables_2d(i,j)%clu
lairef = variables_2d(i,j)%lai
vtyperef = variables_2d(i,j)%vtype
ffracref = variables_2d(i,j)%ffrac
ustref = variables_2d(i,j)%ust
ustref = variables_2d(i,j)%fricv
cszref = variables_2d(i,j)%csz
z0ref = variables_2d(i,j)%z0
z0ref = variables_2d(i,j)%sfcr
molref = variables_2d(i,j)%mol
frpref = variables_2d(i,j)%frp
hgtref = variables_2d(i,j)%href

! ... calculate wind speed from u and v
ubzref = sqrt((uref**2.0) + (vref**2.0))

! ... get scaled canopy model profile and sub-canopy layers
zhc = zk/hcmref
cansublays = floor(hcmref/modres)

! ... check for model vegetation types
if (lu_opt .eq. 0 ) then !VIIRS
if (vtyperef .le. 10 .or. vtyperef .eq. 12) then !VIIRS types
! ... check for valid model vegetation types
if (lu_opt .eq. 0 .or. lu_opt .eq. 1 ) then !VIIRS or MODIS
if (vtyperef .gt. 0 .and. vtyperef .le. 10 .or. vtyperef .eq. 12) then !VIIRS or MODIS types

! ... check for contiguous canopy conditions at each model grid cell
if (hcmref .gt. fch_thresh .and. ffracref .gt. frt_thresh &
Expand Down Expand Up @@ -105,7 +109,8 @@ SUBROUTINE canopy_calcs
call canopy_flameh(flameh_opt, flameh_set, dx_2d(i,j), modres, &
frpref, midflamepoint, flameh_2d(i,j))

if (flameh_2d(i,j) .gt. 0.0) then !only calculate WAF when flameh > 0
if (flameh_2d(i,j) .gt. 0.0 .and. flameh_2d(i,j) .le. hcmref) then
!only calculate WAF when flameh > 0 and <= FH
call canopy_waf(hcmref, lambdars, hgtref, flameh_2d(i,j), &
firetype, d_h, zo_h, canBOT(midflamepoint), &
canTOP(midflamepoint), waf_2d(i,j))
Expand All @@ -121,14 +126,16 @@ SUBROUTINE canopy_calcs

! ... user option to calculate in-canopy eddy photolysis attenuation at height z
if (ifcanphot) then
call canopy_phot(fafraczInt, &
lairef, cluref, cszref, rjcf_3d(i,j,:))
if (cszref .ge. 0.0_rk) then !only calculate if cell isn't dark
call canopy_phot(fafraczInt, &
lairef, cluref, cszref, rjcf_3d(i,j,:))
end if
end if

end if !Contiguous Canopy

else
! write(*,*) 'Warning VIIRS VTYPE ', vtyperef, ' is not supported...continue'
! write(*,*) 'Warning VIIRS/MODIS VTYPE ', vtyperef, ' is not supported...continue'
end if !Vegetation types
else
write(*,*) 'Wrong LU_OPT choice of ', lu_opt, ' in namelist...exiting'
Expand Down Expand Up @@ -159,25 +166,29 @@ SUBROUTINE canopy_calcs
! ... Main loop through model grid cells
do loc=1, nlat*nlon
hcmref = variables(loc)%fh
ubzref = variables(loc)%ws
uref = variables(loc)%ugrd10m
vref = variables(loc)%vgrd10m
cluref = variables(loc)%clu
lairef = variables(loc)%lai
vtyperef = variables(loc)%vtype
ffracref = variables(loc)%ffrac
ustref = variables(loc)%ust
ustref = variables(loc)%fricv
cszref = variables(loc)%csz
z0ref = variables(loc)%z0
z0ref = variables(loc)%sfcr
molref = variables(loc)%mol
frpref = variables(loc)%frp
hgtref = variables(loc)%href

! ... calculate wind speed from u and v
ubzref = sqrt((uref**2.0) + (vref**2.0))

! ... get scaled canopy model profile and sub-canopy layers
zhc = zk/hcmref
cansublays = floor(hcmref/modres)

! ... check for model vegetation types
if (lu_opt .eq. 0 ) then !VIIRS
if (vtyperef .le. 10 .or. vtyperef .eq. 12) then !VIIRS types
! ... check for valid model vegetation types
if (lu_opt .eq. 0 .or. lu_opt .eq. 1 ) then !VIIRS or MODIS
if (vtyperef .gt. 0 .and. vtyperef .le. 10 .or. vtyperef .eq. 12) then !VIIRS or MODIS types

! ... check for contiguous canopy conditions at each model grid cell
if (hcmref .gt. fch_thresh .and. ffracref .gt. frt_thresh &
Expand Down Expand Up @@ -212,7 +223,8 @@ SUBROUTINE canopy_calcs
call canopy_flameh(flameh_opt, flameh_set, dx(loc), modres, &
frpref, midflamepoint, flameh(loc))

if (flameh(loc) .gt. 0.0) then !only calculate WAF when flameh > 0
if (flameh(loc) .gt. 0.0 .and. flameh(loc) .le. hcmref) then
!only calculate WAF when flameh > 0
call canopy_waf(hcmref, lambdars, hgtref, flameh(loc), &
firetype, d_h, zo_h, canBOT(midflamepoint), &
canTOP(midflamepoint), waf(loc))
Expand All @@ -228,14 +240,16 @@ SUBROUTINE canopy_calcs

! ... user option to calculate in-canopy eddy photolysis attenuation at height z
if (ifcanphot) then
call canopy_phot(fafraczInt, &
lairef, cluref, cszref, rjcf(loc, :))
if (cszref .ge. 0.0_rk) then !only calculate if cell isn't dark
call canopy_phot(fafraczInt, &
lairef, cluref, cszref, rjcf(loc, :))
end if
end if

end if !Contiguous Canopy

else
! write(*,*) 'Warning VIIRS VTYPE ', vtyperef, ' is not supported...continue'
! write(*,*) 'Warning VIIRS/MODIS VTYPE ', vtyperef, ' is not supported...continue'
end if !Vegetation types
else
write(*,*) 'Wrong LU_OPT choice of ', lu_opt, ' in namelist...exiting'
Expand Down
9 changes: 6 additions & 3 deletions canopy_canmet_mod.F90
Original file line number Diff line number Diff line change
Expand Up @@ -15,14 +15,15 @@ MODULE canopy_canmet_mod
real(rk) :: lat !latitude of cell/point
real(rk) :: lon !longitude of cell/point
real(rk) :: fh !forest/canopy height (m)
real(rk) :: ws !wind speed at reference height above canopy (m/s)
real(rk) :: ugrd10m !u wind speed at reference height above canopy (m/s)
real(rk) :: vgrd10m !v wind speed at reference height above canopy (m/s)
real(rk) :: clu !clumping index
real(rk) :: lai !leaf area index
integer :: vtype !vegetation type
real(rk) :: ffrac !forest fraction
real(rk) :: ust !friction velocity (u*) (m/s)
real(rk) :: fricv !friction velocity (u*) (m/s)
real(rk) :: csz !cosine of solar zenith angle
real(rk) :: z0 !surface roughness length (m)
real(rk) :: sfcr !surface roughness length (m)
real(rk) :: mol !Monin-Obukhov length (m)
real(rk) :: frp !fire radiative power (MW)
real(rk) :: href !reference height above the canopy (m)
Expand All @@ -34,6 +35,8 @@ MODULE canopy_canmet_mod
real(rk) :: latref !latitude of cell/point
real(rk) :: lonref !longitude of cell/point
real(rk) :: hcmref !Input Canopy Height
real(rk) :: uref !Input above canopy/reference 10-m U wind speed
real(rk) :: vref !Input above canopy/reference 10-m V wind speed
real(rk) :: ubzref !Input above canopy/reference 10-m model wind speed
real(rk) :: cluref !Input canopy clumping index
real(rk) :: lairef !Input leaf area index
Expand Down
Loading