modsauv.F90 26 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622623624625626627628629630631632633634635636637638639640641642643644645646647648649650651652653654655656657658659660661662663
  1. !
  2. ! $Id: modsauv.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_Save
  25. !!
  26. !! Module for the initialization by copy of the grids created by clustering.
  27. !
  28. module Agrif_Save
  29. !
  30. use Agrif_Types
  31. use Agrif_Link
  32. use Agrif_Arrays
  33. use Agrif_Variables
  34. !
  35. implicit none
  36. !
  37. contains
  38. !
  39. !===================================================================================================
  40. ! subroutine Agrif_deallocate_Arrays
  41. !---------------------------------------------------------------------------------------------------
  42. subroutine Agrif_deallocate_Arrays ( var )
  43. !---------------------------------------------------------------------------------------------------
  44. type(Agrif_Variable), pointer :: var
  45. !
  46. if (allocated(var%array1)) deallocate(var%array1)
  47. if (allocated(var%array2)) deallocate(var%array2)
  48. if (allocated(var%array3)) deallocate(var%array3)
  49. if (allocated(var%array4)) deallocate(var%array4)
  50. if (allocated(var%array5)) deallocate(var%array5)
  51. if (allocated(var%array6)) deallocate(var%array6)
  52. !
  53. if (allocated(var%darray1)) deallocate(var%darray1)
  54. if (allocated(var%darray2)) deallocate(var%darray2)
  55. if (allocated(var%darray3)) deallocate(var%darray3)
  56. if (allocated(var%darray4)) deallocate(var%darray4)
  57. if (allocated(var%darray5)) deallocate(var%darray5)
  58. if (allocated(var%darray6)) deallocate(var%darray6)
  59. !
  60. if (allocated(var%sarray1)) deallocate(var%sarray1)
  61. if (allocated(var%sarray2)) deallocate(var%sarray2)
  62. if (allocated(var%sarray3)) deallocate(var%sarray3)
  63. if (allocated(var%sarray4)) deallocate(var%sarray4)
  64. if (allocated(var%sarray5)) deallocate(var%sarray5)
  65. if (allocated(var%sarray6)) deallocate(var%sarray6)
  66. !
  67. !
  68. if (associated(var%oldvalues2D)) deallocate(var%oldvalues2D)
  69. if (allocated(var%posvar)) deallocate(var%posvar)
  70. if (allocated(var%interptab)) deallocate(var%interptab)
  71. if (allocated(var%coords)) deallocate(var%coords)
  72. !---------------------------------------------------------------------------------------------------
  73. end subroutine Agrif_deallocate_Arrays
  74. !---------------------------------------------------------------------------------------------------
  75. subroutine Agrif_deallocate_Arrays_c ( var_c )
  76. !---------------------------------------------------------------------------------------------------
  77. type(Agrif_Variable_c), pointer :: var_c
  78. !
  79. if (allocated(var_c%carray1)) deallocate(var_c%carray1)
  80. if (allocated(var_C%carray2)) deallocate(var_c%carray2)
  81. !
  82. !---------------------------------------------------------------------------------------------------
  83. end subroutine Agrif_deallocate_Arrays_c
  84. !===================================================================================================
  85. !---------------------------------------------------------------------------------------------------
  86. subroutine Agrif_deallocate_Arrays_l ( var_l )
  87. !---------------------------------------------------------------------------------------------------
  88. type(Agrif_Variable_l), pointer :: var_l
  89. !
  90. !
  91. if (allocated(var_l%larray1)) deallocate(var_l%larray1)
  92. if (allocated(var_l%larray2)) deallocate(var_l%larray2)
  93. if (allocated(var_l%larray3)) deallocate(var_l%larray3)
  94. if (allocated(var_l%larray4)) deallocate(var_l%larray4)
  95. if (allocated(var_l%larray5)) deallocate(var_l%larray5)
  96. if (allocated(var_l%larray6)) deallocate(var_l%larray6)
  97. !
  98. !---------------------------------------------------------------------------------------------------
  99. end subroutine Agrif_deallocate_Arrays_l
  100. !---------------------------------------------------------------------------------------------------
  101. subroutine Agrif_deallocate_Arrays_i ( var_i )
  102. !---------------------------------------------------------------------------------------------------
  103. type(Agrif_Variable_i), pointer :: var_i
  104. !
  105. !
  106. if (allocated(var_i%iarray1)) deallocate(var_i%iarray1)
  107. if (allocated(var_i%iarray2)) deallocate(var_i%iarray2)
  108. if (allocated(var_i%iarray3)) deallocate(var_i%iarray3)
  109. if (allocated(var_i%iarray4)) deallocate(var_i%iarray4)
  110. if (allocated(var_i%iarray5)) deallocate(var_i%iarray5)
  111. if (allocated(var_i%iarray6)) deallocate(var_i%iarray6)
  112. !
  113. !---------------------------------------------------------------------------------------------------
  114. end subroutine Agrif_deallocate_Arrays_i
  115. !
  116. !===================================================================================================
  117. ! subroutine Agrif_Free_data_before
  118. !---------------------------------------------------------------------------------------------------
  119. subroutine Agrif_Free_data_before ( Agrif_Gr )
  120. !---------------------------------------------------------------------------------------------------
  121. type(Agrif_Grid), pointer :: Agrif_Gr ! Pointer on the current grid
  122. !
  123. integer :: i
  124. type(Agrif_Variables_List), pointer :: parcours
  125. type(Agrif_Variable), pointer :: var
  126. type(Agrif_Variable_c), pointer :: var_c
  127. type(Agrif_Variable_r), pointer :: var_r
  128. type(Agrif_Variable_l), pointer :: var_l
  129. type(Agrif_Variable_i), pointer :: var_i
  130. !
  131. do i = 1,Agrif_NbVariables(0)
  132. !
  133. var => Agrif_Gr % tabvars(i)
  134. !
  135. if ( .NOT. Agrif_Mygrid % tabvars(i) % restore ) then
  136. call Agrif_deallocate_Arrays(var)
  137. endif
  138. !
  139. if (associated(var%list_interp)) then
  140. call Agrif_Free_list_interp(var%list_interp)
  141. endif
  142. !
  143. enddo
  144. do i=1,Agrif_NbVariables(1)
  145. var_c => Agrif_Gr % tabvars_c(i)
  146. call Agrif_deallocate_Arrays_c(var_c)
  147. enddo
  148. do i=1,Agrif_NbVariables(3)
  149. var_l => Agrif_Gr % tabvars_l(i)
  150. call Agrif_deallocate_Arrays_l(var_l)
  151. enddo
  152. do i=1,Agrif_NbVariables(4)
  153. var_i => Agrif_Gr % tabvars_i(i)
  154. call Agrif_deallocate_Arrays_i(var_i)
  155. enddo
  156. parcours => Agrif_Gr % variables
  157. do i = 1,Agrif_Gr%NbVariables
  158. !
  159. if ( .NOT. parcours%var%root_var % restore ) then
  160. call Agrif_deallocate_Arrays(parcours%var)
  161. endif
  162. !
  163. if (associated(parcours%var%list_interp)) then
  164. call Agrif_Free_list_interp(parcours%var%list_interp)
  165. endif
  166. !
  167. if ( .NOT. parcours%var%root_var % restore ) then
  168. deallocate(parcours%var)
  169. endif
  170. !
  171. parcours => parcours % next
  172. !
  173. enddo
  174. !
  175. if ( Agrif_USE_ONLY_FIXED_GRIDS == 0 ) then
  176. if ( Agrif_Probdim == 1 ) deallocate(Agrif_Gr%tabpoint1D)
  177. if ( Agrif_Probdim == 2 ) deallocate(Agrif_Gr%tabpoint2D)
  178. if ( Agrif_Probdim == 3 ) deallocate(Agrif_Gr%tabpoint3D)
  179. endif
  180. !---------------------------------------------------------------------------------------------------
  181. end subroutine Agrif_Free_data_before
  182. !===================================================================================================
  183. !
  184. !===================================================================================================
  185. ! subroutine Agrif_Free_list_interp
  186. !---------------------------------------------------------------------------------------------------
  187. recursive subroutine Agrif_Free_list_interp ( list_interp )
  188. !---------------------------------------------------------------------------------------------------
  189. type(Agrif_List_Interp_Loc), pointer :: list_interp
  190. !
  191. if (associated(list_interp%suiv)) call Agrif_Free_list_interp(list_interp%suiv)
  192. #if defined AGRIF_MPI
  193. deallocate(list_interp%interp_loc%tab4t)
  194. deallocate(list_interp%interp_loc%memberinall)
  195. deallocate(list_interp%interp_loc%sendtoproc1)
  196. deallocate(list_interp%interp_loc%recvfromproc1)
  197. #endif
  198. deallocate(list_interp%interp_loc)
  199. deallocate(list_interp)
  200. nullify(list_interp)
  201. !---------------------------------------------------------------------------------------------------
  202. end subroutine Agrif_Free_list_interp
  203. !===================================================================================================
  204. !
  205. !===================================================================================================
  206. ! subroutine Agrif_Free_data_after
  207. !---------------------------------------------------------------------------------------------------
  208. subroutine Agrif_Free_data_after ( Agrif_Gr )
  209. !---------------------------------------------------------------------------------------------------
  210. type(Agrif_Grid), pointer :: Agrif_Gr !< Pointer on the current grid
  211. !
  212. integer :: i
  213. type(Agrif_Variables_List), pointer :: parcours, rootparcours
  214. type(Agrif_Variable), pointer :: var
  215. !
  216. do i = 1,Agrif_NbVariables(0)
  217. if ( Agrif_Mygrid % tabvars(i) % restore ) then
  218. var => Agrif_Gr%tabvars(i)
  219. call Agrif_deallocate_Arrays(var)
  220. endif
  221. enddo
  222. !
  223. parcours => Agrif_Gr%variables
  224. rootparcours => Agrif_Mygrid%variables
  225. !
  226. do i = 1,Agrif_Gr%NbVariables
  227. if (rootparcours%var % restore ) then
  228. call Agrif_deallocate_Arrays(parcours%var)
  229. deallocate(parcours%var)
  230. endif
  231. parcours => parcours % next
  232. rootparcours => rootparcours % next
  233. enddo
  234. !---------------------------------------------------------------------------------------------------
  235. end subroutine Agrif_Free_data_after
  236. !===================================================================================================
  237. !
  238. !===================================================================================================
  239. ! subroutine Agrif_CopyFromOld_All
  240. !
  241. !> Called in Agrif_Clustering#Agrif_Init_Hierarchy.
  242. !---------------------------------------------------------------------------------------------------
  243. recursive subroutine Agrif_CopyFromOld_All ( g, oldchildlist )
  244. !---------------------------------------------------------------------------------------------------
  245. type(Agrif_Grid), pointer, intent(inout) :: g !< Pointer on the current grid
  246. type(Agrif_Grid_List), intent(in) :: oldchildlist
  247. !
  248. type(Agrif_PGrid), pointer :: parcours ! Pointer for the recursive procedure
  249. real :: g_eps, eps, oldgrid_eps
  250. integer :: out
  251. integer :: iii
  252. !
  253. out = 0
  254. !
  255. parcours => oldchildlist % first
  256. !
  257. do while (associated(parcours))
  258. !
  259. if ((.NOT. g % fixed) .AND. (parcours % gr %oldgrid)) then
  260. !
  261. g_eps = huge(1.)
  262. oldgrid_eps = huge(1.)
  263. do iii = 1,Agrif_Probdim
  264. g_eps = min(g_eps,g % Agrif_dx(iii))
  265. oldgrid_eps = min(oldgrid_eps, parcours % gr % Agrif_dx(iii))
  266. enddo
  267. !
  268. eps = min(g_eps,oldgrid_eps)/100.
  269. !
  270. do iii = 1,Agrif_Probdim
  271. if (g % Agrif_dx(iii) < (parcours % gr % Agrif_dx(iii) - eps)) then
  272. parcours => parcours % next
  273. out = 1
  274. exit
  275. endif
  276. enddo
  277. !
  278. if ( out == 1 ) cycle
  279. !
  280. call Agrif_CopyFromOld(g,parcours%gr)
  281. !
  282. endif
  283. !
  284. call Agrif_CopyFromOld_All(g,parcours % gr % child_list)
  285. !
  286. parcours => parcours % next
  287. !
  288. enddo
  289. !---------------------------------------------------------------------------------------------------
  290. end subroutine Agrif_CopyFromOld_All
  291. !===================================================================================================
  292. !
  293. !===================================================================================================
  294. ! subroutine Agrif_CopyFromOld
  295. !
  296. !> Call to the Agrif_Copy procedure.
  297. !---------------------------------------------------------------------------------------------------
  298. subroutine Agrif_CopyFromOld ( new_gr, old_gr )
  299. !---------------------------------------------------------------------------------------------------
  300. type(Agrif_Grid), pointer, intent(inout) :: new_gr !< Pointer on the current grid
  301. type(Agrif_Grid), pointer, intent(inout) :: old_gr !< Pointer on an old grid
  302. !
  303. type(Agrif_Variable), pointer :: new_var
  304. type(Agrif_Variable), pointer :: old_var
  305. integer :: i
  306. !
  307. do i = 1,Agrif_NbVariables(0)
  308. if ( Agrif_Mygrid % tabvars(i) % restore ) then
  309. old_var => old_gr % tabvars(i)
  310. new_var => new_gr % tabvars(i)
  311. call Agrif_Copy(new_gr, old_gr, new_var, old_var )
  312. endif
  313. enddo
  314. !---------------------------------------------------------------------------------------------------
  315. end subroutine Agrif_CopyFromOld
  316. !===================================================================================================
  317. !
  318. !===================================================================================================
  319. ! subroutine Agrif_CopyFromOld_AllOneVar
  320. !
  321. !> Called in Agrif_Clustering#Agrif_Init_Hierarchy.
  322. !---------------------------------------------------------------------------------------------------
  323. recursive subroutine Agrif_CopyFromOld_AllOneVar ( g, oldchildlist, indic )
  324. !---------------------------------------------------------------------------------------------------
  325. type(Agrif_Grid), pointer, intent(inout) :: g !< Pointer on the current grid
  326. type(Agrif_Grid_List), intent(in) :: oldchildlist
  327. integer, intent(in) :: indic
  328. !
  329. type(Agrif_PGrid), pointer :: parcours ! Pointer for the recursive procedure
  330. real :: g_eps,eps,oldgrid_eps
  331. integer :: out
  332. integer :: iii
  333. !
  334. out = 0
  335. !
  336. parcours => oldchildlist % first
  337. !
  338. do while (associated(parcours))
  339. !
  340. if ((.NOT. g % fixed) .AND. (parcours % gr %oldgrid)) then
  341. !
  342. g_eps = huge(1.)
  343. oldgrid_eps = huge(1.)
  344. do iii = 1,Agrif_Probdim
  345. g_eps = min(g_eps,g % Agrif_dx(iii))
  346. oldgrid_eps = min(oldgrid_eps,parcours % gr % Agrif_dx(iii))
  347. enddo
  348. !
  349. eps = min(g_eps,oldgrid_eps)/100.
  350. !
  351. do iii = 1,Agrif_Probdim
  352. if (g % Agrif_dx(iii) < (parcours % gr % Agrif_dx(iii) - eps)) then
  353. parcours => parcours % next
  354. out = 1
  355. exit
  356. endif
  357. enddo
  358. if ( out == 1 ) cycle
  359. !
  360. call Agrif_CopyFromOldOneVar(g,parcours%gr,indic)
  361. !
  362. endif
  363. !
  364. call Agrif_CopyFromOld_AllOneVar(g,parcours%gr % child_list,indic)
  365. !
  366. parcours => parcours % next
  367. !
  368. enddo
  369. !---------------------------------------------------------------------------------------------------
  370. end subroutine Agrif_CopyFromOld_AllOneVar
  371. !===================================================================================================
  372. !
  373. !===================================================================================================
  374. ! subroutine Agrif_CopyFromOldOneVar
  375. !
  376. !> Call to Agrif_Copy
  377. !---------------------------------------------------------------------------------------------------
  378. subroutine Agrif_CopyFromOldOneVar ( new_gr, old_gr, indic )
  379. !---------------------------------------------------------------------------------------------------
  380. type(Agrif_Grid), pointer, intent(inout) :: new_gr !< Pointer on the current grid
  381. type(Agrif_Grid), pointer, intent(in) :: old_gr !< Pointer on an old grid
  382. integer, intent(in) :: indic
  383. !
  384. type(Agrif_Variable), pointer :: new_var, old_var
  385. !
  386. new_var => Agrif_Search_Variable(new_gr, -indic)
  387. old_var => Agrif_Search_Variable(old_gr, -indic)
  388. !
  389. call Agrif_array_allocate(new_var,new_var%lb,new_var%ub,new_var%nbdim)
  390. call Agrif_Copy(new_gr, old_gr, new_var,old_var)
  391. !---------------------------------------------------------------------------------------------------
  392. end subroutine Agrif_CopyFromOldOneVar
  393. !===================================================================================================
  394. !
  395. !===================================================================================================
  396. ! subroutine Agrif_Copy
  397. !---------------------------------------------------------------------------------------------------
  398. subroutine Agrif_Copy ( new_gr, old_gr, new_var, old_var )
  399. !---------------------------------------------------------------------------------------------------
  400. type(Agrif_Grid), pointer, intent(in) :: new_gr !< Pointer on the current grid
  401. type(Agrif_Grid), pointer, intent(in) :: old_gr !< Pointer on an old grid
  402. type(Agrif_Variable), pointer, intent(inout) :: new_var
  403. type(Agrif_Variable), pointer, intent(in) :: old_var
  404. !
  405. type(Agrif_Variable), pointer :: root ! Pointer on the variable of the root grid
  406. integer :: nbdim ! Number of dimensions of the current grid
  407. integer, dimension(6) :: pttabnew ! Indexes of the first point in the domain
  408. integer, dimension(6) :: petabnew ! Indexes of the first point in the domain
  409. integer, dimension(6) :: pttabold ! Indexes of the first point in the domain
  410. integer, dimension(6) :: petabold ! Indexes of the first point in the domain
  411. integer, dimension(6) :: nbtabold ! Number of cells in each direction
  412. integer, dimension(6) :: nbtabnew ! Number of cells in each direction
  413. real, dimension(6) :: snew,sold
  414. real, dimension(6) :: dsnew,dsold
  415. real :: eps
  416. integer :: n
  417. !
  418. root => new_var % root_var
  419. nbdim = root % nbdim
  420. !
  421. do n = 1,nbdim
  422. !
  423. select case(root % interptab(n))
  424. !
  425. case('x')
  426. !
  427. pttabnew(n) = root % point(n)
  428. pttabold(n) = root % point(n)
  429. snew(n) = new_gr % Agrif_x(1)
  430. sold(n) = old_gr % Agrif_x(1)
  431. dsnew(n) = new_gr % Agrif_dx(1)
  432. dsold(n) = old_gr % Agrif_dx(1)
  433. !
  434. if (root % posvar(n) == 1) then
  435. petabnew(n) = pttabnew(n) + new_gr % nb(1)
  436. petabold(n) = pttabold(n) + old_gr % nb(1)
  437. else
  438. petabnew(n) = pttabnew(n) + new_gr % nb(1) - 1
  439. petabold(n) = pttabold(n) + old_gr % nb(1) - 1
  440. snew(n) = snew(n) + dsnew(n)/2.
  441. sold(n) = sold(n) + dsold(n)/2.
  442. endif
  443. !
  444. case('y')
  445. !
  446. pttabnew(n) = root % point(n)
  447. pttabold(n) = root % point(n)
  448. snew(n) = new_gr % Agrif_x(2)
  449. sold(n) = old_gr % Agrif_x(2)
  450. dsnew(n) = new_gr % Agrif_dx(2)
  451. dsold(n) = old_gr % Agrif_dx(2)
  452. !
  453. if (root % posvar(n) == 1) then
  454. petabnew(n) = pttabnew(n) + new_gr % nb(2)
  455. petabold(n) = pttabold(n) + old_gr % nb(2)
  456. else
  457. petabnew(n) = pttabnew(n) + new_gr % nb(2) - 1
  458. petabold(n) = pttabold(n) + old_gr % nb(2) - 1
  459. snew(n) = snew(n) + dsnew(n)/2.
  460. sold(n) = sold(n) + dsold(n)/2.
  461. endif
  462. !
  463. case('z')
  464. !
  465. pttabnew(n) = root % point(n)
  466. pttabold(n) = root % point(n)
  467. snew(n) = new_gr % Agrif_x(3)
  468. sold(n) = old_gr % Agrif_x(3)
  469. dsnew(n) = new_gr % Agrif_dx(3)
  470. dsold(n) = old_gr % Agrif_dx(3)
  471. !
  472. if (root % posvar(n) == 1) then
  473. petabnew(n) = pttabnew(n) + new_gr % nb(3)
  474. petabold(n) = pttabold(n) + old_gr % nb(3)
  475. else
  476. petabnew(n) = pttabnew(n) + new_gr % nb(3) - 1
  477. petabold(n) = pttabold(n) + old_gr % nb(3) - 1
  478. snew(n) = snew(n) + dsnew(n)/2.
  479. sold(n) = sold(n) + dsold(n)/2.
  480. endif
  481. !
  482. case('N')
  483. !
  484. call Agrif_get_var_bounds(new_var,pttabnew(n),petabnew(n),n)
  485. !
  486. pttabold(n) = pttabnew(n)
  487. petabold(n) = petabnew(n)
  488. snew(n) = 0.
  489. sold(n) = 0.
  490. dsnew(n) = 1.
  491. dsold(n) = 1.
  492. !
  493. end select
  494. !
  495. enddo
  496. !
  497. do n = 1,nbdim
  498. nbtabnew(n) = petabnew(n) - pttabnew(n)
  499. nbtabold(n) = petabold(n) - pttabold(n)
  500. enddo
  501. !
  502. eps = min(minval(dsnew(1:nbdim)),minval(dsold(1:nbdim))) / 100.
  503. !
  504. do n = 1,nbdim
  505. if (snew(n) > (sold(n)+dsold(n)*nbtabold(n)+eps)) return
  506. if ((snew(n)+dsnew(n)*nbtabnew(n)-eps) < sold(n)) return
  507. enddo
  508. !
  509. call Agrif_CopynD(new_var,old_var,pttabold,petabold,pttabnew,petabnew, &
  510. sold,snew,dsold,dsnew,nbdim)
  511. !---------------------------------------------------------------------------------------------------
  512. end subroutine Agrif_Copy
  513. !===================================================================================================
  514. !
  515. !===================================================================================================
  516. ! subroutine Agrif_CopynD
  517. !
  518. !> Copy from the nD variable Old_Var the nD variable New_Var
  519. !---------------------------------------------------------------------------------------------------
  520. subroutine Agrif_CopynD ( new_var, old_var, pttabold, petabold, pttabnew, petabnew, &
  521. sold, snew, dsold, dsnew, nbdim )
  522. !---------------------------------------------------------------------------------------------------
  523. type(Agrif_Variable), pointer, intent(inout) :: new_var
  524. type(Agrif_Variable), pointer, intent(in) :: old_var
  525. integer, dimension(nbdim), intent(in) :: pttabnew
  526. integer, dimension(nbdim), intent(in) :: petabnew
  527. integer, dimension(nbdim), intent(in) :: pttabold
  528. integer, dimension(nbdim), intent(in) :: petabold
  529. real, dimension(nbdim), intent(in) :: snew, sold
  530. real, dimension(nbdim), intent(in) :: dsnew,dsold
  531. integer, intent(in) :: nbdim
  532. !
  533. integer :: i,j,k,l,m,n,i0,j0,k0,l0,m0,n0
  534. !
  535. real, dimension(nbdim) :: dim_gmin, dim_gmax
  536. real, dimension(nbdim) :: dim_newmin, dim_newmax
  537. real, dimension(nbdim) :: dim_min
  538. integer, dimension(nbdim) :: ind_gmin,ind_newmin, ind_newmax
  539. !
  540. do i = 1,nbdim
  541. !
  542. dim_gmin(i) = sold(i)
  543. dim_gmax(i) = sold(i) + (petabold(i)-pttabold(i)) * dsold(i)
  544. !
  545. dim_newmin(i) = snew(i)
  546. dim_newmax(i) = snew(i) + (petabnew(i)-pttabnew(i)) * dsnew(i)
  547. !
  548. enddo
  549. !
  550. do i = 1,nbdim
  551. if (dim_gmax(i) < dim_newmin(i)) return
  552. if (dim_gmin(i) > dim_newmax(i)) return
  553. enddo
  554. !
  555. do i = 1,nbdim
  556. !
  557. ind_newmin(i) = pttabnew(i) - floor(-(max(dim_gmin(i),dim_newmin(i))-dim_newmin(i))/dsnew(i))
  558. dim_min(i) = snew(i) + (ind_newmin(i)-pttabnew(i))*dsnew(i)
  559. ind_gmin(i) = pttabold(i) + nint((dim_min(i)-dim_gmin(i))/dsold(i))
  560. ind_newmax(i) = pttabnew(i) + int((min(dim_gmax(i),dim_newmax(i))-dim_newmin(i))/dsnew(i))
  561. !
  562. enddo
  563. !
  564. select case (nbdim)
  565. !
  566. case (1)
  567. i0 = ind_gmin(1)
  568. do i = ind_newmin(1),ind_newmax(1)
  569. new_var % array1(i) = old_var % array1(i0)
  570. new_var % restore1D(i) = 1
  571. i0 = i0 + int(dsnew(1)/dsold(1))
  572. enddo
  573. !
  574. case (2)
  575. i0 = ind_gmin(1) ; do i = ind_newmin(1),ind_newmax(1)
  576. j0 = ind_gmin(2) ; do j = ind_newmin(2),ind_newmax(2)
  577. new_var % array2(i,j) = old_var % array2(i0,j0)
  578. new_var % restore2D(i,j) = 1
  579. j0 = j0 + int(dsnew(2)/dsold(2))
  580. enddo
  581. i0 = i0 + int(dsnew(1)/dsold(1))
  582. enddo
  583. !
  584. case (3)
  585. i0 = ind_gmin(1) ; do i = ind_newmin(1),ind_newmax(1)
  586. j0 = ind_gmin(2) ; do j = ind_newmin(2),ind_newmax(2)
  587. k0 = ind_gmin(3) ; do k = ind_newmin(3),ind_newmax(3)
  588. new_var % array3(i,j,k) = old_var % array3(i0,j0,k0)
  589. new_var % restore3D(i,j,k) = 1
  590. k0 = k0 + int(dsnew(3)/dsold(3))
  591. enddo
  592. j0 = j0 + int(dsnew(2)/dsold(2))
  593. enddo
  594. i0 = i0 + int(dsnew(1)/dsold(1))
  595. enddo
  596. !
  597. case (4)
  598. i0 = ind_gmin(1) ; do i = ind_newmin(1),ind_newmax(1)
  599. j0 = ind_gmin(2) ; do j = ind_newmin(2),ind_newmax(2)
  600. k0 = ind_gmin(3) ; do k = ind_newmin(3),ind_newmax(3)
  601. l0 = ind_gmin(4) ; do l = ind_newmin(4),ind_newmax(4)
  602. new_var % array4(i,j,k,l) = old_var % array4(i0,j0,k0,l0)
  603. new_var % restore4D(i,j,k,l) = 1
  604. l0 = l0 + int(dsnew(4)/dsold(4))
  605. enddo
  606. k0 = k0 + int(dsnew(3)/dsold(3))
  607. enddo
  608. j0 = j0 + int(dsnew(2)/dsold(2))
  609. enddo
  610. i0 = i0 + int(dsnew(1)/dsold(1))
  611. enddo
  612. !
  613. case (5)
  614. i0 = ind_gmin(1) ; do i = ind_newmin(1),ind_newmax(1)
  615. j0 = ind_gmin(2) ; do j = ind_newmin(2),ind_newmax(2)
  616. k0 = ind_gmin(3) ; do k = ind_newmin(3),ind_newmax(3)
  617. l0 = ind_gmin(4) ; do l = ind_newmin(4),ind_newmax(4)
  618. m0 = ind_gmin(5) ; do m = ind_newmin(5),ind_newmax(5)
  619. new_var % array5(i,j,k,l,m) = old_var % array5(i0,j0,k0,l0,m0)
  620. new_var % restore5D(i,j,k,l,m) = 1
  621. m0 = m0 + int(dsnew(5)/dsold(5))
  622. enddo
  623. l0 = l0 + int(dsnew(4)/dsold(4))
  624. enddo
  625. k0 = k0 + int(dsnew(3)/dsold(3))
  626. enddo
  627. j0 = j0 + int(dsnew(2)/dsold(2))
  628. enddo
  629. i0 = i0 + int(dsnew(1)/dsold(1))
  630. enddo
  631. !
  632. case (6)
  633. i0 = ind_gmin(1) ; do i = ind_newmin(1),ind_newmax(1)
  634. j0 = ind_gmin(2) ; do j = ind_newmin(2),ind_newmax(2)
  635. k0 = ind_gmin(3) ; do k = ind_newmin(3),ind_newmax(3)
  636. l0 = ind_gmin(4) ; do l = ind_newmin(4),ind_newmax(4)
  637. m0 = ind_gmin(5) ; do m = ind_newmin(5),ind_newmax(5)
  638. n0 = ind_gmin(6) ; do n = ind_newmin(6),ind_newmax(6)
  639. new_var % array6(i,j,k,l,m,n) = old_var % array6(i0,j0,k0,l0,m0,n0)
  640. new_var % restore6D(i,j,k,l,m,n) = 1
  641. n0 = n0 + int(dsnew(6)/dsold(6))
  642. enddo
  643. m0 = m0 + int(dsnew(5)/dsold(5))
  644. enddo
  645. l0 = l0 + int(dsnew(4)/dsold(4))
  646. enddo
  647. k0 = k0 + int(dsnew(3)/dsold(3))
  648. enddo
  649. j0 = j0 + int(dsnew(2)/dsold(2))
  650. enddo
  651. i0 = i0 + int(dsnew(1)/dsold(1))
  652. enddo
  653. !
  654. end select
  655. !---------------------------------------------------------------------------------------------------
  656. end subroutine Agrif_CopynD
  657. !===================================================================================================
  658. !
  659. end module Agrif_Save