source: LMDZ5/trunk/libf/phylmd/undefSTD.F90 @ 3469

Last change on this file since 3469 was 2346, checked in by Ehouarn Millour, 9 years ago

Physics/dynamics separation:

  • remove all references to dimensions.h from physics. nbp_lon (==iim) , nbp_lat (==jjm+1) and nbp_lev (==llm) from mod_grid_phy_lmdz should be used instead.
  • added module regular_lonlat_mod in phy_common to store information about the global (lon-lat) grid cell boundaries and centers.

EM

  • 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: 3.2 KB
RevLine 
[1992]1
[1403]2! $Id: undefSTD.F90 2346 2015-08-21 15:13:46Z ymeurdesoif $
[1418]3
[1992]4SUBROUTINE undefstd(itap, read_climoz)
5  USE netcdf
6  USE dimphy
[2271]7#ifdef CPP_IOIPSL
8  USE phys_state_var_mod
9#endif
[2312]10#ifdef CPP_XIOS
11  USE wxios, ONLY: missing_val
12#endif
[2271]13
[1992]14  IMPLICIT NONE
15  include "clesphys.h"
[2312]16#ifndef CPP_XIOS
[2271]17  REAL :: missing_val
18#endif
[1992]19
20  ! ====================================================================
21
22  ! I. Musat : 09.2004
23
24  ! Calcul * du nombre de pas de temps (FLOAT(ecrit_XXX)-tnondef))
25  ! ou la variable tlevSTD est bien definie (.NE.missing_val),
26  ! et
27  ! * de la somme de tlevSTD => tsumSTD
28
29  ! nout=1 !var. journaliere "day" moyenne sur tous les pas de temps
30  ! ! de la physique
31  ! nout=2 !var. mensuelle "mth" moyennee sur tous les pas de temps
32  ! ! de la physique
33  ! nout=3 !var. mensuelle "NMC" moyennee toutes les ecrit_hf
34
35
36  ! NB: mettre "inst(X)" dans le write_hist*NMC.h !
37  ! ====================================================================
38
39  ! ym#include "dimphy.h"
40  ! variables Input
41
42  ! INTEGER nlevSTD, klevSTD, itap
43  ! PARAMETER(klevSTD=17)
44  INTEGER itap
45  ! REAL dtime
46
47  ! variables locales
48  ! INTEGER i, k, nout, n
49  ! PARAMETER(nout=3) !nout=1 : day; =2 : mth; =3 : NMC
50  INTEGER i, k, n
51  ! REAL freq_calNMC(nout)
52  INTEGER read_climoz
53
54  ! variables Output
55  ! REAL tlevSTD(klon,klevSTD), tsumSTD(klon,klevSTD,nout)
56  ! LOGICAL oknondef(klon,klevSTD,nout)
57  ! REAL tnondef(klon,klevSTD,nout)
58
[2271]59! REAL missing_val
[1992]60
[2271]61! missing_val = nf90_fill_real
62#ifndef CPP_XIOS
63      missing_val=missing_val_nf90
64#endif
[1992]65
66  DO n = 1, nout
67
68
69    ! calcul variables tous les freq_calNMC(n)/dtime pas de temps
70    ! de la physique
71
72    IF (mod(itap,nint(freq_calnmc(n)/dtime))==0) THEN
73      DO k = 1, nlevstd
74        DO i = 1, klon
75          IF (tlevstd(i,k)==missing_val) THEN
76            ! IF(oknondef(i,k,n)) THEN
77            tnondef(i, k, n) = tnondef(i, k, n) + 1.
78            ! ENDIF !oknondef(i,k)
79
80          ELSE IF (tlevstd(i,k)/=missing_val) THEN
81            tsumstd(i, k, n) = tsumstd(i, k, n) + tlevstd(i, k)
82            usumstd(i, k, n) = usumstd(i, k, n) + ulevstd(i, k)
83            vsumstd(i, k, n) = vsumstd(i, k, n) + vlevstd(i, k)
84            wsumstd(i, k, n) = wsumstd(i, k, n) + wlevstd(i, k)
85            phisumstd(i, k, n) = phisumstd(i, k, n) + philevstd(i, k)
86            qsumstd(i, k, n) = qsumstd(i, k, n) + qlevstd(i, k)
87            rhsumstd(i, k, n) = rhsumstd(i, k, n) + rhlevstd(i, k)
88            uvsumstd(i, k, n) = uvsumstd(i, k, n) + uvstd(i, k)
89            vqsumstd(i, k, n) = vqsumstd(i, k, n) + vqstd(i, k)
90            vtsumstd(i, k, n) = vtsumstd(i, k, n) + vtstd(i, k)
91            wqsumstd(i, k, n) = wqsumstd(i, k, n) + wqstd(i, k)
92            vphisumstd(i, k, n) = vphisumstd(i, k, n) + vphistd(i, k)
93            wtsumstd(i, k, n) = wtsumstd(i, k, n) + wtstd(i, k)
94            u2sumstd(i, k, n) = u2sumstd(i, k, n) + u2std(i, k)
95            v2sumstd(i, k, n) = v2sumstd(i, k, n) + v2std(i, k)
96            t2sumstd(i, k, n) = t2sumstd(i, k, n) + t2std(i, k)
97            o3sumstd(i, k, n) = o3sumstd(i, k, n) + o3std(i, k)
98            IF (read_climoz==2) o3daysumstd(i, k, n) = o3daysumstd(i, k, n) + &
99              o3daystd(i, k)
100
101          END IF
102        END DO !i
103      END DO !k
104
105    END IF !MOD(itap,NINT(freq_calNMC(n)/dtime)).EQ.0
106
107  END DO !n
108
109  RETURN
110END SUBROUTINE undefstd
Note: See TracBrowser for help on using the repository browser.