mus_mrtInit_module Module

This module provides the definitions of M and Minv for MRT advection relaxation scheme for all stencils.

The weighted MRT (D3Q27) is based on the following paper Abbas Fakhari, Diogo Bolster, Li-Shi Luo "A weighted multiple-relaxation-time lattice Boltzmann method for multiphase flows and its application to partial coalescence cascades" Journal of Computational Physics, 2017

The MRT (D3Q19) implementation here is taken from:\n J. Toelke, S. Freudiger, and M. Krafczyk, "An adaptive scheme using hierarchical grids for lattice Boltzmann multi-phase flow simulations," Comput. Fluids, vol. 35, pp. 820–830, 2006. \n


Uses

  • module~~mus_mrtinit_module~~UsesGraph module~mus_mrtinit_module mus_mrtInit_module module~env_module env_module module~mus_mrtinit_module->module~env_module module~tem_param_module tem_param_module module~mus_mrtinit_module->module~tem_param_module

Used by

  • module~~mus_mrtinit_module~~UsedByGraph module~mus_mrtinit_module mus_mrtInit_module module~mus_mrt_d3q19_module mus_mrt_d3q19_module module~mus_mrt_d3q19_module->module~mus_mrtinit_module module~mus_mrt_d3q27_module mus_mrt_d3q27_module module~mus_mrt_d3q27_module->module~mus_mrtinit_module module~mus_moments_module mus_moments_module module~mus_moments_module->module~mus_mrtinit_module module~mus_initfluid_module mus_initFluid_module module~mus_initfluid_module->module~mus_mrt_d3q19_module module~mus_initfluid_module->module~mus_mrt_d3q27_module module~mus_initfluidincomp_module mus_initFluidIncomp_module module~mus_initfluidincomp_module->module~mus_mrt_d3q19_module module~mus_initfluidincomp_module->module~mus_mrt_d3q27_module module~mus_scheme_module mus_scheme_module module~mus_scheme_module->module~mus_moments_module module~mus_derivedquantities_module2 mus_derivedQuantities_module2 module~mus_derivedquantities_module2->module~mus_moments_module module~mus_interpolate_verify_module mus_interpolate_verify_module module~mus_interpolate_verify_module->module~mus_moments_module module~mus_d3q27_module mus_d3q27_module module~mus_d3q27_module->module~mus_derivedquantities_module2 module~mus_bc_fluid_noneqexpol_module mus_bc_fluid_nonEqExpol_module module~mus_bc_fluid_noneqexpol_module->module~mus_derivedquantities_module2 module~mus_derquanincomp_module mus_derQuanIncomp_module module~mus_derquanincomp_module->module~mus_derivedquantities_module2 module~mus_flow_module mus_flow_module module~mus_flow_module->module~mus_initfluid_module module~mus_flow_module->module~mus_initfluidincomp_module module~mus_flow_module->module~mus_derivedquantities_module2 module~mus_flow_module->module~mus_interpolate_verify_module module~mus_d2q9_module mus_d2q9_module module~mus_d2q9_module->module~mus_derivedquantities_module2 program~mus_harvesting mus_harvesting program~mus_harvesting->module~mus_scheme_module module~mus_derquan_module mus_derQuan_module module~mus_derquan_module->module~mus_derivedquantities_module2 module~mus_tools_module mus_tools_module module~mus_tools_module->module~mus_scheme_module module~mus_ibm_module mus_IBM_module module~mus_ibm_module->module~mus_derivedquantities_module2 module~mus_hvs_config_module mus_hvs_config_module module~mus_hvs_config_module->module~mus_scheme_module module~mus_d3q19_module mus_d3q19_module module~mus_d3q19_module->module~mus_derivedquantities_module2 module~mus_smagorinsky_module mus_Smagorinsky_module module~mus_smagorinsky_module->module~mus_derivedquantities_module2 module~mus_dynloadbal_module mus_dynLoadBal_module module~mus_dynloadbal_module->module~mus_scheme_module module~mus_interpolate_average_module mus_interpolate_average_module module~mus_interpolate_average_module->module~mus_derivedquantities_module2 module~mus_interpolate_linear_module mus_interpolate_linear_module module~mus_interpolate_linear_module->module~mus_derivedquantities_module2 module~mus_derquanps_module mus_derQuanPS_module module~mus_derquanps_module->module~mus_derivedquantities_module2 module~mus_config_module mus_config_module module~mus_config_module->module~mus_scheme_module module~mus_derquanpoisson_module mus_derQuanPoisson_module module~mus_derquanpoisson_module->module~mus_derivedquantities_module2 module~mus_program_module mus_program_module module~mus_program_module->module~mus_scheme_module module~mus_interpolate_quadratic_module mus_interpolate_quadratic_module module~mus_interpolate_quadratic_module->module~mus_derivedquantities_module2 module~mus_derquanisothermaceq_module mus_derQuanIsothermAcEq_module module~mus_derquanisothermaceq_module->module~mus_derivedquantities_module2

Contents


Variables

TypeVisibilityAttributesNameInitial
real(kind=rk), public, parameter, dimension(19,19):: MMtrD3Q19 =reshape((/1._rk, 1._rk, 1._rk, 1._rk, 1._rk, 1._rk, 1._rk, 1._rk, 1._rk, 1._rk, 1._rk, 1._rk, 1._rk, 1._rk, 1._rk, 1._rk, 1._rk, 1._rk, 1._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 1._rk, 1._rk, 1._rk, 1._rk, 1._rk, 1._rk, 1._rk, 1._rk, 1._rk, 1._rk, 1._rk, 1._rk, -1._rk, -2._rk, -2._rk, -2._rk, -2._rk, -2._rk, -2._rk, 1._rk, 1._rk, 1._rk, 1._rk, 1._rk, 1._rk, 1._rk, 1._rk, 1._rk, 1._rk, 1._rk, 1._rk, 1._rk, -1._rk, 0._rk, 0._rk, 1._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, -1._rk, 1._rk, -1._rk, 1._rk, -1._rk, -1._rk, 1._rk, 1._rk, 0._rk, 2._rk, 0._rk, 0._rk, -2._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, -1._rk, 1._rk, -1._rk, 1._rk, -1._rk, -1._rk, 1._rk, 1._rk, 0._rk, 0._rk, -1._rk, 0._rk, 0._rk, 1._rk, 0._rk, -1._rk, -1._rk, 1._rk, 1._rk, 0._rk, 0._rk, 0._rk, 0._rk, -1._rk, 1._rk, -1._rk, 1._rk, 0._rk, 0._rk, 2._rk, 0._rk, 0._rk, -2._rk, 0._rk, -1._rk, -1._rk, 1._rk, 1._rk, 0._rk, 0._rk, 0._rk, 0._rk, -1._rk, 1._rk, -1._rk, 1._rk, 0._rk, 0._rk, 0._rk, -1._rk, 0._rk, 0._rk, 1._rk, -1._rk, 1._rk, -1._rk, 1._rk, -1._rk, -1._rk, 1._rk, 1._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 2._rk, 0._rk, 0._rk, -2._rk, -1._rk, 1._rk, -1._rk, 1._rk, -1._rk, -1._rk, 1._rk, 1._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 2._rk, -1._rk, -1._rk, 2._rk, -1._rk, -1._rk, -2._rk, -2._rk, -2._rk, -2._rk, 1._rk, 1._rk, 1._rk, 1._rk, 1._rk, 1._rk, 1._rk, 1._rk, 0._rk, -2._rk, 1._rk, 1._rk, -2._rk, 1._rk, 1._rk, -2._rk, -2._rk, -2._rk, -2._rk, 1._rk, 1._rk, 1._rk, 1._rk, 1._rk, 1._rk, 1._rk, 1._rk, 0._rk, 0._rk, 1._rk, -1._rk, 0._rk, 1._rk, -1._rk, 0._rk, 0._rk, 0._rk, 0._rk, -1._rk, -1._rk, -1._rk, -1._rk, 1._rk, 1._rk, 1._rk, 1._rk, 0._rk, 0._rk, -1._rk, 1._rk, 0._rk, -1._rk, 1._rk, 0._rk, 0._rk, 0._rk, 0._rk, -1._rk, -1._rk, -1._rk, -1._rk, 1._rk, 1._rk, 1._rk, 1._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 1._rk, -1._rk, -1._rk, 1._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 1._rk, -1._rk, -1._rk, 1._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 1._rk, -1._rk, -1._rk, 1._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 1._rk, -1._rk, 1._rk, -1._rk, -1._rk, -1._rk, 1._rk, 1._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, -1._rk, -1._rk, 1._rk, 1._rk, 0._rk, 0._rk, 0._rk, 0._rk, 1._rk, -1._rk, 1._rk, -1._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 1._rk, -1._rk, 1._rk, -1._rk, -1._rk, -1._rk, 1._rk, 1._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk/), (/19, 19/), order=(/2, 1/))
real(kind=rk), public, parameter, dimension(19,19):: MMivD3Q19 =reshape((/div1_18, 0._rk, -div1_18, -div1_6, div1_6, 0._rk, 0._rk, 0._rk, 0._rk, div1_12, -div1_12, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, div1_18, 0._rk, -div1_18, 0._rk, 0._rk, -div1_6, div1_6, 0._rk, 0._rk, -div1_24, div1_24, div1_8, -div1_8, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, div1_18, 0._rk, -div1_18, 0._rk, 0._rk, 0._rk, 0._rk, -div1_6, div1_6, -div1_24, div1_24, -div1_8, div1_8, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, div1_18, 0._rk, -div1_18, div1_6, -div1_6, 0._rk, 0._rk, 0._rk, 0._rk, div1_12, -div1_12, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, div1_18, 0._rk, -div1_18, 0._rk, 0._rk, div1_6, -div1_6, 0._rk, 0._rk, -div1_24, div1_24, div1_8, -div1_8, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, div1_18, 0._rk, -div1_18, 0._rk, 0._rk, 0._rk, 0._rk, div1_6, -div1_6, -div1_24, div1_24, -div1_8, div1_8, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, div1_36, div1_24, div1_72, 0._rk, 0._rk, -div1_12, -div1_24, -div1_12, -div1_24, -div1_24, -div1_24, 0._rk, 0._rk, 0._rk, div1_4, 0._rk, 0._rk, -div1_8, div1_8, div1_36, div1_24, div1_72, 0._rk, 0._rk, -div1_12, -div1_24, div1_12, div1_24, -div1_24, -div1_24, 0._rk, 0._rk, 0._rk, -div1_4, 0._rk, 0._rk, -div1_8, -div1_8, div1_36, div1_24, div1_72, 0._rk, 0._rk, div1_12, div1_24, -div1_12, -div1_24, -div1_24, -div1_24, 0._rk, 0._rk, 0._rk, -div1_4, 0._rk, 0._rk, div1_8, div1_8, div1_36, div1_24, div1_72, 0._rk, 0._rk, div1_12, div1_24, div1_12, div1_24, -div1_24, -div1_24, 0._rk, 0._rk, 0._rk, div1_4, 0._rk, 0._rk, div1_8, -div1_8, div1_36, div1_24, div1_72, -div1_12, -div1_24, 0._rk, 0._rk, -div1_12, -div1_24, div1_48, div1_48, -div1_16, -div1_16, 0._rk, 0._rk, div1_4, div1_8, 0._rk, -div1_8, div1_36, div1_24, div1_72, div1_12, div1_24, 0._rk, 0._rk, -div1_12, -div1_24, div1_48, div1_48, -div1_16, -div1_16, 0._rk, 0._rk, -div1_4, -div1_8, 0._rk, -div1_8, div1_36, div1_24, div1_72, -div1_12, -div1_24, 0._rk, 0._rk, div1_12, div1_24, div1_48, div1_48, -div1_16, -div1_16, 0._rk, 0._rk, -div1_4, div1_8, 0._rk, div1_8, div1_36, div1_24, div1_72, div1_12, div1_24, 0._rk, 0._rk, div1_12, div1_24, div1_48, div1_48, -div1_16, -div1_16, 0._rk, 0._rk, div1_4, -div1_8, 0._rk, div1_8, div1_36, div1_24, div1_72, -div1_12, -div1_24, -div1_12, -div1_24, 0._rk, 0._rk, div1_48, div1_48, div1_16, div1_16, div1_4, 0._rk, 0._rk, -div1_8, div1_8, 0._rk, div1_36, div1_24, div1_72, -div1_12, -div1_24, div1_12, div1_24, 0._rk, 0._rk, div1_48, div1_48, div1_16, div1_16, -div1_4, 0._rk, 0._rk, -div1_8, -div1_8, 0._rk, div1_36, div1_24, div1_72, div1_12, div1_24, -div1_12, -div1_24, 0._rk, 0._rk, div1_48, div1_48, div1_16, div1_16, -div1_4, 0._rk, 0._rk, div1_8, div1_8, 0._rk, div1_36, div1_24, div1_72, div1_12, div1_24, div1_12, div1_24, 0._rk, 0._rk, div1_48, div1_48, div1_16, div1_16, div1_4, 0._rk, 0._rk, div1_8, -div1_8, 0._rk, div1_3, -div1_2, div1_6, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk/), (/19, 19/), order=(/2, 1/))
real(kind=rk), public, parameter, dimension(27,27):: WMMtrD3Q27 =reshape((/1._rk, 1._rk, 1._rk, 1._rk, 1._rk, 1._rk, 1._rk, 1._rk, 1._rk, 1._rk, 1._rk, 1._rk, 1._rk, 1._rk, 1._rk, 1._rk, 1._rk, 1._rk, 1._rk, 1._rk, 1._rk, 1._rk, 1._rk, 1._rk, 1._rk, 1._rk, 1._rk, -1._rk, 0._rk, 0._rk, 1._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, -1._rk, 1._rk, -1._rk, 1._rk, -1._rk, -1._rk, 1._rk, 1._rk, -1._rk, -1._rk, -1._rk, -1._rk, 1._rk, 1._rk, 1._rk, 1._rk, 0._rk, 0._rk, -1._rk, 0._rk, 0._rk, 1._rk, 0._rk, -1._rk, -1._rk, 1._rk, 1._rk, 0._rk, 0._rk, 0._rk, 0._rk, -1._rk, 1._rk, -1._rk, 1._rk, -1._rk, -1._rk, 1._rk, 1._rk, -1._rk, -1._rk, 1._rk, 1._rk, 0._rk, 0._rk, 0._rk, -1._rk, 0._rk, 0._rk, 1._rk, -1._rk, 1._rk, -1._rk, 1._rk, -1._rk, -1._rk, 1._rk, 1._rk, 0._rk, 0._rk, 0._rk, 0._rk, -1._rk, 1._rk, -1._rk, 1._rk, -1._rk, 1._rk, -1._rk, 1._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 1._rk, -1._rk, -1._rk, 1._rk, 1._rk, 1._rk, -1._rk, -1._rk, -1._rk, -1._rk, 1._rk, 1._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 1._rk, -1._rk, -1._rk, 1._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 1._rk, -1._rk, -1._rk, 1._rk, 1._rk, -1._rk, -1._rk, 1._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 1._rk, -1._rk, -1._rk, 1._rk, 0._rk, 0._rk, 0._rk, 0._rk, 1._rk, -1._rk, 1._rk, -1._rk, -1._rk, 1._rk, -1._rk, 1._rk, 0._rk, 2._rk, -1._rk, -1._rk, 2._rk, -1._rk, -1._rk, -2._rk, -2._rk, -2._rk, -2._rk, 1._rk, 1._rk, 1._rk, 1._rk, 1._rk, 1._rk, 1._rk, 1._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 1._rk, -1._rk, 0._rk, 1._rk, -1._rk, 0._rk, 0._rk, 0._rk, 0._rk, -1._rk, -1._rk, -1._rk, -1._rk, 1._rk, 1._rk, 1._rk, 1._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 1._rk, 1._rk, 1._rk, 1._rk, 1._rk, 1._rk, 1._rk, 1._rk, 1._rk, 1._rk, 1._rk, 1._rk, 2._rk, 2._rk, 2._rk, 2._rk, 2._rk, 2._rk, 2._rk, 2._rk, -1._rk, 2._rk, 0._rk, 0._rk, -2._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, -1._rk, 1._rk, -1._rk, 1._rk, -1._rk, -1._rk, 1._rk, 1._rk, -4._rk, -4._rk, -4._rk, -4._rk, 4._rk, 4._rk, 4._rk, 4._rk, 0._rk, 0._rk, 2._rk, 0._rk, 0._rk, -2._rk, 0._rk, -1._rk, -1._rk, 1._rk, 1._rk, 0._rk, 0._rk, 0._rk, 0._rk, -1._rk, 1._rk, -1._rk, 1._rk, -4._rk, -4._rk, 4._rk, 4._rk, -4._rk, -4._rk, 4._rk, 4._rk, 0._rk, 0._rk, 0._rk, 2._rk, 0._rk, 0._rk, -2._rk, -1._rk, 1._rk, -1._rk, 1._rk, -1._rk, -1._rk, 1._rk, 1._rk, 0._rk, 0._rk, 0._rk, 0._rk, -4._rk, 4._rk, -4._rk, 4._rk, -4._rk, 4._rk, -4._rk, 4._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 1._rk, -1._rk, 1._rk, -1._rk, -1._rk, -1._rk, 1._rk, 1._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, -1._rk, -1._rk, 1._rk, 1._rk, 0._rk, 0._rk, 0._rk, 0._rk, 1._rk, -1._rk, 1._rk, -1._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 1._rk, -1._rk, 1._rk, -1._rk, -1._rk, -1._rk, 1._rk, 1._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, -1._rk, 1._rk, 1._rk, -1._rk, 1._rk, -1._rk, -1._rk, 1._rk, 0._rk, -1._rk, -1._rk, -1._rk, -1._rk, -1._rk, -1._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 4._rk, 4._rk, 4._rk, 4._rk, 4._rk, 4._rk, 4._rk, 4._rk, 1._rk, -2._rk, 1._rk, 1._rk, -2._rk, 1._rk, 1._rk, -4._rk, -4._rk, -4._rk, -4._rk, 2._rk, 2._rk, 2._rk, 2._rk, 2._rk, 2._rk, 2._rk, 2._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, -1._rk, 1._rk, 0._rk, -1._rk, 1._rk, 0._rk, 0._rk, 0._rk, 0._rk, -2._rk, -2._rk, -2._rk, -2._rk, 2._rk, 2._rk, 2._rk, 2._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, -1._rk, 1._rk, 1._rk, -1._rk, 2._rk, 2._rk, -2._rk, -2._rk, -2._rk, -2._rk, 2._rk, 2._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, -1._rk, 1._rk, 1._rk, -1._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 2._rk, -2._rk, -2._rk, 2._rk, 2._rk, -2._rk, -2._rk, 2._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, -1._rk, 1._rk, 1._rk, -1._rk, 0._rk, 0._rk, 0._rk, 0._rk, 2._rk, -2._rk, 2._rk, -2._rk, -2._rk, 2._rk, -2._rk, 2._rk, 0._rk, -1._rk, 0._rk, 0._rk, 1._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 2._rk, -2._rk, 2._rk, -2._rk, 2._rk, 2._rk, -2._rk, -2._rk, -4._rk, -4._rk, -4._rk, -4._rk, 4._rk, 4._rk, 4._rk, 4._rk, 0._rk, 0._rk, -1._rk, 0._rk, 0._rk, 1._rk, 0._rk, 2._rk, 2._rk, -2._rk, -2._rk, 0._rk, 0._rk, 0._rk, 0._rk, 2._rk, -2._rk, 2._rk, -2._rk, -4._rk, -4._rk, 4._rk, 4._rk, -4._rk, -4._rk, 4._rk, 4._rk, 0._rk, 0._rk, 0._rk, -1._rk, 0._rk, 0._rk, 1._rk, 2._rk, -2._rk, 2._rk, -2._rk, 2._rk, 2._rk, -2._rk, -2._rk, 0._rk, 0._rk, 0._rk, 0._rk, -4._rk, 4._rk, -4._rk, 4._rk, -4._rk, 4._rk, -4._rk, 4._rk, 0._rk, 2._rk, 2._rk, 2._rk, 2._rk, 2._rk, 2._rk, -4._rk, -4._rk, -4._rk, -4._rk, -4._rk, -4._rk, -4._rk, -4._rk, -4._rk, -4._rk, -4._rk, -4._rk, 8._rk, 8._rk, 8._rk, 8._rk, 8._rk, 8._rk, 8._rk, 8._rk, -1._rk/), (/27, 27/), order=(/2, 1/))
real(kind=rk), public, parameter, dimension(27,27):: WMMIvD3Q27 =reshape((/div2_27, -div2_9, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, div1_9, 0._rk, 0._rk, div1_9, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, -div1_18, -div1_18, 0._rk, 0._rk, 0._rk, 0._rk, -div1_18, 0._rk, 0._rk, div1_54, div2_27, 0._rk, -div2_9, 0._rk, 0._rk, 0._rk, 0._rk, -div1_18, div1_6, 0._rk, 0._rk, div1_9, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, -div1_18, div1_36, -div1_12, 0._rk, 0._rk, 0._rk, 0._rk, -div1_18, 0._rk, div1_54, div2_27, 0._rk, 0._rk, -div2_9, 0._rk, 0._rk, 0._rk, -div1_18, -div1_6, 0._rk, 0._rk, 0._rk, div1_9, 0._rk, 0._rk, 0._rk, 0._rk, -div1_18, div1_36, div1_12, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, -div1_18, div1_54, div2_27, div2_9, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, div1_9, 0._rk, 0._rk, -div1_9, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, -div1_18, -div1_18, 0._rk, 0._rk, 0._rk, 0._rk, div1_18, 0._rk, 0._rk, div1_54, div2_27, 0._rk, div2_9, 0._rk, 0._rk, 0._rk, 0._rk, -div1_18, div1_6, 0._rk, 0._rk, -div1_9, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, -div1_18, div1_36, -div1_12, 0._rk, 0._rk, 0._rk, 0._rk, div1_18, 0._rk, div1_54, div2_27, 0._rk, 0._rk, div2_9, 0._rk, 0._rk, 0._rk, -div1_18, -div1_6, 0._rk, 0._rk, 0._rk, -div1_9, 0._rk, 0._rk, 0._rk, 0._rk, -div1_18, div1_36, div1_12, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, div1_18, div1_54, div1_54, 0._rk, -div1_18, -div1_18, 0._rk, div1_6, 0._rk, -div1_36, 0._rk, div1_36, 0._rk, -div1_72, -div1_72, 0._rk, -div1_8, div1_8, 0._rk, 0._rk, -div1_36, 0._rk, 0._rk, -div1_12, 0._rk, 0._rk, div1_36, div1_36, -div1_108, div1_54, 0._rk, -div1_18, div1_18, 0._rk, -div1_6, 0._rk, -div1_36, 0._rk, div1_36, 0._rk, -div1_72, div1_72, 0._rk, -div1_8, -div1_8, 0._rk, 0._rk, -div1_36, 0._rk, 0._rk, div1_12, 0._rk, 0._rk, div1_36, -div1_36, -div1_108, div1_54, 0._rk, div1_18, -div1_18, 0._rk, -div1_6, 0._rk, -div1_36, 0._rk, div1_36, 0._rk, div1_72, -div1_72, 0._rk, div1_8, div1_8, 0._rk, 0._rk, -div1_36, 0._rk, 0._rk, div1_12, 0._rk, 0._rk, -div1_36, div1_36, -div1_108, div1_54, 0._rk, div1_18, div1_18, 0._rk, div1_6, 0._rk, -div1_36, 0._rk, div1_36, 0._rk, div1_72, div1_72, 0._rk, div1_8, -div1_8, 0._rk, 0._rk, -div1_36, 0._rk, 0._rk, -div1_12, 0._rk, 0._rk, -div1_36, -div1_36, -div1_108, div1_54, -div1_18, 0._rk, -div1_18, 0._rk, 0._rk, div1_6, div1_72, -div1_24, div1_36, -div1_72, 0._rk, -div1_72, div1_8, 0._rk, -div1_8, 0._rk, 0._rk, div1_72, -div1_24, 0._rk, 0._rk, -div1_12, div1_36, 0._rk, div1_36, -div1_108, div1_54, div1_18, 0._rk, -div1_18, 0._rk, 0._rk, -div1_6, div1_72, -div1_24, div1_36, div1_72, 0._rk, -div1_72, -div1_8, 0._rk, -div1_8, 0._rk, 0._rk, div1_72, -div1_24, 0._rk, 0._rk, div1_12, -div1_36, 0._rk, div1_36, -div1_108, div1_54, -div1_18, 0._rk, div1_18, 0._rk, 0._rk, -div1_6, div1_72, -div1_24, div1_36, -div1_72, 0._rk, div1_72, div1_8, 0._rk, div1_8, 0._rk, 0._rk, div1_72, -div1_24, 0._rk, 0._rk, div1_12, div1_36, 0._rk, -div1_36, -div1_108, div1_54, div1_18, 0._rk, div1_18, 0._rk, 0._rk, div1_6, div1_72, -div1_24, div1_36, div1_72, 0._rk, div1_72, -div1_8, 0._rk, div1_8, 0._rk, 0._rk, div1_72, -div1_24, 0._rk, 0._rk, -div1_12, -div1_36, 0._rk, -div1_36, -div1_108, div1_54, -div1_18, -div1_18, 0._rk, div1_6, 0._rk, 0._rk, div1_72, div1_24, div1_36, -div1_72, -div1_72, 0._rk, -div1_8, div1_8, 0._rk, 0._rk, 0._rk, div1_72, div1_24, -div1_12, 0._rk, 0._rk, div1_36, div1_36, 0._rk, -div1_108, div1_54, -div1_18, div1_18, 0._rk, -div1_6, 0._rk, 0._rk, div1_72, div1_24, div1_36, -div1_72, div1_72, 0._rk, -div1_8, -div1_8, 0._rk, 0._rk, 0._rk, div1_72, div1_24, div1_12, 0._rk, 0._rk, div1_36, -div1_36, 0._rk, -div1_108, div1_54, div1_18, -div1_18, 0._rk, -div1_6, 0._rk, 0._rk, div1_72, div1_24, div1_36, div1_72, -div1_72, 0._rk, div1_8, div1_8, 0._rk, 0._rk, 0._rk, div1_72, div1_24, div1_12, 0._rk, 0._rk, -div1_36, div1_36, 0._rk, -div1_108, div1_54, div1_18, div1_18, 0._rk, div1_6, 0._rk, 0._rk, div1_72, div1_24, div1_36, div1_72, div1_72, 0._rk, div1_8, -div1_8, 0._rk, 0._rk, 0._rk, div1_72, div1_24, -div1_12, 0._rk, 0._rk, -div1_36, -div1_36, 0._rk, -div1_108, div1_216, -div1_72, -div1_72, -div1_72, div1_24, div1_24, div1_24, 0._rk, 0._rk, div1_72, -div1_72, -div1_72, -div1_72, 0._rk, 0._rk, 0._rk, -div1_8, div1_72, 0._rk, 0._rk, div1_24, div1_24, div1_24, -div1_72, -div1_72, -div1_72, div1_216, div1_216, -div1_72, -div1_72, div1_72, div1_24, -div1_24, -div1_24, 0._rk, 0._rk, div1_72, -div1_72, -div1_72, div1_72, 0._rk, 0._rk, 0._rk, div1_8, div1_72, 0._rk, 0._rk, div1_24, -div1_24, -div1_24, -div1_72, -div1_72, div1_72, div1_216, div1_216, -div1_72, div1_72, -div1_72, -div1_24, -div1_24, div1_24, 0._rk, 0._rk, div1_72, -div1_72, div1_72, -div1_72, 0._rk, 0._rk, 0._rk, div1_8, div1_72, 0._rk, 0._rk, -div1_24, -div1_24, div1_24, -div1_72, div1_72, -div1_72, div1_216, div1_216, -div1_72, div1_72, div1_72, -div1_24, div1_24, -div1_24, 0._rk, 0._rk, div1_72, -div1_72, div1_72, div1_72, 0._rk, 0._rk, 0._rk, -div1_8, div1_72, 0._rk, 0._rk, -div1_24, div1_24, -div1_24, -div1_72, div1_72, div1_72, div1_216, div1_216, div1_72, -div1_72, -div1_72, -div1_24, div1_24, -div1_24, 0._rk, 0._rk, div1_72, div1_72, -div1_72, -div1_72, 0._rk, 0._rk, 0._rk, div1_8, div1_72, 0._rk, 0._rk, -div1_24, div1_24, -div1_24, div1_72, -div1_72, -div1_72, div1_216, div1_216, div1_72, -div1_72, div1_72, -div1_24, -div1_24, div1_24, 0._rk, 0._rk, div1_72, div1_72, -div1_72, div1_72, 0._rk, 0._rk, 0._rk, -div1_8, div1_72, 0._rk, 0._rk, -div1_24, -div1_24, div1_24, div1_72, -div1_72, div1_72, div1_216, div1_216, div1_72, div1_72, -div1_72, div1_24, -div1_24, -div1_24, 0._rk, 0._rk, div1_72, div1_72, div1_72, -div1_72, 0._rk, 0._rk, 0._rk, -div1_8, div1_72, 0._rk, 0._rk, div1_24, -div1_24, -div1_24, div1_72, div1_72, -div1_72, div1_216, div1_216, div1_72, div1_72, div1_72, div1_24, div1_24, div1_24, 0._rk, 0._rk, div1_72, div1_72, div1_72, div1_72, 0._rk, 0._rk, 0._rk, div1_8, div1_72, 0._rk, 0._rk, div1_24, div1_24, div1_24, div1_72, div1_72, div1_72, div1_216, div8_27, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, -div4_9, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, div2_9, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, 0._rk, -div1_27/), (/27, 27/), order=(/2, 1/))

Functions

public function check_mrt_matrix_d3q19() result(test)

Unoptimized explicit implementation

Read more…

Arguments

None

Return Value logical

public function check_mrt_matrix_d3q27() result(test)

Unoptimized explicit implementation

Read more…

Arguments

None

Return Value logical