modvariables.F90 5.8 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146
  1. module Agrif_Variables
  2. !
  3. use Agrif_CurgridFunctions
  4. !
  5. implicit none
  6. !
  7. contains
  8. !
  9. !===================================================================================================
  10. ! subroutine Agrif_Declare_Variable
  11. !
  12. !> Declare a new variable profile
  13. !---------------------------------------------------------------------------------------------------
  14. subroutine Agrif_Declare_Variable ( posvar, firstpoint, raf, lb, ub, varid, torestore )
  15. !---------------------------------------------------------------------------------------------------
  16. integer, dimension(:), intent(in) :: posvar !< position of the variable on the cell
  17. !! (1 for the border of the edge, 2 for the center)
  18. integer, dimension(:), intent(in) :: firstpoint !< index of the first point in the real domain
  19. character(1), dimension(:), intent(in) :: raf !< Array indicating the type of dimension (space or not)
  20. !! for each of them
  21. integer, dimension(:), intent(in) :: lb !< Lower bounds of the array
  22. integer, dimension(:), intent(in) :: ub !< Upper bounds of the array
  23. integer, intent(out) :: varid !< Id number of the newly created variable
  24. logical, optional, intent(in) :: torestore !< Indicates if the array restore is used
  25. !---------------------------------------------------------------------------------------------------
  26. type(Agrif_Variables_List), pointer :: new_varlist
  27. type(Agrif_Variable), pointer :: var
  28. integer :: nbdim, i
  29. logical :: restore
  30. restore = .FALSE.
  31. if ( Agrif_Mygrid % ngridstep /= 0 ) then
  32. if (present(torestore)) restore = torestore
  33. endif
  34. !
  35. nbdim = SIZE(posvar)
  36. !
  37. allocate(new_varlist)
  38. allocate(new_varlist % var)
  39. var => new_varlist % var
  40. allocate(var % posvar(nbdim))
  41. allocate(var % interptab(nbdim))
  42. allocate(var % coords(nbdim))
  43. !
  44. var % nbdim = nbdim
  45. var % interptab = raf(1:nbdim)
  46. var % posvar = posvar(1:nbdim)
  47. var % point(1:nbdim) = firstpoint(1:nbdim)
  48. var % restore = restore
  49. !
  50. do i = 1,nbdim
  51. select case( raf(i) )
  52. case('x') ; var % coords(i) = 1
  53. case('y') ; var % coords(i) = 2
  54. case('z') ; var % coords(i) = 3
  55. case('N') ; var % coords(i) = 0
  56. case default ; var % coords(i) = 0
  57. end select
  58. enddo
  59. !
  60. var % lb(1:nbdim) = lb(1:nbdim)
  61. var % ub(1:nbdim) = ub(1:nbdim)
  62. if ( restore ) then
  63. select case(nbdim)
  64. case(1)
  65. allocate(var % Restore1D(lb(1):ub(1)))
  66. var % Restore1D = 0
  67. case(2)
  68. allocate(var % Restore2D(lb(1):ub(1), &
  69. lb(2):ub(2)))
  70. var % Restore2D = 0
  71. case(3)
  72. allocate(var % Restore3D(lb(1):ub(1), &
  73. lb(2):ub(2), &
  74. lb(3):ub(3)))
  75. var % Restore3D = 0
  76. case(4)
  77. allocate(var % Restore4D(lb(1):ub(1), &
  78. lb(2):ub(2), &
  79. lb(3):ub(3), &
  80. lb(4):ub(4)))
  81. var % Restore4D = 0
  82. case(5)
  83. allocate(var % Restore5D(lb(1):ub(1), &
  84. lb(2):ub(2), &
  85. lb(3):ub(3), &
  86. lb(4):ub(4), &
  87. lb(5):ub(5)))
  88. var % Restore5D = 0
  89. end select
  90. endif
  91. new_varlist % next => Agrif_Curgrid % variables
  92. Agrif_Curgrid % variables => new_varlist
  93. Agrif_Curgrid % Nbvariables = Agrif_Curgrid % Nbvariables + 1
  94. varid = -Agrif_Curgrid % Nbvariables
  95. var % parent_var => Agrif_Search_Variable(Agrif_Curgrid % parent, Agrif_Curgrid % nbvariables)
  96. var % root_var => Agrif_Search_Variable(Agrif_Mygrid, Agrif_Curgrid % nbvariables)
  97. !---------------------------------------------------------------------------------------------------
  98. end subroutine Agrif_Declare_Variable
  99. !===================================================================================================
  100. !
  101. !===================================================================================================
  102. ! function Agrif_Search_Variable
  103. !
  104. !> Returns a pointer to the variable varid for the grid grid.
  105. !---------------------------------------------------------------------------------------------------
  106. function Agrif_Search_Variable ( grid, varid ) result(outvar)
  107. !---------------------------------------------------------------------------------------------------
  108. type(Agrif_Grid), pointer :: grid !< Pointer on the current grid.
  109. integer, intent(in) :: varid !< ID number of the variable we are looking for.
  110. !
  111. type(Agrif_Variable), pointer :: outvar
  112. type(Agrif_Variables_List), pointer :: parcours
  113. integer :: nb, varidinv
  114. !
  115. if ( .not.associated(grid) ) then
  116. outvar => NULL()
  117. return
  118. endif
  119. !
  120. parcours => grid % variables
  121. if (.not. associated(parcours)) then ! can occur on the grand mother grid
  122. outvar => NULL() ! during the first call by agrif_mygrid
  123. return
  124. endif
  125. varidinv = 1 + grid % nbvariables - varid
  126. do nb = 1,varidinv-1
  127. parcours => parcours % next
  128. enddo
  129. outvar => parcours % var
  130. !---------------------------------------------------------------------------------------------------
  131. end function Agrif_Search_variable
  132. !===================================================================================================
  133. !
  134. end module Agrif_Variables