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