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