53 #ifndef ROL_POISSONINVERSION_HPP
54 #define ROL_POISSONINVERSION_HPP
60 #include "Teuchos_LAPACK.hpp"
98 for (
int i = 0; i < this->
nz_; i++) {
100 val += this->
alpha_/2.0 * this->
hz_ * (*zp)[i]*(*zp)[i];
103 val += this->
alpha_ * this->
hz_ * std::sqrt((*zp)[i]*(*zp)[i] + this->
eps_);
106 if ( i < this->nz_-1 ) {
107 val += this->
alpha_ * std::sqrt(std::pow((*zp)[i]-(*zp)[i+1],2.0)+this->
eps_);
126 for (
int i = 0; i < this->
nz_; i++) {
127 (*gp)[i] = this->
alpha_ * this->
hz_ * (*zp)[i]/std::sqrt(std::pow((*zp)[i],2.0)+this->
eps_);
138 for (
int i = 0; i < this->
nz_; i++) {
140 diff = (*zp)[i]-(*zp)[i+1];
141 (*gp)[i] = this->
alpha_ * diff/std::sqrt(std::pow(diff,2.0)+this->
eps_);
143 else if ( i == this->nz_-1 ) {
144 diff = (*zp)[i-1]-(*zp)[i];
145 (*gp)[i] = -this->
alpha_ * diff/std::sqrt(std::pow(diff,2.0)+this->
eps_);
148 diff = (*zp)[i]-(*zp)[i+1];
149 (*gp)[i] = this->
alpha_ * diff/std::sqrt(std::pow(diff,2.0)+this->
eps_);
150 diff = (*zp)[i-1]-(*zp)[i];
151 (*gp)[i] -= this->
alpha_ * diff/std::sqrt(std::pow(diff,2.0)+this->
eps_);
172 for (
int i = 0; i < this->
nz_; i++) {
173 (*hvp)[i] = this->
alpha_*this->
hz_*(*vp)[i]*this->
eps_/std::pow(std::pow((*zp)[i],2.0)+this->
eps_,3.0/2.0);
188 for (
int i = 0; i < this->
nz_; i++) {
190 diff1 = (*zp)[i]-(*zp)[i+1];
191 diff1 = this->
eps_/std::pow(std::pow(diff1,2.0)+this->
eps_,3.0/2.0);
192 (*hvp)[i] = this->
alpha_* ((*vp)[i]*diff1 - (*vp)[i+1]*diff1);
194 else if ( i == this->nz_-1 ) {
195 diff2 = (*zp)[i-1]-(*zp)[i];
196 diff2 = this->
eps_/std::pow(std::pow(diff2,2.0)+this->
eps_,3.0/2.0);
197 (*hvp)[i] = this->
alpha_* (-(*vp)[i-1]*diff2 + (*vp)[i]*diff2);
200 diff1 = (*zp)[i]-(*zp)[i+1];
201 diff1 = this->
eps_/std::pow(std::pow(diff1,2.0)+this->
eps_,3.0/2.0);
202 diff2 = (*zp)[i-1]-(*zp)[i];
203 diff2 = this->
eps_/std::pow(std::pow(diff2,2.0)+this->
eps_,3.0/2.0);
204 (*hvp)[i] = this->
alpha_* (-(*vp)[i-1]*diff2 + (*vp)[i]*(diff1 + diff2) - (*vp)[i+1]*diff1);
219 for (
int i = 0; i < this->
nu_; i++) {
221 (*Mfp)[i] = this->
hu_/6.0*(4.0*(*fp)[i] + (*fp)[i+1]);
223 else if ( i == this->nu_-1 ) {
224 (*Mfp)[i] = this->
hu_/6.0*((*fp)[i-1] + 4.0*(*fp)[i]);
227 (*Mfp)[i] = this->
hu_/6.0*((*fp)[i-1] + 4.0*(*fp)[i] + (*fp)[i+1]);
244 std::vector<Real> d(this->
nu_,1.0);
245 std::vector<Real> o(this->
nu_-1,1.0);
246 for (
int i = 0; i < this->
nu_; i++ ) {
247 d[i] = (std::exp((*zp)[i]) + std::exp((*zp)[i+1]))/this->
hu_;
248 if ( i < this->nu_-1 ) {
249 o[i] *= -std::exp((*zp)[i+1])/this->
hu_;
254 Teuchos::LAPACK<int,Real> lp;
258 lp.PTTRF(this->nu_,&d[0],&o[0],&info);
259 lp.PTTRS(this->nu_,nhrs,&d[0],&o[0],&(*bp)[0],ldb,&info);
282 for (
int i = 0; i < this->
nu_; i++) {
284 (*Bdp)[i] = 1.0/this->
hu_*( std::exp((*zp)[i])*(*up)[i]*(*dp)[i]
285 + std::exp((*zp)[i+1])*((*up)[i]-(*up)[i+1])*(*dp)[i+1] );
287 else if ( i == this->nu_-1 ) {
288 (*Bdp)[i] = 1.0/this->
hu_*( std::exp((*zp)[i])*((*up)[i]-(*up)[i-1])*(*dp)[i]
289 + std::exp((*zp)[i+1])*(*up)[i]*(*dp)[i+1] );
292 (*Bdp)[i] = 1.0/this->
hu_*( std::exp((*zp)[i])*((*up)[i]-(*up)[i-1])*(*dp)[i]
293 + std::exp((*zp)[i+1])*((*up)[i]-(*up)[i+1])*(*dp)[i+1] );
313 for (
int i = 0; i < this->
nz_; i++) {
315 (*Bdp)[i] = std::exp((*zp)[i])/this->
hu_*(*up)[i]*(*dp)[i];
317 else if ( i == this->nz_-1 ) {
318 (*Bdp)[i] = std::exp((*zp)[i])/this->
hu_*(*up)[i-1]*(*dp)[i-1];
321 (*Bdp)[i] = std::exp((*zp)[i])/this->
hu_*( ((*up)[i]-(*up)[i-1])*((*dp)[i]-(*dp)[i-1]) );
344 for (
int i = 0; i < this->
nz_; i++) {
346 (*Bdp)[i] = (*vp)[i]*std::exp((*zp)[i])/this->
hu_*(*up)[i]*(*dp)[i];
348 else if ( i == this->nz_-1 ) {
349 (*Bdp)[i] = (*vp)[i]*std::exp((*zp)[i])/this->
hu_*(*up)[i-1]*(*dp)[i-1];
352 (*Bdp)[i] = (*vp)[i]*std::exp((*zp)[i])/this->
hu_*( ((*up)[i]-(*up)[i-1])*((*dp)[i]-(*dp)[i-1]) );
362 Teuchos::RCP<std::vector<Real> > bp = Teuchos::rcp(
new std::vector<Real> (this->
nu_, 0.0) );
363 for (
int i = 0; i < this->
nu_; i++ ) {
364 if ( (Real)(i+1)*this->
hu_ < 0.5 ) {
365 (*bp)[i] = 2.0*k1*this->
hu_;
368 (*bp)[i] = (k1+k2)*this->
hu_;
370 else if ( (Real)(i+1)*this->
hu_ > 0.5 ) {
371 (*bp)[i] = 2.0*k2*this->
hu_;
389 for (
int i = 0; i < this->
nu_; i++) {
392 StdVector<Real> Mres( Teuchos::rcp(
new std::vector<Real>(this->nu_,0.0) ) );
428 for (
int i = 0; i < this->
nu_; i++) {
431 StdVector<Real> Mres( Teuchos::rcp(
new std::vector<Real>(this->nu_,0.0) ) );
508 Teuchos::SerialDenseMatrix<int, Real> H(dim, dim);
510 H = computeDenseHessian<Real>(obj, x);
513 std::vector<std::vector<Real> > eigenvals = computeEigenvalues<Real>(H);
514 std::sort((eigenvals[0]).begin(), (eigenvals[0]).end());
517 Real inertia = (eigenvals[0])[0];
518 Real correction = 0.0;
519 if ( inertia <= 0.0 ) {
520 correction = 2.0*std::abs(inertia);
521 if ( inertia == 0.0 ) {
523 while ( eigenvals[0][cnt] == 0.0 ) {
526 correction = 0.5*eigenvals[0][cnt];
527 if ( cnt == dim-1 ) {
531 for (
int i=0; i<dim; i++) {
532 H(i,i) += correction;
537 Teuchos::SerialDenseMatrix<int, Real> invH = computeInverse<Real>(H);
540 Teuchos::SerialDenseVector<int, Real> hv_teuchos(Teuchos::View, &((*hvp)[0]), dim);
541 const Teuchos::SerialDenseVector<int, Real> v_teuchos(Teuchos::View, &(vp[0]), dim);
542 hv_teuchos.multiply(Teuchos::NO_TRANS, Teuchos::NO_TRANS, 1.0, invH, v_teuchos, 0.0);
564 for (
int i=0; i<n; i++) {
568 Real h = 1.0/((Real)n+1);
572 for(
int i=0; i<n; i++ ) {
574 if ( pt >= 0.0 && pt < 0.5 ) {
575 (*xp)[i] = std::log(k1);
577 else if ( pt >= 0.5 && pt < 1.0 ) {
578 (*xp)[i] = std::log(k2);
Provides the interface to evaluate objective functions.
void apply_transposed_linearized_control_operator(Vector< Real > &Bd, const Vector< Real > &z, const Vector< Real > &d, const Vector< Real > &u)
virtual void scale(const Real alpha)=0
Compute where .
void getPoissonInversion(Teuchos::RCP< Objective< Real > > &obj, Vector< Real > &x0, Vector< Real > &x)
void solve_state_sensitivity_equation(Vector< Real > &w, const Vector< Real > &v, const Vector< Real > &u, const Vector< Real > &z)
virtual void plus(const Vector &x)=0
Compute , where .
virtual void axpy(const Real alpha, const Vector &x)
Compute where .
void gradient(Vector< Real > &g, const Vector< Real > &z, Real &tol)
Compute gradient.
Real dot(const Vector< Real > &x) const
Compute where .
virtual void hessVec(Vector< Real > &hv, const Vector< Real > &v, const Vector< Real > &x, Real &tol)
Apply Hessian approximation to vector.
Contains definitions for helper functions in ROL.
Defines the linear algebra or vector space interface.
void solve_poisson(Vector< Real > &u, const Vector< Real > &z, Vector< Real > &b)
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)
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)
Real evaluate_target(Real x)
Provides the std::vector implementation of the ROL::Vector interface.
void invHessVec(Vector< Real > &hv, const Vector< Real > &v, const Vector< Real > &x, Real &tol)
Apply inverse Hessian approximation to vector.
void solve_state_equation(Vector< Real > &u, const Vector< Real > &z)
void reg_hessVec(Vector< Real > &hv, const Vector< Real > &v, const Vector< Real > &z)
Teuchos::RCP< std::vector< Real > > downCast(Vector< Real > &x)
void apply_mass(Vector< Real > &Mf, const Vector< Real > &f)
void solve_adjoint_equation(Vector< Real > &p, const Vector< Real > &u, const Vector< Real > &z)
Objective_PoissonInversion(int nz=32, Real alpha=1.e-4)
void apply_linearized_control_operator(Vector< Real > &Bd, const Vector< Real > &z, const Vector< Real > &d, const Vector< Real > &u)
Teuchos::RCP< const std::vector< Real > > constDownCast(const Vector< Real > &x)
Real value(const Vector< Real > &z, Real &tol)
Compute value.
virtual void set(const Vector &x)
Set where .
Real reg_value(const Vector< Real > &z)
void reg_gradient(Vector< Real > &g, const Vector< Real > &z)
Poisson material inversion.
static const double ROL_EPSILON
Platform-dependent machine epsilon.