polyCub/ 0000755 0001762 0000144 00000000000 13427056412 011673 5 ustar ligges users polyCub/inst/ 0000755 0001762 0000144 00000000000 13427003213 012637 5 ustar ligges users polyCub/inst/CITATION 0000644 0001762 0000144 00000000642 13426603244 014007 0 ustar ligges users citHeader("To cite", sQuote("polyCub"), "in publications, please refer to:")
bibentry(
bibtype = "Article", key = "R:polyCub",
author = "Sebastian Meyer",
title = "{polyCub}: An {R} package for Integration over Polygons",
journal = "Journal of Open Source Software",
issn = "2475-9066",
year = "2019",
volume = "4",
number = "34",
pages = "1056",
doi = "10.21105/joss.01056"
)
polyCub/inst/doc/ 0000755 0001762 0000144 00000000000 13427003213 013404 5 ustar ligges users polyCub/inst/doc/polyCub.R 0000644 0001762 0000144 00000004504 13427003207 015152 0 ustar ligges users ## ------------------------------------------------------------------------
library("polyCub")
## ----example-f-----------------------------------------------------------
f <- function (s, sigma = 5)
{
exp(-rowSums(s^2)/2/sigma^2) / (2*pi*sigma^2)
}
## ----example-polygon-----------------------------------------------------
hexagon <- list(
list(x = c(7.33, 7.33, 3, -1.33, -1.33, 3),
y = c(-0.5, 4.5, 7, 4.5, -0.5, -3))
)
## ----example, fig.width = 3, fig.height = 2.5----------------------------
plotpolyf(hexagon, f, xlim = c(-8,8), ylim = c(-8,8))
## ----product-Gauss, echo = -1, fig.show = "hold"-------------------------
par(mar = c(3,3,1,2))
polyCub.SV(hexagon, f, nGQ = 3, plot = TRUE)
## ------------------------------------------------------------------------
nrow(polyCub.SV(hexagon, f = NULL, nGQ = 3)[[1]]$nodes)
## ------------------------------------------------------------------------
polyCub.SVn <- function (polyregion, f, ..., nGQ = 20) {
nw <- polyCub.SV(polyregion, f = NULL, ..., nGQ = nGQ)
## nw is a list with one element per polygon of 'polyregion'
res <- sapply(nw, function (x)
c(result = sum(x$weights * f(x$nodes, ...)), nEval = nrow(x$nodes)))
structure(sum(res["result",]), nEval = sum(res["nEval",]))
}
polyCub.SVn(hexagon, f, nGQ = 3)
## ------------------------------------------------------------------------
for (nGQ in c(1:5, 10, 20)) {
result <- polyCub.SVn(hexagon, f, nGQ = nGQ)
cat(sprintf("nGQ = %2i: %.12f (n=%i)\n", nGQ, result, attr(result, "nEval")))
}
## ---- message = FALSE----------------------------------------------------
library("spatstat")
hexagon.owin <- owin(poly = hexagon)
## ----midpoint, echo = -1, fig.show = "hold"------------------------------
par(mar = c(3,3,1,3), xaxs = "i", yaxs = "i")
polyCub.midpoint(hexagon.owin, f, eps = 0.5, plot = TRUE)
## ------------------------------------------------------------------------
intrfr <- function (R, sigma = 5)
{
(1 - exp(-R^2/2/sigma^2))/2/pi
}
## ------------------------------------------------------------------------
polyCub.iso(hexagon, intrfr = intrfr, center = c(0,0))
## ------------------------------------------------------------------------
gpclibPermit() # accept gpclib license (prohibits commercial use)
polyCub.exact.Gauss(hexagon.owin, mean = c(0,0), Sigma = 5^2*diag(2))
polyCub/inst/doc/polyCub.html 0000644 0001762 0000144 00000254610 13427003213 015717 0 ustar ligges users
Getting started with polyCub
Getting started with polyCub
Sebastian Meyer
2019-02-07

The R package polyCub implements cubature (numerical integration) over polygonal domains. It solves the problem of integrating a continuously differentiable function f(x,y) over simple closed polygons.
For the special case of a rectangular domain along the axes, the package cubature is more appropriate (cf. CRAN Task View: Numerical Mathematics
).
Polygon representations
The integration domain is described by a polygonal boundary (or multiple polygons, including holes). Various R packages for spatial data analysis provide classes for polygons. The implementations differ in vertex order (which direction represents a hole) and if the first vertex is repeated.
All of polyCub’s cubature methods understand
Internally, polyCub uses its auxiliary xylist()
function to extract a plain list of lists of vertex coordinates from these classes, such that vertices are ordered anticlockwise and the first vertex is not repeated (i.e., the "owin"
convention).
Cubature methods
The following cubature methods are implemented in polyCub:
polyCub.SV()
: Product Gauss cubature
polyCub.midpoint()
: Two-dimensional midpoint rule
polyCub.iso()
: Adaptive cubature for radially symmetric functions \(f(x,y) = f_r(\lVert(x-x_0,y-y_0)\rVert)\)
polyCub.exact.Gauss()
: Accurate (but slow) integration of the bivariate Gaussian density
The following section details and illustrates the different cubature methods.
Illustrations
We consider the integration of a function f(x,y) which all of the above cubature methods can handle: an isotropic zero-mean Gaussian density. polyCub expects the integrand f to take a two-column coordinate matrix as its first argument (as opposed to separate arguments for the x and y coordinates), so:
We use a simple hexagon as polygonal integration domain, here specified via an "xylist"
of vertex coordinates:
hexagon <- list(
list(x = c(7.33, 7.33, 3, -1.33, -1.33, 3),
y = c(-0.5, 4.5, 7, 4.5, -0.5, -3))
)
An image of the function and the integration domain can be produced using polyCub’s rudimentary (but convenient) plotting utility:

1. Product Gauss cubature: polyCub.SV()
The polyCub package provides an R-interfaced C-translation of “polygauss: Matlab code for Gauss-like cubature over polygons” (Sommariva and Vianello, 2013, http://www.math.unipd.it/~alvise/software.html), an algorithm described in Sommariva and Vianello (2007, BIT Numerical Mathematics, https://doi.org/10.1007/s10543-007-0131-2). The cubature rule is based on Green’s integration formula and incorporates appropriately transformed weights and nodes of one-dimensional Gauss-Legendre quadrature in both dimensions, thus the name “product Gauss cubature”. It is exact for all bivariate polynomials if the number of cubature nodes is sufficiently large (depending on the degree of the polynomial).
For the above example, a reasonable approximation is already obtained with degree nGQ = 3
of the one-dimensional Gauss-Legendre quadrature:

The involved nodes (displayed in the figure above) and weights can be extracted by calling polyCub.SV()
with f = NULL
, e.g., to determine the number of nodes:
For illustration, we create a variant of polyCub.SV()
, which returns the number of function evaluations as an attribute:
polyCub.SVn <- function (polyregion, f, ..., nGQ = 20) {
nw <- polyCub.SV(polyregion, f = NULL, ..., nGQ = nGQ)
## nw is a list with one element per polygon of 'polyregion'
res <- sapply(nw, function (x)
c(result = sum(x$weights * f(x$nodes, ...)), nEval = nrow(x$nodes)))
structure(sum(res["result",]), nEval = sum(res["nEval",]))
}
polyCub.SVn(hexagon, f, nGQ = 3)
#> [1] 0.2741456
#> attr(,"nEval")
#> [1] 72
We can use this function to investigate how the accuracy of the approximation depends on the degree nGQ
and the associated number of cubature nodes:
2. Two-dimensional midpoint rule: polyCub.midpoint()
The two-dimensional midpoint rule in polyCub is a simple wrapper around as.im.function()
and integral.im()
from package spatstat. In other words, the polygon is represented by a binary pixel image and the integral is approximated as the sum of (pixel area * f(pixel midpoint)) over all pixels whose midpoint is part of the polygon.
To use polyCub.midpoint()
, we need to convert our polygon to spatstat’s “owin” class:
Using a pixel size of eps = 0.5
(here yielding 270 pixels), we obtain:

3. Adaptive cubature for isotropic functions: polyCub.iso()
A radially symmetric function can be expressed in terms of the distance r from its point of symmetry: f(r). If the antiderivative of r times f(r), called intrfr()
, is analytically available, Green’s theorem leads us to a cubature rule which only needs one-dimensional numerical integration. More specifically, intrfr()
will be integrate()
d along the edges of the polygon. The mathematical details are given in Meyer and Held (2014, The Annals of Applied Statistics, https://doi.org/10.1214/14-AOAS743, Supplement B, Section 2.4).
For the bivariate Gaussian density f
defined above, the integral from 0 to R of r*f(r)
is analytically available as:
With this information, we can apply the cubature rule as follows:
Note that we do not even need the original function f
.
If intrfr()
is missing, it can be approximated numerically using integrate()
for r*f(r)
as well, but the overall integration will then be much less efficient than product Gauss cubature.
Package polyCub exposes a C-version of polyCub.iso()
for use by other R packages (notably surveillance) via LinkingTo: polyCub
and #include <polyCubAPI.h>
. This requires the intrfr()
function to be implemented in C as well. See https://github.com/bastistician/polyCub/blob/master/tests/testthat/polyiso_powerlaw.c for an example.
4. Integration of the bivariate Gaussian density: polyCub.exact.Gauss()
Abramowitz and Stegun (1972, Section 26.9, Example 9) offer a formula for the integral of the bivariate Gaussian density over a triangle with one vertex at the origin. This formula can be used after triangulation of the polygonal domain (polyCub currently uses tristrip()
from the gpclib package). The core of the formula is an integral of the bivariate Gaussian density with zero mean, unit variance and some correlation over an infinite rectangle [h, Inf] x [0, Inf], which can be computed accurately using pmvnorm()
from the mvtnorm package.
For the above example, we obtain:
The required triangulation as well as the numerous calls of pmvnorm()
make this integration algorithm quiet cumbersome. For large-scale integration tasks, it is thus advisable to resort to the general-purpose product Gauss cubature rule polyCub.SV()
.
Note: polyCub provides an auxiliary function circleCub.Gauss()
to calculate the integral of an isotropic Gaussian density over a circular domain (which requires nothing more than a single call of pchisq()
).
Benchmark
We use the last result from polyCub.exact.Gauss()
as a reference value and tune the number of cubature nodes in polyCub.SV()
and polyCub.midpoint()
until the absolute error is below 10^-8. This leads to nGQ = 4
for product Gauss cubature and a 1200 x 1200 pixel image for the midpoint rule. For polyCub.iso()
, we keep the default tolerance levels of integrate()
. For comparison, we also run polyCub.iso()
without the analytically derived intrfr
function, which leads to a double-integrate
approximation.
The median runtimes [ms] of the different cubature methods are given below.
benchmark <- microbenchmark::microbenchmark(
SV = polyCub.SV(hexagon.owin, f, nGQ = 4),
midpoint = polyCub.midpoint(hexagon.owin, f, dimyx = 1200),
iso = polyCub.iso(hexagon.owin, intrfr = intrfr, center = c(0,0)),
iso_double_approx = polyCub.iso(hexagon.owin, f, center = c(0,0)),
exact = polyCub.exact.Gauss(hexagon.owin, mean = c(0,0), Sigma = 5^2*diag(2)),
times = 6,
check = function (values) all(abs(unlist(values) - 0.274144773813434) < 1e-8))
SV |
0.13 |
midpoint |
312.72 |
iso |
0.32 |
iso_double_approx |
4.76 |
exact |
7.81 |
The general-purpose SV-method is the clear winner of this small competition. A disadvantage of that method is that the number of cubature nodes needs to be tuned manually. This also holds for the midpoint rule, which is by far the slowest option. In contrast, the “iso”-method for radially symmetric functions is based on R’s integrate()
function, which implements automatic tolerance levels. Furthermore, the “iso”-method can also be used with “spiky” integrands, such as a heavy-tailed power-law kernel \(f(r) = (r+1)^{-2}\).
polyCub/inst/doc/polyCub.Rmd 0000644 0001762 0000144 00000025573 13426777177 015533 0 ustar ligges users ---
title: "Getting started with polyCub"
author: "Sebastian Meyer"
date: "`r Sys.Date()`"
output: rmarkdown::html_vignette
vignette: >
%\VignetteIndexEntry{Getting started with polyCub}
%\VignetteEngine{knitr::rmarkdown}
---
```{R setup, include = FALSE, purl = FALSE}
knitr::opts_chunk$set(collapse = TRUE, comment = "#>")
```
The R package **polyCub** implements *cubature* (numerical integration)
over *polygonal* domains.
It solves the problem of integrating a continuously differentiable
function f(x,y) over simple closed polygons.
For the special case of a rectangular domain along the axes, the package
[**cubature**](https://CRAN.R-project.org/package=cubature)
is more appropriate (cf.
[`CRAN Task View: Numerical Mathematics`](https://CRAN.R-project.org/view=NumericalMathematics)).
## Polygon representations
The integration domain is described by a polygonal boundary
(or multiple polygons, including holes).
Various R packages for spatial data analysis provide classes for polygons.
The implementations differ in vertex order (which direction represents a hole)
and if the first vertex is repeated.
All of **polyCub**'s cubature methods understand
* `"owin"` from package [**spatstat**](https://CRAN.R-project.org/package=spatstat),
* `"gpc.poly"` from [**rgeos**](https://CRAN.R-project.org/package=rgeos)
(or [**gpclib**](https://CRAN.R-project.org/package=gpclib)), and
* `"SpatialPolygons"` from package [**sp**](https://CRAN.R-project.org/package=sp).
Internally, **polyCub** uses its auxiliary `xylist()` function to extract
a plain list of lists of vertex coordinates from these classes,
such that vertices are ordered anticlockwise and the first vertex is not
repeated (i.e., the `"owin"` convention).
## Cubature methods
The following cubature methods are implemented in **polyCub**:
1. `polyCub.SV()`: Product Gauss cubature
2. `polyCub.midpoint()`: Two-dimensional midpoint rule
3. `polyCub.iso()`: Adaptive cubature for radially symmetric functions
$f(x,y) = f_r(\lVert(x-x_0,y-y_0)\rVert)$
4. `polyCub.exact.Gauss()`: Accurate (but slow) integration of the
bivariate Gaussian density
The following section details and illustrates the different cubature methods.
## Illustrations
```{R}
library("polyCub")
```
We consider the integration of a function f(x,y) which all of the above
cubature methods can handle: an isotropic zero-mean Gaussian density.
**polyCub** expects the integrand f to take a two-column
coordinate matrix as its first argument (as opposed to separate arguments
for the x and y coordinates), so:
```{R example-f}
f <- function (s, sigma = 5)
{
exp(-rowSums(s^2)/2/sigma^2) / (2*pi*sigma^2)
}
```
We use a simple hexagon as polygonal integration domain,
here specified via an `"xylist"` of vertex coordinates:
```{R example-polygon}
hexagon <- list(
list(x = c(7.33, 7.33, 3, -1.33, -1.33, 3),
y = c(-0.5, 4.5, 7, 4.5, -0.5, -3))
)
```
An image of the function and the integration domain can be produced using
**polyCub**'s rudimentary (but convenient) plotting utility:
```{R example, fig.width = 3, fig.height = 2.5}
plotpolyf(hexagon, f, xlim = c(-8,8), ylim = c(-8,8))
```
### 1. Product Gauss cubature: `polyCub.SV()`
The **polyCub** package provides an R-interfaced C-translation of
"polygauss: Matlab code for Gauss-like cubature over polygons"
(Sommariva and Vianello, 2013, ),
an algorithm described in Sommariva and Vianello (2007,
*BIT Numerical Mathematics*, ).
The cubature rule is based on Green's integration formula and incorporates
appropriately transformed weights and nodes of one-dimensional
Gauss-Legendre quadrature in both dimensions,
thus the name "product Gauss cubature".
It is exact for all bivariate polynomials if the number of cubature nodes
is sufficiently large (depending on the degree of the polynomial).
For the above example, a reasonable approximation is already obtained
with degree `nGQ = 3` of the one-dimensional Gauss-Legendre quadrature:
```{R product-Gauss, echo = -1, fig.show = "hold"}
par(mar = c(3,3,1,2))
polyCub.SV(hexagon, f, nGQ = 3, plot = TRUE)
```
The involved nodes (displayed in the figure above) and weights can be
extracted by calling `polyCub.SV()` with `f = NULL`, e.g., to determine
the number of nodes:
```{R}
nrow(polyCub.SV(hexagon, f = NULL, nGQ = 3)[[1]]$nodes)
```
For illustration, we create a variant of `polyCub.SV()`,
which returns the number of function evaluations as an attribute:
```{R}
polyCub.SVn <- function (polyregion, f, ..., nGQ = 20) {
nw <- polyCub.SV(polyregion, f = NULL, ..., nGQ = nGQ)
## nw is a list with one element per polygon of 'polyregion'
res <- sapply(nw, function (x)
c(result = sum(x$weights * f(x$nodes, ...)), nEval = nrow(x$nodes)))
structure(sum(res["result",]), nEval = sum(res["nEval",]))
}
polyCub.SVn(hexagon, f, nGQ = 3)
```
We can use this function to investigate how the accuracy of the
approximation depends on the degree `nGQ` and the associated number of
cubature nodes:
```{R}
for (nGQ in c(1:5, 10, 20)) {
result <- polyCub.SVn(hexagon, f, nGQ = nGQ)
cat(sprintf("nGQ = %2i: %.12f (n=%i)\n", nGQ, result, attr(result, "nEval")))
}
```
### 2. Two-dimensional midpoint rule: `polyCub.midpoint()`
The two-dimensional midpoint rule in **polyCub** is a simple wrapper
around `as.im.function()` and `integral.im()` from package **spatstat**.
In other words, the polygon is represented by a binary pixel image and
the integral is approximated as the sum of (pixel area * f(pixel midpoint))
over all pixels whose midpoint is part of the polygon.
To use `polyCub.midpoint()`, we need to convert our polygon to
**spatstat**'s "owin" class:
```{R, message = FALSE}
library("spatstat")
hexagon.owin <- owin(poly = hexagon)
```
Using a pixel size of `eps = 0.5` (here yielding 270 pixels), we obtain:
```{R midpoint, echo = -1, fig.show = "hold"}
par(mar = c(3,3,1,3), xaxs = "i", yaxs = "i")
polyCub.midpoint(hexagon.owin, f, eps = 0.5, plot = TRUE)
```
### 3. Adaptive cubature for *isotropic* functions: `polyCub.iso()`
A radially symmetric function can be expressed in terms of
the distance r from its point of symmetry: f(r).
If the antiderivative of r times f(r), called `intrfr()`, is
analytically available, Green's theorem leads us to a cubature rule
which only needs *one-dimensional* numerical integration.
More specifically, `intrfr()` will be `integrate()`d along the edges of
the polygon. The mathematical details are given in
Meyer and Held (2014, *The Annals of Applied Statistics*,
, Supplement B, Section 2.4).
For the bivariate Gaussian density `f` defined above,
the integral from 0 to R of `r*f(r)` is analytically available as:
```{R}
intrfr <- function (R, sigma = 5)
{
(1 - exp(-R^2/2/sigma^2))/2/pi
}
```
With this information, we can apply the cubature rule as follows:
```{R}
polyCub.iso(hexagon, intrfr = intrfr, center = c(0,0))
```
Note that we do not even need the original function `f`.
If `intrfr()` is missing, it can be approximated numerically using
`integrate()` for `r*f(r)` as well, but the overall integration will then
be much less efficient than product Gauss cubature.
Package **polyCub** exposes a C-version of `polyCub.iso()`
for use by other R packages (notably
[**surveillance**](https://CRAN.R-project.org/package=surveillance)) via
`LinkingTo: polyCub` and `#include `.
This requires the `intrfr()` function to be implemented in C as well. See
for an example.
### 4. Integration of the *bivariate Gaussian density*: `polyCub.exact.Gauss()`
Abramowitz and Stegun (1972, Section 26.9, Example 9) offer a formula for
the integral of the bivariate Gaussian density over a triangle with one
vertex at the origin. This formula can be used after triangulation of
the polygonal domain (**polyCub** currently uses `tristrip()`
from the [**gpclib**](https://CRAN.R-project.org/package=gpclib) package).
The core of the formula is an integral of the bivariate Gaussian density
with zero mean, unit variance and some correlation over an infinite rectangle
[h, Inf] x [0, Inf], which can be computed accurately using `pmvnorm()`
from the [**mvtnorm**](https://CRAN.R-project.org/package=mvtnorm) package.
For the above example, we obtain:
```{R}
gpclibPermit() # accept gpclib license (prohibits commercial use)
polyCub.exact.Gauss(hexagon.owin, mean = c(0,0), Sigma = 5^2*diag(2))
```
The required triangulation as well as the numerous calls of `pmvnorm()`
make this integration algorithm quiet cumbersome. For large-scale
integration tasks, it is thus advisable to resort to the general-purpose
product Gauss cubature rule `polyCub.SV()`.
Note: **polyCub** provides an auxiliary function `circleCub.Gauss()` to
calculate the integral of an *isotropic* Gaussian density over a *circular*
domain (which requires nothing more than a single call of `pchisq()`).
## Benchmark
We use the last result from `polyCub.exact.Gauss()` as a reference value and
tune the number of cubature nodes in `polyCub.SV()` and `polyCub.midpoint()`
until the absolute error is below 10^-8.
This leads to `nGQ = 4` for product Gauss cubature
and a 1200 x 1200 pixel image for the midpoint rule.
For `polyCub.iso()`, we keep the default tolerance levels of `integrate()`.
For comparison, we also run `polyCub.iso()` without the analytically derived
`intrfr` function, which leads to a double-`integrate` approximation.
The median runtimes [ms] of the different cubature methods are given below.
```{r benchmark, purl = FALSE, eval = identical(Sys.getenv("NOT_CRAN"), "true")}
benchmark <- microbenchmark::microbenchmark(
SV = polyCub.SV(hexagon.owin, f, nGQ = 4),
midpoint = polyCub.midpoint(hexagon.owin, f, dimyx = 1200),
iso = polyCub.iso(hexagon.owin, intrfr = intrfr, center = c(0,0)),
iso_double_approx = polyCub.iso(hexagon.owin, f, center = c(0,0)),
exact = polyCub.exact.Gauss(hexagon.owin, mean = c(0,0), Sigma = 5^2*diag(2)),
times = 6,
check = function (values) all(abs(unlist(values) - 0.274144773813434) < 1e-8))
```
```{r, purl = FALSE, eval = FALSE}
summary(benchmark, unit = "ms")[c("expr", "median")]
```
```{r, purl = FALSE, echo = FALSE, eval = identical(Sys.getenv("NOT_CRAN"), "true")}
knitr::kable(summary(benchmark, unit = "ms")[c("expr", "median")], digits = 2)
```
The general-purpose SV-method is the clear winner of this small competition.
A disadvantage of that method is that the number of cubature nodes needs to be
tuned manually. This also holds for the midpoint rule, which is by far the
slowest option. In contrast, the "iso"-method for radially symmetric functions
is based on R's `integrate()` function, which implements automatic tolerance
levels. Furthermore, the "iso"-method can also be used with "spiky" integrands,
such as a heavy-tailed power-law kernel $f(r) = (r+1)^{-2}$.
polyCub/inst/include/ 0000755 0001762 0000144 00000000000 13163463332 014273 5 ustar ligges users polyCub/inst/include/polyCubAPI.h 0000644 0001762 0000144 00000003257 13163463332 016422 0 ustar ligges users /*******************************************************************************
* Header file with wrapper functions for the C-routines provided by polyCub
*
* Copyright (C) 2017 Sebastian Meyer
*
* This file is part of the R package "polyCub",
* free software under the terms of the GNU General Public License, version 2,
* a copy of which is available at https://www.R-project.org/Licenses/.
******************************************************************************/
#include // NULL
#include // SEXP
#include // R_GetCCallable
typedef double (*intrfr_fn) (double, double*);
void polyCub_iso(
double *x, double *y, // vertex coordinates (open)
int *L, // number of vertices
intrfr_fn intrfr, // F(R)
double *pars, // parameters for F(R)
double *center_x, double *center_y, // center of isotropy
int *subdivisions, double *epsabs, double *epsrel, // Rdqags options
int *stop_on_error, // !=0 means to stop at first ier > 0
double *value, double *abserr, int *neval) // results
{
static void(*fun)(double*,double*,int*,intrfr_fn,double*,double*,double*,
int*,double*,double*,int*,double*,double*,int*) = NULL;
if (fun == NULL)
fun = (void(*)(double*,double*,int*,intrfr_fn,double*,double*,double*,
int*,double*,double*,int*,double*,double*,int*))
R_GetCCallable("polyCub", "polyiso");
fun(x, y, L, intrfr, pars, center_x, center_y,
subdivisions, epsabs, epsrel, stop_on_error,
value, abserr, neval);
return;
}
polyCub/tests/ 0000755 0001762 0000144 00000000000 13357377573 013055 5 ustar ligges users polyCub/tests/testthat/ 0000755 0001762 0000144 00000000000 13357377573 014715 5 ustar ligges users polyCub/tests/testthat/polyiso_powerlaw.c 0000644 0001762 0000144 00000003447 13163463332 020467 0 ustar ligges users /*******************************************************************************
* Example of using the C-routine "polyCub_iso", see also test-polyiso.R
*
* Copyright (C) 2015,2017 Sebastian Meyer
*
* This file is part of the R package "polyCub",
* free software under the terms of the GNU General Public License, version 2,
* a copy of which is available at https://www.R-project.org/Licenses/.
******************************************************************************/
#include
#include
// F(R) example
static double intrfr_powerlaw(double R, double *logpars)
{
double sigma = exp(logpars[0]);
double d = exp(logpars[1]);
if (d == 1.0) {
return R - sigma * log(R/sigma + 1);
} else if (d == 2.0) {
return log(R/sigma + 1) - R/(R+sigma);
} else {
return (R*pow(R+sigma,1-d) - (pow(R+sigma,2-d) - pow(sigma,2-d))/(2-d)) / (1-d);
}
}
// function to be called from R
void C_polyiso_powerlaw(
double *x, double *y, // vertex coordinates (open)
int *L, // number of vertices
//intrfr_fn intrfr, // F(R)
double *pars, // parameters for F(R)
double *center_x, double *center_y, // center of isotropy
int *subdivisions, double *epsabs, double *epsrel, // Rdqags options
int *stop_on_error, // !=0 means to stop at first ier > 0
double *value, double *abserr, int *neval) // results
{
polyCub_iso(x, y, L,
intrfr_powerlaw,
pars, center_x, center_y,
subdivisions, epsabs, epsrel,
stop_on_error,
value, abserr, neval);
return;
}
polyCub/tests/testthat/test-NWGL.R 0000644 0001762 0000144 00000001674 13357377573 016574 0 ustar ligges users context("Validation of cached Gauss-Legendre nodes/weights")
test_that("statmod::gauss.quad() still gives the same result", {
new.NWGL <- lapply(
X = seq_len(61L),
FUN = function (n)
unname(statmod::gauss.quad(n = n, kind = "legendre"))
)
expect_equal(new.NWGL, .NWGL, check.attributes = FALSE)
})
test_that("polyCub.SV() can fetch nodes and weights from 'statmod'", {
diamond <- list(list(x = c(1,2,1,0), y = c(1,2,3,2)))
nw <- polyCub.SV(diamond, f = NULL, nGQ = 83)
expect_is(nw, "list")
})
test_that("polyCub.SV() can reduce nodes with zero weight", {
rectangle <- list(list(x = c(-1,1,1,-1), y = c(1,1,2,2)))
##nw0 <- polyCub.SV(rectangle, f = NULL, nGQ = 3, engine = "C")[[1]] # 0s
nw <- polyCub.SV(rectangle, f = NULL, nGQ = 3, engine = "C+reduce")[[1]]
expect_true(all(nw$weights != 0))
##f <- function (s) 1 # => calculate area (= 2)
expect_equal(sum(nw$weights), 2)
})
polyCub/tests/testthat/test-regression.R 0000644 0001762 0000144 00000001163 13165175161 020160 0 ustar ligges users context("Regression tests")
hexagon <- list(
list(x = c(7.33, 7.33, 3, -1.33, -1.33, 3),
y = c(-0.5, 4.5, 7, 4.5, -0.5, -3))
)
f <- function (s) (rowSums(s^2)+1)^-2
##plotpolyf(hexagon, f)
test_that("isotropic cubature can handle control list for integrate()", {
## previosly, passing control arguments did not work
int1 <- polyCub.iso(hexagon, f, center=c(0,0), control=list(rel.tol=1e-3))
int2 <- polyCub.iso(hexagon, f, center=c(0,0), control=list(rel.tol=1e-8))
## results are almost but not identical
expect_equal(int1, int2, tolerance = 1e-3)
expect_false(identical(int1, int2))
})
polyCub/tests/testthat/test-polyCub.R 0000644 0001762 0000144 00000004415 13357136712 017422 0 ustar ligges users context("Correctness of cubature methods")
### set up test case
## bivariate, isotropic Gaussian density
f <- function (s, mean, sd)
dnorm(s[,1], mean=mean[1], sd=sd) * dnorm(s[,2], mean=mean[2], sd=sd)
## circular domain represented by a polygon
r <- 5
center <- c(3,2)
npoly <- 128
disc.owin <- spatstat::disc(radius=r, centre=center, npoly=npoly)
## parameters for f
m <- c(1,1)
sd <- 3
## target value of the integral over the _polygonal_ circle
intExact <- 0.65844436
## taken from exact.Gauss cubature
test_that("gpclibCheck() fails without prior license agreement", {
if (gpclibPermitStatus())
skip("gpclib license has already been accepted")
expect_error(polyCub:::gpclibCheck())
})
if (requireNamespace("mvtnorm") && gpclibPermit()) {
## run this conditionally since gpclib might not be available on all
## platforms (as pointed out by Uwe Ligges, 2014-04-20)
test_that("polyCub.exact.Gauss returns validated result", {
int <- polyCub.exact.Gauss(disc.owin, mean=m, Sigma=sd^2*diag(2))
expect_equal(int, intExact, tolerance=1e-8, check.attributes=FALSE)
})
}
### perform the tests (check against each other)
test_that("polyCub.exact.Gauss and circleCub.Gauss give similar results", {
## exact value of the integral over the _real_ circle
intExact_circle <- circleCub.Gauss(center=center, r=r, mean=m, sd=sd)
## how well this fits with the exact integral over a polyonal approximation
## of the circle depends of course on 'npoly'
expect_equal(intExact, intExact_circle,
tolerance=0.001, check.attributes=FALSE)
})
test_that("midpoint-cubature is correct", {
int <- polyCub.midpoint(disc.owin, f, mean=m, sd=sd, dimyx=500)
expect_equal(int, intExact, tolerance=0.001, check.attributes=FALSE)
})
test_that("SV-cubature is correct", {
intC <- polyCub.SV(disc.owin, f, mean=m, sd=sd, nGQ=3, engine="C")
intR <- polyCub.SV(disc.owin, f, mean=m, sd=sd, nGQ=3, engine="R")
expect_equal(intC, intR)
expect_equal(intC, intExact, tolerance=0.0001, check.attributes=FALSE)
})
test_that("isotropic cubature is correct", {
## using a numerical approximation of intrfr
int0 <- polyCub.iso(disc.owin, f, mean=m, sd=sd, center=m)
expect_equal(int0, intExact, check.attributes=FALSE)
})
polyCub/tests/testthat/test-polyiso.R 0000644 0001762 0000144 00000010135 13357131525 017474 0 ustar ligges users context("polyCub_iso C-routine (API)")
## CAVE (as of R-3.4.0 with testthat 1.0.2):
## During R CMD check, tools:::.runPackageTests() sets R_TESTS=startup.Rs,
## a file which is created in the parent directory "tests", see
## file.path(R.home("share"), "R", "tests-startup.R")
## for its contents. However, testthat tests are run with the working directory
## set to here, so auxiliary R sessions initiated here would fail when trying
## to source() the R_TESTS file on startup, see the system Rprofile file
## file.path(R.home("library"), "base", "R", "Rprofile")
## for what happens. Solution: unset R_TESTS (or set to "") for sub-R processes.
## function to call an R CMD with environment variables
## 'env' specified as a named character vector
Rcmd <- function (args, env = character(), ...) {
stopifnot(is.vector(env, mode = "character"),
!is.null(names(env)))
if (.Platform$OS.type == "windows") {
if (length(env)) {
## the 'env' argument of system2() is not supported on Windows
setenv <- function (envs) {
old <- Sys.getenv(names(envs), unset = NA, names = TRUE)
set <- !is.na(envs)
if (any(set)) do.call(Sys.setenv, as.list(envs[set]))
if (any(!set)) Sys.unsetenv(names(envs)[!set])
invisible(old)
}
oldenv <- setenv(env)
on.exit(setenv(oldenv))
}
system2(command = file.path(R.home("bin"), "Rcmd.exe"),
args = args, ...)
} else {
system2(command = file.path(R.home("bin"), "R"),
args = c("CMD", args),
env = paste(names(env), env, sep = "="), ...)
}
}
message("compiling polyiso_powerlaw.c using R CMD SHLIB")
shlib_error <- Rcmd(
args = c("SHLIB", "--clean", "polyiso_powerlaw.c"),
env = c("PKG_CPPFLAGS" = paste0(
"-I", system.file("include", package="polyCub")
),
"R_TESTS" = "")
)
if (shlib_error)
skip("failed to build the shared object/DLL for the polyCub_iso example")
## load shared object/DLL
myDLL <- paste0("polyiso_powerlaw", .Platform$dynlib.ext)
loadNamespace("polyCub")
dyn.load(myDLL)
## R function calling C_polyiso_powerlaw
polyiso_powerlaw <- function (xypoly, logpars, center,
subdivisions = 100L,
rel.tol = .Machine$double.eps^0.25,
abs.tol = rel.tol, stop.on.error = TRUE)
{
.C("C_polyiso_powerlaw",
as.double(xypoly$x), as.double(xypoly$y), as.integer(length(xypoly$x)),
as.double(logpars),
as.double(center[1L]), as.double(center[2L]),
as.integer(subdivisions), as.double(abs.tol), as.double(rel.tol),
as.integer(stop.on.error),
value = double(1L), abserr = double(1L), neval = integer(1L)
)[c("value", "abserr", "neval")]
}
## example polygon and function parameters
diamond <- list(x = c(1,2,1,0), y = c(1,2,3,2))
logpars <- log(c(0.5, 1))
center <- c(0.5,2.5) # lies on an edge (to cover that case as well)
(res <- polyiso_powerlaw(xypoly = diamond,
logpars = logpars,
center = center))
## compare with R implementation
intrfr.powerlaw <- function (R, logpars)
{
sigma <- exp(logpars[[1L]])
d <- exp(logpars[[2L]])
if (d == 1) {
R - sigma * log(R/sigma + 1)
} else if (d == 2) {
log(R/sigma + 1) - R/(R+sigma)
} else {
(R*(R+sigma)^(1-d) - ((R+sigma)^(2-d) - sigma^(2-d))/(2-d)) / (1-d)
}
}
(orig <- polyCub:::polyCub1.iso(poly = diamond,
intrfr = intrfr.powerlaw,
logpars = logpars,
center = center))
test_that("C and R implementations give equal results", {
expect_equal(res$value, orig[1L])
expect_equal(res$abserr, orig[2L])
})
## microbenchmark::microbenchmark(
## polyCub:::polyCub1.iso(diamond, intrfr.powerlaw, logpars, center=center),
## polyiso_powerlaw(diamond, logpars, center=center),
## times = 1000)
## ## 140 mus vs. 35 mus
dyn.unload(myDLL)
file.remove(myDLL)
polyCub/tests/test-all.R 0000644 0001762 0000144 00000000136 13111116012 014664 0 ustar ligges users if (require("testthat") && packageVersion("testthat") >= "0.9") {
test_check("polyCub")
}
polyCub/src/ 0000755 0001762 0000144 00000000000 13427003213 012451 5 ustar ligges users polyCub/src/polyCub.SV.h 0000644 0001762 0000144 00000001601 13427003213 014564 0 ustar ligges users /*******************************************************************************
* Header file of polyCub.SV.c
*
* Copyright (C) 2017 Sebastian Meyer
*
* This file is part of the R package "polyCub",
* free software under the terms of the GNU General Public License, version 2,
* a copy of which is available at https://www.R-project.org/Licenses/.
******************************************************************************/
void C_polygauss(
double *x, double *y, // vertex coordinates (open) of a polygon
double *s_M, double *w_M, // nodes & weights of Gauss-Legendre quadrature
double *s_N, double *w_N, // of degree M=N+1 and N, respectively
double *alpha, // base-line
int *L, int *M, int *N, // L: number of edges/vertices
// result: nodes and weights of length (<=) M*N per edge
double *nodes_x, double *nodes_y, double *weights);
polyCub/src/polyCub.SV.c 0000644 0001762 0000144 00000005716 13427003213 014572 0 ustar ligges users /*******************************************************************************
* C-version of .polygauss.side()
*
* Copyright (C) 2014,2017 Sebastian Meyer
*
* This file is part of the R package "polyCub",
* free software under the terms of the GNU General Public License, version 2,
* a copy of which is available at https://www.R-project.org/Licenses/.
******************************************************************************/
#include "polyCub.SV.h"
static void C_polygauss_side(
double *x1, double *y1, double *x2, double *y2,
double *s_loc, double *w_loc, double *s_N, double *w_N,
double *alpha,
int *loc, int *N, // lengths (loc is M=N+1 or N)
// *loc * *N nodes and weights will be computed
double *nodes_x, double *nodes_y, double *weights)
{
double half_pt_x = (*x1 + *x2) / 2.0;
double half_length_x = (*x2 - *x1) / 2.0;
double half_pt_y = (*y1 + *y2) / 2.0;
double half_length_y = (*y2 - *y1) / 2.0;
double x_gauss_side, y_gauss_side, scaling_fact_minus;
int idx;
for (int i = 0; i < *loc; i++) {
// GAUSSIAN POINTS ON THE SIDE
x_gauss_side = half_pt_x + half_length_x * s_loc[i];
y_gauss_side = half_pt_y + half_length_y * s_loc[i];
scaling_fact_minus = (x_gauss_side - *alpha) / 2.0;
// COMPUTE NODES AND WEIGHTS
for (int j = 0; j < *N; j++) {
idx = j * *loc + i; // use same order as in R implementation
nodes_x[idx] = *alpha + scaling_fact_minus * (s_N[j] + 1.0);
nodes_y[idx] = y_gauss_side;
weights[idx] = half_length_y*scaling_fact_minus * w_loc[i] * w_N[j];
}
}
}
/***
* Function to be called from R to loop over all polygon edges,
* calling the above C_polygauss_side() for each
***/
void C_polygauss(
double *x, double *y, // vertex coordinates (open) of a polygon
double *s_M, double *w_M, // nodes & weights of Gauss-Legendre quadrature
double *s_N, double *w_N, // of degree M=N+1 and N, respectively
double *alpha, // base-line
int *L, int *M, int *N, // L: number of edges/vertices
// result: nodes and weights of length (<=) M*N per edge
double *nodes_x, double *nodes_y, double *weights)
{
int idxTo, idxBlock;
double x1, y1, x2, y2;
for (int i = 0; i < *L; i++) {
x1 = x[i]; y1 = y[i];
if (i == *L-1) idxTo = 0; else idxTo = i+1;
x2 = x[idxTo]; y2 = y[idxTo];
// if edge is on base-line or is orthogonal to it -> skip
if ((x1 == *alpha && x2 == *alpha) || (y2 == y1))
continue;
idxBlock = i * *M * *N; // start index of nodes of edge i
if (x2 == x1)
// side is parallel to base-line -> use degree N in both dimensions
C_polygauss_side(&x1, &y1, &x2, &y2,
s_N, w_N, s_N, w_N, alpha,
N, N,
nodes_x + idxBlock, nodes_y + idxBlock, weights + idxBlock);
else
// use degrees M and N, respectively
C_polygauss_side(&x1, &y1, &x2, &y2,
s_M, w_M, s_N, w_N, alpha,
M, N,
nodes_x + idxBlock, nodes_y + idxBlock, weights + idxBlock);
}
}
polyCub/src/polyCub.iso.h 0000644 0001762 0000144 00000001767 13427003213 015043 0 ustar ligges users /*******************************************************************************
* Header file of polyCub.iso.c
*
* Copyright (C) 2017 Sebastian Meyer
*
* This file is part of the R package "polyCub",
* free software under the terms of the GNU General Public License, version 2,
* a copy of which is available at https://www.R-project.org/Licenses/.
******************************************************************************/
typedef double (*intrfr_fn) (double, double*);
void polyiso(
double *x, double *y, // vertex coordinates (open)
int *L, // number of vertices
intrfr_fn intrfr, // F(R)
double *pars, // parameters for F(R)
double *center_x, double *center_y, // center of isotropy
int *subdivisions, double *epsabs, double *epsrel, // Rdqags options
int *stop_on_error, // !=0 means to stop at first ier > 0
double *value, double *abserr, int *neval); // results
polyCub/src/polyCub.iso.c 0000644 0001762 0000144 00000011714 13427003213 015027 0 ustar ligges users /*******************************************************************************
* C-version of polyCub1.iso()
*
* Copyright (C) 2015,2017 Sebastian Meyer
*
* This file is part of the R package "polyCub",
* free software under the terms of the GNU General Public License, version 2,
* a copy of which is available at https://www.R-project.org/Licenses/.
******************************************************************************/
/* The corresponding math is derived in Supplement B (Section 2.4) of
* Meyer and Held (2014): "Power-law models for infectious disease spread."
* The Annals of Applied Statistics, 8(3), 1612-1639.
* https://doi.org/10.1214/14-AOAS743SUPPB
*/
#include // R_FINITE, otherwise math.h would suffice
#include // error
#include // R_alloc
#include // Rprintf
#include // Rdqags
// header file defines the intrfr_fn type
#include "polyCub.iso.h"
// integrand for the edge (x0,y0) -> (x1,y1), see Equation 7
static double lineIntegrand(
double t,
double x0, double y0, double x1, double y1,
intrfr_fn intrfr, double *pars)
{
double num = y1*x0 - x1*y0; // numerator term
// point on the edge corresponding to t
double px = x0 + t*(x1-x0);
double py = y0 + t*(y1-y0);
double norm2 = px*px + py*py;
// evaluate F(R) = int_0^R r*f(r) dr at R=||(px,py)||
double inti = intrfr(sqrt(norm2), pars);
if (!R_FINITE(inti))
error("non-finite intrfr value at R=%f", sqrt(norm2));
return num*inti/norm2;
}
// set of parameters for line integration (passed via the *ex argument)
typedef struct {
double x0, y0, x1, y1;
intrfr_fn intrfr;
double *pars;
} Params;
// vectorized lineIntegrand for use with Rdqags
static void myintegr_fn(double *x, int n, void *ex)
{
Params *param = (Params *) ex;
for(int i = 0; i < n; i++) {
x[i] = lineIntegrand(x[i],
param->x0, param->y0, param->x1, param->y1,
param->intrfr, param->pars);
}
return;
}
// calculate line integral for one edge (x0,y0) -> (x1,y1)
// using Gauss-Kronrod quadrature via Rdqags as declared in ,
// implemented in R/src/appl/integrate.c,
// and used in R/src/library/stats/src/integrate.c
static void polyiso_side(
double x0, double y0, double x1, double y1, // 2 vertices
intrfr_fn intrfr, double *pars, // F(R)
int subdivisions, double *epsabs, double *epsrel, // control
double *result, double *abserr, int *neval, int *ier) // results
{
double num = y1*x0 - x1*y0; // numerator in lineIntegrand
// for any point p on the edge
if (num == 0.0) { // 'center' is part of this polygon edge
*result = 0.0;
*abserr = 0.0;
//*last = 0;
*neval = 0;
*ier = 0;
return;
}
// set of parameters for lineIntegrand
Params param = {x0, y0, x1, y1, intrfr, pars};
// prepare for Rdqags
double lower = 0.0;
double upper = 1.0;
int lenw = 4 * subdivisions;
int last; // unused
int *iwork = (int *) R_alloc((size_t) subdivisions, sizeof(int));
double *work = (double *) R_alloc((size_t) lenw, sizeof(double));
Rdqags(myintegr_fn, ¶m, &lower, &upper,
epsabs, epsrel,
result, abserr, neval, ier, // results
&subdivisions, &lenw, &last, iwork, work);
return;
}
// line integration along the edges of a polygon
void polyiso(
double *x, double *y, // vertex coordinates (open)
int *L, // number of vertices
intrfr_fn intrfr, // F(R)
double *pars, // parameters for F(R)
double *center_x, double *center_y, // center of isotropy
int *subdivisions, double *epsabs, double *epsrel, // Rdqags options
int *stop_on_error, // !=0 means to stop at first ier > 0
double *value, double *abserr, int *neval) // results
{
// auxiliary variables
double resulti, abserri;
int nevali, ieri;
double x0, y0, x1, y1;
int idxTo;
// initialize result at 0 (do += for each polygon edge);
*value = 0.0;
*abserr = 0.0;
*neval = 0;
for (int i = 0; i < *L; i++) {
x0 = x[i] - *center_x;
y0 = y[i] - *center_y;
idxTo = (i == *L-1) ? 0 : i+1;
x1 = x[idxTo] - *center_x;
y1 = y[idxTo] - *center_y;
polyiso_side(x0, y0, x1, y1,
intrfr, pars,
*subdivisions, epsabs, epsrel,
&resulti, &abserri, &nevali, &ieri);
if (ieri > 0) {
if (*stop_on_error == 0) {
Rprintf("abnormal termination of integration routine (%i)\n", ieri);
} else {
error("abnormal termination of integration routine (%i)\n", ieri);
}
}
*value += resulti;
*abserr += abserri;
*neval += nevali;
}
return;
}
polyCub/src/init.c 0000644 0001762 0000144 00000002362 13427003213 013563 0 ustar ligges users /*******************************************************************************
* Registering native routines (entry points in compiled code)
*
* Copyright (C) 2017,2019 Sebastian Meyer
*
* This file is part of the R package "polyCub",
* free software under the terms of the GNU General Public License, version 2,
* a copy of which is available at https://www.R-project.org/Licenses/.
******************************************************************************/
#include // for NULL
#include // for SEXP types
#include
#include "polyCub.SV.h"
#include "polyCub.iso.h"
// types array (could be omitted)
static R_NativePrimitiveArgType C_polygauss_t[] = {
REALSXP, REALSXP, REALSXP, REALSXP, REALSXP, REALSXP, REALSXP,
/*L, M, N:*/ INTSXP, INTSXP, INTSXP,
/*results:*/ REALSXP, REALSXP, REALSXP
};
static const R_CMethodDef CEntries[] = {
{"C_polygauss", (DL_FUNC) &C_polygauss, 13, C_polygauss_t},
{NULL, NULL, 0, NULL}
};
void R_init_polyCub(DllInfo *dll)
{
R_registerRoutines(dll, CEntries, NULL, NULL, NULL);
R_useDynamicSymbols(dll, FALSE);
//R_forceSymbols(dll, TRUE); // would require R >= 3.0.0
R_RegisterCCallable("polyCub", "polyiso", (DL_FUNC) &polyiso);
}
polyCub/NAMESPACE 0000644 0001762 0000144 00000003024 13357664114 013117 0 ustar ligges users # Generated by roxygen2: do not edit by hand
S3method(xylist,Polygon)
S3method(xylist,Polygons)
S3method(xylist,SpatialPolygons)
S3method(xylist,default)
S3method(xylist,gpc.poly)
S3method(xylist,owin)
export(.polyCub.iso)
export(as.owin.Polygon)
export(as.owin.Polygons)
export(as.owin.SpatialPolygons)
export(as.owin.gpc.poly)
export(checkintrfr)
export(circleCub.Gauss)
export(gpc2owin)
export(gpclibPermit)
export(gpclibPermitStatus)
export(owin2gpc)
export(plotpolyf)
export(polyCub)
export(polyCub.SV)
export(polyCub.exact.Gauss)
export(polyCub.iso)
export(polyCub.midpoint)
export(xylist)
exportMethods(coerce)
if(getRversion() >= "3.6.0") { # delayed registration
S3method(spatstat::as.owin, SpatialPolygons)
S3method(spatstat::as.owin, Polygons)
S3method(spatstat::as.owin, Polygon)
}
if(getRversion() >= "3.6.0") { # delayed registration
S3method(spatstat::as.owin, gpc.poly)
}
import(methods)
importClassesFrom(sp,Polygon)
importClassesFrom(sp,Polygons)
importClassesFrom(sp,SpatialPolygons)
importClassesFrom(sp,owin)
importFrom(grDevices,extendrange)
importFrom(grDevices,gray)
importFrom(grDevices,heat.colors)
importFrom(grDevices,xy.coords)
importFrom(graphics,image)
importFrom(graphics,lines)
importFrom(graphics,points)
importFrom(graphics,polygon)
importFrom(sp,Polygons)
importFrom(sp,SpatialPolygons)
importFrom(sp,coordinates)
importFrom(sp,plot)
importFrom(stats,cov2cor)
importFrom(stats,dist)
importFrom(stats,integrate)
importFrom(stats,pchisq)
importFrom(stats,pnorm)
useDynLib(polyCub, .registration = TRUE)
polyCub/NEWS.md 0000644 0001762 0000144 00000017327 13426776565 013024 0 ustar ligges users polyCub 0.7.1 (2019-02-07)
==========================
* Added a *getting started* `vignette("polyCub")`
(suggested by @wrathematics in
[openjournals/joss-reviews#1056](https://github.com/openjournals/joss-reviews/issues/1056)).
* fix minor compiler warning about missing `types` field in `R_CMethodDef`
(@wrathematics, [#1](https://github.com/bastistician/polyCub/issues/1)).
polyCub 0.7.0 (2018-10-11)
==========================
* Package **polyCub** no longer attaches package
[**sp**](https://CRAN.R-project.org/package=sp)
(moved from "Depends" to "Imports").
* The R code of the examples is no longer installed by default.
Use the `--example` flag of R CMD INSTALL to achieve that.
* The README now exemplifies the four different cubature rules.
polyCub 0.6.1 (2017-10-02)
==========================
* The exported C-function `polyCub_iso()` ...
* did not handle its `stop_on_error` argument correctly
(it would always stop on error).
* now detects non-finite `intrfr` function values and gives an
informative error message (rather than just reporting "abnormal
termination of integration routine").
* Package **polyCub** no longer strictly depends on package
[**spatstat**](https://CRAN.R-project.org/package=spatstat).
It is only required for `polyCub.midpoint()` and for polygon input of
class `"owin"`.
polyCub 0.6.0 (2017-05-24)
==========================
* Added full C-implementation of `polyCub.iso()`, which is exposed as
`"polyCub_iso"` for use by other R packages (notably future versions of
[**surveillance**](https://CRAN.R-project.org/package=surveillance))
via `LinkingTo: polyCub` and `#include `.
* Accommodate CRAN checks:
add missing import from **graphics**,
register native routines and disable symbol search
polyCub 0.5-2 (2015-02-25)
==========================
* `polyCub.midpoint()` works directly with input polygons of classes
`"gpc.poly"` and `"SpatialPolygons"`, since package **polyCub** now
registers corresponding `as.owin`-methods.
* `polyCub.exact.Gauss()` did not work if the `tristrip` of the
transformed input polygon contained degenerate triangles (spotted by
Ignacio Quintero).
* Line integration in `polyCub.iso()` could break due to division by zero
if the `center` point was part of the polygon boundary.
polyCub 0.5-1 (2014-10-24)
==========================
* Nodes and weights for `polyCub.SV()` were only cached up to `nGQ=59`,
not 60 as announced in version 0.5-0. Fixed that which also makes
examples truly run without **statmod**.
* In `polyCub.SV()`, the new special setting `f=NULL` means to only
compute nodes and weights.
* Internal changes to the `"gpc.poly"` converters to accommodate
[**spatstat**](https://CRAN.R-project.org/package=spatstat) 1.39-0.
polyCub 0.5-0 (2014-05-07)
==========================
* `polyCub.SV()` gained an argument `engine` to choose among available
implementations. The new and faster C-implementation is the default.
There should not be any numerical differences in the result of the
cubature.
* Package [**statmod**](https://CRAN.R-project.org/package=statmod) is no
longer strictly required (imported). Nodes and weights for
Gauss-Legendre quadrature in `polyCub.SV()` are now cached in the
**polyCub** package up to `nGQ=60`. **statmod**`::gauss.quad` is only
queried for a higher number of nodes.
polyCub 0.4-3 (2014-03-14)
==========================
* `polyCub.iso()` ...
* could not handle additional arguments for `integrate()` given in the
`control` list.
* now applies the `control` arguments also to the numerical
approximation of `intrfr`.
* The `checkintrfr()` function is exported and documented.
* Added a CITATION file.
polyCub 0.4-2 (2014-02-12)
==========================
* `plotpolyf()` ...
* gained an additional argument `print.args`, an optional list of
arguments passed to `print.trellis()` if `use.lattice=TRUE`.
* passed a *data frame* of coordinates to `f` instead of a matrix as
documented.
polyCub 0.4-1 (2013-12-05)
==========================
* This version solely fixes a missing NAMESPACE import to make
package **polyCub** again compatible with older versions of
[**spatstat**](https://CRAN.R-project.org/package=spatstat) (< 1.33-0).
polyCub 0.4-0 (2013-11-19)
==========================
INFRASTRUCTURE
--------------
* [**rgeos**](https://CRAN.R-project.org/package=rgeos) (and therefore the
GEOS library) is no longer strictly required (moved from "Imports" to
"Suggests").
* Added `coerce`-methods from `"Polygons"` (or `"SpatialPolygons"` or
`"Polygon"`) to `"owin"` (`as(..., "owin")`).
* S4-style `coerce`-methods between `"gpc.poly"` and `"Polygons"`/`"owin"`
have been removed from the package (since we no longer import the formal
class `"gpc.poly"` from **gpclib** or **rgeos**). However, there are two
new functions `gpc2owin` and `owin2gpc` similar to those dropped from
[**spatstat**](https://CRAN.R-project.org/package=spatstat) since
version 1.34-0.
* Moved `discpoly()` back to
[**surveillance**](https://CRAN.R-project.org/package=surveillance)
since it is only used there.
* The latter two changes cause
[**surveillance**](https://CRAN.R-project.org/package=surveillance)
version 1.6-0 to be incompatible with this new version of **polyCub**.
Appropriate modifications have been made in the new version 1.7-0 of
**surveillance**.
SPEED-UP `polyCub.SV()`
-----------------------
* thorough optimization of `polyCub.SV()`-related code resulted in about
27% speed-up:
* use `mapply()` instead of a `for`-loop
* avoid `cbind()`
* use `tcrossprod()`
* less object copying
MINOR CHANGES
-------------
* `xylist()` is now exported. It simply extracts polygon coordinates from
various spatial classes (with same unifying intention as `xy.coords()`).
* A `polyregion` of class `"SpatialPolygons"` of length more than 1 now
works in `polyCub`-methods.
* Use aspect ratio of 1 in `plotpolyf()`.
polyCub 0.3-1 (2013-08-22)
==========================
* This version solely fixes a few typos and a technical note from `R CMD
check` in the current R development version (also import packages into
the NAMESPACE which are listed in the "Depends" field).
polyCub 0.3-0 (2013-07-06)
==========================
* New cubature method `polyCub.iso()` specific to isotropic functions
(thanks to Emil Hedevang for the basic idea).
* New function `plotpolyf()` to plot a polygonal domain on top of an image
of a bivariate function.
* The package now depends on R >= 2.15.0 (for `.rowSums()`).
* The package no longer registers `"owin"` as an S4-class since we depend
on the **sp** package which does the job. This avoids a spurious warning
(in `.simpleDuplicateClass()`) upon package installation.
* In `discpoly()`, the argument `r` has been renamed to `radius`. This is
backward compatible by partial argument matching in old code.
polyCub 0.2-0 (2013-05-09)
==========================
* This is the initial version of the **polyCub** package mainly built on
functions previously maintained within the
[**surveillance**](https://CRAN.R-project.org/package=surveillance)
package. These methods for cubature of polygonal domains have been
outsourced into this separate **polyCub** package since they are of
general use for other packages as well.
* The **polyCub** package has more documentation and tests, avoids the use
of [**gpclib**](https://CRAN.R-project.org/package=gpclib) as far as
possible (using [**rgeos**](https://CRAN.R-project.org/package=rgeos)
instead), and solves a compatibility issue with package
[**maptools**](https://CRAN.R-project.org/package=maptools) (use
`setClass("owin")` instead of `setOldClass("owin")`).
polyCub/R/ 0000755 0001762 0000144 00000000000 13424543542 012076 5 ustar ligges users polyCub/R/polyCub.SV.R 0000644 0001762 0000144 00000034127 13424543542 014174 0 ustar ligges users ################################################################################
### polyCub.SV: Product Gauss Cubature over Polygonal Domains
###
### Copyright (C) 2009-2014,2017-2018 Sebastian Meyer
###
### This file is part of the R package "polyCub",
### free software under the terms of the GNU General Public License, version 2,
### a copy of which is available at https://www.R-project.org/Licenses/.
################################################################################
##' Product Gauss Cubature over Polygonal Domains
##'
##' Product Gauss cubature over polygons as proposed by
##' Sommariva and Vianello (2007).
##'
##' @inheritParams plotpolyf
##' @param f a two-dimensional real-valued function to be integrated over
##' \code{polyregion} (or \code{NULL} to only compute nodes and weights).
##' As its first argument it must take a coordinate matrix, i.e., a
##' numeric matrix with two columns, and it must return a numeric vector of
##' length the number of coordinates.
##' @param nGQ degree of the one-dimensional Gauss-Legendre quadrature rule
##' (default: 20) as implemented in function \code{\link[statmod]{gauss.quad}}
##' of package \pkg{statmod}. Nodes and weights up to \code{nGQ=60} are cached
##' in \pkg{polyCub}, for larger degrees \pkg{statmod} is required.
##' @param alpha base-line of the (rotated) polygon at \eqn{x = \alpha} (see
##' Sommariva and Vianello (2007) for an explication). If \code{NULL} (default),
##' the midpoint of the x-range of each polygon is chosen if no \code{rotation}
##' is performed, and otherwise the \eqn{x}-coordinate of the rotated point
##' \code{"P"} (see \code{rotation}). If \code{f} has its maximum value at the
##' origin \eqn{(0,0)}, e.g., the bivariate Gaussian density with zero mean,
##' \code{alpha = 0} is a reasonable choice.
##' @param rotation logical (default: \code{FALSE}) or a list of points
##' \code{"P"} and \code{"Q"} describing the preferred direction. If
##' \code{TRUE}, the polygon is rotated according to the vertices \code{"P"} and
##' \code{"Q"}, which are farthest apart (see Sommariva and Vianello, 2007). For
##' convex polygons, this rotation guarantees that all nodes fall inside the
##' polygon.
##' @param engine character string specifying the implementation to use.
##' Up to \pkg{polyCub} version 0.4-3, the two-dimensional nodes and weights
##' were computed by \R functions and these are still available by setting
##' \code{engine = "R"}.
##' The new C-implementation is now the default (\code{engine = "C"}) and
##' requires approximately 30\% less computation time.\cr
##' The special setting \code{engine = "C+reduce"} will discard redundant nodes
##' at (0,0) with zero weight resulting from edges on the base-line
##' \eqn{x = \alpha} or orthogonal to it.
##' This extra cleaning is only worth its cost for computationally intensive
##' functions \code{f} over polygons which really have some edges on the
##' baseline or parallel to the x-axis. Note that the old \R
##' implementation does not have such unset zero nodes and weights.
##' @param plot logical indicating if an illustrative plot of the numerical
##' integration should be produced.
##' @return The approximated value of the integral of \code{f} over
##' \code{polyregion}.\cr
##' In the case \code{f = NULL}, only the computed nodes and weights are
##' returned in a list of length the number of polygons of \code{polyregion},
##' where each component is a list with \code{nodes} (a numeric matrix with
##' two columns), \code{weights} (a numeric vector of length
##' \code{nrow(nodes)}), the rotation \code{angle}, and \code{alpha}.
##' @author Sebastian Meyer\cr
##' These R and C implementations of product Gauss cubature are based on the
##' original \acronym{MATLAB} implementation \code{polygauss} by Sommariva and
##' Vianello (2007), which is available under the GNU GPL (>=2) license from
##' \url{http://www.math.unipd.it/~alvise/software.html}.
##' @references
##' Sommariva, A. and Vianello, M. (2007):
##' Product Gauss cubature over polygons based on Green's integration formula.
##' \emph{BIT Numerical Mathematics}, \bold{47} (2), 441-453.\cr
##' DOI-Link: \url{https://doi.org/10.1007/s10543-007-0131-2}
##' @keywords math spatial
##' @family polyCub-methods
##' @importFrom graphics points
##' @example examples/setting.R
##' @example examples/polyCub.SV.R
##' @export
polyCub.SV <- function (polyregion, f, ...,
nGQ = 20, alpha = NULL, rotation = FALSE, engine = "C",
plot = FALSE)
{
polys <- xylist(polyregion) # transform to something like "owin$bdry"
# which means anticlockwise vertex order with
# first vertex not repeated
stopifnot(isScalar(nGQ), nGQ > 0,
is.null(alpha) || (isScalar(alpha) && !is.na(alpha)))
## COMPUTE NODES AND WEIGHTS OF 1D GAUSS QUADRATURE RULE.
## DEGREE "N" (as requested) (ORDER GAUSS PRIMITIVE)
nw_N <- gauss.quad(nGQ)
## DEGREE "M" = N+1 (ORDER GAUSS INTEGRATION)
nw_M <- gauss.quad(nGQ + 1)
## Special case f=NULL: compute and return nodes and weights only
if (is.null(f)) {
return(lapply(X = polys, FUN = polygauss, nw_MN = c(nw_M, nw_N),
alpha = alpha, rotation = rotation, engine = engine))
}
## Cubature over every single polygon of the "polys" list
f <- match.fun(f)
int1 <- function (poly) {
nw <- polygauss(poly, c(nw_M, nw_N), alpha, rotation, engine)
fvals <- f(nw$nodes, ...)
cubature_val <- sum(nw$weights * fvals)
## if (!isTRUE(all.equal(0, cubature_val))) {
## if ((1 - 2 * as.numeric(poly$hole)) * sign(cubature_val) == -1)
## warning("wrong sign if positive integral")
## }
cubature_val
}
respolys <- vapply(X=polys, FUN=int1, FUN.VALUE=0, USE.NAMES=FALSE)
int <- sum(respolys)
### ILLUSTRATION ###
if (plot) {
plotpolyf(polys, f, ..., use.lattice=FALSE)
for (i in seq_along(polys)) {
nw <- polygauss(polys[[i]], c(nw_M, nw_N), alpha, rotation, engine)
points(nw$nodes, cex=0.6, pch = i) #, col=1+(nw$weights<=0)
}
}
###################
int
}
## this wrapper provides a partially memoized version of
## unname(statmod::gauss.quad(n, kind="legendre"))
gauss.quad <- function (n)
{
if (n <= 61) { # results cached in R/sysdata.rda
.NWGL[[n]]
} else if (requireNamespace("statmod")) {
unname(statmod::gauss.quad(n = n, kind = "legendre"))
} else {
stop("package ", sQuote("statmod"), " is required for nGQ > 60")
}
}
##' Calculate 2D Nodes and Weights of the Product Gauss Cubature
##'
##' @param xy list with elements \code{"x"} and \code{"y"} containing the
##' polygon vertices in \emph{anticlockwise} order (otherwise the result of the
##' cubature will have a negative sign) with first vertex not repeated at the
##' end (like \code{owin.object$bdry}).
##' @param nw_MN unnamed list of nodes and weights of one-dimensional Gauss
##' quadrature rules of degrees \eqn{N} and \eqn{M=N+1} (as returned by
##' \code{\link[statmod]{gauss.quad}}): \code{list(s_M, w_M, s_N, w_N)}.
##' @inherit polyCub.SV params references
##' @keywords internal
##' @useDynLib polyCub, .registration = TRUE
polygauss <- function (xy, nw_MN, alpha = NULL, rotation = FALSE, engine = "C")
{
## POLYGON ROTATION
xyrot <- if (identical(FALSE, rotation)) {
if (is.null(alpha)) { # choose midpoint of x-range
xrange <- range(xy[["x"]])
alpha <- (xrange[1L] + xrange[2L]) / 2
}
angle <- 0
xy[c("x", "y")]
} else {
## convert to coordinate matrix
xy <- cbind(xy[["x"]], xy[["y"]], deparse.level=0)
## determine P and Q
if (identical(TRUE, rotation)) { # automatic choice of rotation angle
## such that for a convex polygon all nodes fall inside the polygon
QP <- vertexpairmaxdist(xy)
Q <- QP[1L,,drop=TRUE]
P <- QP[2L,,drop=TRUE]
} else if (is.list(rotation)) { # predefined rotation
P <- rotation$P
Q <- rotation$Q
stopifnot(is.vector(P, mode="numeric") && length(P) == 2L,
is.vector(Q, mode="numeric") && length(Q) == 2L)
stopifnot(any(P != Q))
rotation <- TRUE
} else {
stop("'rotation' must be logical or a list of points ",
"\"P\" and \"Q\"")
}
rotmat <- rotmatPQ(P,Q)
angle <- attr(rotmat, "angle")
if (is.null(alpha)) {
Prot <- rotmat %*% P
alpha <- Prot[1]
}
xyrot <- xy %*% t(rotmat) # = t(rotmat %*% t(xy))
## convert back to list
list(x = xyrot[,1L,drop=TRUE], y = xyrot[,2L,drop=TRUE])
}
## number of vertices
L <- length(xyrot[[1L]])
## COMPUTE 2D NODES AND WEIGHTS.
if (engine == "R") {
toIdx <- c(seq.int(2, L), 1L)
nwlist <- mapply(.polygauss.side,
xyrot[[1L]], xyrot[[2L]],
xyrot[[1L]][toIdx], xyrot[[2L]][toIdx],
MoreArgs = c(nw_MN, alpha),
SIMPLIFY = FALSE, USE.NAMES = FALSE)
nodes <- c(lapply(nwlist, "[[", 1L),
lapply(nwlist, "[[", 2L),
recursive=TRUE)
dim(nodes) <- c(length(nodes)/2, 2L)
weights <- unlist(lapply(nwlist, "[[", 3L),
recursive=FALSE, use.names=FALSE)
} else { # use C-implementation
## degrees of cubature and vector template for results
M <- length(nw_MN[[1L]])
N <- length(nw_MN[[3L]])
zerovec <- double(L*M*N)
## rock'n'roll
nwlist <- .C(C_polygauss,
as.double(xyrot[[1L]]), as.double(xyrot[[2L]]),
as.double(nw_MN[[1L]]), as.double(nw_MN[[2L]]),
as.double(nw_MN[[3L]]), as.double(nw_MN[[4L]]),
as.double(alpha),
as.integer(L), as.integer(M), as.integer(N),
x = zerovec, y = zerovec, w = zerovec)[c("x", "y", "w")]
nodes <- cbind(nwlist[[1L]], nwlist[[2L]], deparse.level=0)
weights <- nwlist[[3L]]
## remove unset nodes from edges on baseline or orthogonal to it
## (note that the R implementation does not return such redundant nodes)
if (engine == "C+reduce" && any(unset <- weights == 0)) {
nodes <- nodes[!unset,]
weights <- weights[!unset]
}
}
## back-transform rotated nodes by t(t(rotmat) %*% t(nodes))
## (inverse of rotation matrix is its transpose)
list(nodes = if (rotation) nodes %*% rotmat else nodes,
weights = weights, angle = angle, alpha = alpha)
}
## The working horse .polygauss.side below is an R translation
## of the original MATLAB implementation by Sommariva and Vianello (2007).
.polygauss.side <- function (x1, y1, x2, y2, s_loc, w_loc, s_N, w_N, alpha)
{
if ((x1 == alpha && x2 == alpha) || (y2 == y1))
## side lies on base-line or is orthogonal to it -> skip
return(NULL)
if (x2 == x1) { # side is parallel to base-line => degree N
s_loc <- s_N
w_loc <- w_N
}
half_pt_x <- (x1+x2)/2
half_length_x <- (x2-x1)/2
half_pt_y <- (y1+y2)/2
half_length_y <- (y2-y1)/2
## GAUSSIAN POINTS ON THE SIDE.
x_gauss_side <- half_pt_x + half_length_x * s_loc
y_gauss_side <- half_pt_y + half_length_y * s_loc
scaling_fact_minus <- (x_gauss_side - alpha) / 2
## construct nodes and weights: x and y coordinates ARE STORED IN MATRICES.
## A COUPLE WITH THE SAME INDEX IS A POINT, i.e. P_i=(x(k),y(k)).
## Return in an unnamed list of nodes_x, nodes_y, weights
## (there is no need for c(nodes_x) and c(weights))
list(
alpha + tcrossprod(scaling_fact_minus, s_N + 1), # degree_loc x N
rep.int(y_gauss_side, length(s_N)), # length: degree_loc*N
tcrossprod(half_length_y*scaling_fact_minus*w_loc, w_N) # degree_loc x N
)
}
## NOTE: The above .polygauss.side() function is already efficient R code.
## Passing via C only at this deep level (see below) turned out to be
## slower than staying with R! However, stepping into C already for
## looping over the edges in polygauss() improves the speed.
## ## @useDynLib polyCub C_polygauss_side
## .polygauss.side <- function (x1, y1, x2, y2, s_M, w_M, s_N, w_N, alpha)
## {
## if ((x1 == alpha && x2 == alpha) || (y2 == y1))
## ## side lies on base-line or is orthogonal to it -> skip
## return(NULL)
##
## parallel2baseline <- x2 == x1 # side is parallel to base-line => degree N
## M <- length(s_M)
## N <- length(s_N)
## loc <- if (parallel2baseline) N else M
## zerovec <- double(loc * N)
## .C(C_polygauss_side,
## as.double(x1), as.double(y1), as.double(x2), as.double(y2),
## as.double(if (parallel2baseline) s_N else s_M),
## as.double(if (parallel2baseline) w_N else w_M),
## as.double(s_N), as.double(w_N), as.double(alpha),
## as.integer(loc), as.integer(N),
## x = zerovec, y = zerovec, w = zerovec)[c("x", "y", "w")]
## }
##' @importFrom stats dist
vertexpairmaxdist <- function (xy)
{
## compute euclidean distance matrix
distances <- dist(xy)
size <- attr(distances, "Size")
## select two points with maximum distance
maxdistidx <- which.max(distances)
lowertri <- seq_along(distances) == maxdistidx
mat <- matrix(FALSE, size, size)
mat[lower.tri(mat)] <- lowertri
QPidx <- which(mat, arr.ind=TRUE, useNames=FALSE)[1L,]
xy[QPidx,]
}
rotmatPQ <- function (P, Q)
{
direction_axis <- (Q-P) / vecnorm(Q-P)
## determine rotation angle [radian]
rot_angle_x <- acos(direction_axis[1L])
rot_angle_y <- acos(direction_axis[2L])
rot_angle <- if (rot_angle_y <= pi/2) {
if (rot_angle_x <= pi/2) -rot_angle_y else rot_angle_y
} else {
if (rot_angle_x <= pi/2) pi-rot_angle_y else rot_angle_y
}
## rotation matrix
rot_matrix <- diag(cos(rot_angle), nrow=2L)
rot_matrix[2:3] <- c(-1,1) * sin(rot_angle) # clockwise rotation
structure(rot_matrix, angle=rot_angle)
}
polyCub/R/polyCub.iso.R 0000644 0001762 0000144 00000025217 13357335027 014440 0 ustar ligges users ################################################################################
### polyCub.iso: Cubature of Isotropic Functions over Polygonal Domains
###
### Copyright (C) 2013-2018 Sebastian Meyer
###
### This file is part of the R package "polyCub",
### free software under the terms of the GNU General Public License, version 2,
### a copy of which is available at https://www.R-project.org/Licenses/.
################################################################################
#' Cubature of Isotropic Functions over Polygonal Domains
#'
#' \code{polyCub.iso} numerically integrates a radially symmetric function
#' \eqn{f(x,y) = f_r(||(x,y)-\boldsymbol{\mu}||)}{f(x,y) = f_r(||(x,y)-\mu||)},
#' with \eqn{\mu} being the center of isotropy, over a polygonal domain.
#' It internally approximates a line integral along the polygon boundary using
#' \code{\link{integrate}}. The integrand requires the antiderivative of
#' \eqn{r f_r(r)}), which should be supplied as argument \code{intrfr}
#' (\code{f} itself is only required if \code{check.intrfr=TRUE}).
#' The two-dimensional integration problem thereby reduces to an efficient
#' adaptive quadrature in one dimension.
#' If \code{intrfr} is not available analytically, \code{polyCub.iso} can use a
#' numerical approximation (meaning \code{integrate} within \code{integrate}),
#' but the general-purpose cubature method \code{\link{polyCub.SV}} might be
#' more efficient in this case.
#' See Meyer and Held (2014, Supplement B, Section 2.4) for mathematical
#' details.
#'
#' @inheritParams plotpolyf
#' @param intrfr a \code{function(R, ...)}, which implements the (analytical)
#' antiderivative of \eqn{r f_r(r)} from 0 to \code{R}. The first argument
#' must be vectorized but not necessarily named \code{R}.\cr
#' If \code{intrfr} is missing, it will be approximated numerically via
#' \code{\link{integrate}(function(r, ...)
#' r * f(cbind(x0 + r, y0), ...), 0, R, ..., control=control)},
#' where \code{c(x0, y0)} is the \code{center} of isotropy.
#' Note that \code{f} will \emph{not} be checked for isotropy.
#' @param ... further arguments for \code{f} or \code{intrfr}.
#' @param center numeric vector of length 2, the center of isotropy.
#' @param control list of arguments passed to \code{\link{integrate}}, the
#' quadrature rule used for the line integral along the polygon boundary.
#' @param check.intrfr logical (or numeric vector) indicating if
#' (for which \code{r}'s) the supplied \code{intrfr} function should be
#' checked against a numeric approximation. This check requires \code{f}
#' to be specified. If \code{TRUE}, the set of test
#' \code{r}'s defaults to a \code{\link{seq}} of length 20 from 1 to
#' the maximum absolute x or y coordinate of any edge of the \code{polyregion}.
#' @param plot logical indicating if an image of the function should be plotted
#' together with the polygonal domain, i.e.,
#' \code{\link{plotpolyf}(polyregion, f, \dots)}.
#' @return The approximate integral of the isotropic function
#' \code{f} over \code{polyregion}.\cr
#' If the \code{intrfr} function is provided (which is assumed to be exact), an
#' upper bound for the absolute integration error is appended to the result as
#' attribute \code{"abs.error"}. It equals the sum of the absolute errors
#' reported by all \code{\link{integrate}} calls
#' (there is one for each edge of \code{polyregion}).
#' @author Sebastian Meyer
#'
#' The basic mathematical formulation of this efficient integration for radially
#' symmetric functions was ascertained with great support by
#' Emil Hedevang (2013), Dept. of Mathematics, Aarhus University, Denmark.
#' @references
#' Hedevang, E. (2013). Personal communication at the Summer School on Topics in
#' Space-Time Modeling and Inference (May 2013, Aalborg, Denmark).
#'
#' Meyer, S. and Held, L. (2014).
#' Power-law models for infectious disease spread.
#' \emph{The Annals of Applied Statistics}, \bold{8} (3), 1612-1639.\cr
#' DOI-Link: \url{https://doi.org/10.1214/14-AOAS743},
#' \href{https://arxiv.org/abs/1308.5115}{arXiv:1308.5115}
#' @seealso
#' \code{system.file("include", "polyCubAPI.h", package = "polyCub")}
#' for a full C-implementation of this cubature method (for a \emph{single}
#' polygon). The corresponding C-routine \code{polyCub_iso} can be used by
#' other \R packages, notably \pkg{surveillance}, via \samp{LinkingTo: polyCub}
#' (in the \file{DESCRIPTION}) and \samp{#include } (in suitable
#' \file{/src} files). Note that the \code{intrfr} function must then also be
#' supplied as a C-routine. An example can be found in the package tests.
#' @keywords math spatial
#' @family polyCub-methods
#' @example examples/polyCub.iso.R
#' @importFrom stats integrate
#' @export
polyCub.iso <- function (polyregion, f, intrfr, ..., center,
control = list(), check.intrfr = FALSE, plot = FALSE)
{
polys <- xylist(polyregion) # transform to something like "owin$bdry"
# which means anticlockwise vertex order with
# first vertex not repeated
getError <- !missing(intrfr) # can't estimate error of double approximation
center <- as.vector(center, mode = "numeric")
stopifnot(length(center) == 2L, is.finite(center))
## check 'intrfr' function
rs <- if (isTRUE(check.intrfr)) {
seq(1, max(abs(unlist(lapply(polys, "[", c("x","y"))))), length.out=20L)
} else if (identical(check.intrfr, FALSE)) {
numeric(0L)
} else {
check.intrfr
}
intrfr <- checkintrfr(intrfr, f, ..., center=center, control=control, rs=rs)
## plot polygon and function image
if (plot) plotpolyf(polys, f, ...)
## do the cubature over all polygons of the 'polys' list
.polyCub.iso(polys, intrfr, ..., center=center,
control=control, .witherror=getError)
}
##' Check the Integral of \eqn{r f_r(r)}
##'
##' This function is auxiliary to \code{\link{polyCub.iso}}.
##' The (analytical) integral of \eqn{r f_r(r)} from 0 to \eqn{R} is checked
##' against a numeric approximation using \code{\link{integrate}} for various
##' values of the upper bound \eqn{R}. A warning is issued if inconsistencies
##' are found.
##'
##' @inheritParams polyCub.iso
##' @param rs numeric vector of upper bounds for which to check the validity of
##' \code{intrfr}. If it has length 0 (default), no checks are performed.
##' @param tolerance of \code{\link{all.equal.numeric}} when comparing
##' \code{intrfr} results with numerical integration. Defaults to the
##' relative tolerance used for \code{integrate}.
##' @return The \code{intrfr} function. If it was not supplied, its quadrature
##' version using \code{integrate} is returned.
##' @importFrom stats integrate
##' @export
checkintrfr <- function (intrfr, f, ..., center, control = list(),
rs = numeric(0L), tolerance = control$rel.tol)
{
doCheck <- length(rs) > 0L
if (!missing(f)) {
f <- match.fun(f)
rfr <- function (r, ...)
r * f(cbind(center[1L]+r, center[2L], deparse.level=0L), ...)
quadrfr1 <- function (R, ...) integrate(rfr, 0, R, ...)$value
if (length(control))
body(quadrfr1)[[2L]] <- as.call(c(as.list(body(quadrfr1)[[2L]]),
control))
quadrfr <- function (R, ...)
vapply(X = R, FUN = quadrfr1, FUN.VALUE = 0, ..., USE.NAMES = FALSE)
if (missing(intrfr)) {
return(quadrfr)
} else if (doCheck) {
cat("Checking 'intrfr' against a numeric approximation ... ")
stopifnot(is.vector(rs, mode="numeric"))
if (is.null(tolerance))
tolerance <- eval(formals(integrate)$rel.tol)
ana <- intrfr(rs, ...)
num <- quadrfr(rs, ...)
if (!isTRUE(comp <- all.equal(num, ana, tolerance=tolerance))) {
cat("\n->", comp, "\n")
warning("'intrfr' might be incorrect: ", comp)
} else cat("OK\n")
}
} else if (doCheck) {
stop("numerical verification of 'intrfr' requires 'f'")
}
match.fun(intrfr)
}
##' @description
##' \code{.polyCub.iso} is a \dQuote{bare-bone} version of \code{polyCub.iso}.
##' @rdname polyCub.iso
##' @param polys something like \code{owin$bdry}, but see \code{\link{xylist}}.
##' @param .witherror logical indicating if an upper bound for the absolute
##' integration error should be attached as an attribute to the result?
##' @export
.polyCub.iso <- function (polys, intrfr, ..., center,
control = list(), .witherror = FALSE)
{
ints <- lapply(polys, polyCub1.iso,
intrfr, ..., center=center,
control=control, .witherror=.witherror)
if (.witherror) {
res <- sum(vapply(X=ints, FUN="[", FUN.VALUE=0, 1L, USE.NAMES=FALSE))
attr(res, "abs.error") <-
sum(vapply(X=ints, FUN="[", FUN.VALUE=0, 2L, USE.NAMES=FALSE))
res
} else {
sum(unlist(ints, recursive=FALSE, use.names=FALSE))
}
}
## cubature method for a single polygon
polyCub1.iso <- function (poly, intrfr, ..., center,
control = list(), .witherror = TRUE)
{
xy <- cbind(poly[["x"]], poly[["y"]], deparse.level=0L)
nedges <- nrow(xy)
intedges <- erredges <- numeric(nedges)
for (i in seq_len(nedges)) {
v0 <- xy[i, ] - center
v1 <- xy[if (i==nedges) 1L else i+1L, ] - center
int <- lineInt(v0, v1, intrfr, ..., control=control)
intedges[i] <- int$value
erredges[i] <- int$abs.error
}
int <- sum(intedges)
## if (!is.null(poly$hole) && !isTRUE(all.equal(0, int))) {
## if ((1 - 2 * as.numeric(poly$hole)) * sign(int) == -1)
## warning("wrong sign if positive integral")
## }
if (.witherror) {
c(int, sum(erredges))
} else {
int
}
}
## line integral for one edge
##' @importFrom stats integrate
lineInt <- function (v0, v1, intrfr, ..., control = list())
{
d <- v1 - v0
num <- v1[2L]*v0[1L] - v1[1L]*v0[2L] # = d[2]*p[,1] - d[1]*p[,2]
# for any point p on the edge
if (num == 0) { # i.e., if 'center' is part of this polygon edge
return(list(value = 0, abs.error = 0))
}
integrand <- function (t) {
## get the points on the edge corresponding to t
p <- cbind(v0[1L] + t*d[1L], v0[2L] + t*d[2L], deparse.level=0L)
norm2 <- .rowSums(p^2, length(t), 2L)
ints <- intrfr(sqrt(norm2), ...)
##ints[is.infinite(ints)] <- 1e300
num * ints / norm2
}
if (length(control)) { # use slower do.call()-construct
do.call("integrate", c(list(integrand, 0, 1), control))
} else {
integrate(integrand, 0, 1)
}
}
polyCub/R/polyCub.R 0000644 0001762 0000144 00000004405 13424543542 013641 0 ustar ligges users ################################################################################
### polyCub: Wrapper Function for the Various Cubature Methods
###
### Copyright (C) 2009-2013,2019 Sebastian Meyer
###
### This file is part of the R package "polyCub",
### free software under the terms of the GNU General Public License, version 2,
### a copy of which is available at https://www.R-project.org/Licenses/.
################################################################################
#' Wrapper Function for the Various Cubature Methods
#'
#' The wrapper function \code{polyCub} can be used to call specific cubature
#' methods via its \code{method} argument. It calls \code{\link{polyCub.SV}}
#' by default, which implements general-purpose product Gauss cubature.
#'
#' @inheritParams plotpolyf
#' @param f a two-dimensional real-valued function to be integrated over
#' \code{polyregion}. As its first argument it must take a coordinate matrix,
#' i.e., a numeric matrix with two columns, and it must return a numeric vector
#' of length the number of coordinates.\cr
#' For the \code{"exact.Gauss"} \code{method},
#' \code{f} is ignored since it is specific to the bivariate normal density.
#' @param method choose one of the implemented cubature methods (partial
#' argument matching is applied), see \code{help("\link{polyCub-package}")}
#' for an overview. Defaults to using product Gauss cubature
#' implemented in \code{\link{polyCub.SV}}.
#' @param ... arguments of \code{f} or of the specific \code{method}.
#' @param plot logical indicating if an illustrative plot of the numerical
#' integration should be produced.
#' @return The approximated integral of \code{f} over \code{polyregion}.
#' @seealso Details and examples in the \code{vignette("polyCub")}
#' and on the method-specific help pages.
#' @family polyCub-methods
#' @keywords math spatial
#' @export
polyCub <- function (polyregion, f,
method = c("SV", "midpoint", "iso", "exact.Gauss"), ...,
plot = FALSE)
{
method <- match.arg(method)
cl <- match.call()
cl$method <- NULL
cl[[1]] <- call("::", as.name("polyCub"),
as.name(paste("polyCub", method, sep=".")))
if (method == "exact.Gauss") cl$f <- NULL
int <- eval(cl, parent.frame())
int
}
polyCub/R/sysdata.rda 0000644 0001762 0000144 00000044644 12375637762 014266 0 ustar ligges users 7zXZ i"6 ! X|Ih] )TW"nRʟ)'dz$&}T1}} IL9d]hϐf Qcv\ DÿNրGsFcw}8=ɥ""(U=b4Yq:>q x6J6LUB4x}0&[g"QCBb9,xgeދ|`$ٔ'Nz#3<@?4kw,as