modarrays.F90 39 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622623624625626627628629630631632633634635636637638639640641642643644645646647648649650651652653654655656657658659660661662663664665666667668669670671672673674675676677678679680681682683684685686687688689690691692693694695696697698699700701702703704705706707708709710711712713714715716717718719720721722723724725726727728729730731732733734735736737738739740741742743744745746747748749750751752753754755756757758759760761762763764765766767768769770771772773774775776777778779780781782783784785786787788789790791792793794795796797798799800801802803804805806807808809810811812813814815816817818819820821822823824825826
  1. !
  2. ! $Id: modarrays.F90 4779 2014-09-19 14:21:37Z rblod $
  3. !
  4. ! AGRIF (Adaptive Grid Refinement In Fortran)
  5. !
  6. ! Copyright (C) 2003 Laurent Debreu (Laurent.Debreu@imag.fr)
  7. ! Christophe Vouland (Christophe.Vouland@imag.fr)
  8. !
  9. ! This program is free software; you can redistribute it and/or modify
  10. ! it under the terms of the GNU General Public License as published by
  11. ! the Free Software Foundation; either version 2 of the License, or
  12. ! (at your option) any later version.
  13. !
  14. ! This program is distributed in the hope that it will be useful,
  15. ! but WITHOUT ANY WARRANTY; without even the implied warranty of
  16. ! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
  17. ! GNU General Public License for more details.
  18. !
  19. ! You should have received a copy of the GNU General Public License
  20. ! along with this program; if not, write to the Free Software
  21. ! Foundation, Inc., 59 Temple Place- Suite 330, Boston, MA 02111-1307, USA.
  22. !
  23. !
  24. !> Module Agrif_Arrays
  25. !
  26. module Agrif_Arrays
  27. !
  28. use Agrif_Types
  29. use Agrif_Grids
  30. !
  31. implicit none
  32. !
  33. #if defined AGRIF_MPI
  34. interface
  35. subroutine Agrif_InvLoc ( indloc, proc_id, dir, indglob )
  36. integer, intent(in) :: indloc !< local index
  37. integer, intent(in) :: proc_id !< rank of the proc calling this function
  38. integer, intent(in) :: dir !< direction of the index
  39. integer, intent(out) :: indglob !< global index
  40. end subroutine Agrif_InvLoc
  41. end interface
  42. private :: Agrif_InvLoc
  43. #endif
  44. !
  45. contains
  46. !
  47. !===================================================================================================
  48. ! subroutine Agrif_Childbounds
  49. !
  50. !> Computes the global indices of the child grid
  51. !---------------------------------------------------------------------------------------------------
  52. subroutine Agrif_Childbounds ( nbdim, &
  53. lb_var, ub_var, &
  54. lb_tab, ub_tab, &
  55. proc_id, &
  56. coords, &
  57. lb_tab_true, ub_tab_true, memberin )
  58. !---------------------------------------------------------------------------------------------------
  59. integer, intent(in) :: nbdim !< Number of dimensions
  60. integer, dimension(nbdim), intent(in) :: lb_var !< Local lower boundary on the current processor
  61. integer, dimension(nbdim), intent(in) :: ub_var !< Local upper boundary on the current processor
  62. integer, dimension(nbdim), intent(in) :: lb_tab !< Global lower boundary of the variable
  63. integer, dimension(nbdim), intent(in) :: ub_tab !< Global upper boundary of the variable
  64. integer, intent(in) :: proc_id !< Current processor
  65. integer, dimension(nbdim), intent(in) :: coords
  66. integer, dimension(nbdim), intent(out) :: lb_tab_true !< Global value of lb_var on the current processor
  67. integer, dimension(nbdim), intent(out) :: ub_tab_true !< Global value of ub_var on the current processor
  68. logical, intent(out) :: memberin
  69. !
  70. integer :: i, coord_i
  71. integer :: lb_glob_index, ub_glob_index ! Lower and upper global indices
  72. !
  73. do i = 1, nbdim
  74. !
  75. coord_i = coords(i)
  76. !
  77. #if defined AGRIF_MPI
  78. call Agrif_InvLoc( lb_var(i), proc_id, coord_i, lb_glob_index )
  79. call Agrif_InvLoc( ub_var(i), proc_id, coord_i, ub_glob_index )
  80. #else
  81. lb_glob_index = lb_var(i)
  82. ub_glob_index = ub_var(i)
  83. #endif
  84. lb_tab_true(i) = max(lb_tab(i), lb_glob_index)
  85. ub_tab_true(i) = min(ub_tab(i), ub_glob_index)
  86. enddo
  87. !
  88. memberin = .true.
  89. do i = 1,nbdim
  90. if (ub_tab_true(i) < lb_tab_true(i)) then
  91. memberin = .false.
  92. exit
  93. endif
  94. enddo
  95. !---------------------------------------------------------------------------------------------------
  96. end subroutine Agrif_Childbounds
  97. !===================================================================================================
  98. !
  99. !===================================================================================================
  100. subroutine Agrif_get_var_global_bounds( var, lubglob, nbdim )
  101. !---------------------------------------------------------------------------------------------------
  102. type(Agrif_Variable), intent(in) :: var
  103. integer, dimension(nbdim,2), intent(out) :: lubglob
  104. integer, intent(in) :: nbdim
  105. !
  106. #if defined AGRIF_MPI
  107. include 'mpif.h'
  108. integer, dimension(nbdim) :: lb, ub
  109. integer, dimension(nbdim,2) :: iminmaxg
  110. integer :: i, code, coord_i
  111. #endif
  112. !
  113. #if !defined AGRIF_MPI
  114. call Agrif_get_var_bounds_array(var, lubglob(:,1), lubglob(:,2), nbdim)
  115. #else
  116. call Agrif_get_var_bounds_array(var, lb, ub, nbdim)
  117. do i = 1,nbdim
  118. coord_i = var % root_var % coords(i)
  119. call Agrif_InvLoc( lb(i), Agrif_Procrank, coord_i, iminmaxg(i,1) )
  120. call Agrif_InvLoc( ub(i), Agrif_Procrank, coord_i, iminmaxg(i,2) )
  121. enddo
  122. !
  123. iminmaxg(1:nbdim,2) = - iminmaxg(1:nbdim,2)
  124. call MPI_ALLREDUCE(iminmaxg, lubglob, 2*nbdim, MPI_INTEGER, MPI_MIN, Agrif_mpi_comm, code)
  125. lubglob(1:nbdim,2) = - lubglob(1:nbdim,2)
  126. #endif
  127. !---------------------------------------------------------------------------------------------------
  128. end subroutine Agrif_get_var_global_bounds
  129. !===================================================================================================
  130. !
  131. !===================================================================================================
  132. ! subroutine Agrif_get_var_bounds
  133. !
  134. !> Gets the lower and the upper boundaries of a variable, for one particular direction.
  135. !---------------------------------------------------------------------------------------------------
  136. subroutine Agrif_get_var_bounds ( variable, lower, upper, index )
  137. !---------------------------------------------------------------------------------------------------
  138. type(Agrif_Variable), intent(in) :: variable !< Variable for which we want to extract boundaries
  139. integer, intent(out) :: lower !< Lower bound
  140. integer, intent(out) :: upper !< Upper bound
  141. integer, intent(in) :: index !< Direction for wich we want to know the boundaries
  142. !
  143. lower = variable % lb(index)
  144. upper = variable % ub(index)
  145. !---------------------------------------------------------------------------------------------------
  146. end subroutine Agrif_get_var_bounds
  147. !===================================================================================================
  148. !
  149. !===================================================================================================
  150. ! subroutine Agrif_get_var_bounds_array
  151. !
  152. !> Gets the lower and the upper boundaries of a table.
  153. !---------------------------------------------------------------------------------------------------
  154. subroutine Agrif_get_var_bounds_array ( variable, lower, upper, nbdim )
  155. !---------------------------------------------------------------------------------------------------
  156. type(Agrif_Variable), intent(in) :: variable !< Variable for which we want to extract boundaries
  157. integer, dimension(nbdim), intent(out) :: lower !< Lower bounds array
  158. integer, dimension(nbdim), intent(out) :: upper !< Upper bounds array
  159. integer, intent(in) :: nbdim !< Numer of dimensions of the variable
  160. !
  161. lower = variable % lb(1:nbdim)
  162. upper = variable % ub(1:nbdim)
  163. !---------------------------------------------------------------------------------------------------
  164. end subroutine Agrif_get_var_bounds_array
  165. !===================================================================================================
  166. !
  167. !===================================================================================================
  168. ! subroutine Agrif_array_allocate
  169. !
  170. !> Allocates data array in \b variable, according to the required dimension.
  171. !---------------------------------------------------------------------------------------------------
  172. subroutine Agrif_array_allocate ( variable, lb, ub, nbdim )
  173. !---------------------------------------------------------------------------------------------------
  174. type(Agrif_Variable), intent(inout) :: variable !< Variable struct for allocation
  175. integer, dimension(nbdim), intent(in) :: lb !< Lower bound
  176. integer, dimension(nbdim), intent(in) :: ub !< Upper bound
  177. integer, intent(in) :: nbdim !< Dimension of the array
  178. !
  179. select case (nbdim)
  180. case (1) ; allocate(variable%array1(lb(1):ub(1)))
  181. case (2) ; allocate(variable%array2(lb(1):ub(1),lb(2):ub(2)))
  182. case (3) ; allocate(variable%array3(lb(1):ub(1),lb(2):ub(2),lb(3):ub(3)))
  183. case (4) ; allocate(variable%array4(lb(1):ub(1),lb(2):ub(2),lb(3):ub(3),lb(4):ub(4)))
  184. case (5) ; allocate(variable%array5(lb(1):ub(1),lb(2):ub(2),lb(3):ub(3),lb(4):ub(4),lb(5):ub(5)))
  185. case (6) ; allocate(variable%array6(lb(1):ub(1),lb(2):ub(2),lb(3):ub(3),lb(4):ub(4),lb(5):ub(5),lb(6):ub(6)))
  186. end select
  187. !---------------------------------------------------------------------------------------------------
  188. end subroutine Agrif_array_allocate
  189. !===================================================================================================
  190. !
  191. !===================================================================================================
  192. ! subroutine Agrif_array_deallocate
  193. !
  194. !> Dellocates data array in \b variable, according to the required dimension.
  195. !---------------------------------------------------------------------------------------------------
  196. subroutine Agrif_array_deallocate ( variable, nbdim )
  197. !---------------------------------------------------------------------------------------------------
  198. type(Agrif_Variable), intent(inout) :: variable !< Variable struct for deallocation
  199. integer, intent(in) :: nbdim !< Dimension of the array
  200. !
  201. select case (nbdim)
  202. case (1) ; deallocate(variable%array1)
  203. case (2) ; deallocate(variable%array2)
  204. case (3) ; deallocate(variable%array3)
  205. case (4) ; deallocate(variable%array4)
  206. case (5) ; deallocate(variable%array5)
  207. case (6) ; deallocate(variable%array6)
  208. end select
  209. !---------------------------------------------------------------------------------------------------
  210. end subroutine Agrif_array_deallocate
  211. !===================================================================================================
  212. !
  213. !===================================================================================================
  214. ! subroutine Agrif_var_set_array_tozero
  215. !
  216. !> Reset the value of the data array in \b variable, according to the required dimension.
  217. !---------------------------------------------------------------------------------------------------
  218. subroutine Agrif_var_set_array_tozero ( variable, nbdim )
  219. !---------------------------------------------------------------------------------------------------
  220. type(Agrif_Variable), intent(inout) :: variable !< Variable
  221. integer, intent(in) :: nbdim !< Dimension of the array you want to reset
  222. !
  223. select case (nbdim)
  224. case (1) ; call Agrif_set_array_tozero_1D(variable%array1)
  225. case (2) ; call Agrif_set_array_tozero_2D(variable%array2)
  226. case (3) ; call Agrif_set_array_tozero_3D(variable%array3)
  227. case (4) ; call Agrif_set_array_tozero_4D(variable%array4)
  228. case (5) ; call Agrif_set_array_tozero_5D(variable%array5)
  229. case (6) ; call Agrif_set_array_tozero_6D(variable%array6)
  230. end select
  231. !---------------------------------------------------------------------------------------------------
  232. contains
  233. !---------------------------------------------------------------------------------------------------
  234. subroutine Agrif_set_array_tozero_1D ( array )
  235. real, dimension(:), intent(out) :: array
  236. array = 0.
  237. end subroutine Agrif_set_array_tozero_1D
  238. !
  239. subroutine Agrif_set_array_tozero_2D ( array )
  240. real, dimension(:,:), intent(out) :: array
  241. array = 0.
  242. end subroutine Agrif_set_array_tozero_2D
  243. !
  244. subroutine Agrif_set_array_tozero_3D ( array )
  245. real, dimension(:,:,:), intent(out) :: array
  246. array = 0.
  247. end subroutine Agrif_set_array_tozero_3D
  248. !
  249. subroutine Agrif_set_array_tozero_4D ( array )
  250. real, dimension(:,:,:,:), intent(out) :: array
  251. array = 0.
  252. end subroutine Agrif_set_array_tozero_4D
  253. !
  254. subroutine Agrif_set_array_tozero_5D ( array )
  255. real, dimension(:,:,:,:,:), intent(out) :: array
  256. array = 0.
  257. end subroutine Agrif_set_array_tozero_5D
  258. !
  259. subroutine Agrif_set_array_tozero_6D ( array )
  260. real, dimension(:,:,:,:,:,:), intent(out) :: array
  261. array = 0.
  262. end subroutine Agrif_set_array_tozero_6D
  263. !---------------------------------------------------------------------------------------------------
  264. end subroutine Agrif_var_set_array_tozero
  265. !===================================================================================================
  266. !
  267. !===================================================================================================
  268. ! subroutine Agrif_var_copy_array
  269. !
  270. !> Copy a part of data array from var2 to var1
  271. !---------------------------------------------------------------------------------------------------
  272. subroutine Agrif_var_copy_array ( var1, inf1, sup1, var2, inf2, sup2, nbdim )
  273. !---------------------------------------------------------------------------------------------------
  274. type(Agrif_Variable), intent(inout) :: var1 !< Modified variable
  275. integer, dimension(nbdim), intent(in) :: inf1 !< Lower boundary for var1
  276. integer, dimension(nbdim), intent(in) :: sup1 !< Upper boundary for var1
  277. type(Agrif_Variable), intent(in) :: var2 !< Input variable
  278. integer, dimension(nbdim), intent(in) :: inf2 !< Lower boundary for var2
  279. integer, dimension(nbdim), intent(in) :: sup2 !< Upper boundary for var2
  280. integer, intent(in) :: nbdim !< Dimension of the array
  281. !
  282. select case (nbdim)
  283. case (1) ; var1%array1(inf1(1):sup1(1)) = var2%array1(inf2(1):sup2(1))
  284. case (2) ; call Agrif_copy_array_2d( var1%array2, var2%array2, &
  285. lbound(var1%array2), lbound(var2%array2), inf1, sup1, inf2, sup2 )
  286. case (3) ; call Agrif_copy_array_3d( var1%array3, var2%array3, &
  287. lbound(var1%array3), lbound(var2%array3), inf1, sup1, inf2, sup2 )
  288. case (4) ; call Agrif_copy_array_4d( var1%array4, var2%array4, &
  289. lbound(var1%array4), lbound(var2%array4), inf1, sup1, inf2, sup2 )
  290. case (5) ; var1%array5(inf1(1):sup1(1), &
  291. inf1(2):sup1(2), &
  292. inf1(3):sup1(3), &
  293. inf1(4):sup1(4), &
  294. inf1(5):sup1(5)) = var2%array5(inf2(1):sup2(1), &
  295. inf2(2):sup2(2), &
  296. inf2(3):sup2(3), &
  297. inf2(4):sup2(4), &
  298. inf2(5):sup2(5))
  299. case (6) ; var1%array6(inf1(1):sup1(1), &
  300. inf1(2):sup1(2), &
  301. inf1(3):sup1(3), &
  302. inf1(4):sup1(4), &
  303. inf1(5):sup1(5), &
  304. inf1(6):sup1(6)) = var2%array6(inf2(1):sup2(1), &
  305. inf2(2):sup2(2), &
  306. inf2(3):sup2(3), &
  307. inf2(4):sup2(4), &
  308. inf2(5):sup2(5), &
  309. inf2(6):sup2(6))
  310. end select
  311. !---------------------------------------------------------------------------------------------------
  312. contains
  313. !---------------------------------------------------------------------------------------------------
  314. subroutine Agrif_copy_array_2d ( tabout, tabin, l, m, inf1, sup1, inf2, sup2 )
  315. integer, dimension(2), intent(in) :: l, m
  316. integer, dimension(2), intent(in) :: inf1, sup1
  317. integer, dimension(2), intent(in) :: inf2, sup2
  318. real, dimension(l(1):,l(2):), intent(out) :: tabout
  319. real, dimension(m(1):,m(2):), intent(in) :: tabin
  320. tabout(inf1(1):sup1(1), &
  321. inf1(2):sup1(2)) = tabin(inf2(1):sup2(1), &
  322. inf2(2):sup2(2))
  323. end subroutine Agrif_copy_array_2d
  324. !
  325. subroutine Agrif_copy_array_3d ( tabout, tabin, l, m, inf1, sup1, inf2, sup2 )
  326. integer, dimension(3), intent(in) :: l, m
  327. integer, dimension(3), intent(in) :: inf1, sup1
  328. integer, dimension(3), intent(in) :: inf2,sup2
  329. real, dimension(l(1):,l(2):,l(3):), intent(out) :: tabout
  330. real, dimension(m(1):,m(2):,m(3):), intent(in) :: tabin
  331. tabout(inf1(1):sup1(1), &
  332. inf1(2):sup1(2), &
  333. inf1(3):sup1(3)) = tabin(inf2(1):sup2(1), &
  334. inf2(2):sup2(2), &
  335. inf2(3):sup2(3))
  336. end subroutine Agrif_copy_array_3d
  337. !
  338. subroutine Agrif_copy_array_4d ( tabout, tabin, l, m, inf1, sup1, inf2, sup2 )
  339. integer, dimension(4), intent(in) :: l, m
  340. integer, dimension(4), intent(in) :: inf1, sup1
  341. integer, dimension(4), intent(in) :: inf2, sup2
  342. real, dimension(l(1):,l(2):,l(3):,l(4):), intent(out) :: tabout
  343. real, dimension(m(1):,m(2):,m(3):,m(4):), intent(in) :: tabin
  344. tabout(inf1(1):sup1(1), &
  345. inf1(2):sup1(2), &
  346. inf1(3):sup1(3), &
  347. inf1(4):sup1(4)) = tabin(inf2(1):sup2(1), &
  348. inf2(2):sup2(2), &
  349. inf2(3):sup2(3), &
  350. inf2(4):sup2(4))
  351. end subroutine Agrif_copy_array_4d
  352. !---------------------------------------------------------------------------------------------------
  353. end subroutine Agrif_var_copy_array
  354. !===================================================================================================
  355. !
  356. !===================================================================================================
  357. ! subroutine Agrif_var_full_copy_array
  358. !
  359. !> Copy the full data array from var2 to var1
  360. !---------------------------------------------------------------------------------------------------
  361. subroutine Agrif_var_full_copy_array ( var1, var2, nbdim )
  362. !---------------------------------------------------------------------------------------------------
  363. type(Agrif_Variable), intent(inout) :: var1 !< Modified variable
  364. type(Agrif_Variable), intent(in) :: var2 !< Input variable
  365. integer, intent(in) :: nbdim !< Dimension of the array
  366. !
  367. select case (nbdim)
  368. case (1) ; var1 % array1 = var2 % array1
  369. case (2) ; var1 % array2 = var2 % array2
  370. case (3) ; var1 % array3 = var2 % array3
  371. case (4) ; var1 % array4 = var2 % array4
  372. case (5) ; var1 % array5 = var2 % array5
  373. case (6) ; var1 % array6 = var2 % array6
  374. end select
  375. !---------------------------------------------------------------------------------------------------
  376. end subroutine Agrif_var_full_copy_array
  377. !===================================================================================================
  378. !
  379. !===================================================================================================
  380. ! subroutine GiveAgrif_SpecialValueToTab_mpi
  381. !
  382. !> Copy \b value in data array \b var2 where it is present in \b var1.
  383. !---------------------------------------------------------------------------------------------------
  384. subroutine GiveAgrif_SpecialValueToTab_mpi ( var1, var2, bounds, value, nbdim )
  385. !---------------------------------------------------------------------------------------------------
  386. type(Agrif_Variable), intent(in) :: var1 !< Modified variable
  387. type(Agrif_Variable), intent(inout) :: var2 !< Input variable
  388. integer, dimension(:,:,:), intent(in) :: bounds !< Bound for both arrays
  389. real, intent(in) :: value !< Special value
  390. integer, intent(in) :: nbdim !< Dimension of the array
  391. !
  392. select case (nbdim)
  393. case (1)
  394. where (var1 % array1(bounds(1,1,2):bounds(1,2,2)) == value )
  395. var2 % array1(bounds(1,1,1):bounds(1,2,1)) = value
  396. end where
  397. case (2)
  398. where (var1 % array2(bounds(1,1,2):bounds(1,2,2), &
  399. bounds(2,1,2):bounds(2,2,2)) == value)
  400. var2 % array2(bounds(1,1,1):bounds(1,2,1), &
  401. bounds(2,1,1):bounds(2,2,1)) = value
  402. end where
  403. case (3)
  404. where (var1 % array3(bounds(1,1,2):bounds(1,2,2), &
  405. bounds(2,1,2):bounds(2,2,2), &
  406. bounds(3,1,2):bounds(3,2,2)) == value)
  407. var2 % array3(bounds(1,1,1):bounds(1,2,1), &
  408. bounds(2,1,1):bounds(2,2,1), &
  409. bounds(3,1,1):bounds(3,2,1)) = value
  410. end where
  411. case (4)
  412. where (var1 % array4(bounds(1,1,2):bounds(1,2,2), &
  413. bounds(2,1,2):bounds(2,2,2), &
  414. bounds(3,1,2):bounds(3,2,2), &
  415. bounds(4,1,2):bounds(4,2,2)) == value)
  416. var2 % array4(bounds(1,1,1):bounds(1,2,1), &
  417. bounds(2,1,1):bounds(2,2,1), &
  418. bounds(3,1,1):bounds(3,2,1), &
  419. bounds(4,1,1):bounds(4,2,1)) = value
  420. end where
  421. case (5)
  422. where (var1 % array5(bounds(1,1,2):bounds(1,2,2), &
  423. bounds(2,1,2):bounds(2,2,2), &
  424. bounds(3,1,2):bounds(3,2,2), &
  425. bounds(4,1,2):bounds(4,2,2), &
  426. bounds(5,1,2):bounds(5,2,2)) == value)
  427. var2 % array5(bounds(1,1,1):bounds(1,2,1), &
  428. bounds(2,1,1):bounds(2,2,1), &
  429. bounds(3,1,1):bounds(3,2,1), &
  430. bounds(4,1,1):bounds(4,2,1), &
  431. bounds(5,1,1):bounds(5,2,1)) = value
  432. end where
  433. case (6)
  434. where (var1 % array6(bounds(1,1,2):bounds(1,2,2), &
  435. bounds(2,1,2):bounds(2,2,2), &
  436. bounds(3,1,2):bounds(3,2,2), &
  437. bounds(4,1,2):bounds(4,2,2), &
  438. bounds(5,1,2):bounds(5,2,2), &
  439. bounds(6,1,2):bounds(6,2,2)) == value)
  440. var2 % array6(bounds(1,1,1):bounds(1,2,1), &
  441. bounds(2,1,1):bounds(2,2,1), &
  442. bounds(3,1,1):bounds(3,2,1), &
  443. bounds(4,1,1):bounds(4,2,1), &
  444. bounds(5,1,1):bounds(5,2,1), &
  445. bounds(6,1,1):bounds(6,2,1)) = value
  446. end where
  447. end select
  448. !---------------------------------------------------------------------------------------------------
  449. end subroutine GiveAgrif_SpecialValueToTab_mpi
  450. !===================================================================================================
  451. !
  452. ! no more used ???
  453. #if 0
  454. !===================================================================================================
  455. ! subroutine GiveAgrif_SpecialValueToTab
  456. !---------------------------------------------------------------------------------------------------
  457. subroutine GiveAgrif_SpecialValueToTab ( var1, var2, &
  458. lower, upper, Value, nbdim )
  459. !---------------------------------------------------------------------------------------------------
  460. TYPE(Agrif_Variable), pointer :: var1
  461. TYPE(Agrif_Variable), pointer :: var2
  462. INTEGER, intent(in) :: nbdim
  463. INTEGER, DIMENSION(nbdim), intent(in) :: lower, upper
  464. REAL, intent(in) :: Value
  465. !
  466. select case (nbdim)
  467. case (1)
  468. where (var1 % array1( lower(1):upper(1)) == Value)
  469. var2 % array1(lower(1):upper(1)) = Value
  470. end where
  471. case (2)
  472. where (var1 % array2( lower(1):upper(1), &
  473. lower(2):upper(2)) == Value)
  474. var2 % array2(lower(1):upper(1), &
  475. lower(2):upper(2)) = Value
  476. end where
  477. case (3)
  478. where (var1 % array3( lower(1):upper(1), &
  479. lower(2):upper(2), &
  480. lower(3):upper(3)) == Value)
  481. var2 % array3(lower(1):upper(1), &
  482. lower(2):upper(2), &
  483. lower(3):upper(3)) = Value
  484. end where
  485. case (4)
  486. where (var1 % array4( lower(1):upper(1), &
  487. lower(2):upper(2), &
  488. lower(3):upper(3), &
  489. lower(4):upper(4)) == Value)
  490. var2 % array4(lower(1):upper(1), &
  491. lower(2):upper(2), &
  492. lower(3):upper(3), &
  493. lower(4):upper(4)) = Value
  494. end where
  495. case (5)
  496. where (var1 % array5( lower(1):upper(1), &
  497. lower(2):upper(2), &
  498. lower(3):upper(3), &
  499. lower(4):upper(4), &
  500. lower(5):upper(5)) == Value)
  501. var2 % array5(lower(1):upper(1), &
  502. lower(2):upper(2), &
  503. lower(3):upper(3), &
  504. lower(4):upper(4), &
  505. lower(5):upper(5)) = Value
  506. end where
  507. case (6)
  508. where (var1 % array6( lower(1):upper(1), &
  509. lower(2):upper(2), &
  510. lower(3):upper(3), &
  511. lower(4):upper(4), &
  512. lower(5):upper(5), &
  513. lower(6):upper(6)) == Value)
  514. var2 % array6(lower(1):upper(1), &
  515. lower(2):upper(2), &
  516. lower(3):upper(3), &
  517. lower(4):upper(4), &
  518. lower(5):upper(5), &
  519. lower(6):upper(6)) = Value
  520. end where
  521. end select
  522. !---------------------------------------------------------------------------------------------------
  523. end subroutine GiveAgrif_SpecialValueToTab
  524. !===================================================================================================
  525. #endif
  526. !
  527. #if defined AGRIF_MPI
  528. !===================================================================================================
  529. ! subroutine Agrif_var_replace_value
  530. !
  531. !> Replace \b value by \var2 content in \var1 data array.
  532. !---------------------------------------------------------------------------------------------------
  533. subroutine Agrif_var_replace_value ( var1, var2, lower, upper, value, nbdim )
  534. !---------------------------------------------------------------------------------------------------
  535. type(Agrif_Variable), intent(inout) :: var1 !< Modified variable
  536. type(Agrif_Variable), intent(in) :: var2 !< Input variable
  537. integer, dimension(nbdim), intent(in) :: lower !< Lower bound
  538. integer, dimension(nbdim), intent(in) :: upper !< Upper bound
  539. real, intent(in) :: value !< Special value
  540. integer, intent(in) :: nbdim !< Dimension of the array
  541. !
  542. integer :: i,j,k,l,m,n
  543. !
  544. select case (nbdim)
  545. case (1)
  546. do i = lower(1),upper(1)
  547. if (var1%array1(i) == value) then
  548. var1%array1(i) = var2%array1(i)
  549. endif
  550. enddo
  551. case (2)
  552. do j = lower(2),upper(2)
  553. do i = lower(1),upper(1)
  554. if (var1%array2(i,j) == value) then
  555. var1%array2(i,j) = var2%array2(i,j)
  556. endif
  557. enddo
  558. enddo
  559. case (3)
  560. do k = lower(3),upper(3)
  561. do j = lower(2),upper(2)
  562. do i = lower(1),upper(1)
  563. if (var1%array3(i,j,k) == value) then
  564. var1%array3(i,j,k) = var2%array3(i,j,k)
  565. endif
  566. enddo
  567. enddo
  568. enddo
  569. case (4)
  570. do l = lower(4),upper(4)
  571. do k = lower(3),upper(3)
  572. do j = lower(2),upper(2)
  573. do i = lower(1),upper(1)
  574. if (var1%array4(i,j,k,l) == value) then
  575. var1%array4(i,j,k,l) = var2%array4(i,j,k,l)
  576. endif
  577. enddo
  578. enddo
  579. enddo
  580. enddo
  581. case (5)
  582. do m = lower(5),upper(5)
  583. do l = lower(4),upper(4)
  584. do k = lower(3),upper(3)
  585. do j = lower(2),upper(2)
  586. do i = lower(1),upper(1)
  587. if (var1%array5(i,j,k,l,m) == value) then
  588. var1%array5(i,j,k,l,m) = var2%array5(i,j,k,l,m)
  589. endif
  590. enddo
  591. enddo
  592. enddo
  593. enddo
  594. enddo
  595. case (6)
  596. do n = lower(6),upper(6)
  597. do m = lower(5),upper(5)
  598. do l = lower(4),upper(4)
  599. do k = lower(3),upper(3)
  600. do j = lower(2),upper(2)
  601. do i = lower(1),upper(1)
  602. if (var1%array6(i,j,k,l,m,n) == value) then
  603. var1%array6(i,j,k,l,m,n) = var2%array6(i,j,k,l,m,n)
  604. endif
  605. enddo
  606. enddo
  607. enddo
  608. enddo
  609. enddo
  610. enddo
  611. end select
  612. !---------------------------------------------------------------------------------------------------
  613. end subroutine Agrif_var_replace_value
  614. !===================================================================================================
  615. #endif
  616. !
  617. !===================================================================================================
  618. ! subroutine PreProcessToInterpOrUpdate
  619. !---------------------------------------------------------------------------------------------------
  620. subroutine PreProcessToInterpOrUpdate ( parent, child, &
  621. nb_child, ub_child, &
  622. lb_child, lb_parent, &
  623. s_child, s_parent, &
  624. ds_child, ds_parent, nbdim, interp )
  625. !---------------------------------------------------------------------------------------------------
  626. type(Agrif_Variable), pointer, intent(in) :: parent !< Variable on the parent grid
  627. type(Agrif_Variable), pointer, intent(in) :: child !< Variable on the child grid
  628. integer, dimension(6), intent(out) :: nb_child !< Number of cells on the child grid
  629. integer, dimension(6), intent(out) :: ub_child !< Upper bound on the child grid
  630. integer, dimension(6), intent(out) :: lb_child !< Lower bound on the child grid
  631. integer, dimension(6), intent(out) :: lb_parent !< Lower bound on the parent grid
  632. real, dimension(6), intent(out) :: s_child !< Child grid position (s_root = 0)
  633. real, dimension(6), intent(out) :: s_parent !< Parent grid position (s_root = 0)
  634. real, dimension(6), intent(out) :: ds_child !< Child grid dx (ds_root = 1)
  635. real, dimension(6), intent(out) :: ds_parent !< Parent grid dx (ds_root = 1)
  636. integer, intent(out) :: nbdim !< Number of dimensions
  637. logical, intent(in) :: interp !< .true. if preprocess for interpolation, \n
  638. !! .false. if preprocess for update
  639. !
  640. type(Agrif_Variable), pointer :: root_var
  641. type(Agrif_Grid), pointer :: Agrif_Child_Gr
  642. type(Agrif_Grid), pointer :: Agrif_Parent_Gr
  643. integer :: n
  644. !
  645. Agrif_Child_Gr => Agrif_Curgrid
  646. Agrif_Parent_Gr => Agrif_Curgrid % parent
  647. !
  648. root_var => child % root_var
  649. !
  650. ! Number of dimensions of the current grid
  651. nbdim = root_var % nbdim
  652. !
  653. do n = 1,nbdim
  654. !
  655. ! Value of interptab(n) can be either x,y,z or N for a no space dimension
  656. select case(root_var % interptab(n))
  657. !
  658. case('x')
  659. !
  660. lb_child(n) = root_var % point(n)
  661. lb_parent(n) = root_var % point(n)
  662. nb_child(n) = Agrif_Child_Gr % nb(1)
  663. s_child(n) = Agrif_Child_Gr % Agrif_x(1)
  664. s_parent(n) = Agrif_Parent_Gr % Agrif_x(1)
  665. ds_child(n) = Agrif_Child_Gr % Agrif_dx(1)
  666. ds_parent(n) = Agrif_Parent_Gr % Agrif_dx(1)
  667. !
  668. if ( root_var % posvar(n) == 1 ) then
  669. ub_child(n) = lb_child(n) + Agrif_Child_Gr % nb(1)
  670. else
  671. ub_child(n) = lb_child(n) + Agrif_Child_Gr % nb(1) - 1
  672. s_child(n) = s_child(n) + 0.5*ds_child(n)
  673. s_parent(n) = s_parent(n) + 0.5*ds_parent(n)
  674. endif
  675. !
  676. case('y')
  677. !
  678. lb_child(n) = root_var % point(n)
  679. lb_parent(n) = root_var % point(n)
  680. nb_child(n) = Agrif_Child_Gr % nb(2)
  681. s_child(n) = Agrif_Child_Gr % Agrif_x(2)
  682. s_parent(n) = Agrif_Parent_Gr % Agrif_x(2)
  683. ds_child(n) = Agrif_Child_Gr % Agrif_dx(2)
  684. ds_parent(n) = Agrif_Parent_Gr % Agrif_dx(2)
  685. !
  686. if (root_var % posvar(n)==1) then
  687. ub_child(n) = lb_child(n) + Agrif_Child_Gr % nb(2)
  688. else
  689. ub_child(n) = lb_child(n) + Agrif_Child_Gr % nb(2) - 1
  690. s_child(n) = s_child(n) + 0.5*ds_child(n)
  691. s_parent(n) = s_parent(n) + 0.5*ds_parent(n)
  692. endif
  693. !
  694. case('z')
  695. !
  696. lb_child(n) = root_var % point(n)
  697. lb_parent(n) = root_var % point(n)
  698. nb_child(n) = Agrif_Child_Gr % nb(3)
  699. s_child(n) = Agrif_Child_Gr % Agrif_x(3)
  700. s_parent(n) = Agrif_Parent_Gr % Agrif_x(3)
  701. ds_child(n) = Agrif_Child_Gr % Agrif_dx(3)
  702. ds_parent(n) = Agrif_Parent_Gr % Agrif_dx(3)
  703. !
  704. if (root_var % posvar(n)==1) then
  705. ub_child(n) = lb_child(n) + Agrif_Child_Gr % nb(3)
  706. else
  707. ub_child(n) = lb_child(n) + Agrif_Child_Gr % nb(3) - 1
  708. s_child(n) = s_child(n) + 0.5*ds_child(n)
  709. s_parent(n) = s_parent(n) + 0.5*ds_parent(n)
  710. endif
  711. !
  712. case('N') ! No space dimension
  713. !
  714. ! The next coefficients are calculated in order to do a simple copy of
  715. ! values of the grid variable when the interpolation routine is
  716. ! called for this dimension.
  717. !
  718. if (interp) then
  719. call Agrif_get_var_bounds(parent, lb_child(n), ub_child(n), n)
  720. nb_child(n) = parent % ub(n) - parent % lb(n)
  721. else
  722. call Agrif_get_var_bounds(child, lb_child(n), ub_child(n), n)
  723. nb_child(n) = child % ub(n) - child % lb(n)
  724. endif
  725. !
  726. ! No interpolation but only a copy of the values of the grid variable
  727. lb_parent(n) = lb_child(n)
  728. s_child(n) = 0.
  729. s_parent(n) = 0.
  730. ds_child(n) = 1.
  731. ds_parent(n) = 1.
  732. !
  733. end select
  734. !
  735. enddo
  736. !---------------------------------------------------------------------------------------------------
  737. end subroutine PreProcessToInterpOrUpdate
  738. !===================================================================================================
  739. !
  740. #if defined AGRIF_MPI
  741. !===================================================================================================
  742. ! subroutine Agrif_GetLocalBoundaries
  743. !---------------------------------------------------------------------------------------------------
  744. subroutine Agrif_GetLocalBoundaries ( tab1, tab2, coord, lb, ub, deb, fin )
  745. !---------------------------------------------------------------------------------------------------
  746. integer, intent(in) :: tab1
  747. integer, intent(in) :: tab2
  748. integer, intent(in) :: coord
  749. integer, intent(in) :: lb, ub
  750. integer, intent(out) :: deb, fin
  751. !
  752. integer :: imin, imax
  753. integer :: i1, i2
  754. !
  755. call Agrif_InvLoc(lb, AGRIF_ProcRank, coord, imin)
  756. call Agrif_InvLoc(ub, AGRIF_ProcRank, coord, imax)
  757. !
  758. if ( imin > tab2 ) then
  759. i1 = imax - imin
  760. else
  761. i1 = max(tab1 - imin,0)
  762. endif
  763. !
  764. if (imax < tab1) then
  765. i2 = -(imax - imin)
  766. else
  767. i2 = min(tab2 - imax,0)
  768. endif
  769. !
  770. deb = lb + i1
  771. fin = ub + i2
  772. !---------------------------------------------------------------------------------------------------
  773. end subroutine Agrif_GetLocalBoundaries
  774. !===================================================================================================
  775. !
  776. !===================================================================================================
  777. ! subroutine Agrif_GlobalToLocalBounds
  778. !
  779. !> For a global index located on the current processor, tabarray gives the corresponding local index
  780. !---------------------------------------------------------------------------------------------------
  781. subroutine Agrif_GlobalToLocalBounds ( locbounds, lb_var, ub_var, lb_glob, ub_glob, &
  782. coords, nbdim, rank, member )
  783. !---------------------------------------------------------------------------------------------------
  784. integer, dimension(nbdim,2,2), intent(out) :: locbounds !< Local values of \b lb_glob and \b ub_glob
  785. integer, dimension(nbdim), intent(in) :: lb_var !< Local lower boundary on the current processor
  786. integer, dimension(nbdim), intent(in) :: ub_var !< Local upper boundary on the current processor
  787. integer, dimension(nbdim), intent(in) :: lb_glob !< Global lower boundary
  788. integer, dimension(nbdim), intent(in) :: ub_glob !< Global upper boundary
  789. integer, dimension(nbdim), intent(in) :: coords
  790. integer, intent(in) :: nbdim !< Dimension of the array
  791. integer, intent(in) :: rank !< Rank of the processor
  792. logical, intent(out) :: member
  793. !
  794. integer :: i, i1, k
  795. integer :: nbloc(nbdim)
  796. !
  797. locbounds(:,1,:) = HUGE(1)
  798. locbounds(:,2,:) = -HUGE(1)
  799. !
  800. nbloc = 0
  801. !
  802. do i = 1,nbdim
  803. !
  804. call Agrif_InvLoc(lb_var(i), rank, coords(i), i1)
  805. !
  806. do k = lb_glob(i)+lb_var(i)-i1,ub_glob(i)+lb_var(i)-i1
  807. !
  808. if ( (k >= lb_var(i)) .AND. (k <= ub_var(i)) ) then
  809. nbloc(i) = 1
  810. locbounds(i,1,1) = min(locbounds(i,1,1),k-lb_var(i)+i1)
  811. locbounds(i,2,1) = max(locbounds(i,2,1),k-lb_var(i)+i1)
  812. locbounds(i,1,2) = min(locbounds(i,1,2),k)
  813. locbounds(i,2,2) = max(locbounds(i,2,2),k)
  814. endif
  815. enddo
  816. enddo
  817. member = ( sum(nbloc) == nbdim )
  818. !---------------------------------------------------------------------------------------------------
  819. end subroutine Agrif_GlobalToLocalBounds
  820. !===================================================================================================
  821. #endif
  822. !
  823. end module Agrif_Arrays