Octopus
current_to_mxll_field.F90
Go to the documentation of this file.
1!! Copyright (C) 2022 F. Bonafé
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#include "global.h"
19
21 use debug_oct_m
23 use global_oct_m
24 use grid_oct_m
26 use, intrinsic :: iso_fortran_env
28 use mesh_oct_m
31 use parser_oct_m
35 use space_oct_m
38 use unit_oct_m
41
42 implicit none
43
44 private
45 integer, parameter :: COS2 = 1
46 public :: &
48
57 private
58 logical, public :: grid_based_partner = .true.
59 complex(real64), allocatable :: rs_current(:, :)
60
62 integer, public :: partner_np = 0
63 real(real64), allocatable, public :: partner_charge(:)
64 real(real64), allocatable, public :: partner_pos(:,:)
65 real(real64), allocatable, public :: partner_vel(:,:)
66
67 integer :: reg_type
68 real(real64) :: reg_width
69 real(real64), allocatable :: reg(:)
70
71 type(space_t), pointer :: system_space => null()
72 type(lattice_vectors_t), pointer :: system_latt => null()
73
74
75 contains
76 procedure :: init => current_to_mxll_field_init
77
78 procedure :: init_space_latt => current_to_mxll_field_init_space_latt
79
80 procedure :: do_mapping => current_to_mxll_field_do_mapping
81
82 procedure :: calculate_energy => current_to_mxll_field_calculate_energy
83
86
87
90 end interface current_to_mxll_field_t
91
92contains
93
94 subroutine current_to_mxll_field_init(this, gr, ndim)
95 class(current_to_mxll_field_t), intent(inout) :: this
96 type(grid_t), target, intent(in) :: gr
97 integer, intent(in) :: ndim
98
100
101 call field_transfer_init(this, gr, ndim)
102 safe_allocate(this%rs_current(1:gr%np, 1:ndim))
103
105 end subroutine current_to_mxll_field_init
106
107 ! ---------------------------------------------------------
108 function current_to_mxll_field_constructor(partner) result(this)
109 class(interaction_partner_t), target, intent(inout) :: partner
110 class(current_to_mxll_field_t), pointer :: this
111
114 allocate(this)
115
116 this%label = "current_to_mxll_field"
117 this%partner => partner
118
119 ! TODO: For point-wise systems, this should be different, as we
120 ! need there the charge and velocity to do the charge regularization.
121 ! Once each partner has a basis class, we can do a derived type
122 ! TODO: look for user input option of the regularization function for point partner
123 this%couplings_from_partner = [current]
124
125 !%Variable RegularizationFunction
126 !%Type integer
127 !%Default COS2
128 !%Section Maxwell
129 !%Description
130 !% The current arising from charged point particles must be mapped onto the Maxwell
131 !% propagation grid. This requires a smearing or regularization function $\phi(\mathbf{r})$ attached to
132 !% each particle position $\mathbf{r}_i$ with user defined cutoff width, $\sigma$
133 !%Option COS2 1
134 !% $\phi(r)=\text{cos}^2(\frac{\pi}{2}\frac{|\mathbf{r}-\mathbf{r}_i|}{\sigma})$
135 !% if $|\mahtbf{r}-\mathbf{r}_i|<\sigma$, and 0 otherwise.
136 !%End
137
138 call parse_variable(partner%namespace, 'RegularizationFunction', cos2, this%reg_type)
139 if (.not. varinfo_valid_option('RegularizationFunction', this%reg_type)) &
140 call messages_input_error(partner%namespace, 'RegularizationFunction')
141
142 !%Variable RegularizationFunctionWidth
143 !%Type float
144 !%Default 2
145 !%Section Maxwell
146 !%Description
147 !% The current arising from charged point particles must be mapped onto the Maxwell
148 !% propagation grid. This requires a smearing or regularization function $\phi(\mathbf{r})$ attached to
149 !% each particle position $\mathbf{r}_i$ with user defined cutoff width, $\sigma$.
150 !% Default 2 bohrradii
151 !%End
153 call parse_variable(partner%namespace, 'RegularizationFunctionWidth', 2.0_real64, &
154 this%reg_width, units_inp%length)
156 this%intra_interaction = .false.
162 subroutine current_to_mxll_field_init_space_latt(this, space, latt)
163 class(current_to_mxll_field_t), intent(inout) :: this
164 class(space_t), target, intent(in) :: space
165 type(lattice_vectors_t), target, intent(in) :: latt
166
168
169 this%system_space => space
170 this%system_latt => latt
174
175 ! ---------------------------------------------------------
176 subroutine current_to_mxll_field_finalize(this)
177 type(current_to_mxll_field_t), intent(inout) :: this
178
180
181 call this%end()
182 nullify(this%system_space)
183 nullify(this%system_latt)
184 safe_deallocate_a(this%rs_current)
185
188
189 ! ---------------------------------------------------------
190 subroutine current_to_mxll_field_do_mapping(this)
191 class(current_to_mxll_field_t), intent(inout) :: this
192
193 integer :: part_ind, i_dim, ip
194 type(submesh_t) :: submesh
195 real(real64) :: norm, time
196
197 push_sub_with_profile(current_to_mxll_field_do_mapping)
198
199 if(this%grid_based_partner) then
200 call this%regridding%do_transfer(this%system_field, this%partner_field)
201 else
202 this%system_field = m_zero
203 do part_ind = 1, this%partner_np
204 call submesh_init(submesh, this%system_space, this%system_gr, &
205 this%system_latt, this%partner_pos(:,part_ind), this%reg_width)
206
207 safe_allocate(this%reg(1:submesh%np))
208 if (this%reg_type == cos2) then
209 do ip = 1, submesh%np
210 this%reg(ip) = cos( (m_pi/m_two) * (submesh%r(ip)/this%reg_width) )**2
211 end do
212 end if
213
214 norm = dsm_integrate(this%system_gr, submesh, this%reg)
215
216 do i_dim = 1, this%system_space%dim
217 call submesh_add_to_mesh(submesh, this%reg, this%system_field(:,i_dim), &
218 this%partner_charge(part_ind)*this%partner_vel(i_dim,part_ind)/norm)
219 end do
220 safe_deallocate_a(this%reg)
221
222 call submesh_end(submesh)
223 end do
224
225 end if
226 ! add the RS current to the interpolation
227 call build_rs_current_state(this%system_field, this%system_gr, this%rs_current)
228 assert(this%interpolation_initialized)
229 time = this%partner%quantities(this%couplings_from_partner(1))%iteration%value()
230 call this%interpolation%add_time(time, this%rs_current)
231
232 pop_sub_with_profile(current_to_mxll_field_do_mapping)
234
235 ! ---------------------------------------------------------
237 class(current_to_mxll_field_t), intent(inout) :: this
238
240
241 ! interaction energy is zero, since it is only re-gridding the quantities of one system
242 ! on the mesh of the other
243 this%energy = m_zero
244
247
249
250!! Local Variables:
251!! mode: f90
252!! coding: utf-8
253!! End:
double cos(double __x) __attribute__((__nothrow__
subroutine current_to_mxll_field_init_space_latt(this, space, latt)
subroutine current_to_mxll_field_finalize(this)
class(current_to_mxll_field_t) function, pointer current_to_mxll_field_constructor(partner)
subroutine current_to_mxll_field_do_mapping(this)
subroutine current_to_mxll_field_init(this, gr, ndim)
subroutine current_to_mxll_field_calculate_energy(this)
This module implements the field transfer.
subroutine, public field_transfer_init(this, gr, ndim)
the system field is allocated and initialized to 0
real(real64), parameter, public m_two
Definition: global.F90:189
real(real64), parameter, public m_zero
Definition: global.F90:187
real(real64), parameter, public m_pi
some mathematical constants
Definition: global.F90:185
This module implements the underlying real-space grid.
Definition: grid.F90:117
This module defines classes and functions for interaction partners.
This module defines the meshes, which are used in Octopus.
Definition: mesh.F90:118
subroutine, public messages_input_error(namespace, var, details, row, column)
Definition: messages.F90:723
This module defines the quantity_t class and the IDs for quantities, which can be exposed by a system...
Definition: quantity.F90:137
integer, parameter, public current
Definition: quantity.F90:146
Implementation details for regridding.
Definition: regridding.F90:170
subroutine, public build_rs_current_state(current_state, mesh, rs_current_state, ep_field, np)
real(real64) function, public dsm_integrate(mesh, sm, ff, reduce)
Definition: submesh.F90:1163
subroutine, public submesh_end(this)
Definition: submesh.F90:737
subroutine, public submesh_init(this, space, mesh, latt, center, rc)
Definition: submesh.F90:282
brief This module defines the class unit_t which is used by the unit_systems_oct_m module.
Definition: unit.F90:132
This module defines the unit system, used for input and output.
type(unit_system_t), public units_inp
the units systems for reading and writing
Class to transfer a current to a Maxwell field.
class defining the field_transfer interaction
A submesh is a type of mesh, used for the projectors in the pseudopotentials It contains points on a ...
Definition: submesh.F90:177
int true(void)