/[lmdze]/trunk/dyn3d/dynetat0.f
ViewVC logotype

Diff of /trunk/dyn3d/dynetat0.f

Parent Directory Parent Directory | Revision Log Revision Log | View Patch Patch

revision 18 by guez, Thu Aug 7 12:29:13 2008 UTC revision 40 by guez, Tue Feb 22 13:49:36 2011 UTC
# Line 1  Line 1 
1  module dynetat0_m  module dynetat0_m
2    
   ! This module is clean: no C preprocessor directive, no include line.  
   
3    IMPLICIT NONE    IMPLICIT NONE
4    
5  contains    INTEGER day_ini
6    
7    SUBROUTINE dynetat0(vcov, ucov, teta, q, masse, ps, phis, time)  contains
8    
9      ! From dynetat0.F, version 1.2 2004/06/22 11:45:30    SUBROUTINE dynetat0(vcov, ucov, teta, q, masse, ps, phis, time_0)
10    
11      ! Authors:  P. Le Van, L. Fairhead      ! From dynetat0.F, version 1.2, 2004/06/22 11:45:30
12      ! Objet : lecture de l'état initial      ! Authors: P. Le Van, L. Fairhead
13        ! This procedure reads the initial state of the atmosphere.
14    
15      use dimens_m, only: iim, jjm, llm, nqmx      use comconst, only: im, dtvr, jm, lllm
     use comconst, only: im, cpp, dtvr, g, kappa, jm, lllm, omeg, rad  
16      use comvert, only: pa      use comvert, only: pa
     use logic, only: fxyhypb, ysinus  
17      use comgeom, only: rlonu, rlatu, rlonv, rlatv, cu_2d, cv_2d, aire_2d      use comgeom, only: rlonu, rlatu, rlonv, rlatv, cu_2d, cv_2d, aire_2d
18      use serre, only: clon, clat, grossismy, grossismx      use dimens_m, only: iim, jjm, llm, nqmx
     use temps, only: day_ref, day_ini, itau_dyn, annee_ref  
19      use ener, only: etot0, ang0, ptot0, stot0, ztot0      use ener, only: etot0, ang0, ptot0, stot0, ztot0
20      use iniadvtrac_m, only: tname      use iniadvtrac_m, only: tname
21      use netcdf95, only: nf95_open, nf95_inq_varid, handle_err, NF95_CLOSE      use logic, only: fxyhypb, ysinus
22      use netcdf, only: NF90_NOWRITE, NF90_GET_VAR, NF90_NOERR      use serre, only: clon, clat, grossismy, grossismx
23      use numer_rec, only: assert      use netcdf95, only: NF95_GET_VAR, nf95_open, nf95_inq_varid, NF95_CLOSE
24        use netcdf, only: NF90_NOWRITE, NF90_NOERR
25        use nr_util, only: assert
26        use temps, only: day_ref, itau_dyn, annee_ref
27    
28      !   Arguments:      ! Arguments:
29      REAL, intent(out):: vcov(: , :), ucov(:, :), teta(:, :)      REAL, intent(out):: vcov(: , :), ucov(:, :), teta(:, :)
30      REAL, intent(out):: q(:, :, :), masse(:, :)      REAL, intent(out):: q(:, :, :, :) ! (iim + 1, jjm + 1, llm, nqmx)
31      REAL, intent(out):: ps(:) ! in Pa      REAL, intent(out):: masse(:, :)
32      REAL, intent(out):: phis(:, :)      REAL, intent(out):: ps(:, :) ! (iim + 1, jjm + 1) in Pa
33      REAL, intent(out):: time      REAL, intent(out):: phis(:, :) ! (iim + 1, jjm + 1)
34        REAL, intent(out):: time_0
35      !   Variables  
36      INTEGER length, iq      ! Variables
37      PARAMETER (length = 100)      INTEGER iq
38      REAL tab_cntrl(length) ! tableau des parametres du run      INTEGER, PARAMETER:: length = 100
39      INTEGER ierr, nid, nvarid      REAL tab_cntrl(length) ! tableau des paramètres du run
40        INTEGER ierr, ncid, varid
41    
42      !-----------------------------------------------------------------------      !-----------------------------------------------------------------------
43    
44      print *, "Call sequence information: dynetat0"      print *, "Call sequence information: dynetat0"
45    
46      call assert(size(vcov, 1) == (iim + 1) * jjm, "dynetat0 vcov 1")      call assert(size(vcov, 1) == (iim + 1) * jjm, "dynetat0 vcov 1")
47      call assert((/size(ucov, 1), size(teta, 1), size(q, 1), size(masse, 1), &      call assert((/size(ucov, 1), size(teta, 1), size(masse, 1)/) &
48           size(ps)/) == (iim + 1) * (jjm + 1), "dynetat0 (iim + 1) * (jjm + 1)")           == (iim + 1) * (jjm + 1), "dynetat0 (iim + 1) * (jjm + 1)")
49      call assert(shape(phis) == (/iim + 1, jjm + 1/), "dynetat0 phis")      call assert((/size(ps, 1), size(phis, 1), size(q, 1)/) == iim + 1, &
50      call assert((/size(vcov, 2), size(ucov, 2), size(teta, 2), size(q, 2), &           "dynetat0 iim")
51        call assert((/size(ps, 2), size(phis, 2), size(q, 2)/) == jjm + 1, &
52             "dynetat0 jjm")
53        call assert((/size(vcov, 2), size(ucov, 2), size(teta, 2), size(q, 3), &
54           size(masse, 2)/) == llm, "dynetat0 llm")           size(masse, 2)/) == llm, "dynetat0 llm")
55      call assert(size(q, 3) == nqmx, "dynetat0 q 3")      call assert(size(q, 4) == nqmx, "dynetat0 q nqmx")
56    
57      ! Fichier état initial :      ! Fichier état initial :
58      call nf95_open("start.nc", NF90_NOWRITE, nid)      call nf95_open("start.nc", NF90_NOWRITE, ncid)
59    
60      call nf95_inq_varid(nid, "controle", nvarid)      call nf95_inq_varid(ncid, "controle", varid)
61      ierr = NF90_GET_VAR(nid, nvarid, tab_cntrl)      call NF95_GET_VAR(ncid, varid, tab_cntrl)
     call handle_err("dynetat0, controle", ierr, nid)  
   
     im         = int(tab_cntrl(1))  
     jm         = int(tab_cntrl(2))  
     lllm       = int(tab_cntrl(3))  
     day_ref    = int(tab_cntrl(4))  
     annee_ref  = int(tab_cntrl(5))  
     omeg       = tab_cntrl(7)  
     dtvr       = tab_cntrl(12)  
     etot0      = tab_cntrl(13)  
     ptot0      = tab_cntrl(14)  
     ztot0      = tab_cntrl(15)  
     stot0      = tab_cntrl(16)  
     ang0       = tab_cntrl(17)  
     pa         = tab_cntrl(18)  
     clon       = tab_cntrl(20)  
     clat       = tab_cntrl(21)  
     grossismx  = tab_cntrl(22)  
     grossismy  = tab_cntrl(23)  
   
     IF (tab_cntrl(24) == 1.)  THEN  
        fxyhypb  = .TRUE.  
     ELSE  
        fxyhypb = .FALSE.  
        ysinus  = .FALSE.  
        IF (tab_cntrl(27) == 1.) ysinus = .TRUE.  
     ENDIF  
62    
63      day_ini = tab_cntrl(30)      im = int(tab_cntrl(1))
64        jm = int(tab_cntrl(2))
65        lllm = int(tab_cntrl(3))
66        call assert(im == iim, "dynetat0 im iim")
67        call assert(jm == jjm, "dynetat0 jm jjm")
68        call assert(lllm == llm, "dynetat0 lllm llm")
69    
70        day_ref = int(tab_cntrl(4))
71        annee_ref = int(tab_cntrl(5))
72        dtvr = tab_cntrl(12)
73        etot0 = tab_cntrl(13)
74        ptot0 = tab_cntrl(14)
75        ztot0 = tab_cntrl(15)
76        stot0 = tab_cntrl(16)
77        ang0 = tab_cntrl(17)
78        pa = tab_cntrl(18)
79        clon = tab_cntrl(20)
80        clat = tab_cntrl(21)
81        grossismx = tab_cntrl(22)
82        grossismy = tab_cntrl(23)
83        fxyhypb = tab_cntrl(24) == 1.
84        if (.not. fxyhypb) ysinus = tab_cntrl(27) == 1.
85      itau_dyn = tab_cntrl(31)      itau_dyn = tab_cntrl(31)
86    
87      PRINT *, 'rad = ', rad      call NF95_INQ_VARID (ncid, "rlonu", varid)
88      PRINT *, 'omeg = ', omeg      call NF95_GET_VAR(ncid, varid, rlonu)
89      PRINT *, 'g = ', g  
90      PRINT *, 'cpp = ', cpp      call NF95_INQ_VARID (ncid, "rlatu", varid)
91      PRINT *, 'kappa = ', kappa      call NF95_GET_VAR(ncid, varid, rlatu)
92    
93      IF (im /= iim)  THEN      call NF95_INQ_VARID (ncid, "rlonv", varid)
94         PRINT 1, im, iim      call NF95_GET_VAR(ncid, varid, rlonv)
95         STOP 1  
96      ELSE  IF (jm /= jjm)  THEN      call NF95_INQ_VARID (ncid, "rlatv", varid)
97         PRINT 2, jm, jjm      call NF95_GET_VAR(ncid, varid, rlatv)
98         STOP 1  
99      ELSE  IF (lllm /= llm)  THEN      call NF95_INQ_VARID (ncid, "cu", varid)
100         PRINT 3, lllm, llm      call NF95_GET_VAR(ncid, varid, cu_2d)
101         STOP 1  
102      ENDIF      call NF95_INQ_VARID (ncid, "cv", varid)
103        call NF95_GET_VAR(ncid, varid, cv_2d)
104      call NF95_INQ_VARID (nid, "rlonu", nvarid)  
105      ierr = NF90_GET_VAR(nid, nvarid, rlonu)      call NF95_INQ_VARID (ncid, "aire", varid)
106      call handle_err("dynetat0, rlonu", ierr, nid)      call NF95_GET_VAR(ncid, varid, aire_2d)
107    
108      call NF95_INQ_VARID (nid, "rlatu", nvarid)      call NF95_INQ_VARID (ncid, "phisinit", varid)
109      ierr = NF90_GET_VAR(nid, nvarid, rlatu)      call NF95_GET_VAR(ncid, varid, phis)
110      call handle_err("dynetat0, rlatu", ierr, nid)  
111        call NF95_INQ_VARID (ncid, "temps", varid)
112      call NF95_INQ_VARID (nid, "rlonv", nvarid)      call NF95_GET_VAR(ncid, varid, time_0)
113      ierr = NF90_GET_VAR(nid, nvarid, rlonv)  
114      call handle_err("dynetat0, rlonv", ierr, nid)      day_ini = tab_cntrl(30) + INT(time_0)
115        time_0 = time_0 - INT(time_0)
116      call NF95_INQ_VARID (nid, "rlatv", nvarid)      ! {0 <= time0 < 1}
117      ierr = NF90_GET_VAR(nid, nvarid, rlatv)  
118      call handle_err("dynetat0, rlatv", ierr, nid)      call NF95_INQ_VARID (ncid, "ucov", varid)
119        call NF95_GET_VAR(ncid, varid, ucov, count_nc=(/iim + 1, jjm + 1, llm, 1/))
120      call NF95_INQ_VARID (nid, "cu", nvarid)  
121      ierr = NF90_GET_VAR(nid, nvarid, cu_2d)      call NF95_INQ_VARID (ncid, "vcov", varid)
122      call handle_err("dynetat0, cu", ierr, nid)      call NF95_GET_VAR(ncid, varid, vcov, count_nc=(/iim + 1, jjm, llm, 1/))
123    
124      call NF95_INQ_VARID (nid, "cv", nvarid)      call NF95_INQ_VARID (ncid, "teta", varid)
125      ierr = NF90_GET_VAR(nid, nvarid, cv_2d)      call NF95_GET_VAR(ncid, varid, teta, count_nc=(/iim + 1, jjm + 1, llm, 1/))
     call handle_err("dynetat0, cv", ierr, nid)  
   
     call NF95_INQ_VARID (nid, "aire", nvarid)  
     ierr = NF90_GET_VAR(nid, nvarid, aire_2d)  
     call handle_err("dynetat0, aire", ierr, nid)  
   
     call NF95_INQ_VARID (nid, "phisinit", nvarid)  
     ierr = NF90_GET_VAR(nid, nvarid, phis)  
     call handle_err("dynetat0, phisinit", ierr, nid)  
   
     call NF95_INQ_VARID (nid, "temps", nvarid)  
     ierr = NF90_GET_VAR(nid, nvarid, time)  
     call handle_err("dynetat0, temps", ierr, nid)  
   
     call NF95_INQ_VARID (nid, "ucov", nvarid)  
     ierr = NF90_GET_VAR(nid, nvarid, ucov, count=(/iim + 1, jjm + 1, llm/))  
     call handle_err("dynetat0, ucov", ierr, nid)  
   
     call NF95_INQ_VARID (nid, "vcov", nvarid)  
     ierr = NF90_GET_VAR(nid, nvarid, vcov, count=(/iim + 1, jjm, llm/))  
     call handle_err("dynetat0, vcov", ierr, nid)  
   
     call NF95_INQ_VARID (nid, "teta", nvarid)  
     ierr = NF90_GET_VAR(nid, nvarid, teta, count=(/iim + 1, jjm + 1, llm/))  
     call handle_err("dynetat0, teta", ierr, nid)  
126    
127      DO iq = 1, nqmx      DO iq = 1, nqmx
128         call NF95_INQ_VARID(nid, tname(iq), nvarid, ierr)         call NF95_INQ_VARID(ncid, tname(iq), varid, ierr)
129         IF (ierr  /=  NF90_NOERR) THEN         IF (ierr /= NF90_NOERR) THEN
130            PRINT *, 'dynetat0: le champ "' // tname(iq) // '" est absent, ' // &            PRINT *, 'dynetat0: "' // tname(iq) // '" not found, ' // &
131                 "il est donc initialisé à zéro."                 "setting it to zero..."
132            q(:, :, iq) = 0.            q(:, :, :, iq) = 0.
133         ELSE         ELSE
134            ierr = NF90_GET_VAR(nid, nvarid, q(:, :, iq), &            call NF95_GET_VAR(ncid, varid, q(:, :, :, iq))
                count=(/iim + 1, jjm + 1, llm/))  
           call handle_err("dynetat0, " // tname(iq), ierr, nid)  
135         ENDIF         ENDIF
136      ENDDO      ENDDO
137    
138      call NF95_INQ_VARID (nid, "masse", nvarid)      call NF95_INQ_VARID (ncid, "masse", varid)
139      ierr = NF90_GET_VAR(nid, nvarid, masse, count=(/iim + 1, jjm + 1, llm/))      call NF95_GET_VAR(ncid, varid, masse, count_nc=(/iim + 1, jjm + 1, llm, 1/))
140      call handle_err("dynetat0, masse", ierr, nid)  
141        call NF95_INQ_VARID (ncid, "ps", varid)
142      call NF95_INQ_VARID (nid, "ps", nvarid)      call NF95_GET_VAR(ncid, varid, ps)
143      ierr = NF90_GET_VAR(nid, nvarid, ps, count=(/iim + 1, jjm + 1/))  
144      call handle_err("dynetat0, ps", ierr, nid)      call NF95_CLOSE(ncid)
   
     call NF95_CLOSE(nid)  
   
     day_ini=day_ini+INT(time)  
     time=time-INT(time)  
   
 1   FORMAT(//10x, 'la valeur de im =', i4, 2x, &  
          'lue sur le fichier de demarrage est differente de la valeur ' &  
          // 'parametree iim =', i4//)  
 2   FORMAT(//10x, 'la valeur de jm =', i4, 2x, &  
          'lue sur le fichier de demarrage est differente de la valeur ' &  
          // 'parametree jjm =', i4//)  
 3   FORMAT(//10x, 'la valeur de lmax =', i4, 2x, &  
          'lue sur le fichier demarrage est differente de la valeur ' &  
          // 'parametree llm =', i4//)  
145    
146    END SUBROUTINE dynetat0    END SUBROUTINE dynetat0
147    

Legend:
Removed from v.18  
changed lines
  Added in v.40

  ViewVC Help
Powered by ViewVC 1.1.21