Implementation of the Moment-of-Fluid (MOF) method. More...
Namespaces | |
module | fields_mof |
Module containing the fields related to the Moment-of-Fluid method. | |
module | variables_mof |
Module containing the parameters of the Moment-of-Fluid method. | |
module | variables_vof_plic |
Module containing the parameters of the Moment-of-Fluid method. | |
Classes | |
type | type_mof_phase_extension::t_mof_filament |
Filament structure for MOF extension. More... | |
type | type_mof_phase_extension::t_mof_phase_extension |
Definition of the MOF extension of the t_phase_geometry type. More... | |
Functions | |
subroutine, public | mod_mof2d_analytical_reconstruction::mof2d_symmetric_analytical_reconstruction (polygon, cell_centroid, ref_centroid1, ref_centroid2, ref_volume1, ref_volume2, normal, residual1, residual2, polygon1, polygon2) |
Analytic MOF symmetric reconstruction for two materials in a rectangular cell. | |
subroutine, public | mod_mof2d_analytical_reconstruction::mof2d_analytical_reconstruction (polygon, ref_centroid, ref_volume, normal, residual, polygon_full, polygon_empty) |
Analytic MOF reconstruction for two materials in a rectangular cell. | |
subroutine, public | mod_mof2d_backward_advection::mof2d_backward_advection (mof_phases, boundary_condition, velocity1, velocity2, time_step) |
Advection of the MOF phases using the backward advection algorithm. | |
double precision elemental function | mod_mof2d_coordinate_system_wrapper::cg2w_polygon_volume (polygon) |
Wrapper function for coordinate system abstraction. | |
pure subroutine | mod_mof2d_coordinate_system_wrapper::cg2w_polygon_centroid (polygon, centroid, volume) |
Wrapper function for coordinate system abstraction. | |
pure subroutine | mod_mof2d_coordinate_system_wrapper::cg2w_flood_polygon (polygon, normal, volume, s0, s1, polygon_full, polygon_empty) |
Wrapper function for coordinate system abstraction. | |
subroutine | mod_mof2d_global_reconstruction::mof2d_global_reconstruction (phases) |
Reconstruct the interface on every cells in 2D. | |
pure subroutine, public | mod_mof2d_minimization_reconstruction::mof2d_symmetric_minimization_reconstruction (polygon, ref_centroid1, ref_centroid2, ref_volume1, ref_volume2, angle, normal, residual1, residual2, polygon1, polygon2) |
Find the MOF linear reconstruction in a convex cell using a minimization algorithm. | |
pure subroutine, public | mod_mof2d_minimization_reconstruction::mof2d_minimization_reconstruction (polygon, ref_centroid, ref_volume, angle, normal, residual, polygon_full, polygon_empty) |
Find the MOF linear reconstruction in a convex cell using a minimization algorithm. | |
subroutine, public | mod_mof_mpi_exchange_filaments::mof_mpi_exchange_filaments (mof_phases) |
Exchange the filament structure between neighbor processors. | |
subroutine, public | mod_mof2d_one_cell_reconstruction::mof2d_one_cell_reconstruction (cell, ref_volume, ref_centroid, polygon_list, best_cost) |
Reconstruct the optimal location of the materials in one cell. | |
pure subroutine, public | mod_mof3d_analytic_centroid::mof3d_compute_analytic_gradient (s_angles, ref_centroid1, ref_centroid2, ref_volume, c, objective, gradient) |
Compute the centroid and the gradient of the objective function in rectangular hexahedral cell. | |
pure subroutine, public | mod_mof3d_analytic_centroid::mof3d_transform_angles_to_reference (transformation, orig_angles, ref_angles) |
Transform the spherical angles in the original configuration to the reference configuration. | |
pure subroutine, public | mod_mof3d_analytic_centroid::mof3d_transform_angles_to_original (transformation, ref_angles, orig_angles) |
Transform the spherical angles in the reference configuration to the original configuration. | |
pure subroutine, public | mod_mof3d_analytic_centroid::mof3d_compute_residual_analytic (s_angles, ref_centroid1, ref_centroid2, ref_volume, c, residual, jacobian) |
Compute the residual and its jacobian in a rectangular hexahedral cell. | |
pure subroutine, public | mod_mof3d_bfgs::mof3d_bfgs (cell, cell_volume, cell_centroid, ref_centroid1, ref_centroid2, ref_volume, use_analytic, angles, normal, stat, residual) |
Broyden-Fletcher-Goldfarb-Shanno method to minimize the objective function of MOF. | |
pure character(len=:) function, allocatable, public | mod_mof3d_bfgs::mof3d_get_stop_criterion (stat) |
Convert the stop criterion code to a human readable string. | |
subroutine | mod_mof3d_global_reconstruction::mof3d_global_reconstruction (mof_phases) |
Reconstruct the interface on every cells in 3D. | |
pure subroutine, public | mod_mof3d_gradient::mof3d_compute_gradient (polyhedron, cell_volume, cell_centroid, angles, ref_volume, ref_centroid1, ref_centroid2, use_analytic, objective, gradient) |
Compute the gradient of the objective function of MOF. | |
pure subroutine, public | mod_mof3d_gradient::mof3d_compute_residual (polyhedron, cell_volume, cell_centroid, angles, ref_volume, ref_centroid1, ref_centroid2, use_analytic, residual, jacobian) |
Compute the residual of MOF and its jacobian. | |
pure subroutine | mod_mof3d_initial_angles::mof3d_tetra_compute_initial_angles (p0, p1, p2, p3, ref_centroid1, ref_centroid2, ref_volume, best_angles) |
Compute the initial angles for the minimization algorithm for tetrahedrons. | |
pure subroutine | mod_mof3d_initial_angles::mof3d_hexa_compute_initial_angles (c, ref_centroid1, ref_centroid2, ref_volume, best_angles) |
Compute the initial angles for the minimization for hexahedrons. | |
pure subroutine, public | mod_mof3d_one_cell_reconstruction::mof3d_one_cell_reconstruction (cell, ref_volume, ref_centroid, polyhedron_list, best_cost) |
Reconstruct the optimal location of the materials in one cell. | |
pure subroutine, public | mod_mof3d_tetra_analytic_centroid::mof3d_tetra_compute_analytic_gradient (p0, p1, p2, p3, angles, ref_centroid1, ref_centroid2, ref_volume, objective, gradient) |
Compute the centroid and the gradient of the objective function in a tetrahedral cell. | |
pure subroutine, public | mod_mof3d_tetra_analytic_centroid::mof3d_tetra_compute_residual (p0, p1, p2, p3, angles, ref_centroid1, ref_centroid2, ref_volume, residual, jacobian) |
Compute the centroid and the derivatives in a tetrahedral cell. | |
pure subroutine, public | mod_mof3d_tetra_analytic_centroid::mof3d_tetra_compute_transformation (p0, p1, p2, p3, transformation) |
Compute the matrix to transform the reference tetrahedron to the original (transformed) tetrahedron. | |
pure subroutine, public | mod_mof3d_tetra_analytic_centroid::mof3d_tetra_transform_angles_to_reference (transformation, orig_angles, ref_angles) |
Transform the spherical angles in the original configuration to the reference configuration. | |
pure subroutine, public | mod_mof3d_tetra_analytic_centroid::mof3d_tetra_transform_angles_to_original (transformation, ref_angles, orig_angles) |
Transform the spherical angles in the reference configuration to the original configuration. | |
pure subroutine, public | mod_mof3d_tetra_analytic_centroid::mof3d_tetra_compute_analytic_gradient_reference (angles, ref_centroid, volume, objective, gradient) |
Compute the centroid and the gradient of the objective function in the reference tetrahedron. | |
subroutine | mod_mof_advect_phases::mof_advect_phases () |
Interface to advect phases using the Moment-of-Fluid representation. | |
subroutine | mod_mof_apply_boundary_conditions::mof_apply_boundary_conditions (mof_phases, mof_bc) |
Apply boundary conditions the MOF phases. | |
subroutine | mod_mof_compute_cfl::mof_compute_cfl (nb_time_step, time_step) |
Compute the time step and the number of time sub-step according the MOF CFL. | |
Moment-of-Fluid (MOF) is a semi-Lagrangian method to advect interfaces in multiphase flow simulation. This method is the apex of the VOF methods using piecewise linear interface reconstruction. Namely, MOF can manage the reconstruction of a large number of interfaces with automatic ordering using only a one-cell stencil (compared to the 9 (27 in 3D) point stencil of the PLIC method).
MOF is a sharp interface method to represent n immiscible fluids on a mesh; n can be greater than 2. Each cell contains the volume fraction and the centroid of all the fluids.
Consider a cell \(\Omega\) and a partition \(\{\omega_1,\cdots,\omega_n\}\) of this cell. The volume of a given domain \(\omega\) is denoted \(|\omega|=M_0(\omega)\). It is given by
\[ M_0(\omega) = \int_\omega dx \]
The first momentum of a given domain \(\omega\) is denoted \(M_1(\omega)\). It is given by
\[ M_1(\omega) = \int_\omega x dx \]
The centroid of a given domain \(\omega\) is denoted \(x_c(\omega)\). It is given by
\[ x_c(\omega) = M_1(\omega)/M_0(\omega) \]
From the volume fractions and the centroids we are able to reconstruct a piecewise linear interface between the n fluids.
Consider a cell \(\Omega\) containing two sub-domains \(\omega_1\) and \(\omega_2\). The reconstruction consists in finding two sub-domains \(\omega_1^\ell\) and \(\omega_2^\ell\) such as:
This reconstruction is unique (see Dyadechko & Shashkov 2005) and can be generalized to n fluids (see Dyadechko & Shashkov 2008).
This implementation of MOF uses a backward advection scheme to advect the fluids as proposed by Dyadechko & Shashkov 2005. The description of the algorithm is available in the subroutine mof2d_backward_advection.
pure subroutine mod_mof2d_coordinate_system_wrapper::cg2w_flood_polygon | ( | type(t_polygon), intent(in) | polygon, |
double precision, dimension(2), intent(in) | normal, | ||
double precision, intent(in) | volume, | ||
double precision, dimension(2), intent(out) | s0, | ||
double precision, dimension(2), intent(out) | s1, | ||
type(t_polygon), intent(out) | polygon_full, | ||
type(t_polygon), intent(out) | polygon_empty ) |
[in] | polygon | Any convex polygon whose points are defined in counterclockwise order |
[in] | normal | The direction of flood given by a unit vector |
[in] | volume | The required volume fraction of the flooded area |
[out] | s0,s1 | The extreme points of the interface line segment |
[out] | polygon_full | Part of the polygon full of fluid |
[out] | polygon_empty | Part of the polygon without fluid |
pure subroutine mod_mof2d_coordinate_system_wrapper::cg2w_polygon_centroid | ( | type(t_polygon), intent(in) | polygon, |
double precision, dimension(2), intent(out) | centroid, | ||
double precision, intent(out), optional | volume ) |
[in] | polygon | Any polygon |
[out] | centroid | Centroid of the polygon |
[out] | volume | Volume of the polygon |
double precision elemental function mod_mof2d_coordinate_system_wrapper::cg2w_polygon_volume | ( | type(t_polygon), intent(in) | polygon | ) |
[in] | polygon | Any polygon |
subroutine, public mod_mof2d_analytical_reconstruction::mof2d_analytical_reconstruction | ( | type(t_polygon), intent(in) | polygon, |
double precision, dimension(2), intent(in) | ref_centroid, | ||
double precision, intent(in) | ref_volume, | ||
double precision, dimension(2), intent(out) | normal, | ||
double precision, intent(out) | residual, | ||
type(t_polygon), intent(out) | polygon_full, | ||
type(t_polygon), intent(out) | polygon_empty ) |
The volume fraction is expected to be in range [0,1/2]
The polygon
must be defined in this order:
!! p4 p3 !! ┌──────────┐ !! │ │ !! │ │c2 ▲ y !! │ c1 │ ┃ !! └──────────┘ ┗━━▶ x !! p1 p2 O !!
subroutine, public mod_mof2d_backward_advection::mof2d_backward_advection | ( | type(t_phase_geometry), dimension(:), intent(inout) | mof_phases, |
type(t_boundary_condition), dimension(:), intent(in) | boundary_condition, | ||
type(t_face_field), intent(in) | velocity1, | ||
type(t_face_field), intent(in) | velocity2, | ||
double precision | time_step ) |
This is a second order advection scheme. It requires the velocity at the beginning of the time step and at the end of the time step.
Algorithm description:
!! ┌──────────┐ ┌──────────┐ ┌──────────┐ ┌──────────┐ !! │ │ │ │ │ │ │ │ !! │ Phase 1 │ │ Phase 2 │ │ Phase 3 │ … │ Phase n │ !! │ │ │ │ │ │ │ │ !! └──────────┘ └──────────┘ └──────────┘ └──────────┘ !! ╰───────────────┴──────┬──────┴───────────────╯ !! │ !! ╓─────────────────────╖ !! ╭────────╢ Begin loop on cells ║ !! │ ╙─────────────────────╜ !! │ │ !! │ ┌────────────────────────────────┐ !! │ │ Backward advection of the cell │ !! │ └────────────────────────────────┘ !! │ │ !! │ ┌────────────────────────────────────┐ !! │ │ Intersect the surrounding polygons │ !! │ └────────────────────────────────────┘ !! │ │ !! │ ┌────────────────────────────────────┐ !! │ │ Forward advection of the centroids │ !! │ └────────────────────────────────────┘ !! │ │ !! │ ┌───────────────────┐ !! │ │ Volume correction │ !! │ └───────────────────┘ !! │ │ !! │ ╓───────────────────╖ !! ╰─────────╢ End loop on cells ║ !! ╙───────────────────╜ !! │ !! ┌────────────────┐ !! │ Reconstruction │ !! └────────────────┘ !!
[in,out] | mof_phases | list of MOF phases. |
[in] | boundary_condition | boundary conditions. |
[in] | velocity1 | velocity at the beginning of the time step. |
[in] | velocity2 | velocity at the end of the time step. |
[in] | time_step | time step. |
subroutine mod_mof2d_global_reconstruction::mof2d_global_reconstruction | ( | type(t_phase_geometry), dimension(:), intent(inout) | phases | ) |
Note: the 'is_mixed' flag must be set before reconstruction.
pure subroutine, public mod_mof2d_minimization_reconstruction::mof2d_minimization_reconstruction | ( | type(t_polygon), intent(in) | polygon, |
double precision, dimension(2), intent(in) | ref_centroid, | ||
double precision, intent(in) | ref_volume, | ||
double precision, intent(inout) | angle, | ||
double precision, dimension(2), intent(out) | normal, | ||
double precision, intent(out) | residual, | ||
type(t_polygon), intent(out) | polygon_full, | ||
type(t_polygon), intent(out) | polygon_empty ) |
[in] | polygon | convex polygon to split into full and empty part |
[in] | ref_centroid | reference centroid |
[in] | ref_volume | reference volume |
[in,out] | angle | angle first guess. Returns the final angle |
[out] | normal | return the interface normal |
[out] | residual | return the value of the objective function |
[out] | polygon_full | contains the full part of the initial polygon |
[out] | polygon_empty | contains the empty part of the initial polygon |
subroutine, public mod_mof2d_one_cell_reconstruction::mof2d_one_cell_reconstruction | ( | type(t_polygon), intent(in) | cell, |
double precision, dimension(:), intent(in) | ref_volume, | ||
double precision, dimension(:,:), intent(in) | ref_centroid, | ||
type(t_polygon), dimension(:), intent(out) | polygon_list, | ||
double precision, intent(out) | best_cost ) |
cell
must be defined the same way: !! 4┌──────┐3 !! │ │ !! 1└──────┘2 !!
[in] | cell | Convex polyhedral cell. |
[in] | ref_volume | Array of reference volumes. |
[in] | ref_centroid | Array of reference centroids. |
[out] | polygon_list | Array that contains the reconstructed polygons. |
[out] | best_cost | Value of the optimal cost. |
subroutine, public mod_mof2d_analytical_reconstruction::mof2d_symmetric_analytical_reconstruction | ( | type(t_polygon), intent(in) | polygon, |
double precision, dimension(2), intent(in) | cell_centroid, | ||
double precision, dimension(2), intent(in) | ref_centroid1, | ||
double precision, dimension(2), intent(in) | ref_centroid2, | ||
double precision, intent(in) | ref_volume1, | ||
double precision, intent(in) | ref_volume2, | ||
double precision, dimension(2), intent(out) | normal, | ||
double precision, intent(out) | residual1, | ||
double precision, intent(out) | residual2, | ||
type(t_polygon), intent(out) | polygon1, | ||
type(t_polygon), intent(out) | polygon2 ) |
The volume fraction of the first material is expected to be in range [0,1/2].
The polygon
must be defined in this order:
!! p4 p3 !! ┌──────────┐ !! │ │ !! │ │c2 ▲ y !! │ c1 │ ┃ !! └──────────┘ ┗━━▶ x !! p1 p2 O !!
pure subroutine, public mod_mof2d_minimization_reconstruction::mof2d_symmetric_minimization_reconstruction | ( | type(t_polygon), intent(in) | polygon, |
double precision, dimension(2), intent(in) | ref_centroid1, | ||
double precision, dimension(2), intent(in) | ref_centroid2, | ||
double precision, intent(in) | ref_volume1, | ||
double precision, intent(in) | ref_volume2, | ||
double precision, intent(inout) | angle, | ||
double precision, dimension(2), intent(out) | normal, | ||
double precision, intent(out) | residual1, | ||
double precision, intent(out) | residual2, | ||
type(t_polygon), intent(out) | polygon1, | ||
type(t_polygon), intent(out) | polygon2 ) |
[in] | polygon | convex polygon to split into full and empty part |
[in] | ref_centroid1 | reference centroid of material 1 |
[in] | ref_centroid2 | reference centroid of material 2 |
[in] | ref_volume1 | reference volume of material 1 |
[in] | ref_volume2 | reference volume of material 2 |
[in,out] | angle | angle first guess. Returns the final angle |
[out] | normal | return the interface normal |
[out] | residual1 | return the squared distance of the reference centroid of material 1 to its reference centroid |
[out] | residual2 | return the squared distance of the reference centroid of material 2 to its reference centroid |
[out] | polygon1 | contains the polygon of material 1 |
[out] | polygon2 | contains the polygon of material 2 |
pure subroutine, public mod_mof3d_bfgs::mof3d_bfgs | ( | type(t_polyhedron), intent(in) | cell, |
double precision, intent(in) | cell_volume, | ||
double precision, dimension(3), intent(in) | cell_centroid, | ||
double precision, dimension(3), intent(in) | ref_centroid1, | ||
double precision, dimension(3), intent(in) | ref_centroid2, | ||
double precision, intent(in) | ref_volume, | ||
logical, intent(in) | use_analytic, | ||
double precision, dimension(2), intent(inout) | angles, | ||
double precision, dimension(3), intent(out) | normal, | ||
integer, dimension(3), intent(out) | stat, | ||
double precision, dimension(2), intent(out) | residual ) |
Reference:
[in] | cell | Convex polyhedron. |
[in] | cell_volume | Cell volume. |
[in] | cell_centroid | Cell centroid. |
[in] | ref_centroid1 | Coordinates of the reference centroid of the material 1. |
[in] | ref_centroid2 | Coordinates of the reference centroid of the material 2. |
[in] | ref_volume | Reference volume of the material 1. |
[in] | use_analytic | Use analytic formulas. |
[in,out] | angles | Initial guess and final angles. |
[out] | normal | Unit normal of the reconstructed interface (material 1 → material 2). |
[out] | stat | Gradient calls. stat(1): in BFGS. stat(2): in line search. stat(3): exit status. |
[out] | residual | Residuals. residual(1): derivative, residual(2): objective function. |
pure subroutine, public mod_mof3d_analytic_centroid::mof3d_compute_analytic_gradient | ( | double precision, dimension(2), intent(in) | s_angles, |
double precision, dimension(3), intent(in) | ref_centroid1, | ||
double precision, dimension(3), intent(in) | ref_centroid2, | ||
double precision, intent(in) | ref_volume, | ||
double precision, dimension(3), intent(in) | c, | ||
double precision, intent(out) | objective, | ||
double precision, dimension(2), intent(out) | gradient ) |
[in] | s_angles | Spherical angles in the unit cube. |
[in] | ref_centroid1 | Coordinates of the reference centroid of the material 1. |
[in] | ref_centroid2 | Coordinates of the reference centroid of the material 2. |
[in] | ref_volume | Reference volume. |
[in] | c | Dimensions of the cell. |
[out] | objective | Objective function. |
[out] | gradient | Gradient of the objective function. |
pure subroutine, public mod_mof3d_gradient::mof3d_compute_gradient | ( | type(t_polyhedron), intent(in) | polyhedron, |
double precision, intent(in) | cell_volume, | ||
double precision, dimension(3), intent(in) | cell_centroid, | ||
double precision, dimension(2), intent(in) | angles, | ||
double precision, intent(in) | ref_volume, | ||
double precision, dimension(3), intent(in) | ref_centroid1, | ||
double precision, dimension(3), intent(in) | ref_centroid2, | ||
logical, intent(in) | use_analytic, | ||
double precision, intent(out) | objective, | ||
double precision, dimension(2), intent(out) | gradient ) |
[in] | polyhedron | Convex polyhedron. |
[in] | cell_volume | Volume of the cell. |
[in] | cell_centroid | Coordinates of the centroid of the cell. |
[in] | angles | Spherical angles. |
[in] | ref_volume | Reference volume of the material 1. |
[in] | ref_centroid1 | Coordinates of the reference centroid of the material 1. |
[in] | ref_centroid2 | Coordinates of the reference centroid of the material 2. |
[in] | use_analytic | Use analytic gradient. |
[out] | objective | Objective function. |
[out] | gradient | Coordinates of the gradient. |
pure subroutine, public mod_mof3d_gradient::mof3d_compute_residual | ( | type(t_polyhedron), intent(in) | polyhedron, |
double precision, intent(in) | cell_volume, | ||
double precision, dimension(3), intent(in) | cell_centroid, | ||
double precision, dimension(2), intent(in) | angles, | ||
double precision, intent(in) | ref_volume, | ||
double precision, dimension(3), intent(in) | ref_centroid1, | ||
double precision, dimension(3), intent(in) | ref_centroid2, | ||
logical, intent(in) | use_analytic, | ||
double precision, dimension(:), intent(out) | residual, | ||
double precision, dimension(:,:), intent(out) | jacobian ) |
[in] | polyhedron | Convex polyhedron. |
[in] | cell_volume | Volume of the cell. |
[in] | cell_centroid | Coordinates of the centroid of the cell. |
[in] | angles | Spherical angles. |
[in] | ref_volume | Reference volume of the material 1. |
[in] | ref_centroid1 | Coordinates of the reference centroid of the material 1. |
[in] | ref_centroid2 | Coordinates of the reference centroid of the material 2. |
[in] | use_analytic | Use analytic gradient. |
[out] | residual | Residual function. |
[out] | jacobian | Jacobian of the residual. |
pure subroutine, public mod_mof3d_analytic_centroid::mof3d_compute_residual_analytic | ( | double precision, dimension(2), intent(in) | s_angles, |
double precision, dimension(3), intent(in) | ref_centroid1, | ||
double precision, dimension(3), intent(in) | ref_centroid2, | ||
double precision, intent(in) | ref_volume, | ||
double precision, dimension(3), intent(in) | c, | ||
double precision, dimension(:), intent(out) | residual, | ||
double precision, dimension(:,:), intent(out) | jacobian ) |
[in] | s_angles | Spherical angles in the unit cube. |
[in] | ref_centroid1 | Coordinates of the reference centroid of the material 1. |
[in] | ref_centroid2 | Coordinates of the reference centroid of the material 2. |
[in] | ref_volume | Reference volume of the material 1. |
[in] | c | Dimensions of the cell. |
[out] | residual | Residual function. |
[out] | jacobian | Jacobian of the residual function. |
pure character(len=:) function, allocatable, public mod_mof3d_bfgs::mof3d_get_stop_criterion | ( | integer, dimension(3), intent(in) | stat | ) |
[in] | stat | Statistic array from the BFGS. |
pure subroutine mod_mof3d_initial_angles::mof3d_hexa_compute_initial_angles | ( | double precision, dimension(3), intent(in) | c, |
double precision, dimension(3), intent(in) | ref_centroid1, | ||
double precision, dimension(3), intent(in) | ref_centroid2, | ||
double precision, intent(in) | ref_volume, | ||
double precision, dimension(2), intent(out) | best_angles ) |
[in] | c | Dimensions of the hexahedron. |
[in] | ref_centroid1 | Coordinates of the first reference centroid. |
[in] | ref_centroid2 | Coordinates of the second reference centroid. |
[in] | ref_volume | Reference volume. |
[out] | best_angles | Initial angles. |
pure subroutine, public mod_mof3d_one_cell_reconstruction::mof3d_one_cell_reconstruction | ( | type(t_polyhedron), intent(in) | cell, |
double precision, dimension(:), intent(in) | ref_volume, | ||
double precision, dimension(:,:), intent(in) | ref_centroid, | ||
type(t_polyhedron), dimension(:), intent(out) | polyhedron_list, | ||
double precision, intent(out) | best_cost ) |
[in] | cell | Convex polyhedral cell. |
[in] | ref_volume | Array of reference volumes. |
[in] | ref_centroid | Array of reference centroids. |
[out] | polyhedron_list | Array that contains the reconstructed polyedrons. |
[out] | best_cost | Value of the optimal cost. |
pure subroutine, public mod_mof3d_tetra_analytic_centroid::mof3d_tetra_compute_analytic_gradient | ( | double precision, dimension(3), intent(in) | p0, |
double precision, dimension(3), intent(in) | p1, | ||
double precision, dimension(3), intent(in) | p2, | ||
double precision, dimension(3), intent(in) | p3, | ||
double precision, dimension(2), intent(in) | angles, | ||
double precision, dimension(3), intent(in) | ref_centroid1, | ||
double precision, dimension(3), intent(in) | ref_centroid2, | ||
double precision, intent(in) | ref_volume, | ||
double precision, intent(out) | objective, | ||
double precision, dimension(2), intent(out) | gradient ) |
[in] | p0,p1,p2,p3 | Vertices of the tetrahedron. |
[in] | angles | Spherical angles. |
[in] | ref_centroid1 | Coordinates of the reference centroid of the material 1. |
[in] | ref_centroid2 | Coordinates of the reference centroid of the material 2. |
[in] | ref_volume | Reference volume. |
[out] | objective | Objective function. |
[out] | gradient | Gradient of the objective function. |
pure subroutine, public mod_mof3d_tetra_analytic_centroid::mof3d_tetra_compute_analytic_gradient_reference | ( | double precision, dimension(2), intent(in) | angles, |
double precision, dimension(3), intent(in) | ref_centroid, | ||
double precision, intent(in) | volume, | ||
double precision, intent(out) | objective, | ||
double precision, dimension(2), intent(out) | gradient ) |
[in] | angles | Spherical angles. |
[in] | ref_centroid | Coordinates of the reference centroid. |
[in] | volume | Reference volume. |
[out] | objective | Objective function. |
[out] | gradient | Gradient of the objective function. |
pure subroutine mod_mof3d_initial_angles::mof3d_tetra_compute_initial_angles | ( | double precision, dimension(3), intent(in) | p0, |
double precision, dimension(3), intent(in) | p1, | ||
double precision, dimension(3), intent(in) | p2, | ||
double precision, dimension(3), intent(in) | p3, | ||
double precision, dimension(3), intent(in) | ref_centroid1, | ||
double precision, dimension(3), intent(in) | ref_centroid2, | ||
double precision, intent(in) | ref_volume, | ||
double precision, dimension(2), intent(out) | best_angles ) |
[in] | p0,p1,p2,p3 | Vertices of the tetrahedron. |
[in] | ref_centroid1 | Coordinates of the first reference centroid. |
[in] | ref_centroid2 | Coordinates of the second reference centroid. |
[in] | ref_volume | Reference volume. |
[out] | best_angles | Initial angles. |
pure subroutine, public mod_mof3d_tetra_analytic_centroid::mof3d_tetra_compute_residual | ( | double precision, dimension(3), intent(in) | p0, |
double precision, dimension(3), intent(in) | p1, | ||
double precision, dimension(3), intent(in) | p2, | ||
double precision, dimension(3), intent(in) | p3, | ||
double precision, dimension(2), intent(in) | angles, | ||
double precision, dimension(3), intent(in) | ref_centroid1, | ||
double precision, dimension(3), intent(in) | ref_centroid2, | ||
double precision, intent(in) | ref_volume, | ||
double precision, dimension(:), intent(out) | residual, | ||
double precision, dimension(:,:), intent(out) | jacobian ) |
[in] | p0,p1,p2,p3 | vertices of the tetrahedron. |
[in] | angles | Spherical angles. |
[in] | ref_centroid1 | Coordinates of the reference centroid of the material 1. |
[in] | ref_centroid2 | Coordinates of the reference centroid of the material 2. |
[in] | ref_volume | Reference volume. |
[out] | residual | Residual function. |
[out] | jacobian | Partial derivatives of the residual. |
pure subroutine, public mod_mof3d_tetra_analytic_centroid::mof3d_tetra_compute_transformation | ( | double precision, dimension(3), intent(in) | p0, |
double precision, dimension(3), intent(in) | p1, | ||
double precision, dimension(3), intent(in) | p2, | ||
double precision, dimension(3), intent(in) | p3, | ||
double precision, dimension(3,3), intent(out) | transformation ) |
The transformation is given by X(ξ) = Aξ + b. This routine returns the matrix A providing the vertices of the original (transformed) tetrahedron.
[in] | p0,p1,p2,p3 | Vertices of the tetrahedron. |
[out] | transformation | Tranformation matrix. |
pure subroutine, public mod_mof3d_tetra_analytic_centroid::mof3d_tetra_transform_angles_to_original | ( | double precision, dimension(3,3), intent(in) | transformation, |
double precision, dimension(2), intent(in) | ref_angles, | ||
double precision, dimension(2), intent(out) | orig_angles ) |
[in] | transformation | Tranformation matrix. |
[in] | ref_angles | Angles in the reference configuration. |
[out] | orig_angles | Angles in the original configuration. |
pure subroutine, public mod_mof3d_tetra_analytic_centroid::mof3d_tetra_transform_angles_to_reference | ( | double precision, dimension(3,3), intent(in) | transformation, |
double precision, dimension(2), intent(in) | orig_angles, | ||
double precision, dimension(2), intent(out) | ref_angles ) |
[in] | transformation | Tranformation matrix. |
[in] | orig_angles | Angles in the original configuration. |
[out] | ref_angles | Angles in the reference configuration. |
pure subroutine, public mod_mof3d_analytic_centroid::mof3d_transform_angles_to_original | ( | double precision, dimension(3), intent(in) | transformation, |
double precision, dimension(2), intent(in) | ref_angles, | ||
double precision, dimension(2), intent(out) | orig_angles ) |
[in] | transformation | Tranformation matrix. |
[in] | ref_angles | Angles in the reference configuration. |
[out] | orig_angles | Angles in the original configuration. |
pure subroutine, public mod_mof3d_analytic_centroid::mof3d_transform_angles_to_reference | ( | double precision, dimension(3), intent(in) | transformation, |
double precision, dimension(2), intent(in) | orig_angles, | ||
double precision, dimension(2), intent(out) | ref_angles ) |
[in] | transformation | Tranformation matrix. |
[in] | orig_angles | Angles in the original configuration. |
[out] | ref_angles | Angles in the reference configuration. |
subroutine mod_mof_compute_cfl::mof_compute_cfl | ( | integer, intent(out) | nb_time_step, |
double precision, intent(out) | time_step ) |
[out] | nb_time_step | Number of time sub-steps |
[out] | time_step | Time step |
subroutine, public mod_mof_mpi_exchange_filaments::mof_mpi_exchange_filaments | ( | type(t_phase_geometry), dimension(:), intent(inout) | mof_phases | ) |
[in,out] | mof_phases | Array of phases |