MUQ  0.4.3
GradientPiece.cpp
Go to the documentation of this file.
2 
3 using namespace muq::Modeling;
4 
5 
6 
7 
8 GradientPiece::GradientPiece(std::shared_ptr<ModPiece> const& basePieceIn,
9  unsigned int const outWrtIn,
10  unsigned int const inWrtIn) : ModPiece(GetInputSizes(basePieceIn,outWrtIn),
11  GetOutputSizes(basePieceIn, inWrtIn)),
12  basePiece(basePieceIn),
13  outWrt(outWrtIn),
14  inWrt(inWrtIn)
15 {
16 }
17 
18 
20 {
21  ref_vector<Eigen::VectorXd> baseInputs(input.begin(), input.end()-1);
22 
23  outputs.resize(1);
24  outputs.at(0) = basePiece->Gradient(outWrt, inWrt, baseInputs, input.at(input.size()-1).get());
25 }
26 
27 void GradientPiece::ApplyJacobianImpl(unsigned int const outputDimWrt,
28  unsigned int const inputDimWrt,
29  ref_vector<Eigen::VectorXd> const& input,
30  Eigen::VectorXd const& vec)
31 {
32  assert(outputDimWrt==0);
33 
34  ref_vector<Eigen::VectorXd> baseInputs(input.begin(), input.end()-1);
35  jacobianAction = basePiece->ApplyHessian(outWrt, inWrt, inputDimWrt, baseInputs, input.at(input.size()-1).get(), vec);
36 }
37 
38 
39 Eigen::VectorXi GradientPiece::GetInputSizes(std::shared_ptr<ModPiece> const& basePiece,
40  unsigned int const outWrt)
41 {
42  unsigned int numInputs = basePiece->inputSizes.size()+1;
43 
44  Eigen::VectorXi inputSizes(numInputs);
45  inputSizes.head(numInputs-1) = basePiece->inputSizes;
46  inputSizes(numInputs-1) = basePiece->outputSizes(outWrt); // extra input for adjoint variable
47 
48  return inputSizes;
49 }
50 
51 Eigen::VectorXi GradientPiece::GetOutputSizes(std::shared_ptr<ModPiece> const& basePiece,
52  unsigned int const inWrt)
53 {
54  Eigen::VectorXi outputSizes(1);
55  outputSizes(0) = basePiece->inputSizes(inWrt);
56  return outputSizes;
57 }
static Eigen::VectorXi GetOutputSizes(std::shared_ptr< ModPiece > const &basePiece, unsigned int const inWrt)
virtual void EvaluateImpl(ref_vector< Eigen::VectorXd > const &input) override
static Eigen::VectorXi GetInputSizes(std::shared_ptr< ModPiece > const &basePiece, unsigned int const outWrt)
virtual void ApplyJacobianImpl(unsigned int const outputDimWrt, unsigned int const inputDimWrt, ref_vector< Eigen::VectorXd > const &input, Eigen::VectorXd const &vec) override
GradientPiece(std::shared_ptr< ModPiece > const &basePieceIn, unsigned int const outWrt, unsigned int const inWrt)
const unsigned int inWrt
Definition: GradientPiece.h:37
const unsigned int outWrt
Definition: GradientPiece.h:37
std::shared_ptr< ModPiece > basePiece
Definition: GradientPiece.h:35
Provides an abstract interface for defining vector-valued model components.
Definition: ModPiece.h:148
const Eigen::VectorXi inputSizes
Definition: ModPiece.h:469
std::vector< Eigen::VectorXd > outputs
Definition: ModPiece.h:503
const Eigen::VectorXi outputSizes
Definition: ModPiece.h:472
Eigen::VectorXd jacobianAction
Definition: ModPiece.h:505
int numInputs
The number of inputs.
Definition: WorkPiece.h:501
std::vector< std::reference_wrapper< const T > > ref_vector
A vector of references to something ...
Definition: WorkPiece.h:37