MUQ  0.4.3
SplitVector.cpp
Go to the documentation of this file.
2 
3 using namespace muq::Modeling;
4 
5 SplitVector::SplitVector(std::vector<int> const& ind,
6  std::vector<int> const& size,
7  unsigned int const insize) : SplitVector(Eigen::Map<const Eigen::VectorXi>(&ind[0],ind.size()),
8  Eigen::Map<const Eigen::VectorXi>(&size[0],size.size()),
9  insize){};
10 
11 SplitVector::SplitVector(Eigen::VectorXi const& ind,
12  Eigen::VectorXi const& size,
13  unsigned int const insize) : ModPiece(Eigen::VectorXi::Constant(1, insize), size),
14  ind(ind),
15  size(size) {
16  assert(ind.size()==size.size());
17  assert(size.sum()<=insize);
18  assert(ind.maxCoeff()<insize);
19 }
20 
22  const Eigen::VectorXd& in = inputs[0];
23 
24  outputs.resize(ind.size());
25  for( unsigned int i=0; i<ind.size(); ++i ) {
26  outputs[i] = in.segment(ind(i), size(i)).eval();
27  }
28 }
29 
30 void SplitVector::JacobianImpl(unsigned int const outwrt, unsigned int const inwrt, ref_vector<Eigen::VectorXd> const& inputs) {
31  assert(inwrt==0);
32  jacobian = Eigen::MatrixXd::Zero(size(outwrt), inputSizes(0));
33  jacobian.block(0, ind(outwrt), size(outwrt), size(outwrt)) += Eigen::MatrixXd::Identity(size(outwrt), size(outwrt)).eval();
34 }
35 
36 void SplitVector::GradientImpl(unsigned int const outwrt, unsigned int const inwrt, ref_vector<Eigen::VectorXd> const& inputs, Eigen::VectorXd const& sens) {
37  assert(inwrt==0);
38  assert(sens.size()==size(outwrt));
39  gradient = Eigen::VectorXd::Zero(inputSizes(0));
40  gradient.segment(ind(outwrt), size(outwrt)) += sens;
41 }
42 
43 void SplitVector::ApplyJacobianImpl(unsigned int const outwrt, unsigned int const inwrt, ref_vector<Eigen::VectorXd> const& inputs, Eigen::VectorXd const& targ) {
44  assert(inwrt==0);
45  assert(targ.size()==inputSizes(0));
46  jacobianAction = targ.segment(ind(outwrt), size(outwrt));
47 }
48 
49 void SplitVector::ApplyHessianImpl(unsigned int const outwrt, unsigned int const inwrt1,unsigned int const inwrt2, ref_vector<Eigen::VectorXd> const& inputs, Eigen::VectorXd const& sens, Eigen::VectorXd const& targ) {
50 
51  assert(inwrt1==0);
52  assert(inwrt2<2);
53 
54  hessAction = Eigen::VectorXd::Zero(inputSizes(0));
55 
56  // If hessian wrt sensitivity...
57  if(inwrt2==1){
58  hessAction.segment(ind(outwrt), size(outwrt)) = Eigen::VectorXd::Ones(size(outwrt));
59  }
60 }
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
Eigen::VectorXd hessAction
Definition: ModPiece.h:507
Eigen::VectorXd jacobianAction
Definition: ModPiece.h:505
Provides a mechanism for splitting a vector into multiple pieces.
Definition: SplitVector.h:42
SplitVector(std::vector< int > const &ind, std::vector< int > const &size, unsigned int const insize)
Definition: SplitVector.cpp:5
virtual void JacobianImpl(unsigned int const outwrt, unsigned int const inwrt, ref_vector< Eigen::VectorXd > const &inputs) override
Definition: SplitVector.cpp:30
const Eigen::VectorXi size
Definition: SplitVector.h:80
const Eigen::VectorXi ind
Definition: SplitVector.h:78
virtual void EvaluateImpl(ref_vector< Eigen::VectorXd > const &inputs) override
Definition: SplitVector.cpp:21
virtual void GradientImpl(unsigned int const outwrt, unsigned int const inwrt, ref_vector< Eigen::VectorXd > const &inputs, Eigen::VectorXd const &sens) override
Definition: SplitVector.cpp:36
virtual void ApplyJacobianImpl(unsigned int const outwrt, unsigned int const inwrt, ref_vector< Eigen::VectorXd > const &inputs, Eigen::VectorXd const &targ) override
Definition: SplitVector.cpp:43
virtual void ApplyHessianImpl(unsigned int const outwrt, unsigned int const inwrt1, unsigned int const inwrt2, ref_vector< Eigen::VectorXd > const &inputs, Eigen::VectorXd const &sens, Eigen::VectorXd const &targ) override
Definition: SplitVector.cpp:49
std::vector< std::reference_wrapper< const T > > ref_vector
A vector of references to something ...
Definition: WorkPiece.h:37