Octopus
fourier_space.F90
Go to the documentation of this file.
1!! Copyright (C) 2002-2006 M. Marques, A. Castro, A. Rubio, G. Bertsch
2!!
3!! This program is free software; you can redistribute it and/or modify
4!! it under the terms of the GNU General Public License as published by
5!! the Free Software Foundation; either version 2, or (at your option)
6!! any later version.
7!!
8!! This program is distributed in the hope that it will be useful,
9!! but WITHOUT ANY WARRANTY; without even the implied warranty of
10!! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
11!! GNU General Public License for more details.
12!!
13!! You should have received a copy of the GNU General Public License
14!! along with this program; if not, write to the Free Software
15!! Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA
16!! 02110-1301, USA.
17!!
18
19#include "global.h"
20
22 use accel_oct_m
23 use cube_oct_m
25 use debug_oct_m
26 use global_oct_m
27 use, intrinsic :: iso_fortran_env
28 use math_oct_m
31 use fft_oct_m
32#ifdef HAVE_OPENMP
33 use omp_lib
34#endif
35#ifdef HAVE_PFFT
36 use pfft_oct_m
37#endif
39 use types_oct_m
40
41 implicit none
42 private
43 public :: &
56
58 private
59 real(real64), allocatable :: dop(:, :, :)
60 complex(real64), allocatable :: zop(:, :, :)
61 logical :: in_device_memory = .false.
62 type(accel_mem_t) :: op_buffer
63 logical :: real_op
64
65 !Parameters used to generate this kernel
66 real(real64), allocatable, public :: qq(:)
67 real(real64), public :: singularity = m_zero
68 real(real64), public :: mu = m_zero
69 real(real64), public :: alpha = m_zero
70 real(real64), public :: beta = m_zero
71 end type fourier_space_op_t
72
73contains
74
75 ! ---------------------------------------------------------
79 subroutine cube_function_alloc_fs(cube, cf, force_alloc)
80 type(cube_t), target, intent(in) :: cube
81 type(cube_function_t), intent(inout) :: cf
82 logical, optional, intent(in) :: force_alloc
83
84 integer :: n1, n2, n3
85 logical :: is_allocated
86
88
89 assert(.not. associated(cf%fs))
90 assert(allocated(cube%fft))
91
92 cf%forced_alloc = optional_default(force_alloc, .false.)
93
94 n1 = max(1, cube%fs_n(1))
95 n2 = max(1, cube%fs_n(2))
96 n3 = max(1, cube%fs_n(3))
97
98 is_allocated = .false.
99
100 select case (cube%fft%library)
101 case (fftlib_pfft)
102 if (.not. cf%forced_alloc) then
103 is_allocated = .true.
104 if (any(cube%fs_n(1:3) == 0)) then
105 cf%fs => cube%fft%fs_data(1:1,1:1,1:1)
106 else
107 cf%fs => cube%fft%fs_data(1:n3,1:n1,1:n2)
108 end if
109 else ! force allocate transposed with PFFT
110 is_allocated = .true.
111 safe_allocate(cf%fs(1:n3, 1:n1, 1:n2))
112 end if
113 case (fftlib_accel)
114 if (cf%in_device_memory) then
115 is_allocated = .true.
116 call accel_create_buffer(cf%fourier_space_buffer, accel_mem_read_write, type_cmplx, product(cube%fs_n(1:3)))
117 end if
118
119 case (fftlib_fftw)
120 if (.not. cf%forced_alloc) then
121 is_allocated = .true.
122 cf%fs => cube%fft%fs_data(1:cube%fs_n(1), 1:cube%fs_n(2), 1:cube%fs_n(3))
123 end if
124 end select
125
126 if (.not. is_allocated) then
127 safe_allocate(cf%fs(1:cube%fs_n(1), 1:cube%fs_n(2), 1:cube%fs_n(3)))
128 end if
129
131 end subroutine cube_function_alloc_fs
132
133
134 ! ---------------------------------------------------------
136 subroutine cube_function_free_fs(cube, cf)
137 type(cube_t), intent(in) :: cube
138 type(cube_function_t), intent(inout) :: cf
139
140 logical :: deallocated
141
142 push_sub(cube_function_free_fs)
143
144 assert(allocated(cube%fft))
145
146 deallocated = .false.
147
148 select case (cube%fft%library)
149 case (fftlib_pfft)
150 if (.not. cf%forced_alloc) then
151 deallocated = .true.
152 nullify(cf%fs)
153 end if
155 if (cf%in_device_memory) then
156 deallocated = .true.
157 call accel_release_buffer(cf%fourier_space_buffer)
158 end if
160 if (.not. cf%forced_alloc) then
161 deallocated = .true.
162 nullify(cf%fs)
163 end if
164 end select
165
166 if (.not. deallocated) then
167 assert(associated(cf%fs))
168 safe_deallocate_p(cf%fs)
169 end if
170
171 pop_sub(cube_function_free_fs)
173
174
175 ! ---------------------------------------------------------
176 subroutine fourier_space_op_end(this)
177 type(fourier_space_op_t), intent(inout) :: this
178
179 push_sub(fourier_space_op_end)
180
181 if (this%in_device_memory) then
182 call accel_release_buffer(this%op_buffer)
183 this%in_device_memory = .false.
184 end if
185 safe_deallocate_a(this%dop)
186 safe_deallocate_a(this%zop)
187 safe_deallocate_a(this%qq)
188
189 pop_sub(fourier_space_op_end)
190 end subroutine fourier_space_op_end
191
192#include "undef.F90"
193#include "real.F90"
194#include "fourier_space_inc.F90"
195
196#include "undef.F90"
197#include "complex.F90"
198#include "fourier_space_inc.F90"
199
200end module fourier_space_oct_m
201
202
203!! Local Variables:
204!! mode: f90
205!! coding: utf-8
206!! End:
integer, parameter, public accel_mem_read_write
Definition: accel.F90:183
subroutine, public accel_release_buffer(this)
Definition: accel.F90:1246
Fast Fourier Transform module. This module provides a single interface that works with different FFT ...
Definition: fft.F90:118
integer, parameter, public fftlib_accel
Definition: fft.F90:182
integer, parameter, public fftlib_pfft
Definition: fft.F90:182
integer, parameter, public fftlib_fftw
Definition: fft.F90:182
subroutine, public zfourier_space_op_apply(this, cube, cf)
Applies a multiplication factor to the Fourier space grid. This is a local function.
subroutine, public dfourier_space_op_apply(this, cube, cf)
Applies a multiplication factor to the Fourier space grid. This is a local function.
subroutine, public fourier_space_op_end(this)
subroutine, public dfourier_space_op_init(this, cube, op, in_device)
subroutine, public zcube_function_fs2rs(cube, cf)
subroutine, public dcube_function_fs2rs(cube, cf)
subroutine, public cube_function_free_fs(cube, cf)
Deallocates the Fourier space grid.
subroutine, public zfourier_space_op_init(this, cube, op, in_device)
subroutine, public zcube_function_rs2fs(cube, cf)
The following routines convert the function between real space and Fourier space Note that the dimens...
subroutine, public cube_function_alloc_fs(cube, cf, force_alloc)
Allocates locally the Fourier space grid, if PFFT library is not used. Otherwise, it assigns the PFFT...
subroutine, public dcube_function_rs2fs(cube, cf)
The following routines convert the function between real space and Fourier space Note that the dimens...
real(real64), parameter, public m_zero
Definition: global.F90:187
This module is intended to contain "only mathematical" functions and procedures.
Definition: math.F90:115
The low level module to work with the PFFT library. http:
Definition: pfft.F90:164
type(type_t), public type_cmplx
Definition: types.F90:134
int true(void)