source: LMDZ4/branches/LMDZ4_AR5/libf/phylmd/moy_undefSTD.F @ 2532

Last change on this file since 2532 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.0 KB
Line 
1!
2! $Header$
3!
4      SUBROUTINE moy_undefSTD(itap,freq_outNMC,freq_moyNMC)
5      USE netcdf
6      USE dimphy
7      USE phys_state_var_mod ! Variables sauvegardees de la physique
8      IMPLICIT none
9c
10c====================================================================
11c
12c I. Musat : 09.2004
13c
14c Moyenne - a des frequences differentes - des valeurs bien definies
15c         (.NE.missing_val) des variables interpolees a un niveau de
16c         pression.
17c 1) les variables de type "day" (nout=1) ou "mth" (nout=2) sont sommees
18c    tous les pas de temps de la physique
19c
20c 2) les variables de type "NMC" (nout=3) sont calculees a partir
21c    des valeurs instantannees toutes les 6 heures
22c
23c
24c NB: mettre "inst(X)" dans le write_hist*NMC.h !
25c====================================================================
26cym#include "dimensions.h"
27cym      integer jjmp1
28cym      parameter (jjmp1=jjm+1-1/jjm)
29cym#include "dimphy.h"
30c
31c
32c variables Input
33c     INTEGER nlevSTD, klevSTD, itap
34c     PARAMETER(klevSTD=17)
35      INTEGER itap
36c
37c variables locales
38c     INTEGER i, k, nout, n
39c     PARAMETER(nout=3) !nout=1 day/nout=2 mth/nout=3 NMC
40      INTEGER i, k, n
41c     REAL dtime, freq_outNMC(nout), freq_moyNMC(nout)
42      REAL freq_outNMC(nout), freq_moyNMC(nout)
43c
44c variables Output
45c     REAL tnondef(klon,klevSTD,nout)
46c     REAL tsumSTD(klon,klevSTD,nout)
47c
48      REAL missing_val
49c
50      missing_val=nf90_fill_real
51c
52      DO n=1, nout
53c
54c calcul 1 fois par jour
55c
56       IF(MOD(itap,NINT(freq_outNMC(n)/dtime)).EQ.0) THEN
57c
58        DO k=1, nlevSTD
59         DO i=1, klon
60          IF(tnondef(i,k,n).NE.(freq_moyNMC(n))) THEN
61           tsumSTD(i,k,n)=tsumSTD(i,k,n)/
62     $     (freq_moyNMC(n)-tnondef(i,k,n))
63cIM BEG
64          usumSTD(i,k,n)=usumSTD(i,k,n)/
65     $     (freq_moyNMC(n)-tnondef(i,k,n))
66          vsumSTD(i,k,n)=vsumSTD(i,k,n)/
67     $     (freq_moyNMC(n)-tnondef(i,k,n))
68          wsumSTD(i,k,n)=wsumSTD(i,k,n)/
69     $     (freq_moyNMC(n)-tnondef(i,k,n))
70          phisumSTD(i,k,n)=phisumSTD(i,k,n)/
71     $     (freq_moyNMC(n)-tnondef(i,k,n))
72          qsumSTD(i,k,n)=qsumSTD(i,k,n)/
73     $     (freq_moyNMC(n)-tnondef(i,k,n))
74          rhsumSTD(i,k,n)=rhsumSTD(i,k,n)/
75     $     (freq_moyNMC(n)-tnondef(i,k,n))
76          uvsumSTD(i,k,n)=uvsumSTD(i,k,n)/
77     $     (freq_moyNMC(n)-tnondef(i,k,n))
78          vqsumSTD(i,k,n)=vqsumSTD(i,k,n)/
79     $     (freq_moyNMC(n)-tnondef(i,k,n))
80          vTsumSTD(i,k,n)=vTsumSTD(i,k,n)/
81     $     (freq_moyNMC(n)-tnondef(i,k,n))
82          wqsumSTD(i,k,n)=wqsumSTD(i,k,n)/
83     $     (freq_moyNMC(n)-tnondef(i,k,n))
84          vphisumSTD(i,k,n)=vphisumSTD(i,k,n)/
85     $     (freq_moyNMC(n)-tnondef(i,k,n))
86          wTsumSTD(i,k,n)=wTsumSTD(i,k,n)/
87     $     (freq_moyNMC(n)-tnondef(i,k,n))
88          u2sumSTD(i,k,n)=u2sumSTD(i,k,n)/
89     $     (freq_moyNMC(n)-tnondef(i,k,n))
90          v2sumSTD(i,k,n)=v2sumSTD(i,k,n)/
91     $     (freq_moyNMC(n)-tnondef(i,k,n))
92          T2sumSTD(i,k,n)=T2sumSTD(i,k,n)/
93     $     (freq_moyNMC(n)-tnondef(i,k,n))
94          O3sumSTD(i,k,n)=O3sumSTD(i,k,n)/
95     $     (freq_moyNMC(n)-tnondef(i,k,n))
96          O3daysumSTD(i,k,n)=O3daysumSTD(i,k,n)/
97     $     (freq_moyNMC(n)-tnondef(i,k,n))
98cIM END
99          ELSE
100           tsumSTD(i,k,n)=missing_val
101           usumSTD(i,k,n)=missing_val
102           vsumSTD(i,k,n)=missing_val
103           wsumSTD(i,k,n)=missing_val
104           phisumSTD(i,k,n)=missing_val
105           qsumSTD(i,k,n)=missing_val
106           rhsumSTD(i,k,n)=missing_val
107           uvsumSTD(i,k,n)=missing_val
108           vqsumSTD(i,k,n)=missing_val
109           vTsumSTD(i,k,n)=missing_val
110           wqsumSTD(i,k,n)=missing_val
111           vphisumSTD(i,k,n)=missing_val
112           wTsumSTD(i,k,n)=missing_val
113           u2sumSTD(i,k,n)=missing_val
114           v2sumSTD(i,k,n)=missing_val
115           T2sumSTD(i,k,n)=missing_val
116           O3sumSTD(i,k,n)=missing_val
117           O3daysumSTD(i,k,n)=missing_val
118          ENDIF !tnondef(i,k,n).NE.(freq_moyNMC(n))
119         ENDDO !i
120        ENDDO !k
121       ENDIF !MOD(itap,NINT(freq_outNMC(n)/dtime)).EQ.0
122c
123      ENDDO !n
124c
125      RETURN
126      END 
Note: See TracBrowser for help on using the repository browser.