KTH framework for Nek5000 toolboxes; testing version  0.0.1
pstat3D.f
Go to the documentation of this file.
1 
6 !=======================================================================
10  subroutine pstat3d_register()
11  implicit none
12 
13  include 'SIZE'
14  include 'INPUT'
15  include 'FRAMELP'
16  include 'PSTAT3D'
17 
18  ! local variables
19  integer lpmid, il
20  real ltim
21  character*2 str
22 
23  ! functions
24  real dnekclock
25 !-----------------------------------------------------------------------
26  ! timing
27  ltim = dnekclock()
28 
29  ! check if the current module was already registered
30  call mntr_mod_is_name_reg(lpmid,pstat_name)
31  if (lpmid.gt.0) then
32  call mntr_warn(lpmid,
33  $ 'module ['//trim(pstat_name)//'] already registered')
34  return
35  endif
36 
37  ! find parent module
38  call mntr_mod_is_name_reg(lpmid,'FRAME')
39  if (lpmid.le.0) then
40  lpmid = 1
41  call mntr_abort(lpmid,
42  $ 'parent module ['//'FRAME'//'] not registered')
43  endif
44 
45  ! register module
46  call mntr_mod_reg(pstat_id,lpmid,pstat_name,
47  $ 'Post processing for 3D statistics')
48 
49  ! register timers
50  call mntr_tmr_is_name_reg(lpmid,'FRM_TOT')
51  ! total time
52  call mntr_tmr_reg(pstat_tmr_tot_id,lpmid,pstat_id,
53  $ 'PSTAT_TOT','Pstat total time',.false.)
54  lpmid = pstat_tmr_tot_id
55  ! initialisation
56  call mntr_tmr_reg(pstat_tmr_ini_id,lpmid,pstat_id,
57  $ 'PSTAT_INI','Pstat initialisation time',.true.)
58  ! averaging
59  call mntr_tmr_reg(pstat_tmr_avg_id,lpmid,pstat_id,
60  $ 'PSTAT_AVG','Pstat averaging time',.true.)
61  ! new field calculation
62  call mntr_tmr_reg(pstat_tmr_new_id,lpmid,pstat_id,
63  $ 'PSTAT_NEW','Pstat new field calculation time',.true.)
64  ! interpolation
65  call mntr_tmr_reg(pstat_tmr_int_id,lpmid,pstat_id,
66  $ 'PSTAT_INT','Pstat interpolation time',.true.)
67 
68  ! register and set active section
69  call rprm_sec_reg(pstat_sec_id,pstat_id,'_'//adjustl(pstat_name),
70  $ 'Runtime paramere section for pstat module')
71  call rprm_sec_set_act(.true.,pstat_sec_id)
72 
73  ! register parameters
74  call rprm_rp_reg(pstat_nfile_id,pstat_sec_id,'STS_NFILE',
75  $ 'Number of stat files',rpar_int,1,0.0,.false.,' ')
76  call rprm_rp_reg(pstat_stime_id,pstat_sec_id,'STS_STIME',
77  $ 'Statistics starting time',rpar_real,1,0.0,.false.,' ')
78  call rprm_rp_reg(pstat_nstep_id,pstat_sec_id,'STS_NSTEP',
79  $ 'Number of steps between averaging (in sts file)',
80  $ rpar_int,10,0.0,.false.,' ')
81 
82  ! set initialisation flag
83  pstat_ifinit=.false.
84 
85  ! timing
86  ltim = dnekclock() - ltim
87  call mntr_tmr_add(pstat_tmr_tot_id,1,ltim)
88 
89  return
90  end subroutine
91 !=======================================================================
95  subroutine pstat3d_init()
96  implicit none
97 
98  include 'SIZE'
99  include 'FRAMELP'
100  include 'PSTAT3D'
101 
102  ! local variables
103  integer itmp, il
104  real rtmp, ltim
105  logical ltmp
106  character*20 ctmp
107  integer tmp_swfield(pstat_svar)
108 
109  ! functions
110  real dnekclock
111 !-----------------------------------------------------------------------
112 ! check if the module was already initialised
113  if (pstat_ifinit) then
114  call mntr_warn(pstat_id,
115  $ 'module ['//trim(pstat_name)//'] already initiaised.')
116  return
117  endif
118 
119  ! timing
120  ltim = dnekclock()
121 
122  ! get runtime parameters
123  call rprm_rp_get(itmp,rtmp,ltmp,ctmp,pstat_nfile_id,rpar_int)
124  pstat_nfile = abs(itmp)
125  call rprm_rp_get(itmp,rtmp,ltmp,ctmp,pstat_stime_id,rpar_real)
126  pstat_stime = rtmp
127  call rprm_rp_get(itmp,rtmp,ltmp,ctmp,pstat_nstep_id,rpar_int)
128  pstat_nstep = abs(itmp)
129 
130  ! set field swapping array
131  do il = 1,26
132  tmp_swfield(il) = il
133  enddo
134  do il = 27,32
135  tmp_swfield(il) = il+1
136  enddo
137  tmp_swfield(33) = 27
138  tmp_swfield(34) = 38
139  do il = 35,38
140  tmp_swfield(il) = il-1
141  enddo
142  do il = 39,44
143  tmp_swfield(il) = il
144  enddo
145  call icopy(pstat_swfield,tmp_swfield,pstat_svar)
146 
147  ! everything is initialised
148  pstat_ifinit=.true.
149 
150  ! timing
151  ltim = dnekclock() - ltim
152  call mntr_tmr_add(pstat_tmr_ini_id,1,ltim)
153 
154  return
155  end subroutine
156 !=======================================================================
160  logical function pstat3d_is_initialised()
161  implicit none
162 
163  include 'SIZE'
164  include 'PSTAT3D'
165 !-----------------------------------------------------------------------
166  pstat3d_is_initialised = pstat_ifinit
167 
168  return
169  end function
170 !=======================================================================
173  subroutine pstat3d_main
174  implicit none
175 
176  include 'SIZE'
177  include 'FRAMELP'
178  include 'PSTAT3D'
179 
180  ! local variables
181  integer il
182 
183  ! simple timing
184  real ltim
185 
186  ! functions
187  real dnekclock
188 !-----------------------------------------------------------------------
189  ! read and average fields
190  ltim = dnekclock()
191  call mntr_log(pstat_id,lp_inf,'Field averaging')
192  call pstat3d_sts_avg
193  ltim = dnekclock() - ltim
194  call mntr_tmr_add(pstat_tmr_avg_id,1,ltim)
195 
196  ! calculate new fields
197  ltim = dnekclock()
198  call mntr_log(pstat_id,lp_inf,'New field calculation')
199  call pstat3d_nfield
200  ltim = dnekclock() - ltim
201  call mntr_tmr_add(pstat_tmr_new_id,1,ltim)
202 
203  ! interpolate into the set of points
204  ltim = dnekclock()
205  call mntr_log(pstat_id,lp_inf,'Point interpolation')
206  call pstat3d_interp
207  ltim = dnekclock() - ltim
208  call mntr_tmr_add(pstat_tmr_int_id,1,ltim)
209 
210  return
211  end subroutine
212 !=======================================================================
215  subroutine pstat3d_sts_avg
216  implicit none
217 
218  include 'SIZE'
219  include 'INPUT'
220  include 'RESTART'
221  include 'TSTEP'
222  include 'SOLN'
223  include 'FRAMELP'
224  include 'PSTAT3D'
225 
226  ! local variables
227  integer il, jl ! loop index
228  integer ierr ! error mark
229  integer nvec ! single field length
230  character*3 prefix ! file prefix
231  character*2 fins ! number of a file in a section
232  character*132 fname ! file name
233  character*132 bname ! base name
234  character*6 str ! file number
235  integer nps1,nps0,npsr ! number of passive scalars in the file
236  integer istepr ! number of time step in restart files
237  real ltime, dtime ! simulation time and time update
238  real rtmp
239 
240 !-----------------------------------------------------------------------
241  ! no regular mesh; important for file name generation
242  ifreguo = .false.
243 
244  call io_init
245 
246  ierr=0
247  ! open files on i/o nodes
248  ! get base name (SESSION)
249  bname = trim(adjustl(session))
250 
251  ! mark variables to be read
252  ifgetx=.false.
253  ifgetz=.false.
254  ifgetu=.true.
255  ifgetw=.true.
256  ifgetp=.false.
257  ifgett=.true.
258  do il=1, npscal
259  ifgtps(il)=.false.
260  enddo
261  ! mark variables for saving
262  ifxyo = .true.
263  ifvo = .true.
264  ifpo = .false.
265  ifto = .true.
266  do il=1, npscal
267  ifpsco(il)= .false.
268  enddo
269 
270  ! initilise vectors
271  call rzero(pstat_ruavg,lx1**ldim*lelt*pstat_svar)
272  nvec = lx1*ly1*lz1*nelt
273 
274  ! loop over stat files
275  ! to be able to save averaged fields in the order of s** files
276  ! I average each s** file separately
277  do jl=1, pstat_finset
278  ! initial time and step count
279  ltime = pstat_stime
280  istepr = 0
281 
282  ! get prefix for reading
283  write(fins,'(i2.2)') jl
284  prefix = 's'//trim(fins)
285 
286  ! loop over time samples
287  do il = 1,pstat_nfile
288 
289  ! get file name
290  call io_mfo_fname(fname,bname,prefix,ierr)
291  write(str,'(i5.5)') il
292  fname = trim(fname)//trim(str)
293 
294  fid0 = 0
295  call addfid(fname,fid0)
296 
297  ! add directory name
298  fname = 'DATA/'//trim(fname)
299 
300  !call load_fld(fname)
301  call mfi(fname,il)
302 
303  ! calculate interval and update time
304  dtime = timer - ltime
305 
306  ! accumulate fileds
307  call add2s2(pstat_ruavg(1,1,pstat_swfield(1,jl)),
308  $ vx,dtime,nvec)
309  call add2s2(pstat_ruavg(1,1,pstat_swfield(2,jl)),
310  $ vy,dtime,nvec)
311  call add2s2(pstat_ruavg(1,1,pstat_swfield(3,jl)),
312  $ vz,dtime,nvec)
313  call add2s2(pstat_ruavg(1,1,pstat_swfield(4,jl)),
314  $ t,dtime,nvec)
315 ! enddo
316 
317  ! sum number of time steps
318  istepr = istepr + istpr
319 
320  ! update time
321  ltime = timer
322 
323  end do
324 
325  ! divide by time span
326  if (ltime.ne.pstat_stime) then
327  rtmp = 1.0/(ltime-pstat_stime)
328  else
329  rtmp = 1.0
330  endif
331  do il = 1,4
332  call cmult(pstat_ruavg(1,1,pstat_swfield(il,jl)),rtmp,nvec)
333  enddo
334  ! save the averaged filed
335  ! get prefix for writing
336  write(fins,'(i2.2)') jl
337  prefix = 't'//trim(fins)
338  call outpost(pstat_ruavg(1,1,pstat_swfield(1,jl)),
339  $ pstat_ruavg(1,1,pstat_swfield(2,jl)),
340  $ pstat_ruavg(1,1,pstat_swfield(3,jl)),pr,
341  $ pstat_ruavg(1,1,pstat_swfield(4,jl)),prefix)
342  enddo
343 
344  ! save data for file header
345  pstat_etime = ltime
346  pstat_istepr = istepr
347 
348  return
349  end subroutine
350 !=======================================================================
353  subroutine pstat3d_nfield
354  implicit none
355 
356  include 'SIZE'
357  include 'INPUT'
358  include 'SOLN'
359  include 'FRAMELP'
360  include 'PSTAT3D'
361 
362  ! local variables
363  integer il, jl ! loop index
364  integer nvec ! single field length
365  integer itmp
366  real rtmp, rho
367 !-----------------------------------------------------------------------
368  ! update derivative arrays
369  !call geom_reset(1)
370 
371  ! initilise vectors
372  call rzero(pstat_runew,lx1**ldim*lelt*pstat_dvar)
373  call rzero(pstat_rutmp,lx1**ldim*lelt*pstat_tvar)
374  nvec = lx1*ly1*lz1*nelt
375  rho = 1.0
376 
377  ! mean variables U, V, W, P already calculated and stored under
378  ! U (avg 1), V (avg 2), W (avg 3), P (avg 4)
379  ! Reynolds-stress tensor diagonal terms and pressure RMS
380  ! U*U (tmp 1), V*V (tmp 2), W*W (tmp 3), P*P (tmp 4)
381  ! uu (avg 5), vv (avg 6), ww (avg 7), pp (avg 8)
382  do il = 1, 4
383  call col3(pstat_rutmp(1,1,il),pstat_ruavg(1,1,il),
384  $ pstat_ruavg(1,1,il),nvec)
385  call sub2(pstat_ruavg(1,1,4+il),pstat_rutmp(1,1,il),nvec)
386  enddo
387 
388  ! Reynolds-stress tensor off diagonal terms
389  call subcol3(pstat_ruavg(1,1,9),pstat_ruavg(1,1,1), ! uv (avg 9)
390  $ pstat_ruavg(1,1,2),nvec)
391  call subcol3(pstat_ruavg(1,1,10),pstat_ruavg(1,1,2), ! vw (avg 10)
392  $ pstat_ruavg(1,1,3),nvec)
393  call subcol3(pstat_ruavg(1,1,11),pstat_ruavg(1,1,1), ! uw (avg 11)
394  $ pstat_ruavg(1,1,3),nvec)
395 
396  ! pressure skewness
397  call col3(pstat_rutmp(1,1,5),pstat_rutmp(1,1,4), ! P*P*P (tmp 5)
398  $ pstat_ruavg(1,1,4),nvec)
399  call sub2(pstat_ruavg(1,1,27),pstat_rutmp(1,1,5),nvec) ! ppp (avg 27)
400  rtmp = -3.0
401  call admcol3(pstat_ruavg(1,1,27),pstat_ruavg(1,1,4),
402  $ pstat_ruavg(1,1,8),rtmp,nvec)
403 
404  ! pressure flattness
405  rtmp = -4.0
406  call admcol3(pstat_ruavg(1,1,38),pstat_ruavg(1,1,4), ! pppp (avg 38)
407  $ pstat_ruavg(1,1,27),rtmp,nvec)
408  rtmp = -6.0
409  call admcol3(pstat_ruavg(1,1,38),pstat_ruavg(1,1,8),
410  $ pstat_rutmp(1,1,4),rtmp,nvec)
411  call subcol3(pstat_ruavg(1,1,38),pstat_ruavg(1,1,4),
412  $ pstat_rutmp(1,1,5),nvec)
413 
414  ! Skewness tensor
415  ! diagonal terms
416  ! uuu (avg 24), vvv (avg 25), www (avg 26)
417  rtmp = -3.0
418  do il= 1, 3
419  call admcol3(pstat_ruavg(1,1,23+il),pstat_ruavg(1,1,il),
420  $ pstat_ruavg(1,1,4+il),rtmp,nvec)
421  call subcol3(pstat_ruavg(1,1,23+il),pstat_ruavg(1,1,il),
422  $ pstat_rutmp(1,1,il),nvec)
423  enddo
424 
425  ! off diagonal terms
426  rtmp = -2.0
427  call admcol3(pstat_ruavg(1,1,28),pstat_ruavg(1,1,1), ! uuv (avg 28)
428  $ pstat_ruavg(1,1,9),rtmp,nvec)
429  call subcol3(pstat_ruavg(1,1,28),pstat_ruavg(1,1,2),
430  $ pstat_ruavg(1,1,5),nvec)
431  call subcol3(pstat_ruavg(1,1,28),pstat_ruavg(1,1,2),
432  $ pstat_rutmp(1,1,1),nvec)
433  call admcol3(pstat_ruavg(1,1,29),pstat_ruavg(1,1,1), ! uuw (avg 29)
434  $ pstat_ruavg(1,1,11),rtmp,nvec)
435  call subcol3(pstat_ruavg(1,1,29),pstat_ruavg(1,1,3),
436  $ pstat_ruavg(1,1,5),nvec)
437  call subcol3(pstat_ruavg(1,1,29),pstat_ruavg(1,1,3),
438  $ pstat_rutmp(1,1,1),nvec)
439  call admcol3(pstat_ruavg(1,1,30),pstat_ruavg(1,1,2), ! uvv (avg 30)
440  $ pstat_ruavg(1,1,9),rtmp,nvec)
441  call subcol3(pstat_ruavg(1,1,30),pstat_ruavg(1,1,1),
442  $ pstat_ruavg(1,1,6),nvec)
443  call subcol3(pstat_ruavg(1,1,30),pstat_ruavg(1,1,1),
444  $ pstat_rutmp(1,1,2),nvec)
445  call admcol3(pstat_ruavg(1,1,31),pstat_ruavg(1,1,2), ! vvw (avg 31)
446  $ pstat_ruavg(1,1,10),rtmp,nvec)
447  call subcol3(pstat_ruavg(1,1,31),pstat_ruavg(1,1,3),
448  $ pstat_ruavg(1,1,6),nvec)
449  call subcol3(pstat_ruavg(1,1,31),pstat_ruavg(1,1,3),
450  $ pstat_rutmp(1,1,2),nvec)
451  call admcol3(pstat_ruavg(1,1,32),pstat_ruavg(1,1,3), ! uww (avg 32)
452  $ pstat_ruavg(1,1,11),rtmp,nvec)
453  call subcol3(pstat_ruavg(1,1,32),pstat_ruavg(1,1,1),
454  $ pstat_ruavg(1,1,7),nvec)
455  call subcol3(pstat_ruavg(1,1,32),pstat_ruavg(1,1,1),
456  $ pstat_rutmp(1,1,3),nvec)
457  call admcol3(pstat_ruavg(1,1,33),pstat_ruavg(1,1,3), ! vww (avg 33)
458  $ pstat_ruavg(1,1,10),rtmp,nvec)
459  call subcol3(pstat_ruavg(1,1,33),pstat_ruavg(1,1,2),
460  $ pstat_ruavg(1,1,7),nvec)
461  call subcol3(pstat_ruavg(1,1,33),pstat_ruavg(1,1,2),
462  $ pstat_rutmp(1,1,3),nvec)
463  call subcol3(pstat_ruavg(1,1,34),pstat_ruavg(1,1,1), ! uvw (avg 34)
464  $ pstat_ruavg(1,1,10),nvec)
465  call subcol3(pstat_ruavg(1,1,34),pstat_ruavg(1,1,2),
466  $ pstat_ruavg(1,1,11),nvec)
467  call subcol3(pstat_ruavg(1,1,34),pstat_ruavg(1,1,3),
468  $ pstat_ruavg(1,1,9),nvec)
469  call subcol4(pstat_ruavg(1,1,34),pstat_ruavg(1,1,1),
470  $ pstat_ruavg(1,1,2),pstat_ruavg(1,1,3),nvec)
471 
472  ! Velocity gradient tensor
473  ! dUdx (new 1), dUdy (new 2), dUdz (new 3),
474  ! dVdx (new 4), dVdy (new 5), dVdz (new 6),
475  ! dWdx (new 7), dWdy (new 8), dWdz (new 9)
476  do il=1,3
477  itmp = (il-1)*3
478  call gradm1(pstat_runew(1,1,itmp+1),pstat_runew(1,1,itmp+2),
479  $ pstat_runew(1,1,itmp+3),pstat_ruavg(1,1,il))
480  enddo
481 
482  ! Dissipation tensor
483  rtmp = -2.0*param(2)
484  ! diagonal terms
485  ! Dxx (avg 39), Dyy (avg 40), Dzz (avg 41)
486  do il=1,3
487  itmp = (il-1)*3
488  do jl = 1,3
489  call subcol3(pstat_ruavg(1,1,38+il),
490  $ pstat_runew(1,1,itmp+jl),pstat_runew(1,1,itmp+jl),
491  $ nvec)
492  enddo
493  call cmult(pstat_ruavg(1,1,38+il),rtmp,nvec)
494  enddo
495  ! off diagonal terms
496  do jl = 1,3
497  call subcol3(pstat_ruavg(1,1,42),pstat_runew(1,1,jl), ! Dxy (avg 42)
498  $ pstat_runew(1,1,3+jl),nvec)
499  enddo
500  call cmult(pstat_ruavg(1,1,42),rtmp,nvec)
501  do jl = 1,3
502  call subcol3(pstat_ruavg(1,1,43),pstat_runew(1,1,jl), ! Dxz (avg 43)
503  $ pstat_runew(1,1,6+jl),nvec)
504  enddo
505  call cmult(pstat_ruavg(1,1,43),rtmp,nvec)
506  do jl = 1,3
507  call subcol3(pstat_ruavg(1,1,44),pstat_runew(1,1,3+jl), ! Dyz (avg 44)
508  $ pstat_runew(1,1,6+jl),nvec)
509  enddo
510  call cmult(pstat_ruavg(1,1,44),rtmp,nvec)
511 
512  ! Derivatives of the Reynolds-stress tensor
513  ! duudx (tmp 2), duudy (tmp 3), duudz (tmp 4),
514  ! dvvdx (tmp 5), dvvdy (tmp 6), dvvdz (tmp 7),
515  ! dwwdx (tmp 8), dwwdy (tmp 9), dwwdz (tmp 10),
516  ! duvdx (tmp 11), duvdy (tmp 12), duvdz (tmp 13),
517  ! dvwdx (tmp 14), dvwdy (tmp 15), dvwdz (tmp 16),
518  ! duwdx (tmp 17), duwdy (tmp 18), duwdz (tmp 19)
519  do il = 1,3
520  itmp = (il-1)*3
521  call gradm1(pstat_rutmp(1,1,itmp+2),pstat_rutmp(1,1,itmp+3),
522  $ pstat_rutmp(1,1,itmp+4),pstat_ruavg(1,1,4+il))
523  itmp = (il+2)*3
524  call gradm1(pstat_rutmp(1,1,itmp+2),pstat_rutmp(1,1,itmp+3),
525  $ pstat_rutmp(1,1,itmp+4),pstat_ruavg(1,1,8+il))
526  enddo
527 
528  ! Mean convection tensor
529  ! Cxx (new 10), Cyy (new 11), Czz (new 12), Cxy (new 13), Cyz (new 14), Cxz (new 15)
530  do il = 1,6
531  itmp = (il-1)*3
532  call vdot3 (pstat_runew(1,1,9+il),
533  $ pstat_ruavg(1,1,1),pstat_ruavg(1,1,2),pstat_ruavg(1,1,3),
534  $ pstat_rutmp(1,1,itmp+2),pstat_rutmp(1,1,itmp+3),
535  $ pstat_rutmp(1,1,itmp+4),nvec)
536  enddo
537 
538  ! Second derivatives of the Reynolds-stress tensor components
539  call gradm1(pstat_rutmp(1,1,1),pstat_rutmp(1,1,20), ! d2uudx2 (tmp 1)
540  $ pstat_rutmp(1,1,21),pstat_rutmp(1,1,2))
541  call gradm1(pstat_rutmp(1,1,20),pstat_rutmp(1,1,2), ! d2uudy2 (tmp 2)
542  $ pstat_rutmp(1,1,21),pstat_rutmp(1,1,3))
543  call gradm1(pstat_rutmp(1,1,20),pstat_rutmp(1,1,21), ! d2uudz2 (tmp 3)
544  $ pstat_rutmp(1,1,3),pstat_rutmp(1,1,4))
545  call gradm1(pstat_rutmp(1,1,4),pstat_rutmp(1,1,20), ! d2vvdx2 (tmp 4)
546  $ pstat_rutmp(1,1,21),pstat_rutmp(1,1,5))
547  call gradm1(pstat_rutmp(1,1,20),pstat_rutmp(1,1,5), ! d2vvdy2 (tmp 5)
548  $ pstat_rutmp(1,1,21),pstat_rutmp(1,1,6))
549  call gradm1(pstat_rutmp(1,1,20),pstat_rutmp(1,1,21), ! d2vvdz2 (tmp 6)
550  $ pstat_rutmp(1,1,6),pstat_rutmp(1,1,7))
551  call gradm1(pstat_rutmp(1,1,7),pstat_rutmp(1,1,20), ! d2wwdx2 (tmp 7)
552  $ pstat_rutmp(1,1,21),pstat_rutmp(1,1,8))
553  call gradm1(pstat_rutmp(1,1,20),pstat_rutmp(1,1,8), ! d2wwdy2 (tmp 8)
554  $ pstat_rutmp(1,1,21),pstat_rutmp(1,1,9))
555  call gradm1(pstat_rutmp(1,1,20),pstat_rutmp(1,1,21), ! d2wwdz2 (tmp 9)
556  $ pstat_rutmp(1,1,9),pstat_rutmp(1,1,10))
557  call gradm1(pstat_rutmp(1,1,10),pstat_rutmp(1,1,20), ! d2uvdx2 (tmp 10)
558  $ pstat_rutmp(1,1,21),pstat_rutmp(1,1,11))
559  call gradm1(pstat_rutmp(1,1,20),pstat_rutmp(1,1,11), ! d2uvdy2 (tmp 11)
560  $ pstat_rutmp(1,1,21),pstat_rutmp(1,1,12))
561  call gradm1(pstat_rutmp(1,1,20),pstat_rutmp(1,1,21), ! d2uvdz2 (tmp 12)
562  $ pstat_rutmp(1,1,12),pstat_rutmp(1,1,13))
563  call gradm1(pstat_rutmp(1,1,13),pstat_rutmp(1,1,20), ! d2vwdx2 (tmp 13)
564  $ pstat_rutmp(1,1,21),pstat_rutmp(1,1,14))
565  call gradm1(pstat_rutmp(1,1,20),pstat_rutmp(1,1,14), ! d2vwdy2 (tmp 14)
566  $ pstat_rutmp(1,1,21),pstat_rutmp(1,1,15))
567  call gradm1(pstat_rutmp(1,1,20),pstat_rutmp(1,1,21), ! d2vwdz2 (tmp 15)
568  $ pstat_rutmp(1,1,15),pstat_rutmp(1,1,16))
569  call gradm1(pstat_rutmp(1,1,16),pstat_rutmp(1,1,20), ! d2uwdx2 (tmp 16)
570  $ pstat_rutmp(1,1,21),pstat_rutmp(1,1,17))
571  call gradm1(pstat_rutmp(1,1,20),pstat_rutmp(1,1,17), ! d2uwdy2 (tmp 17)
572  $ pstat_rutmp(1,1,21),pstat_rutmp(1,1,18))
573  call gradm1(pstat_rutmp(1,1,20),pstat_rutmp(1,1,21), ! d2uwdz2 (tmp 18)
574  $ pstat_rutmp(1,1,18),pstat_rutmp(1,1,19))
575 
576  ! Viscous diffusion tensor
577  ! VDxx (new 16), VDyy (new 17), VDzz (new 18), VDxy (new 19), VDyz (new 20), VDxz (new 21)
578  rtmp = param(2)
579  do il = 1,6
580  itmp = (il-1)*3
581  call add4(pstat_runew(1,1,15+il),pstat_rutmp(1,1,itmp+1),
582  $ pstat_rutmp(1,1,itmp+2),pstat_rutmp(1,1,itmp+3),nvec)
583  call cmult(pstat_runew(1,1,15+il),rtmp,nvec)
584  enddo
585 
586  ! Derivatives of the triple-product terms
587  call gradm1(pstat_rutmp(1,1,1),pstat_rutmp(1,1,20), ! duuudx (tmp 1)
588  $ pstat_rutmp(1,1,21),pstat_ruavg(1,1,24))
589  call gradm1(pstat_rutmp(1,1,4),pstat_rutmp(1,1,11), ! duvvdx (tmp 4), duvvdy (tmp 11)
590  $ pstat_rutmp(1,1,20),pstat_ruavg(1,1,30))
591  call gradm1(pstat_rutmp(1,1,7),pstat_rutmp(1,1,20), ! duwwdx (tmp 7), duwwdz (tmp 18)
592  $ pstat_rutmp(1,1,15),pstat_ruavg(1,1,32))
593  call gradm1(pstat_rutmp(1,1,10),pstat_rutmp(1,1,2), ! duuvdx (tmp 10), duuvdy (tmp 2)
594  $ pstat_rutmp(1,1,20),pstat_ruavg(1,1,28))
595  call gradm1(pstat_rutmp(1,1,13),pstat_rutmp(1,1,20), ! duuwdx (tmp 16), duuwdz (tmp 3)
596  $ pstat_rutmp(1,1,3),pstat_ruavg(1,1,29))
597  call gradm1(pstat_rutmp(1,1,16),pstat_rutmp(1,1,14), ! duvwdx (tmp 13), duvwdy (tmp 17), duvwdz (tmp 12)
598  $ pstat_rutmp(1,1,12),pstat_ruavg(1,1,34))
599  call gradm1(pstat_rutmp(1,1,20),pstat_rutmp(1,1,5), ! dvvvdy (tmp 5)
600  $ pstat_rutmp(1,1,21),pstat_ruavg(1,1,25))
601  call gradm1(pstat_rutmp(1,1,20),pstat_rutmp(1,1,8), ! dvwwdy (tmp 8), dvwwdz (tmp 15)
602  $ pstat_rutmp(1,1,18),pstat_ruavg(1,1,33))
603  call gradm1(pstat_rutmp(1,1,20),pstat_rutmp(1,1,17), ! dvvwdy (tmp 14), dvvwdz (tmp 6)
604  $ pstat_rutmp(1,1,6),pstat_ruavg(1,1,31))
605  call gradm1(pstat_rutmp(1,1,20),pstat_rutmp(1,1,21), ! dwwwdz (tmp 9)
606  $ pstat_rutmp(1,1,9),pstat_ruavg(1,1,26))
607 
608  ! Turbulent transport tensor
609  ! Txx (new 22), Tyy (new 23), Tzz (new 24), Txy (new 25), Tyz (new 26), Txz (new 27)
610  do il = 1,6
611  itmp = (il-1)*3
612  call add4(pstat_runew(1,1,21+il),pstat_rutmp(1,1,itmp+1),
613  $ pstat_rutmp(1,1,itmp+2),pstat_rutmp(1,1,itmp+3),nvec)
614  call chsign(pstat_runew(1,1,21+il),nvec)
615  enddo
616 
617  ! Pressure strain tensor
618  ! P*dUdx (tmp 13), P*dUdy (tmp 14), P*dUdz (tmp 15), P*dVdx (tmp 16), P*dVdy (tmp 17), P*dVdz (tmp 18),
619  ! P*dWdx (tmp 19), P*dWdy (tmp 20), P*dWdz (tmp 21)
620  ! pdudx (tmp 1), pdudy (tmp 2), pdudz (tmp 3), pdvdx (tmp 4), pdvdy (tmp 5), pdvdz (tmp 6),
621  ! pdwdx (tmp 7), pdwdy (tmp 8), pdwdz (tmp 9)
622  do il = 1,9
623  call col3(pstat_rutmp(1,1,12+il),pstat_ruavg(1,1,4), ! P*d[UVW]d[xyz]
624  $ pstat_runew(1,1,il),nvec)
625  call sub3(pstat_rutmp(1,1,il),pstat_ruavg(1,1,14+il),
626  $ pstat_rutmp(1,1,12+il),nvec)
627  enddo
628 
629  ! diagonal terms
630  ! PSxx (avg 18), PSyy (avg 19), PSzz (avg 20)
631  rtmp = -2.0/rho
632  do il = 1,3
633  itmp = (il-1)*4
634  call copy(pstat_ruavg(1,1,17+il),pstat_rutmp(1,1,itmp+1),nvec)
635  call cmult(pstat_ruavg(1,1,17+il),rtmp,nvec)
636  enddo
637 
638  ! off diagonal terms
639  rtmp = -1.0/rho
640  call add3(pstat_ruavg(1,1,21),pstat_rutmp(1,1,2), ! PSxy (avg 21)
641  $ pstat_rutmp(1,1,4),nvec)
642  call cmult(pstat_ruavg(1,1,21),rtmp,nvec)
643  call add3(pstat_ruavg(1,1,22),pstat_rutmp(1,1,3), ! PSxz (avg 22)
644  $ pstat_rutmp(1,1,7),nvec)
645  call cmult(pstat_ruavg(1,1,22),rtmp,nvec)
646  call add3(pstat_ruavg(1,1,23),pstat_rutmp(1,1,6), ! PSyz (avg 23)
647  $ pstat_rutmp(1,1,8),nvec)
648  call cmult(pstat_ruavg(1,1,23),rtmp,nvec)
649 
650  ! Derivatives of the pressure-velocity products
651  ! dpudx (tmp 1), dpudy (tmp 2), dpudz (tmp 3), dpvdx (tmp 4), dpvdy (tmp 5), dpvdz (tmp 6),
652  ! dpwdx (tmp 7), dpwdy (tmp 8), dpwdz (tmp 9)
653  do il = 1,3
654  itmp = (il-1)*3
655  call gradm1(pstat_rutmp(1,1,itmp+1),pstat_rutmp(1,1,itmp+2),
656  $ pstat_rutmp(1,1,itmp+3),pstat_ruavg(1,1,11+il))
657  enddo
658 
659  ! Derivatives of the mean pressure field
660  call gradm1(pstat_rutmp(1,1,10),pstat_rutmp(1,1,11), ! dPdx (tmp 10), dPdy (tmp 11), dPdz (tmp 12)
661  $ pstat_rutmp(1,1,12),pstat_ruavg(1,1,4))
662  ! save pressure gradient for output
663  do il = 1, 3
664  call copy(pstat_pgrad(1,1,il),pstat_rutmp(1,1,9+il),nvec)
665  end do
666 
667  ! Pressure transport tensor
668  ! update dp[uvw]d[xyz]
669  do il = 1,9
670  call sub2(pstat_rutmp(1,1,il),pstat_rutmp(1,1,12+il),nvec)
671  enddo
672  do il = 1,3
673  itmp = (il-1)*3
674  do jl = 1,3
675  call subcol3(pstat_rutmp(1,1,itmp+jl),pstat_ruavg(1,1,il),
676  $ pstat_rutmp(1,1,9+jl),nvec)
677  enddo
678  enddo
679 
680  ! diagonal terms
681  ! PTxx (avg 12), PTyy (avg 13) PTzz (avg 14)
682  rtmp = -2.0/rho
683  do il=1,3
684  itmp = (il-1)*4
685  call copy(pstat_ruavg(1,1,11+il),pstat_rutmp(1,1,itmp+1),nvec)
686  call cmult(pstat_ruavg(1,1,11+il),rtmp,nvec)
687  enddo
688  ! off diagonal terms
689  rtmp = -1.0/rho
690  call add3(pstat_ruavg(1,1,15),pstat_rutmp(1,1,2), ! PTxy (avg 15)
691  $ pstat_rutmp(1,1,4),nvec)
692  call cmult(pstat_ruavg(1,1,15),rtmp,nvec)
693  call add3(pstat_ruavg(1,1,16),pstat_rutmp(1,1,3), ! PTxz (avg 16)
694  $ pstat_rutmp(1,1,7),nvec)
695  call cmult(pstat_ruavg(1,1,16),rtmp,nvec)
696  call add3(pstat_ruavg(1,1,17),pstat_rutmp(1,1,6), ! PTyz (avg 17)
697  $ pstat_rutmp(1,1,8),nvec)
698  call cmult(pstat_ruavg(1,1,17),rtmp,nvec)
699 
700  ! Production tensor
701  rtmp = -2.0
702  call vdot3(pstat_rutmp(1,1,1), ! Pxx (tmp 1)
703  $ pstat_ruavg(1,1,5),pstat_ruavg(1,1,9),pstat_ruavg(1,1,11),
704  $ pstat_runew(1,1,1),pstat_runew(1,1,2),pstat_runew(1,1,3),
705  $ nvec)
706  call cmult(pstat_rutmp(1,1,1),rtmp,nvec)
707  call vdot3(pstat_rutmp(1,1,2), ! Pyy (tmp 2)
708  $ pstat_ruavg(1,1,9),pstat_ruavg(1,1,6),pstat_ruavg(1,1,10),
709  $ pstat_runew(1,1,4),pstat_runew(1,1,5),pstat_runew(1,1,6),
710  $ nvec)
711  call cmult(pstat_rutmp(1,1,2),rtmp,nvec)
712  call vdot3(pstat_rutmp(1,1,3), ! Pzz (tmp 3)
713  $ pstat_ruavg(1,1,11),pstat_ruavg(1,1,10),pstat_ruavg(1,1,7),
714  $ pstat_runew(1,1,7),pstat_runew(1,1,8),pstat_runew(1,1,9),
715  $ nvec)
716  call cmult(pstat_rutmp(1,1,3),rtmp,nvec)
717  call vdot3(pstat_rutmp(1,1,4), ! Pxy (tmp 4)
718  $ pstat_ruavg(1,1,9),pstat_ruavg(1,1,6),pstat_ruavg(1,1,10),
719  $ pstat_runew(1,1,1),pstat_runew(1,1,2),pstat_runew(1,1,3),
720  $ nvec)
721  call vdot3(pstat_rutmp(1,1,5),
722  $ pstat_ruavg(1,1,5),pstat_ruavg(1,1,9),pstat_ruavg(1,1,11),
723  $ pstat_runew(1,1,4),pstat_runew(1,1,5),pstat_runew(1,1,6),
724  $ nvec)
725  call add2(pstat_rutmp(1,1,4),pstat_rutmp(1,1,5),nvec)
726  call chsign(pstat_rutmp(1,1,4),nvec)
727  call vdot3(pstat_rutmp(1,1,5), ! Pxz (tmp 5)
728  $ pstat_ruavg(1,1,11),pstat_ruavg(1,1,10),pstat_ruavg(1,1,7),
729  $ pstat_runew(1,1,1),pstat_runew(1,1,2),pstat_runew(1,1,3),
730  $ nvec)
731  call vdot3(pstat_rutmp(1,1,6),
732  $ pstat_ruavg(1,1,5),pstat_ruavg(1,1,9),pstat_ruavg(1,1,11),
733  $ pstat_runew(1,1,7),pstat_runew(1,1,8),pstat_runew(1,1,9),
734  $ nvec)
735  call add2(pstat_rutmp(1,1,5),pstat_rutmp(1,1,6),nvec)
736  call chsign(pstat_rutmp(1,1,5),nvec)
737  call vdot3(pstat_rutmp(1,1,6), ! Pyz (tmp 6)
738  $ pstat_ruavg(1,1,11),pstat_ruavg(1,1,10),pstat_ruavg(1,1,7),
739  $ pstat_runew(1,1,4),pstat_runew(1,1,5),pstat_runew(1,1,6),
740  $ nvec)
741  call vdot3(pstat_rutmp(1,1,7),
742  $ pstat_ruavg(1,1,9),pstat_ruavg(1,1,6),pstat_ruavg(1,1,10),
743  $ pstat_runew(1,1,7),pstat_runew(1,1,8),pstat_runew(1,1,9),
744  $ nvec)
745  call add2(pstat_rutmp(1,1,6),pstat_rutmp(1,1,7),nvec)
746  call chsign(pstat_rutmp(1,1,6),nvec)
747 
748  ! Velocity-pressure-gradient tensor
749  ! Pixx (tmp 7), Piyy (tmp 8), Pizz (tmp 9), Pixy (tmp 10), Pixz (tmp 11), Piyz (tmp 12)
750  do il=1,6
751  call sub3(pstat_rutmp(1,1,6+il),pstat_ruavg(1,1,11+il),
752  $ pstat_ruavg(1,1,17+il),nvec)
753  enddo
754 
755  ! TKE budget
756  rtmp = 0.5
757  call add4(pstat_rutmp(1,1,13),pstat_rutmp(1,1,1), ! Pk (tmp 13)
758  $ pstat_rutmp(1,1,2),pstat_rutmp(1,1,3),nvec)
759  call cmult(pstat_rutmp(1,1,13),rtmp,nvec)
760  call add4(pstat_rutmp(1,1,14),pstat_ruavg(1,1,39), ! Dk (tmp 14)
761  $ pstat_ruavg(1,1,40),pstat_ruavg(1,1,41),nvec)
762  call cmult(pstat_rutmp(1,1,14),rtmp,nvec)
763  call add4(pstat_rutmp(1,1,15),pstat_runew(1,1,22), ! Tk (tmp 15)
764  $ pstat_runew(1,1,23),pstat_runew(1,1,24),nvec)
765  call cmult(pstat_rutmp(1,1,15),rtmp,nvec)
766  call add4(pstat_rutmp(1,1,16),pstat_runew(1,1,16), ! VDk (tmp 16)
767  $ pstat_runew(1,1,17),pstat_runew(1,1,18),nvec)
768  call cmult(pstat_rutmp(1,1,16),rtmp,nvec)
769  call add4(pstat_rutmp(1,1,17),pstat_rutmp(1,1,7), ! Pik (tmp 17)
770  $ pstat_rutmp(1,1,8),pstat_rutmp(1,1,9),nvec)
771  call cmult(pstat_rutmp(1,1,17),rtmp,nvec)
772  call add4(pstat_rutmp(1,1,18),pstat_runew(1,1,10), ! Ck (tmp 18)
773  $ pstat_runew(1,1,11),pstat_runew(1,1,12),nvec)
774  call cmult(pstat_rutmp(1,1,18),rtmp,nvec)
775  call add4(pstat_rutmp(1,1,19),pstat_rutmp(1,1,13), ! Resk (tmp 19)
776  $ pstat_rutmp(1,1,14),pstat_rutmp(1,1,15),nvec)
777  call add3(pstat_rutmp(1,1,20),pstat_rutmp(1,1,16),
778  $ pstat_rutmp(1,1,17),nvec)
779  call sub2(pstat_rutmp(1,1,20),pstat_rutmp(1,1,18),nvec)
780  call add2(pstat_rutmp(1,1,19),pstat_rutmp(1,1,20),nvec)
781 
782  ! write down fields
783  call pstat3d_mfo()
784 
785  return
786  end subroutine
787 !=======================================================================
790  subroutine pstat3d_interp
791  implicit none
792 
793  include 'SIZE'
794  include 'GEOM'
795  include 'FRAMELP'
796  include 'PSTAT3D'
797 
798  ! global data structures
799  integer mid,mp,nekcomm,nekgroup,nekreal
800  common /nekmpi/ mid,mp,nekcomm,nekgroup,nekreal
801 
802  ! local variables
803  integer ifpts ! findpts flag
804  real tol ! interpolation tolerance
805  integer nt
806  integer npt_max ! max communication set
807  integer nxf, nyf, nzf ! fine mesh for bb-test
808  real bb_t ! relative size to expand bounding boxes by
809  real toldist
810  parameter(toldist = 5e-6)
811 
812  integer rcode(lhis), proc(lhis),elid(lhis)
813  real dist(lhis), rst(ldim*lhis)
814 
815  integer nfail, npass
816 
817  integer il ! loop index
818  !integer nvec ! single field length
819 
820  ! functions
821  integer iglsum
822 
823 !#define DEBUG
824 #ifdef DEBUG
825  character*3 str1, str2
826  integer iunit, ierr, jl
827  ! call number
828  integer icalld
829  save icalld
830  data icalld /0/
831 #endif
832 !-----------------------------------------------------------------------
833  ! read point position
834  call pstat3d_mfi_interp
835 
836  ! initialise interpolation tool
837  tol = 5e-13
838  nt = lx1*ly1*lz1*lelt
839  npt_max = 256
840  nxf = 2*lx1
841  nyf = 2*ly1
842  nzf = 2*lz1
843  bb_t = 0.01
844 
845  ! start interpolation tool on given mesh
846  call fgslib_findpts_setup(ifpts,nekcomm,mp,ldim,xm1,ym1,zm1,
847  & lx1,ly1,lz1,nelt,nxf,nyf,nzf,bb_t,nt,nt,npt_max,tol)
848 
849  ! identify points
850  call fgslib_findpts(ifpts,rcode,1,proc,1,elid,1,rst,ldim,dist,1,
851  & pstat_int_pts(1,1),ldim,pstat_int_pts(2,1),ldim,
852  & pstat_int_pts(ldim,1),ldim,pstat_npt)
853 
854  ! find problems with interpolation
855  nfail = 0
856  do il = 1,pstat_npt
857  ! check return code
858  if (rcode(il).eq.1) then
859  if (sqrt(dist(il)).gt.toldist) nfail = nfail + 1
860  elseif(rcode(il).eq.2) then
861  nfail = nfail + 1
862  endif
863  enddo
864  nfail = iglsum(nfail,1)
865 
866 #ifdef DEBUG
867  ! for testing
868  ! to output refinement
869  icalld = icalld+1
870  call io_file_freeid(iunit, ierr)
871  write(str1,'(i3.3)') nid
872  write(str2,'(i3.3)') icalld
873  open(unit=iunit,file='INTfpts.txt'//str1//'i'//str2)
874 
875  write(iunit,*) pstat_nptot, pstat_npt, nfail
876  do il=1,pstat_npt
877  write(iunit,*) il, proc(il), elid(il), rcode(il), dist(il),
878  $ (rst(jl+(il-1)*ldim),jl=1,ldim)
879  enddo
880 
881  close(iunit)
882 #endif
883 
884  if (nfail.gt.0) call mntr_abort(pstat_id,
885  $ 'pstat_interp: Points not mapped')
886 
887  ! Interpolate averaged fields
888  do il=1,pstat_svar
889  call fgslib_findpts_eval(ifpts,pstat_int_avg(1,il),1,
890  & rcode,1,proc,1,elid,1,rst,ndim,pstat_npt,
891  & pstat_ruavg(1,1,il))
892  enddo
893 
894 #ifdef DEBUG
895  ! for testing
896  ! to output refinement
897  call io_file_freeid(iunit, ierr)
898  write(str1,'(i3.3)') nid
899  write(str2,'(i3.3)') icalld
900  open(unit=iunit,file='INTavg.txt'//str1//'i'//str2)
901 
902  write(iunit,*) pstat_nptot, pstat_npt
903  do il=1,pstat_npt
904  write(iunit,*) il, (pstat_int_avg(il,jl),jl=1,4)
905  enddo
906 
907  close(iunit)
908 #endif
909 
910  ! Interpolate tmp fields
911  do il=1,pstat_tvar
912  call fgslib_findpts_eval(ifpts,pstat_int_tmp(1,il),1,
913  & rcode,1,proc,1,elid,1,rst,ndim,pstat_npt,
914  & pstat_rutmp(1,1,il))
915  enddo
916 
917 #ifdef DEBUG
918  ! for testing
919  ! to output refinement
920  call io_file_freeid(iunit, ierr)
921  write(str1,'(i3.3)') nid
922  write(str2,'(i3.3)') icalld
923  open(unit=iunit,file='INTtmp.txt'//str1//'i'//str2)
924 
925  write(iunit,*) pstat_nptot, pstat_npt
926  do il=1,pstat_npt
927  write(iunit,*) il, (pstat_int_tmp(il,jl),jl=1,4)
928  enddo
929 
930  close(iunit)
931 #endif
932 
933  ! Interpolate new fields
934  do il=1,pstat_dvar
935  call fgslib_findpts_eval(ifpts,pstat_int_new(1,il),1,
936  & rcode,1,proc,1,elid,1,rst,ndim,pstat_npt,
937  & pstat_runew(1,1,il))
938  enddo
939 
940 #ifdef DEBUG
941  ! for testing
942  ! to output refinement
943  call io_file_freeid(iunit, ierr)
944  write(str1,'(i3.3)') nid
945  write(str2,'(i3.3)') icalld
946  open(unit=iunit,file='INTnew.txt'//str1//'i'//str2)
947 
948  write(iunit,*) pstat_nptot, pstat_npt
949  do il=1,pstat_npt
950  write(iunit,*) il, (pstat_int_new(il,jl),jl=1,4)
951  enddo
952 
953  close(iunit)
954 #endif
955 
956  ! Interpolate pressure gradient
957  do il=1,ldim
958  call fgslib_findpts_eval(ifpts,pstat_int_pgr(1,il),1,
959  & rcode,1,proc,1,elid,1,rst,ndim,pstat_npt,
960  & pstat_pgrad(1,1,il))
961  enddo
962 
963  ! finalise interpolation tool
964  call fgslib_findpts_free(ifpts)
965 
966  ! write down interpolated values
967  call pstat3d_mfo_interp
968 
969 #undef DEBUG
970  return
971  end subroutine
972 !=======================================================================
subroutine io_file_freeid(iunit, ierr)
Get free file unit number and store max unit value.
Definition: io_tools.f:47
subroutine io_mfo_fname(fname, bname, prefix, ierr)
Generate file name according to nek rulles without opening the file.
Definition: io_tools.f:109
subroutine mntr_tmr_is_name_reg(mid, mname)
Check if timer name is registered and return its id.
Definition: mntrtmr.f:146
subroutine mntr_warn(mid, logs)
Write warning message.
Definition: mntrlog.f:803
subroutine mntr_tmr_add(mid, icount, time)
Check if timer id is registered. This operation is performed locally.
Definition: mntrtmr.f:237
subroutine mntr_mod_is_name_reg(mid, mname)
Check if module name is registered and return its id.
Definition: mntrlog.f:459
subroutine mntr_abort(mid, logs)
Abort simulation.
Definition: mntrlog.f:837
subroutine mntr_log(mid, priority, logs)
Write log message.
Definition: mntrlog.f:600
subroutine mntr_mod_reg(mid, pmid, mname, mdscr)
Register new module.
Definition: mntrlog.f:346
subroutine mntr_tmr_reg(mid, pmid, modid, mname, mdscr, ifsum)
Register new timer.
Definition: mntrtmr.f:16
subroutine pstat3d_nfield
Calculate new fileds.
Definition: pstat3D.f:354
subroutine pstat3d_mfi_interp
Read interpolation points position and redistribute them.
Definition: pstat3D_IO.f:83
subroutine pstat3d_mfo_interp
Geather data and write it down.
Definition: pstat3D_IO.f:276
subroutine pstat3d_init()
Initilise pstat module.
Definition: pstat3D.f:96
subroutine pstat3d_interp
Interpolate int the set of points.
Definition: pstat3D.f:791
subroutine pstat3d_register()
Register post processing statistics module.
Definition: pstat3D.f:11
subroutine pstat3d_main
Main interface of pstat module.
Definition: pstat3D.f:174
subroutine pstat3d_sts_avg
Read in fields and average them.
Definition: pstat3D.f:216
logical function pstat3d_is_initialised()
Check if module was initialised.
Definition: pstat3D.f:161
subroutine pstat3d_mfo
Write field data data to the file.
Definition: pstat3D_IO.f:10
subroutine rprm_rp_get(ipval, rpval, lpval, cpval, rpid, ptype)
Get runtime parameter form active section. This operation is performed locally.
Definition: rprm.f:883
subroutine rprm_rp_reg(rpid, mid, pname, pdscr, ptype, ipval, rpval, lpval, cpval)
Register new runtime parameter.
Definition: rprm.f:484
subroutine rprm_sec_set_act(ifact, rpid)
Set section's activation flag. Master value is broadcasted.
Definition: rprm.f:422
subroutine rprm_sec_reg(rpid, mid, pname, pdscr)
Register new parameter section.
Definition: rprm.f:165
subroutine mfi(fname_in, ifile)
Definition: ic.f:2355
subroutine addfid(fname, fid)
Definition: ic.f:2512
subroutine col3(a, b, c, n)
Definition: math.f:598
subroutine sub2(a, b, n)
Definition: math.f:164
subroutine icopy(a, b, n)
Definition: math.f:289
subroutine admcol3(a, b, c, d, n)
Definition: math.f:1556
subroutine add2(a, b, n)
Definition: math.f:622
subroutine add2s2(a, b, c1, n)
Definition: math.f:690
subroutine copy(a, b, n)
Definition: math.f:260
subroutine add3(a, b, c, n)
Definition: math.f:644
subroutine add4(a, b, c, d, n)
Definition: math.f:728
subroutine subcol3(a, b, c, n)
Definition: math.f:186
subroutine subcol4(a, b, c, d, n)
Definition: math.f:197
subroutine sub3(a, b, c, n)
Definition: math.f:175
subroutine cmult(a, const, n)
Definition: math.f:315
subroutine vdot3(dot, u1, u2, u3, v1, v2, v3, n)
Definition: math.f:462
subroutine rzero(a, n)
Definition: math.f:208
subroutine chsign(a, n)
Definition: math.f:305
subroutine gradm1(ux, uy, uz, u)
Definition: navier5.f:463
subroutine outpost(v1, v2, v3, vp, vt, name3)
Definition: prepost.f:1378
subroutine io_init
Definition: prepost.f:1024