Changeset 27 for trunk/src/readreceptors.f90
- Timestamp:
- Jul 3, 2014, 2:55:50 PM (10 years ago)
- File:
-
- 1 edited
Legend:
- Unmodified
- Added
- Removed
-
trunk/src/readreceptors.f90
r4 r27 27 27 ! * 28 28 ! Author: A. Stohl * 29 ! *30 29 ! 1 August 1996 * 30 ! HSO, 14 August 2013 31 ! Added optional namelist input 31 32 ! * 32 33 !***************************************************************************** … … 51 52 character(len=16) :: receptor 52 53 54 integer :: readerror 55 real :: lon,lat ! for namelist input, lon/lat are used instead of x,y 56 integer,parameter :: unitreceptorout=2 57 58 ! declare namelist 59 namelist /receptors/ & 60 receptor, lon, lat 61 62 lon=-999.9 53 63 54 64 ! For backward runs, do not allow receptor output. Thus, set number of receptors to zero … … 64 74 !************************************************************ 65 75 66 open(unitreceptor,file=path(1)(1:length(1))//'RECEPTORS', & 67 status='old',err=999) 76 open(unitreceptor,file=path(1)(1:length(1))//'RECEPTORS',form='formatted',status='old',err=999) 68 77 69 call skplin(5,unitreceptor) 78 ! try namelist input 79 read(unitreceptor,receptors,iostat=readerror) 70 80 81 ! prepare namelist output if requested 82 if (nmlout.eqv..true.) then 83 open(unitreceptorout,file=path(2)(1:length(2))//'RECEPTORS.namelist',access='append',status='new',err=1000) 84 endif 71 85 72 ! Read the names and coordinates of the receptors 73 !************************************************ 86 if ((lon.lt.-900).or.(readerror.ne.0)) then 74 87 75 j=0 76 100 j=j+1 88 close(unitreceptor) 89 open(unitreceptor,file=path(1)(1:length(1))//'RECEPTORS',status='old',err=999) 90 call skplin(5,unitreceptor) 91 92 ! Read the names and coordinates of the receptors 93 !************************************************ 94 95 j=0 96 100 j=j+1 77 97 read(unitreceptor,*,end=99) 78 98 read(unitreceptor,*,end=99) … … 89 109 endif 90 110 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 #### '111 write(*,*) ' #### FLEXPART MODEL ERROR! TOO MANY RECEPTOR #### ' 112 write(*,*) ' #### POINTS ARE GIVEN. #### ' 113 write(*,*) ' #### MAXIMUM NUMBER IS ',maxreceptor,' #### ' 114 write(*,*) ' #### PLEASE MAKE CHANGES IN FILE RECEPTORS #### ' 95 115 endif 96 116 receptorname(j)=receptor … … 100 120 ym=r_earth*dy/180.*pi 101 121 receptorarea(j)=xm*ym 122 123 ! write receptors file in namelist format to output directory if requested 124 if (nmlout.eqv..true.) then 125 lon=x 126 lat=y 127 write(unitreceptorout,nml=receptors) 128 endif 129 102 130 goto 100 103 131 104 99 numreceptor=j-1 105 close(unitreceptor) 132 99 numreceptor=j-1 133 134 else ! continue with namelist input 135 136 j=0 137 do while (readerror.eq.0) 138 j=j+1 139 lon=-999.9 140 read(unitreceptor,receptors,iostat=readerror) 141 if ((lon.lt.-900).or.(readerror.ne.0)) then 142 readerror=1 143 else 144 if (j.gt.maxreceptor) then 145 write(*,*) ' #### FLEXPART MODEL ERROR! TOO MANY RECEPTOR #### ' 146 write(*,*) ' #### POINTS ARE GIVEN. #### ' 147 write(*,*) ' #### MAXIMUM NUMBER IS ',maxreceptor,' #### ' 148 write(*,*) ' #### PLEASE MAKE CHANGES IN FILE RECEPTORS #### ' 149 endif 150 receptorname(j)=receptor 151 xreceptor(j)=(lon-xlon0)/dx ! transform to grid coordinates 152 yreceptor(j)=(lat-ylat0)/dy 153 xm=r_earth*cos(lat*pi/180.)*dx/180.*pi 154 ym=r_earth*dy/180.*pi 155 receptorarea(j)=xm*ym 156 endif 157 158 ! write receptors file in namelist format to output directory if requested 159 if (nmlout.eqv..true.) then 160 write(unitreceptorout,nml=receptors) 161 endif 162 163 end do 164 numreceptor=j-1 165 close(unitreceptor) 166 167 endif 168 169 if (nmlout.eqv..true.) then 170 close(unitreceptorout) 171 endif 172 106 173 return 107 174 108 175 109 999 176 999 write(*,*) ' #### FLEXPART MODEL ERROR! FILE "RECEPTORS" #### ' 110 177 write(*,*) ' #### CANNOT BE OPENED IN THE DIRECTORY #### ' 111 178 write(*,'(a)') path(1)(1:length(1)) 112 179 stop 113 180 181 1000 write(*,*) ' #### FLEXPART MODEL ERROR! FILE "RECEPTORS" #### ' 182 write(*,*) ' #### CANNOT BE OPENED IN THE DIRECTORY #### ' 183 write(*,'(a)') path(2)(1:length(2)) 184 stop 185 114 186 end subroutine readreceptors
Note: See TracChangeset
for help on using the changeset viewer.