Implementing injection with current density #48

Merged
JorgeGonz merged 1 commit from feature/Amps_per_m2 into development 2024-03-29 13:09:30 +01:00
8 changed files with 133 additions and 77 deletions
Showing only changes of commit d86b3a3417 - Show all commits

Implementing injection with current density

WARNING: This current denstiy will be multiplied by the reference
length, no the surface area that is being used for injection!

New units in the injection of particles 'Am2' to inject a density
current. Manual has been modified accordingly.

Reference parameters are now also printed in the case folder.
Jorge Gonzalez 2024-03-28 09:45:46 +01:00

Binary file not shown.

View file

@ -72,7 +72,7 @@
The \Gls{fpakc} is a simulation tool that models species in plasma (ions, electrons and neutrals) following the trajectories of macro-particles as they move and interact between them and the boundaries of the domain. The \Gls{fpakc} is a simulation tool that models species in plasma (ions, electrons and neutrals) following the trajectories of macro-particles as they move and interact between them and the boundaries of the domain.
Particles properties are scattered into a finite element mesh in 1, 2 or three dimensions, with the possibility to choose different geometries. Particles properties are scattered into a finite element mesh in 1, 2 or three dimensions, with the possibility to choose different geometries.
The official repository can be found at: \url{https://gitlab.com/JorgeGonz/fpakc.git}. The official repository can be found at: \url{https://gitlab.com/JorgeGonz/fpakc.git}.
The code is currently in very early steps of development and further improvements are expected very soon. The code is currently in the very early steps of development and further refinements are expected very soon.
%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%% %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
\section{Main Guidelines} \section{Main Guidelines}
@ -86,11 +86,11 @@
\item \acrshort{fpakc} is coded in a \textit{understandable} way. \item \acrshort{fpakc} is coded in a \textit{understandable} way.
This means that the code is required to be written in a clear way that is easy to understand and maintain. This means that the code is required to be written in a clear way that is easy to understand and maintain.
Variables and procedure names need to be self-understanding. Variables and procedure names need to be self-understanding.
This ease the process of fixing bugs and improving the codes by a large team of developers. This eases the process of fixing bugs and improving the codes by a large team of developers.
For more information, please refer to the \acrshort{fpakc} Coding Style document. For more information, please refer to the \acrshort{fpakc} Coding Style document.
\item \acrshort{fpakc} requires to be ease to use. \item \acrshort{fpakc} requires being ease to use.
Input files are required to be in a \textit{human} format, meaning that the different options can be easily understood without constant reference to the user guide. Input files are required to be in a \textit{human} format, meaning that the different options can be easily understood without constant reference to the user guide.
\acrshort{fpakc} is aimed to be used in a wide range of applications and by a variety of scientist: from very established ones to newcomers to the field and also students. \acrshort{fpakc} is aimed to be used in a wide range of applications and by various scientists: from well-established ones to newcomers to the field and also students.
\end{enumerate} \end{enumerate}
These are foundation stones of \acrshort{fpakc} and its development and should always be followed, at least for the releases in the official repository. These are foundation stones of \acrshort{fpakc} and its development and should always be followed, at least for the releases in the official repository.
@ -105,16 +105,16 @@
\section{The Particle Method} \section{The Particle Method}
\Gls{fpakc} uses macro-particles to simulate the dynamics of different plasma species (mainly ions, electrons and neutrals). \Gls{fpakc} uses macro-particles to simulate the dynamics of different plasma species (mainly ions, electrons and neutrals).
These macro-particles could represent a large amount of real particles. These macro-particles could represent a large amount of real particles.
For now own, macro-particles will be referred as just particles by abusing of language. For now own, macro-particles will be referred as just particles by abuse of language.
During the initiation phase, the input and mesh file(s) are reading. During the initiation phase, the input and mesh file(s) are reading.
If an initial distribution for a species is specified in the input file, particles to match that distribution are loaded into the cells. If an initial distribution for a species is specified in the input file, particles to match that distribution are loaded into the cells.
The general steps performed in each iteration are: The general steps performed in each iteration are:
\begin{enumerate} \begin{enumerate}
\item Firstly, new particles are introduced into the domain as specified in the input file. \item Firstly, new particles are introduced into the domain as specified in the input file.
\item Particles are then pushed accounting for possible acceleration by external forces. \item Particles are then pushed, accounting for possible acceleration by external forces.
During this process, if a particle changes cell it is found using the connectivity between elements. During this process, if a particle changes cell, it is found using the connectivity between elements.
If a particle encounters a boundary instead a new cell, the interaction between the boundary and the wall is computed. If a particle encounters a boundary instead a new cell, the interaction between the boundary and the wall are computed.
A particle may abandon the computational domain and is no longer accounted for. A particle may abandon the computational domain and is no longer accounted for.
\item Next, collisions for the particles inside each cell are carried out. \item Next, collisions for the particles inside each cell are carried out.
This may include different collision processes for each particle. This may include different collision processes for each particle.
@ -124,10 +124,10 @@
\item Finally, particle properties are scattered among the mesh nodes. \item Finally, particle properties are scattered among the mesh nodes.
These properties are density, momentum and the stress tensor. These properties are density, momentum and the stress tensor.
\item If requested, the electromagnetic field is computed. \item If requested, the electromagnetic field is computed.
\item If the number of iteration requires writing output files, it is done after all steps for the particles is completed. \item If the number of iteration requires writing output files, it is done after all steps for the particles are completed.
\end{enumerate} \end{enumerate}
\Gls{fpakc} has the capability to configure all the behavior of the simulation via the input file. \Gls{fpakc} has the capability to configure all the behaviour of the simulation via the input file.
Parameters as injection, the kind of pusher used for each species, boundary conditions or collisions are user-input parameters and will be described in Chap.~\ref{ch:input_file}. Parameters as injection, the kind of pusher used for each species, boundary conditions or collisions are user-input parameters and will be described in Chap.~\ref{ch:input_file}.
%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%% %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
@ -168,8 +168,8 @@
%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%% %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
\section{Find new cell} \section{Find new cell}
Once the position and velocity of the particle are updated, the new cell that contains the particle is searched. Once the position and velocity of the particle are updated, the new cell that contains the particle is searched.
This is done by a neighbor search, starting from the previous cell containing the particle. This is done by a neighbour search, starting from the previous cell containing the particle.
In the process of finding the new cell, it is possible that a particle encounters a boundary. In the process of finding the new cell, a particle might encounter a boundary.
When the particle interacts with the boundary, the particle may continue its life in the simulation or might be eliminated from it. When the particle interacts with the boundary, the particle may continue its life in the simulation or might be eliminated from it.
Once that the new cell is found or that the particle life has been terminated, the pushing is complete. Once that the new cell is found or that the particle life has been terminated, the pushing is complete.
If a secondary mesh is used for the Monte-Carlo Collision method, the new cell in that mesh in which the particle reside is also found by the same method, although no interaction with the boundaries is accounted for this step. If a secondary mesh is used for the Monte-Carlo Collision method, the new cell in that mesh in which the particle reside is also found by the same method, although no interaction with the boundaries is accounted for this step.
@ -178,7 +178,7 @@
\section{Variable Weighting Scheme\label{sec:weightingScheme}} \section{Variable Weighting Scheme\label{sec:weightingScheme}}
One of the issues in particle simulations, specially for axial-symmetrical cases, is that due to the disparate volume of cells, specially close to the axis, the statistics in some cells is usually poor. One of the issues in particle simulations, specially for axial-symmetrical cases, is that due to the disparate volume of cells, specially close to the axis, the statistics in some cells is usually poor.
To try to fix that, the possibility to include a Variable Weighting Scheme in the simulations is available in \Gls{fpakc}. To try to fix that, the possibility to include a Variable Weighting Scheme in the simulations is available in \Gls{fpakc}.
These schemes detect when a particle change cells and split it if necessary to improve statistics. These schemes detect when a particle changes cells and split it if necessary to improve statistics.
The use of a Variable Weighting Scheme is defined by the user in the input file. The use of a Variable Weighting Scheme is defined by the user in the input file.
Beware that this can increase the number of particles in the simulation and increase computational time. Beware that this can increase the number of particles in the simulation and increase computational time.
@ -189,16 +189,16 @@
\Gls{fpakc} distinguish between two types of interactions: \acrfull{mcc} and \acrfull{cs}. \Gls{fpakc} distinguish between two types of interactions: \acrfull{mcc} and \acrfull{cs}.
\acrshort{mcc} refers to the process in which two particles interact in short range. \acrshort{mcc} refers to the process in which two particles interact in short range.
These processes include, but are not limited to: elastic collisions, ionization/recombination, charge-exchange, excitation/de-excitation\ldots These processes include, but are not limited to: elastic collisions, ionization/recombination, charge-exchange, excitation/de-excitation\ldots
A secondary mesh, with cell sizes in the range of the mean-free path, can be used for this type of collisions. A secondary mesh, with cell sizes in the range of the mean-free path, can be used for this type of collision.
\acrshort{cs} refers to the large range interaction that a charged species suffer do to the charge of other particles. \acrshort{cs} refers to the large range interaction that a charged species suffer do to the charge of other particles.
The interactions between the different species is defined by the user. The interactions between the different species is defined by the user.
%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%% %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
\subsection{\acrlong{mcc}} \subsection{\acrlong{mcc}}
For each cell the maximum number of collisions between particle is computed. For each cell, the maximum number of collisions between particle is computed.
For each collision, a random pair of particles is chosen. For each collision, a random pair of particles is chosen.
A loop over all possible collisions for the pair of particles chosen is performed. A loop over all possible collisions for the pair of particles chosen is performed.
If a random number is above the probability of collision for that specific type, the collision take place. If a random number is above the probability of collision for that specific type, the collision takes place.
If not, the next type for the particle pair is checked. If not, the next type for the particle pair is checked.
Below are described the type of collision process implemented in \acrshort{fpakc}: Below are described the type of collision process implemented in \acrshort{fpakc}:
@ -219,7 +219,7 @@
\item Recombination. \item Recombination.
When an electron and an ion interact, there is a possibility for them to be recombined into a neutral particle. When an electron and an ion interact, there is a possibility for them to be recombined into a neutral particle.
The photon emitted by this process is not modelled yet. The photons emitted by this process are not modelled yet.
\end{itemize} \end{itemize}
%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%% %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
@ -238,7 +238,7 @@
Once that the pushing is complete, the array of particles that remain inside the domain is copied to a new array. Once that the pushing is complete, the array of particles that remain inside the domain is copied to a new array.
The new array containing only the particles inside the domain will be the one used in the next steps. The new array containing only the particles inside the domain will be the one used in the next steps.
In this section, particles are assigned to the list of particles inside each individual cell. In this section, particles are assigned to the list of particles inside each individual cell.
Unfortunately, this is done right now without parallelisation and is very CPU consuming. Unfortunately, this is done right now without parallelization and is very CPU consuming.
%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%% %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
\section{Probing}\label{sec:probing} \section{Probing}\label{sec:probing}
@ -250,21 +250,21 @@
The user can decide the grid width and the number of points in each direction. The user can decide the grid width and the number of points in each direction.
The distribution function will be calculated and wrote with a time step decided by the user. The distribution function will be calculated and wrote with a time step decided by the user.
If a particle velocity resides outside of the velocity grid (in any direction), it wont be added to the tally of the distribution function. If a particle velocity resides outside the velocity grid (in any direction), it will not be added to the tally of the distribution function.
Due to the limitation of only taking into account particles in the cell, and not neighbour particles, two probes for the same species at different positions but in the same cell will output the same results. Due to the limitation of only considering particles in the cell, and not neighbour particles, two probes for the same species at different positions but in the same cell will output the same results.
A more advance method taking into account distance between the particles and the probe position as well as particles in neighbour cells could be implemented to improve the statistics of the distribution function. A more advance method considering distance between the particles and the probe position as well as particles in neighbour cells could be implemented to improve the statistics of the distribution function.
%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%% %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
\section{Scattering} \section{Scattering}
The properties of each particle are deposited in the nodes from the containing cell. The properties of each particle are deposited in the nodes from the containing cell.
This process depend on the cell type, but in general, each node receive a proportional part of the particle properties as a function of the particle position inside the cell. This process depends on the cell type, but in general, each node receives a proportional part of the particle properties as a function of the particle position inside the cell.
Figure \ref{fig:scatteringQuad} shows how a particle at a generic position $p(x_1, x_2)$ inside the cell is scattered to the four nodes. The figure \ref{fig:scatteringQuad} shows how a particle at a generic position $p(x_1, x_2)$ inside the cell is scattered to the four nodes.
\begin{wrapfigure}{l}{0.4\textwidth} \begin{wrapfigure}{l}{0.4\textwidth}
\centering \centering
\includegraphics{figures/scatteringQuad} \includegraphics{figures/scatteringQuad}
\caption{\label{fig:scatteringQuad}Example of how a particle is weighted in a quadrilateral cell.} \caption{\label{fig:scatteringQuad}Example of how a particle is weighted in a quadrilateral cell.}
\end{wrapfigure} \end{wrapfigure}
Each node receives a proportional part of the area formed by dividing the cell in for rectangles using as an additional vertex the particle position. Each node receives a proportional part of the area formed by dividing the cell in for rectangles, using as an additional vertex the particle position.
These properties are dimensionless, but they are converted to the correct units once the output is printed. These properties are dimensionless, but they are converted to the correct units once the output is printed.
%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%% %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
@ -273,11 +273,11 @@
%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%% %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
\section{Average scheme} \section{Average scheme}
Particle-in-cell codes has an intrinsic statistical noise associated with them. Particle-in-cell codes have an intrinsic statistical noise associated with them.
Although this can be reduced by increasing the number of particles, this also increases the CPU requirements of the case. Although this can be reduced by increasing the number of particles, this also increases the CPU requirements of the case.
It is quite common that most cases reach a quasi-steady state after a number of iterations and time-average results can be obtained after to improve analysis, plotting and restarting the case using these time-average results as new species backgrounds. It is quite common that most cases reach a quasi-steady state after a number of iterations and time-average results can be obtained after to improve analysis, plotting and restarting the case using these time-average results as new species backgrounds.
Although this is possible to do once the simulation is finished with post-processing tools, this is limited to the amount of iterations printed. Although this is possible to do once the simulation is finished with post-processing tools, this is limited to the number of iterations printed.
\Gls{fpakc} implements a simple average scheme that, after a start time provided by the user, scores a mean and standard deviation of all the main species properties, and the electromagnetic field. \Gls{fpakc} implements a simple average scheme that, after a start time provided by the user, scores a mean and standard deviation of all the main species properties, and the electromagnetic field.
This scheme is based on the Welford's online algorithm~\cite{welford1962note}. This scheme is based on the Welford's online algorithm~\cite{welford1962note}.
@ -286,7 +286,7 @@
%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%% %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
\chapter{Installation} \chapter{Installation}
\section{Required Packages} \section{Required Packages}
In order to properly compile \gls{fpakc}, the following packages are required. To properly compile \gls{fpakc}, the following packages are required.
\subsection{Gfortran} \subsection{Gfortran}
The \Gls{opensource} free compiler \Gls{gfortran}\cite{gfortranURL} from GCC is the basic way to compile \acrshort{fpakc}. The \Gls{opensource} free compiler \Gls{gfortran}\cite{gfortranURL} from GCC is the basic way to compile \acrshort{fpakc}.
It is distributed with all GNU/Linux distributions. It is distributed with all GNU/Linux distributions.
@ -369,7 +369,7 @@ make
%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%% %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
\section{Case file} \section{Case file}
The required format for the case file is \Gls{json}. The required format for the case file is \Gls{json}.
\Gls{json} is a case-sensitive format, so input must be written with the correct capitalisation. \Gls{json} is a case-sensitive format, so input must be written with the correct capitalization.
The basic structure and options available for the case file are explained below. The basic structure and options available for the case file are explained below.
The order of the objects and variables is irrelevant, but the structure needs to be maintained. The order of the objects and variables is irrelevant, but the structure needs to be maintained.
%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%% %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
@ -380,9 +380,9 @@ make
\item \textbf{path}: Character. \item \textbf{path}: Character.
Path for the output files. This path is also used to locate the mesh input file. Path for the output files. This path is also used to locate the mesh input file.
\item \textbf{folder}: Character. \item \textbf{folder}: Character.
Base name of the folder in wich output files are placed. Base name of the folder in which output files are placed.
The date and time is appended to this name. The date and time is appended to this name.
If none is provided, only the date and time is writted as the folder name. If none is provided, only the date and time is written as the folder name.
\item \textbf{triggerOutput}: Integer. \item \textbf{triggerOutput}: Integer.
Determines the number of iterations between writing output files for macroscopic quantities. Determines the number of iterations between writing output files for macroscopic quantities.
\item \textbf{cpuTime}: Logical. \item \textbf{cpuTime}: Logical.
@ -515,7 +515,7 @@ make
\item \textbf{absorption}: Particle is eliminated from the domain. \item \textbf{absorption}: Particle is eliminated from the domain.
The particle is first moved into the edge and its properties are scattered among the edge nodes. The particle is first moved into the edge and its properties are scattered among the edge nodes.
\item \textbf{transparent}: Particle abandon the numerical domain. \item \textbf{transparent}: Particle abandon the numerical domain.
\item \textbf{wallTemperature}: Reflective wall with cosntant temperature that exchange heat with particles. \item \textbf{wallTemperature}: Reflective wall with constant temperature that exchange heat with particles.
Required parameters are: Required parameters are:
\begin{itemize} \begin{itemize}
\item \textbf{temperature}: Real. \item \textbf{temperature}: Real.
@ -526,8 +526,8 @@ make
Specific heat capacity of the material. Specific heat capacity of the material.
\end{itemize} \end{itemize}
\item \textbf{ionization}: Per each particle crossing the surface with this type of boundary, a number of ionization events are calculated. \item \textbf{ionization}: Per each particle crossing the surface with this type of boundary, a number of ionization events are calculated.
A pair of ion-electron is generated for each ionization event taking as a reference a neutral background. A pair of ion-electron is generated for each ionization event, taking as a reference a neutral background.
Secondary electron is taken as same type as incident particle. The secondary electron is taken as the same type as the incident particle.
The available input is: The available input is:
\begin{itemize} \begin{itemize}
\item \textbf{neutral}: Object. \item \textbf{neutral}: Object.
@ -540,7 +540,7 @@ make
\item \textbf{mass}: Real. \item \textbf{mass}: Real.
Units in $\unit{kg}$. Units in $\unit{kg}$.
Mass of neutral species. Mass of neutral species.
If missing, the mass of the ion is ussed If missing, the mass of the ion is used
\item \textbf{density}: Real. \item \textbf{density}: Real.
Units in $\unit{m^{-3}}$. Units in $\unit{m^{-3}}$.
Density of neutral background. Density of neutral background.
@ -558,18 +558,18 @@ make
\end{itemize} \end{itemize}
\item \textbf{effectiveTime}: Real. \item \textbf{effectiveTime}: Real.
Units in $\unit{s}$. Units in $\unit{s}$.
As the particle is no longer simulated once it crossed the boundary, this time represent the effective time in which the particle produces ionization processes in the neutral background. As the particle is no longer simulated once it crossed the boundary, this time represents the effective time in which the particle produces ionization processes in the neutral background.
Required parameter. Required parameter.
\item \textbf{energyThreashold}: Real. \item \textbf{energyThreashold}: Real.
Units in $\unit{eV}$. Units in $\unit{eV}$.
Ionization energy threshold for the simulated process. Ionization energy threshold for the simulated process.
Required parameter. Required parameter.
\item \textbf{crossSection}: Character. \item \textbf{crossSection}: Character.
Complete path to the cross section data for the ionization process. Complete path to the cross-section data for the ionization process.
\end{itemize} \end{itemize}
\item \textbf{axis}: Identifies the symmetry axis for 2D cylindrical simulations. \item \textbf{axis}: Identifies the symmetry axis for 2D cylindrical simulations.
If for some reason a particle interact with this axis, it is reflected. If , for some reason, a particle interacts with this axis, it is reflected.
\end{itemize} \end{itemize}
\end{itemize} \end{itemize}
@ -596,7 +596,7 @@ make
%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%% %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
\subsection{inject} \subsection{inject}
The array \textbf{inject} specifies the injection of particles from different surfaces. The array \textbf{inject} specifies the injection of particles from different surfaces.
The injection of particles need to be associated to a physicalSurface in the mesh file. The injection of particles needs to be associated to a physicalSurface in the mesh file.
Multiple injections can be associated to the same surface. Multiple injections can be associated to the same surface.
\begin{itemize} \begin{itemize}
\item \textbf{name}: Character. \item \textbf{name}: Character.
@ -610,7 +610,9 @@ make
Available values are: Available values are:
\begin{itemize} \begin{itemize}
\item \textbf{A}: Ampere. \item \textbf{A}: Ampere.
\item \textbf{sccm}: Standard cubic centimeter. \item \textbf{Am2}: Ampere per square meter.
This value will be multiplied by the square of the reference length to convert to A, it will not consider the surface area of injection.
\item \textbf{sccm}: Standard cubic centimetre.
\item \textbf{part/s}: Particles (real) per second. \item \textbf{part/s}: Particles (real) per second.
\end{itemize} \end{itemize}
\item \textbf{v}: Real. \item \textbf{v}: Real.
@ -627,7 +629,7 @@ make
\begin{itemize} \begin{itemize}
\item \textbf{Maxwellian}: Maxwellian distribution of temperature \textbf{T} and mean \textbf{v} times the value of \textbf{n} in the specified direction. \item \textbf{Maxwellian}: Maxwellian distribution of temperature \textbf{T} and mean \textbf{v} times the value of \textbf{n} in the specified direction.
\item \textbf{Half-Maxwellian}: Half-Maxwellian distribution of temperature \textbf{T} and mean \textbf{v} times the value of \textbf{n} in the specified direction. \item \textbf{Half-Maxwellian}: Half-Maxwellian distribution of temperature \textbf{T} and mean \textbf{v} times the value of \textbf{n} in the specified direction.
Only takes into account the positive part of the half-Maxwellian. Only considers the positive part of the half-Maxwellian.
\item \textbf{Delta}: Dirac's delta distribution function. All particles are injected with velocity \textbf{v} times the value of \textbf{n} in the specified direction. \item \textbf{Delta}: Dirac's delta distribution function. All particles are injected with velocity \textbf{v} times the value of \textbf{n} in the specified direction.
\end{itemize} \end{itemize}
\item \textbf{T}: Real. \item \textbf{T}: Real.
@ -651,7 +653,7 @@ make
\item \textbf{radius}: Real. \item \textbf{radius}: Real.
Reference atomic radius in $\unit{m}$. Reference atomic radius in $\unit{m}$.
\item \textbf{crossSection}: Real. \item \textbf{crossSection}: Real.
Reference cross section in $\unit{m^2}$. Reference cross-section in $\unit{m^2}$.
If this value is present, radius is ignored. If this value is present, radius is ignored.
\end{itemize} \end{itemize}
@ -677,8 +679,8 @@ make
Indicates the type of pusher used for each species: Indicates the type of pusher used for each species:
\begin{itemize} \begin{itemize}
\item \textbf{Neutral}: Pushes a particle without any external force. \item \textbf{Neutral}: Pushes a particle without any external force.
\item \textbf{Electrostatic}: Pushes a particle including the effect of the electrostatic field. \item \textbf{Electrostatic}: Pushes a particle, including the effect of the electrostatic field.
\item \textbf{Electromagnetic}: Pushes particles accounting for the electromagnetic field. \item \textbf{Electromagnetic}: Pushes a particle, accounting for the electromagnetic field.
\end{itemize} \end{itemize}
\item \textbf{WeightingScheme}: Character. \item \textbf{WeightingScheme}: Character.
Indicates the variable weighting scheme to be used in the simulation. Indicates the variable weighting scheme to be used in the simulation.
@ -726,11 +728,11 @@ make
\end{itemize} \end{itemize}
%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%% %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
\subsection{interactions}\label{ssec:input_interactions} \subsection{interactions}\label{ssec:input_interactions}
This object determine the different interactions among species. This object determines the different interactions among species.
Acceptable values are: Acceptable values are:
\begin{itemize} \begin{itemize}
\item \textbf{folderCollisions}: Character. \item \textbf{folderCollisions}: Character.
Indicates the path to in which the cross section tables are allocated. Indicates the path to in which the cross-section tables are allocated.
\item \textbf{meshCollisions}: Character. \item \textbf{meshCollisions}: Character.
Determines a specific mesh for \acrshort{mcc} processes. Determines a specific mesh for \acrshort{mcc} processes.
The file needs to be located in the folder \textbf{output.folder}. The file needs to be located in the folder \textbf{output.folder}.
@ -757,7 +759,7 @@ make
Accepted values are \textbf{elastic}, \textbf{chargeExchange}, \textbf{ionization} and \textbf{recombination}. Accepted values are \textbf{elastic}, \textbf{chargeExchange}, \textbf{ionization} and \textbf{recombination}.
Please refer to Sec.~\ref{ssec:collisions} for a description of the different collision types. Please refer to Sec.~\ref{ssec:collisions} for a description of the different collision types.
\item \textbf{crossSection}: Character. \item \textbf{crossSection}: Character.
File in \textbf{interactions.folderCollisions} that contains the cross section data as a 1D table of relative energy (in $\unit{eV}$) and cross section (in $\unit{m^-2}$). File in \textbf{interactions.folderCollisions} that contains the cross-section data as a 1D table of relative energy (in $\unit{eV}$) and cross-section (in $\unit{m^-2}$).
\item \textbf{energyThreshold}: Real. \item \textbf{energyThreshold}: Real.
Energy threshold of the collisional process in $\unit{eV}$. Energy threshold of the collisional process in $\unit{eV}$.
Only valid for \textbf{ionization} and \textbf{recombination} processes. Only valid for \textbf{ionization} and \textbf{recombination} processes.
@ -778,7 +780,7 @@ make
\begin{itemize} \begin{itemize}
\item \textbf{species\_i}, \textbf{species\_j}: Character. \item \textbf{species\_i}, \textbf{species\_j}: Character.
Define the two species involved in the collision processes. Define the two species involved in the collision processes.
Order is indiferent. Order is indifferent.
\end{itemize} \end{itemize}
\end{itemize} \end{itemize}
@ -804,9 +806,9 @@ make
%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%% %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
\section{1D Emissive Cathode (1D\_Cathode)} \section{1D Emissive Cathode (1D\_Cathode)}
Emission from a 1D cathode in both, cartesian and radial coordinates. Emission from a 1D cathode in both, Cartesian and radial coordinates.
Both cases insert the same amount of electrons from the minimum coordinate and have the same boundary conditions for particles and the electrostatic field. Both cases insert the same number of electrons from the minimum coordinate and have the same boundary conditions for particles and the electrostatic field.
This case is useful to ilustrate hoy \acrshort{fpakc} can deal with different geometries by just modifying some parameters in the input file. This case is useful to illustrate how \acrshort{fpakc} can deal with different geometries by just modifying some parameters in the input file.
The same mesh file (\lstinline|mesh.msh|) is used for both cases. The same mesh file (\lstinline|mesh.msh|) is used for both cases.
%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%% %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%

View file

@ -27,6 +27,10 @@ PROGRAM fpakc
!Reads the json configuration file !Reads the json configuration file
CALL readConfig(inputFile) CALL readConfig(inputFile)
!Create output folder and initial files
CALL initOutput(inputFile)
!Do '0' iteration !Do '0' iteration
t = tInitial t = tInitial

View file

@ -84,20 +84,6 @@ MODULE moduleInput
CALL readParallel(config) CALL readParallel(config)
CALL checkStatus(config, "readParallel") CALL checkStatus(config, "readParallel")
!If everything is correct, creates the output folder
CALL EXECUTE_COMMAND_LINE('mkdir ' // path // folder )
!Copies input file to output folder
CALL EXECUTE_COMMAND_LINE('cp ' // inputFile // ' ' // path // folder)
!Copies particle mesh
IF (mesh%dimen > 0) THEN
CALL EXECUTE_COMMAND_LINE('cp ' // pathMeshParticle // ' ' // path // folder)
IF (doubleMesh) THEN
CALL EXECUTE_COMMAND_LINE('cp ' // pathMeshColl // ' ' // path // folder)
END IF
END IF
END SUBROUTINE readConfig END SUBROUTINE readConfig
!Checks the status of the JSON case file and, if failed, exits the execution. !Checks the status of the JSON case file and, if failed, exits the execution.
@ -813,8 +799,8 @@ MODULE moduleInput
REAL(8), DIMENSION(:), ALLOCATABLE:: v0 REAL(8), DIMENSION(:), ALLOCATABLE:: v0
REAL(8):: effTime REAL(8):: effTime
REAL(8):: eThreshold !Energy threshold REAL(8):: eThreshold !Energy threshold
INTEGER:: speciesID INTEGER:: speciesID, electronSecondaryID
CHARACTER(:), ALLOCATABLE:: speciesName, crossSection CHARACTER(:), ALLOCATABLE:: speciesName, crossSection, electronSecondary
LOGICAL:: found LOGICAL:: found
INTEGER:: nTypes INTEGER:: nTypes
@ -869,9 +855,18 @@ MODULE moduleInput
CALL config%get(object // '.crossSection', crossSection, found) CALL config%get(object // '.crossSection', crossSection, found)
IF (.NOT. found) CALL criticalError("missing parameter 'crossSection' for neutrals in ionization", 'readBoundary') IF (.NOT. found) CALL criticalError("missing parameter 'crossSection' for neutrals in ionization", 'readBoundary')
CALL config%get(object // '.electronSecondary', electronSecondary, found)
electronSecondaryID = speciesName2Index(electronSecondary)
IF (found) THEN
CALL initIonization(boundary(i)%bTypes(s)%obj, species(s)%obj%m, m0, n0, v0, T0, &
speciesID, effTime, crossSection, eThreshold,electronSecondaryID)
ELSE
CALL initIonization(boundary(i)%bTypes(s)%obj, species(s)%obj%m, m0, n0, v0, T0, & CALL initIonization(boundary(i)%bTypes(s)%obj, species(s)%obj%m, m0, n0, v0, T0, &
speciesID, effTime, crossSection, eThreshold) speciesID, effTime, crossSection, eThreshold)
END IF
CASE('wallTemperature') CASE('wallTemperature')
CALL config%get(object // '.temperature', Tw, found) CALL config%get(object // '.temperature', Tw, found)
IF (.NOT. found) CALL criticalError("temperature not found for wallTemperature boundary type", 'readBoundary') IF (.NOT. found) CALL criticalError("temperature not found for wallTemperature boundary type", 'readBoundary')
@ -1380,5 +1375,37 @@ MODULE moduleInput
END SUBROUTINE readParallel END SUBROUTINE readParallel
SUBROUTINE initOutput(inputFile)
USE moduleRefParam
USE moduleMesh, ONLY: mesh, doubleMesh, pathMeshParticle, pathMeshColl
USE moduleOutput, ONLY: path, folder
IMPLICIT NONE
CHARACTER(:), ALLOCATABLE, INTENT(in):: inputFile
INTEGER:: fileReference = 30
!If everything is correct, creates the output folder
CALL EXECUTE_COMMAND_LINE('mkdir ' // path // folder )
!Copies input file to output folder
CALL EXECUTE_COMMAND_LINE('cp ' // inputFile // ' ' // path // folder)
!Copies particle mesh
IF (mesh%dimen > 0) THEN
CALL EXECUTE_COMMAND_LINE('cp ' // pathMeshParticle // ' ' // path // folder)
IF (doubleMesh) THEN
CALL EXECUTE_COMMAND_LINE('cp ' // pathMeshColl // ' ' // path // folder)
END IF
END IF
! Write commit of fpakc
CALL SYSTEM('git rev-parse HEAD > ' // path // folder // '/' // 'fpakc_commit.txt')
! Write file with reference values
OPEN (fileReference, file=path // folder // '/' // 'reference.txt')
WRITE(fileReference, "(7(1X,A20))") 'L_ref', 'v_ref', 'ti_ref', 'Vol_ref', 'EF_ref', 'Volt_ref', 'B_ref'
WRITE(fileReference, "(7(1X,ES20.6E3))") L_ref, v_ref, ti_ref, Vol_ref, EF_ref, Volt_ref, B_ref
CLOSE(fileReference)
END SUBROUTINE initOutput
END MODULE moduleInput END MODULE moduleInput

View file

@ -147,7 +147,13 @@ MODULE moduleMeshBoundary
ALLOCATE(newElectron) ALLOCATE(newElectron)
ALLOCATE(newIon) ALLOCATE(newIon)
IF (ASSOCIATED(bound%electronSecondary)) THEN
newElectron%species => bound%electronSecondary
ELSE
newElectron%species => part%species newElectron%species => part%species
END IF
newIon%species => bound%species newIon%species => bound%species
newElectron%v = v0 + (1.D0 + bound%deltaV*v0/NORM2(v0)) newElectron%v = v0 + (1.D0 + bound%deltaV*v0/NORM2(v0))

View file

@ -38,6 +38,7 @@ MODULE moduleBoundary
TYPE, PUBLIC, EXTENDS(boundaryGeneric):: boundaryIonization TYPE, PUBLIC, EXTENDS(boundaryGeneric):: boundaryIonization
REAL(8):: m0, n0, v0(1:3), vTh !Properties of background neutrals. REAL(8):: m0, n0, v0(1:3), vTh !Properties of background neutrals.
CLASS(speciesGeneric), POINTER:: species !Ion species CLASS(speciesGeneric), POINTER:: species !Ion species
CLASS(speciesCharged), POINTER:: electronSecondary !Pointer to species considerer as secondary electron
TYPE(table1D):: crossSection TYPE(table1D):: crossSection
REAL(8):: effectiveTime REAL(8):: effectiveTime
REAL(8):: eThreshold REAL(8):: eThreshold
@ -103,17 +104,19 @@ MODULE moduleBoundary
END SUBROUTINE initWallTemperature END SUBROUTINE initWallTemperature
SUBROUTINE initIonization(boundary, me, m0, n0, v0, T0, speciesID, effTime, crossSection, eThreshold) SUBROUTINE initIonization(boundary, me, m0, n0, v0, T0, ion, effTime, crossSection, eThreshold, electronSecondary)
USE moduleRefParam USE moduleRefParam
USE moduleSpecies USE moduleSpecies
USE moduleCaseParam USE moduleCaseParam
USE moduleConstParam USE moduleConstParam
USE moduleErrors
IMPLICIT NONE IMPLICIT NONE
CLASS(boundaryGeneric), ALLOCATABLE, INTENT(out):: boundary CLASS(boundaryGeneric), ALLOCATABLE, INTENT(out):: boundary
REAL(8), INTENT(in):: me !Electron mass REAL(8), INTENT(in):: me !Electron mass
REAL(8), INTENT(in):: m0, n0, v0(1:3), T0 !Neutral properties REAL(8), INTENT(in):: m0, n0, v0(1:3), T0 !Neutral properties
INTEGER:: speciesID INTEGER, INTENT(in):: ion
INTEGER, OPTIONAL, INTENT(in):: electronSecondary
REAL(8):: effTime REAL(8):: effTime
CHARACTER(:), ALLOCATABLE, INTENT(in):: crossSection CHARACTER(:), ALLOCATABLE, INTENT(in):: crossSection
REAL(8), INTENT(in):: eThreshold REAL(8), INTENT(in):: eThreshold
@ -126,7 +129,22 @@ MODULE moduleBoundary
boundary%n0 = n0 * Vol_ref boundary%n0 = n0 * Vol_ref
boundary%v0 = v0 / v_ref boundary%v0 = v0 / v_ref
boundary%vTh = DSQRT(kb*T0/m0)/v_ref boundary%vTh = DSQRT(kb*T0/m0)/v_ref
boundary%species => species(speciesID)%obj boundary%species => species(ion)%obj
IF (PRESENT(electronSecondary)) THEN
SELECT TYPE(sp => species(electronSecondary)%obj)
TYPE IS(speciesCharged)
boundary%electronSecondary => sp
CLASS DEFAULT
CALL criticalError("Species " // sp%name // " chosen for " // &
"secondary electron is not a charged species", 'initIonization')
END SELECT
ELSE
boundary%electronSecondary => NULL()
END IF
boundary%effectiveTime = effTime / ti_ref boundary%effectiveTime = effTime / ti_ref
CALL boundary%crossSection%init(crossSection) CALL boundary%crossSection%init(crossSection)
CALL boundary%crossSection%convert(eV2J/(m_ref*v_ref**2), 1.D0/L_ref**2) CALL boundary%crossSection%convert(eV2J/(m_ref*v_ref**2), 1.D0/L_ref**2)

View file

@ -110,6 +110,10 @@ MODULE moduleInject
!Input current in Ampers !Input current in Ampers
self%nParticles = INT(flow*tauInject*ti_ref/(qe*species(sp)%obj%weight)) self%nParticles = INT(flow*tauInject*ti_ref/(qe*species(sp)%obj%weight))
CASE ("Am2")
!Input current in Ampers per square meter
self%nParticles = INT(flow*tauInject*ti_ref*L_ref**2/(qe*species(sp)%obj%weight))
CASE ("part/s") CASE ("part/s")
!Input current in Ampers !Input current in Ampers
self%nParticles = INT(flow*tauInject*ti_ref/species(sp)%obj%weight) self%nParticles = INT(flow*tauInject*ti_ref/species(sp)%obj%weight)

View file

@ -517,11 +517,6 @@ MODULE moduleSolver
INTEGER, INTENT(in):: t INTEGER, INTENT(in):: t
IF (t == tInitial) THEN
CALL SYSTEM('git rev-parse HEAD > ' // path // folder // '/' // 'fpack_commit.txt')
END IF
CALL outputProbes(t) CALL outputProbes(t)
counterOutput = counterOutput + 1 counterOutput = counterOutput + 1