ICFERST  22-06
Reservoir simulator based on DCVFEM, Dynamic Mesh optimisation and Surface-based modelling
cv_advection Module Reference

This module contains all the tools to assemble and solve the equations and fields associated with the CV mesh, i.e. transport equation, continuity equation, DCVFE gradient matrix and the Laplacian system for the zeta-potential. More...

Data Types

interface  dg_derivs_all
 calculates derivatives of vector fields More...
 
interface  pack_loc_all
 Packs field together to be later on used for high order computations. More...
 

Functions/Subroutines

subroutine cv_assemb (state, packed_state, final_phase, Mdims, CV_GIdims, CV_funs, Mspars, ndgln, Mdisopt, Mmat, upwnd, tracer, velocity, density, multi_absorp, DIAG_SCALE_PRES, DIAG_SCALE_PRES_COUP, INV_B, DEN_ALL, DENOLD_ALL, CV_DISOPT, CV_DG_VEL_INT_OPT, DT, CV_THETA, CV_BETA, SUF_SIG_DIAGTEN_BC, DERIV, CV_P, SOURCT_ALL, ABSORBT_ALL, VOLFRA_PORE, GETCV_DISC, GETCT, IGOT_T2, IGOT_THETA_FLUX, GET_THETA_FLUX, USE_THETA_FLUX, THETA_FLUX, ONE_M_THETA_FLUX, THETA_FLUX_J, ONE_M_THETA_FLUX_J, THETA_GDIFF, MEAN_PORE_CV, MASS_MN_PRES, THERMAL, got_free_surf, MASS_SUF, MASS_ELE_TRANSP, TDIFFUSION, saturation, VAD_parameter, Phase_with_Pc, Courant_number, Permeability_tensor_field, calculate_mass_delta, eles_with_pipe, pipes_aux, porous_heat_coef, porous_heat_coef_old, outfluxes, solving_compositional, nonlinear_iteration, assemble_collapsed_to_one_phase)
 This subroutines generates the transport equation for a cv field. It also can generate the Continuity equation if GETCT = .true. and also generate the gradient matrix of the momentum equation for the DCVFE method if the option is activated. More...
 
pure subroutine onvdlim_ano_many (NFIELD, TDLIM, TDCEN, INCOME, ETDNEW_PELE, ETDNEW_PELEOT, XI_LIMIT, TUPWIN, TUPWI2, DENOIN, CTILIN, DENOOU, CTILOU, FTILIN, FTILOU)
 This sub calculates the limited face values TDADJ(1...SNGI) from the central difference face values TDCEN(1...SNGI) using a NVD shceme. INCOME(1...SNGI)=1 for incomming to element ELE else =0. LIBETA is the flux limiting parameter. TDMAX(PELE)=maximum of the surrounding 6 element values of element PELE. TDMIN(PELE)=minimum of the surrounding 6 element values of element PELE. PELEOT=element at other side of current face. ELEOT2=element at other side of the element ELEOTH. ELESID=element next to oposing current face. DENOIN, CTILIN, DENOOU, CTILOU, FTILIN, FTILOU => memory The elements are arranged in this order: ELEOT2,ELE, PELEOT, ELESID. More...
 
subroutine is_field_constant (IGOT_T_CONST, IGOT_T_CONST_VALUE, T_ALL, CV_NONODS)
 Checks whether a field is constant or not. Although I think this check is terrible, no alternatives have provided the same functionality so far. More...
 
subroutine pack_loc (LOC_F, T_ALL, NPHASE, IPT, IGOT_T_PACK)
 Copies memory into the same large array to then perform a projection from CV to FE. More...
 
subroutine unpack_loc (LOC_F, T_ALL, NPHASE, IPT, IGOT_T_PACK, IGOT_T_CONST, IGOT_T_CONST_VALUE)
 If PACK then UNpack loc_f into T_ALL as long at IGOT_T==1 and STORE and not already in storage. More...
 
subroutine pack_or_unpack_loc (LOC_F, T_ALL, NPHASE, NFIELD, IPT, PACK, STORE, IGOT_T)
 If PACK then pack T_ALL into LOC_F as long at IGOT_T==1 and STORE and not already in storage. More...
 
subroutine pack_loc_all1 (LOC_F, field1, oldfield1, field2, oldfield2, field3, oldfield3, IGOT_T_PACK, use_volume_frac_T2, nfield)
 This subroutine is for fields that have already size final_phase - 1 Checks if the fields are constant or not, stored in IGOT_T_PACK, based on that introduces the field into LOC_F to later on apply the limiters on all the fields at once. More...
 
subroutine pack_loc_all2 (LOC_F, field1, oldfield1, field2, oldfield2, field3, oldfield3, IGOT_T_PACK, use_volume_frac_T2, start_phase, final_phase, nodi)
 This subrotuine is for fields that are bigger than final_phase - start_phase Checks if the fields are constant or not, stored in IGOT_T_PACK, based on that introduces the field into LOC_F to later on apply the limiters on all the fields at once. More...
 
subroutine pack_loc_all3 (LOC_F, field1, oldfield1, field2, oldfield2, field3, oldfield3, IGOT_T_PACK, use_volume_frac_T2, start_phase, final_phase, nodi)
 This subrotuine is for fields that are bigger than final_phase - start_phase Checks if the fields are constant or not, stored in IGOT_T_PACK, based on that introduces the field into LOC_F to later on apply the limiters on all the fields at once. This one is for integer fields. More...
 
subroutine i_pack_loc (LOC_F, T_ALL, NPHASE, IPT, IGOT_T_PACK)
 If PACK then pack T_ALL into LOC_F as long at IGOT_T==1 and STORE and not already in storage. More...
 
subroutine unpack_loc_all (LOC_F, field1, oldfield1, field2, oldfield2, field3, oldfield3, IGOT_T_PACK, IGOT_T_CONST, IGOT_T_CONST_VALUE, use_volume_frac_T2, nfield)
 If PACK then UNpack loc_f into T_ALL as long at IGOT_T==1 and STORE and not already in storage. Checks if the fields are constant or not, stored in IGOT_T_PACK, based on that introduces the LOC_F into the field or use a constant value. This is after the limiters have been applied. More...
 
integer function cv_count_faces (Mdims, CV_ELE_TYPE, CV_GIdims)
 This subroutine counts then number of faces in the control volume space. More...
 
subroutine find_other_side (CV_OTHER_LOC, CV_NLOC, U_OTHER_LOC, U_NLOC, MAT_OTHER_LOC, INTEGRAT_AT_GI, X_NLOC, XU_NLOC, X_NDGLN, XU_NDGLN, CV_SNLOC, CVFEM_ON_FACE, X_SHARE, ELE, ELE2, FINELE, COLELE, DISTCONTINUOUS_METHOD)
 We are on the boundary or next to another element. Determine CV_OTHER_LOC, U_OTHER_LOC. CVFEM_ON_FACE(CV_KLOC,GI)=.TRUE. if CV_KLOC is on the face that GI is centred on. Look for these nodes on the other elements. ELE2=0 also when we are between elements but are trying to integrate across the middle of a CV. More...
 
subroutine proj_cv_to_fem (packed_state, fempsi, psi, Mdims, CV_GIdims, CV_funs, Mspars, ndgln, igetct, X, mass_ele, mass_mn_pres, tracer, psi_ave, psi_int, activate_limiters)
 Calls to generate the transport equation for the saturation. Embeded an FPI with backtracking method is uncluded Subroutine description: (1) determine FEMT (finite element wise) etc from T (control volume wise) (2) (optional) calculate psi_int (area) and psi_ave (barycentre) over each CV. More...
 
subroutine dg_derivs_all1 (FEMT, FEMTOLD, DTX_ELE, DTOLDX_ELE, NDIM, NPHASE, NCOMP, TOTELE, CV_NDGLN, XCV_NDGLN, X_NLOC, X_NDGLN, CV_NGI, CV_NLOC, CVWEIGHT, N, NLX, NLY, NLZ, X_N, X_NLX, X_NLY, X_NLZ, X_NONODS, X, Y, Z, NFACE, FACE_ELE, CV_SLOCLIST, X_SLOCLIST, CV_SNLOC, X_SNLOC, WIC_T_BC, SUF_T_BC, SBCVNGI, SBCVFEN, SBWEIGH, X_SBCVFEN, X_SBCVFENSLX, X_SBCVFENSLY, get_gradU, state, P0Mesh)
 calculates derivatives of vector fields More...
 
subroutine dg_derivs_all2 (FEMT, FEMTOLD, DTX_ELE, DTOLDX_ELE, NDIM, NPHASE, TOTELE, CV_NDGLN, XCV_NDGLN, X_NLOC, X_NDGLN, CV_NGI, CV_NLOC, CVWEIGHT, N, NLX, NLY, NLZ, X_N, X_NLX, X_NLY, X_NLZ, X_NONODS, X, Y, Z, NFACE, FACE_ELE, CV_SLOCLIST, X_SLOCLIST, CV_SNLOC, X_SNLOC, WIC_T_BC, SUF_T_BC, SBCVNGI, SBCVFEN, SBWEIGH, X_SBCVFEN, X_SBCVFENSLX, X_SBCVFENSLY, P0Mesh)
 Computes the derivatives of vector fields. More...
 
subroutine diffus_cal_coeff (DIFF_COEF_DIVDX, DIFF_COEFOLD_DIVDX, CV_NLOC, MAT_NLOC, NPHASE, MAT_NDGLN, SMATFEN, SCVFEN, GI, NDIM, TDIFFUSION, HDC, T_CV_NODJ, T_CV_NODI, TOLD_CV_NODJ, TOLD_CV_NODI, ELE, ELE2, CVNORMX_ALL, LOC_DTX_ELE_ALL, LOC_DTOLDX_ELE_ALL, LOC2_DTX_ELE_ALL, LOC2_DTOLDX_ELE_ALL, LOC_WIC_T_BC, CV_OTHER_LOC, MAT_OTHER_LOC, CV_SNLOC, CV_SLOC2LOC, on_domain_boundary, between_elements, has_anisotropic_diffusivity)
 This sub calculates the effective diffusion coefficientd DIFF_COEF_DIVDX, DIFF_COEFOLD_DIVDX based on a non-linear method and a non-oscillating scheme. It requires the derivatives of the field obtained using DG_DERIVS_ALL DG_DERIVS_ALL. More...
 
subroutine linear_high_diffus_cal_coeff_stress_or_tensor (STRESS_IJ_ELE_EXT, S_INV_NNX_MAT12, STRESS_FORM, STRESS_FORM_STAB, ZERO_OR_TWO_THIRDS, U_SNLOC, U_NLOC, CV_SNLOC, NPHASE, SBUFEN_REVERSED, SBCVFEN_REVERSED, SDETWEI, SBCVNGI, NDIM, SLOC_UDIFFUSION, SLOC_UDIFFUSION_VOL, SLOC2_UDIFFUSION, SLOC2_UDIFFUSION_VOL, DIFF_GI_ADDED, ON_BOUNDARY, SNORMXN_ALL)
 This sub calculates the effective diffusion coefficientd STRESS_IJ_ELE_EXT it only works for between element contributions. based on a high order scheme. The matrix S_INV_NNX_MAT12 is used to calculate the rows of the matrix with STRESS_IJ_ELE_EXT. This implements the stress and tensor form of diffusion and calculates a jump conidition. which is in DIFF_COEF_DIVDX, DIFF_COEFOLD_DIVDX The coefficient are in N_DOT_DKDU, N_DOT_DKDUOLD. look at the manual DG treatment of viscocity. More...
 
subroutine calc_stress_ten (STRESS_IJ, ZERO_OR_TWO_THIRDS, NDIM, UFENX_ILOC, UFENX_JLOC, TEN_XX, TEN_VOL)
 determine stress form of viscocity... More...
 
subroutine calc_stress_ten_reduce (STRESS_IJ, ZERO_OR_TWO_THIRDS, NDIM, FEN_TEN_XX, FEN_TEN_VOL, UFENX_JLOC)
 determine stress form of viscocity... More...
 
subroutine cal_lim_vol_adjust (TMIN_STORE, TMIN, T, TMIN_NOD, RESET_STORE, MASS_CV, CV_NODI_IPHA, CV_NODJ_IPHA, IPHASE, CV_NONODS, INCOME)
 Adjust TMIN to take into account different sized CV's. if RESET_STORE then reset TMIN to orginal values. More...
 
subroutine dgsimplnorm_all (NLOC, SNLOC, NDIM, XL_ALL, XSL_ALL, NORMX_ALL)
 Form approximate surface normal (NORMX_ALL(1),NORMX_ALL(2),NORMX_ALL(3)) More...
 
subroutine dgsimplnorm (ELE, SILOC2ILOC, NLOC, SNLOC, XONDGL, X, Y, Z, NORMX, NORMY, NORMZ)
 Form approximate surface normal (NORMX,NORMY,NORMZ) More...
 
subroutine isotropic_limiter_all (T_ALL, TOLD_ALL, T2_ALL, T2OLD_ALL, DEN_ALL, DENOLD_ALL, IGOT_T2, NPHASE, CV_NONODS, nsmall_colm, SMALL_CENTRM, SMALL_FINDRM, SMALL_COLM, STOTEL, CV_SNLOC, CV_SNDGLN, SUF_T_BC_ALL, SUF_T2_BC_ALL, SUF_D_BC_ALL, WIC_T_BC_ALL, WIC_T2_BC_ALL, WIC_D_BC_ALL, MASS_CV, TOLDUPWIND_MAT_ALL, DENOLDUPWIND_MAT_ALL, T2OLDUPWIND_MAT_ALL, TUPWIND_MAT_ALL, DENUPWIND_MAT_ALL, T2UPWIND_MAT_ALL)
 Computes the limited values at the interface, not as generic as the anisotropic one and never used actually... More...
 
subroutine calc_sele (ELE, ELE3, SELE, CV_SILOC, CV_ILOC, U_SLOC2LOC, CV_SLOC2LOC, FACE_ELE, gi, CV_funs, Mdims, CV_GIdims, CV_NDGLN, U_NDGLN, CV_SNDGLN, U_SNDGLN)
 Calculate surface element, surface control volume: SELE, CV_SILOC, U_SLOC2LOC, CV_SLOC2LOC for a face on the boundary of the domain. More...
 
subroutine put_in_ct_rhs (GET_C_IN_CV_ADVDIF_AND_CALC_C_CV, ct_rhs_phase_cv_nodi, ct_rhs_phase_cv_nodj, final_phase, Mdims, CV_funs, ndgln, Mmat, GI, between_elements, on_domain_boundary, ELE, ELE2, SELE, HDC, MASS_ELE, JCOUNT_KLOC, JCOUNT_KLOC2, ICOUNT_KLOC, ICOUNT_KLOC2, C_JCOUNT_KLOC, C_JCOUNT_KLOC2, C_ICOUNT_KLOC, C_ICOUNT_KLOC2, U_OTHER_LOC, U_SLOC2LOC, CV_SLOC2LOC, SCVDETWEI, CVNORMX_ALL, DEN_ALL, CV_NODI, CV_NODJ, WIC_U_BC_ALL, WIC_P_BC_ALL, SUF_P_BC_ALL, UGI_COEF_ELE_ALL, UGI_COEF_ELE2_ALL, NDOTQ, NDOTQOLD, LIMT, LIMDT, LIMDTOLD, LIMT_HAT, NDOTQ_HAT, FTHETA_T2, ONE_M_FTHETA_T2OLD, FTHETA_T2_J, ONE_M_FTHETA_T2OLD_J, integrate_other_side_and_not_boundary, loc_u, THETA_VEL, UDGI_IMP_ALL, RCON, RCON_J, NDOTQ_IMP, X_ALL, SUF_D_BC_ALL, gravty)
 This subroutine caculates the discretised cty eqn acting on the velocities i.e. MmatCT, MmatCT_RHS It also computes the gradient matrix using the DCVFE method. More...
 
subroutine calc_anisotrop_lim (Mmat, T_ALL, TOLD_ALL, DEN_ALL, DENOLD_ALL, T2_ALL, T2OLD_ALL, FEMT_ALL, FEMTOLD_ALL, FEMDEN_ALL, FEMDENOLD_ALL, FEMT2_ALL, FEMT2OLD_ALL, USE_FEMT, TUPWIND_MAT_ALL, TOLDUPWIND_MAT_ALL, DENUPWIND_MAT_ALL, DENOLDUPWIND_MAT_ALL, T2UPWIND_MAT_ALL, T2OLDUPWIND_MAT_ALL, IGOT_T2, NPHASE, CV_NONODS, CV_NLOC, TOTELE, CV_NDGLN, SMALL_FINDRM, SMALL_CENTRM, SMALL_COLM, NSMALL_COLM, X_NDGLN, X_NONODS, NDIM, X_ALL, XC_CV_ALL, use_reflect)
 For the anisotropic limiting scheme we find the upwind values by interpolation using the subroutine FINPTS or IFINPTS; the upwind value for each node pair is stored in the matrices TUPWIND AND. More...
 
subroutine triloccords (Xp, Yp, Zp, N1, N2, N3, N4, X1, Y1, Z1, X2, Y2, Z2, X3, Y3, Z3, X4, Y4, Z4)
 
subroutine triloccords2d (Xp, Yp, N1, N2, N3, X1, Y1, X2, Y2, X3, Y3)
 
logical function shock_front_in_ele (ele, Mdims, sat, ndgln, Imble_frac)
 Detects whether the element has a shockfront or not. More...
 
subroutine generate_laplacian_system (Mdims, packed_state, ndgln, Mmat, Mspars, CV_funs, CV_GIdims, Sigma_field, Solution, K_fields, F_fields, intface_val_type)
 : In this method we assemble and solve the Laplacian system using at least P1 elements The equation solved is the following: Div sigma Grad X = - SUM (Div K Grad F) with Neuman BCs = 0 where K and F are passed down as a vector. Therefore for n entries the SUM will be performed over n fields Example: F = (3, nphase, cv_nonods) would include three terms in the RHS and the same for K If harmonic average then we perform the harmonic average of sigma and K IMPORTANT: This subroutine requires the PHsparsity to be generated Note that this method solves considering FE fields. If using CV you may incur in an small error. More...
 
subroutine scvdetnx (Mdims, ndgln, X_ALL, CV_funs, CV_GIdims, on_domain_boundary, between_elements, ELE, GI, SCVDETWEI, CVNORMX_ALL, XC_ALL, X_NOD, X_NODJ)
 

Variables

integer, parameter wic_t_bc_dirichlet = 1
 
integer, parameter wic_t_bc_robin = 2
 
integer, parameter wic_t_bc_diri_adv_and_robin = 3
 
integer, parameter wic_d_bc_dirichlet = 1
 
integer, parameter wic_u_bc_dirichlet = 1
 
integer, parameter wic_u_bc_robin = 2
 
integer, parameter wic_u_bc_diri_adv_and_robin = 3
 
integer, parameter wic_u_bc_dirichlet_inout = 2
 
integer, parameter wic_p_bc_dirichlet = 1
 
integer, parameter wic_p_bc_free = 2
 

Detailed Description

This module contains all the tools to assemble and solve the equations and fields associated with the CV mesh, i.e. transport equation, continuity equation, DCVFE gradient matrix and the Laplacian system for the zeta-potential.

Function/Subroutine Documentation

◆ cal_lim_vol_adjust()

subroutine cv_advection::cal_lim_vol_adjust ( real, intent(inout)  TMIN_STORE,
real, dimension( : ), intent(inout)  TMIN,
real, dimension( : ), intent(in)  T,
integer, dimension( : ), intent(in)  TMIN_NOD,
logical, intent(in)  RESET_STORE,
real, dimension( : ), intent(in)  MASS_CV,
integer, intent(in)  CV_NODI_IPHA,
integer, intent(in)  CV_NODJ_IPHA,
integer, intent(in)  IPHASE,
integer, intent(in)  CV_NONODS,
real, intent(in)  INCOME 
)

Adjust TMIN to take into account different sized CV's. if RESET_STORE then reset TMIN to orginal values.

Here is the caller graph for this function:

◆ calc_anisotrop_lim()

subroutine cv_advection::calc_anisotrop_lim ( type (multi_matrices), intent(inout)  Mmat,
real, dimension( :, : ), intent(in)  T_ALL,
real, dimension( :, : ), intent(in)  TOLD_ALL,
real, dimension( :, : ), intent(in)  DEN_ALL,
real, dimension( :, : ), intent(in)  DENOLD_ALL,
real, dimension( :,:), intent(in), pointer  T2_ALL,
real, dimension( :,:), intent(in), pointer  T2OLD_ALL,
real, dimension( :, :), intent(in)  FEMT_ALL,
real, dimension( :, :), intent(in)  FEMTOLD_ALL,
real, dimension( :, :), intent(in)  FEMDEN_ALL,
real, dimension( :, :), intent(in)  FEMDENOLD_ALL,
real, dimension( :, :), intent(in)  FEMT2_ALL,
real, dimension( :, :), intent(in)  FEMT2OLD_ALL,
logical, intent(in)  USE_FEMT,
real, dimension( :, : ), intent(inout)  TUPWIND_MAT_ALL,
real, dimension( :, : ), intent(inout)  TOLDUPWIND_MAT_ALL,
real, dimension( :, : ), intent(inout)  DENUPWIND_MAT_ALL,
real, dimension( :, : ), intent(inout)  DENOLDUPWIND_MAT_ALL,
real, dimension( :, : ), intent(inout)  T2UPWIND_MAT_ALL,
real, dimension( :, : ), intent(inout)  T2OLDUPWIND_MAT_ALL,
integer, intent(in)  IGOT_T2,
integer, intent(in)  NPHASE,
integer, intent(in)  CV_NONODS,
integer, intent(in)  CV_NLOC,
integer, intent(in)  TOTELE,
integer, dimension( : ), intent(in)  CV_NDGLN,
integer, dimension( : ), intent(in)  SMALL_FINDRM,
integer, dimension( : ), intent(in)  SMALL_CENTRM,
integer, dimension( : ), intent(in)  SMALL_COLM,
integer, intent(in)  NSMALL_COLM,
integer, dimension(: ), intent(in)  X_NDGLN,
integer, intent(in)  X_NONODS,
integer, intent(in)  NDIM,
real, dimension( :,: ), intent(in)  X_ALL,
real, dimension( ndim, cv_nonods), intent(in)  XC_CV_ALL,
logical, intent(in)  use_reflect 
)

For the anisotropic limiting scheme we find the upwind values by interpolation using the subroutine FINPTS or IFINPTS; the upwind value for each node pair is stored in the matrices TUPWIND AND.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ calc_sele()

subroutine cv_advection::calc_sele ( integer, intent(in)  ELE,
integer, intent(inout)  ELE3,
integer, intent(inout)  SELE,
integer, intent(inout)  CV_SILOC,
integer, intent(in)  CV_ILOC,
integer, dimension( : ), intent(inout)  U_SLOC2LOC,
integer, dimension( : ), intent(inout)  CV_SLOC2LOC,
integer, dimension( :, : ), intent(in)  FACE_ELE,
integer, intent(in)  gi,
type(multi_shape_funs), intent(in)  CV_funs,
type(multi_dimensions), intent(in)  Mdims,
type(multi_gi_dimensions), intent(in)  CV_GIdims,
integer, dimension( : ), intent(in)  CV_NDGLN,
integer, dimension( : ), intent(in)  U_NDGLN,
integer, dimension( : ), intent(in)  CV_SNDGLN,
integer, dimension( : ), intent(in)  U_SNDGLN 
)

Calculate surface element, surface control volume: SELE, CV_SILOC, U_SLOC2LOC, CV_SLOC2LOC for a face on the boundary of the domain.

Here is the caller graph for this function:

◆ calc_stress_ten()

subroutine cv_advection::calc_stress_ten ( real, dimension( :, : ), intent(inout)  STRESS_IJ,
real, intent(in)  ZERO_OR_TWO_THIRDS,
integer, intent(in)  NDIM,
real, dimension( : ), intent(in)  UFENX_ILOC,
real, dimension( : ), intent(in)  UFENX_JLOC,
real, dimension( :,: ), intent(in)  TEN_XX,
real, intent(in)  TEN_VOL 
)

determine stress form of viscocity...

Here is the caller graph for this function:

◆ calc_stress_ten_reduce()

subroutine cv_advection::calc_stress_ten_reduce ( real, dimension( ndim, ndim ), intent(inout)  STRESS_IJ,
real, intent(in)  ZERO_OR_TWO_THIRDS,
integer, intent(in)  NDIM,
real, dimension( ndim, ndim ), intent(in)  FEN_TEN_XX,
real, dimension( ndim ), intent(in)  FEN_TEN_VOL,
real, dimension( ndim ), intent(in)  UFENX_JLOC 
)

determine stress form of viscocity...

Here is the caller graph for this function:

◆ cv_assemb()

subroutine cv_advection::cv_assemb ( type( state_type ), dimension( : ), intent(inout)  state,
type( state_type ), intent(inout)  packed_state,
integer, intent(in)  final_phase,
type(multi_dimensions), intent(in)  Mdims,
type(multi_gi_dimensions), intent(in)  CV_GIdims,
type(multi_shape_funs), intent(inout)  CV_funs,
type(multi_sparsities), intent(in)  Mspars,
type(multi_ndgln), intent(in)  ndgln,
type (multi_discretization_opts Mdisopt,
type (multi_matrices), intent(inout)  Mmat,
type (porous_adv_coefs), intent(inout)  upwnd,
type(tensor_field), intent(inout), target  tracer,
type(tensor_field), intent(in)  velocity,
type(tensor_field), intent(in), target  density,
type(multi_absorption), intent(inout)  multi_absorp,
real, dimension( :, : ), intent(inout), allocatable  DIAG_SCALE_PRES,
real, dimension( :, :, : ), intent(inout), allocatable  DIAG_SCALE_PRES_COUP,
real, dimension( :, :, : ), intent(inout), allocatable  INV_B,
real, dimension( :, : ), intent(inout), target  DEN_ALL,
real, dimension( :, : ), intent(inout)  DENOLD_ALL,
integer, intent(in)  CV_DISOPT,
integer, intent(in)  CV_DG_VEL_INT_OPT,
real, intent(in)  DT,
real, intent(in)  CV_THETA,
real, intent(in)  CV_BETA,
real, dimension( :, : ), intent(in)  SUF_SIG_DIAGTEN_BC,
real, dimension( :, : ), intent(in)  DERIV,
real, dimension( :, :, : ), intent(in)  CV_P,
real, dimension( :, : ), intent(in)  SOURCT_ALL,
real, dimension( :, :, : ), intent(in), pointer  ABSORBT_ALL,
real, dimension( :, : ), intent(in)  VOLFRA_PORE,
logical, intent(in)  GETCV_DISC,
logical, intent(in)  GETCT,
integer, intent(in)  IGOT_T2,
integer, intent(in)  IGOT_THETA_FLUX,
logical, intent(in)  GET_THETA_FLUX,
logical, intent(in)  USE_THETA_FLUX,
real, dimension( :, : ), intent(inout), optional  THETA_FLUX,
real, dimension( :, : ), intent(inout), optional  ONE_M_THETA_FLUX,
real, dimension( :, : ), intent(inout), optional  THETA_FLUX_J,
real, dimension( :, : ), intent(inout), optional  ONE_M_THETA_FLUX_J,
real, dimension( :, : ), intent(inout)  THETA_GDIFF,
real, dimension( :, : ), intent(inout)  MEAN_PORE_CV,
real, dimension( : ), intent(inout)  MASS_MN_PRES,
logical, intent(in)  THERMAL,
logical, intent(in)  got_free_surf,
real, dimension( : ), intent(inout)  MASS_SUF,
real, dimension( : ), intent(inout), optional  MASS_ELE_TRANSP,
real, dimension( :, :, :, : ), intent(in), optional  TDIFFUSION,
type(tensor_field), intent(in), optional, target  saturation,
real, dimension(:), intent(in), optional  VAD_parameter,
integer, intent(in), optional  Phase_with_Pc,
real, dimension(:), intent(inout), optional  Courant_number,
type( tensor_field ), intent(in), optional, pointer  Permeability_tensor_field,
real, dimension(:,:), optional  calculate_mass_delta,
type(pipe_coords), dimension(:), intent(in), optional  eles_with_pipe,
type (multi_pipe_package), intent(in)  pipes_aux,
real, dimension( :), intent(in), optional  porous_heat_coef,
real, dimension( :), intent(in), optional  porous_heat_coef_old,
type (multi_outfluxes), intent(inout), optional  outfluxes,
logical, intent(in), optional  solving_compositional,
integer, intent(in), optional  nonlinear_iteration,
logical, intent(in), optional  assemble_collapsed_to_one_phase 
)

This subroutines generates the transport equation for a cv field. It also can generate the Continuity equation if GETCT = .true. and also generate the gradient matrix of the momentum equation for the DCVFE method if the option is activated.

Author
Chris Pain, Pablo Salinas
Parameters
stateLinked list containing all the fields defined in diamond and considered by Fluidity
packed_stateLinked list containing all the fields used by IC-FERST, memory partially shared with state
final_phaseThis is the final phase to be assembled, in this way we can assemble from phase 1 to final_phase not necessarily being for all the phases
MdimsData type storing all the dimensions describing the mesh, fields, nodes, etc
CV_GIdimsGauss integration numbers for CV fields
FE_GIdimsGauss integration numbers for FE fields
CV_funsShape functions for the CV mesh
FE_funsShape functions for the FE mesh
MsparsSparsity of the matrices
ndglnGlobal to local variables
MdisoptDiscretisation options
MmatMatrices for ICFERST
upwndSigmas to compute the fluxes at the interphase for porous media
tracerTracer considered for the transport equation
densityDensity of the field
velocityVelocity of the field
multi_absorpAbsoprtion of associated with the transport field
CV_DISOPT,CV_DG_VEL_INT_OPT,IGOT_THETA_FLUX.More Discretisation options
IGOT_T2.True (1) if solving for a tracer, false otherwise
DIAG_SCALE_PRESDiagonal scaling of (distributed) pressure matrix (used to treat pressure implicitly)
DIAG_SCALE_PRES_COUPDiagonal scaling of (distributed) pressure matrix (for wells)
INV_BCoupling term of the wells
MASS_MN_PRES??
MASS_SUF??
DEN_ALLDensity of the field, different memory to the input field density, used to apply the Boussinesq approximation
DENOLD_ALLDensity of the field, different memory to the input field density, used to apply the Boussinesq approximation
THETA_GDIFF! (nphase,Mdimscv_nonods)
THETA_FLUX,ONE_M_THETA_FLUX,THETA_FLUX_J,ONE_M_THETA_FLUX_J
DT,CV_THETA,CV_BETATime step size and discretisation options
SUF_SIG_DIAGTEN_BC.Like upwnd but for the boundary
DERIVDerivative of the density against the pressure (nphase,Mdimscv_nonods)
CV_P! Control volume pressure, useless for the DCVFEM (1,Mdimsnpres,Mdimscv_nonods)
SOURCT_ALLSource term of the tracer equation
ABSORBT_ALLAbsorption to be used here
VOLFRA_POREPorosity field (Mdimsnpres,Mdimstotele)
GETCV_DISCobtain the transport equation
GETCTobtain the continuity equation
GET_THETA_FLUX,USE_THETA_FLUX,got_free_surf???
THERMALtrue if solving for heat transport
MEAN_PORE_CVPorosity defined control volume wise
MASS_ELE_TRANSPMass of the elements
saturationPhaseVolumeFraction field
TDIFFUSIONDiffusion associated with the tracer field
Phase_with_PcField that defines the capillary pressure, i.e. non-wetting phase
VAD_parameterVanishing artificial diffusion parameter
Courant_numberObvious ins't it?
Permeability_tensor_fieldPermeability field
calculate_mass_deltaVariable used to control the mass conservation of the system
eles_with_pipeElements that have a pipe
pipes_auxInformation required to define wells
porous_heat_coef,porous_heat_coef_oldincludes an average of porous and fluid heat coefficients
solving_compositionalSolving for composition if true
assemble_collapsed_to_one_phaseCollapses phases and solves for one single temperature. When there is thermal equilibrium
outfluxesContains all the fields required to compute the outfluxes of the model and create the outfluxes.csv file. Computed when assembling the continuity equation
nonlinear_iterationCurrent non-linear iteration
Here is the call graph for this function:
Here is the caller graph for this function:

◆ cv_count_faces()

integer function cv_advection::cv_count_faces ( type(multi_dimensions), intent(in)  Mdims,
integer, intent(in)  CV_ELE_TYPE,
type( multi_gi_dimensions ), intent(in), optional  CV_GIdims 
)

This subroutine counts then number of faces in the control volume space.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ dg_derivs_all1()

subroutine cv_advection::dg_derivs_all1 ( real, dimension( :, :, : ), intent(in)  FEMT,
real, dimension( :, :, : ), intent(in)  FEMTOLD,
real, dimension( :, :, :, :, : ), intent(inout)  DTX_ELE,
real, dimension( :, :, :, :, : ), intent(inout)  DTOLDX_ELE,
integer, intent(in)  NDIM,
integer, intent(in)  NPHASE,
integer, intent(in)  NCOMP,
integer, intent(in)  TOTELE,
integer, dimension( : ), intent(in)  CV_NDGLN,
integer, dimension( : ), intent(in)  XCV_NDGLN,
integer, intent(in)  X_NLOC,
integer, dimension( : ), intent(in)  X_NDGLN,
integer, intent(in)  CV_NGI,
integer, intent(in)  CV_NLOC,
real, dimension( : ), intent(in)  CVWEIGHT,
real, dimension( :, : ), intent(in)  N,
real, dimension( :, : ), intent(in)  NLX,
real, dimension( :, : ), intent(in)  NLY,
real, dimension( :, : ), intent(in)  NLZ,
real, dimension( :, : ), intent(in)  X_N,
real, dimension( :, : ), intent(in)  X_NLX,
real, dimension( :, : ), intent(in)  X_NLY,
real, dimension( :, : ), intent(in)  X_NLZ,
integer, intent(in)  X_NONODS,
real, dimension( : ), intent(in)  X,
real, dimension( : ), intent(in)  Y,
real, dimension( : ), intent(in)  Z,
integer, intent(in)  NFACE,
integer, dimension( :, : ), intent(in)  FACE_ELE,
integer, dimension( :, : ), intent(in)  CV_SLOCLIST,
integer, dimension( :, : ), intent(in)  X_SLOCLIST,
integer, intent(in)  CV_SNLOC,
integer, intent(in)  X_SNLOC,
integer, dimension( :, :, : ), intent(in)  WIC_T_BC,
real, dimension( :, :, : ), intent(in)  SUF_T_BC,
integer, intent(in)  SBCVNGI,
real, dimension( :, : ), intent(in)  SBCVFEN,
real, dimension( : ), intent(in)  SBWEIGH,
real, dimension( :, : ), intent(in)  X_SBCVFEN,
real, dimension( :, : ), intent(in)  X_SBCVFENSLX,
real, dimension( :, : ), intent(in)  X_SBCVFENSLY,
logical, intent(in)  get_gradU,
type( state_type), dimension( : ), intent(inout)  state,
logical, optional  P0Mesh 
)

calculates derivatives of vector fields

Here is the call graph for this function:

◆ dg_derivs_all2()

subroutine cv_advection::dg_derivs_all2 ( real, dimension( :, : ), intent(in)  FEMT,
real, dimension( :, : ), intent(in)  FEMTOLD,
real, dimension( :, :, :, : ), intent(inout)  DTX_ELE,
real, dimension( :, :, :, : ), intent(inout)  DTOLDX_ELE,
integer, intent(in)  NDIM,
integer, intent(in)  NPHASE,
integer, intent(in)  TOTELE,
integer, dimension( : ), intent(in)  CV_NDGLN,
integer, dimension( : ), intent(in)  XCV_NDGLN,
integer, intent(in)  X_NLOC,
integer, dimension( : ), intent(in)  X_NDGLN,
integer, intent(in)  CV_NGI,
integer, intent(in)  CV_NLOC,
real, dimension( : ), intent(in)  CVWEIGHT,
real, dimension( :, : ), intent(in)  N,
real, dimension( :, : ), intent(in)  NLX,
real, dimension( :, : ), intent(in)  NLY,
real, dimension( :, : ), intent(in)  NLZ,
real, dimension( :, : ), intent(in)  X_N,
real, dimension( :, : ), intent(in)  X_NLX,
real, dimension( :, : ), intent(in)  X_NLY,
real, dimension( :, : ), intent(in)  X_NLZ,
integer, intent(in)  X_NONODS,
real, dimension( : ), intent(in)  X,
real, dimension( : ), intent(in)  Y,
real, dimension( : ), intent(in)  Z,
integer, intent(in)  NFACE,
integer, dimension( :, : ), intent(in)  FACE_ELE,
integer, dimension( :, : ), intent(in)  CV_SLOCLIST,
integer, dimension( :, : ), intent(in)  X_SLOCLIST,
integer, intent(in)  CV_SNLOC,
integer, intent(in)  X_SNLOC,
integer, dimension( :, :, : ), intent(in)  WIC_T_BC,
real, dimension( :, :, : ), intent(in)  SUF_T_BC,
integer, intent(in)  SBCVNGI,
real, dimension( :, : ), intent(in)  SBCVFEN,
real, dimension( : ), intent(in)  SBWEIGH,
real, dimension( :, : ), intent(in)  X_SBCVFEN,
real, dimension( :, : ), intent(in)  X_SBCVFENSLX,
real, dimension( :, : ), intent(in)  X_SBCVFENSLY,
logical, optional  P0Mesh 
)

Computes the derivatives of vector fields.

Here is the call graph for this function:

◆ dgsimplnorm()

subroutine cv_advection::dgsimplnorm ( integer, intent(in)  ELE,
integer, dimension( : ), intent(in)  SILOC2ILOC,
integer, intent(in)  NLOC,
integer, intent(in)  SNLOC,
integer, dimension( : ), intent(in)  XONDGL,
real, dimension( : ), intent(in)  X,
real, dimension( : ), intent(in)  Y,
real, dimension( : ), intent(in)  Z,
real, intent(inout)  NORMX,
real, intent(inout)  NORMY,
real, intent(inout)  NORMZ 
)

Form approximate surface normal (NORMX,NORMY,NORMZ)

Here is the caller graph for this function:

◆ dgsimplnorm_all()

subroutine cv_advection::dgsimplnorm_all ( integer, intent(in)  NLOC,
integer, intent(in)  SNLOC,
integer, intent(in)  NDIM,
real, dimension( ndim, nloc ), intent(in)  XL_ALL,
real, dimension( ndim, snloc ), intent(in)  XSL_ALL,
real, dimension( ndim ), intent(inout)  NORMX_ALL 
)

Form approximate surface normal (NORMX_ALL(1),NORMX_ALL(2),NORMX_ALL(3))

Here is the caller graph for this function:

◆ diffus_cal_coeff()

subroutine cv_advection::diffus_cal_coeff ( real, dimension( nphase ), intent(inout)  DIFF_COEF_DIVDX,
real, dimension( nphase ), intent(inout)  DIFF_COEFOLD_DIVDX,
integer, intent(in)  CV_NLOC,
integer, intent(in)  MAT_NLOC,
integer, intent(in)  NPHASE,
integer, dimension( : ), intent(in)  MAT_NDGLN,
real, dimension( :, : ), intent(in)  SMATFEN,
real, dimension( :, : ), intent(in)  SCVFEN,
integer, intent(in)  GI,
integer, intent(in)  NDIM,
real, dimension( :, :, :, : ), intent(in)  TDIFFUSION,
real, intent(in)  HDC,
real, dimension( nphase ), intent(in)  T_CV_NODJ,
real, dimension( nphase ), intent(in)  T_CV_NODI,
real, dimension( nphase ), intent(in)  TOLD_CV_NODJ,
real, dimension( nphase ), intent(in)  TOLD_CV_NODI,
integer, intent(in)  ELE,
integer, intent(in)  ELE2,
real, dimension( : ), intent(in)  CVNORMX_ALL,
real, dimension( ndim, nphase, cv_nloc ), intent(in)  LOC_DTX_ELE_ALL,
real, dimension( ndim, nphase, cv_nloc ), intent(in)  LOC_DTOLDX_ELE_ALL,
real, dimension( ndim, nphase, cv_nloc ), intent(in)  LOC2_DTX_ELE_ALL,
real, dimension( ndim, nphase, cv_nloc ), intent(in)  LOC2_DTOLDX_ELE_ALL,
integer, dimension( nphase ), intent(in)  LOC_WIC_T_BC,
integer, dimension( : ), intent(in)  CV_OTHER_LOC,
integer, dimension( : ), intent(in)  MAT_OTHER_LOC,
integer, intent(in)  CV_SNLOC,
integer, dimension( : ), intent(in)  CV_SLOC2LOC,
logical, intent(in)  on_domain_boundary,
logical, intent(in)  between_elements,
logical, intent(in)  has_anisotropic_diffusivity 
)

This sub calculates the effective diffusion coefficientd DIFF_COEF_DIVDX, DIFF_COEFOLD_DIVDX based on a non-linear method and a non-oscillating scheme. It requires the derivatives of the field obtained using DG_DERIVS_ALL DG_DERIVS_ALL.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ find_other_side()

subroutine cv_advection::find_other_side ( integer, dimension( : ), intent(inout)  CV_OTHER_LOC,
integer, intent(in)  CV_NLOC,
integer, dimension( : ), intent(inout)  U_OTHER_LOC,
integer, intent(in)  U_NLOC,
integer, dimension( : ), intent(inout)  MAT_OTHER_LOC,
logical, intent(inout)  INTEGRAT_AT_GI,
integer, intent(in)  X_NLOC,
integer, intent(in)  XU_NLOC,
integer, dimension( : ), intent(in)  X_NDGLN,
integer, dimension( : ), intent(in)  XU_NDGLN,
integer, intent(in)  CV_SNLOC,
logical, dimension( : ), intent(in)  CVFEM_ON_FACE,
logical, dimension( : ), intent(inout)  X_SHARE,
integer, intent(in)  ELE,
integer, intent(inout)  ELE2,
integer, dimension( : ), intent(in)  FINELE,
integer, dimension( : ), intent(in)  COLELE,
logical, intent(in)  DISTCONTINUOUS_METHOD 
)

We are on the boundary or next to another element. Determine CV_OTHER_LOC, U_OTHER_LOC. CVFEM_ON_FACE(CV_KLOC,GI)=.TRUE. if CV_KLOC is on the face that GI is centred on. Look for these nodes on the other elements. ELE2=0 also when we are between elements but are trying to integrate across the middle of a CV.

Here is the caller graph for this function:

◆ generate_laplacian_system()

subroutine cv_advection::generate_laplacian_system ( type(multi_dimensions), intent(in)  Mdims,
type( state_type ), intent(inout)  packed_state,
type(multi_ndgln), intent(in)  ndgln,
type (multi_matrices), intent(inout)  Mmat,
type(multi_sparsities), intent(in)  Mspars,
type(multi_shape_funs), intent(inout)  CV_funs,
type(multi_gi_dimensions), intent(in)  CV_GIdims,
real, dimension(:,:), intent(in)  Sigma_field,
type( scalar_field ), intent(inout)  Solution,
real, dimension(:,:,:), intent(in)  K_fields,
real, dimension(:,:,:), intent(in)  F_fields,
integer, intent(in)  intface_val_type 
)

: In this method we assemble and solve the Laplacian system using at least P1 elements The equation solved is the following: Div sigma Grad X = - SUM (Div K Grad F) with Neuman BCs = 0 where K and F are passed down as a vector. Therefore for n entries the SUM will be performed over n fields Example: F = (3, nphase, cv_nonods) would include three terms in the RHS and the same for K If harmonic average then we perform the harmonic average of sigma and K IMPORTANT: This subroutine requires the PHsparsity to be generated Note that this method solves considering FE fields. If using CV you may incur in an small error.

Parameters
MdimsData type storing all the dimensions describing the mesh, fields, nodes, etc
packed_stateLinked list containing all the fields used by IC-FERST, memory partially shared with state
ndglnGlobal to local variables
MmatMatrices for ICFERST
MsparsSparsity of the matrices
CV_funsShape functions for the CV mesh
CV_GIdimsGauss integration numbers for CV fields
Sigma_fieldCoefficient multipliying the left-hand side term
SolutionSolution field of type scalar_field, required to ensure consistency
K_fieldsCoefficients multipliying the RHS terms
F_fieldsRHS fierlds
intface_val_type0 = no interpolation; 1 Harmonic mean; !20 for SP solver, harmonic mean considering charge; negative normal mean
Here is the call graph for this function:

◆ i_pack_loc()

subroutine cv_advection::i_pack_loc ( integer, dimension(:), intent(inout)  LOC_F,
integer, dimension(nphase), intent(in)  T_ALL,
integer, intent(in)  NPHASE,
integer, intent(inout)  IPT,
logical, dimension(nphase), intent(in)  IGOT_T_PACK 
)

If PACK then pack T_ALL into LOC_F as long at IGOT_T==1 and STORE and not already in storage.

◆ is_field_constant()

subroutine cv_advection::is_field_constant ( logical  IGOT_T_CONST,
real  IGOT_T_CONST_VALUE,
real, dimension(cv_nonods)  T_ALL,
integer  CV_NONODS 
)

Checks whether a field is constant or not. Although I think this check is terrible, no alternatives have provided the same functionality so far.

◆ isotropic_limiter_all()

subroutine cv_advection::isotropic_limiter_all ( real, dimension( :, : ), intent(in)  T_ALL,
real, dimension( :, : ), intent(in)  TOLD_ALL,
real, dimension( :, : ), intent(in)  T2_ALL,
real, dimension( :, : ), intent(in)  T2OLD_ALL,
real, dimension( :, : ), intent(in)  DEN_ALL,
real, dimension( :, : ), intent(in)  DENOLD_ALL,
integer, intent(in)  IGOT_T2,
integer, intent(in)  NPHASE,
integer, intent(in)  CV_NONODS,
integer, intent(in)  nsmall_colm,
integer, dimension( : ), intent(in)  SMALL_CENTRM,
integer, dimension( : ), intent(in)  SMALL_FINDRM,
integer, dimension( : ), intent(in)  SMALL_COLM,
integer, intent(in)  STOTEL,
integer, intent(in)  CV_SNLOC,
integer, dimension( : ), intent(in)  CV_SNDGLN,
real, dimension( :, :, : ), intent(in), pointer  SUF_T_BC_ALL,
real, dimension( :, :, : ), intent(in), pointer  SUF_T2_BC_ALL,
real, dimension( :, :, : ), intent(in), pointer  SUF_D_BC_ALL,
integer, dimension( :,:, : ), intent(in)  WIC_T_BC_ALL,
integer, dimension( :,:, : ), intent(in)  WIC_T2_BC_ALL,
integer, dimension( :,:, : ), intent(in)  WIC_D_BC_ALL,
real, dimension( : ), intent(in)  MASS_CV,
real, dimension( :, : ), intent(inout)  TOLDUPWIND_MAT_ALL,
real, dimension( :, : ), intent(inout)  DENOLDUPWIND_MAT_ALL,
real, dimension( :, : ), intent(inout)  T2OLDUPWIND_MAT_ALL,
real, dimension( :, : ), intent(inout)  TUPWIND_MAT_ALL,
real, dimension( :, : ), intent(inout)  DENUPWIND_MAT_ALL,
real, dimension( :, : ), intent(inout)  T2UPWIND_MAT_ALL 
)

Computes the limited values at the interface, not as generic as the anisotropic one and never used actually...

Here is the call graph for this function:
Here is the caller graph for this function:

◆ linear_high_diffus_cal_coeff_stress_or_tensor()

subroutine cv_advection::linear_high_diffus_cal_coeff_stress_or_tensor ( real, dimension( ndim, ndim, nphase, u_snloc, 2*u_nloc ), intent(inout)  STRESS_IJ_ELE_EXT,
real, dimension( ndim, u_snloc, 2*u_nloc ), intent(inout)  S_INV_NNX_MAT12,
logical, intent(in)  STRESS_FORM,
logical, intent(in)  STRESS_FORM_STAB,
real, intent(in)  ZERO_OR_TWO_THIRDS,
integer, intent(in)  U_SNLOC,
integer, intent(in)  U_NLOC,
integer, intent(in)  CV_SNLOC,
integer, intent(in)  NPHASE,
real, dimension( sbcvngi, u_snloc ), intent(in)  SBUFEN_REVERSED,
real, dimension( sbcvngi, cv_snloc ), intent(in)  SBCVFEN_REVERSED,
real, dimension( sbcvngi ), intent(in)  SDETWEI,
integer, intent(in)  SBCVNGI,
integer, intent(in)  NDIM,
real, dimension( ndim,ndim,nphase,cv_snloc ), intent(in)  SLOC_UDIFFUSION,
real, dimension( nphase,cv_snloc ), intent(in)  SLOC_UDIFFUSION_VOL,
real, dimension( ndim,ndim,nphase,cv_snloc ), intent(in)  SLOC2_UDIFFUSION,
real, dimension( nphase,cv_snloc ), intent(in)  SLOC2_UDIFFUSION_VOL,
real, dimension( ndim, ndim,ndim, nphase, sbcvngi), intent(in)  DIFF_GI_ADDED,
logical, intent(in)  ON_BOUNDARY,
real, dimension( ndim, sbcvngi ), intent(in)  SNORMXN_ALL 
)

This sub calculates the effective diffusion coefficientd STRESS_IJ_ELE_EXT it only works for between element contributions. based on a high order scheme. The matrix S_INV_NNX_MAT12 is used to calculate the rows of the matrix with STRESS_IJ_ELE_EXT. This implements the stress and tensor form of diffusion and calculates a jump conidition. which is in DIFF_COEF_DIVDX, DIFF_COEFOLD_DIVDX The coefficient are in N_DOT_DKDU, N_DOT_DKDUOLD. look at the manual DG treatment of viscocity.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ onvdlim_ano_many()

pure subroutine cv_advection::onvdlim_ano_many ( integer, intent(in)  NFIELD,
real, dimension( nfield ), intent(inout)  TDLIM,
real, dimension( nfield ), intent(in)  TDCEN,
real, dimension( nfield ), intent(in)  INCOME,
real, dimension( nfield ), intent(in)  ETDNEW_PELE,
real, dimension( nfield ), intent(in)  ETDNEW_PELEOT,
real, dimension( nfield ), intent(in)  XI_LIMIT,
real, dimension( nfield ), intent(in)  TUPWIN,
real, dimension( nfield ), intent(in)  TUPWI2,
real, dimension( nfield ), intent(inout)  DENOIN,
real, dimension( nfield ), intent(inout)  CTILIN,
real, dimension( nfield ), intent(inout)  DENOOU,
real, dimension( nfield ), intent(inout)  CTILOU,
real, dimension( nfield ), intent(inout)  FTILIN,
real, dimension( nfield ), intent(inout)  FTILOU 
)

This sub calculates the limited face values TDADJ(1...SNGI) from the central difference face values TDCEN(1...SNGI) using a NVD shceme. INCOME(1...SNGI)=1 for incomming to element ELE else =0. LIBETA is the flux limiting parameter. TDMAX(PELE)=maximum of the surrounding 6 element values of element PELE. TDMIN(PELE)=minimum of the surrounding 6 element values of element PELE. PELEOT=element at other side of current face. ELEOT2=element at other side of the element ELEOTH. ELESID=element next to oposing current face. DENOIN, CTILIN, DENOOU, CTILOU, FTILIN, FTILOU => memory The elements are arranged in this order: ELEOT2,ELE, PELEOT, ELESID.

This sub finds the neighbouring elements. Suppose that this is the face IFACE.

| ELEOT2 | ELEOTH | ELE | ELESID |

TAIN THALF TAOUT

TEXTIN TEXOUT

Here is the call graph for this function:
Here is the caller graph for this function:

◆ pack_loc()

subroutine cv_advection::pack_loc ( real, dimension(:), intent(inout)  LOC_F,
real, dimension(:), intent(in)  T_ALL,
integer, intent(in)  NPHASE,
integer, intent(inout)  IPT,
logical, dimension(:), intent(in)  IGOT_T_PACK 
)

Copies memory into the same large array to then perform a projection from CV to FE.

Here is the caller graph for this function:

◆ pack_loc_all1()

subroutine cv_advection::pack_loc_all1 ( real, dimension(:), intent(inout)  LOC_F,
real, dimension(:), intent(in)  field1,
real, dimension(:), intent(in)  oldfield1,
real, dimension(:), intent(in)  field2,
real, dimension(:), intent(in)  oldfield2,
real, dimension(:), intent(in)  field3,
real, dimension(:), intent(in)  oldfield3,
logical, dimension(:,:), intent(in)  IGOT_T_PACK,
logical, intent(in)  use_volume_frac_T2,
integer, intent(in)  nfield 
)

This subroutine is for fields that have already size final_phase - 1 Checks if the fields are constant or not, stored in IGOT_T_PACK, based on that introduces the field into LOC_F to later on apply the limiters on all the fields at once.

Here is the call graph for this function:

◆ pack_loc_all2()

subroutine cv_advection::pack_loc_all2 ( real, dimension(:), intent(inout)  LOC_F,
real, dimension(:,:), intent(in)  field1,
real, dimension(:,:), intent(in)  oldfield1,
real, dimension(:,:), intent(in)  field2,
real, dimension(:,:), intent(in)  oldfield2,
real, dimension(:,:), intent(in)  field3,
real, dimension(:,:), intent(in)  oldfield3,
logical, dimension(:,:), intent(in)  IGOT_T_PACK,
logical, intent(in)  use_volume_frac_T2,
integer, intent(in)  start_phase,
integer, intent(in)  final_phase,
integer, intent(in)  nodi 
)

This subrotuine is for fields that are bigger than final_phase - start_phase Checks if the fields are constant or not, stored in IGOT_T_PACK, based on that introduces the field into LOC_F to later on apply the limiters on all the fields at once.

Here is the call graph for this function:

◆ pack_loc_all3()

subroutine cv_advection::pack_loc_all3 ( integer, dimension(:), intent(inout)  LOC_F,
integer, dimension(:,:), intent(in)  field1,
integer, dimension(:,:), intent(in)  oldfield1,
integer, dimension(:,:), intent(in)  field2,
integer, dimension(:,:), intent(in)  oldfield2,
integer, dimension(:,:), intent(in)  field3,
integer, dimension(:,:), intent(in)  oldfield3,
logical, dimension(:,:), intent(in)  IGOT_T_PACK,
logical, intent(in)  use_volume_frac_T2,
integer, intent(in)  start_phase,
integer, intent(in)  final_phase,
integer, intent(in)  nodi 
)

This subrotuine is for fields that are bigger than final_phase - start_phase Checks if the fields are constant or not, stored in IGOT_T_PACK, based on that introduces the field into LOC_F to later on apply the limiters on all the fields at once. This one is for integer fields.

Here is the call graph for this function:

◆ pack_or_unpack_loc()

subroutine cv_advection::pack_or_unpack_loc ( real, dimension(nfield), intent(inout)  LOC_F,
real, dimension(nphase), intent(inout)  T_ALL,
integer, intent(in)  NPHASE,
integer, intent(in)  NFIELD,
integer, intent(inout)  IPT,
logical, intent(in)  PACK,
logical, intent(in)  STORE,
integer, intent(in)  IGOT_T 
)

If PACK then pack T_ALL into LOC_F as long at IGOT_T==1 and STORE and not already in storage.

◆ proj_cv_to_fem()

subroutine cv_advection::proj_cv_to_fem ( type(state_type)  packed_state,
type(tensor_field_pointer), dimension(:), intent(inout)  fempsi,
type(tensor_field_pointer), dimension(:), intent(inout)  psi,
type(multi_dimensions), intent(in)  Mdims,
type(multi_gi_dimensions), intent(in)  CV_GIdims,
type(multi_shape_funs), intent(inout)  CV_funs,
type(multi_sparsities), intent(in)  Mspars,
type(multi_ndgln), intent(in)  ndgln,
integer, intent(in)  igetct,
real, dimension(:,:), intent(in)  X,
real, dimension(:), intent(inout)  mass_ele,
real, dimension(:), intent(inout)  mass_mn_pres,
type(tensor_field), intent(in)  tracer,
type(vector_field_pointer), dimension(:), intent(inout)  psi_ave,
type(vector_field_pointer), dimension(:), intent(inout)  psi_int,
logical, intent(in), optional  activate_limiters 
)

Calls to generate the transport equation for the saturation. Embeded an FPI with backtracking method is uncluded Subroutine description: (1) determine FEMT (finite element wise) etc from T (control volume wise) (2) (optional) calculate psi_int (area) and psi_ave (barycentre) over each CV.

Author
Chris Pain, Pablo Salinas
Parameters
packed_state! local state data
fempsi! finite element field data
psi! finite volume field data
Mdims! dimension data
CV_GIdims! gauss integer dimension data
CV_funs! control volume shape function data
Mspars! sparsity data
ndgln! global numbering data
igetct! whether to get CT matrix
X! coordinates of the elements
mass_ele! finite element mass
mass_mn_pres! ??
tracerfield to be projected
activate_limitersare the limiters on?
psi_int! control volume area
psi_ave! control volume barycentre
Here is the call graph for this function:
Here is the caller graph for this function:

◆ put_in_ct_rhs()

subroutine cv_advection::put_in_ct_rhs ( logical, intent(in)  GET_C_IN_CV_ADVDIF_AND_CALC_C_CV,
real, dimension( : ), intent(inout)  ct_rhs_phase_cv_nodi,
real, dimension( : ), intent(inout)  ct_rhs_phase_cv_nodj,
integer, intent(in)  final_phase,
type(multi_dimensions), intent(in)  Mdims,
type(multi_shape_funs), intent(in)  CV_funs,
type(multi_ndgln), intent(in)  ndgln,
type (multi_matrices), intent(inout)  Mmat,
integer, intent(in)  GI,
logical, intent(in)  between_elements,
logical, intent(in)  on_domain_boundary,
integer, intent(in)  ELE,
integer, intent(in)  ELE2,
integer, intent(in)  SELE,
real, intent(in)  HDC,
real, dimension( : ), intent(in)  MASS_ELE,
integer, dimension( : ), intent(in)  JCOUNT_KLOC,
integer, dimension( : ), intent(in)  JCOUNT_KLOC2,
integer, dimension( : ), intent(in)  ICOUNT_KLOC,
integer, dimension( : ), intent(in)  ICOUNT_KLOC2,
integer, dimension( : ), intent(in)  C_JCOUNT_KLOC,
integer, dimension( : ), intent(in)  C_JCOUNT_KLOC2,
integer, dimension( : ), intent(in)  C_ICOUNT_KLOC,
integer, dimension( : ), intent(in)  C_ICOUNT_KLOC2,
integer, dimension( : ), intent(in)  U_OTHER_LOC,
integer, dimension( : ), intent(in)  U_SLOC2LOC,
integer, dimension( : ), intent(in)  CV_SLOC2LOC,
real, dimension( : ), intent(in)  SCVDETWEI,
real, dimension( :, : ), intent(in)  CVNORMX_ALL,
real, dimension( :, : ), intent(in)  DEN_ALL,
integer, intent(in)  CV_NODI,
integer, intent(in)  CV_NODJ,
integer, dimension(:,:,:)  WIC_U_BC_ALL,
integer, dimension(:,:,:)  WIC_P_BC_ALL,
real, dimension( :, :, : ), intent(inout)  SUF_P_BC_ALL,
real, dimension( :, :, : ), intent(in)  UGI_COEF_ELE_ALL,
real, dimension( :, :, : ), intent(in)  UGI_COEF_ELE2_ALL,
real, dimension( : ), intent(in)  NDOTQ,
real, dimension( : ), intent(in)  NDOTQOLD,
real, dimension( : ), intent(in)  LIMT,
real, dimension( : ), intent(in)  LIMDT,
real, dimension( : ), intent(in)  LIMDTOLD,
real, dimension( : ), intent(in)  LIMT_HAT,
real, intent(in)  NDOTQ_HAT,
real, dimension( : ), intent(in)  FTHETA_T2,
real, dimension( : ), intent(in)  ONE_M_FTHETA_T2OLD,
real, dimension( : ), intent(in)  FTHETA_T2_J,
real, dimension( : ), intent(in)  ONE_M_FTHETA_T2OLD_J,
logical, intent(in)  integrate_other_side_and_not_boundary,
real, dimension( :, :, : ), intent(in)  loc_u,
real, dimension( : ), intent(in)  THETA_VEL,
real, dimension( mdims%ndim, final_phase ), intent(inout)  UDGI_IMP_ALL,
real, dimension( final_phase ), intent(inout)  RCON,
real, dimension( final_phase ), intent(inout)  RCON_J,
real, dimension( final_phase ), intent(inout)  NDOTQ_IMP,
real, dimension(:,:), pointer  X_ALL,
real, dimension(:,:,: ), pointer  SUF_D_BC_ALL,
real  gravty 
)

This subroutine caculates the discretised cty eqn acting on the velocities i.e. MmatCT, MmatCT_RHS It also computes the gradient matrix using the DCVFE method.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ scvdetnx()

subroutine cv_advection::scvdetnx ( type(multi_dimensions), intent(in)  Mdims,
type(multi_ndgln), intent(in)  ndgln,
real, dimension(:,:)  X_ALL,
type(multi_shape_funs), intent(in)  CV_funs,
type(multi_gi_dimensions), intent(in)  CV_GIdims,
logical, intent(in)  on_domain_boundary,
logical, intent(in)  between_elements,
integer, intent(in)  ELE,
integer, intent(in)  GI,
real, dimension( : ), intent(inout)  SCVDETWEI,
real, dimension( mdims%ndim, cv_gidims%scvngi ), intent(inout)  CVNORMX_ALL,
real, dimension( mdims%ndim ), intent(in)  XC_ALL,
integer, intent(in)  X_NOD,
integer, intent(in)  X_NODJ 
)

  • this subroutine computed the surface area at the Gi point (SCVDETWEI)
  • this subroutine calculates the control volume (CV)
  • CVNORMX, CVNORMY, CVNORMZ normals at the Gaussian
  • integration points GI. NODI = the current global
  • node number for the co-ordinates.
  • (XC,YC,ZC) is the centre of CV NODI

- date last modified : 25/08/2020

Here is the call graph for this function:
Here is the caller graph for this function:

◆ shock_front_in_ele()

logical function cv_advection::shock_front_in_ele ( integer  ele,
type(multi_dimensions), intent(in)  Mdims,
real, dimension(:,:), intent(in)  sat,
type(multi_ndgln), intent(in)  ndgln,
real, dimension(:), intent(in)  Imble_frac 
)

Detects whether the element has a shockfront or not.

Here is the caller graph for this function:

◆ triloccords()

subroutine cv_advection::triloccords ( real  Xp,
real  Yp,
real  Zp,
real  N1,
real  N2,
real  N3,
real  N4,
real  X1,
real  Y1,
real  Z1,
real  X2,
real  Y2,
real  Z2,
real  X3,
real  Y3,
real  Z3,
real  X4,
real  Y4,
real  Z4 
)
Here is the call graph for this function:
Here is the caller graph for this function:

◆ triloccords2d()

subroutine cv_advection::triloccords2d ( real  Xp,
real  Yp,
real  N1,
real  N2,
real  N3,
real  X1,
real  Y1,
real  X2,
real  Y2,
real  X3,
real  Y3 
)
Here is the call graph for this function:
Here is the caller graph for this function:

◆ unpack_loc()

subroutine cv_advection::unpack_loc ( real, dimension(:), intent(in)  LOC_F,
real, dimension(nphase), intent(inout)  T_ALL,
integer, intent(in)  NPHASE,
integer, intent(inout)  IPT,
logical, dimension(nphase), intent(in)  IGOT_T_PACK,
logical, dimension(nphase), intent(in)  IGOT_T_CONST,
real, dimension(nphase), intent(in)  IGOT_T_CONST_VALUE 
)

If PACK then UNpack loc_f into T_ALL as long at IGOT_T==1 and STORE and not already in storage.

Here is the caller graph for this function:

◆ unpack_loc_all()

subroutine cv_advection::unpack_loc_all ( real, dimension(:), intent(in)  LOC_F,
real, dimension(:), intent(inout)  field1,
real, dimension(:), intent(inout)  oldfield1,
real, dimension(:), intent(inout)  field2,
real, dimension(:), intent(inout)  oldfield2,
real, dimension(:), intent(inout)  field3,
real, dimension(:), intent(inout)  oldfield3,
logical, dimension(:,:), intent(in)  IGOT_T_PACK,
logical, dimension(:,:), intent(in)  IGOT_T_CONST,
real, dimension(:,:), intent(in)  IGOT_T_CONST_VALUE,
logical, intent(in)  use_volume_frac_T2,
integer, intent(in)  nfield 
)

If PACK then UNpack loc_f into T_ALL as long at IGOT_T==1 and STORE and not already in storage. Checks if the fields are constant or not, stored in IGOT_T_PACK, based on that introduces the LOC_F into the field or use a constant value. This is after the limiters have been applied.

Here is the call graph for this function:
Here is the caller graph for this function:

Variable Documentation

◆ wic_d_bc_dirichlet

integer, parameter cv_advection::wic_d_bc_dirichlet = 1

◆ wic_p_bc_dirichlet

integer, parameter cv_advection::wic_p_bc_dirichlet = 1

◆ wic_p_bc_free

integer, parameter cv_advection::wic_p_bc_free = 2

◆ wic_t_bc_diri_adv_and_robin

integer, parameter cv_advection::wic_t_bc_diri_adv_and_robin = 3

◆ wic_t_bc_dirichlet

integer, parameter cv_advection::wic_t_bc_dirichlet = 1

◆ wic_t_bc_robin

integer, parameter cv_advection::wic_t_bc_robin = 2

◆ wic_u_bc_diri_adv_and_robin

integer, parameter cv_advection::wic_u_bc_diri_adv_and_robin = 3

◆ wic_u_bc_dirichlet

integer, parameter cv_advection::wic_u_bc_dirichlet = 1

◆ wic_u_bc_dirichlet_inout

integer, parameter cv_advection::wic_u_bc_dirichlet_inout = 2

◆ wic_u_bc_robin

integer, parameter cv_advection::wic_u_bc_robin = 2