source: LMDZ6/branches/Amaury_dev/libf/phylmd/lmdz_moy_undefstd.f90 @ 5220

Last change on this file since 5220 was 5160, checked in by abarral, 3 months ago

Put .h into modules

  • 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: 5.4 KB
Line 
1MODULE lmdz_moy_undefstd
2  IMPLICIT NONE; PRIVATE
3  PUBLIC moy_undefstd
4CONTAINS
5  SUBROUTINE moy_undefstd(itap, itapm1)
6    USE netcdf, ONLY: nf90_fill_real
7    USE dimphy
8    USE phys_state_var_mod
9    USE lmdz_wxios, ONLY: missing_val_xios => missing_val, using_xios
10
11    USE phys_cal_mod, ONLY: mth_len
12    USE lmdz_clesphys
13
14    IMPLICIT NONE
15    REAL :: missing_val
16
17    ! ====================================================================
18
19    ! I. Musat : 09.2004
20
21    ! Moyenne - a des frequences differentes - des valeurs bien definies
22    ! (.NE.missing_val) des variables interpolees a un niveau de
23    ! pression.
24    ! 1) les variables de type "day" (nout=1) ou "mth" (nout=2) sont sommees
25    ! tous les pas de temps de la physique
26
27    ! 2) les variables de type "NMC" (nout=3) sont calculees a partir
28    ! des valeurs instantannees toutes les 6 heures
29
30
31    ! NB: mettre "inst(X)" dans le write_hist*NMC.h !
32    ! ====================================================================
33
34
35    ! variables Input
36    ! INTEGER nlevSTD, klevSTD, itap
37    ! PARAMETER(klevSTD=17)
38    INTEGER itap, itapm1
39
40    ! variables locales
41    ! INTEGER i, k, nout, n
42    ! PARAMETER(nout=3) !nout=1 day/nout=2 mth/nout=3 NMC
43    INTEGER i, k, n
44    ! REAL dtime, freq_outNMC(nout), freq_moyNMC(nout)
45    ! REAL freq_outNMC(nout), freq_calNMC(nout)
46    REAL freq_moynmc(nout)
47
48    ! variables Output
49    ! REAL tnondef(klon,klevSTD,nout)
50    ! REAL tsumSTD(klon,klevSTD,nout)
51
52    REAL un_jour
53    PARAMETER (un_jour = 86400.)
54    ! REAL missing_val
55
56    ! missing_val = nf90_fill_real
57    IF (using_xios) THEN
58      missing_val = missing_val_xios
59    ELSE
60      missing_val = missing_val_nf90
61    ENDIF
62
63    DO n = 1, nout
64      IF (freq_outnmc(n)<0) THEN
65        freq_moynmc(n) = (mth_len * un_jour) / freq_calnmc(n)
66        ! PRINT*,'moy_undefSTD n freq_out freq_moy =',
67        ! $n,freq_moyNMC(n)
68      ELSE
69        freq_moynmc(n) = freq_outnmc(n) / freq_calnmc(n)
70      END IF
71
72      ! calcul 1 fois pas mois, 1 fois par jour ou toutes les 6h
73
74      IF (n==1 .AND. itap==itapm1 .OR. n>1 .AND. mod(itap, nint(freq_outnmc(n) / &
75              phys_tstep))==0) THEN
76
77        ! PRINT*,'moy_undefSTD n itap itapm1',n,itap,itapm1
78
79        DO k = 1, nlevstd
80          DO i = 1, klon
81            IF (tnondef(i, k, n)/=(freq_moynmc(n))) THEN
82              tsumstd(i, k, n) = tsumstd(i, k, n) / (freq_moynmc(n) - tnondef(i, k, n &
83                      ))
84              usumstd(i, k, n) = usumstd(i, k, n) / (freq_moynmc(n) - tnondef(i, k, n &
85                      ))
86              vsumstd(i, k, n) = vsumstd(i, k, n) / (freq_moynmc(n) - tnondef(i, k, n &
87                      ))
88              wsumstd(i, k, n) = wsumstd(i, k, n) / (freq_moynmc(n) - tnondef(i, k, n &
89                      ))
90              phisumstd(i, k, n) = phisumstd(i, k, n) / &
91                      (freq_moynmc(n) - tnondef(i, k, n))
92              qsumstd(i, k, n) = qsumstd(i, k, n) / (freq_moynmc(n) - tnondef(i, k, n &
93                      ))
94              rhsumstd(i, k, n) = rhsumstd(i, k, n) / (freq_moynmc(n) - tnondef(i, k &
95                      , n))
96              uvsumstd(i, k, n) = uvsumstd(i, k, n) / (freq_moynmc(n) - tnondef(i, k &
97                      , n))
98              vqsumstd(i, k, n) = vqsumstd(i, k, n) / (freq_moynmc(n) - tnondef(i, k &
99                      , n))
100              vtsumstd(i, k, n) = vtsumstd(i, k, n) / (freq_moynmc(n) - tnondef(i, k &
101                      , n))
102              wqsumstd(i, k, n) = wqsumstd(i, k, n) / (freq_moynmc(n) - tnondef(i, k &
103                      , n))
104              vphisumstd(i, k, n) = vphisumstd(i, k, n) / &
105                      (freq_moynmc(n) - tnondef(i, k, n))
106              wtsumstd(i, k, n) = wtsumstd(i, k, n) / (freq_moynmc(n) - tnondef(i, k &
107                      , n))
108              u2sumstd(i, k, n) = u2sumstd(i, k, n) / (freq_moynmc(n) - tnondef(i, k &
109                      , n))
110              v2sumstd(i, k, n) = v2sumstd(i, k, n) / (freq_moynmc(n) - tnondef(i, k &
111                      , n))
112              t2sumstd(i, k, n) = t2sumstd(i, k, n) / (freq_moynmc(n) - tnondef(i, k &
113                      , n))
114              o3sumstd(i, k, n) = o3sumstd(i, k, n) / (freq_moynmc(n) - tnondef(i, k &
115                      , n))
116              o3daysumstd(i, k, n) = o3daysumstd(i, k, n) / &
117                      (freq_moynmc(n) - tnondef(i, k, n))
118            ELSE
119              tsumstd(i, k, n) = missing_val
120              usumstd(i, k, n) = missing_val
121              vsumstd(i, k, n) = missing_val
122              wsumstd(i, k, n) = missing_val
123              phisumstd(i, k, n) = missing_val
124              qsumstd(i, k, n) = missing_val
125              rhsumstd(i, k, n) = missing_val
126              uvsumstd(i, k, n) = missing_val
127              vqsumstd(i, k, n) = missing_val
128              vtsumstd(i, k, n) = missing_val
129              wqsumstd(i, k, n) = missing_val
130              vphisumstd(i, k, n) = missing_val
131              wtsumstd(i, k, n) = missing_val
132              u2sumstd(i, k, n) = missing_val
133              v2sumstd(i, k, n) = missing_val
134              t2sumstd(i, k, n) = missing_val
135              o3sumstd(i, k, n) = missing_val
136              o3daysumstd(i, k, n) = missing_val
137            END IF !tnondef(i,k,n).NE.(freq_moyNMC(n))
138          END DO !i
139        END DO !k
140      END IF !MOD(itap,NINT(freq_outNMC(n)/phys_tstep)).EQ.0
141
142    END DO !n
143
144  END SUBROUTINE moy_undefstd
145END MODULE lmdz_moy_undefstd
Note: See TracBrowser for help on using the repository browser.