MUQ  0.4.3
CombineVectors.cpp
Go to the documentation of this file.
2 
3 using namespace muq::Modeling;
4 
5 CombineVectors::CombineVectors(Eigen::VectorXi const& inputSizes) : ModPiece(inputSizes, Eigen::VectorXi::Constant(1, inputSizes.sum())) {}
6 
8  Eigen::VectorXd out(outputSizes(0));
9  unsigned int ind = 0;
10  for( unsigned int i=0; i<inputs.size(); ++i ) {
11  out.segment(ind, inputs[i].get().size()) = inputs[i].get();
12  ind += inputs[i].get().size();
13  }
14 
15  outputs.resize(1);
16  outputs[0] = out;
17 }
18 
19 void CombineVectors::JacobianImpl(unsigned int const outwrt, unsigned int const inwrt, ref_vector<Eigen::VectorXd> const& inputs) {
20  assert(outwrt==0);
21 
22  jacobian = Eigen::MatrixXd::Zero(outputSizes(0), inputSizes(inwrt));
23  const unsigned int ind = inputSizes.segment(0, inwrt).sum();
24  jacobian.block(ind, 0, inputSizes(inwrt), inputSizes(inwrt)) = Eigen::MatrixXd::Identity(inputSizes(inwrt), inputSizes(inwrt));
25 }
26 
27 void CombineVectors::GradientImpl(unsigned int const outwrt, unsigned int const inwrt, ref_vector<Eigen::VectorXd> const& inputs, Eigen::VectorXd const& sens) {
28  assert(outwrt==0);
29 
30  const unsigned int ind = inputSizes.segment(0, inwrt).sum();
31  gradient = sens.segment(ind, inputSizes(inwrt));
32 }
33 
34 void CombineVectors::ApplyJacobianImpl(unsigned int const outwrt, unsigned int const inwrt, ref_vector<Eigen::VectorXd> const& inputs, Eigen::VectorXd const& targ) {
35  assert(outwrt==0);
36 
37  const unsigned int ind = inputSizes.segment(0, inwrt).sum();
38  jacobianAction = Eigen::VectorXd::Zero(outputSizes(0));
39  jacobianAction.segment(ind, inputSizes(inwrt)) = targ;
40 }
virtual void GradientImpl(unsigned int const outwrt, unsigned int const inwrt, ref_vector< Eigen::VectorXd > const &inputs, Eigen::VectorXd const &sens) override
virtual void ApplyJacobianImpl(unsigned int const outwrt, unsigned int const inwrt, ref_vector< Eigen::VectorXd > const &inputs, Eigen::VectorXd const &targ) override
virtual void JacobianImpl(unsigned int const outwrt, unsigned int const inwrt, ref_vector< Eigen::VectorXd > const &inputs) override
CombineVectors(Eigen::VectorXi const &inputSizes)
virtual void EvaluateImpl(ref_vector< Eigen::VectorXd > const &inputs) override
Provides an abstract interface for defining vector-valued model components.
Definition: ModPiece.h:148
Eigen::MatrixXd jacobian
Definition: ModPiece.h:506
const Eigen::VectorXi inputSizes
Definition: ModPiece.h:469
std::vector< Eigen::VectorXd > outputs
Definition: ModPiece.h:503
Eigen::VectorXd gradient
Definition: ModPiece.h:504
const Eigen::VectorXi outputSizes
Definition: ModPiece.h:472
Eigen::VectorXd jacobianAction
Definition: ModPiece.h:505
std::vector< std::reference_wrapper< const T > > ref_vector
A vector of references to something ...
Definition: WorkPiece.h:37
auto get(const nlohmann::detail::iteration_proxy_value< IteratorType > &i) -> decltype(i.key())
Definition: json.h:3956