[yt-users] AMR grid creation in spherical coordinates
Kacper Kowalik
xarthisius.kk at gmail.com
Wed May 11 08:59:17 PDT 2016
On 05/11/2016 07:36 AM, Andrea Negri wrote:
> Hi all,
>
> I'm trying to build an AMR in yt in spherical coordinates. I was able
> to it in the following code. But I have a couple of questions.
>
> There is something I don't understand: except for the level 0 grid,
> the right_edge arrays seems to be useless...
> In addition, as you can notice the number of grids in theta direction
> increased as a factor of 5. Is there a way to preserve the number of
> grid in theta (and in phi) from level 0 to level 1?
Hi Andrea,
unfortunately "load_amr_grids" has the same deficiency as GDF format, it
won't handle non uniform mesh refinement.
However, I was able to read you data as a semi-structured mesh:
import yt
import h5py as h5
import numpy as np
dens = None
xgrid = np.empty((0,))
ygrid = np.empty((0,))
zgrid = np.empty((0,))
with h5.File('amr_zeus0.h5', 'r') as h5f:
for grid in h5f["data"]:
LE = h5f["data"][grid].attrs["left_edge"][:]
RE = h5f["data"][grid].attrs["right_edge"][:]
nb = h5f["data"][grid].attrs["n_b"][:]
x = np.linspace(LE[0], RE[0], nb[0], endpoint=False)
y = np.linspace(LE[1], RE[1], nb[1], endpoint=False)
z = np.linspace(LE[2], RE[2], nb[2], endpoint=False)
temp = h5f["data"][grid]["density"][:]
if dens is None:
dens = temp
else:
dens = np.concatenate((dens, temp))
xgrid = np.concatenate((xgrid, x))
ygrid = y # They're always the same so it doesn't matter which
zgrid = z # one we pick
DRE = h5f["simulation_parameters"].attrs["domain_right_edge"][:]
xgrid = np.concatenate((xgrid, [DRE[0]]))
ygrid = np.concatenate((ygrid, [DRE[1]]))
zgrid = np.concatenate((zgrid, [DRE[2]]))
xgrid /= yt.YTQuantity(1, 'kpc').in_cgs().d
coordinates, connectivity = yt.hexahedral_connectivity(xgrid, ygrid, zgrid)
bbox = np.array([[np.min(xgrid), np.max(xgrid)],
[np.min(ygrid), np.max(ygrid)],
[np.min(zgrid), np.max(zgrid)]])
data = {"density" : dens}
ds = yt.load_hexahedral_mesh(data, connectivity, coordinates, 1.0,
bbox=bbox)
slc = yt.SlicePlot(ds, 'y', "density")
slc.save()
It currently uses gdf file created using your conversion script as I'm
more familiar with where to look for things like grid coordinates, raw
data etc. However, I'm pretty sure you can read the required data
straight from zeus hdf5 files.
Cheers,
Kacper
>
>
>
> import numpy as np
> import yt
>
> refine_by=5
> top_grid_dim = [100, 10, 2]
> n1=100
> n2=10
> n3=2
>
> grid_data = [
> dict(left_edge = [0.0, 0.0, 0.0],
> right_edge = [1.0, np.pi, np.pi*2.],
> level = 0,
> dimensions = [n1, n2, n3]),
> dict(left_edge = [0., 0., 0.],
> right_edge = [0.5, np.pi, np.pi*2.],
> level = 1,
> dimensions = [n1, n2*refine_by, n3*refine_by]),
> ]
>
>
> for g in grid_data:
> g["density"] = (np.random.random(g["dimensions"]), "g/cm**3")
>
> ds = yt.load_amr_grids(grid_data, top_grid_dim ,
> geometry='spherical', refine_by=refine_by, length_unit='kpc')
>
> ds.print_stats()
>
> yt.SlicePlot(ds, 'phi', 'density').save()
> yt.SlicePlot(ds, 'theta', 'density').save()
> yt.SlicePlot(ds, 'r', 'density').save()
> _______________________________________________
> yt-users mailing list
> yt-users at lists.spacepope.org
> http://lists.spacepope.org/listinfo.cgi/yt-users-spacepope.org
>
-------------- next part --------------
A non-text attachment was scrubbed...
Name: signature.asc
Type: application/pgp-signature
Size: 819 bytes
Desc: OpenPGP digital signature
URL: <http://lists.spacepope.org/pipermail/yt-users-spacepope.org/attachments/20160511/e9cb7a55/attachment.sig>
More information about the yt-users
mailing list