Fast: Flexible Aerodynamic Solver Technolgy
Preamble
Fast is a common/front module for all Fast series solvers.
Fast is only available for use with the pyTree interface. You must import the module:
import Fast.PyTree as Fast
List of functions
– Actions
| Set numeric data dictionary in bases. | |
| Set numeric data dictionary in zones. | |
| Load a single tree. | |
| Load all updated trees. | |
| Save a single tree. | |
| Save all updated trees. | 
Contents
Actions
- Fast.PyTree.setNum2Base(a, numb)
- Set the numb dictionary to bases. Exists also as in place version (_setNum2Base) that modifies a and returns None. - Parameters:
- a (Base, pyTree) – input data 
- numb (dictionary) – numerics for base 
 
 - the keys of numb dictionary are: - ‘temporal_scheme’: numerical scheme for time integration. Possible values are - ‘explicit’ (RK3 scheme, see p49 http://publications.onera.fr/exl-php/util/documents/accede_document.php) 
- ‘implicit’ (BDF2 (second order Gear scheme for global time stepping) or BDF1 (first order Euler scheme for local time stepping)) 
- ‘implicit_local’ (see p107 http://publications.onera.fr/exl-php/docs/ILS_DOC/227155/DOC356618_s1.pdf) 
- default value is ‘implicit’ 
 
- ‘ss_iteration’: Newton iterations for implicit schemes. - default value is 30 
 
- ‘modulo_verif’: computation period for cfl (RK3 or BDF2), Newton convergence (all), and predictor estimation (implicit_local only) - default value is 200 
 
 - Example of use: - # - setNum2Base (pyTree) - import Converter.PyTree as C import Converter.Internal as Internal import Fast.PyTree as Fast a = C.newPyTree(['Base']) numb = { 'temporal_scheme': 'explicit', 'ss_iteration': 30, 'modulo_verif':20 } Fast._setNum2Base(a, numb) Internal.printTree(a) 
- Fast.PyTree.setNum2Zones(a, numz)
- Set the numz dictionary to zones. Exists also as in place version (_setNum2Zones) that modifies a and returns None. - Parameters:
- a (Zone, Base, pyTree) – input data 
- numz (dictionary) – input data 
 
 - the keys of numz dictionary are: - ‘scheme’: numerical scheme for convective flux reconstruction. Possible values are - ‘ausmpred’ (AUSM+(P) reduced scheme, see p49 https://tel.archives-ouvertes.fr/pastel-00834850/document) 
- ‘senseur’ (AUSM+(P)-based hybrid centered/decentered sensor scheme for DNS/LES, see p50 https://tel.archives-ouvertes.fr/pastel-00834850/document) 
- ‘roe’ (classical Roe scheme) 
- default value is ‘ausmpred’ 
 
- ‘slope’: Slope reconstruction method. Possible values are - ‘o1’ (first order MUSCL reconstrution, only valid for roe scheme) 
- ‘minmod’ (second order MUSCL reconstrution with minmod limiter, only valid for roe scheme) 
- ‘o3’ (third order MUSCL reconstrution, see p50 https://tel.archives-ouvertes.fr/pastel-00834850/document) 
- ‘o3sc’ (third order MUSCL reconstrution with slope limiter, see - Lugrin scheme)
- ‘o5’ (fifth order MUSCL reconstrution) 
- ‘o5sc’ (fifth order MUSCL reconstrution with slope limiter) 
- default value is ‘o3’ 
 
- ‘senseurType’: Sensor mode for the ‘senseur’ scheme. Possible values are - 0 (correction for speed only) 
- 1 (correction for speed, density and pressure) 
 
- ‘motion’: Motion mode for moving walls. Possible values are - ‘none’ (no motion) 
- ‘rigid’ (rigid motion defined in Fast, see p47 https://tel.archives-ouvertes.fr/tel-01011273/document) 
- ‘rigid_ext’ (rigid motion defined externally by RigidMotion) 
- ‘deformation’ (ALE with deformation) 
- default value is ‘none’ 
 
- ‘time_step’: time step value - default value is 1e-4 
 
- ‘time_step_nature’: time step nature. Possible values are - ‘global’ 
- ‘local’ 
- default value is ‘global’ 
 
- ‘cfl’: Courant-Friedrichs-Lewy condition (only for time_step_nature’=’local’) - default value is 1 
 
- ‘epsi_newton’: Newton stopping criteria based on Loo norm - default value is 0.1 
 
- ‘inj1_newton_tol’: Newton tolerence for BCinj1 inflow condition - default value is 1e-5 
 
- ‘inj1_newton_nit’: Newton iterations for BCinj1 inflow condition - default value is 10 
 
- ‘psiroe’: Harten correction parameter - default value is 0.1 
 
- ‘prandtltb’: turbulent Prandtl number (only active for ‘model’=’NSTurbulent’) - default value is 0.92 
 
- ‘ransmodel’: RANS turbulence models. Possible values are - ‘SA’ (Standard Spalart-Allmaras model) 
- ‘SA_comp’ (SA with mixing layer compressibility correction, see https://turbmodels.larc.nasa.gov/spalart.html) 
- default value is ‘SA’ 
 
- ‘SA_add_RotCorr’: add rotation correction to the SA model (SA-R) - 0 
- 1 (toggle rotation correction) 
- default value is 0 
 
- ‘SA_add_LowRe’: add low Reynolds correction to the SA model (SA-LRe) - 0 
- 1 (toggle low-Re correction) 
- default value is 0 
 
- ‘DES’: toggle the turbulence models for DES simulations. Possible values are - ‘none’ (SA computation) 
- ‘zdes1’ (ZDES mode 1, see https://link.springer.com/content/pdf/10.1007%2Fs00162-011-0240-z.pdf) 
- ‘zdes1_w’ (ZDES mode 1 by Chauvet) 
- ‘zdes2’ (ZDES mode 2) 
- ‘zdes2_w’ (ZDES mode 2 by Chauvet) 
- ‘zdes3’ (ZDES mode 3, see p118 https://tel.archives-ouvertes.fr/tel-01365361/document) 
- default value is ‘none’ 
 
- ‘DES_debug’: toggle data extraction for DES post-analysis. Possible values are - 0 
- 1 (save delta and fd functions in the FlowSolution#Centers node) 
- default value is ‘none’ 
 
- ‘sgsmodel’: LES turbulence models. Possible values are - ‘Miles’ (Monotonically integrated LES, for which ViscosityEddy==LaminarViscosity) 
- ‘smsm’ (Selective mixed scale model, see Lenormand et al https://doi.org/10.1002/(SICI)1097-0363(20000229)32:4%3C369::AID-FLD943%3E3.0.CO;2-6) 
- default value is ‘Miles’ 
 
- ‘extract_res’: toggle residual extraction during simulations. Possible values are - 0 (don’t save) 
- 1 (save div(F_Euler-F_viscous) in the FlowSolution#Centers node) 
- 2 (save dqdt + div(F_Euler-F_viscous) in the FlowSolution#Centers node) 
- default value is 0 
 
- ‘source’: toggle a source term for CFD simulations. Possible values are - 0 
- 1 (read a source terme in the FlowSolution#Centers node. The conservative variables centers:Density_src, centers:MomentumX_src, centers:MomentumY_src, centers:MomentumZ_src and centers:EnergyStagnationDensity_src are used.) 
- default value is 0 
 
- ‘ratiom’: maximum cut-off for mut/mu - default value is 10000 
 
 - Example of use: - # - setNum2Zones (pyTree) - import Converter.PyTree as C import Converter.Internal as Internal import Fast.PyTree as Fast import Generator.PyTree as G a = G.cart((0,0,0), (1,1,1), (10,10,10) ) t = C.newPyTree(['Base', a]) numz = { 'scheme': 'ausmpred', 'time_step_nature': 'global', 'time_step': 0.001 } Fast._setNum2Zones(t, numz) Internal.printTree(t) 
- Fast.PyTree.loadTree(fileName='t.cgns', split='single', graph=False)
- Load a tree from a file (solution tree t, connectivity tree tc, or statistics tree ts). The file can be a single CGNS file (e.g. ‘t.cgns’), or can be split with one file per processor (e.g. ‘t/t_1.cgns’, ‘t/t_2.cgns’, etc.). - When executed in sequential mode, the function returns a complete tree. - When executed in mpi mode, the function returns a partial tree on each processor. - Warning - in MPI mode, the tree must already be distributed with the correct number of processors (with ‘proc’ nodes within the zones). - If graph is true, create and return the communication graph for Chimera and abutting transfers (for MPI use only). The dictionary is of the form of: graph={‘graphID’:[], ‘graphIBC’:[], ‘procDict’:[], ‘procList’:[]}. - Warning - to create the communication graph, the loaded tree must be the connectivity tree (tc.cgns), not the solution tree (t.cgns). - Parameters:
- fileName (string) – name of the file to load 
- split (string) – file format: ‘single’ or ‘multiple’ 
- graph (boolean) – true for returning communication graph 
 
- Returns:
- t or tc or ts or (tc, graph) 
 - # - loadFile (pyTree) - import Fast.PyTree as Fast import Converter.PyTree as C import Generator.PyTree as G import Converter.Mpi as Cmpi # Save a file a = G.cart((0,0,0),(1,1,1),(11,11,11)) a = Cmpi.setProc(a,0) b = G.cart((10,0,0),(1,1,1),(11,11,11)) b = Cmpi.setProc(b,1) t = C.newPyTree(['Base',a,b]) Fast.saveTree(t, fileName='t.cgns', split='single') # Load it back Fast.loadTree(fileName='t.cgns', split='single') - Note - since Fast 4.1, loadFile and loadTree have been merged. 
- Fast.PyTree.saveTree(t, fileName='restart.cgns', split='single', compress=0)
- Save a tree to a file (solution tree t, connectivity tree tc, or statistics tree ts). The file can be a single CGNS file (e.g. ‘t.cgns’), or can be split with one file per processor (e.g. ‘t/t_1.cgns’, ‘t/t_2.cgns’, etc.). The tree is also cleared of unnecessary flow fields and temporary work arrays. - Warning - in MPI mode, the tree must already be distributed with the correct number of processors (with ‘proc’ nodes within the zones). - Before saving the file, the flow fields can be compressed using the compress parameter. - Parameters:
- t (tree) – tree to save 
- fileName (string) – name of the file to load 
- split (string) – file format: ‘single’ or ‘multiple’ 
- compress (integer) – compression parameter: 0 (none), 1 (cart only), or 2 (all) 
 
 - # - saveFile (pyTree) - import Fast.PyTree as Fast import Converter.PyTree as C import Generator.PyTree as G import Converter.Mpi as Cmpi a = G.cart((0,0,0),(1,1,1),(11,11,11)) a = Cmpi.setProc(a,0) b = G.cart((10,0,0),(1,1,1),(11,11,11)) b = Cmpi.setProc(b,1) t = C.newPyTree(['Base',a,b]) Fast.saveTree(t, fileName='t.cgns', split='single') - Note - since Fast 4.1, saveFile and saveTree have been merged. 
- Fast.PyTree.load(fileName='t.cgns', fileNameC=None, fileNameS=None, split='single')
- Load all trees at once (solution tree t, connectivity tree tc, and/or statistics tree ts), using multiple calls of the Fast.PyTree.loadTree function. The files can be single CGNS files (e.g. ‘t.cgns’), or can be split with one file per processor (e.g. ‘t/t_1.cgns’, ‘t/t_2.cgns’, etc.). - Warning - in MPI mode, the trees must already be distributed with the correct number of processors (with ‘proc’ nodes within the zones). - If fileNameC is not None, the connectivity tree is loaded and the communication graph is automatically created. - If filenameS is not None, the statistics tree is loaded. - Parameters:
- fileName (string) – name of the solution file to load 
- fileNameC (string) – name of the connectivity file to load 
- fileNameS (string) – name of the statistics file to load 
- split (string) – file format: ‘single’ or ‘multiple’ 
 
- Returns:
- (t, tc, ts, graph) 
 - # - load (pyTree) - import Fast.PyTree as Fast import Converter.PyTree as C import Generator.PyTree as G import Converter.Internal as Internal import Connector.PyTree as X import Converter.Mpi as Cmpi # Create case if Cmpi.rank == 0: a = G.cart((0,0,0),(1,1,1),(11,11,11)) a = Cmpi.setProc(a,0) b = G.cart((10,0,0),(1,1,1),(11,11,11)) b = Cmpi.setProc(b,1) t = C.newPyTree(['Base',a,b]) t = X.connectMatch(t, dim=3) t = Internal.addGhostCells(t,t,2,adaptBCs=0) C.convertPyTree2File(t, 't.cgns') tc = C.node2Center(t) tc = X.setInterpData(t, tc, loc='centers', storage='inverse') C.convertPyTree2File(tc, 'tc.cgns') Cmpi.barrier() t,tc,ts,graph = Fast.load('t.cgns', 'tc.cgns', split='single') #print(t) #print(tc) #print(ts) #print('procDict = ', graph['procDict']) #print('graphID = ', graph['graphID']) 
- Fast.PyTree.save(t, fileName='restart.cgns', tc=None, fileNameC='tc_restart.cgns', ts=None, fileNameS='tstat_restart.cgns', split='single', compress=0)
- Save all trees at once (solution tree t, connectivity tree tc, and/or statistics tree ts), using multiple calls of the Fast.PyTree.saveTree function. The files can be single CGNS files (e.g. ‘t.cgns’), or can be split with one file per processor (e.g. ‘t/t_1.cgns’, ‘t/t_2.cgns’, etc.). The trees are also cleared of unnecessary flow fields and temporary work arrays. - Warning - in MPI mode, the trees must already be distributed with the correct number of processors (with ‘proc’ nodes within the zones). - If tc is not None, the connectivity tree is saved. - If ts is not None, the statistics tree is saved. - Before saving the file, the flow fields can be compressed using the compress parameter. - Parameters:
- t (tree) – updated solution tree to save 
- fileName (string) – name of the solution file 
- tc (tree) – updated connectivity tree to save 
- fileNameC (string) – name of the connectivity file 
- ts (tree) – updated statistics tree to save 
- fileNameS (string) – name of the statistics file 
- split (string) – file format: ‘single’ or ‘multiple’ 
- compress (integer) – compression parameter: 0 (none), 1 (cart only), or 2 (all) 
 
 - Example of use: - # - save (pyTree) - import Fast.PyTree as Fast import Converter.PyTree as C import Generator.PyTree as G import Converter.Mpi as Cmpi a = G.cart((0,0,0),(1,1,1),(11,11,11)) a = Cmpi.setProc(a,0) b = G.cart((10,0,0),(1,1,1),(11,11,11)) b = Cmpi.setProc(b,1) t = C.newPyTree(['Base',a,b]) tc = C.node2Center(t) Fast.save(t, fileName='t.cgns', tc=tc, fileNameC='tc.cgns', split='single') 
