source: trunk/LMDZ.GENERIC/libf/phystd/writeg1d.F @ 2613

Last change on this file since 2613 was 1524, checked in by emillour, 9 years ago

All GCMS:
More updates to enforce dynamics/physics separation:

get rid of references to "temps_mod" from physics packages;
make a "time_phylmdz_mod.F90" module to store that
information and fill it via "iniphysiq".

EM

File size: 6.8 KB
RevLine 
[135]1      SUBROUTINE writeg1d(ngrid,nx,x,nom,titre)
[1397]2      USE comg1d_mod, ONLY: g1d_nomfich,g1d_premier,g1d_unitfich,
3     &  g1d_irec,g1d_nvar,g1d_nomvar,g1d_titrevar,g1d_dimvar,g1d_nlayer
[135]4      IMPLICIT NONE
5
6c.......................................................................
7c
8c  ecriture de x pour GRADS-1D
9c
10c  in :
11c         * ngrid      ---> pour controler que l'on est bien en 1D
12c         * nx         ---> taille du vecteur a stocker
13c                             "1" pour une variable de surface
14c                             "nlayer" pour une variable de centre de couche
15c                             "nlayer+1" pour une variable d'interface
16c         * x          ---> variable a stocker
17c         * nom        ---> nom "pour grads"
18c         * titre      ---> titre "pour grads"
19c
20c.......................................................................
21c
22
23c
24c.......................................................................
25c  declaration des arguments
26c
27      INTEGER ngrid,nx,i
28      REAL*4 xr4(1000)
29      REAL x(nx)
30      CHARACTER*(*) nom,titre
31c
32c  declaration des arguments
33c.......................................................................
34c  declaration des variables locales
35c
36      INTEGER ilayer,ivar
37      LOGICAL test
38c
39c  declaration des variables locales
40c.......................................................................
41c  controle 1D
42c
43c     print*,'ngrid=',ngrid
44      IF (ngrid.NE.1) return
45c
46c  controle 1D
47c.......................................................................
48c  copy => force en reel 4 pour l'ecriture dans grads1d.dat
49
50      do i=1,nx
51        xr4(i) = x(i)
52      enddo
53
54c  copy => force en reel 4 pour l'ecriture dans grads1d.dat
55c.......................................................................
56c  ouverture du fichier au premier appel
57
58
59      g1d_nomfich='g1d.dat'
60
61      IF (g1d_premier) THEN
62        OPEN (g1d_unitfich,FILE=g1d_nomfich
63     &       ,FORM='unformatted',ACCESS='direct',RECL=4)
64        g1d_irec=0
65        g1d_nvar=0
66        g1d_premier=.false.
67      ENDIF
68
69c  ouverture du fichier au premier appel
70c.......................................................................
71c  pour l'ecriture du fichier ctl
72
73      test=.true.
74      DO ivar=1,g1d_nvar
75        IF (nom.EQ.g1d_nomvar(ivar)) test=.false.
76        IF (nx .GT. 1000) then
77          print*,'ERROR:  nx > 1000 dans writeg1d.F'
78          print*,'Changer la dimension de xr4'
79          call exit(1)
80        ENDIF
81      ENDDO
82      IF (test) THEN
83        g1d_nvar=g1d_nvar+1
84        g1d_nomvar(g1d_nvar)=nom
85        g1d_titrevar(g1d_nvar)=titre
86        IF (nx.EQ.1) THEN
87           g1d_dimvar(g1d_nvar)=0
88        ELSEIF (nx.EQ.g1d_nlayer) THEN
89           g1d_dimvar(g1d_nvar)=g1d_nlayer
90        ELSEIF (nx.EQ.g1d_nlayer+1) THEN
91           g1d_dimvar(g1d_nvar)=g1d_nlayer+1
92        ELSE
93           PRINT *,'._. probleme de dimension dans GRADS-1D ._.'
94           print*,'NX = ',nx
95           print*,'g1d_nlayer = ',g1d_nlayer
96        ENDIF
97      ENDIF
98c
99c  pour l'ecriture du fichier ctl
100c.......................................................................
101c  ecriture
102c
103      IF (nx.EQ.1) THEN
104        g1d_irec=g1d_irec+1
105        WRITE(g1d_unitfich,REC=g1d_irec) xr4(1)
106      ELSE
107        DO ilayer=1,g1d_nlayer
108          g1d_irec=g1d_irec+1
109          WRITE(g1d_unitfich,REC=g1d_irec) xr4(ilayer)
110        ENDDO
111      ENDIF
112c
113c  ecriture
114c.......................................................................
115c
11610001 CONTINUE
117c
118c.......................................................................
119c
120      RETURN
121      END
122
123
124
125
126
127c *********************************************************************
128c *********************************************************************
129
130      SUBROUTINE endg1d(ngrid,nlayer,zlayer,ndt)
[1524]131      USE time_phylmdz_mod, ONLY: dtphys, daysec
[1397]132      USE comg1d_mod, ONLY: g1d_nomfich,g1d_unitfich,g1d_nvar,
133     &  g1d_nomvar,g1d_titrevar,g1d_dimvar,g1d_nlayer,g1d_unitctl,
134     &  g1d_nomctl,saveG1D
[135]135      IMPLICIT NONE
136c.......................................................................
137c
138c  ecriture du fichier de controle pour GRADS-1D
139c
140c  in :
141c         * ngrid      ---> pour controler que l'on est bien en 1D
142c         * nlayer     ---> nombre de couches
143c         * zlayer     ---> altitude au centre de chaque couche (km)
144c         * ndt        ---> nombre de pas de temps
145c
146c.......................................................................
147c
148
149
150c
151c.......................................................................
152c  declaration des arguments
153c
154      INTEGER ngrid,nlayer
155      REAL zlayer(nlayer)
156      INTEGER ndt
157c
158c  declaration des arguments
159c.......................................................................
160c  declaration des variables locales
161c
162      INTEGER ivar,ilayer
163c
164
165
166!      integer saveG1D
167
168c  declaration des variables locales
169c.......................................................................
170c  contole 1D
171c
172      IF (ngrid.NE.1) GOTO 10001
173c
174c  contole 1D
175c.......................................................................
176c
177      IF (nlayer.ne.g1d_nlayer)
178     &PRINT *,'._. probleme de dimension dans GRADS-1D (endg1d.F) '
179c
180c.......................................................................
181c
182      CLOSE (g1d_unitfich)
183c
184c.......................................................................
185
186
187      OPEN (g1d_unitctl,FILE=g1d_nomctl,FORM='formatted',RECL=4*100)
188      WRITE (g1d_unitctl,'(a4,2x,a1,a20)') 'DSET','^',g1d_nomfich
189      WRITE (g1d_unitctl,'(a5,2x,a20)') 'UNDEF ','1.E+30'
190      WRITE (g1d_unitctl,'(a11)') 'FORMAT YREV'
191      WRITE (g1d_unitctl,'(a5,2x,a30)') 'TITLE ','champs 1D'
192      WRITE (g1d_unitctl,'(a5,i4,a20)') 'XDEF ',1,' LINEAR 0 1'
193      WRITE (g1d_unitctl,'(a5,i4,a20)') 'YDEF ',1,' LINEAR 0 1'
194      WRITE (g1d_unitctl,'(a5,i4,a20)') 'ZDEF ',g1d_nlayer,' LEVELS'
195      WRITE (g1d_unitctl,'(5(1x,f13.5))')
196     &      (zlayer(ilayer),ilayer=1,g1d_nlayer)
197
198c     Writing true timestep in g1d.ctl (in planet "minutes"= sol/(60*24)
199!      ivar =min( max(int(1440.*dtphys/daysec +0.5),1) , 99)   
200!      WRITE (g1d_unitctl,'(a4,2x,i10,a19,i2,a2)')
201!     &      'TDEF ',ndt,' LINEAR 01JAN2000 ', ivar,'MN '
202
203      ivar =min( max(int(1440.*dtphys*saveG1D/daysec +0.5),1) , 99)
204      ! not sure ivar is right, but it doesnt matter
205      WRITE (g1d_unitctl,'(a4,2x,i10,a19,i2,a2)')
206     &      'TDEF ',ndt/saveG1D,' LINEAR 01JAN2000 ', ivar,'MN '
207
208      WRITE (g1d_unitctl,'(a5,i5)') 'VARS ',g1d_nvar
209      DO ivar=1,g1d_nvar
210      WRITE (g1d_unitctl,'(a9,3x,i4,i3,1x,a39)')
211     &       g1d_nomvar(ivar),g1d_dimvar(ivar),99,g1d_titrevar(ivar)
212      ENDDO
213      WRITE (g1d_unitctl,'(a7)') 'ENDVARS'
214      CLOSE (g1d_unitctl)
215c
216c.......................................................................
217c
21810001 CONTINUE
219c
220c.......................................................................
221c
222      RETURN
223      END
Note: See TracBrowser for help on using the repository browser.