source: LMDZ6/branches/contrails/libf/dyn3d_common/inigrads.f90 @ 5445

Last change on this file since 5445 was 5297, checked in by abarral, 8 weeks 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: 1.7 KB
Line 
1!
2! $Header$
3!
4subroutine inigrads(if,im &
5        ,x,fx,xmin,xmax,jm,y,ymin,ymax,fy,lm,z,fz &
6        ,dt,file,titlel)
7  USE gradsdef_mod_h
8  implicit none
9
10  integer :: if,im,jm,lm,i,j,l
11  real :: x(im),y(jm),z(lm),fx,fy,fz,dt
12  real :: xmin,xmax,ymin,ymax
13
14  character(len=*),intent(in) :: file
15  character(len=*),intent(in) :: titlel
16
17  ! data unit/66,32,34,36,38,40,42,44,46,48/
18  integer :: nf
19  save nf
20  data nf/0/
21
22  unit(1)=66
23  unit(2)=32
24  unit(3)=34
25  unit(4)=36
26  unit(5)=38
27  unit(6)=40
28  unit(7)=42
29  unit(8)=44
30  unit(9)=46
31
32  if (if.le.nf) stop 'verifier les appels a inigrads'
33
34  print*,'Entree dans inigrads'
35
36  nf=if
37  title(if)=titlel
38  ivar(if)=0
39
40  fichier(if)=trim(file)
41
42  firsttime(if)=.true.
43  dtime(if)=dt
44
45  iid(if)=1
46  ifd(if)=im
47  imd(if)=im
48  do i=1,im
49     xd(i,if)=x(i)*fx
50     if(xd(i,if).lt.xmin) iid(if)=i+1
51     if(xd(i,if).le.xmax) ifd(if)=i
52  enddo
53  print*,'On stoke du point ',iid(if),'  a ',ifd(if),' en x'
54
55  jid(if)=1
56  jfd(if)=jm
57  jmd(if)=jm
58  do j=1,jm
59     yd(j,if)=y(j)*fy
60     if(yd(j,if).gt.ymax) jid(if)=j+1
61     if(yd(j,if).ge.ymin) jfd(if)=j
62  enddo
63  print*,'On stoke du point ',jid(if),'  a ',jfd(if),' en y'
64
65  print*,'Open de dat'
66  print*,'file=',file
67  print*,'fichier(if)=',fichier(if)
68
69  print*,4*(ifd(if)-iid(if))*(jfd(if)-jid(if))
70  print*,trim(file)//'.dat'
71
72  OPEN (unit(if)+1,FILE=trim(file)//'.dat' &
73        ,FORM='unformatted', &
74        ACCESS='direct' &
75        ,RECL=4*(ifd(if)-iid(if)+1)*(jfd(if)-jid(if)+1))
76
77  print*,'Open de dat ok'
78
79  lmd(if)=lm
80  do l=1,lm
81     zd(l,if)=z(l)*fz
82  enddo
83
84  irec(if)=0
85
86  print*,if,imd(if),jmd(if),lmd(if)
87  print*,'if,imd(if),jmd(if),lmd(if)'
88
89  return
90end subroutine inigrads
Note: See TracBrowser for help on using the repository browser.