|
| 1 | +#!/usr/bin/env python |
| 2 | + |
| 3 | +## \file run.py |
| 4 | +# \brief Deforming bump in channel. |
| 5 | +# \version 7.5.1 "Blackbird" |
| 6 | +# |
| 7 | +# SU2 Project Website: https://su2code.github.io |
| 8 | +# |
| 9 | +# The SU2 Project is maintained by the SU2 Foundation |
| 10 | +# (http://su2foundation.org) |
| 11 | +# |
| 12 | +# Copyright 2012-2023, SU2 Contributors (cf. AUTHORS.md) |
| 13 | +# |
| 14 | +# SU2 is free software; you can redistribute it and/or |
| 15 | +# modify it under the terms of the GNU Lesser General Public |
| 16 | +# License as published by the Free Software Foundation; either |
| 17 | +# version 2.1 of the License, or (at your option) any later version. |
| 18 | +# |
| 19 | +# SU2 is distributed in the hope that it will be useful, |
| 20 | +# but WITHOUT ANY WARRANTY; without even the implied warranty of |
| 21 | +# MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU |
| 22 | +# Lesser General Public License for more details. |
| 23 | +# |
| 24 | +# You should have received a copy of the GNU Lesser General Public |
| 25 | +# License along with SU2. If not, see <http://www.gnu.org/licenses/>. |
| 26 | + |
| 27 | +import pysu2 |
| 28 | +from mpi4py import MPI |
| 29 | +import numpy as np |
| 30 | + |
| 31 | +def main(): |
| 32 | + comm = MPI.COMM_WORLD |
| 33 | + rank = comm.Get_rank() |
| 34 | + |
| 35 | + # Initialize the corresponding driver of SU2, this includes solver preprocessing. |
| 36 | + try: |
| 37 | + SU2Driver = pysu2.CSinglezoneDriver('config.cfg', 1, comm) |
| 38 | + except TypeError as exception: |
| 39 | + print('A TypeError occured in pysu2.CDriver : ', exception) |
| 40 | + raise |
| 41 | + |
| 42 | + # Get the ID of the marker we want to deform. |
| 43 | + AllMarkerIDs = SU2Driver.GetMarkerIndices() |
| 44 | + MarkerName = 'interface' |
| 45 | + MarkerID = AllMarkerIDs[MarkerName] if MarkerName in AllMarkerIDs else -1 |
| 46 | + |
| 47 | + # Number of vertices on the specified marker (per rank). |
| 48 | + nVertex = SU2Driver.GetNumberMarkerNodes(MarkerID) if MarkerID >= 0 else 0 |
| 49 | + |
| 50 | + # Retrieve some control parameters from the driver. |
| 51 | + deltaT = SU2Driver.GetUnsteady_TimeStep() |
| 52 | + TimeIter = SU2Driver.GetTime_Iter() |
| 53 | + nTimeIter = SU2Driver.GetnTimeIter() |
| 54 | + time = TimeIter * deltaT |
| 55 | + |
| 56 | + # Extract the initial position of each node on the moving marker. |
| 57 | + CoordX = np.zeros(nVertex) |
| 58 | + CoordY = np.zeros(nVertex) |
| 59 | + for iVertex in range(nVertex): |
| 60 | + iPoint = SU2Driver.GetMarkerNode(MarkerID, iVertex) |
| 61 | + CoordX[iVertex], CoordY[iVertex] = SU2Driver.GetInitialCoordinates(iPoint) |
| 62 | + |
| 63 | + if rank == 0: |
| 64 | + print("\n------------------------------ Begin Solver -----------------------------\n") |
| 65 | + |
| 66 | + # The time loop is defined in Python so that we have acces to SU2 functionalities at each time step. |
| 67 | + while (TimeIter < nTimeIter): |
| 68 | + # Apply the surface deformation. |
| 69 | + for iVertex in range(nVertex): |
| 70 | + dy = np.real(DeformFunction(CoordX[iVertex] - 0.9, time)) |
| 71 | + SU2Driver.SetMarkerDisplacements(MarkerID, iVertex, (0.0, dy)) |
| 72 | + |
| 73 | + # Time iteration preprocessing. |
| 74 | + SU2Driver.Preprocess(TimeIter) |
| 75 | + |
| 76 | + # Run one time iteration (e.g. dual-time). |
| 77 | + SU2Driver.Run() |
| 78 | + SU2Driver.Postprocess() |
| 79 | + SU2Driver.Update() |
| 80 | + |
| 81 | + # Monitor the solver and output solution to file if required. |
| 82 | + stopCalc = SU2Driver.Monitor(TimeIter) |
| 83 | + SU2Driver.Output(TimeIter) |
| 84 | + |
| 85 | + if (stopCalc == True): |
| 86 | + break |
| 87 | + |
| 88 | + # Update control parameters |
| 89 | + TimeIter += 1 |
| 90 | + time += deltaT |
| 91 | + |
| 92 | + # Postprocess the solver and exit cleanly |
| 93 | + SU2Driver.Postprocessing() |
| 94 | + |
| 95 | + |
| 96 | +# Imposed deformation |
| 97 | +def DeformFunction(x,t): |
| 98 | + A1 = 0.016 |
| 99 | + L = 1.0 |
| 100 | + k1 = 4.730/L |
| 101 | + return A1*(np.cosh(k1*x) - np.cos(k1*x) + ((np.cosh(k1*L)-np.cos(k1*L))*(np.sin(k1*x)-np.sinh(k1*x)))/(np.sinh(k1*L) - np.sin(k1*L)))*(np.exp(1j*t*2) + np.exp(-1j*t*2)) |
| 102 | + |
| 103 | + |
| 104 | +if __name__ == '__main__': |
| 105 | + main() |
| 106 | + |
0 commit comments