All Classes Files Functions Variables Pages
module_snapshot_volume.f90
Go to the documentation of this file.
1 !=====================================================================================================================================
2 ! EFISPEC3D !
3 ! (Elements FInis SPECtraux 3D) !
4 ! !
5 ! http://efispec.free.fr !
6 ! !
7 ! !
8 ! This file is part of EFISPEC3D !
9 ! Please refer to http://efispec.free.fr if you use it or part of it !
10 ! !
11 ! !
12 !1 ---> French License: CeCILL V2 !
13 ! !
14 ! Copyright BRGM 2009 contributeurs : Florent DE MARTIN !
15 ! David MICHEA !
16 ! Philippe THIERRY !
17 ! !
18 ! Contact: f.demartin at brgm.fr !
19 ! !
20 ! Ce logiciel est un programme informatique servant a resoudre l'equation du !
21 ! mouvement en trois dimensions via une methode des elements finis spectraux. !
22 ! !
23 ! Ce logiciel est regi par la licence CeCILL soumise au droit francais et !
24 ! respectant les principes de diffusion des logiciels libres. Vous pouvez !
25 ! utiliser, modifier et/ou redistribuer ce programme sous les conditions de la !
26 ! licence CeCILL telle que diffusee par le CEA, le CNRS et l'INRIA sur le site !
27 ! "http://www.cecill.info". !
28 ! !
29 ! En contrepartie de l'accessibilite au code source et des droits de copie, de !
30 ! modification et de redistribution accordes par cette licence, il n'est offert !
31 ! aux utilisateurs qu'une garantie limitee. Pour les memes raisons, seule une !
32 ! responsabilite restreinte pese sur l'auteur du programme, le titulaire des !
33 ! droits patrimoniaux et les concedants successifs. !
34 ! !
35 ! A cet egard l'attention de l'utilisateur est attiree sur les risques associes !
36 ! au chargement, a l'utilisation, a la modification et/ou au developpement et a !
37 ! la reproduction du logiciel par l'utilisateur etant donne sa specificite de !
38 ! logiciel libre, qui peut le rendre complexe a manipuler et qui le reserve donc !
39 ! a des developpeurs et des professionnels avertis possedant des connaissances !
40 ! informatiques approfondies. Les utilisateurs sont donc invites a charger et !
41 ! tester l'adequation du logiciel a leurs besoins dans des conditions permettant !
42 ! d'assurer la securite de leurs systemes et ou de leurs donnees et, plus !
43 ! generalement, a l'utiliser et l'exploiter dans les memes conditions de !
44 ! securite. !
45 ! !
46 ! Le fait que vous puissiez acceder a cet en-tete signifie que vous avez pris !
47 ! connaissance de la licence CeCILL et que vous en avez accepte les termes. !
48 ! !
49 ! !
50 !2 ---> International license: GNU GPL V3 !
51 ! !
52 ! EFISPEC3D is a computer program that solves the three-dimensional equations of !
53 ! motion using a finite spectral-element method. !
54 ! !
55 ! Copyright (C) 2009 Florent DE MARTIN !
56 ! !
57 ! Contact: f.demartin at brgm.fr !
58 ! !
59 ! This program is free software: you can redistribute it and/or modify it under !
60 ! the terms of the GNU General Public License as published by the Free Software !
61 ! Foundation, either version 3 of the License, or (at your option) any later !
62 ! version. !
63 ! !
64 ! This program is distributed in the hope that it will be useful, but WITHOUT ANY !
65 ! WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR A !
66 ! PARTICULAR PURPOSE. See the GNU General Public License for more details. !
67 ! !
68 ! You should have received a copy of the GNU General Public License along with !
69 ! this program. If not, see http://www.gnu.org/licenses/. !
70 ! !
71 ! !
72 !3 ---> Third party libraries !
73 ! !
74 ! EFISPEC3D uses the following of third party libraries: !
75 ! !
76 ! --> METIS 5.1.0 !
77 ! see http://glaros.dtc.umn.edu/gkhome/metis/metis/overview !
78 ! !
79 ! --> Lib_VTK_IO !
80 ! see S. Zaghi's website: https://github.com/szaghi/Lib_VTK_IO !
81 ! !
82 ! --> INTERP_LINEAR !
83 ! see J. Burkardt website: http://people.sc.fsu.edu/~jburkardt/ !
84 ! !
85 ! --> SAC !
86 ! http://ds.iris.edu/files/sac-manual/ !
87 ! !
88 ! --> EXODUS II !
89 ! http://sourceforge.net/projects/exodusii/ !
90 ! !
91 ! --> NETCDF !
92 ! http://www.unidata.ucar.edu/software/netcdf/ !
93 ! !
94 ! --> HDF5 !
95 ! http://www.hdfgroup.org/HDF5/ !
96 ! !
97 ! Some of these libraries are located in directory lib and pre-compiled !
98 ! with intel compiler for x86_64 platform. !
99 ! !
100 ! !
101 !4 ---> Related Articles !
102 ! !
103 ! De Martin, F., Matsushima, M., Kawase, H. (BSSA, 2013) !
104 ! Impact of geometric effects on near-surface Green's functions !
105 ! doi:10.1785/0120130039 !
106 ! !
107 ! Aochi, H., Ducellier, A., Dupros, F., Delatre, M., Ulrich, T., De Martin, F., !
108 ! Yoshimi, M., (Pure Appl. Geophys. 2013) !
109 ! Finite Difference Simulations of Seismic Wave Propagation for the 2007 Mw 6.6 !
110 ! Niigata-ken Chuetsu-Oki Earthquake: Validity of Models and Reliable !
111 ! Input Ground Motion in the Near-Field !
112 ! doi:10.1007/s00024-011-0429-5 !
113 ! !
114 ! De Martin, F. (BSSA, 2011) !
115 ! Verification of a Spectral-Element Method Code for the Southern California !
116 ! Earthquake Center LOH.3 Viscoelastic Case !
117 ! doi:10.1785/0120100305 !
118 ! !
119 !=====================================================================================================================================
120 
124 
128 
129  implicit none
130 
131  private
132 
133  public :: init_snapshot_volume
134  private :: init_vtk_mesh
135  private :: get_efi_hexa
136  private :: init_vtk_node
137  private :: init_vtk_node_numbering
138  private :: init_vtk_hexa_connectivity
139  public :: write_snapshot_volume
140  public :: write_snapshot_volume_vtk
141  public :: write_collection_vtk_vol
142  private :: write_medium
143  private :: collection_vtk_volume
144 
145  contains
146 
147 !
148 !
151 !***********************************************************************************************************************************************************************************
152  subroutine init_snapshot_volume()
153 !***********************************************************************************************************************************************************************************
154 
155  use mod_global_variables, only:&
156  cg_prefix&
157  ,cg_myrank&
158  ,ig_hexa_gll_glonum&
159  ,ig_nhexa&
160  ,ig_hexa_snapshot&
161  ,rg_snapshot_volume_xmin&
162  ,rg_snapshot_volume_xmax&
163  ,rg_snapshot_volume_ymin&
164  ,rg_snapshot_volume_ymax&
165  ,rg_snapshot_volume_zmin&
166  ,rg_snapshot_volume_zmax&
167  ,ig_vtk_hexa_conn_snapshot&
168  ,ig_vtk_node_gll_glonum_snapshot&
169  ,rg_vtk_node_x_snapshot&
170  ,rg_vtk_node_y_snapshot&
171  ,rg_vtk_node_z_snapshot&
172  ,ig_vtk_cell_type_snapshot&
173  ,ig_vtk_offset_snapshot&
174  ,ig_vtk_nhexa_snapshot&
175  ,ig_vtk_nnode_snapshot
176 
177  use mod_init_memory
178 
179  use mod_vtk_io
180 
181  implicit none
182 
183  call get_efi_hexa(ig_hexa_snapshot &
184  ,ig_nhexa &
185  ,rg_snapshot_volume_xmin &
186  ,rg_snapshot_volume_xmax &
187  ,rg_snapshot_volume_ymin &
188  ,rg_snapshot_volume_ymax &
189  ,rg_snapshot_volume_zmin &
190  ,rg_snapshot_volume_zmax )
191 
192  call init_vtk_mesh(ig_hexa_snapshot &
193  ,ig_hexa_gll_glonum &
194  ,ig_vtk_hexa_conn_snapshot &
195  ,rg_vtk_node_x_snapshot &
196  ,rg_vtk_node_y_snapshot &
197  ,rg_vtk_node_z_snapshot &
198  ,ig_vtk_node_gll_glonum_snapshot &
199  ,ig_vtk_cell_type_snapshot &
200  ,ig_vtk_offset_snapshot &
201  ,ig_vtk_nhexa_snapshot &
202  ,ig_vtk_nnode_snapshot )
203 
204  deallocate(ig_hexa_snapshot)
205 
206 ! call write_medium()
207 
208  return
209 
210 !***********************************************************************************************************************************************************************************
211  end subroutine init_snapshot_volume
212 !***********************************************************************************************************************************************************************************
213 
214 !
215 !
227 !***********************************************************************************************************************************************************************************
228  subroutine get_efi_hexa(efi_hexa,efi_nhexa,xmin,xmax,ymin,ymax,zmin,zmax)
229 !***********************************************************************************************************************************************************************************
230 
232 
233  use mod_init_memory
234 
235  implicit none
236 
237  integer, intent(out), dimension(:), allocatable :: efi_hexa
238  integer, intent( in) :: efi_nhexa
239  real , intent( in) :: xmin
240  real , intent( in) :: xmax
241  real , intent( in) :: ymin
242  real , intent( in) :: ymax
243  real , intent( in) :: zmin
244  real , intent( in) :: zmax
245 
246  real :: x
247  real :: y
248  real :: z
249 
250  integer :: ios
251  integer :: ihexa
252  integer :: nhexa
253 
254  logical , dimension(efi_nhexa) :: is_snapshot
255 
256 !
257 !---->first pass to get the number of EFISPEC's hexahedron elements located inside the volume used for snapshot. The centroid of hexahedron elements is used as reference for their location.
258 
259  nhexa = 0
260  is_snapshot(:) = .false.
261 
262  do ihexa = 1,efi_nhexa
263 
264  call compute_hexa_point_coord(ihexa,0.0,0.0,0.0,x,y,z)
265 
266  if ( x >= xmin .and. x <= xmax .and.&
267  y >= ymin .and. y <= ymax .and.&
268  z >= zmin .and. z <= zmax ) then
269 
270  is_snapshot(ihexa) = .true.
271  nhexa = nhexa + 1
272 
273  endif
274 
275  enddo
276 
277 !
278 !---->second pass to store hexahedron elements
279 
280  ios = init_array_int(efi_hexa,nhexa,"efi_hexa")
281 
282  nhexa = 0
283 
284  do ihexa = 1,efi_nhexa
285 
286  if (is_snapshot(ihexa)) then
287 
288  nhexa = nhexa + 1
289 
290  efi_hexa(nhexa) = ihexa
291 
292  endif
293 
294  enddo
295 
296  return
297 
298 !***********************************************************************************************************************************************************************************
299  end subroutine get_efi_hexa
300 !***********************************************************************************************************************************************************************************
301 
302 
303 !
304 !
327 !***********************************************************************************************************************************************************************************
328  subroutine init_vtk_mesh(efi_hexa,efi_hexa_gll_glonum,vtk_hexa_conn,vtk_node_x,vtk_node_y,vtk_node_z,vtk_node_gll_glonum,vtk_cell_type,vtk_offset,vtk_nhexa,vtk_nnode)
329 !***********************************************************************************************************************************************************************************
330 
331  use mod_global_variables, only :&
332  ig_ngll&
333  ,rg_gll_coordinate&
334  ,ig_ngll_total&
335  ,ig_nhexa
336 
338 
339  use mod_init_memory
340 
341  implicit none
342 
343  integer , intent( in), dimension(:) :: efi_hexa
344  integer , intent( in), dimension(IG_NGLL,IG_NGLL,IG_NGLL,ig_nhexa) :: efi_hexa_gll_glonum
345  integer , intent(out), dimension(:), allocatable :: vtk_hexa_conn
346  integer , intent(out), dimension(:), allocatable :: vtk_node_gll_glonum
347  integer(kind=1), intent(out), dimension(:), allocatable :: vtk_cell_type
348  integer , intent(out), dimension(:), allocatable :: vtk_offset
349  real , intent(out), dimension(:), allocatable :: vtk_node_x
350  real , intent(out), dimension(:), allocatable :: vtk_node_y
351  real , intent(out), dimension(:), allocatable :: vtk_node_z
352  integer , intent(out) :: vtk_nhexa
353  integer , intent(out) :: vtk_nnode
354 
355  integer , dimension(ig_ngll_total) :: vtk_node
356  integer :: efi_nhexa
357  integer :: ios
358  integer :: ihexa
359 
360  efi_nhexa = size(efi_hexa)
361 
362  !
363  !
364  !**********************************************************
365  !->find which GLL nodes must be used to split hexa by GLL
366  !**********************************************************
367  call init_vtk_node(efi_hexa,efi_nhexa,efi_hexa_gll_glonum,vtk_node,ig_ngll_total,vtk_nnode)
368 
369  !
370  !
371  !********************************************************************
372  !->renumber array vtk_node
373  !********************************************************************
374  call init_vtk_node_numbering(vtk_node,vtk_node_gll_glonum,vtk_nnode,ig_ngll_total)
375 
376  !
377  !
378  !*************************************************
379  !->create hexa connectivity array 'vtk_hexa_conn'
380  !*************************************************
381  call init_vtk_hexa_connectivity(efi_hexa,efi_nhexa,vtk_node,ig_ngll_total,efi_hexa_gll_glonum,vtk_hexa_conn,vtk_nhexa)
382 
383  !
384  !
385  !********************************************
386  !->create hexa coordinates
387  !********************************************
388  ios = init_array_real(vtk_node_x,vtk_nnode,"vtk_node_x")
389  ios = init_array_real(vtk_node_y,vtk_nnode,"vtk_node_y")
390  ios = init_array_real(vtk_node_z,vtk_nnode,"vtk_node_z")
391 
392  call get_node_gll_value(rg_gll_coordinate(1,:),vtk_node_gll_glonum,vtk_node_x)
393  call get_node_gll_value(rg_gll_coordinate(2,:),vtk_node_gll_glonum,vtk_node_y)
394  call get_node_gll_value(rg_gll_coordinate(3,:),vtk_node_gll_glonum,vtk_node_z)
395 
396  !
397  !
398  !********************************************
399  !->create hexa cell type
400  !********************************************
401  allocate(vtk_cell_type(vtk_nhexa))
402 
403  do ihexa = 1,vtk_nhexa
404  vtk_cell_type(ihexa) = 12
405  enddo
406 
407  !
408  !
409  !********************************************
410  !->create hexa offset
411  !********************************************
412  ios = init_array_int(vtk_offset,vtk_nhexa,"vtk_offset")
413 
414  do ihexa = 1,vtk_nhexa
415  vtk_offset(ihexa) = ihexa*8
416  enddo
417 
418  return
419 
420 !***********************************************************************************************************************************************************************************
421  end subroutine init_vtk_mesh
422 !***********************************************************************************************************************************************************************************
423 
424 !
425 !
436 !***********************************************************************************************************************************************************************************
437  subroutine init_vtk_node(hexa,nhexa,efi_hexa_gll_glonum,gll,ngll,ngll_unique)
438 !***********************************************************************************************************************************************************************************
439 
440  use mod_global_variables, only :&
441  ig_ngll&
442  ,ig_nhexa
443 
444  implicit none
445 
446  integer, intent( in) :: nhexa
447  integer, intent( in) :: ngll
448  integer, intent( in), dimension(nhexa) :: hexa
449  integer, intent( in), dimension(IG_NGLL,IG_NGLL,IG_NGLL,ig_nhexa) :: efi_hexa_gll_glonum
450  integer, intent(out), dimension(ngll) :: gll
451  integer, intent(out) :: ngll_unique
452 
453  integer :: myhexa
454  integer :: ihexa
455  integer :: igll
456  integer :: k
457  integer :: l
458  integer :: m
459 
460 
461  do igll = 1,ngll
462 
463  gll(igll) = 0
464 
465  enddo
466 
467 
468  do ihexa = 1,nhexa
469 
470  myhexa = hexa(ihexa)
471 
472  do k = 1,ig_ngll
473  do l = 1,ig_ngll
474  do m = 1,ig_ngll
475 
476  igll = efi_hexa_gll_glonum(m,l,k,myhexa)
477 
478  gll(igll) = 1
479 
480  enddo
481  enddo
482  enddo
483 
484  enddo
485 
486  ngll_unique = 0
487 
488  do igll = 1,ngll
489 
490  ngll_unique = ngll_unique + gll(igll)
491 
492  enddo
493 
494  return
495 
496 !***********************************************************************************************************************************************************************************
497  end subroutine init_vtk_node
498 !***********************************************************************************************************************************************************************************
499 
500 
501 !
502 !
511 !***********************************************************************************************************************************************************************************
512  subroutine init_vtk_node_numbering(node,gll,nnode,ngll)
513 !***********************************************************************************************************************************************************************************
514 
515  use mod_global_variables, only :&
516  ig_ngll&
517  ,ig_nhexa
518 
519  use mod_init_memory
520 
521  implicit none
522 
523  integer, intent( in) :: ngll
524  integer, intent( in) :: nnode
525  integer, intent(out), dimension(ngll) :: node
526  integer, intent(out), dimension(:), allocatable :: gll
527 
528  integer :: ios
529  integer :: igll
530  integer :: jgll
531 
532  ios = init_array_int(gll,nnode,"gll")
533 
534  jgll = 0
535 
536  do igll = 1,ngll
537 
538  if (node(igll) == 1) then
539 
540  jgll = jgll + 1
541 
542  node(igll) = jgll - 1
543 
544  gll(jgll) = igll
545 
546  endif
547 
548  enddo
549 
550  return
551 
552 !***********************************************************************************************************************************************************************************
553  end subroutine init_vtk_node_numbering
554 !***********************************************************************************************************************************************************************************
555 
556 !
557 !
569 !***********************************************************************************************************************************************************************************
570  subroutine init_vtk_hexa_connectivity(hexa,efi_nhexa,node,ngll,efi_hexa_gll_glonum,vtk_hexa,vtk_nhexa)
571 !***********************************************************************************************************************************************************************************
572 
573  use mod_global_variables, only :&
574  ig_ngll&
575  ,ig_ndof&
576  ,ig_nhexa
577 
578  use mod_init_memory
579 
580  implicit none
581 
582  integer, intent( in) :: efi_nhexa
583  integer, intent( in) :: ngll
584  integer, intent( in), dimension(ngll) :: node
585  integer, intent( in), dimension(efi_nhexa) :: hexa
586  integer, intent( in), dimension(IG_NGLL,IG_NGLL,IG_NGLL,ig_nhexa) :: efi_hexa_gll_glonum
587  integer, intent(out), dimension(:), allocatable :: vtk_hexa
588  integer, intent(out) :: vtk_nhexa
589 
590  integer, dimension(8) :: mygll
591  integer :: ihexa
592  integer :: jhexa
593  integer :: myhexa
594  integer :: k
595  integer :: l
596  integer :: m
597  integer :: ios
598 
599  vtk_nhexa = efi_nhexa*(ig_ngll-1)**ig_ndof
600  ios = init_array_int(vtk_hexa,vtk_nhexa*8,"vtk_hexa")
601  jhexa = 0
602 
603  do ihexa = 1,efi_nhexa
604 
605  myhexa = hexa(ihexa)
606 
607  do k = 1,ig_ngll-1
608  do l = 1,ig_ngll-1
609  do m = 1,ig_ngll-1
610 
611  jhexa = jhexa + 1
612 
613  mygll(1) = efi_hexa_gll_glonum(m ,l ,k ,myhexa)
614  mygll(2) = efi_hexa_gll_glonum(m+1,l ,k ,myhexa)
615  mygll(3) = efi_hexa_gll_glonum(m+1,l+1,k ,myhexa)
616  mygll(4) = efi_hexa_gll_glonum(m ,l+1,k ,myhexa)
617  mygll(5) = efi_hexa_gll_glonum(m ,l ,k+1,myhexa)
618  mygll(6) = efi_hexa_gll_glonum(m+1,l ,k+1,myhexa)
619  mygll(7) = efi_hexa_gll_glonum(m+1,l+1,k+1,myhexa)
620  mygll(8) = efi_hexa_gll_glonum(m ,l+1,k+1,myhexa)
621 
622  vtk_hexa((jhexa-1)*8+1) = node(mygll(1))
623  vtk_hexa((jhexa-1)*8+2) = node(mygll(2))
624  vtk_hexa((jhexa-1)*8+3) = node(mygll(3))
625  vtk_hexa((jhexa-1)*8+4) = node(mygll(4))
626  vtk_hexa((jhexa-1)*8+5) = node(mygll(5))
627  vtk_hexa((jhexa-1)*8+6) = node(mygll(6))
628  vtk_hexa((jhexa-1)*8+7) = node(mygll(7))
629  vtk_hexa((jhexa-1)*8+8) = node(mygll(8))
630 
631  enddo
632  enddo
633  enddo
634 
635  enddo
636 
637  return
638 
639 !***********************************************************************************************************************************************************************************
640  end subroutine init_vtk_hexa_connectivity
641 !***********************************************************************************************************************************************************************************
642 
643 
644 !
645 !
648 !***********************************************************************************************************************************************************************************
650 !***********************************************************************************************************************************************************************************
651 
652  use mod_global_variables, only :&
653  ig_idt&
654  ,cg_prefix&
655  ,cg_myrank&
656  ,rg_gll_displacement&
657  ,rg_gll_velocity&
658  ,rg_gll_acceleration&
659  ,ig_vtk_nhexa_snapshot&
660  ,ig_snapshot_volume_saving_incr&
661  ,lg_snapshot_volume_displacement&
662  ,lg_snapshot_volume_velocity&
663  ,lg_snapshot_volume_acceleration
664 
665  implicit none
666 
667  character(len=255) :: fname
668  character(len= 6) :: csnapshot
669 
670 
671  write(csnapshot,'(i6.6)') ig_idt
672 
673 !
674 !
675 !---->snapshot displacement
676  if ( lg_snapshot_volume_displacement .and. (mod(ig_idt-1,ig_snapshot_volume_saving_incr) == 0) .and. (ig_vtk_nhexa_snapshot > 0) ) then
677 
678  fname = trim(cg_prefix)//".volume.snapshot.uxyz."//trim(csnapshot)//".cpu."//trim(cg_myrank)//".vtu"
679 
680  call write_snapshot_volume_vtk(fname,rg_gll_displacement(1,:),"ux",rg_gll_displacement(2,:),"uy",rg_gll_displacement(3,:),"uz")
681 
682  endif
683 
684 !
685 !
686 !---->snapshot velocity
687  if ( lg_snapshot_volume_velocity .and. (mod(ig_idt-1,ig_snapshot_volume_saving_incr) == 0) .and. (ig_vtk_nhexa_snapshot > 0) ) then
688 
689  fname = trim(cg_prefix)//".volume.snapshot.vxyz."//trim(csnapshot)//".cpu."//trim(cg_myrank)//".vtu"
690 
691  call write_snapshot_volume_vtk(fname,rg_gll_velocity(1,:),"vx",rg_gll_velocity(2,:),"vy",rg_gll_velocity(3,:),"vz")
692 
693  endif
694 
695 !
696 !
697 !---->snapshot acceleration
698  if ( lg_snapshot_volume_acceleration .and. (mod(ig_idt-1,ig_snapshot_volume_saving_incr) == 0) .and. (ig_vtk_nhexa_snapshot > 0) ) then
699 
700  fname = trim(cg_prefix)//".volume.snapshot.axyz."//trim(csnapshot)//".cpu."//trim(cg_myrank)//".vtu"
701 
702  call write_snapshot_volume_vtk(fname,rg_gll_acceleration(1,:),"ax",rg_gll_acceleration(2,:),"ay",rg_gll_acceleration(3,:),"az")
703 
704  endif
705 
706  return
707 
708 !***********************************************************************************************************************************************************************************
709  end subroutine write_snapshot_volume
710 !***********************************************************************************************************************************************************************************
711 
712 
713 !
714 !
724 !***********************************************************************************************************************************************************************************
725  subroutine write_snapshot_volume_vtk(fname,gll_var_x,name_var_x,gll_var_y,name_var_y,gll_var_z,name_var_z)
726 !***********************************************************************************************************************************************************************************
727 
728  use mod_global_variables, only:&
729  cg_prefix&
730  ,cg_myrank&
731  ,ig_ndof&
732  ,ig_ngll_total&
733  ,rg_gll_displacement&
734  ,rg_gll_velocity&
735  ,rg_gll_acceleration&
736  ,ig_snapshot_volume_saving_incr&
737  ,ig_vtk_hexa_conn_snapshot&
738  ,ig_vtk_node_gll_glonum_snapshot&
739  ,rg_vtk_node_x_snapshot&
740  ,rg_vtk_node_y_snapshot&
741  ,rg_vtk_node_z_snapshot&
742  ,ig_vtk_cell_type_snapshot&
743  ,ig_vtk_offset_snapshot&
744  ,ig_vtk_nhexa_snapshot&
745  ,ig_vtk_nnode_snapshot
746 
747  use mod_gll_value
748 
749  use mod_vtk_io
750 
751  implicit none
752 
753  character(len=255) , intent(in) :: fname
754  real , optional, intent(in), dimension(ig_ngll_total) :: gll_var_x
755  real , optional, intent(in), dimension(ig_ngll_total) :: gll_var_y
756  real , optional, intent(in), dimension(ig_ngll_total) :: gll_var_z
757 
758  character(len=*) , optional, intent(in) :: name_var_x
759  character(len=*) , optional, intent(in) :: name_var_y
760  character(len=*) , optional, intent(in) :: name_var_z
761 
762  real, dimension(ig_vtk_nnode_snapshot) :: vtk_var_x
763  real, dimension(ig_vtk_nnode_snapshot) :: vtk_var_y
764  real, dimension(ig_vtk_nnode_snapshot) :: vtk_var_z
765 
766  integer :: ios
767 
768 !
769 !---->transfer GLL nodes values to VTK nodes
770  if (present(gll_var_x)) call get_node_gll_value(gll_var_x,ig_vtk_node_gll_glonum_snapshot,vtk_var_x)
771  if (present(gll_var_y)) call get_node_gll_value(gll_var_y,ig_vtk_node_gll_glonum_snapshot,vtk_var_y)
772  if (present(gll_var_z)) call get_node_gll_value(gll_var_z,ig_vtk_node_gll_glonum_snapshot,vtk_var_z)
773 
774 !
775 !---->initialize VTK file
776 
777  ios = vtk_ini_xml( &
778  output_format = 'BINARY' &
779  ,filename = trim(fname) &
780  ,mesh_topology = 'UnstructuredGrid' )
781 
782  ios = vtk_geo_xml( &
783  nn = ig_vtk_nnode_snapshot &
784  ,nc = ig_vtk_nhexa_snapshot &
785  ,x = rg_vtk_node_x_snapshot &
786  ,y = rg_vtk_node_y_snapshot &
787  ,z = rg_vtk_node_z_snapshot )
788 
789  ios = vtk_con_xml( &
790  nc = ig_vtk_nhexa_snapshot &
791  ,connect = ig_vtk_hexa_conn_snapshot &
792  ,offset = ig_vtk_offset_snapshot &
793  ,cell_type = ig_vtk_cell_type_snapshot )
794 
795  ios = vtk_dat_xml( &
796  var_location = 'NODE' &
797  ,var_block_action = 'OPEN' )
798 
799 !
800 !---->write variable values to VTK XML file
801  if (present(gll_var_x)) then
802 
803  ios = vtk_var_xml( &
804  nc_nn = ig_vtk_nnode_snapshot &
805  ,varname = name_var_x &
806  ,var = vtk_var_x )
807 
808  endif
809 
810  if (present(gll_var_y)) then
811 
812  ios = vtk_var_xml( &
813  nc_nn = ig_vtk_nnode_snapshot &
814  ,varname = name_var_y &
815  ,var = vtk_var_y )
816 
817  endif
818 
819  if (present(gll_var_z)) then
820 
821  ios = vtk_var_xml( &
822  nc_nn = ig_vtk_nnode_snapshot &
823  ,varname = name_var_z &
824  ,var = vtk_var_z )
825 
826  endif
827 
828 !
829 !---->close VTK XML file
830  ios = vtk_dat_xml( &
831  var_location = 'NODE' &
832  ,var_block_action = 'CLOSE')
833 
834  ios = vtk_geo_xml()
835 
836  ios = vtk_end_xml()
837 
838  return
839 
840 !***********************************************************************************************************************************************************************************
841  end subroutine write_snapshot_volume_vtk
842 !***********************************************************************************************************************************************************************************
843 
844 !
845 !
848 !***********************************************************************************************************************************************************************************
849  subroutine write_medium()
850 !***********************************************************************************************************************************************************************************
851 
852  use mod_global_variables, only :&
853  ig_ngll_total&
854  ,ig_nhexa&
855  ,ig_ngll&
856  ,cg_prefix&
857  ,cg_myrank&
858  ,ig_vtk_nhexa_snapshot&
859  ,ig_hexa_gll_glonum&
860  ,rg_hexa_gll_rhovs2&
861  ,rg_hexa_gll_rhovp2&
862  ,rg_hexa_gll_rho
863 
864  implicit none
865 
866  real, dimension(ig_ngll_total) :: gll_vs
867  real, dimension(ig_ngll_total) :: gll_vp
868  real, dimension(ig_ngll_total) :: gll_rho
869 
870  integer :: ihexa
871  integer :: igll
872  integer :: k
873  integer :: l
874  integer :: m
875 
876  character(len=255) :: fname
877 
878  do ihexa = 1,ig_nhexa
879 
880  do k = 1,ig_ngll
881  do l = 1,ig_ngll
882  do m = 1,ig_ngll
883 
884  igll = ig_hexa_gll_glonum(m,l,k,ihexa)
885 
886  gll_vs(igll) = sqrt(rg_hexa_gll_rhovs2(m,l,k,ihexa)/rg_hexa_gll_rho(m,l,k,ihexa))
887  gll_vp(igll) = sqrt(rg_hexa_gll_rhovp2(m,l,k,ihexa)/rg_hexa_gll_rho(m,l,k,ihexa))
888  gll_rho(igll) = rg_hexa_gll_rho(m,l,k,ihexa)
889 
890  enddo
891  enddo
892  enddo
893 
894  enddo
895 
896 
897  if (ig_vtk_nhexa_snapshot > 0) then
898 
899  fname = trim(cg_prefix)//".medium.cpu."//trim(cg_myrank)//".vtu"
900 
901  call write_snapshot_volume_vtk(fname,gll_vs,"vs",gll_vp,"vp",gll_rho,"rho")
902 
903  endif
904 
905  return
906 
907 !***********************************************************************************************************************************************************************************
908  end subroutine write_medium
909 !***********************************************************************************************************************************************************************************
910 
911 
912 !
913 !
918 !***********************************************************************************************************************************************************************************
920 !***********************************************************************************************************************************************************************************
921 
922  use mod_global_variables, only :&
923  lg_snapshot_volume_displacement&
924  ,lg_snapshot_volume_velocity&
925  ,lg_snapshot_volume_acceleration
926 
927  implicit none
928 
929  if (lg_snapshot_volume_displacement) then
930  call collection_vtk_volume("volume.snapshot.uxyz")
931  endif
932 
933  if (lg_snapshot_volume_velocity) then
934  call collection_vtk_volume("volume.snapshot.vxyz")
935  endif
936 
937  if (lg_snapshot_volume_acceleration) then
938  call collection_vtk_volume("volume.snapshot.axyz")
939  endif
940 
941  return
942 !***********************************************************************************************************************************************************************************
943  end subroutine write_collection_vtk_vol
944 !***********************************************************************************************************************************************************************************
945 
946 
947 !
948 !
952 !***********************************************************************************************************************************************************************************
953  subroutine collection_vtk_volume(varname)
954 !***********************************************************************************************************************************************************************************
955 
956  use mpi
957 
958  use mod_global_variables, only :&
959  get_newunit&
960  ,cg_prefix&
961  ,rg_dt&
962  ,ig_ndt&
963  ,ig_ncpu&
964  ,ig_myrank&
965  ,ig_snapshot_volume_saving_incr&
966  ,ig_vtk_nhexa_snapshot
967 
968  implicit none
969 
970  character(len=*), intent(in) :: varname
971 
972  real :: time
973 
974  integer, dimension(ig_ncpu) :: vtk_nhexa
975  integer :: icpu
976  integer :: istep
977  integer :: myunit
978  integer :: ios
979 
980  character(len=255) :: fname
981  character(len=6 ) :: csnapshot
982  character(len=6 ) :: crank
983 
984 
985 !
986 !---->cpuO gets the number of hexa in VTK's mesh of other cpus
987  call mpi_gather(ig_vtk_nhexa_snapshot,1,mpi_integer,vtk_nhexa,1,mpi_integer,0,mpi_comm_world,ios)
988 
989  if (ig_myrank == 0) then
990 
991  open(unit=get_newunit(myunit),file=trim(cg_prefix)//".collection."//trim(varname)//".pvd")
992 
993  write(unit=myunit,fmt='(a)') "<?xml version=""1.0""?>"
994  write(unit=myunit,fmt='(a)') "<VTKFile type=""Collection"" version=""0.1"" byte_order=""BigEndian"">"
995  write(unit=myunit,fmt='(a)') " <Collection>"
996 
997  do istep = 1,ig_ndt,ig_snapshot_volume_saving_incr
998 
999  write(csnapshot,'(I6.6)') istep
1000 
1001  time = real(istep-1)*rg_dt
1002 
1003  do icpu = 1,ig_ncpu
1004 
1005  if (vtk_nhexa(icpu) > 0) then
1006 
1007  write(crank,'(I6.6)') icpu-1
1008 
1009  fname = trim(cg_prefix)//"."//trim(varname)//"."//trim(csnapshot)//".cpu."//trim(crank)//".vtu"
1010 
1011  write(unit=myunit,fmt='(a,E14.7,3a)') " <DataSet timestep=""",time,""" group="""" part=""0"" file=""",trim(fname),"""/>"
1012 
1013  endif
1014 
1015  enddo
1016 
1017  enddo
1018 
1019  write(unit=myunit,fmt='(a)') " </Collection>"
1020  write(unit=myunit,fmt='(a)') "</VTKFile>"
1021 
1022  close(myunit)
1023 
1024  endif
1025 
1026  return
1027 !***********************************************************************************************************************************************************************************
1028  end subroutine collection_vtk_volume
1029 !***********************************************************************************************************************************************************************************
1030 
1031 end module mod_snapshot_volume
This module contains subroutines to allocate arrays and to compute an approximation of the total RAM ...
subroutine, private init_vtk_node_numbering(node, gll, nnode, ngll)
This subroutine generates i) VTK nodes&#39; numbering from 0 to n-1 and ii) indirection array from VTK&#39;s ...
integer(i4p) function, public vtk_end_xml()
The VTK_END_XML function finalizes opened files.
integer(i4p) function, public vtk_con_xml(NC, connect, offset, cell_type)
The VTK_CON_XML function must be called when unstructured grid topology is used. It saves the connect...
subroutine, private write_medium()
This subroutine prepares the arrays for writing elastic properties of the medium. ...
integer(i4p) function, public vtk_ini_xml(output_format, filename, mesh_topology, nx1, nx2, ny1, ny2, nz1, nz2)
The VTK_INI_XML function is used for initializing file. This function must be the first to be called...
subroutine, private init_vtk_hexa_connectivity(hexa, efi_nhexa, node, ngll, efi_hexa_gll_glonum, vtk_hexa, vtk_nhexa)
This subroutine creates the connectivity array for defining VTK hexahedron elements from VTK geometri...
This module contains subroutines to compute and write snapshots of a mesh composed by hexahedron elem...
overloading of VTK_GEO_XML
This module defines all global variables of EFISPEC3D. Scalar variables are initialized directly in t...
Interface init_array_real to redirect allocation to n-rank arrays.
subroutine, public write_snapshot_volume_vtk(fname, gll_var_x, name_var_x, gll_var_y, name_var_y, gll_var_z, name_var_z)
This subroutine writes volume snapshot in a VTK XML file that contains a maximum of three vectors...
Interface for redirection to subroutines mod_gll_value::get_node_gll_value_x or mod_gll_value::get_no...
integer(i4p) function, public vtk_dat_xml(var_location, var_block_action)
The VTK_DAT_XML function opens or closes CellData/PointData structures.
This module contains subroutines to compute global -coordinates of a given point in the physical doma...
subroutine, public compute_hexa_point_coord(ihexa, xi, eta, zeta, x, y, z)
subroutine to compute -coordinates of a point knowing to which hexahedron element it belongs and its ...
This module programmed by S. Zaghi (https://github.com/szaghi/Lib_VTK_IO) contains functions to write...
subroutine, public write_collection_vtk_vol()
This subroutine selects which ParaView collection files should be written (displacement, velocity and/or acceleration collections) depending on value of variables mod_global_variables::lg_snapshot_volume_displacement, mod_global_variables::lg_snapshot_volume_velocity and mod_global_variables::lg_snapshot_volume_acceleration.
subroutine, private collection_vtk_volume(varname)
This subroutine write Paraview collection file *.pvd of VTK XML volume snapshot files *...
subroutine, private get_efi_hexa(efi_hexa, efi_nhexa, xmin, xmax, ymin, ymax, zmin, zmax)
This subroutine finds and stores the number of EFISPEC&#39;s hexahedron elements used for generating VTK ...
This module contains subroutines to get GLL nodes values from global GLL nodes numbering.
subroutine, public init_snapshot_volume()
This subroutine intializes the VTK mesh used for volume snapshots by calling &#39;get_efi_hexa&#39; and &#39;init...
subroutine, private init_vtk_mesh(efi_hexa, efi_hexa_gll_glonum, vtk_hexa_conn, vtk_node_x, vtk_node_y, vtk_node_z, vtk_node_gll_glonum, vtk_cell_type, vtk_offset, vtk_nhexa, vtk_nnode)
This subroutines manages the generation of a VTK mesh from spectral elements GLL nodes.
Interface init_array_int to redirect allocation to n-rank arrays.
overloading of VTK_VAR_XML
subroutine, private init_vtk_node(hexa, nhexa, efi_hexa_gll_glonum, gll, ngll, ngll_unique)
This subroutine finds which GLL nodes to use for generating a VTK mesh.
subroutine, public write_snapshot_volume()
This subroutine manages output format for volume snapshot. For now, only VTK format is available...