MPI-AMRVAC 3.1
The MPI - Adaptive Mesh Refinement - Versatile Advection Code (development version)
Loading...
Searching...
No Matches
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
9contains
10!=============================================================================
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
18integer :: idir
19logical, save :: firstcollapse=.true.
20!-----------------------------------------------------------------------------
21
22if (firstcollapse) then
23 firstcollapse=.false.
24end if
25
26do idir=1, ndim
27 if (collapse(idir)) call put_collapse(idir)
28end do
29
31end subroutine write_collapsed
32!=============================================================================
33subroutine put_collapse(dir)
37integer, intent(in) :: dir
38! .. local ..
39double precision,dimension(0:nw+nwauxio) :: normconv
40integer :: jgrid, igrid, Morton_no
41!-----------------------------------------------------------------------------
42if(.not.slab) call mpistop("collapse only for slab cartesian cases")
43
44call allocate_collapsed(dir)
45
46jgrid=0
47do 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)
53end do
54
55! Reduce to head-node:
56if (mype==0) then
57 call mpi_reduce(mpi_in_place,collapseddata,size(collapseddata),mpi_double_precision,mpi_sum,0,icomm,ierrmpi)
58else
59 call mpi_reduce(collapseddata,collapseddata,size(collapseddata),mpi_double_precision,mpi_sum,0,icomm,ierrmpi)
60end if
61call mpi_barrier(icomm, ierrmpi)
62
63{^nooned
64select case(collapse_type)
65case('vti')
66 call output_collapsed_vti(dir,normconv)
67case('csv')
68 call output_collapsed_csv(dir,normconv)
69case('default')
70 call mpistop("Unknown filetype for collapsed views")
71end select
72 }{^ifoned
73call mpistop("sorry, 1D collapsed output routine not yet implemented (should be easy)...")
74}
75
76! If we need the subnodes later, remove deallocation here:
77do jgrid=1,morton_stop(mype)-morton_start(mype)+1
78 call dealloc_subnode(jgrid)
79end do
80deallocate(collapseddata)
81
82end subroutine put_collapse
83!=============================================================================
84subroutine output_collapsed_csv(dir,normconv)
87integer, intent(in) :: dir
88double precision,dimension(0:nw+nwauxio),intent(in):: normconv
89
90{^nooned
91double precision :: dxdim^DM, xdim^LIM^DM
92integer :: ix^DM
93}
94integer :: iw
95integer, dimension(ndim) :: myshape
96logical :: fileopen
97character(len=name_len) :: wnamei(1:nw+nwauxio),xandwnamei(1:ndim+nw+nwauxio)
98character(len=1024) :: filename, outfilehead, line
99!-----------------------------------------------------------------------------
100
101if (mype==0) then
102
103! Get coordinates:
104{^ifthreed
105select case(dir)
106case (1)
107 dxdim1 = dx(2,collapselevel)
108 dxdim2 = dx(3,collapselevel)
109 xdim^lim1=xprob^lim2;
110 xdim^lim2=xprob^lim3;
111case (2)
112 dxdim1 = dx(1,collapselevel)
113 dxdim2 = dx(3,collapselevel)
114 xdim^lim1=xprob^lim1;
115 xdim^lim2=xprob^lim3;
116case (3)
117 dxdim1 = dx(1,collapselevel)
118 dxdim2 = dx(2,collapselevel)
119 xdim^lim1=xprob^lim1;
120 xdim^lim2=xprob^lim2;
121case 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")
127end select
128}
129{^iftwod
130select case(dir)
131case (1)
132 dxdim1 = dx(2,collapselevel)
133 xdim^lim1=xprob^lim2;
134case (2)
135 dxdim1 = dx(1,collapselevel)
136 xdim^lim1=xprob^lim1;
137case default
138 dxdim1 = 1
139 xdim^lim1=xprob^lim2;
140 call mpistop("slice direction not clear in output_collapsed_csv")
141end 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}
169close(unitcollapse)
170
171end if
172end subroutine output_collapsed_csv
173!=============================================================================
174subroutine output_collapsed_vti(dir,normconv)
177integer, intent(in) :: dir
178double precision,dimension(0:nw+nwauxio),intent(in):: normconv
179
180{^nooned
181double precision :: dxdim^DM, xdim^LIM^DM
182integer :: ix^DM
183}
184double precision :: origin(1:3), spacing(1:3)
185integer*8 :: offset
186integer :: wholeExtent(1:6), size_single, length, size_length
187integer :: iw
188integer, dimension(ndim) :: myshape
189logical :: fileopen
190character(len=1024) :: filename, outfilehead, line
191character(len=name_len) :: wnamei(1:nw+nwauxio),xandwnamei(1:ndim+nw+nwauxio)
192character :: buf
193!-----------------------------------------------------------------------------
194if (mype==0) then
195offset=0
196size_single=4
197size_length=4
198! Get coordinates:
199{^ifthreed
200select case(dir)
201case (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;
206case (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;
211case (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;
216case 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")
222end select
223}
224{^iftwod
225select case(dir)
226case (1)
227 dxdim1 = dx(2,1)*2.0d0**(1-collapselevel)
228 xdim^lim1=xprob^lim2;
229case (2)
230 dxdim1 = dx(1,1)*2.0d0**(1-collapselevel)
231 xdim^lim1=xprob^lim1;
232case default
233 dxdim1 = 1
234 xdim^lim1=xprob^lim2;
235 call mpistop("slice direction not clear in output_collapsed_vti")
236end select
237}
238
239origin=0
240spacing=zero
241wholeextent=0
242myshape = shape(collapseddata)
243{^ifoned
244length = size_single
245}
246{^nooned
247length = ^dm&myshape(^dm)*
248length = 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:
262call getheadernames(wnamei,xandwnamei,outfilehead)
263
264! generate xml header
265write(unitcollapse,'(a)')'<?xml version="1.0"?>'
266write(unitcollapse,'(a)',advance='no') '<VTKFile type="ImageData"'
267if(type_endian==1)then
268 write(unitcollapse,'(a)')' version="0.1" byte_order="LittleEndian">'
269else
270 write(unitcollapse,'(a)')' version="0.1" byte_order="BigEndian">'
271endif
272! following corresponds to uniform cartesian grid
273write(unitcollapse,'(a,3(1pe14.6),a,6(i10),a,3(1pe14.6),a)')' <ImageData Origin="',&
274 origin,'" WholeExtent="',wholeextent,'" Spacing="',spacing,'">'
275write(unitcollapse,'(a)')'<FieldData>'
276write(unitcollapse,'(2a)')'<DataArray type="Float32" Name="TIME" ',&
277 'NumberOfTuples="1" format="ascii">'
279write(unitcollapse,'(a)')'</DataArray>'
280write(unitcollapse,'(a)')'</FieldData>'
281
282! we write one VTK PIECE
283write(unitcollapse,'(a,6(i10),a)') &
284 '<Piece Extent="',wholeextent,'">'
285write(unitcollapse,'(a)')'<CellData>'
286
287do 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
295enddo
296
297write(unitcollapse,'(a)')'</CellData>'
298write(unitcollapse,'(a)')'</Piece>'
299write(unitcollapse,'(a)')'</ImageData>'
300write(unitcollapse,'(a)')'<AppendedData encoding="raw">'
301close(unitcollapse)
302
303open(unitcollapse,file=filename,access='stream',form='unformatted',position='append')
304buf='_'
305write(unitcollapse) trim(buf)
306
307do 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 }
322enddo
323
324close(unitcollapse)
325open(unitcollapse,file=filename,status='unknown',form='formatted',position='append')
326write(unitcollapse,'(a)')'</AppendedData>'
327write(unitcollapse,'(a)')'</VTKFile>'
328close(unitcollapse)
329
330end if
331
332end subroutine output_collapsed_vti
333!=============================================================================
334subroutine allocate_collapsed(dir)
336integer, intent(in) :: dir
337integer :: dim^D
338integer :: nx^D
339!-----------------------------------------------------------------------------
340! allocate array for the collapsed data:
341! number of cells per grid.
342nx^d=ixmhi^d-ixmlo^d+1;
343{dim^d=ng^d(1)*2**(collapselevel-1)*nx^d\}
344{^ifthreed
345select case(dir)
346case (1)
347 allocate(collapseddata(&
348 dim2,dim3,nw+nwauxio))
349case (2)
350 allocate(collapseddata(&
351 dim1,dim3,nw+nwauxio))
352case (3)
353 allocate(collapseddata(&
354 dim1,dim2,nw+nwauxio))
355case default
356 call mpistop("slice direction not clear in allocate_collapsed")
357end select
358}
359{^iftwod
360select case(dir)
361case (1)
362 allocate(collapseddata(&
363 dim2,nw+nwauxio))
364case (2)
365 allocate(collapseddata(&
366 dim1,nw+nwauxio))
367case default
368 call mpistop("slice direction not clear in allocate_collapsed")
369end select
370}
371{^ifoned
372 allocate(collapseddata(nw+nwauxio))
373}
374collapseddata = zero
375
376end subroutine allocate_collapsed
377!=============================================================================
378subroutine integrate_subnode(igrid,jgrid,dir)
381integer, intent(in) :: igrid, jgrid, dir
382! .. local ..
383type(tree_node_ptr) :: tree
384integer :: nx^D
385integer :: ig^D, level, ig^Dtargetmin, ig^Dtargetmax
386integer :: ix^Dtarget^LIM, idim^Dtarget^LIM
387integer :: ixMdim^LLIM^D, ix^Dorig, ix^D
388{^nooned
389integer :: igdim^DM
390}
391!-----------------------------------------------------------------------------
392tree%node => igrid_to_node(igrid, mype)%node
393{^d& ig^d = tree%node%ig^d; }
394level = tree%node%level
395! number of cells per grid.
396nx^d=ixmhi^d-ixmlo^d+1;
397
398! assumes uniform cartesian grid with factor 2 refinement per dimension
399{^d&
400if (level > collapselevel) then
401 ig^dtargetmin = int(dble(ig^d-1)/2.0d0**(level-collapselevel))+1
402 ig^dtargetmax = ig^dtargetmin
403else 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
406else
407 ig^dtargetmin = ig^d
408 ig^dtargetmax = ig^d
409end if
410ix^dtargetmin = nx^d*(ig^dtargetmin-1)+1
411ix^dtargetmax = nx^d*ig^dtargetmax
412\}
413
414{^ifthreed
415select case(dir)
416case (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;
423case (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;
430case (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;
437case default
438 call mpistop("slice direction not clear in integrate_subnode")
439end select
440
441if (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
449else
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\}
455end if
456}
457{^iftwod
458select case(dir)
459case (1)
460 igdim1=(ig2-1)*nx2
461 idim1target^lim=ix2target^lim;
462 ixmdim^llim1=ixm^llim2;
463case (2)
464 igdim1=(ig1-1)*nx1
465 idim1target^lim=ix1target^lim;
466 ixmdim^llim1=ixm^llim1;
467case default
468 call mpistop("slice direction not clear in integrate_subnode")
469end select
470
471if (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
477else
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\}
483end if
484}
485{^ifoned
486collapseddata(1:nw+nwauxio) = collapseddata(1:nw+nwauxio) + ps_sub(jgrid)%w(1:nw+nwauxio)
487}
488
489end subroutine integrate_subnode
490!=============================================================================
491subroutine collapse_subnode(igrid,jgrid,dir,normconv)
496
497integer, intent(in) :: igrid, jgrid, dir
498double precision,dimension(0:nw+nwauxio),intent(out) :: normconv
499! .. local ..
500double precision :: dx^D
501double precision, dimension(ixMlo^D-1:ixMhi^D,ndim) :: xC_TMP, xC
502double precision, dimension(ixMlo^D:ixMhi^D,ndim) :: xCC_TMP, xCC
503double precision, dimension(ixMlo^D-1:ixMhi^D,nw+nwauxio) :: wC_TMP
504double precision, dimension(ixMlo^D:ixMhi^D,nw+nwauxio) :: wCC_TMP
505integer :: ixC^L, ixCC^L
506integer :: ix, iw
507!-----------------------------------------------------------------------------
508ps_sub(jgrid)%w=zero
509dx^d=rnode(rpdx^d_,igrid);
510
511call calc_x(igrid,xc,xcc)
512call calc_grid(unitslice,igrid,xc,xcc,xc_tmp,xcc_tmp,wc_tmp,wcc_tmp,normconv,&
513 ixc^l,ixcc^l,.true.)
514
515{^ifthreed
516select case (dir)
517case (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)
536case (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)
555case (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)
574case default
575 print*, 'subnode, dir: ', dir
576 call mpistop("slice direction not clear in collapse_subnode")
577end select
578}
579{^iftwod
580select case (dir)
581case (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)
600case (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)
619case default
620 call mpistop("slice direction not clear in collapse_subnode")
621end select
622}
623{^ifoned
624if(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
628else
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
634end if
635ps_sub(jgrid)%x(1:ndim) = ps(igrid)%x(ixmlo1,1:ndim)
636}
637
638end subroutine collapse_subnode
639!=============================================================================
640end module mod_collapse
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
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)
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.
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...
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:55
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