/[lmdze]/trunk/libf/dyn3d/dynetat0.f90
ViewVC logotype

Diff of /trunk/libf/dyn3d/dynetat0.f90

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

revision 37 by guez, Thu Dec 2 17:11:04 2010 UTC revision 38 by guez, Thu Jan 6 17:52:19 2011 UTC
# Line 8  contains Line 8  contains
8    
9    SUBROUTINE dynetat0(vcov, ucov, teta, q, masse, ps, phis, time_0)    SUBROUTINE dynetat0(vcov, ucov, teta, q, masse, ps, phis, time_0)
10    
11      ! From dynetat0.F, version 1.2 2004/06/22 11:45:30      ! From dynetat0.F, version 1.2, 2004/06/22 11:45:30
12      ! Authors:  P. Le Van, L. Fairhead      ! Authors: P. Le Van, L. Fairhead
13      ! Objet : lecture de l'état initial      ! Objet : lecture de l'état initial
14    
15      use dimens_m, only: iim, jjm, llm, nqmx      use comconst, only: im, dtvr, jm, lllm, omeg
     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, 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 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      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(:, :, :), masse(:, :)
31      REAL, intent(out):: ps(:) ! in Pa      REAL, intent(out):: ps(:) ! in Pa
32      REAL, intent(out):: phis(:, :)      REAL, intent(out):: phis(:, :)
33      REAL, intent(out):: time_0      REAL, intent(out):: time_0
34    
35      !   Variables      ! Variables
36      INTEGER length, iq      INTEGER iq
37      PARAMETER (length = 100)      INTEGER, PARAMETER:: length = 100
38      REAL tab_cntrl(length) ! tableau des parametres du run      REAL tab_cntrl(length) ! tableau des paramètres du run
39      INTEGER ierr, ncid, nvarid      INTEGER ierr, ncid, varid
40    
41      !-----------------------------------------------------------------------      !-----------------------------------------------------------------------
42    
# Line 53  contains Line 53  contains
53      ! Fichier état initial :      ! Fichier état initial :
54      call nf95_open("start.nc", NF90_NOWRITE, ncid)      call nf95_open("start.nc", NF90_NOWRITE, ncid)
55    
56      call nf95_inq_varid(ncid, "controle", nvarid)      call nf95_inq_varid(ncid, "controle", varid)
57      ierr = NF90_GET_VAR(ncid, nvarid, tab_cntrl)      call NF95_GET_VAR(ncid, varid, tab_cntrl)
     call handle_err("dynetat0, controle", ierr, ncid)  
   
     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  
58    
59      day_ini = tab_cntrl(30)      im = int(tab_cntrl(1))
60      itau_dyn = tab_cntrl(31)      jm = int(tab_cntrl(2))
61        lllm = int(tab_cntrl(3))
62        call assert(im == iim, "dynetat0 im iim")
63        call assert(jm == jjm, "dynetat0 jm jjm")
64        call assert(lllm == llm, "dynetat0 lllm llm")
65    
66        day_ref = int(tab_cntrl(4))
67        annee_ref = int(tab_cntrl(5))
68    
69      PRINT *, 'rad = ', rad      omeg = tab_cntrl(7)
70      PRINT *, 'omeg = ', omeg      PRINT *, 'omeg = ', omeg
71      PRINT *, 'g = ', g  
72      PRINT *, 'cpp = ', cpp      dtvr = tab_cntrl(12)
73      PRINT *, 'kappa = ', kappa      etot0 = tab_cntrl(13)
74        ptot0 = tab_cntrl(14)
75      IF (im /= iim)  THEN      ztot0 = tab_cntrl(15)
76         PRINT 1, im, iim      stot0 = tab_cntrl(16)
77         STOP 1      ang0 = tab_cntrl(17)
78      ELSE  IF (jm /= jjm)  THEN      pa = tab_cntrl(18)
79         PRINT 2, jm, jjm      clon = tab_cntrl(20)
80         STOP 1      clat = tab_cntrl(21)
81      ELSE  IF (lllm /= llm)  THEN      grossismx = tab_cntrl(22)
82         PRINT 3, lllm, llm      grossismy = tab_cntrl(23)
83         STOP 1      fxyhypb = tab_cntrl(24) == 1.
84      ENDIF      if (.not. fxyhypb) ysinus = tab_cntrl(27) == 1.
85        itau_dyn = tab_cntrl(31)
86      call NF95_INQ_VARID (ncid, "rlonu", nvarid)  
87      ierr = NF90_GET_VAR(ncid, nvarid, rlonu)      call NF95_INQ_VARID (ncid, "rlonu", varid)
88      call handle_err("dynetat0, rlonu", ierr, ncid)      call NF95_GET_VAR(ncid, varid, rlonu)
89    
90      call NF95_INQ_VARID (ncid, "rlatu", nvarid)      call NF95_INQ_VARID (ncid, "rlatu", varid)
91      ierr = NF90_GET_VAR(ncid, nvarid, rlatu)      call NF95_GET_VAR(ncid, varid, rlatu)
92      call handle_err("dynetat0, rlatu", ierr, ncid)  
93        call NF95_INQ_VARID (ncid, "rlonv", varid)
94      call NF95_INQ_VARID (ncid, "rlonv", nvarid)      call NF95_GET_VAR(ncid, varid, rlonv)
95      ierr = NF90_GET_VAR(ncid, nvarid, rlonv)  
96      call handle_err("dynetat0, rlonv", ierr, ncid)      call NF95_INQ_VARID (ncid, "rlatv", varid)
97        call NF95_GET_VAR(ncid, varid, rlatv)
98      call NF95_INQ_VARID (ncid, "rlatv", nvarid)  
99      ierr = NF90_GET_VAR(ncid, nvarid, rlatv)      call NF95_INQ_VARID (ncid, "cu", varid)
100      call handle_err("dynetat0, rlatv", ierr, ncid)      call NF95_GET_VAR(ncid, varid, cu_2d)
101    
102      call NF95_INQ_VARID (ncid, "cu", nvarid)      call NF95_INQ_VARID (ncid, "cv", varid)
103      ierr = NF90_GET_VAR(ncid, nvarid, cu_2d)      call NF95_GET_VAR(ncid, varid, cv_2d)
104      call handle_err("dynetat0, cu", ierr, ncid)  
105        call NF95_INQ_VARID (ncid, "aire", varid)
106      call NF95_INQ_VARID (ncid, "cv", nvarid)      call NF95_GET_VAR(ncid, varid, aire_2d)
107      ierr = NF90_GET_VAR(ncid, nvarid, cv_2d)  
108      call handle_err("dynetat0, cv", ierr, ncid)      call NF95_INQ_VARID (ncid, "phisinit", varid)
109        call NF95_GET_VAR(ncid, varid, phis)
110      call NF95_INQ_VARID (ncid, "aire", nvarid)  
111      ierr = NF90_GET_VAR(ncid, nvarid, aire_2d)      call NF95_INQ_VARID (ncid, "temps", varid)
112      call handle_err("dynetat0, aire", ierr, ncid)      call NF95_GET_VAR(ncid, varid, time_0)
113    
114      call NF95_INQ_VARID (ncid, "phisinit", nvarid)      day_ini = tab_cntrl(30) + INT(time_0)
115      ierr = NF90_GET_VAR(ncid, nvarid, phis)      time_0 = time_0 - INT(time_0)
116      call handle_err("dynetat0, phisinit", ierr, ncid)      ! {0 <= time0 < 1}
117    
118      call NF95_INQ_VARID (ncid, "temps", nvarid)      call NF95_INQ_VARID (ncid, "ucov", varid)
119      ierr = NF90_GET_VAR(ncid, nvarid, time_0)      call NF95_GET_VAR(ncid, varid, ucov, count_nc=(/iim + 1, jjm + 1, llm/))
120      call handle_err("dynetat0, temps", ierr, ncid)  
121        call NF95_INQ_VARID (ncid, "vcov", varid)
122      call NF95_INQ_VARID (ncid, "ucov", nvarid)      call NF95_GET_VAR(ncid, varid, vcov, count_nc=(/iim + 1, jjm, llm/))
123      ierr = NF90_GET_VAR(ncid, nvarid, ucov, count=(/iim + 1, jjm + 1, llm/))  
124      call handle_err("dynetat0, ucov", ierr, ncid)      call NF95_INQ_VARID (ncid, "teta", varid)
125        call NF95_GET_VAR(ncid, varid, teta, count_nc=(/iim + 1, jjm + 1, llm/))
     call NF95_INQ_VARID (ncid, "vcov", nvarid)  
     ierr = NF90_GET_VAR(ncid, nvarid, vcov, count=(/iim + 1, jjm, llm/))  
     call handle_err("dynetat0, vcov", ierr, ncid)  
   
     call NF95_INQ_VARID (ncid, "teta", nvarid)  
     ierr = NF90_GET_VAR(ncid, nvarid, teta, count=(/iim + 1, jjm + 1, llm/))  
     call handle_err("dynetat0, teta", ierr, ncid)  
126    
127      DO iq = 1, nqmx      DO iq = 1, nqmx
128         call NF95_INQ_VARID(ncid, 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(ncid, nvarid, q(:, :, iq), &            call NF95_GET_VAR(ncid, varid, q(:, :, iq), &
135                 count=(/iim + 1, jjm + 1, llm/))                 count_nc=(/iim + 1, jjm + 1, llm/))
           call handle_err("dynetat0, " // tname(iq), ierr, ncid)  
136         ENDIF         ENDIF
137      ENDDO      ENDDO
138    
139      call NF95_INQ_VARID (ncid, "masse", nvarid)      call NF95_INQ_VARID (ncid, "masse", varid)
140      ierr = NF90_GET_VAR(ncid, nvarid, masse, count=(/iim + 1, jjm + 1, llm/))      call NF95_GET_VAR(ncid, varid, masse, count_nc=(/iim + 1, jjm + 1, llm/))
     call handle_err("dynetat0, masse", ierr, ncid)  
   
     call NF95_INQ_VARID (ncid, "ps", nvarid)  
     ierr = NF90_GET_VAR(ncid, nvarid, ps, count=(/iim + 1, jjm + 1/))  
     call handle_err("dynetat0, ps", ierr, ncid)  
141    
142      call NF95_CLOSE(ncid)      call NF95_INQ_VARID (ncid, "ps", varid)
143        call NF95_GET_VAR(ncid, varid, ps, count_nc=(/iim + 1, jjm + 1/))
144    
145      day_ini=day_ini+INT(time_0)      call NF95_CLOSE(ncid)
     time_0=time_0-INT(time_0)  
     ! {0 <= time0 < 1}  
   
 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//)  
146    
147    END SUBROUTINE dynetat0    END SUBROUTINE dynetat0
148    

Legend:
Removed from v.37  
changed lines
  Added in v.38

  ViewVC Help
Powered by ViewVC 1.1.21