RegridderManager.F90 Source File


This file depends on

sourcefile~~regriddermanager.f90~~EfferentGraph sourcefile~regriddermanager.f90 RegridderManager.F90 sourcefile~errorhandling.f90 ErrorHandling.F90 sourcefile~regriddermanager.f90->sourcefile~errorhandling.f90 sourcefile~esmfregridderfactory.f90 EsmfRegridderFactory.F90 sourcefile~regriddermanager.f90->sourcefile~esmfregridderfactory.f90 sourcefile~nullregridder.f90 NullRegridder.F90 sourcefile~regriddermanager.f90->sourcefile~nullregridder.f90 sourcefile~regridderfactory.f90 RegridderFactory.F90 sourcefile~regriddermanager.f90->sourcefile~regridderfactory.f90 sourcefile~regridderfactoryvector.f90 RegridderFactoryVector.F90 sourcefile~regriddermanager.f90->sourcefile~regridderfactoryvector.f90 sourcefile~regridderspec.f90~2 RegridderSpec.F90 sourcefile~regriddermanager.f90->sourcefile~regridderspec.f90~2 sourcefile~regridderspecvector.f90 RegridderSpecVector.F90 sourcefile~regriddermanager.f90->sourcefile~regridderspecvector.f90 sourcefile~regriddervector.f90 RegridderVector.F90 sourcefile~regriddermanager.f90->sourcefile~regriddervector.f90 sourcefile~mapl_throw.f90 MAPL_Throw.F90 sourcefile~errorhandling.f90->sourcefile~mapl_throw.f90 sourcefile~esmfregridderfactory.f90->sourcefile~errorhandling.f90 sourcefile~esmfregridderfactory.f90->sourcefile~nullregridder.f90 sourcefile~esmfregridderfactory.f90->sourcefile~regridderfactory.f90 sourcefile~esmfregridderfactory.f90->sourcefile~regridderspec.f90~2 sourcefile~regridderparam.f90 RegridderParam.F90 sourcefile~esmfregridderfactory.f90->sourcefile~regridderparam.f90 sourcefile~routehandlemanager.f90 RoutehandleManager.F90 sourcefile~esmfregridderfactory.f90->sourcefile~routehandlemanager.f90 sourcefile~routehandleparam.f90 RoutehandleParam.F90 sourcefile~esmfregridderfactory.f90->sourcefile~routehandleparam.f90 sourcefile~nullregridder.f90->sourcefile~errorhandling.f90 sourcefile~regridderfactoryvector.f90->sourcefile~regridderfactory.f90 sourcefile~geom_mgr.f90 geom_mgr.F90 sourcefile~regridderspec.f90~2->sourcefile~geom_mgr.f90 sourcefile~regridderspec.f90~2->sourcefile~regridderparam.f90 sourcefile~regridderspecvector.f90->sourcefile~regridderspec.f90~2 sourcefile~geomspec.f90 GeomSpec.F90 sourcefile~geom_mgr.f90->sourcefile~geomspec.f90 sourcefile~geomutilities.f90 GeomUtilities.F90 sourcefile~geom_mgr.f90->sourcefile~geomutilities.f90 sourcefile~routehandlemanager.f90->sourcefile~errorhandling.f90 sourcefile~routehandlespec.f90 RoutehandleSpec.F90 sourcefile~routehandlemanager.f90->sourcefile~routehandlespec.f90 sourcefile~routehandlespecvector.f90 RoutehandleSpecVector.F90 sourcefile~routehandlemanager.f90->sourcefile~routehandlespecvector.f90 sourcefile~routehandlevector.f90 RoutehandleVector.F90 sourcefile~routehandlemanager.f90->sourcefile~routehandlevector.f90 sourcefile~routehandleparam.f90->sourcefile~errorhandling.f90 sourcefile~routehandleparam.f90->sourcefile~geom_mgr.f90 sourcefile~geomutilities.f90->sourcefile~errorhandling.f90 sourcefile~routehandlespec.f90->sourcefile~errorhandling.f90 sourcefile~routehandlespec.f90->sourcefile~geom_mgr.f90 sourcefile~routehandlespec.f90->sourcefile~routehandleparam.f90 sourcefile~routehandlespecvector.f90->sourcefile~routehandlespec.f90

Files dependent on this one

sourcefile~~regriddermanager.f90~~AfferentGraph sourcefile~regriddermanager.f90 RegridderManager.F90 sourcefile~regridder_mgr.f90 regridder_mgr.F90 sourcefile~regridder_mgr.f90->sourcefile~regriddermanager.f90 sourcefile~regridaction.f90 RegridAction.F90 sourcefile~regridaction.f90->sourcefile~regridder_mgr.f90 sourcefile~test_regriddermanager.pf Test_RegridderManager.pf sourcefile~test_regriddermanager.pf->sourcefile~regridder_mgr.f90 sourcefile~test_routehandlemanager.pf Test_RouteHandleManager.pf sourcefile~test_routehandlemanager.pf->sourcefile~regridder_mgr.f90 sourcefile~fieldspec.f90~2 FieldSpec.F90 sourcefile~fieldspec.f90~2->sourcefile~regridaction.f90 sourcefile~bracketspec.f90 BracketSpec.F90 sourcefile~bracketspec.f90->sourcefile~fieldspec.f90~2 sourcefile~make_itemspec.f90 make_itemSpec.F90 sourcefile~make_itemspec.f90->sourcefile~fieldspec.f90~2 sourcefile~modelverticalgrid.f90 ModelVerticalGrid.F90 sourcefile~modelverticalgrid.f90->sourcefile~fieldspec.f90~2 sourcefile~test_addfieldspec.pf Test_AddFieldSpec.pf sourcefile~test_addfieldspec.pf->sourcefile~fieldspec.f90~2 sourcefile~test_bracketspec.pf Test_BracketSpec.pf sourcefile~test_bracketspec.pf->sourcefile~fieldspec.f90~2 sourcefile~test_fieldinfo.pf Test_FieldInfo.pf sourcefile~test_fieldinfo.pf->sourcefile~fieldspec.f90~2 sourcefile~test_fieldspec.pf Test_FieldSpec.pf sourcefile~test_fieldspec.pf->sourcefile~fieldspec.f90~2

Source Code

#include "MAPL_Generic.h"
module mapl3g_RegridderManager

   use mapl3g_RegridderSpec
   use mapl3g_Regridder
   use mapl3g_NullRegridder
   use mapl3g_RegridderFactory

   use mapl3g_RegridderFactoryVector
   use mapl3g_RegridderSpecVector
   use mapl3g_RegridderVector
   use mapl3g_EsmfRegridderFactory

   use mapl_ErrorHandlingMod
   implicit none
   private

   public :: RegridderManager
   public :: regridder_manager ! singleton
   public :: get_regridder_manager

   type :: RegridderManager
      private
      type(RegridderFactoryVector) :: factories
      ! Next two vectors grow together
      type(RegridderSpecVector) :: specs
      type(RegridderVector) :: regridders
   contains
      procedure :: get_regridder
      procedure :: add_factory
      procedure :: make_regridder
      procedure :: add_regridder
      procedure :: delete_regridder
   end type RegridderManager

   interface RegridderManager
      procedure new_RegridderManager
   end interface RegridderManager

   type(RegridderManager), target, protected :: regridder_manager

contains

   function new_RegridderManager() result(mgr)
      type(RegridderManager) :: mgr

      ! Load default factories

      call mgr%add_factory(EsmfRegridderFactory())
!!$      call mgr%add_factory(horzHorzFluxRegridderFactory())

   end function new_RegridderManager


   ! TODO - do we need an RC here for duplicate name?
   subroutine add_factory(this, factory)
      class(RegridderManager), intent(inout) :: this
      class(RegridderFactory), intent(in) :: factory
      call this%factories%push_back(factory)
   end subroutine add_factory


   subroutine add_regridder(this, spec, regriddr)
      class(RegridderManager), intent(inout) :: this
      class(RegridderSpec), intent(in) :: spec
      class(Regridder), intent(in) :: regriddr

      call this%specs%push_back(spec)
      call this%regridders%push_back(regriddr)
     
   end subroutine add_regridder

   subroutine delete_regridder(this, spec, rc)
      class(RegridderManager), target, intent(inout) :: this
      class(RegridderSpec), intent(in) :: spec
      integer, optional, intent(out) :: rc

      integer :: status
      type(RegridderSpecVectorIterator) :: spec_iter
      type(RegridderVectorIterator) :: regridder_iter

      associate (specs => this%specs, regridders => this%regridders)
        associate (b => specs%begin(), e => specs%end())

          spec_iter = find(b, e, spec)
          _ASSERT(spec_iter /= e, 'spec not found in RegridderManager.')

          regridder_iter = regridders%begin() + (spec_iter - b)
          regridder_iter = regridders%erase(regridder_iter)

          spec_iter = specs%erase(spec_iter)

        end associate
      end associate

      _RETURN(_SUCCESS)
   end subroutine delete_regridder

   function get_regridder(this, spec, rc) result(regriddr)
      class(Regridder), pointer :: regriddr
      class(RegridderManager), target, intent(inout) :: this
      class(RegridderSpec), intent(in) :: spec
      integer, optional, intent(out) :: rc

      integer :: status
      class(Regridder), allocatable :: tmp_regridder

      associate (b => this%specs%begin(), e => this%specs%end())
        associate (iter => find(b, e, spec))
          if (iter /= e) then
             regriddr => this%regridders%of((iter-b+1))
             _RETURN(_SUCCESS)
          end if

          tmp_regridder = this%make_regridder(spec, _RC)
          call this%add_regridder(spec, tmp_regridder)
          regriddr => this%regridders%back()

        end associate
      end associate

      _RETURN(_SUCCESS)
   end function get_regridder

   function make_regridder(this, spec, rc) result(regriddr)
      class(Regridder), allocatable :: regriddr
      class(RegridderManager), target, intent(in) :: this
      class(RegridderSpec), intent(in) :: spec
      integer, optional, intent(out) :: rc

      integer :: status
      integer :: i
      class(RegridderFactory), pointer :: factory

      regriddr = NULL_REGRIDDER
      do i = 1, this%factories%size()
         factory => this%factories%of(i)
         if (factory%supports(spec%get_param())) then
            deallocate(regriddr) ! workaround for gfortran 12.3
            regriddr = factory%make_regridder(spec, _RC)
            _RETURN(_SUCCESS)
         end if
      end do

      _FAIL('No factory found to make regridder for spec.')
   end function make_regridder

   function get_regridder_manager() result(regridder_mgr)
      type(RegridderManager), pointer :: regridder_mgr
      logical :: init = .false.

      if (.not. init) then
         regridder_manager = RegridderManager()
         init = .true.
      end if

      regridder_mgr => regridder_manager
         

   end function get_regridder_manager

end module mapl3g_RegridderManager