VisIt/instrumentCode: psopen_visit.F90

File psopen_visit.F90, 67.3 KB (added by Jens Henrik Goebbert, 8 years ago)
Line 
1!----------------------------------------------------------------------------
2! This file is part of psOpen
3!
4! Version 2.2
5!
6! Copyright (C) 2011-2012 Jens Henrik Goebbert <jens.henrik.goebbert()rwth-aachen.de>
7!
8! psOpen is free software; you can redistribute it and/or modify
9! it under the terms of the GNU General Public License as published by
10! the Free Software Foundation; either version 2 of the License, or
11! (at your option) any later version.
12
13! This program is distributed in the hope that it will be useful,
14! but WITHOUT ANY WARRANTY; without even the implied warranty of
15! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
16! GNU General Public License for more details.
17
18! You should have received a copy of the GNU General Public License
19! along with this program; if not, write to the Free Software
20! Foundation, Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA
21!
22! Please consider mentioning the copyright owners in any publication
23! if results of psOpen are used.
24!----------------------------------------------------------------------------
25
26!========================================
27!> @addtogroup _visit_
28!! @{
29!!
30!> @file psopen_visit.F90
31!! @brief file of the visit module
32!! @details
33!! @author Jens Henrik Goebbert
34!!
35!! @}
36!========================================
37
38!!!! TODO: check https://software.intel.com/en-us/forums/intel-visual-fortran-compiler-for-windows/topic/514614 !!!
39
40#include "psopen_defines.inc"
41
42#define VISIT_SUCCESS_OR_FAILURE(err_val) \
43 if (err_val /= VISIT_OKAY) then; \
44 if (vis_mpi_proc_id == 0) write (*,*) ' FAILED with err=', err_val; \
45 return; \
46 else; \
47 if (vis_mpi_proc_id == 0) write (*,*) ' SUCCESS'; \
48 end if
49#define WRITE_RANK0_NOADVANCE(fmt_val,val,verbose_val) if(vis_mpi_proc_id == 0 .and. verbose_val <= vis_verbose_level) write(*,fmt=fmt_val,advance='no') val
50#define WRITE_RANK0(fmt_val,val,verbose_val) if(vis_mpi_proc_id == 0 .and. verbose_val <= vis_verbose_level) write(*,fmt=fmt_val) val
51#define WRITE_RANK0_ARGS2(fmt_val,val1,val2,verbose_val) if(vis_mpi_proc_id == 0 .and. verbose_val <= vis_verbose_level) write(*,fmt=fmt_val) val1,val2
52#define WRITE_RANK0_ARGS3(fmt_val,val1,val2,val3,verbose_val) if(vis_mpi_proc_id == 0 .and. verbose_val <= vis_verbose_level) write(*,fmt=fmt_val) val1,val2,val3
53
54module psopen_visit
55#ifdef USE_MPI_MODULE
56 use mpi
57 implicit none
58#else
59 implicit none
60 include 'mpif.h'
61#endif
62 include "visitfortransimV2interface.inc"
63 save
64 private
65
66 ! verbose level <= 1: show debug messages in init_visit()
67 ! verbose level <= 2: show debug messages in visit_checkstatus()
68 ! verbose level <= 3: -- not used --
69 ! verbose level <= 4: show debug messages in callback functions
70 ! verbose level <= 5: show all debug messages
71 integer, parameter, public :: vis_verbose_level = 2
72
73 integer, public :: do_visit = 0
74 integer, public :: visit_port = -1
75
76 integer, public :: vis_mpi_proc_id = -1
77 integer, public :: vis_mpi_nproc = -1
78 integer, public :: vis_mpi_comm = -1
79
80 integer, parameter, public :: SIMMODE_RUN = 1
81 integer, parameter, public :: SIMMODE_PAUSE = 2
82 integer, parameter, public :: SIMMODE_STEP = 3
83 integer, public :: sim_mode = SIMMODE_PAUSE
84
85 integer, parameter, public :: VISMODE_DISCONNECTED = 0
86 integer, parameter, public :: VISMODE_UPDATE_FAST = 1
87 integer, parameter, public :: VISMODE_UPDATE_NEVER = 2
88 integer, parameter, public :: VISMODE_UPDATE_ONCE = 3
89 integer, parameter, public :: VISMODE_UPDATE_SYNC = 4
90 integer, public :: vis_mode = VISMODE_DISCONNECTED
91
92 integer, public :: conn_freq = 1
93 integer, public :: sync_mode = -1
94 integer, public :: sync_count = -1
95
96 integer, parameter, public :: SIMCMD_NONE = 0
97 integer, parameter, public :: SIMCMD_DUMP = 1
98 integer, parameter, public :: SIMCMD_STAT = 2
99 integer, parameter, public :: SIMCMD_SYNC = 3
100 integer, public :: sim_cmd = SIMCMD_NONE
101
102 integer, parameter, public :: VISIT_COMMAND_PROCESS = 0
103 integer, parameter, public :: VISIT_COMMAND_SUCCESS = 1
104 integer, parameter, public :: VISIT_COMMAND_FAILURE = 2
105
106 logical :: visit_tracefile_onlyrank0 = .true.
107 character (len=1024) :: visit_tracefile = '' ! trace file off by default - else set 'visit_trace.out'
108 character (len=1024) :: visit_path = VISIT_F77NULLSTRING ! empty = take visit from $PATH
109 character (len=1024) :: visit_args = VISIT_F77NULLSTRING
110 !
111 ! "visit_args" will be complemented by additional args passed via the host profile (Launch Profiles -> Additional Arguments) from the VisIt GUI.
112 ! ---------------------------
113 ! -pid
114 ! this option adds process ids to debug log files - logs are not overwritten by different processes writing to a single file
115 ! -debug [level][d]
116 ! level 1-5
117 ! d "d" for additional __FILE__ and __LINE__ information in log files
118 ! -reverse_launch
119 ! -fixed-buffer-sockets
120 ! this option enables VisIt to connect to machines which expect fixed size send/recv buffers (like default settings of BlueGene)
121 !
122
123 character (len=1024) :: visit_simbase = 'psOpen'
124 character (len=1024) :: visit_simcomment = 'psOpen Simulation'
125 character (len=1024) :: visit_simpath = VISIT_F77NULLSTRING ! directory path to where simulation was started
126 character (len=1024) :: visit_siminput = VISIT_F77NULLSTRING ! path+name to simulation input file
127 character (len=1024) :: visit_simgui = VISIT_F77NULLSTRING ! xml user interface file, which VisIt can use to create a custom user interface for controlling the simulation
128
129 integer, public :: minRealIndex (3), maxRealIndex (3)
130 integer, public :: minGhosts (3), maxGhosts (3)
131
132 real, dimension (:, :, :, :), allocatable, target :: visit_vardata_vec3df_mem
133 real, dimension (:, :, :, :), pointer, public :: visit_vardata_vec3df
134
135 real, dimension (:, :, :), allocatable, target :: visit_vardata_sca3df_mem
136 real, dimension (:, :, :), pointer, public :: visit_vardata_sca3df
137
138 integer, public :: sim_iloop = -1
139 double precision, public :: sim_time = -1.d0
140
141 integer :: sim_iloop_last = -1
142 integer :: visit_cmd_serie = 0
143
144 public :: init_visit, visit_checkstatus
145 public :: visit_varupdate_sca3df, visit_varupdate_vec3df
146 public :: visit_getdata_sca3df, visit_getdata_vec3df, visit_metadata_vec3d
147
148contains
149
150 ! --------------------------------------
151 !
152 ! initialize VisIt plugin
153 !
154 ! --------------------------------------
155 subroutine init_visit (ierr)
156 use psopen_data
157 implicit none
158
159 ! function args
160 integer :: ierr
161
162 ! other vars
163 integer :: env_len_in, env_len = 1024
164 character (len=1), dimension(:), allocatable :: env
165 character (len=1024) :: visit_tracefile_wid
166
167 ierr = 0
168
169 vis_mpi_proc_id = mpi_proc_id
170 vis_mpi_nproc = mpi_nproc
171 vis_mpi_comm = mpi_mycomm
172
173 ! check if memory arrays are enabled in psopen_data.F90
174 if(.not. associated(cwork10)) then; PRTERR1('work10 not enabled in psopen_data_F90'); end if
175 if(.not. associated(cwork11)) then; PRTERR1('work11 not enabled in psopen_data_F90'); end if
176 if(.not. associated(cwork12)) then; PRTERR1('work12 not enabled in psopen_data_F90'); end if
177
178 if(.not. associated(cwork20)) then; PRTERR1('work20 not enabled in psopen_data_F90'); end if
179 if(.not. associated(cwork21)) then; PRTERR1('work21 not enabled in psopen_data_F90'); end if
180 if(.not. associated(cwork22)) then; PRTERR1('work22 not enabled in psopen_data_F90'); end if
181
182 if(.not. associated(cwork30)) then; PRTERR1('work30 not enabled in psopen_data_F90'); end if
183 if(.not. associated(cwork31)) then; PRTERR1('work31 not enabled in psopen_data_F90'); end if
184 if(.not. associated(cwork32)) then; PRTERR1('work32 not enabled in psopen_data_F90'); end if
185
186 if(.not. associated(cwork100)) then; PRTERR1('work100 not enabled in psopen_data_F90'); end if
187
188 ! tell libsim whether the simulation is parallel
189 if (vis_mpi_nproc > 1) then
190
191 ! By default, the Libsim runtime library will create a duplicate of the MPI_COMM_WORLD communicator unless you specify an alternate communicator.
192 ! visitsetmpicommunicator() was introduced in VisIt 2.9.2, but does not work before VisIt 2.10.0
193 !WRITE_RANK0_NOADVANCE('(a)', 'VisIt: set MPI communicator ...', 1)
194 !ierr = visitsetmpicommunicator(vis_mpi_comm); CHKWARN1(ierr, VISIT_OKAY)
195 !VISIT_SUCCESS_OR_FAILURE(ierr)
196
197 WRITE_RANK0_NOADVANCE('(a)', 'VisIt: set parallel mode ...', 1 )
198 ierr = visItSetParallel(1); CHKWARN1(ierr, VISIT_OKAY)
199 VISIT_SUCCESS_OR_FAILURE(ierr)
200 end if
201
202 WRITE_RANK0_NOADVANCE('(a)', 'VisIt: set parallel rank ...', 1 )
203 ierr = visItSetParallelRank(vis_mpi_proc_id); CHKWARN1(ierr, VISIT_OKAY)
204 VISIT_SUCCESS_OR_FAILURE(ierr)
205
206 ! dump trace of visit function calls for better debugging
207 if (len_trim(visit_tracefile) > 0) then
208 WRITE_RANK0_NOADVANCE('(a)', 'VisIt: open trace file ...', 1 )
209 if(visit_tracefile_onlyrank0) then
210 if (vis_mpi_proc_id == 0) then
211 ierr = visItOpenTraceFile(visit_tracefile,len_trim(visit_tracefile)); CHKWARN1(ierr, VISIT_OKAY)
212 VISIT_SUCCESS_OR_FAILURE(ierr)
213 end if
214 else
215 write(visit_tracefile_wid, fmt='(a,a,i9.9)') trim(visit_tracefile), '_', vis_mpi_proc_id
216 ierr = visItOpenTraceFile(visit_tracefile_wid,len_trim(visit_tracefile_wid)); CHKWARN1(ierr, VISIT_OKAY)
217 VISIT_SUCCESS_OR_FAILURE(ierr)
218 end if
219 end if
220
221 ! set visits binary directory
222 if (len_trim(visit_path) > 0) then
223 WRITE_RANK0_NOADVANCE('(a)', 'VisIt: set directory ...', 1 )
224 ierr = visItSetDirectory(visit_path,len_trim(visit_path)); CHKWARN1(ierr, VISIT_OKAY)
225 VISIT_SUCCESS_OR_FAILURE(ierr)
226 end if
227
228! ! set port on which visit is listening
229! if (visit_port > 0) then
230! visit_args
231! end if
232
233 ! parse command line arguments to visit
234 if (len_trim(visit_args) > 0) then
235 WRITE_RANK0_NOADVANCE('(a)', 'VisIt: set options ...', 1 )
236 ierr = visItSetOptions(visit_args,len_trim(visit_args)); CHKWARN1(ierr, VISIT_OKAY)
237 VISIT_SUCCESS_OR_FAILURE(ierr)
238 end if
239
240 ! prepare visit environment for loading plugins
241 WRITE_RANK0_NOADVANCE('(a)', 'VisIt: setup enviroment ...', 1 )
242
243 !---- best way to setup environment if VisIt >2.10.0 installed (or an older,but patched version)
244 ! get VisIt environment from rank 0 ...
245 ! (read the output of command "visit -env -engine")
246 ! ATTENTION: this needs VisIt 2.10.1 at minimum (or visitgetenv() in visitfortransimV2interface.c to be patched)
247 if (vis_mpi_proc_id == 0) then
248 env_len_in = env_len+1
249 ! loop until string is large enough to hold all characters
250 do while (env_len+1 == env_len_in) ! see comment at visItGetEnv(..) for env_len+1
251 env_len = env_len *2
252 env_len_in = env_len
253 !WRITE_RANK0_ARGS2('(a,i6)', 'env_len_in: ', env_len_in, 1)
254 if(allocated(env)) deallocate(env)
255 allocate(character(len=1) :: env(env_len), stat=ierr); CHKERRQ0 (ierr)
256 ! if max. characters of env are too small to hold full string, visItGetEnv() returns env_len, which is max. characters of env minus 1 (for "\0")
257 ierr = visItGetEnv(env,env_len); CHKWARN1(ierr, VISIT_OKAY)
258 !WRITE_RANK0_ARGS2('(a,i6)', 'env_len_out: ', env_len, 1)
259 end do
260 VISIT_SUCCESS_OR_FAILURE(ierr)
261 end if
262 call MPI_BCAST(env_len, 1, MPI_INTEGER, 0, vis_mpi_comm, ierr); CHKERRQ1(ierr, MPI_SUCCESS)
263 if(vis_mpi_proc_id /= 0) then
264 allocate(character(len=1) :: env(env_len), stat=ierr); CHKERRQ0 (ierr)
265 end if
266 ! ... and pass the VisIt environment to all other processors collectively
267 ierr = visItSetupEnv2(env,env_len); CHKWARN1(ierr, VISIT_OKAY)
268 VISIT_SUCCESS_OR_FAILURE(ierr)
269 deallocate(env)
270
271! !---- else if VisIt <=2.10.0 installed
272! ierr = visItSetupEnv(); CHKWARN1(ierr, VISIT_OKAY)
273! VISIT_SUCCESS_OR_FAILURE(ierr)
274
275 ! initialize libsim library and write .sim2 file to visit_simpath
276 if (vis_mpi_proc_id == 0) then
277 WRITE_RANK0_NOADVANCE('(a)', 'VisIt: initialization ...', 1 )
278 ierr = visItInitializeSim(visit_simbase, len_trim(visit_simbase), &
279 visit_simcomment, len_trim(visit_simcomment), &
280 visit_simpath, len_trim(visit_simpath), &
281 visit_siminput, len_trim(visit_siminput), &
282 visit_simgui, len_trim(visit_simgui), &
283 VISIT_F77NULLSTRING, VISIT_F77NULLSTRINGLEN); CHKWARN1(ierr, VISIT_OKAY)
284 VISIT_SUCCESS_OR_FAILURE(ierr)
285 end if
286
287 ! calc subdomain size of this mpi process
288 ! (values stored at x)
289 !
290 ! 0 1 2 3 4 5 6 7 8 9 10
291 ! |~x~:~x~|-x-:-x-:-x-:-x-:-x-:-x-|-x-:-x-|
292 ! ^ ^
293 ! minRealIndex maxRealIndex
294
295 ! default for no ghostcell
296 minGhosts(:) = 0
297 maxGhosts(:) = 0
298 minRealIndex(:) = 0
299 maxRealIndex(1) = isize(1)
300 maxRealIndex(2) = isize(2)
301 maxRealIndex(3) = isize(3)
302
303 ! set ghostcells (except on mesh-borders)
304 if(istart(1) /= 1) then
305 minGhosts(1) = gsize
306 minRealIndex(1) = minGhosts(1)
307 maxRealIndex(1) = minGhosts(1) +isize(1)
308 end if
309 if(istart(2) /= 1) then
310 minGhosts(2) = gsize
311 minRealIndex(2) = minGhosts(2)
312 maxRealIndex(2) = minGhosts(2) +isize(2)
313 end if
314 if(istart(3) /= 1) then
315 minGhosts(3) = gsize
316 minRealIndex(3) = minGhosts(3)
317 maxRealIndex(3) = minGhosts(3) +isize(3)
318 end if
319 if(iend(1) /= n1) maxGhosts(1) = gsize
320 if(iend(2) /= n2) maxGhosts(2) = gsize
321 if(iend(3) /= n3) maxGhosts(3) = gsize
322
323 ! allocate memory for variables
324 allocate(visit_vardata_vec3df_mem( 3, 0:(maxRealIndex(1)+maxGhosts(1)), &
325 0:(maxRealIndex(2)+maxGhosts(2)), &
326 0:(maxRealIndex(3)+maxGhosts(3))), stat=ierr); CHKERRQ0(ierr)
327 visit_vardata_vec3df => visit_vardata_vec3df_mem
328
329 allocate(visit_vardata_sca3df_mem( 0:(maxRealIndex(1)+maxGhosts(1)), &
330 0:(maxRealIndex(2)+maxGhosts(2)), &
331 0:(maxRealIndex(3)+maxGhosts(3))), stat=ierr); CHKERRQ0(ierr)
332 visit_vardata_sca3df => visit_vardata_sca3df_mem
333
334 end subroutine init_visit
335
336 ! --------------------------------------
337 !
338 ! check connection to VisIT server
339 !
340 ! --------------------------------------
341 recursive integer function visit_checkstatus(sim_iloop_in, sim_time_in)
342 implicit none
343 ! function args
344 integer, intent(in) :: sim_iloop_in
345 double precision, intent(in) :: sim_time_in
346
347 ! other vars
348 integer :: ierr
349 integer :: blocking_call
350 integer :: visitState
351 integer :: visitResult
352
353 ierr = 0
354 visit_checkstatus = 0
355
356 ! -- QUICK RETURN, IF NOTHING TO DO --------------------------------------
357
358 ! do we want to check for new status now?
359 ! - if simulation is doing nothing, we always check for new status
360 ! - if visit_checkstatus() is called recursively, we always check for new status
361 ! - if simulation is doing something, we check with the frequency of 'conn_freq'
362 if(sim_iloop_in > sim_iloop_last .and. &
363 sim_iloop_in < sim_iloop_last +conn_freq) then
364 return
365 end if
366 sim_iloop = sim_iloop_in
367 sim_time = sim_time_in
368
369 ! if single step than switch to pause
370 if (sim_mode == SIMMODE_STEP) then
371 sim_mode = SIMMODE_PAUSE
372 return ! go back to solver once
373 end if
374
375 ! -- UPDATE PLOTS, IF REQUIRED --------------------------------------------
376
377 ! update time and plots
378 if( (vis_mode /= VISMODE_DISCONNECTED) .and. & ! if we are connected
379 (vis_mode /= VISMODE_UPDATE_NEVER) .and. & ! if update enabled
380 (sim_iloop /= sim_iloop_last) ) then ! if new iteration reached
381 if( (vis_mode == VISMODE_UPDATE_FAST) .or. & ! if always update
382 (vis_mode == VISMODE_UPDATE_ONCE) .or. & ! if one single update requested
383 (vis_mode == VISMODE_UPDATE_SYNC) ) then ! if sync mode
384 ierr = visItTimestepChanged(); CHKERRQ1(ierr,VISIT_OKAY)
385 ierr = visItUpdatePlots(); CHKERRQ1(ierr,VISIT_OKAY) ! get meta data ...
386 end if
387 if(vis_mode == VISMODE_UPDATE_ONCE) then
388 vis_mode = VISMODE_UPDATE_NEVER ! go back to 'update never'
389 end if
390 end if
391 sim_iloop_last = sim_iloop
392
393 ! -- CHECK FOR NEW COMMANDS, IF REQUIRED ------------------------------------
394
395 ! detect input from proc 0 and broadcast that input to all others
396 if (vis_mpi_proc_id == 0) then
397 blocking_call = 0
398
399 ! sync simulation and visualization update ?
400 if (vis_mode == VISMODE_UPDATE_SYNC) then
401 ! wait for new commands only if called from simulation loop directly
402 if(visit_cmd_serie == 0) then
403 ! only pause after a certain number of command series
404 if(sync_count >= 0) then
405 WRITE_RANK0('(a)', 'VisIt: wait for input/sync ...', 2)
406 blocking_call = 1
407 end if
408 end if
409 ! pause simulation ?
410 else if (sim_mode == SIMMODE_PAUSE) then
411 if(visit_cmd_serie == 0) then
412 WRITE_RANK0('(a)', 'VisIt: wait for input ...', 2)
413 end if
414 blocking_call = 1
415 end if
416
417 ! check for new commands in queue
418 visitState = VisItDetectInputWithTimeout(blocking_call, 0, -1); ! blocking?, no console input
419 end if
420 call MPI_BCAST(visitState, 1, MPI_INTEGER, 0, vis_mpi_comm, ierr); CHKERRQ1(ierr, MPI_SUCCESS)
421 visit_cmd_serie = visit_cmd_serie +1
422
423 ! -- PROCESS RETURNED VISIT STATE ------------------------------------------
424
425 ! VisItDetectInput() returns with "no new input" = cmd queue empty
426 ! (visitState == 0 is never reached in blocking mode!)
427 if (visitState == 0) then
428 WRITE_RANK0('(a)', 'VisIt: nothing to do', 5)
429 visit_checkstatus = 0
430 visit_cmd_serie = 0 ! command series has ended
431
432 ! return to simulation, if number of command series has been reached
433 if(sync_count <= 0) then
434 sync_count = sync_mode ! reset sync_count
435 WRITE_RANK0('(a)', 'VisIt: return to simulation', 5)
436 return
437 else
438 sync_count = sync_count-1
439 end if
440
441 ! VisItDetectInput() returns with "VisIt is trying to connect" -> open connection/stop tstep
442 else if(visitState == 1) then
443
444 WRITE_RANK0_NOADVANCE('(a)', 'VisIt: try to connect to client ... ', 2)
445 visit_checkstatus = 1
446
447 visitResult = visitAttemptConnection()
448 if (visitResult == VISIT_OKAY) then
449 WRITE_RANK0('(a)', 'connected and listening', 2)
450 vis_mode = VISMODE_UPDATE_FAST
451 conn_freq= 1
452 else
453 WRITE_RANK0('(a)', 'connection attempt FAILED', 2)
454 vis_mode = VISMODE_DISCONNECTED
455 end if
456
457 ! VisItDetectInput() returns with "VisIt wants to tell the engine something"
458 else if (visitState == 2) then
459
460 if(visit_cmd_serie == 1) then
461 WRITE_RANK0('(a)', 'VisIt: process command', 2)
462 else
463 WRITE_RANK0_ARGS2('(a,i6)', 'VisIt process command ', visit_cmd_serie, 5)
464 end if
465 visit_checkstatus = 2
466
467 ! process VisIt command
468 visitResult = visit_processCommand ()
469
470 ! disconnect on an error or closed connection
471 if (visitResult /= VISIT_OKAY) then
472
473 visitResult = VisItDisconnect ()
474 if (visitResult == VISIT_OKAY) then
475 WRITE_RANK0('(a)', 'disconnected.', 2)
476 else
477 WRITE_RANK0('(a)', 'disconnect FAILED!', 2)
478 end if
479
480 sim_mode = SIMMODE_RUN
481 vis_mode = VISMODE_DISCONNECTED
482 conn_freq= 1
483 end if
484
485 ! VisItDetectInput() returns with "error" -> cleanup/disconnect
486 else if (visitState < 0) then
487
488 WRITE_RANK0_NOADVANCE('(a)', 'VisIt: cannot recover from error ... ', 2)
489 visit_checkstatus = -1
490
491 visitResult = VisItDisconnect()
492 if (visitResult == VISIT_OKAY) then
493 WRITE_RANK0('(a)', 'disconnected.', 2)
494 else
495 WRITE_RANK0('(a)', 'disconnect FAILED!', 2)
496 end if
497
498 vis_mode = VISMODE_DISCONNECTED
499 sim_mode = SIMMODE_RUN
500
501 end if
502
503 ! check for other waiting cmds until all processed
504 visit_checkstatus = visit_checkstatus(sim_iloop, sim_time)
505 return
506
507 end function visit_checkstatus
508
509 ! --------------------------------------
510 !
511 ! check connection to VisIT server
512 !
513 ! --------------------------------------
514 function visit_processCommand ()
515 implicit none
516
517 ! function args
518 integer :: visit_processCommand
519
520 ! other vars
521 integer :: ierr, ret, success, command
522
523 visit_processCommand = VISIT_ERROR
524
525 if (vis_mpi_proc_id == 0) then
526
527 ! get commands from visit client and broadcast them
528 ! to other mpi_procs using visitSlaveProcessCallback()
529 success = visitProcessEngineCommand ()
530 if (success == VISIT_OKAY) then
531 command = VISIT_COMMAND_SUCCESS
532 ret = VISIT_OKAY
533 else
534 command = VISIT_COMMAND_FAILURE
535 ret = VISIT_ERROR
536 end if
537
538 ! broadcast success or failure
539 ! info: the corresponding broadcast for "process"
540 ! is deep inside of visitProcessEngineCommand()
541 call MPI_BCAST (command, 1, MPI_INTEGER, 0, vis_mpi_comm, ierr); CHKERRQ1(ierr, MPI_SUCCESS)
542
543 else
544
545 command = VISIT_COMMAND_PROCESS
546 do while (command == VISIT_COMMAND_PROCESS)
547 call MPI_BCAST (command, 1, MPI_INTEGER, 0, vis_mpi_comm, ierr); CHKERRQ1(ierr, MPI_SUCCESS)
548
549 if (command == VISIT_COMMAND_PROCESS) then
550 success = visitProcessEngineCommand ()
551 else if (command == VISIT_COMMAND_SUCCESS) then
552 ret = VISIT_OKAY
553 else ! command == VISIT_COMMAND_FAILURE
554 ret = VISIT_ERROR
555 end if
556 end do
557
558 end if
559
560 visit_processCommand = ret
561 end function visit_processCommand
562
563 ! --------------------------------------
564 !
565 ! update variables as float/real(4) (write from sim-memory to visit-memory)
566 !
567 ! --------------------------------------
568 subroutine visit_varupdate_sca3df (simmem, visitmem, ghostmem)
569 use psopen_data
570#ifdef USE_P3DFFT
571 use p3dfft_ghosts
572#endif
573#ifdef USE_NB3DFFT
574 use nb3dfft_ghosts
575#endif
576 implicit none
577
578 ! function vars
579 real(kind=MP), dimension (:, :, :), pointer, intent (in) :: simmem
580 real, dimension (:, :, :), pointer, intent (out) :: visitmem
581 real(kind=MP), dimension (:, :, :), pointer :: ghostmem
582
583 ! other vars
584 integer :: i, j, k
585 double precision :: ghost_time
586 integer (kind=8) :: memsize1d
587
588 memsize1d = int(msize(1)+2*gsize,8) *int(msize(2)+2*gsize,8) *int(msize(3)+2*gsize,8)
589
590 if(gsize>0) then
591 WRITE_RANK0('(a)', 'VisIt: update scalar data (incl. ghostcells)', 5)
592 ! TODO: do a direct copy from simmem first and then
593 ! copy only ghostcell related data to ghostmem and
594 ! loop only over these ghostcells afterwards
595 ghostmem (:, :, :) = simmem (:, :, :)
596 call ghosts_update (ghostmem, 3, ghost_time)
597 do k = 0, maxRealIndex (3) + maxGhosts (3)
598 do j = 0, maxRealIndex (2) + maxGhosts (2)
599 do i = 0, maxRealIndex (1) + maxGhosts (1)
600! visitmem(i,j,k) = ghosts_ijk_debug(ghostmem, istart(1)-minGhosts(1)+i, &
601! istart(2)-minGhosts(2)+j, &
602! istart(3)-minGhosts(3)+k, memsize1d)
603 visitmem(i,j,k) = ghosts_ijk(ghostmem, istart(1)-minGhosts(1)+i, &
604 istart(2)-minGhosts(2)+j, &
605 istart(3)-minGhosts(3)+k)
606 end do
607 end do
608 end do
609 else
610 WRITE_RANK0('(a)', 'VisIt: update scalar data', 5)
611 do k = 0, maxRealIndex(3)
612 do j = 0, maxRealIndex(2)
613 do i = 0, maxRealIndex(1)
614 visitmem(i,j,k) = simmem( istart(1) +i, &
615 istart(2) +j, &
616 istart(3) +k)
617 end do
618 end do
619 end do
620 end if
621
622 end subroutine visit_varupdate_sca3df
623
624 ! --------------------------------------
625 !
626 ! update variables as float/real(4) (write from sim-memory to visit-memory)
627 !
628 ! --------------------------------------
629 subroutine visit_varupdate_vec3df (simmem1, simmem2, simmem3, visitmem, ghostmem)
630 use psopen_data
631#ifdef USE_P3DFFT
632 use p3dfft_ghosts
633#endif
634#ifdef USE_NB3DFFT
635 use nb3dfft_ghosts
636#endif
637 implicit none
638
639 ! function vars
640 real(kind=MP), dimension (:, :, :), pointer, intent (in) :: simmem1, simmem2, simmem3
641 real, dimension (:, :, :, :), pointer, intent (out) :: visitmem
642 real(kind=MP), dimension (:, :, :), pointer :: ghostmem
643
644 ! other vars
645 integer :: i, j, k
646 double precision :: ghost_time
647 integer (kind=8) :: memsize1d
648
649 memsize1d = int(msize(1)+2*gsize,8) *int(msize(2)+2*gsize,8) *int(msize(3)+2*gsize,8)
650
651 WRITE_RANK0('(a)', 'VisIt: update vector data', 5)
652
653 if(gsize > 0) then
654 ! TODO: do a direct copy from simmem first and then
655 ! copy only ghostcell related data to ghostmem and
656 ! loop only over these ghostcells afterwards
657 ghostmem (:, :, :) = simmem1 (:, :, :)
658 call ghosts_update (ghostmem, 3, ghost_time)
659 do k = 0, maxRealIndex (3) + maxGhosts (3)
660 do j = 0, maxRealIndex (2) + maxGhosts (2)
661 do i = 0, maxRealIndex (1) + maxGhosts (1)
662 visitmem(1,i,j,k) = ghosts_ijk(ghostmem, istart(1)-minGhosts(1)+i, &
663 istart(2)-minGhosts(2)+j, &
664 istart(3)-minGhosts(3)+k)
665
666! visitmem(1,i,j,k) = ghosts_ijk_debug(ghostmem, istart(1)-minGhosts(1)+i, &
667! istart(2)-minGhosts(2)+j, &
668! istart(3)-minGhosts(3)+k, memsize1d)
669! visitmem(1,i,j,k) = istart(3)-minGhosts(3)+k
670! visitmem(1,i,j,k) = 0.0
671! if( (istart(1)-minGhosts(1)+i) == 3 .and. &
672! (istart(2)-minGhosts(2)+j) == 1 .and. &
673! (istart(3)-minGhosts(3)+k) == 1) then
674! visitmem(1,i,j,k) = 1.0
675! end if
676! if(visitmem(1,i,j,k) < -1.0 .or. visitmem(1,i,j,k) > 1.0.or. visitmem(1,i,j,k) /= visitmem(1,i,j,k)) then
677! WRITE_RANK0_ARGS6('(a,i6,f,3(i6))', WRITE_ARGS9('VALUE 1-ERROR: ',vis_mpi_proc_id, visitmem(2,i,j,k), i,j,k), 5)
678! endif
679 end do
680 end do
681 end do
682 else
683 do k = 0, maxRealIndex(3)
684 do j = 0, maxRealIndex(2)
685 do i = 0, maxRealIndex(1)
686 visitmem(1,i,j,k) = simmem1( istart(1) +i, &
687 istart(2) +j, &
688 istart(3) +k)
689 end do
690 end do
691 end do
692 end if
693
694 if(gsize > 0) then
695 ! TODO: do a direct copy from simmem first and then
696 ! copy only ghostcell related data to ghostmem and
697 ! loop only over these ghostcells afterwards
698 ghostmem (:, :, :) = simmem2 (:, :, :)
699 call ghosts_update (ghostmem, 3, ghost_time)
700 do k = 0, maxRealIndex (3) + maxGhosts (3)
701 do j = 0, maxRealIndex (2) + maxGhosts (2)
702 do i = 0, maxRealIndex (1) + maxGhosts (1)
703 visitmem(2,i,j,k) = ghosts_ijk(ghostmem, istart(1)-minGhosts(1)+i, &
704 istart(2)-minGhosts(2)+j, &
705 istart(3)-minGhosts(3)+k)
706
707! visitmem(2,i,j,k) = ghosts_ijk_debug(ghostmem, istart(1)-minGhosts(1)+i, &
708! istart(2)-minGhosts(2)+j, &
709! istart(3)-minGhosts(3)+k, memsize1d)
710! visitmem(2,i,j,k) = istart(3)-minGhosts(3)+k
711! visitmem(2,i,j,k) = 0.0
712! if( (istart(1)-minGhosts(1)+i) == 1 .and. &
713! (istart(2)-minGhosts(2)+j) == 5 .and. &
714! (istart(3)-minGhosts(3)+k) == 1) then
715! visitmem(2,i,j,k) = 1.0
716! end if
717! if(visitmem(2,i,j,k) < -1.0 .or. visitmem(2,i,j,k) > 1.0 .or. visitmem(2,i,j,k) /= visitmem(2,i,j,k)) then
718! WRITE_RANK0_ARGS6('(a,i6,f,3(i6))', WRITE_ARGS9('VALUE 2-ERROR: ',vis_mpi_proc_id, visitmem(2,i,j,k), i,j,k), 5)
719! end if
720 end do
721 end do
722 end do
723 else
724 do k = 0, maxRealIndex(3)
725 do j = 0, maxRealIndex(2)
726 do i = 0, maxRealIndex(1)
727 visitmem(2,i,j,k) = simmem2( istart(1) +i, &
728 istart(2) +j, &
729 istart(3) +k)
730 end do
731 end do
732 end do
733 end if
734
735 if(gsize > 0) then
736 ! TODO: do a direct copy from simmem first and then
737 ! copy only ghostcell related data to ghostmem and
738 ! loop only over these ghostcells afterwards
739 ghostmem (:, :, :) = simmem3 (:, :, :)
740 call ghosts_update (ghostmem, 3, ghost_time)
741 do k = 0, maxRealIndex (3) + maxGhosts (3)
742 do j = 0, maxRealIndex (2) + maxGhosts (2)
743 do i = 0, maxRealIndex (1) + maxGhosts (1)
744 visitmem(3,i,j,k) = ghosts_ijk(ghostmem, istart(1)-minGhosts(1)+i, &
745 istart(2)-minGhosts(2)+j, &
746 istart(3)-minGhosts(3)+k)
747! visitmem(3,i,j,k) = ghosts_ijk_debug(ghostmem, istart(1)-minGhosts(1)+i, &
748! istart(2)-minGhosts(2)+j, &
749! istart(3)-minGhosts(3)+k, memsize1d)
750! visitmem(3,i,j,k) = istart(3)-minGhosts(3)+k
751! visitmem(3,i,j,k) = 0.0
752! if( (istart(1)-minGhosts(1)+i) == 1 .and. &
753! (istart(2)-minGhosts(2)+j) == 1 .and. &
754! (istart(3)-minGhosts(3)+k) == 8) then
755! visitmem(3,i,j,k) = 1.0
756! end if
757! if(visitmem(3,i,j,k) < -1.0 .or. visitmem(3,i,j,k) > 1.0.or. visitmem(3,i,j,k) /= visitmem(3,i,j,k)) then
758! WRITE_RANK0_ARGS6('(a,i6,f,3(i6))', WRITE_ARGS9('VALUE 3-ERROR: ',vis_mpi_proc_id, visitmem(3,i,j,k), i,j,k), 5)
759! end if
760 end do
761 end do
762 end do
763 else
764 do k = 0, maxRealIndex(3)
765 do j = 0, maxRealIndex(2)
766 do i = 0, maxRealIndex(1)
767 visitmem(3,i,j,k) = simmem3( istart(1) +i, &
768 istart(2) +j, &
769 istart(3) +k)
770 end do
771 end do
772 end do
773 end if
774
775 end subroutine visit_varupdate_vec3df
776
777 ! --------------------------------------
778 !
779 ! add mesh variables to meta data
780 !
781 ! --------------------------------------
782 subroutine visit_metadata_vec3d(md, mname, vname, ierr)
783 implicit none
784
785 ! function vars
786 integer, intent(inout) :: md
787 character(len=*) :: mname, vname
788 integer, intent(out) :: ierr
789
790 ! other vars
791 integer :: vmd
792 character(len=100) :: meshname, vecname
793
794 ! to pass it to C/VisIt, we have to copy character to stack, that somehow fixes some problems
795 meshname = repeat(' ',len(meshname))
796 vecname = repeat(' ',len(vecname))
797 meshname = trim(mname)
798 vecname = trim(vname)
799
800 ierr = 1
801
802 !meshname = 'mesh3d'
803 !vecname = 'velocity'
804 !WRITE_RANK0_ARGS3('(a,a,i6)', 'meshn : .',trim(meshname),'.', len_trim(meshname), 5)
805 !WRITE_RANK0_ARGS3('(a,a,i6)', 'vecna : .', trim(vecname),'.', len_trim(vecname), 5)
806 if (visitmdvaralloc(vmd) == VISIT_OKAY) then
807 ierr = visitmdVarSetName(vmd, vecname, len_trim(vecname)); CHKERRQ1(ierr,VISIT_OKAY)
808 ierr = visitmdVarSetMeshName(vmd, meshname, len_trim(meshname)); CHKERRQ1(ierr,VISIT_OKAY)
809 ierr = visitmdVarSetCentering(vmd, VISIT_VARCENTERING_NODE); CHKERRQ1(ierr,VISIT_OKAY)
810 ierr = visitmdVarSetType(vmd, VISIT_VARTYPE_VECTOR); CHKERRQ1(ierr,VISIT_OKAY)
811
812 ierr = visitmdSimAddVariable(md, vmd); CHKERRQ1(ierr,VISIT_OKAY)
813
814 ierr = 0
815 end if
816
817 end subroutine visit_metadata_vec3d
818
819 ! --------------------------------------
820 !
821 ! add mesh scalar-variable to meta data
822 !
823 ! --------------------------------------
824 integer function visit_getdata_sca3df(meshname, vecname, simmem, ghostmem, ierr)
825 use psopen_data ! using MP
826 implicit none
827
828 ! function vars
829 character(len=*) :: meshname, vecname
830 real(kind=MP), dimension (:, :, :), pointer, intent (in) :: simmem
831 real(kind=MP), dimension (:, :, :), pointer :: ghostmem
832 integer, intent(out) :: ierr
833
834 ! function
835 integer :: nTuples
836
837 ierr = 1
838 visit_getdata_sca3df = VISIT_INVALID_HANDLE
839
840 if (visitvardataalloc(visit_getdata_sca3df) == VISIT_OKAY) then
841 WRITE_RANK0('(a)', 'VisIt: get scalar data', 5)
842
843 call visit_varupdate_sca3df (simmem, visit_vardata_sca3df, ghostmem)
844
845 ! set variable velocity with 3 doubles for each node
846 nTuples = size (visit_vardata_sca3df, 1) * size (visit_vardata_sca3df, 2) * size (visit_vardata_sca3df, 3)
847 ierr = visitvardatasetf (visit_getdata_sca3df, VISIT_OWNER_COPY, 1, nTuples, visit_vardata_sca3df); CHKERRQ1(ierr,VISIT_OKAY)
848
849 ierr = 0
850
851 else
852 visit_getdata_sca3df = VISIT_INVALID_HANDLE
853 WRITE_RANK0('(a)', 'VisIt: get scalar data FAILED - visitVarDataAlloc() failed', 5)
854 end if
855
856 end function visit_getdata_sca3df
857
858 ! --------------------------------------
859 !
860 ! add mesh vector-variables to meta data
861 !
862 ! --------------------------------------
863 integer function visit_getdata_vec3df(meshname, vecname, simmem1, simmem2, simmem3, ghostmem, ierr)
864 use psopen_data ! using MP
865 implicit none
866
867 ! function vars
868 character(len=*) :: meshname, vecname
869 real(kind=MP), dimension (:, :, :), pointer, intent (in) :: simmem1, simmem2, simmem3
870 real(kind=MP), dimension (:, :, :), pointer :: ghostmem
871 integer, intent(out) :: ierr
872
873 ! function
874 integer :: nTuples, nComponents
875
876 ierr = 1
877 visit_getdata_vec3df = VISIT_INVALID_HANDLE
878
879 if (visitvardataalloc(visit_getdata_vec3df) == VISIT_OKAY) then
880 WRITE_RANK0('(a)', 'VisIt: get vector data', 5)
881
882 call visit_varupdate_vec3df (simmem1, simmem2, simmem3, visit_vardata_vec3df, ghostmem)
883
884 ! set variable velocity with 3 doubles for each node
885 nTuples = size (visit_vardata_vec3df, 2) * size (visit_vardata_vec3df, 3) * size (visit_vardata_vec3df, 4)
886 nComponents = size (visit_vardata_vec3df, 1)
887 ierr = visitvardatasetf (visit_getdata_vec3df, VISIT_OWNER_COPY, nComponents, nTuples, visit_vardata_vec3df); CHKERRQ1(ierr,VISIT_OKAY)
888
889 ierr = 0
890
891 else
892 visit_getdata_vec3df = VISIT_INVALID_HANDLE
893 WRITE_RANK0('(a)', 'VisIt: get vector data FAILED - visitVarDataAlloc() failed', 5)
894 end if
895
896 end function visit_getdata_vec3df
897
898end module psopen_visit
899
900!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
901!
902! These functions must be defined to satisfy the visitfortransimV2interface lib.
903!
904!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
905
906!--------------------------------------------------------------------
907! visitgetmetadata()
908!--------------------------------------------------------------------
909integer function visitgetmetadata ()
910 use psopen_visit
911 implicit none
912 include "visitfortransimV2interface.inc"
913
914 ! other vars
915 integer :: md, m1, emd, cmd, ierr
916
917 character(len=100) :: meshname, vecname, scaname
918 integer :: vmd
919
920 visitgetmetadata = VISIT_INVALID_HANDLE
921
922 WRITE_RANK0('(a)', 'VisIt: get metadata ', 4)
923
924 ! set meta data for simulation
925 if (visitmdsimalloc(md) == VISIT_OKAY) then
926
927 ierr = visitmdSimSetCycleTime (md, sim_iloop, sim_time)
928 if (sim_mode == SIMMODE_RUN) then
929 ierr = visitmdsimsetmode(md, VISIT_SIMMODE_RUNNING); CHKERRQ1(ierr,VISIT_OKAY)
930 else
931 ierr = visitmdsimsetmode(md, VISIT_SIMMODE_STOPPED); CHKERRQ1(ierr,VISIT_OKAY)
932 end if
933
934 ! Add meshes to meta data
935 if (visitmdmeshalloc(m1) == VISIT_OKAY) then
936 ierr = visitmdMeshSetName(m1, "mesh3d", 6); CHKERRQ1(ierr,VISIT_OKAY)
937 ierr = visitmdMeshSetMeshType(m1,VISIT_MESHTYPE_RECTILINEAR); CHKERRQ1(ierr,VISIT_OKAY)
938 ierr = visitmdMeshSetTopologicalDim(m1, 3); CHKERRQ1(ierr,VISIT_OKAY)
939 ierr = visitmdMeshSetSpatialDim(m1, 3); CHKERRQ1(ierr,VISIT_OKAY)
940 ierr = visitmdMeshSetNumDomains(m1, vis_mpi_nproc); CHKERRQ1(ierr,VISIT_OKAY)
941
942 !ierr = visitmdMeshSetXUnits(m1, 'pt.', 3); CHKERRQ1(ierr,VISIT_OKAY)
943 !ierr = visitmdMeshSetYUnits(m1, 'pt.', 3); CHKERRQ1(ierr,VISIT_OKAY)
944 !ierr = visitmdMeshSetZUnits(m1, 'pt.', 3); CHKERRQ1(ierr,VISIT_OKAY)
945
946 ierr = visitmdMeshSetXLabel(m1, 'X', 1); CHKERRQ1(ierr,VISIT_OKAY)
947 ierr = visitmdMeshSetYLabel(m1, 'Y', 1); CHKERRQ1(ierr,VISIT_OKAY)
948 ierr = visitmdMeshSetZLabel(m1, 'Z', 1); CHKERRQ1(ierr,VISIT_OKAY)
949
950 ierr = visitmdsimaddmesh (md, m1)
951 end if
952
953 ! Add mesh variables to meta data
954! call visit_metadata_vec3d(md, 'velocity', 'mesh3d', ierr)
955 meshname = 'mesh3d'
956 vecname = 'velocity'
957 if (visitmdvaralloc(vmd) == VISIT_OKAY) then
958 ierr = visitmdVarSetName(vmd, trim(vecname), len_trim(vecname)); CHKERRQ1(ierr,VISIT_OKAY)
959 ierr = visitmdVarSetMeshName(vmd, trim(meshname), len_trim(meshname)); CHKERRQ1(ierr,VISIT_OKAY)
960 ierr = visitmdVarSetCentering(vmd, VISIT_VARCENTERING_NODE); CHKERRQ1(ierr,VISIT_OKAY)
961 ierr = visitmdVarSetType(vmd, VISIT_VARTYPE_VECTOR); CHKERRQ1(ierr,VISIT_OKAY)
962
963 ierr = visitmdSimAddVariable(md, vmd); CHKERRQ1(ierr,VISIT_OKAY)
964
965 ierr = 0
966 end if
967 ! Add mesh variables to meta data
968! call visit_metadata_vec3d(md, 'velocity_fluc', 'mesh3d', ierr)
969 meshname = 'mesh3d'
970 vecname = 'velocity_fluc'
971 if (visitmdvaralloc(vmd) == VISIT_OKAY) then
972 ierr = visitmdVarSetName(vmd, trim(vecname), len_trim(vecname)); CHKERRQ1(ierr,VISIT_OKAY)
973 ierr = visitmdVarSetMeshName(vmd, trim(meshname), len_trim(meshname)); CHKERRQ1(ierr,VISIT_OKAY)
974 ierr = visitmdVarSetCentering(vmd, VISIT_VARCENTERING_NODE); CHKERRQ1(ierr,VISIT_OKAY)
975 ierr = visitmdVarSetType(vmd, VISIT_VARTYPE_VECTOR); CHKERRQ1(ierr,VISIT_OKAY)
976
977 ierr = visitmdSimAddVariable(md, vmd); CHKERRQ1(ierr,VISIT_OKAY)
978
979 ierr = 0
980 end if
981
982 ! Add mesh variables to meta data
983! call visit_metadata_vec3d(md, 'vorticity','mesh3d', ierr)
984 meshname = 'mesh3d'
985 vecname = 'vorticity'
986 if (visitmdvaralloc(vmd) == VISIT_OKAY) then
987 ierr = visitmdVarSetName(vmd, trim(vecname), len_trim(vecname)); CHKERRQ1(ierr,VISIT_OKAY)
988 ierr = visitmdVarSetMeshName(vmd, trim(meshname), len_trim(meshname)); CHKERRQ1(ierr,VISIT_OKAY)
989 ierr = visitmdVarSetCentering(vmd, VISIT_VARCENTERING_NODE); CHKERRQ1(ierr,VISIT_OKAY)
990 ierr = visitmdVarSetType(vmd, VISIT_VARTYPE_VECTOR); CHKERRQ1(ierr,VISIT_OKAY)
991
992 ierr = visitmdSimAddVariable(md, vmd); CHKERRQ1(ierr,VISIT_OKAY)
993
994 ierr = 0
995 end if
996
997 ! Add mesh variables to meta data
998 !call visit_metadata_sca3d(md, 'kturb','mesh3d', ierr)
999 meshname = 'mesh3d'
1000 scaname = 'k'
1001 if (visitmdvaralloc(vmd) == VISIT_OKAY) then
1002 ierr = visitmdVarSetName(vmd, trim(scaname), len_trim(scaname)); CHKERRQ1(ierr,VISIT_OKAY)
1003 ierr = visitmdVarSetMeshName(vmd, trim(meshname), len_trim(meshname)); CHKERRQ1(ierr,VISIT_OKAY)
1004 ierr = visitmdVarSetCentering(vmd, VISIT_VARCENTERING_NODE); CHKERRQ1(ierr,VISIT_OKAY)
1005 ierr = visitmdVarSetType(vmd, VISIT_VARTYPE_SCALAR); CHKERRQ1(ierr,VISIT_OKAY)
1006
1007 ierr = visitmdSimAddVariable(md, vmd); CHKERRQ1(ierr,VISIT_OKAY)
1008
1009 ierr = 0
1010 end if
1011
1012 ! Add mesh variables to meta data
1013 !call visit_metadata_sca3d(md, 'kturb','mesh3d', ierr)
1014 meshname = 'mesh3d'
1015 scaname = 'k_fluc'
1016 if (visitmdvaralloc(vmd) == VISIT_OKAY) then
1017 ierr = visitmdVarSetName(vmd, trim(scaname), len_trim(scaname)); CHKERRQ1(ierr,VISIT_OKAY)
1018 ierr = visitmdVarSetMeshName(vmd, trim(meshname), len_trim(meshname)); CHKERRQ1(ierr,VISIT_OKAY)
1019 ierr = visitmdVarSetCentering(vmd, VISIT_VARCENTERING_NODE); CHKERRQ1(ierr,VISIT_OKAY)
1020 ierr = visitmdVarSetType(vmd, VISIT_VARTYPE_SCALAR); CHKERRQ1(ierr,VISIT_OKAY)
1021
1022 ierr = visitmdSimAddVariable(md, vmd); CHKERRQ1(ierr,VISIT_OKAY)
1023
1024 ierr = 0
1025 end if
1026
1027 ! Add mesh variables to meta data
1028 !call visit_metadata_sca3d(md, 'duds','mesh3d', ierr)
1029 meshname = 'mesh3d'
1030 scaname = 'duds'
1031 if (visitmdvaralloc(vmd) == VISIT_OKAY) then
1032 ierr = visitmdVarSetName(vmd, trim(scaname), len_trim(scaname)); CHKERRQ1(ierr,VISIT_OKAY)
1033 ierr = visitmdVarSetMeshName(vmd, trim(meshname), len_trim(meshname)); CHKERRQ1(ierr,VISIT_OKAY)
1034 ierr = visitmdVarSetCentering(vmd, VISIT_VARCENTERING_NODE); CHKERRQ1(ierr,VISIT_OKAY)
1035 ierr = visitmdVarSetType(vmd, VISIT_VARTYPE_SCALAR); CHKERRQ1(ierr,VISIT_OKAY)
1036
1037 ierr = visitmdSimAddVariable(md, vmd); CHKERRQ1(ierr,VISIT_OKAY)
1038
1039 ierr = 0
1040 end if
1041
1042 ! Add mesh variables to meta data
1043 !call visit_metadata_sca3d(md, 'duds','mesh3d', ierr)
1044 meshname = 'mesh3d'
1045 scaname = 'duds_fluc'
1046 if (visitmdvaralloc(vmd) == VISIT_OKAY) then
1047 ierr = visitmdVarSetName(vmd, trim(scaname), len_trim(scaname)); CHKERRQ1(ierr,VISIT_OKAY)
1048 ierr = visitmdVarSetMeshName(vmd, trim(meshname), len_trim(meshname)); CHKERRQ1(ierr,VISIT_OKAY)
1049 ierr = visitmdVarSetCentering(vmd, VISIT_VARCENTERING_NODE); CHKERRQ1(ierr,VISIT_OKAY)
1050 ierr = visitmdVarSetType(vmd, VISIT_VARTYPE_SCALAR); CHKERRQ1(ierr,VISIT_OKAY)
1051
1052 ierr = visitmdSimAddVariable(md, vmd); CHKERRQ1(ierr,VISIT_OKAY)
1053
1054 ierr = 0
1055 end if
1056
1057 ! Add mesh variables to meta data
1058 !call visit_metadata_sca3d(md, 'lambda2','mesh3d', ierr)
1059 meshname = 'mesh3d'
1060 scaname = 'lambda2'
1061 if (visitmdvaralloc(vmd) == VISIT_OKAY) then
1062 ierr = visitmdVarSetName(vmd, trim(scaname), len_trim(scaname)); CHKERRQ1(ierr,VISIT_OKAY)
1063 ierr = visitmdVarSetMeshName(vmd, trim(meshname), len_trim(meshname)); CHKERRQ1(ierr,VISIT_OKAY)
1064 ierr = visitmdVarSetCentering(vmd, VISIT_VARCENTERING_NODE); CHKERRQ1(ierr,VISIT_OKAY)
1065 ierr = visitmdVarSetType(vmd, VISIT_VARTYPE_SCALAR); CHKERRQ1(ierr,VISIT_OKAY)
1066
1067 ierr = visitmdSimAddVariable(md, vmd); CHKERRQ1(ierr,VISIT_OKAY)
1068
1069 ierr = 0
1070 end if
1071
1072 ! Add expressions to meta data
1073 if (visitmdexpralloc(emd) == VISIT_OKAY) then
1074 ierr = visitmdExprSetName(emd, "nid", 3); CHKERRQ1(ierr,VISIT_OKAY)
1075 ierr = visitmdExprSetDefinition(emd, "nodeid(mesh3d)", 14); CHKERRQ1(ierr,VISIT_OKAY)
1076 ierr = visitmdExprSetType(emd, VISIT_VARTYPE_SCALAR); CHKERRQ1(ierr,VISIT_OKAY)
1077
1078 ierr = visitmdSimAddExpression(md, emd); CHKERRQ1(ierr,VISIT_OKAY)
1079 end if
1080
1081 ! Add simulation commands to meta data
1082
1083 ! set simulation mode
1084 if (visitmdcmdalloc(cmd) == VISIT_OKAY) then
1085 ierr = visitmdcmdsetname(cmd, "sim. run", 8); CHKERRQ1(ierr,VISIT_OKAY)
1086 ierr = visitmdsimaddgenericcommand(md, cmd); CHKERRQ1(ierr,VISIT_OKAY)
1087 end if
1088 if (visitmdcmdalloc(cmd) == VISIT_OKAY) then
1089 ierr = visitmdcmdsetname(cmd, "sim. pause", 10); CHKERRQ1(ierr,VISIT_OKAY)
1090 ierr = visitmdsimaddgenericcommand(md, cmd); CHKERRQ1(ierr,VISIT_OKAY)
1091 end if
1092 if (visitmdcmdalloc(cmd) == VISIT_OKAY) then
1093 ierr = visitmdcmdsetname(cmd, "sim. step", 9); CHKERRQ1(ierr,VISIT_OKAY)
1094 ierr = visitmdsimaddgenericcommand(md, cmd); CHKERRQ1(ierr,VISIT_OKAY)
1095 end if
1096
1097 ! set connection mode
1098 if (visitmdcmdalloc(cmd) == VISIT_OKAY) then
1099 ierr = visitmdcmdsetname(cmd, "conn. freq. 1", 13); CHKERRQ1(ierr,VISIT_OKAY)
1100 ierr = visitmdsimaddgenericcommand(md, cmd); CHKERRQ1(ierr,VISIT_OKAY)
1101 end if
1102 if (visitmdcmdalloc(cmd) == VISIT_OKAY) then
1103 ierr = visitmdcmdsetname(cmd, "conn. freq. 5", 13); CHKERRQ1(ierr,VISIT_OKAY)
1104 ierr = visitmdsimaddgenericcommand(md, cmd); CHKERRQ1(ierr,VISIT_OKAY)
1105 end if
1106 if (visitmdcmdalloc(cmd) == VISIT_OKAY) then
1107 ierr = visitmdcmdsetname(cmd, "conn. freq. 10", 14); CHKERRQ1(ierr,VISIT_OKAY)
1108 ierr = visitmdsimaddgenericcommand(md, cmd); CHKERRQ1(ierr,VISIT_OKAY)
1109 end if
1110
1111 ! set visualization mode
1112 if (visitmdcmdalloc(cmd) == VISIT_OKAY) then
1113 ierr = visitmdcmdsetname(cmd, "update fast", 11); CHKERRQ1(ierr,VISIT_OKAY)
1114 ierr = visitmdsimaddgenericcommand(md, cmd); CHKERRQ1(ierr,VISIT_OKAY)
1115 end if
1116 if (visitmdcmdalloc(cmd) == VISIT_OKAY) then
1117 ierr = visitmdcmdsetname(cmd, "update never", 12); CHKERRQ1(ierr,VISIT_OKAY)
1118 ierr = visitmdsimaddgenericcommand(md, cmd); CHKERRQ1(ierr,VISIT_OKAY)
1119 end if
1120 if (visitmdcmdalloc(cmd) == VISIT_OKAY) then
1121 ierr = visitmdcmdsetname(cmd, "update once", 11); CHKERRQ1(ierr,VISIT_OKAY)
1122 ierr = visitmdsimaddgenericcommand(md, cmd); CHKERRQ1(ierr,VISIT_OKAY)
1123 end if
1124
1125 if (visitmdcmdalloc(cmd) == VISIT_OKAY) then
1126 ierr = visitmdcmdsetname(cmd, "update sync 1", 13); CHKERRQ1(ierr,VISIT_OKAY)
1127 ierr = visitmdsimaddgenericcommand(md, cmd); CHKERRQ1(ierr,VISIT_OKAY)
1128 end if
1129 if (visitmdcmdalloc(cmd) == VISIT_OKAY) then
1130 ierr = visitmdcmdsetname(cmd, "update sync (movie)", 19); CHKERRQ1(ierr,VISIT_OKAY)
1131 ierr = visitmdsimaddgenericcommand(md, cmd); CHKERRQ1(ierr,VISIT_OKAY)
1132 end if
1133 if (visitmdcmdalloc(cmd) == VISIT_OKAY) then
1134 ierr = visitmdcmdsetname(cmd, "---", 3); CHKERRQ1(ierr,VISIT_OKAY)
1135 ierr = visitmdsimaddgenericcommand(md, cmd); CHKERRQ1(ierr,VISIT_OKAY)
1136 end if
1137
1138 ! set special commands
1139 if (visitmdcmdalloc(cmd) == VISIT_OKAY) then
1140 ierr = visitmdcmdsetname(cmd, "dump", 4); CHKERRQ1(ierr,VISIT_OKAY)
1141 ierr = visitmdsimaddgenericcommand(md, cmd); CHKERRQ1(ierr,VISIT_OKAY)
1142 end if
1143 if (visitmdcmdalloc(cmd) == VISIT_OKAY) then
1144 ierr = visitmdcmdsetname(cmd, "stat", 4); CHKERRQ1(ierr,VISIT_OKAY)
1145 ierr = visitmdsimaddgenericcommand(md, cmd); CHKERRQ1(ierr,VISIT_OKAY)
1146 end if
1147 end if
1148
1149 visitgetmetadata = md
1150
1151end function visitgetmetadata
1152
1153!-----------------------------------------------------------------
1154! visitcommandcallback
1155!-----------------------------------------------------------------
1156subroutine visitcommandcallback (cmd, lcmd, args, largs)
1157 use psopen_visit
1158 implicit none
1159 include "visitfortransimV2interface.inc"
1160
1161 ! function args
1162 integer :: lcmd, largs
1163 character (len=1) :: cmd (lcmd)
1164 character (len=1) :: args (largs)
1165
1166 ! other vars
1167 integer :: result
1168
1169 WRITE_RANK0_NOADVANCE('(a)', 'VisIt: callback command ... ', 4)
1170
1171 ! do not get confused if disconnection and command call overlaps
1172 if(vis_mode == VISMODE_DISCONNECTED) return
1173
1174 ! Handle the commands that we define in visitgetmetadata.
1175
1176 ! set simulation mode
1177 if (visitstrcmp(cmd, lcmd, "sim. run", 8) == 0) then
1178 WRITE_RANK0('(a)', 'sim run', 4)
1179 sim_mode = SIMMODE_RUN
1180 vis_mode = VISMODE_UPDATE_FAST
1181 !conn_freq= do not change
1182 else if (visitstrcmp(cmd, lcmd, "sim. pause", 10) == 0) then
1183 WRITE_RANK0('(a)', 'sim pause', 4)
1184 sim_mode = SIMMODE_PAUSE
1185 vis_mode = VISMODE_UPDATE_FAST
1186 !conn_freq= do not change
1187 else if (visitstrcmp(cmd, lcmd, "sim. step", 9) == 0) then
1188 WRITE_RANK0('(a)', 'sim step', 4)
1189 sim_mode = SIMMODE_STEP
1190 vis_mode = VISMODE_UPDATE_FAST
1191 !conn_freq= do not change
1192
1193 ! set connection frequency
1194 else if (visitstrcmp(cmd, lcmd, "conn. freq. 1", 13) == 0) then
1195 WRITE_RANK0('(a)', 'conn. freq. 1', 4)
1196 !sim_mode = do not change
1197 !vis_mode = do not change
1198 conn_freq = 1
1199 else if (visitstrcmp(cmd, lcmd, "conn. freq. 5", 13) == 0) then
1200 WRITE_RANK0('(a)', 'conn. freq. 5', 4)
1201 !sim_mode = do not change
1202 !vis_mode = do not change
1203 conn_freq = 5
1204 else if (visitstrcmp(cmd, lcmd, "conn. freq. 10", 14) == 0) then
1205 WRITE_RANK0('(a)', 'conn. freq. 10', 4)
1206 !sim_mode = do not change
1207 !vis_mode = do not change
1208 conn_freq = 10
1209
1210 ! set visualization mode
1211 else if (visitstrcmp(cmd, lcmd, "update fast", 11) == 0) then
1212 WRITE_RANK0('(a)', 'update fast', 4)
1213 !sim_mode = do not change
1214 vis_mode = VISMODE_UPDATE_FAST
1215 !conn_freq = do not change
1216 else if (visitstrcmp(cmd, lcmd, "update never", 12) == 0) then
1217 WRITE_RANK0('(a)', 'update never', 4)
1218 !sim_mode = do not change
1219 vis_mode = VISMODE_UPDATE_NEVER
1220 !conn_freq = do not change
1221 else if (visitstrcmp(cmd, lcmd, "update once", 11) == 0) then
1222 WRITE_RANK0('(a)', 'update once', 4)
1223 !sim_mode = do not change
1224 vis_mode = VISMODE_UPDATE_ONCE
1225 !conn_freq = do not change
1226
1227 ! for sync we have to make sure, to exit recursive call of visit_checkstatus()
1228 ! - disable blocking call to detect input
1229 ! - enable blocking call after visit_checkstatus() is called from simulation loop again
1230 else if (visitstrcmp(cmd, lcmd, "update sync 1", 13) == 0) then
1231 WRITE_RANK0('(a)', 'update sync 1', 4)
1232 vis_mode = VISMODE_UPDATE_SYNC
1233 sync_mode = 1
1234 sync_count = sync_mode
1235 else if (visitstrcmp(cmd, lcmd, "update sync (movie)", 19) == 0) then
1236 WRITE_RANK0('(a)', 'update sync (movie)', 4)
1237 vis_mode = VISMODE_UPDATE_SYNC
1238 sync_mode = 4
1239 sync_count = sync_mode
1240 else if (visitstrcmp(cmd, lcmd, "---", 3) == 0) then
1241 WRITE_RANK0('(a)', '---', 4)
1242 sync_mode = 0
1243 sync_count = sync_mode
1244
1245 else if (visitstrcmp(cmd, lcmd, "dump", 4) == 0) then
1246 WRITE_RANK0('(a)', 'dump', 4)
1247 sim_cmd = SIMCMD_DUMP
1248 else if (visitstrcmp(cmd, lcmd, "stat", 4) == 0) then
1249 WRITE_RANK0('(a)', 'stat', 4)
1250 sim_cmd = SIMCMD_STAT
1251 else
1252 WRITE_RANK0_ARGS2(*, 'unknown = ',cmd, 4)
1253 end if
1254
1255end subroutine visitcommandcallback
1256
1257!---------------------------------------------------------------------------
1258! visitgetmesh
1259!---------------------------------------------------------------------------
1260integer function visitgetmesh (domain, name, lname)
1261 use psopen_data
1262 use psopen_visit
1263 implicit none
1264 include "visitfortransimV2interface.inc"
1265
1266 ! function args
1267 integer domain, lname
1268 character (len=1) name (lname)
1269
1270 ! other vars
1271 integer :: h, hx, hy, hz, ierr
1272 real :: cx (0:isize(1)+2*gsize)
1273 real :: cy (0:isize(2)+2*gsize)
1274 real :: cz (0:isize(3)+2*gsize)
1275 integer :: i, j, k
1276 integer :: baseIndex (3)
1277 integer :: xfac, yfac, zfac
1278 real(kind=MP) :: xshift, yshift, zshift
1279
1280 WRITE_RANK0_ARGS2(*, 'VisIt: get mesh ',name, 4)
1281
1282 !------------------------------------------------------!
1283 ! domain = sub-domain(1) + sub-domain(2) !
1284 ! (values stored at x = cell-centered) !
1285 ! |-x-:-x-:-x-:-x-:-x-:-x-:-x-:-x-| !
1286 ! !
1287 ! 0 1 2 3 4 5 (cell ids) !
1288 ! |-x-:-x-:-x-:-x-|~x~:~x~| sub-domain 1 !
1289 ! 0 1 2 3 4 5 6 (coordinates) !
1290 ! ^ ^ !
1291 ! minRealIndex maxRealIndex !
1292 ! !
1293 ! 0 1 2 3 4 5 (cell ids) !
1294 ! |~x~:~x~|-x-:-x-:-x-:-x-| sub-domain 2 !
1295 ! 0 1 2 3 4 5 6 (coordinates) !
1296 ! ^ ^ !
1297 ! minRealIndex maxRealIndex !
1298 ! !
1299 ! :-x-: grid cell !
1300 ! :~x~: ghost cell !
1301 !------------------------------------------------------!
1302
1303 visitgetmesh = VISIT_INVALID_HANDLE
1304 h = VISIT_INVALID_HANDLE
1305 hx = VISIT_INVALID_HANDLE
1306 hy = VISIT_INVALID_HANDLE
1307 hz = VISIT_INVALID_HANDLE
1308
1309 ! set base index
1310 ! baseIndex is an offset in X,Y,Z that will be added to your mesh’s zone numbers and node numbers (=cell number) when VisIt displays information about your mesh.
1311 ! You can leave these values set at zero. However, when you want to create a multi-domain mesh that has global zone and node numbers (=cell number),
1312 ! you should set the values for baseIndex
1313 baseIndex (1) = istart (1) - 1
1314 baseIndex (2) = istart (2) - 1
1315 baseIndex (3) = istart (3) - 1
1316
1317 !WRITE_RANK0('(a)', 'CX corr---------------------', 4)
1318 ! set to ascending coords and start with 0.0
1319 xfac = 1
1320 xshift = 0._MP
1321 if(coordX(1) > coordX(2)) xfac = -1
1322 if(coordX(1) /= 0._MP) xshift = -coordX(1)*xfac
1323 ! calc coords shown in VisIt
1324 do i = 0, maxRealIndex(1) + maxGhosts(1)
1325 cx (i) = (coordX (istart(1)-minGhosts(1)+i) *xfac) +xshift
1326 !if(vis_mpi_proc_id == 0) write(*,*) i, cx(i)
1327 end do
1328
1329 !WRITE_RANK0('(a)', 'CY corr---------------------', 4)
1330 ! set to ascending coords and start with 0.0
1331 yfac = 1
1332 yshift = 0._MP
1333 if(coordY(1) > coordY(2)) yfac = -1
1334 if(coordY(1) /= 0._MP) yshift = -coordY(1)*yfac
1335 ! calc coords shown in VisIt
1336 do j = 0, maxRealIndex(2) + maxGhosts(2)
1337 cy (j) = (coordY (istart(2)-minGhosts(2)+j) *yfac) +yshift
1338 ! if(vis_mpi_proc_id == 0) write(*,*) j, cy(j)
1339 end do
1340
1341 !WRITE_RANK0('(a)', 'CZ corr---------------------', 4)
1342 ! set to ascending coords and start with 0.0
1343 zfac = 1
1344 zshift = 0._MP
1345 if(coordZ(1) > coordZ(2)) zfac = -1
1346 if(coordZ(1) /= 0._MP) zshift = -coordZ(1)*zfac
1347 ! calc coords shown in VisIt
1348 do k = 0, maxRealIndex(3) + maxGhosts(3)
1349 cz (k) = (coordZ (istart(3)-minGhosts(3)+k) *zfac) +zshift
1350 ! if(vis_mpi_proc_id == 0) write(*,*) k, cz(k)
1351 end do
1352
1353 if (visitstrcmp(name, lname, "mesh3d", 6) == 0) then
1354
1355 ! allocate a rectilinear mesh
1356 if (visitrectmeshalloc(h) == VISIT_OKAY) then
1357 WRITE_RANK0('(a)', 'VisIt: mesh3d set coords', 4)
1358
1359 ! set the coordinates
1360 ierr = visitvardataalloc(hx); CHKERRQ1(ierr,VISIT_OKAY)
1361 ierr = visitvardataalloc(hy); CHKERRQ1(ierr,VISIT_OKAY)
1362 ierr = visitvardataalloc(hz); CHKERRQ1(ierr,VISIT_OKAY)
1363 ! attention: visitvardatasetf(<handle>, <flag>, <no.components per tuple>, <no.tuples>, <memory>)
1364 !!! NOT <start in array>,<end in array> !!!
1365 ierr = visitvardatasetf(hx,VISIT_OWNER_COPY,1,maxRealIndex(1) +maxGhosts(1) +1,cx); CHKERRQ1(ierr,VISIT_OKAY)
1366 ierr = visitvardatasetf(hy,VISIT_OWNER_COPY,1,maxRealIndex(2) +maxGhosts(2) +1,cy); CHKERRQ1(ierr,VISIT_OKAY)
1367 ierr = visitvardatasetf(hz,VISIT_OWNER_COPY,1,maxRealIndex(3) +maxGhosts(3) +1,cz); CHKERRQ1(ierr,VISIT_OKAY)
1368
1369 ierr = visitrectmeshsetcoordsxyz(h, hx, hy, hz); CHKERRQ1(ierr,VISIT_OKAY)
1370
1371 ! set the real indices of the mesh of this mpi process (which tells VisIt about ghostcells)
1372 ierr = visitrectmeshsetrealindices(h, minRealIndex, maxRealIndex); CHKERRQ1(ierr,VISIT_OKAY)
1373
1374 ! set the base index of the mesh of this mpi process
1375 ierr = visitrectmeshsetbaseindex(h, baseIndex); CHKERRQ1(ierr,VISIT_OKAY)
1376
1377 else
1378 WRITE_RANK0('(a)', 'VisIt: mesh3d set coords ... FAILED', 4)
1379 end if
1380
1381 end if
1382
1383 visitgetmesh = h
1384
1385end function visitgetmesh
1386
1387!---------------------------------------------------------------------------
1388! visitgetvariable
1389!---------------------------------------------------------------------------
1390integer function visitgetvariable (domain, name, lname)
1391 use psopen_data
1392 use psopen_visit
1393 use psopen_fft
1394 use psopen_utils_flow
1395 use psopen_utils_stream
1396 use psopen_utils_vortex
1397 implicit none
1398 include "visitfortransimV2interface.inc"
1399
1400 ! function args
1401 integer domain, lname
1402 character (len=1) name (lname)
1403
1404 ! other vars
1405 integer :: h
1406 integer :: k, ierr
1407
1408 WRITE_RANK0_ARGS2(*, 'VisIt: get variable ',name, 4)
1409
1410 h = VISIT_INVALID_HANDLE
1411 visitgetvariable = VISIT_INVALID_HANDLE
1412
1413 ! velocity
1414 if (visitstrcmp(name, lname, "velocity", 8) == 0) then
1415 call fft_c2r(cu_arch, rwork10)
1416 call fft_c2r(cv_arch, rwork11)
1417 call fft_c2r(cw_arch, rwork12)
1418 visitgetvariable = visit_getdata_vec3df('mesh3d', 'velocity', rwork10, rwork11, rwork12, rwork100, ierr)
1419
1420 ! velocity fuctuations
1421 else &
1422 if (visitstrcmp(name, lname, "velocity_fluc", 13) == 0) then
1423 call fft_c2r(cu_arch, rwork10)
1424 call fft_c2r(cv_arch, rwork11)
1425 call fft_c2r(cw_arch, rwork12)
1426 if(.not. read_u_zprofile) then
1427 call calc_zlayer_mean(rwork10, u_zprofile, ierr)
1428 endif
1429 do k=istart(3), iend(3)
1430 rwork10(:,:,k) = rwork10(:,:,k) -u_zprofile(k)
1431 end do
1432 visitgetvariable = visit_getdata_vec3df('mesh3d', 'velocity_fluc', rwork10, rwork11, rwork12, rwork100, ierr)
1433
1434 ! vorticity
1435 else &
1436 if (visitstrcmp(name, lname, "vorticity", 9) == 0) then
1437 call calc_vorticity( cu_arch, cv_arch, cw_arch, cwork10, cwork11, cwork12, cwork100)
1438 call fft_c2r(cwork10, rwork10)
1439 call fft_c2r(cwork11, rwork11)
1440 call fft_c2r(cwork12, rwork12)
1441 visitgetvariable = visit_getdata_vec3df('mesh3d', 'vorticity', rwork10, rwork11, rwork12, rwork100, ierr)
1442
1443 ! vorticity fluctuations
1444 else &
1445 if (visitstrcmp(name, lname, "vorticity_fluc", 14) == 0) then
1446 call calc_vorticity( cu_arch, cv_arch, cw_arch, cwork10, cwork11, cwork12, cwork100)
1447 call fft_c2r(cwork10, rwork10)
1448 call fft_c2r(cwork11, rwork11)
1449 call fft_c2r(cwork12, rwork12)
1450 visitgetvariable = visit_getdata_vec3df('mesh3d', 'vorticity_fluc', rwork10, rwork11, rwork12, rwork100, ierr)
1451
1452 ! kinetic energy
1453 else &
1454 if (visitstrcmp(name, lname, "k", 1) == 0) then
1455 call fft_c2r(cu_arch, rwork10)
1456 call fft_c2r(cv_arch, rwork11)
1457 call fft_c2r(cw_arch, rwork12)
1458 call calc_kturb( rwork10, rwork11, rwork12, rwork20, ierr)
1459 visitgetvariable = visit_getdata_sca3df('mesh3d', 'k', rwork20, rwork100, ierr)
1460
1461 ! kinetic energy fluctuations
1462 else &
1463 if (visitstrcmp(name, lname, "k_fluc", 6) == 0) then
1464 call fft_c2r(cu_arch, rwork10)
1465 call fft_c2r(cv_arch, rwork11)
1466 call fft_c2r(cw_arch, rwork12)
1467 if(.not. read_u_zprofile) then
1468 call calc_zlayer_mean(rwork10, u_zprofile, ierr)
1469 endif
1470 do k=istart(3), iend(3)
1471 rwork10(:,:,k) = rwork10(:,:,k) -u_zprofile(k)
1472 end do
1473 call calc_kturb( rwork10, rwork11, rwork12, rwork20, ierr)
1474 visitgetvariable = visit_getdata_sca3df('mesh3d', 'k_fluc', rwork20, rwork100, ierr)
1475
1476 ! du/ds
1477 else &
1478 if (visitstrcmp(name, lname, "duds", 4) == 0) then
1479 call fft_c2r(cu_arch, rwork10)
1480 call fft_c2r(cv_arch, rwork11)
1481 call fft_c2r(cw_arch, rwork12)
1482 call calc_streamline_duds( rwork10, rwork11, rwork12, rwork32, &
1483 rwork20, rwork21, rwork22, &
1484 cwork30, cwork31, rwork100, ierr)
1485 visitgetvariable = visit_getdata_sca3df('mesh3d', 'duds', rwork32, rwork100, ierr)
1486
1487 ! du/ds fluctuations
1488 else &
1489 if (visitstrcmp(name, lname, "duds_fluc", 9) == 0) then
1490 call fft_c2r(cu_arch, rwork10)
1491 call fft_c2r(cv_arch, rwork11)
1492 call fft_c2r(cw_arch, rwork12)
1493 if(.not. read_u_zprofile) then
1494 call calc_zlayer_mean(rwork10, u_zprofile, ierr)
1495 endif
1496 do k=istart(3), iend(3)
1497 rwork10(:,:,k) = rwork10(:,:,k) -u_zprofile(k)
1498 end do
1499 call calc_streamline_duds( rwork10, rwork11, rwork12, rwork32, &
1500 rwork20, rwork21, rwork22, &
1501 cwork30, cwork31, rwork100, ierr)
1502 visitgetvariable = visit_getdata_sca3df('mesh3d', 'duds_fluc', rwork32, rwork100, ierr)
1503
1504 ! lambda2
1505 else &
1506 if (visitstrcmp(name, lname, "lambda2", 7) == 0) then
1507 call calc_vortex_lambda123(cu_arch, cv_arch, cw_arch, &
1508 cwork10, cwork11, cwork12, cwork20, cwork21, cwork22, &
1509 rwork10, rwork11, rwork12, rwork20, rwork21, rwork22, &
1510 rwork100, rwork101, rwork102, rwork110, rwork111, rwork112, &
1511 cwork120, cwork121, cwork122, &
1512 rwork120, rwork121, rwork122, & ! eigenvalues 1,2,3
1513 ierr)
1514 visitgetvariable = visit_getdata_sca3df('mesh3d', 'lambda2', rwork121, rwork100, ierr)
1515
1516 else
1517 WRITE_RANK0_ARGS2(*, 'VisIt: variable UNKNOWN - ',name, 4)
1518 endif
1519
1520 WRITE_RANK0('(a)', 'VisIt: get variable ... FINISHED', 4)
1521
1522end function visitgetvariable
1523
1524!---------------------------------------------------------------------------
1525! visitgetmixedvariable
1526!---------------------------------------------------------------------------
1527integer function visitgetmixedvariable (domain, name, lname)
1528 use psopen_visit
1529 implicit none
1530 include "visitfortransimV2interface.inc"
1531
1532 ! function args
1533 integer domain, lname
1534 character (len=1) name (lname)
1535
1536 ! other vars
1537 integer :: h
1538 integer :: k, ierr
1539
1540 WRITE_RANK0_ARGS2(*, 'VisIt: get mixed variable ',name, 4)
1541
1542 h = VISIT_INVALID_HANDLE
1543 visitgetmixedvariable = VISIT_INVALID_HANDLE
1544
1545 WRITE_RANK0_ARGS2(*, 'VisIt: mixed variable UNKNOWN - ',name, 4)
1546
1547 WRITE_RANK0('(a)', 'VisIt: get mixed variable ... FINISHED', 4)
1548
1549end function visitgetmixedvariable
1550
1551!---------------------------------------------------------------------------
1552! visitbroadcastintfunction
1553!---------------------------------------------------------------------------
1554integer function visitbroadcastintfunction (value, sender)
1555 use psopen_visit
1556#ifdef USE_MPI_MODULE
1557 use mpi
1558 implicit none
1559#else
1560 implicit none
1561 include 'mpif.h'
1562#endif
1563 include "visitfortransimV2interface.inc"
1564
1565 ! function vars
1566 integer :: value, sender
1567
1568 ! other vars
1569 integer :: ierr
1570
1571 call MPI_BCAST(value, 1, MPI_INTEGER, sender, vis_mpi_comm, ierr); CHKERRQ1(ierr, MPI_SUCCESS)
1572 visitbroadcastintfunction = ierr
1573
1574end function visitbroadcastintfunction
1575
1576!---------------------------------------------------------------------------
1577! visitbroadcaststringfunction
1578!---------------------------------------------------------------------------
1579integer function visitbroadcaststringfunction (str, lstr, sender)
1580 use psopen_visit
1581#ifdef USE_MPI_MODULE
1582 use mpi
1583 implicit none
1584#else
1585 implicit none
1586 include 'mpif.h'
1587#endif
1588 include "visitfortransimV2interface.inc"
1589
1590 ! function vars
1591 integer :: lstr, sender
1592 character (len=1) :: str (lstr)
1593
1594 ! other vars
1595 integer :: ierr
1596
1597 call MPI_BCAST(str, lstr, MPI_CHARACTER, sender, vis_mpi_comm, ierr); CHKERRQ1(ierr, MPI_SUCCESS)
1598 visitbroadcaststringfunction = ierr
1599
1600end function visitbroadcaststringfunction
1601
1602!---------------------------------------------------------------------------
1603! visitslaveprocesscallback
1604!---------------------------------------------------------------------------
1605subroutine visitslaveprocesscallback ()
1606 use psopen_visit
1607#ifdef USE_MPI_MODULE
1608 use mpi
1609 implicit none
1610#else
1611 implicit none
1612 include 'mpif.h'
1613#endif
1614
1615 ! other vars
1616 integer :: c, ierr
1617
1618 c = VISIT_COMMAND_PROCESS
1619 call MPI_BCAST(c, 1, MPI_INTEGER, 0, vis_mpi_comm, ierr); CHKERRQ1(ierr, MPI_SUCCESS)
1620end subroutine visitslaveprocesscallback
1621
1622!---------------------------------------------------------------------------
1623! visitactivatetimestep
1624!---------------------------------------------------------------------------
1625integer function visitactivatetimestep ()
1626 implicit none
1627 include "visitfortransimV2interface.inc"
1628 visitactivatetimestep = VISIT_OKAY
1629end function visitactivatetimestep
1630
1631!---------------------------------------------------------------------------
1632! visitgetcurve
1633!---------------------------------------------------------------------------
1634integer function visitgetcurve (name, lname)
1635 implicit none
1636 character * 8 name
1637 integer lname
1638 include "visitfortransimV2interface.inc"
1639 visitgetcurve = VISIT_INVALID_HANDLE
1640end function visitgetcurve
1641
1642!---------------------------------------------------------------------------
1643! visitgetdomainlist
1644!---------------------------------------------------------------------------
1645integer function visitgetdomainlist (name, lname)
1646 use psopen_visit
1647 implicit none
1648 include "visitfortransimV2interface.inc"
1649
1650 ! function vars
1651 integer :: lname
1652 character (len=1) :: name (lname)
1653
1654 ! other vars
1655 integer :: h, dl, ierr
1656
1657 WRITE_RANK0_ARGS2(*, 'VisIt: get domain list ',name, 4)
1658
1659 visitgetdomainlist = VISIT_INVALID_HANDLE
1660 h = VISIT_INVALID_HANDLE
1661 dl = VISIT_INVALID_HANDLE
1662
1663 ! Tell VisIt that there are as many domains as processors and this
1664 ! processor just has one of them.
1665 h = VISIT_INVALID_HANDLE
1666 if (visitdomainlistalloc(h) == VISIT_OKAY) then
1667 if (visitvardataalloc(dl) == VISIT_OKAY) then
1668 ierr = visitvardataseti(dl, VISIT_OWNER_SIM, 1, 1, vis_mpi_proc_id); CHKERRQ1(ierr,VISIT_OKAY)
1669 ierr = visitdomainlistsetdomains(h, vis_mpi_nproc, dl); CHKERRQ1(ierr,VISIT_OKAY)
1670 end if
1671 end if
1672 visitgetdomainlist = h
1673
1674end function visitgetdomainlist
1675
1676!---------------------------------------------------------------------------
1677! visitgetdomainbounds
1678!---------------------------------------------------------------------------
1679integer function visitgetdomainbounds (name, lname)
1680 implicit none
1681 character * 8 name
1682 integer lname
1683 include "visitfortransimV2interface.inc"
1684 visitgetdomainbounds = VISIT_INVALID_HANDLE
1685end function visitgetdomainbounds
1686
1687!---------------------------------------------------------------------------
1688! visitgetdomainnesting
1689!---------------------------------------------------------------------------
1690integer function visitgetdomainnesting (name, lname)
1691 implicit none
1692 character * 8 name
1693 integer lname
1694 include "visitfortransimV2interface.inc"
1695 visitgetdomainnesting = VISIT_INVALID_HANDLE
1696end function visitgetdomainnesting
1697
1698!---------------------------------------------------------------------------
1699! visitgetmaterial
1700!---------------------------------------------------------------------------
1701integer function visitgetmaterial (domain, name, lname)
1702 implicit none
1703 character * 8 name
1704 integer domain, lname
1705 include "visitfortransimV2interface.inc"
1706 visitgetmaterial = VISIT_INVALID_HANDLE
1707end function visitgetmaterial