source: LMDZ5/branches/LMDZ5_SPLA/libf/phylmd/moy_undefSTD.F90 @ 5441

Last change on this file since 5441 was 1992, checked in by lguez, 11 years ago

Converted to free source form files in libf/phylmd which were still in
fixed source form. The conversion was done using the polish mode of
the NAG Fortran Compiler.

In addition to converting to free source form, the processing of the
files also:

-- indented the code (including comments);

-- set Fortran keywords to uppercase, and set all other identifiers
to lower case;

-- added qualifiers to end statements (for example "end subroutine
conflx", instead of "end");

-- changed the terminating statements of all DO loops so that each
loop ends with an ENDDO statement (instead of a labeled continue).

-- replaced #include by include.

  • 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.7 KB
Line 
1
2! $Id: moy_undefSTD.F90 1992 2014-03-05 13:19:12Z fhourdin $
3
4SUBROUTINE 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"
11
12  ! ====================================================================
13
14  ! I. Musat : 09.2004
15
16  ! Moyenne - a des frequences differentes - des valeurs bien definies
17  ! (.NE.missing_val) des variables interpolees a un niveau de
18  ! pression.
19  ! 1) les variables de type "day" (nout=1) ou "mth" (nout=2) sont sommees
20  ! tous les pas de temps de la physique
21
22  ! 2) les variables de type "NMC" (nout=3) sont calculees a partir
23  ! des valeurs instantannees toutes les 6 heures
24
25
26  ! NB: mettre "inst(X)" dans le write_hist*NMC.h !
27  ! ====================================================================
28  ! ym#include "dimensions.h"
29  ! ym      integer jjmp1
30  ! ym      parameter (jjmp1=jjm+1-1/jjm)
31  ! ym#include "dimphy.h"
32
33
34  ! variables Input
35  ! INTEGER nlevSTD, klevSTD, itap
36  ! PARAMETER(klevSTD=17)
37  INTEGER itap, itapm1
38
39  ! variables locales
40  ! INTEGER i, k, nout, n
41  ! PARAMETER(nout=3) !nout=1 day/nout=2 mth/nout=3 NMC
42  INTEGER i, k, n
43  ! REAL dtime, freq_outNMC(nout), freq_moyNMC(nout)
44  ! REAL freq_outNMC(nout), freq_calNMC(nout)
45  REAL freq_moynmc(nout)
46
47  ! variables Output
48  ! REAL tnondef(klon,klevSTD,nout)
49  ! REAL tsumSTD(klon,klevSTD,nout)
50
51  REAL un_jour
52  PARAMETER (un_jour=86400.)
53  REAL missing_val
54
55  missing_val = nf90_fill_real
56
57  DO n = 1, nout
58    IF (freq_outnmc(n)<0) THEN
59      freq_moynmc(n) = (mth_len*un_jour)/freq_calnmc(n)
60      ! print*,'moy_undefSTD n freq_out freq_moy =',
61      ! $n,freq_moyNMC(n)
62    ELSE
63      freq_moynmc(n) = freq_outnmc(n)/freq_calnmc(n)
64    END IF
65
66    ! calcul 1 fois pas mois, 1 fois par jour ou toutes les 6h
67
68    IF (n==1 .AND. itap==itapm1 .OR. n>1 .AND. mod(itap,nint(freq_outnmc(n)/ &
69        dtime))==0) THEN
70
71      ! print*,'moy_undefSTD n itap itapm1',n,itap,itapm1
72
73      DO k = 1, nlevstd
74        DO i = 1, klon
75          IF (tnondef(i,k,n)/=(freq_moynmc(n))) THEN
76            tsumstd(i, k, n) = tsumstd(i, k, n)/(freq_moynmc(n)-tnondef(i,k,n &
77              ))
78            usumstd(i, k, n) = usumstd(i, k, n)/(freq_moynmc(n)-tnondef(i,k,n &
79              ))
80            vsumstd(i, k, n) = vsumstd(i, k, n)/(freq_moynmc(n)-tnondef(i,k,n &
81              ))
82            wsumstd(i, k, n) = wsumstd(i, k, n)/(freq_moynmc(n)-tnondef(i,k,n &
83              ))
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)/(freq_moynmc(n)-tnondef(i,k,n &
87              ))
88            rhsumstd(i, k, n) = rhsumstd(i, k, n)/(freq_moynmc(n)-tnondef(i,k &
89              ,n))
90            uvsumstd(i, k, n) = uvsumstd(i, k, n)/(freq_moynmc(n)-tnondef(i,k &
91              ,n))
92            vqsumstd(i, k, n) = vqsumstd(i, k, n)/(freq_moynmc(n)-tnondef(i,k &
93              ,n))
94            vtsumstd(i, k, n) = vtsumstd(i, k, n)/(freq_moynmc(n)-tnondef(i,k &
95              ,n))
96            wqsumstd(i, k, n) = wqsumstd(i, k, n)/(freq_moynmc(n)-tnondef(i,k &
97              ,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)/(freq_moynmc(n)-tnondef(i,k &
101              ,n))
102            u2sumstd(i, k, n) = u2sumstd(i, k, n)/(freq_moynmc(n)-tnondef(i,k &
103              ,n))
104            v2sumstd(i, k, n) = v2sumstd(i, k, n)/(freq_moynmc(n)-tnondef(i,k &
105              ,n))
106            t2sumstd(i, k, n) = t2sumstd(i, k, n)/(freq_moynmc(n)-tnondef(i,k &
107              ,n))
108            o3sumstd(i, k, n) = o3sumstd(i, k, n)/(freq_moynmc(n)-tnondef(i,k &
109              ,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          END IF !tnondef(i,k,n).NE.(freq_moyNMC(n))
132        END DO !i
133      END DO !k
134    END IF !MOD(itap,NINT(freq_outNMC(n)/dtime)).EQ.0
135
136  END DO !n
137
138  RETURN
139END SUBROUTINE moy_undefstd
Note: See TracBrowser for help on using the repository browser.