source: flexpart.git/src/readreceptors.f90 @ 3481cc1

10.4.1_peseiGFS_025bugfixes+enhancementsdevrelease-10release-10.4.1scaling-bug
Last change on this file since 3481cc1 was 3481cc1, checked in by Ignacio Pisso <ip@…>, 4 years ago

move license from headers to a different file

  • Property mode set to 100644
File size: 5.6 KB
Line 
1subroutine readreceptors
2
3  !*****************************************************************************
4  !                                                                            *
5  !     This routine reads the user specifications for the receptor points.    *
6  !                                                                            *
7  !     Author: A. Stohl                                                       *
8  !     1 August 1996                                                          *
9  !     HSO, 14 August 2013
10  !     Added optional namelist input
11  !                                                                            *
12  !*****************************************************************************
13  !                                                                            *
14  ! Variables:                                                                 *
15  ! receptorarea(maxreceptor)  area of dx*dy at location of receptor           *
16  ! receptorname(maxreceptor)  names of receptors                              *
17  ! xreceptor,yreceptor  coordinates of receptor points                        *
18  !                                                                            *
19  ! Constants:                                                                 *
20  ! unitreceptor         unit connected to file RECEPTORS                      *
21  !                                                                            *
22  !*****************************************************************************
23
24  use par_mod
25  use com_mod
26
27  implicit none
28
29  integer :: j
30  real :: x,y,xm,ym
31  character(len=16) :: receptor
32
33  integer :: readerror
34  real :: lon,lat   ! for namelist input, lon/lat are used instead of x,y
35  integer,parameter :: unitreceptorout=2
36
37  ! declare namelist
38  namelist /receptors/ &
39    receptor, lon, lat
40
41  lon=-999.9
42
43  ! For backward runs, do not allow receptor output. Thus, set number of receptors to zero
44  !*****************************************************************************
45
46  if (ldirect.lt.0) then
47    numreceptor=0
48    return
49  endif
50
51
52  ! Open the RECEPTORS file and read output grid specifications
53  !************************************************************
54
55  open(unitreceptor,file=path(1)(1:length(1))//'RECEPTORS',form='formatted',status='old',err=999)
56
57  ! try namelist input
58  read(unitreceptor,receptors,iostat=readerror)
59
60  ! prepare namelist output if requested
61  if (nmlout.and.lroot) then
62    open(unitreceptorout,file=path(2)(1:length(2))//'RECEPTORS.namelist',&
63         &access='append',status='replace',err=1000)
64  endif
65
66  if ((lon.lt.-900).or.(readerror.ne.0)) then
67
68    close(unitreceptor)
69    open(unitreceptor,file=path(1)(1:length(1))//'RECEPTORS',status='old',err=999)
70    call skplin(5,unitreceptor)
71
72    ! Read the names and coordinates of the receptors
73    !************************************************
74
75    j=0
76100 j=j+1
77    read(unitreceptor,*,end=99)
78    read(unitreceptor,*,end=99)
79    read(unitreceptor,*,end=99)
80    read(unitreceptor,'(4x,a16)',end=99) receptor
81    call skplin(3,unitreceptor)
82    read(unitreceptor,'(4x,f11.4)',end=99) x
83    call skplin(3,unitreceptor)
84    read(unitreceptor,'(4x,f11.4)',end=99) y
85    if ((x.eq.0.).and.(y.eq.0.).and. &
86         (receptor.eq.'                ')) then
87      j=j-1
88      goto 100
89    endif
90    if (j.gt.maxreceptor) then
91      write(*,*) ' #### FLEXPART MODEL ERROR! TOO MANY RECEPTOR #### '
92      write(*,*) ' #### POINTS ARE GIVEN.                       #### '
93      write(*,*) ' #### MAXIMUM NUMBER IS ',maxreceptor,'       #### '
94      write(*,*) ' #### PLEASE MAKE CHANGES IN FILE RECEPTORS   #### '
95    endif
96    receptorname(j)=receptor
97    xreceptor(j)=(x-xlon0)/dx       ! transform to grid coordinates
98    yreceptor(j)=(y-ylat0)/dy
99    xm=r_earth*cos(y*pi/180.)*dx/180.*pi
100    ym=r_earth*dy/180.*pi
101    receptorarea(j)=xm*ym
102
103    ! write receptors file in namelist format to output directory if requested
104    if (nmlout.and.lroot) then
105      lon=x
106      lat=y
107      write(unitreceptorout,nml=receptors)
108    endif
109
110    goto 100
111
11299  numreceptor=j-1
113
114  else ! continue with namelist input
115
116    j=0
117    do while (readerror.eq.0)
118      j=j+1
119      lon=-999.9
120      read(unitreceptor,receptors,iostat=readerror)
121      if ((lon.lt.-900).or.(readerror.ne.0)) then
122        readerror=1
123      else
124        if (j.gt.maxreceptor) then
125          write(*,*) ' #### FLEXPART MODEL ERROR! TOO MANY RECEPTOR #### '
126          write(*,*) ' #### POINTS ARE GIVEN.                       #### '
127          write(*,*) ' #### MAXIMUM NUMBER IS ',maxreceptor,'       #### '
128          write(*,*) ' #### PLEASE MAKE CHANGES IN FILE RECEPTORS   #### '
129        endif
130        receptorname(j)=receptor
131        xreceptor(j)=(lon-xlon0)/dx       ! transform to grid coordinates
132        yreceptor(j)=(lat-ylat0)/dy
133        xm=r_earth*cos(lat*pi/180.)*dx/180.*pi
134        ym=r_earth*dy/180.*pi
135        receptorarea(j)=xm*ym
136      endif
137
138      ! write receptors file in namelist format to output directory if requested
139      if (nmlout.and.lroot) then
140        write(unitreceptorout,nml=receptors)
141      endif
142
143    end do
144    numreceptor=j-1
145    close(unitreceptor)
146
147  endif
148
149  if (nmlout.and.lroot) then
150    close(unitreceptorout)
151  endif
152
153  return
154
155
156999 write(*,*) ' #### FLEXPART MODEL ERROR! FILE "RECEPTORS"  #### '
157  write(*,*) ' #### CANNOT BE OPENED IN THE DIRECTORY       #### '
158  write(*,'(a)') path(1)(1:length(1))
159  stop
160
1611000 write(*,*) ' #### FLEXPART MODEL ERROR! FILE "RECEPTORS"    #### '
162  write(*,*) ' #### CANNOT BE OPENED IN THE DIRECTORY       #### '
163  write(*,'(a)') path(2)(1:length(2))
164  stop
165
166end subroutine readreceptors
Note: See TracBrowser for help on using the repository browser.
hosted by ZAMG