Octopus
field_transfer.F90
Go to the documentation of this file.
1!! Copyright (C) 2023 S. Ohlmann
2!!
3!! This Source Code Form is subject to the terms of the Mozilla Public
4!! License, v. 2.0. If a copy of the MPL was not distributed with this
5!! file, You can obtain one at https://mozilla.org/MPL/2.0/.
6!!
7#include "global.h"
8
48
50 use debug_oct_m
51 use global_oct_m
52 use grid_oct_m
55 use mesh_oct_m
64
65 implicit none
66
67 private
68 public :: &
71
76 type, extends(interaction_t), abstract :: field_transfer_t
77 private
78 real(real64), allocatable, public :: partner_field(:,:)
79 real(real64), allocatable, public :: system_field(:,:)
80 type(grid_t), pointer, public :: system_gr => null()
81 type(regridding_t), pointer, public :: regridding => null()
82 type(time_interpolation_t), pointer, public :: interpolation => null()
83 integer, public :: ndim
84 logical, public :: interpolation_initialized = .false.
85
86 contains
87 procedure :: init => field_transfer_init
88
89 procedure :: init_from_partner => field_transfer_init_from_partner
90
91 procedure :: init_interpolation => field_transfer_init_interpolation
92
93 procedure :: calculate => field_transfer_calculate
94
95 procedure :: do_mapping => field_transfer_do_mapping
96
97 procedure :: dfield_transfer_interpolate, zfield_transfer_interpolate
98 generic :: interpolate => dfield_transfer_interpolate, zfield_transfer_interpolate
99
100 procedure :: calculate_energy => field_transfer_calculate_energy
101
102 procedure :: read_restart => field_transfer_read_restart
103
104 procedure :: write_restart => field_transfer_write_restart
105
106 procedure :: end => field_transfer_end
107
108 end type field_transfer_t
109
110contains
111
113 subroutine field_transfer_init(this, gr, ndim)
114 class(field_transfer_t), intent(inout) :: this
115 type(grid_t), target, intent(in) :: gr
116 integer, intent(in) :: ndim
117
118 push_sub(field_transfer_init)
119
120 this%system_gr => gr
121 this%ndim = ndim
122 safe_allocate(this%system_field(1:gr%np, 1:ndim))
123 this%system_field(:,:) = m_zero
124
125 pop_sub(field_transfer_init)
126 end subroutine field_transfer_init
127
129 subroutine field_transfer_init_from_partner(this, partner_gr, partner_space, partner_namespace)
130 class(field_transfer_t), intent(inout) :: this
131 type(grid_t), intent(in) :: partner_gr
132 class(space_t), intent(in) :: partner_space
133 type(namespace_t), intent(in) :: partner_namespace
134
136
137 safe_allocate(this%partner_field(1:partner_gr%np, 1:this%ndim))
138 this%partner_field(:,:) = m_zero
139 this%regridding => regridding_t(this%system_gr, partner_gr, partner_space, partner_namespace)
140
143
147 subroutine field_transfer_init_interpolation(this, depth, label, cmplx)
148 class(field_transfer_t), intent(inout) :: this
149 integer, intent(in) :: depth
150 character(len=*), intent(in) :: label
151 logical, optional, intent(in) :: cmplx
152
154
155 this%interpolation => time_interpolation_t(this%system_gr%np, this%ndim, depth, &
156 optional_default(cmplx, .false.), trim(label) // '-' // trim(this%partner%namespace%get()))
157 this%interpolation_initialized = .true.
158
161
162 subroutine field_transfer_end(this)
163 class(field_transfer_t), intent(inout) :: this
164
165 push_sub(field_transfer_end)
166
167 call interaction_end(this)
168 safe_deallocate_a(this%partner_field)
169 safe_deallocate_a(this%system_field)
170 safe_deallocate_p(this%regridding)
171 safe_deallocate_p(this%interpolation)
172
174 end subroutine field_transfer_end
177 class(field_transfer_t), intent(inout) :: this
180
181 ! nothing to do here
182 ! the mapping is done when the quantity is copied to the interaction
183 ! for this interaction, we do not need input from the system, so the
184 ! values do only change when the partner quantity changes
185
187 end subroutine field_transfer_calculate
198 class(field_transfer_t), intent(inout) :: this
199 real(real64) :: time
200
201 push_sub_with_profile(field_transfer_do_mapping)
202
203 assert(this%interpolation_initialized)
204
205 if (allocated(this%partner_field)) then
206 call this%regridding%do_transfer(this%system_field, this%partner_field)
207 end if
209 associate(coupling => this%partner%quantities%get(this%couplings_from_partner(1)))
210 time = coupling%iteration%value()
211 end associate
212 call this%interpolation%add_time(time, this%system_field)
213
214 pop_sub_with_profile(field_transfer_do_mapping)
215 end subroutine field_transfer_do_mapping
216
218 subroutine dfield_transfer_interpolate(this, time, field)
219 class(field_transfer_t), intent(in) :: this
220 real(real64), intent(in) :: time
221 real(real64), contiguous, intent(out) :: field(:, :)
222
223
224 push_sub_with_profile(dfield_transfer_interpolate)
225
226 call this%interpolation%interpolate(time, field)
227
228 pop_sub_with_profile(dfield_transfer_interpolate)
229 end subroutine dfield_transfer_interpolate
230
232 subroutine zfield_transfer_interpolate(this, time, field)
233 class(field_transfer_t), intent(in) :: this
234 real(real64), intent(in) :: time
235 complex(real64), contiguous, intent(out) :: field(:, :)
236
237
238 push_sub_with_profile(zfield_transfer_interpolate)
239
240 call this%interpolation%interpolate(time, field)
241
242 pop_sub_with_profile(zfield_transfer_interpolate)
243 end subroutine zfield_transfer_interpolate
244
245 subroutine field_transfer_calculate_energy(this)
246 class(field_transfer_t), intent(inout) :: this
247
249
250 ! interaction energy is zero, since it is only re-gridding the quantities of one system
251 ! on the mesh of the other
252 this%energy = m_zero
253
256
257 subroutine field_transfer_read_restart(this, mesh, space, restart, err)
258 class(field_transfer_t), intent(inout) :: this
259 class(mesh_t), intent(in) :: mesh
260 class(space_t), intent(in) :: space
261 type(restart_t), intent(in) :: restart
262 integer, intent(out) :: err
263
264 push_sub_with_profile(field_transfer_read_restart)
265
266 call this%interpolation%read_restart(mesh, space, restart, err)
267
268 pop_sub_with_profile(field_transfer_read_restart)
269 end subroutine field_transfer_read_restart
270
271 subroutine field_transfer_write_restart(this, mesh, space, restart, err)
272 class(field_transfer_t), intent(inout) :: this
273 class(mesh_t), intent(in) :: mesh
274 class(space_t), intent(in) :: space
275 type(restart_t), intent(in) :: restart
276 integer, intent(out) :: err
277
278 push_sub_with_profile(field_transfer_write_restart)
279
280 call this%interpolation%write_restart(mesh, space, restart, err)
281
282 pop_sub_with_profile(field_transfer_write_restart)
283 end subroutine field_transfer_write_restart
284end module field_transfer_oct_m
285
286!! Local Variables:
287!! mode: f90
288!! coding: utf-8
289!! End:
This module implements the field transfer.
subroutine field_transfer_end(this)
subroutine field_transfer_read_restart(this, mesh, space, restart, err)
subroutine zfield_transfer_interpolate(this, time, field)
return the interpolated field for a given time
subroutine field_transfer_write_restart(this, mesh, space, restart, err)
subroutine field_transfer_do_mapping(this)
perform the regridding and add the system field to the time interpolator using the time of the quanti...
subroutine field_transfer_init_from_partner(this, partner_gr, partner_space, partner_namespace)
the partner field is allocated and initialized to 0; moreover the regridding structure is initialized
subroutine field_transfer_init_interpolation(this, depth, label, cmplx)
the time interpolation is initialized; it needs to know the depth which is usually given by the order...
subroutine, public field_transfer_init(this, gr, ndim)
the system field is allocated and initialized to 0
subroutine dfield_transfer_interpolate(this, time, field)
return the interpolated field for a given time
subroutine field_transfer_calculate_energy(this)
subroutine field_transfer_calculate(this)
real(real64), parameter, public m_zero
Definition: global.F90:191
This module implements the underlying real-space grid.
Definition: grid.F90:119
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
This module defines the quantity_t class and the IDs for quantities, which can be exposed by a system...
Definition: quantity.F90:140
Implementation details for regridding.
Definition: regridding.F90:172
class defining the field_transfer interaction
abstract interaction class
contains the information of the meshes and provides the transfer functions
Definition: regridding.F90:198
int true(void)