Xyce  6.1
N_NLS_NOX_AugmentLinSys_OPStart.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_NLS_NOX_AugmentLinSys_OPStart.C,v $
27 //
28 // Purpose : Algorithm for augmenting the Jacobian for pseudo
29 // transient solves.
30 //
31 // Special Notes :
32 //
33 // Creator : Dave Shirley, PSSI
34 //
35 // Creation Date : 05/08/06
36 //
37 // Revision Information:
38 // ---------------------
39 //
40 // Revision Number: $Revision: 1.35 $
41 //
42 // Revision Date : $Date: 2015/04/08 19:18:29 $
43 //
44 // Current Owner : $Author: tvrusso $
45 //-------------------------------------------------------------------------
46 
47 #include <Xyce_config.h>
48 
49 
50 // ---------- Standard Includes ----------
51 
52 #include <N_UTL_fwd.h>
53 
54 // ---------- Xyce Includes ----------
55 
56 
57 #include "N_LAS_Vector.h"
58 #include "N_LAS_Matrix.h"
59 #include "N_ERH_ErrorMgr.h"
61 
62 namespace Xyce {
63 namespace Nonlinear {
64 namespace N_NLS_NOX {
65 
66 //-----------------------------------------------------------------------------
67 // Function : AugmentLinSysOPStart::AugmentLinSysOPStart
68 // Purpose : constructor
69 // Special Notes :
70 // Scope : public
71 // Creator : Dave Shirley, PSSI
72 // Creation Date :
73 //-----------------------------------------------------------------------------
75  IO::InitialConditionsData::NodeNamePairMap & op_in,
76  const NodeNameMap & allNodes_in, N_PDS_Comm * pdsCommPtr
77  ) : op_ (op_in),
78  allNodes_(allNodes_in),
79  residualPtr_ (0),
80  solutionPtr_ (0),
81  pdsCommPtr_ (pdsCommPtr),
82  rSize_ (0),
83  skipSet (false)
84 {
85 }
86 
87 //-----------------------------------------------------------------------------
88 // Function : AugmentLinSysOPStart::AugmentLinSysOPStart
89 // Purpose : destructor
90 // Special Notes :
91 // Scope : public
92 // Creator : Dave Shirley, PSSI
93 // Creation Date :
94 //-----------------------------------------------------------------------------
96 {
97 }
98 
99 //-----------------------------------------------------------------------------
100 // Function : AugmentLinSysOPStart::augmentResidual
101 // Purpose : augments the residual to support DCOP restart.
102 //
103 // Special Notes : erkeite: This function doesn't seem to have much purpose.
104 // It saves pointers to the residual and solution vectors, so
105 // that it can use them in the augmentJacobian function, below.
106 //
107 // erkeite: The operations that need to be done on the
108 // residual depend on the jacobian structure. In
109 // particular, they depend on what is happening with
110 // voltage sources. Naively applying the 1 on the diagonal
111 // and 0 on the residual will result in singular matrices
112 // if the nodes in question are already being set by a
113 // voltage source.
114 //
115 // Scope : public
116 // Creator : Dave Shirley, PSSI
117 // Creation Date :
118 //-----------------------------------------------------------------------------
120  (const Linear::Vector * solution, Linear::Vector * residual_vector)
121 {
122  residualPtr_ = residual_vector;
123  solutionPtr_ = solution;
124 
125 #ifdef Xyce_PARALLEL_MPI
126  pmap_ = residualPtr_->pmap();
127 #endif
128 
129  return;
130 }
131 
132 //-----------------------------------------------------------------------------
133 // Function : AugmentLinSysOPStart::augmentJacobian
134 // Purpose :
135 // Special Notes :
136 // Scope : public
137 // Creator : Dave Shirley, PSSI
138 // Creation Date :
139 //-----------------------------------------------------------------------------
140 void AugmentLinSysOPStart::augmentJacobian(Linear::Matrix * jacobian)
141 {
142  IO::InitialConditionsData::NodeNamePairMap::iterator op_i;
143  IO::InitialConditionsData::NodeNamePairMap::iterator op_end = op_.end();
144  int i, row, rowLen, global_row, numRows;
145  std::vector<int> col;
146  std::vector<double> val;
147  std::vector<double> resTmp;
148  bool diag;
149  int GID;
150  std::map<int,std::string> rowOut;
151 
152  // erkeite: Determine which rows have diagonal elements already.
153  // A lot of this code is needed to avoid causing the matrix to
154  // be singular, which is easy to do when dealing with rows/cols
155  // associated with Vsrc devices.
156  numRows = jacobian->getLocalNumRows();
157  if (!skipSet)
158  {
159  for (row=0 ; row<numRows ; ++row)
160  {
161 #ifdef Xyce_PARALLEL_MPI
162  global_row = pmap_->localToGlobalIndex(row);
163 #else
164  global_row = row;
165 #endif
166  rowLen = jacobian->getRowLength(global_row);
167  //rowLen = jacobian->getLocalRowLength(row);
168  col.resize(rowLen);
169  val.resize(rowLen);
170  jacobian->getRowCopy(global_row, rowLen, rowLen, &val[0], &col[0]);
171  diag = false;
172  for (i=0 ; i<rowLen ; ++i)
173  {
174  GID = col[i];
175  if (GID == global_row)
176  {
177  if (val[i] != 0)
178  {
179  diag = true;
180  }
181  }
182  }
183  if (!diag)
184  {
185  skipLID.insert(row);
186  skipGID.insert(global_row);
187  }
188  }
189 #ifdef Xyce_PARALLEL_MPI
190  int numG = skipGID.size();
191  int numG_tot;
192  pdsCommPtr_->sumAll(&numG, &numG_tot, 1);
193  int myPos;
194  pdsCommPtr_->scanSum(&numG, &myPos, 1);
195  myPos -= numG;
196  std::vector<int> buf(numG_tot,0);
197  std::vector<int> buf_g(numG_tot,0);
198  std::set<int>::iterator skip_i=skipGID.begin();
199  std::set<int>::iterator skip_end=skipGID.end();
200  for ( ; skip_i != skip_end ; ++skip_i)
201  buf[myPos++] = *skip_i;
202  pdsCommPtr_->sumAll (&buf[0], &buf_g[0], numG_tot);
203  skipGID.clear();
204  for (i=0 ; i < numG_tot ; i++)
205  skipGID.insert(buf_g[i]);
206 #endif
207  skipSet = true;
208  }
209 
210  op_i = op_.begin();
211  std::set<int>::iterator skipLIDEnd = skipLID.end();
212  std::set<int>::iterator skipGIDEnd = skipGID.end();
213  for ( ; op_i != op_end ; ++op_i)
214  {
215  row = (*op_i).second.first;
216  if (skipLID.find(row) == skipLIDEnd)
217  {
218 #ifdef Xyce_PARALLEL_MPI
219  global_row = pmap_->localToGlobalIndex(row);
220 #else
221  global_row = row;
222 #endif
223  rowLen = jacobian->getRowLength(global_row);
224  col.resize(rowLen);
225  val.resize(rowLen);
226  jacobian->getRowCopy(global_row, rowLen, rowLen, &val[0], &col[0]); // from the A-matrix
227  // zero out all residual elements that are not in the skip list.
228  Linear::Vector & residual = (*residualPtr_);
229  residual[row] = 0;
230  for (i=0 ; i<rowLen ; i++)
231  {
232  GID = col[i];
233  // check if this column is in the skip list. It will be in the skip list if
234  // the diagonal element, (GID,GID) doensn't exist.
235  if (skipGID.find(GID) == skipGIDEnd)
236  {
237  if (GID == global_row) // in other words, if this is the diagonal, then set to "1".
238  {
239  val[i] = 1;
240  }
241  else // otherwise, set all non-diagonal matrix entries to zero.
242  {
243  val[i] = 0;
244  }
245  }
246  else // if (GID,GID) does not exist, then residual element is modified
247  {
248  residual[row] += val[i] *
249  (const_cast<Linear::Vector*>(solutionPtr_))->getElementByGlobalIndex(GID);
250  }
251  }
252  jacobian->shirleyPutRow(global_row, rowLen, &val[0], &col[0]);
253  }
254  }
255 
256  return;
257 }
258 
259 }}}
Pure virtual class to augment a linear system.
void augmentResidual(const Xyce::Linear::Vector *solution, Xyce::Linear::Vector *residual_vector)
Augments the Residual.
void augmentJacobian(Xyce::Linear::Matrix *jacobian)
Augments the Jacobian.
Xyce::IO::InitialConditionsData::NodeNamePairMap & op_
map of specified variables
AugmentLinSysOPStart(Xyce::IO::InitialConditionsData::NodeNamePairMap &, const Xyce::NodeNameMap &, N_PDS_Comm *)