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