Octopus
cube_function.F90
Go to the documentation of this file.
1!! Copyright (C) 2002-2011 M. Marques, A. Castro, A. Rubio, G. Bertsch, M. Oliveira
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 comm_oct_m
24 use cube_oct_m
25 use debug_oct_m
26 use fft_oct_m
27 use global_oct_m
28 use math_oct_m
29 use mesh_oct_m
33 use mpi_oct_m
40 use types_oct_m
41
42 implicit none
43 private
44 public :: &
66
68 ! Components are public by default
69 real(real64), pointer, contiguous :: dRS(:, :, :) => null()
70 complex(real64), pointer, contiguous :: zRS(:, :, :) => null()
71 complex(real64), pointer :: FS(:, :, :) => null()
72 logical :: forced_alloc = .false.
73 logical :: in_device_memory = .false.
74 type(accel_mem_t) :: real_space_buffer
75 type(accel_mem_t) :: fourier_space_buffer
76 end type cube_function_t
77
78contains
79
80#include "undef.F90"
81#include "real.F90"
82#include "cube_function_inc.F90"
83
84#include "undef.F90"
85#include "complex.F90"
86#include "cube_function_inc.F90"
87
88end module cube_function_oct_m
89
90
91!! Local Variables:
92!! mode: f90
93!! coding: utf-8
94!! End:
real(real64) function, public dcube_function_surface_average(cube, cf)
This function calculates the surface average of any function.
subroutine, public zmesh_to_cube(mesh, mf, cube, cf)
The next two subroutines convert a function between the normal mesh and the cube.
complex(real64) function, public zcube_function_surface_average(cube, cf)
This function calculates the surface average of any function.
subroutine, public dmesh_to_cube(mesh, mf, cube, cf)
The next two subroutines convert a function between the normal mesh and the cube.
subroutine, public dcube_to_submesh(cube, cf, sm, mf)
subroutine, public dcube_to_mesh(cube, cf, mesh, mf)
subroutine, public dcube_function_allgather(cube, cf, cf_local, order, gatherfs)
subroutine, public zcube_to_mesh_parallel(cube, cf, mesh, mf, map)
subroutine, public zcube_to_mesh(cube, cf, mesh, mf)
subroutine, public dcube_function_alloc_rs(cube, cf, in_device, force_alloc)
Allocates locally the real space grid, if PFFT library is not used. Otherwise, it assigns the PFFT re...
subroutine, public zcube_to_submesh(cube, cf, sm, mf)
subroutine, public zmesh_to_cube_parallel(mesh, mf, cube, cf, map)
The next two subroutines convert a function between the normal mesh and the cube in parallel.
subroutine, public dcube_function_free_rs(cube, cf)
Deallocates the real space grid.
subroutine, public dsubmesh_to_cube(sm, mf, cube, cf)
The next two subroutines convert a function between a submesh and the cube.
subroutine, public zsubmesh_to_cube(sm, mf, cube, cf)
The next two subroutines convert a function between a submesh and the cube.
subroutine, public zcube_function_allgather(cube, cf, cf_local, order, gatherfs)
subroutine, public dmesh_to_cube_parallel(mesh, mf, cube, cf, map)
The next two subroutines convert a function between the normal mesh and the cube in parallel.
subroutine, public zcube_function_free_rs(cube, cf)
Deallocates the real space grid.
subroutine, public dcube_to_mesh_parallel(cube, cf, mesh, mf, map)
subroutine, public zcube_function_alloc_rs(cube, cf, in_device, force_alloc)
Allocates locally the real space grid, if PFFT library is not used. Otherwise, it assigns the PFFT re...
Fast Fourier Transform module. This module provides a single interface that works with different FFT ...
Definition: fft.F90:118
This module is intended to contain "only mathematical" functions and procedures.
Definition: math.F90:115
This module defines the meshes, which are used in Octopus.
Definition: mesh.F90:118
Some general things and nomenclature:
Definition: par_vec.F90:171
This module is an helper to perform ring-pattern communications among all states.