Two-dimentional, special relativistic, electromagnetic particle-in-cell simulation code for general puposes in space and astrophysical plasmas.
- Solves the Vlasov-Maxwell equations by the particle-in-cell method
- Buneman-Boris method for the equation of motion of particles.
- Implicit FDTD scheme for the Maxwell equation (Hoshino, ApJ, 2013)
- Esirkepov's charge conservation scheme for the current deposit with the 2nd-order shape function (Esirkepov, CPC, 2001)
- Written in Fortran 90/95
- Hybrid parallelization by MPI and OpenMP
- 1D domain decomposition in the x direction
- SIMD optimization and efficient cache usage
- MPI-IO raw data output with JSON-based metadata
- Python scripts for HDF5 format convertor and quicklook
-
Fortran compiler
- We prepare setting files for Makefile for different compilers, including Intel Fortran, GCC-Fortran, Fujitsu compiler
- Because of the requirement of JSON-Fortran library used in this code (https://github.com/jacobwilliams/json-fortran), GCC-Fortran's version must be greater than 4.9.
-
MPI library
- We have tested with MPICH and Open MPI
- The code works with the vender's MPI library on Fujitsu's supercomputer systems
-
Python [OPTIONAL]
- The code generates raw binary data and corresponding JSON-based metadata files. A Python script in each physics directory converts them to HDF5 files as a post process
- A sample python script is prepared for quick look of the results
$ git clone git@github.com:WumingCode/WumingPIC2D.git
WumingPIC2D
├── Makefile
│
├── README.md
│
├── common
│ └── common files of PIC algorithms
│
├── common.mk
│
├── compiler-fujitsu.mk
│
├── compiler-gcc.mk
│
├── compiler-intel.mk
│
├── include
│ └── directory for module files
│
├── lib
│ └── directory for common library
│
├── proj
│ ├── shock
│ │ └── collsion-less shock simulation setup files and scripts for post process
│ └── weibel
│ └── Weibel instability simulation setup files and scripts for post process
│
├── python
│ └── json2hdf5.py - A python script to convert JSON files to HDF5 metadata
│
└── utils
└── utility files for MPI-IO and JSON output
-
Move to the installed directory.
$ cd ./WumingPIC2D
-
Copy one of
comiler-*.mk
files depending on your compiler environment tocompiler.mk
.
For instance, copycompiler-gcc.mk
if you are using gfortran.$ cp compiler-gcc.mk compiler.mk
-
Make a common library.
Compile via$ make
and make sure
libwuming*.a
are genererated in the library directorylib/
.
You are now ready for executing a specific physics problem.
Following physics problem setups are available at present.
Go to one of the physics problem directories proj/*
and make an executable main.out
.
For instance,
$ cd proj/weibel
$ make
will create an executable for Weibel instability. This will read parameters from a configuration file in JSON format. You may copy a sample configuration file config_sample.json
to config.json
:
$ cp config_sample.json config.json
and edit it as you like. By default, running the code via
$ mpiexec -np 4 ./main.out
will try to read config.json
.
If you want, you may explicitly specify the filename with a command line argument:
$ mpiexec -np 4 ./main.out somethingelse.json
in which case sometingelse.json
will be read.
Configuration parameters that are common for all physics problems are as follows.
config
verbose
Print verbose messages if >= 1.datadir
Data directory to which all the output will be saved.max_elapsed
Maximum elapsed time. The code will stop automatically when the elapsed
time goes beyond this limit. A snapshot will be saved for restart.max_it
Maximum number of iteration step.intvl_ptcl
Interval of time step for entire particle data output.intvl_mom
Interval of time step for for moment data output.intvl_orb
Interval of time step for tracer particle data output.restart_file
Snapshot file to be read during initialzation. If this is specified,
the code will start from this state. Otherwise, it will start from the
initial condition. Note that this parameter will be overwritten by the
code automatically when it finishes.
For specific problem-dependent parameters, please refer README.md
in each of physics problem directories.
The code will produce a lot of files in a specified directory. The output should
appear in pairs of *.json
and *.raw
files. A .json
file in JSON format
describes meta data and how the actual simulation data are stored in a .raw
file, which contains raw data in binary format.
The JSON files can be processed to generate HDF5 format files for data analysis
via a script json2hdf5.py
which is located in python/
directory.
For instance in the working directory,
$ python ../../python/json2hdf5.py *.json
will process all JSON files in the current directory and generate HDF5 format files for each JSON format file. You can use the generated HDF5 files for data analysis.
When the code running time goes beyond a given elapsed time limit, it will save a snapshot data and stop running. By specifying a snapshot in the configuration file, you may restart the run.
By default, the configuration file will be overwritten by the code so that you can restart the run next time with the exactly same command.
For instance, if you run the code via
$ mpiexec -np 4 ./main.out
and the elapsed time limit is reached, it will overwrite the configuration file
config.json
to properly set restart_file
option.
In the next time, you may run again via
$ mpiexec -np 4 ./main.out
then the previous snapshot data will be read automatically.
Join Slack workspace via https://join.slack.com/t/wumingpic/shared_invite/zt-xlm8cixg-NOV33dyorO1Whc4~FcVJ0g .
WumingPIC2D code uses
- JSON-Fortran API for reading/writing JSON files from Fortran.
- Amano's MPI-IO, JSON, HDF5 utitlity files
WumingPIC2D code is distributed under the MIT license.