uniform_gradient_field_mod Module

Provides an extension of the uniform field type which also appears to have a uniform gradient. This was written to allow some of the same code used when solving the plume in a Runge-Kutta solver, where I pass in uniform fields rather than ones using a Chebyshev grid.


Uses

  • module~~uniform_gradient_field_mod~~UsesGraph module~uniform_gradient_field_mod uniform_gradient_field_mod utils_mod utils_mod module~uniform_gradient_field_mod->utils_mod factual_mod factual_mod module~uniform_gradient_field_mod->factual_mod iso_fortran_env iso_fortran_env module~uniform_gradient_field_mod->iso_fortran_env

Used by

  • module~~uniform_gradient_field_mod~~UsedByGraph module~uniform_gradient_field_mod uniform_gradient_field_mod module~upstream_plume_mod upstream_plume_mod module~upstream_plume_mod->module~uniform_gradient_field_mod module~asymmetric_plume_mod asymmetric_plume_mod module~asymmetric_plume_mod->module~upstream_plume_mod module~static_plume_mod static_plume_mod module~static_plume_mod->module~upstream_plume_mod module~plume_mod plume_mod module~plume_mod->module~upstream_plume_mod

Contents


Interfaces

public interface uniform_gradient_field

  • private function constructor(val, grad) result(this)

    Author
    Chris MacMackin
    Date
    July 2017

    Creates a new scalar field with a uniform value across all of space but a non-zero gradient.

    Arguments

    Type IntentOptional AttributesName
    real(kind=r8), intent(in) :: val

    The value of the field

    real(kind=r8), intent(in), dimension(:):: grad

    An array in which the ith element contains the gradient in the _i_th direction. Directions corresponding to values of i greater than the size of the array are taken to have a gradient of zero.

    Return Value type(uniform_gradient_field)

    A scalar field initated based on teh arguments of this function.


Derived Types

type, public, extends(uniform_scalar_field) :: uniform_gradient_field

A type of uniform field which also has a uniform gradient. Of course, this is impossible in practice, but it can be useful for tricking certain routines into working properly. Ideally a whole new derived type would be created which just holds the value and gradient at a single point, but the emphasis is on getting something quickly. Note that the gradient is not propagated across operations--the result of all overloaded operators is just a normal uniform field with no gradient.

Components

TypeVisibility AttributesNameInitial
real(kind=r8), public, dimension(:), allocatable:: grad

The values of the gradient in each direction.

Constructor

private function constructor(val, grad)

Creates a new scalar field with a uniform value across all of space but a non-zero gradient.

Type-Bound Procedures

procedure, public :: d_dx => uniform_gradient_d_dx

procedure, private :: gradient => uniform_gradient_gradient

procedure, private :: is_equal => uniform_gradient_is_equal

Checks fields are equal within a tolerance

procedure, private :: assign_field => uniform_gradient_assign


Functions

private function constructor(val, grad) result(this)

Author
Chris MacMackin
Date
July 2017

Creates a new scalar field with a uniform value across all of space but a non-zero gradient.

Arguments

Type IntentOptional AttributesName
real(kind=r8), intent(in) :: val

The value of the field

real(kind=r8), intent(in), dimension(:):: grad

An array in which the ith element contains the gradient in the _i_th direction. Directions corresponding to values of i greater than the size of the array are taken to have a gradient of zero.

Return Value type(uniform_gradient_field)

A scalar field initated based on teh arguments of this function.

private function uniform_gradient_d_dx(this, dir, order) result(res)

Author
Chris MacMackin
Date
July 2017

Arguments

Type IntentOptional AttributesName
class(uniform_gradient_field), intent(in) :: this
integer, intent(in) :: dir

Direction in which to differentiate

integer, intent(in), optional :: order

Order of the derivative, default = 1

Return Value class(scalar_field), pointer

private function uniform_gradient_gradient(this) result(res)

Author
Chris MacMackin
Date
July 2017

Arguments

Type IntentOptional AttributesName
class(uniform_gradient_field), intent(in) :: this

Return Value class(vector_field), pointer

The result of this operation

private function uniform_gradient_is_equal(this, rhs) result(iseq)

Author
Chris MacMackin
Date
July 2017

Evaluates whether two scalar fields are equal within a tolerance, specified by set_tol.

Arguments

Type IntentOptional AttributesName
class(uniform_gradient_field), intent(in) :: this
class(scalar_field), intent(in) :: rhs

Return Value logical


Subroutines

private impure elemental subroutine uniform_gradient_assign(this, rhs)

Author
Chris MacMackin
Date
July 2017

Arguments

Type IntentOptional AttributesName
class(uniform_gradient_field), intent(inout) :: this
class(scalar_field), intent(in) :: rhs