source: LMDZ4/trunk/libf/phylmd/moy_undefSTD.F @ 1352

Last change on this file since 1352 was 1352, 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.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
58c      print*,'n freq_ini=',n,itap,freq_outNMC(n),missing_val
59c
60        DO k=1, nlevSTD
61         DO i=1, klon
62          IF(tnondef(i,k,n).NE.(freq_moyNMC(n))) THEN
63           tsumSTD(i,k,n)=tsumSTD(i,k,n)/
64     $     (freq_moyNMC(n)-tnondef(i,k,n))
65cIM BEG
66          usumSTD(i,k,n)=usumSTD(i,k,n)/
67     $     (freq_moyNMC(n)-tnondef(i,k,n))
68          vsumSTD(i,k,n)=vsumSTD(i,k,n)/
69     $     (freq_moyNMC(n)-tnondef(i,k,n))
70          wsumSTD(i,k,n)=wsumSTD(i,k,n)/
71     $     (freq_moyNMC(n)-tnondef(i,k,n))
72          phisumSTD(i,k,n)=phisumSTD(i,k,n)/
73     $     (freq_moyNMC(n)-tnondef(i,k,n))
74          qsumSTD(i,k,n)=qsumSTD(i,k,n)/
75     $     (freq_moyNMC(n)-tnondef(i,k,n))
76          rhsumSTD(i,k,n)=rhsumSTD(i,k,n)/
77     $     (freq_moyNMC(n)-tnondef(i,k,n))
78          uvsumSTD(i,k,n)=uvsumSTD(i,k,n)/
79     $     (freq_moyNMC(n)-tnondef(i,k,n))
80          vqsumSTD(i,k,n)=vqsumSTD(i,k,n)/
81     $     (freq_moyNMC(n)-tnondef(i,k,n))
82          vTsumSTD(i,k,n)=vTsumSTD(i,k,n)/
83     $     (freq_moyNMC(n)-tnondef(i,k,n))
84          wqsumSTD(i,k,n)=wqsumSTD(i,k,n)/
85     $     (freq_moyNMC(n)-tnondef(i,k,n))
86          vphisumSTD(i,k,n)=vphisumSTD(i,k,n)/
87     $     (freq_moyNMC(n)-tnondef(i,k,n))
88          wTsumSTD(i,k,n)=wTsumSTD(i,k,n)/
89     $     (freq_moyNMC(n)-tnondef(i,k,n))
90          u2sumSTD(i,k,n)=u2sumSTD(i,k,n)/
91     $     (freq_moyNMC(n)-tnondef(i,k,n))
92          v2sumSTD(i,k,n)=v2sumSTD(i,k,n)/
93     $     (freq_moyNMC(n)-tnondef(i,k,n))
94          T2sumSTD(i,k,n)=T2sumSTD(i,k,n)/
95     $     (freq_moyNMC(n)-tnondef(i,k,n))
96cIM END
97          ELSE
98           tsumSTD(i,k,n)=missing_val
99           usumSTD(i,k,n)=missing_val
100           vsumSTD(i,k,n)=missing_val
101           wsumSTD(i,k,n)=missing_val
102           phisumSTD(i,k,n)=missing_val
103           qsumSTD(i,k,n)=missing_val
104           rhsumSTD(i,k,n)=missing_val
105           uvsumSTD(i,k,n)=missing_val
106           vqsumSTD(i,k,n)=missing_val
107           vTsumSTD(i,k,n)=missing_val
108           wqsumSTD(i,k,n)=missing_val
109           vphisumSTD(i,k,n)=missing_val
110           wTsumSTD(i,k,n)=missing_val
111           u2sumSTD(i,k,n)=missing_val
112           v2sumSTD(i,k,n)=missing_val
113           T2sumSTD(i,k,n)=missing_val
114          ENDIF !tnondef(i,k,n).NE.(freq_moyNMC(n))
115c
116c          if(n.EQ.1.AND.i.EQ.4513.AND.k.EQ.17) THEN
117           if(n.EQ.1.AND.i.EQ.1128.AND.k.EQ.17) THEN
118            print*,'itap rlon rlat tlevSTD=',itap,rlon(i),rlat(i),
119     $tlevSTD(i,k)
120           endif
121c
122         ENDDO !i
123        ENDDO !k
124       ENDIF !MOD(itap,NINT(freq_outNMC(n)/dtime)).EQ.0
125c
126      ENDDO !n
127c
128      RETURN
129      END 
Note: See TracBrowser for help on using the repository browser.