source: LMDZ6/trunk/libf/phylmd/ecrad.v1.5.1/setup_aerosol_optics_lmdz_m.F90 @ 5322

Last change on this file since 5322 was 4534, checked in by lguez, 18 months ago

Use netcdf95 constant

File size: 4.8 KB
Line 
1module setup_aerosol_optics_lmdz_m
2
3  implicit none
4
5contains
6
7  subroutine setup_aerosol_optics_lmdz(ao, file_name)
8
9    ! Read aerosol optical properties. Note differences with
10    ! "radiation_aerosol_optics_data::setup_aerosol_optics":
11
12    ! -- The input NetCDF file is not flat, it contains NetCDF groups.
13
14    ! -- We do not define ao%ssa_mono_phobic, ao%g_mono_phobic,
15    ! ao%lidar_ratio_mono_phobic, ao%ssa_mono_philic,
16    ! ao%g_mono_philic, ao%lidar_ratio_mono_philic. They are not in
17    ! the input NetCDF file and they are not used by ECRad.
18
19    ! -- We do not define ao%description_phobic_str and
20    ! ao%description_philic_str. We just leave the initialization
21    ! value, which is a blank.
22
23    ! -- We have to cshift the shortwave fields because the the
24    ! shortwave bands are in ascending order in the NetCDF file while
25    ! they are not in ECRad.
26
27    use radiation_aerosol_optics_data, only: aerosol_optics_type, &
28         IAerosolClassUndefined
29    use netcdf95, only: nf95_open, nf95_inq_grp_full_ncid, nf95_close, &
30         nf95_inq_dimid, nf95_inq_varid, nf95_inquire_dimension, &
31         nf95_get_var, nf95_gw_var, nf95_nowrite
32
33    type(aerosol_optics_type), intent(out):: ao
34
35    character(len=*), intent(in):: file_name
36    ! NetCDF file containing the aerosol optics data
37
38    ! Local:
39    integer ncid, grpid, dimid, varid
40
41    !-----------------------------------------------------------------------
42
43    ao%use_hydrophilic = .true.
44    ao%use_monochromatic = .true.
45    call nf95_open(file_name, nf95_nowrite, ncid)
46    call nf95_inq_grp_full_ncid(ncid, "Hydrophilic", grpid)
47    call nf95_inq_dimid(grpid, "hur", dimid)
48    call nf95_inquire_dimension(grpid, dimid, nclen = ao%nrh)
49    allocate(ao%rh_lower(ao%nrh))
50    call nf95_inq_varid(grpid, "hur_bounds", varid)
51    call nf95_get_var(grpid, varid, ao%rh_lower, count_nc = [1, ao%nrh])
52
53    ! Hydrophilic/LW_bands:
54    call nf95_inq_grp_full_ncid(ncid, "Hydrophilic/LW_bands", grpid)
55    call nf95_inq_varid(grpid, "asymmetry", varid)
56    call nf95_gw_var(grpid, varid, ao%g_lw_philic)
57    call nf95_inq_varid(grpid, "single_scat_alb", varid)
58    call nf95_gw_var(grpid, varid, ao%ssa_lw_philic)
59    call nf95_inq_varid(grpid, "mass_ext", varid)
60    call nf95_gw_var(grpid, varid, ao%mass_ext_lw_philic)
61
62    ! Hydrophilic/SW_bands:
63    call nf95_inq_grp_full_ncid(ncid, "Hydrophilic/SW_bands", grpid)
64    call nf95_inq_varid(grpid, "asymmetry", varid)
65    call nf95_gw_var(grpid, varid, ao%g_sw_philic)
66    ao%g_sw_philic = cshift(ao%g_sw_philic, 1)
67    call nf95_inq_varid(grpid, "single_scat_alb", varid)
68    call nf95_gw_var(grpid, varid, ao%ssa_sw_philic)
69    ao%ssa_sw_philic = cshift(ao%ssa_sw_philic, 1)
70    call nf95_inq_varid(grpid, "mass_ext", varid)
71    call nf95_gw_var(grpid, varid, ao%mass_ext_sw_philic)
72    ao%mass_ext_sw_philic = cshift(ao%mass_ext_sw_philic, 1)
73
74    ! Hydrophilic/Monochromatic:
75    call nf95_inq_grp_full_ncid(ncid, "Hydrophilic/Monochromatic", grpid)
76    call nf95_inq_varid(grpid, "mass_ext", varid)
77    call nf95_gw_var(grpid, varid, ao%mass_ext_mono_philic)
78
79    ! Hydrophobic/LW_bands:
80    call nf95_inq_grp_full_ncid(ncid, "Hydrophobic/LW_bands", grpid)
81    call nf95_inq_varid(grpid, "asymmetry", varid)
82    call nf95_gw_var(grpid, varid, ao%g_lw_phobic)
83    call nf95_inq_varid(grpid, "single_scat_alb", varid)
84    call nf95_gw_var(grpid, varid, ao%ssa_lw_phobic)
85    call nf95_inq_varid(grpid, "mass_ext", varid)
86    call nf95_gw_var(grpid, varid, ao%mass_ext_lw_phobic)
87
88    ! Hydrophobic/SW_bands:
89    call nf95_inq_grp_full_ncid(ncid, "Hydrophobic/SW_bands", grpid)
90    call nf95_inq_varid(grpid, "asymmetry", varid)
91    call nf95_gw_var(grpid, varid, ao%g_sw_phobic)
92    ao%g_sw_phobic = cshift(ao%g_sw_phobic, 1)
93    call nf95_inq_varid(grpid, "single_scat_alb", varid)
94    call nf95_gw_var(grpid, varid, ao%ssa_sw_phobic)
95    ao%ssa_sw_phobic = cshift(ao%ssa_sw_phobic, 1)
96    call nf95_inq_varid(grpid, "mass_ext", varid)
97    call nf95_gw_var(grpid, varid, ao%mass_ext_sw_phobic)
98    ao%mass_ext_sw_phobic = cshift(ao%mass_ext_sw_phobic, 1)
99
100    ! Hydrophobic/Monochromatic:
101    call nf95_inq_grp_full_ncid(ncid, "Hydrophobic/Monochromatic", grpid)
102    call nf95_inq_varid(grpid, "mass_ext", varid)
103    call nf95_gw_var(grpid, varid, ao%mass_ext_mono_phobic)
104
105    call nf95_close(ncid)
106
107    ! Get array sizes
108    ao%n_bands_lw = size(ao%mass_ext_lw_phobic, 1)
109    ao%n_bands_sw = size(ao%mass_ext_sw_phobic, 1)
110    ao%n_mono_wl = size(ao%mass_ext_mono_phobic, 1)
111    ao%n_type_phobic = size(ao%mass_ext_lw_phobic, 2)
112    ao%n_type_philic = size(ao%mass_ext_lw_philic, 3)
113
114    ! Allocate memory for mapping arrays
115    ao%ntype = ao%n_type_phobic + ao%n_type_philic
116    allocate(ao%iclass(ao%ntype))
117    allocate(ao%itype(ao%ntype))
118
119    ao%iclass = IAerosolClassUndefined
120    ao%itype  = 0
121
122  end subroutine setup_aerosol_optics_lmdz
123
124end module setup_aerosol_optics_lmdz_m
Note: See TracBrowser for help on using the repository browser.