version 0.6.0

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. More...
 
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. More...
 
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. More...
 
double precision elemental function mod_mof2d_coordinate_system_wrapper::cg2w_polygon_volume (polygon)
 Wrapper function for coordinate system abstraction. More...
 
pure subroutine mod_mof2d_coordinate_system_wrapper::cg2w_polygon_centroid (polygon, centroid, volume)
 Wrapper function for coordinate system abstraction. More...
 
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. More...
 
subroutine mod_mof2d_global_reconstruction::mof2d_global_reconstruction (phases)
 Reconstruct the interface on every cells in 2D. More...
 
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. More...
 
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. More...
 
subroutine, public mod_mof_mpi_exchange_filaments::mof_mpi_exchange_filaments (mof_phases)
 Exchange the filament structure between neighbor processors. More...
 
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. More...
 
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. More...
 
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. More...
 
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. More...
 
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. More...
 
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. More...
 
pure character(len=:) function, allocatable, public mod_mof3d_bfgs::mof3d_get_stop_criterion (stat)
 Convert the stop criterion code to a human readable string. More...
 
pure subroutine mod_mof3d_gauss_newton::mof3d_gauss_newton2 (cell, cell_volume, cell_centroid, ref_centroid1, ref_centroid2, ref_volume, use_analytic, angles, normal, stat, residual)
 Gauss-Νewton method to minimize the objective function of MOF. More...
 
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. More...
 
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. More...
 
pure subroutine mod_mof3d_gradient::compute_residual_geometric (polyhedron, cell_volume, cell_centroid, angles, ref_volume, ref_centroid1, ref_centroid2, residual, jacobian)
 Compute the residual and the jacobian using geometric approaches. More...
 
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. More...
 
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. More...
 
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. More...
 
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. More...
 
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. More...
 
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. More...
 
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. More...
 
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. More...
 
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. More...
 
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. More...
 

Detailed Description

Warning
This method is currently implemented in 2D but not in 3D.

Introduction

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).

Description of the representation

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) \]

Reconstruction

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).

Advection

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.

References

Function Documentation

◆ cg2w_flood_polygon()

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 
)
Parameters
[in]polygonAny convex polygon whose points are defined in counterclockwise order
[in]normalThe direction of flood given by a unit vector
[in]volumeThe required volume fraction of the flooded area
[out]s0,s1The extreme points of the interface line segment
[out]polygon_fullPart of the polygon full of fluid
[out]polygon_emptyPart of the polygon without fluid

◆ cg2w_polygon_centroid()

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 
)
Parameters
[in]polygonAny polygon
[out]centroidCentroid of the polygon
[out]volumeVolume of the polygon

◆ cg2w_polygon_volume()

double precision elemental function mod_mof2d_coordinate_system_wrapper::cg2w_polygon_volume ( type(t_polygon), intent(in)  polygon)
Parameters
[in]polygonAny polygon

◆ compute_residual_geometric()

pure subroutine mod_mof3d_gradient::compute_residual_geometric ( 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,
double precision, dimension(:), intent(out)  residual,
double precision, dimension(:,:), intent(out)  jacobian 
)
Parameters
[in]polyhedronConvex polyhedron.
[in]cell_volumeVolume of the cell.
[in]cell_centroidCoordinates of the centroid of the cell.
[in]anglesSpherical angles.
[in]ref_volumeReference volume of the material 1.
[in]ref_centroid1Coordinates of the reference centroid of the material 1.
[in]ref_centroid2Coordinates of the reference centroid of the material 2.
[out]residualResidual function.

◆ mof2d_analytical_reconstruction()

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

◆ mof2d_backward_advection()

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 │
                     └────────────────┘
Parameters
[in,out]mof_phaseslist of MOF phases.
[in]boundary_conditionboundary conditions.
[in]velocity1velocity at the beginning of the time step.
[in]velocity2velocity at the end of the time step.
[in]time_steptime step.

◆ mof2d_global_reconstruction()

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.

◆ mof2d_minimization_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 
)
Parameters
[in]polygonconvex polygon to split into full and empty part
[in]ref_centroidreference centroid
[in]ref_volumereference volume
[in,out]angleangle first guess. Returns the final angle
[out]normalreturn the interface normal
[out]residualreturn the value of the objective function
[out]polygon_fullcontains the full part of the initial polygon
[out]polygon_emptycontains the empty part of the initial polygon

◆ mof2d_one_cell_reconstruction()

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 
)
Warning
This routine depends on mod_mof2d_analytical_reconstruction::mof2d_analytical_reconstruction, cell must be defined the same way:
   4┌──────┐3
    │      │
   1└──────┘2
Parameters
[in]cellConvex polyhedral cell.
[in]ref_volumeArray of reference volumes.
[in]ref_centroidArray of reference centroids.
[out]polygon_listArray that contains the reconstructed polygons.
[out]best_costValue of the optimal cost.

◆ mof2d_symmetric_analytical_reconstruction()

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

◆ mof2d_symmetric_minimization_reconstruction()

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 
)
Parameters
[in]polygonconvex polygon to split into full and empty part
[in]ref_centroid1reference centroid of material 1
[in]ref_centroid2reference centroid of material 2
[in]ref_volume1reference volume of material 1
[in]ref_volume2reference volume of material 2
[in,out]angleangle first guess. Returns the final angle
[out]normalreturn the interface normal
[out]residual1return the squared distance of the reference centroid of material 1 to its reference centroid
[out]residual2return the squared distance of the reference centroid of material 2 to its reference centroid
[out]polygon1contains the polygon of material 1
[out]polygon2contains the polygon of material 2

◆ mof3d_bfgs()

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:

  • R. Fletcher, Practical methods of optimization, John Wiley & Sons, 1987.
Parameters
[in]cellConvex polyhedron.
[in]cell_volumeCell volume.
[in]cell_centroidCell centroid.
[in]ref_centroid1Coordinates of the reference centroid of the material 1.
[in]ref_centroid2Coordinates of the reference centroid of the material 2.
[in]ref_volumeReference volume of the material 1.
[in]use_analyticUse analytic formulas.
[in,out]anglesInitial guess and final angles.
[out]normalUnit normal of the reconstructed interface (material 1 → material 2).
[out]statGradient calls. stat(1): in BFGS. stat(2): in line search. stat(3): exit status.
[out]residualResiduals. residual(1): derivative, residual(2): objective function.

◆ mof3d_compute_analytic_gradient()

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 
)
Parameters
[in]s_anglesSpherical angles in the unit cube.
[in]ref_centroid1Coordinates of the reference centroid of the material 1.
[in]ref_centroid2Coordinates of the reference centroid of the material 2.
[in]ref_volumeReference volume.
[in]cDimensions of the cell.
[out]objectiveObjective function.
[out]gradientGradient of the objective function.

◆ mof3d_compute_gradient()

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 
)
Parameters
[in]polyhedronConvex polyhedron.
[in]cell_volumeVolume of the cell.
[in]cell_centroidCoordinates of the centroid of the cell.
[in]anglesSpherical angles.
[in]ref_volumeReference volume of the material 1.
[in]ref_centroid1Coordinates of the reference centroid of the material 1.
[in]ref_centroid2Coordinates of the reference centroid of the material 2.
[in]use_analyticUse analytic gradient.
[out]objectiveObjective function.
[out]gradientCoordinates of the gradient.

◆ mof3d_compute_residual()

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 
)
Parameters
[in]polyhedronConvex polyhedron.
[in]cell_volumeVolume of the cell.
[in]cell_centroidCoordinates of the centroid of the cell.
[in]anglesSpherical angles.
[in]ref_volumeReference volume of the material 1.
[in]ref_centroid1Coordinates of the reference centroid of the material 1.
[in]ref_centroid2Coordinates of the reference centroid of the material 2.
[in]use_analyticUse analytic gradient.
[out]residualResidual function.
[out]jacobianJacobian of the residual.

◆ mof3d_compute_residual_analytic()

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 
)
Parameters
[in]s_anglesSpherical angles in the unit cube.
[in]ref_centroid1Coordinates of the reference centroid of the material 1.
[in]ref_centroid2Coordinates of the reference centroid of the material 2.
[in]ref_volumeReference volume of the material 1.
[in]cDimensions of the cell.
[out]residualResidual function.
[out]jacobianJacobian of the residual function.

◆ mof3d_gauss_newton2()

pure subroutine mod_mof3d_gauss_newton::mof3d_gauss_newton2 ( 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, intent(out)  stat,
double precision, dimension(2), intent(out)  residual 
)

Reference:

  • R. Fletcher, Practical methods of optimization, John Wiley & Sons, 1987.
Parameters
[in]cellConvex polyhedron.
[in]cell_volumeCell volume.
[in]cell_centroidCell centroid.
[in]ref_centroid1Coordinates of the reference centroid of the material 1.
[in]ref_centroid2Coordinates of the reference centroid of the material 2.
[in]ref_volumeReference volume of the material 1.
[in]use_analyticUse analytic formulas.
[in,out]anglesInitial guess and final angles.
[out]normalUnit normal of the reconstructed interface (material 1 → material 2).
[out]statGradient calls.
[out]residualResiduals. residual(1): derivative, residual(2): objective function.

◆ mof3d_get_stop_criterion()

pure character(len=:) function, allocatable, public mod_mof3d_bfgs::mof3d_get_stop_criterion ( integer, dimension(3), intent(in)  stat)
Parameters
[in]statStatistic array from the BFGS.
Returns
criterion: String containing the description of the stop criterion.

◆ mof3d_hexa_compute_initial_angles()

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 
)
Parameters
[in]cDimensions of the hexahedron.
[in]ref_centroid1Coordinates of the first reference centroid.
[in]ref_centroid2Coordinates of the second reference centroid.
[in]ref_volumeReference volume.
[out]best_anglesInitial angles.

◆ mof3d_one_cell_reconstruction()

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 
)
Parameters
[in]cellConvex polyhedral cell.
[in]ref_volumeArray of reference volumes.
[in]ref_centroidArray of reference centroids.
[out]polyhedron_listArray that contains the reconstructed polyedrons.
[out]best_costValue of the optimal cost.

◆ mof3d_tetra_compute_analytic_gradient()

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 
)
Parameters
[in]p0,p1,p2,p3Vertices of the tetrahedron.
[in]anglesSpherical angles.
[in]ref_centroid1Coordinates of the reference centroid of the material 1.
[in]ref_centroid2Coordinates of the reference centroid of the material 2.
[in]ref_volumeReference volume.
[out]objectiveObjective function.
[out]gradientGradient of the objective function.

◆ mof3d_tetra_compute_analytic_gradient_reference()

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 
)
Parameters
[in]anglesSpherical angles.
[in]ref_centroidCoordinates of the reference centroid.
[in]ref_volumeReference volume.
[out]objectiveObjective function.
[out]gradientGradient of the objective function.

◆ mof3d_tetra_compute_initial_angles()

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 
)
Parameters
[in]p0,p1,p2,p3Vertices of the tetrahedron.
[in]ref_centroid1Coordinates of the first reference centroid.
[in]ref_centroid2Coordinates of the second reference centroid.
[in]ref_volumeReference volume.
[out]best_anglesInitial angles.

◆ mof3d_tetra_compute_residual()

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 
)
Parameters
[in]p0,p1,p2,p3vertices of the tetrahedron.
[in]anglesSpherical angles.
[in]ref_centroid1Coordinates of the reference centroid of the material 1.
[in]ref_centroid2Coordinates of the reference centroid of the material 2.
[in]ref_volumeReference volume.
[out]residualResidual function.
[out]jacobianPartial derivatives of the residual.

◆ mof3d_tetra_compute_transformation()

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.

Parameters
[in]p0,p1,p2,p3Vertices of the tetrahedron.
[out]transformationTranformation matrix.

◆ mof3d_tetra_transform_angles_to_original()

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 
)
Parameters
[in]transformationTranformation matrix.
[in]ref_anglesAngles in the reference configuration.
[out]orig_anglesAngles in the original configuration.

◆ mof3d_tetra_transform_angles_to_reference()

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 
)
Parameters
[in]transformationTranformation matrix.
[in]orig_anglesAngles in the original configuration.
[out]ref_anglesAngles in the reference configuration.

◆ mof3d_transform_angles_to_original()

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 
)
Parameters
[in]transformationTranformation matrix.
[in]ref_anglesAngles in the reference configuration.
[out]orig_anglesAngles in the original configuration.

◆ mof3d_transform_angles_to_reference()

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 
)
Parameters
[in]transformationTranformation matrix.
[in]orig_anglesAngles in the original configuration.
[out]ref_anglesAngles in the reference configuration.

◆ mof_compute_cfl()

subroutine mod_mof_compute_cfl::mof_compute_cfl ( integer, intent(out)  nb_time_step,
double precision, intent(out)  time_step 
)
Parameters
[out]nb_time_stepNumber of time sub-steps
[out]time_stepTime step

◆ mof_mpi_exchange_filaments()

subroutine, public mod_mof_mpi_exchange_filaments::mof_mpi_exchange_filaments ( type(t_phase_geometry), dimension(:), intent(inout)  mof_phases)
Parameters
[in,out]mof_phasesArray of phases