Xyce
6.1
Main Page
Related Pages
Modules
Namespaces
Classes
Files
File List
File Members
All
Classes
Namespaces
Files
Functions
Variables
Typedefs
Enumerations
Enumerator
Friends
Macros
Groups
Pages
N_TIA_NoTimeIntegration.h
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.h,v $
27
//
28
// Purpose : This file defines the classes for the "no time integration"
29
// method.
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.32 $
41
//
42
// Revision Date : $Date: 2014/02/24 23:49:26 $
43
//
44
// Current Owner : $Author: tvrusso $
45
//-----------------------------------------------------------------------------
46
47
#ifndef Xyce_N_TIA_NO_TIME_INTEGRATION_H
48
#define Xyce_N_TIA_NO_TIME_INTEGRATION_H
49
50
// ---------- Standard Includes ----------
51
#ifdef HAVE_CMATH
52
#include <cmath>
53
#else
54
#include <math.h>
55
#endif
56
57
// ---------- Xyce Includes ----------
58
#include <
N_TIA_DataStore.h
>
59
#include <
N_TIA_StepErrorControl.h
>
60
#include <
N_TIA_TimeIntegrationMethods.h
>
61
62
#define MATRIX_FAILSAFE 1
63
#define NUM_LIMIT 1.0e-20
64
65
//-----------------------------------------------------------------------------
66
// Class : N_TIA_NoTimeIntegration
67
// Purpose : Class objects for use during DC analysis (applies only when
68
// all time derivatives are set to 0) (derived from
69
// N_TIA_TimeIntegrationMethod)
70
// Special Notes :
71
// Creator : Buddy Watts, SNL
72
// Creation Date : 6/01/00
73
//-----------------------------------------------------------------------------
74
class
N_TIA_NoTimeIntegration
:
public
N_TIA_TimeIntegrationMethod
75
{
76
public
:
77
78
// Destructor
79
~N_TIA_NoTimeIntegration
();
80
81
// Predict solution at next time point (No Integration).
82
virtual
void
obtainPredictor
() {
ds
.
usePreviousSolAsPredictor
(); }
83
84
// Evaluate the predictor derivative formula (No Integration).
85
virtual
void
obtainPredictorDeriv
() { }
86
87
// Evaluate the corrector derivative formula (No Integration).
88
virtual
void
obtainCorrectorDeriv
();
89
90
// Compute an estimate of the error in the integration step (No Integration).
91
virtual
void
updateDerivsBlock
(
const
std::list< index_pair > & solGIDList,
92
const
std::list< index_pair > & staGIDList)
93
{ }
94
95
#ifdef MATRIX_FAILSAFE
96
virtual
double
partialTimeDeriv
() {
return
NUM_LIMIT
; }
97
98
#else
99
virtual
double
partialTimeDeriv
() {
return
0.0; }
100
101
#endif
102
103
// Compute an estimate of the error in the integration step (No Integration).
104
virtual
double
computeErrorEstimate
() {
return
0.0; }
105
106
// Interpolate solution approximation at prescribed time point (No
107
// Integration).
108
virtual
bool
interpolateSolution
(
double
timepoint,
109
N_LAS_Vector * tmpSolVectorPtr, std::vector<N_LAS_Vector*> & historyVec )
110
{
return
false
;};
111
112
// Computes the step adjustment (No Integration).
113
virtual
double
computeExpoStepAdjust
(
double
stepadjust) {
return
0.0; }
114
115
// Gets the time-integration order (No Integration).
116
virtual
int
getOrder
() {
return
0; }
117
virtual
int
getUsedOrder
() {
return
0; }
118
119
virtual
void
getInitialQnorm
(
N_TIA_TwoLevelError
& tle);
120
virtual
void
setupTwoLevelError
(
N_TIA_TwoLevelError
& tle);
121
122
// Time-integration factory.
123
static
N_TIA_TimeIntegrationMethod
*
factory
(
N_TIA_TIAParams
& tiaP,
124
N_TIA_StepErrorControl
& secTmp ,
125
N_TIA_DataStore
& dsTmp);
126
127
// 03/08/04 erkeite: New functions necessary new-DAE:
128
// Evaluate residual for nonlinear solver
129
void
obtainResidual
();
130
131
// Evaluate Jacobian for nonlinear solver
132
void
obtainJacobian
();
133
134
// Apply Jacobian for nonlinear solver
135
void
applyJacobian
(
const
N_LAS_Vector& input, N_LAS_Vector& result);
136
137
private
:
138
// Copied from N_TIA_BackwardDifferentiation15.h
139
double
alphas
;
// $\alpha_s$ fixed-leading coefficient of this BDF method
140
141
// Constructor (private).
142
N_TIA_NoTimeIntegration
(
N_TIA_TIAParams
& tiaP,
143
N_TIA_StepErrorControl
& secTmp,
144
N_TIA_DataStore
& dsTmp);
145
};
146
147
//-----------------------------------------------------------------------------
148
// Function : N_TIA_NoTimeIntegration::factory
149
// Purpose :
150
// Special Notes :
151
// Scope : public
152
// Creator : Eric Keiter, SNL, Parallel Computational Sciences
153
// Creation Date : 6/08/01
154
//-----------------------------------------------------------------------------
155
inline
N_TIA_TimeIntegrationMethod
*
N_TIA_NoTimeIntegration::factory
156
(
N_TIA_TIAParams
& tiaP,
157
N_TIA_StepErrorControl
& secTmp ,
158
N_TIA_DataStore
& dsTmp)
159
{
160
N_TIA_TimeIntegrationMethod
* timPtr =
161
new
N_TIA_NoTimeIntegration
(tiaP,secTmp,dsTmp);
162
return
timPtr;
163
}
164
165
#endif // Xyce_N_TIA_NO_TIME_INTEGRATION_H
src
TimeIntegrationPKG
include
N_TIA_NoTimeIntegration.h
Generated on Mon Mar 24 2014 10:54:40 for Xyce by
1.8.3.1