KTH framework for Nek5000 toolboxes; testing version  0.0.1
visit.f
Go to the documentation of this file.
1 c---------------------------------------------------------------------
2 c VisIt Simulation Code for Nek5000
3 c
4 c Provide a connection to VisIt simulation code. You will be able
5 c to connect to VisIt and drive the simulation though VisIt or
6 c have sim code send data to VisIt.
7 c---------------------------------------------------------------------
8 
9 c---------------------------------------------------------------------
10 c visit_init
11 c
12 c Initialize VisIt. This will create the sim file needed by Visit
13 c and create the initial VisIt setup.
14 c---------------------------------------------------------------------
15  subroutine visit_init
16  implicit none
17  include "visitfortransimV2interface.inc"
18  include "mpif.h"
19 c // local variables
20  integer err
21 c // SIMSTATE common block
22  integer runflag, endflag
23  common /simstate/ runflag, endflag
24  save /simstate/
25 c // PARALLEL state common block
26  integer par_rank, par_size
27  common /parallel/ par_rank, par_size
28  save /parallel/
29 
30 #ifdef VISIT_STOP
31 c // The sim will wait for VisIt to connect after first step.
32  runflag = 0
33 #else
34 c // Default is to run sim and VisIt can connect any time.
35  runflag = 1
36 #endif
37  endflag = 0
38 
39 c // Determine the rank and size of this MPI task so we can tell
40 c // VisIt's libsim about it.
41  call mpi_comm_rank(mpi_comm_world, par_rank, err)
42  call mpi_comm_size(mpi_comm_world, par_size, err)
43  if(par_size.gt.1) then
44  err = visitsetparallel(1)
45  endif
46  err = visitsetparallelrank(par_rank)
47 
48  call simulationarguments()
49 c // TODO: look for the visitsetupenv2 function.
50 c // Has better scaling, but has not been release for fortran.
51  err = visitsetupenv()
52 c // Have the master process create the sim file.
53  if(par_rank.eq.0) then
54  err = visitinitializesim("nek5000", 7,
55  . "Nek5000 Simulation", 18,
56  . "/no/useful/path", 15,
57  . visit_f77nullstring, visit_f77nullstringlen,
58  . visit_f77nullstring, visit_f77nullstringlen,
59  . visit_f77nullstring, visit_f77nullstringlen)
60  endif
61 
62  end
63 
64 c---------------------------------------------------------------------
65 c simulationarguments
66 c This routine handles command line arguments
67 c -dir <VisIt directory>
68 c -options <VisIt Options>
69 c -trace <VisIt trace file>
70 c---------------------------------------------------------------------
71  subroutine simulationarguments()
72  implicit none
73  character (len=80) str
74  integer err, i, N, len
75  integer visitsetoptions, visitsetdirectory, visitopentracefile
76 
77  n = iargc()
78  i = 1
79  len = 80
80 5 if (i.le.n) then
81  call getarg(i, str)
82  if(str.eq."-dir") then
83  call getarg(i+1, str)
84  err = visitsetdirectory(str, len)
85  i = i + 1
86  elseif(str.eq."-options") then
87  call getarg(i+1, str)
88  err = visitsetoptions(str, len)
89  i = i + 1
90  elseif(str.eq."-trace") then
91  call getarg(i+1, str)
92  err = visitopentracefile(str, len)
93  i = i + 1
94  endif
95  i = i + 1
96  goto 5
97  endif
98  end
99 
100 c---------------------------------------------------------------------
101 c processvisitcommand
102 c---------------------------------------------------------------------
103  integer function processvisitcommand()
104  implicit none
105  include "mpif.h"
106  include "visitfortransimV2interface.inc"
107 c // PARALLEL state common block
108  integer par_rank, par_size
109  common /parallel/ par_rank, par_size
110  integer command, e, doloop, success, ret
111  integer visit_command_process
112  integer visit_command_success
113  integer visit_command_failure
114  parameter(visit_command_process = 0)
115  parameter(visit_command_success = 1)
116  parameter(visit_command_failure = 2)
117 
118  if(par_rank.eq.0) then
119  success = visitprocessenginecommand()
120 
121  if(success.gt.0) then
122  command = visit_command_success
123  ret = 1
124  else
125  command = visit_command_failure
126  ret = 0
127  endif
128 
129  call mpi_bcast(command,1,mpi_integer,0,mpi_comm_world,e)
130  else
131  doloop = 1
132 2345 call mpi_bcast(command,1,mpi_integer,0,mpi_comm_world,e)
133  if(command.eq.visit_command_process) then
134  success = visitprocessenginecommand()
135  elseif(command.eq.visit_command_success) then
136  ret = 1
137  doloop = 0
138  else
139  ret = 0
140  doloop = 0
141  endif
142 
143  if(doloop.ne.0) then
144  goto 2345
145  endif
146  endif
147  processvisitcommand = ret
148  end
149 
150 c---------------------------------------------------------------------
151 c visit_check
152 c---------------------------------------------------------------------
153  subroutine visit_check()
154  implicit none
155  include "mpif.h"
156  include "visitfortransimV2interface.inc"
157 c // functions
158  integer processvisitcommand
159 c // local variables
160  integer visitstate, result, blocking, ierr
161 c // SIMSTATE common block
162  integer runflag, endflag
163  common /simstate/ runflag, endflag
164 c // PARALLEL state common block
165  integer par_rank, par_size
166  common /parallel/ par_rank, par_size
167 
168 c // If at the end of sim, lets not force an update.
169  if(endflag.eq.0) then
170 c // Check if we are connected to VisIt
171  result = visitisconnected()
172  if(result.eq.1) then
173 c // Tell VisIt that the timestep changed
174  result = visittimestepchanged()
175 c // Tell VisIt to update its plots
176  result = visitupdateplots()
177  endif
178  endif
179 
180 c write (6, 2000)
181  2000 format('VisIt Check!')
182  do 10
183 c // If we are running don't block
184  if(runflag.eq.1) then
185  blocking = 0
186  else
187  blocking = 1
188  endif
189 
190 c // Detect input from VisIt on processor 0 and then broadcast
191 c // the results of that input to all processors.
192  if(par_rank.eq.0) then
193  visitstate = visitdetectinput(blocking, -1)
194  endif
195  call mpi_bcast(visitstate,1,mpi_integer,0,mpi_comm_world,ierr)
196 
197  if (visitstate.eq.0) then
198 c // Okay - Process time step.
199  goto 1234
200  elseif (visitstate.eq.1) then
201 c // Attempt to Connect VisIt
202  ierr = runflag
203  runflag = 0
204  result = visitattemptconnection()
205  if (result.eq.1) then
206  write (6, 2001)
207  2001 format('VisIt connected!')
208  else
209  write (6, 2002)
210  2002 format('VisIt did not connected!')
211  endif
212  flush( 6 )
213  runflag = ierr
214  elseif (visitstate.eq.2) then
215 c // Engine socket input
216 c ierr = runflag
217 c runflag = 0
218  if (processvisitcommand().eq.0) then
219  result = visitdisconnect()
220 c // If VisIt is disconnect lets run.
221  runflag = 1
222  endif
223 c // Check if user wants to exit sim.
224  if(runflag.eq.2) then
225  goto 1234
226  endif
227  elseif (visitstate.eq.3) then
228 c // Console socket input
229  elseif (visitstate.lt.0) then
230 c // Error
231  goto 1234
232  endif
233 10 continue
234 1234 end
235 
236 c---------------------------------------------------------------------
237 c visit_end
238 c---------------------------------------------------------------------
239  subroutine visit_end()
240  implicit none
241  include "visitfortransimV2interface.inc"
242 c // local variables
243  integer result
244 c // SIMSTATE common block
245  integer runflag, endflag
246  common /simstate/ runflag, endflag
247 
248 c // Check if we are connected to VisIt
249  result = visitisconnected()
250  if(result.eq.1) then
251 c // Let VisIt exit the sim.
252 
253 c // This will tell the visit_check function we are at the end.
254  endflag = 1
255  runflag = 0
256 
257  do 10
258  call visit_check()
259 
260 c // User asked to finish.
261  if (endflag.eq.2) then
262  EXIT
263  endif
264 10 continue
265  endif
266 
267  end
268 
269 cccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc
270 c
271 c These functions must be defined to satisfy the
272 c visitfortransimV2interface lib.
273 c
274 cccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc
275 
276 c---------------------------------------------------------------------
277 c visitcommandcallback
278 c Handle User defined functions from VisIt user interface.
279 c---------------------------------------------------------------------
280  subroutine visitcommandcallback (cmd, lcmd, args, largs)
281  implicit none
282  character*8 cmd, args
283  integer lcmd, largs
284  include "visitfortransimV2interface.inc"
285 c // SIMSTATE common block
286  integer runflag, endflag
287  common /simstate/ runflag, endflag
288 
289 c // Handle the commands that we define in visitgetmetadata.
290  if(visitstrcmp(cmd, lcmd, "stop", 4).eq.0) then
291  runflag = 0
292  elseif(visitstrcmp(cmd, lcmd, "step", 4).eq.0) then
293  runflag = 2
294  elseif(visitstrcmp(cmd, lcmd, "run", 3).eq.0) then
295  runflag = 1
296  elseif(visitstrcmp(cmd, lcmd, "exit", 4).eq.0) then
297  call exitt()
298  elseif(visitstrcmp(cmd, lcmd, "finish", 6).eq.0) then
299  if(endflag.eq.1) then
300  endflag = 2
301  runflag = 2
302  endif
303  endif
304  end
305 
306 c---------------------------------------------------------------------
307 c visitbroadcastintfunction
308 c---------------------------------------------------------------------
309  integer function visitbroadcastintfunction(value, sender)
310  implicit none
311  include "mpif.h"
312  integer value, sender, ierr
313  call mpi_bcast(value,1,mpi_integer,sender,mpi_comm_world,ierr)
315  end
316 
317 c---------------------------------------------------------------------
318 c visitbroadcaststringfunction
319 c---------------------------------------------------------------------
320  integer function visitbroadcaststringfunction(str, lstr, sender)
321  implicit none
322  include "mpif.h"
323  character*8 str
324  integer lstr, sender, ierr
325  call mpi_bcast(str,lstr,mpi_character,sender,mpi_comm_world,ierr)
327  end
328 
329 c---------------------------------------------------------------------
330 c visitslaveprocesscallback
331 c---------------------------------------------------------------------
333  implicit none
334  include "mpif.h"
335  integer c, ierr, VISIT_COMMAND_PROCESS
336  parameter(visit_command_process = 0)
337  c = visit_command_process
338  call mpi_bcast(c, 1, mpi_integer, 0, mpi_comm_world, ierr)
339  end
340 
341 c---------------------------------------------------------------------
342 c visitactivatetimestep
343 c---------------------------------------------------------------------
344  integer function visitactivatetimestep()
345  implicit none
346  include "visitfortransimV2interface.inc"
347  visitactivatetimestep = visit_okay
348  end
349 
350 c---------------------------------------------------------------------
351 c visitgetmetadata
352 c This function tells VisIt about all of the data that this
353 c simulation will output. Mesh, variables, expressions and commands.
354 c---------------------------------------------------------------------
355  integer function visitgetmetadata()
356  include 'SIZE'
357  include 'TSTEP'
358  include 'INPUT'
359  include 'PARALLEL'
360  include "visitfortransimV2interface.inc"
361 c // SIMSTATE common block
362  integer runflag, endflag
363  common /simstate/ runflag, endflag
364 c // PARALLEL state common block
365  integer par_rank, par_size
366  common /parallel/ par_rank, par_size
367 c // Local variables
368  integer md, mmd, vmd, cmd, emd, err, idim, k
369  character*3 seqstring
370 
371  if(visitmdsimalloc(md).eq.visit_okay) then
372  err = visitmdsimsetcycletime(md, istep, time)
373  if(runflag.eq.1) then
374  err = visitmdsimsetmode(md, visit_simmode_running)
375  else
376  err = visitmdsimsetmode(md, visit_simmode_stopped)
377  endif
378 
379 c // TODO: ask why this changes after first save?
380 c write(6,*) 'VISIT: IFXYO',IFXYO
381 c if(IFXYO) then
382 c // Add a mesh
383  if(if3d) then
384  idim = 3
385  else
386  idim = 2
387  endif
388 
389  if(visitmdmeshalloc(mmd).eq.visit_okay) then
390  err = visitmdmeshsetname(mmd, "mesh", 4)
391  err = visitmdmeshsetmeshtype(mmd,
392  . visit_meshtype_curvilinear)
393  err = visitmdmeshsettopologicaldim(mmd, idim)
394  err = visitmdmeshsetspatialdim(mmd, idim)
395  err = visitmdmeshsetnumdomains(mmd, nelgt)
396  err = visitmdmeshsetdomaintitle(mmd, "Domains", 7)
397  err = visitmdmeshsetdomainpiecename(mmd, "domain", 6)
398 c err = visitmdmeshsetxunits(mmd, "cm", 2)
399 c err = visitmdmeshsetyunits(mmd, "cm", 2)
400  err = visitmdmeshsetxlabel(mmd, "X-Axis", 6)
401  err = visitmdmeshsetylabel(mmd, "Y-Axis", 6)
402  if(if3d) then
403 c err = visitmdmeshsetzunits(mmd, "cm", 2)
404  err = visitmdmeshsetzlabel(mmd, "Z-Axis", 6)
405  endif
406  err = visitmdsimaddmesh(md, mmd)
407  endif
408 c endif
409 
410  if(ifvo) then
411 c // Add a X velocity variable on the mesh.
412  if(visitmdvaralloc(vmd).eq.visit_okay) then
413  err = visitmdvarsetname(vmd, "x_velocity", 10)
414  err = visitmdvarsetmeshname(vmd, "mesh", 4)
415 c err = visitmdvarsetunits(vmd, "cm", 2)
416  err = visitmdvarsetcentering(vmd, visit_varcentering_node)
417  err = visitmdvarsettype(vmd, visit_vartype_scalar)
418 c err = visitmdvarsettreatasascii(vmd, 0)
419  err = visitmdsimaddvariable(md, vmd)
420  endif
421 
422 c // Add a Y velocity variable on the mesh.
423  if(visitmdvaralloc(vmd).eq.visit_okay) then
424  err = visitmdvarsetname(vmd, "y_velocity", 10)
425  err = visitmdvarsetmeshname(vmd, "mesh", 4)
426 c err = visitmdvarsetunits(vmd, "cm", 2)
427  err = visitmdvarsetcentering(vmd, visit_varcentering_node)
428  err = visitmdvarsettype(vmd, visit_vartype_scalar)
429 c err = visitmdvarsettreatasascii(vmd, 0)
430  err = visitmdsimaddvariable(md, vmd)
431  endif
432 
433  if(if3d) then
434 c // Add a Z velocity variable on the mesh.
435  if(visitmdvaralloc(vmd).eq.visit_okay) then
436  err = visitmdvarsetname(vmd, "z_velocity", 10)
437  err = visitmdvarsetmeshname(vmd, "mesh", 4)
438 c err = visitmdvarsetunits(vmd, "cm", 2)
439  err = visitmdvarsetcentering(vmd, visit_varcentering_node)
440  err = visitmdvarsettype(vmd, visit_vartype_scalar)
441 c err = visitmdvarsettreatasascii(vmd, 0)
442  err = visitmdsimaddvariable(md, vmd)
443  endif
444  endif
445 
446 c // Add a velocity expression.
447  if(visitmdexpralloc(emd).eq.visit_okay) then
448  err = visitmdexprsetname(emd, "velocity", 8)
449  if(if3d) then
450  err = visitmdexprsetdefinition(emd,
451  . "{x_velocity, y_velocity, z_velocity}", 36)
452  else
453  err = visitmdexprsetdefinition(emd,
454  . "{x_velocity, y_velocity}", 24)
455  endif
456  err = visitmdexprsettype(emd, visit_vartype_vector)
457 
458  err = visitmdsimaddexpression(md, emd)
459  endif
460 
461 c // Add a velocity magnitude expression.
462  if(visitmdexpralloc(emd).eq.visit_okay) then
463  err = visitmdexprsetname(emd, "velocity_mag", 12)
464  err = visitmdexprsetdefinition(emd,
465  . "magnitude(velocity)", 19)
466  err = visitmdexprsettype(emd, visit_vartype_scalar)
467 
468  err = visitmdsimaddexpression(md, emd)
469  endif
470  endif
471 
472  if(ifpo) then
473 c // Add a pressure variable on the mesh.
474  if(visitmdvaralloc(vmd).eq.visit_okay) then
475  err = visitmdvarsetname(vmd, "pressure", 8)
476  err = visitmdvarsetmeshname(vmd, "mesh", 4)
477 c err = visitmdvarsetunits(vmd, "cm", 2)
478  err = visitmdvarsetcentering(vmd, visit_varcentering_node)
479  err = visitmdvarsettype(vmd, visit_vartype_scalar)
480 c err = visitmdvarsettreatasascii(vmd, 0)
481  err = visitmdsimaddvariable(md, vmd)
482  endif
483  endif
484 
485  if(ifto) then
486 c // Add a temperature variable on the mesh.
487  if(visitmdvaralloc(vmd).eq.visit_okay) then
488  err = visitmdvarsetname(vmd, "temperature", 11)
489  err = visitmdvarsetmeshname(vmd, "mesh", 4)
490 c err = visitmdvarsetunits(vmd, "cm", 2)
491  err = visitmdvarsetcentering(vmd, visit_varcentering_node)
492  err = visitmdvarsettype(vmd, visit_vartype_scalar)
493 c err = visitmdvarsettreatasascii(vmd, 0)
494  err = visitmdsimaddvariable(md, vmd)
495  endif
496  endif
497 
498  do k=1,ldimt-1
499 c // Add a user defined variable on the mesh.
500  if(ifpsco(k)) then
501  if(visitmdvaralloc(vmd).eq.visit_okay) then
502  write (seqstring,'(I0)') k
503  err = visitmdvarsetname(vmd, "s"//trim(seqstring),
504  . 1+len_trim(seqstring))
505  err = visitmdvarsetmeshname(vmd, "mesh", 4)
506 c err = visitmdvarsetunits(vmd, "cm", 2)
507  err = visitmdvarsetcentering(vmd, visit_varcentering_node)
508  err = visitmdvarsettype(vmd, visit_vartype_scalar)
509 c err = visitmdvarsettreatasascii(vmd, 0)
510  err = visitmdsimaddvariable(md, vmd)
511  endif
512  endif
513  enddo
514 
515 
516 c // Add simulation commands
517  err = visitmdcmdalloc(cmd)
518  if(err.eq.visit_okay) then
519  err = visitmdcmdsetname(cmd, "stop", 4)
520  err = visitmdsimaddgenericcommand(md, cmd)
521  endif
522  err = visitmdcmdalloc(cmd)
523  if(err.eq.visit_okay) then
524  err = visitmdcmdsetname(cmd, "step", 4)
525  err = visitmdsimaddgenericcommand(md, cmd)
526  endif
527  err = visitmdcmdalloc(cmd)
528  if(err.eq.visit_okay) then
529  err = visitmdcmdsetname(cmd, "run", 3)
530  err = visitmdsimaddgenericcommand(md, cmd)
531  endif
532  err = visitmdcmdalloc(cmd)
533  if(err.eq.visit_okay) then
534  err = visitmdcmdsetname(cmd, "exit", 4)
535  err = visitmdsimaddgenericcommand(md, cmd)
536  endif
537  err = visitmdcmdalloc(cmd)
538  if(err.eq.visit_okay) then
539  err = visitmdcmdsetname(cmd, "finish", 6)
540  err = visitmdsimaddgenericcommand(md, cmd)
541  endif
542  endif
543  visitgetmetadata = md
544  end
545 
546 c---------------------------------------------------------------------
547 c visitgetmesh
548 c Use this function to return mesh data to VisIt.
549 c---------------------------------------------------------------------
550  integer function visitgetmesh(domain, name, lname)
551  include 'SIZE'
552  include 'TOTAL'
553  character*8 name
554  integer domain, lname
555  include "visitfortransimV2interface.inc"
556 
557 c // local variables
558  integer vh, x, y, z, err, dl, dmdims(3)
559 
560  vh = visit_invalid_handle
561 
562 c // Fortran starting index 1, but VisIt is 0
563 c // Also need to get local domain index
564  domain = gllel(domain + 1)
565 
566  if(visitstrcmp(name, lname, "mesh", 4).eq.0) then
567  if(visitcurvmeshalloc(vh).eq.visit_okay) then
568  err = visitvardataalloc(x)
569  err = visitvardataalloc(y)
570  if(if3d) err = visitvardataalloc(z)
571 
572  dl = lx1 * ly1 * lz1
573  err = visitvardatasetd(x, visit_owner_sim, 1, dl,
574  . xm1(1,1,1,domain))
575  err = visitvardatasetd(y, visit_owner_sim, 1, dl,
576  . ym1(1,1,1,domain))
577  if(if3d) then
578  err = visitvardatasetd(z, visit_owner_sim, 1, dl,
579  . zm1(1,1,1,domain))
580  endif
581 
582  dmdims(1) = lx1
583  dmdims(2) = ly1
584  dmdims(3) = lz1
585  if(if3d) then
586  err = visitcurvmeshsetcoordsxyz(vh, dmdims, x, y, z)
587  else
588  err = visitcurvmeshsetcoordsxy(vh, dmdims, x, y)
589  endif
590  endif
591  endif
592 
593  visitgetmesh = vh
594  end
595 
596 c---------------------------------------------------------------------
597 c visitgetmaterial
598 c Use this function to return material data to VisIt.
599 c---------------------------------------------------------------------
600  integer function visitgetmaterial(domain, name, lname)
601  implicit none
602  character*8 name
603  integer domain, lname
604  include "visitfortransimV2interface.inc"
605  visitgetmaterial = visit_invalid_handle
606  end
607 
608 c---------------------------------------------------------------------
609 c visitgetvariable
610 c Use this function to return variable data to VisIt.
611 c---------------------------------------------------------------------
612  integer function visitgetvariable(gdomain, name, lname)
613 c implicit none
614  include 'SIZE'
615  include 'TOTAL'
616  character*8 name
617  integer gdomain, lname
618  include "visitfortransimV2interface.inc"
619 c // local vars
620  integer h, nvals, err, domain, k
621 
622  nvals = lx1 * ly1 * lz1
623 
624  h = visit_invalid_handle
625 
626 c // Fortran starting index 1, but VisIt is 0
627 c // Also need to get local domain index
628  gdomain = gdomain + 1
629  domain = gllel(gdomain)
630 
631  if(visitvardataalloc(h).eq.visit_okay) then
632  if(visitstrcmp(name, lname, "temperature", 11).eq.0) then
633  err = visitvardatasetd(h, visit_owner_sim, 1, nvals,
634  . t(1,1,1,domain,1))
635  elseif(visitstrcmp(name, lname, "pressure", 8).eq.0) then
636  if(gdomain.le.nelgv) then
637  err = visitvardatasetd(h, visit_owner_sim, 1, nvals,
638  . pr(1,1,1,domain))
639  endif
640  elseif(visitstrcmp(name, lname, "x_velocity", 10).eq.0) then
641  if(gdomain.le.nelgv) then
642  err = visitvardatasetd(h, visit_owner_sim, 1, nvals,
643  . vx(1,1,1,domain))
644  endif
645  elseif(visitstrcmp(name, lname, "y_velocity", 10).eq.0) then
646  if(gdomain.le.nelgv) then
647  err = visitvardatasetd(h, visit_owner_sim, 1, nvals,
648  . vy(1,1,1,domain))
649  endif
650  elseif(visitstrcmp(name, lname, "z_velocity", 10).eq.0) then
651  if(gdomain.le.nelgv) then
652  err = visitvardatasetd(h, visit_owner_sim, 1, nvals,
653  . vz(1,1,1,domain))
654  endif
655  elseif(visitstrcmp("s", 1, name, 1).eq.0) then
656 c // Handle the user define variables.
657  read( name(2:lname), '(i10)' ) k
658  if(ifpsco(k)) then
659  k = k + 1
660  err = visitvardatasetd(h, visit_owner_sim, 1, nvals,
661  . t(1,1,1,domain,k))
662  endif
663  endif
664  endif
665 
666  visitgetvariable = h
667  end
668 
669 c---------------------------------------------------------------------
670 c visitgetcurve
671 c Use this function to return curve data to VisIt.
672 c---------------------------------------------------------------------
673  integer function visitgetcurve(handle, name, lname)
674  implicit none
675  character*8 name
676  integer handle, lname
677  include "visitfortransimV2interface.inc"
678  visitgetcurve = visit_invalid_handle
679  end
680 
681 c---------------------------------------------------------------------
682 c visitgetdomainlist
683 c This function returns a list of domains owned by this process.
684 c---------------------------------------------------------------------
685  integer function visitgetdomainlist()
686  include 'SIZE'
687  include 'PARALLEL'
688  include "visitfortransimV2interface.inc"
689 c // local vars
690  integer h, dl, err
691 
692  h = visit_invalid_handle
693 
694  if(visitdomainlistalloc(h).eq.visit_okay) then
695  if(visitvardataalloc(dl).eq.visit_okay) then
696 c // Hack to work around the index difference between
697 c // fortran and VisIt C. Temp change and make copy.
698  do i=1,nelt
699  lglel(i) = lglel(i) - 1
700  enddo
701  err = visitvardataseti(dl, visit_owner_copy,1,nelt,lglel)
702  err = visitdomainlistsetdomains(h, nelgt, dl)
703 c // restore correct fortran values.
704  do i=1,nelt
705  lglel(i) = lglel(i) + 1
706  enddo
707  endif
708  endif
709 
711  end
712 
713 c---------------------------------------------------------------------
714 c visitgetdomainbounds
715 c This function allows VisIt to create ghost zones between domains.
716 c---------------------------------------------------------------------
717  integer function visitgetdomainbounds(name, lname)
718  implicit none
719  character*8 name
720  integer lname
721  include "visitfortransimV2interface.inc"
722  visitgetdomainbounds = visit_invalid_handle
723  end
724 
725 c---------------------------------------------------------------------
726 c visitgetdomainnesting
727 c This is used to tell VisIt how AMR patches are nested.
728 c---------------------------------------------------------------------
729  integer function visitgetdomainnesting(name, lname)
730  implicit none
731  character*8 name
732  integer lname
733  include "visitfortransimV2interface.inc"
734  visitgetdomainnesting = visit_invalid_handle
735  end
736 
737 c---------------------------------------------------------------------
738 c visitgetmixedvariable
739 c This is needed to run with newer version of Visit (2.12.0 tested).
740 c TODO: The comment should be changed
741 c---------------------------------------------------------------------
742  integer function visitgetmixedvariable(domain, name, lname)
743  implicit none
744  character*8 name
745  integer domain, lname
746  include "visitfortransimV2interface.inc"
747  visitgetmixedvariable = visit_invalid_handle
748  end
void exitt()
Definition: comm_mpi.f:604
integer function gllel(ieg)
Definition: dprocmap.f:183
subroutine mpi_comm_rank(comm, me, ierror)
Definition: mpi_dummy.f:389
subroutine mpi_comm_size(comm, nprocs, ierror)
Definition: mpi_dummy.f:410
subroutine mpi_bcast(data, n, datatype, node, comm, ierror)
Definition: mpi_dummy.f:213
integer function visitgetvariable(gdomain, name, lname)
Definition: visit.f:613
subroutine visit_check()
Definition: visit.f:154
integer function visitactivatetimestep()
Definition: visit.f:345
integer function visitgetcurve(handle, name, lname)
Definition: visit.f:674
integer function visitgetmaterial(domain, name, lname)
Definition: visit.f:601
subroutine simulationarguments()
Definition: visit.f:72
integer function processvisitcommand()
Definition: visit.f:104
integer function visitgetdomainbounds(name, lname)
Definition: visit.f:718
integer function visitgetmesh(domain, name, lname)
Definition: visit.f:551
integer function visitgetmetadata()
Definition: visit.f:356
integer function visitgetdomainlist()
Definition: visit.f:686
integer function visitgetdomainnesting(name, lname)
Definition: visit.f:730
integer function visitgetmixedvariable(domain, name, lname)
Definition: visit.f:743
subroutine visitcommandcallback(cmd, lcmd, args, largs)
Definition: visit.f:281
subroutine visit_end()
Definition: visit.f:240
integer function visitbroadcaststringfunction(str, lstr, sender)
Definition: visit.f:321
subroutine visitslaveprocesscallback()
Definition: visit.f:333
subroutine visit_init
Definition: visit.f:16
integer function visitbroadcastintfunction(value, sender)
Definition: visit.f:310