"""Create mesh and time boundary for the Okushiri island validation
"""

from anuga.pmesh.mesh import *
from anuga.pmesh.mesh_interface import create_mesh_from_regions
from anuga.coordinate_transforms.geo_reference import Geo_reference
from anuga.geospatial_data import Geospatial_data
from anuga.config import netcdf_float

from Scientific.IO.NetCDF import NetCDFFile

import project

def prepare_bathymetry(filename):
    """Convert benchmark 2 bathymetry to NetCDF pts file.
    This is a 'throw-away' code taylor made for files like
    'Benchmark_2_bathymetry.txt' from the LWRU2 benchmark
    """

    print 'Creating', filename
    
    # Read the ascii (.txt) version of this file,
    # make it comma separated and invert the bathymetry
    # (Below mean sea level should be negative)
    infile = open(filename[:-4] + '.txt')

    points = []
    attribute = []
    for line in infile.readlines()[1:]: #Skip first line (the header)
        fields = line.strip().split()

        x = float(fields[0])
        y = float(fields[1])
        z = -float(fields[2]) # Bathymetry is inverted in original file
        
        points.append([x,y])
        attribute.append(z)
    infile.close()

    # Convert to geospatial data and store as NetCDF
    G = Geospatial_data(data_points=points,
                        attributes=attribute)
    G.export_points_file(filename)
    


def prepare_timeboundary(filename):
    """Convert benchmark 2 time series to NetCDF tms file.
    This is a 'throw-away' code taylor made for files like
    'Benchmark_2_input.txt' from the LWRU2 benchmark
    """

    print 'Creating', filename

    # Read the ascii (.txt) version of this file
    fid = open(filename[:-4] + '.txt')

    # Skip first line
    line = fid.readline()

    # Read remaining lines
    lines = fid.readlines()
    fid.close()


    N = len(lines)
    T = num.zeros(N, num.float)  #Time
    Q = num.zeros(N, num.float)  #Values

    for i, line in enumerate(lines):
        fields = line.split()

        T[i] = float(fields[0])
        Q[i] = float(fields[1])


    # Create tms NetCDF file

    fid = NetCDFFile(filename, 'w')
    fid.institution = 'Geoscience Australia'
    fid.description = 'Input wave for Benchmark 2'
    fid.starttime = 0.0
    fid.createDimension('number_of_timesteps', len(T))
    fid.createVariable('time', netcdf_float, ('number_of_timesteps',))
    fid.variables['time'][:] = T

    fid.createVariable('stage', netcdf_float, ('number_of_timesteps',))
    fid.variables['stage'][:] = Q[:]

    fid.createVariable('xmomentum', netcdf_float, ('number_of_timesteps',))
    fid.variables['xmomentum'][:] = 0.0

    fid.createVariable('ymomentum', netcdf_float, ('number_of_timesteps',))
    fid.variables['ymomentum'][:] = 0.0

    fid.close()


#-------------------------------------------------------------
if __name__ == "__main__":


    # Prepare bathymetry
    prepare_bathymetry(project.bathymetry_filename)

    # Prepare time boundary
    prepare_timeboundary(project.boundary_filename)


    #--------------------------------------------------------------------------
    # Create the triangular mesh based on overall clipping polygon with a
    # tagged
    # boundary and interior regions defined in project.py along with
    # resolutions (maximal area of per triangle) for each polygon
    #--------------------------------------------------------------------------


    base_resolution = 1 # Use this to coarsen or refine entire mesh. 

    # Basic geometry and bounding polygon
    xleft   = 0
    xright  = 5.448
    ybottom = 0
    ytop    = 3.402

    point_sw = [xleft, ybottom]
    point_se = [xright, ybottom]
    point_nw = [xleft, ytop]    
    point_ne = [xright, ytop]

    bounding_polygon = [point_se,
                        point_ne,
                        point_nw,
                        point_sw]
    

    # Localised refined area for gulleys
    xl = 4.8
    xr = 5.3
    yb = 1.6
    yt = 2.3
    p0 = [xl, yb]
    p1 = [xr, yb]
    p2 = [xr, yt]
    p3 = [xl, yt]
    
    gulleys = [p0, p1, p2, p3]
    

    # Island area and drawdown region (original)
    #island_0 = [xleft + 2*(xright-xleft)/3+1.2, ytop-0.5]
    #island_1 = [xleft + 2*(xright-xleft)/3+0.5, ybottom + 2*(ytop-ybottom)/3]
    #island_2 = [xleft + (xright-xleft)/2+0.3, ybottom + 2*(ytop-ybottom)/3-0.3]
    #island_3 = [xleft + (xright-xleft)/2+0.3, ybottom + (ytop-ybottom)/3+0.3]
    #island_4 = [xleft + 2*(xright-xleft)/3+0.4, ybottom + (ytop-ybottom)/3-0.3]
    #island_5 = [xleft + 2*(xright-xleft)/3+1.2, ybottom+0.2]
    #island_6 = [xl-.01, yb]  #OK
    #island_7 = [xl-.01, yt]  #OK      


    # Island area and drawdown region
    island_0 = [xleft + 2*(xright-xleft)/3+1.2, ytop-0.5]
    island_1 = [xleft + 2*(xright-xleft)/3+0.5, ybottom + 2*(ytop-ybottom)/3]
    island_2 = [xleft + (xright-xleft)/2+0.4, ybottom + 2*(ytop-ybottom)/3-0.3]
    island_3 = [xleft + (xright-xleft)/2+0.4, ybottom + (ytop-ybottom)/3+0.3]
    island_4 = [xleft + 2*(xright-xleft)/3+0.4, ybottom + (ytop-ybottom)/3-0.3]
    island_5 = [xleft + 2*(xright-xleft)/3+1.2, ybottom+0.8]
    island_6 = [xl-.01, yb]  # Keep right edge just off the gulleys
    island_7 = [xl-.01, yt]
 
    island = [island_0, island_1, island_2,
              island_3, island_4, island_5,
              island_6, island_7]


    # Region spanning half right hand side of domain just inside boundary (org)
    #rhs_nw = [xleft + (xright-xleft)/3+1, ytop-0.02]
    #rhs_sw = [xleft + (xright-xleft)/3+1, ybottom+0.02]
    #rhs_se = [xright-0.02, ybottom+0.02]
    #rhs_ne = [xright-0.02, ytop-0.02]

    # Region spanning half right hand side of domain just inside boundary
    rhs_nw = [xleft + (xright-xleft)/3+1, ytop-1.4]
    rhs_sw = [xleft + (xright-xleft)/3+1, ybottom+0.5]
    rhs_se = [xright-0.1, ybottom+0.2]
    rhs_ne = [xright-0.1, ytop-0.2]        

    rhs_region = [rhs_nw, rhs_ne, rhs_se, rhs_sw]

    
    # Interior regions and creation of mesh
    interior_regions = [[rhs_region, 0.0005],
                        [island, 0.0002*base_resolution],
                        [gulleys, 0.00002*base_resolution]]    

    meshname = project.mesh_filename + '.msh'
    m = create_mesh_from_regions(bounding_polygon,
                                 boundary_tags={'wall': [0, 1, 3],
                                                'wave': [2]},     
                                 maximum_triangle_area=0.1*base_resolution,
                                 interior_regions=interior_regions,
                                 filename=project.mesh_filename,
                                 verbose=True)

