new_NS_Basis.F90 Source File


This file depends on

sourcefile~~new_ns_basis.f90~~EfferentGraph sourcefile~new_ns_basis.f90 new_NS_Basis.F90 sourcefile~base_base.f90 Base_Base.F90 sourcefile~new_ns_basis.f90->sourcefile~base_base.f90 sourcefile~vectorbasis.f90 VectorBasis.F90 sourcefile~new_ns_basis.f90->sourcefile~vectorbasis.f90 sourcefile~constants.f90 Constants.F90 sourcefile~base_base.f90->sourcefile~constants.f90 sourcefile~keywordenforcer.f90 KeywordEnforcer.F90 sourcefile~base_base.f90->sourcefile~keywordenforcer.f90 sourcefile~mapl_range.f90 MAPL_Range.F90 sourcefile~base_base.f90->sourcefile~mapl_range.f90 sourcefile~maplgrid.f90 MaplGrid.F90 sourcefile~base_base.f90->sourcefile~maplgrid.f90 sourcefile~errorhandling.f90 ErrorHandling.F90 sourcefile~vectorbasis.f90->sourcefile~errorhandling.f90 sourcefile~fieldblas.f90 FieldBLAS.F90 sourcefile~vectorbasis.f90->sourcefile~fieldblas.f90 sourcefile~fieldpointerutilities.f90 FieldPointerUtilities.F90 sourcefile~vectorbasis.f90->sourcefile~fieldpointerutilities.f90 sourcefile~internalconstants.f90 InternalConstants.F90 sourcefile~constants.f90->sourcefile~internalconstants.f90 sourcefile~mathconstants.f90 MathConstants.F90 sourcefile~constants.f90->sourcefile~mathconstants.f90 sourcefile~physicalconstants.f90 PhysicalConstants.F90 sourcefile~constants.f90->sourcefile~physicalconstants.f90 sourcefile~mapl_throw.f90 MAPL_Throw.F90 sourcefile~errorhandling.f90->sourcefile~mapl_throw.f90 sourcefile~fieldblas.f90->sourcefile~fieldpointerutilities.f90 sourcefile~mapl_exceptionhandling.f90 MAPL_ExceptionHandling.F90 sourcefile~fieldblas.f90->sourcefile~mapl_exceptionhandling.f90 sourcefile~fieldpointerutilities.f90->sourcefile~mapl_exceptionhandling.f90 sourcefile~mapl_range.f90->sourcefile~mapl_exceptionhandling.f90 sourcefile~maplgrid.f90->sourcefile~constants.f90 sourcefile~maplgrid.f90->sourcefile~errorhandling.f90 sourcefile~maplgrid.f90->sourcefile~keywordenforcer.f90 sourcefile~mapl_sort.f90 MAPL_Sort.F90 sourcefile~maplgrid.f90->sourcefile~mapl_sort.f90 sourcefile~pflogger_stub.f90 pflogger_stub.F90 sourcefile~maplgrid.f90->sourcefile~pflogger_stub.f90 sourcefile~mapl_exceptionhandling.f90->sourcefile~errorhandling.f90 sourcefile~mapl_exceptionhandling.f90->sourcefile~mapl_throw.f90 sourcefile~mapl_sort.f90->sourcefile~mapl_exceptionhandling.f90 sourcefile~pfl_keywordenforcer.f90 PFL_KeywordEnforcer.F90 sourcefile~pflogger_stub.f90->sourcefile~pfl_keywordenforcer.f90 sourcefile~wraparray.f90 WrapArray.F90 sourcefile~pflogger_stub.f90->sourcefile~wraparray.f90 sourcefile~physicalconstants.f90->sourcefile~mathconstants.f90

Source Code

#include "MAPL_ErrLog.h"

submodule (mapl3g_VectorBasis) new_NS_Basis_smod
   use mapl_base, only: MAPL_GridGetCorners
contains


   module function new_NS_Basis(geom, rc) result(basis)
      type(VectorBasis) :: basis
      type(ESMF_Geom), intent(inout) :: geom
      integer, optional, intent(out) :: rc

      integer :: status
      real(kind=ESMF_KIND_R8), pointer :: longitudes(:)
      real(kind=ESMF_KIND_R8), pointer :: latitudes(:)

      call create_fields(basis%elements, geom, _RC)
      call MAPL_GeomGetCoords(geom, longitudes, latitudes, _RC)
      call fill_fields(basis, longitudes, latitudes, _RC)

      _RETURN(ESMF_SUCCESS)

   contains

      subroutine fill_fields(basis, longitudes, latitudes, rc)
         type(VectorBasis), intent(inout) :: basis
         real(kind=ESMF_KIND_R8), intent(in) :: longitudes(:)
         real(kind=ESMF_KIND_R8), intent(in) :: latitudes(:)
         integer, optional, intent(out) :: rc

         integer :: status
         type(Ptr_1d) :: x(NI,NJ)
         integer :: i, j, n
         real(kind=ESMF_KIND_R8) :: local_basis(NI,NJ)

         do j = 1, NJ
            do i = 1, NI
               call assign_fptr(basis%elements(i,j), x(i,j)%ptr, _RC)
            end do
         end do

         do n = 1, size(x(1,1)%ptr)
            local_basis = fill_element(longitudes(i), latitudes(i))

            do j = 1, NJ
               do i = 1, NI
                  x(i,j)%ptr(n) = local_basis(i,j)
               end do
            end do

         end do

         _RETURN(ESMF_SUCCESS)
      end subroutine fill_fields

      pure function fill_element(longitude, latitude) result(x)
         real(kind=ESMF_KIND_R8) :: x(NI,NJ)
         real(kind=ESMF_KIND_R8), intent(in) :: longitude
         real(kind=ESMF_KIND_R8), intent(in) :: latitude

         x(:,1) = [ -sin(longitude), cos(longitude), 0._ESMF_KIND_R8 ]
         x(:,2) = [ -sin(latitude)*cos(longitude), -sin(latitude)*sin(longitude), cos(latitude) ]

      end function fill_element

   end function new_NS_Basis

end submodule new_NS_Basis_smod