19#ifndef ROL_POISSONINVERSION_HPP
20#define ROL_POISSONINVERSION_HPP
26#include "ROL_LAPACK.hpp"
27#include "ROL_LinearAlgebra.hpp"
41 typedef typename vector::size_type
uint;
79 for (
uint i = 0; i <
nz_; i++) {
81 val +=
alpha_/2.0 *
hz_ * (*zp)[i]*(*zp)[i];
88 val +=
alpha_ * std::sqrt(std::pow((*zp)[i]-(*zp)[i+1],2.0)+
eps_);
103 ROL::Ptr<const vector> zp =
getVector(z);
106 for (
uint i = 0; i <
nz_; i++) {
107 (*gp)[i] =
alpha_ *
hz_ * (*zp)[i]/std::sqrt(std::pow((*zp)[i],2.0)+
eps_);
111 ROL::Ptr<const vector> zp =
getVector(z);
115 for (
uint i = 0; i <
nz_; i++) {
117 diff = (*zp)[i]-(*zp)[i+1];
118 (*gp)[i] =
alpha_ * diff/std::sqrt(std::pow(diff,2.0)+
eps_);
120 else if ( i ==
nz_-1 ) {
121 diff = (*zp)[i-1]-(*zp)[i];
122 (*gp)[i] = -
alpha_ * diff/std::sqrt(std::pow(diff,2.0)+
eps_);
125 diff = (*zp)[i]-(*zp)[i+1];
126 (*gp)[i] =
alpha_ * diff/std::sqrt(std::pow(diff,2.0)+
eps_);
127 diff = (*zp)[i-1]-(*zp)[i];
128 (*gp)[i] -=
alpha_ * diff/std::sqrt(std::pow(diff,2.0)+
eps_);
143 ROL::Ptr<const vector> zp =
getVector(z);
144 ROL::Ptr<const vector> vp =
getVector(v);
147 for (
uint i = 0; i <
nz_; i++) {
148 (*hvp)[i] =
alpha_*
hz_*(*vp)[i]*
eps_/std::pow(std::pow((*zp)[i],2.0)+
eps_,3.0/2.0);
152 ROL::Ptr<const vector> zp =
getVector(z);
153 ROL::Ptr<const vector> vp =
getVector(v);
158 for (
uint i = 0; i <
nz_; i++) {
160 diff1 = (*zp)[i]-(*zp)[i+1];
161 diff1 =
eps_/std::pow(std::pow(diff1,2.0)+
eps_,3.0/2.0);
162 (*hvp)[i] =
alpha_* ((*vp)[i]*diff1 - (*vp)[i+1]*diff1);
164 else if ( i ==
nz_-1 ) {
165 diff2 = (*zp)[i-1]-(*zp)[i];
166 diff2 =
eps_/std::pow(std::pow(diff2,2.0)+
eps_,3.0/2.0);
167 (*hvp)[i] =
alpha_* (-(*vp)[i-1]*diff2 + (*vp)[i]*diff2);
170 diff1 = (*zp)[i]-(*zp)[i+1];
171 diff1 =
eps_/std::pow(std::pow(diff1,2.0)+
eps_,3.0/2.0);
172 diff2 = (*zp)[i-1]-(*zp)[i];
173 diff2 =
eps_/std::pow(std::pow(diff2,2.0)+
eps_,3.0/2.0);
174 (*hvp)[i] =
alpha_* (-(*vp)[i-1]*diff2 + (*vp)[i]*(diff1 + diff2) - (*vp)[i+1]*diff1);
184 ROL::Ptr<const vector> fp =
getVector(f);
187 for (
uint i = 0; i <
nu_; i++) {
189 (*Mfp)[i] =
hu_/6.0*(4.0*(*fp)[i] + (*fp)[i+1]);
191 else if ( i ==
nu_-1 ) {
192 (*Mfp)[i] =
hu_/6.0*((*fp)[i-1] + 4.0*(*fp)[i]);
195 (*Mfp)[i] =
hu_/6.0*((*fp)[i-1] + 4.0*(*fp)[i] + (*fp)[i+1]);
203 ROL::Ptr<const vector> zp =
getVector(z);
210 for (
uint i = 0; i <
nu_; i++ ) {
211 d[i] = (std::exp((*zp)[i]) + std::exp((*zp)[i+1]))/
hu_;
213 o[i] *= -std::exp((*zp)[i+1])/
hu_;
218 ROL::LAPACK<int,Real> lp;
222 lp.PTTRF(
nu_,&d[0],&o[0],&info);
223 lp.PTTRS(
nu_,nhrs,&d[0],&o[0],&(*bp)[0],ldb,&info);
236 ROL::Ptr<const vector> zp =
getVector(z);
237 ROL::Ptr<const vector> up =
getVector(u);
238 ROL::Ptr<const vector> dp =
getVector(d);
241 for (
uint i = 0; i <
nu_; i++) {
243 (*Bdp)[i] = 1.0/
hu_*( std::exp((*zp)[i])*(*up)[i]*(*dp)[i]
244 + std::exp((*zp)[i+1])*((*up)[i]-(*up)[i+1])*(*dp)[i+1] );
246 else if ( i ==
nu_-1 ) {
247 (*Bdp)[i] = 1.0/
hu_*( std::exp((*zp)[i])*((*up)[i]-(*up)[i-1])*(*dp)[i]
248 + std::exp((*zp)[i+1])*(*up)[i]*(*dp)[i+1] );
251 (*Bdp)[i] = 1.0/
hu_*( std::exp((*zp)[i])*((*up)[i]-(*up)[i-1])*(*dp)[i]
252 + std::exp((*zp)[i+1])*((*up)[i]-(*up)[i+1])*(*dp)[i+1] );
261 ROL::Ptr<const vector> zp =
getVector(z);
262 ROL::Ptr<const vector> up =
getVector(u);
263 ROL::Ptr<const vector> dp =
getVector(d);
266 for (
uint i = 0; i <
nz_; i++) {
268 (*Bdp)[i] = std::exp((*zp)[i])/
hu_*(*up)[i]*(*dp)[i];
270 else if ( i ==
nz_-1 ) {
271 (*Bdp)[i] = std::exp((*zp)[i])/
hu_*(*up)[i-1]*(*dp)[i-1];
274 (*Bdp)[i] = std::exp((*zp)[i])/
hu_*( ((*up)[i]-(*up)[i-1])*((*dp)[i]-(*dp)[i-1]) );
282 ROL::Ptr<const vector> zp =
getVector(z);
283 ROL::Ptr<const vector> vp =
getVector(v);
284 ROL::Ptr<const vector> up =
getVector(u);
285 ROL::Ptr<const vector> dp =
getVector(d);
288 for (
uint i = 0; i <
nz_; i++) {
290 (*Bdp)[i] = (*vp)[i]*std::exp((*zp)[i])/
hu_*(*up)[i]*(*dp)[i];
292 else if ( i ==
nz_-1 ) {
293 (*Bdp)[i] = (*vp)[i]*std::exp((*zp)[i])/
hu_*(*up)[i-1]*(*dp)[i-1];
296 (*Bdp)[i] = (*vp)[i]*std::exp((*zp)[i])/
hu_*( ((*up)[i]-(*up)[i-1])*((*dp)[i]-(*dp)[i-1]) );
310 ROL::Ptr<vector> bp = ROL::makePtr<vector>(
nu_, 0.0);
311 for (
uint i = 0; i <
nu_; i++ ) {
312 if ( (Real)(i+1)*
hu_ < 0.5 ) {
313 (*bp)[i] = 2.0*k1*
hu_;
315 else if ( std::abs((Real)(i+1)*
hu_ - 0.5) < ROL_EPSILON<Real>() ) {
316 (*bp)[i] = (k1+k2)*
hu_;
318 else if ( (Real)(i+1)*
hu_ > 0.5 ) {
319 (*bp)[i] = 2.0*k2*
hu_;
333 ROL::Ptr<const vector> up =
getVector(u);
334 ROL::Ptr<vector> rp = ROL::makePtr<vector>(
nu_,0.0);
337 for (
uint i = 0; i <
nu_; i++) {
349 SV b( ROL::makePtr<vector>(
nu_,0.0) );
359 SV res( ROL::makePtr<vector>(
nu_,0.0) );
361 SV res1( ROL::makePtr<vector>(
nu_,0.0) );
374 ROL::Ptr<vector> up = ROL::makePtr<vector>(
nu_,0.0);
380 ROL::Ptr<vector> rp = ROL::makePtr<vector>(
nu_,0.0);
383 for (
uint i = 0; i <
nu_; i++) {
387 ROL::Ptr<V> Mres = res.
clone();
389 return 0.5*Mres->dot(res) +
reg_value(z);
397 SV u( ROL::makePtr<vector>(
nu_,0.0) );
401 SV p( ROL::makePtr<std::vector<Real>>(
nu_,0.0) );
408 SV g_reg( ROL::makePtr<vector>(
nz_,0.0) );
420 SV u( ROL::makePtr<vector>(
nu_,0.0) );
424 SV p( ROL::makePtr<vector>(
nu_,0.0) );
428 SV w( ROL::makePtr<vector>(
nu_,0.0) );
432 SV q( ROL::makePtr<vector>(
nu_,0.0) );
439 SV tmp( ROL::makePtr<vector>(
nz_,0.0) );
449 SV hv_reg( ROL::makePtr<vector>(
nz_,0.0) );
467 int dim =
static_cast<int>(vp.size());
471 ROL::LA::Matrix<Real> H(
dim,
dim);
473 H = computeDenseHessian<Real>(obj, x);
476 std::vector<vector> eigenvals = computeEigenvalues<Real>(H);
477 std::sort((eigenvals[0]).begin(), (eigenvals[0]).end());
480 Real inertia = (eigenvals[0])[0];
481 Real correction = 0.0;
482 if ( inertia <= 0.0 ) {
483 correction = 2.0*std::abs(inertia);
484 if ( inertia == 0.0 ) {
486 while ( eigenvals[0][cnt] == 0.0 ) {
489 correction = 0.5*eigenvals[0][cnt];
490 if ( cnt ==
dim-1 ) {
494 for (
int i=0; i<
dim; i++) {
495 H(i,i) += correction;
500 ROL::LA::Matrix<Real> invH = computeInverse<Real>(H);
503 ROL::LA::Vector<Real> hv_la(ROL::LA::View, &((*hvp)[0]),
dim);
504 ROL::LA::Vector<Real> v_la(ROL::LA::View, &(vp[0]),
dim);
507 for (
int i = 0; i <
dim; i++) {
509 for (
int j = 0; j <
dim; j++) {
510 hv_la(i) += invH(i,j) * v_la(j);
526 return ROL::makePtr<Objective_PoissonInversion<Real>>(n,1.e-6);
533 ROL::Ptr<std::vector<Real> > x0p = ROL::makePtr<std::vector<Real>>(n,0.0);
534 for (
int i = 0; i < n; i++ ) {
537 return ROL::makePtr<StdVector<Real>>(x0p);
544 ROL::Ptr<std::vector<Real> > xp = ROL::makePtr<std::vector<Real>>(n,0.0);
545 Real h = 1.0/((Real)n+1), pt = 0.0, k1 = 1.0, k2 = 2.0;
546 for(
int li = 0; li < n; li++ ) {
548 if ( pt >= 0.0 && pt < 0.5 ) {
549 (*xp)[li] = std::log(k1);
551 else if ( pt >= 0.5 && pt < 1.0 ) {
552 (*xp)[li] = std::log(k2);
555 return ROL::makePtr<StdVector<Real>>(xp);
Contains definitions for helper functions in ROL.
Contains definitions of test objective functions.
Provides the interface to evaluate objective functions.
virtual void hessVec(Vector< Real > &hv, const Vector< Real > &v, const Vector< Real > &x, Real &tol)
Apply Hessian approximation to vector.
Provides the ROL::Vector interface for scalar values, to be used, for example, with scalar constraint...
void axpy(const Real alpha, const Vector< Real > &x)
Compute where .
virtual Ptr< Vector< Real > > clone() const
Clone to make a new (uninitialized) vector.
Defines the linear algebra or vector space interface.
virtual void set(const Vector &x)
Set where .
virtual void scale(const Real alpha)=0
Compute where .
virtual void plus(const Vector &x)=0
Compute , where .
virtual void axpy(const Real alpha, const Vector &x)
Compute where .
Poisson material inversion.
void solve_adjoint_equation(Vector< Real > &p, const Vector< Real > &u, const Vector< Real > &z)
void solve_poisson(Vector< Real > &u, const Vector< Real > &z, Vector< Real > &b)
std::vector< Real > vector
void solve_adjoint_sensitivity_equation(Vector< Real > &q, const Vector< Real > &w, const Vector< Real > &v, const Vector< Real > &p, const Vector< Real > &u, const Vector< Real > &z)
Real reg_value(const Vector< Real > &z)
void apply_mass(Vector< Real > &Mf, const Vector< Real > &f)
void reg_hessVec(Vector< Real > &hv, const Vector< Real > &v, const Vector< Real > &z)
ROL::Ptr< const vector > getVector(const V &x)
void apply_transposed_linearized_control_operator(Vector< Real > &Bd, const Vector< Real > &z, const Vector< Real > &d, const Vector< Real > &u)
void solve_state_sensitivity_equation(Vector< Real > &w, const Vector< Real > &v, const Vector< Real > &u, const Vector< Real > &z)
void gradient(Vector< Real > &g, const Vector< Real > &z, Real &tol)
Compute gradient.
void solve_state_equation(Vector< Real > &u, const Vector< Real > &z)
void apply_linearized_control_operator(Vector< Real > &Bd, const Vector< Real > &z, const Vector< Real > &d, const Vector< Real > &u)
Objective_PoissonInversion(uint nz=32, Real alpha=1.e-4)
ROL::Ptr< vector > getVector(V &x)
Real value(const Vector< Real > &z, Real &tol)
Compute value.
void reg_gradient(Vector< Real > &g, const Vector< Real > &z)
void apply_transposed_linearized_control_operator_2(Vector< Real > &Bd, const Vector< Real > &z, const Vector< Real > &v, const Vector< Real > &d, const Vector< Real > &u)
void invHessVec(Vector< Real > &hv, const Vector< Real > &v, const Vector< Real > &x, Real &tol)
Apply inverse Hessian approximation to vector.
Real evaluate_target(Real x)
Ptr< Vector< Real > > getInitialGuess(void) const
getPoissonInversion(void)
Ptr< Objective< Real > > getObjective(void) const
Ptr< Vector< Real > > getSolution(const int i=0) const