Report Typos and Errors    
Semi-Lagrangian Library
Modular library for kinetic and gyrokinetic simulations of plasmas in fusion energy devices.
sll_m_pic_visu_parallel.F90
Go to the documentation of this file.
1 !**************************************************************
2 ! Copyright INRIA
3 ! Authors :
4 ! CALVI project team
5 !
6 ! This code SeLaLib (for Semi-Lagrangian-Library)
7 ! is a parallel library for simulating the plasma turbulence
8 ! in a tokamak.
9 !
10 ! This software is governed by the CeCILL-B license
11 ! under French law and abiding by the rules of distribution
12 ! of free software. You can use, modify and redistribute
13 ! the software under the terms of the CeCILL-B license as
14 ! circulated by CEA, CNRS and INRIA at the following URL
15 ! "http://www.cecill.info".
16 !**************************************************************
17 
21 
23 
24 !+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
25 #include "sll_working_precision.h"
26 
27  use sll_m_collective, only: &
31 
32  use sll_m_pic_visu, only: &
34 
35  use sll_m_utilities, only: &
37 
38  use sll_m_xdmf, only: &
40 
41  use mpi, only: &
42  mpi_sum
43 
44  implicit none
45 
46  public :: &
48 
49  private
50 !+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
51 
52 contains
53 
56  subroutine sll_s_distribution_xdmf_coll(plot_name, x, v, w, &
57  xmin, xmax, nx, vmin, vmax, nv, iplot, collective, root_rank)
58 
59  character(len=*), intent(in) :: plot_name
60  sll_real64, dimension(:), intent(in) :: x
61  sll_real64, dimension(:), intent(in) :: v
62  sll_real64, dimension(:), intent(in) :: w
63  sll_real64, intent(in) :: xmin
64  sll_real64, intent(in) :: xmax
65  sll_int32, intent(in) :: nx
66  sll_real64, intent(in) :: vmin
67  sll_real64, intent(in) :: vmax
68  sll_int32, intent(in) :: nv
69  sll_int32, intent(in) :: iplot
70  type(sll_t_collective_t), pointer :: collective
71  sll_int32, optional, intent(in) :: root_rank
72 
73  sll_real64, dimension(nx, nv) :: df_local
74  sll_real64, dimension(nx, nv) :: df
75  sll_real64 :: delta_x
76  sll_real64 :: delta_v
77  character(len=4) :: fin
78 
79  call sll_s_int2string(iplot, fin)
80 
81  delta_x = (xmax - xmin)/real(nx - 1, f64)
82  delta_v = (vmax - vmin)/real(nv - 1, f64)
83 
84  call sll_s_compute_df_cic(x, v, w, xmin, xmax, nx, vmin, vmax, nv, df_local)
85 
86  if (present(root_rank)) then
87 
88  !Write the result only to node with root_rank
89  call sll_o_collective_reduce(collective, df_local, nx*nv, mpi_sum, root_rank, df)
90 
91  if (sll_f_get_collective_rank(collective) == root_rank) then
92  call sll_s_xdmf_corect2d_nodes(plot_name//'_'//fin, df, "density", xmin, delta_x, vmin, delta_v)
93  end if
94 
95  else
96 
97  !Write the result to all nodes, and use parallel io
98 
99  end if
100 
101  end subroutine sll_s_distribution_xdmf_coll
102 
103 end module sll_m_pic_visu_parallel
Reduces values on all processes to a single value.
Parallelizing facility.
integer(kind=i32) function, public sll_f_get_collective_rank(col)
Determines the rank of the calling process in the communicator.
This module provides some routines for plotting during PIC simulations with MPI. It extends sll_m_pic...
subroutine, public sll_s_distribution_xdmf_coll(plot_name, x, v, w, xmin, xmax, nx, vmin, vmax, nv, iplot, collective, root_rank)
VisIt readable output for particles density Data file format could be XML, HDF5 or Binary (not fully ...
This module provides some routines for plotting during PIC simulations.
subroutine, public sll_s_compute_df_cic(xp, yp, wp, xmin, xmax, nx, ymin, ymax, ny, df)
Compute grid field from particles distribution with the CIC scheme (Cloud In.
Some common numerical utilities.
subroutine, public sll_s_int2string(istep, cstep)
Convert an integer < 9999 to a 4 characters string.
Implements the functions to write xdmf file plotable by VisIt.
Definition: sll_m_xdmf.F90:24
subroutine, public sll_s_xdmf_corect2d_nodes(file_name, array, array_name, eta1_min, delta_eta1, eta2_min, delta_eta2, file_format, iplot, time)
Subroutine to write a 2D array in xdmf format The field is describe on a cartesian mesh Axis are perp...
Definition: sll_m_xdmf.F90:265
Wrapper around the communicator.
    Report Typos and Errors