Octopus
nlcc.F90
Go to the documentation of this file.
1!! Copyright (C) 2002-2006 M. Marques, A. Castro, A. Rubio, G. Bertsch
2!! Copyright (C) 2021 Nicolas Tancogne-Dejean
3!!
4!! This program is free software; you can redistribute it and/or modify
5!! it under the terms of the GNU General Public License as published by
6!! the Free Software Foundation; either version 2, or (at your option)
7!! any later version.
8!!
9!! This program is distributed in the hope that it will be useful,
10!! but WITHOUT ANY WARRANTY; without even the implied warranty of
11!! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
12!! GNU General Public License for more details.
13!!
14!! You should have received a copy of the GNU General Public License
15!! along with this program; if not, write to the Free Software
16!! Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA
17!! 02110-1301, USA.
18!!
19
20#include "global.h"
21
22module nlcc_oct_m
23 use atom_oct_m
24 use comm_oct_m
25 use debug_oct_m
28 use global_oct_m
31 use ions_oct_m
34 use mesh_oct_m
39 use ps_oct_m
42 use space_oct_m
47
48 implicit none
49
50 private
51 public :: &
52 nlcc_t
53
54 type, extends(density_interaction_t) :: nlcc_t
55 private
56
57 type(mesh_t), pointer :: mesh => null()
58 type(space_t), pointer :: space => null()
59
60 ! This information should be copied, and not obtained thanks to pointer
61 ! This is a temporary change here
62 type(distributed_t), pointer, public :: atoms_dist => null()
63 type(atom_t), pointer, public :: atom(:) => null()
64 real(real64), pointer, public :: pos(:,:) => null()
65 type(lattice_vectors_t), pointer :: latt => null()
66
67 contains
68 procedure :: init => nlcc_init
69 procedure :: bind => nlcc_bind
70 procedure :: calculate => nlcc_calculate
71 procedure :: calculate_energy => nlcc_calculate_energy
72 procedure :: end => nlcc_end
73 final :: nlcc_finalize
74 end type nlcc_t
75
76 interface nlcc_t
77 module procedure nlcc_constructor
78 end interface nlcc_t
79
80contains
81
82 ! ---------------------------------------------------------
83 function nlcc_constructor(partner) result(this)
84 class(interaction_partner_t), target, intent(inout) :: partner
85 class(nlcc_t), pointer :: this
86
87 push_sub(nlcc_constructor)
88
89 allocate(this)
90
91 this%label = "nlcc"
92
93 this%partner => partner
94
95 pop_sub(nlcc_constructor)
96 end function nlcc_constructor
97
98 ! ---------------------------------------------------------
99 subroutine nlcc_init(this, mesh, ions)
100 class(nlcc_t), intent(inout) :: this
101 class(mesh_t), target, intent(in) :: mesh
102 type(ions_t), target, intent(in) :: ions
103
104
105 push_sub(nlcc_init)
106
107 this%mesh => mesh
108
109 safe_allocate(this%density(1:mesh%np,1))
110
111 this%atoms_dist => ions%atoms_dist
112 this%atom => ions%atom
113 this%space => ions%space
114 this%pos => ions%pos
115 this%latt => ions%latt
116
117 pop_sub(nlcc_init)
118 end subroutine nlcc_init
119
120 ! ---------------------------------------------------------
123 subroutine nlcc_bind(this, ions)
124 class(nlcc_t), intent(inout) :: this
125 type(ions_t), target, intent(in) :: ions
126
127 push_sub(nlcc_bind)
128
129 this%atoms_dist => ions%atoms_dist
130 this%atom => ions%atom
131 this%space => ions%space
132 this%pos => ions%pos
133 this%latt => ions%latt
134
135 pop_sub(nlcc_bind)
136 end subroutine nlcc_bind
137
138 ! ---------------------------------------------------------
139 subroutine nlcc_calculate(this)
140 class(nlcc_t), intent(inout) :: this
141
142 integer :: ia
143
144 push_sub(nlcc_calculate)
145
146 call profiling_in("NLCC_INTER")
147
148 this%density(:,:) = m_zero
150 do ia = this%atoms_dist%start, this%atoms_dist%end
151 if (this%atom(ia)%species%is_ps_with_nlcc()) then
152 call species_get_nlcc(this%atom(ia)%species, this%space, this%latt, this%pos(:, ia), this%mesh, &
153 this%density(:,1), accumulate=.true.)
154 end if
155 end do
156
157 ! reduce over atoms if required
158 if (this%atoms_dist%parallel) then
159 call comm_allreduce(this%atoms_dist%mpi_grp, this%density(:,1))
160 end if
161
162 call profiling_out("NLCC_INTER")
165 end subroutine nlcc_calculate
167 ! ---------------------------------------------------------
168 subroutine nlcc_end(this)
169 class(nlcc_t), intent(inout) :: this
170
171 push_sub(nlcc_end)
172
173 safe_deallocate_a(this%density)
174 nullify(this%mesh)
175
176 call interaction_end(this)
177
178 pop_sub(nlcc_end)
179 end subroutine nlcc_end
180
181
182 ! ---------------------------------------------------------
183 subroutine nlcc_finalize(this)
184 type(nlcc_t), intent(inout) :: this
185
186 push_sub(nlcc_finalize)
187
188 call nlcc_end(this)
189
190 pop_sub(nlcc_finalize)
191 end subroutine nlcc_finalize
192
193 ! ---------------------------------------------------------
194 subroutine nlcc_calculate_energy(this)
195 class(nlcc_t), intent(inout) :: this
196
197 push_sub(nlcc_calculate_energy)
198
199 !NLCC has no energy by itself. The energy is obtained through exchange and correlation
200
201 this%energy = m_zero
202
203 pop_sub(nlcc_calculate_energy)
204 end subroutine nlcc_calculate_energy
205
206end module nlcc_oct_m
207
208!! Local Variables:
209!! mode: f90
210!! coding: utf-8
211!! End:
real(real64), parameter, public m_zero
Definition: global.F90:191
This module defines the abstract interaction_t class, and some auxiliary classes for interactions.
subroutine, public interaction_end(this)
This module defines classes and functions for interaction partners.
This module defines the meshes, which are used in Octopus.
Definition: mesh.F90:120
subroutine nlcc_end(this)
Definition: nlcc.F90:264
subroutine nlcc_calculate(this)
Definition: nlcc.F90:235
subroutine nlcc_finalize(this)
Definition: nlcc.F90:279
subroutine nlcc_init(this, mesh, ions)
Definition: nlcc.F90:195
class(nlcc_t) function, pointer nlcc_constructor(partner)
Definition: nlcc.F90:179
subroutine nlcc_calculate_energy(this)
Definition: nlcc.F90:290
subroutine nlcc_bind(this, ions)
Rebind NLCC helper pointers after copying the ionic geometry. The movers in NLCC reference ionic arra...
Definition: nlcc.F90:219
subroutine, public profiling_out(label)
Increment out counter and sum up difference between entry and exit time.
Definition: profiling.F90:631
subroutine, public profiling_in(label, exclude)
Increment in counter and save entry time.
Definition: profiling.F90:554
Definition: ps.F90:116
This module defines the quantity_t class and the IDs for quantities, which can be exposed by a system...
Definition: quantity.F90:140
subroutine, public species_get_nlcc(species, space, latt, pos, mesh, rho_core, accumulate)
int true(void)