Changeset 9062


Ignore:
Timestamp:
Mar 14, 2014, 11:31:16 AM (11 years ago)
Author:
davies
Message:

Allowing plot_utils.Make_Geotif to take an xyz array

File:
1 edited

Legend:

Unmodified
Added
Removed
  • trunk/anuga_core/source/anuga/utilities/plot_utils.py

    r9052 r9062  
    583583##################################################################################
    584584
    585 def Make_Geotif(swwFile, output_quantities=['depth'],
     585def Make_Geotif(swwFile=None,
     586             output_quantities=['depth'],
    586587             myTimeStep=1, CellSize=5.0,
    587588             lower_left=None, upper_right=None,
     
    594595             verbose=False):
    595596    """
    596         Make a georeferenced tif by nearest-neighbour interpolation of sww file outputs.
     597        Make a georeferenced tif by nearest-neighbour interpolation of sww file outputs (or a 3-column array with xyz Points)
    597598
    598599        You must supply projection information as either a proj4string or an integer EPSG_CODE (but not both!)
    599600
    600         INPUTS: swwFile -- name of sww file
     601        INPUTS: swwFile -- name of sww file, OR a 3-column array with x/y/z
     602                    points. In the latter case x and y are assumed to be in georeferenced
     603                    coordinates.  The output raster will contain 'z', and will have a name-tag
     604                    based on the name in 'output_quantities'.
    601605                output_quantities -- list of quantitiies to plot, e.g. ['depth', 'velocity', 'stage','elevation','depthIntegratedVelocity']
    602606                myTimeStep -- list containing time-index of swwFile to plot (e.g. [1, 10, 32] ) or 'last', or 'max', or 'all'
     
    614618               
    615619    """
     620
     621    #import pdb
     622    #pdb.set_trace()
     623
    616624    try:
    617625        import gdal
     
    624632        from matplotlib import nxutils
    625633    except:
    626         raise Exception, 'Required modules not installed for mkgeotif'
    627 
     634        raise Exception, 'Required modules not installed for Make_Geotif'
     635
     636
     637    # Check whether swwFile is an array, and if so, redefine various inputs to
     638    # make the code work
     639    if(type(swwFile)==scipy.ndarray):
     640        import copy
     641        xyzPoints=copy.copy(swwFile)
     642        swwFile=None
    628643
    629644    if(((EPSG_CODE is None) & (proj4string is None) )|
     
    638653        pass
    639654
    640     # Read in ANUGA outputs
    641     # FIXME: It would be good to support reading of data subsets
    642     if(verbose):
    643         print 'Reading sww File ...'
    644     p=util.get_output(swwFile,min_allowed_height)
    645     p2=util.get_centroids(p,velocity_extrapolation)
    646     swwIn=scipy.io.netcdf_file(swwFile)
    647     xllcorner=swwIn.xllcorner
    648     yllcorner=swwIn.yllcorner
    649 
    650     if(myTimeStep=='all'):
    651         myTimeStep=range(len(p2.time))
    652     # Ensure myTimeStep is a list
    653     if type(myTimeStep)!=list:
    654         myTimeStep=[myTimeStep]
    655 
    656     if(verbose):
    657         print 'Extracting required data ...'
    658     # Get ANUGA points
    659     swwX=p2.x+xllcorner
    660     swwY=p2.y+yllcorner
     655    if(swwFile is not None):
     656        # Read in ANUGA outputs
     657        # FIXME: It would be good to support reading of data subsets
     658        if(verbose):
     659            print 'Reading sww File ...'
     660        p=util.get_output(swwFile,min_allowed_height)
     661        p2=util.get_centroids(p,velocity_extrapolation)
     662        swwIn=scipy.io.netcdf_file(swwFile)
     663        xllcorner=swwIn.xllcorner
     664        yllcorner=swwIn.yllcorner
     665
     666        if(myTimeStep=='all'):
     667            myTimeStep=range(len(p2.time))
     668        # Ensure myTimeStep is a list
     669        if type(myTimeStep)!=list:
     670            myTimeStep=[myTimeStep]
     671
     672        if(verbose):
     673            print 'Extracting required data ...'
     674        # Get ANUGA points
     675        swwX=p2.x+xllcorner
     676        swwY=p2.y+yllcorner
     677    else:
     678        # Get the point data from the 3-column array
     679        if(xyzPoints.shape[1]!=3):
     680            raise Exception, 'If an array is passed, it must have exactly 3 columns'
     681        if(len(output_quantities)!=1):
     682            raise Exception, 'Can only have 1 output quantity when passing an array'
     683        swwX=xyzPoints[:,0]
     684        swwY=xyzPoints[:,1]
     685        myTimeStep=['pointData']
    661686
    662687    # Grid for meshing
     
    681706    swwXY=scipy.array([swwX[:],swwY[:]]).transpose()
    682707    # Get index of nearest point
    683     index_qFun=scipy.interpolate.NearestNDInterpolator(swwXY,scipy.arange(len(swwX),dtype='int32').transpose())
     708    index_qFun=scipy.interpolate.NearestNDInterpolator(swwXY,scipy.arange(len(swwX),dtype='int64').transpose())
    684709    gridXY_array=scipy.array([scipy.concatenate(gridX),scipy.concatenate(gridY)]).transpose()
    685710    gridqInd=index_qFun(gridXY_array)
     
    702727       
    703728
    704             if(myTS!='max'):
     729            #if(myTS!='max'):
     730            if(type(myTS)=='int'):
    705731                if(output_quantity=='stage'):
    706732                    gridq=p2.stage[myTS,:][gridqInd]
     
    715741                    gridq=p2.elev[gridqInd]
    716742                timestepString=str(round(p2.time[myTS]))
    717             else:
     743            elif (myTS=='max'):
    718744                if(output_quantity=='stage'):
    719745                    gridq=p2.stage.max(axis=0)[gridqInd]
     
    728754                    gridq=p2.elev[gridqInd]
    729755                timestepString='max'
     756            elif(myTS=='pointData'):
     757                gridq=xyzPoints[:,2][gridqInd]
     758
    730759
    731760            if(bounding_polygon is not None):
     
    734763
    735764            # Make name for output file
    736             output_name=output_dir+'/'+os.path.splitext(os.path.basename(swwFile))[0] + '_'+\
    737                         output_quantity+'_'+timestepString+\
    738                         '_'+str(myTS)+'.tif'
     765            if(myTS!='pointData'):
     766                output_name=output_dir+'/'+os.path.splitext(os.path.basename(swwFile))[0] + '_'+\
     767                            output_quantity+'_'+timestepString+\
     768                            '_'+str(myTS)+'.tif'
     769            else:
     770                output_name=output_dir+'/'+'PointData_'+output_quantity+'.tif'
    739771
    740772            if(verbose):
Note: See TracChangeset for help on using the changeset viewer.