Report Typos and Errors    
Semi-Lagrangian Library
Modular library for kinetic and gyrokinetic simulations of plasmas in fusion energy devices.
sll_m_poisson_1d_hmf.F90
Go to the documentation of this file.
1 !**************************************************************
2 ! Copyright INRIA
3 !
4 ! This code SeLaLib (for Semi-Lagrangian-Library)
5 ! is a parallel library for simulating the plasma turbulence
6 ! in a tokamak.
7 !
8 ! This software is governed by the CeCILL-B license
9 ! under French law and abiding by the rules of distribution
10 ! of free software. You can use, modify and redistribute
11 ! the software under the terms of the CeCILL-B license as
12 ! circulated by CEA, CNRS and INRIA at the following URL
13 ! "http://www.cecill.info".
14 !**************************************************************
15 
28 
29 !+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
30 #include "sll_assert.h"
31 #include "sll_memory.h"
32 #include "sll_working_precision.h"
33 
36  use sll_m_fft, only: sll_t_fft, &
42 
43  implicit none
44 
45  public :: &
50 
51  private
52 
55 
56  sll_int32 :: nc_eta1
57  sll_real64 :: eta1_min
58  sll_real64 :: eta1_max
59  sll_real64, pointer :: tmp(:)
60  sll_comp64, pointer :: rhok(:)
61  type(sll_t_fft) :: fw
62  type(sll_t_fft) :: bw
63 
64  contains
65 
66  procedure, pass(self) :: init => sll_s_poisson_1d_hmf_init
67  procedure, pass(self) :: free => sll_s_free_poisson_1d_hmf
68  procedure, pass(poisson) :: compute_phi_from_rho => compute_phi_from_rho_1d_hmf
69  procedure, pass(poisson) :: compute_e_from_rho => compute_e_from_rho_1d_hmf
70 
71  end type sll_t_poisson_1d_hmf
72 
73 contains
74 
76  subroutine sll_s_poisson_1d_hmf_init(self, eta1_min, eta1_max, nc_eta1, error)
77 
78  class(sll_t_poisson_1d_hmf), intent(out) :: self
79  sll_int32, intent(in) :: nc_eta1
80  sll_int32, intent(out) :: error
81  sll_real64, intent(in) :: eta1_min
82  sll_real64, intent(in) :: eta1_max
83 
84  error = 0
85  ! geometry
86  self%nc_eta1 = nc_eta1
87  self%eta1_min = eta1_min
88  self%eta1_max = eta1_max
89 
90  sll_allocate(self%rhok(1:nc_eta1), error)
91  sll_allocate(self%tmp(1:nc_eta1), error)
92 
93  call sll_s_fft_init_r2c_1d(self%fw, nc_eta1, self%tmp, self%rhok)
94  call sll_s_fft_init_c2r_1d(self%bw, nc_eta1, self%rhok, self%tmp)
95 
96  end subroutine sll_s_poisson_1d_hmf_init
97 
99  subroutine sll_s_solve_poisson_1d_hmf(self, field, rhs)
100 
101  class(sll_t_poisson_1d_hmf), intent(inout) :: self
102  sll_real64, dimension(:), intent(out) :: field
103  sll_real64, dimension(:), intent(in) :: rhs
104 
105  sll_assert(size(rhs) == self%nc_eta1 + 1)
106 
107  self%tmp = rhs(1:self%nc_eta1)
108  call sll_s_fft_exec_r2c_1d(self%fw, self%tmp, self%rhok)
109 
110  self%rhok(2) = self%rhok(2)*sll_p_pi*sll_p_i1 ! We keep only one mode
111  self%rhok(1) = sll_p_i0
112  self%rhok(3:) = sll_p_i0
113 
114  call sll_s_fft_exec_c2r_1d(self%bw, self%rhok, self%tmp)
115 
116  field(1:self%nc_eta1) = self%tmp/real(self%nc_eta1, f64)
117  ! complete last term by periodicity
118  field(self%nc_eta1 + 1) = field(1)
119 
120  end subroutine sll_s_solve_poisson_1d_hmf
121 
123  subroutine compute_phi_from_rho_1d_hmf(poisson, phi, rho)
124  class(sll_t_poisson_1d_hmf) :: poisson
125  sll_real64, dimension(:), intent(in) :: rho
126  sll_real64, dimension(:), intent(out) :: phi
127 
128  stop ' compute phi from rho is not implemented '
129 
130  end subroutine compute_phi_from_rho_1d_hmf
131 
133  subroutine compute_e_from_rho_1d_hmf(poisson, E, rho)
134  class(sll_t_poisson_1d_hmf) :: poisson
135  sll_real64, dimension(:), intent(in) :: rho
136  sll_real64, dimension(:), intent(out) :: e
137 
138  call sll_s_solve_poisson_1d_hmf(poisson, e, rho)
139 
140  end subroutine compute_e_from_rho_1d_hmf
141 
143  subroutine sll_s_free_poisson_1d_hmf(self, error)
144 
145  class(sll_t_poisson_1d_hmf) :: self
146  sll_int32, intent(inout) :: error
147 
148  deallocate (self%rhok, stat=error)
149  deallocate (self%tmp, stat=error)
150 
151  end subroutine sll_s_free_poisson_1d_hmf
152 
153 end module sll_m_poisson_1d_hmf
Fortran module where set some physical and mathematical constants.
real(kind=f64), parameter, public sll_p_pi
complex(kind=f64), parameter, public sll_p_i1
complex(kind=f64), parameter, public sll_p_i0
FFT interface for FFTW.
subroutine, public sll_s_fft_exec_r2c_1d(plan, array_in, array_out)
Compute fast Fourier transform in real to complex mode.
subroutine, public sll_s_fft_free(plan)
Delete a plan.
subroutine, public sll_s_fft_init_r2c_1d(plan, nx, array_in, array_out, normalized, aligned, optimization)
Create new 1d real to complex plan for forward FFT.
subroutine, public sll_s_fft_exec_c2r_1d(plan, array_in, array_out)
Compute fast Fourier transform in complex to real mode.
subroutine, public sll_s_fft_init_c2r_1d(plan, nx, array_in, array_out, normalized, aligned, optimization)
Create new 1d complex to real plan for backward FFT.
Module interface to solve Poisson equation in 1D.
Module to solve Poisson equation for the HMF model.
subroutine, public sll_s_solve_poisson_1d_hmf(self, field, rhs)
Solve the 1d equation for the Vlasov-HMF model.
subroutine, public sll_s_poisson_1d_hmf_init(self, eta1_min, eta1_max, nc_eta1, error)
Initialize the poisson 1d hmf solver.
subroutine, public sll_s_free_poisson_1d_hmf(self, error)
Free the poisson 1d hmf solver.
subroutine compute_e_from_rho_1d_hmf(poisson, E, rho)
solves
subroutine compute_phi_from_rho_1d_hmf(poisson, phi, rho)
solves
Type for FFT plan in SeLaLib.
Implementation of the poisson 1d solver for the Vlasov-HMF model.
    Report Typos and Errors