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