.. DO NOT EDIT.
.. THIS FILE WAS AUTOMATICALLY GENERATED BY SPHINX-GALLERY.
.. TO MAKE CHANGES, EDIT THE SOURCE PYTHON FILE:
.. "examples/00-fluent/modeling_ablation.py"
.. LINE NUMBERS ARE GIVEN BELOW.

.. only:: html

    .. note::
        :class: sphx-glr-download-link-note

        :ref:`Go to the end <sphx_glr_download_examples_00-fluent_modeling_ablation.py>`
        to download the full example code.

.. rst-class:: sphx-glr-example-title

.. _sphx_glr_examples_00-fluent_modeling_ablation.py:

.. _modeling_ablation:

Modeling Ablation
-------------------------------------------

.. GENERATED FROM PYTHON SOURCE LINES 30-59

Objective
=====================================================================================

Ablation is an effective treatment used to protect an atmospheric reentry vehicle from
the damaging effects of external high temperatures caused by shock wave and viscous
heating. The ablative material is chipped away due to surface reactions that remove a
significant amount of heat and keep the vehicle surface temperature below the melting
point. In this tutorial, Fluent ablation model is demonstrated for a reendtry vehicle
geometry simplified as a 3D wedge.

This tutorial demonstrates how to do the following:

* Define boundary conditions for a high-speed flow.
* Set up the ablation model to model effects of a moving boundary due to ablation.
* Initiate and solve the transient simulation using the density-based solver.

Problem Description:
====================

The geometry of the 3D wedge considered in this tutorial is shown in following figure.
The air flow passes around a nose of a re-entry vehicle operating under high speed
conditions. The inlet air has a temperature of 4500 K, a gauge pressure of 13500 Pa,
and a Mach number of 3. The domain is bounded above and below by symmetry planes
(displayed in yellow). As the ablative coating chips away, the surface of the wall
moves. The moving of the surface is modeled using dynamic meshes. The surface moving
rate is estimated by Vieille's empirical law:

where r is the surface moving rate, p is the absolute pressure, and A and n are model
parameters. In the considered case, A = 5 and n = 0.1.

.. GENERATED FROM PYTHON SOURCE LINES 62-65

.. math::

   r = A \cdot p^n

.. GENERATED FROM PYTHON SOURCE LINES 68-71

.. image:: ../../_static/ablation-problem-schematic.png
   :align: center
   :alt: Problem Schematic

.. GENERATED FROM PYTHON SOURCE LINES 75-77

Import required libraries/modules
==================================================================================

.. GENERATED FROM PYTHON SOURCE LINES 77-82

.. code-block:: Python


    import ansys.fluent.core as pyfluent
    from ansys.fluent.core import examples
    from ansys.fluent.visualization.pyvista import Graphics


.. GENERATED FROM PYTHON SOURCE LINES 83-85

Download example file
==================================================================================

.. GENERATED FROM PYTHON SOURCE LINES 85-90

.. code-block:: Python

    import_filename = examples.download_file(
        "ablation.msh.h5",
        "pyfluent/examples/Ablation-tutorial",
    )


.. GENERATED FROM PYTHON SOURCE LINES 91-93

Fluent Solution Setup
==================================================================================

.. GENERATED FROM PYTHON SOURCE LINES 93-98

.. code-block:: Python


    from ansys.fluent.visualization import set_config  # noqa: E402

    set_config(blocking=True, set_view_on_display="isometric")


.. GENERATED FROM PYTHON SOURCE LINES 99-101

Launch Fluent session with solver mode and print Fluent version
==================================================================================

.. GENERATED FROM PYTHON SOURCE LINES 101-110

.. code-block:: Python


    solver = pyfluent.launch_fluent(
        product_version="25.1.0",
        dimension=3,
        precision="double",
        processor_count=4,
    )
    print(solver.get_fluent_version())


.. GENERATED FROM PYTHON SOURCE LINES 111-113

Import mesh
==================================================================================

.. GENERATED FROM PYTHON SOURCE LINES 113-116

.. code-block:: Python


    solver.file.read_case(file_name=import_filename)


.. GENERATED FROM PYTHON SOURCE LINES 117-119

Define models
==================================================================================

.. GENERATED FROM PYTHON SOURCE LINES 119-126

.. code-block:: Python


    solver.setup.general.solver.type = "density-based-implicit"
    solver.setup.general.solver.time = "unsteady-1st-order"
    solver.setup.general.operating_conditions.operating_pressure = 0.0
    solver.setup.models.energy = {"enabled": True}
    solver.setup.models.ablation = {"enabled": True}


.. GENERATED FROM PYTHON SOURCE LINES 127-129

Define material
=================================================================

.. GENERATED FROM PYTHON SOURCE LINES 129-136

.. code-block:: Python


    solver.setup.materials.fluid["air"] = {
        "density": {"option": "ideal-gas"},
        "molecular_weight": {"value": 28.966, "option": "constant"},
    }
    solver.setup.materials.fluid["air"] = {"density": {"option": "ideal-gas"}}


.. GENERATED FROM PYTHON SOURCE LINES 137-139

Define boundary conditions
==========================

.. GENERATED FROM PYTHON SOURCE LINES 139-159

.. code-block:: Python


    solver.setup.boundary_conditions.set_zone_type(
        zone_list=["inlet"], new_type="pressure-far-field"
    )
    solver.setup.boundary_conditions.pressure_far_field["inlet"].momentum.gauge_pressure = (
        13500
    )
    solver.setup.boundary_conditions.pressure_far_field["inlet"].momentum.mach_number = 3
    solver.setup.boundary_conditions.pressure_far_field["inlet"].thermal.temperature = 4500
    solver.setup.boundary_conditions.pressure_far_field[
        "inlet"
    ].turbulence.turbulent_intensity = 0.001

    solver.setup.boundary_conditions.pressure_outlet["outlet"].momentum.gauge_pressure = (
        13500
    )
    solver.setup.boundary_conditions.pressure_outlet[
        "outlet"
    ].momentum.prevent_reverse_flow = True


.. GENERATED FROM PYTHON SOURCE LINES 160-166

Ablation boundary condition (Vielles Model)
++++++++++++++++++++++++++++++++++++++++++++
Once you have specified the ablation boundary conditions for the wall,
Ansys Fluent automatically enables the Dynamic Mesh model with the Smoothing and
Remeshing options, #creates the wall-ablation dynamic mesh zone, and configure
appropriate dynamic mesh settings.

.. GENERATED FROM PYTHON SOURCE LINES 166-173

.. code-block:: Python


    solver.setup.boundary_conditions.wall[
        "wall_ablation"
    ].ablation.ablation_select_model = "Vielle's Model"
    solver.setup.boundary_conditions.wall["wall_ablation"].ablation.ablation_vielle_a = 5
    solver.setup.boundary_conditions.wall["wall_ablation"].ablation.ablation_vielle_n = 0.1


.. GENERATED FROM PYTHON SOURCE LINES 174-176

Define dynamic mesh controls
============================

.. GENERATED FROM PYTHON SOURCE LINES 176-258

.. code-block:: Python


    solver.tui.define.dynamic_mesh.dynamic_mesh("yes")
    solver.tui.define.dynamic_mesh.zones.create(
        "interior--flow",
        "deforming",
        "faceted",
        "no",
        "no",
        "yes",
        "no",
        "yes",
        "yes",
        "no",
        "yes",
    )
    solver.tui.define.dynamic_mesh.zones.create(
        "outlet",
        "deforming",
        "faceted",
        "no",
        "yes",
        "no",
        "yes",
        "yes",
        "coefficient-based",
        "0.1",
        "yes",
    )
    solver.tui.define.dynamic_mesh.zones.create(
        "symm1",
        "deforming",
        "plane",
        "0",
        "-0.04",
        "0",
        "0",
        "-1",
        "0",
        "no",
        "yes",
        "no",
        "yes",
        "yes",
        "coefficient-based",
        "0.1",
        "yes",
    )
    solver.tui.define.dynamic_mesh.zones.create(
        "symm2",
        "deforming",
        "plane",
        "0",
        "0.04",
        "0",
        "0",
        "1",
        "0",
        "no",
        "yes",
        "no",
        "yes",
        "yes",
        "coefficient-based",
        "0.1",
        "yes",
    )
    solver.tui.define.dynamic_mesh.zones.create(
        "wall_ablation",
        "user-defined",
        "**ablation**",
        "no",
        "no",
        "189",
        "constant",
        "0",
        "yes",
        "yes",
        "0.7",
        "no",
        "no",
    )


.. GENERATED FROM PYTHON SOURCE LINES 259-261

Define solver settings
=======================

.. GENERATED FROM PYTHON SOURCE LINES 261-266

.. code-block:: Python


    solver.setup.general.solver.time = "unsteady-2nd-order"
    solver.solution.controls.limits = {"max_temperature": 25000}
    solver.solution.monitor.residual.equations["energy"].absolute_criteria = 1e-06


.. GENERATED FROM PYTHON SOURCE LINES 267-269

Create report definitions
==========================

.. GENERATED FROM PYTHON SOURCE LINES 269-328

.. code-block:: Python


    solver.solution.report_definitions.drag["drag_force_x"] = {}
    solver.solution.report_definitions.drag["drag_force_x"].zones = ["wall_ablation"]

    solver.solution.monitor.report_plots.create(name="drag_force_x")
    solver.solution.monitor.report_plots["drag_force_x"].report_defs = "drag_force_x"
    solver.tui.solve.report_plots.axes(
        "drag_force_x", "numbers", "float", "4", "exponential", "2", "q"
    )

    solver.solution.monitor.report_files.create(name="drag_force_x")
    solver.solution.monitor.report_files["drag_force_x"] = {
        "report_defs": ["drag_force_x"],
        "file_name": r"drag_force_x.out",
    }

    solver.solution.report_definitions.surface["pressure_avg_abl_wall"] = {}
    solver.solution.report_definitions.surface["pressure_avg_abl_wall"].report_type = (
        "surface-areaavg"
    )
    solver.solution.report_definitions.surface["pressure_avg_abl_wall"].field = "pressure"
    solver.solution.report_definitions.surface["pressure_avg_abl_wall"].surface_names = [
        "wall_ablation"
    ]

    solver.solution.monitor.report_plots.create(name="pressure_avg_abl_wall")
    solver.solution.monitor.report_plots["pressure_avg_abl_wall"].report_defs = (
        "pressure_avg_abl_wall"
    )
    solver.tui.solve.report_plots.axes(
        "pressure_avg_abl_wall", "numbers", "float", "4", "exponential", "2", "q"
    )

    solver.solution.monitor.report_files.create(name="pressure_avg_abl_wall")
    solver.solution.monitor.report_files["pressure_avg_abl_wall"] = {
        "report_defs": ["pressure_avg_abl_wall"],
        "file_name": r"pressure_avg_abl_wall.out",
    }

    solver.solution.report_definitions.surface["recede_point"] = {}
    solver.solution.report_definitions.surface["recede_point"].report_type = (
        "surface-vertexmax"
    )
    solver.solution.report_definitions.surface["recede_point"].field = "z-coordinate"
    solver.solution.report_definitions.surface["recede_point"].surface_names = [
        "wall_ablation"
    ]

    solver.solution.monitor.report_plots.create(name="recede_point")
    solver.solution.monitor.report_plots["recede_point"].report_defs = "recede_point"
    solver.tui.solve.report_plots.axes(
        "recede_point", "numbers", "float", "4", "exponential", "2", "q"
    )
    solver.solution.monitor.report_files.create(name="recede_point")
    solver.solution.monitor.report_files["recede_point"] = {
        "report_defs": ["recede_point"],
        "file_name": r"recede_point.out",
    }


.. GENERATED FROM PYTHON SOURCE LINES 329-331

Initialize and Save case
========================

.. GENERATED FROM PYTHON SOURCE LINES 331-339

.. code-block:: Python


    solver.tui.solve.initialize.compute_defaults.pressure_far_field("inlet")
    solver.solution.initialization.initialization_type = "standard"
    solver.solution.initialization.standard_initialize()
    solver.solution.run_calculation.transient_controls.time_step_size = 1e-06

    solver.file.write(file_type="case", file_name="ablation.cas.h5")


.. GENERATED FROM PYTHON SOURCE LINES 340-343

Run the calculation
===================
Note: It may take about half an hour to finish the calculation.

.. GENERATED FROM PYTHON SOURCE LINES 343-348

.. code-block:: Python


    solver.solution.run_calculation.dual_time_iterate(
        time_step_count=100, max_iter_per_step=20
    )


.. GENERATED FROM PYTHON SOURCE LINES 349-352

Save simulation data
====================
Write case and data files

.. GENERATED FROM PYTHON SOURCE LINES 352-354

.. code-block:: Python

    solver.file.write(file_type="case-data", file_name="ablation_Solved.cas.h5")


.. GENERATED FROM PYTHON SOURCE LINES 355-357

Post Processing
==================================================================================

.. GENERATED FROM PYTHON SOURCE LINES 359-361

Display plots
=============

.. GENERATED FROM PYTHON SOURCE LINES 363-366

.. image:: ../../_static/ablation-residual.png
   :align: center
   :alt: residual

.. GENERATED FROM PYTHON SOURCE LINES 368-369

Scaled residual plot

.. GENERATED FROM PYTHON SOURCE LINES 371-374

.. image:: ../../_static/ablation-drag_force_x.png
   :align: center
   :alt: Drag force in x direction

.. GENERATED FROM PYTHON SOURCE LINES 376-377

History of the drag force on the ablation wall

.. GENERATED FROM PYTHON SOURCE LINES 379-382

.. image:: ../../_static/ablation-avg_pressure.png
   :align: center
   :alt: Average pressure on the ablation wall

.. GENERATED FROM PYTHON SOURCE LINES 384-385

History of the averaged pressure on the ablation wall

.. GENERATED FROM PYTHON SOURCE LINES 387-390

.. image:: ../../_static/ablation-recede_point.png
   :align: center
   :alt: Recede point (albation)

.. GENERATED FROM PYTHON SOURCE LINES 392-393

Recede point (deformation due to ablation)

.. GENERATED FROM PYTHON SOURCE LINES 395-398

Display contour
================
Following contours are displayed in the Fluent GUI:

.. GENERATED FROM PYTHON SOURCE LINES 398-416

.. code-block:: Python


    solver.results.surfaces.plane_surface.create(name="mid_plane")
    solver.results.surfaces.plane_surface["mid_plane"].method = "zx-plane"

    solver.results.graphics.contour.create(name="contour_pressure")
    solver.results.graphics.contour["contour_pressure"] = {
        "field": "pressure",
        "surfaces_list": ["mid_plane"],
    }
    solver.results.graphics.contour.display(object_name="contour_pressure")

    solver.results.graphics.contour.create(name="contour_mach")
    solver.results.graphics.contour["contour_mach"] = {
        "field": "mach-number",
        "surfaces_list": ["mid_plane"],
    }
    solver.results.graphics.contour.display(object_name="contour_mach")


.. GENERATED FROM PYTHON SOURCE LINES 417-420

Post processing with PyVista (3D visualization)
===============================================
Following graphics is displayed in the a new window/notebook

.. GENERATED FROM PYTHON SOURCE LINES 420-426

.. code-block:: Python


    graphics_session1 = Graphics(solver)
    contour1 = graphics_session1.Contours["contour-1"]
    contour1.field = "pressure"
    contour1.surfaces_list = ["mid_plane"]
    contour1.display()

.. GENERATED FROM PYTHON SOURCE LINES 427-430

.. image:: ../../_static/ablation-pressure.png
   :align: center
   :alt: Static Pressure Contour

.. GENERATED FROM PYTHON SOURCE LINES 432-433

Static Pressure Contour

.. GENERATED FROM PYTHON SOURCE LINES 433-440

.. code-block:: Python


    contour1.field = "mach-number"
    contour1.range.option = "auto-range-off"
    contour1.range.auto_range_off.minimum = 0.5
    contour1.range.auto_range_off.maximum = 3.0
    contour1.display()


.. GENERATED FROM PYTHON SOURCE LINES 441-444

.. image:: ../../_static/ablation-mach-number.png
   :align: center
   :alt: Mach Number Contour

.. GENERATED FROM PYTHON SOURCE LINES 446-447

Mach Number Contour

.. GENERATED FROM PYTHON SOURCE LINES 449-451

Close the session
==================================================================================

.. GENERATED FROM PYTHON SOURCE LINES 451-454

.. code-block:: Python


    solver.exit()



.. _sphx_glr_download_examples_00-fluent_modeling_ablation.py:

.. only:: html

  .. container:: sphx-glr-footer sphx-glr-footer-example

    .. container:: sphx-glr-download sphx-glr-download-jupyter

      :download:`Download Jupyter notebook: modeling_ablation.ipynb <modeling_ablation.ipynb>`

    .. container:: sphx-glr-download sphx-glr-download-python

      :download:`Download Python source code: modeling_ablation.py <modeling_ablation.py>`

    .. container:: sphx-glr-download sphx-glr-download-zip

      :download:`Download zipped: modeling_ablation.zip <modeling_ablation.zip>`


.. only:: html

 .. rst-class:: sphx-glr-signature

    `Gallery generated by Sphinx-Gallery <https://sphinx-gallery.github.io>`_