|
| 1 | +--- |
| 2 | +title: Static Fluid-Structure Interaction (FSI) via Python-wrapper |
| 3 | +permalink: /tutorials/Static_FSI_Python/ |
| 4 | +written_by: rsanfer |
| 5 | +for_version: 7.0.3 |
| 6 | +revised_by: rsanfer |
| 7 | +revision_date: 2020-03-04 |
| 8 | +revised_version: 7.0.3 |
| 9 | +solver: Multiphysics |
| 10 | +requires: pysu2 library |
| 11 | +complexity: advanced |
| 12 | +follows: Static_FSI |
| 13 | +userguide: Build-SU2-Linux-MacOS |
| 14 | +--- |
| 15 | + |
| 16 | +### Goals |
| 17 | + |
| 18 | +This tutorial uses SU2's python wrapper and its native solvers for incompressible flow and solid mechanics to solve a steady-state Fluid-Structure Interaction problem. This document will cover: |
| 19 | +- Setting up a python script to run SU2 |
| 20 | +- Extracting flow loads and structural displacements from two different python instances of SU2 |
| 21 | +- Coupling the two instances |
| 22 | + |
| 23 | +In this tutorial, we will solve the exact same problem as for the [Static FSI](../Static_FSI/) tutorial. Please take a moment to familiarize yourself with that tutorial in case you haven't yet. |
| 24 | + |
| 25 | +### Resources |
| 26 | + |
| 27 | +You can find the resources for this tutorial in the folder [python_fsi](https://github.com/su2code/Tutorials/tree/feature_python_fsi/multiphysics/python_fsi) of the [Tutorials repository](https://github.com/su2code/Tutorials). There is a [python script](https://github.com/su2code/Tutorials/tree/feature_python_fsi/multiphysics/python_fsi/run_fsi_primal.py) and two sub-config files for the [flow](https://github.com/su2code/Tutorials/tree/feature_python_fsi/multiphysics/python_fsi/config_channel.cfg) and [structural](https://github.com/su2code/Tutorials/tree/feature_python_fsi/multiphysics/python_fsi/config_cantilever.cfg) subproblems. |
| 28 | + |
| 29 | +Moreover, you will need two mesh files for the [flow domain](https://github.com/su2code/Tutorials/tree/feature_python_fsi/multiphysics/python_fsi/mesh_channel.su2) and the [cantilever](https://github.com/su2code/Tutorials/tree/feature_python_fsi/multiphysics/python_fsi/mesh_cantilever.su2). |
| 30 | + |
| 31 | +### Background |
| 32 | + |
| 33 | +For this tutorial, you will need to use advanced features of SU2, in particular the python-wrapped version of the code which needs to be built from source. |
| 34 | + |
| 35 | +This tutorial has been tested on a linux system with the following specs |
| 36 | + |
| 37 | +``` |
| 38 | +Linux kernel 5.3.18-1-MANJARO |
| 39 | +GCC compilers version 9.2.0 |
| 40 | +Open MPI version 4.0.2 |
| 41 | +Python version 3.8.1 |
| 42 | +SWIG version 4.0.1 |
| 43 | +``` |
| 44 | + |
| 45 | +compiling the code from source using the following meson settings |
| 46 | + |
| 47 | +``` |
| 48 | +./meson.py build -Dwith-mpi=enabled -Denable-autodiff=true -Denable-pywrapper=true |
| 49 | +``` |
| 50 | + |
| 51 | +It requires an adequate setup of the system and a correct linkage of both the ```mpi4py``` and ```swig``` libraries. For questions, updates, or notes on the usability of this tutorials on different systems or configurations, please use the comment section below. |
| 52 | + |
| 53 | +#### Mesh Description |
| 54 | + |
| 55 | +The cantilever is discretized using 1000 4-node quad elements with linear interpolation. The fluid domain is discretized using a finite volume mesh with 7912 nodes and 15292 triangular volumes. The wet interface is matching between the fluid and structural domains. |
| 56 | + |
| 57 | +#### Configuration File Options |
| 58 | + |
| 59 | +We reuse the [flow](https://github.com/su2code/Tutorials/tree/feature_python_fsi/multiphysics/python_fsi/config_channel.cfg) and [structural](https://github.com/su2code/Tutorials/tree/feature_python_fsi/multiphysics/python_fsi/config_cantilever.cfg) config files from the [Static FSI](../Static_FSI/), with only very minor modifications. |
| 60 | + |
| 61 | +First of all, SU2 will see each instance as a single-zone problem. Therefore, it is necessary to set the output appropriately in order to prevent overwriting files. For the flow config |
| 62 | + |
| 63 | +``` |
| 64 | +OUTPUT_FILES = (RESTART, PARAVIEW) |
| 65 | +SOLUTION_FILENAME = solution_fsi_steady_0 |
| 66 | +RESTART_FILENAME = restart_fsi_steady_0 |
| 67 | +VOLUME_FILENAME = fsi_steady_0 |
| 68 | +``` |
| 69 | + |
| 70 | +and for the structural config |
| 71 | + |
| 72 | +``` |
| 73 | +OUTPUT_FILES = (RESTART, PARAVIEW) |
| 74 | +SOLUTION_FILENAME = solution_fsi_steady_1 |
| 75 | +RESTART_FILENAME = restart_fsi_steady_1 |
| 76 | +VOLUME_FILENAME = fsi_steady_1 |
| 77 | +``` |
| 78 | + |
| 79 | +Next, we choose the output that we require to the history files. In this tutorial, we are interested on the drag coefficient of the cantilever once deformed. Therefore, for the flow config, we need to include the ```AERO_COEFF``` |
| 80 | + |
| 81 | +``` |
| 82 | +HISTORY_OUTPUT = ITER, RMS_RES, AERO_COEFF |
| 83 | +CONV_FILENAME= history_0 |
| 84 | +``` |
| 85 | + |
| 86 | +while for the structural config it's enough with the residuals |
| 87 | + |
| 88 | +``` |
| 89 | +HISTORY_OUTPUT = ITER, RMS_RES |
| 90 | +CONV_FILENAME= history_1 |
| 91 | +``` |
| 92 | + |
| 93 | +Finally, it is necessary to indicate the structural solver that the boundary solution will be used to update a fluid field. This is done by setting the config options |
| 94 | + |
| 95 | +``` |
| 96 | +MARKER_FLUID_LOAD = ( feabound ) |
| 97 | +MARKER_DEFORM_MESH = ( feabound ) |
| 98 | +``` |
| 99 | + |
| 100 | +on the structural config file. |
| 101 | + |
| 102 | +#### Applying coupling conditions to the individual domains |
| 103 | + |
| 104 | +The key part of this tutorial is the [python script](https://github.com/su2code/Tutorials/tree/feature_python_fsi/multiphysics/python_fsi/run_fsi_primal.py) to run the FSI problem. Please take a moment to evaluate its contents, as we will go through some of its most important aspects. |
| 105 | + |
| 106 | +First, we will need to import the SU2 python library along with mpi4py, using |
| 107 | + |
| 108 | +``` |
| 109 | +import pysu2 |
| 110 | +from mpi4py import MPI |
| 111 | +``` |
| 112 | + |
| 113 | +We define the names of the config files required by SU2 using |
| 114 | + |
| 115 | +``` |
| 116 | +flow_filename = "config_channel.cfg" |
| 117 | +fea_filename = "config_cantilever.cfg" |
| 118 | +``` |
| 119 | + |
| 120 | +We will exemplify the initialization of SU2 using the flow domain. First, we create a driver object using |
| 121 | + |
| 122 | +``` |
| 123 | +FlowDriver = pysu2.CSinglezoneDriver(flow_filename, 1, comm); |
| 124 | +``` |
| 125 | + |
| 126 | +which is analogous to the SU2 driver in the C++ executable. We use a ```comm``` that is imported from ```mpi4py``` |
| 127 | + |
| 128 | +``` |
| 129 | +comm = MPI.COMM_WORLD |
| 130 | +``` |
| 131 | + |
| 132 | +and the ```flow_filename``` previously defined. The ```1``` corresponds to the number of zones (remember, we consider each instance of SU2 to be independent). |
| 133 | + |
| 134 | +Next, we need to identify the FSI boundary. We retrieve the list of all the flow boundary tags and their associated indices using |
| 135 | + |
| 136 | +``` |
| 137 | +FlowMarkerList = FlowDriver.GetAllBoundaryMarkersTag() |
| 138 | +FlowMarkerIDs = FlowDriver.GetAllBoundaryMarkers() |
| 139 | +``` |
| 140 | + |
| 141 | +For the boundary ```FlowMarkerName = 'flowbound' ```, we identify the corresponding index |
| 142 | + |
| 143 | +``` |
| 144 | +if FlowMarkerName in FlowMarkerList and FlowMarkerName in FlowMarkerIDs.keys(): |
| 145 | + FlowMarkerID = FlowMarkerIDs[FlowMarkerName] |
| 146 | +``` |
| 147 | + |
| 148 | +and, finally, the number of vertices on the boundary is |
| 149 | + |
| 150 | +``` |
| 151 | +nVertex_Marker_Flow = FlowDriver.GetNumberVertices(FlowMarkerID) |
| 152 | +``` |
| 153 | + |
| 154 | +An analogous process can be followed for the FSI boundary on the structural domain. In this case, the boundaries are matching, so the number of vertices is the same. |
| 155 | + |
| 156 | +We start the FSI loop and limit it to the same number of iterations required for convergence in the [Static FSI](../Static_FSI/) tutorial |
| 157 | + |
| 158 | +``` |
| 159 | +for i in range(17): |
| 160 | +``` |
| 161 | + |
| 162 | +First, the flow solution is run |
| 163 | + |
| 164 | +``` |
| 165 | + FlowDriver.ResetConvergence() |
| 166 | + FlowDriver.Preprocess(0) |
| 167 | + FlowDriver.Run() |
| 168 | + FlowDriver.Postprocess() |
| 169 | + stopCalc = FlowDriver.Monitor(0) |
| 170 | +``` |
| 171 | + |
| 172 | +and the flow loads are recovered using |
| 173 | + |
| 174 | +``` |
| 175 | + flow_loads=[] |
| 176 | + for iVertex in range(nVertex_Marker_Flow): |
| 177 | + vertexLoad = FlowDriver.GetFlowLoad(FlowMarkerID, iVertex) |
| 178 | + flow_loads.append(vertexLoad) |
| 179 | +``` |
| 180 | + |
| 181 | +The ```flow_loads``` array now contains the loads in all the vertices of the flow FSI interface. Now, we need to set the flow loads to the FEA nodes. By construction for this case, the vertex IDs are matching for both meshes except for vertex 0 and 1, which are inverted. Therefore, we set the flow loads on the structural domain using |
| 182 | + |
| 183 | +``` |
| 184 | + FEADriver.SetFEA_Loads(FEAMarkerID, 0, flow_loads[1][0], flow_loads[1][1], flow_loads[1][2]) |
| 185 | + FEADriver.SetFEA_Loads(FEAMarkerID, 1, flow_loads[0][0], flow_loads[0][1], flow_loads[0][2]) |
| 186 | + for iVertex in range(2, nVertex_Marker_FEA): |
| 187 | + FEADriver.SetFEA_Loads(FEAMarkerID, iVertex, flow_loads[iVertex][0], flow_loads[iVertex][1], flow_loads[iVertex][2]) |
| 188 | +``` |
| 189 | + |
| 190 | +You can ensure the vertices are coincidental by the coordinates of the nodes using ```FlowDriver.GetVertexCoordX(FlowMarkerID, iVertex)``` and ```FlowDriver.GetVertexCoordY(FlowMarkerID, iVertex)``` for the flow domain, and ```FEADriver.GetVertexCoordX(FEAMarkerID, iVertex)``` and ```FEADriver.GetVertexCoordY(FEAMarkerID, iVertex)``` for the structural domain. |
| 191 | + |
| 192 | +Next, the structural simulation is run with |
| 193 | + |
| 194 | +``` |
| 195 | + FEADriver.ResetConvergence() |
| 196 | + FEADriver.Preprocess(0) |
| 197 | + FEADriver.Run() |
| 198 | + FEADriver.Postprocess() |
| 199 | + stopCalc = FEADriver.Monitor(0) |
| 200 | +``` |
| 201 | + |
| 202 | +and the structural displacements at the ```feabound``` interface are retrieved using |
| 203 | + |
| 204 | +``` |
| 205 | + fea_disp=[] |
| 206 | + for iVertex in range(nVertex_Marker_FEA): |
| 207 | + vertexDisp = FEADriver.GetFEA_Displacements(FEAMarkerID, iVertex) |
| 208 | + fea_disp.append(vertexDisp) |
| 209 | +``` |
| 210 | + |
| 211 | +Finally, these boundary displacements are imposed to the flow mesh |
| 212 | + |
| 213 | +``` |
| 214 | + FlowDriver.SetMeshDisplacement(FlowMarkerID, 0, fea_disp[1][0], fea_disp[1][1], fea_disp[1][2]) |
| 215 | + FlowDriver.SetMeshDisplacement(FlowMarkerID, 1, fea_disp[0][0], fea_disp[0][1], fea_disp[0][2]) |
| 216 | + for iVertex in range(2, nVertex_Marker_FEA): |
| 217 | + FlowDriver.SetMeshDisplacement(FlowMarkerID, iVertex, fea_disp[iVertex][0], fea_disp[iVertex][1], fea_disp[iVertex][2]) |
| 218 | +``` |
| 219 | + |
| 220 | +Once the loop is completed, it only remains to write the solution of each domain to file using |
| 221 | + |
| 222 | +``` |
| 223 | +FlowDriver.Output(0) |
| 224 | +FEADriver.Output(0) |
| 225 | +``` |
| 226 | + |
| 227 | +### Running SU2 |
| 228 | + |
| 229 | +Follow the links provided to download the [python script](https://github.com/su2code/Tutorials/tree/feature_python_fsi/multiphysics/python_fsi/run_fsi_primal.py) and the two sub-config files for the [flow](https://github.com/su2code/Tutorials/tree/feature_python_fsi/multiphysics/python_fsi/config_channel.cfg) and [structural](https://github.com/su2code/Tutorials/tree/feature_python_fsi/multiphysics/python_fsi/config_cantilever.cfg) subproblems. |
| 230 | + |
| 231 | +Also, you will need the two mesh files for the [flow domain](https://github.com/su2code/Tutorials/tree/feature_python_fsi/multiphysics/python_fsi/mesh_channel.su2) and the [cantilever](https://github.com/su2code/Tutorials/tree/feature_python_fsi/multiphysics/python_fsi/mesh_cantilever.su2). |
| 232 | + |
| 233 | +Execute the code using Python |
| 234 | + |
| 235 | +``` |
| 236 | +$ python run_fsi_primal.py |
| 237 | +``` |
| 238 | + |
| 239 | +The convergence history of each individual domain will be printed to screen, starting with the flow domain |
| 240 | + |
| 241 | +``` |
| 242 | ++---------------------------------------------------+ |
| 243 | +| Inner_Iter| rms[P]| rms[U]| rms[V]| |
| 244 | ++---------------------------------------------------+ |
| 245 | +| 0| -4.799717| -19.870540| -32.000000| |
| 246 | +| 10| -5.421646| -4.859176| -5.558343| |
| 247 | +| 20| -6.200459| -5.584293| -6.120761| |
| 248 | +| 30| -6.852188| -6.205093| -6.958953| |
| 249 | +| 40| -7.502751| -6.854100| -8.034484| |
| 250 | +| 50| -8.316190| -7.662838| -8.492639| |
| 251 | +| 60| -9.153994| -8.509430| -9.696622| |
| 252 | +| 70| -10.119472| -9.498203| -10.117162| |
| 253 | +| 80| -11.035958| -10.448171| -11.143732| |
| 254 | +| 90| -11.613994| -10.982943| -11.842820| |
| 255 | +| 100| -12.264740| -11.616077| -12.556728| |
| 256 | +| 104| -12.684964| -12.046580| -12.778620| |
| 257 | +``` |
| 258 | + |
| 259 | +followed by the structural domain |
| 260 | + |
| 261 | +``` |
| 262 | ++----------------------------------------------------------------+ |
| 263 | +| Inner_Iter| rms[U]| rms[R]| rms[E]| VonMises| |
| 264 | ++----------------------------------------------------------------+ |
| 265 | +| 0| -1.121561| -2.706864| -4.419338| 2.2876e+03| |
| 266 | +| 1| -1.708123| 1.269382| -1.716845| 2.2429e+03| |
| 267 | +| 2| -2.583705| 0.631041| -3.175890| 2.2237e+03| |
| 268 | +| 3| -2.508726| -0.795411| -5.686295| 2.0701e+03| |
| 269 | +| 4| -2.751613| -1.306797| -6.550284| 2.0544e+03| |
| 270 | +| 5| -3.147799| -1.271416| -6.914834| 2.0364e+03| |
| 271 | +| 6| -2.962716| -2.664273| -7.800357| 2.0171e+03| |
| 272 | +| 7| -4.083257| -2.080066| -8.553973| 2.0155e+03| |
| 273 | +| 8| -4.522566| -4.400710| -11.005428| 2.0149e+03| |
| 274 | +| 9| -7.252524| -5.240218| -14.876499| 2.0149e+03| |
| 275 | +| 10| -10.834476| -10.711888| -23.634504| 2.0149e+03| |
| 276 | +``` |
| 277 | + |
| 278 | +After 17 iterations, both the flow and structural fields will be successfully converged, |
| 279 | + |
| 280 | +``` |
| 281 | ++---------------------------------------------------+ |
| 282 | +| Inner_Iter| rms[P]| rms[U]| rms[V]| |
| 283 | ++---------------------------------------------------+ |
| 284 | +| 0| -15.808168| -15.358431| -15.769282| |
| 285 | +| 10| -16.714759| -16.074145| -17.174991| |
| 286 | +
|
| 287 | ++----------------------------------------------------------------+ |
| 288 | +| Inner_Iter| rms[U]| rms[R]| rms[E]| VonMises| |
| 289 | ++----------------------------------------------------------------+ |
| 290 | +| 0| -11.937375| -11.637725| -25.999458| 1.8531e+03| |
| 291 | +| 1| -15.794016| -11.636749| -28.541279| 1.8531e+03| |
| 292 | +| 2| -16.165350| -11.621837| -28.528692| 1.8531e+03| |
| 293 | +| 3| -16.056047| -11.608302| -28.544055| 1.8531e+03| |
| 294 | +| 4| -16.156344| -11.627675| -28.565975| 1.8531e+03| |
| 295 | +| 5| -16.204699| -11.633976| -28.581277| 1.8531e+03| |
| 296 | +``` |
| 297 | + |
| 298 | +Which correspond to the exact same residuals as for the last iteration on the [Static FSI](../Static_FSI/) tutorial. This can be tested by running the latter with |
| 299 | + |
| 300 | +``` |
| 301 | +WRT_ZONE_CONV = YES |
| 302 | +``` |
| 303 | + |
| 304 | +which yields for its last iteration |
| 305 | + |
| 306 | +``` |
| 307 | ++----------------------------------------------------------------+ |
| 308 | +| Zone 0 (Incomp. Fluid) | |
| 309 | ++----------------------------------------------------------------+ |
| 310 | +| Outer_Iter| Inner_Iter| rms[P]| rms[U]| rms[V]| |
| 311 | ++----------------------------------------------------------------+ |
| 312 | +| 16| 0| -15.808168| -15.358431| -15.769282| |
| 313 | +| 16| 10| -16.714759| -16.074145| -17.174991| |
| 314 | ++-----------------------------------------------------------------------------+ |
| 315 | +| Zone 1 (Structure) | |
| 316 | ++-----------------------------------------------------------------------------+ |
| 317 | +| Outer_Iter| Inner_Iter| rms[U]| rms[R]| rms[E]| VonMises| |
| 318 | ++-----------------------------------------------------------------------------+ |
| 319 | +| 16| 0| -11.937375| -11.637725| -25.999458| 1.8531e+03| |
| 320 | +| 16| 1| -15.794016| -11.636749| -28.541279| 1.8531e+03| |
| 321 | +| 16| 2| -16.165350| -11.621837| -28.528692| 1.8531e+03| |
| 322 | +| 16| 3| -16.056047| -11.608302| -28.544055| 1.8531e+03| |
| 323 | +| 16| 4| -16.156344| -11.627675| -28.565975| 1.8531e+03| |
| 324 | +| 16| 5| -16.204699| -11.633976| -28.581277| 1.8531e+03| |
| 325 | ++----------------------------------------------------------------+ |
| 326 | +| Multizone Summary | |
| 327 | ++----------------------------------------------------------------+ |
| 328 | +| Outer_Iter| avg[bgs][0]| avg[bgs][1]|MinVolume[0]|DeformIter[0| |
| 329 | ++----------------------------------------------------------------+ |
| 330 | +| 16| -11.228939| -9.130581| 8.7338e-10| 44| |
| 331 | +``` |
| 332 | + |
| 333 | +The displacement field on the structural domain and the velocity field on the flow domain obtained in ```fsi_steady_1.vtu```_and ```fsi_steady_0.vtu``` respectively are shown below: |
| 334 | + |
| 335 | + |
| 336 | + |
| 337 | + |
| 338 | + |
| 339 | +which is, as expected, identical to the result of the [Static FSI](../Static_FSI/) tutorial. The convergence of the drag coefficient is as follows: |
| 340 | + |
| 341 | + |
| 342 | + |
| 343 | +The drag coefficient will be used as the objective function for the adjoint computation in the [next tutorial](../Adjoint_FSI_Python/). |
| 344 | + |
| 345 | +### Attribution |
| 346 | + |
| 347 | +If you are using this content for your research, please kindly cite the following reference in your derived works: |
| 348 | + |
| 349 | +Sanchez, R. _et al._ (2016), ASSESSMENT OF THE FLUID-STRUCTURE INTERACTION CAPABILITIES FOR AERONAUTICAL APPLICATIONS OF THE OPEN-SOURCE SOLVER SU2, _Proceedings of the VII European Congress on Computational Methods in Applied Sciences and Engineering_. DOI: [10.7712/100016.1903.6597](https://doi.org/10.7712/100016.1903.6597) |
| 350 | + |
| 351 | +<dl> |
| 352 | +This work is licensed under a <a rel="license" href="http://creativecommons.org/licenses/by/4.0/">Creative Commons Attribution 4.0 International License</a> |
| 353 | +<br /> |
| 354 | +<a rel="license" href="http://creativecommons.org/licenses/by/4.0/"><img alt="Creative Commons Licence" style="border-width:0" src="https://i.creativecommons.org/l/by/4.0/88x31.png" /></a> |
| 355 | +</dl> |
| 356 | + |
0 commit comments