xref: /phasta/phSolver/compressible/timedata.f90 (revision 3e4d567892ef7e8d6a7306f3925c3cf144cef766)
1      !--------------------------------------------------------
2      ! Initialize the Probe Point Arrays and write the Header
3      ! -------------------------------------------------------
4      subroutine initProbePoints()
5        !Tests if the probe point file xyzts.dat exists, loads probe
6        !point locations, initializes a number of arrays used by
7        !timedataC, and writes the initial header for the output file.
8        !
9        ! Rewritten by:             Nicholas Mati  2014-04-18
10        ! Revision history:
11        !  2014-04-18   Code moved from itrdrv to here
12
13        use timedataC
14        use mkdir_mod
15        include "common.h"
16        include "mpif.h"
17
18        logical :: exVarts
19
20        !Test if xyzts.dat exists and broadcast the result.
21        if(myrank.eq.master) inquire(file='xyzts.dat',exist=exts)
22        if(numpe .gt. 1) then
23          call MPI_BARRIER(MPI_COMM_WORLD, ierr)
24          call MPI_BCAST(exts,1,MPI_INTEGER,master,MPI_COMM_WORLD,ierr)
25        endif
26
27        if(.not. exts) return
28        call readProbePoints
29
30        allocate (statptts(ntspts,2))
31        allocate (parptts( ntspts,nsd))
32        allocate (varts(   ntspts,ndof))
33
34        statptts(:,:) = 0
35        parptts(:,:) = zero
36        varts(:,:) = zero
37        ivartsbuff = 0
38        vartsResetbuffer = .false.
39
40        allocate (ivarts(    ntspts*ndof))
41        allocate (ivartsg(   ntspts*ndof))
42        allocate (vartssoln( ntspts*ndof))
43        allocate (vartssolng(ntspts*ndof))
44        allocate (vartsbuff( ntspts,ndof,nbuff))
45        allocate (vartsbuffstep(nbuff))
46
47        !test if the varts folder exists. If it doesn't create it.
48        if(myrank .eq. master) then
49          inquire(file="./varts/.", exist=exVarts)
50          if(.not. exVarts) then
51            call mkdir("varts")
52          endif
53        endif
54
55!       initProbePoints = exts
56!     end function
57      end subroutine
58
59
60      !------------------------
61      ! Read Probe Point Input
62      !------------------------
63      subroutine readProbePoints
64        ! Original Code written by: ??             ????-??-??
65        ! Rewritten by:             Nicholas Mati  2014-04-18
66        ! Revision history:
67        !  2014-04-18   Rewritten code moved from itrdrv to here.
68        !
69        !Reads the file xyzts.dat for probe point locations, write
70        !frequency, tolerance, ... The file is expected to have the
71        !form:
72        ! ntspts freq tolpt iterat nbuff
73        ! x1 y1 z1
74        ! x2 y2 z2
75        ! ...
76        ! xN yN zN
77        !
78        ! ... where ntspts is the number of probe points and freq is the
79        ! number of steps to take before flushing data. If nbuff is
80        ! zero, the number of time steps between restarts, ntout, is
81        ! used.
82
83        use timedataC
84        include "common.h"
85        include "mpif.h"
86
87        if(myrank.eq.master) then
88          open(unit=626,file='xyzts.dat',status='old')
89          read(626,*) ntspts, freq, tolpt, iterat, nbuff
90        endif
91
92        !Broadcase out the header of xyzts.dat. These should probably
93        !be combined into two calls, but this is quick and dirty.
94        if(numpe .gt. 1) then
95          call MPI_BARRIER(MPI_COMM_WORLD, ierr)
96          call MPI_Bcast(ntspts, 1, MPI_INTEGER,          master,
97     &                              MPI_COMM_WORLD,       ierr)
98          call MPI_Bcast(freq,   1, MPI_INTEGER,          master,
99     &                              MPI_COMM_WORLD,       ierr)
100          call MPI_Bcast(tolpt,  1, MPI_DOUBLE_PRECISION, master,
101     &                              MPI_COMM_WORLD,       ierr)
102          call MPI_Bcast(iterat, 1, MPI_INTEGER,          master,
103     &                              MPI_COMM_WORLD,       ierr)
104          call MPI_Bcast(nbuff,  1, MPI_INTEGER,          master,
105     &                              MPI_COMM_WORLD,       ierr)
106        endif
107
108        allocate (ptts(    ntspts,nsd))
109
110        !Read probe point coordinates and broadcast to the rest of the
111        !processors
112        if(myrank .eq. master) then
113          do jj=1,ntspts        ! read coordinate data where solution desired
114             read(626,*) ptts(jj,1), ptts(jj,2), ptts(jj,3)
115           enddo
116         close(626)
117        endif
118
119        if(numpe .gt. 1) then
120          call MPI_BARRIER(MPI_COMM_WORLD, ierr)
121          call MPI_Bcast(ptts, ntspts*nsd, MPI_DOUBLE_PRECISION,
122     &                          master,     MPI_COMM_WORLD,       ierr)
123        endif
124
125        if (nbuff .eq. 0)
126     &    nbuff=ntout
127      end subroutine
128
129
130      !-----------------------------
131      ! Write the Header varts file
132      !-----------------------------
133      subroutine TD_writeHeader(fvarts)
134        !Creates the file fvarts and writes the data header.
135        !fvarts:    Name The file to create
136
137        use timedataC
138        include "common.h"
139
140        character(len=*) fvarts
141
142        !Open the output varts file and write the header
143         if (myrank .eq. master) then
144
145           !fvarts='varts/varts'
146           !fvarts=trim(fvarts)//trim(cname2(lstep))
147           !fvarts=trim(fvarts)//'.dat'
148           open(unit=1001, file=fvarts, status='unknown')
149
150           !Write the time step
151           write(1001, *) "Time Step: ", Delt(1)
152           write(1001, *)
153
154           !Write the probe locations to varts.ts.dat so that post
155           !processing tools actually know what point goes where.
156           !From experience, it's difficult to keep this straight.
157           write(1001, *)
158     &                 "Probe ID   x              y              z"
159           do jj = 1, ntspts
160             write(1001, "(I5, T12, 3(F16.12))") jj, ptts(jj,1:3)
161           enddo
162           write(1001, *)
163
164           !write the output format string. This can't be hard
165           !coded because ntspts is not known in advance.
166           write(vartsIOFrmtStr, '("(I8, ", I4, "(E15.7e2))")')
167     &           ndof*ntspts
168
169           !Header to delinieate the probe locations with the data.
170           write(1001, *) "Probe Data:"
171           close(unit=1001)
172         endif  ! if(myrank .eq. master)
173      end subroutine
174
175
176
177      !------------------------
178      ! Accumulate Probe Data
179      !------------------------
180      subroutine TD_bufferData()
181
182        use timedataC
183        include "common.h"
184        include "mpif.h"
185
186        integer :: icheck, istp, nstp
187
188        if (mod(lstep,freq).eq.0) then
189          if(vartsResetBuffer) then
190            ivartsBuff = 0
191            vartsResetBuffer = .false.
192          endif
193
194          !------------------------
195          !Merge Data across parts
196          !------------------------
197          if (numpe > 1) then
198            !load the contents of varts into vartssoln
199            do jj = 1, ntspts
200               vartssoln((jj-1)*ndof+1:jj*ndof)=varts(jj,:)
201               ivarts=zero
202            enddo
203
204            !mark which points have been computed on this processor
205            do k=1,ndof*ntspts
206               if(vartssoln(k).ne.zero) ivarts(k)=1
207            enddo
208
209            !merge the solution
210            call MPI_REDUCE(vartssoln, vartssolng, ndof*ntspts,
211     &           MPI_DOUBLE_PRECISION, MPI_SUM, master,
212     &           MPI_COMM_WORLD, ierr)
213
214            call MPI_REDUCE(ivarts, ivartsg, ndof*ntspts,
215     &           MPI_INTEGER, MPI_SUM, master,
216     &           MPI_COMM_WORLD, ierr)
217
218             !if the probe point happened to span multiple partitions,
219             !divide the sum by the number of contributing partitions.
220             if (myrank.eq.master) then
221               do jj = 1, ntspts
222                 indxvarts = (jj-1)*ndof
223                 do k=1,ndof
224                   if(ivartsg(indxvarts+k).ne.0) then ! none of the vartssoln(parts) were non zero
225                      varts(jj,k) =
226     &                    vartssolng(indxvarts+k) / ivartsg(indxvarts+k)
227                   endif
228                 enddo !loop over states / DoF
229               enddo !loop over probe points
230             endif !only on master
231          endif !only if numpe > 1
232
233          ivartsBuff = ivartsBuff + 1
234          if (myrank.eq.master) then
235            if(ivartsBuff .gt. nbuff) then
236              write(*,*) "WARNING: vartsbuff has overflowed"
237              ivartsBuff = nbuff
238            endif
239
240            vartsBuffStep(ivartsBuff) = lstep
241            do jj = 1, ntspts
242              vartsbuff(jj,1:ndof, ivartsBuff) = varts(jj,1:ndof)
243            enddo
244          endif
245        endif
246
247        varts(:,:) = zero
248
249      end subroutine
250
251
252      !------------
253      ! Write Data
254      !------------
255      subroutine TD_writeData(fvarts, forceFlush)
256        !writes the probe point data accumulated durring calls to
257        !TD_bufferData(). Note that actual file IO only occurs when the
258        !buffer is full or when DT_writeData is called with forceFlush
259        !set to true. Also note that TD_writeHeader must be called prior
260        !to calling DT_writeData.
261        use timedataC
262        include "common.h"
263
264        character(len=*) :: fvarts
265        logical :: forceFlush
266!        logical, optional :: forceflush
267        logical :: flush
268        integer :: k, ibuf
269
270        if (myrank.eq.master) then
271
272          !if provided, use the default value passed in to determine
273          !wheather to flush the buffer
274!          if(present(forceFlush)) then   !optional version breaks the
275            flush = forceFlush            !compiler on Bluegene?
276!          else
277!            flush = .false.    !set the default value
278!          endif
279
280          !make sure incomplete buffers get purged at the end of a run
281          !regardless of the default.
282!         if(ivartsBuff .eq. nbuff) flush = .true.
283          if(mod(lstep, nbuff) .eq. 0) flush = .true.
284          if(vartsResetBuffer) flush = .false.  !Prevent repeated calls without updating
285                                                          !the buffer from writting multiple times.
286
287          if(flush) then   !flush the buffer to disc
288            open(unit=1001, file      = fvarts,   status = "old",
289     &                  position = "append", action = "write")
290            do ibuf = 1,ivartsBuff
291              write(1001, vartsIOFrmtStr)
292     &                    vartsBuffStep(ibuf),                  !write the time step in the first column.
293     &                   ((vartsbuff(jj,k,ibuf),  k=1, ndof)   !loop over the variables that you actually want to output.
294     &                                         , jj=1, ntspts) !loop over probe points
295            enddo
296
297            close(1001)
298
299            vartsResetBuffer = .true.
300!            ivartsBuff = 0      !need to reset ivartsBuff because
301!                                !writeDate can be called consecutively
302          endif !only dump when buffer full
303        endif !only on master
304
305!              call flush(1001)
306!              call fsync(1001)
307
308               !Code for writting one file per probe point
309!              do jj = 1, ntspts !loop through probe points
310!                ifile = 1000+jj
311!                do ibuf=1,nbuff
312!                  write(ifile,555) lstep-1 -nbuff+ibuf,
313!     &               (vartsbuff(jj,k,ibuf) , k=1, ndof)
314!!     &              vartsbuff(jj,:,ibuf)
315!
316!                enddo ! buff empty
317!
318!                call flush(ifile)
319!              enddo ! jj ntspts
320
321
322!         varts(:,:) = zero ! reset the array for next step   !MOVED FOR Mach Control
323! 555     format(i6,6(2x,E12.5e2))
324
325      end subroutine
326
327
328      subroutine TD_finalize()
329       use timedataC
330
331        deallocate(ivarts)
332        deallocate(ivartsg)
333        deallocate(vartssoln)
334        deallocate(vartssolng)
335        deallocate(vartsbuff)
336        deallocate(vartsbuffstep)
337
338        deallocate(ptts)
339        deallocate(varts)
340      end subroutine
341
342
343      !---------------------
344      ! allocate the arrays
345      !---------------------
346      subroutine sTD
347        !Allocates the arrays statptts, ptts, parptts, and varts for use
348        !in itrdrv and ??
349        !Subroutine is Depricated.
350
351       use timedataC
352        include "common.h"
353
354        allocate (statptts(ntspts,2))
355        allocate (ptts(ntspts,nsd))
356        allocate (parptts(ntspts,nsd))
357        allocate (varts(ntspts,ndof))
358
359        return
360      end
361
362      !-------------------
363      ! delete the arrays
364      !-------------------
365      subroutine dTD
366        !Deallocates ptts and varts
367       use timedataC
368
369        deallocate (ptts)
370        deallocate (varts)
371
372        return
373      end
374