Octopus
time_interpolation.F90
Go to the documentation of this file.
1!! Copyright (C) 2002-2006 M. Marques, A. Castro, A. Rubio, G. Bertsch
2!! Copyright (C) 2023 S. Ohlmann
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
23 use debug_oct_m
24 use global_oct_m
26 use math_oct_m
27 use mesh_oct_m
32 use space_oct_m
33
34 implicit none
35
36 private
37 public :: &
39
41 private
42 complex(real64), allocatable :: zfield(:, :, :)
43 real(real64), allocatable :: dfield(:, :, :)
44 real(real64), allocatable :: times(:)
45 integer :: max_depth
46 integer :: depth
47 integer :: np
48 integer :: dim
49 logical :: cmplx
50 character(len=MAX_PATH_LEN) :: label
51 contains
52 procedure :: dtime_interpolation_add_time, ztime_interpolation_add_time
53 generic :: add_time => dtime_interpolation_add_time, ztime_interpolation_add_time
54 procedure :: dtime_interpolation_interpolate, ztime_interpolation_interpolate
55 generic :: interpolate => dtime_interpolation_interpolate, ztime_interpolation_interpolate
56 procedure :: read_restart => time_interpolation_read_restart
57 procedure :: write_restart => time_interpolation_write_restart
60
61 interface time_interpolation_t
62 procedure time_interpolation_constructor
63 end interface time_interpolation_t
64
65contains
66 function time_interpolation_constructor(np, dim, depth, cmplx, label) result(this)
67 integer, intent(in) :: np
68 integer, intent(in) :: dim
69 integer, intent(in) :: depth
70 logical, intent(in) :: cmplx
71 character(len=*), intent(in) :: label
72 class(time_interpolation_t), pointer :: this
73
75
76 safe_allocate(this)
77
78 this%np = np
79 this%dim = dim
80 this%max_depth = depth
81 this%depth = 0
82 this%cmplx = cmplx
83 this%label = trim(label)
84
85 if (this%cmplx) then
86 safe_allocate(this%zfield(1:np, 1:dim, 1:depth))
87 else
88 safe_allocate(this%dfield(1:np, 1:dim, 1:depth))
89 end if
90 safe_allocate(this%times(1:depth))
91
94
95 subroutine time_interpolation_finalize(this)
96 type(time_interpolation_t), intent(inout) :: this
97
99
100 if (this%cmplx) then
101 safe_deallocate_a(this%zfield)
102 else
103 safe_deallocate_a(this%dfield)
104 end if
105 safe_deallocate_a(this%times)
106
108 end subroutine time_interpolation_finalize
109
110 subroutine time_interpolation_write_restart(this, mesh, space, restart, err)
111 class(time_interpolation_t), intent(in) :: this
112 class(mesh_t), intent(in) :: mesh
113 class(space_t), intent(in) :: space
114 type(restart_t), intent(in) :: restart
115 integer, intent(out) :: err
116
117 integer :: itime, idim, err_restart, iunit
118 character(len=MAX_PATH_LEN) :: filename, lines(1)
119
121
122 ! write fields
123 err = 0
124 do itime = 1, this%depth
125 do idim = 1, this%dim
126 write(filename, '(a1,i2.2,a1,i3.3)') '_', itime, '_', idim
127 filename = "field_" // trim(this%label) // trim(filename)
128 if (this%cmplx) then
129 call restart_write_mesh_function(restart, space, filename, mesh, &
130 this%zfield(1:this%np, idim, itime), err_restart)
131 else
132 call restart_write_mesh_function(restart, space, filename, mesh, &
133 this%dfield(1:this%np, idim, itime), err_restart)
134 end if
135 if (err_restart /= 0) err = err + 1
136 end do
137 end do
139 ! write times
140 call drestart_write_binary(restart, "field_times_"//trim(this%label), this%depth, this%times(1:this%depth), err_restart)
141 if (err_restart /= 0) err = err + 1
143 ! write depth
144 iunit = restart_open(restart, "field_"//trim(this%label))
145 write(lines(1), '(i2.2)') this%depth
146 call restart_write(restart, iunit, lines, 1, err_restart)
147 if (err_restart /= 0) err = err + 1
148 call restart_close(restart, iunit)
152
153 subroutine time_interpolation_read_restart(this, mesh, space, restart, err)
154 class(time_interpolation_t), intent(inout) :: this
155 class(mesh_t), intent(in) :: mesh
156 class(space_t), intent(in) :: space
157 type(restart_t), intent(in) :: restart
158 integer, intent(out) :: err
160 integer :: itime, idim, err_restart, iunit
161 character(len=MAX_PATH_LEN) :: filename, lines(1)
162
164
165 ! read depth
166 iunit = restart_open(restart, "field_"//trim(this%label))
167 call restart_read(restart, iunit, lines, 1, err_restart)
168 if (err_restart /= 0) then
169 err = err + 1
171 return
172 else
173 read(lines(1), '(i2.2)') this%depth
174 end if
175 call restart_close(restart, iunit)
176
177 ! read fields
178 err = 0
179 do itime = 1, this%depth
180 do idim = 1, this%dim
181 write(filename, '(a1,i2.2,a1,i3.3)') '_', itime, '_', idim
182 filename = "field_" // trim(this%label) // trim(filename)
183 if (this%cmplx) then
184 call restart_read_mesh_function(restart, space, filename, mesh, &
185 this%zfield(1:this%np, idim, itime), err_restart)
186 else
187 call restart_read_mesh_function(restart, space, filename, mesh, &
188 this%dfield(1:this%np, idim, itime), err_restart)
189 end if
190 if (err_restart /= 0) err = err + 1
191 end do
192 end do
193
194 ! read times
195 call drestart_read_binary(restart, "field_times_"//trim(this%label), this%depth, this%times(1:this%depth), err_restart)
196 if (err_restart /= 0) err = err + 1
197
200
201#include "real.F90"
202#include "time_interpolation_inc.F90"
203#include "undef.F90"
204
205#include "complex.F90"
206#include "time_interpolation_inc.F90"
207#include "undef.F90"
208
210
211!! Local Variables:
212!! mode: f90
213!! coding: utf-8
214!! End:
This is the common interface to a simple-minded polynomical interpolation procedure (simple use of th...
Definition: math.F90:186
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
subroutine, public restart_read(restart, iunit, lines, nlines, ierr)
Definition: restart.F90:933
subroutine, public restart_close(restart, iunit)
Close a file previously opened with restart_open.
Definition: restart.F90:950
subroutine, public restart_write(restart, iunit, lines, nlines, ierr)
Definition: restart.F90:906
integer function, public restart_open(restart, filename, status, position, silent)
Open file "filename" found inside the current restart directory. Depending on the type of restart,...
Definition: restart.F90:859
subroutine time_interpolation_read_restart(this, mesh, space, restart, err)
subroutine time_interpolation_write_restart(this, mesh, space, restart, err)
subroutine time_interpolation_finalize(this)
class(time_interpolation_t) function, pointer time_interpolation_constructor(np, dim, depth, cmplx, label)
Describes mesh distribution to nodes.
Definition: mesh.F90:186