metric_graph {MetricGraph} | R Documentation |
Metric graph
Description
Class representing a general metric graph.
Details
A graph object created from vertex and edge matrices, or from an
sp::SpatialLines
object where each line is representing and edge. For more details,
see the vignette:
vignette("metric_graph", package = "MetricGraph")
Value
Object of R6Class
for creating metric graphs.
Public fields
V
Matrix with positions in Euclidean space of the vertices of the graph.
nV
The number of vertices.
E
Matrix with the edges of the graph, where each row represents an edge,
E[i,1]
is the vertex at the start of the ith edge andE[i,2]
is the vertex at the end of the edge.nE
The number of edges.
edge_lengths
Vector with the lengths of the edges in the graph.
C
Constraint matrix used to set Kirchhoff constraints.
CoB
Change-of-basis object used for Kirchhoff constraints.
PtV
Vector with the indices of the vertices which are observation locations.
mesh
Mesh object used for plotting.
edges
The coordinates of the edges in the graph.
vertices
The coordinates of the vertices in the graph, along with several attributes.
geo_dist
Geodesic distances between the vertices in the graph.
res_dist
Resistance distances between the observation locations.
Laplacian
The weighted graph Laplacian of the vertices in the graph. The weights are given by the edge lengths.
characteristics
List with various characteristics of the graph.
Methods
Public methods
Method new()
Create a new metric_graph
object.
Usage
metric_graph$new( edges = NULL, V = NULL, E = NULL, vertex_unit = NULL, length_unit = vertex_unit, edge_weights = 1, kirchhoff_weights = NULL, longlat = FALSE, crs = NULL, proj4string = NULL, which_longlat = "sp", project = FALSE, project_data = FALSE, which_projection = "Winkel tripel", tolerance = list(vertex_vertex = 0.001, vertex_edge = 0.001, edge_edge = 0), check_connected = TRUE, remove_deg2 = FALSE, merge_close_vertices = TRUE, factor_merge_close_vertices = 1, remove_circles = TRUE, verbose = 1, lines = deprecated() )
Arguments
edges
A list containing coordinates as
m x 2
matrices (that is, ofmatrix
type) or m x 2 data frames (data.frame
type) of sequence of points connected by straightlines. Alternatively, you can also prove an object of typeSpatialLinesDataFrame
orSpatialLines
(fromsp
package) orMULTILINESTRING
(fromsf
package).V
n x 2 matrix with Euclidean coordinates of the n vertices.
E
m x 2 matrix where each row represents one of the m edges.
vertex_unit
The unit in which the vertices are specified. The options are 'degrees' (the great circle distance in km), 'km', 'm' and 'miles'. The default is
NULL
, which means no unit. However, if you setlength_unit
, you need to setvertex_unit
.length_unit
The unit in which the lengths will be computed. The options are 'km', 'm' and 'miles'. The default is
vertex_unit
. Observe that ifvertex_unit
isNULL
,length_unit
can only beNULL
. Ifvertex_unit
is 'degrees', then the default value forlength_unit
is 'km'.edge_weights
Either a number, a numerical vector with length given by the number of edges, providing the edge weights, or a
data.frame
with the number of rows being equal to the number of edges, where each row gives a vector of weights to its corresponding edge. Can be changed by using theset_edge_weights()
method.kirchhoff_weights
If non-null, the name (or number) of the column of
edge_weights
that contain the Kirchhoff weights. Must be equal to 1 (orTRUE
) in caseedge_weights
is a single number and those are the Kirchhoff weights.longlat
If
TRUE
, then it is assumed that the coordinates are given. in Longitude/Latitude and that distances should be computed in meters. IfTRUE
it takes precedence oververtex_unit
andlength_unit
, and is equivalent tovertex_unit = 'degrees'
andlength_unit = 'm'
.crs
Coordinate reference system to be used in case
longlat
is set toTRUE
andwhich_longlat
issf
. Object of class crs. The default issf::st_crs(4326)
.proj4string
Projection string of class CRS-class to be used in case
longlat
is set toTRUE
andwhich_longlat
issp
. The default issp::CRS("+proj=longlat +datum=WGS84")
.which_longlat
Compute the distance using which package? The options are
sp
andsf
. The default issp
.project
If
longlat
isTRUE
should a projection be used to compute the distances to be used for the tolerances (seetolerance
below)? The default isFALSE
. WhenTRUE
, the construction of the graph is faster.project_data
If
longlat
isTRUE
should the vertices be project to planar coordinates? The default isFALSE
. WhenTRUE
, the construction of the graph is faster.which_projection
Which projection should be used in case
project
isTRUE
? The options areRobinson
,Winkel tripel
or a proj4string. The default isWinkel tripel
.tolerance
List that provides tolerances during the construction of the graph:
-
vertex_vertex
Vertices that are closer than this number are merged (default = 1e-7). -
vertex_edge
If a vertex at the end of one edge is closer than this number to another edge, this vertex is connected to that edge (default = 1e-7). Previouslyvertex_line
, which is now deprecated. -
edge_edge
If two edges at some point are closer than this number, a new vertex is added at that point and the two edges are connected (default = 0). -
vertex_line
, Deprecated. Usevertex_edge
instead. -
line_line
, Deprecated. Useedge_edge
instead.
In case
longlat = TRUE
, the tolerances are given inlength_unit
.-
check_connected
If
TRUE
, it is checked whether the graph is connected and a warning is given if this is not the case.remove_deg2
Set to
TRUE
to remove all vertices of degree 2 in the initialization. Default isFALSE
.merge_close_vertices
should an additional step to merge close vertices be done?
factor_merge_close_vertices
Which factor to be multiplied by tolerance
vertex_vertex
when merging close vertices at the additional step?remove_circles
All circlular edges with a length smaller than this number are removed. If
TRUE
, thevertex_vertex
tolerance will be used. IfFALSE
, no circles will be removed.verbose
Print progress of graph creation. There are 3 levels of verbose, level 0, 1 and 2. In level 0, no messages are printed. In level 1, only messages regarding important steps are printed. Finally, in level 2, messages detailing all the steps are printed. The default is 1.
lines
Details
A graph object can be initialized in two ways. The first method
is to specify V and E. In this case, all edges are assumed to be straight
lines. The second option is to specify the graph via the lines
input.
In this case, the vertices are set by the end points of the lines.
Thus, if two lines are intersecting somewhere else, this will not be
viewed as a vertex.
Returns
A metric_graph
object.
Method set_edge_weights()
Sets the edge weights
Usage
metric_graph$set_edge_weights( weights = rep(1, self$nE), kirchhoff_weights = NULL )
Arguments
weights
Either a number, a numerical vector with length given by the number of edges, providing the edge weights, or a
data.frame
with the number of rows being equal to the number of edges, where each row gives a vector of weights to its corresponding edge.kirchhoff_weights
If non-null, the name (or number) of the column of
weights
that contain the Kirchhoff weights. Must be equal to 1 (orTRUE
) in caseweights
is a single number and those are the Kirchhoff weights.
Returns
No return value. Called for its side effects.
Method get_edge_weights()
Gets the edge weights
Usage
metric_graph$get_edge_weights(data.frame = FALSE, tibble = TRUE)
Arguments
data.frame
If the edge weights are given as vectors, should the result be returned as a data.frame?
tibble
Should the edge weights be returned as tibble?
Returns
A vector or data.frame
containing the edge weights.
Method get_vertices_incomp_dir()
Gets vertices with incompatible directions
Usage
metric_graph$get_vertices_incomp_dir()
Returns
A vector containing the vertices with incompatible directions.
Method summary()
Prints a summary of various informations of the graph
Usage
metric_graph$summary( messages = FALSE, compute_characteristics = TRUE, check_euclidean = TRUE, check_distance_consistency = TRUE )
Arguments
messages
Should message explaining how to build the results be given for missing quantities?
compute_characteristics
Should the characteristics of the graph be computed?
check_euclidean
Check if the graph has Euclidean edges?
check_distance_consistency
Check the distance consistency assumption?
Returns
No return value. Called for its side effects.
Method print()
Prints various characteristics of the graph
Usage
metric_graph$print()
Returns
No return value. Called for its side effects.
Method compute_characteristics()
Computes various characteristics of the graph
Usage
metric_graph$compute_characteristics(check_euclidean = FALSE)
Arguments
check_euclidean
Also check if the graph has Euclidean edges? This essentially means that the distance consistency check will also be perfomed. If the graph does not have Euclidean edges due to another reason rather than the distance consistency, then it will already be indicated that the graph does not have Euclidean edges.
Returns
No return value. Called for its side effects. The computed characteristics
are stored in the characteristics
element of the metric_graph
object.
Method check_euclidean()
Check if the graph has Euclidean edges.
Usage
metric_graph$check_euclidean()
Returns
Returns TRUE
if the graph has Euclidean edges, or FALSE
otherwise.
The result is stored in the characteristics
element of the metric_graph
object.
The result is displayed when the graph is printed.
Method check_distance_consistency()
Checks distance consistency of the graph.
Usage
metric_graph$check_distance_consistency()
Returns
No return value.
The result is stored in the characteristics
element of the metric_graph
object.
The result is displayed when the graph is printed.
Method compute_geodist()
Computes shortest path distances between the vertices in the graph
Usage
metric_graph$compute_geodist(full = FALSE, obs = TRUE, group = NULL)
Arguments
full
Should the geodesic distances be computed for all the available locations? If
FALSE
, it will be computed separately for the locations of each group.obs
Should the geodesic distances be computed at the observation locations?
group
Vector or list containing which groups to compute the distance for. If
NULL
, it will be computed for all groups.
Returns
No return value. Called for its side effects. The computed geodesic
distances are stored in the geo_dist
element of the metric_graph
object.
Method compute_geodist_PtE()
Computes shortest path distances between the vertices in the graph.
Usage
metric_graph$compute_geodist_PtE( PtE, normalized = TRUE, include_vertices = TRUE )
Arguments
PtE
Points to compute the metric for.
normalized
are the locations in PtE in normalized distance?
include_vertices
Should the original vertices be included in the distance matrix?
Returns
A matrix containing the geodesic distances.
Method compute_geodist_mesh()
Computes shortest path distances between the vertices in the mesh.
Usage
metric_graph$compute_geodist_mesh()
Returns
No return value. Called for its side effects. The geodesic distances
on the mesh are stored in mesh$geo_dist
in the metric_graph
object.
Method compute_resdist()
Computes the resistance distance between the observation locations.
Usage
metric_graph$compute_resdist( full = FALSE, obs = TRUE, group = NULL, check_euclidean = FALSE, include_vertices = FALSE )
Arguments
full
Should the resistance distances be computed for all the available locations. If
FALSE
, it will be computed separately for the locations of each group.obs
Should the resistance distances be computed at the observation locations?
group
Vector or list containing which groups to compute the distance for. If
NULL
, it will be computed for all groups.check_euclidean
Check if the graph used to compute the resistance distance has Euclidean edges? The graph used to compute the resistance distance has the observation locations as vertices.
include_vertices
Should the vertices of the graph be also included in the resulting matrix when using
FULL=TRUE
?
Returns
No return value. Called for its side effects. The geodesic distances
are stored in the res_dist
element of the metric_graph
object.
Method compute_resdist_PtE()
Computes the resistance distance between the observation locations.
Usage
metric_graph$compute_resdist_PtE( PtE, normalized = TRUE, include_vertices = FALSE, check_euclidean = FALSE )
Arguments
PtE
Points to compute the metric for.
normalized
Are the locations in PtE in normalized distance?
include_vertices
Should the original vertices be included in the Laplacian matrix?
check_euclidean
Check if the graph used to compute the resistance distance has Euclidean edges? The graph used to compute the resistance distance has the observation locations as vertices.
Returns
A matrix containing the resistance distances.
Method get_degrees()
Returns the degrees of the vertices in the metric graph.
Usage
metric_graph$get_degrees(which = "degree")
Arguments
which
If "degree", returns the degree of the vertex. If "indegree", returns the indegree, and if "outdegree", it returns the outdegree.
Returns
A vector containing the degrees of the vertices.
Method compute_PtE_edges()
Computes the relative positions of the coordinates of the edges and save it as an attribute to each edge. This improves the quality of plots obtained by the plot_function()
method, however it might be costly to compute.
Usage
metric_graph$compute_PtE_edges()
Returns
No return value, called for its side effects.
Method compute_resdist_mesh()
Computes the resistance metric between the vertices in the mesh.
Usage
metric_graph$compute_resdist_mesh()
Returns
No return value. Called for its side effects. The geodesic distances
on the mesh are stored in the mesh$res_dist
element in the metric_graph
object.
Method compute_laplacian()
Computes the weigthed graph Laplacian for the graph.
Usage
metric_graph$compute_laplacian(full = FALSE, obs = TRUE, group = NULL)
Arguments
full
Should the resistance distances be computed for all the available locations. If
FALSE
, it will be computed separately for the locations of each group.obs
Should the resistance distances be computed at the observation locations? It will only compute for locations in which there is at least one observations that is not NA.
group
Vector or list containing which groups to compute the Laplacian for. If
NULL
, it will be computed for all groups.
Returns
No reutrn value. Called for its side effects. The Laplacian is stored
in the Laplacian
element in the metric_graph
object.
Method prune_vertices()
Removes vertices of degree 2 from the metric graph.
Usage
metric_graph$prune_vertices(check_weights = TRUE, verbose = FALSE)
Arguments
check_weights
If
TRUE
will only prune edges with different weights.verbose
Print progress of pruning. There are 3 levels of verbose, level 0, 1 and 2. In level 0, no messages are printed. In level 1, only messages regarding important steps are printed. Finally, in level 2, messages detailing all the steps are printed. The default is 1.
Details
Vertices of degree 2 are removed as long as the corresponding edges that would be merged are compatible in terms of direction.
Returns
No return value. Called for its side effects.
Method get_groups()
Gets the groups from the data.
Usage
metric_graph$get_groups(get_cols = FALSE)
Arguments
get_cols
Should the names of the columns that created the group variable be returned?
Returns
A vector containing the available groups in the internal data.
Method get_PtE()
Gets PtE from the data.
Usage
metric_graph$get_PtE()
Arguments
group
For which group, should the PtE be returned?
NULL
means that all PtEs available will be returned.include_group
Should the group be included as a column? If
TRUE
, the PtEs for each group will be concatenated, otherwise a single matrix containing the unique PtEs will be returned.
Returns
A matrix with two columns, where the first column contains the edge number and the second column contains the distance on edge of the observation locations.
Method get_edge_lengths()
Gets the edge lengths with the corresponding unit.
Usage
metric_graph$get_edge_lengths(unit = NULL)
Arguments
unit
If non-NULL, changes from
length_unit
from the graph construction tounit
.
Returns
a vector with the length unit (if the graph was constructed with a length unit).
Method get_locations()
Gets the spatial locations from the data.
Usage
metric_graph$get_locations()
Returns
A data.frame
object with observation locations. If longlat = TRUE
, the column names are lon and lat, otherwise the column names are x and y.
Method observation_to_vertex()
Adds observation locations as vertices in the graph.
Usage
metric_graph$observation_to_vertex(tolerance = 1e-15, mesh_warning = TRUE)
Arguments
tolerance
Observations locations are merged to a single vertex if they are closer than this number (given in relative edge distance between 0 and 1). The default is
1e-15
.mesh_warning
Display a warning if the graph structure change and the metric graph has a mesh object.
share_weights
Should the same weight be shared among the split edges? If
FALSE
, the weights will be removed, and a common weight given by 1 will be given.
Returns
No return value. Called for its side effects.
Method edgeweight_to_data()
Turns edge weights into data on the metric graph
Usage
metric_graph$edgeweight_to_data( loc = NULL, mesh = FALSE, data_loc = FALSE, weight_col = NULL, add = TRUE, data_coords = c("PtE", "spatial"), normalized = FALSE, tibble = TRUE, verbose = 1, suppress_warnings = FALSE, return = FALSE )
Arguments
loc
A
matrix
ordata.frame
with two columns containing the locations to generate the data from the edge weights. Ifdata_coords
is 'spatial', the first column must be the x-coordinate of the data, and the second column must be the y-coordinate. Ifdata_coords
is 'PtE', the first column must be the edge number and the second column must be the distance on edge.mesh
Should the data be generated to the mesh locations? In this case, the
loc
argument will be ignored. Observe that the metric graph must have a mesh built for one to use this option. CAUTION: To add edgeweight to data to both the data locations and mesh locations, please, add at the data locations first, then to mesh locations.data_loc
Should the data be generated to the data locations? In this case, the
loc
argument will be ignored. Observe that the metric graph must have data for one to use this option. CAUTION: To add edgeweight to data to both the data locations and mesh locations, please, add at the data locations first, then to mesh locations.weight_col
Which columns of the edge weights should be turned into data? If
NULL
, all columns will be turned into data.add
Should the data generated be added to the metric graph internal data?
data_coords
To be used only if
mesh
isFALSE
. It decides which coordinate system to use. IfPtE
, the user must provideedge_number
anddistance_on_edge
, otherwise ifspatial
, the user must providecoord_x
andcoord_y
.normalized
if TRUE, then the distances in
distance_on_edge
are assumed to be normalized to (0,1). Default FALSE.tibble
Should the data be returned as a
tidyr::tibble
?verbose
Print progress of the steps when adding observations. There are 3 levels of verbose, level 0, 1 and 2. In level 0, no messages are printed. In level 1, only messages regarding important steps are printed. Finally, in level 2, messages detailing all the steps are printed. The default is 1.
suppress_warnings
Suppress warnings related to duplicated observations?
return
Should the data be returned? If
return_removed
isTRUE
, only the removed locations will be return (if there is any).
Method get_mesh_locations()
Returns a list or a matrix with the mesh locations.
Usage
metric_graph$get_mesh_locations(bru = FALSE, loc = NULL, normalized = TRUE)
Arguments
bru
Should an 'inlabru'-friendly list be returned?
loc
If
bru
is set toTRUE
, the name of the location variable. The default name is 'loc'.normalized
If TRUE, then the distances in
distance_on_edge
are assumed to be normalized to (0,1). Default TRUE.
Returns
A list or a matrix containing the mesh locations.
Method clear_observations()
Clear all observations from the metric_graph
object.
Usage
metric_graph$clear_observations()
Returns
No return value. Called for its side effects.
Method process_data()
Process data to the metric graph data format.
Usage
metric_graph$process_data( data = NULL, edge_number = "edge_number", distance_on_edge = "distance_on_edge", coord_x = "coord_x", coord_y = "coord_y", data_coords = c("PtE", "spatial"), group = NULL, group_sep = ".", normalized = FALSE, tibble = TRUE, duplicated_strategy = "closest", include_distance_to_graph = TRUE, only_return_removed = FALSE, tolerance = max(self$edge_lengths)/2, verbose = FALSE, suppress_warnings = FALSE, Spoints = lifecycle::deprecated() )
Arguments
data
A
data.frame
or named list containing the observations. In case of groups, the data.frames for the groups should be stacked vertically, with a column indicating the index of the group. Ifdata
is notNULL
, it takes priority over any eventual data inSpoints
.edge_number
Column (or entry on the list) of the
data
that contains the edge numbers. If not supplied, the column with name "edge_number" will be chosen. Will not be used ifSpoints
is notNULL
.distance_on_edge
Column (or entry on the list) of the
data
that contains the edge numbers. If not supplied, the column with name "distance_on_edge" will be chosen. Will not be used ifSpoints
is notNULL
.coord_x
Column (or entry on the list) of the
data
that contains the x coordinate. If not supplied, the column with name "coord_x" will be chosen. Will not be used ifSpoints
is notNULL
or ifdata_coords
isPtE
.coord_y
Column (or entry on the list) of the
data
that contains the y coordinate. If not supplied, the column with name "coord_x" will be chosen. Will not be used ifSpoints
is notNULL
or ifdata_coords
isPtE
.data_coords
It decides which coordinate system to use. If
PtE
, the user must provideedge_number
anddistance_on_edge
, otherwise ifspatial
, the user must providecoord_x
andcoord_y
. The optioneuclidean
is . Usespatial
instead.group
Vector. If the data is grouped (for example measured at different time points), this argument specifies the columns (or entries on the list) in which the group variables are stored. It will be stored as a single column
.group
with the combined entries.group_sep
separator character for creating the new group variable when grouping two or more variables.
normalized
if TRUE, then the distances in
distance_on_edge
are assumed to be normalized to (0,1). Default FALSE.tibble
Should the data be returned as a
tidyr::tibble
?duplicated_strategy
Which strategy to handle observations on the same location on the metric graph (that is, if there are two or more observations projected at the same location). The options are 'closest' and 'jitter'. If 'closest', only the closest observation will be used. If 'jitter', a small perturbation will be performed on the projected observation location. The default is 'closest'.
include_distance_to_graph
When
data_coord
is 'spatial', should the distance of the observations to the graph be included as a column?only_return_removed
Should the removed data (if it exists) when using 'closest'
duplicated_strategy
be returned instead of the processed data?tolerance
Parameter to control a warning when adding observations. If the distance of some location and the closest point on the graph is greater than the tolerance, the function will display a warning. This helps detecting mistakes on the input locations when adding new data.
verbose
If
TRUE
, report steps and times.suppress_warnings
Suppress warnings related to duplicated observations?
Spoints
Returns
No return value. Called for its side effects. The observations are
stored in the data
element of the metric_graph
object.
Method add_observations()
Add observations to the metric graph.
Usage
metric_graph$add_observations( data = NULL, edge_number = "edge_number", distance_on_edge = "distance_on_edge", coord_x = "coord_x", coord_y = "coord_y", data_coords = c("PtE", "spatial"), group = NULL, group_sep = ".", normalized = FALSE, clear_obs = FALSE, tibble = FALSE, tolerance = max(self$edge_lengths)/2, duplicated_strategy = "closest", include_distance_to_graph = TRUE, return_removed = TRUE, verbose = 1, suppress_warnings = FALSE, Spoints = lifecycle::deprecated() )
Arguments
data
A
data.frame
or named list containing the observations. In case of groups, the data.frames for the groups should be stacked vertically, with a column indicating the index of the group.data
can also be ansf
object or aSpatialPointsDataFrame
object. in which casedata_coords
will automatically be spatial, and there is no need to specify thecoord_x
orcoord_y
arguments.edge_number
Column (or entry on the list) of the
data
that contains the edge numbers. If not supplied, the column with name "edge_number" will be chosen. Will not be used ifSpoints
is notNULL
.distance_on_edge
Column (or entry on the list) of the
data
that contains the edge numbers. If not supplied, the column with name "distance_on_edge" will be chosen. Will not be used ifSpoints
is notNULL
.coord_x
Column (or entry on the list) of the
data
that contains the x coordinate. If not supplied, the column with name "coord_x" will be chosen. Will not be used ifSpoints
is notNULL
or ifdata_coords
isPtE
.coord_y
Column (or entry on the list) of the
data
that contains the y coordinate. If not supplied, the column with name "coord_x" will be chosen. Will not be used ifSpoints
is notNULL
or ifdata_coords
isPtE
.data_coords
It decides which coordinate system to use. If
PtE
, the user must provideedge_number
anddistance_on_edge
, otherwise ifspatial
, the user must providecoord_x
andcoord_y
. The optioneuclidean
is . Usespatial
instead.group
Vector. If the data is grouped (for example measured at different time points), this argument specifies the columns (or entries on the list) in which the group variables are stored. It will be stored as a single column
.group
with the combined entries.group_sep
separator character for creating the new group variable when grouping two or more variables.
normalized
if TRUE, then the distances in
distance_on_edge
are assumed to be normalized to (0,1). Default FALSE.clear_obs
Should the existing observations be removed before adding the data?
tibble
Should the data be returned as a
tidyr::tibble
?tolerance
Parameter to control a warning when adding observations. If the distance of some location and the closest point on the graph is greater than the tolerance, the function will display a warning. This helps detecting mistakes on the input locations when adding new data.
duplicated_strategy
Which strategy to handle observations on the same location on the metric graph (that is, if there are two or more observations projected at the same location). The options are 'closest' and 'jitter'. If 'closest', only the closest observation will be used. If 'jitter', a small perturbation will be performed on the projected observation location. The default is 'closest'.
include_distance_to_graph
When
data_coord
is 'spatial', should the distance of the observations to the graph be included as a column?return_removed
Should the removed data (if it exists) when using 'closest'
duplicated_strategy
be returned?verbose
Print progress of the steps when adding observations. There are 3 levels of verbose, level 0, 1 and 2. In level 0, no messages are printed. In level 1, only messages regarding important steps are printed. Finally, in level 2, messages detailing all the steps are printed. The default is 1.
suppress_warnings
Suppress warnings related to duplicated observations?
Spoints
Returns
No return value. Called for its side effects. The observations are
stored in the data
element of the metric_graph
object.
Method mutate()
Use dplyr::mutate
function on the internal metric graph data object.
Usage
metric_graph$mutate(..., .drop_na = FALSE, .drop_all_na = TRUE)
Arguments
...
Arguments to be passed to
dplyr::mutate()
..drop_na
Should the rows with at least one NA for one of the columns be removed? DEFAULT is
FALSE
..drop_all_na
Should the rows with all variables being NA be removed? DEFAULT is
TRUE
.
Details
A wrapper to use dplyr::mutate()
within the internal metric graph data object.
Returns
A tidyr::tibble
object containing the resulting data list after the mutate.
Method drop_na()
Use tidyr::drop_na()
function on the internal metric graph data object.
Usage
metric_graph$drop_na(...)
Arguments
...
Arguments to be passed to
tidyr::drop_na()
.
Details
A wrapper to use dplyr::drop_na()
within the internal metric graph data object.
Returns
A tidyr::tibble
object containing the resulting data list after the drop_na.
Method select()
Use dplyr::select
function on the internal metric graph data object.
Usage
metric_graph$select(..., .drop_na = FALSE, .drop_all_na = TRUE)
Arguments
...
Arguments to be passed to
dplyr::select()
..drop_na
Should the rows with at least one NA for one of the columns be removed? DEFAULT is
FALSE
..drop_all_na
Should the rows with all variables being NA be removed? DEFAULT is
TRUE
.
Details
A wrapper to use dplyr::select()
within the internal metric graph data object. Observe that it is a bit different from directly using dplyr::select()
since it does not allow to remove the internal positions that are needed for the metric_graph methods to work.
Returns
A tidyr::tibble
object containing the resulting data list after the selection.
Method filter()
Use dplyr::filter
function on the internal metric graph data object.
Usage
metric_graph$filter(..., .drop_na = FALSE, .drop_all_na = TRUE)
Arguments
...
Arguments to be passed to
dplyr::filter()
..drop_na
Should the rows with at least one NA for one of the columns be removed? DEFAULT is
FALSE
..drop_all_na
Should the rows with all variables being NA be removed? DEFAULT is
TRUE
.
Details
A wrapper to use dplyr::filter()
within the internal metric graph data object.
Returns
A tidyr::tibble
object containing the resulting data list after the filter.
Method summarise()
Use dplyr::summarise
function on the internal metric graph data object grouped by the spatial locations and the internal group variable.
Usage
metric_graph$summarise( ..., .include_graph_groups = FALSE, .groups = NULL, .drop_na = FALSE, .drop_all_na = TRUE )
Arguments
...
Arguments to be passed to
dplyr::summarise()
..include_graph_groups
Should the internal graph groups be included in the grouping variables? The default is
FALSE
. This means that, when summarising, the data will be grouped by the internal group variable together with the spatial locations..groups
A vector of strings containing the names of the columns to be additionally grouped, when computing the summaries. The default is
NULL
..drop_na
Should the rows with at least one NA for one of the columns be removed? DEFAULT is
FALSE
..drop_all_na
Should the rows with all variables being NA be removed? DEFAULT is
TRUE
.
Details
A wrapper to use dplyr::summarise()
within the internal metric graph data object grouped by manually inserted groups (optional), the internal group variable (optional) and the spatial locations. Observe that if the integral group variable was not used as a grouping variable for the summarise, a new column, called .group
, will be added, with the same value 1 for all rows.
Returns
A tidyr::tibble
object containing the resulting data list after the summarise.
Method get_data()
Return the internal data with the option to filter by groups.
Usage
metric_graph$get_data( group = NULL, tibble = TRUE, drop_na = FALSE, drop_all_na = TRUE )
Arguments
group
A vector contaning which groups should be returned? The default is
NULL
, which gives the result for the all groups.tibble
Should the data be returned as a
tidyr::tibble
?drop_na
Should the rows with at least one NA for one of the columns be removed? DEFAULT is
FALSE
.drop_all_na
Should the rows with all variables being NA be removed? DEFAULT is
TRUE
.
Method buildDirectionalConstraints()
Build directional ODE constraint matrix from edges.
Usage
metric_graph$buildDirectionalConstraints(alpha = 1)
Arguments
alpha
how many derivatives the processes has
Details
Currently not implemented for circles (edges that start and end in the same vertex)
Returns
No return value. Called for its side effects.
Method buildC()
Build Kirchoff constraint matrix from edges.
Usage
metric_graph$buildC(alpha = 2, edge_constraint = FALSE)
Arguments
alpha
the type of constraint (currently only supports 2)
edge_constraint
if TRUE, add constraints on vertices of degree 1
Details
Currently not implemented for circles (edges that start and end in the same vertex)
Returns
No return value. Called for its side effects.
Method build_mesh()
Builds mesh object for graph.
Usage
metric_graph$build_mesh( h = NULL, n = NULL, continuous = TRUE, continuous.outs = FALSE, continuous.deg2 = FALSE )
Arguments
h
Maximum distance between mesh nodes (should be provided if n is not provided).
n
Maximum number of nodes per edge (should be provided if h is not provided).
continuous
If
TRUE
(default), the mesh contains only one node per vertex. IfFALSE
, each vertex v is split into deg(v) disconnected nodes to allow for the creation of discontinuities at the vertices.continuous.outs
If
continuous = FALSE
andcontinuous.outs = TRUE
, continuity is assumed for the outgoing edges from each vertex.continuous.deg2
If
TRUE
, continuity is assumed at degree 2 vertices.
Details
The mesh is a list with the objects:
-
PtE
The mesh locations excluding the original vertices; -
V
The verties of the mesh; -
E
The edges of the mesh; -
n_e
The number of vertices in the mesh per original edge in the graph; -
h_e
The mesh width per edge in the graph; -
ind
The indices of the vertices in the mesh; -
VtE
All mesh locations including the original vertices.
Returns
No return value. Called for its side effects. The mesh is stored in
the mesh
element of the metric_graph
object.
Method compute_fem()
Build mass and stiffness matrices for given mesh object.
Usage
metric_graph$compute_fem(petrov = FALSE)
Arguments
petrov
Compute Petrov-Galerkin matrices? (default
FALSE
). These are defined asCpet_{ij} = <\phi_i, \psi_j>
andGpet_{ij} = <d\phi_i, \psi_j>
, where\psi_{i}
are piecewise constant basis functions on the edges of the mesh.
Details
The function builds: The matrix C
which is the mass matrix with
elements C_{ij} = <\phi_i, \phi_j>
, the matrix G
which is the stiffness
matrix with elements G_{ij} = <d\phi_i, d\phi_j>
, the matrix B
with
elements B_{ij} = <d\phi_i, \phi_j>
, the matrix D
with elements
D_{ij} = \sum_{v\in V}\phi_i(v)\phi_j(v)
, and the vector with weights
<\phi_i, 1>
.
Returns
No return value. Called for its side effects. The finite element
matrices C
, G
and B
are stored in the mesh
element in the
metric_graph
object. If petrov=TRUE
, the corresponding Petrov-Galerkin
matrices are stored in Cpet
and Gpet
.
Method mesh_A()
Deprecated - Computes observation matrix for mesh.
in favour of metric_graph$fem_basis()
.
Usage
metric_graph$mesh_A(PtE)
Arguments
PtE
Locations given as (edge number in graph, normalized location on edge)
Details
For n locations and a mesh with m nodes, A
is an n x m matrix with
elements A_{ij} = \phi_j(s_i)
.
Returns
The observation matrix.
Method fem_basis()
Computes observation matrix for mesh.
Usage
metric_graph$fem_basis(PtE)
Arguments
PtE
Locations given as (edge number in graph, normalized location on edge)
Details
For n locations and a mesh with m nodes, A
is an n x m matrix with
elements A_{ij} = \phi_j(s_i)
.
Returns
The observation matrix.
Method VtEfirst()
Find one edge corresponding to each vertex.
Usage
metric_graph$VtEfirst()
Returns
A nV x 2 matrix the first element of the i
th row is the edge
number corresponding to the i
th vertex and the second value is 0
if the vertex is at the start of the edge and 1 if the vertex
is at the end of the edge.
Method plot()
Plots the metric graph.
Usage
metric_graph$plot( data = NULL, newdata = NULL, group = 1, plotly = FALSE, interactive = FALSE, vertex_size = 3, vertex_color = "black", edge_width = 0.3, edge_color = "black", data_size = 1, support_width = 0.5, support_color = "gray", mesh = FALSE, X = NULL, X_loc = NULL, p = NULL, degree = FALSE, direction = FALSE, edge_weight = NULL, edge_width_weight = NULL, scale_color_main = ggplot2::scale_color_viridis_c(option = "D"), scale_color_weights = ggplot2::scale_color_viridis_c(option = "C"), scale_color_degree = ggplot2::scale_color_viridis_d(option = "D"), add_new_scale_weights = TRUE, ... )
Arguments
data
Which column of the data to plot? If
NULL
, no data will be plotted.newdata
A dataset of class
metric_graph_data
, obtained by anyget_data()
,mutate()
,filter()
,summarise()
,drop_na()
methods of metric graphs, see the vignette on data manipulation for more details.group
If there are groups, which group to plot? If
group
is a number, it will be the index of the group as stored internally. Ifgroup
is a character, then the group will be chosen by its name.plotly
Use plot_ly for 3D plot (default
FALSE
). This option requires the 'plotly' package.interactive
Only works for 2d plots. If
TRUE
, an interactive plot will be displayed. Unfortunately,interactive
is not compatible withedge_weight
ifadd_new_scale_weights
is TRUE.vertex_size
Size of the vertices.
vertex_color
Color of vertices.
edge_width
Line width for edges. If
edge_width_weight
is notNULL
, this determines the maximum edge width.edge_color
Color of edges.
data_size
Size of markers for data.
support_width
For 3D plot, width of support lines.
support_color
For 3D plot, color of support lines.
mesh
Plot the mesh locations?
X
Additional values to plot.
X_loc
Locations of the additional values in the format (edge, normalized distance on edge).
p
Existing objects obtained from 'ggplot2' or 'plotly' to add the graph to
degree
Show the degrees of the vertices?
direction
Show the direction of the edges?
edge_weight
Which column from edge weights to plot? If
NULL
edge weights are not plotted. To plot the edge weights when the metric graphedge_weights
is a vector instead of adata.frame
, simply set to 1.edge_weight
is only available for 2d plots. For 3d plots with edge weights, please use theplot_function()
method.edge_width_weight
Which column from edge weights to determine the edges widths? If
NULL
edge width will be determined fromedge_width
.scale_color_main
Color scale for the data to be plotted.
scale_color_weights
Color scale for the edge weights. Will only be used if
add_new_scale_weights
is TRUE.scale_color_degree
Color scale for the degrees.
add_new_scale_weights
Should a new color scale for the edge weights be created?
...
Additional arguments to pass to
ggplot()
orplot_ly()
Returns
A plot_ly
(if plotly = TRUE
) or ggplot
object.
Method plot_connections()
Plots the connections in the graph
Usage
metric_graph$plot_connections()
Returns
No return value. Called for its side effects.
Method is_tree()
Checks if the graph is a tree (without considering directions)
Usage
metric_graph$is_tree()
Returns
TRUE if the graph is a tree and FALSE otherwise.
Method plot_function()
Plots continuous function on the graph.
Usage
metric_graph$plot_function( data = NULL, newdata = NULL, group = 1, X = NULL, plotly = FALSE, improve_plot = FALSE, continuous = TRUE, edge_weight = NULL, vertex_size = 5, vertex_color = "black", edge_width = 1, edge_color = "black", line_width = NULL, line_color = "rgb(0,0,200)", scale_color = ggplot2::scale_color_viridis_c(option = "d"), support_width = 0.5, support_color = "gray", p = NULL, ... )
Arguments
data
Which column of the data to plot? If
NULL
, no data will be plotted.newdata
A dataset of class
metric_graph_data
, obtained by anyget_data()
,mutate()
,filter()
,summarise()
,drop_na()
methods of metric graphs, see the vignette on data manipulation for more details.group
If there are groups, which group to plot? If
group
is a number, it will be the index of the group as stored internally. Ifgroup
is a character, then the group will be chosen by its name.X
A vector with values for the function evaluated at the mesh in the graph
plotly
If
TRUE
, then the plot is shown in 3D. This option requires the package 'plotly'.improve_plot
Should the original edge coordinates be added to the data with linearly interpolated values to improve the plot?
continuous
Should continuity be assumed when the plot uses
newdata
?edge_weight
Which column from edge weights to plot? If
NULL
edge weights are not plotted. To plot the edge weights when the metric graphedge_weights
is a vector instead of adata.frame
, simply set to 1.vertex_size
Size of the vertices.
vertex_color
Color of vertices.
edge_width
Width for edges.
edge_color
For 3D plot, color of edges.
line_width
For 3D plot, line width of the function curve.
line_color
Color of the function curve.
scale_color
Color scale to be used for data and weights.
support_width
For 3D plot, width of support lines.
support_color
For 3D plot, color of support lines.
p
Previous plot to which the new plot should be added.
...
Additional arguments for
ggplot()
orplot_ly()
Returns
Either a ggplot
(if plotly = FALSE
) or a plot_ly
object.
Method plot_movie()
Plots a movie of a continuous function evolving on the graph.
Usage
metric_graph$plot_movie( X, plotly = TRUE, vertex_size = 5, vertex_color = "black", edge_width = 1, edge_color = "black", line_width = NULL, line_color = "rgb(0,0,200)", ... )
Arguments
X
A m x T matrix where the ith column represents the function at the ith time, evaluated at the mesh locations.
plotly
If
TRUE
, then plot is shown in 3D. This option requires the package 'plotly'.vertex_size
Size of the vertices.
vertex_color
Color of vertices.
edge_width
Width for edges.
edge_color
For 3D plot, color of edges.
line_width
For 3D plot, line width of the function curve.
line_color
Color of the function curve.
...
Additional arguments for ggplot or plot_ly.
Returns
Either a ggplot
(if plotly=FALSE
) or a plot_ly
object.
Method add_mesh_observations()
Add observations on mesh to the object.
Usage
metric_graph$add_mesh_observations(data = NULL, group = NULL)
Arguments
data
A
data.frame
or named list containing the observations. In case of groups, the data.frames for the groups should be stacked vertically, with a column indicating the index of the group. Ifdata_frame
is notNULL
, it takes priority over any eventual data inSpoints
.group
If the data_frame contains groups, one must provide the column in which the group indices are stored.
Returns
No return value. Called for its side effects. The observations are
stored in the data
element in the metric_graph
object.
Method get_initial_graph()
Returns a copy of the initial metric graph.
Usage
metric_graph$get_initial_graph()
Returns
A metric_graph
object.
Method coordinates()
Convert between locations on the graph and Euclidean coordinates.
Usage
metric_graph$coordinates(PtE = NULL, XY = NULL, normalized = TRUE)
Arguments
PtE
Matrix with locations on the graph (edge number and normalized position on the edge).
XY
Matrix with locations in Euclidean space
normalized
If
TRUE
, it is assumed that the positions inPtE
are normalized to (0,1), and the object returned ifXY
is specified contains normalized locations.
Returns
If PtE
is specified, then a matrix with Euclidean coordinates of
the locations is returned. If XY
is provided, then a matrix with the
closest locations on the graph is returned.
Gets the edge weights
data.frame If the edge weights are given as vectors, should the result be returned as a data.frame?
A vector or data.frame
containing the edge weights.
data List containing data on the metric graph.
Method clone()
The objects of this class are cloneable with this method.
Usage
metric_graph$clone(deep = FALSE)
Arguments
deep
Whether to make a deep clone.
Examples
edge1 <- rbind(c(0, 0), c(2, 0))
edge2 <- rbind(c(2, 0), c(1, 1))
edge3 <- rbind(c(1, 1), c(0, 0))
edges <- list(edge1, edge2, edge3)
graph <- metric_graph$new(edges)
graph$plot()