SABER
bump_interpolation_mod Module Reference

Bump Interpolation module. More...

Data Types

type  bump_interpolator
 

Functions/Subroutines

subroutine bint_init (self, config, comm, in_funcspace, out_funcspace, masks)
 Linked list implementation. More...
 
subroutine bint_driver (self, mpl, rng, nam, geom)
 Initialize bump to perform interpolation. More...
 
subroutine bint_apply (self, infields, outfields)
 Apply interpolation. More...
 
subroutine apply_interp (self, infield, outfield)
 Subroutine: apply_interp Purpose: low-level routine to apply the interpolation to a single field on a single level. More...
 
subroutine bint_apply_ad (self, fields_outgrid, fields_ingrid)
 Apply interpolator operator adjoint. More...
 
subroutine apply_interp_ad (self, fld_outgrid, fld_ingrid)
 Subroutine: apply_interp_ad Purpose: low-level routine to apply the adjoint of the interpolation operator to a single field on a single level. More...
 
subroutine bint_deallocate_outgrid (self)
 Release memory (partial) by deallocating output grid. More...
 
subroutine bint_delete (self)
 Release all memory. More...
 
subroutine dummy (bump)
 Dummy finalization. More...
 

Variables

integer, parameter max_string = 1024
 
type(registry_t), public bump_interpolator_registry
 Registry for bump_interpolator objects. More...
 

Detailed Description

Bump Interpolation module.

Date
Jan, 2020: Created by M. Miesch (JCSDA/UCAR) based on previously exising type_obsop.F90 written by Benjamin Menetrier (JCSDA/UCAR, CERFACS, METEO-FRANCE, IRIT)

Function/Subroutine Documentation

◆ apply_interp()

subroutine bump_interpolation_mod::apply_interp ( class(bump_interpolator), intent(inout)  self,
real(kind_real), dimension(self%bump%geom%nc0a,self%nlev), intent(in)  infield,
real(kind_real), dimension(self%nout_local,self%nlev), intent(out)  outfield 
)
private

Subroutine: apply_interp Purpose: low-level routine to apply the interpolation to a single field on a single level.

Parameters
[in]infieldinput field
[out]outfieldoutput field

Definition at line 480 of file bump_interpolation_mod.F90.

◆ apply_interp_ad()

subroutine bump_interpolation_mod::apply_interp_ad ( class(bump_interpolator), intent(inout)  self,
real(kind_real), dimension(self%nout_local,self%nlev), intent(in)  fld_outgrid,
real(kind_real), dimension(self%bump%geom%nc0a,self%nlev), intent(out)  fld_ingrid 
)
private

Subroutine: apply_interp_ad Purpose: low-level routine to apply the adjoint of the interpolation operator to a single field on a single level.

Parameters
[in]mplbump MPI data
[in]geombump geometry data
[in]fld_outgridfield on output grid
[out]fld_ingridfield on input grid

Definition at line 609 of file bump_interpolation_mod.F90.

◆ bint_apply()

subroutine bump_interpolation_mod::bint_apply ( class(bump_interpolator), intent(inout)  self,
type(fieldset_type), intent(in)  infields,
type(fieldset_type), intent(inout)  outfields 
)
private

Apply interpolation.

Parameters
[in]instate= input fields represented as an atlas fieldset, created from a functionspace
[out]outstate= output fields represented as an atlas fieldset. If the fields that constitudte the fieldset are not already allocated by the caller, then they will be created and allocated by this method. So, the user can optionally pass this routine an empty output fieldset.

Definition at line 399 of file bump_interpolation_mod.F90.

◆ bint_apply_ad()

subroutine bump_interpolation_mod::bint_apply_ad ( class(bump_interpolator), intent(inout)  self,
type(fieldset_type), intent(in)  fields_outgrid,
type(fieldset_type), intent(inout)  fields_ingrid 
)
private

Apply interpolator operator adjoint.

Parameters
[in]fields_outgrid= These are the fields on the second grid, i.e. the target grid of the original interpolation. For the adjoint, these fields are treated as an input.
[in,out]fields_ingrid= These are the fields defined on the first grid, i.e. the source grid of the original interpolation. For the adjoint, these are treated as an output. The caller can optionally pass this arguement as an empty fieldset and the routine will create and allocate each component of the fieldset. Or, if the field components of the fieldset are already allocated by the caller, then this routine will merely replace the field values with the result of the computation.

Definition at line 523 of file bump_interpolation_mod.F90.

◆ bint_deallocate_outgrid()

subroutine bump_interpolation_mod::bint_deallocate_outgrid ( class(bump_interpolator), intent(inout)  self)
private

Release memory (partial) by deallocating output grid.

Definition at line 637 of file bump_interpolation_mod.F90.

◆ bint_delete()

subroutine bump_interpolation_mod::bint_delete ( class(bump_interpolator), intent(inout)  self)
private

Release all memory.

Definition at line 649 of file bump_interpolation_mod.F90.

◆ bint_driver()

subroutine bump_interpolation_mod::bint_driver ( class(bump_interpolator), intent(inout)  self,
type(mpl_type), intent(inout)  mpl,
type(rng_type), intent(inout)  rng,
type(nam_type), intent(in)  nam,
type(geom_type), intent(in)  geom 
)

Initialize bump to perform interpolation.

Parameters
[in,out]mpl< MPI data
[in,out]rng< Random number generator
[in]nam< Namelist
[in]geom< Geometry

Definition at line 267 of file bump_interpolation_mod.F90.

◆ bint_init()

subroutine bump_interpolation_mod::bint_init ( class(bump_interpolator), intent(inout)  self,
class(fckit_configuration), intent(in)  config,
type(fckit_mpi_comm), intent(in)  comm,
class(atlas_functionspace), intent(in)  in_funcspace,
class(atlas_functionspace), intent(in)  out_funcspace,
type(fieldset_type), intent(in), optional  masks 
)
private

Linked list implementation.

Initialize interpolation object

The input and output fields are atlas_FieldSet objects that are assumed to be created from atlas functionspaces. So, they have the grid and mesh information built in.

Parameters
[in]config= configuration
[in]in_afs= This is the input grid, rendered as an atlas functionspace, so it includes information about the mesh (parallel connectivity) as well.
[in]out_afs= This is the output grid, also represented as an altas_functionspace.
[in]masks= This contains metadata needed for the interpolation. It is rendered as an atlas FieldSet with the following named fields. Each of these named fields is optional; if omitted default values will be provided
  • area : cell area
  • vunit : vertical unit
  • gmask : geometry mask
  • smask : sampling mask

Definition at line 133 of file bump_interpolation_mod.F90.

◆ dummy()

subroutine bump_interpolation_mod::dummy ( type(bump_interpolator), intent(inout)  bump)
private

Dummy finalization.

Definition at line 664 of file bump_interpolation_mod.F90.

Here is the caller graph for this function:

Variable Documentation

◆ bump_interpolator_registry

type(registry_t), public bump_interpolation_mod::bump_interpolator_registry

Registry for bump_interpolator objects.

Linked list interface - defines registry_t type Global registry

Definition at line 99 of file bump_interpolation_mod.F90.

◆ max_string

integer, parameter bump_interpolation_mod::max_string = 1024
private

Definition at line 43 of file bump_interpolation_mod.F90.