source: LMDZ5/branches/LMDZ5V1.0-dev/libf/phylmd/cv3a_uncompress.F @ 2185

Last change on this file since 2185 was 1398, checked in by musat, 14 years ago

Last corrections for CMIP5:

  • Add O3 at standard level files histmthNMC.nc
  • Add positive attribute "down" for vertical axes for all output files
  • Replace "inst" by "ave" for hist*NMC.nc files to have time_counter and bounds for time axis (Marie-Alice's hint)
  • Correct units for vertical axes : mb instead of hPa
  • Add mass flux at the bottom of clouds
  • Comment non initialized variables (s_capCL, s_oliqCL, s_cteiCL, s_trmb1, s_trmb2, s_trmb3) for the output files
  • Geopotential field phy850, phi700, phi500, etc are modified to "geopotential height and are called z850, z700, z500, etc
  • Meaning of specific humidity outputs - ovapinit and ovap - were interchanged
  • Fields albs, albslw become alb1, alb2 in output files
  • Correct title for rugs_* fields
  • Correct units for pbase and ptop are Pa (not mb)
  • Correct ndayrain field

FH/JLD/JYG/MAF/IM

  • Property svn:eol-style set to native
  • Property svn:keywords set to Author Date Id Revision
File size: 4.2 KB
Line 
1      SUBROUTINE cv3a_uncompress(nloc,len,ncum,nd,ntra,idcum
2     :         ,iflag,kbas,ktop
3     :         ,precip,cbmf,sig,w0,ptop2
4     :         ,ft,fq,fu,fv,ftra
5     :         ,Ma,mip,Vprecip,upwd,dnwd,dnwd0
6     :         ,qcondc,wd,cape,cin
7     :         ,tvp
8     :         ,ftd,fqd
9     :         ,Plim1,Plim2,asupmax,supmax0
10     :         ,asupmaxmin
11     o         ,iflag1,kbas1,ktop1
12     :         ,precip1,cbmf1,sig1,w01,ptop21
13     :         ,ft1,fq1,fu1,fv1,ftra1
14     :         ,Ma1,mip1,Vprecip1,upwd1,dnwd1,dnwd01
15     :         ,qcondc1,wd1,cape1,cin1
16     :         ,tvp1
17     :         ,ftd1,fqd1
18     :         ,Plim11,Plim21,asupmax1,supmax01
19     :         ,asupmaxmin1     )
20***************************************************************
21*                                                             *
22* CV3A_UNCOMPRESS                                             *
23*                                                             *
24*                                                             *
25* written by   : Sandrine Bony-Lena , 17/05/2003, 11.22.15    *
26* modified by  : Jean-Yves Grandpeix, 23/06/2003, 10.36.17    *
27***************************************************************
28*
29      implicit none
30
31#include "cv3param.h"
32
33c inputs:
34      integer nloc, len, ncum, nd, ntra
35      integer idcum(nloc)
36      integer iflag(nloc),kbas(nloc),ktop(nloc)
37      real precip(nloc),cbmf(nloc)
38      real sig(nloc,nd), w0(nloc,nd),ptop2(nloc)
39      real ft(nloc,nd), fq(nloc,nd), fu(nloc,nd), fv(nloc,nd)
40      real ftra(nloc,nd,ntra)
41      real Ma(nloc,nd),mip(nloc,nd),Vprecip(nloc,nd+1)
42      real upwd(nloc,nd),dnwd(nloc,nd),dnwd0(nloc,nd)
43      real qcondc(nloc,nd)
44      real wd(nloc),cape(nloc),cin(nloc)
45      real tvp(nloc,nd)
46      real ftd(nloc,nd), fqd(nloc,nd)
47      real Plim1(nloc),Plim2(nloc)
48      real asupmax(nloc,nd),supmax0(nloc)
49      real asupmaxmin(nloc)
50
51c outputs:
52      integer iflag1(len),kbas1(len),ktop1(len)
53      real precip1(len),cbmf1(len)
54      real sig1(len,nd), w01(len,nd),ptop21(len)
55      real ft1(len,nd), fq1(len,nd), fu1(len,nd), fv1(len,nd)
56      real ftra1(len,nd,ntra)
57      real Ma1(len,nd),mip1(len,nd),Vprecip1(len,nd+1)
58      real upwd1(len,nd),dnwd1(len,nd),dnwd01(len,nd)
59      real qcondc1(len,nd)
60      real wd1(len),cape1(len),cin1(len)
61      real tvp1(len,nd)
62      real ftd1(len,nd), fqd1(len,nd)
63      real Plim11(len),Plim21(len)
64      real asupmax1(len,nd),supmax01(len)
65      real asupmaxmin1(len)
66c
67c local variables:
68      integer i,k,j,k1,k2
69
70        do 2000 i=1,ncum
71         ptop21(idcum(i))=ptop2(i)
72         precip1(idcum(i))=precip(i)
73         cbmf1(idcum(i))=cbmf(i)
74         iflag1(idcum(i))=iflag(i)
75         kbas1(idcum(i))=kbas(i)
76         ktop1(idcum(i))=ktop(i)
77         wd1(idcum(i))=wd(i)
78         cape1(idcum(i))=cape(i)
79         cin1(idcum(i))=cin(i)
80         Plim11(idcum(i))=Plim1(i)
81         Plim21(idcum(i))=Plim2(i)
82         supmax01(idcum(i))=supmax0(i)
83         asupmaxmin1(idcum(i))=asupmaxmin(i)
84 2000   continue
85
86        do 2020 k=1,nd
87          do 2010 i=1,ncum
88            sig1(idcum(i),k)=sig(i,k)
89            w01(idcum(i),k)=w0(i,k)
90            ft1(idcum(i),k)=ft(i,k)
91            fq1(idcum(i),k)=fq(i,k)
92            fu1(idcum(i),k)=fu(i,k)
93            fv1(idcum(i),k)=fv(i,k)
94            Ma1(idcum(i),k)=Ma(i,k)
95            mip1(idcum(i),k)=mip(i,k)
96            Vprecip1(idcum(i),k)=Vprecip(i,k)
97            upwd1(idcum(i),k)=upwd(i,k)
98            dnwd1(idcum(i),k)=dnwd(i,k)
99            dnwd01(idcum(i),k)=dnwd0(i,k)
100            qcondc1(idcum(i),k)=qcondc(i,k)
101            tvp1(idcum(i),k)=tvp(i,k)
102            ftd1(idcum(i),k)=ftd(i,k)
103            fqd1(idcum(i),k)=fqd(i,k)
104            asupmax1(idcum(i),k)=asupmax(i,k)
105 2010     continue
106 2020   continue
107
108        do 2040 i=1,ncum
109          sig1(idcum(i),nd)=sig(i,nd)
1102040    continue
111
112
113        do 2100 j=1,ntra
114c oct3         do 2110 k=1,nl
115         do 2110 k=1,nd ! oct3
116          do 2120 i=1,ncum
117            ftra1(idcum(i),k,j)=ftra(i,k,j)
118 2120     continue
119 2110    continue
120 2100   continue
121c
122c        do 2220 k2=1,nd
123c         do 2210 k1=1,nd
124c          do 2200 i=1,ncum
125c            ment1(idcum(i),k1,k2) = ment(i,k1,k2)
126c            sij1(idcum(i),k1,k2) = sij(i,k1,k2)
127c2200      enddo
128c2210     enddo
129c2220    enddo
130
131       RETURN
132      END
133
Note: See TracBrowser for help on using the repository browser.