Xyce  6.1
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
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-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_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.46.2.1 $
41 //
42 // Revision Date : $Date: 2014/08/19 16:28:08 $
43 //
44 // Current Owner : $Author: erkeite $
45 //-----------------------------------------------------------------------------
46 
47 #include <Xyce_config.h>
48 
49 
50 // ---------- Standard Includes ----------
51 #include <N_UTL_Misc.h>
52 
53 #ifdef HAVE_IOSTREAM
54 #include <iostream>
55 #else
56 #include <iostream.h>
57 #endif
58 
59 // ---------- Xyce Includes ----------
62 #include <N_TIA_StepErrorControl.h>
63 #include <N_TIA_DataStore.h>
64 #include <N_TIA_TIAParams.h>
65 #include <N_TIA_TwoLevelError.h>
66 #include <N_ERH_ErrorMgr.h>
67 
68 #include <N_LAS_Vector.h>
69 #include <N_LAS_Matrix.h>
70 #include <N_LAS_System.h>
71 
72 //Need Epetra Stuff for now, since directly manipulating Epetra Objs for
73 //obtainJacobian. This will be abstracted later
74 #include <Epetra_Import.h>
75 #include <Epetra_CrsMatrix.h>
76 #include <Epetra_Map.h>
77 
78 //-----------------------------------------------------------------------------
79 // Function : N_TIA_NoTimeIntegration::N_TIA_NoTimeIntegration
80 // Purpose : constructor
81 // Special Notes :
82 // Scope : public
83 // Creator : Buddy Watts, SNL
84 // Creation Date : 6/01/00
85 //-----------------------------------------------------------------------------
87  (N_TIA_TIAParams & tiaP,
88  N_TIA_StepErrorControl & secTmp,
89  N_TIA_DataStore & dsTmp)
90 
91  : N_TIA_TimeIntegrationMethod(tiaP,secTmp,dsTmp)
92 {
93  leadingCoeff = 1.0;
94  alphas = -1.0;
95  return;
96 }
97 
98 //-----------------------------------------------------------------------------
99 // Function : N_TIA_NoTimeIntegration::~N_TIA_NoTimeIntegration
100 // Purpose : destructor
101 // Special Notes :
102 // Scope : public
103 // Creator : Buddy Watts, SNL
104 // Creation Date : 6/01/00
105 //-----------------------------------------------------------------------------
107 
108 
109 //-----------------------------------------------------------------------------
110 // Function : N_TIA_NoTimeIntegration::obtainCorrectorDeriv
111 // Purpose : Evaluate the Corrector Derivative Formula
112 // Special Notes : For "no integration" the derivatives should always be zero.
113 // Scope : public
114 // Creator : Eric Keiter, SNL
115 // Creation Date : 01/10/01
116 //-----------------------------------------------------------------------------
118 {
119  ds.nextSolutionDerivPtr->putScalar(0.0);
120  ds.nextStateDerivPtr->putScalar(0.0);
121  ds.nextStoreLeadCurrQDerivPtr->putScalar(0.0);
122 }
123 
124 
125 //-----------------------------------------------------------------------------
126 // Function : N_TIA_NoTimeIntegration::obtainResidual
127 //
128 // Purpose : This function returns the residual for the steady state
129 // case.
130 //
131 // Special Notes : For "no integration" the derivatives should always be zero,
132 // so don't add in the dqdt term.
133 //
134 // This function is only called in the new-DAE case.
135 //
136 // Scope : public
137 // Creator : Eric Keiter, SNL
138 // Creation Date : 03/09/04
139 //-----------------------------------------------------------------------------
141 {
142  ds.RHSVectorPtr->linearCombo
143  (0.0,*ds.RHSVectorPtr,+1.0,*ds.daeFVectorPtr);
144 #ifdef SEPARATE_F_AND_B
145  ds.RHSVectorPtr->linearCombo
146  (1.0,*ds.RHSVectorPtr,-1.0,*ds.daeBVectorPtr);
147 #endif
148 
149  // since the nonlinear solver is expecting a -f, scale by -1.0:
150  ds.RHSVectorPtr->scale(-1.0);
151 
152  // if voltage limiting is on, add it in:
153  if (ds.limiterFlag)
154  {
155  (ds.RHSVectorPtr)->daxpy(
156  *(ds.RHSVectorPtr), +1.0, *(ds.dFdxdVpVectorPtr));
157  }
158 }
159 
160 //-----------------------------------------------------------------------------
161 // Function : N_TIA_NoTimeIntegration::obtainSensitivityResiduals
162 //
163 // Purpose : This function returns the residual for the steady state
164 // case.
165 //
166 // Special Notes : For "no integration" the derivatives should always be zero,
167 // so don't add in the dqdt term.
168 //
169 // This function is only called in the new-DAE case.
170 //
171 // Scope : public
172 // Creator : Eric Keiter, SNL
173 // Creation Date :
174 //-----------------------------------------------------------------------------
176 {
177  int numParams = ds.sensRHSPtrVector.size();
178 
179  for (int ip=0; ip<numParams;++ip)
180  {
181  N_LAS_Vector & RHSVec = *(ds.sensRHSPtrVector[ip]);
182  N_LAS_Vector & dfdpVec = *(ds.nextDfdpPtrVector[ip]);
183  N_LAS_Vector & dbdpVec = *(ds.nextDbdpPtrVector[ip]);
184 
185  RHSVec.linearCombo(0.0,RHSVec,+1.0,dfdpVec);
186 #ifdef SEPARATE_F_AND_B
187  RHSVec.linearCombo(1.0,RHSVec,-1.0,dbdpVec);
188 #endif
189 
190  // since the nonlinear solver is expecting a -f, scale by -1.0:
191  RHSVec.scale(-1.0);
192  }
193 }
194 
195 //-----------------------------------------------------------------------------
196 // Function : N_TIA_NoTimeIntegration::obtainJacobian
197 //
198 // Purpose : Returns the full Jacobian matrix for the steady state
199 // case.
200 //
201 // Special Notes : For "no integration" the derivatives should always be zero.
202 // However, to prevent singular matrices, the dQdt term is
203 // added in anyway, with an assumed very large time step.
204 //
205 // There may be a better way to handle this - the assumed
206 // large time step leads to a very small alpha/dt, which
207 // could (maybe?) have an adverse effect on the matrix
208 // conditioning.
209 //
210 // The singular matrix problem will happen for the case of
211 // capacitors in serial. Any voltage node which is *only*
212 // connected to capacitors will not really be part of the
213 // system of equations. Ideally nodes like this would just
214 // be removed from the system altogether. No current is
215 // passing through them, and they are really just floating
216 // in space.
217 //
218 // For the time being, however, what we do is put some
219 // bogus C*alpha/dt terms into the Jacobian. If dt is
220 // large, then the C*alpha/dt is very small, and doesn't
221 // significantly change any Jacobian entries, except for
222 // zero entries.
223 //
224 // This function is only called in the new-DAE case.
225 //
226 // Scope : public
227 // Creator : Eric Keiter, SNL
228 // Creation Date : 03/09/04
229 //-----------------------------------------------------------------------------
231 {
232  N_LAS_Matrix & dQdx = *(ds.dQdxMatrixPtr);
233  N_LAS_Matrix & dFdx = *(ds.dFdxMatrixPtr);
234  N_LAS_Matrix & Jac = *(ds.JMatrixPtr);
235 
236  Jac.linearCombo( 1.0e-20, dQdx, 1.0, dFdx );
237 }
238 
239 //-----------------------------------------------------------------------------
240 // Function : N_TIA_NoTimeIntegration::applyJacobian
241 //
242 // Purpose : Applies the Jacobian operator for the steady state
243 // case.
244 //
245 // Special Notes :
246 // This function is only called in the new-DAE HB (matrix-free) case.
247 //
248 // Scope : public
249 // Creator : Todd Coffey, Ting Mei
250 // Creation Date : 07/29/08
251 //-----------------------------------------------------------------------------
252 void N_TIA_NoTimeIntegration::applyJacobian (const N_LAS_Vector& input, N_LAS_Vector& result)
253 {
254  N_LAS_Vector & dQdxV = *(ds.dQdxVecVectorPtr);
255  N_LAS_Vector & dFdxV = *(ds.dFdxVecVectorPtr);
256  result.linearCombo( 1.0e-20, dQdxV, 1.0, dFdxV );
257 }
258 
259 //-----------------------------------------------------------------------------
260 // Function : N_TIA_NoTimeIntegration::getInitialQnorm
261 // Purpose : Needed by 2-level solves.
262 // Special Notes :
263 // Scope : public
264 // Creator : Eric Keiter, SNL, Parallel Computational Sciences
265 // Creation Date : 03/18/07
266 //-----------------------------------------------------------------------------
268 {
269  tle.q1HistorySum = 0.0;
270 }
271 
272 //-----------------------------------------------------------------------------
273 // Function : N_TIA_NoTimeIntegration::setupTwoLevelError
274 // Purpose : Needed by 2-level solves.
275 // Special Notes :
276 // Scope : public
277 // Creator : Eric Keiter, SNL, Parallel Computational Sciences
278 // Creation Date : 03/15/07
279 //-----------------------------------------------------------------------------
281 {
282  tle.xErrorSum = 0.0;
283  tle.qErrorSum = 0.0;
284  tle.xErrorSum_m1 = 0.0;
285  tle.xErrorSum_m2 = 0.0;
286  tle.innerSize = ds.globalLength ();
287 }
288