MPI-AMRVAC  3.0
The MPI - Adaptive Mesh Refinement - Versatile Advection Code
mod_slice.t
Go to the documentation of this file.
1 !> Writes D-1 slice, can do so in various formats, depending on slice_type
2 module mod_slice
4  implicit none
5 
6  !> Maximum number of slices
7  integer, parameter :: nslicemax=1000
8 
9  !> Slice coordinates, see @ref slices.md
10  double precision :: slicecoord(nslicemax)
11 
12  !> the file number of slices
13  integer :: slicenext
14 
15  !> Number of slices to output
16  integer :: nslices
17 
18  !> The slice direction for each slice
19  integer :: slicedir(nslicemax)
20 
21  !> choose data type of slice: vtu, vtuCC, dat, or csv
22  character(len=std_len) :: slice_type
23 
24  !> tag for MPI message
25  integer, private :: itag
26 
27 contains
28 
29  subroutine write_slice
31  ! Writes a D-1 slice
32  ! by Oliver Porth
33  ! 22.Nov 2011
34  integer :: islice
35 
36  do islice=1,nslices
37  call put_slice(slicedir(islice),slicecoord(islice))
38  end do
39 
41  end subroutine write_slice
42 
43  subroutine put_slice(dir,xslice)
46  ! Writes a D-1 slice
47  ! For ONED simulations, output will be appended to one csv-file per slice
48  ! slices are sensitive to the saveprim switch and
49  ! can contain auxiliary io variables (nwauxio)
50  ! Thus csv or vtu(CC)-files with primitive variables are obtained.
51  ! when slice_type='dat', we save a D-1 dat file for potential restarts
52  ! by Oliver Porth
53  ! 22.Nov 2011
54  integer, intent(in) :: dir
55  double precision, intent(in) :: xslice
56  ! .. local ..
57  integer :: Njgrid, jgrid
58  integer, dimension(ndim-1) :: ixsubGlo, ixsubGhi
59  integer, dimension(ndim-1) :: ixsubMlo, ixsubMhi
60  integer :: size_subblock_io, nx^D, slice_fh, nwexpand
61  integer :: type_subblock_io, type_subblockC_io, type_subblock_x_io, type_subblockC_x_io
62  integer, dimension(ndim) :: sizes, subsizes, start
63  double precision,dimension(0:nw+nwauxio) :: normconv
64 
65  ! Preamble:
66  nx^d=ixmhi^d-ixmlo^d+1;
67  slice_fh=unitslice
68 
69  if (ndim==1) then
70  nwexpand = nwauxio
71  else
72  if(slice_type/='dat')then
73  nwexpand = nwauxio
74  else
75  nwexpand = 0
76  endif
77  end if
78 
79  ! Do a last consistency check:
80  select case(dir)
81  {case(^d)
82  if(xslice<xprobmin^d.or.xslice>xprobmax^d) &
83  call mpistop("slice out of bounds")
84  \}
85  end select
86 
87  ! Traverse the forest and fill nodes:
88  call select_slice(dir,xslice,.false.,slice_fh,normconv)
89 
90  ! Create the MPI-datatype and select indices:
91  {^ifthreed
92  select case(dir)
93  case (1)
94  sizes(1) = ixghi2; sizes(2) = ixghi3;
95  ixsubglo(1) = ixglo2; ixsubglo(2) = ixglo3;
96  ixsubghi(1) = ixghi2; ixsubghi(2) = ixghi3;
97  subsizes(1)=nx2;subsizes(2)=nx3;
98  start(1)=ixmlo2-1;start(2)=ixmlo3-1;
99  size_subblock_io=nx2*nx3*(nw+nwexpand)*size_double
100  case (2)
101  sizes(1) = ixghi1; sizes(2) = ixghi3;
102  ixsubglo(1) = ixglo1; ixsubglo(2) = ixglo3;
103  ixsubghi(1) = ixghi1; ixsubghi(2) = ixghi3;
104  subsizes(1)=nx1;subsizes(2)=nx3;
105  start(1)=ixmlo1-1;start(2)=ixmlo3-1;
106  size_subblock_io=nx1*nx3*(nw+nwexpand)*size_double
107  case (3)
108  ixsubglo(1) = ixglo1; ixsubglo(2) = ixglo2;
109  ixsubghi(1) = ixghi1; ixsubghi(2) = ixghi2;
110  sizes(1) = ixghi1; sizes(2) = ixghi2;
111  subsizes(1)=nx1;subsizes(2)=nx2;
112  start(1)=ixmlo1-1;start(2)=ixmlo2-1;
113  size_subblock_io=nx1*nx2*(nw+nwexpand)*size_double
114  case default
115  call mpistop("slice direction not clear in put_slice")
116  end select
117  }
118  {^iftwod
119  select case(dir)
120  case (1)
121  ixsubglo(1) = ixglo2; ixsubghi(1) = ixghi2;
122  sizes(1) = ixghi2
123  subsizes(1)=nx2
124  start(1)=ixmlo2-1
125  size_subblock_io=nx2*(nw+nwexpand)*size_double
126  case (2)
127  ixsubglo(1) = ixglo1; ixsubghi(1) = ixghi1;
128  sizes(1) = ixghi1
129  subsizes(1)=nx1
130  start(1)=ixmlo1-1
131  size_subblock_io=nx1*(nw+nwexpand)*size_double
132  case default
133  call mpistop("slice direction not clear in put_slice")
134  end select
135  }
136  {^ifoned
137  size_subblock_io=(nw+nwexpand)*size_double
138  }
139 
140  {^nooned
141  {^de&ixsubmlo(^de-1) = ixsubglo(^de-1)+nghostcells;}
142  {^de&ixsubmhi(^de-1) = ixsubghi(^de-1)-nghostcells;}
143  }
144 
145  sizes(ndim)=nw+nwexpand
146  subsizes(ndim)=nw+nwexpand
147  start(ndim)=0
148 
149  ! Types for center variables:
150  call mpi_type_create_subarray(ndim,sizes,subsizes,start, &
151  mpi_order_fortran,mpi_double_precision, &
152  type_subblock_io,ierrmpi)
153  call mpi_type_commit(type_subblock_io,ierrmpi)
154 
155  sizes(ndim)=^nd
156  subsizes(ndim)=^nd
157  start(ndim)=0
158  call mpi_type_create_subarray(ndim,sizes,subsizes,start, &
159  mpi_order_fortran,mpi_double_precision, &
160  type_subblock_x_io,ierrmpi)
161  call mpi_type_commit(type_subblock_x_io,ierrmpi)
162 
163 
164  ! Types for corner variables:
165  subsizes(1:ndim-1) = subsizes(1:ndim-1) + 1
166  start(1:ndim-1) = start(1:ndim-1) - 1
167  sizes(ndim)=nw+nwexpand
168  subsizes(ndim)=nw+nwexpand
169  start(ndim)=0
170 
171  call mpi_type_create_subarray(ndim,sizes,subsizes,start, &
172  mpi_order_fortran,mpi_double_precision, &
173  type_subblockc_io,ierrmpi)
174  call mpi_type_commit(type_subblockc_io,ierrmpi)
175 
176  sizes(ndim)=^nd
177  subsizes(ndim)=^nd
178  start(ndim)=0
179  call mpi_type_create_subarray(ndim,sizes,subsizes,start, &
180  mpi_order_fortran,mpi_double_precision, &
181  type_subblockc_x_io,ierrmpi)
182  call mpi_type_commit(type_subblockc_x_io,ierrmpi)
183 
184 
185  ! local number of sub-grids:
187 
188  ! Now output using various schemes:
189  if (ndim==1) then
190  call put_slice_zerod
191  else
192  select case(slice_type)
193  case ('csv')
194  call put_slice_csv
195  case ('dat')
196  call put_slice_dat
197  case ('vtu', 'vtuCC')
198  call put_slice_vtu
199  end select
200  end if
201 
202  do jgrid=1,njgrid
203  call dealloc_subnode(jgrid)
204  end do
205 
206  call mpi_type_free(type_subblock_io,ierrmpi)
207  call mpi_type_free(type_subblock_x_io,ierrmpi)
208  call mpi_type_free(type_subblockc_io,ierrmpi)
209  call mpi_type_free(type_subblockc_x_io,ierrmpi)
210 
211 
212  contains
213 
214  subroutine put_slice_vtu
215 
216  use mod_calculate_xw
217  character(len=1024) :: filename, xlabel
218  character(len=79) :: xxlabel
219  logical :: fileopen
220  character(len=name_len) :: wnamei(1:nw+nwauxio),xandwnamei(1:ndim+nw+nwauxio)
221  character(len=1024) :: outfilehead
222  integer :: status(MPI_STATUS_SIZE), ipe
223 
224  if (mype==0) then
225 
226  inquire(slice_fh,opened=fileopen)
227  if(.not.fileopen)then
228  ! generate filename:
229  write(xlabel,"(D9.2)")xslice
230  xxlabel=trim(xlabel)
231  if(xslice>=zero)then
232  write(xxlabel(1:1),"(a)") "+"
233  endif
234  write(filename,"(a,i1.1,a,i4.4,a)") trim(base_filename)//'_d',dir,'_x'//trim(xxlabel)//'_n',slicenext,'.vtu'
235  open(slice_fh,file=filename,status='unknown',form='formatted')
236  end if
237  ! get and write the header:
238  call getheadernames(wnamei,xandwnamei,outfilehead)
239  ! generate xml header
240  write(slice_fh,'(a)')'<?xml version="1.0"?>'
241  write(slice_fh,'(a)',advance='no') '<VTKFile type="UnstructuredGrid"'
242  if(type_endian==1)then
243  write(slice_fh,'(a)')' version="0.1" byte_order="LittleEndian">'
244  else
245  write(slice_fh,'(a)')' version="0.1" byte_order="BigEndian">'
246  endif
247  write(slice_fh,'(a)')' <UnstructuredGrid>'
248  write(slice_fh,'(a)')'<FieldData>'
249  write(slice_fh,'(2a)')'<DataArray type="Float32" Name="TIME" ',&
250  'NumberOfTuples="1" format="ascii">'
251  write(slice_fh,*) real(global_time*time_convert_factor)
252  write(slice_fh,'(a)')'</DataArray>'
253  write(slice_fh,'(a)')'</FieldData>'
254 
255  ! write to file:
256  do jgrid=1, njgrid
257  call write_slice_vtk(jgrid,slice_fh,wnamei)
258  end do
259 
260  ! create a recv buffer using allocate, will be deallocated at the end of the routine:
261  call alloc_subnode(njgrid+1,dir,nwauxio)
262 
263  end if
264 
265  ! Also communicate the normconv array since processor zero might not have it yet:
266  if (npe>1) then
267  do ipe=1,npe-1
268  do jgrid=1,morton_sub_stop(ipe)-morton_sub_start(ipe)+1
269  itag=morton_sub_start(ipe)+jgrid-1
270  itag=itag*5
271  if (ipe == mype ) then
272  call mpi_send(ps_sub(jgrid)%x,1,type_subblock_x_io,0,itag,icomm,ierrmpi)
273  call mpi_send(ps_sub(jgrid)%w,1,type_subblock_io,0,itag+1,icomm,ierrmpi)
274  call mpi_send(ps_sub(jgrid)%xC,1,type_subblockc_x_io,0,itag+2,icomm,ierrmpi)
275  call mpi_send(ps_sub(jgrid)%wC,1,type_subblockc_io,0,itag+3,icomm,ierrmpi)
276  call mpi_send(normconv,nw+nwauxio+1,mpi_double_precision,0,itag+4,icomm,ierrmpi)
277  end if
278  if (mype == 0) then
279  call mpi_recv(ps_sub(njgrid+1)%x,1,type_subblock_x_io,ipe,itag,icomm,status,ierrmpi)
280  call mpi_recv(ps_sub(njgrid+1)%w,1,type_subblock_io,ipe,itag+1,icomm,status,ierrmpi)
281  call mpi_recv(ps_sub(njgrid+1)%xC,1,type_subblockc_x_io,ipe,itag+2,icomm,status,ierrmpi)
282  call mpi_recv(ps_sub(njgrid+1)%wC,1,type_subblockc_io,ipe,itag+3,icomm,status,ierrmpi)
283  call mpi_recv(normconv,nw+nwauxio+1,mpi_double_precision,ipe,&
284  itag+4,icomm,status,ierrmpi)
285  call write_slice_vtk(njgrid+1,slice_fh,wnamei)
286  end if
287  end do
288  end do
289  endif
290 
291  if (mype==0) then
292 
293  write(slice_fh,'(a)')'</UnstructuredGrid>'
294  write(slice_fh,'(a)')'</VTKFile>'
295  close(slice_fh)
296  call dealloc_subnode(njgrid+1)
297 
298  end if
299 
300  end subroutine put_slice_vtu
301 
302  subroutine write_slice_vtk(jgrid,slice_fh,wnamei)
303 
304  ! this only works for 2D and 3D, 1D reduction (line to point) not allowed
305  integer, intent(in) :: jgrid, slice_fh
306  character(len=name_len), intent(in) :: wnamei(1:nw+nwauxio)
307  ! This remainder part only for more than 1D, but nesting with NOONED gives problems
308  {#IFNDEF D1
309  ! .. local ..
310  integer :: ixC^L, ixCC^L, nc, np, iw
311  integer :: nx^DM, nxC^DM, icell, ix^DM
312  double precision :: x_VTK(1:3)
313  integer :: VTK_type
314  double precision, parameter :: minvalue = 1.0d-99, maxvalue = 1.0d+99
315 
316  {^dm&ixccmin^dm = ixsubmlo(^dm);}
317  {^dm&ixccmax^dm = ixsubmhi(^dm);}
318  {^dm&ixcmin^dm = ixsubmlo(^dm)-1;}
319  {^dm&ixcmax^dm = ixsubmhi(^dm);}
320 
321  nx^dm=ixccmax^dm-ixccmin^dm+1;
322  nxc^dm=nx^dm+1;
323  nc={nx^dm*} ! Number of cells per subgrid
324  np={nxc^dm*} ! Number of corner points per subgrid
325 
326  ! we write out every grid as one VTK PIECE
327  write(slice_fh,'(a,i7,a,i7,a)') &
328  '<Piece NumberOfPoints="',np,'" NumberOfCells="',nc,'">'
329 
330  !==============================
331  ! celldata or pointdata?
332  !==============================
333  select case(slice_type)
334 
335  case('vtuCC') ! celldata
336  write(slice_fh,'(a)')'<CellData>'
337  do iw=1,nw+nwauxio
338  if(iw<=nw) then
339  if(.not.w_write(iw)) cycle
340  endif
341  write(slice_fh,'(a,a,a)')&
342  '<DataArray type="Float64" Name="',trim(wnamei(iw)),'" format="ascii">'
343  write(slice_fh,'(200(1pe14.6))') {^dm&(|}roundoff_minmax(ps_sub(jgrid)%w(ix^dm,iw)*normconv(iw),minvalue,maxvalue),{ix^dm=ixccmin^dm,ixccmax^dm)}
344  write(slice_fh,'(a)')'</DataArray>'
345  enddo
346  write(slice_fh,'(a)')'</CellData>'
347 
348 
349  case('vtu') ! pointdata
350  write(slice_fh,'(a)')'<PointData>'
351  do iw=1,nw+nwauxio
352  if(iw<=nw) then
353  if(.not.w_write(iw)) cycle
354  endif
355  write(slice_fh,'(a,a,a)')&
356  '<DataArray type="Float64" Name="',trim(wnamei(iw)),'" format="ascii">'
357  write(slice_fh,'(200(1pe14.6))') {^dm&(|}roundoff_minmax(ps_sub(jgrid)%wC(ix^dm,iw)*normconv(iw),minvalue,maxvalue),{ix^dm=ixcmin^dm,ixcmax^dm)}
358  write(slice_fh,'(a)')'</DataArray>'
359  enddo
360  write(slice_fh,'(a)')'</PointData>'
361 
362 
363  end select
364  !==============================
365  ! Done: celldata or pointdata?
366  !==============================
367 
368  !==============================
369  ! output Cornerpoints
370  !==============================
371  write(slice_fh,'(a)')'<Points>'
372  write(slice_fh,'(a)')'<DataArray type="Float32" NumberOfComponents="3" format="ascii">'
373  ! write cell corner coordinates in a backward dimensional loop, always 3D output
374  {do ix^dmb=ixcmin^dmb,ixcmax^dmb \}
375  x_vtk(1:3)=zero;
376  x_vtk(1:ndim)=ps_sub(jgrid)%xC(ix^dm,1:ndim)*normconv(0);
377  write(slice_fh,'(3(1pe14.6))') x_vtk
378  {^dm&end do \}
379  write(slice_fh,'(a)')'</DataArray>'
380  write(slice_fh,'(a)')'</Points>'
381  !==============================
382  ! Done: output Cornerpoints
383  !==============================
384 
385  !==============================
386  ! cell Metainformation
387  !==============================
388  write(slice_fh,'(a)')'<Cells>'
389 
390  ! connectivity part
391  write(slice_fh,'(a)')'<DataArray type="Int32" Name="connectivity" format="ascii">'
392 
393  {^dm& do ix^dmb=1,nx^dmb\}
394  {^iftwod
395  write(slice_fh,'(2(i7,1x))')ix1-1,ix1
396  }{^ifthreed
397  write(slice_fh,'(4(i7,1x))')(ix2-1)*nxc1+ix1-1, &
398  (ix2-1)*nxc1+ix1,ix2*nxc1+ix1-1,ix2*nxc1+ix1
399  }{^dm& end do\}
400 
401  write(slice_fh,'(a)')'</DataArray>'
402 
403  ! offsets data array
404  write(slice_fh,'(a)')'<DataArray type="Int32" Name="offsets" format="ascii">'
405  do icell=1,nc
406  write(slice_fh,'(i7)') icell*(2**(^nd-1))
407  end do
408  write(slice_fh,'(a)')'</DataArray>'
409 
410  ! VTK cell type data array
411  write(slice_fh,'(a)')'<DataArray type="Int32" Name="types" format="ascii">'
412  ! VTK_LINE=3; VTK_PIXEL=8; VTK_VOXEL=11 -> vtk-syntax
413  {^iftwod vtk_type=3 \}
414  {^ifthreed vtk_type=8 \}
415  do icell=1,nc
416  write(slice_fh,'(i2)') vtk_type
417  enddo
418  write(slice_fh,'(a)')'</DataArray>'
419 
420  write(slice_fh,'(a)')'</Cells>'
421  !==============================
422  ! Done: cell Metainformation
423  !==============================
424  write(slice_fh,'(a)')'</Piece>'
425 
426  }
427  end subroutine write_slice_vtk
428 
429  subroutine put_slice_csv
430 
431  use mod_calculate_xw
432  character(len=1024) :: filename, xlabel
433  character(len=79) :: xxlabel
434  logical :: fileopen
435  character(len=name_len) :: wnamei(1:nw+nwauxio),xandwnamei(1:ndim+nw+nwauxio)
436  character(len=1024) :: outfilehead
437  integer :: iw, ipe, itag
438  character(len=1024) :: line
439  integer :: status(MPI_STATUS_SIZE)
440 
441  if (mype==0) then
442  inquire(slice_fh,opened=fileopen)
443  if(.not.fileopen)then
444  ! generate filename:
445  write(xlabel,"(D9.2)")xslice
446  xxlabel=trim(xlabel)
447  if(xslice>=zero)then
448  write(xxlabel(1:1),"(a)") "+"
449  endif
450  write(filename,"(a,i1.1,a,i4.4,a)") trim(base_filename)//'_d',dir,'_x'//trim(xxlabel)//'_n',slicenext,'.csv'
451  open(slice_fh,file=filename,status='unknown',form='formatted')
452  end if
453  ! get and write the header:
454  call getheadernames(wnamei,xandwnamei,outfilehead)
455  line=''
456  do iw=1,ndim+nw+nwauxio-1
457  line = trim(line)//trim(xandwnamei(iw))//', '
458  end do
459  line = trim(line)//trim(xandwnamei(ndim+nw+nwauxio))
460  write(slice_fh,'(a)')trim(line)
461  ! create a recv buffer using allocate, will be deallocated at the end of the routine:
462  call alloc_subnode(njgrid+1,dir,nwauxio)
463 
464  ! write to file:
465  do jgrid=1, njgrid
466  call put_slice_line(jgrid,slice_fh)
467  end do
468  end if
469 
470  ! Also communicate the normconv array since processor zero might not have it yet:
471  if (npe>1) then
472  do ipe=1,npe-1
473  do jgrid=1,morton_sub_stop(ipe)-morton_sub_start(ipe)+1
474  itag=morton_sub_start(ipe)+jgrid-1
475  if (ipe == mype ) then
476  call mpi_send(ps_sub(jgrid)%x,1,type_subblock_x_io,0,itag,icomm,ierrmpi)
477  call mpi_send(ps_sub(jgrid)%w,1,type_subblock_io,0,itag,icomm,ierrmpi)
478  call mpi_send(normconv,nw+nwauxio+1,mpi_double_precision,0,itag,icomm,ierrmpi)
479  end if
480  if (mype == 0) then
481  call mpi_recv(ps_sub(njgrid+1)%x,1,type_subblock_x_io,ipe,itag,icomm,status,ierrmpi)
482  call mpi_recv(ps_sub(njgrid+1)%w,1,type_subblock_io,ipe,itag,icomm,status,ierrmpi)
483  call mpi_recv(normconv,nw+nwauxio+1,mpi_double_precision,ipe,&
484  itag,icomm,status,ierrmpi)
485  call put_slice_line(njgrid+1,slice_fh)
486  end if
487  end do
488  end do
489  endif
490 
491  if (mype==0) then
492  close(slice_fh)
493  call dealloc_subnode(njgrid+1)
494  end if
495 
496  end subroutine put_slice_csv
497 
498  subroutine put_slice_line(jout,file_handle)
499  integer, intent(in) :: jout, file_handle
500  ! .. local ..
501  character(len=1024) ::line, data
502  integer :: ix^D,idir,iw
503  double precision, parameter :: minvalue = 1.0d-99, maxvalue = 1.0d+99
504 
505  {^ifthreed
506  do ix2=ixsubmlo(2),ixsubmhi(2)
507  do ix1=ixsubmlo(1),ixsubmhi(1)
508  \}
509  {^iftwod
510  do ix1=ixsubmlo(1),ixsubmhi(1)
511  }
512  ! Format the line:
513  line = ''
514  do idir=1,ndim
515  {^ifthreed
516  write(data,"(es14.6)")roundoff_minmax(ps_sub(jout)%x(ix1,ix2,idir),minvalue,maxvalue)
517  }
518  {^iftwod
519  write(data,"(es14.6)")roundoff_minmax(ps_sub(jout)%x(ix1,idir),minvalue,maxvalue)
520  }
521  {^ifoned
522  write(data,"(es14.6)")roundoff_minmax(ps_sub(jout)%x(idir),minvalue,maxvalue)
523  }
524 
525  line = trim(line)//trim(data)//', '
526  end do
527  do iw = 1,nw+nwauxio-1
528  {^ifthreed
529  write(data,"(es14.6)")roundoff_minmax(ps_sub(jout)%w(ix1,ix2,iw)*normconv(iw),minvalue,maxvalue)
530  }
531  {^iftwod
532  write(data,"(es14.6)")roundoff_minmax(ps_sub(jout)%w(ix1,iw)*normconv(iw),minvalue,maxvalue)
533  }
534  {^ifoned
535  write(data,"(es14.6)")roundoff_minmax(ps_sub(jout)%w(iw)*normconv(iw),minvalue,maxvalue)
536  }
537  line = trim(line)//trim(data)//', '
538  end do
539  {^ifthreed
540  write(data,"(es14.6)")roundoff_minmax(ps_sub(jout)%w(ix1,ix2,nw+nwauxio)*normconv(nw+nwauxio),minvalue,maxvalue)
541  }
542  {^iftwod
543  write(data,"(es14.6)")roundoff_minmax(ps_sub(jout)%w(ix1,nw+nwauxio)*normconv(nw+nwauxio),minvalue,maxvalue)
544  }
545  line = trim(line)//trim(data)
546  write(file_handle,'(a)')trim(line)
547  {^iftwod
548  end do
549  }
550  {^ifthreed
551  end do
552  end do
553  \}
554 
555  end subroutine put_slice_line
556 
557  subroutine put_slice_dat
558 
559  integer, dimension(max_blocks) :: iorequest
560  integer, dimension(MPI_STATUS_SIZE,max_blocks) :: iostatus
561  integer(kind=MPI_OFFSET_KIND) :: offset
562  integer :: nsubleafs
563  character(len=1024) :: filename, xlabel
564  character(len=79) :: xxlabel
565  integer :: amode, status(MPI_STATUS_SIZE), iwrite
566 
567  nsubleafs=morton_sub_stop(npe-1)
568  ! generate filename
569  write(xlabel,"(D9.2)")xslice
570  xxlabel=trim(xlabel)
571  if(xslice>=zero)then
572  write(xxlabel(1:1),"(a)") "+"
573  endif
574  write(filename,"(a,i1.1,a,i4.4,a)") trim(base_filename)//'_d',dir,'_x'//trim(xxlabel)//'_n',slicenext,'.dat'
575 
576  if(mype==0) then
577  open(unit=slice_fh,file=filename,status='replace')
578  close(unit=slice_fh)
579  end if
580 
581  amode=ior(mpi_mode_create,mpi_mode_wronly)
582  call mpi_file_open(icomm,filename,amode,mpi_info_null,slice_fh,ierrmpi)
583  iorequest=mpi_request_null
584  iwrite=0
585 
586  do jgrid=1,njgrid
587  iwrite=iwrite+1
588  offset=int(size_subblock_io,kind=mpi_offset_kind) &
589  *int(morton_sub_start(mype)+jgrid-2,kind=mpi_offset_kind)
590  call mpi_file_iwrite_at(slice_fh,offset,ps_sub(jgrid)%w,1,type_subblock_io, &
591  iorequest(iwrite),ierrmpi)
592  end do
593 
594  if (iwrite>0) call mpi_waitall(iwrite,iorequest,iostatus,ierrmpi)
595  call mpi_barrier(icomm, ierrmpi)
596  call mpi_file_close(slice_fh,ierrmpi)
597 
598  if (mype==0) then
599  amode=ior(mpi_mode_append,mpi_mode_wronly)
600  call mpi_file_open(mpi_comm_self,filename,amode,mpi_info_null, &
601  slice_fh,ierrmpi)
602 
603  call select_slice(dir,xslice,.true.,slice_fh,normconv)
604 
605  {call mpi_file_write(slice_fh,subsizes(^de-1),1,mpi_integer,status,ierrmpi)\}
606 ! call MPI_FILE_WRITE(slice_fh,eqpar,neqpar+nspecialpar, &
607 ! MPI_DOUBLE_PRECISION,status,ierrmpi)
608  call mpi_file_write(slice_fh,nsubleafs,1,mpi_integer,status,ierrmpi)
609  call mpi_file_write(slice_fh,levmax_sub,1,mpi_integer,status,ierrmpi)
610  call mpi_file_write(slice_fh,ndim-1,1,mpi_integer,status,ierrmpi)
611  call mpi_file_write(slice_fh,ndir,1,mpi_integer,status,ierrmpi)
612  call mpi_file_write(slice_fh,nw,1,mpi_integer,status,ierrmpi)
613 ! call MPI_FILE_WRITE(slice_fh,neqpar+nspecialpar,1,MPI_INTEGER,status,ierrmpi)
614  call mpi_file_write(slice_fh,it,1,mpi_integer,status,ierrmpi)
615  call mpi_file_write(slice_fh,global_time,1,mpi_double_precision,status,ierrmpi)
616 
617  call mpi_file_close(slice_fh,ierrmpi)
618  end if
619 
620  end subroutine put_slice_dat
621 
622  subroutine put_slice_zerod
623 
624  use mod_calculate_xw
625  integer:: iw
626  character(len=name_len) :: wnamei(1:nw+nwauxio),xandwnamei(1:ndim+nw+nwauxio)
627  character(len=1024) :: outfilehead
628  logical, save :: opened=.false.
629  character(len=1024) ::line, data
630  character(len=1024) :: filename, xlabel
631  character(len=79) :: xxlabel
632  integer :: amode, iwrite, status(MPI_STATUS_SIZE)
633 
634  {^ifoned
635  ! generate filename:
636  write(xlabel,"(D9.2)")xslice
637  xxlabel=trim(xlabel)
638  if(xslice>=zero)then
639  write(xxlabel(1:1),"(a)") "+"
640  endif
641  write(filename,"(a,i1.1,a,i4.4,a)") trim(base_filename)//'_d',dir,'_x'//trim(xxlabel)//'.csv'
642 
643  ! Open for header:
644  if (.not. opened .and. mype==0) then
645  amode=ior(mpi_mode_create,mpi_mode_wronly)
646  amode=ior(amode,mpi_mode_append)
647  call mpi_file_open(mpi_comm_self,filename,amode, &
648  mpi_info_null,slice_fh,ierrmpi)
649  ! Create the header:
650  call getheadernames(wnamei,xandwnamei,outfilehead)
651  line = 'time'
652  do iw=1,nw+nwauxio+ndim
653  line = trim(line)//', '//trim(xandwnamei(iw))
654  enddo
655  ! Write header:
656  call mpi_file_write(slice_fh,line,len_trim(line), &
657  mpi_character,status,ierrmpi)
658  call mpi_file_write(slice_fh,achar(10),1,mpi_character,status,ierrmpi)
659  call mpi_file_close(slice_fh, ierrmpi)
660  opened=.true.
661  endif ! opened
662  ! Now we let the processor, which holds the data, write - therefore barrier.
663  call mpi_barrier(icomm,ierrmpi)
664 
665  ! There should be only one processor holding the point:
666  if (njgrid > 0) then
667  amode=ior(mpi_mode_create,mpi_mode_wronly)
668  amode=ior(amode,mpi_mode_append)
669  call mpi_file_open(mpi_comm_self,filename,amode, &
670  mpi_info_null,slice_fh,ierrmpi)
671 
672  ! Format the line:
673  write(data,"(es14.6)")global_time
674  line = trim(data)
675  write(data,"(es14.6)")ps_sub(1)%x
676  line = trim(line)//', '//trim(data)
677  do iw = 1,nw+nwauxio
678  write(data,"(es14.6)")ps_sub(1)%w(iw)*normconv(iw)
679  line = trim(line)//', '//trim(data)
680  end do
681  !
682  call mpi_file_write(slice_fh,trim(line),len_trim(line), &
683  mpi_character,status,ierrmpi)
684  call mpi_file_write(slice_fh,achar(10),1,mpi_character,status,ierrmpi)
685  call mpi_file_close(slice_fh, ierrmpi)
686  end if ! Njgrid > 0
687  }
688  end subroutine put_slice_zerod
689 
690  end subroutine put_slice
691 
692  subroutine select_slice(dir,xslice,writeonly,file_handle,normconv)
695  integer, intent(in) :: dir
696  double precision, intent(in) :: xslice
697  integer, intent(in) :: file_handle
698  logical, intent(in) :: writeonly
699  double precision,dimension(0:nw+nwauxio),intent(out) :: normconv
700  ! .. local ..
701  integer :: ig^D, jgrid, slice_fh, ipe, mylevmax
702  integer, dimension(nlevelshi) :: igslice
703 
704  jgrid = 0
705  mylevmax = 0
706 
707  ! Find the global slice index for every level:
708  call get_igslice(dir,xslice,igslice)
709 
710  ! Traverse forest to find grids indicating slice:
711  {^ifthreed
712  select case(dir)
713  case (1)
714  ig1 = igslice(1)
715  do ig3=1,ng3(1)
716  do ig2=1,ng2(1)
717  call traverse_slice(tree_root(ig^d))
718  end do
719  end do
720  case (2)
721  ig2 = igslice(1)
722  do ig3=1,ng3(1)
723  do ig1=1,ng1(1)
724  call traverse_slice(tree_root(ig^d))
725  end do
726  end do
727  case (3)
728  ig3 = igslice(1)
729  do ig2=1,ng2(1)
730  do ig1=1,ng1(1)
731  call traverse_slice(tree_root(ig^d))
732  end do
733  end do
734  case default
735  call mpistop("slice direction not clear in select_slice")
736  end select
737  }
738  {^iftwod
739  select case(dir)
740  case (1)
741  ig1 = igslice(1)
742  do ig2=1,ng2(1)
743  call traverse_slice(tree_root(ig^d))
744  end do
745  case (2)
746  ig2 = igslice(1)
747  do ig1=1,ng1(1)
748  call traverse_slice(tree_root(ig^d))
749  end do
750  case default
751  call mpistop("slice direction not clear in select_slice")
752  end select
753  }
754  {^ifoned
755  ig1 = igslice(1)
756  call traverse_slice(tree_root(ig^d))
757  }
758 
759  if (.not.writeonly) then
760  ! Synchronize the levmax_sub for output (only rank 0 needs it):
761  levmax_sub = mylevmax
762  call mpi_allreduce(mpi_in_place,levmax_sub,1,mpi_integer,mpi_max,icomm,ierrmpi)
763 
764  ! Communicate the subgrid indices according to new Morton sub-sfc:
765  morton_sub_start(:) = 1
766  do ipe=0,npe-1
767  call mpi_gather(jgrid,1,mpi_integer,morton_sub_stop,1,mpi_integer,ipe,icomm,ierrmpi)
768  end do
769 
770  do ipe = 0, npe-2
772  morton_sub_stop(ipe+1) = morton_sub_stop(ipe)+morton_sub_stop(ipe+1)
773  end do
774  end if
775 
776  contains
777 
778  recursive subroutine traverse_slice(tree)
779  implicit none
780  type(tree_node_ptr) :: tree
781  integer :: ic^d
782  integer, dimension(MPI_STATUS_SIZE) :: status
783 
784  if (writeonly) then
785  call mpi_file_write(file_handle,tree%node%leaf,1,mpi_logical, &
786  status,ierrmpi)
787  end if
788 
789  if (tree%node%leaf) then
790  if (tree%node%ipe == mype.and..not.writeonly) then
791  mylevmax = max(mylevmax,tree%node%level)
792  call fill_subnode(tree%node%igrid,tree%node%active,jgrid,dir,xslice,normconv)
793  end if
794  return
795  end if
796  ! We are out for leaves now, continue for branches
797 
798  ! Get the correct child:
799  select case (dir)
800  {case (^d)
801  ic^d = igslice(tree%node%level+1) - 2 * tree%node%ig^d + 2
802  \}
803  case default
804  call mpistop("slice direction not clear in traverse_slice")
805  end select
806 
807  ! Recursively descend into the correct child branch:
808  {^ifthreed
809  select case(dir)
810  case (1)
811  do ic3=1,2
812  do ic2=1,2
813  call traverse_slice(tree%node%child(ic^d))
814  end do
815  end do
816  case (2)
817  do ic3=1,2
818  do ic1=1,2
819  call traverse_slice(tree%node%child(ic^d))
820  end do
821  end do
822  case (3)
823  do ic2=1,2
824  do ic1=1,2
825  call traverse_slice(tree%node%child(ic^d))
826  end do
827  end do
828  case default
829  call mpistop("slice direction not clear in traverse_slice")
830  end select
831  }
832  {^iftwod
833  select case(dir)
834  case (1)
835  do ic2=1,2
836  call traverse_slice(tree%node%child(ic^d))
837  end do
838  case (2)
839  do ic1=1,2
840  call traverse_slice(tree%node%child(ic^d))
841  end do
842  case default
843  call mpistop("slice direction not clear in traverse_slice")
844  end select
845  }
846  {^ifoned
847  call traverse_slice(tree%node%child(ic^d))
848  }
849 
850  end subroutine traverse_slice
851 
852  end subroutine select_slice
853 
854  subroutine fill_subnode(igrid,active,jgrid,dir,xslice,normconv)
856  use mod_calculate_xw
857  integer, intent(in) :: igrid, dir
858  integer, intent(inout) :: jgrid
859  logical, intent(in) :: active
860  double precision, intent(in) :: xslice
861  double precision,dimension(0:nw+nwauxio),intent(out) :: normconv
862  ! .. local ..
863  double precision, dimension(ixMlo^D-1:ixMhi^D,ndim) :: xC_TMP, xC
864  double precision, dimension(ixMlo^D:ixMhi^D,ndim) :: xCC_TMP, xCC
865  double precision, dimension(ixMlo^D-1:ixMhi^D,nw+nwauxio) :: wC_TMP
866  double precision, dimension(ixMlo^D:ixMhi^D,nw+nwauxio) :: wCC_TMP
867  double precision, dimension(ixG^T) :: x_save
868  integer :: ixslice, nwexpand, ixC^L, ixCC^L
869  logical :: mask(ixG^T)
870 
871  mask=.false.
872  ! increase grid-count:
873  jgrid=jgrid+1
874  ! Allocate subdim solution array:
875  if (ndim==1) then
876  nwexpand = nwauxio
877  else
878  if(slice_type/='dat')then
879  nwexpand = nwauxio
880  else
881  nwexpand = 0
882  endif
883  end if
884  call alloc_subnode(jgrid,dir,nwexpand)
885  call fill_subnode_info(igrid,jgrid,dir)
886 
887  mask(ixm^t)=.true.
888 
889  ! Now hunt for the index closest to the slice:
890  {^ifoned
891  ixslice = minloc(dabs(xslice-ps(igrid)%x(:,dir)),1,mask(:))
892  }
893  {^iftwod
894  select case (dir)
895  case (1)
896  ixslice = minloc(dabs(xslice-ps(igrid)%x(:,ixmlo2,dir)),1,mask(:,ixmlo2))
897  case (2)
898  ixslice = minloc(dabs(xslice-ps(igrid)%x(ixmlo1,:,dir)),1,mask(ixmlo1,:))
899  case default
900  call mpistop("slice direction not clear in fill_subnode")
901  end select
902  }
903  {^ifthreed
904  select case (dir)
905  case (1)
906  ixslice = minloc(dabs(xslice-ps(igrid)%x(:,ixmlo2,ixmlo3,dir)),1,mask(:,ixmlo2,ixmlo3))
907  case (2)
908  ixslice = minloc(dabs(xslice-ps(igrid)%x(ixmlo1,:,ixmlo3,dir)),1,mask(ixmlo1,:,ixmlo3))
909  case (3)
910  ixslice = minloc(dabs(xslice-ps(igrid)%x(ixmlo1,ixmlo2,:,dir)),1,mask(ixmlo1,ixmlo2,:))
911  case default
912  call mpistop("slice direction not clear in fill_subnode")
913  end select
914  }
915 
916  call calc_x(igrid,xc,xcc)
917  ! Set the coordinate to be exactly on the slice:
918  xc(:^d&,dir) = xslice
919  xcc(:^d&,dir) = xslice
920  call calc_grid(unitslice,igrid,xc,xcc,xc_tmp,xcc_tmp,wc_tmp,wcc_tmp,normconv,&
921  ixc^l,ixcc^l,.true.)
922  ! CC stands for CellCenter
923  ! C stands for Corner
924 
925  ! Fill the subdimensional solution and position array:
926  {^ifoned
927  ps_sub(jgrid)%w(1:nw+nwexpand) = wcc_tmp(ixslice,1:nw+nwexpand)
928  ps_sub(jgrid)%x(1:ndim) = xcc_tmp(ixslice,1:ndim)
929  ps_sub(jgrid)%wC(1:nw+nwexpand) = wc_tmp(ixslice,1:nw+nwexpand)
930  ps_sub(jgrid)%xC(1:ndim) = xc_tmp(ixslice,1:ndim)
931  }
932  {^iftwod
933  select case (dir)
934  case (1)
935  ps_sub(jgrid)%w(ixccmin2:ixccmax2,1:nw+nwexpand) &
936  = wcc_tmp(ixslice,ixccmin2:ixccmax2,1:nw+nwexpand)
937  ps_sub(jgrid)%x(ixccmin2:ixccmax2,1:ndim) &
938  = xcc_tmp(ixslice,ixccmin2:ixccmax2,1:ndim)
939  ps_sub(jgrid)%wC(ixcmin2:ixcmax2,1:nw+nwexpand) &
940  = wc_tmp(ixslice,ixcmin2:ixcmax2,1:nw+nwexpand)
941  ps_sub(jgrid)%xC(ixcmin2:ixcmax2,1:ndim) &
942  = xc_tmp(ixslice,ixcmin2:ixcmax2,1:ndim)
943  case (2)
944  ps_sub(jgrid)%w(ixccmin1:ixccmax1,1:nw+nwexpand) &
945  = wcc_tmp(ixccmin1:ixccmax1,ixslice,1:nw+nwexpand)
946  ps_sub(jgrid)%x(ixccmin1:ixccmax1,1:ndim) &
947  = xcc_tmp(ixccmin1:ixccmax1,ixslice,1:ndim)
948  ps_sub(jgrid)%wC(ixcmin1:ixcmax1,1:nw+nwexpand) &
949  = wc_tmp(ixcmin1:ixcmax1,ixslice,1:nw+nwexpand)
950  ps_sub(jgrid)%xC(ixcmin1:ixcmax1,1:ndim) &
951  = xc_tmp(ixcmin1:ixcmax1,ixslice,1:ndim)
952  case default
953  call mpistop("slice direction not clear in fill_subnode")
954  end select
955  }
956  {^ifthreed
957  select case (dir)
958  case (1)
959  ps_sub(jgrid)%w(ixccmin2:ixccmax2,ixccmin3:ixccmax3,1:nw+nwexpand) = &
960  wcc_tmp(ixslice,ixccmin2:ixccmax2,ixccmin3:ixccmax3,1:nw+nwexpand)
961  ps_sub(jgrid)%x(ixccmin2:ixccmax2,ixccmin3:ixccmax3,1:ndim) = &
962  xcc_tmp(ixslice,ixccmin2:ixccmax2,ixccmin3:ixccmax3,1:ndim)
963  ps_sub(jgrid)%wC(ixcmin2:ixcmax2,ixcmin3:ixcmax3,1:nw+nwexpand) = &
964  wc_tmp(ixslice,ixcmin2:ixcmax2,ixcmin3:ixcmax3,1:nw+nwexpand)
965  ps_sub(jgrid)%xC(ixcmin2:ixcmax2,ixcmin3:ixcmax3,1:ndim) = &
966  xc_tmp(ixslice,ixcmin2:ixcmax2,ixcmin3:ixcmax3,1:ndim)
967  case (2)
968  ps_sub(jgrid)%w(ixccmin1:ixccmax1,ixccmin3:ixccmax3,1:nw+nwexpand) = &
969  wcc_tmp(ixccmin1:ixccmax1,ixslice,ixccmin3:ixccmax3,1:nw+nwexpand)
970  ps_sub(jgrid)%x(ixccmin1:ixccmax1,ixccmin3:ixccmax3,1:ndim) = &
971  xcc_tmp(ixccmin1:ixccmax1,ixslice,ixccmin3:ixccmax3,1:ndim)
972  ps_sub(jgrid)%wC(ixcmin1:ixcmax1,ixcmin3:ixcmax3,1:nw+nwexpand) = &
973  wc_tmp(ixcmin1:ixcmax1,ixslice,ixcmin3:ixcmax3,1:nw+nwexpand)
974  ps_sub(jgrid)%xC(ixcmin1:ixcmax1,ixcmin3:ixcmax3,1:ndim) = &
975  xc_tmp(ixcmin1:ixcmax1,ixslice,ixcmin3:ixcmax3,1:ndim)
976  case (3)
977  ps_sub(jgrid)%w(ixccmin1:ixccmax1,ixccmin2:ixccmax2,1:nw+nwexpand) = &
978  wcc_tmp(ixccmin1:ixccmax1,ixccmin2:ixccmax2,ixslice,1:nw+nwexpand)
979  ps_sub(jgrid)%x(ixccmin1:ixccmax1,ixccmin2:ixccmax2,1:ndim) = &
980  xcc_tmp(ixccmin1:ixccmax1,ixccmin2:ixccmax2,ixslice,1:ndim)
981  ps_sub(jgrid)%wC(ixcmin1:ixcmax1,ixcmin2:ixcmax2,1:nw+nwexpand) = &
982  wc_tmp(ixcmin1:ixcmax1,ixcmin2:ixcmax2,ixslice,1:nw+nwexpand)
983  ps_sub(jgrid)%xC(ixcmin1:ixcmax1,ixcmin2:ixcmax2,1:ndim) = &
984  xc_tmp(ixcmin1:ixcmax1,ixcmin2:ixcmax2,ixslice,1:ndim)
985  case default
986  call mpistop("slice direction not clear in fill_subnode")
987  end select
988  }
989 
990  end subroutine fill_subnode
991 
992  subroutine alloc_subnode(jgrid,dir,nwexpand)
994  integer, intent(in) :: jgrid, dir, nwexpand
995 
996  ! take care, what comes out is not necessarily a right handed system!
997  {^ifoned
998  allocate(ps_sub(jgrid)%w(1:nw+nwexpand),ps_sub(jgrid)%x(1:ndim))
999  allocate(ps_sub(jgrid)%wC(1:nw+nwexpand),ps_sub(jgrid)%xC(1:ndim))
1000  }
1001  {^iftwod
1002  select case (dir)
1003  case (1)
1004  allocate(ps_sub(jgrid)%w(ixglo2:ixghi2,1:nw+nwexpand),&
1005  ps_sub(jgrid)%x(ixglo2:ixghi2,1:ndim), &
1006  ps_sub(jgrid)%wC(ixglo2:ixghi2,1:nw+nwexpand),&
1007  ps_sub(jgrid)%xC(ixglo2:ixghi2,1:ndim))
1008  case (2)
1009  allocate(ps_sub(jgrid)%w(ixglo1:ixghi1,1:nw+nwexpand),&
1010  ps_sub(jgrid)%x(ixglo1:ixghi1,1:ndim), &
1011  ps_sub(jgrid)%wC(ixglo1:ixghi1,1:nw+nwexpand),&
1012  ps_sub(jgrid)%xC(ixglo1:ixghi1,1:ndim))
1013  case default
1014  call mpistop("slice direction not clear in alloc_subnode")
1015  end select
1016  }
1017  {^ifthreed
1018  select case (dir)
1019  case (1)
1020  allocate(ps_sub(jgrid)%w(ixglo2:ixghi2,ixglo3:ixghi3,1:nw+nwexpand),&
1021  ps_sub(jgrid)%x(ixglo2:ixghi2,ixglo3:ixghi3,1:ndim), &
1022  ps_sub(jgrid)%wC(ixglo2:ixghi2,ixglo3:ixghi3,1:nw+nwexpand),&
1023  ps_sub(jgrid)%xC(ixglo2:ixghi2,ixglo3:ixghi3,1:ndim))
1024  case (2)
1025  allocate(ps_sub(jgrid)%w(ixglo1:ixghi1,ixglo3:ixghi3,1:nw+nwexpand),&
1026  ps_sub(jgrid)%x(ixglo1:ixghi1,ixglo3:ixghi3,1:ndim), &
1027  ps_sub(jgrid)%wC(ixglo1:ixghi1,ixglo3:ixghi3,1:nw+nwexpand),&
1028  ps_sub(jgrid)%xC(ixglo1:ixghi1,ixglo3:ixghi3,1:ndim))
1029  case (3)
1030  allocate(ps_sub(jgrid)%w(ixglo1:ixghi1,ixglo2:ixghi2,1:nw+nwexpand),&
1031  ps_sub(jgrid)%x(ixglo1:ixghi1,ixglo2:ixghi2,1:ndim), &
1032  ps_sub(jgrid)%wC(ixglo1:ixghi1,ixglo2:ixghi2,1:nw+nwexpand),&
1033  ps_sub(jgrid)%xC(ixglo1:ixghi1,ixglo2:ixghi2,1:ndim))
1034  case default
1035  call mpistop("slice direction not clear in alloc_subnode")
1036  end select
1037  }
1038  end subroutine alloc_subnode
1039 
1040  subroutine dealloc_subnode(jgrid)
1042  integer, intent(in) :: jgrid
1043 
1044  if (jgrid==0) then
1045  call mpistop("trying to delete a non-existing grid in dealloc_subnode")
1046  end if
1047 
1048  deallocate(ps_sub(jgrid)%w,ps_sub(jgrid)%x,ps_sub(jgrid)%wC,ps_sub(jgrid)%xC)
1049 
1050  ! reset the global node info:
1051  node_sub(:,jgrid)=0
1052  rnode_sub(:,jgrid)=zero
1053 
1054  end subroutine dealloc_subnode
1055 
1056  subroutine fill_subnode_info(igrid,jgrid,dir)
1058  integer, intent(in) :: igrid,jgrid,dir
1059 
1060  node_sub(plevel_,jgrid)=node(plevel_,igrid)
1061  {^ifthreed
1062  select case(dir)
1063  case (1)
1064  node_sub(pig1_,jgrid)=node(pig2_,igrid)
1065  node_sub(pig2_,jgrid)=node(pig3_,igrid)
1066  rnode_sub(rpdx1_,jgrid)=rnode(rpdx2_,igrid)
1067  rnode_sub(rpdx2_,jgrid)=rnode(rpdx3_,igrid)
1068  rnode_sub(rpxmin1_,jgrid)=rnode(rpxmin2_,igrid)
1069  rnode_sub(rpxmin2_,jgrid)=rnode(rpxmin3_,igrid)
1070  case (2)
1071  node_sub(pig1_,jgrid)=node(pig1_,igrid)
1072  node_sub(pig2_,jgrid)=node(pig3_,igrid)
1073  rnode_sub(rpdx1_,jgrid)=rnode(rpdx1_,igrid)
1074  rnode_sub(rpdx2_,jgrid)=rnode(rpdx3_,igrid)
1075  rnode_sub(rpxmin1_,jgrid)=rnode(rpxmin1_,igrid)
1076  rnode_sub(rpxmin2_,jgrid)=rnode(rpxmin3_,igrid)
1077  case (3)
1078  node_sub(pig1_,jgrid)=node(pig1_,igrid)
1079  node_sub(pig2_,jgrid)=node(pig2_,igrid)
1080  rnode_sub(rpdx1_,jgrid)=rnode(rpdx1_,igrid)
1081  rnode_sub(rpdx2_,jgrid)=rnode(rpdx2_,igrid)
1082  rnode_sub(rpxmin1_,jgrid)=rnode(rpxmin1_,igrid)
1083  rnode_sub(rpxmin2_,jgrid)=rnode(rpxmin2_,igrid)
1084  case default
1085  call mpistop("slice direction not clear in fill_subnode_info")
1086  end select
1087  }
1088  {^iftwod
1089  select case(dir)
1090  case (1)
1091  node_sub(pig1_,jgrid)=node(pig2_,igrid)
1092  rnode_sub(rpdx1_,jgrid)=rnode(rpdx2_,igrid)
1093  rnode_sub(rpxmin1_,jgrid)=rnode(rpxmin2_,igrid)
1094  case (2)
1095  node_sub(pig1_,jgrid)=node(pig1_,igrid)
1096  rnode_sub(rpdx1_,jgrid)=rnode(rpdx1_,igrid)
1097  rnode_sub(rpxmin1_,jgrid)=rnode(rpxmin1_,igrid)
1098  case default
1099  call mpistop("slice direction not clear in fill_subnode_info")
1100  end select
1101  }
1102  {^ifoned
1103  }
1104 
1105  end subroutine fill_subnode_info
1106 
1107  subroutine get_igslice(dir,x,igslice)
1109  integer, intent(in) :: dir
1110  double precision, intent(in) :: x
1111  integer, dimension(nlevelshi), intent(out) :: igslice
1112  ! .. local ..
1113  integer level
1114 
1115  if (x.ne.x) &
1116  call mpistop("get_igslice: your slice position is NaN!")
1117 
1118  select case (dir)
1119  {case (^d)
1120  do level = 1, refine_max_level
1121  igslice(level) = int((x-xprobmin^d)/dg^d(level))+1
1122  ! Gets out of domain when x==xprobmax^D, not caught by put_slice, so limit:
1123  if (x>=xprobmax^d) igslice(level) = int((xprobmax^d-xprobmin^d)/dg^d(level))
1124  ! This is already caught by control in put_slice, but anyways:
1125  if (x<=xprobmin^d) igslice(level) = 1
1126  end do\}
1127  case default
1128  call mpistop("slice direction not clear in get_igslice")
1129  end select
1130  end subroutine get_igslice
1131 
1132  double precision function roundoff_minmax(val,minval,maxval)
1133  implicit none
1134  double precision,intent(in) :: val, minval, maxval
1135 
1136  roundoff_minmax = val
1137 
1138  if (abs(roundoff_minmax) .le. minval) then
1139  roundoff_minmax = 0.0d0
1140  end if
1141 
1142  if (roundoff_minmax .gt. maxval) then
1143  roundoff_minmax = maxval
1144  else if (roundoff_minmax .lt. -maxval) then
1145  roundoff_minmax = -maxval
1146  end if
1147 
1148  ! Replace NaN with maxval (e.g. Paraview chokes on ASCII NaN):
1150 
1151  end function roundoff_minmax
1152 
1153 end module mod_slice
1154 
subroutine mpistop(message)
Exit MPI-AMRVAC with an error message.
Definition: comm_lib.t:194
recursive subroutine traverse_slice(tree)
Definition: mod_slice.t:779
subroutine put_slice_dat
Definition: mod_slice.t:558
subroutine put_slice_zerod
Definition: mod_slice.t:623
subroutine put_slice_csv
Definition: mod_slice.t:430
subroutine write_slice_vtk(jgrid, slice_fh, wnamei)
Definition: mod_slice.t:303
subroutine put_slice_line(jout, file_handle)
Definition: mod_slice.t:499
subroutine put_slice_vtu
Definition: mod_slice.t:215
Module with basic data types used in amrvac.
integer, parameter size_double
Size (in bytes) of double precision real.
Handles computations for coordinates and variables in output.
subroutine calc_grid(qunit, igrid, xC, xCC, xC_TMP, xCC_TMP, wC_TMP, wCC_TMP, normconv, ixCL, ixCCL, first)
Compute both corner as well as cell-centered values for output.
subroutine getheadernames(wnamei, xandwnamei, outfilehead)
get all variables names
subroutine calc_x(igrid, xC, xCC)
computes cell corner (xC) and cell center (xCC) coordinates
Module with basic grid data structures.
Definition: mod_forest.t:2
integer, dimension(:), allocatable, save morton_sub_start
Definition: mod_forest.t:67
integer, dimension(:), allocatable, save morton_sub_stop
Definition: mod_forest.t:67
type(tree_node_ptr), dimension(:^d &), allocatable, save tree_root
Pointers to the coarse grid.
Definition: mod_forest.t:29
This module contains definitions of global parameters and variables and some generic functions/subrou...
integer, parameter unitslice
double precision global_time
The global simulation time.
double precision, dimension(:), allocatable dg
extent of grid blocks in domain per dimension, in array over levels
integer, parameter ndim
Number of spatial dimensions for grid variables.
double precision time_convert_factor
Conversion factor for time unit.
integer icomm
The MPI communicator.
integer mype
The rank of the current MPI task.
integer, parameter plevel_
integer, dimension(:), allocatable, parameter d
integer ixm
the mesh range (within a block with ghost cells)
integer ierrmpi
A global MPI error return code.
integer npe
The number of MPI tasks.
integer nwauxio
Number of auxiliary variables that are only included in output.
logical, dimension(:), allocatable w_write
double precision, dimension(:,:), allocatable rnode
Corner coordinates.
integer nghostcells
Number of ghost cells surrounding a grid.
double precision, dimension(:,:), allocatable rnode_sub
character(len=std_len) base_filename
Base file name for simulation output, which will be followed by a number.
integer refine_max_level
Maximal number of AMR levels.
integer, dimension(:,:), allocatable node
integer, dimension(:,:), allocatable node_sub
Writes D-1 slice, can do so in various formats, depending on slice_type.
Definition: mod_slice.t:2
subroutine fill_subnode_info(igrid, jgrid, dir)
Definition: mod_slice.t:1057
subroutine get_igslice(dir, x, igslice)
Definition: mod_slice.t:1108
subroutine fill_subnode(igrid, active, jgrid, dir, xslice, normconv)
Definition: mod_slice.t:855
integer slicenext
the file number of slices
Definition: mod_slice.t:13
integer, parameter nslicemax
Maximum number of slices.
Definition: mod_slice.t:7
subroutine alloc_subnode(jgrid, dir, nwexpand)
Definition: mod_slice.t:993
subroutine put_slice(dir, xslice)
Definition: mod_slice.t:44
subroutine write_slice
Definition: mod_slice.t:30
character(len=std_len) slice_type
choose data type of slice: vtu, vtuCC, dat, or csv
Definition: mod_slice.t:22
integer, dimension(nslicemax) slicedir
The slice direction for each slice.
Definition: mod_slice.t:19
double precision function roundoff_minmax(val, minval, maxval)
Definition: mod_slice.t:1133
integer nslices
Number of slices to output.
Definition: mod_slice.t:16
double precision, dimension(nslicemax) slicecoord
Slice coordinates, see Slice output.
Definition: mod_slice.t:10
subroutine select_slice(dir, xslice, writeonly, file_handle, normconv)
Definition: mod_slice.t:693
subroutine dealloc_subnode(jgrid)
Definition: mod_slice.t:1041
Pointer to a tree_node.
Definition: mod_forest.t:6