Compute the isentropic Mach number on a blade

Description

Compute the 3D isentropic Mach number field on a blade’s skin.

Parameters

  • base: Base

    The input base corresponds to the blade’s skin.

  • TreatmentPtris: TreatmentPtris

    An object defined upstream with at least in_vars and out_vars options that compute the meridional coordinate and the isentropic pressure variable.

    For more details, please refer to the corresponding TreatmentPtris<convention>.

  • h_H_name: str, default = ‘CoordinateReducedHeight’

    The coordinate reduced height’s name. Basically, ‘CoordinateReducedHeight’ from the antares.treatment.turbomachine.TreatmenthH.TreatmenthH.

  • h_H_range: list(float, float), default = None

    Allows to clip the h_H_name’s range of values. Might be useful to avoid technological effects. Either let it unfilled then the geometry bounds the variable h_H_name. Or a list of float clip it : [minimum_h_H_name, maximum _h_H_name].

  • number_of_heights: int, default = 101

    Spatial discretisation of the height, based on the h_H_name. In other words, number of profiles the blade will be composed of.

  • d_D_range: str, default = ‘0_1’

    Defines the curvilinear reduced coordinates’ range of value. Either ranged from -1 to 1 or from 0 to 1.

  • begin_unwrap: str, default = ‘LE’

    Convention : ‘LE’ (leading edge) or ‘TE’ (trailing edge). In other words, either sorted the blade from the LE to the TE or vice versa. Therefore, if begin_unwrap is set to ‘LE’, the curvilinear coordinate’s zero is at the leading edge.

  • communicator: antares.parallel.controller.ParallelController, default = PARA

    Communicator to communicate with a subset of the original group of processes at once. The default communicator is the singleton PARA created at the import of antares.

Preconditions

The treatment must be applied on a 3D Zone corresponding to a blade skin and its fillet. The Base must have only one Zone and one Instant. Moreover, the treatment requires a reduced height coordinate (h_H_name) in its Instant, for instance, computed with antares.treatment.turbomachine.TreatmenthH.TreatmenthH. Last but not least, a TreatmentPtris<convention> has to be defined upstream with options ‘in_vars’ and ‘out_vars’.

Postconditions

The input base is returned, extended with the variables:

  • Mis

    The isentropic Mach number with a reference pressure based on the given TreatmentPtris.

  • d

    The curvilinear coordinate around each profil discretizing the 3D blade.

  • d_D

    More precisely, out_vars[0] variable’s name defined in the TreatmentPtris<convention>: d_D by default but may be changed by the user. This variable corresponds to the reduced form of the previous one.

Example

import antares

# Fill Ptris information. For more details,
# refer to the corresponding treatment.
myt_Ptris = antares.Treatment('Ptris'+'convention_name')
myt_Ptris['in_vars'] = ['x', 'y', 'z', 'r', 'CoordinateReducedMeridional',
                        'gamma', 'Cp', 'omega', 'Psta', 'Tsta']
myt_Ptris['out_vars'] = ['d_D','Ptris']

# Parameters the Mis calculation on the blade
myt = antares.Treatment('MisOnBlade')
myt['base'] = blade_skin
myt['number_of_heights'] = 101
myt['d_D_range'] = '0_1'
myt['TreatmentPtris'] = myt_Ptris
blade_skin_result = myt.execute()

Main functions

class antares.treatment.turbomachine.TreatmentMisOnBlade.TreatmentMisOnBlade
execute()

Compute the isentropic Mach number on a 3D blade.

Returns:

The 3D blade with the isentropic Mach number computed

Return type:

Base

Example

"""
This script may ne run in parallel with
mpirun -np 16 python misonblade.py

There is no data available for this example.
"""
import os
import numpy as np
import antares
from antares.parallel.controller import PARA

output = 'OUTPUT'
if not os.path.isdir(output):
    os.makedirs(output)

reader = antares.Reader('hdf_cgns')
reader['filename'] = 'data_with_many_zones.cgns'
reader['rename_vars'] = False
base = reader.read()

hub_pts = np.genfromtxt('hub_meridional_line.dat', unpack=False)
shroud_pts = np.genfromtxt('shroud_meridional_line.dat', unpack=False)

t = antares.Treatment('hH')
t['base'] = base
t['families'] = ['ROW1']  # family (volume zone) where to compute hH
t['hub_pts'] = hub_pts
t['shroud_pts'] = shroud_pts
t['extension'] = 1.5
t['number_of_heights'] = 21
t['dmax'] = 1e-05
t['precision'] = 1e-05
t['coordinates'] = ['CoordinateX', 'CoordinateY', 'CoordinateZ']
base, paramgrid = t.execute()

# get the base corresponding to a blade given by a family name
# in parallel, some processors may not hold a part of the blade
try:
    skin_blade = base[base.families['BLADE1']]
except KeyError:
    skin_blade = antares.Base()

myt = antares.Treatment('merge')
myt['base'] = skin_blade
myt['duplicates_detection'] = True
myt['tolerance_decimals'] = 13
skin_blade = myt.execute()

# rename the zone
if skin_blade:
    skin_blade['zone%d'%PARA.rank] = skin_blade[0]
    del skin_blade[0]

skin_blade.set_computer_model('constant_gamma')
skin_blade.compute('gamma')
skin_blade.compute('CP')
skin_blade.compute('R')
skin_blade.compute('Pressure')
skin_blade.compute('Temperature')
wrot = 0#17188.7 * 2*np.pi/60.
skin_blade.compute('omega = {}'.format(wrot))

# compute the reference isentropic Ptr
myt_Ptris = antares.Treatment('ptris<convention>')
myt_Ptris['in_vars'] = ['CoordinateX', 'CoordinateY', 'CoordinateZ', 'R', 'CoordinateReducedMeridional',
                        'gamma', 'CP', 'omega', 'Pressure', 'Temperature']
myt_Ptris['out_vars'] = ['d_D','Ptris']

# Parameters the Mis calculation on the blade
myt = antares.Treatment('MisOnBlade')
myt['base'] = skin_blade
myt['number_of_heights'] = 101
myt['d_D_range'] = '0_1'
myt['TreatmentPtris'] = myt_Ptris
blade_skin_result = myt.execute()

writer = antares.Writer('hdf_cgns')
writer['filename'] = os.path.join(output, 'blade_skin_para_{}p.cgns'.format(PARA.size))
writer['coordinates'] = ['CoordinateX', 'CoordinateY', 'CoordinateZ']
writer['base'] = blade_skin_result
writer.dump()