source: LMDZ6/trunk/libf/dyn3d/wrgrads.F90 @ 5300

Last change on this file since 5300 was 5297, checked in by abarral, 2 days ago

Turn gradsdef.h coefils.h into a module

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