Xyce  6.1
N_ANP_ROL_DC_Optimization.h
Go to the documentation of this file.
1 //-----------------------------------------------------------------------------
2 // Copyright Notice
3 //
4 // Copyright 2002 Sandia Corporation. Under the terms
5 // of Contract DE-AC04-94AL85000 with Sandia Corporation, the U.S.
6 // Government retains certain rights in this software.
7 //
8 // Xyce(TM) Parallel Electrical Simulator
9 // Copyright (C) 2002-2014 Sandia Corporation
10 //
11 // This program is free software: you can redistribute it and/or modify
12 // it under the terms of the GNU General Public License as published by
13 // the Free Software Foundation, either version 3 of the License, or
14 // (at your option) any later version.
15 //
16 // This program is distributed in the hope that it will be useful,
17 // but WITHOUT ANY WARRANTY; without even the implied warranty of
18 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
19 // GNU General Public License for more details.
20 //
21 // You should have received a copy of the GNU General Public License
22 // along with this program. If not, see <http://www.gnu.org/licenses/>.
23 //-----------------------------------------------------------------------------
24 
25 //-----------------------------------------------------------------------------
26 // Filename : $RCSfile: N_ANP_ROL_DC_Optimization.h,v $
27 // Purpose :
28 // Special Notes :
29 // Creator :
30 // Creation Date :
31 //
32 // Revision Information:
33 // ---------------------
34 // Revision Number: $Revision: 1.5 $
35 // Revision Date : $Date: 2015/10/27 19:24:39 $
36 // Current Owner : $Author: tvrusso $
37 //-----------------------------------------------------------------------------
38 
39 #ifndef Xyce_N_ANP_ROL_DC_Optimizationh
40 #define Xyce_N_ANP_ROL_DC_Optimizationh
41 
42 #ifndef FD_HESSIAN
43 #define FD_HESSIAN 0 // 0 to set all hessvec to zero (Gauss-Newton Hessian); else FD Hessian
44 #endif
45 
46 #ifndef IDENTITY_PRECONDITIONER
47 #define IDENTITY_PRECONDITIONER 1
48 #endif
49 
50 #ifdef Xyce_ROL
51 
52 // ---------- ROL Includes ------------//
53 
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"
61 #include "ROL_XyceVector.hpp"
62 
63 #include <iostream>
64 #include <fstream>
65 #include <string>
66 #include <math.h>
67 
68 // ------------ Xyce Includes -----------//
69 
70 #include <Xyce_config.h>
71 #include <N_ANP_AnalysisManager.h>
72 #include <N_ANP_DCSweep.h>
73 #include <N_ANP_OutputMgrAdapter.h>
74 
75 // TT: perhaps not all of the following are needed
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>
88 #include <N_NLS_Manager.h>
89 #include <N_NLS_Sensitivity.h>
90 #include <N_NLS_fwd.h>
91 #include <N_NLS_NonLinearSolver.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>
96 #include <N_TIA_DataStore.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>
104 
105 
106 namespace Xyce {
107 namespace Analysis {
108 
109 enum sensDiffMode
110 {
111  SENS_FWD,
112  SENS_REV,
113  SENS_CNT,
115 };
116 
117 // **************** EqualityConstraint class *************************************//
118 
119 template<class Real>
120 class EqualityConstraint_ROL_DC : virtual public ::ROL::ParametrizedEqualityConstraint_SimOpt<Real>
121 {
122 private:
123  Real nu_, nc_, nz_;
124  Nonlinear::NonLinearSolver & nls_;
125  Linear::System & linearSystem_;
126  Loader::NonlinearEquationLoader & nEqLoader_;
127 
128  bool forceFD_;
129  Linear::Vector * origFVectorPtr_;
130  Linear::Vector * pertFVectorPtr_;
131  Linear::Vector * origQVectorPtr_;
132  Linear::Vector * pertQVectorPtr_;
133  Linear::Vector * origBVectorPtr_;
134  Linear::Vector * pertBVectorPtr_;
135 
136  void setControlParams(const ::ROL::Vector<Real> &z)
137  {
138  Teuchos::RCP<const std::vector<Real> > zp = (Teuchos::dyn_cast<::ROL::StdVector<Real> >(const_cast<::ROL::Vector<Real> &>(z))).getVector();
139  // Set parameters to values given by vector z
140  std::string paramName;
141  Real paramValue;
142  for (int i=0;i<rolSweep_.numParams_;i++)
143  {
144  paramName = rolSweep_.paramNameVec_[i];
145  paramValue = (*zp)[i];
146  //std::cout << "Setting " << paramName << " to " << paramValue << std::endl;
147  nEqLoader_.setParam(paramName,paramValue);
148  }
149  }
150 
151  void setUncertainParams()
152  {
153  //Loader::NonlinearEquationLoader & nEqLoader = analysisManager_.getNonlinearEquationLoader();
154  // Set uncertain parameters
155  const std::vector<Real> uparam = ::ROL::ParametrizedEqualityConstraint_SimOpt<Real>::getParameter();
156  std::string paramName;
157  Real paramValue;
158  for (int i=0;i<rolSweep_.uncertainParams_.size();i++)
159  {
160  paramName = rolSweep_.uncertainParams_[i];
161  paramValue = uparam[i];
162  //std::cout << "Setting " << paramName << " to " << paramValue << std::endl;
163  nEqLoader_.setParam(paramName,paramValue);
164  }
165  }
166 
167 protected:
168  ROL & rolSweep_;
169  AnalysisManager & analysisManager_;
170 
171 public:
172  EqualityConstraint_ROL_DC(
173  Real nu, // # of solution variables
174  Real nc, // # of constraint equations
175  Real nz, // # of optimization variables
176  AnalysisManager & analysisManager,
177  Loader::NonlinearEquationLoader & loader,
178  Nonlinear::NonLinearSolver & nls,
179  Linear::System & linearSystem,
180  ROL & rolSweep)
181  : nu_(nu),
182  nc_(nc),
183  nz_(nz),
184  forceFD_(false),
185  analysisManager_(analysisManager),
186  nEqLoader_(loader),
187  rolSweep_(rolSweep),
188  nls_(nls),
189  linearSystem_(linearSystem),
190  origFVectorPtr_(0),
191  pertFVectorPtr_(0),
192  origQVectorPtr_(0),
193  pertQVectorPtr_(0),
194  origBVectorPtr_(0),
195  pertBVectorPtr_(0)
196  {}
197 
198  void value(::ROL::Vector<Real> &c, const ::ROL::Vector<Real> &u, const ::ROL::Vector<Real> &z, Real &tol)
199  {
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();
202 
203  setControlParams(z);
204  setUncertainParams();
205 
206  bool success;
207 
208  bool vstatus = nEqLoader_.getVoltageLimiterStatus();
209  nEqLoader_.setVoltageLimiterStatus(false);// turns off voltage limiting
210 
211  for (int i=0;i<nc_;i++)
212  {
213  //std::string par = "V2";
214  //std::cout << (*(*up)[i])[0][1] << std::endl;
215  //nEqLoader_.setParam(par,(*(*up)[i])[0][0]);
216  rolSweep_.setSweepValue(i); // update source term; note: more direct approach would be to use setParam
217 
218  //analysisManager_.getDataStore()->setZeroHistory(); // not necessary
219 
220  // Set next solution vector pointer to the values given in vector u[i]
221  // Linear::Vector * temp = (*up)[i].getRawPtr();
222  // success = analysisManager_.getDataStore()->setNextSolVectorPtr(temp); // this causes segfault
223 
224  // Using newly added function in DataStore class
225  success = analysisManager_.getDataStore()->setNextSolVectorPtr(*(*up)[i]);
226 
227  // Load rhs
228  nEqLoader_.loadRHS();
229 
230  // Get rhs
231  *((*cp)[i]) = *(linearSystem_.getRHSVector());
232  //(*cp)[i]->printPetraObject(Xyce::dout());
233  }
234  nEqLoader_.setVoltageLimiterStatus(vstatus);
235 
236 
237  }
238 
239  void solve(::ROL::Vector<Real> &u, const ::ROL::Vector<Real> &z, Real &tol)
240  {
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());
242 
243  bool success;
244  setControlParams(z);
245  setUncertainParams();
246 
247  // Perform the sweep to obtain vector of solutions
248  success = rolSweep_.doInit() && rolSweep_.doLoopProcess() && rolSweep_.doFinish();
249 
250  for (int i=0;i<nc_;i++)
251  {
252  *((*up)[i]) = *(rolSweep_.solutionPtrVector_[i]);
253  //(*up)[i]->printPetraObject(Xyce::dout());
254  }
255 
256  }
257 
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)
259  {
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();
263 
264  bool success;
265  setControlParams(z);
266  setUncertainParams();
267 
268  bool vstatus = nEqLoader_.getVoltageLimiterStatus();
269  nEqLoader_.setVoltageLimiterStatus(false);// turns off voltage limiting
270 
271  for (int i=0;i<nc_;i++)
272  {
273  // Set next solution vector pointer to the values given in vector u[i]
274  // Linear::Vector * temp = (*up)[i].getRawPtr();
275  // success = analysisManager_.getDataStore()->setNextSolVectorPtr(temp); //segfault
276 
277  // Using newly added function in DataStore class
278  success = analysisManager_.getDataStore()->setNextSolVectorPtr(*(*up)[i]);
279 
280  //rolSweep_.setSweepValue(i);// not needed here
281  // Load rhs: appears to be necessary
282  nEqLoader_.loadRHS();
283 
284  // Apply Jacobian: not working (HB analysis?)
285  // success = nls_.applyJacobian( dynamic_cast< const Linear::Vector& > (const_cast< const Linear::MultiVector& >( *((*vp)[i]) ) ), dynamic_cast< Linear::Vector &>( *((*jvp)[i]) ) );
286 
287  // Load Jacobian
288  success = nEqLoader_.loadJacobian();
289 
290  // Get Jacobian
291  Linear::Matrix & Jac = *(linearSystem_.getJacobianMatrix());
292  Jac.scale(-1.0);
293 
294  // Apply Jacobian
295  Jac.matvec(false, const_cast< const Linear::Vector& >( *((*vp)[i]) ), *((*jvp)[i]) );
296 
297  // // Test
298  // std::cout << "Jac_1 Success = " << success << std::endl;
299  // (*jvp)[i]->printPetraObject(Xyce::dout());
300 
301  }
302  nEqLoader_.setVoltageLimiterStatus(vstatus);
303 
304 
305  }
306 
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)
308  {
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();
312 
313  bool success;
314  setControlParams(z);
315  setUncertainParams();
316 
317  bool vstatus = nEqLoader_.getVoltageLimiterStatus();
318  nEqLoader_.setVoltageLimiterStatus(false);// turns off voltage limiting
319 
320  for (int i=0;i<nc_;i++)
321  {
322  // Set next solution vector pointer to the values given in vector u[i]
323  // Linear::Vector * temp = (*up)[i].getRawPtr();
324  // success = analysisManager_.getDataStore()->setNextSolVectorPtr(temp); //segfault
325 
326  // Using newly added function in DataStore class
327  success = analysisManager_.getDataStore()->setNextSolVectorPtr(*(*up)[i]);
328 
329  //rolSweep_.setSweepValue(i);// not needed here
330  // Load rhs: appears to be necessary
331  nEqLoader_.loadRHS();
332 
333  // Load Jacobian
334  success = nEqLoader_.loadJacobian();
335 
336  // Get Jacobian
337  Linear::Matrix & Jac = *(linearSystem_.getJacobianMatrix());
338  Jac.scale(-1.0);
339 
340  // Apply Jacobian transpose
341  // Jac.setUseTranspose(true);
342  Jac.matvec(true, const_cast< const Linear::Vector& >( *((*vp)[i]) ), *((*jvp)[i]) );
343 
344  }
345  nEqLoader_.setVoltageLimiterStatus(vstatus);
346 
347 
348  }
349 
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)
351  {
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();
355 
356  bool success;
357  setControlParams(z);
358  setUncertainParams();
359 
360  bool vstatus = nEqLoader_.getVoltageLimiterStatus();
361  nEqLoader_.setVoltageLimiterStatus(false);// turns off voltage limiting
362 
363  for (int i=0;i<nc_;i++)
364  {
365  // Set next solution vector pointer to the values given in vector u[i]
366  // Linear::Vector * temp = (*up)[i].getRawPtr();
367  // success = analysisManager_.getDataStore()->setNextSolVectorPtr(temp); //segfault
368 
369  // Using newly added function in DataStore class
370  success = analysisManager_.getDataStore()->setNextSolVectorPtr(*(*up)[i]);
371 
372  //rolSweep_.setSweepValue(i);// not needed here
373  // Load rhs: appears to be necessary
374  nEqLoader_.loadRHS();
375 
376  // Load Jacobian
377  success = nEqLoader_.loadJacobian();
378 
379  // save rhs and Newton vectors
380  Linear::Vector * savedRHSVectorPtr_ = nls_.rhsVectorPtr_;
381  Linear::Vector * savedNewtonVectorPtr_ = nls_.NewtonVectorPtr_;
382  //savedRHSVectorPtr_->update(1.0, *(nls_.rhsVectorPtr_),0.0);
383  //savedNewtonVectorPtr_->update(1.0, *(nls_.NewtonVectorPtr_),0.0);
384 
385  // Linear::Vector * temp1 = (*vp)[i].getRawPtr();
386  // nls_.rhsVectorPtr_->update(-1.0, *temp1, 0.0);
387  *(nls_.rhsVectorPtr_) = *(*vp)[i];
388  nls_.rhsVectorPtr_->scale(-1.0); // solver expects negative rhs
389 
390  int status = nls_.lasSolverPtr_->solve(false);
391  if (status!=0)
392  {
393  std::string msg("Sensitivity::solveAdjoint. Solver failed\n");
394  N_ERH_ErrorMgr::report (N_ERH_ErrorMgr::DEV_FATAL, msg);
395  }
396 
397  //(*jvp)[i]->update(1.0, *(nls_.NewtonVectorPtr_),0.0);
398  *((*jvp)[i]) = *(nls_.NewtonVectorPtr_);
399 
400  // Restore the RHS and Newton vectors.
401  nls_.rhsVectorPtr_->update(1.0, *(savedRHSVectorPtr_),0.0);
402  nls_.NewtonVectorPtr_->update(1.0, *(savedNewtonVectorPtr_),0.0);
403 
404  }
405  nEqLoader_.setVoltageLimiterStatus(vstatus);
406 
407  }
408 
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)
410  {
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();
415 
416  bool success;
417  setControlParams(z);
418  setUncertainParams();
419 
420  bool vstatus = nEqLoader_.getVoltageLimiterStatus();
421  nEqLoader_.setVoltageLimiterStatus(false);// turns off voltage limiting
422 
423  for (int i=0;i<nc_;i++)
424  {
425  // Set next solution vector pointer to the values given in vector u[i]
426  // Linear::Vector * temp = (*up)[i].getRawPtr();
427  // success = analysisManager_.getDataStore()->setNextSolVectorPtr(temp); //segfault
428 
429  // Using newly added function in DataStore class
430  success = analysisManager_.getDataStore()->setNextSolVectorPtr(*(*up)[i]);
431 
432  //rolSweep_.setSweepValue(i);// not needed here
433  // Load rhs: appears to be necessary
434  nEqLoader_.loadRHS();
435 
436  // Load Jacobian
437  success = nEqLoader_.loadJacobian();
438 
439  // save rhs and Newton vectors
440  Linear::Vector * savedRHSVectorPtr_ = nls_.rhsVectorPtr_;
441  Linear::Vector * savedNewtonVectorPtr_ = nls_.NewtonVectorPtr_;
442  //savedRHSVectorPtr_->update(1.0, *(nls_.rhsVectorPtr_),0.0);
443  //savedNewtonVectorPtr_->update(1.0, *(nls_.NewtonVectorPtr_),0.0);
444 
445  // Linear::Vector * temp1 = (*vp)[i].getRawPtr();
446  // nls_.rhsVectorPtr_->update(-1.0, *temp1, 0.0);
447  *(nls_.rhsVectorPtr_) = *(*vp)[i];
448  nls_.rhsVectorPtr_->scale(-1.0);
449 
450  int status = nls_.lasSolverPtr_->solveTranspose(false);// some reusefactors
451  if (status!=0)
452  {
453  std::string msg("Sensitivity::solveAdjoint. Solver failed\n");
454  N_ERH_ErrorMgr::report (N_ERH_ErrorMgr::DEV_FATAL, msg);
455  }
456 
457  //(*jvp)[i]->update(1.0, *(nls_.NewtonVectorPtr_),0.0);
458  *((*jvp)[i]) = *(nls_.NewtonVectorPtr_);
459 
460  // Restore the RHS and Newton vectors.
461  nls_.rhsVectorPtr_->update(1.0, *(savedRHSVectorPtr_),0.0);
462  nls_.NewtonVectorPtr_->update(1.0, *(savedNewtonVectorPtr_),0.0);
463 
464  }
465  nEqLoader_.setVoltageLimiterStatus(true);
466 
467  }
468 
469 
470 
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)
472 {
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();
476 
477  bool success;
478  jv.zero();
479  setControlParams(z);
480  setUncertainParams();
481 
482  bool vstatus = nEqLoader_.getVoltageLimiterStatus();
483  nEqLoader_.setVoltageLimiterStatus(false);// turns off voltage limiting
484 
485  // loop over constraint equations
486  for (int k=0;k<nc_;k++)
487  {
488  rolSweep_.setSweepValue(k); // update source term; Important!
489 
490  // Note: perhaps all source terms need be updated
491 
492  // Set next solution vector pointer to the values given in vector u[i]
493  // Linear::Vector * temp = (*up)[k].getRawPtr();
494  // success = analysisManager_.getDataStore()->setNextSolVectorPtr(temp); //segfault
495 
496  // Using newly added function in DataStore class
497  success = analysisManager_.getDataStore()->setNextSolVectorPtr(*(*up)[k]);
498 
499  // Load rhs
500  success = nEqLoader_.loadRHS();
501 
502  int iparam;
503  std::string msg;
504 
505  analysisManager_.getDataStore()->paramOrigVals_.clear();
506 
507  // it is necessary to load the Jacobian here to make sure we have the most
508  // up-to-date matrix. The Jacobian is not loaded for the final
509  // evaluation of the residual in the Newton solve.
510  success = nEqLoader_.loadJacobian ();
511 
512  // Loop over the vector of parameters. For each parameter, find the
513  // device entity (a model or an instance) which corresponds to it, and
514  // perform the finite difference calculation.
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 )
520  {
521  std::string paramName(*iterParam);
522  // std::cout << "Computing derivative wrt " << rolSweep_.paramNameVec_[iparam] << std::endl;
523  // get the original value of this parameter, to be used for scaling and/or
524  // numerical derivatives later.
525  double paramOrig = 0.0;
526  bool found = nEqLoader_.getParamAndReduce(analysisManager_.getPDSManager()->getPDSComm()->comm(),paramName, paramOrig);
527 
528  if (!found)
529  {
530  std::string msg("EqualityConstraint::applyJacobian_2: cannot find parameter ");
531  msg += paramName;
532  N_ERH_ErrorMgr::report (N_ERH_ErrorMgr::DEV_FATAL, msg);
533  }
534 
535  analysisManager_.getDataStore()->paramOrigVals_.push_back(paramOrig);
536 
537  // check if derivative is available analytically. If not, take FD.
538  bool analyticAvailable = false;
539  if (!forceFD_)
540  {
541  analyticAvailable = nEqLoader_.analyticSensitivitiesAvailable (paramName);
542  }
543  if (analyticAvailable)
544  {
545  //std::cout << "Analytic derivs available" << std::endl;
546  std::vector<double> dfdpVec;
547  std::vector<double> dqdpVec;
548  std::vector<double> dbdpVec;
549 
550  std::vector<int> FindicesVec;
551  std::vector<int> QindicesVec;
552  std::vector<int> BindicesVec;
553 
554  nEqLoader_.getAnalyticSensitivities(paramName,dfdpVec,dqdpVec,dbdpVec,FindicesVec, QindicesVec, BindicesVec);
555 
556  rolSweep_.mydfdpPtrVector_[iparam]->putScalar(0.0);
557  rolSweep_.mydqdpPtrVector_[iparam]->putScalar(0.0);
558  rolSweep_.mydbdpPtrVector_[iparam]->putScalar(0.0);
559 
560  int Fsize=FindicesVec.size();
561  for (int i=0;i<Fsize;++i)
562  {
563  Linear::Vector & dfdpRef = *(rolSweep_.mydfdpPtrVector_[iparam]);
564  dfdpRef[FindicesVec[i]] += dfdpVec[i];
565  }
566 
567  int Qsize=QindicesVec.size();
568  for (int i=0;i<Qsize;++i)
569  {
570  Linear::Vector & dqdpRef = *(rolSweep_.mydqdpPtrVector_[iparam]);
571  dqdpRef[QindicesVec[i]] += dqdpVec[i];
572  }
573 
574  int Bsize=BindicesVec.size();
575  for (int i=0;i<Bsize;++i)
576  {
577  Linear::Vector & dbdpRef = *(rolSweep_.mydbdpPtrVector_[iparam]);
578  dbdpRef[BindicesVec[i]] += dbdpVec[i];
579  }
580 
581  rolSweep_.mydfdpPtrVector_[iparam]->fillComplete();
582  rolSweep_.mydqdpPtrVector_[iparam]->fillComplete();
583  rolSweep_.mydbdpPtrVector_[iparam]->fillComplete();
584 
585  if (DEBUG_NONLINEAR && isActive(Diag::SENS_SOLVER))
586  {
587  Xyce::dout() << *iterParam << ": ";
588  Xyce::dout().setf(std::ios::scientific);
589  Xyce::dout() << std::endl;
590 
591  for (int k1 = 0; k1 < analysisManager_.getDataStore()->solutionSize; ++k1)
592  {
593  Xyce::dout()
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]
597  <<std::endl;
598  }
599  }
600  }
601  else
602  {
603  // Xyce::dout() << "Analytical derivatives NOT available" << std::endl;
604 
605  if (DEBUG_NONLINEAR && isActive(Diag::SENS_SOLVER))
606  {
607  Xyce::dout() << std::endl << " Calculating numerical df/dp, dq/dp and db/dp for: ";
608  Xyce::dout() << *iterParam << std::endl;
609  }
610 
611  // save a copy of the DAE vectors
612  origFVectorPtr_ = linearSystem_.builder().createVector();
613  pertFVectorPtr_ = linearSystem_.builder().createVector();
614 
615  origQVectorPtr_ = linearSystem_.builder().createVector();
616  pertQVectorPtr_ = linearSystem_.builder().createVector();
617 
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);
623 
624  // now perturb the value of this parameter.
625  double sqrtEta_ = pow(10,int(log10(fabs(paramOrig))));; // TT: this should be parameter specific (see Sensitivity class)
626  double dp = sqrtEta_ * 1.e-5 * (1.0 + fabs(paramOrig));
627  double paramPerturbed = paramOrig;
628 
629  int difference=SENS_FWD; // TT: this should be optional
630  if (difference==SENS_FWD)
631  {
632  paramPerturbed += dp;
633  }
634  else if (difference==SENS_REV)
635  {
636  paramPerturbed -= dp;
637  }
638  else if (difference==SENS_CNT)
639  {
640  static std::string tmp = "difference=central not supported.\n";
641  N_ERH_ErrorMgr::report( N_ERH_ErrorMgr::USR_FATAL_0, tmp);
642  }
643  else
644  {
645  static std::string tmp = "difference not recognized!\n";
646  N_ERH_ErrorMgr::report( N_ERH_ErrorMgr::USR_FATAL_0, tmp);
647  }
648 
649  if (DEBUG_NONLINEAR && isActive(Diag::SENS_SOLVER))
650  {
651  int maxParamStringSize_ = 16; // TT: not sure what this should be
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
656  <<std::endl;
657  }
658 
659  nEqLoader_.setParam (paramName, paramPerturbed);
660 
661  // Now that the parameter has been perturbed,
662  // calculate the numerical derivative.
663 
664  rolSweep_.setSweepValue(k);// needed?
665  // Load F,Q and B.
666  nEqLoader_.loadRHS();
667 
668  // save the perturbed DAE vectors
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);
672 
673  // calculate the df/dp vector.
674  double rdp=1/dp;
675  rolSweep_.mydfdpPtrVector_[iparam]->putScalar(0.0);
676  rolSweep_.mydfdpPtrVector_[iparam]->addVec (+1.0, *(pertFVectorPtr_)); //+Fperturb
677  rolSweep_.mydfdpPtrVector_[iparam]->addVec (-1.0, *(origFVectorPtr_)); //-Forig
678  rolSweep_.mydfdpPtrVector_[iparam]->scale(rdp);
679 
680  // calculate the dq/dp vector.
681  rolSweep_.mydqdpPtrVector_[iparam]->putScalar(0.0);
682  rolSweep_.mydqdpPtrVector_[iparam]->addVec (+1.0, *(pertQVectorPtr_)); //+Fperturb
683  rolSweep_.mydqdpPtrVector_[iparam]->addVec (-1.0, *(origQVectorPtr_)); //-Forig
684  rolSweep_.mydqdpPtrVector_[iparam]->scale(rdp);
685 
686  // calculate the db/dp vector.
687  rolSweep_.mydbdpPtrVector_[iparam]->putScalar(0.0);
688  rolSweep_.mydbdpPtrVector_[iparam]->addVec (+1.0, *(pertBVectorPtr_)); //+Fperturb
689  rolSweep_.mydbdpPtrVector_[iparam]->addVec (-1.0, *(origBVectorPtr_)); //-Forig
690  rolSweep_.mydbdpPtrVector_[iparam]->scale(rdp);
691 
692  if (DEBUG_NONLINEAR && isActive(Diag::SENS_SOLVER))
693  {
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;
697 
698  for (int k1 = 0; k1 < analysisManager_.getDataStore()->solutionSize; ++k1)
699  {
700 
701  Xyce::dout()
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]
705  <<std::endl;
706  }
707 
708  Xyce::dout() << std::endl;
709  for (int k1 = 0; k1 < analysisManager_.getDataStore()->solutionSize; ++k1)
710  {
711  Xyce::dout()
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]
715  <<std::endl;
716  }
717 
718  Xyce::dout() << std::endl ;
719  for (int k1 = 0; k1 < analysisManager_.getDataStore()->solutionSize; ++k1)
720  {
721  Xyce::dout()
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]
725  <<std::endl;
726 
727  }
728 
729  // std::ostringstream filename;
730  // filename << netlistFilename_ << "_dfdp";
731  // filename << std::setw(3) << std::setfill('0') << iparam;
732  // filename << ".txt";
733  // dfdpPtrVector_[iparam]->writeToFile(const_cast<char *>(filename.str().c_str()));
734 
735  // filename.str("");
736  // filename << netlistFilename_ << "_fpert";
737  // filename << std::setw(3) << std::setfill('0') << iparam;
738  // filename << ".txt";
739  // pertFVectorPtr_->writeToFile(const_cast<char *>(filename.str().c_str()));
740 
741  // filename.str("");
742  // filename << netlistFilename_ << "_dqdp";
743  // filename << std::setw(3) << std::setfill('0') << iparam;
744  // filename << ".txt";
745  // dqdpPtrVector_[iparam]->writeToFile(const_cast<char *>(filename.str().c_str()));
746 
747  // filename.str("");
748  // filename << netlistFilename_ << "_qpert";
749  // filename << std::setw(3) << std::setfill('0') << iparam;
750  // filename << ".txt";
751  // pertQVectorPtr_->writeToFile(const_cast<char *>(filename.str().c_str()));
752 
753  // filename.str("");
754  // filename << netlistFilename_ << "_dbdp";
755  // filename << std::setw(3) << std::setfill('0') << iparam;
756  // filename << ".txt";
757  // dbdpPtrVector_[iparam]->writeToFile(const_cast<char *>(filename.str().c_str()));
758 
759  // filename.str("");
760  // filename << netlistFilename_ << "_bpert";
761  // filename << std::setw(3) << std::setfill('0') << iparam;
762  // filename << ".txt";
763  // pertBVectorPtr_->writeToFile(const_cast<char *>(filename.str().c_str()));
764  }
765 
766  // now reset the parameter and rhs to previous values.
767  nEqLoader_.setParam (paramName, paramOrig);
768 
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);
772  }
773 
774  // Now collect dfdp, dqdp, dbdp into dFdp (as done in obtainSensitivityResiduals in NoTimeIntegration, that is we ignore dqdp)
775 
776  Linear::Vector & sensRHS = *(rolSweep_.mysensRHSPtrVector_[iparam]);
777  Linear::Vector & dfdpRef = *(rolSweep_.mydfdpPtrVector_[iparam]);
778  Linear::Vector & dbdpRef = *(rolSweep_.mydbdpPtrVector_[iparam]);
779 
780  sensRHS.linearCombo(0.0,sensRHS,+1.0,dfdpRef);
781  sensRHS.linearCombo(1.0,sensRHS,-1.0,dbdpRef);
782 
783  sensRHS.scale(-1.0);
784 
785  // multiply by v[iparam]
786  sensRHS.scale( (*vp)[iparam] );
787 
788  // add to jvp[k]
789  (*jvp)[k]->addVec(1.0,sensRHS);
790 
791  }// end of parameters loop
792 
793  }// end of constraints loop
794  nEqLoader_.setVoltageLimiterStatus(vstatus);
795 
796 }
797 
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)
799  {
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());
803 
804  bool success;
805  jv.zero();
806  setControlParams(z);
807  setUncertainParams();
808 
809  bool vstatus = nEqLoader_.getVoltageLimiterStatus();
810  nEqLoader_.setVoltageLimiterStatus(false);// turns off voltage limiting
811 
812  // loop over constraint equations
813  for (int k=0;k<nc_;k++)
814  {
815  rolSweep_.setSweepValue(k); // update source term; Important!
816 
817  // Set next solution vector pointer to the values given in vector u[i]
818  // Linear::Vector * temp = (*up)[k].getRawPtr();
819  // success = analysisManager_.getDataStore()->setNextSolVectorPtr(temp); //segfault
820 
821  // Using newly added function in DataStore class
822  success = analysisManager_.getDataStore()->setNextSolVectorPtr(*(*up)[k]);
823 
824  // Load rhs
825  success = nEqLoader_.loadRHS();
826 
827  int iparam;
828  std::string msg;
829 
830  analysisManager_.getDataStore()->paramOrigVals_.clear();
831 
832  // it is necessary to load the Jacobian here to make sure we have the most
833  // up-to-date matrix. The Jacobian is not loaded for the final
834  // evaluation of the residual in the Newton solve.
835  nEqLoader_.loadJacobian ();
836 
837  // Loop over the vector of parameters. For each parameter, find the
838  // device entity (a model or an instance) which corresponds to it, and
839  // perform the finite difference calculation.
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 )
845  {
846  std::string paramName(*iterParam);
847  // std::cout << "Computing derivative wrt " << rolSweep_.paramNameVec_[iparam] << std::endl;
848  // get the original value of this parameter, to be used for scaling and/or
849  // numerical derivatives later.
850  double paramOrig = 0.0;
851  bool found = nEqLoader_.getParamAndReduce(analysisManager_.getPDSManager()->getPDSComm()->comm(),paramName, paramOrig);
852 
853  if (!found)
854  {
855  std::string msg("EqualityConstraint::applyJacobian_2: cannot find parameter ");
856  msg += paramName;
857  N_ERH_ErrorMgr::report (N_ERH_ErrorMgr::DEV_FATAL, msg);
858  }
859 
860  analysisManager_.getDataStore()->paramOrigVals_.push_back(paramOrig);
861 
862  // check if derivative is available analytically. If not, take FD.
863  bool analyticAvailable = false;
864  if (!forceFD_)
865  {
866  analyticAvailable = nEqLoader_.analyticSensitivitiesAvailable (paramName);
867  }
868  if (analyticAvailable)
869  {
870  std::vector<double> dfdpVec;
871  std::vector<double> dqdpVec;
872  std::vector<double> dbdpVec;
873 
874  std::vector<int> FindicesVec;
875  std::vector<int> QindicesVec;
876  std::vector<int> BindicesVec;
877 
878  nEqLoader_.getAnalyticSensitivities(paramName,dfdpVec,dqdpVec,dbdpVec,FindicesVec, QindicesVec, BindicesVec);
879 
880  rolSweep_.mydfdpPtrVector_[iparam]->putScalar(0.0);
881  rolSweep_.mydqdpPtrVector_[iparam]->putScalar(0.0);
882  rolSweep_.mydbdpPtrVector_[iparam]->putScalar(0.0);
883 
884  int Fsize=FindicesVec.size();
885  for (int i=0;i<Fsize;++i)
886  {
887  Linear::Vector & dfdpRef = *(rolSweep_.mydfdpPtrVector_[iparam]);
888  dfdpRef[FindicesVec[i]] += dfdpVec[i];
889  }
890 
891  int Qsize=QindicesVec.size();
892  for (int i=0;i<Qsize;++i)
893  {
894  Linear::Vector & dqdpRef = *(rolSweep_.mydqdpPtrVector_[iparam]);
895  dqdpRef[QindicesVec[i]] += dqdpVec[i];
896  }
897 
898  int Bsize=BindicesVec.size();
899  for (int i=0;i<Bsize;++i)
900  {
901  Linear::Vector & dbdpRef = *(rolSweep_.mydbdpPtrVector_[iparam]);
902  dbdpRef[BindicesVec[i]] += dbdpVec[i];
903  }
904 
905  rolSweep_.mydfdpPtrVector_[iparam]->fillComplete();
906  rolSweep_.mydqdpPtrVector_[iparam]->fillComplete();
907  rolSweep_.mydbdpPtrVector_[iparam]->fillComplete();
908 
909  if (DEBUG_NONLINEAR && isActive(Diag::SENS_SOLVER))
910  {
911  Xyce::dout() << *iterParam << ": ";
912  Xyce::dout().setf(std::ios::scientific);
913  Xyce::dout() << std::endl;
914 
915  for (int k1 = 0; k1 < analysisManager_.getDataStore()->solutionSize; ++k1)
916  {
917  Xyce::dout()
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]
921  <<std::endl;
922  }
923  }
924  }
925  else
926  {
927  // std::cout << "Analytical derivatives NOT available" << std::endl;
928 
929  if (DEBUG_NONLINEAR && isActive(Diag::SENS_SOLVER))
930  {
931  Xyce::dout() << std::endl << " Calculating numerical df/dp, dq/dp and db/dp for: ";
932  Xyce::dout() << *iterParam << std::endl;
933  }
934 
935  // save a copy of the DAE vectors
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);
939 
940  // now perturb the value of this parameter.
941  double sqrtEta_ = pow(10,int(log10(fabs(paramOrig))));; // TT: this should be parameter specific
942  double dp = sqrtEta_ * 1.e-4 * (1.0 + fabs(paramOrig));
943  double paramPerturbed = paramOrig;
944 
945  int difference=SENS_FWD; // TT: this should be optional
946  if (difference==SENS_FWD)
947  {
948  paramPerturbed += dp;
949  }
950  else if (difference==SENS_REV)
951  {
952  paramPerturbed -= dp;
953  }
954  else if (difference==SENS_CNT)
955  {
956  static std::string tmp = "difference=central not supported.\n";
957  N_ERH_ErrorMgr::report( N_ERH_ErrorMgr::USR_FATAL_0, tmp);
958  }
959  else
960  {
961  static std::string tmp = "difference not recognized!\n";
962  N_ERH_ErrorMgr::report( N_ERH_ErrorMgr::USR_FATAL_0, tmp);
963  }
964 
965  if (DEBUG_NONLINEAR && isActive(Diag::SENS_SOLVER))
966  {
967  int maxParamStringSize_ = 16; // TT: not sure what this should be
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
972  <<std::endl;
973  }
974 
975  nEqLoader_.setParam (paramName, paramPerturbed);
976 
977  // Now that the parameter has been perturbed,
978  // calculate the numerical derivative.
979 
980  rolSweep_.setSweepValue(k);// needed?
981  // Load F,Q and B.
982  nEqLoader_.loadRHS();
983 
984  // save the perturbed DAE vectors
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);
988 
989  // calculate the df/dp vector.
990  double rdp=1/dp;
991  rolSweep_.mydfdpPtrVector_[iparam]->putScalar(0.0);
992  rolSweep_.mydfdpPtrVector_[iparam]->addVec (+1.0, *(pertFVectorPtr_)); //+Fperturb
993  rolSweep_.mydfdpPtrVector_[iparam]->addVec (-1.0, *(origFVectorPtr_)); //-Forig
994  rolSweep_.mydfdpPtrVector_[iparam]->scale(rdp);
995 
996  // calculate the dq/dp vector.
997  rolSweep_.mydqdpPtrVector_[iparam]->putScalar(0.0);
998  rolSweep_.mydqdpPtrVector_[iparam]->addVec (+1.0, *(pertQVectorPtr_)); //+Fperturb
999  rolSweep_.mydqdpPtrVector_[iparam]->addVec (-1.0, *(origQVectorPtr_)); //-Forig
1000  rolSweep_.mydqdpPtrVector_[iparam]->scale(rdp);
1001 
1002  // calculate the db/dp vector.
1003  rolSweep_.mydbdpPtrVector_[iparam]->putScalar(0.0);
1004  rolSweep_.mydbdpPtrVector_[iparam]->addVec (+1.0, *(pertBVectorPtr_)); //+Fperturb
1005  rolSweep_.mydbdpPtrVector_[iparam]->addVec (-1.0, *(origBVectorPtr_)); //-Forig
1006  rolSweep_.mydbdpPtrVector_[iparam]->scale(rdp);
1007 
1008  if (DEBUG_NONLINEAR && isActive(Diag::SENS_SOLVER))
1009  {
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;
1013 
1014  for (int k1 = 0; k1 < analysisManager_.getDataStore()->solutionSize; ++k1)
1015  {
1016 
1017  Xyce::dout()
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]
1021  <<std::endl;
1022  }
1023 
1024  Xyce::dout() << std::endl;
1025  for (int k1 = 0; k1 < analysisManager_.getDataStore()->solutionSize; ++k1)
1026  {
1027  Xyce::dout()
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]
1031  <<std::endl;
1032  }
1033 
1034  Xyce::dout() << std::endl ;
1035  for (int k1 = 0; k1 < analysisManager_.getDataStore()->solutionSize; ++k1)
1036  {
1037  Xyce::dout()
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]
1041  <<std::endl;
1042 
1043  }
1044 
1045  // std::ostringstream filename;
1046  // filename << netlistFilename_ << "_dfdp";
1047  // filename << std::setw(3) << std::setfill('0') << iparam;
1048  // filename << ".txt";
1049  // dfdpPtrVector_[iparam]->writeToFile(const_cast<char *>(filename.str().c_str()));
1050 
1051  // filename.str("");
1052  // filename << netlistFilename_ << "_fpert";
1053  // filename << std::setw(3) << std::setfill('0') << iparam;
1054  // filename << ".txt";
1055  // pertFVectorPtr_->writeToFile(const_cast<char *>(filename.str().c_str()));
1056 
1057  // filename.str("");
1058  // filename << netlistFilename_ << "_dqdp";
1059  // filename << std::setw(3) << std::setfill('0') << iparam;
1060  // filename << ".txt";
1061  // dqdpPtrVector_[iparam]->writeToFile(const_cast<char *>(filename.str().c_str()));
1062 
1063  // filename.str("");
1064  // filename << netlistFilename_ << "_qpert";
1065  // filename << std::setw(3) << std::setfill('0') << iparam;
1066  // filename << ".txt";
1067  // pertQVectorPtr_->writeToFile(const_cast<char *>(filename.str().c_str()));
1068 
1069  // filename.str("");
1070  // filename << netlistFilename_ << "_dbdp";
1071  // filename << std::setw(3) << std::setfill('0') << iparam;
1072  // filename << ".txt";
1073  // dbdpPtrVector_[iparam]->writeToFile(const_cast<char *>(filename.str().c_str()));
1074 
1075  // filename.str("");
1076  // filename << netlistFilename_ << "_bpert";
1077  // filename << std::setw(3) << std::setfill('0') << iparam;
1078  // filename << ".txt";
1079  // pertBVectorPtr_->writeToFile(const_cast<char *>(filename.str().c_str()));
1080  }
1081 
1082  // now reset the parameter and rhs to previous values.
1083  nEqLoader_.setParam (paramName, paramOrig);
1084 
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);
1088  }
1089 
1090  // Now collect dfdp, dqdp, dbdp into dFdp (as done in obtainSensitivityResiduals in NoTimeIntegration, that is we ignore dqdp)
1091 
1092  Linear::Vector & sensRHS = *(rolSweep_.mysensRHSPtrVector_[iparam]);
1093  Linear::Vector & dfdpRef = *(rolSweep_.mydfdpPtrVector_[iparam]);
1094  Linear::Vector & dbdpRef = *(rolSweep_.mydbdpPtrVector_[iparam]);
1095 
1096  sensRHS.linearCombo(0.0,sensRHS,+1.0,dfdpRef);
1097  sensRHS.linearCombo(1.0,sensRHS,-1.0,dbdpRef);
1098 
1099  sensRHS.scale(-1.0);
1100 
1101  // Take dot product with k-th vector in v
1102  Real dot = sensRHS.dotProduct(*((*vp)[k]));
1103 
1104  // Add to jvp
1105  (*jvp)[iparam] += dot;
1106 
1107  }// end of parameters loop
1108 
1109  }// end of constraints loop
1110  nEqLoader_.setVoltageLimiterStatus(vstatus);
1111 
1112  }
1113 
1114 #if FD_HESSIAN==0
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)
1116  {
1117  ahwv.zero();
1118  }
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)
1120  {
1121  ahwv.zero();
1122  }
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)
1124  {
1125  ahwv.zero();
1126  }
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)
1128  {
1129  ahwv.zero();
1130  }
1131 #endif
1132 
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)
1135  {
1136  pv.set(v.dual());
1137  }
1138 #endif
1139 
1140 };
1141 
1142 template<class Real>
1143 class Objective_DC_L2Norm : public ::ROL::ParametrizedObjective_SimOpt<Real>
1144 {
1145 private:
1146  Real alpha_; // Penalty Parameter
1147  int ns_,nz_;
1148  Teuchos::RCP<std::vector<Real> > Imeas_;
1149  int solindex_;
1150 
1151 public:
1152 
1153  Objective_DC_L2Norm(Real alpha, int ns, int nz)
1154  {
1155  ns_ = ns; // number of sources (same as number of constraints)
1156  nz_ = nz; // number of optimization variables
1157  alpha_ = alpha;
1158  solindex_ = 0; // 0 for diode_forTimur, 1 for diopde, 2 for diode_direct, 5 for vbic (collector current)
1159  Imeas_ = Teuchos::rcp(new std::vector<Real>(ns_,0.0));
1160  Real temp;
1161  std::ifstream measurements("measurementsDC.txt");
1162  std::cout << "Reading measurements" << std::endl;
1163  if (measurements.is_open())
1164  {
1165  for (int i=0;i<ns_;i++)
1166  {
1167  measurements >> temp;
1168  measurements >> temp;
1169  (*Imeas_)[i] = temp;
1170  }
1171  }
1172  else
1173  {
1174  std::cout << "Could not open measurements file" << std::endl;
1175  return;
1176  }
1177  measurements.close();
1178 
1179  }
1180 
1181  Real value( const ::ROL::Vector<Real> &u, const ::ROL::Vector<Real> &z, Real &tol )
1182  {
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();
1185  Real val = 0.0;
1186  Real current;
1187  // Here we need to know which index corresponds to the state variable we need (e.g., current);
1188  // solution index (solindex) is specified in constructor
1189  for (int i=0;i<ns_;i++)
1190  {
1191  current = (*(*up)[i])[solindex_]; // vector
1192  val += (current - (*Imeas_)[i])*(current - (*Imeas_)[i]);
1193  }
1194  return 0.5*val;
1195  }
1196 
1197  void gradient_1( ::ROL::Vector<Real> &g, const ::ROL::Vector<Real> &u, const ::ROL::Vector<Real> &z, Real &tol )
1198  {
1199  // Unwrap g
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());
1201  // Unwrap x
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();
1205  // COMPUTE GRADIENT WRT U
1206  g.zero();
1207  Real current;
1208  for (int i=0; i<ns_; i++)
1209  {
1210  current = (*(*up)[i])[solindex_]; // vector
1211  (*(*gup)[i])[solindex_] = (current-(*Imeas_)[i]);
1212  }
1213 
1214  }
1215 
1216  void gradient_2( ::ROL::Vector<Real> &g, const ::ROL::Vector<Real> &u, const ::ROL::Vector<Real> &z, Real &tol )
1217  {
1218  // Unwrap g
1219  Teuchos::RCP<std::vector<Real> > gzp = Teuchos::rcp_const_cast<std::vector<Real> >((Teuchos::dyn_cast<const ::ROL::StdVector<Real> >(g)).getVector());
1220  // Unwrap x
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();
1223  // COMPUTE GRADIENT WRT Z
1224  for (int i=0; i<nz_; i++)
1225  {
1226  (*gzp)[i] = 0.0;
1227  }
1228  }
1229 
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 )
1232  {
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());
1234  // Unwrap v
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();
1236 
1237  hv.zero();
1238  // COMPUTE GRADIENT WRT U
1239  for(int i=0; i<ns_; i++)
1240  {
1241  (*(*hvup)[i])[solindex_] = (*(*vup)[i])[solindex_]; // picking the same component of vup[i] that corresponds to current in up[i]
1242  }
1243  }
1244 
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 )
1247  {
1248  hv.zero();
1249  }
1250 
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 )
1253  {
1254  hv.zero();
1255  }
1256 
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 )
1259  {
1260  hv.zero();
1261  }
1262 
1263 };
1264 
1265 template<class Real>
1266 class Objective_DC_AMP : virtual public ::ROL::ParametrizedObjective_SimOpt<Real>
1267 {
1268 private:
1269  int ns_,nz_,first_,middle_,last_;
1270  int solindex_;
1271 
1272 public:
1273 
1274  Objective_DC_AMP(int ns, int nz)
1275  {
1276  ns_ = ns; // number of sources (same as number of constraints)
1277  nz_ = nz; // number of optimization variables
1278  first_ = 0;
1279  middle_ = ns_/2;
1280  last_ = ns_-1;
1281  solindex_ = 3; // V(3)
1282  }
1283 
1284  Real value( const ::ROL::Vector<Real> &u, const ::ROL::Vector<Real> &z, Real &tol )
1285  {
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();
1288  Real val = 0.0;
1289  Real V3,t;
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++)
1296  {
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 );
1300  }
1301  return 0.5*val;
1302  }
1303 
1304  void gradient_1( ::ROL::Vector<Real> &g, const ::ROL::Vector<Real> &u, const ::ROL::Vector<Real> &z, Real &tol )
1305  {
1306  // Unwrap g
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());
1308  // Unwrap x
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();
1312  // COMPUTE GRADIENT WRT U
1313  g.zero();
1314  Real V3,t;
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++)
1323  {
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 );
1330  }
1331  }
1332 
1333  void gradient_2(::ROL::Vector<Real> &g, const ::ROL::Vector<Real> &u, const ::ROL::Vector<Real> &z, Real &tol )
1334  {
1335  // Unwrap g
1336  Teuchos::RCP<std::vector<Real> > gzp = Teuchos::rcp_const_cast<std::vector<Real> >((Teuchos::dyn_cast<const ::ROL::StdVector<Real> >(g)).getVector());
1337  // Unwrap x
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();
1340  // COMPUTE GRADIENT WRT Z
1341  for (int i=0; i<nz_; i++)
1342  {
1343  (*gzp)[i] = 0.0;
1344  }
1345  }
1346 
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 )
1349  {
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());
1351  // Unwrap v
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();
1354  hv.zero();
1355  Real V3,t;
1356  Real V3_0 = (*(*vup)[first_ ])[solindex_];
1357  Real V3_10 = (*(*vup)[middle_])[solindex_];
1358  Real V3_20 = (*(*vup)[last_ ])[solindex_];
1359  Real sum1 = 0.0;
1360  for(int i=1;i<ns_-1;i++)
1361  {
1362  sum1 += std::pow((-1.0 + i*2.0/(ns_-1)),2);
1363  }
1364  sum1 /= 4.0;
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++)
1369  {
1370  t = -1.0 + i*2.0/(ns_-1);
1371  V3 = (*(*vup)[i])[solindex_];
1372  (*(*hvup)[i])[solindex_] += V3
1373  + (-t/2.0)*V3_0
1374  - V3_10
1375  + ( t/2.0)*V3_20;
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;
1379  }
1380  }
1381 
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 )
1384  {
1385  hv.zero();
1386  }
1387 
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 )
1390  {
1391  hv.zero();
1392  }
1393 
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 )
1396  {
1397  hv.zero();
1398  }
1399 
1400 };
1401 
1402 template<class Real>
1403 class Penalty_DC_AMP : virtual public ::ROL::ParametrizedObjective_SimOpt<Real>
1404 {
1405 private:
1406  int ptype_; // penalty type
1407  Real alpha_; // penalty parameter
1408  Real epsilon_; // penalty parameter
1409  Real ampl_; // amplification parameter
1410  int ns_,nz_,first_,middle_,last_;
1411  int solindex_;
1412 
1413 public:
1414 
1415  Penalty_DC_AMP(int ptype, Real alpha, Real ampl, int ns, int nz)
1416  {
1417  ns_ = ns; // number of sources (same as number of constraints)
1418  nz_ = nz; // number of optimization variables
1419  first_ = 0;
1420  middle_ = ns_/2;
1421  last_ = ns_-1;
1422  ptype_ = ptype;
1423  alpha_ = alpha;
1424  solindex_ = 3; // V(3)
1425  epsilon_ = 10.0;
1426  ampl_ = ampl;
1427  }
1428 
1429  void switch_ptype(int ptype)
1430  {
1431  ptype_ = ptype;
1432  }
1433 
1434  Real value( const ::ROL::Vector<Real> &u, const ::ROL::Vector<Real> &z, Real &tol )
1435  {
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();
1438  Real val = 0;
1439  Real A = ((*(*up)[first_])[solindex_] - (*(*up)[last_])[solindex_] )/2.0; // negative gain
1440  if(ptype_==1)
1441  val += alpha_*(0.5+0.5*tanh(epsilon_*(ampl_ - A*A)));
1442  if(ptype_==2)
1443  val += alpha_*( (ampl_ - A*A > 0) ? ampl_-A*A : 0 );
1444  if(ptype_==3)
1445  val += alpha_*std::pow(( (ampl_ - A*A > 0) ? ampl_-A*A : 0 ),2);
1446  if(ptype_==4)
1447  val += alpha_*std::pow(( (A - ampl_ > 0) ? A - ampl_ : 0 ),2); // Moreau-Yosida penalty
1448  if(ptype_==5)
1449  val += alpha_*2.*A; // maximize gain
1450  return 0.5*val;
1451  }
1452 
1453  void gradient_1( ::ROL::Vector<Real> &g, const ::ROL::Vector<Real> &u, const ::ROL::Vector<Real> &z, Real &tol )
1454  {
1455  // Unwrap g
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());
1457  // Unwrap x
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();
1461  // COMPUTE GRADIENT WRT U
1462  g.zero();
1463  Real A = ((*(*up)[first_])[solindex_] - (*(*up)[last_])[solindex_] )/2.0;
1464  if(ptype_==1)
1465  {
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);
1468  }
1469 
1470  if(ptype_==2)
1471  {
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 );
1474  }
1475 
1476  if(ptype_==3)
1477  {
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 );
1480  }
1481 
1482  if(ptype_==4)
1483  {
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 );
1486  }
1487 
1488  if(ptype_==5)
1489  {
1490  (*(*gup)[first_])[solindex_] = 0.5*alpha_;
1491  (*(*gup)[last_ ])[solindex_] = -0.5*alpha_;
1492  }
1493  }
1494 
1495  void gradient_2( ::ROL::Vector<Real> &g, const ::ROL::Vector<Real> &u, const ::ROL::Vector<Real> &z, Real &tol )
1496  {
1497  // Unwrap g
1498  Teuchos::RCP<std::vector<Real> > gzp = Teuchos::rcp_const_cast<std::vector<Real> >((Teuchos::dyn_cast<const ::ROL::StdVector<Real> >(g)).getVector());
1499  // Unwrap x
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();
1502  // COMPUTE GRADIENT WRT Z
1503  for (int i=0; i<nz_; i++)
1504  {
1505  (*gzp)[i] = 0.0;
1506  }
1507  }
1508 
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 )
1511  {
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());
1513  // Unwrap v
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();
1516 
1517  hv.zero();
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;
1521 
1522  if(ptype_==1)
1523  {
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);
1526  }
1527 
1528  if(ptype_==2)
1529  {
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);
1532  }
1533 
1534  if(ptype_==3)
1535  {
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);
1538  }
1539 
1540  if(ptype_==4)
1541  {
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);
1544  }
1545 
1546  }
1547 
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 )
1550  {
1551  hv.zero();
1552  }
1553 
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 )
1556  {
1557  hv.zero();
1558  }
1559 
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 )
1562  {
1563  hv.zero();
1564  }
1565 
1566 };
1567 
1568 template<class Real>
1569 class SumObjective : public virtual ::ROL::Objective<Real>
1570 {
1571 
1572 private:
1573  std::vector< Teuchos::RCP< ::ROL::Objective<Real> > > objVec_;
1574  std::vector<bool> types_;
1575 
1576 public:
1577 
1578  SumObjective(std::vector< Teuchos::RCP< ::ROL::Objective<Real> > > & objVec,
1579  std::vector<bool> & types)
1580  : objVec_(objVec), types_(types)
1581  {}
1582 
1583  virtual void update( const ::ROL::Vector<Real> &x, bool flag = true, int iter = -1 )
1584  {
1585  for(int i=0;i<objVec_.size();i++)
1586  {
1587  if(types_[i]==0)
1588  {
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); // CVaR vector
1591  }
1592  else if(types_[i]==1)
1593  {
1594  objVec_[i]->update(x,flag,iter);
1595  }
1596  }
1597  }
1598 
1599  virtual Real value( const ::ROL::Vector<Real> &x, Real &tol )
1600  {
1601  Real val = 0;
1602  for(int i=0;i<objVec_.size();i++)
1603  {
1604  if(types_[i]==0)
1605  {
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); // CVaR vector
1608  }
1609  else if(types_[i]==1)
1610  {
1611  val += objVec_[i]->value(x,tol);
1612  }
1613  }
1614  return val;
1615  }
1616 
1617  virtual void gradient( ::ROL::Vector<Real> &g, const ::ROL::Vector<Real> &x, Real &tol )
1618  {
1619  g.zero();
1620  Teuchos::RCP< ::ROL::Vector<Real> > temp = g.clone();
1621  for(int i=0;i<objVec_.size();i++)
1622  {
1623  if(types_[i]==0)
1624  {
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); // CVaR vector
1628  g.plus(tempc);
1629  }
1630  else if(types_[i]==1)
1631  {
1632  objVec_[i]->gradient(*temp,x,tol);
1633  g.plus(*temp);
1634  }
1635  }
1636  }
1637 
1638  virtual void hessVec( ::ROL::Vector<Real> &hv, const ::ROL::Vector<Real> &v, const ::ROL::Vector<Real> &x, Real &tol )
1639  {
1640  hv.zero();
1641  Teuchos::RCP< ::ROL::Vector<Real> > temp = hv.clone();
1642  for(int i=0;i<objVec_.size();i++)
1643  {
1644  if(types_[i]==0)
1645  {
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); // CVaR vector
1650  hv.plus(tempc);
1651  }
1652  else if(types_[i]==1)
1653  {
1654  objVec_[i]->hessVec(*temp,v,x,tol);
1655  hv.plus(*temp);
1656  }
1657  }
1658  }
1659 
1660 };
1661 
1662 template<class Real>
1663 class BoundConstraint_ROL_DC : public ::ROL::BoundConstraint<Real>
1664 {
1665 
1666 private:
1667  /// Vector of lower bounds
1668  std::vector<Real> x_lo_;
1669  /// Vector of upper bounds
1670  std::vector<Real> x_up_;
1671  /// Half of the minimum distance between upper and lower bounds
1672  Real min_diff_;
1673  /// Scaling for the epsilon margin
1674  Real scale_;
1675  int n;
1676 
1677 public:
1678  BoundConstraint_ROL_DC( Real scale, std::vector<Real> lo_bounds, std::vector<Real> up_bounds )
1679  {
1680  //assume lo_bounds and up_bounds are same size
1681  n = lo_bounds.size();
1682  for (int i=0;i<n;i++)
1683  {
1684  x_lo_.push_back(lo_bounds[i]);
1685  }
1686  for (int i=0;i<n;i++)
1687  {
1688  x_up_.push_back(up_bounds[i]);
1689  }
1690 
1691  scale_ = scale;
1692  min_diff_ = 1.e+200;
1693  for(int i=0;i<n;i++)
1694  {
1695  min_diff_ = std::min(x_up_[i]-x_lo_[i],min_diff_);
1696  }
1697  min_diff_ *= 0.5;
1698  }
1699 
1700  void project( ::ROL::Vector<Real> &x )
1701  {
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++)
1705  {
1706  (*ex)[i] = std::max(x_lo_[i],std::min(x_up_[i],(*ex)[i]));
1707  }
1708  }
1709  bool isFeasible( const ::ROL::Vector<Real> &x )
1710  {
1711  Teuchos::RCP<const std::vector<Real> > ex =
1712  (Teuchos::dyn_cast<::ROL::StdVector<Real> >(const_cast<::ROL::Vector<Real> &>(x))).getVector();
1713  bool bs = true;
1714  for (int i=0;i<n;i++)
1715  {
1716  bs = bs && (((*ex)[i] >= this->x_lo_[i])?true:false);
1717  bs = bs && (((*ex)[i] <= this->x_up_[i])?true:false);
1718  }
1719  return bs;
1720  }
1721 
1722  void pruneLowerActive(::ROL::Vector<Real> &v, const ::ROL::Vector<Real> &x, Real eps)
1723  {
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++ )
1730  {
1731  if ( ((*ex)[i] <= this->x_lo_[i]+epsn) )
1732  {
1733  (*ev)[i] = 0.0;
1734  }
1735  }
1736  }
1737 
1738  void pruneUpperActive(::ROL::Vector<Real> &v, const ::ROL::Vector<Real> &x, Real eps)
1739  {
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++ )
1746  {
1747  if ( ((*ex)[i] >= this->x_up_[i]-epsn) )
1748  {
1749  (*ev)[i] = 0.0;
1750  }
1751  }
1752  }
1753 
1754  void pruneActive(::ROL::Vector<Real> &v, const ::ROL::Vector<Real> &x, Real eps)
1755  {
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++ )
1762  {
1763  if ( ((*ex)[i] <= this->x_lo_[i]+epsn) ||
1764  ((*ex)[i] >= this->x_up_[i]-epsn) )
1765  {
1766  (*ev)[i] = 0.0;
1767  }
1768  }
1769  }
1770 
1771  void pruneLowerActive(::ROL::Vector<Real> &v, const ::ROL::Vector<Real> &g, const ::ROL::Vector<Real> &x, Real eps)
1772  {
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++ )
1781  {
1782  if ( ((*ex)[i] <= this->x_lo_[i]+epsn && (*eg)[i] > 0.0) )
1783  {
1784  (*ev)[i] = 0.0;
1785  }
1786  }
1787  }
1788 
1789  void pruneUpperActive(::ROL::Vector<Real> &v, const ::ROL::Vector<Real> &g, const ::ROL::Vector<Real> &x, Real eps)
1790  {
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++ )
1799  {
1800  if ( ((*ex)[i] >= this->x_up_[i]-epsn && (*eg)[i] < 0.0) )
1801  {
1802  (*ev)[i] = 0.0;
1803  }
1804  }
1805  }
1806 
1807  void pruneActive(::ROL::Vector<Real> &v, const ::ROL::Vector<Real> &g, const ::ROL::Vector<Real> &x, Real eps)
1808  {
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++ )
1817  {
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) )
1820  {
1821  (*ev)[i] = 0.0;
1822  }
1823  }
1824  }
1825 
1826  void setVectorToUpperBound( ::ROL::Vector<Real> &u )
1827  {
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) );
1831  u.set(*up);
1832  }
1833 
1834  void setVectorToLowerBound( ::ROL::Vector<Real> &l )
1835  {
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) );
1839  l.set(*lp);
1840  }
1841 
1842 };
1843 
1844 } // namespace Analysis
1845 } // namespace Xyce`
1846 
1847 #endif // Xyce_ROL
1848 
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.
Definition: N_DEV_Pars.h:1224
AnalysisManager & analysisManager_
Definition: N_ANP_AC.C:960
Linear::System & linearSystem_
Definition: N_ANP_AC.C:961