source: LMDZ5/trunk/libf/phylmd/moy_undefSTD.F @ 1916

Last change on this file since 1916 was 1915, checked in by musat, 11 years ago

D'autres oublis pour les sorties annuelles
IM

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