Xyce  6.1
N_TIA_NoTimeIntegration.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_TIA_NoTimeIntegration.C,v $
27 //
28 // Purpose : This file contains the functions which define the
29 // time integration methods classes.
30 //
31 // Special Notes :
32 //
33 // Creator : Eric Keiter
34 //
35 // Creation Date : 7/21/00
36 //
37 // Revision Information:
38 // ---------------------
39 //
40 // Revision Number: $Revision: 1.64 $
41 //
42 // Revision Date : $Date: 2015/09/03 22:11:40 $
43 //
44 // Current Owner : $Author: rlschie $
45 //-----------------------------------------------------------------------------
46 
47 #include <Xyce_config.h>
48 
49 // ---------- Standard Includes ----------
50 #include <iostream>
51 
52 // ---------- Xyce Includes ----------
53 #include <N_ERH_ErrorMgr.h>
54 #include <N_IO_InitialConditions.h>
55 #include <N_LAS_Matrix.h>
56 #include <N_LAS_System.h>
57 #include <N_LAS_Vector.h>
58 #include <N_TIA_DataStore.h>
60 #include <N_TIA_StepErrorControl.h>
61 #include <N_TIA_TIAParams.h>
63 #include <N_TIA_TwoLevelError.h>
64 #include <N_ANP_OutputMgrAdapter.h>
65 #include <N_UTL_Diagnostic.h>
66 #include <N_UTL_FeatureTest.h>
67 
68 //Need Epetra Stuff for now, since directly manipulating Epetra Objs for
69 //obtainJacobian. This will be abstracted later
70 #include <Epetra_Import.h>
71 #include <Epetra_CrsMatrix.h>
72 #include <Epetra_Map.h>
73 
74 namespace Xyce {
75 namespace TimeIntg {
76 
77 const char *
79 
80 //-----------------------------------------------------------------------------
81 // Function : NoTimeIntegration::factory
82 // Purpose :
83 // Special Notes :
84 // Scope : public
85 // Creator : Eric Keiter, SNL, Parallel Computational Sciences
86 // Creation Date : 6/08/01
87 //-----------------------------------------------------------------------------
88 TimeIntegrationMethod *
90  const TIAParams & tia_params,
91  StepErrorControl & step_error_control,
92  DataStore & data_store)
93 {
94  return new NoTimeIntegration(tia_params, step_error_control, data_store);
95 }
96 
97 //-----------------------------------------------------------------------------
98 // Function : NoTimeIntegration::NoTimeIntegration
99 // Purpose : constructor
100 // Special Notes :
101 // Scope : public
102 // Creator : Buddy Watts, SNL
103 // Creation Date : 6/01/00
104 //-----------------------------------------------------------------------------
106  const TIAParams & tia_params,
107  StepErrorControl & secTmp,
108  DataStore & dsTmp)
110  leadingCoeff(1.0),
111  sec(secTmp),
112  ds(dsTmp)
113 {
114  leadingCoeff = 1.0;
115  alphas = -1.0;
116  return;
117 }
118 
119 //-----------------------------------------------------------------------------
120 // Function : NoTimeIntegration::~NoTimeIntegration
121 // Purpose : destructor
122 // Special Notes :
123 // Scope : public
124 // Creator : Buddy Watts, SNL
125 // Creation Date : 6/01/00
126 //-----------------------------------------------------------------------------
128 
129 
130 void
132  const TIAParams & tia_params)
133 {
136 }
137 
138 void
140  const TIAParams & tia_params)
141 {}
142 
143 //-----------------------------------------------------------------------------
144 // Function : NoTimeIntegration::obtainCorrectorDeriv
145 // Purpose : Evaluate the Corrector Derivative Formula
146 // Special Notes : For "no integration" the derivatives should always be zero.
147 // Scope : public
148 // Creator : Eric Keiter, SNL
149 // Creation Date : 01/10/01
150 //-----------------------------------------------------------------------------
152 {
153  ds.nextSolutionDerivPtr->putScalar(0.0);
154  ds.nextStateDerivPtr->putScalar(0.0);
155  ds.nextStoreLeadCurrQDerivPtr->putScalar(0.0);
156 }
157 
158 
159 //-----------------------------------------------------------------------------
160 // Function : NoTimeIntegration::obtainResidual
161 //
162 // Purpose : This function returns the residual for the steady state
163 // case.
164 //
165 // Special Notes : For "no integration" the derivatives should always be zero,
166 // so don't add in the dqdt term.
167 //
168 // This function is only called in the new-DAE case.
169 //
170 // Scope : public
171 // Creator : Eric Keiter, SNL
172 // Creation Date : 03/09/04
173 //-----------------------------------------------------------------------------
175 {
176  ds.RHSVectorPtr->linearCombo
177  (0.0,*ds.RHSVectorPtr,+1.0,*ds.daeFVectorPtr);
178  ds.RHSVectorPtr->linearCombo
179  (1.0,*ds.RHSVectorPtr,-1.0,*ds.daeBVectorPtr);
180 
181  // since the nonlinear solver is expecting a -f, scale by -1.0:
182  ds.RHSVectorPtr->scale(-1.0);
183 
184  // if voltage limiting is on, add it in:
185  if (ds.limiterFlag)
186  {
187  (ds.RHSVectorPtr)->daxpy(
188  *(ds.RHSVectorPtr), +1.0, *(ds.dFdxdVpVectorPtr));
189  }
190 }
191 
192 //-----------------------------------------------------------------------------
193 // Function : NoTimeIntegration::obtainSensitivityResiduals
194 //
195 // Purpose : This function returns the residual for the steady state
196 // case.
197 //
198 // Special Notes : For "no integration" the derivatives should always be zero,
199 // so don't add in the dqdt term.
200 //
201 // This function is only called in the new-DAE case.
202 //
203 // Scope : public
204 // Creator : Eric Keiter, SNL
205 // Creation Date :
206 //-----------------------------------------------------------------------------
208 {
209  int numParams = ds.sensRHSPtrVector.size();
210 
211  for (int ip=0; ip<numParams;++ip)
212  {
213  Linear::Vector & RHSVec = *(ds.sensRHSPtrVector[ip]);
214  Linear::Vector & dfdpVec = *(ds.nextDfdpPtrVector[ip]);
215  Linear::Vector & dbdpVec = *(ds.nextDbdpPtrVector[ip]);
216 
217  RHSVec.linearCombo(0.0,RHSVec,+1.0,dfdpVec);
218  RHSVec.linearCombo(1.0,RHSVec,-1.0,dbdpVec);
219 
220  // since the nonlinear solver is expecting a -f, scale by -1.0:
221  RHSVec.scale(-1.0);
222  }
223 }
224 
225 //-----------------------------------------------------------------------------
226 // Function : NoTimeIntegration::obtainJacobian
227 //
228 // Purpose : Returns the full Jacobian matrix for the steady state
229 // case.
230 //
231 // Special Notes : For "no integration" the derivatives should always be zero.
232 // However, to prevent singular matrices, the dQdt term is
233 // added in anyway, with an assumed very large time step.
234 //
235 // There may be a better way to handle this - the assumed
236 // large time step leads to a very small alpha/dt, which
237 // could (maybe?) have an adverse effect on the matrix
238 // conditioning.
239 //
240 // The singular matrix problem will happen for the case of
241 // capacitors in serial. Any voltage node which is *only*
242 // connected to capacitors will not really be part of the
243 // system of equations. Ideally nodes like this would just
244 // be removed from the system altogether. No current is
245 // passing through them, and they are really just floating
246 // in space.
247 //
248 // For the time being, however, what we do is put some
249 // bogus C*alpha/dt terms into the Jacobian. If dt is
250 // large, then the C*alpha/dt is very small, and doesn't
251 // significantly change any Jacobian entries, except for
252 // zero entries.
253 //
254 // This function is only called in the new-DAE case.
255 //
256 // Scope : public
257 // Creator : Eric Keiter, SNL
258 // Creation Date : 03/09/04
259 //-----------------------------------------------------------------------------
261 {
262  Linear::Matrix & dQdx = *(ds.dQdxMatrixPtr);
263  Linear::Matrix & dFdx = *(ds.dFdxMatrixPtr);
264  Linear::Matrix & Jac = *(ds.JMatrixPtr);
265 
266  Jac.linearCombo( 1.0e-20, dQdx, 1.0, dFdx );
267 }
268 
269 //-----------------------------------------------------------------------------
270 // Function : NoTimeIntegration::applyJacobian
271 //
272 // Purpose : Applies the Jacobian operator for the steady state
273 // case.
274 //
275 // Special Notes :
276 // This function is only called in the new-DAE HB (matrix-free) case.
277 //
278 // Scope : public
279 // Creator : Todd Coffey, Ting Mei
280 // Creation Date : 07/29/08
281 //-----------------------------------------------------------------------------
282 void NoTimeIntegration::applyJacobian (const Linear::Vector& input, Linear::Vector& result)
283 {
284  Linear::Vector & dQdxV = *(ds.dQdxVecVectorPtr);
285  Linear::Vector & dFdxV = *(ds.dFdxVecVectorPtr);
286  result.linearCombo( 1.0e-20, dQdxV, 1.0, dFdxV );
287 }
288 
289 //-----------------------------------------------------------------------------
290 // Function : NoTimeIntegration::getInitialQnorm
291 // Purpose : Needed by 2-level solves.
292 // Special Notes :
293 // Scope : public
294 // Creator : Eric Keiter, SNL, Parallel Computational Sciences
295 // Creation Date : 03/18/07
296 //-----------------------------------------------------------------------------
298 {
299  tle.q1HistorySum = 0.0;
300 }
301 
302 //-----------------------------------------------------------------------------
303 // Function : NoTimeIntegration::setupTwoLevelError
304 // Purpose : Needed by 2-level solves.
305 // Special Notes :
306 // Scope : public
307 // Creator : Eric Keiter, SNL, Parallel Computational Sciences
308 // Creation Date : 03/15/07
309 //-----------------------------------------------------------------------------
311 {
312  tle.xErrorSum = 0.0;
313  tle.qErrorSum = 0.0;
314  tle.xErrorSum_m1 = 0.0;
315  tle.xErrorSum_m2 = 0.0;
316  tle.innerSize = ds.globalLength ();
317 }
318 
320  Analysis::OutputMgrAdapter & outputManagerAdapter,
321  const TIAParams & tia_params,
322  const double time,
323  Linear::Vector * solnVecPtr,
324  const bool doNotInterpolate,
325  const std::vector<double> & outputInterpolationTimes,
326  bool skipPrintLineOutput )
327 {
328  if (DEBUG_TIME && isActive(Diag::TIME_OUTPUT))
329  Xyce::dout() << "Calling conventional outputs!" << std::endl;
330 
331  outputManagerAdapter.tranOutput( time, *solnVecPtr,
333  ds.objectiveVec_,
336  skipPrintLineOutput);
337 
338  return true;
339 }
340 
341 //-----------------------------------------------------------------------------
342 // Function : TimeIntegrationMethod::saveOutputSolution
343 // Purpose :
344 // Special Notes : For the old method functions (old-DAE/ODE) no interpolation
345 // is possible, so this function calls directly through to the
346 // output manager.
347 // Scope : public
348 // Creator : Eric Keiter, SNL, 1437
349 // Creation Date : 10/25/07
350 //-----------------------------------------------------------------------------
352  Parallel::Machine comm,
353  IO::InitialConditionsManager & initial_conditions_manager,
354  const NodeNameMap & node_name_map,
355  const TIAParams & tia_params,
356  Linear::Vector * solnVecPtr,
357  const double saveTime,
358  const bool doNotInterpolate)
359 {
360  // outputManagerAdapter.outputDCOP( *(solnVecPtr) );
361  initial_conditions_manager.outputDCOP(comm, node_name_map, *solnVecPtr);
362 
363  return true;
364 }
365 
366 } // namespace TimeIntg
367 } // namespace Xyce
Linear::Vector * tmpLeadDeltaVPtr
void applyJacobian(const Linear::Vector &input, Linear::Vector &result)
StepErrorControl & sec
Reference to step-error control object.
bool saveOutputSolution(Parallel::Machine comm, IO::InitialConditionsManager &initial_conditions_manager, const NodeNameMap &node_name_map, const TIAParams &tia_params, Linear::Vector *solnVecPtr, const double saveTime, const bool doNotInterpolate)
Pure virtual class to augment a linear system.
Linear::Vector * nextSolutionDerivPtr
std::vector< double > dOdpVec_
std::vector< double > scaled_dOdpVec_
void completeStep(const TIAParams &tia_params)
Linear::Vector * currStorePtr
std::vector< double > scaled_dOdpAdjVec_
std::vector< Linear::Vector * > nextDbdpPtrVector
std::vector< Linear::Vector * > nextDfdpPtrVector
bool printOutputSolution(Analysis::OutputMgrAdapter &outputManagerAdapter, const TIAParams &tia_params, const double time, Linear::Vector *solnVecPtr, const bool doNotInterpolate, const std::vector< double > &outputInterpolationTimes, bool skipPrintLineOutput)
Linear::Vector * currStatePtr
std::vector< double > objectiveVec_
double alphas
$$ fixed-leading coefficient of this BDF method
Linear::Vector * tmpLeadCurrentVectorPtr
NoTimeIntegration(const TIAParams &tiaP, StepErrorControl &secTmp, DataStore &dsTmp)
Linear::Vector * daeFVectorPtr
double leadingCoeff
Time-integration method leading coefficient value.
Linear::Vector * nextStoreLeadCurrQDerivPtr
Linear::Matrix * dFdxMatrixPtr
Linear::Matrix * JMatrixPtr
Linear::Vector * dFdxdVpVectorPtr
void rejectStep(const TIAParams &tia_params)
Linear::Vector * RHSVectorPtr
DataStore & ds
Reference to the TIA data-store object.
std::vector< double > dOdpAdjVec_
Linear::Vector * dQdxVecVectorPtr
Linear::Vector * tmpLeadCurrentQDerivVectorPtr
void getTwoLevelError(TwoLevelError &tle) const
Linear::Vector * dFdxVecVectorPtr
Linear::Vector * nextStateDerivPtr
std::vector< Linear::Vector * > sensRHSPtrVector
Linear::Matrix * dQdxMatrixPtr
void tranOutput(double time, Linear::Vector &currSolutionPtr, Linear::Vector &stateVecPtr, Linear::Vector &storeVecPtr, Linear::Vector &lead_current_vector, Linear::Vector &junction_voltage_vector, Linear::Vector &lead_current_dqdt_vector, std::vector< double > &objectiveVec_, std::vector< double > &dOdpVec_, std::vector< double > &dOdpAdjVec_, std::vector< double > &scaled_dOdpVec_, std::vector< double > &scaled_dOdpAdjVec_, bool skipPrintLineOutput=false)
void getInitialQnorm(TwoLevelError &tle) const
Linear::Vector * daeBVectorPtr
static TimeIntegrationMethod * factory(const TIAParams &tia_params, StepErrorControl &step_error_control, DataStore &data_store)