Skip to content
New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

Updated options documentation #33

Merged
merged 8 commits into from
Nov 12, 2020
Merged
Show file tree
Hide file tree
Changes from all commits
Commits
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
31 changes: 31 additions & 0 deletions doc/BC.rst
Original file line number Diff line number Diff line change
@@ -0,0 +1,31 @@
.. _pyhyp_BC:

Using the BC option
===================

This page describes how to use the ``BC`` option to specify boundary conditions at boundary edges of the surface mesh.

Here is an example of a dictionary that can be used with ``BC``::

"BC":{1:{"iLow":"ySymm"}, 3:{"jHigh":"splay", "iHigh":"xConst"}}

Each entry in the dictionary has a key and at least one nested key-value pair.
The key is the 1-based block number, the nested key is the boundary edge specification, and the value is the boundary condition.
For the first entry, these are ``1``, ``iLow``, and ``ySymm``, respectively.

The 1-based block number and boundary edge specification for a boundary edge can be determined using Tecplot:

#. Load the surface mesh file into Tecplot.
#. Use the Probe tool to select a point on the boundary edge of interest.
Use Ctrl+click to snap on to a boundary edge.
#. Select Zone/Cell Info in the toolbar on the right.
#. The number shown after Zone is the 1-based block number.
#. The edge specification depends on the values of I and J.
The edge is ``iLow`` if I=1, ``iHigh`` if I=I-Max, ``jLow`` if J=1, or ``jHigh`` if J=J-Max.
Only one of these is true for any boundary edge.

The supported boundary conditions are:

* ``splay`` for free edges
* ``xSymm``, ``ySymm``, ``zSymm`` for symmetry planes
* ``xConst``, ``yConst``, ``zConst``, ``xyConst``, ``yzConst``, ``xzConst`` for constant planes
1 change: 1 addition & 0 deletions doc/index.rst
Original file line number Diff line number Diff line change
Expand Up @@ -32,5 +32,6 @@ Contents:
install
plot3d
cgns
BC
options
API
250 changes: 159 additions & 91 deletions doc/options.rst
Original file line number Diff line number Diff line change
Expand Up @@ -5,94 +5,162 @@ Options

Here are the options currently available in pyHyp.

=================== ========== =======================================================================================
Parameter Type Description
=================== ========== =======================================================================================
``inputFile`` ``char`` Name of the file that contains the surface mesh.
This is a file that has been generated in an external meshing program, typically ICEMCFD.

``fileType`` ``char`` Type of the input file. Use either ``Plot3d`` or ``CGNS``.

``N`` ``int`` Number of grid levels to march.
This determines the grid dimension in the off-wall direction.
Typically this should be a "multi-grid" friendly number.

``s0`` ``float`` Initial off-wall (normal) spacing of grid.
This is taken to be constant across the entire geometry.
The units are consistent with the rest of the geometry.

``rMin`` ``float`` Relative distance in the normal direction to march.
It is specified as a multiple of the radius of the sphere enclosing the initial surface geometry.
If symmetry is specified, the full mirrored geometry is used to compute the sphere's radius.
Most wing geometries will have ``rMin`` between 10 and 20, that is the farfield boundary is 20 spans away from the geometry.

``cMax`` ``float`` The maximum permissible ratio of marching direction length to the any other in-plane edge.
This parameter effectively operates as a CFL-type limit.
If a step would require a step which would result in a ratio ``c`` greater than ``cMax``, the step is automatically split internally to respect this user-supplied limit.
Typical values of ``cMax`` are around 6-8.
Increased robustness can be achieved at the expense of computational cost by lowering ``cMax``.

``nonLinear`` ``float`` Use the nonlinear formulation.
This is experimental and not currently recommended and may not work at all.

``slExp`` ``float`` Exponent for the :math:`S_l` computation.
The :math:`S_l` value serves the same purpose as found in Chan et al. but the computation is different.
The :math:`S_l` computation in Chan is given as :math:`\sqrt{\frac{N-1}{l-1}}` for :math:`l > 2`.

``ps0`` ``float`` Initial pseudo offwall spacing.
This spacing **must** be less than or equal to ``s0``.
This is actual spacing the hyperbolic scheme uses.
The solver may take many pseudo steps before the first real grid level at ``s0``.

``pGridRatio`` ``float`` The ratio between successive levels in the pseudo grid.
This will be typically somewhere between ~1.05 for large grids to 1.2 for small grids.
This number is **not** the actual grid spacing of the final grid; that spacing ratio is computed and displayed at the beginning of a calculation.
The ``pGridRatio`` **must** be smaller than that number.

``epsE`` ``float`` The explict smoothing parameter.
See the :ref:`Theory<pyhyp_theory>` section for more information.
Typical values are approximately 1.0. Increasing the explicit smoothing may result in a smoother grid, at the expense of orhtogonality.
If the geometry is very sharp corners, too much explicit smoothing will cause the solver to rapidly "soften" the corner and the grid will fold back on itself.
In concave corners, additional smoothing will prevent lines from crossing (avoiding negative cells).

``epsI`` ``float`` Implicit smoothing parameter.
See the :ref:`Theory<pyhyp_theory>` section for more information.
Typical values are from 2.0 to 6.0.
Generally increasing the implicit coefficient results in a more stable solution procedure.
Usually this value should be twice the explicit smoothing parameter.

``theta`` ``float`` Kinsley-Barth coefficient See the :ref:`Theory<pyhyp_theory>` section for more information.
Only a single theta value is used for both directions.
Typical values are ~2.0 to ~4.0.

``volCoef`` ``float`` Coefficient used in point-Jacobi local volume smoothing algorithm.
Typically this value is 0.16 and need not be modified.
Use more ``volSmoothIter`` for stronger local smoothing.

``volBlend`` ``float`` The global volume blending coefficient.
See the :ref:`Theory<pyhyp_theory>` section for more information.
This value will typically be very small, especially if you widely varying cell sizes.
Typically values are from ~0 to 0.001.
Default is 0.0001.

``volSmoothIter`` ``int`` The number of point-Jacobi local volume smoothing iterations to perform.
Typical values are ~5 to ~25.
Default is 10.

``kspRelTol`` ``float`` Tolerance for the solution of the linear system at each iteration.
Typically :math:`1\times 10^{-8}` is sufficient.
Very difficult cases may benefit from a tighter convergence tolerance.

``kspMaxIts`` ``int`` Maximum number of iterations to perform for each step.
Default is 500 which should be sufficient for most cases.

``preConLag`` ``int`` Lag the update of the preconditioner by this number of iterations.
The default value of 10 will typically not need to be changed.

``kspSubspaceSize`` ``int`` Size of the ksp subspace.
Default is 50.
Very large and difficult problems may befefit from a larger subspace size.

``writeMetrics`` ``bool`` Flag to write the mesh gradients to the solution file.
This option should only be used for debugging purposes.
=================== ========== =======================================================================================
.. list-table::
:widths: 5 5 90
:header-rows: 1

* - Parameter
- Type
- Description

* - ``inputFile``
- ``str``
- Name of the file that contains the surface mesh.
This is a file that has been generated in an external meshing program, typically ICEMCFD.

* - ``fileType``
- ``str``
- Type of the input file.
Use either ``Plot3d`` or ``CGNS``.

* - ``unattachedEdgesAreSymmetry``
- ``bool``
- Automatically applies symmetry boundary conditions to any edges that do not interface with another block.
This option works in many cases but does not work for all surface meshes.
If you encounter negative volumes near the symmetry plane, try explicitly setting the symmetry boundary conditions using the ``BC`` option.

* - ``outerFaceBC``
- ``str``
- Specifies the boundary condition at the outermost face of the extruded mesh.
Use either ``farfield`` or ``overset``.

* - ``BC``
- ``dict``
- Specifies boundary condition information for specific block edges. See :ref:`here<pyhyp_BC>` for details.

* - ``families``
- ``str`` / ``dict``
- Name given to wall surfaces.
If a dictionary is submitted, each wall patch can be named separately.
This can help with applying operations to specific wall patches.

* - ``N``
- ``int``
- Number of grid levels to march.
This determines the grid dimension in the off-wall direction.
Typically this should be a "multi-grid" friendly number.

* - ``s0``
- ``float``
- Initial off-wall (normal) spacing of grid.
This is taken to be constant across the entire geometry.
The units are consistent with the rest of the geometry.

* - ``rMin``
- ``float``
- Relative distance in the normal direction to march.
It is specified as a multiple of the radius of the sphere enclosing the initial surface geometry.
If symmetry is specified, the full mirrored geometry is used to compute the sphere's radius.
Most wing geometries will have ``rMin`` between 10 and 20, that is the farfield boundary is 20 spans away from the geometry.

* - ``cMax``
- ``float``
- The maximum permissible ratio of marching direction length to the any other in-plane edge.
This parameter effectively operates as a CFL-type limit.
If a step would require a step which would result in a ratio ``c`` greater than ``cMax``, the step is automatically split internally to respect this user-supplied limit.
Typical values of ``cMax`` are around 6-8.
Increased robustness can be achieved at the expense of computational cost by lowering ``cMax``.

* - ``nonLinear``
- ``float``
- Use the nonlinear formulation.
This is experimental and not currently recommended and may not work at all.

* - ``slExp``
- ``float``
- Exponent for the :math:`S_l` computation.
The :math:`S_l` value serves the same purpose as found in Chan et al. but the computation is different.
The :math:`S_l` computation in Chan is given as :math:`\sqrt{\frac{N-1}{l-1}}` for :math:`l > 2`.

* - ``ps0``
- ``float``
- Initial pseudo offwall spacing.
This spacing **must** be less than or equal to ``s0``.
This is actual spacing the hyperbolic scheme uses.
The solver may take many pseudo steps before the first real grid level at ``s0``.

* - ``pGridRatio``
- ``float``
- The ratio between successive levels in the pseudo grid.
This will be typically somewhere between ~1.05 for large grids to 1.2 for small grids.
This number is **not** the actual grid spacing of the final grid; that spacing ratio is computed and displayed at the beginning of a calculation.
The ``pGridRatio`` **must** be smaller than that number.

* - ``epsE``
- ``float``
- The explict smoothing parameter.
See the :ref:`Theory<pyhyp_theory>` section for more information.
Typical values are approximately 1.0. Increasing the explicit smoothing may result in a smoother grid, at the expense of orhtogonality.
If the geometry is very sharp corners, too much explicit smoothing will cause the solver to rapidly "soften" the corner and the grid will fold back on itself.
In concave corners, additional smoothing will prevent lines from crossing (avoiding negative cells).

* - ``epsI``
- ``float``
- Implicit smoothing parameter.
See the :ref:`Theory<pyhyp_theory>` section for more information.
Typical values are from 2.0 to 6.0.
Generally increasing the implicit coefficient results in a more stable solution procedure.
Usually this value should be twice the explicit smoothing parameter.

* - ``theta``
- ``float``
- Kinsley-Barth coefficient See the :ref:`Theory<pyhyp_theory>` section for more information.
Only a single theta value is used for both directions.
Typical values are ~2.0 to ~4.0.

* - ``volCoef``
- ``float``
- Coefficient used in point-Jacobi local volume smoothing algorithm.
Typically this value is 0.16 and need not be modified.
Use more ``volSmoothIter`` for stronger local smoothing.

* - ``volBlend``
- ``float``
- The global volume blending coefficient.
See the :ref:`Theory<pyhyp_theory>` section for more information.
This value will typically be very small, especially if you widely varying cell sizes.
Typically values are from ~0 to 0.001.
Default is 0.0001.

* - ``volSmoothIter``
- ``int``
- The number of point-Jacobi local volume smoothing iterations to perform.
Typical values are ~5 to ~25.
Default is 10.

* - ``kspRelTol``
- ``float``
- Tolerance for the solution of the linear system at each iteration.
Typically :math:`1\times 10^{-8}` is sufficient.
Very difficult cases may benefit from a tighter convergence tolerance.

* - ``kspMaxIts``
- ``int``
- Maximum number of iterations to perform for each step.
Default is 500 which should be sufficient for most cases.

* - ``preConLag``
- ``int``
- Lag the update of the preconditioner by this number of iterations.
The default value of 10 will typically not need to be changed.

* - ``kspSubspaceSize``
- ``int``
- Size of the ksp subspace.
Default is 50.
Very large and difficult problems may befefit from a larger subspace size.

* - ``writeMetrics``
- ``bool``
- Flag to write the mesh gradients to the solution file.
This option should only be used for debugging purposes.
9 changes: 5 additions & 4 deletions pyhyp/pyHyp.py
Original file line number Diff line number Diff line change
Expand Up @@ -619,14 +619,15 @@ def __init__(self, comm=None, options=None, debug=False):
)
BCToSet = BCs[blkBC][edgeKey].lower()
if BCToSet.lower() not in BCMap.keys():

raise Error(
"Boundary condition specification unknown. Must be one of: "
"'splay', 'BCXSymm', 'BCYSymm', 'BCZSymm', "
"'BCXConst', 'BCYConst', 'BCZConst, 'BCXYConst, "
"'BCYZConst or BCXZConst'. %s" % helpStr
"'splay', 'xSymm', 'ySymm', 'zSymm', "
"'xConst', 'yConst', 'zConst, 'xyConst, "
"'yzConst or xzConst'. %s" % helpStr
)

fBCs[edgeMap[lKey], blkBC - 1] = BCMap[BCToSet]
fBCs[edgeMap[lKey], blkBC - 1] = BCMap[BCToSet]

# Set the boundary condition information into fortran
self.hyp.hypinput.bcs = fBCs
Expand Down