Xyce  6.1
N_ANP_ModelEvaluator.C
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-2015 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_ModelEvaluator.C,v $
27 // Purpose : This file contains the functions which define the
28 // EpetraExt::ModelEvaluator interface for Xyce
29 // Special Notes :
30 // Creator : Todd Coffey, 1414
31 // Creation Date : early 2009
32 //
33 // Revision Information:
34 // ---------------------
35 // Revision Number: $Revision: 1.32 $
36 // Revision Date : $Date: 2015/04/08 19:18:20 $
37 // Current Owner : $Author: tvrusso $
38 //-----------------------------------------------------------------------------
39 
40 #include <Xyce_config.h>
41 
42 #include <N_ANP_ModelEvaluator.h>
43 #include <N_CIR_Xyce.h>
44 #include <N_LAS_Vector.h>
45 #include <N_LAS_BlockVector.h>
46 #include <N_LAS_BlockMatrix.h>
47 #include <N_LAS_BlockSystemHelpers.h>
48 #include <N_PDS_ParMap.h>
49 
50 #include <Epetra_Map.h>
51 #include <Epetra_Vector.h>
52 #include <Epetra_CrsGraph.h>
53 #include <Epetra_CrsMatrix.h>
54 #include <Epetra_Operator.h>
55 
56 using Teuchos::rcp;
57 using Teuchos::is_null;
58 
59 namespace Xyce {
60 namespace Analysis {
61 
62 RCP<Linear::BlockVector> convertEpetraToNLASBlockVectorView(
63  const RCP<const Epetra_Vector>& vec,
64  const RCP<N_PDS_ParMap> & map
65  )
66 {
67  // Add check that the size of vec is the same as 2*(size of map)
68  int Nv = vec->GlobalLength();
69  int Nm = map->numGlobalEntities();
70  TEUCHOS_TEST_FOR_EXCEPTION( Nv != 2*Nm, std::logic_error,
71  "Error! The size of the vector must be twice the size of the map!"
72  );
73  Epetra_Vector* nonconst_vec = const_cast<Epetra_Vector*>(&*vec);
74  RCP<Linear::BlockVector> nlasBlockVec = rcp(new Linear::BlockVector( nonconst_vec, map, 2, false ));
75  return nlasBlockVec;
76 }
77 
78 
79 RCP<Linear::Vector> convertEpetraToNLASVectorView(
80  const RCP<const Epetra_Vector>& vec
81  )
82 {
83  Epetra_Vector* nonconst_vec = const_cast<Epetra_Vector*>(&*vec);
84  RCP<Linear::Vector> nlasVec = rcp(new Linear::Vector( nonconst_vec, false ));
85  return nlasVec;
86 }
87 
88 
90  :isInitialized_(false)
91 {
92 }
93 
94 
95 void ModelEvaluator::initialize(int iargs, char* cargs[])
96 {
97  xycePtr_ = rcp(new Circuit::Simulator());
98  xycePtr_->initialize(iargs,cargs);
99  xycePtr_->initializeTransientModel();
100  this->setupInOutArgs_();
101  this->setupMapsAndGraphs_(); // must be called after xycePtr_->initialize
102  //voltLimFVector_ = rcp(new Linear::Vector(*x_map_,*x_ognd_map_));
103  //voltLimQVector_ = rcp(new Linear::Vector(*x_map_,*x_ognd_map_));
104 
105  x_gnd_ = rcp(new Linear::Vector(*x_map_,*x_ognd_map_));
106  z_gnd_ = rcp(new Linear::Vector(*x_map_,*x_ognd_map_));
107  xdot_gnd_ = rcp(new Linear::Vector(*x_map_,*x_ognd_map_));
108  zdot_gnd_ = rcp(new Linear::Vector(*x_map_,*x_ognd_map_));
109  f_0_gnd_ = rcp(new Linear::Vector(*x_map_,*x_ognd_map_));
110  f_0_gnd_.release(); // Xyce deletes this pointer in TimeIntg::DataStore
111  f_1_gnd_ = rcp(new Linear::Vector(*x_map_,*x_ognd_map_));
112  f_1_gnd_.release(); // Xyce deletes this pointer in TimeIntg::DataStore
113  eVec_ = rcp(new Linear::Vector(*x_map_,*x_ognd_map_));
114  eVec_->putScalar(1.0);
115 
116  dQdx_gnd_matrix_ = rcp(new Linear::Matrix( &*dQdx_ognd_graph_, &*dQdx_graph_ ));
117  dQdx_gnd_matrix_.release(); // Xyce deletes this pointer in TimeIntg::DataStore
118  dFdx_gnd_matrix_ = rcp(new Linear::Matrix( &*dFdx_ognd_graph_, &*dFdx_graph_ ));
119  dFdx_gnd_matrix_.release(); // Xyce deletes this pointer in TimeIntg::DataStore
120 
121  isInitialized_ = true;
122 }
123 
124 
126 {
127  //xycePtr_->finalize();
128 }
129 
130 
131 std::vector<std::string> ModelEvaluator::getVariableNames()
132 {
133  return xycePtr_->getVariableNames();
134 }
135 
136 
138 {
139  Np_ = 2;
140  Ng_ = 3;
141  EpetraExt::ModelEvaluator::InArgsSetup inArgs;
142  inArgs.setSupports(IN_ARG_t,true);
143  inArgs.setSupports(IN_ARG_x,true);
144  inArgs.setSupports(IN_ARG_x_dot,true);
145  inArgs.setSupports(IN_ARG_alpha,true);
146  inArgs.setSupports(IN_ARG_beta,true);
147  inArgs.set_Np(Np_);
148  inArgs_ = inArgs;
149 
150  EpetraExt::ModelEvaluator::OutArgsSetup outArgs;
151  outArgs.setSupports(OUT_ARG_f,true);
152  outArgs.setSupports(OUT_ARG_W,true);
153  outArgs.set_Np_Ng(Np_,Ng_);
154  outArgs_ = outArgs;
155 }
156 
157 
159 {
161 
162  // create block map for conversion between Xyce and Trilinos
163  int Size = 2; // two blocks
164  int BlockSize = x_map_->numGlobalEntities();
165  blockMap_ = Linear::createBlockParMap( Size, *x_map_ );
166 
167  // Augment the base graph to include the diagonal entries.
169  rcp(new Epetra_CrsGraph(
170  Copy,
171  *(x_map_->petraBlockMap()),
172  0
173  )
174  );
175  // Copy graph non-zeros from dQdx and dFdx graphs:
176  // W = [ dFdx, I ]
177  // [ dQdx, I ]
178  int MaxIndices = dFdx_graph_->MaxNumIndices();
179  std::vector<int> Indices;
180  Indices.resize(MaxIndices+1);
181  int NumIndices;
182  int BaseRow;
183  for( int j = 0; j < BlockSize; ++j )
184  {
185  BaseRow = x_map_->localToGlobalIndex(j);
186  dFdx_graph_->ExtractGlobalRowCopy( BaseRow, MaxIndices, NumIndices, &Indices[0] );
187  // Here we add the diagonal:
188  Indices[NumIndices] = BaseRow;
189  NumIndices++;
190  dFdx_graph_with_diagonal_->InsertGlobalIndices( BaseRow, NumIndices, &Indices[0] );
191  }
192  dFdx_graph_with_diagonal_->TransformToLocal();
193 }
194 
195 
196 EpetraExt::ModelEvaluator::InArgs ModelEvaluator::createInArgs() const
197 {
198  return inArgs_;
199 }
200 
201 
202 EpetraExt::ModelEvaluator::OutArgs ModelEvaluator::createOutArgs() const
203 {
204  return outArgs_;
205 }
206 
207 
208 RCP<const Epetra_Map> ModelEvaluator::get_x_map() const
209 {
210  return rcp( blockMap_->petraMap(), false );
211 }
212 
213 
214 RCP<const Epetra_Map> ModelEvaluator::get_f_map() const
215 {
216  return rcp( blockMap_->petraMap(), false );
217 }
218 
219 
220 RCP<const Epetra_Map> ModelEvaluator::get_p_map(int p) const
221 {
222  TEUCHOS_ASSERT( ((0 <= p) && (p < Np_)) );
223  return Teuchos::rcpFromRef(*s_map_->petraMap());
224 }
225 
226 
227 RCP<const Epetra_Map> ModelEvaluator::get_g_map(int p) const
228 {
229  TEUCHOS_ASSERT( ((0 <= p) && (p < Ng_)) );
230  if (p == 0) { // s
231  return Teuchos::rcpFromRef(*s_map_->petraMap());
232  }
233  else if (p == 1) { // voltLimQ
234  return rcp( blockMap_->petraMap(), false );
235  }
236  else if (p == 2) { // voltLim F
237  return rcp( blockMap_->petraMap(), false );
238  }
239  return Teuchos::null; // Should never get here.
240 }
241 
242 
243 RCP<const Epetra_Map> ModelEvaluator::get_small_x_map() const
244 {
245  return Teuchos::rcpFromRef(*x_map_->petraMap());
246 }
247 
248 
249 void ModelEvaluator::evalModel( const InArgs& inArgs, const OutArgs& outArgs ) const
250 {
251  TEUCHOS_TEST_FOR_EXCEPTION( !isInitialized_, std::logic_error,
252  "Error! Please call initialize before evalModel"
253  );
254  // Get data out of inArgs
255  double t = inArgs.get_t();
256  RCP<const Epetra_Vector> x_in = inArgs.get_x().assert_not_null();
257  RCP<const Epetra_Vector> xdot_in = inArgs.get_x_dot().assert_not_null();
258  RCP<Epetra_Vector> s_out = outArgs.get_g(0);
259  RCP<Epetra_Vector> store_out = outArgs.get_g(0); // ???
260 
261  RCP<Linear::BlockVector> xyce_x = convertEpetraToNLASBlockVectorView(x_in,x_map_);
262  RCP<Linear::BlockVector> xyce_xdot = convertEpetraToNLASBlockVectorView(xdot_in,x_map_);
263 
264  // xyce_x = (x,z)
265  Linear::Vector & x = xyce_x->block(0);
266  Linear::Vector & z = xyce_x->block(1);
267  Linear::Vector & xdot = xyce_xdot->block(0);
268  Linear::Vector & zdot = xyce_xdot->block(1);
269 
270  // The base Linear::Vector map has an extra element at the end for ground and
271  // we're choosing for the moment to copy into this type of data structure
272  // rather than fix our blockMap.
273  x_gnd_->update(1.0,x,0.0); // x_gnd_ = x
274  z_gnd_->update(1.0,z,0.0); // z_gnd_ = z
275  xdot_gnd_->update(1.0,xdot,0.0); // xdot_gnd_ = xdot
276  zdot_gnd_->update(1.0,zdot,0.0); // zdot_gnd_ = zdot
277 
278  Linear::Vector & x_ref = *x_gnd_;
279  Linear::Vector & z_ref = *z_gnd_;
280  //Linear::Vector & xdot_ref = *xdot_gnd_;
281  Linear::Vector & zdot_ref = *zdot_gnd_;
282 
283  if (!Teuchos::is_null(s_out)) {
284  // We're loading the state vectors
285  RCP<Linear::Vector> xyce_s_out = convertEpetraToNLASVectorView(s_out);
286  Linear::Vector & s_out = *xyce_s_out;
287 
288  RCP<Linear::Vector> xyce_store_out = convertEpetraToNLASVectorView(store_out);
289  Linear::Vector & store_out = *xyce_store_out;
290 
291  xycePtr_->evalTransientModelState(
292  t,
293  & x_ref,
294  & s_out,
295  & store_out
296  );
297 
298  } else {
299  // We're loading f and W
300  double alpha = inArgs.get_alpha();
301  double beta = inArgs.get_beta();
302  RCP<const Epetra_Vector> s_in = inArgs.get_p(0).assert_not_null();
303  RCP<const Epetra_Vector> sdot_in = inArgs.get_p(1).assert_not_null();
304 
305  RCP<Epetra_Vector> f_out = outArgs.get_f().assert_not_null();
306  RCP<Epetra_Operator> W_out = outArgs.get_W().assert_not_null();
307  RCP<Linear::BlockMatrix> bMat;
308  {
309  // Pull out dFdx and dQdx matrices:
310  RCP<Epetra_CrsMatrix> crsMat = Teuchos::rcp_dynamic_cast<Epetra_CrsMatrix>(W_out,true);
311  std::string label = "Linear::BlockMatrix";
312  bMat = Teuchos::get_extra_data<RCP<Linear::BlockMatrix> >(crsMat,label);
313  //Linear::Matrix& dFdx_tmp = bMat->block(0,0);
314  //Linear::Matrix& dQdx_tmp = bMat->block(1,0);
315 
316  // Copy dFdx_tmp into dFdx_gnd_matrix_
317  dFdx_gnd_matrix_->put(0.0);
318  // Commented out
319  //dFdx_gnd_matrix_->add(dFdx_tmp);
320 
321  // Copy dQdx_tmp into dQdx_gnd_matrix_
322  dQdx_gnd_matrix_->put(0.0);
323  // Commented out
324  //dQdx_gnd_matrix_->add(dQdx_tmp);
325 
326  }
327  RCP<Epetra_Vector> voltLimQVector = outArgs.get_g(1);
328  RCP<Epetra_Vector> voltLimFVector = outArgs.get_g(2);
329 
330  RCP<Linear::Vector> xyce_s = convertEpetraToNLASVectorView(s_in);
331  RCP<Linear::Vector> xyce_sdot = convertEpetraToNLASVectorView(sdot_in);
332  RCP<Linear::BlockVector> xyce_f = convertEpetraToNLASBlockVectorView(f_out,x_map_);
333  RCP<Linear::Vector> xyce_store = convertEpetraToNLASVectorView(s_in);
334  // 08/05/09 tscoffe: This is a hack to avoid having to deal with voltage
335  // limiting for the moment.
336  RCP<Linear::Vector> xyce_voltLimQ;
337  if (Teuchos::is_null(voltLimQVector)) {
338  if (Teuchos::is_null(tempVoltLimQVector_)) {
339  tempVoltLimQVector_ = rcp(new Linear::Vector(*x_map_,*x_ognd_map_));
340  }
341  xyce_voltLimQ = tempVoltLimQVector_;
342  }
343  else {
344  xyce_voltLimQ = convertEpetraToNLASVectorView(voltLimQVector);
345  }
346  RCP<Linear::Vector> xyce_voltLimF;
347  if (Teuchos::is_null(voltLimFVector)) {
348  if (Teuchos::is_null(tempVoltLimFVector_)) {
349  tempVoltLimFVector_ = rcp(new Linear::Vector(*x_map_,*x_ognd_map_));
350  }
351  xyce_voltLimF = tempVoltLimFVector_;
352  }
353  else {
354  xyce_voltLimF = convertEpetraToNLASVectorView(voltLimFVector);
355  }
356  Linear::Vector & s = *xyce_s;
357  Linear::Vector & store = *xyce_store;
358  Linear::Vector & sdot = *xyce_sdot;
359  Linear::Vector & f_0 = xyce_f->block(0);
360  Linear::Vector & f_1 = xyce_f->block(1);
361  //f_0_gnd_->update(1.0,f_0,0.0); // f_0_gnd_ = f_0
362  //f_1_gnd_->update(1.0,f_1,0.0); // f_1_gnd_ = f_1
363  f_0_gnd_->putScalar(0.0);
364  f_1_gnd_->putScalar(0.0);
365  Linear::Vector & f_0_ref = *f_0_gnd_;
366  Linear::Vector & f_1_ref = *f_1_gnd_;
367  Linear::Vector & voltLimQ_ref = *xyce_voltLimQ;
368  Linear::Vector & voltLimF_ref = *xyce_voltLimF;
369 
370  // Eval model through AnalysisManager
371  // F(t,x,xdot,s,sdot) = [ zdot + f, q - z ]^T
372  xycePtr_->evalTransientModel(
373  t,
374  & x_ref,
375  & x_ref,
376  & x_ref,
377  & s,
378  & s,
379  & s,
380  & sdot,
381  & store,
382  & store,
383  & store,
384  & store, // this should e storLeadCurrQ !
385  & store, // this should next lead current !
386  & store, // this should current lead current!
387  & store, // this should last lead current !
388  & store, // this should lead current Q component !
389  & store, // this should next junction voltage !
390  & store, // this should current junction voltage !
391  & store, // this should last junction voltage !
392  & f_1_ref, // f_1 = q
393  & f_0_ref, // f_0 = f
394  & f_0_ref, // THIS WILL NOT WORK. It needs to be the B-vector.
395  & voltLimF_ref,
396  & voltLimQ_ref,
397  &* dQdx_gnd_matrix_,
399  );
400 
401  // Assemble outbound residual:
402  // f_0 = zdot + f
403  // f_1 = q - z
404  f_0_ref.update(1.0,zdot_ref,+1.0);
405  f_1_ref.update(-1.0,z_ref,1.0);
406 
407  // Now we copy the data back into our client's Epetra_Vectors
408  f_0.update(1.0,f_0_ref,0.0); // f_0 = f_0_ref
409  f_1.update(1.0,f_1_ref,0.0); // f_1 = f_1_ref
410 
411  {
412  // Assemble the outbound Jacobian:
413  // W = [ beta*dFdx_gnd_matrix_, alpha*I ]
414  // [ beta*dQdx_gnd_matrix_, -beta*I ]
415  Linear::Matrix & ul = bMat->block(0,0);
416  ul.put(0.0);
417  ul.add(*dFdx_gnd_matrix_);
418  ul.scale(beta);
419  Linear::Matrix & ll = bMat->block(1,0);
420  ll.put(0.0);
421  ll.add(*dQdx_gnd_matrix_);
422  ll.scale(beta);
423  Linear::Matrix & ur = bMat->block(0,1);
424  Linear::Matrix & lr = bMat->block(1,1);
425  ur.put(0.0);
426  lr.put(0.0);
427  ur.replaceDiagonal(*eVec_);
428  ur.scale(alpha);
429  lr.replaceDiagonal(*eVec_);
430  lr.scale(-beta);
431 
432  // Copy all data from blocks to global matrix if necessary.
433  bMat->assembleGlobalMatrix();
434  }
435  }
436 }
437 
438 
439 RCP<Epetra_Operator> ModelEvaluator::create_W() const
440 {
441  // Define which blocks of the blockMatrix are non-zero (ALL):
442  int Size = 2;
443  std::vector< std::vector<int> > Cols;
444  Cols.resize(Size);
445  for (int i=0 ; i<Size ; ++i) {
446  Cols[i].resize(Size);
447  for (int j=0 ; j<Size ; ++j) {
448  Cols[i][j] = j;
449  }
450  }
451 
452  int MaxGID = x_map_->maxGlobalEntity();
453  int offset=1;
454  while ( offset <= MaxGID ) offset *= 10;
455 
456  RCP<Epetra_CrsGraph> blockGraph = Linear::createBlockGraph( offset, Cols, *blockMap_, *dFdx_graph_with_diagonal_ );
457 
458  // Create the Linear::BlockMatrix:
459  RCP<Linear::BlockMatrix> bMat = rcp( new Linear::BlockMatrix( Size, offset, Cols, *blockGraph, *dFdx_graph_with_diagonal_ ) );
460  bMat.release(); // will leak memory now.
461 
462  // Get the underlying Epetra_CrsMatrix out of this
463  Epetra_CrsMatrix & crsMat = bMat->epetraObj();
464  // We attach the RCP<Linear::BlockMatrix> to the RCP object for the Epetra_Operator, so its not deleted.
465  // NOTE: Since the blocks of bMat are NOT views of the underlying Epetra_CrsMatrix, one needs to be careful.
466  RCP<Epetra_CrsMatrix> outMatrix = rcp(&crsMat,false); // do not delete the underlying epetraObj()
467  std::string label = "Linear::BlockMatrix";
468  Teuchos::set_extra_data( bMat, label, Teuchos::outArg(outMatrix) );
469 
470  return outMatrix;
471 }
472 
473 
475 {
476  return isInitialized_;
477 }
478 
479 } // namespace Analysis
480 } // namespace Xyce
EpetraExt::ModelEvaluator::InArgs createInArgs() const
Teuchos::RCP< const Epetra_Map > get_p_map(int p) const
RCP< Linear::Vector > convertEpetraToNLASVectorView(const RCP< const Epetra_Vector > &vec)
RCP< Epetra_CrsGraph > dFdx_graph_
EpetraExt::ModelEvaluator::OutArgs outArgs_
Pure virtual class to augment a linear system.
Teuchos::RCP< Epetra_Operator > create_W() const
EpetraExt::ModelEvaluator::InArgs inArgs_
RCP< Epetra_CrsGraph > dQdx_graph_
RCP< Circuit::Simulator > xycePtr_
void initialize(int iargs, char *cargs[])
RCP< Linear::Matrix > dQdx_gnd_matrix_
std::vector< std::string > getVariableNames()
RCP< Linear::Vector > tempVoltLimFVector_
Teuchos::RCP< const Epetra_Map > get_f_map() const
RCP< Linear::Vector > tempVoltLimQVector_
Teuchos::RCP< const Epetra_Map > get_small_x_map() const
RCP< Epetra_CrsGraph > dFdx_ognd_graph_
RCP< Epetra_CrsGraph > dQdx_ognd_graph_
EpetraExt::ModelEvaluator::OutArgs createOutArgs() const
Teuchos::RCP< const Epetra_Map > get_g_map(int p) const
void evalModel(const InArgs &inArgs, const OutArgs &outArgs) const
RCP< Linear::BlockVector > convertEpetraToNLASBlockVectorView(const RCP< const Epetra_Vector > &vec, const RCP< Epetra_Map > &map)
Teuchos::RCP< const Epetra_Map > get_x_map() const
RCP< Epetra_CrsGraph > dFdx_graph_with_diagonal_
RCP< Linear::Matrix > dFdx_gnd_matrix_