Functions
Alternatively to using the GUI, one has access to all the functions after installing the package.
API Functions
rtmsim.start_rtmsim
— Methodstart_rtmsim(inputfilename)
Reads the text input file and calls the solver with the read parameters.
Arguments:
- inputfilename :: String
The complete set of input parameters can be accessed in the input file. The following paragraph shows an example for such an input file:
1 #i_model
meshfiles/mesh_permeameter1_foursets.bdf #meshfilename
200 #tmax
1.01325e5 1.225 1.4 0.06 #p_ref rho_ref gamma mu_resin_val
1.35e5 1.0e5 #p_a_val p_init_val
3e-3 0.7 3e-10 1 1 0 0 #t_val porosity_val K_val alpha_val refdir1_val refdir2_val refdir3_val
3e-3 0.7 3e-10 1 1 0 0 #t1_val porosity1_val K1_val alpha1_val refdir11_val refdir21_val refdir31_val
3e-3 0.7 3e-10 1 1 0 0 #t2_val porosity2_val K2_val alpha2_val refdir12_val refdir22_val refdir32_val
3e-3 0.7 3e-10 1 1 0 0 #t3_val porosity3_val K3_val alpha3_val refdir13_val refdir23_val refdir33_val
3e-3 0.7 3e-10 1 1 0 0 #t4_val porosity4_val K4_val alpha4_val refdir14_val refdir24_val refdir34_val
1 0 0 0 #patchtype1val patchtype2val patchtype3val patchtype4val
0 results.jld2 #i_restart restartfilename
0 0.01 #i_interactive r_p
16 #n_pics
Meaning of the variables:
i_model
: Identifier for physical model (Default value is 1)meshfilename
: Mesh filename.tmax
: Maximum simulation time.p_ref rho_ref gamma mu_resin_val
: Parameters for the adiabatic equation of state and dynamic viscosity of resin used in the Darcy term.p_a_val p_init_val
: Absolut pressure value for injection port and for initial cavity pressure.t_val porosity_val K_val alpha_val refdir1_val refdir2_val refdir3_val
: Properties of the cells in the main preform: The vector(refdir1_val,refdir2_val,refdir3_val)
is projected onto the cell in order to define the first principal cell direction. The second principal cell direction is perpendicular to the first one in the plane spanned by the cell nodes. The principal cell directions are used as the principal permeabilty directions. The cell properties are defined by the thicknesst_val
, the porosityporosity_val
, the permeabilityK_val
in the first principal cell direction, the permeablityalpha_val
in the second principal direction.t1_val porosity1_val K1_val alpha1_val refdir11_val refdir21_val refdir31_val
etc.: Properties for up to four additional cell regions if preform.patchtype1val patchtype2val patchtype3val patchtype4val
: These regions are used to specify the location of the pressure boundary conditions and to specify regions with different permeability, porosity and thickness properties (e.g. for different part thickness and layup or for race tracking which are regions with very high permeability typically at the boundary of the preforms). Vents need not be specified. Parameterspatchtype1val
define the patch type. Numerical values 0, 1, 2 and 3 are allowed with the following interpretation:- 0 .. the patch is ignored
- 1 .. the patch represents an inlet gate, where the specified injection pressure level applies
- 2 .. the patch specifies a preform region
- 3 .. the patch represents a vent, where the specified initial pressure level applies
i_restart restartfilename
: Start with new simulation if0
or continue previous simulation if1
from specified filei_interactive r_p
: Select the inlet ports graphically if i_interactive equal to1
and inlet ports have specified radiusn_pics
: Number of intermediate output files, supposed to be a multiple of4
Entries are separated by one blank.
Unit test:
MODULE_ROOT=splitdir(splitdir(pathof(rtmsim))[1])[1]; inputfilename=joinpath(MODULE_ROOT,"inputfiles","input.txt"); rtmsim.start_rtmsim(inputfilename);
rtmsim.rtmsim_rev1
— Methodrtmsim_rev1(param)
RTMsim solver with the following main steps:
- Simulation initialization
- Read mesh file and prepare patches
- Find neighbouring cells
- Assign parameters to cells
- Create local cell coordinate systems
- Calculate initial time step
- Array initialization
- Define simulation time and intermediate output times
- Boundary conditions
- (Optional initialization if
i_model=2,3,..
) - Time evolution (for loops over all indices inside a while loop for time evolution)
- Calculation of correction factors for cell thickness, porosity, permeability, viscosity
- Pressure gradient calculation
- Numerical flux function calculation
- Update of rho, u, v, gamma and p according to conservation laws and equation of state
- Boundary conditions
- Prepare arrays for next time step
- Saving of intermediate data
- (Opional time marching etc. for i_model=2,3,...)
- Calculation of adaptive time step
Arguments: Data structure param
with
i_model :: Int64 =1
meshfilename :: String ="input.txt"
tmax :: Float64 =200
p_ref :: Float64 =1.01325e5
rho_ref :: Float64 =1.225
gamma :: Float64 =1.4
mu_resin_val :: Float64 =0.06
p_a_val :: Float64 =1.35e5
p_init_val :: Float64 =1.0e5
t_val :: Float64 =3e-3
porosity_val :: Float64 =0.7
K_val :: Float64 =3e-10
alpha_val :: Float64 =1.0
refdir1_val :: Float64 =1.0
refdir2_val :: Float64 =0.0
refdir3_val :: Float64 =0.0
t1_val :: Float64 =3e-3
porosity1_val :: Float64 =0.7
K1_val :: Float64 =3e-10
alpha1_val :: Float64 =1.0
refdir11_val :: Float64 =1.0
refdir21_val :: Float64 =0.0
refdir31_val :: Float64 =0.0
t2_val :: Float64 =3e-3
porosity2_val :: Float64 =0.7
K2_val :: Float64 =3e-10
alpha2_val :: Float64 =1.0
refdir12_val :: Float64 =1.0
refdir22_val :: Float64 =0.0
refdir32_val :: Float64 =0.0
t3_val :: Float64 =3e-3
porosity3_val :: Float64 =0.7
K3_val :: Float64 =3e-10
alpha3_val :: Float64 =1.0
refdir13_val :: Float64 =0.0
refdir23_val :: Float64 =0.0
refdir33_val :: Float64 =0.0
t4_val :: Float64 =3e-3
porosity4_val :: Float64 =0.7
K4_val :: Float64 =3e-10
alpha4_val :: Float64 =1.0
refdir14_val :: Float64 =1.0
refdir24_val :: Float64 =0.0
refdir34_val :: Float64 =0.0
patchtype1val :: Int64 =1
patchtype2val :: Int64 =0
patchtype3val :: Int64 =0
patchtype4val :: Int64 =0
i_restart :: Int64 =0
restartfilename :: String ="results.jld2"
i_interactive :: Int64 =0
r_p :: Float64 =0.01
n_pics :: Int64 =16
Meaning of the variables:
i_model
: Identifier for physical model (Default value is 1)meshfilename
: Mesh filename.tmax
: Maximum simulation time.p_ref rho_ref gamma mu_resin_val
: Parameters for the adiabatic equation of state and dynamic viscosity of resin used in the Darcy term.p_a_val p_init_val
: Absolut pressure value for injection port and for initial cavity pressure.t_val porosity_val K_val alpha_val refdir1_val refdir2_val refdir3_val
: Properties of the cells in the main preform: The vector(refdir1_val,refdir2_val,refdir3_val)
is projected onto the cell in order to define the first principal cell direction. The second principal cell direction is perpendicular to the first one in the plane spanned by the cell nodes. The principal cell directions are used as the principal permeabilty directions. The cell properties are defined by the thicknesst_val
, the porosityporosity_val
, the permeabilityK_val
in the first principal cell direction, the permeablityalpha_val
in the second principal direction.t1_val porosity1_val K1_val alpha1_val refdir11_val refdir21_val refdir31_val
etc.: Properties for up to four additional cell regions if preform.patchtype1val patchtype2val patchtype3val patchtype4val
: These regions are used to specify the location of the pressure boundary conditions and to specify regions with different permeability, porosity and thickness properties (e.g. for different part thickness and layup or for race tracking which are regions with very high permeability typically at the boundary of the preforms). Vents need not be specified. Parameterspatchtype1val
define the patch type. Numerical values 0, 1, 2 and 3 are allowed with the following interpretation:- 0 .. the patch is ignored
- 1 .. the patch represents an inlet gate, where the specified injection pressure level applies
- 2 .. the patch specifies a preform region
- 3 .. the patch represents a vent, where the specified initial pressure level applies
i_restart restartfilename
: Start with new simulation if0
or continue previous simulation if1
from specified filei_interactive r_p
: Select the inlet ports graphically if i_interactive equal to1
and inlet ports have specified radiusn_pics
: Number of intermediate output files, supposed to be a multiple of4
Unit tests:
MODULE_ROOT=splitdir(splitdir(pathof(rtmsim))[1])[1]; meshfilename=joinpath(MODULE_ROOT,"meshfiles","mesh_permeameter1_foursets.bdf"); param=rtmsim.input_vals(1,meshfilename,200, 101325,1.225,1.4,0.06, 1.35e5,1.00e5, 3e-3,0.7,3e-10,1,1,0,0, 3e-3,0.7,3e-10,1,1,0,0, 3e-3,0.7,3e-11,1,1,0,0, 3e-3,0.7,3e-11,1,1,0,0, 3e-3,0.7,3e-9,1,1,0,0, 1,2,2,2,0,"results.jld2",0,0.01,16); rtmsim.rtmsim_rev1(param);
Addtional unit tests:
MODULE_ROOT=splitdir(splitdir(pathof(rtmsim))[1])[1]; meshfilename=joinpath(MODULE_ROOT,"meshfiles","mesh_permeameter1_foursets.bdf"); param=rtmsim.input_vals(1,meshfilename,200, 101325,1.225,1.4,0.06, 1.35e5,1.00e5, 3e-3,0.7,3e-10,1,1,0,0, 3e-3,0.7,3e-10,1,1,0,0, 3e-3,0.7,3e-11,1,1,0,0, 3e-3,0.7,3e-11,1,1,0,0, 3e-3,0.7,3e-9,1,1,0,0, 1,0,0,0, 0,"results.jld2",0,0.01,16); rtmsim.rtmsim_rev1(param);
for starting a simulation with one pressure inlet port (sets 2, 3 and 4 are not used and consequently the preform parameters are ignored; since set 1 is a pressure inlet, also the parameters for set 1 are ignored and the only relevant parameter for the specified set is the pressure difference between injection and initial cavity pressure)
MODULE_ROOT=splitdir(splitdir(pathof(rtmsim))[1])[1]; meshfilename=joinpath(MODULE_ROOT,"meshfiles","mesh_permeameter1_foursets.bdf"); param=rtmsim.input_vals(1,meshfilename,200, 101325,1.225,1.4,0.06, 1.35e5,1.00e5, 3e-3,0.7,3e-10,1,1,0,0, 3e-3,0.7,3e-10,1,1,0,0, 3e-3,0.7,3e-11,1,1,0,0, 3e-3,0.7,3e-11,1,1,0,0, 3e-3,0.7,3e-9,1,1,0,0, 1,2,2,2, 0,"results.jld2",0,0.01,16); rtmsim.rtmsim_rev1(param);
for starting a simulation with different patches and race tracking
MODULE_ROOT=splitdir(splitdir(pathof(rtmsim))[1])[1]; meshfilename=joinpath(MODULE_ROOT,"meshfiles","mesh_permeameter1_foursets.bdf"); param=rtmsim.input_vals(1,meshfilename,200, 101325,1.225,1.4,0.06, 1.35e5,1.00e5, 3e-3,0.7,3e-10,1,1,0,0, 3e-3,0.7,3e-10,1,1,0,0, 3e-3,0.7,3e-11,1,1,0,0, 3e-3,0.7,3e-11,1,1,0,0, 3e-3,0.7,3e-9,1,1,0,0, 1,2,2,2, 1,"results.jld2",0,0.01,16); rtmsim.rtmsim_rev1(param);
for continuing the previous simulation
rtmsim.plot_mesh
— Methodfunction plot_mesh(meshfilename,i_mode)
Create mesh plot with cells with i_mode==1
and create mesh plots with cell center nodes with i_mode==2
for manual selection of inlet ports.
Arguments:
- meshfilename :: String
- i_mode :: Int
Unit test:
MODULE_ROOT=splitdir(splitdir(pathof(rtmsim))[1])[1]; meshfilename=joinpath(MODULE_ROOT,"meshfiles","mesh_permeameter1_foursets.bdf"); rtmsim.plot_mesh(meshfilename,1);
Additional unit tests:
MODULE_ROOT=splitdir(splitdir(pathof(rtmsim))[1])[1]; meshfilename=joinpath(MODULE_ROOT,"meshfiles","mesh_permeameter1_foursets.bdf"); rtmsim.plot_mesh(meshfilename,2);
for the manual selection of inlet ports with left mouse button click while key p is pressed
MODULE_ROOT=splitdir(splitdir(pathof(rtmsim))[1])[1]; meshfilename=joinpath(MODULE_ROOT,"meshfiles","mesh_permeameter1_foursets.bdf"); param=rtmsim.input_vals(1,meshfilename,200, 0.35e5,1.205,1.4,0.06, 0.35e5,0.00e5, 3e-3,0.7,3e-10,1,1,0,0, 3e-3,0.7,3e-10,1,1,0,0, 3e-3,0.7,3e-10,1,1,0,0, 3e-3,0.7,3e-10,1,1,0,0, 3e-3,0.7,3e-10,1,1,0,0, 0,0,0,0, 0,"results.jld2",1,0.01,16); rtmsim.rtmsim_rev1(param);
for starting only with the interactively selected inlet ports
rtmsim.plot_sets
— Methodfunction plot_sets(meshfilename)
Create a plot with the up to four cell sets defined in the mesh file.
Unit test:
MODULE_ROOT=splitdir(splitdir(pathof(rtmsim))[1])[1]; meshfilename=joinpath(MODULE_ROOT,"meshfiles","mesh_permeameter1_foursets.bdf"); rtmsim.plot_sets(meshfilename);
rtmsim.read_nastran_mesh
— Methodfunction read_nastran_mesh(input_struct)
Read file in Nastran format with fixed length (8 digits), nodes (GRIDS
) defined in global coordinate system.
Arguments:
meshfilename :: String
paramset, paramset1, paramset2, paramset3, paramset3 :: Vector{Float}
patchtype1val,patchtype1val1,patchtype1val2,patchtype1val3,patchtype1val4 :: Int
i_interactive :: Int
r_p :: Float
Unit test:
MODULE_ROOT=splitdir(splitdir(pathof(rtmsim))[1])[1]; meshfilename=joinpath(MODULE_ROOT,"meshfiles","mesh_permeameter1_foursets.bdf"); paramset=[0.5,0.3,3e-10,1.0,1.0,0.0,0.0];paramset1=paramset;paramset2=paramset;paramset3=paramset;paramset4=paramset;patchtype1val=-1;patchtype2val=-1;patchtype3val=-1;patchtype4val=-1;i_interactive=0;r_p=0.01; input_struct=rtmsim.input_args_read_mesh(meshfilename,paramset,paramset1,paramset2,paramset3,paramset4,patchtype1val,patchtype2val,patchtype3val,patchtype4val,i_interactive,r_p); return_struct=rtmsim.read_mesh(input_struct);
rtmsim.plot_results
— Methodfunction plot_results(resultsfilename)
Create contour plots of the filling factor and the pressure after loading a results file.
Arguments:
- resultsfilename :: String
Unit test:
WORK_DIR=pwd(); resultsfilename=joinpath(WORK_DIR,"results.jld2"); rtmsim.plot_results(resultsfilename);
rtmsim.plot_overview
— Methodfunction plot_overview(n_out,n_pics)
Create filling factor contour plots. n_out
is the index of the last output file, if n_out==-1
the output file with the highest index is chosen. Consider the last n_pics
for creating the contour plots at four equidistant time intervals, if n_pics==-1
all available output files are considered.
Arguments:
- n_out :: Int
- n_pics :: Int
Unit test:
rtmsim.plot_overview(-1,-1)
rtmsim.plot_filling
— Methodfunction plot_filling(n_out,n_pics)
Create a window showing the filling factor contour plot at a selected time instance. Selection is with slider bar. n_out
is the index of the last output file, if n_out==-1
the output file with the highest index is chosen. Consider the last n_pics
for creating the contour plots at four equidistant time intervals, if n_pics==-1
all available output files are considered.
Arguments:
- n_out :: Int
- n_pics :: Int
Unit test:
rtmsim.plot_filling(-1,-1)
rtmsim.gui
— Methodfunction gui()
Opens the GUI.
No arguments.
Unit test:
rtmsim.gui();
to open the GUI.
rtmsim.numerical_gradient
— Methodfunction numerical_gradient(input_struct)
Calculates the pressure gradient from the cell values of the neighbouring cells.
i_method=1
.. Least square solution to determine gradienti_method=2
.. Least square solution to determine gradient with limiteri_method=3
.. RUntime optimized least square solution to determine gradient
Arguments: Data structure with
i_method :: Int
ind :: Int
p_old :: Vector{Float}
cellneighoursarray :: Array{Float,2}
cellcentertocellcenterx, cellcentertocellcentery :: Array{Float,2}
Return: Data structure with
dpdx :: Float
dpdy :: Float
rtmsim.numerical_flux_function
— Methodfunction numerical_flux_function(input_struct)
Evaluates the numerical flux functions at the cell boundaries. -i_method==1
.. first order upwinding
Arguments: Data structure with
i_method :: Int
vars_P, vars_A :: 4-element Vector{Float}
meshparameters :: 3-element Vector{Float}
Return: Data structure with
F_rho_num_add :: Float
F_u_num_add :: Float
F_v_num_add :: Float
F_gamma_num_add :: Float
F_gamma_num1_add :: Float
Unit tests:
input_struct=rtmsim.input_args_flux(1,[1.0; 1.2; 0.0; 0.0],[1.0; 0.4; 0.0; 0.0],[1.0;0.0;1.0]);return_struct=rtmsim.numerical_flux_function(input_struct);println("F_rho="*string(return_struct.F_rho_num_add)*", F_u="*string(return_struct.F_u_num_add))
gives a comparison with the one-dimensional numerical flux for the continuity and momentum equations. The result for
i_mehtod=1
isF_rho=0.8, F_u=0.96
. The first is calculated without upwindingF_rho=rho*u=0.5*(1.0+1.0)*0.5*(1.2+0.4)=0.8
, the second is calculated without upwinding in the first factor and with upwinding in the second factorF_u=(rho*u)*(u)=0.8*1.2=0.96
. Mass density and velocity are[rho,u]=[1.0,1.2]
in the considered and[rho,u]=[1.0,0.4]
in the neighbouring cells.y
velocity and fluid fractiongamma
are not considered in this test case.input_struct=rtmsim.input_args_flux(1,[1.225; 1.2; 0.4; 0.9],[1.0; 0.4; 1.2; 0.1],[1/sqrt(2);1/sqrt(2);1.0]);return_struct=rtmsim.numerical_flux_function(input_struct);println("F_rho="*string(return_struct.F_rho_num_add)*", F_u="*string(return_struct.F_u_num_add)*", F_v="*string(return_struct.F_v_num_add)*", F_gamma="*string(return_struct.F_gamma_num_add)*", F_gamma_num1="*string(return_struct.F_gamma_num1_add))
is a general test case with
[rho,u,v,gamma]=[1.225,1.2,0.4,0.9]
in the considered cell,[rho,u,v,gamma]=[1.225,1.2,0.4,0.9]
in the neighbouring cell, outwards pointing cell normal vector[1/sqrt(2),1/sqrt(2)]
and boundary face area1.0
. The result fori_mehtod=1
isF_rho=1.2586500705120547, F_u=1.5103800846144655, F_v=0.5034600282048219, F_gamma=1.0182337649086284, F_gamma_num1=1.131370849898476
.
rtmsim.numerical_flux_function_boundary
— Methodnumerical_flux_function_boundary(input_struct)
Evaluates the numerical flux functions at the cell boundaries to pressure inlet or outlet.
i_method==1
.. first order upwinding
Arguments: Data structure with
i_method :: Int
vars_P, vars_A ::
4-elementVector{Float}
meshparameters ::
3-elementVector{Float}
n_dot_u :: Float
Return: Data structure with
F_rho_num_add :: Float
F_u_num_add :: Float
F_v_num_add :: Float
F_gamma_num_add :: Float
F_gamma_num1_add :: Float
Unit tests:
input_struct=rtmsim.input_args_flux_boundary(1,[1.0; 1.2; 0.0; 0.0],[1.0; 0.0; 0.0; 0.0],[1.0;0.0;1.0],1.2);return_struct=rtmsim.numerical_flux_function(input_struct);println("F_rho="*string(return_struct.F_rho_num_add)*", F_u="*string(return_struct.F_u_num_add))
gives a comparison with the one-dimensional inflow from a cell with mass density
rho=1.0
with velocityu=1.2
at the cell boundary. Mass density and velocity in the considered cell arerho=1.0
andu=1.2
.y
velocity and fluid fractiongamma
are not considered in this test case. The result fori_mehtod=1
isF_rho=0.6, F_u=0.72
.input_struct=rtmsim.input_args_flux_boundary(1,[1.225; 1.2; 0.4; 0.9],[1.0; 0.4; 1.2; 0.1],[1/sqrt(2);1/sqrt(2);1.0],-0.8);return_struct=rtmsim.numerical_flux_function_boundary(input_struct);println("F_rho="*string(return_struct.F_rho_num_add)*", F_u="*string(return_struct.F_u_num_add)*", F_v="*string(return_struct.F_v_num_add)*", F_gamma="*string(return_struct.F_gamma_num_add)*", F_gamma_num1="*string(return_struct.F_gamma_num1_add))
is a general case. The result for
i_mehtod=1
isF_rho=-0.8900000000000001, F_u=-0.3560000000000001, F_v=-1.068, F_gamma=-0.08000000000000002, F_gamma_num1=-0.8
.
rtmsim.create_coordinate_systems
— Methodfunction create_coordinate_systems(input_struct)
Define the local cell coordinate system and the transformation matrix from the local cell coordinate system from the neighbouring cell to the local cell coordinate system of the considered cell.
Arguments: Data structure with
N :: Int
cellgridid :: Array{Int,2}
gridx :: Vector{Float}
gridy :: Vector{Float}
gridz :: Vector{Float}
cellcenterx :: Vector{Float}
cellcentery :: Vector{Float}
cellcenterz :: Vector{Float}
faces :: Array{Int,2}
cellneighboursarray :: Array{Int,2}
celldirection :: Array{Float,2}
cellthickness :: Vector{Float}
maxnumberofneighbours :: Int
Meaning of the arguments:
N
: Number of cellscellgridid
: The i-th line contains the three IDs of the nodes which form the cellgridx
,gridy
,gridz
: i-th component of these vectors contain the x, y and z coordinates of node with ID icellcenterx
,cellcentery
,cellcenterz
: i-th component of these vectors contain the x, y and z coordinates of geometric cell centers of cell with ID ifaces
: Array with three columns. The three entries n1, n2, c1 in a cell are indices and describe a line between nodes with IDs n1 and n2 and this boundary belongs to cell with ID c1.cellneighboursarray
: The i-th line contains the indices of the neighbouring cells of cell with ID i. Number of columns is given by maxnumberofneighbours. Array is initialized with -9 and only the positive entries are considered.celldirection
: The i-th line contains the x, y and z coordinates of the unit normal vector which is projected on the cell to define the cell coordinate system for cell with ID i.cellthickness
: The i-th line contains the thickness of cell with ID i.maxnumberofneighbours
: Number of columns of arraycellneighboursarray
. Default value is10
. If more cell neighbours in the mesh, an error occurs and this value must be increased.
Return: Data structure with
cellvolume :: Vector{Float}
cellcentertocellcenterx :: Array{Float,2}
cellcentertocellcentery :: Array{Float,2}
T11 :: Array{Float,2}
T12 :: Array{Float,2}
T21 :: Array{Float,2}
T22 :: Array{Float,2}
cellfacenormalx :: Array{Float,2}
cellfacenormaly :: Array{Float,2}
cellfacearea :: Array{Float,2}
rtmsim.assign_parameters
— Methodfunction assign_parameters(input_struct)
Assign properties to cells.
Arguments: Data structure with
i_interactive :: Int
celltype :: Vector{Int}
patchparameters0 :: Vector{Float}
patchparameters1 :: Vector{Float}
patchparameters2 :: Vector{Float}
patchparameters3 :: Vector{Float}
patchparameters4 :: Vector{Float}
patchtype1val :: Int
patchtype2val :: Int
patchtype3val :: Int
patchtype4val :: Int
patchids1 :: Vector{Int}
patchids2 :: Vector{Int}
patchids3 :: Vector{Int}
patchids4 :: Vector{Int}
inletpatchids :: Vector{Int}
mu_resin_val :: Float
N :: Int
Meaning of the arguments:
i_interactive
: 1..if pressure inlet cells are selected manually, 0..elsecelltype[i]
: Describes cell type. -1..pressure inlet, -2..pressure outlet, -3..cell with wall boundary, 1..interior cellpatchparameters0
,patchparameters1
,patchparameters2
,patchparameters3
,patchparameters4
:: 7-element vector with parameters for the main preform and the four sets if used as patch. The seven elements are cellporosity,cellthickness,cellpermeability,cellalpha,celldirection[1],celldirection[2],celldirection[3].patchtype1val
,patchtype2val
,patchtype3val
,patchtype4val
:: Type of patch. 1..patchids1
,patchids2
,patchids3
,patchids4
:: Type of set. 0..ignored, 1..pressure inlet, 2..patch, 3..pressure outletinletpatchids
:: Vector with the IDs of the cells which are inlet cells.mu_resin_val
:: Kinematic viscosity value.N
: Number of cells.
Return: Data structure with
cellthickness :: Vector{Float64}
cellporosity :: Vector{Float64}
cellpermeability :: Vector{Float64}
cellalpha :: Vector{Float64}
celldirection :: Array{Float64,2}
cellviscosity :: Vector{Float64}
celltype :: Vector{Int64}
rtmsim.create_faces
— Methodfunction create_faces(input_struct)
Find the set with the IDs of the neighbouring cells and identify wall cells.
Arguments: Data structure with
cellgridid :: Array{Int,2}
N :: Int
maxnumberofneighbours :: Int
Meaning of the arguments:
cellgridid
: The i-th line contains the three IDs of the nodes which form the cellN
: Number of cells.maxnumberofneighbours
: Number of columns of arraycellneighboursarray
. Default value is10
. If more cell neighbours in the mesh, an error occurs and this value must be increased.
Return: Data structure with
faces :: Array{Int,2}
cellneighboursarray :: Array{Int,2}
celltype :: Vector{Int}
Additional features
The source code is prepared for the following extensions:
- Import mesh file in different format. Selection is based on the extension of the mesh file.
- Input parameter
i_model
(for iso-thermal RTM=1
) is used for adding additional functionalities. E.g. adding temperature and degree-of-cure equations with variable resin viscosity or for VARI with variable porosity, permeability and cavity thickness. - Parameter
i_method
in the functions for numerical differentiation and flux functions can be used to implement different numerical schemes. E.g. gradient limiter or second-order upwinding.
If you are interested, please have a look at the contribution item in the community standards.