def initialise_field(fs, name, fname, fpath='.', outputdir=None, op=CoupledOptions(), **kwargs): """ Initialise bathymetry field with results from a previous simulation. :arg fs: field will live in this finite element space. :arg name: name used internally for field. :arg fname: file name (without '.h5' extension). :kwarg fpath: directory to read the data from. :kwarg op: :class:`Options` parameter object. :kwarg index_str: optional five digit string. :kwarg delete: toggle deletion of the file. """ delete = kwargs.get('delete', False) index_str = kwargs.get('index_str', None) fname = get_filename(fname, index_str) with timed_stage('initialising {:s}'.format(name)): f = Function(fs, name=name) with DumbCheckpoint(os.path.join(fpath, fname), mode=FILE_READ) as chk: chk.load(f) # Optionally delete the HDF5 file if delete: os.remove(os.path.join(fpath, fname) + '.h5') # Plot to PVD if outputdir is not None and op.plot_pvd: File(os.path.join(outputdir, '_'.join([name, 'imported.pvd']))).write(f) return f
def test_bathymetry_io(family): op = CoupledOptions(bathymetry_family=family) # Create some bathymetry field mesh = get_mesh(2) x, y = SpatialCoordinate(mesh) P1 = get_function_space(mesh, 'scalar', family.upper(), 1) b = interpolate(x + 1, P1) # Save it to file plexname = 'myplex' index_str = '00000' fpath = create_directory('tmp') export_bathymetry(b, fpath, plexname=plexname, op=op, index_str=index_str) # Read from the file and check consistency # TODO: Make consistent by reading from file b_new = initialise_bathymetry(mesh, fpath, outputdir=fpath, op=op, index_str=index_str, delete=True) assert np.allclose(b.dat.data, b_new.dat.data) # Clean up os.remove(os.path.join(fpath, plexname + '.h5')) os.remove(os.path.join(fpath, 'bathymetry_out.pvd')) os.remove(os.path.join(fpath, 'bathymetry_out_0.vtu')) os.remove(os.path.join(fpath, 'bathymetry_imported.pvd')) os.remove(os.path.join(fpath, 'bathymetry_imported_0.vtu')) os.rmdir(fpath)
def export_field(f, name, fname, fpath='.', plexname='myplex', op=CoupledOptions(), index_str=None): """ Export some field to be used in a subsequent simulation. :arg f: field (Firedrake :class:`Function`) to be stored. :arg name: name used internally for field. :arg fname: filename to save the data to. :kwarg fpath: directory to save the data to. :kwarg plexname: file name to be used for the DMPlex data file. :kwarg op: :class:`Options` parameter object. :kwarg index_str: optional five digit string. """ if not os.path.exists(fpath): os.makedirs(fpath) op.print_debug( "I/O: Exporting {:s} for subsequent simulation".format(name)) # Create checkpoint to HDF5 fname = get_filename(fname, index_str) with DumbCheckpoint(os.path.join(fpath, fname), mode=FILE_CREATE) as chk: chk.store(f, name=name) # Plot to PVD if op.plot_pvd: File(os.path.join(fpath, '_'.join([name, 'out.pvd']))).write(f) # Save mesh to DMPlex format if plexname is not None: save_mesh(f.function_space().mesh(), plexname, fpath)
def initialise_bathymetry(mesh, fpath, op=CoupledOptions(), **kwargs): """ Initialise bathymetry field with results from a previous simulation. :arg mesh: field will be defined in finite element space on this mesh. :arg fpath: directory to read the data from. :kwarg op: :class:`Options` parameter object. """ # TODO: Would be nice to have consistency: here mesh is an arg but below it is read from file fs = FunctionSpace(mesh, op.bathymetry_family.upper(), 1) return initialise_field(fs, 'bathymetry', 'bathymetry', fpath, **kwargs)
def export_hydrodynamics(uv, elev, fpath='.', plexname='myplex', op=CoupledOptions(), **kwargs): """ Export velocity and elevation to be used in a subsequent simulation :arg uv: velocity field to be stored. :arg elev: elevation field to be stored. :kwarg fpath: directory to save the data to. :kwarg plexname: file name to be used for the DMPlex data file. :kwarg op: :class:`Options` parameter object. :kwarg index_str: optional five digit string. """ index_str = kwargs.get('index_str', None) if not os.path.exists(fpath): os.makedirs(fpath) op.print_debug("I/O: Exporting fields for subsequent simulation") # Check consistency of meshes mesh = elev.function_space().mesh() assert mesh == uv.function_space().mesh() # Export velocity name = "velocity" fname = get_filename(name, index_str) with DumbCheckpoint(os.path.join(fpath, fname), mode=FILE_CREATE) as chk: chk.store(uv, name=name) # Export elevation name = "elevation" fname = get_filename(name, index_str) with DumbCheckpoint(os.path.join(fpath, fname), mode=FILE_CREATE) as chk: chk.store(elev, name=name) # Plot to .pvd if op.plot_pvd: uv_proj = Function(VectorFunctionSpace(mesh, "CG", 1), name="Initial velocity") uv_proj.project(uv) File(os.path.join(fpath, 'velocity_out.pvd')).write(uv_proj) elev_proj = Function(FunctionSpace(mesh, "CG", 1), name="Initial elevation") elev_proj.project(elev) File(os.path.join(fpath, 'elevation_out.pvd')).write(elev_proj) # Export mesh if plexname is not None: save_mesh(mesh, plexname, fpath)
def test_hydro_io(pair): element_pair = { 'dg-dg': (('DG', 1), ('DG', 1)), 'dg-cg': (('DG', 1), ('CG', 2)), 'cg-cg': (('CG', 2), ('CG', 1)), }[pair] op = CoupledOptions(family=pair) # Create some hydrodynamics fields mesh = get_mesh(2) x, y = SpatialCoordinate(mesh) U = get_function_space(mesh, 'vector', *element_pair[0]) H = get_function_space(mesh, 'scalar', *element_pair[1]) uv = interpolate(as_vector([x, y]), U) elev = interpolate(x * y, H) # Save them to file fnames = ('velocity', 'elevation') plexname = 'myplex' index_str = '00000' fpath = create_directory('tmp') export_hydrodynamics(uv, elev, fpath, plexname=plexname, op=op, index_str=index_str) # Read from the file and check consistency uv_new, elev_new = initialise_hydrodynamics(fpath, outputdir=fpath, op=op, index_str=index_str, delete=True) assert np.allclose(uv.dat.data, uv_new.dat.data) assert np.allclose(elev.dat.data, elev_new.dat.data) # Clean up for fname in fnames: for extension in ('.pvd', '_0.vtu'): os.remove( os.path.join(fpath, '{:s}_out{:s}'.format(fname, extension))) os.remove( os.path.join(fpath, '{:s}_imported{:s}'.format(fname, extension))) os.rmdir(fpath)
def initialise_hydrodynamics(inputdir='.', outputdir=None, op=CoupledOptions(), **kwargs): """ Initialise velocity and elevation with results from a previous simulation. :kwarg inputdir: directory to read the data from. :kwarg outputdir: directory to optionally plot the data in .pvd format. :kwarg op: :class:`Options` parameter object. :kwarg delete: toggle deletion of the file. :kwarg plexname: file name used for the DMPlex data file. :kwarg variant: relates to distribution of quadrature nodes in an element. """ plexname = kwargs.get('plexname', 'myplex') variant = kwargs.get('variant', 'equispaced') delete = kwargs.get('delete', False) index_str = kwargs.get('index_str', None) # Get finite element if op.family == 'dg-dg': uv_element = VectorElement("DG", triangle, 1) if variant is None: elev_element = FiniteElement("DG", triangle, 1) else: elev_element = FiniteElement("DG", triangle, 1, variant=variant) elif op.family == 'dg-cg': uv_element = VectorElement("DG", triangle, 1) if variant is None: elev_element = FiniteElement("CG", triangle, 2) else: elev_element = FiniteElement("CG", triangle, 2, variant=variant) elif op.family == 'cg-cg': uv_element = VectorElement("CG", triangle, 2) if variant is None: elev_element = FiniteElement("CG", triangle, 1) else: elev_element = FiniteElement("CG", triangle, 1, variant=variant) # Load mesh with timed_stage('mesh'): mesh = op.default_mesh if plexname is None else load_mesh( plexname, inputdir, delete=delete) # Load velocity name = "velocity" fname = get_filename(name, index_str) U = FunctionSpace(mesh, uv_element) with timed_stage('initialising {:s}'.format(name)): with DumbCheckpoint(os.path.join(inputdir, fname), mode=FILE_READ) as chk: uv_init = Function(U, name=name) chk.load(uv_init) # Optionally delete the velocity HDF5 file if delete: os.remove(os.path.join(inputdir, fname) + '.h5') # Load elevation name = "elevation" fname = get_filename(name, index_str) H = FunctionSpace(mesh, elev_element) with timed_stage('initialising {:s}'.format(name)): with DumbCheckpoint(os.path.join(inputdir, fname), mode=FILE_READ) as chk: elev_init = Function(H, name=name) chk.load(elev_init) # Optionally delete the elevation HDF5 file if delete: os.remove(os.path.join(inputdir, fname) + '.h5') # Plot to .pvd if outputdir is not None and op.plot_pvd: uv_proj = Function(VectorFunctionSpace(mesh, "CG", 1), name="Initial velocity") uv_proj.project(uv_init) File(os.path.join(outputdir, "velocity_imported.pvd")).write(uv_proj) elev_proj = Function(FunctionSpace(mesh, "CG", 1), name="Initial elevation") elev_proj.project(elev_init) File(os.path.join(outputdir, "elevation_imported.pvd")).write(elev_proj) return uv_init, elev_init