forked from dnowacki-usgs/swanmod
-
Notifications
You must be signed in to change notification settings - Fork 0
/
SwanFindObstacles.f90
145 lines (145 loc) · 4.94 KB
/
SwanFindObstacles.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
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
subroutine SwanFindObstacles ( cross )
!
! --|-----------------------------------------------------------|--
! | Delft University of Technology |
! | Faculty of Civil Engineering and Geosciences |
! | Environmental Fluid Mechanics Section |
! | P.O. Box 5048, 2600 GA Delft, The Netherlands |
! | |
! | Programmer: Marcel Zijlema |
! --|-----------------------------------------------------------|--
!
!
! SWAN (Simulating WAves Nearshore); a third generation wave model
! Copyright (C) 1993-2015 Delft University of Technology
!
! This program 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 2 of
! the License, or (at your option) any later version.
!
! This program 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.
!
! A copy of the GNU General Public License is available at
! http://www.gnu.org/copyleft/gpl.html#SEC3
! or by writing to the Free Software Foundation, Inc.,
! 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
!
!
! Authors
!
! 40.80: Marcel Zijlema
!
! Updates
!
! 40.80, March 2008: New subroutine
!
! Purpose
!
! Search for obstacles in computational grid and store them
!
! Method
!
! for each face an obstacle is found when they crossed each other
!
! Modules used
!
use ocpcomm4
use swcomm3
use m_obsta
use SwanGriddata
use SwanGridobjects
!
implicit none
!
! Argument variables
!
integer, dimension(nfaces), intent(out) :: cross ! contains sequence number of obstacles for each face
! where they crossing or zero if no crossing
!
! Local variables
!
integer, save :: ient = 0 ! number of entries in this subroutine
integer :: iface ! loop counter over faces
integer :: j ! loop counter
integer :: k ! loop counter
integer :: vb ! vertex of begin of present face
integer :: ve ! vertex of end of present face
!
real, dimension(2) :: xobs ! x-coordinate of obstacle point
real, dimension(2) :: xv ! x-coordinate of vertex of face
real, dimension(2) :: yobs ! y-coordinate of obstacle point
real, dimension(2) :: yv ! y-coordinate of vertex of face
!
logical :: SwanCrossObstacle ! indicate whether a face cross an obstacle
!
type(OBSTDAT), pointer :: cobst ! pointer to a considered obstacle
!
type(facetype), dimension(:), pointer :: face ! datastructure for faces with their attributes
type(verttype), dimension(:), pointer :: vert ! datastructure for vertices with their attributes
!
! Structure
!
! Description of the pseudo code
!
! Source text
!
if (ltrace) call strace (ient,'SwanFindObstacles')
!
! point to vertex and face objects
!
vert => gridobject%vert_grid
face => gridobject%face_grid
!
! initialize array cross
!
cross = 0
!
! go to first obstacle
!
cobst => FOBSTAC
!
do j = 1, NUMOBS
!
xobs(1) = cobst%XCRP(1)
yobs(1) = cobst%YCRP(1)
!
do k = 2, cobst%NCRPTS ! number of corner points of considered obstacle
!
xobs(2) = cobst%XCRP(k)
yobs(2) = cobst%YCRP(k)
!
! loop over faces
!
faceloop : do iface = 1, nfaces
!
if ( face(iface)%atti(FMARKER) == 1 ) cycle faceloop ! boundary face
!
vb = face(iface)%atti(FACEV1)
ve = face(iface)%atti(FACEV2)
!
xv(1) = vert(vb)%attr(VERTX)
yv(1) = vert(vb)%attr(VERTY)
xv(2) = vert(ve)%attr(VERTX)
yv(2) = vert(ve)%attr(VERTY)
!
if ( SwanCrossObstacle( xv, yv, xobs, yobs ) ) cross(iface) = j
!
enddo faceloop
!
xobs(1) = xobs(2)
yobs(1) = yobs(2)
!
enddo
!
! go to next obstacle, if present
!
if (.not.associated(cobst%NEXTOBST)) exit
cobst => cobst%NEXTOBST
!
enddo
!
end subroutine SwanFindObstacles