atl_maxwell_flux_module Module

module that holds all routines to calculate the flux for hyperbolic Maxwell equations.


Uses

  • module~~atl_maxwell_flux_module~~UsesGraph module~atl_maxwell_flux_module atl_maxwell_flux_module module~ply_poly_project_module ply_poly_project_module module~atl_maxwell_flux_module->module~ply_poly_project_module module~env_module env_module module~atl_maxwell_flux_module->module~env_module module~ply_oversample_module ply_oversample_module module~atl_maxwell_flux_module->module~ply_oversample_module module~ply_poly_project_module->module~env_module module~ply_fxt_module ply_fxt_module module~ply_poly_project_module->module~ply_fxt_module module~tem_logging_module tem_logging_module module~ply_poly_project_module->module~tem_logging_module module~ply_prj_header_module ply_prj_header_module module~ply_poly_project_module->module~ply_prj_header_module module~ply_dof_module ply_dof_module module~ply_poly_project_module->module~ply_dof_module module~ply_legfpt_3d_module ply_legFpt_3D_module module~ply_poly_project_module->module~ply_legfpt_3d_module module~tem_aux_module tem_aux_module module~ply_poly_project_module->module~tem_aux_module module~ply_legfpt_module ply_legFpt_module module~ply_poly_project_module->module~ply_legfpt_module module~ply_nodes_header_module ply_nodes_header_module module~ply_poly_project_module->module~ply_nodes_header_module module~tem_tools_module tem_tools_module module~ply_poly_project_module->module~tem_tools_module module~ply_dynarray_project_module ply_dynarray_project_module module~ply_poly_project_module->module~ply_dynarray_project_module module~ply_l2p_module ply_l2p_module module~ply_poly_project_module->module~ply_l2p_module module~ply_nodes_module ply_nodes_module module~ply_poly_project_module->module~ply_nodes_module module~ply_legfpt_2d_module ply_legFpt_2D_module module~ply_poly_project_module->module~ply_legfpt_2d_module module~ply_oversample_module->module~ply_poly_project_module module~ply_oversample_module->module~env_module module~ply_oversample_module->module~ply_dof_module module~ply_fxt_module->module~env_module fxt_fwrap fxt_fwrap module~ply_fxt_module->fxt_fwrap module~ply_fxt_header_module ply_fxt_header_module module~ply_fxt_module->module~ply_fxt_header_module module~ply_prj_header_module->module~env_module module~ply_prj_header_module->module~tem_logging_module module~ply_prj_header_module->module~tem_aux_module module~ply_prj_header_module->module~tem_tools_module module~aotus_module aotus_module module~ply_prj_header_module->module~aotus_module module~aot_out_module aot_out_module module~ply_prj_header_module->module~aot_out_module module~ply_prj_header_module->module~ply_fxt_header_module fftw_wrap fftw_wrap module~ply_prj_header_module->fftw_wrap module~ply_fpt_header_module ply_fpt_header_module module~ply_prj_header_module->module~ply_fpt_header_module module~ply_l2p_header_module ply_l2p_header_module module~ply_prj_header_module->module~ply_l2p_header_module module~ply_dof_module->module~env_module module~ply_legfpt_3d_module->module~env_module module~ply_legfpt_3d_module->module~ply_legfpt_module module~ply_legfpt_3d_module->fftw_wrap iso_c_binding iso_c_binding module~ply_legfpt_3d_module->iso_c_binding module~ply_legfpt_module->module~env_module module~ply_polybaseexc_module ply_polyBaseExc_module module~ply_legfpt_module->module~ply_polybaseexc_module module~ply_legfpt_module->fftw_wrap module~tem_compileconf_module tem_compileconf_module module~ply_legfpt_module->module~tem_compileconf_module module~ply_legfpt_module->iso_c_binding module~ply_legfpt_module->module~ply_fpt_header_module module~ply_nodes_header_module->module~env_module module~ply_dynarray_project_module->module~env_module module~ply_dynarray_project_module->module~tem_logging_module module~ply_dynarray_project_module->module~ply_prj_header_module module~ply_dynarray_project_module->module~aotus_module module~ply_l2p_module->module~env_module module~ply_l2p_module->module~tem_logging_module module~ply_l2p_module->module~tem_aux_module module~ply_l2p_module->module~tem_compileconf_module module~ply_modg_basis_module ply_modg_basis_module module~ply_l2p_module->module~ply_modg_basis_module module~ply_space_integration_module ply_space_integration_module module~ply_l2p_module->module~ply_space_integration_module module~ply_lagrange_module ply_lagrange_module module~ply_l2p_module->module~ply_lagrange_module module~ply_nodeset_module ply_nodeset_module module~ply_l2p_module->module~ply_nodeset_module module~ply_l2p_module->module~ply_l2p_header_module module~ply_nodes_module->module~env_module module~ply_nodes_module->module~tem_aux_module module~ply_nodes_module->module~ply_nodes_header_module module~ply_nodes_module->module~aotus_module module~ply_nodes_module->fftw_wrap module~ply_nodes_module->module~ply_nodeset_module module~ply_legfpt_2d_module->module~env_module module~ply_legfpt_2d_module->module~ply_legfpt_module module~ply_legfpt_2d_module->fftw_wrap module~ply_legfpt_2d_module->iso_c_binding

Used by

  • module~~atl_maxwell_flux_module~~UsedByGraph module~atl_maxwell_flux_module atl_maxwell_flux_module module~atl_modg_maxwelldivcor_kernel_module atl_modg_maxwellDivCor_kernel_module module~atl_modg_maxwelldivcor_kernel_module->module~atl_maxwell_flux_module module~atl_modg_maxwell_kernel_module atl_modg_maxwell_kernel_module module~atl_modg_maxwell_kernel_module->module~atl_maxwell_flux_module proc~compute_rhs_cubes_modg compute_rhs_cubes_modg proc~compute_rhs_cubes_modg->module~atl_modg_maxwelldivcor_kernel_module proc~compute_rhs_cubes_modg->module~atl_modg_maxwell_kernel_module

Contents


Interfaces

public interface atl_maxwell_flux

Interface for fluxes of pure Maxwell equations.

  • private subroutine maxwell_flux_cube(left, right, left_mu, left_epsi, right_mu, right_epsi, flux)

    Subroutine to calculate the flux for pure Maxwell equations without any divergence cleaning on the reference cubic face.

    Read more…

    Arguments

    TypeIntentOptionalAttributesName
    real(kind=rk), intent(in) :: left(6)

    Left state vector (as conservative variables). The order of this vector has to be \f$ (D_x, D_y, D_z, B_1, B_2, B_3) \f$ where E and B denoted electric field vetor and magnetic field (also called magnetic induction) vector.

    real(kind=rk), intent(in) :: right(6)

    Right state vector (as conservative variables). The order of this vector has to be (D_x, D_y, D_z, B_1, B_2, B_3) where E and B denoted the electric field vetor and magnetic field (also called magnetic induction) vector.

    real(kind=rk), intent(in) :: left_mu

    The magnetic permeability of the left element.

    real(kind=rk), intent(in) :: left_epsi

    The electric permitivity of the left element.

    real(kind=rk), intent(in) :: right_mu

    The magnetic permeability of the right element.

    real(kind=rk), intent(in) :: right_epsi

    The electric permitivity of the right element.

    real(kind=rk), intent(out) :: flux(6)

    The flux between left and right cell. The order of this vector is the same as the input arguments.

  • private subroutine maxwell_flux_cube_vec(nTotalFaces, nSides, nFaceDofs, faceRep, faceFlux, leftPos, rightPos, var, material_left, material_right)

    calculate flux of pure maxwell equation directly on the

    Read more…

    Arguments

    TypeIntentOptionalAttributesName
    integer, intent(in) :: nTotalFaces
    integer, intent(in) :: nSides
    integer, intent(in) :: nFaceDofs
    real(kind=rk), intent(in) :: faceRep(nTotalFaces,nFaceDofs,6,2)
    real(kind=rk), intent(inout) :: faceFlux(nTotalFaces,nFaceDofs,6,2)
    integer, intent(in) :: leftPos(nSides)
    integer, intent(in) :: rightPos(nsides)
    integer, intent(in) :: var(6)
    real(kind=rk), intent(in) :: material_left(nSides,1,2)
    real(kind=rk), intent(in) :: material_right(nSides,1,2)
  • private subroutine maxwell_flux_nonconst_cube_vec(nTotalFaces, nSides, nFaceDofs, faceRep, faceFlux, leftPos, rightPos, var, material_left, material_right, poly_proj, modalCoeffs, pntVal, nodalNumFlux, numFluxBuffer)

    calculate flux of pure maxwell equation directly on the

    Read more…

    Arguments

    TypeIntentOptionalAttributesName
    integer, intent(in) :: nTotalFaces
    integer, intent(in) :: nSides
    integer, intent(in) :: nFaceDofs
    real(kind=rk), intent(in) :: faceRep(nTotalFaces,nFaceDofs,6,2)

    The modal representation on the faces, left and right trace.

    real(kind=rk), intent(inout) :: faceFlux(nTotalFaces,nFaceDofs,6,2)

    The fluxes for all faces, for left and right elements.

    integer, intent(in) :: leftPos(nSides)

    Positions for the left and right elements of all faces

    integer, intent(in) :: rightPos(nsides)

    Positions for the left and right elements of all faces

    integer, intent(in) :: var(6)

    Variable rotation indices

    real(kind=rk), intent(in) :: material_left(nSides,nFaceDofs,2)

    Material parameters for the left faces.

    real(kind=rk), intent(in) :: material_right(nSides,nFaceDofs,2)

    Material parameters for the right faces.

    FPT to convert between nodes and modes.

    type(ply_poly_project_type) :: poly_proj

    Data for projection method

    real(kind=rk), intent(inout) :: modalCoeffs(:,:,:)

    Working array for the left and right modal coefficients

    real(kind=rk), intent(inout) :: pntVal(:,:,:)

    Working array for the left and right point values

    real(kind=rk), intent(inout) :: nodalNumFlux(:,:)

    Working array for the nodal flux

    real(kind=rk), intent(inout) :: numFluxBuffer(:,:)

    Working array for the modal numerical flux

public interface atl_maxwell_hc_flux

Interface for fluxes of Maxwell equations with hyperbolic divergence cleaning.

  • private subroutine maxwell_hc_flux_cube(left, right, mat_left, mat_right, flux)

    Subroutine to calculate the flux for pure Maxwell equations with

    Read more…

    Arguments

    TypeIntentOptionalAttributesName
    real(kind=rk), intent(in) :: left(8)

    Left state vector (as conservative variables). The order of this vector has to be \f$ (D_x, D_y, D_z, B_1, B_2, B_3, phi, psi) \f$ where E and B denoted electric field vetor and magnetic field (also called magnetic induction) vector.

    real(kind=rk), intent(in) :: right(8)

    Right state vector (as conservative variables). The order of this vector has to be (D_x, D_y, D_z, B_1, B_2, B_3, phi, psi) where E and B denoted the electric field vetor and magnetic field (also called magnetic induction) vector.

    real(kind=rk), intent(in) :: mat_left(4)

    Material for the left face

    real(kind=rk), intent(in) :: mat_right(4)

    Material for the right face !> The magnetic permeability of the left element. real(kind=rk), intent(in) :: left_mu !> The electric permitivity of the left element. real(kind=rk), intent(in) :: left_epsi !> Parameter for the magnetic correction on the left element. real(kind=rk), intent(in) :: left_gam !> Parameter for the electric correction on the left element. real(kind=rk), intent(in) :: left_chi !> The magnetic permeability of the right element. real(kind=rk), intent(in) :: right_mu !> The electric permitivity of the right element. real(kind=rk), intent(in) :: right_epsi !> Parameter for the magnetic correction on the right element. real(kind=rk), intent(in) :: right_gam !> Parameter for the electric correction on the right element. real(kind=rk), intent(in) :: right_chi

    real(kind=rk), intent(out) :: flux(8)

    The flux between left and right cell. The order of this vector is the same as the input arguments.

    JZ: Old implementation of the flux. There must be an error somewhere. I have to check my solution for the Riemann problem agian. So, I replaced the flux by a simple Lax-Friedrich type flux.

    real(kind=rk) :: left_speedOfLight, right_speedOfLight real(kind=rk) :: inv_denom_mu, inv_denom_epsi ! --------------------------------------------------------------------------

    ! The speed of light in the left and right element left_speedOfLight = 1.0_rk / sqrt( left_mu * left_epsi ) right_speedOfLight = 1.0_rk / sqrt( right_mu * right_epsi )

    ! The inverse of the denominators inv_denom_mu = 1.0_rk / ((-1.0_rk)left_speedOfLightleft_mu - right_speedOfLightright_mu) inv_denom_epsi = 1.0_rk / (left_speedOfLightleft_epsi + right_speedOfLight*right_epsi)

    ! D_x flux(1) = inv_denom_mu * ( & & (-1.0_rk)left_chileft(1)/left_epsi & & - (-1.0_rk)right_chiright(1)/right_epsi & & ) & & + ( & & left_chileft_chileft(7) & & + right_chiright_chiright(7) & & ) / ( left_epsileft_mu + right_epsiright_mu ) ! B_x flux(4) = inv_denom_epsi * ( & & (-1.0_rk)left_gamleft(4)/left_mu & & - (-1.0_rk)right_gamright(4)/right_mu) & & + ( & & left_gamleft_gamleft(8) & & + right_gamright_gamright(8) & & ) / ( left_epsileft_mu + right_epsiright_mu )

    ! the flux for phi (electric correction) flux(7) = ( & & left_speedOfLightleft(1) & & + right_speedOfLightright(1) & & + left_speedOfLightleft_chileft(7) & & - right_speedOfLightright_chiright(7) & & ) / (left_speedOfLight + right_speedOfLight)

    ! the flux for psi (magnetic correction) flux(8) = ( & & left_speedOfLightleft(4) & & + right_speedOfLightright(4) & & + left_speedOfLightleft_gamleft(8) & & - right_speedOfLightright_gamright(8) & & ) / (left_speedOfLight + right_speedOfLight)

    ! the flux for D_y flux(2) = ( & & ( (-1.0_rkleft(2) / left_epsi) & & - (-1.0_rkright(2) / right_epsi) ) & & - ( left_speedOfLight * left(6) & & + right_speedOfLight * right(6) )) ! the flux for B_z flux(6) = ( & & ( left_speedOfLight * left(2) & & + right_speedOfLight * right(2) ) & & + ( ( left(6) / left_mu ) & & - ( right(6) / right_mu ) )) ! the flux for D_z flux(3) = ( & & ( ( -1.0_rk * left(3) / left_epsi ) & & - ( -1.0_rk * right(3) / right_epsi ) ) & & + ( left_speedOfLight * left(5) & & + right_speedOfLight * right(5) ) & & )

    ! the flux for B_y flux(5) = ( & & ( -1.0_rk * left_speedOfLight * left(3) & & - right_speedOfLight * right(3) ) & & + ( ( left(5) / left_mu ) & & - ( right(5) / right_mu) ) & & )

    ! Normalize the calculated fluxes
    flux(2:3) = inv_denom_mu * flux(2:3)
    flux(5:6) = inv_denom_epsi * flux(5:6)
    
  • private subroutine maxwell_hc_flux_cube_vec(nTotalFaces, nSides, nFaceDofs, faceRep, faceFlux, leftPos, rightPos, var, material_left, material_right)

    calculate flux of maxwell equation with hyperbolic divergence

    Read more…

    Arguments

    TypeIntentOptionalAttributesName
    integer, intent(in) :: nTotalFaces
    integer, intent(in) :: nSides
    integer, intent(in) :: nFaceDofs
    real(kind=rk), intent(in) :: faceRep(nTotalFaces,nFaceDofs,8,2)
    real(kind=rk), intent(inout) :: faceFlux(nTotalFaces,nFaceDofs,8,2)
    integer, intent(in) :: leftPos(nSides)
    integer, intent(in) :: rightPos(nsides)
    integer, intent(in) :: var(8)
    real(kind=rk), intent(in) :: material_left(nSides,1,4)
    real(kind=rk), intent(in) :: material_right(nSides,1,4)

    integer :: iSide, left, right, iDof real(kind=rk) :: left_mu, right_mu real(kind=rk) :: left_epsi, right_epsi real(kind=rk) :: left_gam, right_gam real(kind=rk) :: left_chi, right_chi real(kind=rk) :: left_speedOfLight, right_speedOfLight real(kind=rk) :: inv_denom_mu, inv_denom_epsi ! --------------------------------------------------------------------------

    ! flux for D_X, B_Z, D_y and B_z do iDof = 1, nFaceDofs

    do iSide = 1, nSides

     ! The position of the left and right element in the state
     ! vector.
     left = leftPos(iSide)
     right = rightPos(iSide)
    
     ! The material parameters of the left and right element
     left_mu = material_left(iSide,1,1)
     left_epsi = material_left(iSide,1,2)
     left_gam = material_left(iSide,1,3)
     left_chi = material_left(iSide,1,4)
     left_speedOfLight = 1.0_rk / sqrt( left_mu * left_epsi )
     right_mu = material_right(iSide,1,1)
     right_epsi = material_right(iSide,1,2)
     right_gam = material_right(iSide,1,3)
     right_chi = material_right(iSide,1,4)
     right_speedOfLight = 1.0_rk / sqrt( right_mu * right_epsi )
    
     ! The inverse of the denominators
     inv_denom_mu = 1.0_rk / ((-1.0_rk)*left_speedOfLight*left_mu - right_speedOfLight*right_mu)
     inv_denom_epsi = 1.0_rk / (left_speedOfLight*left_epsi + right_speedOfLight*right_epsi)
    
     ! flux for D_x
     faceFlux(left,iDof,var(1),2) = inv_denom_mu * (                                   &
                                    &       (-1.0_rk)*left_chi*faceRep(left,iDof,var(1),2)/left_epsi    &
                                    &     - (-1.0_rk)*right_chi*faceRep(right,iDof,var(1),1)/right_epsi &
                                    &              )                                   &
                                    & + (                                              &
                                    &        left_chi*left_chi*faceRep(left,iDof,var(7),2) &
                                    &      + right_chi*right_chi*faceRep(right,iDof,var(7),1) &
                                    &   ) / ( left_epsi*left_mu + right_epsi*right_mu  )
    
     ! flux for B_x
     faceFlux(left,iDof,var(4),2) = inv_denom_epsi * (                                 &
                                    &        (-1.0_rk)*left_gam*faceRep(left,iDof,var(4),2)/left_mu &
                                    &      - (-1.0_rk)*right_gam*faceRep(right,iDof,var(4),1)/right_mu) &
                                    & + (                                              &
                                    &       left_gam*left_gam*faceRep(left,iDof,var(8),2) &
                                    &     + right_gam*right_gam*faceRep(right,iDof,var(8),1) &
                                    &   ) / ( left_epsi*left_mu + right_epsi*right_mu  )
    
     ! flux for phi (electric correction)
     faceFlux(left,iDof,var(7),2) = (                                            &
                                    &    left_speedOfLight*faceRep(left,iDof,var(1),2) &
                                    &  + right_speedOfLight*faceRep(right,iDof,var(1),1) &
                                    &  + left_speedOfLight*left_chi*faceRep(left,iDof,var(7),2) &
                                    &  - right_speedOfLight*right_chi*faceRep(right,iDof,var(7),1)   &
                                    & ) / (left_speedOfLight + right_speedOfLight)
    
     ! flux for psi (magnetic correction)
     faceFlux(left,iDof,var(8),2) = (                                            &
                                    &    left_speedOfLight*faceRep(left,iDof,var(4),2)               &
                                    &  + right_speedOfLight*faceRep(right,iDof,var(4),1)             &
                                    &  + left_speedOfLight*left_gam*faceRep(left,iDof,var(8),2)      &
                                    &  - right_speedOfLight*right_gam*faceRep(right,iDof,var(8),1)   &
                                    & ) / (left_speedOfLight + right_speedOfLight)
    
     ! the flux for D_y
     faceFlux(left,iDof,var(2),2) = ( &
       &                  ( ((-1.0_rk)*faceRep(left,iDof,var(2),2) / left_epsi)        &
       &                    - ((-1.0_rk)*faceRep(right,iDof,var(2),1) / right_epsi)  ) &
       &               -  ( left_speedOfLight * faceRep(left,iDof,var(6),2)          &
       &                    + right_speedOfLight * faceRep(right,iDof,var(6),1) ))
    
    
     ! the flux for B_z
     faceFlux(left,iDof,var(6),2) = ( &
       &                  ( left_speedOfLight * faceRep(left,iDof,var(2),2)     &
       &                  + right_speedOfLight * faceRep(right,iDof,var(2),1) )&
       &                + ( ( faceRep(left,iDof,var(6),2) / left_mu )     &
       &                  - ( faceRep(right,iDof,var(6),1) / right_mu ) ))
    
     ! the flux for D_z
     faceFlux(left,iDof,var(3),2) = ( &
       &                  ( ( (-1.0_rk) * faceRep(left,iDof,var(3),2) / left_epsi )     &
       &                  - ( (-1.0_rk) * faceRep(right,iDof,var(3),1) / right_epsi ) ) &
       &                + ( left_speedOfLight * faceRep(left,iDof,var(5),2)           &
       &                  + right_speedOfLight * faceRep(right,iDof,var(5),1) )       &
       &                            )
    
     ! the flux for B_y
     faceFlux(left,iDof,var(5),2) = ( &
       &                  ( (-1.0_rk) * left_speedOfLight * faceRep(left,iDof,var(3),2) &
       &                  - right_speedOfLight * faceRep(right,iDof,var(3),1) )       &
       &                + ( ( faceRep(left,iDof,var(5),2) / left_mu )                 &
       &                  - ( faceRep(right,iDof,var(5),1) / right_mu) )              &
       &                            )
    
     ! Normalize the calculated fluxes
     faceFlux(left,iDof,var(2),2) = inv_denom_mu * faceFlux(left,iDof,var(2),2)
     faceFlux(left,iDof,var(6),2) = inv_denom_epsi * faceFlux(left,iDof,var(6),2)
     faceFlux(left,iDof,var(3),2) = inv_denom_mu * faceFlux(left,iDof,var(3),2)
     faceFlux(left,iDof,var(5),2) = inv_denom_epsi * faceFlux(left,iDof,var(5),2)
    
     ! Assign the same flux for both adjacent elements
     faceFlux(right,iDof,:,1) = faceFlux(left,iDof,:,2)
    

    end do

    end do

  • private subroutine maxwell_hc_flux_nonconst_cube_vec(nTotalFaces, nSides, nFaceDofs, faceRep, faceFlux, leftPos, rightPos, var, material_left, material_right, poly_proj, left_modalCoeffs, right_modalCoeffs, left_pntVal, right_pntVal, nodalNumFlux, numFluxBuffer)

    calculate flux of maxwell equation with hyperbolic divergence

    Read more…

    Arguments

    TypeIntentOptionalAttributesName
    integer, intent(in) :: nTotalFaces
    integer, intent(in) :: nSides
    integer, intent(in) :: nFaceDofs
    real(kind=rk), intent(in) :: faceRep(nTotalFaces,nFaceDofs,8,2)

    The modal representation on the faces, left and right trace.

    real(kind=rk), intent(inout) :: faceFlux(nTotalFaces,nFaceDofs,8,2)

    The fluxes for all faces, for left and right elements.

    integer, intent(in) :: leftPos(nSides)

    Positions for the left and right elements of all faces

    integer, intent(in) :: rightPos(nsides)

    Positions for the left and right elements of all faces

    integer, intent(in) :: var(8)

    Variable rotation indices

    real(kind=rk), intent(in) :: material_left(nSides,nFaceDofs,4)

    Material parameters for the left faces.

    real(kind=rk), intent(in) :: material_right(nSides,nFaceDofs,4)

    Material parameters for the right faces.

    type(ply_poly_project_type) :: poly_proj

    Data for projection method !> Working array for the left and right modal coefficients real(kind=rk), intent(inout) :: left_modalCoeffs((fpt%nQuadPoints)2,8) real(kind=rk), intent(inout) :: right_modalCoeffs((fpt%nQuadPoints)2,8) !> Working array for the left and right point values real(kind=rk), intent(inout) :: left_pntVal((fpt%nQuadPoints)2,8) real(kind=rk), intent(inout) :: right_pntVal((fpt%nQuadPoints)2,8) !> Working array for the nodal flux real(kind=rk), intent(inout) :: nodalNumFlux((fpt%nQuadPoints)2,8) !> Working array for the modal numerical flux real(kind=rk), intent(inout) :: numFluxBuffer((fpt%nQuadPoints)2,8)

    real(kind=rk), intent(inout), allocatable:: left_modalCoeffs(:,:)

    Working array for the left and right modal coefficients

    real(kind=rk), intent(inout), allocatable:: right_modalCoeffs(:,:)
    real(kind=rk), intent(inout), allocatable:: left_pntVal(:,:)

    Working array for the left and right point values

    real(kind=rk), intent(inout), allocatable:: right_pntVal(:,:)
    real(kind=rk), intent(inout), allocatable:: nodalNumFlux(:,:)

    Working array for the nodal flux

    real(kind=rk), intent(inout), allocatable:: numFluxBuffer(:,:)

    Working array for the modal numerical flux


Functions

public function atl_physFluxMaxwellDivCor(state, material) result(flux)

Function for physical flux of the Maxwell equations in terms of D and B.

Arguments

TypeIntentOptionalAttributesName
real(kind=rk), intent(in) :: state(8)

State to compute the fluxes from (D,B).

real(kind=rk), intent(in) :: material(4)

Material parameters (mu, epsilon) the flux calculation

Return Value real(kind=rk)(8)

The resulting flux in x direction


Subroutines

private subroutine maxwell_flux_cube(left, right, left_mu, left_epsi, right_mu, right_epsi, flux)

Subroutine to calculate the flux for pure Maxwell equations without any divergence cleaning on the reference cubic face.

Read more…

Arguments

TypeIntentOptionalAttributesName
real(kind=rk), intent(in) :: left(6)

Left state vector (as conservative variables). The order of this vector has to be \f$ (D_x, D_y, D_z, B_1, B_2, B_3) \f$ where E and B denoted electric field vetor and magnetic field (also called magnetic induction) vector.

real(kind=rk), intent(in) :: right(6)

Right state vector (as conservative variables). The order of this vector has to be (D_x, D_y, D_z, B_1, B_2, B_3) where E and B denoted the electric field vetor and magnetic field (also called magnetic induction) vector.

real(kind=rk), intent(in) :: left_mu

The magnetic permeability of the left element.

real(kind=rk), intent(in) :: left_epsi

The electric permitivity of the left element.

real(kind=rk), intent(in) :: right_mu

The magnetic permeability of the right element.

real(kind=rk), intent(in) :: right_epsi

The electric permitivity of the right element.

real(kind=rk), intent(out) :: flux(6)

The flux between left and right cell. The order of this vector is the same as the input arguments.

private subroutine maxwell_flux_nonconst_cube_vec(nTotalFaces, nSides, nFaceDofs, faceRep, faceFlux, leftPos, rightPos, var, material_left, material_right, poly_proj, modalCoeffs, pntVal, nodalNumFlux, numFluxBuffer)

calculate flux of pure maxwell equation directly on the

Read more…

Arguments

TypeIntentOptionalAttributesName
integer, intent(in) :: nTotalFaces
integer, intent(in) :: nSides
integer, intent(in) :: nFaceDofs
real(kind=rk), intent(in) :: faceRep(nTotalFaces,nFaceDofs,6,2)

The modal representation on the faces, left and right trace.

real(kind=rk), intent(inout) :: faceFlux(nTotalFaces,nFaceDofs,6,2)

The fluxes for all faces, for left and right elements.

integer, intent(in) :: leftPos(nSides)

Positions for the left and right elements of all faces

integer, intent(in) :: rightPos(nsides)

Positions for the left and right elements of all faces

integer, intent(in) :: var(6)

Variable rotation indices

real(kind=rk), intent(in) :: material_left(nSides,nFaceDofs,2)

Material parameters for the left faces.

real(kind=rk), intent(in) :: material_right(nSides,nFaceDofs,2)

Material parameters for the right faces.

FPT to convert between nodes and modes.

type(ply_poly_project_type) :: poly_proj

Data for projection method

real(kind=rk), intent(inout) :: modalCoeffs(:,:,:)

Working array for the left and right modal coefficients

real(kind=rk), intent(inout) :: pntVal(:,:,:)

Working array for the left and right point values

real(kind=rk), intent(inout) :: nodalNumFlux(:,:)

Working array for the nodal flux

real(kind=rk), intent(inout) :: numFluxBuffer(:,:)

Working array for the modal numerical flux

private subroutine maxwell_flux_cube_vec(nTotalFaces, nSides, nFaceDofs, faceRep, faceFlux, leftPos, rightPos, var, material_left, material_right)

calculate flux of pure maxwell equation directly on the

Read more…

Arguments

TypeIntentOptionalAttributesName
integer, intent(in) :: nTotalFaces
integer, intent(in) :: nSides
integer, intent(in) :: nFaceDofs
real(kind=rk), intent(in) :: faceRep(nTotalFaces,nFaceDofs,6,2)
real(kind=rk), intent(inout) :: faceFlux(nTotalFaces,nFaceDofs,6,2)
integer, intent(in) :: leftPos(nSides)
integer, intent(in) :: rightPos(nsides)
integer, intent(in) :: var(6)
real(kind=rk), intent(in) :: material_left(nSides,1,2)
real(kind=rk), intent(in) :: material_right(nSides,1,2)

private subroutine maxwell_hc_flux_cube(left, right, mat_left, mat_right, flux)

Subroutine to calculate the flux for pure Maxwell equations with

Read more…

Arguments

TypeIntentOptionalAttributesName
real(kind=rk), intent(in) :: left(8)

Left state vector (as conservative variables). The order of this vector has to be \f$ (D_x, D_y, D_z, B_1, B_2, B_3, phi, psi) \f$ where E and B denoted electric field vetor and magnetic field (also called magnetic induction) vector.

real(kind=rk), intent(in) :: right(8)

Right state vector (as conservative variables). The order of this vector has to be (D_x, D_y, D_z, B_1, B_2, B_3, phi, psi) where E and B denoted the electric field vetor and magnetic field (also called magnetic induction) vector.

real(kind=rk), intent(in) :: mat_left(4)

Material for the left face

real(kind=rk), intent(in) :: mat_right(4)

Material for the right face !> The magnetic permeability of the left element. real(kind=rk), intent(in) :: left_mu !> The electric permitivity of the left element. real(kind=rk), intent(in) :: left_epsi !> Parameter for the magnetic correction on the left element. real(kind=rk), intent(in) :: left_gam !> Parameter for the electric correction on the left element. real(kind=rk), intent(in) :: left_chi !> The magnetic permeability of the right element. real(kind=rk), intent(in) :: right_mu !> The electric permitivity of the right element. real(kind=rk), intent(in) :: right_epsi !> Parameter for the magnetic correction on the right element. real(kind=rk), intent(in) :: right_gam !> Parameter for the electric correction on the right element. real(kind=rk), intent(in) :: right_chi

real(kind=rk), intent(out) :: flux(8)

The flux between left and right cell. The order of this vector is the same as the input arguments.

JZ: Old implementation of the flux. There must be an error somewhere. I have to check my solution for the Riemann problem agian. So, I replaced the flux by a simple Lax-Friedrich type flux.

real(kind=rk) :: left_speedOfLight, right_speedOfLight real(kind=rk) :: inv_denom_mu, inv_denom_epsi ! --------------------------------------------------------------------------

! The speed of light in the left and right element left_speedOfLight = 1.0_rk / sqrt( left_mu * left_epsi ) right_speedOfLight = 1.0_rk / sqrt( right_mu * right_epsi )

! The inverse of the denominators inv_denom_mu = 1.0_rk / ((-1.0_rk)left_speedOfLightleft_mu - right_speedOfLightright_mu) inv_denom_epsi = 1.0_rk / (left_speedOfLightleft_epsi + right_speedOfLight*right_epsi)

! D_x flux(1) = inv_denom_mu * ( & & (-1.0_rk)left_chileft(1)/left_epsi & & - (-1.0_rk)right_chiright(1)/right_epsi & & ) & & + ( & & left_chileft_chileft(7) & & + right_chiright_chiright(7) & & ) / ( left_epsileft_mu + right_epsiright_mu ) ! B_x flux(4) = inv_denom_epsi * ( & & (-1.0_rk)left_gamleft(4)/left_mu & & - (-1.0_rk)right_gamright(4)/right_mu) & & + ( & & left_gamleft_gamleft(8) & & + right_gamright_gamright(8) & & ) / ( left_epsileft_mu + right_epsiright_mu )

! the flux for phi (electric correction) flux(7) = ( & & left_speedOfLightleft(1) & & + right_speedOfLightright(1) & & + left_speedOfLightleft_chileft(7) & & - right_speedOfLightright_chiright(7) & & ) / (left_speedOfLight + right_speedOfLight)

! the flux for psi (magnetic correction) flux(8) = ( & & left_speedOfLightleft(4) & & + right_speedOfLightright(4) & & + left_speedOfLightleft_gamleft(8) & & - right_speedOfLightright_gamright(8) & & ) / (left_speedOfLight + right_speedOfLight)

! the flux for D_y flux(2) = ( & & ( (-1.0_rkleft(2) / left_epsi) & & - (-1.0_rkright(2) / right_epsi) ) & & - ( left_speedOfLight * left(6) & & + right_speedOfLight * right(6) )) ! the flux for B_z flux(6) = ( & & ( left_speedOfLight * left(2) & & + right_speedOfLight * right(2) ) & & + ( ( left(6) / left_mu ) & & - ( right(6) / right_mu ) )) ! the flux for D_z flux(3) = ( & & ( ( -1.0_rk * left(3) / left_epsi ) & & - ( -1.0_rk * right(3) / right_epsi ) ) & & + ( left_speedOfLight * left(5) & & + right_speedOfLight * right(5) ) & & )

! the flux for B_y flux(5) = ( & & ( -1.0_rk * left_speedOfLight * left(3) & & - right_speedOfLight * right(3) ) & & + ( ( left(5) / left_mu ) & & - ( right(5) / right_mu) ) & & )

! Normalize the calculated fluxes
flux(2:3) = inv_denom_mu * flux(2:3)
flux(5:6) = inv_denom_epsi * flux(5:6)

private subroutine maxwell_hc_flux_cube_vec(nTotalFaces, nSides, nFaceDofs, faceRep, faceFlux, leftPos, rightPos, var, material_left, material_right)

calculate flux of maxwell equation with hyperbolic divergence

Read more…

Arguments

TypeIntentOptionalAttributesName
integer, intent(in) :: nTotalFaces
integer, intent(in) :: nSides
integer, intent(in) :: nFaceDofs
real(kind=rk), intent(in) :: faceRep(nTotalFaces,nFaceDofs,8,2)
real(kind=rk), intent(inout) :: faceFlux(nTotalFaces,nFaceDofs,8,2)
integer, intent(in) :: leftPos(nSides)
integer, intent(in) :: rightPos(nsides)
integer, intent(in) :: var(8)
real(kind=rk), intent(in) :: material_left(nSides,1,4)
real(kind=rk), intent(in) :: material_right(nSides,1,4)

integer :: iSide, left, right, iDof real(kind=rk) :: left_mu, right_mu real(kind=rk) :: left_epsi, right_epsi real(kind=rk) :: left_gam, right_gam real(kind=rk) :: left_chi, right_chi real(kind=rk) :: left_speedOfLight, right_speedOfLight real(kind=rk) :: inv_denom_mu, inv_denom_epsi ! --------------------------------------------------------------------------

! flux for D_X, B_Z, D_y and B_z do iDof = 1, nFaceDofs

do iSide = 1, nSides

 ! The position of the left and right element in the state
 ! vector.
 left = leftPos(iSide)
 right = rightPos(iSide)

 ! The material parameters of the left and right element
 left_mu = material_left(iSide,1,1)
 left_epsi = material_left(iSide,1,2)
 left_gam = material_left(iSide,1,3)
 left_chi = material_left(iSide,1,4)
 left_speedOfLight = 1.0_rk / sqrt( left_mu * left_epsi )
 right_mu = material_right(iSide,1,1)
 right_epsi = material_right(iSide,1,2)
 right_gam = material_right(iSide,1,3)
 right_chi = material_right(iSide,1,4)
 right_speedOfLight = 1.0_rk / sqrt( right_mu * right_epsi )

 ! The inverse of the denominators
 inv_denom_mu = 1.0_rk / ((-1.0_rk)*left_speedOfLight*left_mu - right_speedOfLight*right_mu)
 inv_denom_epsi = 1.0_rk / (left_speedOfLight*left_epsi + right_speedOfLight*right_epsi)

 ! flux for D_x
 faceFlux(left,iDof,var(1),2) = inv_denom_mu * (                                   &
                                &       (-1.0_rk)*left_chi*faceRep(left,iDof,var(1),2)/left_epsi    &
                                &     - (-1.0_rk)*right_chi*faceRep(right,iDof,var(1),1)/right_epsi &
                                &              )                                   &
                                & + (                                              &
                                &        left_chi*left_chi*faceRep(left,iDof,var(7),2) &
                                &      + right_chi*right_chi*faceRep(right,iDof,var(7),1) &
                                &   ) / ( left_epsi*left_mu + right_epsi*right_mu  )

 ! flux for B_x
 faceFlux(left,iDof,var(4),2) = inv_denom_epsi * (                                 &
                                &        (-1.0_rk)*left_gam*faceRep(left,iDof,var(4),2)/left_mu &
                                &      - (-1.0_rk)*right_gam*faceRep(right,iDof,var(4),1)/right_mu) &
                                & + (                                              &
                                &       left_gam*left_gam*faceRep(left,iDof,var(8),2) &
                                &     + right_gam*right_gam*faceRep(right,iDof,var(8),1) &
                                &   ) / ( left_epsi*left_mu + right_epsi*right_mu  )

 ! flux for phi (electric correction)
 faceFlux(left,iDof,var(7),2) = (                                            &
                                &    left_speedOfLight*faceRep(left,iDof,var(1),2) &
                                &  + right_speedOfLight*faceRep(right,iDof,var(1),1) &
                                &  + left_speedOfLight*left_chi*faceRep(left,iDof,var(7),2) &
                                &  - right_speedOfLight*right_chi*faceRep(right,iDof,var(7),1)   &
                                & ) / (left_speedOfLight + right_speedOfLight)

 ! flux for psi (magnetic correction)
 faceFlux(left,iDof,var(8),2) = (                                            &
                                &    left_speedOfLight*faceRep(left,iDof,var(4),2)               &
                                &  + right_speedOfLight*faceRep(right,iDof,var(4),1)             &
                                &  + left_speedOfLight*left_gam*faceRep(left,iDof,var(8),2)      &
                                &  - right_speedOfLight*right_gam*faceRep(right,iDof,var(8),1)   &
                                & ) / (left_speedOfLight + right_speedOfLight)

 ! the flux for D_y
 faceFlux(left,iDof,var(2),2) = ( &
   &                  ( ((-1.0_rk)*faceRep(left,iDof,var(2),2) / left_epsi)        &
   &                    - ((-1.0_rk)*faceRep(right,iDof,var(2),1) / right_epsi)  ) &
   &               -  ( left_speedOfLight * faceRep(left,iDof,var(6),2)          &
   &                    + right_speedOfLight * faceRep(right,iDof,var(6),1) ))


 ! the flux for B_z
 faceFlux(left,iDof,var(6),2) = ( &
   &                  ( left_speedOfLight * faceRep(left,iDof,var(2),2)     &
   &                  + right_speedOfLight * faceRep(right,iDof,var(2),1) )&
   &                + ( ( faceRep(left,iDof,var(6),2) / left_mu )     &
   &                  - ( faceRep(right,iDof,var(6),1) / right_mu ) ))

 ! the flux for D_z
 faceFlux(left,iDof,var(3),2) = ( &
   &                  ( ( (-1.0_rk) * faceRep(left,iDof,var(3),2) / left_epsi )     &
   &                  - ( (-1.0_rk) * faceRep(right,iDof,var(3),1) / right_epsi ) ) &
   &                + ( left_speedOfLight * faceRep(left,iDof,var(5),2)           &
   &                  + right_speedOfLight * faceRep(right,iDof,var(5),1) )       &
   &                            )

 ! the flux for B_y
 faceFlux(left,iDof,var(5),2) = ( &
   &                  ( (-1.0_rk) * left_speedOfLight * faceRep(left,iDof,var(3),2) &
   &                  - right_speedOfLight * faceRep(right,iDof,var(3),1) )       &
   &                + ( ( faceRep(left,iDof,var(5),2) / left_mu )                 &
   &                  - ( faceRep(right,iDof,var(5),1) / right_mu) )              &
   &                            )

 ! Normalize the calculated fluxes
 faceFlux(left,iDof,var(2),2) = inv_denom_mu * faceFlux(left,iDof,var(2),2)
 faceFlux(left,iDof,var(6),2) = inv_denom_epsi * faceFlux(left,iDof,var(6),2)
 faceFlux(left,iDof,var(3),2) = inv_denom_mu * faceFlux(left,iDof,var(3),2)
 faceFlux(left,iDof,var(5),2) = inv_denom_epsi * faceFlux(left,iDof,var(5),2)

 ! Assign the same flux for both adjacent elements
 faceFlux(right,iDof,:,1) = faceFlux(left,iDof,:,2)

end do

end do

private subroutine maxwell_hc_flux_nonconst_cube_vec(nTotalFaces, nSides, nFaceDofs, faceRep, faceFlux, leftPos, rightPos, var, material_left, material_right, poly_proj, left_modalCoeffs, right_modalCoeffs, left_pntVal, right_pntVal, nodalNumFlux, numFluxBuffer)

calculate flux of maxwell equation with hyperbolic divergence

Read more…

Arguments

TypeIntentOptionalAttributesName
integer, intent(in) :: nTotalFaces
integer, intent(in) :: nSides
integer, intent(in) :: nFaceDofs
real(kind=rk), intent(in) :: faceRep(nTotalFaces,nFaceDofs,8,2)

The modal representation on the faces, left and right trace.

real(kind=rk), intent(inout) :: faceFlux(nTotalFaces,nFaceDofs,8,2)

The fluxes for all faces, for left and right elements.

integer, intent(in) :: leftPos(nSides)

Positions for the left and right elements of all faces

integer, intent(in) :: rightPos(nsides)

Positions for the left and right elements of all faces

integer, intent(in) :: var(8)

Variable rotation indices

real(kind=rk), intent(in) :: material_left(nSides,nFaceDofs,4)

Material parameters for the left faces.

real(kind=rk), intent(in) :: material_right(nSides,nFaceDofs,4)

Material parameters for the right faces.

type(ply_poly_project_type) :: poly_proj

Data for projection method !> Working array for the left and right modal coefficients real(kind=rk), intent(inout) :: left_modalCoeffs((fpt%nQuadPoints)2,8) real(kind=rk), intent(inout) :: right_modalCoeffs((fpt%nQuadPoints)2,8) !> Working array for the left and right point values real(kind=rk), intent(inout) :: left_pntVal((fpt%nQuadPoints)2,8) real(kind=rk), intent(inout) :: right_pntVal((fpt%nQuadPoints)2,8) !> Working array for the nodal flux real(kind=rk), intent(inout) :: nodalNumFlux((fpt%nQuadPoints)2,8) !> Working array for the modal numerical flux real(kind=rk), intent(inout) :: numFluxBuffer((fpt%nQuadPoints)2,8)

real(kind=rk), intent(inout), allocatable:: left_modalCoeffs(:,:)

Working array for the left and right modal coefficients

real(kind=rk), intent(inout), allocatable:: right_modalCoeffs(:,:)
real(kind=rk), intent(inout), allocatable:: left_pntVal(:,:)

Working array for the left and right point values

real(kind=rk), intent(inout), allocatable:: right_pntVal(:,:)
real(kind=rk), intent(inout), allocatable:: nodalNumFlux(:,:)

Working array for the nodal flux

real(kind=rk), intent(inout), allocatable:: numFluxBuffer(:,:)

Working array for the modal numerical flux