19 logical,
save :: firstcollapse=.true.
22 if (firstcollapse)
then
37 integer,
intent(in) :: dir
39 integer :: jgrid, igrid, Morton_no
40 double precision,
dimension(0:nw+nwauxio) :: normconv
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
89 character(len=1024) :: filename, outfilehead, line
92 character(len=name_len) :: wnamei(1:nw+nwauxio),xandwnamei(1:ndim+nw+nwauxio)
93 integer,
dimension(ndim) :: myshape
96 double precision :: dxdim^DM, xdim^LIM^DM
108 xdim^lim1=
xprob^lim2;
109 xdim^lim2=
xprob^lim3;
113 xdim^lim1=
xprob^lim1;
114 xdim^lim2=
xprob^lim3;
118 xdim^lim1=
xprob^lim1;
119 xdim^lim2=
xprob^lim2;
123 xdim^lim1=
xprob^lim2;
124 xdim^lim2=
xprob^lim3;
125 call mpistop(
"slice direction not clear in output_collapsed_csv")
132 xdim^lim1=
xprob^lim2;
135 xdim^lim1=
xprob^lim1;
138 xdim^lim1=
xprob^lim2;
139 call mpistop(
"slice direction not clear in output_collapsed_csv")
144 if(.not.fileopen)
then
146 write(filename,
"(a,i1.1,a,i1.1,a,i4.4,a)") trim(
base_filename) //
'_d', &
148 open(
unitcollapse,file=filename,status=
'unknown',form=
'formatted')
153 do iw=1,ndim+nw+nwauxio-1
154 if (iw .eq. dir) cycle
155 line = trim(line)//trim(xandwnamei(iw))//
', '
157 line = trim(line)//trim(xandwnamei(ndim+nw+nwauxio))
159 myshape = shape(collapseddata)
161 {^dm&
do ix^dmb = 1,myshape(^dmb)\}
164 {^dm& dxdim^dm*dble(ix^dm)+xdimmin^dm}, &
165 (collapseddata(ix^dm,iw)*normconv(iw),iw=1,nw+nwauxio)
176 integer,
intent(in) :: dir
177 double precision,
dimension(0:nw+nwauxio),
intent(in):: normconv
178 character(len=1024) :: filename, outfilehead, line
181 character(len=name_len) :: wnamei(1:nw+nwauxio),xandwnamei(1:ndim+nw+nwauxio)
182 integer,
dimension(ndim) :: myshape
185 double precision :: dxdim^DM, xdim^LIM^DM
187 double precision :: origin(1:3), spacing(1:3)
188 integer :: wholeExtent(1:6), size_single, length, size_length
202 xdim^lim1=
xprob^lim2;
203 xdim^lim2=
xprob^lim3;
207 xdim^lim1=
xprob^lim1;
208 xdim^lim2=
xprob^lim3;
212 xdim^lim1=
xprob^lim1;
213 xdim^lim2=
xprob^lim2;
217 xdim^lim1=
xprob^lim2;
218 xdim^lim2=
xprob^lim3;
219 call mpistop(
"slice direction not clear in output_collapsed_vti")
226 xdim^lim1=
xprob^lim2;
229 xdim^lim1=
xprob^lim1;
232 xdim^lim1=
xprob^lim2;
233 call mpistop(
"slice direction not clear in output_collapsed_vti")
240 myshape = shape(collapseddata)
245 length = ^dm&myshape(^dm)*
246 length = length*size_single
247 {^dm&wholeextent(^dm*2)=myshape(^dm); }
248 {^dm&spacing(^dm) = dxdim^dm; }
249 {^dm&origin(^
d) = xdimmin^dm; }
253 if(.not.fileopen)
then
255 write(filename,
"(a,i1.1,a,i1.1,a,i4.4,a)") trim(
base_filename)//
'_d',dir,&
257 open(
unitcollapse,file=filename,status=
'unknown',form=
'formatted')
264 write(
unitcollapse,
'(a)',advance=
'no')
'<VTKFile type="ImageData"'
266 write(
unitcollapse,
'(a)')
' version="0.1" byte_order="LittleEndian">'
268 write(
unitcollapse,
'(a)')
' version="0.1" byte_order="BigEndian">'
271 write(
unitcollapse,
'(a,3(1pe14.6),a,6(i10),a,3(1pe14.6),a)')
' <ImageData Origin="',&
272 origin,
'" WholeExtent="',wholeextent,
'" Spacing="',spacing,
'">'
274 write(
unitcollapse,
'(2a)')
'<DataArray type="Float32" Name="TIME" ',&
275 'NumberOfTuples="1" format="ascii">'
282 '<Piece Extent="',wholeextent,
'">'
290 '<DataArray type="Float32" Name="',trim(wnamei(iw)),&
291 '" format="appended" offset="',offset,
'"/>'
292 offset = offset + length + size_length
298 write(
unitcollapse,
'(a)')
'<AppendedData encoding="raw">'
301 open(
unitcollapse,file=filename,access=
'stream',form=
'unformatted',position=
'append')
311 write(
unitcollapse) real(collapseddata(iw)*normconv(iw))
314 write(
unitcollapse) (real(collapseddata(ix1,iw)*normconv(iw)),ix1=1,myshape(1))
317 write(
unitcollapse) ((real(collapseddata(ix1,ix2,iw)*normconv(iw)),ix1=1,&
318 myshape(1)),ix2=1,myshape(2))
323 open(
unitcollapse,file=filename,status=
'unknown',form=
'formatted',position=
'append')
334 integer,
intent(in) :: dir
340 nx^d=ixmhi^d-ixmlo^d+1;
345 allocate(collapseddata(&
348 allocate(collapseddata(&
351 allocate(collapseddata(&
354 call mpistop(
"slice direction not clear in allocate_collapsed")
360 allocate(collapseddata(&
363 allocate(collapseddata(&
366 call mpistop(
"slice direction not clear in allocate_collapsed")
370 allocate(collapseddata(nw+
nwauxio))
379 integer,
intent(in) :: igrid, jgrid, dir
383 integer :: ig^D, level, ig^Dtargetmin, ig^Dtargetmax
384 integer :: ix^Dtarget^LIM, idim^Dtarget^LIM
385 integer :: ixMdim^LLIM^D, ix^Dorig, ix^D
391 {^d& ig^d = tree%node%ig^d; }
392 level = tree%node%level
394 nx^d=ixmhi^d-ixmlo^d+1;
399 ig^dtargetmin = int(dble(ig^d-1)/2.0d0**(level-
collapselevel))+1
400 ig^dtargetmax = ig^dtargetmin
408 ix^dtargetmin = nx^d*(ig^dtargetmin-1)+1
409 ix^dtargetmax = nx^d*ig^dtargetmax
417 idim1target^lim=ix2target^lim;
418 idim2target^lim=ix3target^lim;
419 ixmdim^llim1=
ixm^llim2;
420 ixmdim^llim2=
ixm^llim3;
424 idim1target^lim=ix1target^lim;
425 idim2target^lim=ix3target^lim;
426 ixmdim^llim1=
ixm^llim1;
427 ixmdim^llim2=
ixm^llim3;
431 idim1target^lim=ix1target^lim;
432 idim2target^lim=ix2target^lim;
433 ixmdim^llim1=
ixm^llim1;
434 ixmdim^llim2=
ixm^llim2;
436 call mpistop(
"slice direction not clear in integrate_subnode")
440 do ix1orig = ixmdimlo1,ixmdimhi1
441 do ix2orig = ixmdimlo2,ixmdimhi2
443 collapseddata(ix1,ix2,1:nw+
nwauxio) = collapseddata(ix1,ix2,1:nw+
nwauxio) &
449 do ix^dm = idim^dmtargetmin,idim^dmtargetmax\}
451 collapseddata(ix1,ix2,1:nw+
nwauxio) = collapseddata(ix1,ix2,1:nw+
nwauxio) + ps_sub(jgrid)%w(ix1orig,ix2orig,1:nw+
nwauxio)
459 idim1target^lim=ix2target^lim;
460 ixmdim^llim1=
ixm^llim2;
463 idim1target^lim=ix1target^lim;
464 ixmdim^llim1=
ixm^llim1;
466 call mpistop(
"slice direction not clear in integrate_subnode")
470 do ix1orig = ixmdimlo1,ixmdimhi1
472 collapseddata(ix1,1:nw+
nwauxio) = collapseddata(ix1,1:nw+
nwauxio) &
477 do ix^dm = idim^dmtargetmin,idim^dmtargetmax\}
479 collapseddata(ix1,1:nw+
nwauxio) = collapseddata(ix1,1:nw+
nwauxio) + ps_sub(jgrid)%w(ix1orig,1:nw+
nwauxio)
495 integer,
intent(in) :: igrid, jgrid, dir
496 double precision,
dimension(0:nw+nwauxio),
intent(out) :: normconv
499 double precision :: dx^D
500 double precision,
dimension(ixMlo^D-1:ixMhi^D,ndim) :: xC_TMP, xC
501 double precision,
dimension(ixMlo^D:ixMhi^D,ndim) :: xCC_TMP, xCC
502 double precision,
dimension(ixMlo^D-1:ixMhi^D,nw+nwauxio) :: wC_TMP
503 double precision,
dimension(ixMlo^D:ixMhi^D,nw+nwauxio) :: wCC_TMP
504 integer :: ixC^L, ixCC^L
518 ps_sub(jgrid)%w(ixmlo2:ixmhi2,ixmlo3:ixmhi3,1:nw+
nwauxio) = &
519 ps_sub(jgrid)%w(ixmlo2:ixmhi2,ixmlo3:ixmhi3,1:nw+
nwauxio) &
520 + wcc_tmp(ix,ixmlo2:ixmhi2,ixmlo3:ixmhi3,1:nw+
nwauxio) * dx1
525 ps_sub(jgrid)%w(ixmlo2:ixmhi2,ixmlo3:ixmhi3,iw) = &
526 ps_sub(jgrid)%w(ixmlo2:ixmhi2,ixmlo3:ixmhi3,iw) &
527 + wcc_tmp(ix,ixmlo2:ixmhi2,ixmlo3:ixmhi3,iw) * &
528 ps(igrid)%dx(ix,ixmlo2:ixmhi2,ixmlo3:ixmhi3,1)
532 ps_sub(jgrid)%x(ixmlo2:ixmhi2,ixmlo3:ixmhi3,1:
ndim) = &
533 ps(igrid)%x(ixmlo1,ixmlo2:ixmhi2,ixmlo3:ixmhi3,1:
ndim)
537 ps_sub(jgrid)%w(ixmlo1:ixmhi1,ixmlo3:ixmhi3,1:nw+
nwauxio) = &
538 ps_sub(jgrid)%w(ixmlo1:ixmhi1,ixmlo3:ixmhi3,1:nw+
nwauxio) &
539 + wcc_tmp(ixmlo1:ixmhi1,ix,ixmlo3:ixmhi3,1:nw+
nwauxio) * dx2
544 ps_sub(jgrid)%w(ixmlo1:ixmhi1,ixmlo3:ixmhi3,iw) = &
545 ps_sub(jgrid)%w(ixmlo1:ixmhi1,ixmlo3:ixmhi3,iw) &
546 + wcc_tmp(ixmlo1:ixmhi1,ix,ixmlo3:ixmhi3,iw) * &
547 ps(igrid)%dx(ixmlo1:ixmhi1,ix,ixmlo3:ixmhi3,2)
551 ps_sub(jgrid)%x(ixmlo1:ixmhi1,ixmlo3:ixmhi3,1:
ndim) = &
552 ps(igrid)%x(ixmlo1:ixmhi1,ixmlo2,ixmlo3:ixmhi3,1:
ndim)
556 ps_sub(jgrid)%w(ixmlo1:ixmhi1,ixmlo2:ixmhi2,1:nw+
nwauxio) = &
557 ps_sub(jgrid)%w(ixmlo1:ixmhi1,ixmlo2:ixmhi2,1:nw+
nwauxio) &
558 + wcc_tmp(ixmlo1:ixmhi1,ixmlo2:ixmhi2,ix,1:nw+
nwauxio) * dx3
563 ps_sub(jgrid)%w(ixmlo1:ixmhi1,ixmlo2:ixmhi2,iw) = &
564 ps_sub(jgrid)%w(ixmlo1:ixmhi1,ixmlo2:ixmhi2,iw) &
565 + wcc_tmp(ixmlo1:ixmhi1,ixmlo2:ixmhi2,ix,iw) * &
566 ps(igrid)%dx(ixmlo1:ixmhi1,ixmlo2:ixmhi2,ix,3)
570 ps_sub(jgrid)%x(ixmlo1:ixmhi1,ixmlo2:ixmhi2,1:
ndim) = &
571 ps(igrid)%x(ixmlo1:ixmhi1,ixmlo2:ixmhi2,ixmlo3,1:
ndim)
573 print*,
'subnode, dir: ', dir
574 call mpistop(
"slice direction not clear in collapse_subnode")
582 ps_sub(jgrid)%w(ixmlo2:ixmhi2,1:nw+
nwauxio) = &
583 ps_sub(jgrid)%w(ixmlo2:ixmhi2,1:nw+
nwauxio) &
584 + wcc_tmp(ix,ixmlo2:ixmhi2,1:nw+
nwauxio) * dx1
589 ps_sub(jgrid)%w(ixmlo2:ixmhi2,iw) = &
590 ps_sub(jgrid)%w(ixmlo2:ixmhi2,iw) &
591 + wcc_tmp(ix,ixmlo2:ixmhi2,iw) * &
592 ps(igrid)%dx(ix,ixmlo2:ixmhi2,1)
596 ps_sub(jgrid)%x(ixmlo2:ixmhi2,1:
ndim) = &
597 ps(igrid)%x(ixmlo1,ixmlo2:ixmhi2,1:
ndim)
601 ps_sub(jgrid)%w(ixmlo1:ixmhi1,1:nw+
nwauxio) = &
602 ps_sub(jgrid)%w(ixmlo1:ixmhi1,1:nw+
nwauxio) &
603 + wcc_tmp(ixmlo1:ixmhi1,ix,1:nw+
nwauxio) * dx2
608 ps_sub(jgrid)%w(ixmlo1:ixmhi1,iw) = &
609 ps_sub(jgrid)%w(ixmlo1:ixmhi1,iw) &
610 + wcc_tmp(ixmlo1:ixmhi1,ix,iw) * &
611 ps(igrid)%dx(ixmlo1:ixmhi1,ix,2)
615 ps_sub(jgrid)%x(ixmlo1:ixmhi1,1:
ndim) = &
616 ps(igrid)%x(ixmlo1:ixmhi1,ixmlo2,1:
ndim)
618 call mpistop(
"slice direction not clear in collapse_subnode")
629 ps_sub(jgrid)%w(iw) = ps_sub(jgrid)%w(iw) + wcc_tmp(ix,iw) * ps(igrid)%dx(ix,1)
633 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.
integer, 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
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