Octopus
energy_mxll.F90
Go to the documentation of this file.
1!! Copyright (C) 2020 F. Bonafe, R. Jestaedt, H. Appel, N. Tancogne-Dejean
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 debug_oct_m
23 use mesh_oct_m
24 use global_oct_m
27
28 implicit none
29
30 private
31
32 public :: &
35
36
37 type energy_mxll_t
38 ! Components are public by default
39 ! Energies
40 real(real64) :: energy = m_zero
41 real(real64) :: boundaries = m_zero
42 real(real64) :: e_energy = m_zero
43 real(real64) :: b_energy = m_zero
44 real(real64) :: energy_plane_waves = m_zero
45 real(real64) :: e_energy_plane_waves = m_zero
46 real(real64) :: b_energy_plane_waves = m_zero
47
48 real(real64) :: energy_trans = m_zero
49 real(real64) :: energy_long = m_zero
50 real(real64) :: e_energy_trans = m_zero
51 real(real64) :: b_energy_trans = m_zero
52 real(real64) :: energy_incident_waves = m_zero
53 end type energy_mxll_t
54
55contains
56
57 subroutine energy_density_calc(mesh, st, rs_field, energy_dens, e_energy_dens, b_energy_dens, plane_waves_check, &
58 rs_field_plane_waves, energy_dens_plane_waves)
59 class(mesh_t), intent(in) :: mesh
60 type(states_mxll_t), intent(in) :: st
61 complex(real64), intent(in) :: rs_field(:,:)
62 real(real64), intent(inout) :: energy_dens(:)
63 real(real64), intent(inout) :: e_energy_dens(:)
64 real(real64), intent(inout) :: b_energy_dens(:)
65 logical, optional, intent(in) :: plane_waves_check
66 complex(real64), optional, intent(in) :: rs_field_plane_waves(:,:)
67 real(real64), optional, intent(inout) :: energy_dens_plane_waves(:)
68
69 integer :: idim, ip
70
71 push_sub(energy_density_calc)
72
73 call profiling_in('ENERGY_DENSITY_CALC')
74
75 e_energy_dens(:) = m_zero
76 b_energy_dens(:) = m_zero
77 do ip = 1, mesh%np
78 do idim = 1, st%dim
79 e_energy_dens(ip) = e_energy_dens(ip) + real(rs_field(ip,idim))**2
80 b_energy_dens(ip) = b_energy_dens(ip) + aimag(rs_field(ip,idim))**2
81 end do
82 energy_dens(ip) = e_energy_dens(ip) + b_energy_dens(ip)
83 end do
84
85 if (present(rs_field_plane_waves) .and. present(energy_dens_plane_waves) .and. &
86 optional_default(plane_waves_check, .false.)) then
87 energy_dens_plane_waves(:) = m_zero
88 do ip = 1, mesh%np
89 do idim = 1, st%dim
90 energy_dens_plane_waves(ip) = energy_dens_plane_waves(ip) &
91 + real(conjg(rs_field_plane_waves(ip,idim)) * rs_field_plane_waves(ip,idim))
92 end do
93 end do
94 end if
95
96 call profiling_out('ENERGY_DENSITY_CALC')
97
98 pop_sub(energy_density_calc)
99 end subroutine energy_density_calc
100end module energy_mxll_oct_m
101
102!! Local Variables:
103!! mode: f90
104!! coding: utf-8
105!! End:
subroutine, public energy_density_calc(mesh, st, rs_field, energy_dens, e_energy_dens, b_energy_dens, plane_waves_check, rs_field_plane_waves, energy_dens_plane_waves)
real(real64), parameter, public m_zero
Definition: global.F90:187
This module defines the meshes, which are used in Octopus.
Definition: mesh.F90:118
subroutine, public profiling_out(label)
Increment out counter and sum up difference between entry and exit time.
Definition: profiling.F90:623
subroutine, public profiling_in(label, exclude)
Increment in counter and save entry time.
Definition: profiling.F90:552