Skip to content

Commit 085dee7

Browse files
committed
Moves dist calculation to fortran code
1 parent c517a49 commit 085dee7

File tree

4 files changed

+197
-13
lines changed

4 files changed

+197
-13
lines changed

R/twdtw_reduce_time.R

+53-13
Original file line numberDiff line numberDiff line change
@@ -82,21 +82,28 @@ twdtwReduceTime = function(x,
8282
py <- py[,names(px),drop=FALSE]
8383

8484
# Compute local cost matrix
85-
cm <- proxy::dist(py, px, method = dist.method)
85+
#cm <- proxy::dist(py, px, method = dist.method)
8686

87-
if(!is.null(weight.fun)){
88-
# Get day of the year for pattern and time series
89-
doyy <- as.numeric(format(ty, "%j"))
90-
doyx <- as.numeric(format(tx, "%j"))
91-
92-
# Compute time-weght matrix
93-
w <- .g(proxy::dist(doyy, doyx, method = dist.method))
94-
95-
# Apply time-weight to local cost matrix
96-
cm <- weight.fun(cm, w)
97-
}
87+
# Get day of the year for pattern and time series
88+
doyy <- as.numeric(format(ty, "%j"))
89+
doyx <- as.numeric(format(tx, "%j"))
90+
91+
# if(!is.null(weight.fun)){
92+
# # Get day of the year for pattern and time series
93+
# doyy <- as.numeric(format(ty, "%j"))
94+
# doyx <- as.numeric(format(tx, "%j"))
95+
#
96+
# # Compute time-weght matrix
97+
# w <- .g(proxy::dist(doyy, doyx, method = dist.method))
98+
#
99+
# # Apply time-weight to local cost matrix
100+
# cm <- weight.fun(cm, w)
101+
# }
98102
# Compute accumulated DTW cost matrix
99-
internals <- .computecost(cm = cm, step.matrix = step.matrix)
103+
#internals <- dtwSat:::.computecost(cm = cm, step.matrix = step.matrix)
104+
xm = na.omit(cbind(doyx, as.matrix(px)))
105+
ym = na.omit(cbind(doyy, as.matrix(py)))
106+
internals = .computecost_fast(xm, ym, step.matrix)
100107

101108
# Find all low cost candidates
102109
a <- internals$startingMatrix[internals$N,1:internals$M]
@@ -147,3 +154,36 @@ twdtwReduceTime = function(x,
147154
return(out)
148155
}
149156

157+
# @useDynLib dtwSat computecost
158+
.computecost_fast = function(xm, ym, step.matrix){
159+
160+
# cm = rbind(0, cm)
161+
n = nrow(ym)
162+
m = nrow(xm)
163+
d = ncol(ym)
164+
165+
if(is.loaded("computecostfast", PACKAGE = "dtwSat", type = "Fortran")){
166+
out = .Fortran(computecostfast,
167+
XM = matrix(as.double(xm), m, d),
168+
YM = matrix(as.double(ym), n, d),
169+
CM = matrix(as.double(0), n+1, m),
170+
DM = matrix(as.integer(0), n+1, m),
171+
VM = matrix(as.integer(0), n+1, m),
172+
SM = matrix(as.integer(step.matrix), nrow(step.matrix), ncol(step.matrix)),
173+
N = as.integer(n),
174+
M = as.integer(m),
175+
D = as.integer(d),
176+
NS = as.integer(nrow(step.matrix)))
177+
} else {
178+
stop("Fortran computecostfast lib is not loaded")
179+
}
180+
# sqrt(sum((ym[1,-1] - xm[1,-1])*(ym[1,-1] - xm[1,-1])))
181+
res = list()
182+
res$costMatrix = out$CM[-1,]
183+
res$directionMatrix = out$DM[-1,]
184+
res$startingMatrix = out$VM[-1,]
185+
res$stepPattern = step.matrix
186+
res$N = n
187+
res$M = m
188+
res
189+
}

src/computecost_fast.f

+96
Original file line numberDiff line numberDiff line change
@@ -0,0 +1,96 @@
1+
CCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCC
2+
C C
3+
C (c) Victor Maus <[email protected]> C
4+
C Institute for Geoinformatics (IFGI) C
5+
C University of Muenster (WWU), Germany C
6+
C C
7+
C Earth System Science Center (CCST) C
8+
C National Institute for Space Research (INPE), Brazil C
9+
C C
10+
C C
11+
C Efficient computation of DTW cost matrix - 2015-10-16 C
12+
C C
13+
C This function was adpted from the C function 'computeCM' C
14+
C implemented in the R package 'dtw' by Toni Giorgino. C
15+
CCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCC
16+
C
17+
C XM - matrix with the time series (N,D)
18+
C YM - matrix with the temporal profile (M,D)
19+
C CM - Output cumulative cost matrix
20+
C DM - Direction matrix
21+
C VM - Starting points matrix
22+
C SM - Matrix of step patterns
23+
C N - Number of rows in CM, DM, and VM - time series
24+
C M - Number of columns CM, DM, and VM - temporal profile
25+
C D - Number of spectral dimensions including time in XM and YM
26+
C NS - Number of rows in SM
27+
SUBROUTINE computecostfast(XM, YM, CM, DM, VM, SM, N, M, D, NS)
28+
C I/O Variables
29+
INTEGER N, M, D, NS, SM(NS,4), DM(N+1,M), VM(N+1,M)
30+
DOUBLE PRECISION XM(M,D), YM(N,D), CM(N+1,M)
31+
C Internals
32+
DOUBLE PRECISION W, CP(NS), VMIN
33+
INTEGER I, J, IL(NS), JL(NS), K, PK, KMIN, ZERO, ONE
34+
PARAMETER(ZERO=0,ONE=1)
35+
REAL NAN, INF
36+
NAN = ZERO
37+
NAN = NAN / NAN
38+
INF = HUGE(ZERO)
39+
VM(1,1) = 1
40+
C Initialize the firt row and col of the matrices
41+
DO 21 I = 2, N+1
42+
CM(I,1) = CM(I-1,1) + distance(YM, XM, N, M, D, I-1, 1)
43+
C WRITE (*,*) 'The distance ',I,',1 is ', CM(I,1)
44+
DM(I,1) = 3
45+
VM(I,1) = 1
46+
21 CONTINUE
47+
DO 31 J = 2, M
48+
CM(2,J) = CM(2,J-1) + distance(YM, XM, N, M, D, 1, J)
49+
C WRITE (*,*) 'The distance 2,',J,' is ', CM(2,J)
50+
DM(1,J) = 2
51+
VM(1,J) = J
52+
31 CONTINUE
53+
C Compute cumulative cost matrix
54+
DO 32 J = 2, M
55+
DO 22 I = 2, N+1
56+
C Calculate local distance
57+
CM(I,J) = distance(YM, XM, N, M, D, I-1, J)
58+
C Initialize list of step cost
59+
DO 10 K = 1, NS
60+
CP(K) = NAN
61+
10 CONTINUE
62+
DO 11 K = 1, NS
63+
PK = SM(K,1)
64+
IL(K) = I - SM(K,2)
65+
JL(K) = J - SM(K,3)
66+
IF ((IL(K).GT.ZERO).AND.(JL(K).GT.ZERO)) THEN
67+
W = SM(K,4)
68+
IF (W.EQ.-ONE) THEN
69+
CP(PK) = CM(IL(K),JL(K))
70+
ELSE
71+
CP(PK) = CP(PK) + CM(IL(K),JL(K))*W
72+
ENDIF
73+
ENDIF
74+
11 CONTINUE
75+
KMIN = -ONE
76+
VMIN = INF
77+
ILMIN = -ONE
78+
JLMIN = -ONE
79+
DO 12 K = 1, NS
80+
PK = SM(K,1)
81+
IF (CP(PK).EQ.CP(PK).AND.CP(PK).LT.VMIN) THEN
82+
KMIN = PK
83+
VMIN = CP(PK)
84+
ILMIN = IL(K)
85+
JLMIN = JL(K)
86+
ENDIF
87+
12 CONTINUE
88+
IF (KMIN.GT.-ONE) THEN
89+
CM(I,J) = VMIN
90+
DM(I,J) = KMIN
91+
VM(I,J) = VM(ILMIN, JLMIN)
92+
ENDIF
93+
22 CONTINUE
94+
32 CONTINUE
95+
99 CONTINUE
96+
END

src/distance.f

+46
Original file line numberDiff line numberDiff line change
@@ -0,0 +1,46 @@
1+
CCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCC
2+
C C
3+
C (c) Victor Maus <[email protected]> C
4+
C Institute for Geoinformatics (IFGI) C
5+
C University of Muenster (WWU), Germany C
6+
C C
7+
C Earth System Science Center (CCST) C
8+
C National Institute for Space Research (INPE), Brazil C
9+
C C
10+
C C
11+
C Computation allapsed time - 2016-01-22 C
12+
C C
13+
CCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCCC
14+
C
15+
C XM - matrix with the time series (N,D)
16+
C YM - matrix with the temporal profile (M,D)
17+
C N - Number of rows in CM, DM, and VM - time series
18+
C M - Number of columns CM, DM, and VM - temporal profile
19+
C D - Number of spectral dimensions including time in XM and YM
20+
C I - Single point in the time series to calculate the local distance
21+
C J - Single point in the temporal profile to calculate the local distance
22+
REAL FUNCTION distance(YM, XM, N, M, D, I, J)
23+
INTEGER N, M, D, I, J, K
24+
DOUBLE PRECISION XM(M,D), YM(N,D), PC, TD, BD, CD, HPC
25+
PARAMETER(ZERO=0,ONE=1)
26+
REAL NAN, INF
27+
NAN = ZERO
28+
NAN = NAN / NAN
29+
PC = 366
30+
HPC = PC/2
31+
distance = NAN
32+
C Compute ellapsed time difference
33+
TD = SQRT((YM(I,1) - XM(J,1)) * (YM(I,1) - XM(J,1)))
34+
IF (TD.GT.HPC) THEN
35+
TD = ABS(PC - TD)
36+
ENDIF
37+
CD = ZERO
38+
DO 30 K = 2, D
39+
BD = YM(I,K) - XM(J,K)
40+
CD = CD + (BD * BD)
41+
30 CONTINUE
42+
C WRITE (*,*) 'The value of X J', J, 'is', XM(J,D)
43+
distance = SQRT(CD) + 1.0 / (1.0 + EXP(-0.1 * (TD - 50.0)))
44+
RETURN
45+
END
46+

src/init.c

+2
Original file line numberDiff line numberDiff line change
@@ -10,12 +10,14 @@
1010
/* .Fortran calls */
1111
extern void F77_NAME(bestmatches)(void *, void *, void *, void *, void *, void *, void *, void *, void *, void *, void *);
1212
extern void F77_NAME(computecost)(void *, void *, void *, void *, void *, void *, void *);
13+
extern void F77_NAME(computecostfast)(void *, void *, void *, void *, void *, void *, void *, void *, void *, void *);
1314
extern void F77_NAME(g)(void *, void *, void *, void *);
1415
extern void F77_NAME(tracepath)(void *, void *, void *, void *, void *, void *, void *, void *, void *, void *, void *);
1516

1617
static const R_FortranMethodDef FortranEntries[] = {
1718
{"bestmatches", (DL_FUNC) &F77_NAME(bestmatches), 11},
1819
{"computecost", (DL_FUNC) &F77_NAME(computecost), 7},
20+
{"computecostfast", (DL_FUNC) &F77_NAME(computecostfast), 10},
1921
{"g", (DL_FUNC) &F77_NAME(g), 4},
2022
{"tracepath", (DL_FUNC) &F77_NAME(tracepath), 11},
2123
{NULL, NULL, 0}

0 commit comments

Comments
 (0)