[816] | 1 | program tmc |
---|
| 2 | |
---|
| 3 | ! SL 01/2010: |
---|
| 4 | ! This program reads 4D (lon-lat-alt-time) fields recast in log P coordinates |
---|
| 5 | ! |
---|
| 6 | ! it computes angular momentum transport from high-frequency outputs: |
---|
| 7 | ! |
---|
| 8 | ! totvang -- 2D -- Meridional transport of angular momentum, total (m3 s-2) |
---|
| 9 | ! totwang -- 2D -- Vertical transport of angular momentum, total (m3 s-2) |
---|
| 10 | ! mmcvang -- 2D -- Meridional transport of angular momentum, by MMC (m3 s-2) |
---|
| 11 | ! mmcwang -- 2D -- Vertical transport of angular momentum, by MMC (m3 s-2) |
---|
| 12 | ! trsvang -- 2D -- Meridional transport of angular momentum, transients (m3 s-2) |
---|
| 13 | ! trswang -- 2D -- Vertical transport of angular momentum, transients (m3 s-2) |
---|
| 14 | ! stnvang -- 2D -- Meridional transport of angular momentum, stationaries (m3 s-2) |
---|
| 15 | ! stnwang -- 2D -- Vertical transport of angular momentum, stationaries (m3 s-2) |
---|
| 16 | ! dmass -- 2D -- Mass in each cell (dmassmeanbar) |
---|
| 17 | ! |
---|
| 18 | ! Minimal requirements and dependencies: |
---|
| 19 | ! The dataset must include the following data: |
---|
| 20 | ! - pressure vertical coordinate |
---|
| 21 | ! - surface pressure |
---|
| 22 | ! - zonal, meridional and vertical (Pa/s) winds |
---|
| 23 | ! - altitude above areoid |
---|
| 24 | ! |
---|
| 25 | ! Convention: qbar <=> zonal average / qstar = q - qbar |
---|
| 26 | ! qmean <=> temporal average / qprim = q - qmean |
---|
| 27 | ! |
---|
| 28 | ! Therefore: ((qv)mean)bar (total) |
---|
| 29 | ! = qmeanbar * vmeanbar (mmc) |
---|
| 30 | ! + (qstarmean * vstarmean)bar (stn) |
---|
| 31 | ! + (qprim * vprim)meanbar (trs) |
---|
| 32 | ! |
---|
| 33 | !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! |
---|
| 34 | |
---|
| 35 | implicit none |
---|
| 36 | |
---|
| 37 | include "netcdf.inc" ! NetCDF definitions |
---|
| 38 | |
---|
| 39 | character (len=128) :: infile ! input file name (name_P.nc) |
---|
| 40 | character (len=128) :: outfile ! output file name |
---|
| 41 | |
---|
| 42 | character (len=64) :: text ! to store some text |
---|
| 43 | integer infid ! NetCDF input file ID |
---|
| 44 | integer outfid ! NetCDF output file ID |
---|
| 45 | integer lon_dimid,lat_dimid,alt_dimid,time_dimid ! NetCDF dimension IDs |
---|
| 46 | integer lon_varid,lat_varid,alt_varid,time_varid |
---|
| 47 | integer,dimension(2) :: datashape2d ! shape of 3D datasets |
---|
| 48 | integer,dimension(3) :: datashape3d ! shape of 3D datasets |
---|
| 49 | integer,dimension(4) :: datashape4d ! shape of 4D datasets |
---|
| 50 | |
---|
| 51 | real :: miss_val=9.99e+33 ! special "missing value" to specify missing data |
---|
| 52 | real,parameter :: miss_val_def=9.99e+33 ! default value for "missing value" |
---|
| 53 | real :: pi |
---|
| 54 | real,dimension(:),allocatable :: lon ! longitude |
---|
| 55 | integer lonlength ! # of grid points along longitude |
---|
| 56 | real,dimension(:),allocatable :: lat ! latitude |
---|
| 57 | real,dimension(:),allocatable :: latrad ! latitude in radian |
---|
| 58 | integer latlength ! # of grid points along latitude |
---|
| 59 | real,dimension(:),allocatable :: plev ! Pressure levels (Pa) |
---|
| 60 | integer altlength ! # of grid point along altitude (of input datasets) |
---|
| 61 | real,dimension(:),allocatable :: time ! time |
---|
| 62 | integer timelength ! # of points along time |
---|
| 63 | real,dimension(:,:,:),allocatable :: ps ! surface pressure |
---|
| 64 | real,dimension(:,:,:,:),allocatable :: vitu ! zonal wind (in m/s) |
---|
| 65 | real,dimension(:,:,:,:),allocatable :: vitv ! meridional wind (in m/s) |
---|
| 66 | real,dimension(:,:,:,:),allocatable :: vitw ! vertical wind (in Pa/s, then converted in m/s) |
---|
| 67 | real,dimension(:,:,:,:),allocatable :: za ! above areoid levels (m) |
---|
| 68 | |
---|
| 69 | !!! output variables |
---|
| 70 | real,dimension(:,:),allocatable :: totvang ! merid transport of ang momentum |
---|
| 71 | real,dimension(:,:),allocatable :: totwang ! verti transport of ang momentum |
---|
| 72 | real,dimension(:,:),allocatable :: mmcvang ! merid transport of ang momentum |
---|
| 73 | real,dimension(:,:),allocatable :: mmcwang ! verti transport of ang momentum |
---|
| 74 | real,dimension(:,:),allocatable :: trsvang ! merid transport of ang momentum |
---|
| 75 | real,dimension(:,:),allocatable :: trswang ! verti transport of ang momentum |
---|
| 76 | real,dimension(:,:),allocatable :: stnvang ! merid transport of ang momentum |
---|
| 77 | real,dimension(:,:),allocatable :: stnwang ! verti transport of ang momentum |
---|
| 78 | real,dimension(:,:),allocatable :: dmassmeanbar |
---|
| 79 | |
---|
| 80 | ! local variables |
---|
| 81 | real :: deltalat,deltalon ! lat and lon intervals in radians |
---|
| 82 | real,dimension(:,:,:,:),allocatable :: rayon ! distance to center (m) |
---|
| 83 | real,dimension(:,:,:,:),allocatable :: grav ! gravity field (m s-2) |
---|
| 84 | real,dimension(:,:,:,:),allocatable :: dmass ! mass in cell (kg) |
---|
| 85 | real,dimension(:,:,:),allocatable :: dmassmean |
---|
| 86 | |
---|
| 87 | real,dimension(:,:,:,:),allocatable :: osam ! planetary rotation specific ang (m2/s) |
---|
| 88 | real,dimension(:,:,:,:),allocatable :: rsam ! zonal wind specific ang (m2/s) |
---|
| 89 | real,dimension(:,:,:,:),allocatable :: tsam ! total specific ang (m2/s) |
---|
| 90 | |
---|
| 91 | real,dimension(:,:,:,:),allocatable :: vang ! v * specific ang (m3/s) |
---|
| 92 | real,dimension(:,:,:,:),allocatable :: wang ! w * specific ang (m3/s) |
---|
| 93 | real,dimension(:,:,:,:),allocatable :: vstar |
---|
| 94 | real,dimension(:,:,:,:),allocatable :: wstar |
---|
| 95 | real,dimension(:,:,:,:),allocatable :: angstar |
---|
| 96 | real,dimension(:,:,:,:),allocatable :: vprim |
---|
| 97 | real,dimension(:,:,:,:),allocatable :: wprim |
---|
| 98 | real,dimension(:,:,:,:),allocatable :: angprim |
---|
| 99 | real,dimension(:,:,:,:),allocatable :: angprimvprim |
---|
| 100 | real,dimension(:,:,:,:),allocatable :: angprimwprim |
---|
| 101 | |
---|
| 102 | ! lon,lat,alt |
---|
| 103 | real,dimension(:,:,:),allocatable :: vangmean |
---|
| 104 | real,dimension(:,:,:),allocatable :: wangmean |
---|
| 105 | real,dimension(:,:,:),allocatable :: vmean |
---|
| 106 | real,dimension(:,:,:),allocatable :: wmean |
---|
| 107 | real,dimension(:,:,:),allocatable :: angmean |
---|
| 108 | real,dimension(:,:,:),allocatable :: angprimvprimmean |
---|
| 109 | real,dimension(:,:,:),allocatable :: angprimwprimmean |
---|
| 110 | real,dimension(:,:,:),allocatable :: vstarmean |
---|
| 111 | real,dimension(:,:,:),allocatable :: wstarmean |
---|
| 112 | real,dimension(:,:,:),allocatable :: angstarmean |
---|
| 113 | real,dimension(:,:,:),allocatable :: angstarmeanvstarmean |
---|
| 114 | real,dimension(:,:,:),allocatable :: angstarmeanwstarmean |
---|
| 115 | real,dimension(:,:,:),allocatable :: v3d ! intermediate for vbar |
---|
| 116 | real,dimension(:,:,:),allocatable :: w3d ! intermediate for wbar |
---|
| 117 | real,dimension(:,:,:),allocatable :: ang3d ! intermediate for angbar |
---|
| 118 | |
---|
| 119 | ! lat,alt,time |
---|
| 120 | real,dimension(:,:,:),allocatable :: vbar |
---|
| 121 | real,dimension(:,:,:),allocatable :: wbar |
---|
| 122 | real,dimension(:,:,:),allocatable :: angbar |
---|
| 123 | |
---|
| 124 | !lat,alt |
---|
| 125 | real,dimension(:,:),allocatable :: vmeanbar |
---|
| 126 | real,dimension(:,:),allocatable :: wmeanbar |
---|
| 127 | real,dimension(:,:),allocatable :: angmeanbar |
---|
| 128 | real,dimension(:,:),allocatable :: vbar2d ! intermediate for vbar |
---|
| 129 | real,dimension(:,:),allocatable :: wbar2d ! intermediate for wbar |
---|
| 130 | real,dimension(:,:),allocatable :: angbar2d ! intermediate for qbar |
---|
| 131 | |
---|
| 132 | integer ierr,ierr1,ierr2 ! NetCDF routines return codes |
---|
| 133 | integer i,j,ilon,ilat,ilev,itim ! for loops |
---|
| 134 | logical :: lmdflag |
---|
| 135 | |
---|
| 136 | include "planet.h" |
---|
| 137 | |
---|
| 138 | !=============================================================================== |
---|
| 139 | ! 1. Input parameters |
---|
| 140 | !=============================================================================== |
---|
| 141 | |
---|
| 142 | pi = 2.*asin(1.) |
---|
| 143 | |
---|
| 144 | write(*,*) "" |
---|
| 145 | write(*,*) "You are working on the atmosphere of ",planet |
---|
| 146 | |
---|
| 147 | !=============================================================================== |
---|
| 148 | ! 1.1 Input file |
---|
| 149 | !=============================================================================== |
---|
| 150 | |
---|
| 151 | write(*,*) "" |
---|
| 152 | write(*,*) "Program valid for files with pressure axis (*_P.nc)" |
---|
| 153 | write(*,*) "Enter input file name:" |
---|
| 154 | |
---|
| 155 | read(*,'(a128)') infile |
---|
| 156 | write(*,*) "" |
---|
| 157 | |
---|
| 158 | ! open input file |
---|
| 159 | |
---|
| 160 | ierr = NF_OPEN(infile,NF_NOWRITE,infid) |
---|
| 161 | if (ierr.ne.NF_NOERR) then |
---|
| 162 | write(*,*) 'ERROR: Pb opening file ',trim(infile) |
---|
| 163 | stop "" |
---|
| 164 | endif |
---|
| 165 | |
---|
| 166 | !=============================================================================== |
---|
| 167 | ! 1.2 Get grids in lon,lat,alt(pressure),time |
---|
| 168 | !=============================================================================== |
---|
| 169 | |
---|
| 170 | call get_iddim(infid,lat_varid,latlength,lon_varid,lonlength,& |
---|
| 171 | alt_varid,altlength,time_varid,timelength,lmdflag ) |
---|
| 172 | |
---|
| 173 | allocate(lon(lonlength)) |
---|
| 174 | ierr=NF_GET_VAR_REAL(infid,lon_varid,lon) |
---|
| 175 | if (ierr.ne.NF_NOERR) stop "Error: Failed reading longitude" |
---|
| 176 | |
---|
| 177 | allocate(lat(latlength)) |
---|
| 178 | ierr=NF_GET_VAR_REAL(infid,lat_varid,lat) |
---|
| 179 | if (ierr.ne.NF_NOERR) stop "Error: Failed reading lat" |
---|
| 180 | |
---|
| 181 | allocate(latrad(latlength)) |
---|
| 182 | latrad = lat*pi/180. |
---|
| 183 | |
---|
| 184 | ! Lat, lon pressure intervals |
---|
| 185 | deltalat = abs(latrad(2)-latrad(1)) |
---|
| 186 | deltalon = abs(lon(2)-lon(1))*pi/180. |
---|
| 187 | |
---|
| 188 | allocate(plev(altlength)) |
---|
| 189 | ierr=NF_GET_VAR_REAL(infid,alt_varid,plev) |
---|
| 190 | if (ierr.ne.NF_NOERR) stop "Error: Failed reading altitude (ie pressure levels)" |
---|
| 191 | |
---|
| 192 | allocate(time(timelength)) |
---|
| 193 | ierr=NF_GET_VAR_REAL(infid,time_varid,time) |
---|
| 194 | if (ierr.ne.NF_NOERR) stop "Error: Failed reading time" |
---|
| 195 | |
---|
| 196 | !=============================================================================== |
---|
| 197 | ! 1.3 Get output file name |
---|
| 198 | !=============================================================================== |
---|
| 199 | write(*,*) "" |
---|
| 200 | !write(*,*) "Enter output file name" |
---|
| 201 | !read(*,*) outfile |
---|
| 202 | outfile=infile(1:len_trim(infile)-3)//"_TMC.nc" |
---|
| 203 | write(*,*) "Output file name is: "//trim(outfile) |
---|
| 204 | |
---|
| 205 | |
---|
| 206 | |
---|
| 207 | !=============================================================================== |
---|
| 208 | ! 2.1 Store needed fields |
---|
| 209 | !=============================================================================== |
---|
| 210 | |
---|
| 211 | !=============================================================================== |
---|
| 212 | ! 2.1.1 Surface pressure |
---|
| 213 | !=============================================================================== |
---|
| 214 | allocate(ps(lonlength,latlength,timelength)) |
---|
| 215 | |
---|
| 216 | text="ps" |
---|
| 217 | call get_var3d(infid,lonlength,latlength,timelength,text,ps,ierr1,ierr2) |
---|
| 218 | if (ierr1.ne.NF_NOERR) then |
---|
| 219 | write(*,*) " looking for psol instead... " |
---|
| 220 | text="psol" |
---|
| 221 | call get_var3d(infid,lonlength,latlength,timelength,text,ps,ierr1,ierr2) |
---|
| 222 | if (ierr1.ne.NF_NOERR) stop "Error: Failed to get psol ID" |
---|
| 223 | endif |
---|
| 224 | if (ierr2.ne.NF_NOERR) stop "Error: Failed reading surface pressure" |
---|
| 225 | |
---|
| 226 | !=============================================================================== |
---|
| 227 | ! 2.1.3 Winds |
---|
| 228 | !=============================================================================== |
---|
| 229 | allocate(vitu(lonlength,latlength,altlength,timelength)) |
---|
| 230 | allocate(vitv(lonlength,latlength,altlength,timelength)) |
---|
| 231 | allocate(vitw(lonlength,latlength,altlength,timelength)) |
---|
| 232 | |
---|
| 233 | ! zonal wind vitu (in m/s) |
---|
| 234 | text="vitu" |
---|
| 235 | call get_var4d(infid,lonlength,latlength,altlength,timelength,text,vitu,ierr1,ierr2) |
---|
| 236 | if (ierr1.ne.NF_NOERR) stop "Error: Failed to get vitu ID" |
---|
| 237 | if (ierr2.ne.NF_NOERR) stop "Error: Failed reading zonal wind" |
---|
| 238 | |
---|
| 239 | ! meridional wind vitv (in m/s) |
---|
| 240 | text="vitv" |
---|
| 241 | call get_var4d(infid,lonlength,latlength,altlength,timelength,text,vitv,ierr1,ierr2) |
---|
| 242 | if (ierr1.ne.NF_NOERR) stop "Error: Failed to get vitv ID" |
---|
| 243 | if (ierr2.ne.NF_NOERR) stop "Error: Failed reading meridional wind" |
---|
| 244 | |
---|
| 245 | ! vertical wind vitw (in Pa/s) |
---|
| 246 | text="vitw" |
---|
| 247 | call get_var4d(infid,lonlength,latlength,altlength,timelength,text,vitw,ierr1,ierr2) |
---|
| 248 | if (ierr1.ne.NF_NOERR) stop "Error: Failed to get vitw ID" |
---|
| 249 | if (ierr2.ne.NF_NOERR) stop "Error: Failed reading vertical wind" |
---|
| 250 | |
---|
| 251 | !=============================================================================== |
---|
| 252 | ! 2.1.4 Altitude above areoide |
---|
| 253 | !=============================================================================== |
---|
| 254 | ! Only needed if g(z) on Titan... |
---|
| 255 | |
---|
| 256 | ! allocate(za(lonlength,latlength,altlength,timelength)) |
---|
| 257 | |
---|
| 258 | ! text="zareoid" |
---|
| 259 | ! call get_var4d(infid,lonlength,latlength,altlength,timelength,text,za,ierr1,ierr2) |
---|
| 260 | ! if (ierr1.ne.NF_NOERR) stop "Error: Failed to get za ID" |
---|
| 261 | ! if (ierr2.ne.NF_NOERR) stop "Error: Failed reading zareoid" |
---|
| 262 | |
---|
| 263 | !=============================================================================== |
---|
| 264 | ! 2.2 Computations |
---|
| 265 | !=============================================================================== |
---|
| 266 | |
---|
| 267 | !=============================================================================== |
---|
| 268 | ! 2.2.2 Mass in cells |
---|
| 269 | !=============================================================================== |
---|
| 270 | allocate(rayon(lonlength,latlength,altlength,timelength)) |
---|
| 271 | allocate( grav(lonlength,latlength,altlength,timelength)) |
---|
| 272 | allocate(dmass(lonlength,latlength,altlength,timelength)) |
---|
| 273 | allocate(dmassmean(lonlength,latlength,altlength)) |
---|
| 274 | allocate(dmassmeanbar(latlength,altlength)) |
---|
| 275 | |
---|
| 276 | do itim=1,timelength |
---|
| 277 | do ilon=1,lonlength |
---|
| 278 | do ilat=1,latlength |
---|
| 279 | do ilev=1,altlength |
---|
| 280 | ! Need to be consistent with GCM computations |
---|
| 281 | ! if (za(ilon,ilat,ilev,itim).ne.miss_val) then |
---|
| 282 | rayon(ilon,ilat,ilev,itim) = a0 |
---|
| 283 | ! rayon(ilon,ilat,ilev,itim) = za(ilon,ilat,ilev,itim) + a0 |
---|
| 284 | grav(ilon,ilat,ilev,itim) = g0*a0*a0 & |
---|
| 285 | /(rayon(ilon,ilat,ilev,itim)*rayon(ilon,ilat,ilev,itim)) |
---|
| 286 | ! else |
---|
| 287 | ! rayon(ilon,ilat,ilev,itim) = miss_val |
---|
| 288 | ! grav(ilon,ilat,ilev,itim) = miss_val |
---|
| 289 | ! endif |
---|
| 290 | enddo |
---|
| 291 | enddo |
---|
| 292 | enddo |
---|
| 293 | enddo ! timelength |
---|
| 294 | |
---|
| 295 | call cellmass(infid,latlength,lonlength,altlength,timelength, & |
---|
| 296 | lmdflag,deltalon,deltalat,latrad,plev,ps,grav,rayon, & |
---|
| 297 | dmass ) |
---|
| 298 | |
---|
| 299 | call moytim(lonlength,latlength,altlength,timelength,miss_val,dmass,dmassmean) |
---|
| 300 | call moyzon2(lonlength,latlength,altlength,miss_val,dmassmean,dmassmeanbar) |
---|
| 301 | |
---|
| 302 | print*,"OK dmass" |
---|
| 303 | |
---|
| 304 | !=============================================================================== |
---|
| 305 | ! 2.2.3 Specific angular momentum |
---|
| 306 | !=============================================================================== |
---|
| 307 | allocate(osam(lonlength,latlength,altlength,timelength)) |
---|
| 308 | allocate(rsam(lonlength,latlength,altlength,timelength)) |
---|
| 309 | allocate(tsam(lonlength,latlength,altlength,timelength)) |
---|
| 310 | |
---|
| 311 | do itim=1,timelength |
---|
| 312 | do ilon=1,lonlength |
---|
| 313 | do ilat=1,latlength |
---|
| 314 | do ilev=1,altlength |
---|
| 315 | if (rayon(ilon,ilat,ilev,itim).ne.miss_val) then |
---|
| 316 | osam(ilon,ilat,ilev,itim) = & |
---|
| 317 | rayon(ilon,ilat,ilev,itim)*rayon(ilon,ilat,ilev,itim) & |
---|
| 318 | * cos(latrad(ilat))*cos(latrad(ilat)) & |
---|
| 319 | * omega |
---|
| 320 | rsam(ilon,ilat,ilev,itim) = vitu(ilon,ilat,ilev,itim) & |
---|
| 321 | * rayon(ilon,ilat,ilev,itim)*cos(latrad(ilat)) |
---|
| 322 | tsam(ilon,ilat,ilev,itim) = osam(ilon,ilat,ilev,itim)& |
---|
| 323 | + rsam(ilon,ilat,ilev,itim) |
---|
| 324 | else |
---|
| 325 | osam(ilon,ilat,ilev,itim) = miss_val |
---|
| 326 | rsam(ilon,ilat,ilev,itim) = miss_val |
---|
| 327 | tsam(ilon,ilat,ilev,itim) = miss_val |
---|
| 328 | endif |
---|
| 329 | enddo |
---|
| 330 | enddo |
---|
| 331 | enddo |
---|
| 332 | enddo ! timelength |
---|
| 333 | |
---|
| 334 | print*,"debut tprt ang" |
---|
| 335 | |
---|
| 336 | !=============================================================================== |
---|
| 337 | ! 2.2.4 Angular momentum transport |
---|
| 338 | !=============================================================================== |
---|
| 339 | ! allocations |
---|
| 340 | !------------- |
---|
| 341 | allocate(vang(lonlength,latlength,altlength,timelength)) |
---|
| 342 | allocate(wang(lonlength,latlength,altlength,timelength)) |
---|
| 343 | allocate( vstar(lonlength,latlength,altlength,timelength)) |
---|
| 344 | allocate( wstar(lonlength,latlength,altlength,timelength)) |
---|
| 345 | allocate(angstar(lonlength,latlength,altlength,timelength)) |
---|
| 346 | allocate( vprim(lonlength,latlength,altlength,timelength)) |
---|
| 347 | allocate( wprim(lonlength,latlength,altlength,timelength)) |
---|
| 348 | allocate(angprim(lonlength,latlength,altlength,timelength)) |
---|
| 349 | allocate(angprimvprim(lonlength,latlength,altlength,timelength)) |
---|
| 350 | allocate(angprimwprim(lonlength,latlength,altlength,timelength)) |
---|
| 351 | |
---|
| 352 | ! lon,lat,alt |
---|
| 353 | allocate(vangmean(lonlength,latlength,altlength)) |
---|
| 354 | allocate(wangmean(lonlength,latlength,altlength)) |
---|
| 355 | allocate( vmean(lonlength,latlength,altlength)) |
---|
| 356 | allocate( wmean(lonlength,latlength,altlength)) |
---|
| 357 | allocate(angmean(lonlength,latlength,altlength)) |
---|
| 358 | allocate(angprimvprimmean(lonlength,latlength,altlength)) |
---|
| 359 | allocate(angprimwprimmean(lonlength,latlength,altlength)) |
---|
| 360 | allocate( vstarmean(lonlength,latlength,altlength)) |
---|
| 361 | allocate( wstarmean(lonlength,latlength,altlength)) |
---|
| 362 | allocate(angstarmean(lonlength,latlength,altlength)) |
---|
| 363 | allocate(angstarmeanvstarmean(lonlength,latlength,altlength)) |
---|
| 364 | allocate(angstarmeanwstarmean(lonlength,latlength,altlength)) |
---|
| 365 | allocate( v3d(lonlength,latlength,altlength)) |
---|
| 366 | allocate( w3d(lonlength,latlength,altlength)) |
---|
| 367 | allocate(ang3d(lonlength,latlength,altlength)) |
---|
| 368 | |
---|
| 369 | ! lat,alt,time |
---|
| 370 | allocate( vbar(latlength,altlength,timelength)) |
---|
| 371 | allocate( wbar(latlength,altlength,timelength)) |
---|
| 372 | allocate(angbar(latlength,altlength,timelength)) |
---|
| 373 | |
---|
| 374 | !lat,alt |
---|
| 375 | allocate( vmeanbar(latlength,altlength)) |
---|
| 376 | allocate( wmeanbar(latlength,altlength)) |
---|
| 377 | allocate(angmeanbar(latlength,altlength)) |
---|
| 378 | allocate( vbar2d(latlength,altlength)) |
---|
| 379 | allocate( wbar2d(latlength,altlength)) |
---|
| 380 | allocate(angbar2d(latlength,altlength)) |
---|
| 381 | |
---|
| 382 | allocate(totvang(latlength,altlength)) |
---|
| 383 | allocate(totwang(latlength,altlength)) |
---|
| 384 | allocate(mmcvang(latlength,altlength)) |
---|
| 385 | allocate(mmcwang(latlength,altlength)) |
---|
| 386 | allocate(trsvang(latlength,altlength)) |
---|
| 387 | allocate(trswang(latlength,altlength)) |
---|
| 388 | allocate(stnvang(latlength,altlength)) |
---|
| 389 | allocate(stnwang(latlength,altlength)) |
---|
| 390 | |
---|
| 391 | ! intermediates |
---|
| 392 | !----------------- |
---|
| 393 | |
---|
| 394 | do itim=1,timelength |
---|
| 395 | v3d(:,:,:) = vitv(:,:,:,itim) |
---|
| 396 | w3d(:,:,:) = vitw(:,:,:,itim) |
---|
| 397 | ang3d(:,:,:) = tsam(:,:,:,itim) |
---|
| 398 | call moyzon2(lonlength,latlength,altlength,miss_val, v3d, vbar2d) |
---|
| 399 | call moyzon2(lonlength,latlength,altlength,miss_val, w3d, wbar2d) |
---|
| 400 | call moyzon2(lonlength,latlength,altlength,miss_val,ang3d,angbar2d) |
---|
| 401 | vbar(:,:,itim) = vbar2d(:,:) |
---|
| 402 | wbar(:,:,itim) = wbar2d(:,:) |
---|
| 403 | angbar(:,:,itim) = angbar2d(:,:) |
---|
| 404 | enddo ! timelength |
---|
| 405 | |
---|
| 406 | do ilon=1,lonlength |
---|
| 407 | do ilat=1,latlength |
---|
| 408 | do ilev=1,altlength |
---|
| 409 | do itim=1,timelength |
---|
| 410 | if ((vitv(ilon,ilat,ilev,itim).ne.miss_val).and. & |
---|
| 411 | (vbar(ilat,ilev,itim) .ne.miss_val)) then |
---|
| 412 | vstar(ilon,ilat,ilev,itim) = vitv(ilon,ilat,ilev,itim)- vbar(ilat,ilev,itim) |
---|
| 413 | else |
---|
| 414 | vstar(ilon,ilat,ilev,itim) = miss_val |
---|
| 415 | endif |
---|
| 416 | if ((vitw(ilon,ilat,ilev,itim).ne.miss_val).and. & |
---|
| 417 | (wbar(ilat,ilev,itim) .ne.miss_val)) then |
---|
| 418 | wstar(ilon,ilat,ilev,itim) = vitw(ilon,ilat,ilev,itim)- wbar(ilat,ilev,itim) |
---|
| 419 | else |
---|
| 420 | wstar(ilon,ilat,ilev,itim) = miss_val |
---|
| 421 | endif |
---|
| 422 | if (( tsam(ilon,ilat,ilev,itim).ne.miss_val).and. & |
---|
| 423 | (angbar(ilat,ilev,itim) .ne.miss_val)) then |
---|
| 424 | angstar(ilon,ilat,ilev,itim) = tsam(ilon,ilat,ilev,itim)-angbar(ilat,ilev,itim) |
---|
| 425 | else |
---|
| 426 | angstar(ilon,ilat,ilev,itim) = miss_val |
---|
| 427 | endif |
---|
| 428 | enddo |
---|
| 429 | enddo |
---|
| 430 | enddo |
---|
| 431 | enddo ! lonlength |
---|
| 432 | call moytim(lonlength,latlength,altlength,timelength,miss_val, vstar, vstarmean) |
---|
| 433 | call moytim(lonlength,latlength,altlength,timelength,miss_val, wstar, wstarmean) |
---|
| 434 | call moytim(lonlength,latlength,altlength,timelength,miss_val,angstar,angstarmean) |
---|
| 435 | do ilon=1,lonlength |
---|
| 436 | do ilat=1,latlength |
---|
| 437 | do ilev=1,altlength |
---|
| 438 | if ((angstarmean(ilon,ilat,ilev).ne.miss_val).and. & |
---|
| 439 | ( vstarmean(ilon,ilat,ilev).ne.miss_val)) then |
---|
| 440 | angstarmeanvstarmean(ilon,ilat,ilev) = angstarmean(ilon,ilat,ilev)*vstarmean(ilon,ilat,ilev) |
---|
| 441 | else |
---|
| 442 | angstarmeanvstarmean(ilon,ilat,ilev) = miss_val |
---|
| 443 | endif |
---|
| 444 | if ((angstarmean(ilon,ilat,ilev).ne.miss_val).and. & |
---|
| 445 | ( wstarmean(ilon,ilat,ilev).ne.miss_val)) then |
---|
| 446 | angstarmeanwstarmean(ilon,ilat,ilev) = angstarmean(ilon,ilat,ilev)*wstarmean(ilon,ilat,ilev) |
---|
| 447 | else |
---|
| 448 | angstarmeanwstarmean(ilon,ilat,ilev) = miss_val |
---|
| 449 | endif |
---|
| 450 | enddo |
---|
| 451 | enddo |
---|
| 452 | enddo ! lonlength |
---|
| 453 | |
---|
| 454 | call moytim(lonlength,latlength,altlength,timelength,miss_val,vitv, vmean) |
---|
| 455 | call moytim(lonlength,latlength,altlength,timelength,miss_val,vitw, wmean) |
---|
| 456 | call moytim(lonlength,latlength,altlength,timelength,miss_val,tsam,angmean) |
---|
| 457 | call moyzon2(lonlength,latlength,altlength,miss_val, vmean, vmeanbar) |
---|
| 458 | call moyzon2(lonlength,latlength,altlength,miss_val, wmean, wmeanbar) |
---|
| 459 | call moyzon2(lonlength,latlength,altlength,miss_val,angmean,angmeanbar) |
---|
| 460 | |
---|
| 461 | do ilon=1,lonlength |
---|
| 462 | do ilat=1,latlength |
---|
| 463 | do ilev=1,altlength |
---|
| 464 | do itim=1,timelength |
---|
| 465 | if ((vitv(ilon,ilat,ilev,itim).ne.miss_val).and. & |
---|
| 466 | (tsam(ilon,ilat,ilev,itim).ne.miss_val)) then |
---|
| 467 | vang(ilon,ilat,ilev,itim) = vitv(ilon,ilat,ilev,itim)*tsam(ilon,ilat,ilev,itim) |
---|
| 468 | else |
---|
| 469 | vang(ilon,ilat,ilev,itim) = miss_val |
---|
| 470 | endif |
---|
| 471 | if ((vitw(ilon,ilat,ilev,itim).ne.miss_val).and. & |
---|
| 472 | (tsam(ilon,ilat,ilev,itim).ne.miss_val)) then |
---|
| 473 | wang(ilon,ilat,ilev,itim) = vitw(ilon,ilat,ilev,itim)*tsam(ilon,ilat,ilev,itim) |
---|
| 474 | else |
---|
| 475 | wang(ilon,ilat,ilev,itim) = miss_val |
---|
| 476 | endif |
---|
| 477 | enddo |
---|
| 478 | enddo |
---|
| 479 | enddo |
---|
| 480 | enddo ! lonlength |
---|
| 481 | call moytim(lonlength,latlength,altlength,timelength,miss_val,vang,vangmean) |
---|
| 482 | call moytim(lonlength,latlength,altlength,timelength,miss_val,wang,wangmean) |
---|
| 483 | |
---|
| 484 | do ilon=1,lonlength |
---|
| 485 | do ilat=1,latlength |
---|
| 486 | do ilev=1,altlength |
---|
| 487 | do itim=1,timelength |
---|
| 488 | if ((vitv(ilon,ilat,ilev,itim).ne.miss_val).and. & |
---|
| 489 | (vmean(ilon,ilat,ilev) .ne.miss_val)) then |
---|
| 490 | vprim(ilon,ilat,ilev,itim) = vitv(ilon,ilat,ilev,itim)- vmean(ilon,ilat,ilev) |
---|
| 491 | else |
---|
| 492 | vprim(ilon,ilat,ilev,itim) = miss_val |
---|
| 493 | endif |
---|
| 494 | if ((vitw(ilon,ilat,ilev,itim).ne.miss_val).and. & |
---|
| 495 | (wmean(ilon,ilat,ilev) .ne.miss_val)) then |
---|
| 496 | wprim(ilon,ilat,ilev,itim) = vitw(ilon,ilat,ilev,itim)- wmean(ilon,ilat,ilev) |
---|
| 497 | else |
---|
| 498 | wprim(ilon,ilat,ilev,itim) = miss_val |
---|
| 499 | endif |
---|
| 500 | if ((tsam(ilon,ilat,ilev,itim).ne.miss_val).and. & |
---|
| 501 | (angmean(ilon,ilat,ilev) .ne.miss_val)) then |
---|
| 502 | angprim(ilon,ilat,ilev,itim) = tsam(ilon,ilat,ilev,itim)-angmean(ilon,ilat,ilev) |
---|
| 503 | else |
---|
| 504 | angprim(ilon,ilat,ilev,itim) = miss_val |
---|
| 505 | endif |
---|
| 506 | if ((angprim(ilon,ilat,ilev,itim).ne.miss_val).and. & |
---|
| 507 | ( vprim(ilon,ilat,ilev,itim).ne.miss_val)) then |
---|
| 508 | angprimvprim(ilon,ilat,ilev,itim) = angprim(ilon,ilat,ilev,itim)*vprim(ilon,ilat,ilev,itim) |
---|
| 509 | else |
---|
| 510 | angprimvprim(ilon,ilat,ilev,itim) = miss_val |
---|
| 511 | endif |
---|
| 512 | if ((angprim(ilon,ilat,ilev,itim).ne.miss_val).and. & |
---|
| 513 | ( wprim(ilon,ilat,ilev,itim).ne.miss_val)) then |
---|
| 514 | angprimwprim(ilon,ilat,ilev,itim) = angprim(ilon,ilat,ilev,itim)*wprim(ilon,ilat,ilev,itim) |
---|
| 515 | else |
---|
| 516 | angprimwprim(ilon,ilat,ilev,itim) = miss_val |
---|
| 517 | endif |
---|
| 518 | enddo |
---|
| 519 | enddo |
---|
| 520 | enddo |
---|
| 521 | enddo ! lonlength |
---|
| 522 | call moytim(lonlength,latlength,altlength,timelength,miss_val,& |
---|
| 523 | angprimvprim,angprimvprimmean) |
---|
| 524 | call moytim(lonlength,latlength,altlength,timelength,miss_val,& |
---|
| 525 | angprimwprim,angprimwprimmean) |
---|
| 526 | |
---|
| 527 | ! ang transport terms |
---|
| 528 | !---------------------- |
---|
| 529 | |
---|
| 530 | call moyzon2(lonlength,latlength,altlength,miss_val,vangmean,totvang) |
---|
| 531 | call moyzon2(lonlength,latlength,altlength,miss_val,wangmean,totwang) |
---|
| 532 | |
---|
| 533 | do ilat=1,latlength |
---|
| 534 | do ilev=1,altlength |
---|
| 535 | if ((angmeanbar(ilat,ilev).ne.miss_val).and. & |
---|
| 536 | ( vmeanbar(ilat,ilev).ne.miss_val)) then |
---|
| 537 | mmcvang(ilat,ilev) = angmeanbar(ilat,ilev)*vmeanbar(ilat,ilev) |
---|
| 538 | else |
---|
| 539 | mmcvang(ilat,ilev) = miss_val |
---|
| 540 | endif |
---|
| 541 | if ((angmeanbar(ilat,ilev).ne.miss_val).and. & |
---|
| 542 | ( wmeanbar(ilat,ilev).ne.miss_val)) then |
---|
| 543 | mmcwang(ilat,ilev) = angmeanbar(ilat,ilev)*wmeanbar(ilat,ilev) |
---|
| 544 | else |
---|
| 545 | mmcwang(ilat,ilev) = miss_val |
---|
| 546 | endif |
---|
| 547 | enddo |
---|
| 548 | enddo |
---|
| 549 | |
---|
| 550 | call moyzon2(lonlength,latlength,altlength,miss_val,angprimvprimmean,trsvang) |
---|
| 551 | call moyzon2(lonlength,latlength,altlength,miss_val,angprimwprimmean,trswang) |
---|
| 552 | |
---|
| 553 | call moyzon2(lonlength,latlength,altlength,miss_val,angstarmeanvstarmean,stnvang) |
---|
| 554 | call moyzon2(lonlength,latlength,altlength,miss_val,angstarmeanwstarmean,stnwang) |
---|
| 555 | |
---|
| 556 | |
---|
| 557 | print*,"End of computations" |
---|
| 558 | |
---|
| 559 | !=============================================================================== |
---|
| 560 | ! 3. Create output file |
---|
| 561 | !=============================================================================== |
---|
| 562 | |
---|
| 563 | ! Create output file |
---|
| 564 | ierr=NF_CREATE(outfile,NF_CLOBBER,outfid) |
---|
| 565 | if (ierr.ne.NF_NOERR) then |
---|
| 566 | write(*,*)"Error: could not create file ",outfile |
---|
| 567 | stop |
---|
| 568 | endif |
---|
| 569 | |
---|
| 570 | !=============================================================================== |
---|
| 571 | ! 3.1. Define and write dimensions |
---|
| 572 | !=============================================================================== |
---|
| 573 | |
---|
| 574 | call write_dim(outfid,lonlength,latlength,altlength,timelength, & |
---|
| 575 | lon,lat,plev,time,lon_dimid,lat_dimid,alt_dimid,time_dimid) |
---|
| 576 | |
---|
| 577 | !=============================================================================== |
---|
| 578 | ! 3.2. Define and write variables |
---|
| 579 | !=============================================================================== |
---|
| 580 | |
---|
| 581 | datashape2d(1)=lat_dimid |
---|
| 582 | datashape2d(2)=alt_dimid |
---|
| 583 | |
---|
| 584 | call write_var2d(outfid,datashape2d,latlength,altlength,& |
---|
| 585 | "totvang ", "tot hor trpt of ang ","m3 s-2 ",miss_val,& |
---|
| 586 | totvang ) |
---|
| 587 | |
---|
| 588 | call write_var2d(outfid,datashape2d,latlength,altlength,& |
---|
| 589 | "totwang ", "tot ver trpt of ang ","m3 s-2 ",miss_val,& |
---|
| 590 | totwang ) |
---|
| 591 | |
---|
| 592 | call write_var2d(outfid,datashape2d,latlength,altlength,& |
---|
| 593 | "mmcvang ", "MMC hor trpt of ang ","m3 s-2 ",miss_val,& |
---|
| 594 | mmcvang ) |
---|
| 595 | |
---|
| 596 | call write_var2d(outfid,datashape2d,latlength,altlength,& |
---|
| 597 | "mmcwang ", "MMC ver trpt of ang ","m3 s-2 ",miss_val,& |
---|
| 598 | mmcwang ) |
---|
| 599 | |
---|
| 600 | call write_var2d(outfid,datashape2d,latlength,altlength,& |
---|
| 601 | "trsvang ", "trs hor trpt of ang ","m3 s-2 ",miss_val,& |
---|
| 602 | trsvang ) |
---|
| 603 | |
---|
| 604 | call write_var2d(outfid,datashape2d,latlength,altlength,& |
---|
| 605 | "trswang ", "trs ver trpt of ang ","m3 s-2 ",miss_val,& |
---|
| 606 | trswang ) |
---|
| 607 | |
---|
| 608 | call write_var2d(outfid,datashape2d,latlength,altlength,& |
---|
| 609 | "stnvang ", "stn hor trpt of ang ","m3 s-2 ",miss_val,& |
---|
| 610 | stnvang ) |
---|
| 611 | |
---|
| 612 | call write_var2d(outfid,datashape2d,latlength,altlength,& |
---|
| 613 | "stnwang ", "stn ver trpt of ang ","m3 s-2 ",miss_val,& |
---|
| 614 | stnwang ) |
---|
| 615 | |
---|
| 616 | call write_var2d(outfid,datashape2d,latlength,altlength,& |
---|
| 617 | "dmass ", "mass in each cell ","kg ",miss_val,& |
---|
| 618 | dmassmeanbar ) |
---|
| 619 | |
---|
| 620 | |
---|
| 621 | !!!! Close output file |
---|
| 622 | ierr=NF_CLOSE(outfid) |
---|
| 623 | if (ierr.ne.NF_NOERR) write(*,*) 'Error, failed to close output file ',outfile |
---|
| 624 | |
---|
| 625 | |
---|
| 626 | end program |
---|