pktools  2.6.4
Processing Kernel for geospatial data
pkkalman.cc
1 /**********************************************************************
2 pkkalman.cc: produce kalman filtered raster time series
3 Copyright (C) 2008-2014 Pieter Kempeneers
4 
5 This file is part of pktools
6 
7 pktools is free software: you can redistribute it and/or modify
8 it under the terms of the GNU General Public License as published by
9 the Free Software Foundation, either version 3 of the License, or
10 (at your option) any later version.
11 
12 pktools is distributed in the hope that it will be useful,
13 but WITHOUT ANY WARRANTY; without even the implied warranty of
14 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
15 GNU General Public License for more details.
16 
17 You should have received a copy of the GNU General Public License
18 along with pktools. If not, see <http://www.gnu.org/licenses/>.
19 ***********************************************************************/
20 #include <sstream>
21 #include <vector>
22 #include <algorithm>
23 #include "base/Optionpk.h"
24 #include "base/Vector2d.h"
25 #include "imageclasses/ImgReaderGdal.h"
26 #include "imageclasses/ImgWriterGdal.h"
27 #include "algorithms/StatFactory.h"
28 #include "algorithms/ImgRegression.h"
29 
30 /******************************************************************************/
76 using namespace std;
77 /*------------------
78  Main procedure
79  ----------------*/
80 int main(int argc,char **argv) {
81  Optionpk<string> direction_opt("dir","direction","direction to run model (forward|backward|smooth)","forward");
82  Optionpk<string> model_opt("mod","model","model input datasets, e.g., MODIS (use: -mod model1 -mod model2 etc.)");
83  Optionpk<string> observation_opt("obs","observation","observation input datasets, e.g., landsat (use: -obs obs1 -obs obs2 etc.)");
84  Optionpk<int> tmodel_opt("tmod","tmodel","time sequence of model input. Sequence must have exact same length as model input. Leave empty to have default sequence 0,1,2,etc.");
85  Optionpk<int> tobservation_opt("tobs","tobservation","time sequence of observation input. Sequence must have exact same length as observation input)");
86  Optionpk<string> projection_opt("a_srs", "a_srs", "Override the projection for the output file (leave blank to copy from input file, use epsg:3035 to use European projection and force to European grid");
87  Optionpk<string> outputfw_opt("ofw", "outputfw", "Output raster dataset for forward model");
88  Optionpk<string> outputbw_opt("obw", "outputbw", "Output raster dataset for backward model");
89  Optionpk<string> outputfb_opt("ofb", "outputfb", "Output raster dataset for smooth model");
90  Optionpk<double> modnodata_opt("modnodata", "modnodata", "invalid value for model input", 0);
91  Optionpk<double> obsnodata_opt("obsnodata", "obsnodata", "invalid value for observation input", 0);
92  Optionpk<double> modoffset_opt("modoffset", "modoffset", "offset used to read model input dataset (value=offset+scale*readValue)");
93  Optionpk<double> obsoffset_opt("obsoffset", "obsoffset", "offset used to read observation input dataset (value=offset+scale*readValue)");
94  Optionpk<double> modscale_opt("modscale", "modscale", "scale used to read model input dataset (value=offset+scale*readValue)");
95  Optionpk<double> obsscale_opt("obsscale", "obsscale", "scale used to read observation input dataset (value=offset+scale*readValue)");
96  Optionpk<double> eps_opt("eps", "eps", "epsilon for non zero division", 0.00001);
97  Optionpk<double> uncertModel_opt("um", "uncertmodel", "Multiply this value with std dev of first model image to obtain uncertainty of model",2);
98  Optionpk<double> uncertObs_opt("uo", "uncertobs", "Uncertainty of valid observations",0);
99  Optionpk<double> weight_opt("w", "weight", "Penalize outliers in measurement via weights. Use first weight to penalize small measurements wrt model and second weight to penalize large measurements wrt model");
100  Optionpk<double> deltaObs_opt("dobs", "deltaobs", "Lower and upper thresholds for relative pixel differences (in percentage): (observation-model)/model. For instance to force the observation within a +/- 10 % interval, use: -dobs -10 -dobs 10 (equivalent to -dobs 10). Leave empty to always update on observation");
101  Optionpk<double> uncertNodata_opt("unodata", "uncertnodata", "Uncertainty in case of no-data values in observation", 10000);
102  Optionpk<double> regTime_opt("rt", "regtime", "Weight for regression in time series", 1.0);
103  Optionpk<double> regSensor_opt("rs", "regsensor", "Weight for regression model - measurement (model - observation).");
104  Optionpk<int> down_opt("down", "down", "Downsampling factor for reading model data to calculate regression");
105  Optionpk<float> threshold_opt("th", "threshold", "threshold for selecting samples (randomly). Provide probability in percentage (>0) or absolute (<0).", 0);
106  Optionpk<int> minreg_opt("minreg", "minreg", "Minimum number of pixels to take into account for regression", 5, 2);
107  // Optionpk<bool> regObs_opt("regobs", "regobs", "Perform regression between modeled and observed value",false);
108  // Optionpk<double> checkDiff_opt("diff", "diff", "Flag observation as invalid if difference with model is above uncertainty",false);
109  Optionpk<unsigned short> window_opt("win", "window", "window size for calculating regression (use 0 for global)", 0);
110  // Optionpk<string> mask_opt("m", "mask", "Use the first band of the specified file as a validity mask. Nodata values can be set with the option msknodata.");
111  // Optionpk<double> msknodata_opt("msknodata", "msknodata", "Mask value(s) not to consider for filtering. First value will be set in output image.", 0);
112 
113  Optionpk<string> oformat_opt("of", "oformat", "Output image format (see also gdal_translate). Empty string: inherit from input image","GTiff",2);
114  Optionpk<string> option_opt("co", "co", "Creation option for output file. Multiple options can be specified.");
115  Optionpk<short> verbose_opt("v", "verbose", "verbose mode when positive", 0);
116 
117  bool doProcess;//stop process when program was invoked with help option (-h --help)
118  try{
119  doProcess=direction_opt.retrieveOption(argc,argv);
120  model_opt.retrieveOption(argc,argv);
121  observation_opt.retrieveOption(argc,argv);
122  tmodel_opt.retrieveOption(argc,argv);
123  tobservation_opt.retrieveOption(argc,argv);
124  projection_opt.retrieveOption(argc,argv);
125  outputfw_opt.retrieveOption(argc,argv);
126  outputbw_opt.retrieveOption(argc,argv);
127  outputfb_opt.retrieveOption(argc,argv);
128  modnodata_opt.retrieveOption(argc,argv);
129  obsnodata_opt.retrieveOption(argc,argv);
130  modoffset_opt.retrieveOption(argc,argv);
131  modscale_opt.retrieveOption(argc,argv);
132  obsoffset_opt.retrieveOption(argc,argv);
133  obsscale_opt.retrieveOption(argc,argv);
134  eps_opt.retrieveOption(argc,argv);
135  uncertModel_opt.retrieveOption(argc,argv);
136  uncertObs_opt.retrieveOption(argc,argv);
137  weight_opt.retrieveOption(argc,argv);
138  deltaObs_opt.retrieveOption(argc,argv);
139  uncertNodata_opt.retrieveOption(argc,argv);
140  regTime_opt.retrieveOption(argc,argv);
141  regSensor_opt.retrieveOption(argc,argv);
142  down_opt.retrieveOption(argc,argv);
143  threshold_opt.retrieveOption(argc,argv);
144  minreg_opt.retrieveOption(argc,argv);
145  // regObs_opt.retrieveOption(argc,argv);
146  // checkDiff_opt.retrieveOption(argc,argv);
147  window_opt.retrieveOption(argc,argv);
148  // mask_opt.retrieveOption(argc,argv);
149  // msknodata_opt.retrieveOption(argc,argv);
150  oformat_opt.retrieveOption(argc,argv);
151  option_opt.retrieveOption(argc,argv);
152  verbose_opt.retrieveOption(argc,argv);
153  }
154  catch(string predefinedString){
155  std::cout << predefinedString << std::endl;
156  exit(0);
157  }
158  if(!doProcess){
159  std::cerr << "short option -h shows basic options only, use long option --help to show all options" << std::endl;
160  exit(0);//help was invoked, stop processing
161  }
162 
163  if(deltaObs_opt.size()==1){
164  if(deltaObs_opt[0]<=0)
165  deltaObs_opt.push_back(-deltaObs_opt[0]);
166  else
167  deltaObs_opt.insert(deltaObs_opt.begin(),-deltaObs_opt[0]);
168  }
169  if(weight_opt.size()==1){
170  weight_opt.push_back(weight_opt[0]);
171  }
172 
173  if(down_opt.empty()){
174  std::cerr << "short option -h shows basic options only, use long option --help to show all options" << std::endl;
175  exit(0);//help was invoked, stop processing
176  }
177  try{
178  ostringstream errorStream;
179  if(model_opt.size()<2){
180  errorStream << "Error: no model dataset selected, use option -mod" << endl;
181  throw(errorStream.str());
182  }
183  if(observation_opt.size()<1){
184  errorStream << "Error: no observation dataset selected, use option -obs" << endl;
185  throw(errorStream.str());
186  }
187  if(direction_opt[0]=="smooth"){
188  if(outputfw_opt.empty()){
189  errorStream << "Error: output forward datasets is not provided, use option -ofw" << endl;
190  throw(errorStream.str());
191  }
192  if(outputbw_opt.empty()){
193  errorStream << "Error: output backward datasets is not provided, use option -obw" << endl;
194  throw(errorStream.str());
195  }
196  if(outputfb_opt.empty()){
197  errorStream << "Error: output smooth datasets is not provided, use option -ofb" << endl;
198  throw(errorStream.str());
199  }
200  }
201  else{
202  if(direction_opt[0]=="forward"&&outputfw_opt.empty()){
203  errorStream << "Error: output forward datasets is not provided, use option -ofw" << endl;
204  throw(errorStream.str());
205  }
206  else if(direction_opt[0]=="backward"&&outputbw_opt.empty()){
207  errorStream << "Error: output backward datasets is not provided, use option -obw" << endl;
208  throw(errorStream.str());
209  }
210 
211  if(model_opt.size()<observation_opt.size()){
212  errorStream << "Error: sequence of models should be larger than observations" << endl;
213  throw(errorStream.str());
214  }
215  if(tmodel_opt.size()!=model_opt.size()){
216  if(tmodel_opt.empty())
217  cout << "Warning: time sequence is not provided, self generating time sequence from 0 to " << model_opt.size() << endl;
218  else
219  cout << "Warning: time sequence provided (" << tmodel_opt.size() << ") does not match number of model raster datasets (" << model_opt.size() << ")" << endl;
220  tmodel_opt.clear();
221  for(int tindex=0;tindex<model_opt.size();++tindex)
222  tmodel_opt.push_back(tindex);
223  }
224  if(tobservation_opt.size()!=observation_opt.size()){
225  errorStream << "Error: time sequence for observation must match size of observation dataset" << endl;
226  throw(errorStream.str());
227  }
228  }
229  }
230  catch(string errorString){
231  std::cout << errorString << std::endl;
232  exit(1);
233  }
234 
236  stat.setNoDataValues(modnodata_opt);
238  // vector<ImgReaderGdal> imgReaderModel(model_opt.size());
239  // vector<ImgReaderGdal> imgReaderObs(observation_opt.size());
240  ImgReaderGdal imgReaderModel1;
241  ImgReaderGdal imgReaderModel2;
242  ImgReaderGdal imgReaderObs;
243  ImgWriterGdal imgWriterEst;
244 
245  imgReaderObs.open(observation_opt[0]);
246 
247  int ncol=imgReaderObs.nrOfCol();
248  int nrow=imgReaderObs.nrOfRow();
249  if(projection_opt.empty())
250  projection_opt.push_back(imgReaderObs.getProjection());
251  double geotransform[6];
252  imgReaderObs.getGeoTransform(geotransform);
253 
254  string imageType=imgReaderObs.getImageType();
255  if(oformat_opt.size())//default
256  imageType=oformat_opt[0];
257  if(option_opt.findSubstring("INTERLEAVE=")==option_opt.end()){
258  string theInterleave="INTERLEAVE=";
259  theInterleave+=imgReaderObs.getInterleave();
260  option_opt.push_back(theInterleave);
261  }
262 
263  if(down_opt.empty()){
264  imgReaderModel1.open(model_opt[0]);
265  double resModel=imgReaderModel1.getDeltaX();
266  double resObs=imgReaderObs.getDeltaX();
267  int down=static_cast<int>(ceil(resModel/resObs));
268  if(!(down%2))
269  down+=1;
270  down_opt.push_back(down);
271  imgReaderModel1.close();
272  }
273  imgReaderObs.close();
274 
275  if(regSensor_opt.empty())
276  regSensor_opt.push_back(1.0/down_opt[0]);
277  //hiero
278  // ImgReaderGdal maskReader;
279  // double colMask=0;
280  // double rowMask=0;
281 
282  // if(mask_opt.size()){
283  // try{
284  // if(verbose_opt[0]>=1)
285  // std::cout << "opening mask image file " << mask_opt[0] << std::endl;
286  // maskReader.open(mask_opt[0]);
287  // maskReader.setNoData(msknodata_opt);
288  // }
289  // catch(string error){
290  // cerr << error << std::endl;
291  // exit(2);
292  // }
293  // catch(...){
294  // cerr << "error catched" << std::endl;
295  // exit(1);
296  // }
297  // }
298 
299  int obsindex=0;
300 
301  const char* pszMessage;
302  void* pProgressArg=NULL;
303  GDALProgressFunc pfnProgress=GDALTermProgress;
304  double progress=0;
305 
306  imgreg.setDown(down_opt[0]);
307  imgreg.setThreshold(threshold_opt[0]);
308 
309  double c0modGlobal=0;//time regression coefficient c0 (offset) calculated on entire image
310  double c1modGlobal=1;//time regression coefficient c1 (scale) calculated on entire image
311  double c0mod=0;//time regression coefficient c0 (offset) calculated on local window
312  double c1mod=1;//time regression coefficient c1 (scale) calculated on local window
313 
314  double c0obs=0;//offset
315  double c1obs=1;//scale
316  double errObs=uncertNodata_opt[0];//start with high initial value in case we do not have first observation at time 0
317 
318  vector<int> relobsindex;
319  // cout << tmodel_opt << endl;
320  // cout << tobservation_opt << endl;
321 
322  for(int tindex=0;tindex<tobservation_opt.size();++tindex){
323  vector<int>::iterator modit;
324  modit=upper_bound(tmodel_opt.begin(),tmodel_opt.end(),tobservation_opt[tindex]);
325  int relpos=modit-tmodel_opt.begin()-1;
326  assert(relpos>=0);//todo: for now, we assume model is available at time before first measurement
327  relobsindex.push_back(relpos);
328  if(verbose_opt[0])
329  cout << "observation " << tindex << ": " << "relative position in model time series is " << relpos << ", date of observation is (tobservation_opt[tindex]): " << tobservation_opt[tindex] << ", relobsindex.back(): " << relobsindex.back() << ", filename observation: " << observation_opt[tindex] << ", filename of corresponding model: " << model_opt[relpos] << endl;
330  // if(verbose_opt[0])
331  // cout << "tobservation_opt[tindex] " << tobservation_opt[tindex] << " " << relobsindex.back() << endl;
332  }
333 
334  int ndigit=log(1.0*tmodel_opt.back())/log(10.0)+1;
335 
336  double geox=0;
337  double geoy=0;
338 
339  if(find(direction_opt.begin(),direction_opt.end(),"forward")!=direction_opt.end()){
341  cout << "Running forward model" << endl;
342  obsindex=0;
343  //initialization
344  string output;
345  if(outputfw_opt.size()==model_opt.size()){
346  output=outputfw_opt[0];
347  }
348  else{
349  ostringstream outputstream;
350  outputstream << outputfw_opt[0] << "_";
351  outputstream << setfill('0') << setw(ndigit) << tmodel_opt[0];
352  outputstream << ".tif";
353  //test
354  // outputstream << outputfw_opt[0] << "_" << tmodel_opt[0] << ".tif";
355  output=outputstream.str();
356  }
357  if(verbose_opt[0])
358  cout << "Opening image " << output << " for writing " << endl;
359 
360  imgWriterEst.open(output,ncol,nrow,2,GDT_Float32,imageType,option_opt);
361  imgWriterEst.setProjectionProj4(projection_opt[0]);
362  imgWriterEst.setGeoTransform(geotransform);
363  imgWriterEst.GDALSetNoDataValue(obsnodata_opt[0]);
364 
365  if(verbose_opt[0]){
366  cout << "processing time " << tmodel_opt[0] << endl;
367  if(obsindex<relobsindex.size())
368  cout << "next observation " << tmodel_opt[relobsindex[obsindex]] << endl;
369  else
370  cout << "There is no next observation" << endl;
371  }
372 
373  try{
374  imgReaderModel1.open(model_opt[0]);
375  imgReaderModel1.setNoData(modnodata_opt);
376  if(modoffset_opt.size())
377  imgReaderModel1.setOffset(modoffset_opt[0]);
378  if(modscale_opt.size())
379  imgReaderModel1.setScale(modscale_opt[0]);
380  }
381  catch(string errorString){
382  cerr << errorString << endl;
383  }
384  catch(...){
385  cerr << "Error opening file " << model_opt[0] << endl;
386  }
387 
388  //calculate standard deviation of image to serve as model uncertainty
389  GDALRasterBand* rasterBand;
390  rasterBand=imgReaderModel1.getRasterBand(0);
391  double minValue, maxValue, meanValue, stdDev;
392  void* pProgressData;
393  rasterBand->ComputeStatistics(0,&minValue,&maxValue,&meanValue,&stdDev,pfnProgress,pProgressData);
394  double modRow=0;
395  double modCol=0;
396  if(relobsindex[0]>0){//initialize output_opt[0] as model[0]
397  //write first model as output
398  if(verbose_opt[0])
399  cout << "write first model as output" << endl;
400  for(int irow=0;irow<nrow;++irow){
401  vector<double> estReadBuffer;
402  vector<double> estWriteBuffer(ncol);
403  vector<double> uncertWriteBuffer(ncol);
404  // vector<double> lineMask;
405  imgWriterEst.image2geo(0,irow,geox,geoy);
406  imgReaderModel1.geo2image(geox,geoy,modCol,modRow);
407  if(modRow<0||modRow>=imgReaderModel1.nrOfRow()){
408  cerr << "Error: geo coordinates (" << geox << "," << geoy << ") not covered in model image " << imgReaderModel1.getFileName() << endl;
409  assert(modRow>=0&&modRow<imgReaderModel1.nrOfRow());
410  }
411  try{
412  imgReaderModel1.readData(estReadBuffer,GDT_Float64,modRow);
413  //simple nearest neighbor
414  //stat.nearUp(estReadBuffer,estWriteBuffer);
415 
416  // double oldRowMask=-1;//keep track of row mask to optimize number of line readings
417  for(int icol=0;icol<ncol;++icol){
418  imgWriterEst.image2geo(icol,irow,geox,geoy);
419  imgReaderModel1.geo2image(geox,geoy,modCol,modRow);
420  double modValue=estReadBuffer[modCol];
421  if(imgReaderModel1.isNoData(modValue)){
422  estWriteBuffer[icol]=obsnodata_opt[0];
423  uncertWriteBuffer[icol]=uncertNodata_opt[0];
424  }
425  else{
426  //todo: should take into account regression model-obs...
427  estWriteBuffer[icol]=modValue;
428  uncertWriteBuffer[icol]=uncertModel_opt[0]*stdDev;
429  }
430  }
431  imgWriterEst.writeData(estWriteBuffer,GDT_Float64,irow,0);
432  imgWriterEst.writeData(uncertWriteBuffer,GDT_Float64,irow,1);
433  }
434  catch(string errorString){
435  cerr << errorString << endl;
436  }
437  catch(...){
438  cerr << "Error writing file " << imgWriterEst.getFileName() << endl;
439  }
440  }
441  }
442  else{//we have a measurement
443  if(verbose_opt[0])
444  cout << "we have a measurement at initial time" << endl;
445  imgReaderObs.open(observation_opt[0]);
446  imgReaderObs.getGeoTransform(geotransform);
447  imgReaderObs.setNoData(obsnodata_opt);
448  if(obsoffset_opt.size())
449  imgReaderObs.setOffset(obsoffset_opt[0]);
450  if(obsscale_opt.size())
451  imgReaderObs.setScale(obsscale_opt[0]);
452 
453  if(regSensor_opt[0]>0){
454  errObs=regSensor_opt[0]*imgreg.getRMSE(imgReaderModel1,imgReaderObs,c0obs,c1obs,0,0,verbose_opt[0]);
455  if(errObs<0){
456  c0obs=0;
457  c1obs=1;
458  }
459  }
460  else{
461  c0obs=0;
462  c1obs=1;
463  errObs=0;
464  }
465  if(verbose_opt[0])
466  cout << "c0obs, c1obs: " << c0obs << ", " << c1obs << endl;
467 
468  for(int irow=0;irow<nrow;++irow){
469  vector<double> estReadBuffer;
470  imgWriterEst.image2geo(0,irow,geox,geoy);
471  imgReaderModel1.geo2image(geox,geoy,modCol,modRow);
472  assert(modRow>=0&&modRow<imgReaderModel1.nrOfRow());
473  imgReaderModel1.readData(estReadBuffer,GDT_Float64,modRow);
474  vector<double> obsLineBuffer;
475  vector<double> estWriteBuffer(ncol);
476  vector<double> uncertWriteBuffer(ncol);
477  vector<double> uncertObsLineBuffer;
478  // vector<double> lineMask;
479  // imgReaderObs.readData(estWriteBuffer,GDT_Float64,irow,0);
480  imgReaderObs.readData(obsLineBuffer,GDT_Float64,irow,0);
481 
482  if(imgReaderObs.nrOfBand()>1)
483  imgReaderObs.readData(uncertObsLineBuffer,GDT_Float64,irow,1);
484  // double oldRowMask=-1;//keep track of row mask to optimize number of line readings
485  for(int icol=0;icol<ncol;++icol){
486  imgWriterEst.image2geo(icol,irow,geox,geoy);
487  imgReaderModel1.geo2image(geox,geoy,modCol,modRow);
488  assert(modRow>=0&&modRow<imgReaderModel1.nrOfRow());
489  double modValue=estReadBuffer[modCol];
490  if(imgReaderModel1.isNoData(modValue)){//model is nodata: retain observation
491  estWriteBuffer[icol]=obsLineBuffer[icol];
492  if(imgReaderObs.isNoData(obsLineBuffer[icol])){
493  estWriteBuffer[icol]=obsnodata_opt[0];
494  uncertWriteBuffer[icol]=uncertNodata_opt[0];
495  }
496  else if(uncertObsLineBuffer.size()>icol)
497  uncertWriteBuffer[icol]=uncertObsLineBuffer[icol];
498  else
499  uncertWriteBuffer[icol]=uncertObs_opt[0];
500  }
501  else{//model is valid: calculate estimate from model
502  double errMod=uncertModel_opt[0]*stdDev;
503  errMod*=regTime_opt[0];
504  // double certNorm=(errMod*errMod+errObs*errObs);
505  // double certMod=errObs*errObs/certNorm;
506  // double certObs=errMod*errMod/certNorm;
507  // double regTime=0;
508  // double regSensor=(c0obs+c1obs*modValue)*certMod;
509  // estWriteBuffer[icol]=regTime+regSensor;
510  estWriteBuffer[icol]=modValue;
511  double totalUncertainty=errMod;
512  // if(errMod<eps_opt[0])
513  // totalUncertainty=errObs;
514  // else if(errObs<eps_opt[0])
515  // totalUncertainty=errMod;
516  // else{
517  // totalUncertainty=1.0/errMod/errMod+1/errObs/errObs;
518  // totalUncertainty=sqrt(1.0/totalUncertainty);
519  // }
520  uncertWriteBuffer[icol]=totalUncertainty;//in case observation is not valid
521  }
522  // uncertWriteBuffer[icol]+=uncertReadBuffer[icol];
523  //measurement update
524  if(!imgReaderObs.isNoData(obsLineBuffer[icol])){
525  double kalmanGain=1;
526  double uncertObs=uncertObs_opt[0];
527  if(uncertObsLineBuffer.size()>icol)
528  uncertObs=uncertObsLineBuffer[icol];
529  else if(weight_opt.size()||deltaObs_opt.size()){
530  vector<double> obsWindowBuffer;//buffer for observation to calculate average corresponding to model pixel
531  int minCol=(icol>down_opt[0]/2) ? icol-down_opt[0]/2 : 0;
532  int maxCol=(icol+down_opt[0]/2<imgReaderObs.nrOfCol()) ? icol+down_opt[0]/2 : imgReaderObs.nrOfCol()-1;
533  int minRow=(irow>down_opt[0]/2) ? irow-down_opt[0]/2 : 0;
534  int maxRow=(irow+down_opt[0]/2<imgReaderObs.nrOfRow()) ? irow+down_opt[0]/2 : imgReaderObs.nrOfRow()-1;
535 
536  imgReaderObs.readDataBlock(obsWindowBuffer,GDT_Float64,minCol,maxCol,minRow,maxRow,0);
537  statfactory::StatFactory statobs;
538  statobs.setNoDataValues(obsnodata_opt);
539  double obsMeanValue=(statobs.mean(obsWindowBuffer)-c0obs)/c1obs;
540  double difference=obsMeanValue-modValue;
541  if(modValue){
542  double relativeDifference=difference/modValue;
543  if(deltaObs_opt.size()){
544  assert(deltaObs_opt.size()>1);
545  if(100*relativeDifference<deltaObs_opt[0])//lower bound
546  kalmanGain=0;
547  else if(100*relativeDifference>deltaObs_opt[1])//upper bound
548  kalmanGain=0;
549  }
550  else if(weight_opt.size()){
551  assert(weight_opt.size()>1);
552  if(obsMeanValue<modValue)
553  uncertObs=weight_opt[0]*relativeDifference;
554  else if(obsMeanValue>modValue)
555  uncertObs=weight_opt[1]*relativeDifference;
556  }
557  }
558  if(uncertObs<=0)
559  uncertObs=0;
560  if(verbose_opt[0]>1)
561  cout << "obsMeanValue:" << obsMeanValue << ", modValue: " << modValue << endl;
562  }
563  if(kalmanGain>0){
564  if((uncertWriteBuffer[icol]+uncertObs)>eps_opt[0])
565  kalmanGain=uncertWriteBuffer[icol]/(uncertWriteBuffer[icol]+uncertObs);
566  }
567  assert(kalmanGain<=1);
568  estWriteBuffer[icol]+=kalmanGain*(obsLineBuffer[icol]-estWriteBuffer[icol]);
569  uncertWriteBuffer[icol]*=(1-kalmanGain);
570  }
571  }
572  imgWriterEst.writeData(estWriteBuffer,GDT_Float64,irow,0);
573  imgWriterEst.writeData(uncertWriteBuffer,GDT_Float64,irow,1);
574  }
575  imgReaderObs.close();
576  ++obsindex;
577  }
578  imgReaderModel1.close();
579  imgWriterEst.close();
580 
581  for(int modindex=1;modindex<model_opt.size();++modindex){
582  if(verbose_opt[0]){
583  cout << "processing time " << tmodel_opt[modindex] << endl;
584  if(obsindex<relobsindex.size())
585  cout << "next observation " << tmodel_opt[relobsindex[obsindex]] << endl;
586  else
587  cout << "There is no next observation" << endl;
588  }
589  string output;
590  if(outputfw_opt.size()==model_opt.size())
591  output=outputfw_opt[modindex];
592  else{
593  ostringstream outputstream;
594  outputstream << outputfw_opt[0] << "_";
595  outputstream << setfill('0') << setw(ndigit) << tmodel_opt[modindex];
596  outputstream << ".tif";
597  // outputstream << outputfw_opt[0] << "_" << tmodel_opt[modindex] << ".tif";
598  output=outputstream.str();
599  }
600 
601  //two band output band0=estimation, band1=uncertainty
602  imgWriterEst.open(output,ncol,nrow,2,GDT_Float32,imageType,option_opt);
603  imgWriterEst.setProjectionProj4(projection_opt[0]);
604  imgWriterEst.setGeoTransform(geotransform);
605  imgWriterEst.GDALSetNoDataValue(obsnodata_opt[0]);
606 
607  //calculate regression between two subsequence model inputs
608  imgReaderModel1.open(model_opt[modindex-1]);
609  imgReaderModel1.setNoData(modnodata_opt);
610  if(modoffset_opt.size())
611  imgReaderModel1.setOffset(modoffset_opt[0]);
612  if(modscale_opt.size())
613  imgReaderModel1.setScale(modscale_opt[0]);
614  imgReaderModel2.open(model_opt[modindex]);
615  imgReaderModel2.setNoData(modnodata_opt);
616  if(modoffset_opt.size())
617  imgReaderModel2.setOffset(modoffset_opt[0]);
618  if(modscale_opt.size())
619  imgReaderModel2.setScale(modscale_opt[0]);
620  //calculate regression
621  //we could re-use the points from second image from last run, but
622  //to keep it general, we must redo it (overlap might have changed)
623 
624  pfnProgress(progress,pszMessage,pProgressArg);
625 
626  if(verbose_opt[0])
627  cout << "Calculating regression for " << imgReaderModel1.getFileName() << " " << imgReaderModel2.getFileName() << endl;
628 
629  double errMod=imgreg.getRMSE(imgReaderModel1,imgReaderModel2,c0modGlobal,c1modGlobal,0,0);
630  errMod*=regTime_opt[0];
631 
632  // double errMod=imgreg.getRMSE(imgReaderModel1,imgReaderModel2,c0modGlobal,c1modGlobal,verbose_opt[0]);
633  if(verbose_opt[0])
634  cout << "c0modGlobal, c1modGlobal: " << c0modGlobal << ", " << c1modGlobal << endl;
635 
636  bool update=false;
637  if(obsindex<relobsindex.size()){
638  update=(relobsindex[obsindex]==modindex);
639  }
640  if(update){
641  if(verbose_opt[0])
642  cout << "***update " << relobsindex[obsindex] << " = " << modindex << " " << observation_opt[obsindex] << " ***" << endl;
643 
644  imgReaderObs.open(observation_opt[obsindex]);
645  imgReaderObs.getGeoTransform(geotransform);
646  imgReaderObs.setNoData(obsnodata_opt);
647  if(obsoffset_opt.size())
648  imgReaderObs.setOffset(obsoffset_opt[0]);
649  if(obsscale_opt.size())
650  imgReaderObs.setScale(obsscale_opt[0]);
651  //calculate regression between model and observation
652  if(verbose_opt[0])
653  cout << "Calculating regression for " << imgReaderModel2.getFileName() << " " << imgReaderObs.getFileName() << endl;
654  if(regSensor_opt[0]>0){
655  errObs=regSensor_opt[0]*imgreg.getRMSE(imgReaderModel2,imgReaderObs,c0obs,c1obs,0,0,verbose_opt[0]);
656  if(errObs<0){
657  c0obs=0;
658  c1obs=1;
659  }
660  }
661  else{
662  c0obs=0;
663  c1obs=1;
664  errObs=0;
665  }
666  if(verbose_opt[0])
667  cout << "c0obs, c1obs, errObs: " << c0obs << ", " << c1obs << ", " << errObs << endl;
668  }
669  //prediction (also to fill cloudy pixels in update mode)
670  string input;
671  if(outputfw_opt.size()==model_opt.size())
672  input=outputfw_opt[modindex-1];
673  else{
674  ostringstream outputstream;
675  outputstream << outputfw_opt[0] << "_";
676  outputstream << setfill('0') << setw(ndigit) << tmodel_opt[modindex-1];
677  outputstream << ".tif";
678  // outputstream << outputfw_opt[0] << "_" << tmodel_opt[modindex-1] << ".tif";
679  input=outputstream.str();
680  }
681  if(verbose_opt[0])
682  cout << "opening " << input << endl;
683  ImgReaderGdal imgReaderEst(input);
684  imgReaderEst.setNoData(obsnodata_opt);
685  if(obsoffset_opt.size())
686  imgReaderEst.setOffset(obsoffset_opt[0]);
687  if(obsscale_opt.size())
688  imgReaderEst.setScale(obsscale_opt[0]);
689 
690  vector< vector<double> > obsLineVector(down_opt[0]);
691  vector<double> obsLineBuffer;
692  vector<double> obsWindowBuffer;//buffer for observation to calculate average corresponding to model pixel
693  vector<double> model1LineBuffer;
694  vector<double> model2LineBuffer;
695  vector<double> model1buffer;//buffer for model 1 to calculate time regression based on window
696  vector<double> model2buffer;//buffer for model 2 to calculate time regression based on window
697  vector<double> uncertObsLineBuffer;
698  vector<double> estReadBuffer;
699  // vector<double> estWindowBuffer;//buffer for estimate to calculate average corresponding to model pixel
700  vector<double> uncertReadBuffer;
701  vector<double> estWriteBuffer(ncol);
702  vector<double> uncertWriteBuffer(ncol);
703  // vector<double> lineMask;
704 
705  //initialize obsLineVector if update
706  if(update){
707  if(verbose_opt[0])
708  cout << "initialize obsLineVector" << endl;
709  assert(down_opt[0]%2);//window size must be odd
710  for(int iline=-down_opt[0]/2;iline<down_opt[0]/2+1;++iline){
711  if(iline<0)//replicate line 0
712  imgReaderObs.readData(obsLineVector[iline+down_opt[0]/2],GDT_Float64,0,0);
713  else
714  imgReaderObs.readData(obsLineVector[iline+down_opt[0]/2],GDT_Float64,iline,0);
715  }
716  }
717  for(int irow=0;irow<imgWriterEst.nrOfRow();++irow){
718  //do not read from imgReaderObs, because we read entire window for each pixel...
719  imgReaderEst.readData(estReadBuffer,GDT_Float64,irow,0);
720  imgReaderEst.readData(uncertReadBuffer,GDT_Float64,irow,1);
721  //read model2 in case current estimate is nodata
722  imgReaderEst.image2geo(0,irow,geox,geoy);
723  imgReaderModel2.geo2image(geox,geoy,modCol,modRow);
724  assert(modRow>=0&&modRow<imgReaderModel2.nrOfRow());
725  imgReaderModel2.readData(model2LineBuffer,GDT_Float64,modRow,0);
726 
727  imgReaderModel1.geo2image(geox,geoy,modCol,modRow);
728  assert(modRow>=0&&modRow<imgReaderModel1.nrOfRow());
729  imgReaderModel1.readData(model1LineBuffer,GDT_Float64,modRow,0);
730 
731  if(update){
732  int maxRow=(irow+down_opt[0]/2<imgReaderEst.nrOfRow()) ? irow+down_opt[0]/2 : imgReaderEst.nrOfRow()-1;
733  obsLineVector.erase(obsLineVector.begin());
734  imgReaderObs.readData(obsLineBuffer,GDT_Float64,maxRow,0);
735  obsLineVector.push_back(obsLineBuffer);
736  obsLineBuffer=obsLineVector[down_opt[0]/2];
737  // imgReaderObs.readData(obsLineBuffer,GDT_Float64,irow,0);
738  if(imgReaderObs.nrOfBand()>1)
739  imgReaderObs.readData(uncertObsLineBuffer,GDT_Float64,irow,1);
740  }
741  // double oldRowMask=-1;//keep track of row mask to optimize number of line readings
742  for(int icol=0;icol<imgWriterEst.nrOfCol();++icol){
743  imgReaderEst.image2geo(icol,irow,geox,geoy);
744  int minCol=(icol>down_opt[0]/2) ? icol-down_opt[0]/2 : 0;
745  int maxCol=(icol+down_opt[0]/2<imgReaderEst.nrOfCol()) ? icol+down_opt[0]/2 : imgReaderEst.nrOfCol()-1;
746  int minRow=(irow>down_opt[0]/2) ? irow-down_opt[0]/2 : 0;
747  int maxRow=(irow+down_opt[0]/2<imgReaderEst.nrOfRow()) ? irow+down_opt[0]/2 : imgReaderEst.nrOfRow()-1;
748  if(update){
749  obsWindowBuffer.clear();
750  for(int iline=0;iline<obsLineVector.size();++iline){
751  for(int isample=minCol;isample<=maxCol;++isample){
752  assert(isample<obsLineVector[iline].size());
753  obsWindowBuffer.push_back(obsLineVector[iline][isample]);
754  }
755  }
756  // imgReaderObs.readDataBlock(obsWindowBuffer,GDT_Float64,minCol,maxCol,minRow,maxRow,0);
757  }
758  double estValue=estReadBuffer[icol];
759  double estMeanValue=0;//stat.mean(estWindowBuffer);
760  double nvalid=0;
761  //time update
762  imgReaderModel2.geo2image(geox,geoy,modCol,modRow);
763  assert(modRow>=0&&modRow<imgReaderModel2.nrOfRow());
764  double modValue=model2LineBuffer[modCol];
765  bool estNodata=imgReaderEst.isNoData(estValue);//validity of current estimate
766  if(estNodata){
767  //we have not found any valid data yet, better here to take the current model value if valid
768  if(imgReaderModel2.isNoData(modValue)){//if both estimate and model are no-data, set obs to nodata
769  estWriteBuffer[icol]=obsnodata_opt[0];
770  uncertWriteBuffer[icol]=uncertNodata_opt[0];
771  }
772  else{
773  estWriteBuffer[icol]=modValue;
774  uncertWriteBuffer[icol]=uncertModel_opt[0]*stdDev;
775  }
776  }
777  else{
778  if(window_opt[0]>0){
779  try{
780  // imgReaderModel2.geo2image(geox,geoy,modCol,modRow);//did that already
781  minCol=(modCol>window_opt[0]/2) ? modCol-window_opt[0]/2 : 0;
782  maxCol=(modCol+window_opt[0]/2<imgReaderModel2.nrOfCol()) ? modCol+window_opt[0]/2 : imgReaderModel2.nrOfCol()-1;
783  minRow=(modRow>window_opt[0]/2) ? modRow-window_opt[0]/2 : 0;
784  maxRow=(modRow+window_opt[0]/2<imgReaderModel2.nrOfRow()) ? modRow+window_opt[0]/2 : imgReaderModel2.nrOfRow()-1;
785  imgReaderModel2.readDataBlock(model2buffer,GDT_Float64,minCol,maxCol,minRow,maxRow,0);
786 
787  imgReaderModel1.geo2image(geox,geoy,modCol,modRow);
788  assert(modRow>=0&&modRow<imgReaderModel1.nrOfRow());
789  minCol=(modCol>window_opt[0]/2) ? modCol-window_opt[0]/2 : 0;
790  maxCol=(modCol+window_opt[0]/2<imgReaderModel1.nrOfCol()) ? modCol+window_opt[0]/2 : imgReaderModel1.nrOfCol()-1;
791  minRow=(modRow>window_opt[0]/2) ? modRow-window_opt[0]/2 : 0;
792  maxRow=(modRow+window_opt[0]/2<imgReaderModel1.nrOfRow()) ? modRow+window_opt[0]/2 : imgReaderModel1.nrOfRow()-1;
793  imgReaderModel1.readDataBlock(model1buffer,GDT_Float64,minCol,maxCol,minRow,maxRow,0);
794  // imgReaderEst.image2geo(icol,irow,geox,geoy);
795  }
796  catch(string errorString){
797  cerr << "Error reading data block for " << minCol << "-" << maxCol << ", " << minRow << "-" << maxRow << endl;
798  }
799  //erase no-data from buffer
800  vector<double>::iterator it1=model1buffer.begin();
801  vector<double>::iterator it2=model2buffer.begin();
802  while(it1!=model1buffer.end()&&it2!=model2buffer.end()){
803  //erase nodata
804  bool modNodata=false;
805  modNodata=modNodata||imgReaderModel1.isNoData(*it1);
806  modNodata=modNodata||imgReaderModel2.isNoData(*it2);
807  if(modNodata){
808  model1buffer.erase(it1);
809  model2buffer.erase(it2);
810  }
811  else{
812  ++it1;
813  ++it2;
814  }
815  }
816  if(model1buffer.size()>minreg_opt[0]&&model2buffer.size()>minreg_opt[0]){
817  errMod=stat.linear_regression_err(model1buffer,model2buffer,c0mod,c1mod);
818  errMod*=regTime_opt[0];
819  }
820  else{//use global regression...
821  c0mod=c0modGlobal;
822  c1mod=c1modGlobal;
823  }
824  }
825  else{
826  c0mod=c0modGlobal;
827  c1mod=c1modGlobal;
828  }
829  double certNorm=(errMod*errMod+errObs*errObs);
830  double certMod=errObs*errObs/certNorm;
831  double certObs=errMod*errMod/certNorm;
832  double regTime=(c0mod+c1mod*estValue)*certObs;
833  double regSensor=(c0obs+c1obs*modValue)*certMod;
834  // double regSensor=(c0obs+c1obs*estValue)*certObs;
835  estWriteBuffer[icol]=regTime+regSensor;
836  //test
837  // if(regTime<regSensor){
838  // cout << "regTime = (" << c0mod << "+" << c1mod << "*" << estValue << ")*" << certObs << " = " << regTime << endl;
839  // cout << "regSensor = (" << c0obs << "+" << c1obs << "*" << modValue << ")*" << certMod << " = " << regSensor << endl;
840  // assert(regTime+regSensor>0);
841  // assert(regTime+regSensor<=1);
842  // }
843 
844  double totalUncertainty=0;
845  if(errMod<eps_opt[0])
846  totalUncertainty=errObs;
847  else if(errObs<eps_opt[0])
848  totalUncertainty=errMod;
849  else{
850  totalUncertainty=1.0/errMod/errMod+1/errObs/errObs;
851  totalUncertainty=sqrt(1.0/totalUncertainty);
852  }
853  uncertWriteBuffer[icol]=totalUncertainty+uncertReadBuffer[icol];
854  }
855  //measurement update
856  if(update&&!imgReaderObs.isNoData(obsLineBuffer[icol])){
857  double kalmanGain=1;
858  double uncertObs=uncertObs_opt[0];
859  if(uncertObsLineBuffer.size()>icol)
860  uncertObs=uncertObsLineBuffer[icol];
861  else if(weight_opt.size()>1||deltaObs_opt.size()){
862  statfactory::StatFactory statobs;
863  statobs.setNoDataValues(obsnodata_opt);
864  double obsMeanValue=(statobs.mean(obsWindowBuffer)-c0obs)/c1obs;
865  double difference=obsMeanValue-modValue;
866  if(modValue){
867  double relativeDifference=difference/modValue;
868  if(deltaObs_opt.size()){
869  assert(deltaObs_opt.size()>1);
870  if(100*relativeDifference<deltaObs_opt[0])//lower bound
871  kalmanGain=0;
872  else if(100*relativeDifference>deltaObs_opt[1])//upper bound
873  kalmanGain=0;
874  }
875  else if(weight_opt.size()){
876  assert(weight_opt.size()>1);
877  if(obsMeanValue<modValue)
878  uncertObs=weight_opt[0]*relativeDifference;
879  else if(obsMeanValue>modValue)
880  uncertObs=weight_opt[1]*relativeDifference;
881  }
882  }
883  if(uncertObs<=0)
884  uncertObs=0;
885  if(verbose_opt[0]>1)
886  cout << "obsMeanValue:" << obsMeanValue << ", modValue: " << modValue << endl;
887  }
888  if(kalmanGain>0){
889  if((uncertWriteBuffer[icol]+uncertObs)>eps_opt[0])
890  kalmanGain=uncertWriteBuffer[icol]/(uncertWriteBuffer[icol]+uncertObs);
891  }
892  assert(kalmanGain<=1);
893  estWriteBuffer[icol]+=kalmanGain*(obsLineBuffer[icol]-estWriteBuffer[icol]);
894  uncertWriteBuffer[icol]*=(1-kalmanGain);
895  }
896  }
897  imgWriterEst.writeData(estWriteBuffer,GDT_Float64,irow,0);
898  imgWriterEst.writeData(uncertWriteBuffer,GDT_Float64,irow,1);
899  progress=static_cast<float>((irow+1.0)/imgWriterEst.nrOfRow());
900  pfnProgress(progress,pszMessage,pProgressArg);
901  }
902 
903  imgWriterEst.close();
904  imgReaderEst.close();
905 
906  if(update){
907  imgReaderObs.close();
908  ++obsindex;
909  }
910  imgReaderModel1.close();
911  imgReaderModel2.close();
912  }
913  }
914  if(find(direction_opt.begin(),direction_opt.end(),"backward")!=direction_opt.end()){
916  cout << "Running backward model" << endl;
917  obsindex=relobsindex.size()-1;
918  //initialization
919  string output;
920  if(outputbw_opt.size()==model_opt.size())
921  output=outputbw_opt.back();
922  else{
923  ostringstream outputstream;
924  outputstream << outputbw_opt[0] << "_";
925  outputstream << setfill('0') << setw(ndigit) << tmodel_opt.back();
926  outputstream << ".tif";
927  // outputstream << outputbw_opt[0] << "_" << tmodel_opt.back() << ".tif";
928  output=outputstream.str();
929  }
930  if(verbose_opt[0])
931  cout << "Opening image " << output << " for writing " << endl;
932  imgWriterEst.open(output,ncol,nrow,2,GDT_Float32,imageType,option_opt);
933  imgWriterEst.setProjectionProj4(projection_opt[0]);
934  imgWriterEst.setGeoTransform(geotransform);
935  imgWriterEst.GDALSetNoDataValue(obsnodata_opt[0]);
936 
937  if(verbose_opt[0]){
938  cout << "processing time " << tmodel_opt.back() << endl;
939  if(obsindex<relobsindex.size())
940  cout << "next observation " << tmodel_opt[relobsindex[obsindex]] << endl;
941  else
942  cout << "There is no next observation" << endl;
943  }
944 
945  try{
946  imgReaderModel1.open(model_opt.back());
947  imgReaderModel1.setNoData(modnodata_opt);
948  if(modoffset_opt.size())
949  imgReaderModel1.setOffset(modoffset_opt[0]);
950  if(modscale_opt.size())
951  imgReaderModel1.setScale(modscale_opt[0]);
952  }
953  catch(string errorString){
954  cerr << errorString << endl;
955  }
956  catch(...){
957  cerr << "Error opening file " << model_opt[0] << endl;
958  }
959 
960  //calculate standard deviation of image to serve as model uncertainty
961  GDALRasterBand* rasterBand;
962  rasterBand=imgReaderModel1.getRasterBand(0);
963  double minValue, maxValue, meanValue, stdDev;
964  void* pProgressData;
965  rasterBand->ComputeStatistics(0,&minValue,&maxValue,&meanValue,&stdDev,pfnProgress,pProgressData);
966  double modRow=0;
967  double modCol=0;
968  if(relobsindex.back()<model_opt.size()-1){//initialize output_opt.back() as model[0]
969  //write last model as output
970  if(verbose_opt[0])
971  cout << "write last model as output" << endl;
972  for(int irow=0;irow<nrow;++irow){
973  vector<double> estReadBuffer;
974  vector<double> estWriteBuffer(ncol);
975  vector<double> uncertWriteBuffer(ncol);
976  // vector<double> lineMask;
977  imgWriterEst.image2geo(0,irow,geox,geoy);
978  imgReaderModel1.geo2image(geox,geoy,modCol,modRow);
979  assert(modRow>=0&&modRow<imgReaderModel1.nrOfRow());
980  try{
981  imgReaderModel1.readData(estReadBuffer,GDT_Float64,modRow);
982  //simple nearest neighbor
983  //stat.nearUp(estReadBuffer,estWriteBuffer);
984 
985  // double oldRowMask=-1;//keep track of row mask to optimize number of line readings
986  for(int icol=0;icol<imgWriterEst.nrOfCol();++icol){
987  imgWriterEst.image2geo(icol,irow,geox,geoy);
988  imgReaderModel1.geo2image(geox,geoy,modCol,modRow);
989  double modValue=estReadBuffer[modCol];
990  if(imgReaderModel1.isNoData(modValue)){
991  estWriteBuffer[icol]=obsnodata_opt[0];
992  uncertWriteBuffer[icol]=uncertNodata_opt[0];
993  }
994  else{
995  estWriteBuffer[icol]=modValue;
996  uncertWriteBuffer[icol]=uncertModel_opt[0]*stdDev;
997  }
998  }
999  imgWriterEst.writeData(estWriteBuffer,GDT_Float64,irow,0);
1000  imgWriterEst.writeData(uncertWriteBuffer,GDT_Float64,irow,1);
1001  }
1002  catch(string errorString){
1003  cerr << errorString << endl;
1004  }
1005  catch(...){
1006  cerr << "Error writing file " << imgWriterEst.getFileName() << endl;
1007  }
1008  }
1009  }
1010  else{//we have an measurement at end time
1011  if(verbose_opt[0])
1012  cout << "we have an measurement at end time" << endl;
1013  imgReaderObs.open(observation_opt.back());
1014  imgReaderObs.getGeoTransform(geotransform);
1015  imgReaderObs.setNoData(obsnodata_opt);
1016  if(obsoffset_opt.size())
1017  imgReaderObs.setOffset(obsoffset_opt[0]);
1018  if(obsscale_opt.size())
1019  imgReaderObs.setScale(obsscale_opt[0]);
1020 
1021  if(regSensor_opt[0]>0){
1022  errObs=regSensor_opt[0]*imgreg.getRMSE(imgReaderModel1,imgReaderObs,c0obs,c1obs,0,0,verbose_opt[0]);
1023  if(errObs<0){
1024  c0obs=0;
1025  c1obs=1;
1026  }
1027  }
1028  else{
1029  c0obs=0;
1030  c1obs=1;
1031  errObs=0;
1032  }
1033  if(verbose_opt[0])
1034  cout << "c0obs, c1obs: " << c0obs << ", " << c1obs << endl;
1035 
1036  for(int irow=0;irow<nrow;++irow){
1037  vector<double> estReadBuffer;
1038  imgWriterEst.image2geo(0,irow,geox,geoy);
1039  imgReaderModel1.geo2image(geox,geoy,modCol,modRow);
1040  assert(modRow>=0&&modRow<imgReaderModel1.nrOfRow());
1041  imgReaderModel1.readData(estReadBuffer,GDT_Float64,modRow);
1042  vector<double> obsLineBuffer;
1043  vector<double> estWriteBuffer(ncol);
1044  vector<double> uncertWriteBuffer(ncol);
1045  vector<double> uncertObsLineBuffer;
1046  // vector<double> lineMask;
1047  // imgReaderObs.readData(estWriteBuffer,GDT_Float64,irow,0);
1048  imgReaderObs.readData(obsLineBuffer,GDT_Float64,irow,0);
1049 
1050  if(imgReaderObs.nrOfBand()>1)
1051  imgReaderObs.readData(uncertObsLineBuffer,GDT_Float64,irow,1);
1052  // double oldRowMask=-1;//keep track of row mask to optimize number of line readings
1053  for(int icol=0;icol<imgWriterEst.nrOfCol();++icol){
1054  imgWriterEst.image2geo(icol,irow,geox,geoy);
1055  imgReaderModel1.geo2image(geox,geoy,modCol,modRow);
1056  assert(modRow>=0&&modRow<imgReaderModel1.nrOfRow());
1057  double modValue=estReadBuffer[modCol];
1058  if(imgReaderModel1.isNoData(modValue)){//model is nodata: retain observation
1059  estWriteBuffer[icol]=obsLineBuffer[icol];
1060  if(imgReaderObs.isNoData(obsLineBuffer[icol])){
1061  estWriteBuffer[icol]=obsnodata_opt[0];
1062  uncertWriteBuffer[icol]=uncertNodata_opt[0];
1063  }
1064  else if(uncertObsLineBuffer.size()>icol)
1065  uncertWriteBuffer[icol]=uncertObsLineBuffer[icol];
1066  else
1067  uncertWriteBuffer[icol]=uncertObs_opt[0];
1068  }
1069  else{//model is valid: calculate estimate from model
1070  double errMod=uncertModel_opt[0]*stdDev;
1071  errMod*=regTime_opt[0];
1072  // double certNorm=(errMod*errMod+errObs*errObs);
1073  // double certMod=errObs*errObs/certNorm;
1074  // double certObs=errMod*errMod/certNorm;
1075  // double regTime=0;
1076  // double regSensor=(c0obs+c1obs*modValue)*certMod;
1077  // estWriteBuffer[icol]=regTime+regSensor;
1078  estWriteBuffer[icol]=modValue;
1079  double totalUncertainty=errMod;
1080  // if(errMod<eps_opt[0])
1081  // totalUncertainty=errObs;
1082  // else if(errObs<eps_opt[0])
1083  // totalUncertainty=errMod;
1084  // else{
1085  // totalUncertainty=1.0/errMod/errMod+1/errObs/errObs;
1086  // totalUncertainty=sqrt(1.0/totalUncertainty);
1087  // }
1088  uncertWriteBuffer[icol]=totalUncertainty;//in case observation is not valid
1089  }
1090  //measurement update
1091  if(!imgReaderObs.isNoData(obsLineBuffer[icol])){
1092  double kalmanGain=1;
1093  double uncertObs=uncertObs_opt[0];
1094  if(uncertObsLineBuffer.size()>icol)
1095  uncertObs=uncertObsLineBuffer[icol];
1096  else if(weight_opt.size()>1||deltaObs_opt.size()){
1097  vector<double> obsWindowBuffer;//buffer for observation to calculate average corresponding to model pixel
1098  int minCol=(icol>down_opt[0]/2) ? icol-down_opt[0]/2 : 0;
1099  int maxCol=(icol+down_opt[0]/2<imgReaderObs.nrOfCol()) ? icol+down_opt[0]/2 : imgReaderObs.nrOfCol()-1;
1100  int minRow=(irow>down_opt[0]/2) ? irow-down_opt[0]/2 : 0;
1101  int maxRow=(irow+down_opt[0]/2<imgReaderObs.nrOfRow()) ? irow+down_opt[0]/2 : imgReaderObs.nrOfRow()-1;
1102 
1103  imgReaderObs.readDataBlock(obsWindowBuffer,GDT_Float64,minCol,maxCol,minRow,maxRow,0);
1104  statfactory::StatFactory statobs;
1105  statobs.setNoDataValues(obsnodata_opt);
1106  double obsMeanValue=(statobs.mean(obsWindowBuffer)-c0obs)/c1obs;
1107  double difference=obsMeanValue-modValue;
1108  if(modValue){
1109  double relativeDifference=difference/modValue;
1110  if(deltaObs_opt.size()){
1111  assert(deltaObs_opt.size()>1);
1112  if(100*relativeDifference<deltaObs_opt[0])//lower bound
1113  kalmanGain=0;
1114  else if(100*relativeDifference>deltaObs_opt[1])//upper bound
1115  kalmanGain=0;
1116  }
1117  else if(weight_opt.size()){
1118  assert(weight_opt.size()>1);
1119  if(obsMeanValue<modValue)
1120  uncertObs=weight_opt[0]*relativeDifference;
1121  else if(obsMeanValue>modValue)
1122  uncertObs=weight_opt[1]*relativeDifference;
1123  }
1124  }
1125  if(uncertObs<=0)
1126  uncertObs=0;
1127  if(verbose_opt[0]>1)
1128  cout << "obsMeanValue:" << obsMeanValue << ", modValue: " << modValue << endl;
1129  }
1130  if(kalmanGain>0){
1131  if((uncertWriteBuffer[icol]+uncertObs)>eps_opt[0])
1132  kalmanGain=uncertWriteBuffer[icol]/(uncertWriteBuffer[icol]+uncertObs);
1133  }
1134  assert(kalmanGain<=1);
1135  estWriteBuffer[icol]+=kalmanGain*(obsLineBuffer[icol]-estWriteBuffer[icol]);
1136  uncertWriteBuffer[icol]*=(1-kalmanGain);
1137  }
1138  }
1139  imgWriterEst.writeData(estWriteBuffer,GDT_Float64,irow,0);
1140  imgWriterEst.writeData(uncertWriteBuffer,GDT_Float64,irow,1);
1141  }
1142  imgReaderObs.close();
1143  --obsindex;
1144  }
1145  imgReaderModel1.close();
1146  imgWriterEst.close();
1147 
1148  for(int modindex=model_opt.size()-2;modindex>=0;--modindex){
1149  if(verbose_opt[0]){
1150  cout << "processing time " << tmodel_opt[modindex] << endl;
1151  if(obsindex<relobsindex.size())
1152  cout << "next observation " << tmodel_opt[relobsindex[obsindex]] << endl;
1153  else
1154  cout << "There is no next observation" << endl;
1155  }
1156  string output;
1157  if(outputbw_opt.size()==model_opt.size())
1158  output=outputbw_opt[modindex];
1159  else{
1160  ostringstream outputstream;
1161  outputstream << outputbw_opt[0] << "_";
1162  outputstream << setfill('0') << setw(ndigit) << tmodel_opt[modindex];
1163  outputstream << ".tif";
1164  // outputstream << outputbw_opt[0] << "_" << tmodel_opt[modindex] << ".tif";
1165  output=outputstream.str();
1166  }
1167 
1168  //two band output band0=estimation, band1=uncertainty
1169  imgWriterEst.open(output,ncol,nrow,2,GDT_Float32,imageType,option_opt);
1170  imgWriterEst.setProjectionProj4(projection_opt[0]);
1171  imgWriterEst.setGeoTransform(geotransform);
1172  imgWriterEst.GDALSetNoDataValue(obsnodata_opt[0]);
1173 
1174  //calculate regression between two subsequence model inputs
1175  imgReaderModel1.open(model_opt[modindex+1]);
1176  imgReaderModel1.setNoData(modnodata_opt);
1177  if(modoffset_opt.size())
1178  imgReaderModel1.setOffset(modoffset_opt[0]);
1179  if(modscale_opt.size())
1180  imgReaderModel1.setScale(modscale_opt[0]);
1181  imgReaderModel2.open(model_opt[modindex]);
1182  imgReaderModel2.setNoData(modnodata_opt);
1183  if(modoffset_opt.size())
1184  imgReaderModel2.setOffset(modoffset_opt[0]);
1185  if(modscale_opt.size())
1186  imgReaderModel2.setScale(modscale_opt[0]);
1187  //calculate regression
1188  //we could re-use the points from second image from last run, but
1189  //to keep it general, we must redo it (overlap might have changed)
1190 
1191  pfnProgress(progress,pszMessage,pProgressArg);
1192 
1193  if(verbose_opt[0])
1194  cout << "Calculating regression for " << imgReaderModel1.getFileName() << " " << imgReaderModel2.getFileName() << endl;
1195 
1196  double errMod=imgreg.getRMSE(imgReaderModel1,imgReaderModel2,c0modGlobal,c1modGlobal,0,0);
1197  errMod*=regTime_opt[0];
1198 
1199  // double errMod=imgreg.getRMSE(imgReaderModel1,imgReaderModel2,c0modGlobal,c1modGlobal,verbose_opt[0]);
1200  if(verbose_opt[0])
1201  cout << "c0modGlobal, c1modGlobal: " << c0modGlobal << ", " << c1modGlobal << endl;
1202 
1203  bool update=false;
1204  if(obsindex<relobsindex.size()){
1205  update=(relobsindex[obsindex]==modindex);
1206  }
1207  if(update){
1208  if(verbose_opt[0])
1209  cout << "***update " << relobsindex[obsindex] << " = " << modindex << " " << observation_opt[obsindex] << " ***" << endl;
1210 
1211  imgReaderObs.open(observation_opt[obsindex]);
1212  imgReaderObs.getGeoTransform(geotransform);
1213  imgReaderObs.setNoData(obsnodata_opt);
1214  if(obsoffset_opt.size())
1215  imgReaderObs.setOffset(obsoffset_opt[0]);
1216  if(obsscale_opt.size())
1217  imgReaderObs.setScale(obsscale_opt[0]);
1218  //calculate regression between model and observation
1219  if(verbose_opt[0])
1220  cout << "Calculating regression for " << imgReaderModel2.getFileName() << " " << imgReaderObs.getFileName() << endl;
1221  if(regSensor_opt[0]>0){
1222  errObs=regSensor_opt[0]*imgreg.getRMSE(imgReaderModel2,imgReaderObs,c0obs,c1obs,0,0,verbose_opt[0]);
1223  if(errObs<0){
1224  c0obs=0;
1225  c1obs=1;
1226  }
1227  }
1228  else{
1229  c0obs=0;
1230  c1obs=1;
1231  errObs=0;
1232  }
1233  if(verbose_opt[0])
1234  cout << "c0obs, c1obs: " << c0obs << ", " << c1obs << endl;
1235  }
1236  //prediction (also to fill cloudy pixels in update mode)
1237  string input;
1238  if(outputbw_opt.size()==model_opt.size())
1239  input=outputbw_opt[modindex+1];
1240  else{
1241  ostringstream outputstream;
1242  outputstream << outputbw_opt[0] << "_";
1243  outputstream << setfill('0') << setw(ndigit) << tmodel_opt[modindex+1];
1244  outputstream << ".tif";
1245  // outputstream << outputbw_opt[0] << "_" << tmodel_opt[modindex+1] << ".tif";
1246  input=outputstream.str();
1247  }
1248  ImgReaderGdal imgReaderEst(input);
1249  imgReaderEst.setNoData(obsnodata_opt);
1250  if(obsoffset_opt.size())
1251  imgReaderEst.setOffset(obsoffset_opt[0]);
1252  if(obsscale_opt.size())
1253  imgReaderEst.setScale(obsscale_opt[0]);
1254 
1255  vector< vector<double> > obsLineVector(down_opt[0]);
1256  vector<double> obsLineBuffer;
1257  vector<double> obsWindowBuffer;//buffer for observation to calculate average corresponding to model pixel
1258  vector<double> model1LineBuffer;
1259  vector<double> model2LineBuffer;
1260  vector<double> model1buffer;//buffer for model 1 to calculate time regression based on window
1261  vector<double> model2buffer;//buffer for model 2 to calculate time regression based on window
1262  vector<double> uncertObsLineBuffer;
1263  vector<double> estReadBuffer;
1264  // vector<double> estWindowBuffer;//buffer for estimate to calculate average corresponding to model pixel
1265  vector<double> uncertReadBuffer;
1266  vector<double> estWriteBuffer(ncol);
1267  vector<double> uncertWriteBuffer(ncol);
1268  // vector<double> lineMask;
1269 
1270  //initialize obsLineVector
1271  if(update){
1272  assert(down_opt[0]%2);//window size must be odd
1273  for(int iline=-down_opt[0]/2;iline<down_opt[0]/2+1;++iline){
1274  if(iline<0)//replicate line 0
1275  imgReaderObs.readData(obsLineVector[iline+down_opt[0]/2],GDT_Float64,0,0);
1276  else
1277  imgReaderObs.readData(obsLineVector[iline+down_opt[0]/2],GDT_Float64,iline,0);
1278  }
1279  }
1280  for(int irow=0;irow<imgWriterEst.nrOfRow();++irow){
1281  assert(irow<imgReaderEst.nrOfRow());
1282  //do not read from imgReaderObs, because we read entire window for each pixel...
1283  imgReaderEst.readData(estReadBuffer,GDT_Float64,irow,0);
1284  imgReaderEst.readData(uncertReadBuffer,GDT_Float64,irow,1);
1285  //read model2 in case current estimate is nodata
1286  imgReaderEst.image2geo(0,irow,geox,geoy);
1287  imgReaderModel2.geo2image(geox,geoy,modCol,modRow);
1288  assert(modRow>=0&&modRow<imgReaderModel2.nrOfRow());
1289  imgReaderModel2.readData(model2LineBuffer,GDT_Float64,modRow,0);
1290 
1291  imgReaderModel1.geo2image(geox,geoy,modCol,modRow);
1292  assert(modRow>=0&&modRow<imgReaderModel1.nrOfRow());
1293  imgReaderModel1.readData(model1LineBuffer,GDT_Float64,modRow,0);
1294 
1295  if(update){
1296  int maxRow=(irow+down_opt[0]/2<imgReaderEst.nrOfRow()) ? irow+down_opt[0]/2 : imgReaderEst.nrOfRow()-1;
1297  obsLineVector.erase(obsLineVector.begin());
1298  imgReaderObs.readData(obsLineBuffer,GDT_Float64,maxRow,0);
1299  obsLineVector.push_back(obsLineBuffer);
1300  obsLineBuffer=obsLineVector[down_opt[0]/2];
1301  // imgReaderObs.readData(obsLineBuffer,GDT_Float64,irow,0);
1302  if(imgReaderObs.nrOfBand()>1)
1303  imgReaderObs.readData(uncertObsLineBuffer,GDT_Float64,irow,1);
1304  }
1305  // double oldRowMask=-1;//keep track of row mask to optimize number of line readings
1306  for(int icol=0;icol<imgWriterEst.nrOfCol();++icol){
1307  imgReaderEst.image2geo(icol,irow,geox,geoy);
1308  int minCol=(icol>down_opt[0]/2) ? icol-down_opt[0]/2 : 0;
1309  int maxCol=(icol+down_opt[0]/2<imgReaderEst.nrOfCol()) ? icol+down_opt[0]/2 : imgReaderEst.nrOfCol()-1;
1310  int minRow=(irow>down_opt[0]/2) ? irow-down_opt[0]/2 : 0;
1311  int maxRow=(irow+down_opt[0]/2<imgReaderEst.nrOfRow()) ? irow+down_opt[0]/2 : imgReaderEst.nrOfRow()-1;
1312  if(update){
1313  obsWindowBuffer.clear();
1314  for(int iline=0;iline<obsLineVector.size();++iline){
1315  for(int isample=minCol;isample<=maxCol;++isample){
1316  assert(isample<obsLineVector[iline].size());
1317  obsWindowBuffer.push_back(obsLineVector[iline][isample]);
1318  }
1319  }
1320  // imgReaderObs.readDataBlock(obsWindowBuffer,GDT_Float64,minCol,maxCol,minRow,maxRow,0);
1321  }
1322  double estValue=estReadBuffer[icol];
1323  double estMeanValue=0;//stat.mean(estWindowBuffer);
1324  double nvalid=0;
1325  //time update
1326  imgReaderModel2.geo2image(geox,geoy,modCol,modRow);
1327  assert(modRow>=0&&modRow<imgReaderModel2.nrOfRow());
1328  double modValue=model2LineBuffer[modCol];
1329  bool estNodata=imgReaderEst.isNoData(estValue);
1330  if(estNodata){
1331  //we have not found any valid data yet, better here to take the current model value if valid
1332  if(imgReaderModel2.isNoData(modValue)){//if both estimate and model are no-data, set obs to nodata
1333  estWriteBuffer[icol]=obsnodata_opt[0];
1334  uncertWriteBuffer[icol]=uncertNodata_opt[0];
1335  }
1336  else{
1337  estWriteBuffer[icol]=modValue;
1338  uncertWriteBuffer[icol]=uncertModel_opt[0]*stdDev;
1339  }
1340  }
1341  else{
1342  if(window_opt[0]>0){
1343  try{
1344  // imgReaderModel2.geo2image(geox,geoy,modCol,modRow);//did that already
1345  minCol=(modCol>window_opt[0]/2) ? modCol-window_opt[0]/2 : 0;
1346  maxCol=(modCol+window_opt[0]/2<imgReaderModel2.nrOfCol()) ? modCol+window_opt[0]/2 : imgReaderModel2.nrOfCol()-1;
1347  minRow=(modRow>window_opt[0]/2) ? modRow-window_opt[0]/2 : 0;
1348  maxRow=(modRow+window_opt[0]/2<imgReaderModel2.nrOfRow()) ? modRow+window_opt[0]/2 : imgReaderModel2.nrOfRow()-1;
1349  imgReaderModel2.readDataBlock(model2buffer,GDT_Float64,minCol,maxCol,minRow,maxRow,0);
1350 
1351  imgReaderModel1.geo2image(geox,geoy,modCol,modRow);
1352  assert(modRow>=0&&modRow<imgReaderModel1.nrOfRow());
1353  minCol=(modCol>window_opt[0]/2) ? modCol-window_opt[0]/2 : 0;
1354  maxCol=(modCol+window_opt[0]/2<imgReaderModel1.nrOfCol()) ? modCol+window_opt[0]/2 : imgReaderModel1.nrOfCol()-1;
1355  minRow=(modRow>window_opt[0]/2) ? modRow-window_opt[0]/2 : 0;
1356  maxRow=(modRow+window_opt[0]/2<imgReaderModel1.nrOfRow()) ? modRow+window_opt[0]/2 : imgReaderModel1.nrOfRow()-1;
1357  imgReaderModel1.readDataBlock(model1buffer,GDT_Float64,minCol,maxCol,minRow,maxRow,0);
1358  // imgReaderEst.image2geo(icol,irow,geox,geoy);
1359  }
1360  catch(string errorString){
1361  cerr << "Error reading data block for " << minCol << "-" << maxCol << ", " << minRow << "-" << maxRow << endl;
1362  }
1363  //erase no-data from buffer
1364  vector<double>::iterator it1=model1buffer.begin();
1365  vector<double>::iterator it2=model2buffer.begin();
1366  while(it1!=model1buffer.end()&&it2!=model2buffer.end()){
1367  //erase nodata
1368  bool modNodata=false;
1369  modNodata=modNodata||imgReaderModel1.isNoData(*it1);
1370  modNodata=modNodata||imgReaderModel2.isNoData(*it2);
1371  if(modNodata){
1372  model1buffer.erase(it1);
1373  model2buffer.erase(it2);
1374  }
1375  else{
1376  ++it1;
1377  ++it2;
1378  }
1379  }
1380  if(model1buffer.size()>minreg_opt[0]&&model2buffer.size()>minreg_opt[0]){
1381  errMod=stat.linear_regression_err(model1buffer,model2buffer,c0mod,c1mod);
1382  errMod*=regTime_opt[0];
1383  }
1384  else{//use global regression...
1385  c0mod=c0modGlobal;
1386  c1mod=c1modGlobal;
1387  }
1388  }
1389  else{
1390  c0mod=c0modGlobal;
1391  c1mod=c1modGlobal;
1392  }
1393  double certNorm=(errMod*errMod+errObs*errObs);
1394  double certMod=errObs*errObs/certNorm;
1395  double certObs=errMod*errMod/certNorm;
1396  double regTime=(c0mod+c1mod*estValue)*certObs;
1397 
1398  // double regSensor=(c0obs+c1obs*estValue)*certObs;
1399  double regSensor=(c0obs+c1obs*modValue)*certMod;
1400  estWriteBuffer[icol]=regTime+regSensor;
1401  double totalUncertainty=0;
1402  if(errMod<eps_opt[0])
1403  totalUncertainty=errObs;
1404  else if(errObs<eps_opt[0])
1405  totalUncertainty=errMod;
1406  else{
1407  totalUncertainty=1.0/errMod/errMod+1/errObs/errObs;
1408  totalUncertainty=sqrt(1.0/totalUncertainty);
1409  }
1410  uncertWriteBuffer[icol]=totalUncertainty+uncertReadBuffer[icol];
1411  }
1412  //measurement update
1413  if(update&&!imgReaderObs.isNoData(obsLineBuffer[icol])){
1414  double kalmanGain=1;
1415  double uncertObs=uncertObs_opt[0];
1416  if(uncertObsLineBuffer.size()>icol)
1417  uncertObs=uncertObsLineBuffer[icol];
1418  else if(weight_opt.size()>1||deltaObs_opt.size()){
1419  statfactory::StatFactory statobs;
1420  statobs.setNoDataValues(obsnodata_opt);
1421  double obsMeanValue=(statobs.mean(obsWindowBuffer)-c0obs)/c1obs;
1422  double difference=obsMeanValue-modValue;
1423  if(modValue){
1424  double relativeDifference=difference/modValue;
1425  if(deltaObs_opt.size()){
1426  assert(deltaObs_opt.size()>1);
1427  if(100*relativeDifference<deltaObs_opt[0])//lower bound
1428  kalmanGain=0;
1429  else if(100*relativeDifference>deltaObs_opt[1])//upper bound
1430  kalmanGain=0;
1431  }
1432  else if(weight_opt.size()){
1433  assert(weight_opt.size()>1);
1434  if(obsMeanValue<modValue)
1435  uncertObs=weight_opt[0]*relativeDifference;
1436  else if(obsMeanValue>modValue)
1437  uncertObs=weight_opt[1]*relativeDifference;
1438  }
1439  }
1440  if(uncertObs<=0)
1441  uncertObs=0;
1442  if(verbose_opt[0]>1)
1443  cout << "obsMeanValue:" << obsMeanValue << ", modValue: " << modValue << endl;
1444  }
1445  if(kalmanGain>0){
1446  if((uncertWriteBuffer[icol]+uncertObs)>eps_opt[0])
1447  kalmanGain=uncertWriteBuffer[icol]/(uncertWriteBuffer[icol]+uncertObs);
1448  }
1449  assert(kalmanGain<=1);
1450  estWriteBuffer[icol]+=kalmanGain*(obsLineBuffer[icol]-estWriteBuffer[icol]);
1451  uncertWriteBuffer[icol]*=(1-kalmanGain);
1452  }
1453  }
1454  imgWriterEst.writeData(estWriteBuffer,GDT_Float64,irow,0);
1455  imgWriterEst.writeData(uncertWriteBuffer,GDT_Float64,irow,1);
1456  progress=static_cast<float>((irow+1.0)/imgWriterEst.nrOfRow());
1457  pfnProgress(progress,pszMessage,pProgressArg);
1458  }
1459 
1460  imgWriterEst.close();
1461  imgReaderEst.close();
1462 
1463  if(update){
1464  imgReaderObs.close();
1465  --obsindex;
1466  }
1467  imgReaderModel1.close();
1468  imgReaderModel2.close();
1469  }
1470  }
1471  if(find(direction_opt.begin(),direction_opt.end(),"smooth")!=direction_opt.end()){
1473  cout << "Running smooth model" << endl;
1474  obsindex=0;
1475  for(int modindex=0;modindex<model_opt.size();++modindex){
1476  if(verbose_opt[0]){
1477  cout << "processing time " << tmodel_opt[modindex] << endl;
1478  if(obsindex<relobsindex.size())
1479  cout << "next observation " << tmodel_opt[relobsindex[obsindex]] << endl;
1480  else
1481  cout << "There is no next observation" << endl;
1482  }
1483  string output;
1484  if(outputfb_opt.size()==model_opt.size())
1485  output=outputfb_opt[modindex];
1486  else{
1487  ostringstream outputstream;
1488  outputstream << outputfb_opt[0] << "_";
1489  outputstream << setfill('0') << setw(ndigit) << tmodel_opt[modindex];
1490  outputstream << ".tif";
1491  // outputstream << outputfb_opt[0] << "_" << tmodel_opt[modindex] << ".tif";
1492  output=outputstream.str();
1493  }
1494 
1495  //two band output band0=estimation, band1=uncertainty
1496  imgWriterEst.open(output,ncol,nrow,2,GDT_Float32,imageType,option_opt);
1497  imgWriterEst.setProjectionProj4(projection_opt[0]);
1498  imgWriterEst.setGeoTransform(geotransform);
1499  imgWriterEst.GDALSetNoDataValue(obsnodata_opt[0]);
1500 
1501  //open forward and backward estimates
1502  //we assume forward in model and backward in observation...
1503 
1504  string inputfw;
1505  string inputbw;
1506  if(outputfw_opt.size()==model_opt.size())
1507  inputfw=outputfw_opt[modindex];
1508  else{
1509  ostringstream outputstream;
1510  outputstream << outputfw_opt[0] << "_";
1511  outputstream << setfill('0') << setw(ndigit) << tmodel_opt[modindex];
1512  outputstream << ".tif";
1513  // outputstream << outputfw_opt[0] << "_" << tmodel_opt[modindex] << ".tif";
1514  inputfw=outputstream.str();
1515  }
1516  if(outputbw_opt.size()==model_opt.size())
1517  inputbw=outputbw_opt[modindex];
1518  else{
1519  ostringstream outputstream;
1520  outputstream << outputbw_opt[0] << "_";
1521  outputstream << setfill('0') << setw(ndigit) << tmodel_opt[modindex];
1522  outputstream << ".tif";
1523  // outputstream << outputbw_opt[0] << "_" << tmodel_opt[modindex] << ".tif";
1524  inputbw=outputstream.str();
1525  }
1526  ImgReaderGdal imgReaderForward(inputfw);
1527  ImgReaderGdal imgReaderBackward(inputbw);
1528  imgReaderForward.setNoData(obsnodata_opt);
1529  if(obsoffset_opt.size())
1530  imgReaderForward.setOffset(obsoffset_opt[0]);
1531  if(obsscale_opt.size())
1532  imgReaderForward.setScale(obsscale_opt[0]);
1533  imgReaderBackward.setNoData(obsnodata_opt);
1534  if(obsoffset_opt.size())
1535  imgReaderBackward.setOffset(obsoffset_opt[0]);
1536  if(obsscale_opt.size())
1537  imgReaderBackward.setScale(obsscale_opt[0]);
1538 
1539  vector<double> estForwardBuffer;
1540  vector<double> estBackwardBuffer;
1541  vector<double> uncertObsLineBuffer;
1542  vector<double> uncertForwardBuffer;
1543  vector<double> uncertBackwardBuffer;
1544  vector<double> uncertReadBuffer;
1545  vector<double> estWriteBuffer(ncol);
1546  vector<double> uncertWriteBuffer(ncol);
1547  // vector<double> lineMask;
1548 
1549  bool update=false;
1550  if(obsindex<relobsindex.size()){
1551  update=(relobsindex[obsindex]==modindex);
1552  }
1553 
1554  if(update){
1555  if(verbose_opt[0])
1556  cout << "***update " << relobsindex[obsindex] << " = " << modindex << " " << observation_opt[obsindex] << " ***" << endl;
1557  imgReaderObs.open(observation_opt[obsindex]);
1558  imgReaderObs.getGeoTransform(geotransform);
1559  imgReaderObs.setNoData(obsnodata_opt);
1560  if(obsoffset_opt.size())
1561  imgReaderObs.setOffset(obsoffset_opt[0]);
1562  if(obsscale_opt.size())
1563  imgReaderObs.setScale(obsscale_opt[0]);
1564  //calculate regression between model and observation
1565  }
1566 
1567  pfnProgress(progress,pszMessage,pProgressArg);
1568 
1569  for(int irow=0;irow<imgWriterEst.nrOfRow();++irow){
1570  assert(irow<imgReaderForward.nrOfRow());
1571  assert(irow<imgReaderBackward.nrOfRow());
1572  imgReaderForward.readData(estForwardBuffer,GDT_Float64,irow,0);
1573  imgReaderBackward.readData(estBackwardBuffer,GDT_Float64,irow,0);
1574  imgReaderForward.readData(uncertForwardBuffer,GDT_Float64,irow,1);
1575  imgReaderBackward.readData(uncertBackwardBuffer,GDT_Float64,irow,1);
1576 
1577  if(update){
1578  imgReaderObs.readData(estWriteBuffer,GDT_Float64,irow,0);
1579  if(imgReaderObs.nrOfBand()>1)
1580  imgReaderObs.readData(uncertObsLineBuffer,GDT_Float64,irow,1);
1581  }
1582 
1583  // double oldRowMask=-1;//keep track of row mask to optimize number of line readings
1584  for(int icol=0;icol<imgWriterEst.nrOfCol();++icol){
1585  imgWriterEst.image2geo(icol,irow,geox,geoy);
1586  double A=estForwardBuffer[icol];
1587  double B=estBackwardBuffer[icol];
1588  double C=uncertForwardBuffer[icol]*uncertForwardBuffer[icol];
1589  double D=uncertBackwardBuffer[icol]*uncertBackwardBuffer[icol];
1590  double uncertObs=uncertObs_opt[0];
1591 
1592  // if(update){//check for nodata in observation
1593  // if(imgReaderObs.isNoData(estWriteBuffer[icol]))
1594  // uncertObs=uncertNodata_opt[0];
1595  // else if(uncertObsLineBuffer.size()>icol)
1596  // uncertObs=uncertObsLineBuffer[icol];
1597  // }
1598 
1599  double noemer=(C+D);
1600  //todo: consistently check for division by zero...
1601  if(imgReaderForward.isNoData(A)&&imgReaderBackward.isNoData(B)){
1602  estWriteBuffer[icol]=obsnodata_opt[0];
1603  uncertWriteBuffer[icol]=uncertNodata_opt[0];
1604  }
1605  else if(imgReaderForward.isNoData(A)){
1606  estWriteBuffer[icol]=B;
1607  uncertWriteBuffer[icol]=uncertBackwardBuffer[icol];
1608  }
1609  else if(imgReaderForward.isNoData(B)){
1610  estWriteBuffer[icol]=A;
1611  uncertWriteBuffer[icol]=uncertForwardBuffer[icol];
1612  }
1613  else{
1614  if(noemer<eps_opt[0]){//simple average if both uncertainties are ~>0
1615  estWriteBuffer[icol]=0.5*(A+B);
1616  uncertWriteBuffer[icol]=uncertObs;
1617  }
1618  else{
1619  estWriteBuffer[icol]=(A*D+B*C)/noemer;
1620  double P=0;
1621  if(C>eps_opt[0])
1622  P+=1.0/C;
1623  if(D>eps_opt[0])
1624  P+=1.0/D;
1625  if(uncertObs*uncertObs>eps_opt[0])
1626  P-=1.0/uncertObs/uncertObs;
1627  if(P>eps_opt[0])
1628  P=sqrt(1.0/P);
1629  else
1630  P=0;
1631  uncertWriteBuffer[icol]=P;
1632  }
1633  }
1634  }
1635  imgWriterEst.writeData(estWriteBuffer,GDT_Float64,irow,0);
1636  imgWriterEst.writeData(uncertWriteBuffer,GDT_Float64,irow,1);
1637  progress=static_cast<float>((irow+1.0)/imgWriterEst.nrOfRow());
1638  pfnProgress(progress,pszMessage,pProgressArg);
1639  }
1640 
1641  imgWriterEst.close();
1642  imgReaderForward.close();
1643  imgReaderBackward.close();
1644  if(update){
1645  imgReaderObs.close();
1646  ++obsindex;
1647  }
1648  }
1649  }
1650  // if(mask_opt.size())
1651  // maskReader.close();
1652 }
1653