source: LMDZ5/branches/LMDZ5V1.0-dev/libf/dyn3dpar/wrgrads.F @ 2865

Last change on this file since 2865 was 774, checked in by Laurent Fairhead, 17 years ago

Suite du merge entre la version et la HEAD: quelques modifications de
Yann sur le

LF

  • Property svn:eol-style set to native
  • Property svn:keywords set to Author Date Id Revision
File size: 3.5 KB
Line 
1!
2! $Header$
3!
4      subroutine wrgrads(if,nl,field,name,titlevar)
5      implicit none
6
7c   Declarations
8c    if indice du fichier
9c    nl nombre de couches
10c    field   champ
11c    name    petit nom
12c    titlevar   Titre
13
14#include "gradsdef.h"
15
16c   arguments
17      integer if,nl
18      real field(imx*jmx*lmx)
19      character*10 name,file
20      character*10 titlevar
21
22c   local
23
24      integer im,jm,lm,i,j,l,lnblnk,iv,iii,iji,iif,ijf
25
26      logical writectl
27
28
29      writectl=.false.
30
31      print*,if,iid(if),jid(if),ifd(if),jfd(if)
32      iii=iid(if)
33      iji=jid(if)
34      iif=ifd(if)
35      ijf=jfd(if)
36      im=iif-iii+1
37      jm=ijf-iji+1
38      lm=lmd(if)
39
40      print*,'im,jm,lm,name,firsttime(if)'
41      print*,im,jm,lm,name,firsttime(if)
42
43      if(firsttime(if)) then
44         if(name.eq.var(1,if)) then
45            firsttime(if)=.false.
46            ivar(if)=1
47         print*,'fin de l initialiation de l ecriture du fichier'
48         print*,file
49           print*,'fichier no: ',if
50           print*,'unit ',unit(if)
51           print*,'nvar  ',nvar(if)
52           print*,'vars ',(var(iv,if),iv=1,nvar(if))
53         else
54            ivar(if)=ivar(if)+1
55            nvar(if)=ivar(if)
56            var(ivar(if),if)=name
57            tvar(ivar(if),if)=titlevar(1:lnblnk(titlevar))
58            nld(ivar(if),if)=nl
59            print*,'initialisation ecriture de ',var(ivar(if),if)
60            print*,'if ivar(if) nld ',if,ivar(if),nld(ivar(if),if)
61         endif
62         writectl=.true.
63         itime(if)=1
64      else
65         ivar(if)=mod(ivar(if),nvar(if))+1
66         if (ivar(if).eq.nvar(if)) then
67            writectl=.true.
68            itime(if)=itime(if)+1
69         endif
70
71         if(var(ivar(if),if).ne.name) then
72           print*,'Il faut stoker la meme succession de champs a chaque'
73           print*,'pas de temps'
74           print*,'fichier no: ',if
75           print*,'unit ',unit(if)
76           print*,'nvar  ',nvar(if)
77           print*,'vars ',(var(iv,if),iv=1,nvar(if))
78
79           stop
80         endif
81      endif
82
83      print*,'ivar(if),nvar(if),var(ivar(if),if),writectl'
84      print*,ivar(if),nvar(if),var(ivar(if),if),writectl
85      do l=1,nl
86         irec(if)=irec(if)+1
87c        print*,'Ecrit rec=',irec(if),iii,iif,iji,ijf,
88c    s (l-1)*imd(if)*jmd(if)+(iji-1)*imd(if)+iii
89c    s ,(l-1)*imd(if)*jmd(if)+(ijf-1)*imd(if)+iif
90         write(unit(if)+1,rec=irec(if))
91     s   ((field((l-1)*imd(if)*jmd(if)+(j-1)*imd(if)+i)
92     s   ,i=iii,iif),j=iji,ijf)
93      enddo
94      if (writectl) then
95
96      file=fichier(if)
97c   WARNING! on reecrase le fichier .ctl a chaque ecriture
98      open(unit(if),file=file(1:lnblnk(file))//'.ctl'
99     &         ,form='formatted',status='unknown')
100      write(unit(if),'(a5,1x,a40)')
101     &       'DSET ','^'//file(1:lnblnk(file))//'.dat'
102
103      write(unit(if),'(a12)') 'UNDEF 1.0E30'
104      write(unit(if),'(a5,1x,a40)') 'TITLE ',title(if)
105      call formcoord(unit(if),im,xd(iii,if),1.,.false.,'XDEF')
106      call formcoord(unit(if),jm,yd(iji,if),1.,.true.,'YDEF')
107      call formcoord(unit(if),lm,zd(1,if),1.,.false.,'ZDEF')
108      write(unit(if),'(a4,i10,a30)')
109     &       'TDEF ',itime(if),' LINEAR 02JAN1987 1MO '
110      write(unit(if),'(a4,2x,i5)') 'VARS',nvar(if)
111      do iv=1,nvar(if)
112c        print*,'if,var(iv,if),nld(iv,if),nld(iv,if)-1/nld(iv,if)'
113c        print*,if,var(iv,if),nld(iv,if),nld(iv,if)-1/nld(iv,if)
114         write(unit(if),1000) var(iv,if),nld(iv,if)-1/nld(iv,if)
115     &     ,99,tvar(iv,if)
116      enddo
117      write(unit(if),'(a7)') 'ENDVARS'
118c
1191000  format(a5,3x,i4,i3,1x,a39)
120
121      close(unit(if))
122
123      endif ! writectl
124
125      return
126
127      END
128
Note: See TracBrowser for help on using the repository browser.