#!/usr/bin/env python3 # -*- coding: utf-8 -*- # Copyright 2022 University of Liège # # Licensed under the Apache License, Version 2.0 (the "License"); # you may not use this file except in compliance with the License. # You may obtain a copy of the License at # # http://www.apache.org/licenses/LICENSE-2.0 # # Unless required by applicable law or agreed to in writing, software # distributed under the License is distributed on an "AS IS" BASIS, # WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied. # See the License for the specific language governing permissions and # limitations under the License. # @author Paul Dechamps # @date 2022 # Test the blast implementation. The test case is a compressible attached transitional flow. # Tested functionalities; # - Time integration. # - Two time-marching techniques (agressive and safe). # - Transition routines. # - Forced transition. # - Compressible flow routines. # - Laminar and turbulent flow. # # Untested functionalities. # - Separation routines. # - Multiple failure case at one iteration. # - Response to unconverged inviscid solver. # Imports. import blast.utils as viscUtils import numpy as np import os.path from fwk.wutils import parseargs import fwk from fwk.testing import * from fwk.coloring import ccolors import math def cfgInviscid(nthrds, verb): # Parameters return { # Options 'scenario' : 'aerodynamic', 'task' : 'analysis', 'Threads' : nthrds, # number of threads 'Verb' : verb, # verbosity # Model (geometry or mesh) 'File' : os.path.abspath(os.path.join(os.path.abspath(__file__), '../../../models/dart/rae_3.geo')), # Input file containing the model 'Pars' : {'spn' : 1.0, 'lgt' : 6.0, 'wdt' : 3.0, 'hgt' : 6.0, 'msN' : 0.02, 'msF' : 1}, # parameters for input file model 'Dim' : 3, # problem dimension 'Format' : 'vtk', # save format (vtk or gmsh) # Markers 'Fluid' : 'field', # name of physical group containing the fluid 'Farfield' : ['upstream', 'farfield', 'downstream'], # LIST of names of physical groups containing the farfield boundaries (upstream/downstream should be first/last element) 'Wings' : ['wing'], # LIST of names of physical groups containing the lifting surface boundary 'Wakes' : ['wake'], # LIST of names of physical group containing the wake 'WakeTips' : ['wakeTip'], # LIST of names of physical group containing the edge of the wake 'Tes' : ['te'], # LIST of names of physical group containing the trailing edge 'Symmetry': 'symmetry', # name of physical group containing the symmetry BC 'dbc' : True, 'Upstream' : 'upstream', # Freestream 'M_inf' : 0.8, # freestream Mach number 'AoA' : 0.0, # freestream angle of attack # Geometry 'S_ref' : 1.0, # reference surface length 'c_ref' : 1.0, # reference chord length 'x_ref' : 0.0, # reference point for moment computation (x) 'y_ref' : 0.0, # reference point for moment computation (y) 'z_ref' : 0.0, # reference point for moment computation (z) # Numerical 'LSolver' : 'GMRES', # inner solver (Pardiso, MUMPS or GMRES) 'G_fill' : 2, # fill-in factor for GMRES preconditioner 'G_tol' : 1e-5, # tolerance for GMRES 'G_restart' : 50, # restart for GMRES 'Rel_tol' : 1e-6, # relative tolerance on solver residual 'Abs_tol' : 1e-8, # absolute tolerance on solver residual 'Max_it' : 50 # solver maximum number of iterations } def cfgBlast(verb): return { 'Re' : 1e7, # Freestream Reynolds number 'Minf' : 0.8, # Freestream Mach number (used for the computation of the time step only) 'CFL0' : 1, # Inital CFL number of the calculation 'sections' : np.linspace(0.01, 0.95, 3), 'writeSections': [0.2, 0.4, 0.6, 0.8, 1.0], 'Sym':[0.], 'span': 1., 'interpolator': 'Rbf', 'rbftype': 'linear', 'smoothing': 1e-8, 'degree': 0, 'neighbors': 10, 'saveTag': 4, 'Verb': verb, # Verbosity level of the solver 'couplIter': 5, # Maximum number of iterations 'couplTol' : 5e-2, # Tolerance of the VII methodology 'iterPrint': 1, # int, number of iterations between outputs 'resetInv' : False, # bool, flag to reset the inviscid calculation at every iteration. 'xtrF' : [0., 0.], # Forced transition location 'nDim' : 3 } def main(): # Timer. tms = fwk.Timers() tms['total'].start() args = parseargs() icfg = cfgInviscid(args.k, args.verb) vcfg = cfgBlast(args.verb) parsViscous = {'spn': icfg['Pars']['spn'], 'lgt': icfg['Pars']['lgt'], 'nLe': 25, 'nMid': 50, 'nTe': 10, 'nSpan': 60, 'nWake': 25, 'progLe': 1.1, 'progMid': 1.0, 'progTe': 1.0, 'progSpan': 1.0, 'progWake': 1.15} vMeshFile = os.path.abspath(os.path.join(os.path.abspath(__file__), '../../../models/dart/rae_3_visc.geo')) vMsh = viscUtils.mesh(vMeshFile, parsViscous) vcfg['vMsh'] = vMsh tms['pre'].start() coupler, isol, vsol = viscUtils.initBlast(icfg, vcfg) tms['pre'].stop() print(ccolors.ANSI_BLUE + 'PySolving...' + ccolors.ANSI_RESET) tms['solver'].start() aeroCoeffs = coupler.run(write=False) tms['solver'].stop() # Display results. print(ccolors.ANSI_BLUE + 'PyRes...' + ccolors.ANSI_RESET) print(' Re M alpha Cl Cd Cd_wake Cdp Cdf Cm') print('{0:6.1f}e6 {1:8.2f} {2:8.1f} {3:8.4f} {4:8.4f} {5:8.4f} {6:8.4f} {7:8.4f} {8:8.4f}'.format(vcfg['Re']/1e6, isol.getMinf(), isol.getAoA()*180/math.pi, isol.getCl(), vsol.Cdf + isol.getCd(), vsol.Cdt, vsol.Cdp, vsol.Cdf, isol.getCm())) # Write results to file. vSolution = viscUtils.getSolution(isol.sec, write=True, toW='all') # Save results isol.save(sfx='_viscous') tms['total'].stop() print(ccolors.ANSI_BLUE + 'PyTiming...' + ccolors.ANSI_RESET) print('CPU statistics') print(tms) print('SOLVERS statistics') print(coupler.tms) # Test solution print(ccolors.ANSI_BLUE + 'PyTesting...' + ccolors.ANSI_RESET) tests = CTests() tests.add(CTest('Cl', isol.getCl(), 0.135, 5e-2)) tests.add(CTest('Cd wake', vsol.Cdt, 0.0074, 1e-3, forceabs=True)) tests.add(CTest('Cd integral', vsol.Cdf + isol.getCd(), 0.0140, 1e-3, forceabs=True)) tests.add(CTest('Cdf', vsol.Cdf, 0.0061, 1e-3, forceabs=True)) tests.add(CTest('Iterations', len(aeroCoeffs), 3, 0, forceabs=True)) tests.run() # eof print('') if __name__ == "__main__": main()