def calc_surface_range(self):
        # ~~> Dependancies
        args = subset_variables_slf("FREE SURFACE", self.slf.varnames)[0]
        # ~~> New variable name
        calcs = {'vars': [["SURFACE RANGE   ", "M               "]]}

        # ~~> Initial value for new variable
        def init(vrs, ivars, t_0, t_i):
            return [
                np.array(vrs[ivars[0]], copy=True),
                np.array(vrs[ivars[0]], copy=True)
            ]

        calcs.update({'init': (init, args)})

        # ~~> Calculation for new variable
        def calc(vrs, ivars, t_0, t_i, vari):
            return [
                np.minimum(vrs[ivars[0]], vari[0]),
                np.maximum(vrs[ivars[0]], vari[1])
            ]

        calcs.update({'calc': (calc, args)})

        # ~~> Conclusion step for new variable
        def stop(t_0, t_i, vari):
            return [vari[1] - vari[0]]

        calcs.update({'stop': stop})
        # ~~> Store
        self.calcs.append(calcs)
    def calc_residual_velocity(self):
        # ~~> Dependancies
        args = subset_variables_slf("VELOCITY U;VELOCITY V",
                                    self.slf.varnames)[0]
        # ~~> New variable name
        calcs = {
            'vars': [["RESIDUAL U      ", "M/S             "],
                     ["RESIDUAL V      ", "M/S             "]]
        }

        # ~~> Initial value for new variable
        def init(vrs, ivars, t_0, t_i):
            return [
                np.zeros(self.slf.npoin3, dtype=np.float32),
                np.zeros(self.slf.npoin3, dtype=np.float32)
            ]

        calcs.update({'init': (init, args)})

        # ~~> Calculation for new variable
        def calc(vrs, ivars, t_0, t_i, vari):
            vari[0] += vars[ivars[0]]
            vari[1] += vars[ivars[1]]
            return vari

        calcs.update({'calc': (calc, args)})

        # ~~> Conclusion step for new variable
        def stop(t_0, t_i, vari):
            return vari / (t_i - t_0)

        calcs.update({'stop': stop})
        # ~~> Store
        self.calcs.append(calcs)
    def calc_maximum_speed(self):
        # ~~> Dependancies
        args = subset_variables_slf("VELOCITY U;VELOCITY V",
                                    self.slf.varnames)[0]
        # ~~> New variable name
        calcs = {'vars': [["MAXIMUM SPEED   ", "M/S             "]]}

        # ~~> Initial value for new variable
        def init(vrs, ivars, t_0, t_i):
            return [np.zeros(self.slf.npoin3, dtype=np.float32)]

        calcs.update({'init': (init, args)})

        # ~~> Calculation for new variable
        def calc(vrs, ivars, t_0, t_i, vari):
            tmp = np.power(vrs[ivars[0]], 2) + np.power(vrs[ivars[1]], 2)
            return [np.maximum(vari[0], np.power(tmp, (1.0 / 2.0)))]

        calcs.update({'calc': (calc, args)})

        # ~~> Conclusion step for new variable
        def stop(t_0, t_i, vari):
            return vari

        calcs.update({'stop': stop})
        # ~~> Store
        self.calcs.append(calcs)
 def update_vars(self, vrs):
     allvars = self.slf.varnames
     allvars.extend(self.slf.cldnames)
     allunits = self.slf.varunits
     allunits.extend(self.slf.cldunits)
     self.slf.varindex = subset_variables_slf(vrs, allvars)[0]
     nbv1 = 0
     varnames = []
     varunits = []
     nbv2 = 0
     cldnames = []
     cldunits = []
     for i in self.slf.varindex:
         if i <= self.slf.nbv1:
             nbv1 += 1
             varnames.append(allvars[i])
             varunits.append(allunits[i])
         else:
             nbv2 += 1
             varnames.append(allvars[i])
             varunits.append(allunits[i])
     self.slf.nbv1 = nbv1
     self.slf.varnames = varnames
     self.slf.varunits = varunits
     self.slf.nbv2 = nbv2
     self.slf.cldnames = cldnames
     self.slf.cldunits = cldunits
예제 #5
0
    def calc_water_depth(self):
        self.slf.varnames.append("WATER DEPTH     ")
        self.slf.varunits.append("M               ")
        self.slf.nbv1 += 1
        args = subset_variables_slf("FREE SURFACE;BOTTOM",
                                    self.slf.varnames)[0]

        def calc(vrs, ivars):
            return [vrs[ivars[0]] - vrs[ivars[1]]]

        self.calcs.append([calc, args])
예제 #6
0
    def calc_kinetic_energy(self):
        self.slf.varnames.append("KINETIC ENERGY  ")
        self.slf.varunits.append("?               ")
        self.slf.nbv1 += 1
        args = subset_variables_slf("VELOCITY U;VELOCITY V",
                                    self.slf.varnames)[0]

        def calc(vrs, ivars):
            var1 = np.square(vrs[ivars[0]])
            var2 = np.square(vrs[ivars[1]])
            return [np.power(var1 + var2, (3.0 / 2.0))]

        self.calcs.append([calc, args])
    def calc_peak_time_modulo_m2(self):
        # ~~> Dependancies
        args = subset_variables_slf("FREE SURFACE", self.slf.varnames)[0]
        # ~~> New variable name
        calcs = {
            'vars': [["TIME OF PEAK    ", "H               "],
                     ["SURFACE AT PEAK ", "M               "]]
        }

        # ~~> Initial value for new variable
        def init(vrs, ivars, t_0, t_i):
            return [
                np.zeros(self.slf.npoin3, dtype=np.float32),
                np.zeros(self.slf.npoin3, dtype=np.float32)
            ]

        calcs.update({'init': (init, args)})

        # ~~> Calculation for new variable
        def calc(vrs, ivars, t_0, t_i, vari):
            for ipoin in range(self.slf.npoin3):
                if vrs[ivars[0]][ipoin] > vari[1][ipoin]:
                    # Modulo M2 (this could be an argument)
                    vari[0][ipoin] = ((t_i - t_0) / 3600.0) % 12.42
                    vari[1][ipoin] = vrs[ivars[0]][ipoin]
            return vari

        calcs.update({'calc': (calc, args)})

        # ~~> Conclusion step for new variable
        def stop(t_0, t_i, vari):
            return vari

        calcs.update({'stop': stop})
        # ~~> Store
        self.calcs.append(calcs)
예제 #8
0
def main():
    """ Main function of convertToIni """

    # <<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
    # ~~~~ Dependencies towards other modules ~~~~~~~~~~~~~~~~~~~~~~~~~~
    from argparse import ArgumentParser, RawDescriptionHelpFormatter
    from data_manip.formats.selafin import Selafin
    from data_manip.extraction.parser_selafin import subset_variables_slf, \
        get_value_history_slf
    from pretel.meshes import xys_locate_mesh

    # <<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
    # ~~~~ Reads config file ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
    print('\n\nInterpreting command line options\n' + '~' * 72 + '\n')
    parser = ArgumentParser(
        formatter_class=RawDescriptionHelpFormatter,
        description=('''\n
A script to map 2D or 3D outter model results stored in a SELAFIN file,\
 onto the
   one frame of contained SELAFIN file of your choosing (your MESH).
      '''),
        usage=' (--help for help)\n---------\n       =>  '
        '%(prog)s  geo-mesh.slf in-result.slf out-init.slf \n---------')
    parser.add_argument("args", default='', nargs=3)
    options = parser.parse_args()

    # <<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
    # ~~~~ slf new mesh ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
    geo_file = options.args[0]
    if not path.exists(geo_file):
        raise TelemacException(
            '... the provided geo_file does not seem to exist: {}'
            '\n\n'.format(geo_file))


# Find corresponding (x,y) in corresponding new mesh
    print('   +> getting hold of the GEO file and of its bathymetry')
    geo = Selafin(geo_file)
    xys = np.vstack((geo.meshx, geo.meshy)).T
    _ = geo.get_variables_at(0,
                             subset_variables_slf("BOTTOM: ",
                                                  geo.varnames)[0])[0]

    # <<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
    # ~~~~ slf existing res ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
    slf_file = options.args[1]
    if not path.exists(slf_file):
        raise TelemacException(
            '... the provided geo_file does not seem to exist: {}'
            '\n\n'.format(slf_file))

    slf = Selafin(slf_file)
    slf.set_kd_tree()
    slf.set_mpl_tri()

    print('   +> support extraction')
    # Extract triangles and weights in 2D
    support2d = []
    ibar = 0
    pbar = ProgressBar(maxval=len(xys)).start()
    for xyi in xys:
        support2d.append(
            xys_locate_mesh(xyi, slf.ikle2, slf.meshx, slf.meshy, slf.tree,
                            slf.neighbours))
        ibar += 1
        pbar.update(ibar)
    pbar.finish()
    # Extract support in 3D
    support3d = list(zip(support2d, len(xys) * [range(slf.nplan)]))

    # <<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
    # ~~~~ writes INI header ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~

    ini_file = options.args[2]
    ini = Selafin('')
    ini.fole = {}
    ini.fole.update({'hook': open(ini_file, 'wb')})
    ini.fole.update({'name': ini_file})
    ini.fole.update({'endian': ">"})  # big endian
    ini.fole.update({'float': ('f', 4)})  # single precision

    # Meta data and variable names
    ini.title = ''
    ini.nbv1 = 5
    # /!\ ELEVATION has to be the first variable
    # (for possible vertical re-interpolation within TELEMAC)
    ini.varnames = [
        'ELEVATION Z     ', 'VELOCITY U      ', 'VELOCITY V      ',
        'SALINITY        ', 'TEMPERATURE     '
    ]
    ini.varunits = [
        'M               ', 'M/S             ', 'M/S             ',
        '                ', '                '
    ]
    ini.nvar = ini.nbv1
    ini.varindex = range(ini.nvar)

    # sizes and mesh connectivity
    ini.nplan = slf.nplan
    ini.ndp2 = 3
    ini.ndp3 = 6
    ini.npoin2 = geo.npoin2
    ini.npoin3 = geo.npoin2 * ini.nplan
    ini.nelem2 = geo.nelem2
    ini.nelem3 = ini.nelem2 * (ini.nplan - 1)

    print('   +> setting connectivity')
    ini.ikle3 = \
        np.repeat(geo.npoin2*np.arange(ini.nplan-1),
                  geo.nelem2*ini.ndp3)\
        .reshape((geo.nelem2*(ini.nplan-1), ini.ndp3)) + \
        np.tile(np.add(np.tile(geo.ikle2, 2),
                np.repeat(geo.npoin2*np.arange(2), geo.ndp2)),
                (ini.nplan-1, 1))
    ini.ipob3 = np.ravel(
        np.add(
            np.repeat(geo.ipob2, ini.nplan).reshape((geo.npoin2, ini.nplan)),
            geo.npoin2 * np.arange(ini.nplan)).T)
    ini.iparam = [0, 0, 0, 0, 0, 0, ini.nplan, 0, 0, 0]

    # Mesh coordinates
    ini.meshx = geo.meshx
    ini.meshy = geo.meshy

    print('   +> writing header')
    # Write header
    ini.append_header_slf()

    # <<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
    # ~~~~ writes INI core ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~

    print('   +> setting variables')
    ini.tags['times'] = slf.tags['times']
    # VARIABLE extraction
    vrs = subset_variables_slf(
        "ELEVATION Z: ;VELOCITY U: ;VELOCITY V: "
        ";SALINITY: ;TEMPERATURE: ", slf.varnames)

    # Read / Write data for first time step
    zeros = np.zeros((ini.npoin3, 1), dtype=np.float)

    print('   +> extracting variables')
    data = get_value_history_slf(slf.file, slf.tags, [0], support3d, slf.nvar,
                                 slf.npoin3, slf.nplan, vrs)

    # special case for TEMPERATURE and SALINITY
    data[3] = np.maximum(data[3], zeros)
    data[4] = np.maximum(data[4], zeros)
    print('   +> correcting variables')
    # duplicate values below bottom
    data = np.reshape(
        np.transpose(
            np.reshape(np.ravel(data), (ini.nvar, ini.npoin2, ini.nplan)),
            (0, 2, 1)), (ini.nvar, ini.npoin3))
    print('   +> writing variables')
    ini.append_core_time_slf(0)
    ini.append_core_vars_slf(data)

    # Close ini_file
    ini.fole['hook'].close()

    # <<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
    # ~~~~ Jenkins' success message ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
    print('\n\nMy work is done\n\n')

    sys.exit(0)
예제 #9
0
def generate_bnd(cli_file,
                 geo_file,
                 slf_file,
                 bnd_file,
                 varnames,
                 varunits,
                 showbar=True):
    """
    @param cli_file
    @param geo_file
    @param slf_file
    @param bnd_file
    @param varnames
    @param varunits
    @param showbar (boolean) If True display a showbar for the progress
    """

    # <<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
    # ~~~~ cli+slf new mesh ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
    if not path.exists(cli_file):
        raise TelemacException(\
             '... the provided cli_file does not seem to exist:'
             ' {}\n\n'.format(cli_file))
    if not path.exists(geo_file):
        raise TelemacException(\
            '... the provided geo_file does not seem to exist: '
            '{}\n\n'.format(geo_file))

    if len(varnames) != len(varunits):
        raise TelemacException(\
          'Not the same number of variables and units\nvarnames: {}\nvarunits: {}'
          '{}\n\n'.format(varnames, varunits))

    # Read the new CLI file to get boundary node numbers
    print('   +> getting hold of the Conlim file and of its liquid boundaries')
    cli = Conlim(cli_file)
    # Keeping only open boundary nodes
    bor = np.extract(cli.bor['lih'] != 2, cli.bor['n'])

    # Find corresponding (x,y) in corresponding new mesh
    print('   +> getting hold of the GEO file and of its bathymetry')
    geo = Selafin(geo_file)
    xys = np.vstack((geo.meshx[bor - 1], geo.meshy[bor - 1])).T
    _ = geo.get_variables_at(0,\
                  subset_variables_slf("BOTTOM: ", geo.varnames)[0])[0]

    # <<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
    # ~~~~ slf existing res ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
    if not path.exists(slf_file):
        raise TelemacException(\
               '... the provided slf_file does not seem to exist: '
               '{}\n\n'.format(slf_file))
    slf = Selafin(slf_file)
    slf.set_kd_tree()
    slf.set_mpl_tri()

    print('   +> support extraction')
    # Extract triangles and weigths in 2D
    support2d = []
    if showbar:
        ibar = 0
        pbar = ProgressBar(maxval=len(xys)).start()
    for xyi in xys:
        support2d.append(
            xys_locate_mesh(xyi, slf.ikle2, slf.meshx, slf.meshy, slf.tree,
                            slf.neighbours))
        if showbar:
            ibar += 1
            pbar.update(ibar)
    if showbar:
        pbar.finish()
    # Extract support in 3D
    support3d = list(zip(support2d, len(xys) * [range(slf.nplan)]))

    # <<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
    # ~~~~ writes BND header ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~

    bnd = Selafin('')
    bnd.fole = {}
    bnd.fole.update({'hook': open(bnd_file, 'wb')})
    bnd.fole.update({'name': bnd_file})
    bnd.fole.update({'endian': ">"})  # big endian
    bnd.fole.update({'float': ('f', 4)})  # single precision

    # Meta data and variable names
    bnd.title = ''
    bnd.nbv1 = len(varnames)
    # /!\ ELEVATION has to be the first variable
    # (for possible vertical re-interpolation within TELEMAC)

    bnd.varnames = []
    bnd.varunits = []
    for var, unit in zip(varnames, varunits):
        new_var = var + (16 - len(var)) * " "
        new_unit = unit + (16 - len(unit)) * " "
        bnd.varnames.append(new_var)
        bnd.varunits.append(new_unit)

    bnd.nvar = bnd.nbv1
    bnd.varindex = range(bnd.nvar)

    # Sizes and mesh connectivity
    bnd.nplan = slf.nplan
    # Number of nodes per boundary element  (ndp2 in 2D and ndp3 in 3D)
    bnd.ndp2 = 2
    bnd.ndp3 = 4
    bnd.npoin2 = len(bor)
    bnd.npoin3 = bnd.npoin2 * slf.nplan
    bnd.iparam = [0, 0, 0, 0, 0, 0, bnd.nplan, 0, 0, 1]
    bnd.ipob2 = bor  # /!\ note that ipobo keeps the original numbering
    print('   +> masking and setting connectivity')
    # Set the array that only includes elements of geo.ikle2
    # with at least two nodes in bor
    array_1d = np.in1d(geo.ikle2, np.sort(bor - 1))
    mask = geo.ikle2[np.where(
        np.sum(array_1d.reshape(geo.nelem2, geo.ndp2), axis=1) == 2)]
    # this ikle2 keeps the original numbering
    ikle2 = np.ravel(mask)[np.in1d(mask,
                                   np.sort(bor - 1))].reshape(len(mask), 2)
    # ~~> re-numbering ikle2 as a local connectivity matrix
    knolg, _ = np.unique(np.ravel(ikle2), return_index=True)
    knogl = dict(zip(knolg, range(len(knolg))))
    bnd.ikle2 = -np.ones_like(ikle2, dtype=np.int)
    for k in range(len(ikle2)):
        # /!\ bnd.ikle2 has a local numbering, fit to the boundary elements
        bnd.ikle2[k] = [knogl[ikle2[k][0]], knogl[ikle2[k][1]]]
    # Last few numbers
    bnd.nelem2 = len(bnd.ikle2)
    if slf.nplan > 1:
        bnd.nelem3 = bnd.nelem2 * (slf.nplan - 1)
    else:
        bnd.nelem3 = bnd.nelem2
        bnd.ndp3 = bnd.ndp2
    # 3D structures
    if slf.nplan > 1:
        bnd.ipob3 = np.ravel(np.add(np.repeat(bnd.ipob2, slf.nplan)\
                                      .reshape((bnd.npoin2, slf.nplan)),
                                    bnd.npoin2*np.arange(slf.nplan)).T)
        bnd.ikle3 = \
            np.repeat(bnd.npoin2*np.arange(slf.nplan-1),
                      bnd.nelem2*bnd.ndp3)\
              .reshape((bnd.nelem2*(slf.nplan-1), bnd.ndp3)) + \
            np.tile(np.add(np.tile(bnd.ikle2, 2),
                           np.repeat(bnd.npoin2*np.arange(2), bnd.ndp2)),
                    (slf.nplan-1, 1))
    else:
        bnd.ipob3 = bnd.ipob2
        bnd.ikle3 = bnd.ikle2
    # Mesh coordinates
    bnd.meshx = geo.meshx[bor - 1]
    bnd.meshy = geo.meshy[bor - 1]

    print('   +> writing header')
    # Write header
    bnd.append_header_slf()

    # <<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
    # ~~~~ writes BND core ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~

    print('   +> setting variables')
    # TIME and DATE extraction
    bnd.datetime = slf.datetime
    bnd.tags['times'] = slf.tags['times']
    # VARIABLE extraction
    list_var = varnames[0] + ": "
    for var in varnames[1:]:
        list_var += ";" + var + ": "

    vrs = subset_variables_slf(list_var, slf.varnames)

    # Read / Write data, one time step at a time to support large files
    print('   +> reading / writing variables')
    if showbar:
        pbar = ProgressBar(maxval=len(slf.tags['times'])).start()
    zeros = np.zeros((bnd.npoin3, 1), dtype=np.float)
    for itime in range(len(slf.tags['times'])):
        data = get_value_history_slf(slf.file, slf.tags, [itime], support3d,
                                     slf.nvar, slf.npoin3, slf.nplan, vrs)
        data = np.reshape(
            np.transpose(
                np.reshape(np.ravel(data), (bnd.nvar, bnd.npoin2, bnd.nplan)),
                (0, 2, 1)), (bnd.nvar, bnd.npoin3))
        bnd.append_core_time_slf(itime)
        bnd.append_core_vars_slf(data)
        if showbar:
            pbar.update(itime)
    if showbar:
        pbar.finish()

    # Close bnd_file
    bnd.fole['hook'].close()
예제 #10
0
def generate_atm(geo_file, slf_file, atm_file, ll2utm):

# <<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
# ~~~~ slf new mesh ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
    if not path.exists(geo_file):
        raise TelemacException(\
            '... the provided geo_file does not '
            'seem to exist: {}\n\n'.format(geo_file))

# Find corresponding (x,y) in corresponding new mesh
    print('   +> getting hold of the GEO file')
    geo = Selafin(geo_file)
    if ll2utm is not None:
        zone = int(ll2utm[:-1])
        zone_letter = ll2utm[-1]
        x, y = to_latlon(geo.meshx, geo.meshy, zone, zone_letter)
    else:
        x = geo.meshx
        y = geo.meshy
    xys = np.vstack((x, y)).T

# <<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
# ~~~~ slf existing res ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
    if not path.exists(slf_file):
        raise TelemacException(\
                '... the provided slf_file does not seem to exist: '
                '{}\n\n'.format(slf_file))
    slf = Selafin(slf_file)
    slf.set_kd_tree()
    slf.set_mpl_tri()

    print('   +> support extraction')
    # Extract triangles and weights in 2D
    support2d = []
    ibar = 0
    pbar = ProgressBar(maxval=len(xys)).start()
    for xyi in xys:
        support2d.append(xys_locate_mesh(xyi, slf.ikle2, slf.meshx, slf.meshy,
                                         slf.tree, slf.neighbours))
        ibar += 1
        pbar.update(ibar)
    pbar.finish()
    # Extract support in 3D
    support3d = list(zip(support2d, len(xys)*[range(slf.nplan)]))

# <<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
# ~~~~ writes ATM header ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~

    atm = Selafin('')
    atm.fole = {}
    atm.fole.update({'hook':open(atm_file, 'wb')})
    atm.fole.update({'name':atm_file})
    atm.fole.update({'endian':">"})     # big endian
    atm.fole.update({'float':('f', 4)})  # single precision

    # Meta data and variable names
    atm.title = ''
    atm.varnames = []
    atm.varunits = []
    if 'WIND VELOCITY U ' in slf.varnames:
        atm.varnames.append('WIND VELOCITY U ')
        atm.varunits.append('M/S             ')
    if 'WIND VELOCITY V ' in slf.varnames:
        atm.varnames.append('WIND VELOCITY V ')
        atm.varunits.append('M/S             ')
    if 'SURFACE PRESSURE' in slf.varnames:
        atm.varnames.append('SURFACE PRESSURE')
        atm.varunits.append('UI              ')
    if 'AIR TEMPERATURE ' in slf.varnames:
        atm.varnames.append('AIR TEMPERATURE ')
        atm.varunits.append('DEGREES         ')
    if not atm.varnames:
        raise TelemacException(
            'There are no meteorological variables to convert!')
    atm.nbv1 = len(atm.varnames)
    atm.nvar = atm.nbv1
    atm.varindex = range(atm.nvar)

    # Sizes and mesh connectivity
    atm.nplan = slf.nplan          # it should be 2d but why the heack not ...
    atm.ndp2 = slf.ndp2
    atm.ndp3 = slf.ndp3
    atm.npoin2 = geo.npoin2
    atm.npoin3 = geo.npoin2*atm.nplan
    atm.nelem2 = geo.nelem2

    print('   +> setting connectivity')
    if atm.nplan > 1:
        atm.nelem3 = geo.nelem2*(atm.nplan-1)
        atm.ikle2 = geo.ikle2
        atm.ikle3 = \
            np.repeat(geo.npoin2*np.arange(atm.nplan-1),
                      geo.nelem2*atm.ndp3)\
              .reshape((geo.nelem2*(atm.nplan-1), atm.ndp3)) + \
            np.tile(np.add(np.tile(geo.ikle2, 2),
                           np.repeat(geo.npoin2*np.arange(2),
                                     geo.ndp2)),
                    (atm.nplan-1, 1))
        atm.ipob2 = geo.ipob2
        atm.ipob3 = np.ravel(np.add(np.repeat(geo.ipob2, atm.nplan)\
                                      .reshape((geo.npoin2, atm.nplan)),
                                    geo.npoin2*np.arange(atm.nplan)).T)
    else:
        atm.nelem3 = geo.nelem2
        atm.ikle2 = geo.ikle2
        atm.ikle3 = geo.ikle3
        atm.ipob2 = geo.ipob2
        atm.ipob3 = geo.ipob3
    atm.iparam = [0, 0, 0, 0, 0, 0, 0, np.count_nonzero(atm.ipob2), 0, 1]

    # Mesh coordinates
    atm.meshx = geo.meshx
    atm.meshy = geo.meshy

    print('   +> writing header')
    # Write header
    atm.datetime = slf.datetime
    atm.append_header_slf()

# <<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
# ~~~~ writes ATM core ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~

    print('   +> setting variables')
    # TIME and DATE extraction
    atm.tags['times'] = slf.tags['times']
    # VARIABLE extraction
    vrs = subset_variables_slf(';'.join([var+': ' for var in atm.varnames]),
                               slf.varnames)

    # Read / Write data, one time step at a time to support large files
    pbar = ProgressBar(maxval=len(slf.tags['times'])).start()
    for time in range(len(slf.tags['times'])):

        data = get_value_history_slf(slf.file, slf.tags, [time], support3d,
                                     slf.nvar, slf.npoin3, slf.nplan, vrs)
        # special cases ?
        atm.append_core_time_slf(time)
        atm.append_core_vars_slf(np.reshape(np.transpose(\
                  np.reshape(np.ravel(data),
                             (atm.nvar, atm.npoin2, atm.nplan)),
                  (0, 2, 1)),
                                            (atm.nvar, atm.npoin3)))
        pbar.update(time)
    pbar.finish()

    # Close atm_file
    atm.fole['hook'].close()
def main():
    """ Main function of convertToSPG """
    # <<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
    # ~~~~ Dependencies towards other modules ~~~~~~~~~~~~~~~~~~~~~~~~~~
    from argparse import ArgumentParser, RawDescriptionHelpFormatter
    from data_manip.formats.selafin import Selafin
    from data_manip.extraction.parser_selafin import subset_variables_slf, \
                                       get_value_history_slf
    from pretel.meshes import xys_locate_mesh

    # <<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
    # ~~~~ Reads config file ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
    print('\n\nInterpreting command line options\n' + '~' * 72 + '\n')
    parser = ArgumentParser(\
        formatter_class=RawDescriptionHelpFormatter,
        description=('''\n
Creates a binary sponge file from a global model (SLF form)
    by interpolating on a given GEO model domain. The time series in
    the SPG file are extracted only at the nodes where the SPONGE mask
    value is more than 0.5.
        '''),
        usage=' (--help for help)\n---------\n       =>  '\
              '%(prog)s  geo-mask.slf in-result.slf out-sponge.slf \n---------')
    parser.add_argument("args", default='', nargs=3)
    options = parser.parse_args()

    # <<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
    # ~~~~ slf new mesh ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
    geo_file = options.args[0]
    if not path.exists(geo_file):
        raise TelemacException(\
                'Could not find geo_file: {}\n\n'.format(geo_file))

    # Read the new GEO file and its SPONGE mask variable
    print('   +> getting hold of the GEO file and of its SPONGE mask')
    geo = Selafin(geo_file)
    _, spg = geo.get_variables_at(0, subset_variables_slf(\
              "BOTTOM: ;SPONGE mask: ", geo.varnames)[0])[0:2]
    print('   +> extracting the masked elements')
    # Keeping only masked nodes
    array_1d = np.in1d(geo.ikle2, np.sort(np.where(spg > 0.)[0]))
    mask = geo.ikle2[np.where(
        np.sum(array_1d.reshape(geo.nelem2, geo.ndp2), axis=1) > 0)]
    bor = np.unique(mask) + 1

    # Find corresponding (x,y) for the mask
    print('   +> getting hold of the GEO file and of its bathymetry')
    xys = np.vstack((geo.meshx[bor - 1], geo.meshy[bor - 1])).T

    # <<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
    # ~~~~ slf existing res ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
    slf_file = options.args[1]
    if not path.exists(geo_file):
        raise TelemacException(\
                'Could not find slf_file: {}\n\n'.format(slf_file))
    slf = Selafin(slf_file)
    slf.set_kd_tree()
    slf.set_mpl_tri()

    print('   +> support extraction')
    # Extract triangles and weigths in 2D
    support2d = []
    ibar = 0
    pbar = ProgressBar(maxval=len(xys)).start()
    for xyi in xys:
        support2d.append(
            xys_locate_mesh(xyi, slf.ikle2, slf.meshx, slf.meshy, slf.tree,
                            slf.neighbours))
        ibar += 1
        pbar.update(ibar)
    pbar.finish()
    # Extract support in 3D
    support3d = zip(support2d, len(xys) * [range(slf.nplan)])

    # <<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
    # ~~~~ writes BND header ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~

    bnd_file = options.args[2]
    bnd = Selafin('')
    bnd.fole = {}
    bnd.fole.update({'hook': open(bnd_file, 'wb')})
    bnd.fole.update({'name': bnd_file})
    bnd.fole.update({'endian': ">"})  # big endian
    bnd.fole.update({'float': ('f', 4)})  # single precision

    # Meta data and variable names
    bnd.title = ''
    bnd.nbv1 = 5
    # /!\ ELEVATION has to be the first variable
    # (for possible vertical re-interpolation within TELEMAC)
    bnd.varnames = ['ELEVATION Z     ', \
                    'VELOCITY U      ', 'VELOCITY V      ', \
                    'SALINITY        ', 'TEMPERATURE     ']
    bnd.varunits = ['M               ', \
                    'M/S             ', 'M/S             ', \
                    '                ', '                ']
    bnd.nvar = bnd.nbv1
    bnd.varindex = range(bnd.nvar)

    # Sizes and mesh connectivity
    bnd.nplan = slf.nplan
    bnd.ndp2 = 3
    bnd.ndp3 = 6
    bnd.npoin2 = len(bor)
    bnd.npoin3 = bnd.npoin2 * slf.nplan
    bnd.iparam = [0, 0, 0, 0, 0, 0, bnd.nplan, 0, 0, 0]
    bnd.ipob2 = bor  # /!\ Note that ipobO keeps the original numbering
    print('   +> masking and setting connectivity')
    # Set the array that only includes elements of geo.ikle2
    # with at least two nodes in bor
    ikle2 = mask
    # ~~> re-numbering ikle2 as a local connectivity matrix
    knolg, _ = np.unique(np.ravel(ikle2), return_index=True)
    knogl = dict(zip(knolg, range(len(knolg))))
    bnd.ikle2 = -np.ones_like(ikle2, dtype=np.int)
    for k in range(len(ikle2)):
        # /!\ bnd.ikle2 has a local numbering, fit to the boundary elements
        bnd.ikle2[k] = [
            knogl[ikle2[k][0]], knogl[ikle2[k][1]], knogl[ikle2[k][2]]
        ]
    # Last few numbers
    bnd.nelem2 = len(bnd.ikle2)
    if slf.nplan > 1:
        bnd.nelem3 = bnd.nelem2 * (slf.nplan - 1)
    else:
        bnd.nelem3 = bnd.nelem2
    # 3D structures
    if slf.nplan > 1:
        bnd.ipob3 = np.ravel(np.add(np.repeat(bnd.ipob2, slf.nplan)\
                                              .reshape((bnd.npoin2, slf.nplan)),
                                    bnd.npoin2*np.arange(slf.nplan)).T)
        bnd.ikle3 = \
            np.repeat(bnd.npoin2*np.arange(slf.nplan-1),
                      bnd.nelem2*bnd.ndp3)\
              .reshape((bnd.nelem2*(slf.nplan-1), bnd.ndp3)) + \
            np.tile(np.add(np.tile(bnd.ikle2, 2),
                           np.repeat(bnd.npoin2*np.arange(2), bnd.ndp2)),
                    (slf.nplan-1, 1))
    else:
        bnd.ipob3 = bnd.ipob2
        bnd.ikle3 = bnd.ikle2
    # Mesh coordinates
    bnd.meshx = geo.meshx[bor - 1]
    bnd.meshy = geo.meshy[bor - 1]

    print('   +> writing header')
    # Write header
    bnd.append_header_slf()

    # <<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
    # ~~~~ writes BND core ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~

    print('   +> setting variables')
    # TIME and DATE extraction
    bnd.tags['times'] = slf.tags['times']
    # VARIABLE extraction
    vrs = subset_variables_slf("ELEVATION Z: ;VELOCITY U: ;VELOCITY V: "\
                                      ";SALINITY: ;TEMPERATURE: ", slf.varnames)

    # Read / Write data, one time step at a time to support large files
    print('   +> reading / writing variables')
    pbar = ProgressBar(maxval=len(slf.tags['times'])).start()
    zeros = np.zeros((bnd.npoin3, 1), dtype=np.float)
    for time in range(len(slf.tags['times'])):
        data = get_value_history_slf(slf.file, slf.tags, [time], support3d,
                                     slf.nvar, slf.npoin3, slf.nplan, vrs)
        # special case for TEMPERATURE and SALINITY
        data[3] = np.maximum(data[3], zeros)
        data[4] = np.maximum(data[4], zeros)
        data = np.reshape(
            np.transpose(
                np.reshape(np.ravel(data), (bnd.nvar, bnd.npoin2, bnd.nplan)),
                (0, 2, 1)), (bnd.nvar, bnd.npoin3))
        bnd.append_core_time_slf(time)
        bnd.append_core_vars_slf(data)
        pbar.update(time)
    pbar.finish()

    # Close bnd_file
    bnd.fole['hook'].close()

    # <<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<
    # ~~~~ Jenkins' success message ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
    print('\n\nMy work is done\n\n')

    sys.exit(0)