MPI-AMRVAC 3.1
The MPI - Adaptive Mesh Refinement - Versatile Advection Code (development version)
Loading...
Searching...
No Matches
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
4 use mod_comm_lib, only: mpistop
5 implicit none
6
7 !> Slice coordinates, see @ref slices.md
8 double precision :: slicecoord(1000)
9
10 !> Maximum number of slices
11 integer, parameter :: nslicemax=1000
12
13 !> Number of slices to output
14 integer :: nslices
15
16 !> The slice direction for each slice
17 integer :: slicedir(nslicemax)
18
19 !> tag for MPI message
20 integer, private :: itag
21
22 !> choose data type of slice: vtu, vtuCC, dat, or csv
23 character(len=std_len) :: slice_type
24
25 !> igrid and ipe of the leaf blocks on the slice, with the first (1,1:2) stores the nleafs and nparents
26 integer, allocatable :: sfc_sub(:,:)
27
28contains
29
30 subroutine write_slice
32 ! Writes a D-1 slice
33 ! by Oliver Porth
34 ! 22.Nov 2011
35 integer :: islice
36
37 do islice=1,nslices
38 call put_slice(slicedir(islice),slicecoord(islice))
39 end do
40
42 end subroutine write_slice
43
44 subroutine put_slice(dir,xslice)
47 ! Writes a D-1 slice
48 ! For ONED simulations, output will be appended to one csv-file per slice
49 ! slices are sensitive to the saveprim switch and
50 ! can contain auxiliary io variables (nwauxio)
51 ! Thus csv or vtu(CC)-files with primitive variables are obtained.
52 ! when slice_type='dat', we save a D-1 dat file for potential restarts
53 ! by Oliver Porth
54 ! 22.Nov 2011
55 integer, intent(in) :: dir
56 double precision, intent(in) :: xslice
57 ! .. local ..
58 double precision,dimension(0:nw+nwauxio) :: normconv
59 integer :: Njgrid, jgrid
60 integer, dimension(ndim-1) :: ixsubGlo, ixsubGhi
61 integer, dimension(ndim-1) :: ixsubMlo, ixsubMhi
62 integer :: size_subblock_io, nx^D, slice_fh, nwexpand
63 integer :: type_subblock_io, type_subblockC_io, type_subblock_x_io, type_subblockC_x_io
64 integer, dimension(ndim) :: sizes, subsizes, start
65
66 double precision, dimension(ndim-1) :: xprobminsub, xprobmaxsub
67 integer, dimension(ndim-1) :: block_nx_sub, domain_nx_sub
68 logical, dimension(ndim-1) :: periodBsub
69
70 ! Preamble:
71 nx^d=ixmhi^d-ixmlo^d+1;
72 slice_fh=unitslice
73
74 if (ndim==1) then
75 nwexpand = nwauxio
76 else
77 if(slice_type/='dat')then
78 nwexpand = nwauxio
79 else
80 nwexpand = 0
81 endif
82 end if
83
84 ! Do a last consistency check:
85 select case(dir)
86 {case(^d)
87 if(xslice<xprobmin^d.or.xslice>xprobmax^d) &
88 call mpistop("slice out of bounds")
89 \}
90 end select
91
92 ! Traverse the forest and fill nodes:
93 call select_slice(dir,xslice,.false.,slice_fh,normconv)
94
95 ! Create the MPI-datatype and select indices:
96 {^ifthreed
97 select case(dir)
98 case (1)
99 sizes(1) = ixghi2; sizes(2) = ixghi3;
100 ixsubglo(1) = ixglo2; ixsubglo(2) = ixglo3;
101 ixsubghi(1) = ixghi2; ixsubghi(2) = ixghi3;
102 subsizes(1)=nx2;subsizes(2)=nx3;
103 start(1)=ixmlo2-1;start(2)=ixmlo3-1;
104 size_subblock_io=nx2*nx3*(nw+nwexpand)*size_double
105 xprobminsub(1)=xprobmin2; xprobmaxsub(1)=xprobmax2;
106 xprobminsub(2)=xprobmin3; xprobmaxsub(2)=xprobmax3;
107 domain_nx_sub(1)=domain_nx2; domain_nx_sub(2)=domain_nx3;
108 block_nx_sub(1)=block_nx2; block_nx_sub(2)=block_nx3;
109 periodbsub(1)=periodb(2); periodbsub(2)=periodb(3);
110 case (2)
111 sizes(1) = ixghi1; sizes(2) = ixghi3;
112 ixsubglo(1) = ixglo1; ixsubglo(2) = ixglo3;
113 ixsubghi(1) = ixghi1; ixsubghi(2) = ixghi3;
114 subsizes(1)=nx1;subsizes(2)=nx3;
115 start(1)=ixmlo1-1;start(2)=ixmlo3-1;
116 size_subblock_io=nx1*nx3*(nw+nwexpand)*size_double
117 xprobminsub(1)=xprobmin1; xprobmaxsub(1)=xprobmax1;
118 xprobminsub(2)=xprobmin3; xprobmaxsub(2)=xprobmax3;
119 domain_nx_sub(1)=domain_nx1; domain_nx_sub(2)=domain_nx3;
120 block_nx_sub(1)=block_nx1; block_nx_sub(2)=block_nx3;
121 periodbsub(1)=periodb(1); periodbsub(2)=periodb(3);
122 case (3)
123 ixsubglo(1) = ixglo1; ixsubglo(2) = ixglo2;
124 ixsubghi(1) = ixghi1; ixsubghi(2) = ixghi2;
125 sizes(1) = ixghi1; sizes(2) = ixghi2;
126 subsizes(1)=nx1;subsizes(2)=nx2;
127 start(1)=ixmlo1-1;start(2)=ixmlo2-1;
128 size_subblock_io=nx1*nx2*(nw+nwexpand)*size_double
129 xprobminsub(1)=xprobmin1; xprobmaxsub(1)=xprobmax1;
130 xprobminsub(2)=xprobmin2; xprobmaxsub(2)=xprobmax2;
131 domain_nx_sub(1)=domain_nx1; domain_nx_sub(2)=domain_nx2;
132 block_nx_sub(1)=block_nx1; block_nx_sub(2)=block_nx2;
133 periodbsub(1)=periodb(1); periodbsub(2)=periodb(2);
134 case default
135 call mpistop("slice direction not clear in put_slice")
136 end select
137 }
138 {^iftwod
139 select case(dir)
140 case (1)
141 ixsubglo(1) = ixglo2; ixsubghi(1) = ixghi2;
142 sizes(1) = ixghi2
143 subsizes(1)=nx2
144 start(1)=ixmlo2-1
145 size_subblock_io=nx2*(nw+nwexpand)*size_double
146 xprobminsub(1)=xprobmin2; xprobmaxsub(1)=xprobmax2;
147 domain_nx_sub(1)=domain_nx2;
148 block_nx_sub(1)=block_nx2;
149 periodbsub(1)=periodb(2);
150 case (2)
151 ixsubglo(1) = ixglo1; ixsubghi(1) = ixghi1;
152 sizes(1) = ixghi1
153 subsizes(1)=nx1
154 start(1)=ixmlo1-1
155 size_subblock_io=nx1*(nw+nwexpand)*size_double
156 xprobminsub(1)=xprobmin1; xprobmaxsub(1)=xprobmax1;
157 domain_nx_sub(1)=domain_nx1;
158 block_nx_sub(1)=block_nx1;
159 periodbsub(1)=periodb(1);
160 case default
161 call mpistop("slice direction not clear in put_slice")
162 end select
163 }
164 {^ifoned
165 size_subblock_io=(nw+nwexpand)*size_double
166 }
167
168 {^nooned
169 {^de&ixsubmlo(^de-1) = ixsubglo(^de-1)+nghostcells;}
170 {^de&ixsubmhi(^de-1) = ixsubghi(^de-1)-nghostcells;}
171 }
172
173 sizes(ndim)=nw+nwexpand
174 subsizes(ndim)=nw+nwexpand
175 start(ndim)=0
176
177 ! Types for center variables:
178 call mpi_type_create_subarray(ndim,sizes,subsizes,start, &
179 mpi_order_fortran,mpi_double_precision, &
180 type_subblock_io,ierrmpi)
181 call mpi_type_commit(type_subblock_io,ierrmpi)
182
183 sizes(ndim)=^nd
184 subsizes(ndim)=^nd
185 start(ndim)=0
186 call mpi_type_create_subarray(ndim,sizes,subsizes,start, &
187 mpi_order_fortran,mpi_double_precision, &
188 type_subblock_x_io,ierrmpi)
189 call mpi_type_commit(type_subblock_x_io,ierrmpi)
190
191
192 ! Types for corner variables:
193 subsizes(1:ndim-1) = subsizes(1:ndim-1) + 1
194 start(1:ndim-1) = start(1:ndim-1) - 1
195 sizes(ndim)=nw+nwexpand
196 subsizes(ndim)=nw+nwexpand
197 start(ndim)=0
198
199 call mpi_type_create_subarray(ndim,sizes,subsizes,start, &
200 mpi_order_fortran,mpi_double_precision, &
201 type_subblockc_io,ierrmpi)
202 call mpi_type_commit(type_subblockc_io,ierrmpi)
203
204 sizes(ndim)=^nd
205 subsizes(ndim)=^nd
206 start(ndim)=0
207 call mpi_type_create_subarray(ndim,sizes,subsizes,start, &
208 mpi_order_fortran,mpi_double_precision, &
209 type_subblockc_x_io,ierrmpi)
210 call mpi_type_commit(type_subblockc_x_io,ierrmpi)
211
212
213 ! local number of sub-grids:
215
216 ! Now output using various schemes:
217 if (ndim==1) then
218 call put_slice_zerod
219 else
220 select case(slice_type)
221 case ('csv')
222 call put_slice_csv
223 case ('dat')
224 call put_slice_dat
225 case ('vtu', 'vtuCC')
226 call put_slice_vtu
227 end select
228 end if
229
230 do jgrid=1,njgrid
231 call dealloc_subnode(jgrid)
232 end do
233
234 call mpi_type_free(type_subblock_io,ierrmpi)
235 call mpi_type_free(type_subblock_x_io,ierrmpi)
236 call mpi_type_free(type_subblockc_io,ierrmpi)
237 call mpi_type_free(type_subblockc_x_io,ierrmpi)
238
239 if (allocated(sfc_sub)) deallocate(sfc_sub)
240
241 contains
242
243 subroutine put_slice_vtu
244
246 integer :: status(MPI_STATUS_SIZE), ipe
247 logical :: fileopen
248 character(len=1024) :: filename, xlabel
249 character(len=79) :: xxlabel
250 character(len=name_len) :: wnamei(1:nw+nwauxio),xandwnamei(1:ndim+nw+nwauxio)
251 character(len=1024) :: outfilehead
252
253 if (mype==0) then
254
255 inquire(slice_fh,opened=fileopen)
256 if(.not.fileopen)then
257 ! generate filename:
258 write(xlabel,"(D9.2)")xslice
259 xxlabel=trim(xlabel)
260 if(xslice>=zero)then
261 write(xxlabel(1:1),"(a)") "+"
262 endif
263 write(filename,"(a,i1.1,a,i4.4,a)") trim(base_filename)//'_d',dir,'_x'//trim(xxlabel)//'_n',slicenext,'.vtu'
264 open(slice_fh,file=filename,status='unknown',form='formatted')
265 end if
266 ! get and write the header:
267 call getheadernames(wnamei,xandwnamei,outfilehead)
268 ! generate xml header
269 write(slice_fh,'(a)')'<?xml version="1.0"?>'
270 write(slice_fh,'(a)',advance='no') '<VTKFile type="UnstructuredGrid"'
271 if(type_endian==1)then
272 write(slice_fh,'(a)')' version="0.1" byte_order="LittleEndian">'
273 else
274 write(slice_fh,'(a)')' version="0.1" byte_order="BigEndian">'
275 endif
276 write(slice_fh,'(a)')' <UnstructuredGrid>'
277 write(slice_fh,'(a)')'<FieldData>'
278 write(slice_fh,'(2a)')'<DataArray type="Float32" Name="TIME" ',&
279 'NumberOfTuples="1" format="ascii">'
280 write(slice_fh,*) real(global_time*time_convert_factor)
281 write(slice_fh,'(a)')'</DataArray>'
282 write(slice_fh,'(a)')'</FieldData>'
283
284 ! write to file:
285 do jgrid=1, njgrid
286 call write_slice_vtk(jgrid,slice_fh,wnamei)
287 end do
288
289 ! create a recv buffer using allocate, will be deallocated at the end of the routine:
290 call alloc_subnode(njgrid+1,dir,nwauxio)
291
292 end if
293
294 ! Also communicate the normconv array since processor zero might not have it yet:
295 if (npe>1) then
296 do ipe=1,npe-1
297 do jgrid=1,morton_sub_stop(ipe)-morton_sub_start(ipe)+1
298 itag=morton_sub_start(ipe)+jgrid-1
299 itag=itag*5
300 if (ipe == mype ) then
301 call mpi_send(ps_sub(jgrid)%x,1,type_subblock_x_io,0,itag,icomm,ierrmpi)
302 call mpi_send(ps_sub(jgrid)%w,1,type_subblock_io,0,itag+1,icomm,ierrmpi)
303 call mpi_send(ps_sub(jgrid)%xC,1,type_subblockc_x_io,0,itag+2,icomm,ierrmpi)
304 call mpi_send(ps_sub(jgrid)%wC,1,type_subblockc_io,0,itag+3,icomm,ierrmpi)
305 call mpi_send(normconv,nw+nwauxio+1,mpi_double_precision,0,itag+4,icomm,ierrmpi)
306 end if
307 if (mype == 0) then
308 call mpi_recv(ps_sub(njgrid+1)%x,1,type_subblock_x_io,ipe,itag,icomm,status,ierrmpi)
309 call mpi_recv(ps_sub(njgrid+1)%w,1,type_subblock_io,ipe,itag+1,icomm,status,ierrmpi)
310 call mpi_recv(ps_sub(njgrid+1)%xC,1,type_subblockc_x_io,ipe,itag+2,icomm,status,ierrmpi)
311 call mpi_recv(ps_sub(njgrid+1)%wC,1,type_subblockc_io,ipe,itag+3,icomm,status,ierrmpi)
312 call mpi_recv(normconv,nw+nwauxio+1,mpi_double_precision,ipe,&
313 itag+4,icomm,status,ierrmpi)
314 call write_slice_vtk(njgrid+1,slice_fh,wnamei)
315 end if
316 end do
317 end do
318 endif
319
320 if (mype==0) then
321
322 write(slice_fh,'(a)')'</UnstructuredGrid>'
323 write(slice_fh,'(a)')'</VTKFile>'
324 close(slice_fh)
325 call dealloc_subnode(njgrid+1)
326
327 end if
328
329 end subroutine put_slice_vtu
330
331 subroutine write_slice_vtk(jgrid,slice_fh,wnamei)
332
333 ! this only works for 2D and 3D, 1D reduction (line to point) not allowed
334 integer, intent(in) :: jgrid, slice_fh
335 character(len=name_len), intent(in) :: wnamei(1:nw+nwauxio)
336 ! This remainder part only for more than 1D, but nesting with NOONED gives problems
337 {#IFNDEF D1
338 ! .. local ..
339 double precision :: x_VTK(1:3)
340 double precision, parameter :: minvalue = 1.0d-99, maxvalue = 1.0d+99
341 integer :: ixC^L, ixCC^L, nc, np, iw
342 integer :: nx^DM, nxC^DM, icell, ix^DM
343 integer :: VTK_type
344
345 {^dm&ixccmin^dm = ixsubmlo(^dm);}
346 {^dm&ixccmax^dm = ixsubmhi(^dm);}
347 {^dm&ixcmin^dm = ixsubmlo(^dm)-1;}
348 {^dm&ixcmax^dm = ixsubmhi(^dm);}
349
350 nx^dm=ixccmax^dm-ixccmin^dm+1;
351 nxc^dm=nx^dm+1;
352 nc={nx^dm*} ! Number of cells per subgrid
353 np={nxc^dm*} ! Number of corner points per subgrid
354
355 ! we write out every grid as one VTK PIECE
356 write(slice_fh,'(a,i7,a,i7,a)') &
357 '<Piece NumberOfPoints="',np,'" NumberOfCells="',nc,'">'
358
359 !==============================
360 ! celldata or pointdata?
361 !==============================
362 select case(slice_type)
363
364 case('vtuCC') ! celldata
365 write(slice_fh,'(a)')'<CellData>'
366 do iw=1,nw+nwauxio
367 if(iw<=nw) then
368 if(.not.w_write(iw)) cycle
369 endif
370 write(slice_fh,'(a,a,a)')&
371 '<DataArray type="Float64" Name="',trim(wnamei(iw)),'" format="ascii">'
372 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)}
373 write(slice_fh,'(a)')'</DataArray>'
374 enddo
375 write(slice_fh,'(a)')'</CellData>'
376
377
378 case('vtu') ! pointdata
379 write(slice_fh,'(a)')'<PointData>'
380 do iw=1,nw+nwauxio
381 if(iw<=nw) then
382 if(.not.w_write(iw)) cycle
383 endif
384 write(slice_fh,'(a,a,a)')&
385 '<DataArray type="Float64" Name="',trim(wnamei(iw)),'" format="ascii">'
386 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)}
387 write(slice_fh,'(a)')'</DataArray>'
388 enddo
389 write(slice_fh,'(a)')'</PointData>'
390
391
392 end select
393 !==============================
394 ! Done: celldata or pointdata?
395 !==============================
396
397 !==============================
398 ! output Cornerpoints
399 !==============================
400 write(slice_fh,'(a)')'<Points>'
401 write(slice_fh,'(a)')'<DataArray type="Float32" NumberOfComponents="3" format="ascii">'
402 ! write cell corner coordinates in a backward dimensional loop, always 3D output
403 {do ix^dmb=ixcmin^dmb,ixcmax^dmb \}
404 x_vtk(1:3)=zero;
405 x_vtk(1:ndim)=ps_sub(jgrid)%xC(ix^dm,1:ndim)*normconv(0);
406 write(slice_fh,'(3(1pe14.6))') x_vtk
407 {^dm&end do \}
408 write(slice_fh,'(a)')'</DataArray>'
409 write(slice_fh,'(a)')'</Points>'
410 !==============================
411 ! Done: output Cornerpoints
412 !==============================
413
414 !==============================
415 ! cell Metainformation
416 !==============================
417 write(slice_fh,'(a)')'<Cells>'
418
419 ! connectivity part
420 write(slice_fh,'(a)')'<DataArray type="Int32" Name="connectivity" format="ascii">'
421
422 {^dm& do ix^dmb=1,nx^dmb\}
423 {^iftwod
424 write(slice_fh,'(2(i7,1x))')ix1-1,ix1
425 }{^ifthreed
426 write(slice_fh,'(4(i7,1x))')(ix2-1)*nxc1+ix1-1, &
427 (ix2-1)*nxc1+ix1,ix2*nxc1+ix1-1,ix2*nxc1+ix1
428 }{^dm& end do\}
429
430 write(slice_fh,'(a)')'</DataArray>'
431
432 ! offsets data array
433 write(slice_fh,'(a)')'<DataArray type="Int32" Name="offsets" format="ascii">'
434 do icell=1,nc
435 write(slice_fh,'(i7)') icell*(2**(^nd-1))
436 end do
437 write(slice_fh,'(a)')'</DataArray>'
438
439 ! VTK cell type data array
440 write(slice_fh,'(a)')'<DataArray type="Int32" Name="types" format="ascii">'
441 ! VTK_LINE=3; VTK_PIXEL=8; VTK_VOXEL=11 -> vtk-syntax
442 {^iftwod vtk_type=3 \}
443 {^ifthreed vtk_type=8 \}
444 do icell=1,nc
445 write(slice_fh,'(i2)') vtk_type
446 enddo
447 write(slice_fh,'(a)')'</DataArray>'
448
449 write(slice_fh,'(a)')'</Cells>'
450 !==============================
451 ! Done: cell Metainformation
452 !==============================
453 write(slice_fh,'(a)')'</Piece>'
454
455 }
456 end subroutine write_slice_vtk
457
458 subroutine put_slice_csv
459
461 integer :: iw, ipe, itag
462 integer :: status(MPI_STATUS_SIZE)
463 logical :: fileopen
464 character(len=1024) :: filename, xlabel
465 character(len=79) :: xxlabel
466 character(len=name_len) :: wnamei(1:nw+nwauxio),xandwnamei(1:ndim+nw+nwauxio)
467 character(len=1024) :: outfilehead
468 character(len=1024) :: line
469
470 if (mype==0) then
471 inquire(slice_fh,opened=fileopen)
472 if(.not.fileopen)then
473 ! generate filename:
474 write(xlabel,"(D9.2)")xslice
475 xxlabel=trim(xlabel)
476 if(xslice>=zero)then
477 write(xxlabel(1:1),"(a)") "+"
478 endif
479 write(filename,"(a,i1.1,a,i4.4,a)") trim(base_filename)//'_d',dir,'_x'//trim(xxlabel)//'_n',slicenext,'.csv'
480 open(slice_fh,file=filename,status='unknown',form='formatted')
481 end if
482 ! get and write the header:
483 call getheadernames(wnamei,xandwnamei,outfilehead)
484 line=''
485 do iw=1,ndim+nw+nwauxio-1
486 line = trim(line)//trim(xandwnamei(iw))//', '
487 end do
488 line = trim(line)//trim(xandwnamei(ndim+nw+nwauxio))
489 write(slice_fh,'(a)')trim(line)
490 ! create a recv buffer using allocate, will be deallocated at the end of the routine:
491 call alloc_subnode(njgrid+1,dir,nwauxio)
492
493 ! write to file:
494 do jgrid=1, njgrid
495 call put_slice_line(jgrid,slice_fh)
496 end do
497 end if
498
499 ! Also communicate the normconv array since processor zero might not have it yet:
500 if (npe>1) then
501 do ipe=1,npe-1
502 do jgrid=1,morton_sub_stop(ipe)-morton_sub_start(ipe)+1
503 itag=morton_sub_start(ipe)+jgrid-1
504 if (ipe == mype ) then
505 call mpi_send(ps_sub(jgrid)%x,1,type_subblock_x_io,0,itag,icomm,ierrmpi)
506 call mpi_send(ps_sub(jgrid)%w,1,type_subblock_io,0,itag,icomm,ierrmpi)
507 call mpi_send(normconv,nw+nwauxio+1,mpi_double_precision,0,itag,icomm,ierrmpi)
508 end if
509 if (mype == 0) then
510 call mpi_recv(ps_sub(njgrid+1)%x,1,type_subblock_x_io,ipe,itag,icomm,status,ierrmpi)
511 call mpi_recv(ps_sub(njgrid+1)%w,1,type_subblock_io,ipe,itag,icomm,status,ierrmpi)
512 call mpi_recv(normconv,nw+nwauxio+1,mpi_double_precision,ipe,&
513 itag,icomm,status,ierrmpi)
514 call put_slice_line(njgrid+1,slice_fh)
515 end if
516 end do
517 end do
518 endif
519
520 if (mype==0) then
521 close(slice_fh)
522 call dealloc_subnode(njgrid+1)
523 end if
524
525 end subroutine put_slice_csv
526
527 subroutine put_slice_line(jout,file_handle)
528 integer, intent(in) :: jout, file_handle
529 ! .. local ..
530 double precision, parameter :: minvalue = 1.0d-99, maxvalue = 1.0d+99
531 integer :: ix^D,idir,iw
532 character(len=1024) ::line, data
533
534 {^ifthreed
535 do ix2=ixsubmlo(2),ixsubmhi(2)
536 do ix1=ixsubmlo(1),ixsubmhi(1)
537 \}
538 {^iftwod
539 do ix1=ixsubmlo(1),ixsubmhi(1)
540 }
541 ! Format the line:
542 line = ''
543 do idir=1,ndim
544 {^ifthreed
545 write(data,"(es14.6)")roundoff_minmax(ps_sub(jout)%x(ix1,ix2,idir),minvalue,maxvalue)
546 }
547 {^iftwod
548 write(data,"(es14.6)")roundoff_minmax(ps_sub(jout)%x(ix1,idir),minvalue,maxvalue)
549 }
550 {^ifoned
551 write(data,"(es14.6)")roundoff_minmax(ps_sub(jout)%x(idir),minvalue,maxvalue)
552 }
553
554 line = trim(line)//trim(data)//', '
555 end do
556 do iw = 1,nw+nwauxio-1
557 {^ifthreed
558 write(data,"(es14.6)")roundoff_minmax(ps_sub(jout)%w(ix1,ix2,iw)*normconv(iw),minvalue,maxvalue)
559 }
560 {^iftwod
561 write(data,"(es14.6)")roundoff_minmax(ps_sub(jout)%w(ix1,iw)*normconv(iw),minvalue,maxvalue)
562 }
563 {^ifoned
564 write(data,"(es14.6)")roundoff_minmax(ps_sub(jout)%w(iw)*normconv(iw),minvalue,maxvalue)
565 }
566 line = trim(line)//trim(data)//', '
567 end do
568 {^ifthreed
569 write(data,"(es14.6)")roundoff_minmax(ps_sub(jout)%w(ix1,ix2,nw+nwauxio)*normconv(nw+nwauxio),minvalue,maxvalue)
570 }
571 {^iftwod
572 write(data,"(es14.6)")roundoff_minmax(ps_sub(jout)%w(ix1,nw+nwauxio)*normconv(nw+nwauxio),minvalue,maxvalue)
573 }
574 line = trim(line)//trim(data)
575 write(file_handle,'(a)')trim(line)
576 {^iftwod
577 end do
578 }
579 {^ifthreed
580 end do
581 end do
582 \}
583
584 end subroutine put_slice_line
585
586 subroutine put_slice_dat
587
591
592 integer, dimension(max_blocks*2) :: iorequest
593 integer, dimension(MPI_STATUS_SIZE,max_blocks) :: iostatus
594 integer(kind=MPI_OFFSET_KIND) :: offset, offset1
595 integer :: nsubleafs
596 integer :: amode, status(MPI_STATUS_SIZE), iwrite
597 character(len=1024) :: filename, xlabel
598 character(len=79) :: xxlabel
599 character(len=name_len) :: dname
600
601 integer :: nsubparents, iw, ileaf, inode
602 integer(kind=MPI_OFFSET_KIND) :: block_offset, tree_offset
603 integer, allocatable :: block_ig_sub(:,:)
604 integer, allocatable :: block_lvl_sub(:)
605 integer(kind=MPI_OFFSET_KIND), allocatable :: block_offset_sub(:)
606 logical, allocatable :: isleaf_sub(:)
607 type(tree_node), pointer :: pnode
608
609 nsubleafs=morton_sub_stop(npe-1)
610
611 ! only on processor 0, prepare all the data for writing header
612 if (mype==0) then
613
614 if (nsubleafs /= sfc_sub(0,1)) error stop "nsubleafs /= sfc_sub(0,1)"
615 nsubparents = sfc_sub(0,2)
616 allocate(block_ig_sub(ndim-1,nsubleafs))
617 allocate(block_lvl_sub(nsubleafs))
618 allocate(block_offset_sub(nsubleafs))
619 allocate(isleaf_sub(nsubleafs+nsubparents))
620
621 block_offset = int(0,kind=mpi_offset_kind)
622 ileaf = 0
623 ! get block_ig_sub, block_lvl_sub, block_offset_sub
624 do inode=1,nsubleafs+nsubparents
625 if (sfc_sub(inode,1) == 0) then
626 isleaf_sub(inode) = .false.
627 else
628 isleaf_sub(inode) = .true.
629 pnode => igrid_to_node(sfc_sub(inode,1),sfc_sub(inode,2))%node
630 ileaf = ileaf + 1
631 {^ifthreed
632 select case (dir)
633 case (1)
634 block_ig_sub(:,ileaf) = [ pnode%ig2, pnode%ig3 ]
635 case (2)
636 block_ig_sub(:,ileaf) = [ pnode%ig1, pnode%ig3 ]
637 case (3)
638 block_ig_sub(:,ileaf) = [ pnode%ig1, pnode%ig2 ]
639 end select
640 }
641 {^iftwod
642 select case (dir)
643 case (1)
644 block_ig_sub(:,ileaf) = [ pnode%ig2 ]
645 case (2)
646 block_ig_sub(:,ileaf) = [ pnode%ig1 ]
647 end select
648 }
649 block_lvl_sub(ileaf) = pnode%level
650 ! will not save ghost cells, start from offset_block_data
651 block_offset_sub(ileaf) = block_offset
652 block_offset = block_offset + int(2*(ndim-1)*size_int,kind=mpi_offset_kind)
653 block_offset = block_offset + int(product(block_nx_sub)*nw*size_double,kind=mpi_offset_kind)
654 end if
655 end do
656 end if
657
658 ! generate filename
659 write(xlabel,"(D9.2)")xslice
660 xxlabel=trim(xlabel)
661 if(xslice>=zero)then
662 write(xxlabel(1:1),"(a)") "+"
663 endif
664 write(filename,"(a,i1.1,a,i4.4,a)") trim(base_filename)//'_d',dir,'_x'//trim(xxlabel)//'_n',slicenext,'.dat'
665
666 if(mype==0) then
667 open(unit=slice_fh,file=filename,status='replace')
668 close(unit=slice_fh)
669 end if
670
671 ! header writing
672 if (mype==0) then
673 ! call select_slice(dir,xslice,.true.,slice_fh,normconv)
674 amode=ior(mpi_mode_append,mpi_mode_wronly)
675 call mpi_file_open(mpi_comm_self,filename,amode,mpi_info_null, &
676 slice_fh,ierrmpi)
677
678 call mpi_file_write(slice_fh, version_number, 1, mpi_integer, status, ierrmpi)
679 ! two offsets for tree and block data
680 call mpi_file_write(slice_fh, 0, 1, mpi_integer, status, ierrmpi)
681 call mpi_file_write(slice_fh, 0, 1, mpi_integer, status, ierrmpi)
682 call mpi_file_write(slice_fh, nw, 1, mpi_integer, status, ierrmpi)
683 call mpi_file_write(slice_fh, ndir, 1, mpi_integer, status, ierrmpi)
684 call mpi_file_write(slice_fh, ndim-1, 1, mpi_integer, status, ierrmpi)
685 call mpi_file_write(slice_fh, levmax_sub, 1, mpi_integer, status, ierrmpi)
686 call mpi_file_write(slice_fh, nsubleafs, 1, mpi_integer, status, ierrmpi)
687 call mpi_file_write(slice_fh, nsubparents, 1, mpi_integer, status, ierrmpi)
688 call mpi_file_write(slice_fh, it, 1, mpi_integer, status, ierrmpi)
689 ! Note: It is nice when this double has an even number of 4 byte
690 ! integers before it (for alignment)
691 call mpi_file_write(slice_fh, global_time, 1, mpi_double_precision, status, ierrmpi)
692
693 call mpi_file_write(slice_fh, xprobminsub, ndim-1, &
694 mpi_double_precision, status, ierrmpi)
695 call mpi_file_write(slice_fh, xprobmaxsub, ndim-1, &
696 mpi_double_precision, status, ierrmpi)
697 call mpi_file_write(slice_fh, domain_nx_sub, ndim-1, &
698 mpi_integer, status, ierrmpi)
699 call mpi_file_write(slice_fh, block_nx_sub, ndim-1, &
700 mpi_integer, status, ierrmpi)
701
702 ! Periodicity (assume all variables are periodic if one is)
703 call mpi_file_write(slice_fh, periodbsub, ndim-1, mpi_logical, status, ierrmpi)
704
705 ! Geometry
706 call mpi_file_write(slice_fh, geometry_name(1:name_len), &
707 name_len, mpi_character, status, ierrmpi)
708
709 ! Write stagger grid mark
710 call mpi_file_write(slice_fh, stagger_grid, 1, mpi_logical, status, ierrmpi)
711
712 do iw = 1, nw
713 ! using directly trim(adjustl((cons_wnames(iw)))) in MPI_FILE_WRITE call
714 ! does not work, there will be trailing characters
715 dname = trim(adjustl((cons_wnames(iw))))
716 call mpi_file_write(slice_fh, dname, name_len, mpi_character, status, ierrmpi)
717 end do
718
719 ! Physics related information
720 call mpi_file_write(slice_fh, physics_type, name_len, mpi_character, status, ierrmpi)
721
722 ! Format:
723 ! integer :: n_par
724 ! double precision :: values(n_par)
725 ! character(n_par * name_len) :: names
726 call phys_write_info(slice_fh)
727
728 ! Write snapshotnext etc., which is useful for restarting.
729 ! Note we add one, since snapshotnext is updated *after* this procedure
730 if(pass_wall_time.or.save_now) then
731 call mpi_file_write(slice_fh, snapshotnext, 1, mpi_integer, status, ierrmpi)
732 else
733 call mpi_file_write(slice_fh, snapshotnext+1, 1, mpi_integer, status, ierrmpi)
734 end if
735 call mpi_file_write(slice_fh, slicenext, 1, mpi_integer, status, ierrmpi)
736 call mpi_file_write(slice_fh, collapsenext, 1, mpi_integer, status, ierrmpi)
737
738 call mpi_file_get_position(slice_fh, tree_offset, ierrmpi)
739 call mpi_file_write(slice_fh, isleaf_sub, nsubleafs+nsubparents, mpi_logical, status, ierrmpi)
740 call mpi_file_write(slice_fh, block_lvl_sub, nsubleafs, mpi_integer, status, ierrmpi)
741 call mpi_file_write(slice_fh, block_ig_sub, size(block_ig_sub), mpi_integer, status, ierrmpi)
742
743 ! add the block_offset(above) to block_offset_sub
744 call mpi_file_get_position(slice_fh, block_offset, ierrmpi)
745 block_offset_sub(:) = block_offset_sub(:) + block_offset + &
746 int(nsubleafs*2*size_int, kind=mpi_offset_kind)
747
748 call mpi_file_write(slice_fh, block_offset_sub, nsubleafs, mpi_offset, status, ierrmpi)
749
750 ! get the right block_offset
751 call mpi_file_get_position(slice_fh, block_offset, ierrmpi)
752 if (block_offset - tree_offset /= &
753 (nsubleafs + nsubparents) * size_logical + &
754 nsubleafs * ((ndim-1+1) * size_int + 2 * size_int)) then
755 print *, "Warning: MPI_OFFSET type /= 8 bytes"
756 print *, "This *could* cause problems when reading .dat files"
757 end if
758
759 ! come back to write offsets
760 call mpi_file_seek(slice_fh, int(size_int,kind=mpi_offset_kind), mpi_seek_set, ierrmpi)
761 call mpi_file_write(slice_fh, int(tree_offset), 1, mpi_integer, status, ierrmpi)
762 call mpi_file_write(slice_fh, int(block_offset), 1, mpi_integer, status, ierrmpi)
763
764 call mpi_file_close(slice_fh,ierrmpi)
765 end if
766
767 ! Broadcast nsubleafs and nsubparents to all processors
768 call mpi_bcast(nsubleafs, 1, mpi_integer, 0, icomm, ierrmpi)
769 call mpi_bcast(nsubparents, 1, mpi_integer, 0, icomm, ierrmpi)
770
771 ! Allocate block_offset_sub on all processors (if not already allocated on processor 0)
772 if (mype /= 0) then
773 allocate(block_offset_sub(nsubleafs))
774 end if
775
776 ! Broadcast block_offset_sub from processor 0 to all processors
777 call mpi_bcast(block_offset_sub, nsubleafs, mpi_offset, 0, icomm, ierrmpi)
778
779 amode=ior(mpi_mode_create,mpi_mode_wronly)
780 call mpi_file_open(icomm,filename,amode,mpi_info_null,slice_fh,ierrmpi)
781 iorequest=mpi_request_null
782 iwrite=0
783
784 do jgrid=1,njgrid
785 iwrite=iwrite+1
786 offset=block_offset_sub(int(morton_sub_start(mype)+jgrid-1))
787 call mpi_file_iwrite_at(slice_fh,offset,0,2*(ndim-1),mpi_integer,iorequest(iwrite),ierrmpi)
788 iwrite=iwrite+1
789 offset1=offset+int(2*(ndim-1)*size_int,kind=mpi_offset_kind)
790 call mpi_file_iwrite_at(slice_fh,offset1,ps_sub(jgrid)%w,1,type_subblock_io, &
791 iorequest(iwrite),ierrmpi)
792 end do
793
794 if (iwrite>0) call mpi_waitall(iwrite,iorequest,iostatus,ierrmpi)
795 call mpi_barrier(icomm, ierrmpi)
796 call mpi_file_close(slice_fh,ierrmpi)
797
798 deallocate(block_offset_sub)
799 if (mype == 0) then
800 deallocate(block_ig_sub)
801 deallocate(block_lvl_sub)
802 deallocate(isleaf_sub)
803 end if
804
805 end subroutine put_slice_dat
806
807 subroutine put_slice_zerod
808
810 integer:: iw
811 integer :: amode, iwrite, status(MPI_STATUS_SIZE)
812 logical, save :: opened=.false.
813 character(len=name_len) :: wnamei(1:nw+nwauxio),xandwnamei(1:ndim+nw+nwauxio)
814 character(len=1024) :: outfilehead
815 character(len=1024) ::line, data
816 character(len=1024) :: filename, xlabel
817 character(len=79) :: xxlabel
818
819 {^ifoned
820 ! generate filename:
821 write(xlabel,"(D9.2)")xslice
822 xxlabel=trim(xlabel)
823 if(xslice>=zero)then
824 write(xxlabel(1:1),"(a)") "+"
825 endif
826 write(filename,"(a,i1.1,a,i4.4,a)") trim(base_filename)//'_d',dir,'_x'//trim(xxlabel)//'.csv'
827
828 ! Open for header:
829 if (.not. opened .and. mype==0) then
830 amode=ior(mpi_mode_create,mpi_mode_wronly)
831 amode=ior(amode,mpi_mode_append)
832 call mpi_file_open(mpi_comm_self,filename,amode, &
833 mpi_info_null,slice_fh,ierrmpi)
834 ! Create the header:
835 call getheadernames(wnamei,xandwnamei,outfilehead)
836 line = 'time'
837 do iw=1,nw+nwauxio+ndim
838 line = trim(line)//', '//trim(xandwnamei(iw))
839 enddo
840 ! Write header:
841 call mpi_file_write(slice_fh,line,len_trim(line), &
842 mpi_character,status,ierrmpi)
843 call mpi_file_write(slice_fh,achar(10),1,mpi_character,status,ierrmpi)
844 call mpi_file_close(slice_fh, ierrmpi)
845 opened=.true.
846 endif ! opened
847 ! Now we let the processor, which holds the data, write - therefore barrier.
848 call mpi_barrier(icomm,ierrmpi)
849
850 ! There should be only one processor holding the point:
851 if (njgrid > 0) then
852 amode=ior(mpi_mode_create,mpi_mode_wronly)
853 amode=ior(amode,mpi_mode_append)
854 call mpi_file_open(mpi_comm_self,filename,amode, &
855 mpi_info_null,slice_fh,ierrmpi)
856
857 ! Format the line:
858 write(data,"(es14.6)")global_time
859 line = trim(data)
860 write(data,"(es14.6)")ps_sub(1)%x
861 line = trim(line)//', '//trim(data)
862 do iw = 1,nw+nwauxio
863 write(data,"(es14.6)")ps_sub(1)%w(iw)*normconv(iw)
864 line = trim(line)//', '//trim(data)
865 end do
866 !
867 call mpi_file_write(slice_fh,trim(line),len_trim(line), &
868 mpi_character,status,ierrmpi)
869 call mpi_file_write(slice_fh,achar(10),1,mpi_character,status,ierrmpi)
870 call mpi_file_close(slice_fh, ierrmpi)
871 end if ! Njgrid > 0
872 }
873 end subroutine put_slice_zerod
874
875 end subroutine put_slice
876
877 subroutine select_slice(dir,xslice,writeonly,file_handle,normconv)
880 integer, intent(in) :: dir
881 double precision, intent(in) :: xslice
882 integer, intent(in) :: file_handle
883 logical, intent(in) :: writeonly
884 double precision,dimension(0:nw+nwauxio),intent(out) :: normconv
885 ! .. local ..
886 integer :: ig^D, jgrid, slice_fh, ipe, mylevmax
887 integer, dimension(nlevelshi) :: igslice
888 integer, allocatable :: iglevel1_sub_sfc(:)
889 integer :: ig^L, ngsub1, jglevel1, igsorted, kgrid, ngmax
890
891 kgrid = 0
892 jgrid = 0
893 mylevmax = 0
894 jglevel1 = 0
895
896 ! Find the global slice index for every level:
897 call get_igslice(dir,xslice,igslice)
898 igmin^d=1;igmax^d=ng^d(1);
899
900 ! Traverse forest to find grids indicating slice:
901 select case(dir)
902 {case(^d)
903 igmin^d=igslice(1);
904 igmax^d=igslice(1);
905 \}
906 case default
907 call mpistop("slice direction not clear in select_slice")
908 end select
909
910 ! total number of level 1 blocks in the slice
911 ngsub1={(igmax^d-igmin^d+1)|*}
912 allocate(iglevel1_sub_sfc(ngsub1))
913 if (.not. allocated(sfc_sub)) then
914 {^ifthreed
915 ngmax = ngsub1*(4**refine_max_level-1)
916 \}
917 {^iftwod
918 ngmax = ngsub1*(2**refine_max_level-1)
919 \}
920 {^ifoned
921 ngmax = ngsub1
922 \}
923 allocate(sfc_sub(0:ngmax,1:3)) ! the first three are nleafs, nparents and isleaf (0 or 1)
924 sfc_sub = 0
925 end if
926
927 ! get the sfc index of the level 1 blocks in the slice
928 {do ig^db=igmin^db,igmax^db\}
929 jglevel1=jglevel1+1
930 iglevel1_sub_sfc(jglevel1)=iglevel1_sfc(ig^d)
931 {end do \}
932
933 ! Sort the array from small to big using quicksort
934 call quicksort_integer(iglevel1_sub_sfc, 1, ngsub1)
935
936 ! traverse level 1 blocks depending on the sub morton order
937 do igsorted=1,ngsub1
938 {ig^d=sfc_iglevel1(^d,iglevel1_sub_sfc(igsorted))\}
939 call traverse_slice(tree_root(ig^d))
940 end do
941 deallocate(iglevel1_sub_sfc)
942
943 if (.not.writeonly) then
944 ! Synchronize the levmax_sub for output (only rank 0 needs it):
945 levmax_sub = mylevmax
946 call mpi_allreduce(mpi_in_place,levmax_sub,1,mpi_integer,mpi_max,icomm,ierrmpi)
947
948 ! Communicate the subgrid indices according to new Morton sub-sfc:
949 morton_sub_start(:) = 1
950 do ipe=0,npe-1
951 call mpi_gather(jgrid,1,mpi_integer,morton_sub_stop,1,mpi_integer,ipe,icomm,ierrmpi)
952 end do
953
954 do ipe = 0, npe-2
957 end do
958 end if
959
960 contains
961
962 recursive subroutine traverse_slice(tree)
963 implicit none
964 type(tree_node_ptr) :: tree
965 integer :: ic^d
966 integer, dimension(MPI_STATUS_SIZE) :: status
967
968 ! if (writeonly) then
969 ! call MPI_FILE_WRITE(file_handle,tree%node%leaf,1,MPI_LOGICAL, &
970 ! status,ierrmpi)
971 ! end if
972 if (mype==0) then
973 kgrid = kgrid + 1
974 sfc_sub(kgrid,1) = tree%node%igrid
975 sfc_sub(kgrid,2) = tree%node%ipe
976 if (tree%node%leaf) then
977 sfc_sub(0,1) = sfc_sub(0,1) + 1
978 sfc_sub(kgrid,3) = 1
979 else
980 sfc_sub(0,2) = sfc_sub(0,2) + 1
981 sfc_sub(kgrid,3) = 0
982 end if
983 end if
984
985 if (tree%node%leaf) then
986 if (tree%node%ipe == mype.and..not.writeonly) then
987 mylevmax = max(mylevmax,tree%node%level)
988 call fill_subnode(tree%node%igrid,tree%node%active,jgrid,dir,xslice,normconv)
989 end if
990 return
991 end if
992 ! We are out for leaves now, continue for branches
993
994 ! Get the correct child:
995 select case (dir)
996 {case (^d)
997 ic^d = igslice(tree%node%level+1) - 2 * tree%node%ig^d + 2
998 \}
999 case default
1000 call mpistop("slice direction not clear in traverse_slice")
1001 end select
1002
1003 ! Recursively descend into the correct child branch:
1004 {^ifthreed
1005 select case(dir)
1006 case (1)
1007 do ic3=1,2
1008 do ic2=1,2
1009 call traverse_slice(tree%node%child(ic^d))
1010 end do
1011 end do
1012 case (2)
1013 do ic3=1,2
1014 do ic1=1,2
1015 call traverse_slice(tree%node%child(ic^d))
1016 end do
1017 end do
1018 case (3)
1019 do ic2=1,2
1020 do ic1=1,2
1021 call traverse_slice(tree%node%child(ic^d))
1022 end do
1023 end do
1024 case default
1025 call mpistop("slice direction not clear in traverse_slice")
1026 end select
1027 }
1028 {^iftwod
1029 select case(dir)
1030 case (1)
1031 do ic2=1,2
1032 call traverse_slice(tree%node%child(ic^d))
1033 end do
1034 case (2)
1035 do ic1=1,2
1036 call traverse_slice(tree%node%child(ic^d))
1037 end do
1038 case default
1039 call mpistop("slice direction not clear in traverse_slice")
1040 end select
1041 }
1042 {^ifoned
1043 call traverse_slice(tree%node%child(ic^d))
1044 }
1045
1046 end subroutine traverse_slice
1047
1048 end subroutine select_slice
1049
1050 subroutine fill_subnode(igrid,active,jgrid,dir,xslice,normconv)
1053 integer, intent(in) :: igrid, dir
1054 integer, intent(inout) :: jgrid
1055 logical, intent(in) :: active
1056 double precision, intent(in) :: xslice
1057 double precision,dimension(0:nw+nwauxio),intent(out) :: normconv
1058 ! .. local ..
1059 double precision, dimension(ixMlo^D-1:ixMhi^D,ndim) :: xC_TMP, xC
1060 double precision, dimension(ixMlo^D:ixMhi^D,ndim) :: xCC_TMP, xCC
1061 double precision, dimension(ixMlo^D-1:ixMhi^D,nw+nwauxio) :: wC_TMP
1062 double precision, dimension(ixMlo^D:ixMhi^D,nw+nwauxio) :: wCC_TMP
1063 double precision, dimension(ixG^T) :: x_save
1064 integer :: ixslice, nwexpand, ixC^L, ixCC^L
1065 logical :: mask(ixG^T)
1066
1067 mask=.false.
1068 ! increase grid-count:
1069 jgrid=jgrid+1
1070 ! Allocate subdim solution array:
1071 if (ndim==1) then
1072 nwexpand = nwauxio
1073 else
1074 if(slice_type/='dat')then
1075 nwexpand = nwauxio
1076 else
1077 nwexpand = 0
1078 endif
1079 end if
1080 call alloc_subnode(jgrid,dir,nwexpand)
1081 call fill_subnode_info(igrid,jgrid,dir)
1082
1083 mask(ixm^t)=.true.
1084
1085 ! Now hunt for the index closest to the slice:
1086 {^ifoned
1087 ixslice = minloc(dabs(xslice-ps(igrid)%x(:,dir)),1,mask(:))
1088 }
1089 {^iftwod
1090 select case (dir)
1091 case (1)
1092 ixslice = minloc(dabs(xslice-ps(igrid)%x(:,ixmlo2,dir)),1,mask(:,ixmlo2))
1093 case (2)
1094 ixslice = minloc(dabs(xslice-ps(igrid)%x(ixmlo1,:,dir)),1,mask(ixmlo1,:))
1095 case default
1096 call mpistop("slice direction not clear in fill_subnode")
1097 end select
1098 }
1099 {^ifthreed
1100 select case (dir)
1101 case (1)
1102 ixslice = minloc(dabs(xslice-ps(igrid)%x(:,ixmlo2,ixmlo3,dir)),1,mask(:,ixmlo2,ixmlo3))
1103 case (2)
1104 ixslice = minloc(dabs(xslice-ps(igrid)%x(ixmlo1,:,ixmlo3,dir)),1,mask(ixmlo1,:,ixmlo3))
1105 case (3)
1106 ixslice = minloc(dabs(xslice-ps(igrid)%x(ixmlo1,ixmlo2,:,dir)),1,mask(ixmlo1,ixmlo2,:))
1107 case default
1108 call mpistop("slice direction not clear in fill_subnode")
1109 end select
1110 }
1111
1112 call calc_x(igrid,xc,xcc)
1113 ! Set the coordinate to be exactly on the slice:
1114 xc(:^d&,dir) = xslice
1115 xcc(:^d&,dir) = xslice
1116 call calc_grid(unitslice,igrid,xc,xcc,xc_tmp,xcc_tmp,wc_tmp,wcc_tmp,normconv,&
1117 ixc^l,ixcc^l,.true.)
1118 ! CC stands for CellCenter
1119 ! C stands for Corner
1120
1121 ! Fill the subdimensional solution and position array:
1122 {^ifoned
1123 ps_sub(jgrid)%w(1:nw+nwexpand) = wcc_tmp(ixslice,1:nw+nwexpand)
1124 ps_sub(jgrid)%x(1:ndim) = xcc_tmp(ixslice,1:ndim)
1125 ps_sub(jgrid)%wC(1:nw+nwexpand) = wc_tmp(ixslice,1:nw+nwexpand)
1126 ps_sub(jgrid)%xC(1:ndim) = xc_tmp(ixslice,1:ndim)
1127 }
1128 {^iftwod
1129 select case (dir)
1130 case (1)
1131 ps_sub(jgrid)%w(ixccmin2:ixccmax2,1:nw+nwexpand) &
1132 = wcc_tmp(ixslice,ixccmin2:ixccmax2,1:nw+nwexpand)
1133 ps_sub(jgrid)%x(ixccmin2:ixccmax2,1:ndim) &
1134 = xcc_tmp(ixslice,ixccmin2:ixccmax2,1:ndim)
1135 ps_sub(jgrid)%wC(ixcmin2:ixcmax2,1:nw+nwexpand) &
1136 = wc_tmp(ixslice,ixcmin2:ixcmax2,1:nw+nwexpand)
1137 ps_sub(jgrid)%xC(ixcmin2:ixcmax2,1:ndim) &
1138 = xc_tmp(ixslice,ixcmin2:ixcmax2,1:ndim)
1139 case (2)
1140 ps_sub(jgrid)%w(ixccmin1:ixccmax1,1:nw+nwexpand) &
1141 = wcc_tmp(ixccmin1:ixccmax1,ixslice,1:nw+nwexpand)
1142 ps_sub(jgrid)%x(ixccmin1:ixccmax1,1:ndim) &
1143 = xcc_tmp(ixccmin1:ixccmax1,ixslice,1:ndim)
1144 ps_sub(jgrid)%wC(ixcmin1:ixcmax1,1:nw+nwexpand) &
1145 = wc_tmp(ixcmin1:ixcmax1,ixslice,1:nw+nwexpand)
1146 ps_sub(jgrid)%xC(ixcmin1:ixcmax1,1:ndim) &
1147 = xc_tmp(ixcmin1:ixcmax1,ixslice,1:ndim)
1148 case default
1149 call mpistop("slice direction not clear in fill_subnode")
1150 end select
1151 }
1152 {^ifthreed
1153 select case (dir)
1154 case (1)
1155 ps_sub(jgrid)%w(ixccmin2:ixccmax2,ixccmin3:ixccmax3,1:nw+nwexpand) = &
1156 wcc_tmp(ixslice,ixccmin2:ixccmax2,ixccmin3:ixccmax3,1:nw+nwexpand)
1157 ps_sub(jgrid)%x(ixccmin2:ixccmax2,ixccmin3:ixccmax3,1:ndim) = &
1158 xcc_tmp(ixslice,ixccmin2:ixccmax2,ixccmin3:ixccmax3,1:ndim)
1159 ps_sub(jgrid)%wC(ixcmin2:ixcmax2,ixcmin3:ixcmax3,1:nw+nwexpand) = &
1160 wc_tmp(ixslice,ixcmin2:ixcmax2,ixcmin3:ixcmax3,1:nw+nwexpand)
1161 ps_sub(jgrid)%xC(ixcmin2:ixcmax2,ixcmin3:ixcmax3,1:ndim) = &
1162 xc_tmp(ixslice,ixcmin2:ixcmax2,ixcmin3:ixcmax3,1:ndim)
1163 case (2)
1164 ps_sub(jgrid)%w(ixccmin1:ixccmax1,ixccmin3:ixccmax3,1:nw+nwexpand) = &
1165 wcc_tmp(ixccmin1:ixccmax1,ixslice,ixccmin3:ixccmax3,1:nw+nwexpand)
1166 ps_sub(jgrid)%x(ixccmin1:ixccmax1,ixccmin3:ixccmax3,1:ndim) = &
1167 xcc_tmp(ixccmin1:ixccmax1,ixslice,ixccmin3:ixccmax3,1:ndim)
1168 ps_sub(jgrid)%wC(ixcmin1:ixcmax1,ixcmin3:ixcmax3,1:nw+nwexpand) = &
1169 wc_tmp(ixcmin1:ixcmax1,ixslice,ixcmin3:ixcmax3,1:nw+nwexpand)
1170 ps_sub(jgrid)%xC(ixcmin1:ixcmax1,ixcmin3:ixcmax3,1:ndim) = &
1171 xc_tmp(ixcmin1:ixcmax1,ixslice,ixcmin3:ixcmax3,1:ndim)
1172 case (3)
1173 ps_sub(jgrid)%w(ixccmin1:ixccmax1,ixccmin2:ixccmax2,1:nw+nwexpand) = &
1174 wcc_tmp(ixccmin1:ixccmax1,ixccmin2:ixccmax2,ixslice,1:nw+nwexpand)
1175 ps_sub(jgrid)%x(ixccmin1:ixccmax1,ixccmin2:ixccmax2,1:ndim) = &
1176 xcc_tmp(ixccmin1:ixccmax1,ixccmin2:ixccmax2,ixslice,1:ndim)
1177 ps_sub(jgrid)%wC(ixcmin1:ixcmax1,ixcmin2:ixcmax2,1:nw+nwexpand) = &
1178 wc_tmp(ixcmin1:ixcmax1,ixcmin2:ixcmax2,ixslice,1:nw+nwexpand)
1179 ps_sub(jgrid)%xC(ixcmin1:ixcmax1,ixcmin2:ixcmax2,1:ndim) = &
1180 xc_tmp(ixcmin1:ixcmax1,ixcmin2:ixcmax2,ixslice,1:ndim)
1181 case default
1182 call mpistop("slice direction not clear in fill_subnode")
1183 end select
1184 }
1185
1186 end subroutine fill_subnode
1187
1188 subroutine alloc_subnode(jgrid,dir,nwexpand)
1190 integer, intent(in) :: jgrid, dir, nwexpand
1191
1192 ! take care, what comes out is not necessarily a right handed system!
1193 {^ifoned
1194 allocate(ps_sub(jgrid)%w(1:nw+nwexpand),ps_sub(jgrid)%x(1:ndim))
1195 allocate(ps_sub(jgrid)%wC(1:nw+nwexpand),ps_sub(jgrid)%xC(1:ndim))
1196 }
1197 {^iftwod
1198 select case (dir)
1199 case (1)
1200 allocate(ps_sub(jgrid)%w(ixglo2:ixghi2,1:nw+nwexpand),&
1201 ps_sub(jgrid)%x(ixglo2:ixghi2,1:ndim), &
1202 ps_sub(jgrid)%wC(ixglo2:ixghi2,1:nw+nwexpand),&
1203 ps_sub(jgrid)%xC(ixglo2:ixghi2,1:ndim))
1204 case (2)
1205 allocate(ps_sub(jgrid)%w(ixglo1:ixghi1,1:nw+nwexpand),&
1206 ps_sub(jgrid)%x(ixglo1:ixghi1,1:ndim), &
1207 ps_sub(jgrid)%wC(ixglo1:ixghi1,1:nw+nwexpand),&
1208 ps_sub(jgrid)%xC(ixglo1:ixghi1,1:ndim))
1209 case default
1210 call mpistop("slice direction not clear in alloc_subnode")
1211 end select
1212 }
1213 {^ifthreed
1214 select case (dir)
1215 case (1)
1216 allocate(ps_sub(jgrid)%w(ixglo2:ixghi2,ixglo3:ixghi3,1:nw+nwexpand),&
1217 ps_sub(jgrid)%x(ixglo2:ixghi2,ixglo3:ixghi3,1:ndim), &
1218 ps_sub(jgrid)%wC(ixglo2:ixghi2,ixglo3:ixghi3,1:nw+nwexpand),&
1219 ps_sub(jgrid)%xC(ixglo2:ixghi2,ixglo3:ixghi3,1:ndim))
1220 case (2)
1221 allocate(ps_sub(jgrid)%w(ixglo1:ixghi1,ixglo3:ixghi3,1:nw+nwexpand),&
1222 ps_sub(jgrid)%x(ixglo1:ixghi1,ixglo3:ixghi3,1:ndim), &
1223 ps_sub(jgrid)%wC(ixglo1:ixghi1,ixglo3:ixghi3,1:nw+nwexpand),&
1224 ps_sub(jgrid)%xC(ixglo1:ixghi1,ixglo3:ixghi3,1:ndim))
1225 case (3)
1226 allocate(ps_sub(jgrid)%w(ixglo1:ixghi1,ixglo2:ixghi2,1:nw+nwexpand),&
1227 ps_sub(jgrid)%x(ixglo1:ixghi1,ixglo2:ixghi2,1:ndim), &
1228 ps_sub(jgrid)%wC(ixglo1:ixghi1,ixglo2:ixghi2,1:nw+nwexpand),&
1229 ps_sub(jgrid)%xC(ixglo1:ixghi1,ixglo2:ixghi2,1:ndim))
1230 case default
1231 call mpistop("slice direction not clear in alloc_subnode")
1232 end select
1233 }
1234 end subroutine alloc_subnode
1235
1236 subroutine dealloc_subnode(jgrid)
1238 integer, intent(in) :: jgrid
1239
1240 if (jgrid==0) then
1241 call mpistop("trying to delete a non-existing grid in dealloc_subnode")
1242 end if
1243
1244 deallocate(ps_sub(jgrid)%w,ps_sub(jgrid)%x,ps_sub(jgrid)%wC,ps_sub(jgrid)%xC)
1245
1246 ! reset the global node info:
1247 node_sub(:,jgrid)=0
1248 rnode_sub(:,jgrid)=zero
1249
1250 end subroutine dealloc_subnode
1251
1252 subroutine fill_subnode_info(igrid,jgrid,dir)
1254 integer, intent(in) :: igrid,jgrid,dir
1255
1256 node_sub(plevel_,jgrid)=node(plevel_,igrid)
1257 {^ifthreed
1258 select case(dir)
1259 case (1)
1260 node_sub(pig1_,jgrid)=node(pig2_,igrid)
1261 node_sub(pig2_,jgrid)=node(pig3_,igrid)
1262 rnode_sub(rpdx1_,jgrid)=rnode(rpdx2_,igrid)
1263 rnode_sub(rpdx2_,jgrid)=rnode(rpdx3_,igrid)
1264 rnode_sub(rpxmin1_,jgrid)=rnode(rpxmin2_,igrid)
1265 rnode_sub(rpxmin2_,jgrid)=rnode(rpxmin3_,igrid)
1266 case (2)
1267 node_sub(pig1_,jgrid)=node(pig1_,igrid)
1268 node_sub(pig2_,jgrid)=node(pig3_,igrid)
1269 rnode_sub(rpdx1_,jgrid)=rnode(rpdx1_,igrid)
1270 rnode_sub(rpdx2_,jgrid)=rnode(rpdx3_,igrid)
1271 rnode_sub(rpxmin1_,jgrid)=rnode(rpxmin1_,igrid)
1272 rnode_sub(rpxmin2_,jgrid)=rnode(rpxmin3_,igrid)
1273 case (3)
1274 node_sub(pig1_,jgrid)=node(pig1_,igrid)
1275 node_sub(pig2_,jgrid)=node(pig2_,igrid)
1276 rnode_sub(rpdx1_,jgrid)=rnode(rpdx1_,igrid)
1277 rnode_sub(rpdx2_,jgrid)=rnode(rpdx2_,igrid)
1278 rnode_sub(rpxmin1_,jgrid)=rnode(rpxmin1_,igrid)
1279 rnode_sub(rpxmin2_,jgrid)=rnode(rpxmin2_,igrid)
1280 case default
1281 call mpistop("slice direction not clear in fill_subnode_info")
1282 end select
1283 }
1284 {^iftwod
1285 select case(dir)
1286 case (1)
1287 node_sub(pig1_,jgrid)=node(pig2_,igrid)
1288 rnode_sub(rpdx1_,jgrid)=rnode(rpdx2_,igrid)
1289 rnode_sub(rpxmin1_,jgrid)=rnode(rpxmin2_,igrid)
1290 case (2)
1291 node_sub(pig1_,jgrid)=node(pig1_,igrid)
1292 rnode_sub(rpdx1_,jgrid)=rnode(rpdx1_,igrid)
1293 rnode_sub(rpxmin1_,jgrid)=rnode(rpxmin1_,igrid)
1294 case default
1295 call mpistop("slice direction not clear in fill_subnode_info")
1296 end select
1297 }
1298 {^ifoned
1299 }
1300
1301 end subroutine fill_subnode_info
1302
1303 subroutine get_igslice(dir,x,igslice)
1305 integer, intent(in) :: dir
1306 double precision, intent(in) :: x
1307 integer, dimension(nlevelshi), intent(out) :: igslice
1308 ! .. local ..
1309 double precision :: distance
1310 integer :: level
1311 !double precision :: xsgrid(ndim,nlevelshi),qs(ndim),xmgrid(ndim,3),xnew,xlgrid(ndim,3)
1312 !integer :: nbefore
1313
1314 if (x.ne.x) &
1315 call mpistop("get_igslice: your slice position is NaN!")
1316
1317 select case (dir)
1318 {case (^d)
1319 select case (stretch_type(^d))
1320 case (stretch_none) ! This dimension is unstretched
1321 if (x<=xprobmin^d) then
1322 igslice=1
1323 else if (x>=xprobmax^d) then
1324 do level=1,refine_max_level
1325 igslice(level)=ng^d(level)
1326 end do
1327 else
1328 distance=x-xprobmin^d
1329 do level = 1, refine_max_level
1330 igslice(level) = ceiling(distance/dg^d(level))
1331 end do
1332 end if
1333
1334 case (stretch_uni) ! Uniform stretching
1335
1336 if (x<=xprobmin^d) then
1337 igslice=1
1338 else if (x>=xprobmax^d) then
1339 do level=1,refine_max_level
1340 igslice(level)=ng^d(level)
1341 end do
1342 else
1343 distance=x-xprobmin^d
1344 do level=1,refine_max_level
1345 igslice(level)= ceiling(dlog(distance/dxfirst(level,^d)*(qstretch(level,^d)-1.d0)+1.d0)&
1346 /dlog(qstretch(level,^d))/dble(block_nx^d))
1347 end do
1348 end if
1349
1350 case (stretch_symm) ! Symmetric stretching
1351 ! symmetric stretch about 0.5*(xprobmin+xprobmax)
1352 if (x<=xprobmin^d) then
1353 igslice=1
1354 else if (x>=xprobmax^d) then
1355 do level=1,refine_max_level
1356 igslice(level)=ng^d(level)
1357 end do
1358 else if(x<xprobmin^d+xstretch^d) then
1359 distance=xprobmin^d+xstretch^d-x
1360 ! stretch to left from xprobmin+xstretch
1361 do level=1,refine_max_level
1362 igslice(level) = nstretchedblocks(level,^d)/2-int(dlog(distance*(qstretch(level,^d)-one)/&
1363 dxfirst(level,^d)+one)/dlog(qstretch(level,^d))/dble(block_nx^d))
1364 end do
1365 else if(x>xprobmax^d-xstretch^d) then
1366 distance=x-xprobmax^d+xstretch^d
1367 ! stretch to right from xprobmax-xstretch
1368 do level=1,refine_max_level
1369 igslice(level) = ceiling(dlog(distance*(qstretch(level,^d)-one)/&
1370 dxfirst(level,^d)+one)/dlog(qstretch(level,^d))/dble(block_nx^d))+ng^d(level)-nstretchedblocks(level,^d)/2
1371 end do
1372 else
1373 ! possible non-stretched central part
1374 distance=x-xprobmin^d-xstretch^d
1375 do level=1,refine_max_level
1376 igslice(level)=nstretchedblocks(level,^d)/2+ceiling(distance/dg^d(level))
1377 end do
1378 end if
1379 !xsgrid(^D,1)=half*(xprobmax^D-xprobmin^D)&
1380 ! /(half*domain_nx^D-half*nstretchedblocks_baselevel(^D)*block_nx^D &
1381 ! +(one-qstretch_baselevel(^D)**(half*nstretchedblocks_baselevel(^D)*block_nx^D)) &
1382 ! /(one-qstretch_baselevel(^D)))
1383 !xmgrid(^D,1)=xsgrid(^D,1)*(domain_nx^D-nstretchedblocks_baselevel(^D)*block_nx^D)
1384 !xlgrid(^D,1)=xsgrid(^D,1)
1385
1386 !if (x .ge. xmgrid(^D,1)*half) then
1387 ! xnew=x-xmgrid(^D,1)*half
1388 ! do level=1,refine_max_level
1389 ! qs(^D)=qstretch_baselevel(^D)**(one/2.d0**(level-1))
1390 ! if (level .gt. 1) xsgrid(^D,level)=xsgrid(^D,level-1)/(one+qs(^D))
1391 ! nbefore=(domain_nx^D/block_nx^D-nstretchedblocks_baselevel(^D)/2)*2**(level-1)
1392 ! igslice(level)=floor((dlog(one-xnew/xsgrid(^D,level)*(one-qs(^D)))/dlog(qs(^D))-1.e-16)/block_nx^D)+1+nbefore
1393 ! if (x>=xprobmax^D) igslice(level)=domain_nx^D*2**(level-1)
1394 ! end do
1395 !else if(x .ge. -xmgrid(^D,1)*half) then
1396 ! xnew=x+xmgrid(^D,1)*half
1397 ! do level=1,refine_max_level
1398 ! nbefore=nstretchedblocks_baselevel(^D)/2*2**(level-1)
1399 ! if (level .gt. 1) xlgrid(^D,level)=xlgrid(^D,level-1)/2.d0
1400 ! igslice(level)=floor(((xnew-1.e-16)/xlgrid(^D,level))/block_nx^D)+1+nbefore
1401 ! end do
1402 !else
1403 ! xnew=xmgrid(^D,1)*half-x
1404 ! do level=1,refine_max_level
1405 ! qs(^D)=qstretch_baselevel(^D)**(one/2.d0**(level-1))
1406 ! if (level .gt. 1) xsgrid(^D,level)=xsgrid(^D,level-1)/(one+qs(^D))
1407 ! igslice(level)=nstretchedblocks_baselevel(^D)/2*2**(level-1) &
1408 ! -floor((dlog(one-xnew/xsgrid(^D,level)*(one-qs(^D)))/dlog(qs(^D)))/block_nx^D)
1409 ! if (x<=xprobmin^D) igslice(level)=1
1410 ! end do
1411 !end if
1412
1413 case default
1414 call mpistop("stretch type not supported by get_igslice")
1415 end select
1416 \}
1417 case default
1418 call mpistop("slice direction not clear in get_igslice")
1419 end select
1420
1421 end subroutine get_igslice
1422
1423 double precision function roundoff_minmax(val,minval,maxval)
1424 implicit none
1425 double precision,intent(in) :: val, minval, maxval
1426
1427 roundoff_minmax = val
1428
1429 if (abs(roundoff_minmax) .le. minval) then
1430 roundoff_minmax = 0.0d0
1431 end if
1432
1433 if (roundoff_minmax .gt. maxval) then
1434 roundoff_minmax = maxval
1435 else if (roundoff_minmax .lt. -maxval) then
1436 roundoff_minmax = -maxval
1437 end if
1438
1439 ! Replace NaN with maxval (e.g. Paraview chokes on ASCII NaN):
1441
1442 end function roundoff_minmax
1443
1444 recursive subroutine quicksort_integer(arr, left, right)
1445 implicit none
1446 integer, intent(inout) :: arr(:)
1447 integer, intent(in) :: left, right
1448 ! .. local ..
1449 integer :: pivot, i, j, temp
1450
1451 if (left < right) then
1452 ! Choose pivot (using middle element)
1453 pivot = arr((left + right) / 2)
1454 i = left - 1
1455 j = right + 1
1456
1457 do
1458 ! Find elements to swap
1459 do
1460 i = i + 1
1461 if (arr(i) >= pivot) exit
1462 end do
1463
1464 do
1465 j = j - 1
1466 if (arr(j) <= pivot) exit
1467 end do
1468
1469 if (i >= j) exit
1470
1471 ! Swap elements
1472 temp = arr(i)
1473 arr(i) = arr(j)
1474 arr(j) = temp
1475 end do
1476
1477 ! Recursively sort subarrays
1478 call quicksort_integer(arr, left, j)
1479 call quicksort_integer(arr, j + 1, right)
1480 end if
1481
1482 end subroutine quicksort_integer
1483
1484
1485end module mod_slice
1486
recursive subroutine traverse_slice(tree)
Definition mod_slice.t:963
subroutine put_slice_vtu
Definition mod_slice.t:244
Module with basic data types used in amrvac.
integer, parameter size_double
Size (in bytes) of double precision real.
integer, parameter name_len
Default length for names (of e.g. variables)
integer, parameter size_logical
Size (in bytes) of default logical.
integer, parameter size_int
Size (in bytes) of default integer.
Handles computations for coordinates and variables in output.
subroutine getheadernames(wnamei, xandwnamei, outfilehead)
get all variables names
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 calc_x(igrid, xc, xcc)
computes cell corner (xC) and cell center (xCC) coordinates
subroutine, public mpistop(message)
Exit MPI-AMRVAC with an error message.
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
integer, dimension(:^d &), allocatable, save iglevel1_sfc
iglevel1_sfc(ig^D) gives the Morton number for grid ig^D
Definition mod_forest.t:50
integer, dimension(:,:), allocatable, save sfc_iglevel1
Space filling curve for level 1 grid. sfc_iglevel1(^D, MN) gives ig^D (the spatial index of the grid)
Definition mod_forest.t:47
type(tree_node_ptr), dimension(:^d &), allocatable, save tree_root
Pointers to the coarse grid.
Definition mod_forest.t:29
type(tree_node_ptr), dimension(:,:), allocatable, save igrid_to_node
Array to go from an [igrid, ipe] index to a node pointer.
Definition mod_forest.t:32
This module contains definitions of global parameters and variables and some generic functions/subrou...
integer, parameter stretch_uni
Unidirectional stretching from a side.
double precision global_time
The global simulation time.
double precision xstretch
physical extent of stretched border in symmetric stretching
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, dimension(:,:), allocatable nstretchedblocks
(even) number of (symmetrically) stretched blocks per level and dimension
integer icomm
The MPI communicator.
integer, dimension(:), allocatable ng
number of grid blocks in domain per dimension, in array over levels
integer mype
The rank of the current MPI task.
integer block_nx
number of cells for each dimension in grid block excluding ghostcells
double precision, dimension(:), allocatable, parameter d
integer slicenext
IO: slice output number/label.
integer ixm
the mesh range of a physical block without ghost cells
integer ierrmpi
A global MPI error return code.
double precision, dimension(:,:), allocatable qstretch
Stretching factors and first cell size for each AMR level and dimension.
integer npe
The number of MPI tasks.
integer nwauxio
Number of auxiliary variables that are only included in output.
integer, parameter stretch_none
No stretching.
logical, dimension(:), allocatable w_write
if true write the w variable in output
logical, dimension(ndim) periodb
True for dimensions with periodic boundaries.
double precision, dimension(:,:), allocatable rnode
Corner coordinates.
integer nghostcells
Number of ghost cells surrounding a grid.
integer, dimension(ndim) stretch_type
What kind of stretching is used per dimension.
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, parameter stretch_symm
Symmetric stretching around the center.
double precision, dimension(:,:), allocatable dxfirst
integer, dimension(:,:), allocatable node
integer, dimension(:,:), allocatable node_sub
logical save_now
whether a manually inserted snapshot is saved
integer, parameter version_number
Version number of the .dat file output.
This module defines the procedures of a physics module. It contains function pointers for the various...
Definition mod_physics.t:4
procedure(sub_write_info), pointer phys_write_info
Definition mod_physics.t:77
character(len=name_len) physics_type
String describing the physics type of the simulation.
Definition mod_physics.t:50
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:1253
subroutine get_igslice(dir, x, igslice)
Definition mod_slice.t:1304
subroutine fill_subnode(igrid, active, jgrid, dir, xslice, normconv)
Definition mod_slice.t:1051
integer, parameter nslicemax
Maximum number of slices.
Definition mod_slice.t:11
subroutine alloc_subnode(jgrid, dir, nwexpand)
Definition mod_slice.t:1189
subroutine put_slice(dir, xslice)
Definition mod_slice.t:45
subroutine write_slice
Definition mod_slice.t:31
integer, dimension(:,:), allocatable sfc_sub
igrid and ipe of the leaf blocks on the slice, with the first (1,1:2) stores the nleafs and nparents
Definition mod_slice.t:26
double precision, dimension(1000) slicecoord
Slice coordinates, see Slice output.
Definition mod_slice.t:8
character(len=std_len) slice_type
choose data type of slice: vtu, vtuCC, dat, or csv
Definition mod_slice.t:23
integer, dimension(nslicemax) slicedir
The slice direction for each slice.
Definition mod_slice.t:17
double precision function roundoff_minmax(val, minval, maxval)
Definition mod_slice.t:1424
integer nslices
Number of slices to output.
Definition mod_slice.t:14
subroutine select_slice(dir, xslice, writeonly, file_handle, normconv)
Definition mod_slice.t:878
subroutine dealloc_subnode(jgrid)
Definition mod_slice.t:1237
recursive subroutine quicksort_integer(arr, left, right)
Definition mod_slice.t:1445
Pointer to a tree_node.
Definition mod_forest.t:6
The data structure that contains information about a tree node/grid block.
Definition mod_forest.t:11