Functions available to any gadget3 model

Details

g3_env is the top-level environment that any gadget3 model uses, populated with utility functions.

NB: Several functions have _vec variants. Due to TMB limitations these should be used when you have a vector not scalar input.

ADREPORT

TMB's ADREPORT function. See sdreport documentation

as_integer

C++ compatible equivalent to as.integer

as.numeric

R as.numeric or TMB asDouble

assert_msg

C++/R function that ensures expression is true, or stops model.

assert_msg(x > 0, "x must be positive")

avoid_zero / avoid_zero_vec

Adds small value to input to ensure output is never zero

bounded / bounded_vec

Ensures x is within limits a & b.

bounded_vec(x, 100, 1000)

g3_matrix_vec

Apply matrix transformation tf to vector vec, return resultant vector.

g3_matrix_vec(tf, vec)

lgamma_vec

Vector equivalent of lgamma

logspace_add / logspace_add_vec

TMB's logspace_add, essentially a differentiable version of pmax.

logspace_minmax_vec

Differentiable equivalent of pmax(pmin(vec, upper), lower). scale influences the sharpness of inflection points at lower & upper, should be ~1e5, depending on ranges of input values.

logspace_minmax_vec(vec, lower, upper, scale)

normalize_vec

Divide vector a by it's sum, i.e. so it now sums to 1

nvl

Return first non-null argument. NB: No C++ implementation.

pow_vec

Vector equivalent of ^

ratio_add_vec

Sum orig_vec & new_vec according to ratio of orig_amount & new_amount

nonconform_add / nonconform_mult / nonconform_div / nonconform_div_avz

Scalar sum/multiply/divide 2 non-conforming arrays, by treating the latter as a vector and repeating as many times as required. Results will be structured identically to the first array.

nonconform_div_avz(x, y) is equivalent to nonconform_div(x, avoid_zero_vec(y))

REPORT

TMB's REPORT function.

REprintf

Equivalent of RCpp REprintf

Rprintf

Equivalent of RCpp Rprintf

Examples

## avoid_zero / avoid_zero_vec
g3_eval(quote( c( avoid_zero(0), avoid_zero(10) ) ))
#> [1] 6.931472e-04 1.000000e+01
g3_eval(quote( avoid_zero_vec(0:5) ))
#> [1] 0.0006931472 1.0000000000 2.0000000000 3.0000000000 4.0000000000
#> [6] 5.0000000000

## bounded / bounded_vec
curve(g3_eval(quote( bounded(x, 100, 200) ), x = x), -100, 100)


## logspace_add
curve(g3_eval(quote( logspace_add(x, 10) ), x = x), 0, 40)


## logspace_minmax_vec
curve(g3_eval(quote( logspace_minmax_vec(x, 30, 60, scale = 1e5)), x = x), 0, 100)


## normalize_vec
g3_eval(quote( normalize_vec(c( 4, 4, 8, 2 )) ))
#> [1] 0.2222222 0.2222222 0.4444444 0.1111111

## nonconform_mult
g3_eval(quote( nonconform_mult(
    array(seq(0, 4*5*6), dim = c(4,5,6)),
    c(1e1, 1e2, 1e3, 1e4)) ))
#> , , 1
#> 
#>       [,1]  [,2]   [,3]   [,4]   [,5]
#> [1,]     0    40     80    120    160
#> [2,]   100   500    900   1300   1700
#> [3,]  2000  6000  10000  14000  18000
#> [4,] 30000 70000 110000 150000 190000
#> 
#> , , 2
#> 
#>        [,1]   [,2]   [,3]   [,4]   [,5]
#> [1,]    200    240    280    320    360
#> [2,]   2100   2500   2900   3300   3700
#> [3,]  22000  26000  30000  34000  38000
#> [4,] 230000 270000 310000 350000 390000
#> 
#> , , 3
#> 
#>        [,1]   [,2]   [,3]   [,4]   [,5]
#> [1,]    400    440    480    520    560
#> [2,]   4100   4500   4900   5300   5700
#> [3,]  42000  46000  50000  54000  58000
#> [4,] 430000 470000 510000 550000 590000
#> 
#> , , 4
#> 
#>        [,1]   [,2]   [,3]   [,4]   [,5]
#> [1,]    600    640    680    720    760
#> [2,]   6100   6500   6900   7300   7700
#> [3,]  62000  66000  70000  74000  78000
#> [4,] 630000 670000 710000 750000 790000
#> 
#> , , 5
#> 
#>        [,1]   [,2]   [,3]   [,4]   [,5]
#> [1,]    800    840    880    920    960
#> [2,]   8100   8500   8900   9300   9700
#> [3,]  82000  86000  90000  94000  98000
#> [4,] 830000 870000 910000 950000 990000
#> 
#> , , 6
#> 
#>         [,1]    [,2]    [,3]    [,4]    [,5]
#> [1,]    1000    1040    1080    1120    1160
#> [2,]   10100   10500   10900   11300   11700
#> [3,]  102000  106000  110000  114000  118000
#> [4,] 1030000 1070000 1110000 1150000 1190000
#> 

## nonconform_div_avz
g3_eval(quote( nonconform_div_avz(
    array(seq(0, 4*5*6), dim = c(4,5,6)),
    c(1e1, 1e2, 0, 1e4)) ))
#> , , 1
#> 
#>           [,1]      [,2]       [,3]       [,4]       [,5]
#> [1,]    0.0000    0.4000     0.8000     1.2000     1.6000
#> [2,]    0.0100    0.0500     0.0900     0.1300     0.1700
#> [3,] 2885.3901 8656.1702 14426.9504 20197.7306 25968.5107
#> [4,]    0.0003    0.0007     0.0011     0.0015     0.0019
#> 
#> , , 2
#> 
#>            [,1]       [,2]       [,3]       [,4]       [,5]
#> [1,]     2.0000     2.4000     2.8000     3.2000     3.6000
#> [2,]     0.2100     0.2500     0.2900     0.3300     0.3700
#> [3,] 31739.2909 37510.0711 43280.8512 49051.6314 54822.4116
#> [4,]     0.0023     0.0027     0.0031     0.0035     0.0039
#> 
#> , , 3
#> 
#>            [,1]       [,2]       [,3]       [,4]       [,5]
#> [1,]     4.0000     4.4000     4.8000     5.2000     5.6000
#> [2,]     0.4100     0.4500     0.4900     0.5300     0.5700
#> [3,] 60593.1917 66363.9719 72134.7520 77905.5322 83676.3124
#> [4,]     0.0043     0.0047     0.0051     0.0055     0.0059
#> 
#> , , 4
#> 
#>            [,1]       [,2]        [,3]        [,4]        [,5]
#> [1,]     6.0000     6.4000      6.8000      7.2000      7.6000
#> [2,]     0.6100     0.6500      0.6900      0.7300      0.7700
#> [3,] 89447.0925 95217.8727 100988.6529 106759.4330 112530.2132
#> [4,]     0.0063     0.0067      0.0071      0.0075      0.0079
#> 
#> , , 5
#> 
#>             [,1]        [,2]        [,3]        [,4]        [,5]
#> [1,]      8.0000      8.4000      8.8000      9.2000      9.6000
#> [2,]      0.8100      0.8500      0.8900      0.9300      0.9700
#> [3,] 118300.9934 124071.7735 129842.5537 135613.3338 141384.1140
#> [4,]      0.0083      0.0087      0.0091      0.0095      0.0099
#> 
#> , , 6
#> 
#>             [,1]        [,2]        [,3]        [,4]        [,5]
#> [1,]     10.0000     10.4000     10.8000     11.2000     11.6000
#> [2,]      1.0100      1.0500      1.0900      1.1300      1.1700
#> [3,] 147154.8942 152925.6743 158696.4545 164467.2347 170238.0148
#> [4,]      0.0103      0.0107      0.0111      0.0115      0.0119
#> 
g3_eval(quote( nonconform_div(
    array(seq(0, 4*5*6), dim = c(4,5,6)),
    avoid_zero(c(1e1, 1e2, 0, 1e4))) ))
#> , , 1
#> 
#>           [,1]      [,2]       [,3]       [,4]       [,5]
#> [1,]    0.0000    0.4000     0.8000     1.2000     1.6000
#> [2,]    0.0100    0.0500     0.0900     0.1300     0.1700
#> [3,] 2885.3901 8656.1702 14426.9504 20197.7306 25968.5107
#> [4,]    0.0003    0.0007     0.0011     0.0015     0.0019
#> 
#> , , 2
#> 
#>            [,1]       [,2]       [,3]       [,4]       [,5]
#> [1,]     2.0000     2.4000     2.8000     3.2000     3.6000
#> [2,]     0.2100     0.2500     0.2900     0.3300     0.3700
#> [3,] 31739.2909 37510.0711 43280.8512 49051.6314 54822.4116
#> [4,]     0.0023     0.0027     0.0031     0.0035     0.0039
#> 
#> , , 3
#> 
#>            [,1]       [,2]       [,3]       [,4]       [,5]
#> [1,]     4.0000     4.4000     4.8000     5.2000     5.6000
#> [2,]     0.4100     0.4500     0.4900     0.5300     0.5700
#> [3,] 60593.1917 66363.9719 72134.7520 77905.5322 83676.3124
#> [4,]     0.0043     0.0047     0.0051     0.0055     0.0059
#> 
#> , , 4
#> 
#>            [,1]       [,2]        [,3]        [,4]        [,5]
#> [1,]     6.0000     6.4000      6.8000      7.2000      7.6000
#> [2,]     0.6100     0.6500      0.6900      0.7300      0.7700
#> [3,] 89447.0925 95217.8727 100988.6529 106759.4330 112530.2132
#> [4,]     0.0063     0.0067      0.0071      0.0075      0.0079
#> 
#> , , 5
#> 
#>             [,1]        [,2]        [,3]        [,4]        [,5]
#> [1,]      8.0000      8.4000      8.8000      9.2000      9.6000
#> [2,]      0.8100      0.8500      0.8900      0.9300      0.9700
#> [3,] 118300.9934 124071.7735 129842.5537 135613.3338 141384.1140
#> [4,]      0.0083      0.0087      0.0091      0.0095      0.0099
#> 
#> , , 6
#> 
#>             [,1]        [,2]        [,3]        [,4]        [,5]
#> [1,]     10.0000     10.4000     10.8000     11.2000     11.6000
#> [2,]      1.0100      1.0500      1.0900      1.1300      1.1700
#> [3,] 147154.8942 152925.6743 158696.4545 164467.2347 170238.0148
#> [4,]      0.0103      0.0107      0.0111      0.0115      0.0119
#>