ROL
ROL_DiodeCircuit.hpp
Go to the documentation of this file.
1 #ifndef ROL_DIODECIRCUIT_HPP
2 #define ROL_DIODECIRCUIT_HPP
3 
4 #include "ROL_StdVector.hpp"
5 #include "ROL_Objective.hpp"
7 
8 #include <iostream>
9 #include <fstream>
10 #include <string>
11 
18 namespace ROL {
19 namespace ZOO {
20 
32  template<class Real>
33  class Objective_DiodeCircuit : public Objective<Real> {
34  private:
36  Real Vth_;
38  Teuchos::RCP<std::vector<Real> > Imeas_;
40  Teuchos::RCP<std::vector<Real> > Vsrc_;
42  bool lambertw_;
44  Real noise_;
49 
50  public:
51 
59  Objective_DiodeCircuit(Real Vth, Real Vsrc_min, Real Vsrc_max, Real Vsrc_step, Real true_Is, Real true_Rs, bool lambertw, Real noise, bool use_adjoint, int use_hessvec){
61  Vth_ = Vth;
62  use_adjoint_ = use_adjoint;
63  use_hessvec_ = use_hessvec;
64  int n = (Vsrc_max-Vsrc_min)/Vsrc_step + 1;
65  Vsrc_ = Teuchos::rcp(new std::vector<Real>(n,0.0));
66  Imeas_ = Teuchos::rcp(new std::vector<Real>(n,0.0));
67  std::ofstream output ("Measurements.dat");
68  Real left = 0.0; Real right = 1.0;
69  if(lambertw_){
70  // Using Lambert-W function
71  std::cout << "Generating data using Lambert-W function." << std::endl;
72  for(int i=0;i<n;i++){
73  (*Vsrc_)[i] = Vsrc_min+i*Vsrc_step;
74  (*Imeas_)[i] = lambertWCurrent(true_Is,true_Rs,(*Vsrc_)[i]);
75  if(noise>0.0){
76  (*Imeas_)[i] += noise*pow(10,int(log10((*Imeas_)[i])))*(( (Real)rand() / (Real)RAND_MAX ) * (right - left) + left);
77  }
78  // Write generated data into file
79  if(output.is_open()){
80  output << std::setprecision(8) << std::scientific << (*Vsrc_)[i] << " " << (*Imeas_)[i] << "\n";
81  }
82  }
83  }
84  else{
85  // Using Newton's method
86  std::cout << "Generating data using Newton's method." << std::endl;
87  for(int i=0;i<n;i++){
88  (*Vsrc_)[i] = Vsrc_min+i*Vsrc_step;
89  Real I0 = 1.e-12; // initial guess for Newton
90  (*Imeas_)[i] = Newton(I0,Vsrc_min+i*Vsrc_step,true_Is,true_Rs);
91  if(noise>0.0){
92  (*Imeas_)[i] += noise*pow(10,int(log10((*Imeas_)[i])))*(( (Real)rand() / (Real)RAND_MAX ) * (right - left) + left);
93  }
94  // Write generated data into file
95  if(output.is_open()){
96  output << std::setprecision(8) << std::scientific << (*Vsrc_)[i] << " " << (*Imeas_)[i] << "\n";
97  }
98  }
99  }
100 
101  output.close();
102  }
103 
111  Objective_DiodeCircuit(Real Vth, std::ifstream& input_file, bool lambertw, Real noise, bool use_adjoint, int use_hessvec){
112  lambertw_ = lambertw;
113  Vth_ = Vth;
114  use_adjoint_ = use_adjoint;
115  use_hessvec_ = use_hessvec;
116  std::string line;
117  int dim = 0;
118  for(int k=0;std::getline(input_file,line);++k){dim=k;} // count number of lines
119  input_file.clear(); // reset to beginning of file
120  input_file.seekg(0,std::ios::beg);
121  Vsrc_ = Teuchos::rcp(new std::vector<Real>(dim,0.0));
122  Imeas_ = Teuchos::rcp(new std::vector<Real>(dim,0.0));
123  double Vsrc, Imeas;
124  std::cout << "Using input file to generate data." << "\n";
125  for(int i=0;i<dim;i++){
126  input_file >> Vsrc;
127  input_file >> Imeas;
128  (*Vsrc_)[i] = Vsrc;
129  (*Imeas_)[i] = Imeas;
130  }
131  input_file.close();
132  }
133 
135  void set_method(bool lambertw){
137  }
138 
149  Real diode(const Real I, const Real Vsrc, const Real Is, const Real Rs){
150  return I-Is*(exp((Vsrc-I*Rs)/Vth_)-1);
151  }
152 
154  Real diodeI(const Real I, const Real Vsrc, const Real Is, const Real Rs){
155  return 1+Is*exp((Vsrc-I*Rs)/Vth_)*(Rs/Vth_);
156  }
157 
159  Real diodeIs(const Real I, const Real Vsrc, const Real Is, const Real Rs){
160  return 1-exp((Vsrc-I*Rs)/Vth_);
161  }
162 
164  Real diodeRs(const Real I, const Real Vsrc, const Real Is, const Real Rs){
165  return Is*exp((Vsrc-I*Rs)/Vth_)*(I/Vth_);
166  }
167 
169  Real diodeII(const Real I, const Real Vsrc, const Real Is, const Real Rs){
170  return -Is*exp((Vsrc-I*Rs)/Vth_)*(Rs/Vth_)*(Rs/Vth_);
171  }
172 
174  Real diodeIIs(const Real I, const Real Vsrc, const Real Is, const Real Rs){
175  return exp((Vsrc-I*Rs)/Vth_)*(Rs/Vth_);
176  }
177 
179  Real diodeIRs(const Real I, const Real Vsrc, const Real Is, const Real Rs){
180  return (Is/Vth_)*exp((Vsrc-I*Rs)/Vth_)*(1-(I*Rs)/Vth_);
181  }
182 
184  Real diodeIsIs(const Real I, const Real Vsrc, const Real Is, const Real Rs){
185  return 0;
186  }
187 
189  Real diodeIsRs(const Real I, const Real Vsrc, const Real Is, const Real Rs){
190  return exp((Vsrc-I*Rs)/Vth_)*(I/Vth_);
191  }
192 
194  Real diodeRsRs(const Real I, const Real Vsrc, const Real Is, const Real Rs){
195  return -Is*exp((Vsrc-I*Rs)/Vth_)*(I/Vth_)*(I/Vth_);
196  }
197 
205  Real Newton(const Real I, const Real Vsrc, const Real Is, const Real Rs){
206  double EPS = 1.e-16;
207  double TOL = 1.e-13;
208  int MAXIT = 200;
209  Real IN = I;
210  Real fval = diode(IN,Vsrc,Is,Rs);
211  Real dfval = 0.0;
212  Real IN_tmp = 0.0;
213  Real fval_tmp = 0.0;
214  Real alpha = 1.0;
215  for ( int i = 0; i < MAXIT; i++ ) {
216  if ( std::fabs(fval) < TOL ) {
217  // std::cout << "converged with |fval| = " << std::fabs(fval) << " and TOL = " << TOL << "\n";
218  break;
219  }
220  dfval = diodeI(IN,Vsrc,Is,Rs);
221  if( std::fabs(dfval) < EPS ){
222  std::cout << "denominator is too small" << std::endl;
223  break;
224  }
225 
226  alpha = 1.0;
227  IN_tmp = IN - alpha*fval/dfval;
228  fval_tmp = diode(IN_tmp,Vsrc,Is,Rs);
229  while ( std::fabs(fval_tmp) >= (1.0-1.e-4*alpha)*std::fabs(fval) ) {
230  alpha /= 2.0;
231  IN_tmp = IN - alpha*fval/dfval;
232  fval_tmp = diode(IN_tmp,Vsrc,Is,Rs);
233  if ( alpha < std::sqrt(EPS) ) {
234  // std::cout << "Step tolerance met\n";
235  break;
236  }
237  }
238  IN = IN_tmp;
239  fval = fval_tmp;
240  // if ( i == MAXIT-1){
241  // std::cout << "did not converge " << std::fabs(fval) << "\n";
242  // }
243  }
244  return IN;
245  }
246 
278  void lambertw(Real x, Real &w, int &ierr, Real &xi){
279  int i=0, maxit = 10;
280  const Real turnpt = -exp(-1.), c1 = 1.5, c2 = .75;
281  Real r, r2, r3, s, mach_eps, relerr = 1., diff;
282  mach_eps = 2.e-15; // float:2e-7
283  ierr = 0;
284 
285  if( x > c1){
286  w = c2*log(x);
287  xi = log( x/ w) - w;
288  }
289  else{
290  if( x >= 0.0){
291  w = x;
292  if( x == 0. ) return;
293  if( x < (1-c2) ) w = x*(1.-x + c1*x*x);
294  xi = - w;
295  }
296  else{
297  if( x >= turnpt){
298  if( x > -0.2 ){
299  w = x*(1.0-x + c1*x*x);
300  xi = log(1.0-x + c1*x*x) - w;
301  }
302  else{
303  diff = x-turnpt;
304  if( diff < 0.0 ) diff = -diff;
305  w = -1 + sqrt(2.0*exp(1.))*sqrt(x-turnpt);
306  if( diff == 0.0 ) return;
307  xi = log( x/ w) - w;
308  }
309  }
310  else{
311  ierr = 1; // x is not in the domain.
312  w = -1.0;
313  return;
314  }
315  }
316  }
317 
318  while( relerr > mach_eps && i<maxit){
319  r = xi/(w+1.0); //singularity at w=-1
320  r2 = r*r;
321  r3 = r2*r;
322  s = 6.*(w+1.0)*(w+1.0);
323  w = w * ( 1.0 + r + r2/(2.0*( w+1.0)) - (2. * w -1.0)*r3/s );
324  if( w * x < 0.0 ) w = -w;
325  xi = log( x/ w) - w;
326 
327  if( x>1.0 ){
328  relerr = xi / w;
329  }
330  else{
331  relerr = xi;
332  }
333  if(relerr < 0.0 ) relerr = -relerr;
334  ++i;
335  }
336  if( i == maxit ) ierr = 2;
337  }
338 
351  Real lambertWCurrent(Real Is, Real Rs, Real Vsrc){
352  Real arg1 = (Vsrc + Is*Rs)/Vth_;
353  Real evd = exp(arg1);
354  Real lambWArg = Is*Rs*evd/Vth_;
355  Real lambWReturn;
356  int ierr;
357  Real lambWError;
358  lambertw(lambWArg, lambWReturn, ierr, lambWError);
359  if(ierr == 1){std::cout << "LambertW error: argument is not in the domain" << std::endl; return -1.0;}
360  if(ierr == 2){std::cout << "LambertW error: BUG!" << std::endl;}
361  Real Id = -Is+Vth_*(lambWReturn)/Rs;
362  //Real Gd = lambWReturn / ((1 + lambWReturn)*RS);
363  return Id;
364  }
365 
366 
369  Teuchos::RCP<std::vector<Real> > Ip =
370  Teuchos::rcp_const_cast<std::vector<Real> >((Teuchos::dyn_cast<StdVector<Real> >(I)).getVector());
371  Teuchos::RCP<const std::vector<Real> > Sp =
372  (Teuchos::dyn_cast<StdVector<Real> >(const_cast<Vector<Real> &>(S))).getVector();
373 
374  int n = Ip->size();
375 
376  if(lambertw_){
377  // Using Lambert-W function
378  Real lambval;
379  for(int i=0;i<n;i++){
380  lambval = lambertWCurrent((*Sp)[0],(*Sp)[1],(*Vsrc_)[i]);
381  (*Ip)[i] = lambval;
382  }
383  }
384  else{
385  // Using Newton's method
386  Real I0 = 1.e-12; // Initial guess for Newton
387  for(int i=0;i<n;i++){
388  (*Ip)[i] = Newton(I0,(*Vsrc_)[i],(*Sp)[0],(*Sp)[1]);
389  }
390  }
391  }
392 
393 
401  Real value(const Vector<Real> &S, Real &tol){
402  Teuchos::RCP<const std::vector<Real> > Sp =
403  (Teuchos::dyn_cast<StdVector<Real> >(const_cast<Vector<Real> &>(S))).getVector();
404  int n = Imeas_->size();
405  StdVector<Real> I( Teuchos::rcp( new std::vector<Real>(n,0.0) ) );
406  Teuchos::RCP<std::vector<Real> > Ip =
407  Teuchos::rcp_const_cast<std::vector<Real> >((Teuchos::dyn_cast<StdVector<Real> >(I)).getVector());
408 
409  // Solve state equation
410  this->solve_circuit(I,S);
411  Real val = 0;
412 
413  for(int i=0;i<n;i++){
414  val += ((*Ip)[i]-(*Imeas_)[i])*((*Ip)[i]-(*Imeas_)[i]);
415  }
416  return val/2.0;
417  }
418 
426  void solve_adjoint(Vector<Real> &lambda, const Vector<Real> &I, const Vector<Real> &S){
427 
428  Teuchos::RCP<std::vector<Real> > lambdap =
429  Teuchos::rcp_const_cast<std::vector<Real> >((Teuchos::dyn_cast<StdVector<Real> >(lambda)).getVector());
430 
431  Teuchos::RCP<const std::vector<Real> > Ip =
432  (Teuchos::dyn_cast<StdVector<Real> >(const_cast<Vector<Real> &>(I))).getVector();
433  Teuchos::RCP<const std::vector<Real> > Sp =
434  (Teuchos::dyn_cast<StdVector<Real> >(const_cast<Vector<Real> &>(S))).getVector();
435 
436  int n = Ip->size();
437  for(int i=0;i<n;i++){
438  (*lambdap)[i] = ((*Imeas_)[i]-(*Ip)[i])/diodeI((*Ip)[i],(*Vsrc_)[i],(*Sp)[0],(*Sp)[1]);
439  }
440  }
441 
450 
451  Teuchos::RCP<std::vector<Real> > sensp =
452  Teuchos::rcp_const_cast<std::vector<Real> >((Teuchos::dyn_cast<StdVector<Real> >(sens)).getVector());
453  Teuchos::RCP<const std::vector<Real> > Ip =
454  (Teuchos::dyn_cast<StdVector<Real> >(const_cast<Vector<Real> &>(I))).getVector();
455  Teuchos::RCP<const std::vector<Real> > Sp =
456  (Teuchos::dyn_cast<StdVector<Real> >(const_cast<Vector<Real> &>(S))).getVector();
457 
458  int n = Ip->size();
459  for(int i=0;i<n;i++){
460  (*sensp)[i] = -diodeIs((*Ip)[i],(*Vsrc_)[i],(*Sp)[0],(*Sp)[1])/diodeI((*Ip)[i],(*Vsrc_)[i],(*Sp)[0],(*Sp)[1]);
461  }
462  }
463 
472 
473  Teuchos::RCP<std::vector<Real> > sensp =
474  Teuchos::rcp_const_cast<std::vector<Real> >((Teuchos::dyn_cast<StdVector<Real> >(sens)).getVector());
475  Teuchos::RCP<const std::vector<Real> > Ip =
476  (Teuchos::dyn_cast<StdVector<Real> >(const_cast<Vector<Real> &>(I))).getVector();
477  Teuchos::RCP<const std::vector<Real> > Sp =
478  (Teuchos::dyn_cast<StdVector<Real> >(const_cast<Vector<Real> &>(S))).getVector();
479 
480  int n = Ip->size();
481  for(int i=0;i<n;i++){
482  (*sensp)[i] = -diodeRs((*Ip)[i],(*Vsrc_)[i],(*Sp)[0],(*Sp)[1])/diodeI((*Ip)[i],(*Vsrc_)[i],(*Sp)[0],(*Sp)[1]);
483  }
484  }
485 
487  void gradient(Vector<Real> &g, const Vector<Real> &S, Real &tol){
488 
489  Teuchos::RCP<std::vector<Real> > gp =
490  Teuchos::rcp_const_cast<std::vector<Real> >((Teuchos::dyn_cast<StdVector<Real> >(g)).getVector());
491  Teuchos::RCP<const std::vector<Real> > Sp =
492  (Teuchos::dyn_cast<StdVector<Real> >(const_cast<Vector<Real> &>(S))).getVector();
493 
494  int n = Imeas_->size();
495 
496  StdVector<Real> I( Teuchos::rcp( new std::vector<Real>(n,0.0) ) );
497  Teuchos::RCP<std::vector<Real> > Ip =
498  Teuchos::rcp_const_cast<std::vector<Real> >((Teuchos::dyn_cast<StdVector<Real> >(I)).getVector());
499 
500  // Solve state equation
501  this->solve_circuit(I,S);
502 
503  if(use_adjoint_){
504  // Compute the gradient of the reduced objective function using adjoint computation
505  StdVector<Real> lambda( Teuchos::rcp( new std::vector<Real>(n,0.0) ) );
506  Teuchos::RCP<std::vector<Real> > lambdap =
507  Teuchos::rcp_const_cast<std::vector<Real> >((Teuchos::dyn_cast<StdVector<Real> >(lambda)).getVector());
508 
509  // Solve adjoint equation
510  this->solve_adjoint(lambda,I,S);
511 
512  // Compute gradient
513  (*gp)[0] = 0.0; (*gp)[1] = 0.0;
514  for(int i=0;i<n;i++){
515  (*gp)[0] += diodeIs((*Ip)[i],(*Vsrc_)[i],(*Sp)[0],(*Sp)[1])*(*lambdap)[i];
516  (*gp)[1] += diodeRs((*Ip)[i],(*Vsrc_)[i],(*Sp)[0],(*Sp)[1])*(*lambdap)[i];
517  }
518  }
519  else{
520  // Compute the gradient of the reduced objective function using sensitivities
521  StdVector<Real> sensIs( Teuchos::rcp( new std::vector<Real>(n,0.0) ) );
522  StdVector<Real> sensRs( Teuchos::rcp( new std::vector<Real>(n,0.0) ) );
523  // Solve sensitivity equations
524  this->solve_sensitivity_Is(sensIs,I,S);
525  this->solve_sensitivity_Rs(sensRs,I,S);
526 
527  Teuchos::RCP<std::vector<Real> > sensIsp = Teuchos::rcp_const_cast<std::vector<Real> >((Teuchos::dyn_cast<StdVector<Real> >(sensIs)).getVector());
528  Teuchos::RCP<std::vector<Real> > sensRsp = Teuchos::rcp_const_cast<std::vector<Real> >((Teuchos::dyn_cast<StdVector<Real> >(sensRs)).getVector());
529 
530  // Write sensitivities into file
531  std::ofstream output ("Sensitivities.dat");
532  for(int k=0;k<n;k++){
533  if(output.is_open()){
534  output << std::scientific << (*sensIsp)[k] << " " << (*sensRsp)[k] << "\n";
535  }
536  }
537  output.close();
538  // Compute gradient
539  (*gp)[0] = 0.0; (*gp)[1] = 0.0;
540  for(int i=0;i<n;i++){
541  (*gp)[0] += ((*Ip)[i]-(*Imeas_)[i])*(*sensIsp)[i];
542  (*gp)[1] += ((*Ip)[i]-(*Imeas_)[i])*(*sensRsp)[i];
543  }
544  }
545  }
546 
547 
548 
556  void hessVec( Vector<Real> &hv, const Vector<Real> &v, const Vector<Real> &S, Real &tol ){
557  if(use_hessvec_==0){
558  // Use finite-difference approximation
559  // Modification of parent class function that takes into accout different scale of components
560  Teuchos::RCP<const std::vector<Real> > vp =
561  (Teuchos::dyn_cast<StdVector<Real> >(const_cast<Vector<Real> &>(v))).getVector();
562  Teuchos::RCP<const std::vector<Real> > Sp =
563  (Teuchos::dyn_cast<StdVector<Real> >(const_cast<Vector<Real> &>(S))).getVector();
564  Real gtol = std::sqrt(ROL_EPSILON);
565 
566  // Get Step Length
567  Real h = std::max(1.0,S.norm()/v.norm())*tol;
568  //Real h = 2.0/(v.norm()*v.norm())*tol;
569 
570  // Find the scale of componenets of S
571  Real Is_scale = pow( 10,int( log10( (*Sp)[0] ) ) );
572  Real Rs_scale = pow( 10,int( log10( (*Sp)[1] ) ) );
573  // Apply scaling
574  Real h1 = Is_scale*h;
575  Real h2 = Rs_scale*h;
576 
577  // Compute Gradient at S
578  Teuchos::RCP<Vector<Real> > g = S.clone();
579  this->gradient(*g,S,gtol);
580 
581  // Compute New Step S + h*v
582  Teuchos::RCP<std::vector<Real> > Snewp = Teuchos::rcp( new std::vector<Real> (2, 0.0) );
583  ROL::StdVector<Real> Snew(Snewp);
584  (*Snewp)[0] = (*Sp)[0] + h1*(*vp)[0];
585  (*Snewp)[1] = (*Sp)[1] + h2*(*vp)[1];
586 
587  // Compute Gradient at x + h*v
588  hv.zero();
589  this->gradient(hv,Snew,gtol);
590 
591  // Compute Newton Quotient
592  hv.axpy(-1.0,*g);
593  hv.scale(1.0/std::sqrt(h1*h1+h2*h2));
594  }
595  else if(use_hessvec_==1){
596  Teuchos::RCP<std::vector<Real> > hvp =
597  Teuchos::rcp_const_cast<std::vector<Real> >((Teuchos::dyn_cast<StdVector<Real> >(hv)).getVector());
598  Teuchos::RCP<const std::vector<Real> > vp =
599  (Teuchos::dyn_cast<StdVector<Real> >(const_cast<Vector<Real> &>(v))).getVector();
600  Teuchos::RCP<const std::vector<Real> > Sp =
601  (Teuchos::dyn_cast<StdVector<Real> >(const_cast<Vector<Real> &>(S))).getVector();
602 
603  int n = Imeas_->size();
604 
605  StdVector<Real> I( Teuchos::rcp( new std::vector<Real>(n,0.0) ) );
606  Teuchos::RCP<std::vector<Real> > Ip =
607  Teuchos::rcp_const_cast<std::vector<Real> >((Teuchos::dyn_cast<StdVector<Real> >(I)).getVector());
608 
609  // Solve state equation
610  this->solve_circuit(I,S);
611 
612 
613  StdVector<Real> lambda( Teuchos::rcp( new std::vector<Real>(n,0.0) ) );
614  Teuchos::RCP<std::vector<Real> > lambdap =
615  Teuchos::rcp_const_cast<std::vector<Real> >((Teuchos::dyn_cast<StdVector<Real> >(lambda)).getVector());
616 
617  // Solve adjoint equation
618  this->solve_adjoint(lambda,I,S);
619 
620  StdVector<Real> w( Teuchos::rcp( new std::vector<Real>(n,0.0) ) );
621  Teuchos::RCP<std::vector<Real> > wp =
622  Teuchos::rcp_const_cast<std::vector<Real> >((Teuchos::dyn_cast<StdVector<Real> >(w)).getVector());
623 
624  // Solve state sensitivity equation
625  for(int i=0;i<n;i++){
626  (*wp)[i] = ( (*vp)[0] * diodeIs( (*Ip)[i],(*Vsrc_)[i],(*Sp)[0],(*Sp)[1] ) + (*vp)[1] * diodeRs( (*Ip)[i],(*Vsrc_)[i],(*Sp)[0],(*Sp)[1] ) ) / diodeI((*Ip)[i],(*Vsrc_)[i],(*Sp)[0],(*Sp)[1]);
627  }
628 
629  StdVector<Real> p( Teuchos::rcp( new std::vector<Real>(n,0.0) ) );
630  Teuchos::RCP<std::vector<Real> > pp =
631  Teuchos::rcp_const_cast<std::vector<Real> >((Teuchos::dyn_cast<StdVector<Real> >(p)).getVector());
632 
633  // Solve for p
634  for(int j=0;j<n;j++){
635  (*pp)[j] = ( (*wp)[j] + (*lambdap)[j] * diodeII( (*Ip)[j],(*Vsrc_)[j],(*Sp)[0],(*Sp)[1] ) * (*wp)[j] - (*lambdap)[j] * diodeIIs( (*Ip)[j],(*Vsrc_)[j],(*Sp)[0],(*Sp)[1] ) * (*vp)[0] - (*lambdap)[j] * diodeIRs( (*Ip)[j],(*Vsrc_)[j],(*Sp)[0],(*Sp)[1] ) * (*vp)[1] ) / diodeI( (*Ip)[j],(*Vsrc_)[j],(*Sp)[0],(*Sp)[1] );
636  }
637 
638  // Assemble Hessian-vector product
639  (*hvp)[0] = 0.0;(*hvp)[1] = 0.0;
640  for(int k=0;k<n;k++){
641  (*hvp)[0] += diodeIs( (*Ip)[k],(*Vsrc_)[k],(*Sp)[0],(*Sp)[1] ) * (*pp)[k] - (*lambdap)[k] * (*wp)[k] * diodeIIs( (*Ip)[k],(*Vsrc_)[k],(*Sp)[0],(*Sp)[1] ) + (*lambdap)[k] * (*vp)[0] * diodeIsIs( (*Ip)[k],(*Vsrc_)[k],(*Sp)[0],(*Sp)[1] ) + (*lambdap)[k] * (*vp)[1] * diodeIsRs( (*Ip)[k],(*Vsrc_)[k],(*Sp)[0],(*Sp)[1] );
642  (*hvp)[1] += diodeRs( (*Ip)[k],(*Vsrc_)[k],(*Sp)[0],(*Sp)[1] ) * (*pp)[k] - (*lambdap)[k] * (*wp)[k] * diodeIRs( (*Ip)[k],(*Vsrc_)[k],(*Sp)[0],(*Sp)[1] ) + (*lambdap)[k] * (*vp)[0] * diodeIsRs( (*Ip)[k],(*Vsrc_)[k],(*Sp)[0],(*Sp)[1] ) + (*lambdap)[k] * (*vp)[1] * diodeRsRs( (*Ip)[k],(*Vsrc_)[k],(*Sp)[0],(*Sp)[1] );
643  }
644  }
645  else if(use_hessvec_==2){
646  //Gauss-Newton approximation
647  Teuchos::RCP<std::vector<Real> > hvp =
648  Teuchos::rcp_const_cast<std::vector<Real> >((Teuchos::dyn_cast<StdVector<Real> >(hv)).getVector());
649  Teuchos::RCP<const std::vector<Real> > vp =
650  (Teuchos::dyn_cast<StdVector<Real> >(const_cast<Vector<Real> &>(v))).getVector();
651  Teuchos::RCP<const std::vector<Real> > Sp =
652  (Teuchos::dyn_cast<StdVector<Real> >(const_cast<Vector<Real> &>(S))).getVector();
653 
654  int n = Imeas_->size();
655 
656  StdVector<Real> I( Teuchos::rcp( new std::vector<Real>(n,0.0) ) );
657  Teuchos::RCP<std::vector<Real> > Ip =
658  Teuchos::rcp_const_cast<std::vector<Real> >((Teuchos::dyn_cast<StdVector<Real> >(I)).getVector());
659 
660  // Solve state equation
661  this->solve_circuit(I,S);
662 
663  // Compute sensitivities
664  StdVector<Real> sensIs( Teuchos::rcp( new std::vector<Real>(n,0.0) ) );
665  StdVector<Real> sensRs( Teuchos::rcp( new std::vector<Real>(n,0.0) ) );
666  // Solve sensitivity equations
667  this->solve_sensitivity_Is(sensIs,I,S);
668  this->solve_sensitivity_Rs(sensRs,I,S);
669  Teuchos::RCP<std::vector<Real> > sensIsp = Teuchos::rcp_const_cast<std::vector<Real> >((Teuchos::dyn_cast<StdVector<Real> >(sensIs)).getVector());
670  Teuchos::RCP<std::vector<Real> > sensRsp = Teuchos::rcp_const_cast<std::vector<Real> >((Teuchos::dyn_cast<StdVector<Real> >(sensRs)).getVector());
671 
672  // Compute approximate Hessian
673  Real H11 = 0.0; Real H12 = 0.0; Real H22 = 0.0;
674  for(int k=0;k<n;k++){
675  H11 += (*sensIsp)[k]*(*sensIsp)[k];
676  H12 += (*sensIsp)[k]*(*sensRsp)[k];
677  H22 += (*sensRsp)[k]*(*sensRsp)[k];
678  }
679 
680  // Compute approximate Hessian-times-vector
681  (*hvp)[0] = H11*(*vp)[0] + H12*(*vp)[1];
682  (*hvp)[1] = H12*(*vp)[0] + H22*(*vp)[1];
683  }
684  else{
685  this->ROL::Objective<Real>::hessVec( hv, v, S, tol ); // Use parent class function
686  }
687  }
688 
700  void generate_plot(Real Is_lo, Real Is_up, Real Is_step, Real Rs_lo, Real Rs_up, Real Rs_step){
701  Teuchos::RCP<std::vector<double> > S_rcp = Teuchos::rcp(new std::vector<double>(2,0.0) );
702  ROL::StdVector<double> S(S_rcp);
703  std::ofstream output ("Objective.dat");
704 
705  Real Is = 0.0;
706  Real Rs = 0.0;
707  Real val = 0.0;
708  Real tol = 1.e-16;
709  int n = (Is_up-Is_lo)/Is_step + 1;
710  int m = (Rs_up-Rs_lo)/Rs_step + 1;
711  for(int i=0;i<n;i++){
712  Is = Is_lo + i*Is_step;
713  for(int j=0;j<m;j++){
714  Rs = Rs_lo + j*Rs_step;
715  (*S_rcp)[0] = Is;
716  (*S_rcp)[1] = Rs;
717  val = this->value(S,tol);
718  if(output.is_open()){
719  output << std::scientific << Is << " " << Rs << " " << val << "\n";
720  }
721  }
722  }
723  output.close();
724  }
725 
726 
727  };
728 
735  template<class Real>
737  private:
739  std::vector<Real> x_lo_;
741  std::vector<Real> x_up_;
743  Real min_diff_;
745  Real scale_;
746  public:
747  BoundConstraint_DiodeCircuit( Real scale, Real lo_Is, Real up_Is, Real lo_Rs, Real up_Rs ){
748  x_lo_.push_back(lo_Is);
749  x_lo_.push_back(lo_Rs);
750 
751  x_up_.push_back(up_Is);
752  x_up_.push_back(up_Rs);
753 
754  scale_ = scale;
755  min_diff_ = 0.5*std::min(x_up_[0]-x_lo_[0],x_up_[1]-x_lo_[1]);
756  }
757  void project( Vector<Real> &x ) {
758  Teuchos::RCP<std::vector<Real> > ex =
759  Teuchos::rcp_const_cast<std::vector<Real> >((Teuchos::dyn_cast<StdVector<Real> >(x)).getVector());
760  (*ex)[0] = std::max(x_lo_[0],std::min(x_up_[0],(*ex)[0]));
761  (*ex)[1] = std::max(x_lo_[1],std::min(x_up_[1],(*ex)[1]));
762  }
763  bool isFeasible( const Vector<Real> &x ) {
764  Teuchos::RCP<const std::vector<Real> > ex =
765  (Teuchos::dyn_cast<StdVector<Real> >(const_cast<Vector<Real> &>(x))).getVector();
766  return ((*ex)[0] >= this->x_lo_[0] && (*ex)[1] >= this->x_lo_[1] &&
767  (*ex)[0] <= this->x_up_[0] && (*ex)[1] <= this->x_up_[1]);
768  }
769  void pruneLowerActive(Vector<Real> &v, const Vector<Real> &x, Real eps) {
770  Teuchos::RCP<const std::vector<Real> > ex =
771  (Teuchos::dyn_cast<StdVector<Real> >(const_cast<Vector<Real> &>(x))).getVector();
772  Teuchos::RCP<std::vector<Real> > ev =
773  Teuchos::rcp_const_cast<std::vector<Real> >((Teuchos::dyn_cast<StdVector<Real> >(v)).getVector());
774  Real epsn = std::min(this->scale_*eps,this->min_diff_);
775  //epsn *= this->scale_;
776  for ( int i = 0; i < 2; i++ ) {
777  if ( ((*ex)[i] <= this->x_lo_[i]+epsn) ) {
778  (*ev)[i] = 0.0;
779  }
780  }
781  }
782  void pruneUpperActive(Vector<Real> &v, const Vector<Real> &x, Real eps) {
783  Teuchos::RCP<const std::vector<Real> > ex =
784  (Teuchos::dyn_cast<StdVector<Real> >(const_cast<Vector<Real> &>(x))).getVector();
785  Teuchos::RCP<std::vector<Real> > ev =
786  Teuchos::rcp_const_cast<std::vector<Real> >((Teuchos::dyn_cast<StdVector<Real> >(v)).getVector());
787  Real epsn = std::min(this->scale_*eps,this->min_diff_);
788  //epsn *= this->scale_;
789  for ( int i = 0; i < 2; i++ ) {
790  if ( ((*ex)[i] >= this->x_up_[i]-epsn) ) {
791  (*ev)[i] = 0.0;
792  }
793  }
794  }
795  void pruneActive(Vector<Real> &v, const Vector<Real> &x, Real eps) {
796  Teuchos::RCP<const std::vector<Real> > ex =
797  (Teuchos::dyn_cast<StdVector<Real> >(const_cast<Vector<Real> &>(x))).getVector();
798  Teuchos::RCP<std::vector<Real> > ev =
799  Teuchos::rcp_const_cast<std::vector<Real> >((Teuchos::dyn_cast<StdVector<Real> >(v)).getVector());
800  Real epsn = std::min(this->scale_*eps,this->min_diff_);
801  //epsn *= this->scale_;
802  for ( int i = 0; i < 2; i++ ) {
803  if ( ((*ex)[i] <= this->x_lo_[i]+epsn) ||
804  ((*ex)[i] >= this->x_up_[i]-epsn) ) {
805  (*ev)[i] = 0.0;
806  }
807  }
808  }
809  void pruneLowerActive(Vector<Real> &v, const Vector<Real> &g, const Vector<Real> &x, Real eps) {
810  Teuchos::RCP<const std::vector<Real> > ex =
811  (Teuchos::dyn_cast<StdVector<Real> >(const_cast<Vector<Real> &>(x))).getVector();
812  Teuchos::RCP<const std::vector<Real> > eg =
813  (Teuchos::dyn_cast<StdVector<Real> >(const_cast<Vector<Real> &>(g))).getVector();
814  Teuchos::RCP<std::vector<Real> > ev =
815  Teuchos::rcp_const_cast<std::vector<Real> >((Teuchos::dyn_cast<StdVector<Real> >(v)).getVector());
816  Real epsn = std::min(this->scale_*eps,this->min_diff_);
817  //epsn *= this->scale_;
818  for ( int i = 0; i < 2; i++ ) {
819  if ( ((*ex)[i] <= this->x_lo_[i]+epsn && (*eg)[i] > 0.0) ) {
820  (*ev)[i] = 0.0;
821  }
822  }
823  }
824  void pruneUpperActive(Vector<Real> &v, const Vector<Real> &g, const Vector<Real> &x, Real eps) {
825  Teuchos::RCP<const std::vector<Real> > ex =
826  (Teuchos::dyn_cast<StdVector<Real> >(const_cast<Vector<Real> &>(x))).getVector();
827  Teuchos::RCP<const std::vector<Real> > eg =
828  (Teuchos::dyn_cast<StdVector<Real> >(const_cast<Vector<Real> &>(g))).getVector();
829  Teuchos::RCP<std::vector<Real> > ev =
830  Teuchos::rcp_const_cast<std::vector<Real> >((Teuchos::dyn_cast<StdVector<Real> >(v)).getVector());
831  Real epsn = std::min(this->scale_*eps,this->min_diff_);
832  //epsn *= this->scale_;
833  for ( int i = 0; i < 2; i++ ) {
834  if ( ((*ex)[i] >= this->x_up_[i]-epsn && (*eg)[i] < 0.0) ) {
835  (*ev)[i] = 0.0;
836  }
837  }
838  }
839  void pruneActive(Vector<Real> &v, const Vector<Real> &g, const Vector<Real> &x, Real eps) {
840  Teuchos::RCP<const std::vector<Real> > ex =
841  (Teuchos::dyn_cast<StdVector<Real> >(const_cast<Vector<Real> &>(x))).getVector();
842  Teuchos::RCP<const std::vector<Real> > eg =
843  (Teuchos::dyn_cast<StdVector<Real> >(const_cast<Vector<Real> &>(g))).getVector();
844  Teuchos::RCP<std::vector<Real> > ev =
845  Teuchos::rcp_const_cast<std::vector<Real> >((Teuchos::dyn_cast<StdVector<Real> >(v)).getVector());
846  Real epsn = std::min(this->scale_*eps,this->min_diff_);
847  //epsn *= this->scale_;
848  for ( int i = 0; i < 2; i++ ) {
849  if ( ((*ex)[i] <= this->x_lo_[i]+epsn && (*eg)[i] > 0.0) ||
850  ((*ex)[i] >= this->x_up_[i]-epsn && (*eg)[i] < 0.0) ) {
851  (*ev)[i] = 0.0;
852  }
853  }
854  }
855 
857  Teuchos::RCP<std::vector<Real> > us = Teuchos::rcp( new std::vector<Real>(2,0.0) );
858  us->assign(this->x_up_.begin(),this->x_up_.end());
859  Teuchos::RCP<ROL::Vector<Real> > up = Teuchos::rcp( new ROL::StdVector<Real>(us) );
860  u.set(*up);
861  }
862 
864  Teuchos::RCP<std::vector<Real> > ls = Teuchos::rcp( new std::vector<Real>(2,0.0) );
865  ls->assign(this->x_lo_.begin(),this->x_lo_.end());
866  Teuchos::RCP<ROL::Vector<Real> > lp = Teuchos::rcp( new ROL::StdVector<Real>(ls) );
867  l.set(*lp);
868  }
869  };
870 
871 
872 
873  // template<class Real>
874  // void getDiodeCircuit( Teuchos::RCP<Objective<Real> > &obj, Vector<Real> &x0, Vector<Real> &x ) {
875  // // Cast Initial Guess and Solution Vectors
876  // Teuchos::RCP<std::vector<Real> > x0p =
877  // Teuchos::rcp_const_cast<std::vector<Real> >((Teuchos::dyn_cast<StdVector<Real> >(x0)).getVector());
878  // Teuchos::RCP<std::vector<Real> > xp =
879  // Teuchos::rcp_const_cast<std::vector<Real> >((Teuchos::dyn_cast<StdVector<Real> >(x)).getVector());
880 
881  // int n = xp->size();
882 
883  // // Resize Vectors
884  // n = 2;
885  // x0p->resize(n);
886  // xp->resize(n);
887 
888  // // Instantiate Objective Function
889  // obj = Teuchos::rcp( new Objective_DiodeCircuit<Real> (0.02585,0.0,1.0,1.e-2));
890  // //ROL::Objective_DiodeCircuit<Real> obj(0.02585,0.0,1.0,1.e-2);
891 
892  // // Get Initial Guess
893  // (*x0p)[0] = 1.e-13;
894  // (*x0p)[1] = 0.2;
895 
896  // // Get Solution
897  // (*xp)[0] = 1.e-12;
898  // (*xp)[1] = 0.25;
899 
900  // }
901 
902 
903 } //end namespace ZOO
904 } //end namespace ROL
905 
906 #endif
Provides the interface to evaluate objective functions.
bool use_adjoint_
If true, use adjoint gradient computation, else compute gradient using sensitivities.
Real diodeIsRs(const Real I, const Real Vsrc, const Real Is, const Real Rs)
Second derivative of diode equation wrt Is and Rs.
Real noise_
Percentage of noise to add to measurements; if 0.0 - no noise.
Real diodeIsIs(const Real I, const Real Vsrc, const Real Is, const Real Rs)
Second derivative of diode equation wrt Is^2.
void pruneUpperActive(Vector< Real > &v, const Vector< Real > &x, Real eps)
Set variables to zero if they correspond to the upper -active set.
virtual void scale(const Real alpha)=0
Compute where .
Real min_diff_
Half of the minimum distance between upper and lower bounds.
virtual void axpy(const Real alpha, const Vector &x)
Compute where .
Definition: ROL_Vector.hpp:141
Objective_DiodeCircuit(Real Vth, std::ifstream &input_file, bool lambertw, Real noise, bool use_adjoint, int use_hessvec)
A constructor using data from given file.
Real scale_
Scaling for the epsilon margin.
BoundConstraint_DiodeCircuit(Real scale, Real lo_Is, Real up_Is, Real lo_Rs, Real up_Rs)
Real lambertWCurrent(Real Is, Real Rs, Real Vsrc)
Find currents using Lambert-W function.
int use_hessvec_
0 - use FD(with scaling), 1 - use exact implementation (with second order derivatives), 2 - use Gauss-Newton approximation (first order derivatives only)
void pruneActive(Vector< Real > &v, const Vector< Real > &x, Real eps)
Set variables to zero if they correspond to the -active set.
virtual void hessVec(Vector< Real > &hv, const Vector< Real > &v, const Vector< Real > &x, Real &tol)
Apply Hessian approximation to vector.
void solve_sensitivity_Is(Vector< Real > &sens, const Vector< Real > &I, const Vector< Real > &S)
Solve the sensitivity equation wrt Is.
void solve_adjoint(Vector< Real > &lambda, const Vector< Real > &I, const Vector< Real > &S)
Solve the adjoint equation.
void gradient(Vector< Real > &g, const Vector< Real > &S, Real &tol)
Compute the gradient of the reduced objective function either using adjoint or using sensitivities...
Teuchos::RCP< std::vector< Real > > Imeas_
Vector of measured currents in DC analysis (data)
Real Vth_
Thermal voltage (constant)
virtual Teuchos::RCP< Vector > clone() const =0
Clone to make a new (uninitialized) vector.
void setVectorToUpperBound(ROL::Vector< Real > &u)
Set the input vector to the upper bound.
void generate_plot(Real Is_lo, Real Is_up, Real Is_step, Real Rs_lo, Real Rs_up, Real Rs_step)
Generate data to plot objective function.
virtual void zero()
Set to zero vector.
Definition: ROL_Vector.hpp:155
Defines the linear algebra or vector space interface.
Definition: ROL_Vector.hpp:72
Bound constraints on optimization parameters.
Teuchos::RCP< std::vector< Real > > Vsrc_
Vector of source voltages in DC analysis (input)
void pruneLowerActive(Vector< Real > &v, const Vector< Real > &x, Real eps)
Set variables to zero if they correspond to the lower -active set.
The diode circuit problem.
Real value(const Vector< Real > &S, Real &tol)
Evaluate objective function.
void pruneLowerActive(Vector< Real > &v, const Vector< Real > &g, const Vector< Real > &x, Real eps)
Set variables to zero if they correspond to the lower -binding set.
Real diodeI(const Real I, const Real Vsrc, const Real Is, const Real Rs)
Derivative of diode equation wrt I.
Provides the std::vector implementation of the ROL::Vector interface.
std::vector< Real > x_lo_
Vector of lower bounds.
void pruneActive(Vector< Real > &v, const Vector< Real > &g, const Vector< Real > &x, Real eps)
Set variables to zero if they correspond to the -binding set.
bool isFeasible(const Vector< Real > &x)
Check if the vector, v, is feasible.
void lambertw(Real x, Real &w, int &ierr, Real &xi)
Lambert-W function for diodes.
Real Newton(const Real I, const Real Vsrc, const Real Is, const Real Rs)
Newton's method with line search.
void hessVec(Vector< Real > &hv, const Vector< Real > &v, const Vector< Real > &S, Real &tol)
Compute the Hessian-vector product of the reduced objective function.
Provides the interface to apply upper and lower bound constraints.
void solve_circuit(Vector< Real > &I, const Vector< Real > &S)
Solve circuit given optimization parameters Is and Rs.
void project(Vector< Real > &x)
Project optimization variables onto the bounds.
std::vector< Real > x_up_
Vector of upper bounds.
bool lambertw_
If true, use Lambert-W function to solve circuit, else use Newton's method.
virtual void set(const Vector &x)
Set where .
Definition: ROL_Vector.hpp:194
virtual Real norm() const =0
Returns where .
Real diode(const Real I, const Real Vsrc, const Real Is, const Real Rs)
Diode equation.
Real diodeIs(const Real I, const Real Vsrc, const Real Is, const Real Rs)
Derivative of diode equation wrt Is.
Real diodeRsRs(const Real I, const Real Vsrc, const Real Is, const Real Rs)
Second derivative of diode equation wrt Rs^2.
Objective_DiodeCircuit(Real Vth, Real Vsrc_min, Real Vsrc_max, Real Vsrc_step, Real true_Is, Real true_Rs, bool lambertw, Real noise, bool use_adjoint, int use_hessvec)
A constructor generating data.
void solve_sensitivity_Rs(Vector< Real > &sens, const Vector< Real > &I, const Vector< Real > &S)
Solve the sensitivity equation wrt Rs.
void setVectorToLowerBound(ROL::Vector< Real > &l)
Set the input vector to the lower bound.
void pruneUpperActive(Vector< Real > &v, const Vector< Real > &g, const Vector< Real > &x, Real eps)
Set variables to zero if they correspond to the upper -binding set.
Real diodeIIs(const Real I, const Real Vsrc, const Real Is, const Real Rs)
Second derivative of diode equation wrt I and Is.
Real diodeIRs(const Real I, const Real Vsrc, const Real Is, const Real Rs)
Second derivative of diode equation wrt I and Rs.
Real diodeII(const Real I, const Real Vsrc, const Real Is, const Real Rs)
Second derivative of diode equation wrt I^2.
static const double ROL_EPSILON
Platform-dependent machine epsilon.
Definition: ROL_Types.hpp:115
Real diodeRs(const Real I, const Real Vsrc, const Real Is, const Real Rs)
Derivative of diode equation wrt Rs.
void set_method(bool lambertw)
Change the method for solving the circuit if needed.