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