/* * \file Observations.cpp * \brief: Implementation of Observations class, derived from DataSet class. */ /*Headers: {{{*/ #ifdef HAVE_CONFIG_H #include #else #error "Cannot compile with HAVE_CONFIG_H symbol! run configure first!" #endif #include #include #include #include #include "../Options/Options.h" #include "./Observations.h" #include "./Observation.h" #include "../../datastructures/datastructures.h" #include "../../shared/shared.h" #include "./Quadtree.h" #include "./Variogram.h" #include "../../toolkits/toolkits.h" using namespace std; /*}}}*/ /*Object constructors and destructor*/ /*FUNCTION Observations::Observations(){{{*/ Observations::Observations(){ this->quadtree = NULL; return; } /*}}}*/ /*FUNCTION Observations::Observations(IssmPDouble* observations_list,IssmPDouble* x,IssmPDouble* y,int n,Options* options){{{*/ Observations::Observations(IssmPDouble* observations_list,IssmPDouble* x,IssmPDouble* y,int n,Options* options){ /*Intermediaries*/ int i,maxdepth,level,counter,index; int xi,yi; IssmPDouble xmin,xmax,ymin,ymax; IssmPDouble offset,minlength,minspacing,mintrimming,maxtrimming; Observation *observation = NULL; /*Check that observations is not empty*/ if(n==0) _error_("No observation found"); /*Get extrema*/ xmin=x[0]; ymin=y[0]; xmax=x[0]; ymax=y[0]; for(i=1;iGet(&mintrimming,"mintrimming",-1.e+21); options->Get(&maxtrimming,"maxtrimming",+1.e+21); options->Get(&minspacing,"minspacing",0.01); if(minspacing<=0) _error_("minspacing must > 0"); /*Get Minimum box size*/ if(options->GetOption("boxlength")){ options->Get(&minlength,"boxlength"); if(minlength<=0)_error_("boxlength should be a positive number"); maxdepth=reCast(log(max(xmax-xmin,ymax-ymin)/minlength +1)/log(2.0)); } else{ maxdepth = 30; minlength=max(xmax-xmin,ymax-ymin)/IssmPDouble((1L<quadtree = new Quadtree(xmin,xmax,ymin,ymax,maxdepth); /*Add observations one by one*/ counter = 0; for(i=0;imaxtrimming) continue; if(observations_list[i]quadtree->ClosestObs(&index,x[i],y[i]); if(index>=0){ observation=dynamic_cast(this->GetObjectByOffset(index)); if(pow(observation->x-x[i],2)+pow(observation->y-y[i],2) < minspacing) continue; } this->quadtree->IntergerCoordinates(&xi,&yi,x[i],y[i]); this->quadtree->QuadtreeDepth2(&level,xi,yi); if((int)level <= maxdepth){ observation = new Observation(x[i],y[i],xi,yi,counter++,observations_list[i]); this->quadtree->Add(observation); this->AddObject(observation); } else{ /*We need to average with the current observations*/ this->quadtree->AddAndAverage(x[i],y[i],observations_list[i]); } } _printf0_("done\n"); _printf0_("Initial number of observations: " << n << "\n"); _printf0_(" Final number of observations: " << this->quadtree->NbObs << "\n"); } /*}}}*/ /*FUNCTION Observations::~Observations(){{{*/ Observations::~Observations(){ delete quadtree; return; } /*}}}*/ /*Methods*/ /*FUNCTION Observations::ClosestObservation{{{*/ void Observations::ClosestObservation(IssmPDouble *px,IssmPDouble *py,IssmPDouble *pobs,IssmPDouble x_interp,IssmPDouble y_interp,IssmPDouble radius){ /*Output and Intermediaries*/ int nobs,i,index; IssmPDouble hmin,h2,hmin2,radius2; int *indices = NULL; Observation *observation = NULL; /*If radius is not provided or is 0, return all observations*/ if(radius==0) radius=this->quadtree->root->length; /*First, find closest point in Quadtree (fast but might not be the true closest obs)*/ this->quadtree->ClosestObs(&index,x_interp,y_interp); if(index>=0){ observation=dynamic_cast(this->GetObjectByOffset(index)); hmin = sqrt((observation->x-x_interp)*(observation->x-x_interp) + (observation->y-y_interp)*(observation->y-y_interp)); if(hminquadtree->RangeSearch(&indices,&nobs,x_interp,y_interp,radius); for (i=0;i(this->GetObjectByOffset(indices[i])); h2 = (observation->x-x_interp)*(observation->x-x_interp) + (observation->y-y_interp)*(observation->y-y_interp); if(i==0){ hmin2 = h2; index = indices[i]; } else{ if(h2=0){ observation=dynamic_cast(this->GetObjectByOffset(index)); *px=observation->x; *py=observation->y; *pobs=observation->value; } else{ *px=UNDEF; *py=UNDEF; *pobs=UNDEF; } xDelete(indices); }/*}}}*/ /*FUNCTION Observations::Distances{{{*/ void Observations::Distances(IssmPDouble* distances,IssmPDouble *x,IssmPDouble *y,int n,IssmPDouble radius){ IssmPDouble xi,yi,obs; for(int i=0;iClosestObservation(&xi,&yi,&obs,x[i],y[i],radius); if(xi==UNDEF && yi==UNDEF) distances[i]=UNDEF; else distances[i]=sqrt( (x[i]-xi)*(x[i]-xi) + (y[i]-yi)*(y[i]-yi) ); } }/*}}}*/ /*FUNCTION Observations::ObservationList(IssmPDouble **px,IssmPDouble **py,IssmPDouble **pobs,int* pnobs,IssmPDouble x_interp,IssmPDouble y_interp,IssmPDouble radius,int maxdata){{{*/ void Observations::ObservationList(IssmPDouble **px,IssmPDouble **py,IssmPDouble **pobs,int* pnobs,IssmPDouble x_interp,IssmPDouble y_interp,IssmPDouble radius,int maxdata){ /*Output and Intermediaries*/ bool stop; int nobs,tempnobs,i,j,k,n,counter; IssmPDouble h2,radius2; int *indices = NULL; int *tempindices = NULL; IssmPDouble *dists = NULL; IssmPDouble *x = NULL; IssmPDouble *y = NULL; IssmPDouble *obs = NULL; Observation *observation = NULL; /*If radius is not provided or is 0, return all observations*/ if(radius==0) radius=this->quadtree->root->length; /*Compute radius square*/ radius2 = radius*radius; /*Find all observations that are in radius*/ this->quadtree->RangeSearch(&tempindices,&tempnobs,x_interp,y_interp,radius); if(tempnobs){ indices = xNew(tempnobs); dists = xNew(tempnobs); } nobs = 0; for (i=0;i(this->GetObjectByOffset(tempindices[i])); h2 = (observation->x-x_interp)*(observation->x-x_interp) + (observation->y-y_interp)*(observation->y-y_interp); if(nobs==maxdata && h2>radius2) continue; if(nobs<=maxdata){ indices[nobs] = tempindices[i]; dists[nobs] = h2; nobs++; } if(nobs==1) continue; /*Sort all dists up to now*/ n=nobs-1; stop = false; for(k=0;k(dists); xDelete(tempindices); if(nobs){ /*Allocate vectors*/ x = xNew(nobs); y = xNew(nobs); obs = xNew(nobs); /*Loop over all observations and fill in x, y and obs*/ for (i=0;i(this->GetObjectByOffset(indices[i])); observation->WriteXYObs(&x[i],&y[i],&obs[i]); } } /*Assign output pointer*/ xDelete(indices); *px=x; *py=y; *pobs=obs; *pnobs=nobs; }/*}}}*/ /*FUNCTION Observations::ObservationList(IssmPDouble **px,IssmPDouble **py,IssmPDouble **pobs,int* pnobs){{{*/ void Observations::ObservationList(IssmPDouble **px,IssmPDouble **py,IssmPDouble **pobs,int* pnobs){ /*Output and Intermediaries*/ int nobs; IssmPDouble *x = NULL; IssmPDouble *y = NULL; IssmPDouble *obs = NULL; Observation *observation = NULL; nobs = this->Size(); if(nobs){ x = xNew(nobs); y = xNew(nobs); obs = xNew(nobs); for(int i=0;iSize();i++){ observation=dynamic_cast(this->GetObjectByOffset(i)); observation->WriteXYObs(&x[i],&y[i],&obs[i]); } } /*Assign output pointer*/ *px=x; *py=y; *pobs=obs; *pnobs=nobs; }/*}}}*/ /*FUNCTION Observations::InterpolationIDW{{{*/ void Observations::InterpolationIDW(IssmPDouble *pprediction,IssmPDouble x_interp,IssmPDouble y_interp,IssmPDouble radius,int mindata,int maxdata,IssmPDouble power){ /*Intermediaries*/ int i,n_obs; IssmPDouble prediction; IssmPDouble numerator,denominator,h,weight; IssmPDouble *x = NULL; IssmPDouble *y = NULL; IssmPDouble *obs = NULL; /*Some checks*/ _assert_(maxdata>0); _assert_(pprediction); _assert_(power>0); /*If radius is not provided or is 0, return all observations*/ if(radius==0) radius=this->quadtree->root->length; /*Get list of observations for current point*/ this->ObservationList(&x,&y,&obs,&n_obs,x_interp,y_interp,radius,maxdata); /*If we have less observations than mindata, return UNDEF*/ if(n_obs(x); xDelete(y); xDelete(obs); }/*}}}*/ /*FUNCTION Observations::InterpolationKriging{{{*/ void Observations::InterpolationKriging(IssmPDouble *pprediction,IssmPDouble *perror,IssmPDouble x_interp,IssmPDouble y_interp,IssmPDouble radius,int mindata,int maxdata,Variogram* variogram){ /*Intermediaries*/ int i,j,n_obs; IssmPDouble prediction,error; IssmPDouble numerator,denominator,ratio; IssmPDouble *x = NULL; IssmPDouble *y = NULL; IssmPDouble *obs = NULL; IssmPDouble *Gamma = NULL; IssmPDouble *GinvG0 = NULL; IssmPDouble *Ginv1 = NULL; IssmPDouble *GinvZ = NULL; IssmPDouble *gamma0 = NULL; IssmPDouble *ones = NULL; /*Some checks*/ _assert_(mindata>0 && maxdata>0); _assert_(pprediction && perror); /*If radius is not provided or is 0, return all observations*/ if(radius==0) radius=this->quadtree->root->length; /*Get list of observations for current point*/ this->ObservationList(&x,&y,&obs,&n_obs,x_interp,y_interp,radius,maxdata); /*If we have less observations than mindata, return UNDEF*/ if(n_obs(n_obs*n_obs); gamma0 = xNew(n_obs); ones = xNew(n_obs); /*First: Create semivariogram matrix for observations*/ for(i=0;iSemiVariogram(x[i]-x[j],y[i]-y[j]); Gamma[i*n_obs+j] = variogram->Covariance(x[i]-x[j],y[i]-y[j]); Gamma[j*n_obs+i] = Gamma[i*n_obs+j]; } } for(i=0;iSemiVariogram(x[i]-x_interp,y[i]-y_interp); for(i=0;iCovariance(x[i]-x_interp,y[i]-y_interp); /*Solve the three linear systems*/ #if _HAVE_GSL_ DenseGslSolve(&GinvG0,Gamma,gamma0,n_obs); // Gamma^-1 gamma0 DenseGslSolve(&Ginv1, Gamma,ones,n_obs); // Gamma^-1 ones DenseGslSolve(&GinvZ, Gamma,obs,n_obs); // Gamma^-1 Z #else _error_("GSL is required"); #endif /*Prepare predictor*/ numerator=-1.; denominator=0.; for(i=0;i(x); xDelete(y); xDelete(obs); xDelete(Gamma); xDelete(gamma0); xDelete(ones); xDelete(GinvG0); xDelete(Ginv1); xDelete(GinvZ); }/*}}}*/ /*FUNCTION Observations::InterpolationNearestNeighbor{{{*/ void Observations::InterpolationNearestNeighbor(IssmPDouble *pprediction,IssmPDouble x_interp,IssmPDouble y_interp,IssmPDouble radius){ /*Intermediaries*/ IssmPDouble x,y,obs; /*Get clostest observation*/ this->ClosestObservation(&x,&y,&obs,x_interp,y_interp,radius); /*Assign output pointer*/ *pprediction = obs; }/*}}}*/ /*FUNCTION Observations::InterpolationV4{{{*/ void Observations::InterpolationV4(IssmPDouble *pprediction,IssmPDouble x_interp,IssmPDouble y_interp,IssmPDouble radius,int mindata,int maxdata){ /* Reference: David T. Sandwell, Biharmonic spline interpolation of GEOS-3 * and SEASAT altimeter data, Geophysical Research Letters, 2, 139-142, * 1987. Describes interpolation using value or gradient of value in any * dimension.*/ /*Intermediaries*/ int i,j,n_obs; IssmPDouble prediction,h; IssmPDouble *x = NULL; IssmPDouble *y = NULL; IssmPDouble *obs = NULL; IssmPDouble *Green = NULL; IssmPDouble *weights = NULL; IssmPDouble *g = NULL; /*Some checks*/ _assert_(maxdata>0); _assert_(pprediction); /*If radius is not provided or is 0, return all observations*/ if(radius==0) radius=this->quadtree->root->length; /*Get list of observations for current point*/ this->ObservationList(&x,&y,&obs,&n_obs,x_interp,y_interp,radius,maxdata); /*If we have less observations than mindata, return UNDEF*/ if(n_obs(n_obs*n_obs); g = xNew(n_obs); /*First: distance vector*/ for(i=0;i0){ g[i] = h*h*(log(h)-1.); } else{ g[i] = 0.; } } /*Build Green function matrix*/ for(i=0;i0){ Green[j*n_obs+i] = h*h*(log(h)-1.); } else{ Green[j*n_obs+i] = 0.; } Green[i*n_obs+j] = Green[j*n_obs+i]; } } /*Compute weights*/ #if _HAVE_GSL_ DenseGslSolve(&weights,Green,obs,n_obs); // Green^-1 obs #else _error_("GSL is required"); #endif /*Interpolate*/ prediction = 0; for(i=0;i(x); xDelete(y); xDelete(obs); xDelete(Green); xDelete(g); xDelete(weights); }/*}}}*/ /*FUNCTION Observations::QuadtreeColoring{{{*/ void Observations::QuadtreeColoring(IssmPDouble* A,IssmPDouble *x,IssmPDouble *y,int n){ int xi,yi,level; for(int i=0;iquadtree->IntergerCoordinates(&xi,&yi,x[i],y[i]); this->quadtree->QuadtreeDepth(&level,xi,yi); A[i]=(IssmPDouble)level; } }/*}}}*/ /*FUNCTION Observations::Variomap{{{*/ void Observations::Variomap(IssmPDouble* gamma,IssmPDouble *x,int n){ /*Output and Intermediaries*/ int i,j,k; IssmPDouble distance; Observation *observation1 = NULL; Observation *observation2 = NULL; IssmPDouble *counter = xNew(n); for(j=0;jSize();i++){ observation1=dynamic_cast(this->GetObjectByOffset(i)); for(j=i+1;jSize();j++){ observation2=dynamic_cast(this->GetObjectByOffset(j)); distance=sqrt(pow(observation1->x - observation2->x,2) + pow(observation1->y - observation2->y,2)); if(distance>x[n-1]) continue; int index = int(distance/(x[1]-x[0])); if(index>n-1) index = n-1; if(index<0) index = 0; gamma[index] += 1./2.*pow(observation1->value - observation2->value,2); counter[index] += 1.; } } /*Normalize semivariogram*/ gamma[0]=0.; for(k=0;k(counter); }/*}}}*/