-
Notifications
You must be signed in to change notification settings - Fork 0
/
coordtrafo.f90
126 lines (106 loc) · 5.48 KB
/
coordtrafo.f90
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
!***********************************************************************
!* Copyright 2012,2013 *
!* Jerome Brioude, Delia Arnold, Andreas Stohl, Wayne Angevine, *
!* John Burkhart, Massimo Cassiani, Adam Dingwell, Richard C Easter, Sabine Eckhardt,*
!* Stephanie Evan, Jerome D Fast, Don Morton, Ignacio Pisso, *
!* Petra Seibert, Gerard Wotawa, Caroline Forster, Harald Sodemann, *
!* *
!* This file is part of FLEXPART WRF *
!* *
!* FLEXPART is free software: you can redistribute it and/or modify *
!* it under the terms of the GNU General Public License as published by*
!* the Free Software Foundation, either version 3 of the License, or *
!* (at your option) any later version. *
!* *
!* FLEXPART is distributed in the hope that it will be useful, *
!* but WITHOUT ANY WARRANTY; without even the implied warranty of *
!* MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
!* GNU General Public License for more details. *
!* *
!* You should have received a copy of the GNU General Public License *
!* along with FLEXPART. If not, see <http://www.gnu.org/licenses/>. *
!***********************************************************************
subroutine coordtrafo
!**********************************************************************
! *
! Note: This is the FLEXPART_WRF version of subroutine coordtrafo. *
! *
! FLEXPART MODEL SUBROUTINE COORDTRAFO *
! *
!**********************************************************************
! *
! AUTHOR: G. WOTAWA *
! DATE: 1994-02-07 *
! LAST UPDATE: 1996-05-18 A. STOHL *
! *
! Dec 2005, R. Easter - changed names of "*lon0*" & "*lat0*" variables*
! July 2012, J Brioude: modification following flexpart 9 *
!**********************************************************************
! *
! DESCRIPTION: This subroutine transforms x and y coordinates of *
! particle release points to grid coordinates. *
! *
!**********************************************************************
use point_mod
use par_mod
use com_mod
implicit none
integer :: i,j,k
if (numpoint.eq.0) goto 30
! TRANSFORM X- AND Y- COORDINATES OF STARTING POINTS TO GRID COORDINATES
!***********************************************************************
do i=1,numpoint
xpoint1(i)=(xpoint1(i)-xmet0)/dx
xpoint2(i)=(xpoint2(i)-xmet0)/dx
ypoint1(i)=(ypoint1(i)-ymet0)/dy
ypoint2(i)=(ypoint2(i)-ymet0)/dy
end do
15 continue
! CHECK IF RELEASE POINTS ARE WITHIN DOMAIN
!******************************************
do i=1,numpoint
if ((ypoint1(i).lt.1.e-6).or.(ypoint1(i).ge.real(nymin1)-1.e-6) &
.or.(ypoint2(i).lt.1.e-6).or.(ypoint2(i).ge.real(nymin1)-1.e-6) &
.or.(xpoint1(i).lt.1.e-6).or.(xpoint1(i).ge.real(nxmin1)-1.e-6) &
.or.(xpoint2(i).lt.1.e-6).or.(xpoint2(i).ge.real(nxmin1)-1.e-6)) then
! if ((ypoint1(i).lt.1.e-6).or.(ypoint1(i).ge.real(nymin1)-1.e-6) &
! .or.(ypoint2(i).lt.1.e-6).or.(ypoint2(i).ge.real(nymin1)-1.e-6) &
! .or.((.not.xglobal).and.((xpoint1(i).lt.1.e-6).or. &
! (xpoint1(i).ge.real(nxmin1)-1.e-6).or.(xpoint2(i).lt.1.e-6).or. &
! (xpoint2(i).ge.real(nxmin1)-1.e-6)))) then
write(*,*) ' NOTICE: RELEASE POINT OUT OF DOMAIN DETECTED.'
write(*,*) ' IT IS REMOVED NOW ... '
write(*,*) ' DEBUG INFO:'
write(*,*) ' Point number:',i
write(*,*) ' xpoint1,ypoint1:',xpoint1,ypoint1
write(*,*) ' xpoint2,ypoint2:',xpoint2,ypoint2
! write(*,*) ' COMMENT: ',compoint(i)
if (i.lt.numpoint) then
do j=i+1,numpoint
xpoint1(j-1)=xpoint1(j)
ypoint1(j-1)=ypoint1(j)
xpoint2(j-1)=xpoint2(j)
ypoint2(j-1)=ypoint2(j)
zpoint1(j-1)=zpoint1(j)
zpoint2(j-1)=zpoint2(j)
npart(j-1)=npart(j)
kindz(j-1)=kindz(j)
ireleasestart(j-1)=ireleasestart(j)
ireleaseend(j-1)=ireleaseend(j)
if (j.le.2000) compoint(j-1)=compoint(j)
do k=1,nspec
xmass(j-1,k)=xmass(j,k)
end do
enddo
endif
numpoint=numpoint-1
if (numpoint.gt.0) goto 15
endif
end do
30 if (numpoint.eq.0) then
write(*,*) ' FLEXPART MODEL SUBROUTINE COORDTRAFO: ERROR ! '
write(*,*) ' NO PARTICLE RELEASES ARE DEFINED!'
write(*,*) ' CHECK FILE RELEASES...'
stop
endif
end subroutine coordtrafo