SPEED
MAKE_MONITOR_FILES.f90
Go to the documentation of this file.
1! Copyright (C) 2012 The SPEED FOUNDATION
2! Author: Ilario Mazzieri
3!
4! This file is part of SPEED.
5!
6! SPEED is free software; you can redistribute it and/or modify it
7! under the terms of the GNU Affero General Public License as
8! published by the Free Software Foundation, either version 3 of the
9! License, or (at your option) any later version.
10!
11! SPEED is distributed in the hope that it will be useful, but
12! WITHOUT ANY WARRANTY; without even the implied warranty of
13! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
14! Affero General Public License for more details.
15!
16! You should have received a copy of the GNU Affero General Public License
17! along with SPEED. If not, see <http://www.gnu.org/licenses/>.
18
33
34 subroutine make_monitor_files(nmonitlst,elem_mlst,local_el_num, ne_loc,&
35 count_monitor, filename, monitor_file, mpi_id, &
36 option_out_var)
37
38 implicit none
39
40 integer*4, intent(inout) :: count_monitor
41 integer*4 :: imon, ielem, ne_loc, ie, mpi_id, unit_monitor, i
42 integer*4 :: nmonitlst
43
44 integer*4, dimension(nmonitlst) :: elem_mlst
45 integer*4, dimension(ne_loc) :: local_el_num
46 integer*4, dimension (6) :: option_out_var
47 integer*4, dimension(:), allocatable :: monitor_index
48
49 character*70 :: filename, monitor_file, filename_new
50 character*5 :: filesuffix
51
52! write(*,*) 'Making MONITOR files'
53
54 count_monitor = 0
55 do imon = 1, nmonitlst
56 ielem = elem_mlst(imon) !ie = index hexahedra containing monitor
57 call get_indloc_from_indglo(local_el_num, ne_loc, ielem, ie)
58 if (ie .ne. 0) count_monitor = count_monitor + 1
59 enddo
60
61 allocate(monitor_index(count_monitor)); count_monitor = 0
62
63 do imon = 1, nmonitlst
64 ielem = elem_mlst(imon) !ie = index hexahedra containing monitor
65 call get_indloc_from_indglo(local_el_num, ne_loc, ielem, ie)
66 if (ie .ne. 0) then
67 count_monitor = count_monitor + 1
68 monitor_index(count_monitor) = imon !node_mlst(imon)
69 endif
70 enddo
71
72 unit_monitor = 40 + mpi_id
73
74 if (count_monitor .ne. 0) then
75
76 ! Loop over variables
77 do i = 1,6
78
79 filesuffix = ''
80 if (i .eq. 1 .and. option_out_var(i) .eq. 1 ) filesuffix = '.D'
81 if (i .eq. 2 .and. option_out_var(i) .eq. 1 ) filesuffix = '.V'
82 if (i .eq. 3 .and. option_out_var(i) .eq. 1 ) filesuffix = '.A'
83 if (i .eq. 4 .and. option_out_var(i) .eq. 1 ) filesuffix = '.S'
84 if (i .eq. 5 .and. option_out_var(i) .eq. 1 ) filesuffix = '.E'
85
86 ! Skip writing if necessary
87 if (option_out_var(i) .eq. 1) then
88 ! Create output filenames for data output
89
90! write(*,'(A,I0,A)') 'mpi_id ', mpi_id, ' suffix ' // filesuffix
91 write(filename, '(A,I5.5,A5)') 'MONITOR', mpi_id, filesuffix
92! write(*,*) 'Variable = ', i
93! write(*,*) 'Output filename = ', filename
94! write(*,*) 'Monitor directory: ', monitor_file
95! write(*,*) 'length ', len_trim(monitor_file)
96
97 ! Prepend directory name
98! if(len_trim(monitor_file) + len_trim(filename) + 1 .le. 255) then
99 if(len_trim(monitor_file) .ne. 70) then
100 filename_new = monitor_file(1:len_trim(monitor_file)) // '/' // filename
101 else
102 filename_new = filename
103 endif
104
105 write(*,*) 'Creating file ', trim(filename_new)
106 open(unit_monitor,file=trim(filename_new)) !open MONITORXXXXX.D, etc...
107 close(unit_monitor)
108 endif
109
110 enddo
111
112 write(filename, '(A,I5.5,A5)') 'MONITOR', mpi_id, '.INFO'
113! write(*,*) 'INFO file: ', filename
114
115 if(len_trim(monitor_file) .ne. 70) then
116! if(len_trim(monitor_file) + len_trim(filename) + 1 .le. 255) then
117 filename_new = monitor_file(1:len_trim(monitor_file)) // '/' // filename
118 else
119 filename_new = filename
120 endif
121
122 write(*,*) 'Creating file ', trim(filename_new)
123 open(unit_monitor,file=trim(filename_new))
124 write(unit_monitor,*) count_monitor, monitor_index
125 close(unit_monitor)
126
127 endif
128
129 ! write(*,*)
130 deallocate(monitor_index)
131
132 end subroutine make_monitor_files
133
134
135
136
subroutine get_indloc_from_indglo(local_el, nel_loc, ie, ic)
Returns local id from global id.