39 #ifndef Xyce_N_ANP_ROL_DC_Optimizationh
40 #define Xyce_N_ANP_ROL_DC_Optimizationh
43 #define FD_HESSIAN 0 // 0 to set all hessvec to zero (Gauss-Newton Hessian); else FD Hessian
46 #ifndef IDENTITY_PRECONDITIONER
47 #define IDENTITY_PRECONDITIONER 1
54 #include "ROL_StdVector.hpp"
55 #include "ROL_CVaRVector.hpp"
56 #include "ROL_Vector_SimOpt.hpp"
57 #include "ROL_EqualityConstraint_SimOpt.hpp"
58 #include "ROL_ParametrizedEqualityConstraint_SimOpt.hpp"
59 #include "ROL_Objective_SimOpt.hpp"
60 #include "ROL_ParametrizedObjective_SimOpt.hpp"
70 #include <Xyce_config.h>
76 #include <N_ERH_ErrorMgr.h>
77 #include <N_IO_CircuitBlock.h>
78 #include <N_IO_OptionBlock.h>
79 #include <N_IO_PkgOptionsMgr.h>
80 #include <N_IO_SpiceSeparatedFieldTool.h>
81 #include <N_LAS_Builder.h>
82 #include <N_LAS_Matrix.h>
83 #include <N_LAS_Problem.h>
84 #include <N_LAS_Solver.h>
85 #include <N_LAS_System.h>
86 #include <N_LAS_Vector.h>
92 #include <N_PDS_Comm.h>
93 #include <N_PDS_MPI.h>
94 #include <N_PDS_Manager.h>
95 #include <N_PDS_Serial.h>
97 #include <N_TOP_Topology.h>
98 #include <N_UTL_Algorithm.h>
99 #include <N_UTL_Diagnostic.h>
100 #include <N_UTL_Expression.h>
101 #include <N_UTL_ExtendedString.h>
102 #include <N_UTL_FeatureTest.h>
103 #include <N_UTL_OptionBlock.h>
120 class EqualityConstraint_ROL_DC :
virtual public ::ROL::ParametrizedEqualityConstraint_SimOpt<Real>
124 Nonlinear::NonLinearSolver & nls_;
126 Loader::NonlinearEquationLoader & nEqLoader_;
129 Linear::Vector * origFVectorPtr_;
130 Linear::Vector * pertFVectorPtr_;
131 Linear::Vector * origQVectorPtr_;
132 Linear::Vector * pertQVectorPtr_;
133 Linear::Vector * origBVectorPtr_;
134 Linear::Vector * pertBVectorPtr_;
136 void setControlParams(const ::ROL::Vector<Real> &z)
138 Teuchos::RCP<const std::vector<Real> > zp = (Teuchos::dyn_cast<::ROL::StdVector<Real> >(
const_cast<::ROL::Vector<Real> &
>(z))).getVector();
140 std::string paramName;
142 for (
int i=0;i<rolSweep_.numParams_;i++)
144 paramName = rolSweep_.paramNameVec_[i];
145 paramValue = (*zp)[i];
147 nEqLoader_.setParam(paramName,paramValue);
151 void setUncertainParams()
155 const std::vector<Real> uparam = ::ROL::ParametrizedEqualityConstraint_SimOpt<Real>::getParameter();
156 std::string paramName;
158 for (
int i=0;i<rolSweep_.uncertainParams_.size();i++)
160 paramName = rolSweep_.uncertainParams_[i];
161 paramValue = uparam[i];
163 nEqLoader_.setParam(paramName,paramValue);
172 EqualityConstraint_ROL_DC(
176 AnalysisManager & analysisManager,
177 Loader::NonlinearEquationLoader & loader,
178 Nonlinear::NonLinearSolver & nls,
179 Linear::System & linearSystem,
198 void value(::ROL::Vector<Real> &c, const ::ROL::Vector<Real> &u, const ::ROL::Vector<Real> &z, Real &tol)
200 Teuchos::RCP<std::vector<Teuchos::RCP<Linear::Vector> > > cp = Teuchos::rcp_const_cast<std::vector<Teuchos::RCP<Linear::Vector> > >((Teuchos::dyn_cast<Linear::ROL_XyceVector<Real> >(c)).getVector());
201 Teuchos::RCP<const std::vector<Teuchos::RCP<Linear::Vector> > > up = (Teuchos::dyn_cast< Linear::ROL_XyceVector<Real> >(
const_cast< ::ROL::Vector<Real> &
>(u))).getVector();
204 setUncertainParams();
208 bool vstatus = nEqLoader_.getVoltageLimiterStatus();
209 nEqLoader_.setVoltageLimiterStatus(
false);
211 for (
int i=0;i<nc_;i++)
216 rolSweep_.setSweepValue(i);
225 success = analysisManager_.getDataStore()->setNextSolVectorPtr(*(*up)[i]);
228 nEqLoader_.loadRHS();
231 *((*cp)[i]) = *(linearSystem_.getRHSVector());
234 nEqLoader_.setVoltageLimiterStatus(vstatus);
239 void solve(::ROL::Vector<Real> &u, const ::ROL::Vector<Real> &z, Real &tol)
241 Teuchos::RCP<std::vector<Teuchos::RCP<Linear::Vector> > > up = Teuchos::rcp_const_cast<std::vector<Teuchos::RCP<Linear::Vector> > >((Teuchos::dyn_cast<Linear::ROL_XyceVector<Real> >(u)).getVector());
245 setUncertainParams();
248 success = rolSweep_.doInit() && rolSweep_.doLoopProcess() && rolSweep_.doFinish();
250 for (
int i=0;i<nc_;i++)
252 *((*up)[i]) = *(rolSweep_.solutionPtrVector_[i]);
258 void applyJacobian_1(::ROL::Vector<Real> &jv, const ::ROL::Vector<Real> &v, const ::ROL::Vector<Real> &u, const ::ROL::Vector<Real> &z, Real &tol)
260 Teuchos::RCP<std::vector<Teuchos::RCP<Linear::Vector> > > jvp = Teuchos::rcp_const_cast<std::vector<Teuchos::RCP<Linear::Vector> > >((Teuchos::dyn_cast<Linear::ROL_XyceVector<Real> >(jv)).getVector());
261 Teuchos::RCP<const std::vector<Teuchos::RCP<Linear::Vector> > > vp = (Teuchos::dyn_cast< Linear::ROL_XyceVector<Real> >(
const_cast< ::ROL::Vector<Real> &
>(v))).getVector();
262 Teuchos::RCP<const std::vector<Teuchos::RCP<Linear::Vector> > > up = (Teuchos::dyn_cast< Linear::ROL_XyceVector<Real> >(
const_cast< ::ROL::Vector<Real> &
>(u))).getVector();
266 setUncertainParams();
268 bool vstatus = nEqLoader_.getVoltageLimiterStatus();
269 nEqLoader_.setVoltageLimiterStatus(
false);
271 for (
int i=0;i<nc_;i++)
278 success = analysisManager_.getDataStore()->setNextSolVectorPtr(*(*up)[i]);
282 nEqLoader_.loadRHS();
288 success = nEqLoader_.loadJacobian();
291 Linear::Matrix & Jac = *(linearSystem_.getJacobianMatrix());
295 Jac.matvec(
false, const_cast< const Linear::Vector& >( *((*vp)[i]) ), *((*jvp)[i]) );
302 nEqLoader_.setVoltageLimiterStatus(vstatus);
307 void applyAdjointJacobian_1(::ROL::Vector<Real> &jv, const ::ROL::Vector<Real> &v, const ::ROL::Vector<Real> &u, const ::ROL::Vector<Real> &z, Real &tol)
309 Teuchos::RCP<std::vector<Teuchos::RCP<Linear::Vector> > > jvp = Teuchos::rcp_const_cast<std::vector<Teuchos::RCP<Linear::Vector> > >((Teuchos::dyn_cast<Linear::ROL_XyceVector<Real> >(jv)).getVector());
310 Teuchos::RCP<const std::vector<Teuchos::RCP<Linear::Vector> > > vp = (Teuchos::dyn_cast< Linear::ROL_XyceVector<Real> >(
const_cast< ::ROL::Vector<Real> &
>(v))).getVector();
311 Teuchos::RCP<const std::vector<Teuchos::RCP<Linear::Vector> > > up = (Teuchos::dyn_cast< Linear::ROL_XyceVector<Real> >(
const_cast< ::ROL::Vector<Real> &
>(u))).getVector();
315 setUncertainParams();
317 bool vstatus = nEqLoader_.getVoltageLimiterStatus();
318 nEqLoader_.setVoltageLimiterStatus(
false);
320 for (
int i=0;i<nc_;i++)
327 success = analysisManager_.getDataStore()->setNextSolVectorPtr(*(*up)[i]);
331 nEqLoader_.loadRHS();
334 success = nEqLoader_.loadJacobian();
337 Linear::Matrix & Jac = *(linearSystem_.getJacobianMatrix());
342 Jac.matvec(
true, const_cast< const Linear::Vector& >( *((*vp)[i]) ), *((*jvp)[i]) );
345 nEqLoader_.setVoltageLimiterStatus(vstatus);
350 void applyInverseJacobian_1(::ROL::Vector<Real> &jv, const ::ROL::Vector<Real> &v, const ::ROL::Vector<Real> &u, const ::ROL::Vector<Real> &z, Real &tol)
352 Teuchos::RCP<std::vector<Teuchos::RCP<Linear::Vector> > > jvp = Teuchos::rcp_const_cast<std::vector<Teuchos::RCP<Linear::Vector> > >((Teuchos::dyn_cast<Linear::ROL_XyceVector<Real> >(jv)).getVector());
353 Teuchos::RCP<const std::vector<Teuchos::RCP<Linear::Vector> > > vp = (Teuchos::dyn_cast< Linear::ROL_XyceVector<Real> >(
const_cast< ::ROL::Vector<Real> &
>(v))).getVector();
354 Teuchos::RCP<const std::vector<Teuchos::RCP<Linear::Vector> > > up = (Teuchos::dyn_cast< Linear::ROL_XyceVector<Real> >(
const_cast< ::ROL::Vector<Real> &
>(u))).getVector();
358 setUncertainParams();
360 bool vstatus = nEqLoader_.getVoltageLimiterStatus();
361 nEqLoader_.setVoltageLimiterStatus(
false);
363 for (
int i=0;i<nc_;i++)
370 success = analysisManager_.getDataStore()->setNextSolVectorPtr(*(*up)[i]);
374 nEqLoader_.loadRHS();
377 success = nEqLoader_.loadJacobian();
380 Linear::Vector * savedRHSVectorPtr_ = nls_.rhsVectorPtr_;
381 Linear::Vector * savedNewtonVectorPtr_ = nls_.NewtonVectorPtr_;
387 *(nls_.rhsVectorPtr_) = *(*vp)[i];
388 nls_.rhsVectorPtr_->scale(-1.0);
390 int status = nls_.lasSolverPtr_->solve(
false);
393 std::string msg(
"Sensitivity::solveAdjoint. Solver failed\n");
394 N_ERH_ErrorMgr::report (N_ERH_ErrorMgr::DEV_FATAL, msg);
398 *((*jvp)[i]) = *(nls_.NewtonVectorPtr_);
401 nls_.rhsVectorPtr_->update(1.0, *(savedRHSVectorPtr_),0.0);
402 nls_.NewtonVectorPtr_->update(1.0, *(savedNewtonVectorPtr_),0.0);
405 nEqLoader_.setVoltageLimiterStatus(vstatus);
409 void applyInverseAdjointJacobian_1(::ROL::Vector<Real> &jv, const ::ROL::Vector<Real> &v, const ::ROL::Vector<Real> &u, const ::ROL::Vector<Real> &z, Real &tol)
411 Teuchos::RCP<std::vector<Teuchos::RCP<Linear::Vector> > > jvp = Teuchos::rcp_const_cast<std::vector<Teuchos::RCP<Linear::Vector> > >((Teuchos::dyn_cast<Linear::ROL_XyceVector<Real> >(jv)).getVector());
412 Teuchos::RCP<const std::vector<Teuchos::RCP<Linear::Vector> > > vp = (Teuchos::dyn_cast< Linear::ROL_XyceVector<Real> >(
const_cast< ::ROL::Vector<Real> &
>(v))).getVector();
413 Teuchos::RCP<const std::vector<Teuchos::RCP<Linear::Vector> > > up = (Teuchos::dyn_cast< Linear::ROL_XyceVector<Real> >(
const_cast< ::ROL::Vector<Real> &
>(u))).getVector();
414 Teuchos::RCP<const std::vector<Real> > zp = (Teuchos::dyn_cast<::ROL::StdVector<Real> >(
const_cast<::ROL::Vector<Real> &
>(z))).getVector();
418 setUncertainParams();
420 bool vstatus = nEqLoader_.getVoltageLimiterStatus();
421 nEqLoader_.setVoltageLimiterStatus(
false);
423 for (
int i=0;i<nc_;i++)
430 success = analysisManager_.getDataStore()->setNextSolVectorPtr(*(*up)[i]);
434 nEqLoader_.loadRHS();
437 success = nEqLoader_.loadJacobian();
440 Linear::Vector * savedRHSVectorPtr_ = nls_.rhsVectorPtr_;
441 Linear::Vector * savedNewtonVectorPtr_ = nls_.NewtonVectorPtr_;
447 *(nls_.rhsVectorPtr_) = *(*vp)[i];
448 nls_.rhsVectorPtr_->scale(-1.0);
450 int status = nls_.lasSolverPtr_->solveTranspose(
false);
453 std::string msg(
"Sensitivity::solveAdjoint. Solver failed\n");
454 N_ERH_ErrorMgr::report (N_ERH_ErrorMgr::DEV_FATAL, msg);
458 *((*jvp)[i]) = *(nls_.NewtonVectorPtr_);
461 nls_.rhsVectorPtr_->update(1.0, *(savedRHSVectorPtr_),0.0);
462 nls_.NewtonVectorPtr_->update(1.0, *(savedNewtonVectorPtr_),0.0);
465 nEqLoader_.setVoltageLimiterStatus(
true);
471 void applyJacobian_2(::ROL::Vector<Real> &jv, const ::ROL::Vector<Real> &v, const ::ROL::Vector<Real> &u, const ::ROL::Vector<Real> &z, Real &tol)
473 Teuchos::RCP<std::vector<Teuchos::RCP<Linear::Vector> > > jvp = Teuchos::rcp_const_cast<std::vector<Teuchos::RCP<Linear::Vector> > >((Teuchos::dyn_cast<Linear::ROL_XyceVector<Real> >(jv)).getVector());
474 Teuchos::RCP<const std::vector<Teuchos::RCP<Linear::Vector> > > up = (Teuchos::dyn_cast< Linear::ROL_XyceVector<Real> >(
const_cast< ::ROL::Vector<Real> &
>(u))).getVector();
475 Teuchos::RCP<const std::vector<Real> > vp = (Teuchos::dyn_cast<::ROL::StdVector<Real> >(
const_cast<::ROL::Vector<Real> &
>(v))).getVector();
480 setUncertainParams();
482 bool vstatus = nEqLoader_.getVoltageLimiterStatus();
483 nEqLoader_.setVoltageLimiterStatus(
false);
486 for (
int k=0;k<nc_;k++)
488 rolSweep_.setSweepValue(k);
497 success = analysisManager_.getDataStore()->setNextSolVectorPtr(*(*up)[k]);
500 success = nEqLoader_.loadRHS();
505 analysisManager_.getDataStore()->paramOrigVals_.clear();
510 success = nEqLoader_.loadJacobian ();
515 std::vector<std::string>::iterator firstParam = rolSweep_.paramNameVec_.begin ();
516 std::vector<std::string>::iterator lastParam = rolSweep_.paramNameVec_.end ();
517 std::vector<std::string>::iterator iterParam;
518 for ( iterParam=firstParam, iparam=0;
519 iterParam!=lastParam; ++iterParam, ++iparam )
521 std::string paramName(*iterParam);
525 double paramOrig = 0.0;
526 bool found = nEqLoader_.getParamAndReduce(analysisManager_.getPDSManager()->getPDSComm()->comm(),paramName, paramOrig);
530 std::string msg(
"EqualityConstraint::applyJacobian_2: cannot find parameter ");
532 N_ERH_ErrorMgr::report (N_ERH_ErrorMgr::DEV_FATAL, msg);
535 analysisManager_.getDataStore()->paramOrigVals_.push_back(paramOrig);
538 bool analyticAvailable =
false;
541 analyticAvailable = nEqLoader_.analyticSensitivitiesAvailable (paramName);
543 if (analyticAvailable)
546 std::vector<double> dfdpVec;
547 std::vector<double> dqdpVec;
548 std::vector<double> dbdpVec;
550 std::vector<int> FindicesVec;
551 std::vector<int> QindicesVec;
552 std::vector<int> BindicesVec;
554 nEqLoader_.getAnalyticSensitivities(paramName,dfdpVec,dqdpVec,dbdpVec,FindicesVec, QindicesVec, BindicesVec);
556 rolSweep_.mydfdpPtrVector_[iparam]->putScalar(0.0);
557 rolSweep_.mydqdpPtrVector_[iparam]->putScalar(0.0);
558 rolSweep_.mydbdpPtrVector_[iparam]->putScalar(0.0);
560 int Fsize=FindicesVec.size();
561 for (
int i=0;i<Fsize;++i)
563 Linear::Vector & dfdpRef = *(rolSweep_.mydfdpPtrVector_[iparam]);
564 dfdpRef[FindicesVec[i]] += dfdpVec[i];
567 int Qsize=QindicesVec.size();
568 for (
int i=0;i<Qsize;++i)
570 Linear::Vector & dqdpRef = *(rolSweep_.mydqdpPtrVector_[iparam]);
571 dqdpRef[QindicesVec[i]] += dqdpVec[i];
574 int Bsize=BindicesVec.size();
575 for (
int i=0;i<Bsize;++i)
577 Linear::Vector & dbdpRef = *(rolSweep_.mydbdpPtrVector_[iparam]);
578 dbdpRef[BindicesVec[i]] += dbdpVec[i];
581 rolSweep_.mydfdpPtrVector_[iparam]->fillComplete();
582 rolSweep_.mydqdpPtrVector_[iparam]->fillComplete();
583 rolSweep_.mydbdpPtrVector_[iparam]->fillComplete();
585 if (DEBUG_NONLINEAR && isActive(Diag::SENS_SOLVER))
587 Xyce::dout() << *iterParam <<
": ";
588 Xyce::dout().setf(std::ios::scientific);
589 Xyce::dout() << std::endl;
591 for (
int k1 = 0; k1 < analysisManager_.getDataStore()->solutionSize; ++k1)
594 <<
"dfdp["<<std::setw(3)<<k1<<
"]= "<<std::setw(15)<<std::scientific<<std::setprecision(8)<<(*(rolSweep_.mydfdpPtrVector_[iparam]))[k1]
595 <<
" dqdp["<<std::setw(3)<<k1<<
"]= "<<std::setw(15)<<std::scientific<<std::setprecision(8)<<(*(rolSweep_.mydqdpPtrVector_[iparam]))[k1]
596 <<
" dbdp["<<std::setw(3)<<k1<<
"]= "<<std::setw(15)<<std::scientific<<std::setprecision(8)<<(*(rolSweep_.mydbdpPtrVector_[iparam]))[k1]
605 if (DEBUG_NONLINEAR && isActive(Diag::SENS_SOLVER))
607 Xyce::dout() << std::endl <<
" Calculating numerical df/dp, dq/dp and db/dp for: ";
608 Xyce::dout() << *iterParam << std::endl;
612 origFVectorPtr_ = linearSystem_.builder().createVector();
613 pertFVectorPtr_ = linearSystem_.builder().createVector();
615 origQVectorPtr_ = linearSystem_.builder().createVector();
616 pertQVectorPtr_ = linearSystem_.builder().createVector();
618 origBVectorPtr_ = linearSystem_.builder().createVector();
619 pertBVectorPtr_ = linearSystem_.builder().createVector();
620 origFVectorPtr_->update(1.0, *(analysisManager_.getDataStore()->daeFVectorPtr), 0.0);
621 origQVectorPtr_->update(1.0, *(analysisManager_.getDataStore()->daeQVectorPtr), 0.0);
622 origBVectorPtr_->update(1.0, *(analysisManager_.getDataStore()->daeBVectorPtr), 0.0);
625 double sqrtEta_ = pow(10,
int(log10(fabs(paramOrig))));;
626 double dp = sqrtEta_ * 1.e-5 * (1.0 + fabs(paramOrig));
627 double paramPerturbed = paramOrig;
632 paramPerturbed += dp;
636 paramPerturbed -= dp;
640 static std::string tmp =
"difference=central not supported.\n";
641 N_ERH_ErrorMgr::report( N_ERH_ErrorMgr::USR_FATAL_0, tmp);
645 static std::string tmp =
"difference not recognized!\n";
646 N_ERH_ErrorMgr::report( N_ERH_ErrorMgr::USR_FATAL_0, tmp);
649 if (DEBUG_NONLINEAR && isActive(Diag::SENS_SOLVER))
651 int maxParamStringSize_ = 16;
652 Xyce::dout() << std::setw(maxParamStringSize_)<< *iterParam
653 <<
" dp = " << std::setw(11)<< std::scientific<< std::setprecision(4) << dp
654 <<
" original value = " << std::setw(16)<< std::scientific<< std::setprecision(9) << paramOrig
655 <<
" modified value = " << std::setw(16)<< std::scientific<< std::setprecision(9) << paramPerturbed
659 nEqLoader_.setParam (paramName, paramPerturbed);
664 rolSweep_.setSweepValue(k);
666 nEqLoader_.loadRHS();
669 pertFVectorPtr_->update(1.0, *(analysisManager_.getDataStore()->daeFVectorPtr), 0.0);
670 pertQVectorPtr_->update(1.0, *(analysisManager_.getDataStore()->daeQVectorPtr), 0.0);
671 pertBVectorPtr_->update(1.0, *(analysisManager_.getDataStore()->daeBVectorPtr), 0.0);
675 rolSweep_.mydfdpPtrVector_[iparam]->putScalar(0.0);
676 rolSweep_.mydfdpPtrVector_[iparam]->addVec (+1.0, *(pertFVectorPtr_));
677 rolSweep_.mydfdpPtrVector_[iparam]->addVec (-1.0, *(origFVectorPtr_));
678 rolSweep_.mydfdpPtrVector_[iparam]->scale(rdp);
681 rolSweep_.mydqdpPtrVector_[iparam]->putScalar(0.0);
682 rolSweep_.mydqdpPtrVector_[iparam]->addVec (+1.0, *(pertQVectorPtr_));
683 rolSweep_.mydqdpPtrVector_[iparam]->addVec (-1.0, *(origQVectorPtr_));
684 rolSweep_.mydqdpPtrVector_[iparam]->scale(rdp);
687 rolSweep_.mydbdpPtrVector_[iparam]->putScalar(0.0);
688 rolSweep_.mydbdpPtrVector_[iparam]->addVec (+1.0, *(pertBVectorPtr_));
689 rolSweep_.mydbdpPtrVector_[iparam]->addVec (-1.0, *(origBVectorPtr_));
690 rolSweep_.mydbdpPtrVector_[iparam]->scale(rdp);
692 if (DEBUG_NONLINEAR && isActive(Diag::SENS_SOLVER))
694 Xyce::dout() << *iterParam <<
": ";
695 Xyce::dout().width(15); Xyce::dout().precision(7); Xyce::dout().setf(std::ios::scientific);
696 Xyce::dout() <<
"deviceSens_dp = " << dp << std::endl;
698 for (
int k1 = 0; k1 < analysisManager_.getDataStore()->solutionSize; ++k1)
702 <<
"fpert["<<std::setw(3)<<k1<<
"]= "<<std::setw(15)<<std::scientific<<std::setprecision(8)<<(*(pertFVectorPtr_))[k1]
703 <<
" forig["<<std::setw(3)<<k1<<
"]= "<<std::setw(15)<<std::scientific<<std::setprecision(8)<<(*(origFVectorPtr_))[k1]
704 <<
" dfdp ["<<std::setw(3)<<k1<<
"]= "<<std::setw(15)<<std::scientific<<std::setprecision(8)<<(*(rolSweep_.mydfdpPtrVector_[iparam]))[k1]
708 Xyce::dout() << std::endl;
709 for (
int k1 = 0; k1 < analysisManager_.getDataStore()->solutionSize; ++k1)
712 <<
"qpert["<<std::setw(3)<<k1<<
"]= "<<std::setw(15)<<std::scientific<<std::setprecision(8)<<(*(pertQVectorPtr_))[k1]
713 <<
" qorig["<<std::setw(3)<<k1<<
"]= "<<std::setw(15)<<std::scientific<<std::setprecision(8)<<(*(origQVectorPtr_))[k1]
714 <<
" dqdp ["<<std::setw(3)<<k1<<
"]= "<<std::setw(15)<<std::scientific<<std::setprecision(8)<<(*(rolSweep_.mydqdpPtrVector_[iparam]))[k1]
718 Xyce::dout() << std::endl ;
719 for (
int k1 = 0; k1 < analysisManager_.getDataStore()->solutionSize; ++k1)
722 <<
"bpert["<<std::setw(3)<<k1<<
"]= "<<std::setw(15)<<std::scientific<<std::setprecision(8)<<(*(pertBVectorPtr_))[k1]
723 <<
" borig["<<std::setw(3)<<k1<<
"]= "<<std::setw(15)<<std::scientific<<std::setprecision(8)<<(*(origBVectorPtr_))[k1]
724 <<
" dbdp ["<<std::setw(3)<<k1<<
"]= "<<std::setw(15)<<std::scientific<<std::setprecision(8)<<(*(rolSweep_.mydbdpPtrVector_[iparam]))[k1]
767 nEqLoader_.setParam (paramName, paramOrig);
769 analysisManager_.getDataStore()->daeFVectorPtr->update(1.0, *(origFVectorPtr_), 0.0);
770 analysisManager_.getDataStore()->daeQVectorPtr->update(1.0, *(origQVectorPtr_), 0.0);
771 analysisManager_.getDataStore()->daeBVectorPtr->update(1.0, *(origBVectorPtr_), 0.0);
776 Linear::Vector & sensRHS = *(rolSweep_.mysensRHSPtrVector_[iparam]);
777 Linear::Vector & dfdpRef = *(rolSweep_.mydfdpPtrVector_[iparam]);
778 Linear::Vector & dbdpRef = *(rolSweep_.mydbdpPtrVector_[iparam]);
780 sensRHS.linearCombo(0.0,sensRHS,+1.0,dfdpRef);
781 sensRHS.linearCombo(1.0,sensRHS,-1.0,dbdpRef);
786 sensRHS.scale( (*vp)[iparam] );
789 (*jvp)[k]->addVec(1.0,sensRHS);
794 nEqLoader_.setVoltageLimiterStatus(vstatus);
798 void applyAdjointJacobian_2(::ROL::Vector<Real> &jv, const ::ROL::Vector<Real> &v, const ::ROL::Vector<Real> &u, const ::ROL::Vector<Real> &z, Real &tol)
800 Teuchos::RCP<const std::vector<Teuchos::RCP<Linear::Vector> > > vp = (Teuchos::dyn_cast< Linear::ROL_XyceVector<Real> >(
const_cast< ::ROL::Vector<Real> &
>(v))).getVector();
801 Teuchos::RCP<const std::vector<Teuchos::RCP<Linear::Vector> > > up = (Teuchos::dyn_cast< Linear::ROL_XyceVector<Real> >(
const_cast< ::ROL::Vector<Real> &
>(u))).getVector();
802 Teuchos::RCP<std::vector<Real> > jvp = Teuchos::rcp_const_cast<std::vector<Real> >((Teuchos::dyn_cast<::ROL::StdVector<Real> >(jv)).getVector());
807 setUncertainParams();
809 bool vstatus = nEqLoader_.getVoltageLimiterStatus();
810 nEqLoader_.setVoltageLimiterStatus(
false);
813 for (
int k=0;k<nc_;k++)
815 rolSweep_.setSweepValue(k);
822 success = analysisManager_.getDataStore()->setNextSolVectorPtr(*(*up)[k]);
825 success = nEqLoader_.loadRHS();
830 analysisManager_.getDataStore()->paramOrigVals_.clear();
835 nEqLoader_.loadJacobian ();
840 std::vector<std::string>::iterator firstParam = rolSweep_.paramNameVec_.begin ();
841 std::vector<std::string>::iterator lastParam = rolSweep_.paramNameVec_.end ();
842 std::vector<std::string>::iterator iterParam;
843 for ( iterParam=firstParam, iparam=0;
844 iterParam!=lastParam; ++iterParam, ++iparam )
846 std::string paramName(*iterParam);
850 double paramOrig = 0.0;
851 bool found = nEqLoader_.getParamAndReduce(analysisManager_.getPDSManager()->getPDSComm()->comm(),paramName, paramOrig);
855 std::string msg(
"EqualityConstraint::applyJacobian_2: cannot find parameter ");
857 N_ERH_ErrorMgr::report (N_ERH_ErrorMgr::DEV_FATAL, msg);
860 analysisManager_.getDataStore()->paramOrigVals_.push_back(paramOrig);
863 bool analyticAvailable =
false;
866 analyticAvailable = nEqLoader_.analyticSensitivitiesAvailable (paramName);
868 if (analyticAvailable)
870 std::vector<double> dfdpVec;
871 std::vector<double> dqdpVec;
872 std::vector<double> dbdpVec;
874 std::vector<int> FindicesVec;
875 std::vector<int> QindicesVec;
876 std::vector<int> BindicesVec;
878 nEqLoader_.getAnalyticSensitivities(paramName,dfdpVec,dqdpVec,dbdpVec,FindicesVec, QindicesVec, BindicesVec);
880 rolSweep_.mydfdpPtrVector_[iparam]->putScalar(0.0);
881 rolSweep_.mydqdpPtrVector_[iparam]->putScalar(0.0);
882 rolSweep_.mydbdpPtrVector_[iparam]->putScalar(0.0);
884 int Fsize=FindicesVec.size();
885 for (
int i=0;i<Fsize;++i)
887 Linear::Vector & dfdpRef = *(rolSweep_.mydfdpPtrVector_[iparam]);
888 dfdpRef[FindicesVec[i]] += dfdpVec[i];
891 int Qsize=QindicesVec.size();
892 for (
int i=0;i<Qsize;++i)
894 Linear::Vector & dqdpRef = *(rolSweep_.mydqdpPtrVector_[iparam]);
895 dqdpRef[QindicesVec[i]] += dqdpVec[i];
898 int Bsize=BindicesVec.size();
899 for (
int i=0;i<Bsize;++i)
901 Linear::Vector & dbdpRef = *(rolSweep_.mydbdpPtrVector_[iparam]);
902 dbdpRef[BindicesVec[i]] += dbdpVec[i];
905 rolSweep_.mydfdpPtrVector_[iparam]->fillComplete();
906 rolSweep_.mydqdpPtrVector_[iparam]->fillComplete();
907 rolSweep_.mydbdpPtrVector_[iparam]->fillComplete();
909 if (DEBUG_NONLINEAR && isActive(Diag::SENS_SOLVER))
911 Xyce::dout() << *iterParam <<
": ";
912 Xyce::dout().setf(std::ios::scientific);
913 Xyce::dout() << std::endl;
915 for (
int k1 = 0; k1 < analysisManager_.getDataStore()->solutionSize; ++k1)
918 <<
"dfdp["<<std::setw(3)<<k1<<
"]= "<<std::setw(15)<<std::scientific<<std::setprecision(8)<<(*(rolSweep_.mydfdpPtrVector_[iparam]))[k1]
919 <<
" dqdp["<<std::setw(3)<<k1<<
"]= "<<std::setw(15)<<std::scientific<<std::setprecision(8)<<(*(rolSweep_.mydqdpPtrVector_[iparam]))[k1]
920 <<
" dbdp["<<std::setw(3)<<k1<<
"]= "<<std::setw(15)<<std::scientific<<std::setprecision(8)<<(*(rolSweep_.mydbdpPtrVector_[iparam]))[k1]
929 if (DEBUG_NONLINEAR && isActive(Diag::SENS_SOLVER))
931 Xyce::dout() << std::endl <<
" Calculating numerical df/dp, dq/dp and db/dp for: ";
932 Xyce::dout() << *iterParam << std::endl;
936 origFVectorPtr_->update(1.0, *(analysisManager_.getDataStore()->daeFVectorPtr), 0.0);
937 origQVectorPtr_->update(1.0, *(analysisManager_.getDataStore()->daeQVectorPtr), 0.0);
938 origBVectorPtr_->update(1.0, *(analysisManager_.getDataStore()->daeBVectorPtr), 0.0);
941 double sqrtEta_ = pow(10,
int(log10(fabs(paramOrig))));;
942 double dp = sqrtEta_ * 1.e-4 * (1.0 + fabs(paramOrig));
943 double paramPerturbed = paramOrig;
948 paramPerturbed += dp;
952 paramPerturbed -= dp;
956 static std::string tmp =
"difference=central not supported.\n";
957 N_ERH_ErrorMgr::report( N_ERH_ErrorMgr::USR_FATAL_0, tmp);
961 static std::string tmp =
"difference not recognized!\n";
962 N_ERH_ErrorMgr::report( N_ERH_ErrorMgr::USR_FATAL_0, tmp);
965 if (DEBUG_NONLINEAR && isActive(Diag::SENS_SOLVER))
967 int maxParamStringSize_ = 16;
968 Xyce::dout() << std::setw(maxParamStringSize_)<< *iterParam
969 <<
" dp = " << std::setw(11)<< std::scientific<< std::setprecision(4) << dp
970 <<
" original value = " << std::setw(16)<< std::scientific<< std::setprecision(9) << paramOrig
971 <<
" modified value = " << std::setw(16)<< std::scientific<< std::setprecision(9) << paramPerturbed
975 nEqLoader_.setParam (paramName, paramPerturbed);
980 rolSweep_.setSweepValue(k);
982 nEqLoader_.loadRHS();
985 pertFVectorPtr_->update(1.0, *(analysisManager_.getDataStore()->daeFVectorPtr), 0.0);
986 pertQVectorPtr_->update(1.0, *(analysisManager_.getDataStore()->daeQVectorPtr), 0.0);
987 pertBVectorPtr_->update(1.0, *(analysisManager_.getDataStore()->daeBVectorPtr), 0.0);
991 rolSweep_.mydfdpPtrVector_[iparam]->putScalar(0.0);
992 rolSweep_.mydfdpPtrVector_[iparam]->addVec (+1.0, *(pertFVectorPtr_));
993 rolSweep_.mydfdpPtrVector_[iparam]->addVec (-1.0, *(origFVectorPtr_));
994 rolSweep_.mydfdpPtrVector_[iparam]->scale(rdp);
997 rolSweep_.mydqdpPtrVector_[iparam]->putScalar(0.0);
998 rolSweep_.mydqdpPtrVector_[iparam]->addVec (+1.0, *(pertQVectorPtr_));
999 rolSweep_.mydqdpPtrVector_[iparam]->addVec (-1.0, *(origQVectorPtr_));
1000 rolSweep_.mydqdpPtrVector_[iparam]->scale(rdp);
1003 rolSweep_.mydbdpPtrVector_[iparam]->putScalar(0.0);
1004 rolSweep_.mydbdpPtrVector_[iparam]->addVec (+1.0, *(pertBVectorPtr_));
1005 rolSweep_.mydbdpPtrVector_[iparam]->addVec (-1.0, *(origBVectorPtr_));
1006 rolSweep_.mydbdpPtrVector_[iparam]->scale(rdp);
1008 if (DEBUG_NONLINEAR && isActive(Diag::SENS_SOLVER))
1010 Xyce::dout() << *iterParam <<
": ";
1011 Xyce::dout().width(15); Xyce::dout().precision(7); Xyce::dout().setf(std::ios::scientific);
1012 Xyce::dout() <<
"deviceSens_dp = " << dp << std::endl;
1014 for (
int k1 = 0; k1 < analysisManager_.getDataStore()->solutionSize; ++k1)
1018 <<
"fpert["<<std::setw(3)<<k1<<
"]= "<<std::setw(15)<<std::scientific<<std::setprecision(8)<<(*(pertFVectorPtr_))[k1]
1019 <<
" forig["<<std::setw(3)<<k1<<
"]= "<<std::setw(15)<<std::scientific<<std::setprecision(8)<<(*(origFVectorPtr_))[k1]
1020 <<
" dfdp ["<<std::setw(3)<<k1<<
"]= "<<std::setw(15)<<std::scientific<<std::setprecision(8)<<(*(rolSweep_.mydfdpPtrVector_[iparam]))[k1]
1024 Xyce::dout() << std::endl;
1025 for (
int k1 = 0; k1 < analysisManager_.getDataStore()->solutionSize; ++k1)
1028 <<
"qpert["<<std::setw(3)<<k1<<
"]= "<<std::setw(15)<<std::scientific<<std::setprecision(8)<<(*(pertQVectorPtr_))[k1]
1029 <<
" qorig["<<std::setw(3)<<k1<<
"]= "<<std::setw(15)<<std::scientific<<std::setprecision(8)<<(*(origQVectorPtr_))[k1]
1030 <<
" dqdp ["<<std::setw(3)<<k1<<
"]= "<<std::setw(15)<<std::scientific<<std::setprecision(8)<<(*(rolSweep_.mydqdpPtrVector_[iparam]))[k1]
1034 Xyce::dout() << std::endl ;
1035 for (
int k1 = 0; k1 < analysisManager_.getDataStore()->solutionSize; ++k1)
1038 <<
"bpert["<<std::setw(3)<<k1<<
"]= "<<std::setw(15)<<std::scientific<<std::setprecision(8)<<(*(pertBVectorPtr_))[k1]
1039 <<
" borig["<<std::setw(3)<<k1<<
"]= "<<std::setw(15)<<std::scientific<<std::setprecision(8)<<(*(origBVectorPtr_))[k1]
1040 <<
" dbdp ["<<std::setw(3)<<k1<<
"]= "<<std::setw(15)<<std::scientific<<std::setprecision(8)<<(*(rolSweep_.mydbdpPtrVector_[iparam]))[k1]
1083 nEqLoader_.setParam (paramName, paramOrig);
1085 analysisManager_.getDataStore()->daeFVectorPtr->update(1.0, *(origFVectorPtr_), 0.0);
1086 analysisManager_.getDataStore()->daeQVectorPtr->update(1.0, *(origQVectorPtr_), 0.0);
1087 analysisManager_.getDataStore()->daeBVectorPtr->update(1.0, *(origBVectorPtr_), 0.0);
1092 Linear::Vector & sensRHS = *(rolSweep_.mysensRHSPtrVector_[iparam]);
1093 Linear::Vector & dfdpRef = *(rolSweep_.mydfdpPtrVector_[iparam]);
1094 Linear::Vector & dbdpRef = *(rolSweep_.mydbdpPtrVector_[iparam]);
1096 sensRHS.linearCombo(0.0,sensRHS,+1.0,dfdpRef);
1097 sensRHS.linearCombo(1.0,sensRHS,-1.0,dbdpRef);
1099 sensRHS.scale(-1.0);
1102 Real dot = sensRHS.dotProduct(*((*vp)[k]));
1105 (*jvp)[iparam] += dot;
1110 nEqLoader_.setVoltageLimiterStatus(vstatus);
1115 void applyAdjointHessian_11(::ROL::Vector<Real> &ahwv, const ::ROL::Vector<Real> &w, const ::ROL::Vector<Real> &v, const ::ROL::Vector<Real> &u, const ::ROL::Vector<Real> &z, Real &tol)
1119 void applyAdjointHessian_12(::ROL::Vector<Real> &ahwv, const ::ROL::Vector<Real> &w, const ::ROL::Vector<Real> &v, const ::ROL::Vector<Real> &u, const ::ROL::Vector<Real> &z, Real &tol)
1123 void applyAdjointHessian_21(::ROL::Vector<Real> &ahwv, const ::ROL::Vector<Real> &w, const ::ROL::Vector<Real> &v, const ::ROL::Vector<Real> &u, const ::ROL::Vector<Real> &z, Real &tol)
1127 void applyAdjointHessian_22(::ROL::Vector<Real> &ahwv, const ::ROL::Vector<Real> &w, const ::ROL::Vector<Real> &v, const ::ROL::Vector<Real> &u, const ::ROL::Vector<Real> &z, Real &tol)
1133 #if IDENTITY_PRECONDITIONER==1
1134 void applyPreconditioner(::ROL::Vector<Real> &pv, const ::ROL::Vector<Real> &v, const ::ROL::Vector<Real> &x, const::ROL::Vector<Real> &g, Real &tol)
1142 template<
class Real>
1143 class Objective_DC_L2Norm :
public ::ROL::ParametrizedObjective_SimOpt<Real>
1148 Teuchos::RCP<std::vector<Real> > Imeas_;
1153 Objective_DC_L2Norm(Real alpha,
int ns,
int nz)
1159 Imeas_ = Teuchos::rcp(
new std::vector<Real>(ns_,0.0));
1161 std::ifstream measurements(
"measurementsDC.txt");
1162 std::cout <<
"Reading measurements" << std::endl;
1163 if (measurements.is_open())
1165 for (
int i=0;i<ns_;i++)
1167 measurements >> temp;
1168 measurements >> temp;
1169 (*Imeas_)[i] = temp;
1174 std::cout <<
"Could not open measurements file" << std::endl;
1177 measurements.close();
1181 Real
value( const ::ROL::Vector<Real> &u, const ::ROL::Vector<Real> &z, Real &tol )
1183 Teuchos::RCP<const std::vector<Teuchos::RCP<Linear::Vector> > > up = (Teuchos::dyn_cast< Linear::ROL_XyceVector<Real> >(
const_cast< ::ROL::Vector<Real> &
>(u))).getVector();
1184 Teuchos::RCP<const std::vector<Real> > zp = (Teuchos::dyn_cast< ::ROL::StdVector<Real> >(
const_cast< ::ROL::Vector<Real> &
>(z))).getVector();
1189 for (
int i=0;i<ns_;i++)
1191 current = (*(*up)[i])[solindex_];
1192 val += (current - (*Imeas_)[i])*(current - (*Imeas_)[i]);
1197 void gradient_1( ::ROL::Vector<Real> &g, const ::ROL::Vector<Real> &u, const ::ROL::Vector<Real> &z, Real &tol )
1200 Teuchos::RCP< std::vector<Teuchos::RCP<Linear::Vector> > > gup = Teuchos::rcp_const_cast< std::vector<Teuchos::RCP<Linear::Vector> > >((Teuchos::dyn_cast<Linear::ROL_XyceVector<Real> >(g)).getVector());
1202 Teuchos::RCP<const std::vector<Teuchos::RCP<Linear::Vector> > > up = (Teuchos::dyn_cast< Linear::ROL_XyceVector<Real> >(
const_cast< ::ROL::Vector<Real> &
>(u))).getVector();
1203 Teuchos::RCP<const std::vector<Real> > zp =
1204 (Teuchos::dyn_cast< ::ROL::StdVector<Real> >(
const_cast< ::ROL::Vector<Real> &
>(z))).getVector();
1208 for (
int i=0; i<ns_; i++)
1210 current = (*(*up)[i])[solindex_];
1211 (*(*gup)[i])[solindex_] = (current-(*Imeas_)[i]);
1216 void gradient_2( ::ROL::Vector<Real> &g, const ::ROL::Vector<Real> &u, const ::ROL::Vector<Real> &z, Real &tol )
1219 Teuchos::RCP<std::vector<Real> > gzp = Teuchos::rcp_const_cast<std::vector<Real> >((Teuchos::dyn_cast<const ::ROL::StdVector<Real> >(g)).getVector());
1221 Teuchos::RCP<const std::vector<Teuchos::RCP<Linear::Vector> > > up = (Teuchos::dyn_cast< Linear::ROL_XyceVector<Real> >(
const_cast< ::ROL::Vector<Real> &
>(u))).getVector();
1222 Teuchos::RCP<const std::vector<Real> > zp = (Teuchos::dyn_cast< ::ROL::StdVector<Real> >(
const_cast< ::ROL::Vector<Real> &
>(z))).getVector();
1224 for (
int i=0; i<nz_; i++)
1230 void hessVec_11( ::ROL::Vector<Real> &hv, const ::ROL::Vector<Real> &v,
1231 const ::ROL::Vector<Real> &u, const ::ROL::Vector<Real> &z, Real &tol )
1233 Teuchos::RCP< std::vector<Teuchos::RCP<Linear::Vector> > > hvup = Teuchos::rcp_const_cast< std::vector<Teuchos::RCP<Linear::Vector> > >((Teuchos::dyn_cast<Linear::ROL_XyceVector<Real> >(hv)).getVector());
1235 Teuchos::RCP<const std::vector<Teuchos::RCP<Linear::Vector> > > vup = (Teuchos::dyn_cast< Linear::ROL_XyceVector<Real> >(
const_cast< ::ROL::Vector<Real> &
>(v))).getVector();
1239 for(
int i=0; i<ns_; i++)
1241 (*(*hvup)[i])[solindex_] = (*(*vup)[i])[solindex_];
1245 void hessVec_12( ::ROL::Vector<Real> &hv, const ::ROL::Vector<Real> &v,
1246 const ::ROL::Vector<Real> &u, const ::ROL::Vector<Real> &z, Real &tol )
1251 void hessVec_21( ::ROL::Vector<Real> &hv, const ::ROL::Vector<Real> &v,
1252 const ::ROL::Vector<Real> &u, const ::ROL::Vector<Real> &z, Real &tol )
1257 void hessVec_22( ::ROL::Vector<Real> &hv, const ::ROL::Vector<Real> &v,
1258 const ::ROL::Vector<Real> &u, const ::ROL::Vector<Real> &z, Real &tol )
1265 template<
class Real>
1266 class Objective_DC_AMP :
virtual public ::ROL::ParametrizedObjective_SimOpt<Real>
1269 int ns_,nz_,first_,middle_,last_;
1274 Objective_DC_AMP(
int ns,
int nz)
1284 Real
value( const ::ROL::Vector<Real> &u, const ::ROL::Vector<Real> &z, Real &tol )
1286 Teuchos::RCP<const std::vector<Teuchos::RCP<Linear::Vector> > > up = (Teuchos::dyn_cast< Linear::ROL_XyceVector<Real> >(
const_cast< ::ROL::Vector<Real> &
>(u))).getVector();
1287 Teuchos::RCP<const std::vector<Real> > zp = (Teuchos::dyn_cast< ::ROL::StdVector<Real> >(
const_cast< ::ROL::Vector<Real> &
>(z))).getVector();
1290 Real V3_0 = (*(*up)[first_ ])[solindex_];
1291 Real V3_10 = (*(*up)[middle_])[solindex_];
1292 Real V3_20 = (*(*up)[last_ ])[solindex_];
1293 Real A = (V3_0 - V3_20)/2.0;
1294 std::cout <<
"A = " << A <<
' ' << (*zp)[0] <<
' ' << (*zp)[1] << std::endl;
1295 for (
int i=0;i<ns_;i++)
1297 V3 = (*(*up)[i])[solindex_];
1298 t = -1.0 + i*(2.0/(ns_-1));
1299 val += ( A*t - V3 + V3_10 )*( A*t - V3 + V3_10 );
1304 void gradient_1( ::ROL::Vector<Real> &g, const ::ROL::Vector<Real> &u, const ::ROL::Vector<Real> &z, Real &tol )
1307 Teuchos::RCP< std::vector<Teuchos::RCP<Linear::Vector> > > gup = Teuchos::rcp_const_cast< std::vector<Teuchos::RCP<Linear::Vector> > >((Teuchos::dyn_cast<Linear::ROL_XyceVector<Real> >(g)).getVector());
1309 Teuchos::RCP<const std::vector<Teuchos::RCP<Linear::Vector> > > up = (Teuchos::dyn_cast< Linear::ROL_XyceVector<Real> >(
const_cast< ::ROL::Vector<Real> &
>(u))).getVector();
1310 Teuchos::RCP<const std::vector<Real> > zp =
1311 (Teuchos::dyn_cast< ::ROL::StdVector<Real> >(
const_cast< ::ROL::Vector<Real> &
>(z))).getVector();
1315 Real V3_0 = (*(*up)[first_ ])[solindex_];
1316 Real V3_10 = (*(*up)[middle_])[solindex_];
1317 Real V3_20 = (*(*up)[last_ ])[solindex_];
1318 Real A = (V3_0 - V3_20)/2.0;
1319 (*(*gup)[first_ ])[solindex_] += ( -A - V3_0 + V3_10 )*(-3.0/2.0) + ( A - V3_20 + V3_10 )/2.0;
1320 (*(*gup)[last_ ])[solindex_] += ( A - V3_20 + V3_10 )*(-3.0/2.0) + ( -A - V3_0 + V3_10 )/2.0;
1321 (*(*gup)[middle_])[solindex_] += 2*V3_10 - V3_0 - V3_20;
1322 for (
int i=1; i<ns_-1; i++)
1324 V3 = (*(*up)[i])[solindex_];
1325 t = -1.0 + i*(2.0/(ns_-1));
1326 (*(*gup)[i ])[solindex_] -= ( A*t - V3 + V3_10 );
1327 (*(*gup)[first_ ])[solindex_] += ( A*t - V3 + V3_10 )*(t/2.0);
1328 (*(*gup)[last_ ])[solindex_] -= ( A*t - V3 + V3_10 )*(t/2.0);
1329 (*(*gup)[middle_])[solindex_] += ( A*t - V3 + V3_10 );
1333 void gradient_2(::ROL::Vector<Real> &g, const ::ROL::Vector<Real> &u, const ::ROL::Vector<Real> &z, Real &tol )
1336 Teuchos::RCP<std::vector<Real> > gzp = Teuchos::rcp_const_cast<std::vector<Real> >((Teuchos::dyn_cast<const ::ROL::StdVector<Real> >(g)).getVector());
1338 Teuchos::RCP<const std::vector<Teuchos::RCP<Linear::Vector> > > up = (Teuchos::dyn_cast< Linear::ROL_XyceVector<Real> >(
const_cast< ::ROL::Vector<Real> &
>(u))).getVector();
1339 Teuchos::RCP<const std::vector<Real> > zp = (Teuchos::dyn_cast< ::ROL::StdVector<Real> >(
const_cast< ::ROL::Vector<Real> &
>(z))).getVector();
1341 for (
int i=0; i<nz_; i++)
1347 void hessVec_11( ::ROL::Vector<Real> &hv, const ::ROL::Vector<Real> &v,
1348 const ::ROL::Vector<Real> &u, const ::ROL::Vector<Real> &z, Real &tol )
1350 Teuchos::RCP< std::vector<Teuchos::RCP<Linear::Vector> > > hvup = Teuchos::rcp_const_cast< std::vector<Teuchos::RCP<Linear::Vector> > >((Teuchos::dyn_cast<Linear::ROL_XyceVector<Real> >(hv)).getVector());
1352 Teuchos::RCP<const std::vector<Teuchos::RCP<Linear::Vector> > > up = (Teuchos::dyn_cast< Linear::ROL_XyceVector<Real> >(
const_cast< ::ROL::Vector<Real> &
>(u))).getVector();
1353 Teuchos::RCP<const std::vector<Teuchos::RCP<Linear::Vector> > > vup = (Teuchos::dyn_cast< Linear::ROL_XyceVector<Real> >(
const_cast< ::ROL::Vector<Real> &
>(v))).getVector();
1356 Real V3_0 = (*(*vup)[first_ ])[solindex_];
1357 Real V3_10 = (*(*vup)[middle_])[solindex_];
1358 Real V3_20 = (*(*vup)[last_ ])[solindex_];
1360 for(
int i=1;i<ns_-1;i++)
1362 sum1 += std::pow((-1.0 + i*2.0/(ns_-1)),2);
1365 (*(*hvup)[first_ ])[solindex_] += (sum1+5.0/2.0)*V3_0 - V3_10 - (sum1+3.0/2.0)*V3_20;
1366 (*(*hvup)[middle_])[solindex_] = - V3_0 - V3_20;
1367 (*(*hvup)[last_ ])[solindex_] += -(sum1+3.0/2.0)*V3_0 - V3_10 + (sum1+5.0/2.0)*V3_20;
1368 for(
int i=1; i<ns_-1; i++)
1370 t = -1.0 + i*2.0/(ns_-1);
1371 V3 = (*(*vup)[i])[solindex_];
1372 (*(*hvup)[i])[solindex_] += V3
1376 (*(*hvup)[first_ ])[solindex_] -= (t/2.0)*V3;
1377 (*(*hvup)[middle_])[solindex_] += (i==middle_) ? (ns_-1)*V3 : -V3;
1378 (*(*hvup)[last_ ])[solindex_] += (t/2.0)*V3;
1382 void hessVec_12(::ROL::Vector<Real> &hv, const ::ROL::Vector<Real> &v,
1383 const ::ROL::Vector<Real> &u, const ::ROL::Vector<Real> &z, Real &tol )
1388 void hessVec_21( ::ROL::Vector<Real> &hv, const ::ROL::Vector<Real> &v,
1389 const ::ROL::Vector<Real> &u, const ::ROL::Vector<Real> &z, Real &tol )
1394 void hessVec_22( ::ROL::Vector<Real> &hv, const ::ROL::Vector<Real> &v,
1395 const ::ROL::Vector<Real> &u, const ::ROL::Vector<Real> &z, Real &tol )
1402 template<
class Real>
1403 class Penalty_DC_AMP :
virtual public ::ROL::ParametrizedObjective_SimOpt<Real>
1410 int ns_,nz_,first_,middle_,last_;
1415 Penalty_DC_AMP(
int ptype, Real alpha, Real ampl,
int ns,
int nz)
1429 void switch_ptype(
int ptype)
1434 Real
value( const ::ROL::Vector<Real> &u, const ::ROL::Vector<Real> &z, Real &tol )
1436 Teuchos::RCP<const std::vector<Teuchos::RCP<Linear::Vector> > > up = (Teuchos::dyn_cast< Linear::ROL_XyceVector<Real> >(
const_cast< ::ROL::Vector<Real> &
>(u))).getVector();
1437 Teuchos::RCP<const std::vector<Real> > zp = (Teuchos::dyn_cast< ::ROL::StdVector<Real> >(
const_cast< ::ROL::Vector<Real> &
>(z))).getVector();
1439 Real A = ((*(*up)[first_])[solindex_] - (*(*up)[last_])[solindex_] )/2.0;
1441 val += alpha_*(0.5+0.5*tanh(epsilon_*(ampl_ - A*A)));
1443 val += alpha_*( (ampl_ - A*A > 0) ? ampl_-A*A : 0 );
1445 val += alpha_*std::pow(( (ampl_ - A*A > 0) ? ampl_-A*A : 0 ),2);
1447 val += alpha_*std::pow(( (A - ampl_ > 0) ? A - ampl_ : 0 ),2);
1453 void gradient_1( ::ROL::Vector<Real> &g, const ::ROL::Vector<Real> &u, const ::ROL::Vector<Real> &z, Real &tol )
1456 Teuchos::RCP< std::vector<Teuchos::RCP<Linear::Vector> > > gup = Teuchos::rcp_const_cast< std::vector<Teuchos::RCP<Linear::Vector> > >((Teuchos::dyn_cast<Linear::ROL_XyceVector<Real> >(g)).getVector());
1458 Teuchos::RCP<const std::vector<Teuchos::RCP<Linear::Vector> > > up = (Teuchos::dyn_cast< Linear::ROL_XyceVector<Real> >(
const_cast< ::ROL::Vector<Real> &
>(u))).getVector();
1459 Teuchos::RCP<const std::vector<Real> > zp =
1460 (Teuchos::dyn_cast< ::ROL::StdVector<Real> >(
const_cast< ::ROL::Vector<Real> &
>(z))).getVector();
1463 Real A = ((*(*up)[first_])[solindex_] - (*(*up)[last_])[solindex_] )/2.0;
1466 (*(*gup)[first_])[solindex_] -= 0.5*(0.5*alpha_*epsilon_*A)/std::pow(cosh(epsilon_*(ampl_-A*A)),2);
1467 (*(*gup)[last_ ])[solindex_] += 0.5*(0.5*alpha_*epsilon_*A)/std::pow(cosh(epsilon_*(ampl_-A*A)),2);
1472 (*(*gup)[first_])[solindex_] -= 0.5*alpha_*( (ampl_ - A*A > 0) ? A : 0 );
1473 (*(*gup)[last_ ])[solindex_] += 0.5*alpha_*( (ampl_ - A*A > 0) ? A : 0 );
1478 (*(*gup)[first_])[solindex_] -= 0.5*alpha_*( (ampl_ - A*A > 0) ? 8.0*A - 2.0*A*A*A : 0 );
1479 (*(*gup)[last_ ])[solindex_] += 0.5*alpha_*( (ampl_ - A*A > 0) ? 8.0*A - 2.0*A*A*A : 0 );
1484 (*(*gup)[first_])[solindex_] += alpha_*( (A - ampl_ > 0) ? 0.5*(A - ampl_) : 0 );
1485 (*(*gup)[last_ ])[solindex_] -= alpha_*( (A - ampl_ > 0) ? 0.5*(A - ampl_) : 0 );
1490 (*(*gup)[first_])[solindex_] = 0.5*alpha_;
1491 (*(*gup)[last_ ])[solindex_] = -0.5*alpha_;
1495 void gradient_2( ::ROL::Vector<Real> &g, const ::ROL::Vector<Real> &u, const ::ROL::Vector<Real> &z, Real &tol )
1498 Teuchos::RCP<std::vector<Real> > gzp = Teuchos::rcp_const_cast<std::vector<Real> >((Teuchos::dyn_cast<const ::ROL::StdVector<Real> >(g)).getVector());
1500 Teuchos::RCP<const std::vector<Teuchos::RCP<Linear::Vector> > > up = (Teuchos::dyn_cast< Linear::ROL_XyceVector<Real> >(
const_cast< ::ROL::Vector<Real> &
>(u))).getVector();
1501 Teuchos::RCP<const std::vector<Real> > zp = (Teuchos::dyn_cast< ::ROL::StdVector<Real> >(
const_cast< ::ROL::Vector<Real> &
>(z))).getVector();
1503 for (
int i=0; i<nz_; i++)
1509 void hessVec_11( ::ROL::Vector<Real> &hv, const ::ROL::Vector<Real> &v,
1510 const ::ROL::Vector<Real> &u, const ::ROL::Vector<Real> &z, Real &tol )
1512 Teuchos::RCP< std::vector<Teuchos::RCP<Linear::Vector> > > hvup = Teuchos::rcp_const_cast< std::vector<Teuchos::RCP<Linear::Vector> > >((Teuchos::dyn_cast<Linear::ROL_XyceVector<Real> >(hv)).getVector());
1514 Teuchos::RCP<const std::vector<Teuchos::RCP<Linear::Vector> > > up = (Teuchos::dyn_cast< Linear::ROL_XyceVector<Real> >(
const_cast< ::ROL::Vector<Real> &
>(u))).getVector();
1515 Teuchos::RCP<const std::vector<Teuchos::RCP<Linear::Vector> > > vup = (Teuchos::dyn_cast< Linear::ROL_XyceVector<Real> >(
const_cast< ::ROL::Vector<Real> &
>(v))).getVector();
1518 Real V3_0 = (*(*vup)[first_ ])[solindex_];
1519 Real V3_20 = (*(*vup)[last_ ])[solindex_];
1520 Real A = ((*(*up)[first_])[solindex_] - (*(*up)[last_])[solindex_] )/2.0;
1524 (*(*hvup)[first_])[solindex_] -= 0.5*alpha_*( (sinh(epsilon_*(ampl_-A*A))/std::pow(cosh(epsilon_*(ampl_-A*A)),3.0))*std::pow(epsilon_*A,2.0) + epsilon_/(std::pow(cosh(epsilon_*(ampl_-A*A)),2.0)*4.0) )*(V3_0 - V3_20);
1525 (*(*hvup)[last_ ])[solindex_] += 0.5*alpha_*( (sinh(epsilon_*(ampl_-A*A))/std::pow(cosh(epsilon_*(ampl_-A*A)),3.0))*std::pow(epsilon_*A,2.0) + epsilon_/(std::pow(cosh(epsilon_*(ampl_-A*A)),2.0)*4.0) )*(V3_0 - V3_20);
1530 (*(*hvup)[first_])[solindex_] += 0.5*alpha_*( (ampl_ - A*A > 0) ? -0.5 : 0 )*(V3_0 - V3_20);
1531 (*(*hvup)[last_ ])[solindex_] -= 0.5*alpha_*( (ampl_ - A*A > 0) ? -0.5 : 0 )*(V3_0 - V3_20);
1536 (*(*hvup)[first_])[solindex_] += 0.5*alpha_*( (ampl_ - A*A > 0) ? 3.0*A*A - 4.0 : 0 )*(V3_0 - V3_20);
1537 (*(*hvup)[last_])[solindex_] -= 0.5*alpha_*( (ampl_ - A*A > 0) ? 3.0*A*A - 4.0 : 0 )*(V3_0 - V3_20);
1542 (*(*hvup)[first_])[solindex_] += alpha_*( (A - ampl_ > 0) ? 0.25 : 0 )*(V3_0 - V3_20);
1543 (*(*hvup)[last_ ])[solindex_] += alpha_*( (A - ampl_ > 0) ? 0.25 : 0 )*(V3_20 - V3_0);
1548 void hessVec_12( ::ROL::Vector<Real> &hv, const ::ROL::Vector<Real> &v,
1549 const ::ROL::Vector<Real> &u, const ::ROL::Vector<Real> &z, Real &tol )
1554 void hessVec_21( ::ROL::Vector<Real> &hv, const ::ROL::Vector<Real> &v,
1555 const ::ROL::Vector<Real> &u, const ::ROL::Vector<Real> &z, Real &tol )
1560 void hessVec_22( ::ROL::Vector<Real> &hv, const ::ROL::Vector<Real> &v,
1561 const ::ROL::Vector<Real> &u, const ::ROL::Vector<Real> &z, Real &tol )
1568 template<
class Real>
1569 class SumObjective :
public virtual ::ROL::Objective<Real>
1573 std::vector< Teuchos::RCP< ::ROL::Objective<Real> > > objVec_;
1574 std::vector<bool> types_;
1578 SumObjective(std::vector< Teuchos::RCP< ::ROL::Objective<Real> > > & objVec,
1579 std::vector<bool> & types)
1580 : objVec_(objVec), types_(types)
1583 virtual void update( const ::ROL::Vector<Real> &x,
bool flag =
true,
int iter = -1 )
1585 for(
int i=0;i<objVec_.size();i++)
1589 ::ROL::CVaRVector<Real> &xc = Teuchos::dyn_cast< ::ROL::CVaRVector<Real> >(
const_cast< ::ROL::Vector<Real> &
>(x));
1590 objVec_[i]->update(*(xc.getVector()),flag,iter);
1592 else if(types_[i]==1)
1594 objVec_[i]->update(x,flag,iter);
1599 virtual Real
value( const ::ROL::Vector<Real> &x, Real &tol )
1602 for(
int i=0;i<objVec_.size();i++)
1606 ::ROL::CVaRVector<Real> &xc = Teuchos::dyn_cast< ::ROL::CVaRVector<Real> >(
const_cast< ::ROL::Vector<Real> &
>(x));
1607 val += objVec_[i]->value(*xc.getVector(),tol);
1609 else if(types_[i]==1)
1611 val += objVec_[i]->value(x,tol);
1617 virtual void gradient( ::ROL::Vector<Real> &g, const ::ROL::Vector<Real> &x, Real &tol )
1620 Teuchos::RCP< ::ROL::Vector<Real> > temp = g.clone();
1621 for(
int i=0;i<objVec_.size();i++)
1625 ::ROL::CVaRVector<Real> &tempc = Teuchos::dyn_cast< ::ROL::CVaRVector<Real> >(*temp);
1626 ::ROL::CVaRVector<Real> &xc = Teuchos::dyn_cast< ::ROL::CVaRVector<Real> >(
const_cast< ::ROL::Vector<Real> &
>(x));
1627 objVec_[i]->gradient(
const_cast< ::ROL::Vector<Real> &
>(*tempc.getVector()),*xc.getVector(),tol);
1630 else if(types_[i]==1)
1632 objVec_[i]->gradient(*temp,x,tol);
1638 virtual void hessVec( ::ROL::Vector<Real> &hv, const ::ROL::Vector<Real> &v, const ::ROL::Vector<Real> &x, Real &tol )
1641 Teuchos::RCP< ::ROL::Vector<Real> > temp = hv.clone();
1642 for(
int i=0;i<objVec_.size();i++)
1646 ::ROL::CVaRVector<Real> &tempc = Teuchos::dyn_cast< ::ROL::CVaRVector<Real> >(*temp);
1647 ::ROL::CVaRVector<Real> &xc = Teuchos::dyn_cast< ::ROL::CVaRVector<Real> >(
const_cast< ::ROL::Vector<Real> &
>(x));
1648 ::ROL::CVaRVector<Real> &vc = Teuchos::dyn_cast< ::ROL::CVaRVector<Real> >(
const_cast< ::ROL::Vector<Real> &
>(v));
1649 objVec_[i]->hessVec(
const_cast< ::ROL::Vector<Real> &
>(*tempc.getVector()),*vc.getVector(),*xc.getVector(),tol);
1652 else if(types_[i]==1)
1654 objVec_[i]->hessVec(*temp,v,x,tol);
1662 template<
class Real>
1663 class BoundConstraint_ROL_DC :
public ::ROL::BoundConstraint<Real>
1668 std::vector<Real> x_lo_;
1670 std::vector<Real> x_up_;
1678 BoundConstraint_ROL_DC( Real scale, std::vector<Real> lo_bounds, std::vector<Real> up_bounds )
1681 n = lo_bounds.size();
1682 for (
int i=0;i<n;i++)
1684 x_lo_.push_back(lo_bounds[i]);
1686 for (
int i=0;i<n;i++)
1688 x_up_.push_back(up_bounds[i]);
1692 min_diff_ = 1.e+200;
1693 for(
int i=0;i<n;i++)
1695 min_diff_ = std::min(x_up_[i]-x_lo_[i],min_diff_);
1700 void project( ::ROL::Vector<Real> &x )
1702 Teuchos::RCP<std::vector<Real> > ex =
1703 Teuchos::rcp_const_cast<std::vector<Real> >((Teuchos::dyn_cast<::ROL::StdVector<Real> >(x)).getVector());
1704 for (
int i=0;i<n;i++)
1706 (*ex)[i] = std::max(x_lo_[i],std::min(x_up_[i],(*ex)[i]));
1709 bool isFeasible( const ::ROL::Vector<Real> &x )
1711 Teuchos::RCP<const std::vector<Real> > ex =
1712 (Teuchos::dyn_cast<::ROL::StdVector<Real> >(
const_cast<::ROL::Vector<Real> &
>(x))).getVector();
1714 for (
int i=0;i<n;i++)
1716 bs = bs && (((*ex)[i] >= this->x_lo_[i])?
true:
false);
1717 bs = bs && (((*ex)[i] <= this->x_up_[i])?
true:
false);
1722 void pruneLowerActive(::ROL::Vector<Real> &v, const ::ROL::Vector<Real> &x, Real eps)
1724 Teuchos::RCP<const std::vector<Real> > ex =
1725 (Teuchos::dyn_cast<::ROL::StdVector<Real> >(
const_cast<::ROL::Vector<Real> &
>(x))).getVector();
1726 Teuchos::RCP<std::vector<Real> > ev =
1727 Teuchos::rcp_const_cast<std::vector<Real> >((Teuchos::dyn_cast<::ROL::StdVector<Real> >(v)).getVector());
1728 Real epsn = std::min(this->scale_*eps,this->min_diff_);
1729 for (
int i = 0; i < n; i++ )
1731 if ( ((*ex)[i] <= this->x_lo_[i]+epsn) )
1738 void pruneUpperActive(::ROL::Vector<Real> &v, const ::ROL::Vector<Real> &x, Real eps)
1740 Teuchos::RCP<const std::vector<Real> > ex =
1741 (Teuchos::dyn_cast<::ROL::StdVector<Real> >(
const_cast<::ROL::Vector<Real> &
>(x))).getVector();
1742 Teuchos::RCP<std::vector<Real> > ev =
1743 Teuchos::rcp_const_cast<std::vector<Real> >((Teuchos::dyn_cast<::ROL::StdVector<Real> >(v)).getVector());
1744 Real epsn = std::min(this->scale_*eps,this->min_diff_);
1745 for (
int i = 0; i < n; i++ )
1747 if ( ((*ex)[i] >= this->x_up_[i]-epsn) )
1754 void pruneActive(::ROL::Vector<Real> &v, const ::ROL::Vector<Real> &x, Real eps)
1756 Teuchos::RCP<const std::vector<Real> > ex =
1757 (Teuchos::dyn_cast<::ROL::StdVector<Real> >(
const_cast<::ROL::Vector<Real> &
>(x))).getVector();
1758 Teuchos::RCP<std::vector<Real> > ev =
1759 Teuchos::rcp_const_cast<std::vector<Real> >((Teuchos::dyn_cast<::ROL::StdVector<Real> >(v)).getVector());
1760 Real epsn = std::min(this->scale_*eps,this->min_diff_);
1761 for (
int i = 0; i < n; i++ )
1763 if ( ((*ex)[i] <= this->x_lo_[i]+epsn) ||
1764 ((*ex)[i] >= this->x_up_[i]-epsn) )
1771 void pruneLowerActive(::ROL::Vector<Real> &v, const ::ROL::Vector<Real> &g, const ::ROL::Vector<Real> &x, Real eps)
1773 Teuchos::RCP<const std::vector<Real> > ex =
1774 (Teuchos::dyn_cast<::ROL::StdVector<Real> >(
const_cast<::ROL::Vector<Real> &
>(x))).getVector();
1775 Teuchos::RCP<const std::vector<Real> > eg =
1776 (Teuchos::dyn_cast<::ROL::StdVector<Real> >(
const_cast<::ROL::Vector<Real> &
>(g))).getVector();
1777 Teuchos::RCP<std::vector<Real> > ev =
1778 Teuchos::rcp_const_cast<std::vector<Real> >((Teuchos::dyn_cast<::ROL::StdVector<Real> >(v)).getVector());
1779 Real epsn = std::min(this->scale_*eps,this->min_diff_);
1780 for (
int i = 0; i < n; i++ )
1782 if ( ((*ex)[i] <= this->x_lo_[i]+epsn && (*eg)[i] > 0.0) )
1789 void pruneUpperActive(::ROL::Vector<Real> &v, const ::ROL::Vector<Real> &g, const ::ROL::Vector<Real> &x, Real eps)
1791 Teuchos::RCP<const std::vector<Real> > ex =
1792 (Teuchos::dyn_cast<::ROL::StdVector<Real> >(
const_cast<::ROL::Vector<Real> &
>(x))).getVector();
1793 Teuchos::RCP<const std::vector<Real> > eg =
1794 (Teuchos::dyn_cast<::ROL::StdVector<Real> >(
const_cast<::ROL::Vector<Real> &
>(g))).getVector();
1795 Teuchos::RCP<std::vector<Real> > ev =
1796 Teuchos::rcp_const_cast<std::vector<Real> >((Teuchos::dyn_cast<::ROL::StdVector<Real> >(v)).getVector());
1797 Real epsn = std::min(this->scale_*eps,this->min_diff_);
1798 for (
int i = 0; i < n; i++ )
1800 if ( ((*ex)[i] >= this->x_up_[i]-epsn && (*eg)[i] < 0.0) )
1807 void pruneActive(::ROL::Vector<Real> &v, const ::ROL::Vector<Real> &g, const ::ROL::Vector<Real> &x, Real eps)
1809 Teuchos::RCP<const std::vector<Real> > ex =
1810 (Teuchos::dyn_cast<::ROL::StdVector<Real> >(
const_cast<::ROL::Vector<Real> &
>(x))).getVector();
1811 Teuchos::RCP<const std::vector<Real> > eg =
1812 (Teuchos::dyn_cast<::ROL::StdVector<Real> >(
const_cast<::ROL::Vector<Real> &
>(g))).getVector();
1813 Teuchos::RCP<std::vector<Real> > ev =
1814 Teuchos::rcp_const_cast<std::vector<Real> >((Teuchos::dyn_cast<::ROL::StdVector<Real> >(v)).getVector());
1815 Real epsn = std::min(this->scale_*eps,this->min_diff_);
1816 for (
int i = 0; i < n; i++ )
1818 if ( ((*ex)[i] <= this->x_lo_[i]+epsn && (*eg)[i] > 0.0) ||
1819 ((*ex)[i] >= this->x_up_[i]-epsn && (*eg)[i] < 0.0) )
1826 void setVectorToUpperBound( ::ROL::Vector<Real> &u )
1828 Teuchos::RCP<std::vector<Real> > us = Teuchos::rcp(
new std::vector<Real>(2,0.0) );
1829 us->assign(this->x_up_.begin(),this->x_up_.end());
1830 Teuchos::RCP<::ROL::Vector<Real> > up = Teuchos::rcp( new ::ROL::StdVector<Real>(us) );
1834 void setVectorToLowerBound( ::ROL::Vector<Real> &l )
1836 Teuchos::RCP<std::vector<Real> > ls = Teuchos::rcp(
new std::vector<Real>(2,0.0) );
1837 ls->assign(this->x_lo_.begin(),this->x_lo_.end());
1838 Teuchos::RCP<::ROL::Vector<Real> > lp = Teuchos::rcp( new ::ROL::StdVector<Real>(ls) );
1849 #endif // Xyce_N_ANP_ROL_DC_Optimization
Pure virtual class to augment a linear system.
const T & value(const ParameterBase &entity, const Descriptor &descriptor)
Returns the value of the parameter for the entity.
AnalysisManager & analysisManager_
Linear::System & linearSystem_