geometry.F90 coverage: 44.44 %func 25.50 %block
1) module Geometry_module
2)
3) use PFLOTRAN_Constants_module
4)
5) implicit none
6)
7) private
8)
9) #include "petsc/finclude/petscsys.h"
10)
11) type, public :: point3d_type
12) PetscReal :: x
13) PetscReal :: y
14) PetscReal :: z
15) end type point3d_type
16)
17) type, public :: polygonal_volume_type
18) type(point3d_type), pointer :: xy_coordinates(:)
19) type(point3d_type), pointer :: xz_coordinates(:)
20) type(point3d_type), pointer :: yz_coordinates(:)
21) end type polygonal_volume_type
22)
23) interface GeometryPointInPolygon
24) module procedure GeometryPointInPolygon1
25) module procedure GeometryPointInPolygon2
26) end interface GeometryPointInPolygon
27)
28) public :: GeometryCreatePolygonalVolume, &
29) GeometryReadCoordinates, &
30) GeometryReadCoordinate, &
31) GeometryPointInPolygonalVolume, &
32) GeometryCopyCoordinates, &
33) GeometryDestroyPolygonalVolume, &
34) GeometryPointInPolygon
35)
36) contains
37)
38) ! ************************************************************************** !
39)
40) function GeometryCreatePolygonalVolume()
41) !
42) ! Creates a polygonal volume. I.e. a volume
43) ! defined by 3 polygons, on each plane of the
44) ! principle coordinate system x,y,z
45) !
46) ! Author: Glenn Hammond
47) ! Date: 01/12/12
48) !
49)
50) implicit none
51)
52) type(polygonal_volume_type), pointer :: GeometryCreatePolygonalVolume
53)
54) type(polygonal_volume_type), pointer :: polygonal_volume
55)
56) allocate(polygonal_volume)
57) nullify(polygonal_volume%xy_coordinates)
58) nullify(polygonal_volume%xz_coordinates)
59) nullify(polygonal_volume%yz_coordinates)
60)
61) GeometryCreatePolygonalVolume => polygonal_volume
62)
63) end function GeometryCreatePolygonalVolume
64)
65) ! ************************************************************************** !
66)
67) subroutine GeometryReadCoordinates(input,option,region_name,coordinates)
68) !
69) ! Reads a list of coordinates
70) !
71) ! Author: Glenn Hammond
72) ! Date: 01/12/12
73) !
74)
75) use Input_Aux_module
76) use Option_module
77)
78) implicit none
79)
80) type(input_type), pointer :: input
81) type(option_type) :: option
82) character(len=MAXWORDLENGTH) :: region_name
83) type(point3d_type), pointer :: coordinates(:)
84)
85) PetscInt :: icount
86) PetscInt, parameter :: max_num_coordinates = 100
87) type(point3d_type) :: temp_coordinates(max_num_coordinates)
88)
89) icount = 0
90) do
91) call InputReadPflotranString(input,option)
92) call InputReadStringErrorMsg(input,option,'REGION')
93) if (InputCheckExit(input,option)) exit
94) icount = icount + 1
95) if (icount > max_num_coordinates) then
96) write(option%io_buffer, &
97) '(''Number of coordinates in region '',a, &
98) &'' exceeds limit of '',i3)') region_name, max_num_coordinates
99) option%io_buffer = trim(option%io_buffer) // &
100) ' Increase size of PetscInt, parameter :: max_num_coordinates ' // &
101) ' in RegionReadCoordinates()'
102) call printErrMsg(option)
103) endif
104) call GeometryReadCoordinate(input,option,temp_coordinates(icount),'REGION')
105) enddo
106) allocate(coordinates(icount))
107) do icount = 1, size(coordinates)
108) coordinates(icount)%x = temp_coordinates(icount)%x
109) coordinates(icount)%y = temp_coordinates(icount)%y
110) coordinates(icount)%z = temp_coordinates(icount)%z
111) enddo
112)
113) end subroutine GeometryReadCoordinates
114)
115) ! ************************************************************************** !
116)
117) subroutine GeometryReadCoordinate(input,option,coordinate,error_string)
118) !
119) ! Reads a list of coordinates
120) !
121) ! Author: Glenn Hammond
122) ! Date: 01/12/12
123) !
124)
125) use Input_Aux_module
126) use Option_module
127)
128) implicit none
129)
130) type(input_type) :: input
131) type(option_type) :: option
132) character(len=*) :: error_string
133) type(point3d_type) :: coordinate
134)
135) call InputReadDouble(input,option,coordinate%x)
136) call InputErrorMsg(input,option,'x-coordinate',error_string)
137) call InputReadDouble(input,option,coordinate%y)
138) call InputErrorMsg(input,option,'y-coordinate',error_string)
139) call InputReadDouble(input,option,coordinate%z)
140) call InputErrorMsg(input,option,'z-coordinate',error_string)
141)
142) end subroutine GeometryReadCoordinate
143)
144) ! ************************************************************************** !
145)
146) function GeometryPointInPolygonalVolume(x,y,z,polygonal_volume,option)
147) !
148) ! Determines whether a point in xyz space is
149) ! within a polygonal volume defined by polygons
150) !
151) ! Author: Glenn Hammond
152) ! Date: 01/12/12
153) !
154)
155) use Option_module
156)
157) implicit none
158)
159) PetscReal :: x, y, z
160) type(polygonal_volume_type) :: polygonal_volume
161) type(option_type) :: option
162)
163) PetscBool :: xy, xz, yz
164) PetscBool :: GeometryPointInPolygonalVolume
165)
166) GeometryPointInPolygonalVolume = PETSC_FALSE
167)
168) ! if a polygon is not defined in a particular direction, it is assumed that
169) ! the point is within the "infinite" polygon
170)
171) ! XY plane
172) if (associated(polygonal_volume%xy_coordinates)) then
173) if (.not.GeometryPointInPolygon(x,y,z,Z_DIRECTION, &
174) polygonal_volume%xy_coordinates)) then
175) return
176) endif
177) endif
178)
179) ! XZ plane
180) if (associated(polygonal_volume%xz_coordinates)) then
181) if (.not.GeometryPointInPolygon(x,y,z,Y_DIRECTION, &
182) polygonal_volume%xz_coordinates)) then
183) return
184) endif
185) endif
186)
187) ! YZ plane
188) if (associated(polygonal_volume%yz_coordinates)) then
189) if (.not.GeometryPointInPolygon(x,y,z,X_DIRECTION, &
190) polygonal_volume%yz_coordinates)) then
191) return
192) endif
193) endif
194)
195) ! if point is not within any of polygons above, the function will return
196) ! prior to this point.
197) GeometryPointInPolygonalVolume = PETSC_TRUE
198)
199) end function GeometryPointInPolygonalVolume
200)
201) ! ************************************************************************** !
202)
203) function GeometryPointInPolygon1(x,y,z,axis,coordinates)
204) !
205) ! Determines whether a point in xyz space is within
206) ! a 2d polygon based on coordinate object
207) !
208) ! Author: Glenn Hammond
209) ! Date: 01/12/12
210) !
211)
212) implicit none
213)
214) PetscReal :: x, y, z
215) PetscInt :: axis
216) type(point3d_type) :: coordinates(:)
217)
218) PetscInt, parameter :: max_num_coordinates = 100
219) PetscReal :: xx, yy
220) PetscInt :: i, num_coordinates
221) PetscReal :: xx_array(max_num_coordinates), yy_array(max_num_coordinates)
222)
223) PetscBool :: GeometryPointInPolygon1
224)
225) GeometryPointInPolygon1 = PETSC_FALSE
226)
227) num_coordinates = size(coordinates)
228) select case(axis)
229) case(Z_DIRECTION)
230) xx = x
231) yy = y
232) do i = 1, num_coordinates
233) xx_array(i) = coordinates(i)%x
234) yy_array(i) = coordinates(i)%y
235) enddo
236) case(Y_DIRECTION)
237) xx = x
238) yy = z
239) do i = 1, num_coordinates
240) xx_array(i) = coordinates(i)%x
241) yy_array(i) = coordinates(i)%z
242) enddo
243) case(X_DIRECTION)
244) xx = y
245) yy = z
246) do i = 1, num_coordinates
247) xx_array(i) = coordinates(i)%y
248) yy_array(i) = coordinates(i)%z
249) enddo
250) end select
251) if (num_coordinates == 2) then
252) GeometryPointInPolygon1 = &
253) GeometryPointInRectangle(xx,yy,xx_array,yy_array)
254) else
255) GeometryPointInPolygon1 = &
256) GeometryPointInPolygon2(xx,yy,xx_array,yy_array,num_coordinates)
257) endif
258)
259) end function GeometryPointInPolygon1
260)
261) ! ************************************************************************** !
262)
263) function GeometryPointInPolygon2(x,y,x_array,y_array,num_coordinates)
264) !
265) ! Determines whether a point in xy space is within
266) ! a polygon
267) !
268) ! Author: Glenn Hammond
269) ! Date: 01/12/12
270) !
271)
272) implicit none
273)
274) PetscReal :: x, y
275) PetscReal :: x_array(:), y_array(:)
276) PetscInt :: num_coordinates
277)
278) PetscBool :: GeometryPointInPolygon2
279) PetscInt :: i, j
280)
281) GeometryPointInPolygon2 = PETSC_FALSE
282)
283) j = 1
284) do i = 1, num_coordinates
285) j = i + 1
286) if (i == num_coordinates) j = 1
287) if ((y_array(i) < y .and. y_array(j) >= y) .or. &
288) (y_array(j) < y .and. y_array(i) >= y)) then
289) if ((x_array(i) + &
290) (y-y_array(i))/(y_array(j)-y_array(i))*(x_array(j)-x_array(i))) &
291) < x) then
292) GeometryPointInPolygon2 = .not.GeometryPointInPolygon2
293) endif
294) endif
295) enddo
296)
297) end function GeometryPointInPolygon2
298)
299) ! ************************************************************************** !
300)
301) function GeometryPointInRectangle(x,y,x_array,y_array)
302) !
303) ! Determines whether a point in xy space is within
304) ! a box
305) !
306) ! Author: Glenn Hammond
307) ! Date: 01/12/12
308) !
309)
310) implicit none
311)
312) PetscReal :: x, y
313) PetscReal :: x_array(:), y_array(:)
314)
315) PetscBool :: GeometryPointInRectangle
316)
317) PetscReal :: x1, y1
318) PetscReal :: x2, y2
319)
320) GeometryPointInRectangle = PETSC_FALSE
321)
322) ! only using first 2 values in each array
323)
324) x1 = minval(x_array(1:2))
325) y1 = minval(y_array(1:2))
326) x2 = maxval(x_array(1:2))
327) y2 = maxval(y_array(1:2))
328)
329)
330) if (x >= x1 .and. x <= x2 .and. y >= y1 .and. y <= y2) &
331) GeometryPointInRectangle = PETSC_TRUE
332)
333) end function GeometryPointInRectangle
334)
335) ! ************************************************************************** !
336)
337) subroutine GeometryCopyCoordinates(coordinates_in,coordinates_out)
338) !
339) ! Deallocates a polygonal volume object
340) !
341) ! Author: Glenn Hammond
342) ! Date: 01/12/12
343) !
344)
345) implicit none
346)
347) type(point3d_type) :: coordinates_in(:)
348) type(point3d_type), pointer :: coordinates_out(:)
349)
350) PetscInt :: num_coordinates
351) PetscInt :: i
352)
353) num_coordinates = size(coordinates_in)
354) allocate(coordinates_out(num_coordinates))
355) do i = 1, num_coordinates
356) coordinates_out(i)%x = coordinates_in(i)%x
357) coordinates_out(i)%y = coordinates_in(i)%y
358) coordinates_out(i)%z = coordinates_in(i)%z
359) enddo
360)
361) end subroutine GeometryCopyCoordinates
362)
363) ! ************************************************************************** !
364)
365) subroutine GeometryDestroyPolygonalVolume(polygonal_volume)
366) !
367) ! Deallocates a polygonal volume object
368) !
369) ! Author: Glenn Hammond
370) ! Date: 11/01/09
371) !
372)
373) implicit none
374)
375) type(polygonal_volume_type), pointer :: polygonal_volume
376)
377) if (.not.associated(polygonal_volume)) return
378)
379) if (associated(polygonal_volume%xy_coordinates)) &
380) deallocate(polygonal_volume%xy_coordinates)
381) nullify(polygonal_volume%xy_coordinates)
382) if (associated(polygonal_volume%xz_coordinates)) &
383) deallocate(polygonal_volume%xz_coordinates)
384) nullify(polygonal_volume%xz_coordinates)
385) if (associated(polygonal_volume%yz_coordinates)) &
386) deallocate(polygonal_volume%yz_coordinates)
387) nullify(polygonal_volume%yz_coordinates)
388)
389) deallocate(polygonal_volume)
390) nullify(polygonal_volume)
391)
392) end subroutine GeometryDestroyPolygonalVolume
393)
394) end module Geometry_module