17#pragma warning(disable:4503)
29#include "viennagrid/config/default_configs.hpp"
67template<
typename DeviceType>
70 typedef typename DeviceType::segment_type SegmentType;
73 SegmentType
const & gate_contact = device.segment(1);
74 SegmentType
const & source_contact = device.segment(2);
75 SegmentType
const & gate_oxide = device.segment(3);
76 SegmentType
const & drain_contact = device.segment(4);
77 SegmentType
const & source = device.segment(5);
78 SegmentType
const & drain = device.segment(6);
79 SegmentType
const & body = device.segment(7);
80 SegmentType
const & body_contact = device.segment(8);
83 std::cout <<
"* init_device(): Setting materials..." << std::endl;
101 std::cout <<
"* init_device(): Setting doping..." << std::endl;
102 device.set_doping_n(1e24, source);
103 device.set_doping_p(1e8, source);
105 device.set_doping_n(1e24, drain);
106 device.set_doping_p(1e8, drain);
108 device.set_doping_n(1e17, body);
109 device.set_doping_p(1e15, body);
115 device.set_contact_potential(0.8, gate_contact);
116 device.set_contact_potential(0.0, source_contact);
117 device.set_contact_potential(1.0, drain_contact);
118 device.set_contact_potential(0.0, body_contact);
130int main(
int argc,
char **argv)
138 typedef DeviceType::segment_type SegmentType;
139 std::string LSolver(
"petsc_parallel_linear_solver");
141 PetscInitialize(&argc, &argv, (
char *) 0,
"");
142 MPI_Comm_size(MPI_COMM_WORLD, &size);
143 opt = getopt(argc, argv,
"seS");
145 case 's': s = 1;h = 0 ;
break;
146 case 'e': s = 0;h = 1 ;
break;
147 case 'S': s = 0;h = 0 ;
break;
164 std::cout <<
"* main(): Creating and scaling device..." << std::endl;
166 device.load_mesh(
"../examples/data/mosfet840.mesh");
173 std::cout <<
"* main(): Initializing device..." << std::endl;
175 if(s == 1 && size != 1) device.refine(size) ;
187 std::cout <<
"* main(): Creating DD simulator..." << std::endl;
222 std::cout <<
"* main(): Launching DD simulator..." << std::endl;
239 SegmentType
const & drain_contact = device.segment(4);
240 SegmentType
const & body = device.segment(6);
267 <<
"* main(): Setting up first-order SHE (semi-self-consistent using 40 Gummel iterations)..."
314 std::cout <<
"* main(): Computing first-order SHE..." << std::endl;
318 dd_simulator.potential());
320 dd_simulator.electron_density());
322 dd_simulator.hole_density());
339 <<
"* main(): Writing energy distribution function from first-order SHE result..."
370 <<
"* main(): Results can now be viewed with your favorite VTK viewer (e.g. ParaView)."
373 <<
"* main(): Don't forget to scale the z-axis by about a factor of 1e12 when examining the distribution function."
375 std::cout << std::endl;
376 std::cout <<
"*********************************************************"
378 std::cout <<
"* ViennaSHE finished successfully *"
380 std::cout <<
"*********************************************************"
The main SHE configuration class. To be adjusted by the user for his/her needs.
long max_expansion_order() const
Returns the current maximum expansion order.
bool with_holes() const
Returns true if holes are considered in the simulation.
nonlinear_solver_config_type & nonlinear_solver()
Returns the configuration object for the nonlinear solver.
void set_electron_equation(equation_id equ_id)
linear_solver_config_type & linear_solver()
Returns the configuration object for the linear solver.
double energy_spacing() const
Returns the uniform energy spacing of discrete energies.
void set_hole_equation(equation_id equ_id)
bool with_electrons() const
Returns true if electrons are considered in the simulation.
viennashe::she::scatter_config & scattering()
Returns the configuration object for scattering.
Defines the physical properties of a device, e.g. doping. This is the implementation for 2d and highe...
ionized_impurity_scattering_parameters const & ionized_impurity() const
Returns the parameters for ionized impurity scattering. Const-version.
acoustic_phonon_scattering_parameters const & acoustic_phonon() const
Returns the parameters for acoustic phonon scattering. Const-version.
optical_phonon_scattering_parameters const & optical_phonon() const
Returns the parameters for optical phonon scattering. Const-version.
Class for self-consistent SHE simulations.
void set(long solver_id)
Sets a new linear solver. Provide IDs defined in.
std::size_t max_iters() const
Returns the maximum number of iterative linear solver iterations.
void setArgv(char **argv)
std::size_t max_iters() const
Returns the maximum number of nonlinear solver iterations.
double damping() const
Returns the damping for the nonlinear solver.
void setThreshold(long threshold)
Convenience header, which includes all core functionality available in ViennaSHE.
void init_device(DeviceType &device, double Vg_init, double Nd, double Na)
Initalizes the MOS 1D device for testing.
std::string hole_density()
std::string electron_density()
std::string preamble()
Prints the ViennaSHE preamble (header). Used in all the examples as well as the standalone-applicatio...
A class referring to hafnium dioxide.
A class referring to silicon and providing certain material parameters Note that this is the default ...
static const double q
Elementary charge.