Xyce  6.1
N_DEV_DiodePDE_DAE.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 //-------------------------------------------------------------------------
27 // Filename : $RCSfile: N_DEV_DiodePDE_DAE.C,v $
28 //
29 // Purpose : One dimensional PDE device, new-DAE functions.
30 //
31 // Special Notes :
32 //
33 // Creator : Eric R. Keiter, SNL, Parallel Computational Sciences
34 //
35 // Creation Date : 05/13/05
36 //
37 // Revision Information:
38 // ---------------------
39 //
40 // Revision Number: $Revision: 1.39.2.1 $
41 //
42 // Revision Date : $Date: 2015/04/02 18:20:12 $
43 //
44 // Current Owner : $Author: tvrusso $
45 //-------------------------------------------------------------------------
46 
47 #include <Xyce_config.h>
48 
49 
50 // ---------- Standard Includes ----------
51 #include <iostream>
52 #include <N_UTL_Math.h>
53 #include <cstdio>
54 
55 // ---------- Xyce Includes ----------
56 #include <N_DEV_DiodePDE.h>
57 #include <N_DEV_ExternData.h>
58 #include <N_DEV_SolverState.h>
59 #include <N_DEV_DeviceOptions.h>
60 #include <N_DEV_MatrixLoadData.h>
61 
62 #include <N_LAS_Vector.h>
63 #include <N_LAS_Matrix.h>
64 #include <N_LAS_System.h>
65 #include <N_LAS_Builder.h>
66 
67 // default number of mesh points:
68 #define NUM_MESH_POINTS 11
69 
70 // default maximum number of nonzero entries in a matrix row
71 #define MAX_COLS_PER_ROW 10
72 
73 namespace Xyce {
74 namespace Device {
75 namespace DiodePDE {
76 
77 //-----------------------------------------------------------------------------
78 // Function : Instance::loadDAEFVector
79 // Purpose :
80 // Special Notes :
81 // Scope : public
82 // Creator : Eric Keiter, SNL, Parallel Computational Sciences
83 // Creation Date : 5/18/05
84 //-----------------------------------------------------------------------------
86 {
87  bool bsuccess = true;
88  bool bs1;
89 
90  if (DEBUG_DEVICE && isActive(Diag::DEVICE_PARAMETERS) && getSolverState().debugTimeFlag)
91  {
92  Xyce::dout() << "loadDAEFVector: doubleDCOPStep="<<getSolverState().doubleDCOPStep;
93  if (getSolverState().dcopFlag) Xyce::dout() << " DCOP load" << std::endl;
94  else Xyce::dout() << " Transient load" << std::endl;
95  }
96 
97  if ((getSolverState().dcopFlag) && getSolverState().doubleDCOPStep == 0)
98  {
99  equationSet = 0;
100  bs1 = loadDAEFNonlinPoisson ();
101  }
102  else
103  {
104  equationSet = 1;
105 
106  if (getSolverState().twoLevelNewtonCouplingMode==Nonlinear::INNER_PROBLEM ||
107  getSolverState().twoLevelNewtonCouplingMode==Nonlinear::FULL_PROBLEM)
108  {
109  bs1 = loadDAEFDDFormulation ();
110  }
111  else if (getSolverState().twoLevelNewtonCouplingMode==Nonlinear::OUTER_PROBLEM)
112  {
114  }
115  else
116  {
117  std::string msg = "Instance::loadDAEFVector."
118  "Invalid coupling Mode";
119  N_ERH_ErrorMgr::report ( N_ERH_ErrorMgr::DEV_FATAL,msg);
120  }
121  }
122 
123  return (bsuccess && bs1);
124 }
125 
126 //-----------------------------------------------------------------------------
127 // Function : Instance::loadDAEFNonlinPoisson
128 // Purpose : Loads the F-vector the nonlinear poisson calculation.
129 //
130 // Special Notes : This should be identical to the original loadRHS function,
131 // for the nonlinear poisson. All of that goes into the
132 // FVector, but with the opposite sign.
133 //
134 // Scope : public
135 // Creator : Eric Keiter, SNL, Parallel Computational Sciences
136 // Creation Date : 5/18/05
137 //-----------------------------------------------------------------------------
139 {
140  double * fvec = extData.daeFVectorRawPtr;
141 
142  bool bsuccess = loadVecNLPoisson ( fvec );
143 
144  return bsuccess;
145 }
146 
147 //-----------------------------------------------------------------------------
148 // Function : Instance::loadDAEFDDFormulation
149 // Purpose : This function should be called from the loadDAEFVector
150 // function when solving the drift-diffusion equations.
151 // Special Notes :
152 //
153 // Scope : private
154 // Creator : Eric Keiter, SNL, Parallel Computational Sciences
155 // Creation Date : 5/18/05
156 //-----------------------------------------------------------------------------
158 {
159  double * fvec = extData.daeFVectorRawPtr;
160 
161  bool bsuccess = loadVecDDForm ( fvec );
162 
163  return bsuccess;
164 }
165 
166 //-----------------------------------------------------------------------------
167 // Function : Instance::loadDAEFExtractedConductance
168 // Purpose :
169 // Special Notes :
170 // Scope : private
171 // Creator : Eric Keiter, SNL, Parallel Computational Sciences
172 // Creation Date : 05/18/05
173 //-----------------------------------------------------------------------------
175 {
176  double coef;
177 
178  Linear::Vector & fVec = *(extData.daeFVectorPtr) ;
179 
180  // KCL equations for the various connecting terminals:
181  int iRow = 0;
182  int iBC;
183  for (iBC=0;iBC<bcVec.size();++iBC, ++iRow)
184  {
185  coef = bcVec[iBC].currentSum;
186 
187  if (iBC==1) coef *= -1.0;
188 
189  double voltLimFac = 0.0;
191  {
192  int iCol;
193  for (iCol=0; iCol < numElectrodes ; ++iCol)
194  {
195  if (bcVec[iCol].gid == -1) continue;
196 
197  double vdiff = bcVec[iCol].Vckt_final - bcVec[iCol].Vckt_orig;
198 
199  vdiff *= scalingVars.V0;
200 
201  voltLimFac += vdiff * condVec[iRow][iCol];
202  }
203  }
204 
205  fVec[bcVec[iBC].lid] += -coef + voltLimFac;
206  }
207 
208  // Load in the zeros...
209  for (int i=0;i<NX;++i)
210  {
211  fVec[ li_Vrowarray[i] ] = 0.0;
212  fVec[ li_Nrowarray[i] ] = 0.0;
213  fVec[ li_Prowarray[i] ] = 0.0;
214  }
215 
216  return true;
217 }
218 
219 //-----------------------------------------------------------------------------
220 // Function : Instance::loadDAEQVector
221 // Purpose :
222 // Special Notes :
223 // Scope : public
224 // Creator : Eric Keiter, SNL, Parallel Computational Sciences
225 // Creation Date : 5/18/05
226 //-----------------------------------------------------------------------------
228 {
229  bool bsuccess = true;
230  bool bs1;
231 
232  Linear::Vector & qvec = *(extData.daeQVectorPtr);
233 
234  if ((getSolverState().dcopFlag) && getSolverState().doubleDCOPStep == 0)
235  {
236  equationSet = 0;
237 
238  bs1 = true; // no-op
239  }
240  else
241  {
242  equationSet = 1;
243 
246  {
247  bs1 = loadDAEQDDFormulation ();
248  }
249  else if (getSolverState().twoLevelNewtonCouplingMode==Nonlinear::OUTER_PROBLEM)
250  {
252  }
253  else
254  {
255  std::string msg = "Instance::loadDAEQVector."
256  "Invalid coupling Mode";
257  N_ERH_ErrorMgr::report ( N_ERH_ErrorMgr::DEV_FATAL,msg);
258  }
259  }
260 
261  return (bsuccess && bs1);
262 }
263 
264 
265 //-----------------------------------------------------------------------------
266 // Function : Instance::loadDAEQDDFormulation
267 // Purpose :
268 // Special Notes :
269 // Scope : private
270 // Creator : Eric Keiter, SNL, Parallel Computational Sciences
271 // Creation Date : 05/18/05
272 //-----------------------------------------------------------------------------
274 {
275  Linear::Vector & qvec = *(extData.daeQVectorPtr);
276  int i;
277 
278  // mesh points for the PDE problem:
279  // only do the interior mesh points - boundaries have no dndt terms.
280  for (i=1;i<NX-1;++i)
281  {
282  qvec[ li_Nrowarray[i] ] = -nnVec[i]*scalingVars.t0;
283  qvec[ li_Prowarray[i] ] = -npVec[i]*scalingVars.t0;
284  } // ip_iter, row loop...
285 
286  return true;
287 }
288 
289 //-----------------------------------------------------------------------------
290 // Function : Instance::loadDAEQExtractedConductance
291 // Purpose :
292 // Special Notes :
293 // Scope : private
294 // Creator : Eric Keiter, SNL, Parallel Computational Sciences
295 // Creation Date : 05/18/05
296 //-----------------------------------------------------------------------------
298 {
299  return true;
300 }
301 
302 //-----------------------------------------------------------------------------
303 // Function : Instance::loadDAEdFdx
304 // Purpose : This function performs an analytic Jacobian matrix load for
305 // the diode-pde class, for the case of solving a nonlinear
306 // poisson equation.
307 // Special Notes :
308 //
309 // Scope : private
310 // Creator : Eric Keiter, SNL, Parallel Computational Sciences
311 // Creation Date : 05/18/05
312 //-----------------------------------------------------------------------------
314 {
315  bool bsuccess;
316 
317  if ((getSolverState().dcopFlag) && getSolverState().doubleDCOPStep == 0)
318  {
319  bsuccess = loadDAEdFdxNonlinPoisson ();
320  }
321  else
322  {
323  if (getSolverState().twoLevelNewtonCouplingMode==Nonlinear::INNER_PROBLEM ||
324  getSolverState().twoLevelNewtonCouplingMode==Nonlinear::FULL_PROBLEM)
325  {
326  bsuccess = loadDAEdFdxDDFormulation ();
327  }
328  else if (getSolverState().twoLevelNewtonCouplingMode==Nonlinear::OUTER_PROBLEM)
329  {
330  bsuccess = loadDAEdFdxExtractedConductance ();
331  }
332  else
333  {
334  Report::DevelFatal().in("Instance::loadDAEdFdx") << "Invalid coupling Mode" << numElectrodes;
335  }
336  }
337 
338  return bsuccess;
339 }
340 
341 //-----------------------------------------------------------------------------
342 // Function : Instance::loadDAEdFdxNonlinPoisson
343 // Purpose : This function performs an analytic Jacobian matrix load for
344 // the diode-pde class, for the case of solving a nonlinear
345 // poisson equation.
346 // Special Notes :
347 //
348 // Scope : private
349 // Creator : Eric Keiter, SNL, Parallel Computational Sciences
350 // Creation Date : 05/18/05
351 //-----------------------------------------------------------------------------
353 {
355 }
356 
357 //-----------------------------------------------------------------------------
358 // Function : Instance::loadDAEdFdxDDFormulation
359 // Purpose :
360 // Special Notes :
361 // Scope : private
362 // Creator : Eric Keiter, SNL, Parallel Computational Sciences
363 // Creation Date : 05/18/05
364 //-----------------------------------------------------------------------------
366 {
367  return loadMatDDForm ( *(extData.dFdxMatrixPtr));
368 }
369 
370 //-----------------------------------------------------------------------------
371 // Function : Instance::loadDAEdFdxExtractedConductance
372 // Purpose :
373 // Special Notes :
374 // Scope : private
375 // Creator : Eric Keiter, SNL, Parallel Computational Sciences
376 // Creation Date : 05/18/05
377 //-----------------------------------------------------------------------------
379 {
380  bool bsuccess = true;
381  bool bs1 = true;
382 
383  Linear::Matrix & mat = *(extData.dFdxMatrixPtr);
384 
385  int iRow, iCol;
386  int colsLID[10];
387  double valsLID[10];
388  int i;
389  for (i=0;i<10;++i)
390  {
391  colsLID[i] = -1;
392  valsLID[i] = 0.0;
393  }
394 
395  // first put 1's on the diagonals of all the mesh-rows:
396  for (i=0;i<NX;++i)
397  {
398  valsLID[0] = 1.0;
399 
400  iRow = li_Vrowarray[i];
401  if (iRow != -1)
402  {
403  colsLID[0] = iRow;
404  bs1 = mat.putLocalRow(iRow,1, valsLID, colsLID);
405  bsuccess = bsuccess && bs1;
406  }
407 
408  iRow = li_Nrowarray[i];
409  if (iRow != -1)
410  {
411  colsLID[0] = iRow;
412  bs1 = mat.putLocalRow(iRow,1, valsLID, colsLID);
413  bsuccess = bsuccess && bs1;
414  }
415 
416  iRow = li_Prowarray[i];
417  if (iRow != -1)
418  {
419  colsLID[0] = iRow;
420  bs1 = mat.putLocalRow(iRow,1, valsLID, colsLID);
421  bsuccess = bsuccess && bs1;
422  }
423  }
424 
425  // now load the equivalent conductances.
426  for (iRow=0; iRow < numElectrodes ; ++iRow)
427  {
428 
429  int iRowLID = bcVec[iRow].gid;
430 
431  if (iRowLID == -1) continue;
432 
433  int count = 0;
434 
435  for (iCol=0; iCol < numElectrodes ; ++iCol)
436  {
437  if (bcVec[iCol].gid == -1) continue;
438 
439  colsLID[count] = bcVec[iCol].gid;
440  valsLID[count] = condVec[iRow][iCol];
441  ++count;
442  }
443 
444  bs1 = mat.sumIntoLocalRow (iRowLID, count, valsLID, colsLID);
445  bsuccess = bsuccess && bs1;
446  }
447 
448  if (DEBUG_DEVICE && isActive(Diag::DEVICE_PARAMETERS))
449  {
450  int iE1, iE2;
451  char tmpString[128]; for (int i=0;i<128;++i) tmpString[i] = 0;
452 
453  Xyce::dout() << std::endl;
454  sprintf(tmpString,"ConArray:\t "); Xyce::dout() << std::string(tmpString);
455  for (iE2 = 0; iE2 < numElectrodes; ++ iE2)
456  {
457  sprintf(tmpString,"\t%14s",bcVec[iE2].eName.c_str()); Xyce::dout() << std::string(tmpString);
458  }
459  Xyce::dout() << std::endl;
460 
461  for (iE1 = 0; iE1 < numElectrodes; ++iE1)
462  {
463  sprintf(tmpString,"ConArray:\t%14s",bcVec[iE1].eName.c_str()); Xyce::dout() << std::string(tmpString);
464  for (iE2 = 0; iE2 < numElectrodes; ++ iE2)
465  {
466  sprintf(tmpString,"\t%14.4e",condVec[iE1][iE2]); Xyce::dout() << std::string(tmpString);
467  }
468  Xyce::dout() << std::endl;
469  }
470  Xyce::dout() << std::endl;
471  }
472 
473  return bsuccess;
474 }
475 
476 //-----------------------------------------------------------------------------
477 // Function : Instance::loadDAEdQdx
478 // Purpose :
479 // Special Notes :
480 // Scope : private
481 // Creator : Eric Keiter, SNL, Parallel Computational Sciences
482 // Creation Date : 05/18/05
483 //-----------------------------------------------------------------------------
485 {
486  bool bsuccess;
487  Linear::Matrix & dQdxMat = *(extData.dQdxMatrixPtr);
488 
489  if ((getSolverState().dcopFlag) && getSolverState().doubleDCOPStep == 0)
490  {
491  bsuccess = true; // no-op
492  }
493  else
494  {
497  {
498  bsuccess = loadDAEdQdxDDFormulation ();
499  }
500  else if (getSolverState().twoLevelNewtonCouplingMode==Nonlinear::OUTER_PROBLEM)
501  {
502  bsuccess = loadDAEdQdxExtractedConductance ();
503  }
504  else
505  {
506  Report::DevelFatal().in("Instance::loadDAEdQdx") << "Invalid coupling Mode " << numElectrodes;
507  }
508  }
509 
510  return bsuccess;
511 }
512 
513 //-----------------------------------------------------------------------------
514 // Function : Instance::loadDAEdQdxDDFormulation
515 // Purpose :
516 // Special Notes :
517 // Scope : private
518 // Creator : Eric Keiter, SNL, Parallel Computational Sciences
519 // Creation Date : 05/18/05
520 //-----------------------------------------------------------------------------
522 {
523  bool bsuccess = true;
524 
525  Linear::Matrix & dQdxMat = *(extData.dQdxMatrixPtr);
526 
527  // load the rows associated with the PDE mesh:
528  // use direct matrix access:
529  for (int i=1;i<NX-1;++i)
530  {
531  int li_Nrow = li_Nrowarray[i];
532  int li_Prow = li_Prowarray[i];
533 
534  // Note::these terms need to be the SAME sign as the terms
535  // in loadDAEQVector!
536 
537  // electron continuity row:
538  // derivative w.r.t. nnVec[i ]:
539  dQdxMat[li_Nrow][li_Ncolarray[i][1]] = -scalingVars.t0;
540 
541  // hole continuity row:
542  // derivative w.r.t. npVec[i ]:
543  dQdxMat[li_Prow][li_Pcolarray[i][1]] = -scalingVars.t0;
544  }
545 
546  return bsuccess;
547 }
548 
549 //-----------------------------------------------------------------------------
550 // Function : Instance::loadDAEdQdxExtractedConductance
551 // Purpose :
552 // Special Notes :
553 // Scope : private
554 // Creator : Eric Keiter, SNL, Parallel Computational Sciences
555 // Creation Date : 05/18/05
556 //-----------------------------------------------------------------------------
558 {
559  bool bsuccess = true;
560  return bsuccess;
561 }
562 
563 } // namespace DiodePDE
564 } // namespace Device
565 } // namespace Xyce
Linear::Vector * daeQVectorPtr
Pure virtual class to augment a linear system.
bool loadMatDDForm(Linear::Matrix &mat)
Nonlinear::TwoLevelNewtonMode twoLevelNewtonCouplingMode
bool loadMatNLPoisson(Linear::Matrix &mat)
const DeviceOptions & getDeviceOptions() const
Linear::Matrix * dFdxMatrixPtr
std::vector< std::vector< int > > li_Pcolarray
const SolverState & getSolverState() const
std::vector< std::vector< double > > condVec
Linear::Vector * daeFVectorPtr
std::vector< std::vector< int > > li_Ncolarray
Linear::Matrix * dQdxMatrixPtr