source: trunk/LMDZ.VENUS/Tools/dmass.F90 @ 1243

Last change on this file since 1243 was 816, checked in by slebonnois, 12 years ago

SL: tools for postprocessing (Veznus and Titan); see DOC/documentation/vt-tools.pdf

File size: 5.4 KB
RevLine 
[816]1subroutine cellmass(infid,latlength,lonlength,altlength,timelength, &
2                     lmdflag,deltalon,deltalat,latrad,plev,ps,grav,rayon, &
3                     dmass )
4
5implicit none
6
7include "netcdf.inc" ! NetCDF definitions
8
9! arguments
10integer infid ! NetCDF input file ID
11integer lonlength ! # of grid points along longitude
12integer latlength ! # of grid points along latitude
13integer altlength ! # of grid point along altitude (of input datasets)
14integer timelength ! # of points along time
15logical lmdflag ! true=LMD, false=CAM
16real :: deltalat,deltalon ! lat and lon intervals in radians
17real,dimension(latlength) :: latrad ! latitudes in radians
18real,dimension(altlength) :: plev ! Pressure levels (Pa)
19real,dimension(lonlength,latlength,timelength),intent(in) :: ps ! surface pressure
20real,dimension(lonlength,latlength,altlength,timelength),intent(in) :: grav,rayon
21real,dimension(lonlength,latlength,altlength,timelength),intent(out) :: dmass
22
23!local
24character (len=64) :: text ! to store some text
25real,dimension(:),allocatable :: aps,bps ! hybrid vertical coordinates
26real,dimension(:,:,:,:),allocatable :: press ! GCM atmospheric pressure
27real,dimension(:),allocatable :: pp ! Pressure levels (Pa)
28real,dimension(:,:,:,:),allocatable :: deltap ! pressure thickness of each layer (Pa)
29integer ierr,ierr1,ierr2 ! NetCDF routines return codes
30integer i,j,ilon,ilat,ilev,itim ! for loops
31real :: tmpp ! temporary pressure
32real :: p0 ! GCM reference pressure
33integer tmpvarid
34real :: miss_val=9.99e+33 ! special "missing value" to specify missing data
35real,parameter :: miss_val_def=9.99e+33 ! default value for "missing value"
36
37include "planet.h"
38
39!===============================================================================
40! 2.2.1 Pressure intervals
41!===============================================================================
42allocate(deltap(lonlength,latlength,altlength,timelength))
43allocate(press(lonlength,latlength,altlength,timelength))
44allocate(pp(altlength))
45
46if(lmdflag) then  ! LMD vs CAM
47  text="pres"
48  call get_var4d(infid,lonlength,latlength,altlength,timelength,text,press,ierr1,ierr2)
49  if (ierr1.ne.NF_NOERR) then
50! assume we are using _P files
51do itim=1,timelength
52 do ilon=1,lonlength
53  do ilat=1,latlength
54     press(ilon,ilat,:,itim) = plev
55  enddo
56 enddo
57enddo   
58  else
59    if (ierr2.ne.NF_NOERR) stop "Error: Failed reading pres"
60  endif
61
62else  ! LMD vs CAM
63  ierr=NF_INQ_VARID(infid,"P0",tmpvarid)
64  if (ierr.ne.NF_NOERR) then
65    write(*,*) "Failed to get P0 ID, used psref=",psref
66    p0=psref
67  else
68    ierr=NF_GET_VAR_REAL(infid,tmpvarid,p0)
69    if (ierr.ne.NF_NOERR) then
70      write(*,*) "Failed reading reference pressure, used psref=",psref
71      p0=psref
72    endif
73  endif
74  ! hybrid coordinate aps/hyam
75  text="hyam"
76  ierr=NF_INQ_VARID(infid,text,tmpvarid)
77  if (ierr.ne.NF_NOERR) then
78    stop "Error: Could not find aps/hyam coordinates"
79  else
80    allocate(aps(altlength))
81    ierr=NF_GET_VAR_REAL(infid,tmpvarid,aps)
82    if (ierr.ne.NF_NOERR) then
83      stop "Error: Failed reading aps/hyam"
84    endif
85    call reverselev(altlength,aps)
86  endif
87
88  ! hybrid coordinate hybm
89  text="hybm"
90  ierr=NF_INQ_VARID(infid,text,tmpvarid)
91  if (ierr.ne.NF_NOERR) then
92    stop "Error: Failed to get bps/hybm ID"
93  else
94    allocate(bps(altlength))
95    ierr=NF_GET_VAR_REAL(infid,tmpvarid,bps)
96    if (ierr.ne.NF_NOERR) then
97      stop "Error: Failed reading bps/hybm"
98    endif
99    call reverselev(altlength,bps)
100  endif
101
102  aps=aps*p0
103
104  do itim=1,timelength
105    do ilev=1,altlength
106      do ilat=1,latlength
107        do ilon=1,lonlength
108          press(ilon,ilat,ilev,itim)=aps(ilev)+bps(ilev)*ps(ilon,ilat,itim)
109        enddo
110      enddo
111    enddo
112  enddo
113
114endif  ! LMD vs CAM
115
116do itim=1,timelength
117
118do ilon=1,lonlength
119 do ilat=1,latlength
120  tmpp=ps(ilon,ilat,itim)
121  pp=press(ilon,ilat,:,itim)
122if (tmpp.ge.pp(1)) then
123  deltap(ilon,ilat,1,itim)= tmpp - exp((log(pp(1))+log(pp(2)))/2.) ! initialization, rectified with ps below
124else
125  deltap(ilon,ilat,1,itim)= miss_val
126endif
127do ilev=2,altlength-1
128  deltap(ilon,ilat,ilev,itim)= exp((log(pp(ilev-1))+log(pp(ilev)))/2.)-&
129                     exp((log(pp(ilev))+log(pp(ilev+1)))/2.)
130enddo
131deltap(ilon,ilat,altlength,itim)=exp((log(pp(altlength-1))+log(pp(altlength)))/2.)
132 enddo
133enddo
134
135do ilon=1,lonlength
136 do ilat=1,latlength
137  tmpp=ps(ilon,ilat,itim)
138  pp=press(ilon,ilat,:,itim)
139  do ilev=altlength,2,-1
140    if ((pp(ilev).le.tmpp).and.(pp(ilev-1).gt.tmpp)) then
141      deltap(ilon,ilat,ilev,itim)= tmpp - &
142               exp((log(pp(ilev))+log(pp(ilev+1)))/2.)
143    elseif (pp(ilev).gt.tmpp) then
144      deltap(ilon,ilat,ilev,itim)=miss_val
145    endif
146  enddo
147 enddo
148enddo
149
150enddo ! timelength
151
152!===============================================================================
153! 2.2.2 Mass in cells
154!===============================================================================
155
156do itim=1,timelength
157
158do ilon=1,lonlength
159 do ilat=1,latlength
160  do ilev=1,altlength
161    if ((deltap(ilon,ilat,ilev,itim).ne.miss_val) &
162   .and. (rayon(ilon,ilat,ilev,itim).ne.miss_val) &
163   .and.  (grav(ilon,ilat,ilev,itim).ne.miss_val)) then
164     dmass(ilon,ilat,ilev,itim) = &
165               rayon(ilon,ilat,ilev,itim)*cos(latrad(ilat))*deltalon &
166             * rayon(ilon,ilat,ilev,itim)*deltalat &
167             * deltap(ilon,ilat,ilev,itim)/grav(ilon,ilat,ilev,itim)
168    else
169     dmass(ilon,ilat,ilev,itim) = miss_val
170    endif
171  enddo
172 enddo
173enddo
174
175enddo ! timelength
176
177end subroutine cellmass
Note: See TracBrowser for help on using the repository browser.