Update ChangeLog
[gss-tcad.git] / src / include / em_fem.h
blobf0356e2bfeb7af3b40cec2375ef0e716b38ffe26
1 /*****************************************************************************/
2 /* 8888888 88888888 88888888 */
3 /* 8 8 8 */
4 /* 8 8 8 */
5 /* 8 88888888 88888888 */
6 /* 8 8888 8 8 */
7 /* 8 8 8 8 */
8 /* 888888 888888888 888888888 */
9 /* */
10 /* A Two-Dimensional General Purpose Semiconductor Simulator. */
11 /* */
12 /* GSS 0.4x */
13 /* Last update: Jan 19, 2006 */
14 /* */
15 /* Gong Ding */
16 /* gdiso@ustc.edu */
17 /* NINT, No.69 P.O.Box, Xi'an City, China */
18 /* */
19 /*****************************************************************************/
21 #ifndef _em_fem_h_
22 #define _em_fem_h_
24 #include "mathfunc.h"
25 #include "petscksp.h"
26 #include "bsolver.h"
27 #include "mesh.h"
28 #include "cmdbuf.h"
29 #include <vector>
30 using namespace std;
32 class EM_FEM_Solver : virtual public BSolver
34 private:
35 int N; // the scale if matrix
36 int total_nodes; // node number of all zones
37 PetscScalar v; //insert variable
38 PetscInt its; // iteration number
39 PetscInt *rows,count; //variables for add dirichlet boundary condition
40 PetscReal norm;
42 PetscErrorCode ierr;
43 KSPConvergedReason reason;
44 KSP solver;
45 PC prec;
46 //vector<int> zofs; // the offset value of each zone in solution vector
48 int eleNode[3]; //store nodes' global number of each element
49 PetscScalar eleA[3]; //assistant variables
50 PetscScalar eleB[3];
51 PetscScalar eleC[3];
52 PetscScalar eleMatrix[3][3]; //element matrix
53 PetscScalar eleV[3]; //element right vector
55 PetscScalar RealAlphaX; //only used for PML zone
56 PetscScalar RealAlphaY;
57 PetscScalar RealBeta;
58 PetscScalar ImageAlphaX;
59 PetscScalar ImageAlphaY;
60 PetscScalar ImageBeta;
61 PetscScalar Electri0; //only used for PML zone
62 PetscScalar Magne0;
63 PetscScalar TransVar;
64 PetscScalar CurrentX;
65 PetscScalar CurrentY;
67 PetscScalar RealPermi; //Real part of permitivity
68 PetscScalar ImagePermi; //Image part ```````
69 PetscScalar InitPhi;
70 PetscScalar PhiDiff;
71 PetscScalar SourceAlpha;
73 public:
74 Vec x,b; // solution and right vectors
75 //Vectors to save Ex,Ey,Hx,Hy in complex form
76 vector<PetscScalar> TMEZ,TMHX,TMHY,TEHZ,TEEX,TEEY;
77 PetscScalar *TM; //Needle of TM mode for output the result
78 PetscScalar *TE;
79 Mat A; // matrix
80 PetscScalar XLeft; //minimun value of x in the whole model
81 PetscScalar XRight;
82 PetscScalar YTop;
83 PetscScalar YBottom;
84 PetscScalar VacXMax; //Find max X in vacuum zone
85 PetscScalar VacXMin;
86 PetscScalar VacYMax;
87 PetscScalar VacYMin; //These 4 variable is used to find 4 corners of vacuum zone
88 PetscScalar lamda; //wave length of incidence light
89 PetscScalar power; //intensity of incidence light
90 PetscScalar WeightOfTE; //weight of TE mode in incidence light
91 PetscScalar WeightOfTM;
92 PetscScalar OptAngle; //incident angle of light --degree
93 PetscScalar OptTheta; //incident angle of light --radian
94 PetscScalar eta; //quantum efficiency of semiconductor material
95 PetscScalar WaveVector0; //wave vector in vacuum
96 PetscScalar SourceTM; //Source term for TM mode
97 PetscScalar SourceTE;
98 PetscScalar DeltaX; //variable to scale grid size
99 int m; //constant used in PML
100 PetscScalar R;
101 PetscScalar EpsInVac; //Variables used when changing power to field
102 PetscScalar MuInVac;
103 bool append;
104 public:
105 int init_solver(SolveDefine &sv);
106 int do_solve(SolveDefine &sv);
107 int TM_solution_update(SolveDefine &sv);
108 int TE_solution_update(SolveDefine &sv);
109 int field_calculation(SolveDefine &sv);
110 int generation_cal(SolveDefine &sv);
111 int destroy_solver(SolveDefine &sv);
112 EM_FEM_Solver();
115 #endif