MPI-AMRVAC  3.1
The MPI - Adaptive Mesh Refinement - Versatile Advection Code (development version)
mod_collapse.t
Go to the documentation of this file.
1 !> Collapses D-dimensional output to D-1 view by line-of-sight integration
2 ! for now: only for non-stretched cartesian cases, with LOS along coordinate
3 ! writes out as either csv or vti file (collapse_type)
6  use mod_comm_lib, only: mpistop
7  implicit none
8 
9 contains
10 !=============================================================================
11 subroutine write_collapsed
13 ! Writes a collapsed view of the data integrated over one grid-direction.
14 ! E.g. column density maps.
15 ! Uses flat interpolation throughout.
16 ! by Oliver Porth
17 ! 6.Nov 2013
18 integer :: idir
19 logical, save :: firstcollapse=.true.
20 !-----------------------------------------------------------------------------
21 
22 if (firstcollapse) then
23  firstcollapse=.false.
24 end if
25 
26 do idir=1, ndim
27  if (collapse(idir)) call put_collapse(idir)
28 end do
29 
31 end subroutine write_collapsed
32 !=============================================================================
33 subroutine put_collapse(dir)
37 integer, intent(in) :: dir
38 ! .. local ..
39 double precision,dimension(0:nw+nwauxio) :: normconv
40 integer :: jgrid, igrid, Morton_no
41 !-----------------------------------------------------------------------------
42 if(.not.slab) call mpistop("collapse only for slab cartesian cases")
43 
44 call allocate_collapsed(dir)
45 
46 jgrid=0
47 do morton_no=morton_start(mype),morton_stop(mype)
48  igrid=sfc_to_igrid(morton_no)
49  jgrid=jgrid+1
50  call alloc_subnode(jgrid,dir,nwauxio)
51  call collapse_subnode(igrid,jgrid,dir,normconv)
52  call integrate_subnode(igrid,jgrid,dir)
53 end do
54 
55 ! Reduce to head-node:
56 if (mype==0) then
57  call mpi_reduce(mpi_in_place,collapseddata,size(collapseddata),mpi_double_precision,mpi_sum,0,icomm,ierrmpi)
58 else
59  call mpi_reduce(collapseddata,collapseddata,size(collapseddata),mpi_double_precision,mpi_sum,0,icomm,ierrmpi)
60 end if
61 call mpi_barrier(icomm, ierrmpi)
62 
63 {^nooned
64 select case(collapse_type)
65 case('vti')
66  call output_collapsed_vti(dir,normconv)
67 case('csv')
68  call output_collapsed_csv(dir,normconv)
69 case('default')
70  call mpistop("Unknown filetype for collapsed views")
71 end select
72  }{^ifoned
73 call mpistop("sorry, 1D collapsed output routine not yet implemented (should be easy)...")
74 }
75 
76 ! If we need the subnodes later, remove deallocation here:
77 do jgrid=1,morton_stop(mype)-morton_start(mype)+1
78  call dealloc_subnode(jgrid)
79 end do
80 deallocate(collapseddata)
81 
82 end subroutine put_collapse
83 !=============================================================================
84 subroutine output_collapsed_csv(dir,normconv)
87 integer, intent(in) :: dir
88 double precision,dimension(0:nw+nwauxio),intent(in):: normconv
89 
90 {^nooned
91 double precision :: dxdim^DM, xdim^LIM^DM
92 integer :: ix^DM
93 }
94 integer :: iw
95 integer, dimension(ndim) :: myshape
96 logical :: fileopen
97 character(len=name_len) :: wnamei(1:nw+nwauxio),xandwnamei(1:ndim+nw+nwauxio)
98 character(len=1024) :: filename, outfilehead, line
99 !-----------------------------------------------------------------------------
100 
101 if (mype==0) then
102 
103 ! Get coordinates:
104 {^ifthreed
105 select case(dir)
106 case (1)
107  dxdim1 = dx(2,collapselevel)
108  dxdim2 = dx(3,collapselevel)
109  xdim^lim1=xprob^lim2;
110  xdim^lim2=xprob^lim3;
111 case (2)
112  dxdim1 = dx(1,collapselevel)
113  dxdim2 = dx(3,collapselevel)
114  xdim^lim1=xprob^lim1;
115  xdim^lim2=xprob^lim3;
116 case (3)
117  dxdim1 = dx(1,collapselevel)
118  dxdim2 = dx(2,collapselevel)
119  xdim^lim1=xprob^lim1;
120  xdim^lim2=xprob^lim2;
121 case default
122  dxdim1 = 1
123  dxdim2 = 1
124  xdim^lim1=xprob^lim2;
125  xdim^lim2=xprob^lim3;
126  call mpistop("slice direction not clear in output_collapsed_csv")
127 end select
128 }
129 {^iftwod
130 select case(dir)
131 case (1)
132  dxdim1 = dx(2,collapselevel)
133  xdim^lim1=xprob^lim2;
134 case (2)
135  dxdim1 = dx(1,collapselevel)
136  xdim^lim1=xprob^lim1;
137 case default
138  dxdim1 = 1
139  xdim^lim1=xprob^lim2;
140  call mpistop("slice direction not clear in output_collapsed_csv")
141 end select
142 }
143 
144  inquire(unitcollapse,opened=fileopen)
145  if(.not.fileopen)then
146  ! generate filename:
147  write(filename,"(a,i1.1,a,i1.1,a,i4.4,a)") trim(base_filename) // '_d', &
148  dir,'_l',collapselevel,'_n',collapsenext,'.csv'
149  open(unitcollapse,file=filename,status='unknown',form='formatted')
150  end if
151  ! get and write the header:
152  call getheadernames(wnamei,xandwnamei,outfilehead)
153  line=''
154  do iw=1,ndim+nw+nwauxio-1
155  if (iw .eq. dir) cycle
156  line = trim(line)//trim(xandwnamei(iw))//', '
157  end do
158  line = trim(line)//trim(xandwnamei(ndim+nw+nwauxio))
159  write(unitcollapse,'(a)')trim(line)
160  myshape = shape(collapseddata)
161 {^nooned
162 {^dm& do ix^dmb = 1,myshape(^dmb)\}
163  ! following assumes uniform cartesian grid
164  write(unitcollapse,'(200(1pe20.12))') &
165  {^dm& dxdim^dm*dble(ix^dm)+xdimmin^dm}, &
166  (collapseddata(ix^dm,iw)*normconv(iw),iw=1,nw+nwauxio)
167 {^dm& enddo\}
168 }
169 close(unitcollapse)
170 
171 end if
172 end subroutine output_collapsed_csv
173 !=============================================================================
174 subroutine output_collapsed_vti(dir,normconv)
177 integer, intent(in) :: dir
178 double precision,dimension(0:nw+nwauxio),intent(in):: normconv
179 
180 {^nooned
181 double precision :: dxdim^DM, xdim^LIM^DM
182 integer :: ix^DM
183 }
184 double precision :: origin(1:3), spacing(1:3)
185 integer*8 :: offset
186 integer :: wholeExtent(1:6), size_single, length, size_length
187 integer :: iw
188 integer, dimension(ndim) :: myshape
189 logical :: fileopen
190 character(len=1024) :: filename, outfilehead, line
191 character(len=name_len) :: wnamei(1:nw+nwauxio),xandwnamei(1:ndim+nw+nwauxio)
192 character :: buf
193 !-----------------------------------------------------------------------------
194 if (mype==0) then
195 offset=0
196 size_single=4
197 size_length=4
198 ! Get coordinates:
199 {^ifthreed
200 select case(dir)
201 case (1)
202  dxdim1 = dx(2,1)*2.0d0**(1-collapselevel)
203  dxdim2 = dx(3,1)*2.0d0**(1-collapselevel)
204  xdim^lim1=xprob^lim2;
205  xdim^lim2=xprob^lim3;
206 case (2)
207  dxdim1 = dx(1,1)*2.0d0**(1-collapselevel)
208  dxdim2 = dx(3,1)*2.0d0**(1-collapselevel)
209  xdim^lim1=xprob^lim1;
210  xdim^lim2=xprob^lim3;
211 case (3)
212  dxdim1 = dx(1,1)*2.0d0**(1-collapselevel)
213  dxdim2 = dx(2,1)*2.0d0**(1-collapselevel)
214  xdim^lim1=xprob^lim1;
215  xdim^lim2=xprob^lim2;
216 case default
217  dxdim1 = 1
218  dxdim2 = 1
219  xdim^lim1=xprob^lim2;
220  xdim^lim2=xprob^lim3;
221  call mpistop("slice direction not clear in output_collapsed_vti")
222 end select
223 }
224 {^iftwod
225 select case(dir)
226 case (1)
227  dxdim1 = dx(2,1)*2.0d0**(1-collapselevel)
228  xdim^lim1=xprob^lim2;
229 case (2)
230  dxdim1 = dx(1,1)*2.0d0**(1-collapselevel)
231  xdim^lim1=xprob^lim1;
232 case default
233  dxdim1 = 1
234  xdim^lim1=xprob^lim2;
235  call mpistop("slice direction not clear in output_collapsed_vti")
236 end select
237 }
238 
239 origin=0
240 spacing=zero
241 wholeextent=0
242 myshape = shape(collapseddata)
243 {^ifoned
244 length = size_single
245 }
246 {^nooned
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; }
252 }
253 
254  inquire(unitcollapse,opened=fileopen)
255  if(.not.fileopen)then
256  ! generate filename:
257  write(filename,"(a,i1.1,a,i1.1,a,i4.4,a)") trim(base_filename)//'_d',dir,&
258  '_l',collapselevel,'_n',collapsenext,'.vti'
259  open(unitcollapse,file=filename,status='unknown',form='formatted')
260  end if
261 ! get the header:
262 call getheadernames(wnamei,xandwnamei,outfilehead)
263 
264 ! generate xml header
265 write(unitcollapse,'(a)')'<?xml version="1.0"?>'
266 write(unitcollapse,'(a)',advance='no') '<VTKFile type="ImageData"'
267 if(type_endian==1)then
268  write(unitcollapse,'(a)')' version="0.1" byte_order="LittleEndian">'
269 else
270  write(unitcollapse,'(a)')' version="0.1" byte_order="BigEndian">'
271 endif
272 ! following corresponds to uniform cartesian grid
273 write(unitcollapse,'(a,3(1pe14.6),a,6(i10),a,3(1pe14.6),a)')' <ImageData Origin="',&
274  origin,'" WholeExtent="',wholeextent,'" Spacing="',spacing,'">'
275 write(unitcollapse,'(a)')'<FieldData>'
276 write(unitcollapse,'(2a)')'<DataArray type="Float32" Name="TIME" ',&
277  'NumberOfTuples="1" format="ascii">'
279 write(unitcollapse,'(a)')'</DataArray>'
280 write(unitcollapse,'(a)')'</FieldData>'
281 
282 ! we write one VTK PIECE
283 write(unitcollapse,'(a,6(i10),a)') &
284  '<Piece Extent="',wholeextent,'">'
285 write(unitcollapse,'(a)')'<CellData>'
286 
287 do iw=1,nw+nwauxio
288  if(iw<=nw) then
289  if(.not.w_write(iw)) cycle
290  endif
291  write(unitcollapse,'(a,a,a,i16,a)')&
292  '<DataArray type="Float32" Name="',trim(wnamei(iw)),&
293  '" format="appended" offset="',offset,'"/>'
294  offset = offset + length + size_length
295 enddo
296 
297 write(unitcollapse,'(a)')'</CellData>'
298 write(unitcollapse,'(a)')'</Piece>'
299 write(unitcollapse,'(a)')'</ImageData>'
300 write(unitcollapse,'(a)')'<AppendedData encoding="raw">'
301 close(unitcollapse)
302 
303 open(unitcollapse,file=filename,access='stream',form='unformatted',position='append')
304 buf='_'
305 write(unitcollapse) trim(buf)
306 
307 do iw=1,nw+nwauxio
308  if(iw<=nw) then
309  if(.not.w_write(iw)) cycle
310  endif
311  write(unitcollapse) length
312  {^ifoned
313  write(unitcollapse) real(collapseddata(iw)*normconv(iw))
314  }
315  {^iftwod
316  write(unitcollapse) (real(collapseddata(ix1,iw)*normconv(iw)),ix1=1,myshape(1))
317  }
318  {^ifthreed
319  write(unitcollapse) ((real(collapseddata(ix1,ix2,iw)*normconv(iw)),ix1=1,&
320  myshape(1)),ix2=1,myshape(2))
321  }
322 enddo
323 
324 close(unitcollapse)
325 open(unitcollapse,file=filename,status='unknown',form='formatted',position='append')
326 write(unitcollapse,'(a)')'</AppendedData>'
327 write(unitcollapse,'(a)')'</VTKFile>'
328 close(unitcollapse)
329 
330 end if
331 
332 end subroutine output_collapsed_vti
333 !=============================================================================
334 subroutine allocate_collapsed(dir)
336 integer, intent(in) :: dir
337 integer :: dim^D
338 integer :: nx^D
339 !-----------------------------------------------------------------------------
340 ! allocate array for the collapsed data:
341 ! number of cells per grid.
342 nx^d=ixmhi^d-ixmlo^d+1;
343 {dim^d=ng^d(1)*2**(collapselevel-1)*nx^d\}
344 {^ifthreed
345 select case(dir)
346 case (1)
347  allocate(collapseddata(&
348  dim2,dim3,nw+nwauxio))
349 case (2)
350  allocate(collapseddata(&
351  dim1,dim3,nw+nwauxio))
352 case (3)
353  allocate(collapseddata(&
354  dim1,dim2,nw+nwauxio))
355 case default
356  call mpistop("slice direction not clear in allocate_collapsed")
357 end select
358 }
359 {^iftwod
360 select case(dir)
361 case (1)
362  allocate(collapseddata(&
363  dim2,nw+nwauxio))
364 case (2)
365  allocate(collapseddata(&
366  dim1,nw+nwauxio))
367 case default
368  call mpistop("slice direction not clear in allocate_collapsed")
369 end select
370 }
371 {^ifoned
372  allocate(collapseddata(nw+nwauxio))
373 }
374 collapseddata = zero
375 
376 end subroutine allocate_collapsed
377 !=============================================================================
378 subroutine integrate_subnode(igrid,jgrid,dir)
381 integer, intent(in) :: igrid, jgrid, dir
382 ! .. local ..
383 type(tree_node_ptr) :: tree
384 integer :: nx^D
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
388 {^nooned
389 integer :: igdim^DM
390 }
391 !-----------------------------------------------------------------------------
392 tree%node => igrid_to_node(igrid, mype)%node
393 {^d& ig^d = tree%node%ig^d; }
394 level = tree%node%level
395 ! number of cells per grid.
396 nx^d=ixmhi^d-ixmlo^d+1;
397 
398 ! assumes uniform cartesian grid with factor 2 refinement per dimension
399 {^d&
400 if (level > collapselevel) then
401  ig^dtargetmin = int(dble(ig^d-1)/2.0d0**(level-collapselevel))+1
402  ig^dtargetmax = ig^dtargetmin
403 else if (level < collapselevel) then
404  ig^dtargetmin = int(2.0d0**(collapselevel-level))*(ig^d-1)+1
405  ig^dtargetmax = int(2.0d0**(collapselevel-level))*ig^d
406 else
407  ig^dtargetmin = ig^d
408  ig^dtargetmax = ig^d
409 end if
410 ix^dtargetmin = nx^d*(ig^dtargetmin-1)+1
411 ix^dtargetmax = nx^d*ig^dtargetmax
412 \}
413 
414 {^ifthreed
415 select case(dir)
416 case (1)
417  igdim1=(ig2-1)*nx2
418  igdim2=(ig3-1)*nx3
419  idim1target^lim=ix2target^lim;
420  idim2target^lim=ix3target^lim;
421  ixmdim^llim1=ixm^llim2;
422  ixmdim^llim2=ixm^llim3;
423 case (2)
424  igdim1=(ig1-1)*nx1
425  igdim2=(ig3-1)*nx3
426  idim1target^lim=ix1target^lim;
427  idim2target^lim=ix3target^lim;
428  ixmdim^llim1=ixm^llim1;
429  ixmdim^llim2=ixm^llim3;
430 case (3)
431  igdim1=(ig1-1)*nx1
432  igdim2=(ig2-1)*nx2
433  idim1target^lim=ix1target^lim;
434  idim2target^lim=ix2target^lim;
435  ixmdim^llim1=ixm^llim1;
436  ixmdim^llim2=ixm^llim2;
437 case default
438  call mpistop("slice direction not clear in integrate_subnode")
439 end select
440 
441 if (level >= collapselevel) then
442  do ix1orig = ixmdimlo1,ixmdimhi1
443  do ix2orig = ixmdimlo2,ixmdimhi2
444 {^dm& ix^dm = int(dble(ix^dmorig-nghostcells+igdim^dm-1)*2.0d0**(collapselevel-level))+1\}
445  collapseddata(ix1,ix2,1:nw+nwauxio) = collapseddata(ix1,ix2,1:nw+nwauxio) &
446  + ps_sub(jgrid)%w(ix1orig,ix2orig,1:nw+nwauxio) / 2.0d0**(2*(level-collapselevel))
447  end do
448  end do
449 else
450 {^dm&
451  do ix^dm = idim^dmtargetmin,idim^dmtargetmax\}
452  {^dm& ix^dmorig = int(dble(ix^dm-idim^dmtargetmin)/2.0d0**(collapselevel-level))+1+nghostcells \}
453  collapseddata(ix1,ix2,1:nw+nwauxio) = collapseddata(ix1,ix2,1:nw+nwauxio) + ps_sub(jgrid)%w(ix1orig,ix2orig,1:nw+nwauxio)
454  {^dm& enddo\}
455 end if
456 }
457 {^iftwod
458 select case(dir)
459 case (1)
460  igdim1=(ig2-1)*nx2
461  idim1target^lim=ix2target^lim;
462  ixmdim^llim1=ixm^llim2;
463 case (2)
464  igdim1=(ig1-1)*nx1
465  idim1target^lim=ix1target^lim;
466  ixmdim^llim1=ixm^llim1;
467 case default
468  call mpistop("slice direction not clear in integrate_subnode")
469 end select
470 
471 if (level >= collapselevel) then
472  do ix1orig = ixmdimlo1,ixmdimhi1
473 {^dm& ix^dm = int(dble(ix^dmorig-nghostcells+igdim^dm-1)*2.0d0**(collapselevel-level))+1\}
474  collapseddata(ix1,1:nw+nwauxio) = collapseddata(ix1,1:nw+nwauxio) &
475  + ps_sub(jgrid)%w(ix1orig,1:nw+nwauxio) / 2.0d0**(level-collapselevel)
476  end do
477 else
478 {^dm&
479  do ix^dm = idim^dmtargetmin,idim^dmtargetmax\}
480  {^dm& ix^dmorig = int(dble(ix^dm-idim^dmtargetmin)/2.0d0**(collapselevel-level))+1+nghostcells \}
481  collapseddata(ix1,1:nw+nwauxio) = collapseddata(ix1,1:nw+nwauxio) + ps_sub(jgrid)%w(ix1orig,1:nw+nwauxio)
482  {^dm& enddo\}
483 end if
484 }
485 {^ifoned
486 collapseddata(1:nw+nwauxio) = collapseddata(1:nw+nwauxio) + ps_sub(jgrid)%w(1:nw+nwauxio)
487 }
488 
489 end subroutine integrate_subnode
490 !=============================================================================
491 subroutine collapse_subnode(igrid,jgrid,dir,normconv)
494  use mod_physics, only: phys_to_primitive
495  use mod_calculate_xw
496 
497 integer, intent(in) :: igrid, jgrid, dir
498 double precision,dimension(0:nw+nwauxio),intent(out) :: normconv
499 ! .. local ..
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
506 integer :: ix, iw
507 !-----------------------------------------------------------------------------
508 ps_sub(jgrid)%w=zero
509 dx^d=rnode(rpdx^d_,igrid);
510 
511 call calc_x(igrid,xc,xcc)
512 call calc_grid(unitslice,igrid,xc,xcc,xc_tmp,xcc_tmp,wc_tmp,wcc_tmp,normconv,&
513  ixc^l,ixcc^l,.true.)
514 
515 {^ifthreed
516 select case (dir)
517 case (1)
518  if(slab_uniform) then
519  do ix=ixmlo1,ixmhi1
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
523  end do
524  else
525  do iw=1,nw+nwauxio
526  do ix=ixmlo1,ixmhi1
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)
531  end do
532  end do
533  end if
534  ps_sub(jgrid)%x(ixmlo2:ixmhi2,ixmlo3:ixmhi3,1:ndim) = &
535  ps(igrid)%x(ixmlo1,ixmlo2:ixmhi2,ixmlo3:ixmhi3,1:ndim)
536 case (2)
537  if(slab_uniform) then
538  do ix=ixmlo2,ixmhi2
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
542  end do
543  else
544  do iw=1,nw+nwauxio
545  do ix=ixmlo2,ixmhi2
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)
550  end do
551  end do
552  endif
553  ps_sub(jgrid)%x(ixmlo1:ixmhi1,ixmlo3:ixmhi3,1:ndim) = &
554  ps(igrid)%x(ixmlo1:ixmhi1,ixmlo2,ixmlo3:ixmhi3,1:ndim)
555 case (3)
556  if(slab_uniform) then
557  do ix=ixmlo3,ixmhi3
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
561  end do
562  else
563  do iw=1,nw+nwauxio
564  do ix=ixmlo3,ixmhi3
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)
569  end do
570  end do
571  endif
572  ps_sub(jgrid)%x(ixmlo1:ixmhi1,ixmlo2:ixmhi2,1:ndim) = &
573  ps(igrid)%x(ixmlo1:ixmhi1,ixmlo2:ixmhi2,ixmlo3,1:ndim)
574 case default
575  print*, 'subnode, dir: ', dir
576  call mpistop("slice direction not clear in collapse_subnode")
577 end select
578 }
579 {^iftwod
580 select case (dir)
581 case (1)
582  if(slab_uniform) then
583  do ix=ixmlo1,ixmhi1
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
587  end do
588  else
589  do iw=1,nw+nwauxio
590  do ix=ixmlo1,ixmhi1
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)
595  end do
596  end do
597  end if
598  ps_sub(jgrid)%x(ixmlo2:ixmhi2,1:ndim) = &
599  ps(igrid)%x(ixmlo1,ixmlo2:ixmhi2,1:ndim)
600 case (2)
601  if(slab_uniform) then
602  do ix=ixmlo2,ixmhi2
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
606  end do
607  else
608  do iw=1,nw+nwauxio
609  do ix=ixmlo2,ixmhi2
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)
614  end do
615  end do
616  end if
617  ps_sub(jgrid)%x(ixmlo1:ixmhi1,1:ndim) = &
618  ps(igrid)%x(ixmlo1:ixmhi1,ixmlo2,1:ndim)
619 case default
620  call mpistop("slice direction not clear in collapse_subnode")
621 end select
622 }
623 {^ifoned
624 if(slab_uniform) then
625  do ix=ixmlo1,ixmhi1
626  ps_sub(jgrid)%w(1:nw+nwauxio) = ps_sub(jgrid)%w(1:nw+nwauxio) + wcc_tmp(ix,1:nw+nwauxio) * dx1
627  end do
628 else
629  do iw=1,nw+nwauxio
630  do ix=ixmlo1,ixmhi1
631  ps_sub(jgrid)%w(iw) = ps_sub(jgrid)%w(iw) + wcc_tmp(ix,iw) * ps(igrid)%dx(ix,1)
632  end do
633  end do
634 end if
635 ps_sub(jgrid)%x(1:ndim) = ps(igrid)%x(ixmlo1,1:ndim)
636 }
637 
638 end subroutine collapse_subnode
639 !=============================================================================
640 end module mod_collapse
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.
Definition: mod_collapse.t:4
subroutine collapse_subnode(igrid, jgrid, dir, normconv)
Definition: mod_collapse.t:492
subroutine output_collapsed_csv(dir, normconv)
Definition: mod_collapse.t:85
subroutine allocate_collapsed(dir)
Definition: mod_collapse.t:335
subroutine put_collapse(dir)
Definition: mod_collapse.t:34
subroutine output_collapsed_vti(dir, normconv)
Definition: mod_collapse.t:175
subroutine write_collapsed
Definition: mod_collapse.t:12
subroutine integrate_subnode(igrid, jgrid, dir)
Definition: mod_collapse.t:379
subroutine, public mpistop(message)
Exit MPI-AMRVAC with an error message.
Definition: mod_comm_lib.t:208
Module with basic grid data structures.
Definition: mod_forest.t:2
integer, dimension(:), allocatable, save sfc_to_igrid
Go from a Morton number to an igrid index (for a single processor)
Definition: mod_forest.t:53
integer, dimension(:), allocatable, save morton_start
First Morton number per processor.
Definition: mod_forest.t:62
integer, dimension(:), allocatable, save morton_stop
Last Morton number per processor.
Definition: mod_forest.t:65
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 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...
Definition: mod_physics.t:4
procedure(sub_convert), pointer phys_to_primitive
Definition: mod_physics.t:59
Writes D-1 slice, can do so in various formats, depending on slice_type.
Definition: mod_slice.t:2
subroutine alloc_subnode(jgrid, dir, nwexpand)
Definition: mod_slice.t:994
subroutine dealloc_subnode(jgrid)
Definition: mod_slice.t:1042
Module with all the methods that users can customize in AMRVAC.
procedure(aux_output), pointer usr_aux_output
Pointer to a tree_node.
Definition: mod_forest.t:6