film Module

FE solution of the Reynolds equation

Description of the film module

This is the main module for the finite element solution of the Reynolds equation.

Definition of FE_film

FE_FILM is a data structure containing a FE_MESH and some data describing a lubrication problem

Solution procedure

This module can be used to create a FE_FILM, assemble the FE_FIM and solve it. Some post-treatements like fluxes, forces are available.


Uses

  • module~~film~~UsesGraph module~film film module~surfile surfile module~film->module~surfile module~data_arch data_arch module~film->module~data_arch module~solver solver module~film->module~solver module~num_param num_param module~film->module~num_param module~elements elements module~film->module~elements module~fluid_law fluid_law module~film->module~fluid_law module~data_film_hd data_film_hd module~film->module~data_film_hd module~mesh mesh module~film->module~mesh omp_lib omp_lib module~film->omp_lib module~surfile->module~data_arch module~sort_arrays sort_arrays module~surfile->module~sort_arrays iso_c_binding iso_c_binding module~surfile->iso_c_binding iso_fortran_env iso_fortran_env module~data_arch->iso_fortran_env module~solver->module~data_arch hsl_ma48_double hsl_ma48_double module~solver->hsl_ma48_double module~solver->iso_fortran_env module~solver->module~sort_arrays module~mumps_wrapper mumps_wrapper module~solver->module~mumps_wrapper module~gen_param gen_param module~solver->module~gen_param module~mumfpack mumfpack module~solver->module~mumfpack module~solver->iso_c_binding module~sulu_wrapper sulu_wrapper module~solver->module~sulu_wrapper module~num_param->module~data_arch module~num_param->iso_fortran_env module~elements->module~data_arch module~fluid_law->module~data_arch module~data_film_hd->module~data_arch module~data_film_hd->module~fluid_law module~mesh->module~data_arch module~sort_arrays->module~data_arch module~gen_param->module~data_arch module~mumfpack->iso_c_binding module~sulu_wrapper->iso_c_binding

Used by

  • module~~film~~UsedByGraph module~film film module~ms_film ms_film module~ms_film->module~film module~test_musst test_musst module~test_musst->module~film module~test_musst->module~ms_film module~inout_files inout_files module~test_musst->module~inout_files module~inout_files->module~film module~inout_files->module~ms_film program~main main program~main->module~test_musst

Contents


Variables

TypeVisibility AttributesNameInitial
integer(kind=I4), public, parameter:: H1_N =1

code for bottom surface height

integer(kind=I4), public, parameter:: H2_N =2

code for top surface height

integer(kind=I4), public, parameter:: H_N =3

code for film thickness

integer(kind=I4), public, parameter:: P_N =4

code for film absolute pressure

integer(kind=I4), public, parameter:: RHO_N =5

code for fluid density

integer(kind=I4), public, parameter:: T_N =6

code for film absolute temperature

integer(kind=I4), public, parameter:: DRHODP_N =7

code for film compressibility

integer(kind=I4), public, parameter:: MU_N =8

code for fluid viscosity

integer(kind=I4), public, parameter:: VX_N =9

code for surface velocity along . Should be modified (surfaces 1 and 2)

integer(kind=I4), public, parameter:: VY_N =10

code for surface velocity along . Should be modified (surfaces 1 and 2)

integer(kind=I4), public, parameter:: HG_C =1

code for groove depth on the stationnary surface , cell variable

integer(kind=I4), public, parameter:: PEK_C =2

code for Peclet number along

integer(kind=I4), public, parameter:: PEE_C =3

code for Peclet number along

integer(kind=I4), public, parameter:: REY =1

code for imposed pressure at the boundary

integer(kind=I4), public, parameter:: ASS =1

code for assembly of the system

integer(kind=I4), public, parameter:: NO_ASS =0

code for no assemble, computation of the residual only

integer(kind=I4), public, parameter:: NO_BC =-1

code for computation of the residual everywhere, even at the boundaries (fluxes computations)

type(SCALE_SURF), private :: scal_tmp

object SCALE_SURF

logical(kind=I4), public, parameter:: BC_SPLINE =.false.

instead of linearly interpolating the boundary pressures, it can be done in a smoother way. NOT TO BE USED YET.


Derived Types

type, public :: PRC_TAB

PRC_TAB stores some precomputed coefficients for the finite element matrices

Components

TypeVisibility AttributesNameInitial
integer(kind=I4), public :: ng

Gauss point number along a direction

real(kind=R8), public, dimension(:), allocatable:: pg

point coordinates in a direction

real(kind=R8), public, dimension(:), allocatable:: wg

point weight

real(kind=R8), public, dimension(:,:,:), allocatable:: vni4

for each node, shape function at Gauss points

real(kind=R8), public, dimension(:,:,:), allocatable:: vni4x

for each node, shape function derivative at Gauss points

real(kind=R8), public, dimension(:,:,:), allocatable:: vni4y

for each node, shape function derivative at Gauss points

real(kind=R8), public, dimension(:,:,:), allocatable:: vni4d

for each node, upwind shape function at Gauss points

real(kind=R8), public, dimension(:,:,:), allocatable:: vcal

14 values calculated at the Gauss points

type, public :: FE_FILM

FE_FILM stores the whole stuff related to a film: the nodal variables, the mesh, etc.

Components

TypeVisibility AttributesNameInitial
integer(kind=I4), public :: n_vn

number of nodal variables

integer(kind=I4), public :: n_vc

number of variables on cells

type(FE_MESH), public :: m

mesh of the film

type(DATA_FILM), public :: data_f

data of the problem

type(PRC_TAB), public :: prc

precomputation

type(NUM_PAR), public :: num_p

numerical param for iterative solution

real(kind=R8), public, dimension(:,:), allocatable:: vn

nodal variables table

real(kind=R8), public, dimension(:,:), allocatable:: vc

cell variables table

integer(kind=I4), public, dimension(:,:), allocatable:: bc

boundary nodes code

character(len=20), public, dimension(:), allocatable:: vn_name

nodal variable names

character(len=20), public, dimension(:), allocatable:: vc_name

cell variable names

Type-Bound Procedures

procedure, public :: fx

force computation along

procedure, public :: fy

force computation along

procedure, public :: fz

force computation along


Functions

private function fz(fe_f)

Read more…

Arguments

Type IntentOptional AttributesName
class(FE_FILM), intent(inout) :: fe_f

fluid type

Return Value real(kind=R8)

private function fx(fe_f)

Read more…

Arguments

Type IntentOptional AttributesName
class(FE_FILM), intent(inout) :: fe_f

FE film

Return Value real(kind=R8)

private function fy(fe_f)

Read more…

Arguments

Type IntentOptional AttributesName
class(FE_FILM), intent(inout) :: fe_f

FE film

Return Value real(kind=R8)


Subroutines

public subroutine create_rect_FE_film(data_f, num_p, fe_f)

Read more…

Arguments

Type IntentOptional AttributesName
type(DATA_FILM), intent(inout) :: data_f

data of the film

type(NUM_PAR), intent(in) :: num_p

numerical param for iterative solution

type(FE_FILM), intent(inout) :: fe_f

FE film data type

public subroutine solve_FE_film(fe_f, mat, bc, flag_ass)

Read more…

Arguments

Type IntentOptional AttributesName
type(FE_FILM), intent(inout) :: fe_f

FE film data type

type(MAT_SOLV), intent(inout) :: mat

matrices for solving

real(kind=R8), intent(in), dimension(MAX_NNC):: bc

boundary conditions at the corners

logical(kind=I4), intent(in), optional :: flag_ass

optional flag for assembly

public subroutine apply_bc_FE_film(fe_f, bc)

Read more…

Arguments

Type IntentOptional AttributesName
type(FE_FILM), intent(inout) :: fe_f

FE film data type

real(kind=R8), intent(in), dimension(MAX_NNC):: bc

boundary conditions at the corners

public subroutine apply_bc_FE_film_simple(fe_f, bc)

Read more…

Arguments

Type IntentOptional AttributesName
type(FE_FILM), intent(inout) :: fe_f

FE film data type

real(kind=R8), intent(in), dimension(MAX_NNC):: bc

boundary conditions at the corners

public subroutine assembly_FE_film_reynolds(fe_f, mat, ass_c)

Read more…

Arguments

Type IntentOptional AttributesName
type(FE_FILM), intent(inout) :: fe_f

FE_film data type

type(MAT_SOLV), intent(inout) :: mat

matrices for solving

integer(kind=I4), intent(in) :: ass_c

assembly type

public subroutine assemble_in_mat_sol(mat, num, nelt, nline, tind, m_elt, compt)

Read more…

Arguments

Type IntentOptional AttributesName
type(MAT_SOLV), intent(inout) :: mat

mat_solv type

integer(kind=I4), intent(in) :: num

element number

integer(kind=I4), intent(in) :: nelt

size of elemental matrix

integer(kind=I4), intent(in) :: nline

number of lines

integer(kind=I4), intent(in), dimension(nelt):: tind

index table of elemental matrix

real(kind=R8), intent(in), dimension(nelt, nelt):: m_elt

elemental matrix

integer(kind=I4), intent(inout), dimension(2):: compt

number to index the position in the solver matrix

private subroutine elementary_assembly_FE_film_reynolds(fe_f, ke_ij, be_i, ind_e, e, ass_c)

Read more…

Arguments

Type IntentOptional AttributesName
type(FE_FILM), intent(inout) :: fe_f

FE film

real(kind=R8), intent(out), dimension(MAX_NNE, MAX_NNE):: ke_ij

elementary matrix

real(kind=R8), intent(out), dimension(MAX_NNE):: be_i

elementary rhs member

integer(kind=I4), intent(out), dimension(MAX_NNE):: ind_e

elementary index member

integer(kind=I4), intent(in) :: e

element number

integer(kind=I4), intent(in) :: ass_c

assembly type

public subroutine elementary_full_domain_FE_film_reynolds(fe_f, mat, ke_ij, be_i, ind_e)

Read more…

Arguments

Type IntentOptional AttributesName
type(FE_FILM), intent(inout) :: fe_f

FE film

type(MAT_SOLV), intent(inout) :: mat

solver type matrices

real(kind=R8), intent(out), dimension(MAX_NNC, MAX_NNC):: ke_ij

elementary matrix

real(kind=R8), intent(inout), dimension(MAX_NNC ):: be_i

elementary rhs member

integer(kind=I4), intent(out), dimension(MAX_NNC ):: ind_e

elementary index member

public subroutine init_prc_tab(fe_f)

Read more…

Arguments

Type IntentOptional AttributesName
type(FE_FILM), intent(inout) :: fe_f

private subroutine compute_prc_tables_reynolds_supg(fe_f, e)

Read more…

Arguments

Type IntentOptional AttributesName
type(FE_FILM), intent(inout) :: fe_f

FE film

integer(kind=I4), intent(in) :: e

element number

private subroutine length_width_elem(spdx, spdy, x, y, length, width)

Read more…

Arguments

Type IntentOptional AttributesName
real(kind=R8), intent(in) :: spdx

fluid velocity along axis

real(kind=R8), intent(in) :: spdy

fluid velocity along axis

real(kind=R8), intent(in), dimension(4):: x

corner abscissae

real(kind=R8), intent(in), dimension(4):: y

corner ordinates

real(kind=R8), intent(out) :: length

fluid element length

real(kind=R8), intent(out) :: width

fluid element width

public subroutine compute_corner_fluxes(fe_f, mat, bf)

Read more…

Arguments

Type IntentOptional AttributesName
type(FE_FILM), intent(inout) :: fe_f

FE film

type(MAT_SOLV), intent(inout) :: mat

solver type matrices

real(kind=R8), intent(out), dimension(MAX_NNC):: bf

table of fluxes at the corner

public subroutine save_fe_field(fe_f, file_name, code, nodal)

Read more…

Arguments

Type IntentOptional AttributesName
type(FE_FILM), intent(in) :: fe_f
character(len=*), intent(in) :: file_name
integer(kind=I4), intent(in) :: code
logical(kind=I4), intent(in) :: nodal

if false : cell value, if true : nodal value