[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 * |
---|
[3405994] | 12 | ! HSO, 14 August 2013: Added optional namelist input * |
---|
| 13 | ! PS, 2020-05-28: correct bug in nml input, code cosmetics * |
---|
[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 | |
---|
[3405994] | 36 | integer :: ierr |
---|
[b4d29ce] | 37 | real :: lon,lat ! for namelist input, lon/lat are used instead of x,y |
---|
[3405994] | 38 | integer,parameter :: iunitreceptorout=2 |
---|
[b4d29ce] | 39 | |
---|
| 40 | ! declare namelist |
---|
[3405994] | 41 | namelist /receptors/ receptor, lon, lat |
---|
[e200b7a] | 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 | |
---|
[3405994] | 51 | ! prepare namelist output if requested |
---|
| 52 | if (nmlout .and. lroot) & |
---|
| 53 | open(iunitreceptorout,file=path(2)(1:length(2))//'RECEPTORS.namelist', & |
---|
| 54 | status='replace',err=1000) |
---|
[e200b7a] | 55 | |
---|
[3405994] | 56 | ! Open the RECEPTORS file and read output grid specifications |
---|
| 57 | !************************************************************ |
---|
| 58 | ! try namelist input |
---|
| 59 | open(unitreceptor,file=path(1)(1:length(1))//'RECEPTORS',status='old',err=999) |
---|
| 60 | read(unitreceptor,receptors,iostat=ierr) |
---|
| 61 | close(unitreceptor) |
---|
[e200b7a] | 62 | |
---|
[3405994] | 63 | open(unitreceptor,file=path(1)(1:length(1))//'RECEPTORS') |
---|
[e200b7a] | 64 | |
---|
[3405994] | 65 | if (ierr.ne.0) then ! not namelist |
---|
| 66 | |
---|
[b4d29ce] | 67 | call skplin(5,unitreceptor) |
---|
[e200b7a] | 68 | |
---|
[b4d29ce] | 69 | ! Read the names and coordinates of the receptors |
---|
| 70 | !************************************************ |
---|
| 71 | |
---|
| 72 | j=0 |
---|
| 73 | 100 j=j+1 |
---|
[3405994] | 74 | read(unitreceptor,*,end=99) |
---|
| 75 | read(unitreceptor,*,end=99) |
---|
| 76 | read(unitreceptor,*,end=99) |
---|
| 77 | read(unitreceptor,'(4x,a16)',end=99) receptor |
---|
| 78 | call skplin(3,unitreceptor) |
---|
| 79 | read(unitreceptor,'(4x,f11.4)',end=99) x |
---|
| 80 | call skplin(3,unitreceptor) |
---|
| 81 | read(unitreceptor,'(4x,f11.4)',end=99) y |
---|
| 82 | if (x.eq.0. .and. y.eq.0. .and. receptor.eq.' ') then |
---|
| 83 | j=j-1 |
---|
| 84 | goto 100 |
---|
| 85 | endif |
---|
| 86 | if (j.gt.maxreceptor) goto 998 ! ERROR - STOP |
---|
| 87 | receptorname(j)=receptor |
---|
| 88 | xreceptor(j)=(x-xlon0)/dx ! transform to grid coordinates |
---|
| 89 | yreceptor(j)=(y-ylat0)/dy |
---|
| 90 | xm=r_earth*cos(y*pi/180.)*dx/180.*pi |
---|
| 91 | ym=r_earth*dy/180.*pi |
---|
| 92 | receptorarea(j)=xm*ym |
---|
| 93 | |
---|
| 94 | ! write receptors file in namelist format to output directory if requested |
---|
| 95 | if (nmlout .and. lroot) then |
---|
| 96 | lon=x |
---|
| 97 | lat=y |
---|
| 98 | write(iunitreceptorout,nml=receptors) |
---|
| 99 | endif |
---|
| 100 | |
---|
| 101 | goto 100 ! read next |
---|
| 102 | |
---|
| 103 | 99 numreceptor=j-1 ! read all |
---|
[b4d29ce] | 104 | |
---|
| 105 | else ! continue with namelist input |
---|
| 106 | |
---|
| 107 | j=0 |
---|
[3405994] | 108 | do while (ierr.eq.0) |
---|
[b4d29ce] | 109 | j=j+1 |
---|
[3405994] | 110 | read(unitreceptor,receptors,iostat=ierr) |
---|
| 111 | if (ierr.eq.0) then |
---|
| 112 | if (j.gt.maxreceptor) goto 998 ! ERROR - STOP |
---|
[b4d29ce] | 113 | receptorname(j)=receptor |
---|
| 114 | xreceptor(j)=(lon-xlon0)/dx ! transform to grid coordinates |
---|
| 115 | yreceptor(j)=(lat-ylat0)/dy |
---|
| 116 | xm=r_earth*cos(lat*pi/180.)*dx/180.*pi |
---|
| 117 | ym=r_earth*dy/180.*pi |
---|
| 118 | receptorarea(j)=xm*ym |
---|
| 119 | ! write receptors file in namelist format to output directory if requested |
---|
[3405994] | 120 | if (nmlout.and.lroot) & |
---|
| 121 | write(iunitreceptorout,nml=receptors) |
---|
[b4d29ce] | 122 | endif |
---|
| 123 | end do |
---|
| 124 | numreceptor=j-1 |
---|
| 125 | close(unitreceptor) |
---|
| 126 | |
---|
[3405994] | 127 | endif ! end reading nml input |
---|
[b4d29ce] | 128 | |
---|
[3405994] | 129 | if (nmlout.and.lroot) & |
---|
| 130 | close(iunitreceptorout) |
---|
[b4d29ce] | 131 | |
---|
[e200b7a] | 132 | return |
---|
| 133 | |
---|
[3405994] | 134 | 998 continue |
---|
| 135 | write(*,*) ' #### FLEXPART MODEL ERROR! TOO MANY RECEPTOR #### ' |
---|
| 136 | write(*,*) ' #### POINTS ARE GIVEN. #### ' |
---|
| 137 | write(*,*) ' #### MAXIMUM NUMBER IS ',maxreceptor,' #### ' |
---|
| 138 | write(*,*) ' #### PLEASE MAKE CHANGES IN FILE RECEPTORS #### ' |
---|
| 139 | stop |
---|
[e200b7a] | 140 | |
---|
[3405994] | 141 | 999 continue |
---|
| 142 | write(*,*) ' #### FLEXPART MODEL ERROR! FILE "RECEPTORS" #### ' |
---|
[e200b7a] | 143 | write(*,*) ' #### CANNOT BE OPENED IN THE DIRECTORY #### ' |
---|
| 144 | write(*,'(a)') path(1)(1:length(1)) |
---|
| 145 | stop |
---|
| 146 | |
---|
[3405994] | 147 | 1000 continue |
---|
| 148 | write(*,*) ' #### FLEXPART MODEL ERROR! FILE "RECEPTORS" #### ' |
---|
[b4d29ce] | 149 | write(*,*) ' #### CANNOT BE OPENED IN THE DIRECTORY #### ' |
---|
| 150 | write(*,'(a)') path(2)(1:length(2)) |
---|
| 151 | stop |
---|
| 152 | |
---|
[e200b7a] | 153 | end subroutine readreceptors |
---|