source: trunk/LMDZ.MARS/libf/phymars/writeg1d.F @ 171

Last change on this file since 171 was 141, checked in by emillour, 14 years ago

Mars: EM

minor bug correction in writeg1d (JYC)
add "-check" to debug option in makegcm_ifort

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