0.6.2-dev0
FORCES
FORtran lib for Comp. Env. Sys.
Loading...
Searching...
No Matches
mo_poly Module Reference

Polygon calculations. More...

Data Types

interface  areapoly
 Area of polygon. More...
 
interface  center_of_mass
 Center of mass of polygon. More...
 
interface  inpoly
 Determination point of polygon. More...
 
interface  mod_pole
 Modify polygon so it covers pole correctly. More...
 
interface  mod_shift
 Shifts the (longitude) value 180 degrees. More...
 
interface  orientpoly
 Check orientation of polygon. More...
 

Functions/Subroutines

real(sp) function areapoly_sp (coord)
 Area of polygon.
 
real(sp) function, dimension(2) center_of_mass_sp (coord)
 Center of mass of polygon.
 
subroutine inpoly_sp (p, coord, erg)
 Determination point of polygon.
 
integer(i4) function orientpoly_sp (coord)
 Check orientation of polygon.
 
real(sp) function, dimension(:,:), allocatable mod_pole_sp (coord, meridian_arg)
 Modify polygon so it covers pole correctly.
 
elemental real(sp) function mod_shift_sp (x_coord, meridian_arg)
 Shifts the (longitude) value 180 degrees.
 
real(dp) function areapoly_dp (coord)
 Area of polygon.
 
real(dp) function, dimension(2) center_of_mass_dp (coord)
 Center of mass of polygon.
 
subroutine inpoly_dp (p, coord, erg)
 Determination point of polygon.
 
integer(i4) function orientpoly_dp (coord)
 Check orientation of polygon.
 
real(dp) function, dimension(:,:), allocatable mod_pole_dp (coord, meridian_arg)
 Modify polygon so it covers pole correctly.
 
elemental real(dp) function mod_shift_dp (x_coord, meridian_arg)
 Shifts the (longitude) value 180 degrees.
 

Detailed Description

Polygon calculations.

This module determines wether a 2D point lies inside, outside, or on the vertice/edge of a 2D polygon and is part of the UFZ CHS Fortran library.

Changelog
  • Juliane Mai, July 2012
    • written
  • Maren Goehler, July 2012
    • area & center of mass
Author
Juliane Mai
Date
Jul 2012

Function/Subroutine Documentation

◆ areapoly_dp()

real(dp) function mo_poly::areapoly_dp ( real(dp), dimension(:,:), intent(in)  coord)
private

Area of polygon.

Function for computing the area of a polygon (2D, convex or not). Function for computing the area of a polygon The polygon can be convex or not.

The method is only applicable for 2D polygons and points.

Example

polygon(:,1) = (/ 1.0_dp,2.0_dp,2.0_dp,1.0_dp /)
polygon(:,2) = (/ 1.0_dp,1.0_dp,2.0_dp,2.0_dp /)
area = areapoly( polygon )
! --> area = 1

See also example in test directory.

Literature

  1. http://de.wikipedia.org/wiki/Geometrischer_Schwerpunkt
Returns
Area of polygon
Parameters
[in]coordcoordinates of the polygon in question

Definition at line 487 of file mo_poly.f90.

Referenced by mo_poly::center_of_mass::center_of_mass_dp().

Here is the caller graph for this function:

◆ areapoly_sp()

real(sp) function mo_poly::areapoly_sp ( real(sp), dimension(:,:), intent(in)  coord)
private

Area of polygon.

Function for computing the area of a polygon (2D, convex or not). Function for computing the area of a polygon The polygon can be convex or not.

The method is only applicable for 2D polygons and points.

Example

polygon(:,1) = (/ 1.0_dp,2.0_dp,2.0_dp,1.0_dp /)
polygon(:,2) = (/ 1.0_dp,1.0_dp,2.0_dp,2.0_dp /)
area = areapoly( polygon )
! --> area = 1

See also example in test directory.

Literature

  1. http://de.wikipedia.org/wiki/Geometrischer_Schwerpunkt
Returns
Area of polygon
Parameters
[in]coordcoordinates of the polygon in question

Definition at line 112 of file mo_poly.f90.

Referenced by mo_poly::center_of_mass::center_of_mass_sp().

Here is the caller graph for this function:

◆ center_of_mass_dp()

real(dp) function, dimension(2) mo_poly::center_of_mass_dp ( real(dp), dimension(:,:), intent(in)  coord)
private

Center of mass of polygon.

Function for computing the center of mass of a polygon (2D, convex or not).

\[ A = \sum_{i}(x_{i}y_{i+1}-x_{i+1}y_{i}) \]

\[ x_s = \frac{1}{6A} \sum_i(x_i+x_{i+1})(x_iy_{i+1}-x_{i+1}y_i) \]

\[ y_s = \frac{1}{6A} \sum_i(y_i+y_{i+1})(x_iy_{i+1}-x_{i+1}y_i) \]

Example

polygon(:,1) = (/ 1.0_dp,2.0_dp,2.0_dp,1.0_dp /)
polygon(:,2) = (/ 1.0_dp,1.0_dp,2.0_dp,2.0_dp /)
com = center_of_mass( polygon )
! --> com = (/1.5_dp, 1.5_dp/)

See also example in test directory

Literature

  1. http://de.wikipedia.org/wiki/Geometrischer_Schwerpunkt
Returns
Center of mass of polygon.
Parameters
[in]coordcoordinates of polygon in question

Definition at line 543 of file mo_poly.f90.

◆ center_of_mass_sp()

real(sp) function, dimension(2) mo_poly::center_of_mass_sp ( real(sp), dimension(:,:), intent(in)  coord)
private

Center of mass of polygon.

Function for computing the center of mass of a polygon (2D, convex or not).

\[ A = \sum_{i}(x_{i}y_{i+1}-x_{i+1}y_{i}) \]

\[ x_s = \frac{1}{6A} \sum_i(x_i+x_{i+1})(x_iy_{i+1}-x_{i+1}y_i) \]

\[ y_s = \frac{1}{6A} \sum_i(y_i+y_{i+1})(x_iy_{i+1}-x_{i+1}y_i) \]

Example

polygon(:,1) = (/ 1.0_dp,2.0_dp,2.0_dp,1.0_dp /)
polygon(:,2) = (/ 1.0_dp,1.0_dp,2.0_dp,2.0_dp /)
com = center_of_mass( polygon )
! --> com = (/1.5_dp, 1.5_dp/)

See also example in test directory

Literature

  1. http://de.wikipedia.org/wiki/Geometrischer_Schwerpunkt
Returns
Center of mass of polygon.
Parameters
[in]coordcoordinates of polygon in question

Definition at line 168 of file mo_poly.f90.

◆ inpoly_dp()

subroutine mo_poly::inpoly_dp ( real(dp), dimension(2), intent(in)  p,
real(dp), dimension(:, :), intent(in)  coord,
integer(i4), intent(out)  erg 
)
private

Determination point of polygon.

Determines whether a 2D point is inside, outside or on vertex of a polygon (2D, convex or not).

The original version of the source code (pnpoly) was implemented by W. Randolph Franklin. It was insufficiently assigning vertex/edge points.

Example

polygon(:,1) = (/ 1.0_dp,2.0_dp,2.0_dp,1.0_dp /)
polygon(:,2) = (/ 1.0_dp,1.0_dp,2.0_dp,2.0_dp /)
point = (/ 1.5, 1.5 /)
call inpoly( point, polygon, inside )
! --> inside = 1 ... point is inside the polygon

See also example in test directory

Literature

  1. https://wrf.ecse.rpi.edu/Research/Short_Notes/pnpoly.html
Returns
Whether point is inside (=1), outside (=-1) or on a vertex/edge of the polygon (=0)
Parameters
[in]ppoint in question
[in]coordcoordinates of the polygon
[out]ergresult: inside: erg = 1 outside: erg = -1 on vertex/edge: erg = 0

Definition at line 607 of file mo_poly.f90.

◆ inpoly_sp()

subroutine mo_poly::inpoly_sp ( real(sp), dimension(2), intent(in)  p,
real(sp), dimension(:, :), intent(in)  coord,
integer(i4), intent(out)  erg 
)
private

Determination point of polygon.

Determines whether a 2D point is inside, outside or on vertex of a polygon (2D, convex or not).

The original version of the source code (pnpoly) was implemented by W. Randolph Franklin. It was insufficiently assigning vertex/edge points.

Example

polygon(:,1) = (/ 1.0_dp,2.0_dp,2.0_dp,1.0_dp /)
polygon(:,2) = (/ 1.0_dp,1.0_dp,2.0_dp,2.0_dp /)
point = (/ 1.5, 1.5 /)
call inpoly( point, polygon, inside )
! --> inside = 1 ... point is inside the polygon

See also example in test directory

Literature

  1. https://wrf.ecse.rpi.edu/Research/Short_Notes/pnpoly.html
Returns
Whether point is inside (=1), outside (=-1) or on a vertex/edge of the polygon (=0)
Parameters
[in]ppoint in question
[in]coordcoordinates of the polygon
[out]ergresult: inside: erg = 1 outside: erg = -1 on vertex/edge: erg = 0

Definition at line 232 of file mo_poly.f90.

◆ mod_pole_dp()

real(dp) function, dimension(:,:), allocatable mo_poly::mod_pole_dp ( real(dp), dimension(:,:), intent(in)  coord,
real(dp), intent(in), optional  meridian_arg 
)
private

Modify polygon so it covers pole correctly.

Modifies a polygon (2D, convex or not) to include pole when passed to inpoly this function is intended to modify a given polygon, so it can be represented in a Cartesian grid the use case is a polygon (e.g. 120,80, -120,80, 0,80 covering the north pole) that is not represented on Cartesian lat-lon grid as polygon. The script inserts additional coordinates, so the pole is covered (e.g. 180,80, 180,90, -180,90, -180,80) See test cases for examples.

Returns
modified coordinates
Parameters
[in]coordcoordinates of the polygon in question
[in]meridian_argmeridian that represents discontinuity, defaults to 180.0

Definition at line 728 of file mo_poly.f90.

◆ mod_pole_sp()

real(sp) function, dimension(:,:), allocatable mo_poly::mod_pole_sp ( real(sp), dimension(:,:), intent(in)  coord,
real(sp), intent(in), optional  meridian_arg 
)
private

Modify polygon so it covers pole correctly.

Modifies a polygon (2D, convex or not) to include pole when passed to inpoly this function is intended to modify a given polygon, so it can be represented in a Cartesian grid the use case is a polygon (e.g. 120,80, -120,80, 0,80 covering the north pole) that is not represented on Cartesian lat-lon grid as polygon. The script inserts additional coordinates, so the pole is covered (e.g. 180,80, 180,90, -180,90, -180,80) See test cases for examples.

Returns
modified coordinates
Parameters
[in]coordcoordinates of the polygon in question
[in]meridian_argmeridian that represents discontinuity, defaults to 180.0

Definition at line 353 of file mo_poly.f90.

◆ mod_shift_dp()

elemental real(dp) function mo_poly::mod_shift_dp ( real(dp), intent(in)  x_coord,
real(dp), intent(in), optional  meridian_arg 
)
private

Shifts the (longitude) value 180 degrees.

Modify a coordinate value

Returns
Shifted value
Parameters
[in]x_coordcoordinates of the polygon in question
[in]meridian_argmeridian that represents discontinuity, defaults to 180.0

Definition at line 818 of file mo_poly.f90.

◆ mod_shift_sp()

elemental real(sp) function mo_poly::mod_shift_sp ( real(sp), intent(in)  x_coord,
real(sp), intent(in), optional  meridian_arg 
)
private

Shifts the (longitude) value 180 degrees.

Modify a coordinate value

Returns
Shifted value
Parameters
[in]x_coordcoordinates of the polygon in question
[in]meridian_argmeridian that represents discontinuity, defaults to 180.0

Definition at line 443 of file mo_poly.f90.

◆ orientpoly_dp()

integer(i4) function mo_poly::orientpoly_dp ( real(dp), dimension(:,:), intent(in)  coord)
private

Check orientation of polygon.

Function for checking the orientation of a polygon (2D, convex or not).

Returns
Integer indicating orientation (counter-clockwise: -1, clockwise: 1, none: 0)
Parameters
[in]coordcoordinates of the polygon in question

result: clockwise: orientpoly = 1 counter-clockwise: orientpoly = -1 flat line or similar irregular shape: orientpoly = 0

Definition at line 691 of file mo_poly.f90.

◆ orientpoly_sp()

integer(i4) function mo_poly::orientpoly_sp ( real(sp), dimension(:,:), intent(in)  coord)
private

Check orientation of polygon.

Function for checking the orientation of a polygon (2D, convex or not).

Returns
Integer indicating orientation (counter-clockwise: -1, clockwise: 1, none: 0)
Parameters
[in]coordcoordinates of the polygon in question

result: clockwise: orientpoly = 1 counter-clockwise: orientpoly = -1 flat line or similar irregular shape: orientpoly = 0

Definition at line 316 of file mo_poly.f90.