11 #ifndef OOPS_ASSIMILATION_SADDLEPOINTLMPMATRIX_H_
12 #define OOPS_ASSIMILATION_SADDLEPOINTLMPMATRIX_H_
14 #include <Eigen/Dense>
19 #include <boost/noncopyable.hpp>
54 template<
typename MODEL,
typename OBS>
65 void setup(
const std::vector<SPVector_> &,
const std::vector<SPVector_> &);
78 void Fmultiply(Eigen::VectorXd &, Eigen::VectorXd &)
const;
94 template<
typename MODEL,
typename OBS>
102 template<
typename MODEL,
typename OBS>
104 const std::vector<SPVector_> & pqVEC) {
109 for (
unsigned ii = 0; ii < xyVEC.size(); ++ii) {
110 xyVEC_.push_back(xyVEC[ii]);
111 pqVEC_.push_back(pqVEC[ii]);
114 nvec_ = xyVEC_.size();
117 spvecinit_.reset(
new SPVector_(xyVEC_[0]));
118 (*spvecinit_).zero();
124 for (
unsigned jv = 0; jv < nvec_; ++jv) {
125 this->Pinitmultiply(xyVEC_[jv], ww);
128 RpVEC_.push_back(yy.
lambda());
129 RqVEC_.push_back(yy.
dx());
133 ZMat_.resize(nvec_, nvec_);
135 for (
unsigned jj = 0; jj < nvec_; ++jj) {
136 for (
unsigned ii = 0; ii < nvec_; ++ii) {
137 ZMat_(ii, jj) = calpha * dot_product(RpVEC_[ii], xyVEC_[jj].lambda());
144 Ik = Eigen::MatrixXd::Identity(2*nvec_, 2*nvec_);
147 FMat_.resize(2*nvec_, 2*nvec_);
148 Eigen::VectorXd output(2*nvec_);
149 Eigen::VectorXd input(2*nvec_);
150 for (
int ii = 0; ii < 2*nvec_; ii++) {
152 Fmultiply(input, output);
153 FMat_.col(ii) = output;
159 template<
typename MODEL,
typename OBS>
185 j_.jb().multiplyB(x.
lambda().
dx(), diag.
dx());
186 for (
unsigned jj = 0; jj < j_.nterms(); ++jj) {
211 template<
typename MODEL,
typename OBS>
222 Eigen::VectorXd rhs(2*nvec_);
223 Eigen::VectorXd finvx(2*nvec_);
224 this->Gmultiply(z, rhs);
225 finvx = FMat_.colPivHouseholderQr().solve(rhs);
226 this->Rmultiply(finvx, svw);
233 template<
typename MODEL,
typename OBS>
235 Eigen::VectorXd & z)
const {
236 Eigen::VectorXd v1(nvec_);
237 Eigen::VectorXd v2(nvec_);
238 Eigen::VectorXd ax(nvec_);
239 Eigen::VectorXd axt(nvec_);
240 Eigen::MatrixXd ZMatt(nvec_, nvec_);
241 ZMatt = ZMat_.transpose();
244 for (
unsigned jj = 0; jj < nvec_; ++jj) {
245 v1(jj) = dot_product(RpVEC_[jj], x.
lambda());
246 v2(jj) = dot_product(RqVEC_[jj], x.
dx());
249 ax = ZMat_.colPivHouseholderQr().solve(v1);
250 axt = ZMatt.colPivHouseholderQr().solve(v2);
255 template<
typename MODEL,
typename OBS>
262 for (
unsigned ii = 0; ii < nvec_; ++ii) {
263 ww.
axpy(x(ii + nvec_), RpVEC_[ii]);
264 zz.
axpy(x(ii), RqVEC_[ii]);
271 template<
typename MODEL,
typename OBS>
273 Eigen::VectorXd & fx)
const {
277 Eigen::VectorXd ww(2*nvec_);
280 this->Rmultiply(x, tmp);
282 this->Gmultiply(tmp2, ww);
void zero()
Linear algebra operators.
void axpy(const double, const ControlIncrement &)
Weak Constraint 4D-Var Cost Function.
Container of dual space vectors for all terms of the cost function.
void append(std::unique_ptr< GeneralizedDepartures > &&)
std::shared_ptr< const GeneralizedDepartures > getv(const unsigned) const
void axpy(const double, const DualVector &)
The preconditioner for the saddle-point minimizer.
CostFunction< MODEL, OBS > CostFct_
std::vector< SPVector_ > xyVEC_
void Fmultiply(Eigen::VectorXd &, Eigen::VectorXd &) const
std::vector< SPVector_ > pqVEC_
ControlIncrement< MODEL, OBS > CtrlInc_
std::vector< CtrlInc_ > RqVEC_
std::vector< LagVector_ > RpVEC_
SaddlePointVector< MODEL, OBS > SPVector_
DualVector< MODEL, OBS > LagVector_
void multiply(const SPVector_ &, SPVector_ &) const
SaddlePointLMPMatrix(const CostFct_ &j)
void Pinitmultiply(const SPVector_ &, SPVector_ &) const
void setup(const std::vector< SPVector_ > &, const std::vector< SPVector_ > &)
void Rmultiply(Eigen::VectorXd &, SPVector_ &) const
void Gmultiply(const SPVector_ &, Eigen::VectorXd &) const
CostFctWeak< MODEL, OBS > CostFctWeak_
std::unique_ptr< SPVector_ > spvecinit_
The preconditioner for the saddle-point minimizer.
void multiply(const SPVector_ &, SPVector_ &) const
Control vector for the saddle point formulation.
const Multipliers_ & lambda() const
Accessor method to get the lambda_ component.
const CtrlInc_ & dx() const
Accessor method to get the dx_ component.
The namespace for the main oops code.