19 logical,
save :: firstcollapse=.true.
22 if (firstcollapse)
then
37 integer,
intent(in) :: dir
39 double precision,
dimension(0:nw+nwauxio) :: normconv
40 integer :: jgrid, igrid, Morton_no
42 if(.not.
slab)
call mpistop(
"collapse only for slab cartesian cases")
57 call mpi_reduce(mpi_in_place,collapseddata,
size(collapseddata),mpi_double_precision,mpi_sum,0,
icomm,
ierrmpi)
59 call mpi_reduce(collapseddata,collapseddata,
size(collapseddata),mpi_double_precision,mpi_sum,0,
icomm,
ierrmpi)
70 call mpistop(
"Unknown filetype for collapsed views")
73 call mpistop(
"sorry, 1D collapsed output routine not yet implemented (should be easy)...")
80 deallocate(collapseddata)
87 integer,
intent(in) :: dir
88 double precision,
dimension(0:nw+nwauxio),
intent(in):: normconv
91 double precision :: dxdim^DM, xdim^LIM^DM
95 integer,
dimension(ndim) :: myshape
97 character(len=name_len) :: wnamei(1:nw+nwauxio),xandwnamei(1:ndim+nw+nwauxio)
98 character(len=1024) :: filename, outfilehead, line
109 xdim^lim1=
xprob^lim2;
110 xdim^lim2=
xprob^lim3;
114 xdim^lim1=
xprob^lim1;
115 xdim^lim2=
xprob^lim3;
119 xdim^lim1=
xprob^lim1;
120 xdim^lim2=
xprob^lim2;
124 xdim^lim1=
xprob^lim2;
125 xdim^lim2=
xprob^lim3;
126 call mpistop(
"slice direction not clear in output_collapsed_csv")
133 xdim^lim1=
xprob^lim2;
136 xdim^lim1=
xprob^lim1;
139 xdim^lim1=
xprob^lim2;
140 call mpistop(
"slice direction not clear in output_collapsed_csv")
145 if(.not.fileopen)
then
147 write(filename,
"(a,i1.1,a,i1.1,a,i4.4,a)") trim(
base_filename) //
'_d', &
149 open(
unitcollapse,file=filename,status=
'unknown',form=
'formatted')
154 do iw=1,ndim+nw+nwauxio-1
155 if (iw .eq. dir) cycle
156 line = trim(line)//trim(xandwnamei(iw))//
', '
158 line = trim(line)//trim(xandwnamei(ndim+nw+nwauxio))
160 myshape = shape(collapseddata)
162 {^dm&
do ix^dmb = 1,myshape(^dmb)\}
165 {^dm& dxdim^dm*dble(ix^dm)+xdimmin^dm}, &
166 (collapseddata(ix^dm,iw)*normconv(iw),iw=1,nw+nwauxio)
177 integer,
intent(in) :: dir
178 double precision,
dimension(0:nw+nwauxio),
intent(in):: normconv
181 double precision :: dxdim^DM, xdim^LIM^DM
184 double precision :: origin(1:3), spacing(1:3)
186 integer :: wholeExtent(1:6), size_single, length, size_length
188 integer,
dimension(ndim) :: myshape
190 character(len=1024) :: filename, outfilehead, line
191 character(len=name_len) :: wnamei(1:nw+nwauxio),xandwnamei(1:ndim+nw+nwauxio)
204 xdim^lim1=
xprob^lim2;
205 xdim^lim2=
xprob^lim3;
209 xdim^lim1=
xprob^lim1;
210 xdim^lim2=
xprob^lim3;
214 xdim^lim1=
xprob^lim1;
215 xdim^lim2=
xprob^lim2;
219 xdim^lim1=
xprob^lim2;
220 xdim^lim2=
xprob^lim3;
221 call mpistop(
"slice direction not clear in output_collapsed_vti")
228 xdim^lim1=
xprob^lim2;
231 xdim^lim1=
xprob^lim1;
234 xdim^lim1=
xprob^lim2;
235 call mpistop(
"slice direction not clear in output_collapsed_vti")
242 myshape = shape(collapseddata)
247 length = ^dm&myshape(^dm)*
248 length = length*size_single
249 {^dm&wholeextent(^dm*2)=myshape(^dm); }
250 {^dm&spacing(^dm) = dxdim^dm; }
251 {^dm&origin(^
d) = xdimmin^dm; }
255 if(.not.fileopen)
then
257 write(filename,
"(a,i1.1,a,i1.1,a,i4.4,a)") trim(
base_filename)//
'_d',dir,&
259 open(
unitcollapse,file=filename,status=
'unknown',form=
'formatted')
266 write(
unitcollapse,
'(a)',advance=
'no')
'<VTKFile type="ImageData"'
268 write(
unitcollapse,
'(a)')
' version="0.1" byte_order="LittleEndian">'
270 write(
unitcollapse,
'(a)')
' version="0.1" byte_order="BigEndian">'
273 write(
unitcollapse,
'(a,3(1pe14.6),a,6(i10),a,3(1pe14.6),a)')
' <ImageData Origin="',&
274 origin,
'" WholeExtent="',wholeextent,
'" Spacing="',spacing,
'">'
276 write(
unitcollapse,
'(2a)')
'<DataArray type="Float32" Name="TIME" ',&
277 'NumberOfTuples="1" format="ascii">'
284 '<Piece Extent="',wholeextent,
'">'
292 '<DataArray type="Float32" Name="',trim(wnamei(iw)),&
293 '" format="appended" offset="',offset,
'"/>'
294 offset = offset + length + size_length
300 write(
unitcollapse,
'(a)')
'<AppendedData encoding="raw">'
303 open(
unitcollapse,file=filename,access=
'stream',form=
'unformatted',position=
'append')
313 write(
unitcollapse) real(collapseddata(iw)*normconv(iw))
316 write(
unitcollapse) (real(collapseddata(ix1,iw)*normconv(iw)),ix1=1,myshape(1))
319 write(
unitcollapse) ((real(collapseddata(ix1,ix2,iw)*normconv(iw)),ix1=1,&
320 myshape(1)),ix2=1,myshape(2))
325 open(
unitcollapse,file=filename,status=
'unknown',form=
'formatted',position=
'append')
336 integer,
intent(in) :: dir
342 nx^d=ixmhi^d-ixmlo^d+1;
347 allocate(collapseddata(&
350 allocate(collapseddata(&
353 allocate(collapseddata(&
356 call mpistop(
"slice direction not clear in allocate_collapsed")
362 allocate(collapseddata(&
365 allocate(collapseddata(&
368 call mpistop(
"slice direction not clear in allocate_collapsed")
372 allocate(collapseddata(nw+
nwauxio))
381 integer,
intent(in) :: igrid, jgrid, dir
385 integer :: ig^D, level, ig^Dtargetmin, ig^Dtargetmax
386 integer :: ix^Dtarget^LIM, idim^Dtarget^LIM
387 integer :: ixMdim^LLIM^D, ix^Dorig, ix^D
393 {^d& ig^d = tree%node%ig^d; }
394 level = tree%node%level
396 nx^d=ixmhi^d-ixmlo^d+1;
401 ig^dtargetmin = int(dble(ig^d-1)/2.0d0**(level-
collapselevel))+1
402 ig^dtargetmax = ig^dtargetmin
410 ix^dtargetmin = nx^d*(ig^dtargetmin-1)+1
411 ix^dtargetmax = nx^d*ig^dtargetmax
419 idim1target^lim=ix2target^lim;
420 idim2target^lim=ix3target^lim;
421 ixmdim^llim1=
ixm^llim2;
422 ixmdim^llim2=
ixm^llim3;
426 idim1target^lim=ix1target^lim;
427 idim2target^lim=ix3target^lim;
428 ixmdim^llim1=
ixm^llim1;
429 ixmdim^llim2=
ixm^llim3;
433 idim1target^lim=ix1target^lim;
434 idim2target^lim=ix2target^lim;
435 ixmdim^llim1=
ixm^llim1;
436 ixmdim^llim2=
ixm^llim2;
438 call mpistop(
"slice direction not clear in integrate_subnode")
442 do ix1orig = ixmdimlo1,ixmdimhi1
443 do ix2orig = ixmdimlo2,ixmdimhi2
445 collapseddata(ix1,ix2,1:nw+
nwauxio) = collapseddata(ix1,ix2,1:nw+
nwauxio) &
451 do ix^dm = idim^dmtargetmin,idim^dmtargetmax\}
453 collapseddata(ix1,ix2,1:nw+
nwauxio) = collapseddata(ix1,ix2,1:nw+
nwauxio) + ps_sub(jgrid)%w(ix1orig,ix2orig,1:nw+
nwauxio)
461 idim1target^lim=ix2target^lim;
462 ixmdim^llim1=
ixm^llim2;
465 idim1target^lim=ix1target^lim;
466 ixmdim^llim1=
ixm^llim1;
468 call mpistop(
"slice direction not clear in integrate_subnode")
472 do ix1orig = ixmdimlo1,ixmdimhi1
474 collapseddata(ix1,1:nw+
nwauxio) = collapseddata(ix1,1:nw+
nwauxio) &
479 do ix^dm = idim^dmtargetmin,idim^dmtargetmax\}
481 collapseddata(ix1,1:nw+
nwauxio) = collapseddata(ix1,1:nw+
nwauxio) + ps_sub(jgrid)%w(ix1orig,1:nw+
nwauxio)
497 integer,
intent(in) :: igrid, jgrid, dir
498 double precision,
dimension(0:nw+nwauxio),
intent(out) :: normconv
500 double precision :: dx^D
501 double precision,
dimension(ixMlo^D-1:ixMhi^D,ndim) :: xC_TMP, xC
502 double precision,
dimension(ixMlo^D:ixMhi^D,ndim) :: xCC_TMP, xCC
503 double precision,
dimension(ixMlo^D-1:ixMhi^D,nw+nwauxio) :: wC_TMP
504 double precision,
dimension(ixMlo^D:ixMhi^D,nw+nwauxio) :: wCC_TMP
505 integer :: ixC^L, ixCC^L
520 ps_sub(jgrid)%w(ixmlo2:ixmhi2,ixmlo3:ixmhi3,1:nw+
nwauxio) = &
521 ps_sub(jgrid)%w(ixmlo2:ixmhi2,ixmlo3:ixmhi3,1:nw+
nwauxio) &
522 + wcc_tmp(ix,ixmlo2:ixmhi2,ixmlo3:ixmhi3,1:nw+
nwauxio) * dx1
527 ps_sub(jgrid)%w(ixmlo2:ixmhi2,ixmlo3:ixmhi3,iw) = &
528 ps_sub(jgrid)%w(ixmlo2:ixmhi2,ixmlo3:ixmhi3,iw) &
529 + wcc_tmp(ix,ixmlo2:ixmhi2,ixmlo3:ixmhi3,iw) * &
530 ps(igrid)%dx(ix,ixmlo2:ixmhi2,ixmlo3:ixmhi3,1)
534 ps_sub(jgrid)%x(ixmlo2:ixmhi2,ixmlo3:ixmhi3,1:
ndim) = &
535 ps(igrid)%x(ixmlo1,ixmlo2:ixmhi2,ixmlo3:ixmhi3,1:
ndim)
539 ps_sub(jgrid)%w(ixmlo1:ixmhi1,ixmlo3:ixmhi3,1:nw+
nwauxio) = &
540 ps_sub(jgrid)%w(ixmlo1:ixmhi1,ixmlo3:ixmhi3,1:nw+
nwauxio) &
541 + wcc_tmp(ixmlo1:ixmhi1,ix,ixmlo3:ixmhi3,1:nw+
nwauxio) * dx2
546 ps_sub(jgrid)%w(ixmlo1:ixmhi1,ixmlo3:ixmhi3,iw) = &
547 ps_sub(jgrid)%w(ixmlo1:ixmhi1,ixmlo3:ixmhi3,iw) &
548 + wcc_tmp(ixmlo1:ixmhi1,ix,ixmlo3:ixmhi3,iw) * &
549 ps(igrid)%dx(ixmlo1:ixmhi1,ix,ixmlo3:ixmhi3,2)
553 ps_sub(jgrid)%x(ixmlo1:ixmhi1,ixmlo3:ixmhi3,1:
ndim) = &
554 ps(igrid)%x(ixmlo1:ixmhi1,ixmlo2,ixmlo3:ixmhi3,1:
ndim)
558 ps_sub(jgrid)%w(ixmlo1:ixmhi1,ixmlo2:ixmhi2,1:nw+
nwauxio) = &
559 ps_sub(jgrid)%w(ixmlo1:ixmhi1,ixmlo2:ixmhi2,1:nw+
nwauxio) &
560 + wcc_tmp(ixmlo1:ixmhi1,ixmlo2:ixmhi2,ix,1:nw+
nwauxio) * dx3
565 ps_sub(jgrid)%w(ixmlo1:ixmhi1,ixmlo2:ixmhi2,iw) = &
566 ps_sub(jgrid)%w(ixmlo1:ixmhi1,ixmlo2:ixmhi2,iw) &
567 + wcc_tmp(ixmlo1:ixmhi1,ixmlo2:ixmhi2,ix,iw) * &
568 ps(igrid)%dx(ixmlo1:ixmhi1,ixmlo2:ixmhi2,ix,3)
572 ps_sub(jgrid)%x(ixmlo1:ixmhi1,ixmlo2:ixmhi2,1:
ndim) = &
573 ps(igrid)%x(ixmlo1:ixmhi1,ixmlo2:ixmhi2,ixmlo3,1:
ndim)
575 print*,
'subnode, dir: ', dir
576 call mpistop(
"slice direction not clear in collapse_subnode")
584 ps_sub(jgrid)%w(ixmlo2:ixmhi2,1:nw+
nwauxio) = &
585 ps_sub(jgrid)%w(ixmlo2:ixmhi2,1:nw+
nwauxio) &
586 + wcc_tmp(ix,ixmlo2:ixmhi2,1:nw+
nwauxio) * dx1
591 ps_sub(jgrid)%w(ixmlo2:ixmhi2,iw) = &
592 ps_sub(jgrid)%w(ixmlo2:ixmhi2,iw) &
593 + wcc_tmp(ix,ixmlo2:ixmhi2,iw) * &
594 ps(igrid)%dx(ix,ixmlo2:ixmhi2,1)
598 ps_sub(jgrid)%x(ixmlo2:ixmhi2,1:
ndim) = &
599 ps(igrid)%x(ixmlo1,ixmlo2:ixmhi2,1:
ndim)
603 ps_sub(jgrid)%w(ixmlo1:ixmhi1,1:nw+
nwauxio) = &
604 ps_sub(jgrid)%w(ixmlo1:ixmhi1,1:nw+
nwauxio) &
605 + wcc_tmp(ixmlo1:ixmhi1,ix,1:nw+
nwauxio) * dx2
610 ps_sub(jgrid)%w(ixmlo1:ixmhi1,iw) = &
611 ps_sub(jgrid)%w(ixmlo1:ixmhi1,iw) &
612 + wcc_tmp(ixmlo1:ixmhi1,ix,iw) * &
613 ps(igrid)%dx(ixmlo1:ixmhi1,ix,2)
617 ps_sub(jgrid)%x(ixmlo1:ixmhi1,1:
ndim) = &
618 ps(igrid)%x(ixmlo1:ixmhi1,ixmlo2,1:
ndim)
620 call mpistop(
"slice direction not clear in collapse_subnode")
631 ps_sub(jgrid)%w(iw) = ps_sub(jgrid)%w(iw) + wcc_tmp(ix,iw) * ps(igrid)%dx(ix,1)
635 ps_sub(jgrid)%x(1:
ndim) = ps(igrid)%x(ixmlo1,1:
ndim)
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
Collapses D-dimensional output to D-1 view by line-of-sight integration.
subroutine collapse_subnode(igrid, jgrid, dir, normconv)
subroutine output_collapsed_csv(dir, normconv)
subroutine allocate_collapsed(dir)
subroutine put_collapse(dir)
subroutine output_collapsed_vti(dir, normconv)
subroutine write_collapsed
subroutine integrate_subnode(igrid, jgrid, dir)
subroutine, public mpistop(message)
Exit MPI-AMRVAC with an error message.
Module with basic grid data structures.
integer, dimension(:), allocatable, save sfc_to_igrid
Go from a Morton number to an igrid index (for a single processor)
integer, dimension(:), allocatable, save morton_start
First Morton number per processor.
integer, dimension(:), allocatable, save morton_stop
Last Morton number per processor.
type(tree_node_ptr), dimension(:,:), allocatable, save igrid_to_node
Array to go from an [igrid, ipe] index to a node pointer.
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 xprob
minimum and maximum domain boundaries for each dimension
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, dimension(:), allocatable ng
number of grid blocks in domain per dimension, in array over levels
integer mype
The rank of the current MPI task.
logical, dimension(ndim) collapse
If collapse(DIM) is true, generate output integrated over DIM.
double precision, dimension(:), allocatable, parameter d
integer ixm
the mesh range of a physical block without ghost cells
integer ierrmpi
A global MPI error return code.
character(len=std_len) collapse_type
logical slab
Cartesian geometry or not.
integer nwauxio
Number of auxiliary variables that are only included in output.
logical, dimension(:), allocatable w_write
if true write the w variable in output
double precision, dimension(:,:), allocatable rnode
Corner coordinates.
integer, parameter unitcollapse
double precision, dimension(:,:), allocatable dx
integer nghostcells
Number of ghost cells surrounding a grid.
integer collapselevel
The level at which to produce line-integrated / collapsed output.
logical slab_uniform
uniform Cartesian geometry or not (stretched Cartesian)
character(len=std_len) base_filename
Base file name for simulation output, which will be followed by a number.
This module defines the procedures of a physics module. It contains function pointers for the various...
procedure(sub_convert), pointer phys_to_primitive
Writes D-1 slice, can do so in various formats, depending on slice_type.
subroutine alloc_subnode(jgrid, dir, nwexpand)
subroutine dealloc_subnode(jgrid)
Module with all the methods that users can customize in AMRVAC.
procedure(aux_output), pointer usr_aux_output