modinit.F90 8.0 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209
  1. ! AGRIF (Adaptive Grid Refinement In Fortran)
  2. !
  3. ! Copyright (C) 2003 Laurent Debreu (Laurent.Debreu@imag.fr)
  4. ! Christophe Vouland (Christophe.Vouland@imag.fr)
  5. !
  6. ! This program is free software; you can redistribute it and/or modify
  7. ! it under the terms of the GNU General Public License as published by
  8. ! the Free Software Foundation; either version 2 of the License, or
  9. ! (at your option) any later version.
  10. !
  11. ! This program is distributed in the hope that it will be useful,
  12. ! but WITHOUT ANY WARRANTY; without even the implied warranty of
  13. ! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
  14. ! GNU General Public License for more details.
  15. !
  16. ! You should have received a copy of the GNU General Public License
  17. ! along with this program; if not, write to the Free Software
  18. ! Foundation, Inc., 59 Temple Place- Suite 330, Boston, MA 02111-1307, USA.
  19. !
  20. !
  21. !> Module Agrif_Init.
  22. !>
  23. !> Several operations on the variables of the current grid (creation, instanciation, ...)
  24. !! used during the creation of the grid hierarchy and during the time integration.
  25. !
  26. module Agrif_Init
  27. !
  28. use Agrif_Grids
  29. use Agrif_Link
  30. use Agrif_Mpp
  31. !
  32. implicit none
  33. !
  34. contains
  35. !
  36. !===================================================================================================
  37. ! subroutine Agrif_Allocation
  38. !
  39. !> Allocates the arrays containing the values of the variables of the current grd.
  40. !---------------------------------------------------------------------------------------------------
  41. subroutine Agrif_Allocation ( Agrif_Gr, procname )
  42. !---------------------------------------------------------------------------------------------------
  43. type(Agrif_Grid), pointer :: Agrif_Gr !< Pointer on the current grid
  44. procedure(alloc_proc), optional :: procname !< Allocation procedure (Default: Agrif_Allocationcalls)
  45. !
  46. if ( present(procname) ) then
  47. call procname(Agrif_Gr)
  48. else
  49. call Agrif_Allocationcalls(Agrif_Gr)
  50. endif
  51. Agrif_Gr % allocation_is_done = .true.
  52. !
  53. if ( Agrif_USE_ONLY_FIXED_GRIDS == 0 ) then
  54. !
  55. if ( Agrif_Probdim == 1 ) allocate( Agrif_Gr%tabpoint1D(Agrif_Gr%nb(1)+1) )
  56. if ( Agrif_Probdim == 2 ) allocate( Agrif_Gr%tabpoint2D(Agrif_Gr%nb(1)+1, &
  57. Agrif_Gr%nb(2)+1) )
  58. if ( Agrif_Probdim == 3 ) allocate( Agrif_Gr%tabpoint3D(Agrif_Gr%nb(1)+1, &
  59. Agrif_Gr%nb(2)+1, &
  60. Agrif_Gr%nb(3)+1) )
  61. endif
  62. !---------------------------------------------------------------------------------------------------
  63. end subroutine Agrif_Allocation
  64. !===================================================================================================
  65. !
  66. !===================================================================================================
  67. ! subroutine Agrif_Instance
  68. !
  69. !> Make the pointer Agrif_Types::Agrif_Curgrid point to Agrif_Gr
  70. !---------------------------------------------------------------------------------------------------
  71. subroutine Agrif_Instance ( Agrif_Gr )
  72. !---------------------------------------------------------------------------------------------------
  73. type(Agrif_Grid), pointer :: Agrif_Gr !< Pointer on the current grid
  74. !
  75. Agrif_Curgrid => Agrif_Gr
  76. Agrif_tabvars => Agrif_Curgrid % tabvars
  77. Agrif_tabvars_c => Agrif_Curgrid % tabvars_c
  78. Agrif_tabvars_r => Agrif_Curgrid % tabvars_r
  79. Agrif_tabvars_l => Agrif_Curgrid % tabvars_l
  80. Agrif_tabvars_i => Agrif_Curgrid % tabvars_i
  81. !
  82. #if defined AGRIF_MPI
  83. if ( Agrif_Gr % communicator /= -1 ) then
  84. call Agrif_MPI_switch_comm( Agrif_Gr % communicator )
  85. endif
  86. #endif
  87. !
  88. call Agrif_Get_numberofcells(Agrif_Gr)
  89. call Agrif_InitWorkSpace()
  90. !---------------------------------------------------------------------------------------------------
  91. end subroutine Agrif_Instance
  92. !===================================================================================================
  93. !
  94. !===================================================================================================
  95. ! subroutine Agrif_initialisations
  96. !---------------------------------------------------------------------------------------------------
  97. subroutine Agrif_initialisations ( Agrif_Gr )
  98. !---------------------------------------------------------------------------------------------------
  99. type(Agrif_Grid), pointer :: Agrif_Gr !< Pointer on the current grid
  100. !
  101. integer :: i
  102. type(Agrif_Variable), pointer :: var => NULL()
  103. type(Agrif_Variable_c), pointer :: var_c => NULL()
  104. type(Agrif_Variable_r), pointer :: var_r => NULL()
  105. type(Agrif_Variable_l), pointer :: var_l => NULL()
  106. type(Agrif_Variable_i), pointer :: var_i => NULL()
  107. !
  108. do i = 1,Agrif_NbVariables(0)
  109. !
  110. var => Agrif_Gr % tabvars(i)
  111. var % nbdim = 0
  112. !
  113. if (allocated(var%array1)) then
  114. var % nbdim = 1
  115. var % lb(1:1) = lbound(var%array1)
  116. var % ub(1:1) = ubound(var%array1)
  117. endif
  118. if (allocated(var%array2)) then
  119. var % nbdim = 2
  120. var % lb(1:2) = lbound(var%array2)
  121. var % ub(1:2) = ubound(var%array2)
  122. endif
  123. if (allocated(var%array3)) then
  124. var % nbdim = 3
  125. var % lb(1:3) = lbound(var%array3)
  126. var % ub(1:3) = ubound(var%array3)
  127. endif
  128. if (allocated(var%array4)) then
  129. var % nbdim = 4
  130. var % lb(1:4) = lbound(var%array4)
  131. var % ub(1:4) = ubound(var%array4)
  132. endif
  133. if (allocated(var%array5)) then
  134. var % nbdim = 5
  135. var % lb(1:5) = lbound(var%array5)
  136. var % ub(1:5) = ubound(var%array5)
  137. endif
  138. if (allocated(var%array6)) then
  139. var % nbdim = 6
  140. var % lb(1:6) = lbound(var%array6)
  141. var % ub(1:6) = ubound(var%array6)
  142. endif
  143. !
  144. if (allocated(var%darray1)) var % nbdim = 1
  145. if (allocated(var%darray2)) var % nbdim = 2
  146. if (allocated(var%darray3)) var % nbdim = 3
  147. if (allocated(var%darray4)) var % nbdim = 4
  148. if (allocated(var%darray5)) var % nbdim = 5
  149. if (allocated(var%darray6)) var % nbdim = 6
  150. !
  151. if (allocated(var%sarray1)) var % nbdim = 1
  152. if (allocated(var%sarray2)) var % nbdim = 2
  153. if (allocated(var%sarray3)) var % nbdim = 3
  154. if (allocated(var%sarray4)) var % nbdim = 4
  155. if (allocated(var%sarray5)) var % nbdim = 5
  156. if (allocated(var%sarray6)) var % nbdim = 6
  157. !
  158. enddo
  159. do i = 1,Agrif_NbVariables(1)
  160. !
  161. var_c => Agrif_Gr % tabvars_c(i)
  162. var_c % nbdim = 0
  163. !
  164. if (allocated(var_c%carray1)) var_c % nbdim = 1
  165. if (allocated(var_c%carray2)) var_c % nbdim = 2
  166. !
  167. enddo
  168. do i = 1,Agrif_NbVariables(2)
  169. !
  170. var_r => Agrif_Gr % tabvars_r(i)
  171. var_r % nbdim = 0
  172. !
  173. enddo
  174. do i = 1,Agrif_NbVariables(3)
  175. !
  176. var_l => Agrif_Gr % tabvars_l(i)
  177. var_l % nbdim = 0
  178. !
  179. if (allocated(var_l%larray1)) var_l % nbdim = 1
  180. if (allocated(var_l%larray2)) var_l % nbdim = 2
  181. if (allocated(var_l%larray3)) var_l % nbdim = 3
  182. if (allocated(var_l%larray4)) var_l % nbdim = 4
  183. if (allocated(var_l%larray5)) var_l % nbdim = 5
  184. if (allocated(var_l%larray6)) var_l % nbdim = 6
  185. !
  186. enddo
  187. do i = 1,Agrif_NbVariables(4)
  188. !
  189. var_i => Agrif_Gr % tabvars_i(i)
  190. var_i % nbdim = 0
  191. !
  192. if (allocated(var_i%iarray1)) var_i % nbdim = 1
  193. if (allocated(var_i%iarray2)) var_i % nbdim = 2
  194. if (allocated(var_i%iarray3)) var_i % nbdim = 3
  195. if (allocated(var_i%iarray4)) var_i % nbdim = 4
  196. if (allocated(var_i%iarray5)) var_i % nbdim = 5
  197. if (allocated(var_i%iarray6)) var_i % nbdim = 6
  198. !
  199. enddo
  200. !---------------------------------------------------------------------------------------------------
  201. end subroutine Agrif_initialisations
  202. !===================================================================================================
  203. !
  204. end module Agrif_Init