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