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

Last change on this file since 3439 was 1454, checked in by slebonnois, 9 years ago

SL: corrections in Titan and Venus tools

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