modbcfunction.F90 30 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622623624625626627628629630631632633634635636637638639640641642643644645646647648649650651652653654655656657658659660661662663664665666667668669670671672673674675676677678679680
  1. !
  2. ! $Id: modbcfunction.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. !> Module Agrif_BcFunction.
  24. !
  25. module Agrif_BcFunction
  26. !
  27. ! Modules used:
  28. !
  29. use Agrif_Boundary
  30. use Agrif_Update
  31. use Agrif_Save
  32. !
  33. implicit none
  34. !
  35. interface Agrif_Set_Parent
  36. module procedure Agrif_Set_Parent_int, &
  37. Agrif_Set_Parent_real4, &
  38. Agrif_Set_Parent_real8
  39. end interface
  40. !
  41. interface Agrif_Save_Forrestore
  42. module procedure Agrif_Save_Forrestore0d, &
  43. Agrif_Save_Forrestore2d, &
  44. Agrif_Save_Forrestore3d, &
  45. Agrif_Save_Forrestore4d
  46. end interface
  47. !
  48. contains
  49. !
  50. !===================================================================================================
  51. ! subroutine Agrif_Set_parent_int
  52. !
  53. !> To set the TYPE of the variable
  54. !---------------------------------------------------------------------------------------------------
  55. subroutine Agrif_Set_parent_int(tabvarsindic,value)
  56. !---------------------------------------------------------------------------------------------------
  57. integer, intent(in) :: tabvarsindic !< indice of the variable in tabvars
  58. integer, intent(in) :: value !< input value
  59. !
  60. Agrif_Curgrid % parent % tabvars_i(tabvarsindic) % iarray0 = value
  61. !---------------------------------------------------------------------------------------------------
  62. end subroutine Agrif_Set_parent_int
  63. !===================================================================================================
  64. !
  65. !===================================================================================================
  66. ! subroutine Agrif_Set_parent_real4
  67. !---------------------------------------------------------------------------------------------------
  68. !> To set the TYPE of the variable
  69. !---------------------------------------------------------------------------------------------------
  70. subroutine Agrif_Set_parent_real4 ( tabvarsindic, value )
  71. !---------------------------------------------------------------------------------------------------
  72. integer, intent(in) :: tabvarsindic !< indice of the variable in tabvars
  73. real(kind=4),intent(in) :: value !< input value
  74. !
  75. Agrif_Curgrid % parent % tabvars_r(tabvarsindic) % array0 = value
  76. Agrif_Curgrid % parent % tabvars_r(tabvarsindic) % sarray0 = value
  77. !---------------------------------------------------------------------------------------------------
  78. end subroutine Agrif_Set_parent_real4
  79. !===================================================================================================
  80. !
  81. !===================================================================================================
  82. ! subroutine Agrif_Set_parent_real8
  83. !---------------------------------------------------------------------------------------------------
  84. !> To set the TYPE of the variable
  85. !---------------------------------------------------------------------------------------------------
  86. subroutine Agrif_Set_parent_real8 ( tabvarsindic, value )
  87. !---------------------------------------------------------------------------------------------------
  88. integer, intent(in) :: tabvarsindic !< indice of the variable in tabvars
  89. real(kind=8),intent(in) :: value !< input value
  90. !
  91. Agrif_Curgrid % parent % tabvars_r(tabvarsindic) % darray0 = value
  92. !---------------------------------------------------------------------------------------------------
  93. end subroutine Agrif_Set_parent_real8
  94. !===================================================================================================
  95. !
  96. !===================================================================================================
  97. ! subroutine Agrif_Set_bc
  98. !---------------------------------------------------------------------------------------------------
  99. subroutine Agrif_Set_bc ( tabvarsindic, bcinfsup, Interpolationshouldbemade )
  100. !---------------------------------------------------------------------------------------------------
  101. integer, intent(in) :: tabvarsindic !< indice of the variable in tabvars
  102. integer, dimension(2), intent(in) :: bcinfsup !< bcinfsup
  103. logical, optional, intent(in) :: Interpolationshouldbemade !< interpolation should be made
  104. !
  105. integer :: indic ! indice of the variable in tabvars
  106. type(Agrif_Variable), pointer :: var
  107. !
  108. indic = Agrif_Curgrid % tabvars_i(tabvarsindic) % iarray0
  109. !
  110. if (indic <= 0) then
  111. var => Agrif_Search_Variable(Agrif_Curgrid,-indic)
  112. else
  113. print*,"Agrif_Set_bc : warning indic >= 0 !!!"
  114. var => Agrif_Curgrid % tabvars(indic)
  115. endif
  116. if (.not.associated(var)) return ! Grand mother grid case
  117. !
  118. if ( Agrif_Curgrid % fixedrank /= 0 ) then
  119. if ( .not.associated(var % oldvalues2D) ) then
  120. allocate(var % oldvalues2D(2,1))
  121. var % interpIndex = -1
  122. var % oldvalues2D = 0.
  123. endif
  124. if ( present(Interpolationshouldbemade) ) then
  125. var % Interpolationshouldbemade = Interpolationshouldbemade
  126. endif
  127. endif
  128. !
  129. var % bcinf = bcinfsup(1)
  130. var % bcsup = bcinfsup(2)
  131. !---------------------------------------------------------------------------------------------------
  132. end subroutine Agrif_Set_bc
  133. !===================================================================================================
  134. !
  135. !===================================================================================================
  136. ! subroutine Agrif_Set_interp
  137. !---------------------------------------------------------------------------------------------------
  138. subroutine Agrif_Set_interp ( tabvarsindic, interp, interp1, interp2, interp3 , interp4)
  139. !---------------------------------------------------------------------------------------------------
  140. integer, intent(in) :: tabvarsindic !< indice of the variable in tabvars
  141. integer, optional, intent(in) :: interp, interp1, interp2, interp3, interp4
  142. !
  143. integer :: indic ! indice of the variable in tabvars
  144. type(Agrif_Variable), pointer :: var
  145. !
  146. indic = Agrif_Curgrid % tabvars_i(tabvarsindic) % iarray0
  147. !
  148. if (indic <= 0) then
  149. var => Agrif_Search_Variable(Agrif_Mygrid,-indic)
  150. else
  151. print*,"Agrif_Set_interp : warning indic >= 0 !!!"
  152. var => Agrif_Mygrid % tabvars(indic)
  153. endif
  154. !
  155. var % type_interp = Agrif_Constant
  156. !
  157. if (present(interp)) var % type_interp = interp
  158. if (present(interp1)) var % type_interp(1) = interp1
  159. if (present(interp2)) var % type_interp(2) = interp2
  160. if (present(interp3)) var % type_interp(3) = interp3
  161. if (present(interp4)) var % type_interp(4) = interp4
  162. !---------------------------------------------------------------------------------------------------
  163. end subroutine Agrif_Set_interp
  164. !===================================================================================================
  165. !
  166. !===================================================================================================
  167. ! subroutine Agrif_Set_bcinterp
  168. !---------------------------------------------------------------------------------------------------
  169. subroutine Agrif_Set_bcinterp ( tabvarsindic, interp, interp1, interp2, interp3, interp4, &
  170. interp11, interp12, interp21, interp22 )
  171. !---------------------------------------------------------------------------------------------------
  172. INTEGER, intent(in) :: tabvarsindic !< indice of the variable in tabvars
  173. INTEGER, OPTIONAL, intent(in) :: interp, interp1, interp2, interp3, interp4
  174. INTEGER, OPTIONAL, intent(in) :: interp11, interp12, interp21, interp22
  175. !
  176. INTEGER :: indic ! indice of the variable in tabvars
  177. TYPE(Agrif_Variable), pointer :: var
  178. !
  179. indic = Agrif_Curgrid%tabvars_i(tabvarsindic)%iarray0
  180. !
  181. if (indic <= 0) then
  182. var => Agrif_Search_Variable(Agrif_Mygrid,-indic)
  183. else
  184. print*,"Agrif_Set_bcinterp : warning indic >= 0 !!!"
  185. var => Agrif_Mygrid % tabvars(indic)
  186. endif
  187. !
  188. var % type_interp_bc = Agrif_Constant
  189. !
  190. if (present(interp)) var % type_interp_bc = interp
  191. if (present(interp1)) var % type_interp_bc(:,1) = interp1
  192. if (present(interp11)) var % type_interp_bc(1,1) = interp11
  193. if (present(interp12)) var % type_interp_bc(1,2) = interp12
  194. if (present(interp2)) var % type_interp_bc(:,2) = interp2
  195. if (present(interp21)) var % type_interp_bc(2,1) = interp21
  196. if (present(interp22)) var % type_interp_bc(2,2) = interp22
  197. if (present(interp3)) var % type_interp_bc(:,3) = interp3
  198. if (present(interp4)) var % type_interp_bc(:,4) = interp4
  199. !---------------------------------------------------------------------------------------------------
  200. end subroutine Agrif_Set_bcinterp
  201. !===================================================================================================
  202. !
  203. !===================================================================================================
  204. ! subroutine Agrif_Set_UpdateType
  205. !---------------------------------------------------------------------------------------------------
  206. subroutine Agrif_Set_UpdateType ( tabvarsindic, update, update1, update2, &
  207. update3, update4, update5 )
  208. !---------------------------------------------------------------------------------------------------
  209. INTEGER, intent(in) :: tabvarsindic !< indice of the variable in tabvars
  210. INTEGER, OPTIONAL, intent(in) :: update, update1, update2, update3, update4, update5
  211. !
  212. INTEGER :: indic ! indice of the variable in tabvars
  213. type(Agrif_Variable), pointer :: root_var
  214. !
  215. indic = Agrif_Curgrid%tabvars_i(tabvarsindic)%iarray0
  216. !
  217. if (indic <= 0) then
  218. root_var => Agrif_Search_Variable(Agrif_Mygrid,-indic)
  219. else
  220. print*,"Agrif_Set_UpdateType : warning indic >= 0 !!!"
  221. root_var => Agrif_Mygrid % tabvars(indic)
  222. endif
  223. !
  224. root_var % type_update = Agrif_Update_Copy
  225. if (present(update)) root_var % type_update = update
  226. if (present(update1)) root_var % type_update(1) = update1
  227. if (present(update2)) root_var % type_update(2) = update2
  228. if (present(update3)) root_var % type_update(3) = update3
  229. if (present(update4)) root_var % type_update(4) = update4
  230. if (present(update5)) root_var % type_update(5) = update5
  231. !---------------------------------------------------------------------------------------------------
  232. end subroutine Agrif_Set_UpdateType
  233. !===================================================================================================
  234. !
  235. !===================================================================================================
  236. ! subroutine Agrif_Set_restore
  237. !---------------------------------------------------------------------------------------------------
  238. subroutine Agrif_Set_restore ( tabvarsindic )
  239. !---------------------------------------------------------------------------------------------------
  240. INTEGER, intent(in) :: tabvarsindic !< indice of the variable in tabvars
  241. !
  242. INTEGER :: indic ! indice of the variable in tabvars
  243. !
  244. indic = Agrif_Curgrid%tabvars_i(tabvarsindic)%iarray0
  245. !
  246. Agrif_Mygrid%tabvars(indic) % restore = .TRUE.
  247. !---------------------------------------------------------------------------------------------------
  248. end subroutine Agrif_Set_restore
  249. !===================================================================================================
  250. !
  251. !===================================================================================================
  252. ! subroutine Agrif_Init_variable
  253. !---------------------------------------------------------------------------------------------------
  254. subroutine Agrif_Init_variable ( tabvarsindic, procname )
  255. !---------------------------------------------------------------------------------------------------
  256. INTEGER, intent(in) :: tabvarsindic !< indice of the variable in tabvars
  257. procedure() :: procname !< Data recovery procedure
  258. !
  259. if ( Agrif_Curgrid%level <= 0 ) return
  260. !
  261. call Agrif_Interp_variable(tabvarsindic, procname)
  262. call Agrif_Bc_variable(tabvarsindic, procname, 1.)
  263. !---------------------------------------------------------------------------------------------------
  264. end subroutine Agrif_Init_variable
  265. !===================================================================================================
  266. !
  267. !===================================================================================================
  268. ! subroutine Agrif_Bc_variable
  269. !---------------------------------------------------------------------------------------------------
  270. subroutine Agrif_Bc_variable ( tabvarsindic, procname, calledweight )
  271. !---------------------------------------------------------------------------------------------------
  272. integer, intent(in) :: tabvarsindic !< indice of the variable in tabvars
  273. procedure() :: procname
  274. real, optional, intent(in) :: calledweight
  275. !
  276. real :: weight
  277. logical :: pweight
  278. integer :: indic
  279. integer :: nbdim
  280. type(Agrif_Variable), pointer :: root_var
  281. type(Agrif_Variable), pointer :: parent_var
  282. type(Agrif_Variable), pointer :: child_var
  283. type(Agrif_Variable), pointer :: child_tmp ! Temporary variable on the child grid
  284. !
  285. if ( Agrif_Curgrid%level <= 0 ) return
  286. !
  287. indic = Agrif_Curgrid%tabvars_i(tabvarsindic)%iarray0
  288. !
  289. if ( present(calledweight) ) then
  290. weight = calledweight
  291. pweight = .true.
  292. else
  293. weight = 0.
  294. pweight = .false.
  295. endif
  296. !
  297. if (indic <= 0) then
  298. child_var => Agrif_Search_Variable(Agrif_Curgrid,-indic)
  299. parent_var => child_var % parent_var
  300. root_var => child_var % root_var
  301. else
  302. print*,"Agrif_Bc_variable : warning indic >= 0 !!!"
  303. child_var => Agrif_Curgrid % tabvars(indic)
  304. parent_var => Agrif_Curgrid % parent % tabvars(indic)
  305. root_var => Agrif_Mygrid % tabvars(indic)
  306. endif
  307. !
  308. nbdim = root_var % nbdim
  309. !
  310. select case( nbdim )
  311. case(1)
  312. allocate(parray1(child_var%lb(1):child_var%ub(1)))
  313. case(2)
  314. allocate(parray2(child_var%lb(1):child_var%ub(1), &
  315. child_var%lb(2):child_var%ub(2) ))
  316. case(3)
  317. allocate(parray3(child_var%lb(1):child_var%ub(1), &
  318. child_var%lb(2):child_var%ub(2), &
  319. child_var%lb(3):child_var%ub(3) ))
  320. case(4)
  321. allocate(parray4(child_var%lb(1):child_var%ub(1), &
  322. child_var%lb(2):child_var%ub(2), &
  323. child_var%lb(3):child_var%ub(3), &
  324. child_var%lb(4):child_var%ub(4) ))
  325. case(5)
  326. allocate(parray5(child_var%lb(1):child_var%ub(1), &
  327. child_var%lb(2):child_var%ub(2), &
  328. child_var%lb(3):child_var%ub(3), &
  329. child_var%lb(4):child_var%ub(4), &
  330. child_var%lb(5):child_var%ub(5) ))
  331. case(6)
  332. allocate(parray6(child_var%lb(1):child_var%ub(1), &
  333. child_var%lb(2):child_var%ub(2), &
  334. child_var%lb(3):child_var%ub(3), &
  335. child_var%lb(4):child_var%ub(4), &
  336. child_var%lb(5):child_var%ub(5), &
  337. child_var%lb(6):child_var%ub(6) ))
  338. end select
  339. !
  340. ! Create temporary child variable
  341. allocate(child_tmp)
  342. !
  343. child_tmp % root_var => root_var
  344. child_tmp % oldvalues2D => child_var % oldvalues2D
  345. !
  346. ! Index indicating if a space interpolation is necessary
  347. child_tmp % interpIndex = child_var % interpIndex
  348. child_tmp % list_interp => child_var % list_interp
  349. child_tmp % Interpolationshouldbemade = child_var % Interpolationshouldbemade
  350. !
  351. child_tmp % point = child_var % point
  352. child_tmp % lb = child_var % lb
  353. child_tmp % ub = child_var % ub
  354. !
  355. child_tmp % bcinf = child_var % bcinf
  356. child_tmp % bcsup = child_var % bcsup
  357. !
  358. child_tmp % childarray = child_var % childarray
  359. child_tmp % memberin = child_var % memberin
  360. !
  361. call Agrif_CorrectVariable(parent_var, child_tmp, pweight, weight, procname)
  362. !
  363. child_var % childarray = child_tmp % childarray
  364. child_var % memberin = child_tmp % memberin
  365. !
  366. child_var % oldvalues2D => child_tmp % oldvalues2D
  367. child_var % list_interp => child_tmp % list_interp
  368. !
  369. child_var % interpIndex = child_tmp % interpIndex
  370. !
  371. deallocate(child_tmp)
  372. !
  373. select case( nbdim )
  374. case(1); deallocate(parray1)
  375. case(2); deallocate(parray2)
  376. case(3); deallocate(parray3)
  377. case(4); deallocate(parray4)
  378. case(5); deallocate(parray5)
  379. case(6); deallocate(parray6)
  380. end select
  381. !---------------------------------------------------------------------------------------------------
  382. end subroutine Agrif_Bc_variable
  383. !===================================================================================================
  384. !
  385. !===================================================================================================
  386. ! subroutine Agrif_Interp_variable
  387. !---------------------------------------------------------------------------------------------------
  388. subroutine Agrif_Interp_variable ( tabvarsindic, procname )
  389. !---------------------------------------------------------------------------------------------------
  390. integer, intent(in) :: tabvarsindic !< indice of the variable in tabvars
  391. procedure() :: procname !< Data recovery procedure
  392. !
  393. integer :: nbdim
  394. integer :: indic ! indice of the variable in tabvars
  395. logical :: torestore
  396. type(Agrif_Variable), pointer :: root_var
  397. type(Agrif_Variable), pointer :: parent_var ! Variable on the parent grid
  398. type(Agrif_Variable), pointer :: child_var ! Variable on the parent grid
  399. type(Agrif_Variable), pointer :: child_tmp ! Temporary variable on the child grid
  400. !
  401. if ( Agrif_Curgrid%level <= 0 ) return
  402. !
  403. indic = Agrif_Curgrid%tabvars_i(tabvarsindic)%iarray0
  404. !
  405. if (indic <= 0) then
  406. child_var => Agrif_Search_Variable(Agrif_Curgrid,-indic)
  407. parent_var => child_var % parent_var
  408. root_var => child_var % root_var
  409. else
  410. print*,"Agrif_Interp_variable : warning indic >= 0 !!!"
  411. child_var => Agrif_Curgrid % tabvars(indic)
  412. parent_var => Agrif_Curgrid % parent % tabvars(indic)
  413. root_var => Agrif_Mygrid % tabvars(indic)
  414. endif
  415. !
  416. nbdim = root_var % nbdim
  417. torestore = root_var % restore
  418. !
  419. allocate(child_tmp)
  420. !
  421. child_tmp % root_var => root_var
  422. child_tmp % nbdim = root_var % nbdim
  423. child_tmp % point = child_var % point
  424. child_tmp % lb = child_var % lb
  425. child_tmp % ub = child_var % ub
  426. child_tmp % interpIndex = child_var % interpIndex
  427. child_tmp % list_interp => child_var % list_interp
  428. child_tmp % Interpolationshouldbemade = child_var % Interpolationshouldbemade
  429. !
  430. if ( torestore ) then
  431. select case( nbdim )
  432. case(1)
  433. parray1 = child_var % array1
  434. child_tmp % restore1D => child_var % restore1D
  435. case(2)
  436. parray2 = child_var % array2
  437. child_tmp % restore2D => child_var % restore2D
  438. case(3)
  439. parray3 = child_var % array3
  440. child_tmp % restore3D => child_var % restore3D
  441. case(4)
  442. parray4 = child_var % array4
  443. child_tmp % restore4D => child_var % restore4D
  444. case(5)
  445. parray5 = child_var % array5
  446. child_tmp % restore5D => child_var % restore5D
  447. case(6)
  448. parray6 = child_var % array6
  449. child_tmp % restore6D => child_var % restore6D
  450. end select
  451. endif
  452. !
  453. call Agrif_InterpVariable(parent_var, child_tmp, torestore, procname)
  454. !
  455. child_var % list_interp => child_tmp % list_interp
  456. !
  457. deallocate(child_tmp)
  458. !---------------------------------------------------------------------------------------------------
  459. end subroutine Agrif_Interp_variable
  460. !===================================================================================================
  461. !
  462. !===================================================================================================
  463. ! subroutine Agrif_Update_Variable
  464. !---------------------------------------------------------------------------------------------------
  465. subroutine Agrif_Update_Variable ( tabvarsindic, procname, &
  466. locupdate, locupdate1, locupdate2, locupdate3, locupdate4 )
  467. !---------------------------------------------------------------------------------------------------
  468. integer, intent(in) :: tabvarsindic !< Indice of the variable in tabvars
  469. procedure() :: procname !< Data recovery procedure
  470. integer, dimension(2), intent(in), optional :: locupdate
  471. integer, dimension(2), intent(in), optional :: locupdate1
  472. integer, dimension(2), intent(in), optional :: locupdate2
  473. integer, dimension(2), intent(in), optional :: locupdate3
  474. integer, dimension(2), intent(in), optional :: locupdate4
  475. !---------------------------------------------------------------------------------------------------
  476. integer :: indic
  477. integer :: nbdim
  478. integer, dimension(6) :: updateinf ! First positions where interpolations are calculated
  479. integer, dimension(6) :: updatesup ! Last positions where interpolations are calculated
  480. type(Agrif_Variable), pointer :: root_var
  481. type(Agrif_Variable), pointer :: parent_var
  482. type(Agrif_Variable), pointer :: child_var
  483. !
  484. if ( Agrif_Root() .AND. (.not.agrif_coarse) ) return
  485. if (agrif_curgrid%grand_mother_grid) return
  486. !
  487. indic = Agrif_Curgrid%tabvars_i(tabvarsindic)%iarray0
  488. !
  489. if (indic <= 0) then
  490. child_var => Agrif_Search_Variable(Agrif_Curgrid, -indic)
  491. parent_var => child_var % parent_var
  492. if (.not.associated(parent_var)) then
  493. ! can occur during the first update of Agrif_Coarsegrid (if any)
  494. parent_var => Agrif_Search_Variable(Agrif_Curgrid % parent, -indic)
  495. child_var % parent_var => parent_var
  496. endif
  497. root_var => child_var % root_var
  498. else
  499. print*,"Agrif_Update_Variable : warning indic >= 0 !!!"
  500. root_var => Agrif_Mygrid % tabvars(indic)
  501. child_var => Agrif_Curgrid % tabvars(indic)
  502. parent_var => Agrif_Curgrid % parent % tabvars(indic)
  503. endif
  504. !
  505. nbdim = root_var % nbdim
  506. !
  507. updateinf = -99
  508. updatesup = -99
  509. !
  510. if ( present(locupdate) ) then
  511. updateinf(1:nbdim) = locupdate(1)
  512. updatesup(1:nbdim) = locupdate(2)
  513. endif
  514. !
  515. if ( present(locupdate1) ) then
  516. updateinf(1) = locupdate1(1)
  517. updatesup(1) = locupdate1(2)
  518. endif
  519. !
  520. if ( present(locupdate2) ) then
  521. updateinf(2) = locupdate2(1)
  522. updatesup(2) = locupdate2(2)
  523. endif
  524. if ( present(locupdate3) ) then
  525. updateinf(3) = locupdate3(1)
  526. updatesup(3) = locupdate3(2)
  527. endif
  528. if ( present(locupdate4) ) then
  529. updateinf(4) = locupdate4(1)
  530. updatesup(4) = locupdate4(2)
  531. endif
  532. !
  533. call Agrif_UpdateVariable( parent_var, child_var, updateinf, updatesup, procname )
  534. !---------------------------------------------------------------------------------------------------
  535. end subroutine Agrif_Update_Variable
  536. !===================================================================================================
  537. !
  538. !===================================================================================================
  539. ! subroutine Agrif_Save_ForRestore0D
  540. !---------------------------------------------------------------------------------------------------
  541. subroutine Agrif_Save_ForRestore0D ( tabvarsindic0, tabvarsindic )
  542. !---------------------------------------------------------------------------------------------------
  543. integer, intent(in) :: tabvarsindic0, tabvarsindic
  544. !
  545. type(Agrif_Variable), pointer :: root_var, save_var
  546. integer :: nbdim
  547. !
  548. root_var => Agrif_Mygrid % tabvars(tabvarsindic0)
  549. save_var => Agrif_Curgrid % tabvars(tabvarsindic0)
  550. nbdim = root_var % nbdim
  551. !
  552. select case(nbdim)
  553. case(2); call Agrif_Save_ForRestore2D(save_var % array2, tabvarsindic)
  554. case(3); call Agrif_Save_ForRestore3D(save_var % array3, tabvarsindic)
  555. case(4); call Agrif_Save_ForRestore4D(save_var % array4, tabvarsindic)
  556. end select
  557. !---------------------------------------------------------------------------------------------------
  558. end subroutine Agrif_Save_ForRestore0D
  559. !===================================================================================================
  560. !
  561. !===================================================================================================
  562. ! subroutine Agrif_Save_ForRestore2D
  563. !---------------------------------------------------------------------------------------------------
  564. subroutine Agrif_Save_ForRestore2D ( q, tabvarsindic )
  565. !---------------------------------------------------------------------------------------------------
  566. real, dimension(:,:), intent(in) :: q
  567. integer, intent(in) :: tabvarsindic
  568. !
  569. type(Agrif_Variable), pointer :: root_var, save_var
  570. integer :: indic
  571. !
  572. indic = tabvarsindic
  573. if (tabvarsindic >= 0) then
  574. if (Agrif_Curgrid%tabvars_i(tabvarsindic)%nbdim == 0) then
  575. indic = Agrif_Curgrid%tabvars_i(tabvarsindic)%iarray0
  576. endif
  577. endif
  578. !
  579. if (indic <= 0) then
  580. save_var => Agrif_Search_Variable(Agrif_Curgrid,-indic)
  581. root_var => Agrif_Search_Variable(Agrif_Mygrid,-indic)
  582. else
  583. save_var => Agrif_Curgrid % tabvars(indic)
  584. root_var => Agrif_Mygrid % tabvars(indic)
  585. endif
  586. !
  587. if ( .not.allocated(save_var%array2) ) then
  588. allocate(save_var%array2(save_var%lb(1):save_var%ub(1), &
  589. save_var%lb(2):save_var%ub(2)))
  590. endif
  591. !
  592. save_var % array2 = q
  593. root_var % restore = .true.
  594. !---------------------------------------------------------------------------------------------------
  595. end subroutine Agrif_Save_ForRestore2D
  596. !===================================================================================================
  597. !
  598. !===================================================================================================
  599. ! subroutine Agrif_Save_ForRestore3D
  600. !---------------------------------------------------------------------------------------------------
  601. subroutine Agrif_Save_ForRestore3D ( q, tabvarsindic )
  602. !---------------------------------------------------------------------------------------------------
  603. real, dimension(:,:,:), intent(in) :: q
  604. integer, intent(in) :: tabvarsindic
  605. !
  606. type(Agrif_Variable), pointer :: root_var, save_var
  607. integer :: indic
  608. !
  609. indic = tabvarsindic
  610. if (tabvarsindic >= 0) then
  611. if (Agrif_Curgrid%tabvars_i(tabvarsindic)%nbdim == 0) then
  612. indic = Agrif_Curgrid%tabvars_i(tabvarsindic)%iarray0
  613. endif
  614. endif
  615. !
  616. if (indic <= 0) then
  617. save_var => Agrif_Search_Variable(Agrif_Curgrid,-indic)
  618. root_var => Agrif_Search_Variable(Agrif_Mygrid,-indic)
  619. else
  620. save_var => Agrif_Curgrid % tabvars(indic)
  621. root_var => Agrif_Mygrid % tabvars(indic)
  622. endif
  623. !
  624. if ( .not.allocated(save_var%array3) ) then
  625. allocate(save_var%array3(save_var%lb(1):save_var%ub(1), &
  626. save_var%lb(2):save_var%ub(2), &
  627. save_var%lb(3):save_var%ub(3)))
  628. endif
  629. !
  630. save_var % array3 = q
  631. root_var % restore = .true.
  632. !---------------------------------------------------------------------------------------------------
  633. end subroutine Agrif_Save_ForRestore3D
  634. !===================================================================================================
  635. !
  636. !===================================================================================================
  637. ! subroutine Agrif_Save_ForRestore4D
  638. !---------------------------------------------------------------------------------------------------
  639. subroutine Agrif_Save_ForRestore4D ( q, tabvarsindic )
  640. !---------------------------------------------------------------------------------------------------
  641. real, dimension(:,:,:,:), intent(in) :: q
  642. integer, intent(in) :: tabvarsindic
  643. !
  644. type(Agrif_Variable), pointer :: root_var, save_var
  645. integer :: indic
  646. !
  647. indic = tabvarsindic
  648. if (tabvarsindic >= 0) then
  649. if (Agrif_Curgrid%tabvars_i(tabvarsindic)%nbdim == 0) then
  650. indic = Agrif_Curgrid%tabvars_i(tabvarsindic)%iarray0
  651. endif
  652. endif
  653. !
  654. if (indic <= 0) then
  655. save_var => Agrif_Search_Variable(Agrif_Curgrid,-indic)
  656. root_var => Agrif_Search_Variable(Agrif_Mygrid,-indic)
  657. else
  658. save_var => Agrif_Curgrid % tabvars(indic)
  659. root_var => Agrif_Mygrid % tabvars(indic)
  660. endif
  661. !
  662. if (.not.allocated(save_var%array4)) then
  663. allocate(save_var%array4(save_var%lb(1):save_var%ub(1),&
  664. save_var%lb(2):save_var%ub(2),&
  665. save_var%lb(3):save_var%ub(3),&
  666. save_var%lb(4):save_var%ub(4)))
  667. endif
  668. !
  669. save_var % array4 = q
  670. root_var % restore = .true.
  671. !---------------------------------------------------------------------------------------------------
  672. end subroutine Agrif_Save_ForRestore4D
  673. !===================================================================================================
  674. !
  675. end module Agrif_BcFunction