The Open FUSION Toolkit 1.0.0-8905cc5
Modeling tools for plasma and fusion research and engineering
|
In this example we demonstrate how to compute eigenvalues/eigenvectors and run a time-domain simulation for a large ThinCurr model using HODLR matrix compression.
pip
or other standard methods.To load the ThinCurr python module we need to tell python where to the module is located. This can be done either through the PYTHONPATH
environment variable or within a script using sys.path.append()
as below, where we look for the environement variable OFT_ROOTPATH
to provide the path to where the OpenFUSIONToolkit is installed (/Applications/OFT
for binaries on macOS).
We now create a OFT_env instance for execution using four threads and a ThinCurr instance that utilizes that execution environment.
Once created, we setup the model from an existing HDF5 and XML mesh definition using setup_model().
We also initialize I/O for this model using setup_io() to enable output of plotting files for 3D visualization in VisIt, Paraview, or using pyvista below.
#---------------------------------------------- Open FUSION Toolkit Initialized Development branch: v1_beta6 Revision id: 681e857 Parallelization Info: # of MPI tasks = 1 # of NUMA nodes = 1 # of OpenMP threads = 4 Fortran input file = /var/folders/52/n5qxh27n4w19qxzqygz2btbw0000gn/T/oft_64966/oftpyin XML input file = none Integer Precisions = 4 8 Float Precisions = 4 8 16 Complex Precisions = 4 8 LA backend = native #---------------------------------------------- Creating thin-wall model Orientation depth = 26500 Loading V(t) driver coils Loading I(t) driver coils # of points = 22580 # of edges = 67150 # of cells = 44560 # of holes = 11 # of Vcoils = 0 # of closures = 0 # of Icoils = 1 Building holes Loading region resistivity: 1 1.2570E-05
To visualize the results we define three Mirnov sensors spanning one port period, just outside the torus. While every sensor in ThinCurr is a flux loop, helper classes for defining specific types of common sensors (eg. circular loop or Mirnovs) are available in OpenFUSIONToolkit.ThinCurr.sensor, which also contains the save_sensors() function for saving this information to a ThinCurr-compatible file format.
After defining the sensors we use compute_Msensor() to setup the sensors and compute mutual matrices between the sensors and the model (Msensor
) and the sensors and Icoils (Msc
).
Loading floop information: # of floops = 3 Building element->sensor inductance matrix Building coil->sensor inductance matrix
With the model setup, we can now compute the self-inductance matrix using HODLR. When HODLR is used the result is a pointer to the Fortran operator, which is stored at tw_torus.Lmat_hodlr. As in any other case, by default, the resistivity matrix is not moved to python as it is sparse and converting to dense representation would require an increase in memory. These matrices correspond to the \(\textrm{L}\) and \(\textrm{R}\) matrices for the physical system
\(\textrm{L} \frac{\partial I}{\partial t} + \textrm{R} I = V\)
Building coil<->element inductance matrices Time = 0s Building coil<->coil inductance matrix
Partitioning grid for block low rank compressed operators nBlocks = 32 Avg block size = 686 # of SVD = 167 # of ACA = 161 Building block low rank inductance operator Building hole and Vcoil columns Building diagonal blocks 10% 20% 30% 40% 50% 60% 70% 80% 90% Building off-diagonal blocks using ACA+ 10% 20% 30% 40% 50% 60% 70% 80% 90% Compression ratio: 6.1% ( 2.94E+07/ 4.82E+08) Time = 1m 17s Saving HODLR matrix to file: HOLDR_L.save
Building resistivity matrix
With \(\textrm{L}\) and \(\textrm{R}\) matrices we can now compute the eigenvalues and eigenvectors of the system \(\textrm{L} I = \lambda \textrm{R} I\), where the eigenvalues \(\lambda = \tau_{L/R}\) are the decay time-constants of the current distribution corresponding to each eigenvector.
Starting eigenvalue solve Time = 3.2132510000000001 Eigenvalues 3.8803529812470830E-002 2.3122405888831986E-002 2.3121319178033790E-002 2.1444708809391269E-002 2.0605567088517945E-002
The resulting currents can be saved for plotting using tw_torus.save_current(). Here we save each of the five eigenvectors for visualization. Once all fields have been saved for plotting tw_torus.build_XDMF() to generate the XDMF descriptor files for plotting with VisIt of Paraview. This method also returns a XDMF_plot_file object, which can be used to read and interact with plot data in Python (see below).
Removing old Xdmf files Removed 23 files Creating output files: oft_xdmf.XXXX.h5 Found Group: thincurr Found Mesh: icoils # of blocks: 1 Found Mesh: smesh # of blocks: 1
Now we plot the current vectors on the plate showing the longest-lived eddy current structure, which corresponds to a large circulation on the plate. The XDMF_plot_file class provides functionality to work with the data stored in OFT plot files, including methods to generate information for 3D plotting in Python using pyvista.
Plotting data is always associated with a specific mesh, which is itself associated with a physics group. In this case ThinCurr
is the physics group and the data we are interested in is stored on the surface mesh smesh
. The ‘plot_data['ThinCurr’]['smesh']` is a XDMF_plot_mesh object with further functionality for accessing data.
To plot the first eigenvector we first get the pyvista plotting grid using get_pyvista_grid() and then retrieve the vertex-centered field (J_01_v
) using get_field()
With the model fully defined we can now use run_td() to perform a time-domain simulation. In this case we simulate 40 ms using a timestep of 0.2 ms (200 steps). We also specify the current in the single I-coil defined in the XML input file as a function of time (coil_currs
), where the first column specifies time points in ascending order and the remaining columns specify coil currents at each time point.
Starting simulation Timestep 10 2.00000009E-03 22.7324562 32 Timestep 20 4.00000019E-03 44.1192398 18 Timestep 30 6.00000005E-03 42.2026978 32 Timestep 40 8.00000038E-03 39.9049339 32 Timestep 50 9.99999978E-03 37.7692719 32 Timestep 60 1.20000001E-02 35.7698593 32 Timestep 70 1.40000004E-02 33.8907852 32 Timestep 80 1.60000008E-02 32.1206093 32 Timestep 90 1.79999992E-02 30.4503765 32 Timestep 100 1.99999996E-02 28.8726788 32 Timestep 110 2.19999999E-02 27.3811455 32 Timestep 120 2.40000002E-02 25.9701900 32 Timestep 130 2.60000005E-02 24.6347961 32 Timestep 140 2.80000009E-02 23.3704300 32 Timestep 150 2.99999993E-02 22.1729183 32 Timestep 160 3.20000015E-02 21.0384312 32 Timestep 170 3.40000018E-02 19.9634094 32 Timestep 180 3.59999985E-02 18.9445438 32 Timestep 190 3.79999988E-02 17.9787312 32 Timestep 200 3.99999991E-02 17.0630760 32
Building block low rank magnetic field operator Building hole and Vcoil columns Building diagonal blocks 10% 20% 30% 40% 50% 60% 70% 80% 90% Building off-diagonal blocks using ACA+ 10% 20% 30% 40% 50% 60% 70% 80% 90% Compression ratio: 7.1% ( 1.06E+08/ 1.49E+09) Time = 3m 4s Saving HODLR matrix from file: HODLR_B.save
Post-processing simulation Removing old Xdmf files Removed 2 files Creating output files: oft_xdmf.XXXX.h5 Found Group: thincurr Found Mesh: icoils # of blocks: 1 Found Mesh: smesh # of blocks: 1
OFT History file: floops.hist Number of fields = 4 Number of entries = 201 Fields: time: Simulation time [s] (d1) Bz_1: No description (d1) Bz_2: No description (d1) Bz_3: No description (d1)