-
Notifications
You must be signed in to change notification settings - Fork 0
/
polyCub.Rmd
320 lines (248 loc) · 11.8 KB
/
polyCub.Rmd
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
---
title: "Getting started with polyCub"
author: "Sebastian Meyer"
date: "`r Sys.Date()`"
output:
markdown::html_format:
options:
smartypants: false
---
<!--
%\VignetteIndexEntry{Getting started with polyCub}
%\VignetteEngine{knitr::knitr}
%\VignetteDepends{polyCub, spatstat.geom, markdown}
%% ^ implied by the engine
-->
```{R setup, include = FALSE, purl = FALSE}
knitr::opts_chunk$set(collapse = TRUE, comment = "#>",
fig.width = 4, fig.height = 3.5)
## maintainer-mode options for building the vignette:
if (USE_GPCLIB <- identical(Sys.getenv("R_GPCLIBPERMIT"), "true")) {
if (!requireNamespace("gpclib")) # unavailable in --as-cran checks
USE_GPCLIB <- FALSE
}
if (identical(Sys.getenv("NOT_CRAN"), "true")) {
stopifnot(nzchar(Sys.which("pngquant"))) # to notice when building
knitr::knit_hooks$set(pngquant = knitr::hook_pngquant)
knitr::opts_chunk$set(pngquant = "--speed=1")
DO_BENCHMARK <- USE_GPCLIB
} else {
DO_BENCHMARK <- FALSE
}
if (DO_BENCHMARK)
stopifnot(requireNamespace("microbenchmark"))
```
<img src="../man/figures/logo.png" align="right" alt="" width="120" />
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.geom**](https://CRAN.R-project.org/package=spatstat.geom),
* `"gpc.poly"` from [**gpclib**](https://github.com/rdpeng/gpclib/),
* `"SpatialPolygons"` from package [**sp**](https://CRAN.R-project.org/package=sp), and
* `"(MULTI)POLYGON"` from package [**sf**](https://CRAN.R-project.org/package=sf).
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 (for exterior boundaries)
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()` (*currently disabled*): Accurate 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}
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, <https://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:
```{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.geom**.
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.geom**'s "owin" class:
```{R, message = FALSE}
library("spatstat.geom")
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*,
<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:
```{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 <polyCubAPI.h>`.
This requires the `intrfr()` function to be implemented in C as well. See
<https://github.com/bastistician/polyCub/blob/master/tests/polyiso_powerlaw.c>
for an example.
### 4. Integration of the *bivariate Gaussian density*: `polyCub.exact.Gauss()`
*This cubature method is currently disabled*
([#2](https://github.com/bastistician/polyCub/issues/2)).
It requires polygon triangulation originally performed using `tristrip()`
from the [**gpclib**](https://github.com/rdpeng/gpclib/) package;
unfortunately, it has become unavailable from mainstream repositories.
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.
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 obtained:
```{R, purl = FALSE, eval = USE_GPCLIB}
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 quite 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 = DO_BENCHMARK}
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 = 9,
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 = DO_BENCHMARK}
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}$.