File: simLocal.spatialProcess.R

package info (click to toggle)
r-cran-fields 16.3.1-1
  • links: PTS, VCS
  • area: main
  • in suites: forky, sid, trixie
  • size: 4,972 kB
  • sloc: fortran: 1,021; ansic: 288; sh: 35; makefile: 2
file content (232 lines) | stat: -rw-r--r-- 7,905 bytes parent folder | download
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
#
# fields  is a package for analysis of spatial data written for
# the R software environment.
# Copyright (C) 2024 Colorado School of Mines
# 1500 Illinois St., Golden, CO 80401
# Contact: Douglas Nychka,  douglasnychka@gmail.com,
#
# This program is free software; you can redistribute it and/or modify
# it under the terms of the GNU General Public License as published by
# the Free Software Foundation; either version 2 of the License, or
# (at your option) any later version.
# This program is distributed in the hope that it will be useful,
# but WITHOUT ANY WARRANTY; without even the implied warranty of
# MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
# GNU General Public License for more details.
#
# You should have received a copy of the GNU General Public License
# along with the R software environment if not, write to the Free Software
# Foundation, Inc., 51 Franklin St, Fifth Floor, Boston, MA  02110-1301  USA
# or see http://www.r-project.org/Licenses/GPL-2
##END HEADER
"simLocal.spatialProcess" <- function(mKrigObject,  
    predictionGridList = NULL,
    simulationGridList = NULL, 
        gridRefinement = 1, 
                    np = 2, 
                     M = 1,
                    nx = 80,
                    ny = 80, 
               verbose = FALSE,
                 delta = NULL, giveWarnings=TRUE,
                  fast = FALSE, 
                NNSize = 5, 
                          ...)
#
# NOTE throughout $x is first dimension of the grid in a gridList but also 
# $x in mKrig object is the _matrix_ of locations
# 
    {
  
   nObs<- nrow( mKrigObject$x)
   sDimension<- ncol(mKrigObject$x)
   if ( sDimension > 2) {
        stop("conditional simulation only implemented for 1 and 2 dimensions")
   }
   
   if( sDimension == 1 & fast ){
     stop("fast prediction not implemented in 1 D")
   }
   
    # create prediction set of points based on what is passed
    # and if the grid is not specified 
   if (is.null(predictionGridList)) {
     # these adjustments insure there are enough grid
     # points beyond the range of the locations.
     # Put xr[1] in the middle of the  npth grid box
     # and xr[2] in to the  nx - np
     predictionGridList<- makePredictionGridList(
       mKrigObject=mKrigObject, 
       nx=nx, 
       ny=ny, 
       np=np
     )
   }
   
   nx <- length(predictionGridList$x)
   ny <- ifelse(sDimension>=2, length(predictionGridList$y) ,1 )
   
# 
# check that predictionGrid is equally spaced
# this is needed because of the fast simulation algorithm
    checkPredictGrid( predictionGridList) 
#
#
   if (is.null(simulationGridList)) {
     simulationGridList<- makeSimulationGrid(
       predictionGridList,
       gridRefinement)
   }
#
# #
# # round off the grids so that they match to 8 digits
# # that way prediction grid is  precisely a subset of
# # simulation grid
         predictionGridList$x<- signif(predictionGridList$x, 8)
         simulationGridList$x<- signif(simulationGridList$x, 8)
         if( sDimension ==2){
         predictionGridList$y<- signif(predictionGridList$y, 8)
         simulationGridList$y<- signif(simulationGridList$y, 8)
         }
         indexSubset<-  list( x=match(predictionGridList$x,
                                      simulationGridList$x))
         # # shortcut to avoid if statement for predicted in for
         # # loop
         # indexSubset$y = rep(1, length( indexSubset$x) )
         if( sDimension ==2){
         indexSubset$y=match(predictionGridList$y,
                                      simulationGridList$y)
         }
         else{
           indexSubset$y=1
         }
         
#  
# core covariance parameters from spatial model   
    tau <-    mKrigObject$summary["tau"]
    sigma2 <- mKrigObject$summary["sigma2"]
    aRange<-  mKrigObject$summary["aRange"]
    Covariance <- mKrigObject$args$Covariance
# wipe out some extraneous components that are not used by the Covariance
# function.
    covArgs0 <- mKrigObject$args
    covArgs0$Covariance<- NULL 
    covArgs0$distMat <- NULL
    covArgs0$onlyUpper<- NULL
    covArgs0$aRange<- NULL
    
#
# set up various arrays for reuse during the simulation
    nObs <- nrow(mKrigObject$x)
#  
    
    timeCESetup<- system.time(
    # set up object for simulating on a grid using circulant embedding
    CEObject<- circulantEmbeddingSetup(simulationGridList,
                                   cov.function = mKrigObject$cov.function,
                                       cov.args = mKrigObject$args,
                                          delta = delta )
    )[3]
    
    #
    if (verbose) {
        cat("dim of full circulant matrix ", CEObject$M, 
            fill = TRUE)
    }
#
# weights crucial to fast off grid simulation
#
      timeOffGridSetup <- system.time(
        offGridObject <- offGridWeights(
          mKrigObject$x,
          simulationGridList,
          mKrigObject,
          np = np,
          giveWarnings = giveWarnings
        )
      )[3]
   
    #
    # find conditional mean field from initial fit
      hHat <- predictSurface(mKrigObject,
                           gridList = predictionGridList,
                           fast=fast, 
                           NNSize= NNSize,
                            ...)$z
      
      sdNugget<- tau* sqrt(1/mKrigObject$weights)
#      
# setup output array to hold ensemble
#  in 1D case ny=1 
#    
    out <- array(NA, c( nx, ny, M))
    t1<-t2<- t3<- rep( NA, M)
    
##########################################################################################
### begin the big loop
##########################################################################################
    
     for (k in 1:M) {
       if( k%%10 ==0 ){
        cat(k, " ")
       }
        # simulate full field
        t1[k]<- system.time(
        hTrue<- as.matrix(sqrt(sigma2) * circulantEmbedding(CEObject))
        )[3]
        
        # NOTE: fixed part of model (null space) does not need to be simulated
        # because the  estimator is unbiased for this part.
        # the variability is still captured because the fixed part
        # is still estimated as part of the predict step below
        #
        t2[k]<- system.time(
        hData <- offGridObject$B%*%c(hTrue) + 
                    (offGridObject$SE)%*%rnorm(nObs) 
              )[3]
        ySynthetic <- hData + sdNugget*rnorm(nObs)
        #
        # predict at grid using these data
        # and subtract from synthetic 'true' value
        #
        t3[k]<-system.time(
            spatialError <- predictSurface.mKrig(mKrigObject,
                                 gridList = predictionGridList,
                                 ynew = ySynthetic,
                                 fast=fast, 
                                 NNSize= NNSize,
                                 giveWarnings = FALSE,
                                 ...)$z
         
             )[3]
 
        
        # add the error to the actual estimate  (conditional mean)
        # subset  hTrue to the prediction grid
        # note for 1D $y is 1. 
       
        out[,, k] <- hHat + (spatialError -  
                               hTrue[indexSubset$x,indexSubset$y])
       
     }
    
    cat(" ", fill=TRUE)
    
    return(list(x = predictionGridList$x,
                y = predictionGridList$y,
                z = out, 
                hHat= hHat,
                timing=c( CESetup=timeCESetup,
                          OffSetup=timeOffGridSetup,
                                  CE = median(t1), 
                             OffGrid = median(t2),
                                mKrig = median(t3)
                          ),
                gridRefinement=gridRefinement,
                M= CEObject$M,
                #simulationGridList= simulationGridList,
                timingFull = cbind( t1, t2,t3),
             call = match.call())
           )
}