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