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