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