2D AMRCLAW
outval.f
Go to the documentation of this file.
1 c
3 c =======================================================================
4  subroutine outval(val,nvar,mitot,mjtot,mptr,outgrd,naux,aux)
5 c =======================================================================
6 c
7  use amr_module
8  implicit double precision (a-h,o-z)
9 
10  dimension val(nvar,mitot,mjtot)
11  dimension aux(naux,mitot,mjtot)
12  logical outgrd
13 
14 
15 c ::::::::::::::::::::::OUTVAL :::::::::::::::::::::::::::::::
16 c print solution and aux. variables to output.
17 c if cell outside domain, don't print soln. value - nothing
18 c currently in ghost cells.
19 c ::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::::
20 c
21  if (.not. outgrd) go to 99
22  level = node(nestlevel,mptr)
23  hx = hxposs(level)
24  hy = hyposs(level)
25  cornx = rnode(cornxlo,mptr) - nghost*hx
26  corny = rnode(cornylo,mptr) - nghost*hy
27 c
28  do 25 j=nghost+1,mjtot-nghost
29  do 20 i=nghost+1,mitot-nghost
30 
31  x = cornx + hx*(dble(i)-.5d0)
32  y = corny + hy*(dble(j)-.5d0)
33  write(outunit,107) x,y,i,j,(val(ivar,i,j),ivar=1,nvar)
34  107 format(2hx=,f6.3,2hy=,f6.3,3h,i=,i3,3h,j=,i3,' a=',
35  * e25.15)
36 c * 5(e9.3,1x))
37  if (naux.gt.0) write(outunit,108) (aux(iaux,i,j),iaux=1,naux)
38  108 format(1x,'aux = ',7(e10.3,1x))
39 
40  20 continue
41  25 continue
42 
43  99 return
44  end
integer, parameter cornxlo
x-coordinate of the left border of this grid
Definition: amr_module.f90:143
real(kind=8), dimension(maxlv) hyposs
Definition: amr_module.f90:193
real(kind=8), dimension(maxlv) hxposs
Definition: amr_module.f90:193
integer, dimension(nsize, maxgr) node
Definition: amr_module.f90:198
integer, parameter nestlevel
AMR level of the grid.
Definition: amr_module.f90:44
real(kind=8), dimension(rsize, maxgr) rnode
Definition: amr_module.f90:193
integer, parameter outunit
Definition: amr_module.f90:290
subroutine outval(val, nvar, mitot, mjtot, mptr, outgrd, naux, aux)
print solution and aux.
Definition: outval.f:5
integer, parameter cornylo
y-coordinate of the lower border of this grid
Definition: amr_module.f90:145
integer nghost
Definition: amr_module.f90:232
The module contains the definition of a "node descriptor" as well as other global variables used duri...
Definition: amr_module.f90:21