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

Last change on this file since 5258 was 5246, checked in by abarral, 6 weeks ago

Convert fixed-form to free-form sources .F -> .{f,F}90
(WIP: some .F remain, will be handled in subsequent commits)

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