MUQ  0.4.3
ObservationInformation.cpp
Go to the documentation of this file.
2 
3 using namespace muq::Approximation;
4 
5 
6 void ObservationInformation::FillSelfCov(std::shared_ptr<KernelBase> kernel,
7  Eigen::Ref<Eigen::MatrixXd> covBlock)
8 {
9  covBlock = H->Apply( H->Apply( BuildBaseCovariance(kernel) ).transpose() ) + obsCov;
10 }
11 
12 void ObservationInformation::FillCrossCov(Eigen::Ref<const Eigen::VectorXd> const& otherLoc,
13  std::shared_ptr<KernelBase> kernel,
14  Eigen::Ref<Eigen::MatrixXd> covBlock)
15 {
16  covBlock = H->Apply( BuildBaseCovariance(otherLoc, kernel) );
17 }
18 
19 void ObservationInformation::FillCrossCov(std::shared_ptr<ObservationInformation> otherObs,
20  std::shared_ptr<KernelBase> kernel,
21  Eigen::Ref<Eigen::MatrixXd> covBlock)
22 {
23  covBlock = otherObs->H->Apply( H->Apply( BuildBaseCovariance(otherObs, kernel) ).transpose() ).transpose();
24 };
25 
26 Eigen::MatrixXd ObservationInformation::BuildBaseCovariance(Eigen::Ref<const Eigen::VectorXd> const& otherLoc,
27  std::shared_ptr<KernelBase> kernel)
28 {
29  return kernel->BuildCovariance(loc, otherLoc);
30 }
31 
32 Eigen::MatrixXd ObservationInformation::BuildBaseCovariance(std::shared_ptr<KernelBase> kernel)
33 {
34  return kernel->BuildCovariance(loc,loc);
35 }
36 
37 Eigen::MatrixXd ObservationInformation::BuildBaseCovariance(std::shared_ptr<ObservationInformation> otherObs,
38  std::shared_ptr<KernelBase> kernel)
39 {
40  // Check to see if the other observation is a derivative observation
41  auto derivObs = std::dynamic_pointer_cast<DerivativeObservation>(otherObs);
42  if(derivObs){
43  return derivObs->BuildBaseCovariance(shared_from_this(), kernel).transpose();
44  }else{
45  return kernel->BuildCovariance(loc, otherObs->loc);
46  }
47 }
48 
49 
50 
51 Eigen::MatrixXd DerivativeObservation::BuildBaseCovariance(Eigen::Ref<const Eigen::VectorXd> const& otherLoc,
52  std::shared_ptr<KernelBase> kernel)
53 {
54  Eigen::MatrixXd output(derivCoords.size() * kernel->coDim, kernel->coDim);
55  for(int i=0; i<derivCoords.size(); ++i)
56  output.block(i*derivCoords.size(),0,kernel->coDim,kernel->coDim) = kernel->GetPosDerivative(otherLoc, loc, derivCoords.at(i));
57 
58  return output;
59 }
60 
61 Eigen::MatrixXd DerivativeObservation::BuildBaseCovariance(std::shared_ptr<KernelBase> kernel)
62 {
63  Eigen::MatrixXd output(derivCoords.size() * kernel->coDim, derivCoords.size() * kernel->coDim);
64  for(int i=0; i<derivCoords.size(); ++i){
65  for(int j=0; j<derivCoords.size(); ++j){
66 
67  std::vector<int> allWrts = derivCoords.at(i);
68  allWrts.insert(allWrts.end(), derivCoords.at(j).begin(), derivCoords.at(j).end());
69  output.block(i*derivCoords.size(),j*derivCoords.size(),kernel->coDim,kernel->coDim) = kernel->GetPosDerivative(loc, loc, allWrts);
70  }
71  }
72 
73  return output;
74 
75 }
76 
77 Eigen::MatrixXd DerivativeObservation::BuildBaseCovariance(std::shared_ptr<ObservationInformation> otherObs,
78  std::shared_ptr<KernelBase> kernel)
79 {
80 
81  // Check to see if the other observation is a derivative observation
82  auto derivObs = std::dynamic_pointer_cast<DerivativeObservation>(otherObs);
83  if(derivObs){
84 
85  Eigen::MatrixXd output(derivCoords.size() * kernel->coDim, derivObs->derivCoords.size() * kernel->coDim);
86 
87  for(int i=0; i<derivCoords.size(); ++i){
88  for(int j=0; j<derivObs->derivCoords.size(); ++j){
89 
90  std::vector<int> allWrts = derivCoords.at(i);
91  allWrts.insert(allWrts.end(), derivObs->derivCoords.at(j).begin(), derivObs->derivCoords.at(j).end());
92  output.block(i*derivCoords.size(),j*derivObs->derivCoords.size(),kernel->coDim,kernel->coDim) = kernel->GetPosDerivative(loc, otherObs->loc, allWrts);
93  }
94  }
95 
96  return output;
97 
98  }else{
99  return BuildBaseCovariance(otherObs->loc, kernel);
100  }
101 }
virtual Eigen::MatrixXd BuildBaseCovariance(Eigen::Ref< const Eigen::VectorXd > const &otherObs, std::shared_ptr< KernelBase > kernel) override
std::vector< std::vector< int > > derivCoords
virtual Eigen::MatrixXd BuildBaseCovariance(Eigen::Ref< const Eigen::VectorXd > const &otherObs, std::shared_ptr< KernelBase > kernel)
virtual void FillCrossCov(Eigen::Ref< const Eigen::VectorXd > const &otherLoc, std::shared_ptr< KernelBase > kernel, Eigen::Ref< Eigen::MatrixXd > covBlock)
std::shared_ptr< muq::Modeling::LinearOperator > H
virtual void FillSelfCov(std::shared_ptr< KernelBase > kernel, Eigen::Ref< Eigen::MatrixXd > covBlock)