mod_domain Module

module~~mod_domain~~UsesGraph module~mod_domain mod_domain json_module json_module json_module->module~mod_domain module~mod_spectrum mod_spectrum json_module->module~mod_spectrum module~mod_spectrum->module~mod_domain datetime_module datetime_module datetime_module->module~mod_domain datetime_module->module~mod_spectrum module~mod_grid mod_grid module~mod_grid->module~mod_domain module~mod_const mod_const module~mod_const->module~mod_domain module~mod_const->module~mod_spectrum module~mod_precision mod_precision module~mod_precision->module~mod_domain module~mod_precision->module~mod_spectrum module~mod_precision->module~mod_grid module~mod_precision->module~mod_const module~mod_utility mod_utility module~mod_precision->module~mod_utility module~mod_utility->module~mod_spectrum module~mod_utility->module~mod_grid iso_c_binding iso_c_binding iso_c_binding->module~mod_spectrum iso_fortran_env iso_fortran_env iso_fortran_env->module~mod_precision
Help

Used By

module~~mod_domain~~UsedByGraph module~mod_domain mod_domain module~mod_time_integration mod_time_integration module~mod_domain->module~mod_time_integration
Help


Interfaces

public interface domain_type

  • private function constructor(grid, spectrum, shallow_water_mode) result(domain)

    Constructor function for the domain object.

    Arguments

    Type IntentOptional AttributesName
    type(grid_type), intent(in) :: grid

    Input grid instance

    type(spectrum_type), intent(in) :: spectrum

    Input spectrum instance

    logical, intent(in), optional :: shallow_water_mode

    Logical switch to enable shallow water solver

    Return Value type(domain_type)


Derived Types

type, public :: domain_type

Components

TypeVisibility AttributesNameInitial
character(len=:), private, allocatable:: type_name
type(grid_type), private :: grid
type(spectrum_type), private, dimension(:,:), allocatable:: spectrum
logical, private :: shallow_water_mode
type(datetime), private :: start_time

Simulation start time

type(datetime), private :: end_time

Simulation end time

type(timedelta), private :: time_step

Time step [s]

real(kind=rk), private, dimension(:,:), allocatable:: dx
real(kind=rk), private, dimension(:,:), allocatable:: dy
real(kind=rk), private, dimension(:,:), allocatable:: u
real(kind=rk), private, dimension(:,:), allocatable:: v
real(kind=rk), private, dimension(:,:), allocatable:: eta
real(kind=rk), private, dimension(:,:), allocatable:: depth
integer(kind=ik), private, dimension(2):: lb
integer(kind=ik), private, dimension(2):: ub
integer(kind=ik), private :: nfreqs
integer(kind=ik), private :: ndirs

Constructor

private function constructor(grid, spectrum, shallow_water_mode)

Constructor function for the domain object.

Type-Bound Procedures

procedure, public, pass(self) :: frequencyMoment
procedure, public, pass(self) :: getCurrent_u
procedure, public, pass(self) :: getCurrent_v
procedure, public, pass(self) :: getGravity
procedure, public, pass(self) :: getGrid
procedure, public, pass(self) :: getGridSpacingXWithHalo
procedure, public, pass(self) :: getGridSpacingYWithHalo
procedure, public, pass(self) :: getDepth
procedure, public, pass(self) :: getElevation
procedure, public, pass(self) :: getFrequency
procedure, public, pass(self) :: getDirections
procedure, public, pass(self) :: getLowerBounds
procedure, public, pass(self) :: getUpperBounds
procedure, public, pass(self) :: getSpectrum
procedure, public, pass(self) :: getSpectrumArray
procedure, public, pass(self) :: getPhaseSpeed
procedure, public, pass(self) :: getGroupSpeed
procedure, public, pass(self) :: getSurfaceTension
procedure, public, pass(self) :: getAirDensity
procedure, public, pass(self) :: getWaterDensity
procedure, public, pass(self) :: isAllocated
procedure, public, pass(self) :: meanPeriod
procedure, public, pass(self) :: meanPeriodZeroCrossing
procedure, public, pass(self) :: setDepth
procedure, public, pass(self) :: setElevation
procedure, public, pass(self) :: setGravity
procedure, public, pass(self) :: setSurfaceTension
procedure, public, pass(self) :: setAirDensity
procedure, public, pass(self) :: setWaterDensity
procedure, public, pass(self) :: significantWaveHeight
procedure, public, pass(self) :: wavenumberMoment
procedure, public, pass(self) :: writeJSON
procedure, private, pass(self) :: advect1dRank1
procedure, private, pass(self) :: advect1dRank2
procedure, private, pass(self) :: advect2dRank2
procedure, private, pass(self) :: assign_spectrum_array_1d
procedure, private, pass(self) :: assign_spectrum_array_2d
procedure, private, pass(self) :: domain_add_domain
procedure, private, pass(self) :: domain_add_real
procedure, private, pass(self) :: domain_sub_domain
procedure, private, pass(self) :: domain_sub_real
procedure, private, pass(self) :: domain_mult_domain
procedure, private, pass(self) :: domain_mult_real
procedure, private, pass(self) :: domain_div_domain
procedure, private, pass(self) :: domain_div_real
procedure, private, pass(self) :: domain_unary_minus
procedure, private, pass(self) :: real_add_domain
procedure, private, pass(self) :: real_sub_domain
procedure, private, pass(self) :: real_mult_domain
procedure, private, pass(self) :: real_div_domain
procedure, private, pass(self) :: eq
procedure, private, pass(self) :: neq
procedure, private, pass(self) :: setSpectrum1d
procedure, private, pass(self) :: setSpectrum2d
procedure, private, pass(self) :: setSpectrumArray1d1d
procedure, private, pass(self) :: setSpectrumArray1d2d
procedure, private, pass(self) :: setSpectrumArray2d2d
generic, public :: advect => advect1dRank1, advect1dRank2, advect2dRank2
generic, public :: setSpectrum => setSpectrum1d, setSpectrum2d
generic, public :: setSpectrumArray => setSpectrumArray1d1d, setSpectrumArray1d2d, setSpectrumArray2d2d
generic, public :: assignment(=) => assign_spectrum_array_1d, assign_spectrum_array_2d
generic, public :: operator(+) => domain_add_domain, domain_add_real, real_add_domain
generic, public :: operator(-) => domain_sub_domain, domain_sub_real, domain_unary_minus, real_sub_domain
generic, public :: operator(*) => domain_mult_domain, domain_mult_real, real_mult_domain
generic, public :: operator(/) => domain_div_domain, domain_div_real, real_div_domain
generic, public :: operator(==) => eq
generic, public :: operator(/=) => neq

Functions

private function constructor(grid, spectrum, shallow_water_mode) result(domain)

Constructor function for the domain object.

Arguments

Type IntentOptional AttributesName
type(grid_type), intent(in) :: grid

Input grid instance

type(spectrum_type), intent(in) :: spectrum

Input spectrum instance

logical, intent(in), optional :: shallow_water_mode

Logical switch to enable shallow water solver

Return Value type(domain_type)

private pure function advect1dRank1(self, advection_method, halowidth, directional_type) result(adv)

Computes the advective tendency for the domain instance given the desired advection method as an input function and the number of halo cells. This function works only in cases where ndirs == 1.

Arguments

Type IntentOptional AttributesName
class(domain_type), intent(in) :: self

domain instance

public pure function advection_method(f, u, dx) result(tendency)
Arguments
Type IntentOptional AttributesName
real(kind=rk), intent(in), dimension(:,:):: f
real(kind=rk), intent(in), dimension(:,:):: u
real(kind=rk), intent(in), dimension(:):: dx
Return Value real(kind=rk), dimension(:,:), allocatable
integer(kind=ik), intent(in) :: halowidth

number of halo cells to use in the advection method

integer(kind=ik), intent(in), dimension(:):: directional_type

A global constant that helps resolve the interface of this specific prodedure

Return Value type(domain_type)

private pure function advect1dRank2(self, advection_method, halowidth, directional_type) result(adv)

Computes the advective tendency for the domain instance given the desired advection method as an input function and the number of halo cells. This function works both when ndirs == 1 (omnidirectional) and when ndirs > 1 (directional).

Arguments

Type IntentOptional AttributesName
class(domain_type), intent(in) :: self

domain instance

public pure function advection_method(f, u, dx) result(tendency)
Arguments
Type IntentOptional AttributesName
real(kind=rk), intent(in), dimension(:,:,:):: f
real(kind=rk), intent(in), dimension(:,:,:):: u
real(kind=rk), intent(in), dimension(:):: dx
Return Value real(kind=rk), dimension(:,:,:), allocatable
integer(kind=ik), intent(in) :: halowidth

number of halo cells to use in the advection method

integer(kind=ik), intent(in), dimension(:,:):: directional_type

A global constant that helps resolve the interface of this specific prodedure

Return Value type(domain_type)

private pure function advect2dRank2(self, advection_method, halowidth) result(adv)

Computes the advective tendency for the domain instance given the desired advection method as an input function and the number of halo cells. This function works both when ndirs == 1 (omnidirectional) and when ndirs > 1 (directional).

Arguments

Type IntentOptional AttributesName
class(domain_type), intent(in) :: self

domain instance

public pure function advection_method(f, u, v, dx, dy) result(tendency)
Arguments
Type IntentOptional AttributesName
real(kind=rk), intent(in), dimension(:,:,:,:):: f
real(kind=rk), intent(in), dimension(:,:,:,:):: u
real(kind=rk), intent(in), dimension(:,:,:,:):: v
real(kind=rk), intent(in), dimension(:,:):: dx
real(kind=rk), intent(in), dimension(:,:):: dy
Return Value real(kind=rk), dimension(:,:,:,:), allocatable
integer(kind=ik), intent(in), dimension(:):: halowidth

number of halo cells to use in the advection method

Return Value type(domain_type)

private pure elemental function isAllocated(self)

Returns the allocation status of the domains sub-components.

Arguments

Type IntentOptional AttributesName
class(domain_type), intent(in) :: self

domain instance

Return Value logical

private pure elemental function eq(self, d2)

Logical equality comparison function. Overloads the /= operator.

Arguments

Type IntentOptional AttributesName
class(domain_type), intent(in) :: self

l.h.s. domain instance

class(domain_type), intent(in) :: d2

r.h.s. domain instance

Return Value logical

private pure elemental function neq(self, d2)

Logical inequality comparison function. Overloads the /= operator.

Arguments

Type IntentOptional AttributesName
class(domain_type), intent(in) :: self

l.h.s. domain instance

class(domain_type), intent(in) :: d2

r.h.s. domain instance

Return Value logical

private pure elemental function domain_add_domain(self, d2) result(domain)

Returns a sum of two domain instances.

Arguments

Type IntentOptional AttributesName
class(domain_type), intent(in) :: self

l.h.s. domain instance

class(domain_type), intent(in) :: d2

r.h.s. domain instance

Return Value type(domain_type)

private pure elemental function domain_sub_domain(self, d2) result(domain)

Returns a difference between two domain instances.

Arguments

Type IntentOptional AttributesName
class(domain_type), intent(in) :: self

l.h.s. domain instance

class(domain_type), intent(in) :: d2

r.h.s. domain instance

Return Value type(domain_type)

private pure elemental function domain_unary_minus(self) result(domain)

Returns a negative domain instances.

Arguments

Type IntentOptional AttributesName
class(domain_type), intent(in) :: self

domain instance

Return Value type(domain_type)

private pure elemental function domain_mult_domain(self, d2) result(domain)

Returns a product of two domain instances.

Arguments

Type IntentOptional AttributesName
class(domain_type), intent(in) :: self

l.h.s. domain instance

class(domain_type), intent(in) :: d2

r.h.s. domain instance

Return Value type(domain_type)

private pure elemental function domain_div_domain(self, d2) result(domain)

Returns a division of two domain instances.

Arguments

Type IntentOptional AttributesName
class(domain_type), intent(in) :: self

l.h.s. domain instance

class(domain_type), intent(in) :: d2

r.h.s. domain instance

Return Value type(domain_type)

private pure elemental function domain_add_real(self, a) result(domain)

Returns a sum of a domain instance and a real number.

Arguments

Type IntentOptional AttributesName
class(domain_type), intent(in) :: self

l.h.s. domain instance

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

r.h.s. real number

Return Value type(domain_type)

private pure elemental function domain_sub_real(self, a) result(domain)

Returns a difference between a domain instance and a real number.

Arguments

Type IntentOptional AttributesName
class(domain_type), intent(in) :: self

l.h.s. domain instance

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

r.h.s. real number

Return Value type(domain_type)

private pure elemental function domain_mult_real(self, a) result(domain)

Returns a product of a domain instance and a real number.

Arguments

Type IntentOptional AttributesName
class(domain_type), intent(in) :: self

l.h.s. domain instance

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

r.h.s. real number

Return Value type(domain_type)

private pure elemental function domain_div_real(self, a) result(domain)

Returns a division of a domain instance and a real number.

Arguments

Type IntentOptional AttributesName
class(domain_type), intent(in) :: self

l.h.s. domain instance

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

r.h.s. real number

Return Value type(domain_type)

private pure elemental function real_add_domain(a, self) result(domain)

Returns a sum of a real number and a domain instance.

Arguments

Type IntentOptional AttributesName
real(kind=rk), intent(in) :: a

l.h.s. real number

class(domain_type), intent(in) :: self

r.h.s. domain instance

Return Value type(domain_type)

private pure elemental function real_sub_domain(a, self) result(domain)

Returns a difference between a real number and a domain instance.

Arguments

Type IntentOptional AttributesName
real(kind=rk), intent(in) :: a

l.h.s. real number

class(domain_type), intent(in) :: self

r.h.s. domain instance

Return Value type(domain_type)

private pure elemental function real_mult_domain(a, self) result(domain)

Returns a product of a real number and a domain instance.

Arguments

Type IntentOptional AttributesName
real(kind=rk), intent(in) :: a

l.h.s. real number

class(domain_type), intent(in) :: self

r.h.s. domain instance

Return Value type(domain_type)

private pure elemental function real_div_domain(a, self) result(domain)

Returns a product of a real number and a domain instance.

Arguments

Type IntentOptional AttributesName
real(kind=rk), intent(in) :: a

l.h.s. real number

class(domain_type), intent(in) :: self

r.h.s. domain instance

Return Value type(domain_type)

private pure function getCurrent_u(self) result(u)

Returns the 3-d array with values of Eulerian velocity (mean current) in x-direction [m/s].

Arguments

Type IntentOptional AttributesName
class(domain_type), intent(in) :: self

Domain instance

Return Value real(kind=rk), dimension(:,:,:), allocatable

Eulerian u-velocity [m/s]

private pure function getCurrent_v(self) result(v)

Returns the 3-d array with values of Eulerian velocity (mean current) in y-direction [m/s].

Arguments

Type IntentOptional AttributesName
class(domain_type), intent(in) :: self

Domain instance

Return Value real(kind=rk), dimension(:,:,:), allocatable

Eulerian v-velocity [m/s]

private pure function getDepth(self) result(depth)

Returns the mean water depth [m] array.

Arguments

Type IntentOptional AttributesName
class(domain_type), intent(in) :: self

Domain instance

Return Value real(kind=rk), dimension(:,:), allocatable

Mean water depth [m]

private pure function getElevation(self) result(elevation)

Returns the mean water elevation [m] array.

Arguments

Type IntentOptional AttributesName
class(domain_type), intent(in) :: self

Domain instance

Return Value real(kind=rk), dimension(:,:), allocatable

Mean water elevation [m]

private pure function getFrequency(self) result(frequency)

Returns the frequency [Hz] array.

Arguments

Type IntentOptional AttributesName
class(domain_type), intent(in) :: self

Domain instance

Return Value real(kind=rk), dimension(:), allocatable

Frequency [Hz]

private pure function getDirections(self) result(directions)

Returns the spectral direction bins [rad].

Arguments

Type IntentOptional AttributesName
class(domain_type), intent(in) :: self

Domain instance

Return Value real(kind=rk), dimension(:), allocatable

Directions [rad]

private pure function getGravity(self) result(grav)

Returns the gravitational acceleration [m/s^2] array.

Arguments

Type IntentOptional AttributesName
class(domain_type), intent(in) :: self

Domain instance

Return Value real(kind=rk), dimension(:,:), allocatable

Gravitational acceleration [m/s^2]

private pure function getGrid(self) result(grid)

Returns the grid instance that is the component of the domain.

Arguments

Type IntentOptional AttributesName
class(domain_type), intent(in) :: self

Domain instance

Return Value type(grid_type)

Grid instance component

private pure function getLowerBounds(self) result(lb)

Returns the lower bounds of the domain instance.

Arguments

Type IntentOptional AttributesName
class(domain_type), intent(in) :: self

Domain instance

Return Value integer(kind=ik), dimension(2)

Lower bound indices

private pure function getUpperBounds(self) result(ub)

Returns the upper bounds of the domain instance.

Arguments

Type IntentOptional AttributesName
class(domain_type), intent(in) :: self

Domain instance

Return Value integer(kind=ik), dimension(2)

Upper bound indices

private pure function getSpectrum(self) result(spectrum)

Returns the array of spectrum instances.

Arguments

Type IntentOptional AttributesName
class(domain_type), intent(in) :: self

Domain instance

Return Value type(spectrum_type), dimension(:,:), allocatable

Array of spectrum instances

private pure function getSpectrumArray(self, halowidth, periodic) result(spectrum_array)

Returns a 4-dimensional spectrum array, where the first two dimensions are frequency and directional dimensions and the second two are spatial x and y dimensions.

Arguments

Type IntentOptional AttributesName
class(domain_type), intent(in) :: self

Domain instance

integer(kind=ik), intent(in), dimension(2):: halowidth

Integers indicating how many cells to allocate for halo points

logical, intent(in) :: periodic

If .true., halo cells will be updated with values corresponding to periodic boundary conditions

Return Value real(kind=rk), dimension(:,:,:,:), allocatable

Spectrum array

private pure function getPhaseSpeed(self) result(cp)

Returns a 3-d array with phase speed values [m/s].

Arguments

Type IntentOptional AttributesName
class(domain_type), intent(in) :: self

Domain instance

Return Value real(kind=rk), dimension(:,:,:), allocatable

Phase speed [m/s] array

private pure function getGroupSpeed(self, halowidth, periodic) result(cg)

Returns a 3-d array with group speed values [m/s].

Arguments

Type IntentOptional AttributesName
class(domain_type), intent(in) :: self

Domain instance

integer(kind=ik), intent(in), dimension(2):: halowidth

Integers indicating how many cells to allocate for halo points

logical, intent(in) :: periodic

If .true., halo cells will be updated with values corresponding to periodic boundary conditions

Return Value real(kind=rk), dimension(:,:,:), allocatable

Group speed [m/s] array

private pure function getGridSpacingXWithHalo(self, halowidth, periodic) result(dx)

Returns grid spacing array in x-direction including halo cells.

Arguments

Type IntentOptional AttributesName
class(domain_type), intent(in) :: self

Domain instance

integer(kind=ik), intent(in), dimension(2):: halowidth

Integer width of halo region

logical, intent(in) :: periodic

If .true., halo cells will be updated with values corresponding to periodic boundary conditions

Return Value real(kind=rk), dimension(:,:), allocatable

Grid spacing in x [m]

private pure function getGridSpacingYWithHalo(self, halowidth, periodic) result(dy)

Returns grid spacing array in y-direction including halo cells.

Arguments

Type IntentOptional AttributesName
class(domain_type), intent(in) :: self

Domain instance

integer(kind=ik), intent(in), dimension(2):: halowidth

Integer width of halo region

logical, intent(in) :: periodic

If .true., halo cells will be updated with values corresponding to periodic boundary conditions

Return Value real(kind=rk), dimension(:,:), allocatable

Grid spacing in y [m]

private pure function getSurfaceTension(self) result(surface_tension)

Returns the surface tension [N/m].

Arguments

Type IntentOptional AttributesName
class(domain_type), intent(in) :: self

Domain instance

Return Value real(kind=rk), dimension(:,:), allocatable

Surface tension [N/m]

private pure function getAirDensity(self) result(air_density)

Returns the air density [kg/m^3].

Arguments

Type IntentOptional AttributesName
class(domain_type), intent(in) :: self

Domain instance

Return Value real(kind=rk), dimension(:,:), allocatable

Air density [kg/m^3]

private pure function getWaterDensity(self) result(water_density)

Returns the water density [kg/m^3].

Arguments

Type IntentOptional AttributesName
class(domain_type), intent(in) :: self

Domain instance

Return Value real(kind=rk), dimension(:,:), allocatable

Water density [kg/m^3]

private pure function frequencyMoment(self, n) result(moment)

Returns the spectral frequency moment of order n.

Arguments

Type IntentOptional AttributesName
class(domain_type), intent(in) :: self

Spectrum instance

integer(kind=ik), intent(in) :: n

Order

Return Value real(kind=rk), dimension(:,:), allocatable

private pure function wavenumberMoment(self, n) result(moment)

Returns the spectral frequency moment of order n.

Arguments

Type IntentOptional AttributesName
class(domain_type), intent(in) :: self

Spectrum instance

integer(kind=ik), intent(in) :: n

Order

Return Value real(kind=rk), dimension(:,:), allocatable

private pure function meanPeriod(self)

Returns the mean wave period [s] for the whole domain.

Arguments

Type IntentOptional AttributesName
class(domain_type), intent(in) :: self

Domain instance

Return Value real(kind=rk), dimension(:,:), allocatable

Mean period [s] array

private pure function meanPeriodZeroCrossing(self)

Returns the zero-crossing mean wave period [s] for the whole domain.

Arguments

Type IntentOptional AttributesName
class(domain_type), intent(in) :: self

Domain instance

Return Value real(kind=rk), dimension(:,:), allocatable

Mean period [s] array

private pure function significantWaveHeight(self) result(hs)

Returns the significant wave height [m] for the whole domain.

Arguments

Type IntentOptional AttributesName
class(domain_type), intent(in) :: self

Domain instance

Return Value real(kind=rk), dimension(:,:), allocatable

Significant wave height [m] array


Subroutines

private pure subroutine assign_spectrum_array_1d(self, spectrum_array)

Assigns a 1-d array of spectrum instances to a domain instance. This procedure overloads the assignment ('=') operator.

Arguments

Type IntentOptional AttributesName
class(domain_type), intent(inout) :: self

l.h.s. domain instance

class(spectrum_type), intent(in), dimension(:):: spectrum_array

r.h.s. array of spectrum instances

private pure subroutine assign_spectrum_array_2d(self, spectrum_array)

Assigns a 2-d array of spectrum instances to a domain instance. This procedure overloads the assignment ('=') operator.

Arguments

Type IntentOptional AttributesName
class(domain_type), intent(inout) :: self

l.h.s. domain instance

class(spectrum_type), intent(in), dimension(:,:):: spectrum_array

r.h.s. array of spectrum instances

private pure subroutine setDepth(self, depth)

Sets the mean water depth [m].

Arguments

Type IntentOptional AttributesName
class(domain_type), intent(inout) :: self

Domain instance

real(kind=rk), intent(in), dimension(:,:):: depth

Mean water depth [m]

private pure subroutine setElevation(self, elevation)

Sets the mean water elevation [m].

Arguments

Type IntentOptional AttributesName
class(domain_type), intent(inout) :: self

Domain instance

real(kind=rk), intent(in), dimension(:,:):: elevation

Mean water elevation [m]

private pure subroutine setGravity(self, grav)

Sets the gravitational acceleration [m/s^2].

Arguments

Type IntentOptional AttributesName
class(domain_type), intent(inout) :: self

Domain instance

real(kind=rk), intent(in), dimension(:,:):: grav

Gravitational acceleration [m/s^2]

private pure subroutine setSurfaceTension(self, surface_tension)

Sets the surface tension [N/m^2].

Arguments

Type IntentOptional AttributesName
class(domain_type), intent(inout) :: self

Domain instance

real(kind=rk), intent(in), dimension(:,:):: surface_tension

Surface tension [N/m^2]

private pure subroutine setAirDensity(self, air_density)

Sets the air density [kg/m^3].

Arguments

Type IntentOptional AttributesName
class(domain_type), intent(inout) :: self

Domain instance

real(kind=rk), intent(in), dimension(:,:):: air_density

Air density [kg/m^3]

private pure subroutine setWaterDensity(self, water_density)

Sets the water density [kg/m^3].

Arguments

Type IntentOptional AttributesName
class(domain_type), intent(inout) :: self

Domain instance

real(kind=rk), intent(in), dimension(:,:):: water_density

Water density [kg/m^3]

private pure subroutine setSpectrum1d(self, spectrum)

Sets the 1-d spectrum array. This procedure is overloaded by the generic procedure setSpectrum.

Arguments

Type IntentOptional AttributesName
class(domain_type), intent(inout) :: self

Domain instance

type(spectrum_type), intent(in), dimension(:):: spectrum

Input 1-d array of spectrum object instances

private pure subroutine setSpectrum2d(self, spectrum)

Sets the 2-d spectrum array. This procedure is overloaded by the generic procedure setSpectrum.

Arguments

Type IntentOptional AttributesName
class(domain_type), intent(inout) :: self

Domain instance

type(spectrum_type), intent(in), dimension(:,:):: spectrum

Input 2-d array of spectrum object instances

private pure subroutine setSpectrumArray1d1d(self, spectrum_array)

Sets the spectrum instances based on input spectrum array. This implementation is for omnidirectional spectrum in 1-d space (1d-1d)

Arguments

Type IntentOptional AttributesName
class(domain_type), intent(inout) :: self

Domain instance

real(kind=rk), intent(in), dimension(:,:):: spectrum_array

Spectrum array

private pure subroutine setSpectrumArray1d2d(self, spectrum_array)

Sets the spectrum instances based on input spectrum array. This implementation is for setting 1-d spectrum into 2-d physical space of 2-d spectrum into 1-d physical space.

Arguments

Type IntentOptional AttributesName
class(domain_type), intent(inout) :: self

Domain instance

real(kind=rk), intent(in), dimension(:,:,:):: spectrum_array

Spectrum array

private pure subroutine setSpectrumArray2d2d(self, spectrum_array)

Sets the spectrum instances based on input spectrum array. This implementation is for directional spectrum in 2-d space (2d-2d)

Arguments

Type IntentOptional AttributesName
class(domain_type), intent(inout) :: self

Domain instance

real(kind=rk), intent(in), dimension(:,:,:,:):: spectrum_array

Spectrum array

private subroutine writeJSON(self, filename, minify)

Writes a spectrum instance to a JSON file.

Arguments

Type IntentOptional AttributesName
class(domain_type), intent(in) :: self

domain instance

character(len=*), intent(in) :: filename

JSON file name

logical, intent(in) :: minify

Logical switch to minify the JSON file