[14996] | 1 | /*
|
---|
| 2 | * \file Observations.cpp
|
---|
| 3 | * \brief: Implementation of Observations class, derived from DataSet class.
|
---|
| 4 | */
|
---|
| 5 |
|
---|
| 6 | /*Headers: {{{*/
|
---|
| 7 | #ifdef HAVE_CONFIG_H
|
---|
| 8 | #include <config.h>
|
---|
| 9 | #else
|
---|
| 10 | #error "Cannot compile with HAVE_CONFIG_H symbol! run configure first!"
|
---|
| 11 | #endif
|
---|
| 12 |
|
---|
| 13 | #include <vector>
|
---|
| 14 | #include <functional>
|
---|
| 15 | #include <algorithm>
|
---|
| 16 | #include <iostream>
|
---|
| 17 |
|
---|
[15012] | 18 | #include "../Options/Options.h"
|
---|
[14996] | 19 | #include "./Observations.h"
|
---|
| 20 | #include "./Observation.h"
|
---|
[15067] | 21 | #include "../../datastructures/datastructures.h"
|
---|
[14996] | 22 | #include "../../shared/shared.h"
|
---|
| 23 |
|
---|
| 24 | #include "./Quadtree.h"
|
---|
[19105] | 25 | #include "./Covertree.h"
|
---|
[14996] | 26 | #include "./Variogram.h"
|
---|
[14997] | 27 | #include "../../toolkits/toolkits.h"
|
---|
[14996] | 28 |
|
---|
| 29 | using namespace std;
|
---|
| 30 | /*}}}*/
|
---|
| 31 |
|
---|
| 32 | /*Object constructors and destructor*/
|
---|
[18301] | 33 | Observations::Observations(){/*{{{*/
|
---|
[19105] | 34 | this->treetype = 0;
|
---|
| 35 | this->quadtree = NULL;
|
---|
| 36 | this->covertree = NULL;
|
---|
[14996] | 37 | return;
|
---|
| 38 | }
|
---|
| 39 | /*}}}*/
|
---|
[18301] | 40 | Observations::Observations(IssmPDouble* observations_list,IssmPDouble* x,IssmPDouble* y,int n,Options* options){/*{{{*/
|
---|
[14996] | 41 |
|
---|
[19105] | 42 | /*Check that there are observations*/
|
---|
| 43 | if(n<=0) _error_("No observation found");
|
---|
| 44 |
|
---|
| 45 | /*Get tree type (FIXME)*/
|
---|
[25836] | 46 | double dtree = 0.;
|
---|
[19105] | 47 | options->Get(&dtree,"treetype",1.);
|
---|
| 48 | this->treetype = reCast<int>(dtree);
|
---|
| 49 | switch(this->treetype){
|
---|
| 50 | case 1:
|
---|
| 51 | this->covertree = NULL;
|
---|
| 52 | this->InitQuadtree(observations_list,x,y,n,options);
|
---|
| 53 | break;
|
---|
| 54 | case 2:
|
---|
| 55 | this->quadtree = NULL;
|
---|
| 56 | this->InitCovertree(observations_list,x,y,n,options);
|
---|
| 57 | break;
|
---|
| 58 | default:
|
---|
| 59 | _error_("Tree type "<<this->treetype<<" not supported yet (1: quadtree, 2: covertree)");
|
---|
| 60 | }
|
---|
| 61 | }
|
---|
| 62 | /*}}}*/
|
---|
| 63 | Observations::~Observations(){/*{{{*/
|
---|
| 64 | switch(this->treetype){
|
---|
| 65 | case 1:
|
---|
| 66 | delete this->quadtree;
|
---|
| 67 | break;
|
---|
| 68 | case 2:
|
---|
| 69 | delete this->covertree;
|
---|
| 70 | break;
|
---|
| 71 | default:
|
---|
[24686] | 72 | _printf_("Tree type "<<this->treetype<<" not supported yet (1: quadtree, 2: covertree)");
|
---|
[19105] | 73 | }
|
---|
| 74 | return;
|
---|
| 75 | }
|
---|
| 76 | /*}}}*/
|
---|
| 77 |
|
---|
| 78 | /*Initialize data structures*/
|
---|
[21341] | 79 | void Observations::InitCovertree(IssmPDouble* observations_list,IssmPDouble* x,IssmPDouble* y,int n,Options* options){/*{{{*/
|
---|
| 80 |
|
---|
| 81 | /*Intermediaries*/
|
---|
| 82 | IssmPDouble minspacing,mintrimming,maxtrimming;
|
---|
| 83 |
|
---|
| 84 | /*Checks*/
|
---|
| 85 | _assert_(n);
|
---|
| 86 |
|
---|
| 87 | /*Get trimming limits*/
|
---|
| 88 | options->Get(&mintrimming,"mintrimming",-1.e+21);
|
---|
| 89 | options->Get(&maxtrimming,"maxtrimming",+1.e+21);
|
---|
| 90 | options->Get(&minspacing,"minspacing",0.01);
|
---|
| 91 | if(minspacing<=0) _error_("minspacing must > 0");
|
---|
| 92 |
|
---|
| 93 | /*Get maximum distance between 2 points
|
---|
| 94 | * maxDist should be the maximum distance that any two points
|
---|
| 95 | * can have between each other. IE p.distance(q) < maxDist for all
|
---|
| 96 | * p,q that you will ever try to insert. The cover tree may be invalid
|
---|
| 97 | * if an inaccurate maxDist is given.*/
|
---|
| 98 | IssmPDouble xmin = x[0];
|
---|
| 99 | IssmPDouble xmax = x[0];
|
---|
| 100 | IssmPDouble ymin = y[0];
|
---|
| 101 | IssmPDouble ymax = y[0];
|
---|
| 102 | for(int i=1;i<n;i++){
|
---|
| 103 | if(x[i]<xmin) xmin=x[i];
|
---|
| 104 | if(x[i]>xmax) xmax=x[i];
|
---|
| 105 | if(y[i]<ymin) ymin=y[i];
|
---|
| 106 | if(y[i]>ymax) ymax=y[i];
|
---|
| 107 | }
|
---|
| 108 | IssmPDouble maxDist = sqrt(pow(xmax-xmin,2)+pow(ymax-ymin,2));
|
---|
| 109 | IssmPDouble base = 2.;
|
---|
| 110 | int maxdepth = ceilf(log(maxDist)/log(base));
|
---|
| 111 |
|
---|
| 112 | _printf0_("Generating covertree with a maximum depth " << maxdepth <<"... ");
|
---|
| 113 | this->covertree=new Covertree(maxdepth);
|
---|
| 114 |
|
---|
| 115 | for(int i=0;i<n;i++){
|
---|
| 116 |
|
---|
| 117 | /*First check limits*/
|
---|
| 118 | if(observations_list[i]>maxtrimming) continue;
|
---|
| 119 | if(observations_list[i]<mintrimming) continue;
|
---|
| 120 |
|
---|
| 121 | /*Second, check that this observation is not too close from another one*/
|
---|
| 122 | Observation newobs = Observation(x[i],y[i],observations_list[i]);
|
---|
| 123 | if(i>0 && this->covertree->getRoot()){
|
---|
| 124 | /*Get closest obs and see if it is too close*/
|
---|
| 125 | std::vector<Observation> kNN=(this->covertree->kNearestNeighbors(newobs,1));
|
---|
| 126 | Observation oldobs = (*kNN.begin());
|
---|
| 127 | if(oldobs.distance(newobs)<minspacing) continue;
|
---|
| 128 | }
|
---|
| 129 |
|
---|
| 130 | this->covertree->insert(newobs);
|
---|
| 131 | }
|
---|
| 132 | _printf0_("done\n");
|
---|
| 133 | }
|
---|
| 134 | /*}}}*/
|
---|
[19105] | 135 | void Observations::InitQuadtree(IssmPDouble* observations_list,IssmPDouble* x,IssmPDouble* y,int n,Options* options){/*{{{*/
|
---|
| 136 |
|
---|
[14996] | 137 | /*Intermediaries*/
|
---|
| 138 | int i,maxdepth,level,counter,index;
|
---|
| 139 | int xi,yi;
|
---|
| 140 | IssmPDouble xmin,xmax,ymin,ymax;
|
---|
| 141 | IssmPDouble offset,minlength,minspacing,mintrimming,maxtrimming;
|
---|
| 142 | Observation *observation = NULL;
|
---|
| 143 |
|
---|
[19105] | 144 | /*Checks*/
|
---|
| 145 | _assert_(n);
|
---|
[14996] | 146 |
|
---|
| 147 | /*Get extrema*/
|
---|
| 148 | xmin=x[0]; ymin=y[0];
|
---|
| 149 | xmax=x[0]; ymax=y[0];
|
---|
| 150 | for(i=1;i<n;i++){
|
---|
| 151 | xmin=min(xmin,x[i]); ymin=min(ymin,y[i]);
|
---|
| 152 | xmax=max(xmax,x[i]); ymax=max(ymax,y[i]);
|
---|
| 153 | }
|
---|
| 154 | offset=0.05*(xmax-xmin); xmin-=offset; xmax+=offset;
|
---|
| 155 | offset=0.05*(ymax-ymin); ymin-=offset; ymax+=offset;
|
---|
| 156 |
|
---|
| 157 | /*Get trimming limits*/
|
---|
| 158 | options->Get(&mintrimming,"mintrimming",-1.e+21);
|
---|
| 159 | options->Get(&maxtrimming,"maxtrimming",+1.e+21);
|
---|
| 160 | options->Get(&minspacing,"minspacing",0.01);
|
---|
| 161 | if(minspacing<=0) _error_("minspacing must > 0");
|
---|
| 162 |
|
---|
| 163 | /*Get Minimum box size*/
|
---|
| 164 | if(options->GetOption("boxlength")){
|
---|
| 165 | options->Get(&minlength,"boxlength");
|
---|
| 166 | if(minlength<=0)_error_("boxlength should be a positive number");
|
---|
| 167 | maxdepth=reCast<int,IssmPDouble>(log(max(xmax-xmin,ymax-ymin)/minlength +1)/log(2.0));
|
---|
| 168 | }
|
---|
| 169 | else{
|
---|
| 170 | maxdepth = 30;
|
---|
| 171 | minlength=max(xmax-xmin,ymax-ymin)/IssmPDouble((1L<<maxdepth)-1);
|
---|
| 172 | }
|
---|
| 173 |
|
---|
| 174 | /*Initialize Quadtree*/
|
---|
[15100] | 175 | _printf0_("Generating quadtree with a maximum box size " << minlength << " (depth=" << maxdepth << ")... ");
|
---|
[14996] | 176 | this->quadtree = new Quadtree(xmin,xmax,ymin,ymax,maxdepth);
|
---|
| 177 |
|
---|
| 178 | /*Add observations one by one*/
|
---|
| 179 | counter = 0;
|
---|
| 180 | for(i=0;i<n;i++){
|
---|
| 181 |
|
---|
| 182 | /*First check limits*/
|
---|
| 183 | if(observations_list[i]>maxtrimming) continue;
|
---|
| 184 | if(observations_list[i]<mintrimming) continue;
|
---|
| 185 |
|
---|
[20500] | 186 | /*Second, check that this observation is not too close from another one*/
|
---|
[14996] | 187 | this->quadtree->ClosestObs(&index,x[i],y[i]);
|
---|
| 188 | if(index>=0){
|
---|
[19105] | 189 | observation=xDynamicCast<Observation*>(this->GetObjectByOffset(index));
|
---|
[14996] | 190 | if(pow(observation->x-x[i],2)+pow(observation->y-y[i],2) < minspacing) continue;
|
---|
| 191 | }
|
---|
| 192 |
|
---|
| 193 | this->quadtree->IntergerCoordinates(&xi,&yi,x[i],y[i]);
|
---|
| 194 | this->quadtree->QuadtreeDepth2(&level,xi,yi);
|
---|
| 195 | if((int)level <= maxdepth){
|
---|
| 196 | observation = new Observation(x[i],y[i],xi,yi,counter++,observations_list[i]);
|
---|
| 197 | this->quadtree->Add(observation);
|
---|
| 198 | this->AddObject(observation);
|
---|
| 199 | }
|
---|
| 200 | else{
|
---|
| 201 | /*We need to average with the current observations*/
|
---|
| 202 | this->quadtree->AddAndAverage(x[i],y[i],observations_list[i]);
|
---|
| 203 | }
|
---|
| 204 | }
|
---|
[15104] | 205 | _printf0_("done\n");
|
---|
[15100] | 206 | _printf0_("Initial number of observations: " << n << "\n");
|
---|
| 207 | _printf0_(" Final number of observations: " << this->quadtree->NbObs << "\n");
|
---|
[14996] | 208 | }
|
---|
| 209 | /*}}}*/
|
---|
[19105] | 210 |
|
---|
[14996] | 211 | /*Methods*/
|
---|
[18301] | 212 | void Observations::ClosestObservation(IssmPDouble *px,IssmPDouble *py,IssmPDouble *pobs,IssmPDouble x_interp,IssmPDouble y_interp,IssmPDouble radius){/*{{{*/
|
---|
[14996] | 213 |
|
---|
[20500] | 214 | switch(this->treetype){
|
---|
| 215 | case 1:
|
---|
| 216 | this->ClosestObservationQuadtree(px,py,pobs,x_interp,y_interp,radius);
|
---|
| 217 | break;
|
---|
| 218 | case 2:
|
---|
| 219 | this->ClosestObservationCovertree(px,py,pobs,x_interp,y_interp,radius);
|
---|
| 220 | break;
|
---|
| 221 | default:
|
---|
| 222 | _error_("Tree type "<<this->treetype<<" not supported yet (1: quadtree, 2: covertree)");
|
---|
| 223 | }
|
---|
| 224 |
|
---|
| 225 | }/*}}}*/
|
---|
[21341] | 226 | void Observations::ClosestObservationCovertree(IssmPDouble *px,IssmPDouble *py,IssmPDouble *pobs,IssmPDouble x_interp,IssmPDouble y_interp,IssmPDouble radius){/*{{{*/
|
---|
| 227 |
|
---|
| 228 | IssmPDouble hmin = UNDEF;
|
---|
| 229 |
|
---|
| 230 | if(this->covertree->getRoot()){
|
---|
| 231 | /*Get closest obs and see if it is too close*/
|
---|
| 232 | Observation newobs = Observation(x_interp,y_interp,0.);
|
---|
| 233 | std::vector<Observation> kNN=(this->covertree->kNearestNeighbors(newobs,1));
|
---|
| 234 | Observation observation = (*kNN.begin());
|
---|
| 235 | hmin = observation.distance(newobs);
|
---|
| 236 | if(hmin<=radius){
|
---|
| 237 | *px = observation.x;
|
---|
| 238 | *py = observation.y;
|
---|
| 239 | *pobs = observation.value;
|
---|
| 240 | return;
|
---|
| 241 | }
|
---|
| 242 | }
|
---|
| 243 |
|
---|
| 244 | *px = UNDEF;
|
---|
| 245 | *py = UNDEF;
|
---|
| 246 | *pobs = UNDEF;
|
---|
| 247 | }/*}}}*/
|
---|
[20500] | 248 | void Observations::ClosestObservationQuadtree(IssmPDouble *px,IssmPDouble *py,IssmPDouble *pobs,IssmPDouble x_interp,IssmPDouble y_interp,IssmPDouble radius){/*{{{*/
|
---|
| 249 |
|
---|
[14996] | 250 | /*Output and Intermediaries*/
|
---|
| 251 | int nobs,i,index;
|
---|
[16560] | 252 | IssmPDouble hmin,h2,hmin2;
|
---|
[14996] | 253 | int *indices = NULL;
|
---|
| 254 | Observation *observation = NULL;
|
---|
| 255 |
|
---|
| 256 | /*If radius is not provided or is 0, return all observations*/
|
---|
| 257 | if(radius==0) radius=this->quadtree->root->length;
|
---|
| 258 |
|
---|
[22758] | 259 | /*For CPPcheck*/
|
---|
| 260 | hmin = 2*radius;
|
---|
| 261 |
|
---|
[14996] | 262 | /*First, find closest point in Quadtree (fast but might not be the true closest obs)*/
|
---|
| 263 | this->quadtree->ClosestObs(&index,x_interp,y_interp);
|
---|
| 264 | if(index>=0){
|
---|
[19105] | 265 | observation=xDynamicCast<Observation*>(this->GetObjectByOffset(index));
|
---|
[14996] | 266 | hmin = sqrt((observation->x-x_interp)*(observation->x-x_interp) + (observation->y-y_interp)*(observation->y-y_interp));
|
---|
| 267 | if(hmin<radius) radius=hmin;
|
---|
| 268 | }
|
---|
| 269 |
|
---|
| 270 | /*Find all observations that are in radius*/
|
---|
| 271 | this->quadtree->RangeSearch(&indices,&nobs,x_interp,y_interp,radius);
|
---|
| 272 | for (i=0;i<nobs;i++){
|
---|
[19105] | 273 | observation=xDynamicCast<Observation*>(this->GetObjectByOffset(indices[i]));
|
---|
[14996] | 274 | h2 = (observation->x-x_interp)*(observation->x-x_interp) + (observation->y-y_interp)*(observation->y-y_interp);
|
---|
| 275 | if(i==0){
|
---|
| 276 | hmin2 = h2;
|
---|
| 277 | index = indices[i];
|
---|
| 278 | }
|
---|
| 279 | else{
|
---|
| 280 | if(h2<hmin2){
|
---|
| 281 | hmin2 = h2;
|
---|
| 282 | index = indices[i];
|
---|
| 283 | }
|
---|
| 284 | }
|
---|
[25836] | 285 | }
|
---|
[14996] | 286 |
|
---|
| 287 | /*Assign output pointer*/
|
---|
[17806] | 288 | if(nobs || hmin==radius){
|
---|
[19105] | 289 | observation=xDynamicCast<Observation*>(this->GetObjectByOffset(index));
|
---|
[17806] | 290 | *px = observation->x;
|
---|
| 291 | *py = observation->y;
|
---|
| 292 | *pobs = observation->value;
|
---|
[14996] | 293 | }
|
---|
| 294 | else{
|
---|
[17806] | 295 | *px = UNDEF;
|
---|
| 296 | *py = UNDEF;
|
---|
| 297 | *pobs = UNDEF;
|
---|
[14996] | 298 | }
|
---|
| 299 | xDelete<int>(indices);
|
---|
| 300 |
|
---|
| 301 | }/*}}}*/
|
---|
[18301] | 302 | void Observations::Distances(IssmPDouble* distances,IssmPDouble *x,IssmPDouble *y,int n,IssmPDouble radius){/*{{{*/
|
---|
[14996] | 303 |
|
---|
| 304 | IssmPDouble xi,yi,obs;
|
---|
| 305 |
|
---|
| 306 | for(int i=0;i<n;i++){
|
---|
| 307 | this->ClosestObservation(&xi,&yi,&obs,x[i],y[i],radius);
|
---|
[17806] | 308 | if(xi==UNDEF && yi==UNDEF){
|
---|
[14996] | 309 | distances[i]=UNDEF;
|
---|
[17806] | 310 | }
|
---|
| 311 | else{
|
---|
[14996] | 312 | distances[i]=sqrt( (x[i]-xi)*(x[i]-xi) + (y[i]-yi)*(y[i]-yi) );
|
---|
[17806] | 313 | }
|
---|
[14996] | 314 | }
|
---|
| 315 | }/*}}}*/
|
---|
[20500] | 316 | void Observations::ObservationList(IssmPDouble **px,IssmPDouble **py,IssmPDouble **pobs,int* pnobs){/*{{{*/
|
---|
| 317 |
|
---|
| 318 | /*Output and Intermediaries*/
|
---|
| 319 | int nobs;
|
---|
| 320 | IssmPDouble *x = NULL;
|
---|
| 321 | IssmPDouble *y = NULL;
|
---|
| 322 | IssmPDouble *obs = NULL;
|
---|
| 323 | Observation *observation = NULL;
|
---|
| 324 |
|
---|
| 325 | nobs = this->Size();
|
---|
| 326 |
|
---|
| 327 | if(nobs){
|
---|
| 328 | x = xNew<IssmPDouble>(nobs);
|
---|
| 329 | y = xNew<IssmPDouble>(nobs);
|
---|
| 330 | obs = xNew<IssmPDouble>(nobs);
|
---|
[25836] | 331 | int i=0;
|
---|
| 332 | for(Object* & object: this->objects){
|
---|
| 333 | observation=xDynamicCast<Observation*>(object);
|
---|
[20500] | 334 | observation->WriteXYObs(&x[i],&y[i],&obs[i]);
|
---|
[25836] | 335 | i++;
|
---|
[20500] | 336 | }
|
---|
| 337 | }
|
---|
| 338 |
|
---|
| 339 | /*Assign output pointer*/
|
---|
| 340 | *px=x;
|
---|
| 341 | *py=y;
|
---|
| 342 | *pobs=obs;
|
---|
| 343 | *pnobs=nobs;
|
---|
| 344 | }/*}}}*/
|
---|
[18301] | 345 | void Observations::ObservationList(IssmPDouble **px,IssmPDouble **py,IssmPDouble **pobs,int* pnobs,IssmPDouble x_interp,IssmPDouble y_interp,IssmPDouble radius,int maxdata){/*{{{*/
|
---|
[14996] | 346 |
|
---|
[20500] | 347 | switch(this->treetype){
|
---|
| 348 | case 1:
|
---|
| 349 | this->ObservationListQuadtree(px,py,pobs,pnobs,x_interp,y_interp,radius,maxdata);
|
---|
| 350 | break;
|
---|
| 351 | case 2:
|
---|
| 352 | this->ObservationListCovertree(px,py,pobs,pnobs,x_interp,y_interp,radius,maxdata);
|
---|
| 353 | break;
|
---|
| 354 | default:
|
---|
| 355 | _error_("Tree type "<<this->treetype<<" not supported yet (1: quadtree, 2: covertree)");
|
---|
[19105] | 356 | }
|
---|
[20500] | 357 | }/*}}}*/
|
---|
[21341] | 358 | void Observations::ObservationListCovertree(double **px,double **py,double **pobs,int* pnobs,double x_interp,double y_interp,double radius,int maxdata){/*{{{*/
|
---|
| 359 |
|
---|
| 360 | double *x = NULL;
|
---|
| 361 | double *y = NULL;
|
---|
| 362 | double *obs = NULL;
|
---|
| 363 | Observation observation=Observation(x_interp,y_interp,0.);
|
---|
| 364 | std::vector<Observation> kNN;
|
---|
| 365 |
|
---|
| 366 | kNN=(this->covertree->kNearestNeighbors(observation, maxdata));
|
---|
| 367 | //cout << "kNN's size: " << kNN.size() << " (maxdata = " <<maxdata<<")"<<endl;
|
---|
| 368 |
|
---|
| 369 | //kNN is sort from closest to farthest neighbor
|
---|
| 370 | //searches for the first neighbor that is out of radius
|
---|
| 371 | //deletes and resizes the kNN vector
|
---|
| 372 | vector<Observation>::iterator it;
|
---|
| 373 | if(radius>0.){
|
---|
| 374 | for (it = kNN.begin(); it != kNN.end(); ++it) {
|
---|
| 375 | //(*it).print();
|
---|
| 376 | //cout << "\n" << (*it).distance(observation) << endl;
|
---|
| 377 | if ((*it).distance(observation) > radius) {
|
---|
| 378 | break;
|
---|
| 379 | }
|
---|
| 380 | }
|
---|
| 381 | kNN.erase(it, kNN.end());
|
---|
| 382 | }
|
---|
| 383 |
|
---|
| 384 | /*Allocate vectors*/
|
---|
| 385 | x = new double[kNN.size()];
|
---|
| 386 | y = new double[kNN.size()];
|
---|
| 387 | obs = new double[kNN.size()];
|
---|
| 388 |
|
---|
| 389 | /*Loop over all observations and fill in x, y and obs*/
|
---|
| 390 | int i = 0;
|
---|
| 391 | for(it = kNN.begin(); it != kNN.end(); ++it) {
|
---|
| 392 | (*it).WriteXYObs((*it), &x[i], &y[i], &obs[i]);
|
---|
| 393 | i++;
|
---|
| 394 | }
|
---|
| 395 |
|
---|
| 396 | *px=x;
|
---|
| 397 | *py=y;
|
---|
| 398 | *pobs=obs;
|
---|
| 399 | *pnobs = kNN.size();
|
---|
| 400 | }/*}}}*/
|
---|
[20500] | 401 | void Observations::ObservationListQuadtree(IssmPDouble **px,IssmPDouble **py,IssmPDouble **pobs,int* pnobs,IssmPDouble x_interp,IssmPDouble y_interp,IssmPDouble radius,int maxdata){/*{{{*/
|
---|
[19105] | 402 |
|
---|
[14996] | 403 | /*Output and Intermediaries*/
|
---|
| 404 | bool stop;
|
---|
| 405 | int nobs,tempnobs,i,j,k,n,counter;
|
---|
| 406 | IssmPDouble h2,radius2;
|
---|
| 407 | int *indices = NULL;
|
---|
| 408 | int *tempindices = NULL;
|
---|
| 409 | IssmPDouble *dists = NULL;
|
---|
| 410 | IssmPDouble *x = NULL;
|
---|
| 411 | IssmPDouble *y = NULL;
|
---|
| 412 | IssmPDouble *obs = NULL;
|
---|
| 413 | Observation *observation = NULL;
|
---|
| 414 |
|
---|
| 415 | /*If radius is not provided or is 0, return all observations*/
|
---|
[17806] | 416 | if(radius==0.) radius=this->quadtree->root->length*2.;
|
---|
[14996] | 417 |
|
---|
| 418 | /*Compute radius square*/
|
---|
| 419 | radius2 = radius*radius;
|
---|
| 420 |
|
---|
| 421 | /*Find all observations that are in radius*/
|
---|
| 422 | this->quadtree->RangeSearch(&tempindices,&tempnobs,x_interp,y_interp,radius);
|
---|
| 423 | if(tempnobs){
|
---|
| 424 | indices = xNew<int>(tempnobs);
|
---|
| 425 | dists = xNew<IssmPDouble>(tempnobs);
|
---|
| 426 | }
|
---|
| 427 | nobs = 0;
|
---|
[17806] | 428 | for(i=0;i<tempnobs;i++){
|
---|
[19105] | 429 | observation=xDynamicCast<Observation*>(this->GetObjectByOffset(tempindices[i]));
|
---|
[14996] | 430 | h2 = (observation->x-x_interp)*(observation->x-x_interp) + (observation->y-y_interp)*(observation->y-y_interp);
|
---|
| 431 |
|
---|
| 432 | if(nobs==maxdata && h2>radius2) continue;
|
---|
[17806] | 433 | if(nobs<maxdata){
|
---|
[14996] | 434 | indices[nobs] = tempindices[i];
|
---|
| 435 | dists[nobs] = h2;
|
---|
| 436 | nobs++;
|
---|
| 437 | }
|
---|
| 438 | if(nobs==1) continue;
|
---|
| 439 |
|
---|
| 440 | /*Sort all dists up to now*/
|
---|
| 441 | n=nobs-1;
|
---|
| 442 | stop = false;
|
---|
| 443 | for(k=0;k<n-1;k++){
|
---|
| 444 | if(h2<dists[k]){
|
---|
| 445 | counter=1;
|
---|
| 446 | for(int jj=k;jj<n;jj++){
|
---|
| 447 | j = n-counter;
|
---|
| 448 | dists[j+1] = dists[j];
|
---|
| 449 | indices[j+1] = indices[j];
|
---|
| 450 | counter++;
|
---|
| 451 | }
|
---|
| 452 | dists[k] = h2;
|
---|
| 453 | indices[k] = tempindices[i];
|
---|
| 454 | stop = true;
|
---|
| 455 | break;
|
---|
| 456 | }
|
---|
| 457 | if(stop) break;
|
---|
| 458 | }
|
---|
[25836] | 459 | }
|
---|
[14996] | 460 | xDelete<IssmPDouble>(dists);
|
---|
| 461 | xDelete<int>(tempindices);
|
---|
| 462 |
|
---|
| 463 | if(nobs){
|
---|
| 464 | /*Allocate vectors*/
|
---|
| 465 | x = xNew<IssmPDouble>(nobs);
|
---|
| 466 | y = xNew<IssmPDouble>(nobs);
|
---|
| 467 | obs = xNew<IssmPDouble>(nobs);
|
---|
| 468 |
|
---|
| 469 | /*Loop over all observations and fill in x, y and obs*/
|
---|
[17806] | 470 | for(i=0;i<nobs;i++){
|
---|
[19105] | 471 | observation=xDynamicCast<Observation*>(this->GetObjectByOffset(indices[i]));
|
---|
[14996] | 472 | observation->WriteXYObs(&x[i],&y[i],&obs[i]);
|
---|
| 473 | }
|
---|
| 474 | }
|
---|
| 475 |
|
---|
| 476 | /*Assign output pointer*/
|
---|
| 477 | xDelete<int>(indices);
|
---|
| 478 | *px=x;
|
---|
| 479 | *py=y;
|
---|
| 480 | *pobs=obs;
|
---|
| 481 | *pnobs=nobs;
|
---|
| 482 | }/*}}}*/
|
---|
[18301] | 483 | void Observations::InterpolationIDW(IssmPDouble *pprediction,IssmPDouble x_interp,IssmPDouble y_interp,IssmPDouble radius,int mindata,int maxdata,IssmPDouble power){/*{{{*/
|
---|
[14996] | 484 |
|
---|
| 485 | /*Intermediaries*/
|
---|
| 486 | int i,n_obs;
|
---|
| 487 | IssmPDouble prediction;
|
---|
| 488 | IssmPDouble numerator,denominator,h,weight;
|
---|
| 489 | IssmPDouble *x = NULL;
|
---|
| 490 | IssmPDouble *y = NULL;
|
---|
| 491 | IssmPDouble *obs = NULL;
|
---|
| 492 |
|
---|
| 493 | /*Some checks*/
|
---|
| 494 | _assert_(maxdata>0);
|
---|
| 495 | _assert_(pprediction);
|
---|
| 496 | _assert_(power>0);
|
---|
| 497 |
|
---|
| 498 | /*Get list of observations for current point*/
|
---|
| 499 | this->ObservationList(&x,&y,&obs,&n_obs,x_interp,y_interp,radius,maxdata);
|
---|
| 500 |
|
---|
| 501 | /*If we have less observations than mindata, return UNDEF*/
|
---|
| 502 | if(n_obs<mindata){
|
---|
[25836] | 503 | prediction = UNDEF;
|
---|
[14996] | 504 | }
|
---|
| 505 | else{
|
---|
| 506 | numerator = 0.;
|
---|
| 507 | denominator = 0.;
|
---|
| 508 | for(i=0;i<n_obs;i++){
|
---|
| 509 | h = sqrt( (x[i]-x_interp)*(x[i]-x_interp) + (y[i]-y_interp)*(y[i]-y_interp));
|
---|
| 510 | if (h<0.0000001){
|
---|
| 511 | numerator = obs[i];
|
---|
| 512 | denominator = 1.;
|
---|
| 513 | break;
|
---|
| 514 | }
|
---|
| 515 | weight = 1./pow(h,power);
|
---|
| 516 | numerator += weight*obs[i];
|
---|
| 517 | denominator += weight;
|
---|
| 518 | }
|
---|
[25836] | 519 | prediction = numerator/denominator;
|
---|
[14996] | 520 | }
|
---|
| 521 |
|
---|
| 522 | /*clean-up*/
|
---|
| 523 | *pprediction = prediction;
|
---|
| 524 | xDelete<IssmPDouble>(x);
|
---|
| 525 | xDelete<IssmPDouble>(y);
|
---|
| 526 | xDelete<IssmPDouble>(obs);
|
---|
| 527 | }/*}}}*/
|
---|
[18301] | 528 | void Observations::InterpolationKriging(IssmPDouble *pprediction,IssmPDouble *perror,IssmPDouble x_interp,IssmPDouble y_interp,IssmPDouble radius,int mindata,int maxdata,Variogram* variogram){/*{{{*/
|
---|
[14996] | 529 |
|
---|
| 530 | /*Intermediaries*/
|
---|
| 531 | int i,j,n_obs;
|
---|
| 532 | IssmPDouble prediction,error;
|
---|
[17806] | 533 | IssmPDouble *x = NULL;
|
---|
| 534 | IssmPDouble *y = NULL;
|
---|
| 535 | IssmPDouble *obs = NULL;
|
---|
| 536 | IssmPDouble *Lambda = NULL;
|
---|
[14996] | 537 |
|
---|
| 538 | /*Some checks*/
|
---|
| 539 | _assert_(mindata>0 && maxdata>0);
|
---|
| 540 | _assert_(pprediction && perror);
|
---|
| 541 |
|
---|
| 542 | /*Get list of observations for current point*/
|
---|
| 543 | this->ObservationList(&x,&y,&obs,&n_obs,x_interp,y_interp,radius,maxdata);
|
---|
| 544 |
|
---|
| 545 | /*If we have less observations than mindata, return UNDEF*/
|
---|
| 546 | if(n_obs<mindata){
|
---|
[25836] | 547 | *pprediction = -999.0;
|
---|
| 548 | *perror = -999.0;
|
---|
[14996] | 549 | return;
|
---|
| 550 | }
|
---|
| 551 |
|
---|
| 552 | /*Allocate intermediary matrix and vectors*/
|
---|
[17806] | 553 | IssmPDouble* A = xNew<IssmPDouble>((n_obs+1)*(n_obs+1));
|
---|
| 554 | IssmPDouble* B = xNew<IssmPDouble>(n_obs+1);
|
---|
[14996] | 555 |
|
---|
[25836] | 556 | IssmPDouble unbias = variogram->Covariance(0.,0.);
|
---|
[14996] | 557 | /*First: Create semivariogram matrix for observations*/
|
---|
| 558 | for(i=0;i<n_obs;i++){
|
---|
[20500] | 559 | //printf("%g %g ==> %g\n",x[i],y[i],sqrt(pow(x[i]-x_interp,2)+pow(y[i]-y_interp,2)));
|
---|
[14996] | 560 | for(j=0;j<=i;j++){
|
---|
[17806] | 561 | A[i*(n_obs+1)+j] = variogram->Covariance(x[i]-x[j],y[i]-y[j]);
|
---|
| 562 | A[j*(n_obs+1)+i] = A[i*(n_obs+1)+j];
|
---|
[14996] | 563 | }
|
---|
[17806] | 564 | A[i*(n_obs+1)+n_obs] = unbias;
|
---|
| 565 | //A[i*(n_obs+1)+n_obs] = 1.;
|
---|
[14996] | 566 | }
|
---|
[17806] | 567 | for(i=0;i<n_obs;i++) A[n_obs*(n_obs+1)+i]=unbias;
|
---|
| 568 | //for(i=0;i<n_obs;i++) A[n_obs*(n_obs+1)+i]=1.;
|
---|
| 569 | A[n_obs*(n_obs+1)+n_obs] = 0.;
|
---|
[14996] | 570 |
|
---|
| 571 | /*Get semivariogram vector associated to this location*/
|
---|
[17806] | 572 | for(i=0;i<n_obs;i++) B[i] = variogram->Covariance(x[i]-x_interp,y[i]-y_interp);
|
---|
| 573 | B[n_obs] = unbias;
|
---|
| 574 | //B[n_obs] = 1.;
|
---|
[14996] | 575 |
|
---|
| 576 | /*Solve the three linear systems*/
|
---|
| 577 | #if _HAVE_GSL_
|
---|
[17806] | 578 | DenseGslSolve(&Lambda,A,B,n_obs+1); // Gamma^-1 Z
|
---|
[14996] | 579 | #else
|
---|
| 580 | _error_("GSL is required");
|
---|
| 581 | #endif
|
---|
| 582 |
|
---|
[17806] | 583 | /*Compute predictor*/
|
---|
[14996] | 584 | prediction = 0.;
|
---|
[17806] | 585 | for(i=0;i<n_obs;i++) prediction += Lambda[i]*obs[i];
|
---|
[14996] | 586 |
|
---|
[17806] | 587 | /*Compute error (GSLIB p15 eq II.14)*/
|
---|
| 588 | error = variogram->Covariance(0.,0.)*(1. - Lambda[n_obs]);;
|
---|
| 589 | for(i=0;i<n_obs;i++) error += -Lambda[i]*B[i];
|
---|
| 590 |
|
---|
[14996] | 591 | /*clean-up*/
|
---|
| 592 | *pprediction = prediction;
|
---|
| 593 | *perror = error;
|
---|
| 594 | xDelete<IssmPDouble>(x);
|
---|
| 595 | xDelete<IssmPDouble>(y);
|
---|
| 596 | xDelete<IssmPDouble>(obs);
|
---|
[17806] | 597 | xDelete<IssmPDouble>(A);
|
---|
| 598 | xDelete<IssmPDouble>(B);
|
---|
| 599 | xDelete<IssmPDouble>(Lambda);
|
---|
[14996] | 600 | }/*}}}*/
|
---|
[18301] | 601 | void Observations::InterpolationNearestNeighbor(IssmPDouble *pprediction,IssmPDouble x_interp,IssmPDouble y_interp,IssmPDouble radius){/*{{{*/
|
---|
[14996] | 602 |
|
---|
| 603 | /*Intermediaries*/
|
---|
| 604 | IssmPDouble x,y,obs;
|
---|
| 605 |
|
---|
| 606 | /*Get clostest observation*/
|
---|
| 607 | this->ClosestObservation(&x,&y,&obs,x_interp,y_interp,radius);
|
---|
| 608 |
|
---|
| 609 | /*Assign output pointer*/
|
---|
| 610 | *pprediction = obs;
|
---|
| 611 | }/*}}}*/
|
---|
[18301] | 612 | void Observations::InterpolationV4(IssmPDouble *pprediction,IssmPDouble x_interp,IssmPDouble y_interp,IssmPDouble radius,int mindata,int maxdata){/*{{{*/
|
---|
[14996] | 613 | /* Reference: David T. Sandwell, Biharmonic spline interpolation of GEOS-3
|
---|
| 614 | * and SEASAT altimeter data, Geophysical Research Letters, 2, 139-142,
|
---|
| 615 | * 1987. Describes interpolation using value or gradient of value in any
|
---|
| 616 | * dimension.*/
|
---|
| 617 |
|
---|
| 618 | /*Intermediaries*/
|
---|
| 619 | int i,j,n_obs;
|
---|
| 620 | IssmPDouble prediction,h;
|
---|
| 621 | IssmPDouble *x = NULL;
|
---|
| 622 | IssmPDouble *y = NULL;
|
---|
| 623 | IssmPDouble *obs = NULL;
|
---|
| 624 | IssmPDouble *Green = NULL;
|
---|
| 625 | IssmPDouble *weights = NULL;
|
---|
| 626 | IssmPDouble *g = NULL;
|
---|
| 627 |
|
---|
| 628 | /*Some checks*/
|
---|
| 629 | _assert_(maxdata>0);
|
---|
| 630 | _assert_(pprediction);
|
---|
| 631 |
|
---|
| 632 | /*Get list of observations for current point*/
|
---|
| 633 | this->ObservationList(&x,&y,&obs,&n_obs,x_interp,y_interp,radius,maxdata);
|
---|
| 634 |
|
---|
| 635 | /*If we have less observations than mindata, return UNDEF*/
|
---|
| 636 | if(n_obs<mindata || n_obs<2){
|
---|
[25836] | 637 | prediction = UNDEF;
|
---|
[14996] | 638 | }
|
---|
| 639 | else{
|
---|
| 640 |
|
---|
| 641 | /*Allocate intermediary matrix and vectors*/
|
---|
| 642 | Green = xNew<IssmPDouble>(n_obs*n_obs);
|
---|
| 643 | g = xNew<IssmPDouble>(n_obs);
|
---|
| 644 |
|
---|
| 645 | /*First: distance vector*/
|
---|
| 646 | for(i=0;i<n_obs;i++){
|
---|
| 647 | h = sqrt( (x[i]-x_interp)*(x[i]-x_interp) + (y[i]-y_interp)*(y[i]-y_interp) );
|
---|
| 648 | if(h>0){
|
---|
| 649 | g[i] = h*h*(log(h)-1.);
|
---|
| 650 | }
|
---|
| 651 | else{
|
---|
| 652 | g[i] = 0.;
|
---|
| 653 | }
|
---|
| 654 | }
|
---|
| 655 |
|
---|
| 656 | /*Build Green function matrix*/
|
---|
| 657 | for(i=0;i<n_obs;i++){
|
---|
| 658 | for(j=0;j<=i;j++){
|
---|
| 659 | h = sqrt( (x[i]-x[j])*(x[i]-x[j]) + (y[i]-y[j])*(y[i]-y[j]) );
|
---|
| 660 | if(h>0){
|
---|
| 661 | Green[j*n_obs+i] = h*h*(log(h)-1.);
|
---|
| 662 | }
|
---|
| 663 | else{
|
---|
| 664 | Green[j*n_obs+i] = 0.;
|
---|
| 665 | }
|
---|
| 666 | Green[i*n_obs+j] = Green[j*n_obs+i];
|
---|
| 667 | }
|
---|
[17806] | 668 | /*Zero diagonal (should be done already, but just in case)*/
|
---|
| 669 | Green[i*n_obs+i] = 0.;
|
---|
[14996] | 670 | }
|
---|
| 671 |
|
---|
| 672 | /*Compute weights*/
|
---|
| 673 | #if _HAVE_GSL_
|
---|
| 674 | DenseGslSolve(&weights,Green,obs,n_obs); // Green^-1 obs
|
---|
| 675 | #else
|
---|
| 676 | _error_("GSL is required");
|
---|
| 677 | #endif
|
---|
| 678 |
|
---|
| 679 | /*Interpolate*/
|
---|
| 680 | prediction = 0;
|
---|
| 681 | for(i=0;i<n_obs;i++) prediction += weights[i]*g[i];
|
---|
| 682 |
|
---|
| 683 | }
|
---|
| 684 |
|
---|
| 685 | /*clean-up*/
|
---|
| 686 | *pprediction = prediction;
|
---|
| 687 | xDelete<IssmPDouble>(x);
|
---|
| 688 | xDelete<IssmPDouble>(y);
|
---|
| 689 | xDelete<IssmPDouble>(obs);
|
---|
| 690 | xDelete<IssmPDouble>(Green);
|
---|
| 691 | xDelete<IssmPDouble>(g);
|
---|
| 692 | xDelete<IssmPDouble>(weights);
|
---|
| 693 | }/*}}}*/
|
---|
[18301] | 694 | void Observations::QuadtreeColoring(IssmPDouble* A,IssmPDouble *x,IssmPDouble *y,int n){/*{{{*/
|
---|
[14996] | 695 |
|
---|
[19105] | 696 | if(this->treetype!=1) _error_("Tree type is not quadtree");
|
---|
[14996] | 697 | int xi,yi,level;
|
---|
| 698 |
|
---|
| 699 | for(int i=0;i<n;i++){
|
---|
| 700 | this->quadtree->IntergerCoordinates(&xi,&yi,x[i],y[i]);
|
---|
| 701 | this->quadtree->QuadtreeDepth(&level,xi,yi);
|
---|
| 702 | A[i]=(IssmPDouble)level;
|
---|
| 703 | }
|
---|
| 704 |
|
---|
| 705 | }/*}}}*/
|
---|
[18301] | 706 | void Observations::Variomap(IssmPDouble* gamma,IssmPDouble *x,int n){/*{{{*/
|
---|
[14996] | 707 |
|
---|
| 708 | /*Output and Intermediaries*/
|
---|
| 709 | int i,j,k;
|
---|
| 710 | IssmPDouble distance;
|
---|
| 711 | Observation *observation1 = NULL;
|
---|
| 712 | Observation *observation2 = NULL;
|
---|
| 713 |
|
---|
| 714 | IssmPDouble *counter = xNew<IssmPDouble>(n);
|
---|
| 715 | for(j=0;j<n;j++) counter[j] = 0.0;
|
---|
| 716 | for(j=0;j<n;j++) gamma[j] = 0.0;
|
---|
| 717 |
|
---|
[25836] | 718 | for(Object* & object : this->objects){
|
---|
| 719 | observation1=xDynamicCast<Observation*>(object);
|
---|
[14996] | 720 |
|
---|
[25836] | 721 | for(Object* & object : this->objects){
|
---|
| 722 | observation2=xDynamicCast<Observation*>(object);
|
---|
[14996] | 723 |
|
---|
[16137] | 724 | distance=sqrt(pow(observation1->x - observation2->x,2) + pow(observation1->y - observation2->y,2));
|
---|
[14996] | 725 | if(distance>x[n-1]) continue;
|
---|
| 726 |
|
---|
| 727 | int index = int(distance/(x[1]-x[0]));
|
---|
| 728 | if(index>n-1) index = n-1;
|
---|
| 729 | if(index<0) index = 0;
|
---|
| 730 |
|
---|
[16137] | 731 | gamma[index] += 1./2.*pow(observation1->value - observation2->value,2);
|
---|
[14996] | 732 | counter[index] += 1.;
|
---|
| 733 | }
|
---|
| 734 | }
|
---|
| 735 |
|
---|
| 736 | /*Normalize semivariogram*/
|
---|
| 737 | gamma[0]=0.;
|
---|
| 738 | for(k=0;k<n;k++){
|
---|
| 739 | if(counter[k]) gamma[k] = gamma[k]/counter[k];
|
---|
| 740 | }
|
---|
| 741 |
|
---|
| 742 | /*Assign output pointer*/
|
---|
| 743 | xDelete<IssmPDouble>(counter);
|
---|
| 744 | }/*}}}*/
|
---|