Neko 0.9.99
A portable framework for high-order spectral element flow simulations
Loading...
Searching...
No Matches
average_fields_in_time.f90
Go to the documentation of this file.
1
4 use neko
5 implicit none
6
7 character(len=NEKO_FNAME_LEN) :: inputchar, mesh_fname, fld_fname, output_fname
8 type(file_t) :: fld_file, part_file, output_file
9 real(kind=rp) :: start_time
10 type(fld_file_data_t) :: fld_data, fld_data_avg
11 integer :: argc, i
12
13 argc = command_argument_count()
14
15 if ((argc .lt. 3) .or. (argc .gt. 3)) then
16 if (pe_rank .eq. 0) then
17 write(*,*) 'Usage: ./average_fields_in_time field_series_name.fld start_time output_name.fld'
18 write(*,*) 'Example command: ./average_fields_in_time mean_field104.fld 103.2 mean_field_avg.fld'
19 write(*,*) 'Computes the average field over the fld files described in mean_field104.nek5000'
20 write(*,*) 'The start time is the time at which the first file startsto collect stats'
21 write(*,*) 'The files need to be aranged chronological order.'
22 write(*,*) 'The average field is then stored in a fld series, i.e. output_name.nek5000 and output_name.f00000'
23 end if
24 stop
25 end if
26
27 call neko_init
28
29 call get_command_argument(1, inputchar)
30 read(inputchar, *) fld_fname
31 fld_file = file_t(trim(fld_fname))
32 call get_command_argument(2, inputchar)
33 read(inputchar, *) start_time
34 call get_command_argument(3, inputchar)
35 read(inputchar, *) output_fname
36
37 call fld_data%init()
38 call fld_data_avg%init()
39
40 call fld_file%read(fld_data_avg)
41
42 call fld_data_avg%scale(fld_data_avg%time-start_time)
43 if (pe_rank .eq. 0) write(*,*) fld_data_avg%nelv, fld_data_avg%n_scalars
44
45 do i = 1, fld_data_avg%meta_nsamples-1
46 if (pe_rank .eq. 0) write(*,*) 'Reading file:', i+1
47 call fld_file%read(fld_data)
48 call fld_data%scale(fld_data%time-fld_data_avg%time)
49 call fld_data_avg%add(fld_data)
50
51 if (pe_rank .eq. 0) write(*,*) 'dt', fld_data%time - fld_data_avg%time
52 fld_data_avg%time = fld_data%time
53 end do
54 call fld_data_avg%scale(1.0_rp/(fld_data_avg%time-start_time))
55
56 output_file = file_t(trim(output_fname))
57
58
59
60 if (pe_rank .eq. 0) write(*,*) 'Writing file: ', trim(output_fname)
61 call output_file%write(fld_data_avg, fld_data_avg%time)
62 if (pe_rank .eq. 0) write(*,*) 'Done'
63
64 call neko_finalize
65
66end program average_fields_in_time
program average_fields_in_time
Program to sum up averaged fields computed for statistics and mean field Martin Karp 27/01-23.
NEKTON fld file format.
Definition fld_file.f90:35
Master module.
Definition neko.f90:34
subroutine neko_finalize(c)
Definition neko.f90:303
subroutine neko_init(c)
Definition neko.f90:130