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