Changeset 17340
- Timestamp:
- 02/21/14 15:31:47 (11 years ago)
- File:
-
- 1 edited
Legend:
- Unmodified
- Added
- Removed
-
issm/trunk-jpl/src/c/classes/kriging/Observations.cpp
r17329 r17340 354 354 IssmPDouble prediction,error; 355 355 IssmPDouble numerator,denominator,ratio; 356 IssmPDouble *x = NULL; 357 IssmPDouble *y = NULL; 358 IssmPDouble *obs = NULL; 359 IssmPDouble *Gamma = NULL; 360 IssmPDouble *GinvG0 = NULL; 361 IssmPDouble *Ginv1 = NULL; 362 IssmPDouble *GinvZ = NULL; 363 IssmPDouble *gamma0 = NULL; 364 IssmPDouble *ones = NULL; 356 IssmPDouble *x = NULL; 357 IssmPDouble *y = NULL; 358 IssmPDouble *obs = NULL; 359 IssmPDouble *Lambda = NULL; 365 360 366 361 /*Some checks*/ … … 382 377 383 378 /*Allocate intermediary matrix and vectors*/ 384 Gamma = xNew<IssmPDouble>(n_obs*n_obs); 385 gamma0 = xNew<IssmPDouble>(n_obs); 386 ones = xNew<IssmPDouble>(n_obs); 379 IssmPDouble* A = xNew<IssmPDouble>((n_obs+1)*(n_obs+1)); 380 IssmPDouble* B = xNew<IssmPDouble>(n_obs+1); 387 381 388 382 /*First: Create semivariogram matrix for observations*/ 389 383 for(i=0;i<n_obs;i++){ 390 384 for(j=0;j<=i;j++){ 391 //Gamma[i*n_obs+j] = variogram->SemiVariogram(x[i]-x[j],y[i]-y[j]); 392 Gamma[i*n_obs+j] = variogram->Covariance(x[i]-x[j],y[i]-y[j]); 393 Gamma[j*n_obs+i] = Gamma[i*n_obs+j]; 394 } 395 } 396 for(i=0;i<n_obs;i++) ones[i]=1; 385 A[i*(n_obs+1)+j] = variogram->Covariance(x[i]-x[j],y[i]-y[j]); 386 A[j*(n_obs+1)+i] = A[i*(n_obs+1)+j]; 387 } 388 A[i*(n_obs+1)+n_obs] = 1.; 389 } 390 for(i=0;i<n_obs;i++) A[n_obs*(n_obs+1)+i]=1.; 391 A[n_obs*(n_obs+1)+n_obs] = 0.; 397 392 398 393 /*Get semivariogram vector associated to this location*/ 399 //for(i=0;i<n_obs;i++) gamma0[i] = variogram->SemiVariogram(x[i]-x_interp,y[i]-y_interp);400 for(i=0;i<n_obs;i++) gamma0[i] = variogram->Covariance(x[i]-x_interp,y[i]-y_interp);394 for(i=0;i<n_obs;i++) B[i] = variogram->Covariance(x[i]-x_interp,y[i]-y_interp); 395 B[n_obs] = 1.; 401 396 402 397 /*Solve the three linear systems*/ 403 398 #if _HAVE_GSL_ 404 DenseGslTripleSolve(&GinvG0,&Ginv1,&GinvZ,Gamma,gamma0,ones,obs,n_obs); 405 //DenseGslSolve(&GinvG0,Gamma,gamma0,n_obs); // Gamma^-1 gamma0 406 //DenseGslSolve(&Ginv1, Gamma,ones,n_obs); // Gamma^-1 ones 407 //DenseGslSolve(&GinvZ, Gamma,obs,n_obs); // Gamma^-1 Z 399 DenseGslSolve(&Lambda,A,B,n_obs+1); // Gamma^-1 Z 408 400 #else 409 401 _error_("GSL is required"); 410 402 #endif 411 403 412 /*Prepare predictor*/ 413 numerator=-1.; denominator=0.; 414 for(i=0;i<n_obs;i++) numerator +=GinvG0[i]; 415 for(i=0;i<n_obs;i++) denominator+=Ginv1[i]; 416 ratio=numerator/denominator; 417 404 /*Compute predictor*/ 418 405 prediction = 0.; 419 error = - numerator*numerator/denominator; 420 for(i=0;i<n_obs;i++) prediction += (gamma0[i]-ratio)*GinvZ[i]; 421 for(i=0;i<n_obs;i++) error += gamma0[i]*GinvG0[i]; 406 for(i=0;i<n_obs;i++) prediction += Lambda[i]*obs[i]; 407 408 /*Compute error (GSLIB p15 eq II.14)*/ 409 error = variogram->Covariance(0.,0.)*(1. - Lambda[n_obs]);; 410 for(i=0;i<n_obs;i++) error += -Lambda[i]*B[i]; 422 411 423 412 /*clean-up*/ … … 427 416 xDelete<IssmPDouble>(y); 428 417 xDelete<IssmPDouble>(obs); 429 xDelete<IssmPDouble>(Gamma); 430 xDelete<IssmPDouble>(gamma0); 431 xDelete<IssmPDouble>(ones); 432 xDelete<IssmPDouble>(GinvG0); 433 xDelete<IssmPDouble>(Ginv1); 434 xDelete<IssmPDouble>(GinvZ); 418 xDelete<IssmPDouble>(A); 419 xDelete<IssmPDouble>(B); 420 xDelete<IssmPDouble>(Lambda); 435 421 436 422 }/*}}}*/
Note:
See TracChangeset
for help on using the changeset viewer.