1 /*****************************************************************************/
2 /* 8888888 88888888 88888888 */
5 /* 8 88888888 88888888 */
8 /* 888888 888888888 888888888 */
10 /* A Two-Dimensional General Purpose Semiconductor Simulator. */
13 /* Last update: Nov 23, 2005 */
17 /* NINT, No.69 P.O.Box, Xi'an City, China */
19 /*****************************************************************************/
23 void ISZone::F1Q_ddm_inner(int i
,PetscScalar
*x
,PetscScalar
*f
, ODE_Formula
&ODE_F
, vector
<int> & zofs
)
25 const VoronoiCell
* pcell
= &pzone
->davcell
[i
];
26 PetscScalar Vi
= x
[zofs
[zone_index
]+i
];
27 PetscScalar div_grad_P
= 0;
28 for(int j
=0;j
<pcell
->nb_num
;j
++)
30 int nb
= pcell
->nb_array
[j
];
31 PetscScalar Vr
= x
[zofs
[zone_index
]+nb
]; //potential of nb node
32 div_grad_P
+= pcell
->elen
[j
]/pcell
->ilen
[j
]*(Vr
-Vi
)/pcell
->area
;
34 f
[zofs
[zone_index
]+i
] = div_grad_P
;
38 void ISZone::F1Q_ddm_gatebc(int i
,PetscScalar
*x
,PetscScalar
*f
, ODE_Formula
&ODE_F
, vector
<int> & zofs
,DABC
&bc
)
41 int size
= pzone
->davcell
.size();
42 const VoronoiCell
* pcell
= &pzone
->davcell
[i
];
43 GateBC
*pbc
= dynamic_cast<GateBC
* >(bc
.Get_pointer(pcell
->bc_index
-1));
45 for(int j
=0;j
<electrode
.size();j
++)
47 if(electrode
[j
]==pcell
->bc_index
-1) {gate_equ
=j
;break;}
49 f
[zofs
[zone_index
]+i
] = x
[zofs
[zone_index
]+i
] - x
[zofs
[zone_index
]+equ_num
*size
+gate_equ
] + pbc
->WorkFunction
;
53 void ISZone::F1Q_ddm_semiconductor_insulator_interface(int i
,PetscScalar
*x
,PetscScalar
*f
, ODE_Formula
&ODE_F
, vector
<int> & zofs
, DABC
&bc
,SMCZone
*pz
, int n
)
55 f
[zofs
[zone_index
]+i
] = x
[zofs
[zone_index
]+i
] - x
[zofs
[pz
->zone_index
]+5*n
+0];
58 void ISZone::F1Q_ddm_electrode_insulator_interface(int i
,PetscScalar
*x
,PetscScalar
*f
, ODE_Formula
&ODE_F
, vector
<int> & zofs
, DABC
&bc
,ElZone
*pz
, int n
)
60 f
[zofs
[zone_index
]+i
] = x
[zofs
[zone_index
]+i
] - x
[zofs
[pz
->zone_index
]+n
];
63 void ISZone::F1Q_ddm_insulator_insulator_interface(int i
,PetscScalar
*x
,PetscScalar
*f
, ODE_Formula
&ODE_F
, vector
<int> & zofs
, DABC
&bc
,ISZone
*pz
, int n
)
65 if(zone_index
> pz
->pzone
->zone_index
)
67 f
[zofs
[zone_index
]+i
] = x
[zofs
[zone_index
]+i
] - x
[zofs
[pz
->zone_index
]+n
];
70 const VoronoiCell
* pcell
= &pzone
->davcell
[i
];
71 PetscScalar Vi
= x
[zofs
[zone_index
]+i
];
72 PetscScalar grad_P
= 0;
73 for(int j
=0;j
<pcell
->nb_num
;j
++)
75 int nb
= pcell
->nb_array
[j
];
76 PetscScalar Vr
= x
[zofs
[zone_index
]+nb
]; //potential of nb node
77 grad_P
+= aux
[i
].eps
*pcell
->elen
[j
]/pcell
->ilen
[j
]*(Vr
-Vi
);
79 const VoronoiCell
* ncell
= pz
->pzone
->davcell
.GetPointer(n
);
80 for(int j
=0;j
<ncell
->nb_num
;j
++)
82 int nb
= ncell
->nb_array
[j
];
83 PetscScalar Vr
= x
[zofs
[pz
->zone_index
]+nb
]; //potential of nb node
84 grad_P
+= pz
->aux
[n
].eps
*ncell
->elen
[j
]/ncell
->ilen
[j
]*(Vr
-Vi
);
86 f
[zofs
[zone_index
]+i
] = grad_P
;
89 void ISZone::F1Q_ddm_chargebc(int i
,PetscScalar
*x
,PetscScalar
*f
, ODE_Formula
&ODE_F
, vector
<int> & zofs
,DABC
&bc
)
92 int size
= pzone
->davcell
.size();
93 const VoronoiCell
* pcell
= &pzone
->davcell
[i
];
94 ChargedContactBC
*pbc
= dynamic_cast<ChargedContactBC
* >(bc
.Get_pointer(pcell
->bc_index
-1));
96 for(int j
=0;j
<electrode
.size();j
++)
98 if(electrode
[j
]==pcell
->bc_index
-1) {charge_equ
=j
;break;}
100 f
[zofs
[zone_index
]+i
] = x
[zofs
[zone_index
]+i
] - x
[zofs
[zone_index
]+equ_num
*size
+charge_equ
];
104 void ISZone::F1Q_gate_electrode(int i
,PetscScalar
*x
,PetscScalar
*f
, ODE_Formula
&ODE_F
, vector
<int> & zofs
, DABC
&bc
, PetscScalar DeviceDepth
)
107 int size
= pzone
->davcell
.size();
108 int bc_index
= electrode
[i
];
109 GateBC
*pbc
= dynamic_cast <GateBC
* > (bc
.Get_pointer(bc_index
));
110 PetscScalar current
=0;
111 for(int j
=0;j
<bc
[bc_index
].psegment
->node_array
.size();j
++)
113 int node
=bc
[bc_index
].psegment
->node_array
[j
];
114 const VoronoiCell
* pcell
= pzone
->davcell
.GetPointer(node
);
115 PetscScalar Vi
= x
[zofs
[zone_index
]+node
]; //potential of node i
117 for(int k
=0;k
<pcell
->nb_num
;k
++)
119 int nb
= pcell
->nb_array
[k
];
120 PetscScalar Vj
= x
[zofs
[zone_index
]+nb
]; //potential of nb node
121 //displacement current
122 current
+= DeviceDepth
*pcell
->elen
[k
]*aux
[node
].eps
*((Vi
-Vj
)-(fs
[node
].P
-fs
[nb
].P
))/pcell
->ilen
[k
]/ODE_F
.dt
;
125 if(pbc
->electrode_type
==VoltageBC
)
127 PetscScalar V
= pbc
->Vapp
;
128 PetscScalar R
= pbc
->R
;
129 PetscScalar C
= pbc
->C
;
130 PetscScalar L
= pbc
->L
;
131 PetscScalar In
= pbc
->current
;
132 PetscScalar Icn
= pbc
->cap_current
;
133 PetscScalar Pn
= pbc
->potential
;
134 PetscScalar P
= x
[zofs
[zone_index
]+equ_num
*size
+i
];
135 f
[zofs
[zone_index
]+equ_num
*size
+i
]=(L
/ODE_F
.dt
+R
)*current
-V
+(1+(L
/ODE_F
.dt
+R
)*C
/ODE_F
.dt
)*P
136 -(L
/ODE_F
.dt
+R
)*C
/ODE_F
.dt
*Pn
-L
/ODE_F
.dt
*(In
+Icn
);
137 pbc
->Set_Current_new(current
);
139 pbc
->Set_Potential_new(x
[zofs
[zone_index
]+equ_num
*size
+i
]);
143 void ISZone::F1Q_charge_electrode(int i
,PetscScalar
*x
,PetscScalar
*f
, ODE_Formula
&ODE_F
, vector
<int> & zofs
, DABC
&bc
)
146 int size
= pzone
->davcell
.size();
147 int bc_index
= electrode
[i
];
148 ChargedContactBC
*pbc
= dynamic_cast <ChargedContactBC
* > (bc
.Get_pointer(bc_index
));
149 PetscScalar grad_P
=0;
150 for(int j
=0;j
<bc
[bc_index
].psegment
->node_array
.size();j
++)
152 int node
=bc
[bc_index
].psegment
->node_array
[j
];
153 const VoronoiCell
* pcell
= pzone
->davcell
.GetPointer(node
);
154 PetscScalar Vi
= x
[zofs
[zone_index
]+node
]; //potential of node i
155 PetscScalar bc_len
= (pcell
->ilen
[0]+pcell
->ilen
[pcell
->nb_num
-1])/2;
156 for(int k
=0;k
<pcell
->nb_num
;k
++)
158 int nb
= pcell
->nb_array
[k
];
159 PetscScalar Vj
= x
[zofs
[zone_index
]+nb
]; //potential of nb node
160 grad_P
+= aux
[node
].eps
*(Vj
-Vi
)/pcell
->ilen
[k
]*pcell
->elen
[k
];
163 f
[zofs
[zone_index
]+equ_num
*size
+i
]=grad_P
+pbc
->QF
;
166 //------------------------------------------------------------------------------------------------
168 void ISZone::J1Q_ddm_inner(int i
,PetscScalar
*x
,Mat
*jac
, ODE_Formula
&ODE_F
, vector
<int> & zofs
)
170 const VoronoiCell
* pcell
= &pzone
->davcell
[i
];
171 PetscScalar d_div_grad_P
= 0;
172 for(int j
=0;j
<pcell
->nb_num
;j
++)
174 int nb
= pcell
->nb_array
[j
];
175 PetscScalar value
= pcell
->elen
[j
]/pcell
->ilen
[j
]/pcell
->area
;
176 MatSetValue(*jac
,zofs
[zone_index
]+i
,zofs
[zone_index
]+nb
,value
,INSERT_VALUES
);
177 d_div_grad_P
+= -pcell
->elen
[j
]/pcell
->ilen
[j
]/pcell
->area
;
179 MatSetValue(*jac
,zofs
[zone_index
]+i
,zofs
[zone_index
]+i
,d_div_grad_P
,INSERT_VALUES
);
183 void ISZone::J1Q_ddm_gatebc(int i
,PetscScalar
*x
,Mat
*jac
, ODE_Formula
&ODE_F
, vector
<int> & zofs
,DABC
&bc
)
185 MatSetValue(*jac
,zofs
[zone_index
]+i
,zofs
[zone_index
]+i
,1.0,INSERT_VALUES
);
188 int size
= pzone
->davcell
.size();
189 const VoronoiCell
* pcell
= pzone
->davcell
.GetPointer(i
);
191 for(int j
=0;j
<electrode
.size();j
++)
192 if(electrode
[j
]==pcell
->bc_index
-1)
194 MatSetValue(*jac
,zofs
[zone_index
]+i
,zofs
[zone_index
]+equ_num
*size
+gate_equ
,-1.0,INSERT_VALUES
);
198 void ISZone::J1Q_ddm_semiconductor_insulator_interface(int i
,PetscScalar
*x
,Mat
*jac
, ODE_Formula
&ODE_F
, vector
<int> & zofs
, DABC
&bc
,
201 MatSetValue(*jac
,zofs
[zone_index
]+i
,zofs
[zone_index
]+i
,1.0,INSERT_VALUES
);
202 MatSetValue(*jac
,zofs
[zone_index
]+i
,zofs
[pz
->zone_index
]+5*n
+0,-1.0,INSERT_VALUES
);
206 void ISZone::J1Q_ddm_electrode_insulator_interface(int i
,PetscScalar
*x
,Mat
*jac
, ODE_Formula
&ODE_F
, vector
<int> & zofs
, DABC
&bc
,
209 MatSetValue(*jac
,zofs
[zone_index
]+i
,zofs
[zone_index
]+i
,1.0,INSERT_VALUES
);
210 MatSetValue(*jac
,zofs
[zone_index
]+i
,zofs
[pz
->zone_index
]+n
,-1.0,INSERT_VALUES
);
213 void ISZone::J1Q_ddm_insulator_insulator_interface(int i
,PetscScalar
*x
,Mat
*jac
, ODE_Formula
&ODE_F
, vector
<int> & zofs
, DABC
&bc
,
216 if(zone_index
> pz
->pzone
->zone_index
)
218 MatSetValue(*jac
,zofs
[zone_index
]+i
,zofs
[zone_index
]+i
,1.0,INSERT_VALUES
);
219 MatSetValue(*jac
,zofs
[zone_index
]+i
,zofs
[pz
->zone_index
]+n
,-1.0,INSERT_VALUES
);
222 const VoronoiCell
* pcell
= &pzone
->davcell
[i
];
223 PetscScalar d_div_grad_P
= 0;
224 for(int j
=0;j
<pcell
->nb_num
;j
++)
226 int nb
= pcell
->nb_array
[j
];
227 PetscScalar value
= aux
[i
].eps
*pcell
->elen
[j
]/pcell
->ilen
[j
];
228 MatSetValue(*jac
,zofs
[zone_index
]+i
,zofs
[zone_index
]+nb
,value
,INSERT_VALUES
);
229 d_div_grad_P
+= -aux
[i
].eps
*pcell
->elen
[j
]/pcell
->ilen
[j
];
232 const VoronoiCell
* ncell
= pz
->pzone
->davcell
.GetPointer(n
);
233 for(int j
=0;j
<ncell
->nb_num
;j
++)
235 int nb
= ncell
->nb_array
[j
];
236 PetscScalar value
= pz
->aux
[n
].eps
*ncell
->elen
[j
]/ncell
->ilen
[j
];
237 MatSetValue(*jac
,zofs
[zone_index
]+i
,zofs
[pz
->zone_index
]+nb
,value
,INSERT_VALUES
);
238 d_div_grad_P
+= -pz
->aux
[n
].eps
*ncell
->elen
[j
]/ncell
->ilen
[j
];
240 MatSetValue(*jac
,zofs
[zone_index
]+i
,zofs
[zone_index
]+i
,d_div_grad_P
,INSERT_VALUES
);
243 void ISZone::J1Q_ddm_chargebc(int i
,PetscScalar
*x
,Mat
*jac
, ODE_Formula
&ODE_F
, vector
<int> & zofs
,DABC
&bc
)
246 int size
= pzone
->davcell
.size();
247 const VoronoiCell
* pcell
= pzone
->davcell
.GetPointer(i
);
249 for(int j
=0;j
<electrode
.size();j
++)
250 if(electrode
[j
]==pcell
->bc_index
-1)
251 {charge_equ
=j
;break;}
252 MatSetValue(*jac
,zofs
[zone_index
]+i
,zofs
[zone_index
]+i
,1.0,INSERT_VALUES
);
253 MatSetValue(*jac
,zofs
[zone_index
]+i
,zofs
[zone_index
]+equ_num
*size
+charge_equ
,-1.0,INSERT_VALUES
);
257 void ISZone::J1Q_gate_electrode(int i
,PetscScalar
*x
,Mat
*jac
, ODE_Formula
&ODE_F
, vector
<int> & zofs
, DABC
&bc
, PetscScalar DeviceDepth
)
260 int size
= pzone
->davcell
.size();
261 int bc_index
= electrode
[i
];
262 GateBC
*pbc
= dynamic_cast <GateBC
* > (bc
.Get_pointer(bc_index
));
264 int matrix_row
= zofs
[zone_index
]+equ_num
*size
+i
;
265 PetscScalar R
= pbc
->R
;
266 PetscScalar C
= pbc
->C
;
267 PetscScalar L
= pbc
->L
;
268 PetscScalar e
= mt
->e
;
270 MatAssemblyBegin(*jac
,MAT_FLUSH_ASSEMBLY
);
271 MatAssemblyEnd(*jac
,MAT_FLUSH_ASSEMBLY
);
272 for(int j
=0;j
<bc
[bc_index
].psegment
->node_array
.size();j
++)
274 int node
=bc
[bc_index
].psegment
->node_array
[j
];
275 const VoronoiCell
* pcell
= pzone
->davcell
.GetPointer(node
);
276 for(int k
=0;k
<pcell
->nb_num
;k
++)
278 int nb
= pcell
->nb_array
[k
];
279 //for displacement current
280 PetscScalar dJ_dVi
= aux
[node
].eps
/pcell
->ilen
[k
]/ODE_F
.dt
*pcell
->elen
[k
];
281 PetscScalar dJ_dVr
= -aux
[node
].eps
/pcell
->ilen
[k
]/ODE_F
.dt
*pcell
->elen
[k
];
282 if(pbc
->electrode_type
==VoltageBC
)
284 MatSetValue(*jac
,matrix_row
,zofs
[zone_index
]+node
,DeviceDepth
*dJ_dVi
*(L
/ODE_F
.dt
+R
),ADD_VALUES
);
285 MatSetValue(*jac
,matrix_row
,zofs
[zone_index
]+nb
,DeviceDepth
*dJ_dVr
*(L
/ODE_F
.dt
+R
),ADD_VALUES
);
289 MatAssemblyBegin(*jac
,MAT_FLUSH_ASSEMBLY
);
290 MatAssemblyEnd(*jac
,MAT_FLUSH_ASSEMBLY
);
291 if(pbc
->electrode_type
==VoltageBC
)
292 MatSetValue(*jac
,matrix_row
,matrix_row
,1+(L
/ODE_F
.dt
+R
)*C
/ODE_F
.dt
,INSERT_VALUES
); //dJ/dP
297 void ISZone::J1Q_charge_electrode(int i
,PetscScalar
*x
,Mat
*jac
, ODE_Formula
&ODE_F
, vector
<int> & zofs
, DABC
&bc
)
300 int size
= pzone
->davcell
.size();
301 int bc_index
= electrode
[i
];
302 ChargedContactBC
*pbc
= dynamic_cast <ChargedContactBC
* > (bc
.Get_pointer(bc_index
));
303 int matrix_row
= zofs
[zone_index
]+equ_num
*size
+i
;
305 MatAssemblyBegin(*jac
,MAT_FLUSH_ASSEMBLY
);
306 MatAssemblyEnd(*jac
,MAT_FLUSH_ASSEMBLY
);
307 for(int j
=0;j
<bc
[bc_index
].psegment
->node_array
.size();j
++)
309 int node
=bc
[bc_index
].psegment
->node_array
[j
];
310 const VoronoiCell
* pcell
= pzone
->davcell
.GetPointer(node
);
311 PetscScalar dP_dVi
=0;
312 for(int k
=0;k
<pcell
->nb_num
;k
++)
314 int nb
= pcell
->nb_array
[k
];
315 PetscScalar dP_dVj
=aux
[node
].eps
/pcell
->ilen
[k
]*pcell
->elen
[k
];
316 dP_dVi
+=-aux
[node
].eps
/pcell
->ilen
[k
]*pcell
->elen
[k
];
317 MatSetValue(*jac
,matrix_row
,zofs
[zone_index
]+nb
,dP_dVj
,ADD_VALUES
);
319 MatSetValue(*jac
,matrix_row
,zofs
[zone_index
]+node
,dP_dVi
,ADD_VALUES
);
321 MatAssemblyBegin(*jac
,MAT_FLUSH_ASSEMBLY
);
322 MatAssemblyEnd(*jac
,MAT_FLUSH_ASSEMBLY
);
326 //----------------------------------------------------------------------------------------------
327 void ISZone::F1Q_efield_update(PetscScalar
*x
,vector
<int> & zofs
, DABC
&bc
,vector
<BZoneData
*>zonedata
)
329 PetscScalar a
=0,b
=0,c
=0,P_x
=0,P_y
=0,w
=0;
330 //calculate Ex Ey with least-squares gradient construction
331 for(int i
=0; i
<pzone
->davcell
.size();i
++)
333 a
=0,b
=0,c
=0,P_x
=0,P_y
=0,w
=0;
334 const VoronoiCell
*pcell
= pzone
->davcell
.GetPointer(i
);
335 for(int k
=0;k
<pcell
->nb_num
;k
++)
337 const VoronoiCell
*ncell
= pzone
->davcell
.GetPointer(pcell
->nb_array
[k
]);
338 PetscScalar dx
= ncell
->x
- pcell
->x
;
339 PetscScalar dy
= ncell
->y
- pcell
->y
;
340 PetscScalar dP
= x
[zofs
[zone_index
]+pcell
->nb_array
[k
]]- x
[zofs
[zone_index
]+i
];
341 w
=1.0/sqrt(dx
*dx
+dy
*dy
);
348 if(pcell
->bc_index
&&bc
[pcell
->bc_index
-1].BCType
==InsulatorInterface
)
350 InsulatorInterfaceBC
*pbc
;
351 pbc
= dynamic_cast<InsulatorInterfaceBC
*>(bc
.Get_pointer(pcell
->bc_index
-1));
352 int n_zone
= pbc
->pinterface
->Find_neighbor_zone_index(zone_index
);
353 int n_node
= pbc
->pinterface
->Find_neighbor_node_index(zone_index
,i
);
354 SMCZone
* pz
= dynamic_cast<SMCZone
*>(zonedata
[n_zone
]);
355 const VoronoiCell
* dcell
= pz
->pzone
->davcell
.GetPointer(n_node
);
356 for(int k
=1;k
<dcell
->nb_num
-1;k
++)
358 const VoronoiCell
*ncell
= pz
->pzone
->davcell
.GetPointer(dcell
->nb_array
[k
]);
359 PetscScalar dx
= ncell
->x
- dcell
->x
;
360 PetscScalar dy
= ncell
->y
- dcell
->y
;
361 PetscScalar dP
= x
[zofs
[n_zone
]+5*dcell
->nb_array
[k
]+0]- x
[zofs
[n_zone
]+5*n_node
+0];
362 w
=1.0/sqrt(dx
*dx
+dy
*dy
);
370 aux
[i
].Ex
= -(c
*P_x
-b
*P_y
)/(a
*c
-b
*b
);
371 aux
[i
].Ey
= -(a
*P_y
-b
*P_x
)/(a
*c
-b
*b
);