calculus Module

we do not want the full-grid sizes (lx1,lx2,lx3) in scope since we routinely need to do subgrid derivatives ROUTINES BELOW DO NOT ACCOUNT FOR METRIC FACTORS... As such they need to really be renamed to avoid confusion (they aren't curvilinear derivatives)


Uses

  • module~~calculus~~UsesGraph module~calculus calculus module~phys_consts phys_consts module~calculus->module~phys_consts module~mesh mesh module~calculus->module~mesh iso_fortran_env iso_fortran_env module~phys_consts->iso_fortran_env module~mesh->module~phys_consts

Used by

  • module~~calculus~~UsedByGraph module~calculus calculus module~sources sources module~sources->module~calculus module~potential_comm potential_comm module~potential_comm->module~calculus module~potential_mumps potential_mumps module~potential_comm->module~potential_mumps module~potential2d potential2d module~potential2d->module~calculus module~potential2d->module~potential_mumps module~gradient gradient module~gradient->module~calculus module~div div module~div->module~calculus module~multifluid multifluid module~multifluid->module~calculus module~multifluid->module~sources module~potential_mumps->module~calculus module~integral integral module~integral->module~calculus module~potential_root potential_root module~potential_root->module~potential_comm module~sources_mpi sources_mpi module~sources_mpi->module~sources program~gemini3d Gemini3D program~gemini3d->module~potential_comm program~gemini3d->module~multifluid module~potential_worker potential_worker module~potential_worker->module~potential_comm

Contents


Interfaces

public interface grad3D1

public interface grad3D2

public interface grad3D3

public interface div3D

public interface grad2D2

  • private function grad2D2_curv_23(f, x, lbnd, ubnd)

    Arguments

    Type IntentOptional AttributesName
    real(kind=wp), intent(in), dimension(:,:):: f
    type(curvmesh), intent(in) :: x
    integer, intent(in) :: lbnd
    integer, intent(in) :: ubnd

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

public interface grad2D3

  • private function grad2D3_curv_23(f, x, lbnd, ubnd)

    Arguments

    Type IntentOptional AttributesName
    real(kind=wp), intent(in), dimension(:,:):: f
    type(curvmesh), intent(in) :: x
    integer, intent(in) :: lbnd
    integer, intent(in) :: ubnd

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

public interface grad2D1_curv_alt

  • private function grad2D1_curv_alt_23(f, x, lbnd, ubnd)

    Arguments

    Type IntentOptional AttributesName
    real(kind=wp), intent(in), dimension(:,:):: f
    type(curvmesh), intent(in) :: x
    integer, intent(in) :: lbnd
    integer, intent(in) :: ubnd

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

public interface integral3D1

private interface integral2D1

private interface integral2D2

interface

  • private module function grad3D1_curv_3(f, x, lbnd1, ubnd1, lbnd2, ubnd2, lbnd3, ubnd3)

    Arguments

    Type IntentOptional AttributesName
    real(kind=wp), intent(in), dimension(:,:,:):: f
    type(curvmesh), intent(in) :: x
    integer, intent(in) :: lbnd1
    integer, intent(in) :: ubnd1
    integer, intent(in) :: lbnd2
    integer, intent(in) :: ubnd2
    integer, intent(in) :: lbnd3
    integer, intent(in) :: ubnd3

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

interface

  • private module function grad3D1_curv_23(f, x, lbnd1, ubnd1, lbnd2, ubnd2, lbnd3, ubnd3)

    Arguments

    Type IntentOptional AttributesName
    real(kind=wp), intent(in), dimension(:,:,:):: f
    type(curvmesh), intent(in) :: x
    integer, intent(in) :: lbnd1
    integer, intent(in) :: ubnd1
    integer, intent(in) :: lbnd2
    integer, intent(in) :: ubnd2
    integer, intent(in) :: lbnd3
    integer, intent(in) :: ubnd3

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

interface

  • private module function grad3D2_curv_3(f, x, lbnd1, ubnd1, lbnd2, ubnd2, lbnd3, ubnd3)

    Arguments

    Type IntentOptional AttributesName
    real(kind=wp), intent(in), dimension(:,:,:):: f
    type(curvmesh), intent(in) :: x
    integer, intent(in) :: lbnd1
    integer, intent(in) :: ubnd1
    integer, intent(in) :: lbnd2
    integer, intent(in) :: ubnd2
    integer, intent(in) :: lbnd3
    integer, intent(in) :: ubnd3

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

interface

  • private module function grad3D2_curv_23(f, x, lbnd1, ubnd1, lbnd2, ubnd2, lbnd3, ubnd3)

    Arguments

    Type IntentOptional AttributesName
    real(kind=wp), intent(in), dimension(:,:,:):: f
    type(curvmesh), intent(in) :: x
    integer, intent(in) :: lbnd1
    integer, intent(in) :: ubnd1
    integer, intent(in) :: lbnd2
    integer, intent(in) :: ubnd2
    integer, intent(in) :: lbnd3
    integer, intent(in) :: ubnd3

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

interface

  • private module function grad3D3_curv_3(f, x, lbnd1, ubnd1, lbnd2, ubnd2, lbnd3, ubnd3)

    Arguments

    Type IntentOptional AttributesName
    real(kind=wp), intent(in), dimension(:,:,:):: f
    type(curvmesh), intent(in) :: x
    integer, intent(in) :: lbnd1
    integer, intent(in) :: ubnd1
    integer, intent(in) :: lbnd2
    integer, intent(in) :: ubnd2
    integer, intent(in) :: lbnd3
    integer, intent(in) :: ubnd3

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

interface

  • private module function grad3D3_curv_23(f, x, lbnd1, ubnd1, lbnd2, ubnd2, lbnd3, ubnd3)

    Arguments

    Type IntentOptional AttributesName
    real(kind=wp), intent(in), dimension(:,:,:):: f
    type(curvmesh), intent(in) :: x
    integer, intent(in) :: lbnd1
    integer, intent(in) :: ubnd1
    integer, intent(in) :: lbnd2
    integer, intent(in) :: ubnd2
    integer, intent(in) :: lbnd3
    integer, intent(in) :: ubnd3

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

interface

  • private pure module function integral3D1_curv(f, x, lbnd, ubnd)

    Arguments

    Type IntentOptional AttributesName
    real(kind=wp), intent(in), dimension(:,:,:):: f
    type(curvmesh), intent(in) :: x
    integer, intent(in) :: lbnd
    integer, intent(in) :: ubnd

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

interface

  • public pure module function integral3D1_curv_alt(f, x, lbnd, ubnd)

    Arguments

    Type IntentOptional AttributesName
    real(kind=wp), intent(in), dimension(:,:,:):: f
    type(curvmesh), intent(in) :: x
    integer, intent(in) :: lbnd
    integer, intent(in) :: ubnd

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

interface

  • private pure module function integral2D1_curv(f, x, lbnd, ubnd)

    Arguments

    Type IntentOptional AttributesName
    real(kind=wp), intent(in), dimension(:,:):: f
    type(curvmesh), intent(in) :: x
    integer, intent(in) :: lbnd
    integer, intent(in) :: ubnd

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

interface

  • private pure module function integral2D1_curv_alt(f, x, lbnd, ubnd)

    Arguments

    Type IntentOptional AttributesName
    real(kind=wp), intent(in), dimension(:,:):: f
    type(curvmesh), intent(in) :: x
    integer, intent(in) :: lbnd
    integer, intent(in) :: ubnd

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

interface

  • private pure module function integral2D2_curv(f, x, lbnd, ubnd)

    Arguments

    Type IntentOptional AttributesName
    real(kind=wp), intent(in), dimension(:,:):: f
    type(curvmesh), intent(in) :: x
    integer, intent(in) :: lbnd
    integer, intent(in) :: ubnd

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

interface

  • private pure module function integral2D2_curv_alt(f, x, lbnd, ubnd)

    Arguments

    Type IntentOptional AttributesName
    real(kind=wp), intent(in), dimension(:,:):: f
    type(curvmesh), intent(in) :: x
    integer, intent(in) :: lbnd
    integer, intent(in) :: ubnd

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

interface

  • private module function div3D_curv_23(f1, f2, f3, x, lbnd1, ubnd1, lbnd2, ubnd2, lbnd3, ubnd3)

    Arguments

    Type IntentOptional AttributesName
    real(kind=wp), intent(in), dimension(:,:,:):: f1

    regardless of what has been passed into the function here, we assume these start at 1

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

    regardless of what has been passed into the function here, we assume these start at 1

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

    regardless of what has been passed into the function here, we assume these start at 1

    type(curvmesh), intent(in) :: x
    integer, intent(in) :: lbnd1
    integer, intent(in) :: ubnd1
    integer, intent(in) :: lbnd2
    integer, intent(in) :: ubnd2
    integer, intent(in) :: lbnd3
    integer, intent(in) :: ubnd3

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

interface

  • private module function div3D_curv_3(f1, f2, f3, x, lbnd1, ubnd1, lbnd2, ubnd2, lbnd3, ubnd3)

    Arguments

    Type IntentOptional AttributesName
    real(kind=wp), intent(in), dimension(:,:,:):: f1
    real(kind=wp), intent(in), dimension(:,:,:):: f2
    real(kind=wp), intent(in), dimension(:,:,:):: f3
    type(curvmesh), intent(in) :: x
    integer, intent(in) :: lbnd1
    integer, intent(in) :: ubnd1
    integer, intent(in) :: lbnd2
    integer, intent(in) :: ubnd2
    integer, intent(in) :: lbnd3
    integer, intent(in) :: ubnd3

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


Functions

public pure function ETD_uncoupled(f, P, L, dt)

Arguments

Type IntentOptional AttributesName
real(kind=wp), intent(in), dimension(:,:,:):: f
real(kind=wp), intent(in), dimension(:,:,:):: P
real(kind=wp), intent(in), dimension(:,:,:):: L
real(kind=wp), intent(in) :: dt

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

private pure function grad2D1_curv(f, x, lbnd, ubnd)

Arguments

Type IntentOptional AttributesName
real(kind=wp), intent(in), dimension(:,:):: f
type(curvmesh), intent(in) :: x
integer, intent(in) :: lbnd
integer, intent(in) :: ubnd

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

private pure function grad2D1_curv_alt_3(f, x, lbnd, ubnd)

Arguments

Type IntentOptional AttributesName
real(kind=wp), intent(in), dimension(:,:):: f
type(curvmesh), intent(in) :: x
integer, intent(in) :: lbnd
integer, intent(in) :: ubnd

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

private function grad2D1_curv_alt_23(f, x, lbnd, ubnd)

Arguments

Type IntentOptional AttributesName
real(kind=wp), intent(in), dimension(:,:):: f
type(curvmesh), intent(in) :: x
integer, intent(in) :: lbnd
integer, intent(in) :: ubnd

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

private pure function grad2D2_curv(f, x, lbnd, ubnd)

Arguments

Type IntentOptional AttributesName
real(kind=wp), intent(in), dimension(:,:):: f
type(curvmesh), intent(in) :: x
integer, intent(in) :: lbnd
integer, intent(in) :: ubnd

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

private function grad2D2_curv_23(f, x, lbnd, ubnd)

Arguments

Type IntentOptional AttributesName
real(kind=wp), intent(in), dimension(:,:):: f
type(curvmesh), intent(in) :: x
integer, intent(in) :: lbnd
integer, intent(in) :: ubnd

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

private pure function grad2D3_curv(f, x, lbnd, ubnd)

Arguments

Type IntentOptional AttributesName
real(kind=wp), intent(in), dimension(:,:):: f
type(curvmesh), intent(in) :: x
integer, intent(in) :: lbnd
integer, intent(in) :: ubnd

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

private function grad2D3_curv_23(f, x, lbnd, ubnd)

Arguments

Type IntentOptional AttributesName
real(kind=wp), intent(in), dimension(:,:):: f
type(curvmesh), intent(in) :: x
integer, intent(in) :: lbnd
integer, intent(in) :: ubnd

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

private pure function grad3D3_curv_periodic(f, x, lbnd, ubnd)

this assumes that the backward difference for hte first cell has been set the the same as the forward difference for the final cell, i.e. x%dx3all(lbnd)==x%dx3all(ubnd+1). In general when doing periodic grids it is probably best to hard code all of the differences outside the domain to be equal (or to use uniform meshes)

Arguments

Type IntentOptional AttributesName
real(kind=wp), intent(in), dimension(:,:,:):: f
type(curvmesh), intent(in) :: x
integer, intent(in) :: lbnd
integer, intent(in) :: ubnd

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

public pure function grad2D3_curv_periodic(f, x, lbnd, ubnd)

this assumes that the backward difference for hte first cell has been set the same as the forward difference for the final cell, i.e. x%dx3all(lbnd)==x%dx3all(ubnd+1). In general when doing periodic grids it is probably best to hard code all of the differences outside the domain to be equal (or to use uniform meshes)

Arguments

Type IntentOptional AttributesName
real(kind=wp), intent(in), dimension(:,:):: f
type(curvmesh), intent(in) :: x
integer, intent(in) :: lbnd
integer, intent(in) :: ubnd

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

public pure function chapman_a(alt, nmax, alt0, H)

Arguments

Type IntentOptional AttributesName
real(kind=wp), intent(in), dimension(:,:,:):: alt
real(kind=wp), intent(in) :: nmax
real(kind=wp), intent(in) :: alt0
real(kind=wp), intent(in), dimension(:,:,:):: H

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