Changeset 1483 for LMDZ4/branches/LMDZ4_AR5/libf/phylmd/readaerosol.F90
- Timestamp:
- Feb 9, 2011, 4:21:55 PM (14 years ago)
- File:
-
- 1 edited
Legend:
- Unmodified
- Added
- Removed
-
LMDZ4/branches/LMDZ4_AR5/libf/phylmd/readaerosol.F90
r1321 r1483 7 7 CONTAINS 8 8 9 SUBROUTINE readaerosol(name_aero, type, iyr_in, klev_src, pt_ap, pt_b, pt_out, psurf, load)9 SUBROUTINE readaerosol(name_aero, type, filename, iyr_in, klev_src, pt_ap, pt_b, pt_out, psurf, load) 10 10 11 11 !**************************************************************************************** … … 27 27 ! Input arguments 28 28 CHARACTER(len=7), INTENT(IN) :: name_aero 29 CHARACTER(len=*), INTENT(IN) :: type ! correspond to aer_type in clesphys.h 29 CHARACTER(len=*), INTENT(IN) :: type ! actuel, annuel, scenario or preind 30 CHARACTER(len=8), INTENT(IN) :: filename 30 31 INTEGER, INTENT(IN) :: iyr_in 31 32 … … 58 59 ! get_aero_fromfile returns pt_out allocated and initialized with data for 12 month 59 60 ! pt_out has dimensions (klon, klev_src, 12) 60 CALL get_aero_fromfile(name_aero, cyear, klev_src, pt_ap, pt_b, p0, pt_out, psurf, load)61 CALL get_aero_fromfile(name_aero, cyear, filename, klev_src, pt_ap, pt_b, p0, pt_out, psurf, load) 61 62 62 63 … … 67 68 ! get_aero_fromfile returns pt_out allocated and initialized with data for 12 month 68 69 ! pt_out has dimensions (klon, klev_src, 12) 69 CALL get_aero_fromfile(name_aero, cyear, klev_src, pt_ap, pt_b, p0, pt_out, psurf, load)70 CALL get_aero_fromfile(name_aero, cyear, filename, klev_src, pt_ap, pt_b, p0, pt_out, psurf, load) 70 71 71 72 ELSE IF (type == 'annuel') THEN … … 76 77 ! get_aero_fromfile returns pt_out allocated and initialized with data for nbr_tsteps month 77 78 ! pt_out has dimensions (klon, klev_src, 12) 78 CALL get_aero_fromfile(name_aero, cyear, klev_src, pt_ap, pt_b, p0, pt_out, psurf, load)79 CALL get_aero_fromfile(name_aero, cyear, filename, klev_src, pt_ap, pt_b, p0, pt_out, psurf, load) 79 80 80 81 ELSE IF (type == 'scenario') THEN … … 86 87 ! get_aero_fromfile returns pt_out allocated and initialized with data for 12 month 87 88 ! pt_out has dimensions (klon, klev_src, 12) 88 CALL get_aero_fromfile(name_aero, cyear, klev_src, pt_ap, pt_b, p0, pt_out, psurf, load)89 CALL get_aero_fromfile(name_aero, cyear, filename, klev_src, pt_ap, pt_b, p0, pt_out, psurf, load) 89 90 90 91 ELSE IF (iyr_in .GE. 2100) THEN … … 93 94 ! get_aero_fromfile returns pt_out allocated and initialized with data for 12 month 94 95 ! pt_out has dimensions (klon, klev_src, 12) 95 CALL get_aero_fromfile(name_aero, cyear, klev_src, pt_ap, pt_b, p0, pt_out, psurf, load)96 CALL get_aero_fromfile(name_aero, cyear, filename, klev_src, pt_ap, pt_b, p0, pt_out, psurf, load) 96 97 97 98 ELSE … … 113 114 ! get_aero_fromfile returns pt_out allocated and initialized with data for 12 month 114 115 ! pt_out has dimensions (klon, klev_src, 12) 115 CALL get_aero_fromfile(name_aero, cyear, klev_src, pt_ap, pt_b, p0, pt_out, psurf, load)116 CALL get_aero_fromfile(name_aero, cyear, filename, klev_src, pt_ap, pt_b, p0, pt_out, psurf, load) 116 117 117 118 ! If to read two decades: … … 125 126 ! get_aero_fromfile returns pt_2 allocated and initialized with data for 12 month 126 127 ! pt_2 has dimensions (klon, klev_src, 12) 127 CALL get_aero_fromfile(name_aero, cyear, klev_src2, pt_ap, pt_b, p0, pt_2, psurf2, load2)128 CALL get_aero_fromfile(name_aero, cyear, filename, klev_src2, pt_ap, pt_b, p0, pt_2, psurf2, load2) 128 129 ! Test for same number of vertical levels 129 130 IF (klev_src /= klev_src2) THEN … … 160 161 161 162 ELSE 162 WRITE(lunout,*)'This option is not implemented : aer_type = ', type 163 WRITE(lunout,*)'This option is not implemented : aer_type = ', type,' name_aero=',name_aero 163 164 CALL abort_gcm('readaerosol','Error : aer_type parameter not accepted',1) 164 165 END IF ! type … … 168 169 169 170 170 SUBROUTINE get_aero_fromfile(varname, cyr, klev_src, pt_ap, pt_b, p0, pt_year, psurf_out, load_out)171 SUBROUTINE get_aero_fromfile(varname, cyr, filename, klev_src, pt_ap, pt_b, p0, pt_year, psurf_out, load_out) 171 172 !**************************************************************************************** 172 173 ! Read 12 month aerosol from file and distribute to local process on physical grid. … … 200 201 CHARACTER(len=7), INTENT(IN) :: varname 201 202 CHARACTER(len=4), INTENT(IN) :: cyr 203 CHARACTER(len=8), INTENT(IN) :: filename 202 204 203 205 ! Output arguments … … 213 215 ! Local variables 214 216 CHARACTER(len=30) :: fname 215 CHARACTER(len=8) :: filename='aerosols'216 217 CHARACTER(len=30) :: cvar 217 218 INTEGER :: ncid, dimid, varid … … 242 243 ! 1) Open file 243 244 !**************************************************************************************** 244 fname = filename//cyr//'.nc' 245 ! Add suffix to filename 246 fname = trim(filename)//cyr//'.nc' 245 247 246 WRITE(lunout,*) 'reading ', TRIM(fname)248 WRITE(lunout,*) 'reading variable ',TRIM(varname),' in file ', TRIM(fname) 247 249 CALL check_err( nf90_open(TRIM(fname), NF90_NOWRITE, ncid) ) 248 250 … … 283 285 CALL abort_gcm('get_aero_fromfile', 'latitudes do not correspond between file and model',1) 284 286 END IF 285 286 ! 1.5) Check number of month in file opened287 !288 !**************************************************************************************************289 ierr = nf90_inq_dimid(ncid, 'TIME',dimid)290 CALL check_err( nf90_inquire_dimension(ncid, dimid, len = nbr_tsteps) )291 ! IF (nbr_tsteps /= 12 .AND. nbr_tsteps /= 14) THEN292 IF (nbr_tsteps /= 12 ) THEN293 CALL abort_gcm('get_aero_fromfile', 'not the right number of months in aerosol file read (should be 12 for the moment)',1)294 ENDIF295 287 296 288 … … 335 327 336 328 IF (new_file) THEN 329 330 ! **) Check number of month in file opened 331 !************************************************************************************************** 332 ierr = nf90_inq_dimid(ncid, 'TIME',dimid) 333 CALL check_err( nf90_inquire_dimension(ncid, dimid, len = nbr_tsteps) ) 334 ! IF (nbr_tsteps /= 12 .AND. nbr_tsteps /= 14) THEN 335 IF (nbr_tsteps /= 12 ) THEN 336 CALL abort_gcm('get_aero_fromfile', & 337 'not the right number of months in aerosol file read (should be 12 for the moment)',1) 338 ENDIF 337 339 338 340 ! ++) Read the aerosol concentration month by month and concatenate to total variable varyear
Note: See TracChangeset
for help on using the changeset viewer.