diffusion Module

This module sets up the ionospheric diffusion problem and then passes it off to the parabolic solvers.

sizes are not imported in case we want to do subgrid diffusion


Uses

  • module~~diffusion~~UsesGraph module~diffusion diffusion module~phys_consts phys_consts module~diffusion->module~phys_consts module~pdeparabolic PDEparabolic module~diffusion->module~pdeparabolic module~grid grid module~diffusion->module~grid module~mesh mesh module~diffusion->module~mesh iso_fortran_env iso_fortran_env module~phys_consts->iso_fortran_env module~pdeparabolic->module~phys_consts module~vendor_lapack95 vendor_lapack95 module~pdeparabolic->module~vendor_lapack95 module~grid->module~phys_consts module~grid->module~mesh module~mpimod mpimod module~grid->module~mpimod module~reader reader module~grid->module~reader module~grid->iso_fortran_env module~mesh->module~phys_consts module~mpimod->module~phys_consts module~mpimod->iso_fortran_env mpi mpi module~mpimod->mpi module~reader->module~phys_consts module~reader->iso_fortran_env module~vendor_lapack95->iso_fortran_env

Used by

  • module~~diffusion~~UsedByGraph module~diffusion diffusion module~multifluid multifluid module~multifluid->module~diffusion program~gemini3d Gemini3D program~gemini3d->module~multifluid

Contents


Interfaces

public interface TRBDF23D

  • private function TRBDF23D_curv(f, A, B, C, D, E, dt, x)

    SOLVE A 3D SEQUENCE OF 1D DIFFUSION PROBLEMS. GHOST CELLS ARE ACCOMODATED AS THEY PROVIDE A CONVENIENT MEMORY SPACE FOR BOUNDARY CONDITIONS.

    Read more…

    Arguments

    Type IntentOptional AttributesName
    real(kind=wp), intent(in), dimension(-1:,-1:,-1:):: f

    expected to include ghosts

    real(kind=wp), intent(in), dimension(:,:,:):: A

    trimmed to grid size

    real(kind=wp), intent(in), dimension(:,:,:):: B

    trimmed to grid size

    real(kind=wp), intent(in), dimension(:,:,:):: C

    trimmed to grid size

    real(kind=wp), intent(in), dimension(:,:,:):: D

    trimmed to grid size

    real(kind=wp), intent(in), dimension(:,:,:):: E

    trimmed to grid size

    real(kind=wp), intent(in) :: dt
    type(curvmesh), intent(in) :: x

    Return Value real(kind=wp), dimension(-1:size(f,1)-2,-1:size(f,2)-2,-1:size(f,3)-2)

public interface backEuler3D

  • private function backEuler3D_curv(f, A, B, C, D, E, dt, x)

    SOLVE A 3D SEQUENCE OF 1D DIFFUSION PROBLEMS. GHOST CELLS ARE ACCOMODATED AS THEY PROVIDE A CONVENIENT MEMORY SPACE FOR BOUNDARY CONDITIONS.

    Read more…

    Arguments

    Type IntentOptional AttributesName
    real(kind=wp), intent(in), dimension(-1:,-1:,-1:):: f
    real(kind=wp), intent(in), dimension(:,:,:):: A
    real(kind=wp), intent(in), dimension(:,:,:):: B
    real(kind=wp), intent(in), dimension(:,:,:):: C
    real(kind=wp), intent(in), dimension(:,:,:):: D
    real(kind=wp), intent(in), dimension(:,:,:):: E
    real(kind=wp), intent(in) :: dt
    type(curvmesh), intent(in) :: x

    Return Value real(kind=wp), dimension(-1:size(f,1)-2,-1:size(f,2)-2,-1:size(f,3)-2)


Functions

private function backEuler3D_curv(f, A, B, C, D, E, dt, x)

SOLVE A 3D SEQUENCE OF 1D DIFFUSION PROBLEMS. GHOST CELLS ARE ACCOMODATED AS THEY PROVIDE A CONVENIENT MEMORY SPACE FOR BOUNDARY CONDITIONS.

Read more…

Arguments

Type IntentOptional AttributesName
real(kind=wp), intent(in), dimension(-1:,-1:,-1:):: f
real(kind=wp), intent(in), dimension(:,:,:):: A
real(kind=wp), intent(in), dimension(:,:,:):: B
real(kind=wp), intent(in), dimension(:,:,:):: C
real(kind=wp), intent(in), dimension(:,:,:):: D
real(kind=wp), intent(in), dimension(:,:,:):: E
real(kind=wp), intent(in) :: dt
type(curvmesh), intent(in) :: x

Return Value real(kind=wp), dimension(-1:size(f,1)-2,-1:size(f,2)-2,-1:size(f,3)-2)

private function TRBDF23D_curv(f, A, B, C, D, E, dt, x)

SOLVE A 3D SEQUENCE OF 1D DIFFUSION PROBLEMS. GHOST CELLS ARE ACCOMODATED AS THEY PROVIDE A CONVENIENT MEMORY SPACE FOR BOUNDARY CONDITIONS.

Read more…

Arguments

Type IntentOptional AttributesName
real(kind=wp), intent(in), dimension(-1:,-1:,-1:):: f

expected to include ghosts

real(kind=wp), intent(in), dimension(:,:,:):: A

trimmed to grid size

real(kind=wp), intent(in), dimension(:,:,:):: B

trimmed to grid size

real(kind=wp), intent(in), dimension(:,:,:):: C

trimmed to grid size

real(kind=wp), intent(in), dimension(:,:,:):: D

trimmed to grid size

real(kind=wp), intent(in), dimension(:,:,:):: E

trimmed to grid size

real(kind=wp), intent(in) :: dt
type(curvmesh), intent(in) :: x

Return Value real(kind=wp), dimension(-1:size(f,1)-2,-1:size(f,2)-2,-1:size(f,3)-2)


Subroutines

public pure subroutine diffusion_prep(isp, x, lambda, betacoeff, ns, T, A, B, C, D, E, Tn, Teinf)

COMPUTE COEFFICIENTS IN DIFFUSION EQUATIONS AND LOAD UP GHOST CELLS

Read more…

Arguments

Type IntentOptional AttributesName
integer, intent(in) :: isp
type(curvmesh), intent(in) :: x
real(kind=wp), intent(in), dimension(:,:,:):: lambda
real(kind=wp), intent(in), dimension(:,:,:):: betacoeff
real(kind=wp), intent(in), dimension(-1:,-1:,-1:):: ns
real(kind=wp), intent(inout), dimension(-1:,-1:,-1:):: T
real(kind=wp), intent(out), dimension(1:size(ns,1)-4,1:size(ns,2)-4,1:size(ns,3)-4):: A
real(kind=wp), intent(out), dimension(1:size(ns,1)-4,1:size(ns,2)-4,1:size(ns,3)-4):: B
real(kind=wp), intent(out), dimension(1:size(ns,1)-4,1:size(ns,2)-4,1:size(ns,3)-4):: C
real(kind=wp), intent(out), dimension(1:size(ns,1)-4,1:size(ns,2)-4,1:size(ns,3)-4):: D
real(kind=wp), intent(out), dimension(1:size(ns,1)-4,1:size(ns,2)-4,1:size(ns,3)-4):: E
real(kind=wp), intent(in), dimension(:,:,:):: Tn
real(kind=wp), intent(in) :: Teinf