source: LMDZ4/branches/LMDZ4V5.0-dev/libf/phylmd/moy_undefSTD.F @ 1350

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

Add 3 output files for standard pressure levels AR5 exercice and flags
to manage their computation and output frequencies
histhfNMC.nc with 3 standard pressure levels
histdayNMC.nc with 8 (or may have 17) standard pressure levels
histmthNMC.nc with 17 standard pressure levels
Add 3 flags in the physiq.def file: freq_calNMC(3), freq_outNMC(3) and lev_histdayNMC
freq_calNMC(3) : computation frequency of variables on standard pressure levels

and by default has fallowing values (in fact physics' time step dtime)

freq_calNMC(1)=900.
freq_calNMC(2)=900.
freq_calNMC(3)=900.
freq_outNMC(3) : output frequency of variables on standard pressure levels

with following default values

freq_out(1) = 2592000. (30 days)
freq_out(2) = 86400. (1 day)
freq_out(3) = 21600. (6 hours)
lev_histdayNMC is 8 by default but may be switched to 17 (if we need more levels for a particular run)
IM

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